KTH framework for Nek5000 toolboxes; testing version  0.0.1
io_tools.f
Go to the documentation of this file.
1 
6 !=======================================================================
10  subroutine io_register()
11  implicit none
12 
13  include 'FRAMELP'
14  include 'IOTOOLD'
15 
16  ! local variables
17  integer lpmid
18 !-----------------------------------------------------------------------
19  ! check if the current module was already registered
20  call mntr_mod_is_name_reg(lpmid,io_name)
21  if (lpmid.gt.0) then
22  call mntr_warn(lpmid,
23  $ 'module ['//trim(io_name)//'] already registered')
24  return
25  endif
26 
27  ! find parent module
28  call mntr_mod_is_name_reg(lpmid,'FRAME')
29  if (lpmid.le.0) then
30  lpmid = 1
31  call mntr_abort(lpmid,
32  $ 'Parent module ['//'FRAME'//'] not registered')
33  endif
34 
35  ! register module
36  call mntr_mod_reg(io_id,lpmid,io_name,'I/O TOOLS')
37 
38  return
39  end subroutine
40 !=======================================================================
46  subroutine io_file_freeid(iunit, ierr)
47  implicit none
48 
49  include 'FRAMELP'
50  include 'IOTOOLD'
51 
52  ! argument list
53  integer iunit
54  integer ierr
55 
56  ! local variables
57  logical ifcnnd ! is unit connected
58 !-----------------------------------------------------------------------
59  ! initialise variables
60  ierr=0
61  iunit = io_iunit_min
62 
63  do
64  inquire(unit=iunit,opened=ifcnnd,iostat=ierr)
65  if(ifcnnd) then
66  iunit = iunit +1
67  else
68  exit
69  endif
70  enddo
71 
72  if (iunit.gt.io_iunit_max) io_iunit_max = iunit
73 
74  return
75  end subroutine
76 !=======================================================================
80  subroutine io_file_close()
81  implicit none
82 
83  include 'FRAMELP'
84  include 'IOTOOLD'
85 
86  ! local variables
87  integer iunit, ierr
88  logical ifcnnd ! is unit connected
89 !-----------------------------------------------------------------------
90  do iunit = io_iunit_min, io_iunit_max
91  inquire(unit=iunit,opened=ifcnnd,iostat=ierr)
92  if(ifcnnd) close(iunit)
93  enddo
94  io_iunit_max = io_iunit_min
95 
96  return
97  end subroutine
98 !=======================================================================
108  subroutine io_mfo_fname(fname,bname,prefix,ierr)
109  implicit none
110 
111  include 'SIZE'
112  include 'INPUT' ! IFREGUO, IFMPIIO
113  include 'RESTART' ! NFILEO
114  include 'FRAMELP'
115  include 'IOTOOLD'
116 
117 
118  ! argument list
119  character*132 fname, bname
120  character*3 prefix
121  integer ierr
122 
123  ! local variables
124  integer ndigit, itmp
125  real rfileo
126 
127  character*(*) six
128  parameter(six='??????')
129 !-----------------------------------------------------------------------
130  ! initialise variables
131  ierr = 0
132  fname = ''
133 
134  ! numbe or IO nodes
135  if (ifmpiio) then
136  rfileo = 1
137  else
138  rfileo = nfileo
139  endif
140  ndigit = log10(rfileo) + 1
141 
142  ! Add directory
143  if (ifdiro) fname = 'A'//six(1:ndigit)//'/'
144 
145  ! Add prefix
146  if (prefix(1:1).ne.' '.and.prefix(2:2).ne.' '
147  $ .and.prefix(3:3).ne.' ')
148  $ fname = trim(fname)//trim(adjustl(prefix))
149 
150  ! Add SESSION
151  fname = trim(fname)//trim(adjustl(bname))
152 
153  if (ifreguo) fname = trim(fname)//'_reg'
154 
155  ! test string length
156  itmp = len_trim(fname)
157  if (itmp.eq.0) then
158  call mntr_error(io_id,'io_mfo_fname; zero lenght fname.')
159  ierr = 1
160  return
161  elseif ((itmp+ndigit+2+5).gt.132) then
162  call mntr_error(io_id,'io_mfo_fname; fname too long.')
163  ierr = 2
164  return
165  endif
166 
167  ! Add file-id holder and .f appendix
168  fname = trim(fname)//six(1:ndigit)//'.f'
169 
170  return
171  end subroutine
172 !=======================================================================
181  subroutine io_mbyte_open(hname,ierr)
182  implicit none
183 
184  include 'SIZE' ! NID
185  include 'INPUT' ! ifmpiio
186  include 'RESTART' ! fid0, pid0, ifh_mbyte
187  include 'FRAMELP'
188  include 'IOTOOLD'
189 
190  ! argumnt list
191  integer fid, ierr
192  character*132 hname
193 
194  ! local variables
195  character*132 fname
196  integer itmp
197 !-----------------------------------------------------------------------
198  ! initialise variables
199  ierr = 0
200  ! work on local copy
201  fname = trim(adjustl(hname))
202 
203  ! test string length
204  itmp = len_trim(fname)
205  if (itmp.eq.0) then
206  call mntr_error(io_id,'io_mbyte_open; zero lenght fname.')
207  ierr = 1
208  return
209  endif
210 
211  ! add file number
212  call addfid(fname,fid0)
213 
214  call mntr_log(io_id,lp_ess,'Opening file: '//trim(fname))
215  if(ifmpiio) then
216  call byte_open_mpi(fname,ifh_mbyte,.false.,ierr)
217  else
218  ! add ending character; required by C
219  fname = trim(fname)//char(0)
220  call byte_open(fname,ierr)
221  endif
222 
223  return
224  end subroutine
225 !=======================================================================
231  subroutine io_mbyte_close(ierr)
232  implicit none
233 
234  include 'SIZE' ! NID
235  include 'INPUT' ! ifmpiio
236  include 'RESTART' ! pid0, ifh_mbyte
237 
238  ! argumnt list
239  integer ierr
240 
241  ! local variables
242  character*132 fname
243  integer itmp
244 !-----------------------------------------------------------------------
245  ! initialise variables
246  ierr = 0
247 
248  ! close the file
249  if (nid.eq.pid0) then
250  if(ifmpiio) then
251  call byte_close_mpi(ifh_mbyte,ierr)
252  else
253  call byte_close(ierr)
254  endif
255  endif
256 
257  return
258  end subroutine
259 !=======================================================================
271  subroutine io_mfov(offs,lvx,lvy,lvz,lnx,lny,lnz,
272  $ lnel,lnelg,lndim)
273  implicit none
274 
275  include 'SIZE'
276  include 'INPUT'
277  include 'RESTART'
278  include 'FRAMELP'
279  include 'IOTOOLD'
280 
281  ! argumnt list
282  integer*8 offs
283  integer lnx,lny,lnz,lnel,lnelg,lndim
284  real lvx(lnx,lny,lnz,lnel), lvy(lnx,lny,lnz,lnel),
285  $ lvz(lnx,lny,lnz,lnel)
286 
287  ! local variables
288  integer*8 loffs
289  integer il, ik, itmp
290 
291  real rvx(lxo*lxo*(1 + (ldim-2)*(lxo-1))*lelt),
292  $ rvy(lxo*lxo*(1 + (ldim-2)*(lxo-1))*lelt),
293  $ rvz(lxo*lxo*(1 + (ldim-2)*(lxo-1))*lelt)
294  common /scruz/ rvx, rvy, rvz
295 !-----------------------------------------------------------------------
296  if (ifreguo) then
297  ! check size of mapping space
298  if (nrg.gt.lxo) then
299  call mntr_warn(io_id,
300  $ 'io_mfov; nrg too large, reset to lxo!')
301  nrg = lxo
302  endif
303 
304  ! map to regular mesh
305  ! this code works with square element only
306  itmp = nrg**lndim
307  if (lndim.eq.2) then
308  ik=1
309  do il=1,lnel
310  call map2reg_2di_e(rvx(ik),nrg,lvx(1,1,1,il),lnx)
311  ik = ik + itmp
312  enddo
313  ik=1
314  do il=1,lnel
315  call map2reg_2di_e(rvy(ik),nrg,lvy(1,1,1,il),lnx)
316  ik = ik + itmp
317  enddo
318  else
319  ik = 1
320  do il=1,lnel
321  call map2reg_3di_e(rvx(ik),nrg,lvx(1,1,1,il),lnx)
322  ik = ik + itmp
323  enddo
324  ik = 1
325  do il=1,lnel
326  call map2reg_3di_e(rvy(ik),nrg,lvy(1,1,1,il),lnx)
327  ik = ik + itmp
328  enddo
329  ik = 1
330  do il=1,lnel
331  call map2reg_3di_e(rvz(ik),nrg,lvz(1,1,1,il),lnx)
332  ik = ik + itmp
333  enddo
334  endif
335 
336  ! shift offset taking onto account elements on processes with smaller id
337  itmp = 1 + (lndim-2)*(nrg-1)
338  ! to ensure proper integer prolongation
339  loffs = offs + int(nelb,8)*int(lndim*wdsizo*nrg*nrg*itmp,8)
340  call byte_set_view(loffs,ifh_mbyte)
341 
342  ! write vector
343  call mfo_outv(rvx,rvy,rvz,lnel,nrg,nrg,itmp)
344 
345  ! update offset
346  offs = offs + int(lnelg,8)*int(lndim*wdsizo*nrg*nrg*itmp,8)
347  else
348  ! shift offset taking onto account elements on processes with smaller id
349  ! to ensure proper integer prolongation
350  loffs = offs + int(nelb,8)*int(lndim*wdsizo*lnx*lny*lnz,8)
351  call byte_set_view(loffs,ifh_mbyte)
352 
353  ! write vector
354  call mfo_outv(lvx,lvy,lvz,lnel,lnx,lny,lnz)
355 
356  ! update offset
357  offs = offs + int(lnelg,8)*int(lndim*wdsizo*lnx*lny*lnz,8)
358  endif
359 
360  return
361  end subroutine
362 !=======================================================================
374  subroutine io_mfos(offs,lvs,lnx,lny,lnz,lnel,lnelg,lndim)
375  implicit none
376 
377  include 'SIZE'
378  include 'INPUT'
379  include 'RESTART'
380  include 'FRAMELP'
381  include 'IOTOOLD'
382 
383  ! argumnt list
384  integer*8 offs
385  integer lnx,lny,lnz,lnel,lnelg,lndim
386  real lvs(lnx,lny,lnz,lnel)
387 
388  ! local variables
389  integer*8 loffs
390  integer il, ik, itmp
391 
392  real rvs(lxo*lxo*(1 + (ldim-2)*(lxo-1))*lelt)
393  common /scruz/ rvs
394 !-----------------------------------------------------------------------
395  if (ifreguo) then
396  ! check size of mapping space
397  if (nrg.gt.lxo) then
398  call mntr_warn(io_id,
399  $ 'io_mfos; nrg too large, reset to lxo!')
400  nrg = lxo
401  endif
402 
403  ! map to regular mesh
404  ! this code works with square element only
405  itmp = nrg**lndim
406  if (lndim.eq.2) then
407  ik=1
408  do il=1,lnel
409  call map2reg_2di_e(rvs(ik),nrg,lvs(1,1,1,il),lnx)
410  ik = ik + itmp
411  enddo
412  else
413  ik = 1
414  do il=1,lnel
415  call map2reg_3di_e(rvs(ik),nrg,lvs(1,1,1,il),lnx)
416  ik = ik + itmp
417  enddo
418  endif
419 
420  ! shift offset taking onto account elements on processes with smaller id
421  itmp = 1 + (lndim-2)*(nrg-1)
422  ! to ensure proper integer prolongation
423  loffs = offs + int(nelb,8)*int(wdsizo*nrg*nrg*itmp,8)
424  call byte_set_view(loffs,ifh_mbyte)
425 
426  ! write vector
427  call mfo_outs(rvs,lnel,nrg,nrg,itmp)
428 
429  ! update offset
430  offs = offs + int(lnelg,8)*int(wdsizo*nrg*nrg*itmp,8)
431  else
432  ! shift offset taking onto account elements on processes with smaller id
433  ! to ensure proper integer prolongation
434  loffs = offs + int(nelb,8)*int(wdsizo*lnx*lny*lnz,8)
435  call byte_set_view(loffs,ifh_mbyte)
436 
437  ! write vector
438  call mfo_outs(lvs,lnel,lnx,lny,lnz)
439 
440  ! update offset
441  offs = offs + int(lnelg,8)*int(wdsizo*lnx*lny*lnz,8)
442  endif
443 
444  return
445  end subroutine
446 !=======================================================================
457  subroutine io_mfiv(offs,uf,vf,wf,lnx,lny,lnz,lnel,ifskip)
458  implicit none
459  include 'SIZE'
460  include 'INPUT'
461  include 'PARALLEL'
462  include 'RESTART'
463  include 'FRAMELP'
464  include 'IOTOOLD'
465 
466  ! argument list
467  integer*8 offs
468  integer lnx,lny,lnz,lnel
469  real uf(lnx*lny*lnz,lnel),vf(lnx*lny*lnz,lnel),
470  $ wf(lnx*lny*lnz,lnel)
471  logical ifskip
472 
473  ! local variables
474  integer lndim, nxyzr, nxyzw, nxyzv, mlen
475  integer num_recv, num_avail, nread, nelrr
476  integer el, il, kl, ll, ierr
477  integer ei, eg, jnid, jeln
478  integer msg_id(lelt)
479  integer*8 i8tmp
480 
481  ! read buffer
482  integer lrbs
483  parameter(lrbs=20*lx1*ly1*lz1*lelt)
484  real*4 w2(lrbs)
485  common /vrthov/ w2
486 
487  ! communication buffer
488  integer lwk
489  parameter(lwk = 14*lx1*ly1*lz1*lelt)
490  real*4 wk(lwk)
491  common /scrns/ wk
492 
493  ! functions
494  integer irecv, iglmax
495 !-----------------------------------------------------------------------
496  call nekgsync() ! clear outstanding message queues.
497 
498  ! check element size
499  if ((nxr.ne.lnx).or.(nyr.ne.lny).or.(nzr.ne.lnz)) then
500  call mntr_abort(io_id,'io_mfiv, wrong element size')
501  endif
502 
503  if (lnz.gt.1) then
504  lndim = 3
505  else
506  lndim = 2
507  endif
508 
509  nxyzr = lndim*lnx*lny*lnz ! element size
510  mlen = nxyzr*wdsizr ! message length
511  if (wdsizr.eq.8) then
512  nxyzw = 2*nxyzr
513  else
514  nxyzw = nxyzr
515  endif
516 
517  ! shift offset
518  i8tmp = offs + int(nelbr,8)*int(mlen,8)
519  call byte_set_view(i8tmp,ifh_mbyte)
520 
521  ! check message buffer
522  num_recv = mlen
523  num_avail = lwk*4
524  call lim_chk(num_recv,num_avail,' ',' ','io_mfiv a')
525 
526  ! setup read buffer
527  if (nid.eq.pid0r) then
528  i8tmp = int(nxyzw,8)*int(nelr,8)
529  nread = i8tmp/int(lrbs,8)
530  if (mod(i8tmp,int(lrbs,8)).ne.0) nread = nread + 1
531  if(ifmpiio) nread = iglmax(nread,1) ! needed because of collective read
532  nelrr = nelr/nread
533  endif
534  call bcast(nelrr,4)
535  call lim_chk(nxyzw*nelrr,lrbs,' ',' ','io_mfiv b')
536 
537  ! reset error flag
538  ierr = 0
539 
540  if (ifskip) then ! just read deata and do not use it (avoid communication)
541 
542  if (nid.eq.pid0r) then ! only i/o nodes will read
543  ! read blocks of size nelrr
544  kl = 0
545  do il = 1,nread
546  if (il.eq.nread) then ! clean-up
547  nelrr = nelr - (nread-1)*nelrr
548  if (nelrr.lt.0) nelrr = 0
549  endif
550 
551  if (ierr.eq.0) then
552  if (ifmpiio) then
553  call byte_read_mpi(w2,nxyzw*nelrr,-1,ifh_mbyte,ierr)
554  else
555  call byte_read (w2,nxyzw*nelrr,ierr)
556  endif
557  endif
558  enddo
559  endif
560 
561  call nekgsync()
562 
563  else ! read and use data
564  ! pre-post recieves
565  if (np.gt.1) then
566  ll = 1
567  do el=1,nelt
568  msg_id(el) = irecv(el,wk(ll),mlen)
569  ll = ll+nxyzw
570  enddo
571  endif
572 
573  if (nid.eq.pid0r.and.np.gt.1) then ! only i/o nodes will read
574  ! read blocks of size nelrr
575  kl = 0
576  do il = 1,nread
577  if (il.eq.nread) then ! clean-up
578  nelrr = nelr - (nread-1)*nelrr
579  if (nelrr.lt.0) nelrr = 0
580  endif
581 
582  if (ierr.eq.0) then
583  if (ifmpiio) then
584  call byte_read_mpi(w2,nxyzw*nelrr,-1,ifh_mbyte,ierr)
585  else
586  call byte_read (w2,nxyzw*nelrr,ierr)
587  endif
588  endif
589 
590  ! distribute data across target processors
591  ll = 1
592  do el = kl+1,kl+nelrr
593  jnid = gllnid(er(el)) ! where is er(e) now?
594  jeln = gllel(er(el))
595  if(ierr.ne.0) call rzero(w2(ll),mlen)
596  call csend(jeln,w2(ll),mlen,jnid,0) ! blocking send
597  ll = ll+nxyzw
598  enddo
599  kl = kl + nelrr
600  enddo
601  elseif (np.eq.1) then
602  if (ifmpiio) then
603  call byte_read_mpi(wk,nxyzw*nelr,-1,ifh_mbyte,ierr)
604  else
605  call byte_read(wk,nxyzw*nelr,ierr)
606  endif
607  endif
608 
609  ! distinguish between vector lengths
610  nxyzv = lnx*lny*lnz
611  if (wdsizr.eq.8) then
612  nxyzw = 2*nxyzv
613  else
614  nxyzw = nxyzv
615  endif
616 
617  ll = 1
618  do el=1,nelt
619  if (np.gt.1) then
620  call msgwait(msg_id(el))
621  ei = el
622  elseif(np.eq.1) then
623  ei = er(el)
624  endif
625  if (if_byte_sw) then
626  if (wdsizr.eq.8) then
627  call byte_reverse8(wk(ll),nxyzr*2,ierr)
628  else
629  call byte_reverse(wk(ll),nxyzr,ierr)
630  endif
631  endif
632  ! copy data
633  if (wdsizr.eq.4) then
634  call copy4r(uf(1,ei),wk(ll),nxyzv)
635  call copy4r(vf(1,ei),wk(ll + nxyzw),nxyzv)
636  if (lndim.eq.3)
637  $ call copy4r(wf(1,ei),wk(ll + 2*nxyzw),nxyzv)
638  else
639  call copy(uf(1,ei),wk(ll),nxyzv)
640  call copy(vf(1,ei),wk(ll + nxyzw),nxyzv)
641  if (lndim.eq.3)
642  $ call copy(wf(1,ei),wk(ll + 2*nxyzw),nxyzv)
643  endif
644  ll = ll+ldim*nxyzw
645  enddo
646 
647  endif
648 
649  ! update offset
650  offs = offs + int(nelgr,8)*int(mlen,8)
651 
652  call err_chk(ierr,'Error reading restart data,in io_mfiv.$')
653  return
654  end subroutine
655 !=======================================================================
666  subroutine io_mfis(offs,uf,lnx,lny,lnz,lnel,ifskip)
667  implicit none
668  include 'SIZE'
669  include 'INPUT'
670  include 'PARALLEL'
671  include 'RESTART'
672  include 'FRAMELP'
673  include 'IOTOOLD'
674 
675  ! argument list
676  integer*8 offs
677  integer lnx,lny,lnz,lnel
678  real uf(lnx*lny*lnz,lnel)
679  logical ifskip
680 
681  ! local variables
682  integer nxyzr, nxyzw, mlen
683  integer num_recv, num_avail, nread, nelrr
684  integer el, il, kl, ll, ierr
685  integer ei, eg, jnid, jeln
686  integer msg_id(lelt)
687  integer*8 i8tmp
688 
689  ! read buffer
690  integer lrbs
691  parameter(lrbs=20*lx1*ly1*lz1*lelt)
692  real*4 w2(lrbs)
693  common /vrthov/ w2
694 
695  ! communication buffer
696  integer lwk
697  parameter(lwk = 14*lx1*ly1*lz1*lelt)
698  real*4 wk(lwk)
699  common /scrns/ wk
700 
701  ! functions
702  integer irecv, iglmax
703 !-----------------------------------------------------------------------
704  call nekgsync() ! clear outstanding message queues.
705 
706  ! check element size
707  if ((nxr.ne.lnx).or.(nyr.ne.lny).or.(nzr.ne.lnz)) then
708  call mntr_abort(io_id,'io_mfis, wrong element size')
709  endif
710 
711  nxyzr = lnx*lny*lnz ! element size
712  mlen = nxyzr*wdsizr ! message length
713  if (wdsizr.eq.8) then
714  nxyzw = 2*nxyzr
715  else
716  nxyzw = nxyzr
717  endif
718 
719  ! shift offset
720  i8tmp = offs + int(nelbr,8)*int(mlen,8)
721  call byte_set_view(i8tmp,ifh_mbyte)
722 
723  ! check message buffer
724  num_recv = mlen
725  num_avail = lwk*4
726  call lim_chk(num_recv,num_avail,' ',' ','io_mfis a')
727 
728  ! setup read buffer
729  if (nid.eq.pid0r) then
730  i8tmp = int(nxyzw,8)*int(nelr,8)
731  nread = i8tmp/int(lrbs,8)
732  if (mod(i8tmp,int(lrbs,8)).ne.0) nread = nread + 1
733  if(ifmpiio) nread = iglmax(nread,1) ! needed because of collective read
734  nelrr = nelr/nread
735  endif
736  call bcast(nelrr,4)
737  call lim_chk(nxyzw*nelrr,lrbs,' ',' ','io_mfis b')
738 
739  ! reset error flag
740  ierr = 0
741 
742  if (ifskip) then ! just read deata and do not use it (avoid communication)
743 
744  if (nid.eq.pid0r) then ! only i/o nodes will read
745  ! read blocks of size nelrr
746  kl = 0
747  do il = 1,nread
748  if (il.eq.nread) then ! clean-up
749  nelrr = nelr - (nread-1)*nelrr
750  if (nelrr.lt.0) nelrr = 0
751  endif
752 
753  if (ierr.eq.0) then
754  if (ifmpiio) then
755  call byte_read_mpi(w2,nxyzw*nelrr,-1,ifh_mbyte,ierr)
756  else
757  call byte_read (w2,nxyzw*nelrr,ierr)
758  endif
759  endif
760  enddo
761  endif
762 
763  call nekgsync()
764 
765  else ! read and use data
766  ! pre-post recieves
767  if (np.gt.1) then
768  ll = 1
769  do el=1,nelt
770  msg_id(el) = irecv(el,wk(ll),mlen)
771  ll = ll+nxyzw
772  enddo
773  endif
774 
775  if (nid.eq.pid0r.and.np.gt.1) then ! only i/o nodes will read
776  ! read blocks of size nelrr
777  kl = 0
778  do il = 1,nread
779  if (il.eq.nread) then ! clean-up
780  nelrr = nelr - (nread-1)*nelrr
781  if (nelrr.lt.0) nelrr = 0
782  endif
783 
784  if (ierr.eq.0) then
785  if (ifmpiio) then
786  call byte_read_mpi(w2,nxyzw*nelrr,-1,ifh_mbyte,ierr)
787  else
788  call byte_read (w2,nxyzw*nelrr,ierr)
789  endif
790  endif
791 
792  ! distribute data across target processors
793  ll = 1
794  do el = kl+1,kl+nelrr
795  jnid = gllnid(er(el)) ! where is er(e) now?
796  jeln = gllel(er(el))
797  if(ierr.ne.0) call rzero(w2(ll),mlen)
798  call csend(jeln,w2(ll),mlen,jnid,0) ! blocking send
799  ll = ll+nxyzw
800  enddo
801  kl = kl + nelrr
802  enddo
803  elseif (np.eq.1) then
804  if (ifmpiio) then
805  call byte_read_mpi(wk,nxyzw*nelr,-1,ifh_mbyte,ierr)
806  else
807  call byte_read(wk,nxyzw*nelr,ierr)
808  endif
809  endif
810 
811  ll = 1
812  do el=1,nelt
813  if (np.gt.1) then
814  call msgwait(msg_id(el))
815  ei = el
816  elseif(np.eq.1) then
817  ei = er(el)
818  endif
819  if (if_byte_sw) then
820  if (wdsizr.eq.8) then
821  call byte_reverse8(wk(ll),nxyzw,ierr)
822  else
823  call byte_reverse(wk(ll),nxyzw,ierr)
824  endif
825  endif
826  ! copy data
827  if (wdsizr.eq.4) then
828  call copy4r(uf(1,ei),wk(ll),nxyzr)
829  else
830  call copy (uf(1,ei),wk(ll),nxyzr)
831  endif
832  ll = ll+nxyzw
833  enddo
834 
835  endif
836 
837  ! update offset
838  offs = offs + int(nelgr,8)*int(mlen,8)
839 
840  call err_chk(ierr,'Error reading restart data,in io_mfis.$')
841  return
842  end subroutine
843 !=======================================================================
#define byte_reverse8
Definition: byte.c:34
#define byte_open
Definition: byte.c:35
#define byte_reverse
Definition: byte.c:33
#define byte_close
Definition: byte.c:36
#define byte_read
Definition: byte.c:38
subroutine byte_open_mpi(fnamei, mpi_fh, ifro, ierr)
Definition: byte_mpi.f:13
subroutine byte_read_mpi(buf, icount, iorank, mpi_fh, ierr)
Definition: byte_mpi.f:48
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 msgwait(imsg)
Definition: comm_mpi.f:489
subroutine nekgsync()
Definition: comm_mpi.f:502
subroutine csend(mtype, buf, len, jnid, jpid)
Definition: comm_mpi.f:303
subroutine bcast(buf, len)
Definition: comm_mpi.f:431
subroutine err_chk(ierr, string)
Definition: comm_mpi.f:558
subroutine lim_chk(n, m, avar5, lvar5, sub_name10)
Definition: convect.f:454
integer function gllel(ieg)
Definition: dprocmap.f:183
integer function gllnid(ieg)
Definition: dprocmap.f:161
subroutine io_mfos(offs, lvs, lnx, lny, lnz, lnel, lnelg, lndim)
Write single scalar to the file.
Definition: io_tools.f:375
subroutine io_mfov(offs, lvx, lvy, lvz, lnx, lny, lnz, lnel, lnelg, lndim)
Write single vector to the file.
Definition: io_tools.f:273
subroutine io_file_close()
Close all opened files up to sotred max unit numer.
Definition: io_tools.f:81
subroutine io_mbyte_open(hname, ierr)
Open field file.
Definition: io_tools.f:182
subroutine io_file_freeid(iunit, ierr)
Get free file unit number and store max unit value.
Definition: io_tools.f:47
subroutine io_register()
Register io tool module.
Definition: io_tools.f:11
subroutine io_mfo_fname(fname, bname, prefix, ierr)
Generate file name according to nek rulles without opening the file.
Definition: io_tools.f:109
subroutine io_mbyte_close(ierr)
Close field file.
Definition: io_tools.f:232
subroutine io_mfis(offs, uf, lnx, lny, lnz, lnel, ifskip)
Read scalar filed from the file.
Definition: io_tools.f:667
subroutine io_mfiv(offs, uf, vf, wf, lnx, lny, lnz, lnel, ifskip)
Read vector filed from the file.
Definition: io_tools.f:458
subroutine mntr_warn(mid, logs)
Write warning message.
Definition: mntrlog.f:803
subroutine mntr_mod_is_name_reg(mid, mname)
Check if module name is registered and return its id.
Definition: mntrlog.f:459
subroutine mntr_abort(mid, logs)
Abort simulation.
Definition: mntrlog.f:837
subroutine mntr_log(mid, priority, logs)
Write log message.
Definition: mntrlog.f:600
subroutine mntr_error(mid, logs)
Write error message.
Definition: mntrlog.f:820
subroutine mntr_mod_reg(mid, pmid, mname, mdscr)
Register new module.
Definition: mntrlog.f:346
subroutine addfid(fname, fid)
Definition: ic.f:2512
subroutine copy(a, b, n)
Definition: math.f:260
subroutine rzero(a, n)
Definition: math.f:208
subroutine map2reg_3di_e(uf, n, uc, m)
Definition: postpro.f:594
subroutine map2reg_2di_e(uf, n, uc, m)
Definition: postpro.f:565
subroutine mfo_outv(u, v, w, nel, mx, my, mz)
Definition: prepost.f:1726
subroutine copy4r(a, b, n)
Definition: prepost.f:568
subroutine mfo_outs(u, nel, mx, my, mz)
Definition: prepost.f:1633