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

Go to the source code of this file.

Functions/Subroutines

subroutine zmumps_elt_distrib (n, nelt, na_elt8, comm, myid, slavef, ielptr_loc8, relptr_loc8, eltvar_loc, eltval_loc, lintarr, ldblarr, keep, keep8, maxelt_size, frtptr, frtelt, a, la, fils, id, root)
subroutine zmumps_elt_fill_buf (elnodes, elval, sizei, sizer, dest, nbuf, nbrecords, bufi, bufr, comm)
subroutine zmumps_maxelt_size (eltptr, nelt, maxelt_size)
subroutine zmumps_scale_element (n, sizei, sizer, eltvar, eltval, seltval, lseltval, rowsca, colsca, k50)

Function/Subroutine Documentation

◆ zmumps_elt_distrib()

subroutine zmumps_elt_distrib ( integer n,
integer nelt,
integer(8) na_elt8,
integer comm,
integer myid,
integer slavef,
integer(8), dimension( nelt + 1 ), intent(in) ielptr_loc8,
integer(8), dimension( nelt + 1 ), intent(inout) relptr_loc8,
integer, dimension( lintarr ) eltvar_loc,
complex(kind=8), dimension( ldblarr ) eltval_loc,
integer(8), intent(in) lintarr,
integer(8), intent(in) ldblarr,
integer, dimension(500) keep,
integer(8), dimension(150) keep8,
integer maxelt_size,
integer, dimension( n+1 ) frtptr,
integer, dimension( nelt ) frtelt,
complex(kind=8), dimension( la ) a,
integer(8), intent(in) la,
integer, dimension ( n ) fils,
type(zmumps_struc) id,
type(zmumps_root_struc) root )

Definition at line 14 of file zfac_distrib_ELT.F.

