OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
spmd_tri25vox.F File Reference
#include "implicit_f.inc"
#include "spmd.inc"
#include "com01_c.inc"
#include "com04_c.inc"
#include "task_c.inc"
#include "timeri_c.inc"
#include "sms_c.inc"
#include "i25edge_c.inc"
#include "assert.inc"

Go to the source code of this file.

Functions/Subroutines

subroutine spmd_tri25vox (nsv, nsn, x, v, ms, bminmal, weight, stifn, nin, isendto, ircvfrom, iad_elem, fr_elem, nsnr, igap, gap_s, itab, kinet, ifq, inacti, nsnfiold, intth, ieles, areas, temp, num_imp, nodnx_sms, gap_s_l, ityp, irtlm, i24_time_s, i24_frfi, i24_pene_old, i24_stif_old, nbinflg, ilev, i24_icont_i, intfric, ipartfrics, itied, ivis2, if_adh, ledge, nedge, lndedge, stfm, nedge_local, gape, gap_e_l, stfe, edg_bisector, vtx_bisector, admsr, irect, ebinflg, mvoisin, iedge, icodt, iskew, ipartfric_e, e2s_nod_normal, istif_msdt, stifmsdt_s, stifmsdt_edg, ifsub_carea, intarean)

Function/Subroutine Documentation

◆ spmd_tri25vox()

subroutine spmd_tri25vox ( integer, dimension(*) nsv,
integer nsn,
x,
v,
ms,
bminmal,
integer, dimension(*) weight,
stifn,
integer nin,
integer, dimension(ninter+1,*) isendto,
integer, dimension(ninter+1,*) ircvfrom,
integer, dimension(2,*) iad_elem,
integer, dimension(*) fr_elem,
integer nsnr,
integer igap,
gap_s,
integer, dimension(*) itab,
integer, dimension(*) kinet,
integer ifq,
integer inacti,
integer, dimension(*) nsnfiold,
integer intth,
integer, dimension(*) ieles,
areas,
temp,
integer num_imp,
integer, dimension(*) nodnx_sms,
gap_s_l,
integer ityp,
integer, dimension(*) irtlm,
i24_time_s,
i24_frfi,
i24_pene_old,
i24_stif_old,
integer, dimension(*) nbinflg,
integer ilev,
integer, dimension(*) i24_icont_i,
integer intfric,
integer, dimension(*) ipartfrics,
integer itied,
integer ivis2,
integer, dimension(*) if_adh,
integer, dimension(lndedge,nedge) ledge,
integer nedge,
integer lndedge,
stfm,
integer, intent(in) nedge_local,
gape,
gap_e_l,
stfe,
real*4, dimension(3,4,*) edg_bisector,
real*4, dimension(3,2,*) vtx_bisector,
integer, dimension(4,*) admsr,
integer, dimension(4,*) irect,
integer, dimension(*), intent(in) ebinflg,
integer, dimension(4,*), intent(in) mvoisin,
integer, intent(in) iedge,
integer, dimension(*), intent(in) icodt,
integer, dimension(*), intent(in) iskew,
integer, dimension(*) ipartfric_e,
real*4, dimension(3,*) e2s_nod_normal,
integer, intent(in) istif_msdt,
dimension(nsn), intent(in) stifmsdt_s,
dimension(nedge), intent(in) stifmsdt_edg,
integer, intent(in) ifsub_carea,
dimension(numnod), intent(in) intarean )

Definition at line 39 of file spmd_tri25vox.F.

