xref: /phasta/phSolver/common/readnblk.f (revision 02ce253ef78063ddb6a3537928184e2070a2a441)
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
3259599516SKenneth E. Jansen      include "common.h"
33*02ce253eSCameron Smith
349a6935e5SKenneth E. Jansen      real*8, target, allocatable :: xread(:,:), qread(:,:), acread(:,:)
359a6935e5SKenneth E. Jansen      real*8, target, allocatable :: uread(:,:)
369a6935e5SKenneth E. Jansen      real*8, target, allocatable :: BCinpread(:,:)
379a6935e5SKenneth E. Jansen      integer, target, allocatable :: iperread(:), iBCtmpread(:)
389a6935e5SKenneth E. Jansen      integer, target, allocatable :: ilworkread(:), nBCread(:)
3959599516SKenneth E. Jansen      character*10 cname2
4059599516SKenneth E. Jansen      character*8 mach2
4159599516SKenneth E. Jansen      character*30 fmt1
4259599516SKenneth E. Jansen      character*255 fname1,fnamer,fnamelr
4359599516SKenneth E. Jansen      character*255 warning
4459599516SKenneth E. Jansen      integer igeomBAK, ibndc, irstin, ierr
459a6935e5SKenneth E. Jansen      integer, target :: intfromfile(50) ! integers read from headers
4659599516SKenneth E. Jansen      integer :: descriptor, descriptorG, GPID, color, nfiles, nfields
4759599516SKenneth E. Jansen      integer ::  numparts, nppf
4859599516SKenneth E. Jansen      integer :: ierr_io, numprocs, itmp, itmp2
49ade0e30fSCameron Smith      integer :: ignored
5059599516SKenneth E. Jansen      character*255 fname2, temp2
5159599516SKenneth E. Jansen      character*64 temp1
52e5afe575SCameron Smith      type(c_ptr) :: handle
53e5afe575SCameron Smith      character(len=1024) :: dataInt, dataDbl
54e5afe575SCameron Smith      dataInt = c_char_'integer'//c_null_char
55e5afe575SCameron Smith      dataDbl = c_char_'double'//c_null_char
5659599516SKenneth E. Jansenc
5759599516SKenneth E. Jansenc.... determine the step number to start with
5859599516SKenneth E. Jansenc
5959599516SKenneth E. Jansen      open(unit=72,file='numstart.dat',status='old')
6059599516SKenneth E. Jansen      read(72,*) irstart
6159599516SKenneth E. Jansen      close(72)
6259599516SKenneth E. Jansen      lstep=irstart ! in case restart files have no fields
6359599516SKenneth E. Jansen
6459599516SKenneth E. Jansen      fname1='geombc.dat'
6559599516SKenneth E. Jansen      fname1= trim(fname1)  // cname2(myrank+1)
6659599516SKenneth E. Jansen
6759599516SKenneth E. Jansen      itmp=1
6859599516SKenneth E. Jansen      if (irstart .gt. 0) itmp = int(log10(float(irstart)))+1
6959599516SKenneth E. Jansen      write (fmt1,"('(''restart.'',i',i1,',1x)')") itmp
7059599516SKenneth E. Jansen      write (fnamer,fmt1) irstart
7159599516SKenneth E. Jansen      fnamer = trim(fnamer) // cname2(myrank+1)
7259599516SKenneth E. Jansenc
7359599516SKenneth E. Jansenc.... open input files
7459599516SKenneth E. Jansenc.... input the geometry parameters
7559599516SKenneth E. Jansenc
7659599516SKenneth E. Jansen      nfiles = nsynciofiles
7759599516SKenneth E. Jansen      numparts = numpe !This is the common settings. Beware if you try to compute several parts per process
7859599516SKenneth E. Jansen
7959599516SKenneth E. Jansen      itwo=2
8059599516SKenneth E. Jansen      ione=1
8159599516SKenneth E. Jansen      ieleven=11
8259599516SKenneth E. Jansen      itmp = int(log10(float(myrank+1)))+1
8359599516SKenneth E. Jansen
842efdc748SKenneth E. Jansen      if(nsynciofiles.gt.0) then
852efdc748SKenneth E. Jansen        call phio_openfile_read(c_char_'geombc-dat.' // char(0),
862efdc748SKenneth E. Jansen     &                          nfiles, fhandle);
872efdc748SKenneth E. Jansen      else
882efdc748SKenneth E. Jansen        call phio_openfile_read(c_char_'geombc.dat.' // char(0),
892efdc748SKenneth E. Jansen     &                          nfiles, fhandle);
902efdc748SKenneth E. Jansen      endif
9159599516SKenneth E. Jansen
92d5d2f64dSCameron Smith      call phio_readheader(fhandle,c_char_'number of nodes' // char(0),
93e5afe575SCameron Smith     & c_loc(numnp),ione, dataInt, iotype)
9459599516SKenneth E. Jansen
95d5d2f64dSCameron Smith      call phio_readheader(fhandle,c_char_'number of modes' // char(0),
96e5afe575SCameron Smith     & c_loc(nshg),ione, dataInt, iotype)
9759599516SKenneth E. Jansen
98d5d2f64dSCameron Smith      call phio_readheader(fhandle,
99e5afe575SCameron Smith     &  c_char_'number of interior elements' // char(0),
100e5afe575SCameron Smith     &  c_loc(numel),ione, dataInt, iotype)
10159599516SKenneth E. Jansen
102d5d2f64dSCameron Smith      call phio_readheader(fhandle,
103e5afe575SCameron Smith     &  c_char_'number of boundary elements' // char(0),
104e5afe575SCameron Smith     &  c_loc(numelb),ione, dataInt, iotype)
10559599516SKenneth E. Jansen
106d5d2f64dSCameron Smith      call phio_readheader(fhandle,
107e5afe575SCameron Smith     &  c_char_'maximum number of element nodes' // char(0),
108e5afe575SCameron Smith     &  c_loc(nen),ione, dataInt, iotype)
10959599516SKenneth E. Jansen
110d5d2f64dSCameron Smith      call phio_readheader(fhandle,
111e5afe575SCameron Smith     &  c_char_'number of interior tpblocks' // char(0),
112e5afe575SCameron Smith     &  c_loc(nelblk),ione, dataInt, iotype)
11359599516SKenneth E. Jansen
114d5d2f64dSCameron Smith      call phio_readheader(fhandle,
115e5afe575SCameron Smith     & c_char_'number of boundary tpblocks' // char(0),
116e5afe575SCameron Smith     & c_loc(nelblb),ione, dataInt, iotype)
11759599516SKenneth E. Jansen
118d5d2f64dSCameron Smith      call phio_readheader(fhandle,
119e5afe575SCameron Smith     & c_char_'number of nodes with Dirichlet BCs' // char(0),
120e5afe575SCameron Smith     & c_loc(numpbc),ione, dataInt, iotype)
12159599516SKenneth E. Jansen
122d5d2f64dSCameron Smith      call phio_readheader(fhandle,
123e5afe575SCameron Smith     & c_char_'number of shape functions' // char(0),
124e5afe575SCameron Smith     & c_loc(ntopsh),ione, dataInt, iotype)
12559599516SKenneth E. Jansenc
12659599516SKenneth E. Jansenc.... calculate the maximum number of boundary element nodes
12759599516SKenneth E. Jansenc
12859599516SKenneth E. Jansen      nenb = 0
12959599516SKenneth E. Jansen      do i = 1, melCat
13059599516SKenneth E. Jansen         if (nen .eq. nenCat(i,nsd)) nenb = max(nenCat(i,nsd-1), nenb)
13159599516SKenneth E. Jansen      enddo
13259599516SKenneth E. Jansen      if (myrank == master) then
13359599516SKenneth E. Jansen         if (nenb .eq. 0) call error ('input   ','nen     ',nen)
13459599516SKenneth E. Jansen      endif
13559599516SKenneth E. Jansenc
13659599516SKenneth E. Jansenc.... setup some useful constants
13759599516SKenneth E. Jansenc
13859599516SKenneth E. Jansen      I3nsd  = nsd / 3          ! nsd=3 integer flag
13959599516SKenneth E. Jansen      E3nsd  = float(I3nsd)     ! nsd=3 real    flag
14059599516SKenneth E. Jansen      if(matflg(1,1).lt.0) then
14159599516SKenneth E. Jansen         nflow = nsd + 1
14259599516SKenneth E. Jansen      else
14359599516SKenneth E. Jansen         nflow = nsd + 2
14459599516SKenneth E. Jansen      endif
14559599516SKenneth E. Jansen      ndof   = nsd + 2
14659599516SKenneth E. Jansen      nsclr=impl(1)/100
14759599516SKenneth E. Jansen      ndof=ndof+nsclr           ! number of sclr transport equations to solve
14859599516SKenneth E. Jansen
14959599516SKenneth E. Jansen      ndofBC = ndof + I3nsd     ! dimension of BC array
15059599516SKenneth E. Jansen      ndiBCB = 2                ! dimension of iBCB array
15159599516SKenneth E. Jansen      ndBCB  = ndof + 1         ! dimension of BCB array
15259599516SKenneth E. Jansen      nsymdf = (ndof*(ndof + 1)) / 2 ! symm. d.o.f.'s
15359599516SKenneth E. Jansenc
15459599516SKenneth E. Jansenc.... ----------------------> Communication tasks <--------------------
15559599516SKenneth E. Jansenc
15659599516SKenneth E. Jansen      if(numpe > 1) then
157d5d2f64dSCameron Smith         call phio_readheader(fhandle,
158e5afe575SCameron Smith     &    c_char_'size of ilwork array' // char(0),
159e5afe575SCameron Smith     &    c_loc(nlwork),ione, dataInt, iotype)
16059599516SKenneth E. Jansen
161d5d2f64dSCameron Smith         call phio_readheader(fhandle,
162e5afe575SCameron Smith     &    c_char_'ilwork' //char(0),
163e5afe575SCameron Smith     &    c_loc(nlwork),ione, dataInt, iotype)
16459599516SKenneth E. Jansen
16559599516SKenneth E. Jansen         allocate( point2ilwork(nlwork) )
16659599516SKenneth E. Jansen         allocate( ilworkread(nlwork) )
167d5d2f64dSCameron Smith         call phio_readdatablock(fhandle, c_char_'ilwork' // char(0),
168bc62cfd4SCameron Smith     &      c_loc(ilworkread), nlwork, dataInt , iotype)
16959599516SKenneth E. Jansen
17059599516SKenneth E. Jansen         point2ilwork = ilworkread
17159599516SKenneth E. Jansen         call ctypes (point2ilwork)
17259599516SKenneth E. Jansen      else
17359599516SKenneth E. Jansen           nlwork=1
17459599516SKenneth E. Jansen           allocate( point2ilwork(1))
17559599516SKenneth E. Jansen           nshg0 = nshg
17659599516SKenneth E. Jansen      endif
17759599516SKenneth E. Jansen
17859599516SKenneth E. Jansen      itwo=2
17959599516SKenneth E. Jansen
180d5d2f64dSCameron Smith      call phio_readheader(fhandle,
181e5afe575SCameron Smith     & c_char_'co-ordinates' // char(0),
182e5afe575SCameron Smith     & c_loc(intfromfile),itwo, dataDbl, iotype)
18359599516SKenneth E. Jansen      numnp=intfromfile(1)
18459599516SKenneth E. Jansen      allocate( point2x(numnp,nsd) )
18559599516SKenneth E. Jansen      allocate( xread(numnp,nsd) )
18659599516SKenneth E. Jansen      ixsiz=numnp*nsd
187d5d2f64dSCameron Smith      call phio_readdatablock(fhandle,
188bc62cfd4SCameron Smith     & c_char_'co-ordinates' // char(0),
189bc62cfd4SCameron Smith     & c_loc(xread),ixsiz, dataDbl, iotype)
19059599516SKenneth E. Jansen      point2x = xread
19159599516SKenneth E. Jansenc
19259599516SKenneth E. Jansenc.... read in and block out the connectivity
19359599516SKenneth E. Jansenc
19459599516SKenneth E. Jansen      call genblk (IBKSIZ)
19559599516SKenneth E. Jansenc
19659599516SKenneth E. Jansenc.... read the boundary condition mapping array
19759599516SKenneth E. Jansenc
19859599516SKenneth E. Jansen      ione=1
199d5d2f64dSCameron Smith      call phio_readheader(fhandle,
200e5afe575SCameron Smith     & c_char_'bc mapping array' // char(0),
201e5afe575SCameron Smith     & c_loc(nshg),ione, dataInt, iotype)
20259599516SKenneth E. Jansen
20359599516SKenneth E. Jansen      allocate( nBC(nshg) )
20459599516SKenneth E. Jansen
20559599516SKenneth E. Jansen      allocate( nBCread(nshg) )
20659599516SKenneth E. Jansen
207d5d2f64dSCameron Smith      call phio_readdatablock(fhandle,
208bc62cfd4SCameron Smith     & c_char_'bc mapping array' // char(0),
209bc62cfd4SCameron Smith     & c_loc(nBCread), nshg, dataInt, iotype)
21059599516SKenneth E. Jansen
21159599516SKenneth E. Jansen      nBC=nBCread
21259599516SKenneth E. Jansenc
21359599516SKenneth E. Jansenc.... read the temporary iBC array
21459599516SKenneth E. Jansenc
21559599516SKenneth E. Jansen      ione=1
216d5d2f64dSCameron Smith      call phio_readheader(fhandle,
217e5afe575SCameron Smith     & c_char_'bc codes array' // char(0),
218e5afe575SCameron Smith     & c_loc(numpbc),ione, dataInt, iotype)
21959599516SKenneth E. Jansen
22059599516SKenneth E. Jansen      if ( numpbc > 0 ) then
22159599516SKenneth E. Jansen        allocate( iBCtmp(numpbc) )
22259599516SKenneth E. Jansen        allocate( iBCtmpread(numpbc) )
22359599516SKenneth E. Jansen      else
22459599516SKenneth E. Jansen        allocate( iBCtmp(1) )
22559599516SKenneth E. Jansen        allocate( iBCtmpread(1) )
22659599516SKenneth E. Jansen      endif
227d5d2f64dSCameron Smith      call phio_readdatablock(fhandle,
228bc62cfd4SCameron Smith     & c_char_'bc codes array' // char(0),
229bc62cfd4SCameron Smith     & c_loc(iBCtmpread), numpbc, dataInt, iotype)
23059599516SKenneth E. Jansen
23159599516SKenneth E. Jansen      if ( numpbc > 0 ) then
23259599516SKenneth E. Jansen         iBCtmp=iBCtmpread
23359599516SKenneth E. Jansen      else  ! sometimes a partition has no BC's
23459599516SKenneth E. Jansen         deallocate( iBCtmpread)
23559599516SKenneth E. Jansen         iBCtmp=0
23659599516SKenneth E. Jansen      endif
23759599516SKenneth E. Jansenc
23859599516SKenneth E. Jansenc.... read boundary condition data
23959599516SKenneth E. Jansenc
24059599516SKenneth E. Jansen      ione=1
241d5d2f64dSCameron Smith      call phio_readheader(fhandle,
242e5afe575SCameron Smith     & c_char_'boundary condition array' // char(0),
243e5afe575SCameron Smith     & c_loc(intfromfile),ione, dataDbl, iotype)
24459599516SKenneth E. Jansen
24559599516SKenneth E. Jansen      if ( numpbc > 0 ) then
24659599516SKenneth E. Jansen         allocate( BCinp(numpbc,ndof+7) )
24759599516SKenneth E. Jansen         nsecondrank=intfromfile(1)/numpbc
24859599516SKenneth E. Jansen         allocate( BCinpread(numpbc,nsecondrank) )
24959599516SKenneth E. Jansen         iBCinpsiz=intfromfile(1)
25059599516SKenneth E. Jansen      else
25159599516SKenneth E. Jansen         allocate( BCinp(1,ndof+7) )
25259599516SKenneth E. Jansen         allocate( BCinpread(0,0) ) !dummy
25359599516SKenneth E. Jansen         iBCinpsiz=intfromfile(1)
25459599516SKenneth E. Jansen      endif
25559599516SKenneth E. Jansen
256d5d2f64dSCameron Smith      call phio_readdatablock(fhandle,
257bc62cfd4SCameron Smith     & c_char_'boundary condition array' // char(0),
258bc62cfd4SCameron Smith     & c_loc(BCinpread), iBCinpsiz, dataDbl, iotype)
25959599516SKenneth E. Jansen
26059599516SKenneth E. Jansen      if ( numpbc > 0 ) then
26159599516SKenneth E. Jansen         BCinp(:,1:(ndof+7))=BCinpread(:,1:(ndof+7))
26259599516SKenneth E. Jansen      else  ! sometimes a partition has no BC's
26359599516SKenneth E. Jansen         deallocate(BCinpread)
26459599516SKenneth E. Jansen         BCinp=0
26559599516SKenneth E. Jansen      endif
26659599516SKenneth E. Jansenc
26759599516SKenneth E. Jansenc.... read periodic boundary conditions
26859599516SKenneth E. Jansenc
26959599516SKenneth E. Jansen      ione=1
270d5d2f64dSCameron Smith      call phio_readheader(fhandle,
271e5afe575SCameron Smith     & c_char_'periodic masters array' // char(0),
272e5afe575SCameron Smith     & c_loc(nshg), ione, dataInt, iotype)
27359599516SKenneth E. Jansen      allocate( point2iper(nshg) )
27459599516SKenneth E. Jansen      allocate( iperread(nshg) )
275d5d2f64dSCameron Smith      call phio_readdatablock(fhandle,
276bc62cfd4SCameron Smith     & c_char_'periodic masters array' // char(0),
277bc62cfd4SCameron Smith     & c_loc(iperread), nshg, dataInt, iotype)
27859599516SKenneth E. Jansen      point2iper=iperread
27959599516SKenneth E. Jansenc
28059599516SKenneth E. Jansenc.... generate the boundary element blocks
28159599516SKenneth E. Jansenc
28259599516SKenneth E. Jansen      call genbkb (ibksiz)
28359599516SKenneth E. Jansenc
28459599516SKenneth E. Jansenc  Read in the nsons and ifath arrays if needed
28559599516SKenneth E. Jansenc
28659599516SKenneth E. Jansenc  There is a fundamental shift in the meaning of ifath based on whether
28759599516SKenneth E. Jansenc  there exist homogenous directions in the flow.
28859599516SKenneth E. Jansenc
28959599516SKenneth E. Jansenc  HOMOGENOUS DIRECTIONS EXIST:  Here nfath is the number of inhomogenous
29059599516SKenneth E. Jansenc  points in the TOTAL mesh.  That is to say that each partition keeps a
29159599516SKenneth E. Jansenc  link to  ALL inhomogenous points.  This link is furthermore not to the
29259599516SKenneth E. Jansenc  sms numbering but to the original structured grid numbering.  These
29359599516SKenneth E. Jansenc  inhomogenous points are thought of as fathers, with their sons being all
29459599516SKenneth E. Jansenc  the points in the homogenous directions that have this father's
29559599516SKenneth E. Jansenc  inhomogeneity.  The array ifath takes as an arguement the sms numbering
29659599516SKenneth E. Jansenc  and returns as a result the father.
29759599516SKenneth E. Jansenc
29859599516SKenneth E. Jansenc  In this case nsons is the number of sons that each father has and ifath
29959599516SKenneth E. Jansenc  is an array which tells the
30059599516SKenneth E. Jansenc
30159599516SKenneth E. Jansenc  NO HOMOGENOUS DIRECTIONS.  In this case the mesh would grow to rapidly
30259599516SKenneth E. Jansenc  if we followed the above strategy since every partition would index its
30359599516SKenneth E. Jansenc  points to the ENTIRE mesh.  Furthermore, there would never be a need
30459599516SKenneth E. Jansenc  to average to a node off processor since there is no spatial averaging.
30559599516SKenneth E. Jansenc  Therefore, to properly account for this case we must recognize it and
30659599516SKenneth E. Jansenc  inerrupt certain actions (i.e. assembly of the average across partitions).
30759599516SKenneth E. Jansenc  This case is easily identified by noting that maxval(nsons) =1 (i.e. no
30859599516SKenneth E. Jansenc  father has any sons).  Reiterating to be clear, in this case ifath does
30959599516SKenneth E. Jansenc  not point to a global numbering but instead just points to itself.
31059599516SKenneth E. Jansenc
31159599516SKenneth E. Jansen      nfath=1  ! some architectures choke on a zero or undeclared
31259599516SKenneth E. Jansen                 ! dimension variable.  This sets it to a safe, small value.
31359599516SKenneth E. Jansen      if(((iLES .lt. 20) .and. (iLES.gt.0))
31459599516SKenneth E. Jansen     &                   .or. (itwmod.gt.0)  ) then ! don't forget same
31559599516SKenneth E. Jansen                                                    ! conditional in proces.f
31659599516SKenneth E. Jansen                                                    ! needed for alloc
31759599516SKenneth E. Jansen         ione=1
31859599516SKenneth E. Jansen         if(nohomog.gt.0) then
319d5d2f64dSCameron Smith            call phio_readheader(fhandle,
320e5afe575SCameron Smith     &       c_char_'number of father-nodes' // char(0),
321e5afe575SCameron Smith     &       c_loc(nfath), ione, dataInt, iotype)
32259599516SKenneth E. Jansen
323d5d2f64dSCameron Smith            call phio_readheader(fhandle,
324e5afe575SCameron Smith     &       c_char_'number of son-nodes for each father' // char(0),
325e5afe575SCameron Smith     &       c_loc(nfath), ione, dataInt, iotype)
32659599516SKenneth E. Jansen
32759599516SKenneth E. Jansen            allocate (point2nsons(nfath))
32859599516SKenneth E. Jansen
329d5d2f64dSCameron Smith            call phio_readdatablock(fhandle,
330bc62cfd4SCameron Smith     &       c_char_'number of son-nodes for each father' // char(0),
331bc62cfd4SCameron Smith     &       c_loc(point2nsons),nfath, dataInt, iotype)
33259599516SKenneth E. Jansen
333d5d2f64dSCameron Smith            call phio_readheader(fhandle,
334e5afe575SCameron Smith     &       c_char_'keyword ifath' // char(0),
335e5afe575SCameron Smith     &       c_loc(nshg), ione, dataInt, iotype);
33659599516SKenneth E. Jansen
33759599516SKenneth E. Jansen            allocate (point2ifath(nshg))
33859599516SKenneth E. Jansen
339d5d2f64dSCameron Smith            call phio_readdatablock(fhandle,
340bc62cfd4SCameron Smith     &       c_char_'keyword ifath' // char(0),
341bc62cfd4SCameron Smith     &       c_loc(point2ifath), nshg, dataInt, iotype)
34259599516SKenneth E. Jansen
34359599516SKenneth E. Jansen            nsonmax=maxval(point2nsons)
34459599516SKenneth E. Jansen         else  ! this is the case where there is no homogeneity
34559599516SKenneth E. Jansen               ! therefore ever node is a father (too itself).  sonfath
34659599516SKenneth E. Jansen               ! (a routine in NSpre) will set this up but this gives
34759599516SKenneth E. Jansen               ! you an option to avoid that.
34859599516SKenneth E. Jansen            nfath=nshg
34959599516SKenneth E. Jansen            allocate (point2nsons(nfath))
35059599516SKenneth E. Jansen            point2nsons=1
35159599516SKenneth E. Jansen            allocate (point2ifath(nshg))
35259599516SKenneth E. Jansen            do i=1,nshg
35359599516SKenneth E. Jansen               point2ifath(i)=i
35459599516SKenneth E. Jansen            enddo
35559599516SKenneth E. Jansen            nsonmax=1
35659599516SKenneth E. Jansen         endif
35759599516SKenneth E. Jansen      else
35859599516SKenneth E. Jansen         allocate (point2nsons(1))
35959599516SKenneth E. Jansen         allocate (point2ifath(1))
36059599516SKenneth E. Jansen      endif
36159599516SKenneth E. Jansen
362d5d2f64dSCameron Smith      call phio_closefile_read(fhandle);
36359599516SKenneth E. Jansenc.... Read restart files
3642efdc748SKenneth E. Jansen      if(nsynciofiles.gt.0) then
36536adee64SCameron Smith        fnamer = c_char_"restart-dat."//c_null_char
3662efdc748SKenneth E. Jansen      else
3672efdc748SKenneth E. Jansen        fnamer = c_char_"restart."//c_null_char
3682efdc748SKenneth E. Jansen      endif
3692efdc748SKenneth E. Jansen
37036adee64SCameron Smith      call phio_appendStep(fnamer, irstart)
371d5d2f64dSCameron Smith      call phio_openfile_read(fnamer, nfiles, fhandle)
37259599516SKenneth E. Jansen
37359599516SKenneth E. Jansen      ithree=3
37459599516SKenneth E. Jansen
37559599516SKenneth E. Jansen      itmp = int(log10(float(myrank+1)))+1
37659599516SKenneth E. Jansen
37759599516SKenneth E. Jansen      intfromfile=0
378d5d2f64dSCameron Smith      call phio_readheader(fhandle,
379e5afe575SCameron Smith     & c_char_'solution' // char(0),
380e5afe575SCameron Smith     & c_loc(intfromfile), ithree, dataInt, iotype)
38159599516SKenneth E. Jansenc
38259599516SKenneth E. Jansenc.... read the values of primitive variables into q
38359599516SKenneth E. Jansenc
38459599516SKenneth E. Jansen      allocate( qold(nshg,ndof) )
38559599516SKenneth E. Jansen      if(intfromfile(1).ne.0) then
38659599516SKenneth E. Jansen         nshg2=intfromfile(1)
38759599516SKenneth E. Jansen         ndof2=intfromfile(2)
38859599516SKenneth E. Jansen         lstep=intfromfile(3)
38959599516SKenneth E. Jansen         if(ndof2.ne.ndof) then
39059599516SKenneth E. Jansen
39159599516SKenneth E. Jansen         endif
39259599516SKenneth E. Jansen        if (nshg2 .ne. nshg)
39359599516SKenneth E. Jansen     &        call error ('restar  ', 'nshg   ', nshg)
39459599516SKenneth E. Jansen         allocate( qread(nshg,ndof2) )
39559599516SKenneth E. Jansen         iqsiz=nshg*ndof2
396d5d2f64dSCameron Smith         call phio_readdatablock(fhandle,
397bc62cfd4SCameron Smith     &    c_char_'solution' // char(0),
398bc62cfd4SCameron Smith     &    c_loc(qread),iqsiz, dataDbl,iotype)
39959599516SKenneth E. Jansen         qold(:,1:ndof)=qread(:,1:ndof)
40059599516SKenneth E. Jansen         deallocate(qread)
40159599516SKenneth E. Jansen      else
40259599516SKenneth E. Jansen         if (myrank.eq.master) then
40359599516SKenneth E. Jansen            if (matflg(1,1).eq.0) then ! compressible
40459599516SKenneth E. Jansen               warning='Solution is set to zero (with p and T to one)'
40559599516SKenneth E. Jansen            else
40659599516SKenneth E. Jansen               warning='Solution is set to zero'
40759599516SKenneth E. Jansen            endif
40859599516SKenneth E. Jansen            write(*,*) warning
40959599516SKenneth E. Jansen         endif
41059599516SKenneth E. Jansen         qold=zero
41159599516SKenneth E. Jansen         if (matflg(1,1).eq.0) then ! compressible
41259599516SKenneth E. Jansen            qold(:,1)=one ! avoid zero pressure
41359599516SKenneth E. Jansen            qold(:,nflow)=one ! avoid zero temperature
41459599516SKenneth E. Jansen         endif
41559599516SKenneth E. Jansen      endif
41659599516SKenneth E. Jansen
41759599516SKenneth E. Jansen      intfromfile=0
418d5d2f64dSCameron Smith      call phio_readheader(fhandle,
419e5afe575SCameron Smith     & c_char_'time derivative of solution' // char(0),
420e5afe575SCameron Smith     & c_loc(intfromfile), ithree, dataInt, iotype)
42159599516SKenneth E. Jansen      allocate( acold(nshg,ndof) )
42259599516SKenneth E. Jansen      if(intfromfile(1).ne.0) then
42359599516SKenneth E. Jansen         nshg2=intfromfile(1)
42459599516SKenneth E. Jansen         ndof2=intfromfile(2)
42559599516SKenneth E. Jansen         lstep=intfromfile(3)
42659599516SKenneth E. Jansen
42759599516SKenneth E. Jansen         if (nshg2 .ne. nshg)
42859599516SKenneth E. Jansen     &        call error ('restar  ', 'nshg   ', nshg)
42959599516SKenneth E. Jansen         allocate( acread(nshg,ndof2) )
43059599516SKenneth E. Jansen         acread=zero
43159599516SKenneth E. Jansen         iacsiz=nshg*ndof2
432d5d2f64dSCameron Smith         call phio_readdatablock(fhandle,
433bc62cfd4SCameron Smith     &    c_char_'time derivative of solution' // char(0),
434bc62cfd4SCameron Smith     &    c_loc(acread), iacsiz, dataDbl,iotype)
43559599516SKenneth E. Jansen         acold(:,1:ndof)=acread(:,1:ndof)
43659599516SKenneth E. Jansen         deallocate(acread)
43759599516SKenneth E. Jansen      else
43859599516SKenneth E. Jansen         if (myrank.eq.master) then
43959599516SKenneth E. Jansen            warning='Time derivative of solution is set to zero (SAFE)'
44059599516SKenneth E. Jansen            write(*,*) warning
44159599516SKenneth E. Jansen         endif
44259599516SKenneth E. Jansen         acold=zero
44359599516SKenneth E. Jansen      endif
44459599516SKenneth E. Jansencc
44559599516SKenneth E. Jansencc.... read the header and check it against the run data
44659599516SKenneth E. Jansencc
44759599516SKenneth E. Jansen      if (ideformwall.eq.1) then
44859599516SKenneth E. Jansen
44959599516SKenneth E. Jansen          intfromfile=0
450d5d2f64dSCameron Smith          call phio_readheader(fhandle,
451e5afe575SCameron Smith     &     c_char_'displacement' // char(0),
452e5afe575SCameron Smith     &     c_loc(intfromfile), ithree, dataInt, iotype)
45359599516SKenneth E. Jansen
45459599516SKenneth E. Jansen         nshg2=intfromfile(1)
45559599516SKenneth E. Jansen         ndisp=intfromfile(2)
45659599516SKenneth E. Jansen         lstep=intfromfile(3)
45759599516SKenneth E. Jansen         if(ndisp.ne.nsd) then
45859599516SKenneth E. Jansen            warning='WARNING ndisp not equal nsd'
45959599516SKenneth E. Jansen            write(*,*) warning , ndisp
46059599516SKenneth E. Jansen         endif
46159599516SKenneth E. Jansen         if (nshg2 .ne. nshg)
46259599516SKenneth E. Jansen     &        call error ('restar  ', 'nshg   ', nshg)
46359599516SKenneth E. Jansenc
46459599516SKenneth E. Jansenc.... read the values of primitive variables into uold
46559599516SKenneth E. Jansenc
46659599516SKenneth E. Jansen         allocate( uold(nshg,nsd) )
46759599516SKenneth E. Jansen         allocate( uread(nshg,nsd) )
46859599516SKenneth E. Jansen
46959599516SKenneth E. Jansen         iusiz=nshg*nsd
47059599516SKenneth E. Jansen
471d5d2f64dSCameron Smith         call phio_readdatablock(fhandle,
472bc62cfd4SCameron Smith     &    c_char_'displacement' // char(0),
473bc62cfd4SCameron Smith     &    c_loc(uread), iusiz, dataDbl, iotype)
47459599516SKenneth E. Jansen
47559599516SKenneth E. Jansen         uold(:,1:nsd)=uread(:,1:nsd)
47659599516SKenneth E. Jansen       else
47759599516SKenneth E. Jansen         allocate( uold(nshg,nsd) )
47859599516SKenneth E. Jansen         uold(:,1:nsd) = zero
47959599516SKenneth E. Jansen       endif
48059599516SKenneth E. Jansenc
48159599516SKenneth E. Jansenc.... close c-binary files
48259599516SKenneth E. Jansenc
483d5d2f64dSCameron Smith      call phio_closefile_read(fhandle)
48459599516SKenneth E. Jansen
48559599516SKenneth E. Jansen      deallocate(xread)
48659599516SKenneth E. Jansen      if ( numpbc > 0 )  then
48759599516SKenneth E. Jansen         deallocate(bcinpread)
48859599516SKenneth E. Jansen         deallocate(ibctmpread)
48959599516SKenneth E. Jansen      endif
49059599516SKenneth E. Jansen      deallocate(iperread)
49159599516SKenneth E. Jansen      if(numpe.gt.1)
49259599516SKenneth E. Jansen     &     deallocate(ilworkread)
49359599516SKenneth E. Jansen      deallocate(nbcread)
49459599516SKenneth E. Jansen
49559599516SKenneth E. Jansen      return
49659599516SKenneth E. Jansen 994  call error ('input   ','opening ', igeomBAK)
49759599516SKenneth E. Jansen 995  call error ('input   ','opening ', igeomBAK)
49859599516SKenneth E. Jansen 997  call error ('input   ','end file', igeomBAK)
49959599516SKenneth E. Jansen 998  call error ('input   ','end file', igeomBAK)
50059599516SKenneth E. Jansen      end
50159599516SKenneth E. Jansenc
50259599516SKenneth E. Jansenc No longer called but kept around in case....
50359599516SKenneth E. Jansenc
50459599516SKenneth E. Jansen      subroutine genpzero(iBC)
50559599516SKenneth E. Jansen
50659599516SKenneth E. Jansen      use pointer_data
50759599516SKenneth E. Jansen      include "common.h"
50859599516SKenneth E. Jansen      integer iBC(nshg)
50959599516SKenneth E. Jansenc
51059599516SKenneth E. Jansenc....  check to see if any of the nodes have a dirichlet pressure
51159599516SKenneth E. Jansenc
51259599516SKenneth E. Jansen      pzero=1
51359599516SKenneth E. Jansen      if (any(btest(iBC,2))) pzero=0
51459599516SKenneth E. Jansen      do iblk = 1, nelblb
51559599516SKenneth E. Jansen         npro = lcblkb(1,iblk+1)-lcblkb(1,iblk)
51659599516SKenneth E. Jansen         do i=1, npro
51759599516SKenneth E. Jansen            iBCB1=miBCB(iblk)%p(i,1)
51859599516SKenneth E. Jansenc
51959599516SKenneth E. Jansenc.... check to see if any of the nodes have a Neumann pressure
52059599516SKenneth E. Jansenc     but not periodic (note that
52159599516SKenneth E. Jansenc
52259599516SKenneth E. Jansen            if(btest(iBCB1,1)) pzero=0
52359599516SKenneth E. Jansen         enddo
52459599516SKenneth E. Jansenc
52559599516SKenneth E. Jansenc.... share results with other processors
52659599516SKenneth E. Jansenc
52759599516SKenneth E. Jansen         pzl=pzero
52859599516SKenneth E. Jansen         if (numpe .gt. 1)
52959599516SKenneth E. Jansen     &        call MPI_ALLREDUCE (pzl, pzero, 1,
53059599516SKenneth E. Jansen     &        MPI_DOUBLE_PRECISION,MPI_MIN, MPI_COMM_WORLD,ierr)
53159599516SKenneth E. Jansen      enddo
53259599516SKenneth E. Jansen      return
53359599516SKenneth E. Jansen      end
534