OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
cfac_asm.F
Go to the documentation of this file.
1C
2C This file is part of MUMPS 5.5.1, released
3C on Tue Jul 12 13:17:24 UTC 2022
4C
5C
6C Copyright 1991-2022 CERFACS, CNRS, ENS Lyon, INP Toulouse, Inria,
7C Mumps Technologies, University of Bordeaux.
8C
9C This version of MUMPS is provided to you free of charge. It is
10C released under the CeCILL-C license
11C (see doc/CeCILL-C_V1-en.txt, doc/CeCILL-C_V1-fr.txt, and
12C https://cecill.info/licences/Licence_CeCILL-C_V1-en.html)
13C
14 SUBROUTINE cmumps_asm_slave_master(N, INODE, IW, LIW, A, LA,
15 & ISON, NBROWS, NBCOLS, ROWLIST,
16 & VALSON, PTLUST_S, PTRAST, STEP, PIMASTER,
17 & OPASSW, IWPOSCB, MYID, KEEP,KEEP8, IS_ofType5or6,
18 & LDA_VALSON )
19 USE cmumps_load
20 IMPLICIT NONE
21 INTEGER KEEP(500)
22 INTEGER(8) KEEP8(150)
23 INTEGER(8) :: LA
24 INTEGER N,LIW,MYID
25 INTEGER INODE,ISON, IWPOSCB
26 INTEGER NBROWS, NBCOLS, LDA_VALSON
27 INTEGER(8) :: PTRAST(KEEP(28))
28 INTEGER IW(LIW), STEP(N), PIMASTER(KEEP(28)),
29 & ptlust_s(keep(28)), rowlist(nbrows)
30 COMPLEX A(LA), VALSON(LDA_VALSON,NBROWS)
31 DOUBLE PRECISION OPASSW
32 LOGICAL, INTENT(IN) :: IS_ofType5or6
33 INTEGER(8) :: POSELT, POSEL1, APOS, JJ2
34 INTEGER HF,HS, NSLAVES, NFRONT, NASS1,
35 & ioldps, istchk, lstk, nslson,nelim,
36 & npivs,ncols,j1,jj,jj1,nrows,
37 & ldafs_pere, ibeg, diag
38 include 'mumps_headers.h'
39 LOGICAL SAME_PROC
40 IOLDPS = ptlust_s(step(inode))
41 poselt = ptrast(step(inode))
42 nfront = iw(ioldps+keep(ixsz))
43 nass1 = iabs(iw(ioldps + 2+keep(ixsz)))
44 nslaves= iw(ioldps+5+keep(ixsz))
45 IF (keep(50).EQ.0) THEN
46 ldafs_pere = nfront
47 ELSE
48 IF ( nslaves .eq. 0 ) THEN
49 ldafs_pere = nfront
50 ELSE
51 ldafs_pere = nass1
52 ENDIF
53 ENDIF
54 hf = 6 + nslaves + keep(ixsz)
55 posel1 = poselt - int(ldafs_pere,8)
56 istchk = pimaster(step(ison))
57 lstk = iw(istchk+keep(ixsz))
58 nslson = iw(istchk + 5+keep(ixsz))
59 hs = 6 + nslson + keep(ixsz)
60 opassw = opassw + dble(nbrows*nbcols)
61 nelim = iw(istchk + 1+keep(ixsz))
62 npivs = iw(istchk + 3+keep(ixsz))
63 IF (npivs.LT.0) npivs = 0
64 ncols = npivs + lstk
65 same_proc = (istchk.LT.iwposcb)
66 IF (same_proc) THEN
67 nrows = ncols
68 ELSE
69 nrows = iw(istchk+2+keep(ixsz))
70 ENDIF
71 j1 = istchk + nrows + hs + npivs
72 IF (keep(50).EQ.0) THEN
73 IF (is_oftype5or6) THEN
74 apos = posel1 + int(rowlist(1),8) * int(ldafs_pere,8)
75 DO jj = 1, nbrows
76 DO jj1 = 1, nbcols
77 jj2 = apos + int(jj1-1,8)
78 a(jj2)=a(jj2)+valson(jj1,jj)
79 ENDDO
80 apos = apos + int(ldafs_pere,8)
81 ENDDO
82 ELSE
83 DO 170 jj = 1, nbrows
84 apos = posel1 + int(rowlist(jj),8) * int(ldafs_pere,8)
85 DO 160 jj1 = 1, nbcols
86 jj2 = apos + int(iw(j1 + jj1 - 1) - 1,8)
87 a(jj2) = a(jj2) + valson(jj1,jj)
88 160 CONTINUE
89 170 CONTINUE
90 ENDIF
91 ELSE
92 IF (is_oftype5or6) THEN
93 apos = posel1 + int(rowlist(1),8) * int(ldafs_pere,8)
94 diag = rowlist(1)
95 DO jj = 1, nbrows
96 DO jj1 = 1, diag
97 jj2 = apos+int(jj1-1,8)
98 a(jj2) = a(jj2) + valson(jj1,jj)
99 ENDDO
100 diag = diag+1
101 apos = apos + int(ldafs_pere,8)
102 ENDDO
103 ELSE
104 DO jj = 1, nbrows
105 IF (rowlist(jj).LE.nass1.and..NOT.is_oftype5or6) THEN
106 apos = posel1 + int(rowlist(jj) - 1,8)
107 DO jj1 = 1, nelim
108 jj2 = apos + int(iw(j1+jj1-1),8)*int(ldafs_pere,8)
109 a(jj2) = a(jj2) + valson(jj1,jj)
110 ENDDO
111 ibeg = nelim+1
112 ELSE
113 ibeg = 1
114 ENDIF
115 apos = posel1 + int(rowlist(jj),8) * int(ldafs_pere,8)
116 DO jj1 = ibeg, nbcols
117 IF (rowlist(jj).LT.iw(j1 + jj1 - 1)) EXIT
118 jj2 = apos + int(iw(j1 + jj1 - 1) - 1,8)
119 a(jj2) = a(jj2) + valson(jj1,jj)
120 ENDDO
121 ENDDO
122 ENDIF
123 ENDIF
124 RETURN
125 END SUBROUTINE cmumps_asm_slave_master
127 & (n, inode, iw, liw, a, la,
128 & nbrows, nbcols,
129 & opassw, opeliw, step, ptrist, ptrast, itloc,
130 & rhs_mumps, fils, ptrarw, ptraiw, intarr, dblarr,
131 & icntl, keep,keep8, myid, lrgroups)
133 IMPLICIT NONE
134 INTEGER N,LIW
135 INTEGER(8) :: LA
136 INTEGER KEEP(500), ICNTL(60)
137 INTEGER(8) KEEP8(150)
138 INTEGER INODE, MYID
139 INTEGER NBROWS, NBCOLS
140 INTEGER(8) :: PTRAST(KEEP(28))
141 INTEGER IW(LIW), ITLOC(N+KEEP(253)), STEP(N),
142 & ptrist(keep(28)), fils(n)
143 INTEGER(8), INTENT(IN) :: PTRARW(N), PTRAIW(N)
144 COMPLEX :: RHS_MUMPS(KEEP(255))
145 COMPLEX :: A(LA)
146 INTEGER :: INTARR(KEEP8(27))
147 COMPLEX :: DBLARR(KEEP8(26))
148 DOUBLE PRECISION OPASSW, OPELIW
149 INTEGER, INTENT(IN) :: LRGROUPS(N)
150 INTEGER(8) :: POSELT
151 COMPLEX, DIMENSION(:), POINTER :: A_PTR
152 INTEGER(8) :: LA_PTR
153 INTEGER IOLDPS, NBCOLF, NBROWF, NSLAVES, HF,
154 & k1,k2,k,j,jpos,nass
155 COMPLEX ZERO
156 parameter( zero = (0.0e0,0.0e0) )
157 include 'mumps_headers.h'
158 ioldps = ptrist(step(inode))
159 CALL cmumps_dm_set_dynptr( iw(ioldps+xxs), a, la,
160 & ptrast(step(inode)), iw(ioldps+xxd), iw(ioldps+xxr),
161 & a_ptr, poselt, la_ptr )
162 nbcolf = iw(ioldps+keep(ixsz))
163 nbrowf = iw(ioldps+2+keep(ixsz))
164 nass = iw(ioldps+1+keep(ixsz))
165 nslaves = iw(ioldps+5+keep(ixsz))
166 hf = 6 + nslaves + keep(ixsz)
167 IF (nass.LT.0) THEN
168 nass = -nass
169 iw(ioldps+1+keep(ixsz)) = nass
170 CALL cmumps_asm_slave_arrowheads(inode, n, iw, liw,
171 & ioldps, a_ptr(poselt), la_ptr, 1_8, keep, keep8,
172 & itloc, fils, ptraiw, ptrarw, intarr, dblarr,
173 & keep8(27), keep8(26),
174 & rhs_mumps, lrgroups)
175 ENDIF
176 IF (nbrows.GT.0) THEN
177 k1 = ioldps + hf + nbrowf
178 k2 = k1 + nbcolf - 1
179 jpos = 1
180 DO k = k1, k2
181 j = iw(k)
182 itloc(j) = jpos
183 jpos = jpos + 1
184 ENDDO
185 ENDIF
186 RETURN
187 END SUBROUTINE cmumps_asm_slave_to_slave_init
189 & (n, inode, iw, liw, nbrows, step, ptrist,
190 & itloc, rhs_mumps, keep,keep8)
191 IMPLICIT NONE
192 INTEGER N, LIW
193 INTEGER KEEP(500)
194 INTEGER(8) KEEP8(150)
195 INTEGER INODE
196 INTEGER NBROWS
197 INTEGER IW(LIW), ITLOC(N+KEEP(253)), STEP(N),
198 & ptrist(keep(28))
199 COMPLEX :: RHS_MUMPS(KEEP(255))
200 include 'mumps_headers.h'
201 INTEGER IOLDPS, NBCOLF, NBROWF, NSLAVES, HF,
202 & K1,K2,K,J
203 ioldps = ptrist(step(inode))
204 nbcolf = iw(ioldps+keep(ixsz))
205 nbrowf = iw(ioldps+2+keep(ixsz))
206 nslaves = iw(ioldps+5+keep(ixsz))
207 hf = 6 + nslaves+keep(ixsz)
208 IF (nbrows.GT.0) THEN
209 k1 = ioldps + hf + nbrowf
210 k2 = k1 + nbcolf - 1
211 DO k = k1, k2
212 j = iw(k)
213 itloc(j) = 0
214 ENDDO
215 ENDIF
216 RETURN
217 END SUBROUTINE cmumps_asm_slave_to_slave_end
218 SUBROUTINE cmumps_asm_slave_to_slave(N, INODE, IW, LIW, A, LA,
219 & NBROWS, NBCOLS, ROWLIST, COLLIST, VALSON,
220 & OPASSW, OPELIW, STEP, PTRIST, PTRAST, ITLOC,
221 & RHS_MUMPS, FILS,
222 & ICNTL, KEEP,KEEP8, MYID, IS_ofType5or6, LDA_VALSON)
224 IMPLICIT NONE
225 INTEGER N,LIW
226 INTEGER(8) :: LA
227 INTEGER KEEP(500), ICNTL(60)
228 INTEGER(8) KEEP8(150)
229 INTEGER INODE, MYID
230 LOGICAL, intent(in) :: IS_ofType5or6
231 INTEGER NBROWS, NBCOLS, LDA_VALSON
232 INTEGER ROWLIST(NBROWS), COLLIST(NBCOLS)
233 INTEGER IW(LIW), ITLOC(N+KEEP(253)), STEP(N),
234 & ptrist(keep(28)), fils(n)
235 COMPLEX :: RHS_MUMPS(KEEP(255))
236 INTEGER(8) :: PTRAST(KEEP(28))
237 COMPLEX A(LA), VALSON(LDA_VALSON,NBROWS)
238 DOUBLE PRECISION OPASSW, OPELIW
239 INTEGER(8) :: POSEL1, POSELT, APOS, K8
240 INTEGER IOLDPS, NBCOLF, NBROWF, NSLAVES, HF,
241 & i,j,nass,idiag
242 COMPLEX, POINTER, DIMENSION(:) :: A_PTR
243 INTEGER(8) :: LA_PTR
244 include 'mumps_headers.h'
245 ioldps = ptrist(step(inode))
246 CALL cmumps_dm_set_dynptr( iw(ioldps+xxs), a, la,
247 & ptrast(step(inode)), iw(ioldps+xxd), iw(ioldps+xxr),
248 & a_ptr, poselt, la_ptr )
249 nbcolf = iw(ioldps+keep(ixsz))
250 nbrowf = iw(ioldps+2+keep(ixsz))
251 nass = iw(ioldps+1+keep(ixsz))
252 IF ( nbrows .GT. nbrowf ) THEN
253 WRITE(*,*) ' ERR: ERROR : NBROWS > NBROWF'
254 WRITE(*,*) ' ERR: INODE =', inode
255 WRITE(*,*) ' ERR: NBROW=',nbrows,'NBROWF=',nbrowf
256 WRITE(*,*) ' ERR: ROW_LIST=', rowlist
257 WRITE(*,*) ' ERR: NBCOLF/NASS=', nbcolf, nass
258 CALL mumps_abort()
259 END IF
260 nslaves = iw(ioldps+5+keep(ixsz))
261 hf = 6 + nslaves+keep(ixsz)
262 IF (nbrows.GT.0) THEN
263 posel1 = poselt - int(nbcolf,8)
264 IF (keep(50).EQ.0) THEN
265 IF (is_oftype5or6) THEN
266 apos = posel1 + int(rowlist(1),8) * int(nbcolf,8)
267 DO i=1, nbrows
268 DO j = 1, nbcols
269 a_ptr(apos+int(j-1,8)) = a_ptr( apos+int(j-1,8)) +
270 & valson(j,i)
271 ENDDO
272 apos = apos + int(nbcolf,8)
273 END DO
274 ELSE
275 DO i=1,nbrows
276 apos = posel1 + int(rowlist(i),8) * int(nbcolf,8)
277 DO j=1,nbcols
278 k8 = apos + int(itloc(collist(j)),8) - 1_8
279 a_ptr(k8) = a_ptr(k8) + valson(j,i)
280 ENDDO
281 ENDDO
282 ENDIF
283 ELSE
284 IF (is_oftype5or6) THEN
285 apos = posel1 + int(rowlist(1),8) * int(nbcolf,8)
286 & + int((nbrows-1),8)*int(nbcolf,8)
287 idiag = 0
288 DO i=nbrows,1,-1
289 DO j=1,nbcols-idiag
290 k8 = apos+int(j-1,8)
291 a_ptr(k8) = a_ptr(k8) + valson(j,i)
292 ENDDO
293 apos = apos - int(nbcolf,8)
294 idiag = idiag + 1
295 ENDDO
296 ELSE
297 DO i=1,nbrows
298 apos = posel1 + int(rowlist(i),8) * int(nbcolf,8)
299 DO j=1,nbcols
300 IF (itloc(collist(j)) .EQ. 0) THEN
301 EXIT
302 ENDIF
303 k8 = apos + int(itloc(collist(j)),8) - 1_8
304 a_ptr(k8) = a_ptr(k8) + valson(j,i)
305 ENDDO
306 ENDDO
307 ENDIF
308 ENDIF
309 opassw = opassw + dble(nbrows*nbcols)
310 ENDIF
311 RETURN
312 END SUBROUTINE cmumps_asm_slave_to_slave
313 SUBROUTINE cmumps_ldlt_asm_niv12_ip( A, LA,
314 & IAFATH, NFRONT, NASS1,
315 & IACB, NCOLS, LCB,
316 & IW, NROWS, NELIM, ETATASS,
317 & CB_IS_COMPRESSED )
318 IMPLICIT NONE
319 INTEGER NFRONT, NASS1
320 INTEGER(8) :: LA
321 INTEGER NCOLS, NROWS, NELIM
322 INTEGER(8) :: LCB
323 COMPLEX A( LA )
324 INTEGER(8) :: IAFATH, IACB
325 INTEGER IW( NCOLS )
326 INTEGER ETATASS
327 LOGICAL CB_IS_COMPRESSED
328 COMPLEX ZERO
329 PARAMETER( ZERO = (0.0e0,0.0e0) )
330 INTEGER I, J
331 INTEGER(8) :: APOS, POSELT
332 INTEGER(8) :: IPOSCB, IBEGCBROW, IENDFRONT
333 LOGICAL RESET_TO_ZERO, RISK_OF_SAME_POS,
334 & risk_of_same_pos_this_line
335 iendfront = iafath+int(nfront,8)*int(nfront,8)-1_8
336 iposcb=1_8
337 reset_to_zero = iacb .LT. iendfront + 1_8
338 risk_of_same_pos = iacb + lcb .EQ. iendfront + 1_8
339 risk_of_same_pos_this_line = .false.
340 DO i=1, nrows
341 poselt = int(iw(i)-1,8) * int(nfront,8)
342 IF (.NOT. cb_is_compressed ) THEN
343 iposcb = 1_8 + int(i - 1,8) * int(ncols,8)
344 IF (iacb+iposcb-1_8 .GE. iendfront + 1_8) THEN
345 reset_to_zero = .false.
346 ENDIF
347 ENDIF
348 IF ( risk_of_same_pos ) THEN
349 IF (i.EQ.nrows .OR. .NOT. cb_is_compressed) THEN
350 IF ( iafath + poselt + int(iw(i)-1,8) .EQ.
351 & iacb+iposcb+int(i-1-1,8)) THEN
352 risk_of_same_pos_this_line = .true.
353 ENDIF
354 ENDIF
355 ENDIF
356 IF (reset_to_zero) THEN
357 IF ( risk_of_same_pos_this_line ) THEN
358 DO j=1, i
359 apos = poselt + int(iw( j ),8)
360 IF (iafath + apos - 1_8.NE. iacb+iposcb-1_8) THEN
361 a(iafath+ apos -1_8) = a(iacb+iposcb-1_8)
362 a(iacb+iposcb-1_8) = zero
363 ENDIF
364 iposcb = iposcb + 1_8
365 ENDDO
366 ELSE
367#if defined(__ve__)
368!NEC$ IVDEP
369#endif
370 DO j=1, i
371 apos = poselt + int(iw( j ),8)
372 a(iafath+ apos -1_8) = a(iacb+iposcb-1_8)
373 a(iacb+iposcb-1_8) = zero
374 iposcb = iposcb + 1_8
375 ENDDO
376 ENDIF
377 ELSE
378#if defined(__ve__)
379!NEC$ IVDEP
380#endif
381 DO j=1, i
382 apos = poselt + int(iw( j ),8)
383 a(iafath+ apos -1_8) = a(iacb+iposcb-1_8)
384 iposcb = iposcb + 1_8
385 ENDDO
386 ENDIF
387 IF (.NOT. cb_is_compressed ) THEN
388 ibegcbrow = iacb+iposcb-1_8
389 IF ( ibegcbrow .LE. iendfront ) THEN
390 a(ibegcbrow:ibegcbrow+int(ncols-i,8)-1_8)=zero
391 ENDIF
392 ENDIF
393 IF (iacb+iposcb-1_8 .GE. iendfront + 1_8) THEN
394 reset_to_zero = .false.
395 ENDIF
396 ENDDO
397 RETURN
398 END SUBROUTINE cmumps_ldlt_asm_niv12_ip
399 SUBROUTINE cmumps_ldlt_asm_niv12( A, LA, SON_A,
400 & IAFATH, NFRONT, NASS1,
401 & NCOLS, LCB,
402 & IW, NROWS, NELIM, ETATASS,
403 & CB_IS_COMPRESSED
404!$ & , K360
405 & )
406 IMPLICIT NONE
407 INTEGER NFRONT, NASS1
408 INTEGER(8) :: LA
409 INTEGER NCOLS, NROWS, NELIM
410 INTEGER(8) :: LCB
411 COMPLEX A( LA )
412 COMPLEX SON_A( LCB )
413 INTEGER(8) :: IAFATH
414 INTEGER IW( NCOLS )
415 INTEGER ETATASS
416 LOGICAL CB_IS_COMPRESSED
417!$ INTEGER, INTENT(in):: K360
418 COMPLEX ZERO
419 PARAMETER( ZERO = (0.0e0,0.0e0) )
420 INTEGER I, J
421 INTEGER(8) :: APOS, POSELT
422 INTEGER(8) :: IPOSCB
423!$ LOGICAL :: OMP_FLAG
424 IF ((etatass.EQ.0) .OR. (etatass.EQ.1)) THEN
425 iposcb = 1_8
426#if defined(__ve__)
427!NEC$ IVDEP
428#endif
429 DO i = 1, nelim
430 poselt = int( iw( i ) - 1, 8 ) * int(nfront, 8)
431 IF (.NOT. cb_is_compressed) THEN
432 iposcb = 1_8 + int( i - 1, 8 ) * int(ncols,8)
433 ENDIF
434#if defined(__ve__)
435!NEC$ IVDEP
436#endif
437 DO j = 1, i
438 apos = poselt + int(iw( j ),8)
439 a(iafath+ apos -1_8) = a(iafath+ apos -1_8)
440 & + son_a(iposcb)
441 iposcb = iposcb + 1_8
442 END DO
443 END DO
444 ENDIF
445 IF ((etatass.EQ.0).OR.(etatass.EQ.1)) THEN
446!$ OMP_FLAG = (NROWS-NELIM).GE.K360
447!$OMP PARALLEL DO PRIVATE(IPOSCB, POSELT, J, APOS) IF (OMP_FLAG)
448 DO i = nelim + 1, nrows
449 IF (cb_is_compressed) THEN
450 iposcb = (int(i,8) * int(i-1,8)) / 2_8 + 1_8
451 ELSE
452 iposcb = int(i-1,8) * int(ncols,8) + 1_8
453 ENDIF
454 poselt = int(iw( i ),8)
455 IF (poselt.LE. int(nass1,8)) THEN
456#if defined(__ve__)
457!NEC$ IVDEP
458#endif
459 DO j = 1, nelim
460 apos = poselt + int( iw( j ) - 1, 8 ) * int(nfront,8)
461 a(iafath+apos-1_8) = a(iafath+apos-1_8) +
462 & son_a(iposcb)
463 iposcb = iposcb + 1_8
464 END DO
465 ELSE
466 poselt = int( iw( i ) - 1, 8 ) * int(nfront, 8)
467#if defined(__ve__)
468!NEC$ IVDEP
469#endif
470 DO j = 1, nelim
471 apos = poselt + int(iw( j ), 8)
472 a(iafath+apos-1_8) = a(iafath+apos-1_8)
473 & + son_a(iposcb)
474 iposcb = iposcb + 1_8
475 END DO
476 ENDIF
477 IF (etatass.EQ.1) THEN
478 poselt = int( iw( i ) - 1, 8 ) * int(nfront, 8)
479 DO j = nelim + 1, i
480 IF (iw(j).GT.nass1) EXIT
481 apos = poselt + int(iw( j ), 8)
482 a(iafath+apos-1_8) = a(iafath+apos-1_8)
483 & + son_a(iposcb)
484 iposcb = iposcb +1_8
485 END DO
486 ELSE
487 poselt = int( iw( i ) - 1, 8 ) * int(nfront, 8)
488#if defined(__ve__)
489!NEC$ IVDEP
490#endif
491 DO j = nelim + 1, i
492 apos = poselt + int(iw( j ), 8)
493 a(iafath+apos-1_8) = a(iafath+apos-1_8)
494 & + son_a(iposcb)
495 iposcb = iposcb + 1_8
496 END DO
497 ENDIF
498 END DO
499!$OMP END PARALLEL DO
500 ELSE
501 DO i= nrows, nelim+1, -1
502 IF (cb_is_compressed) THEN
503 iposcb = (int(i,8)*int(i+1,8))/2_8
504 ELSE
505 iposcb = int(i-1,8) * int(ncols,8) + int(i,8)
506 ENDIF
507 poselt = int(iw( i ),8)
508 IF (poselt.LE.int(nass1,8)) EXIT
509 poselt = int( iw( i ) - 1, 8 ) * int(nfront, 8)
510 DO j=i,nelim+1, -1
511 IF (iw(j).LE.nass1) EXIT
512 apos = poselt + int(iw( j ), 8)
513 a(iafath+apos-1_8) = a(iafath+apos-1_8)
514 & + son_a(iposcb)
515 iposcb = iposcb - 1_8
516 ENDDO
517 ENDDO
518 ENDIF
519 RETURN
520 END SUBROUTINE cmumps_ldlt_asm_niv12
521 SUBROUTINE cmumps_restore_indices(N, ISON, INODE, IWPOSCB,
522 & PIMASTER, PTLUST_S, IW, LIW, STEP, KEEP,KEEP8)
523 IMPLICIT NONE
524 INTEGER N, ISON, INODE, IWPOSCB
525 INTEGER KEEP(500), STEP(N)
526 INTEGER(8) KEEP8(150)
527 INTEGER PIMASTER(KEEP(28)), PTLUST_S(KEEP(28))
528 INTEGER LIW
529 INTEGER IW(LIW)
530 INTEGER ISTCHK, LSTK, NSLSON, HS, NROWS, NCOLS, NPIVS, NELIM
531 INTEGER IOLDPS, NFRONT, NSLAVES, ICT11, HF
532 INTEGER J1, J2, J3, JJ, JPOS
533 LOGICAL SAME_PROC
534 include 'mumps_headers.h'
535 istchk = pimaster(step(ison))
536 lstk = iw(istchk+keep(ixsz))
537 nslson = iw(istchk+5+keep(ixsz))
538 hs = 6 + nslson + keep(ixsz)
539 nelim = iw(istchk + 1+keep(ixsz))
540 npivs = iw(istchk + 3+keep(ixsz))
541 ncols = npivs + lstk
542 IF ( npivs < 0 ) npivs = 0
543 same_proc = istchk < iwposcb
544 IF (same_proc) THEN
545 nrows = ncols
546 ELSE
547 nrows = iw(istchk+2+keep(ixsz))
548 ENDIF
549 j1 = istchk + nrows + hs + npivs
550 IF (keep(50).NE.0) THEN
551 j2 = j1 + lstk - 1
552 DO jj = j1, j2
553 iw(jj) = iw(jj - nrows)
554 ENDDO
555 ELSE
556 j2 = j1 + lstk - 1
557 j3 = j1 + nelim
558 DO jj = j3, j2
559 iw(jj) = iw(jj - nrows)
560 ENDDO
561 IF (nelim .NE. 0) THEN
562 ioldps = ptlust_s(step(inode))
563 nfront = iw(ioldps+keep(ixsz))
564 nslaves= iw(ioldps+5+keep(ixsz))
565 hf = 6 + nslaves+keep(ixsz)
566 ict11 = ioldps + hf - 1 + nfront
567 j3 = j3 - 1
568 DO 190 jj = j1, j3
569 jpos = iw(jj) + ict11
570 iw(jj) = iw(jpos)
571 190 CONTINUE
572 ENDIF
573 ENDIF
574 RETURN
575 END SUBROUTINE cmumps_restore_indices
576 SUBROUTINE cmumps_asm_max(
577 & N, INODE, IW, LIW, A, LA,
578 & ISON, NBCOLS,
579 & VALSON, PTLUST_S, PTRAST, STEP, PIMASTER,
580 & OPASSW, IWPOSCB,MYID, KEEP,KEEP8 )
581 USE cmumps_load
582 IMPLICIT NONE
583 INTEGER KEEP(500)
584 INTEGER(8) KEEP8(150)
585 INTEGER(8) :: LA
586 INTEGER N,LIW,MYID
587 INTEGER INODE,ISON,IWPOSCB
588 INTEGER NBCOLS
589 INTEGER IW(LIW), STEP(N),
590 & pimaster(keep(28)),
591 & ptlust_s(keep(28))
592 INTEGER(8) PTRAST(KEEP(28))
593 COMPLEX A(LA)
594 REAL VALSON(NBCOLS)
595 DOUBLE PRECISION OPASSW
596 INTEGER HF,HS, NSLAVES, NASS1,
597 & ioldps, istchk,
598 & lstk, nslson,nelim,npivs,ncols, j1,
599 & jj1,nrows
600 INTEGER(8) POSELT, APOS, JJ2
601 INCLUDE 'mumps_headers.h'
602 LOGICAL SAME_PROC
603 INTRINSIC real
604 ioldps = ptlust_s(step(inode))
605 poselt = ptrast(step(inode))
606 nass1 = iabs(iw(ioldps + 2 + keep(ixsz)))
607 nslaves= iw(ioldps+5 + keep(ixsz))
608 hf = 6 + nslaves + keep(ixsz)
609 istchk = pimaster(step(ison))
610 lstk = iw(istchk + keep(ixsz))
611 nslson = iw(istchk + 5 + keep(ixsz))
612 hs = 6 + nslson + keep(ixsz)
613 nelim = iw(istchk + 1 + keep(ixsz))
614 npivs = iw(istchk + 3 + keep(ixsz))
615 IF (npivs.LT.0) npivs = 0
616 ncols = npivs + lstk
617 same_proc = (istchk.LT.iwposcb)
618 IF (same_proc) THEN
619 nrows = ncols
620 ELSE
621 nrows = iw(istchk+2 + keep(ixsz))
622 ENDIF
623 j1 = istchk + nrows + hs + npivs
624 apos = poselt + int(nass1,8)*int(nass1,8) - 1_8
625 DO jj1 = 1, nbcols
626 jj2 = apos+int(iw(j1 + jj1 - 1),8)
627 IF(real(a(jj2)) .LT. valson(jj1)) THEN
628 a(jj2) = cmplx(valson(jj1),kind=kind(a))
629 ENDIF
630 ENDDO
631 RETURN
632 END SUBROUTINE cmumps_asm_max
633 SUBROUTINE cmumps_asm_slave_arrowheads(INODE, N, IW, LIW, IOLDPS,
634 & A, LA, POSELT, KEEP, KEEP8,
635 & ITLOC, FILS, PTRAIW, PTRARW, INTARR, DBLARR,
636 & LINTARR, LDBLARR, RHS_MUMPS, LRGROUPS)
637!$ USE OMP_LIB
638 USE cmumps_ana_lr, ONLY : get_cut
639 USE cmumps_lr_core, ONLY : max_cluster
641 IMPLICIT NONE
642 INTEGER, intent(in) :: N, LIW, IOLDPS, INODE
643 INTEGER(8), intent(in) :: LA, POSELT
644 INTEGER(8), intent(in) :: LINTARR, LDBLARR
645 INTEGER, intent(in) :: IW(LIW)
646 INTEGER, intent(in) :: KEEP(500)
647 INTEGER(8), intent(in) :: KEEP8(150)
648 INTEGER, intent(inout) :: ITLOC(N+KEEP(253))
649 COMPLEX, intent(inout) :: A(LA)
650 COMPLEX, intent(in) :: RHS_MUMPS(KEEP(255))
651 COMPLEX, intent(in) :: DBLARR(LDBLARR)
652 INTEGER, intent(in) :: INTARR(LINTARR)
653 INTEGER, intent(in) :: FILS(N)
654 INTEGER(8), intent(in) :: PTRAIW(N), PTRARW(N)
655 INTEGER, INTENT(IN) :: LRGROUPS(N)
656!$ INTEGER :: NOMP
657!$ INTEGER(8) :: CHUNK8
658!$ INTEGER :: CHUNK
659 include 'mumps_headers.h'
660 INTEGER :: HF, NBROWF, NBCOLF, NASS, NSLAVES
661 INTEGER :: ILOC, J, K, K1, K2, JPOS, IJROW
662 INTEGER :: IN
663 INTEGER(8) :: J18, J28, JJ8, JK8
664 INTEGER(8) :: APOS, ICT12
665 INTEGER(8) :: AINPUT8
666 INTEGER, POINTER, DIMENSION(:) :: BEGS_BLR_LS
667 INTEGER :: NB_BLR_LS, NPARTSCB, NPARTSASS, MAXI_CLUSTER,
668 & ibcksz2, minsize, topdiag
669 INTEGER(8) :: JJ3
670 INTEGER :: K1RHS, K2RHS, JFirstRHS
671 COMPLEX ZERO
672 parameter( zero = (0.0e0,0.0e0) )
673 nbcolf = iw(ioldps+keep(ixsz))
674 nbrowf = iw(ioldps+2+keep(ixsz))
675 nass = iw(ioldps+1+keep(ixsz))
676 nslaves= iw(ioldps+5 + keep(ixsz))
677 hf = 6 + nslaves + keep(ixsz)
678!$ NOMP = OMP_GET_MAX_THREADS()
679 IF (keep(50) .EQ. 0 .OR. nbrowf .LT. keep(63)) THEN
680!$ CHUNK8 = int(KEEP(361),8)
681!$OMP PARALLEL DO PRIVATE(JJ8) SCHEDULE(STATIC, CHUNK8)
682!$OMP& IF (int(NBROWF,8)*int(NBCOLF,8) > int(KEEP(361),8)
683!$OMP& .AND. NOMP .GT. 1)
684 DO jj8=poselt, poselt+int(nbrowf,8)*int(nbcolf,8)-1_8
685 a(jj8) = zero
686 ENDDO
687!$OMP END PARALLEL DO
688 ELSE
689 topdiag = 0
690 IF (iw(ioldps+xxlr).GE.1) THEN
691 CALL get_cut(iw(ioldps+hf:ioldps+hf+nbrowf-1), 0,
692 & nbrowf, lrgroups, npartscb,
693 & npartsass, begs_blr_ls)
694 nb_blr_ls = npartscb
695 call max_cluster(begs_blr_ls,nb_blr_ls+1,maxi_cluster)
696 DEALLOCATE(begs_blr_ls)
697 CALL compute_blr_vcs(keep(472), ibcksz2, keep(488), nass)
698 minsize = int(ibcksz2 / 2)
699 topdiag = max(2*minsize + maxi_cluster-1, topdiag)
700 ENDIF
701!$ CHUNK = max( KEEP(360)/2,
702!$ & ((NBROWF+NOMP-1)/NOMP +2) / 3 )
703!$OMP PARALLEL DO PRIVATE(APOS,JJ3,JJ8) SCHEDULE(STATIC,CHUNK)
704!$OMP& IF (NBROWF .GT. KEEP(360) .AND. NOMP .GT. 1)
705 DO jj8 = 0_8, int(nbrowf-1,8)
706 apos = poselt+ jj8*int(nbcolf,8)
707 jj3 = min( int(nbcolf,8) - 1_8,
708 & jj8 + int(nbcolf-nbrowf,8) + topdiag )
709 a(apos: apos+jj3) = zero
710 ENDDO
711!$OMP END PARALLEL DO
712 ENDIF
713 k1 = ioldps + hf + nbrowf
714 k2 = k1 + nass - 1
715 jpos = 1
716 DO k = k1, k2
717 j = iw(k)
718 itloc(j) = -jpos
719 jpos = jpos + 1
720 ENDDO
721 k1 = ioldps + hf
722 k2 = k1 + nbrowf - 1
723 jpos = 1
724 IF ((keep(253).GT.0).AND.(keep(50).NE.0)) THEN
725 k1rhs = 0
726 k2rhs = -1
727 DO k = k1, k2
728 j = iw(k)
729 itloc(j) = jpos
730 IF ((k1rhs.EQ.0).AND.(j.GT.n)) THEN
731 k1rhs = k
732 jfirstrhs=j-n
733 ENDIF
734 jpos = jpos + 1
735 ENDDO
736 IF (k1rhs.GT.0) k2rhs=k2
737 IF ( k2rhs.GE.k1rhs ) THEN
738 in = inode
739 DO WHILE (in.GT.0)
740 ijrow = -itloc(in)
741 DO k = k1rhs, k2rhs
742 j = iw(k)
743 iloc = itloc(j)
744 apos = poselt+int(iloc-1,8)*int(nbcolf,8) +
745 & int(ijrow-1,8)
746 a(apos) = a(apos) + rhs_mumps(
747 & (jfirstrhs+(k-k1rhs)-1)*keep(254)+in)
748 ENDDO
749 in = fils(in)
750 ENDDO
751 ENDIF
752 ELSE
753 DO k = k1, k2
754 j = iw(k)
755 itloc(j) = jpos
756 jpos = jpos + 1
757 ENDDO
758 ENDIF
759 in = inode
760 DO WHILE (in.GT.0)
761 ainput8 = ptrarw(in)
762 jk8 = ptraiw(in)
763 jj8 = jk8 + 1_8
764 j18 = jj8 + 1_8
765 j28 = j18 + intarr(jk8)
766 ijrow = -itloc(intarr(j18))
767 ict12 = poselt +int(- nbcolf + ijrow - 1,8)
768 DO jj8= j18,j28
769 iloc = itloc(intarr(jj8))
770 IF (iloc.GT.0) THEN
771 apos = ict12 + int(iloc,8)*int(nbcolf,8)
772 a(apos) = a(apos) + dblarr(ainput8)
773 ENDIF
774 ainput8 = ainput8 + 1_8
775 ENDDO
776 in = fils(in)
777 ENDDO
778 k1 = ioldps + hf
779 k2 = k1 + nbrowf + nass - 1
780 DO k = k1, k2
781 j = iw(k)
782 itloc(j) = 0
783 ENDDO
784 RETURN
785 END SUBROUTINE cmumps_asm_slave_arrowheads
786 SUBROUTINE cmumps_set_parpivt1 ( INODE, NFRONT, NASS1, KEEP,
787 & LR_ACTIVATED, PARPIV_T1)
788 IMPLICIT NONE
789 INTEGER, intent(in) :: INODE, NFRONT, NASS1, KEEP(500)
790 LOGICAL, intent(in) :: LR_ACTIVATED
791 INTEGER, intent(out) :: PARPIV_T1
792 INTEGER :: NCB
793 LOGICAL, EXTERNAL :: CMUMPS_IS_TRSM_LARGE_ENOUGH,
795 parpiv_t1 = keep(269)
796 IF (parpiv_t1.EQ.-3) THEN
797 parpiv_t1 = 0
798 ENDIF
799 IF (parpiv_t1.EQ.77) THEN
800 parpiv_t1 = 0
801 ENDIF
802 IF (parpiv_t1.EQ.0) RETURN
803 IF ( (parpiv_t1.EQ.-2).AND.lr_activated ) THEN
804 parpiv_t1 = 1
805 ENDIF
806 ncb = nfront-nass1
807 IF (parpiv_t1.EQ.-2) THEN
808 IF (
809 & ( cmumps_is_trsm_large_enough( nass1, ncb
810 & )
811 & )
812 & .OR.
813 & ( cmumps_is_gemm_large_enough( ncb, ncb, nass1
814 & )
815 & )
816 & ) THEN
817 parpiv_t1 = 1
818 ELSE
819 parpiv_t1 = 0
820 ENDIF
821 ENDIF
822 IF (ncb.EQ.keep(253)) THEN
823 parpiv_t1 = 0
824 ENDIF
825 RETURN
826 END SUBROUTINE cmumps_set_parpivt1
828 & ( m, n
829 & )
830 IMPLICIT NONE
831 INTEGER, INTENT(in) :: m, n
832 DOUBLE PRECISION :: ai
833 INTEGER, PARAMETER :: thres_ai = 400
834 ai = ( dble(m)*dble(n) ) /
835 & ( dble(m)/dble(2) + dble(2)*dble(n) )
836 cmumps_is_trsm_large_enough = (ai.GE.dble(thres_ai))
837 RETURN
838 END FUNCTION cmumps_is_trsm_large_enough
840 & ( m, n, k
841 & )
842 IMPLICIT NONE
843 INTEGER, INTENT(in) :: m, n, k
844 DOUBLE PRECISION :: ai
845 INTEGER, PARAMETER :: thres_ai = 400
846 ai = ( dble(2)*dble(m)*dble(n)*dble(k) ) /
847 & ( dble(m)*dble(n) + dble(m)*dble(k) + dble(k)*dble(n) )
848 cmumps_is_gemm_large_enough = (ai.GE.dble(thres_ai))
849 RETURN
850 END FUNCTION cmumps_is_gemm_large_enough
851 SUBROUTINE cmumps_parpivt1_set_max ( INODE,
852 & A, LAELL8, KEEP, NFRONT,
853 & NASS1, NVSCHUR_K253, NB_POSTPONED)
854 IMPLICIT NONE
855 INTEGER(8), intent(in) :: LAELL8
856 INTEGER, intent(in) :: INODE
857 INTEGER, intent(in) :: KEEP(500), NFRONT, NASS1,
858 & NVSCHUR_K253
859 INTEGER, intent(in) :: NB_POSTPONED
860 COMPLEX, intent(inout) :: A(LAELL8)
861 INTEGER(8) :: APOSMAX, APOS, NASS1_8, NFRONT_8
862 INTEGER :: I, J, NCB
863 COMPLEX :: ZERO
864 REAL :: RMAX
865 parameter( zero = (0.0e0,0.0e0) )
866 nass1_8 = int(nass1, 8)
867 nfront_8 = int(nfront, 8)
868 ncb = nfront-nass1-nvschur_k253
869 IF ((ncb.EQ.0).AND.(nvschur_k253.EQ.0)) CALL mumps_abort()
870 aposmax = laell8 - nass1_8 + 1_8
871 a(aposmax:aposmax+nass1_8-1_8)= zero
872 IF (ncb.EQ.0) RETURN
873 IF (keep(50).EQ.2) THEN
874 apos = 1_8 + (nass1_8*nfront_8)
875 DO i = 1, ncb
876 DO j = 1, nass1
877 rmax = real(a(aposmax+int(j,8)-1_8))
878 rmax = max(rmax, abs(a(apos+int(j,8)-1_8)))
879 a(aposmax+int(j,8)-1_8) = cmplx(rmax,kind=kind(a))
880 ENDDO
881 apos = apos+nfront_8
882 ENDDO
883 ELSE
884 apos = 1_8 + nass1_8
885 DO i = 1, nass1
886 rmax = real(a(aposmax+int(i,8)-1_8))
887 DO j = 1, ncb
888 rmax = max(rmax, abs(a(apos+int(j,8)-1)))
889 ENDDO
890 a(aposmax+int(i,8)-1_8) = cmplx(rmax,kind=kind(a))
891 apos = apos+nfront_8
892 ENDDO
893 ENDIF
894 CALL cmumps_update_parpiv_entries ( inode,
895 & keep, a(aposmax), nass1, nb_postponed)
896 RETURN
897 END SUBROUTINE cmumps_parpivt1_set_max
899 & KEEP, PARPIV, LPARPIV,
900 & NB_POSTPONED)
901 IMPLICIT NONE
902 INTEGER, intent(in) :: INODE, LPARPIV, KEEP(500)
903 COMPLEX, intent(inout):: PARPIV(LPARPIV)
904 INTEGER, intent(in) :: NB_POSTPONED
905 INTEGER :: I
906 REAL :: EPS, RMIN, RZERO, RTMP
907 REAL :: RMAX
908 LOGICAL :: UPDATE_PARPIV
909 parameter( rzero = 0.0e0 )
910 update_parpiv=.false.
911 rmin = huge(rzero)
912 rmax = rzero
913 eps = sqrt(epsilon(rzero))*0.01e0
914 DO i = 1, lparpiv
915 rtmp = real(parpiv(i))
916 IF (rtmp.GT.rzero) THEN
917 rmin = min(rmin, rtmp)
918 ELSE
919 update_parpiv=.true.
920 ENDIF
921 IF (rtmp.LE.eps) update_parpiv=.true.
922 rmax= max(rmax,real(parpiv(i)))
923 ENDDO
924 IF (update_parpiv) THEN
925 IF (rmin.LT.huge(rmin)) THEN
926 rmax= min(rmax, eps)
927 DO i = 1, lparpiv-nb_postponed
928 rtmp = real(parpiv(i))
929 IF (rtmp.LE.eps) THEN
930 parpiv(i) = cmplx(-rmax, kind=kind(parpiv))
931 ENDIF
932 ENDDO
933 IF (nb_postponed.GT.0) THEN
934 DO i=lparpiv-nb_postponed+1, lparpiv
935 rtmp = real(parpiv(i))
936 IF (rtmp.LE.eps) THEN
937 parpiv(i) = cmplx(-rmax, kind=kind(parpiv))
938 ENDIF
939 ENDDO
940 ENDIF
941 ENDIF
942 ENDIF
943 RETURN
944 END SUBROUTINE cmumps_update_parpiv_entries
946 & (n, inode, iw, liw, a, la, keep, perm,
947 & ioldps, poselt,
948 & nfront, nass1, lr_activated, parpiv_t1,
949 & nb_postponed)
951 IMPLICIT NONE
952 INTEGER, intent(in) :: N, INODE, LIW, IOLDPS,
953 & nfront, nass1, nb_postponed
954 INTEGER(8), intent(in) :: LA, POSELT
955 INTEGER, intent(in) :: IW (LIW), PERM(N), KEEP(500)
956 LOGICAL, intent(in) :: LR_ACTIVATED
957 COMPLEX, intent(inout) :: A(LA)
958 INTEGER, intent(inout) :: PARPIV_T1
959 INTEGER :: NVSCHUR_K253, IROW_L
960 INTEGER(8) :: LAELL8, NFRONT8
961 include 'mumps_headers.h'
962 IF (parpiv_t1.EQ.-999) THEN
963 CALL cmumps_set_parpivt1 ( inode, nfront, nass1, keep,
964 & lr_activated, parpiv_t1)
965 ELSE IF ((parpiv_t1.NE.0.AND.parpiv_t1.NE.1)) THEN
966 parpiv_t1 = 0
967 ENDIF
968 IF (parpiv_t1.NE.0) THEN
969 IF ((keep(114).EQ.1) .AND. (keep(116).GT.0) ) THEN
970 irow_l = ioldps+6+keep(ixsz)+nass1
972 & n,
973 & nfront-nass1,
974 & keep(116),
975 & iw(irow_l), perm,
976 & nvschur_k253 )
977 ELSE
978 nvschur_k253 = keep(253)
979 ENDIF
980 nfront8 = int(nfront,8)
981 laell8 = nfront8 * nfront8 + int(nass1,8)
982 CALL cmumps_parpivt1_set_max ( inode,
983 & a(poselt), laell8, keep,
984 & nfront, nass1, nvschur_k253,
985 & nb_postponed )
986 ENDIF
987 RETURN
float cmplx[2]
Definition pblas.h:136
#define mumps_abort
Definition VE_Metis.h:25
subroutine cmumps_ldlt_asm_niv12_ip(a, la, iafath, nfront, nass1, iacb, ncols, lcb, iw, nrows, nelim, etatass, cb_is_compressed)
Definition cfac_asm.F:318
subroutine cmumps_asm_slave_to_slave_end(n, inode, iw, liw, nbrows, step, ptrist, itloc, rhs_mumps, keep, keep8)
Definition cfac_asm.F:191
subroutine cmumps_asm_slave_arrowheads(inode, n, iw, liw, ioldps, a, la, poselt, keep, keep8, itloc, fils, ptraiw, ptrarw, intarr, dblarr, lintarr, ldblarr, rhs_mumps, lrgroups)
Definition cfac_asm.F:637
subroutine cmumps_restore_indices(n, ison, inode, iwposcb, pimaster, ptlust_s, iw, liw, step, keep, keep8)
Definition cfac_asm.F:523
subroutine cmumps_asm_slave_to_slave_init(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)
Definition cfac_asm.F:132
subroutine cmumps_ldlt_asm_niv12(a, la, son_a, iafath, nfront, nass1, ncols, lcb, iw, nrows, nelim, etatass, cb_is_compressed)
Definition cfac_asm.F:406
subroutine cmumps_asm_slave_to_slave(n, inode, iw, liw, a, la, nbrows, nbcols, rowlist, collist, valson, opassw, opeliw, step, ptrist, ptrast, itloc, rhs_mumps, fils, icntl, keep, keep8, myid, is_oftype5or6, lda_valson)
Definition cfac_asm.F:223
subroutine cmumps_asm_max(n, inode, iw, liw, a, la, ison, nbcols, valson, ptlust_s, ptrast, step, pimaster, opassw, iwposcb, myid, keep, keep8)
Definition cfac_asm.F:581
subroutine cmumps_update_parpiv_entries(inode, keep, parpiv, lparpiv, nb_postponed)
Definition cfac_asm.F:901
logical function cmumps_is_gemm_large_enough(m, n, k)
Definition cfac_asm.F:842
subroutine cmumps_parpivt1_set_nvschur_max(n, inode, iw, liw, a, la, keep, perm, ioldps, poselt, nfront, nass1, lr_activated, parpiv_t1, nb_postponed)
Definition cfac_asm.F:950
subroutine cmumps_set_parpivt1(inode, nfront, nass1, keep, lr_activated, parpiv_t1)
Definition cfac_asm.F:788
subroutine cmumps_parpivt1_set_max(inode, a, laell8, keep, nfront, nass1, nvschur_k253, nb_postponed)
Definition cfac_asm.F:854
subroutine cmumps_asm_slave_master(n, inode, iw, liw, a, la, ison, nbrows, nbcols, rowlist, valson, ptlust_s, ptrast, step, pimaster, opassw, iwposcb, myid, keep, keep8, is_oftype5or6, lda_valson)
Definition cfac_asm.F:19
logical function cmumps_is_trsm_large_enough(m, n)
Definition cfac_asm.F:830
#define min(a, b)
Definition macros.h:20
#define max(a, b)
Definition macros.h:21
subroutine get_cut(iwr, nass, ncb, lrgroups, npartscb, npartsass, cut)
Definition cana_lr.F:25
subroutine cmumps_dm_set_dynptr(cb_state, a, la, pamaster_or_ptrast, ixxd, ixxr, son_a, iachk, recsize)
subroutine cmumps_get_size_schur_in_front(n, ncb, size_schur, row_indices, perm, nvschur)
subroutine max_cluster(cut, cut_size, maxi_cluster)
Definition clr_core.F:1304
subroutine compute_blr_vcs(k472, ibcksz, maxsize, nass)
Definition lr_common.F:18