OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
cmumps_ana_aux_m Module Reference

Functions/Subroutines

subroutine cmumps_ana_f (n, nz8, irn, icn, liwalloc, ikeep1, ikeep2, ikeep3, iord, nfsiz, fils, frere, listvar_schur, size_schur, icntl, info, keep, keep8, nslaves, piv, cntl4, colsca, rowsca, norig_arg, sizeofblocks, gcomp_provided_in, gcomp)
subroutine cmumps_ana_n_dist (id, ptrar)
subroutine cmumps_ana_o (n, nz, mtrans, perm, ikeepalloc, idirn, idjcn, ida, idrowsca, idcolsca, work2, keep, icntl, info, infog)

Function/Subroutine Documentation

◆ cmumps_ana_f()

subroutine cmumps_ana_aux_m::cmumps_ana_f ( integer, intent(in) n,
integer(8), intent(in) nz8,
integer, dimension(:), pointer irn,
integer, dimension(:), pointer icn,
integer(8), intent(in) liwalloc,
integer, dimension(:), intent(inout) ikeep1,
integer, dimension(:), intent(inout) ikeep2,
integer, dimension(:), intent(inout) ikeep3,
integer, intent(inout) iord,
integer, dimension(:), intent(out) nfsiz,
integer, dimension(:), intent(out) fils,
integer, dimension(:), intent(out) frere,
integer, dimension(:), intent(in) listvar_schur,
integer, intent(in) size_schur,
integer, dimension(60), intent(in) icntl,
integer, dimension(80), intent(inout) info,
integer, dimension(500), intent(inout) keep,
integer(8), dimension(150), intent(inout) keep8,
integer, intent(in) nslaves,
integer, dimension(:), intent(inout) piv,
real cntl4,
real, dimension(:), pointer colsca,
real, dimension(:), pointer rowsca,
integer, intent(in), optional norig_arg,
integer, dimension(n), intent(in), optional sizeofblocks,
logical, intent(in), optional gcomp_provided_in,
type(compact_graph_t), optional gcomp )

Definition at line 22 of file cana_aux.F.

