xref: /phasta/phSolver/common/readnblk.f (revision e81a6dc1f924882a3384f8e849eaba3a3575cfa5)
159599516SKenneth E. Jansenc  readnblk.f (pronounce "Reed and Block Dot Eff") contains:
259599516SKenneth E. Jansenc
359599516SKenneth E. Jansenc    module readarrays ("Red Arrays") -- contains the arrays that
459599516SKenneth E. Jansenc     are read in from binary files but not immediately blocked
559599516SKenneth E. Jansenc     through pointers.
659599516SKenneth E. Jansenc
759599516SKenneth E. Jansenc    subroutine readnblk ("Reed and Block") -- allocates space for
859599516SKenneth E. Jansenc     and reads data to be contained in module readarrays.  Reads
959599516SKenneth E. Jansenc     all remaining data and blocks them with pointers.
1059599516SKenneth E. Jansenc
1159599516SKenneth E. Jansen      module readarrays
1259599516SKenneth E. Jansen
1359599516SKenneth E. Jansen      real*8, allocatable :: point2x(:,:)
1459599516SKenneth E. Jansen      real*8, allocatable :: qold(:,:)
1559599516SKenneth E. Jansen      real*8, allocatable :: uold(:,:)
1659599516SKenneth E. Jansen      real*8, allocatable :: acold(:,:)
1759599516SKenneth E. Jansen      integer, allocatable :: iBCtmp(:)
1859599516SKenneth E. Jansen      real*8, allocatable :: BCinp(:,:)
1959599516SKenneth E. Jansen
2059599516SKenneth E. Jansen      integer, allocatable :: point2ilwork(:)
2159599516SKenneth E. Jansen      integer, allocatable :: nBC(:)
2259599516SKenneth E. Jansen      integer, allocatable :: point2iper(:)
239a6935e5SKenneth E. Jansen      integer, target, allocatable :: point2ifath(:)
249a6935e5SKenneth E. Jansen      integer, target, allocatable :: point2nsons(:)
2559599516SKenneth E. Jansen
2659599516SKenneth E. Jansen      end module
2759599516SKenneth E. Jansen
2859599516SKenneth E. Jansen      subroutine readnblk
29e5afe575SCameron Smith      use iso_c_binding
3059599516SKenneth E. Jansen      use readarrays
31e5afe575SCameron Smith      use phio
32ab645d52SCameron Smith      use syncio
33ab645d52SCameron Smith      use posixio
34ab645d52SCameron Smith      use streamio
3559599516SKenneth E. Jansen      include "common.h"
3602ce253eSCameron Smith
379a6935e5SKenneth E. Jansen      real*8, target, allocatable :: xread(:,:), qread(:,:), acread(:,:)
389a6935e5SKenneth E. Jansen      real*8, target, allocatable :: uread(:,:)
399a6935e5SKenneth E. Jansen      real*8, target, allocatable :: BCinpread(:,:)
409a6935e5SKenneth E. Jansen      integer, target, allocatable :: iperread(:), iBCtmpread(:)
419a6935e5SKenneth E. Jansen      integer, target, allocatable :: ilworkread(:), nBCread(:)
4259599516SKenneth E. Jansen      character*10 cname2
4359599516SKenneth E. Jansen      character*8 mach2
4459599516SKenneth E. Jansen      character*30 fmt1
4559599516SKenneth E. Jansen      character*255 fname1,fnamer,fnamelr
4659599516SKenneth E. Jansen      character*255 warning
4759599516SKenneth E. Jansen      integer igeomBAK, ibndc, irstin, ierr
489a6935e5SKenneth E. Jansen      integer, target :: intfromfile(50) ! integers read from headers
4959599516SKenneth E. Jansen      integer :: descriptor, descriptorG, GPID, color, nfiles, nfields
5059599516SKenneth E. Jansen      integer ::  numparts, nppf
5159599516SKenneth E. Jansen      integer :: ierr_io, numprocs, itmp, itmp2
52ade0e30fSCameron Smith      integer :: ignored
53d7abaf6cSCameron Smith      integer :: fileFmt
5459599516SKenneth E. Jansen      character*255 fname2, temp2
5559599516SKenneth E. Jansen      character*64 temp1
56e5afe575SCameron Smith      type(c_ptr) :: handle
57e5afe575SCameron Smith      character(len=1024) :: dataInt, dataDbl
58e5afe575SCameron Smith      dataInt = c_char_'integer'//c_null_char
59e5afe575SCameron Smith      dataDbl = c_char_'double'//c_null_char
6059599516SKenneth E. Jansenc
6159599516SKenneth E. Jansenc.... determine the step number to start with
6259599516SKenneth E. Jansenc
6359599516SKenneth E. Jansen      open(unit=72,file='numstart.dat',status='old')
6459599516SKenneth E. Jansen      read(72,*) irstart
6559599516SKenneth E. Jansen      close(72)
6659599516SKenneth E. Jansen      lstep=irstart ! in case restart files have no fields
6759599516SKenneth E. Jansen
6859599516SKenneth E. Jansen      fname1='geombc.dat'
6959599516SKenneth E. Jansen      fname1= trim(fname1)  // cname2(myrank+1)
7059599516SKenneth E. Jansen
7159599516SKenneth E. Jansen      itmp=1
7259599516SKenneth E. Jansen      if (irstart .gt. 0) itmp = int(log10(float(irstart)))+1
7359599516SKenneth E. Jansen      write (fmt1,"('(''restart.'',i',i1,',1x)')") itmp
7459599516SKenneth E. Jansen      write (fnamer,fmt1) irstart
7559599516SKenneth E. Jansen      fnamer = trim(fnamer) // cname2(myrank+1)
7659599516SKenneth E. Jansenc
7759599516SKenneth E. Jansenc.... open input files
7859599516SKenneth E. Jansenc.... input the geometry parameters
7959599516SKenneth E. Jansenc
8059599516SKenneth E. Jansen      nfiles = nsynciofiles
8159599516SKenneth E. Jansen      numparts = numpe !This is the common settings. Beware if you try to compute several parts per process
8259599516SKenneth E. Jansen
8359599516SKenneth E. Jansen      itwo=2
8459599516SKenneth E. Jansen      ione=1
8559599516SKenneth E. Jansen      ieleven=11
8659599516SKenneth E. Jansen      itmp = int(log10(float(myrank+1)))+1
8759599516SKenneth E. Jansen
88ab645d52SCameron Smith      if( nsynciofiles .eq. -1 ) then
89d7abaf6cSCameron Smith        fileFmt = PHIO_STREAM
90ab645d52SCameron Smith        call streamio_setup(grstream, fhandle)
91ab645d52SCameron Smith      else if( nsynciofiles .eq. 0 ) then
92ab645d52SCameron Smith        call posixio_setup(fhandle, c_char_'r')
93ab645d52SCameron Smith      else if( nsynciofiles .gt. 1 ) then
94ab645d52SCameron Smith        call syncio_setup_read(nsynciofiles, fhandle)
952efdc748SKenneth E. Jansen      end if
96d7abaf6cSCameron Smith      call phio_constructName(fileFmt,
97d7abaf6cSCameron Smith     &        c_char_'geombc' // char(0), fname1)
98d7abaf6cSCameron Smith      call phio_openfile(fname1, fhandle);
9959599516SKenneth E. Jansen
100d5d2f64dSCameron Smith      call phio_readheader(fhandle,c_char_'number of nodes' // char(0),
101e5afe575SCameron Smith     & c_loc(numnp),ione, dataInt, iotype)
10259599516SKenneth E. Jansen
103d5d2f64dSCameron Smith      call phio_readheader(fhandle,c_char_'number of modes' // char(0),
104e5afe575SCameron Smith     & c_loc(nshg),ione, dataInt, iotype)
10559599516SKenneth E. Jansen
106d5d2f64dSCameron Smith      call phio_readheader(fhandle,
107e5afe575SCameron Smith     &  c_char_'number of interior elements' // char(0),
108e5afe575SCameron Smith     &  c_loc(numel),ione, dataInt, iotype)
10959599516SKenneth E. Jansen
110d5d2f64dSCameron Smith      call phio_readheader(fhandle,
111e5afe575SCameron Smith     &  c_char_'number of boundary elements' // char(0),
112e5afe575SCameron Smith     &  c_loc(numelb),ione, dataInt, iotype)
11359599516SKenneth E. Jansen
114d5d2f64dSCameron Smith      call phio_readheader(fhandle,
115e5afe575SCameron Smith     &  c_char_'maximum number of element nodes' // char(0),
116e5afe575SCameron Smith     &  c_loc(nen),ione, dataInt, iotype)
11759599516SKenneth E. Jansen
118d5d2f64dSCameron Smith      call phio_readheader(fhandle,
119e5afe575SCameron Smith     &  c_char_'number of interior tpblocks' // char(0),
120e5afe575SCameron Smith     &  c_loc(nelblk),ione, dataInt, iotype)
12159599516SKenneth E. Jansen
122d5d2f64dSCameron Smith      call phio_readheader(fhandle,
123e5afe575SCameron Smith     & c_char_'number of boundary tpblocks' // char(0),
124e5afe575SCameron Smith     & c_loc(nelblb),ione, dataInt, iotype)
12559599516SKenneth E. Jansen
126d5d2f64dSCameron Smith      call phio_readheader(fhandle,
127e5afe575SCameron Smith     & c_char_'number of nodes with Dirichlet BCs' // char(0),
128e5afe575SCameron Smith     & c_loc(numpbc),ione, dataInt, iotype)
12959599516SKenneth E. Jansen
130d5d2f64dSCameron Smith      call phio_readheader(fhandle,
131e5afe575SCameron Smith     & c_char_'number of shape functions' // char(0),
132e5afe575SCameron Smith     & c_loc(ntopsh),ione, dataInt, iotype)
13359599516SKenneth E. Jansenc
13459599516SKenneth E. Jansenc.... calculate the maximum number of boundary element nodes
13559599516SKenneth E. Jansenc
13659599516SKenneth E. Jansen      nenb = 0
13759599516SKenneth E. Jansen      do i = 1, melCat
13859599516SKenneth E. Jansen         if (nen .eq. nenCat(i,nsd)) nenb = max(nenCat(i,nsd-1), nenb)
13959599516SKenneth E. Jansen      enddo
14059599516SKenneth E. Jansen      if (myrank == master) then
14159599516SKenneth E. Jansen         if (nenb .eq. 0) call error ('input   ','nen     ',nen)
14259599516SKenneth E. Jansen      endif
14359599516SKenneth E. Jansenc
14459599516SKenneth E. Jansenc.... setup some useful constants
14559599516SKenneth E. Jansenc
14659599516SKenneth E. Jansen      I3nsd  = nsd / 3          ! nsd=3 integer flag
14759599516SKenneth E. Jansen      E3nsd  = float(I3nsd)     ! nsd=3 real    flag
14859599516SKenneth E. Jansen      if(matflg(1,1).lt.0) then
14959599516SKenneth E. Jansen         nflow = nsd + 1
15059599516SKenneth E. Jansen      else
15159599516SKenneth E. Jansen         nflow = nsd + 2
15259599516SKenneth E. Jansen      endif
15359599516SKenneth E. Jansen      ndof   = nsd + 2
15459599516SKenneth E. Jansen      nsclr=impl(1)/100
15559599516SKenneth E. Jansen      ndof=ndof+nsclr           ! number of sclr transport equations to solve
15659599516SKenneth E. Jansen
15759599516SKenneth E. Jansen      ndofBC = ndof + I3nsd     ! dimension of BC array
15859599516SKenneth E. Jansen      ndiBCB = 2                ! dimension of iBCB array
15959599516SKenneth E. Jansen      ndBCB  = ndof + 1         ! dimension of BCB array
16059599516SKenneth E. Jansen      nsymdf = (ndof*(ndof + 1)) / 2 ! symm. d.o.f.'s
16159599516SKenneth E. Jansenc
16259599516SKenneth E. Jansenc.... ----------------------> Communication tasks <--------------------
16359599516SKenneth E. Jansenc
16459599516SKenneth E. Jansen      if(numpe > 1) then
165d5d2f64dSCameron Smith         call phio_readheader(fhandle,
166e5afe575SCameron Smith     &    c_char_'size of ilwork array' // char(0),
167e5afe575SCameron Smith     &    c_loc(nlwork),ione, dataInt, iotype)
16859599516SKenneth E. Jansen
169d5d2f64dSCameron Smith         call phio_readheader(fhandle,
170e5afe575SCameron Smith     &    c_char_'ilwork' //char(0),
171e5afe575SCameron Smith     &    c_loc(nlwork),ione, dataInt, iotype)
17259599516SKenneth E. Jansen
17359599516SKenneth E. Jansen         allocate( point2ilwork(nlwork) )
17459599516SKenneth E. Jansen         allocate( ilworkread(nlwork) )
175d5d2f64dSCameron Smith         call phio_readdatablock(fhandle, c_char_'ilwork' // char(0),
176bc62cfd4SCameron Smith     &      c_loc(ilworkread), nlwork, dataInt , iotype)
17759599516SKenneth E. Jansen
17859599516SKenneth E. Jansen         point2ilwork = ilworkread
17959599516SKenneth E. Jansen         call ctypes (point2ilwork)
18059599516SKenneth E. Jansen      else
18159599516SKenneth E. Jansen           nlwork=1
18259599516SKenneth E. Jansen           allocate( point2ilwork(1))
18359599516SKenneth E. Jansen           nshg0 = nshg
18459599516SKenneth E. Jansen      endif
18559599516SKenneth E. Jansen
18659599516SKenneth E. Jansen      itwo=2
18759599516SKenneth E. Jansen
188d5d2f64dSCameron Smith      call phio_readheader(fhandle,
189e5afe575SCameron Smith     & c_char_'co-ordinates' // char(0),
190e5afe575SCameron Smith     & c_loc(intfromfile),itwo, dataDbl, iotype)
19159599516SKenneth E. Jansen      numnp=intfromfile(1)
19259599516SKenneth E. Jansen      allocate( point2x(numnp,nsd) )
19359599516SKenneth E. Jansen      allocate( xread(numnp,nsd) )
19459599516SKenneth E. Jansen      ixsiz=numnp*nsd
195d5d2f64dSCameron Smith      call phio_readdatablock(fhandle,
196bc62cfd4SCameron Smith     & c_char_'co-ordinates' // char(0),
197bc62cfd4SCameron Smith     & c_loc(xread),ixsiz, dataDbl, iotype)
19859599516SKenneth E. Jansen      point2x = xread
19959599516SKenneth E. Jansenc
20059599516SKenneth E. Jansenc.... read in and block out the connectivity
20159599516SKenneth E. Jansenc
20259599516SKenneth E. Jansen      call genblk (IBKSIZ)
20359599516SKenneth E. Jansenc
20459599516SKenneth E. Jansenc.... read the boundary condition mapping array
20559599516SKenneth E. Jansenc
20659599516SKenneth E. Jansen      ione=1
207d5d2f64dSCameron Smith      call phio_readheader(fhandle,
208e5afe575SCameron Smith     & c_char_'bc mapping array' // char(0),
209e5afe575SCameron Smith     & c_loc(nshg),ione, dataInt, iotype)
21059599516SKenneth E. Jansen
21159599516SKenneth E. Jansen      allocate( nBC(nshg) )
21259599516SKenneth E. Jansen
21359599516SKenneth E. Jansen      allocate( nBCread(nshg) )
21459599516SKenneth E. Jansen
215d5d2f64dSCameron Smith      call phio_readdatablock(fhandle,
216bc62cfd4SCameron Smith     & c_char_'bc mapping array' // char(0),
217bc62cfd4SCameron Smith     & c_loc(nBCread), nshg, dataInt, iotype)
21859599516SKenneth E. Jansen
21959599516SKenneth E. Jansen      nBC=nBCread
22059599516SKenneth E. Jansenc
22159599516SKenneth E. Jansenc.... read the temporary iBC array
22259599516SKenneth E. Jansenc
22359599516SKenneth E. Jansen      ione=1
224d5d2f64dSCameron Smith      call phio_readheader(fhandle,
225e5afe575SCameron Smith     & c_char_'bc codes array' // char(0),
226e5afe575SCameron Smith     & c_loc(numpbc),ione, dataInt, iotype)
22759599516SKenneth E. Jansen
22859599516SKenneth E. Jansen      if ( numpbc > 0 ) then
22959599516SKenneth E. Jansen        allocate( iBCtmp(numpbc) )
23059599516SKenneth E. Jansen        allocate( iBCtmpread(numpbc) )
23159599516SKenneth E. Jansen      else
23259599516SKenneth E. Jansen        allocate( iBCtmp(1) )
23359599516SKenneth E. Jansen        allocate( iBCtmpread(1) )
23459599516SKenneth E. Jansen      endif
235d5d2f64dSCameron Smith      call phio_readdatablock(fhandle,
236bc62cfd4SCameron Smith     & c_char_'bc codes array' // char(0),
237bc62cfd4SCameron Smith     & c_loc(iBCtmpread), numpbc, dataInt, iotype)
23859599516SKenneth E. Jansen
23959599516SKenneth E. Jansen      if ( numpbc > 0 ) then
24059599516SKenneth E. Jansen         iBCtmp=iBCtmpread
24159599516SKenneth E. Jansen      else  ! sometimes a partition has no BC's
24259599516SKenneth E. Jansen         deallocate( iBCtmpread)
24359599516SKenneth E. Jansen         iBCtmp=0
24459599516SKenneth E. Jansen      endif
24559599516SKenneth E. Jansenc
24659599516SKenneth E. Jansenc.... read boundary condition data
24759599516SKenneth E. Jansenc
24859599516SKenneth E. Jansen      ione=1
249d5d2f64dSCameron Smith      call phio_readheader(fhandle,
250e5afe575SCameron Smith     & c_char_'boundary condition array' // char(0),
251e5afe575SCameron Smith     & c_loc(intfromfile),ione, dataDbl, iotype)
25259599516SKenneth E. Jansen
25359599516SKenneth E. Jansen      if ( numpbc > 0 ) then
25459599516SKenneth E. Jansen         allocate( BCinp(numpbc,ndof+7) )
25559599516SKenneth E. Jansen         nsecondrank=intfromfile(1)/numpbc
25659599516SKenneth E. Jansen         allocate( BCinpread(numpbc,nsecondrank) )
25759599516SKenneth E. Jansen         iBCinpsiz=intfromfile(1)
25859599516SKenneth E. Jansen      else
25959599516SKenneth E. Jansen         allocate( BCinp(1,ndof+7) )
26059599516SKenneth E. Jansen         allocate( BCinpread(0,0) ) !dummy
26159599516SKenneth E. Jansen         iBCinpsiz=intfromfile(1)
26259599516SKenneth E. Jansen      endif
26359599516SKenneth E. Jansen
264d5d2f64dSCameron Smith      call phio_readdatablock(fhandle,
265bc62cfd4SCameron Smith     & c_char_'boundary condition array' // char(0),
266bc62cfd4SCameron Smith     & c_loc(BCinpread), iBCinpsiz, dataDbl, iotype)
26759599516SKenneth E. Jansen
26859599516SKenneth E. Jansen      if ( numpbc > 0 ) then
26959599516SKenneth E. Jansen         BCinp(:,1:(ndof+7))=BCinpread(:,1:(ndof+7))
27059599516SKenneth E. Jansen      else  ! sometimes a partition has no BC's
27159599516SKenneth E. Jansen         deallocate(BCinpread)
27259599516SKenneth E. Jansen         BCinp=0
27359599516SKenneth E. Jansen      endif
27459599516SKenneth E. Jansenc
27559599516SKenneth E. Jansenc.... read periodic boundary conditions
27659599516SKenneth E. Jansenc
27759599516SKenneth E. Jansen      ione=1
278d5d2f64dSCameron Smith      call phio_readheader(fhandle,
279e5afe575SCameron Smith     & c_char_'periodic masters array' // char(0),
280e5afe575SCameron Smith     & c_loc(nshg), ione, dataInt, iotype)
28159599516SKenneth E. Jansen      allocate( point2iper(nshg) )
28259599516SKenneth E. Jansen      allocate( iperread(nshg) )
283d5d2f64dSCameron Smith      call phio_readdatablock(fhandle,
284bc62cfd4SCameron Smith     & c_char_'periodic masters array' // char(0),
285bc62cfd4SCameron Smith     & c_loc(iperread), nshg, dataInt, iotype)
28659599516SKenneth E. Jansen      point2iper=iperread
28759599516SKenneth E. Jansenc
28859599516SKenneth E. Jansenc.... generate the boundary element blocks
28959599516SKenneth E. Jansenc
29059599516SKenneth E. Jansen      call genbkb (ibksiz)
29159599516SKenneth E. Jansenc
29259599516SKenneth E. Jansenc  Read in the nsons and ifath arrays if needed
29359599516SKenneth E. Jansenc
29459599516SKenneth E. Jansenc  There is a fundamental shift in the meaning of ifath based on whether
29559599516SKenneth E. Jansenc  there exist homogenous directions in the flow.
29659599516SKenneth E. Jansenc
29759599516SKenneth E. Jansenc  HOMOGENOUS DIRECTIONS EXIST:  Here nfath is the number of inhomogenous
29859599516SKenneth E. Jansenc  points in the TOTAL mesh.  That is to say that each partition keeps a
29959599516SKenneth E. Jansenc  link to  ALL inhomogenous points.  This link is furthermore not to the
30059599516SKenneth E. Jansenc  sms numbering but to the original structured grid numbering.  These
30159599516SKenneth E. Jansenc  inhomogenous points are thought of as fathers, with their sons being all
30259599516SKenneth E. Jansenc  the points in the homogenous directions that have this father's
30359599516SKenneth E. Jansenc  inhomogeneity.  The array ifath takes as an arguement the sms numbering
30459599516SKenneth E. Jansenc  and returns as a result the father.
30559599516SKenneth E. Jansenc
30659599516SKenneth E. Jansenc  In this case nsons is the number of sons that each father has and ifath
30759599516SKenneth E. Jansenc  is an array which tells the
30859599516SKenneth E. Jansenc
30959599516SKenneth E. Jansenc  NO HOMOGENOUS DIRECTIONS.  In this case the mesh would grow to rapidly
31059599516SKenneth E. Jansenc  if we followed the above strategy since every partition would index its
31159599516SKenneth E. Jansenc  points to the ENTIRE mesh.  Furthermore, there would never be a need
31259599516SKenneth E. Jansenc  to average to a node off processor since there is no spatial averaging.
31359599516SKenneth E. Jansenc  Therefore, to properly account for this case we must recognize it and
31459599516SKenneth E. Jansenc  inerrupt certain actions (i.e. assembly of the average across partitions).
31559599516SKenneth E. Jansenc  This case is easily identified by noting that maxval(nsons) =1 (i.e. no
31659599516SKenneth E. Jansenc  father has any sons).  Reiterating to be clear, in this case ifath does
31759599516SKenneth E. Jansenc  not point to a global numbering but instead just points to itself.
31859599516SKenneth E. Jansenc
31959599516SKenneth E. Jansen      nfath=1  ! some architectures choke on a zero or undeclared
32059599516SKenneth E. Jansen                 ! dimension variable.  This sets it to a safe, small value.
32159599516SKenneth E. Jansen      if(((iLES .lt. 20) .and. (iLES.gt.0))
32259599516SKenneth E. Jansen     &                   .or. (itwmod.gt.0)  ) then ! don't forget same
32359599516SKenneth E. Jansen                                                    ! conditional in proces.f
32459599516SKenneth E. Jansen                                                    ! needed for alloc
32559599516SKenneth E. Jansen         ione=1
32659599516SKenneth E. Jansen         if(nohomog.gt.0) then
327d5d2f64dSCameron Smith            call phio_readheader(fhandle,
328e5afe575SCameron Smith     &       c_char_'number of father-nodes' // char(0),
329e5afe575SCameron Smith     &       c_loc(nfath), ione, dataInt, iotype)
33059599516SKenneth E. Jansen
331d5d2f64dSCameron Smith            call phio_readheader(fhandle,
332e5afe575SCameron Smith     &       c_char_'number of son-nodes for each father' // char(0),
333e5afe575SCameron Smith     &       c_loc(nfath), ione, dataInt, iotype)
33459599516SKenneth E. Jansen
33559599516SKenneth E. Jansen            allocate (point2nsons(nfath))
33659599516SKenneth E. Jansen
337d5d2f64dSCameron Smith            call phio_readdatablock(fhandle,
338bc62cfd4SCameron Smith     &       c_char_'number of son-nodes for each father' // char(0),
339bc62cfd4SCameron Smith     &       c_loc(point2nsons),nfath, dataInt, iotype)
34059599516SKenneth E. Jansen
341d5d2f64dSCameron Smith            call phio_readheader(fhandle,
342e5afe575SCameron Smith     &       c_char_'keyword ifath' // char(0),
343e5afe575SCameron Smith     &       c_loc(nshg), ione, dataInt, iotype);
34459599516SKenneth E. Jansen
34559599516SKenneth E. Jansen            allocate (point2ifath(nshg))
34659599516SKenneth E. Jansen
347d5d2f64dSCameron Smith            call phio_readdatablock(fhandle,
348bc62cfd4SCameron Smith     &       c_char_'keyword ifath' // char(0),
349bc62cfd4SCameron Smith     &       c_loc(point2ifath), nshg, dataInt, iotype)
35059599516SKenneth E. Jansen
35159599516SKenneth E. Jansen            nsonmax=maxval(point2nsons)
35259599516SKenneth E. Jansen         else  ! this is the case where there is no homogeneity
35359599516SKenneth E. Jansen               ! therefore ever node is a father (too itself).  sonfath
35459599516SKenneth E. Jansen               ! (a routine in NSpre) will set this up but this gives
35559599516SKenneth E. Jansen               ! you an option to avoid that.
35659599516SKenneth E. Jansen            nfath=nshg
35759599516SKenneth E. Jansen            allocate (point2nsons(nfath))
35859599516SKenneth E. Jansen            point2nsons=1
35959599516SKenneth E. Jansen            allocate (point2ifath(nshg))
36059599516SKenneth E. Jansen            do i=1,nshg
36159599516SKenneth E. Jansen               point2ifath(i)=i
36259599516SKenneth E. Jansen            enddo
36359599516SKenneth E. Jansen            nsonmax=1
36459599516SKenneth E. Jansen         endif
36559599516SKenneth E. Jansen      else
36659599516SKenneth E. Jansen         allocate (point2nsons(1))
36759599516SKenneth E. Jansen         allocate (point2ifath(1))
36859599516SKenneth E. Jansen      endif
36959599516SKenneth E. Jansen
370d7abaf6cSCameron Smith      call phio_closefile(fhandle);
37159599516SKenneth E. Jansenc.... Read restart files
3722efdc748SKenneth E. Jansen      if(nsynciofiles.gt.0) then
37336adee64SCameron Smith        fnamer = c_char_"restart-dat."//c_null_char
3742efdc748SKenneth E. Jansen      else
3752efdc748SKenneth E. Jansen        fnamer = c_char_"restart."//c_null_char
3762efdc748SKenneth E. Jansen      endif
3772efdc748SKenneth E. Jansen
378d7abaf6cSCameron Smith      if( nsynciofiles .eq. -1 ) then
379d7abaf6cSCameron Smith        fileFmt = PHIO_STREAM
380d7abaf6cSCameron Smith        call streamio_setup(grstream, fhandle)
381d7abaf6cSCameron Smith      else if( nsynciofiles .eq. 0 ) then
382d7abaf6cSCameron Smith        call posixio_setup(fhandle, c_char_'r')
383d7abaf6cSCameron Smith      else if( nsynciofiles .gt. 1 ) then
384d7abaf6cSCameron Smith        call syncio_setup_read(nsynciofiles, fhandle)
385d7abaf6cSCameron Smith      end if
386d7abaf6cSCameron Smith      call phio_constructName(fileFmt,
387d7abaf6cSCameron Smith     &        c_char_'restart' // char(0), fnamer)
388*e81a6dc1SCameron Smith      call phio_appendInt(fnamer, irstart)
389d7abaf6cSCameron Smith      call phio_openfile(fnamer, fhandle);
39059599516SKenneth E. Jansen
39159599516SKenneth E. Jansen      ithree=3
39259599516SKenneth E. Jansen
39359599516SKenneth E. Jansen      itmp = int(log10(float(myrank+1)))+1
39459599516SKenneth E. Jansen
39559599516SKenneth E. Jansen      intfromfile=0
396d5d2f64dSCameron Smith      call phio_readheader(fhandle,
397e5afe575SCameron Smith     & c_char_'solution' // char(0),
398e5afe575SCameron Smith     & c_loc(intfromfile), ithree, dataInt, iotype)
39959599516SKenneth E. Jansenc
40059599516SKenneth E. Jansenc.... read the values of primitive variables into q
40159599516SKenneth E. Jansenc
40259599516SKenneth E. Jansen      allocate( qold(nshg,ndof) )
40359599516SKenneth E. Jansen      if(intfromfile(1).ne.0) then
40459599516SKenneth E. Jansen         nshg2=intfromfile(1)
40559599516SKenneth E. Jansen         ndof2=intfromfile(2)
40659599516SKenneth E. Jansen         lstep=intfromfile(3)
40759599516SKenneth E. Jansen         if(ndof2.ne.ndof) then
40859599516SKenneth E. Jansen
40959599516SKenneth E. Jansen         endif
41059599516SKenneth E. Jansen        if (nshg2 .ne. nshg)
41159599516SKenneth E. Jansen     &        call error ('restar  ', 'nshg   ', nshg)
41259599516SKenneth E. Jansen         allocate( qread(nshg,ndof2) )
41359599516SKenneth E. Jansen         iqsiz=nshg*ndof2
414d5d2f64dSCameron Smith         call phio_readdatablock(fhandle,
415bc62cfd4SCameron Smith     &    c_char_'solution' // char(0),
416bc62cfd4SCameron Smith     &    c_loc(qread),iqsiz, dataDbl,iotype)
41759599516SKenneth E. Jansen         qold(:,1:ndof)=qread(:,1:ndof)
41859599516SKenneth E. Jansen         deallocate(qread)
41959599516SKenneth E. Jansen      else
42059599516SKenneth E. Jansen         if (myrank.eq.master) then
42159599516SKenneth E. Jansen            if (matflg(1,1).eq.0) then ! compressible
42259599516SKenneth E. Jansen               warning='Solution is set to zero (with p and T to one)'
42359599516SKenneth E. Jansen            else
42459599516SKenneth E. Jansen               warning='Solution is set to zero'
42559599516SKenneth E. Jansen            endif
42659599516SKenneth E. Jansen            write(*,*) warning
42759599516SKenneth E. Jansen         endif
42859599516SKenneth E. Jansen         qold=zero
42959599516SKenneth E. Jansen         if (matflg(1,1).eq.0) then ! compressible
43059599516SKenneth E. Jansen            qold(:,1)=one ! avoid zero pressure
43159599516SKenneth E. Jansen            qold(:,nflow)=one ! avoid zero temperature
43259599516SKenneth E. Jansen         endif
43359599516SKenneth E. Jansen      endif
43459599516SKenneth E. Jansen
43559599516SKenneth E. Jansen      intfromfile=0
436d5d2f64dSCameron Smith      call phio_readheader(fhandle,
437e5afe575SCameron Smith     & c_char_'time derivative of solution' // char(0),
438e5afe575SCameron Smith     & c_loc(intfromfile), ithree, dataInt, iotype)
43959599516SKenneth E. Jansen      allocate( acold(nshg,ndof) )
44059599516SKenneth E. Jansen      if(intfromfile(1).ne.0) then
44159599516SKenneth E. Jansen         nshg2=intfromfile(1)
44259599516SKenneth E. Jansen         ndof2=intfromfile(2)
44359599516SKenneth E. Jansen         lstep=intfromfile(3)
44459599516SKenneth E. Jansen
44559599516SKenneth E. Jansen         if (nshg2 .ne. nshg)
44659599516SKenneth E. Jansen     &        call error ('restar  ', 'nshg   ', nshg)
44759599516SKenneth E. Jansen         allocate( acread(nshg,ndof2) )
44859599516SKenneth E. Jansen         acread=zero
44959599516SKenneth E. Jansen         iacsiz=nshg*ndof2
450d5d2f64dSCameron Smith         call phio_readdatablock(fhandle,
451bc62cfd4SCameron Smith     &    c_char_'time derivative of solution' // char(0),
452bc62cfd4SCameron Smith     &    c_loc(acread), iacsiz, dataDbl,iotype)
45359599516SKenneth E. Jansen         acold(:,1:ndof)=acread(:,1:ndof)
45459599516SKenneth E. Jansen         deallocate(acread)
45559599516SKenneth E. Jansen      else
45659599516SKenneth E. Jansen         if (myrank.eq.master) then
45759599516SKenneth E. Jansen            warning='Time derivative of solution is set to zero (SAFE)'
45859599516SKenneth E. Jansen            write(*,*) warning
45959599516SKenneth E. Jansen         endif
46059599516SKenneth E. Jansen         acold=zero
46159599516SKenneth E. Jansen      endif
46259599516SKenneth E. Jansencc
46359599516SKenneth E. Jansencc.... read the header and check it against the run data
46459599516SKenneth E. Jansencc
46559599516SKenneth E. Jansen      if (ideformwall.eq.1) then
46659599516SKenneth E. Jansen
46759599516SKenneth E. Jansen          intfromfile=0
468d5d2f64dSCameron Smith          call phio_readheader(fhandle,
469e5afe575SCameron Smith     &     c_char_'displacement' // char(0),
470e5afe575SCameron Smith     &     c_loc(intfromfile), ithree, dataInt, iotype)
47159599516SKenneth E. Jansen
47259599516SKenneth E. Jansen         nshg2=intfromfile(1)
47359599516SKenneth E. Jansen         ndisp=intfromfile(2)
47459599516SKenneth E. Jansen         lstep=intfromfile(3)
47559599516SKenneth E. Jansen         if(ndisp.ne.nsd) then
47659599516SKenneth E. Jansen            warning='WARNING ndisp not equal nsd'
47759599516SKenneth E. Jansen            write(*,*) warning , ndisp
47859599516SKenneth E. Jansen         endif
47959599516SKenneth E. Jansen         if (nshg2 .ne. nshg)
48059599516SKenneth E. Jansen     &        call error ('restar  ', 'nshg   ', nshg)
48159599516SKenneth E. Jansenc
48259599516SKenneth E. Jansenc.... read the values of primitive variables into uold
48359599516SKenneth E. Jansenc
48459599516SKenneth E. Jansen         allocate( uold(nshg,nsd) )
48559599516SKenneth E. Jansen         allocate( uread(nshg,nsd) )
48659599516SKenneth E. Jansen
48759599516SKenneth E. Jansen         iusiz=nshg*nsd
48859599516SKenneth E. Jansen
489d5d2f64dSCameron Smith         call phio_readdatablock(fhandle,
490bc62cfd4SCameron Smith     &    c_char_'displacement' // char(0),
491bc62cfd4SCameron Smith     &    c_loc(uread), iusiz, dataDbl, iotype)
49259599516SKenneth E. Jansen
49359599516SKenneth E. Jansen         uold(:,1:nsd)=uread(:,1:nsd)
49459599516SKenneth E. Jansen       else
49559599516SKenneth E. Jansen         allocate( uold(nshg,nsd) )
49659599516SKenneth E. Jansen         uold(:,1:nsd) = zero
49759599516SKenneth E. Jansen       endif
49859599516SKenneth E. Jansenc
49959599516SKenneth E. Jansenc.... close c-binary files
50059599516SKenneth E. Jansenc
501d7abaf6cSCameron Smith      call phio_closefile(fhandle)
50259599516SKenneth E. Jansen
50359599516SKenneth E. Jansen      deallocate(xread)
50459599516SKenneth E. Jansen      if ( numpbc > 0 )  then
50559599516SKenneth E. Jansen         deallocate(bcinpread)
50659599516SKenneth E. Jansen         deallocate(ibctmpread)
50759599516SKenneth E. Jansen      endif
50859599516SKenneth E. Jansen      deallocate(iperread)
50959599516SKenneth E. Jansen      if(numpe.gt.1)
51059599516SKenneth E. Jansen     &     deallocate(ilworkread)
51159599516SKenneth E. Jansen      deallocate(nbcread)
51259599516SKenneth E. Jansen
51359599516SKenneth E. Jansen      return
51459599516SKenneth E. Jansen 994  call error ('input   ','opening ', igeomBAK)
51559599516SKenneth E. Jansen 995  call error ('input   ','opening ', igeomBAK)
51659599516SKenneth E. Jansen 997  call error ('input   ','end file', igeomBAK)
51759599516SKenneth E. Jansen 998  call error ('input   ','end file', igeomBAK)
51859599516SKenneth E. Jansen      end
51959599516SKenneth E. Jansenc
52059599516SKenneth E. Jansenc No longer called but kept around in case....
52159599516SKenneth E. Jansenc
52259599516SKenneth E. Jansen      subroutine genpzero(iBC)
52359599516SKenneth E. Jansen
52459599516SKenneth E. Jansen      use pointer_data
52559599516SKenneth E. Jansen      include "common.h"
52659599516SKenneth E. Jansen      integer iBC(nshg)
52759599516SKenneth E. Jansenc
52859599516SKenneth E. Jansenc....  check to see if any of the nodes have a dirichlet pressure
52959599516SKenneth E. Jansenc
53059599516SKenneth E. Jansen      pzero=1
53159599516SKenneth E. Jansen      if (any(btest(iBC,2))) pzero=0
53259599516SKenneth E. Jansen      do iblk = 1, nelblb
53359599516SKenneth E. Jansen         npro = lcblkb(1,iblk+1)-lcblkb(1,iblk)
53459599516SKenneth E. Jansen         do i=1, npro
53559599516SKenneth E. Jansen            iBCB1=miBCB(iblk)%p(i,1)
53659599516SKenneth E. Jansenc
53759599516SKenneth E. Jansenc.... check to see if any of the nodes have a Neumann pressure
53859599516SKenneth E. Jansenc     but not periodic (note that
53959599516SKenneth E. Jansenc
54059599516SKenneth E. Jansen            if(btest(iBCB1,1)) pzero=0
54159599516SKenneth E. Jansen         enddo
54259599516SKenneth E. Jansenc
54359599516SKenneth E. Jansenc.... share results with other processors
54459599516SKenneth E. Jansenc
54559599516SKenneth E. Jansen         pzl=pzero
54659599516SKenneth E. Jansen         if (numpe .gt. 1)
54759599516SKenneth E. Jansen     &        call MPI_ALLREDUCE (pzl, pzero, 1,
54859599516SKenneth E. Jansen     &        MPI_DOUBLE_PRECISION,MPI_MIN, MPI_COMM_WORLD,ierr)
54959599516SKenneth E. Jansen      enddo
55059599516SKenneth E. Jansen      return
55159599516SKenneth E. Jansen      end
552