OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
zfac_asm_ELT.F File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine zmumps_elt_asm_s_2_s_init (nelt, frt_ptr, frt_elt, n, inode, iw, liw, a, la, nbrows, nbcols, opassw, opeliw, step, ptrist, ptrast, itloc, rhs_mumps, fils, ptrarw, ptraiw, intarr, dblarr, icntl, keep, keep8, myid, lrgroups)
subroutine zmumps_asm_slave_elements (inode, n, nelt, iw, liw, ioldps, a, la, poselt, keep, keep8, itloc, fils, ptraiw, ptrarw, intarr, dblarr, lintarr, ldblarr, frt_ptr, frt_elt, rhs_mumps, lrgroups)

Function/Subroutine Documentation

◆ zmumps_asm_slave_elements()

subroutine zmumps_asm_slave_elements ( integer, intent(in) inode,
integer, intent(in) n,
integer, intent(in) nelt,
integer, dimension(liw), intent(in) iw,
integer, intent(in) liw,
integer, intent(in) ioldps,
complex(kind=8), dimension(la), intent(inout) a,
integer(8), intent(in) la,
integer(8), intent(in) poselt,
integer, dimension(500), intent(in) keep,
integer(8), dimension(150), intent(in) keep8,
integer, dimension(n+keep(253)), intent(inout) itloc,
integer, dimension(n), intent(in) fils,
integer(8), dimension(nelt+1), intent(in) ptraiw,
integer(8), dimension(nelt+1), intent(in) ptrarw,
integer, dimension(lintarr), intent(in) intarr,
complex(kind=8), dimension(ldblarr), intent(in) dblarr,
integer(8), intent(in) lintarr,
integer(8), intent(in) ldblarr,
integer, dimension(n+1), intent(in) frt_ptr,
integer, dimension(nelt), intent(in) frt_elt,
complex(kind=8), dimension(keep(255)), intent(in) rhs_mumps,
integer, dimension(n), intent(in) lrgroups )

Definition at line 79 of file zfac_asm_ELT.F.

