KTH framework for Nek5000 toolboxes; testing version  0.0.1
reader_rea.f
Go to the documentation of this file.
1 c-----------------------------------------------------------------------
2  subroutine rdparam
3 C
4 C .Read in parameters supplied by preprocessor and
5 C (eventually) echo check.
6 C
7 C .Broadcast run parameters to all processors
8 C
9  include 'SIZE'
10  include 'INPUT'
11  include 'PARALLEL'
12  include 'CTIMER'
13 
14  character*132 string(100)
15 
16  vnekton = 3 ! dummy not really used anymore
17 
18  optlevel = 1! fixed for now
19  loglevel = 1! fixed for now
20 
21  IF(nid.EQ.0) THEN
22  READ(9,*,err=400)
23  READ(9,*,err=400)
24  READ(9,*,err=400) ldimr
25  READ(9,*,err=400) nparam
26  DO 20 i=1,nparam
27  READ(9,*,err=400) param(i)
28  20 CONTINUE
29  ENDIF
30  call bcast(ldimr,isize)
31  call bcast(nparam,isize)
32  call bcast(param ,200*wdsize)
33 
34  npscal=int(param(23))
35  npscl1=npscal+1
36  npscl2=npscal+2
37 
38  IF (npscl1.GT.ldimt) THEN
39  if(nid.eq.0) then
40  WRITE(6,21) ldimt,npscl1
41  21 FORMAT(//,2x,'Error: This NEKTON Solver has been compiled'
42  $ /,2x,' for',i4,' passive scalars. This run'
43  $ /,2x,' requires that LDIMT be set to',i4,'.')
44  endif
45  call exitt
46  ENDIF
47 
48 c Use same tolerances for all fields
49  restol(0) = param(22) ! mesh
50  restol(1) = param(22)
51  do i=1,ldimt
52  restol(1+i) = param(22)
53  enddo
54  call bcast(restol, (ldimt1+1)*wdsize)
55 
56 c
57 c Read in the passive scalar conduct and rhocp's:
58 c
59 c fluid
60 c .viscosity is PARAM(2)
61 c .if it is negative, it indicates that Re has been input
62 c .therefore, redefine PARAM(2) = -1.0/PARAM(2)
63 c
64  if(param(2) .lt.0.0) param(2) = -1.0/param(2)
65  if(param(8) .lt.0.0) param(8) = -1.0/param(8)
66  if(param(29).lt.0.0) param(29) = -1.0/param(29)
67 C
68  cpfld(1,1)=param(2)
69  cpfld(1,2)=param(1)
70 C temperature
71  cpfld(2,1)=param(8)
72  cpfld(2,2)=param(7)
73  cpfld(2,3)=param(9)
74 c
75 c passive scalars
76 c
77  IF(nid.EQ.0) THEN
78  READ(9,*,err=400) nskip
79  IF (nskip.GT.0 .AND. npscal.GT.0) THEN
80  READ(9,*,err=400)(cpfld(i,1),i=3,npscl2)
81  IF(npscl2.LT.9)READ(9,*)
82  READ(9,*,err=400)(cpfld(i,2),i=3,npscl2)
83  IF(npscl2.LT.9)READ(9,*)
84  do i=3,npscl2
85  if (cpfld(i,1).lt.0) cpfld(i,1) = -1./cpfld(i,1)
86  if (cpfld(i,2).lt.0) cpfld(i,2) = -1./cpfld(i,2)
87  enddo
88  ELSE
89  DO 25 i=1,nskip
90  READ(9,*,err=500)
91  25 CONTINUE
92  ENDIF
93  ENDIF
94  call bcast(cpfld,wdsize*ldimt1*3)
95 
96 C
97 C Read logical equation type descriptors....
98 C
99  iftmsh(0) = .false.
100  ifprojfld(0) = .false.
101  ifdgfld(0) = .false.
102  do i=1,npscl2
103  iftmsh(i) = .false.
104  ifadvc(i) = .false.
105  ifdgfld(i) = .false.
106  iffilter(i) = .false.
107  ifdiff(i) = .true.
108  ifdeal(i) = .true. ! still depends on param(99)
109  ifprojfld(i) = .false.
110  if (param(94).gt.0) ifprojfld(i) = .true.
111  enddo
112 
113  do i=1,npscl1
114  idpss(i) = 0 ! use Helmholtz for passive scalars
115  enddo
116 
117  ifflow = .false.
118  ifheat = .false.
119  iftran = .false.
120  ifaxis = .false.
121  ifaziv = .false.
122  ifstrs = .false.
123  iflomach = .false.
124  ifmodel = .false.
125  ifkeps = .false.
126  ifmvbd = .false.
127  ifchar = .false.
128  ifdg = .false.
129  ifanls = .false.
130  ifmoab = .false.
131  ifcoup = .false.
132  ifvcoup = .false.
133  ifmhd = .false.
134  ifessr = .false.
135  iftmsh(0) = .false.
136  ifuservp = .false.
137  ifcons = .false. ! Use conservation form?
138  ifusermv = .false.
139  ifcyclic = .false.
140  ifsync = .false.
141  ifexplvis = .false.
142  ifschclob = .false.
143 c IFSPLIT = .false.
144 
145  ifbase = .true.
146  ifpert = .false.
147 
148  ifreguo = .false. ! by default we dump the data based on the GLL mesh
149 
150  ifrich = .false.
151 
152  IF(nid.EQ.0) READ(9,*,err=500) nlogic
153  call bcast(nlogic,isize)
154  IF(nlogic.GT.100) THEN
155  if(nid.eq.0)
156  $ write(6,*) 'ABORT: Too many logical switches', nlogic
157  call exitt
158  ENDIF
159 
160  if(nid.eq.0) READ(9,'(A132)',err=500) (string(i),i=1,nlogic)
161  call bcast(string,100*132*csize)
162 
163  do i = 1,nlogic
164  call capit(string(i),132)
165  if (indx1(string(i),'IFTMSH' ,6).gt.0) then
166  read(string(i),*,err=490) (iftmsh(ii),ii=0,npscl2)
167  elseif (indx1(string(i),'IFNAV' ,5).gt.0 .and.
168  & indx1(string(i),'IFADVC' ,6).gt.0) then
169  read(string(i),*,err=490) (ifadvc(ii),ii=1,npscl2)
170  elseif (indx1(string(i),'IFADVC' ,6).gt.0) then
171  read(string(i),*,err=490) (ifadvc(ii),ii=1,npscl2)
172  elseif (indx1(string(i),'IFFLOW' ,6).gt.0) then
173  read(string(i),*) ifflow
174  elseif (indx1(string(i),'IFHEAT' ,6).gt.0) then
175  read(string(i),*) ifheat
176  elseif (indx1(string(i),'IFTRAN' ,6).gt.0) then
177  read(string(i),*) iftran
178  elseif (indx1(string(i),'IFAXIS' ,6).gt.0) then
179  read(string(i),*) ifaxis
180  elseif (indx1(string(i),'IFAZIV' ,6).gt.0) then
181  read(string(i),*) ifaziv
182  elseif (indx1(string(i),'IFSTRS' ,6).gt.0) then
183  read(string(i),*) ifstrs
184  elseif (indx1(string(i),'IFLO' ,4).gt.0) then
185  read(string(i),*) iflomach
186  elseif (indx1(string(i),'IFMGRID',7).gt.0) then
187 c read(string(i),*) IFMGRID
188  elseif (indx1(string(i),'IFKEPS' ,6).gt.0) then
189  read(string(i),*) ifkeps
190  elseif (indx1(string(i),'IFMODEL',7).gt.0) then
191  read(string(i),*) ifmodel
192  elseif (indx1(string(i),'IFMVBD' ,6).gt.0) then
193  read(string(i),*) ifmvbd
194  elseif (indx1(string(i),'IFCHAR' ,6).gt.0) then
195  read(string(i),*) ifchar
196  elseif (indx1(string(i),'IFDG' ,4).gt.0) then
197  read(string(i),*) ifdg
198  elseif (indx1(string(i),'IFANLS' ,6).gt.0) then
199  read(string(i),*) ifanls
200  elseif (indx1(string(i),'IFCOUP' ,6).gt.0) then
201  read(string(i),*) ifcoup
202  elseif (indx1(string(i),'IFVCOUP' ,7).gt.0) then
203  read(string(i),*) ifvcoup
204  elseif (indx1(string(i),'IFMHD' ,5).gt.0) then
205  read(string(i),*) ifmhd
206  elseif (indx1(string(i),'IFCONS' ,6).gt.0) then
207  read(string(i),*) ifcons
208  elseif (indx1(string(i),'IFUSERVP',8).gt.0) then
209  read(string(i),*) ifuservp
210  elseif (indx1(string(i),'IFUSERMV',8).gt.0) then
211  read(string(i),*) ifusermv
212  elseif (indx1(string(i),'IFCYCLIC',8).gt.0) then
213  read(string(i),*) ifcyclic
214  elseif (indx1(string(i),'IFPERT' ,6).gt.0) then
215  read(string(i),*) ifpert
216  elseif (indx1(string(i),'IFBASE' ,6).gt.0) then
217  read(string(i),*) ifbase
218  elseif (indx1(string(i),'IFSYNC' ,6).gt.0) then
219  read(string(i),*) ifsync
220  elseif (indx1(string(i),'IFSCHCLOB',9).gt.0) then
221  read(string(i),*) ifschclob
222  elseif (indx1(string(i),'IFSPLIT' ,7).gt.0) then
223 c read(string,*) IFSPLIT
224  else
225  if(nid.eq.0) then
226  write(6,'(1X,2A)') 'ABORT: Unknown logical flag', string
227  write(6,'(30(A,/))')
228  & ' Available logical flags:',
229  & ' IFTMSH' ,
230  & ' IFADVC' ,
231  & ' IFFLOW' ,
232  & ' IFHEAT' ,
233  & ' IFTRAN' ,
234  & ' IFAXIS' ,
235  & ' IFCYCLIC' ,
236  & ' IFSTRS' ,
237  & ' IFLOMACH' ,
238  & ' IFMGRID' ,
239  & ' IFKEPS' ,
240  & ' IFMVBD' ,
241  & ' IFCHAR' ,
242  & ' IFDG' ,
243  & ' IFANLS' ,
244  & ' IFUSERVP' ,
245  & ' IFUSERMV' ,
246  & ' IFSYNC' ,
247  & ' IFCYCLIC' ,
248  & ' IFSPLIT' ,
249  & ' IFEXPLVIS',
250  & ' IFCONS' ,
251  & ' IFCOUP' ,
252  & ' IFVCOUP'
253  endif
254  call exitt
255  endif
256  490 continue
257  enddo
258 
259  ifmgrid = .false.
260  if (ifsplit) ifmgrid = .true.
261 
262  if (ifaxis.and..not.ifsplit) then ! Use standard Schwarz/PCG solver
263  ifmgrid = .false.
264  param(42) = 1.00000 ! p042 0=gmres/1=pcg
265  param(43) = 1.00000 ! p043 0=semg/1=schwarz
266  param(44) = 1.00000 ! p044 0=E-based/1=A-based prec.
267  endif
268 
269  if (param(29).ne.0.) ifmhd = .true.
270  if (ifmhd) ifessr = .true.
271  if (ifmhd) npscl1 = npscl1 + 1
272  if (param(30).gt.0) ifuservp = .true.
273  if (param(31).ne.0.) ifpert = .true.
274  if (param(31).lt.0.) ifbase = .false. ! don't time adv base flow
275  npert = abs(param(31))
276 
277  IF (npscl1.GT.ldimt .AND. ifmhd) THEN
278  if(nid.eq.0) then
279  WRITE(6,22) ldimt,npscl1
280  22 FORMAT(/s,2x,'Error: This NEKTON Solver has been compiled'
281  $ /,2x,' for',i4,' passive scalars. A MHD run'
282  $ /,2x,' requires that LDIMT be set to',i4,'.')
283  endif
284  call exitt
285  ENDIF
286 
287  if (ifmvbd) then
288  if (lx1.ne.lx1m.or.ly1.ne.ly1m.or.lz1.ne.lz1m)
289  $ call exitti('Need lx1m=lx1 etc. in SIZE . $',lx1m)
290  endif
291 
292  ifldmhd = npscal + 3
293  if (ifmhd) then
294  cpfld(ifldmhd,1) = param(29) ! magnetic viscosity
295  cpfld(ifldmhd,2) = param( 1) ! magnetic rho same as for fluid
296  endif
297 C
298 C Set up default time dependent coefficients - NSTEPS,DT.
299 C
300  if (.not.iftran) then
301  if (ifflow.and.ifsplit) then
302  iftran=.true.
303  else
304  param(11) = 1.0
305  param(12) = 1.0
306  param(19) = 0.0
307  endif
308  endif
309 C
310 C Do some checks
311 C
312  IF(ldimr.NE.ldim)THEN
313  IF(nid.EQ.0) THEN
314  WRITE(6,10) ldim,ldimr
315  10 FORMAT(//,2x,'ERROR: This NEKTON Solver has been compiled'
316  $ /,2x,' for spatial dimension equal to',i2,'.'
317  $ /,2x,' The data file has dimension',i2,'.')
318  ENDIF
319  call exitt
320  ENDIF
321  IF (ldim.EQ.3) if3d=.true.
322  IF (ldim.NE.3) if3d=.false.
323 
324  if (if3d) then
325  if (ly1.ne.lx1.or.lz1.ne.lx1) then
326  if (nid.eq.0) write(6,13) lx1,ly1,lz1
327  13 format('ERROR: lx1,ly1,lz1:',3i5,' must be equal for 3D')
328  call exitt
329  endif
330  if (ly2.ne.lx2.or.lz2.ne.lx2) then
331  if (nid.eq.0) write(6,14) lx2,ly2,lz2
332  14 format('ERROR: lx2,ly2,lz2:',3i5,' must be equal for 3D')
333  call exitt
334  endif
335  else
336  if (ly1.ne.lx1.or.lz1.ne.1) then
337  if (nid.eq.0) write(6,12) lx1,ly1,lz1
338  12 format('ERROR: ',3i5,' must have lx1=ly1; lz1=1, for 2D')
339  call exitt
340  endif
341  if (ly2.ne.lx2.or.lz2.ne.1) then
342  if (nid.eq.0) write(6,11) lx2,ly2,lz2
343  11 format('ERROR: ',3i5,' must have lx2=ly2; lz2=1, for 2D')
344  call exitt
345  endif
346  endif
347 
348  if (lgmres.lt.5 .and. param(42).eq.0) then
349  if(nid.eq.0) write(6,*)
350  $ 'WARNING: lgmres might be too low!'
351  endif
352 
353 
354  if (ifsplit) then
355  if (lx1.ne.lx2) then
356  if (nid.eq.0) write(6,43) lx1,lx2
357  43 format('ERROR: lx1,lx2:',2i4,' must be equal for IFSPLIT=T')
358  call exitt
359  endif
360  else
361  if (lx2.lt.lx1-2) then
362  if (nid.eq.0) write(6,44) lx1,lx2
363  44 format('ERROR: lx1,lx2:',2i4,' lx2 must be lx-2 for IFSPLIT=F')
364  call exitt
365  endif
366  endif
367 
368  if (param(40).eq.3 .and. .not.ifsplit) then
369  call exitti
370  $ ('ERROR: Selected preconditioner requires lx2=lx1$',lx2)
371  endif
372 
373  if (ifcvode) then
374  if(nid.eq.0) write(6,*)
375  $ 'ABORT: Using CVODE requires .par file!'
376  call exitt
377  endif
378 
379  if (ifsplit .and. ifuservp .and. .not.ifstrs) then
380  if(nid.eq.0) write(6,*)
381  $ 'Enable stress formulation to support PN/PN and IFUSERVP=T'
382  ifstrs = .true.
383  endif
384 
385  if (ifcyclic .and. .not.ifstrs) then
386  if(nid.eq.0) write(6,*)
387  $ 'Enable stress formulation to support cyclic BC'
388  ifstrs = .true.
389  endif
390 
391  ktest = (lx1-lx1m) + (ly1-ly1m) + (lz1-lz1m)
392  if (ifstrs.and.ktest.ne.0) then
393  if(nid.eq.0) write(6,*)
394  $ 'ABORT: Stress formulation requires lx1m=lx1, etc. in SIZE'
395  call exitt
396  endif
397 
398 c if (ifsplit .and. ifstrs) then
399 c if(nid.eq.0) write(6,*)
400 c $ 'ABORT: Stress formulation in Pn-Pn is not supported'
401 c call exitt
402 c endif
403 
404  if (ifsplit .and. ifmhd) then
405  if(nid.eq.0) write(6,*)
406  $ 'ABORT: MHD in Pn-Pn is not supported'
407  call exitt
408  endif
409 
410  if (ifneknekc.and.(nelgv.ne.nelgt)) call exitti(
411  $ 'ABORT: nek-nek not supported w/ conj. ht transfer$',1)
412 
413  if (ifchar.and.(nelgv.ne.nelgt)) call exitti(
414  $ 'ABORT: IFCHAR curr. not supported w/ conj. ht transfer$',nelgv)
415 
416  if (ifmhd .and. lbx1.ne.lx1) then
417  if(nid.eq.0) write(6,*)
418  $ 'ABORT: For MHD, need lbx1=lx1, etc.; Change SIZE '
419  call exitt
420  endif
421 
422  if (ifpert .and. lpx1.ne.lx1) then
423  if(nid.eq.0) write(6,*)
424  $ 'ABORT: For Lyapunov, need lpx1=lx1, etc.; Change SIZE '
425  endif
426 
427  if (if3d) ifaxis = .false.
428 
429  if (iflomach .and. .not.ifsplit) then
430  if(nid.eq.0) write(6,*)
431  $ 'ABORT: For lowMach, need lx2=lx1, etc.; Change SIZE '
432  call exitt
433  endif
434 
435  if (iflomach .and. .not.ifheat) then
436  if(nid.eq.0) write(6,*)
437  $ 'ABORT For lowMach, need ifheat=true; Change IFHEAT'
438  call exitt
439  endif
440 
441  if (ifmhd) ifchar = .false. ! For now, at least.
442 
443 c set dealiasing handling
444  if (param(99).lt.0) then
445  param(99) = -1 ! No dealiasing
446  else
447  param(99) = 4 ! default
448  if (ifaxis) param(99) = 3 ! For now, at least.
449  if (ifmvbd) param(99) = 3 ! For now, at least.
450  endif
451 
452  if (ifchar .and. param(99).lt.0) then
453  if (nid.eq.0) write(6,*)
454  & 'ABORT: Characteristic scheme needs dealiasing!'
455  call exitt
456  endif
457 
458  if (.not.ifsplit .and. ifaxis .and. ifstrs) then
459  if (nid.eq.0) write(6,*)
460  $ 'ABORT: Axisymetric and stress formulation not supported ' //
461  $ 'for PN/PN-2$'
462  call exitt
463  endif
464 
465  if (param(99).gt.-1 .and. (lxd.lt.lx1 .or. lyd.lt.ly1 .or.
466  & lzd.lt.lz1)) then
467  if(nid.eq.0) write(6,*)
468  & 'ABORT: Dealiasing space too small; Check lxd,lyd,lzd in SIZE '
469  call exitt
470  endif
471 
472 c SET PRESSURE SOLVER DEFAULTS, ADJUSTED IN USR FILE ONLY
473  param(41) = 0 ! use additive SEMG
474  ! 1 use hybrid SEMG (not yet working... but coming soon!)
475  param(42) = 0 ! use GMRES for iterative solver w/ nonsymmetric weighting
476  ! 1 use PCG for iterative solver, do not use weighting
477  param(43) = 0 ! use additive multilevel scheme (requires param(42).eq.0)
478  ! 1 use original 2 level scheme
479  param(44) = 0 ! base top-level additive Schwarz on restrictions of E
480  ! 1 base top-level additive Schwarz on restrictions of A
481 
482 c SET DEFAULT TO 6, ADJUSTED IN USR FILE ONLY
483  param(66) = 6
484  param(67) = 6
485 
486  param(59) = 1 ! No fast operator eval, ADJUSTED IN USR FILE ONLY
487  param(33) = 0
488 
489  fem_amg_param(1) = 0
490  crs_param(1) = 0
491 
492  filtertype = 0
493  if (param(103).gt.0) then
494  filtertype = 1
495  call ltrue(iffilter,size(iffilter))
496  endif
497 
498  return
499 
500 C
501 C Error handling:
502 C
503  400 CONTINUE
504  if(nid.eq.0) WRITE(6,401)
505  401 FORMAT(2x,'ERROR READING PARAMETER DATA'
506  $ ,/,2x,'ABORTING IN ROUTINE RDPARAM.')
507  call exitt
508 C
509  500 CONTINUE
510  if(nid.eq.0) WRITE(6,501)
511  501 FORMAT(2x,'ERROR READING LOGICAL DATA'
512  $ ,/,2x,'ABORTING IN ROUTINE RDPARAM.')
513  call exitt
514 C
515  return
516  END
517 c-----------------------------------------------------------------------
518  subroutine rdmesh
519 C
520 C .Read number of elements
521 C
522 C .Construct sequential element-processor partition according
523 C to number of elements and processors
524 C
525 C .Selectively read mesh (defined by element vertices, and group numbers)
526 C on each processor
527 C
528  include 'SIZE'
529  include 'INPUT'
530  include 'PARALLEL'
531  character*1 adum
532  real dum(4)
533 
534 
535 c Read elemental mesh data, formatted
536  iffmtin = .true.
537 
538  nsides=ldim*2
539  DO 40 ieg=1,nelgt
540  IF (gllnid(ieg).EQ.nid) THEN
541  iel=gllel(ieg)
542 
543  igroup(iel) = 0
544 c read(9,30,err=31,end=600) igroup(iel)
545 c 30 format(43x,i5)
546  read(9,*,err=31,end=600) adum
547  31 continue
548 
549 C Read Corner data
550  IF(ldim.EQ.2)THEN
551  READ(9,*,err=500,END=600) (XC(IC,IEL),IC=1,4)
552  READ(9,*,err=500,END=600) (YC(IC,IEL),IC=1,4)
553  call rzero (zc(1 ,iel) ,4)
554  ELSE IF(ldim.EQ.3)THEN
555  READ(9,*,err=500,END=600) (XC(IC,IEL),IC=1,4)
556  READ(9,*,err=500,END=600) (YC(IC,IEL),IC=1,4)
557  READ(9,*,err=500,END=600) (ZC(IC,IEL),IC=1,4)
558  READ(9,*,err=500,END=600) (XC(IC,IEL),IC=5,8)
559  READ(9,*,err=500,END=600) (YC(IC,IEL),IC=5,8)
560  READ(9,*,err=500,END=600) (ZC(IC,IEL),IC=5,8)
561  ENDIF
562  ELSE
563 C Skip over this data for element NOT on this processor
564  READ(9,41,err=500,END=600) adum
565 C Read Corner data
566  IF(ldim.EQ.2)THEN
567  READ(9,41,err=500,END=600) adum
568  READ(9,41,err=500,END=600) adum
569  ELSE IF(ldim.EQ.3)THEN
570  READ(9,41,err=500,END=600) adum
571  READ(9,41,err=500,END=600) adum
572  READ(9,41,err=500,END=600) adum
573  READ(9,41,err=500,END=600) adum
574  READ(9,41,err=500,END=600) adum
575  READ(9,41,err=500,END=600) adum
576  ENDIF
577  ENDIF
578  40 CONTINUE
579  41 FORMAT(a1)
580 C
581 C End of mesh read.
582 C
583  return
584 C
585 C Error handling:
586 C
587  400 CONTINUE
588  if(nid.eq.0) WRITE(6,401)
589  401 FORMAT(2x,'ERROR READING SCALE FACTORS, CHECK READ FILE'
590  $ ,/,2x,'ABORTING IN ROUTINE RDMESH.')
591  call exitt
592 
593  500 CONTINUE
594  if(nid.eq.0) WRITE(6,501) ieg
595  501 FORMAT(2x,'ERROR READING MESH DATA NEAR ELEMENT',i12
596  $ ,/,2x,'ABORTING IN ROUTINE RDMESH.')
597  call exitt
598 
599  600 CONTINUE
600  if(nid.eq.0) WRITE(6,601) ieg
601  601 FORMAT(2x,'ERROR 2 READING MESH DATA NEAR ELEMENT',i12
602  $ ,/,2x,'ABORTING IN ROUTINE RDMESH.')
603  call exitt
604 
605  return
606  end
607 c-----------------------------------------------------------------------
608  subroutine rdcurve
609 C
610 C .Read curve side data
611 C
612 C .Disperse curve side data to all processors according
613 C to sequential partition scheme
614 C
615 C
616  include 'SIZE'
617  include 'INPUT'
618  include 'PARALLEL'
619  CHARACTER*1 ANS
620 C
621 C
622 C
623  IF (iffmtin) THEN
624 C
625 C Read formatted curve side data
626 C
627  READ(9,*)
628  READ(9,*)ncurve
629  CALL rzero(curve ,72*lelt)
630  CALL blank(ccurve,12*lelt)
631  IF (ncurve.GT.0) THEN
632  DO 50 icurve=1,ncurve
633  IF (nelgt.LT.1000) THEN
634  READ(9,60,err=500,END=500) IEDG,IEG,R1,R2,R3,R4,R5,ans
635  ELSEIF (nelgt.LT.1 000 000) THEN
636  READ(9,61,err=500,END=500) IEDG,IEG,R1,R2,R3,R4,R5,ans
637  ELSE
638  READ(9,62,err=500,END=500) IEDG,IEG,R1,R2,R3,R4,R5,ans
639  ENDIF
640  60 FORMAT(i3,i3 ,5g14.6,1x,a1)
641  61 FORMAT(i2,i6 ,5g14.6,1x,a1)
642  62 FORMAT(i2,i12,5g14.6,1x,a1)
643 
644  IF (gllnid(ieg).EQ.nid) THEN
645  iel=gllel(ieg)
646  curve(1,iedg,iel)=r1
647  curve(2,iedg,iel)=r2
648  curve(3,iedg,iel)=r3
649  curve(4,iedg,iel)=r4
650  curve(5,iedg,iel)=r5
651  ccurve( iedg,iel)=ans
652  ENDIF
653  50 CONTINUE
654  ENDIF
655  return
656 C
657 C Error handling:
658 C
659  500 CONTINUE
660  if(nid.eq.0) WRITE(6,501)
661  501 FORMAT(2x,'ERROR READING CURVE SIDE DATA'
662  $ ,/,2x,'ABORTING IN ROUTINE RDCURVE.')
663  call exitt
664  return
665 C
666  ELSE
667 C
668 C Read unformatted curve side data
669 C
670  READ(8) ncurve
671  CALL rzero(curve ,72*lelt)
672  CALL blank(ccurve,12*lelt)
673  IF (ncurve.GT.0) THEN
674  DO 1050 icurve=1,ncurve
675  READ(8,err=1500,END=1500) IEDG,IEG,R1,R2,R3,R4,R5,ans
676  IF (gllnid(ieg).EQ.nid) THEN
677  iel=gllel(ieg)
678  curve(1,iedg,iel)=r1
679  curve(2,iedg,iel)=r2
680  curve(3,iedg,iel)=r3
681  curve(4,iedg,iel)=r4
682  curve(5,iedg,iel)=r5
683  ccurve( iedg,iel)=ans
684  ENDIF
685  1050 CONTINUE
686  ENDIF
687  return
688 C
689 C Error handling:
690 C
691  1500 CONTINUE
692  if(nid.eq.0) WRITE(6,1501)
693  1501 FORMAT(2x,'ERROR READING unformatted CURVE SIDE DATA'
694  $ ,/,2x,'ABORTING IN ROUTINE RDCURVE.')
695  call exitt
696 C
697  return
698  ENDIF
699  END
700 c-----------------------------------------------------------------------
701  subroutine rdbdry
702 C
703 C .Read Boundary Conditions (and connectivity data)
704 C
705 C .Disperse boundary condition data to all processors
706 C according to sequential partition scheme
707 C
708  include 'SIZE'
709  include 'INPUT'
710  include 'PARALLEL'
711  include 'SCRCT'
712  CHARACTER CBC1*1,CBC3*3,CHTEMP*1,CHTMP3*3
713  equivalence(chtemp,chtmp3)
714  character*132 string
715 C
716 C Set up TEMPORARY value for NFIELD - NFLDT
717 C
718  nfldt = 1
719  IF (ifheat) nfldt=2+npscal
720  if (ifmhd ) nfldt=2+npscal+1
721  nbcs = nfldt
722  ibcs = 2
723  IF (ifflow) ibcs = 1
724  nsides = 2*ldim
725 C
726 C Read boundary conditions for all fields
727 C
728  lcbc=18*lelt*(ldimt1 + 1)
729  lrbc=30*lelt*(ldimt1 + 1)
730  CALL rzero(bc ,lrbc)
731  CALL blank(cbc,lcbc)
732 C
733 C-----------------------------------------------------------------
734 C Formatted Reads
735 C-----------------------------------------------------------------
736 C
737  IF (iffmtin) THEN
738 C
739  READ(9,*,err=500,END=500) ! ***** BOUNDARY CONDITIONS *****
740  ibcnew = 1
741  DO 100 ifield=ibcnew,nbcs ! DO 100 IFIELD=IBCS,NBCS
742  nel=nelgt
743  if (.not.iftmsh(ifield)) nel = nelgv
744 C Fluid and/or thermal
745  read(9,81) string ! ***** FLUID BOUNDARY CONDITIONS *****
746  call capit(string,132)
747 
748 c write(6,*) 'reading BC:',ifield,ibcs,nbcs
749 c write(6,81) string
750 c if1 = indx1(string,'NO ',3)
751 c write(6,*) if1,' if NO. quit.',ifield,ibcs,nbcs
752 c write(6,*) ifield,iftmsh(ifield),nel,' iftmsh'
753 c call exitt
754 
755 
756  if (indx1(string,'NO ',3).eq.0) then ! we have acitve bc info
757 C
758  IF(vnekton .LE. 2.52) nbcrea = 3
759  IF(vnekton .GE. 2.55) nbcrea = 5
760 C
761  DO 80 ieg=1,nel
762  DO 80 iside=1,nsides
763  IF (gllnid(ieg).EQ.nid) THEN
764  iel=gllel(ieg)
765  IF (nelgt.LT.1000) THEN
766  READ(9,50,err=500,END=500)
767  $ chtemp,
768  $ cbc(iside,iel,ifield),id1,id2,
769  $ (bc(ii,iside,iel,ifield),ii=1,nbcrea)
770 c write(6,50)
771 c $ CHTEMP,
772 c $ CBC(ISIDE,IEL,IFIELD),ID1,ID2,
773 c $ (BC(II,ISIDE,IEL,IFIELD),II=1,NBCREA)
774  50 FORMAT(a1,a3,2i3,5g14.6)
775  ELSEIF (nelgt.LT.100 000) THEN
776  READ(9,51,err=500,END=500)
777  $ chtemp,
778  $ cbc(iside,iel,ifield),id1,id2,
779  $ (bc(ii,iside,iel,ifield),ii=1,nbcrea)
780  51 FORMAT(a1,a3,i5,i1,5g14.6)
781  ELSEIF (nelgt.LT.1 000 000) THEN
782  READ(9,52,err=500,END=500)
783  $ chtemp,
784  $ cbc(iside,iel,ifield),id1,
785  $ (bc(ii,iside,iel,ifield),ii=1,nbcrea)
786  52 FORMAT(a1,a3,i6,5g14.6)
787  ELSE
788  READ(9,53,err=500,END=500)
789  $ chtemp,
790  $ cbc(iside,iel,ifield),id1,
791  $ (bc(ii,iside,iel,ifield),ii=1,nbcrea)
792  53 FORMAT(a1,a3,i12,5g18.11)
793  ENDIF
794 C Mesh B.C.'s in 1st column of 1st field
795  IF (chtemp.NE.' ') cbc(iside,iel,0)(1:1)= chtemp
796 C check for fortran function as denoted by lower case bc's:
797  cbc1=cbc(iside,iel,ifield)
798  cbc3=cbc(iside,iel,ifield)
799  icbc1=ichar(cbc1)
800 c IF (ICBC1.GE.97.AND.ICBC1.LE.122) THEN
801 c IF(CBC3(3:3).NE.'i')NLINES=BC(1,ISIDE,IEL,IFIELD)
802 c IF(CBC3(3:3).EQ.'i')NLINES=BC(4,ISIDE,IEL,IFIELD)
803 c DO 60 I=1,NLINES
804 c 60 READ(9,*,ERR=500,END=500)
805 c ENDIF
806  ELSE
807  READ(9,*,err=500,END=500) cbc1 ! dummy read, pff 4/28/05
808  ENDIF
809  80 CONTINUE
810  endif
811  81 format(a132)
812  100 CONTINUE
813 C
814 C END OF BC READ
815 C
816 C Check for dummy line: "NO THERMAL B.C.'S"
817  IF (nfldt.EQ.1) READ(9,*,err=500,END=500)
818 C
819  return
820 C
821 C Error handling:
822 C
823  500 CONTINUE
824  if(nid.eq.0) WRITE(6,501) ifield,ieg
825  501 FORMAT(2x,'ERROR READING BOUNDARY CONDITIONS FOR FIELD',i4,i12
826  $ ,/,2x,'ABORTING IN ROUTINE RDBDRY.')
827  call exitt
828  return
829 C
830 C
831  ELSE
832 C
833 C-----------------------------------------------------------------
834 C UNformatted Reads
835 C-----------------------------------------------------------------
836 C
837 c READ(8,ERR=500,END=500)
838  DO 1100 ifield=ibcs,nbcs
839  nel=nelgt
840 C Fluid and/or thermal
841  nbcrea = 5
842 C
843  DO 1080 ieg=1,nel
844  DO 1080 iside=1,nsides
845  IF (gllnid(ieg).EQ.nid) THEN
846  iel=gllel(ieg)
847  READ(8,err=1500,END=1500)
848  $ chtmp3,
849  $ cbc(iside,iel,ifield),id1,id2,
850  $ (bc(ii,iside,iel,ifield),ii=1,nbcrea)
851 C
852 C Mesh B.C.'s in 1st column of 1st field
853  IF (chtemp.NE.' ') cbc(iside,iel,0)(1:1)= chtemp
854 C check for fortran function as denoted by lower case bc's:
855  ELSE
856  iel=1
857  READ(8,err=1500,END=1500) CHTMP3,
858  $ cbcs(iside,iel),id1,id2,(bcs(ii,iside,iel),ii=1,nbcrea)
859 C check for fortran function as denoted by lower case bcs:
860  ENDIF
861  1080 CONTINUE
862  1100 CONTINUE
863 C
864 C END OF BC READ
865 C
866  return
867 C
868 C Error handling:
869 C
870  1500 CONTINUE
871  if(nid.eq.0) WRITE(6,1501) ifield,ieg
872  1501 FORMAT(2x,'ERROR READING BOUNDARY CONDITIONS FOR FIELD',i4,i12
873  $ ,/,2x,'(unformatted) ABORTING IN ROUTINE RDBDRY.')
874  call exitt
875  ENDIF
876 C
877  return
878  END
879 c-----------------------------------------------------------------------
880  subroutine rdicdf
881 C
882 C .Read Initial Conditions / Drive Force
883 C
884 C .Broadcast ICFILE to all processors
885 C
886  include 'SIZE'
887  include 'INPUT'
888  include 'PARALLEL'
889 
890  character*132 line
891  logical ifgtil
892 
893  ierr = 0
894 
895  if (nid.eq.0) then ! Read names of restart files
896 
897  call blank(initc,15*132)
898  read (9,80,err=200,end=200) line
899  call capit(line,132)
900  if (indx1(line,'RESTART',7).ne.0) then
901  if (.not.ifgtil(nskip,line)) goto 200
902 C read(line,*,err=200,end=200) nskip
903  do 50 i=1,nskip
904  read(9,80,err=200,end=200) initc(i)
905  50 continue
906  read(9,80,err=200,end=200) line
907  endif
908  80 format(a132)
909 
910  if (.not.ifgtil(nskip,line)) goto 200
911 
912 C Read initial conditions
913  do 100 i=1,nskip
914  read(9,80,err=200,end=200) line
915  100 continue
916 
917 C Read drive force data
918  read(9,*,err=200,end=200)
919  read(9,*,err=200,end=200) nskip
920  do 110 i=1,nskip
921  read(9,80,err=200,end=200) line
922  110 continue
923  endif
924 
925  ierr = iglmax(ierr,1)
926  if (ierr.eq.0) then
927  call bcast(initc,15*132*csize)
928  return
929  else
930  goto 210
931  endif
932 c
933 c Error handling:
934 c
935  200 ierr = 1
936  ierr = iglmax(ierr,1)
937 
938  210 continue
939  if (nid.eq.0) write(6,300)
940  300 format(2x,'Error reading initial condition/drive force data'
941  $ ,/,2x,'aborting in routine rdicdf.')
942  call exitti('rdicdf error$',ierr)
943 
944  return
945  end
946 c-----------------------------------------------------------------------
947  subroutine rdmatp
948 C
949 C .Read materials property data
950 C
951 C .Disperse material properties to all processors according
952 C to sequential partition scheme
953 C
954  include 'SIZE'
955  include 'INPUT'
956  include 'PARALLEL'
957 
958  CHARACTER*132 LINE
959 C
960  CALL izero(matype,16*ldimt1)
961  CALL rzero(cpgrp ,48*ldimt1)
962 C
963 C Read material property data
964 C
965  IF(nid.EQ.0) THEN
966  READ(9,*,err=200,END=200)
967  READ(9,*,err=200,END=200) nskip
968  READ(9,*,err=200,END=200) npacks
969  DO 100 iig=1,npacks
970  ifvps=.true.
971  READ(9,*)igrp,ifld,itype
972  matype(igrp,ifld)=itype
973  DO 100 iprop=1,3
974  IF(itype.EQ.1) READ(9,* ) cpgrp(igrp,ifld,iprop)
975  IF(itype.EQ.2) READ(9,80) line
976  80 FORMAT(a132)
977  100 CONTINUE
978  ENDIF
979 
980  CALL bcast(matype,16*ldimt1*isize)
981  CALL bcast(cpgrp ,48*ldimt1*wdsize)
982 
983  return
984 C
985 C Error handling:
986 C
987  200 CONTINUE
988  if(nid.eq.0) WRITE(6,201)
989  201 FORMAT(2x,'ERROR READING MATERIAL PROPERTIES DATA'
990  $ ,/,2x,'ABORTING IN ROUTINE RDMATP.')
991  call exitt
992 C
993  return
994  END
995 c-----------------------------------------------------------------------
996  subroutine rdhist
997 C
998 C .Read history data
999 C
1000 C .Broadcast to all processors
1001 C
1002  include 'SIZE'
1003  include 'INPUT'
1004  include 'PARALLEL'
1005 C
1006  ierr=0
1007  if(nid.eq.0) then
1008 c read history data
1009  read (9,*)
1010  read (9,*,err=200,end=200) nhis
1011  do i = 1,nhis
1012  read (9,*)
1013  enddo
1014  endif
1015 
1016  return
1017 C
1018 C Error handling:
1019 C
1020  200 CONTINUE
1021  if(nid.eq.0) WRITE(6,201)
1022  201 FORMAT(2x,'ERROR READING HISTORY DATA'
1023  $ ,/,2x,'ABORTING IN ROUTINE RDHIST.')
1024  call exitt
1025 C
1026  return
1027  END
1028 c-----------------------------------------------------------------------
1029  subroutine rdout
1030 C
1031 C .Read output specs
1032 C
1033 C .Broadcast to all processors
1034 C
1035  include 'SIZE'
1036  include 'INPUT'
1037  include 'PARALLEL'
1038 
1039  logical lbuf(5+ldimt1)
1040 
1041  call lfalse(lbuf,5+ldimt1)
1042  iflag = 0 ! Check for valid ipsco read
1043 
1044  IF(nid.EQ.0) THEN ! Read output specs
1045 
1046  READ(9,*,err=200,END=200)
1047  READ(9,*,err=200,END=200) nouts
1048  READ(9,*,err=200,END=200) ifxyo
1049  READ(9,*,err=200,END=200) ifvo
1050  READ(9,*,err=200,END=200) ifpo
1051  READ(9,*,err=200,END=200) ifto
1052  READ(9,*,err=200,END=200) ifbo ! IFTGO
1053 
1054  lbuf(1) = ifxyo
1055  lbuf(2) = ifvo
1056  lbuf(3) = ifpo
1057  lbuf(4) = ifto
1058  lbuf(5) = ifbo
1059 
1060  k = 5
1061 
1062  call lfalse(ifpsco,ldimt1)
1063  read(9,*,err=200,end=200) ipsco
1064  if (ipsco.gt.0) then
1065  if (ipsco.gt.ldimt1) then ! Invalid ifpsco read
1066  iflag = 1
1067  else
1068  do i=1,ipsco
1069  read(9,*,err=200,end=200) ifpsco(i)
1070  k = k+1
1071  lbuf(k) = ifpsco(i)
1072  enddo
1073  endif
1074  endif
1075 
1076  endif
1077 
1078 
1079  iflag = iglmax(iflag,1) ! Check for valid ipsco read
1080  if (iflag.gt.0) call exitti ! Invalid ifpsco read
1081  $ ('Error in rdout. Increase ldimt1 in SIZE to$',ipsco)
1082 
1083  k = 5+ldimt1
1084  call bcast(lbuf ,lsize*k)
1085  call bcast(ipsco,isize )
1086 
1087  ifxyo = lbuf(1)
1088  ifvo = lbuf(2)
1089  ifpo = lbuf(3)
1090  ifto = lbuf(4)
1091  ifbo = lbuf(5)
1092 
1093  k = 5
1094  do i=1,ipsco
1095  k = k+1
1096  ifpsco(i) = lbuf(k)
1097  enddo
1098 
1099  return
1100 
1101 C
1102 C Error handling:
1103 C
1104  200 CONTINUE
1105  WRITE(6,201)
1106  201 FORMAT(2x,'ERROR READING OUTPUT SPECIFICATION DATA'
1107  $ ,/,2x,'ABORTING IN ROUTINE RDOUT.')
1108  call exitt
1109 C
1110  return
1111  END
1112 c-----------------------------------------------------------------------
1113  subroutine rdobj
1114 C
1115 C .Read objects
1116 C
1117 C .Broadcast to all processors
1118 C
1119  include 'SIZE'
1120  include 'INPUT'
1121  include 'PARALLEL'
1122 C
1123 C Default if no data is read No Objects
1124 C
1125  ierr=0
1126  IF(nid.EQ.0) THEN
1127  nobj=0
1128  READ(9,*,err=200,END=200)
1129  READ(9,*,err=200,END=200) nobj
1130 
1131  IF(nobj.GT.maxobj) ierr=1
1132 
1133  if(ierr.eq.0) then
1134  DO 10 iobj = 1,nobj
1135  READ(9,*,err=200,END=200) NMEMBER(IOBJ)
1136  IF(nmember(iobj).GT.maxmbr)THEN
1137  print*,'ERROR: Too many members in object ',iobj
1138  ierr=2
1139  ENDIF
1140  if(ierr.eq.0) then
1141  DO 5 member=1,nmember(iobj)
1142  READ(9,*,err=200,END=200) OBJECT(IOBJ,MEMBER,1),
1143  $ object(iobj,member,2)
1144  5 CONTINUE
1145  endif
1146  10 CONTINUE
1147  write(6,*) nobj,' objects found'
1148  $ ,(nmember(k),k=1,nobj)
1149  endif
1150  endif
1151  call err_chk(ierr,'ERROR, too many objects:$')
1152 
1153  call bcast(nobj ,isize)
1154  call bcast(nmember,maxobj*isize)
1155  call bcast(object ,maxobj*maxmbr*2*isize)
1156 
1157 
1158  return
1159 C
1160 C Error handling: For old versions, default to no objects
1161 C
1162  200 CONTINUE
1163  nobj=0
1164 
1165  return
1166  END
1167 
subroutine lfalse(IFA, N)
Definition: bdry.f:1813
void exitt()
Definition: comm_mpi.f:604
subroutine exitti(stringi, idata)
Definition: comm_mpi.f:535
subroutine bcast(buf, len)
Definition: comm_mpi.f:431
subroutine err_chk(ierr, string)
Definition: comm_mpi.f:558
integer function gllel(ieg)
Definition: dprocmap.f:183
integer function gllnid(ieg)
Definition: dprocmap.f:161
logical function ifgtil(IVALUE, LINE)
Definition: ic.f:1699
subroutine capit(lettrs, n)
Definition: ic.f:1648
integer function indx1(S1, S2, L2)
Definition: ic.f:1357
subroutine blank(A, N)
Definition: math.f:19
function iglmax(a, n)
Definition: math.f:913
subroutine izero(a, n)
Definition: math.f:215
subroutine ltrue(a, n)
Definition: math.f:237
subroutine rzero(a, n)
Definition: math.f:208
subroutine rdout
Definition: reader_rea.f:1030
subroutine rdmesh
Definition: reader_rea.f:519
subroutine rdparam
Definition: reader_rea.f:3
subroutine rdcurve
Definition: reader_rea.f:609
subroutine rdbdry
Definition: reader_rea.f:702
subroutine rdicdf
Definition: reader_rea.f:881
subroutine rdobj
Definition: reader_rea.f:1114
subroutine rdhist
Definition: reader_rea.f:997
subroutine rdmatp
Definition: reader_rea.f:948