34 IMPLICIT NONE
35 INTEGER, INTENT(IN) :: N, SIZE_SCHUR, NSLAVES
36 INTEGER(8), INTENT(IN) :: NZ8
37 INTEGER(8), INTENT(IN) :: LIWALLOC
38 INTEGER, INTENT(in) :: LISTVAR_SCHUR(:)
39 INTEGER, POINTER :: IRN(:), ICN(:)
40 INTEGER, INTENT(IN) :: ICNTL(60)
41 INTEGER, INTENT(INOUT) :: IORD
42 INTEGER, INTENT(INOUT) :: INFO(80), KEEP(500)
43 INTEGER(8), INTENT(INOUT) :: KEEP8(150)
44 INTEGER, INTENT(OUT) :: NFSIZ(:), FILS(:), FRERE(:)
45 INTEGER, INTENT(INOUT) :: PIV(:)
46 INTEGER, INTENT(INOUT) :: IKEEP1(:), IKEEP2(:), IKEEP3(:)
47 REAL :: CNTL4
48 REAL, POINTER :: COLSCA(:), ROWSCA(:)
49#if defined(metis) || defined(parmetis) || defined(metis4) || defined(parmetis3)
50 INTEGER, INTENT(IN) :: METIS_OPTIONS(40)
51#endif
52 INTEGER, INTENT(IN), OPTIONAL :: NORIG_ARG
53 INTEGER, INTENT(IN), OPTIONAL :: SIZEOFBLOCKS(N)
54 LOGICAL, INTENT(IN), OPTIONAL :: GCOMP_PROVIDED_IN
55 TYPE(COMPACT_GRAPH_T), OPTIONAL :: GCOMP
56 INTEGER, DIMENSION(:), ALLOCATABLE, TARGET :: IWALLOC
57 INTEGER, DIMENSION(:), POINTER :: IW
58 INTEGER(8), DIMENSION(:), ALLOCATABLE, TARGET :: IPEALLOC
59 INTEGER(8), DIMENSION(:), POINTER :: IPE
60 INTEGER(8), DIMENSION(:), ALLOCATABLE :: IPQ8
61 INTEGER, DIMENSION(:,:), ALLOCATABLE :: PTRAR
62 INTEGER, DIMENSION(:), ALLOCATABLE :: PARENT
63 INTEGER, DIMENSION(:), ALLOCATABLE :: IWL1
64 INTEGER NBBUCK
65 INTEGER, DIMENSION(:), ALLOCATABLE :: WTEMP
66 INTEGER IERR
67 INTEGER I, K, NCMPA, IN, IFSON
68 INTEGER(8) :: J8, I8
69 INTEGER :: NORIG
70 INTEGER(8) :: IFIRST, ILAST
71 INTEGER(8) IWFR8
72 INTEGER NEMIN, LP, MP, LDIAG, ITEMP, symmetry
73 INTEGER NBQD, AvgDens
74 LOGICAL PROK, COMPRESS_SCHUR, LPOK, COMPUTE_PERM
75#if defined(metis) || defined(parmetis) || defined(metis4) || defined(parmetis3)
76#if defined(metis4) || defined(parmetis3)
77 INTEGER NUMFLAG
78#endif
79 INTEGER METIS_IDX_SIZE
80 INTEGER OPT_METIS_SIZE
81#endif
82#if defined(scotch) || defined(ptscotch)
83 INTEGER :: SCOTCH_INT_SIZE
84#endif
85#if defined(pord)
86 INTEGER :: PORD_INT_SIZE
87#endif
88 REAL, DIMENSION(:), ALLOCATABLE :: COLSCA_TEMP
89 INTEGER THRESH, IVersion
90 LOGICAL AGG6
91 INTEGER MINSYM
92 parameter(minsym=50)
93 INTEGER(8) :: K79REF
94 parameter(k79ref=12000000_8)
95 INTEGER, PARAMETER :: LIDUMMY = 1
96 INTEGER :: IDUMMY(1)
97 INTEGER MTRANS, COMPRESS,NCMP,IERROR,J,JPERM,NCST
98 INTEGER TOTEL
99#if defined(pord)
100 INTEGER TOTW
101#endif
102 INTEGER WEIGHTUSED
103#if defined(scotch) || defined(ptscotch)
104 INTEGER WEIGHTREQUESTED
105#endif
106 LOGICAL SCOTCH_SYMBOLIC
107 LOGICAL IDENT,SPLITROOT
108 LOGICAL FREE_CENTRALIZED_MATRIX
109 LOGICAL GCOMP_PROVIDED
110 LOGICAL INPLACE64_GRAPH_COPY, INPLACE64_RESTORE_GRAPH
111 INTEGER(8) :: LIW8, NZG8
112 DOUBLE PRECISION TIMEB
113 EXTERNAL mumps_ana_h, cmumps_ana_j,
116#if defined(OLDDFS)
117 EXTERNAL cmumps_ana_l
118#endif
119 EXTERNAL cmumps_gnew_schur
122 IF (liwalloc.GT.0_8) THEN
123 ALLOCATE( iwalloc(liwalloc), stat = ierr )
124 IF ( ierr .GT. 0 ) THEN
125 info( 1 ) = -7
126 CALL mumps_set_ierror(liwalloc,info(2))
127 GOTO 90
128 ENDIF
129 ENDIF
130 ALLOCATE( iwl1(n), stat = ierr )
131 IF ( ierr .GT. 0 ) THEN
132 info( 1 ) = -7
133 info( 2 ) = n
134 GOTO 90
135 ENDIF
136 ALLOCATE( ipealloc(n+1), stat = ierr )
137 IF ( ierr .GT. 0 ) THEN
138 info( 1 ) = -7
139 info( 2 ) = (n+1)*keep(10)
140 GOTO 90
141 ENDIF
142 ALLOCATE( ptrar(n,3), stat = ierr )
143 IF ( ierr .GT. 0 ) THEN
144 info( 1 ) = -7
145 info( 2 ) = 3*n
146 GOTO 90
147 ENDIF
148 scotch_symbolic=(keep(270).EQ.0)
149 symmetry = info(8)
150 nbqd = 0
151 gcomp_provided=.false.
152 weightused = 0
153 norig = n
154 IF (present(norig_arg)) THEN
155 norig=norig_arg
156 ENDIF
157 IF (present(gcomp_provided_in))
158 & gcomp_provided = gcomp_provided_in
159 IF (gcomp_provided.AND.(.NOT. present(gcomp))) THEN
160 info(1) = -900
161 WRITE(6,*) " INTERNAL ERROR in MUMPS(ANA_F) ",
162 & gcomp_provided_in, present(gcomp)
163 info(2) = 1
164 RETURN
165 ENDIF
166 IF ( (liwalloc.EQ.0_8).AND.(.not.gcomp_provided)) THEN
167 info(1) = -900
168 WRITE(6,*) " INTERNAL ERROR in MUMPS(ANA_F) ",
169 & "LIWALLOC, GCOMP_PROVIDED=", liwalloc, gcomp_provided
170 info(2) = 2
171 RETURN
172 ENDIF
173 IF (gcomp_provided) THEN
174 nzg8 = gcomp%NZG
175 liw8 = nzg8 + int(gcomp%NG,8)+1_8
176 iw => gcomp%ADJ(1:liw8)
177 ipe => gcomp%IPE(1:gcomp%NG+1)
178 DO i=1,gcomp%NG
179 ptrar(i,2) = int(ipe(i+1)-ipe(i))
180 ENDDO
181 ELSE
182 liw8 = liwalloc
183 nzg8 = nz8
184 iw => iwalloc(1:liw8)
185 ipe => ipealloc(1:n+1)
186 ENDIF
187 lp = icntl(1)
188 mp = icntl(3)
189 lpok = ((lp.GT.0).AND.(icntl(4).GE.1))
190 prok = ((mp.GT.0).AND.(icntl(4).GE.2))
191 ldiag = icntl(4)
192 compress_schur = .false.
193 IF (prok) THEN
194 IF (present(gcomp)) THEN
195 WRITE(mp,'(A,I10,A,I13,A)') " Processing a graph of size:", n
196 & ," with ", gcomp%NZG, " edges"
197 ELSE
198 WRITE(mp,'(A,I10)') " Processing a graph of size:", n
199 ENDIF
200 ENDIF
201 IF (gcomp_provided) THEN
202 free_centralized_matrix = .false.
203 ELSE
204 free_centralized_matrix = (
205 & (keep(54).EQ.3).AND.
206 & (keep(494).EQ.0).AND.
207 & (keep(106).NE.3)
208 & )
209 ENDIF
210 inplace64_graph_copy = .false.
211 inplace64_restore_graph = .true.
212 IF (keep(1).LT.0) keep(1) = 0
213 nemin = keep(1)
214 IF (ldiag.GT.2 .AND. mp.GT.0) THEN
215 IF (present(sizeofblocks)) THEN
216 k = min(10,gcomp%NG)
217 IF (ldiag.EQ.4) k = gcomp%NG
218 WRITE (mp,99909) n, nzg8, info(1)
219 i8= 0_8
220 WRITE(mp,'(A)') " Graph adjacency "
221 DO j=1, k
222 ifirst = gcomp%IPE(j)
223 ilast= min(gcomp%IPE(j+1)-1,gcomp%IPE(j)+k-1)
224 write(mp,'(A,I10)') " .... node/column:", j
225 write(mp,'(8X,10I9)')
226 & (gcomp%ADJ(i8),i8=ifirst,ilast)
227 ENDDO
228 ELSE
229 j8 = min(nzg8, 10_8)
230 IF (ldiag .EQ.4) j8 = nzg8
231 WRITE (mp,99999) n, nzg8, liw8, info(1)
232 IF (j8.GT.0_8) WRITE (mp,99998) (irn(i8),icn(i8),i8=1_8,j8)
233 ENDIF
234 k = min0(10,n)
235 IF (ldiag.EQ.4) k = n
236 IF (iord.EQ.1 .AND. k.GT.0) THEN
237 WRITE (mp,99997) (ikeep1(i),i=1,k)
238 ENDIF
239 ENDIF
240 ncmp = n
241 IF (keep(60).NE.0) THEN
242 IF ((size_schur.LE.0 ).OR.
243 & (size_schur.GE.n) ) GOTO 90
244 ENDIF
245#if defined(metis) || defined(parmetis) || defined(metis4) || defined(parmetis3)
246 IF ( ( keep(60).NE.0).AND.(size_schur.GT.0)
247 & .AND.
248 & ((iord.EQ.7).OR.(iord.EQ.5))
249 & )THEN
250 compress_schur=.true.
251 ncmp = n-size_schur
252 ALLOCATE(ipq8(n),stat=ierr)
253 IF ( ierr .GT. 0 ) THEN
254 info( 1 ) = -7
255 info( 2 ) = n*keep(10)
256 ENDIF
257 CALL cmumps_gnew_schur(n,ncmp,nz8,irn(1), icn(1), iw(1), liw8,
258 & ipe(1), ptrar(1,2),
259 & ipq8, iwl1, iwfr8, keep8(126), keep8(127),
260 & info(1), info(2), icntl, symmetry,
261 & keep(50), nbqd, avgdens,
262 & keep(264), keep(265),
263 & listvar_schur(1), size_schur, frere(1), fils(1),
264 & inplace64_graph_copy)
265 info(8) = symmetry
266 inplace64_graph_copy = inplace64_graph_copy.AND.
267 & (.NOT.free_centralized_matrix)
268 DEALLOCATE(ipq8)
269 iord = 5
270 keep(95) = 1
271 nbqd = 0
272 ELSE
273#endif
274 IF (gcomp_provided) THEN
275 iwfr8 = gcomp%NZG+1_8
276 ELSE
277 ALLOCATE(ipq8(n),stat=ierr)
278 IF ( ierr .GT. 0 ) THEN
279 info( 1 ) = -7
280 info( 2 ) = n*keep(10)
281 ENDIF
282 IF ( prok ) THEN
283 CALL mumps_secdeb( timeb )
284 ENDIF
285 CALL cmumps_ana_gnew(n,nz8,irn(1), icn(1), iw(1), liw8,
286 & ipe(1), ptrar(1,2),
287 & ipq8, iwl1, iwfr8, keep8(126), keep8(127),
288 & info(1), info(2), icntl, symmetry,
289 & keep(50), nbqd, avgdens, keep(264), keep(265),
290 & .true., inplace64_graph_copy)
291 info(8) = symmetry
292 inplace64_graph_copy = inplace64_graph_copy.AND.
293 & (.NOT.free_centralized_matrix)
294 DEALLOCATE(ipq8)
295 ENDIF
296#if defined(metis) || defined(parmetis) || defined(metis4) || defined(parmetis3)
297 ENDIF
298#endif
299 IF(nbqd .GT. 0) THEN
300 IF( keep(50) .EQ. 2 .AND. icntl(12) .EQ. 0 ) THEN
301 IF(keep(95) .NE. 1) THEN
302 IF ( prok )
303 & WRITE( mp,*)
304 & 'Compressed/constrained ordering set OFF'
305 keep(95) = 1
306 ENDIF
307 ENDIF
308 ENDIF
309 IF ( (keep(60).NE.0) .AND. (iord.GT.1) .AND.
310 & .NOT. compress_schur ) THEN
311 iord = 0
312 ENDIF
313 IF ( (keep(50).EQ.2)
314 & .AND. (keep(95) .EQ. 3)
315 & .AND. (iord .EQ. 7) ) THEN
316 iord = 2
317 ENDIF
318 CALL mumps_set_ordering( norig, keep,
319 & keep(50), nslaves, iord,
320 & nbqd, avgdens,
321 & prok, mp )
322 IF(keep(50) .EQ. 2) THEN
323 IF(keep(95) .EQ. 3 .AND. iord .NE. 2) THEN
324 IF (prok) WRITE(mp,*)
325 & 'WARNING: CMUMPS_ANA_F constrained ordering not '//
326 & ' available with selected ordering. Move to' //
327 & ' compressed ordering.'
328 keep(95) = 2
329 ENDIF
330 ELSE
331 keep(95) = 1
332 ENDIF
333 mtrans = keep(23)
334 compress = keep(95) - 1
335 IF(compress .GT. 0 .AND. keep(52) .EQ. -2) THEN
336 IF(cntl4 .GE. 0.0e0) THEN
337 IF (keep(1).LE.8) THEN
338 nemin = 16
339 ELSE
340 nemin = 2*keep(1)
341 ENDIF
342 ENDIF
343 ENDIF
344 IF(mtrans .GT. 0 .AND. keep(50) .EQ. 2) THEN
345 keep(23) = 0
346 ENDIF
347 IF (compress .EQ. 2) THEN
348 IF (iord.NE.2) THEN
349 WRITE(*,*) "IORD not compatible with COMPRESS:",
350 & iord, compress
351 CALL mumps_abort()
352 ENDIF
354 & n,piv(1),frere(1),fils(1),nfsiz(1),ikeep1(1),
355 & ncst,keep,keep8, rowsca(1)
356 & )
357 ENDIF
358 IF ( iord .NE. 1 ) THEN
359 IF (compress .GE. 1) THEN
360 ALLOCATE(ipq8(n),stat=ierr)
361 IF ( ierr .GT. 0 ) THEN
362 info( 1 ) = -7
363 info( 2 ) = n*keep(10)
364 ENDIF
366 & n, nz8, irn(1), icn(1), piv(1),
367 & ncmp, iw(1), liw8, ipe(1), ptrar(1,2), ipq8,
368 & iwl1, fils(1), iwfr8,
369 & ierror, keep, keep8, icntl, inplace64_graph_copy)
370 DEALLOCATE(ipq8)
371 symmetry = 100
372 ENDIF
373 IF ( (symmetry.LT.minsym).AND.(keep(50).EQ.0) ) THEN
374 IF(keep(23) .EQ. 7 ) THEN
375 keep(23) = -5
376 GOTO 90
377 ELSE IF(keep(23) .EQ. -9876543) THEN
378 ident = .true.
379 keep(23) = 5
380 IF (prok) WRITE(mp,'(A)')
381 & ' ... Apply column permutation (already computed)'
382 DO j=1,n
383 jperm = piv(j)
384 fils(jperm) = j
385 IF (jperm.NE.j) ident = .false.
386 ENDDO
387 IF (.NOT.ident) THEN
388 DO j8=1_8,nz8
389 j = icn(j8)
390 IF ((j.LE.0).OR.(j.GT.n)) cycle
391 icn(j8) = fils(j)
392 ENDDO
393 ALLOCATE(colsca_temp(n), stat=ierr)
394 IF ( ierr > 0 ) THEN
395 info( 1 ) = -7
396 info( 2 ) = n
397 GOTO 90
398 ENDIF
399 DO j = 1, n
400 colsca_temp(j)=colsca(j)
401 ENDDO
402 DO j=1, n
403 colsca(fils(j))=colsca_temp(j)
404 ENDDO
405 DEALLOCATE(colsca_temp)
406 IF (prok)
407 & WRITE(mp,'(/A)')
408 & ' WARNING input matrix data modified'
409 ALLOCATE(ipq8(n),stat=ierr)
410 IF ( ierr .GT. 0 ) THEN
411 info( 1 ) = -7
412 info( 2 ) = n*keep(10)
413 ENDIF
414 CALL cmumps_ana_gnew
415 & (n,nz8,irn(1), icn(1), iw(1), liw8,
416 & ipe(1), ptrar(1,2),
417 & ipq8, iwl1, iwfr8, keep8(126), keep8(127),
418 & info(1), info(2), icntl, symmetry, keep(50),
419 & nbqd, avgdens, keep(264), keep(265),
420 & .true.,inplace64_graph_copy)
421 info(8) = symmetry
422 DEALLOCATE(ipq8)
423 ncmp = n
424 ELSE
425 keep(23) = 0
426 ENDIF
427 ENDIF
428 ELSE IF (keep(23) .EQ. 7 .OR. keep(23) .EQ. -9876543 ) THEN
429 IF (prok) WRITE(mp,'(A)')
430 & ' ... No column permutation'
431 keep(23) = 0
432 ENDIF
433 ENDIF
434 IF (free_centralized_matrix
435 & .AND.compress.EQ.0.AND.(.NOT.compress_schur)) THEN
436 deallocate(irn)
437 NULLIFY(irn)
438 deallocate(icn)
439 NULLIFY(icn)
440 ENDIF
441 inplace64_restore_graph =
442 & inplace64_restore_graph.AND.(compress.NE.1)
443 ALLOCATE( parent( n ), stat = ierr )
444 IF ( ierr .GT. 0 ) THEN
445 info( 1 ) = -7
446 info( 2 ) = n
447 GOTO 90
448 ENDIF
449 IF (iord.NE.1 .AND. iord.NE.5) THEN
450 IF ( keep(60) .NE. 0 ) THEN
451 iord = 0
452 ENDIF
453 IF (prok) THEN
454 IF (iord.EQ.2) THEN
455 WRITE(mp,'(A)') ' Ordering based on AMF '
456#if defined(scotch) || defined(ptscotch)
457 ELSE IF (iord.EQ.3) THEN
458 WRITE(mp,'(A)') ' Ordering based on SCOTCH '
459#endif
460#if defined(pord)
461 ELSE IF (iord.EQ.4) THEN
462 WRITE(mp,'(A)') ' Ordering based on PORD '
463#endif
464 ELSE IF (iord.EQ.6) THEN
465 WRITE(mp,'(A)') ' Ordering based on QAMD '
466 ELSE
467 WRITE(mp,'(A)') ' Ordering based on AMD '
468 ENDIF
469 ENDIF
470 IF ( prok ) THEN
471 CALL mumps_secdeb( timeb )
472 ENDIF
473 IF ( keep(60) .NE. 0 ) THEN
474 CALL mumps_hamd(n, liw8, ipe(1), iwfr8, ptrar(1,2), iw(1),
475 & iwl1, ikeep1(1),
476 & ikeep2(1), ncmpa, fils(1), ikeep3(1),
477 & ptrar, ptrar(1,3),
478 & parent,
479 & listvar_schur(1), size_schur)
480 IF (keep(60)==1) THEN
481 keep(20) = listvar_schur(1)
482 ELSE
483 keep(38) = listvar_schur(1)
484 ENDIF
485 ELSE
486 IF ( .false. ) THEN
487#if defined(pord)
488 ELSEIF (iord .EQ. 4) THEN
489 CALL mumps_pord_intsize(pord_int_size)
490 totw = n
491 IF ( (compress .EQ. 1)
492 & .OR.
493 & ( (norig.NE.n).AND.present(sizeofblocks) )
494 & ) THEN
495 IF (compress .EQ. 1) THEN
496 DO i=1,keep(93)/2
497 iwl1(i) = 2
498 ENDDO
499 DO i=1+keep(93)/2,ncmp
500 iwl1(i) = 1
501 ENDDO
502 ELSE IF
503 & ( (norig.NE.n).AND.present(sizeofblocks) ) THEN
504 totw = norig
505 DO i= 1, n
506 iwl1(i) = sizeofblocks(i)
507 ENDDO
508 ENDIF
509 IF (pord_int_size .EQ. 64) THEN
510 CALL mumps_pordf_wnd_mixedto64(ncmp, iwfr8-1_8,
511 & ipe, iw,
512 & iwl1, ncmpa, totw, parent,
513 & info(1), lp, lpok, keep(10),
514 & inplace64_graph_copy
515 & )
516 ELSE IF (pord_int_size .EQ. 32) THEN
517 CALL mumps_pordf_wnd_mixedto32(ncmp, iwfr8-1_8,
518 & ipe, iw,
519 & iwl1, ncmpa, totw, parent,
520 & info(1), lp, lpok, keep(10))
521 ELSE
522 WRITE(*,*)
523 & "Internal error in PORD wrappers, PORD_INT_SIZE=",
524 & pord_int_size
525 CALL mumps_abort()
526 ENDIF
527 IF ( ncmpa .NE. 0 ) THEN
528 write(6,*) ' Out PORD, NCMPA=', ncmpa
529 info( 1 ) = -9999
530 info( 2 ) = 4
531 GOTO 90
532 ENDIF
533 IF (info(1) .LT.0) GOTO 90
534 IF (compress.EQ.1) THEN
535 CALL cmumps_get_elim_tree(ncmp,parent,iwl1,fils(1))
536 CALL cmumps_get_perm_from_pe(ncmp,parent,ikeep1(1),
537 & frere(1),ptrar(1,1))
538 DO i=1,ncmp
539 ikeep2(ikeep1(i))=i
540 ENDDO
541 ENDIF
542 ELSE
543 IF (pord_int_size.EQ.64) THEN
544 CALL mumps_pordf_mixedto64(ncmp, iwfr8-1_8, ipe,
545 & iw,
546 & iwl1, ncmpa, parent,
547 & info(1), lp, lpok, keep(10),
548 & inplace64_graph_copy
549 & )
550 ELSE IF (pord_int_size.EQ.32) THEN
551 CALL mumps_pordf_mixedto32(ncmp, iwfr8-1_8, ipe,
552 & iw,
553 & iwl1, ncmpa, parent,
554 & info(1), lp, lpok, keep(10))
555 ELSE
556 WRITE(*,*)
557 & "Internal error in PORD wrappers, PORD_INT_SIZE=",
558 & pord_int_size
559 CALL mumps_abort()
560 ENDIF
561 ENDIF
562 IF ( ncmpa .NE. 0 ) THEN
563 write(6,*) ' Out PORD, NCMPA=', ncmpa
564 info( 1 ) = -9999
565 info( 2 ) = 4
566 GOTO 90
567 ENDIF
568 IF (info(1) .LT. 0) GOTO 90
569#endif
570#if defined(scotch) || defined(ptscotch)
571 ELSEIF (iord .EQ. 3) THEN
572 CALL mumps_scotch_intsize(scotch_int_size)
573 IF ( (compress .EQ. 1)
574 & .OR.
575 & ( (norig.NE.n).AND.present(sizeofblocks) )
576 & ) THEN
577 weightrequested=1
578 IF (compress .EQ. 1) THEN
579 DO i=1,keep(93)/2
580 iwl1(i) = 2
581 ENDDO
582 DO i=1+keep(93)/2,ncmp
583 iwl1(i) = 1
584 ENDDO
585 ELSE IF
586 & ( (norig.NE.n).AND.present(sizeofblocks) ) THEN
587 DO i= 1, n
588 iwl1(i) = sizeofblocks(i)
589 ENDDO
590 ENDIF
591 ELSE
592 weightrequested = 0
593 DO i= 1, n
594 iwl1(i) = 1
595 ENDDO
596 ENDIF
597 IF (scotch_int_size.EQ.32) THEN
598 IF (keep(10).EQ.1) THEN
599 info(1) = -52
600 info(2) = 2
601 ELSE
602 CALL mumps_scotch_mixedto32(ncmp,
603 & iwfr8-1_8, ipe,
604 & parent, iwfr8,
605 & ptrar(1,2), iw, iwl1, ikeep1,
606 & ikeep2, ncmpa, info, lp, lpok,
607 & weightused, weightrequested, scotch_symbolic)
608 ENDIF
609 ELSE IF (scotch_int_size.EQ.64) THEN
610 CALL mumps_scotch_mixedto64(ncmp,
611 & iwfr8-1_8, ipe,
612 & parent, iwfr8,
613 & ptrar(1,2), iw, iwl1, ikeep1,
614 & ikeep2, ncmpa, info, lp, lpok, keep(10),
615 & inplace64_graph_copy,
616 & weightused, weightrequested, scotch_symbolic)
617 ELSE
618 WRITE(*,*)
619 & "Internal error in SCOTCH wrappers, SCOTCH_INT_SIZE=",
620 & scotch_int_size
621 CALL mumps_abort()
622 ENDIF
623 IF (info(1) .LT. 0) GOTO 90
624 IF (.NOT. scotch_symbolic) THEN
625 IF ( compress .EQ. 1 ) THEN
626 CALL cmumps_expand_permutation(n,ncmp,keep(94),
627 & keep(93),piv(1),ikeep1(1),ikeep2(1))
628 compress = -1
629 ENDIF
630 ELSE IF ( (compress .EQ. 1)
631 & .OR.
632 & ( (norig.NE.n).AND.present(sizeofblocks).AND.
633 & (weightused.EQ.0) )
634 & ) THEN
635 CALL cmumps_get_elim_tree(ncmp,parent,iwl1,fils(1))
636 CALL cmumps_get_perm_from_pe(ncmp,parent,ikeep1(1),
637 & frere(1),ptrar(1,1))
638 DO i=1,ncmp
639 ikeep2(ikeep1(i))=i
640 ENDDO
641 ENDIF
642#endif
643 ELSEIF (iord .EQ. 2) THEN
644 nbbuck = 2*n
645 compute_perm=.false.
646 IF(compress .GE. 1) THEN
647 compute_perm=.true.
648 DO i=1,keep(93)/2
649 iwl1(i) = 2
650 ENDDO
651 DO i=1+keep(93)/2,ncmp
652 iwl1(i) = 1
653 ENDDO
654 totel = keep(93)+keep(94)
655 ELSE
656 iwl1(1) = -1
657 totel = n
658 ENDIF
659 IF (present(sizeofblocks)) THEN
660 IF (compress.GE.1) THEN
661 CALL mumps_abort()
662 ENDIF
663 nbbuck = max(nbbuck, norig-n)
664 nbbuck = max(nbbuck, 2*norig)
665 ncmp = n
666 totel = norig
667 DO i= 1, n
668 iwl1(i) = sizeofblocks(i)
669 ENDDO
670 ENDIF
671 ALLOCATE( wtemp( 0: nbbuck + 1), stat = ierr )
672 IF ( ierr .GT. 0 ) THEN
673 info( 1 ) = -7
674 info( 2 ) = nbbuck+2
675 GOTO 90
676 ENDIF
677 IF(compress .LE. 1) THEN
678 CALL mumps_hamf4
679 & (totel, ncmp, compute_perm, nbbuck, liw8, ipe(1),
680 & iwfr8, ptrar(1,2),
681 & iw(1), iwl1, ikeep1(1), ikeep2(1), ncmpa, fils(1),
682 & ikeep3(1), ptrar, ptrar(1,3), wtemp, parent(1))
683 ELSE
684 IF(prok) WRITE(mp,'(A)')
685 & ' Constrained Ordering based on AMF'
686 CALL mumps_cst_amf(ncmp, nbbuck, liw8, ipe(1),
687 & iwfr8, ptrar(1,2),
688 & iw(1), iwl1, ikeep1(1), ikeep2(1), ncmpa, fils(1),
689 & ikeep3(1), ptrar, ptrar(1,3), wtemp,
690 & nfsiz(1), frere(1), parent(1))
691 ENDIF
692 DEALLOCATE(wtemp)
693 ELSEIF (iord .EQ. 6) THEN
694 ALLOCATE( wtemp( n ), stat = ierr )
695 IF ( ierr .GT. 0 ) THEN
696 info( 1 ) = -7
697 info( 2 ) = n
698 GOTO 90
699 ENDIF
700 thresh = 1
701 iversion = 2
702 compute_perm=.false.
703 IF(compress .EQ. 1) THEN
704 compute_perm=.true.
705 DO i=1,keep(93)/2
706 iwl1(i) = 2
707 ENDDO
708 DO i=1+keep(93)/2,ncmp
709 iwl1(i) = 1
710 ENDDO
711 totel = keep(93)+keep(94)
712 ELSE
713 iwl1(1) = -1
714 totel = n
715 ENDIF
716 IF (present(sizeofblocks)) THEN
717 IF (compress.EQ.1) THEN
718 CALL mumps_abort()
719 ENDIF
720 ncmp = n
721 totel = norig
722 DO i= 1, n
723 iwl1(i) = sizeofblocks(i)
724 ENDDO
725 ENDIF
726 CALL mumps_qamd
727 & (totel,compute_perm,iversion, thresh, wtemp,
728 & ncmp, liw8, ipe(1), iwfr8, ptrar(1,2), iw(1),
729 & iwl1, ikeep1(1), ikeep2(1), ncmpa, fils(1),
730 & ikeep3(1), ptrar, ptrar(1,3), parent(1))
731 DEALLOCATE(wtemp)
732 ELSE
733 compute_perm=.false.
734 IF(compress .EQ. 1) THEN
735 compute_perm=.true.
736 DO i=1,keep(93)/2
737 iwl1(i) = 2
738 ENDDO
739 DO i=1+keep(93)/2,ncmp
740 iwl1(i) = 1
741 ENDDO
742 totel = keep(93)+keep(94)
743 ELSE
744 iwl1(1) = -1
745 totel = n
746 ENDIF
747 IF (present(sizeofblocks)) THEN
748 IF (compress.EQ.1) THEN
749 CALL mumps_abort()
750 ENDIF
751 ncmp = n
752 totel = norig
753 DO i= 1, n
754 iwl1(i) = sizeofblocks(i)
755 ENDDO
756 ENDIF
757 CALL mumps_ana_h(totel, compute_perm,
758 & ncmp, liw8, ipe(1), iwfr8, ptrar(1,2),
759 & iw(1), iwl1, ikeep1(1), ikeep2(1), ncmpa, fils(1),
760 & ikeep3(1), ptrar, ptrar(1,3), parent(1))
761 ENDIF
762 ENDIF
763 IF(compress .GE. 1) THEN
764 CALL cmumps_expand_permutation(n,ncmp,keep(94),keep(93),
765 & piv(1),ikeep1(1),ikeep2(1))
766 compress = -1
767 ENDIF
768 IF ( prok ) THEN
769 CALL mumps_secfin( timeb )
770#if defined(scotch) || defined(ptscotch)
771 IF (iord.EQ.3) THEN
772 WRITE( mp, '(A,F12.4)' )
773 & ' ELAPSED TIME SPENT IN SCOTCH reordering =', timeb
774 ENDIF
775#endif
776 ENDIF
777 ENDIF
778#if defined(metis) || defined(parmetis) || defined(metis4) || defined(parmetis3)
779 IF (iord.EQ.5) THEN
780 IF (prok) THEN
781 WRITE(mp,'(A)') ' Ordering based on METIS'
782 ENDIF
783 IF ( prok ) THEN
784 CALL mumps_secdeb( timeb )
785 ENDIF
786 CALL mumps_metis_idxsize(metis_idx_size)
787 IF (keep(10).EQ.1.AND.metis_idx_size.NE.64) THEN
788 info(1) = -52
789 info(2) = 1
790 GOTO 90
791 ENDIF
792#if defined(metis4) || defined(parmetis3)
793 numflag = 1
794 opt_metis_size = 8
795#else
796 opt_metis_size = 40
797#endif
798 IF (compress .EQ. 1) THEN
799 DO i=1,keep(93)/2
800 frere(i) = 2
801 ENDDO
802 DO i=keep(93)/2+1,ncmp
803 frere(i) = 1
804 ENDDO
805#if defined(metis4) || defined(parmetis3)
806 IF (metis_idx_size .EQ.32) THEN
807 CALL mumps_metis_nodewnd_mixedto32(
808 & ncmp, ipe, iw, frere,
809 & numflag, metis_options(1), opt_metis_size,
810 & ikeep2, ikeep1, info(1), lp, lpok )
811 ELSE IF (metis_idx_size .EQ.64) THEN
812 CALL mumps_metis_nodewnd_mixedto64(
813 & ncmp, ipe, iw, frere,
814 & numflag, metis_options(1), opt_metis_size,
815 & ikeep2, ikeep1, info(1), lp, lpok, keep(10),
816 & inplace64_graph_copy )
817 ELSE
818 WRITE(*,*)
819 & "Internal error in METIS wrappers, METIS_IDX_SIZE=",
820 & metis_idx_size
821 CALL mumps_abort()
822 ENDIF
823 ELSE
824 IF ((norig.NE.n).AND.present(sizeofblocks)) THEN
825 DO i=1, n
826 frere(i) = sizeofblocks(i)
827 ENDDO
828 IF (metis_idx_size .EQ.32) THEN
829 CALL mumps_metis_nodewnd_mixedto32(
830 & ncmp, ipe, iw, frere,
831 & numflag, metis_options(1), opt_metis_size,
832 & ikeep2, ikeep1, info(1), lp, lpok )
833 ELSE IF (metis_idx_size .EQ.64) THEN
834 CALL mumps_metis_nodewnd_mixedto64(
835 & ncmp, ipe, iw, frere,
836 & numflag, metis_options(1), opt_metis_size,
837 & ikeep2, ikeep1, info(1), lp, lpok, keep(10),
838 & inplace64_graph_copy )
839 ELSE
840 WRITE(*,*)
841 & "Internal error in METIS wrappers, METIS_IDX_SIZE=",
842 & metis_idx_size
843 CALL mumps_abort()
844 ENDIF
845 ELSE
846 IF (metis_idx_size .EQ.32) THEN
847 CALL mumps_metis_nodend_mixedto32(
848 & ncmp, ipe, iw, numflag,
849 & metis_options(1), opt_metis_size,
850 & ikeep2, ikeep1, info(1), lp, lpok )
851 ELSE IF (metis_idx_size .EQ.64) THEN
852 CALL mumps_metis_nodend_mixedto64(
853 & ncmp, ipe, iw, numflag,
854 & metis_options(1), opt_metis_size,
855 & ikeep2, ikeep1, info(1), lp,lpok,keep(10),
856 & liw8, inplace64_graph_copy,
857 & inplace64_restore_graph)
858 ELSE
859 WRITE(*,*)
860 & "Internal error in METIS wrappers, METIS_IDX_SIZE=",
861 & metis_idx_size
862 CALL mumps_abort()
863 ENDIF
864 ENDIF
865 ENDIF
866#else
867 ELSE
868 IF (present(sizeofblocks)) THEN
869 DO i=1,n
870 frere(i) = sizeofblocks(i)
871 ENDDO
872 ELSE
873 DO i=1,ncmp
874 frere(i) = 1
875 ENDDO
876 ENDIF
877 ENDIF
878 IF (metis_idx_size .EQ. 32) THEN
879 CALL mumps_metis_nodend_mixedto32(
880 & ncmp, ipe, iw, frere,
881 & metis_options(1), opt_metis_size,
882 & ikeep2, ikeep1, info(1), lp, lpok )
883 ELSE IF (metis_idx_size .EQ. 64) THEN
884 CALL mumps_metis_nodend_mixedto64(
885 & ncmp, ipe, iw, frere,
886 & metis_options(1), opt_metis_size,
887 & ikeep2, ikeep1, info(1), lp,lpok,keep(10),
888 & liw8, inplace64_graph_copy,
889 & inplace64_restore_graph)
890 ELSE
891 IF (lpok) WRITE(lp,*)
892 & "Internal error in METIS wrappers, METIS_IDX_SIZE=",
893 & metis_idx_size
894 CALL mumps_abort()
895 ENDIF
896#endif
897 IF (info(1) .LT.0) GOTO 90
898 IF ( prok ) THEN
899 CALL mumps_secfin( timeb )
900 WRITE( mp, '(A,F12.4)' )
901 & ' ELAPSED TIME SPENT IN METIS reordering =', timeb
902 ENDIF
903 IF ( compress_schur ) THEN
905 & n, ncmp, ikeep1(1),ikeep2(1),
906 & listvar_schur(1), size_schur, fils(1))
907 compress = -1
908 ENDIF
909 IF (compress .EQ. 1) THEN
910 CALL cmumps_expand_permutation(n,ncmp,keep(94),
911 & keep(93),piv(1),ikeep1(1),ikeep2(1))
912 compress = -1
913 ENDIF
914 ENDIF
915#endif
916 IF (prok) THEN
917 IF (iord.EQ.1) THEN
918 WRITE(mp,'(A)') ' Ordering given is used'
919 ENDIF
920 ENDIF
921 IF (iord.EQ.1 .OR. iord.EQ.5 .OR. compress.EQ.-1
922 & .OR. ( (iord.EQ.3).AND.(.NOT.scotch_symbolic) )
923 & .OR.
924 & ( (norig.NE.n).AND.present(sizeofblocks) .AND.(iord.EQ.3)
925 & .AND. (weightused.EQ.0)
926 & )
927 & ) THEN
928 IF ((keep(106).EQ.1).OR.(keep(106).EQ.2).OR.(keep(106).EQ.4)
929 & .OR.(keep(60).NE.0)) THEN
930 IF ( compress .EQ. -1 ) THEN
931 ALLOCATE(ipq8(n),stat=ierr)
932 IF ( ierr .GT. 0 ) THEN
933 info( 1 ) = -7
934 info( 2 ) = n*keep(10)
935 ENDIF
936 CALL cmumps_ana_gnew(n,nz8,irn(1),icn(1),iw(1),liw8,
937 & ipe(1), ptrar(1,2),
938 & ipq8, iwl1, iwfr8, keep8(126), keep8(127),
939 & info(1), info(2), icntl, symmetry, keep(50),
940 & nbqd, avgdens, keep(264),keep(265), .true.,
941 & inplace64_graph_copy)
942 DEALLOCATE(ipq8)
943 ENDIF
944 compress = 0
945 IF (keep(106).EQ.2) THEN
946 IF (prok) THEN
947 WRITE(mp,*) " SYMBOLIC based on column counts "
948 ENDIF
949 IF (present(sizeofblocks)) THEN
950 DO i=1, n
951 frere(i) = sizeofblocks(i)
952 ENDDO
953 ELSE
954 frere(1) = -1
955 ENDIF
956 CALL mumps_wrap_ginp94 (
957 & n, ipe(1), iw(1), iwfr8,
958 & ikeep1(1),
959 & frere(1),
960 & keep(60), listvar_schur(1), size_schur,
961 & keep(378),
962 & iwl1, parent,
963 & ikeep2(1), ikeep3(1), nfsiz(1),
964 & ptrar(1,1), ptrar(1,2), ptrar(1,3),
965 & info )
966 IF (info(1).LT.0) GOTO 90
967 ELSE IF ((keep(106).EQ.4).AND.(keep(60).EQ.0).AND.
968 & (.NOT.present(sizeofblocks) .OR. (norig.EQ.n))
969 & ) THEN
970 WRITE(mp,*) " Undefined option for ICNTL(58) "
971 info(1)= -99998
972 GOTO 90
973 ELSE
974 ALLOCATE( wtemp( 2*n ), stat = ierr )
975 IF ( ierr .GT. 0 ) THEN
976 info( 1 ) = -7
977 info( 2 ) = 2*n
978 GOTO 90
979 ENDIF
980 thresh = -1
981 IF (keep(60) == 0) THEN
982 itemp = 0
983 ELSE
984 itemp = size_schur
985 ENDIF
986 agg6 =.false.
987 IF (present(sizeofblocks)) THEN
988 DO i=1, n
989 iwl1(i) = sizeofblocks(i)
990 ENDDO
991 totel = norig
992 ELSE
993 iwl1(1) = -1
994 totel = n
995 ENDIF
996 CALL mumps_symqamd(thresh, wtemp,
997 & n, totel, liw8, ipe(1), iwfr8, ptrar(1,2), iw(1),
998 & iwl1(1), wtemp(n+1),
999 & ikeep2(1), ncmpa, fils(1), ikeep3(1), ptrar,
1000 & ptrar(1,3),ikeep1(1), listvar_schur(1), itemp,
1001 & agg6, parent)
1002 DEALLOCATE(wtemp)
1003 ENDIF
1004 ELSE
1005 CALL cmumps_ana_j(n, nz8, irn(1), icn(1), ikeep1(1), iw(1),
1006 & liw8, ipe(1),
1007 & ptrar(1,2), iwl1, iwfr8,
1008 & info(1),info(2), mp)
1009 IF (keep(60) .EQ. 0) THEN
1010 itemp = 0
1011 ELSE
1012 itemp = size_schur
1013 ENDIF
1014 CALL cmumps_ana_k(n, ipe(1), iw(1), liw8, iwfr8, ikeep1(1),
1015 & ikeep2(1), iwl1,
1016 & ptrar, ncmpa, itemp, parent)
1017 ENDIF
1018 ENDIF
1019 IF (keep(60) .NE. 0) THEN
1020 IF (keep(60)==1) THEN
1021 keep(20) = listvar_schur(1)
1022 ELSE
1023 keep(38) = listvar_schur(1)
1024 ENDIF
1025 ENDIF
1026#if defined(OLDDFS)
1027 CALL cmumps_ana_l
1028 & (n, parent, iwl1, ikeep1(1), ikeep2(1), ikeep3(1),
1029 & nfsiz, info(6), fils(1), frere(1), ptrar(1,3),
1030 & nemin, keep(60))
1031#else
1032 IF (allocated(ipealloc)) DEALLOCATE(ipealloc)
1033 ALLOCATE(wtemp(n), stat=ierr)
1034 IF ( ierr .GT. 0 ) THEN
1035 info( 1 ) = -7
1036 info( 2 ) = n
1037 GOTO 90
1038 ENDIF
1039 IF (present(sizeofblocks)) THEN
1040 CALL cmumps_ana_lnew
1041 & (n, parent, iwl1, ikeep1(1), ikeep2(1), ikeep3(1),
1042 & nfsiz(1), ptrar, info(6), fils(1), frere(1),
1043 & ptrar(1,3), nemin, wtemp, keep(60),
1044 & keep(20),keep(38),ptrar(1,2),keep(104),iw(1),keep(50),
1045 & icntl(13), keep(37), keep(197), nslaves, keep(250).EQ.1
1046 & , .true. , sizeofblocks, n
1047 & )
1048 ELSE
1049 CALL cmumps_ana_lnew
1050 & (n, parent, iwl1, ikeep1(1), ikeep2(1), ikeep3(1),
1051 & nfsiz(1), ptrar, info(6), fils(1), frere(1),
1052 & ptrar(1,3), nemin, wtemp, keep(60),
1053 & keep(20),keep(38),ptrar(1,2),keep(104),iw(1),keep(50),
1054 & icntl(13), keep(37), keep(197), nslaves, keep(250).EQ.1
1055 & , .false., idummy, lidummy )
1056 ENDIF
1057 DEALLOCATE(wtemp)
1058#endif
1059 IF (keep(60).NE.0) THEN
1060 IF (keep(60)==1) THEN
1061 in = keep(20)
1062 ELSE
1063 in = keep(38)
1064 ENDIF
1065 DO WHILE (in.GT.0)
1066 in = fils(in)
1067 END DO
1068 ifson = -in
1069 IF (keep(60)==1) THEN
1070 in = keep(20)
1071 ELSE
1072 in = keep(38)
1073 ENDIF
1074 DO i=2,size_schur
1075 fils(in) = listvar_schur(i)
1076 in = fils(in)
1077 frere(in) = n+1
1078 ENDDO
1079 fils(in) = -ifson
1080 ENDIF
1081 CALL cmumps_ana_m(ikeep2(1),
1082 & ptrar(1,3), info(6),
1083 & info(5), keep(2), keep(50),
1084 & keep8(101), keep(108), keep(5),
1085 & keep(6), keep(226), keep(253))
1086 keep(59) = info(5)
1087 IF ( keep(53) .NE. 0 ) THEN
1088 CALL mumps_make1root( n, frere(1), fils(1), nfsiz(1),
1089 & keep(20) )
1090 END IF
1091 IF ( (keep(48) == 4 .AND. keep8(21).GT.0_8)
1092 & .OR.
1093 & (keep(48)==5 .AND. keep8(21) .GT. 0_8 )
1094 & .OR.
1095 & (keep(24).NE.0.AND.keep8(21).GT.0_8) ) THEN
1096 CALL cmumps_set_k821_surface(keep8(21), keep(2),
1097 & keep(48), keep(50), nslaves)
1098 END IF
1099 IF (keep(210).LT.0.OR.keep(210).GT.2) THEN
1100 keep(210)=0
1101 ENDIF
1102 IF (keep(210).EQ.0.AND.keep(201).GT.0) THEN
1103 keep(210)=1
1104 ENDIF
1105 IF (keep(210).EQ.0.AND.keep(201).EQ.0) THEN
1106 keep(210)=2
1107 ENDIF
1108 IF (keep(210).EQ.2) THEN
1109 keep8(79)=huge(keep8(79))
1110 ENDIF
1111 IF (keep(210).EQ.1.AND.keep8(79).LE.0_8) THEN
1112 keep8(79)=k79ref * int(nslaves,8)
1113 ENDIF
1114 IF ( (keep(79).EQ.0).OR.(keep(79).EQ.2).OR.
1115 & (keep(79).EQ.3).OR.(keep(79).EQ.5).OR.
1116 & (keep(79).EQ.6)
1117 & ) THEN
1118 IF (keep(210).EQ.1) THEN
1119 splitroot = .false.
1120 IF ( keep(62).GE.1) THEN
1121 iwl1(1) = -1
1122 IF (present(sizeofblocks)) THEN
1123 DO i= 1, n
1124 iwl1(i) = sizeofblocks(i)
1125 ENDDO
1126 ENDIF
1127 CALL cmumps_cutnodes(n, frere(1), fils(1), nfsiz(1),
1128 & iwl1(1), n, info(6),
1129 & nslaves, keep,keep8, splitroot,
1130 & mp, ldiag, info(1), info(2))
1131 IF (info(1).LT.0) GOTO 90
1132 IF (prok) THEN
1133 WRITE(mp,*) " Number of split nodes in pre-splitting=",
1134 & keep(61)
1135 ENDIF
1136 ENDIF
1137 ENDIF
1138 ENDIF
1139 splitroot = ((icntl(13).GT.0 .AND. nslaves.GT.icntl(13)) .OR.
1140 & icntl(13).EQ.-1 )
1141 IF (keep(53) .NE. 0) THEN
1142 splitroot = .true.
1143 ENDIF
1144 splitroot = (splitroot.AND.( (keep(60).EQ.0) ))
1145 IF (splitroot) THEN
1146 iwl1(1) = -1
1147 IF (present(sizeofblocks)) THEN
1148 DO i= 1, n
1149 iwl1(i) = sizeofblocks(i)
1150 ENDDO
1151 ENDIF
1152 CALL cmumps_cutnodes(n, frere(1), fils(1), nfsiz(1),
1153 & iwl1(1), n, info(6),
1154 & nslaves, keep,keep8, splitroot,
1155 & mp, ldiag, info(1), info(2))
1156 IF (info(1).LT.0) GOTO 90
1157 IF ( keep(53) .NE. 0 ) THEN
1158 CALL mumps_make1root( n, frere(1), fils(1), nfsiz(1),
1159 & keep(20) )
1160 ENDIF
1161 ENDIF
1162 IF (ldiag.GT.2 .AND. mp.GT.0) THEN
1163 k = min0(10,n)
1164 IF (ldiag.EQ.4) k = n
1165 IF (k.GT.0) WRITE (mp,99987) (nfsiz(i),i=1,k)
1166 IF (k.GT.0) WRITE (mp,99989) (fils(i),i=1,k)
1167 IF (k.GT.0) WRITE (mp,99988) (frere(i),i=1,k)
1168 ENDIF
1169 GO TO 90
1170 90 CONTINUE
1171 IF (info(1) .NE. 0) THEN
1172 IF ((lp.GT.0).AND.(icntl(4).GE.1))
1173 & WRITE (lp,99996) info(1), info(2)
1174 ENDIF
1175 IF (allocated(iwalloc)) DEALLOCATE(iwalloc)
1176 IF (allocated(iwl1)) DEALLOCATE(iwl1)
1177 IF (allocated(ipealloc)) DEALLOCATE(ipealloc)
1178 IF (allocated(ptrar)) DEALLOCATE(ptrar)
1179 IF (allocated(parent)) DEALLOCATE(parent)
1180 RETURN
118199999 FORMAT (/'Entering ordering phase with ...'/
1182 & ' N NNZ LIW INFO(1)'/,
1183 & 6x, i10, i11, i12, i10)
118499998 FORMAT ('Matrix entries: IRN() ICN()'/
1185 & (i12, i9, i12, i9, i12, i9))
118699909 FORMAT (/'Entering ordering phase with graph dimensions ...'/
1187 & ' |V| |E| INFO(1)'/,
1188 & 10x, i10, i13, i10)
118999997 FORMAT ('IKEEP1(.)=', 10i8/(12x, 10i8))
119099996 FORMAT
1191 & (/'** Error/warning return ** from Analysis * INFO(1:2)= ',
1192 & (i3, i16))
119399989 FORMAT ('FILS (.) =', 10i9/(11x, 10i9))
119499988 FORMAT ('FRERE(.) =', 10i9/(11x, 10i9))
119599987 FORMAT ('NFSIZ(.) =', 10i9/(11x, 10i9))
#define mumps_abort
Definition VE_Metis.h:25
subroutine mumps_wrap_ginp94(n, ipe, iw, liw8, perm, sizeofblocks, keep60, listvar_schur, size_schur, keep378, colcount, parent, porder, iwtmp1, iwtmp2, iwtmp3, iwtmp4, iwtmp5, info)
Definition ana_AMDMF.F:757
subroutine mumps_qamd(totel, compute_perm, iversion, thresh, ndense, n, iwlen, pe, pfree, len, iw, nv, elen, last, ncmpa, degree, head, next, w, parent)
subroutine mumps_ana_h(totel, compute_perm, n, iwlen, pe, pfree, len, iw, nv, elen, last, ncmpa, degree, head, next, w, parent)
subroutine mumps_symqamd(thresh, ndense, n, totel, iwlen, pe, pfree, len, iw, nv, elen, last, ncmpa, degree, head, next, w, perm, listvar_schur, size_schur, agg6, parent)
subroutine mumps_cst_amf(n, nbbuck, iwlen, pe, pfree, len, iw, nv, elen, last, ncmpa, degree, wf, next, w, head, constraint, theson, parent)
subroutine mumps_hamd(n, iwlen, pe, pfree, len, iw, nv, elen, last, ncmpa, degree, head, next, w, parent, listvar_schur, size_schur)
subroutine mumps_hamf4(norig, n, compute_perm, nbbuck, iwlen, pe, pfree, len, iw, nv, elen, last, ncmpa, degree, wf, next, w, head, parent)
subroutine mumps_set_ordering(n, keep, sym, nprocs, iord, nbqd, avgdens, prok, mp)
subroutine cmumps_ldlt_compress(n, nz, irn, icn, piv, ncmp, iw, lw, ipe, len, iq, flag, icmp, iwfr, ierror, keep, keep8, icntl, inplace64_graph_copy)
subroutine cmumps_set_constraints(n, piv, frere, fils, nfsiz, ikeep, ncst, keep, keep8, rowsca)
subroutine cmumps_get_perm_from_pe(n, pe, invperm, nfils, work)
subroutine cmumps_get_elim_tree(n, pe, nv, work)
subroutine cmumps_gnew_schur(na, n, nz, irn, icn, iw, lw, ipe, len, iq, flag, iwfr, nrorm, niorm, iflag, ierror, icntl, symmetry, sym, nbqd, avgdens, keep264, keep265, listvar_schur, size_schur, atoao, aotoa, inplace64_graph_copy)
subroutine cmumps_expand_permutation(n, ncmp, n11, n22, piv, invperm, perm)
subroutine cmumps_expand_perm_schur(na, ncmp, invperm, perm, listvar_schur, size_schur, aotoa)
subroutine cmumps_ana_gnew(n, nz, irn, icn, iw, lw, ipe, len, iq, flag, iwfr, nrorm, niorm, iflag, ierror, icntl, symmetry, sym, nbqd, avgdens, keep264, keep265, printstat, inplace64_graph_copy)
Definition cana_aux.F:3231
subroutine cmumps_ana_k(n, ipe, iw, lw, iwfr, ips, ipv, nv, flag, ncmpa, size_schur, parent)
Definition cana_aux.F:1931
subroutine cmumps_set_k821_surface(keep821, keep2, keep48, keep50, nslaves)
Definition cana_aux.F:3554
subroutine cmumps_ana_j(n, nz, irn, icn, perm, iw, lw, ipe, iq, flag, iwfr, iflag, ierror, mp)
Definition cana_aux.F:2032
subroutine cmumps_cutnodes(n, frere, fils, nfsiz, sizeofblocks, lsizeofblocks, nsteps, nslaves, keep, keep8, splitroot, mp, ldiag, info1, info2)
Definition cana_aux.F:2919
subroutine cmumps_ana_m(ne, nd, nsteps, maxfr, maxelim, k50, sizefac_tot, maxnpiv, k5, k6, panel_size, k253)
Definition cana_aux.F:2780
subroutine cmumps_ana_lnew(n, ipe, nv, ips, ne, na, nfsiz, node, nsteps, fils, frere, nd, nemin, subord, keep60, keep20, keep38, namalg, namalgmax, cumul, keep50, icntl13, keep37, keep197, nslaves, allow_amalg_tiny_nodes, blkon, sizeofblocks, lsizeofblocks)
Definition cana_aux.F:2412
#define min(a, b)
Definition macros.h:20
#define max(a, b)
Definition macros.h:21
subroutine mumps_secfin(t)
subroutine mumps_set_ierror(size8, ierror)
subroutine mumps_secdeb(t)
subroutine mumps_make1root(n, frere, fils, nfsiz, theroot)