83!$ USE OMP_LIB
84 USE zmumps_ana_lr, ONLY : get_cut
85 USE zmumps_lr_core, ONLY : max_cluster
87 IMPLICIT NONE
88 INTEGER, intent(in) :: N, NELT, LIW, IOLDPS, INODE
89 INTEGER(8), intent(in) :: LA, POSELT, LINTARR, LDBLARR
90 INTEGER, intent(in) :: IW(LIW)
91 INTEGER, intent(in) :: KEEP(500)
92 INTEGER(8), intent(in) :: KEEP8(150)
93 INTEGER, intent(inout) :: ITLOC(N+KEEP(253))
94 COMPLEX(kind=8), intent(inout) :: A(LA)
95 COMPLEX(kind=8), intent(in) :: RHS_MUMPS(KEEP(255))
96 INTEGER, intent(in) :: INTARR(LINTARR)
97 COMPLEX(kind=8), intent(in) :: DBLARR(LDBLARR)
98 INTEGER, intent(in) :: FRT_PTR(N+1), FRT_ELT(NELT)
99 INTEGER, intent(in) :: FILS(N)
100 INTEGER(8), intent(in) :: PTRAIW(NELT+1), PTRARW(NELT+1)
101 INTEGER, INTENT(IN) :: LRGROUPS(N)
102!$ INTEGER :: NOMP
103!$ INTEGER(8) :: CHUNK8
104!$ INTEGER :: CHUNK
105 include 'mumps_headers.h'
106 INTEGER :: HF, NBROWF, NBCOLF, NASS, NSLAVES
107 INTEGER :: ILOC, IELL, ELTI, ELBEG, NUMELT
108 INTEGER(8) :: SIZE_ELTI8
109 INTEGER :: I, J, K, K1, K2
110 INTEGER :: IPOS, IPOS1, IPOS2, JPOS, IJROW
111 INTEGER :: IN
112 INTEGER(8) :: II8, JJ8, J18, J28
113 INTEGER(8) :: AINPUT8
114 INTEGER(8) :: AII8
115 INTEGER(8) :: APOS, APOS2, ICT12
116 INTEGER, POINTER, DIMENSION(:) :: BEGS_BLR_LS
117 INTEGER :: NB_BLR_LS, NPARTSCB, NPARTSASS, MAXI_CLUSTER,
118 & IBCKSZ2, MINSIZE, TOPDIAG
119 INTEGER(8) :: JJ3
120 INTEGER :: K1RHS, K2RHS, JFirstRHS
121 COMPLEX(kind=8) ZERO
122 parameter( zero = (0.0d0,0.0d0) )
123 nbcolf = iw(ioldps+keep(ixsz))
124 nbrowf = iw(ioldps+2+keep(ixsz))
125 nass = iw(ioldps+1+keep(ixsz))
126 nslaves= iw(ioldps+5 + keep(ixsz))
127 hf = 6 + nslaves + keep(ixsz)
128!$ NOMP = OMP_GET_MAX_THREADS()
129 IF (keep(50) .EQ. 0 .OR. nbrowf .LT. keep(63)) THEN
130!$ CHUNK8 = int(KEEP(361),8)
131!$OMP PARALLEL DO PRIVATE(JJ8) SCHEDULE(STATIC, CHUNK8)
132!$OMP& IF (int(NBROWF,8)*int(NBCOLF,8) > int(KEEP(361),8)
133!$OMP& .AND. NOMP .GT. 1)
134 DO jj8=poselt, poselt+int(nbrowf,8)*int(nbcolf,8)-1_8
135 a(jj8) = zero
136 ENDDO
137!$OMP END PARALLEL DO
138 ELSE
139 topdiag = 0
140 IF (iw(ioldps+xxlr).GE.1) THEN
141 CALL get_cut(iw(ioldps+hf:ioldps+hf+nbrowf-1), 0,
142 & nbrowf, lrgroups, npartscb,
143 & npartsass, begs_blr_ls)
144 nb_blr_ls = npartscb
145 call max_cluster(begs_blr_ls,nb_blr_ls+1,maxi_cluster)
146 DEALLOCATE(begs_blr_ls)
147 CALL compute_blr_vcs(keep(472), ibcksz2, keep(488), nass)
148 minsize = int(ibcksz2 / 2)
149 topdiag = max(2*minsize + maxi_cluster-1, topdiag)
150 ENDIF
151!$ CHUNK = max( KEEP(360)/2,
152!$ & ((NBROWF+NOMP-1)/NOMP +2) / 3 )
153!$OMP PARALLEL DO PRIVATE(APOS,JJ3,JJ8) SCHEDULE(STATIC,CHUNK)
154!$OMP& IF (NBROWF .GT. KEEP(360) .AND. NOMP .GT. 1)
155 DO jj8 = 0_8, int(nbrowf-1,8)
156 apos = poselt+ jj8*int(nbcolf,8)
157 jj3 = min( int(nbcolf,8) - 1_8,
158 & jj8 + int(nbcolf-nbrowf,8) + topdiag )
159 a(apos: apos+jj3) = zero
160 ENDDO
161!$OMP END PARALLEL DO
162 ENDIF
163 k1 = ioldps + hf + nbrowf
164 k2 = k1 + nbcolf - 1
165 jpos = 1
166 DO k = k1, k2
167 j = iw(k)
168 itloc(j) = -jpos
169 jpos = jpos + 1
170 END DO
171 k1 = ioldps + hf
172 k2 = k1 + nbrowf - 1
173 jpos = 1
174 IF ((keep(253).GT.0).AND.(keep(50).NE.0)) THEN
175 k1rhs = 0
176 k2rhs = -1
177 DO k = k1, k2
178 j = iw(k)
179 itloc(j) = -itloc(j)*nbcolf + jpos
180 IF ((k1rhs.EQ.0).AND.(j.GT.n)) THEN
181 k1rhs = k
182 jfirstrhs=j-n
183 ENDIF
184 jpos = jpos + 1
185 ENDDO
186 IF (k1rhs.GT.0) k2rhs=k2
187 IF ( k2rhs.GE.k1rhs ) THEN
188 in = inode
189 DO WHILE (in.GT.0)
190 ijrow = -itloc(in)
191 DO k = k1rhs, k2rhs
192 j = iw(k)
193 i = itloc(j)
194 iloc = mod(i,nbcolf)
195 apos = poselt+int(iloc-1,8)*int(nbcolf,8) +
196 & int(ijrow-1,8)
197 a(apos) = a(apos) + rhs_mumps(
198 & (jfirstrhs+(k-k1rhs)-1)*keep(254)+ in)
199 ENDDO
200 in = fils(in)
201 ENDDO
202 ENDIF
203 ELSE
204 DO k = k1, k2
205 j = iw(k)
206 itloc(j) = -itloc(j)*nbcolf + jpos
207 jpos = jpos + 1
208 END DO
209 ENDIF
210 elbeg = frt_ptr(inode)
211 numelt = frt_ptr(inode+1) - elbeg
212 DO iell=elbeg,elbeg+numelt-1
213 elti = frt_elt(iell)
214 j18= ptraiw(elti)
215 j28= ptraiw(elti+1)-1_8
216 aii8 = ptrarw(elti)
217 size_elti8 = j28 - j18 + 1_8
218 DO ii8=j18,j28
219 i = itloc(intarr(ii8))
220 IF (keep(50).EQ.0) THEN
221 IF (i.LE.0) cycle
222 ainput8 = aii8 + ii8 - j18
223 ipos = mod(i,nbcolf)
224 ict12 = poselt + int(ipos-1,8) * int(nbcolf,8)
225 DO jj8 = j18, j28
226 jpos = itloc(intarr(jj8))
227 IF (jpos.LE.0) THEN
228 jpos = -jpos
229 ELSE
230 jpos = jpos/nbcolf
231 END IF
232 apos2 = ict12 + int(jpos - 1,8)
233 a(apos2) = a(apos2) + dblarr(ainput8)
234 ainput8 = ainput8 + size_elti8
235 END DO
236 ELSE
237 IF ( i .EQ. 0 ) THEN
238 aii8 = aii8 + j28 - ii8 + 1_8
239 cycle
240 ENDIF
241 IF ( i .LE. 0 ) THEN
242 ipos1 = -i
243 ipos2 = 0
244 ELSE
245 ipos1 = i/nbcolf
246 ipos2 = mod(i,nbcolf)
247 END IF
248 ict12 = poselt + int(ipos2-1,8)*int(nbcolf,8)
249 DO jj8=ii8,j28
250 aii8 = aii8 + 1_8
251 j = itloc(intarr(jj8))
252 IF ( j .EQ. 0 ) cycle
253 IF ( ipos2.EQ.0 .AND. j.LE.0) cycle
254 IF ( j .LE. 0 ) THEN
255 jpos = -j
256 ELSE
257 jpos = j/nbcolf
258 END IF
259 IF ( (ipos1.GE.jpos) .AND. (ipos2.GT.0) ) THEN
260 apos2 = ict12 + int(jpos - 1,8)
261 a(apos2) = a(apos2) + dblarr(aii8-1_8)
262 END IF
263 IF ( (ipos1.LT.jpos) .AND. (j.GT.0) ) THEN
264 ipos = mod(j,nbcolf)
265 jpos = ipos1
266 apos2 = poselt + int(ipos-1,8)*int(nbcolf,8)
267 & + int(jpos - 1,8)
268 a(apos2) = a(apos2) + dblarr(aii8-1_8)
269 END IF
270 END DO
271 END IF
272 END DO
273 END DO
274 k1 = ioldps + hf + nbrowf
275 k2 = k1 + nbcolf - 1
276 DO k = k1, k2
277 j = iw(k)
278 itloc(j) = 0
279 END DO
#define min(a, b)
Definition macros.h:20
#define max(a, b)
Definition macros.h:21
subroutine compute_blr_vcs(k472, ibcksz, maxsize, nass)
Definition lr_common.F:18
subroutine get_cut(iwr, nass, ncb, lrgroups, npartscb, npartsass, cut)
Definition zana_lr.F:25
subroutine max_cluster(cut, cut_size, maxi_cluster)
Definition zlr_core.F:1304