24 IMPLICIT NONE
25 INTEGER N, NELT
26 INTEGER(8) :: NA_ELT8
27 INTEGER COMM, MYID, SLAVEF, MAXELT_SIZE, MSGLEN
28 INTEGER(8), intent(IN) :: LA
29 INTEGER FRTPTR( N+1 )
30 INTEGER FRTELT( NELT ), FILS ( N )
31 INTEGER KEEP(500)
32 INTEGER(8) KEEP8(150)
33 INTEGER(8), INTENT(IN) :: IELPTR_LOC8( NELT + 1 )
34 INTEGER(8), INTENT(INOUT) :: RELPTR_LOC8( NELT + 1 )
35 INTEGER(8), INTENT(IN) :: LINTARR, LDBLARR
36 INTEGER ELTVAR_LOC( LINTARR )
37 COMPLEX(kind=8) ELTVAL_LOC( LDBLARR )
38 COMPLEX(kind=8) A( LA )
39 TYPE(ZMUMPS_STRUC) :: id
40 TYPE(ZMUMPS_ROOT_STRUC) :: root
41 INTEGER numroc
42 EXTERNAL numroc
43 include 'mpif.h'
44 include 'mumps_tags.h'
45 INTEGER :: IERR_MPI
46 INTEGER :: STATUS(MPI_STATUS_SIZE)
47 INTEGER MSGTAG
48 INTEGER allocok
49 INTEGER I, DEST, MAXELT_REAL_SIZE, MPG, IEL, SIZEI, SIZER
50 INTEGER NBRECORDS, NBUF
51 INTEGER(8) :: RECV_IELTPTR8
52 INTEGER(8) :: RECV_RELTPTR8
53 INTEGER INODE
54 INTEGER(8) :: IELTPTR8, RELTPTR8
55 LOGICAL FINI, PROKG, I_AM_SLAVE, EARLYT3ROOTINS
56 INTEGER(8) :: PTR_ROOT
57 INTEGER LOCAL_M, LOCAL_N, LP, IBEG, IGLOB, JGLOB
58 INTEGER ARROW_ROOT
59 INTEGER IELT, J, NB_REC, IREC
60 INTEGER(8) :: K8, IVALPTR8
61 INTEGER ILOCROOT, JLOCROOT, IPOSROOT, JPOSROOT, IPTR
62 INTEGER JCOL_GRID, IROW_GRID
63 INTEGER NBELROOT
64 INTEGER MASTER
65 parameter( master = 0 )
66 COMPLEX(kind=8) VAL
67 COMPLEX(kind=8) ZERO
68 parameter( zero = (0.0d0,0.0d0) )
69 INTEGER, DIMENSION( :, : ), ALLOCATABLE :: BUFI
70 COMPLEX(kind=8), DIMENSION( :, : ), ALLOCATABLE :: BUFR
71 COMPLEX(kind=8), DIMENSION( : ), ALLOCATABLE :: TEMP_ELT_R
72 INTEGER, DIMENSION( : ), ALLOCATABLE :: TEMP_ELT_I
73 INTEGER(8), DIMENSION( : ), ALLOCATABLE :: ELROOTPOS8
74 INTEGER, DIMENSION( : ), ALLOCATABLE, TARGET :: RG2LALLOC
75 INTEGER, DIMENSION( : ), POINTER :: RG2L
76 mpg = id%ICNTL(3)
77 lp = id%ICNTL(1)
78 i_am_slave = ( keep(46) .eq. 1 .or. myid .ne.master )
79 prokg = ( mpg > 0 .and. myid .eq. master )
80 prokg = (prokg.AND.(id%ICNTL(4).GE.2))
81 keep(49) = 0
82 arrow_root = 0
83 earlyt3rootins = keep(200) .EQ.0
84 & .OR. ( keep(200) .LT. 0 .AND. keep(400) .EQ. 0 )
85 IF ( myid .eq. master ) THEN
86 IF ( keep(46) .eq. 0 ) THEN
87 nbuf = slavef
88 ELSE
89 nbuf = slavef - 1
90 END IF
91 nbrecords = keep(39)
92 IF (na_elt8 < int(nbrecords,8)) THEN
93 nbrecords = int(na_elt8)
94 ENDIF
95 IF ( keep(50) .eq. 0 ) THEN
96 maxelt_real_size = maxelt_size * maxelt_size
97 ELSE
98 maxelt_real_size = maxelt_size * (maxelt_size+1)/2
99 END IF
100 IF ( maxelt_real_size .GT. keep(39) ) THEN
101 nbrecords = maxelt_real_size
102 IF ( mpg .GT. 0 ) THEN
103 WRITE(mpg,*)
104 & ' ** Warning : For element distrib NBRECORDS set to ',
105 & maxelt_real_size,' because one element is large'
106 END IF
107 END IF
108 ALLOCATE( bufi( 2*nbrecords+1, nbuf ), stat=allocok )
109 IF ( allocok .gt. 0 ) THEN
110 id%INFO(1) = -13
111 id%INFO(2) = 2*nbrecords + 1
112 GOTO 100
113 END IF
114 ALLOCATE( bufr( nbrecords+1, nbuf ), stat=allocok )
115 IF ( allocok .gt. 0 ) THEN
116 id%INFO(1) = -13
117 id%INFO(2) = nbrecords + 1
118 GOTO 100
119 END IF
120 IF ( keep(52) .ne. 0 ) THEN
121 ALLOCATE( temp_elt_r( maxelt_real_size ), stat =allocok )
122 IF ( allocok .gt. 0 ) THEN
123 id%INFO(1) = -13
124 id%INFO(2) = maxelt_real_size
125 GOTO 100
126 END IF
127 END IF
128 ALLOCATE( temp_elt_i( maxelt_size ), stat=allocok )
129 IF ( allocok .gt. 0 ) THEN
130 id%INFO(1) = -13
131 id%INFO(2) = maxelt_size
132 GOTO 100
133 END IF
134 IF ( keep(38) .ne. 0 ) THEN
135 nbelroot = frtptr(keep(38)+1)-frtptr(keep(38))
136 IF ( earlyt3rootins ) THEN
137 ALLOCATE( elrootpos8( max(nbelroot,1) ),
138 & stat = allocok )
139 IF ( allocok .gt. 0 ) THEN
140 id%INFO(1) = -13
141 id%INFO(2) = nbelroot
142 GOTO 100
143 END IF
144 ENDIF
145 IF (keep(46) .eq. 0 ) THEN
146 ALLOCATE( rg2lalloc( n ), stat = allocok )
147 IF ( allocok .gt. 0 ) THEN
148 id%INFO(1) = -13
149 id%INFO(2) = n
150 GOTO 100
151 END IF
152 inode = keep(38)
153 i = 1
154 DO WHILE ( inode .GT. 0 )
155 rg2lalloc( inode ) = i
156 inode = fils( inode )
157 i = i + 1
158 END DO
159 rg2l => rg2lalloc
160 ELSE
161 rg2l => root%RG2L_ROW
162 END IF
163 END IF
164 DO i = 1, nbuf
165 bufi( 1, i ) = 0
166 bufr( 1, i ) = zero
167 END DO
168 END IF
169 100 CONTINUE
170 CALL mumps_propinfo( id%ICNTL(1), id%INFO(1), comm, myid )
171 IF ( id%INFO(1) .LT. 0 ) RETURN
172 CALL mpi_bcast( nbrecords, 1, mpi_integer, master,
173 & comm, ierr_mpi )
174 recv_ieltptr8 = 1_8
175 recv_reltptr8 = 1_8
176 IF ( myid .eq. master ) THEN
177 nbelroot = 0
178 reltptr8 = 1_8
179 relptr_loc8(1) = 1
180 DO iel = 1, nelt
181 ieltptr8 = int(id%ELTPTR( iel ),8)
182 sizei = int(int(id%ELTPTR( iel + 1 ),8) - ieltptr8)
183 IF ( keep( 50 ) .eq. 0 ) THEN
184 sizer = sizei * sizei
185 ELSE
186 sizer = sizei * ( sizei + 1 ) / 2
187 END IF
188 dest = id%ELTPROC( iel )
189 IF ( dest .eq. -2 ) THEN
190 nbelroot = nbelroot + 1
191 frtelt( frtptr(keep(38)) + nbelroot - 1 ) = iel
192 elrootpos8( nbelroot ) = reltptr8
193 GOTO 200
194 END IF
195 IF ( dest .ge. 0 .and. keep(46) .eq. 0 ) dest = dest + 1
196 IF ( keep(52) .ne. 0 ) THEN
197 CALL zmumps_scale_element( n, sizei, sizer,
198 & id%ELTVAR( ieltptr8 ), id%A_ELT( reltptr8 ),
199 & temp_elt_r(1), maxelt_real_size,
200 & id%ROWSCA(1), id%COLSCA(1), keep(50) )
201 END IF
202 IF ( dest .eq. 0 .or. ( dest .eq. -1 .and. keep(46) .ne. 0 ) )
203 & THEN
204 eltvar_loc( recv_ieltptr8: recv_ieltptr8 + sizei - 1 )
205 & = id%ELTVAR( ieltptr8: ieltptr8 + sizei - 1 )
206 recv_ieltptr8 = recv_ieltptr8 + sizei
207 IF ( keep(52) .ne. 0 ) THEN
208 eltval_loc( recv_reltptr8: recv_reltptr8 + sizer - 1)
209 & = temp_elt_r( 1: sizer )
210 recv_reltptr8 = recv_reltptr8 + sizer
211 END IF
212 END IF
213 IF ( dest .NE. 0 .AND. dest. ne. -3 ) THEN
214 IF ( keep(52) .eq. 0 ) THEN
216 & id%ELTVAR(ieltptr8),
217 & id%A_ELT (reltptr8),
218 & sizei, sizer,
219 &
220 & dest, nbuf, nbrecords,
221 & bufi, bufr, comm )
222 ELSE
224 & id%ELTVAR(ieltptr8),
225 & temp_elt_r( 1 ),
226 & sizei, sizer,
227 &
228 & dest, nbuf, nbrecords,
229 & bufi, bufr, comm )
230 END IF
231 END IF
232 200 CONTINUE
233 reltptr8 = reltptr8 + sizer
234 IF ( keep(46) .eq. 0 .OR. keep(52) .eq. 0 ) THEN
235 relptr_loc8( iel + 1 ) = reltptr8
236 ELSE
237 relptr_loc8( iel + 1 ) = recv_reltptr8
238 ENDIF
239 END DO
240 IF ( keep(46) .eq. 0 .OR. keep(52) .eq. 0 ) THEN
241 keep8(26) = reltptr8 - 1_8
242 ELSE
243 keep8(26) = recv_reltptr8 - 1_8
244 ENDIF
245 IF ( reltptr8 - 1_8 .NE. na_elt8 ) THEN
246 WRITE(*,*) " ** Internal error in ZMUMPS_ELT_DISTRIB",
247 & reltptr8 - 1_8, na_elt8
248 CALL mumps_abort()
249 END IF
250 dest = -2
251 ieltptr8 = 1_8
252 reltptr8 = 1_8
253 sizei = 1
254 sizer = 1
256 & id%ELTVAR(ieltptr8),
257 & id%A_ELT (reltptr8),
258 & sizei, sizer,
259 &
260 & dest, nbuf, nbrecords,
261 & bufi, bufr, comm )
262 IF ( keep(52) .NE. 0 ) DEALLOCATE( temp_elt_r )
263 ELSE
264 fini = ( recv_ieltptr8 .eq. ielptr_loc8( nelt+1 )
265 & .and. recv_reltptr8 .eq. relptr_loc8( nelt+1 ) )
266 DO WHILE ( .not. fini )
267 CALL mpi_probe( master, mpi_any_tag,
268 & comm, status, ierr_mpi )
269 msgtag = status( mpi_tag )
270 SELECT CASE ( msgtag )
271 CASE( elt_int )
272 CALL mpi_get_count( status, mpi_integer,
273 & msglen, ierr_mpi )
274 CALL mpi_recv( eltvar_loc( recv_ieltptr8 ), msglen,
275 & mpi_integer, master, elt_int,
276 & comm, status, ierr_mpi )
277 recv_ieltptr8 = recv_ieltptr8 + msglen
278 CASE( elt_real )
279 CALL mpi_get_count( status, mpi_double_complex,
280 & msglen, ierr_mpi )
281 CALL mpi_recv( eltval_loc( recv_reltptr8 ), msglen,
282 & mpi_double_complex, master, elt_real,
283 & comm, status, ierr_mpi )
284 recv_reltptr8 = recv_reltptr8 + msglen
285 END SELECT
286 fini = ( recv_ieltptr8 .eq. ielptr_loc8( nelt+1 )
287 & .and. recv_reltptr8 .eq. relptr_loc8( nelt+1 ) )
288 END DO
289 END IF
290 IF ( keep(38) .NE. 0 .AND. earlyt3rootins ) THEN
291 IF ( i_am_slave .and. root%yes ) THEN
292 CALL zmumps_get_root_info(root,
293 & local_m, local_n, ptr_root, la)
294 CALL zmumps_set_root_to_zero(root, keep, a, la)
295 END IF
296 IF ( myid .NE. master ) THEN
297 ALLOCATE( bufi( nbrecords * 2 + 1, 1 ), stat = allocok )
298 IF ( allocok .GT. 0 ) THEN
299 id%INFO(1) = -13
300 id%INFO(2) = nbrecords * 2 + 1
301 GOTO 250
302 END IF
303 ALLOCATE( bufr( nbrecords, 1 ) , stat = allocok )
304 IF ( allocok .GT. 0 ) THEN
305 id%INFO(1) = -13
306 id%INFO(2) = nbrecords
307 END IF
308 END IF
309 250 CONTINUE
310 CALL mumps_propinfo( id%ICNTL(1), id%INFO(1), comm, myid )
311 IF ( id%INFO(1) .LT. 0 ) RETURN
312 IF ( myid .eq. master ) THEN
313 DO iptr = frtptr(keep(38)), frtptr(keep(38)+1) - 1
314 ielt = frtelt( iptr )
315 sizei = id%ELTPTR( ielt + 1 ) - id%ELTPTR( ielt )
316 DO i = 1, sizei
317 temp_elt_i( i ) = rg2l
318 & ( id%ELTVAR( id%ELTPTR(ielt) + i - 1 ) )
319 END DO
320 ivalptr8 = elrootpos8( iptr - frtptr(keep(38)) + 1 ) - 1
321 k8 = 1_8
322 DO j = 1, sizei
323 jglob = id%ELTVAR( id%ELTPTR( ielt ) + j - 1 )
324 IF ( keep(50).eq. 0 ) THEN
325 ibeg = 1
326 ELSE
327 ibeg = j
328 END IF
329 DO i = ibeg, sizei
330 iglob = id%ELTVAR( id%ELTPTR( ielt ) + i - 1 )
331 IF ( keep(52) .eq. 0 ) THEN
332 val = id%A_ELT( ivalptr8 + k8 )
333 ELSE
334 val = id%A_ELT( ivalptr8 + k8 ) *
335 & id%ROWSCA( iglob ) * id%COLSCA( jglob )
336 END IF
337 IF ( keep(50).eq.0 ) THEN
338 iposroot = temp_elt_i( i )
339 jposroot = temp_elt_i( j )
340 ELSE
341 IF ( temp_elt_i(i) .GT. temp_elt_i(j) ) THEN
342 iposroot = temp_elt_i(i)
343 jposroot = temp_elt_i(j)
344 ELSE
345 iposroot = temp_elt_i(j)
346 jposroot = temp_elt_i(i)
347 END IF
348 END IF
349 irow_grid = mod( ( iposroot - 1 )/root%MBLOCK,
350 & root%NPROW )
351 jcol_grid = mod( ( jposroot - 1 )/root%NBLOCK,
352 & root%NPCOL )
353 IF ( keep(46) .eq. 0 ) THEN
354 dest = irow_grid * root%NPCOL + jcol_grid + 1
355 ELSE
356 dest = irow_grid * root%NPCOL + jcol_grid
357 END IF
358 IF ( dest .eq. master ) THEN
359 ilocroot = root%MBLOCK * ( ( iposroot - 1 ) /
360 & ( root%MBLOCK * root%NPROW ) )
361 & + mod( iposroot - 1, root%MBLOCK ) + 1
362 jlocroot = root%NBLOCK * ( ( jposroot - 1 ) /
363 & ( root%NBLOCK * root%NPCOL ) )
364 & + mod( jposroot - 1, root%NBLOCK ) + 1
365 arrow_root = arrow_root + 1
366 IF (keep(60)==0) THEN
367 a( ptr_root
368 & + int(jlocroot - 1,8) * int(local_m,8)
369 & + int(ilocroot - 1,8) )
370 & = a( ptr_root
371 & + int(jlocroot - 1,8) * int(local_m,8)
372 & + int(ilocroot - 1,8) )
373 & + val
374 ELSE
375 root%SCHUR_POINTER( int(jlocroot-1,8)
376 & * int(root%SCHUR_LLD,8)
377 & + int(ilocroot,8) )
378 & = root%SCHUR_POINTER( int(jlocroot-1,8)
379 & * int(root%SCHUR_LLD,8)
380 & + int(ilocroot,8) )
381 & + val
382 ENDIF
383 ELSE
385 & iposroot, jposroot, val, dest, bufi, bufr, nbrecords,
386 & nbuf, lp, comm )
387 END IF
388 k8 = k8 + 1_8
389 END DO
390 END DO
391 END DO
393 & bufi, bufr, nbrecords,
394 & nbuf, lp, comm, keep(46) )
395 ELSE
396 fini = .false.
397 DO WHILE ( .not. fini )
398 CALL mpi_recv( bufi(1,1), 2*nbrecords+1,
399 & mpi_integer, master,
400 & arrowhead,
401 & comm, status, ierr_mpi )
402 nb_rec = bufi(1,1)
403 arrow_root = arrow_root + nb_rec
404 IF (nb_rec.LE.0) THEN
405 fini = .true.
406 nb_rec = -nb_rec
407 ENDIF
408 IF (nb_rec.EQ.0) EXIT
409 CALL mpi_recv( bufr(1,1), nbrecords, mpi_double_complex,
410 & master, arrowhead,
411 & comm, status, ierr_mpi )
412 DO irec = 1, nb_rec
413 iposroot = bufi( irec * 2, 1 )
414 jposroot = bufi( irec * 2 + 1, 1 )
415 val = bufr( irec, 1 )
416 ilocroot = root%MBLOCK * ( ( iposroot - 1 ) /
417 & ( root%MBLOCK * root%NPROW ) )
418 & + mod( iposroot - 1, root%MBLOCK ) + 1
419 jlocroot = root%NBLOCK * ( ( jposroot - 1 ) /
420 & ( root%NBLOCK * root%NPCOL ) )
421 & + mod( jposroot - 1, root%NBLOCK ) + 1
422 IF (keep(60).eq.0) THEN
423 a( ptr_root + int(jlocroot-1,8) * int(local_m,8)
424 & + int(ilocroot-1,8))
425 & = a( ptr_root + int(jlocroot-1,8) * int(local_m,8)
426 & + int(ilocroot-1,8))
427 & + val
428 ELSE
429 root%SCHUR_POINTER(int(jlocroot-1,8)
430 & * int(root%SCHUR_LLD,8)
431 & + int(ilocroot,8) )
432 & = root%SCHUR_POINTER( int(jlocroot - 1,8)
433 & * int(root%SCHUR_LLD,8)
434 & + int(ilocroot,8))
435 & + val
436 ENDIF
437 END DO
438 END DO
439 DEALLOCATE( bufi )
440 DEALLOCATE( bufr )
441 END IF
442 END IF
443 IF ( myid .eq. master ) THEN
444 DEALLOCATE( bufi )
445 DEALLOCATE( bufr )
446 IF (allocated(elrootpos8)) DEALLOCATE(elrootpos8)
447 IF (keep(38).ne.0) THEN
448 IF (keep(46) .eq. 0 ) THEN
449 DEALLOCATE(rg2lalloc)
450 ENDIF
451 ENDIF
452 DEALLOCATE( temp_elt_i )
453 END IF
454 keep(49) = arrow_root
455 RETURN
#define mumps_abort
Definition VE_Metis.h:25
subroutine mumps_propinfo(icntl, info, comm, id)
#define max(a, b)
Definition macros.h:21
subroutine mpi_recv(buf, cnt, datatype, source, tag, comm, status, ierr)
Definition mpi.f:461
subroutine mpi_get_count(status, datatype, cnt, ierr)
Definition mpi.f:296
subroutine mpi_probe(source, tag, comm, status, ierr)
Definition mpi.f:449
subroutine mpi_bcast(buffer, cnt, datatype, root, comm, ierr)
Definition mpi.f:205
integer function numroc(n, nb, iproc, isrcproc, nprocs)
Definition mpi.f:786
initmumps id
subroutine zmumps_set_root_to_zero(root, keep, a, la)
subroutine zmumps_arrow_finish_send_buf(bufi, bufr, nbrecords, nbufs, lp, comm, type_parall)
subroutine zmumps_get_root_info(root, local_m, local_n, ptr_root, la)
subroutine zmumps_arrow_fill_send_buf_elt(isend_shr, jsend_shr, val_shr, dest_shr, bufi, bufr, nbrecords, nbufs, lp, comm)
subroutine zmumps_scale_element(n, sizei, sizer, eltvar, eltval, seltval, lseltval, rowsca, colsca, k50)
subroutine zmumps_elt_fill_buf(elnodes, elval, sizei, sizer, dest, nbuf, nbrecords, bufi, bufr, comm)

