FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
hecmw_jadm_33.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_33
16  public :: hecmw_jad_finalize_33
18  public :: hecmw_jad_matvec_33
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 :: wp1(:), wp2(:), wp3(:)
27  integer(kind=kint) :: INITIALIZED = 0
28 
29 contains
30 
31  subroutine hecmw_jad_init_33(hecMAT)
32  type(hecmwst_matrix) :: hecmat
33  allocate(wp1(hecmat%NP), wp2(hecmat%NP), wp3(hecmat%NP))
34  allocate(ajad((hecmat%NPL+hecmat%NPU)*9))
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_33
41 
42  subroutine hecmw_jad_finalize_33()
43  deallocate(ajad)
44  deallocate(jajad)
45  deallocate(jadord)
46  deallocate(iajad)
47  deallocate(wp1,wp2,wp3)
48  initialized = 0
49  end subroutine hecmw_jad_finalize_33
50 
52  integer(kind=kint) :: hecmw_jad_is_initialized_33
53  hecmw_jad_is_initialized_33 = initialized
54  end function hecmw_jad_is_initialized_33
55 
56  subroutine hecmw_jad_matvec_33(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
65  real(kind=kreal) :: x1, x2, x3
66 
67  start_time= hecmw_wtime()
68  call hecmw_update_r (hecmesh, x, hecmat%NP, 3)
69  end_time= hecmw_wtime()
70  commtime = commtime + end_time - start_time
71 
72  d => hecmat%D
73 
74  !$OMP PARALLEL PRIVATE(i)
75  !$OMP DO
76  do i= 1, hecmat%N
77  x1= x(3*i-2)
78  x2= x(3*i-1)
79  x3= x(3*i )
80  y(3*i -2)= d(9*i-8)*x1 + d(9*i-7)*x2 + d(9*i-6)*x3
81  y(3*i -1)= d(9*i-5)*x1 + d(9*i-4)*x2 + d(9*i-3)*x3
82  y(3*i )= d(9*i-2)*x1 + d(9*i-1)*x2 + d(9*i )*x3
83  enddo
84  !$OMP END DO
85  !$OMP END PARALLEL
86 
87  call matjad(hecmat%N, mjad, iajad, jajad, ajad, jadord, x, y, wp1, wp2, wp3)
88  end subroutine hecmw_jad_matvec_33
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
102  integer(kind = kint) :: i, j, js, je, in, jc
103  integer(kind = kint), allocatable :: len(:), lenz(:), jadreord(:)
104 
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(jc*9-8:jc*9) = hecmat%AL(9*j-8:9*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(jc*9-8:jc*9) = hecmat%AU(9*j-8:9*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, MJAD, IAJAD, JAJAD, AJAD, JADORD, X, Y, W1, W2, W3)
171  integer(kind=kint) :: N, MJAD
172  integer(kind=kint) :: IAJAD(*), JAJAD(*), JADORD(*)
173  real(kind=kreal) :: ajad(*), x(*), y(*), w1(*), w2(*), w3(*)
174 
175  integer(kind=kint) :: I, K, NZ, IXX
176  real(kind=kreal) :: x1, x2, x3
177 
178  !$OMP PARALLEL PRIVATE(I,K,X1,X2,X3,IXX)
179  !$OMP DO
180  do i=1,n
181  w1(i)=0.d0
182  w2(i)=0.d0
183  w3(i)=0.d0
184  enddo
185  !$OMP END DO
186 
187  do nz=1,mjad
188  !$OMP DO
189  do k=iajad(nz),iajad(nz+1)-1
190  x1=x(jajad(k)*3-2)
191  x2=x(jajad(k)*3-1)
192  x3=x(jajad(k)*3 )
193  ixx = k-iajad(nz)+1
194  w1(ixx)=w1(ixx) + ajad(k*9-8)*x1 + ajad(k*9-7)*x2 + ajad(k*9-6)*x3
195  w2(ixx)=w2(ixx) + ajad(k*9-5)*x1 + ajad(k*9-4)*x2 + ajad(k*9-3)*x3
196  w3(ixx)=w3(ixx) + ajad(k*9-2)*x1 + ajad(k*9-1)*x2 + ajad(k*9-0)*x3
197  enddo
198  !$OMP END DO
199  enddo
200 
201  !$OMP DO
202  do i=1,n
203  y(3*i-2)=y(3*i-2)+w1(jadord(i))
204  y(3*i-1)=y(3*i-1)+w2(jadord(i))
205  y(3*i )=y(3*i )+w3(jadord(i))
206  enddo
207  !$OMP END DO
208  !$OMP END PARALLEL
209  end subroutine matjad
210 
211 end module hecmw_jad_type_33
hecmw_jad_type_33::matjad
subroutine matjad(N, MJAD, IAJAD, JAJAD, AJAD, JADORD, X, Y, W1, W2, W3)
Definition: hecmw_jadm_33.f90:170
hecmw_util::hecmw_wtime
real(kind=kreal) function hecmw_wtime()
Definition: hecmw_util_f.F90:549
hecmw_jad_type_33::hecmw_jad_matvec_33
subroutine, public hecmw_jad_matvec_33(hecMESH, hecMAT, X, Y, COMMtime)
Definition: hecmw_jadm_33.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_jad_type_33::hecmw_jad_is_initialized_33
integer(kind=kint) function, public hecmw_jad_is_initialized_33()
Definition: hecmw_jadm_33.f90:52
hecmw_util::kreal
integer(kind=4), parameter kreal
Definition: hecmw_util_f.F90:16
hecmw_jad_type_33::hecmw_jad_finalize_33
subroutine, public hecmw_jad_finalize_33()
Definition: hecmw_jadm_33.f90:43
hecmw_jad_type_33::hecmw_jad_init_33
subroutine, public hecmw_jad_init_33(hecMAT)
Definition: hecmw_jadm_33.f90:32
hecmw_jad_type_33
Jagged Diagonal Matrix storage for vector processors. Original code was provided by JAMSTEC.
Definition: hecmw_jadm_33.f90:8
m_hecmw_comm_f::hecmw_update_r
subroutine hecmw_update_r(hecMESH, val, n, m)
Definition: hecmw_comm_f.F90:683
hecmw_util::hecmwst_matrix
Definition: hecmw_util_f.F90:444