◆ cmumps_ana_n_dist()

subroutine cmumps_ana_aux_m::cmumps_ana_n_dist ( type(cmumps_struc), intent(inout), target id,
integer(8), dimension(:), intent(out), target ptrar )

Definition at line 1197 of file cana_aux.F.

1198 USE cmumps_struc_def, ONLY : cmumps_struc
1199 IMPLICIT NONE
1200 include 'mpif.h'
1201 TYPE(CMUMPS_STRUC), INTENT(INOUT), TARGET :: id
1202 INTEGER(8), INTENT(OUT), TARGET :: PTRAR(:)
1203 INTEGER :: IERR, allocok
1204 INTEGER :: IOLD, JOLD, INEW, JNEW
1205 INTEGER(8) :: K, INZ
1206 INTEGER, POINTER :: IIRN(:), IJCN(:)
1207 INTEGER(8), POINTER :: IWORK1(:), IWORK2(:)
1208 LOGICAL :: IDO
1209 IF(id%KEEP(54) .EQ. 3) THEN
1210 iirn => id%IRN_loc
1211 ijcn => id%JCN_loc
1212 inz = id%KEEP8(29)
1213 iwork1 => ptrar(id%N+1:id%N+id%N)
1214 allocate(iwork2(id%N),stat=allocok)
1215 IF (allocok > 0 ) THEN
1216 id%INFO(1) = -7
1217 id%INFO(2) = id%N
1218 RETURN
1219 ENDIF
1220 ido = .true.
1221 ELSE
1222 iirn => id%IRN
1223 ijcn => id%JCN
1224 inz = id%KEEP8(28)
1225 iwork1 => ptrar(1:id%N)
1226 iwork2 => ptrar(id%N+1:id%N+id%N)
1227 ido = id%MYID .EQ. 0
1228 END IF
1229 DO 50 iold=1,id%N
1230 iwork1(iold) = 0_8
1231 iwork2(iold) = 0_8
1232 50 CONTINUE
1233 IF(ido) THEN
1234 DO 70 k=1_8,inz
1235 iold = iirn(k)
1236 jold = ijcn(k)
1237 IF ( (iold.GT.id%N).OR.(jold.GT.id%N).OR.(iold.LT.1)
1238 & .OR.(jold.LT.1) ) GOTO 70
1239 IF (iold.NE.jold) THEN
1240 inew = id%SYM_PERM(iold)
1241 jnew = id%SYM_PERM(jold)
1242 IF ( id%KEEP( 50 ) .EQ. 0 ) THEN
1243 IF (inew.LT.jnew) THEN
1244 iwork2(iold) = iwork2(iold) + 1_8
1245 ELSE
1246 iwork1(jold) = iwork1(jold) + 1_8
1247 ENDIF
1248 ELSE
1249 IF ( inew .LT. jnew ) THEN
1250 iwork1( iold ) = iwork1( iold ) + 1_8
1251 ELSE
1252 iwork1( jold ) = iwork1( jold ) + 1_8
1253 END IF
1254 ENDIF
1255 ENDIF
1256 70 CONTINUE
1257 END IF
1258 IF (id%KEEP(54) .EQ. 3) THEN
1259 CALL mpi_allreduce(iwork1(1), ptrar(1), id%N,
1260 & mpi_integer8, mpi_sum, id%COMM, ierr )
1261 CALL mpi_allreduce(iwork2(1), ptrar(id%N+1), id%N,
1262 & mpi_integer8, mpi_sum, id%COMM, ierr )
1263 deallocate(iwork2)
1264 ELSE
1265 CALL mpi_bcast( ptrar(1), 2*id%N, mpi_integer8,
1266 & 0, id%COMM, ierr )
1267 END IF
1268 RETURN
subroutine mpi_allreduce(sendbuf, recvbuf, cnt, datatype, operation, comm, ierr)
Definition mpi.f:103
subroutine mpi_bcast(buffer, cnt, datatype, root, comm, ierr)
Definition mpi.f:205
initmumps id

