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