My Project
Data Types | Functions/Subroutines | Variables
probes Module Reference

Data Types

interface  assignment(=)
 
interface  myprobe
 
type  probe_obj
 

Functions/Subroutines

subroutine init_nml_probe
 
subroutine alloc_probe (PROBE, N)
 
subroutine assign_probe (A, B)
 
subroutine myprobe_vec (PROBE, VEC)
 
subroutine myprobe_arr (PROBE, ARR)
 
subroutine set_probes (P_ON, NP, FNM)
 
subroutine open_probes
 
subroutine dump_probe_data
 
subroutine probe_store (PROBE)
 

Variables

character(len=80) probe_interval
 
integer probe_location
 
integer, dimension(2) probe_levels
 
character(len=80) probe_title
 
character(len=80) probe_description
 
character(len=80) probe_variable
 
character(len=80) probe_var_name
 
type(probe_obj), dimension(:), allocatable glb_probe
 
type(probe_obj), dimension(:), allocatable lcl_probe
 
logical probe_on
 
integer global_probes
 
integer local_probes
 
character input_fname
 

Function/Subroutine Documentation

◆ alloc_probe()

subroutine probes::alloc_probe ( type(probe_obj), dimension(:), allocatable  PROBE,
integer, intent(in)  N 
)

Definition at line 174 of file mod_probe.f90.

174  IMPLICIT NONE
175  INTEGER,INTENT(IN) :: N
176  TYPE(PROBE_OBJ),Allocatable :: PROBE(:)
177  INTEGER :: STATUS,I
178 
179  ALLOCATE(probe(n),stat=status)
180  IF(status/=0) CALL fatal_error("MOD_PROBE: COULD NOT ALLOCATE PROBE TYPE")
181 
182  ! INITIALIZE TYPE DATA
183  probe%MINE= .false.
184  probe%O_INT= zerotime
185  probe%O_NEXT= zerotime
186  probe%O_NUM=0
187  probe%D_LOC=0
188  probe%D_LOC_GL=0
189  probe%K_ONE=0
190  probe%K_TWO=0
191  probe%xloc=0.0_sp
192  probe%yloc=0.0_sp
193  probe%lonloc=0.0_sp
194  probe%latloc=0.0_sp
195  probe%dpth=0.0_sp
196  probe%D_TIT ='none'
197  probe%D_DES ='none'
198  probe%VAR ='none'
199  probe%VNAME ='none'
200  probe%FILENAME ='none'
201  DO i=1,n
202  NULLIFY(probe(i)%SCL)
203  NULLIFY(probe(i)%VEC)
204  END DO
205 
type(time) zerotime
Definition: mod_main.f90:830
integer n
Definition: mod_main.f90:55
Here is the call graph for this function:
Here is the caller graph for this function:

◆ assign_probe()

subroutine probes::assign_probe ( type(probe_obj), intent(out)  A,
type(probe_obj), intent(in)  B 
)

Definition at line 210 of file mod_probe.f90.

210  IMPLICIT NONE
211  TYPE(PROBE_OBJ), INTENT(OUT) ::A
212  TYPE(PROBE_OBJ), INTENT(IN) ::B
213 
214  a%MINE = b%MINE
215  a%O_INT = b%O_INT
216  a%O_NEXT = b%O_NEXT
217  a%O_NUM = b%O_NUM
218  a%D_LOC = b%D_LOC
219  a%D_LOC_GL = b%D_LOC_GL
220  a%K_ONE = b%K_ONE
221  a%K_TWO = b%K_TWO
222  a%xloc = b%xloc
223  a%yloc = b%yloc
224  a%lonloc = b%lonloc
225  a%latloc = b%latloc
226  a%dpth = b%dpth
227  a%D_TIT = b%D_TIT
228  a%D_DES = b%D_DES
229  a%VAR = b%VAR
230  a%VNAME = b%VNAME
231  a%FILENAME = b%FILENAME
232 
233  a%VEC =>b%VEC
234  a%SCL =>b%SCL
235 

◆ dump_probe_data()

subroutine probes::dump_probe_data ( )

Definition at line 761 of file mod_probe.f90.

