My Project
Public Member Functions | List of all members
probes::myprobe Interface Reference

Public Member Functions

subroutine myprobe_vec (PROBE, VEC)
 
subroutine myprobe_arr (PROBE, ARR)
 

Detailed Description

Definition at line 114 of file mod_probe.f90.

Member Function/Subroutine Documentation

◆ myprobe_arr()

subroutine probes::myprobe::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

◆ myprobe_vec()

subroutine probes::myprobe::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

The documentation for this interface was generated from the following file: