OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
tools_common.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 mumps_find_unit(IUNIT)
15 IMPLICIT NONE
16 INTEGER :: IUNIT
17 INTEGER, PARAMETER :: UNIT_MIN = 10
18 INTEGER, PARAMETER :: UNIT_MAX = 500
19 INTEGER :: I
20 LOGICAL :: BUSY
21 iunit = -1
22 DO i = unit_min, unit_max
23 INQUIRE(unit=i, opened=busy)
24 IF ( .NOT. busy ) THEN
25 iunit = i
26 RETURN
27 END IF
28 ENDDO
29 RETURN
30 END SUBROUTINE mumps_find_unit
31 SUBROUTINE mumps_make1root( N, FRERE, FILS, NFSIZ, THEROOT )
32 IMPLICIT NONE
33 INTEGER, intent( in ) :: N
34 INTEGER, intent( in ) :: NFSIZ( N )
35 INTEGER, intent( inout ) :: FRERE( N ), FILS( N )
36 INTEGER, intent( out ) :: THEROOT
37 INTEGER INODE, IROOT, IFILS, IN, IROOTLAST, SIZE
38 iroot = -9999
39 SIZE = 0
40 DO inode = 1, n
41 IF ( frere( inode ) .EQ. 0 ) THEN
42 IF ( nfsiz( inode ) .GT. SIZE ) THEN
43 SIZE = nfsiz( inode )
44 iroot = inode
45 END IF
46 ENDIF
47 END DO
48 in = iroot
49 DO WHILE ( fils( in ) .GT. 0 )
50 in = fils( in )
51 END DO
52 irootlast = in
53 ifils = - fils( in )
54 DO inode = 1, n
55 IF ( frere( inode ) .eq. 0 .and. inode .ne. iroot ) THEN
56 IF ( ifils .eq. 0 ) THEN
57 fils( irootlast ) = - inode
58 frere( inode ) = -iroot
59 ifils = inode
60 ELSE
61 frere( inode ) = -fils( irootlast )
62 fils( irootlast ) = - inode
63 END IF
64 END IF
65 END DO
66 theroot = iroot
67 RETURN
68 END SUBROUTINE mumps_make1root
69 INTEGER FUNCTION mumps_encode_tpn_iproc(TPN,IPROC,K199)
70 INTEGER, INTENT(IN) :: TPN, iproc, k199
71 IF (k199 < 0) THEN
72 mumps_encode_tpn_iproc = iproc + ishft(tpn+1, 24)
73 ELSE
74 mumps_encode_tpn_iproc = (tpn-1)*k199+iproc+1
75 ENDIF
76 RETURN
77 END FUNCTION mumps_encode_tpn_iproc
78 INTEGER FUNCTION mumps_typenode_rough(PROCINFO_INODE, K199)
79 IMPLICIT NONE
80 INTEGER k199
81 INTEGER procinfo_inode
82 IF (k199 < 0) THEN
83 mumps_typenode_rough = ishft(procinfo_inode,-24) - 1
84 ELSE
85 mumps_typenode_rough = (procinfo_inode-1+2*k199)/k199 - 1
86 ENDIF
87 RETURN
88 END FUNCTION mumps_typenode_rough
89 INTEGER FUNCTION mumps_typenode(PROCINFO_INODE, K199)
90 IMPLICIT NONE
91 INTEGER k199
92 INTEGER procinfo_inode, tpn
93 IF (k199 < 0) THEN
94 tpn = ishft(procinfo_inode,-24) - 1
95 IF (tpn < 1 ) THEN
96 tpn = 1
97 ELSE IF (tpn.GE.4) THEN
98 tpn = 2
99 ENDIF
100 ELSE
101 IF (procinfo_inode <= k199 ) THEN
102 tpn = 1
103 ELSE
104 tpn = (procinfo_inode-1+2*k199)/k199 - 1
105 IF ( tpn .LT. 1 ) tpn = 1
106 IF (tpn.EQ.4.OR.tpn.EQ.5.OR.tpn.EQ.6) tpn = 2
107 END IF
108 END IF
109 mumps_typenode = tpn
110 RETURN
111 END FUNCTION mumps_typenode
112 SUBROUTINE mumps_typeandprocnode( TPN,
113 & MUMPS_PROCNODE, PROCINFO_INODE, K199 )
114 INTEGER, INTENT(IN) :: K199, PROCINFO_INODE
115 INTEGER, intent(out) :: TPN, MUMPS_PROCNODE
116 IF (k199 < 0 ) THEN
117 mumps_procnode=iand(procinfo_inode,
118#if defined(MUMPS_F2003)
119 & int(b"111111111111111111111111"))
120#else
121 & 16777215)
122#endif
123 tpn = ishft(procinfo_inode,-24) - 1
124 IF (tpn < 1 ) THEN
125 tpn = 1
126 ELSE IF (tpn.GE.4) THEN
127 tpn = 2
128 ENDIF
129 ELSE
130 IF (k199 == 1) THEN
131 mumps_procnode = 0
132 IF (procinfo_inode <= k199) THEN
133 tpn = 1
134 ELSE
135 tpn = 3
136 ENDIF
137 ELSE
138 tpn = (procinfo_inode-1+2*k199)/k199-1
139 mumps_procnode = (procinfo_inode-1+2*k199)-
140 & (tpn+1)*k199
141 IF (tpn .LT. 1) THEN
142 tpn = 1
143 ELSE IF (tpn .ge. 4) THEN
144 tpn = 2
145 ENDIF
146 ENDIF
147 ENDIF
148 RETURN
149 END SUBROUTINE mumps_typeandprocnode
150 INTEGER FUNCTION mumps_procnode(PROCINFO_INODE, K199)
151 IMPLICIT NONE
152 INTEGER k199
153 INTEGER procinfo_inode
154 IF ( k199 < 0 ) THEN
155 mumps_procnode=iand(procinfo_inode,
156#if defined(MUMPS_F2003)
157 & int(b"111111111111111111111111"))
158#else
159 & 16777215 )
160#endif
161 ELSE
162 IF (k199 == 1) THEN
164 ELSE
165 mumps_procnode=mod(2*k199+procinfo_inode-1,k199)
166 END IF
167 ENDIF
168 RETURN
169 END FUNCTION mumps_procnode
170 INTEGER FUNCTION mumps_typesplit (PROCINFO_INODE, K199)
171 IMPLICIT NONE
172 INTEGER, intent(in) :: k199
173 INTEGER procinfo_inode, tpn
174 IF (k199 < 0) THEN
175 tpn = ishft(procinfo_inode,-24) - 1
176 IF (tpn < 1 ) tpn = 1
177 ELSE
178 IF (procinfo_inode <= k199 ) THEN
179 tpn = 1
180 ELSE
181 tpn = (procinfo_inode-1+2*k199)/k199 - 1
182 IF ( tpn .LT. 1 ) tpn = 1
183 ENDIF
184 ENDIF
185 mumps_typesplit = tpn
186 RETURN
187 END FUNCTION mumps_typesplit
188 LOGICAL FUNCTION mumps_rootssarbr( PROCINFO_INODE, K199 )
189 IMPLICIT NONE
190 INTEGER k199
191 INTEGER tpn, procinfo_inode
192 IF (k199 < 0) THEN
193 tpn = ishft(procinfo_inode,-24) - 1
194 ELSE
195 tpn = (procinfo_inode-1+2*k199)/k199 - 1
196 ENDIF
197 mumps_rootssarbr = ( tpn .eq. 0 )
198 RETURN
199 END FUNCTION mumps_rootssarbr
200 LOGICAL FUNCTION mumps_inssarbr( PROCINFO_INODE, K199 )
201 IMPLICIT NONE
202 INTEGER k199
203 INTEGER tpn, procinfo_inode
204 IF (k199 < 0) THEN
205 tpn = ishft(procinfo_inode,-24) - 1
206 ELSE
207 tpn = (procinfo_inode-1+k199+k199)/k199 - 1
208 ENDIF
209 mumps_inssarbr = ( tpn .eq. -1 )
210 RETURN
211 END FUNCTION mumps_inssarbr
212 LOGICAL FUNCTION mumps_in_or_root_ssarbr
213 & ( procinfo_inode, k199 )
214 IMPLICIT NONE
215 INTEGER k199
216 INTEGER tpn, procinfo_inode
217 IF (k199 < 0) THEN
218 tpn = ishft(procinfo_inode,-24) - 1
219 ELSE
220 tpn = (procinfo_inode-1+k199+k199)/k199 - 1
221 ENDIF
223 & ( tpn .eq. -1 .OR. tpn .eq. 0 )
224 RETURN
225 END FUNCTION mumps_in_or_root_ssarbr
227 & SSARBR, INODE, DAD, N,
228 & KEEP28,
229 & STEP, PROCNODE_STEPS, K199)
230 IMPLICIT NONE
231 INTEGER, INTENT(IN) :: N, KEEP28, K199, INODE
232 INTEGER, INTENT(IN) :: DAD(KEEP28), PROCNODE_STEPS(KEEP28)
233 INTEGER, INTENT(IN) :: STEP(N)
234 LOGICAL, INTENT(OUT) :: SSARBR
235 INTEGER :: DADINODE, TYPEDAD
236 LOGICAL, EXTERNAL :: MUMPS_INSSARBR
237 INTEGER, EXTERNAL :: MUMPS_TYPENODE
238 ssarbr = .false.
239 dadinode = dad(step(inode))
240 IF (dadinode .NE. 0) THEN
241 typedad = mumps_typenode(procnode_steps(step(dadinode)),
242 & k199)
243 IF (typedad.EQ.1) THEN
244 ssarbr=mumps_inssarbr(procnode_steps(step(dadinode)),
245 & k199)
246 ENDIF
247 ENDIF
248 RETURN
249 END SUBROUTINE mumps_set_ssarbr_dad
250 LOGICAL FUNCTION mumps_i_am_candidate( MYID, SLAVEF, INODE,
251 & NMB_PAR2, ISTEP_TO_INIV2 , K71, STEP, N,
252 & CANDIDATES, KEEP24 )
253 IMPLICIT NONE
254 INTEGER myid, slavef, inode, nmb_par2, keep24, i
255 INTEGER k71, n
256 INTEGER istep_to_iniv2 ( k71 ), step ( n )
257 INTEGER candidates(slavef+1, max(nmb_par2,1))
258 INTEGER ncand, posinode
259 mumps_i_am_candidate = .false.
260 IF (keep24 .eq. 0) RETURN
261 posinode = istep_to_iniv2( step(inode) )
262 ncand = candidates( slavef+1, posinode )
263 DO i = 1, ncand
264 IF (myid .EQ. candidates( i, posinode ))
265 & mumps_i_am_candidate = .true.
266 END DO
267 RETURN
268 END FUNCTION mumps_i_am_candidate
269 SUBROUTINE mumps_secdeb(T)
270 DOUBLE PRECISION T
271 DOUBLE PRECISION MPI_WTIME
272 EXTERNAL mpi_wtime
273 t=mpi_wtime()
274 RETURN
275 END SUBROUTINE mumps_secdeb
276 SUBROUTINE mumps_secfin(T)
277 DOUBLE PRECISION T
278 DOUBLE PRECISION MPI_WTIME
279 EXTERNAL mpi_wtime
280 t=mpi_wtime()-t
281 RETURN
282 END SUBROUTINE mumps_secfin
283 SUBROUTINE mumps_sort_doubles( N, VAL, ID )
284 INTEGER N
285 INTEGER ID( N )
286 DOUBLE PRECISION VAL( N )
287 INTEGER I, ISWAP
288 DOUBLE PRECISION SWAP
289 LOGICAL DONE
290 done = .false.
291 DO WHILE ( .NOT. done )
292 done = .true.
293 DO i = 1, n - 1
294 IF ( val( i ) .GT. val( i + 1 ) ) THEN
295 done = .false.
296 iswap = id( i )
297 id( i ) = id( i + 1 )
298 id( i + 1 ) = iswap
299 swap = val( i )
300 val( i ) = val( i + 1 )
301 val( i + 1 ) = swap
302 END IF
303 END DO
304 END DO
305 RETURN
306 END SUBROUTINE mumps_sort_doubles
307 SUBROUTINE mumps_sort_doubles_dec( N, VAL, ID )
308 INTEGER N
309 INTEGER ID( N )
310 DOUBLE PRECISION VAL( N )
311 INTEGER I, ISWAP
312 DOUBLE PRECISION SWAP
313 LOGICAL DONE
314 done = .false.
315 DO WHILE ( .NOT. done )
316 done = .true.
317 DO i = 1, n - 1
318 IF ( val( i ) .LT. val( i + 1 ) ) THEN
319 done = .false.
320 iswap = id( i )
321 id( i ) = id( i + 1 )
322 id( i + 1 ) = iswap
323 swap = val( i )
324 val( i ) = val( i + 1 )
325 val( i + 1 ) = swap
326 END IF
327 END DO
328 END DO
329 RETURN
330 END SUBROUTINE mumps_sort_doubles_dec
331#if defined (PESSL)
332 SUBROUTINE descinit( DESC, M, N, MB, NB, IRSRC, ICSRC, ICTXT,
333 & LLD, INFO )
334 INTEGER ICSRC, ICTXT, INFO, IRSRC, LLD, M, MB, N, NB
335 INTEGER DESC( * )
336 INTEGER BLOCK_CYCLIC_2D, CSRC_, CTXT_, DLEN_, DTYPE_,
337 & lld_, mb_, m_, nb_, n_, rsrc_
338# if defined(DESC8)
339 parameter( dlen_ = 8, dtype_ = 1,
340 & ctxt_ = 7, m_ = 1, n_ = 2, mb_ = 3, nb_ = 4,
341 & rsrc_ = 5, csrc_ = 6, lld_ = 8 )
342# else
343 parameter( block_cyclic_2d = 1, dlen_ = 9, dtype_ = 1,
344 & ctxt_ = 2, m_ = 3, n_ = 4, mb_ = 5, nb_ = 6,
345 & rsrc_ = 7, csrc_ = 8, lld_ = 9 )
346# endif
347 INTEGER MYCOL, MYROW, NPCOL, NPROW
348 EXTERNAL blacs_gridinfo, PXERBLA
349 INTEGER NUMROC
350 EXTERNAL numroc
351 INTRINSIC max, min
352 CALL blacs_gridinfo( ictxt, nprow, npcol, myrow, mycol )
353 info = 0
354 IF( m.LT.0 ) THEN
355 info = -2
356 ELSE IF( n.LT.0 ) THEN
357 info = -3
358 ELSE IF( mb.LT.1 ) THEN
359 info = -4
360 ELSE IF( nb.LT.1 ) THEN
361 info = -5
362 ELSE IF( irsrc.LT.0 .OR. irsrc.GE.nprow ) THEN
363 info = -6
364 ELSE IF( icsrc.LT.0 .OR. icsrc.GE.npcol ) THEN
365 info = -7
366 ELSE IF( nprow.EQ.-1 ) THEN
367 info = -8
368 ELSE IF( lld.LT.max( 1, numroc( m, mb, myrow, irsrc,
369 & nprow ) ) ) THEN
370 info = -9
371 END IF
372 IF( info.NE.0 )
373 & CALL pxerbla( ictxt, 'DESCINIT', -info )
374# ifndef DESC8
375 desc( dtype_ ) = block_cyclic_2d
376# endif
377 desc( m_ ) = max( 0, m )
378 desc( n_ ) = max( 0, n )
379 desc( mb_ ) = max( 1, mb )
380 desc( nb_ ) = max( 1, nb )
381 desc( rsrc_ ) = max( 0, min( irsrc, nprow-1 ) )
382 desc( csrc_ ) = max( 0, min( icsrc, npcol-1 ) )
383 desc( ctxt_ ) = ictxt
384 desc( lld_ ) = max( lld, max( 1, numroc( desc( m_ ), desc( mb_ ),
385 & myrow, desc( rsrc_ ), nprow ) ) )
386 RETURN
387 END SUBROUTINE descinit
388 SUBROUTINE pxerbla( ICTXT, SRNAME, INFO )
389 INTEGER ICTXT, INFO
390 CHARACTER*(*) SRNAME
391 INTEGER MYCOL, MYROW, NPCOL, NPROW
392 EXTERNAL blacs_gridinfo
393 CALL blacs_gridinfo( ictxt, nprow, npcol, myrow, mycol )
394 WRITE( *, fmt = 9999 ) myrow, mycol, srname, info
395 9999 FORMAT( '{', i5, ',', i5, '}: On entry to ', a,
396 & ' parameter number', i4, ' had an illegal value' )
397 END SUBROUTINE pxerbla
398#endif
399 SUBROUTINE mumps_mem_centralize(MYID, COMM, INFO, INFOG, IRANK)
400 IMPLICIT NONE
401 INTEGER MYID, COMM, IRANK, INFO, INFOG(2)
402 include 'mpif.h'
403 INTEGER IERR_MPI, MASTER
404#if defined(WORKAROUNDINTELILP64MPI2INTEGER)
405 INTEGER(4) :: TEMP1(2),TEMP2(2)
406#else
407 INTEGER :: TEMP1(2),TEMP2(2)
408#endif
409 parameter( master = 0 )
410 CALL mpi_reduce( info, infog(1), 1, mpi_integer,
411 & mpi_max, master, comm, ierr_mpi )
412 CALL mpi_reduce( info, infog(2), 1, mpi_integer,
413 & mpi_sum, master, comm, ierr_mpi )
414 temp1(1) = info
415 temp1(2) = myid
416 CALL mpi_reduce( temp1, temp2, 1, mpi_2integer,
417 & mpi_maxloc, master, comm, ierr_mpi )
418 IF ( myid.eq. master ) THEN
419 IF ( infog(1) .ne. temp2(1) ) THEN
420 write(*,*) 'Error in MUMPS_MEM_CENTRALIZE'
421 CALL mumps_abort()
422 END IF
423 irank = temp2(2)
424 ELSE
425 irank = -1
426 END IF
427 RETURN
428 END SUBROUTINE mumps_mem_centralize
429 INTEGER FUNCTION mumps_get_pool_length
430 & (max_active_nodes,keep,keep8)
431 IMPLICIT NONE
432 INTEGER max_active_nodes
433 INTEGER keep(500)
434 INTEGER(8) keep8(150)
435 mumps_get_pool_length = max_active_nodes + 1 + 3
436 RETURN
437 END FUNCTION mumps_get_pool_length
439 & nb_prun_roots, Pruned_Roots,
440 & MYROOT, MYID_NODES,
441 & KEEP, KEEP8, STEP, PROCNODE_STEPS,
442 & IPOOL, LPOOL )
443 IMPLICIT NONE
444 INTEGER, INTENT(IN) :: N, MYID_NODES, LPOOL, nb_prun_roots
445 INTEGER KEEP(500)
446 INTEGER(8) KEEP8(150)
447 INTEGER, INTENT(IN) :: STEP(N)
448 INTEGER, INTENT(IN) :: PROCNODE_STEPS(KEEP(28))
449 INTEGER, INTENT(IN) :: Pruned_Roots(nb_prun_roots)
450 INTEGER, INTENT(OUT) :: MYROOT
451 INTEGER, INTENT(OUT) :: IPOOL(LPOOL)
452 INTEGER, EXTERNAL :: MUMPS_PROCNODE
453 INTEGER :: I, INODE
454 myroot = 0
455 DO i = nb_prun_roots, 1, -1
456 inode = pruned_roots(i)
457 IF (mumps_procnode(procnode_steps(step(inode)),
458 & keep(199)) .EQ. myid_nodes) THEN
459 myroot = myroot + 1
460 ipool(myroot) = inode
461 ENDIF
462 END DO
463 RETURN
464 END SUBROUTINE mumps_init_pool_dist_bwd
466 & nb_prun_roots, Pruned_Roots,
467 & MYROOT, MYID_NODES,
468 & KEEP, KEEP8, STEP, PROCNODE_STEPS,
469 & IPOOL, LPOOL, TO_PROCESS )
470 IMPLICIT NONE
471 INTEGER, INTENT(IN) :: N, MYID_NODES, LPOOL, nb_prun_roots
472 INTEGER KEEP(500)
473 INTEGER(8) KEEP8(150)
474 INTEGER, INTENT(IN) :: STEP(N)
475 INTEGER, INTENT(IN) :: PROCNODE_STEPS(KEEP(28))
476 LOGICAL, INTENT(IN) :: TO_PROCESS(KEEP(28))
477 INTEGER, INTENT(IN) :: Pruned_Roots(nb_prun_roots)
478 INTEGER, INTENT(OUT) :: MYROOT
479 INTEGER, INTENT(OUT) :: IPOOL(LPOOL)
480 INTEGER, EXTERNAL :: MUMPS_PROCNODE
481 INTEGER :: I, INODE
482 MYROOT = 0
483 do i = nb_prun_roots, 1, -1
484 inode = pruned_roots(i)
485 IF (mumps_procnode(procnode_steps(step(inode)),
486 & keep(199)) .EQ. myid_nodes) THEN
487 IF ( to_process(step(inode)) ) THEN
488 myroot = myroot + 1
489 ipool(myroot) = inode
490 ENDIF
491 ENDIF
492 END DO
493 RETURN
494 END SUBROUTINE mumps_init_pool_dist_bwd_l0
495 SUBROUTINE mumps_init_pool_dist_na_bwd(N, MYROOT, MYID_NODES,
496 & NA, LNA, KEEP, KEEP8, STEP, PROCNODE_STEPS,
497 & IPOOL, LPOOL )
498 IMPLICIT NONE
499 INTEGER, INTENT(IN) :: N, MYID_NODES, LPOOL, LNA
500 INTEGER KEEP(500)
501 INTEGER(8) KEEP8(150)
502 INTEGER, INTENT(IN) :: STEP(N)
503 INTEGER, INTENT(IN) :: PROCNODE_STEPS(KEEP(28)), NA(LNA)
504 INTEGER, INTENT(OUT) :: IPOOL(LPOOL)
505 INTEGER, INTENT(OUT) :: MYROOT
506 INTEGER, EXTERNAL :: MUMPS_PROCNODE
507 INTEGER :: NBLEAF, NBROOT, I, INODE
508 nbleaf = na(1)
509 nbroot = na(2)
510 myroot = 0
511 DO i = nbroot, 1, -1
512 inode = na(nbleaf+i+2)
513 IF (mumps_procnode(procnode_steps(step(inode)),
514 & keep(199)) .EQ. myid_nodes) THEN
515 myroot = myroot + 1
516 ipool(myroot) = inode
517 ENDIF
518 END DO
519 RETURN
520 END SUBROUTINE mumps_init_pool_dist_na_bwd
521 SUBROUTINE mumps_init_pool_dist_na_bwd_l0(N, MYROOT, MYID_NODES,
522 & NA, LNA, KEEP, KEEP8, STEP, PROCNODE_STEPS,
523 & IPOOL, LPOOL, L0_OMP_MAPPING )
524 IMPLICIT NONE
525 INTEGER, INTENT(IN) :: N, MYID_NODES, LPOOL, LNA
526 INTEGER KEEP(500)
527 INTEGER(8) KEEP8(150)
528 INTEGER, INTENT(IN) :: STEP(N)
529 INTEGER, INTENT(IN) :: PROCNODE_STEPS(KEEP(28)), NA(LNA)
530 INTEGER, INTENT(IN) :: L0_OMP_MAPPING(KEEP(28))
531 INTEGER, INTENT(OUT) :: IPOOL(LPOOL)
532 INTEGER, INTENT(OUT) :: MYROOT
533 INTEGER, EXTERNAL :: MUMPS_PROCNODE
534 INTEGER :: NBLEAF, NBROOT, I, INODE
535 nbleaf = na(1)
536 nbroot = na(2)
537 myroot = 0
538 DO i = nbroot, 1, -1
539 inode = na(nbleaf+i+2)
540 IF (mumps_procnode(procnode_steps(step(inode)),
541 & keep(199)) .EQ. myid_nodes) THEN
542 IF ( l0_omp_mapping(step(inode)).EQ.0 ) THEN
543 myroot = myroot + 1
544 ipool(myroot) = inode
545 ENDIF
546 ENDIF
547 END DO
548 RETURN
549 END SUBROUTINE mumps_init_pool_dist_na_bwd_l0
551 & MYID_NODES,
552 & NA, LNA, KEEP, KEEP8, STEP, PROCNODE_STEPS,
553 & IPOOL, LPOOL, L0_OMP_MAPPING, TO_PROCESS )
554 IMPLICIT NONE
555 INTEGER, INTENT(IN) :: N, MYID_NODES, LPOOL, LNA
556 INTEGER KEEP(500)
557 INTEGER(8) KEEP8(150)
558 INTEGER, INTENT(IN) :: STEP(N)
559 INTEGER, INTENT(IN) :: PROCNODE_STEPS(KEEP(28)), NA(LNA)
560 INTEGER, INTENT(IN) :: L0_OMP_MAPPING(KEEP(28))
561 INTEGER, INTENT(OUT) :: IPOOL(LPOOL)
562 INTEGER, INTENT(OUT) :: MYROOT
563 LOGICAL, INTENT(IN) :: TO_PROCESS( KEEP(28) )
564 INTEGER, EXTERNAL :: MUMPS_PROCNODE
565 INTEGER :: NBLEAF, NBROOT, I, INODE
566 nbleaf = na(1)
567 nbroot = na(2)
568 myroot = 0
569 DO i = nbroot, 1, -1
570 inode = na(nbleaf+i+2)
571 IF (mumps_procnode(procnode_steps(step(inode)),
572 & keep(199)) .EQ. myid_nodes) THEN
573 IF ( l0_omp_mapping(step(inode)).EQ.0 ) THEN
574 IF ( to_process( step(inode) ) ) THEN
575 myroot = myroot + 1
576 ipool(myroot) = inode
577 ENDIF
578 ENDIF
579 ENDIF
580 END DO
581 RETURN
583 SUBROUTINE mumps_init_pool_dist(N, LEAF,
584 & MYID_NODES,
585 & K199, NA, LNA, KEEP,KEEP8, STEP,
586 & PROCNODE_STEPS, IPOOL, LPOOL)
587 IMPLICIT NONE
588 INTEGER N, LEAF, MYID_NODES,
589 & k199, lpool, lna
590 INTEGER KEEP(500)
591 INTEGER(8) KEEP8(150)
592 INTEGER STEP(N)
593 INTEGER PROCNODE_STEPS(KEEP(28)), NA(LNA),
594 & ipool(lpool)
595 INTEGER NBLEAF, INODE, I
596 INTEGER MUMPS_PROCNODE
597 EXTERNAL MUMPS_PROCNODE
598 nbleaf = na(1)
599 leaf = 1
600 DO i = 1, nbleaf
601 inode = na(i+2)
602 IF (mumps_procnode(procnode_steps(step(inode)),keep(199))
603 & .EQ.myid_nodes) THEN
604 ipool(leaf) = inode
605 leaf = leaf + 1
606 ENDIF
607 ENDDO
608 RETURN
609 END SUBROUTINE mumps_init_pool_dist
611 & (n, leaf, myid_nodes,
612 & lleaves, leaves, keep,keep8, step,
613 & procnode_steps, ipool, lpool)
614 IMPLICIT NONE
615 INTEGER N, LEAF, MYID_NODES,
616 & lpool, lleaves
617 INTEGER KEEP(500)
618 INTEGER(8) KEEP8(150)
619 INTEGER STEP(N)
620 INTEGER PROCNODE_STEPS(KEEP(28)), LEAVES(LLEAVES),
621 & ipool(lpool)
622 INTEGER I, INODE
623 INTEGER, EXTERNAL :: MUMPS_PROCNODE
624 leaf = 1
625 DO i = 1, lleaves
626 inode = leaves(i)
627 IF ( mumps_procnode(procnode_steps(step(inode)),keep(199))
628 & .EQ.myid_nodes ) THEN
629 ipool( leaf ) = inode
630 leaf = leaf + 1
631 ENDIF
632 ENDDO
633 RETURN
634 END SUBROUTINE mumps_init_pool_dist_nona
635 SUBROUTINE mumps_init_nroot_dist(N, NBROOT,
636 & NROOT_LOC, MYID_NODES,
637 & SLAVEF, NA, LNA, KEEP, STEP,
638 & PROCNODE_STEPS)
639 IMPLICIT NONE
640 INTEGER, INTENT( OUT ) :: NROOT_LOC
641 INTEGER, INTENT( OUT ) :: NBROOT
642 INTEGER, INTENT( IN ) :: KEEP( 500 )
643 INTEGER, INTENT( IN ) :: SLAVEF
644 INTEGER, INTENT( IN ) :: N
645 INTEGER, INTENT( IN ) :: STEP(N)
646 INTEGER, INTENT( IN ) :: LNA
647 INTEGER, INTENT( IN ) :: NA(LNA)
648 INTEGER, INTENT( IN ) :: PROCNODE_STEPS(KEEP(28))
649 INTEGER, INTENT( IN ) :: MYID_NODES
650 INTEGER, EXTERNAL :: MUMPS_PROCNODE
651 INTEGER :: INODE, I, NBLEAF
652 nbleaf = na(1)
653 nbroot = na(2)
654 nroot_loc = 0
655 DO i = 1, nbroot
656 inode = na(i+2+nbleaf)
657 IF (mumps_procnode(procnode_steps(step(inode)),
658 & keep(199)).EQ.myid_nodes) THEN
659 nroot_loc = nroot_loc + 1
660 END IF
661 ENDDO
662 RETURN
663 END SUBROUTINE mumps_init_nroot_dist
665 & (n, nbrorl, rorl_list,
666 & nrorl_loc, myid_nodes,
667 & slavef, keep, step,
668 & procnode_steps)
669 IMPLICIT NONE
670 INTEGER, INTENT( OUT ) :: NRORL_LOC
671 INTEGER, INTENT( IN ) :: NBRORL
672 INTEGER, INTENT( IN ) :: RORL_LIST(NBRORL)
673 INTEGER, INTENT( IN ) :: KEEP( 500 )
674 INTEGER, INTENT( IN ) :: SLAVEF
675 INTEGER, INTENT( IN ) :: N
676 INTEGER, INTENT( IN ) :: STEP(N)
677 INTEGER, INTENT( IN ) :: PROCNODE_STEPS(KEEP(28))
678 INTEGER, INTENT( IN ) :: MYID_NODES
679 INTEGER I, INODE
680 INTEGER, EXTERNAL :: MUMPS_PROCNODE
681 nrorl_loc = 0
682 DO i = 1, nbrorl
683 inode = rorl_list(i)
684 IF (mumps_procnode(procnode_steps(step(inode)),
685 & keep(199)).EQ.myid_nodes) THEN
686 nrorl_loc = nrorl_loc + 1
687 END IF
688 ENDDO
689 RETURN
690 END SUBROUTINE mumps_nblocal_roots_or_leaves
691 LOGICAL FUNCTION mumps_compare_tab(TAB1,TAB2,LEN1,LEN2)
692 IMPLICIT NONE
693 INTEGER len1 , len2 ,i
694 INTEGER tab1(len1)
695 INTEGER tab2(len2)
696 mumps_compare_tab=.false.
697 IF(len1 .NE. len2) THEN
698 RETURN
699 ENDIF
700 DO i=1 , len1
701 IF(tab1(i) .NE. tab2(i)) THEN
702 RETURN
703 ENDIF
704 ENDDO
705 mumps_compare_tab=.true.
706 RETURN
707 END FUNCTION mumps_compare_tab
708 SUBROUTINE mumps_sort_int( N, VAL, ID )
709 INTEGER N
710 INTEGER ID( N )
711 INTEGER VAL( N )
712 INTEGER I, ISWAP
713 INTEGER SWAP
714 LOGICAL DONE
715 done = .false.
716 DO WHILE ( .NOT. done )
717 done = .true.
718 DO i = 1, n - 1
719 IF ( val( i ) .GT. val( i + 1 ) ) THEN
720 done = .false.
721 iswap = id( i )
722 id( i ) = id( i + 1 )
723 id( i + 1 ) = iswap
724 swap = val( i )
725 val( i ) = val( i + 1 )
726 val( i + 1 ) = swap
727 END IF
728 END DO
729 END DO
730 RETURN
731 END SUBROUTINE mumps_sort_int
732 SUBROUTINE mumps_sort_int_dec( N, VAL, ID )
733 INTEGER N
734 INTEGER ID( N )
735 INTEGER VAL( N )
736 INTEGER I, ISWAP
737 INTEGER SWAP
738 LOGICAL DONE
739 done = .false.
740 DO WHILE ( .NOT. done )
741 done = .true.
742 DO i = 1, n - 1
743 IF ( val( i ) .LT. val( i + 1 ) ) THEN
744 done = .false.
745 iswap = id( i )
746 id( i ) = id( i + 1 )
747 id( i + 1 ) = iswap
748 swap = val( i )
749 val( i ) = val( i + 1 )
750 val( i + 1 ) = swap
751 END IF
752 END DO
753 END DO
754 RETURN
755 END SUBROUTINE mumps_sort_int_dec
756 SUBROUTINE mumps_sort_int8( N, VAL, ID )
757 INTEGER N
758 INTEGER ID( N )
759 INTEGER(8) :: VAL( N )
760 INTEGER I, ISWAP
761 INTEGER(8) SWAP
762 LOGICAL DONE
763 done = .false.
764 DO WHILE ( .NOT. done )
765 done = .true.
766 DO i = 1, n - 1
767 IF ( val( i ) .GT. val( i + 1 ) ) THEN
768 done = .false.
769 iswap = id( i )
770 id( i ) = id( i + 1 )
771 id( i + 1 ) = iswap
772 swap = val( i )
773 val( i ) = val( i + 1 )
774 val( i + 1 ) = swap
775 END IF
776 END DO
777 END DO
778 RETURN
779 END SUBROUTINE mumps_sort_int8
780 SUBROUTINE mumps_abort()
781#if defined(PRINT_BACKTRACE_ON_ABORT)
782#if defined(__INTEL_COMPILER)
783 USE ifcore
784#endif
785 IMPLICIT NONE
786 include 'mpif.h'
787 INTEGER IERR, IERRCODE
788#if defined(__GFORTRAN__)
789 CALL backtrace()
790#endif
791#if defined(__INTEL_COMPILER)
792!$OMP CRITICAL(MUMPS_TRACEBACKQQ)
793 CALL tracebackqq("MUMPS_ABORT calls TRACEBACKQQ:",
794 & user_exit_code=-1)
795!$OMP END CRITICAL(MUMPS_TRACEBACKQQ)
796#endif
797#else
798 IMPLICIT NONE
799 include 'mpif.h'
800 INTEGER IERR, IERRCODE
801#endif
802 ierrcode = -99
803 CALL mpi_abort(mpi_comm_world, ierrcode, ierr)
804 RETURN
805 END SUBROUTINE mumps_abort
806 SUBROUTINE mumps_get_perlu(KEEP12,ICNTL14,
807 & KEEP50,KEEP54,ICNTL6,ICNTL8)
808 IMPLICIT NONE
809 INTEGER, intent(out)::KEEP12
810 INTEGER, intent(in)::ICNTL14,KEEP50,KEEP54,ICNTL6,ICNTL8
811 KEEP12 = icntl14
812 RETURN
813 END SUBROUTINE mumps_get_perlu
814 SUBROUTINE mumps_reducei8( IN, OUT, MPI_OP, ROOT, COMM)
815 IMPLICIT NONE
816 include 'mpif.h'
817 INTEGER ROOT, COMM, MPI_OP
818 INTEGER(8) IN, OUT
819 INTEGER IERR
820 DOUBLE PRECISION DIN, DOUT
821 DIN =dble(in)
822 dout=0.0d0
823 CALL mpi_reduce(din, dout, 1, mpi_double_precision,
824 & mpi_op, root, comm, ierr)
825 out=int(dout,kind=8)
826 RETURN
827 END SUBROUTINE mumps_reducei8
828 SUBROUTINE mumps_allreducei8( IN, OUT, MPI_OP, COMM)
829 IMPLICIT NONE
830 include 'mpif.h'
831 INTEGER COMM, MPI_OP
832 INTEGER(8) IN, OUT
833 INTEGER IERR
834 DOUBLE PRECISION DIN, DOUT
835 din =dble(in)
836 dout=0.0d0
837 CALL mpi_allreduce(din, dout, 1, mpi_double_precision,
838 & mpi_op, comm, ierr)
839 out=int(dout,kind=8)
840 RETURN
841 END SUBROUTINE mumps_allreducei8
842 SUBROUTINE mumps_seti8toi4(I8, I)
843 IMPLICIT NONE
844 INTEGER , INTENT(OUT) :: I
845 INTEGER(8), INTENT(IN) :: I8
846 IF ( i8 .GT. int(huge(i),8) ) THEN
847 i = -int(i8/1000000_8,kind(i))
848 ELSE
849 i = int(i8,kind(i))
850 ENDIF
851 RETURN
852 END SUBROUTINE mumps_seti8toi4
853 SUBROUTINE mumps_abort_on_overflow(I8, STRING)
854 IMPLICIT NONE
855 INTEGER(8), INTENT(IN) :: I8
856 CHARACTER(*), INTENT(IN) :: STRING
857 INTEGER I
858 IF ( i8 .GT. int(huge(i),8)) THEN
859 WRITE(*,*) string
860 CALL mumps_abort()
861 ENDIF
862 RETURN
863 END SUBROUTINE mumps_abort_on_overflow
864 SUBROUTINE mumps_set_ierror( SIZE8, IERROR )
865 INTEGER(8), INTENT(IN) :: SIZE8
866 INTEGER, INTENT(OUT) :: IERROR
867 CALL MUMPS_SETI8TOI4(SIZE8, IERROR)
868 RETURN
869 END SUBROUTINE mumps_set_ierror
870 SUBROUTINE mumps_storei8(I8, INT_ARRAY)
871 IMPLICIT NONE
872 INTEGER(8), intent(in) :: I8
873 INTEGER, intent(out) :: INT_ARRAY(2)
874 INTEGER(kind(0_4)) :: I32
875 INTEGER(8) :: IDIV, IPAR
876 parameter(ipar=int(huge(i32),8))
877 parameter(idiv=ipar+1_8)
878 IF ( i8 .LT. idiv ) THEN
879 int_array(1) = 0
880 int_array(2) = int(i8)
881 ELSE
882 int_array(1) = int(i8 / idiv)
883 int_array(2) = int(mod(i8,idiv))
884 ENDIF
885 RETURN
886 END SUBROUTINE mumps_storei8
887 SUBROUTINE mumps_geti8(I8, INT_ARRAY)
888 IMPLICIT NONE
889 INTEGER(8), intent(out) :: I8
890 INTEGER, intent(in) :: INT_ARRAY(2)
891 INTEGER(kind(0_4)) :: I32
892 INTEGER(8) :: IDIV, IPAR
893 parameter(ipar=int(huge(i32),8))
894 parameter(idiv=ipar+1_8)
895 IF ( int_array(1) .EQ. 0 ) THEN
896 i8=int(int_array(2),8)
897 ELSE
898 i8=int(int_array(1),8)*idiv+int(int_array(2),8)
899 ENDIF
900 RETURN
901 END SUBROUTINE mumps_geti8
902 SUBROUTINE mumps_addi8toarray( INT_ARRAY, I8 )
903 IMPLICIT NONE
904 INTEGER(8), intent(in) :: I8
905 INTEGER, intent(inout) :: INT_ARRAY(2)
906 INTEGER(8) :: I8TMP
907 CALL mumps_geti8(i8tmp, int_array)
908 i8tmp = i8tmp + i8
909 CALL mumps_storei8(i8tmp, int_array)
910 RETURN
911 END SUBROUTINE mumps_addi8toarray
912 SUBROUTINE mumps_subtri8toarray( INT_ARRAY, I8 )
913 IMPLICIT NONE
914 INTEGER(8), intent(in) :: I8
915 INTEGER, intent(inout) :: INT_ARRAY(2)
916 INTEGER(8) :: I8TMP
917 CALL mumps_geti8(i8tmp, int_array)
918 i8tmp = i8tmp - i8
919 CALL mumps_storei8(i8tmp, int_array)
920 RETURN
921 END SUBROUTINE mumps_subtri8toarray
922 FUNCTION mumps_seqana_avail(ICNTL7)
923 LOGICAL :: mumps_seqana_avail
924 INTEGER, INTENT(IN) :: icntl7
925 LOGICAL :: scotch=.false.
926 LOGICAL :: metis =.false.
927#if defined(metis) || defined(parmetis) || defined(metis4) || defined(parmetis3)
928 metis = .true.
929#endif
930#if defined(scotch) || defined(ptscotch)
931 scotch = .true.
932#endif
933 IF ( icntl7 .LT. 0 .OR. icntl7 .GT. 7 ) THEN
934 mumps_seqana_avail = .false.
935 ELSE
936 mumps_seqana_avail = .true.
937 ENDIF
938 IF ( icntl7 .EQ. 5 ) mumps_seqana_avail = metis
939 IF ( icntl7 .EQ. 3 ) mumps_seqana_avail = scotch
940 RETURN
941 END FUNCTION mumps_seqana_avail
942 FUNCTION mumps_parana_avail(WHICH)
943 LOGICAL :: mumps_parana_avail
944 CHARACTER :: which*(*)
945 LOGICAL :: ptscotch=.false., parmetis=.false.
946#if defined(ptscotch)
947 ptscotch = .true.
948#endif
949#if defined(parmetis) || defined(parmetis3)
950 parmetis = .true.
951#endif
952 SELECT CASE(which)
953 CASE('ptscotch','PTSCOTCH')
954 mumps_parana_avail = ptscotch
955 CASE('parmetis','PARMETIS')
956 mumps_parana_avail = parmetis
957 CASE('both','BOTH')
958 mumps_parana_avail = ptscotch .AND. parmetis
959 CASE('any','ANY')
960 mumps_parana_avail = ptscotch .OR. parmetis
961 CASE default
962 write(*,'("Invalid input in MUMPS_PARANA_AVAIL")')
963 END SELECT
964 RETURN
965 END FUNCTION mumps_parana_avail
966 SUBROUTINE mumps_sort_step(N,FRERE,STEP,FILS,
967 & NA,LNA,NE,ND,DAD,LDAD,USE_DAD,
968 & NSTEPS,INFO,LP,
969 & PROCNODE,SLAVEF
970 & )
971 IMPLICIT NONE
972 INTEGER N, NSTEPS, LNA, LP,LDAD
973 INTEGER FRERE(NSTEPS), FILS(N), STEP(N)
974 INTEGER NA(LNA), NE(NSTEPS), ND(NSTEPS)
975 INTEGER DAD(LDAD)
976 LOGICAL USE_DAD
977 INTEGER INFO(80)
978 INTEGER SLAVEF,PROCNODE(NSTEPS)
979 INTEGER POSTORDER,TMP_SWAP
980 INTEGER, DIMENSION (:), ALLOCATABLE :: STEP_TO_NODE
981 INTEGER, DIMENSION (:), ALLOCATABLE :: IPOOL,TNSTK
982 INTEGER I,II,allocok
983 INTEGER NBLEAF,NBROOT,LEAF,IN,INODE,IFATH
984 EXTERNAL MUMPS_TYPENODE
985 INTEGER MUMPS_TYPENODE
986 postorder=1
987 nbleaf = na(1)
988 nbroot = na(2)
989 ALLOCATE( ipool(nbleaf), tnstk(nsteps), stat=allocok )
990 IF (allocok > 0) THEN
991 IF ( lp .GT. 0 )
992 & WRITE(lp,*)'Memory allocation error in MUMPS_SORT_STEP'
993 info(1)=-7
994 info(2)=nsteps
995 RETURN
996 ENDIF
997 DO i=1,nsteps
998 tnstk(i) = ne(i)
999 ENDDO
1000 ALLOCATE(step_to_node(nsteps),stat=allocok)
1001 IF (allocok > 0) THEN
1002 IF ( lp .GT. 0 )
1003 & WRITE(lp,*)'Memory allocation error in
1004 &MUMPS_SORT_STEP'
1005 info(1)=-7
1006 info(2)=nsteps
1007 RETURN
1008 ENDIF
1009 DO i=1,n
1010 IF(step(i).GT.0)THEN
1011 step_to_node(step(i))=i
1012 ENDIF
1013 ENDDO
1014 ipool(1:nbleaf)=na(3:2+nbleaf)
1015 leaf = nbleaf + 1
1016 91 CONTINUE
1017 IF (leaf.NE.1) THEN
1018 leaf = leaf -1
1019 inode = ipool(leaf)
1020 ENDIF
1021 96 CONTINUE
1022 IF (use_dad) THEN
1023 ifath = dad( step(inode) )
1024 ELSE
1025 in = inode
1026 113 in = frere(in)
1027 IF (in.GT.0) GO TO 113
1028 ifath = -in
1029 ENDIF
1030 tmp_swap=frere(step(inode))
1031 frere(step(inode))=frere(postorder)
1032 frere(postorder)=tmp_swap
1033 tmp_swap=nd(step(inode))
1034 nd(step(inode))=nd(postorder)
1035 nd(postorder)=tmp_swap
1036 tmp_swap=ne(step(inode))
1037 ne(step(inode))=ne(postorder)
1038 ne(postorder)=tmp_swap
1039 tmp_swap=procnode(step(inode))
1040 procnode(step(inode))=procnode(postorder)
1041 procnode(postorder)=tmp_swap
1042 IF(use_dad)THEN
1043 tmp_swap=dad(step(inode))
1044 dad(step(inode))=dad(postorder)
1045 dad(postorder)=tmp_swap
1046 ENDIF
1047 tmp_swap=tnstk(step(inode))
1048 tnstk(step(inode))=tnstk(postorder)
1049 tnstk(postorder)=tmp_swap
1050 ii=step_to_node(postorder)
1051 tmp_swap=step(inode)
1052 step(step_to_node(postorder))=tmp_swap
1053 step(inode)=postorder
1054 step_to_node(postorder)=inode
1055 step_to_node(tmp_swap)=ii
1056 in=ii
1057 101 in = fils(in)
1058 IF (in .GT. 0 ) THEN
1059 step(in)=-step(ii)
1060 GOTO 101
1061 ENDIF
1062 in=inode
1063 102 in = fils(in)
1064 IF (in .GT. 0 ) THEN
1065 step(in)=-step(inode)
1066 GOTO 102
1067 ENDIF
1068 postorder = postorder + 1
1069 IF (ifath.EQ.0) THEN
1070 nbroot = nbroot - 1
1071 IF (nbroot.EQ.0) GOTO 116
1072 GOTO 91
1073 ENDIF
1074 tnstk(step(ifath)) = tnstk(step(ifath)) - 1
1075 IF ( tnstk(step(ifath)) .EQ. 0 ) THEN
1076 inode = ifath
1077 GOTO 96
1078 ELSE
1079 GOTO 91
1080 ENDIF
1081 116 CONTINUE
1082 DEALLOCATE(step_to_node)
1083 DEALLOCATE(ipool,tnstk)
1084 RETURN
1085 END SUBROUTINE mumps_sort_step
1086 SUBROUTINE mumps_check_comm_nodes(COMM_NODES, EXIT_FLAG)
1087 IMPLICIT NONE
1088 INTEGER, INTENT(IN) :: COMM_NODES
1089 LOGICAL, INTENT(OUT) :: EXIT_FLAG
1090 include 'mumps_tags.h'
1091 include 'mpif.h'
1092 INTEGER :: STATUS(MPI_STATUS_SIZE), IERR
1093 CALL mpi_iprobe( mpi_any_source, terreur, comm_nodes,
1094 & exit_flag, status, ierr)
1095 RETURN
1096 END SUBROUTINE mumps_check_comm_nodes
1097 SUBROUTINE mumps_get_proc_per_node(K414, MyID, NbProcs, COMM)
1098 IMPLICIT NONE
1099 include 'mpif.h'
1100 INTEGER :: K414, MyID, NbProcs, COMM, ALLOCOK
1101 INTEGER :: ierr,MyNAME_length,MyNAME_length_RCV,i,j
1102 CHARACTER(len=MPI_MAX_PROCESSOR_NAME) :: MyNAME
1103 CHARACTER, dimension(:), allocatable :: MyNAME_TAB,MyNAME_TAB_RCV
1104 logical :: SAME_NAME
1105 call mpi_get_processor_name(myname, myname_length, ierr)
1106 allocate(myname_tab(myname_length), stat=allocok)
1107 IF(allocok.LT.0) THEN
1108 write(*,*) "Allocation error in MUMPS_GET_PROC_PER_NODE"
1109 call mumps_abort()
1110 ENDIF
1111 DO i=1, myname_length
1112 myname_tab(i) = myname(i:i)
1113 ENDDO
1114 k414=0
1115 do i=0, nbprocs-1
1116 if(myid .eq. i) then
1117 myname_length_rcv = myname_length
1118 else
1119 myname_length_rcv = 0
1120 endif
1121 call mpi_bcast(myname_length_rcv,1,mpi_integer,
1122 & i,comm,ierr)
1123 allocate(myname_tab_rcv(myname_length_rcv), stat=allocok)
1124 IF(allocok.LT.0) THEN
1125 write(*,*) "Allocation error in MUMPS_GET_PROC_PER_NODE"
1126 call mumps_abort()
1127 ENDIF
1128 if(myid .eq. i) then
1129 myname_tab_rcv = myname_tab
1130 endif
1131 call mpi_bcast(myname_tab_rcv,myname_length_rcv,mpi_character,
1132 & i,comm,ierr)
1133 same_name=.false.
1134 IF(myname_length .EQ. myname_length_rcv) THEN
1135 DO j=1, myname_length
1136 IF(myname_tab(j) .NE. myname_tab_rcv(j)) THEN
1137 goto 100
1138 ENDIF
1139 ENDDO
1140 same_name=.true.
1141 ENDIF
1142 100 continue
1143 IF(same_name) k414=k414+1
1144 deallocate(myname_tab_rcv)
1145 enddo
1146 deallocate(myname_tab)
1147 END SUBROUTINE mumps_get_proc_per_node
1148 SUBROUTINE mumps_icopy_32to64 (INTAB, SIZETAB, OUTTAB8)
1149 INTEGER, intent(in) :: SIZETAB
1150 INTEGER, intent(in) :: INTAB(SIZETAB)
1151 INTEGER(8), intent(out) :: OUTTAB8(SIZETAB)
1152 INTEGER :: I
1153 DO i=1,sizetab
1154 outtab8(i) = int(intab(i),8)
1155 ENDDO
1156 RETURN
1157 END SUBROUTINE mumps_icopy_32to64
1158 SUBROUTINE mumps_icopy_32to64_64c(INTAB, SIZETAB8, OUTTAB8)
1159 INTEGER(8), intent(in) :: SIZETAB8
1160 INTEGER, intent(in) :: INTAB(SIZETAB8)
1161 INTEGER(8), intent(out) :: OUTTAB8(SIZETAB8)
1162 INTEGER(8) :: I8
1163 LOGICAL :: OMP_FLAG
1164 omp_flag = (sizetab8 .GE.500000_8 )
1165!$OMP PARALLEL DO PRIVATE(I8)
1166!$OMP& IF(OMP_FLAG)
1167 DO i8=1_8, sizetab8
1168 outtab8(i8) = int(intab(i8),8)
1169 ENDDO
1170!$OMP END PARALLEL DO
1171 RETURN
1172 END SUBROUTINE mumps_icopy_32to64_64c
1173 SUBROUTINE mumps_icopy_32to64_64c_ip(IN_OUT_TAB48, SIZETAB)
1174 INTEGER(8), intent(in) :: SIZETAB
1175 INTEGER, intent(inout) :: IN_OUT_TAB48(2*SIZETAB)
1176 CALL mumps_icopy_32to64_64c_ip_rec(in_out_tab48, sizetab)
1177 RETURN
1178 END SUBROUTINE mumps_icopy_32to64_64c_ip
1179 RECURSIVE SUBROUTINE mumps_icopy_32to64_64c_ip_rec(
1180 & IN_OUT_TAB48, SIZETAB)
1181 IMPLICIT NONE
1182 INTEGER(8), intent(in) :: sizetab
1183 INTEGER :: in_out_tab48(2*sizetab)
1184 INTEGER(8) :: ibeg24, ibeg28, size1, size2
1185 IF (sizetab.LE. 1000_8) THEN
1186 CALL mumps_icopy_32to64_64c_ip_c(in_out_tab48(1),
1187 & sizetab)
1188 ELSE
1189 size2 = sizetab / 2
1190 size1 = sizetab - size2
1191 ibeg24 = size1+1
1192 ibeg28 = 2*size1+1_8
1193 CALL mumps_icopy_32to64_64c(in_out_tab48(ibeg24),
1194 & size2, in_out_tab48(ibeg28))
1195 CALL mumps_icopy_32to64_64c_ip_rec(in_out_tab48,
1196 & size1)
1197 ENDIF
1198 RETURN
1199 END SUBROUTINE mumps_icopy_32to64_64c_ip_rec
1200 SUBROUTINE mumps_icopy_64to32(INTAB8, SIZETAB, OUTTAB)
1201 INTEGER, intent(in) :: SIZETAB
1202 INTEGER(8), intent(in) :: INTAB8(SIZETAB)
1203 INTEGER, intent(out) :: OUTTAB(SIZETAB)
1204 INTEGER :: I
1205 DO i=1,sizetab
1206 outtab(i) = int(intab8(i))
1207 ENDDO
1208 RETURN
1209 END SUBROUTINE mumps_icopy_64to32
1210 SUBROUTINE mumps_icopy_64to32_64c (INTAB8, SIZETAB, OUTTAB)
1211 INTEGER(8), intent(in) :: SIZETAB
1212 INTEGER(8), intent(in) :: INTAB8(SIZETAB)
1213 INTEGER, intent(out) :: OUTTAB(SIZETAB)
1214 INTEGER(8) :: I8
1215 DO i8=1_8,sizetab
1216 outtab(i8) = int(intab8(i8))
1217 ENDDO
1218 RETURN
1219 END SUBROUTINE mumps_icopy_64to32_64c
1220 SUBROUTINE mumps_icopy_64to32_64c_ip(IN_OUT_TAB48, SIZETAB)
1221 INTEGER(8), intent(in) :: SIZETAB
1222 INTEGER, intent(inout) :: IN_OUT_TAB48(2*SIZETAB)
1223 CALL mumps_icopy_64to32_64c_ip_rec(in_out_tab48, sizetab)
1224 RETURN
1225 END SUBROUTINE mumps_icopy_64to32_64c_ip
1226 RECURSIVE SUBROUTINE mumps_icopy_64to32_64c_ip_rec(
1227 & IN_OUT_TAB48, SIZETAB)
1228 IMPLICIT NONE
1229 INTEGER(8), intent(in) :: sizetab
1230 INTEGER :: in_out_tab48(2*sizetab)
1231 INTEGER(8) :: ibeg24, ibeg28, size1, size2
1232 IF (sizetab.LE. 1000_8) THEN
1233 CALL mumps_icopy_64to32_64c_ip_c(in_out_tab48(1),
1234 & sizetab)
1235 ELSE
1236 size2 = sizetab / 2
1237 size1 = sizetab - size2
1238 ibeg24 = size1 + 1
1239 ibeg28 = size1 + size1 + 1_8
1240 CALL mumps_icopy_64to32_64c_ip_rec(in_out_tab48,
1241 & size1)
1242 CALL mumps_icopy_64to32_64c(in_out_tab48(ibeg28),
1243 & size2, in_out_tab48(ibeg24))
1244 ENDIF
1245 RETURN
1246 END SUBROUTINE mumps_icopy_64to32_64c_ip_rec
1247 SUBROUTINE mumps_get_nnz_internal( NNZ, NZ, NNZ_i )
1248 INTEGER , INTENT(IN) :: NZ
1249 INTEGER(8), INTENT(IN) :: NNZ
1250 INTEGER(8), INTENT(OUT) :: NNZ_i
1251 IF (nnz > 0_8) THEN
1252 nnz_i = nnz
1253 ELSE
1254 nnz_i = int(nz, 8)
1255 ENDIF
1256 END SUBROUTINE mumps_get_nnz_internal
1258 & N, NSTEPS, STEP, FRERE, FILS,
1259 & NA, LNA, NE, MAXNPIVTREE )
1260 IMPLICIT NONE
1261 INTEGER, intent(in) :: N, NSTEPS, LNA
1262 INTEGER, intent(in) :: FRERE(NSTEPS), FILS(N), STEP(N)
1263 INTEGER, intent(in) :: NA(LNA), NE(NSTEPS)
1264 INTEGER, intent(out) :: MAXNPIVTREE
1265 INTEGER :: IFATH,INODE,ISON
1266 INTEGER :: NPIV,ILEAF,NBLEAF,NBROOT
1267 INTEGER, DIMENSION(:) , ALLOCATABLE :: MAXNPIV
1268 INTEGER :: I, allocok
1269 maxnpivtree = -9999
1270 ALLOCATE ( maxnpiv(nsteps), stat=allocok)
1271 IF (allocok .gt.0) THEN
1272 WRITE(*, *) 'Allocation error in MUMPS_NPIV_CRITICAL_PATH'
1273 & ,nsteps
1274 CALL mumps_abort()
1275 ENDIF
1276 nbleaf = na(1)
1277 nbroot = na(2)
1278 maxnpiv = 0
1279 nbleaf = na(1)
1280 nbroot = na(2)
1281 DO ileaf = 1, nbleaf
1282 inode = na(2+ileaf)
1283 95 CONTINUE
1284 npiv = 0
1285 ison = inode
1286 100 npiv = npiv + 1
1287 ison = fils(ison)
1288 IF (ison .GT. 0 ) GOTO 100
1289 ison = -ison
1290 maxnpiv( step(inode) ) = npiv
1291 DO i = 1, ne(step(inode))
1292 maxnpiv(step(inode)) = max( maxnpiv(step(inode)),
1293 & npiv + maxnpiv(step(ison)) )
1294 ison = frere(step(ison))
1295 ENDDO
1296 ifath = inode
1297 DO WHILE (ifath .GT. 0)
1298 ifath = frere(step(ifath))
1299 ENDDO
1300 ifath = -ifath
1301 IF (ifath.EQ.0) THEN
1302 maxnpivtree = max(maxnpivtree, maxnpiv(step(inode)))
1303 ELSE
1304 IF (frere(step(inode)) .LT. 0) THEN
1305 inode = ifath
1306 GOTO 95
1307 ENDIF
1308 ENDIF
1309 ENDDO
1310 DEALLOCATE( maxnpiv )
1311 RETURN
1312 END SUBROUTINE mumps_npiv_critical_path
1313 SUBROUTINE mumps_ldltpanel_nbtarget( NPIV, NB_TARGET, KEEP )
1314 IMPLICIT NONE
1315 INTEGER, INTENT(IN) :: NPIV
1316 INTEGER, INTENT(IN) :: KEEP(500)
1317 INTEGER, INTENT(OUT) :: NB_TARGET
1318 INTEGER :: NBPANELS, NBCOLMIN, NBPANELSMAX
1319 IF (npiv .EQ. 0) THEN
1320 nb_target = 0
1321 ELSE
1322 nbcolmin = keep(460)
1323 nbpanelsmax = keep(459)
1324 nbpanels = min( (npiv+nbcolmin-1) / nbcolmin, nbpanelsmax )
1325 nb_target = ( npiv+nbpanels-1 ) / nbpanels
1326 ENDIF
1327 RETURN
1328 END SUBROUTINE mumps_ldltpanel_nbtarget
1330 & ( npiv, keep, iw, nb_entries )
1331 IMPLICIT NONE
1332 INTEGER, INTENT(IN) :: NPIV
1333 INTEGER, INTENT(IN) :: KEEP(500), IW(*)
1334 INTEGER(8), INTENT(OUT) :: NB_ENTRIES
1335 INTEGER :: NB_TARGET, NBCOLS_PANEL, NBROWS_PANEL
1336 INTEGER :: ICOL_BEG, ICOL_END, NBPANELS
1337 CALL mumps_ldltpanel_nbtarget( npiv, nb_target, keep )
1338 nb_entries = 0_8
1339 nbrows_panel = npiv
1340 icol_beg = 1
1341 nbpanels = 0
1342 DO WHILE ( icol_beg .LE. npiv )
1343 nbpanels = nbpanels + 1
1344 icol_end = min(nb_target * nbpanels, npiv)
1345 IF (iw(1) .NE. 0) THEN
1346 IF ( iw( icol_end ) < 0 ) THEN
1347 icol_end = icol_end + 1
1348 ENDIF
1349 ENDIF
1350 nbcols_panel = icol_end - icol_beg + 1
1351 nb_entries = nb_entries + int(nbcols_panel,8) *
1352 & int(nbrows_panel,8)
1353 nbrows_panel = nbrows_panel - nbcols_panel
1354 icol_beg = icol_end + 1
1355 ENDDO
1356 RETURN
1357 END SUBROUTINE mumps_ldltpanel_storage
1358 SUBROUTINE mumps_ldltpanel_panelinfos( NPIV, KEEP, IW,
1359 & NB_TARGET, NBPANELS, PANEL_COL, PANEL_POS, PANEL_TABSIZE,
1360 & IGNORE_K459 )
1361 IMPLICIT NONE
1362 INTEGER, INTENT(IN) :: NPIV
1363 INTEGER, INTENT(IN) :: IW( NPIV )
1364 INTEGER, INTENT(IN) :: KEEP(500)
1365 INTEGER, INTENT(IN) :: PANEL_TABSIZE
1366 INTEGER, INTENT(OUT) :: NB_TARGET, NBPANELS
1367 INTEGER, INTENT(OUT) :: PANEL_COL( PANEL_TABSIZE )
1368 INTEGER(8), INTENT(OUT) :: PANEL_POS( PANEL_TABSIZE )
1369 LOGICAL, INTENT(IN) :: IGNORE_K459
1370 INTEGER :: IPANEL, ICOL_END, NBROWS_PANEL, NBCOLS_PANEL
1371 IF ( ignore_k459 ) THEN
1372 nb_target = npiv
1373 ELSE
1374 CALL mumps_ldltpanel_nbtarget( npiv, nb_target, keep )
1375 ENDIF
1376 panel_pos(1) = 1_8
1377 panel_col(1) = 1
1378 nbrows_panel = npiv
1379 nbpanels = 1
1380 IF ( keep(459) .GT. 1 .AND. keep(50) .NE. 0 .AND.
1381 & nb_target.NE.npiv ) THEN
1382 nbpanels = ( npiv + nb_target -1 ) / nb_target
1383 IF ( panel_tabsize .LT. nbpanels + 1 ) THEN
1384 WRITE(*,*) " Internal error in MUMPS_LDLTPANEL_PANELINFOS",
1385 & panel_tabsize, nbpanels
1386 CALL mumps_abort()
1387 ENDIF
1388 DO ipanel=1, nbpanels
1389 icol_end = min(ipanel*nb_target, npiv)
1390 IF ( iw(icol_end) .LT. 0 ) THEN
1391 icol_end = icol_end + 1
1392 ENDIF
1393 nbcols_panel = icol_end - panel_col(ipanel) + 1
1394 panel_pos(ipanel+1) = panel_pos(ipanel) +
1395 & int(nbrows_panel,8)*int(nbcols_panel,8)
1396 panel_col(ipanel+1) = panel_col(ipanel) + nbcols_panel
1397 nbrows_panel = nbrows_panel - nbcols_panel
1398 ENDDO
1399 ELSE
1400 panel_pos(2) = int(npiv,8)*int(npiv,8)+1_8
1401 panel_col(2) = npiv+1
1402 ENDIF
1403 END SUBROUTINE mumps_ldltpanel_panelinfos
1405 & ( npiv, keep, iw, panel_sizes, nbpanels )
1406 IMPLICIT NONE
1407 INTEGER, INTENT(IN) :: NPIV
1408 INTEGER, INTENT(IN) :: KEEP(500), IW(NPIV)
1409 INTEGER(8), INTENT(OUT) :: PANEL_SIZES( KEEP(459) )
1410 INTEGER, INTENT(OUT) :: NBPANELS
1411 INTEGER :: NB_TARGET
1412 INTEGER :: ICOL_BEG, ICOL_END
1413 nbpanels = 0
1414 CALL mumps_ldltpanel_nbtarget( npiv, nb_target, keep )
1415 icol_beg = 1
1416 nbpanels = 0
1417 DO WHILE ( icol_beg .LE. npiv )
1418 nbpanels = nbpanels + 1
1419 icol_end = min(nb_target * nbpanels, npiv)
1420 IF ( iw( icol_end ) < 0 ) THEN
1421 icol_end = icol_end + 1
1422 ENDIF
1423 panel_sizes(nbpanels) = icol_end-icol_beg+1
1424 icol_beg = icol_end + 1
1425 ENDDO
1426 panel_sizes(nbpanels+1:keep(459))=0
1427 RETURN
1428 END SUBROUTINE mumps_ldltpanel_sizes
1430 & ( comm, newcomm, newsize, newrank )
1431 IMPLICIT NONE
1432 include 'mpif.h'
1433 INTEGER, INTENT(IN) :: COMM
1434 INTEGER, INTENT(OUT) :: NEWCOMM, NEWSIZE, NEWRANK
1435 INTEGER :: SMALLEST_ID_ON_SAME_NODE, IPROC, MYID, IERR, NPROCS
1436 INTEGER :: TMPNAME_LENGTH, MYNAME_LENGTH
1437 CHARACTER(len=MPI_MAX_PROCESSOR_NAME) :: MYNAME, TMPNAME
1438 smallest_id_on_same_node = -1
1439 CALL mpi_comm_rank( comm, myid, ierr )
1440 CALL mpi_comm_size( comm, nprocs, ierr )
1441 CALL mpi_get_processor_name(myname, myname_length, ierr )
1442 DO iproc = 0, nprocs - 1
1443 IF (myid .EQ. iproc) THEN
1444 tmpname = myname
1445 tmpname_length = myname_length
1446 ENDIF
1447 CALL mpi_bcast( tmpname_length, 1, mpi_integer,
1448 & iproc, comm, ierr )
1449 CALL mpi_bcast( tmpname, tmpname_length, mpi_character,
1450 & iproc, comm, ierr)
1451 IF (smallest_id_on_same_node .LT. 0) THEN
1452 IF ( tmpname_length .EQ. myname_length ) THEN
1453 IF ( tmpname(1:tmpname_length) .EQ. myname(1:myname_length) )
1454 & THEN
1455 smallest_id_on_same_node = iproc
1456 ENDIF
1457 ENDIF
1458 ENDIF
1459 ENDDO
1460 CALL mpi_comm_split( comm, smallest_id_on_same_node, 0,
1461 & newcomm, ierr )
1462 CALL mpi_comm_rank( newcomm, newrank, ierr )
1463 CALL mpi_comm_size( newcomm, newsize, ierr )
1464 RETURN
1465 END SUBROUTINE mumps_build_arch_node_comm
1466 SUBROUTINE mumps_destroy_arch_node_comm( ARCH_NODE_COMM )
1467 IMPLICIT NONE
1468 INTEGER :: ARCH_NODE_COMM, IERR
1469 INCLUDE 'mpif.h'
1470 CALL mpi_comm_free( arch_node_comm, ierr )
1471 RETURN
1472 END SUBROUTINE mumps_destroy_arch_node_comm
1474 & ( mem_count_allocated, atomic_updates, keep8,
1475 & iflag, ierror, k69upd, k71upd )
1476 IMPLICIT NONE
1477 INTEGER(8), INTENT(IN) :: MEM_COUNT_ALLOCATED
1478 INTEGER(8), INTENT(INOUT) :: KEEP8(150)
1479 LOGICAL, INTENT(IN) :: ATOMIC_UPDATES
1480 INTEGER, INTENT(INOUT) :: IFLAG, IERROR
1481 LOGICAL, INTENT(IN) :: K69UPD
1482 LOGICAL, INTENT(IN) :: K71UPD
1483 INTEGER(8) :: KEEP8TMPCOPY
1484 IF (mem_count_allocated.GT.0) THEN
1485 IF (atomic_updates ) THEN
1486!$OMP ATOMIC CAPTURE
1487 keep8(73) = keep8(73) + mem_count_allocated
1488 keep8tmpcopy = keep8(73)
1489!$OMP END ATOMIC
1490!$OMP ATOMIC UPDATE
1491 keep8(74) = max(keep8(74), keep8tmpcopy)
1492!$OMP END ATOMIC
1493 ELSE
1494 keep8(73) = keep8(73) + mem_count_allocated
1495 keep8tmpcopy = keep8(73)
1496 keep8(74) = max(keep8(74), keep8(73))
1497 ENDIF
1498 IF ( keep8tmpcopy .GT. keep8(75) ) THEN
1499 iflag = -19
1500 CALL mumps_set_ierror(
1501 & (keep8tmpcopy-keep8(75)), ierror)
1502 ENDIF
1503 IF ( k69upd ) THEN
1504 IF ( atomic_updates ) THEN
1505!$OMP ATOMIC CAPTURE
1506 keep8(69) = keep8(69) + mem_count_allocated
1507 keep8tmpcopy = keep8(69)
1508!$OMP END ATOMIC
1509!$OMP ATOMIC UPDATE
1510 keep8(68) = max(keep8(68), keep8tmpcopy)
1511!$OMP END ATOMIC
1512 ELSE
1513 keep8(69) = keep8(69) + mem_count_allocated
1514 keep8(68) = max(keep8(69), keep8(68))
1515 ENDIF
1516 ENDIF
1517 IF ( k71upd ) THEN
1518 IF ( atomic_updates ) THEN
1519!$OMP ATOMIC CAPTURE
1520 keep8(71) = keep8(71) + mem_count_allocated
1521 keep8tmpcopy = keep8(71)
1522!$OMP END ATOMIC
1523!$OMP ATOMIC UPDATE
1524 keep8(70) = max(keep8(70), keep8tmpcopy)
1525!$OMP END ATOMIC
1526 ELSE
1527 keep8(71) = keep8(71) + mem_count_allocated
1528 keep8(70) = max(keep8(71), keep8(70))
1529 ENDIF
1530 ENDIF
1531 ELSE
1532 IF (atomic_updates) THEN
1533!$OMP ATOMIC UPDATE
1534 keep8(73) = keep8(73) + mem_count_allocated
1535!$OMP END ATOMIC
1536 IF ( k69upd ) THEN
1537!$OMP ATOMIC UPDATE
1538 keep8(69) = keep8(69) + mem_count_allocated
1539!$OMP END ATOMIC
1540 ENDIF
1541 IF ( k71upd ) THEN
1542!$OMP ATOMIC UPDATE
1543 keep8(71) = keep8(71) + mem_count_allocated
1544!$OMP END ATOMIC
1545 ENDIF
1546 ELSE
1547 keep8(73) = keep8(73) + mem_count_allocated
1548 IF ( k69upd ) THEN
1549 keep8(69) = keep8(69) + mem_count_allocated
1550 ENDIF
1551 IF ( k71upd ) THEN
1552 keep8(71) = keep8(71) + mem_count_allocated
1553 ENDIF
1554 ENDIF
1555 ENDIF
1556 RETURN
1557 END SUBROUTINE mumps_dm_fac_upd_dyn_memcnts
1559 & SLAVEF,
1560 & KEEP,KEEP8,
1561 & PROCS,
1562 & MEM_DISTRIB, NCB, NFRONT, NSLAVES_NODE,
1563 & TAB_POS, SLAVES_LIST, SIZE_SLAVES_LIST,MYID,INODE,
1564 & TAB_MAXS_ARG,SUP_PROC_ARG,MAX_SURF,NB_ROW_MAX
1565 & )
1566 IMPLICIT NONE
1567 INTEGER, intent(in) :: KEEP(500),SIZE_SLAVES_LIST
1568 INTEGER(8) KEEP8(150)
1569 INTEGER, intent(in) :: SLAVEF, NFRONT, NCB,MYID
1570 INTEGER, intent(in) :: PROCS(SLAVEF+1)
1571 INTEGER(8), intent(in) :: TAB_MAXS_ARG(0:SLAVEF-1)
1572 INTEGER, intent(in) :: SUP_PROC_ARG(2)
1573 INTEGER, intent(in) :: MEM_DISTRIB(0:SLAVEF-1),INODE
1574 INTEGER, intent(out):: SLAVES_LIST(SIZE_SLAVES_LIST)
1575 INTEGER, intent(out):: TAB_POS(SLAVEF+2)
1576 INTEGER, intent(out):: NSLAVES_NODE,NB_ROW_MAX
1577 INTEGER(8), intent(out):: MAX_SURF
1578 LOGICAL :: FORCE_LDLTRegular_NIV2
1579 INTEGER NSLAVES,ACC
1580 INTEGER i,J,NELIM,NB_SUP,K50,NB_ROWS(PROCS(SLAVEF+1))
1581 INTEGER TMP_NROW,X,K
1582 LOGICAL SUP,MEM_CSTR
1583 DOUBLE PRECISION MAX_LOAD,TOTAL_LOAD,VAR,TMP,A,B,C,DELTA,
1584 & load_corr
1585 INTEGER IDWLOAD(SLAVEF)
1586 INTEGER(8) MEM_CONSTRAINT(2)
1587 k50=keep(50)
1588 force_ldltregular_niv2 = .false.
1589 max_surf=0
1590 nb_row_max=0
1591 nelim=nfront-ncb
1592 nb_sup=0
1593 total_load=0.0d0
1594 sup=.false.
1595 IF(sup_proc_arg(1).NE.
1596 & 0)THEN
1597 mem_constraint(1)=tab_maxs_arg(procs(1))
1598 total_load=total_load+dble(sup_proc_arg(1))/100.0d0
1599 nb_sup=nb_sup+1
1600 ENDIF
1601 IF(sup_proc_arg(2).NE.
1602 & 0)THEN
1603 mem_constraint(2)=tab_maxs_arg(procs(procs(slavef+1)))
1604 total_load=total_load+dble(sup_proc_arg(2))/100.0d0
1605 nb_sup=nb_sup+1
1606 ENDIF
1607 total_load=total_load+(procs(slavef+1)-nb_sup)
1608 IF(k50.EQ.0)THEN
1609 max_load=dble( nelim ) * dble( ncb ) +
1610 * dble(ncb) * dble(nelim)*dble(2*nfront-nelim-1)
1611 ELSE
1612 max_load=dble(nelim) * dble( ncb ) *
1613 * dble(nfront+1)
1614 ENDIF
1615 tmp=min(max_load,max_load/total_load)
1616 j=1
1617 DO i=1,procs(slavef+1)
1618 IF((nb_sup.GT.0).AND.(i.EQ.1))THEN
1619 cycle
1620 ELSEIF((nb_sup.EQ.2).AND.(i.EQ.procs(slavef+1)))THEN
1621 cycle
1622 ENDIF
1623 idwload(j)=procs(i)
1624 j=j+1
1625 ENDDO
1626 DO i=1,nb_sup
1627 IF(i.EQ.1)THEN
1628 idwload(j)=procs(1)
1629 ELSE
1630 idwload(j)=procs(procs(slavef+1))
1631 ENDIF
1632 j=j+1
1633 ENDDO
1634 IF ((k50.EQ.0).OR.force_ldltregular_niv2) THEN
1635 acc=0
1636 j=procs(slavef+1)-nb_sup+1
1637 DO i=1,nb_sup
1638 var=dble(sup_proc_arg(i))/100.0d0
1639 tmp_nrow=int(dble(mem_constraint(i))/dble(nfront))
1640 nb_rows(j)=int(max((var*dble(tmp))/
1641 & (dble(nelim)*dble(2*nfront-nelim)),
1642 & dble(1)))
1643 IF(nb_rows(j).GT.tmp_nrow)THEN
1644 nb_rows(j)=tmp_nrow
1645 ENDIF
1646 IF(ncb-acc.LT.nb_rows(j)) THEN
1647 nb_rows(j)=ncb-acc
1648 acc=ncb
1649 EXIT
1650 ENDIF
1651 acc=acc+nb_rows(j)
1652 j=j+1
1653 ENDDO
1654 IF(acc.EQ.ncb)THEN
1655 GOTO 777
1656 ENDIF
1657 DO i=1,procs(slavef+1)-nb_sup
1658 var=1.0d0
1659 tmp_nrow=int((dble(tab_maxs_arg(idwload(i))))/dble(nfront))
1660 nb_rows(i)=int((dble(var)*dble(tmp))/
1661 & (dble(nelim)*dble(2*nfront-nelim)))
1662 IF(nb_rows(i).GT.tmp_nrow)THEN
1663 nb_rows(i)=tmp_nrow
1664 ENDIF
1665 IF(ncb-acc.LT.nb_rows(i)) THEN
1666 nb_rows(i)=ncb-acc
1667 acc=ncb
1668 EXIT
1669 ENDIF
1670 acc=acc+nb_rows(i)
1671 ENDDO
1672 IF(acc.NE.ncb)THEN
1673 IF(procs(slavef+1).EQ.nb_sup)THEN
1674 tmp_nrow=(ncb-acc)/procs(slavef+1)+1
1675 DO i=1,procs(slavef+1)
1676 nb_rows(i)=nb_rows(i)+tmp_nrow
1677 IF(acc+tmp_nrow.GT.ncb)THEN
1678 nb_rows(i)=nb_rows(i)-tmp_nrow+ncb-acc
1679 acc=ncb
1680 EXIT
1681 ENDIF
1682 acc=acc+tmp_nrow
1683 ENDDO
1684 ELSE
1685 tmp_nrow=(ncb-acc)/(procs(slavef+1)-nb_sup)+1
1686 DO i=1,procs(slavef+1)-nb_sup
1687 nb_rows(i)=nb_rows(i)+tmp_nrow
1688 acc=acc+tmp_nrow
1689 IF(acc.GT.ncb) THEN
1690 nb_rows(i)=nb_rows(i)-tmp_nrow+
1691 & (ncb-(acc-tmp_nrow))
1692 EXIT
1693 ENDIF
1694 ENDDO
1695 ENDIF
1696 ENDIF
1697 ELSE
1698 acc=0
1699 i=procs(slavef+1)-nb_sup+1
1700 x=ncb
1701 load_corr=0.0d0
1702 mem_cstr=.false.
1703 DO j=1,nb_sup
1704 var=dble(sup_proc_arg(j))/dble(100)
1705 a=1.0d0
1706 b=dble(x+nelim)
1707 c=-dble(max(mem_constraint(j),0_8))
1708 delta=((b*b)-(4*a*c))
1709 tmp_nrow=int((-b+sqrt(delta))/(2*a))
1710 a=dble(-nelim)
1711 b=dble(nelim)*(dble(-nelim)+dble(2*(x+nelim)+1))
1712 c=-(var*tmp)
1713 delta=(b*b-(4*a*c))
1714 nb_rows(i)=int((-b+sqrt(delta))/(2*a))
1715 IF(nb_rows(i).GT.tmp_nrow)THEN
1716 nb_rows(i)=tmp_nrow
1717 mem_cstr=.true.
1718 ENDIF
1719 IF(acc+nb_rows(i).GT.ncb)THEN
1720 nb_rows(i)=ncb-acc
1721 acc=ncb
1722 x=0
1723 EXIT
1724 ENDIF
1725 x=x-nb_rows(i)
1726 acc=acc+nb_rows(i)
1727 load_corr=load_corr+(dble(nelim) * dble(nb_rows(i)) *
1728 * dble(2*(x+nelim) - nelim - nb_rows(i) + 1))
1729 i=i+1
1730 ENDDO
1731 IF(acc.EQ.ncb)THEN
1732 GOTO 777
1733 ENDIF
1734 IF((procs(slavef+1).NE.nb_sup).AND.mem_cstr)THEN
1735 tmp=(max_load-load_corr)/(procs(slavef+1)-nb_sup)
1736 ENDIF
1737 x=acc
1738 acc=0
1739 DO i=1,procs(slavef+1)-nb_sup
1740 IF (keep(375) .EQ. 1) THEN
1741 var=1.0d0
1742 a=dble(nelim)
1743 b=dble(nelim)*(dble(nelim)+dble(2*acc+1))
1744 c=-(var*tmp)
1745 ELSE
1746 a=1.0d0
1747 b=dble(acc+nelim)
1748 c=-tmp
1749 ENDIF
1750 delta=((b*b)-(4*a*c))
1751 nb_rows(i)=int((-b+sqrt(delta))/(2*a))
1752 IF(ncb-acc-x.LT.nb_rows(i))THEN
1753 nb_rows(i)=ncb-acc-x
1754 acc=ncb-x
1755 EXIT
1756 ENDIF
1757 acc=acc+nb_rows(i)
1758 ENDDO
1759 acc=acc+x
1760 IF(acc.NE.ncb)THEN
1761 IF(procs(slavef+1).EQ.nb_sup)THEN
1762 tmp_nrow=(ncb-acc)/procs(slavef+1)+1
1763 DO i=1,procs(slavef+1)
1764 nb_rows(i)=nb_rows(i)+tmp_nrow
1765 IF(acc+tmp_nrow.GT.ncb)THEN
1766 nb_rows(i)=nb_rows(i)-tmp_nrow+ncb-acc
1767 acc=ncb
1768 EXIT
1769 ENDIF
1770 acc=acc+tmp_nrow
1771 ENDDO
1772 ELSE
1773 nb_rows(procs(slavef+1)-nb_sup)=
1774 & nb_rows(procs(slavef+1)
1775 & -nb_sup)+ncb-acc
1776 ENDIF
1777 ENDIF
1778 ENDIF
1779 777 CONTINUE
1780 nslaves=0
1781 acc=1
1782 j=1
1783 k=1
1784 DO i=1,procs(slavef+1)
1785 IF(nb_rows(i).NE.0)THEN
1786 slaves_list(j)=idwload(i)
1787 tab_pos(j)=acc
1788 acc=acc+nb_rows(i)
1789 nb_row_max=max(nb_row_max,nb_rows(i))
1790 IF(k50.EQ.0)THEN
1791 max_surf=max(int(nb_rows(i),8)*int(ncb,8),int(0,8))
1792 ELSE
1793 max_surf=max(int(nb_rows(i),8)*int(acc,8),int(0,8))
1794 ENDIF
1795 nslaves=nslaves+1
1796 j=j+1
1797 ELSE
1798 slaves_list(procs(slavef+1)-k+1)=idwload(i)
1799 k=k+1
1800 ENDIF
1801 ENDDO
1802 tab_pos(slavef+2) = nslaves
1803 tab_pos(nslaves+1)= ncb+1
1804 nslaves_node=nslaves
1805 END SUBROUTINE mumps_set_parti_regular
#define min(a, b)
Definition macros.h:20
#define max(a, b)
Definition macros.h:21
subroutine mpi_iprobe(source, tag, comm, flag, status, ierr)
Definition mpi.f:360
subroutine mpi_reduce(sendbuf, recvbuf, cnt, datatype, op, root, comm, ierr)
Definition mpi.f:120
subroutine mpi_comm_split(comm, color, key, comm2, ierr)
Definition mpi.f:272
subroutine mpi_allreduce(sendbuf, recvbuf, cnt, datatype, operation, comm, ierr)
Definition mpi.f:103
subroutine mpi_get_processor_name(name, resultlen, ierror)
Definition mpi.f:196
subroutine pxerbla(contxt, srname, info)
Definition mpi.f:1600
subroutine mpi_comm_size(comm, size, ierr)
Definition mpi.f:263
subroutine mpi_bcast(buffer, cnt, datatype, root, comm, ierr)
Definition mpi.f:205
subroutine mpi_comm_free(comm, ierr)
Definition mpi.f:238
subroutine descinit(desc, m, n, mb, nb, irsrc, icsrc, ictxt, lld, info)
Definition mpi.f:777
subroutine mpi_comm_rank(comm, rank, ierr)
Definition mpi.f:254
subroutine mpi_abort(comm, ierrcode, ierr)
Definition mpi.f:153
subroutine blacs_gridinfo(cntxt, nprow, npcol, myrow, mycol)
Definition mpi.f:754
recursive subroutine mumps_icopy_64to32_64c_ip_rec(in_out_tab48, sizetab)
subroutine mumps_nblocal_roots_or_leaves(n, nbrorl, rorl_list, nrorl_loc, myid_nodes, slavef, keep, step, procnode_steps)
subroutine mumps_addi8toarray(int_array, i8)
subroutine mumps_secfin(t)
subroutine mumps_init_pool_dist_na_bwd_l0(n, myroot, myid_nodes, na, lna, keep, keep8, step, procnode_steps, ipool, lpool, l0_omp_mapping)
logical function mumps_in_or_root_ssarbr(procinfo_inode, k199)
subroutine mumps_get_perlu(keep12, icntl14, keep50, keep54, icntl6, icntl8)
subroutine mumps_init_pool_dist(n, leaf, myid_nodes, k199, na, lna, keep, keep8, step, procnode_steps, ipool, lpool)
subroutine mumps_ldltpanel_panelinfos(npiv, keep, iw, nb_target, nbpanels, panel_col, panel_pos, panel_tabsize, ignore_k459)
subroutine mumps_icopy_32to64_64c(intab, sizetab8, outtab8)
subroutine mumps_sort_doubles_dec(n, val, id)
subroutine mumps_allreducei8(in, out, mpi_op, comm)
subroutine mumps_seti8toi4(i8, i)
logical function mumps_rootssarbr(procinfo_inode, k199)
subroutine mumps_check_comm_nodes(comm_nodes, exit_flag)
recursive subroutine mumps_icopy_32to64_64c_ip_rec(in_out_tab48, sizetab)
logical function mumps_parana_avail(which)
subroutine mumps_storei8(i8, int_array)
subroutine mumps_init_pool_dist_bwd_l0(n, nb_prun_roots, pruned_roots, myroot, myid_nodes, keep, keep8, step, procnode_steps, ipool, lpool, to_process)
subroutine mumps_sort_int8(n, val, id)
logical function mumps_inssarbr(procinfo_inode, k199)
subroutine mumps_sort_step(n, frere, step, fils, na, lna, ne, nd, dad, ldad, use_dad, nsteps, info, lp, procnode, slavef)
integer function mumps_get_pool_length(max_active_nodes, keep, keep8)
subroutine mumps_init_pool_dist_bwd(n, nb_prun_roots, pruned_roots, myroot, myid_nodes, keep, keep8, step, procnode_steps, ipool, lpool)
subroutine mumps_sort_int(n, val, id)
subroutine mumps_subtri8toarray(int_array, i8)
subroutine mumps_init_pool_dist_na_bwdl0es(n, myroot, myid_nodes, na, lna, keep, keep8, step, procnode_steps, ipool, lpool, l0_omp_mapping, to_process)
subroutine mumps_icopy_64to32_64c_ip(in_out_tab48, sizetab)
subroutine mumps_ldltpanel_nbtarget(npiv, nb_target, keep)
integer function mumps_typenode(procinfo_inode, k199)
integer function mumps_procnode(procinfo_inode, k199)
subroutine mumps_set_ierror(size8, ierror)
subroutine mumps_icopy_64to32_64c(intab8, sizetab, outtab)
subroutine mumps_init_pool_dist_nona(n, leaf, myid_nodes, lleaves, leaves, keep, keep8, step, procnode_steps, ipool, lpool)
subroutine mumps_set_parti_regular(slavef, keep, keep8, procs, mem_distrib, ncb, nfront, nslaves_node, tab_pos, slaves_list, size_slaves_list, myid, inode, tab_maxs_arg, sup_proc_arg, max_surf, nb_row_max)
subroutine mumps_get_nnz_internal(nnz, nz, nnz_i)
subroutine mumps_init_pool_dist_na_bwd(n, myroot, myid_nodes, na, lna, keep, keep8, step, procnode_steps, ipool, lpool)
integer function mumps_typesplit(procinfo_inode, k199)
logical function mumps_seqana_avail(icntl7)
subroutine mumps_set_ssarbr_dad(ssarbr, inode, dad, n, keep28, step, procnode_steps, k199)
subroutine mumps_reducei8(in, out, mpi_op, root, comm)
subroutine mumps_ldltpanel_sizes(npiv, keep, iw, panel_sizes, nbpanels)
integer function mumps_typenode_rough(procinfo_inode, k199)
subroutine mumps_geti8(i8, int_array)
subroutine mumps_mem_centralize(myid, comm, info, infog, irank)
logical function mumps_compare_tab(tab1, tab2, len1, len2)
subroutine mumps_icopy_32to64_64c_ip(in_out_tab48, sizetab)
subroutine mumps_abort()
subroutine mumps_icopy_64to32(intab8, sizetab, outtab)
subroutine mumps_destroy_arch_node_comm(arch_node_comm)
subroutine mumps_dm_fac_upd_dyn_memcnts(mem_count_allocated, atomic_updates, keep8, iflag, ierror, k69upd, k71upd)
subroutine mumps_find_unit(iunit)
subroutine mumps_abort_on_overflow(i8, string)
subroutine mumps_secdeb(t)
subroutine mumps_init_nroot_dist(n, nbroot, nroot_loc, myid_nodes, slavef, na, lna, keep, step, procnode_steps)
subroutine mumps_get_proc_per_node(k414, myid, nbprocs, comm)
logical function mumps_i_am_candidate(myid, slavef, inode, nmb_par2, istep_to_iniv2, k71, step, n, candidates, keep24)
integer function mumps_encode_tpn_iproc(tpn, iproc, k199)
subroutine mumps_ldltpanel_storage(npiv, keep, iw, nb_entries)
subroutine mumps_icopy_32to64(intab, sizetab, outtab8)
subroutine mumps_typeandprocnode(tpn, mumps_procnode, procinfo_inode, k199)
subroutine mumps_build_arch_node_comm(comm, newcomm, newsize, newrank)
subroutine mumps_sort_doubles(n, val, id)
subroutine mumps_sort_int_dec(n, val, id)
subroutine mumps_make1root(n, frere, fils, nfsiz, theroot)
subroutine mumps_npiv_critical_path(n, nsteps, step, frere, fils, na, lna, ne, maxnpivtree)