◆ zmumps_elt_fill_buf()

subroutine zmumps_elt_fill_buf ( integer, dimension( sizei ) elnodes,
complex(kind=8), dimension( sizer ) elval,
integer sizei,
integer sizer,
integer dest,
integer nbuf,
integer nbrecords,
integer, dimension( 2*nbrecords + 1, nbuf ) bufi,
complex(kind=8), dimension( nbrecords + 1, nbuf ) bufr,
integer comm )

Definition at line 457 of file zfac_distrib_ELT.F.

460 IMPLICIT NONE
461 INTEGER SIZEI, SIZER, DEST, NBUF, NBRECORDS, COMM
462 INTEGER ELNODES( SIZEI ), BUFI( 2*NBRECORDS + 1, NBUF )
463 COMPLEX(kind=8) ELVAL( SIZER ), BUFR( NBRECORDS + 1, NBUF )
464 include 'mumps_tags.h'
465 include 'mpif.h'
466 INTEGER I, IBEG, IEND, IERR_MPI, NBRECR
467 INTEGER NBRECI
468 COMPLEX(kind=8) ZERO
469 parameter( zero = (0.0d0,0.0d0) )
470 IF ( dest .lt. 0 ) THEN
471 ibeg = 1
472 iend = nbuf
473 ELSE
474 ibeg = dest
475 iend = dest
476 END IF
477 DO i = ibeg, iend
478 nbreci = bufi(1,i)
479 IF ( nbreci .ne.0 .and.
480 & ( dest.eq.-2 .or.
481 & nbreci + sizei .GT. 2*nbrecords ) ) THEN
482 CALL mpi_send( bufi(2, i), nbreci, mpi_integer,
483 & i, elt_int, comm, ierr_mpi )
484 bufi(1,i) = 0
485 nbreci = 0
486 END IF
487 nbrecr = int(dble(bufr(1,i))+0.5d0)
488 IF ( nbrecr .ne.0 .and.
489 & ( dest.eq.-2 .or.
490 & nbrecr + sizer .GT. nbrecords ) ) THEN
491 CALL mpi_send( bufr(2, i), nbrecr, mpi_double_complex,
492 & i, elt_real, comm, ierr_mpi )
493 bufr(1,i) = zero
494 nbrecr = 0
495 END IF
496 IF ( dest .ne. -2 ) THEN
497 bufi( 2 + nbreci : 2 + nbreci + sizei - 1, i ) =
498 & elnodes( 1: sizei )
499 bufr( 2 + nbrecr : 2 + nbrecr + sizer - 1, i ) =
500 & elval( 1: sizer )
501 bufi(1,i) = nbreci + sizei
502 bufr(1,i) = cmplx( nbrecr + sizer, kind=kind(bufr) )
503 END IF
504 END DO
505 RETURN
float cmplx[2]
Definition pblas.h:136
subroutine mpi_send(buf, cnt, datatype, dest, tag, comm, ierr)
Definition mpi.f:480

