xref: /phasta/M2NFixBnd/src/readnblk.f (revision 595995161822a203c8467e0e4a253d7bd7d6df32)
1*59599516SKenneth E. Jansenc  readnblk.f (pronounce "Reed and Block Dot Eff") contains:
2*59599516SKenneth E. Jansenc
3*59599516SKenneth E. Jansenc    module readarrays ("Red Arrays") -- contains the arrays that
4*59599516SKenneth E. Jansenc     are read in from binary files but not immediately blocked
5*59599516SKenneth E. Jansenc     through pointers.
6*59599516SKenneth E. Jansenc
7*59599516SKenneth E. Jansenc    subroutine readnblk ("Reed and Block") -- allocates space for
8*59599516SKenneth E. Jansenc     and reads data to be contained in module readarrays.  Reads
9*59599516SKenneth E. Jansenc     all remaining data and blocks them with pointers.
10*59599516SKenneth E. Jansenc
11*59599516SKenneth E. Jansen
12*59599516SKenneth E. Jansen
13*59599516SKenneth E. Jansen      module readarrays
14*59599516SKenneth E. Jansen
15*59599516SKenneth E. Jansen      real*8, allocatable :: point2x(:,:)
16*59599516SKenneth E. Jansen      real*8, allocatable :: qold(:,:)
17*59599516SKenneth E. Jansen      real*8, allocatable :: dwal(:)
18*59599516SKenneth E. Jansen      real*8, allocatable :: errors(:,:)
19*59599516SKenneth E. Jansen      real*8, allocatable :: ybar(:,:)
20*59599516SKenneth E. Jansen      real*8, allocatable :: yphbar(:,:,:)
21*59599516SKenneth E. Jansen      real*8, allocatable :: vort(:,:)
22*59599516SKenneth E. Jansen      real*8, allocatable :: uold(:,:)
23*59599516SKenneth E. Jansen      real*8, allocatable :: acold(:,:)
24*59599516SKenneth E. Jansen      integer, allocatable :: iBCtmp(:)
25*59599516SKenneth E. Jansen      real*8, allocatable :: BCinp(:,:)
26*59599516SKenneth E. Jansen
27*59599516SKenneth E. Jansen      integer, allocatable :: point2ilwork(:)
28*59599516SKenneth E. Jansen      integer, allocatable :: nBC(:)
29*59599516SKenneth E. Jansen      integer, allocatable :: point2iper(:)
30*59599516SKenneth E. Jansen      integer, allocatable :: point2ifath(:)
31*59599516SKenneth E. Jansen      integer, allocatable :: point2nsons(:)
32*59599516SKenneth E. Jansen
33*59599516SKenneth E. Jansen      end module
34*59599516SKenneth E. Jansen
35*59599516SKenneth E. Jansen
36*59599516SKenneth E. Jansen
37*59599516SKenneth E. Jansen      subroutine readnblk
38*59599516SKenneth E. Jansenc
39*59599516SKenneth E. Jansen      use readarrays
40*59599516SKenneth E. Jansen      include "commonM2NFixBnd.h"
41*59599516SKenneth E. Jansen      include "mpif.h"
42*59599516SKenneth E. Jansenc
43*59599516SKenneth E. Jansen      real*8, allocatable :: xread(:,:), qread(:,:), qread1(:)
44*59599516SKenneth E. Jansen      real*8, allocatable :: uread(:,:), acread(:,:)
45*59599516SKenneth E. Jansen      real*8, allocatable :: BCinpread(:,:)
46*59599516SKenneth E. Jansen      real*8 globmax,globmin
47*59599516SKenneth E. Jansen      integer, allocatable :: iperread(:), iBCtmpread(:)
48*59599516SKenneth E. Jansen      integer, allocatable :: ilworkread(:), nBCread(:)
49*59599516SKenneth E. Jansen      character*10 cname2
50*59599516SKenneth E. Jansen      character*30 fmt1
51*59599516SKenneth E. Jansen      character*255 fname1,fnamer,fnamelr
52*59599516SKenneth E. Jansen      character*255 warning
53*59599516SKenneth E. Jansen
54*59599516SKenneth E. Jansen      integer :: descriptor, color, nfiles, nfields
55*59599516SKenneth E. Jansen      integer ::  numparts, nppf
56*59599516SKenneth E. Jansen      character*255 fname2, temp2
57*59599516SKenneth E. Jansen      character*64 temp1
58*59599516SKenneth E. Jansen      integer :: igeom, ibndc, irstin, ierr
59*59599516SKenneth E. Jansen      integer :: ndof, ndoferrors, ndofybar
60*59599516SKenneth E. Jansen      integer :: itmp, itmp2
61*59599516SKenneth E. Jansen      integer :: irstart, irstartmap, iybar
62*59599516SKenneth E. Jansen      integer :: ierror, numphavg, ivort, idwal, idebug, iphavg
63*59599516SKenneth E. Jansen
64*59599516SKenneth E. Jansen      integer intfromfile(50) ! integers read from headers
65*59599516SKenneth E. Jansen      logical exinput
66*59599516SKenneth E. Jansenc
67*59599516SKenneth E. Jansenc
68*59599516SKenneth E. Jansenc.... determine the step number to start with
69*59599516SKenneth E. Jansenc
70*59599516SKenneth E. Jansen!      open(unit=72,file='numstart.dat',status='old')
71*59599516SKenneth E. Jansen!      read(72,*) irstart
72*59599516SKenneth E. Jansen!      close(72)
73*59599516SKenneth E. Jansen
74*59599516SKenneth E. Jansen      if(myrank == 0) then
75*59599516SKenneth E. Jansen        fnamer='M2N_input.dat'
76*59599516SKenneth E. Jansen        fnamer = trim(fnamer) // char(0)
77*59599516SKenneth E. Jansen        inquire(file=fnamer,exist=exinput)
78*59599516SKenneth E. Jansen        if(exinput) then
79*59599516SKenneth E. Jansen          open(unit=72,file=fnamer,status='old')
80*59599516SKenneth E. Jansen          read(72,*) irstart
81*59599516SKenneth E. Jansen          read(72,*) irstartmap
82*59599516SKenneth E. Jansen          read(72,*) iybar
83*59599516SKenneth E. Jansen          read(72,*) ierror
84*59599516SKenneth E. Jansen          read(72,*) numphavg
85*59599516SKenneth E. Jansen          read(72,*) ivort
86*59599516SKenneth E. Jansen          read(72,*) idwal
87*59599516SKenneth E. Jansen          read(72,*) idebug
88*59599516SKenneth E. Jansen          close(72)
89*59599516SKenneth E. Jansen        else
90*59599516SKenneth E. Jansen           write(*,*) 'ERROR: Input file ',
91*59599516SKenneth E. Jansen     &                 trim(fnamer),' does not exist!'
92*59599516SKenneth E. Jansen        endif
93*59599516SKenneth E. Jansen      endif
94*59599516SKenneth E. Jansen      call mpi_barrier(mpi_comm_world, ierr)
95*59599516SKenneth E. Jansen      call mpi_bcast(exinput,1,MPI_LOGICAL,0,mpi_comm_world,ierr)
96*59599516SKenneth E. Jansen      if(.not. exinput) then ! M2NFixBnd_input.dat does not exist. Quit
97*59599516SKenneth E. Jansen!        call mpi_abort(mpi_comm_world,911,ierr)
98*59599516SKenneth E. Jansen!        call mpi_finalize(ierr)
99*59599516SKenneth E. Jansen        return
100*59599516SKenneth E. Jansen      else ! broadcast the information read by rank 0
101*59599516SKenneth E. Jansen        call mpi_bcast(irstart,1,MPI_INTEGER,0,mpi_comm_world,ierr)
102*59599516SKenneth E. Jansen        call mpi_bcast(iybar,1,MPI_INTEGER,0,mpi_comm_world,ierr)
103*59599516SKenneth E. Jansen        call mpi_bcast(ierror,1,MPI_INTEGER,0,mpi_comm_world,ierr)
104*59599516SKenneth E. Jansen        call mpi_bcast(numphavg,1,MPI_INTEGER,0,mpi_comm_world,ierr)
105*59599516SKenneth E. Jansen        call mpi_bcast(ivort,1,MPI_INTEGER,0,mpi_comm_world,ierr)
106*59599516SKenneth E. Jansen        call mpi_bcast(idwal,1,MPI_INTEGER,0,mpi_comm_world,ierr)
107*59599516SKenneth E. Jansen      endif
108*59599516SKenneth E. Jansen
109*59599516SKenneth E. Jansen      if(myrank == 0 ) then
110*59599516SKenneth E. Jansen      ! Print some info
111*59599516SKenneth E. Jansen        write(*,*) 'The solution field is reduced by default'
112*59599516SKenneth E. Jansen
113*59599516SKenneth E. Jansen        if(iybar .gt. 0) then
114*59599516SKenneth E. Jansen          write(*,*) 'The ybar field will also be reduced'
115*59599516SKenneth E. Jansen          iybar = 1 ! security
116*59599516SKenneth E. Jansen        else
117*59599516SKenneth E. Jansen          write(*,*) 'The ybar field will NOT be reduced'
118*59599516SKenneth E. Jansen        endif
119*59599516SKenneth E. Jansen
120*59599516SKenneth E. Jansen        if(ierror .gt. 0) then
121*59599516SKenneth E. Jansen          write(*,*) 'The error field will also be reduced'
122*59599516SKenneth E. Jansen          ierror = 1 ! security
123*59599516SKenneth E. Jansen        else
124*59599516SKenneth E. Jansen          write(*,*) 'The error field will NOT be reduced'
125*59599516SKenneth E. Jansen        endif
126*59599516SKenneth E. Jansen
127*59599516SKenneth E. Jansen        if(numphavg .gt. 0) then
128*59599516SKenneth E. Jansen          write(*,*) 'The phase average fields (',numphavg,
129*59599516SKenneth E. Jansen     &                ') will also be reduced'
130*59599516SKenneth E. Jansen        else
131*59599516SKenneth E. Jansen          write(*,*) 'The phase average fields will NOT be reduced'
132*59599516SKenneth E. Jansen        endif
133*59599516SKenneth E. Jansen
134*59599516SKenneth E. Jansen        if(ivort .gt. 0) then
135*59599516SKenneth E. Jansen          write(*,*) 'The vorticity field will also be reduced'
136*59599516SKenneth E. Jansen          ivort = 1 ! security
137*59599516SKenneth E. Jansen        else
138*59599516SKenneth E. Jansen          write(*,*) 'The vorticity field will NOT be reduced'
139*59599516SKenneth E. Jansen        endif
140*59599516SKenneth E. Jansen
141*59599516SKenneth E. Jansen        if(idwal .gt. 0) then
142*59599516SKenneth E. Jansen          write(*,*) 'The dwal field will also be reduced'
143*59599516SKenneth E. Jansen          idwal = 1 ! security
144*59599516SKenneth E. Jansen        else
145*59599516SKenneth E. Jansen          write(*,*) 'The dwal field will NOT be reduced'
146*59599516SKenneth E. Jansen        endif
147*59599516SKenneth E. Jansen        write(*,*) ''
148*59599516SKenneth E. Jansen      endif
149*59599516SKenneth E. Jansen
150*59599516SKenneth E. Jansen
151*59599516SKenneth E. Jansen      call mpi_barrier(mpi_comm_world, ierr)
152*59599516SKenneth E. Jansen      if(myrank == 0) then
153*59599516SKenneth E. Jansen        write(*,*) 'Reading the geombcRed-dat files'
154*59599516SKenneth E. Jansen      endif
155*59599516SKenneth E. Jansen
156*59599516SKenneth E. Jansen      lstep=irstart ! in case restart files have no fields
157*59599516SKenneth E. Jansen
158*59599516SKenneth E. Jansen      nfiles = nsynciofiles
159*59599516SKenneth E. Jansen      numparts = numpe !This is the common settings. Beware if you try to compute several parts per process
160*59599516SKenneth E. Jansenc
161*59599516SKenneth E. Jansenc.... input the geometry parameters
162*59599516SKenneth E. Jansenc
163*59599516SKenneth E. Jansen
164*59599516SKenneth E. Jansen      color = int(myrank/(numparts/nfiles)) !Should call the color routine in SyncIO here
165*59599516SKenneth E. Jansen      itmp2 = int(log10(float(color+1)))+1
166*59599516SKenneth E. Jansen      write (temp2,"('(''geombcRed-dat.'',i',i1,')')") itmp2
167*59599516SKenneth E. Jansen      write (fnamer,temp2) (color+1)
168*59599516SKenneth E. Jansen      fnamer = trim(fnamer)//char(0)
169*59599516SKenneth E. Jansen
170*59599516SKenneth E. Jansen      ieleven=11
171*59599516SKenneth E. Jansen      ione=1
172*59599516SKenneth E. Jansen
173*59599516SKenneth E. Jansen      itmp = int(log10(float(myrank+1)))+1
174*59599516SKenneth E. Jansen
175*59599516SKenneth E. Jansen      call queryphmpiio(fnamer, nfields, nppf);
176*59599516SKenneth E. Jansen      if (myrank == 0) then
177*59599516SKenneth E. Jansen        write(*,*) 'Number of fields in geombcRed-dat: ',nfields
178*59599516SKenneth E. Jansen        write(*,*) 'Number of parts per file geombcRed-dat: ',nppf
179*59599516SKenneth E. Jansen      endif
180*59599516SKenneth E. Jansen      call initphmpiio( nfields, nppf, nfiles, igeom,
181*59599516SKenneth E. Jansen     & 'read'//char(0))
182*59599516SKenneth E. Jansen      call openfile( fnamer, 'read'//char(0), igeom )
183*59599516SKenneth E. Jansen
184*59599516SKenneth E. Jansen      write (temp1,"('(''number of nodes@'',i',i1,',A1)')") itmp
185*59599516SKenneth E. Jansen      write (fname2,temp1) (myrank+1),'?'
186*59599516SKenneth E. Jansen      call readheader(igeom,fname2//char(0),numnp,ione,
187*59599516SKenneth E. Jansen     & 'integer'//char(0), iotype)
188*59599516SKenneth E. Jansen
189*59599516SKenneth E. Jansen      write (temp1,"('(''number of modes@'',i',i1,',A1)')") itmp
190*59599516SKenneth E. Jansen      write (fname2,temp1) (myrank+1),'?'
191*59599516SKenneth E. Jansen      call readheader(igeom,fname2//char(0),nshg,ione,
192*59599516SKenneth E. Jansen     & 'integer'//char(0), iotype)
193*59599516SKenneth E. Jansen
194*59599516SKenneth E. Jansen      write (temp1,"('(''number of interior elements@'',i',i1,',A1)')")
195*59599516SKenneth E. Jansen     &       itmp
196*59599516SKenneth E. Jansen      write (fname2,temp1) (myrank+1),'?'
197*59599516SKenneth E. Jansen      call readheader(igeom,fname2//char(0),numel,ione,
198*59599516SKenneth E. Jansen     & 'integer'//char(0), iotype)
199*59599516SKenneth E. Jansen
200*59599516SKenneth E. Jansen      write (temp1,"('(''number of boundary elements@'',i',i1,',A1)')")
201*59599516SKenneth E. Jansen     &       itmp
202*59599516SKenneth E. Jansen      write (fname2,temp1) (myrank+1),'?'
203*59599516SKenneth E. Jansen      call readheader(igeom,fname2//char(0),numelb,ione,
204*59599516SKenneth E. Jansen     & 'integer'//char(0),iotype)
205*59599516SKenneth E. Jansen
206*59599516SKenneth E. Jansen      write (temp1,
207*59599516SKenneth E. Jansen     & "('(''maximum number of element nodes@'',i',i1,',A1)')") itmp
208*59599516SKenneth E. Jansen      write (fname2,temp1) (myrank+1),'?'
209*59599516SKenneth E. Jansen      call readheader(igeom,fname2//char(0),nen,ione,
210*59599516SKenneth E. Jansen     &'integer'//char(0),iotype)
211*59599516SKenneth E. Jansen
212*59599516SKenneth E. Jansen      write (temp1,"('(''number of interior tpblocks@'',i',i1,',A1)')")
213*59599516SKenneth E. Jansen     &       itmp
214*59599516SKenneth E. Jansen      write (fname2,temp1) (myrank+1),'?'
215*59599516SKenneth E. Jansen      call readheader(igeom,fname2//char(0),nelblk,ione,
216*59599516SKenneth E. Jansen     & 'integer'//char(0) ,iotype)
217*59599516SKenneth E. Jansen
218*59599516SKenneth E. Jansen      write (temp1,"('(''number of boundary tpblocks@'',i',i1,',A1)')")
219*59599516SKenneth E. Jansen     &       itmp
220*59599516SKenneth E. Jansen      write (fname2,temp1) (myrank+1),'?'
221*59599516SKenneth E. Jansen      call readheader(igeom,fname2//char(0),nelblb,ione,
222*59599516SKenneth E. Jansen     & 'integer'//char(0), iotype)
223*59599516SKenneth E. Jansen
224*59599516SKenneth E. Jansen      write (temp1,
225*59599516SKenneth E. Jansen     & "('(''number of nodes with Dirichlet BCs@'',i',i1,',A1)')") itmp
226*59599516SKenneth E. Jansen      write (fname2,temp1) (myrank+1),'?'
227*59599516SKenneth E. Jansen      call readheader(igeom,fname2//char(0),numpbc,ione,
228*59599516SKenneth E. Jansen     & 'integer'//char(0),iotype)
229*59599516SKenneth E. Jansen
230*59599516SKenneth E. Jansen      write (temp1,"('(''number of shape functions@'',i',i1,',A1)')")
231*59599516SKenneth E. Jansen     &       itmp
232*59599516SKenneth E. Jansen      write (fname2,temp1) (myrank+1),'?'
233*59599516SKenneth E. Jansen      call readheader(igeom,fname2//char(0),ntopsh,ione,
234*59599516SKenneth E. Jansen     & 'integer'//char(0),iotype)
235*59599516SKenneth E. Jansen
236*59599516SKenneth E. Jansen      call closefile( igeom, "read"//char(0) )
237*59599516SKenneth E. Jansen      call finalizephmpiio( igeom )
238*59599516SKenneth E. Jansen
239*59599516SKenneth E. Jansenc
240*59599516SKenneth E. Jansenc.... calculate the maximum number of boundary element nodes
241*59599516SKenneth E. Jansenc
242*59599516SKenneth E. Jansen      nenb = 3 !was initialized to 0 but
243*59599516SKenneth E. Jansen      do i = 1, melCat !melCat is 0 here
244*59599516SKenneth E. Jansen         if (nen .eq. nenCat(i,nsd)) nenb = max(nenCat(i,nsd-1), nenb)
245*59599516SKenneth E. Jansen      enddo
246*59599516SKenneth E. Jansenc
247*59599516SKenneth E. Jansen      if (myrank == master) then
248*59599516SKenneth E. Jansen         if (nenb .eq. 0) call error ('input   '//char(0),
249*59599516SKenneth E. Jansen     &                                'nen     '//char(0),nen)
250*59599516SKenneth E. Jansen      endif
251*59599516SKenneth E. Jansenc
252*59599516SKenneth E. Jansenc.... setup some useful constants
253*59599516SKenneth E. Jansenc
254*59599516SKenneth E. Jansen      I3nsd  = nsd / 3          ! nsd=3 integer flag
255*59599516SKenneth E. Jansen      E3nsd  = float(I3nsd)     ! nsd=3 real    flag
256*59599516SKenneth E. Jansenc
257*59599516SKenneth E. Jansen      if(matflg(1,1).lt.0) then
258*59599516SKenneth E. Jansen         nflow = nsd + 1
259*59599516SKenneth E. Jansen      else
260*59599516SKenneth E. Jansen         nflow = nsd + 2
261*59599516SKenneth E. Jansen      endif
262*59599516SKenneth E. Jansen      ndof   = nsd + 2
263*59599516SKenneth E. Jansen      nsclr=impl(1)/100
264*59599516SKenneth E. Jansen      ndof=ndof+nsclr           ! number of sclr transport equations to solve
265*59599516SKenneth E. Jansen
266*59599516SKenneth E. Jansen      ndofBC = ndof + I3nsd     ! dimension of BC array
267*59599516SKenneth E. Jansen      ndiBCB = 2                ! dimension of iBCB array
268*59599516SKenneth E. Jansen      ndBCB  = ndof + 1         ! dimension of BCB array
269*59599516SKenneth E. Jansenc
270*59599516SKenneth E. Jansen      nsymdf = (ndof*(ndof + 1)) / 2 ! symm. d.o.f.'s
271*59599516SKenneth E. Jansen
272*59599516SKenneth E. Jansen
273*59599516SKenneth E. Jansen
274*59599516SKenneth E. Jansenc
275*59599516SKenneth E. Jansenc.... Read restart files
276*59599516SKenneth E. Jansenc
277*59599516SKenneth E. Jansenc.... read the header and check it against the run data
278*59599516SKenneth E. Jansenc
279*59599516SKenneth E. Jansen      call mpi_barrier(mpi_comm_world, ierr)
280*59599516SKenneth E. Jansen      if(myrank == 0) then
281*59599516SKenneth E. Jansen        write(*,*) 'Reading the RestartRedTmp-dat files'
282*59599516SKenneth E. Jansen      endif
283*59599516SKenneth E. Jansen
284*59599516SKenneth E. Jansen!     Beware in what follows! nshg2 read from the header of the solution,
285*59599516SKenneth E. Jansen!     dwal, errors and ybar can be different from nshg read from the geombc files
286*59599516SKenneth E. Jansen!     This is due to the fact that nshg2 represents the largest vID referenced
287*59599516SKenneth E. Jansen!     by the mapping and not the largest vID present in the mesh.
288*59599516SKenneth E. Jansen!     qold, dwal, errors, ybar must be allocated to nshg and initialized
289*59599516SKenneth E. Jansen!     to a large negative value. The unmapped vertices will be updated
290*59599516SKenneth E. Jansen!     through communication with ilwork.
291*59599516SKenneth E. Jansen
292*59599516SKenneth E. Jansen      itmp=1
293*59599516SKenneth E. Jansen      if (irstart .gt. 0) itmp = int(log10(float(irstart+1)))+1
294*59599516SKenneth E. Jansen
295*59599516SKenneth E. Jansen      write (fmt1,"('(''restartRedTmp-dat.'',i',i1,',1x)')") itmp
296*59599516SKenneth E. Jansen      write (fnamer,fmt1) irstart
297*59599516SKenneth E. Jansen      fnamer = trim(fnamer) // cname2(color+1)
298*59599516SKenneth E. Jansen
299*59599516SKenneth E. Jansen      call queryphmpiio(fnamer//char(0), nfields, nppf);
300*59599516SKenneth E. Jansen      if (myrank == 0) then
301*59599516SKenneth E. Jansen        write(*,*) 'Number of fields in restartRedTmp-dat: ',nfields
302*59599516SKenneth E. Jansen        write(*,*) 'Number of parts per file restartRedTmp-dat: ',nppf
303*59599516SKenneth E. Jansen      endif
304*59599516SKenneth E. Jansen
305*59599516SKenneth E. Jansen      call initphmpiio(nfields,nppf,nfiles,descriptor,
306*59599516SKenneth E. Jansen     & 'read'//char(0))
307*59599516SKenneth E. Jansen      call openfile( fnamer//char(0) ,
308*59599516SKenneth E. Jansen     & 'read'//char(0), descriptor )
309*59599516SKenneth E. Jansen
310*59599516SKenneth E. Jansenc
311*59599516SKenneth E. Jansenc     Read the solution
312*59599516SKenneth E. Jansenc
313*59599516SKenneth E. Jansen      ithree=3
314*59599516SKenneth E. Jansen      itmp = int(log10(float(myrank+1)))+1
315*59599516SKenneth E. Jansen      write (temp1,"('(''solution@'',i',i1,',A1)')") itmp
316*59599516SKenneth E. Jansen      write (fname1,temp1) (myrank+1),'?'
317*59599516SKenneth E. Jansen      fname1 = trim(fname1)
318*59599516SKenneth E. Jansen
319*59599516SKenneth E. Jansen      intfromfile=0
320*59599516SKenneth E. Jansen      call readheader(descriptor,fname1//char(0) ,intfromfile,
321*59599516SKenneth E. Jansen     & ithree,'integer'//char(0), iotype)
322*59599516SKenneth E. Jansen
323*59599516SKenneth E. Jansenc
324*59599516SKenneth E. Jansenc.... read the values of primitive variables into q
325*59599516SKenneth E. Jansenc
326*59599516SKenneth E. Jansen      if(intfromfile(1).ne.0) then
327*59599516SKenneth E. Jansen         nshg2=intfromfile(1)
328*59599516SKenneth E. Jansen         ndof2=intfromfile(2)
329*59599516SKenneth E. Jansen         ndof=ndof2 !This must be the same anyway
330*59599516SKenneth E. Jansen         allocate( qold(nshg,ndof) )
331*59599516SKenneth E. Jansen         qold(:,:) = -9.87654321e32
332*59599516SKenneth E. Jansen         lstep=intfromfile(3)
333*59599516SKenneth E. Jansen         allocate( qread(nshg2,ndof2) )
334*59599516SKenneth E. Jansen
335*59599516SKenneth E. Jansen         if (nshg2 .ne. nshg) then
336*59599516SKenneth E. Jansen           write(*,*) 'nshg from geombc and nshg2 from restart differ'
337*59599516SKenneth E. Jansen     &                //' on rank', myrank, ' :',nshg,nshg2,
338*59599516SKenneth E. Jansen     &               ' - Probably mixing phasta files'
339*59599516SKenneth E. Jansen         endif
340*59599516SKenneth E. Jansen         call mpi_barrier(mpi_comm_world, ierr)
341*59599516SKenneth E. Jansen
342*59599516SKenneth E. Jansen         iqsiz=nshg2*ndof2
343*59599516SKenneth E. Jansen
344*59599516SKenneth E. Jansen         call readdatablock(descriptor,fname1//char(0),qread,iqsiz,
345*59599516SKenneth E. Jansen     &                         'double'//char(0),iotype)
346*59599516SKenneth E. Jansen         qold(1:nshg2,1:ndof2)=qread(1:nshg2,1:ndof2)
347*59599516SKenneth E. Jansen         deallocate(qread)
348*59599516SKenneth E. Jansen      else
349*59599516SKenneth E. Jansen         if (myrank.eq.master) then
350*59599516SKenneth E. Jansen            if (matflg(1,1).eq.0) then ! compressible
351*59599516SKenneth E. Jansen               warning='Solution is set to zero (with p and T to one)'
352*59599516SKenneth E. Jansen            else
353*59599516SKenneth E. Jansen               warning='Solution is set to zero'
354*59599516SKenneth E. Jansen            endif
355*59599516SKenneth E. Jansen            write(*,*) warning// char(0)
356*59599516SKenneth E. Jansen         endif
357*59599516SKenneth E. Jansen         qold=zero
358*59599516SKenneth E. Jansen         if (matflg(1,1).eq.0) then ! compressible
359*59599516SKenneth E. Jansen            qold(:,1)=one ! avoid zero pressure
360*59599516SKenneth E. Jansen            qold(:,nflow)=one ! avoid zero temperature
361*59599516SKenneth E. Jansen         endif
362*59599516SKenneth E. Jansen      endif
363*59599516SKenneth E. Jansen
364*59599516SKenneth E. Jansenc
365*59599516SKenneth E. Jansenc    Read the ybar
366*59599516SKenneth E. Jansenc
367*59599516SKenneth E. Jansen      ndofybar = 0
368*59599516SKenneth E. Jansen      if(iybar==1) then
369*59599516SKenneth E. Jansen        itmp = int(log10(float(myrank+1)))+1
370*59599516SKenneth E. Jansen        write (temp1,"('(''ybar@'',i',i1,',a1)')") itmp
371*59599516SKenneth E. Jansen        write (fname1,temp1) (myrank+1),'?'
372*59599516SKenneth E. Jansen        fname1 = trim(fname1)
373*59599516SKenneth E. Jansen        intfromfile=0
374*59599516SKenneth E. Jansen        call readheader(descriptor,fname1//char(0),intfromfile,
375*59599516SKenneth E. Jansen     &     ithree,'integer'//char(0),iotype)
376*59599516SKenneth E. Jansen        nshg2=intfromfile(1)
377*59599516SKenneth E. Jansen        ndofybar=intfromfile(2)
378*59599516SKenneth E. Jansen        !lstep=intfromfile(3)
379*59599516SKenneth E. Jansen        if(ndofybar .ne. 0) then
380*59599516SKenneth E. Jansen          allocate( ybar(nshg,ndofybar) )
381*59599516SKenneth E. Jansen          ybar(:,:) = -9.87654321e32
382*59599516SKenneth E. Jansen          allocate( qread(nshg2,ndofybar) )
383*59599516SKenneth E. Jansen
384*59599516SKenneth E. Jansen          iqsiz=nshg2*ndofybar
385*59599516SKenneth E. Jansen          call readdatablock(descriptor,fname1//char(0) ,qread,iqsiz,
386*59599516SKenneth E. Jansen     &                   'double'//char(0),iotype)
387*59599516SKenneth E. Jansen          ybar(1:nshg2,1:ndofybar)=qread(1:nshg2,1:ndofybar)
388*59599516SKenneth E. Jansen          deallocate(qread)
389*59599516SKenneth E. Jansen        else
390*59599516SKenneth E. Jansen          write(*,*) 'WARNING: ybar is missing in the restart files'
391*59599516SKenneth E. Jansen          iybar = 0
392*59599516SKenneth E. Jansen        endif
393*59599516SKenneth E. Jansen      endif
394*59599516SKenneth E. Jansen
395*59599516SKenneth E. Jansenc
396*59599516SKenneth E. Jansenc    Read the errors
397*59599516SKenneth E. Jansenc
398*59599516SKenneth E. Jansen      ndoferrors=0
399*59599516SKenneth E. Jansen      if(ierror == 1) then
400*59599516SKenneth E. Jansen        write (temp1,"('(''errors@'',i',i1,',a1)')")
401*59599516SKenneth E. Jansen     &       itmp
402*59599516SKenneth E. Jansen        write (fname1,temp1) (myrank+1),'?'
403*59599516SKenneth E. Jansen        fname1 = trim(fname1)
404*59599516SKenneth E. Jansen        intfromfile=0
405*59599516SKenneth E. Jansen        call readheader(descriptor,fname1//char(0),intfromfile,
406*59599516SKenneth E. Jansen     &     ithree,'integer'//char(0),iotype)
407*59599516SKenneth E. Jansen        nshg2=intfromfile(1)
408*59599516SKenneth E. Jansen        ndoferrors=intfromfile(2)
409*59599516SKenneth E. Jansen        !lstep=intfromfile(3)
410*59599516SKenneth E. Jansen        if(ndoferrors .ne. 0) then
411*59599516SKenneth E. Jansen          allocate( errors(nshg,ndoferrors) )
412*59599516SKenneth E. Jansen          errors(:,:) = -9.87654321e32
413*59599516SKenneth E. Jansen          allocate( qread(nshg2,ndoferrors) )
414*59599516SKenneth E. Jansen          iqsiz=nshg2*ndoferrors
415*59599516SKenneth E. Jansen          call readdatablock(descriptor,fname1//char(0),qread,iqsiz,
416*59599516SKenneth E. Jansen     &                   'double'//char(0),iotype)
417*59599516SKenneth E. Jansen          errors(1:nshg2,1:ndoferrors)=qread(1:nshg2,1:ndoferrors)
418*59599516SKenneth E. Jansen          deallocate(qread)
419*59599516SKenneth E. Jansen        else
420*59599516SKenneth E. Jansen          if(myrank==0) then
421*59599516SKenneth E. Jansen            write(*,*) 'WARNING: errors is missing in the restart files'
422*59599516SKenneth E. Jansen          endif
423*59599516SKenneth E. Jansen          ierror = 0
424*59599516SKenneth E. Jansen        endif
425*59599516SKenneth E. Jansen      endif
426*59599516SKenneth E. Jansen
427*59599516SKenneth E. Jansen
428*59599516SKenneth E. Jansen!
429*59599516SKenneth E. Jansen!     Read the phase_average fields
430*59599516SKenneth E. Jansen!
431*59599516SKenneth E. Jansen      ndofyphbar=0
432*59599516SKenneth E. Jansen      if(numphavg .gt. 0) then
433*59599516SKenneth E. Jansen        do iphavg = 1,numphavg
434*59599516SKenneth E. Jansen          itmp = int(log10(float(myrank+1)))+1
435*59599516SKenneth E. Jansen          itmp2 = int(log10(float(iphavg)))+1
436*59599516SKenneth E. Jansen          write (temp1,
437*59599516SKenneth E. Jansen     &            "('(''phase_average'',i',i1,',''@'',i',i1,',A1)')")
438*59599516SKenneth E. Jansen     &             itmp2, itmp
439*59599516SKenneth E. Jansen          write (fname1,temp1) iphavg,(myrank+1),'?'
440*59599516SKenneth E. Jansen          fname1 = trim(fname1)
441*59599516SKenneth E. Jansen          intfromfile=0
442*59599516SKenneth E. Jansen          call readheader(descriptor,fname1//char(0),intfromfile,
443*59599516SKenneth E. Jansen     &                  ithree,'integer'//char(0),iotype)
444*59599516SKenneth E. Jansen
445*59599516SKenneth E. Jansen          nshg2=intfromfile(1)
446*59599516SKenneth E. Jansen          ndofyphbar=intfromfile(2)
447*59599516SKenneth E. Jansen          !lstep=intfromfile(3)
448*59599516SKenneth E. Jansen
449*59599516SKenneth E. Jansen          if(ndofyphbar.ne.0) then
450*59599516SKenneth E. Jansen            ! Allocate some memory for the first ts only
451*59599516SKenneth E. Jansen            if(iphavg==1) then
452*59599516SKenneth E. Jansen              allocate( yphbar(nshg,ndofyphbar,numphavg) )
453*59599516SKenneth E. Jansen              yphbar(:,:,:) = -9.87654321e32
454*59599516SKenneth E. Jansen            endif
455*59599516SKenneth E. Jansen
456*59599516SKenneth E. Jansen            allocate( qread(nshg2,ndofyphbar) )
457*59599516SKenneth E. Jansen            iqsiz = nshg2*ndofyphbar
458*59599516SKenneth E. Jansen            call readdatablock(descriptor,fname1//char(0),qread,iqsiz,
459*59599516SKenneth E. Jansen     &                     'double'//char(0),iotype)
460*59599516SKenneth E. Jansen            yphbar(1:nshg2,1:ndofyphbar,iphavg) =
461*59599516SKenneth E. Jansen     &                         qread(1:nshg2,1:ndofyphbar)
462*59599516SKenneth E. Jansen            deallocate(qread)
463*59599516SKenneth E. Jansen          else
464*59599516SKenneth E. Jansen            if(myrank==0) then
465*59599516SKenneth E. Jansen              write(*,*)'WARNING: phase_average is missing '//
466*59599516SKenneth E. Jansen     &            'in the restart files'
467*59599516SKenneth E. Jansen            endif
468*59599516SKenneth E. Jansen            numphavg = 0
469*59599516SKenneth E. Jansen            if(iphavg > 1) then
470*59599516SKenneth E. Jansen              deallocate(yphbar)
471*59599516SKenneth E. Jansen            endif
472*59599516SKenneth E. Jansen            exit
473*59599516SKenneth E. Jansen          endif
474*59599516SKenneth E. Jansen        enddo
475*59599516SKenneth E. Jansen      endif
476*59599516SKenneth E. Jansen
477*59599516SKenneth E. Jansen!
478*59599516SKenneth E. Jansen! follow the usual convention for loading the vorticity field
479*59599516SKenneth E. Jansen!
480*59599516SKenneth E. Jansen      ndofvort=0
481*59599516SKenneth E. Jansen      if(ivort == 1) then
482*59599516SKenneth E. Jansen        itmp = int(log10(float(myrank+1)))+1
483*59599516SKenneth E. Jansen        write (temp1,"('(''vorticity@'',i',i1,',a1)')") itmp
484*59599516SKenneth E. Jansen        write (fname1,temp1) (myrank+1),'?'
485*59599516SKenneth E. Jansen        fname1 = trim(fname1)
486*59599516SKenneth E. Jansen        intfromfile=0
487*59599516SKenneth E. Jansen        call readheader(descriptor,fname1//char(0),intfromfile,
488*59599516SKenneth E. Jansen     &     ithree,'integer'//char(0),iotype)
489*59599516SKenneth E. Jansen
490*59599516SKenneth E. Jansen        nshg2=intfromfile(1)
491*59599516SKenneth E. Jansen        ndofvort=intfromfile(2)
492*59599516SKenneth E. Jansen        !lstep=intfromfile(3)
493*59599516SKenneth E. Jansen
494*59599516SKenneth E. Jansen        if(ndofvort .ne. 0) then
495*59599516SKenneth E. Jansen          allocate(vort(nshg,ndofvort))
496*59599516SKenneth E. Jansen          vort(:,:) = -9.87654321e32
497*59599516SKenneth E. Jansen          allocate(qread(nshg2,ndofvort))
498*59599516SKenneth E. Jansen          iqsiz = nshg2*ndofvort
499*59599516SKenneth E. Jansen          call readdatablock(descriptor,fname1//char(0),qread,iqsiz,
500*59599516SKenneth E. Jansen     &                   'double'//char(0),iotype)
501*59599516SKenneth E. Jansen          vort(1:nshg2,1:ndofvort)=qread(1:nshg2,1:ndofvort)
502*59599516SKenneth E. Jansen          deallocate(qread)
503*59599516SKenneth E. Jansen        else
504*59599516SKenneth E. Jansen          if(myrank==0) then
505*59599516SKenneth E. Jansen            write(*,*) 'WARNING: vorticity is missing '//
506*59599516SKenneth E. Jansen     &                 'in the restart files'
507*59599516SKenneth E. Jansen          endif
508*59599516SKenneth E. Jansen          ivort = 0
509*59599516SKenneth E. Jansen        endif
510*59599516SKenneth E. Jansen      endif
511*59599516SKenneth E. Jansen
512*59599516SKenneth E. Jansenc
513*59599516SKenneth E. Jansenc     Read the dwal
514*59599516SKenneth E. Jansenc
515*59599516SKenneth E. Jansen      ndofdwal=0
516*59599516SKenneth E. Jansen      if(idwal==1) then
517*59599516SKenneth E. Jansen        write (temp1,"('(''dwal@'',i',i1,',a1)')") itmp
518*59599516SKenneth E. Jansen        write (fname1,temp1) (myrank+1),'?'
519*59599516SKenneth E. Jansen        fname1 = trim(fname1)
520*59599516SKenneth E. Jansen        intfromfile=0
521*59599516SKenneth E. Jansen        call readheader(descriptor,fname1//char(0),intfromfile,
522*59599516SKenneth E. Jansen     &     ithree,'integer'//char(0),iotype)
523*59599516SKenneth E. Jansen
524*59599516SKenneth E. Jansen        nshg2=intfromfile(1)
525*59599516SKenneth E. Jansen        ndofdwal=intfromfile(2)
526*59599516SKenneth E. Jansen        if(ndofdwal .ne. 0) then
527*59599516SKenneth E. Jansen          if(ndofdwal.ne.1) then
528*59599516SKenneth E. Jansen            warning='WARNING: ndofdwal not equal 1'
529*59599516SKenneth E. Jansen            write(*,*) warning, ndofdwal
530*59599516SKenneth E. Jansen          endif
531*59599516SKenneth E. Jansen          allocate( dwal(nshg) )
532*59599516SKenneth E. Jansen          dwal(:) = -9.87654321e32
533*59599516SKenneth E. Jansen          allocate( qread1(nshg2) )
534*59599516SKenneth E. Jansen          iqsiz=nshg2*1
535*59599516SKenneth E. Jansen          call readdatablock(descriptor,fname1//char(0),qread1,iqsiz,
536*59599516SKenneth E. Jansen     &                   'double'//char(0),iotype)
537*59599516SKenneth E. Jansen          dwal(1:nshg2)=qread1(1:nshg2)
538*59599516SKenneth E. Jansen          deallocate(qread1)
539*59599516SKenneth E. Jansen        else
540*59599516SKenneth E. Jansen          if(myrank==0) then
541*59599516SKenneth E. Jansen            write(*,*) 'WARNING: dwal is missing in the restart files'
542*59599516SKenneth E. Jansen          endif
543*59599516SKenneth E. Jansen          idwal = 0
544*59599516SKenneth E. Jansen        endif
545*59599516SKenneth E. Jansen      endif
546*59599516SKenneth E. Jansen
547*59599516SKenneth E. Jansen!
548*59599516SKenneth E. Jansen!.... close c-binary files
549*59599516SKenneth E. Jansen!
550*59599516SKenneth E. Jansen      call closefile( descriptor, "read"//char(0) )
551*59599516SKenneth E. Jansen      call finalizephmpiio( descriptor )
552*59599516SKenneth E. Jansenc.... ----------------------> Communication tasks <--------------------
553*59599516SKenneth E. Jansenc
554*59599516SKenneth E. Jansen
555*59599516SKenneth E. Jansen      if(numpe > 1) then
556*59599516SKenneth E. Jansen
557*59599516SKenneth E. Jansen        call mpi_barrier(mpi_comm_world, ierr)
558*59599516SKenneth E. Jansen        if(myrank == 0) then
559*59599516SKenneth E. Jansen          write(*,*) 'Reading the geombc-dat files again for ilwork'
560*59599516SKenneth E. Jansen        endif
561*59599516SKenneth E. Jansen
562*59599516SKenneth E. Jansen        color = int(myrank/(numparts/nfiles)) !Should call the color routine in SyncIO here
563*59599516SKenneth E. Jansen        itmp2 = int(log10(float(color+1)))+1
564*59599516SKenneth E. Jansen        write (temp2,"('(''geombcRed-dat.'',i',i1,')')") itmp2
565*59599516SKenneth E. Jansen        write (fnamer,temp2) (color+1)
566*59599516SKenneth E. Jansen        fnamer = trim(fnamer)//char(0)
567*59599516SKenneth E. Jansen
568*59599516SKenneth E. Jansen        ieleven=11
569*59599516SKenneth E. Jansen        ione=1
570*59599516SKenneth E. Jansen
571*59599516SKenneth E. Jansen        itmp = int(log10(float(myrank+1)))+1
572*59599516SKenneth E. Jansen
573*59599516SKenneth E. Jansen        call queryphmpiio(fnamer, nfields, nppf);
574*59599516SKenneth E. Jansen        if (myrank == 0) then
575*59599516SKenneth E. Jansen          write(*,*) 'Number of fields in geombcRed-dat: ',nfields
576*59599516SKenneth E. Jansen          write(*,*) 'Number of parts per file geombcRed-dat: ',nppf
577*59599516SKenneth E. Jansen        endif
578*59599516SKenneth E. Jansen        call initphmpiio( nfields, nppf, nfiles, igeom,
579*59599516SKenneth E. Jansen     &       'read'//char(0))
580*59599516SKenneth E. Jansen        call openfile( fnamer, 'read'//char(0), igeom )
581*59599516SKenneth E. Jansen
582*59599516SKenneth E. Jansen        write (temp1,"('(''size of ilwork array@'',i',i1,',A1)')") itmp
583*59599516SKenneth E. Jansen        write (fname2,temp1) (myrank+1),'?'
584*59599516SKenneth E. Jansen        call readheader(igeom,fname2//char(0),nlwork,ione,
585*59599516SKenneth E. Jansen     &   'integer' // char(0) ,iotype)
586*59599516SKenneth E. Jansen
587*59599516SKenneth E. Jansen        write (temp1,"('(''ilwork@'',i',i1,',A1)')") itmp
588*59599516SKenneth E. Jansen        write (fname2,temp1) (myrank+1),'?'
589*59599516SKenneth E. Jansen        call readheader(igeom,fname2//char(0) ,nlwork,ione,
590*59599516SKenneth E. Jansen     &   'integer'//char(0) , iotype)
591*59599516SKenneth E. Jansen
592*59599516SKenneth E. Jansen        allocate( point2ilwork(nlwork) )
593*59599516SKenneth E. Jansen        allocate( ilworkread(nlwork) )
594*59599516SKenneth E. Jansen        call readdatablock(igeom,fname2//char(0),ilworkread,
595*59599516SKenneth E. Jansen     &                      nlwork,'integer'//char(0) , iotype)
596*59599516SKenneth E. Jansen
597*59599516SKenneth E. Jansen        point2ilwork = ilworkread
598*59599516SKenneth E. Jansen        deallocate(ilworkread)
599*59599516SKenneth E. Jansen
600*59599516SKenneth E. Jansen        call closefile( igeom, "read"//char(0) )
601*59599516SKenneth E. Jansen        call finalizephmpiio( igeom )
602*59599516SKenneth E. Jansen
603*59599516SKenneth E. Jansen        call ctypes (point2ilwork)
604*59599516SKenneth E. Jansen
605*59599516SKenneth E. Jansen      else
606*59599516SKenneth E. Jansen        nlwork=1
607*59599516SKenneth E. Jansen        allocate( point2ilwork(1))
608*59599516SKenneth E. Jansen        nshg0 = nshg
609*59599516SKenneth E. Jansen      endif
610*59599516SKenneth E. Jansen
611*59599516SKenneth E. Jansenc
612*59599516SKenneth E. Jansenc.... -------------------->   communications <-------------------------
613*59599516SKenneth E. Jansenc
614*59599516SKenneth E. Jansen      call mpi_barrier(mpi_comm_world, ierr)  ! make sure every rank is synced here
615*59599516SKenneth E. Jansen      if(myrank == 0) then
616*59599516SKenneth E. Jansen        write(*,*) 'Updating the vertices on the part boundaries'
617*59599516SKenneth E. Jansen      endif
618*59599516SKenneth E. Jansen
619*59599516SKenneth E. Jansen      if (numpe > 1) then
620*59599516SKenneth E. Jansen         ! solution
621*59599516SKenneth E. Jansen          call commuMax (qold, point2ilwork, ndof, 'in '//char(0))
622*59599516SKenneth E. Jansen          call commuMax (qold, point2ilwork, ndof, 'out'//char(0))
623*59599516SKenneth E. Jansen          call mpi_barrier(mpi_comm_world, ierr)  ! make sure everybody is done with ilwork
624*59599516SKenneth E. Jansen
625*59599516SKenneth E. Jansen          ! ybar
626*59599516SKenneth E. Jansen          if(iybar == 1) then
627*59599516SKenneth E. Jansen            call commuMax (ybar, point2ilwork, ndofybar, 'in '//char(0))
628*59599516SKenneth E. Jansen            call commuMax (ybar, point2ilwork, ndofybar, 'out'//char(0))
629*59599516SKenneth E. Jansen            call mpi_barrier(mpi_comm_world, ierr)  ! make sure everybody is done with ilwork
630*59599516SKenneth E. Jansen          endif
631*59599516SKenneth E. Jansen
632*59599516SKenneth E. Jansen         ! errors
633*59599516SKenneth E. Jansen          if(ierror == 1) then
634*59599516SKenneth E. Jansen            call commuMax (errors, point2ilwork, ndoferrors,
635*59599516SKenneth E. Jansen     &                                             'in '//char(0))
636*59599516SKenneth E. Jansen            call commuMax (errors, point2ilwork, ndoferrors,
637*59599516SKenneth E. Jansen     &                                             'out'//char(0))
638*59599516SKenneth E. Jansen            call mpi_barrier(mpi_comm_world, ierr)  ! make sure everybody is done with ilwork
639*59599516SKenneth E. Jansen          endif
640*59599516SKenneth E. Jansen
641*59599516SKenneth E. Jansen          ! phase_average
642*59599516SKenneth E. Jansen          if(numphavg .gt. 0) then
643*59599516SKenneth E. Jansen            do iphavg = 1,numphavg
644*59599516SKenneth E. Jansen              call commuMax (yphbar(:,:,iphavg), point2ilwork,
645*59599516SKenneth E. Jansen     &                                  ndofyphbar, 'in '//char(0))
646*59599516SKenneth E. Jansen              call commuMax (yphbar(:,:,iphavg), point2ilwork,
647*59599516SKenneth E. Jansen     &                                  ndofyphbar, 'out'//char(0))
648*59599516SKenneth E. Jansen              call mpi_barrier(mpi_comm_world, ierr)  ! make sure everybody is done with ilwork
649*59599516SKenneth E. Jansen            enddo
650*59599516SKenneth E. Jansen          endif
651*59599516SKenneth E. Jansen
652*59599516SKenneth E. Jansen          ! vorticity
653*59599516SKenneth E. Jansen          if(ivort == 1) then
654*59599516SKenneth E. Jansen            call commuMax (vort, point2ilwork, ndofvort, 'in '//char(0))
655*59599516SKenneth E. Jansen            call commuMax (vort, point2ilwork, ndofvort, 'out'//char(0))
656*59599516SKenneth E. Jansen            call mpi_barrier(mpi_comm_world, ierr)  ! make sure everybody is done with ilwork
657*59599516SKenneth E. Jansen          endif
658*59599516SKenneth E. Jansen
659*59599516SKenneth E. Jansen         ! dwal
660*59599516SKenneth E. Jansen          if(idwal == 1) then
661*59599516SKenneth E. Jansen            call commuMax (dwal, point2ilwork, 1, 'in '//char(0))
662*59599516SKenneth E. Jansen            call commuMax (dwal, point2ilwork, 1, 'out'//char(0))
663*59599516SKenneth E. Jansen            call mpi_barrier(mpi_comm_world, ierr)  ! make sure everybody is done with ilwork
664*59599516SKenneth E. Jansen          endif
665*59599516SKenneth E. Jansen      endif
666*59599516SKenneth E. Jansen
667*59599516SKenneth E. Jansenc
668*59599516SKenneth E. Jansenc.... -------------------->  Print Min and Max of each field  <-------------------------
669*59599516SKenneth E. Jansenc
670*59599516SKenneth E. Jansen      call mpi_barrier(mpi_comm_world, ierr)
671*59599516SKenneth E. Jansen      if (myrank == 0) then
672*59599516SKenneth E. Jansen        write(*,*) 'Printing min and max of each field component'
673*59599516SKenneth E. Jansen        write(*,*) ''
674*59599516SKenneth E. Jansen      endif
675*59599516SKenneth E. Jansen
676*59599516SKenneth E. Jansen      ! qold
677*59599516SKenneth E. Jansen      do idof=1,ndof
678*59599516SKenneth E. Jansen        call mpi_allreduce(maxval(qold(:,idof)),globmax,1,MPI_DOUBLE,
679*59599516SKenneth E. Jansen     &                                  MPI_MAX, mpi_comm_world, ierr )
680*59599516SKenneth E. Jansen        call mpi_allreduce(minval(qold(:,idof)),globmin,1,MPI_DOUBLE,
681*59599516SKenneth E. Jansen     &                                  MPI_MIN, mpi_comm_world, ierr )
682*59599516SKenneth E. Jansen        if (myrank == 0) then
683*59599516SKenneth E. Jansen          write(*,925) 'max/min qold(',idof,'): ',globmax,globmin
684*59599516SKenneth E. Jansen        endif
685*59599516SKenneth E. Jansen      enddo
686*59599516SKenneth E. Jansen      if (myrank == 0) then
687*59599516SKenneth E. Jansen        write(*,*) ''
688*59599516SKenneth E. Jansen      endif
689*59599516SKenneth E. Jansen
690*59599516SKenneth E. Jansen      ! ybar
691*59599516SKenneth E. Jansen      if(iybar == 1) then
692*59599516SKenneth E. Jansen        do idof=1,ndofybar
693*59599516SKenneth E. Jansen          call mpi_allreduce(maxval(ybar(:,idof)),globmax,1,MPI_DOUBLE,
694*59599516SKenneth E. Jansen     &                                  MPI_MAX, mpi_comm_world, ierr )
695*59599516SKenneth E. Jansen          call mpi_allreduce(minval(ybar(:,idof)),globmin,1,MPI_DOUBLE,
696*59599516SKenneth E. Jansen     &                                  MPI_MIN, mpi_comm_world, ierr )
697*59599516SKenneth E. Jansen          if (myrank == 0) then
698*59599516SKenneth E. Jansen            write(*,925) 'max/min ybar(',idof,'): ',globmax,globmin
699*59599516SKenneth E. Jansen          endif
700*59599516SKenneth E. Jansen        enddo
701*59599516SKenneth E. Jansen        if (myrank == 0) then
702*59599516SKenneth E. Jansen          write(*,*) ''
703*59599516SKenneth E. Jansen        endif
704*59599516SKenneth E. Jansen      endif
705*59599516SKenneth E. Jansen
706*59599516SKenneth E. Jansen      ! errors
707*59599516SKenneth E. Jansen      if(ierror == 1) then
708*59599516SKenneth E. Jansen        do idof=1,ndoferrors
709*59599516SKenneth E. Jansen          call mpi_allreduce(maxval(errors(:,idof)),globmax,1,
710*59599516SKenneth E. Jansen     &                      MPI_DOUBLE,MPI_MAX, mpi_comm_world, ierr )
711*59599516SKenneth E. Jansen          call mpi_allreduce(minval(errors(:,idof)),globmin,1,
712*59599516SKenneth E. Jansen     &                      MPI_DOUBLE,MPI_MIN, mpi_comm_world, ierr )
713*59599516SKenneth E. Jansen          if (myrank == 0) then
714*59599516SKenneth E. Jansen            write(*,925) 'max/min errors(',idof,'): ',globmax,globmin
715*59599516SKenneth E. Jansen          endif
716*59599516SKenneth E. Jansen        enddo
717*59599516SKenneth E. Jansen        if (myrank == 0) then
718*59599516SKenneth E. Jansen          write(*,*) ''
719*59599516SKenneth E. Jansen        endif
720*59599516SKenneth E. Jansen      endif
721*59599516SKenneth E. Jansen
722*59599516SKenneth E. Jansen      ! phase_average
723*59599516SKenneth E. Jansen      if(numphavg .gt. 0) then
724*59599516SKenneth E. Jansen        do iphavg = 1,numphavg
725*59599516SKenneth E. Jansen          do idof=1,ndofyphbar
726*59599516SKenneth E. Jansen            call mpi_allreduce(maxval(yphbar(:,idof,iphavg)), globmax,1,
727*59599516SKenneth E. Jansen     &                     MPI_DOUBLE, MPI_MAX, mpi_comm_world, ierr )
728*59599516SKenneth E. Jansen            call mpi_allreduce(minval(yphbar(:,idof,iphavg)), globmin,1,
729*59599516SKenneth E. Jansen     &                     MPI_DOUBLE, MPI_MIN, mpi_comm_world, ierr )
730*59599516SKenneth E. Jansen            if (myrank == 0) then
731*59599516SKenneth E. Jansen              write(*,926) 'max/min yphbar(',idof,iphavg,'):',
732*59599516SKenneth E. Jansen     &                                globmax,globmin
733*59599516SKenneth E. Jansen            endif
734*59599516SKenneth E. Jansen          enddo
735*59599516SKenneth E. Jansen          if (myrank == 0) then
736*59599516SKenneth E. Jansen            write(*,*) ''
737*59599516SKenneth E. Jansen          endif
738*59599516SKenneth E. Jansen        enddo
739*59599516SKenneth E. Jansen      endif
740*59599516SKenneth E. Jansen
741*59599516SKenneth E. Jansen      ! vorticity
742*59599516SKenneth E. Jansen      if(ivort == 1) then
743*59599516SKenneth E. Jansen        do idof=1,ndofvort
744*59599516SKenneth E. Jansen          call mpi_allreduce(maxval(vort(:,idof)),globmax,1,
745*59599516SKenneth E. Jansen     &                      MPI_DOUBLE,MPI_MAX, mpi_comm_world, ierr )
746*59599516SKenneth E. Jansen          call mpi_allreduce(minval(vort(:,idof)),globmin,1,
747*59599516SKenneth E. Jansen     &                      MPI_DOUBLE,MPI_MIN, mpi_comm_world, ierr )
748*59599516SKenneth E. Jansen          if (myrank == 0) then
749*59599516SKenneth E. Jansen            write(*,925) 'max/min vort(',idof,'): ',globmax,globmin
750*59599516SKenneth E. Jansen          endif
751*59599516SKenneth E. Jansen        enddo
752*59599516SKenneth E. Jansen        if (myrank == 0) then
753*59599516SKenneth E. Jansen          write(*,*) ''
754*59599516SKenneth E. Jansen        endif
755*59599516SKenneth E. Jansen      endif
756*59599516SKenneth E. Jansen
757*59599516SKenneth E. Jansen      ! dwal
758*59599516SKenneth E. Jansen      if(idwal == 1) then
759*59599516SKenneth E. Jansen        call mpi_allreduce(maxval(dwal(:)), globmax, 1, MPI_DOUBLE,
760*59599516SKenneth E. Jansen     &                                 MPI_MAX, mpi_comm_world, ierr )
761*59599516SKenneth E. Jansen        call mpi_allreduce(minval(dwal(:)), globmin, 1, MPI_DOUBLE,
762*59599516SKenneth E. Jansen     &                                 MPI_MIN, mpi_comm_world, ierr )
763*59599516SKenneth E. Jansen        if (myrank == 0) then
764*59599516SKenneth E. Jansen          write(*,925) 'max/min dwal(',1,'):',globmax,globmin
765*59599516SKenneth E. Jansen          write(*,*) ''
766*59599516SKenneth E. Jansen        endif
767*59599516SKenneth E. Jansen      endif
768*59599516SKenneth E. Jansen
769*59599516SKenneth E. Jansen925   format(A,i2,A,2(e24.17,1x))
770*59599516SKenneth E. Jansen926   format(A,i2,i2,A,2(e24.17,1x))
771*59599516SKenneth E. Jansen
772*59599516SKenneth E. Jansenc
773*59599516SKenneth E. Jansenc.... e------------------->  Write data to disks  <-------------------------
774*59599516SKenneth E. Jansenc
775*59599516SKenneth E. Jansen      call mpi_barrier(mpi_comm_world, ierr)
776*59599516SKenneth E. Jansen      if(myrank == 0) then
777*59599516SKenneth E. Jansen        write(*,*) 'Writing the reduced restartRed-dat files'
778*59599516SKenneth E. Jansen      endif
779*59599516SKenneth E. Jansen
780*59599516SKenneth E. Jansen!      call Write_M2NFixBnd(myrank, lstep, nshg,
781*59599516SKenneth E. Jansen!     &          ndof, ndofybar, ndoferrors,
782*59599516SKenneth E. Jansen!     &          qold, ybar, errors, dwal)
783*59599516SKenneth E. Jansen
784*59599516SKenneth E. Jansen      nsynciofieldswriterestart = 1+iybar+ierror+numphavg+ivort+idwal
785*59599516SKenneth E. Jansen      call Write_M2NFixBnd_SolOnly(myrank, lstep, nshg,
786*59599516SKenneth E. Jansen     &          ndof, qold)
787*59599516SKenneth E. Jansen      deallocate(qold)
788*59599516SKenneth E. Jansen
789*59599516SKenneth E. Jansen      ! ybar
790*59599516SKenneth E. Jansen      if(iybar == 1) then
791*59599516SKenneth E. Jansen        call Write_Field(myrank,'a','ybar',4,ybar,'d',
792*59599516SKenneth E. Jansen     &                    nshg,ndofybar,lstep)
793*59599516SKenneth E. Jansen        deallocate(ybar)
794*59599516SKenneth E. Jansen      endif
795*59599516SKenneth E. Jansen
796*59599516SKenneth E. Jansen      ! errors
797*59599516SKenneth E. Jansen      if(ierror == 1) then
798*59599516SKenneth E. Jansen        call Write_Field(myrank,'a','errors',6,errors,'d',
799*59599516SKenneth E. Jansen     &                    nshg,ndoferrors,lstep)
800*59599516SKenneth E. Jansen        deallocate(errors)
801*59599516SKenneth E. Jansen      endif
802*59599516SKenneth E. Jansen
803*59599516SKenneth E. Jansen      ! phase_average
804*59599516SKenneth E. Jansen      if(numphavg .gt. 0) then
805*59599516SKenneth E. Jansen        do iphavg=1,numphavg
806*59599516SKenneth E. Jansen          call write_phavg2(myrank,'a','phase_average',13,
807*59599516SKenneth E. Jansen     &             iphavg,numphavg,yphbar(:,:,iphavg),'d',
808*59599516SKenneth E. Jansen     &             nshg,ndofyphbar,lstep)
809*59599516SKenneth E. Jansen        enddo
810*59599516SKenneth E. Jansen        deallocate(yphbar)
811*59599516SKenneth E. Jansen      endif
812*59599516SKenneth E. Jansen
813*59599516SKenneth E. Jansen      ! vorticity
814*59599516SKenneth E. Jansen      if(ivort == 1) then
815*59599516SKenneth E. Jansen        call Write_Field(myrank,'a','vorticity',9,vort,'d',
816*59599516SKenneth E. Jansen     &                    nshg,ndofvort,lstep)
817*59599516SKenneth E. Jansen        deallocate(vort)
818*59599516SKenneth E. Jansen      endif
819*59599516SKenneth E. Jansen
820*59599516SKenneth E. Jansen      ! dwal
821*59599516SKenneth E. Jansen      if(idwal == 1) then
822*59599516SKenneth E. Jansen        call Write_Field(myrank,'a','dwal',4,dwal,'d',
823*59599516SKenneth E. Jansen     &                    nshg,1,lstep)
824*59599516SKenneth E. Jansen        deallocate(dwal)
825*59599516SKenneth E. Jansen      endif
826*59599516SKenneth E. Jansen
827*59599516SKenneth E. Jansen
828*59599516SKenneth E. Jansen!
829*59599516SKenneth E. Jansen!     Deallocate some remaining memory
830*59599516SKenneth E. Jansen!
831*59599516SKenneth E. Jansen      if(numpe.gt.1) then
832*59599516SKenneth E. Jansen        deallocate(point2ilwork)
833*59599516SKenneth E. Jansen      endif
834*59599516SKenneth E. Jansen
835*59599516SKenneth E. Jansen      return
836*59599516SKenneth E. Jansenc
837*59599516SKenneth E. Jansen 994  call error ('input   ','opening ', igeom)
838*59599516SKenneth E. Jansen 995  call error ('input   ','opening ', igeom)
839*59599516SKenneth E. Jansen 997  call error ('input   ','end file', igeom)
840*59599516SKenneth E. Jansen 998  call error ('input   ','end file', igeom)
841*59599516SKenneth E. Jansenc
842*59599516SKenneth E. Jansen      end
843