KTH framework for Nek5000 toolboxes; testing version  0.0.1
mntrlog.f
Go to the documentation of this file.
1 
6 !=======================================================================
10  subroutine mntr_register_mod(log_thr)
11  implicit none
12 
13  include 'SIZE'
14  include 'FRAMELP'
15  include 'MNTRLOGD'
16  include 'MNTRTMRD'
17 
18  ! argument list
19  integer log_thr
20 
21  ! local variables
22  character*2 str
23  character*200 lstring
24  integer il
25 
26  ! functions
27  integer frame_get_master
28  real dnekclock
29 !-----------------------------------------------------------------------
30  ! simple timing
31  mntr_frame_tmini = dnekclock()
32 
33  ! amd seems to have big problem with block data for arrays
34  mntr_ifinit = .false.
35  mntr_stdl = 1
36  mntr_ifconv = .false.
37  mntr_pid0 = 0
38  mntr_mod_num = 0
39  mntr_mod_mpos = 0
40  do il = 1, mntr_id_max
41  mntr_mod_id(il) = -1
42  mntr_mod_name(il) = mntr_blname
43  end do
44  mntr_tmr_num = 0
45  mntr_tmr_mpos = 0
46  do il = 1, mntr_tmr_id_size
47  mntr_tmr_id(il,1) = -1
48  mntr_tmrv_timer(il,1) = 0.0
49  end do
50  do il = 1, mntr_tmr_id_max
51  mntr_tmr_sum(il) = .false.
52  mntr_tmr_name(il) = mntr_blname
53  end do
54 
55  ! set master node
56  mntr_pid0 = frame_get_master()
57 
58  ! first register framework
59  mntr_frame_id = 1
60  mntr_mod_id(mntr_frame_id) = 0
61  mntr_mod_name(mntr_frame_id) = mntr_frame_name
62  mntr_mod_dscr(mntr_frame_id) = 'Framework backbone'
63  mntr_mod_num = mntr_mod_num + 1
64  mntr_mod_mpos = mntr_mod_mpos + 1
65 
66  ! next monitor
67  mntr_id = 2
68  mntr_mod_id(mntr_id) = mntr_frame_id
69  mntr_mod_name(mntr_id) = mntr_name
70  mntr_mod_dscr(mntr_id) = 'Monitoring module'
71  mntr_mod_num = mntr_mod_num + 1
72  mntr_mod_mpos = mntr_mod_mpos + 1
73 
74  ! set log threshold
75  mntr_lp_def = log_thr
76 
77  ! log changes
78  lstring ='Registered module ['//trim(mntr_mod_name(mntr_frame_id))
79  lstring= trim(lstring)//']: '//trim(mntr_mod_dscr(mntr_frame_id))
80  call mntr_log(mntr_id,lp_inf,trim(lstring))
81 
82  lstring = 'Registered module ['//trim(mntr_mod_name(mntr_id))
83  lstring= trim(lstring)//']: '//trim(mntr_mod_dscr(mntr_id))
84  call mntr_log(mntr_id,lp_inf,trim(lstring))
85 
86  ! register framework timer and get initiaisation time
87  call mntr_tmr_reg(mntr_frame_tmr_id,0,mntr_frame_id,
88  $ 'FRM_TOT','Total elapsed framework time',.false.)
89 
90  write(str,'(I2)') mntr_lp_def
91  call mntr_log(mntr_id,lp_inf,
92  $ 'Initial log threshold set to: '//trim(str))
93 
94  return
95  end subroutine
96 !=======================================================================
99  subroutine mntr_register_par
100  implicit none
101 
102  include 'SIZE'
103  include 'FRAMELP'
104  include 'MNTRLOGD'
105 
106  ! local variables
107  integer rpid,itmp
108  real rtmp
109  logical ltmp
110  character*20 ctmp
111 !-----------------------------------------------------------------------
112  ! register and set active section
113  call rprm_sec_reg(mntr_sec_id,mntr_id,'_'//adjustl(mntr_name),
114  $ 'Runtime parameter section for monitor module')
115  call rprm_sec_set_act(.true.,mntr_sec_id)
116 
117  ! register parameters
118  call rprm_rp_reg(mntr_lp_def_id,mntr_sec_id,'LOGLEVEL',
119  $ 'Logging threshold for toolboxes',rpar_int,mntr_lp_def,
120  $ 0.0,.false.,' ')
121 
122  call rprm_rp_reg(mntr_iftdsc_id,mntr_sec_id,'IFTIMDSCR',
123  $ 'Write timer description in the summary',rpar_log,0,
124  $ 0.0,.false.,' ')
125 
126  call rprm_rp_reg(mntr_wtime_id,mntr_sec_id,'WALLTIME',
127  $ 'Simulation wall time',rpar_str,0,0.0,.false.,'00:00')
128 
129  return
130  end subroutine
131 !=======================================================================
134  subroutine mntr_init
135  implicit none
136 
137  include 'SIZE'
138  include 'FRAMELP'
139  include 'MNTRLOGD'
140 
141  ! local variables
142  integer ierr, nhour, nmin
143  integer itmp
144  real rtmp
145  logical ltmp
146  character*20 ctmp
147  character*2 str
148 !-----------------------------------------------------------------------
149  call rprm_rp_get(itmp,rtmp,ltmp,ctmp,mntr_lp_def_id,rpar_int)
150  mntr_lp_def = itmp
151 
152  call rprm_rp_get(itmp,rtmp,ltmp,ctmp,mntr_iftdsc_id,rpar_log)
153  mntr_iftdsc = ltmp
154 
155  write(str,'(I2)') mntr_lp_def
156  call mntr_log(mntr_id,lp_inf,
157  $ 'Reseting log threshold to: '//trim(str))
158 
159  call rprm_rp_get(itmp,rtmp,ltmp,ctmp,mntr_wtime_id,rpar_str)
160  mntr_wtimes = ctmp
161 
162  ! get wall clock
163  ctmp = trim(adjustl(mntr_wtimes))
164  ! check string format
165  ierr = 0
166  if (ctmp(3:3).ne.':') ierr = 1
167  if (.not.(lge(ctmp(1:1),'0').and.lle(ctmp(1:1),'9'))) ierr = 1
168  if (.not.(lge(ctmp(2:2),'0').and.lle(ctmp(2:2),'9'))) ierr = 1
169  if (.not.(lge(ctmp(4:4),'0').and.lle(ctmp(4:4),'9'))) ierr = 1
170  if (.not.(lge(ctmp(5:5),'0').and.lle(ctmp(5:5),'9'))) ierr = 1
171 
172  if (ierr.eq.0) then
173  read(ctmp(1:2),'(I2)') nhour
174  read(ctmp(4:5),'(I2)') nmin
175  mntr_wtime = 60.0*(nmin +60*nhour)
176  else
177  call mntr_log(mntr_id,lp_inf,'Wrong wall time format')
178  endif
179 
180  ! write summary
182 
183  mntr_ifinit = .true.
184 
185  return
186  end subroutine
187 !=======================================================================
190  subroutine mntr_wclock
191  implicit none
192 
193  include 'SIZE'
194  include 'TSTEP' ! ISTEP, NSTEPS, LASTEP
195  include 'INPUT' ! IFMVBD, IFREGUO
196  include 'PARALLEL' ! WDSIZE
197  include 'CTIMER' ! ETIMES
198  include 'FRAMELP'
199  include 'MNTRLOGD'
200 
201  ! local variables
202  integer il, lstdl
203  real rtmp
204 !-----------------------------------------------------------------------
205  ! double delay step as monitoing routine does not know when checkpointing
206  ! starts and wall clock can be reached in the middle of writnig process
207  lstdl = 2*mntr_stdl+1
208 
209  ! check simulation wall time
210  if (mntr_wtime.gt.0.0) then
211 
212  ! save wall time of the current step
213  do il=lstdl,2,-1
214  mntr_wtstep(il) = mntr_wtstep(il-1)
215  enddo
216  mntr_wtstep(1) = dnekclock() - etimes
217  ! check if simulation is going to exceed wall time, but
218  ! first let read all checkpointing files (necessary for multi file
219  ! checkpointing)
220  if (istep.gt.lstdl) then
221  ! it should be enough for the master to check condition
222  if (nid.eq.mntr_pid0) rtmp = 2.0*mntr_wtstep(1) -
223  $ mntr_wtstep(lstdl)
224  ! broadcast predicted time
225  il = wdsize
226  call bcast(rtmp,il)
227 
228  if (rtmp.gt.mntr_wtime.and.(nsteps-istep).gt.lstdl) then
229  call mntr_log(mntr_id,lp_inf,
230  $ 'Wall clock reached; adjust NSTEPS')
231  nsteps = istep+lstdl
232  endif
233  endif
234  endif
235 
236  ! check convergence flag
237  if (mntr_ifconv.and.(nsteps-istep).gt.lstdl) then
238  call mntr_log(mntr_id,lp_inf,
239  $ 'Simulation converged; adjust NSTEPS')
240  nsteps = istep+lstdl
241  endif
242 
243  ! just to take into account there is istep and kstep,
244  ! and kstep is just a local variable
245  if (istep.ge.nsteps) lastep=1
246 
247  return
248  end subroutine
249 !=======================================================================
253  subroutine mntr_set_step_delay(dstep)
254  implicit none
255 
256  include 'SIZE'
257  include 'FRAMELP'
258  include 'MNTRLOGD'
259 
260  ! argument list
261  integer dstep
262 
263 !-----------------------------------------------------------------------
264  if (dstep.gt.mntr_stdl_max) then
265  call mntr_abort(mntr_id,"Step delay exceeds mntr_stdl_max")
266  else
267  mntr_stdl = max(mntr_stdl,dstep)
268  endif
269 
270  return
271  end subroutine
272 !=======================================================================
276  subroutine mntr_get_step_delay(dstep)
277  implicit none
278 
279  include 'SIZE'
280  include 'FRAMELP'
281  include 'MNTRLOGD'
282 
283  ! argument list
284  integer dstep
285 !-----------------------------------------------------------------------
286  dstep = mntr_stdl
287 
288  return
289  end subroutine
290 !=======================================================================
293  subroutine mntr_set_conv(ifconv)
294  implicit none
295 
296  include 'SIZE'
297  include 'FRAMELP'
298  include 'MNTRLOGD'
299 
300  ! argument list
301  logical ifconv
302 
303 !-----------------------------------------------------------------------
304  mntr_ifconv = ifconv
305 
306  return
307  end subroutine
308 !=======================================================================
312  logical function mntr_is_initialised()
313  implicit none
314 
315  include 'SIZE'
316  include 'FRAMELP'
317  include 'MNTRLOGD'
318 !-----------------------------------------------------------------------
319  mntr_is_initialised = mntr_ifinit
320 
321  return
322  end function
323 !=======================================================================
327  integer function mntr_lp_def_get()
328  implicit none
329 
330  include 'SIZE'
331  include 'FRAMELP'
332  include 'MNTRLOGD'
333 !-----------------------------------------------------------------------
334  mntr_lp_def_get = mntr_lp_def
335 
336  return
337  end function
338 !=======================================================================
345  subroutine mntr_mod_reg(mid,pmid,mname,mdscr)
346  implicit none
347 
348  include 'SIZE'
349  include 'PARALLEL' ! ISIZE
350  include 'FRAMELP'
351  include 'MNTRLOGD'
352 
353  ! argument list
354  integer mid, pmid
355  character*(*) mname, mdscr
356 
357  ! local variables
358  character*10 lname
359  character*132 ldscr
360  integer slen,slena
361 
362  integer il, ipos
363 !-----------------------------------------------------------------------
364  ! check name length
365  slena = len_trim(adjustl(mname))
366  ! remove trailing blanks
367  slen = len_trim(mname) - slena + 1
368  if (slena.gt.mntr_lstl_mnm) then
369  call mntr_log(mntr_id,lp_deb,
370  $ 'too long module name; shortenning')
371  slena = min(slena,mntr_lstl_mnm)
372  endif
373  call blank(lname,mntr_lstl_mnm)
374  lname= mname(slen:slen+slena- 1)
375  call capit(lname,slena)
376 
377  ! check description length
378  slena = len_trim(adjustl(mdscr))
379  ! remove trailing blanks
380  slen = len_trim(mdscr) - slena + 1
381  if (slena.ge.mntr_lstl_mds) then
382  call mntr_log(mntr_id,lp_deb,
383  $ 'too long module description; shortenning')
384  slena = min(slena,mntr_lstl_mnm)
385  endif
386  call blank(ldscr,mntr_lstl_mds)
387  ldscr= mdscr(slen:slen + slena - 1)
388 
389  ! find empty space
390  ipos = 0
391 
392  ! to ensure consistency I do it on master and broadcast result
393  if (nid.eq.mntr_pid0) then
394 
395  ! check if module is already registered
396  do il=1,mntr_mod_mpos
397  if (mntr_mod_id(il).ge.0.and.
398  $ mntr_mod_name(il).eq.lname) then
399  ipos = -il
400  exit
401  endif
402  enddo
403 
404  ! find empty spot
405  if (ipos.eq.0) then
406  do il=1,mntr_id_max
407  if (mntr_mod_id(il).eq.-1) then
408  ipos = il
409  exit
410  endif
411  enddo
412  endif
413  endif
414 
415  ! broadcast mid
416  call bcast(ipos,isize)
417 
418  ! error; no free space found
419  if (ipos.eq.0) then
420  mid = ipos
421  call mntr_abort(mntr_id,
422  $ 'module ['//trim(lname)//'] cannot be registered')
423  ! module already registered
424  elseif (ipos.lt.0) then
425  mid = abs(ipos)
426  call mntr_abort(mntr_id,
427  $ 'Module ['//trim(lname)//'] is already registered')
428  ! new module
429  else
430  mid = ipos
431  ! check if parent module is registered
432  if (pmid.gt.0) then
433  if (mntr_mod_id(pmid).ge.0) then
434  mntr_mod_id(ipos) = pmid
435  else
436  mntr_mod_id(ipos) = 0
437  call mntr_log(mntr_id,lp_inf,
438  $ "Module's ["//trim(lname)//"] parent not registered.")
439  endif
440  else
441  mntr_mod_id(ipos) = 0
442  endif
443  mntr_mod_name(ipos)=lname
444  mntr_mod_dscr(ipos)=ldscr
445  mntr_mod_num = mntr_mod_num + 1
446  if (mntr_mod_mpos.lt.ipos) mntr_mod_mpos = ipos
447  call mntr_log(mntr_id,lp_inf,
448  $ 'Registered module ['//trim(lname)//']: '//trim(ldscr))
449  endif
450 
451  return
452  end subroutine
453 !=======================================================================
458  subroutine mntr_mod_is_name_reg(mid,mname)
459  implicit none
460 
461  include 'SIZE'
462  include 'PARALLEL' ! ISIZE
463  include 'FRAMELP'
464  include 'MNTRLOGD'
465 
466  ! argument list
467  integer mid
468  character*(*) mname
469 
470  ! local variables
471  character*10 lname
472  character*3 str
473  integer slen,slena
474 
475  integer il, ipos
476 !-----------------------------------------------------------------------
477  ! check name length
478  slena = len_trim(adjustl(mname))
479  ! remove trailing blanks
480  slen = len_trim(mname) - slena + 1
481  if (slena.gt.mntr_lstl_mnm) then
482  call mntr_log(mntr_id,lp_deb,
483  $ 'too long module name; shortenning')
484  slena = min(slena,mntr_lstl_mnm)
485  endif
486  call blank(lname,mntr_lstl_mnm)
487  lname= mname(slen:slen+slena- 1)
488  call capit(lname,slena)
489 
490  ! find module
491  ipos = 0
492 
493  ! to ensure consistency I do it on master and broadcast result
494  if (nid.eq.mntr_pid0) then
495  ! check if module is already registered
496  do il=1,mntr_mod_mpos
497  if (mntr_mod_id(il).ge.0.and.
498  $ mntr_mod_name(il).eq.lname) then
499  ipos = il
500  exit
501  endif
502  enddo
503  endif
504 
505  ! broadcast ipos
506  call bcast(ipos,isize)
507 
508  if (ipos.eq.0) then
509  mid = -1
510  call mntr_log(mntr_id,lp_inf,
511  $ 'Module ['//trim(lname)//'] not registered')
512  else
513  mid = ipos
514  write(str,'(I3)') ipos
515  call mntr_log(mntr_id,lp_vrb,
516  $ 'Module ['//trim(lname)//'] registered with mid='//str)
517  endif
518 
519  return
520  end subroutine
521 !=======================================================================
526  logical function mntr_mod_is_id_reg(mid)
527  implicit none
528 
529  include 'SIZE'
530  include 'PARALLEL'
531  include 'FRAMELP'
532  include 'MNTRLOGD'
533 
534  ! argument list
535  integer mid
536 !-----------------------------------------------------------------------
537  mntr_mod_is_id_reg = mntr_mod_id(mid).ge.0
538 
539  return
540  end function
541 !=======================================================================
546  subroutine mntr_mod_get_number(nmod,mmod)
547  implicit none
548 
549  include 'SIZE'
550  include 'FRAMELP'
551  include 'MNTRLOGD'
552 
553  ! argument list
554  integer nmod, mmod
555 !-----------------------------------------------------------------------
556  nmod = mntr_mod_num
557  mmod = mntr_mod_mpos
558 
559  return
560  end subroutine
561 !=======================================================================
567  subroutine mntr_mod_get_info(mname, pmid,mid)
568  implicit none
569 
570  include 'SIZE'
571  include 'FRAMELP'
572  include 'MNTRLOGD'
573 
574  ! argument list
575  character*10 mname
576  integer mid, pmid
577 
578  ! local variables
579  character*5 str
580 !-----------------------------------------------------------------------
581  if (mntr_mod_id(mid).ge.0) then
582  pmid = mntr_mod_id(mid)
583  mname = mntr_mod_name(mid)
584  else
585  mid = -1
586  write(str,'(I3)') mid
587  call mntr_log(mntr_id,lp_vrb,
588  $ 'Module id'//trim(str)//' not registered')
589  endif
590 
591  return
592  end subroutine
593 !=======================================================================
599  subroutine mntr_log(mid,priority,logs)
600  implicit none
601 
602  include 'SIZE'
603  include 'FRAMELP'
604  include 'MNTRLOGD'
605 
606  ! argument list
607  integer mid,priority
608  character*(*) logs
609 
610  ! local variables
611  character*200 llogs
612  character*5 str
613  integer slen, slena
614 !-----------------------------------------------------------------------
615  ! check log priority
616  if (priority.lt.mntr_lp_def) return
617 
618  ! done only by master
619  if (nid.eq.mntr_pid0) then
620 
621  ! check description length
622  slena = len_trim(adjustl(logs))
623  ! remove trailing blanks
624  slen = len_trim(logs) - slena + 1
625  if (slena.ge.mntr_lstl_log) then
626  if (mntr_lp_def.le.lp_deb) write(*,*)' ['//mntr_name//'] ',
627  $ 'too long log string; shortenning'
628  slena = min(slena,mntr_lstl_log)
629  endif
630  call blank(llogs,mntr_lstl_mds)
631  llogs= logs(slen:slen + slena - 1)
632 
633  ! check module id
634  if (mntr_mod_id(mid).ge.0) then
635  ! add module name
636  write(*,*) ' ['//trim(mntr_mod_name(mid))//'] '//trim(llogs)
637  else
638  write(str,'(I3)') mid
639  write(*,*) ' ['//trim(mntr_name)//'] ',
640  $ ' WARNING: module'//trim(str)//' not registered;'
641  write(*,*) 'Log body: '//trim(llogs)
642  endif
643  endif
644 
645  return
646  end subroutine
647 !=======================================================================
654  subroutine mntr_log_local(mid,priority,logs,prid)
655  implicit none
656 
657  include 'SIZE'
658  include 'FRAMELP'
659  include 'MNTRLOGD'
660 
661  ! argument list
662  integer mid,priority, prid
663  character*(*) logs
664 
665  ! local variables
666  character*200 llogs
667  character*5 str
668  integer slen, slena
669 !-----------------------------------------------------------------------
670  ! check log priority
671  if (priority.lt.mntr_lp_def) return
672 
673  ! done only by given process
674 
675  ! check description length
676  slena = len_trim(adjustl(logs))
677  ! remove trailing blanks
678  slen = len_trim(logs) - slena + 1
679  if (slena.ge.mntr_lstl_log) then
680  if (mntr_lp_def.le.lp_deb) write(*,*)' ['//mntr_name//'] ',
681  $ 'too long log string; shortenning'
682  slena = min(slena,mntr_lstl_log)
683  endif
684  call blank(llogs,mntr_lstl_mds)
685  llogs= logs(slen:slen + slena - 1)
686 
687  ! check module id
688  if (mntr_mod_id(mid).ge.0) then
689  ! add module name
690  write(*,*) ' ['//trim(mntr_mod_name(mid))//'] nid= ',prid,
691  $ ' '//trim(llogs)
692  else
693  write(str,'(I3)') mid
694  write(*,*) ' ['//trim(mntr_name)//'] ',
695  $ ' WARNING: module'//trim(str)//' not registered;'
696  write(*,*) 'Log body: nid= ',prid,' '//trim(llogs)
697  endif
698 
699  return
700  end subroutine
701 !=======================================================================
708  subroutine mntr_logi(mid,priority,logs,ivar)
709  implicit none
710 
711  ! argument list
712  integer mid,priority,ivar
713  character*(*) logs
714 
715  ! local variables
716  character*10 str
717 !-----------------------------------------------------------------------
718  write(str,'(I8)') ivar
719  call mntr_log(mid,priority,trim(logs)//' '//trim(str))
720 
721  return
722  end subroutine
723 !=======================================================================
730  subroutine mntr_logr(mid,priority,logs,rvar)
731  implicit none
732 
733  ! argument list
734  integer mid,priority
735  character*(*) logs
736  real rvar
737 
738  ! local variables
739  character*20 str
740 !-----------------------------------------------------------------------
741  write(str,'(E15.8)') rvar
742  call mntr_log(mid,priority,trim(logs)//' '//trim(str))
743 
744  return
745  end subroutine
746 !=======================================================================
754  subroutine mntr_logrv(mid,priority,logs,rvarv,rvarn)
755  implicit none
756 
757  ! argument list
758  integer mid,priority,rvarn
759  character*(*) logs
760  real rvarv(rvarn)
761 
762  ! local variables
763  integer il
764  character*20 str
765 !-----------------------------------------------------------------------
766  call mntr_log(mid,priority,trim(logs))
767  do il = 1, rvarn
768  write(str,'(E15.8)') rvarv(il)
769  call mntr_log(mid,priority,' '//trim(str))
770  end do
771 
772  return
773  end subroutine
774 !=======================================================================
781  subroutine mntr_logl(mid,priority,logs,lvar)
782  implicit none
783 
784  ! argument list
785  integer mid,priority
786  character*(*) logs
787  logical lvar
788 
789  ! local variables
790  character*2 str
791 !-----------------------------------------------------------------------
792  write(str,'(L2)') lvar
793  call mntr_log(mid,priority,trim(logs)//' '//trim(str))
794 
795  return
796  end subroutine
797 !=======================================================================
802  subroutine mntr_warn(mid,logs)
803  implicit none
804 
805  include 'FRAMELP'
806 
807  ! argument list
808  integer mid,priority
809  character*(*) logs
810 !-----------------------------------------------------------------------
811  call mntr_log(mid,lp_inf,'WARNING: '//logs)
812  return
813  end subroutine
814 !=======================================================================
819  subroutine mntr_error(mid,logs)
820  implicit none
821 
822  include 'FRAMELP'
823 
824  ! argument list
825  integer mid
826  character*(*) logs
827 !-----------------------------------------------------------------------
828  call mntr_log(mid,lp_err,'ERROR: '//logs)
829  return
830  end subroutine
831 !=======================================================================
836  subroutine mntr_abort(mid,logs)
837  implicit none
838 
839  include 'FRAMELP'
840 
841  ! argument list
842  integer mid
843  character*(*) logs
844 !-----------------------------------------------------------------------
845  call mntr_log(mid,lp_err,'ABORT: '//logs)
846  call exitt
847  return
848  end subroutine
849 !=======================================================================
855  subroutine mntr_check_abort(mid,ierr,logs)
856  implicit none
857 
858  include 'FRAMELP'
859 
860  ! argument list
861  integer mid,ierr
862  character*(*) logs
863 
864  ! local variables
865  integer imax, imin, itest
866  character*5 str
867  ! functions
868  integer iglmax, iglmin
869 !-----------------------------------------------------------------------
870  imax = iglmax(ierr,1)
871  imin = iglmin(ierr,1)
872  if (imax.gt.0) then
873  itest = imax
874  else
875  itest = imin
876  endif
877 
878  if (itest.ne.0) then
879  write(str,'(I3)') itest
880  call mntr_log(mid,lp_err,
881  $ 'ABORT: '//trim(logs)//' ierr='//trim(str))
882  call exitt
883  endif
884  return
885  end subroutine
886 !=======================================================================
890  implicit none
891 
892  include 'SIZE'
893  include 'FRAMELP'
894  include 'MNTRLOGD'
895 
896  ! local variables
897  integer il, stride
898  parameter(stride=4)
899  integer olist(2,mntr_id_max), ierr
900  character*25 ftm
901  character*3 str
902 !-----------------------------------------------------------------------
903  call mntr_log(mntr_id,lp_prd,
904  $ 'Summary of registered modules')
905 
906  if (nid.eq.mntr_pid0) then
907  ! get ordered list
908  call mntr_mod_get_olist(olist, ierr)
909 
910 
911  if(ierr.eq.0.and.mntr_lp_def.le.lp_prd) then
912  do il=1,mntr_mod_num
913  write(str,'(I3)') stride*(olist(2,il))
914  ftm = '('//trim(str)//'X,"[",A,"] : ",A)'
915  write(*,ftm) mntr_mod_name(olist(1,il)),
916  $ mntr_mod_dscr(olist(1,il))
917  enddo
918  endif
919  endif
920 
921  return
922  end subroutine
923 !=======================================================================
928  subroutine mntr_mod_get_olist(olist,ierr)
929  implicit none
930 
931  include 'SIZE'
932  include 'FRAMELP'
933  include 'MNTRLOGD'
934 
935  ! argument list
936  integer olist(2,mntr_id_max), ierr
937 
938  ! local variables
939  integer ind(mntr_id_max), level, parent, ipos
940  integer slist(2,mntr_id_max), itmp1(2)
941  integer npos, key
942  integer il, jl
943  integer istart, in, itest
944 !-----------------------------------------------------------------------
945  ierr = 0
946 
947  ! sort module index array
948  ! copy data removing possible empty slots
949  npos=0
950  do il=1,mntr_mod_mpos
951  if (mntr_mod_id(il).ge.0) then
952  npos = npos + 1
953  slist(1,npos) = mntr_mod_id(il)
954  slist(2,npos) = il
955  endif
956  enddo
957  if(npos.ne.mntr_mod_num) then
958  ierr = 1
959  call mntr_log_local(mntr_id,lp_inf,
960  $ 'Inconsistent module number; return',mntr_pid0)
961  return
962  endif
963 
964  ! sort with respect to parent id
965  key = 1
966  call ituple_sort(slist,2,npos,key,1,ind,itmp1)
967 
968  ! sort within children of single parent with respect to children id
969  istart = 1
970  itest = slist(1,istart)
971  do il=1,npos
972  if(itest.ne.slist(1,il).or.il.eq.npos) then
973  if (il.eq.npos.and.itest.eq.slist(1,il)) then
974  jl = npos + 1
975  else
976  jl = il
977  endif
978  in = jl - istart
979  if (itest.eq.0.and.in.ne.1) then
980  call mntr_log_local(mntr_id,lp_inf,
981  $ 'Must be single root of the graph; return',mntr_pid0)
982  ierr = 2
983  return
984  endif
985  if (in.gt.1) then
986  key = 2
987  call ituple_sort(slist(1,istart),2,in,key,1,ind,itmp1)
988  endif
989  if (il.ne.npos) then
990  itest = slist(1,il)
991  istart = il
992  endif
993  endif
994  enddo
995 
996  parent = 0
997  level = 0
998  ipos = 1
999  call mntr_build_ord_list(olist,slist,npos,ipos,parent,level)
1000 
1001  return
1002  end subroutine
1003 !=======================================================================
1012  recursive subroutine mntr_build_ord_list(olist,slist,nlist,npos,
1013  $ parent,level)
1014  implicit none
1015 
1016  ! argument list
1017  integer nlist, npos, parent, level
1018  integer olist(2,nlist),slist(2,nlist)
1019 
1020  ! local variables
1021  integer il
1022  integer lparent, llevel
1023 !-----------------------------------------------------------------------
1024  llevel = level + 1
1025  do il=1, nlist
1026  if (slist(1,il).eq.parent) then
1027  slist(1,il) = - parent
1028  lparent = slist(2,il)
1029  olist(1,npos) = lparent
1030  olist(2,npos) = llevel
1031  npos = npos +1
1032  call mntr_build_ord_list(olist,slist,nlist,npos,lparent,
1033  $ llevel)
1034  endif
1035  enddo
1036 
1037  return
1038  end subroutine
1039 !=======================================================================
void exitt()
Definition: comm_mpi.f:604
subroutine bcast(buf, len)
Definition: comm_mpi.f:431
real *8 function dnekclock()
Definition: comm_mpi.f:393
subroutine mntr_logi(mid, priority, logs, ivar)
Write log message adding single integer.
Definition: mntrlog.f:709
subroutine mntr_mod_summary_print()
Print registered modules showing tree structure.
Definition: mntrlog.f:890
recursive subroutine mntr_build_ord_list(olist, slist, nlist, npos, parent, level)
Build ordered list reflecting graph structure.
Definition: mntrlog.f:1014
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_logrv(mid, priority, logs, rvarv, rvarn)
Write log message adding real vector of length n.
Definition: mntrlog.f:755
subroutine mntr_mod_get_olist(olist, ierr)
Provide ordered list of registered modules for printing.
Definition: mntrlog.f:929
integer function mntr_lp_def_get()
Get logging threashold.
Definition: mntrlog.f:328
subroutine mntr_abort(mid, logs)
Abort simulation.
Definition: mntrlog.f:837
subroutine mntr_set_conv(ifconv)
Set convergence flag to shorten simulation.
Definition: mntrlog.f:294
subroutine mntr_log(mid, priority, logs)
Write log message.
Definition: mntrlog.f:600
subroutine mntr_wclock
Monitor simulation wall clock.
Definition: mntrlog.f:191
subroutine mntr_error(mid, logs)
Write error message.
Definition: mntrlog.f:820
subroutine mntr_mod_get_number(nmod, mmod)
Get number of registered modules. This operation is performed locally.
Definition: mntrlog.f:547
subroutine mntr_register_mod(log_thr)
Initialise monitor by registering framework and monitor.
Definition: mntrlog.f:11
subroutine mntr_mod_reg(mid, pmid, mname, mdscr)
Register new module.
Definition: mntrlog.f:346
subroutine mntr_init
Initialise monitor module.
Definition: mntrlog.f:135
subroutine mntr_logl(mid, priority, logs, lvar)
Write log message adding single logical.
Definition: mntrlog.f:782
subroutine mntr_get_step_delay(dstep)
Get step delay.
Definition: mntrlog.f:277
logical function mntr_is_initialised()
Check if module was initialised.
Definition: mntrlog.f:313
subroutine mntr_log_local(mid, priority, logs, prid)
Write log message from given process.
Definition: mntrlog.f:655
subroutine mntr_tmr_reg(mid, pmid, modid, mname, mdscr, ifsum)
Register new timer.
Definition: mntrtmr.f:16
subroutine mntr_check_abort(mid, ierr, logs)
Abort simulation.
Definition: mntrlog.f:856
subroutine mntr_register_par
Register monitor runtime parameters.
Definition: mntrlog.f:100
subroutine mntr_set_step_delay(dstep)
Set number of steps necessary to write proper checkpointing.
Definition: mntrlog.f:254
logical function mntr_mod_is_id_reg(mid)
Check if module id is registered. This operation is performed locally.
Definition: mntrlog.f:527
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
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
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