KTH framework for Nek5000 toolboxes; testing version  0.0.1
rprm.f
Go to the documentation of this file.
1 
6 !=======================================================================
9  subroutine rprm_register
10  implicit none
11 
12  include 'SIZE'
13  include 'FRAMELP'
14  include 'RPRMD'
15 
16  ! local variables
17  integer itmp
18  integer il
19 
20  ! functions
21  integer frame_get_master
22 !-----------------------------------------------------------------------
23  rprm_pid0 = frame_get_master()
24 
25  ! amd seems to have big problem with block data for arrays
26  rprm_ifinit = .false.
27  rprm_pid0 = 0
28  rprm_sec_num = 0
29  rprm_sec_mpos = 0
30  do il = 1, rprm_sec_id_max
31  rprm_sec_id(il) = -1
32  rprm_sec_act(il) = .false.
33  rprm_sec_name(il) = rprm_blname
34  end do
35  rprm_par_num = 0
36  rprm_par_mpos = 0
37  do il = 1, rprm_par_id_size
38  rprm_par_id(il,1) = -1
39  end do
40  do il = 1, rprm_par_id_max
41  rprm_par_name(il) = rprm_blname
42  rprm_parv_int(il) = 0
43  rprm_parv_real(il) = 0.0
44  rprm_parv_log(il) = .false.
45  rprm_parv_str(il) = rprm_blname
46  end do
47 
48  ! check if the current module was already registered
49  call mntr_mod_is_name_reg(itmp,rprm_name)
50  if (itmp.gt.0) then
51  call mntr_warn(itmp,
52  $ 'module ['//trim(rprm_name)//'] already registered')
53  return
54  endif
55 
56  ! find parent module
57  call mntr_mod_is_name_reg(itmp,'FRAME')
58  if (itmp.le.0) then
59  itmp = 1
60  call mntr_abort(itmp,
61  $ 'parent module ['//'FRAME'//'] not registered')
62  endif
63 
64  ! register module
65  call mntr_mod_reg(rprm_id,itmp,rprm_name,'Runtime parameters')
66 
67  ! register and set active section
68  call rprm_sec_reg(rprm_lsec_id,rprm_id,'_'//adjustl(rprm_name),
69  $ 'Runtime parameter section for rprm module')
70  call rprm_sec_set_act(.true.,rprm_lsec_id)
71 
72  ! register parameters
73  call rprm_rp_reg(rprm_ifparf_id,rprm_lsec_id,'PARFWRITE',
74  $ 'Do we write runtime parameter file',rpar_log,0,
75  $ 0.0,.false.,' ')
76 
77  call rprm_rp_reg(rprm_parfnm_id,rprm_lsec_id,'PARFNAME',
78  $ 'Runtime parameter file name for output (without .par)',
79  $ rpar_str,0,0.0,.false.,'outparfile')
80 
81  return
82  end subroutine
83 !=======================================================================
86  subroutine rprm_init
87  implicit none
88 
89  include 'SIZE'
90  include 'FRAMELP'
91  include 'RPRMD'
92 
93  ! local variables
94  integer itmp
95  real rtmp
96  logical ltmp
97  character*20 ctmp
98  integer iunit, ierr
99  character*30 fname
100 !-----------------------------------------------------------------------
101  ! check if the module was already initialised
102  if (rprm_ifinit) then
103  call mntr_warn(rprm_id,
104  $ 'module ['//trim(rprm_name)//'] already initiaised.')
105  return
106  endif
107 
108  ! get runtime parameters
109  call rprm_rp_get(itmp,rtmp,ltmp,ctmp,rprm_ifparf_id,rpar_log)
110  rprm_ifparf = ltmp
111  call rprm_rp_get(itmp,rtmp,ltmp,ctmp,rprm_parfnm_id,rpar_str)
112  rprm_parfnm = ctmp
113 
114  ! write summary
115  iunit = 6
116  call rprm_rp_summary_print(iunit)
117 
118  ! save .par file
119  if (rprm_ifparf) then
120  call io_file_freeid(iunit, ierr)
121  if (ierr.eq.0) then
122  fname=trim(adjustl(rprm_parfnm))//'.par'
123  open(unit=iunit,file=fname,status='new',iostat=ierr)
124  if (ierr.eq.0) then
125  call rprm_rp_summary_print(iunit)
126  close (iunit)
127  else
128  call mntr_log(rprm_id,lp_inf,
129  $ 'ERROR: cannot open output .par file')
130  endif
131  else
132  call mntr_log(rprm_id,lp_inf,
133  $ 'ERROR: cannot allocate iunit for output .par file')
134  endif
135  endif
136 
137  ! everything is initialised
138  rprm_ifinit = .true.
139 
140  return
141  end subroutine
142 !=======================================================================
146  logical function rprm_is_initialised()
147  implicit none
148 
149  include 'SIZE'
150  include 'FRAMELP'
151  include 'RPRMD'
152 !-----------------------------------------------------------------------
153  rprm_is_initialised = rprm_ifinit
154 
155  return
156  end function
157 !=======================================================================
164  subroutine rprm_sec_reg(rpid,mid,pname,pdscr)
165  implicit none
166 
167  include 'SIZE'
168  include 'PARALLEL' ! ISIZE
169  include 'FRAMELP'
170  include 'RPRMD'
171 
172  ! argument list
173  integer rpid, mid
174  character*(*) pname, pdscr
175 
176  ! local variables
177  character*10 mname
178  character*20 lname
179  character*132 ldscr
180  character*200 llog
181  integer slen,slena
182 
183  integer il, ipos
184  integer lval
185 
186  ! functions
187  logical mntr_mod_is_id_reg
188 !-----------------------------------------------------------------------
189  ! check name length
190  slena = len_trim(adjustl(pname))
191  ! remove trailing blanks
192  slen = len_trim(pname) - slena + 1
193  if (slena.gt.rprm_lstl_mnm) then
194  call mntr_log(rprm_id,lp_deb,
195  $ 'too long section name; shortenning')
196  slena = min(slena,rprm_lstl_mnm)
197  endif
198  call blank(lname,rprm_lstl_mnm)
199  lname= pname(slen:slen+slena- 1)
200  call capit(lname,slena)
201 
202  ! check description length
203  slena = len_trim(adjustl(pdscr))
204  ! remove trailing blanks
205  slen = len_trim(pdscr) - slena + 1
206  if (slena.ge.rprm_lstl_mds) then
207  call mntr_log(rprm_id,lp_deb,
208  $ 'too long section description; shortenning')
209  slena = min(slena,rprm_lstl_mnm)
210  endif
211  call blank(ldscr,rprm_lstl_mds)
212  ldscr= pdscr(slen:slen + slena - 1)
213 
214  ! find empty space
215  ipos = 0
216 
217  ! to ensure consistency I do it on master and broadcast result
218  if (nid.eq.rprm_pid0) then
219 
220  ! check if parameter name is already registered
221  do il=1,rprm_sec_mpos
222  if (rprm_sec_id(il).gt.0.and.
223  $ rprm_sec_name(il).eq.lname) then
224  ipos = -il
225  exit
226  endif
227  enddo
228 
229  ! find empty spot
230  if (ipos.eq.0) then
231  do il=1,rprm_sec_id_max
232  if (rprm_sec_id(il).eq.-1) then
233  ipos = il
234  exit
235  endif
236  enddo
237  endif
238  endif
239 
240  ! broadcast ipos
241  call bcast(ipos,isize)
242 
243  ! error; no free space found
244  if (ipos.eq.0) then
245  rpid = ipos
246  call mntr_abort(rprm_id,
247  $ 'Section '//trim(lname)//' cannot be registered')
248  ! section already registered
249  elseif (ipos.lt.0) then
250  rpid = abs(ipos)
251  call mntr_abort(rprm_id,
252  $ 'Section '//trim(lname)//' is already registered')
253  ! new section
254  else
255  rpid = ipos
256  ! check if module is registered
257  if (mntr_mod_is_id_reg(mid)) then
258  rprm_sec_id(ipos) = mid
259  else
260  call mntr_abort(rprm_id,
261  $ "Sections's "//trim(lname)//" module not registered")
262  endif
263  rprm_sec_name(ipos)=lname
264  rprm_sec_dscr(ipos)=ldscr
265  rprm_sec_num = rprm_sec_num + 1
266  if (rprm_sec_mpos.lt.ipos) rprm_sec_mpos = ipos
267 
268  ! logging
269  call mntr_mod_get_info(mname,ipos,mid)
270  llog='Module ['//trim(mname)//'] registered section '
271  llog=trim(llog)//' '//trim(lname)//': '//trim(ldscr)
272  call mntr_log(rprm_id,lp_inf,trim(llog))
273  endif
274 
275  return
276  end subroutine
277 !=======================================================================
283  subroutine rprm_sec_is_name_reg(rpid,mid,pname)
284  implicit none
285 
286  include 'SIZE'
287  include 'PARALLEL' ! ISIZE
288  include 'FRAMELP'
289  include 'RPRMD'
290 
291  ! argument list
292  integer rpid, mid
293  character*(*) pname
294 
295  ! local variables
296  character*3 str
297  character*10 mname
298  character*20 lname
299  character*132 llog
300  integer slen,slena
301 
302  integer il, ipos
303 !-----------------------------------------------------------------------
304  ! check name length
305  slena = len_trim(adjustl(pname))
306  ! remove trailing blanks
307  slen = len_trim(pname) - slena + 1
308  if (slena.gt.rprm_lstl_mnm) then
309  call mntr_log(rprm_id,lp_deb,
310  $ 'too long section name; shortenning')
311  slena = min(slena,rprm_lstl_mnm)
312  endif
313  call blank(lname,rprm_lstl_mnm)
314  lname= pname(slen:slen+slena- 1)
315  call capit(lname,slena)
316 
317  ! find parameter
318  ipos = 0
319 
320  ! to ensure consistency I do it on master and broadcast result
321  if (nid.eq.rprm_pid0) then
322 
323  ! check if parameter name is already registered
324  do il=1,rprm_sec_mpos
325  if (rprm_sec_id(il).gt.0.and.
326  $ rprm_sec_name(il).eq.lname) then
327  ipos = il
328  exit
329  endif
330  enddo
331 
332  endif
333 
334  ! broadcast ipos
335  call bcast(ipos,isize)
336 
337  if (ipos.eq.0) then
338  rpid = -1
339  call mntr_log(rprm_id,lp_inf,
340  $ 'Section '//trim(lname)//' not registered')
341  else
342  rpid = ipos
343  write(str,'(I3)') ipos
344  call mntr_log(rprm_id,lp_vrb,
345  $ 'Section '//trim(lname)//' registered with id = '//trim(str))
346  ! check module
347  if (mid.ne.rprm_sec_id(ipos)) then
348  call mntr_log(rprm_id,lp_inf,
349  $ "Section's "//trim(lname)//" module inconsistent")
350  rpid = -1
351  endif
352  endif
353 
354  return
355  end subroutine
356 !=======================================================================
361  logical function rprm_sec_is_id_reg(rpid)
362  implicit none
363 
364  include 'SIZE'
365  include 'FRAMELP'
366  include 'RPRMD'
367 
368  ! argument list
369  integer rpid
370 !-----------------------------------------------------------------------
371  rprm_sec_is_id_reg = rprm_sec_id(rpid).gt.0
372 
373  return
374  end function
375 !=======================================================================
382  subroutine rprm_sec_get_info(pname,mid,ifact,rpid)
383  implicit none
384 
385  include 'SIZE'
386  include 'FRAMELP'
387  include 'RPRMD'
388 
389  ! argument list
390  integer rpid, mid
391  character*20 pname
392  logical ifact
393 
394  ! local variables
395  character*5 str
396 !-----------------------------------------------------------------------
397  if (rprm_sec_id(rpid).gt.0) then
398  pname = rprm_sec_name(rpid)
399  mid = rprm_sec_id(rpid)
400  ifact = rprm_sec_act(rpid)
401  else
402  write(str,'(I3)') rpid
403  call mntr_log(rprm_id,lp_inf,
404  $ 'Section id'//trim(str)//' not registered')
405  rpid = -1
406  endif
407 
408  return
409  end subroutine
410 !=======================================================================
421  subroutine rprm_sec_set_act(ifact,rpid)
422  implicit none
423 
424  include 'SIZE'
425  include 'PARALLEL'
426  include 'FRAMELP'
427  include 'RPRMD'
428 
429  ! argument list
430  integer rpid
431  logical ifact
432 
433  ! local variables
434  logical lval
435  character*5 str
436 !-----------------------------------------------------------------------
437  if (rprm_sec_id(rpid).gt.0) then
438  ! broadcast pval; to keep consistency
439  lval = ifact
440  call bcast(lval,lsize)
441  rprm_sec_act(rpid) = lval
442  else
443  write(str,'(I3)') rpid
444  call mntr_abort(rprm_id,
445  $ "Section "//trim(str)//" activation error")
446  endif
447 
448  return
449  end subroutine
450 !=======================================================================
455  logical function rprm_sec_is_id_act(rpid)
456  implicit none
457 
458  include 'SIZE'
459  include 'FRAMELP'
460  include 'RPRMD'
461 
462  ! argument list
463  integer rpid
464 !-----------------------------------------------------------------------
465  rprm_sec_is_id_act = rprm_sec_id(rpid).gt.0.and.
466  $ rprm_sec_act(rpid)
467 
468  return
469  end function
470 !=======================================================================
482  subroutine rprm_rp_reg(rpid,mid,pname,pdscr,ptype ,
483  $ ipval, rpval, lpval, cpval)
484  implicit none
485 
486  include 'SIZE'
487  include 'PARALLEL' ! ISIZE
488  include 'FRAMELP'
489  include 'RPRMD'
490 
491 
492  ! argument list
493  integer rpid, mid, ptype, ipval
494  real rpval
495  logical lpval
496  character*(*) pname, pdscr, cpval
497 
498  ! local variables
499  character*10 mname
500  character*20 lname
501  character*132 ldscr
502  character*200 llog
503  integer slen,slena
504 
505  integer il, ipos
506  integer ivall
507  real rvall
508  logical lvall
509  character*20 cvall
510 
511  ! functions
512  logical rprm_sec_is_id_reg
513 !-----------------------------------------------------------------------
514  ! check name length
515  slena = len_trim(adjustl(pname))
516  ! remove trailing blanks
517  slen = len_trim(pname) - slena + 1
518  if (slena.gt.rprm_lstl_mnm) then
519  call mntr_log(rprm_id,lp_deb,
520  $ 'too long parameter name; shortenning')
521  slena = min(slena,rprm_lstl_mnm)
522  endif
523  call blank(lname,rprm_lstl_mnm)
524  lname= pname(slen:slen+slena- 1)
525  call capit(lname,slena)
526 
527  ! check description length
528  slena = len_trim(adjustl(pdscr))
529  ! remove trailing blanks
530  slen = len_trim(pdscr) - slena + 1
531  if (slena.ge.rprm_lstl_mds) then
532  call mntr_log(rprm_id,lp_deb,
533  $ 'too long parameter description; shortenning')
534  slena = min(slena,rprm_lstl_mnm)
535  endif
536  call blank(ldscr,rprm_lstl_mds)
537  ldscr= pdscr(slen:slen + slena - 1)
538 
539  ! find empty space
540  ipos = 0
541 
542  ! to ensure consistency I do it on master and broadcast result
543  if (nid.eq.rprm_pid0) then
544 
545  ! check if parameter name is already registered
546  do il=1,rprm_par_mpos
547  if (rprm_par_id(rprm_par_mark,il).gt.0.and.
548  $ rprm_par_id(rprm_par_mark,il).eq.mid.and.
549  $ rprm_par_name(il).eq.lname) then
550  ipos = -il
551  exit
552  endif
553  enddo
554 
555  ! find empty spot
556  if (ipos.eq.0) then
557  do il=1,rprm_par_id_max
558  if (rprm_par_id(rprm_par_mark,il).eq.-1) then
559  ipos = il
560  exit
561  endif
562  enddo
563  endif
564  endif
565 
566  ! broadcast ipos
567  call bcast(ipos,isize)
568 
569  ! error; no free space found
570  if (ipos.eq.0) then
571  rpid = ipos
572  call mntr_abort(rprm_id,
573  $ 'Parameter '//trim(lname)//' cannot be registered')
574  ! parameter already registered
575  elseif (ipos.lt.0) then
576  rpid = abs(ipos)
577  call mntr_abort(rprm_id,
578  $ 'Parameter '//trim(lname)//' is already registered')
579  ! new parameter
580  else
581  rpid = ipos
582  ! check if section is registered
583  if (rprm_sec_is_id_reg(mid)) then
584  rprm_par_id(rprm_par_mark,ipos) = mid
585  else
586  call mntr_abort(rprm_id,
587  $ "Parameter's "//trim(lname)//" section not registered")
588  endif
589  rprm_par_id(rprm_par_type,ipos) = ptype
590  rprm_par_name(ipos)=lname
591  rprm_par_dscr(ipos)=ldscr
592  rprm_par_num = rprm_par_num + 1
593  if (rprm_par_mpos.lt.ipos) rprm_par_mpos = ipos
594 
595  ! broadcast pval; to keep consistency
596  if (ptype.eq.rpar_int) then
597  ivall = ipval
598  call bcast(ivall,isize)
599  rprm_parv_int(ipos) = ivall
600  elseif (ptype.eq.rpar_real) then
601  rvall = rpval
602  call bcast(rvall,wdsize)
603  rprm_parv_real(ipos) = rvall
604  elseif (ptype.eq.rpar_log) then
605  lvall = lpval
606  call bcast(lvall,lsize)
607  rprm_parv_log(ipos) = lvall
608  elseif (ptype.eq.rpar_str) then
609  ! check value length
610  slena = len_trim(adjustl(cpval))
611  ! remove trailing blanks
612  slen = len_trim(cpval) - slena + 1
613  if (slena.gt.rprm_lstl_mnm) then
614  call mntr_log(rprm_id,lp_deb,
615  $ 'too long parameter default value; shortenning')
616  slena = min(slena,rprm_lstl_mnm)
617  endif
618  call blank(cvall,rprm_lstl_mnm)
619  cvall= cpval(slen:slen+slena- 1)
620  ! broadcast pval; to keep consistency
621  call bcast(cvall,rprm_lstl_mnm*csize)
622  rprm_parv_str(ipos) = cvall
623  else
624  call mntr_abort(rprm_id,
625  $ "Parameter's "//trim(lname)//" wrong type")
626  endif
627 
628  ! logging
629  mname = trim(rprm_sec_name(mid))
630  llog='Section '//trim(mname)//' registered parameter '
631  llog=trim(llog)//' '//trim(lname)//': '//trim(ldscr)
632  call mntr_log(rprm_id,lp_inf,trim(llog))
633  if (ptype.eq.rpar_int) then
634  call mntr_logi(rprm_id,lp_vrb,
635  $ 'Default value '//trim(lname)//' = ',ivall)
636  elseif (ptype.eq.rpar_real) then
637  call mntr_logr(rprm_id,lp_vrb,
638  $ 'Default value '//trim(lname)//' = ',rvall)
639  elseif (ptype.eq.rpar_log) then
640  call mntr_logl(rprm_id,lp_vrb,
641  $ 'Default value '//trim(lname)//' = ',lvall)
642  elseif (ptype.eq.rpar_str) then
643  call mntr_log(rprm_id,lp_vrb,
644  $ 'Default value '//trim(lname)//' = '//trim(cvall))
645  endif
646  endif
647 
648  return
649  end subroutine
650 !=======================================================================
657  subroutine rprm_rp_is_name_reg(rpid,mid,pname,ptype)
658  implicit none
659 
660  include 'SIZE'
661  include 'PARALLEL' ! ISIZE
662  include 'FRAMELP'
663  include 'RPRMD'
664 
665  ! argument list
666  integer rpid, mid, ptype
667  character*(*) pname
668 
669  ! local variables
670  character*3 str
671  character*10 mname
672  character*20 lname
673  character*132 llog
674  integer slen,slena
675 
676  integer il, ipos
677 !-----------------------------------------------------------------------
678  ! check name length
679  slena = len_trim(adjustl(pname))
680  ! remove trailing blanks
681  slen = len_trim(pname) - slena + 1
682  if (slena.gt.rprm_lstl_mnm) then
683  call mntr_log(rprm_id,lp_deb,
684  $ 'too long parameter name; shortenning')
685  slena = min(slena,rprm_lstl_mnm)
686  endif
687  call blank(lname,rprm_lstl_mnm)
688  lname= pname(slen:slen+slena- 1)
689  call capit(lname,slena)
690 
691  ! find parameter
692  ipos = 0
693 
694  ! to ensure consistency I do it on master and broadcast result
695  if (nid.eq.rprm_pid0) then
696 
697  ! check if parameter name is already registered
698  do il=1,rprm_par_mpos
699  if (rprm_par_id(rprm_par_mark,il).gt.0.and.
700  $ rprm_par_name(il).eq.lname) then
701  ipos = il
702  exit
703  endif
704  enddo
705 
706  endif
707 
708  ! broadcast ipos
709  call bcast(ipos,isize)
710 
711  if (ipos.eq.0) then
712  rpid = -1
713  call mntr_log(rprm_id,lp_inf,
714  $ 'Parameter '//trim(lname)//' not registered')
715  else
716  rpid = ipos
717  write(str,'(I3)') ipos
718  call mntr_log(rprm_id,lp_vrb,
719  $ 'Parameter '//trim(lname)//' registered with id = '//trim(str))
720  ! check module
721  if (mid.ne.rprm_par_id(rprm_par_mark,ipos)) then
722  call mntr_log(rprm_id,lp_inf,
723  $ "Parameter's "//trim(lname)//" section inconsistent")
724  rpid = -1
725  endif
726  ! check type
727  if (ptype.ne.rprm_par_id(rprm_par_type,ipos)) then
728  call mntr_log(rprm_id,lp_inf,
729  $ "Parameter's "//trim(lname)//" type inconsistent")
730  rpid = -1
731  endif
732  endif
733 
734  return
735  end subroutine
736 !=======================================================================
742  logical function rprm_rp_is_id_reg(rpid,ptype)
743  implicit none
744 
745  include 'SIZE'
746  include 'FRAMELP'
747  include 'RPRMD'
748 
749  ! argument list
750  integer rpid, ptype
751 !-----------------------------------------------------------------------
752  rprm_rp_is_id_reg = rprm_par_id(rprm_par_mark,rpid).gt.0.and.
753  $ rprm_par_id(rprm_par_type,rpid).eq.ptype
754 
755  return
756  end function
757 !=======================================================================
764  subroutine rprm_rp_get_info(pname,mid,ptype,rpid)
765  implicit none
766 
767  include 'SIZE'
768  include 'FRAMELP'
769  include 'RPRMD'
770 
771  ! argument list
772  integer rpid, mid, ptype
773  character*20 pname
774 
775  ! local variables
776  character*5 str
777 !-----------------------------------------------------------------------
778  if (rprm_par_id(rprm_par_mark,rpid).gt.0) then
779  pname = rprm_par_name(rpid)
780  mid = rprm_par_id(rprm_par_mark,rpid)
781  ptype = rprm_par_id(rprm_par_type,rpid)
782  else
783  write(str,'(I3)') rpid
784  call mntr_log(rprm_id,lp_inf,
785  $ 'Parameter id'//trim(str)//' not registered')
786  rpid = -1
787  endif
788 
789  return
790  end subroutine
791 !=======================================================================
800  subroutine rprm_rp_set(rpid,ptype,ipval,rpval,lpval,cpval)
801  implicit none
802 
803  include 'SIZE'
804  include 'PARALLEL'
805  include 'FRAMELP'
806  include 'RPRMD'
807 
808  ! argument list
809  integer rpid, ptype
810  integer ipval
811  real rpval
812  logical lpval
813  character*(*) cpval
814 
815  ! local variables
816  integer ivall
817  real rvall
818  logical lvall
819  character*20 cvall
820  character*5 str
821  integer slen,slena
822 
823 !-----------------------------------------------------------------------
824  if (rprm_par_id(rprm_par_mark,rpid).gt.0.and.
825  $ rprm_par_id(rprm_par_type,rpid).eq.ptype) then
826  if(rprm_sec_act(rprm_par_id(rprm_par_mark,rpid))) then
827  ! broadcast pval; to keep consistency
828  if (ptype.eq.rpar_int) then
829  ivall = ipval
830  call bcast(ivall,isize)
831  rprm_parv_int(rpid) = ivall
832  elseif (ptype.eq.rpar_real) then
833  rvall = rpval
834  call bcast(rvall,wdsize)
835  rprm_parv_real(rpid) = rvall
836  elseif (ptype.eq.rpar_log) then
837  lvall = lpval
838  call bcast(lvall,lsize)
839  rprm_parv_log(rpid) = lvall
840  elseif (ptype.eq.rpar_str) then
841  ! check value length
842  slena = len_trim(adjustl(cpval))
843  ! remove trailing blanks
844  slen = len_trim(cpval) - slena + 1
845  if (slena.gt.rprm_lstl_mnm) then
846  call mntr_log(rprm_id,lp_deb,
847  $ 'too long parameter value; shortenning')
848  slena = min(slena,rprm_lstl_mnm)
849  endif
850  call blank(cvall,rprm_lstl_mnm)
851  cvall= cpval(slen:slen+slena- 1)
852  ! broadcast pval; to keep consistency
853  call bcast(cvall,rprm_lstl_mnm*csize)
854  rprm_parv_str(rpid) = cvall
855  else
856  write(str,'(I3)') rpid
857  call mntr_abort(rprm_id,
858  $ "Parameter set "//trim(str)//" wrong type")
859  endif
860  else
861  write(str,'(I3)') rpid
862  call mntr_warn(rprm_id,
863  $ "Parameter set "//trim(str)//" section not active")
864  endif
865  else
866  write(str,'(I3)') rpid
867  call mntr_abort(rprm_id,
868  $ "Parameter "//trim(str)//" setting error")
869  endif
870 
871  return
872  end subroutine
873 !=======================================================================
882  subroutine rprm_rp_get(ipval,rpval,lpval,cpval,rpid,ptype)
883  implicit none
884 
885  include 'SIZE'
886  include 'FRAMELP'
887  include 'RPRMD'
888 
889  ! argument list
890  integer rpid, ptype
891  integer ipval
892  real rpval
893  logical lpval
894  character*20 cpval
895 
896  ! local variables
897  character*5 str
898 !-----------------------------------------------------------------------
899  if (rprm_par_id(rprm_par_mark,rpid).gt.0.and.
900  $ rprm_par_id(rprm_par_type,rpid).eq.ptype) then
901  if(rprm_sec_act(rprm_par_id(rprm_par_mark,rpid))) then
902 
903  if (ptype.eq.rpar_int) then
904  ipval = rprm_parv_int(rpid)
905  elseif (ptype.eq.rpar_real) then
906  rpval = rprm_parv_real(rpid)
907  elseif (ptype.eq.rpar_log) then
908  lpval = rprm_parv_log(rpid)
909  elseif (ptype.eq.rpar_str) then
910  cpval = rprm_parv_str(rpid)
911  else
912  write(str,'(I3)') rpid
913  call mntr_abort(rprm_id,
914  $ "Parameter get "//trim(str)//" wrong type")
915  endif
916  else
917  write(str,'(I3)') rpid
918  call mntr_warn(rprm_id,
919  $ "Parameter get "//trim(str)//" section not active")
920  endif
921  else
922  write(str,'(I3)') rpid
923  call mntr_abort(rprm_id,
924  $ "Parameter "//trim(str)//" getting error")
925  endif
926 
927  return
928  end subroutine
929 !=======================================================================
932  subroutine rprm_dict_get()
933  implicit none
934 
935  include 'SIZE'
936  include 'PARALLEL'
937  include 'FRAMELP'
938  include 'RPRMD'
939 
940  ! local variables
941  integer il, jl, kl
942  integer nkey, ifnd, i_out
943  integer nmod, pmid
944  character*132 key, lkey
945  character*1024 val
946  logical ifoundm, ifoundp, ifact
947  integer itmp
948  real rtmp
949 
950  character*20 lname
951  character*132 ldscr
952  character*200 llog
953  integer slen,slena
954 
955 !-----------------------------------------------------------------------
956  ! dictionary exists on master node only
957  if (nid.eq.rprm_pid0) then
958  ! key number in dictionary
959  call finiparser_getdictentries(nkey)
960 
961  do il=1,nkey!rprm_par_mpos
962 
963  ! get a key
964  call finiparser_getpair(key,val,il,ifnd)
965  key = adjustl(key)
966  call capit(key,132)
967 
968  ! find section key belongs to
969  ifoundm=.false.
970  do jl=1,rprm_sec_mpos
971  if (rprm_sec_id(jl).gt.0) then
972  lname=trim(adjustl(rprm_sec_name(jl)))
973  ifnd = index(key,trim(lname))
974  if (ifnd.eq.1) then
975  ! set section to active
976  rprm_sec_act(jl) = .true.
977  ifoundm=.true.
978  ! looking for more than section name
979  if (trim(key).ne.trim(lname)) then
980  ! add variable name
981  ifoundp =.false.
982  do kl=1,rprm_par_mpos
983  if (rprm_par_id(rprm_par_mark,kl).eq.jl) then
984  lkey = trim(lname)//':'//trim(rprm_par_name(kl))
985  if (trim(key).eq.trim(lkey)) then
986  ifoundp=.true.
987  ! read parameter value
988  if (rprm_par_id(rprm_par_type,kl).eq.
989  $ rpar_int) then
990  read(val,*) itmp
991  rprm_parv_int(kl) = itmp
992  elseif (rprm_par_id(rprm_par_type,kl).eq.
993  $ rpar_real) then
994  read(val,*) rtmp
995  rprm_parv_real(kl) = rtmp
996  elseif (rprm_par_id(rprm_par_type,kl).eq.
997  $ rpar_log) then
998  call finiparser_getbool(i_out,trim(lkey),ifnd)
999  if (ifnd.eq.1) then
1000  if (i_out.eq.1) then
1001  rprm_parv_log(kl) = .true.
1002  else
1003  rprm_parv_log(kl) = .false.
1004  endif
1005  else
1006  call mntr_warn(rprm_id,
1007  $ 'Boolean parameter reading error '//trim(key))
1008  endif
1009  elseif (rprm_par_id(rprm_par_type,kl).eq.
1010  $ rpar_str) then
1011  rprm_parv_str(kl) = trim(adjustl(val))
1012  else
1013  call mntr_warn(rprm_id,
1014  $ 'Runtime parameter type missmatch '//trim(key))
1015  endif
1016  exit
1017  endif
1018  endif
1019  enddo
1020  ! is it unknown parameter
1021  if (.not.ifoundp) then
1022  call mntr_warn(rprm_id,
1023  $ 'Unknown runtime parameter '//trim(key))
1024  endif
1025  endif
1026  exit
1027  endif
1028  endif
1029  enddo
1030  if (.not.ifoundm) then
1031  ! possible palce for warning that section not found
1032  endif
1033  enddo
1034  endif
1035 
1036  ! broadcast array data
1037  call bcast(rprm_parv_int,rprm_par_id_max*isize)
1038  call bcast(rprm_parv_real,rprm_par_id_max*wdsize)
1039  call bcast(rprm_parv_log,rprm_par_id_max*lsize)
1040  call bcast(rprm_parv_str,rprm_par_id_max*rprm_lstl_mnm*csize)
1041 
1042  ! broadcast activation lfag
1043  call bcast(rprm_sec_act,rprm_sec_id_max*lsize)
1044 
1045  return
1046  end subroutine
1047 !=======================================================================
1051  subroutine rprm_rp_summary_print(unit)
1052  implicit none
1053 
1054  include 'SIZE'
1055  include 'FRAMELP'
1056  include 'RPRMD'
1057 
1058  ! argument list
1059  integer unit
1060 
1061  ! local variables
1062  integer ind(rprm_par_id_max)
1063  integer offset(2,rprm_par_id_max)
1064  integer slist(2,rprm_par_id_max), itmp1(2)
1065  integer npos, nset, key
1066  integer il, jl
1067  integer istart, in, itest
1068  character*20 str
1069  character*22 sname
1070  character*(*) cmnt
1071  parameter(cmnt='#')
1072 
1073  ! functions
1074  integer mntr_lp_def_get
1075 !-----------------------------------------------------------------------
1076  if (unit.eq.6) then
1077  call mntr_log(rprm_id,lp_prd,
1078  $ 'Summary of registered runtime parameters for active sections')
1079  else
1080  call mntr_log(rprm_id,lp_prd,
1081  $ 'Generated .par file for active sections')
1082  endif
1083 
1084  if (nid.eq.rprm_pid0) then
1085 
1086  ! sort module index array
1087  ! copy data removing possible empty slots
1088  npos=0
1089  do il=1,rprm_par_mpos
1090  in = rprm_par_id(rprm_par_mark,il)
1091  if (in.ge.0.and.rprm_sec_act(in)) then
1092  npos = npos + 1
1093  slist(1,npos) = in
1094  slist(2,npos) = il
1095  endif
1096  enddo
1097 
1098  ! sort with respect to section id
1099  key = 1
1100  call ituple_sort(slist,2,npos,key,1,ind,itmp1)
1101 
1102  ! sort parameters in single section with respect to parameter id
1103  nset = 0
1104  istart = 1
1105  itest = slist(1,istart)
1106  do il=1,npos
1107  if(itest.ne.slist(1,il).or.il.eq.npos) then
1108  if (il.eq.npos.and.itest.eq.slist(1,il)) then
1109  jl = npos + 1
1110  else
1111  jl = il
1112  endif
1113  in = jl - istart
1114  if (in.gt.1) then
1115  key = 2
1116  call ituple_sort(slist(1,istart),2,in,key,1,ind,itmp1)
1117  endif
1118  nset = nset +1
1119  offset(1,nset) = istart
1120  offset(2,nset) = in
1121  if (il.ne.npos) then
1122  itest = slist(1,il)
1123  istart = il
1124  elseif(itest.ne.slist(1,il)) then
1125  nset = nset +1
1126  offset(1,nset) = il
1127  offset(2,nset) = 1
1128  endif
1129  endif
1130  enddo
1131 
1132  if (mntr_lp_def_get().le.lp_prd.or.unit.ne.6) then
1133  if (unit.ne.6) then
1134  write(unit,'(A)') cmnt
1135  write(unit,'(A,A)') cmnt,
1136  $ ' runtime parameter file generated by'
1137  write(unit,'(A,A)') cmnt,' rprm_rp_summary_print'
1138  write(unit,'(A)') cmnt
1139  endif
1140  do il=1,nset
1141  istart = offset(1,il)
1142  in = offset(2,il)
1143  key = rprm_sec_id(slist(1,istart))
1144  sname = '['//trim(rprm_sec_name(slist(1,istart)))//']'
1145  write(unit,'(A)') cmnt
1146  write(unit,'(A,A)') sname,
1147  $ ' '//cmnt//' '//trim(adjustl(rprm_sec_dscr(slist(1,istart))))
1148  do jl = 0, in-1
1149  key = slist(2,istart+jl)
1150  if (rprm_par_id(rprm_par_type,key).eq.
1151  $ rpar_int) then
1152  write(str,'(I8)') rprm_parv_int(key)
1153  elseif (rprm_par_id(rprm_par_type,key).eq.
1154  $ rpar_real) then
1155  write(str,'(E15.8)') rprm_parv_real(key)
1156  elseif (rprm_par_id(rprm_par_type,key).eq.
1157  $ rpar_log) then
1158  if (rprm_parv_log(key)) then
1159  str = 'yes'
1160  else
1161  str = 'no'
1162  endif
1163  elseif (rprm_par_id(rprm_par_type,key).eq.
1164  $ rpar_str) then
1165  str = rprm_parv_str(key)
1166  endif
1167  write(unit,'(A," = ",A,A)')
1168  $ rprm_par_name(key), adjustl(str),
1169  $ ' '//cmnt//' '//trim(adjustl(rprm_par_dscr(key)))
1170  enddo
1171  enddo
1172  if (unit.ne.6) then
1173  write(unit,'(A)') cmnt
1174  write(unit,'(A,A)') cmnt,' end of runtime parameter file'
1175  write(unit,'(A)') cmnt
1176  else
1177  write(unit,'(A1)') ' '
1178  endif
1179  endif
1180  endif
1181 
1182 
1183  return
1184  end subroutine
1185 !=======================================================================
1197  subroutine rprm_check(mod_nkeys, mod_dictkey, mod_n3dkeys,
1198  $ mod_l3dkey, ifsec)
1199  implicit none
1200 
1201  include 'SIZE'
1202  include 'INPUT' ! IF3D
1203  include 'FRAMELP'
1204  include 'RPRMD'
1205 
1206  ! argument list
1207  integer mod_nkeys, mod_n3dkeys, mod_l3dkey(mod_n3dkeys)
1208  character*132 mod_dictkey(mod_nkeys)
1209  logical ifsec
1210 
1211  ! local variables
1212  integer il, jl, ip ! loop index
1213  ! dictionary operations
1214  integer nkey, ifnd, i_out
1215  real d_out
1216  character*132 key, lkey
1217  character*1024 val
1218  logical ifvar, if3dkey
1219 !-----------------------------------------------------------------------
1220  ! check consistency
1221  ! key number in dictionary
1222  call finiparser_getdictentries(nkey)
1223 
1224  ! set marker for finding module's section
1225  ifsec = .false.
1226  do il=1,nkey
1227  ! get a key
1228  call finiparser_getpair(key,val,il,ifnd)
1229  call capit(key,132)
1230 
1231  ! does it belong to current module's section
1232  ifnd = index(key,trim(mod_dictkey(1)))
1233  if (ifnd.eq.1) then
1234  ! section was found, check variable
1235  ifsec = .true.
1236  ifvar = .false.
1237  do ip = mod_nkeys,1,-1
1238  lkey = trim(adjustl(mod_dictkey(1)))
1239  if (ip.gt.1) lkey =trim(adjustl(lkey))//
1240  $ ':'//trim(adjustl(mod_dictkey(ip)))
1241  if(index(key,trim(lkey)).eq.1) then
1242  ifvar = .true.
1243  exit
1244  endif
1245  enddo
1246 
1247  if (ifvar) then
1248  ! check 2D versus 3D
1249  if (.not.if3d) then
1250  if3dkey = .false.
1251  do jl=1,mod_n3dkeys
1252  if (ip.eq.mod_l3dkey(jl)) then
1253  if3dkey = .true.
1254  exit
1255  endif
1256  enddo
1257 
1258  if (if3dkey) then
1259  call mntr_log(rprm_id,lp_inf,
1260  $ 'Module '//trim(mod_dictkey(1)))
1261  call mntr_log(rprm_id,lp_inf,
1262  $ '3D parameter '//trim(key)//' specified for 2D run')
1263  endif
1264  endif
1265  else
1266  ! variable not found
1267  call mntr_log(rprm_id,lp_inf,
1268  $ 'Module '//trim(mod_dictkey(1)))
1269  call mntr_log(rprm_id,lp_inf,
1270  $ 'Unknown runtime parameter: '//trim(key))
1271  endif
1272  endif
1273  enddo
1274 
1275  ! no parameter section; give warning
1276  if (.not.ifsec) then
1277  call mntr_log(rprm_id,lp_inf,'Module '//trim(mod_dictkey(1)))
1278  call mntr_log(rprm_id,lp_inf,
1279  $ 'runtime parameter section not found.')
1280  endif
1281 
1282  return
1283  end subroutine
1284 !=======================================================================
subroutine bcast(buf, len)
Definition: comm_mpi.f:431
subroutine io_file_freeid(iunit, ierr)
Get free file unit number and store max unit value.
Definition: io_tools.f:47
subroutine mntr_logi(mid, priority, logs, ivar)
Write log message adding single integer.
Definition: mntrlog.f:709
subroutine mntr_logr(mid, priority, logs, rvar)
Write log message adding single real.
Definition: mntrlog.f:731
subroutine mntr_warn(mid, logs)
Write warning message.
Definition: mntrlog.f:803
subroutine mntr_mod_is_name_reg(mid, mname)
Check if module name is registered and return its id.
Definition: mntrlog.f:459
subroutine mntr_mod_get_info(mname, pmid, mid)
Get module name an parent id for given module id. This operation is performed locally.
Definition: mntrlog.f:568
subroutine mntr_abort(mid, logs)
Abort simulation.
Definition: mntrlog.f:837
subroutine mntr_log(mid, priority, logs)
Write log message.
Definition: mntrlog.f:600
subroutine mntr_mod_reg(mid, pmid, mname, mdscr)
Register new module.
Definition: mntrlog.f:346
subroutine mntr_logl(mid, priority, logs, lvar)
Write log message adding single logical.
Definition: mntrlog.f:782
subroutine rprm_register
Register runtime parameters database.
Definition: rprm.f:10
subroutine rprm_sec_get_info(pname, mid, ifact, rpid)
Get section info based on its id. This operation is performed locally.
Definition: rprm.f:383
subroutine rprm_rp_is_name_reg(rpid, mid, pname, ptype)
Check if parameter name is registered and return its id. Check flags as well.
Definition: rprm.f:658
subroutine rprm_rp_set(rpid, ptype, ipval, rpval, lpval, cpval)
Set runtime parameter of active section. Master value is broadcasted.
Definition: rprm.f:801
subroutine rprm_rp_summary_print(unit)
Print out summary of registered runtime parameters (active sections only)
Definition: rprm.f:1052
logical function rprm_rp_is_id_reg(rpid, ptype)
Check if parameter id is registered and check type consistency. This operation is performed locally.
Definition: rprm.f:743
subroutine rprm_init
Initialise modules runtime parameters and write summary.
Definition: rprm.f:87
subroutine rprm_dict_get()
Get runtime parameter from nek parser dictionary.
Definition: rprm.f:933
subroutine rprm_rp_get(ipval, rpval, lpval, cpval, rpid, ptype)
Get runtime parameter form active section. This operation is performed locally.
Definition: rprm.f:883
logical function rprm_sec_is_id_act(rpid)
Check if section id is registered and activated. This operation is performed locally.
Definition: rprm.f:456
subroutine rprm_check(mod_nkeys, mod_dictkey, mod_n3dkeys, mod_l3dkey, ifsec)
Check consistency of module's runtime parameters.
Definition: rprm.f:1199
subroutine rprm_sec_is_name_reg(rpid, mid, pname)
Check if section name is registered and return its id. Check mid as well.
Definition: rprm.f:284
subroutine rprm_rp_get_info(pname, mid, ptype, rpid)
Get parameter info based on its id. This operation is performed locally.
Definition: rprm.f:765
logical function rprm_sec_is_id_reg(rpid)
Check if section id is registered. This operation is performed locally.
Definition: rprm.f:362
subroutine rprm_rp_reg(rpid, mid, pname, pdscr, ptype, ipval, rpval, lpval, cpval)
Register new runtime parameter.
Definition: rprm.f:484
subroutine rprm_sec_set_act(ifact, rpid)
Set section's activation flag. Master value is broadcasted.
Definition: rprm.f:422
logical function rprm_is_initialised()
Check if module was initialised.
Definition: rprm.f:147
subroutine rprm_sec_reg(rpid, mid, pname, pdscr)
Register new parameter section.
Definition: rprm.f:165
subroutine capit(lettrs, n)
Definition: ic.f:1648
subroutine blank(A, N)
Definition: math.f:19
subroutine ituple_sort(a, lda, n, key, nkey, ind, aa)
Definition: navier8.f:327