KTH framework for Nek5000 toolboxes; testing version  0.0.1
reader_par.f
Go to the documentation of this file.
1 c-----------------------------------------------------------------------
2  subroutine readat_par
3 C
4 C Read in run parameters from .par file
5 C
6  include 'SIZE'
7  include 'INPUT'
8  include 'RESTART'
9  include 'PARALLEL'
10  include 'CTIMER'
11 c
12  logical ifbswap
13 
14  call setdefaultparam
15 
16  if(nid.eq.0) call par_read(ierr)
17  call bcast(ierr,isize)
18  if(ierr .ne. 0) call exitt
19  call bcastparam
20 
21  call usrdat0
22 
23  call read_re2_hdr(ifbswap, .true.)
24 
25  call chkparam
26 
27  call mapelpr ! read .map file, est. gllnid, etc.
28 
29  call read_re2_data(ifbswap)
30 
31  call nekgsync()
32 
33  return
34  end
35 c-----------------------------------------------------------------------
36  subroutine setdefaultparam
37 C
38  include 'SIZE'
39  include 'INPUT'
40  include 'RESTART'
41  include 'PARALLEL'
42  include 'CTIMER'
43 
44  loglevel = 1
45  optlevel = 1
46 
47  call rzero(param,200)
48  call rzero(uparam,20)
49 
50  param(10) = 0 ! stop at numSteps
51  param(14) = 0 ! iostep
52  param(15) = 0 ! iotime
53 
54  param(21) = 1e-5 ! pressure tolerance
55  param(22) = 1e-7 ! velocity tolerance
56 
57  param(26) = 0.5 ! target Courant number
58  param(27) = 2 ! 2nd order in time
59  param(28) = 0 ! use same torder for mesh solver
60 
61  param(31) = 0 ! zero perturbations
62  param(32) = 0 ! read all BC from .re2
63  param(33) = 0 ! use default field index for BCs to read from .re2
64 
65  param(40) = 0 ! XXT
66 
67  param(41) = 0 ! additive SEMG
68  param(42) = 0 ! GMRES for iterative solver w/ nonsymmetric weighting
69  param(43) = 0 ! additive multilevel scheme (requires param(42).eq.0)
70  param(44) = 0 ! base top-level additive Schwarz on restrictions of E
71 
72  param(47) = 0.4 ! viscosity for mesh elasticity solver
73 
74  param(59) = 1 ! No fast operator eval
75 
76  param(65) = 1 ! just one i/o node
77  param(66) = 6 ! write in binary
78  param(67) = 6 ! read in binary
79  param(93) = mxprev ! number of vectors for projection
80 
81  param(94) = 5 ! projection for helmholz solves (controled by ifprojfld)
82  param(95) = 0 ! projection for pressure solve
83  param(99) = 4 ! enable dealising
84 
85  param(101) = 0 ! no additional modes
86  param(103) = 0 ! filter weight
87  filtertype = 0 ! no filtering
88 
89  param(160) = 1 ! cvode use normal mode
90  param(161) = 2 ! cvode use stiff integration
91  param(162) = 0 ! cvode absolute tolerance
92  param(163) = -1 ! cvode realtive tolerance
93  param(164) = 100 ! cvode don't limit internal dt
94  param(165) = 1 ! cvode increment factor DQJ
95  param(166) = 0 ! cvode use default ratio linear/non-linear tolerances
96  param(167) = 0 ! cvode use no preconditioner
97 
98  restol(0) = param(22)
99  restol(1) = param(22)
100  do i=1,ldimt
101  restol(1+i) = param(22)
102  enddo
103 
104  iftmsh(0) = .false.
105  iftmsh(1) = .false.
106  do i=1,ldimt
107  iftmsh(1+i) = .false.
108  enddo
109 
110  ifxyo = .true.
111  ifvo = .false.
112  ifpo = .false.
113  ifto = .false.
114  do i=1,ldimt1
115  ifpsco(i) = .false.
116  enddo
117 
118  ifadvc(1) = .true.
119  do i=1,ldimt
120  ifadvc(i+1) = .true.
121  enddo
122 
123  iffilter(1) = .false.
124  do i=1,ldimt
125  iffilter(i+1) = .false.
126  enddo
127 
128  ifdgfld(0) = .false.
129  ifdgfld(1) = .false.
130  do i=1,ldimt
131  ifdgfld(i+1) = .false.
132  enddo
133 
134  ifdiff(1) = .true.
135  do i=1,ldimt
136  ifdiff(i+1) = .true.
137  enddo
138 
139  ifdeal(1) = .true.
140  do i=1,ldimt
141  ifdeal(i+1) = .true.
142  enddo
143 
144  do i=1,ldimt
145  idpss(i) = -1
146  enddo
147 
148  meshpartitioner=3 ! HYBRID (RSB+RCB)
149 
150  ifprojfld(0) = .false.
151  ifprojfld(1) = .false.
152  do i=1,ldimt
153  ifprojfld(1+i) = .false.
154  enddo
155 
156  ifflow = .false.
157  ifheat = .false.
158  iftran = .true.
159  ifaxis = .false.
160  ifstrs = .false.
161  iflomach = .false.
162  ifmvbd = .false.
163  ifchar = .false.
164  ifmhd = .false.
165  ifuservp = .false.
166  ifcyclic = .false.
167  ifusermv = .false.
168  ifmgrid = .true.
169  ifessr = .false.
170  ifreguo = .false.
171  ifbase = .true.
172  ifpert = .false.
173  ifaziv = .false.
174  ifmoab = .false.
175  ifcvode = .false.
176 
177  ifdg = .false.
178  ifsync = .false.
179  ifanls = .false.
180  ifcoup = .false.
181  ifvcoup = .false.
182  ifexplvis = .false.
183  ifcons = .false.
184 
185  ifmodel = .false.
186  ifkeps = .false.
187  ifschclob = .false.
188 
189  ifdp0dt = .false.
190  ifreguo = .false. ! dump on the GLL mesh
191 
192  fem_amg_param(1) = 0
193  crs_param(1) = 0
194 
195  call izero(matype,16*ldimt1)
196  call rzero(cpgrp ,48*ldimt1)
197 
198  call blank (initc,15*132)
199 
200  return
201  end
202 c-----------------------------------------------------------------------
203  subroutine par_read(ierr)
204 c
205 c parse .par file and set run parameters
206 c
207 c todo:
208 c - check for invalid values for a given key
209 c - print default value to screen
210 c - separate settings for tol, proj, dealiasing for ps
211 c - mhd support
212 
213  include 'SIZE'
214  include 'INPUT'
215  include 'ADJOINT'
216  include 'RESTART'
217  include 'PARALLEL'
218  include 'CTIMER'
219  include 'TSTEP'
220 
221  character*132 c_out,txt, txt2
222 
223  call finiparser_load(parfle,ierr)
224  if(ierr .ne. 0) return
225 
226  call par_verify(ierr)
227  if(ierr .ne. 0) return
228 
229 c set parameters
230  call finiparser_getdbl(d_out,'general:loglevel',ifnd)
231  if(ifnd .eq. 1) loglevel = d_out
232 
233  call finiparser_getdbl(d_out,'general:optlevel',ifnd)
234  if(ifnd .eq. 1) optlevel = d_out
235 
236  call finiparser_getstring(c_out,'general:stopAt',ifnd)
237  if (ifnd .eq. 1) then
238  call capit(c_out,132)
239  if (index(c_out,'ENDTIME') .eq. 1) then
240  call finiparser_getdbl(d_out,'general:endTime',ifnd)
241  if (ifnd .eq. 1) then
242  param(10) = d_out
243  else
244  write(6,*) 'general:endTime'
245  write(6,*) 'is required for general:stopAt = endTime!'
246  goto 999
247  endif
248  else if (index(c_out,'NUMSTEPS') .eq. 1) then
249  call finiparser_getdbl(d_out,'general:numSteps',ifnd)
250  if (ifnd .eq. 1) then
251  param(11) = d_out
252  else
253  write(6,*) 'general:numSteps'
254  write(6,*) 'is required for general:stopAt = numSteps!'
255  goto 999
256  endif
257  else
258  write(6,*) 'value: ',trim(c_out)
259  write(6,*) 'is invalid for general:stopAt!'
260  goto 999
261  endif
262  else
263  call finiparser_getdbl(d_out,'general:numSteps',ifnd)
264  if (ifnd .eq. 1) then
265  param(11) = d_out
266  endif
267  endif
268 
269  call finiparser_getdbl(d_out,'general:dt',ifnd)
270  if (ifnd .eq. 1) then
271  param(12) = d_out
272  endif
273 
274  param(12) = -1*abs(param(12))
275  call finiparser_getbool(i_out,'general:variableDt',ifnd)
276  if (ifnd .eq. 1) then
277  if (i_out .eq. 1) then
278  param(12) = abs(param(12))
279  call finiparser_getdbl(d_out,'general:targetCFL',ifnd)
280  if (ifnd .eq. 1) then
281  param(26) = d_out
282  else
283  write(6,*) 'general:targetCFL'
284  write(6,*) 'is required for general:variableDt!'
285  goto 999
286  endif
287  endif
288  endif
289 
290  call finiparser_getdbl(d_out,'general:writeInterval',ifnd)
291  if (ifnd .eq. 1) param(15) = d_out
292  call finiparser_getstring(c_out,'general:writeControl',ifnd)
293  if (ifnd .eq. 1) then
294  call capit(c_out,132)
295  if (index(c_out,'RUNTIME') .eq. 1) then
296  param(14) = d_out
297  param(15) = 0
298  else if (index(c_out,'TIMESTEP') .eq. 1) then
299  param(14) = 0
300  param(15) = d_out
301  else if (index(c_out,'ELAPSEDTIME') .eq. 1) then
302  param(14) = d_out
303  param(15) = 0
304  timeioe = 1
305  else
306  write(6,*) 'value: ',trim(c_out)
307  write(6,*) 'is invalid for general:writeControl!'
308  goto 999
309  endif
310  endif
311 
312  call finiparser_getdbl(d_out,'pressure:residualTol',ifnd)
313  if(ifnd .eq. 1) param(21) = d_out
314  call finiparser_getdbl(d_out,'velocity:residualTol',ifnd)
315  if(ifnd .eq. 1) then
316  restol(1) = d_out
317  param(22) = d_out
318  endif
319 
320  call finiparser_find(i_out,'temperature',ifnd)
321  if(ifnd .eq. 1) then
322  ifheat = .true.
323  ifto = .true.
324  idpss(1) = 0 ! Helmholtz is default
325  endif
326 
327  j = 0
328  do i = 1,99
329  write(txt,"('scalar',i2.2)") i
330  call finiparser_find(i_out,txt,ifnd)
331  if (ifnd .eq. 1) j = j + 1
332  enddo
333  if (j.gt.ldimt-1) then
334  write(6,*) 'found more scalars than specified in SIZE!'
335  goto 999
336  endif
337 
338  j = 0
339  do i = 1,ldimt-1
340  write(txt,"('scalar',i2.2)") i
341  call finiparser_find(i_out,txt,ifnd)
342  if (ifnd .eq. 1) then
343  j = j + 1
344  ifpsco(i) = .true.
345  idpss(i+1) = 0 ! Helmholtz is default
346  endif
347  enddo
348  param(23) = j ! number of scalars
349  n = param(23)
350 
351  is = 2
352  if (ifheat) is = 1
353 
354  do i = is,n+1
355 
356  if (i.eq.1) then
357  txt = 'temperature'
358  else
359  write(txt,"('scalar',i2.2)") i-1
360  endif
361 
362  call finiparser_getstring(c_out,trim(txt)//':solver',ifnd)
363  call capit(c_out,132)
364  if(ifnd .eq. 1) then
365  if (index(c_out,'CVODE') .eq. 1) then
366  idpss(i) = 1
367  call finiparser_getdbl(d_out,trim(txt)//':absoluteTol',ifnd)
368  if (ifnd .eq. 1) then
369  atol(i+1) = d_out
370  else
371  write(6,*) trim(txt) // ':absoluteTol'
372  write(6,*) 'is required for ',trim(txt)//':solver = CVODE'
373  goto 999
374  endif
375  else if (index(c_out,'HELM') .eq. 1) then
376  continue
377  else if (index(c_out,'NONE') .eq. 1) then
378  idpss(i) = -1
379  else
380  write(6,*) 'value: ',trim(c_out)
381  write(6,*) 'is invalid for ',trim(txt)//':solver!'
382  goto 999
383  endif
384  else
385  call finiparser_getdbl(d_out,trim(txt)//':residualTol',ifnd)
386  if (ifnd .eq. 1) then
387  restol(i+1) = d_out
388  endif
389  endif
390  call finiparser_getbool(i_out,trim(txt)//':residualProj',ifnd)
391  if (ifnd .eq. 1) then
392  ifprojfld(i+1) = .false.
393  if(i_out .eq. 1) ifprojfld(i+1) = .true.
394  endif
395 
396  enddo
397 
398  call finiparser_getdbl(d_out,'cvode:absoluteTol',ifnd)
399  if(ifnd .eq. 1) param(162) = d_out
400  call finiparser_getdbl(d_out,'cvode:relativeTol',ifnd)
401  if(ifnd .eq. 1) param(163) = d_out
402  call finiparser_getdbl(d_out,'cvode:dtmax',ifnd)
403  if(ifnd .eq. 1) param(164) = d_out
404  call finiparser_getdbl(d_out,'cvode:DQJincrementFactor',ifnd)
405  if(ifnd .eq. 1) param(165) = d_out
406  call finiparser_getdbl(d_out,'cvode:ratioLNLtol',ifnd)
407  if(ifnd .eq. 1) param(166) = d_out
408 
409  call finiparser_getstring(c_out,'cvode:preconditioner',ifnd)
410  call capit(c_out,132)
411  if(ifnd .eq. 1) then
412  if (index(c_out,'USER') .eq. 1) then
413  param(167) = 1
414  else if (index(c_out,'NONE') .eq. 1) then
415  param(167) = 0
416  else
417  write(6,*) 'value: ',trim(c_out)
418  write(6,*) 'is invalid for cvode:preconditioner!'
419  goto 999
420  endif
421  endif
422 
423  j = 0
424  do i = 1,ldimt
425  if (idpss(i).ge.0) j = j + 1
426  enddo
427  if (j .ge. 1) then ! we have to solve for temp and/or ps
428  ifheat = .true.
429  else
430  ifheat = .false.
431  endif
432 
433  call finiparser_getdbl(d_out,'magnetic:viscosity',ifnd)
434  if(ifnd .eq. 1) param(29) = d_out
435  if(param(29).lt.0.0) param(29) = -1.0/param(29)
436 
437  call finiparser_getdbl(d_out,'mesh:numberOfBCFields',ifnd)
438  if(ifnd .eq. 1) param(32) = int(d_out)
439 
440  call finiparser_getdbl(d_out,'mesh:firstBCFieldIndex',ifnd)
441  if(ifnd .eq. 1) param(33) = int(d_out)
442 
443  call finiparser_getstring(c_out,'pressure:solver',ifnd)
444  if (ifnd .eq. 1) then
445  call capit(c_out,132)
446  if (index(c_out,'GMRES') .eq. 1) then
447  param(42) = 0
448  else if (index(c_out,'CGFLEX') .eq. 1) then
449  param(42) = 2
450  else
451  write(6,*) 'value: ',trim(c_out)
452  write(6,*) 'not supported for pressure:solver!'
453  goto 999
454  endif
455  endif
456 
457  call finiparser_getstring(c_out,'pressure:preconditioner',ifnd)
458  if (ifnd .eq. 1) then
459  call capit(c_out,132)
460  if (index(c_out,'SEMG_XXT') .eq. 1) then
461  param(40) = 0
462  else if (index(c_out,'SEMG_AMG_HYPRE') .eq. 1) then
463  param(40) = 2
464  else if (index(c_out,'SEMG_AMG') .eq. 1) then
465  param(40) = 1
466  else if (index(c_out,'FEM_AMG_HYPRE') .eq. 1) then
467  param(40) = 3
468  else
469  write(6,*) 'value: ',trim(c_out)
470  write(6,*) 'is invalid for pressure:preconditioner!'
471  goto 999
472  endif
473  endif
474 
475  call finiparser_getbool(i_out,'general:writeDoublePrecision',ifnd)
476  if(ifnd .eq. 1 .and. i_out .eq. 1) param(63) = 1
477 
478  call finiparser_getdbl(d_out,'general:writeNFiles',ifnd)
479  if(ifnd .eq. 1) param(65) = int(d_out)
480 
481  call finiparser_getbool(i_out,'velocity:residualProj',ifnd)
482  if(ifnd .eq. 1) then
483  ifprojfld(1) = .false.
484  if(i_out .eq. 1) ifprojfld(1) = .true.
485  endif
486 
487  call finiparser_getbool(i_out,'pressure:residualProj',ifnd)
488  if(ifnd .eq. 1) then
489  param(95) = 0
490  if(i_out .eq. 1) param(95) = 5
491  endif
492 
493  call finiparser_getbool(i_out,'general:dealiasing',ifnd)
494  if(ifnd .eq. 1 .and. i_out .eq. 0) param(99) = -1
495 
496 c filtering parameters
497  call finiparser_getstring(c_out,'general:filtering',ifnd)
498  if (ifnd .eq. 1) then
499 c stabilization type: none, explicit or hpfrt
500  call capit(c_out,132)
501  if (index(c_out,'NONE') .eq. 1) then
502  filtertype = 0
503  goto 101
504  else if (index(c_out,'EXPLICIT') .eq. 1) then
505  filtertype = 1
506  call ltrue(iffilter,size(iffilter))
507  else if (index(c_out,'HPFRT') .eq. 1) then
508  filtertype = 2
509  call ltrue(iffilter,size(iffilter))
510  else
511  write(6,*) 'value: ',c_out
512  write(6,*) 'is invalid for general:filtering!'
513  goto 999
514  endif
515  call finiparser_getdbl(d_out,'general:filterWeight',ifnd)
516  if (ifnd .eq. 1) then
517  param(103) = d_out
518  else
519  write(6,*) 'general:filterWeight'
520  write(6,*) 'is required for general:filtering!'
521  goto 999
522  endif
523  call finiparser_getdbl(d_out,'general:filterModes',ifnd)
524  if (ifnd .eq. 1) then
525  param(101) = int(d_out) - 1
526  if (int(param(101)).eq.0) filtertype = 0
527  else
528  call finiparser_getdbl
529  $ (d_out,'general:filterCutoffRatio',ifnd)
530  if (ifnd .eq. 1) then
531  dtmp = anint(lx1*(1.0 - d_out))
532  param(101) = max(dtmp-1,0.0)
533  if (abs(1.0 - d_out).lt.0.01) filtertype = 0
534  else
535  write(6,*) 'general:filterCutoffRatio or filterModes'
536  write(6,*) 'is required for general:filtering!'
537  goto 999
538  endif
539  endif
540  101 continue
541  endif
542 
543  call finiparser_getstring(c_out,'cvode:mode',ifnd)
544  call capit(c_out,132)
545  if (index(c_out,'NORMAL') .eq. 1) param(160) = 1
546  if (index(c_out,'NORMAL_TSTOP' ) .eq. 1) param(160) = 3
547 
548  do i = 1,20
549  call blank(txt,132)
550  write(txt,"('general:userParam',i2.2)") i
551  call finiparser_getdbl(d_out,txt,ifnd)
552  if(ifnd .eq. 1) uparam(i) = d_out
553  enddo
554 
555 c set logical flags
556  call finiparser_getstring(c_out,'general:timeStepper',ifnd)
557  if (ifnd .eq. 1) then
558  call capit(c_out,132)
559 
560  if (index(c_out,'BDF1') .eq. 1) then
561  param(27) = 1
562  else if (index(c_out,'BDF2') .eq. 1) then
563  param(27) = 2
564  else if (index(c_out,'BDF3') .eq. 1) then
565  param(27) = 3
566  else
567  write(6,*) 'value: ',trim(c_out)
568  write(6,*) 'is invalid for general:timeStepper!'
569  goto 999
570  endif
571  endif
572 
573  call finiparser_getstring(c_out,'general:extrapolation',ifnd)
574  if (ifnd .eq. 1) then
575  call capit(c_out,132)
576  if (index(c_out,'OIFS') .eq. 1) then
577  ifchar = .true.
578 
579  call finiparser_getdbl(d_out,'general:targetCFL',ifnd)
580  if (ifnd .eq. 1) then
581  param(26) = d_out
582  else
583  write(6,*) 'general:targetCFL'
584  write(6,*) 'is required for general:extrapolation!'
585  goto 999
586  endif
587  else if (index(c_out,'STANDARD') .eq. 1) then
588  continue
589  else
590  write(6,*) 'value: ',trim(c_out)
591  write(6,*) 'is invalid for general:extrapolation!'
592  goto 999
593  endif
594  endif
595 
596  call finiparser_find(i_out,'velocity',ifnd)
597  if(ifnd .eq. 1) then
598  ifflow = .true.
599  ifvo = .true.
600  ifpo = .true.
601  endif
602 
603  call finiparser_getstring(c_out,'mesh:motion',ifnd)
604  if (ifnd .eq. 1) then
605  call capit(c_out,132)
606  if (index(c_out,'ELASTICITY') .eq. 1) then
607  ifmvbd = .true.
608  call finiparser_getdbl(d_out,'mesh:viscosity',ifnd)
609  if (ifnd .eq. 1) param(47) = d_out
610  call finiparser_getbool(i_out,'mesh:residualProj',ifnd)
611  if (ifnd .eq. 1) then
612  ifprojfld(0) = .false.
613  if(i_out .eq. 1) ifprojfld(0) = .true.
614  endif
615  else if (index(c_out,'USER') .eq. 1) then
616  ifmvbd = .true.
617  ifusermv = .true.
618  else if (index(c_out,'NONE') .eq. 1) then
619  continue
620  else
621  write(6,*) 'value: ',trim(c_out)
622  write(6,*) 'is invalid for mesh:motion!'
623  goto 999
624  endif
625  endif
626  call finiparser_getdbl(d_out,'mesh:residualTol',ifnd)
627  if(ifnd .eq. 1) restol(0) = d_out
628 
629  call finiparser_getbool(i_out,'problemType:axiSymmetry',ifnd)
630  if(ifnd .eq. 1) then
631  ifaxis = .false.
632  if(i_out .eq. 1) ifaxis = .true.
633  endif
634 
635  call finiparser_getbool(i_out,'problemType:swirl',ifnd)
636  if(ifnd .eq. 1) then
637  ifaziv = .false.
638  if(i_out .eq. 1) ifaziv = .true.
639  endif
640 
641  call finiparser_getbool(i_out,'problemType:cyclicBoundaries',ifnd)
642  if(ifnd.eq.1) then
643  ifcyclic = .false.
644  if(i_out .eq. 1) ifcyclic = .true.
645  endif
646 
647  call finiparser_getstring(c_out,'problemType:equation',ifnd)
648  call capit(c_out,132)
649  if (index(c_out,'STEADYSTOKES').eq.1) then
650  iftran = .false.
651  ifadvc(1) = .false.
652  else if (index(c_out,'STOKES').eq.1) then
653  ifadvc(1) = .false.
654  else if (index(c_out,'STEADYHEAT').eq.1) then
655  iftran = .false.
656  ifflow = .false.
657  ifheat = .true.
658  ifadvc(2) = .false.
659  else if (index(c_out,'LOWMACHNS').eq.1) then
660  iflomach = .true.
661  else if (index(c_out,'INCOMPLINNS').eq.1 .or.
662  $ index(c_out,'INCOMPLINADJNS').eq.1) then
663  ifpert = .true.
664  if (index(c_out,'INCOMPLINADJNS').eq.1) ifadj = .true.
665  call finiparser_getdbl
666  $ (d_out,'problemType:numberOfPerturbations',ifnd)
667  if (ifnd .eq. 1) then
668  param(31) = int(d_out)
669  else
670  param(31) = 1
671  endif
672  call finiparser_getbool
673  $ (i_out,'problemType:solveBaseFlow',ifnd)
674  if (ifnd .eq. 1) then
675  ifbase = .false.
676  if(i_out .eq. 1) ifbase = .true.
677  else
678  write(6,*) 'problemType:solveBaseFlow'
679  write(6,*) 'is required for ', trim(c_out)
680  goto 999
681  endif
682  else if (index(c_out,'COMPNS') .eq. 1) then
683 #ifdef CMTNEK
684  continue
685 #else
686  write(6,*) 'value: ',trim(c_out)
687  write(6,*) 'not supported for problemType:equation!'
688  write(6,*) 'Recompile with CMTNEK ...'
689  goto 999
690 #endif
691  else if (index(c_out,'INCOMPMHD') .eq. 1) then
692  write(6,*) 'value: ',trim(c_out)
693  write(6,*) 'not yet supported for problemType:equation!'
694  goto 999
695  endif
696 
697  call finiparser_getbool(i_out,
698  & 'problemType:stressFormulation',ifnd)
699  if(ifnd .eq. 1) then
700  ifstrs = .false.
701  if(i_out .eq. 1) ifstrs = .true.
702  endif
703 
704  call finiparser_getbool(i_out,
705  & 'problemType:variableProperties',ifnd)
706  if(ifnd .eq. 1) then
707  ifuservp = .false.
708  if(i_out .eq. 1) ifuservp = .true.
709  endif
710 
711  call finiparser_getbool(i_out,'problemType:dp0dt',ifnd)
712  if(ifnd .eq. 1) then
713  if(i_out .eq. 1) ifdp0dt = .true.
714  endif
715 
716  call finiparser_getbool(i_out,'cvode:stiff',ifnd)
717  if(ifnd .eq. 1) then
718  param(161) = 1 ! AM
719  if(i_out .eq. 1) param(161) = 2 !BDF
720  endif
721 
722 c set advection
723  call finiparser_getbool(i_out,'velocity:advection',ifnd)
724  if(ifnd .eq. 1) then
725  ifadvc(1) = .false.
726  if(i_out .eq. 1) ifadvc(1) = .true.
727  endif
728 
729  call finiparser_getbool(i_out,'temperature:advection',ifnd)
730  if(ifnd .eq. 1) then
731  ifadvc(2) = .false.
732  if(i_out .eq. 1) ifadvc(2) = .true.
733  endif
734 
735  do i = 1,ldimt-1
736  write(txt,"('scalar',i2.2,a)") i,':advection'
737  call finiparser_getbool(i_out,txt,ifnd)
738  if(ifnd .eq. 1) then
739  ifadvc(i+2) = .false.
740  if(i_out .eq. 1) ifadvc(i+2) = .true.
741  endif
742  enddo
743 
744 c set mesh-field mapping
745  call finiparser_getbool(i_out,'temperature:conjugateHeatTransfer',
746  & ifnd)
747  if(ifnd .eq. 1) then
748  if(i_out .eq. 1) then
749  iftmsh(0) = .true.
750  iftmsh(2) = .true.
751  endif
752  endif
753 
754  do i = 1,ldimt-1
755  write(txt,"('scalar',i2.2,a)") i,':conjugateHeatTransfer'
756  call finiparser_getbool(i_out,txt,ifnd)
757  if(ifnd .eq. 1) then
758  iftmsh(i+2) = .false.
759  if(i_out .eq. 1) iftmsh(i+2) = .true.
760  endif
761  enddo
762 
763 c set output flags
764  call finiparser_getbool(i_out,'mesh:writeToFieldFile',ifnd)
765  if(ifnd .eq. 1) then
766  ifxyo = .false.
767  if(i_out .eq. 1) ifxyo = .true.
768  endif
769 
770  call finiparser_getbool(i_out,'velocity:writeToFieldFile',ifnd)
771  if(ifnd .eq. 1) then
772  ifvo = .false.
773  if(i_out .eq. 1) ifvo = .true.
774  endif
775 
776  call finiparser_getbool(i_out,'pressure:writeToFieldFile',ifnd)
777  if(ifnd .eq. 1) then
778  ifpo = .false.
779  if(i_out .eq. 1) ifpo = .true.
780  endif
781 
782  call finiparser_getbool(i_out,'temperature:writeToFieldFile',ifnd)
783  if(ifnd .eq. 1) then
784  ifto = .false.
785  if(i_out .eq. 1) ifto = .true.
786  endif
787 
788  do i = 1,ldimt-1
789  write(txt,"('scalar',i2.2,a)") i,':writeToFieldFile'
790  call finiparser_getbool(i_out,txt,ifnd)
791  if(ifnd .eq. 1) then
792  ifpsco(i) = .false.
793  if(i_out .eq. 1) ifpsco(i) = .true.
794  endif
795  enddo
796 
797 c set properties
798  call finiparser_getdbl(d_out,'velocity:viscosity',ifnd)
799  if(ifnd .eq. 1) cpfld(1,1) = d_out
800  if (cpfld(1,1) .lt.0.0) cpfld(1,1) = -1.0/cpfld(1,1)
801  call finiparser_getdbl(d_out,'velocity:density',ifnd)
802  if(ifnd .eq. 1) cpfld(1,2) = d_out
803 
804  call finiparser_getdbl(d_out,'temperature:conductivity',ifnd)
805  if(ifnd .eq. 1) cpfld(2,1) = d_out
806  if (cpfld(2,1) .lt.0.0) cpfld(2,1) = -1.0/cpfld(2,1)
807  call finiparser_getdbl(d_out,'temperature:rhoCp',ifnd)
808  if(ifnd .eq. 1) cpfld(2,2) = d_out
809 
810  do i = 1,ldimt-1
811  write(txt,"('scalar',i2.2,a)") i,':diffusivity'
812  call finiparser_getdbl(d_out,txt,ifnd)
813  if(ifnd .eq. 1) cpfld(2+i,1) = d_out
814  if(cpfld(2+i,1) .lt.0.0) cpfld(2+i,1) = -1.0/cpfld(2+i,1)
815  write(txt,"('scalar',i2.2,a)") i,':density'
816  call finiparser_getdbl(d_out,txt,ifnd)
817  if(ifnd .eq. 1) cpfld(2+i,2) = d_out
818  enddo
819 
820 c set restart options
821  call finiparser_findtokens('general:startFrom', ',' , ifnd)
822  do i = 1,min(ifnd,15)
823  call finiparser_gettoken(initc(i),i)
824  if(index(initc(i),'0') .eq. 1) call blank(initc(i),132)
825  enddo
826 
827 c set partitioner options
828  call finiparser_getstring(c_out,'mesh:partitioner',ifnd)
829  call capit(c_out,132)
830  if(index(c_out,'RSB').eq.1) then
831  meshpartitioner=1
832  else if(index(c_out,'RCB').eq.1) then
833  meshpartitioner=2
834  else if (index(c_out,'RCBRSB').eq.1) then
835  meshpartitioner=3
836  else if (index(c_out,'METIS').eq.1) then
837  meshpartitioner=4
838  endif
839 
840 100 if(ierr.eq.0) call finiparser_dump()
841  return
842 
843 c error handling
844  999 continue
845  ierr = 1
846  goto 100
847 
848  end
849 c-----------------------------------------------------------------------
850  subroutine bcastparam
851 C
852 C Broadcast run parameters to all processors
853 C
854  include 'SIZE'
855  include 'INPUT'
856  include 'TSTEP'
857  include 'RESTART'
858  include 'PARALLEL'
859  include 'CTIMER'
860  include 'ADJOINT'
861  include 'CVODE'
862 
863  call bcast(loglevel, isize)
864  call bcast(optlevel, isize)
865 
866  call bcast(param , 200*wdsize)
867  call bcast(uparam, 20*wdsize)
868 
869  call bcast(filtertype, wdsize)
870 
871  call bcast(atol , (ldimt1+1)*wdsize)
872  call bcast(restol, (ldimt1+1)*wdsize)
873 
874  call bcast(ifchar , lsize)
875  call bcast(iftran , lsize)
876  call bcast(ifflow , lsize)
877  call bcast(ifheat , lsize)
878  call bcast(iflomach, lsize)
879  call bcast(ifstrs , lsize)
880  call bcast(ifmvbd , lsize)
881  call bcast(ifusermv, lsize)
882  call bcast(ifdp0dt, lsize)
883  call bcast(ifaxis , lsize)
884  call bcast(ifcyclic, lsize)
885  call bcast(ifmhd , lsize)
886  call bcast(ifuservp, lsize)
887  call bcast(ifpert, lsize)
888  call bcast(ifbase, lsize)
889  call bcast(ifmoab, lsize)
890  call bcast(ifaziv, lsize)
891  call bcast(ifadj , lsize)
892 
893  call bcast(ifadvc , ldimt1*lsize)
894  call bcast(ifdiff , ldimt1*lsize)
895  call bcast(ifdeal , ldimt1*lsize)
896  call bcast(iffilter, ldimt1*lsize)
897 
898  call bcast(idpss , ldimt*isize)
899 
900  call bcast(meshpartitioner,isize)
901 
902  call bcast(iftmsh , (ldimt1+1)*lsize)
903  call bcast(ifprojfld, (ldimt1+1)*lsize)
904 
905  call bcast(cpfld, 3*ldimt1*wdsize)
906 
907  call bcast(ifxyo , lsize)
908  call bcast(ifvo , lsize)
909  call bcast(ifpo , lsize)
910  call bcast(ifto , lsize)
911  call bcast(ifpsco, ldimt1*lsize)
912 
913  call bcast(initc, 15*132*csize)
914 
915  call bcast(timeioe,sizeof(timeioe))
916 
917 c set some internals
918  if (ldim.eq.3) if3d=.true.
919  if (ldim.ne.3) if3d=.false.
920 
921  param(1) = cpfld(1,2)
922  param(2) = cpfld(1,1)
923  param(7) = cpfld(2,2)
924  param(8) = cpfld(2,1)
925 
926  npscal=int(param(23)) ! number of scalar sections
927  npscl1=npscal+1
928  npscl2=npscal+2
929 
930  npert = abs(param(31))
931 
932  if (if3d) ifaxis = .false.
933  if (ifaxis) param(99) = 3 ! For now, at least.
934 
935  ifldmhd = npscal + 3
936  if (ifmhd) then
937  ifessr = .true.
938  ifchar = .false.
939  npscl1 = npscl1 + 1
940  cpfld(ifldmhd,1) = param(29) ! magnetic viscosity
941  cpfld(ifldmhd,2) = param( 1) ! magnetic rho same as for fluid
942  endif
943 
944  if (ifaxis.and..not.ifsplit) then ! Use standard Schwarz/PCG solver
945  ifmgrid = .false.
946  param(42) = 1.00000 ! p042 0=gmres/1=pcg
947  param(43) = 1.00000 ! p043 0=semg/1=schwarz
948  param(44) = 1.00000 ! p044 0=E-based/1=A-based prec.
949  endif
950 
951  if (.not.iftran) then
952  if (ifflow.and.ifsplit) then
953  iftran=.true. ! no steady Stokes support for Pn/Pn
954  else
955  param(11) = 1.0
956  param(12) = 1.0
957  param(19) = 0.0
958  endif
959  endif
960 
961  cv_nfld = 0
962  do i = 1,ldimt
963  if (idpss(i) .gt. 0) then
964  cv_nfld = cv_nfld + 1
965  ifcvfld(i+1) = .true.
966  endif
967  enddo
968  if (cv_nfld.gt.0) ifcvode = .true.
969 
970  return
971  END
972 c-----------------------------------------------------------------------
973  subroutine chkparam
974  include 'SIZE'
975  include 'INPUT'
976  include 'RESTART'
977  include 'PARALLEL'
978  include 'CTIMER'
979 c
980  neltmx=np*lelt
981  nelvmx=np*lelv
982 
983  neltmx=min(neltmx,lelg)
984  nelvmx=min(nelvmx,lelg)
985 
986  nelgt = iglmax(nelgt,1)
987  nelgv = iglmax(nelgv,1)
988 
989  if (npscal+1.gt.ldimt) then
990  if(nid.eq.0) then
991  write(6,21) ldimt,npscal+1
992  21 format(//,2x,'Error: Nek has been compiled'
993  $ /,2x,' for max.',i4,' scalars. This run'
994  $ /,2x,' requires that ldimt be set to',i4,'.')
995  endif
996  call exitt
997  endif
998 
999  if (nelgt.gt.neltmx.or.nelgv.gt.nelvmx) then
1000  if (nid.eq.0) then
1001  lelt_needed = nelgt/np
1002  if (mod(nelgt,np).ne.0) lelt_needed = lelt_needed + 1
1003  write(6,82) lelt,lelg,lelt_needed,np,nelgt
1004  82 format(//,2x,'Error: Nek has has been compiled for'
1005  $ ,/,2x, ' number of elements/proc (lelt):',i12
1006  $ ,/,2x, ' total number of elements (lelg):',i12
1007  $ ,/,2x
1008  $ ,/,2x,'This run requires:'
1009  $ ,/,2x,' lelt >= ',i12,' for np = ',i12
1010  $ ,/,2x,' lelg >= ',i12,/)
1011  endif
1012  call exitt
1013  endif
1014 
1015  if (ifmvbd) then
1016  if (lx1.ne.lx1m.or.ly1.ne.ly1m.or.lz1.ne.lz1m)
1017  $ call exitti
1018  $ ('Error: Mesh motion requires lx1m=lx1 etc. in SIZE . $',lx1m)
1019  endif
1020 
1021  IF(ldimr.NE.ldim) THEN
1022  IF(nid.EQ.0) THEN
1023  WRITE(6,10) ldim,ldimr
1024  10 FORMAT(//,2x,'Error: Nek has been compiled'
1025  $ /,2x,' for spatial dimension equal to',i2,'.'
1026  $ /,2x,' The mesh file has dimension',i2,'.')
1027  ENDIF
1028  call exitt
1029  ENDIF
1030 
1031  if (lpert.lt.npert) then
1032  if(nid.eq.0) write(6,*)
1033  $ 'ERROR: Increase lpert in SIZE to', npert
1034  call exitt
1035  endif
1036 
1037  IF (npscal+1.GT.ldimt .AND. ifmhd) THEN
1038  if(nid.eq.0) then
1039  WRITE(6,22) ldimt,npscal+1
1040  22 FORMAT(/s,2x,'Error: Nek has been compiled'
1041  $ /,2x,' for',i4,' scalars. A MHD run'
1042  $ /,2x,' requires that LDIMT be set to',i4,'.')
1043  endif
1044  call exitt
1045  ENDIF
1046 
1047  if (if3d) then
1048  if (ly1.ne.lx1.or.lz1.ne.lx1) then
1049  if (nid.eq.0) write(6,13) lx1,ly1,lz1
1050  13 format('ERROR: lx1,ly1,lz1:',3i5,' must be equal for 3D')
1051  call exitt
1052  endif
1053  if (ly2.ne.lx2.or.lz2.ne.lx2) then
1054  if (nid.eq.0) write(6,14) lx2,ly2,lz2
1055  14 format('ERROR: lx2,ly2,lz2:',3i5,' must be equal for 3D')
1056  call exitt
1057  endif
1058  else
1059  if (ly1.ne.lx1.or.lz1.ne.1) then
1060  if (nid.eq.0) write(6,12) lx1,ly1,lz1
1061  12 format('ERROR: ',3i5,' must have lx1=ly1; lz1=1, for 2D')
1062  call exitt
1063  endif
1064  if (ly2.ne.lx2.or.lz2.ne.1) then
1065  if (nid.eq.0) write(6,11) lx2,ly2,lz2
1066  11 format('ERROR: ',3i5,' must have lx2=ly2; lz2=1, for 2D')
1067  call exitt
1068  endif
1069  endif
1070 
1071  if (lgmres.lt.5 .and. param(42).eq.0) then
1072  if(nid.eq.0) write(6,*)
1073  $ 'WARNING: lgmres might be too low!'
1074  endif
1075 
1076  if (ifsplit) then
1077  if (lx1.ne.lx2) then
1078  if (nid.eq.0) write(6,43) lx1,lx2
1079  43 format('ERROR: lx1,lx2:',2i4,' must be equal for IFSPLIT=T')
1080  call exitt
1081  endif
1082  else
1083  if (lx2.ne.lx1-2) then
1084  if (nid.eq.0) write(6,44) lx1,lx2
1085  44 format('ERROR: lx1,lx2:',2i4,' lx2 must be lx-2 for IFSPLIT=F')
1086  call exitt
1087  endif
1088  endif
1089 
1090  if (param(40).eq.3 .and. .not.ifsplit) then
1091  call exitti
1092  $ ('ERROR: Selected preconditioner requires lx2=lx1$',lx2)
1093  endif
1094 
1095  if (ifsplit .and. ifuservp .and. .not.ifstrs) then
1096  if(nid.eq.0) write(6,*)
1097  $ 'Enable stress formulation to support PN/PN and IFUSERVP=T'
1098  ifstrs = .true.
1099  endif
1100 
1101  if (ifcyclic .and. .not.ifstrs) then
1102  if(nid.eq.0) write(6,*)
1103  $ 'Enable stress formulation to support cyclic BC'
1104  ifstrs = .true.
1105  endif
1106 
1107  ktest = (lx1-lx1m) + (ly1-ly1m) + (lz1-lz1m)
1108  if (ifstrs.and.ktest.ne.0) then
1109  if(nid.eq.0) write(6,*)
1110  $ 'ERROR: Stress formulation requires lx1m=lx1, etc. in SIZE'
1111  call exitt
1112  endif
1113 
1114  if (ifsplit .and. ifmhd) then
1115  if(nid.eq.0) write(6,*)
1116  $ 'ABORT: No MHD support for Pn-Pn'
1117  call exitt
1118  endif
1119 
1120  if (ifmhd .and. lbx1.ne.lx1) then
1121  if(nid.eq.0) write(6,*)
1122  $ 'ABORT: For MHD, need lbx1=lx1, etc.; Change SIZE '
1123  call exitt
1124  endif
1125 
1126  if (ifpert .and. lpx1.ne.lx1) then
1127  if(nid.eq.0) write(6,*)
1128  $ 'ABORT: For Lyapunov, need lpx1=lx1, etc.; Change SIZE '
1129  endif
1130 
1131  if (ifpert .and. ifsplit) then
1132  if(nid.eq.0) write(6,*)
1133  $ 'ABORT: For Lyapunov, need lx2=lx1-2, etc. in SIZE'
1134  endif
1135 
1136  if (iflomach .and. .not.ifsplit) then
1137  if(nid.eq.0) write(6,*)
1138  $ 'ABORT: For lowMach, need lx2=lx1, etc.; Change SIZE '
1139  call exitt
1140  endif
1141 
1142  if (iflomach .and. .not.ifheat) then
1143  if(nid.eq.0) write(6,*)
1144  $ 'ABORT: For lowMach, need to solve for temperature too!'
1145  call exitt
1146  endif
1147 
1148  if (ifdp0dt .and. .not.iflomach) then
1149  if(nid.eq.0) write(6,*)
1150  $ 'ABORT: Varying p0 requires lowMach! '
1151  call exitt
1152  endif
1153 
1154  if (ifchar .and. param(99).lt.0) then
1155  if (nid.eq.0) write(6,*)
1156  & 'ABORT: Characteristic scheme needs dealiasing!'
1157  call exitt
1158  endif
1159 
1160  if (.not.ifsplit .and. ifaxis .and. ifstrs) then
1161  if (nid.eq.0) write(6,*)
1162  $ 'ABORT: Axisymetric and stress formulation not supported ' //
1163  $ 'for PN/PN-2$'
1164  call exitt
1165  endif
1166 
1167  if (ifchar.and.(nelgv.ne.nelgt)) call exitti(
1168  $ 'ABORT: Characteristics not supported w/ conj. ht transfer$',1)
1169 
1170  if (param(99).gt.-1 .and. (lxd.lt.lx1 .or. lyd.lt.ly1 .or.
1171  & lzd.lt.lz1)) then
1172  if(nid.eq.0) write(6,*)
1173  & 'ABORT: Dealiasing space too small; Check lxd,lyd,lzd in SIZE '
1174  call exitt
1175  endif
1176 
1177  return
1178  end
1179 c-----------------------------------------------------------------------
1180  subroutine par_verify(ierr)
1181 
1182  include 'PARDICT'
1183 
1184 
1185  character*132 key
1186  character*1024 val
1187 
1188  character*132 txt
1189  character*1 tx1(132)
1190  equivalence(tx1,txt)
1191 
1192  ierr = 0
1193 
1194  call finiparser_getdictentries(n)
1195  do i = 1,n
1196  call finiparser_getpair(key,val,i,ifnd)
1197  call capit(key,132)
1198 
1199  is = index(key,'_') ! ignore user keys
1200  if (is.eq.1) goto 10
1201 
1202  do j = 1,pardict_nkeys ! do we find the key in the par-dictionary
1203  if(index(pardictkey(j),key).eq.1) goto 10
1204 
1205  is = index(key,'SCALAR')
1206  if(is .eq. 1) then
1207  call chcopy(txt,key,132)
1208  call chcopy(tx1(is+6),'%%',2)
1209  if(index(pardictkey(j),txt).eq.1) goto 10
1210  endif
1211  enddo
1212  write(6,*) 'ERROR: Par file contains unknown key ', key
1213  ierr = ierr + 1
1214  10 enddo
1215 
1216  return
1217  end
void exitt()
Definition: comm_mpi.f:604
subroutine nekgsync()
Definition: comm_mpi.f:502
subroutine exitti(stringi, idata)
Definition: comm_mpi.f:535
subroutine bcast(buf, len)
Definition: comm_mpi.f:431
#define finiparser_find
Definition: finiparser.c:33
#define finiparser_load
Definition: finiparser.c:28
#define finiparser_dump
Definition: finiparser.c:26
subroutine capit(lettrs, n)
Definition: ic.f:1648
subroutine mapelpr()
Definition: map2.f:3
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 chcopy(a, b, n)
Definition: math.f:281
subroutine rzero(a, n)
Definition: math.f:208
subroutine bcastparam
Definition: reader_par.f:851
subroutine chkparam
Definition: reader_par.f:974
subroutine setdefaultparam
Definition: reader_par.f:37
subroutine par_read(ierr)
Definition: reader_par.f:204
subroutine par_verify(ierr)
Definition: reader_par.f:1181
subroutine readat_par
Definition: reader_par.f:3
subroutine read_re2_hdr(ifbswap, ifverbose)
Definition: reader_re2.f:772
subroutine read_re2_data(ifbswap)
Definition: reader_re2.f:3