OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
spmd_i25front.F
Go to the documentation of this file.
1Copyright> OpenRadioss
2Copyright> Copyright (C) 1986-2025 Altair Engineering Inc.
3Copyright>
4Copyright> This program is free software: you can redistribute it and/or modify
5Copyright> it under the terms of the GNU Affero General Public License as published by
6Copyright> the Free Software Foundation, either version 3 of the License, or
7Copyright> (at your option) any later version.
8Copyright>
9Copyright> This program is distributed in the hope that it will be useful,
10Copyright> but WITHOUT ANY WARRANTY; without even the implied warranty of
11Copyright> MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12Copyright> GNU Affero General Public License for more details.
13Copyright>
14Copyright> You should have received a copy of the GNU Affero General Public License
15Copyright> along with this program. If not, see <https://www.gnu.org/licenses/>.
16Copyright>
17Copyright>
18Copyright> Commercial Alternative: Altair Radioss Software
19Copyright>
20Copyright> As an alternative to this open-source version, Altair also offers Altair Radioss
21Copyright> software under a commercial license. Contact Altair to discuss further if the
22Copyright> commercial version may interest you: https://www.altair.com/radioss/.
23#define TO1D(i,j,k,s1,s2) i+(j-1)*s1+(k-1)*s1*s2
24!||====================================================================
25!|| spmd_i25front_init ../engine/source/mpi/interfaces/spmd_i25front.F
26!||--- called by ------------------------------------------------------
27!|| resol ../engine/source/engine/resol.F
28!||--- calls -----------------------------------------------------
29!|| i25free_bound ../engine/source/interfaces/int25/i25free_bound.F
30!|| spmd_comm_split ../engine/source/mpi/generic/spmd_comm_split.F
31!||--- uses -----------------------------------------------------
32!|| intbufdef_mod ../common_source/modules/interfaces/intbufdef_mod.F90
33!|| message_mod ../engine/share/message_module/message_mod.F
34!|| spmd_comm_world_mod ../engine/source/mpi/spmd_comm_world.F90
35!|| tri7box ../engine/share/modules/tri7box.F
36!||====================================================================
37 SUBROUTINE spmd_i25front_init(ITAB,MAIN_PROC,INTBUF_TAB,IPARI)
38C============================================================================
39C M o d u l e s
40C-----------------------------------------------
41 USE tri7box
42 USE message_mod
43 USE intbufdef_mod
44C-----------------------------------------------
45C I m p l i c i t T y p e s
46C-----------------------------------------------
47 USE spmd_comm_world_mod, ONLY : spmd_comm_world
48#include "implicit_f.inc"
49C-----------------------------------------------
50C M e s s a g e P a s s i n g
51C-----------------------------------------------
52#include "spmd.inc"
53C-----------------------------------------------
54C C o m m o n B l o c k s
55C-----------------------------------------------
56#include "com01_c.inc"
57#include "com04_c.inc"
58#include "param_c.inc"
59#include "task_c.inc"
60#include "assert.inc"
61#include "macro.inc"
62C-----------------------------------------------
63C D u m m y A r g u m e n t s
64C-----------------------------------------------
65 TYPE(intbuf_struct_) INTBUF_TAB(*)
66 INTEGER, INTENT(IN) :: ITAB(NUMNOD),MAIN_PROC(NUMNOD)
67 INTEGER, INTENT(IN) :: IPARI(NPARI,NINTER)
68C-----------------------------------------------
69C L o c a l V a r i a b l e s
70C-----------------------------------------------
71 TYPE(int_pointer), DIMENSION(NINTER,NSPMD) :: SEND_BUF,RECV_BUF
72 INTEGER, DIMENSION(NINTER,NSPMD) :: SIZ_SEND_BUF,SIZ_RECV_BUF
73 INTEGER :: I,J,K,L,II,NIN,ITY,P,UID
74 INTEGER :: IERROR
75 INTEGER :: COLOR,KEY,NEDGE
76 INTEGER :: IED,IE,JE,WGT,NRTM,NSN
77 INTEGER, DIMENSION(:,:), ALLOCATABLE :: TAB
78#ifdef MPI
79 INTEGER :: STATUS(MPI_STATUS_SIZE), MSGTYP, REQ_S(NSPMD),REQ_R(NSPMD)
80#else
81 INTEGER,PARAMETER :: MPI_COMM_NULL = 0
82#endif
83
84 ninter25e = 0
85 DO nin = 1,ninter
86 intbuf_tab(nin)%MPI_COMM = mpi_comm_null
87 intbuf_tab(nin)%RANK = -1
88 intbuf_tab(nin)%NSPMD = 0
89 intbuf_tab(nin)%NB_INTERNAL_EDGES = 0
90 intbuf_tab(nin)%NB_BOUNDARY_EDGES_LOCAL = 0
91 intbuf_tab(nin)%NB_BOUNDARY_EDGES_REMOTE = 0
92 ity = ipari(7,nin)
93 IF( ity == 25 ) THEN
94C Compute free bound
95 nrtm = ipari(4,nin)
96 ALLOCATE(intbuf_tab(nin)%FREE_IRECT_ID(nrtm))
97 CALL i25free_bound(
98 . nrtm,
99 . intbuf_tab(nin)%MVOISIN,
100 . intbuf_tab(nin)%IRECTM,
101 . intbuf_tab(nin)%STFM,
102 . intbuf_tab(nin)%NRTM_FREE,
103 . intbuf_tab(nin)%FREE_IRECT_ID)
104
105C Split communicator
106
107 key = ispmd
108 color = 0
109 nsn = ipari(macro_nsn,nin)
110 nedge = ipari(macro_nedge,nin)
111 IF(ipari(macro_iedge,nin) > 0) ninter25e = ninter25e + 1
112 IF(nedge > 0) THEN
113 color = 1
114
115 ALLOCATE(tab(nedge,3))
116
117 ! Gather iformations
118 DO ie = 1,nedge
119 tab(ie,1) = intbuf_tab(nin)%LEDGE( (ie-1)*nledge + 1 ) !IE
120 tab(ie,2) = intbuf_tab(nin)%LEDGE( (ie-1)*nledge + 3 ) !JE
121 tab(ie,3) = intbuf_tab(nin)%LEDGE( (ie-1)*nledge + 9 ) !weight
122 ENDDO
123
124 ied = count(tab(1:nedge,1) >= 0 .AND. tab(1:nedge,2) >= 0)
125 intbuf_tab(nin)%NB_INTERNAL_EDGES = ied
126
127
128 ied = count(.NOT. (tab(1:nedge,1) >= 0 .AND. tab(1:nedge,2) >= 0)
129 . .AND. tab(1:nedge,3) == 1)
130 intbuf_tab(nin)%NB_BOUNDARY_EDGES_LOCAL = ied
131
132 intbuf_tab(nin)%NB_BOUNDARY_EDGES_REMOTE = nedge
133 . - intbuf_tab(nin)%NB_BOUNDARY_EDGES_LOCAL
134 . - intbuf_tab(nin)%NB_INTERNAL_EDGES
135
136 DEALLOCATE(tab)
137 ENDIF
138#ifdef MPI
139 CALL spmd_comm_split(color,key,intbuf_tab(nin)%MPI_COMM,
140 . intbuf_tab(nin)%RANK, intbuf_tab(nin)%NSPMD)
141#endif
142 ENDIF
143 ENDDO
144
145
146
147#ifdef MPI
148 siz_recv_buf(1:ninter,1:nspmd) = 0
149 siz_send_buf(1:ninter,1:nspmd) = 0
150C sizes of send buffers
151 DO nin = 1,ninter
152 ity = ipari(macro_ity,nin)
153 nsn = ipari(macro_nsn,nin)
154 IF( ity == 25) THEN
155 DO i = 1,nsn
156 p = main_proc(intbuf_tab(nin)%NSV(i))
157 siz_send_buf(nin,p) = siz_send_buf(nin,p) + 1
158 ENDDO
159 ENDIF
160 ENDDO
161
162 CALL mpi_alltoall(siz_send_buf, ninter, mpi_integer,
163 . siz_recv_buf, ninter, mpi_integer,
164 . spmd_comm_world,ierror)
165
166C allocation of buffers
167 DO nin = 1,ninter
168 ity=ipari(macro_ity,nin)
169 nsn = ipari(macro_nsn,nin)
170 IF( ity == 25 ) THEN
171 DO p=1,nspmd
172c WRITE(6,*) ISPMD+1,P, SIZ_SEND_BUF(NIN,P),SIZ_RECV_BUF(NIN,P)
173 ALLOCATE(send_buf(nin,p)%P(siz_send_buf(nin,p)))
174 ALLOCATE(recv_buf(nin,p)%P(siz_recv_buf(nin,p)))
175 ENDDO
176 ENDIF
177 ENDDO
178
179C sizes of send buffers
180 siz_send_buf(1:ninter,1:nspmd) = 0
181 DO nin = 1,ninter
182 ity=ipari(macro_ity,nin)
183 nsn = ipari(macro_nsn,nin)
184 IF( ity == 25) THEN
185 DO i = 1,nsn
186 p = main_proc(intbuf_tab(nin)%NSV(i))
187 siz_send_buf(nin,p) = siz_send_buf(nin,p) + 1
188 ii = siz_send_buf(nin,p)
189 IF(ispmd+1 == p ) THEN
190 send_buf(nin,p)%P(ii) = i
191 ELSE
192 send_buf(nin,p)%P(ii) = itab(intbuf_tab(nin)%NSV(i))
193 ENDIF
194 ENDDO
195 ENDIF
196 ENDDO
197
198 msgtyp = 1000
199 DO nin = 1, ninter
200 DO p = 1, nspmd
201 IF(p /= ispmd + 1 .AND. siz_send_buf(nin,p) > 0) THEN
202 k = siz_send_buf(nin,p)
203c WRITE(6,*) "SEND ->",P,SEND_BUF(NIN,P)%P(1:k)
204 CALL mpi_isend(
205 . send_buf(nin,p)%P(1),k,mpi_integer,it_spmd(p),msgtyp,
206 . spmd_comm_world,req_s(p),ierror)
207 ENDIF
208 IF(p /= ispmd + 1 .AND. siz_recv_buf(nin,p) > 0) THEN
209 k = siz_recv_buf(nin,p)
210 CALL mpi_irecv(
211 . recv_buf(nin,p)%P(1),k,mpi_integer,it_spmd(p),msgtyp,
212 . spmd_comm_world,req_r(p),ierror)
213 ENDIF
214 ENDDO
215c CALL FLUSH(6)
216c CALL MPI_BARRIER(SPMD_COMM_WORLD,IERR)
217
218 DO p = 1,nspmd
219 IF(p /= ispmd + 1 .AND. siz_send_buf(nin,p) > 0) THEN
220 CALL mpi_wait(req_s(p),status,ierror)
221 ENDIF
222 IF(p /= ispmd + 1 .AND. siz_recv_buf(nin,p) > 0) THEN
223 CALL mpi_wait(req_r(p),status,ierror)
224 ENDIF
225 ENDDO
226C ----------------------------------------------------------------
227C RESOLVE UID to LOCAL NODE
228 DO p = 1, nspmd
229c WRITE(6,*) __FILE__,__LINE__, P,SIZ_RECV_BUF(NIN,P)
230 IF(p /= ispmd +1) THEN
231 DO l = 1, siz_recv_buf(nin,p)
232 k = 1
233 uid = recv_buf(nin,p)%P(l)
234c WRITE(6,*) "looking for UID",UID
235 nsn = ipari(macro_nsn,nin)
236 DO WHILE (k<nsn .AND. uid /= itab(intbuf_tab(nin)%NSV(k)))
237 k = k + 1
238 ENDDO
239 IF(itab(intbuf_tab(nin)%NSV(k))==uid) THEN
240c WRITE(6,*) "found at ",K
241 recv_buf(nin,p)%P(l) = k
242 ELSE
243 assert(.false.)
244 ENDIF
245 ENDDO
246 ENDIF
247 ENDDO
248
249 DO p = 1, nspmd
250 IF(p /= ispmd + 1 .AND. siz_recv_buf(nin,p) > 0) THEN
251 k = siz_recv_buf(nin,p)
252 CALL mpi_isend(
253 . recv_buf(nin,p)%P(1),k,mpi_integer,it_spmd(p),msgtyp,
254 . spmd_comm_world,req_s(p),ierror)
255 ENDIF
256 IF(p /= ispmd + 1 .AND. siz_send_buf(nin,p) > 0) THEN
257 k = siz_send_buf(nin,p)
258 CALL mpi_irecv(
259 . send_buf(nin,p)%P(1),k,mpi_integer,it_spmd(p),msgtyp,
260 . spmd_comm_world,req_r(p),ierror)
261 ENDIF
262 ENDDO
263 DO p = 1,nspmd
264 IF(p /= ispmd + 1 .AND. siz_recv_buf(nin,p) > 0) THEN
265 CALL mpi_wait(req_s(p),status,ierror)
266 ENDIF
267 IF(p /= ispmd + 1 .AND. siz_send_buf(nin,p) > 0) THEN
268 CALL mpi_wait(req_r(p),status,ierror)
269 ENDIF
270 ENDDO
271 ENDDO ! NIN
272
273
274 siz_send_buf(1:ninter,1:nspmd) = 0
275 DO nin = 1,ninter
276 ity=ipari(macro_ity,nin)
277 nsn = ipari(macro_nsn,nin)
278 IF( ity == 25) THEN
279 ALLOCATE(intbuf_tab(nin)%NSV_ON_PMAIN(nsn))
280 DO i = 1, nsn
281 p = main_proc(intbuf_tab(nin)%NSV(i))
282 siz_send_buf(nin,p) = siz_send_buf(nin,p) + 1
283 ii = siz_send_buf(nin,p)
284 intbuf_tab(nin)%NSV_ON_PMAIN(i) = send_buf(nin,p)%P(ii)
285 ENDDO
286
287 ENDIF
288 ENDDO
289
290 DO nin = 1,ninter
291 ity=ipari(macro_ity,nin)
292 nsn = ipari(macro_nsn,nin)
293 IF( ity == 25 ) THEN
294 DO p=1,nspmd
295 DEALLOCATE(send_buf(nin,p)%P)
296 DEALLOCATE(recv_buf(nin,p)%P)
297 ENDDO
298 ENDIF
299 ENDDO
300#endif
301 RETURN
302 END
303!||====================================================================
304!|| spmd_i25front_nor ../engine/source/mpi/interfaces/spmd_i25front.F
305!||--- called by ------------------------------------------------------
306!|| inttri ../engine/source/interfaces/intsort/inttri.F
307!||--- calls -----------------------------------------------------
308!||--- uses -----------------------------------------------------
309!|| intbufdef_mod ../common_source/modules/interfaces/intbufdef_mod.F90
310!|| message_mod ../engine/share/message_module/message_mod.F
311!|| mpi_commod ../engine/share/modules/mpi_comm_mod.F
312!|| spmd_comm_world_mod ../engine/source/mpi/spmd_comm_world.F90
313!|| tri25ebox ../engine/share/modules/tri25ebox.F
314!|| tri7box ../engine/share/modules/tri7box.F
315!||====================================================================
316 SUBROUTINE spmd_i25front_nor(IPARI,
317 . INTBUF_TAB,
318 . INTLIST25,
319 . X) ! list of type 25 interfaces
320
321C Update normal (EDG_BISSECTOR and VTX_BISSECTOR)
322C on secnd node that are candidates on other domain
323C out: updated value of EDGE_BISECTOR_FIE and VTX_BISECTOR_FIE
324C and LEDGE (important if rupture) is modified here
325C-----------------------------------------------
326C M o d u l e s
327C-----------------------------------------------
328 USE tri25ebox
329 USE tri7box
330 USE message_mod
331 USE mpi_commod
332 USE intbufdef_mod
333#ifdef WITH_ASSERT
334 use, intrinsic :: iso_fortran_env
335 use, intrinsic :: ieee_arithmetic
336#endif
337C-----------------------------------------------
338C I m p l i c i t T y p e s
339C-----------------------------------------------
340 USE spmd_comm_world_mod, ONLY : spmd_comm_world
341#include "implicit_f.inc"
342C-----------------------------------------------
343C M e s s a g e P a s s i n g
344C-----------------------------------------------
345#include "spmd.inc"
346C-----------------------------------------------
347C C o m m o n B l o c k s
348C-----------------------------------------------
349#include "com01_c.inc"
350#include "com04_c.inc"
351#include "param_c.inc"
352#include "task_c.inc"
353#include "i25edge_c.inc"
354C-----------------------------------------------
355C D u m m y A r g u m e n t s
356C-----------------------------------------------
357 INTEGER :: IPARI(NPARI,NINTER), INTLIST25(*)
358 TYPE(INTBUF_STRUCT_) :: INTBUF_TAB(*)
359 my_real, INTENT(IN) :: X(3,NUMNOD)
360C-----------------------------------------------
361C L o c a l V a r i a b l e s
362C-----------------------------------------------
363#ifdef MPI
364 INTEGER :: NI25
365 INTEGER :: N
366 INTEGER :: NB
367 INTEGER :: IEDGE
368 INTEGER :: P
369 INTEGER :: SEND_SIZE
370 INTEGER :: RECV_SIZE
371 INTEGER :: K,I,I1,I2,IE,JE,IED,L,L0
372 TYPE(mpi_comm_nor_struct) , DIMENSION(NINTER25) :: BUFFERS
373 INTEGER :: IERROR
374 INTEGER :: MSGTYP,MSGOFF
375 INTEGER :: EID
376 INTEGER :: IM,M1,M2
377 INTEGER :: NEDGE_LOCAL
378 INTEGER :: IGLOB,IBEGIN
379 INTEGER :: N1,N2,N3,N4,NN1,NN2,PP,TYPEDG
380 INTEGER :: NRTM
381 REAL *4 :: SP
382#ifdef WITH_ASSERT
383 real (real32) :: nan32
384#endif
385 DATA msgoff/163/
386 INTEGER NTRIA(3,4)
387 DATA NTRIA/1,2,4,2,4,1,0,0,0,4,1,2/
388C-----------------------------------------------
389
390#ifdef MYREAL8
391 INTEGER, PARAMETER :: NB_VALUES = 3 + 12 + 4 + 12 * 2 + 3 * 2!3 vtx bissector + 6 x 2 VTX_BISECTOR
392 ! + LEDGE right / LEDGE left
393#else
394 INTEGER, PARAMETER :: NB_VALUES = 3 + 12 + 4 + 12 + 3 * 2!3 vtx bissector + 6 x 2 VTX_BISECTOR
395 ! + LEDGE right / LEDGE left
396#endif
397
398
399
400C WRITE(6,*) __FILE__,__LINE__
401C===== ALlocations && IRECV
402 DO ni25=1,ninter25
403 n = intlist25(ni25)
404 iedge = ipari(58,n)
405 IF( iedge > 0 ) THEN
406 buffers(ni25)%NBIRECV = 0
407
408 ALLOCATE(buffers(ni25)%SEND_RQ(nspmd))
409 ALLOCATE(buffers(ni25)%RECV_RQ(nspmd))
410 ALLOCATE(buffers(ni25)%IAD_RECV(nspmd+1))
411 ALLOCATE(buffers(ni25)%IAD_SEND(nspmd+1))
412
413 send_size = 0
414 recv_size = 0
415 buffers(ni25)%IAD_SEND(1) = 1
416 buffers(ni25)%IAD_RECV(1) = 1
417 DO p = 1, nspmd
418 nb = nsnsie(n)%P(p)
419 send_size = send_size + nb*nb_values
420 buffers(ni25)%IAD_SEND(p+1) = buffers(ni25)%IAD_SEND(p) + nb
421
422 nb = nsnfie(n)%P(p)
423 recv_size = recv_size + nb*nb_values
424 buffers(ni25)%IAD_RECV(p+1) = buffers(ni25)%IAD_RECV(p) + nb
425 ENDDO
426 ALLOCATE(buffers(ni25)%SEND_BUF(send_size))
427 ALLOCATE(buffers(ni25)%RECV_BUF(recv_size))
428 DO p = 1, nspmd
429 recv_size = nb_values * ( buffers(ni25)%IAD_RECV(p+1)-buffers(ni25)%IAD_RECV(p))
430 buffers(ni25)%RECV_RQ(p) = mpi_request_null
431 IF(recv_size > 0) THEN
432 buffers(ni25)%NBIRECV = buffers(ni25)%NBIRECV + 1
433 msgtyp = msgoff
434 i = buffers(ni25)%IAD_RECV(p)
435 l = (i-1) * nb_values + 1
436 assert(l > 0)
437
438C WRITE(6,*) "RECV FROM",IT_SPMD(P),RECV_SIZE
439C CALL FLUSH(6)
440
441 CALL mpi_irecv(
442 . buffers(ni25)%RECV_BUF(l),
443 . recv_size,
444 . mpi_real4,
445 . it_spmd(p),
446 . msgtyp,
447 . spmd_comm_world,
448 . buffers(ni25)%RECV_RQ(p),
449 . ierror)
450 ENDIF
451 ENDDO
452 ENDIF
453 ENDDO
454
455
456C === Fill buffer and ISEND
457 DO ni25=1,ninter25
458 n = intlist25(ni25)
459 iedge = ipari(58,n)
460 IF( iedge > 0 ) THEN
461 nedge_local = intbuf_tab(n)%NB_INTERNAL_EDGES + intbuf_tab(n)%NB_BOUNDARY_EDGES_LOCAL
462
463 buffers(ni25)%NBISEND = 0
464 DO p = 1, nspmd
465 buffers(ni25)%SEND_RQ(p) = mpi_request_null
466 send_size = ( buffers(ni25)%IAD_SEND(p+1)-buffers(ni25)%IAD_SEND(p)) * nb_values
467 DO i = buffers(ni25)%IAD_SEND(p), buffers(ni25)%IAD_SEND(p+1)-1
468 ied = nsvsie(n)%P(i)
469C ASSERT(IED <= NEDGE)
470 ie = intbuf_tab(n)%LEDGE((ied-1)*nledge+1)
471 je = intbuf_tab(n)%LEDGE((ied-1)*nledge+2)
472
473 debug_e2e(intbuf_tab(n)%LEDGE((ied-1)*nledge+8) == d_es,p-1)
474 m1 = intbuf_tab(n)%LEDGE(5+(ied-1)*nledge)
475 m2 = intbuf_tab(n)%LEDGE(6+(ied-1)*nledge)
476 im = intbuf_tab(n)%LEDGE(10+(ied-1)*nledge)
477C ASSERT( (IM == 1 .OR. IM == -1) )
478 typedg = intbuf_tab(n)%LEDGE((ied-1)*nledge+7)
479 IF(typedg == 1 .AND. ie > 0) THEN
480 nn1 = intbuf_tab(n)%ADMSR(je+(ie-1)*4)
481 nn2 = intbuf_tab(n)%ADMSR(mod(je,4)+1+(ie-1)*4)
482 ELSE
483 nn1 = 0
484 nn2 = 0
485 ENDIF
486
487C C====================================================
488C IF( IM /= 1 .AND. IM /= -1) THEN
489C WRITE(6,*) "IM=",IM
490C CALL FLUSH(6)
491C ASSERT(.FALSE.)
492C ENDIF
493C ASSERT(IE > 0)
494C ASSERT(JE > 0)
495C IF(INTBUF_TAB(N)%IRECTM(JE+(IE-1)*4)==M1 .AND.
496C . INTBUF_TAB(N)%IRECTM(MOD(JE,4)+1+(IE-1)*4) == M2)THEN
497C IM= 1
498C ELSEIF(INTBUF_TAB(N)%IRECTM(JE+(IE-1)*4)==M2 .AND.
499C . INTBUF_TAB(N)%IRECTM(MOD(JE,4)+1+(IE-1)*4)==M1)THEN
500C IM=-1
501C ELSE
502C ASSERT(.FALSE.)
503C ENDIF
504C
505C IF(IM==1)THEN
506C ASSERT(JE+(IE-1)*4 > 0)
507C ASSERT(MOD(JE,4)+1+(IE-1)*4 > 0)
508C ASSERT(JE+(IE-1)*4 <= INTBUF_TAB(N)%S_ADMSR)
509C ASSERT(MOD(JE,4)+1+(IE-1)*4 <= INTBUF_TAB(N)%S_ADMSR)
510C I1=INTBUF_TAB(N)%ADMSR(JE+(IE-1)*4)
511C I2=INTBUF_TAB(N)%ADMSR(MOD(JE,4)+1+(IE-1)*4)
512C ELSE IF(IM==-1) THEN
513C ASSERT(JE+(IE-1)*4 > 0)
514C ASSERT(MOD(JE,4)+1+(IE-1)*4 > 0)
515C ASSERT(JE+(IE-1)*4 <= INTBUF_TAB(N)%S_ADMSR)
516C ASSERT(MOD(JE,4)+1+(IE-1)*4 <= INTBUF_TAB(N)%S_ADMSR)
517C I2=INTBUF_TAB(N)%ADMSR(JE+(IE-1)*4)
518C I1=INTBUF_TAB(N)%ADMSR(MOD(JE,4)+1+(IE-1)*4)
519C ELSE
520C WRITE(6,*) "IM=",IM
521C CALL FLUSH(6)
522C ASSERT(.FALSE.)
523C C STOP
524C END IF
525C C====================================================
526
527C ASSERT(I1 == INTBUF_TAB(N)%LEDGE(11+(IED-1)*NLEDGE))
528C ASSERT(I2 == INTBUF_TAB(N)%LEDGE(12+(IED-1)*NLEDGE))
529
530 nrtm = ipari(4,n)
531
532 printif(4 +(ie-1)*4 > 4 * nrtm,nrtm)
533 printif(4 +(ie-1)*4 > 4 * nrtm,ie)
534
535 IF(ie > 0) THEN ! Used only for solids edges
536 IF(intbuf_tab(n)%IRECTM(3 +(ie-1)*4)
537 . /= intbuf_tab(n)%IRECTM(4 +(ie-1)*4) ) THEN
538 n1 = intbuf_tab(n)%IRECTM( je +(ie-1)*4)
539 n2 = intbuf_tab(n)%IRECTM(mod(je,4) +1+(ie-1)*4)
540 n3 = intbuf_tab(n)%IRECTM(mod(je+1,4)+1+(ie-1)*4)
541 n4 = intbuf_tab(n)%IRECTM(mod(je+2,4)+1+(ie-1)*4)
542 ELSE
543 n1 = intbuf_tab(n)%IRECTM(ntria(1,je)+(ie-1)*4)
544 n2 = intbuf_tab(n)%IRECTM(ntria(2,je)+(ie-1)*4)
545 n3 = intbuf_tab(n)%IRECTM(ntria(3,je)+(ie-1)*4)
546 n4 = n3
547 END IF
548 ENDIF
549
550 i1 = intbuf_tab(n)%LEDGE(11+(ied-1)*nledge)
551 i2 = intbuf_tab(n)%LEDGE(12+(ied-1)*nledge)
552
553 l = (i-1) * nb_values
554
555 IF(ie > 0) THEN
556 buffers(ni25)%SEND_BUF(l+1) = intbuf_tab(n)%EDGE_BISECTOR(to1d(1,je,ie,3,4))
557 buffers(ni25)%SEND_BUF(l+2) = intbuf_tab(n)%EDGE_BISECTOR(to1d(2,je,ie,3,4))
558 buffers(ni25)%SEND_BUF(l+3) = intbuf_tab(n)%EDGE_BISECTOR(to1d(3,je,ie,3,4))
559
560 ELSEIF (ie < 0) THEN
561C Local IRECT of the edge is deleted
562C But ISPMD is still in charge to send that secondary edge.
563 buffers(ni25)%SEND_BUF(l+1) = intbuf_tab(n)%EDGE_BISECTOR(to1d(1,je,-ie,3,4))
564 buffers(ni25)%SEND_BUF(l+2) = intbuf_tab(n)%EDGE_BISECTOR(to1d(2,je,-ie,3,4))
565 buffers(ni25)%SEND_BUF(l+3) = intbuf_tab(n)%EDGE_BISECTOR(to1d(3,je,-ie,3,4))
566 ELSE
567 assert(.false.)
568 ENDIF
569C WRITE(6,*) "Out:",INTBUF_TAB(N)%LEDGE((IED-1)*NLEDGE+8),L
570
571 buffers(ni25)%SEND_BUF(l+4) = intbuf_tab(n)%VTX_BISECTOR(to1d(1,1,i1,3,2))
572 buffers(ni25)%SEND_BUF(l+5) = intbuf_tab(n)%VTX_BISECTOR(to1d(2,1,i1,3,2))
573 buffers(ni25)%SEND_BUF(l+6) = intbuf_tab(n)%VTX_BISECTOR(to1d(3,1,i1,3,2))
574 buffers(ni25)%SEND_BUF(l+7) = intbuf_tab(n)%VTX_BISECTOR(to1d(1,2,i1,3,2))
575 buffers(ni25)%SEND_BUF(l+8) = intbuf_tab(n)%VTX_BISECTOR(to1d(2,2,i1,3,2))
576 buffers(ni25)%SEND_BUF(l+9) = intbuf_tab(n)%VTX_BISECTOR(to1d(3,2,i1,3,2))
577 buffers(ni25)%SEND_BUF(l+10) = intbuf_tab(n)%VTX_BISECTOR(to1d(1,1,i2,3,2))
578 buffers(ni25)%SEND_BUF(l+11) = intbuf_tab(n)%VTX_BISECTOR(to1d(2,1,i2,3,2))
579 buffers(ni25)%SEND_BUF(l+12) = intbuf_tab(n)%VTX_BISECTOR(to1d(3,1,i2,3,2))
580 buffers(ni25)%SEND_BUF(l+13) = intbuf_tab(n)%VTX_BISECTOR(to1d(1,2,i2,3,2))
581 buffers(ni25)%SEND_BUF(l+14) = intbuf_tab(n)%VTX_BISECTOR(to1d(2,2,i2,3,2))
582 buffers(ni25)%SEND_BUF(l+15) = intbuf_tab(n)%VTX_BISECTOR(to1d(3,2,i2,3,2))
583
584
585
586 assert( .NOT. isnan( buffers(ni25)%SEND_BUF(l+5) ))
587 assert( .NOT. isnan( buffers(ni25)%SEND_BUF(l+6) ))
588 assert( .NOT. isnan( buffers(ni25)%SEND_BUF(l+7) ))
589 assert( .NOT. isnan( buffers(ni25)%SEND_BUF(l+8) ))
590 assert( .NOT. isnan( buffers(ni25)%SEND_BUF(l+9) ))
591 assert( .NOT. isnan( buffers(ni25)%SEND_BUF(l+10)))
592 assert( .NOT. isnan( buffers(ni25)%SEND_BUF(l+11)))
593 assert( .NOT. isnan( buffers(ni25)%SEND_BUF(l+12)))
594 assert( .NOT. isnan( buffers(ni25)%SEND_BUF(l+13)))
595 assert( .NOT. isnan( buffers(ni25)%SEND_BUF(l+14)))
596 assert( .NOT. isnan( buffers(ni25)%SEND_BUF(l+15)))
597
598 buffers(ni25)%SEND_BUF(l+16) =
599 . transfer(intbuf_tab(n)%LEDGE((ied-1)*nledge+1),buffers(ni25)%SEND_BUF(l+16))
600 buffers(ni25)%SEND_BUF(l+17) =
601 . transfer(intbuf_tab(n)%LEDGE((ied-1)*nledge+2),buffers(ni25)%SEND_BUF(l+17))
602 buffers(ni25)%SEND_BUF(l+18) =
603 . transfer(intbuf_tab(n)%LEDGE((ied-1)*nledge+3),buffers(ni25)%SEND_BUF(l+18))
604 buffers(ni25)%SEND_BUF(l+19) =
605 . transfer(intbuf_tab(n)%LEDGE((ied-1)*nledge+4),buffers(ni25)%SEND_BUF(l+19))
606
607 eid = intbuf_tab(n)%LEDGE((ied-1)*nledge+ledge_global_id)
608 debug_e2e(eid==d_es, intbuf_tab(n)%LEDGE((ied-1)*nledge+3))
609
610
611 debug_e2e(intbuf_tab(n)%LEDGE((ied-1)*nledge+8) == d_es,x(1,n1))
612 debug_e2e(intbuf_tab(n)%LEDGE((ied-1)*nledge+8) == d_es,x(2,n1))
613 debug_e2e(intbuf_tab(n)%LEDGE((ied-1)*nledge+8) == d_es,x(3,n1))
614 debug_e2e(intbuf_tab(n)%LEDGE((ied-1)*nledge+8) == d_es,x(1,n2))
615 debug_e2e(intbuf_tab(n)%LEDGE((ied-1)*nledge+8) == d_es,x(2,n2))
616 debug_e2e(intbuf_tab(n)%LEDGE((ied-1)*nledge+8) == d_es,x(3,n2))
617 debug_e2e(intbuf_tab(n)%LEDGE((ied-1)*nledge+8) == d_es,x(1,n3))
618 debug_e2e(intbuf_tab(n)%LEDGE((ied-1)*nledge+8) == d_es,x(2,n3))
619 debug_e2e(intbuf_tab(n)%LEDGE((ied-1)*nledge+8) == d_es,x(3,n3))
620 debug_e2e(intbuf_tab(n)%LEDGE((ied-1)*nledge+8) == d_es,x(1,n4))
621 debug_e2e(intbuf_tab(n)%LEDGE((ied-1)*nledge+8) == d_es,x(2,n4))
622 debug_e2e(intbuf_tab(n)%LEDGE((ied-1)*nledge+8) == d_es,x(3,n4))
623
624 IF(ie > 0) THEN
625#ifdef MYREAL8
626 buffers(ni25)%SEND_BUF(l+20:l+20+1)=transfer(x(1,n1),sp,2)
627 buffers(ni25)%SEND_BUF(l+22:l+22+1)=transfer(x(2,n1),sp,2)
628 buffers(ni25)%SEND_BUF(l+24:l+24+1)=transfer(x(3,n1),sp,2)
629 buffers(ni25)%SEND_BUF(l+26:l+26+1)=transfer(x(1,n2),sp,2)
630 buffers(ni25)%SEND_BUF(l+28:l+28+1)=transfer(x(2,n2),sp,2)
631 buffers(ni25)%SEND_BUF(l+30:l+30+1)=transfer(x(3,n2),sp,2)
632 buffers(ni25)%SEND_BUF(l+32:l+32+1)=transfer(x(1,n3),sp,2)
633 buffers(ni25)%SEND_BUF(l+34:l+34+1)=transfer(x(2,n3),sp,2)
634 buffers(ni25)%SEND_BUF(l+36:l+36+1)=transfer(x(3,n3),sp,2)
635 buffers(ni25)%SEND_BUF(l+38:l+38+1)=transfer(x(1,n4),sp,2)
636 buffers(ni25)%SEND_BUF(l+40:l+40+1)=transfer(x(2,n4),sp,2)
637 buffers(ni25)%SEND_BUF(l+42:l+42+1)=transfer(x(3,n4),sp,2)
638 pp = 43
639#else
640 buffers(ni25)%SEND_BUF(l+20) = x(1,n1)
641 buffers(ni25)%SEND_BUF(l+21) = x(2,n1)
642 buffers(ni25)%SEND_BUF(l+22) = x(3,n1)
643 buffers(ni25)%SEND_BUF(l+23) = x(1,n2)
644 buffers(ni25)%SEND_BUF(l+24) = x(2,n2)
645 buffers(ni25)%SEND_BUF(l+25) = x(3,n2)
646 buffers(ni25)%SEND_BUF(l+26) = x(1,n3)
647 buffers(ni25)%SEND_BUF(l+27) = x(2,n3)
648 buffers(ni25)%SEND_BUF(l+28) = x(3,n3)
649 buffers(ni25)%SEND_BUF(l+29) = x(1,n4)
650 buffers(ni25)%SEND_BUF(l+30) = x(2,n4)
651 buffers(ni25)%SEND_BUF(l+31) = x(3,n4)
652 pp = 31
653#endif
654 ELSE
655#ifdef MYREAL8
656 buffers(ni25)%SEND_BUF(l+20) = 0
657 buffers(ni25)%SEND_BUF(l+21) = 0
658 buffers(ni25)%SEND_BUF(l+22) = 0
659 buffers(ni25)%SEND_BUF(l+23) = 0
660 buffers(ni25)%SEND_BUF(l+24) = 0
661 buffers(ni25)%SEND_BUF(l+25) = 0
662 buffers(ni25)%SEND_BUF(l+26) = 0
663 buffers(ni25)%SEND_BUF(l+27) = 0
664 buffers(ni25)%SEND_BUF(l+28) = 0
665 buffers(ni25)%SEND_BUF(l+29) = 0
666 buffers(ni25)%SEND_BUF(l+30) = 0
667 buffers(ni25)%SEND_BUF(l+31) = 0
668 buffers(ni25)%SEND_BUF(l+32) = 0
669 buffers(ni25)%SEND_BUF(l+33) = 0
670 buffers(ni25)%SEND_BUF(l+34) = 0
671 buffers(ni25)%SEND_BUF(l+35) = 0
672 buffers(ni25)%SEND_BUF(l+36) = 0
673 buffers(ni25)%SEND_BUF(l+37) = 0
674 buffers(ni25)%SEND_BUF(l+38) = 0
675 buffers(ni25)%SEND_BUF(l+39) = 0
676 buffers(ni25)%SEND_BUF(l+40) = 0
677 buffers(ni25)%SEND_BUF(l+41) = 0
678 buffers(ni25)%SEND_BUF(l+42) = 0
679 buffers(ni25)%SEND_BUF(l+43) = 0
680
681 pp = 43
682#else
683 buffers(ni25)%SEND_BUF(l+20) = 0
684 buffers(ni25)%SEND_BUF(l+21) = 0
685 buffers(ni25)%SEND_BUF(l+22) = 0
686 buffers(ni25)%SEND_BUF(l+23) = 0
687 buffers(ni25)%SEND_BUF(l+24) = 0
688 buffers(ni25)%SEND_BUF(l+25) = 0
689 buffers(ni25)%SEND_BUF(l+26) = 0
690 buffers(ni25)%SEND_BUF(l+27) = 0
691 buffers(ni25)%SEND_BUF(l+28) = 0
692 buffers(ni25)%SEND_BUF(l+29) = 0
693 buffers(ni25)%SEND_BUF(l+30) = 0
694 buffers(ni25)%SEND_BUF(l+31) = 0
695 pp = 31
696#endif
697 ENDIF
698
699C Send Normal of nodes of second edge : edge solid
700
701 IF(typedg == 1 .AND. nn1 > 0) THEN
702
703 buffers(ni25)%SEND_BUF(l+pp+1) = intbuf_tab(n)%E2S_NOD_NORMAL(3*(nn1-1)+1)
704 buffers(ni25)%SEND_BUF(l+pp+2) = intbuf_tab(n)%E2S_NOD_NORMAL(3*(nn1-1)+2)
705 buffers(ni25)%SEND_BUF(l+pp+3) = intbuf_tab(n)%E2S_NOD_NORMAL(3*(nn1-1)+3)
706
707 buffers(ni25)%SEND_BUF(l+pp+4) = intbuf_tab(n)%E2S_NOD_NORMAL(3*(nn2-1)+1)
708 buffers(ni25)%SEND_BUF(l+pp+5) = intbuf_tab(n)%E2S_NOD_NORMAL(3*(nn2-1)+2)
709 buffers(ni25)%SEND_BUF(l+pp+6) = intbuf_tab(n)%E2S_NOD_NORMAL(3*(nn2-1)+3)
710 ELSE
711 buffers(ni25)%SEND_BUF(l+pp+1) = 0
712 buffers(ni25)%SEND_BUF(l+pp+2) = 0
713 buffers(ni25)%SEND_BUF(l+pp+3) = 0
714
715 buffers(ni25)%SEND_BUF(l+pp+4) = 0
716 buffers(ni25)%SEND_BUF(l+pp+5) = 0
717 buffers(ni25)%SEND_BUF(l+pp+6) = 0
718 ENDIF
719
720C EID = INTBUF_TAB(N)%LEDGE((IED-1)*NLEDGE+8)
721
722 ENDDO
723 IF(send_size > 0) THEN
724 buffers(ni25)%NBISEND = buffers(ni25)%NBISEND + 1
725 msgtyp = msgoff
726 i = buffers(ni25)%IAD_SEND(p)
727 l = (i-1) * nb_values+1
728C WRITE(6,*) "SEND TO",IT_SPMD(P),SEND_SIZE
729C CALL FLUSH(6)
730
731 CALL mpi_isend(
732 . buffers(ni25)%SEND_BUF(l),
733 . send_size,
734 . mpi_real4,
735 . it_spmd(p),
736 . msgtyp,
737 . spmd_comm_world,
738 . buffers(ni25)%SEND_RQ(p),
739 . ierror)
740 ENDIF
741 ENDDO ! PROC
742 ENDIF
743 ENDDO ! NI25
744
745#ifdef WITH_ASSERT
746 nan32 = ieee_value(nan32,ieee_quiet_nan)
747 do ni25=1,ninter25
748 n = intlist25(ni25)
749 iedge = ipari(58,n)
750 if( iedge > 0 ) then
751 do p = 1,nspmd
752 iglob = 0
753 do i =1,nsnfie(n)%p(p)
754 iglob = iglob + 1
755 edg_bisector_fie(n)%p(1,1,iglob) = nan32
756 edg_bisector_fie(n)%p(2,1,iglob) = nan32
757 edg_bisector_fie(n)%p(3,1,iglob) = nan32
758 vtx_bisector_fie(n)%p(1,1,iglob) = nan32
759 vtx_bisector_fie(n)%p(2,1,iglob) = nan32
760 vtx_bisector_fie(n)%p(3,1,iglob) = nan32
761 vtx_bisector_fie(n)%p(1,2,iglob) = nan32
762 vtx_bisector_fie(n)%p(2,2,iglob) = nan32
763 vtx_bisector_fie(n)%p(3,2,iglob) = nan32
764 vtx_bisector_fie(n)%p(1,3,iglob) = nan32
765 vtx_bisector_fie(n)%p(2,3,iglob) = nan32
766 vtx_bisector_fie(n)%p(3,3,iglob) = nan32
767 vtx_bisector_fie(n)%p(1,4,iglob) = nan32
768 vtx_bisector_fie(n)%p(2,4,iglob) = nan32
769 vtx_bisector_fie(n)%p(3,4,iglob) = nan32
770 enddo
771 enddo
772 endif
773 enddo
774#endif
775
776
777 DO ni25=1,ninter25
778 n = intlist25(ni25)
779 iedge = ipari(58,n)
780 IF( iedge > 0 ) THEN
781 DO k = 1,buffers(ni25)%NBIRECV
782C WRITE(6,*) "Go into wait"
783 CALL flush(6)
784 CALL mpi_waitany(nspmd,buffers(ni25)%RECV_RQ,p,mpi_status_ignore,ierror)
785C WRITE(6,*) "RECEIVE FROM",P-1
786 l0 = (buffers(ni25)%IAD_RECV(p) - 1)*nb_values
787 ibegin = 0
788 IF( p > 1) THEN
789 DO l = 1,p-1
790 IF( l - 1 /= ispmd) ibegin = ibegin + nsnfie(n)%P(l)
791 ENDDO
792 ENDIF
793 DO i =1,nsnfie(n)%P(p)
794 l = l0 + (i-1) * nb_values
795 iglob = i + ibegin
796c WRITE(6,*) iglob,"In:",LEDGE_FIE%P(E_GLOBAL_ID,IGLOB),L
797 edg_bisector_fie(n)%P(1,1,iglob) = buffers(ni25)%RECV_BUF(l+1 )
798 edg_bisector_fie(n)%P(2,1,iglob) = buffers(ni25)%RECV_BUF(l+2 )
799 edg_bisector_fie(n)%P(3,1,iglob) = buffers(ni25)%RECV_BUF(l+3 )
800 vtx_bisector_fie(n)%P(1,1,iglob) = buffers(ni25)%RECV_BUF(l+4 )
801 vtx_bisector_fie(n)%P(2,1,iglob) = buffers(ni25)%RECV_BUF(l+5 )
802 vtx_bisector_fie(n)%P(3,1,iglob) = buffers(ni25)%RECV_BUF(l+6 )
803 vtx_bisector_fie(n)%P(1,2,iglob) = buffers(ni25)%RECV_BUF(l+7 )
804 vtx_bisector_fie(n)%P(2,2,iglob) = buffers(ni25)%RECV_BUF(l+8 )
805 vtx_bisector_fie(n)%P(3,2,iglob) = buffers(ni25)%RECV_BUF(l+9 )
806 vtx_bisector_fie(n)%P(1,3,iglob) = buffers(ni25)%RECV_BUF(l+10)
807 vtx_bisector_fie(n)%P(2,3,iglob) = buffers(ni25)%RECV_BUF(l+11)
808 vtx_bisector_fie(n)%P(3,3,iglob) = buffers(ni25)%RECV_BUF(l+12)
809 vtx_bisector_fie(n)%P(1,4,iglob) = buffers(ni25)%RECV_BUF(l+13)
810 vtx_bisector_fie(n)%P(2,4,iglob) = buffers(ni25)%RECV_BUF(l+14)
811 vtx_bisector_fie(n)%P(3,4,iglob) = buffers(ni25)%RECV_BUF(l+15)
812
813
814 assert(.NOT. isnan(buffers(ni25)%RECV_BUF(l+1 )))
815 assert(.NOT. isnan(buffers(ni25)%RECV_BUF(l+2 )))
816 assert(.NOT. isnan(buffers(ni25)%RECV_BUF(l+3 )))
817 assert(.NOT. isnan(buffers(ni25)%RECV_BUF(l+4 )))
818 assert(.NOT. isnan(buffers(ni25)%RECV_BUF(l+5 )))
819 assert(.NOT. isnan(buffers(ni25)%RECV_BUF(l+6 )))
820 assert(.NOT. isnan(buffers(ni25)%RECV_BUF(l+7 )))
821 assert(.NOT. isnan(buffers(ni25)%RECV_BUF(l+8 )))
822 assert(.NOT. isnan(buffers(ni25)%RECV_BUF(l+9 )))
823 assert(.NOT. isnan(buffers(ni25)%RECV_BUF(l+10)))
824 assert(.NOT. isnan(buffers(ni25)%RECV_BUF(l+11)))
825 assert(.NOT. isnan(buffers(ni25)%RECV_BUF(l+12)))
826 assert(.NOT. isnan(buffers(ni25)%RECV_BUF(l+13)))
827 assert(.NOT. isnan(buffers(ni25)%RECV_BUF(l+14)))
828 assert(.NOT. isnan(buffers(ni25)%RECV_BUF(l+15)))
829
830 ledge_fie(n)%P(e_left_seg ,iglob) =
831 . transfer(buffers(ni25)%RECV_BUF(l+16),l0)
832 ledge_fie(n)%P(e_left_id,iglob) =
833 . transfer(buffers(ni25)%RECV_BUF(l+17),l0)
834 ledge_fie(n)%P(e_right_seg ,iglob) =
835 . transfer(buffers(ni25)%RECV_BUF(l+18),l0)
836 ledge_fie(n)%P(e_right_id,iglob) =
837 . transfer(buffers(ni25)%RECV_BUF(l+19),l0)
838
839
840
841 debug_e2e(ledge_fie(n)%P(1,iglob)==d_es, ledge_fie(n)%P(e_right_seg ,iglob))
842
843#ifdef MYREAL8
844 x_seg_fie(n)%P(1,1,iglob) =transfer( buffers(ni25)%RECV_BUF(l+20:l+20+1),one)
845 x_seg_fie(n)%P(2,1,iglob) =transfer( buffers(ni25)%RECV_BUF(l+22:l+22+1),one)
846 x_seg_fie(n)%P(3,1,iglob) =transfer( buffers(ni25)%RECV_BUF(l+24:l+24+1),one)
847 x_seg_fie(n)%P(1,2,iglob) =transfer( buffers(ni25)%RECV_BUF(l+26:l+26+1),one)
848 x_seg_fie(n)%P(2,2,iglob) =transfer( buffers(ni25)%RECV_BUF(l+28:l+28+1),one)
849 x_seg_fie(n)%P(3,2,iglob) =transfer( buffers(ni25)%RECV_BUF(l+30:l+30+1),one)
850 x_seg_fie(n)%P(1,3,iglob) =transfer( buffers(ni25)%RECV_BUF(l+32:l+32+1),one)
851 x_seg_fie(n)%P(2,3,iglob) =transfer( buffers(ni25)%RECV_BUF(l+34:l+34+1),one)
852 x_seg_fie(n)%P(3,3,iglob) =transfer( buffers(ni25)%RECV_BUF(l+36:l+36+1),one)
853 x_seg_fie(n)%P(1,4,iglob) =transfer( buffers(ni25)%RECV_BUF(l+38:l+38+1),one)
854 x_seg_fie(n)%P(2,4,iglob) =transfer( buffers(ni25)%RECV_BUF(l+40:l+40+1),one)
855 x_seg_fie(n)%P(3,4,iglob) =transfer( buffers(ni25)%RECV_BUF(l+42:l+42+1),one)
856 pp = 43
857#else
858 x_seg_fie(n)%P(1,1,iglob) = buffers(ni25)%RECV_BUF(l+20)
859 x_seg_fie(n)%P(2,1,iglob) = buffers(ni25)%RECV_BUF(l+21)
860 x_seg_fie(n)%P(3,1,iglob) = buffers(ni25)%RECV_BUF(l+22)
861 x_seg_fie(n)%P(1,2,iglob) = buffers(ni25)%RECV_BUF(l+23)
862 x_seg_fie(n)%P(2,2,iglob) = buffers(ni25)%RECV_BUF(l+24)
863 x_seg_fie(n)%P(3,2,iglob) = buffers(ni25)%RECV_BUF(l+25)
864 x_seg_fie(n)%P(1,3,iglob) = buffers(ni25)%RECV_BUF(l+26)
865 x_seg_fie(n)%P(2,3,iglob) = buffers(ni25)%RECV_BUF(l+27)
866 x_seg_fie(n)%P(3,3,iglob) = buffers(ni25)%RECV_BUF(l+28)
867 x_seg_fie(n)%P(1,4,iglob) = buffers(ni25)%RECV_BUF(l+29)
868 x_seg_fie(n)%P(2,4,iglob) = buffers(ni25)%RECV_BUF(l+30)
869 x_seg_fie(n)%P(3,4,iglob) = buffers(ni25)%RECV_BUF(l+31)
870 pp = 31
871#endif
872 edg_bisector_fie(n)%P(1,2,iglob) = buffers(ni25)%RECV_BUF(l+pp+1 )
873 edg_bisector_fie(n)%P(2,2,iglob) = buffers(ni25)%RECV_BUF(l+pp+2 )
874 edg_bisector_fie(n)%P(3,2,iglob) = buffers(ni25)%RECV_BUF(l+pp+3 )
875
876 edg_bisector_fie(n)%P(1,3,iglob) = buffers(ni25)%RECV_BUF(l+pp+4 )
877 edg_bisector_fie(n)%P(2,3,iglob) = buffers(ni25)%RECV_BUF(l+pp+5 )
878 edg_bisector_fie(n)%P(3,3,iglob) = buffers(ni25)%RECV_BUF(l+pp+6 )
879
880 ENDDO
881 ENDDO ! RECV
882 ENDIF ! iedge
883 ENDDO ! NINTER25
884
885#ifdef WITH_ASSERT
886C debug mode
887 nan32 = ieee_value(nan32,ieee_quiet_nan)
888 do ni25=1,ninter25
889 n = intlist25(ni25)
890 iedge = ipari(58,n)
891 if( iedge > 0 ) then
892 do p = 1,nspmd
893 iglob = 0
894 do i =1,nsnfie(n)%p(p)
895 iglob = iglob + 1
896 assert(.NOT. ieee_is_nan(edg_bisector_fie(n)%p(1,1,iglob)))
897 assert(.NOT. ieee_is_nan(edg_bisector_fie(n)%p(2,1,iglob)))
898 assert(.NOT. ieee_is_nan(edg_bisector_fie(n)%p(3,1,iglob)))
899 assert(.NOT. ieee_is_nan(vtx_bisector_fie(n)%p(1,1,iglob)))
900 assert(.NOT. ieee_is_nan(vtx_bisector_fie(n)%p(2,1,iglob)))
901 assert(.NOT. ieee_is_nan(vtx_bisector_fie(n)%p(3,1,iglob)))
902 assert(.NOT. ieee_is_nan(vtx_bisector_fie(n)%p(1,2,iglob)))
903 assert(.NOT. ieee_is_nan(vtx_bisector_fie(n)%p(2,2,iglob)))
904 assert(.NOT. ieee_is_nan(vtx_bisector_fie(n)%p(3,2,iglob)))
905 assert(.NOT. ieee_is_nan(vtx_bisector_fie(n)%p(1,3,iglob)))
906 assert(.NOT. ieee_is_nan(vtx_bisector_fie(n)%p(2,3,iglob)))
907 assert(.NOT. ieee_is_nan(vtx_bisector_fie(n)%p(3,3,iglob)))
908 assert(.NOT. ieee_is_nan(vtx_bisector_fie(n)%p(1,4,iglob)))
909 assert(.NOT. ieee_is_nan(vtx_bisector_fie(n)%p(2,4,iglob)))
910 assert(.NOT. ieee_is_nan(vtx_bisector_fie(n)%p(3,4,iglob)))
911 enddo
912 enddo
913 endif
914 enddo
915#endif
916
917
918C --------------- Free Send request SEND wait
919 DO ni25=1,ninter25
920 n = intlist25(ni25)
921 iedge = ipari(58,n)
922 IF( iedge > 0 ) THEN
923 CALL mpi_waitall(nspmd,buffers(ni25)%SEND_RQ,mpi_statuses_ignore,ierror)
924C CALL MPI_WAITALL(NSPMD,BUFFERS(NI25)%RECV_RQ,MPI_STATUSES_IGNORE,IERROR)
925 DEALLOCATE(buffers(ni25)%SEND_BUF)
926 DEALLOCATE(buffers(ni25)%RECV_BUF)
927 DEALLOCATE(buffers(ni25)%SEND_RQ)
928 DEALLOCATE(buffers(ni25)%RECV_RQ)
929 DEALLOCATE(buffers(ni25)%IAD_RECV)
930 DEALLOCATE(buffers(ni25)%IAD_SEND)
931 ENDIF
932 ENDDO
933
934#endif
935 RETURN
936 END
937
subroutine i25free_bound(nrtm, mvoisin, irect, stifm, nrtm_free, free_irect_id)
subroutine mpi_isend(buf, cnt, datatype, dest, tag, comm, ireq, ierr)
Definition mpi.f:382
subroutine mpi_wait(ireq, status, ierr)
Definition mpi.f:525
subroutine mpi_waitall(cnt, array_of_requests, status, ierr)
Definition mpi.f:536
subroutine mpi_alltoall(sendbuf, sendcnt, sendtype, recvbuf, recvcnt, recvtype, comm, ierr)
Definition mpi.f:161
subroutine mpi_waitany(cnt, array_of_requests, index, status, ierr)
Definition mpi.f:549
subroutine mpi_irecv(buf, cnt, datatype, source, tag, comm, ireq, ierr)
Definition mpi.f:372
type(real4_pointer3), dimension(:), allocatable edg_bisector_fie
Definition tri25ebox.F:83
type(real4_pointer3), dimension(:), allocatable vtx_bisector_fie
Definition tri25ebox.F:84
type(real_pointer3), dimension(:), allocatable x_seg_fie
Definition tri25ebox.F:85
type(int_pointer2), dimension(:), allocatable ledge_fie
Definition tri25ebox.F:88
type(int_pointer), dimension(:), allocatable nsnfie
Definition tri7box.F:440
type(int_pointer), dimension(:), allocatable nsnsie
Definition tri7box.F:491
type(int_pointer), dimension(:), allocatable nsvsie
Definition tri7box.F:485
subroutine spmd_comm_split(color, key, subcomm, rank, size_l)
subroutine spmd_i25front_nor(ipari, intbuf_tab, intlist25, x)
subroutine spmd_i25front_init(itab, main_proc, intbuf_tab, ipari)