761 
762 !------------------------------------------------------------------------------|
763 ! WRITE TIME SERIES DATA TO TIME SERIES DATA FILES |
764 !------------------------------------------------------------------------------|
765 
766  USE mod_prec
767  USE all_vars
768  IMPLICIT NONE
769  INTEGER I,K,K1,K2,IUNIT
770  TYPE(TIME) :: OTIME
771 
772  if(dbg_set(dbg_sbr)) write(ipt,*) "START: DUMP_PROBE_DATA"
773 
774 !==============================================================================!
775 ! MAIN LOOP OVER TIME SERIES OUTPUT !
776 !==============================================================================!
777  DO i=1,local_probes
778 
779  !----Return if not on Time Series Write Interval-----------------------------
780  IF( lcl_probe(i)%O_NEXT .GT. inttime) cycle
781 
782  lcl_probe(i)%O_NEXT = inttime + lcl_probe(i)%O_INT
783  !----Open File For Write-----------------------------------------------------
784  iunit = lcl_probe(i)%O_NUM + 100
785  OPEN(unit=iunit,file=lcl_probe(i)%FILENAME,form='FORMATTED',position='APPEND')
786 
787  !----Write Data to File------------------------------------------------------
788  IF(ASSOCIATED(lcl_probe(i)%VEC))THEN
789  k1=lcl_probe(i)%K_one
790  k2=lcl_probe(i)%K_two
791 ! WRITE(IUNIT,*) DAYS(IntTime),(LCL_PROBE(I)%VEC(K),K=K1,K2)
792  WRITE(iunit,'(f15.5,50f9.3)') days(inttime),(lcl_probe(i)%VEC)
793  ELSE IF(ASSOCIATED(lcl_probe(i)%SCL))THEN
794  WRITE(iunit,'(f15.5,50f9.3)') days(inttime),(lcl_probe(i)%SCL)
795  END IF
796 
797  END DO
798 
799  if(dbg_set(dbg_sbr)) write(ipt,*) "END: DUMP_PROBE_DATA"
800  RETURN
type(time) inttime
Definition: mod_main.f90:827
real(dp) function days(MJD)
Definition: mod_time.f90:749
integer ipt
Definition: mod_main.f90:922
Here is the call graph for this function:
Here is the caller graph for this function:

◆ init_nml_probe()

subroutine probes::init_nml_probe ( )

Definition at line 161 of file mod_probe.f90.

161  IMPLICIT NONE
162 
163  probe_interval = 'none'
164  probe_location = -1 !!LOCAL ELEMENT/NODE LOCATION
165  probe_levels = -1 !!FINAL SIGMA LEVEL
166  probe_title = 'none' !!OBJECT TITLE (FOR FILENAMING)
167  probe_description = 'none' !!OBJECT DESCRIPTION
168  probe_variable = 'none' !!VARIABLE TO DUMP
169  probe_var_name = 'none' !!VARIABLE NAMES
170 
Here is the caller graph for this function:

◆ myprobe_arr()