◆ zmumps_elt_asm_s_2_s_init()

subroutine zmumps_elt_asm_s_2_s_init ( integer nelt,
integer, dimension(n+1) frt_ptr,
integer, dimension(nelt) frt_elt,
integer n,
integer inode,
integer, dimension(liw) iw,
integer liw,
complex(kind=8), dimension(la) a,
integer(8) la,
integer nbrows,
integer nbcols,
double precision opassw,
double precision opeliw,
integer, dimension(n) step,
integer, dimension(keep(28)) ptrist,
integer(8), dimension(keep(28)) ptrast,
integer, dimension(n+keep(253)) itloc,
complex(kind=8), dimension(keep(255)) rhs_mumps,
integer, dimension(n) fils,
integer(8), dimension(nelt+1), intent(in) ptrarw,
integer(8), dimension(nelt+1), intent(in) ptraiw,
integer, dimension(keep8(27)) intarr,
complex(kind=8), dimension(keep8(26)) dblarr,
integer, dimension(60) icntl,
integer, dimension(500) keep,
integer(8), dimension(150) keep8,
integer myid,
integer, dimension(n), intent(in) lrgroups )

Definition at line 14 of file zfac_asm_ELT.F.

23 IMPLICIT NONE
24 INTEGER NELT, N,LIW
25 INTEGER(8) :: LA
26 INTEGER KEEP(500), ICNTL(60)
27 INTEGER(8) KEEP8(150)
28 INTEGER INODE, MYID
29 INTEGER NBROWS, NBCOLS
30 INTEGER(8) :: PTRAST(KEEP(28))
31 INTEGER IW(LIW), ITLOC(N+KEEP(253)), STEP(N),
32 & PTRIST(KEEP(28)), FILS(N)
33 INTEGER(8), INTENT(IN) :: PTRARW(NELT+1), PTRAIW(NELT+1)
34 COMPLEX(kind=8) :: RHS_MUMPS(KEEP(255))
35 INTEGER INTARR(KEEP8(27))
36 INTEGER FRT_PTR(N+1), FRT_ELT(NELT)
37 COMPLEX(kind=8) :: A(LA)
38 COMPLEX(kind=8) :: DBLARR(KEEP8(26))
39 DOUBLE PRECISION OPASSW, OPELIW
40 INTEGER, INTENT(IN) :: LRGROUPS(N)
41 INTEGER(8) :: POSELT
42 COMPLEX(kind=8), DIMENSION(:), POINTER :: A_PTR
43 INTEGER(8) :: LA_PTR
44 INTEGER IOLDPS, NBCOLF, NBROWF, NSLAVES, HF,
45 & K1,K2,K,J,JPOS,NASS
46 COMPLEX(kind=8) ZERO
47 parameter( zero = (0.0d0,0.0d0) )
48 include 'mumps_headers.h'
49 ioldps = ptrist(step(inode))
50 CALL zmumps_dm_set_dynptr( iw(ioldps+xxs), a, la,
51 & ptrast(step(inode)), iw(ioldps+xxd), iw(ioldps+xxr),
52 & a_ptr, poselt, la_ptr )
53 nbcolf = iw(ioldps+keep(ixsz))
54 nbrowf = iw(ioldps+2+keep(ixsz))
55 nass = iw(ioldps+1+keep(ixsz))
56 nslaves = iw(ioldps+5+keep(ixsz))
57 hf = 6 + nslaves + keep(ixsz)
58 IF (nass.LT.0) THEN
59 nass = -nass
60 iw(ioldps+1+keep(ixsz)) = nass
61 CALL zmumps_asm_slave_elements( inode, n, nelt, iw, liw,
62 & ioldps, a_ptr(poselt), la_ptr, 1_8, keep, keep8, itloc, fils,
63 & ptraiw, ptrarw,
64 & intarr, dblarr, keep8(27), keep8(26), frt_ptr, frt_elt,
65 & rhs_mumps, lrgroups)
66 ENDIF
67 IF (nbrows.GT.0) THEN
68 k1 = ioldps + hf + nbrowf
69 k2 = k1 + nbcolf - 1
70 jpos = 1
71 DO k = k1, k2
72 j = iw(k)
73 itloc(j) = jpos
74 jpos = jpos + 1
75 ENDDO
76 ENDIF
77 RETURN
subroutine zmumps_dm_set_dynptr(cb_state, a, la, pamaster_or_ptrast, ixxd, ixxr, son_a, iachk, recsize)
subroutine zmumps_asm_slave_elements(inode, n, nelt, iw, liw, ioldps, a, la, poselt, keep, keep8, itloc, fils, ptraiw, ptrarw, intarr, dblarr, lintarr, ldblarr, frt_ptr, frt_elt, rhs_mumps, lrgroups)