xref: /phasta/phSolver/common/readnblk.f (revision 9a6935e58d1bec7fbe804a813bfa0382dfd762c1)
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
1259599516SKenneth E. Jansen
1359599516SKenneth E. Jansen      module readarrays
1459599516SKenneth E. Jansen
1559599516SKenneth E. Jansen      real*8, allocatable :: point2x(:,:)
1659599516SKenneth E. Jansen      real*8, allocatable :: qold(:,:)
1759599516SKenneth E. Jansen      real*8, allocatable :: uold(:,:)
1859599516SKenneth E. Jansen      real*8, allocatable :: acold(:,:)
1959599516SKenneth E. Jansen      integer, allocatable :: iBCtmp(:)
2059599516SKenneth E. Jansen      real*8, allocatable :: BCinp(:,:)
2159599516SKenneth E. Jansen
2259599516SKenneth E. Jansen      integer, allocatable :: point2ilwork(:)
2359599516SKenneth E. Jansen      integer, allocatable :: nBC(:)
2459599516SKenneth E. Jansen      integer, allocatable :: point2iper(:)
25*9a6935e5SKenneth E. Jansen      integer, target, allocatable :: point2ifath(:)
26*9a6935e5SKenneth E. Jansen      integer, target, allocatable :: point2nsons(:)
2759599516SKenneth E. Jansen
2859599516SKenneth E. Jansen      end module
2959599516SKenneth E. Jansen
3059599516SKenneth E. Jansen
3159599516SKenneth E. Jansen
3259599516SKenneth E. Jansen      subroutine readnblk
3359599516SKenneth E. Jansenc
34e5afe575SCameron Smith      use iso_c_binding
3559599516SKenneth E. Jansen      use readarrays
36e5afe575SCameron Smith      use phio
3759599516SKenneth E. Jansen      include "common.h"
3859599516SKenneth E. Jansenc
39*9a6935e5SKenneth E. Jansen      real*8, target, allocatable :: xread(:,:), qread(:,:), acread(:,:)
40*9a6935e5SKenneth E. Jansen      real*8, target, allocatable :: uread(:,:)
41*9a6935e5SKenneth E. Jansen      real*8, target, allocatable :: BCinpread(:,:)
42*9a6935e5SKenneth E. Jansen      integer, target, allocatable :: iperread(:), iBCtmpread(:)
43*9a6935e5SKenneth E. Jansen      integer, target, allocatable :: ilworkread(:), nBCread(:)
4459599516SKenneth E. Jansen      character*10 cname2
4559599516SKenneth E. Jansen      character*8 mach2
4659599516SKenneth E. Jansen!MR CHANGE
4759599516SKenneth E. Jansen!      character*20 fmt1
4859599516SKenneth E. Jansen      character*30 fmt1
4959599516SKenneth E. Jansen!MR CHANGE END
5059599516SKenneth E. Jansen      character*255 fname1,fnamer,fnamelr
5159599516SKenneth E. Jansen      character*255 warning
5259599516SKenneth E. Jansen      integer igeomBAK, ibndc, irstin, ierr
53*9a6935e5SKenneth E. Jansen      integer, target :: intfromfile(50) ! integers read from headers
5459599516SKenneth E. Jansen
5559599516SKenneth E. Jansencccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
5659599516SKenneth E. Jansenccccccccccccccccccccccccccccccccccccccccc New PhastaIO Definition Part ccccccccccccccccccccccccccccccccccccccccc
5759599516SKenneth E. Jansen
5859599516SKenneth E. Jansen      integer :: descriptor, descriptorG, GPID, color, nfiles, nfields
5959599516SKenneth E. Jansen      integer ::  numparts, nppf
6059599516SKenneth E. Jansen      integer :: ierr_io, numprocs, itmp, itmp2
61ade0e30fSCameron Smith      integer :: ignored
6259599516SKenneth E. Jansen      character*255 fname2, temp2
6359599516SKenneth E. Jansen      character*64 temp1
6459599516SKenneth E. Jansen
6559599516SKenneth E. Jansencccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
6659599516SKenneth E. Jansencccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
6759599516SKenneth E. Jansen
68e5afe575SCameron Smith      type(c_ptr) :: handle
69e5afe575SCameron Smith      character(len=1024) :: dataInt, dataDbl
70e5afe575SCameron Smith      dataInt = c_char_'integer'//c_null_char
71e5afe575SCameron Smith      dataDbl = c_char_'double'//c_null_char
7259599516SKenneth E. Jansen
7359599516SKenneth E. Jansenc
7459599516SKenneth E. Jansenc.... determine the step number to start with
7559599516SKenneth E. Jansenc
7659599516SKenneth E. Jansen      open(unit=72,file='numstart.dat',status='old')
7759599516SKenneth E. Jansen      read(72,*) irstart
7859599516SKenneth E. Jansen      close(72)
7959599516SKenneth E. Jansen      lstep=irstart ! in case restart files have no fields
8059599516SKenneth E. Jansen
8159599516SKenneth E. Jansenc
8259599516SKenneth E. Jansen      fname1='geombc.dat'
8359599516SKenneth E. Jansen      fname1= trim(fname1)  // cname2(myrank+1)
8459599516SKenneth E. Jansenc      fnamelr='restart.latest'
8559599516SKenneth E. Jansen
8659599516SKenneth E. Jansen      itmp=1
8759599516SKenneth E. Jansen      if (irstart .gt. 0) itmp = int(log10(float(irstart)))+1
8859599516SKenneth E. Jansen      write (fmt1,"('(''restart.'',i',i1,',1x)')") itmp
8959599516SKenneth E. Jansen      write (fnamer,fmt1) irstart
9059599516SKenneth E. Jansen      fnamer = trim(fnamer) // cname2(myrank+1)
9159599516SKenneth E. Jansenc      fnamelr = trim(fnamelr) // cname2(myrank+1)
9259599516SKenneth E. Jansen
9359599516SKenneth E. Jansenc
9459599516SKenneth E. Jansenc.... open input files
9559599516SKenneth E. Jansenc
9659599516SKenneth E. Jansen
9759599516SKenneth E. Jansen
9859599516SKenneth E. Jansenc      call openfile(  fname1,  'read?', igeomBAK );
9959599516SKenneth E. Jansen
10059599516SKenneth E. Jansen
10159599516SKenneth E. Jansenc
10259599516SKenneth E. Jansenc.... try opening restart.latest.proc before trying restart.stepno.proc
10359599516SKenneth E. Jansenc
10459599516SKenneth E. Jansenc      call openfile(  fnamelr,  'read?', irstin );
10559599516SKenneth E. Jansenc      if ( irstin .eq. 0 )
10659599516SKenneth E. Jansen
10759599516SKenneth E. Jansen!MR CHANGE
10859599516SKenneth E. Jansen!       call openfile( fnamer, 'read?', irstin );
10959599516SKenneth E. Jansen!MR CHANGE END
11059599516SKenneth E. Jansen
11159599516SKenneth E. Jansen! either one will work
11259599516SKenneth E. Jansenc
11359599516SKenneth E. Jansenc.... input the geometry parameters
11459599516SKenneth E. Jansenc
11559599516SKenneth E. Jansen
11659599516SKenneth E. Jansencccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
11759599516SKenneth E. Jansen!MR CHANGE
11859599516SKenneth E. Jansen
11959599516SKenneth E. Jansen!      nfiles = 2
12059599516SKenneth E. Jansen!      nfields = 31
12159599516SKenneth E. Jansen!      numparts = 8
12259599516SKenneth E. Jansen!      nppp = numparts/numpe
12359599516SKenneth E. Jansen!      nppf = numparts/nfiles
12459599516SKenneth E. Jansen
12559599516SKenneth E. Jansen      nfiles = nsynciofiles
12659599516SKenneth E. Jansen!      nfields = nsynciofieldsreadgeombc
12759599516SKenneth E. Jansen      numparts = numpe !This is the common settings. Beware if you try to compute several parts per process
12859599516SKenneth E. Jansen
12959599516SKenneth E. Jansen!      nppp = numparts/numpe
13059599516SKenneth E. Jansen!      nppf = numparts/nfiles
13159599516SKenneth E. Jansen!MR CHANGE END
13259599516SKenneth E. Jansen
13359599516SKenneth E. Jansen      itwo=2
13459599516SKenneth E. Jansen      ione=1
13559599516SKenneth E. Jansen      ieleven=11
13659599516SKenneth E. Jansen      itmp = int(log10(float(myrank+1)))+1
13759599516SKenneth E. Jansen
138d5d2f64dSCameron Smith      call phio_openfile_read(c_char_'geombc-dat.' // char(0), nfiles, fhandle);
13959599516SKenneth E. Jansen
140d5d2f64dSCameron Smith      call phio_readheader(fhandle,c_char_'number of nodes' // char(0),
141e5afe575SCameron Smith     & c_loc(numnp),ione, dataInt, iotype)
14259599516SKenneth E. Jansen
143d5d2f64dSCameron Smith      call phio_readheader(fhandle,c_char_'number of modes' // char(0),
144e5afe575SCameron Smith     & c_loc(nshg),ione, dataInt, iotype)
14559599516SKenneth E. Jansen
146d5d2f64dSCameron Smith      call phio_readheader(fhandle,
147e5afe575SCameron Smith     &  c_char_'number of interior elements' // char(0),
148e5afe575SCameron Smith     &  c_loc(numel),ione, dataInt, iotype)
14959599516SKenneth E. Jansen
150d5d2f64dSCameron Smith      call phio_readheader(fhandle,
151e5afe575SCameron Smith     &  c_char_'number of boundary elements' // char(0),
152e5afe575SCameron Smith     &  c_loc(numelb),ione, dataInt, iotype)
15359599516SKenneth E. Jansen
154d5d2f64dSCameron Smith      call phio_readheader(fhandle,
155e5afe575SCameron Smith     &  c_char_'maximum number of element nodes' // char(0),
156e5afe575SCameron Smith     &  c_loc(nen),ione, dataInt, iotype)
15759599516SKenneth E. Jansen
158d5d2f64dSCameron Smith      call phio_readheader(fhandle,
159e5afe575SCameron Smith     &  c_char_'number of interior tpblocks' // char(0),
160e5afe575SCameron Smith     &  c_loc(nelblk),ione, dataInt, iotype)
16159599516SKenneth E. Jansen
162d5d2f64dSCameron Smith      call phio_readheader(fhandle,
163e5afe575SCameron Smith     & c_char_'number of boundary tpblocks' // char(0),
164e5afe575SCameron Smith     & c_loc(nelblb),ione, dataInt, iotype)
16559599516SKenneth E. Jansen
166d5d2f64dSCameron Smith      call phio_readheader(fhandle,
167e5afe575SCameron Smith     & c_char_'number of nodes with Dirichlet BCs' // char(0),
168e5afe575SCameron Smith     & c_loc(numpbc),ione, dataInt, iotype)
16959599516SKenneth E. Jansen
170d5d2f64dSCameron Smith      call phio_readheader(fhandle,
171e5afe575SCameron Smith     & c_char_'number of shape functions' // char(0),
172e5afe575SCameron Smith     & c_loc(ntopsh),ione, dataInt, iotype)
17359599516SKenneth E. Jansen
17459599516SKenneth E. Jansenc      call closefile( igeom, "read" )
17559599516SKenneth E. Jansenc      call finalizephmpiio( igeom )
17659599516SKenneth E. Jansen
17759599516SKenneth E. Jansen!       if(myrank==0) then
17859599516SKenneth E. Jansen!          print *, "Reading Finished, ********* : ", trim(fname2)
17959599516SKenneth E. Jansen!       endif
18059599516SKenneth E. Jansen
18159599516SKenneth E. Jansen
18259599516SKenneth E. Jansenccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
18359599516SKenneth E. Jansen
18459599516SKenneth E. Jansenc      ieleven=11
18559599516SKenneth E. Jansenc      ione=1
18659599516SKenneth E. Jansenc      fname1='number of nodes?'
18759599516SKenneth E. Jansenc      call readheader(igeomBAK,fname1,numnp,ione,'integer', iotype)
18859599516SKenneth E. Jansenc      fname1='number of modes?'
18959599516SKenneth E. Jansenc      call readheader(igeomBAK,fname1,nshg,ione,'integer', iotype)
19059599516SKenneth E. Jansencc      fname1='number of global modes?'
19159599516SKenneth E. Jansencc      call readheader(igeomBAK,fname1,nshgt,ione,'integer', iotype)
19259599516SKenneth E. Jansenc      fname1='number of interior elements?'
19359599516SKenneth E. Jansenc      call readheader(igeomBAK,fname1,numel,ione,'integer', iotype)
19459599516SKenneth E. Jansenc      fname1='number of boundary elements?'
19559599516SKenneth E. Jansenc      call readheader(igeomBAK,fname1,numelb,ione,'integer', iotype)
19659599516SKenneth E. Jansenc      fname1='maximum number of element nodes?'
19759599516SKenneth E. Jansenc      call readheader(igeomBAK,fname1,nen,ione,'integer', iotype)
19859599516SKenneth E. Jansenc      fname1='number of interior tpblocks?'
19959599516SKenneth E. Jansenc      call readheader(igeomBAK,fname1,nelblk,ione,'integer', iotype)
20059599516SKenneth E. Jansenc      fname1='number of boundary tpblocks?'
20159599516SKenneth E. Jansenc      call readheader(igeomBAK,fname1,nelblb,ione,'integer', iotype)
20259599516SKenneth E. Jansenc      fname1='number of nodes with Dirichlet BCs?'
20359599516SKenneth E. Jansenc      call readheader(igeomBAK,fname1,numpbc,ione,'integer', iotype)
20459599516SKenneth E. Jansenc      fname1='number of shape functions?'
20559599516SKenneth E. Jansenc      call readheader(igeomBAK,fname1,ntopsh,ione,'integer', iotype)
20659599516SKenneth E. Jansen
20759599516SKenneth E. Jansenc
20859599516SKenneth E. Jansenc.... calculate the maximum number of boundary element nodes
20959599516SKenneth E. Jansenc
21059599516SKenneth E. Jansen      nenb = 0
21159599516SKenneth E. Jansen      do i = 1, melCat
21259599516SKenneth E. Jansen         if (nen .eq. nenCat(i,nsd)) nenb = max(nenCat(i,nsd-1), nenb)
21359599516SKenneth E. Jansen      enddo
21459599516SKenneth E. Jansenc
21559599516SKenneth E. Jansen      if (myrank == master) then
21659599516SKenneth E. Jansen         if (nenb .eq. 0) call error ('input   ','nen     ',nen)
21759599516SKenneth E. Jansen      endif
21859599516SKenneth E. Jansenc
21959599516SKenneth E. Jansenc.... setup some useful constants
22059599516SKenneth E. Jansenc
22159599516SKenneth E. Jansen      I3nsd  = nsd / 3          ! nsd=3 integer flag
22259599516SKenneth E. Jansen      E3nsd  = float(I3nsd)     ! nsd=3 real    flag
22359599516SKenneth E. Jansenc
22459599516SKenneth E. Jansen      if(matflg(1,1).lt.0) then
22559599516SKenneth E. Jansen         nflow = nsd + 1
22659599516SKenneth E. Jansen      else
22759599516SKenneth E. Jansen         nflow = nsd + 2
22859599516SKenneth E. Jansen      endif
22959599516SKenneth E. Jansen      ndof   = nsd + 2
23059599516SKenneth E. Jansen      nsclr=impl(1)/100
23159599516SKenneth E. Jansen      ndof=ndof+nsclr           ! number of sclr transport equations to solve
23259599516SKenneth E. Jansen
23359599516SKenneth E. Jansen      ndofBC = ndof + I3nsd     ! dimension of BC array
23459599516SKenneth E. Jansen      ndiBCB = 2                ! dimension of iBCB array
23559599516SKenneth E. Jansen      ndBCB  = ndof + 1         ! dimension of BCB array
23659599516SKenneth E. Jansenc
23759599516SKenneth E. Jansen      nsymdf = (ndof*(ndof + 1)) / 2 ! symm. d.o.f.'s
23859599516SKenneth E. Jansenc
23959599516SKenneth E. Jansenc.... ----------------------> Communication tasks <--------------------
24059599516SKenneth E. Jansenc
24159599516SKenneth E. Jansen
24259599516SKenneth E. Jansencc still read in new
24359599516SKenneth E. Jansen
24459599516SKenneth E. Jansen      if(numpe > 1) then
24559599516SKenneth E. Jansen
24659599516SKenneth E. Jansencc MR CHANGE
247d5d2f64dSCameron Smith         call phio_readheader(fhandle,
248e5afe575SCameron Smith     &    c_char_'size of ilwork array' // char(0),
249e5afe575SCameron Smith     &    c_loc(nlwork),ione, dataInt, iotype)
25059599516SKenneth E. Jansen
251d5d2f64dSCameron Smith         call phio_readheader(fhandle,
252e5afe575SCameron Smith     &    c_char_'ilwork' //char(0),
253e5afe575SCameron Smith     &    c_loc(nlwork),ione, dataInt, iotype)
25459599516SKenneth E. Jansen
25559599516SKenneth E. Jansen         allocate( point2ilwork(nlwork) )
25659599516SKenneth E. Jansen         allocate( ilworkread(nlwork) )
257d5d2f64dSCameron Smith         call phio_readdatablock(fhandle, c_char_'ilwork' // char(0),
258bc62cfd4SCameron Smith     &      c_loc(ilworkread), nlwork, dataInt , iotype)
25959599516SKenneth E. Jansen
26059599516SKenneth E. Jansenc      call closefile( igeom, "read" )
26159599516SKenneth E. Jansenc      call finalizephmpiio( igeom )
26259599516SKenneth E. Jansen
26359599516SKenneth E. Jansenc         fname1='size of ilwork array?'
26459599516SKenneth E. Jansenc         call readheader(igeomBAK,fname1,nlwork,ione,'integer', iotype)
26559599516SKenneth E. Jansen
26659599516SKenneth E. Jansenc         ione=1
26759599516SKenneth E. Jansenc         fname1='ilwork?'
26859599516SKenneth E. Jansenc         call readheader(igeomBAK,fname1,nlwork,ione,'integer', iotype)
26959599516SKenneth E. Jansen
27059599516SKenneth E. Jansenc         allocate( point2ilwork(nlwork) )
27159599516SKenneth E. Jansenc         allocate( ilworkread(nlwork) )
27259599516SKenneth E. Jansenc         call readdatablock(igeomBAK,fname1,ilworkread,
27359599516SKenneth E. Jansenc     &                      nlwork,'integer', iotype)
27459599516SKenneth E. Jansencc MR CHANGE
27559599516SKenneth E. Jansen
27659599516SKenneth E. Jansen
27759599516SKenneth E. Jansen         point2ilwork = ilworkread
27859599516SKenneth E. Jansen         call ctypes (point2ilwork)
27959599516SKenneth E. Jansen      else
28059599516SKenneth E. Jansen           nlwork=1
28159599516SKenneth E. Jansen           allocate( point2ilwork(1))
28259599516SKenneth E. Jansen           nshg0 = nshg
28359599516SKenneth E. Jansen      endif
28459599516SKenneth E. Jansen
28559599516SKenneth E. Jansencccccccccccccccccccccccccccccccccccccccccc
28659599516SKenneth E. Jansen
28759599516SKenneth E. Jansen      itwo=2
28859599516SKenneth E. Jansen
28959599516SKenneth E. Jansenc      print *, "fname2 is", fname2
29059599516SKenneth E. Jansen
29159599516SKenneth E. Jansencc MR CHANGE
29259599516SKenneth E. Jansenc      call initphmpiio(nfields,nppf,nfiles,igeom,'read')
29359599516SKenneth E. Jansenc      call openfile( fnamer, 'read', igeom )
29459599516SKenneth E. JansenCC MR CHANGE
29559599516SKenneth E. Jansen
296d5d2f64dSCameron Smith      call phio_readheader(fhandle,
297e5afe575SCameron Smith     & c_char_'co-ordinates' // char(0),
298e5afe575SCameron Smith     & c_loc(intfromfile),itwo, dataDbl, iotype)
29959599516SKenneth E. Jansen      numnp=intfromfile(1)
30059599516SKenneth E. Jansenc      print *, "read out @@@@@@ is ", numnp
30159599516SKenneth E. Jansen      allocate( point2x(numnp,nsd) )
30259599516SKenneth E. Jansen      allocate( xread(numnp,nsd) )
30359599516SKenneth E. Jansen      ixsiz=numnp*nsd
304d5d2f64dSCameron Smith      call phio_readdatablock(fhandle,
305bc62cfd4SCameron Smith     & c_char_'co-ordinates' // char(0),
306bc62cfd4SCameron Smith     & c_loc(xread),ixsiz, dataDbl, iotype)
30759599516SKenneth E. Jansen      point2x = xread
30859599516SKenneth E. Jansen
30959599516SKenneth E. Jansen!      call closefile( igeom, "read" )
31059599516SKenneth E. Jansen!      call finalizephmpiio( igeom )
31159599516SKenneth E. Jansen
31259599516SKenneth E. Jansen!       print *, "Finished finalize"
31359599516SKenneth E. Jansen
31459599516SKenneth E. Jansenc      deallocate (point2x)
31559599516SKenneth E. Jansenc      deallocate (xread)
31659599516SKenneth E. Jansen
31759599516SKenneth E. Jansencccccccccccccccccccccccccccccccccccccccccc
31859599516SKenneth E. Jansen
31959599516SKenneth E. Jansenc
32059599516SKenneth E. Jansenc.... read the node coordinates
32159599516SKenneth E. Jansenc
32259599516SKenneth E. Jansen
32359599516SKenneth E. Jansenc      itwo=2
32459599516SKenneth E. Jansenc      fname1='co-ordinates?'
32559599516SKenneth E. Jansenc      call readheader(igeomBAK,fname1,intfromfile,itwo, 'double', iotype)
32659599516SKenneth E. Jansenc      numnp=intfromfile(1)
32759599516SKenneth E. Jansencc      nsd=intfromfile(2)
32859599516SKenneth E. Jansenc      allocate( point2x(numnp,nsd) )
32959599516SKenneth E. Jansenc      allocate( xread(numnp,nsd) )
33059599516SKenneth E. Jansenc      ixsiz=numnp*nsd
33159599516SKenneth E. Jansenc      call readdatablock(igeomBAK,fname1,xread,ixsiz, 'double',iotype)
33259599516SKenneth E. Jansenc      point2x = xread
33359599516SKenneth E. Jansen
33459599516SKenneth E. Jansenc
33559599516SKenneth E. Jansenc.... read in and block out the connectivity
33659599516SKenneth E. Jansenc
33759599516SKenneth E. Jansen
33859599516SKenneth E. Jansen! !MR CHANGE
33959599516SKenneth E. Jansen!     This is not necessary but this avoids to have the geombc files opend two times.
34059599516SKenneth E. Jansen!     A better way consists in pasisng the file handler to genblk or make it global or use igeomBAK instead of igeom
34159599516SKenneth E. Jansen!      call closefile( igeom, "read" )
34259599516SKenneth E. Jansen!      call finalizephmpiio( igeom )
34359599516SKenneth E. Jansen! !MR CHANGE END
34459599516SKenneth E. Jansen
34559599516SKenneth E. Jansen      call genblk (IBKSIZ)
34659599516SKenneth E. Jansen
34759599516SKenneth E. Jansen! !MR CHANGE
34859599516SKenneth E. Jansen!      call initphmpiio( nfields, nppf, nfiles, igeom, 'read')
34959599516SKenneth E. Jansen!      call openfile( fnamer, 'read', igeom )
35059599516SKenneth E. Jansen! !MR CHANGE END
35159599516SKenneth E. Jansen
35259599516SKenneth E. Jansenc
35359599516SKenneth E. Jansenc.... read the boundary condition mapping array
35459599516SKenneth E. Jansenc
35559599516SKenneth E. Jansen
35659599516SKenneth E. Jansencc MR CHANGE
35759599516SKenneth E. Jansen!      call initphmpiio(nfields,nppf,nfiles,igeom, 'read')
35859599516SKenneth E. Jansen!      call openfile( fnamer, 'read', igeom )
35959599516SKenneth E. Jansencc MR CHANGE
36059599516SKenneth E. Jansen
36159599516SKenneth E. Jansen      ione=1
362d5d2f64dSCameron Smith      call phio_readheader(fhandle,
363e5afe575SCameron Smith     & c_char_'bc mapping array' // char(0),
364e5afe575SCameron Smith     & c_loc(nshg),ione, dataInt, iotype)
36559599516SKenneth E. Jansen
36659599516SKenneth E. Jansenc      fname1='bc mapping array?'
36759599516SKenneth E. Jansenc      call readheader(igeomBAK,fname1,nshg,
36859599516SKenneth E. Jansenc     &     ione,'integer', iotype)
36959599516SKenneth E. Jansen
37059599516SKenneth E. Jansen      allocate( nBC(nshg) )
37159599516SKenneth E. Jansen
37259599516SKenneth E. Jansen      allocate( nBCread(nshg) )
37359599516SKenneth E. Jansen
37459599516SKenneth E. Jansenc      call readdatablock(igeomBAK,fname1,nBCread,nshg,'integer',iotype)
375d5d2f64dSCameron Smith      call phio_readdatablock(fhandle,
376bc62cfd4SCameron Smith     & c_char_'bc mapping array' // char(0),
377bc62cfd4SCameron Smith     & c_loc(nBCread), nshg, dataInt, iotype)
37859599516SKenneth E. Jansen
37959599516SKenneth E. Jansen      nBC=nBCread
38059599516SKenneth E. Jansen
38159599516SKenneth E. Jansenc
38259599516SKenneth E. Jansenc.... read the temporary iBC array
38359599516SKenneth E. Jansenc
38459599516SKenneth E. Jansen      ione=1
385d5d2f64dSCameron Smith      call phio_readheader(fhandle,
386e5afe575SCameron Smith     & c_char_'bc codes array' // char(0),
387e5afe575SCameron Smith     & c_loc(numpbc),ione, dataInt, iotype)
38859599516SKenneth E. Jansen
38959599516SKenneth E. Jansenc      ione = 1
39059599516SKenneth E. Jansenc      fname1='bc codes array?'
39159599516SKenneth E. Jansenc      call readheader(igeomBAK,fname1,numpbc,
39259599516SKenneth E. Jansenc     &     ione, 'integer', iotype)
39359599516SKenneth E. Jansen
39459599516SKenneth E. Jansen!MR CHANGE
39559599516SKenneth E. Jansen!       if ( numpbc > 0 ) then
39659599516SKenneth E. Jansen!          allocate( iBCtmp(numpbc) )
39759599516SKenneth E. Jansen!          allocate( iBCtmpread(numpbc) )
39859599516SKenneth E. Jansen! c         call readdatablock(igeomBAK,fname1,iBCtmpread,numpbc,
39959599516SKenneth E. Jansen! c     &                      'integer',iotype)
40059599516SKenneth E. Jansen!         call readdatablock(igeom,fname2,iBCtmpread,numpbc,
40159599516SKenneth E. Jansen!      &                      'integer',iotype)
40259599516SKenneth E. Jansen!          iBCtmp=iBCtmpread
40359599516SKenneth E. Jansen!       else  ! sometimes a partition has no BC's
40459599516SKenneth E. Jansen!          allocate( iBCtmp(1) )
40559599516SKenneth E. Jansen!          iBCtmp=0
40659599516SKenneth E. Jansen!       endif
40759599516SKenneth E. Jansen
40859599516SKenneth E. Jansen      if ( numpbc > 0 ) then
40959599516SKenneth E. Jansen        allocate( iBCtmp(numpbc) )
41059599516SKenneth E. Jansen        allocate( iBCtmpread(numpbc) )
41159599516SKenneth E. Jansen      else
41259599516SKenneth E. Jansen        allocate( iBCtmp(1) )
41359599516SKenneth E. Jansen        allocate( iBCtmpread(1) )
41459599516SKenneth E. Jansen      endif
41559599516SKenneth E. Jansenc         call readdatablock(igeomBAK,fname1,iBCtmpread,numpbc,
41659599516SKenneth E. Jansenc     &                      'integer',iotype)
417d5d2f64dSCameron Smith      call phio_readdatablock(fhandle,
418bc62cfd4SCameron Smith     & c_char_'bc codes array' // char(0),
419bc62cfd4SCameron Smith     & c_loc(iBCtmpread), numpbc, dataInt, iotype)
42059599516SKenneth E. Jansen
42159599516SKenneth E. Jansen      if ( numpbc > 0 ) then
42259599516SKenneth E. Jansen         iBCtmp=iBCtmpread
42359599516SKenneth E. Jansen      else  ! sometimes a partition has no BC's
42459599516SKenneth E. Jansen         deallocate( iBCtmpread)
42559599516SKenneth E. Jansen         iBCtmp=0
42659599516SKenneth E. Jansen      endif
42759599516SKenneth E. Jansen!MR CHANGE END
42859599516SKenneth E. Jansen
42959599516SKenneth E. Jansenc
43059599516SKenneth E. Jansenc.... read boundary condition data
43159599516SKenneth E. Jansenc
43259599516SKenneth E. Jansen
43359599516SKenneth E. Jansen      ione=1
43459599516SKenneth E. Jansen
43559599516SKenneth E. Jansenc      ione=1
43659599516SKenneth E. Jansenc      fname1='boundary condition array?'
43759599516SKenneth E. Jansenc      call readheader(igeomBAK,fname1,intfromfile,
43859599516SKenneth E. Jansenc     &     ione, 'double', iotype)
439d5d2f64dSCameron Smith      call phio_readheader(fhandle,
440e5afe575SCameron Smith     & c_char_'boundary condition array' // char(0),
441e5afe575SCameron Smith     & c_loc(intfromfile),ione, dataDbl, iotype)
44259599516SKenneth E. Jansenc here intfromfile(1) contains (ndof+7)*numpbc
44359599516SKenneth E. Jansen!MR CHANGE
44459599516SKenneth E. Jansen!       if ( numpbc > 0 ) then
44559599516SKenneth E. Jansen!          if(intfromfile(1).ne.(ndof+7)*numpbc) then
44659599516SKenneth E. Jansen!            warning='WARNING more data in BCinp than needed: keeping 1st'
44759599516SKenneth E. Jansen!            write(*,*) warning, ndof+7
44859599516SKenneth E. Jansen!          endif
44959599516SKenneth E. Jansen!          allocate( BCinp(numpbc,ndof+7) )
45059599516SKenneth E. Jansen!          nsecondrank=intfromfile(1)/numpbc
45159599516SKenneth E. Jansen!          allocate( BCinpread(numpbc,nsecondrank) )
45259599516SKenneth E. Jansen!          iBCinpsiz=intfromfile(1)
45359599516SKenneth E. Jansen! c         call readdatablock(igeomBAK,fname1,BCinpread,iBCinpsiz,
45459599516SKenneth E. Jansen! c     &                      'double',iotype)
45559599516SKenneth E. Jansen!          call readdatablock(igeom,fname2,BCinpread,iBCinpsiz,
45659599516SKenneth E. Jansen!      &                      'double',iotype)
45759599516SKenneth E. Jansen!          BCinp(:,1:(ndof+7))=BCinpread(:,1:(ndof+7))
45859599516SKenneth E. Jansen!       else  ! sometimes a partition has no BC's
45959599516SKenneth E. Jansen!          allocate( BCinp(1,ndof+7) )
46059599516SKenneth E. Jansen!          BCinp=0
46159599516SKenneth E. Jansen!       endif
46259599516SKenneth E. Jansen
46359599516SKenneth E. Jansen      if ( numpbc > 0 ) then
46459599516SKenneth E. Jansen!         if(intfromfile(1).ne.(ndof+7)*numpbc) then
46559599516SKenneth E. Jansen!           warning='WARNING more data in BCinp than needed: keeping 1st'
46659599516SKenneth E. Jansen!           write(*,*) warning, ndof+7
46759599516SKenneth E. Jansen!         endif
46859599516SKenneth E. Jansen         allocate( BCinp(numpbc,ndof+7) )
46959599516SKenneth E. Jansen         nsecondrank=intfromfile(1)/numpbc
47059599516SKenneth E. Jansen         allocate( BCinpread(numpbc,nsecondrank) )
47159599516SKenneth E. Jansen         iBCinpsiz=intfromfile(1)
47259599516SKenneth E. Jansen      else
47359599516SKenneth E. Jansen         allocate( BCinp(1,ndof+7) )
47459599516SKenneth E. Jansen         allocate( BCinpread(0,0) ) !dummy
47559599516SKenneth E. Jansen         iBCinpsiz=intfromfile(1)
47659599516SKenneth E. Jansen      endif
47759599516SKenneth E. Jansenc         call readdatablock(igeomBAK,fname1,BCinpread,iBCinpsiz,
47859599516SKenneth E. Jansenc     &                      'double',iotype)
47959599516SKenneth E. Jansen
480d5d2f64dSCameron Smith      call phio_readdatablock(fhandle,
481bc62cfd4SCameron Smith     & c_char_'boundary condition array' // char(0),
482bc62cfd4SCameron Smith     & c_loc(BCinpread), iBCinpsiz, dataDbl, iotype)
48359599516SKenneth E. Jansen
48459599516SKenneth E. Jansen      if ( numpbc > 0 ) then
48559599516SKenneth E. Jansen         BCinp(:,1:(ndof+7))=BCinpread(:,1:(ndof+7))
48659599516SKenneth E. Jansen      else  ! sometimes a partition has no BC's
48759599516SKenneth E. Jansen         deallocate(BCinpread)
48859599516SKenneth E. Jansen         BCinp=0
48959599516SKenneth E. Jansen      endif
49059599516SKenneth E. Jansen!MR CHANGE END
49159599516SKenneth E. Jansen
49259599516SKenneth E. Jansenc
49359599516SKenneth E. Jansenc.... read periodic boundary conditions
49459599516SKenneth E. Jansenc
49559599516SKenneth E. Jansen
49659599516SKenneth E. Jansen      ione=1
49759599516SKenneth E. Jansenc      fname1='periodic masters array?'
49859599516SKenneth E. Jansenc      call readheader(igeomBAK,fname1,nshg,
49959599516SKenneth E. Jansenc     &     ione, 'integer', iotype)
500d5d2f64dSCameron Smith      call phio_readheader(fhandle,
501e5afe575SCameron Smith     & c_char_'periodic masters array' // char(0),
502e5afe575SCameron Smith     & c_loc(nshg), ione, dataInt, iotype)
50359599516SKenneth E. Jansen      allocate( point2iper(nshg) )
50459599516SKenneth E. Jansen      allocate( iperread(nshg) )
50559599516SKenneth E. Jansenc      call readdatablock(igeomBAK,fname1,iperread,nshg,
50659599516SKenneth E. Jansenc     &                      'integer',iotype)
507d5d2f64dSCameron Smith      call phio_readdatablock(fhandle,
508bc62cfd4SCameron Smith     & c_char_'periodic masters array' // char(0),
509bc62cfd4SCameron Smith     & c_loc(iperread), nshg, dataInt, iotype)
51059599516SKenneth E. Jansen      point2iper=iperread
51159599516SKenneth E. Jansen
51259599516SKenneth E. Jansen
51359599516SKenneth E. Jansen! !MR CHANGE
51459599516SKenneth E. Jansen!      call closefile( igeom, "read" )
51559599516SKenneth E. Jansen!      call finalizephmpiio( igeom )
51659599516SKenneth E. Jansen! !MR CHANGE END
51759599516SKenneth E. Jansen
51859599516SKenneth E. Jansenc
51959599516SKenneth E. Jansenc.... generate the boundary element blocks
52059599516SKenneth E. Jansenc
52159599516SKenneth E. Jansen      call genbkb (ibksiz)
52259599516SKenneth E. Jansen
52359599516SKenneth E. Jansen
52459599516SKenneth E. Jansen! !MR CHANGE
52559599516SKenneth E. Jansen!       write(*,*) 'HELLO 12 from ', myrank
52659599516SKenneth E. Jansen! !MR CHANGE END
52759599516SKenneth E. Jansen
52859599516SKenneth E. Jansenc
52959599516SKenneth E. Jansenc  Read in the nsons and ifath arrays if needed
53059599516SKenneth E. Jansenc
53159599516SKenneth E. Jansenc  There is a fundamental shift in the meaning of ifath based on whether
53259599516SKenneth E. Jansenc  there exist homogenous directions in the flow.
53359599516SKenneth E. Jansenc
53459599516SKenneth E. Jansenc  HOMOGENOUS DIRECTIONS EXIST:  Here nfath is the number of inhomogenous
53559599516SKenneth E. Jansenc  points in the TOTAL mesh.  That is to say that each partition keeps a
53659599516SKenneth E. Jansenc  link to  ALL inhomogenous points.  This link is furthermore not to the
53759599516SKenneth E. Jansenc  sms numbering but to the original structured grid numbering.  These
53859599516SKenneth E. Jansenc  inhomogenous points are thought of as fathers, with their sons being all
53959599516SKenneth E. Jansenc  the points in the homogenous directions that have this father's
54059599516SKenneth E. Jansenc  inhomogeneity.  The array ifath takes as an arguement the sms numbering
54159599516SKenneth E. Jansenc  and returns as a result the father.
54259599516SKenneth E. Jansenc
54359599516SKenneth E. Jansenc  In this case nsons is the number of sons that each father has and ifath
54459599516SKenneth E. Jansenc  is an array which tells the
54559599516SKenneth E. Jansenc
54659599516SKenneth E. Jansenc  NO HOMOGENOUS DIRECTIONS.  In this case the mesh would grow to rapidly
54759599516SKenneth E. Jansenc  if we followed the above strategy since every partition would index its
54859599516SKenneth E. Jansenc  points to the ENTIRE mesh.  Furthermore, there would never be a need
54959599516SKenneth E. Jansenc  to average to a node off processor since there is no spatial averaging.
55059599516SKenneth E. Jansenc  Therefore, to properly account for this case we must recognize it and
55159599516SKenneth E. Jansenc  inerrupt certain actions (i.e. assembly of the average across partitions).
55259599516SKenneth E. Jansenc  This case is easily identified by noting that maxval(nsons) =1 (i.e. no
55359599516SKenneth E. Jansenc  father has any sons).  Reiterating to be clear, in this case ifath does
55459599516SKenneth E. Jansenc  not point to a global numbering but instead just points to itself.
55559599516SKenneth E. Jansenc
55659599516SKenneth E. Jansen      nfath=1  ! some architectures choke on a zero or undeclared
55759599516SKenneth E. Jansen                 ! dimension variable.  This sets it to a safe, small value.
55859599516SKenneth E. Jansen      if(((iLES .lt. 20) .and. (iLES.gt.0))
55959599516SKenneth E. Jansen     &                   .or. (itwmod.gt.0)  ) then ! don't forget same
56059599516SKenneth E. Jansen                                                    ! conditional in proces.f
56159599516SKenneth E. Jansen
56259599516SKenneth E. Jansenc           read (igeomBAK) nfath  ! nfath already read in input.f,
56359599516SKenneth E. Jansen                                     ! needed for alloc
56459599516SKenneth E. Jansen         ione=1
56559599516SKenneth E. Jansenc         call creadlist(igeomBAK,ione,nfath)
56659599516SKenneth E. Jansenc         fname1='keyword sonfath?'
56759599516SKenneth E. Jansen         if(nohomog.gt.0) then
56859599516SKenneth E. Jansen
56959599516SKenneth E. Jansen
57059599516SKenneth E. Jansen!             fname1='number of father-nodes?'
57159599516SKenneth E. Jansen!             call readheader(igeomBAK,fname1,nfath,ione,'integer', iotype)
57259599516SKenneth E. Jansen
573d5d2f64dSCameron Smith            call phio_readheader(fhandle,
574e5afe575SCameron Smith     &       c_char_'number of father-nodes' // char(0),
575e5afe575SCameron Smith     &       c_loc(nfath), ione, dataInt, iotype)
57659599516SKenneth E. Jansen
57759599516SKenneth E. Jansenc
57859599516SKenneth E. Jansenc     fname1='keyword nsons?'
57959599516SKenneth E. Jansen!             fname1='number of son-nodes for each father?'
58059599516SKenneth E. Jansen!             call readheader(igeomBAK,fname1,nfath,ione,'integer', iotype)
58159599516SKenneth E. Jansen
582d5d2f64dSCameron Smith            call phio_readheader(fhandle,
583e5afe575SCameron Smith     &       c_char_'number of son-nodes for each father' // char(0),
584e5afe575SCameron Smith     &       c_loc(nfath), ione, dataInt, iotype)
58559599516SKenneth E. Jansen
58659599516SKenneth E. Jansen            allocate (point2nsons(nfath))
58759599516SKenneth E. Jansen
58859599516SKenneth E. Jansen!             call readdatablock(igeomBAK,fname1,point2nsons,nfath,
58959599516SKenneth E. Jansen!      &                      'integer',iotype)
590d5d2f64dSCameron Smith            call phio_readdatablock(fhandle,
591bc62cfd4SCameron Smith     &       c_char_'number of son-nodes for each father' // char(0),
592bc62cfd4SCameron Smith     &       c_loc(point2nsons),nfath, dataInt, iotype)
59359599516SKenneth E. Jansen
59459599516SKenneth E. Jansenc
59559599516SKenneth E. Jansen!             fname1='keyword ifath?'
59659599516SKenneth E. Jansen!             call readheader(igeomBAK,fname1,nshg,ione,'integer', iotype)
59759599516SKenneth E. Jansen
598d5d2f64dSCameron Smith            call phio_readheader(fhandle,
599e5afe575SCameron Smith     &       c_char_'keyword ifath' // char(0),
600e5afe575SCameron Smith     &       c_loc(nshg), ione, dataInt, iotype);
60159599516SKenneth E. Jansen
60259599516SKenneth E. Jansen            allocate (point2ifath(nshg))
60359599516SKenneth E. Jansen
60459599516SKenneth E. Jansen!             call readdatablock(igeomBAK,fname1,point2ifath,nshg,
60559599516SKenneth E. Jansen!      &                      'integer',iotype)
606d5d2f64dSCameron Smith            call phio_readdatablock(fhandle,
607bc62cfd4SCameron Smith     &       c_char_'keyword ifath' // char(0),
608bc62cfd4SCameron Smith     &       c_loc(point2ifath), nshg, dataInt, iotype)
60959599516SKenneth E. Jansen
61059599516SKenneth E. Jansenc
61159599516SKenneth E. Jansen            nsonmax=maxval(point2nsons)
61259599516SKenneth E. Jansenc
61359599516SKenneth E. Jansen         else  ! this is the case where there is no homogeneity
61459599516SKenneth E. Jansen               ! therefore ever node is a father (too itself).  sonfath
61559599516SKenneth E. Jansen               ! (a routine in NSpre) will set this up but this gives
61659599516SKenneth E. Jansen               ! you an option to avoid that.
61759599516SKenneth E. Jansen            nfath=nshg
61859599516SKenneth E. Jansen            allocate (point2nsons(nfath))
61959599516SKenneth E. Jansen            point2nsons=1
62059599516SKenneth E. Jansen            allocate (point2ifath(nshg))
62159599516SKenneth E. Jansen            do i=1,nshg
62259599516SKenneth E. Jansen               point2ifath(i)=i
62359599516SKenneth E. Jansen            enddo
62459599516SKenneth E. Jansen            nsonmax=1
62559599516SKenneth E. Jansenc
62659599516SKenneth E. Jansen         endif
62759599516SKenneth E. Jansen      else
62859599516SKenneth E. Jansen         allocate (point2nsons(1))
62959599516SKenneth E. Jansen         allocate (point2ifath(1))
63059599516SKenneth E. Jansen      endif
63159599516SKenneth E. Jansen
632d5d2f64dSCameron Smith      call phio_closefile_read(fhandle);
63359599516SKenneth E. Jansen
63459599516SKenneth E. Jansen! !MR CHANGE
63559599516SKenneth E. Jansen!       write(*,*) 'HELLO 13 from ', myrank
63659599516SKenneth E. Jansen! !MR CHANGE END
63759599516SKenneth E. Jansen
63859599516SKenneth E. Jansenc
63959599516SKenneth E. Jansenc  renumber the master partition for SPEBC
64059599516SKenneth E. Jansenc
64159599516SKenneth E. Jansenc      if((myrank.eq.master).and.(irscale.ge.0)) then
64259599516SKenneth E. Jansenc         call setSPEBC(numnp, nfath, nsonmax)
64359599516SKenneth E. Jansenc         call renum(point2x,point2ifath,point2nsons)
64459599516SKenneth E. Jansenc      endif
64559599516SKenneth E. Jansenc
64659599516SKenneth E. Jansenc.... Read restart files
64759599516SKenneth E. Jansen
64859599516SKenneth E. Jansenc$$$ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
64959599516SKenneth E. Jansenc
65059599516SKenneth E. Jansenc      nfields = 11
65159599516SKenneth E. Jansenc      numparts = 512
65259599516SKenneth E. Jansenc      nppp = numparts/numpe
65359599516SKenneth E. Jansenc      startpart = myrank * nppp +1
65459599516SKenneth E. Jansenc      int endpart = startpart + nppp - 1
65559599516SKenneth E. Jansenc      nppf = numparts/nfiles
65659599516SKenneth E. Jansencc      fnamer = "/users/nliu/PIG4/4-procs_case/restart-file"
65759599516SKenneth E. Jansencc      fnamer="./4-procs_case/restart-file"
65859599516SKenneth E. Jansen
65936adee64SCameron Smith      fnamer = c_char_"restart-dat."//c_null_char
66036adee64SCameron Smith      call phio_appendStep(fnamer, irstart)
661d5d2f64dSCameron Smith      call phio_openfile_read(fnamer, nfiles, fhandle)
66259599516SKenneth E. Jansen
66359599516SKenneth E. Jansen      ithree=3
66459599516SKenneth E. Jansenc      call creadlist(irstin,ithree,nshg2,ndof2,lstep)
66559599516SKenneth E. Jansen
66659599516SKenneth E. Jansen      itmp = int(log10(float(myrank+1)))+1
66759599516SKenneth E. Jansen
66859599516SKenneth E. Jansenc      print *, "Solution is : ", fname1
66959599516SKenneth E. Jansen
67059599516SKenneth E. Jansen      intfromfile=0
671d5d2f64dSCameron Smith      call phio_readheader(fhandle,
672e5afe575SCameron Smith     & c_char_'solution' // char(0),
673e5afe575SCameron Smith     & c_loc(intfromfile), ithree, dataInt, iotype)
67459599516SKenneth E. Jansenc
67559599516SKenneth E. Jansenc.... read the values of primitive variables into q
67659599516SKenneth E. Jansenc
67759599516SKenneth E. Jansen
67859599516SKenneth E. Jansenc      print *, "intfromfile(1) is ", intfromfile(1)
67959599516SKenneth E. Jansenc      print *, "intfromfile(2) is ", intfromfile(2)
68059599516SKenneth E. Jansenc      print *, "intfromfile(3) is ", intfromfile(3)
68159599516SKenneth E. Jansen
68259599516SKenneth E. Jansen      allocate( qold(nshg,ndof) )
68359599516SKenneth E. Jansen      if(intfromfile(1).ne.0) then
68459599516SKenneth E. Jansen         nshg2=intfromfile(1)
68559599516SKenneth E. Jansen         ndof2=intfromfile(2)
68659599516SKenneth E. Jansen         lstep=intfromfile(3)
68759599516SKenneth E. Jansen         if(ndof2.ne.ndof) then
68859599516SKenneth E. Jansen
68959599516SKenneth E. Jansen         endif
69059599516SKenneth E. Jansenc
69159599516SKenneth E. Jansen        if (nshg2 .ne. nshg)
69259599516SKenneth E. Jansen     &        call error ('restar  ', 'nshg   ', nshg)
69359599516SKenneth E. Jansen         allocate( qread(nshg,ndof2) )
69459599516SKenneth E. Jansen         iqsiz=nshg*ndof2
695d5d2f64dSCameron Smith         call phio_readdatablock(fhandle,
696bc62cfd4SCameron Smith     &    c_char_'solution' // char(0),
697bc62cfd4SCameron Smith     &    c_loc(qread),iqsiz, dataDbl,iotype)
69859599516SKenneth E. Jansen         qold(:,1:ndof)=qread(:,1:ndof)
69959599516SKenneth E. Jansen         deallocate(qread)
70059599516SKenneth E. Jansen      else
70159599516SKenneth E. Jansen         if (myrank.eq.master) then
70259599516SKenneth E. Jansen            if (matflg(1,1).eq.0) then ! compressible
70359599516SKenneth E. Jansen               warning='Solution is set to zero (with p and T to one)'
70459599516SKenneth E. Jansen            else
70559599516SKenneth E. Jansen               warning='Solution is set to zero'
70659599516SKenneth E. Jansen            endif
70759599516SKenneth E. Jansen            write(*,*) warning
70859599516SKenneth E. Jansen         endif
70959599516SKenneth E. Jansen         qold=zero
71059599516SKenneth E. Jansen         if (matflg(1,1).eq.0) then ! compressible
71159599516SKenneth E. Jansen            qold(:,1)=one ! avoid zero pressure
71259599516SKenneth E. Jansen            qold(:,nflow)=one ! avoid zero temperature
71359599516SKenneth E. Jansen         endif
71459599516SKenneth E. Jansen      endif
71559599516SKenneth E. Jansen
71659599516SKenneth E. Jansen
71759599516SKenneth E. Jansen! !MR CHANGE
71859599516SKenneth E. Jansen!       write(*,*) 'HELLO 16-8 from ', myrank
71959599516SKenneth E. Jansen! !MR CHANGE END
72059599516SKenneth E. Jansen
72159599516SKenneth E. Jansenc      itmp=1
72259599516SKenneth E. Jansenc      if (myrank .gt. 0) itmp = int(log10(float(myrank)))+1
72359599516SKenneth E. Jansen      intfromfile=0
724d5d2f64dSCameron Smith      call phio_readheader(fhandle,
725e5afe575SCameron Smith     & c_char_'time derivative of solution' // char(0),
726e5afe575SCameron Smith     & c_loc(intfromfile), ithree, dataInt, iotype)
72759599516SKenneth E. Jansen      allocate( acold(nshg,ndof) )
72859599516SKenneth E. Jansen      if(intfromfile(1).ne.0) then
72959599516SKenneth E. Jansen         nshg2=intfromfile(1)
73059599516SKenneth E. Jansen         ndof2=intfromfile(2)
73159599516SKenneth E. Jansen         lstep=intfromfile(3)
73259599516SKenneth E. Jansen
73359599516SKenneth E. Jansenc      print *, "intfromfile(1) is ", intfromfile(1)
73459599516SKenneth E. Jansenc      print *, "intfromfile(2) is ", intfromfile(2)
73559599516SKenneth E. Jansenc      print *, "intfromfile(3) is ", intfromfile(3)
73659599516SKenneth E. Jansen
73759599516SKenneth E. Jansen         if (nshg2 .ne. nshg)
73859599516SKenneth E. Jansen     &        call error ('restar  ', 'nshg   ', nshg)
73959599516SKenneth E. Jansenc
74059599516SKenneth E. Jansen         allocate( acread(nshg,ndof2) )
74159599516SKenneth E. Jansen         acread=zero
74259599516SKenneth E. Jansen         iacsiz=nshg*ndof2
743d5d2f64dSCameron Smith         call phio_readdatablock(fhandle,
744bc62cfd4SCameron Smith     &    c_char_'time derivative of solution' // char(0),
745bc62cfd4SCameron Smith     &    c_loc(acread), iacsiz, dataDbl,iotype)
74659599516SKenneth E. Jansen         acold(:,1:ndof)=acread(:,1:ndof)
74759599516SKenneth E. Jansen         deallocate(acread)
74859599516SKenneth E. Jansen      else
74959599516SKenneth E. Jansen         if (myrank.eq.master) then
75059599516SKenneth E. Jansen            warning='Time derivative of solution is set to zero (SAFE)'
75159599516SKenneth E. Jansen            write(*,*) warning
75259599516SKenneth E. Jansen         endif
75359599516SKenneth E. Jansen         acold=zero
75459599516SKenneth E. Jansen      endif
75559599516SKenneth E. Jansen
75659599516SKenneth E. Jansencccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
75759599516SKenneth E. Jansencc
75859599516SKenneth E. Jansenc
75959599516SKenneth E. Jansencc
76059599516SKenneth E. Jansencc.... read the header and check it against the run data
76159599516SKenneth E. Jansencc
76259599516SKenneth E. Jansen
76359599516SKenneth E. Jansen
76459599516SKenneth E. Jansenc      ithree=3
76559599516SKenneth E. Jansenccc      call creadlist(irstin,ithree,nshg2,ndof2,lstep)
76659599516SKenneth E. Jansenc      fname1='solution?'
76759599516SKenneth E. Jansenc      intfromfile=0
76859599516SKenneth E. Jansenc      call readheader(irstin,fname1,intfromfile,
76959599516SKenneth E. Jansenc     &     ithree,'integer', iotype)
77059599516SKenneth E. Jansencc
77159599516SKenneth E. Jansencc.... read the values of primitive variables into q
77259599516SKenneth E. Jansencc
77359599516SKenneth E. Jansenc      allocate( qold(nshg,ndof) )
77459599516SKenneth E. Jansenc      if(intfromfile(1).ne.0) then
77559599516SKenneth E. Jansenc         nshg2=intfromfile(1)
77659599516SKenneth E. Jansenc         ndof2=intfromfile(2)
77759599516SKenneth E. Jansenc         lstep=intfromfile(3)
77859599516SKenneth E. Jansenc         if(ndof2.ne.ndof) then
77959599516SKenneth E. Jansenc           warning='WARNING more data in restart than needed: keeping 1st '
78059599516SKenneth E. Jansenc           write(*,*) warning , ndof
78159599516SKenneth E. Jansenc         endif
78259599516SKenneth E. Jansencc
78359599516SKenneth E. Jansenc         if (nshg2 .ne. nshg)
78459599516SKenneth E. Jansenc     &        call error ('restar  ', 'nshg   ', nshg)
78559599516SKenneth E. Jansenc         allocate( qread(nshg,ndof2) )
78659599516SKenneth E. Jansen
78759599516SKenneth E. Jansenc         iqsiz=nshg*ndof2
78859599516SKenneth E. Jansenc         call readdatablock(irstin,fname1,qread,iqsiz,
78959599516SKenneth E. Jansenc     &                         'double',iotype)
79059599516SKenneth E. Jansenc         qold(:,1:ndof)=qread(:,1:ndof)
79159599516SKenneth E. Jansenc         deallocate(qread)
79259599516SKenneth E. Jansenc      else
79359599516SKenneth E. Jansenc         if (myrank.eq.master) then
79459599516SKenneth E. Jansenc            if (matflg(1,1).eq.0) then ! compressible
79559599516SKenneth E. Jansenc               warning='Solution is set to zero (with p and T to one)'
79659599516SKenneth E. Jansenc            else
79759599516SKenneth E. Jansenc               warning='Solution is set to zero'
79859599516SKenneth E. Jansenc            endif
79959599516SKenneth E. Jansenc            write(*,*) warning
80059599516SKenneth E. Jansenc         endif
80159599516SKenneth E. Jansenc         qold=zero
80259599516SKenneth E. Jansenc         if (matflg(1,1).eq.0) then ! compressible
80359599516SKenneth E. Jansenc            qold(:,1)=one ! avoid zero pressure
80459599516SKenneth E. Jansenc            qold(:,nflow)=one ! avoid zero temperature
80559599516SKenneth E. Jansenc         endif
80659599516SKenneth E. Jansenc      endif
80759599516SKenneth E. Jansencc
80859599516SKenneth E. Jansenc      fname1='time derivative of solution?'
80959599516SKenneth E. Jansenc      intfromfile=0
81059599516SKenneth E. Jansenc      call readheader(irstin,fname1,intfromfile,
81159599516SKenneth E. Jansenc     &     ithree,'integer', iotype)
81259599516SKenneth E. Jansenc      allocate( acold(nshg,ndof) )
81359599516SKenneth E. Jansenc      if(intfromfile(1).ne.0) then
81459599516SKenneth E. Jansenc         nshg2=intfromfile(1)
81559599516SKenneth E. Jansenc         ndof2=intfromfile(2)
81659599516SKenneth E. Jansenc         lstep=intfromfile(3)
81759599516SKenneth E. Jansenc
81859599516SKenneth E. Jansenc         if (nshg2 .ne. nshg)
81959599516SKenneth E. Jansenc     &        call error ('restar  ', 'nshg   ', nshg)
82059599516SKenneth E. Jansencc
82159599516SKenneth E. Jansenc         allocate( acread(nshg,ndof2) )
82259599516SKenneth E. Jansenc         acread=zero
82359599516SKenneth E. Jansenc
82459599516SKenneth E. Jansenc         iacsiz=nshg*ndof2
82559599516SKenneth E. Jansenc         call readdatablock(irstin,fname1,acread,iacsiz,
82659599516SKenneth E. Jansenc     &                   'double',iotype)
82759599516SKenneth E. Jansenc         acold(:,1:ndof)=acread(:,1:ndof)
82859599516SKenneth E. Jansenc         deallocate(acread)
82959599516SKenneth E. Jansenc      else
83059599516SKenneth E. Jansenc         if (myrank.eq.master) then
83159599516SKenneth E. Jansenc            warning='Time derivative of solution is set to zero (SAFE)'
83259599516SKenneth E. Jansenc            write(*,*) warning
83359599516SKenneth E. Jansenc         endif
83459599516SKenneth E. Jansenc         acold=zero
83559599516SKenneth E. Jansenc      endif
83659599516SKenneth E. Jansen
83759599516SKenneth E. Jansenc      call creadlist(irstin,ithree,nshg2,ndisp,lstep)
83859599516SKenneth E. Jansen
83959599516SKenneth E. Jansen      if (ideformwall.eq.1) then
84059599516SKenneth E. Jansen!          fname1='displacement?'
84159599516SKenneth E. Jansen!          call readheader(irstin,fname1,intfromfile,
84259599516SKenneth E. Jansen!      &        ithree,'integer', iotype)
84359599516SKenneth E. Jansen
84459599516SKenneth E. Jansen          intfromfile=0
845d5d2f64dSCameron Smith          call phio_readheader(fhandle,
846e5afe575SCameron Smith     &     c_char_'displacement' // char(0),
847e5afe575SCameron Smith     &     c_loc(intfromfile), ithree, dataInt, iotype)
84859599516SKenneth E. Jansen
84959599516SKenneth E. Jansen         nshg2=intfromfile(1)
85059599516SKenneth E. Jansen         ndisp=intfromfile(2)
85159599516SKenneth E. Jansen         lstep=intfromfile(3)
85259599516SKenneth E. Jansen         if(ndisp.ne.nsd) then
85359599516SKenneth E. Jansen            warning='WARNING ndisp not equal nsd'
85459599516SKenneth E. Jansen            write(*,*) warning , ndisp
85559599516SKenneth E. Jansen         endif
85659599516SKenneth E. Jansenc
85759599516SKenneth E. Jansen         if (nshg2 .ne. nshg)
85859599516SKenneth E. Jansen     &        call error ('restar  ', 'nshg   ', nshg)
85959599516SKenneth E. Jansenc
86059599516SKenneth E. Jansenc.... read the values of primitive variables into uold
86159599516SKenneth E. Jansenc
86259599516SKenneth E. Jansen
86359599516SKenneth E. Jansen         allocate( uold(nshg,nsd) )
86459599516SKenneth E. Jansen         allocate( uread(nshg,nsd) )
86559599516SKenneth E. Jansen
86659599516SKenneth E. Jansen         iusiz=nshg*nsd
86759599516SKenneth E. Jansen
86859599516SKenneth E. Jansen!          call readdatablock(irstin,fname1,uread,iusiz,
86959599516SKenneth E. Jansen!      &        'double',iotype)
870d5d2f64dSCameron Smith         call phio_readdatablock(fhandle,
871bc62cfd4SCameron Smith     &    c_char_'displacement' // char(0),
872bc62cfd4SCameron Smith     &    c_loc(uread), iusiz, dataDbl, iotype)
87359599516SKenneth E. Jansen
87459599516SKenneth E. Jansen         uold(:,1:nsd)=uread(:,1:nsd)
87559599516SKenneth E. Jansen       else
87659599516SKenneth E. Jansen         allocate( uold(nshg,nsd) )
87759599516SKenneth E. Jansen         uold(:,1:nsd) = zero
87859599516SKenneth E. Jansen       endif
87959599516SKenneth E. Jansen
88059599516SKenneth E. Jansenc
88159599516SKenneth E. Jansenc.... close c-binary files
88259599516SKenneth E. Jansenc
88359599516SKenneth E. Jansen!MR CHANGE
88459599516SKenneth E. Jansen!      call closefile( irstin, "read" )
88559599516SKenneth E. Jansen
886d5d2f64dSCameron Smith      call phio_closefile_read(fhandle)
88759599516SKenneth E. Jansen
88859599516SKenneth E. Jansen!MR CHANGE
88959599516SKenneth E. Jansen!      call closefile( igeomBAK,  "read" )
89059599516SKenneth E. Jansenc
89159599516SKenneth E. Jansen      deallocate(xread)
89259599516SKenneth E. Jansen      if ( numpbc > 0 )  then
89359599516SKenneth E. Jansen         deallocate(bcinpread)
89459599516SKenneth E. Jansen         deallocate(ibctmpread)
89559599516SKenneth E. Jansen      endif
89659599516SKenneth E. Jansen      deallocate(iperread)
89759599516SKenneth E. Jansen      if(numpe.gt.1)
89859599516SKenneth E. Jansen     &     deallocate(ilworkread)
89959599516SKenneth E. Jansen      deallocate(nbcread)
90059599516SKenneth E. Jansen
90159599516SKenneth E. Jansen      return
90259599516SKenneth E. Jansenc
90359599516SKenneth E. Jansen 994  call error ('input   ','opening ', igeomBAK)
90459599516SKenneth E. Jansen 995  call error ('input   ','opening ', igeomBAK)
90559599516SKenneth E. Jansen 997  call error ('input   ','end file', igeomBAK)
90659599516SKenneth E. Jansen 998  call error ('input   ','end file', igeomBAK)
90759599516SKenneth E. Jansenc
90859599516SKenneth E. Jansen      end
90959599516SKenneth E. Jansen
91059599516SKenneth E. Jansenc
91159599516SKenneth E. Jansenc No longer called but kept around in case....
91259599516SKenneth E. Jansenc
91359599516SKenneth E. Jansen      subroutine genpzero(iBC)
91459599516SKenneth E. Jansen
91559599516SKenneth E. Jansen      use pointer_data
91659599516SKenneth E. Jansenc
91759599516SKenneth E. Jansen      include "common.h"
91859599516SKenneth E. Jansen      integer iBC(nshg)
91959599516SKenneth E. Jansenc
92059599516SKenneth E. Jansenc....  check to see if any of the nodes have a dirichlet pressure
92159599516SKenneth E. Jansenc
92259599516SKenneth E. Jansen      pzero=1
92359599516SKenneth E. Jansen      if (any(btest(iBC,2))) pzero=0
92459599516SKenneth E. Jansenc
92559599516SKenneth E. Jansen      do iblk = 1, nelblb
92659599516SKenneth E. Jansen         npro = lcblkb(1,iblk+1)-lcblkb(1,iblk)
92759599516SKenneth E. Jansen         do i=1, npro
92859599516SKenneth E. Jansen            iBCB1=miBCB(iblk)%p(i,1)
92959599516SKenneth E. Jansenc
93059599516SKenneth E. Jansenc.... check to see if any of the nodes have a Neumann pressure
93159599516SKenneth E. Jansenc     but not periodic (note that
93259599516SKenneth E. Jansenc
93359599516SKenneth E. Jansen            if(btest(iBCB1,1)) pzero=0
93459599516SKenneth E. Jansen         enddo
93559599516SKenneth E. Jansenc
93659599516SKenneth E. Jansenc.... share results with other processors
93759599516SKenneth E. Jansenc
93859599516SKenneth E. Jansen         pzl=pzero
93959599516SKenneth E. Jansen         if (numpe .gt. 1)
94059599516SKenneth E. Jansen     &        call MPI_ALLREDUCE (pzl, pzero, 1,
94159599516SKenneth E. Jansen     &        MPI_DOUBLE_PRECISION,MPI_MIN, MPI_COMM_WORLD,ierr)
94259599516SKenneth E. Jansen
94359599516SKenneth E. Jansen      enddo
94459599516SKenneth E. Jansenc
94559599516SKenneth E. Jansenc.... return
94659599516SKenneth E. Jansenc
94759599516SKenneth E. Jansen      return
94859599516SKenneth E. Jansenc
94959599516SKenneth E. Jansen      end
95059599516SKenneth E. Jansen
951