xref: /phasta/phSolver/common/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 :: uold(:,:)
18*59599516SKenneth E. Jansen      real*8, allocatable :: acold(:,:)
19*59599516SKenneth E. Jansen      integer, allocatable :: iBCtmp(:)
20*59599516SKenneth E. Jansen      real*8, allocatable :: BCinp(:,:)
21*59599516SKenneth E. Jansen
22*59599516SKenneth E. Jansen      integer, allocatable :: point2ilwork(:)
23*59599516SKenneth E. Jansen      integer, allocatable :: nBC(:)
24*59599516SKenneth E. Jansen      integer, allocatable :: point2iper(:)
25*59599516SKenneth E. Jansen      integer, allocatable :: point2ifath(:)
26*59599516SKenneth E. Jansen      integer, allocatable :: point2nsons(:)
27*59599516SKenneth E. Jansen
28*59599516SKenneth E. Jansen      end module
29*59599516SKenneth E. Jansen
30*59599516SKenneth E. Jansen
31*59599516SKenneth E. Jansen
32*59599516SKenneth E. Jansen      subroutine readnblk
33*59599516SKenneth E. Jansenc
34*59599516SKenneth E. Jansen      use readarrays
35*59599516SKenneth E. Jansen      include "common.h"
36*59599516SKenneth E. Jansenc
37*59599516SKenneth E. Jansen      real*8, allocatable :: xread(:,:), qread(:,:), acread(:,:)
38*59599516SKenneth E. Jansen      real*8, allocatable :: uread(:,:)
39*59599516SKenneth E. Jansen      real*8, allocatable :: BCinpread(:,:)
40*59599516SKenneth E. Jansen      integer, allocatable :: iperread(:), iBCtmpread(:)
41*59599516SKenneth E. Jansen      integer, allocatable :: ilworkread(:), nBCread(:)
42*59599516SKenneth E. Jansen      character*10 cname2
43*59599516SKenneth E. Jansen      character*8 mach2
44*59599516SKenneth E. Jansen!MR CHANGE
45*59599516SKenneth E. Jansen!      character*20 fmt1
46*59599516SKenneth E. Jansen      character*30 fmt1
47*59599516SKenneth E. Jansen!MR CHANGE END
48*59599516SKenneth E. Jansen      character*255 fname1,fnamer,fnamelr
49*59599516SKenneth E. Jansen      character*255 warning
50*59599516SKenneth E. Jansen      integer igeomBAK, ibndc, irstin, ierr
51*59599516SKenneth E. Jansen      integer intfromfile(50) ! integers read from headers
52*59599516SKenneth E. Jansen
53*59599516SKenneth E. Jansencccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
54*59599516SKenneth E. Jansenccccccccccccccccccccccccccccccccccccccccc New PhastaIO Definition Part ccccccccccccccccccccccccccccccccccccccccc
55*59599516SKenneth E. Jansen
56*59599516SKenneth E. Jansen      integer :: descriptor, descriptorG, GPID, color, nfiles, nfields
57*59599516SKenneth E. Jansen      integer ::  numparts, nppf
58*59599516SKenneth E. Jansen      integer :: ierr_io, numprocs, itmp, itmp2
59*59599516SKenneth E. Jansen      character*255 fname2, temp2
60*59599516SKenneth E. Jansen      character*64 temp1
61*59599516SKenneth E. Jansen
62*59599516SKenneth E. Jansencccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
63*59599516SKenneth E. Jansencccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
64*59599516SKenneth E. Jansen
65*59599516SKenneth E. Jansen
66*59599516SKenneth E. Jansenc
67*59599516SKenneth E. Jansenc.... determine the step number to start with
68*59599516SKenneth E. Jansenc
69*59599516SKenneth E. Jansen      open(unit=72,file='numstart.dat',status='old')
70*59599516SKenneth E. Jansen      read(72,*) irstart
71*59599516SKenneth E. Jansen      close(72)
72*59599516SKenneth E. Jansen      lstep=irstart ! in case restart files have no fields
73*59599516SKenneth E. Jansen
74*59599516SKenneth E. Jansenc
75*59599516SKenneth E. Jansen      fname1='geombc.dat'
76*59599516SKenneth E. Jansen      fname1= trim(fname1)  // cname2(myrank+1)
77*59599516SKenneth E. Jansenc      fnamelr='restart.latest'
78*59599516SKenneth E. Jansen
79*59599516SKenneth E. Jansen      itmp=1
80*59599516SKenneth E. Jansen      if (irstart .gt. 0) itmp = int(log10(float(irstart)))+1
81*59599516SKenneth E. Jansen      write (fmt1,"('(''restart.'',i',i1,',1x)')") itmp
82*59599516SKenneth E. Jansen      write (fnamer,fmt1) irstart
83*59599516SKenneth E. Jansen      fnamer = trim(fnamer) // cname2(myrank+1)
84*59599516SKenneth E. Jansenc      fnamelr = trim(fnamelr) // cname2(myrank+1)
85*59599516SKenneth E. Jansen
86*59599516SKenneth E. Jansenc
87*59599516SKenneth E. Jansenc.... open input files
88*59599516SKenneth E. Jansenc
89*59599516SKenneth E. Jansen
90*59599516SKenneth E. Jansen
91*59599516SKenneth E. Jansenc      call openfile(  fname1,  'read?', igeomBAK );
92*59599516SKenneth E. Jansen
93*59599516SKenneth E. Jansen
94*59599516SKenneth E. Jansenc
95*59599516SKenneth E. Jansenc.... try opening restart.latest.proc before trying restart.stepno.proc
96*59599516SKenneth E. Jansenc
97*59599516SKenneth E. Jansenc      call openfile(  fnamelr,  'read?', irstin );
98*59599516SKenneth E. Jansenc      if ( irstin .eq. 0 )
99*59599516SKenneth E. Jansen
100*59599516SKenneth E. Jansen!MR CHANGE
101*59599516SKenneth E. Jansen!       call openfile( fnamer, 'read?', irstin );
102*59599516SKenneth E. Jansen!MR CHANGE END
103*59599516SKenneth E. Jansen
104*59599516SKenneth E. Jansen! either one will work
105*59599516SKenneth E. Jansenc
106*59599516SKenneth E. Jansenc.... input the geometry parameters
107*59599516SKenneth E. Jansenc
108*59599516SKenneth E. Jansen
109*59599516SKenneth E. Jansencccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
110*59599516SKenneth E. Jansen!MR CHANGE
111*59599516SKenneth E. Jansen
112*59599516SKenneth E. Jansen!      nfiles = 2
113*59599516SKenneth E. Jansen!      nfields = 31
114*59599516SKenneth E. Jansen!      numparts = 8
115*59599516SKenneth E. Jansen!      nppp = numparts/numpe
116*59599516SKenneth E. Jansen!      nppf = numparts/nfiles
117*59599516SKenneth E. Jansen
118*59599516SKenneth E. Jansen      nfiles = nsynciofiles
119*59599516SKenneth E. Jansen!      nfields = nsynciofieldsreadgeombc
120*59599516SKenneth E. Jansen      numparts = numpe !This is the common settings. Beware if you try to compute several parts per process
121*59599516SKenneth E. Jansen
122*59599516SKenneth E. Jansen!      nppp = numparts/numpe
123*59599516SKenneth E. Jansen!      nppf = numparts/nfiles
124*59599516SKenneth E. Jansen!MR CHANGE END
125*59599516SKenneth E. Jansen
126*59599516SKenneth E. Jansenc      fnamer="/home/nliu/develop/PIG4/512-procs_case/geombc-dat"
127*59599516SKenneth E. Jansen
128*59599516SKenneth E. Jansen      color = int(myrank/(numparts/nfiles)) !Should call the color routine in SyncIO here
129*59599516SKenneth E. Jansen      itmp2 = int(log10(float(color+1)))+1
130*59599516SKenneth E. Jansen      write (temp2,"('(''geombc-dat.'',i',i1,')')") itmp2
131*59599516SKenneth E. Jansen      write (fnamer,temp2) (color+1)
132*59599516SKenneth E. Jansen      fnamer = trim(fnamer) // char(0)
133*59599516SKenneth E. Jansen
134*59599516SKenneth E. Jansenc      fnamer="/home/nliu/develop/test-case/512-procs_case/geombc-dat"
135*59599516SKenneth E. Jansen
136*59599516SKenneth E. Jansen      itwo=2
137*59599516SKenneth E. Jansen      ione=1
138*59599516SKenneth E. Jansen      ieleven=11
139*59599516SKenneth E. Jansen      itmp = int(log10(float(myrank+1)))+1
140*59599516SKenneth E. Jansen
141*59599516SKenneth E. Jansen      call queryphmpiio(fnamer, nfields, nppf);
142*59599516SKenneth E. Jansen      if (myrank == 0) then
143*59599516SKenneth E. Jansen        write(*,*) 'Number of fields in geombc-dat: ',nfields
144*59599516SKenneth E. Jansen        write(*,*) 'Number of parts per file geombc-dat: ',nppf
145*59599516SKenneth E. Jansen      endif
146*59599516SKenneth E. Jansen      call initphmpiio( nfields, nppf, nfiles, igeom,
147*59599516SKenneth E. Jansen     & 'read' // char(0))
148*59599516SKenneth E. Jansen      call openfile( fnamer, 'read' // char(0), igeom )
149*59599516SKenneth E. Jansen
150*59599516SKenneth E. Jansen      write (temp1,"('(''number of nodes@'',i',i1,',A1)')") itmp
151*59599516SKenneth E. Jansen      write (fname2,temp1) (myrank+1),'?'
152*59599516SKenneth E. Jansen      call readheader(igeom,fname2 // char(0),numnp,ione,
153*59599516SKenneth E. Jansen     & 'integer' // char(0), iotype)
154*59599516SKenneth E. Jansen
155*59599516SKenneth E. Jansen      write (temp1,"('(''number of modes@'',i',i1,',A1)')") itmp
156*59599516SKenneth E. Jansen      write (fname2,temp1) (myrank+1),'?'
157*59599516SKenneth E. Jansen      call readheader(igeom,fname2 // char(0),nshg,ione,
158*59599516SKenneth E. Jansen     & 'integer' // char(0), iotype)
159*59599516SKenneth E. Jansen
160*59599516SKenneth E. Jansen      write (temp1,"('(''number of interior elements@'',i',i1,',A1)')")
161*59599516SKenneth E. Jansen     &       itmp
162*59599516SKenneth E. Jansen      write (fname2,temp1) (myrank+1),'?'
163*59599516SKenneth E. Jansen      call readheader(igeom,fname2 // char(0),numel,ione,
164*59599516SKenneth E. Jansen     & 'integer' // char(0), iotype)
165*59599516SKenneth E. Jansen
166*59599516SKenneth E. Jansen      write (temp1,"('(''number of boundary elements@'',i',i1,',A1)')")
167*59599516SKenneth E. Jansen     &       itmp
168*59599516SKenneth E. Jansen      write (fname2,temp1) (myrank+1),'?'
169*59599516SKenneth E. Jansen      call readheader(igeom,fname2 // char(0),numelb,ione,
170*59599516SKenneth E. Jansen     & 'integer' // char(0),iotype)
171*59599516SKenneth E. Jansen
172*59599516SKenneth E. Jansen      write (temp1,
173*59599516SKenneth E. Jansen     & "('(''maximum number of element nodes@'',i',i1,',A1)')") itmp
174*59599516SKenneth E. Jansen      write (fname2,temp1) (myrank+1),'?'
175*59599516SKenneth E. Jansen      call readheader(igeom,fname2 // char(0),nen,ione,
176*59599516SKenneth E. Jansen     &'integer' // char(0),iotype)
177*59599516SKenneth E. Jansen
178*59599516SKenneth E. Jansen      write (temp1,"('(''number of interior tpblocks@'',i',i1,',A1)')")
179*59599516SKenneth E. Jansen     &       itmp
180*59599516SKenneth E. Jansen      write (fname2,temp1) (myrank+1),'?'
181*59599516SKenneth E. Jansen      call readheader(igeom,fname2 // char(0),nelblk,ione,
182*59599516SKenneth E. Jansen     & 'integer' // char(0) ,iotype)
183*59599516SKenneth E. Jansen
184*59599516SKenneth E. Jansen      write (temp1,"('(''number of boundary tpblocks@'',i',i1,',A1)')")
185*59599516SKenneth E. Jansen     &       itmp
186*59599516SKenneth E. Jansen      write (fname2,temp1) (myrank+1),'?'
187*59599516SKenneth E. Jansen      call readheader(igeom,fname2 // char(0),nelblb,ione,
188*59599516SKenneth E. Jansen     & 'integer' // char(0), iotype)
189*59599516SKenneth E. Jansen
190*59599516SKenneth E. Jansen      write (temp1,
191*59599516SKenneth E. Jansen     & "('(''number of nodes with Dirichlet BCs@'',i',i1,',A1)')") itmp
192*59599516SKenneth E. Jansen      write (fname2,temp1) (myrank+1),'?'
193*59599516SKenneth E. Jansen      call readheader(igeom,fname2 // char(0),numpbc,ione,
194*59599516SKenneth E. Jansen     & 'integer' // char(0),iotype)
195*59599516SKenneth E. Jansen
196*59599516SKenneth E. Jansen      write (temp1,"('(''number of shape functions@'',i',i1,',A1)')")
197*59599516SKenneth E. Jansen     &       itmp
198*59599516SKenneth E. Jansen      write (fname2,temp1) (myrank+1),'?'
199*59599516SKenneth E. Jansen      call readheader(igeom,fname2 // char(0),ntopsh,ione,
200*59599516SKenneth E. Jansen     & 'integer' // char(0),iotype)
201*59599516SKenneth E. Jansen
202*59599516SKenneth E. Jansenc      call closefile( igeom, "read" )
203*59599516SKenneth E. Jansenc      call finalizephmpiio( igeom )
204*59599516SKenneth E. Jansen
205*59599516SKenneth E. Jansen!       if(myrank==0) then
206*59599516SKenneth E. Jansen!          print *, "Reading Finished, ********* : ", trim(fname2)
207*59599516SKenneth E. Jansen!       endif
208*59599516SKenneth E. Jansen
209*59599516SKenneth E. Jansen
210*59599516SKenneth E. Jansenccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
211*59599516SKenneth E. Jansen
212*59599516SKenneth E. Jansenc      ieleven=11
213*59599516SKenneth E. Jansenc      ione=1
214*59599516SKenneth E. Jansenc      fname1='number of nodes?'
215*59599516SKenneth E. Jansenc      call readheader(igeomBAK,fname1,numnp,ione,'integer', iotype)
216*59599516SKenneth E. Jansenc      fname1='number of modes?'
217*59599516SKenneth E. Jansenc      call readheader(igeomBAK,fname1,nshg,ione,'integer', iotype)
218*59599516SKenneth E. Jansencc      fname1='number of global modes?'
219*59599516SKenneth E. Jansencc      call readheader(igeomBAK,fname1,nshgt,ione,'integer', iotype)
220*59599516SKenneth E. Jansenc      fname1='number of interior elements?'
221*59599516SKenneth E. Jansenc      call readheader(igeomBAK,fname1,numel,ione,'integer', iotype)
222*59599516SKenneth E. Jansenc      fname1='number of boundary elements?'
223*59599516SKenneth E. Jansenc      call readheader(igeomBAK,fname1,numelb,ione,'integer', iotype)
224*59599516SKenneth E. Jansenc      fname1='maximum number of element nodes?'
225*59599516SKenneth E. Jansenc      call readheader(igeomBAK,fname1,nen,ione,'integer', iotype)
226*59599516SKenneth E. Jansenc      fname1='number of interior tpblocks?'
227*59599516SKenneth E. Jansenc      call readheader(igeomBAK,fname1,nelblk,ione,'integer', iotype)
228*59599516SKenneth E. Jansenc      fname1='number of boundary tpblocks?'
229*59599516SKenneth E. Jansenc      call readheader(igeomBAK,fname1,nelblb,ione,'integer', iotype)
230*59599516SKenneth E. Jansenc      fname1='number of nodes with Dirichlet BCs?'
231*59599516SKenneth E. Jansenc      call readheader(igeomBAK,fname1,numpbc,ione,'integer', iotype)
232*59599516SKenneth E. Jansenc      fname1='number of shape functions?'
233*59599516SKenneth E. Jansenc      call readheader(igeomBAK,fname1,ntopsh,ione,'integer', iotype)
234*59599516SKenneth E. Jansen
235*59599516SKenneth E. Jansenc
236*59599516SKenneth E. Jansenc.... calculate the maximum number of boundary element nodes
237*59599516SKenneth E. Jansenc
238*59599516SKenneth E. Jansen      nenb = 0
239*59599516SKenneth E. Jansen      do i = 1, melCat
240*59599516SKenneth E. Jansen         if (nen .eq. nenCat(i,nsd)) nenb = max(nenCat(i,nsd-1), nenb)
241*59599516SKenneth E. Jansen      enddo
242*59599516SKenneth E. Jansenc
243*59599516SKenneth E. Jansen      if (myrank == master) then
244*59599516SKenneth E. Jansen         if (nenb .eq. 0) call error ('input   ','nen     ',nen)
245*59599516SKenneth E. Jansen      endif
246*59599516SKenneth E. Jansenc
247*59599516SKenneth E. Jansenc.... setup some useful constants
248*59599516SKenneth E. Jansenc
249*59599516SKenneth E. Jansen      I3nsd  = nsd / 3          ! nsd=3 integer flag
250*59599516SKenneth E. Jansen      E3nsd  = float(I3nsd)     ! nsd=3 real    flag
251*59599516SKenneth E. Jansenc
252*59599516SKenneth E. Jansen      if(matflg(1,1).lt.0) then
253*59599516SKenneth E. Jansen         nflow = nsd + 1
254*59599516SKenneth E. Jansen      else
255*59599516SKenneth E. Jansen         nflow = nsd + 2
256*59599516SKenneth E. Jansen      endif
257*59599516SKenneth E. Jansen      ndof   = nsd + 2
258*59599516SKenneth E. Jansen      nsclr=impl(1)/100
259*59599516SKenneth E. Jansen      ndof=ndof+nsclr           ! number of sclr transport equations to solve
260*59599516SKenneth E. Jansen
261*59599516SKenneth E. Jansen      ndofBC = ndof + I3nsd     ! dimension of BC array
262*59599516SKenneth E. Jansen      ndiBCB = 2                ! dimension of iBCB array
263*59599516SKenneth E. Jansen      ndBCB  = ndof + 1         ! dimension of BCB array
264*59599516SKenneth E. Jansenc
265*59599516SKenneth E. Jansen      nsymdf = (ndof*(ndof + 1)) / 2 ! symm. d.o.f.'s
266*59599516SKenneth E. Jansenc
267*59599516SKenneth E. Jansenc.... ----------------------> Communication tasks <--------------------
268*59599516SKenneth E. Jansenc
269*59599516SKenneth E. Jansen
270*59599516SKenneth E. Jansencc still read in new
271*59599516SKenneth E. Jansen
272*59599516SKenneth E. Jansen      if(numpe > 1) then
273*59599516SKenneth E. Jansen
274*59599516SKenneth E. Jansencc MR CHANGE
275*59599516SKenneth E. Jansen         write (temp1,"('(''size of ilwork array@'',i',i1,',A1)')") itmp
276*59599516SKenneth E. Jansen         write (fname2,temp1) (myrank+1),'?'
277*59599516SKenneth E. Jansen         call readheader(igeom,fname2 // char(0),nlwork,ione,
278*59599516SKenneth E. Jansen     &   'integer' // char(0) ,iotype)
279*59599516SKenneth E. Jansen
280*59599516SKenneth E. Jansen         write (temp1,"('(''ilwork@'',i',i1,',A1)')") itmp
281*59599516SKenneth E. Jansen         write (fname2,temp1) (myrank+1),'?'
282*59599516SKenneth E. Jansen         call readheader(igeom,fname2 //char(0) ,nlwork,ione,
283*59599516SKenneth E. Jansen     &   'integer' // char(0) , iotype)
284*59599516SKenneth E. Jansen
285*59599516SKenneth E. Jansen         allocate( point2ilwork(nlwork) )
286*59599516SKenneth E. Jansen         allocate( ilworkread(nlwork) )
287*59599516SKenneth E. Jansen         call readdatablock(igeom,fname2 // char(0),ilworkread,
288*59599516SKenneth E. Jansen     &                      nlwork,'integer' // char(0) , iotype)
289*59599516SKenneth E. Jansen
290*59599516SKenneth E. Jansenc      call closefile( igeom, "read" )
291*59599516SKenneth E. Jansenc      call finalizephmpiio( igeom )
292*59599516SKenneth E. Jansen
293*59599516SKenneth E. Jansenc         fname1='size of ilwork array?'
294*59599516SKenneth E. Jansenc         call readheader(igeomBAK,fname1,nlwork,ione,'integer', iotype)
295*59599516SKenneth E. Jansen
296*59599516SKenneth E. Jansenc         ione=1
297*59599516SKenneth E. Jansenc         fname1='ilwork?'
298*59599516SKenneth E. Jansenc         call readheader(igeomBAK,fname1,nlwork,ione,'integer', iotype)
299*59599516SKenneth E. Jansen
300*59599516SKenneth E. Jansenc         allocate( point2ilwork(nlwork) )
301*59599516SKenneth E. Jansenc         allocate( ilworkread(nlwork) )
302*59599516SKenneth E. Jansenc         call readdatablock(igeomBAK,fname1,ilworkread,
303*59599516SKenneth E. Jansenc     &                      nlwork,'integer', iotype)
304*59599516SKenneth E. Jansencc MR CHANGE
305*59599516SKenneth E. Jansen
306*59599516SKenneth E. Jansen
307*59599516SKenneth E. Jansen         point2ilwork = ilworkread
308*59599516SKenneth E. Jansen         call ctypes (point2ilwork)
309*59599516SKenneth E. Jansen      else
310*59599516SKenneth E. Jansen           nlwork=1
311*59599516SKenneth E. Jansen           allocate( point2ilwork(1))
312*59599516SKenneth E. Jansen           nshg0 = nshg
313*59599516SKenneth E. Jansen      endif
314*59599516SKenneth E. Jansen
315*59599516SKenneth E. Jansencccccccccccccccccccccccccccccccccccccccccc
316*59599516SKenneth E. Jansen
317*59599516SKenneth E. Jansen      itwo=2
318*59599516SKenneth E. Jansen      write (temp1,"('(''co-ordinates@'',i',i1,',A1)')") itmp
319*59599516SKenneth E. Jansen      write (fname2,temp1) (myrank+1),'?'
320*59599516SKenneth E. Jansen
321*59599516SKenneth E. Jansenc      print *, "fname2 is", fname2
322*59599516SKenneth E. Jansen
323*59599516SKenneth E. Jansencc MR CHANGE
324*59599516SKenneth E. Jansenc      call initphmpiio(nfields,nppf,nfiles,igeom,'read')
325*59599516SKenneth E. Jansenc      call openfile( fnamer, 'read', igeom )
326*59599516SKenneth E. JansenCC MR CHANGE
327*59599516SKenneth E. Jansen
328*59599516SKenneth E. Jansen      call readheader(igeom,fname2 // char(0),intfromfile,itwo,
329*59599516SKenneth E. Jansen     & 'double' // char(0), iotype)
330*59599516SKenneth E. Jansen      numnp=intfromfile(1)
331*59599516SKenneth E. Jansenc      print *, "read out @@@@@@ is ", numnp
332*59599516SKenneth E. Jansen      allocate( point2x(numnp,nsd) )
333*59599516SKenneth E. Jansen      allocate( xread(numnp,nsd) )
334*59599516SKenneth E. Jansen      ixsiz=numnp*nsd
335*59599516SKenneth E. Jansen      call readdatablock(igeom,fname2 // char(0),xread,ixsiz,
336*59599516SKenneth E. Jansen     & 'double' // char(0), iotype)
337*59599516SKenneth E. Jansen      point2x = xread
338*59599516SKenneth E. Jansen
339*59599516SKenneth E. Jansen!      call closefile( igeom, "read" )
340*59599516SKenneth E. Jansen!      call finalizephmpiio( igeom )
341*59599516SKenneth E. Jansen
342*59599516SKenneth E. Jansen!       print *, "Finished finalize"
343*59599516SKenneth E. Jansen
344*59599516SKenneth E. Jansenc      deallocate (point2x)
345*59599516SKenneth E. Jansenc      deallocate (xread)
346*59599516SKenneth E. Jansen
347*59599516SKenneth E. Jansencccccccccccccccccccccccccccccccccccccccccc
348*59599516SKenneth E. Jansen
349*59599516SKenneth E. Jansenc
350*59599516SKenneth E. Jansenc.... read the node coordinates
351*59599516SKenneth E. Jansenc
352*59599516SKenneth E. Jansen
353*59599516SKenneth E. Jansenc      itwo=2
354*59599516SKenneth E. Jansenc      fname1='co-ordinates?'
355*59599516SKenneth E. Jansenc      call readheader(igeomBAK,fname1,intfromfile,itwo, 'double', iotype)
356*59599516SKenneth E. Jansenc      numnp=intfromfile(1)
357*59599516SKenneth E. Jansencc      nsd=intfromfile(2)
358*59599516SKenneth E. Jansenc      allocate( point2x(numnp,nsd) )
359*59599516SKenneth E. Jansenc      allocate( xread(numnp,nsd) )
360*59599516SKenneth E. Jansenc      ixsiz=numnp*nsd
361*59599516SKenneth E. Jansenc      call readdatablock(igeomBAK,fname1,xread,ixsiz, 'double',iotype)
362*59599516SKenneth E. Jansenc      point2x = xread
363*59599516SKenneth E. Jansen
364*59599516SKenneth E. Jansenc
365*59599516SKenneth E. Jansenc.... read in and block out the connectivity
366*59599516SKenneth E. Jansenc
367*59599516SKenneth E. Jansen
368*59599516SKenneth E. Jansen! !MR CHANGE
369*59599516SKenneth E. Jansen!     This is not necessary but this avoids to have the geombc files opend two times.
370*59599516SKenneth E. Jansen!     A better way consists in pasisng the file handler to genblk or make it global or use igeomBAK instead of igeom
371*59599516SKenneth E. Jansen!      call closefile( igeom, "read" )
372*59599516SKenneth E. Jansen!      call finalizephmpiio( igeom )
373*59599516SKenneth E. Jansen! !MR CHANGE END
374*59599516SKenneth E. Jansen
375*59599516SKenneth E. Jansen      call genblk (IBKSIZ)
376*59599516SKenneth E. Jansen
377*59599516SKenneth E. Jansen! !MR CHANGE
378*59599516SKenneth E. Jansen!      call initphmpiio( nfields, nppf, nfiles, igeom, 'read')
379*59599516SKenneth E. Jansen!      call openfile( fnamer, 'read', igeom )
380*59599516SKenneth E. Jansen! !MR CHANGE END
381*59599516SKenneth E. Jansen
382*59599516SKenneth E. Jansenc
383*59599516SKenneth E. Jansenc.... read the boundary condition mapping array
384*59599516SKenneth E. Jansenc
385*59599516SKenneth E. Jansen
386*59599516SKenneth E. Jansencc MR CHANGE
387*59599516SKenneth E. Jansen!      call initphmpiio(nfields,nppf,nfiles,igeom, 'read')
388*59599516SKenneth E. Jansen!      call openfile( fnamer, 'read', igeom )
389*59599516SKenneth E. Jansencc MR CHANGE
390*59599516SKenneth E. Jansen
391*59599516SKenneth E. Jansen      ione=1
392*59599516SKenneth E. Jansen      write (temp1,"('(''bc mapping array@'',i',i1,',A1)')") itmp
393*59599516SKenneth E. Jansen      write (fname2,temp1) (myrank+1),'?'
394*59599516SKenneth E. Jansen      call readheader(igeom,fname2 // char(0),nshg,ione,
395*59599516SKenneth E. Jansen     & 'integer' // char(0),iotype)
396*59599516SKenneth E. Jansen
397*59599516SKenneth E. Jansenc      fname1='bc mapping array?'
398*59599516SKenneth E. Jansenc      call readheader(igeomBAK,fname1,nshg,
399*59599516SKenneth E. Jansenc     &     ione,'integer', iotype)
400*59599516SKenneth E. Jansen
401*59599516SKenneth E. Jansen      allocate( nBC(nshg) )
402*59599516SKenneth E. Jansen
403*59599516SKenneth E. Jansen      allocate( nBCread(nshg) )
404*59599516SKenneth E. Jansen
405*59599516SKenneth E. Jansenc      call readdatablock(igeomBAK,fname1,nBCread,nshg,'integer',iotype)
406*59599516SKenneth E. Jansen      call readdatablock(igeom,fname2 // char(0),nBCread,nshg,
407*59599516SKenneth E. Jansen     & 'integer' // char(0),iotype)
408*59599516SKenneth E. Jansen
409*59599516SKenneth E. Jansen      nBC=nBCread
410*59599516SKenneth E. Jansen
411*59599516SKenneth E. Jansenc
412*59599516SKenneth E. Jansenc.... read the temporary iBC array
413*59599516SKenneth E. Jansenc
414*59599516SKenneth E. Jansen      ione=1
415*59599516SKenneth E. Jansen      write (temp1,"('(''bc codes array@'',i',i1,',A1)')") itmp
416*59599516SKenneth E. Jansen      write (fname2,temp1) (myrank+1),'?'
417*59599516SKenneth E. Jansen      call readheader(igeom,fname2 // char(0) ,numpbc,ione,
418*59599516SKenneth E. Jansen     & 'integer' // char(0),iotype)
419*59599516SKenneth E. Jansen
420*59599516SKenneth E. Jansenc      ione = 1
421*59599516SKenneth E. Jansenc      fname1='bc codes array?'
422*59599516SKenneth E. Jansenc      call readheader(igeomBAK,fname1,numpbc,
423*59599516SKenneth E. Jansenc     &     ione, 'integer', iotype)
424*59599516SKenneth E. Jansen
425*59599516SKenneth E. Jansen!MR CHANGE
426*59599516SKenneth E. Jansen!       if ( numpbc > 0 ) then
427*59599516SKenneth E. Jansen!          allocate( iBCtmp(numpbc) )
428*59599516SKenneth E. Jansen!          allocate( iBCtmpread(numpbc) )
429*59599516SKenneth E. Jansen! c         call readdatablock(igeomBAK,fname1,iBCtmpread,numpbc,
430*59599516SKenneth E. Jansen! c     &                      'integer',iotype)
431*59599516SKenneth E. Jansen!         call readdatablock(igeom,fname2,iBCtmpread,numpbc,
432*59599516SKenneth E. Jansen!      &                      'integer',iotype)
433*59599516SKenneth E. Jansen!          iBCtmp=iBCtmpread
434*59599516SKenneth E. Jansen!       else  ! sometimes a partition has no BC's
435*59599516SKenneth E. Jansen!          allocate( iBCtmp(1) )
436*59599516SKenneth E. Jansen!          iBCtmp=0
437*59599516SKenneth E. Jansen!       endif
438*59599516SKenneth E. Jansen
439*59599516SKenneth E. Jansen      if ( numpbc > 0 ) then
440*59599516SKenneth E. Jansen        allocate( iBCtmp(numpbc) )
441*59599516SKenneth E. Jansen        allocate( iBCtmpread(numpbc) )
442*59599516SKenneth E. Jansen      else
443*59599516SKenneth E. Jansen        allocate( iBCtmp(1) )
444*59599516SKenneth E. Jansen        allocate( iBCtmpread(1) )
445*59599516SKenneth E. Jansen      endif
446*59599516SKenneth E. Jansenc         call readdatablock(igeomBAK,fname1,iBCtmpread,numpbc,
447*59599516SKenneth E. Jansenc     &                      'integer',iotype)
448*59599516SKenneth E. Jansen      call readdatablock(igeom,fname2 // char(0),iBCtmpread,numpbc,
449*59599516SKenneth E. Jansen     &  'integer' // char(0),iotype)
450*59599516SKenneth E. Jansen
451*59599516SKenneth E. Jansen      if ( numpbc > 0 ) then
452*59599516SKenneth E. Jansen         iBCtmp=iBCtmpread
453*59599516SKenneth E. Jansen      else  ! sometimes a partition has no BC's
454*59599516SKenneth E. Jansen         deallocate( iBCtmpread)
455*59599516SKenneth E. Jansen         iBCtmp=0
456*59599516SKenneth E. Jansen      endif
457*59599516SKenneth E. Jansen!MR CHANGE END
458*59599516SKenneth E. Jansen
459*59599516SKenneth E. Jansenc
460*59599516SKenneth E. Jansenc.... read boundary condition data
461*59599516SKenneth E. Jansenc
462*59599516SKenneth E. Jansen
463*59599516SKenneth E. Jansen      ione=1
464*59599516SKenneth E. Jansen      write (temp1,"('(''boundary condition array@'',i',i1,',A1)')")
465*59599516SKenneth E. Jansen     &     itmp
466*59599516SKenneth E. Jansen      write (fname2,temp1) (myrank+1),'?'
467*59599516SKenneth E. Jansen
468*59599516SKenneth E. Jansenc      ione=1
469*59599516SKenneth E. Jansenc      fname1='boundary condition array?'
470*59599516SKenneth E. Jansenc      call readheader(igeomBAK,fname1,intfromfile,
471*59599516SKenneth E. Jansenc     &     ione, 'double', iotype)
472*59599516SKenneth E. Jansen      call readheader(igeom,fname2 // char(0),intfromfile,
473*59599516SKenneth E. Jansen     &     ione, 'double' // char(0), iotype)
474*59599516SKenneth E. Jansenc here intfromfile(1) contains (ndof+7)*numpbc
475*59599516SKenneth E. Jansen!MR CHANGE
476*59599516SKenneth E. Jansen!       if ( numpbc > 0 ) then
477*59599516SKenneth E. Jansen!          if(intfromfile(1).ne.(ndof+7)*numpbc) then
478*59599516SKenneth E. Jansen!            warning='WARNING more data in BCinp than needed: keeping 1st'
479*59599516SKenneth E. Jansen!            write(*,*) warning, ndof+7
480*59599516SKenneth E. Jansen!          endif
481*59599516SKenneth E. Jansen!          allocate( BCinp(numpbc,ndof+7) )
482*59599516SKenneth E. Jansen!          nsecondrank=intfromfile(1)/numpbc
483*59599516SKenneth E. Jansen!          allocate( BCinpread(numpbc,nsecondrank) )
484*59599516SKenneth E. Jansen!          iBCinpsiz=intfromfile(1)
485*59599516SKenneth E. Jansen! c         call readdatablock(igeomBAK,fname1,BCinpread,iBCinpsiz,
486*59599516SKenneth E. Jansen! c     &                      'double',iotype)
487*59599516SKenneth E. Jansen!          call readdatablock(igeom,fname2,BCinpread,iBCinpsiz,
488*59599516SKenneth E. Jansen!      &                      'double',iotype)
489*59599516SKenneth E. Jansen!          BCinp(:,1:(ndof+7))=BCinpread(:,1:(ndof+7))
490*59599516SKenneth E. Jansen!       else  ! sometimes a partition has no BC's
491*59599516SKenneth E. Jansen!          allocate( BCinp(1,ndof+7) )
492*59599516SKenneth E. Jansen!          BCinp=0
493*59599516SKenneth E. Jansen!       endif
494*59599516SKenneth E. Jansen
495*59599516SKenneth E. Jansen      if ( numpbc > 0 ) then
496*59599516SKenneth E. Jansen!         if(intfromfile(1).ne.(ndof+7)*numpbc) then
497*59599516SKenneth E. Jansen!           warning='WARNING more data in BCinp than needed: keeping 1st'
498*59599516SKenneth E. Jansen!           write(*,*) warning, ndof+7
499*59599516SKenneth E. Jansen!         endif
500*59599516SKenneth E. Jansen         allocate( BCinp(numpbc,ndof+7) )
501*59599516SKenneth E. Jansen         nsecondrank=intfromfile(1)/numpbc
502*59599516SKenneth E. Jansen         allocate( BCinpread(numpbc,nsecondrank) )
503*59599516SKenneth E. Jansen         iBCinpsiz=intfromfile(1)
504*59599516SKenneth E. Jansen      else
505*59599516SKenneth E. Jansen         allocate( BCinp(1,ndof+7) )
506*59599516SKenneth E. Jansen         allocate( BCinpread(0,0) ) !dummy
507*59599516SKenneth E. Jansen         iBCinpsiz=intfromfile(1)
508*59599516SKenneth E. Jansen      endif
509*59599516SKenneth E. Jansenc         call readdatablock(igeomBAK,fname1,BCinpread,iBCinpsiz,
510*59599516SKenneth E. Jansenc     &                      'double',iotype)
511*59599516SKenneth E. Jansen
512*59599516SKenneth E. Jansen      call readdatablock(igeom,fname2 // char(0),BCinpread,iBCinpsiz,
513*59599516SKenneth E. Jansen     &                      'double' // char(0) ,iotype)
514*59599516SKenneth E. Jansen
515*59599516SKenneth E. Jansen      if ( numpbc > 0 ) then
516*59599516SKenneth E. Jansen         BCinp(:,1:(ndof+7))=BCinpread(:,1:(ndof+7))
517*59599516SKenneth E. Jansen      else  ! sometimes a partition has no BC's
518*59599516SKenneth E. Jansen         deallocate(BCinpread)
519*59599516SKenneth E. Jansen         BCinp=0
520*59599516SKenneth E. Jansen      endif
521*59599516SKenneth E. Jansen!MR CHANGE END
522*59599516SKenneth E. Jansen
523*59599516SKenneth E. Jansenc
524*59599516SKenneth E. Jansenc.... read periodic boundary conditions
525*59599516SKenneth E. Jansenc
526*59599516SKenneth E. Jansen
527*59599516SKenneth E. Jansen      ione=1
528*59599516SKenneth E. Jansen      write (temp1,"('(''periodic masters array@'',i',i1,',A1)')") itmp
529*59599516SKenneth E. Jansen      write (fname2,temp1) (myrank+1),'?'
530*59599516SKenneth E. Jansen
531*59599516SKenneth E. Jansenc      fname1='periodic masters array?'
532*59599516SKenneth E. Jansenc      call readheader(igeomBAK,fname1,nshg,
533*59599516SKenneth E. Jansenc     &     ione, 'integer', iotype)
534*59599516SKenneth E. Jansen      call readheader(igeom,fname2 // char(0) ,nshg,
535*59599516SKenneth E. Jansen     &     ione, 'integer' // char(0), iotype)
536*59599516SKenneth E. Jansen      allocate( point2iper(nshg) )
537*59599516SKenneth E. Jansen      allocate( iperread(nshg) )
538*59599516SKenneth E. Jansenc      call readdatablock(igeomBAK,fname1,iperread,nshg,
539*59599516SKenneth E. Jansenc     &                      'integer',iotype)
540*59599516SKenneth E. Jansen      call readdatablock(igeom,fname2 // char(0),iperread,nshg,
541*59599516SKenneth E. Jansen     &                      'integer' // char(0),iotype)
542*59599516SKenneth E. Jansen      point2iper=iperread
543*59599516SKenneth E. Jansen
544*59599516SKenneth E. Jansen
545*59599516SKenneth E. Jansen! !MR CHANGE
546*59599516SKenneth E. Jansen!      call closefile( igeom, "read" )
547*59599516SKenneth E. Jansen!      call finalizephmpiio( igeom )
548*59599516SKenneth E. Jansen! !MR CHANGE END
549*59599516SKenneth E. Jansen
550*59599516SKenneth E. Jansenc
551*59599516SKenneth E. Jansenc.... generate the boundary element blocks
552*59599516SKenneth E. Jansenc
553*59599516SKenneth E. Jansen      call genbkb (ibksiz)
554*59599516SKenneth E. Jansen
555*59599516SKenneth E. Jansen
556*59599516SKenneth E. Jansen! !MR CHANGE
557*59599516SKenneth E. Jansen!       write(*,*) 'HELLO 12 from ', myrank
558*59599516SKenneth E. Jansen! !MR CHANGE END
559*59599516SKenneth E. Jansen
560*59599516SKenneth E. Jansenc
561*59599516SKenneth E. Jansenc  Read in the nsons and ifath arrays if needed
562*59599516SKenneth E. Jansenc
563*59599516SKenneth E. Jansenc  There is a fundamental shift in the meaning of ifath based on whether
564*59599516SKenneth E. Jansenc  there exist homogenous directions in the flow.
565*59599516SKenneth E. Jansenc
566*59599516SKenneth E. Jansenc  HOMOGENOUS DIRECTIONS EXIST:  Here nfath is the number of inhomogenous
567*59599516SKenneth E. Jansenc  points in the TOTAL mesh.  That is to say that each partition keeps a
568*59599516SKenneth E. Jansenc  link to  ALL inhomogenous points.  This link is furthermore not to the
569*59599516SKenneth E. Jansenc  sms numbering but to the original structured grid numbering.  These
570*59599516SKenneth E. Jansenc  inhomogenous points are thought of as fathers, with their sons being all
571*59599516SKenneth E. Jansenc  the points in the homogenous directions that have this father's
572*59599516SKenneth E. Jansenc  inhomogeneity.  The array ifath takes as an arguement the sms numbering
573*59599516SKenneth E. Jansenc  and returns as a result the father.
574*59599516SKenneth E. Jansenc
575*59599516SKenneth E. Jansenc  In this case nsons is the number of sons that each father has and ifath
576*59599516SKenneth E. Jansenc  is an array which tells the
577*59599516SKenneth E. Jansenc
578*59599516SKenneth E. Jansenc  NO HOMOGENOUS DIRECTIONS.  In this case the mesh would grow to rapidly
579*59599516SKenneth E. Jansenc  if we followed the above strategy since every partition would index its
580*59599516SKenneth E. Jansenc  points to the ENTIRE mesh.  Furthermore, there would never be a need
581*59599516SKenneth E. Jansenc  to average to a node off processor since there is no spatial averaging.
582*59599516SKenneth E. Jansenc  Therefore, to properly account for this case we must recognize it and
583*59599516SKenneth E. Jansenc  inerrupt certain actions (i.e. assembly of the average across partitions).
584*59599516SKenneth E. Jansenc  This case is easily identified by noting that maxval(nsons) =1 (i.e. no
585*59599516SKenneth E. Jansenc  father has any sons).  Reiterating to be clear, in this case ifath does
586*59599516SKenneth E. Jansenc  not point to a global numbering but instead just points to itself.
587*59599516SKenneth E. Jansenc
588*59599516SKenneth E. Jansen      nfath=1  ! some architectures choke on a zero or undeclared
589*59599516SKenneth E. Jansen                 ! dimension variable.  This sets it to a safe, small value.
590*59599516SKenneth E. Jansen      if(((iLES .lt. 20) .and. (iLES.gt.0))
591*59599516SKenneth E. Jansen     &                   .or. (itwmod.gt.0)  ) then ! don't forget same
592*59599516SKenneth E. Jansen                                                    ! conditional in proces.f
593*59599516SKenneth E. Jansen
594*59599516SKenneth E. Jansenc           read (igeomBAK) nfath  ! nfath already read in input.f,
595*59599516SKenneth E. Jansen                                     ! needed for alloc
596*59599516SKenneth E. Jansen         ione=1
597*59599516SKenneth E. Jansenc         call creadlist(igeomBAK,ione,nfath)
598*59599516SKenneth E. Jansenc         fname1='keyword sonfath?'
599*59599516SKenneth E. Jansen         if(nohomog.gt.0) then
600*59599516SKenneth E. Jansen
601*59599516SKenneth E. Jansen
602*59599516SKenneth E. Jansen!             fname1='number of father-nodes?'
603*59599516SKenneth E. Jansen!             call readheader(igeomBAK,fname1,nfath,ione,'integer', iotype)
604*59599516SKenneth E. Jansen
605*59599516SKenneth E. Jansen            write (temp1,"('(''number of father-nodes@'',i',i1,',A1)')")
606*59599516SKenneth E. Jansen     &             itmp
607*59599516SKenneth E. Jansen            write (fname2,temp1) (myrank+1),'?'
608*59599516SKenneth E. Jansen            call readheader(igeom,fname2 // char(0),nfath,ione,
609*59599516SKenneth E. Jansen     &      'integer' // char(0), iotype)
610*59599516SKenneth E. Jansen
611*59599516SKenneth E. Jansenc
612*59599516SKenneth E. Jansenc     fname1='keyword nsons?'
613*59599516SKenneth E. Jansen!             fname1='number of son-nodes for each father?'
614*59599516SKenneth E. Jansen!             call readheader(igeomBAK,fname1,nfath,ione,'integer', iotype)
615*59599516SKenneth E. Jansen
616*59599516SKenneth E. Jansen            write (temp1,
617*59599516SKenneth E. Jansen     & "('(''number of son-nodes for each father@'',i',i1,',A1)')") itmp
618*59599516SKenneth E. Jansen            write (fname2,temp1) (myrank+1),'?'
619*59599516SKenneth E. Jansen            call readheader(igeom,fname2 // char(0),nfath,ione,
620*59599516SKenneth E. Jansen     &      'integer' // char(0), iotype)
621*59599516SKenneth E. Jansen
622*59599516SKenneth E. Jansen            allocate (point2nsons(nfath))
623*59599516SKenneth E. Jansen
624*59599516SKenneth E. Jansen!             call readdatablock(igeomBAK,fname1,point2nsons,nfath,
625*59599516SKenneth E. Jansen!      &                      'integer',iotype)
626*59599516SKenneth E. Jansen            call readdatablock(igeom,fname2 // char(0),point2nsons,
627*59599516SKenneth E. Jansen     &                      nfath,'integer' // char(0), iotype)
628*59599516SKenneth E. Jansen
629*59599516SKenneth E. Jansenc
630*59599516SKenneth E. Jansen!             fname1='keyword ifath?'
631*59599516SKenneth E. Jansen!             call readheader(igeomBAK,fname1,nshg,ione,'integer', iotype)
632*59599516SKenneth E. Jansen
633*59599516SKenneth E. Jansen            write (temp1,"('(''keyword ifath@'',i',i1,',A1)')") itmp
634*59599516SKenneth E. Jansen            write (fname2,temp1) (myrank+1),'?'
635*59599516SKenneth E. Jansen            call readheader(igeom,fname2 // char(0),nshg,ione,
636*59599516SKenneth E. Jansen     &      'integer' // char(0), iotype)
637*59599516SKenneth E. Jansen
638*59599516SKenneth E. Jansen            allocate (point2ifath(nshg))
639*59599516SKenneth E. Jansen
640*59599516SKenneth E. Jansen!             call readdatablock(igeomBAK,fname1,point2ifath,nshg,
641*59599516SKenneth E. Jansen!      &                      'integer',iotype)
642*59599516SKenneth E. Jansen            call readdatablock(igeom,fname2 // char(0),point2ifath,
643*59599516SKenneth E. Jansen     &                      nshg,'integer' // char(0) , iotype)
644*59599516SKenneth E. Jansen
645*59599516SKenneth E. Jansenc
646*59599516SKenneth E. Jansen            nsonmax=maxval(point2nsons)
647*59599516SKenneth E. Jansenc
648*59599516SKenneth E. Jansen         else  ! this is the case where there is no homogeneity
649*59599516SKenneth E. Jansen               ! therefore ever node is a father (too itself).  sonfath
650*59599516SKenneth E. Jansen               ! (a routine in NSpre) will set this up but this gives
651*59599516SKenneth E. Jansen               ! you an option to avoid that.
652*59599516SKenneth E. Jansen            nfath=nshg
653*59599516SKenneth E. Jansen            allocate (point2nsons(nfath))
654*59599516SKenneth E. Jansen            point2nsons=1
655*59599516SKenneth E. Jansen            allocate (point2ifath(nshg))
656*59599516SKenneth E. Jansen            do i=1,nshg
657*59599516SKenneth E. Jansen               point2ifath(i)=i
658*59599516SKenneth E. Jansen            enddo
659*59599516SKenneth E. Jansen            nsonmax=1
660*59599516SKenneth E. Jansenc
661*59599516SKenneth E. Jansen         endif
662*59599516SKenneth E. Jansen      else
663*59599516SKenneth E. Jansen         allocate (point2nsons(1))
664*59599516SKenneth E. Jansen         allocate (point2ifath(1))
665*59599516SKenneth E. Jansen      endif
666*59599516SKenneth E. Jansen
667*59599516SKenneth E. Jansen! !MR CHANGE
668*59599516SKenneth E. Jansen      call closefile( igeom, "read" // char(0) )
669*59599516SKenneth E. Jansen      call finalizephmpiio( igeom )
670*59599516SKenneth E. Jansen! !MR CHANGE END
671*59599516SKenneth E. Jansen
672*59599516SKenneth E. Jansen! !MR CHANGE
673*59599516SKenneth E. Jansen!       write(*,*) 'HELLO 13 from ', myrank
674*59599516SKenneth E. Jansen! !MR CHANGE END
675*59599516SKenneth E. Jansen
676*59599516SKenneth E. Jansenc
677*59599516SKenneth E. Jansenc  renumber the master partition for SPEBC
678*59599516SKenneth E. Jansenc
679*59599516SKenneth E. Jansenc      if((myrank.eq.master).and.(irscale.ge.0)) then
680*59599516SKenneth E. Jansenc         call setSPEBC(numnp, nfath, nsonmax)
681*59599516SKenneth E. Jansenc         call renum(point2x,point2ifath,point2nsons)
682*59599516SKenneth E. Jansenc      endif
683*59599516SKenneth E. Jansenc
684*59599516SKenneth E. Jansenc.... Read restart files
685*59599516SKenneth E. Jansen
686*59599516SKenneth E. Jansenc$$$ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
687*59599516SKenneth E. Jansenc
688*59599516SKenneth E. Jansenc      nfields = 11
689*59599516SKenneth E. Jansenc      numparts = 512
690*59599516SKenneth E. Jansenc      nppp = numparts/numpe
691*59599516SKenneth E. Jansenc      startpart = myrank * nppp +1
692*59599516SKenneth E. Jansenc      int endpart = startpart + nppp - 1
693*59599516SKenneth E. Jansenc      nppf = numparts/nfiles
694*59599516SKenneth E. Jansencc      fnamer = "/users/nliu/PIG4/4-procs_case/restart-file"
695*59599516SKenneth E. Jansencc      fnamer="./4-procs_case/restart-file"
696*59599516SKenneth E. Jansen
697*59599516SKenneth E. Jansen      itmp=1
698*59599516SKenneth E. Jansen      if (irstart .gt. 0) itmp = int(log10(float(irstart+1)))+1
699*59599516SKenneth E. Jansen
700*59599516SKenneth E. Jansen      write (fmt1,"('(''restart-dat.'',i',i1,',1x)')") itmp
701*59599516SKenneth E. Jansen
702*59599516SKenneth E. Jansen      write (fnamer,fmt1) irstart
703*59599516SKenneth E. Jansen      fnamer = trim(fnamer) // cname2(color+1)
704*59599516SKenneth E. Jansen
705*59599516SKenneth E. Jansen      call queryphmpiio(fnamer // char(0), nfields, nppf);
706*59599516SKenneth E. Jansen      if (myrank == 0) then
707*59599516SKenneth E. Jansen        write(*,*) 'Number of fields in restart-dat: ',nfields
708*59599516SKenneth E. Jansen        write(*,*) 'Number of parts per file restart-dat: ',nppf
709*59599516SKenneth E. Jansen      endif
710*59599516SKenneth E. Jansen      call initphmpiio(nfields,nppf,nfiles,descriptor,
711*59599516SKenneth E. Jansen     & 'read' // char(0))
712*59599516SKenneth E. Jansen      call openfile( fnamer // char(0) ,
713*59599516SKenneth E. Jansen     & 'read' // char(0), descriptor )
714*59599516SKenneth E. Jansen
715*59599516SKenneth E. Jansen      ithree=3
716*59599516SKenneth E. Jansenc      call creadlist(irstin,ithree,nshg2,ndof2,lstep)
717*59599516SKenneth E. Jansen
718*59599516SKenneth E. Jansen      itmp = int(log10(float(myrank+1)))+1
719*59599516SKenneth E. Jansen      write (temp1,"('(''solution@'',i',i1,',A1)')") itmp
720*59599516SKenneth E. Jansen      write (fname1,temp1) (myrank+1),'?'
721*59599516SKenneth E. Jansen      fname1 = trim(fname1)
722*59599516SKenneth E. Jansen
723*59599516SKenneth E. Jansenc      print *, "Solution is : ", fname1
724*59599516SKenneth E. Jansen
725*59599516SKenneth E. Jansen      intfromfile=0
726*59599516SKenneth E. Jansen      call readheader(descriptor,fname1 // char(0) ,intfromfile,
727*59599516SKenneth E. Jansen     & ithree,'integer' // char(0),
728*59599516SKenneth E. Jansen     &                                                         iotype)
729*59599516SKenneth E. Jansenc
730*59599516SKenneth E. Jansenc.... read the values of primitive variables into q
731*59599516SKenneth E. Jansenc
732*59599516SKenneth E. Jansen
733*59599516SKenneth E. Jansenc      print *, "intfromfile(1) is ", intfromfile(1)
734*59599516SKenneth E. Jansenc      print *, "intfromfile(2) is ", intfromfile(2)
735*59599516SKenneth E. Jansenc      print *, "intfromfile(3) is ", intfromfile(3)
736*59599516SKenneth E. Jansen
737*59599516SKenneth E. Jansen      allocate( qold(nshg,ndof) )
738*59599516SKenneth E. Jansen      if(intfromfile(1).ne.0) then
739*59599516SKenneth E. Jansen         nshg2=intfromfile(1)
740*59599516SKenneth E. Jansen         ndof2=intfromfile(2)
741*59599516SKenneth E. Jansen         lstep=intfromfile(3)
742*59599516SKenneth E. Jansen         if(ndof2.ne.ndof) then
743*59599516SKenneth E. Jansen
744*59599516SKenneth E. Jansen         endif
745*59599516SKenneth E. Jansenc
746*59599516SKenneth E. Jansen        if (nshg2 .ne. nshg)
747*59599516SKenneth E. Jansen     &        call error ('restar  ', 'nshg   ', nshg)
748*59599516SKenneth E. Jansen         allocate( qread(nshg,ndof2) )
749*59599516SKenneth E. Jansen         iqsiz=nshg*ndof2
750*59599516SKenneth E. Jansen         call readdatablock(descriptor,fname1 // char(0),qread,iqsiz,
751*59599516SKenneth E. Jansen     &                         'double' // char(0),iotype)
752*59599516SKenneth E. Jansen         qold(:,1:ndof)=qread(:,1:ndof)
753*59599516SKenneth E. Jansen         deallocate(qread)
754*59599516SKenneth E. Jansen      else
755*59599516SKenneth E. Jansen         if (myrank.eq.master) then
756*59599516SKenneth E. Jansen            if (matflg(1,1).eq.0) then ! compressible
757*59599516SKenneth E. Jansen               warning='Solution is set to zero (with p and T to one)'
758*59599516SKenneth E. Jansen            else
759*59599516SKenneth E. Jansen               warning='Solution is set to zero'
760*59599516SKenneth E. Jansen            endif
761*59599516SKenneth E. Jansen            write(*,*) warning
762*59599516SKenneth E. Jansen         endif
763*59599516SKenneth E. Jansen         qold=zero
764*59599516SKenneth E. Jansen         if (matflg(1,1).eq.0) then ! compressible
765*59599516SKenneth E. Jansen            qold(:,1)=one ! avoid zero pressure
766*59599516SKenneth E. Jansen            qold(:,nflow)=one ! avoid zero temperature
767*59599516SKenneth E. Jansen         endif
768*59599516SKenneth E. Jansen      endif
769*59599516SKenneth E. Jansen
770*59599516SKenneth E. Jansen
771*59599516SKenneth E. Jansen! !MR CHANGE
772*59599516SKenneth E. Jansen!       write(*,*) 'HELLO 16-8 from ', myrank
773*59599516SKenneth E. Jansen! !MR CHANGE END
774*59599516SKenneth E. Jansen
775*59599516SKenneth E. Jansenc      itmp=1
776*59599516SKenneth E. Jansenc      if (myrank .gt. 0) itmp = int(log10(float(myrank)))+1
777*59599516SKenneth E. Jansen      write (temp1,"('(''time derivative of solution@'',i',i1,',A1)')")
778*59599516SKenneth E. Jansen     &       itmp
779*59599516SKenneth E. Jansen      write (fname1,temp1) (myrank+1),'?'
780*59599516SKenneth E. Jansen      fname1 = trim(fname1)
781*59599516SKenneth E. Jansen      intfromfile=0
782*59599516SKenneth E. Jansen      call readheader(descriptor,fname1 // char(0) ,intfromfile,
783*59599516SKenneth E. Jansen     & ithree,'integer' // char(0),
784*59599516SKenneth E. Jansen     &                iotype)
785*59599516SKenneth E. Jansen      allocate( acold(nshg,ndof) )
786*59599516SKenneth E. Jansen      if(intfromfile(1).ne.0) then
787*59599516SKenneth E. Jansen         nshg2=intfromfile(1)
788*59599516SKenneth E. Jansen         ndof2=intfromfile(2)
789*59599516SKenneth E. Jansen         lstep=intfromfile(3)
790*59599516SKenneth E. Jansen
791*59599516SKenneth E. Jansenc      print *, "intfromfile(1) is ", intfromfile(1)
792*59599516SKenneth E. Jansenc      print *, "intfromfile(2) is ", intfromfile(2)
793*59599516SKenneth E. Jansenc      print *, "intfromfile(3) is ", intfromfile(3)
794*59599516SKenneth E. Jansen
795*59599516SKenneth E. Jansen         if (nshg2 .ne. nshg)
796*59599516SKenneth E. Jansen     &        call error ('restar  ', 'nshg   ', nshg)
797*59599516SKenneth E. Jansenc
798*59599516SKenneth E. Jansen         allocate( acread(nshg,ndof2) )
799*59599516SKenneth E. Jansen         acread=zero
800*59599516SKenneth E. Jansen         iacsiz=nshg*ndof2
801*59599516SKenneth E. Jansen         call readdatablock(descriptor,fname1 // char(0),acread,
802*59599516SKenneth E. Jansen     &    iacsiz, 'double' // char(0),iotype)
803*59599516SKenneth E. Jansen         acold(:,1:ndof)=acread(:,1:ndof)
804*59599516SKenneth E. Jansen         deallocate(acread)
805*59599516SKenneth E. Jansen      else
806*59599516SKenneth E. Jansen         if (myrank.eq.master) then
807*59599516SKenneth E. Jansen            warning='Time derivative of solution is set to zero (SAFE)'
808*59599516SKenneth E. Jansen            write(*,*) warning
809*59599516SKenneth E. Jansen         endif
810*59599516SKenneth E. Jansen         acold=zero
811*59599516SKenneth E. Jansen      endif
812*59599516SKenneth E. Jansen
813*59599516SKenneth E. Jansencccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
814*59599516SKenneth E. Jansencc
815*59599516SKenneth E. Jansenc
816*59599516SKenneth E. Jansencc
817*59599516SKenneth E. Jansencc.... read the header and check it against the run data
818*59599516SKenneth E. Jansencc
819*59599516SKenneth E. Jansen
820*59599516SKenneth E. Jansen
821*59599516SKenneth E. Jansenc      ithree=3
822*59599516SKenneth E. Jansenccc      call creadlist(irstin,ithree,nshg2,ndof2,lstep)
823*59599516SKenneth E. Jansenc      fname1='solution?'
824*59599516SKenneth E. Jansenc      intfromfile=0
825*59599516SKenneth E. Jansenc      call readheader(irstin,fname1,intfromfile,
826*59599516SKenneth E. Jansenc     &     ithree,'integer', iotype)
827*59599516SKenneth E. Jansencc
828*59599516SKenneth E. Jansencc.... read the values of primitive variables into q
829*59599516SKenneth E. Jansencc
830*59599516SKenneth E. Jansenc      allocate( qold(nshg,ndof) )
831*59599516SKenneth E. Jansenc      if(intfromfile(1).ne.0) then
832*59599516SKenneth E. Jansenc         nshg2=intfromfile(1)
833*59599516SKenneth E. Jansenc         ndof2=intfromfile(2)
834*59599516SKenneth E. Jansenc         lstep=intfromfile(3)
835*59599516SKenneth E. Jansenc         if(ndof2.ne.ndof) then
836*59599516SKenneth E. Jansenc           warning='WARNING more data in restart than needed: keeping 1st '
837*59599516SKenneth E. Jansenc           write(*,*) warning , ndof
838*59599516SKenneth E. Jansenc         endif
839*59599516SKenneth E. Jansencc
840*59599516SKenneth E. Jansenc         if (nshg2 .ne. nshg)
841*59599516SKenneth E. Jansenc     &        call error ('restar  ', 'nshg   ', nshg)
842*59599516SKenneth E. Jansenc         allocate( qread(nshg,ndof2) )
843*59599516SKenneth E. Jansen
844*59599516SKenneth E. Jansenc         iqsiz=nshg*ndof2
845*59599516SKenneth E. Jansenc         call readdatablock(irstin,fname1,qread,iqsiz,
846*59599516SKenneth E. Jansenc     &                         'double',iotype)
847*59599516SKenneth E. Jansenc         qold(:,1:ndof)=qread(:,1:ndof)
848*59599516SKenneth E. Jansenc         deallocate(qread)
849*59599516SKenneth E. Jansenc      else
850*59599516SKenneth E. Jansenc         if (myrank.eq.master) then
851*59599516SKenneth E. Jansenc            if (matflg(1,1).eq.0) then ! compressible
852*59599516SKenneth E. Jansenc               warning='Solution is set to zero (with p and T to one)'
853*59599516SKenneth E. Jansenc            else
854*59599516SKenneth E. Jansenc               warning='Solution is set to zero'
855*59599516SKenneth E. Jansenc            endif
856*59599516SKenneth E. Jansenc            write(*,*) warning
857*59599516SKenneth E. Jansenc         endif
858*59599516SKenneth E. Jansenc         qold=zero
859*59599516SKenneth E. Jansenc         if (matflg(1,1).eq.0) then ! compressible
860*59599516SKenneth E. Jansenc            qold(:,1)=one ! avoid zero pressure
861*59599516SKenneth E. Jansenc            qold(:,nflow)=one ! avoid zero temperature
862*59599516SKenneth E. Jansenc         endif
863*59599516SKenneth E. Jansenc      endif
864*59599516SKenneth E. Jansencc
865*59599516SKenneth E. Jansenc      fname1='time derivative of solution?'
866*59599516SKenneth E. Jansenc      intfromfile=0
867*59599516SKenneth E. Jansenc      call readheader(irstin,fname1,intfromfile,
868*59599516SKenneth E. Jansenc     &     ithree,'integer', iotype)
869*59599516SKenneth E. Jansenc      allocate( acold(nshg,ndof) )
870*59599516SKenneth E. Jansenc      if(intfromfile(1).ne.0) then
871*59599516SKenneth E. Jansenc         nshg2=intfromfile(1)
872*59599516SKenneth E. Jansenc         ndof2=intfromfile(2)
873*59599516SKenneth E. Jansenc         lstep=intfromfile(3)
874*59599516SKenneth E. Jansenc
875*59599516SKenneth E. Jansenc         if (nshg2 .ne. nshg)
876*59599516SKenneth E. Jansenc     &        call error ('restar  ', 'nshg   ', nshg)
877*59599516SKenneth E. Jansencc
878*59599516SKenneth E. Jansenc         allocate( acread(nshg,ndof2) )
879*59599516SKenneth E. Jansenc         acread=zero
880*59599516SKenneth E. Jansenc
881*59599516SKenneth E. Jansenc         iacsiz=nshg*ndof2
882*59599516SKenneth E. Jansenc         call readdatablock(irstin,fname1,acread,iacsiz,
883*59599516SKenneth E. Jansenc     &                   'double',iotype)
884*59599516SKenneth E. Jansenc         acold(:,1:ndof)=acread(:,1:ndof)
885*59599516SKenneth E. Jansenc         deallocate(acread)
886*59599516SKenneth E. Jansenc      else
887*59599516SKenneth E. Jansenc         if (myrank.eq.master) then
888*59599516SKenneth E. Jansenc            warning='Time derivative of solution is set to zero (SAFE)'
889*59599516SKenneth E. Jansenc            write(*,*) warning
890*59599516SKenneth E. Jansenc         endif
891*59599516SKenneth E. Jansenc         acold=zero
892*59599516SKenneth E. Jansenc      endif
893*59599516SKenneth E. Jansen
894*59599516SKenneth E. Jansenc      call creadlist(irstin,ithree,nshg2,ndisp,lstep)
895*59599516SKenneth E. Jansen
896*59599516SKenneth E. Jansen      if (ideformwall.eq.1) then
897*59599516SKenneth E. Jansen!          fname1='displacement?'
898*59599516SKenneth E. Jansen!          call readheader(irstin,fname1,intfromfile,
899*59599516SKenneth E. Jansen!      &        ithree,'integer', iotype)
900*59599516SKenneth E. Jansen
901*59599516SKenneth E. Jansen          write (temp1,"('(''displacement@'',i',i1,',A1)')")
902*59599516SKenneth E. Jansen     &           itmp
903*59599516SKenneth E. Jansen          write (fname1,temp1) (myrank+1),'?'
904*59599516SKenneth E. Jansen          fname1 = trim(fname1)
905*59599516SKenneth E. Jansen          intfromfile=0
906*59599516SKenneth E. Jansen          call readheader(descriptor,fname1 // char(0),intfromfile,
907*59599516SKenneth E. Jansen     &     ithree,'integer' // char(0),iotype)
908*59599516SKenneth E. Jansen
909*59599516SKenneth E. Jansen         nshg2=intfromfile(1)
910*59599516SKenneth E. Jansen         ndisp=intfromfile(2)
911*59599516SKenneth E. Jansen         lstep=intfromfile(3)
912*59599516SKenneth E. Jansen         if(ndisp.ne.nsd) then
913*59599516SKenneth E. Jansen            warning='WARNING ndisp not equal nsd'
914*59599516SKenneth E. Jansen            write(*,*) warning , ndisp
915*59599516SKenneth E. Jansen         endif
916*59599516SKenneth E. Jansenc
917*59599516SKenneth E. Jansen         if (nshg2 .ne. nshg)
918*59599516SKenneth E. Jansen     &        call error ('restar  ', 'nshg   ', nshg)
919*59599516SKenneth E. Jansenc
920*59599516SKenneth E. Jansenc.... read the values of primitive variables into uold
921*59599516SKenneth E. Jansenc
922*59599516SKenneth E. Jansen
923*59599516SKenneth E. Jansen         allocate( uold(nshg,nsd) )
924*59599516SKenneth E. Jansen         allocate( uread(nshg,nsd) )
925*59599516SKenneth E. Jansen
926*59599516SKenneth E. Jansen         iusiz=nshg*nsd
927*59599516SKenneth E. Jansen
928*59599516SKenneth E. Jansen!          call readdatablock(irstin,fname1,uread,iusiz,
929*59599516SKenneth E. Jansen!      &        'double',iotype)
930*59599516SKenneth E. Jansen         call readdatablock(descriptor,fname1 // char(0) ,uread,iusiz,
931*59599516SKenneth E. Jansen     &                   'double' // char(0),iotype)
932*59599516SKenneth E. Jansen
933*59599516SKenneth E. Jansen         uold(:,1:nsd)=uread(:,1:nsd)
934*59599516SKenneth E. Jansen       else
935*59599516SKenneth E. Jansen         allocate( uold(nshg,nsd) )
936*59599516SKenneth E. Jansen         uold(:,1:nsd) = zero
937*59599516SKenneth E. Jansen       endif
938*59599516SKenneth E. Jansen
939*59599516SKenneth E. Jansenc
940*59599516SKenneth E. Jansenc.... close c-binary files
941*59599516SKenneth E. Jansenc
942*59599516SKenneth E. Jansen!MR CHANGE
943*59599516SKenneth E. Jansen!      call closefile( irstin, "read" )
944*59599516SKenneth E. Jansen
945*59599516SKenneth E. Jansen      call closefile( descriptor, "read" // char(0) )
946*59599516SKenneth E. Jansen      call finalizephmpiio( descriptor )
947*59599516SKenneth E. Jansen
948*59599516SKenneth E. Jansen!MR CHANGE
949*59599516SKenneth E. Jansen!      call closefile( igeomBAK,  "read" )
950*59599516SKenneth E. Jansenc
951*59599516SKenneth E. Jansen      deallocate(xread)
952*59599516SKenneth E. Jansen      if ( numpbc > 0 )  then
953*59599516SKenneth E. Jansen         deallocate(bcinpread)
954*59599516SKenneth E. Jansen         deallocate(ibctmpread)
955*59599516SKenneth E. Jansen      endif
956*59599516SKenneth E. Jansen      deallocate(iperread)
957*59599516SKenneth E. Jansen      if(numpe.gt.1)
958*59599516SKenneth E. Jansen     &     deallocate(ilworkread)
959*59599516SKenneth E. Jansen      deallocate(nbcread)
960*59599516SKenneth E. Jansen
961*59599516SKenneth E. Jansen      return
962*59599516SKenneth E. Jansenc
963*59599516SKenneth E. Jansen 994  call error ('input   ','opening ', igeomBAK)
964*59599516SKenneth E. Jansen 995  call error ('input   ','opening ', igeomBAK)
965*59599516SKenneth E. Jansen 997  call error ('input   ','end file', igeomBAK)
966*59599516SKenneth E. Jansen 998  call error ('input   ','end file', igeomBAK)
967*59599516SKenneth E. Jansenc
968*59599516SKenneth E. Jansen      end
969*59599516SKenneth E. Jansen
970*59599516SKenneth E. Jansenc
971*59599516SKenneth E. Jansenc No longer called but kept around in case....
972*59599516SKenneth E. Jansenc
973*59599516SKenneth E. Jansen      subroutine genpzero(iBC)
974*59599516SKenneth E. Jansen
975*59599516SKenneth E. Jansen      use pointer_data
976*59599516SKenneth E. Jansenc
977*59599516SKenneth E. Jansen      include "common.h"
978*59599516SKenneth E. Jansen      integer iBC(nshg)
979*59599516SKenneth E. Jansenc
980*59599516SKenneth E. Jansenc....  check to see if any of the nodes have a dirichlet pressure
981*59599516SKenneth E. Jansenc
982*59599516SKenneth E. Jansen      pzero=1
983*59599516SKenneth E. Jansen      if (any(btest(iBC,2))) pzero=0
984*59599516SKenneth E. Jansenc
985*59599516SKenneth E. Jansen      do iblk = 1, nelblb
986*59599516SKenneth E. Jansen         npro = lcblkb(1,iblk+1)-lcblkb(1,iblk)
987*59599516SKenneth E. Jansen         do i=1, npro
988*59599516SKenneth E. Jansen            iBCB1=miBCB(iblk)%p(i,1)
989*59599516SKenneth E. Jansenc
990*59599516SKenneth E. Jansenc.... check to see if any of the nodes have a Neumann pressure
991*59599516SKenneth E. Jansenc     but not periodic (note that
992*59599516SKenneth E. Jansenc
993*59599516SKenneth E. Jansen            if(btest(iBCB1,1)) pzero=0
994*59599516SKenneth E. Jansen         enddo
995*59599516SKenneth E. Jansenc
996*59599516SKenneth E. Jansenc.... share results with other processors
997*59599516SKenneth E. Jansenc
998*59599516SKenneth E. Jansen         pzl=pzero
999*59599516SKenneth E. Jansen         if (numpe .gt. 1)
1000*59599516SKenneth E. Jansen     &        call MPI_ALLREDUCE (pzl, pzero, 1,
1001*59599516SKenneth E. Jansen     &        MPI_DOUBLE_PRECISION,MPI_MIN, MPI_COMM_WORLD,ierr)
1002*59599516SKenneth E. Jansen
1003*59599516SKenneth E. Jansen      enddo
1004*59599516SKenneth E. Jansenc
1005*59599516SKenneth E. Jansenc.... return
1006*59599516SKenneth E. Jansenc
1007*59599516SKenneth E. Jansen      return
1008*59599516SKenneth E. Jansenc
1009*59599516SKenneth E. Jansen      end
1010*59599516SKenneth E. Jansen
1011