◆ zmumps_maxelt_size()

subroutine zmumps_maxelt_size ( integer, dimension( nelt + 1 ) eltptr,
integer nelt,
integer maxelt_size )

Definition at line 507 of file zfac_distrib_ELT.F.

508 INTEGER NELT, MAXELT_SIZE
509 INTEGER ELTPTR( NELT + 1 )
510 INTEGER I, S
511 maxelt_size = 0
512 DO i = 1, nelt
513 s = eltptr( i + 1 ) - eltptr( i )
514 maxelt_size = max( s, maxelt_size )
515 END DO
516 RETURN

◆ zmumps_scale_element()

subroutine zmumps_scale_element ( integer n,
integer sizei,
integer sizer,
integer, dimension( sizei ) eltvar,
complex(kind=8), dimension( sizer ) eltval,
complex(kind=8), dimension( lseltval ) seltval,
integer lseltval,
double precision, dimension( n ) rowsca,
double precision, dimension( n ) colsca,
integer k50 )

Definition at line 518 of file zfac_distrib_ELT.F.

522 INTEGER N, SIZEI, SIZER, LSELTVAL, K50
523 INTEGER ELTVAR( SIZEI )
524 COMPLEX(kind=8) ELTVAL( SIZER )
525 COMPLEX(kind=8) SELTVAL( LSELTVAL )
526 DOUBLE PRECISION ROWSCA( N ), COLSCA( N )
527 INTEGER I, J, K
528 k = 1
529 IF ( k50 .eq. 0 ) THEN
530 DO j = 1, sizei
531 DO i = 1, sizei
532 seltval(k) = eltval(k) *
533 & rowsca(eltvar(i)) *
534 & colsca(eltvar(j))
535 k = k + 1
536 END DO
537 END DO
538 ELSE
539 DO j = 1, sizei
540 DO i = j, sizei
541 seltval(k) = eltval(k) *
542 & rowsca(eltvar(i)) *
543 & colsca(eltvar(j))
544 k = k + 1
545 END DO
546 END DO
547 END IF
548 RETURN