subroutine probes::myprobe_arr ( type(probe_obj PROBE,
real(sp), dimension(:,:), allocatable, target  ARR 
)

Definition at line 359 of file mod_probe.f90.

359  USE all_vars, only:xmc,ymc,h1,h,xm,ym,h,lon,lat,lonc,latc
360  IMPLICIT NONE
361  TYPE(PROBE_OBJ) :: PROBE
362  REAL(SP), ALLOCATABLE, TARGET :: ARR(:,:)
363  CHARACTER(LEN=80):: cstr1,cstr2,cstr3
364  INTEGER :: I,J,IBND,PROCMAX
365  !============ENSURE SIGMA RANGE IS TENABLE=============================
366 
367  if(dbg_set(dbg_sbr)) write(ipt,*) "START: MYPROBE_ARR"
368 
369  IF(probe%K_ONE > ubound(arr,2) .or. probe%K_TWO > ubound(arr,2) )THEN
370  CALL fatal_error('ERROR IN PROBE SETUP: PROBE LEVEL RANGE NOT CORRECT FOR VARIABLE: '//trim(probe_variable),&
371  & 'MAKE SURE PROBE LEVELS ARE LESS THAN OR EQUAL TO THE NUMBER OF MODEL LEVELS')
372  END IF
373 
374  IF(probe%K_ONE < 1 .or. probe%K_TWO <1)THEN
375  CALL fatal_error('ERROR IN PROBE SETUP: PROBE LEVEL RANGE NOT CORRECT FOR VARIABLE: '//trim(probe_variable),&
376  & 'MAKE SURE PROBE LEVELS ARE GREATER THAN OR EQUAL TO ONE')
377  END IF
378 
379 
380  IF(probe_levels(1) > probe_levels(2) )THEN
381  CALL fatal_error&
382  &('ERROR IN PROBE SETUP: PROBE LEVEL RANGE NOT CORRECT FOR VARIABLE: '//trim(probe_variable),&
383  & 'THE PROBE LEVEL INTERVAL MUST SPECIFY A VALID RANGE a:b')
384  END IF
385 
386 
387  IF(probe%D_LOC_GL <1)THEN
388  write(cstr1,'(i8)') probe%D_LOC_GL
389  write(cstr2,'(i8)') ngl
390 
391  CALL fatal_error('ERROR IN PROBE SETUP: DATA LOCATION'//trim(cstr1)//' FOR TIME SERIES FILE: '//trim(probe_variable),&
392  & 'IS NOT IN GLOBAL DOMAIN: 1 --> '//trim(cstr2))
393  END IF
394 
395  IF( ubound(arr,1) == nt)THEN
396 
397  IF(probe%D_LOC_GL > ngl)THEN
398 
399  write(cstr1,'(i8)') probe%D_LOC_GL
400  write(cstr2,'(i8)') ngl
401 
402  CALL fatal_error('ERROR IN PROBE SETUP: DATA LOCATION'//trim(cstr1)//' FOR TIME SERIES FILE: '//trim(probe_variable),&
403  & 'IS NOT IN GLOBAL DOMAIN: 1 --> '//trim(cstr2))
404  END IF
405 
406  IF(elid(probe%D_LOC_GL) /= 0) THEN
407 
408  probe%D_LOC = elid(probe%D_LOC_GL)
409  probe%MINE=.true.
410 
411  ! METERS
412  probe%XLOC = xmc(probe%D_LOC)
413  probe%YLOC = ymc(probe%D_LOC)
414  ! SPHERICAL
415  probe%LONLOC = lonc(probe%D_LOC)
416  probe%LATLOC = latc(probe%D_LOC)
417  ! DEPTH
418  probe%DPTH = h1(probe%D_LOC)
419 
420  probe%VEC => arr(probe%D_LOC,probe%K_ONE:probe%K_TWO)
421  END IF
422 
423 
424  ELSE IF ( ubound(arr,1) == mt)THEN
425 
426  IF(probe%D_LOC_GL > mgl)THEN
427  write(cstr1,'(i8)') probe%D_LOC_GL
428  write(cstr2,'(i8)') mgl
429 
430  CALL fatal_error('ERROR IN PROBE SETUP: DATA LOCATION'//trim(cstr1)//' FOR TIME SERIES FILE: '//trim(probe_variable),&
431  & 'IS NOT IN GLOBAL DOMAIN: 1 --> '//trim(cstr2))
432 
433  END IF
434 
435  IF(nlid(probe%D_LOC_GL) == 0) RETURN
436 
437  IF(nlid(probe%D_LOC_GL) > 0) THEN
438 
439  IF(nde_id(nlid(probe%D_LOC_GL)) == 1)THEN !!BOUNDARY NODE
440  DO j=1,nbn
441  IF(bn_lst(j) == probe%D_LOC_GL) ibnd = j
442  END DO
443  !----Choose Processor of Lowest ID to be responsible for node
444  procmax = 10000
445  DO j=1,nprocs
446  IF(bn_ney(ibnd,j)==1) THEN
447  IF(j < procmax) procmax = j
448  END IF
449  END DO
450 
451  IF(procmax == myid) THEN
452  probe%MINE = .true. !!NOT RESPONSIBLE FOR TIME SERIES
453  probe%D_LOC = nlid(probe%D_LOC_GL)
454 
455  ! METERS
456  probe%XLOC = xm(probe%D_LOC)
457  probe%YLOC = ym(probe%D_LOC)
458  ! SPHERICAL
459  probe%LONLOC = lon(probe%D_LOC)
460  probe%LATLOC = lat(probe%D_LOC)
461  ! DEPTH
462  probe%DPTH = h(probe%D_LOC)
463 
464  probe%VEC => arr(probe%D_LOC,probe%K_ONE:probe%K_TWO)
465  END IF
466  ELSE
467  probe%MINE = .true. !!NOT RESPONSIBLE FOR TIME SERIES
468  probe%D_LOC = nlid(probe%D_LOC_GL)
469 
470  ! METERS
471  probe%XLOC = xm(probe%D_LOC)
472  probe%YLOC = ym(probe%D_LOC)
473  ! SPHERICAL
474  probe%LONLOC = lon(probe%D_LOC)
475  probe%LATLOC = lat(probe%D_LOC)
476  ! DEPTH
477  probe%DPTH = h(probe%D_LOC)
478 
479  probe%VEC => arr(probe%D_LOC,probe%K_ONE:probe%K_TWO)
480 
481  END IF
482  END IF
483 
484 
485  ELSE
486  CALL fatal_error('MYPROBE: INVALID VARIABLE SIZE?'&
487  &,'Variable:'//trim(probe%VAR))
488  END IF
489 
490  if(dbg_set(dbg_sbr)) write(ipt,*) "END: MYPROBE_ARR"
491 
integer nbn
Definition: mod_par.f90:75
real(sp), dimension(:), allocatable, target h
Definition: mod_main.f90:1131
integer mt
Definition: mod_main.f90:78
integer, dimension(:), pointer elid
Definition: mod_par.f90:53
integer myid
Definition: mod_main.f90:67
integer, target nprocs
Definition: mod_main.f90:72
real(sp), dimension(:), allocatable, target ymc
Definition: mod_main.f90:994
integer, dimension(:), pointer nde_id
Definition: mod_par.f90:81
real(sp), dimension(:), allocatable, target latc
Definition: mod_main.f90:998
integer, dimension(:), pointer nlid
Definition: mod_par.f90:54
real(sp), dimension(:), allocatable, target xmc
Definition: mod_main.f90:993
real(sp), dimension(:), allocatable, target lonc
Definition: mod_main.f90:997
integer, dimension(:), pointer bn_lst
Definition: mod_par.f90:77
real(sp), dimension(:), allocatable, target xm
Definition: mod_main.f90:991
integer mgl
Definition: mod_main.f90:50
real(sp), dimension(:), allocatable, target lat
Definition: mod_main.f90:996
real(sp), dimension(:), allocatable, target h1
Definition: mod_main.f90:1115
real(sp), dimension(:), allocatable, target lon
Definition: mod_main.f90:995
integer, dimension(:,:), pointer bn_ney
Definition: mod_par.f90:80
integer ipt
Definition: mod_main.f90:922
integer nt
Definition: mod_main.f90:77
integer ngl
Definition: mod_main.f90:49
real(sp), dimension(:), allocatable, target ym
Definition: mod_main.f90:992
Here is the call graph for this function:

◆ myprobe_vec()

subroutine probes::myprobe_vec ( type(probe_obj PROBE,
real(sp), dimension(:), allocatable, target  VEC 
)

Definition at line 239 of file mod_probe.f90.

239  USE all_vars, only:xmc,ymc,h1,h,xm,ym,h,lon,lat,lonc,latc
240  IMPLICIT NONE
241  TYPE(PROBE_OBJ) :: PROBE
242  REAL(SP), ALLOCATABLE, TARGET :: VEC(:)
243  CHARACTER(LEN=80):: cstr1,cstr2,cstr3
244  INTEGER :: I,J,IBND,PROCMAX
245  !============SEE IF DATA POINT IS IN THE GLOBAL DOMAIN=================
246 
247  if(dbg_set(dbg_sbr)) write(ipt,*) "START: MYPROBE_VEC"
248 
249  IF (probe%K_ONE /= -1 .or. probe%K_TWO /= -1) CALL fatal_error&
250  &("ERROR IN PROBE SETUP: PROBE_LEVELS should not be set for vector variables",&
251  & "Do not specify it in the PROBE Namelist object for "//trim(probe_variable))
252 
253  IF(probe%D_LOC_GL <1)THEN
254  write(cstr1,'(i8)') probe%D_LOC_GL
255  write(cstr2,'(i8)') ngl
256 
257  CALL fatal_error('ERROR IN PROBE SETUP: DATA LOCATION'//trim(cstr1)//' FOR TIME SERIES FILE: '//trim(probe_variable),&
258  & 'IS NOT IN GLOBAL DOMAIN: 1 --> '//trim(cstr2))
259  END IF
260 
261  IF( ubound(vec,1) == nt)THEN
262 
263  IF(probe%D_LOC_GL > ngl)THEN
264 
265  write(cstr1,'(i8)') probe%D_LOC_GL
266  write(cstr2,'(i8)') ngl
267 
268  CALL fatal_error('ERROR IN PROBE SETUP: DATA LOCATION'//trim(cstr1)//' FOR TIME SERIES FILE: '//trim(probe_variable),&
269  & 'IS NOT IN GLOBAL DOMAIN: 1 --> '//trim(cstr2))
270  END IF
271 
272  IF(elid(probe%D_LOC_GL) /= 0) THEN
273 
274  probe%D_LOC = elid(probe%D_LOC_GL)
275  probe%MINE=.true.
276 
277  ! METERS
278  probe%XLOC = xmc(probe%D_LOC)
279  probe%YLOC = ymc(probe%D_LOC)
280  ! SPHERICAL
281  probe%LONLOC = lonc(probe%D_LOC)
282  probe%LATLOC = latc(probe%D_LOC)
283  ! DEPTH
284  probe%DPTH = h1(probe%D_LOC)
285 
286  probe%SCL => vec(probe%D_LOC)
287  END IF
288 
289  ELSE IF ( ubound(vec,1) == mt)THEN
290 
291  IF(probe%D_LOC_GL > mgl)THEN
292  write(cstr1,'(i8)') probe%D_LOC_GL
293  write(cstr2,'(i8)') mgl
294 
295  CALL fatal_error('ERROR IN PROBE SETUP: DATA LOCATION'//trim(cstr1)//' FOR TIME SERIES FILE: '//trim(probe_variable),&
296  & 'IS NOT IN GLOBAL DOMAIN: 1 --> '//trim(cstr2))
297 
298  END IF
299 
300  IF(nlid(probe%D_LOC_GL) == 0) RETURN
301 
302  IF(nlid(probe%D_LOC_GL) > 0) THEN
303 
304  IF(nde_id(nlid(probe%D_LOC_GL)) == 1)THEN !!BOUNDARY NODE
305  DO j=1,nbn
306  IF(bn_lst(j) == probe%D_LOC_GL) ibnd = j
307  END DO
308  !----Choose Processor of Lowest ID to be responsible for node
309  procmax = 10000
310  DO j=1,nprocs
311  IF(bn_ney(ibnd,j)==1) THEN
312  IF(j < procmax) procmax = j
313  END IF
314  END DO
315 
316  IF(procmax == myid) THEN
317  probe%MINE = .true. !!NOT RESPONSIBLE FOR TIME SERIES
318  probe%D_LOC = nlid(probe%D_LOC_GL)
319 
320  ! METERS
321  probe%XLOC = xm(probe%D_LOC)
322  probe%YLOC = ym(probe%D_LOC)
323  ! SPHERICAL
324  probe%LONLOC = lon(probe%D_LOC)
325  probe%LATLOC = lat(probe%D_LOC)
326  ! DEPTH
327  probe%DPTH = h(probe%D_LOC)
328 
329  probe%SCL => vec(probe%D_LOC)
330  END IF
331  ELSE
332  probe%MINE = .true. !!NOT RESPONSIBLE FOR TIME SERIES
333  probe%D_LOC = nlid(probe%D_LOC_GL)
334 
335  ! METERS
336  probe%XLOC = xm(probe%D_LOC)
337  probe%YLOC = ym(probe%D_LOC)
338  ! SPHERICAL
339  probe%LONLOC = lon(probe%D_LOC)
340  probe%LATLOC = lat(probe%D_LOC)
341  ! DEPTH
342  probe%DPTH = h(probe%D_LOC)
343 
344  probe%SCL => vec(probe%D_LOC)
345 
346  END IF
347  END IF
348 
349 
350  ELSE
351  CALL fatal_error('MYPROBE: INVALID VARIABLE SIZE (Not equal MT or NT?)'&
352  &,'Variable:'//trim(probe%VAR))
353  END IF
354  if(dbg_set(dbg_sbr)) write(ipt,*) "END: MYPROBE_VEC"
355 
integer nbn
Definition: mod_par.f90:75
real(sp), dimension(:), allocatable, target h
Definition: mod_main.f90:1131
integer mt
Definition: mod_main.f90:78
integer, dimension(:), pointer elid
Definition: mod_par.f90:53
integer myid
Definition: mod_main.f90:67
integer, target nprocs
Definition: mod_main.f90:72
real(sp), dimension(:), allocatable, target ymc
Definition: mod_main.f90:994
integer, dimension(:), pointer nde_id
Definition: mod_par.f90:81
real(sp), dimension(:), allocatable, target latc
Definition: mod_main.f90:998
integer, dimension(:), pointer nlid
Definition: mod_par.f90:54
real(sp), dimension(:), allocatable, target xmc
Definition: mod_main.f90:993
real(sp), dimension(:), allocatable, target lonc
Definition: mod_main.f90:997
integer, dimension(:), pointer bn_lst
Definition: mod_par.f90:77
real(sp), dimension(:), allocatable, target xm
Definition: mod_main.f90:991
integer mgl
Definition: mod_main.f90:50
real(sp), dimension(:), allocatable, target lat
Definition: mod_main.f90:996
real(sp), dimension(:), allocatable, target h1
Definition: mod_main.f90:1115
real(sp), dimension(:), allocatable, target lon
Definition: mod_main.f90:995
integer, dimension(:,:), pointer bn_ney
Definition: mod_par.f90:80
integer ipt
Definition: mod_main.f90:922
integer nt
Definition: mod_main.f90:77
integer ngl
Definition: mod_main.f90:49
real(sp), dimension(:), allocatable, target ym
Definition: mod_main.f90:992
Here is the call graph for this function:

◆ open_probes()

subroutine probes::open_probes ( )

Definition at line 696 of file mod_probe.f90.

696 
697  !------------------------------------------------------------------------------|
698  ! CREATE FILE NAMES AND WRITE HEADER INFORMATION FOR EACH TS OBJECT |
699  !------------------------------------------------------------------------------|
700 
701  USE mod_prec
702  USE all_vars
703  IMPLICIT NONE
704  INTEGER I,IUNIT,ICNT
705  CHARACTER(LEN=120) :: FNAME
706  CHARACTER(LEN=2 ) :: NAC
707  CHARACTER(LEN=3 ) :: APPEND
708  LOGICAL FEXIST
709 
710  if(dbg_set(dbg_sbr)) write(ipt,*) "START: OPEN_PROBES"
711  !
712  !--Open Up Files -> If File Exists Create Secondary File (-01,-02, etc)
713  !
714 
715  DO i=1,local_probes
716 
717  icnt = 0
718  fname = trim(output_dir)//trim(lcl_probe(i)%D_TIT)
719  INQUIRE(file=fname,exist=fexist)
720  DO WHILE(fexist)
721  icnt = icnt + 1
722  IF(icnt .GE. 100) CALL fatal_error&
723  &("Please clean old time seris output in your results directory!")
724  WRITE(nac,'(I2.2)')icnt
725  append = "-"//nac
726  fname = trim(output_dir)//trim(lcl_probe(i)%D_TIT)//trim(append)
727  INQUIRE(file=fname,exist=fexist)
728  END DO
729 
730  iunit = lcl_probe(i)%O_NUM + 100
731  CALL fopen(iunit,fname,'ofr')
732  WRITE(iunit,*)trim(lcl_probe(i)%D_DES)
733  WRITE(iunit,*)trim(lcl_probe(i)%VNAME)
734  WRITE(iunit,*)
735  CALL print_real_time(get_now(),iunit,"MODEL START DATE")
736  WRITE(iunit,*)
737  WRITE(iunit,*)' K1 K2 '
738  WRITE(iunit,'(2(I12,3X))')lcl_probe(i)%K_ONE,lcl_probe(i)%K_TWO
739  WRITE(iunit,*)' X(M) Y(M) DEPTH(M)'
740  WRITE(iunit,'(3(F12.3,3X))')lcl_probe(i)%XLOC,lcl_probe(i)%YLOC,lcl_probe(i)%DPTH
741  WRITE(iunit,*)' LON LAT DEPTH(M)'
742  WRITE(iunit,'(3(F12.3,3X))')lcl_probe(i)%LONLOC,lcl_probe(i)%LATLOC,lcl_probe(i)%DPTH
743  WRITE(iunit,*)
744  WRITE(iunit,*)'DATA FOLLOWS:'
745  WRITE(iunit,*)'Time(days) Data...'
746  CLOSE(iunit)
747  lcl_probe(i)%FILENAME = fname
748  END DO
749 
750  if(dbg_set(dbg_sbr)) write(ipt,*) "END: OPEN_PROBES"
751 
type(time) function get_now()
Definition: mod_time.f90:716
character(len=80) output_dir
Definition: mod_main.f90:184
subroutine print_real_time(mjd, IPT, char, TZONE)
Definition: mod_time.f90:1201
integer ipt
Definition: mod_main.f90:922
Here is the call graph for this function:
Here is the caller graph for this function:

◆ probe_store()

subroutine probes::probe_store ( type(probe_obj PROBE)

Definition at line 813 of file mod_probe.f90.

813 
814  !------------------------------------------------------------------------------|
815  ! PUT TIME SERIES DATA IN TEMPORARY ARRAY |
816  !------------------------------------------------------------------------------|
817 
818  USE mod_prec
819  USE all_vars
820  IMPLICIT NONE
821  TYPE(PROBE_OBJ) :: PROBE
822 
823  if(dbg_set(dbg_sbr)) write(ipt,*) "START: PROBE_STORE"
824 
825  !--Store Data In Temporary Array-------------------------------------------
826  SELECT CASE(trim(probe%VAR))
827 
828  CASE("u")
829  CALL myprobe(probe,u)
830  CASE("v")
831  CALL myprobe(probe,v)
832  CASE("w")
833  CALL myprobe(probe,w)
834  CASE("ww")
835  CALL myprobe(probe,ww)
836  CASE("q2")
837  CALL myprobe(probe,q2)
838  CASE("q2l")
839  CALL myprobe(probe,q2l)
840  CASE("l")
841  CALL myprobe(probe,l)
842  CASE("km")
843  CALL myprobe(probe,km)
844  CASE("kq")
845  CALL myprobe(probe,kq)
846  CASE("kh")
847  CALL myprobe(probe,kh)
848  CASE("t1")
849  CALL myprobe(probe,t1)
850  CASE("s1")
851  CALL myprobe(probe,s1)
852  CASE("rho1")
853  CALL myprobe(probe,rho1)
854  CASE("ua")
855  CALL myprobe(probe,ua)
856  CASE("va")
857  CALL myprobe(probe,va)
858  CASE("el")
859  CALL myprobe(probe,el)
860  CASE DEFAULT
861  CALL fatal_error&
862  &('VARIABLE: '//trim(probe%VAR)//' HAS NOT BEEN SET UP',&
863  & 'FOR TIME SERIES OUTPUT (Did you use CAPITALS by mistake?)',&
864  & 'MODIFY MOD_PROBE TO ADD IT!')
865 
866  END SELECT
867 
868  if(dbg_set(dbg_sbr)) write(ipt,*) "END: PROBE_STORE"
869 
real(sp), dimension(:,:), allocatable, target q2
Definition: mod_main.f90:1290
real(sp), dimension(:,:), allocatable, target km
Definition: mod_main.f90:1293
real(sp), dimension(:), allocatable, target va
Definition: mod_main.f90:1104
real(sp), dimension(:), allocatable, target el
Definition: mod_main.f90:1134
real(sp), dimension(:,:), allocatable, target v
Definition: mod_main.f90:1269
real(sp), dimension(:,:), allocatable, target rho1
Definition: mod_main.f90:1309
real(sp), dimension(:,:), allocatable, target t1
Definition: mod_main.f90:1307
real(sp), dimension(:,:), allocatable, target w
Definition: mod_main.f90:1279
real(sp), dimension(:,:), allocatable, target ww
Definition: mod_main.f90:1280
real(sp), dimension(:,:), allocatable, target q2l
Definition: mod_main.f90:1292
real(sp), dimension(:,:), allocatable, target u
Definition: mod_main.f90:1268
real(sp), dimension(:,:), allocatable, target s1
Definition: mod_main.f90:1308
real(sp), dimension(:), allocatable, target ua
Definition: mod_main.f90:1103
real(sp), dimension(:,:), allocatable, target l
Definition: mod_main.f90:1291
real(sp), dimension(:,:), allocatable, target kh
Definition: mod_main.f90:1294
integer ipt
Definition: mod_main.f90:922
real(sp), dimension(:,:), allocatable, target kq
Definition: mod_main.f90:1295
Here is the call graph for this function:
Here is the caller graph for this function:

◆ set_probes()

subroutine probes::set_probes ( logical, intent(in)  P_ON,
integer, intent(in)  NP,
character(len=*), intent(in)  FNM 
)

Definition at line 497 of file mod_probe.f90.

497  !------------------------------------------------------------------------------|
498  ! READ IN TIME SERIES OBJECTS FROM INPUT |
499  !------------------------------------------------------------------------------|
500  USE mod_par
501  USE lims
502  USE mod_set_time
503 
504  IMPLICIT NONE
505  LOGICAL, INTENT(IN) :: P_ON
506  INTEGER, INTENT(IN) :: NP
507  CHARACTER(LEN=*), INTENT(IN):: FNM
508 
509  CHARACTER(LEN=80):: cstr1,cstr2,cstr3
510  LOGICAL FEXIST,ISLOCAL
511  INTEGER :: I,J,IERR,IOS,STATUS, N_PROBE
512  INTEGER :: PROCMAX, IBND, charnum
513  CHARACTER(LEN=120) :: pathnfile
514  CHARACTER(LEN=4) :: OFLAG
515  TYPE(TIME) :: OTIME
516  INTEGER(ITIME) :: OSTEP
517  if(dbg_set(dbg_sbr)) &
518  & write(ipt,*) "START: SET_PROBES;"
519 
520 
521  probe_on=p_on
522  global_probes=np
523  local_probes=0
524  input_fname=fnm
525 
526  IF (.not. probe_on) THEN
527  if(dbg_set(dbg_log)) &
528  & write(ipt,*) "! Time Series Probes are off"
529  return
530  ELSE
531  if(dbg_set(dbg_log)) then
532  write(ipt,*) "! Time Series Probes are on"
533  write(ipt,*) "! Setting up Probes:"
534  end if
535  END IF
536 
537 
538 
539  charnum = index(probes_file,".nml")
540  if (charnum /= len_trim(probes_file)-3)&
541  & CALL warning("PROBES FILE does not end in .nml", &
542  & trim(probes_file))
543  ! OPEN FILE - try both with appending input dir and without!
544  pathnfile = trim(input_dir)//trim(probes_file)
545  INQUIRE(file=pathnfile,exist=fexist)
546  IF(fexist) THEN
547  Call fopen(probeunit,trim(pathnfile),'cfr')
548  ELSE
549  pathnfile = trim(probes_file)
550  Call fopen(probeunit,trim(pathnfile),'cfr')
551  ! LET FOPEN CALL ERROR IF FILES DOES NOT EXIST
552  END IF
553 
554 
555  CALL alloc_probe(glb_probe,global_probes)
556 
557  !------------------------------------------------------------------------------|
558  ! READ TIME SERIES NAME LIST |
559  !------------------------------------------------------------------------------|
560  n_probe = 0
561  DO
562 
563  CALL init_nml_probe
564 
565  islocal = .false.
566 
567  READ(unit=probeunit, nml=nml_probe,iostat=ios)
568  if (ios /= 0 ) exit
569 
570  n_probe = n_probe +1
571  if (n_probe > global_probes) exit ! To prevent sigsev...
572 
573 
574  !------------------------------------------------------------------------------|
575  ! SET FUNDEMENTAL DATA
576  !------------------------------------------------------------------------------|
577  glb_probe(n_probe)%K_ONE = probe_levels(1)
578  glb_probe(n_probe)%K_TWO = probe_levels(2)
579  glb_probe(n_probe)%D_LOC_GL = probe_location
580 
581  glb_probe(n_probe)%D_TIT =probe_title
582  glb_probe(n_probe)%D_DES =probe_description
583  glb_probe(n_probe)%VAR =probe_variable
584  glb_probe(n_probe)%VNAME =probe_var_name
585 
586  !------------------------------------------------------------------------------|
587  ! SET TIME FOR PROBES
588  !------------------------------------------------------------------------------|
589  glb_probe(n_probe)%O_NEXT = starttime
590  CALL ideal_time_string2time(probe_interval,oflag,otime,ostep)
591  IF (oflag == 'time') THEN ! IF OUTPUT TIME INTERVAL WAS SPECIFIED
592  glb_probe(n_probe)%O_INT = otime
593 
594  ELSE IF(oflag == 'step') THEN ! IF OUTPUT CYCLE INTERVAL WAS SPECIFIED
595  glb_probe(n_probe)%O_INT = imdti * ostep
596 
597  END IF
598 
599  IF (glb_probe(n_probe)%O_INT <= zerotime) CALL fatal_error&
600  &('ERROR IN PROBE SETUP: Time series output interval is less than or equal to zero!')
601 
602  !------------------------------------------------------------------------------|
603  ! POINT STORAGE DATA AND LOCATION TYPE |
604  !------------------------------------------------------------------------------|
605  CALL probe_store(glb_probe(n_probe))
606 
607  IF (glb_probe(n_probe)%MINE) THEN
608  local_probes=local_probes+1
609 
610  glb_probe(n_probe)%O_NUM =n_probe
611 
612 
613  END IF
614 
615  END DO
616 
617  IF(global_probes .ne. n_probe) THEN
618 
619  if(dbg_set(dbg_log)) then
620  write(ipt,*)"Bad NML_PROBE in the Name List!"
621  write(ipt,*)"Specified number of PROBES=",global_probes
622  write(ipt,*)"But Found",n_probe, "; Valid PROBE name list objects.(Printing Last)"
623  write(unit=ipt,nml=nml_probe)
624  end if
625 
626  CALL fatal_error&
627  &('THE NUMBER OF PROBES SPECIFIED IN THE RUN FILE CAN',&
628  &'NOT BE FOUND IN THE PROBE FILE:'//trim(probes_file))
629  END IF
630 
631  CLOSE(probeunit)
632 
633 
634  CALL alloc_probe(lcl_probe,local_probes)
635 
636  n_probe=0
637  DO i = 1,global_probes
638  IF (glb_probe(i)%O_NUM .NE. 0) THEN
639  n_probe=n_probe+1
640  lcl_probe(n_probe)=glb_probe(i)
641  END IF
642  END DO
643  IF(n_probe .NE. local_probes) CALL fatal_error("mod_probe: this should not happen?")
644 
645 
646  DEALLOCATE(glb_probe)
647 
648  !------------------------------------------------------------------------------|
649  ! PRINT STATISTICS ON TIME SERIES OBJECTS TO OUTPUT |
650  !------------------------------------------------------------------------------|
651  IF(global_probes > 0)THEN
652 
653  if(dbg_set(dbg_sbrio)) write(ipt,*) "GLobal probes"&
654  &,global_probes,"Local_probes",local_probes,size(lcl_probe)
655 
656  IF(dbg_set(dbg_log))THEN
657  WRITE(ipt,*)
658  WRITE(ipt,*)'! TIME SERIES OBJECT DATA '
659  WRITE(ipt,*)" OBJ# PROC GLOBAL LOCAL VAR FILENAME"
660  END IF
661 
662 
663  ! THIS IS VERY SLOW - DO NOT USE THIS METHOD INSIDE THE MAIN LOOP!
664  DO j=1,nprocs
665 
666  IF(myid == j) THEN
667  DO i=1,local_probes
668  WRITE(ipt,101)lcl_probe(i)%O_NUM,myid,lcl_probe(i)%D_LOC_GL,lcl_probe(i)%D_LOC, &
669  & trim(lcl_probe(i)%VAR),trim(lcl_probe(i)%D_TIT)
670  END DO
671  END IF
672  END DO
673 
674  END IF
675 
676  !------------------------------------------------------------------------------|
677  ! OPEN UP OUTPUT FILES AND WRITE HEADER INFORMATION |
678  !------------------------------------------------------------------------------|
679  CALL open_probes
680 
681  CALL dump_probe_data
682 
683  if(dbg_set(dbg_sbr)) write(ipt,*) "END: SET_PROBES"
684 
685 101 FORMAT(i5,i5,i8,i8,3x,a6,1x,a30)
type(time) zerotime
Definition: mod_main.f90:830
integer myid
Definition: mod_main.f90:67
integer, target nprocs
Definition: mod_main.f90:72
type(time) starttime
Definition: mod_main.f90:833
integer, parameter probeunit
Definition: mod_main.f90:937
type(time) imdti
Definition: mod_main.f90:848
subroutine ideal_time_string2time(string, flag, ntime, tstep)
character(len=80) input_dir
Definition: mod_main.f90:183
character(len=80) probes_file
Definition: mod_main.f90:796
integer ipt
Definition: mod_main.f90:922
integer ios
Definition: mod_obcs2.f90:81
Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ glb_probe

type(probe_obj), dimension(:), allocatable probes::glb_probe

Definition at line 143 of file mod_probe.f90.

143  TYPE(PROBE_OBJ), Allocatable :: GLB_PROBE(:) !!GLOBAL PROBE OBJECTS

◆ global_probes

integer probes::global_probes

Definition at line 147 of file mod_probe.f90.

147  INTEGER :: GLOBAL_PROBES, LOCAL_PROBES

◆ input_fname

character probes::input_fname

Definition at line 148 of file mod_probe.f90.

148  CHARACTER :: INPUT_FNAME

◆ lcl_probe

type(probe_obj), dimension(:), allocatable probes::lcl_probe

Definition at line 144 of file mod_probe.f90.

144  TYPE(PROBE_OBJ), Allocatable :: LCL_PROBE(:) !!LOCAL PROBE OBJECTS

◆ local_probes

integer probes::local_probes

Definition at line 147 of file mod_probe.f90.

◆ probe_description

character(len=80) probes::probe_description

Definition at line 127 of file mod_probe.f90.

127  CHARACTER(LEN=80) :: PROBE_DESCRIPTION !!OBJECT DESCRIPTION

◆ probe_interval

character(len=80) probes::probe_interval

Definition at line 123 of file mod_probe.f90.

123  CHARACTER(LEN=80) :: PROBE_INTERVAL

◆ probe_levels

integer, dimension(2) probes::probe_levels

Definition at line 125 of file mod_probe.f90.

125  INTEGER :: PROBE_LEVELS(2) !!FINAL SIGMA LEVEL

◆ probe_location

integer probes::probe_location

Definition at line 124 of file mod_probe.f90.

124  INTEGER :: PROBE_LOCATION !!LOCAL ELEMENT/NODE LOCATION

◆ probe_on

logical probes::probe_on

Definition at line 146 of file mod_probe.f90.

146  LOGICAL :: PROBE_ON

◆ probe_title

character(len=80) probes::probe_title

Definition at line 126 of file mod_probe.f90.

126  CHARACTER(LEN=80) :: PROBE_TITLE !!OBJECT TITLE (FOR FILENAMING)

◆ probe_var_name

character(len=80) probes::probe_var_name

Definition at line 129 of file mod_probe.f90.

129  CHARACTER(LEN=80) :: PROBE_VAR_NAME !!VARIABLE NAMES

◆ probe_variable

character(len=80) probes::probe_variable

Definition at line 128 of file mod_probe.f90.

128  CHARACTER(LEN=80) :: PROBE_VARIABLE !!VARIABLE TO DUMP