◆ cmumps_ana_o()

subroutine cmumps_ana_aux_m::cmumps_ana_o ( integer, intent(in) n,
integer(8), intent(in) nz,
integer, intent(inout) mtrans,
integer, dimension(:), intent(out) perm,
integer, dimension(3*n), target ikeepalloc,
integer, dimension(:), pointer idirn,
integer, dimension(:), pointer idjcn,
complex, dimension(:), pointer ida,
real, dimension(:), pointer idrowsca,
real, dimension(:), pointer idcolsca,
integer, dimension(n), target work2,
integer, dimension(500) keep,
integer, dimension(60), intent(in) icntl,
integer, dimension(80), intent(inout) info,
integer, dimension(80), intent(inout) infog )

Definition at line 1270 of file cana_aux.F.

1273 IMPLICIT NONE
1274 INTEGER, INTENT(IN) :: N
1275 INTEGER(8), INTENT(IN) :: NZ
1276 INTEGER, INTENT(OUT) :: PERM(:)
1277 INTEGER, POINTER, DIMENSION(:) :: idIRN, idJCN
1278 COMPLEX, POINTER, DIMENSION(:) :: idA
1279 REAL, POINTER, DIMENSION(:) :: idROWSCA, idCOLSCA
1280 INTEGER, TARGET :: IKEEPALLOC(3*N)
1281 INTEGER, INTENT(INOUT) :: MTRANS
1282 INTEGER :: KEEP(500)
1283 INTEGER, INTENT(IN) :: ICNTL(60)
1284 INTEGER, INTENT(INOUT) :: INFO(80)
1285 INTEGER, INTENT(INOUT) :: INFOG(80)
1286 INTEGER, TARGET :: WORK2(N)
1287 INTEGER :: allocok
1288 INTEGER, ALLOCATABLE, DIMENSION(:) :: IW
1289 REAL, ALLOCATABLE, DIMENSION(:) :: S2
1290 TARGET :: s2
1291 INTEGER ICNTL64(10), INFO64(10)
1292 INTEGER ICNTL_SYM_MWM(10),INFO_SYM_MWM(10)
1293 REAL CNTL64(10)
1294 INTEGER MPRINT,LP, MP
1295 INTEGER JPERM
1296 INTEGER NUMNZ, I, J, JPOS
1297 LOGICAL PROK, IDENT, DUPPLI
1298 INTEGER K50, KER_SIZE, NZER_DIAG, MTRANSLOC,RZ_DIAG
1299 INTEGER(8) :: LIWG
1300 INTEGER(8), DIMENSION(:), ALLOCATABLE :: IPE
1301 INTEGER(8), DIMENSION(:), ALLOCATABLE :: IPQ8
1302 INTEGER :: LSC
1303 INTEGER(8) :: NZTOT, NZREAL, IPIW, LIW, LIWMIN, NZsave,
1304 & K, KPOS, LDW, LDWMIN, IRNW, RSPOS, CSPOS,
1305 & LS2,J8, N8
1306 LOGICAL SCALINGLOC
1307 INTEGER,POINTER,DIMENSION(:) :: ZERODIAG
1308 INTEGER,POINTER,DIMENSION(:) :: STR_KER
1309 INTEGER,POINTER,DIMENSION(:) :: MARKED
1310 INTEGER,POINTER,DIMENSION(:) :: FLAG
1311 INTEGER,POINTER,DIMENSION(:) :: PIV_OUT
1312 REAL THEMIN, THEMAX, COLNORM,MAXDBL, ABSAK
1313 REAL ZERO,TWO,ONE
1314 parameter(zero = 0.0e0,two = 2.0e0,one = 1.0e0)
1315 n8 = int(n,8)
1316 mprint = icntl(3)
1317 lp = icntl(1)
1318 mp = icntl(2)
1319 prok = ((mprint.GT.0).AND.(icntl(4).GE.2))
1320 k50 = keep(50)
1321 scalingloc = .false.
1322 IF(keep(52) .EQ. -2) THEN
1323 IF(.not.associated(ida)) THEN
1324 ELSE
1325 scalingloc = .true.
1326 ENDIF
1327 ELSE IF(keep(52) .EQ. 77) THEN
1328 scalingloc = .true.
1329 IF( mtrans .NE. 5 .AND. mtrans .NE. 6
1330 & .AND. mtrans .NE. 7) THEN
1331 scalingloc = .false.
1332 ENDIF
1333 IF(.not.associated(ida)) THEN
1334 scalingloc = .false.
1335 IF (prok)
1336 & WRITE(mprint,*) 'Analysis: auto scaling OFF because ',
1337 & 'A not provided at analysis '
1338 ENDIF
1339 ENDIF
1340 IF ( (keep(50).EQ.2).AND.(icntl(8).NE.-2).AND.
1341 & (mtrans .EQ. 7 .OR. keep(95) .EQ. 0) ) THEN
1342 zerodiag => ikeepalloc(1:n)
1343 zerodiag = 0
1344 nzer_diag = n
1345 rz_diag = 0
1346 DO k=1,nz
1347 i = idirn(k)
1348 j = idjcn(k)
1349 IF (i.NE.j) cycle
1350 IF ( (j.LE.n).AND.(j.GE.1) ) THEN
1351 IF(zerodiag(i) .EQ. 0) THEN
1352 zerodiag(i) = 1
1353 IF(associated(ida)) THEN
1354 absak= abs(ida(k))
1355 IF(absak .EQ. real(0.0e0)) THEN
1356 rz_diag = rz_diag + 1
1357 ENDIF
1358 ENDIF
1359 nzer_diag = nzer_diag - 1
1360 ENDIF
1361 ENDIF
1362 ENDDO
1363 IF( (nzer_diag+rz_diag) .LT. (n/10) ) THEN
1364 mtrans = 0
1365 keep(95) =1
1366 GOTO 500
1367 ENDIF
1368 ENDIF
1369 IF(scalingloc) THEN
1370 IF (prok) WRITE(mprint,*)
1371 & 'Scaling will be computed during analysis'
1372 ENDIF
1373 IF( mtrans.NE.0 .AND. (.NOT.associated(ida)) ) mtrans=1
1374 mtransloc = mtrans
1375 IF (mtrans.LT.0 .OR. mtrans.GT.7) GO TO 500
1376 IF (k50 .EQ. 0) THEN
1377 IF(.NOT. scalingloc .AND. mtrans .EQ. 7) THEN
1378 GO TO 500
1379 ENDIF
1380 IF(scalingloc) THEN
1381 IF (mtransloc.NE.6) THEN
1382 mtransloc = 5
1383 ENDIF
1384 ENDIF
1385 ELSE
1386 IF (mtrans .EQ. 7) mtransloc = 5
1387 ENDIF
1388 IF(scalingloc .AND. mtransloc .NE. 5 .AND.
1389 & mtransloc .NE. 6 ) THEN
1390 IF (prok) WRITE(mprint,*)
1391 & 'WARNING scaling required: set MTRANS option to 5'
1392 mtransloc = 5
1393 ENDIF
1394 IF (n.EQ.1) THEN
1395 mtrans=0
1396 GO TO 500
1397 ENDIF
1398 IF(k50 .NE. 0) THEN
1399 nztot = 2_8*nz+n8
1400 ELSE
1401 nztot = nz
1402 ENDIF
1403 zerodiag => ikeepalloc(1:n)
1404 str_ker => ikeepalloc(n+1:2*n)
1405 CALL cmumps_mtransi(icntl64,cntl64)
1406 icntl64(1) = icntl(1)
1407 icntl64(2) = icntl(2)
1408 icntl64(3) = icntl(3)
1409 icntl64(4) = -1
1410 IF (icntl(4).EQ.3) icntl64(4) = 0
1411 IF (icntl(4).EQ.4) icntl64(4) = 1
1412 icntl64(5) = -1
1413 IF (prok) THEN
1414 WRITE(mprint,'(A,I3)')
1415 & 'Compute maximum matching (Maximum Transversal):',
1416 & mtransloc
1417 IF (mtransloc.EQ.1)
1418 & WRITE(mprint,'(A,I3)')' ... JOB =',mtransloc
1419 IF (mtransloc.EQ.2)
1420 & WRITE(mprint,'(A,I3,A)')
1421 & ' ... JOB =',mtransloc,': BOTTLENECK THESIS'
1422 IF (mtransloc.EQ.3)
1423 & WRITE(mprint,'(A,I3,A)')
1424 & ' ... JOB =',mtransloc,': BOTTLENECK SIMAX'
1425 IF (mtransloc.EQ.4)
1426 & WRITE(mprint,'(A,I3,A)')
1427 & ' ... JOB =',mtransloc,': MAXIMIZE SUM DIAGONAL'
1428 IF (mtransloc.EQ.5 .OR. mtransloc.EQ.6)
1429 & WRITE(mprint,'(A,I3,A)')
1430 & ' ... JOB =',mtransloc,
1431 & ': MAXIMIZE PRODUCT DIAGONAL AND SCALE'
1432 ENDIF
1433 infog(23) = mtransloc
1434 cntl64(2) = huge(cntl64(2))
1435 irnw = 1
1436 ipiw = irnw + nztot
1437 IF (mtransloc.EQ.1) liwmin = 5_8*n8
1438 IF (mtransloc.EQ.2) liwmin = 3_8*n8
1439 IF (mtransloc.EQ.3) liwmin = 10_8*n8 + nztot
1440 IF (mtransloc.EQ.4) liwmin = 2_8*n8
1441 IF (mtransloc.EQ.5) liwmin = 5_8*n8
1442 IF (mtransloc.EQ.6) liwmin = 5_8*n8 + nztot
1443 liw = liwmin
1444 liwg = liw + nztot
1445 ALLOCATE(iw(liwg), stat=allocok)
1446 IF (allocok .GT. 0 ) THEN
1447 GOTO 410
1448 ENDIF
1449 ALLOCATE( ipq8(n), ipe(n+1), stat = allocok )
1450 IF ( allocok .GT. 0 ) THEN
1451 info( 1 ) = -7
1452 info( 2 ) = (2*n+1)*keep(10)
1453 GOTO 500
1454 ENDIF
1455 IF (mtransloc.EQ.1) THEN
1456 ldwmin = n8+3_8
1457 ENDIF
1458 IF (mtransloc.EQ.2) ldwmin = max( n8+nztot , n8+3_8 )
1459 IF (mtransloc.EQ.3) ldwmin = max( nztot+1_8 , n8+3_8 )
1460 IF (mtransloc.EQ.4) ldwmin = 2_8 * n8 +
1461 & max( nztot , n8+3_8 )
1462 IF (mtransloc.EQ.5) ldwmin = 3_8*n8 + nztot
1463 IF (mtransloc.EQ.6) ldwmin = 4_8*n8 + nztot
1464 ldw = ldwmin
1465 ALLOCATE(s2(ldw), stat=allocok)
1466 IF (allocok .GT. 0 ) THEN
1467 GOTO 430
1468 ENDIF
1469 IF(mtransloc .NE. 1) ldw = ldw-nztot
1470 rspos = nztot
1471 cspos = rspos+n8
1472 nzreal = 0_8
1473 DO 5 j=1,n
1474 ipq8(j) = 0_8
1475 5 CONTINUE
1476 IF(k50 .EQ. 0) THEN
1477 DO 10 k=1,nz
1478 i = idirn(k)
1479 j = idjcn(k)
1480 IF ( (j.LE.n).AND.(j.GE.1).AND.
1481 & (i.LE.n).AND.(i.GE.1) ) THEN
1482 ipq8(j) = ipq8(j) + 1_8
1483 nzreal = nzreal + 1_8
1484 ENDIF
1485 10 CONTINUE
1486 ELSE
1487 zerodiag = 0
1488 nzer_diag = n
1489 rz_diag = 0
1490 DO k=1,nz
1491 i = idirn(k)
1492 j = idjcn(k)
1493 IF ( (j.LE.n).AND.(j.GE.1).AND.
1494 & (i.LE.n).AND.(i.GE.1) ) THEN
1495 ipq8(j) = ipq8(j) + 1_8
1496 nzreal = nzreal + 1_8
1497 IF(i .NE. j) THEN
1498 ipq8(i) = ipq8(i) + 1_8
1499 nzreal = nzreal + 1_8
1500 ELSE
1501 IF (zerodiag(i) .EQ. 0) THEN
1502 zerodiag(i) = 1
1503 IF(associated(ida)) THEN
1504 absak= abs(ida(k))
1505 IF(absak .EQ. real(0.0e0)) THEN
1506 rz_diag = rz_diag + 1
1507 ENDIF
1508 zerodiag(i) = exponent(absak)
1509 if ( zerodiag(i).EQ.0) zerodiag(i)=1
1510 ENDIF
1511 nzer_diag = nzer_diag - 1
1512 ELSE
1513 IF(associated(ida)) THEN
1514 absak= abs(ida(k))
1515 zerodiag(i) = zerodiag(i)+ exponent(absak)
1516 if ( zerodiag(i).EQ.0) zerodiag(i)=1
1517 ENDIF
1518 ENDIF
1519 ENDIF
1520 ENDIF
1521 ENDDO
1522 IF(mtransloc .GE. 4) THEN
1523 DO i =1, n
1524 IF(zerodiag(i) .EQ. 0) THEN
1525 ipq8(i) = ipq8(i) + 1_8
1526 nzreal = nzreal + 1_8
1527 ENDIF
1528 ENDDO
1529 ENDIF
1530 ENDIF
1531 ipe(1) = 1
1532 DO 20 j=1,n
1533 ipe(j+1) = ipe(j)+ipq8(j)
1534 20 CONTINUE
1535 DO 25 j=1, n
1536 ipq8(j ) = ipe(j)
1537 25 CONTINUE
1538 IF(k50 .EQ. 0) THEN
1539 IF (mtransloc.EQ.1) THEN
1540 DO k=1,nz
1541 i = idirn(k)
1542 j = idjcn(k)
1543 IF ( (j.LE.n).AND.(j.GE.1) .AND.
1544 & (i.LE.n).AND.(i.GE.1)) THEN
1545 kpos = ipq8(j)
1546 iw(irnw+kpos-1_8) = i
1547 ipq8(j) = ipq8(j) + 1_8
1548 ENDIF
1549 END DO
1550 ELSE
1551 IF ( .not.associated(ida)) THEN
1552 info(1) = -22
1553 info(2) = 4
1554 GOTO 500
1555 ENDIF
1556 DO k=1,nz
1557 i = idirn(k)
1558 j = idjcn(k)
1559 IF ( (j.LE.n).AND.(j.GE.1) .AND.
1560 & (i.LE.n).AND.(i.GE.1)) THEN
1561 kpos = ipq8(j)
1562 iw(irnw+kpos-1) = i
1563 s2(kpos) = abs(ida(k))
1564 ipq8(j) = ipq8(j) + 1_8
1565 ENDIF
1566 END DO
1567 ENDIF
1568 ELSE
1569 IF (mtransloc.EQ.1) THEN
1570 DO k=1,nz
1571 i = idirn(k)
1572 j = idjcn(k)
1573 IF ( (j.LE.n).AND.(j.GE.1) .AND.
1574 & (i.LE.n).AND.(i.GE.1)) THEN
1575 kpos = ipq8(j)
1576 iw(irnw+kpos-1) = i
1577 ipq8(j) = ipq8(j) + 1_8
1578 IF(i.NE.j) THEN
1579 kpos = ipq8(i)
1580 iw(irnw+kpos-1) = j
1581 ipq8(i) = ipq8(i) + 1_8
1582 ENDIF
1583 ENDIF
1584 ENDDO
1585 ELSE
1586 IF ( .not.associated(ida) ) THEN
1587 info(1) = -22
1588 info(2) = 4
1589 GOTO 500
1590 ENDIF
1591 themax = zero
1592 themin = huge(themin)
1593 DO k=1,nz
1594 i = idirn(k)
1595 j = idjcn(k)
1596 IF ( (j.LE.n).AND.(j.GE.1) .AND.
1597 & (i.LE.n).AND.(i.GE.1)) THEN
1598 kpos = ipq8(j)
1599 iw(irnw+kpos-1_8) = i
1600 s2(kpos) = abs(ida(k))
1601 ipq8(j) = ipq8(j) + 1_8
1602 IF(abs(ida(k)) .GT. themax) THEN
1603 themax = abs(ida(k))
1604 ELSE IF(abs(ida(k)) .LT. themin
1605 & .AND. abs(ida(k)).GT. zero) THEN
1606 themin = abs(ida(k))
1607 ENDIF
1608 IF(i.NE.j) THEN
1609 kpos = ipq8(i)
1610 iw(irnw+kpos-1) = j
1611 s2(kpos) = abs(ida(k))
1612 ipq8(i) = ipq8(i) + 1_8
1613 ENDIF
1614 ENDIF
1615 ENDDO
1616 DO i =1, n
1617 IF(zerodiag(i) .EQ. 0) THEN
1618 kpos = ipq8(i)
1619 iw(irnw+kpos-1) = i
1620 s2(kpos) = zero
1621 ipq8(i) = ipq8(i) + 1_8
1622 ENDIF
1623 ENDDO
1624 IF ( themax .NE. zero ) THEN
1625 cntl64(2) = (log(themax/themin))*(real(n))
1626 & - log(themin) + one
1627 ENDIF
1628 ENDIF
1629 ENDIF
1630 duppli = .false.
1631 nzsave = nzreal
1632 flag => ikeepalloc(2*n+1:3*n)
1633 IF(mtransloc.NE.1) THEN
1634 CALL cmumps_suppress_duppli_val(n,nzreal,ipe(1),iw(irnw),s2,
1635 & perm(1),ipq8(1))
1636 ELSE
1637 CALL cmumps_suppress_duppli_str(n,nzreal,ipe(1),iw(irnw),
1638 & perm(1))
1639 ENDIF
1640 IF(nzreal .NE. nzsave) duppli = .true.
1641 ls2 = nztot
1642 IF ( mtransloc .EQ. 1 ) THEN
1643 ls2 = 1_8
1644 ldw = 1_8
1645 ENDIF
1646 CALL cmumps_mtrans_driver(mtransloc ,n, n, nzreal,
1647 & ipe, iw(irnw), s2(1), ls2,
1648 & numnz, perm(1), liw, iw(ipiw), ldw, s2(ls2+1),
1649 & ipq8,
1650 & icntl64, cntl64, info64, info)
1651 IF (info(1).LT.0) THEN
1652 IF (lp.GT.0 .AND. icntl(4).GE.1)
1653 & WRITE(lp,'(A,I5)')
1654 & ' Not enough memory in MAXTRANS INFO(1)=',info(1)
1655 GOTO 500
1656 ENDIF
1657 IF (info64(1).LT.0) THEN
1658 IF (lp.GT.0 .AND. icntl(4).GE.1)
1659 & WRITE(lp,'(A,I5)')
1660 & ' INTERNAL ERROR in MAXTRANS INFO(1)=',info64(1)
1661 info(1) = -9964
1662 info(2) = info64(1)
1663 GO TO 500
1664 ENDIF
1665 IF (info64(1).GT.0) THEN
1666 IF (mp.GT.0 .AND. icntl(4).GE.2)
1667 & WRITE(mp,'(A,I5)')
1668 & ' WARNING in MAXTRANS INFO(1)=',info64(1)
1669 ENDIF
1670 ker_size = 0
1671 IF(k50 .EQ. 2) THEN
1672 DO i=1,n
1673 IF(zerodiag(i) .EQ. 0) THEN
1674 IF(perm(i) .EQ. i) THEN
1675 ker_size = ker_size + 1
1676 perm(i) = -i
1677 str_ker(ker_size) = i
1678 ENDIF
1679 ENDIF
1680 ENDDO
1681 ENDIF
1682 IF (numnz.LT.n) GO TO 400
1683 IF(k50 .EQ. 0) THEN
1684 ident = .true.
1685 IF (mtrans .EQ. 0 ) GOTO 102
1686 DO 80 j=1,n
1687 jperm = perm(j)
1688 iw(irnw+int(jperm-1,8)) = j
1689 IF (jperm.NE.j) ident = .false.
1690 80 CONTINUE
1691 IF(ident) THEN
1692 mtrans = 0
1693 ELSE
1694 IF(mtrans .EQ. 7) THEN
1695 mtrans = -9876543
1696 GOTO 102
1697 ENDIF
1698 IF (prok) WRITE(mprint,'(A)')
1699 & ' ... Apply column permutation'
1700 DO 100 k=1,nz
1701 j = idjcn(k)
1702 IF ((j.LE.0).OR.(j.GT.n)) GO TO 100
1703 idjcn(k) = iw(irnw+int(j-1,8))
1704 100 CONTINUE
1705 IF (mp.GT.0 .AND. icntl(4).GE.2)
1706 & WRITE(mp,'(/A)')
1707 & ' WARNING input matrix data modified'
1708 ENDIF
1709 102 CONTINUE
1710 IF (scalingloc) THEN
1711 IF ( associated(idcolsca))
1712 & DEALLOCATE( idcolsca )
1713 IF ( associated(idrowsca))
1714 & DEALLOCATE( idrowsca )
1715 ALLOCATE( idcolsca(n), stat=allocok)
1716 IF (allocok .GT.0) THEN
1717 info(1)=-5
1718 info(2)=n
1719 IF ((lp.GE.0).AND.(icntl(4).GE.1)) THEN
1720 WRITE (lp,'(/A)') '** Error in CMUMPS_ANA_O'
1721 WRITE (lp,'(A)')
1722 & '** Failure during allocation of COLSCA'
1723 GOTO 500
1724 ENDIF
1725 ENDIF
1726 ALLOCATE( idrowsca(n), stat=allocok)
1727 IF (allocok .GT.0) THEN
1728 info(1)=-5
1729 info(2)=n
1730 IF ((lp.GE.0).AND.(icntl(4).GE.1)) THEN
1731 WRITE (lp,'(/A)') '** Error in CMUMPS_ANA_O'
1732 WRITE (lp,'(A)')
1733 & '** Failure during allocation of ROWSCA'
1734 GOTO 500
1735 ENDIF
1736 ENDIF
1737 keep(52) = -2
1738 keep(74) = 1
1739 maxdbl = log(huge(maxdbl))
1740 DO j=1,n
1741 IF(s2(rspos+j) .GT. maxdbl) THEN
1742 s2(rspos+j) = zero
1743 ENDIF
1744 IF(s2(cspos+j) .GT. maxdbl) THEN
1745 s2(cspos+j)= zero
1746 ENDIF
1747 ENDDO
1748 DO 105 j=1,n
1749 j8 = int(j,8)
1750 idrowsca(j) = exp(s2(rspos+j8))
1751 IF(idrowsca(j) .EQ. zero) THEN
1752 idrowsca(j) = one
1753 ENDIF
1754 IF ( mtrans .EQ. -9876543 .OR. mtrans.EQ. 0 ) THEN
1755 idcolsca(j)= exp(s2(cspos+j8))
1756 IF(idcolsca(j) .EQ. zero) THEN
1757 idcolsca(j) = one
1758 ENDIF
1759 ELSE
1760 idcolsca(iw(irnw+j8-1_8))= exp(s2(cspos+j8))
1761 IF(idcolsca(iw(irnw+j8-1_8)) .EQ. zero) THEN
1762 idcolsca(iw(irnw+j8-1_8)) = one
1763 ENDIF
1764 ENDIF
1765 105 CONTINUE
1766 ENDIF
1767 ELSE
1768 ident = .false.
1769 IF(scalingloc) THEN
1770 IF ( associated(idcolsca)) DEALLOCATE( idcolsca )
1771 IF ( associated(idrowsca)) DEALLOCATE( idrowsca )
1772 ALLOCATE( idcolsca(n), stat=allocok)
1773 IF (allocok .GT.0) THEN
1774 info(1)=-5
1775 info(2)=n
1776 IF ((lp.GE.0).AND.(icntl(4).GE.1)) THEN
1777 WRITE (lp,'(/A)') '** Error in CMUMPS_ANA_O'
1778 WRITE (lp,'(A)')
1779 & '** Failure during allocation of COLSCA'
1780 GOTO 500
1781 ENDIF
1782 ENDIF
1783 ALLOCATE( idrowsca(n), stat=allocok)
1784 IF (allocok .GT.0) THEN
1785 info(1)=-5
1786 info(2)=n
1787 IF ((lp.GE.0).AND.(icntl(4).GE.1)) THEN
1788 WRITE (lp,'(/A)') '** Error in CMUMPS_ANA_O'
1789 WRITE (lp,'(A)')
1790 & '** Failure during allocation of ROWSCA'
1791 GOTO 500
1792 ENDIF
1793 ENDIF
1794 keep(52) = -2
1795 keep(74) = 1
1796 maxdbl = log(huge(maxdbl))
1797 DO j=1,n
1798 j8 = int(j,8)
1799 IF(s2(rspos+j8)+s2(cspos+j8) .GT. maxdbl) THEN
1800 s2(rspos+j8) = zero
1801 s2(cspos+j8)= zero
1802 ENDIF
1803 ENDDO
1804 DO j=1,n
1805 j8 = int(j,8)
1806 IF(perm(j) .GT. 0) THEN
1807 idrowsca(j) =
1808 & exp((s2(rspos+j8)+s2(cspos+j8))/two)
1809 IF(idrowsca(j) .EQ. zero) THEN
1810 idrowsca(j) = one
1811 ENDIF
1812 idcolsca(j)= idrowsca(j)
1813 ENDIF
1814 ENDDO
1815 DO jpos=1,ker_size
1816 i = str_ker(jpos)
1817 colnorm = zero
1818 DO k = ipe(i),ipe(i+1) - 1
1819 IF ( perm( iw( irnw+k-1_8) ) > 0 ) THEN
1820 colnorm = max(colnorm,s2(j))
1821 ENDIF
1822 ENDDO
1823 colnorm = exp(colnorm)
1824 idrowsca(i) = one / colnorm
1825 idcolsca(i) = idrowsca(i)
1826 ENDDO
1827 ENDIF
1828 IF(mtrans .EQ. 7 .OR. keep(95) .EQ. 0) THEN
1829 IF( (nzer_diag+rz_diag) .LT. (n/10)
1830 & .AND. keep(95) .EQ. 0) THEN
1831 mtrans = 0
1832 keep(95) = 1
1833 GOTO 390
1834 ELSE
1835 IF(keep(95) .EQ. 0) THEN
1836 IF(scalingloc) THEN
1837 keep(95) = 3
1838 ELSE
1839 keep(95) = 2
1840 ENDIF
1841 ENDIF
1842 IF(mtrans .EQ. 7) mtrans = 5
1843 ENDIF
1844 ENDIF
1845 IF(mtrans .EQ. 0) GOTO 390
1846 icntl_sym_mwm = 0
1847 info_sym_mwm = 0
1848 IF(mtrans .EQ. 5 .OR. mtrans .EQ. 6 .OR.
1849 & mtrans .EQ. 7) THEN
1850 icntl_sym_mwm(1) = 0
1851 icntl_sym_mwm(2) = 1
1852 ELSE IF(mtrans .EQ. 4) THEN
1853 icntl_sym_mwm(1) = 2
1854 icntl_sym_mwm(2) = 1
1855 ELSE
1856 icntl_sym_mwm(1) = 0
1857 icntl_sym_mwm(2) = 1
1858 ENDIF
1859 marked => ikeepalloc(n+1:2*n)
1860 flag => ikeepalloc(2*n+1:3*n)
1861 piv_out => work2(1:n)
1862 IF(mtransloc .LT. 4) THEN
1863 lsc = 1
1864 ELSE
1865 lsc = 2*n
1866 ENDIF
1867 CALL cmumps_sym_mwm(
1868 & n, nzreal, ipe, iw(irnw), s2(1),lsc, perm(1),
1869 & zerodiag(1),
1870 & icntl_sym_mwm, s2(lsc+1),marked(1),flag(1),
1871 & piv_out(1), info_sym_mwm)
1872 IF(info_sym_mwm(1) .NE. 0) THEN
1873 WRITE(*,*) '** Error in CMUMPS_ANA_O'
1874 RETURN
1875 ENDIF
1876 IF(info_sym_mwm(3) .EQ. n) THEN
1877 ident = .true.
1878 ELSEIF ( (icntl(12).EQ.0).AND.
1879 & ( (n-info_sym_mwm(4)-info_sym_mwm(3)) .GT. n/10 )
1880 & ) THEN
1881 ident = .true.
1882 keep(95) = 1
1883 ELSE
1884 DO i=1,n
1885 perm(i) = piv_out(i)
1886 ENDDO
1887 ENDIF
1888 keep(93) = info_sym_mwm(4)
1889 keep(94) = info_sym_mwm(3)
1890 IF (ident) mtrans=0
1891 ENDIF
1892 390 IF(mtrans .EQ. 0) THEN
1893 keep(95) = 1
1894 IF (prok) THEN
1895 WRITE (mprint,'(A)')
1896 & ' ... Column permutation not used'
1897 ENDIF
1898 ENDIF
1899 GO TO 500
1900 400 IF ((lp.GE.0).AND.(icntl(4).GE.1))
1901 & WRITE (lp,'(/A)') '** Error: Matrix is structurally singular'
1902 info(1) = -6
1903 info(2) = numnz
1904 GOTO 500
1905 410 IF ((lp.GE.0).AND.(icntl(4).GE.1)) THEN
1906 WRITE (lp,'(/A)') '** Error in CMUMPS_ANA_O'
1907 WRITE (lp,'(A,I14)')
1908 & '** Failure during allocation of INTEGER array of size ',
1909 & liwg
1910 ENDIF
1911 info(1) = -7
1912 CALL mumps_set_ierror(liwg,info(2))
1913 GOTO 500
1914 430 IF ((lp.GE.0).AND.(icntl(4).GE.1)) THEN
1915 WRITE (lp,'(/A)') '** Error in CMUMPS_ANA_O'
1916 WRITE (lp,'(A)') '** Failure during allocation of S2'
1917 ENDIF
1918 info(1) = -5
1919 CALL mumps_set_ierror(ldw,info(2))
1920 500 CONTINUE
1921 IF (allocated(iw)) DEALLOCATE(iw)
1922 IF (allocated(s2)) DEALLOCATE(s2)
1923 IF (allocated(ipe)) DEALLOCATE(ipe)
1924 IF (allocated(ipq8)) DEALLOCATE(ipq8)
1925 RETURN
subroutine cmumps_sym_mwm(n, ne, ip, irn, scaling, lsc, cperm, diag, icntl, weight, marked, flag, piv_out, info)
subroutine cmumps_suppress_duppli_str(n, nz, ip, irn, flag)
Definition cana_aux.F:4008
subroutine cmumps_mtrans_driver(job, m, n, ne, ip, irn, a, la, num, perm, liw, iw, ldw, dw, ipq8, icntl, cntl, info, infomumps)
Definition cana_aux.F:3612
subroutine cmumps_suppress_duppli_val(n, nz, ip, irn, a, flag, posi)
Definition cana_aux.F:3974
subroutine cmumps_mtransi(icntl, cntl)
Definition cana_mtrans.F:39