KTH framework for Nek5000 toolboxes; testing version  0.0.1
stat_IO.f
Go to the documentation of this file.
1 
9 !=======================================================================
12  subroutine stat_mfo()
13  implicit none
14 
15  include 'SIZE'
16  include 'STATD'
17 
18  if (stat_rdim.eq.1) then
19  call stat_mfo_outfld2d()
20  else
21  call stat_mfo_outfld3d()
22  endif
23 
24  return
25  end subroutine
26 !=======================================================================
29  subroutine stat_mfo_outfld3d()
30  implicit none
31 
32  include 'SIZE'
33  include 'INPUT' ! ifpo
34  include 'SOLN' ! pr
35  include 'STATD' !
36 
37  ! local variables
38  logical ifpo_tmp, ifto_tmp, ifxyo_tmp
39 !----------------------------------------------------------------------
40  ifpo_tmp = ifpo
41  ifto_tmp = ifto
42  ifxyo_tmp = ifxyo
43  ifpo=.false.
44  ifto =.true.
45  ifxyo = .true.
46 
47  ! Fields to outpost: <u>t, <v>t, <w>t, <p>t
48  call outpost(stat_ruavg(1,1,1),stat_ruavg(1,1,2),
49  $ stat_ruavg(1,1,3),pr,stat_ruavg(1,1,4),'s01')
50 
51  ! Fields to outpost: <uu>t, <vv>t, <ww>t, <pp>t
52  call outpost(stat_ruavg(1,1,5),stat_ruavg(1,1,6),
53  $ stat_ruavg(1,1,7),pr,stat_ruavg(1,1,8),'s02')
54 
55  ! Fields to outpost: <uv>t, <vw>t,<uw>t, <pu>t
56  call outpost(stat_ruavg(1,1,9),stat_ruavg(1,1,10),
57  $ stat_ruavg(1,1,11),pr,stat_ruavg(1,1,12),'s03')
58 
59  ! Fields to outpost: <pv>t, <pw>t, <pdudx>t, <pdudy>t
60  call outpost(stat_ruavg(1,1,13),stat_ruavg(1,1,14),
61  $ stat_ruavg(1,1,15),pr,stat_ruavg(1,1,16),'s04')
62 
63  ! Fields to outpost: <pdudz>t, <pdvdx>t, <pdvdy>t, <pdvdz>t
64  call outpost(stat_ruavg(1,1,17),stat_ruavg(1,1,18),
65  $ stat_ruavg(1,1,19),pr,stat_ruavg(1,1,20),'s05')
66 
67  !Fields to outpost: <pdwdx>t, <pdwdy>t, <pdwdz>t, <uuu>t
68  call outpost(stat_ruavg(1,1,21),stat_ruavg(1,1,22),
69  $ stat_ruavg(1,1,23),pr,stat_ruavg(1,1,24),'s06')
70 
71  ! Fields to outpost: <vvv>t, <www>t, <uuv>t, <uuw>t
72  call outpost(stat_ruavg(1,1,25),stat_ruavg(1,1,26),
73  $ stat_ruavg(1,1,27),pr,stat_ruavg(1,1,28),'s07')
74 
75  ! Fields to outpost: <vvu>t, <vvw>t, <wwu>t, <wwv>t
76  call outpost(stat_ruavg(1,1,29),stat_ruavg(1,1,30),
77  $ stat_ruavg(1,1,31),pr,stat_ruavg(1,1,32),'s08')
78 
79  ! Fields to outpost: <ppp>t, <pppp>t, <uvw>t, <uuuu>t
80  call outpost(stat_ruavg(1,1,33),stat_ruavg(1,1,34),
81  $ stat_ruavg(1,1,35),pr,stat_ruavg(1,1,36),'s09')
82 
83  ! Fields to outpost: <vvvv>t, <wwww>t, <e11>t, <e22>t
84  call outpost(stat_ruavg(1,1,37),stat_ruavg(1,1,38),
85  $ stat_ruavg(1,1,39),pr,stat_ruavg(1,1,40),'s10')
86 
87  ! Fields to outpost: <e33>t, <e12>t, <e13>t, <e23>t
88  call outpost(stat_ruavg(1,1,41),stat_ruavg(1,1,42),
89  $ stat_ruavg(1,1,43),pr,stat_ruavg(1,1,44),'s11')
90 
91  ifpo=ifpo_tmp
92  ifto=ifto_tmp
93  ifxyo=ifxyo_tmp
94 
95  return
96  end subroutine
97 !=======================================================================
102  subroutine stat_mfo_outfld2d()
103  implicit none
104 
105  include 'SIZE'
106  include 'RESTART'
107  include 'PARALLEL'
108  include 'INPUT'
109  include 'TSTEP'
110  include 'MAP2D'
111  include 'STATD'
112 
113  ! global variablse
114  ! dummy arrays
115  real ur1(lx1,lz1,2*LELT)
116  common /scruz/ ur1
117 
118  ! local variables
119  ! temporary variables to overwrite global values
120  logical ifreguol ! uniform mesh
121  logical ifxyol ! write down mesh
122  integer wdsizol ! store global wdsizo
123  integer nel2DB ! running sum for owned 2D elements
124  integer nelBl ! store global nelB
125 
126  integer il, jl, kl ! loop index
127  integer itmp ! dummy integer
128  integer ierr ! error mark
129  integer nxyzo ! element size
130 
131  character*3 prefix ! file prefix
132 
133  integer*8 offs0, offs ! offset
134  integer*8 stride,strideB ! stride
135 
136  integer ioflds ! fields count
137 
138  real dnbyte ! byte sum
139  real tiostart, tio ! simple timing
140 
141  ! functions
142  integer igl_running_sum
143  real dnekclock_sync, glsum
144 !----------------------------------------------------------------------
145  ! simple timing
146  tiostart=dnekclock_sync()
147 
148  ! intialise I/O
149  ifdiro = .false.
150 
151  ifmpiio = .false.
152  if(abs(param(65)).eq.1 .and. abs(param(66)).eq.6) ifmpiio=.true.
153 #ifdef NOMPIIO
154  ifmpiio = .false.
155 #endif
156 
157  if(ifmpiio) then
158  nfileo = np
159  nproc_o = 1
160  fid0 = 0
161  pid0 = nid
162  pid1 = 0
163  else
164  if(param(65).lt.0) ifdiro = .true. ! p65 < 0 --> multi subdirectories
165  nfileo = abs(param(65))
166  if(nfileo.eq.0) nfileo = 1
167  if(np.lt.nfileo) nfileo=np
168  nproc_o = np / nfileo ! # processors pointing to pid0
169  fid0 = nid/nproc_o ! file id
170  pid0 = nproc_o*fid0 ! my parent i/o node
171  pid1 = min(np-1,pid0+nproc_o-1) ! range of sending procs
172  endif
173 
174  ! save and set global IO variables
175  ! no uniform mesh
176  ifreguol = ifreguo
177  ifreguo = .false.
178 
179  ! save mesh
180  ifxyol = ifxyo
181  ifxyo = .true.
182 
183  ! force double precission
184  wdsizol = wdsizo
185  ! for testing
186  wdsizo = wdsize
187 
188  nrg = lxo
189 
190  ! get number of 2D elements owned by proceesor with smaller nid
191  itmp = map2d_lown
192  nel2db = igl_running_sum(itmp)
193  nel2db = nel2db - map2d_lown
194  ! replace value
195  nelbl = nelb
196  nelb = nel2db
197 
198  ! set element size
199  nxo = stat_nm2
200  nyo = stat_nm3
201  nzo = 1
202  nxyzo = nxo*nyo*nzo
203 
204  ! Save info regarding element ordering; especially important for AMR run
205  call stat_mfo_crd2d
206 
207  ! open files on i/o nodes
208  prefix='sts'
209  ierr=0
210  if (nid.eq.pid0) call mfo_open_files(prefix,ierr)
211 
212  call mntr_check_abort(stat_id,ierr,
213  $ 'Error opening file in stat_mfo_outfld2D.')
214 
215  ! write header, byte key, global ordering
217 
218  ! initial offset: header, test pattern, global ordering
219  offs0 = iheadersize + 4 + isize*int(map2d_gnum,8)
220  offs = offs0
221 
222  ! stride
223  strideb = int(nelb,8)*nxyzo*wdsizo
224  stride = int(map2d_gnum,8)*nxyzo*wdsizo
225 
226  ! count fields
227  ioflds = 0
228 
229  ! write coordinates
230  kl = 0
231  ! copy vector
232  do il=1,map2d_lnum
233  if(map2d_own(il).eq.nid) then
234  call copy(ur1(1,1,2*kl+1),map2d_xm1(1,1,il),nxyzo)
235  call copy(ur1(1,1,2*kl+2),map2d_ym1(1,1,il),nxyzo)
236  kl = kl +1
237  endif
238  enddo
239 
240  ! check consistency
241  ierr = 0
242  if (kl.ne.map2d_lown) ierr=1
243  call mntr_check_abort(stat_id,ierr,'inconsistent map2d_lown 1')
244 
245  ! offset
246  kl = 2*kl
247  offs = offs0 + stride*ioflds + 2*strideb
248  call byte_set_view(offs,ifh_mbyte)
249  call mfo_outs(ur1,kl,nxo,nyo,nzo)
250  ioflds = ioflds + 2
251 
252  ! write fields
253  do jl=1,stat_nvar
254  kl = 0
255  ! copy vector
256  do il=1,map2d_lnum
257  if(map2d_own(il).eq.nid) then
258  kl = kl +1
259  call copy(ur1(1,1,kl),stat_ruavg(1,il,jl),nxyzo)
260  endif
261  enddo
262 
263  ! check consistency
264  ierr = 0
265  if (kl.ne.map2d_lown) ierr=1
266  call mntr_check_abort(stat_id,ierr,'inconsistent map2d_lown 2')
267 
268  ! offset
269  offs = offs0 + stride*ioflds + strideb
270  call byte_set_view(offs,ifh_mbyte)
271  call mfo_outs(ur1,kl,nxo,nyo,nzo)
272  ioflds = ioflds + 1
273  enddo
274 
275  ! write averaging data
277 
278  ! count bytes
279  dnbyte = 1.*ioflds*map2d_lown*wdsizo*nxyzo
280 
281  ierr = 0
282  if (nid.eq.pid0) then
283  if(ifmpiio) then
284  call byte_close_mpi(ifh_mbyte,ierr)
285  else
286  call byte_close(ierr)
287  endif
288  endif
289  call mntr_check_abort(stat_id,ierr,
290  $ 'Error closing file in stat_mfo_outfld2D.')
291 
292  tio = dnekclock_sync()-tiostart
293  if (tio.le.0) tio=1.
294 
295  dnbyte = glsum(dnbyte,1)
296  dnbyte = dnbyte + iheadersize + 4 + isize*map2d_gnum
297  dnbyte = dnbyte/1024/1024
298  if(nio.eq.0) write(6,7) istep,time,dnbyte,dnbyte/tio,
299  & nfileo
300  7 format(/,i9,1pe12.4,' done :: Write checkpoint',/,
301  & 30x,'file size = ',3pg12.2,'MB',/,
302  & 30x,'avg data-throughput = ',0pf7.1,'MB/s',/,
303  & 30x,'io-nodes = ',i5,/)
304 
305  ! set global IO variables back
306  ifreguo = ifreguol
307  ifxyo = ifxyol
308  wdsizo = wdsizol
309  nelb = nelbl
310 
311  return
312  end subroutine
313 !=======================================================================
319  implicit none
320 
321  include 'SIZE'
322  include 'INPUT'
323  include 'RESTART'
324  include 'PARALLEL'
325  include 'TSTEP'
326  include 'SOLN'
327  include 'MAP2D'
328  include 'STATD'
329 
330  ! global variable
331  integer lglist(0:LELT) ! dummy array
332  common /ctmp0/ lglist
333 
334  ! local variables
335  real*4 test_pattern ! byte key
336  integer idum, inelp
337  integer nelo ! number of elements to write
338  integer nfileoo ! number of files to create
339 
340  integer il, jl, kl ! loop index
341  integer mtype ! tag
342 
343  integer ierr ! error mark
344  integer ibsw_out, len
345  integer*8 ioff ! offset
346  logical if_press_mesh ! pessure mesh mark
347 
348  character*132 hdr ! header
349 !-----------------------------------------------------------------------
350  if(ifmpiio) then
351  nfileoo = 1 ! all data into one file
352  nelo = map2d_gnum
353  else
354  nfileoo = nfileo
355  if(nid.eq.pid0) then ! how many elements to dump
356  nelo = map2d_lown
357  do jl = pid0+1,pid1
358  mtype = jl
359  call csend(mtype,idum,isize,jl,0) ! handshake
360  call crecv(mtype,inelp,isize)
361  nelo = nelo + inelp
362  enddo
363  else
364  mtype = nid
365  call crecv(mtype,idum,isize) ! hand-shake
366  call csend(mtype,map2d_lown,isize,pid0,0) ! u4 :=: u8
367  endif
368  endif
369 
370  ! write header
371  ierr = 0
372  if(nid.eq.pid0) then
373  call blank(hdr,132)
374 
375  ! varialbe set
376  call blank(rdcode1,10)
377 
378  ! we save coordinates
379  rdcode1(1)='X'
380  ! and set of fields marked as passive scalars
381  rdcode1(2) = 'S'
382  write(rdcode1(3),'(i1)') stat_nvar/10
383  write(rdcode1(4),'(i1)') stat_nvar-(stat_nvar/10)*10
384 
385  ! no pressure written so pressure format set to false
386  if_press_mesh = .false.
387 
388  write(hdr,1) wdsizo,nxo,nyo,nzo,nelo,map2d_gnum,time,istep,
389  $ fid0, nfileoo, (rdcode1(il),il=1,10),p0th,if_press_mesh
390  1 format('#std',1x,i1,1x,i2,1x,i2,1x,i2,1x,i10,1x,i10,1x,
391  $ e20.13,1x,i9,1x,i6,1x,i6,1x,10a,1pe15.7,1x,l1)
392 
393  ! write test pattern for byte swap
394  test_pattern = 6.54321
395 
396  if(ifmpiio) then
397  ! only rank0 (pid00) will write hdr + test_pattern + time list
398  call byte_write_mpi(hdr,iheadersize/4,pid00,ifh_mbyte,ierr)
399  call byte_write_mpi(test_pattern,1,pid00,ifh_mbyte,ierr)
400  else
401  call byte_write(hdr,iheadersize/4,ierr)
402  call byte_write(test_pattern,1,ierr)
403  endif
404 
405  endif
406 
407  call mntr_check_abort(stat_id,ierr,
408  $ 'Error writing header in stat_mfo_write_hdr2D.')
409 
410  ! write global 2D elements numbering for this group
411  ! copy data
412  lglist(0) = map2d_lown
413  kl = 0
414  do il=1,map2d_lnum
415  if(map2d_own(il).eq.nid) then
416  kl = kl +1
417  lglist(kl) = map2d_gmap(il)
418  endif
419  enddo
420  ! check consistency
421  ierr = 0
422  if (kl.ne.map2d_lown) ierr=1
423 
424  call mntr_check_abort(stat_id,ierr,'inconsistent map2d_lown 3')
425 
426  if(nid.eq.pid0) then
427  if(ifmpiio) then
428  ioff = iheadersize + 4 + int(nelb,8)*isize
429  call byte_set_view (ioff,ifh_mbyte)
430  call byte_write_mpi (lglist(1),lglist(0),-1,ifh_mbyte,ierr)
431  else
432  call byte_write(lglist(1),lglist(0),ierr)
433  endif
434 
435  do jl = pid0+1,pid1
436  mtype = jl
437  call csend(mtype,idum,isize,jl,0) ! handshake
438  len = isize*(lelt+1)
439  call crecv(mtype,lglist,len)
440  if(ierr.eq.0) then
441  if(ifmpiio) then
442  call byte_write_mpi
443  $ (lglist(1),lglist(0),-1,ifh_mbyte,ierr)
444  else
445  call byte_write(lglist(1),lglist(0),ierr)
446  endif
447  endif
448  enddo
449  else
450  mtype = nid
451  call crecv(mtype,idum,isize) ! hand-shake
452 
453  len = isize*(map2d_lown+1)
454  call csend(mtype,lglist,len,pid0,0)
455  endif
456 
457  call mntr_check_abort(stat_id,ierr,
458  $ 'Error writing global nums in stat_mfo_write_hdr2D')
459 
460  return
461  end subroutine
462 !=======================================================================
466  implicit none
467 
468  include 'SIZE'
469  include 'RESTART'
470  include 'MAP2D'
471  include 'STATD'
472 !-----------------------------------------------------------------------
473  ! to be added
474  if(nid.eq.pid0) then
475 
476 
477  endif
478 
479  return
480  end subroutine
481 !=======================================================================
484  subroutine stat_mfo_crd2d
485  implicit none
486 
487  include 'SIZE'
488  include 'INPUT'
489  include 'RESTART'
490  include 'PARALLEL'
491  include 'TSTEP'
492  include 'WZ'
493  include 'MAP2D'
494  include 'STATD'
495 #ifdef AMR
496  include 'AMR'
497 #endif
498 
499  ! local variables
500  character*3 prefix ! file prefix
501 
502  real*4 test_pattern ! byte key
503  integer lglist(0:LELT), itmp(LELT) ! dummy array
504  integer nxl ! number of grid points for interpolation operator
505  parameter(nxl=3)
506  real zgml(nxl), wgtl(nxl) ! gll points and weights for interpolation mesh
507  real iresl(nxl,lx1), itresl(lx1,nxl) ! interpolation operators
508  real rtmp1(nxl,lx1), rtmp2(nxl,nxl) ! temporary arrays for interpolation
509  real rglist(2,0:LELT)
510  integer noutl ! number of bytes to write
511  integer idum, inelp
512  integer nelo ! number of elements to write
513  integer nfileoo ! number of files to create
514 
515  integer il, jl, kl ! loop index
516  integer mtype ! tag
517 
518  integer ierr ! error mark
519  integer ibsw_out, len
520  integer*8 ioff ! offset
521 
522  character*132 hdr ! header
523 !-----------------------------------------------------------------------
524  ! open files on i/o nodes
525  prefix='c2D'
526  ierr=0
527  if (nid.eq.pid0) call mfo_open_files(prefix,ierr)
528 
529  call mntr_check_abort(stat_id,ierr,
530  $ 'Error opening file in stat_mfo_AMRd2D.')
531 
532  ! master-slave communication
533  if(ifmpiio) then
534  nfileoo = 1 ! all data into one file
535  nelo = map2d_gnum
536  else
537  nfileoo = nfileo
538  if(nid.eq.pid0) then ! how many elements to dump
539  nelo = map2d_lown
540  do jl = pid0+1,pid1
541  mtype = jl
542  call csend(mtype,idum,isize,jl,0) ! handshake
543  call crecv(mtype,inelp,isize)
544  nelo = nelo + inelp
545  enddo
546  else
547  mtype = nid
548  call crecv(mtype,idum,isize) ! hand-shake
549  call csend(mtype,map2d_lown,isize,pid0,0) ! u4 :=: u8
550  endif
551  endif
552 
553  ! write header
554  ierr = 0
555  if(nid.eq.pid0) then
556  call blank(hdr,132)
557 
558  call blank(rdcode1,10)
559 
560  write(hdr,1) wdsizo,nelo,map2d_gnum,time,istep,
561  $ fid0, nfileoo
562  1 format('#amr',1x,i2,1x,i10,1x,i10,1x,
563  $ e20.13,1x,i9,1x,i6,1x,i6)
564 
565  ! write test pattern for byte swap
566  test_pattern = 6.54321
567 
568  if(ifmpiio) then
569  ! only rank0 (pid00) will write hdr + test_pattern + time list
570  call byte_write_mpi(hdr,iheadersize/4,pid00,ifh_mbyte,ierr)
571  call byte_write_mpi(test_pattern,1,pid00,ifh_mbyte,ierr)
572  else
573  call byte_write(hdr,iheadersize/4,ierr)
574  call byte_write(test_pattern,1,ierr)
575  endif
576  endif
577 
578  call mntr_check_abort(stat_id,ierr,
579  $ 'Error writing header in stat_mfo_crd2D.')
580 
581  ! write global 2D elements numbering for this group
582  ! copy data
583  lglist(0) = map2d_lown
584  kl = 0
585  do il=1,map2d_lnum
586  if(map2d_own(il).eq.nid) then
587  kl = kl +1
588  lglist(kl) = map2d_gmap(il)
589  endif
590  enddo
591  ! check consistency
592  ierr = 0
593  if (kl.ne.map2d_lown) ierr=1
594 
595  call mntr_check_abort(stat_id,ierr,'inconsistent map2d_lown 4')
596 
597  if(nid.eq.pid0) then
598  noutl = lglist(0)*isize/4
599  if(ifmpiio) then
600  ioff = iheadersize + 4 + int(nelb,8)*isize
601  call byte_set_view (ioff,ifh_mbyte)
602  call byte_write_mpi (lglist(1),noutl,-1,ifh_mbyte,ierr)
603  else
604  call byte_write(lglist(1),noutl,ierr)
605  endif
606 
607  do jl = pid0+1,pid1
608  mtype = jl
609  call csend(mtype,idum,isize,jl,0) ! handshake
610  len = isize*(lelt+1)
611  call crecv(mtype,lglist,len)
612  if(ierr.eq.0) then
613  noutl = lglist(0)*isize/4
614  if(ifmpiio) then
615  call byte_write_mpi(lglist(1),noutl,-1,ifh_mbyte,ierr)
616  else
617  call byte_write(lglist(1),noutl,ierr)
618  endif
619  endif
620  enddo
621  else
622  mtype = nid
623  call crecv(mtype,idum,isize) ! hand-shake
624 
625  len = isize*(map2d_lown+1)
626  call csend(mtype,lglist,len,pid0,0)
627  endif
628 
629  call mntr_check_abort(stat_id,ierr,
630  $ 'Error writing global nums in stat_mfo_crd2D')
631 
632 #ifdef AMR
633  ! write refinement level; this code is not optimal, but I don't do it often
634  do il = 1,nelv
635  itmp(map2d_lmap(il)) = amr_level(il)
636  enddo
637 #else
638  do il = 1,nelt
639  itmp(il) = 0
640  enddo
641 #endif
642  ! copy data
643  lglist(0) = map2d_lown
644  kl = 0
645  do il=1,map2d_lnum
646  if(map2d_own(il).eq.nid) then
647  kl = kl +1
648  lglist(kl) = itmp(il)
649  endif
650  enddo
651 
652  if(nid.eq.pid0) then
653  noutl = lglist(0)*isize/4
654  if(ifmpiio) then
655  ioff = iheadersize + 4 + (int(map2d_gnum,8) + int(nelb,8))
656  $ *isize
657  call byte_set_view (ioff,ifh_mbyte)
658  call byte_write_mpi (lglist(1),noutl,-1,ifh_mbyte,ierr)
659  else
660  call byte_write(lglist(1),noutl,ierr)
661  endif
662 
663  do jl = pid0+1,pid1
664  mtype = jl
665  call csend(mtype,idum,isize,jl,0) ! handshake
666  len = isize*(lelt+1)
667  call crecv(mtype,lglist,len)
668  if(ierr.eq.0) then
669  noutl = lglist(0)*isize/4
670  if(ifmpiio) then
671  call byte_write_mpi(lglist(1),noutl,-1,ifh_mbyte,ierr)
672  else
673  call byte_write(lglist(1),noutl,ierr)
674  endif
675  endif
676  enddo
677  else
678  mtype = nid
679  call crecv(mtype,idum,isize) ! hand-shake
680 
681  len = isize*(map2d_lown+1)
682  call csend(mtype,lglist,len,pid0,0)
683  endif
684 
685  call mntr_check_abort(stat_id,ierr,
686  $ 'Error writing refinement level in stat_mfo_crd2D')
687 
688  ! get 2D elements cell centres
689  ! get interpolation operators
690  call zwgll(zgml,wgtl,nxl)
691  call igllm(iresl,itresl,zgm1,zgml,lx1,nxl,lx1,nxl)
692 
693  ! interpolate mesh and extract element centre
694  rglist(1,0) = real(map2d_lown)
695  kl = 0
696  do il=1,map2d_lnum
697  if(map2d_own(il).eq.nid) then
698  kl = kl +1
699  call mxm(iresl,nxl,map2d_xm1(1,1,il),lx1,rtmp1,lx1)
700  call mxm (rtmp1,nxl,itresl,lx1,rtmp2,nxl)
701  rglist(1,kl) = rtmp2(2,2)
702  call mxm(iresl,nxl,map2d_ym1(1,1,il),lx1,rtmp1,lx1)
703  call mxm (rtmp1,nxl,itresl,lx1,rtmp2,nxl)
704  rglist(2,kl) = rtmp2(2,2)
705  endif
706  enddo
707 
708  ! write down cell centres
709  if(nid.eq.pid0) then
710  noutl = 2*int(rglist(1,0))*wdsizo/4
711  if(ifmpiio) then
712  ioff = iheadersize + 4 + int(map2d_gnum,8)*2*isize
713  $ + int(nelb,8)*2*wdsizo
714  call byte_set_view (ioff,ifh_mbyte)
715  call byte_write_mpi (rglist(1,1),noutl,-1,ifh_mbyte,ierr)
716  else
717  call byte_write(rglist(1,1),noutl,ierr)
718  endif
719 
720  do jl = pid0+1,pid1
721  mtype = jl
722  call csend(mtype,idum,isize,jl,0) ! handshake
723  len = 2*wdsize*(lelt+1)
724  call crecv(mtype,rglist,len)
725  if(ierr.eq.0) then
726  noutl = 2*int(rglist(1,0))*wdsizo/4
727  if(ifmpiio) then
728  call byte_write_mpi
729  $ (rglist(1,1),noutl,-1,ifh_mbyte,ierr)
730  else
731  call byte_write(rglist(1,1),noutl,ierr)
732  endif
733  endif
734  enddo
735  else
736  mtype = nid
737  call crecv(mtype,idum,isize) ! hand-shake
738 
739  len = 2*wdsize*(map2d_lown+1)
740  call csend(mtype,rglist,len,pid0,0)
741  endif
742 
743  call mntr_check_abort(stat_id,ierr,
744  $ 'Error writing element centres in stat_mfo_crd2D')
745 
746  ierr = 0
747  if (nid.eq.pid0) then
748  if(ifmpiio) then
749  call byte_close_mpi(ifh_mbyte,ierr)
750  else
751  call byte_close(ierr)
752  endif
753  endif
754  call mntr_check_abort(stat_id,ierr,
755  $ 'Error closing file in stat_mfo_crd2D.')
756 
757  return
758  end subroutine
759 !=======================================================================
760 
#define byte_write
Definition: byte.c:39
#define byte_close
Definition: byte.c:36
subroutine byte_write_mpi(buf, icount, iorank, mpi_fh, ierr)
Definition: byte_mpi.f:65
subroutine byte_close_mpi(mpi_fh, ierr)
Definition: byte_mpi.f:84
subroutine byte_set_view(ioff_in, mpi_fh)
Definition: byte_mpi.f:97
subroutine crecv(mtype, buf, lenm)
Definition: comm_mpi.f:313
subroutine csend(mtype, buf, len, jnid, jpid)
Definition: comm_mpi.f:303
subroutine mntr_check_abort(mid, ierr, logs)
Abort simulation.
Definition: mntrlog.f:856
subroutine stat_mfo_write_stat2d
Write additional data at the end of the file.
Definition: stat_IO.f:466
subroutine stat_mfo_write_hdr2d
Write hdr, byte key, global ordering of 2D data.
Definition: stat_IO.f:319
subroutine stat_mfo()
Main interface for saving statistics.
Definition: stat_IO.f:13
subroutine stat_mfo_outfld3d()
Statistics muti-file output of 3D data.
Definition: stat_IO.f:30
subroutine stat_mfo_crd2d
Write element centres to the file (in case of AMR level as well)
Definition: stat_IO.f:485
subroutine stat_mfo_outfld2d()
Statistics muti-file output of 2D data.
Definition: stat_IO.f:103
subroutine copy(a, b, n)
Definition: math.f:260
subroutine blank(A, N)
Definition: math.f:19
subroutine mxm(a, n1, b, n2, c, n3)
Definition: mxm_wrapper.f:2
subroutine mfo_open_files(prefix, ierr)
Definition: prepost.f:1071
subroutine mfo_outs(u, nel, mx, my, mz)
Definition: prepost.f:1633
subroutine outpost(v1, v2, v3, vp, vt, name3)
Definition: prepost.f:1378
subroutine zwgll(Z, W, NP)
Definition: speclib.f:108
subroutine igllm(I12, IT12, Z1, Z2, lz1, lz2, ND1, ND2)
Definition: speclib.f:1102