OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
dparrws.F File Reference
#include "implicit_f.inc"
#include "spmd_c.inc"
#include "com01_c.inc"
#include "com04_c.inc"
#include "param_c.inc"
#include "task_c.inc"

Go to the source code of this file.

Functions/Subroutines

subroutine dparrws (nesbw, nstrf, ixc, ixtg, x, nodcut, rwbuf, nprw, nodglob, buf, ixs)

Function/Subroutine Documentation

◆ dparrws()

subroutine dparrws ( integer nesbw,
integer, dimension(*) nstrf,
integer, dimension(nixc,*) ixc,
integer, dimension(nixtg,*) ixtg,
x,
integer nodcut,
rwbuf,
integer, dimension(*) nprw,
integer, dimension(*) nodglob,
integer buf,
integer, dimension(nixs,*) ixs )

Definition at line 31 of file dparrws.F.

34C-----------------------------------------------
35C I m p l i c i t T y p e s
36C-----------------------------------------------
37#include "implicit_f.inc"
38C-----------------------------------------------
39C C o m m o n B l o c k s
40C-----------------------------------------------
41#include "spmd_c.inc"
42#include "com01_c.inc"
43#include "com04_c.inc"
44#include "param_c.inc"
45#include "task_c.inc"
46C-----------------------------------------------
47C D u m m y A r g u m e n t s
48C-----------------------------------------------
49 INTEGER NESBW,NSTRF(*),IXC(NIXC,*),IXTG(NIXTG,*),
50 . NODCUT,NPRW(*), IXS(NIXS,*),BUF,NODGLOB(*)
52 . x(3,*),rwbuf(nrwlp,*)
53C-----------------------------------------------
54C L o c a l V a r i a b l e s
55C-----------------------------------------------
56 INTEGER J, JJ, LEN, I, K, L, KK, K0, K1, K5, K9, N,
57 . N0, N1, N2, N3, N4, N10, NSEG, NSEGC, NSEGTG, ITYP,
58 . UNPACK(15,4), II(8), N5, N6, N7, N8, NSEGS, K3,OW,
59 . WA(6*BUF+4)
61 . xx1, yy1, zz1, xx2, yy2, zz2, xx3, yy3, zz3,
62 . xx4, yy4, zz4, d13, xxc, yyc, zzc
63 INTEGER POWER2(8),IPACK
64 INTEGER :: INDICE
65 INTEGER :: MODE,SIZE_BUFFER_S,SIZE_BUFFER00_R
66 INTEGER, DIMENSION(NSPMD) :: SHIFT_R,NB_ELEM_R
67 INTEGER, DIMENSION(NSECT,3,NSPMD) :: SHIFT_SECT
68 INTEGER, DIMENSION(NSECT+1,3) ::SINDEX
69 INTEGER, DIMENSION(NSECT+1,3,NSPMD) :: RINDEX_PROC
70 INTEGER, DIMENSION(:), ALLOCATABLE :: BUFFER_S
71 INTEGER, DIMENSION(:), ALLOCATABLE :: BUFFER00_R
72
73 DATA power2/1,2,4,8,16,32,64,128/
74C-----------------------------------------------
75 DATA unpack/1,2,1,3,1,2,1,4,1,2,1,3,1,2,1,
76 . 0,0,2,0,3,3,2,0,4,4,2,4,3,3,2,
77 . 0,0,0,0,0,0,3,0,0,0,4,0,4,4,3,
78 . 0,0,0,0,0,0,0,0,0,0,0,0,0,0,4/
79
80! 1st step : count the number of value send by each proc
81! 2nd step : allocate the sended buffer for each proc and the received buffer (proc = 0)
82! 3rd step : initialize and send the sended buffer to proc 0
83! 4th step : write the data (proc=0)
84!
85!
86! Structure of the buffer :
87! *******************
88!
89! sended buffer (local on each proc) :
90! sh = shell ; t = shell3n ; s = solid
91!
92!
93! sh | t | s | sh | t | s | sh | t | s | sh | t | s | sh | t | s | sh | t | s | sh | t | s | sh | t | s | ...
94! \__________/ \__________/ \__________/ \__________/ \__________/ \__________/ \__________/ \__________/ ...
95! nsect = 1 2 3 4 ...
96!
97!
98!
99!
100! received buffer (only on proc 0) :
101!
102! sh | t | s | sh | t | s |... sh | t | s | sh | t | s | sh | t | s |... sh | t | s | ... sh | t | s | sh | t | s | ...
103! \__________/ \__________/ \__________/ \__________/ \__________/ \__________/ ...\__________/ \__________/ ...
104! nsect = 1 2 ... n || 1 2 ... n | || 1 2 ...
105! ||________________________________________||________________________________________| ...||_____________________________
106! | | |
107! proc = 0 1 ... m
108
109! Index arrays : SHIFT_R, SHIFT_SECT, RINDEX_PROC
110! SHIFT_R --> index for the processor
111! SHIFT_SECT --> index for the section for each processor
112! RINDEX_PROC --> index for the element types (shell, shell3n or solid) for each processor and each section
113!
114!
115
116! SHIFT_R(1) SHIFT_R(2) SHIFT_R(m+1)
117!
118! ||________________________________________||________________________________________| ...||_____________________________
119! | | |
120! proc = 0 1 ... m
121
122! SHIFT_SECT :
123! (1,1) (2,1) (n,1) (1,2) ...
124! \__________/ \__________/ \__________/ \__________/ \__________/ \__________/ ...\__________/ \__________/ ...
125! nsect = 1 2 ... n || 1 2 ... n | || 1 2 ...
126!
127! RINDEX_PROC
128! (1) (2) (3)
129! sh | t | s |
130!
131!
132! -----------------------------------------------
133
134! ----------------------
135! count the number of data
136! ----------------------
137 sindex(1:nsect+1,1:3) = 0
138 size_buffer00_r = 0
139 IF(ispmd==0) rindex_proc(1:nsect,1:3,1:nspmd) = 0
140 jj = 0
141
142 IF (nsect>0) THEN
143 k0 = nstrf(25)
144 DO i=1,nsect
145 n0 = numnod + nodcut + i - 1
146 k5=k0+30+nstrf(k0+14)+nstrf(k0+6)
147 1 + 2*nstrf(k0+7) +nstrf(k0+8)*2
148 nsegc = nstrf(k0+9)
149
150! ----------------
151! SHELL
152! ----------------
153 DO j=1,nsegc
154 kk = k5+2*(j-1)
155 n = nstrf(kk)
156 IF(nstrf(kk+1)/=0) THEN
157 jj = jj + 4
158 ENDIF
159 ENDDO
160 sindex(i,1) = jj
161 k9 = k5+2*nstrf(k0+9) +2*nstrf(k0+10)
162 1 +2*nstrf(k0+11)+2*nstrf(k0+12)
163 nsegtg = nstrf(k0+13)
164! ----------------
165! SHELL3N
166! ----------------
167 DO j=1,nsegtg
168 kk = k9+2*(j-1)
169 n = nstrf(kk)
170 IF(nstrf(kk+1)/=0) THEN
171 jj = jj + 4
172 ENDIF
173 ENDDO
174 sindex(i,2) = jj
175
176 k3=k0+30+nstrf(k0+14)+nstrf(k0+6)
177 nsegs=nstrf(k0+7)
178
179! ----------------
180! SOLID
181! ----------------
182 IF(nsegs/=0)THEN
183 jj = jj + 4
184 END IF
185 DO j=1,nsegs
186 kk=k3+2*(j-1)
187 ipack=nstrf(kk+1)
188 IF(ipack/=0)THEN
189 n =nstrf(kk)
190 IF (nspmd == 1) THEN
191 ii(1)=ixs(2,n)-1
192 ii(2)=ixs(3,n)-1
193 ii(3)=ixs(4,n)-1
194 ii(4)=ixs(5,n)-1
195 ii(5)=ixs(6,n)-1
196 ii(6)=ixs(7,n)-1
197 ii(7)=ixs(8,n)-1
198 ii(8)=ixs(9,n)-1
199 ELSE
200 ii(1)=nodglob(ixs(2,n))-1
201 ii(2)=nodglob(ixs(3,n))-1
202 ii(3)=nodglob(ixs(4,n))-1
203 ii(4)=nodglob(ixs(5,n))-1
204 ii(5)=nodglob(ixs(6,n))-1
205 ii(6)=nodglob(ixs(7,n))-1
206 ii(7)=nodglob(ixs(8,n))-1
207 ii(8)=nodglob(ixs(9,n))-1
208 ENDIF
209
210 IF( ii(2)==ii(1).AND.ii(4)==ii(3)
211 . .AND.ii(8)==ii(5).AND.ii(7)==ii(6))THEN
212C tetra4, tetra10
213 n1=mod(ipack/power2(1),2)
214 n2=mod(ipack/power2(3),2)
215 n3=mod(ipack/power2(5),2)
216 n4=mod(ipack/power2(6),2)
217 IF(n1/=0.AND.n2/=0.AND.n3/=0)THEN
218 jj=jj+4
219 END IF
220 IF(n1/=0.AND.n2/=0.AND.n4/=0)THEN
221 jj=jj+4
222 END IF
223 IF(n2/=0.AND.n3/=0.AND.n4/=0)THEN
224 jj=jj+4
225 END IF
226 IF(n3/=0.AND.n1/=0.AND.n4/=0)THEN
227 jj=jj+4
228 END IF
229 ELSE
230C brick, shell16, brick20
231 n1=mod(ipack/power2(1),2)
232 n2=mod(ipack/power2(2),2)
233 n3=mod(ipack/power2(3),2)
234 n4=mod(ipack/power2(4),2)
235 n5=mod(ipack/power2(5),2)
236 n6=mod(ipack/power2(6),2)
237 n7=mod(ipack/power2(7),2)
238 n8=mod(ipack/power2(8),2)
239 IF(n1/=0.AND.n2/=0.AND.n3/=0.AND.n4/=0)THEN
240 jj=jj+4
241 END IF
242 IF(n5/=0.AND.n6/=0.AND.n7/=0.AND.n8/=0)THEN
243 jj=jj+4
244 END IF
245 IF(n1/=0.AND.n5/=0.AND.n6/=0.AND.n2/=0)THEN
246 jj=jj+4
247 END IF
248 IF(n4/=0.AND.n8/=0.AND.n7/=0.AND.n3/=0)THEN
249 jj=jj+4
250 END IF
251 IF(n1/=0.AND.n4/=0.AND.n8/=0.AND.n5/=0)THEN
252 jj=jj+4
253 END IF
254 IF(n2/=0.AND.n3/=0.AND.n7/=0.AND.n6/=0)THEN
255 jj=jj+4
256 END IF
257 END IF
258 END IF
259 END DO
260 k0=nstrf(k0+24)
261 sindex(i,3) = jj
262 ENDDO
263
264 size_buffer_s = jj
265 sindex(nsect+1,1:3) = size_buffer_s
266! -----------------------------------------------
267! ----------------------
268! Allocate the sended buffer
269! ----------------------
270 ALLOCATE( buffer_s(size_buffer_s) )
271 mode = 0
272
273 ALLOCATE( buffer00_r(0) )
274
275 IF(nspmd>1) THEN
276 CALL spmd_gather_wa(mode,size_buffer_s,size_buffer00_r,sindex,rindex_proc,
277 1 buffer_s,buffer00_r,shift_r,nb_elem_r)
278 ELSE
279 size_buffer00_r = size_buffer_s
280 ! shift for the _gatherv comm
281 shift_r(1) = 0
282 rindex_proc(1:nsect,1:3,1) = sindex(1:nsect,1:3)
283 ENDIF
284
285 IF(ispmd==0) THEN
286 DEALLOCATE( buffer00_r )
287 ALLOCATE( buffer00_r(size_buffer00_r) )
288 ENDIF
289
290! ----------------------
291! initialize the sended buffer
292! ----------------------
293 k0 = nstrf(25)
294 jj = 0
295 DO i=1,nsect
296 n0 = numnod + nodcut + i - 1
297 k5=k0+30+nstrf(k0+14)+nstrf(k0+6)
298 1 + 2*nstrf(k0+7) +nstrf(k0+8)*2
299 nsegc = nstrf(k0+9)
300
301 DO j=1,nsegc
302 kk = k5+2*(j-1)
303 n = nstrf(kk)
304 IF(nstrf(kk+1)/=0) THEN
305 n1 = unpack(nstrf(kk+1),1)
306 n2 = unpack(nstrf(kk+1),2)
307 IF(n2==0)THEN
308 n2 = n1
309 n3 = n1
310 ELSE
311 n3 = unpack(nstrf(kk+1),3)
312 IF(n3==0)n3 = n2
313 ENDIF
314 IF (nspmd == 1) THEN
315 buffer_s(jj+1) = n0
316 buffer_s(jj+2) = ixc(1+n1,n)-1
317 buffer_s(jj+3) = ixc(1+n2,n)-1
318 buffer_s(jj+4) = ixc(1+n3,n)-1
319 jj = jj + 4
320 ELSE
321 buffer_s(jj+1) = numnodg + nodcut + i - 1
322 buffer_s(jj+2) = nodglob(ixc(1+n1,n))-1
323 buffer_s(jj+3) = nodglob(ixc(1+n2,n))-1
324 buffer_s(jj+4) = nodglob(ixc(1+n3,n))-1
325 jj = jj + 4
326 ENDIF
327 ENDIF
328 ENDDO
329
330 k9=k5+2*nstrf(k0+9) +2*nstrf(k0+10)
331 1 + 2*nstrf(k0+11)+2*nstrf(k0+12)
332 nsegtg = nstrf(k0+13)
333 DO j=1,nsegtg
334 kk = k9+2*(j-1)
335 n = nstrf(kk)
336 IF(nstrf(kk+1)/=0) THEN
337 n1 = unpack(nstrf(1+kk),1)
338 n2 = unpack(nstrf(1+kk),2)
339 IF(n2==0)THEN
340 n2 = n1
341 n3 = n1
342 ELSE
343 n3 = unpack(nstrf(1+kk),3)
344 IF(n3==0)n3 = n2
345 ENDIF
346 IF (nspmd == 1) THEN
347 buffer_s(jj+1) = n0
348 buffer_s(jj+2) = ixtg(1+n1,n)-1
349 buffer_s(jj+3) = ixtg(1+n2,n)-1
350 buffer_s(jj+4) = ixtg(1+n3,n)-1
351 jj = jj + 4
352 ELSE
353 buffer_s(jj+1) = numnodg + nodcut + i - 1
354 buffer_s(jj+2) = nodglob(ixtg(1+n1,n))-1
355 buffer_s(jj+3) = nodglob(ixtg(1+n2,n))-1
356 buffer_s(jj+4) = nodglob(ixtg(1+n3,n))-1
357 jj = jj + 4
358 ENDIF
359 ENDIF
360 ENDDO
361
362 k3=k0+30+nstrf(k0+14)+nstrf(k0+6)
363 nsegs=nstrf(k0+7)
364
365 IF(nsegs/=0)THEN
366 IF (nspmd == 1) THEN
367 buffer_s(jj+1) = n0
368 buffer_s(jj+2) = n0
369 buffer_s(jj+3) = n0
370 buffer_s(jj+4) = n0
371 jj = jj + 4
372 ELSE
373 buffer_s(jj+1) = numnodg + nodcut + i - 1
374 buffer_s(jj+2) = numnodg + nodcut + i - 1
375 buffer_s(jj+3) = numnodg + nodcut + i - 1
376 buffer_s(jj+4) = numnodg + nodcut + i - 1
377 jj = jj + 4
378 ENDIF
379 END IF
380C
381 DO j=1,nsegs
382 kk=k3+2*(j-1)
383 ipack=nstrf(kk+1)
384 IF(ipack/=0)THEN
385 n =nstrf(kk)
386 IF (nspmd == 1) THEN
387 ii(1)=ixs(2,n)-1
388 ii(2)=ixs(3,n)-1
389 ii(3)=ixs(4,n)-1
390 ii(4)=ixs(5,n)-1
391 ii(5)=ixs(6,n)-1
392 ii(6)=ixs(7,n)-1
393 ii(7)=ixs(8,n)-1
394 ii(8)=ixs(9,n)-1
395 ELSE
396 ii(1)=nodglob(ixs(2,n))-1
397 ii(2)=nodglob(ixs(3,n))-1
398 ii(3)=nodglob(ixs(4,n))-1
399 ii(4)=nodglob(ixs(5,n))-1
400 ii(5)=nodglob(ixs(6,n))-1
401 ii(6)=nodglob(ixs(7,n))-1
402 ii(7)=nodglob(ixs(8,n))-1
403 ii(8)=nodglob(ixs(9,n))-1
404 ENDIF
405
406 IF( ii(2)==ii(1).AND.ii(4)==ii(3)
407 . .AND.ii(8)==ii(5).AND.ii(7)==ii(6))THEN
408C tetra4, tetra10
409 n1=mod(ipack/power2(1),2)
410 n2=mod(ipack/power2(3),2)
411 n3=mod(ipack/power2(5),2)
412 n4=mod(ipack/power2(6),2)
413 IF(n1/=0.AND.n2/=0.AND.n3/=0)THEN
414 buffer_s(jj+1) =ii(1)
415 buffer_s(jj+2) =ii(3)
416 buffer_s(jj+3) =ii(5)
417 buffer_s(jj+4) =ii(5)
418 jj=jj+4
419 END IF
420 IF(n1/=0.AND.n2/=0.AND.n4/=0)THEN
421 buffer_s(jj+1) =ii(1)
422 buffer_s(jj+2) =ii(3)
423 buffer_s(jj+3) =ii(6)
424 buffer_s(jj+4) =ii(6)
425 jj=jj+4
426 END IF
427 IF(n2/=0.AND.n3/=0.AND.n4/=0)THEN
428 buffer_s(jj+1) =ii(3)
429 buffer_s(jj+2) =ii(5)
430 buffer_s(jj+3) =ii(6)
431 buffer_s(jj+4) =ii(6)
432 jj=jj+4
433 END IF
434 IF(n3/=0.AND.n1/=0.AND.n4/=0)THEN
435 buffer_s(jj+1) =ii(5)
436 buffer_s(jj+2) =ii(1)
437 buffer_s(jj+3) =ii(6)
438 buffer_s(jj+4) =ii(6)
439 jj=jj+4
440 END IF
441 ELSE
442C brick, shell16, brick20
443 n1=mod(ipack/power2(1),2)
444 n2=mod(ipack/power2(2),2)
445 n3=mod(ipack/power2(3),2)
446 n4=mod(ipack/power2(4),2)
447 n5=mod(ipack/power2(5),2)
448 n6=mod(ipack/power2(6),2)
449 n7=mod(ipack/power2(7),2)
450 n8=mod(ipack/power2(8),2)
451
452 IF(n1/=0.AND.n2/=0.AND.n3/=0.AND.n4/=0)THEN
453 buffer_s(jj+1) =ii(1)
454 buffer_s(jj+2) =ii(2)
455 buffer_s(jj+3) =ii(3)
456 buffer_s(jj+4) =ii(4)
457 jj=jj+4
458 END IF
459 IF(n5/=0.AND.n6/=0.AND.n7/=0.AND.n8/=0)THEN
460 buffer_s(jj+1) =ii(5)
461 buffer_s(jj+2) =ii(6)
462 buffer_s(jj+3) =ii(7)
463 buffer_s(jj+4) =ii(8)
464 jj=jj+4
465 END IF
466 IF(n1/=0.AND.n5/=0.AND.n6/=0.AND.n2/=0)THEN
467 buffer_s(jj+1) =ii(1)
468 buffer_s(jj+2) =ii(5)
469 buffer_s(jj+3) =ii(6)
470 buffer_s(jj+4) =ii(2)
471 jj=jj+4
472 END IF
473 IF(n4/=0.AND.n8/=0.AND.n7/=0.AND.n3/=0)THEN
474 buffer_s(jj+1) =ii(4)
475 buffer_s(jj+2) =ii(8)
476 buffer_s(jj+3) =ii(7)
477 buffer_s(jj+4) =ii(3)
478 jj=jj+4
479 END IF
480 IF(n1/=0.AND.n4/=0.AND.n8/=0.AND.n5/=0)THEN
481 buffer_s(jj+1) =ii(1)
482 buffer_s(jj+2) =ii(4)
483 buffer_s(jj+3) =ii(8)
484 buffer_s(jj+4) =ii(5)
485 jj=jj+4
486 END IF
487 IF(n2/=0.AND.n3/=0.AND.n7/=0.AND.n6/=0)THEN
488 buffer_s(jj+1) =ii(2)
489 buffer_s(jj+2) =ii(3)
490 buffer_s(jj+3) =ii(7)
491 buffer_s(jj+4) =ii(6)
492 jj=jj+4
493 END IF
494 END IF
495 END IF
496 END DO
497 k0=nstrf(k0+24)
498 ENDDO
499
500! ----------------------
501! send buffer
502! ----------------------
503
504 mode = 1
505 IF(nspmd>1) THEN
506 CALL spmd_gather_wa(mode,size_buffer_s,size_buffer00_r,sindex,rindex_proc,
507 1 buffer_s,buffer00_r,shift_r,nb_elem_r)
508 ELSE
509 buffer00_r(1:size_buffer00_r) = buffer_s(1:size_buffer_s)
510 ENDIF
511 DEALLOCATE( buffer_s )
512 ENDIF
513
514! ----------------------
515! write the received buffer
516! ----------------------
517
518
519
520! received buffer (only on proc 0) :
521!
522! sh | t | s | sh | t | s |... sh | t | s | sh | t | s | sh | t | s |... sh | t | s | ... sh | t | s | sh | t | s | ...
523! \__________/ \__________/ \__________/ \__________/ \__________/ \__________/ ...\__________/ \__________/ ...
524! nsect = 1 2 ... n || 1 2 ... n | || 1 2 ...
525! ||________________________________________||________________________________________| ...||_____________________________
526! | | |
527! proc = 0 1 ... m
528
529
530 IF (ispmd==0.AND.nsect>0) THEN
531 DO i=1,nspmd
532 shift_sect(1,1,i) = 0
533 shift_sect(1,2,i) = rindex_proc(1,1,i)
534 shift_sect(1,3,i) = rindex_proc(1,2,i)
535 DO jj=2,nsect
536 shift_sect(jj,1,i) = rindex_proc(jj-1,3,i)
537 shift_sect(jj,2,i) = rindex_proc(jj,1,i)
538 shift_sect(jj,3,i) = rindex_proc(jj,2,i)
539 ENDDO
540 ENDDO
541
542 DO jj=1,nsect
543 ! ----------------
544 ! SHELL
545 ! ----------------
546 DO i=1,nspmd
547 len = rindex_proc(jj,1,i) - shift_sect(jj,1,i)
548 IF(len>0) THEN
549 indice = 1 + shift_r(i) + shift_sect(jj,1,i)
550 CALL write_i_c(buffer00_r(indice),len)
551 ENDIF
552 ENDDO
553 ! ----------------
554 ! SHELL3N
555 ! ----------------
556 DO i=1,nspmd
557 len = rindex_proc(jj,2,i) - rindex_proc(jj,1,i)
558 IF(len>0) THEN
559 indice = 1 + shift_r(i) + shift_sect(jj,2,i)
560 CALL write_i_c(buffer00_r(indice),len)
561 ENDIF
562 ENDDO
563 ! ----------------
564 ! SOLID
565 ! ----------------
566 DO i=1,nspmd
567 len = rindex_proc(jj,3,i) - rindex_proc(jj,2,i)
568 IF(len>0) THEN
569 indice = 1 + shift_r(i) + shift_sect(jj,3,i)
570 CALL write_i_c(buffer00_r(indice),len)
571 ENDIF
572 ENDDO
573 ENDDO
574 ENDIF
575
576
577
578 IF (ispmd==0) THEN
579 n0 = numnodg + nodcut + nsect
580 n1 = numnodg + nodcut + nsect + nrwall
581
582 DO n=1,nrwall
583 ii(1) = n0
584 ii(2) = n0
585 ii(3) = n0
586 ii(4) = n0
587 CALL write_i_c(ii,4)
588 n0 = n0 + 1
589 k=1
590 n2=n +nrwall
591 n3=n2+nrwall
592 n4=n3+nrwall
593 ityp= nprw(n4)
594
595 IF(iabs(ityp)==1.OR.ityp==4)THEN
596 ii(1) = n1
597 ii(2) = n1 + 3
598 ii(3) = n1 + 2
599 ii(4) = n1 + 1
600 CALL write_i_c(ii,4)
601 n1 = n1 + 4
602 ELSEIF(ityp==2)THEN
603 n10 = n1
604 DO j = 1,23
605 ii(1) = n1
606 ii(2) = n1 + 2
607 ii(3) = n1 + 3
608 ii(4) = n1 + 1
609 CALL write_i_c(ii,4)
610 n1 = n1 + 2
611 ENDDO
612 ii(1) = n1
613 ii(2) = n10
614 ii(3) = n10 + 1
615 ii(4) = n1 + 1
616 CALL write_i_c(ii,4)
617 n1 = n1 + 2
618 ELSEIF(ityp==3)THEN
619 DO i = 1,6
620 DO j = 1,6
621 DO l = 1,6
622 ii(1) = n1
623 ii(2) = n1 + 1
624 ii(3) = n1 + 8
625 ii(4) = n1 + 7
626 CALL write_i_c(ii,4)
627 n1 = n1 + 1
628 ENDDO
629 n1 = n1 + 1
630 ENDDO
631 n1 = n1 + 7
632 ENDDO
633 ENDIF
634 k=k+nprw(n)
635 IF(nprw(n4)==-1)k=k+nint(rwbuf(8,n))
636 ENDDO
637 ENDIF
638
639 IF(ALLOCATED(buffer00_r)) DEALLOCATE( buffer00_r )
640
641 RETURN
#define my_real
Definition cppsort.cpp:32
subroutine spmd_gather_wa(mode, size_buffer_s, size_buffer_r, sindex, rindex_proc, buffer_s, buffer_r, shift_r, nb_elem_r)
void write_i_c(int *w, int *len)