xref: /phasta/phSolver/common/readnblk.f (revision 513954ef803c300cddba2bb96b4a5dac0b93489a)
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
329071d3baSCameron Smith      use phstr
33ab645d52SCameron Smith      use syncio
34ab645d52SCameron Smith      use posixio
35ab645d52SCameron Smith      use streamio
3659599516SKenneth E. Jansen      include "common.h"
3702ce253eSCameron Smith
389a6935e5SKenneth E. Jansen      real*8, target, allocatable :: xread(:,:), qread(:,:), acread(:,:)
399a6935e5SKenneth E. Jansen      real*8, target, allocatable :: uread(:,:)
409a6935e5SKenneth E. Jansen      real*8, target, allocatable :: BCinpread(:,:)
41266200f9SCameron Smith      real*8 :: iotime
429a6935e5SKenneth E. Jansen      integer, target, allocatable :: iperread(:), iBCtmpread(:)
439a6935e5SKenneth E. Jansen      integer, target, allocatable :: ilworkread(:), nBCread(:)
44*513954efSKenneth E. Jansen      integer, target, allocatable :: fncorpread(:)
45*513954efSKenneth E. Jansen      integer fncorpsize
46*513954efSKenneth E. Jansen      character*10 cname2, cname2nd
4759599516SKenneth E. Jansen      character*8 mach2
4859599516SKenneth E. Jansen      character*30 fmt1
4959599516SKenneth E. Jansen      character*255 fname1,fnamer,fnamelr
5059599516SKenneth E. Jansen      character*255 warning
5159599516SKenneth E. Jansen      integer igeomBAK, ibndc, irstin, ierr
529a6935e5SKenneth E. Jansen      integer, target :: intfromfile(50) ! integers read from headers
539f4aafb6SCameron Smith      integer :: descriptor, descriptorG, GPID, color, nfields
5459599516SKenneth E. Jansen      integer ::  numparts, nppf
5559599516SKenneth E. Jansen      integer :: ierr_io, numprocs, itmp, itmp2
56ade0e30fSCameron Smith      integer :: ignored
57d7abaf6cSCameron Smith      integer :: fileFmt
5859599516SKenneth E. Jansen      character*255 fname2, temp2
5959599516SKenneth E. Jansen      character*64 temp1
60e5afe575SCameron Smith      type(c_ptr) :: handle
61e5afe575SCameron Smith      character(len=1024) :: dataInt, dataDbl
62e5afe575SCameron Smith      dataInt = c_char_'integer'//c_null_char
63e5afe575SCameron Smith      dataDbl = c_char_'double'//c_null_char
6459599516SKenneth E. Jansenc
6559599516SKenneth E. Jansenc.... determine the step number to start with
6659599516SKenneth E. Jansenc
6759599516SKenneth E. Jansen      open(unit=72,file='numstart.dat',status='old')
6859599516SKenneth E. Jansen      read(72,*) irstart
6959599516SKenneth E. Jansen      close(72)
7059599516SKenneth E. Jansen      lstep=irstart ! in case restart files have no fields
7159599516SKenneth E. Jansen
7259599516SKenneth E. Jansen      fname1='geombc.dat'
7359599516SKenneth E. Jansen      fname1= trim(fname1)  // cname2(myrank+1)
7459599516SKenneth E. Jansen
7559599516SKenneth E. Jansen      itmp=1
7659599516SKenneth E. Jansen      if (irstart .gt. 0) itmp = int(log10(float(irstart)))+1
7759599516SKenneth E. Jansen      write (fmt1,"('(''restart.'',i',i1,',1x)')") itmp
7859599516SKenneth E. Jansen      write (fnamer,fmt1) irstart
7959599516SKenneth E. Jansen      fnamer = trim(fnamer) // cname2(myrank+1)
8059599516SKenneth E. Jansenc
8159599516SKenneth E. Jansenc.... open input files
8259599516SKenneth E. Jansenc.... input the geometry parameters
8359599516SKenneth E. Jansenc
8459599516SKenneth E. Jansen      numparts = numpe !This is the common settings. Beware if you try to compute several parts per process
8559599516SKenneth E. Jansen
8659599516SKenneth E. Jansen      itwo=2
8759599516SKenneth E. Jansen      ione=1
8859599516SKenneth E. Jansen      ieleven=11
8959599516SKenneth E. Jansen      itmp = int(log10(float(myrank+1)))+1
9059599516SKenneth E. Jansen
91266200f9SCameron Smith      iotime = TMRC()
929f4aafb6SCameron Smith      if( input_mode .eq. -1 ) then
93a486e66cSCameron Smith        call streamio_setup_read(fhandle, geomRestartStream)
949f4aafb6SCameron Smith      else if( input_mode .eq. 0 ) then
95ab645d52SCameron Smith        call posixio_setup(fhandle, c_char_'r')
969f4aafb6SCameron Smith      else if( input_mode .ge. 1 ) then
97ab645d52SCameron Smith        call syncio_setup_read(nsynciofiles, fhandle)
982efdc748SKenneth E. Jansen      end if
99a93de25bSCameron Smith      call phio_constructName(fhandle,
100d7abaf6cSCameron Smith     &        c_char_'geombc' // char(0), fname1)
101d7abaf6cSCameron Smith      call phio_openfile(fname1, fhandle);
10259599516SKenneth E. Jansen
103d5d2f64dSCameron Smith      call phio_readheader(fhandle,c_char_'number of nodes' // char(0),
104e5afe575SCameron Smith     & c_loc(numnp),ione, dataInt, iotype)
10559599516SKenneth E. Jansen
106d5d2f64dSCameron Smith      call phio_readheader(fhandle,c_char_'number of modes' // char(0),
107e5afe575SCameron Smith     & c_loc(nshg),ione, dataInt, iotype)
10859599516SKenneth E. Jansen
109d5d2f64dSCameron Smith      call phio_readheader(fhandle,
110e5afe575SCameron Smith     &  c_char_'number of interior elements' // char(0),
111e5afe575SCameron Smith     &  c_loc(numel),ione, dataInt, iotype)
11259599516SKenneth E. Jansen
113d5d2f64dSCameron Smith      call phio_readheader(fhandle,
114e5afe575SCameron Smith     &  c_char_'number of boundary elements' // char(0),
115e5afe575SCameron Smith     &  c_loc(numelb),ione, dataInt, iotype)
11659599516SKenneth E. Jansen
117d5d2f64dSCameron Smith      call phio_readheader(fhandle,
118e5afe575SCameron Smith     &  c_char_'maximum number of element nodes' // char(0),
119e5afe575SCameron Smith     &  c_loc(nen),ione, dataInt, iotype)
12059599516SKenneth E. Jansen
121d5d2f64dSCameron Smith      call phio_readheader(fhandle,
122e5afe575SCameron Smith     &  c_char_'number of interior tpblocks' // char(0),
123e5afe575SCameron Smith     &  c_loc(nelblk),ione, dataInt, iotype)
12459599516SKenneth E. Jansen
125d5d2f64dSCameron Smith      call phio_readheader(fhandle,
126e5afe575SCameron Smith     & c_char_'number of boundary tpblocks' // char(0),
127e5afe575SCameron Smith     & c_loc(nelblb),ione, dataInt, iotype)
12859599516SKenneth E. Jansen
129d5d2f64dSCameron Smith      call phio_readheader(fhandle,
130e5afe575SCameron Smith     & c_char_'number of nodes with Dirichlet BCs' // char(0),
131e5afe575SCameron Smith     & c_loc(numpbc),ione, dataInt, iotype)
13259599516SKenneth E. Jansen
133d5d2f64dSCameron Smith      call phio_readheader(fhandle,
134e5afe575SCameron Smith     & c_char_'number of shape functions' // char(0),
135e5afe575SCameron Smith     & c_loc(ntopsh),ione, dataInt, iotype)
13659599516SKenneth E. Jansenc
13759599516SKenneth E. Jansenc.... calculate the maximum number of boundary element nodes
13859599516SKenneth E. Jansenc
13959599516SKenneth E. Jansen      nenb = 0
14059599516SKenneth E. Jansen      do i = 1, melCat
14159599516SKenneth E. Jansen         if (nen .eq. nenCat(i,nsd)) nenb = max(nenCat(i,nsd-1), nenb)
14259599516SKenneth E. Jansen      enddo
14359599516SKenneth E. Jansen      if (myrank == master) then
14459599516SKenneth E. Jansen         if (nenb .eq. 0) call error ('input   ','nen     ',nen)
14559599516SKenneth E. Jansen      endif
14659599516SKenneth E. Jansenc
14759599516SKenneth E. Jansenc.... setup some useful constants
14859599516SKenneth E. Jansenc
14959599516SKenneth E. Jansen      I3nsd  = nsd / 3          ! nsd=3 integer flag
15059599516SKenneth E. Jansen      E3nsd  = float(I3nsd)     ! nsd=3 real    flag
15159599516SKenneth E. Jansen      if(matflg(1,1).lt.0) then
15259599516SKenneth E. Jansen         nflow = nsd + 1
15359599516SKenneth E. Jansen      else
15459599516SKenneth E. Jansen         nflow = nsd + 2
15559599516SKenneth E. Jansen      endif
15659599516SKenneth E. Jansen      ndof   = nsd + 2
15759599516SKenneth E. Jansen      nsclr=impl(1)/100
15859599516SKenneth E. Jansen      ndof=ndof+nsclr           ! number of sclr transport equations to solve
15959599516SKenneth E. Jansen
16059599516SKenneth E. Jansen      ndofBC = ndof + I3nsd     ! dimension of BC array
16159599516SKenneth E. Jansen      ndiBCB = 2                ! dimension of iBCB array
16259599516SKenneth E. Jansen      ndBCB  = ndof + 1         ! dimension of BCB array
16359599516SKenneth E. Jansen      nsymdf = (ndof*(ndof + 1)) / 2 ! symm. d.o.f.'s
16459599516SKenneth E. Jansenc
16559599516SKenneth E. Jansenc.... ----------------------> Communication tasks <--------------------
16659599516SKenneth E. Jansenc
16759599516SKenneth E. Jansen      if(numpe > 1) then
168d5d2f64dSCameron Smith         call phio_readheader(fhandle,
169e5afe575SCameron Smith     &    c_char_'size of ilwork array' // char(0),
170e5afe575SCameron Smith     &    c_loc(nlwork),ione, dataInt, iotype)
17159599516SKenneth E. Jansen
172d5d2f64dSCameron Smith         call phio_readheader(fhandle,
173e5afe575SCameron Smith     &    c_char_'ilwork' //char(0),
174e5afe575SCameron Smith     &    c_loc(nlwork),ione, dataInt, iotype)
17559599516SKenneth E. Jansen
17659599516SKenneth E. Jansen         allocate( point2ilwork(nlwork) )
17759599516SKenneth E. Jansen         allocate( ilworkread(nlwork) )
178d5d2f64dSCameron Smith         call phio_readdatablock(fhandle, c_char_'ilwork' // char(0),
179bc62cfd4SCameron Smith     &      c_loc(ilworkread), nlwork, dataInt , iotype)
18059599516SKenneth E. Jansen
18159599516SKenneth E. Jansen         point2ilwork = ilworkread
18259599516SKenneth E. Jansen         call ctypes (point2ilwork)
183*513954efSKenneth E. Jansen
184*513954efSKenneth E. Jansen       if(usingPETSc.eq.1) then
185*513954efSKenneth E. Jansen         fncorpsize = nshg
186*513954efSKenneth E. Jansen         allocate(fncorp(fncorpsize))
187*513954efSKenneth E. Jansen         call gen_ncorp(fncorp, ilworkread, nlwork, fncorpsize)
188*513954efSKenneth E. Jansen!
189*513954efSKenneth E. Jansen! the  following code finds the global range of the owned nodes
190*513954efSKenneth E. Jansen!
191*513954efSKenneth E. Jansen         maxowned=0
192*513954efSKenneth E. Jansen         minowned=maxval(fncorp)
193*513954efSKenneth E. Jansen         do i = 1,nshg
194*513954efSKenneth E. Jansen            if(fncorp(i).gt.0) then  ! don't consider remote copies
195*513954efSKenneth E. Jansen               maxowned=max(maxowned,fncorp(i))
196*513954efSKenneth E. Jansen               minowned=min(minowned,fncorp(i))
197*513954efSKenneth E. Jansen            endif
198*513954efSKenneth E. Jansen         enddo
199*513954efSKenneth E. Jansen!
200*513954efSKenneth E. Jansen!  end of global range code
201*513954efSKenneth E. Jansen!
202*513954efSKenneth E. Jansen         call commuInt(fncorp, point2ilwork, 1, 'out')
203*513954efSKenneth E. Jansen         ncorpsize = fncorpsize
204*513954efSKenneth E. Jansen       endif
20559599516SKenneth E. Jansen      else
20659599516SKenneth E. Jansen           nlwork=1
20759599516SKenneth E. Jansen           allocate( point2ilwork(1))
20859599516SKenneth E. Jansen           nshg0 = nshg
20959599516SKenneth E. Jansen      endif
21059599516SKenneth E. Jansen
21159599516SKenneth E. Jansen      itwo=2
21259599516SKenneth E. Jansen
213d5d2f64dSCameron Smith      call phio_readheader(fhandle,
214e5afe575SCameron Smith     & c_char_'co-ordinates' // char(0),
215e5afe575SCameron Smith     & c_loc(intfromfile),itwo, dataDbl, iotype)
21659599516SKenneth E. Jansen      numnp=intfromfile(1)
21759599516SKenneth E. Jansen      allocate( point2x(numnp,nsd) )
21859599516SKenneth E. Jansen      allocate( xread(numnp,nsd) )
21959599516SKenneth E. Jansen      ixsiz=numnp*nsd
220d5d2f64dSCameron Smith      call phio_readdatablock(fhandle,
221bc62cfd4SCameron Smith     & c_char_'co-ordinates' // char(0),
222bc62cfd4SCameron Smith     & c_loc(xread),ixsiz, dataDbl, iotype)
22359599516SKenneth E. Jansen      point2x = xread
224*513954efSKenneth E. Jansen
225*513954efSKenneth E. Jansenc..............................for Duct
226*513954efSKenneth E. Jansen      if(istretchOutlet.eq.1)then
227*513954efSKenneth E. Jansen
228*513954efSKenneth E. Jansenc...geometry6
229*513954efSKenneth E. Jansen        if(iDuctgeometryType .eq. 6) then
230*513954efSKenneth E. Jansen          xmaxn = 1.276
231*513954efSKenneth E. Jansen          xmaxo = 0.848
232*513954efSKenneth E. Jansen          xmin  = 0.42
233*513954efSKenneth E. Jansenc...geometry8
234*513954efSKenneth E. Jansen        elseif(iDuctgeometryType .eq. 8)then
235*513954efSKenneth E. Jansen          xmaxn=1.6*4.5*0.0254+0.85*1.5
236*513954efSKenneth E. Jansen          xmaxo=1.6*4.5*0.0254+0.85*1.0
237*513954efSKenneth E. Jansen          xmin =1.6*4.5*0.0254+0.85*0.5
238*513954efSKenneth E. Jansen        endif
239*513954efSKenneth E. Jansenc...
240*513954efSKenneth E. Jansen        alpha=(xmaxn-xmaxo)/(xmaxo-xmin)**2
241*513954efSKenneth E. Jansen        where (point2x(:,1) .ge. xmin)
242*513954efSKenneth E. Jansenc..... N=# of current elements from .42 to exit(~40)
243*513954efSKenneth E. Jansenc..... (x_mx-x_mn)/N=.025
244*513954efSKenneth E. Jansenc..... alpha=3    3*.025=.075
245*513954efSKenneth E. Jansen           point2x(:,1)=point2x(:,1)+
246*513954efSKenneth E. Jansen     &     alpha*(point2x(:,1)-xmin)**2
247*513954efSKenneth E. Jansenc..... ftn to stretch x at exit
248*513954efSKenneth E. Jansen        endwhere
249*513954efSKenneth E. Jansen      endif
250*513954efSKenneth E. Jansen
25159599516SKenneth E. Jansenc
25259599516SKenneth E. Jansenc.... read in and block out the connectivity
25359599516SKenneth E. Jansenc
25459599516SKenneth E. Jansen      call genblk (IBKSIZ)
25559599516SKenneth E. Jansenc
25659599516SKenneth E. Jansenc.... read the boundary condition mapping array
25759599516SKenneth E. Jansenc
25859599516SKenneth E. Jansen      ione=1
259d5d2f64dSCameron Smith      call phio_readheader(fhandle,
260e5afe575SCameron Smith     & c_char_'bc mapping array' // char(0),
261e5afe575SCameron Smith     & c_loc(nshg),ione, dataInt, iotype)
26259599516SKenneth E. Jansen
26359599516SKenneth E. Jansen      allocate( nBC(nshg) )
26459599516SKenneth E. Jansen
26559599516SKenneth E. Jansen      allocate( nBCread(nshg) )
26659599516SKenneth E. Jansen
267d5d2f64dSCameron Smith      call phio_readdatablock(fhandle,
268bc62cfd4SCameron Smith     & c_char_'bc mapping array' // char(0),
269bc62cfd4SCameron Smith     & c_loc(nBCread), nshg, dataInt, iotype)
27059599516SKenneth E. Jansen
27159599516SKenneth E. Jansen      nBC=nBCread
27259599516SKenneth E. Jansenc
27359599516SKenneth E. Jansenc.... read the temporary iBC array
27459599516SKenneth E. Jansenc
27559599516SKenneth E. Jansen      ione=1
276d5d2f64dSCameron Smith      call phio_readheader(fhandle,
277e5afe575SCameron Smith     & c_char_'bc codes array' // char(0),
278e5afe575SCameron Smith     & c_loc(numpbc),ione, dataInt, iotype)
27959599516SKenneth E. Jansen
28059599516SKenneth E. Jansen      if ( numpbc > 0 ) then
28159599516SKenneth E. Jansen        allocate( iBCtmp(numpbc) )
28259599516SKenneth E. Jansen        allocate( iBCtmpread(numpbc) )
28359599516SKenneth E. Jansen      else
28459599516SKenneth E. Jansen        allocate( iBCtmp(1) )
28559599516SKenneth E. Jansen        allocate( iBCtmpread(1) )
28659599516SKenneth E. Jansen      endif
287d5d2f64dSCameron Smith      call phio_readdatablock(fhandle,
288bc62cfd4SCameron Smith     & c_char_'bc codes array' // char(0),
289bc62cfd4SCameron Smith     & c_loc(iBCtmpread), numpbc, dataInt, iotype)
29059599516SKenneth E. Jansen
29159599516SKenneth E. Jansen      if ( numpbc > 0 ) then
29259599516SKenneth E. Jansen         iBCtmp=iBCtmpread
29359599516SKenneth E. Jansen      else  ! sometimes a partition has no BC's
29459599516SKenneth E. Jansen         deallocate( iBCtmpread)
29559599516SKenneth E. Jansen         iBCtmp=0
29659599516SKenneth E. Jansen      endif
29759599516SKenneth E. Jansenc
29859599516SKenneth E. Jansenc.... read boundary condition data
29959599516SKenneth E. Jansenc
30059599516SKenneth E. Jansen      ione=1
301d5d2f64dSCameron Smith      call phio_readheader(fhandle,
302e5afe575SCameron Smith     & c_char_'boundary condition array' // char(0),
303e5afe575SCameron Smith     & c_loc(intfromfile),ione, dataDbl, iotype)
30459599516SKenneth E. Jansen
30559599516SKenneth E. Jansen      if ( numpbc > 0 ) then
30659599516SKenneth E. Jansen         allocate( BCinp(numpbc,ndof+7) )
30759599516SKenneth E. Jansen         nsecondrank=intfromfile(1)/numpbc
30859599516SKenneth E. Jansen         allocate( BCinpread(numpbc,nsecondrank) )
30959599516SKenneth E. Jansen         iBCinpsiz=intfromfile(1)
31059599516SKenneth E. Jansen      else
31159599516SKenneth E. Jansen         allocate( BCinp(1,ndof+7) )
31259599516SKenneth E. Jansen         allocate( BCinpread(0,0) ) !dummy
31359599516SKenneth E. Jansen         iBCinpsiz=intfromfile(1)
31459599516SKenneth E. Jansen      endif
31559599516SKenneth E. Jansen
316d5d2f64dSCameron Smith      call phio_readdatablock(fhandle,
317bc62cfd4SCameron Smith     & c_char_'boundary condition array' // char(0),
318bc62cfd4SCameron Smith     & c_loc(BCinpread), iBCinpsiz, dataDbl, iotype)
31959599516SKenneth E. Jansen
32059599516SKenneth E. Jansen      if ( numpbc > 0 ) then
32159599516SKenneth E. Jansen         BCinp(:,1:(ndof+7))=BCinpread(:,1:(ndof+7))
32259599516SKenneth E. Jansen      else  ! sometimes a partition has no BC's
32359599516SKenneth E. Jansen         deallocate(BCinpread)
32459599516SKenneth E. Jansen         BCinp=0
32559599516SKenneth E. Jansen      endif
32659599516SKenneth E. Jansenc
32759599516SKenneth E. Jansenc.... read periodic boundary conditions
32859599516SKenneth E. Jansenc
32959599516SKenneth E. Jansen      ione=1
330d5d2f64dSCameron Smith      call phio_readheader(fhandle,
331e5afe575SCameron Smith     & c_char_'periodic masters array' // char(0),
332e5afe575SCameron Smith     & c_loc(nshg), ione, dataInt, iotype)
33359599516SKenneth E. Jansen      allocate( point2iper(nshg) )
33459599516SKenneth E. Jansen      allocate( iperread(nshg) )
335d5d2f64dSCameron Smith      call phio_readdatablock(fhandle,
336bc62cfd4SCameron Smith     & c_char_'periodic masters array' // char(0),
337bc62cfd4SCameron Smith     & c_loc(iperread), nshg, dataInt, iotype)
33859599516SKenneth E. Jansen      point2iper=iperread
33959599516SKenneth E. Jansenc
34059599516SKenneth E. Jansenc.... generate the boundary element blocks
34159599516SKenneth E. Jansenc
34259599516SKenneth E. Jansen      call genbkb (ibksiz)
34359599516SKenneth E. Jansenc
34459599516SKenneth E. Jansenc  Read in the nsons and ifath arrays if needed
34559599516SKenneth E. Jansenc
34659599516SKenneth E. Jansenc  There is a fundamental shift in the meaning of ifath based on whether
34759599516SKenneth E. Jansenc  there exist homogenous directions in the flow.
34859599516SKenneth E. Jansenc
34959599516SKenneth E. Jansenc  HOMOGENOUS DIRECTIONS EXIST:  Here nfath is the number of inhomogenous
35059599516SKenneth E. Jansenc  points in the TOTAL mesh.  That is to say that each partition keeps a
35159599516SKenneth E. Jansenc  link to  ALL inhomogenous points.  This link is furthermore not to the
35259599516SKenneth E. Jansenc  sms numbering but to the original structured grid numbering.  These
35359599516SKenneth E. Jansenc  inhomogenous points are thought of as fathers, with their sons being all
35459599516SKenneth E. Jansenc  the points in the homogenous directions that have this father's
35559599516SKenneth E. Jansenc  inhomogeneity.  The array ifath takes as an arguement the sms numbering
35659599516SKenneth E. Jansenc  and returns as a result the father.
35759599516SKenneth E. Jansenc
35859599516SKenneth E. Jansenc  In this case nsons is the number of sons that each father has and ifath
35959599516SKenneth E. Jansenc  is an array which tells the
36059599516SKenneth E. Jansenc
36159599516SKenneth E. Jansenc  NO HOMOGENOUS DIRECTIONS.  In this case the mesh would grow to rapidly
36259599516SKenneth E. Jansenc  if we followed the above strategy since every partition would index its
36359599516SKenneth E. Jansenc  points to the ENTIRE mesh.  Furthermore, there would never be a need
36459599516SKenneth E. Jansenc  to average to a node off processor since there is no spatial averaging.
36559599516SKenneth E. Jansenc  Therefore, to properly account for this case we must recognize it and
36659599516SKenneth E. Jansenc  inerrupt certain actions (i.e. assembly of the average across partitions).
36759599516SKenneth E. Jansenc  This case is easily identified by noting that maxval(nsons) =1 (i.e. no
36859599516SKenneth E. Jansenc  father has any sons).  Reiterating to be clear, in this case ifath does
36959599516SKenneth E. Jansenc  not point to a global numbering but instead just points to itself.
37059599516SKenneth E. Jansenc
37159599516SKenneth E. Jansen      nfath=1  ! some architectures choke on a zero or undeclared
37259599516SKenneth E. Jansen                 ! dimension variable.  This sets it to a safe, small value.
37359599516SKenneth E. Jansen      if(((iLES .lt. 20) .and. (iLES.gt.0))
37459599516SKenneth E. Jansen     &                   .or. (itwmod.gt.0)  ) then ! don't forget same
37559599516SKenneth E. Jansen                                                    ! conditional in proces.f
37659599516SKenneth E. Jansen                                                    ! needed for alloc
37759599516SKenneth E. Jansen         ione=1
37859599516SKenneth E. Jansen         if(nohomog.gt.0) then
379d5d2f64dSCameron Smith            call phio_readheader(fhandle,
380e5afe575SCameron Smith     &       c_char_'number of father-nodes' // char(0),
381e5afe575SCameron Smith     &       c_loc(nfath), ione, dataInt, iotype)
38259599516SKenneth E. Jansen
383d5d2f64dSCameron Smith            call phio_readheader(fhandle,
384e5afe575SCameron Smith     &       c_char_'number of son-nodes for each father' // char(0),
385e5afe575SCameron Smith     &       c_loc(nfath), ione, dataInt, iotype)
38659599516SKenneth E. Jansen
38759599516SKenneth E. Jansen            allocate (point2nsons(nfath))
38859599516SKenneth E. Jansen
389d5d2f64dSCameron Smith            call phio_readdatablock(fhandle,
390bc62cfd4SCameron Smith     &       c_char_'number of son-nodes for each father' // char(0),
391bc62cfd4SCameron Smith     &       c_loc(point2nsons),nfath, dataInt, iotype)
39259599516SKenneth E. Jansen
393d5d2f64dSCameron Smith            call phio_readheader(fhandle,
394e5afe575SCameron Smith     &       c_char_'keyword ifath' // char(0),
395e5afe575SCameron Smith     &       c_loc(nshg), ione, dataInt, iotype);
39659599516SKenneth E. Jansen
39759599516SKenneth E. Jansen            allocate (point2ifath(nshg))
39859599516SKenneth E. Jansen
399d5d2f64dSCameron Smith            call phio_readdatablock(fhandle,
400bc62cfd4SCameron Smith     &       c_char_'keyword ifath' // char(0),
401bc62cfd4SCameron Smith     &       c_loc(point2ifath), nshg, dataInt, iotype)
40259599516SKenneth E. Jansen
40359599516SKenneth E. Jansen            nsonmax=maxval(point2nsons)
40459599516SKenneth E. Jansen         else  ! this is the case where there is no homogeneity
40559599516SKenneth E. Jansen               ! therefore ever node is a father (too itself).  sonfath
40659599516SKenneth E. Jansen               ! (a routine in NSpre) will set this up but this gives
40759599516SKenneth E. Jansen               ! you an option to avoid that.
40859599516SKenneth E. Jansen            nfath=nshg
40959599516SKenneth E. Jansen            allocate (point2nsons(nfath))
41059599516SKenneth E. Jansen            point2nsons=1
41159599516SKenneth E. Jansen            allocate (point2ifath(nshg))
41259599516SKenneth E. Jansen            do i=1,nshg
41359599516SKenneth E. Jansen               point2ifath(i)=i
41459599516SKenneth E. Jansen            enddo
41559599516SKenneth E. Jansen            nsonmax=1
41659599516SKenneth E. Jansen         endif
41759599516SKenneth E. Jansen      else
41859599516SKenneth E. Jansen         allocate (point2nsons(1))
41959599516SKenneth E. Jansen         allocate (point2ifath(1))
42059599516SKenneth E. Jansen      endif
42159599516SKenneth E. Jansen
422d7abaf6cSCameron Smith      call phio_closefile(fhandle);
423266200f9SCameron Smith      iotime = TMRC() - iotime
424266200f9SCameron Smith      if (myrank.eq.master) then
425266200f9SCameron Smith        write(*,*) 'time to read geombc (seconds)', iotime
426266200f9SCameron Smith      endif
427266200f9SCameron Smith
42859599516SKenneth E. Jansenc.... Read restart files
429266200f9SCameron Smith      iotime = TMRC()
4309f4aafb6SCameron Smith      if( input_mode .eq. -1 ) then
431a486e66cSCameron Smith        call streamio_setup_read(fhandle, geomRestartStream)
4329f4aafb6SCameron Smith      else if( input_mode .eq. 0 ) then
433d7abaf6cSCameron Smith        call posixio_setup(fhandle, c_char_'r')
4349f4aafb6SCameron Smith      else if( input_mode .ge. 1 ) then
435d7abaf6cSCameron Smith        call syncio_setup_read(nsynciofiles, fhandle)
436d7abaf6cSCameron Smith      end if
437a93de25bSCameron Smith      call phio_constructName(fhandle,
438d7abaf6cSCameron Smith     &        c_char_'restart' // char(0), fnamer)
4399071d3baSCameron Smith      call phstr_appendInt(fnamer, irstart)
4409071d3baSCameron Smith      call phstr_appendStr(fnamer, c_char_'.'//c_null_char)
441d7abaf6cSCameron Smith      call phio_openfile(fnamer, fhandle);
44259599516SKenneth E. Jansen
44359599516SKenneth E. Jansen      ithree=3
44459599516SKenneth E. Jansen
44559599516SKenneth E. Jansen      itmp = int(log10(float(myrank+1)))+1
44659599516SKenneth E. Jansen
44759599516SKenneth E. Jansen      intfromfile=0
448d5d2f64dSCameron Smith      call phio_readheader(fhandle,
449e5afe575SCameron Smith     & c_char_'solution' // char(0),
450e5afe575SCameron Smith     & c_loc(intfromfile), ithree, dataInt, iotype)
45159599516SKenneth E. Jansenc
45259599516SKenneth E. Jansenc.... read the values of primitive variables into q
45359599516SKenneth E. Jansenc
45459599516SKenneth E. Jansen      allocate( qold(nshg,ndof) )
45559599516SKenneth E. Jansen      if(intfromfile(1).ne.0) then
45659599516SKenneth E. Jansen         nshg2=intfromfile(1)
45759599516SKenneth E. Jansen         ndof2=intfromfile(2)
45859599516SKenneth E. Jansen         lstep=intfromfile(3)
45959599516SKenneth E. Jansen         if(ndof2.ne.ndof) then
46059599516SKenneth E. Jansen
46159599516SKenneth E. Jansen         endif
46259599516SKenneth E. Jansen        if (nshg2 .ne. nshg)
46359599516SKenneth E. Jansen     &        call error ('restar  ', 'nshg   ', nshg)
46459599516SKenneth E. Jansen         allocate( qread(nshg,ndof2) )
46559599516SKenneth E. Jansen         iqsiz=nshg*ndof2
466d5d2f64dSCameron Smith         call phio_readdatablock(fhandle,
467bc62cfd4SCameron Smith     &    c_char_'solution' // char(0),
468bc62cfd4SCameron Smith     &    c_loc(qread),iqsiz, dataDbl,iotype)
46959599516SKenneth E. Jansen         qold(:,1:ndof)=qread(:,1:ndof)
47059599516SKenneth E. Jansen         deallocate(qread)
47159599516SKenneth E. Jansen      else
47259599516SKenneth E. Jansen         if (myrank.eq.master) then
47359599516SKenneth E. Jansen            if (matflg(1,1).eq.0) then ! compressible
47459599516SKenneth E. Jansen               warning='Solution is set to zero (with p and T to one)'
47559599516SKenneth E. Jansen            else
47659599516SKenneth E. Jansen               warning='Solution is set to zero'
47759599516SKenneth E. Jansen            endif
47859599516SKenneth E. Jansen            write(*,*) warning
47959599516SKenneth E. Jansen         endif
48059599516SKenneth E. Jansen         qold=zero
48159599516SKenneth E. Jansen         if (matflg(1,1).eq.0) then ! compressible
48259599516SKenneth E. Jansen            qold(:,1)=one ! avoid zero pressure
48359599516SKenneth E. Jansen            qold(:,nflow)=one ! avoid zero temperature
48459599516SKenneth E. Jansen         endif
48559599516SKenneth E. Jansen      endif
48659599516SKenneth E. Jansen
48759599516SKenneth E. Jansen      intfromfile=0
488d5d2f64dSCameron Smith      call phio_readheader(fhandle,
489e5afe575SCameron Smith     & c_char_'time derivative of solution' // char(0),
490e5afe575SCameron Smith     & c_loc(intfromfile), ithree, dataInt, iotype)
49159599516SKenneth E. Jansen      allocate( acold(nshg,ndof) )
49259599516SKenneth E. Jansen      if(intfromfile(1).ne.0) then
49359599516SKenneth E. Jansen         nshg2=intfromfile(1)
49459599516SKenneth E. Jansen         ndof2=intfromfile(2)
49559599516SKenneth E. Jansen         lstep=intfromfile(3)
49659599516SKenneth E. Jansen
49759599516SKenneth E. Jansen         if (nshg2 .ne. nshg)
49859599516SKenneth E. Jansen     &        call error ('restar  ', 'nshg   ', nshg)
49959599516SKenneth E. Jansen         allocate( acread(nshg,ndof2) )
50059599516SKenneth E. Jansen         acread=zero
50159599516SKenneth E. Jansen         iacsiz=nshg*ndof2
502d5d2f64dSCameron Smith         call phio_readdatablock(fhandle,
503bc62cfd4SCameron Smith     &    c_char_'time derivative of solution' // char(0),
504bc62cfd4SCameron Smith     &    c_loc(acread), iacsiz, dataDbl,iotype)
50559599516SKenneth E. Jansen         acold(:,1:ndof)=acread(:,1:ndof)
50659599516SKenneth E. Jansen         deallocate(acread)
50759599516SKenneth E. Jansen      else
50859599516SKenneth E. Jansen         if (myrank.eq.master) then
50959599516SKenneth E. Jansen            warning='Time derivative of solution is set to zero (SAFE)'
51059599516SKenneth E. Jansen            write(*,*) warning
51159599516SKenneth E. Jansen         endif
51259599516SKenneth E. Jansen         acold=zero
51359599516SKenneth E. Jansen      endif
51459599516SKenneth E. Jansencc
51559599516SKenneth E. Jansencc.... read the header and check it against the run data
51659599516SKenneth E. Jansencc
51759599516SKenneth E. Jansen      if (ideformwall.eq.1) then
51859599516SKenneth E. Jansen
51959599516SKenneth E. Jansen          intfromfile=0
520d5d2f64dSCameron Smith          call phio_readheader(fhandle,
521e5afe575SCameron Smith     &     c_char_'displacement' // char(0),
522e5afe575SCameron Smith     &     c_loc(intfromfile), ithree, dataInt, iotype)
52359599516SKenneth E. Jansen
52459599516SKenneth E. Jansen         nshg2=intfromfile(1)
52559599516SKenneth E. Jansen         ndisp=intfromfile(2)
52659599516SKenneth E. Jansen         lstep=intfromfile(3)
52759599516SKenneth E. Jansen         if(ndisp.ne.nsd) then
52859599516SKenneth E. Jansen            warning='WARNING ndisp not equal nsd'
52959599516SKenneth E. Jansen            write(*,*) warning , ndisp
53059599516SKenneth E. Jansen         endif
53159599516SKenneth E. Jansen         if (nshg2 .ne. nshg)
53259599516SKenneth E. Jansen     &        call error ('restar  ', 'nshg   ', nshg)
53359599516SKenneth E. Jansenc
53459599516SKenneth E. Jansenc.... read the values of primitive variables into uold
53559599516SKenneth E. Jansenc
53659599516SKenneth E. Jansen         allocate( uold(nshg,nsd) )
53759599516SKenneth E. Jansen         allocate( uread(nshg,nsd) )
53859599516SKenneth E. Jansen
53959599516SKenneth E. Jansen         iusiz=nshg*nsd
54059599516SKenneth E. Jansen
541d5d2f64dSCameron Smith         call phio_readdatablock(fhandle,
542bc62cfd4SCameron Smith     &    c_char_'displacement' // char(0),
543bc62cfd4SCameron Smith     &    c_loc(uread), iusiz, dataDbl, iotype)
54459599516SKenneth E. Jansen
54559599516SKenneth E. Jansen         uold(:,1:nsd)=uread(:,1:nsd)
54659599516SKenneth E. Jansen       else
54759599516SKenneth E. Jansen         allocate( uold(nshg,nsd) )
54859599516SKenneth E. Jansen         uold(:,1:nsd) = zero
54959599516SKenneth E. Jansen       endif
55059599516SKenneth E. Jansenc
55159599516SKenneth E. Jansenc.... close c-binary files
55259599516SKenneth E. Jansenc
553d7abaf6cSCameron Smith      call phio_closefile(fhandle)
554266200f9SCameron Smith      iotime = TMRC() - iotime
555266200f9SCameron Smith      if (myrank.eq.master) then
556266200f9SCameron Smith        write(*,*) 'time to read restart (seconds)', iotime
557266200f9SCameron Smith      endif
55859599516SKenneth E. Jansen
55959599516SKenneth E. Jansen      deallocate(xread)
56059599516SKenneth E. Jansen      if ( numpbc > 0 )  then
56159599516SKenneth E. Jansen         deallocate(bcinpread)
56259599516SKenneth E. Jansen         deallocate(ibctmpread)
56359599516SKenneth E. Jansen      endif
56459599516SKenneth E. Jansen      deallocate(iperread)
56559599516SKenneth E. Jansen      if(numpe.gt.1)
56659599516SKenneth E. Jansen     &     deallocate(ilworkread)
56759599516SKenneth E. Jansen      deallocate(nbcread)
56859599516SKenneth E. Jansen
56959599516SKenneth E. Jansen      return
57059599516SKenneth E. Jansen 994  call error ('input   ','opening ', igeomBAK)
57159599516SKenneth E. Jansen 995  call error ('input   ','opening ', igeomBAK)
57259599516SKenneth E. Jansen 997  call error ('input   ','end file', igeomBAK)
57359599516SKenneth E. Jansen 998  call error ('input   ','end file', igeomBAK)
57459599516SKenneth E. Jansen      end
57559599516SKenneth E. Jansenc
57659599516SKenneth E. Jansenc No longer called but kept around in case....
57759599516SKenneth E. Jansenc
57859599516SKenneth E. Jansen      subroutine genpzero(iBC)
57959599516SKenneth E. Jansen
58059599516SKenneth E. Jansen      use pointer_data
58159599516SKenneth E. Jansen      include "common.h"
58259599516SKenneth E. Jansen      integer iBC(nshg)
58359599516SKenneth E. Jansenc
58459599516SKenneth E. Jansenc....  check to see if any of the nodes have a dirichlet pressure
58559599516SKenneth E. Jansenc
58659599516SKenneth E. Jansen      pzero=1
58759599516SKenneth E. Jansen      if (any(btest(iBC,2))) pzero=0
58859599516SKenneth E. Jansen      do iblk = 1, nelblb
58959599516SKenneth E. Jansen         npro = lcblkb(1,iblk+1)-lcblkb(1,iblk)
59059599516SKenneth E. Jansen         do i=1, npro
59159599516SKenneth E. Jansen            iBCB1=miBCB(iblk)%p(i,1)
59259599516SKenneth E. Jansenc
59359599516SKenneth E. Jansenc.... check to see if any of the nodes have a Neumann pressure
59459599516SKenneth E. Jansenc     but not periodic (note that
59559599516SKenneth E. Jansenc
59659599516SKenneth E. Jansen            if(btest(iBCB1,1)) pzero=0
59759599516SKenneth E. Jansen         enddo
59859599516SKenneth E. Jansenc
59959599516SKenneth E. Jansenc.... share results with other processors
60059599516SKenneth E. Jansenc
60159599516SKenneth E. Jansen         pzl=pzero
60259599516SKenneth E. Jansen         if (numpe .gt. 1)
60359599516SKenneth E. Jansen     &        call MPI_ALLREDUCE (pzl, pzero, 1,
60459599516SKenneth E. Jansen     &        MPI_DOUBLE_PRECISION,MPI_MIN, MPI_COMM_WORLD,ierr)
60559599516SKenneth E. Jansen      enddo
60659599516SKenneth E. Jansen      return
60759599516SKenneth E. Jansen      end
608