MFIX  2016-1
write_res1_des_mod.f
Go to the documentation of this file.
2 
3  use desmpi
4  use compar, only: mype
5  use compar, only: pe_io
6 
7  use cdist, only: bdist_io
8 
9  use mpi_comm_des, only: des_gather
10  use mpi_comm_des, only: desmpi_gatherv
11 
12  IMPLICIT NONE
13 
14  PRIVATE
15 
16  PUBLIC :: init_write_res_des
17  PUBLIC :: finl_write_res_des
18 
19  PUBLIC :: write_res_des
20  PUBLIC :: write_res_parray
21  PUBLIC :: write_res_carray
22 
23 ! Write scalar and data WITHOUT MPI data collection.
24  INTERFACE write_res_des
25  MODULE PROCEDURE write_res_des_0i, write_res_des_1i
26  MODULE PROCEDURE write_res_des_0d, write_res_des_1d
27  MODULE PROCEDURE write_res_des_0l, write_res_des_1l
28  END INTERFACE
29 
30 ! Write particle array data.
31  INTERFACE write_res_parray
32  MODULE PROCEDURE write_res_parray_1b
33  MODULE PROCEDURE write_res_parray_1i
34  MODULE PROCEDURE write_res_parray_1d
35  MODULE PROCEDURE write_res_parray_1l
36  END INTERFACE
37 
38 ! Write neighbor/collision array data.
39  INTERFACE write_res_carray
40  MODULE PROCEDURE write_res_carray_1i
41  MODULE PROCEDURE write_res_carray_1d
42  MODULE PROCEDURE write_res_carray_1l
43  END INTERFACE
44 
45  INTEGER, PARAMETER :: rdes_unit = 901
46 
47 ! Send/Recv parameters for Particle arrays:
48  INTEGER :: prootcnt, pproccnt
49  INTEGER :: psend
50  INTEGER, allocatable :: pgather(:)
51  INTEGER, allocatable :: pdispls(:)
52 
53 ! Send/Recv parameters for Collision/Neighbor arrays:
54  INTEGER :: crootcnt, cproccnt
55  INTEGER :: csend
56  INTEGER, allocatable :: cgather(:)
57  INTEGER, allocatable :: cdispls(:)
58 
59  CONTAINS
60 
61 !``````````````````````````````````````````````````````````````````````!
62 ! Subroutine: OPEN_RES_DES !
63 ! !
64 ! Purpose: Construct the file name and open the DES RES file. !
65 !``````````````````````````````````````````````````````````````````````!
66  SUBROUTINE open_res_des(BASE)
67 
68  use machine, only: open_n1
69 
70  CHARACTER(len=*), INTENT(IN) :: base
71  CHARACTER(len=32) :: lfname
72 
73  IF(bdist_io) THEN
74  WRITE(lfname,'(A,I4.4,A)') base//'_DES_',mype,'.RES'
75  OPEN(convert='BIG_ENDIAN',unit=rdes_unit, file=lfname, form='UNFORMATTED', &
76  status='UNKNOWN', access='DIRECT', recl=open_n1)
77 
78  ELSEIF(mype == pe_io) THEN
79  WRITE(lfname,'(A,A)') base//'_DES.RES'
80  OPEN(convert='BIG_ENDIAN',unit=rdes_unit, file=lfname, form='UNFORMATTED', &
81  status='UNKNOWN', access='DIRECT', recl=open_n1)
82  ENDIF
83 
84  END SUBROUTINE open_res_des
85 
86 !``````````````````````````````````````````````````````````````````````!
87 ! Subroutine: INIT_WRITE_RES_DES !
88 ! !
89 ! Purpose: Construct the file name and open the DES RES file. !
90 !``````````````````````````````````````````````````````````````````````!
91  SUBROUTINE init_write_res_des(BASE, lVERSION, lNEXT_REC)
92 
93  use compar, only: numpes
94  use mpi_utility, only: global_sum
95 
96  use discretelement, only: pip, ighost_cnt
97  use discretelement, only: neighbors, neighbor_index, neigh_num
98  use functions, only: is_nonexistent
99 
100  CHARACTER(len=*), INTENT(IN) :: BASE
101  DOUBLE PRECISION, INTENT(IN) :: lVERSION
102  INTEGER, INTENT(OUT) :: lNEXT_REC
103 
104 ! Number of real particles on local rank
105  INTEGER :: lPROC
106 ! Total number of real particles.
107  INTEGER :: lGHOST_CNT
108 ! Local gather counts for send/recv
109  INTEGER :: lGatherCnts(0:numpes-1)
110 ! Loop counters
111  INTEGER :: LC1,part
112 
113  CALL open_res_des(base)
114 
115  allocate(pgather(0:numpes-1))
116  allocate(pdispls(0:numpes-1))
117 
118  allocate(cgather(0:numpes-1))
119  allocate(cdispls(0:numpes-1))
120 
121  IF(bdist_io) THEN
122 
123  prootcnt = pip
124  pproccnt = prootcnt
125 
126  lghost_cnt = ighost_cnt
127 
128  crootcnt = neigh_num
129  cproccnt = crootcnt
130  ELSE
131 
132 ! Setup data for particle array data collection:
133  prootcnt = 10
134  pproccnt = pip - ighost_cnt
135 
136 ! Rank 0 gets the total number of gloabl particles.
137  CALL global_sum(pproccnt, prootcnt)
138 
139 ! Serial IO does not store ghost particle data.
140  lghost_cnt = 0
141 
142 ! Construct an array for the Root process that states the number of
143 ! (real) particles on each process.
144  psend = pproccnt
145 
146  lgathercnts = 0
147  lgathercnts(mype) = pproccnt
148 
149  CALL global_sum(lgathercnts, pgather)
150 
151 ! Calculate the displacements for each process in the global array.
152  pdispls(0) = 0
153  DO lproc = 1,numpes-1
154  pdispls(lproc) = pdispls(lproc-1) + pgather(lproc-1)
155  ENDDO
156 
157 ! Setup data for neighbor arrays
158  crootcnt = 10
159 ! Count the number of real neighbors.
160  cproccnt = 0
161  part = 1
162  DO lc1 = 1, neigh_num
163  IF (0 .eq. neighbors(lc1)) EXIT
164  IF (lc1.eq.neighbor_index(part)) THEN
165  part = part + 1
166  ENDIF
167  IF(.NOT.is_nonexistent(part) .AND. .NOT.is_nonexistent(neighbors(lc1))) THEN
168  cproccnt = cproccnt +1
169  ENDIF
170 
171  ENDDO
172 
173 ! Rank 0 gets the total number of global particles.
174  CALL global_sum(cproccnt, crootcnt)
175 
176 ! Construct an array for the Root process that states the number of
177 ! (real) particles on each process.
178  csend = cproccnt
179 
180  lgathercnts = 0
181  lgathercnts(mype) = cproccnt
182 
183  CALL global_sum(lgathercnts, cgather)
184 
185 ! Calculate the displacements for each process in the global array.
186  cdispls(0) = 0
187  DO lproc = 1,numpes-1
188  cdispls(lproc) = cdispls(lproc-1) + cgather(lproc-1)
189  ENDDO
190 
191  ENDIF
192 
193 ! Write out the initial data.
194  lnext_rec = 1
195  CALL write_res_des(lnext_rec, lversion) ! RES file version
196  CALL write_res_des(lnext_rec, prootcnt) ! Number of Particles
197  CALL write_res_des(lnext_rec, lghost_cnt) ! Number of Ghosts
198  CALL write_res_des(lnext_rec, crootcnt) ! Number of neighbors
199 
200  RETURN
201  END SUBROUTINE init_write_res_des
202 
203 !``````````````````````````````````````````````````````````````````````!
204 ! Subroutine: CLOSE_RES_DES !
205 ! !
206 ! Purpose: Close the DES RES file. !
207 !``````````````````````````````````````````````````````````````````````!
208  SUBROUTINE finl_write_res_des
210  IF(bdist_io .OR. mype == pe_io) close(rdes_unit)
211 
212  IF(allocated(dprocbuf)) deallocate(dprocbuf)
213  IF(allocated(drootbuf)) deallocate(drootbuf)
214  IF(allocated(iprocbuf)) deallocate(iprocbuf)
215  IF(allocated(irootbuf)) deallocate(irootbuf)
216 
217  if(allocated(pgather)) deallocate(pgather)
218  if(allocated(pdispls)) deallocate(pdispls)
219 
220  if(allocated(cgather)) deallocate(cgather)
221  if(allocated(cdispls)) deallocate(cdispls)
222 
223  RETURN
224  END SUBROUTINE finl_write_res_des
225 
226 !``````````````````````````````````````````````````````````````````````!
227 ! Subroutine: WRITE_RES_DES_0I !
228 ! !
229 ! Purpose: Write scalar integers to RES file. !
230 !``````````````````````````````````````````````````````````````````````!
231  SUBROUTINE write_res_des_0i(lNEXT_REC, INPUT_I)
232 
233  INTEGER, INTENT(INOUT) :: lNEXT_REC
234  INTEGER, INTENT(IN) :: INPUT_I
235 
236  IF(bdist_io .OR. mype == pe_io) &
237  WRITE(rdes_unit, rec=lnext_rec) input_i
238 
239  lnext_rec = lnext_rec + 1
240 
241  RETURN
242  END SUBROUTINE write_res_des_0i
243 
244 !``````````````````````````````````````````````````````````````````````!
245 ! Subroutine: WRITE_RES_DES_1I !
246 ! !
247 ! Purpose: Write an array of integers to RES file. Note that data is !
248 ! not collected and hsould be on rank 0. !
249 !``````````````````````````````````````````````````````````````````````!
250  SUBROUTINE write_res_des_1i(lNEXT_REC, INPUT_I)
251 
252  INTEGER, INTENT(INOUT) :: lNEXT_REC
253  INTEGER, INTENT(IN) :: INPUT_I(:)
254 
255  INTEGER :: lSIZE
256 
257  lsize = size(input_i)
258 
259  IF(bdist_io .OR. mype == pe_io) &
260  CALL out_bin_512i(rdes_unit, input_i, lsize, lnext_rec)
261 
262  RETURN
263  END SUBROUTINE write_res_des_1i
264 
265 !``````````````````````````````````````````````````````````````````````!
266 ! Subroutine: WRITE_RES_DES_0D !
267 ! !
268 ! Purpose: Write scalar double percision values to RES file. !
269 !``````````````````````````````````````````````````````````````````````!
270  SUBROUTINE write_res_des_0d(lNEXT_REC, INPUT_D)
271 
272  INTEGER, INTENT(INOUT) :: lNEXT_REC
273  DOUBLE PRECISION, INTENT(IN) :: INPUT_D
274 
275  IF(bdist_io .OR. mype == pe_io) &
276  WRITE(rdes_unit, rec=lnext_rec) input_d
277 
278  lnext_rec = lnext_rec + 1
279 
280  RETURN
281  END SUBROUTINE write_res_des_0d
282 
283 !``````````````````````````````````````````````````````````````````````!
284 ! Subroutine: WRITE_RES_DES_1D !
285 ! !
286 ! Purpose: Write an array of double percision values to RES file. Note !
287 ! that data is not collected and should be on rank 0. !
288 !``````````````````````````````````````````````````````````````````````!
289  SUBROUTINE write_res_des_1d(lNEXT_REC, INPUT_D)
290 
291  INTEGER, INTENT(INOUT) :: lNEXT_REC
292  DOUBLE PRECISION, INTENT(IN) :: INPUT_D(:)
293 
294  INTEGER :: lSIZE
295 
296  lsize = size(input_d)
297 
298  IF(bdist_io .OR. mype == pe_io) &
299  CALL out_bin_512(rdes_unit, input_d, lsize, lnext_rec)
300 
301  RETURN
302  END SUBROUTINE write_res_des_1d
303 
304 !``````````````````````````````````````````````````````````````````````!
305 ! Subroutine: WRITE_RES_DES_0L !
306 ! !
307 ! Purpose: Write scalar logical values to RES file. !
308 !``````````````````````````````````````````````````````````````````````!
309  SUBROUTINE write_res_des_0l(lNEXT_REC, INPUT_L)
310 
311  INTEGER, INTENT(INOUT) :: lNEXT_REC
312  LOGICAL, INTENT(IN) :: INPUT_L
313 
314  INTEGER :: INPUT_I
315 
316  input_i = merge(1,0,input_l)
317 
318  IF(bdist_io .OR. mype == pe_io) &
319  WRITE(rdes_unit, rec=lnext_rec) input_i
320 
321  lnext_rec = lnext_rec + 1
322 
323  RETURN
324  END SUBROUTINE write_res_des_0l
325 
326 !``````````````````````````````````````````````````````````````````````!
327 ! Subroutine: WRITE_RES_DES_1L !
328 ! !
329 ! Purpose: Write an array of integers to RES file. Note that data is !
330 ! not collected and hsould be on rank 0. !
331 !``````````````````````````````````````````````````````````````````````!
332  SUBROUTINE write_res_des_1l(lNEXT_REC, INPUT_L)
333 
334  INTEGER, INTENT(INOUT) :: lNEXT_REC
335  LOGICAL, INTENT(IN) :: INPUT_L(:)
336 
337  INTEGER, ALLOCATABLE :: INPUT_I(:)
338 
339  INTEGER :: lSIZE, LC1
340 
341  lsize = size(input_l)
342  ALLOCATE(input_i(lsize))
343 
344  DO lc1=1, lsize
345  input_i(lc1) = merge(1,0,input_l(lc1))
346  ENDDO
347 
348  IF(bdist_io .OR. mype == pe_io) &
349  CALL out_bin_512i(rdes_unit, input_i, lsize, lnext_rec)
350 
351  IF(allocated(input_i)) deallocate(input_i)
352 
353  RETURN
354  END SUBROUTINE write_res_des_1l
355 
356 !``````````````````````````````````````````````````````````````````````!
357 ! Subroutine: WRITE_RES_PARRAY_1B !
358 ! !
359 ! Purpose: Write scalar bytes to RES file. !
360 !``````````````````````````````````````````````````````````````````````!
361  SUBROUTINE write_res_parray_1b(lNEXT_REC, INPUT_B, pLOC2GLB)
362 
363  use desmpi, only: iprocbuf
364  use discretelement, only: max_pip, pip
365  use discretelement, only: iglobal_id
366  use functions, only: is_nonexistent
367 
368  INTEGER, INTENT(INOUT) :: lNEXT_REC
369  INTEGER(KIND=1), INTENT(IN) :: INPUT_B(:)
370  LOGICAL, INTENT(IN), OPTIONAL :: pLOC2GLB
371 
372  INTEGER, ALLOCATABLE :: INPUT_I(:)
373  LOGICAL :: lLOC2GLB
374 ! Loop counters
375  INTEGER :: LC1, LC2
376 
377  lloc2glb = .false.
378  IF(present(ploc2glb)) lloc2glb = ploc2glb
379 
380  allocate(iprocbuf(pproccnt))
381  allocate(irootbuf(prootcnt))
382 
383  allocate(input_i(size(input_b)))
384 
385  input_i(:) = input_b(:)
386 
387  idispls = pdispls
388  igath_sendcnt = psend
389  igathercnts = pgather
390 
391  IF(bdist_io) THEN
392  lc1 = 1
393  IF(lloc2glb) THEN
394  DO lc2 = 1, max_pip
395  IF(lc1 > pip) EXIT
396  IF(is_nonexistent(lc1)) cycle
397  iprocbuf(lc1) = iglobal_id(input_i(lc2))
398  lc1 = lc1 + 1
399  ENDDO
400  ELSE
401  DO lc2 = 1, max_pip
402  IF(lc1 > pip) EXIT
403  IF(is_nonexistent(lc1)) cycle
404  iprocbuf(lc1) = input_i(lc2)
405  lc1 = lc1 + 1
406  ENDDO
407  ENDIF
408  CALL out_bin_512i(rdes_unit, iprocbuf, prootcnt, lnext_rec)
409 
410  ELSE
411  CALL des_gather(input_i, lloc2glb)
412  IF(mype == pe_io) &
413  CALL out_bin_512i(rdes_unit,irootbuf, prootcnt, lnext_rec)
414  ENDIF
415 
416  deallocate(iprocbuf)
417  deallocate(irootbuf)
418  deallocate(input_i)
419 
420  RETURN
421  END SUBROUTINE write_res_parray_1b
422 
423 !``````````````````````````````````````````````````````````````````````!
424 ! Subroutine: WRITE_RES_PARRAY_1I !
425 ! !
426 ! Purpose: Write scalar integers to RES file. !
427 !``````````````````````````````````````````````````````````````````````!
428  SUBROUTINE write_res_parray_1i(lNEXT_REC, INPUT_I, pLOC2GLB)
430  use desmpi, only: iprocbuf
431  use discretelement, only: max_pip, pip
432  use discretelement, only: iglobal_id
433  use functions, only: is_nonexistent
434 
435  INTEGER, INTENT(INOUT) :: lNEXT_REC
436  INTEGER, INTENT(IN) :: INPUT_I(:)
437  LOGICAL, INTENT(IN), OPTIONAL :: pLOC2GLB
438 
439  LOGICAL :: lLOC2GLB
440 ! Loop counters
441  INTEGER :: LC1, LC2
442 
443  lloc2glb = .false.
444  IF(present(ploc2glb)) lloc2glb = ploc2glb
445 
446  allocate(iprocbuf(pproccnt))
447  allocate(irootbuf(prootcnt))
448 
449  idispls = pdispls
450  igath_sendcnt = psend
451  igathercnts = pgather
452 
453  IF(bdist_io) THEN
454  lc1 = 1
455 
456  IF(lloc2glb) THEN
457  DO lc2 = 1, max_pip
458  IF(lc1 > pip) EXIT
459  IF(is_nonexistent(lc1)) cycle
460  iprocbuf(lc1) = iglobal_id(input_i(lc2))
461  lc1 = lc1 + 1
462  ENDDO
463  ELSE
464  DO lc2 = 1, max_pip
465  IF(lc1 > pip) EXIT
466  IF(is_nonexistent(lc1)) cycle
467  iprocbuf(lc1) = input_i(lc2)
468  lc1 = lc1 + 1
469  ENDDO
470  ENDIF
471  CALL out_bin_512i(rdes_unit, iprocbuf, prootcnt, lnext_rec)
472 
473  ELSE
474  CALL des_gather(input_i, lloc2glb)
475  IF(mype == pe_io) &
476  CALL out_bin_512i(rdes_unit,irootbuf, prootcnt, lnext_rec)
477  ENDIF
478 
479  deallocate(iprocbuf)
480  deallocate(irootbuf)
481 
482  RETURN
483  END SUBROUTINE write_res_parray_1i
484 
485 !``````````````````````````````````````````````````````````````````````!
486 ! Subroutine: WRITE_RES_PARRAY_1D !
487 ! !
488 ! Purpose: Write scalar integers to RES file. !
489 !``````````````````````````````````````````````````````````````````````!
490  SUBROUTINE write_res_parray_1d(lNEXT_REC, INPUT_D)
492  use discretelement, only: max_pip, pip
493  use functions, only: is_nonexistent
494 
495  IMPLICIT NONE
496 
497  INTEGER, INTENT(INOUT) :: lNEXT_REC
498  DOUBLE PRECISION, INTENT(IN) :: INPUT_D(:)
499 
500 ! Loop counters
501  INTEGER :: LC1, LC2
502 
503 
504  allocate(dprocbuf(pproccnt))
505  allocate(drootbuf(prootcnt))
506 
507  idispls = pdispls
508  igath_sendcnt = psend
509  igathercnts = pgather
510 
511  IF(bdist_io) THEN
512  lc1 = 1
513  DO lc2 = 1, max_pip
514  IF(lc1 > pip) EXIT
515  IF(is_nonexistent(lc1)) cycle
516  dprocbuf(lc1) = input_d(lc2)
517  lc1 = lc1 + 1
518  ENDDO
519  CALL out_bin_512(rdes_unit, dprocbuf, prootcnt, lnext_rec)
520  ELSE
521  CALL des_gather(input_d)
522  IF(mype == pe_io) &
523  CALL out_bin_512(rdes_unit,drootbuf, prootcnt, lnext_rec)
524  ENDIF
525 
526  deallocate(dprocbuf)
527  deallocate(drootbuf)
528 
529  RETURN
530  END SUBROUTINE write_res_parray_1d
531 
532 !``````````````````````````````````````````````````````````````````````!
533 ! Subroutine: WRITE_RES_PARRAY_1D !
534 ! !
535 ! Purpose: Write scalar integers to RES file. !
536 !``````````````````````````````````````````````````````````````````````!
537  SUBROUTINE write_res_parray_1l(lNEXT_REC, INPUT_L)
539  use desmpi, only: iprocbuf
540  use discretelement, only: max_pip, pip
541  use functions, only: is_nonexistent
542 
543  INTEGER, INTENT(INOUT) :: lNEXT_REC
544  LOGICAL, INTENT(IN) :: INPUT_L(:)
545 
546 ! Loop counters
547  INTEGER :: LC1, LC2
548 
549  allocate(iprocbuf(pproccnt))
550  allocate(irootbuf(prootcnt))
551 
552  idispls = pdispls
553  igath_sendcnt = psend
554  igathercnts = pgather
555 
556  IF(bdist_io) THEN
557  lc1 = 1
558  DO lc2 = 1, max_pip
559  IF(lc1 > pip) EXIT
560  IF(is_nonexistent(lc1)) cycle
561  iprocbuf(lc1) = merge(1,0,input_l(lc2))
562  lc1 = lc1 + 1
563  ENDDO
564  CALL out_bin_512i(rdes_unit, iprocbuf, prootcnt, lnext_rec)
565  ELSE
566  CALL des_gather(input_l)
567  IF(mype == pe_io) &
568  CALL out_bin_512i(rdes_unit,irootbuf, prootcnt, lnext_rec)
569  ENDIF
570 
571  deallocate(iprocbuf)
572  deallocate(irootbuf)
573 
574  RETURN
575  END SUBROUTINE write_res_parray_1l
576 
577 !``````````````````````````````````````````````````````````````````````!
578 ! Subroutine: WRITE_RES_cARRAY_1I !
579 ! !
580 ! Purpose: Write scalar integers to RES file. !
581 !``````````````````````````````````````````````````````````````````````!
582  SUBROUTINE write_res_carray_1i(lNEXT_REC, INPUT_I, pLOC2GLB)
584  use desmpi, only: iprocbuf
585  use discretelement, only: neighbors, neighbor_index, neigh_num
586  USE functions, only: is_nonexistent
587  use discretelement, only: iglobal_id
588 
589  INTEGER, INTENT(INOUT) :: lNEXT_REC
590  INTEGER, INTENT(IN) :: INPUT_I(:)
591  LOGICAL, INTENT(IN), OPTIONAL :: pLOC2GLB
592 
593  LOGICAL :: lLOC2GLB
594 ! Loop counters
595  INTEGER :: LC1, LC2, part
596 
597  lloc2glb = .false.
598  IF(present(ploc2glb)) lloc2glb = ploc2glb
599 
600  allocate(iprocbuf(cproccnt))
601  allocate(irootbuf(crootcnt))
602 
603  idispls = cdispls
604  igath_sendcnt = csend
605  igathercnts = cgather
606 
607  lc2 = 1
608  part = 1
609 
610  DO lc1 = 1, neigh_num
611  IF (0 .eq. neighbors(lc1)) EXIT
612  IF (lc1.eq.neighbor_index(part)) THEN
613  part = part + 1
614  ENDIF
615  IF(.NOT.is_nonexistent(part) .AND. .NOT.is_nonexistent(neighbors(lc1))) THEN
616  IF(lloc2glb) THEN
617  iprocbuf(lc2) = iglobal_id(input_i(lc1))
618  ELSE
619  iprocbuf(lc2) = input_i(lc1)
620  ENDIF
621  lc2 = lc2 + 1
622  ENDIF
623  ENDDO
624 
625  IF(bdist_io) THEN
626  CALL out_bin_512i(rdes_unit, iprocbuf, cproccnt, lnext_rec)
627 
628  ELSE
629  CALL desmpi_gatherv(ptype=1)
630  IF(mype == pe_io) &
631  CALL out_bin_512i(rdes_unit,irootbuf, crootcnt, lnext_rec)
632  ENDIF
633 
634  deallocate(iprocbuf)
635  deallocate(irootbuf)
636 
637  RETURN
638  END SUBROUTINE write_res_carray_1i
639 
640 !``````````````````````````````````````````````````````````````````````!
641 ! Subroutine: WRITE_RES_cARRAY_1D !
642 ! !
643 ! Purpose: Write scalar integers to RES file. !
644 !``````````````````````````````````````````````````````````````````````!
645  SUBROUTINE write_res_carray_1d(lNEXT_REC, INPUT_D)
647  use desmpi, only: dprocbuf ! Local process buffer
648  use desmpi, only: drootbuf ! Root process buffer
649  use discretelement, only: neighbors, neighbor_index, neigh_num
650  USE functions, only: is_nonexistent
651 
652  INTEGER, INTENT(INOUT) :: lNEXT_REC
653  DOUBLE PRECISION, INTENT(IN) :: INPUT_D(:)
654 
655 ! Loop counters
656  INTEGER :: LC1, LC2, part
657 
658  allocate(dprocbuf(cproccnt))
659  allocate(drootbuf(crootcnt))
660 
661  idispls = cdispls
662  igath_sendcnt = csend
663  igathercnts = cgather
664 
665  lc2 = 1
666  part = 1
667  DO lc1 = 1, neigh_num
668  IF (0 .eq. neighbors(lc1)) EXIT
669  IF (lc1.eq.neighbor_index(part)) THEN
670  part = part + 1
671  ENDIF
672  IF(.NOT.is_nonexistent(part) .AND. .NOT.is_nonexistent(neighbors(lc1))) THEN
673  dprocbuf(lc2) = input_d(lc1)
674  lc2 = lc2 + 1
675  ENDIF
676  ENDDO
677 
678  IF(bdist_io) THEN
679  CALL out_bin_512(rdes_unit, dprocbuf, cproccnt, lnext_rec)
680 
681  ELSE
682  CALL desmpi_gatherv(ptype=2)
683  IF(mype == pe_io) &
684  CALL out_bin_512(rdes_unit, drootbuf, crootcnt, lnext_rec)
685  ENDIF
686 
687  deallocate(dprocbuf)
688  deallocate(drootbuf)
689 
690  RETURN
691  END SUBROUTINE write_res_carray_1d
692 
693 
694 !``````````````````````````````````````````````````````````````````````!
695 ! Subroutine: WRITE_RES_cARRAY_1L !
696 ! !
697 ! Purpose: Write scalar integers to RES file. !
698 !``````````````````````````````````````````````````````````````````````!
699  SUBROUTINE write_res_carray_1l(lNEXT_REC, INPUT_L)
701  use desmpi, only: iprocbuf
702  use discretelement, only: neighbors, neighbor_index, neigh_num
703  use functions, only: is_nonexistent
704 
705  INTEGER, INTENT(INOUT) :: lNEXT_REC
706  LOGICAL, INTENT(IN) :: INPUT_L(:)
707 
708 ! Loop counters
709  INTEGER :: LC1, LC2, part
710 
711  allocate(iprocbuf(cproccnt))
712  allocate(irootbuf(crootcnt))
713 
714  idispls = cdispls
715  igath_sendcnt = csend
716  igathercnts = cgather
717 
718 ! Pack the local buffer, skipping data for deleted particles.
719  lc2 = 1
720  part = 1
721  DO lc1 = 1, neigh_num
722  IF (0 .eq. neighbors(lc1)) EXIT
723  IF (lc1.eq.neighbor_index(part)) THEN
724  part = part + 1
725  ENDIF
726  IF(.NOT.is_nonexistent(part) .AND. .NOT.is_nonexistent(neighbors(lc1))) THEN
727  iprocbuf(lc2) = merge(1,0,input_l(lc1))
728  lc2 = lc2 + 1
729  ENDIF
730  ENDDO
731 
732  IF(bdist_io) THEN
733  CALL out_bin_512i(rdes_unit, iprocbuf, cproccnt, lnext_rec)
734 
735  ELSE
736  CALL desmpi_gatherv(ptype=1)
737  IF(mype == pe_io) &
738  CALL out_bin_512i(rdes_unit,irootbuf, crootcnt, lnext_rec)
739  ENDIF
740 
741  deallocate(iprocbuf)
742  deallocate(irootbuf)
743 
744  RETURN
745  END SUBROUTINE write_res_carray_1l
746 
747  END MODULE write_res1_des
integer, dimension(:), allocatable igathercnts
Definition: desmpi_mod.f:48
logical bdist_io
Definition: cdist_mod.f:4
subroutine write_res_carray_1l(lNEXT_REC, INPUT_L)
integer open_n1
Definition: machine_mod.f:5
double precision, dimension(:), allocatable dprocbuf
Definition: desmpi_mod.f:42
subroutine, public init_write_res_des(BASE, lVERSION, lNEXT_REC)
subroutine desmpi_gatherv(ptype, pdebug)
subroutine write_res_parray_1d(lNEXT_REC, INPUT_D)
integer, dimension(:), allocatable irootbuf
Definition: desmpi_mod.f:43
integer numpes
Definition: compar_mod.f:24
integer pe_io
Definition: compar_mod.f:30
subroutine write_res_parray_1i(lNEXT_REC, INPUT_I, pLOC2GLB)
subroutine write_res_parray_1l(lNEXT_REC, INPUT_L)
Definition: cdist_mod.f:2
integer, dimension(:), allocatable iprocbuf
Definition: desmpi_mod.f:44
subroutine out_bin_512(IUNIT, ARRAY, N, NEXT_REC)
Definition: out_bin_512.f:24
subroutine write_res_carray_1d(lNEXT_REC, INPUT_D)
subroutine write_res_carray_1i(lNEXT_REC, INPUT_I, pLOC2GLB)
integer mype
Definition: compar_mod.f:24
integer, dimension(:), allocatable idispls
Definition: desmpi_mod.f:46
integer igath_sendcnt
Definition: desmpi_mod.f:51
double precision, dimension(:), allocatable drootbuf
Definition: desmpi_mod.f:41
subroutine out_bin_512i(IUNIT, ARRAY, N, NEXT_REC)
Definition: out_bin_512i.f:24
subroutine, public finl_write_res_des