OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
ssol_fwd_aux.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 RECURSIVE SUBROUTINE smumps_traiter_message_solve
15 & ( bufr, lbufr, lbufr_bytes,
16 & msgtag, msgsou, myid, slavef, comm,
17 & n, nrhs, ipool, lpool, leaf,
18 & nbfin, nstk_s, iw, liw, a, la, ptrist,
19 & ptrfac, iwcb, liwcb,
20 & wcb, lwcb, poswcb,
21 & pleftwcb, posiwcb,
22 & ptricb,
23 & info, keep, keep8, dkeep, step, procnode_steps,
24 & rhscomp, lrhscomp, posinrhscomp_fwd
25 & , from_pp
26 & )
27 USE smumps_ooc
29 USE smumps_buf
30 IMPLICIT NONE
31 INTEGER lbufr, lbufr_bytes
32 INTEGER msgtag, msgsou, myid, slavef, comm
33 INTEGER liw
34 INTEGER(8), INTENT(IN) :: LA, lwcb
35 INTEGER n, nrhs, lpool, leaf, nbfin, lrhscomp
36 INTEGER liwcb, posiwcb
37 INTEGER(8) :: poswcb, pleftwcb
38 INTEGER info( 80 ), keep( 500)
39 INTEGER(8) keep8(150)
40 REAL, INTENT(INOUT) :: dkeep(230)
41 INTEGER bufr( lbufr )
42 INTEGER ipool( lpool ), nstk_s( n )
43 INTEGER iwcb( liwcb )
44 INTEGER iw( liw )
45 INTEGER ptricb(keep(28)),ptrist(keep(28))
46 INTEGER(8) :: ptrfac(keep(28))
47 INTEGER step(n)
48 INTEGER procnode_steps(keep(28))
49 REAL wcb( lwcb ), a( la )
50 REAL rhscomp( lrhscomp, nrhs )
51 INTEGER, intent(in) :: posinrhscomp_fwd(n)
52 LOGICAL, intent(in) :: from_pp
53 include 'mpif.h'
54 include 'mumps_tags.h'
55 INTEGER(8) :: ptrx, PTRY, ifr8
56 INTEGER ierr, k, jj, jbdeb, jbfin, nrhs_b
57 INTEGER :: iwhdlr, lda_slave
58 INTEGER :: mtype_slave
59 INTEGER finode, fpere, long, ncb, position, ncv, npiv
60 INTEGER pdest, i, iposinrhscomp
61 INTEGER j1
62 INTEGER(8) :: apos
63 LOGICAL dummy
64 LOGICAL flag
65 LOGICAL :: omp_flag
66 EXTERNAL mumps_procnode
67 INTEGER mumps_procnode
68 LOGICAL compress_panel, lr_activated
69 LOGICAL oocwrite_compatible_with_blr
70 REAL alpha, one
71 parameter(one = 1.0e0, alpha=-1.0e0)
72 include 'mumps_headers.h'
73 IF ( msgtag .EQ. racine_solve ) THEN
74 nbfin = nbfin - 1
75 IF ( nbfin .eq. 0 ) GOTO 270
76 ELSE IF (msgtag .EQ. contvec ) THEN
77 position = 0
78 CALL mpi_unpack( bufr, lbufr_bytes, position,
79 & finode, 1, mpi_integer, comm, ierr )
80 CALL mpi_unpack( bufr, lbufr_bytes, position,
81 & fpere, 1, mpi_integer, comm, ierr )
82 CALL mpi_unpack( bufr, lbufr_bytes, position,
83 & ncb, 1, mpi_integer, comm, ierr )
84 CALL mpi_unpack( bufr, lbufr_bytes, position,
85 & jbdeb, 1, mpi_integer, comm, ierr )
86 CALL mpi_unpack( bufr, lbufr_bytes, position,
87 & jbfin, 1, mpi_integer, comm, ierr )
88 CALL mpi_unpack( bufr, lbufr_bytes, position,
89 & long, 1, mpi_integer, comm, ierr )
90 nrhs_b = jbfin-jbdeb+1
91 IF ( ncb .eq. 0 ) THEN
92 ptricb(step(finode)) = -1
93 ELSE
94 IF ( ptricb(step(finode)) .EQ. 0 ) THEN
95 ptricb(step(finode)) = ncb + 1
96 END IF
97 IF ( ( posiwcb - long ) .LT. 0 ) THEN
98 info( 1 ) = -14
99 info( 2 ) = long
100 GOTO 260
101 END IF
102 IF ( poswcb - pleftwcb + 1_8 .LT.
103 & int(long,8) * int(nrhs_b,8)) THEN
104 info( 1 ) = -11
105 CALL mumps_set_ierror(pleftwcb-poswcb-1_8+
106 & int(long,8) * int(nrhs_b,8),
107 & info(2))
108 GOTO 260
109 END IF
110 IF (long .GT. 0) THEN
111 CALL mpi_unpack( bufr, lbufr_bytes, position,
112 & iwcb( 1 ),
113 & long, mpi_integer, comm, ierr )
114 DO k = 1, nrhs_b
115 CALL mpi_unpack( bufr, lbufr_bytes, position,
116 & wcb( pleftwcb ),
117 & long, mpi_real, comm, ierr )
118#if defined(__ve__)
119!NEC$ IVDEP
120#endif
121 DO i = 1, long
122 iposinrhscomp= abs(posinrhscomp_fwd(iwcb(i)))
123 rhscomp(iposinrhscomp,jbdeb+k-1) =
124 & rhscomp(iposinrhscomp,jbdeb+k-1) +
125 & wcb(pleftwcb+i-1)
126 ENDDO
127 END DO
128 ptricb(step(finode)) = ptricb(step(finode)) - long
129 ENDIF
130 END IF
131 IF ( ptricb(step(finode)) == 1 .OR.
132 & ptricb(step(finode)) == -1 ) THEN
133 nstk_s(step(fpere)) = nstk_s(step(fpere)) - 1
134 ptricb(step(finode)) = 0
135 END IF
136 IF ( nstk_s(step(fpere)) .EQ. 0 ) THEN
137 ipool( leaf ) = fpere
138 leaf = leaf + 1
139 IF ( leaf > lpool ) THEN
140 WRITE(*,*)
141 & 'Internal error 1 SMUMPS_TRAITER_MESSAGE_SOLVE',
142 & leaf, lpool
143 CALL mumps_abort()
144 END IF
145 ENDIF
146 ELSEIF ( msgtag .EQ. master2slave ) THEN
147 position = 0
148 CALL mpi_unpack( bufr, lbufr_bytes, position,
149 & finode, 1, mpi_integer, comm, ierr )
150 CALL mpi_unpack( bufr, lbufr_bytes, position,
151 & fpere, 1, mpi_integer, comm, ierr )
152 CALL mpi_unpack( bufr, lbufr_bytes, position,
153 & ncv, 1, mpi_integer, comm, ierr )
154 CALL mpi_unpack( bufr, lbufr_bytes, position,
155 & npiv, 1, mpi_integer, comm, ierr )
156 CALL mpi_unpack( bufr, lbufr_bytes, position,
157 & jbdeb, 1, mpi_integer, comm, ierr )
158 CALL mpi_unpack( bufr, lbufr_bytes, position,
159 & jbfin, 1, mpi_integer, comm, ierr )
160 nrhs_b = jbfin-jbdeb+1
161 ptry = pleftwcb
162 ptrx = pleftwcb + int(ncv,8) * int(nrhs_b,8)
163 pleftwcb = pleftwcb + int(npiv + ncv,8) * int(nrhs_b,8)
164 IF ( poswcb - pleftwcb + 1 .LT. 0 ) THEN
165 info(1) = -11
166 CALL mumps_set_ierror(-poswcb+pleftwcb-1_8,info(2))
167 GO TO 260
168 END IF
169 DO k=1, nrhs_b
170 CALL mpi_unpack( bufr, lbufr_bytes, position,
171 & wcb( ptry + (k-1) * ncv ), ncv,
172 & mpi_real, comm, ierr )
173 ENDDO
174 IF ( npiv .GT. 0 ) THEN
175 DO k=1, nrhs_b
176 CALL mpi_unpack( bufr, lbufr_bytes, position,
177 & wcb( ptrx + (k-1)*npiv ), npiv,
178 & mpi_real, comm, ierr )
179 END DO
180 END IF
181 lr_activated = (iw(ptrist(step(finode))+xxlr).GT.0)
182 compress_panel = (iw(ptrist(step(finode))+xxlr).GE.2)
183 oocwrite_compatible_with_blr =
184 & (.NOT.lr_activated.OR.(.NOT.compress_panel).OR.
185 & (keep(485).EQ.0)
186 & )
187 IF (keep(201).GT.0.AND.oocwrite_compatible_with_blr) THEN
189 & finode,ptrfac,keep,a,la,step,
190 & keep8,n,dummy,ierr)
191 IF(ierr.LT.0)THEN
192 info(1)=ierr
193 info(2)=0
194 GOTO 260
195 ENDIF
196 ENDIF
197 IF ( iw(ptrist(step(finode))+xxlr) .GE. 2 .AND.
198 & keep(485) .EQ. 1 ) THEN
199 iwhdlr = iw(ptrist(step(finode))+xxf)
200 mtype_slave = 1
201 CALL smumps_sol_slave_lr_u( finode, iwhdlr,
202 & -9999,
203 & wcb, lwcb,
204 & npiv, ncv,
205 & ptrx, ptry,
206 & jbdeb, jbfin,
207 & mtype_slave, keep, keep8,
208 & info(1), info(2) )
209 ELSE
210 apos = ptrfac(step(finode))
211 IF (keep(201) .EQ. 1) THEN
212 mtype_slave = 0
213 lda_slave = ncv
214 ELSE
215 mtype_slave = 1
216 lda_slave = npiv
217 ENDIF
219 & ( a, la, apos, npiv,
220 & lda_slave,
221 & ncv,
222 & nrhs_b, wcb, lwcb,
223 & ptrx, npiv,
224 & ptry, ncv,
225 & mtype_slave, keep, one )
226 ENDIF
227 IF ((keep(201).GT.0).AND.oocwrite_compatible_with_blr) THEN
228 CALL smumps_free_factors_for_solve(finode,ptrfac,
229 & keep(28),a,la,.true.,ierr)
230 IF(ierr.LT.0)THEN
231 info(1)=ierr
232 info(2)=0
233 GOTO 260
234 ENDIF
235 ENDIF
236 pleftwcb = pleftwcb - int(npiv,8) * int(nrhs_b,8)
237 pdest = mumps_procnode( procnode_steps(step(fpere)),
238 & keep(199) )
239 IF ( pdest .EQ. myid ) THEN
240 IF ( ptricb(step(finode)) .EQ. 0 ) THEN
241 ncb = iw( ptrist(step(finode)) + 2 + keep(ixsz) )
242 ptricb(step(finode)) = ncb + 1
243 END IF
244 j1 = ptrist(step(finode))+3+keep(ixsz)
245 omp_flag = .false.
246!$ OMP_FLAG = ( JBFIN-JBDEB+1.GE.KEEP(362) .AND.
247!$ & (NCV*(JBFIN-JBDEB+1) .GE. KEEP(363) ) )
248 IF (omp_flag) THEN
249!$OMP PARALLEL DO PRIVATE(I,JJ,IFR8,IPOSINRHSCOMP)
250 DO k=1, nrhs_b
251 ifr8 = ptry+int(k-1,8)*int(ncv,8)
252#if defined(__ve__)
253!NEC$ IVDEP
254#endif
255 DO i = 1,ncv
256 jj = iw(j1+i)
257 iposinrhscomp= abs(posinrhscomp_fwd(jj))
258 rhscomp(iposinrhscomp,jbdeb+k-1)=
259 & rhscomp(iposinrhscomp,jbdeb+k-1)
260 & + wcb(ifr8+int(i-1,8))
261 ENDDO
262 ENDDO
263!$OMP END PARALLEL DO
264 ELSE
265 DO k=1, nrhs_b
266 ifr8 = ptry+int(k-1,8)*int(ncv,8)
267#if defined(__ve__)
268!NEC$ IVDEP
269#endif
270 DO i = 1,ncv
271 jj = iw(j1+i)
272 iposinrhscomp= abs(posinrhscomp_fwd(jj))
273 rhscomp(iposinrhscomp,jbdeb+k-1)=
274 & rhscomp(iposinrhscomp,jbdeb+k-1)
275 & + wcb(ifr8+int(i-1,8))
276 ENDDO
277 ENDDO
278 ENDIF
279 ptricb(step(finode)) = ptricb(step(finode)) - ncv
280 IF ( ptricb( step( finode ) ) == 1 ) THEN
281 nstk_s(step(fpere)) = nstk_s(step(fpere)) - 1
282 ptricb(step(finode)) = 0
283 END IF
284 IF ( nstk_s(step(fpere)) .EQ. 0 ) THEN
285 ipool( leaf ) = fpere
286 leaf = leaf + 1
287 IF ( leaf > lpool ) THEN
288 WRITE(*,*)
289 & 'INTERNAL Error in SMUMPS_TRAITER_MESSAGE_SOLVE',
290 & leaf, lpool
291 CALL mumps_abort()
292 END IF
293 ENDIF
294 ELSE
295 210 CONTINUE
296 CALL smumps_buf_send_vcb( nrhs_b, finode, fpere,
297 & iw(ptrist(step( finode )) + 2 + keep(ixsz) ), ncv,ncv,
298 & iw(ptrist(step(finode))+4+ keep(ixsz) ),
299 & wcb( ptry ), jbdeb, jbfin,
300 & rhscomp, 1, 1, -9999, -9999,
301 & keep, pdest, contvec, comm, ierr )
302 IF ( ierr .EQ. -1 ) THEN
303 CALL smumps_solve_recv_and_treat( .false., flag,
304 & bufr, lbufr, lbufr_bytes,
305 & myid, slavef, comm,
306 & n, nrhs, ipool, lpool, leaf,
307 & nbfin, nstk_s, iw, liw, a, la, ptrist, ptrfac,
308 & iwcb, liwcb,
309 & wcb, lwcb, poswcb, pleftwcb, posiwcb,
310 & ptricb, info, keep,keep8, dkeep, step,
311 & procnode_steps,
312 & rhscomp, lrhscomp, posinrhscomp_fwd
313 & , from_pp
314 & )
315 IF ( info( 1 ) .LT. 0 ) GOTO 270
316 GOTO 210
317 ELSE IF ( ierr .EQ. -2 ) THEN
318 info( 1 ) = -17
319 info( 2 ) = ( ncv + 4 ) * keep( 34 ) +
320 & ncv * keep( 35 )
321 GOTO 260
322 ELSE IF ( ierr .EQ. -3 ) THEN
323 info( 1 ) = -20
324 info( 2 ) = ( ncv + 4 ) * keep( 34 ) +
325 & ncv * keep( 35 )
326 END IF
327 END IF
328 pleftwcb = pleftwcb - int(ncv,8) * int(nrhs_b,8)
329 ELSEIF ( msgtag .EQ. terreur ) THEN
330 info(1) = -001
331 info(2) = msgsou
332 GOTO 270
333 ELSE IF ( (msgtag.EQ.update_load).OR.
334 & (msgtag.EQ.tag_dummy) ) THEN
335 GO TO 270
336 ELSE
337 info(1)=-100
338 info(2)=msgtag
339 GO TO 260
340 ENDIF
341 GO TO 270
342 260 CONTINUE
343 CALL smumps_bdc_error( myid, slavef, comm, keep )
344 270 CONTINUE
345 RETURN
346 END SUBROUTINE smumps_traiter_message_solve
347 SUBROUTINE smumps_solve_node_fwd( INODE,
348 & LASTFSL0STA, LASTFSL0DYN,
349 & BUFR, LBUFR, LBUFR_BYTES,
350 & MYID, SLAVEF, COMM,
351 & N, IPOOL, LPOOL, LEAF,
352 & NBFIN, NSTK_S,
353 & IWCB, LIWCB,
354 & WCB, LWCB, A, LA, IW, LIW,
355 & NRHS, POSWCB, PLEFTWCB, POSIWCB,
356 & PTRICB, PTRIST, PTRFAC, PROCNODE_STEPS,
357 & FILS, STEP, FRERE, DAD,
358 & INFO, KEEP,KEEP8, DKEEP, RHS_ROOT, LRHS_ROOT, MTYPE,
359 & RHSCOMP, LRHSCOMP, POSINRHSCOMP_FWD,
360 &
361 & ISTEP_TO_INIV2, TAB_POS_IN_PERE
362 & , RHS_BOUNDS, LRHS_BOUNDS, DO_NBSPARSE, FROM_PP
363 & , ERROR_WAS_BROADCASTED
364 & )
365 USE smumps_sol_lr
366!$ USE SMUMPS_SOL_L0OMP_M, ONLY: LOCK_FOR_SCATTER
368 USE smumps_ooc
369 USE smumps_buf
370 IMPLICIT NONE
371 INTEGER MTYPE
372 INTEGER, INTENT( IN ) :: INODE, LASTFSL0STA, LASTFSL0DYN
373 INTEGER LBUFR, LBUFR_BYTES
374 INTEGER MYID, SLAVEF, COMM
375 INTEGER LIWCB, LIW, POSIWCB
376 INTEGER(8) :: POSWCB, PLEFTWCB, LWCB
377 INTEGER(8) :: LA
378 INTEGER N, LPOOL, LEAF, NBFIN
379 INTEGER INFO( 80 ), KEEP( 500)
380 INTEGER(8) KEEP8(150)
381 REAL, INTENT(INOUT) :: DKEEP(230)
382 INTEGER BUFR( LBUFR )
383 INTEGER IPOOL( LPOOL ), NSTK_S(KEEP(28))
384 INTEGER IWCB( LIWCB ), IW( LIW )
385 INTEGER NRHS
386 REAL WCB( LWCB )
387 REAL :: A( LA )
388 INTEGER(8) :: LRHS_ROOT
389 REAL RHS_ROOT( LRHS_ROOT )
390 INTEGER PTRICB(KEEP(28)), PTRIST(KEEP(28))
391 INTEGER(8) :: PTRFAC(KEEP(28))
392 INTEGER PROCNODE_STEPS(KEEP(28))
393 INTEGER FILS( N ), STEP( N ), FRERE(KEEP(28)), DAD(KEEP(28))
394 INTEGER ISTEP_TO_INIV2(KEEP(71)),
395 & tab_pos_in_pere(slavef+2,max(1,keep(56)))
396 INTEGER POSINRHSCOMP_FWD(N), LRHSCOMP
397 REAL RHSCOMP(LRHSCOMP, NRHS)
398 LOGICAL, intent(in) :: DO_NBSPARSE
399 INTEGER, intent(in) :: LRHS_BOUNDS
400 INTEGER, intent(in) :: RHS_BOUNDS(LRHS_BOUNDS)
401 LOGICAL, intent(in) :: FROM_PP
402 LOGICAL, intent(out) :: ERROR_WAS_BROADCASTED
404 INTEGER MUMPS_PROCNODE
405 REAL ALPHA,ONE,ZERO
406 parameter(zero=0.0e0, one = 1.0e0, alpha=-1.0e0)
407 INTEGER :: IWHDLR
408 INTEGER JBDEB, JBFIN, NRHS_B
409 INTEGER LDADIAG
410 INTEGER(8) :: APOS, APOS1, IFR8, IFR_ini8
411 INTEGER I, J, K, IPOS, NSLAVES, J1, J2, J3, FPERE, FPERE_MAPPING,
412 & NPIV, NCB, LIELL, JJ, NELIM, IERR
413 INTEGER(8) :: PCB_COURANT, PPIV_COURANT, PPIV_PANEL, PCB_PANEL
414 INTEGER IPOSINRHSCOMP_TMP
415 INTEGER Effective_CB_Size, NUPDATE, ISLAVE, PDEST, FirstIndex
416 LOGICAL FLAG
417 INTEGER :: NUPDATE_NONCRITICAL, IPOSINRHSCOMPLASTFSDYN
418 LOGICAL :: OMP_FLAG
419 include 'mumps_headers.h'
420 INTEGER(8) :: APOSDEB
421 INTEGER TempNROW, TempNCOL, PANEL_SIZE,
422 & jfin, nbj, nupdate_panel,
423 & typef
424 INTEGER LD_WCBPIV
425 INTEGER LD_WCBCB
426 LOGICAL :: LDEQLIELLPANEL
427 LOGICAL :: CBINITZERO
428 INTEGER LDAJ, LDAJ_FIRST_PANEL
429 INTEGER LDAtemp
430 LOGICAL COMPRESS_PANEL, LR_ACTIVATED
431 LOGICAL OOCWRITE_COMPATIBLE_WITH_BLR
432 INTEGER TMP_NBPANELS,
433 & i_pivrptr, i_pivr, ipanel
434 LOGICAL MUST_BE_PERMUTED
435 INTEGER :: SIZEBLOCK, NB, JCourant, NB_LOCK
436 include 'mpif.h'
437 include 'mumps_tags.h'
438 INTEGER DUMMY( 1 )
439 ERROR_WAS_BROADCASTED = .false.
440 dummy(1)=1
441 lr_activated = (iw(ptrist(step(inode))+xxlr).GT.0)
442 compress_panel = (iw(ptrist(step(inode))+xxlr).GE.2)
443 oocwrite_compatible_with_blr =
444 & (.NOT.lr_activated.OR.(.NOT.compress_panel).OR.
445 & (keep(485).EQ.0)
446 & )
447 IF (do_nbsparse) THEN
448 jbdeb= rhs_bounds(2*step(inode)-1)
449 jbfin= rhs_bounds(2*step(inode))
450 ELSE
451 jbdeb = 1
452 jbfin = nrhs
453 ENDIF
454 nrhs_b = jbfin-jbdeb+1
455 IF (do_nbsparse) THEN
456 if (jbdeb.GT.jbfin) then
457 write(6,*) " Internal error 1 in nbsparse :",
458 & jbdeb, jbfin
459 CALL mumps_abort()
460 endif
461 IF (jbdeb.LT.1 .OR. jbdeb.GT.nrhs .or.
462 & jbfin.LT.1 .OR. jbfin.GT.nrhs ) THEN
463 write(6,*) " Internal error 2 in nbsparse :",
464 & jbdeb, jbfin
465 CALL mumps_abort()
466 endif
467 ENDIF
468 IF ( inode .eq. keep( 38 ) .OR. inode .eq.keep( 20 ) ) THEN
469 liell = iw( ptrist( step(inode)) + 3 + keep(ixsz))
470 npiv = liell
471 nelim = 0
472 nslaves = 0
473 ipos = ptrist( step(inode)) + 5 + keep(ixsz)
474 ELSE
475 ipos = ptrist(step(inode)) + 2 + keep(ixsz)
476 liell = iw(ipos-2)+iw(ipos+1)
477 nelim = iw(ipos-1)
478 nslaves = iw( ptrist(step(inode)) + 5 + keep(ixsz) )
479 ipos = ipos + 1
480 npiv = iw(ipos)
481 ipos = ipos + 1
482 IF ((keep(201).GT.0).AND.oocwrite_compatible_with_blr) THEN
484 & inode,ptrfac,keep,a,la,step,
485 & keep8,n,must_be_permuted,ierr)
486 IF(ierr.LT.0)THEN
487 info(1)=ierr
488 info(2)=0
489 error_was_broadcasted = .false.
490 GOTO 270
491 ENDIF
492 IF (keep(201).EQ.1 .AND. keep(50).NE.1) THEN
494 & iw(ipos+1+2*liell+1+nslaves),
495 & must_be_permuted )
496 ENDIF
497 ENDIF
498 nslaves = iw( ptrist(step(inode)) + 5 + keep(ixsz))
499 ipos = ipos + 1 + nslaves
500 END IF
501 IF ( mtype .EQ. 1 .OR. keep(50) .NE. 0 ) THEN
502 j1 = ipos + 1
503 j2 = ipos + liell
504 j3 = ipos + npiv
505 ELSE
506 j1 = ipos + liell + 1
507 j2 = ipos + 2 * liell
508 j3 = ipos + liell + npiv
509 END IF
510 ncb = liell-npiv
511 IF (keep(50).NE.0) THEN
512 IF ( keep(459) .GT. 1 ) THEN
513 ldadiag = -99999
514 ELSE
515 ldadiag = npiv
516 ENDIF
517 ELSE
518 ldadiag = liell
519 ENDIF
520 IF ( inode .eq. keep( 38 ) .OR. inode .eq. keep(20) ) THEN
521 ifr8 = 0_8
522 iposinrhscomp_tmp = posinrhscomp_fwd(iw(j1))
523 ifr_ini8 = ifr8
524 omp_flag = .false.
525!$ OMP_FLAG = ( JBFIN-JBDEB+1.GE.KEEP(362) .AND.
526!$ & (J3-J1+1)*(JBFIN-JBDEB+1) .GE. KEEP(363) )
527 IF (omp_flag) THEN
528!$OMP PARALLEL DO PRIVATE(IFR8,JJ)
529 DO k=jbdeb,jbfin
530 ifr8 = ifr_ini8 + int(k-1,8)*int(npiv,8)
531 DO jj = j1, j3
532 rhs_root(ifr8+int(jj-j1+1,8)) =
533 & rhscomp(iposinrhscomp_tmp+jj-j1,k)
534 ENDDO
535 ENDDO
536!$OMP END PARALLEL DO
537 ELSE
538 DO k=jbdeb,jbfin
539 ifr8 = ifr_ini8 + int(k-1,8)*int(npiv,8)
540 DO jj = j1, j3
541 rhs_root(ifr8+int(jj-j1+1,8)) =
542 & rhscomp(iposinrhscomp_tmp+jj-j1,k)
543 ENDDO
544 ENDDO
545 ENDIF
546 IF ( npiv .LT. liell ) THEN
547 WRITE(*,*) ' Internal error 1 in SMUMPS_SOLVE_NODE_FWD',
548 & npiv, liell
549 CALL mumps_abort()
550 END IF
551 GO TO 270
552 END IF
553 apos = ptrfac(step(inode))
554 IF ( (keep(201).EQ.1).AND.oocwrite_compatible_with_blr ) THEN
555 IF (mtype.EQ.1) THEN
556 IF ((mtype.EQ.1).AND.nslaves.NE.0) THEN
557 tempnrow= npiv+nelim
558 tempncol= npiv
559 ldaj_first_panel=tempnrow
560 ELSE
561 tempnrow= liell
562 tempncol= npiv
563 ldaj_first_panel=tempnrow
564 ENDIF
565 typef=typef_l
566 ELSE
567 tempncol= liell
568 tempnrow= npiv
569 ldaj_first_panel=tempncol
570 typef= typef_u
571 ENDIF
572 panel_size = smumps_ooc_panel_size( ldaj_first_panel )
573 ENDIF
574 ppiv_courant = pleftwcb
575 pleftwcb = pleftwcb + int(liell,8) * int(nrhs_b,8)
576 IF ( poswcb - pleftwcb + 1_8 .LT. 0 ) THEN
577 info(1) = -11
578 CALL mumps_set_ierror(pleftwcb-poswcb-1_8, info(2))
579 error_was_broadcasted = .false.
580 GOTO 270
581 END IF
582 IF (keep(201) .EQ. 1 .AND. oocwrite_compatible_with_blr) THEN
583 ldeqliellpanel = .true.
584 ld_wcbpiv = liell
585 ld_wcbcb = liell
586 pcb_courant = ppiv_courant + npiv
587 ELSE
588 ldeqliellpanel = .false.
589 ld_wcbpiv = npiv
590 ld_wcbcb = ncb
591 pcb_courant = ppiv_courant + int(npiv,8)*int(nrhs_b,8)
592 ENDIF
593 fpere = dad(step(inode))
594 IF ( fpere .NE. 0 ) THEN
595 fpere_mapping = mumps_procnode( procnode_steps(step(fpere)),
596 & keep(199) )
597 ELSE
598 fpere_mapping = -1
599 ENDIF
600 IF ( lastfsl0dyn .LE. n ) THEN
601 cbinitzero = .true.
602 ELSE IF ( fpere_mapping .EQ. myid ) THEN
603 cbinitzero = .true.
604 ELSE
605 cbinitzero = .false.
606 ENDIF
608 & npiv, ncb, liell, cbinitzero, ldeqliellpanel,
609 & rhscomp(1, jbdeb), lrhscomp, nrhs_b,
610 & posinrhscomp_fwd, n,
611 & wcb(ppiv_courant),
612 & iw, liw, j1, j3, j2, keep, dkeep)
613 IF ( npiv .NE. 0 ) THEN
614 IF ((keep(201).EQ.1).AND.oocwrite_compatible_with_blr) THEN
615 aposdeb = apos
616 j = 1
617 ipanel = 0
618 10 CONTINUE
619 ipanel = ipanel + 1
620 jfin = min(j+panel_size-1, npiv)
621 IF (iw(ipos+ liell + jfin) < 0) THEN
622 jfin=jfin+1
623 ENDIF
624 nbj = jfin-j+1
625 ldaj = ldaj_first_panel-j+1
626 IF ( (keep(50).NE.1).AND. must_be_permuted ) THEN
627 CALL smumps_get_ooc_perm_ptr(typef, tmp_nbpanels,
628 & i_pivrptr, i_pivr, ipos+1+2*liell, iw, liw)
629 IF (npiv.EQ.(iw(i_pivrptr+ipanel-1)-1)) THEN
630 must_be_permuted=.false.
631 ELSE
633 & iw( i_pivr+ iw(i_pivrptr+ipanel-1)-
634 & iw(i_pivrptr)),
635 & npiv-iw(i_pivrptr+ipanel-1)+1,
636 & iw(i_pivrptr+ipanel-1)-1,
637 & a(aposdeb),
638 & ldaj, nbj, j-1 )
639 ENDIF
640 ENDIF
641 nupdate_panel = ldaj - nbj
642 ppiv_panel = ppiv_courant+int(j-1,8)
643 pcb_panel = ppiv_panel+int(nbj,8)
644 apos1 = aposdeb+int(nbj,8)
645 IF (mtype.EQ.1) THEN
646#if defined(MUMPS_USE_BLAS2)
647 IF ( nrhs_b == 1 ) THEN
648 CALL strsv( 'L', 'N', 'U', nbj, a(aposdeb), ldaj,
649 & wcb(ppiv_panel), 1 )
650 IF (nupdate_panel.GT.0) THEN
651 CALL sgemv('N', nupdate_panel,nbj,alpha, a(apos1),
652 & ldaj, wcb(ppiv_panel), 1, one,
653 & wcb(pcb_panel), 1)
654 ENDIF
655 ELSE
656#endif
657 CALL strsm( 'L','L','N','U', nbj, nrhs_b, one,
658 & a(aposdeb), ldaj, wcb(ppiv_panel),
659 & liell )
660 IF (nupdate_panel.GT.0) THEN
661 CALL sgemm('N', 'N', nupdate_panel, nrhs_b, nbj,
662 & alpha,
663 & a(apos1), ldaj, wcb(ppiv_panel), liell, one,
664 & wcb(pcb_panel), liell)
665 ENDIF
666#if defined(MUMPS_USE_BLAS2)
667 ENDIF
668#endif
669 ELSE
670#if defined(MUMPS_USE_BLAS2)
671 IF (nrhs_b == 1) THEN
672 CALL strsv( 'L', 'N', 'N', nbj, a(aposdeb), ldaj,
673 & wcb(ppiv_panel), 1 )
674 IF (nupdate_panel.GT.0) THEN
675 CALL sgemv('N',nupdate_panel, nbj, alpha, a(apos1),
676 & ldaj, wcb(ppiv_panel), 1,
677 & one, wcb(pcb_panel), 1 )
678 ENDIF
679 ELSE
680#endif
681 CALL strsm('L','L','N','N',nbj, nrhs_b, one,
682 & a(aposdeb), ldaj, wcb(ppiv_panel),
683 & liell)
684 IF (nupdate_panel.GT.0) THEN
685 CALL sgemm('N', 'N', nupdate_panel, nrhs_b, nbj,
686 & alpha,
687 & a(apos1), ldaj, wcb(ppiv_panel), liell, one,
688 & wcb(pcb_panel), liell)
689 ENDIF
690#if defined(MUMPS_USE_BLAS2)
691 ENDIF
692#endif
693 ENDIF
694 aposdeb = aposdeb+int(ldaj,8)*int(nbj,8)
695 j=jfin+1
696 IF ( j .LE. npiv ) GOTO 10
697 ELSE
698 IF ( iw(ptrist(step(inode))+xxlr) .GE. 2 .AND.
699 & keep(485) .EQ. 1 ) THEN
700 iwhdlr = iw(ptrist(step(inode))+xxf)
702 & inode, n, iwhdlr, npiv, nslaves,
703 & iw, ipos, liw,
704 & liell, wcb, lwcb,
705 & ld_wcbpiv, ld_wcbcb,
706 & ppiv_courant, pcb_courant,
707 & rhscomp, lrhscomp, nrhs,
708 & posinrhscomp_fwd, jbdeb, jbfin,
709 & mtype, keep, keep8, oocwrite_compatible_with_blr,
710 & info(1), info(2) )
711 IF (info(1).LT.0) THEN
712 error_was_broadcasted = .false.
713 GOTO 270
714 ENDIF
715 ELSE IF ( keep(459) .GT. 1 .AND. keep(50) .NE. 0 ) THEN
717 & a, la, apos,
718 & npiv, iw(ipos+liell+1),
719 & nrhs_b, wcb, lwcb, ld_wcbpiv,
720 & ppiv_courant, mtype, keep)
721 ELSE
723 & a, la, apos,
724 & npiv, ldadiag,
725 & nrhs_b, wcb, lwcb, ld_wcbpiv,
726 & ppiv_courant, mtype, keep)
727 ENDIF
728 END IF
729 END IF
730 ncb = liell - npiv
731 IF ( mtype .EQ. 1 ) THEN
732 IF ( nslaves .EQ. 0 .OR. npiv .eq. 0 ) THEN
733 nupdate = ncb
734 ELSE
735 nupdate = nelim
736 END IF
737 IF (keep(459) .GT. 1 .AND. keep(50) .NE. 0) THEN
738 CALL mumps_geti8(apos1, iw(ptrist(step(inode))+xxr))
739 apos1 = apos + apos1 - int(npiv,8)*int(nupdate,8)
740 ELSE
741 apos1 = apos + int(npiv,8) * int(ldadiag,8)
742 ENDIF
743 ELSE
744 apos1 = apos + int(npiv,8)
745 nupdate = ncb
746 END IF
747 IF (keep(201).NE.1) THEN
748 IF ( iw(ptrist(step(inode))+xxlr) .LT. 2 .OR.
749 & keep(485).EQ.0) THEN
750 IF (mtype .EQ. 1) THEN
751 ldatemp = npiv
752 ELSE
753 ldatemp = liell
754 ENDIF
756 & (a, la, apos1,
757 & npiv, ldatemp, nupdate,
758 & nrhs_b, wcb, lwcb, ppiv_courant, ld_wcbpiv,
759 & pcb_courant, ld_wcbcb,
760 & mtype, keep, one)
761 ENDIF
762 END IF
763 IF ( iw(ptrist(step(inode))+xxlr) .LT. 2 .OR.
764 & keep(485).EQ.0) THEN
765 IF (keep(201) .GT. 0 .AND. oocwrite_compatible_with_blr) THEN
767 & inode, n, npiv, liell, nelim, nslaves,
768 & ppiv_courant,
769 & iw, ipos, liw,
770 & a, la, apos,
771 & wcb, lwcb, ld_wcbpiv,
772 & rhscomp, lrhscomp, nrhs,
773 & posinrhscomp_fwd, jbdeb, jbfin,
774 & mtype, keep, oocwrite_compatible_with_blr,
775 & .false.
776 & )
777 ELSE
779 & inode, n, npiv, liell, nelim, nslaves,
780 & ppiv_courant,
781 & iw, ipos, liw,
782 & a, la, apos,
783 & wcb, lwcb, ld_wcbpiv,
784 & rhscomp, lrhscomp, nrhs,
785 & posinrhscomp_fwd, jbdeb, jbfin,
786 & mtype, keep, oocwrite_compatible_with_blr,
787 & .false.
788 & )
789 ENDIF
790 ENDIF
791 IF ((keep(201).EQ.1).AND.oocwrite_compatible_with_blr)
792 &THEN
793 CALL smumps_free_factors_for_solve(inode,ptrfac,keep(28),
794 & a,la,.true.,ierr)
795 IF(ierr.LT.0)THEN
796 info(1)=ierr
797 info(2)=0
798 error_was_broadcasted = .false.
799 GOTO 270
800 ENDIF
801 END IF
802 IF ( fpere .EQ. 0 ) THEN
803 pleftwcb = pleftwcb - int(liell,8) *int(nrhs_b,8)
804 GOTO 270
805 ENDIF
806 IF ( nupdate .NE. 0 .OR. ncb.EQ.0 ) THEN
807 IF (mumps_procnode(procnode_steps(step(fpere)),
808 & keep(199)) .EQ. myid) THEN
809 IF ( ncb .ne. 0 ) THEN
810 ptricb(step(inode)) = ncb + 1
811 nupdate_noncritical = nupdate
812 IF (lastfsl0dyn .LE. n) THEN
813 IF ( lastfsl0dyn .EQ. 0 ) THEN
814 iposinrhscomplastfsdyn = 0
815 ELSE
816 iposinrhscomplastfsdyn =
817 & abs(posinrhscomp_fwd(lastfsl0dyn))
818 ENDIF
819 DO i = 1, nupdate
820 IF ( abs(posinrhscomp_fwd( iw(j3+i) )) .GT.
821 & iposinrhscomplastfsdyn ) THEN
822 IF (abs(step(iw(j3+i))) .GT.
823 & abs(step( lastfsl0sta))
824 & .OR. keep(261) .NE. 1) THEN
825 nupdate_noncritical = i - 1
826 EXIT
827 ENDIF
828 ENDIF
829 ENDDO
830 ENDIF
831 omp_flag = .false.
832!$ OMP_FLAG = ( NRHS_B.GE.KEEP(362) .AND.
833!$ & (NUPDATE*NRHS_B .GE. KEEP(363)) )
834 IF (omp_flag) THEN
835!$OMP PARALLEL DO PRIVATE(I,IFR8,IPOSINRHSCOMP_TMP)
836 DO k = jbdeb, jbfin
837 ifr8 = pcb_courant + int(k-jbdeb,8)*int(ld_wcbcb,8)
838#if defined(__ve__)
839!NEC$ IVDEP
840#endif
841 DO i = 1, nupdate_noncritical
842 iposinrhscomp_tmp =
843 & abs(posinrhscomp_fwd(iw(j3 + i)))
844 rhscomp( iposinrhscomp_tmp, k ) =
845 & rhscomp( iposinrhscomp_tmp, k )
846 & + wcb(ifr8 + int(i-1,8))
847 ENDDO
848 ENDDO
849!$OMP END PARALLEL DO
850 ELSE
851 DO k = jbdeb, jbfin
852 ifr8 = pcb_courant + int(k-jbdeb,8)*int(ld_wcbcb,8)
853#if defined(__ve__)
854!NEC$ IVDEP
855#endif
856 DO i = 1, nupdate_noncritical
857 iposinrhscomp_tmp =
858 & abs(posinrhscomp_fwd(iw(j3 + i)))
859 rhscomp( iposinrhscomp_tmp, k ) =
860 & rhscomp( iposinrhscomp_tmp, k )
861 & + wcb(ifr8 + int(i-1,8))
862 ENDDO
863 ENDDO
864 ENDIF
865 IF ( cbinitzero ) THEN
866 IF ( nupdate .NE. nupdate_noncritical) THEN
867 nb_lock = 1
868 IF (.NOT.do_nbsparse.AND.(keep(400).GT.1)) THEN
869 nb_lock = min(keep(400),nb_lock_max)
870 ENDIF
871 sizeblock = (jbfin-jbdeb+1+nb_lock-1) / nb_lock
872 DO nb = 1, nb_lock
873 jcourant = jbdeb+sizeblock*(nb-1)
874!$ CALL OMP_SET_LOCK(LOCK_FOR_SCATTER(NB))
875 DO k = jcourant, min(jbfin,jcourant+sizeblock-1)
876 ifr8 = pcb_courant + int(k-jbdeb,8)*int(ld_wcbcb,8)
877#if defined(__ve__)
878!NEC$ IVDEP
879#endif
880 DO i = nupdate_noncritical+1, nupdate
881 iposinrhscomp_tmp =
882 & abs(posinrhscomp_fwd(iw(j3 + i)))
883 rhscomp( iposinrhscomp_tmp, k ) =
884 & rhscomp( iposinrhscomp_tmp, k )
885 & + wcb(ifr8 + int(i-1,8))
886 ENDDO
887 ENDDO
888!$ CALL OMP_UNSET_LOCK(LOCK_FOR_SCATTER(NB))
889 ENDDO
890 ENDIF
891 ENDIF
892 ptricb(step( inode )) = ptricb(step( inode )) - nupdate
893 ELSE
894 ptricb(step( inode )) = -1
895 ENDIF
896 ELSE
897 210 CONTINUE
898 CALL smumps_buf_send_vcb( nrhs_b, inode, fpere,
899 & ncb, ld_wcbcb,
900 & nupdate,
901 & iw( j3 + 1 ), wcb( pcb_courant ), jbdeb, jbfin,
902 & rhscomp, 1, 1, -9999, -9999,
903 & keep,
904 & mumps_procnode(procnode_steps(step(fpere)), keep(199)),
905 & contvec,
906 & comm, ierr )
907 IF ( ierr .EQ. -1 ) THEN
908 CALL smumps_solve_recv_and_treat( .false., flag,
909 & bufr, lbufr, lbufr_bytes,
910 & myid, slavef, comm,
911 & n, nrhs, ipool, lpool, leaf,
912 & nbfin, nstk_s, iw, liw, a, la, ptrist, ptrfac,
913 & iwcb, liwcb,
914 & wcb, lwcb, poswcb, pleftwcb, posiwcb,
915 & ptricb, info, keep,keep8, dkeep, step,
916 & procnode_steps,
917 & rhscomp, lrhscomp, posinrhscomp_fwd
918 & , from_pp
919 & )
920 IF ( info( 1 ) .LT. 0 ) THEN
921 error_was_broadcasted = .true.
922 GOTO 270
923 ENDIF
924 GOTO 210
925 ELSE IF ( ierr .EQ. -2 ) THEN
926 info( 1 ) = -17
927 info( 2 ) = nupdate * keep( 35 ) +
928 & ( nupdate + 3 ) * keep( 34 )
929 error_was_broadcasted = .false.
930 GOTO 270
931 ELSE IF ( ierr .EQ. -3 ) THEN
932 info( 1 ) = -20
933 info( 2 ) = nupdate * keep( 35 ) +
934 & ( nupdate + 3 ) * keep( 34 )
935 error_was_broadcasted = .false.
936 GOTO 270
937 END IF
938 ENDIF
939 END IF
940 IF ( nslaves .NE. 0 .AND. mtype .eq. 1
941 & .and. npiv .NE. 0 ) THEN
942 DO islave = 1, nslaves
943 pdest = iw( ptrist(step(inode)) + 5 + islave +keep(ixsz))
945 & keep,keep8, inode, step, n, slavef,
946 & istep_to_iniv2, tab_pos_in_pere,
947 & islave, ncb - nelim,
948 & nslaves,
949 & effective_cb_size, firstindex )
950 222 CONTINUE
951 CALL smumps_buf_send_master2slave( nrhs_b,
952 & inode, fpere,
953 & effective_cb_size, ld_wcbcb, ld_wcbpiv, npiv,
954 & jbdeb, jbfin,
955 & wcb( pcb_courant + nelim + firstindex - 1 ),
956 & wcb( ppiv_courant ),
957 & pdest, comm, keep, ierr )
958 IF ( ierr .EQ. -1 ) THEN
959 CALL smumps_solve_recv_and_treat( .false., flag,
960 & bufr, lbufr, lbufr_bytes,
961 & myid, slavef, comm,
962 & n, nrhs, ipool, lpool, leaf,
963 & nbfin, nstk_s, iw, liw, a, la, ptrist,ptrfac,
964 & iwcb, liwcb,
965 & wcb, lwcb, poswcb, pleftwcb, posiwcb,
966 & ptricb, info, keep,keep8, dkeep, step,
967 & procnode_steps,
968 & rhscomp, lrhscomp, posinrhscomp_fwd
969 & , from_pp
970 & )
971 IF ( info( 1 ) .LT. 0 ) THEN
972 error_was_broadcasted = .true.
973 GOTO 270
974 ENDIF
975 GOTO 222
976 ELSE IF ( ierr .EQ. -2 ) THEN
977 info( 1 ) = -17
978 info( 2 ) = (npiv+effective_cb_size)*nrhs_b*keep(35) +
979 & 6 * keep( 34 )
980 error_was_broadcasted = .false.
981 GOTO 270
982 ELSE IF ( ierr .EQ. -3 ) THEN
983 info( 1 ) = -20
984 info( 2 ) = (npiv+effective_cb_size)*nrhs_b*keep(35) +
985 & 6 * keep( 34 )
986 error_was_broadcasted = .false.
987 GOTO 270
988 END IF
989 END DO
990 END IF
991 pleftwcb = pleftwcb - int(liell,8)*int(nrhs_b,8)
992 270 CONTINUE
993 RETURN
994 END SUBROUTINE smumps_solve_node_fwd
995 RECURSIVE SUBROUTINE smumps_solve_recv_and_treat( BLOQ, FLAG,
996 & BUFR, LBUFR, LBUFR_BYTES,
997 & MYID, SLAVEF, COMM,
998 & N, NRHS, IPOOL, LPOOL, LEAF,
999 & NBFIN, NSTK_S, IW, LIW, A, LA, PTRIST,PTRFAC,
1000 & IWCB, LIWCB,
1001 & WCB, LWCB, POSWCB,
1002 & PLEFTWCB, POSIWCB,
1003 & PTRICB, INFO, KEEP,KEEP8, DKEEP, STEP, PROCNODE_STEPS,
1004 & RHSCOMP, LRHSCOMP, POSINRHSCOMP_FWD
1005 & , FROM_PP
1006 & )
1007 IMPLICIT NONE
1008 LOGICAL bloq
1009 INTEGER lbufr, lbufr_bytes
1010 INTEGER myid, slavef, comm
1011 INTEGER n, nrhs, lpool, leaf, nbfin
1012 INTEGER liwcb, posiwcb
1013 INTEGER(8) :: poswcb, pleftwcb
1014 INTEGER liw
1015 INTEGER(8), INTENT(IN) :: la, lwcb
1016 INTEGER info( 80 ), keep( 500)
1017 INTEGER(8) keep8(150)
1018 REAL, INTENT(INOUT) :: dkeep(230)
1019 INTEGER bufr( lbufr ), ipool(lpool)
1020 INTEGER nstk_s( keep(28) )
1021 INTEGER iwcb( liwcb )
1022 INTEGER iw( liw )
1023 REAL wcb( lwcb ), a( la )
1024 INTEGER ptricb(keep(28)), ptrist(keep(28))
1025 INTEGER(8) :: ptrfac(keep(28))
1026 INTEGER step(n)
1027 INTEGER procnode_steps(keep(28))
1028 LOGICAL flag
1029 INTEGER lrhscomp, posinrhscomp_fwd(n)
1030 REAL rhscomp(lrhscomp,nrhs)
1031 LOGICAL, intent(in) :: from_pp
1032 include 'mpif.h'
1033 include 'mumps_tags.h'
1034 INTEGER :: ierr
1035 INTEGER :: status(mpi_status_size)
1036 INTEGER msgsou, msgtag, msglen
1037 flag = .false.
1038 IF ( bloq ) THEN
1039 flag = .false.
1040 CALL mpi_probe( mpi_any_source, mpi_any_tag,
1041 & comm, status, ierr )
1042 flag = .true.
1043 ELSE
1044 CALL mpi_iprobe( mpi_any_source, mpi_any_tag, comm,
1045 & flag, status, ierr )
1046 END IF
1047 IF ( flag ) THEN
1048 keep(266) = keep(266) -1
1049 msgsou = status( mpi_source )
1050 msgtag = status( mpi_tag )
1051 CALL mpi_get_count( status, mpi_packed, msglen, ierr )
1052 IF ( msglen .GT. lbufr_bytes ) THEN
1053 info(1) = -20
1054 info(2) = msglen
1055 CALL smumps_bdc_error( myid, slavef, comm, keep )
1056 ELSE
1057 CALL mpi_recv( bufr, lbufr_bytes, mpi_packed,
1058 & msgsou, msgtag, comm, status, ierr )
1059 CALL smumps_traiter_message_solve( bufr, lbufr, lbufr_bytes,
1060 & msgtag, msgsou, myid, slavef, comm,
1061 & n, nrhs, ipool, lpool, leaf,
1062 & nbfin, nstk_s, iw, liw, a, la, ptrist, ptrfac,
1063 & iwcb, liwcb,
1064 & wcb, lwcb, poswcb,
1065 & pleftwcb, posiwcb,
1066 & ptricb, info, keep,keep8, dkeep, step,
1067 & procnode_steps,
1068 & rhscomp, lrhscomp, posinrhscomp_fwd
1069 & , from_pp
1070 & )
1071 END IF
1072 END IF
1073 RETURN
1074 END SUBROUTINE smumps_solve_recv_and_treat
1076 & NPIV, NCB, LIELL, CBINITZERO, LDEQLIELLPANEL,
1077 & RHSCOMP, LRHSCOMP, NRHS_B,
1078 & POSINRHSCOMP_FWD, N,
1079 & WCB,
1080 & IW, LIW, J1, J3, J2, KEEP, DKEEP)
1081 IMPLICIT NONE
1082 INTEGER, INTENT( IN ) :: NPIV, NCB, LIELL, N,
1083 & lrhscomp, nrhs_b,
1084 & liw, j1, j2, j3
1085 LOGICAL, INTENT( IN ) :: LDEQLIELLPANEL
1086 LOGICAL, INTENT( IN ) :: CBINITZERO
1087 INTEGER, INTENT( IN ) :: POSINRHSCOMP_FWD( N ), IW( LIW )
1088 REAL, INTENT( INOUT ) :: RHSCOMP( LRHSCOMP, NRHS_B )
1089 REAL, INTENT( OUT ) :: WCB( int(LIELL,8)*
1090 & int(nrhs_b,8) )
1091 INTEGER :: KEEP(500)
1092 REAL :: DKEEP(150)
1093 INTEGER, PARAMETER :: ZERO = 0.0e0
1094 INTEGER(8), PARAMETER :: PPIV_COURANT = 1_8
1095 INTEGER(8) :: PCB_COURANT
1096 INTEGER :: LD_WCBCB, LD_WCBPIV, J, JJ, K, IPOSINRHSCOMP
1097 INTEGER(8) :: IFR8, IFR_ini8
1098 INCLUDE 'mpif.h'
1099 LOGICAL :: OMP_FLAG
1100 IF ( LDEQLIELLPANEL ) THEN
1101 LD_WCBPIV = liell
1102 ld_wcbcb = liell
1103 pcb_courant = ppiv_courant + npiv
1104 ELSE
1105 ld_wcbpiv = npiv
1106 ld_wcbcb = ncb
1107 pcb_courant = ppiv_courant + npiv * nrhs_b
1108 ENDIF
1109 IF ( ldeqliellpanel ) THEN
1110 DO k=1, nrhs_b
1111 ifr8 = ppiv_courant+int(k-1,8)*int(ld_wcbpiv,8)-1_8
1112 iposinrhscomp = posinrhscomp_fwd(iw(j1))
1113 DO jj = j1, j3
1114 ifr8 = ifr8 + 1_8
1115 wcb(ifr8) = rhscomp(iposinrhscomp,k)
1116 iposinrhscomp = iposinrhscomp + 1
1117 ENDDO
1118 IF (ncb.GT.0 .AND. .NOT. cbinitzero) THEN
1119#if defined(__ve__)
1120!NEC$ IVDEP
1121#endif
1122 DO jj = j3+1, j2
1123 j = iw(jj)
1124 ifr8 = ifr8 + 1_8
1125 iposinrhscomp = abs(posinrhscomp_fwd(j))
1126 wcb(ifr8) = rhscomp(iposinrhscomp,k)
1127 rhscomp(iposinrhscomp,k) = zero
1128 ENDDO
1129 ENDIF
1130 ENDDO
1131 ELSE
1132 pcb_courant = ppiv_courant + ld_wcbpiv*nrhs_b
1133 ifr8 = ppiv_courant - 1_8
1134 ifr_ini8 = ifr8
1135 iposinrhscomp = posinrhscomp_fwd(iw(j1))
1136 omp_flag = .false.
1137!$ OMP_FLAG = ( NRHS_B .GE. KEEP(362) .AND.
1138!$ & int(NCB,8)*int(NRHS_B,8) .GE. KEEP(363) )
1139 IF (omp_flag) THEN
1140!$OMP PARALLEL DO PRIVATE(JJ,IFR8)
1141 DO k=1, nrhs_b
1142 ifr8 = ifr_ini8 + int(k-1,8)*int(npiv,8)
1143 DO jj = j1, j3
1144 wcb(ifr8+int(jj-j1+1,8)) =
1145 & rhscomp(iposinrhscomp+jj-j1,k)
1146 ENDDO
1147 ENDDO
1148!$OMP END PARALLEL DO
1149 ELSE
1150 DO k=1, nrhs_b
1151 ifr8 = ifr_ini8 + int(k-1,8)*int(npiv,8)
1152 DO jj = j1, j3
1153 wcb(ifr8+int(jj-j1+1,8)) =
1154 & rhscomp(iposinrhscomp+jj-j1,k)
1155 ENDDO
1156 ENDDO
1157 ENDIF
1158 ifr8 = pcb_courant - 1_8
1159 IF (ncb.GT.0 .AND. .NOT. cbinitzero) THEN
1160 ifr_ini8 = ifr8
1161 omp_flag = .false.
1162!$ OMP_FLAG = ( NRHS_B.GE.KEEP(362) .AND.
1163!$ & NCB*NRHS_B .GE. KEEP(363) )
1164 IF (omp_flag) THEN
1165!$OMP PARALLEL DO PRIVATE (IFR8, JJ, J, IPOSINRHSCOMP)
1166 DO k=1, nrhs_b
1167 ifr8 = ifr_ini8+(k-1)*ncb
1168#if defined(__ve__)
1169!NEC$ IVDEP
1170#endif
1171 DO jj = j3 + 1, j2
1172 j = iw(jj)
1173 iposinrhscomp = abs(posinrhscomp_fwd(j))
1174 wcb(ifr8+int(jj-j3,8)) = rhscomp(iposinrhscomp,k)
1175 rhscomp(iposinrhscomp,k)=zero
1176 ENDDO
1177 ENDDO
1178!$omp END PARALLEL do
1179 ELSE
1180 DO k=1, nrhs_b
1181 ifr8 = ifr_ini8+(k-1)*ncb
1182#if defined(__ve__)
1183!NEC$ IVDEP
1184#endif
1185 DO jj = j3 + 1, j2
1186 j = iw(jj)
1187 iposinrhscomp = abs(posinrhscomp_fwd(j))
1188 wcb(ifr8+int(jj-j3,8)) = rhscomp(iposinrhscomp,k)
1189 rhscomp(iposinrhscomp,k)=zero
1190 ENDDO
1191 ENDDO
1192 ENDIF
1193 ENDIF
1194 ENDIF
1195 IF ( cbinitzero ) THEN
1196 omp_flag = .false.
1197!$ OMP_FLAG = int(NCB,8)*int(NRHS_B,8) .GE. KEEP(363)
1198 IF (omp_flag) THEN
1199!$omp parallel DO collapse(2)
1200 DO k = 1, nrhs_b
1201 DO jj = 1, ncb
1202 wcb(pcb_courant+int(k-1,8)*int(ld_wcbcb,8)+jj-1_8) = zero
1203 ENDDO
1204 ENDDO
1205!$OMP END PARALLEL DO
1206 ELSE
1207 DO k = 1, nrhs_b
1208 DO jj = 1, ncb
1209 wcb(pcb_courant+int(k-1,8)*int(ld_wcbcb,8)+jj-1_8) = zero
1210 ENDDO
1211 ENDDO
1212 ENDIF
1213 ENDIF
1214 RETURN
1215 END SUBROUTINE smumps_rhscomp_to_wcb
#define mumps_abort
Definition VE_Metis.h:25
#define alpha
Definition eval.h:35
subroutine sgemv(trans, m, n, alpha, a, lda, x, incx, beta, y, incy)
SGEMV
Definition sgemv.f:156
subroutine strsv(uplo, trans, diag, n, a, lda, x, incx)
STRSV
Definition strsv.f:149
subroutine strsm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
STRSM
Definition strsm.f:181
subroutine sgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
SGEMM
Definition sgemm.f:187
#define min(a, b)
Definition macros.h:20
#define max(a, b)
Definition macros.h:21
subroutine mpi_recv(buf, cnt, datatype, source, tag, comm, status, ierr)
Definition mpi.f:461
subroutine mpi_iprobe(source, tag, comm, flag, status, ierr)
Definition mpi.f:360
subroutine mpi_get_count(status, datatype, cnt, ierr)
Definition mpi.f:296
subroutine mpi_probe(source, tag, comm, status, ierr)
Definition mpi.f:449
subroutine mpi_unpack(inbuf, insize, position, outbuf, outcnt, datatype, comm, ierr)
Definition mpi.f:514
subroutine mumps_bloc2_get_slave_info(keep, keep8, inode, step, n, slavef, istep_to_iniv2, tab_pos_in_pere islave, ncb, nslaves, size, first_index)
subroutine, public smumps_buf_send_master2slave(nrhs, inode, ifath, eff_cb_size, ld_cb, ld_piv, npiv, jbdeb, jbfin, cb, sol, dest, comm, keep, ierr)
subroutine, public smumps_buf_send_vcb(nrhs_b, node1, node2, ncb, ldw, long, iw, w, jbdeb, jbfin, rhscomp, nrhs, lrhscomp, iposinrhscomp, npiv, keep, dest, tag, comm, ierr)
integer function, public smumps_ooc_panel_size(nnmax)
subroutine smumps_free_factors_for_solve(inode, ptrfac, nsteps, a, la, flag, ierr)
integer, parameter nb_lock_max
Definition ssol_omp_m.F:16
subroutine smumps_sol_fwd_lr_su(inode, n, iwhdlr, npiv_global, nslaves, iw, ipos_init, liw, liell, wcb, lwcb, ld_wcbpiv, ld_wcbcb, ppiv_init, pcb_init, rhscomp, lrhscomp, nrhs, posinrhscomp_fwd, jbdeb, jbfin, mtype, keep, keep8, oocwrite_compatible_with_blr, iflag, ierror)
Definition ssol_lr.F:31
subroutine smumps_sol_slave_lr_u(inode, iwhdlr, npiv_global, wcb, lwcb, ldx, ldy, ptrx_init, ptry_init, jbdeb, jbfin, mtype, keep, keep8, iflag, ierror)
Definition ssol_lr.F:189
subroutine smumps_bdc_error(myid, slavef, comm, keep)
Definition sbcast_int.F:38
subroutine smumps_ooc_pp_check_perm_freed(iw_location, must_be_permuted)
subroutine smumps_get_ooc_perm_ptr(typef, nbpanels, i_pivptr, i_piv, ipos, iw, liw)
subroutine smumps_permute_panel(ipiv, lpiv, ishift, the_panel, nbrow, nbcol, kbeforepanel)
subroutine smumps_sol_ld_and_reload_panel(inode, n, npiv, liell, nelim, nslaves, ppiv_courant, iw, ipos, liw, a, la, apos, wcb, lwcb, ld_wcbpiv, rhscomp, lrhscomp, nrhs, posinrhscomp_fwd, jbdeb, jbfin, mtype, keep, oocwrite_compatible_with_blr, ignore_k459)
Definition ssol_aux.F:1379
subroutine smumps_solve_fwd_trsolve(a, la, apos, npiv, ldadiag, nrhs_b, wcb, lwcb, lda_wcb, ppiv_courant, mtype, keep)
Definition ssol_aux.F:1147
subroutine smumps_sol_ld_and_reload(inode, n, npiv, liell, nelim, nslaves, ppiv_courant, iw, ipos, liw, a, la, apos, wcb, lwcb, ld_wcbpiv, rhscomp, lrhscomp, nrhs, posinrhscomp_fwd, jbdeb, jbfin, mtype, keep, oocwrite_compatible_with_blr, ignore_k459)
Definition ssol_aux.F:1509
subroutine smumps_solve_gemm_update(a, la, apos1, nx, lda, ny, nrhs_b, wcb, lwcb, ptrx, ldx, ptry, ldy, mtype, keep, coef_y)
Definition ssol_aux.F:1326
subroutine smumps_solve_get_ooc_node(inode, ptrfac, keep, a, la, step, keep8, n, must_be_permuted, ierr)
Definition ssol_aux.F:732
subroutine smumps_solve_fwd_panels(a, la, apos, npiv, iw, nrhs_b, wcb, lwcb, lda_wcb, ppiv_courant, mtype, keep)
Definition ssol_aux.F:1225
subroutine smumps_rhscomp_to_wcb(npiv, ncb, liell, cbinitzero, ldeqliellpanel, rhscomp, lrhscomp, nrhs_b, posinrhscomp_fwd, n, wcb, iw, liw, j1, j3, j2, keep, dkeep)
recursive subroutine smumps_solve_recv_and_treat(bloq, flag, bufr, lbufr, lbufr_bytes, myid, slavef, comm, n, nrhs, ipool, lpool, leaf, nbfin, nstk_s, iw, liw, a, la, ptrist, ptrfac, iwcb, liwcb, wcb, lwcb, poswcb, pleftwcb, posiwcb, ptricb, info, keep, keep8, dkeep, step, procnode_steps, rhscomp, lrhscomp, posinrhscomp_fwd, from_pp)
recursive subroutine smumps_traiter_message_solve(bufr, lbufr, lbufr_bytes, msgtag, msgsou, myid, slavef, comm, n, nrhs, ipool, lpool, leaf, nbfin, nstk_s, iw, liw, a, la, ptrist, ptrfac, iwcb, liwcb, wcb, lwcb, poswcb, pleftwcb, posiwcb, ptricb, info, keep, keep8, dkeep, step, procnode_steps, rhscomp, lrhscomp, posinrhscomp_fwd, from_pp)
subroutine smumps_solve_node_fwd(inode, lastfsl0sta, lastfsl0dyn, bufr, lbufr, lbufr_bytes, myid, slavef, comm, n, ipool, lpool, leaf, nbfin, nstk_s, iwcb, liwcb, wcb, lwcb, a, la, iw, liw, nrhs, poswcb, pleftwcb, posiwcb, ptricb, ptrist, ptrfac, procnode_steps, fils, step, frere, dad, info, keep, keep8, dkeep, rhs_root, lrhs_root, mtype, rhscomp, lrhscomp, posinrhscomp_fwd istep_to_iniv2, tab_pos_in_pere, rhs_bounds, lrhs_bounds, do_nbsparse, from_pp, error_was_broadcasted)
integer function mumps_procnode(procinfo_inode, k199)
subroutine mumps_set_ierror(size8, ierror)
subroutine mumps_geti8(i8, int_array)