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 :: uold(:,:) 21*59599516SKenneth E. Jansen real*8, allocatable :: acold(:,:) 22*59599516SKenneth E. Jansen integer, allocatable :: iBCtmp(:) 23*59599516SKenneth E. Jansen real*8, allocatable :: BCinp(:,:) 24*59599516SKenneth E. Jansen 25*59599516SKenneth E. Jansen integer, allocatable :: point2ilwork(:) 26*59599516SKenneth E. Jansen integer, allocatable :: nBC(:) 27*59599516SKenneth E. Jansen integer, allocatable :: point2iper(:) 28*59599516SKenneth E. Jansen integer, allocatable :: point2ifath(:) 29*59599516SKenneth E. Jansen integer, allocatable :: point2nsons(:) 30*59599516SKenneth E. Jansen 31*59599516SKenneth E. Jansen end module 32*59599516SKenneth E. Jansen 33*59599516SKenneth E. Jansen 34*59599516SKenneth E. Jansen 35*59599516SKenneth E. Jansen subroutine readnblk 36*59599516SKenneth E. Jansenc 37*59599516SKenneth E. Jansen use readarrays 38*59599516SKenneth E. Jansen include "commonAcuStat.h" 39*59599516SKenneth E. Jansen include "mpif.h" 40*59599516SKenneth E. Jansenc 41*59599516SKenneth E. Jansen real*8, allocatable :: xread(:,:),qread(:,:) 42*59599516SKenneth E. Jansen real*8, allocatable :: qread1(:) 43*59599516SKenneth E. Jansen real*8, allocatable :: uread(:,:), acread(:,:) 44*59599516SKenneth E. Jansen real*8, allocatable :: BCinpread(:,:) 45*59599516SKenneth E. Jansen real*8 globmax,globmin 46*59599516SKenneth E. Jansen integer, allocatable :: iperread(:), iBCtmpread(:) 47*59599516SKenneth E. Jansen integer, allocatable :: ilworkread(:), nBCread(:) 48*59599516SKenneth E. Jansen character*10 cname2 49*59599516SKenneth E. Jansen character*30 fmt1 50*59599516SKenneth E. Jansen character*255 fname1,fnamer,fnamelr 51*59599516SKenneth E. Jansen character*255 warning 52*59599516SKenneth E. Jansen 53*59599516SKenneth E. Jansen integer :: descriptor, color, nfiles, nfields 54*59599516SKenneth E. Jansen integer :: numparts, nppf 55*59599516SKenneth E. Jansen character*255 fname2, temp2 56*59599516SKenneth E. Jansen character*64 temp1 57*59599516SKenneth E. Jansen integer :: igeom, ibndc, irstin, ierr 58*59599516SKenneth E. Jansen integer :: ndof, ndoferrors, ndofybar, ndofyphbar 59*59599516SKenneth E. Jansen integer :: itmp, itmp2 60*59599516SKenneth E. Jansen 61*59599516SKenneth E. Jansen integer intfromfile(50) ! integers read from headers 62*59599516SKenneth E. Jansen logical exinput 63*59599516SKenneth E. Jansen 64*59599516SKenneth E. Jansen integer :: numts, numphavg, its, isize, sumts, totsumts, debug 65*59599516SKenneth E. Jansen integer :: ierror, idwal, iphavg, ivort, ndofvort 66*59599516SKenneth E. Jansen integer :: descriptorG,GPID 67*59599516SKenneth E. Jansen integer :: timestep, ithree 68*59599516SKenneth E. Jansen integer :: nshg2,ndof2 69*59599516SKenneth E. Jansen integer, allocatable :: tablets(:,:) 70*59599516SKenneth E. Jansen real*8, allocatable :: qybarcumul(:,:),qerrorcumul(:,:) 71*59599516SKenneth E. Jansen real*8, allocatable :: qyphbarcumul(:,:,:), qvorticity(:,:) 72*59599516SKenneth E. Jansen 73*59599516SKenneth E. Jansen debug = 0 74*59599516SKenneth E. Jansen 75*59599516SKenneth E. Jansen numts=-1 76*59599516SKenneth E. Jansen 77*59599516SKenneth E. Jansen if(myrank == 0) then 78*59599516SKenneth E. Jansen fnamer='AcuStat_input.dat' 79*59599516SKenneth E. Jansen fnamer = trim(fnamer) // char(0) 80*59599516SKenneth E. Jansen inquire(file=fnamer,exist=exinput) 81*59599516SKenneth E. Jansen if(exinput) then 82*59599516SKenneth E. Jansen open(unit=24,file=fnamer,status='old') 83*59599516SKenneth E. Jansen read(24,*) numts 84*59599516SKenneth E. Jansen read(24,*) ierror 85*59599516SKenneth E. Jansen read(24,*) numphavg 86*59599516SKenneth E. Jansen read(24,*) ivort 87*59599516SKenneth E. Jansen read(24,*) idwal 88*59599516SKenneth E. Jansen else 89*59599516SKenneth E. Jansen write(*,*) 'ERROR: Input file ', 90*59599516SKenneth E. Jansen & trim(fnamer),' does not exist!' 91*59599516SKenneth E. Jansen endif 92*59599516SKenneth E. Jansen endif 93*59599516SKenneth E. Jansen 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 ! AcuStat_input.dat does not exist. Quit 97*59599516SKenneth E. Jansen return 98*59599516SKenneth E. Jansen endif 99*59599516SKenneth E. Jansen 100*59599516SKenneth E. Jansen call MPI_Bcast(numts,1,MPI_INTEGER,0,MPI_COMM_WORLD, ierr) 101*59599516SKenneth E. Jansen allocate(tablets(numts,2)) ; tablets=-1 102*59599516SKenneth E. Jansen 103*59599516SKenneth E. Jansen if(myrank == 0 ) then 104*59599516SKenneth E. Jansen ! Read the series of time steps 105*59599516SKenneth E. Jansen do its=1,numts 106*59599516SKenneth E. Jansen read(24,*) tablets(its,1), tablets(its,2) 107*59599516SKenneth E. Jansen enddo 108*59599516SKenneth E. Jansen close(24) 109*59599516SKenneth E. Jansen 110*59599516SKenneth E. Jansen ! Print some info 111*59599516SKenneth E. Jansen write(*,*) 'There are ',numts, 'restart files to read for ybar' 112*59599516SKenneth E. Jansen write(*,*) 'The solution field from the last time step will be', 113*59599516SKenneth E. Jansen & ' added' 114*59599516SKenneth E. Jansen 115*59599516SKenneth E. Jansen if(ierror .gt. 0) then 116*59599516SKenneth E. Jansen write(*,*) 'The error field will also be accumulated' 117*59599516SKenneth E. Jansen ierror = 1 ! security 118*59599516SKenneth E. Jansen else 119*59599516SKenneth E. Jansen write(*,*) 'The error field will NOT be accumulated' 120*59599516SKenneth E. Jansen endif 121*59599516SKenneth E. Jansen 122*59599516SKenneth E. Jansen if(numphavg .gt. 0) then 123*59599516SKenneth E. Jansen write(*,*) 'The phase average fields (',numphavg, 124*59599516SKenneth E. Jansen & ') will also be accumulated' 125*59599516SKenneth E. Jansen else 126*59599516SKenneth E. Jansen write(*,*) 'The phase average fields will NOT be accumulated' 127*59599516SKenneth E. Jansen endif 128*59599516SKenneth E. Jansen 129*59599516SKenneth E. Jansen if(ivort .gt. 0) then 130*59599516SKenneth E. Jansen write(*,*) 'The vorticity field will also be added' 131*59599516SKenneth E. Jansen ivort = 1 ! security 132*59599516SKenneth E. Jansen else 133*59599516SKenneth E. Jansen write(*,*) 'The vorticity field will NOT be added' 134*59599516SKenneth E. Jansen endif 135*59599516SKenneth E. Jansen 136*59599516SKenneth E. Jansen if(idwal .gt. 0) then 137*59599516SKenneth E. Jansen write(*,*) 'The dwal field will also be added' 138*59599516SKenneth E. Jansen idwal = 1 ! security 139*59599516SKenneth E. Jansen else 140*59599516SKenneth E. Jansen write(*,*) 'The dwal field will NOT be added' 141*59599516SKenneth E. Jansen endif 142*59599516SKenneth E. Jansen 143*59599516SKenneth E. Jansen endif 144*59599516SKenneth E. Jansen 145*59599516SKenneth E. Jansen isize=2*numts 146*59599516SKenneth E. Jansen call MPI_Bcast(tablets,isize,MPI_INTEGER,0,MPI_COMM_WORLD, ierr) 147*59599516SKenneth E. Jansen call MPI_Bcast(ierror,1,MPI_INTEGER,0,MPI_COMM_WORLD, ierr) 148*59599516SKenneth E. Jansen call MPI_Bcast(numphavg,1,MPI_INTEGER,0,MPI_COMM_WORLD, ierr) 149*59599516SKenneth E. Jansen call MPI_Bcast(idwal,1,MPI_INTEGER,0,MPI_COMM_WORLD, ierr) 150*59599516SKenneth E. Jansen call MPI_Bcast(ivort,1,MPI_INTEGER,0,MPI_COMM_WORLD, ierr) 151*59599516SKenneth E. Jansen 152*59599516SKenneth E. Jansen! write(*,*) 'my rank is: ',myrank,numts 153*59599516SKenneth E. Jansen 154*59599516SKenneth E. Jansen! do irank=0,7 155*59599516SKenneth E. Jansen! if (irank == myrank) then 156*59599516SKenneth E. Jansen! do its=1,numts 157*59599516SKenneth E. Jansen! write(*,*) myrank,tablets(its,1:2) 158*59599516SKenneth E. Jansen! enddo 159*59599516SKenneth E. Jansen! endif 160*59599516SKenneth E. Jansen! call mpi_barrier(MPI_COMM_WORLD,ierr) 161*59599516SKenneth E. Jansen! enddo 162*59599516SKenneth E. Jansen 163*59599516SKenneth E. Jansen nfiles = nsynciofiles 164*59599516SKenneth E. Jansen numparts = numpe 165*59599516SKenneth E. Jansen color = int(myrank/(numparts/nfiles)) 166*59599516SKenneth E. Jansen sumts = 0 167*59599516SKenneth E. Jansen totsumts = 0 168*59599516SKenneth E. Jansen 169*59599516SKenneth E. Jansen ithree=3 170*59599516SKenneth E. Jansen 171*59599516SKenneth E. Jansen! Loop over the restart files to read and accumulate the ybar field 172*59599516SKenneth E. Jansen do its=1,numts 173*59599516SKenneth E. Jansen 174*59599516SKenneth E. Jansen! 175*59599516SKenneth E. Jansen! Read the table 176*59599516SKenneth E. Jansen! 177*59599516SKenneth E. Jansen timestep = tablets(its,1) 178*59599516SKenneth E. Jansen sumts = tablets(its,2) 179*59599516SKenneth E. Jansen totsumts = totsumts + sumts 180*59599516SKenneth E. Jansen 181*59599516SKenneth E. Jansen! 182*59599516SKenneth E. Jansen! Get the right restart files 183*59599516SKenneth E. Jansen! 184*59599516SKenneth E. Jansen itmp=1 185*59599516SKenneth E. Jansen if (timestep .gt. 0) itmp = int(log10(float(timestep+1)))+1 186*59599516SKenneth E. Jansen write (fmt1,"('(''restart-dat.'',i',i1,',1x)')") itmp 187*59599516SKenneth E. Jansen write (fnamer,fmt1) timestep 188*59599516SKenneth E. Jansen fnamer = trim(fnamer) // cname2(color+1) 189*59599516SKenneth E. Jansen 190*59599516SKenneth E. Jansen if(debug==1) then 191*59599516SKenneth E. Jansen write(*,*)myrank,timestep,nfiles,numparts,trim(adjustl(fnamer)) 192*59599516SKenneth E. Jansen endif 193*59599516SKenneth E. Jansen 194*59599516SKenneth E. Jansen! 195*59599516SKenneth E. Jansen! query, init and open 196*59599516SKenneth E. Jansen! 197*59599516SKenneth E. Jansen call queryphmpiio(fnamer//char(0), nfields, nppf); 198*59599516SKenneth E. Jansen if (myrank == 0) then 199*59599516SKenneth E. Jansen write(*,*) '' 200*59599516SKenneth E. Jansen write(*,*) 'Considering time step',timestep,sumts 201*59599516SKenneth E. Jansen write(*,*) 'Number of fields in restart-dat: ',nfields 202*59599516SKenneth E. Jansen write(*,*) 'Number of parts per file restart-dat: ',nppf 203*59599516SKenneth E. Jansen endif 204*59599516SKenneth E. Jansen 205*59599516SKenneth E. Jansen call initphmpiio(nfields,nppf,nfiles,descriptor,'read'//char(0)) 206*59599516SKenneth E. Jansen call openfile( fnamer//char(0), 'read'//char(0), descriptor ) 207*59599516SKenneth E. Jansen 208*59599516SKenneth E. Jansen! 209*59599516SKenneth E. Jansen! Read the ybar field 210*59599516SKenneth E. Jansen! 211*59599516SKenneth E. Jansen itmp = int(log10(float(myrank+1)))+1 212*59599516SKenneth E. Jansen write (temp1,"('(''ybar@'',i',i1,',A1)')") itmp 213*59599516SKenneth E. Jansen write (fname1,temp1) (myrank+1),'?' 214*59599516SKenneth E. Jansen fname1 = trim(fname1) 215*59599516SKenneth E. Jansen 216*59599516SKenneth E. Jansen intfromfile=0 217*59599516SKenneth E. Jansen call readheader(descriptor,fname1//char(0),intfromfile,ithree, 218*59599516SKenneth E. Jansen & 'integer'//char(0),iotype) 219*59599516SKenneth E. Jansen if(intfromfile(1).ne.0) then 220*59599516SKenneth E. Jansen nshg2=intfromfile(1) !Do not use nshg and ndof from common.h here! 221*59599516SKenneth E. Jansen ndofybar=intfromfile(2) 222*59599516SKenneth E. Jansen !lstep2=intfromfile(3) 223*59599516SKenneth E. Jansen else 224*59599516SKenneth E. Jansen write(*,*)'ERROR: Did not find the field requested 225*59599516SKenneth E. Jansen & in AcuStat_input',myrank 226*59599516SKenneth E. Jansen endif 227*59599516SKenneth E. Jansen 228*59599516SKenneth E. Jansen ! Allocate some memory for the first ts only 229*59599516SKenneth E. Jansen if(its==1) then 230*59599516SKenneth E. Jansen allocate( qybarcumul(nshg2,ndofybar) ) ; qybarcumul = 0.d0 231*59599516SKenneth E. Jansen endif 232*59599516SKenneth E. Jansen 233*59599516SKenneth E. Jansen allocate( qread(nshg2,ndofybar) ) ; qread = 0.d0 234*59599516SKenneth E. Jansen iqsiz = nshg2*ndofybar 235*59599516SKenneth E. Jansen call readdatablock(descriptor,fname1//char(0),qread,iqsiz, 236*59599516SKenneth E. Jansen & 'double'//char(0),iotype) 237*59599516SKenneth E. Jansen 238*59599516SKenneth E. Jansen ! Weight correctly the field and accumulate for average 239*59599516SKenneth E. Jansen qybarcumul(:,:) = qybarcumul(:,:) + qread (:,:)*sumts 240*59599516SKenneth E. Jansen deallocate(qread) 241*59599516SKenneth E. Jansen 242*59599516SKenneth E. Jansen! 243*59599516SKenneth E. Jansen! Read the errors 244*59599516SKenneth E. Jansen! 245*59599516SKenneth E. Jansen if(ierror == 1) then 246*59599516SKenneth E. Jansen write (temp1,"('(''errors@'',i',i1,',a1)')") 247*59599516SKenneth E. Jansen & itmp 248*59599516SKenneth E. Jansen write (fname1,temp1) (myrank+1),'?' 249*59599516SKenneth E. Jansen fname1 = trim(fname1) 250*59599516SKenneth E. Jansen intfromfile=0 251*59599516SKenneth E. Jansen call readheader(descriptor,fname1//char(0),intfromfile,ithree, 252*59599516SKenneth E. Jansen & 'integer'//char(0),iotype) 253*59599516SKenneth E. Jansen if(intfromfile(1).ne.0) then 254*59599516SKenneth E. Jansen nshg2=intfromfile(1) 255*59599516SKenneth E. Jansen ndoferrors=intfromfile(2) 256*59599516SKenneth E. Jansen !lstep=intfromfile(3) 257*59599516SKenneth E. Jansen else 258*59599516SKenneth E. Jansen write(*,*)'ERROR: Did not find the field requested 259*59599516SKenneth E. Jansen & in AcuStat_input',myrank 260*59599516SKenneth E. Jansen endif 261*59599516SKenneth E. Jansen 262*59599516SKenneth E. Jansen ! Allocate some memory for the first ts only 263*59599516SKenneth E. Jansen if(its==1) then 264*59599516SKenneth E. Jansen allocate( qerrorcumul(nshg2,ndoferrors) ) ; qerrorcumul = 0.d0 265*59599516SKenneth E. Jansen endif 266*59599516SKenneth E. Jansen 267*59599516SKenneth E. Jansen allocate( qread(nshg2,ndoferrors) ) ; qread = 0.d0 268*59599516SKenneth E. Jansen iqsiz=nshg2*ndoferrors 269*59599516SKenneth E. Jansen call readdatablock(descriptor,fname1//char(0),qread,iqsiz, 270*59599516SKenneth E. Jansen & 'double'//char(0),iotype) 271*59599516SKenneth E. Jansen 272*59599516SKenneth E. Jansen ! Weight correctly the field and accumulate for average 273*59599516SKenneth E. Jansen qerrorcumul(:,:) = qerrorcumul(:,:) + qread(:,:)*sumts 274*59599516SKenneth E. Jansen deallocate(qread) 275*59599516SKenneth E. Jansen endif ! if ierror == 1 276*59599516SKenneth E. Jansen 277*59599516SKenneth E. Jansen! 278*59599516SKenneth E. Jansen! Read the phase_average fields 279*59599516SKenneth E. Jansen! 280*59599516SKenneth E. Jansen if(numphavg .gt. 0) then 281*59599516SKenneth E. Jansen do iphavg = 1,numphavg 282*59599516SKenneth E. Jansen itmp = int(log10(float(myrank+1)))+1 283*59599516SKenneth E. Jansen itmp2 = int(log10(float(iphavg)))+1 284*59599516SKenneth E. Jansen! write (temp1,"('(''phase_average@'',i',i1,',A1)')") itmp 285*59599516SKenneth E. Jansen! write(*,*) temp1 ! ('phase_average@',i1,A1) 286*59599516SKenneth E. Jansen write (temp1, 287*59599516SKenneth E. Jansen & "('(''phase_average'',i',i1,',''@'',i',i1,',A1)')") 288*59599516SKenneth E. Jansen & itmp2, itmp 289*59599516SKenneth E. Jansen! write(*,*) temp1 ! ('phase_average',i1,'@',i1,A1) 290*59599516SKenneth E. Jansen write (fname1,temp1) iphavg,(myrank+1),'?' 291*59599516SKenneth E. Jansen fname1 = trim(fname1) 292*59599516SKenneth E. Jansen if(debug == 1) then 293*59599516SKenneth E. Jansen if(myrank == 0 .or. myrank == numpe-1) then 294*59599516SKenneth E. Jansen write(*,*)'Rank: ',myrank,'- phase: ', 295*59599516SKenneth E. Jansen & iphavg,' - field: ',trim(fname1) 296*59599516SKenneth E. Jansen endif 297*59599516SKenneth E. Jansen endif 298*59599516SKenneth E. Jansen intfromfile=0 299*59599516SKenneth E. Jansen call readheader(descriptor,fname1//char(0),intfromfile, 300*59599516SKenneth E. Jansen & ithree,'integer'//char(0),iotype) 301*59599516SKenneth E. Jansen if(intfromfile(1).ne.0) then 302*59599516SKenneth E. Jansen nshg2=intfromfile(1) !Do not use nshg and ndof from common.h here! 303*59599516SKenneth E. Jansen ndofyphbar=intfromfile(2) 304*59599516SKenneth E. Jansen !lstep2=intfromfile(3) 305*59599516SKenneth E. Jansen else 306*59599516SKenneth E. Jansen write(*,*)'ERROR: Did not find the field requested 307*59599516SKenneth E. Jansen & in AcuStat_input',myrank 308*59599516SKenneth E. Jansen endif 309*59599516SKenneth E. Jansen 310*59599516SKenneth E. Jansen ! Allocate some memory for the first ts only 311*59599516SKenneth E. Jansen if(its==1 .and. iphavg==1) then 312*59599516SKenneth E. Jansen allocate( qyphbarcumul(nshg2,ndofyphbar,numphavg) ) 313*59599516SKenneth E. Jansen qyphbarcumul = 0.d0 314*59599516SKenneth E. Jansen endif 315*59599516SKenneth E. Jansen 316*59599516SKenneth E. Jansen allocate( qread(nshg2,ndofyphbar) ) ; qread = 0.d0 317*59599516SKenneth E. Jansen iqsiz = nshg2*ndofyphbar 318*59599516SKenneth E. Jansen call readdatablock(descriptor,fname1//char(0),qread,iqsiz, 319*59599516SKenneth E. Jansen & 'double'//char(0),iotype) 320*59599516SKenneth E. Jansen 321*59599516SKenneth E. Jansen ! Weight correctly the field and accumulate for average 322*59599516SKenneth E. Jansen qyphbarcumul(:,:,iphavg) = qyphbarcumul(:,:,iphavg) 323*59599516SKenneth E. Jansen & + qread (:,:)*sumts 324*59599516SKenneth E. Jansen deallocate(qread) 325*59599516SKenneth E. Jansen enddo 326*59599516SKenneth E. Jansen endif 327*59599516SKenneth E. Jansen! 328*59599516SKenneth E. Jansen! If last ts, read also the solution, error and dwal field 329*59599516SKenneth E. Jansen! for reduction with M2N or for convenience 330*59599516SKenneth E. Jansen! 331*59599516SKenneth E. Jansen 332*59599516SKenneth E. Jansen if(its == numts) then 333*59599516SKenneth E. Jansen! 334*59599516SKenneth E. Jansen! Read the solution of the last time step 335*59599516SKenneth E. Jansen! 336*59599516SKenneth E. Jansen itmp = int(log10(float(myrank+1)))+1 337*59599516SKenneth E. Jansen write (temp1,"('(''solution@'',i',i1,',A1)')") itmp 338*59599516SKenneth E. Jansen write (fname1,temp1) (myrank+1),'?' 339*59599516SKenneth E. Jansen fname1 = trim(fname1) 340*59599516SKenneth E. Jansen 341*59599516SKenneth E. Jansen intfromfile=0 342*59599516SKenneth E. Jansen call readheader(descriptor,fname1//char(0) ,intfromfile, 343*59599516SKenneth E. Jansen & ithree,'integer'//char(0), iotype) 344*59599516SKenneth E. Jansen 345*59599516SKenneth E. Jansen if(intfromfile(1).ne.0) then 346*59599516SKenneth E. Jansen nshg2=intfromfile(1) 347*59599516SKenneth E. Jansen ndof2=intfromfile(2) 348*59599516SKenneth E. Jansen allocate( qold(nshg2,ndof2) ) 349*59599516SKenneth E. Jansen qold(:,:) = -9.87654321e32 350*59599516SKenneth E. Jansen lstep=intfromfile(3) 351*59599516SKenneth E. Jansen allocate( qread(nshg2,ndof2) ) ; qread = 0.d0 352*59599516SKenneth E. Jansen iqsiz=nshg2*ndof2 353*59599516SKenneth E. Jansen call readdatablock(descriptor,fname1//char(0),qread,iqsiz, 354*59599516SKenneth E. Jansen & 'double'//char(0),iotype) 355*59599516SKenneth E. Jansen qold(1:nshg2,1:ndof2) = qread(1:nshg2,1:ndof2) 356*59599516SKenneth E. Jansen deallocate(qread) 357*59599516SKenneth E. Jansen else 358*59599516SKenneth E. Jansen write(*,*)'ERROR: Did not find the field requested 359*59599516SKenneth E. Jansen & in AcuStat_input',myrank 360*59599516SKenneth E. Jansen endif 361*59599516SKenneth E. Jansen 362*59599516SKenneth E. Jansen! 363*59599516SKenneth E. Jansen! Read the vorticity field from the last time step 364*59599516SKenneth E. Jansen! 365*59599516SKenneth E. Jansen if(ivort == 1) then 366*59599516SKenneth E. Jansen write (temp1,"('(''vorticity@'',i',i1,',a1)')") 367*59599516SKenneth E. Jansen & itmp 368*59599516SKenneth E. Jansen write (fname1,temp1) (myrank+1),'?' 369*59599516SKenneth E. Jansen fname1 = trim(fname1) 370*59599516SKenneth E. Jansen intfromfile=0 371*59599516SKenneth E. Jansen call readheader(descriptor,fname1//char(0),intfromfile, 372*59599516SKenneth E. Jansen & ithree,'integer'//char(0),iotype) 373*59599516SKenneth E. Jansen if(intfromfile(1).ne.0) then 374*59599516SKenneth E. Jansen nshg2=intfromfile(1) 375*59599516SKenneth E. Jansen ndofvort=intfromfile(2) 376*59599516SKenneth E. Jansen allocate( qvorticity(nshg2,ndofvort) ) 377*59599516SKenneth E. Jansen qvorticity(:,:) = -9.87654321e32 378*59599516SKenneth E. Jansen allocate( qread(nshg2,ndofvort) ) ; qread = 0.d0 379*59599516SKenneth E. Jansen iqsiz=nshg2*ndofvort 380*59599516SKenneth E. Jansen call readdatablock(descriptor,fname1//char(0),qread, 381*59599516SKenneth E. Jansen & iqsiz,'double'//char(0),iotype) 382*59599516SKenneth E. Jansen qvorticity(1:nshg2,1:ndofvort) = qread(1:nshg2,1:ndofvort) 383*59599516SKenneth E. Jansen deallocate(qread) 384*59599516SKenneth E. Jansen else 385*59599516SKenneth E. Jansen write(*,*)'ERROR: Did not find the field requested 386*59599516SKenneth E. Jansen & in AcuStat_input',myrank 387*59599516SKenneth E. Jansen endif 388*59599516SKenneth E. Jansen endif 389*59599516SKenneth E. Jansen 390*59599516SKenneth E. Jansen! 391*59599516SKenneth E. Jansen! Read the dwal 392*59599516SKenneth E. Jansen! 393*59599516SKenneth E. Jansen if(idwal == 1) then 394*59599516SKenneth E. Jansen write (temp1,"('(''dwal@'',i',i1,',a1)')") 395*59599516SKenneth E. Jansen & itmp 396*59599516SKenneth E. Jansen write (fname1,temp1) (myrank+1),'?' 397*59599516SKenneth E. Jansen fname1 = trim(fname1) 398*59599516SKenneth E. Jansen intfromfile=0 399*59599516SKenneth E. Jansen call readheader(descriptor,fname1//char(0),intfromfile, 400*59599516SKenneth E. Jansen & ithree,'integer'//char(0),iotype) 401*59599516SKenneth E. Jansen if(intfromfile(1).ne.0) then 402*59599516SKenneth E. Jansen nshg2=intfromfile(1) 403*59599516SKenneth E. Jansen allocate( dwal(nshg2) ) 404*59599516SKenneth E. Jansen dwal(:) = -9.87654321e32 405*59599516SKenneth E. Jansen allocate( qread1(nshg2) ) ; qread1 = 0.d0 406*59599516SKenneth E. Jansen iqsiz=nshg2*1 407*59599516SKenneth E. Jansen call readdatablock(descriptor,fname1//char(0),qread1, 408*59599516SKenneth E. Jansen & iqsiz,'double'//char(0),iotype) 409*59599516SKenneth E. Jansen dwal(1:nshg2) = qread1(1:nshg2) 410*59599516SKenneth E. Jansen deallocate(qread1) 411*59599516SKenneth E. Jansen else 412*59599516SKenneth E. Jansen write(*,*)'ERROR: Did not find the field requested 413*59599516SKenneth E. Jansen & in AcuStat_input',myrank 414*59599516SKenneth E. Jansen endif 415*59599516SKenneth E. Jansen endif !if idwal == 1 416*59599516SKenneth E. Jansen 417*59599516SKenneth E. Jansen 418*59599516SKenneth E. Jansen 419*59599516SKenneth E. Jansen endif ! First time step 420*59599516SKenneth E. Jansen! 421*59599516SKenneth E. Jansen! Close and finalize 422*59599516SKenneth E. Jansen! 423*59599516SKenneth E. Jansen call closefile( descriptor, 'read'//char(0) ) 424*59599516SKenneth E. Jansen call finalizephmpiio( descriptor ) 425*59599516SKenneth E. Jansen 426*59599516SKenneth E. Jansen enddo !numts = number of restart files to read 427*59599516SKenneth E. Jansen 428*59599516SKenneth E. Jansen! 429*59599516SKenneth E. Jansen! Weight correctly the average field 430*59599516SKenneth E. Jansen! 431*59599516SKenneth E. Jansen qybarcumul(:,:) = qybarcumul(:,:)/totsumts 432*59599516SKenneth E. Jansen if (ierror == 1) then 433*59599516SKenneth E. Jansen qerrorcumul(:,:) = qerrorcumul(:,:)/totsumts 434*59599516SKenneth E. Jansen endif 435*59599516SKenneth E. Jansen if(numphavg .gt. 0) then 436*59599516SKenneth E. Jansen do iphavg = 1,numphavg 437*59599516SKenneth E. Jansen qyphbarcumul(:,:,iphavg) = qyphbarcumul(:,:,iphavg)/totsumts 438*59599516SKenneth E. Jansen enddo 439*59599516SKenneth E. Jansen endif 440*59599516SKenneth E. Jansen 441*59599516SKenneth E. Jansen! Write qold, dwal, errors and qybarcumul (= cumulated ybar) 442*59599516SKenneth E. Jansen 443*59599516SKenneth E. Jansen nsynciofieldswriterestart = 2 + ierror + numphavg + ivort + idwal ! 2 for solution and accumulated ybar 444*59599516SKenneth E. Jansen! write(*,*) 'nsynciofieldswriterestart:',nsynciofieldswriterestart 445*59599516SKenneth E. Jansen 446*59599516SKenneth E. Jansen! call Write_AcuStat(myrank, tablets(1,1), tablets(numts,1), 447*59599516SKenneth E. Jansen! & totsumts, nshg2, ndof2, ndoferrors, ndofybar, 448*59599516SKenneth E. Jansen! & qold, dwal, qerrorcumul, qybarcumul) 449*59599516SKenneth E. Jansen 450*59599516SKenneth E. Jansen call mpi_barrier(mpi_comm_world, ierr) 451*59599516SKenneth E. Jansen call Write_AcuStat(myrank, tablets(1,1), tablets(numts,1), 452*59599516SKenneth E. Jansen & totsumts, nshg2, ndof2, qold) 453*59599516SKenneth E. Jansen 454*59599516SKenneth E. Jansen if(ivort == 1) then 455*59599516SKenneth E. Jansen call write_field(myrank,'a','vorticity',9,qvorticity, 456*59599516SKenneth E. Jansen & 'd', nshg2, ndofvort, tablets(numts,1)) 457*59599516SKenneth E. Jansen endif 458*59599516SKenneth E. Jansen 459*59599516SKenneth E. Jansen if (ierror == 1) then 460*59599516SKenneth E. Jansen call write_error(myrank, tablets(numts,1), nshg2, ndoferrors, 461*59599516SKenneth E. Jansen & qerrorcumul) 462*59599516SKenneth E. Jansen endif 463*59599516SKenneth E. Jansen 464*59599516SKenneth E. Jansen call write_field(myrank,'a', 'ybar', 4, 465*59599516SKenneth E. Jansen & qybarcumul, 'd', nshg2, ndofybar, tablets(numts,1)) 466*59599516SKenneth E. Jansen 467*59599516SKenneth E. Jansen if(numphavg .gt. 0) then 468*59599516SKenneth E. Jansen do iphavg = 1,numphavg 469*59599516SKenneth E. Jansen call write_phavg2(myrank,'a','phase_average',13,iphavg, 470*59599516SKenneth E. Jansen & numphavg,qyphbarcumul(:,:,iphavg),'d',nshg2,ndofyphbar, 471*59599516SKenneth E. Jansen & tablets(numts,1)) 472*59599516SKenneth E. Jansen enddo 473*59599516SKenneth E. Jansen endif 474*59599516SKenneth E. Jansen 475*59599516SKenneth E. Jansen if (idwal == 1) then 476*59599516SKenneth E. Jansen call write_field(myrank, 'a', 'dwal', 4, dwal, 477*59599516SKenneth E. Jansen & 'd', nshg2, 1, tablets(numts,1)) 478*59599516SKenneth E. Jansen endif 479*59599516SKenneth E. Jansen 480*59599516SKenneth E. Jansen! 481*59599516SKenneth E. Jansen! Free memory 482*59599516SKenneth E. Jansen! 483*59599516SKenneth E. Jansen deallocate(tablets) 484*59599516SKenneth E. Jansen deallocate(qybarcumul, qold) 485*59599516SKenneth E. Jansen if(ierror == 1) then 486*59599516SKenneth E. Jansen deallocate(qerrorcumul) 487*59599516SKenneth E. Jansen endif 488*59599516SKenneth E. Jansen if(numphavg .gt. 0) then 489*59599516SKenneth E. Jansen deallocate(qyphbarcumul) 490*59599516SKenneth E. Jansen endif 491*59599516SKenneth E. Jansen if(ivort == 1) then 492*59599516SKenneth E. Jansen deallocate(qvorticity) 493*59599516SKenneth E. Jansen endif 494*59599516SKenneth E. Jansen if(idwal == 1) then 495*59599516SKenneth E. Jansen deallocate(dwal) 496*59599516SKenneth E. Jansen endif 497*59599516SKenneth E. Jansen 498*59599516SKenneth E. Jansen return 499*59599516SKenneth E. Jansen 500*59599516SKenneth E. Jansen end 501