55C-----------------------------------------------
56C M o d u l e s
57C-----------------------------------------------
58 USE tri7box
59 USE tri25ebox
60 USE message_mod
61 USE spmd_mod
62C-----------------------------------------------
63C I m p l i c i t T y p e s
64C-----------------------------------------------
65#include "implicit_f.inc"
66C-----------------------------------------------
67C M e s s a g e P a s s i n g
68C-----------------------------------------------
69#include "spmd.inc"
70C-----------------------------------------------
71C C o m m o n B l o c k s
72C-----------------------------------------------
73#include "com01_c.inc"
74#include "com04_c.inc"
75#include "task_c.inc"
76#include "timeri_c.inc"
77#include "sms_c.inc"
78#include "i25edge_c.inc"
79#include "assert.inc"
80C-----------------------------------------------
81C D u m m y A r g u m e n t s
82C-----------------------------------------------
83 INTEGER NIN, NSN, IFQ, INACTI, IGAP,INTTH,NSNR,INTFRIC,
84 . ITIED, IVIS2,
85 . NSNFIOLD(*), NSV(*), WEIGHT(*),
86 . ISENDTO(NINTER+1,*), IRCVFROM(NINTER+1,*),
87 . IAD_ELEM(2,*), FR_ELEM(*), ITAB(*), KINET(*),
88 . IELES(*),NUM_IMP, NODNX_SMS(*),IRTLM(*),ITYP,
89 . NBINFLG(*),ILEV,I24_ICONT_I(*),IPARTFRICS(*),IF_ADH(*),
90 . IPARTFRIC_E(*)
91 INTEGER :: NEDGE, LNDEDGE, LEDGE(LNDEDGE,NEDGE)
92 INTEGER :: ADMSR(4,*),IRECT(4,*)
93 INTEGER, INTENT(IN) :: EBINFLG(*)
94 INTEGER, INTENT(IN) :: NEDGE_LOCAL
95 INTEGER, INTENT(IN) :: MVOISIN(4,*)
96 INTEGER, INTENT(IN) :: IEDGE
97 INTEGER, INTENT(IN) :: ICODT(*)
98 INTEGER, INTENT(IN) :: ISKEW(*)
99 INTEGER, INTENT(IN) :: ISTIF_MSDT, IFSUB_CAREA
100
101C INTEGER :: NSNFIEOLD(*)
102
103 my_real
104 . x(3,*), v(3,*), ms(*), bminmal(*), stifn(*), gap_s(*),
105 . areas(*),temp(*),gap_s_l(*),i24_time_s(*),i24_frfi(6,*),
106 . i24_pene_old(5,*),i24_stif_old(2,*),stfm(*),
107 . gape(*),
108 . gap_e_l(*),
109 . stfe(*)
110 real*4 edg_bisector(3,4,*),vtx_bisector(3,2,*),e2s_nod_normal(3,*)
111 my_real , INTENT(IN) :: stifmsdt_s(nsn), stifmsdt_edg(nedge)
112 my_real , INTENT(IN) :: intarean(numnod)
113
114C-----------------------------------------------
115C L o c a l V a r i a b l e s
116C-----------------------------------------------
117#ifdef MPI
118 INTEGER MSGTYP, I, NOD, LOC_PROC, P, IDEB,
119 . J, L, LEN, NB,
120 . IERROR,REQ_SB(NSPMD),
121 . REQ_RB(NSPMD),KK,NBIRECV,IRINDEXI(NSPMD),
122 . REQ_RD(NSPMD),REQ_SD(NSPMD),REQ_SD2(NSPMD),
123 . REQ_RC(NSPMD),REQ_SC(NSPMD),
124 . INDEXI,ISINDEXI(NSPMD),INDEX(NUMNOD),
125 . NBOX2(2,NSPMD),NBOX(2,NSPMD),
126 . NBX,NBY,NBZ,IX,IY,IZ,
127 . MSGOFF, MSGOFF2, MSGOFF3, MSGOFF4, MSGOFF5,MSGOFF6,
128 . MSGOFF7,
129 . RSIZ, ISIZ, L2, REQ_SD3(NSPMD),REQ_RD2(NSPMD),
130 . REQ_SD4(NSPMD),REQ_RD4(NSPMD),
131 . REQ_SD5(NSPMD),REQ_RD5(NSPMD),
132 . LEN2, RSHIFT, ISHIFT, ND, JDEB, Q, NBB,
133 . NB_EDGE, IDEB_EDGE
134
135 my_real:: xmins,ymins,zmins
136 my_real:: xmaxs,ymaxs,zmaxs
137 INTEGER :: N1,N2 ,NN1,NN2
138 INTEGER :: IX1,IX2,IY1,IY2,IZ1,IZ2
139 INTEGER :: IE,JE,I1,I2
140
141 my_real :: dx,dy,dz
142 my_real :: stf
143
144 DATA msgoff/6000/
145 DATA msgoff2/6001/
146 DATA msgoff3/6002/
147 DATA msgoff4/6003/
148 DATA msgoff5/6004/
149 DATA msgoff6/6006/
150 DATA msgoff7/6007/
151
152 my_real
153 . bminma(6,nspmd),
154 . xmaxb,ymaxb,zmaxb,xminb,yminb,zminb
155
156 TYPE(real_pointer), DIMENSION(NSPMD) :: RBUF
157 TYPE(int_pointer) , DIMENSION(NSPMD) :: IBUF
158 TYPE(int_pointer) , DIMENSION(NSPMD) :: IBUF_EDGE
159 TYPE(real_pointer), DIMENSION(NSPMD) :: RBUF_EDGE
160
161 INTEGER, DIMENSION(:), ALLOCATABLE :: ITAGNSNFI
162 INTEGER, DIMENSION(:), ALLOCATABLE :: INDEX_EDGE
163
164 INTEGER :: NBIRECV_NODE,NBIRECV_EDGE
165 INTEGER :: IAM,JAM,IM,M1,M2
166
167C-----------------------------------------------
168C S o u r c e L i n e s
169C-----------------------------------------------
170C
171C=======================================================================
172C tag of boxes containing facets
173C and creation of candidates
174C=======================================================================
175 loc_proc = ispmd + 1
176 nbx = lrvoxel
177 nby = lrvoxel
178 nbz = lrvoxel
179C
180C Sauvegarde valeur ancienne des nsn frontieres
181C
182 IF(inacti==5.OR.inacti==6.OR.inacti==7.OR.ifq>0
183 . .OR.num_imp>0.OR.itied/=0.OR.ityp==23.OR.ityp==24
184 . .OR.ityp==25) THEN
185 DO p = 1, nspmd
186 nsnfiold(p) = nsnfi(nin)%P(p)
187 IF(iedge > 0) THEN
188 nsnfieold(p) = nsnfie(nin)%P(p)
189 ENDIF
190 END DO
191 END IF
192C
193C minmax box for sorting from i7buce bminma
194C
195 nedge_remote = 0
196 DO p = 1, nspmd
197 nsnfi(nin)%P(p) = 0
198 IF(iedge /= 0) nsnfie(nin)%P(p) = 0
199 ENDDO
200
201 IF(ircvfrom(nin,loc_proc)==0.AND.
202 . isendto(nin,loc_proc)==0) RETURN
203
204 bminma(1,loc_proc) = bminmal(1)
205 bminma(2,loc_proc) = bminmal(2)
206 bminma(3,loc_proc) = bminmal(3)
207 bminma(4,loc_proc) = bminmal(4)
208 bminma(5,loc_proc) = bminmal(5)
209 bminma(6,loc_proc) = bminmal(6)
210C
211C Voxel shipment + min/max box
212C
213 IF(ircvfrom(nin,loc_proc)/=0) THEN
214 DO p = 1, nspmd
215 IF(isendto(nin,p)/=0) THEN
216 IF(p/=loc_proc) THEN
217 msgtyp = msgoff
218 CALL spmd_isend(
219 . crvoxel25(0,0,1,loc_proc),
220 . 2*(lrvoxel25+1)*(lrvoxel25+1),
221 . it_spmd(p),msgtyp,req_sc(p))
222 msgtyp = msgoff2
223 CALL spmd_isend(
224 . bminma(1,loc_proc),6 ,it_spmd(p),msgtyp,req_sb(p))
225 ENDIF
226 ENDIF
227 ENDDO
228 ENDIF
229C
230C Voxel reception + min-max boxes
231C
232 IF(isendto(nin,loc_proc)/=0) THEN
233 nbirecv=0
234 DO p = 1, nspmd
235 IF(ircvfrom(nin,p)/=0) THEN
236 IF(loc_proc/=p) THEN
237 nbirecv=nbirecv+1
238 irindexi(nbirecv)=p
239 msgtyp = msgoff
240 CALL spmd_irecv(
241 . crvoxel25(0,0,1,p),
242 . 2*(lrvoxel+1)*(lrvoxel+1),
243 . it_spmd(p),msgtyp,req_rc(nbirecv))
244 msgtyp = msgoff2
245 CALL spmd_irecv(
246 . bminma(1,p) ,6 ,it_spmd(p),msgtyp,
247 . req_rb(nbirecv))
248 ENDIF
249 ENDIF
250 ENDDO
251 ENDIF
252C
253C sending xrem
254C
255C computation of real and integer sending buffers sizes
256c general case
257 rsiz = 9
258 isiz = 6
259
260 IF(.true.) THEN
261 isiz = isiz + 2
262 ENDIF
263
264c specific cases
265c IGAP=1 or IGAP=2
266 IF(igap==1 .OR. igap==2)THEN
267 rsiz = rsiz + 1
268c IGAP=3
269 ELSEIF(igap==3)THEN
270 rsiz = rsiz + 2
271 ENDIF
272
273C thermic
274 IF(intth > 0 ) THEN
275 rsiz = rsiz + 2
276 isiz = isiz + 1
277 ENDIF
278
279C Interface Adhesion
280 IF(ityp==25.AND.ivis2==-1 ) THEN
281 IF(intth==0) rsiz = rsiz + 1 ! areas
282 isiz = isiz + 2 ! if_adh+ioldnsnfi
283 ENDIF
284
285C Friction
286 IF(intfric > 0 ) THEN
287 isiz = isiz + 1
288 ENDIF
289
290C---Stiffness based on mass and time step
291 IF(istif_msdt > 0) rsiz = rsiz + 1
292C---CAREA output
293 IF(ifsub_carea > 0) rsiz = rsiz + 1
294
295C -- IDTMINS==2
296 IF(idtmins == 2)THEN
297 isiz = isiz + 2
298C -- IDTMINS_INT /= 0
299 ELSEIF(idtmins_int/=0)THEN
300 isiz = isiz + 1
301 END IF
302
303c INT24
304 IF(ityp==24)THEN
305 rsiz = rsiz + 8
306 isiz = isiz + 3
307C-----for NBINFLG
308 IF (ilev==2) isiz = isiz + 1
309
310 ENDIF
311
312c INT25
313 IF(ityp==25)THEN
314 rsiz = rsiz + 3
315 isiz = isiz + 6
316C-----for NBINFLG
317 IF (ilev==2) isiz = isiz + 1
318 ENDIF
319 ideb = 1
320 req_sd4(1:nspmd) = mpi_request_null
321 req_sd5(1:nspmd) = mpi_request_null
322 req_rd(1:nspmd) = mpi_request_null
323 req_rd2(1:nspmd) = mpi_request_null
324 req_rd4(1:nspmd) = mpi_request_null
325 req_rd5(1:nspmd) = mpi_request_null
326
327
328
329 jdeb = 0
330 IF(ityp==25)THEN
331 ALLOCATE(itagnsnfi(numnod),stat=ierror)
332 itagnsnfi(1:numnod) = 0
333 ALLOCATE(index_edge(nedge),stat=ierror)
334 index_edge(1:nedge) = 0
335 END IF
336
337 IF(isendto(nin,loc_proc)/=0) THEN
338 DO kk = 1, nbirecv
339 CALL spmd_waitany(nbirecv,req_rb,indexi)
340 p=irindexi(indexi)
341 CALL spmd_wait(req_rc(indexi))
342C special processing on d.d. keep only the internal nodes
343 DO j = iad_elem(1,p), iad_elem(1,p+1)-1
344 nod = fr_elem(j)
345C weight < 0 temporarily to keep only the non-boundary nodes
346 weight(nod) = weight(nod)*(-1)
347 ENDDO
348C
349 l = ideb
350 nbox(2,p) = 0
351 nb = 0
352 xmaxb = bminma(1,p)
353 ymaxb = bminma(2,p)
354 zmaxb = bminma(3,p)
355 xminb = bminma(4,p)
356 yminb = bminma(5,p)
357 zminb = bminma(6,p)
358C ==================== Secnd nodes =============================
359 DO i=1,nsn
360 nod = nsv(i)
361 IF(weight(nod)==1)THEN
362 IF(stifn(i)>zero)THEN
363 IF(ityp==25.AND.irtlm(4*(i-1)+4)==p)THEN
364 nb = nb + 1
365 index(nb) = i
366 ELSEIF(itied/=0.AND.ityp==7.AND.candf_si(nin)%P(i)/=0) THEN
367 nb = nb + 1
368 index(nb) = i
369 ELSE
370 IF(x(1,nod) < xminb) cycle
371 IF(x(1,nod) > xmaxb) cycle
372 IF(x(2,nod) < yminb) cycle
373 IF(x(2,nod) > ymaxb) cycle
374 IF(x(3,nod) < zminb) cycle
375 IF(x(3,nod) > zmaxb) cycle
376 ix=int(nbx*(x(1,nod)-xminb)/(xmaxb-xminb))
377 IF(ix >= 0 .AND. ix <= nbx) THEN
378 iy=int(nby*(x(2,nod)-yminb)/(ymaxb-yminb))
379 IF(iy >= 0 .AND. iy <= nby) THEN
380 iz=int(nbz*(x(3,nod)-zminb)/(zmaxb-zminb))
381 IF(iz >= 0 .AND. iz <= nbz) THEN
382 IF(btest(crvoxel25(iy,iz,1,p),ix)) THEN
383 nb = nb + 1
384 index(nb) = i
385 ENDIF
386 ENDIF
387 ENDIF
388 ENDIF
389 ENDIF
390 ENDIF
391 ENDIF
392 ENDDO
393 nbox(1,p) = nb
394 DO j = iad_elem(1,p), iad_elem(1,p+1)-1
395 nod = fr_elem(j)
396C reset of weight > 0
397 weight(nod) = weight(nod)*(-1)
398 ENDDO
399
400C ==================== Secnd edges =============================
401 dx=xmaxb-xminb
402 dy=ymaxb-yminb
403 dz=zmaxb-zminb
404 nb_edge = 0
405! loop over edge that are local to ISPMD (ISPMD is main of the edge,
406! even if it is a boundary edge
407 IF(iedge /= 0) THEN
408 DO i=1,nedge_local
409 assert(ledge(9,i) == 1)
410 n1=ledge(5,i)
411 n2=ledge(6,i)
412 assert(n1 > 0)
413 assert(n2 > 0)
414 assert(n1 <= numnod)
415 assert(n2 <= numnod)
416
417 IF(ledge(1,i) > 0) THEN
418C First segment is local
419 stf = stfm(ledge(1,i))
420 ELSEIF (ledge(3,i) > 0) THEN
421C First segment is on the other side of the boundary
422 stf = one
423 IF(mvoisin(ledge(4,i),ledge(3,i)) == 0) stf = 0
424 ELSE !
425
426 ! ISPMD owns a boundary edge, but the local segment is deleted
427 stf = one
428 ENDIF
429 debug_e2e(ledge(8,i) == d_es,p-1)
430 debug_e2e(ledge(8,i) == d_es,stf)
431 debug_e2e(ledge(8,i) == d_es,ledge(7,i))
432
433
434 IF( stf > zero .AND. ledge(7,i) >= 0) THEN
435C
436C GAPE useless here (Redundant with BGAPEMAX on main side) !
437 xmins = min(x(1,n1),x(1,n2))!- GAPE(I)
438 ymins = min(x(2,n1),x(2,n2))!- GAPE(I)
439 zmins = min(x(3,n1),x(3,n2))!- GAPE(I)
440 xmaxs = max(x(1,n1),x(1,n2))!+ GAPE(I)
441 ymaxs = max(x(2,n1),x(2,n2))!+ GAPE(I)
442 zmaxs = max(x(3,n1),x(3,n2))!+ GAPE(I)
443
444 debug_e2e(ledge(8,i) == d_es, xmins)
445 debug_e2e(ledge(8,i) == d_es, ymins)
446 debug_e2e(ledge(8,i) == d_es, zmins)
447 debug_e2e(ledge(8,i) == d_es, xmaxs)
448 debug_e2e(ledge(8,i) == d_es, ymaxs)
449 debug_e2e(ledge(8,i) == d_es, zmaxs)
450
451 ix1=int(nbx*(xmins-xminb)/dx)
452 ix2=int(nbx*(xmaxs-xminb)/dx)
453
454 IF(ix2>=0.AND.ix1<=nbx)THEN
455 iy1=int(nby*(ymins-yminb)/dy)
456 iy2=int(nby*(ymaxs-yminb)/dy)
457
458 IF(iy2>=0.AND.iy1<=nby)THEN
459 iz1=int(nbz*(zmins-zminb)/dz)
460 iz2=int(nbz*(zmaxs-zminb)/dz)
461
462 IF(iz2>=0.AND.iz1<=nbz)THEN
463 ix1=max(ix1,0)
464 ix2=min(ix2,nbx)
465 iy1=max(iy1,0)
466 iy2=min(iy2,nbx)
467 iz1=max(iz1,0)
468 iz2=min(iz2,nbx)
469 DO ix=ix1,ix2
470 DO iy=iy1,iy2
471 DO iz=iz1,iz2
472 IF(btest(crvoxel25(iy,iz,1,p),ix)) THEN
473 nb_edge = nb_edge + 1
474 index_edge(nb_edge) = i
475 debug_e2e(ledge(8,i)==d_es,nb_edge)
476 GOTO 111 !next I
477 END IF
478 END DO
479 END DO
480 END DO
481 ENDIF
482 ENDIF
483 ENDIF
484 111 CONTINUE
485 ENDIF !
486 ENDDO
487 ENDIF ! IEDGE
488
489 nbox(2,p) = nb_edge
490C WRITE(6,*) ISPMD,"sends ",NB_EDGE,"to,",P-1
491 IF(ityp==25)THEN
492 jdeb = 0
493 DO q=1,p-1
494 jdeb = jdeb + nsnsi(nin)%P(q)
495 END DO
496 nbb = nsnsi(nin)%P(p)
497 DO j = 1, nbb
498 nd = nsvsi(nin)%P(jdeb+j)
499 nod= nsv(nd)
500 itagnsnfi(nod)=j
501 END DO
502 END IF
503C
504C Envoi taille msg
505C
506 msgtyp = msgoff3
507 CALL spmd_isend(nbox(1,p),2,it_spmd(p),msgtyp,req_sd(p))
508C
509C Alloc buffer
510C
511 IF( nb_edge > 0) THEN
512 ALLOCATE(ibuf_edge(p)%P(e_ibuf_size*nb_edge))
513 ALLOCATE(rbuf_edge(p)%P(e_rbuf_size*nb_edge))
514
515 l = 0
516 DO j=1,nb_edge
517 i = index_edge(j)
518 assert(i > 0)
519 assert(i <= nedge)
520 ibuf_edge(p)%p(e_global_id + l) = ledge(8,i)
521 ibuf_edge(p)%p(e_left_seg + l) = ledge(1,i)
522 ibuf_edge(p)%p(e_left_id + l) = ledge(2,i)
523 ibuf_edge(p)%p(e_right_seg + l) = ledge(3,i)
524 ibuf_edge(p)%p(e_right_id + l) = ledge(4,i)
525 ibuf_edge(p)%p(e_node1_id + l) = ledge(5,i)
526 ibuf_edge(p)%p(e_node2_id + l) = ledge(6,i)
527 ibuf_edge(p)%p(e_type + l) = ledge(7,i)
528! It is possible that one of the node is not sent
529 ibuf_edge(p)%p(e_node1_globid + l) = itab((ledge(5,i)))
530 ibuf_edge(p)%p(e_node2_globid + l) = itab((ledge(6,i)))
531 ibuf_edge(p)%p(e_local_id + l) = i
532 IF(ilev == 2) THEN
533 ibuf_edge(p)%p(e_ebinflg + l) = ebinflg(i)
534 ELSE
535 ibuf_edge(p)%p(e_ebinflg + l) = 0
536 ENDIF
537 iam= ledge(1,i)
538 jam= ledge(2,i)
539 m1 = ledge(5,i)
540 m2 = ledge(6,i)
541 im = ledge(10,i)
542 ibuf_edge(p)%p(e_im + l) = im
543 IF(idtmins /= 0) THEN
544 IF(idtmins/=2 .AND. idtmins_int == 0) THEN
545 ELSEIF(idtmins==2) THEN
546 ibuf_edge(p)%p(e_nodnx1 + l) = nodnx_sms(m1)
547 ibuf_edge(p)%p(e_nodams1 + l) = m1
548 ibuf_edge(p)%p(e_nodnx2 + l) = nodnx_sms(m2)
549 ibuf_edge(p)%p(e_nodams2 + l) = m2
550 ELSE ! IDTMINS_INT == 0
551 ibuf_edge(p)%p(e_nodnx1 + l) = 0
552 ibuf_edge(p)%p(e_nodams1 + l) = m1
553 ibuf_edge(p)%p(e_nodnx2 + l) = 0
554 ibuf_edge(p)%p(e_nodams2 + l) = m2
555 ENDIF
556 assert(nodnx_sms(m1) >=0)
557 assert(nodnx_sms(m2) >=0)
558 debug_e2e(nodnx_sms(m1) < 0,nodnx_sms(m1))
559 debug_e2e(nodnx_sms(m2) < 0,nodnx_sms(m2))
560 ENDIF ! IDTMINS /= 0
561 IF(intfric > 0) THEN
562 ibuf_edge(p)%p(e_ipartfric_e + l) = ipartfric_e(i)
563 ELSE
564 ibuf_edge(p)%p(e_ipartfric_e + l) = 0
565 ENDIF
566 l = l + e_ibuf_size
567 ENDDO
568
569 l = 0
570 DO j=1,nb_edge
571 i = index_edge(j)
572 rbuf_edge(p)%p(e_x1+ l) = x(1,(ledge(5,i)))
573 rbuf_edge(p)%p(e_y1+ l) = x(2,(ledge(5,i)))
574 rbuf_edge(p)%p(e_z1+ l) = x(3,(ledge(5,i)))
575 rbuf_edge(p)%p(e_x2+ l) = x(1,(ledge(6,i)))
576 rbuf_edge(p)%p(e_y2+ l) = x(2,(ledge(6,i)))
577 rbuf_edge(p)%p(e_z2+ l) = x(3,(ledge(6,i)))
578 rbuf_edge(p)%p(e_vx1+ l) = v(1,(ledge(5,i)))
579 rbuf_edge(p)%p(e_vy1+ l) = v(2,(ledge(5,i)))
580 rbuf_edge(p)%p(e_vz1+ l) = v(3,(ledge(5,i)))
581 rbuf_edge(p)%p(e_vx2+ l) = v(1,(ledge(6,i)))
582 rbuf_edge(p)%p(e_vy2+ l) = v(2,(ledge(6,i)))
583 rbuf_edge(p)%p(e_vz2+ l) = v(3,(ledge(6,i)))
584 rbuf_edge(p)%p(e_ms1+ l) = ms((ledge(5,i)))
585 rbuf_edge(p)%p(e_ms2+ l) = ms((ledge(6,i)))
586 rbuf_edge(p)%p(e_gap+ l) = gape(i)
587 IF(igap == 3) THEN
588 rbuf_edge(p)%p(e_gapl+ l) = gap_e_l(i)
589 ELSE
590 rbuf_edge(p)%p(e_gapl+ l) = 0
591 ENDIF
592 assert(not(isnan( rbuf_edge(p)%p(e_gapl+ l))))
593
594C RBUF_EDGE(P)%p(E_STIFE+ L) = STFM(LEDGE(1,I))
595 rbuf_edge(p)%p(e_stife+ l) = stfe(i)
596 assert(not(isnan(stfe(i))))
597
598C TO DO: single precision
599 l2 = e_edg_bis + l
600
601 ie = abs(ledge(1,i))
602 je = ledge(2,i)
603 iam = ledge(1,i)
604 jam = ledge(2,i)
605 m1 = ledge(5,i)
606 m2 = ledge(6,i)
607 im = ledge(10,i)
608 i1 = ledge(11,i)
609 i2 = ledge(12,i)
610 nn1 = admsr(je,ie)
611 nn2 = admsr(mod(je,4)+1,ie)
612
613
614 rbuf_edge(p)%p(l2:l2+2) = edg_bisector(1:3,je,ie)
615
616 l2 = e_vtx_bis + l
617 rbuf_edge(p)%p(l2:l2+2) = vtx_bisector(1:3,1,i1)
618
619 l2 = l2 + 3
620 rbuf_edge(p)%p(l2:l2+2) = vtx_bisector(1:3,2,i1)
621
622 l2 = l2 + 3
623 rbuf_edge(p)%p(l2:l2+2) = vtx_bisector(1:3,1,i2)
624
625 l2 = l2 + 3
626 rbuf_edge(p)%p(l2:l2+2) = vtx_bisector(1:3,2,i2)
627
628 l2 = l2 + 3
629 rbuf_edge(p)%p(l2:l2+2) = e2s_nod_normal(1:3,nn1)
630
631 l2 = l2 + 3
632 rbuf_edge(p)%p(l2:l2+2) = e2s_nod_normal(1:3,nn2)
633
634 IF(istif_msdt > 0) rbuf_edge(p)%p(e_stife_msdt_fi+ l) = stifmsdt_edg(i)
635
636 l = l + e_rbuf_size
637 ENDDO
638
639c DO J = 1, L-1
640c IF(ISNAN(RBUF_EDGE(P)%p(J))) THEN
641c WRITE(6,*) ISPMD,"NaN found",J,"/",L
642c ENDIF
643C ENDDO
644 ENDIF
645
646 IF (nb > 0) THEN
647 ALLOCATE(rbuf(p)%P(rsiz*nb),stat=ierror)
648 ALLOCATE(ibuf(p)%P(isiz*nb),stat=ierror)
649 l = 0
650 l2= 0
651
652#include "vectorize.inc"
653 DO j = 1, nb
654 i = index(j)
655 nod = nsv(i)
656 rbuf(p)%p(l+1) = x(1,nod)
657 rbuf(p)%p(l+2) = x(2,nod)
658 rbuf(p)%p(l+3) = x(3,nod)
659 rbuf(p)%p(l+4) = v(1,nod)
660 rbuf(p)%p(l+5) = v(2,nod)
661 rbuf(p)%p(l+6) = v(3,nod)
662 rbuf(p)%p(l+7) = ms(nod)
663 rbuf(p)%p(l+8) = stifn(i)
664 ibuf(p)%p(l2+1) = i
665 ibuf(p)%p(l2+2) = itab(nod)
666 ibuf(p)%p(l2+3) = kinet(nod)
667! save specifics IREM and XREM indexes for INT24 sorting
668 ibuf(p)%p(l2+4) = 0 !IGAPXREMP
669 ibuf(p)%p(l2+5) = 0 !I24XREMP
670 ibuf(p)%p(l2+6) = 0 !I24IREMP
671 l = l + rsiz
672 l2 = l2 + isiz
673 END DO
674
675c shift for real variables (prepare for next setting)
676 rshift = 9
677c shift for integer variables (prepare for next setting)
678 ishift = 7
679
680C symmetric plane
681 IF(.true. )THEN
682 l = 0
683#include "vectorize.inc"
684 DO j = 1, nb
685 i = index(j)
686 nod = nsv(i)
687 ibuf(p)%p(l+ishift+0)= icodt(nod)
688 ibuf(p)%p(l+ishift+1)= iskew(nod)
689 l = l + isiz
690 ENDDO
691 ishift = ishift + 2
692 ENDIF
693
694
695
696
697
698c specific cases
699c IGAP=1 or IGAP=2
700 IF(igap==1 .OR. igap==2.OR. igap==5)THEN
701 l = 0
702 igapxremp = rshift
703#include "vectorize.inc"
704 DO j = 1, nb
705 i = index(j)
706 rbuf(p)%p(l+rshift)= gap_s(i)
707 l = l + rsiz
708 ENDDO
709 rshift = rshift + 1
710
711c IGAP=3
712 ELSEIF(igap==3)THEN
713 l = 0
714 igapxremp = rshift
715#include "vectorize.inc"
716 DO j = 1, nb
717 i = index(j)
718 rbuf(p)%p(l+rshift) = gap_s(i)
719 rbuf(p)%p(l+rshift+1)= gap_s_l(i)
720 l = l + rsiz
721 END DO
722 rshift = rshift + 2
723 ENDIF
724
725C thermic
726 IF(intth>0)THEN
727 l = 0
728 l2 = 0
729#include "vectorize.inc"
730 DO j = 1, nb
731 i = index(j)
732 nod = nsv(i)
733 rbuf(p)%p(l+rshift) = temp(nod)
734 rbuf(p)%p(l+rshift+1) = areas(i)
735 ibuf(p)%p(l2+ishift) = ieles(i)
736 l = l + rsiz
737 l2 = l2 + isiz
738 END DO
739 rshift = rshift + 2
740 ishift = ishift + 1
741 ENDIF
742
743C Interface Adhesion
744 IF(ityp==25.AND.ivis2==-1)THEN
745 l = 0
746 l2 = 0
747#include "vectorize.inc"
748 DO j = 1, nb
749 i = index(j)
750 nod = nsv(i)
751 IF(intth==0) rbuf(p)%p(l+rshift) = areas(i)
752 ibuf(p)%p(l2+ishift) = if_adh(i)
753 ibuf(p)%p(l2+ishift+1)=itagnsnfi(nod)
754 IF(intth==0) l = l + rsiz
755 l2 = l2 + isiz
756 END DO
757 IF(intth==0) rshift = rshift + 1
758 ishift = ishift + 2
759 ENDIF
760
761C Friction
762 IF(intfric>0)THEN
763 l2 = 0
764#include "vectorize.inc"
765 DO j = 1, nb
766 i = index(j)
767 ibuf(p)%p(l2+ishift) = ipartfrics(i)
768 l2 = l2 + isiz
769 END DO
770 ishift = ishift + 1
771 ENDIF
772
773 IF(istif_msdt > 0) THEN
774 l = 0
775#include "vectorize.inc"
776 DO j = 1, nb
777 i = index(j)
778 rbuf(p)%p(l+rshift) =stifmsdt_s(i)
779 l = l + rsiz
780 END DO
781 rshift = rshift + 1
782 ENDIF
783
784
785 IF(ifsub_carea > 0) THEN
786 l = 0
787#include "vectorize.inc"
788 DO j = 1, nb
789 i = index(j)
790 nod = nsv(i)
791 rbuf(p)%p(l+rshift) =intarean(nod)
792 l = l + rsiz
793 END DO
794 rshift = rshift + 1
795 ENDIF
796
797C -- IDTMINS==2
798 IF(idtmins==2)THEN
799 l2 = 0
800#include "vectorize.inc"
801 DO j = 1, nb
802 i = index(j)
803 nod = nsv(i)
804 ibuf(p)%p(l2+ishift) = nodnx_sms(nod)
805 ibuf(p)%p(l2+ishift+1)= nod
806 l2 = l2 + isiz
807 END DO
808 ishift = ishift + 2
809
810C -- IDTMINS_INT /= 0
811 ELSEIF(idtmins_int/=0)THEN
812 l2 = 0
813#include "vectorize.inc"
814 DO j = 1, nb
815 i = index(j)
816 nod = nsv(i)
817 ibuf(p)%p(l2+ishift)= nod
818 l2 = l2 + isiz
819 END DO
820 ishift = ishift + 1
821 ENDIF
822
823c INT25
824 IF(ityp==25)THEN
825 l = 0
826 i24xremp = rshift
827#include "vectorize.inc"
828 DO j = 1, nb
829 i = index(j)
830 rbuf(p)%p(l+rshift) =i24_time_s(2*(i-1)+1)
831 rbuf(p)%p(l+rshift+1) =i24_time_s(2*(i-1)+2)
832 rbuf(p)%p(l+rshift+2) =i24_pene_old(5,i) ! used only at time=0
833 l = l + rsiz
834 END DO
835 rshift = rshift + 3
836
837 l2 = 0
838 i24iremp = ishift
839
840#include "vectorize.inc"
841 DO j = 1, nb
842 i = index(j)
843 nod = nsv(i)
844C IRTLM(3,NSN) in TYPE25 / IRTLM(3,-) useless here
845 ibuf(p)%p(l2+ishift) =irtlm(4*(i-1)+1)
846 ibuf(p)%p(l2+ishift+1)=irtlm(4*(i-1)+2)
847C
848C IRTLM(3,I) == local n of the impacted segment is shared but only valid on proc == IRTLM(4,I)
849 ibuf(p)%p(l2+ishift+2)=irtlm(4*(i-1)+3)
850 ibuf(p)%p(l2+ishift+3)=irtlm(4*(i-1)+4)
851 ibuf(p)%p(l2+ishift+4)=i24_icont_i(i)
852 ibuf(p)%p(l2+ishift+5)=itagnsnfi(nod)
853 l2 = l2 + isiz
854 END DO
855 ishift = ishift + 6
856C---pay attention in i25sto.F IREM(I24IREMP+4,N-NSN) is used,
857C----change the shift value when new table was added like IRTLM(3*(I-1)+2)
858 IF (ilev==2) THEN
859 l2 = 0
860#include "vectorize.inc"
861 DO j = 1, nb
862 i = index(j)
863 ibuf(p)%p(l2+ishift)=nbinflg(i)
864 l2 = l2 + isiz
865 END DO
866 END IF
867 ishift = ishift + 1
868
869 END IF !(ITYP==25)
870C
871 !save specifics IREM and XREM indexes for INT24 sorting
872 l2 = 0
873#include "vectorize.inc"
874 DO j = 1, nb
875 i = index(j)
876 nod = nsv(i)
877 !save specifics IREM and XREM indexes for INT24 sorting
878 ibuf(p)%p(l2+4) = igapxremp
879 ibuf(p)%p(l2+5) = i24xremp
880 ibuf(p)%p(l2+6) = i24iremp
881 l2 = l2 + isiz
882 END DO
883 ENDIF ! NB > 0
884
885 IF( nb > 0 ) THEN
886C WRITE(6,*) "Sends",NB,"nodes to",P-1
887 msgtyp = msgoff4
888 CALL spmd_isend(
889 1 rbuf(p)%P(1),nb*rsiz,it_spmd(p),msgtyp,
890 2 req_sd2(p))
891
892 msgtyp = msgoff5
893 CALL spmd_isend(
894 1 ibuf(p)%P(1),nb*isiz,it_spmd(p),msgtyp,
895 2 req_sd3(p))
896 ENDIF
897 IF(nb_edge > 0) THEN
898
899 msgtyp = msgoff6
900 CALL spmd_isend(
901 1 ibuf_edge(p)%P(1),e_ibuf_size*nb_edge ,it_spmd(p),msgtyp,
902 2 req_sd4(p))
903
904 msgtyp = msgoff7
905 CALL spmd_isend(
906 1 rbuf_edge(p)%P(1),e_rbuf_size*nb_edge ,it_spmd(p),msgtyp,
907 2 req_sd5(p))
908 ENDIF ! NB_EDGE > 0
909c ENDIF
910C
911C reset old tag for next P
912 IF(ityp==25)THEN
913C reset
914 nbb = nsnsi(nin)%P(p)
915 DO j = 1, nbb
916 nd = nsvsi(nin)%P(jdeb+j)
917 nod= nsv(nd)
918 itagnsnfi(nod)=0
919 END DO
920 END IF
921 ENDDO
922 ENDIF
923C
924 IF(ityp==25) THEN
925 DEALLOCATE(itagnsnfi)
926 DEALLOCATE(index_edge)
927 ENDIF
928C
929C reception of xrem data
930C
931 nedge_remote = 0
932 IF(ircvfrom(nin,loc_proc)/=0) THEN
933 nsnr = 0
934 l=0
935 DO p = 1, nspmd
936 nsnfi(nin)%P(p) = 0
937 IF(iedge /= 0) nsnfie(nin)%P(p) = 0
938 IF(isendto(nin,p)/=0) THEN
939 IF(loc_proc/=p) THEN
940 msgtyp = msgoff3
941 CALL spmd_recv(nbox2(1,p),2,it_spmd(p),msgtyp)
942 nsnfi(nin)%P(p) = nbox2(1,p)
943
944 IF(iedge /= 0) THEN
945 nedge_remote = nedge_remote + nbox2(2,p)
946 edge_fi(nin)%P(p) = nbox2(2,p)
947 nsnfie(nin)%P(p) = nbox2(2,p)
948 ELSE
949C EDGE_FI(NIN)%P(P) = 0
950C NSNFIE(NIN)%P(P) = 0
951 ENDIF
952
953 IF(nsnfi(nin)%P(p)> 0 .OR. nbox2(2,p) > 0) THEN
954 l=l+1
955 isindexi(l)=p
956 nsnr = nsnr + nsnfi(nin)%P(p)
957 ENDIF
958 ENDIF
959 ENDIF
960 ENDDO
961 nbirecv=l
962C
963C Allocate total size
964C
965
966 IF(nsnr > 0 .OR. nedge_remote > 0 ) THEN
967 ALLOCATE(xrem(rsiz,nsnr),stat=ierror)
968 IF(ierror/=0) THEN
969 CALL ancmsg(msgid=20,anmode=aninfo)
970 CALL arret(2)
971 ENDIF
972 ALLOCATE(irem(isiz,nsnr),stat=ierror)
973 IF(ierror/=0) THEN
974 CALL ancmsg(msgid=20,anmode=aninfo)
975 CALL arret(2)
976 ENDIF
977 IF(iedge /= 0) THEN
978 ALLOCATE(irem_edge(e_ibuf_size,nedge_remote),stat=ierror)
979 IF(ierror/=0) THEN
980 CALL ancmsg(msgid=20,anmode=aninfo)
981 CALL arret(2)
982 ENDIF
983 ALLOCATE(xrem_edge(e_rbuf_size,nedge_remote),stat=ierror)
984 IF(ierror/=0) THEN
985 CALL ancmsg(msgid=20,anmode=aninfo)
986 CALL arret(2)
987 ENDIF
988 ENDIF
989 ideb = 1
990 ideb_edge = 1
991 nbirecv_edge = 0
992 nbirecv_node = 0
993 DO l = 1, nbirecv
994 p = isindexi(l)
995 IF(nsnfi(nin)%P(p) > 0 ) THEN
996 len = nsnfi(nin)%P(p)*rsiz
997 msgtyp = msgoff4
998 nbirecv_node = nbirecv_node + 1
999 CALL spmd_irecv(
1000 1 xrem(1,ideb),len,it_spmd(p),
1001 2 msgtyp,req_rd(nbirecv_node))
1002
1003 len2 = nsnfi(nin)%P(p)*isiz
1004 msgtyp = msgoff5
1005 CALL spmd_irecv(
1006 1 irem(1,ideb),len2,it_spmd(p),
1007 2 msgtyp,req_rd2(nbirecv_node))
1008 ideb = ideb + nsnfi(nin)%P(p)
1009 ENDIF
1010
1011 IF(iedge /= 0) THEN
1012 IF(edge_fi(nin)%P(p) > 0 ) THEN
1013 msgtyp = msgoff6
1014 len2 = edge_fi(nin)%P(p)*e_ibuf_size
1015 nbirecv_edge = nbirecv_edge + 1
1016
1017 CALL spmd_irecv(
1018 1 irem_edge(1,ideb_edge),len2,it_spmd(p),
1019 2 msgtyp,req_rd4(nbirecv_edge))
1020
1021 msgtyp = msgoff7
1022 len2 = edge_fi(nin)%P(p)*e_rbuf_size
1023 CALL spmd_irecv(
1024 1 xrem_edge(1,ideb_edge),len2,it_spmd(p),
1025 2 msgtyp,req_rd5(nbirecv_edge))
1026 ideb_edge = ideb_edge + edge_fi(nin)%P(p)
1027 ENDIF
1028 ENDIF
1029 ENDDO
1030
1031
1032
1033 CALL spmd_waitall(nbirecv_node,req_rd )
1034 CALL spmd_waitall(nbirecv_node,req_rd2)
1035 CALL spmd_waitall(nbirecv_edge,req_rd4)
1036 CALL spmd_waitall(nbirecv_edge,req_rd5)
1037
1038 !set specifics IREM and XREM indexes for INT24 sorting
1039 IF(isiz > 5 .AND. nsnr > 0) THEN
1040 igapxremp = irem(4,1)
1041 i24xremp = irem(5,1)
1042 i24iremp = irem(6,1)
1043 ENDIF
1044 ENDIF
1045 ENDIF
1046C
1047 IF(ircvfrom(nin,loc_proc)/=0) THEN
1048 DO p = 1, nspmd
1049 IF(isendto(nin,p)/=0) THEN
1050 IF(p/=loc_proc) THEN
1051 CALL spmd_wait(req_sb(p))
1052 CALL spmd_wait(req_sc(p))
1053 ENDIF
1054 ENDIF
1055 ENDDO
1056 ENDIF
1057C
1058 IF(isendto(nin,loc_proc)/=0) THEN
1059 DO p = 1, nspmd
1060 IF(ircvfrom(nin,p)/=0) THEN
1061 IF(p/=loc_proc) THEN
1062 CALL spmd_wait(req_sd(p))
1063 IF(nbox(1,p) > 0) THEN
1064 CALL spmd_wait(req_sd2(p))
1065 DEALLOCATE(rbuf(p)%p)
1066 CALL spmd_wait(req_sd3(p))
1067 DEALLOCATE(ibuf(p)%p)
1068 ENDIF
1069 IF(nbox(2,p) > 0) THEN
1070 CALL spmd_wait(req_sd4(p))
1071 DEALLOCATE(ibuf_edge(p)%p)
1072 CALL spmd_wait(req_sd5(p))
1073 DEALLOCATE(rbuf_edge(p)%p)
1074 END IF
1075 ENDIF
1076 ENDIF
1077 ENDDO
1078 ENDIF
1079C
1080#endif
1081 RETURN
#define my_real
Definition cppsort.cpp:32
#define min(a, b)
Definition macros.h:20
#define max(a, b)
Definition macros.h:21
integer nedge_remote
Definition tri25ebox.F:75
integer, dimension(:,:,:,:), allocatable crvoxel25
Definition tri25ebox.F:72
integer, parameter lrvoxel25
Definition tri25ebox.F:71
integer, dimension(:,:), allocatable irem_edge
Definition tri25ebox.F:66
integer, dimension(:), allocatable nsnfieold
Definition tri25ebox.F:97
type(int_pointer), dimension(:), allocatable edge_fi
Definition tri25ebox.F:69
type(int_pointer), dimension(:), allocatable candf_si
Definition tri7box.F:560
type(int_pointer), dimension(:), allocatable nsvsi
Definition tri7box.F:485
integer i24iremp
Definition tri7box.F:423
type(int_pointer), dimension(:), allocatable nsnfie
Definition tri7box.F:440
type(int_pointer), dimension(:), allocatable nsnsi
Definition tri7box.F:491
integer i24xremp
Definition tri7box.F:423
integer igapxremp
Definition tri7box.F:423
integer lrvoxel
Definition tri7box.F:54
type(int_pointer), dimension(:), allocatable nsnfi
Definition tri7box.F:440
integer, dimension(:,:), allocatable irem
Definition tri7box.F:339
subroutine ancmsg(msgid, msgtype, anmode, i1, i2, i3, i4, i5, i6, i7, i8, i9, i10, i11, i12, i13, i14, i15, i16, i17, i18, i19, i20, r1, r2, r3, r4, r5, r6, r7, r8, r9, c1, c2, c3, c4, c5, c6, c7, c8, c9, prmode)
Definition message.F:895
subroutine arret(nn)
Definition arret.F:86