FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
hecmw_jadm_nn.f90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
2 ! Copyright (c) 2019 FrontISTR Commons
3 ! This software is released under the MIT License, see LICENSE.txt
4 !-------------------------------------------------------------------------------
7 
9  use hecmw_util
10  use m_hecmw_comm_f
11  implicit none
12 
13  private
14 
15  public :: hecmw_jad_init_nn
16  public :: hecmw_jad_finalize_nn
18  public :: hecmw_jad_matvec_nn
19 
20  !C---------------------- AU&AL
21  real(kind=kreal), allocatable :: ajad(:)
22  integer(kind=kint), allocatable :: JAJAD(:)
23  integer(kind=kint), allocatable :: JADORD(:)
24  integer(kind=kint), allocatable :: IAJAD(:)
25  integer(kind=kint) :: MJAD
26  real(kind=kreal), allocatable :: wp(:,:)
27  integer(kind=kint) :: INITIALIZED = 0
28 
29 contains
30 
31  subroutine hecmw_jad_init_nn(hecMAT)
32  type(hecmwst_matrix) :: hecmat
33  allocate(wp(hecmat%NDOF,hecmat%NP))
34  allocate(ajad((hecmat%NPL+hecmat%NPU)*hecmat%NDOF*hecmat%NDOF))
35  allocate(jajad(hecmat%NPL+hecmat%NPU))
36  allocate(jadord(hecmat%NP))
37  allocate(iajad(hecmat%NP+1))
38  call repack(hecmat%N, hecmat, mjad, ajad, jajad, iajad, jadord)
39  initialized = 1
40  end subroutine hecmw_jad_init_nn
41 
42  subroutine hecmw_jad_finalize_nn()
43  deallocate(ajad)
44  deallocate(jajad)
45  deallocate(jadord)
46  deallocate(iajad)
47  deallocate(wp)
48  initialized = 0
49  end subroutine hecmw_jad_finalize_nn
50 
52  integer(kind=kint) :: hecmw_jad_is_initialized_nn
53  hecmw_jad_is_initialized_nn = initialized
54  end function hecmw_jad_is_initialized_nn
55 
56  subroutine hecmw_jad_matvec_nn(hecMESH, hecMAT, X, Y, COMMtime)
57  type(hecmwst_local_mesh), intent(in) :: hecmesh
58  type(hecmwst_matrix), intent(in), target :: hecmat
59  real(kind=kreal), intent(in) :: x(:)
60  real(kind=kreal), intent(out) :: y(:)
61  real(kind=kreal), intent(inout) :: commtime
62  real(kind=kreal) :: start_time, end_time
63  real(kind=kreal), pointer :: d(:)
64  integer(kind=kint) :: i,idof,jdof,ndof,ndof2
65 
66  start_time= hecmw_wtime()
67  call hecmw_update_r (hecmesh, x, hecmat%NP,hecmat%NDOF)
68  end_time= hecmw_wtime()
69  commtime = commtime + end_time - start_time
70 
71  d => hecmat%D
72  ndof = hecmat%NDOF
73  ndof2 = ndof*ndof
74  y=0.0d0
75 
76  !$OMP PARALLEL PRIVATE(i, idof, jdof)
77  !$OMP DO
78  do i= 1, hecmat%N
79  do idof=1,hecmat%NDOF
80  do jdof=1,hecmat%NDOF
81  y(ndof*(i-1)+idof) = y(ndof*(i-1)+idof) + d(ndof2*(i-1)+ndof*(idof-1)+jdof)*x(ndof*(i-1)+jdof)
82  end do
83  end do
84  enddo
85  !$OMP END DO
86  !$OMP END PARALLEL
87  call matjad(hecmat%N,hecmat%NDOF, mjad, iajad, jajad, ajad, jadord, x, y, wp)
88  end subroutine hecmw_jad_matvec_nn
89 
90  subroutine repack(N, hecMAT, MJAD, AJAD, JAJAD, IAJAD, JADORD)
91  use hecmw_util
92  !C---------------------------------
93  type (hecmwst_matrix) :: hecmat
94  !C----------------------
95  integer(kind = kint) :: n, mjad
96  real(kind = kreal), dimension(*) :: ajad
97  integer(kind = kint), dimension(*) :: jajad
98  integer(kind = kint), dimension(*) :: iajad
99  integer(kind = kint), dimension(*) :: jadord
100 
101  integer(kind = kint) :: ijad, maxnz, minnz,ndof,ndof2
102  integer(kind = kint) :: i, j, js, je, in, jc
103  integer(kind = kint), allocatable :: len(:), lenz(:), jadreord(:)
104  ndof = hecmat%NDOF;ndof2=ndof*ndof
105  allocate(len(n))
106  allocate(jadreord(n))
107  do i=1,n
108  len(i)= hecmat%indexL(i) - hecmat%indexL(i-1) &
109  & + hecmat%indexU(i) - hecmat%indexU(i-1)
110  end do
111  maxnz=maxval(len(1:n))
112  minnz=minval(len(1:n))
113  mjad =maxnz
114  allocate(lenz(0:mjad))
115  lenz = 0
116  do i=1,n
117  lenz(len(i))=lenz(len(i))+1
118  enddo
119  do i=maxnz-1,minnz,-1
120  lenz(i)=lenz(i)+lenz(i+1)
121  enddo
122  do i=1,n
123  jadord(i)=lenz(len(i))
124  lenz(len(i))=lenz(len(i))-1
125  enddo
126  do i=1,n
127  jadreord(jadord(i))=i
128  enddo
129  do i=1,n
130  lenz(len(jadreord(i)))=i
131  enddo
132  do i=maxnz-1,1,-1
133  lenz(i)=max(lenz(i+1),lenz(i))
134  enddo
135  iajad(1)=1
136  do i=1,maxnz
137  iajad(i+1)=iajad(i)+lenz(i)
138  enddo
139  len=0
140  do i= 1, n
141  ijad=jadord(i)
142  js= hecmat%indexL(i-1) + 1
143  je= hecmat%indexL(i )
144  do j=js,je
145  in = hecmat%itemL(j)
146  len(ijad)=len(ijad)+1
147  jc=iajad(len(ijad))+ijad-1
148  ajad(ndof2*(jc-1)+1:ndof2*(jc)) = hecmat%AL(ndof2*(j-1)+1:ndof2*(j))
149  jajad(jc) = in
150  end do
151  end do
152  do i= 1, n
153  ijad=jadord(i)
154  js= hecmat%indexU(i-1) + 1
155  je= hecmat%indexU(i )
156  do j=js,je
157  in = hecmat%itemU(j)
158  len(ijad)=len(ijad)+1
159  jc=iajad(len(ijad))+ijad-1
160  ajad(ndof2*(jc-1)+1:ndof2*(jc)) = hecmat%AU(ndof2*(j-1)+1:ndof2*(j))
161  jajad(jc) = in
162  end do
163  end do
164  deallocate(len)
165  deallocate(jadreord)
166  deallocate(lenz)
167  end subroutine repack
168 
169  subroutine matjad(N, NDOF, MJAD, IAJAD, JAJAD, AJAD, JADORD, X, Y, W)
171  integer(kind=kint) :: N,NDOF, MJAD,NDOF2
172  integer(kind=kint) :: IAJAD(*), JAJAD(*), JADORD(*)
173  real(kind=kreal) :: ajad(*), x(*), y(*), w(ndof,n)
174  integer(kind=kint) :: I, K, NZ, IXX, idof,jdof
175  ndof2=ndof*ndof
176 
177  w=0.0d0
178 
179  do nz=1,mjad
180  !$OMP PARALLEL PRIVATE(K,IXX,idof,jdof)
181  !$OMP DO
182  do k=iajad(nz),iajad(nz+1)-1
183  ixx = k-iajad(nz)+1
184  do idof = 1, ndof
185  do jdof = 1, ndof
186  w(idof,ixx)=w(idof,ixx)+ajad(ndof2*(k-1)+ndof*(idof-1)+jdof)*x(ndof*(jajad(k)-1)+jdof)
187  end do
188  end do
189  enddo
190  !$OMP END DO
191  !$OMP END PARALLEL
192  enddo
193 
194  !$OMP PARALLEL PRIVATE(I,idof)
195  !$OMP DO
196  do i=1,n
197  do idof = 1, ndof
198  y(ndof*(i-1)+idof)=y(ndof*(i-1)+idof)+w(idof,jadord(i))
199  end do
200  enddo
201  !$OMP END DO
202  !$OMP END PARALLEL
203  end subroutine matjad
204 
205 end module hecmw_jad_type_nn
hecmw_jad_type_nn::hecmw_jad_is_initialized_nn
integer(kind=kint) function, public hecmw_jad_is_initialized_nn()
Definition: hecmw_jadm_nn.f90:52
hecmw_jad_type_nn::matjad
subroutine matjad(N, NDOF, MJAD, IAJAD, JAJAD, AJAD, JADORD, X, Y, W)
Definition: hecmw_jadm_nn.f90:170
hecmw_util::hecmw_wtime
real(kind=kreal) function hecmw_wtime()
Definition: hecmw_util_f.F90:549
hecmw_jad_type_nn::hecmw_jad_matvec_nn
subroutine, public hecmw_jad_matvec_nn(hecMESH, hecMAT, X, Y, COMMtime)
Definition: hecmw_jadm_nn.f90:57
hecmw_util
I/O and Utility.
Definition: hecmw_util_f.F90:7
hecmw_util::hecmwst_local_mesh
Definition: hecmw_util_f.F90:234
m_hecmw_comm_f
Definition: hecmw_comm_f.F90:6
hecmw_util::kreal
integer(kind=4), parameter kreal
Definition: hecmw_util_f.F90:16
hecmw_jad_type_nn
Jagged Diagonal Matrix storage for vector processors. Original code was provided by JAMSTEC.
Definition: hecmw_jadm_nn.f90:8
hecmw_jad_type_nn::hecmw_jad_init_nn
subroutine, public hecmw_jad_init_nn(hecMAT)
Definition: hecmw_jadm_nn.f90:32
m_hecmw_comm_f::hecmw_update_r
subroutine hecmw_update_r(hecMESH, val, n, m)
Definition: hecmw_comm_f.F90:683
hecmw_jad_type_nn::hecmw_jad_finalize_nn
subroutine, public hecmw_jad_finalize_nn()
Definition: hecmw_jadm_nn.f90:43
hecmw_util::hecmwst_matrix
Definition: hecmw_util_f.F90:444