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

Data Types

interface  degrees2meters
 
interface  interp_anodal
 
interface  interp_azonal
 
interface  interp_pnodal
 
interface  interp_pzonal
 
interface  meters2degrees
 
interface  shutdown_check
 

Functions/Subroutines

subroutine initialize_control (NAME)
 
logical function dbg_set (vrb)
 
subroutine dbg_init (IPT_BASE, outtofile)
 
subroutine fatal_error (ER1, ER2, ER3, ER4)
 
subroutine warning (ER1, ER2, ER3, ER4)
 
subroutine pstop
 
subroutine pshutdown
 
subroutine shutdown_check_1d (VAR, MSG)
 
subroutine shutdown_check_2d (VAR, MSG)
 
logical function have_proj (proj_ref)
 
subroutine degrees2meters_scl_flt (LON, LAT, proj_ref, X, Y)
 
subroutine degrees2meters_vec_flt (LON, LAT, proj_ref, X, Y, nsze)
 
subroutine degrees2meters_arr_flt (LON, LAT, proj_ref, X, Y, nsze, msze)
 
subroutine degrees2meters_scl_dbl (LON, LAT, proj_ref, X, Y)
 
subroutine degrees2meters_vec_dbl (LON, LAT, proj_ref, X, Y, nsze)
 
subroutine degrees2meters_arr_dbl (LON, LAT, proj_ref, X, Y, nsze, msze)
 
subroutine meters2degrees_scl_flt (X, Y, proj_ref, LON, LAT)
 
subroutine meters2degrees_vec_flt (X, Y, proj_ref, LON, LAT, nsze)
 
subroutine meters2degrees_arr_flt (X, Y, proj_ref, LON, LAT, nsze, msze)
 
subroutine meters2degrees_scl_dbl (X, Y, proj_ref, LON, LAT)
 
subroutine meters2degrees_vec_dbl (X, Y, proj_ref, LON, LAT, nsze)
 
subroutine meters2degrees_arr_dbl (X, Y, proj_ref, LON, LAT, nsze, msze)
 
real(spa) function interp_anodal_2d_flt (xloc, yloc, i, Field)
 
real(dp) function interp_anodal_2d_dbl (xloc, yloc, i, Field)
 
real(spa) function interp_pnodal_2d_flt (xloc, yloc, i, Field)
 
real(dp) function interp_pnodal_2d_dbl (xloc, yloc, i, Field)
 
real(spa) function interp_anodal_3d_flt (xloc, yloc, sigloc, lvls, i, Field)
 
real(dp) function interp_anodal_3d_dbl (xloc, yloc, sigloc, lvls, i, Field)
 
real(spa) function interp_pnodal_3d_flt (xloc, yloc, sigloc, lvls, i, Field)
 
real(dp) function interp_pnodal_3d_dbl (xloc, yloc, sigloc, lvls, i, Field)
 
real(spa) function interp_azonal_2d_flt (xloc, yloc, i, Field)
 
real(dp) function interp_azonal_2d_dbl (xloc, yloc, i, Field)
 
real(spa) function interp_pzonal_2d_flt (xloc, yloc, i, Field)
 
real(dp) function interp_pzonal_2d_dbl (xloc, yloc, i, Field)
 
real(spa) function interp_azonal_3d_flt (xloc, yloc, sigloc, lvls, i, Field)
 
real(dp) function interp_azonal_3d_dbl (xloc, yloc, sigloc, lvls, i, Field)
 
real(spa) function interp_pzonal_3d_flt (xloc, yloc, sigloc, lvls, i, Field)
 
real(dp) function interp_pzonal_3d_dbl (xloc, yloc, sigloc, lvls, i, Field)
 
integer function find_element_containing (xloc, yloc, GUESS)
 
integer function find_element_containing_robust (xloc, yloc)
 
integer function find_element_containing_quick (xloc, yloc, Guess)
 
logical function isintri (X0, Y0, Xt, Yt)
 
logical function isintriangle (i, x0, y0)
 
logical function sameside (p1, p2, a, b)
 
subroutine grid_neighbor_index (FOUND, IDEX, CNT, ORDER)
 
subroutine write_banner (PAR, NP, ID)
 
subroutine fopen (IUNIT, INSTR, IOPT)
 
integer function open_dat (FNAME, UNIT, PATH)
 
subroutine get_value (LNUM, NUMCHAR, TEXT_LINE, VARNAME, VARTYPE, LOGVAL, STRINGVAL, REALVAL, INTVAL, NVAL)
 
integer function scan_file (UNIT, VNAME, ISCAL, FSCAL, IVEC, FVEC, CVEC, NSZE, CVAL, LVAL)
 
subroutine split_string (instring, delim, outstrings)
 
subroutine path_split (STRING, PATH, FILE, EXTENSION)
 
subroutine test_split_strings
 
real(sp) function limled (a, b, q)
 
real(sp) function limled1 (a, b, alpha)
 
real(sp) function limled2 (a, b, alpha)
 
real(dp) function read_float (ITEM, IERR)
 
real(sp) function read_int (ITEM, IERR)
 
integer function scan_file2 (FNAME, VNAME, ISCAL, FSCAL, IVEC, FVEC, CVEC, NSZE, CVAL, LVAL)
 
integer function scan_file3 (FNAME, VNAME, ISCAL, FSCAL, IVEC, FVEC, CVEC, NSZE, CVAL, LVAL)
 

Variables

integer, parameter sync_tag =3601
 
integer, parameter ext_code =3602
 
integer, parameter wait_code =3603
 
integer restart_code
 
integer nc_code
 
integer ncav_code
 
integer init_code
 
integer nesting_code
 
integer initnest_code
 
integer, parameter dbg_nbr =9
 
integer, parameter dbg_log =0
 
integer, parameter dbg_io =1
 
integer, parameter dbg_scl =2
 
integer, parameter dbg_mpi =3
 
integer, parameter dbg_sbr =4
 
integer, parameter dbg_sbrio =5
 
integer, parameter dbg_vec =6
 
integer, parameter dbg_vrb =7
 
character(len=80) prg_nm
 
integer dbg_lvl
 
logical dbg_par =.false.
 

Function/Subroutine Documentation

◆ dbg_init()

subroutine mod_utils::dbg_init ( integer, intent(in)  IPT_BASE,
logical, intent(in)  outtofile 
)

Definition at line 197 of file mod_utils.f90.

197  use control, only: infofile, prg_name, msr
198  use lims, only: myid
199  implicit none
200  integer, intent(in):: IPT_BASE
201  logical, intent(in):: outtofile
202  character(LEN=3) :: ch3
203  character(len=100) :: debugname
204 
205  if (outtofile .AND. msr) then
206  WRITE(ipt,*)"========================================================================"
207  WRITE(ipt,*)"=== All further standard output goes to the user specified log file ===="
208  WRITE(ipt,*)"=== Any further standard error messages will still print to screen ====="
209  WRITE(ipt,*)"=== LOG FILE NAME: "//trim(infofile)
210  WRITE(ipt,*)"========================================================================"
211  ipt = ipt_base+1
212  CALL fopen(ipt, trim(infofile) ,"ofr")
213 
214  end if
215 
216  if (dbg_par .AND. .NOT. msr) then
217  ipt = ipt_base + myid
218  write(ch3,'(i3.3)') myid
219  debugname=trim(prg_name)//"_DEBUG."&
220  & // trim(adjustl(ch3)) // ".log"
221 
222  CALL fopen(ipt, trim(debugname) ,"ofr")
223 
224  end if
225 
logical msr
Definition: mod_main.f90:101
integer myid
Definition: mod_main.f90:67
character(len=80) infofile
Definition: mod_main.f90:117
character(len=80) prg_name
Definition: mod_main.f90:105
Here is the call graph for this function:
Here is the caller graph for this function:

◆ dbg_set()

logical function mod_utils::dbg_set ( integer, intent(in)  vrb)

Definition at line 182 of file mod_utils.f90.

182  USE control, only: msr
183  IMPLICIT NONE
184  INTEGER, INTENT(IN) :: vrb
185 
186  ! ONLY RETURN TRUE IF LEVEL PASSED IS .LE. THE RUN TIME LEVEL
187  dbg_set =.false.
188  if (vrb <= dbg_lvl) then
189  if ( msr .OR. dbg_par) dbg_set =.true.
190  end if
191 
192  return
logical msr
Definition: mod_main.f90:101

◆ degrees2meters_arr_dbl()

subroutine mod_utils::degrees2meters_arr_dbl ( real(dp), dimension(nsze,msze), intent(in)  LON,
real(dp), dimension(nsze,msze), intent(in)  LAT,
character(len=*), intent(in)  proj_ref,
real(dp), dimension(nsze,msze), intent(out)  X,
real(dp), dimension(nsze,msze), intent(out)  Y,
integer, intent(in)  nsze,
integer, intent(in)  msze 
)

Definition at line 435 of file mod_utils.f90.

435  implicit none
436  integer, intent(in) :: nsze, msze
437  REAL(DP), INTENT(IN) :: LON(nsze,msze), LAT(nsze,msze)
438  character(LEN=*), INTENT(IN) :: proj_ref
439  REAL(DP), INTENT(OUT) :: X(nsze,msze), Y(nsze,msze)
440  x=-1.0_dp
441  y=-1.0_dp

◆ degrees2meters_arr_flt()

subroutine mod_utils::degrees2meters_arr_flt ( real(spa), dimension(nsze,msze), intent(in)  LON,
real(spa), dimension(nsze,msze), intent(in)  LAT,
character(len=*), intent(in)  proj_ref,
real(spa), dimension(nsze,msze), intent(out)  X,
real(spa), dimension(nsze,msze), intent(out)  Y,
integer, intent(in)  nsze,
integer, intent(in)  msze 
)

Definition at line 402 of file mod_utils.f90.

402  implicit none
403  integer, intent(in) :: nsze, msze
404  REAL(SPA), INTENT(IN) :: LON(nsze,msze), LAT(nsze,msze)
405  character(LEN=*), INTENT(IN) :: proj_ref
406  REAL(SPA), intent(out):: X(nsze,msze), Y(nsze,msze)
407 
408  x=-1.0_spa
409  y=-1.0_spa

◆ degrees2meters_scl_dbl()

subroutine mod_utils::degrees2meters_scl_dbl ( real(dp), intent(in)  LON,
real(dp), intent(in)  LAT,
character(len=*), intent(in)  proj_ref,
real(dp), intent(out)  X,
real(dp), intent(out)  Y 
)

Definition at line 413 of file mod_utils.f90.

413 ! USE PROJ4
414  implicit none
415  REAL(DP), INTENT(IN) :: LON, LAT
416  character(LEN=*), INTENT(IN) :: proj_ref
417  REAL(DP), intent(out):: X, Y
418 
419  x=-1.0_dp
420  y=-1.0_dp

◆ degrees2meters_scl_flt()

subroutine mod_utils::degrees2meters_scl_flt ( real(spa), intent(in)  LON,
real(spa), intent(in)  LAT,
character(len=*), intent(in)  proj_ref,
real(spa), intent(out)  X,
real(spa), intent(out)  Y 
)

Definition at line 382 of file mod_utils.f90.

382  implicit none
383  REAL(SPA), INTENT(IN) :: LON, LAT
384  character(LEN=*), INTENT(IN) :: proj_ref
385  REAL(SPA), intent(out):: X, Y
386  x=-1.0_spa
387  y=-1.0_spa

◆ degrees2meters_vec_dbl()

subroutine mod_utils::degrees2meters_vec_dbl ( real(dp), dimension(nsze), intent(in)  LON,
real(dp), dimension(nsze), intent(in)  LAT,
character(len=*), intent(in)  proj_ref,
real(dp), dimension(nsze), intent(out)  X,
real(dp), dimension(nsze), intent(out)  Y,
integer, intent(in)  nsze 
)

Definition at line 424 of file mod_utils.f90.

424 ! USE PROJ4
425  implicit none
426  integer, intent(in) :: nsze
427  REAL(DP), INTENT(IN) :: LON(nsze), LAT(nsze)
428  character(LEN=*), INTENT(IN) :: proj_ref
429  REAL(DP), INTENT(OUT) :: X(nsze), Y(nsze)
430  x=-1.0_dp
431  y=-1.0_dp

◆ degrees2meters_vec_flt()

subroutine mod_utils::degrees2meters_vec_flt ( real(spa), dimension(nsze), intent(in)  LON,
real(spa), dimension(nsze), intent(in)  LAT,
character(len=*), intent(in)  proj_ref,
real(spa), dimension(nsze), intent(out)  X,
real(spa), dimension(nsze), intent(out)  Y,
integer, intent(in)  nsze 
)

Definition at line 391 of file mod_utils.f90.

391  implicit none
392  integer, intent(in) :: nsze
393  REAL(SPA), INTENT(IN) :: LON(nsze), LAT(nsze)
394  character(LEN=*), INTENT(IN) :: proj_ref
395  REAL(SPA), intent(out):: X(nsze), Y(nsze)
396 
397  x=-1.0_spa
398  y=-1.0_spa

◆ fatal_error()

subroutine mod_utils::fatal_error ( character(len=*)  ER1,
character(len=*), optional  ER2,
character(len=*), optional  ER3,
character(len=*), optional  ER4 
)

Definition at line 230 of file mod_utils.f90.

230  USE control, only: prg_name
231  implicit none
232  character(Len=*) :: ER1
233  CHARACTER(LEN=*), OPTIONAL :: ER2
234  CHARACTER(LEN=*), OPTIONAL :: ER3
235  CHARACTER(LEN=*), OPTIONAL :: ER4
236 
237  if(dbg_set(dbg_log)) then
238  write(ipt,*) "%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%"
239  write(ipt,*) trim(prg_name)//" Fatal Error!"
240  write(ipt,*) er1
241  IF(PRESENT(er2)) WRITE(ipt,*) er2
242  IF(PRESENT(er3)) WRITE(ipt,*) er3
243  IF(PRESENT(er4)) WRITE(ipt,*) er4
244  write(ipt,*) "Stopping "//trim(prg_name)
245  write(ipt,*) "%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%"
246  end if
247  call pstop
subroutine pstop
Definition: mod_utils.f90:273
character(len=80) prg_name
Definition: mod_main.f90:105
Here is the call graph for this function:

◆ find_element_containing()

integer function mod_utils::find_element_containing ( real(sp), intent(in)  xloc,
real(sp), intent(in)  yloc,
integer, intent(in), optional  GUESS 
)

Definition at line 1185 of file mod_utils.f90.

1185  !==============================================================================|
1186  ! find home element for points (x,y) |
1187  ! search nearest element to progressively further elements.
1188  !==============================================================================|
1189 
1190  !------------------------------------------------------------------------------|
1191 
1192  use all_vars
1193  implicit none
1194  INTEGER :: EID
1195 
1196  !------------------------------------------------------------------------------|
1197  real(sp),INTENT(IN) :: xloc,yloc
1198  INTEGER, INTENT(IN),OPTIONAL :: GUESS
1199 
1200  IF(PRESENT(guess)) THEN
1201  eid = find_element_containing_quick(xloc,yloc,guess)
1202  IF (eid /= 0) RETURN
1203  END IF
1204 
1205  eid = find_element_containing_robust(xloc,yloc)
1206 
Here is the call graph for this function:

◆ find_element_containing_quick()

integer function mod_utils::find_element_containing_quick ( real(sp), intent(in)  xloc,
real(sp), intent(in)  yloc,
integer, intent(in)  Guess 
)

Definition at line 1271 of file mod_utils.f90.

1271  !==============================================================================|
1272  ! determine which element a location reside in by searching neighboring
1273  ! elements.
1274  !==============================================================================|
1275 
1276  !------------------------------------------------------------------------------|
1277 
1278  use all_vars, only: nv,nbve,ntve,nv
1279  implicit none
1280  INTEGER :: EID
1281 
1282  REAL(SP), INTENT(IN) :: Xloc,Yloc
1283  INTEGER, INTENT(IN) :: GUESS
1284 
1285  integer i,j,k,iney,ncheck
1286  real(sp), dimension(3) :: xlast,ylast,xney,yney
1287 
1288  !==============================================================================|
1289  eid = 0
1290  IF (guess == 0) RETURN
1291 
1292  if(isintriangle(guess,xloc,yloc))then !!particle remains in element
1293  eid = guess
1294  else !!check neighbors
1295  outer: do j=1,3
1296  ncheck = nv(guess,j)
1297  do k=1,ntve(ncheck)
1298  iney = nbve(ncheck,k)
1299  if(isintriangle(iney,xloc,yloc))then
1300  eid = iney
1301  exit outer
1302  end if
1303  end do
1304  end do outer
1305  end if
1306 
1307  return
integer, dimension(:), allocatable, target ntve
Definition: mod_main.f90:1022
integer, dimension(:,:), allocatable, target nv
Definition: mod_main.f90:1018
integer, dimension(:,:), allocatable, target nbve
Definition: mod_main.f90:1034
Here is the call graph for this function:
Here is the caller graph for this function:

◆ find_element_containing_robust()

integer function mod_utils::find_element_containing_robust ( real(sp), intent(in)  xloc,
real(sp), intent(in)  yloc 
)

Definition at line 1215 of file mod_utils.f90.

1215  !==============================================================================|
1216  ! find home element for points (x,y) |
1217  ! search nearest element to progressively further elements.
1218  !==============================================================================|
1219 
1220  !------------------------------------------------------------------------------|
1221 
1222  use all_vars
1223  USE mod_spherical
1224  implicit none
1225  INTEGER :: EID
1226 
1227  !------------------------------------------------------------------------------|
1228  real(sp),INTENT(IN) :: xloc,yloc
1229 
1230  integer i,min_loc
1231  real(sp), dimension(1:nt,1) :: radlist
1232  real(sp), dimension(3) :: xtri,ytri
1233  real(sp) :: radlast
1234  integer :: locij(2), cnt
1235 
1236  !==============================================================================|
1237 
1238  eid = 0
1239 
1240  cnt = 0
1241 ! radlist(1:nt,1) = sqrt((xc(1:nt)-xloc)**2 + (yc(1:nt)-yloc)**2)
1242 
1243 !---------------------------------------------------------------
1244  radlist(1:nt,1) = abs(xc(1:nt)-xloc) + abs(yc(1:nt)-yloc)
1245 
1246  radlast = -1.0_sp
1247  in: do while(cnt < 50)
1248  cnt = cnt+1
1249  locij = minloc(radlist,radlist>radlast)
1250  min_loc = locij(1)
1251  if(min_loc == 0) then
1252  exit in
1253  end if
1254 
1255  if(isintriangle(min_loc,xloc,yloc))then
1256  eid = min_loc
1257  exit in
1258  end if
1259  radlast = radlist(min_loc,1)
1260  end do in
1261 
1262  return
real(sp), dimension(:), allocatable, target yc
Definition: mod_main.f90:1004
real(sp), dimension(:), allocatable, target xc
Definition: mod_main.f90:1003
Here is the call graph for this function:
Here is the caller graph for this function:

◆ fopen()

subroutine mod_utils::fopen ( integer, intent(in)  IUNIT,
character(len=*)  INSTR,
character(len=3), intent(in)  IOPT 
)

Definition at line 1577 of file mod_utils.f90.

1577 !==============================================================================|
1578 ! FOPEN Utility to open non Netcdf files
1579 !==============================================================================|
1580 
1581  IMPLICIT NONE
1582  INTEGER, INTENT(IN) :: IUNIT
1583  CHARACTER(LEN=*) :: INSTR
1584  CHARACTER(LEN=3), INTENT(IN) :: IOPT
1585  CHARACTER(LEN=11) :: FORMSTR
1586  CHARACTER(LEN=7) :: STATSTR
1587  LOGICAL CHECK,FEXIST
1588  CHARACTER(LEN=2) :: cios
1589  integer :: ios
1590 
1591  IF(iopt(1:1) == "c")THEN
1592  statstr = "old"
1593  check = .true.
1594  ELSE IF(iopt(1:1) == "o") THEN
1595  statstr = "unknown"
1596  check = .false.
1597  ELSE
1598  CALL fatal_error("FIRST LETTER IN FOPEN OPTION STRING MUST BE 'c' OR 'o'")
1599  END IF
1600 
1601  IF(iopt(2:2) == "f")THEN
1602  formstr = "formatted"
1603  ELSE IF(iopt(2:2) == "u") THEN
1604  formstr = "unformatted"
1605  ELSE
1606  CALL fatal_error("ERROR PROCESSING FOPEN ON FILE",instr,"2ND LETTER IN FOPEN OPTION STRING MUST BE 'f' OR 'u'")
1607  END IF
1608 
1609  IF(check)THEN
1610  INQUIRE(file=instr,exist=fexist)
1611  IF(.NOT. fexist) CALL fatal_error("FILE "//instr//" NOT FOUND")
1612  END IF
1613 
1614  OPEN(iunit,file=instr,status=trim(statstr),form=trim(formstr),iostat=ios)
1615 
1616 
1617  write(cios,'(i2.2)') ios
1618  if (ios == 9) then
1619  CALL fatal_error("UNABLE TO OPEN THE FILE:",&
1620  & instr, "IOSTAT ERROR#"//cios//"; suggests bad permissions ?")
1621 
1622  elseif (ios ==29) then
1623  CALL fatal_error("UNABLE TO OPEN THE FILE:",&
1624  & instr, "IOSTAT ERROR#"//cios//"; suggests bad directory path ?")
1625 
1626 
1627  elseif (ios /= 0)then
1628  Call fatal_error("UNABLE TO OPEN THE FILE:",instr,"IOSTAT ERROR# &
1629  &"//cios//"; UNNKOWN ERROR ?")
1630  END IF
1631 
1632 
1633  IF(iopt(3:3) == "r") rewind(iunit)
1634 
1635  if(dbg_set(dbg_io)) &
1636  & write(ipt,*) "Opend File: ",instr
1637 
1638 
1639 
integer ios
Definition: mod_obcs2.f90:81
Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_value()

subroutine mod_utils::get_value ( integer, intent(in)  LNUM,
integer, intent(in)  NUMCHAR,
character(len=numchar)  TEXT_LINE,
character(len=40), intent(out)  VARNAME,
character(len=7), intent(out)  VARTYPE,
logical, intent(out)  LOGVAL,
character(len=80), dimension(150), intent(out)  STRINGVAL,
real(dp), dimension(150), intent(inout)  REALVAL,
integer, dimension(150), intent(inout)  INTVAL,
integer, intent(out)  NVAL 
)

Definition at line 1677 of file mod_utils.f90.

1677 
1678 !==============================================================================!
1679  USE mod_prec
1680  IMPLICIT NONE
1681  INTEGER, INTENT(IN) :: LNUM,NUMCHAR
1682  CHARACTER(LEN=NUMCHAR) :: TEXT_LINE
1683  CHARACTER(LEN=40), INTENT(OUT) :: VARNAME
1684  CHARACTER(LEN=7), INTENT(OUT) :: VARTYPE
1685  LOGICAL, INTENT(OUT) :: LOGVAL
1686  CHARACTER(LEN=80), INTENT(OUT) :: STRINGVAL(150)
1687  REAL(DP), INTENT(INOUT) :: REALVAL(150)
1688  INTEGER, INTENT(INOUT) :: INTVAL(150)
1689  INTEGER, INTENT(OUT) :: NVAL
1690 !------------------------------------------------------------------------------!
1691  CHARACTER(LEN=NUMCHAR) :: VARVAL,TEMP,FRAG(200)
1692  CHARACTER(LEN=80) :: TSTRING
1693  CHARACTER(LEN=3) :: ERRSTRING
1694  CHARACTER(LEN=16) :: NUMCHARS
1695  INTEGER LENGTH,EQLOC,LVARVAL,DOTLOC
1696  INTEGER I,J,LOCEX,NP
1697  LOGICAL ONFRAG
1698 
1699 !==============================================================================!
1700  frag = " "
1701  numchars = "0123456789+-Ee. "
1702  vartype = "error"
1703  logval = .false.
1704  length = len_trim(text_line)
1705  WRITE(errstring,"(I3)") lnum
1706  locex = index(text_line,"!")
1707 
1708 !
1709 !-----------------------CHECK FOR BLANK LINE OR COMMENT------------------------!
1710 !
1711  IF(length == 0 .OR. locex==1)THEN
1712  vartype = "no data"
1713  varname = "no data"
1714  RETURN
1715  END IF
1716 
1717 !
1718 !-----------------------CHANGE COMMAS TO BLANKS--------------------------------!
1719 !
1720  DO i=1,length
1721  IF(text_line(i:i) == ",") text_line(i:i) = " "
1722  END DO
1723 !
1724 !-----------------------REMOVING TRAILING COMMENTS-----------------------------!
1725 !
1726  IF(locex /= 0)THEN
1727  temp = text_line(1:locex-1)
1728  text_line = temp
1729  END IF
1730 !
1731 !--------------------ENSURE "=" EXISTS AND DETERMINE LOCATION------------------!
1732 !
1733  eqloc = index(text_line,"=")
1734  IF(eqloc == 0) then
1735  CALL warning(&
1736  &'Could not find correct variable name in the datafile header',&
1737  &'Header comment lines must start with "!", Data lines must contain "="',&
1738  &'DATA LINE '//errstring//': This often occurs if the header variable is missing!')
1739  vartype = "no data"
1740  varname = "no data"
1741  RETURN
1742  END IF
1743 
1744 !
1745 !--------------------SPLIT OFF VARNAME AND VARVAL STRINGS----------------------!
1746 !
1747  varname = text_line(1:eqloc-1)
1748  varval = adjustl(text_line(eqloc+1:length))
1749  lvarval = len_trim(varval)
1750  IF(lvarval == 0) then
1751  CALL warning('IN DATA PARAMETER FILE', &
1752  & 'VARIABLE'//varname//'; LINE'//errstring//' HAS NO ASSOCIATED VALUE')
1753  vartype = "no data"
1754  varname = "no data"
1755  RETURN
1756  END IF
1757 
1758 !-----------------DETERMINE TYPE OF VARVAL-------------------------------------!
1759 !
1760 
1761 !
1762 ! CHECK FOR LOGICAL
1763 !
1764  IF((varval(1:1) == "T" .OR. varval(1:1) == "F") .AND. lvarval == 1)THEN
1765  vartype = "logical"
1766  IF(varval(1:1) == "T") logval = .true.
1767  RETURN
1768  END IF
1769 
1770 !
1771 ! CHECK IF IT IS A STRING (CONTAINS CHARACTERS OTHER THAN 0-9,+,-,e,E,.)
1772 !
1773  DO i=1,lvarval
1774  IF(index(numchars,varval(i:i)) == 0) vartype = "string"
1775  END DO
1776 
1777 !
1778 ! PROCESS STRING (MAY BE MULTIPLE)
1779 !
1780  IF(vartype == "string") THEN
1781  tstring = varval
1782  stringval(1) = tstring
1783  nval = 1
1784  onfrag = .true.
1785  DO i=1,lvarval
1786  IF(varval(i:i) /= " ")THEN
1787  frag(nval) = trim(frag(nval))//varval(i:i)
1788  onfrag = .true.
1789  ELSE
1790  IF(onfrag) nval = nval + 1
1791  onfrag = .false.
1792  END IF
1793  END DO
1794  DO i=1,nval
1795  stringval(i+1) = trim(frag(i))
1796  END DO
1797  RETURN
1798  END IF
1799 
1800 !
1801 ! CHECK IF IT IS A FLOAT
1802 !
1803 
1804  dotloc = index(varval,".")
1805  IF(dotloc /= 0) THEN
1806  vartype = "float"
1807  ELSE
1808  vartype = "integer"
1809  END IF
1810 !
1811 !-----------------FRAGMENT INTO STRINGS FOR MULTIPLE VALUES---------------------!
1812 !
1813  np = 1
1814  onfrag = .true.
1815  DO i=1,lvarval
1816  IF(varval(i:i) /= " ")THEN
1817  frag(np) = trim(frag(np))//varval(i:i)
1818  onfrag = .true.
1819  ELSE
1820  IF(onfrag) np = np + 1
1821  onfrag = .false.
1822  END IF
1823  END DO
1824 !
1825 !-----------------EXTRACT NUMBER(S) FROM CHARACTER STRINGS----------------------!
1826 !
1827 
1828  nval = np
1829  DO i=1,np
1830  temp = trim(frag(i))
1831  IF(vartype == "float") THEN
1832  READ(temp,*)realval(i)
1833  ELSE
1834  READ(temp,*)intval(i)
1835  END IF
1836  END DO
1837 
Here is the call graph for this function:
Here is the caller graph for this function:

◆ grid_neighbor_index()

subroutine mod_utils::grid_neighbor_index ( integer, dimension(:), pointer  FOUND,
integer, dimension(:,:), pointer  IDEX,
integer, dimension(:), pointer  CNT,
integer, dimension(:), pointer  ORDER 
)

Definition at line 1423 of file mod_utils.f90.

1423  USE all_vars
1424  IMPLICIT NONE
1425  INTEGER, POINTER :: FOUND(:)
1426  INTEGER, POINTER :: IDEX(:,:)
1427  INTEGER, POINTER :: CNT(:)
1428  INTEGER, POINTER :: ORDER(:)
1429 
1430  INTEGER :: i,j
1431  INTEGER :: LOOP
1432  INTEGER :: ORD
1433  ! DO NO ALLOCATION HERE - ASSUME ALL VARIALBE ARE ALLOCATED
1434 
1435  ! LOOKING FOR NEIGHBORING NODES
1436  IF (ubound(found,1) == mt) THEN
1437 
1438  loop = 0
1439  ord = 0
1440  DO WHILE(any(found==-1))
1441 
1442  loop = loop +1
1443  IF(loop>m) CALL fatal_error&
1444  &("LOOP COUNT EXCEEDED IN GRID_NEIGHBOR_INDEX")
1445 
1446  DO i=1,m
1447 
1448  IF(found(i) == -1) THEN
1449 
1450  ! look to see if this node has neighbors which are set
1451  DO j=1,ntsn(i)
1452  IF(found(nbsn(i,j))>-1.and.found(nbsn(i,j))<loop )THEN
1453  ! INCREASE THE COUNT FOR AVERAGE
1454  cnt(i)= cnt(i) +1
1455  ! RECORD THE INDEX
1456  idex(i,cnt(i)) = nbsn(i,j)
1457  ! RECORD WHICH LOOP FOUND IT
1458  found(i) = loop
1459 
1460  END IF
1461 
1462 
1463  END DO
1464 
1465  ! RECORD THE ORDER IN WHICH THIS NODE WAS FOUND
1466  IF (cnt(i) >0) THEN
1467  ord = ord +1
1468  order(ord) = i
1469  END IF
1470 
1471  END IF
1472 
1473  END DO
1474 
1475  END DO
1476 
1477  ! LOOKING FOR NEIGHBORING ELEMENTS
1478  ELSEIF (ubound(found,1) == nt) THEN
1479 
1480  loop = 0
1481  ord = 0
1482  DO WHILE(any(found==-1))
1483 
1484  loop = loop +1
1485  IF(loop>n) CALL fatal_error&
1486  &("LOOP COUNT EXCEEDED IN GRID_NEIGHBOR_INDEX")
1487 
1488  DO i=1,n
1489 
1490  IF(found(i) == -1) THEN
1491 
1492  ! look to see if this element has neighbors which are found
1493  DO j=1,3
1494  IF(found(nbe(i,j))>-1 .and.found(nbe(i,j))<loop )THEN
1495  ! INCREASE THE COUNT FOR AVERAGE
1496  cnt(i)= cnt(i) +1
1497  ! RECORD THE INDEX
1498  idex(i,cnt(i)) = nbe(i,j)
1499  ! RECORD WHICH LOOP FOUND IT
1500  found(i) = loop
1501 
1502  END IF
1503 
1504 
1505  END DO
1506 
1507  ! RECORD THE ORDER IN WHICH THIS NODE WAS FOUND
1508  IF (cnt(i) >0) THEN
1509  ord = ord +1
1510  order(ord) = i
1511  END IF
1512 
1513  END IF
1514 
1515  END DO
1516 
1517  END DO
1518 
1519 
1520  ELSE
1521 
1522  CALL fatal_error("PASSED INVALID SIZE TO GRID_NEIGHBOR_INDEX ???")
1523 
1524  END IF
integer, dimension(:), allocatable, target ntsn
Definition: mod_main.f90:1023
integer, dimension(:,:), allocatable, target nbe
Definition: mod_main.f90:1020
integer, dimension(:,:), allocatable, target nbsn
Definition: mod_main.f90:1030
Here is the call graph for this function:
Here is the caller graph for this function:

◆ have_proj()

logical function mod_utils::have_proj ( character(len=*), intent(in)  proj_ref)

Definition at line 376 of file mod_utils.f90.

376  LOGICAL :: RES
377  character(LEN=*), INTENT(IN) :: proj_ref
378  res=.false.
Here is the caller graph for this function:

◆ initialize_control()

subroutine mod_utils::initialize_control ( character(len=*), intent(in)  NAME)

Definition at line 141 of file mod_utils.f90.

141  USE lims
142  USE control
143  IMPLICIT NONE
144  CHARACTER(LEN=*), INTENT(IN) :: NAME
145 
146 !==============================================================================!
147 ! FVCOM VERSION !
148 !==============================================================================!
149 
150  fvcom_version = 'FVCOM_3.0'
151  fvcom_website = 'http://fvcom.smast.umassd.edu, http://codfish.smast.umassd.edu'
152  institution = 'School for Marine Science and Technology'
153 
154  ! Set the IO UNIT value to screen output for now:
155  ipt = 6
156 
157 
158 !==============================================================================!
159 ! SETUP PARALLEL ENVIRONMENT !
160 !==============================================================================!
161  ! DEFAULT SETTINGS
162  serial = .true.
163  par = .false.
164  msr = .true.
165  in_mpi_io_loop = .false.
166  use_mpi_io_mode = .false.
167  myid = 1
168  msrid = 1
169  nprocs = 1
170  prg_name = name
171  mx_nbr_elem = 0
173  ioprocid =-1
174 
175  zerotime%MUSOD = 0
176  zerotime%MJD = 0
177  mpi_fvcom_group = -1
logical serial
Definition: mod_main.f90:100
logical msr
Definition: mod_main.f90:101
type(time) zerotime
Definition: mod_main.f90:830
integer myid
Definition: mod_main.f90:67
logical use_mpi_io_mode
Definition: mod_main.f90:188
logical par
Definition: mod_main.f90:102
integer, target nprocs
Definition: mod_main.f90:72
integer, pointer nprocs_total
Definition: mod_main.f90:71
character(len=80) institution
Definition: mod_main.f90:112
integer mpi_fvcom_group
Definition: mod_main.f90:107
integer ioprocid
Definition: mod_main.f90:69
character(len=80) fvcom_website
Definition: mod_main.f90:113
character(len=80) fvcom_version
Definition: mod_main.f90:111
integer mx_nbr_elem
Definition: mod_main.f90:79
integer msrid
Definition: mod_main.f90:68
integer ipt
Definition: mod_main.f90:922
logical in_mpi_io_loop
Definition: mod_main.f90:104
character(len=80) prg_name
Definition: mod_main.f90:105
Here is the caller graph for this function:

◆ interp_anodal_2d_dbl()

real(dp) function mod_utils::interp_anodal_2d_dbl ( real(dp), intent(in)  xloc,
real(dp), intent(in)  yloc,
integer, intent(in)  i,
real(dp), dimension(:), intent(in), allocatable, target  Field 
)

Definition at line 535 of file mod_utils.f90.

535  IMPLICIT NONE
536  REAL(DP) :: FPT
537  INTEGER, INTENT(IN) :: i ! Cell Number
538  REAL(DP), INTENT(IN):: xloc,yloc
539  REAL(DP), ALLOCATABLE, TARGET, INTENT(IN) ::FIELD(:)
540  REAL(DP), POINTER :: PFIELD(:)
541 
542  pfield => field
543  fpt = interp_pnodal_2d_dbl(xloc,yloc,i,pfield)
544 
Here is the call graph for this function:

◆ interp_anodal_2d_flt()

real(spa) function mod_utils::interp_anodal_2d_flt ( real(spa), intent(in)  xloc,
real(spa), intent(in)  yloc,
integer, intent(in)  i,
real(spa), dimension(:), intent(in), allocatable, target  Field 
)

Definition at line 522 of file mod_utils.f90.

522  IMPLICIT NONE
523  REAL(SPA) :: FPT
524  INTEGER, INTENT(IN) :: i ! Cell Number
525  REAL(SPA), INTENT(IN):: xloc,yloc
526  REAL(SPA), ALLOCATABLE, TARGET, INTENT(IN) ::FIELD(:)
527  REAL(SPA), POINTER :: PFIELD(:)
528 
529  pfield => field
530  fpt = interp_pnodal_2d_flt(xloc,yloc,i,pfield)
531 
Here is the call graph for this function:

◆ interp_anodal_3d_dbl()

real(dp) function mod_utils::interp_anodal_3d_dbl ( real(dp), intent(in)  xloc,
real(dp), intent(in)  yloc,
real(dp), intent(in)  sigloc,
integer, intent(in)  lvls,
integer, intent(in)  i,
real(dp), dimension(:,:), intent(in), allocatable, target  Field 
)

Definition at line 639 of file mod_utils.f90.

639  IMPLICIT NONE
640  REAL(DP) :: FPT
641  INTEGER, INTENT(IN) :: i,lvls ! Cell Number, Number of Levels in data (kb or kbm1)
642  REAL(DP), INTENT(IN):: xloc,yloc,sigloc
643  REAL(DP), ALLOCATABLE, TARGET, INTENT(IN) ::FIELD(:,:)
644  REAL(DP), POINTER :: PFIELD(:,:)
645 
646  pfield => field
647  fpt = interp_pnodal_3d_dbl(xloc,yloc,sigloc,lvls,i,pfield)
648 
Here is the call graph for this function:

◆ interp_anodal_3d_flt()

real(spa) function mod_utils::interp_anodal_3d_flt ( real(spa), intent(in)  xloc,
real(spa), intent(in)  yloc,
real(spa), intent(in)  sigloc,
integer, intent(in)  lvls,
integer, intent(in)  i,
real(spa), dimension(:,:), intent(in), allocatable, target  Field 
)

Definition at line 626 of file mod_utils.f90.

626  IMPLICIT NONE
627  REAL(SPA) :: FPT
628  INTEGER, INTENT(IN) :: i,lvls ! Cell Number, Number of Levels in data (kb or kbm1)
629  REAL(SPA), INTENT(IN):: xloc,yloc,sigloc
630  REAL(SPA), ALLOCATABLE, TARGET, INTENT(IN) ::FIELD(:,:)
631  REAL(SPA), POINTER :: PFIELD(:,:)
632 
633  pfield => field
634  fpt = interp_pnodal_3d_flt(xloc,yloc,sigloc,lvls,i,pfield)
635 
Here is the call graph for this function:

◆ interp_azonal_2d_dbl()

real(dp) function mod_utils::interp_azonal_2d_dbl ( real(dp), intent(in)  xloc,
real(dp), intent(in)  yloc,
integer, intent(in)  i,
real(dp), dimension(:), intent(in), allocatable, target  Field 
)

Definition at line 871 of file mod_utils.f90.

871  IMPLICIT NONE
872  REAL(DP) :: FPT
873  INTEGER, INTENT(IN) :: i ! Cell Number
874  REAL(DP), INTENT(IN):: xloc,yloc
875  REAL(DP), ALLOCATABLE, TARGET, INTENT(IN) ::FIELD(:)
876  REAL(DP), POINTER :: PFIELD(:)
877 
878  pfield => field
879  fpt = interp_pzonal_2d_dbl(xloc,yloc,i,pfield)
880 
Here is the call graph for this function:

◆ interp_azonal_2d_flt()

real(spa) function mod_utils::interp_azonal_2d_flt ( real(spa), intent(in)  xloc,
real(spa), intent(in)  yloc,
integer, intent(in)  i,
real(spa), dimension(:), intent(in), allocatable, target  Field 
)

Definition at line 858 of file mod_utils.f90.

858  IMPLICIT NONE
859  REAL(SPA) :: FPT
860  INTEGER, INTENT(IN) :: i ! Cell Number
861  REAL(SPA), INTENT(IN):: xloc,yloc
862  REAL(SPA), ALLOCATABLE, TARGET, INTENT(IN) ::FIELD(:)
863  REAL(SPA), POINTER :: PFIELD(:)
864 
865  pfield => field
866  fpt = interp_pzonal_2d_flt(xloc,yloc,i,pfield)
867 
Here is the call graph for this function:

◆ interp_azonal_3d_dbl()

real(dp) function mod_utils::interp_azonal_3d_dbl ( real(dp), intent(in)  xloc,
real(dp), intent(in)  yloc,
real(dp), intent(in)  sigloc,
integer, intent(in)  lvls,
integer, intent(in)  i,
real(dp), dimension(:,:), intent(in), allocatable, target  Field 
)

Definition at line 973 of file mod_utils.f90.

973  IMPLICIT NONE
974  REAL(DP) :: FPT
975  INTEGER, INTENT(IN) :: i, lvls ! Cell Number, Number of levels(kb,or kbm1)
976  REAL(DP), INTENT(IN):: xloc,yloc,sigloc
977  REAL(DP), ALLOCATABLE, TARGET, INTENT(IN) ::FIELD(:,:)
978  REAL(DP), POINTER :: PFIELD(:,:)
979 
980  pfield => field
981  fpt = interp_pzonal_3d_dbl(xloc,yloc,sigloc,lvls,i,pfield)
982 
Here is the call graph for this function:

◆ interp_azonal_3d_flt()

real(spa) function mod_utils::interp_azonal_3d_flt ( real(spa), intent(in)  xloc,
real(spa), intent(in)  yloc,
real(spa), intent(in)  sigloc,
integer, intent(in)  lvls,
integer, intent(in)  i,
real(spa), dimension(:,:), intent(in), allocatable, target  Field 
)

Definition at line 960 of file mod_utils.f90.

960  IMPLICIT NONE
961  REAL(SPA) :: FPT
962  INTEGER, INTENT(IN) :: i, lvls ! Cell Number, Number of levels(kb,or kbm1)
963  REAL(SPA), INTENT(IN):: xloc,yloc,sigloc
964  REAL(SPA), ALLOCATABLE, TARGET, INTENT(IN) ::FIELD(:,:)
965  REAL(SPA), POINTER :: PFIELD(:,:)
966 
967  pfield => field
968  fpt = interp_pzonal_3d_flt(xloc,yloc,sigloc,lvls,i,pfield)
969 
Here is the call graph for this function:

◆ interp_pnodal_2d_dbl()

real(dp) function mod_utils::interp_pnodal_2d_dbl ( real(dp), intent(in)  xloc,
real(dp), intent(in)  yloc,
integer, intent(in)  i,
real(dp), dimension(:), intent(in), pointer  Field 
)

Definition at line 585 of file mod_utils.f90.

585  USE all_vars, only : aw0,awx,awy,nv,xc,yc
586  IMPLICIT NONE
587  REAL(DP) :: FPT
588  INTEGER, INTENT(IN) :: i ! Cell Number
589  REAL(DP), INTENT(IN):: xloc,yloc
590  REAL(DP), POINTER, INTENT(IN) :: FIELD(:)
591 
592  REAL(DP):: X0c, Y0c,F0,Fx,Fy
593  INTEGER :: n1,n2,n3
594  REAL(DP) :: dx_sph
595  IF(dbg_set(dbg_sbr)) WRITE(ipt,*) "START: INTERP_PNODAL_2D"
596 
597  !element location (i) and surrounding nodes (n1,n2,n3)
598  n1 = nv(i,1)
599  n2 = nv(i,2)
600  n3 = nv(i,3)
601 
602 !offset from element center
603 !---------------------------------------------------------------
604  !offset from element center
605  x0c = xloc - xc(i)
606  y0c = yloc - yc(i)
607 !---------------------------------------------------------------
608 !---------------------------------------------------------------
609 
610  !linear interpolation of Field
611  f0 = aw0(i,1)*field(n1)+aw0(i,2)*field(n2)+aw0(i,3)*field(n3)
612  fx = awx(i,1)*field(n1)+awx(i,2)*field(n2)+awx(i,3)*field(n3)
613  fy = awy(i,1)*field(n1)+awy(i,2)*field(n2)+awy(i,3)*field(n3)
614  fpt = f0 + fx*x0c + fy*y0c
615 
616 
617  IF(dbg_set(dbg_sbr)) WRITE(ipt,*)"END: INTERP_PNODAL_2D"
real(sp), dimension(:), allocatable, target yc
Definition: mod_main.f90:1004
real(sp), dimension(:,:), allocatable, target awx
Definition: mod_main.f90:1333
real(sp), dimension(:,:), allocatable, target aw0
Definition: mod_main.f90:1335
real(sp), dimension(:,:), allocatable, target awy
Definition: mod_main.f90:1334
integer, dimension(:,:), allocatable, target nv
Definition: mod_main.f90:1018
real(sp), dimension(:), allocatable, target xc
Definition: mod_main.f90:1003
integer ipt
Definition: mod_main.f90:922
Here is the call graph for this function:
Here is the caller graph for this function:

◆ interp_pnodal_2d_flt()

real(spa) function mod_utils::interp_pnodal_2d_flt ( real(spa), intent(in)  xloc,
real(spa), intent(in)  yloc,
integer, intent(in)  i,
real(spa), dimension(:), intent(in), pointer  Field 
)

Definition at line 548 of file mod_utils.f90.

548  USE all_vars, only : aw0,awx,awy,nv,xc,yc
549  IMPLICIT NONE
550  REAL(SPA) :: FPT
551  INTEGER, INTENT(IN) :: i ! Cell Number
552  REAL(SPA), INTENT(IN):: xloc,yloc
553  REAL(SPA), POINTER, INTENT(IN) :: FIELD(:)
554 
555  REAL(SPA):: X0c, Y0c,F0,Fx,Fy
556  INTEGER :: n1,n2,n3
557  REAL(SPA) :: dx_sph
558  IF(dbg_set(dbg_sbr)) WRITE(ipt,*) "START: INTERP_PNODAL_2D"
559 
560  !element location (i) and surrounding nodes (n1,n2,n3)
561  n1 = nv(i,1)
562  n2 = nv(i,2)
563  n3 = nv(i,3)
564 
565 !offset from element center
566 !---------------------------------------------------------------
567 
568  x0c = xloc - xc(i)
569  y0c = yloc - yc(i)
570 !---------------------------------------------------------------
571 !---------------------------------------------------------------
572 
573  !linear interpolation of Field
574  f0 = aw0(i,1)*field(n1)+aw0(i,2)*field(n2)+aw0(i,3)*field(n3)
575  fx = awx(i,1)*field(n1)+awx(i,2)*field(n2)+awx(i,3)*field(n3)
576  fy = awy(i,1)*field(n1)+awy(i,2)*field(n2)+awy(i,3)*field(n3)
577  fpt = f0 + fx*x0c + fy*y0c
578 
579 
580  IF(dbg_set(dbg_sbr)) WRITE(ipt,*)"END: INTERP_PNODAL_2D"
real(sp), dimension(:), allocatable, target yc
Definition: mod_main.f90:1004
real(sp), dimension(:,:), allocatable, target awx
Definition: mod_main.f90:1333
real(sp), dimension(:,:), allocatable, target aw0
Definition: mod_main.f90:1335
real(sp), dimension(:,:), allocatable, target awy
Definition: mod_main.f90:1334
integer, dimension(:,:), allocatable, target nv
Definition: mod_main.f90:1018
real(sp), dimension(:), allocatable, target xc
Definition: mod_main.f90:1003
integer ipt
Definition: mod_main.f90:922
Here is the call graph for this function:
Here is the caller graph for this function:

◆ interp_pnodal_3d_dbl()

real(dp) function mod_utils::interp_pnodal_3d_dbl ( real(dp), intent(in)  xloc,
real(dp), intent(in)  yloc,
real(dp), intent(in)  sigloc,
integer, intent(in)  lvls,
integer, intent(in)  i,
real(dp), dimension(:,:), intent(in), pointer  Field 
)

Definition at line 752 of file mod_utils.f90.

752  USE all_vars, only : aw0,awx,awy,nv,xc,yc,kbm2,kbm1,kb,z1,zz1,dz1,dzz1
753  IMPLICIT NONE
754  REAL(DP) :: FPT
755  INTEGER, INTENT(IN) :: i, lvls ! Cell Number, Number of Levels in data (kb or kbm1)
756  REAL(DP), INTENT(IN):: xloc,yloc,sigloc
757  REAL(DP), POINTER, INTENT(IN) :: FIELD(:,:)
758 
759  REAL(DP):: X0c,Y0c,F0,Fx,Fy,F_LOWER,F_UPPER, alpha,dsig
760  INTEGER :: n1,n2,n3,k1,k2,k
761  REAL(DP) :: dx_sph
762 
763  IF(dbg_set(dbg_sbr)) WRITE(ipt,*) "START: INTERP_PNODAL_3D"
764 
765  !element location (i) and surrounding nodes (n1,n2,n3)
766  n1 = nv(i,1)
767  n2 = nv(i,2)
768  n3 = nv(i,3)
769 
770 !offset from element center
771 !---------------------------------------------------------------
772 
773  x0c = xloc - xc(i)
774  y0c = yloc - yc(i)
775 !---------------------------------------------------------------
776 !---------------------------------------------------------------
777 
778 
779  ! Determine the layer in which the point resides
780  IF(lvls == kbm1) THEN
781 
782  !top
783  if(sigloc >= zz1(i,1))then
784  k1 = 1
785  k2 = 1
786  alpha = -1
787  elseif(sigloc > zz1(i,kbm1)) then!intermediate
788  do k=1,kbm2
789  if(sigloc < zz1(i,k) .and. sigloc >= zz1(i,k+1) )then
790  k1 = k
791  k2 = k+1
792  alpha = (zz1(i,k)-sigloc)/dzz1(i,k)
793  exit
794  endif
795  end do
796  else
797  ! TOP
798  k1 = kbm1
799  k2 = kbm1
800  alpha = -1
801  endif
802  ELSE IF(lvls == kb) THEN
803 
804  !top
805  if(sigloc >= z1(i,1))then
806  k1 = 1
807  k2 = 1
808  alpha = -1
809  elseif(sigloc > z1(i,kb)) then !intermediate
810  do k=1,kbm1
811  if(sigloc < z1(i,k) .and. sigloc >= z1(i,k+1) )then
812  k1 = k
813  k2 = k+1
814  alpha = (z1(i,k)-sigloc)/dz1(i,k)
815  endif
816  end do
817  else
818  !bottom
819  k1 = kbm1
820  k2 = kbm1
821  alpha = -1
822  endif
823 
824  ELSE
825  CALL fatal_error("INTERP_PNODAL_3D: Invalid number of levels passed",&
826  & "(Must be equal to either KB or KBM1")
827  END IF
828 
829  !linear interpolation of Field
830  f0 = aw0(i,1)*field(n1,k1)+aw0(i,2)*field(n2,k1)+aw0(i,3)*field(n3,k1)
831  fx = awx(i,1)*field(n1,k1)+awx(i,2)*field(n2,k1)+awx(i,3)*field(n3,k1)
832  fy = awy(i,1)*field(n1,k1)+awy(i,2)*field(n2,k1)+awy(i,3)*field(n3,k1)
833  f_upper = f0 + fx*x0c + fy*y0c
834 
835  IF(k1 == k2) THEN
836  fpt = f_upper
837 
838  ELSE
839 
840  f0 = aw0(i,1)*field(n1,k2)+aw0(i,2)*field(n2,k2)+aw0(i,3)*field(n3,k2)
841  fx = awx(i,1)*field(n1,k2)+awx(i,2)*field(n2,k2)+awx(i,3)*field(n3,k2)
842  fy = awy(i,1)*field(n1,k2)+awy(i,2)*field(n2,k2)+awy(i,3)*field(n3,k2)
843  f_lower = f0 + fx*x0c + fy*y0c
844 
845  fpt = (alpha)*f_lower + (1.0_dp-alpha)*f_upper
846  END IF
847 
848 
849  IF(dbg_set(dbg_sbr)) WRITE(ipt,*)"END: INTERP_PNODAL_3D"
real(sp), dimension(:), allocatable, target yc
Definition: mod_main.f90:1004
real(sp), dimension(:,:), allocatable, target dzz1
Definition: mod_main.f90:1097
real(sp), dimension(:,:), allocatable, target awx
Definition: mod_main.f90:1333
real(sp), dimension(:,:), allocatable, target aw0
Definition: mod_main.f90:1335
real(sp), dimension(:,:), allocatable, target awy
Definition: mod_main.f90:1334
integer kb
Definition: mod_main.f90:64
integer kbm2
Definition: mod_main.f90:66
integer, dimension(:,:), allocatable, target nv
Definition: mod_main.f90:1018
real(sp), dimension(:,:), allocatable, target zz1
Definition: mod_main.f90:1095
real(sp), dimension(:), allocatable, target xc
Definition: mod_main.f90:1003
real(sp), dimension(:,:), allocatable, target dz1
Definition: mod_main.f90:1096
real(sp), dimension(:,:), allocatable, target z1
Definition: mod_main.f90:1094
integer ipt
Definition: mod_main.f90:922
integer kbm1
Definition: mod_main.f90:65
Here is the call graph for this function:
Here is the caller graph for this function:

◆ interp_pnodal_3d_flt()

real(spa) function mod_utils::interp_pnodal_3d_flt ( real(spa), intent(in)  xloc,
real(spa), intent(in)  yloc,
real(spa), intent(in)  sigloc,
integer, intent(in)  lvls,
integer, intent(in)  i,
real(spa), dimension(:,:), intent(in), pointer  Field 
)

Definition at line 652 of file mod_utils.f90.

652  USE all_vars, only : aw0,awx,awy,nv,xc,yc,kbm2,kbm1,kb,z1,zz1,dz1,dzz1
653  IMPLICIT NONE
654  REAL(SPA) :: FPT
655  INTEGER, INTENT(IN) :: i, lvls ! Cell Number, Number of Levels in data (kb or kbm1)
656  REAL(SPA), INTENT(IN):: xloc,yloc,sigloc
657  REAL(SPA), POINTER, INTENT(IN) :: FIELD(:,:)
658 
659  REAL(SPA):: X0c,Y0c,F0,Fx,Fy,F_LOWER,F_UPPER, alpha,dsig
660  INTEGER :: n1,n2,n3,k1,k2,k
661  REAL(SPA) :: dx_sph
662 
663  IF(dbg_set(dbg_sbr)) WRITE(ipt,*) "START: INTERP_PNODAL_3D"
664 
665  !element location (i) and surrounding nodes (n1,n2,n3)
666  n1 = nv(i,1)
667  n2 = nv(i,2)
668  n3 = nv(i,3)
669 
670 !offset from element center
671 !---------------------------------------------------------------
672 
673  x0c = xloc - xc(i)
674  y0c = yloc - yc(i)
675 !---------------------------------------------------------------
676 !---------------------------------------------------------------
677 
678  ! Determine the layer in which the point resides
679  IF(lvls == kbm1) THEN
680 
681  !top
682  if(sigloc >= zz1(i,1))then
683  k1 = 1
684  k2 = 1
685  alpha = -1
686  elseif(sigloc > zz1(i,kbm1)) then!intermediate
687  do k=1,kbm2
688  if(sigloc < zz1(i,k) .and. sigloc >= zz1(i,k+1) )then
689  k1 = k
690  k2 = k+1
691  alpha = (zz1(i,k)-sigloc)/dzz1(i,k)
692  exit
693  endif
694  end do
695  else
696  ! TOP
697  k1 = kbm1
698  k2 = kbm1
699  alpha = -1
700  endif
701  ELSE IF(lvls == kb) THEN
702 
703  !top
704  if(sigloc >= z1(i,1))then
705  k1 = 1
706  k2 = 1
707  alpha = -1
708  elseif(sigloc > z1(i,kb)) then !intermediate
709  do k=1,kbm1
710  if(sigloc < z1(i,k) .and. sigloc >= z1(i,k+1) )then
711  k1 = k
712  k2 = k+1
713  alpha = (z1(i,k)-sigloc)/dz1(i,k)
714  endif
715  end do
716  else
717  !bottom
718  k1 = kbm1
719  k2 = kbm1
720  alpha = -1
721  endif
722 
723  ELSE
724  CALL fatal_error("INTERP_PNODAL_3D: Invalid number of levels passed",&
725  & "(Must be equal to either KB or KBM1")
726  END IF
727 
728  !linear interpolation of Field
729  f0 = aw0(i,1)*field(n1,k1)+aw0(i,2)*field(n2,k1)+aw0(i,3)*field(n3,k1)
730  fx = awx(i,1)*field(n1,k1)+awx(i,2)*field(n2,k1)+awx(i,3)*field(n3,k1)
731  fy = awy(i,1)*field(n1,k1)+awy(i,2)*field(n2,k1)+awy(i,3)*field(n3,k1)
732  f_upper = f0 + fx*x0c + fy*y0c
733 
734  IF(k1 == k2) THEN
735  fpt = f_upper
736 
737  ELSE
738 
739  f0 = aw0(i,1)*field(n1,k2)+aw0(i,2)*field(n2,k2)+aw0(i,3)*field(n3,k2)
740  fx = awx(i,1)*field(n1,k2)+awx(i,2)*field(n2,k2)+awx(i,3)*field(n3,k2)
741  fy = awy(i,1)*field(n1,k2)+awy(i,2)*field(n2,k2)+awy(i,3)*field(n3,k2)
742  f_lower = f0 + fx*x0c + fy*y0c
743 
744  fpt = (alpha)*f_lower + (1.0_spa-alpha)*f_upper
745  END IF
746 
747 
748  IF(dbg_set(dbg_sbr)) WRITE(ipt,*)"END: INTERP_PNODAL_3D"
real(sp), dimension(:), allocatable, target yc
Definition: mod_main.f90:1004
real(sp), dimension(:,:), allocatable, target dzz1
Definition: mod_main.f90:1097
real(sp), dimension(:,:), allocatable, target awx
Definition: mod_main.f90:1333
real(sp), dimension(:,:), allocatable, target aw0
Definition: mod_main.f90:1335
real(sp), dimension(:,:), allocatable, target awy
Definition: mod_main.f90:1334
integer kb
Definition: mod_main.f90:64
integer kbm2
Definition: mod_main.f90:66
integer, dimension(:,:), allocatable, target nv
Definition: mod_main.f90:1018
real(sp), dimension(:,:), allocatable, target zz1
Definition: mod_main.f90:1095
real(sp), dimension(:), allocatable, target xc
Definition: mod_main.f90:1003
real(sp), dimension(:,:), allocatable, target dz1
Definition: mod_main.f90:1096
real(sp), dimension(:,:), allocatable, target z1
Definition: mod_main.f90:1094
integer ipt
Definition: mod_main.f90:922
integer kbm1
Definition: mod_main.f90:65
Here is the call graph for this function:
Here is the caller graph for this function:

◆ interp_pzonal_2d_dbl()

real(dp) function mod_utils::interp_pzonal_2d_dbl ( real(dp), intent(in)  xloc,
real(dp), intent(in)  yloc,
integer, intent(in)  i,
real(dp), dimension(:), intent(in), pointer  Field 
)

Definition at line 920 of file mod_utils.f90.

920  USE all_vars, only : a1u,a2u,xc,yc,nbe
921  IMPLICIT NONE
922  REAL(DP) :: FPT
923  INTEGER, INTENT(IN) :: i ! Cell Number
924  REAL(DP), INTENT(IN):: xloc,yloc
925  REAL(DP), POINTER, INTENT(IN) :: FIELD(:)
926 
927  REAL(DP):: X0c, Y0c,F0,Fx,Fy
928  INTEGER :: e1,e2,e3
929  REAL(DP) :: dx_sph
930  IF(dbg_set(dbg_sbr)) WRITE(ipt,*) "START: INTERP_PZONAL_2D"
931 
932 !offset from element center
933 !---------------------------------------------------------------
934  !offset from element center
935  x0c = xloc - xc(i)
936  y0c = yloc - yc(i)
937 !---------------------------------------------------------------
938 !---------------------------------------------------------------
939 
940  !Surrounding Element IDs
941  e1 = nbe(i,1)
942  e2 = nbe(i,2)
943  e3 = nbe(i,3)
944 
945  !interpolate Field to the location
946  fx = a1u(i,1)*field(i)+a1u(i,2)*field(e1)+a1u(i,3)*field(e2)+a1u(i,4)*field(e3)
947  fy = a2u(i,1)*field(i)+a2u(i,2)*field(e1)+a2u(i,3)*field(e2)+a2u(i,4)*field(e3)
948  fpt = field(i) + fx*x0c + fy*y0c
949 
950  IF(dbg_set(dbg_sbr)) WRITE(ipt,*) "START: INTERP_PZONAL_2D"
951 
real(sp), dimension(:), allocatable, target yc
Definition: mod_main.f90:1004
real(sp), dimension(:,:), allocatable, target a1u
Definition: mod_main.f90:1331
integer, dimension(:,:), allocatable, target nbe
Definition: mod_main.f90:1020
real(sp), dimension(:,:), allocatable, target a2u
Definition: mod_main.f90:1332
real(sp), dimension(:), allocatable, target xc
Definition: mod_main.f90:1003
integer ipt
Definition: mod_main.f90:922
Here is the call graph for this function:
Here is the caller graph for this function:

◆ interp_pzonal_2d_flt()

real(spa) function mod_utils::interp_pzonal_2d_flt ( real(spa), intent(in)  xloc,
real(spa), intent(in)  yloc,
integer, intent(in)  i,
real(spa), dimension(:), intent(in), pointer  Field 
)

Definition at line 886 of file mod_utils.f90.

886  USE all_vars, only : a1u,a2u,xc,yc,nbe
887  IMPLICIT NONE
888  REAL(SPA) :: FPT
889  INTEGER, INTENT(IN) :: i ! Cell Number
890  REAL(SPA), INTENT(IN):: xloc,yloc
891  REAL(SPA), POINTER, INTENT(IN) :: FIELD(:)
892 
893  REAL(SPA):: X0c, Y0c,F0,Fx,Fy
894  INTEGER :: e1,e2,e3
895  REAL(SPA) :: dx_sph
896  IF(dbg_set(dbg_sbr)) WRITE(ipt,*) "START: INTERP_PZONAL_2D"
897 
898 !offset from element center
899 !---------------------------------------------------------------
900 
901  x0c = xloc - xc(i)
902  y0c = yloc - yc(i)
903 !---------------------------------------------------------------
904 !---------------------------------------------------------------
905  !Surrounding Element IDs
906  e1 = nbe(i,1)
907  e2 = nbe(i,2)
908  e3 = nbe(i,3)
909 
910  !interpolate Field to the location
911  fx = a1u(i,1)*field(i)+a1u(i,2)*field(e1)+a1u(i,3)*field(e2)+a1u(i,4)*field(e3)
912  fy = a2u(i,1)*field(i)+a2u(i,2)*field(e1)+a2u(i,3)*field(e2)+a2u(i,4)*field(e3)
913  fpt = field(i) + fx*x0c + fy*y0c
914 
915  IF(dbg_set(dbg_sbr)) WRITE(ipt,*) "START: INTERP_PZONAL_2D"
916 
real(sp), dimension(:), allocatable, target yc
Definition: mod_main.f90:1004
real(sp), dimension(:,:), allocatable, target a1u
Definition: mod_main.f90:1331
integer, dimension(:,:), allocatable, target nbe
Definition: mod_main.f90:1020
real(sp), dimension(:,:), allocatable, target a2u
Definition: mod_main.f90:1332
real(sp), dimension(:), allocatable, target xc
Definition: mod_main.f90:1003
integer ipt
Definition: mod_main.f90:922
Here is the call graph for this function:
Here is the caller graph for this function:

◆ interp_pzonal_3d_dbl()

real(dp) function mod_utils::interp_pzonal_3d_dbl ( real(dp), intent(in)  xloc,
real(dp), intent(in)  yloc,
real(dp), intent(in)  sigloc,
integer, intent(in)  lvls,
integer, intent(in)  i,
real(dp), dimension(:,:), intent(in), pointer  Field 
)

Definition at line 1083 of file mod_utils.f90.

1083  USE all_vars, only : a1u,a2u,xc,yc,nbe,kbm2,kbm1,kb,z1,dz1,zz1,dzz1
1084  IMPLICIT NONE
1085  REAL(DP) :: FPT
1086  INTEGER, INTENT(IN) :: i,lvls ! Cell Number, Number of levels(kb,or kbm1)
1087  REAL(DP), INTENT(IN):: xloc,yloc,sigloc
1088  REAL(DP), POINTER, INTENT(IN) :: FIELD(:,:)
1089 
1090  REAL(DP):: X0c, Y0c,F0,Fx,Fy,alpha,F_upper,F_lower
1091  INTEGER :: e1,e2,e3,k1,k2,k
1092 
1093  REAL(DP) :: dx_sph
1094  IF(dbg_set(dbg_sbr)) WRITE(ipt,*) "START: INTERP_PZONAL_3D"
1095 
1096 
1097 !offset from element center
1098 !---------------------------------------------------------------
1099  !offset from element center
1100  x0c = xloc - xc(i)
1101  y0c = yloc - yc(i)
1102 !---------------------------------------------------------------
1103 !---------------------------------------------------------------
1104 
1105  !Surrounding Element IDs
1106  e1 = nbe(i,1)
1107  e2 = nbe(i,2)
1108  e3 = nbe(i,3)
1109 
1110 !----determine sigma layers above and below sigloc
1111 
1112  IF(lvls == kbm1) THEN
1113 
1114  k1 = 1
1115  k2 = 1
1116  alpha = -1
1117  if(sigloc < zz1(i,kbm1)) then !!particle near bottom
1118  k1 = kbm1
1119  k2 = kbm1
1120  alpha = -1
1121  else
1122  do k=1,kbm2
1123  if(sigloc < zz1(i,k) .and. sigloc >= zz1(i,k+1) )then
1124  k1 = k
1125  k2 = k+1
1126  alpha = (zz1(i,k)-sigloc)/dzz1(i,k)
1127  exit
1128  endif
1129  end do
1130  end if
1131 
1132  ELSE IF(lvls == kb) THEN
1133 
1134  !surface (default)
1135  k1 = 1
1136  k2 = 1
1137  alpha = -1
1138  !bottom
1139  if(sigloc < z1(i,kb))then
1140  k1 = kb
1141  k2 = kb
1142  alpha = -1
1143  else !intermediate
1144  do k=1,kbm1
1145  if(sigloc < z1(i,k) .and. sigloc >= z1(i,k+1) )then
1146  k1 = k
1147  k2 = k+1
1148  alpha = (z1(i,k)-sigloc)/dz1(i,k)
1149  endif
1150  end do
1151  endif
1152 
1153  ELSE
1154  CALL fatal_error("INTERP_PZONAL_3D: Invalid number of levels passed",&
1155  & "(Must be equal to either KB or KBM1")
1156  END IF
1157 
1158  !interpolate Field to the location
1159  fx = a1u(i,1)*field(i,k1)+a1u(i,2)*field(e1,k1)+a1u(i,3)*field(e2,k1)+a1u(i,4)*field(e3,k1)
1160  fy = a2u(i,1)*field(i,k1)+a2u(i,2)*field(e1,k1)+a2u(i,3)*field(e2,k1)+a2u(i,4)*field(e3,k1)
1161  f_upper = field(i,k1) + fx*x0c + fy*y0c
1162 
1163  IF(k1 == k2) THEN
1164  fpt = f_upper
1165  ELSE
1166 
1167  fx = a1u(i,1)*field(i,k2)+a1u(i,2)*field(e1,k2)+a1u(i,3)*field(e2,k2)+a1u(i,4)*field(e3,k2)
1168  fy = a2u(i,1)*field(i,k2)+a2u(i,2)*field(e1,k2)+a2u(i,3)*field(e2,k2)+a2u(i,4)*field(e3,k2)
1169  f_lower = field(i,k2) + fx*x0c + fy*y0c
1170 
1171  fpt = (alpha)*f_lower + (1.0_dp-alpha)*f_upper
1172  END IF
1173 
1174 
1175  IF(dbg_set(dbg_sbr)) WRITE(ipt,*) "START: INTERP_PZONAL_3D"
1176 
real(sp), dimension(:), allocatable, target yc
Definition: mod_main.f90:1004
real(sp), dimension(:,:), allocatable, target dzz1
Definition: mod_main.f90:1097
real(sp), dimension(:,:), allocatable, target a1u
Definition: mod_main.f90:1331
integer kb
Definition: mod_main.f90:64
integer kbm2
Definition: mod_main.f90:66
integer, dimension(:,:), allocatable, target nbe
Definition: mod_main.f90:1020
real(sp), dimension(:,:), allocatable, target zz1
Definition: mod_main.f90:1095
real(sp), dimension(:,:), allocatable, target a2u
Definition: mod_main.f90:1332
real(sp), dimension(:), allocatable, target xc
Definition: mod_main.f90:1003
real(sp), dimension(:,:), allocatable, target dz1
Definition: mod_main.f90:1096
real(sp), dimension(:,:), allocatable, target z1
Definition: mod_main.f90:1094
integer ipt
Definition: mod_main.f90:922
integer kbm1
Definition: mod_main.f90:65
Here is the call graph for this function:
Here is the caller graph for this function:

◆ interp_pzonal_3d_flt()

real(spa) function mod_utils::interp_pzonal_3d_flt ( real(spa), intent(in)  xloc,
real(spa), intent(in)  yloc,
real(spa), intent(in)  sigloc,
integer, intent(in)  lvls,
integer, intent(in)  i,
real(spa), dimension(:,:), intent(in), pointer  Field 
)

Definition at line 988 of file mod_utils.f90.

988  USE all_vars, only : a1u,a2u,xc,yc,nbe,kbm2,kbm1,kb,z1,dz1,zz1,dzz1
989  IMPLICIT NONE
990  REAL(SPA) :: FPT
991  INTEGER, INTENT(IN) :: i,lvls ! Cell Number, Number of levels(kb,or kbm1)
992  REAL(SPA), INTENT(IN):: xloc,yloc,sigloc
993  REAL(SPA), POINTER, INTENT(IN) :: FIELD(:,:)
994 
995  REAL(SPA):: X0c, Y0c,F0,Fx,Fy,alpha,F_upper,F_lower
996  INTEGER :: e1,e2,e3,k1,k2,k
997  REAL(SPA) :: dx_sph
998  IF(dbg_set(dbg_sbr)) WRITE(ipt,*) "START: INTERP_PZONAL_3D"
999 
1000 !offset from element center
1001 !---------------------------------------------------------------
1002  !offset from element center
1003  x0c = xloc - xc(i)
1004  y0c = yloc - yc(i)
1005 !---------------------------------------------------------------
1006 !---------------------------------------------------------------
1007 
1008  !Surrounding Element IDs
1009  e1 = nbe(i,1)
1010  e2 = nbe(i,2)
1011  e3 = nbe(i,3)
1012 
1013 !----determine sigma layers above and below sigloc
1014 
1015  IF(lvls == kbm1) THEN
1016 
1017  k1 = 1
1018  k2 = 1
1019  alpha = -1
1020  if(sigloc < zz1(i,kbm1)) then !!particle near bottom
1021  k1 = kbm1
1022  k2 = kbm1
1023  alpha = -1
1024  else
1025  do k=1,kbm2
1026  if(sigloc < zz1(i,k) .and. sigloc >= zz1(i,k+1) )then
1027  k1 = k
1028  k2 = k+1
1029  alpha = (zz1(i,k)-sigloc)/dzz1(i,k)
1030  exit
1031  endif
1032  end do
1033  end if
1034 
1035  ELSE IF(lvls == kb) THEN
1036 
1037  !surface (default)
1038  k1 = 1
1039  k2 = 1
1040  alpha = -1
1041  !bottom
1042  if(sigloc < z1(i,kb))then
1043  k1 = kb
1044  k2 = kb
1045  alpha = -1
1046  else !intermediate
1047  do k=1,kbm1
1048  if(sigloc < z1(i,k) .and. sigloc >= z1(i,k+1) )then
1049  k1 = k
1050  k2 = k+1
1051  alpha = (z1(i,k)-sigloc)/dz1(i,k)
1052  endif
1053  end do
1054  endif
1055 
1056  ELSE
1057  CALL fatal_error("INTERP_PZONAL_3D: Invalid number of levels passed",&
1058  & "(Must be equal to either KB or KBM1")
1059  END IF
1060 
1061  !interpolate Field to the location
1062  fx = a1u(i,1)*field(i,k1)+a1u(i,2)*field(e1,k1)+a1u(i,3)*field(e2,k1)+a1u(i,4)*field(e3,k1)
1063  fy = a2u(i,1)*field(i,k1)+a2u(i,2)*field(e1,k1)+a2u(i,3)*field(e2,k1)+a2u(i,4)*field(e3,k1)
1064  f_upper = field(i,k1) + fx*x0c + fy*y0c
1065 
1066  IF(k1 == k2) THEN
1067  fpt = f_upper
1068  ELSE
1069 
1070  fx = a1u(i,1)*field(i,k2)+a1u(i,2)*field(e1,k2)+a1u(i,3)*field(e2,k2)+a1u(i,4)*field(e3,k2)
1071  fy = a2u(i,1)*field(i,k2)+a2u(i,2)*field(e1,k2)+a2u(i,3)*field(e2,k2)+a2u(i,4)*field(e3,k2)
1072  f_lower = field(i,k2) + fx*x0c + fy*y0c
1073 
1074  fpt = (alpha)*f_lower + (1.0-alpha)*f_upper
1075  END IF
1076 
1077 
1078  IF(dbg_set(dbg_sbr)) WRITE(ipt,*) "START: INTERP_PZONAL_3D"
1079 
real(sp), dimension(:), allocatable, target yc
Definition: mod_main.f90:1004
real(sp), dimension(:,:), allocatable, target dzz1
Definition: mod_main.f90:1097
real(sp), dimension(:,:), allocatable, target a1u
Definition: mod_main.f90:1331
integer kb
Definition: mod_main.f90:64
integer kbm2
Definition: mod_main.f90:66
integer, dimension(:,:), allocatable, target nbe
Definition: mod_main.f90:1020
real(sp), dimension(:,:), allocatable, target zz1
Definition: mod_main.f90:1095
real(sp), dimension(:,:), allocatable, target a2u
Definition: mod_main.f90:1332
real(sp), dimension(:), allocatable, target xc
Definition: mod_main.f90:1003
real(sp), dimension(:,:), allocatable, target dz1
Definition: mod_main.f90:1096
real(sp), dimension(:,:), allocatable, target z1
Definition: mod_main.f90:1094
integer ipt
Definition: mod_main.f90:922
integer kbm1
Definition: mod_main.f90:65
Here is the call graph for this function:
Here is the caller graph for this function:

◆ isintri()

logical function mod_utils::isintri ( real(sp), intent(in)  X0,
real(sp), intent(in)  Y0,
real(sp), dimension(3), intent(in)  Xt,
real(sp), dimension(3), intent(in)  Yt 
)

Definition at line 1315 of file mod_utils.f90.

1315  use mod_prec
1316  implicit none
1317  real(sp), intent(in) :: x0,Y0
1318  real(sp), intent(in) :: xt(3),yt(3)
1319  real(sp) :: f1,f2,f3
1320 
1321 !---------------------------------------------------------------
1322 !---------------------------------------------------------------
1323 
1324  isintri = .false.
1325 
1326  if(y0 < minval(yt) .or. y0 > maxval(yt)) then
1327  isintri = .false.
1328  return
1329  endif
1330  if(x0 < minval(xt) .or. x0 > maxval(xt)) then
1331  isintri = .false.
1332  return
1333  endif
1334 
1335  f1 = (y0-yt(1))*(xt(2)-xt(1)) - (x0-xt(1))*(yt(2)-yt(1))
1336  f2 = (y0-yt(3))*(xt(1)-xt(3)) - (x0-xt(3))*(yt(1)-yt(3))
1337  f3 = (y0-yt(2))*(xt(3)-xt(2)) - (x0-xt(2))*(yt(3)-yt(2))
1338  if(f1*f3 >= 0.0_sp .and. f3*f2 >= 0.0_sp) isintri = .true.
1339 !---------------------------------------------------------------
1340 !---------------------------------------------------------------
1341 
1342 
1343  return
1344 
Here is the caller graph for this function:

◆ isintriangle()

logical function mod_utils::isintriangle ( integer, intent(in)  i,
real(sp), intent(in)  x0,
real(sp), intent(in)  y0 
)

Definition at line 1350 of file mod_utils.f90.

1350 !==============================================================================|
1351 ! determine if point (x0,y0) is in triangle defined by nodes (xt(3),yt(3)) |
1352 ! using algorithm used for scene rendering in computer graphics |
1353 ! algorithm works well unless particle happens to lie in a line parallel |
1354 ! to the edge of a triangle. |
1355 ! This can cause problems if you use a regular grid, say for idealized |
1356 ! modelling and you happen to see particles right on edges or parallel to |
1357 ! edges. |
1358 !==============================================================================|
1359 
1360  use mod_prec
1361  use all_vars, only : nv, vx, vy
1362  implicit none
1363  real(sp), intent(in) :: x0,y0
1364  integer, intent(in) :: i
1365  real(sp) :: xt(3),yt(3)
1366  real(sp) :: f1,f2,f3
1367  real(sp) :: x1(2)
1368  real(sp) :: x2(2)
1369  real(sp) :: x3(2)
1370  real(sp) :: p(2)
1371 
1372 !------------------------------------------------------------------------------|
1373 
1374  isintriangle = .false.
1375 
1376  xt = vx(nv(i,1:3))
1377  yt = vy(nv(i,1:3))
1378 
1379  isintriangle = isintri(x0,y0,xt,yt)
1380 
1381 !
1382 ! if(sameside(p,x1,x2,x3).and.sameside(p,x2,x1,x3).and. &
1383 ! sameside(p,x3,x1,x2)) isintriangle = .true.
1384 !!$ if(y0 < minval(yt) .or. y0 > maxval(yt)) then
1385 !!$ isintriangle = .false.
1386 !!$ return
1387 !!$ endif
1388 !!$ if(x0 < minval(xt) .or. x0 > maxval(xt)) then
1389 !!$ isintriangle = .false.
1390 !!$ return
1391 !!$ endif
1392 !!$
1393 !!$ f1 = (y0-yt(1))*(xt(2)-xt(1)) - (x0-xt(1))*(yt(2)-yt(1))
1394 !!$ f2 = (y0-yt(3))*(xt(1)-xt(3)) - (x0-xt(3))*(yt(1)-yt(3))
1395 !!$ f3 = (y0-yt(2))*(xt(3)-xt(2)) - (x0-xt(2))*(yt(3)-yt(2))
1396 !!$ if(f1*f3 >= 0.0_sp .and. f3*f2 >= 0.0_sp) isintriangle = .true.
1397 
1398  return
real(sp), dimension(:), allocatable, target vx
Definition: mod_main.f90:1001
real(sp), dimension(:), allocatable, target vy
Definition: mod_main.f90:1002
integer, dimension(:,:), allocatable, target nv
Definition: mod_main.f90:1018
Here is the call graph for this function:
Here is the caller graph for this function:

◆ limled()

real(sp) function mod_utils::limled ( real(sp)  a,
real(sp)  b,
real(sp)  q 
)

Definition at line 2216 of file mod_utils.f90.

2216 
2217  IMPLICIT NONE
2218 
2219  real(sp) a,b
2220  real(sp) q,R
2221  real(sp) eps
2222  eps = epsilon(eps)
2223 
2224  ! exponent
2225  ! q = 0. !1st order
2226  ! q = 1. !minmod
2227  ! q = 2. !van leer
2228 
2229  r = abs( (a-b)/(abs(a)+abs(b)+eps) )**q
2230  lim = .5*(1-r)*(a+b)
2231 

◆ limled1()

real(sp) function mod_utils::limled1 ( real(sp)  a,
real(sp)  b,
real(sp)  alpha 
)

Definition at line 2235 of file mod_utils.f90.

2235  IMPLICIT NONE
2236 
2237  real(sp) a,b,alpha
2238 
2239  lim = 0.5_sp*(sign(1.,a)+sign(1.,b))*max(min(alpha*abs(a),abs(b)),min(abs(a),alpha*abs(b)))
2240 

◆ limled2()

real(sp) function mod_utils::limled2 ( real(sp)  a,
real(sp)  b,
real(sp)  alpha 
)

Definition at line 2244 of file mod_utils.f90.

2244  IMPLICIT NONE
2245 
2246  real(sp) a,b,alpha
2247 
2248  lim = 0.5_sp*(sign(1.,a)+sign(1.,b))*min(0.5_sp*abs(a+b),alpha*abs(a),alpha*abs(b))
2249 
Here is the caller graph for this function:

◆ meters2degrees_arr_dbl()

subroutine mod_utils::meters2degrees_arr_dbl ( real(dp), dimension(nsze,msze), intent(in)  X,
real(dp), dimension(nsze,msze), intent(in)  Y,
character(len=*), intent(in)  proj_ref,
real(dp), dimension(nsze,msze), intent(out)  LON,
real(dp), dimension(nsze,msze), intent(out)  LAT,
integer, intent(in)  nsze,
integer, intent(in)  msze 
)

Definition at line 504 of file mod_utils.f90.

504  implicit none
505  integer, intent(in) :: nsze, msze
506  REAL(DP), INTENT(OUT) :: LON(nsze,msze), LAT(nsze,msze)
507  character(LEN=*), INTENT(IN) :: proj_ref
508  REAL(DP), intent(IN):: X(nsze,msze), Y(nsze,msze)
509 
510  lon=-1.0_dp
511  lat=-1.0_dp

◆ meters2degrees_arr_flt()

subroutine mod_utils::meters2degrees_arr_flt ( real(spa), dimension(nsze,msze), intent(in)  X,
real(spa), dimension(nsze,msze), intent(in)  Y,
character(len=*), intent(in)  proj_ref,
real(spa), dimension(nsze,msze), intent(out)  LON,
real(spa), dimension(nsze,msze), intent(out)  LAT,
integer, intent(in)  nsze,
integer, intent(in)  msze 
)

Definition at line 470 of file mod_utils.f90.

470  implicit none
471  integer, intent(in) :: nsze, msze
472  REAL(SPA), INTENT(OUT) :: LON(nsze,msze), LAT(nsze,msze)
473  character(LEN=*), INTENT(IN) :: proj_ref
474  REAL(SPA), intent(IN):: X(nsze,msze), Y(nsze,msze)
475 
476  lon=-1.0_spa
477  lat=-1.0_spa

◆ meters2degrees_scl_dbl()

subroutine mod_utils::meters2degrees_scl_dbl ( real(dp), intent(in)  X,
real(dp), intent(in)  Y,
character(len=*), intent(in)  proj_ref,
real(dp), intent(out)  LON,
real(dp), intent(out)  LAT 
)

Definition at line 481 of file mod_utils.f90.

481 ! USE PROJ4
482  implicit none
483  REAL(DP), INTENT(OUT) :: LON, LAT
484  character(LEN=*), INTENT(IN) :: proj_ref
485  REAL(DP), intent(IN):: X, Y
486 
487  lon=-1.0_dp
488  lat=-1.0_dp

◆ meters2degrees_scl_flt()

subroutine mod_utils::meters2degrees_scl_flt ( real(spa), intent(in)  X,
real(spa), intent(in)  Y,
character(len=*), intent(in)  proj_ref,
real(spa), intent(out)  LON,
real(spa), intent(out)  LAT 
)

Definition at line 449 of file mod_utils.f90.

449  implicit none
450  REAL(SPA), INTENT(OUT) :: LON, LAT
451  character(LEN=*), INTENT(IN) :: proj_ref
452  REAL(SPA), intent(IN):: X, Y
453 
454  lon=-1.0_spa
455  lat=-1.0_spa

◆ meters2degrees_vec_dbl()

subroutine mod_utils::meters2degrees_vec_dbl ( real(dp), dimension(nsze), intent(in)  X,
real(dp), dimension(nsze), intent(in)  Y,
character(len=*), intent(in)  proj_ref,
real(dp), dimension(nsze), intent(out)  LON,
real(dp), dimension(nsze), intent(out)  LAT,
integer, intent(in)  nsze 
)

Definition at line 492 of file mod_utils.f90.

492 ! USE PROJ4
493  implicit none
494  integer, intent(in) :: nsze
495  REAL(DP), INTENT(OUT) :: LON(nsze), LAT(nsze)
496  character(LEN=*), INTENT(IN) :: proj_ref
497  REAL(DP), intent(IN):: X(nsze), Y(nsze)
498 
499  lon=-1.0_dp
500  lat=-1.0_dp

◆ meters2degrees_vec_flt()

subroutine mod_utils::meters2degrees_vec_flt ( real(spa), dimension(nsze), intent(in)  X,
real(spa), dimension(nsze), intent(in)  Y,
character(len=*), intent(in)  proj_ref,
real(spa), dimension(nsze), intent(out)  LON,
real(spa), dimension(nsze), intent(out)  LAT,
integer, intent(in)  nsze 
)

Definition at line 459 of file mod_utils.f90.

459  implicit none
460  integer, intent(in) :: nsze
461  REAL(SPA), INTENT(OUT) :: LON(nsze), LAT(nsze)
462  character(LEN=*), INTENT(IN) :: proj_ref
463  REAL(SPA), intent(IN):: X(nsze), Y(nsze)
464 
465  lon=-1.0_spa
466  lat=-1.0_spa

◆ open_dat()

integer function mod_utils::open_dat ( character(len=*)  FNAME,
integer  UNIT,
character(len=*), optional  PATH 
)

Definition at line 1644 of file mod_utils.f90.

1644  USE control, only : input_dir
1645  implicit none
1646  CHARACTER(LEN=*) :: FNAME
1647  INTEGER :: UNIT
1648  CHARACTER(LEN=*), OPTIONAL :: PATH
1649 
1650  CHARACTER(LEN=400) :: PATHNFILE
1651 
1652  open_dat = -1
1653  IF (len_trim(fname) ==0) return
1654 
1655  IF (PRESENT(path)) THEN
1656  IF (len_trim(path) ==0) return
1657 
1658  pathnfile = trim(path)//trim(fname)
1659  ELSE
1660  IF (len_trim(input_dir) ==0) return
1661  pathnfile = trim(input_dir)//trim(fname)
1662  END IF
1663  Call fopen(unit,trim(pathnfile),'cfr')
1664  open_dat = 0
1665 
1666 
character(len=80) input_dir
Definition: mod_main.f90:183
Here is the call graph for this function:

◆ path_split()

subroutine mod_utils::path_split ( character(len=*), intent(in)  STRING,
character(len=*), intent(out)  PATH,
character(len=*), intent(out)  FILE,
character(len=*), intent(out)  EXTENSION 
)

Definition at line 2146 of file mod_utils.f90.

2146  IMPLICIT NONE
2147 
2148  CHARACTER(LEN=*), INTENT(IN) :: STRING
2149  CHARACTER(LEN=*), INTENT(OUT):: PATH,FILE,EXTENSION
2150 
2151  INTEGER :: IDX, LGN, I
2152 
2153  lgn = len_trim(string)
2154 
2155  path = ''
2156  file = ''
2157  extension = ''
2158 
2159  idx = 1
2160  DO i=1,lgn
2161  if(string(i:i) == '/') idx = i
2162  END DO
2163 
2164 
2165  IF(idx>1) THEN
2166  path = string(1:idx)
2167  idx = idx+1
2168  ELSE ! HANDLE THE CASE OF NO PATH
2169  path = './'
2170  END IF
2171 
2172 
2173  DO i = idx,lgn
2174  IF(string(i:i) == '.') THEN
2175  file = string(idx:(i-1))
2176  extension =string(i:lgn)
2177  EXIT
2178  END IF
2179  END DO
2180 
2181  ! HANDLE THE CASE OF NO EXTENSION
2182  IF (len_trim(file) == 0) file = string(idx:lgn)
2183 
2184 

◆ pshutdown()

subroutine mod_utils::pshutdown ( )

Definition at line 280 of file mod_utils.f90.

280  !==============================================================================|
281  stop
Here is the caller graph for this function:

◆ pstop()

subroutine mod_utils::pstop ( )

Definition at line 273 of file mod_utils.f90.

273  !==============================================================================|
274  stop
Here is the caller graph for this function:

◆ read_float()

real(dp) function mod_utils::read_float ( character(len=*), intent(in)  ITEM,
integer, intent(out)  IERR 
)

Definition at line 2253 of file mod_utils.f90.

2253  IMPLICIT NONE
2254  CHARACTER(LEN=*),INTENT(IN) :: ITEM
2255  INTEGER, INTENT(OUT) :: IERR
2256 
2257  LOGICAL :: ISFLOAT
2258  INTEGER :: I
2259 
2260  ierr = -1
2261  isfloat = .false.
2262  read_float = -99999.9
2263  DO i = 1,len_trim(item)
2264  IF (item(i:i) == ".") THEN
2265  isfloat = .true.
2266  cycle
2267  END IF
2268 
2269  IF( lgt("0",item(i:i)) .or. llt("9",item(i:i)) ) THEN
2270  isfloat = .false.
2271  EXIT
2272  END IF
2273 
2274  END DO
2275 
2276  IF(isfloat) THEN
2277  READ(item,*,iostat=ierr) read_float
2278  END IF
2279 

◆ read_int()

real(sp) function mod_utils::read_int ( character(len=*), intent(in)  ITEM,
integer, intent(out)  IERR 
)

Definition at line 2283 of file mod_utils.f90.

2283  IMPLICIT NONE
2284  CHARACTER(LEN=*),INTENT(IN) :: ITEM
2285  INTEGER, INTENT(OUT) :: IERR
2286 
2287  INTEGER :: I
2288 
2289  ierr = -1
2290  read_int = -99999.9
2291  DO i = 1,len_trim(item)
2292 
2293  IF( lgt("0",item(i:i)) .or. llt("9",item(i:i)) ) THEN
2294  RETURN
2295  END IF
2296 
2297  END DO
2298 
2299  READ(item,*,iostat=ierr) read_int
2300 

◆ sameside()

logical function mod_utils::sameside ( real(sp), dimension(2), intent(in)  p1,
real(sp), dimension(2), intent(in)  p2,
real(sp), dimension(2), intent(in)  a,
real(sp), dimension(2), intent(in)  b 
)

Definition at line 1402 of file mod_utils.f90.

1402  real(sp), intent(in) :: p1(2)
1403  real(sp), intent(in) :: p2(2)
1404  real(sp), intent(in) :: a(2)
1405  real(sp), intent(in) :: b(2)
1406  logical value
1407  real(sp) :: cp1,cp2
1408 
1409  cp1 = (b(1)-a(1))*(p1(2)-a(2)) - (b(2)-a(2))*(p1(1)-a(1))
1410  cp2 = (b(1)-a(1))*(p2(2)-a(2)) - (b(2)-a(2))*(p2(1)-a(1))
1411 
1412  value = .false.
1413  if(cp1*cp2 >= 0) value = .true.
1414 
real(kind=dbl_kind), parameter p2
real(kind=dbl_kind), parameter p1

◆ scan_file()

integer function mod_utils::scan_file ( integer, intent(in)  UNIT,
character(len=*)  VNAME,
integer, intent(inout), optional  ISCAL,
real(sp), intent(inout), optional  FSCAL,
integer, dimension(*), intent(inout), optional  IVEC,
real(sp), dimension(*), intent(inout), optional  FVEC,
character(len=80), dimension(*), optional  CVEC,
integer, intent(inout), optional  NSZE,
character(len=80), optional  CVAL,
logical, intent(inout), optional  LVAL 
)

Definition at line 1844 of file mod_utils.f90.

1844 
1845 !==============================================================================|
1846 ! Scan an Input File for a Variable |
1847 ! RETURN VALUE: |
1848 ! 0 = FILE FOUND, VARIABLE VALUE FOUND |
1849 ! -1 = FILE DOES NOT EXIST OR PERMISSIONS ARE INCORRECT |
1850 ! -2 = VARIABLE NOT FOUND OR IMPROPERLY SET |
1851 ! -3 = VARIABLE IS OF DIFFERENT TYPE, CHECK INPUT FILE |
1852 ! -4 = VECTOR PROVIDED BUT DATA IS SCALAR TYPE |
1853 ! -5 = NO DATATYPE DESIRED, EXITING |
1854 ! |
1855 ! REQUIRED INPUT: |
1856 ! UNIT = File UNIT |
1857 ! FSIZE = Length of Filename |
1858 ! |
1859 ! OPTIONAL (MUST PROVIDE ONE) |
1860 ! ISCAL = INTEGER SCALAR |
1861 ! FSCAL = FLOAT SCALAR |
1862 ! CVAL = CHARACTER VARIABLE |
1863 ! LVAL = LOGICAL VARIABLE |
1864 ! IVEC = INTEGER VECTOR ** |
1865 ! FVEC = FLOAT VECTOR ** |
1866 ! CVEC = STRING VECTOR ** (STRINGS OF LENGTH 80) |
1867 ! **NSZE = ARRAY SIZE (MUST BE PROVIDED WITH IVEC/FVEC) |
1868 ! |
1869 !==============================================================================|
1870 
1871  IMPLICIT NONE
1872  INTEGER, INTENT(IN) :: UNIT
1873  CHARACTER(LEN=*) :: VNAME
1874  INTEGER, INTENT(INOUT), OPTIONAL :: ISCAL,IVEC(*)
1875  REAL(SP),INTENT(INOUT), OPTIONAL :: FSCAL,FVEC(*)
1876  CHARACTER(LEN=80), OPTIONAL :: CVAL,CVEC(*)
1877  LOGICAL, INTENT(INOUT), OPTIONAL :: LVAL
1878  INTEGER, INTENT(INOUT), OPTIONAL :: NSZE
1879 
1880 !------------------------------------------------------------------------------|
1881  INTEGER :: SCAN_FILE
1882  REAL(DP) REALVAL(150)
1883  INTEGER INTVAL(150)
1884  CHARACTER(LEN=40 ) :: VARNAME
1885  CHARACTER(LEN=80 ) :: STRINGVAL(150),TITLE
1886  CHARACTER(LEN=400 ) :: INPLINE
1887  CHARACTER(LEN=800) :: TLINE
1888  CHARACTER(LEN=7 ) :: VARTYPE, ENDLINE
1889  CHARACTER(LEN=20 ), DIMENSION(200) :: SET
1890  INTEGER I,NVAL,J,NSET,NLINE,NREP,bgn,nd, LEL
1891  LOGICAL SETYES,ALLSET,CHECK,LOGVAL
1892 
1893 
1894  scan_file = 0
1895 
1896 !==============================================================================|
1897 ! SCAN THE FILE FOR THE VARIABLE NAME |
1898 !==============================================================================|
1899 
1900  ENDLINE='\\'
1901  lel = len_trim(endline)-1
1902 
1903  rewind(unit)
1904 
1905  nset = 0
1906  nline = 0
1907  DO WHILE(.true.)
1908 
1909  if(nline >200) then
1910  CALL warning("Read 200 lines of header with out finding parameters! ")
1911  scan_file=-2
1912  return
1913  end if
1914 
1915  tline(1:len(tline)) = ' '
1916  nrep = 0
1917  nline = nline + 1
1918  READ(unit,'(a)',end=20) inpline
1919  tline = trim(inpline)
1920 
1921 !----PROCESS LINE CONTINUATIONS------------------------------------------------!
1922  DO
1923 ! write(ipt,*) '"'//TRIM(INPLINE)//'"'
1924  i = len_trim(inpline)
1925  IF(i > 1)THEN
1926 
1927  IF( inpline(i-lel:i) == trim(endline))THEN
1928 
1929  nrep = nrep + 1
1930  READ(unit,'(a)',end=20) inpline
1931  nline = nline + 1
1932  bgn = len_trim(tline)+1
1933  nd = bgn +len_trim(inpline)
1934 
1935  tline(bgn:nd) = trim(inpline)
1936  ELSE
1937  EXIT
1938 
1939  END IF
1940  ELSE
1941  EXIT
1942  END IF
1943  END DO
1944 
1945 !----REMOVE LINE CONTINUATION CHARACTER \\-------------------------------------!
1946  IF(nrep > 0)THEN
1947  DO i=2,len_trim(tline)
1948  IF( tline(i-lel:i) == endline) tline(i-lel:i) = ' '
1949 
1950  END DO
1951  END IF
1952 
1953 !----PROCESS THE LINE----------------------------------------------------------!
1954  CALL get_value(nline,len_trim(tline),adjustl(tline),varname,vartype,logval,&
1955  stringval,realval,intval,nval)
1956 
1957 !----IF VARNAME MATCHES, PROCESS VARIABLE AND ERROR-CHECK----------------------!
1958 
1959  IF(trim(varname) == trim(vname))THEN
1960 
1961  IF(PRESENT(iscal))THEN
1962  IF(vartype == 'integer')THEN
1963  iscal = intval(1)
1964  RETURN
1965  ELSE
1966  scan_file = -3
1967  return
1968  END IF
1969  ELSE IF(PRESENT(fscal))THEN
1970  IF(vartype == 'float')THEN
1971  fscal = realval(1)
1972  RETURN
1973  ELSE
1974  scan_file = -3
1975  return
1976  END IF
1977  ELSE IF(PRESENT(cval))THEN
1978  IF(vartype == 'string')THEN
1979  cval = stringval(1)
1980  RETURN
1981  ELSE
1982  scan_file = -3
1983  return
1984  END IF
1985  ELSE IF(PRESENT(lval))THEN
1986  IF(vartype == 'logical')THEN
1987  lval = logval
1988  RETURN
1989  ELSE
1990  scan_file = -3
1991  return
1992  END IF
1993  ELSE IF(PRESENT(ivec))THEN
1994  IF(nval > 1)THEN
1995  IF(vartype == 'integer')THEN
1996  ivec(1:nval) = intval(1:nval)
1997  nsze = nval
1998  RETURN
1999  ELSE
2000  scan_file = -3
2001  return
2002  END IF
2003  ELSE
2004  scan_file = -4
2005  return
2006  END IF
2007  ELSE IF(PRESENT(fvec))THEN
2008  IF(nval > 1)THEN
2009  IF(vartype == 'float')THEN
2010  fvec(1:nval) = realval(1:nval)
2011  nsze = nval
2012  RETURN
2013  ELSE
2014  scan_file = -3
2015  return
2016  END IF
2017  ELSE
2018  scan_file = -4
2019  return
2020  END IF
2021  ELSE IF(PRESENT(cvec))THEN
2022  IF(nval > 0)THEN
2023  IF(vartype == 'string')THEN
2024  cvec(1:nval) = stringval(2:nval+1)
2025  nsze = nval
2026  RETURN
2027  ELSE
2028  scan_file = -3
2029  return
2030  END IF
2031  ELSE
2032  scan_file = -4
2033  return
2034  END IF
2035  ELSE
2036  scan_file = -5
2037  return
2038  END IF
2039  END IF !!VARIABLE IS CORRECT
2040 
2041  END DO !!LOOP OVER INPUT FILE
2042 
2043 20 scan_file = -2
2044 
2045  RETURN
Here is the call graph for this function:

◆ scan_file2()

integer function mod_utils::scan_file2 ( character(len=*)  FNAME,
character(len=*)  VNAME,
integer, intent(inout), optional  ISCAL,
real(sp), intent(inout), optional  FSCAL,
integer, dimension(*), intent(inout), optional  IVEC,
real(sp), dimension(*), intent(inout), optional  FVEC,
character(len=80), dimension(*), optional  CVEC,
integer, intent(inout), optional  NSZE,
character(len=80), optional  CVAL,
logical, intent(inout), optional  LVAL 
)

Definition at line 2304 of file mod_utils.f90.

2304 
2305 !==============================================================================|
2306 ! Scan an Input File for a Variable |
2307 ! RETURN VALUE: |
2308 ! 0 = FILE FOUND, VARIABLE VALUE FOUND |
2309 ! -1 = FILE DOES NOT EXIST OR PERMISSIONS ARE INCORRECT |
2310 ! -2 = VARIABLE NOT FOUND OR IMPROPERLY SET |
2311 ! -3 = VARIABLE IS OF DIFFERENT TYPE, CHECK INPUT FILE |
2312 ! -4 = VECTOR PROVIDED BUT DATA IS SCALAR TYPE |
2313 ! -5 = NO DATATYPE DESIRED, EXITING |
2314 ! |
2315 ! REQUIRED INPUT: |
2316 ! FNAME = File Name |
2317 ! FSIZE = Length of Filename |
2318 ! |
2319 ! OPTIONAL (MUST PROVIDE ONE) |
2320 ! ISCAL = INTEGER SCALAR |
2321 ! FSCAL = FLOAT SCALAR |
2322 ! CVAL = CHARACTER VARIABLE |
2323 ! LVAL = LOGICAL VARIABLE |
2324 ! IVEC = INTEGER VECTOR ** |
2325 ! FVEC = FLOAT VECTOR ** |
2326 ! CVEC = STRING VECTOR ** (STRINGS OF LENGTH 80) |
2327 ! **NSZE = ARRAY SIZE (MUST BE PROVIDED WITH IVEC/FVEC) |
2328 ! |
2329 !==============================================================================|
2330 
2331  USE mod_prec
2332  USE ocpcomm4
2333  IMPLICIT NONE
2334  CHARACTER(LEN=*) :: FNAME,VNAME
2335  INTEGER, INTENT(INOUT), OPTIONAL :: ISCAL,IVEC(*)
2336  REAL(SP),INTENT(INOUT), OPTIONAL :: FSCAL,FVEC(*)
2337  CHARACTER(LEN=80), OPTIONAL :: CVAL,CVEC(*)
2338  LOGICAL, INTENT(INOUT), OPTIONAL :: LVAL
2339  INTEGER, INTENT(INOUT), OPTIONAL :: NSZE
2340 
2341 !------------------------------------------------------------------------------|
2342 
2343  INTEGER :: SCAN_FILE2
2344  REAL(DP) REALVAL(150)
2345  INTEGER INTVAL(150)
2346  CHARACTER(LEN=40 ) :: VARNAME
2347  CHARACTER(LEN=80 ) :: STRINGVAL(150),TITLE
2348  CHARACTER(LEN=80 ) :: INPLINE
2349  CHARACTER(LEN=400) :: TLINE
2350  CHARACTER(LEN=7 ) :: VARTYPE
2351  CHARACTER(LEN=20 ), DIMENSION(200) :: SET
2352  INTEGER I,NVAL,J,NSET,NLINE,NREP
2353  LOGICAL SETYES,ALLSET,CHECK,LOGVAL
2354 
2355 
2356  scan_file2 = 0
2357 !==============================================================================|
2358 ! OPEN THE INPUT FILE |
2359 !==============================================================================|
2360  INQUIRE(file=trim(fname),exist=check)
2361  IF(.NOT.check)THEN
2362  scan_file2 = -1
2363  RETURN
2364  END IF
2365 
2366  OPEN(inputf,file=trim(fname)) ; rewind(inputf)
2367 
2368 !==============================================================================|
2369 ! SCAN THE FILE FOR THE VARIABLE NAME |
2370 !==============================================================================|
2371 
2372  nset = 0
2373  nline = 0
2374  DO WHILE(.true.)
2375  tline(1:len(tline)) = ' '
2376  nrep = 0
2377  nline = nline + 1
2378  READ(inputf,'(a)',end=20) inpline
2379  tline(1:80) = inpline(1:80)
2380 
2381 !----PROCESS LINE CONTINUATIONS------------------------------------------------!
2382  110 CONTINUE
2383  i = len_trim(inpline)
2384  IF(i /= 0)THEN
2385  IF( inpline(i-1:i) == '\\\\')THEN
2386  nrep = nrep + 1
2387  READ(inputf,'(a)',end=20) inpline
2388  nline = nline + 1
2389  tline( nrep*80 + 1 : nrep*80 +80) = inpline(1:80)
2390  GOTO 110
2391  END IF
2392  END IF
2393  IF(nrep > 4)CALL fatal_error("CANNOT HAVE > 4 LINE CONTINUATIONS")
2394 
2395 !----REMOVE LINE CONTINUATION CHARACTER \\-------------------------------------!
2396  IF(nrep > 0)THEN
2397  DO i=2,len_trim(tline)
2398  IF( tline(i-1:i) == '\\\\') tline(i-1:i) = ' '
2399  END DO
2400  END IF
2401 
2402 !----PROCESS THE LINE----------------------------------------------------------!
2403  CALL get_value(nline,len_trim(tline),adjustl(tline),varname,vartype,logval,&
2404  stringval,realval,intval,nval)
2405 
2406 !----IF VARNAME MATCHES, PROCESS VARIABLE AND ERROR-CHECK----------------------!
2407 
2408  IF(trim(varname) == trim(vname))THEN
2409 
2410  IF(PRESENT(iscal))THEN
2411  IF(vartype == 'integer')THEN
2412  iscal = intval(1)
2413  RETURN
2414  ELSE
2415  scan_file2 = -3
2416  END IF
2417  ELSE IF(PRESENT(fscal))THEN
2418  IF(vartype == 'float')THEN
2419  fscal = realval(1)
2420  RETURN
2421  ELSE
2422  scan_file2 = -3
2423  END IF
2424  ELSE IF(PRESENT(cval))THEN
2425  IF(vartype == 'string')THEN
2426  cval = stringval(1)
2427  RETURN
2428  ELSE
2429  scan_file2 = -3
2430  END IF
2431  ELSE IF(PRESENT(lval))THEN
2432  IF(vartype == 'logical')THEN
2433  lval = logval
2434  RETURN
2435  ELSE
2436  scan_file2 = -3
2437  END IF
2438  ELSE IF(PRESENT(ivec))THEN
2439  IF(nval > 1)THEN
2440  IF(vartype == 'integer')THEN
2441  ivec(1:nval) = intval(1:nval)
2442  nsze = nval
2443  RETURN
2444  ELSE
2445  scan_file2 = -3
2446  END IF
2447  ELSE
2448  scan_file2 = -4
2449  END IF
2450  ELSE IF(PRESENT(fvec))THEN
2451  IF(nval > 1)THEN
2452  IF(vartype == 'float')THEN
2453  fvec(1:nval) = realval(1:nval)
2454  nsze = nval
2455  RETURN
2456  ELSE
2457  scan_file2 = -3
2458  END IF
2459  ELSE
2460  scan_file2 = -4
2461  END IF
2462  ELSE IF(PRESENT(cvec))THEN
2463  IF(nval > 0)THEN
2464  IF(vartype == 'string')THEN
2465  cvec(1:nval) = stringval(2:nval+1)
2466  nsze = nval
2467  RETURN
2468  ELSE
2469  scan_file2 = -3
2470  END IF
2471  ELSE
2472  scan_file2 = -4
2473  END IF
2474  ELSE
2475  scan_file2 = -5
2476  END IF
2477  END IF !!VARIABLE IS CORRECT
2478 
2479  END DO !!LOOP OVER INPUT FILE
2480  20 CLOSE(inputf)
2481  scan_file2 = -2
2482 
2483  RETURN
integer inputf
Definition: swmod1.f90:516
Here is the call graph for this function:
Here is the caller graph for this function:

◆ scan_file3()

integer function mod_utils::scan_file3 ( character(len=*)  FNAME,
character(len=*)  VNAME,
integer, intent(inout), optional  ISCAL,
real(sp), intent(inout), optional  FSCAL,
integer, dimension(*), intent(inout), optional  IVEC,
real(sp), dimension(*), intent(inout), optional  FVEC,
character(len=80), dimension(*), optional  CVEC,
integer, intent(inout), optional  NSZE,
character(len=80), optional  CVAL,
logical, intent(inout), optional  LVAL 
)

Definition at line 2488 of file mod_utils.f90.

2488 
2489 !==============================================================================|
2490 ! Scan an Input File for a Variable |
2491 ! RETURN VALUE: |
2492 ! 0 = FILE FOUND, VARIABLE VALUE FOUND |
2493 ! -1 = FILE DOES NOT EXIST OR PERMISSIONS ARE INCORRECT |
2494 ! -2 = VARIABLE NOT FOUND OR IMPROPERLY SET |
2495 ! -3 = VARIABLE IS OF DIFFERENT TYPE, CHECK INPUT FILE |
2496 ! -4 = VECTOR PROVIDED BUT DATA IS SCALAR TYPE |
2497 ! -5 = NO DATATYPE DESIRED, EXITING |
2498 ! |
2499 ! REQUIRED INPUT: |
2500 ! FNAME = File Name |
2501 ! FSIZE = Length of Filename |
2502 ! |
2503 ! OPTIONAL (MUST PROVIDE ONE) |
2504 ! ISCAL = INTEGER SCALAR |
2505 ! FSCAL = FLOAT SCALAR |
2506 ! CVAL = CHARACTER VARIABLE |
2507 ! LVAL = LOGICAL VARIABLE |
2508 ! IVEC = INTEGER VECTOR ** |
2509 ! FVEC = FLOAT VECTOR ** |
2510 ! CVEC = STRING VECTOR ** (STRINGS OF LENGTH 80) |
2511 ! **NSZE = ARRAY SIZE (MUST BE PROVIDED WITH IVEC/FVEC) |
2512 ! |
2513 !==============================================================================|
2514 
2515  USE mod_prec
2516  USE ocpcomm4
2517  IMPLICIT NONE
2518  CHARACTER(LEN=*) :: FNAME,VNAME
2519  INTEGER, INTENT(INOUT), OPTIONAL :: ISCAL,IVEC(*)
2520  REAL(SP),INTENT(INOUT), OPTIONAL :: FSCAL,FVEC(*)
2521  CHARACTER(LEN=80), OPTIONAL :: CVAL,CVEC(*)
2522  LOGICAL, INTENT(INOUT), OPTIONAL :: LVAL
2523  INTEGER, INTENT(INOUT), OPTIONAL :: NSZE
2524 
2525 !------------------------------------------------------------------------------|
2526 
2527  INTEGER :: SCAN_FILE3
2528  REAL(DP) REALVAL(150)
2529  INTEGER INTVAL(150)
2530  CHARACTER(LEN=40 ) :: VARNAME
2531  CHARACTER(LEN=80 ) :: STRINGVAL(150),TITLE
2532  CHARACTER(LEN=80 ) :: INPLINE
2533  CHARACTER(LEN=400) :: TLINE
2534  CHARACTER(LEN=7 ) :: VARTYPE
2535  CHARACTER(LEN=20 ), DIMENSION(200) :: SET
2536  INTEGER I,NVAL,J,NSET,NLINE,NREP
2537  LOGICAL SETYES,ALLSET,CHECK,LOGVAL
2538 
2539 
2540  scan_file3 = 0
2541 !==============================================================================|
2542 ! OPEN THE INPUT FILE |
2543 !==============================================================================|
2544  INQUIRE(file=trim(fname),exist=check)
2545  IF(.NOT.check)THEN
2546  scan_file3 = -1
2547  RETURN
2548  END IF
2549 
2550  OPEN(inputf,file=trim(fname)) ; rewind(inputf)
2551 
2552 !==============================================================================|
2553 ! SCAN THE FILE FOR THE VARIABLE NAME |
2554 !==============================================================================|
2555 
2556  nset = 0
2557  nline = 0
2558  DO WHILE(.true.)
2559  tline(1:len(tline)) = ' '
2560  nrep = 0
2561  nline = nline + 1
2562  READ(inputf,'(a)',end=20) inpline
2563  tline(1:80) = inpline(1:80)
2564 
2565 !----PROCESS LINE CONTINUATIONS------------------------------------------------!
2566  110 CONTINUE
2567  i = len_trim(inpline)
2568  IF(i > 2)THEN
2569  IF( inpline(i-1:i) == '\\')THEN
2570  nrep = nrep + 1
2571  READ(inputf,'(a)',end=20) inpline
2572  nline = nline + 1
2573  tline( nrep*80 + 1 : nrep*80 +80) = inpline(1:80)
2574  GOTO 110
2575  END IF
2576  END IF
2577  IF(nrep > 4)CALL fatal_error("CANNOT HAVE > 4 LINE CONTINUATIONS")
2578 
2579 !----REMOVE LINE CONTINUATION CHARACTER \\-------------------------------------!
2580  IF(nrep > 0)THEN
2581  DO i=2,len_trim(tline)
2582  IF( tline(i-1:i) == '\\') tline(i-1:i) = ' '
2583  END DO
2584  END IF
2585 
2586 
2587 !----PROCESS THE LINE----------------------------------------------------------!
2588  CALL get_value(nline,len_trim(tline),adjustl(tline),varname,vartype,logval,&
2589  stringval,realval,intval,nval)
2590 
2591 !----IF VARNAME MATCHES, PROCESS VARIABLE AND ERROR-CHECK----------------------!
2592 
2593  IF(trim(varname) == trim(vname))THEN
2594 
2595  IF(PRESENT(iscal))THEN
2596  IF(vartype == 'integer')THEN
2597  iscal = intval(1)
2598  CLOSE(inputf)
2599  RETURN
2600  ELSE
2601  scan_file3 = -3
2602  END IF
2603  ELSE IF(PRESENT(fscal))THEN
2604  IF(vartype == 'float')THEN
2605  fscal = realval(1)
2606  CLOSE(inputf)
2607  RETURN
2608  ELSE
2609  scan_file3 = -3
2610  END IF
2611  ELSE IF(PRESENT(cval))THEN
2612  IF(vartype == 'string')THEN
2613  cval = stringval(1)
2614  CLOSE(inputf)
2615  RETURN
2616  ELSE
2617  scan_file3 = -3
2618  END IF
2619  ELSE IF(PRESENT(lval))THEN
2620  IF(vartype == 'logical')THEN
2621  lval = logval
2622  CLOSE(inputf)
2623  RETURN
2624  ELSE
2625  scan_file3 = -3
2626  END IF
2627  ELSE IF(PRESENT(ivec))THEN
2628  IF(nval > 1)THEN
2629  IF(vartype == 'integer')THEN
2630  ivec(1:nval) = intval(1:nval)
2631  nsze = nval
2632  CLOSE(inputf)
2633  RETURN
2634  ELSE
2635  scan_file3 = -3
2636  END IF
2637  ELSE
2638  scan_file3 = -4
2639  END IF
2640  ELSE IF(PRESENT(fvec))THEN
2641  IF(nval > 1)THEN
2642  IF(vartype == 'float')THEN
2643  fvec(1:nval) = realval(1:nval)
2644  nsze = nval
2645  CLOSE(inputf)
2646  RETURN
2647  ELSE
2648  scan_file3 = -3
2649  END IF
2650  ELSE
2651  scan_file3 = -4
2652  END IF
2653  ELSE IF(PRESENT(cvec))THEN
2654  IF(nval > 0)THEN
2655  IF(vartype == 'string')THEN
2656  cvec(1:nval) = stringval(2:nval+1)
2657  nsze = nval
2658  CLOSE(inputf)
2659  RETURN
2660  ELSE
2661  scan_file3 = -3
2662  END IF
2663  ELSE
2664  scan_file3 = -4
2665  END IF
2666  ELSE
2667  scan_file3 = -5
2668  END IF
2669  END IF !!VARIABLE IS CORRECT
2670 
2671  END DO !!LOOP OVER INPUT FILE
2672  20 CLOSE(inputf)
2673  scan_file3 = -2
2674 
2675  RETURN
integer inputf
Definition: swmod1.f90:516
Here is the call graph for this function:
Here is the caller graph for this function:

◆ shutdown_check_1d()

subroutine mod_utils::shutdown_check_1d ( real(sp), dimension(:), intent(in)  VAR,
character(len=*), optional  MSG 
)

Definition at line 289 of file mod_utils.f90.

289 
290  !==============================================================================|
291  USE control
292  IMPLICIT NONE
293  REAL(DP) :: SBUF,RBUF
294  INTEGER :: IERR,I,J
295  REAL(SP), DIMENSION(:), INTENT(IN) :: VAR
296  CHARACTER(LEN=*), OPTIONAL :: MSG
297  !==============================================================================|
298 
299  IF (dbg_set(dbg_sbr)) THEN
300  IF (PRESENT(msg)) THEN
301  WRITE(ipt,*) "START: SHUTDOWN CHECK: "//msg
302  ELSE
303  WRITE(ipt,*) "START: SHUTDOWN CHECK: no msg"
304  END IF
305  END IF
306 
307  !Collect Depth Average to Master Processor
308  sbuf = sum(dble(var(1:ubound(var,1))))
309  rbuf = sbuf
310 
311  !Halt FVCOM if Depth Average = NaN
312  IF(isnan(rbuf))THEN
313  IF (PRESENT(msg)) THEN
314  CALL fatal_error("SHUTDOWN_CHECK FOUND NON FINITE VALUE:",&
315  & msg )
316  ELSE
317  CALL fatal_error('NON FINITE VALUE (DEPTH?) FOUND',&
318  & 'MODEL HAS BECOME UNSTABLE')
319  END IF
320 
321  END IF
322 
323 
324  IF (dbg_set(dbg_sbr)) WRITE(ipt,*) "END: SHUTDOWN CHECK"
325 
326  RETURN
integer ipt
Definition: mod_main.f90:922
Here is the call graph for this function:

◆ shutdown_check_2d()

subroutine mod_utils::shutdown_check_2d ( real(sp), dimension(:,:), intent(in)  VAR,
character(len=*), optional  MSG 
)

Definition at line 330 of file mod_utils.f90.

330 
331  !==============================================================================|
332  USE control
333  IMPLICIT NONE
334  REAL(DP) :: SBUF,RBUF
335  INTEGER :: IERR,I,J
336  REAL(SP), DIMENSION(:,:),INTENT(IN) :: VAR
337  CHARACTER(LEN=*), OPTIONAL :: MSG
338  !==============================================================================|
339 
340  IF (dbg_set(dbg_sbr)) THEN
341  IF (PRESENT(msg)) THEN
342  WRITE(ipt,*) "START: SHUTDOWN CHECK: "//msg
343  ELSE
344  WRITE(ipt,*) "START: SHUTDOWN CHECK: no msg"
345  END IF
346  END IF
347 
348  !Collect Depth Average to Master Processor
349  sbuf = sum(sum(dble(var(1:ubound(var,1),:)),1),1)
350  rbuf = sbuf
351 
352  !Halt FVCOM if Depth Average = NaN
353  IF(isnan(rbuf))THEN
354  IF (PRESENT(msg)) THEN
355  CALL fatal_error("SHUTDOWN_CHECK FOUND NON FINITE VALUE:",&
356  & msg )
357  ELSE
358  CALL fatal_error('NON FINITE VALUE (DEPTH?) FOUND',&
359  & 'MODEL HAS BECOME UNSTABLE')
360  END IF
361 
362  END IF
363 
364 
365  IF (dbg_set(dbg_sbr)) WRITE(ipt,*) "END: SHUTDOWN CHECK"
366 
367  RETURN
integer ipt
Definition: mod_main.f90:922
Here is the call graph for this function:

◆ split_string()

subroutine mod_utils::split_string ( character(len=*), intent(in)  instring,
character, intent(in)  delim,
character(len=*), dimension(:), intent(out), allocatable  outstrings 
)

Definition at line 2050 of file mod_utils.f90.

2050  IMPLICIT NONE
2051  character(len=*), intent(in) :: instring
2052  character, intent(in) :: delim
2053  character(len=*), intent(OUT), ALLOCATABLE :: outstrings(:)
2054  integer :: nlen, i, cnt, prev,next, idx, outlen, lgn
2055 
2056 ! character(len=len(outstrings)), ALLOCATABLE :: out_temp(:)
2057  character(len=len(instring)), ALLOCATABLE :: out_temp(:)
2058 
2059 
2060  ! Get the length of the string
2061  lgn = len_trim(instring)
2062  outlen = len(outstrings)
2063 
2064  ! CHECK FOR DEGERNERATE CASE (EMPTY STRING!)
2065  IF(lgn==0) THEN
2066  ALLOCATE(outstrings(1))
2067  outstrings=""
2068  RETURN
2069  END IF
2070 
2071 
2072  ! Count the number of seperations
2073  cnt = 0
2074  do i = 1,lgn
2075  if(instring(i:i) == delim) cnt=cnt+1
2076  end do
2077  ! If the string is not terminated, count the last entry too...
2078  if(instring(lgn:lgn) /= delim) cnt=cnt+1
2079 
2080  ! Allocate space
2081  ALLOCATE(outstrings(cnt))
2082 
2083  ! Split the string
2084  prev=1
2085  next=0
2086  DO i = 1,cnt
2087 ! write(*,*) "*** '"//trim(instring(prev:lgn))//"'"
2088 
2089  ! Find the first seperation
2090  idx = index(instring(prev:lgn),delim)
2091  if(idx==0) then
2092  ! IF none found, use end of string
2093  idx = lgn+1
2094  else
2095  ! Get the index into the real string length
2096  idx =idx + prev-1
2097  end if
2098 
2099  ! Set that last value to take
2100  next = idx-1
2101 ! write(*,*) I, prev, idx, next
2102 
2103  if(outlen .le. next-prev) Call warning&
2104  ("Insufficent room to split string!")
2105 
2106  ! Copy it into the string array
2107  outstrings(i) = trim(adjustl(instring(prev:next)))
2108 
2109  if(outlen .le. next-prev) Call warning&
2110  ("Insufficent room to split string!","'"//trim(outstrings(i))//"'")
2111 
2112 
2113  ! Set the first character of the next string
2114  prev=idx+1
2115 
2116 ! write(*,*) "! '"//trim(strings(I))//"'"
2117  END DO
2118 
2119  ! REMOVE DEGENERATE CASE FOR 'space' delimiter
2120 
2121  IF(delim == ' ') THEN
2122 
2123  idx = 0
2124  ALLOCATE(out_temp(cnt))
2125  DO i = 1,cnt
2126 
2127  IF(len_trim(outstrings(i)) > 0) THEN
2128  idx = idx + 1
2129  out_temp(idx) = outstrings(i)
2130  END IF
2131 
2132  END DO
2133 
2134  DEALLOCATE(outstrings)
2135  allocate(outstrings(idx))
2136  outstrings=out_temp
2137 
2138 
2139  END IF
2140 
Here is the call graph for this function:
Here is the caller graph for this function:

◆ test_split_strings()

subroutine mod_utils::test_split_strings ( )

Definition at line 2189 of file mod_utils.f90.

2189  implicit none
2190  CHARACTER(len=200) TESTIN
2191  CHARACTER(len=50),allocatable::testout(:)
2192  integer I
2193 
2194 ! TESTIN = "Hello world"
2195 ! TESTIN = "Hello world, Hello world2"
2196 ! TESTIN = "Hello worldmore that fiftymore that fiftymore that fiftymore that fiftymore that fiftymore that fiftymore that fiftymore that fifty,, Hello world2,"
2197  testin = "Hello world, Hello world2, Hello world2, Hello world2, Hello world2"
2198  testin = "Hello world"//achar(10)//" Hello world2"//achar(10)//" Hello world2, Hello world2, Hello world2"
2199 
2200  call split_string(testin,achar(10),testout)
2201 
2202  write(ipt,*) "! "
2203  write(ipt,*) "! TESTING SPLIT STRINGS"
2204  write(ipt,*) "! "
2205 
2206  do i=1,size(testout)
2207  write(ipt,*) "! ",i,"'"//trim(testout(i))//"'"
2208  end do
2209 
integer ipt
Definition: mod_main.f90:922
Here is the call graph for this function:

◆ warning()

subroutine mod_utils::warning ( character(len=*)  ER1,
character(len=*), optional  ER2,
character(len=*), optional  ER3,
character(len=*), optional  ER4 
)

Definition at line 251 of file mod_utils.f90.

251  USE control, only: prg_name
252  implicit none
253  character(Len=*) :: ER1
254  CHARACTER(LEN=*), OPTIONAL :: ER2
255  CHARACTER(LEN=*), OPTIONAL :: ER3
256  CHARACTER(LEN=*), OPTIONAL :: ER4
257 
258  if(dbg_set(dbg_log)) then
259  write(ipt,*) "+++++++++++++++++++++++++++++++++++++++++++++++++"
260  write(ipt,*) trim(prg_name)//" WARNING!"
261  write(ipt,*) er1
262  IF(PRESENT(er2)) WRITE(ipt,*) er2
263  IF(PRESENT(er3)) WRITE(ipt,*) er3
264  IF(PRESENT(er4)) WRITE(ipt,*) er4
265  write(ipt,*) trim(prg_name)//" CONTINUEING"
266  write(ipt,*) "+++++++++++++++++++++++++++++++++++++++++++++++++"
267  end if
268 
character(len=80) prg_name
Definition: mod_main.f90:105
Here is the call graph for this function:
Here is the caller graph for this function:

◆ write_banner()

subroutine mod_utils::write_banner ( logical, intent(in)  PAR,
integer, intent(in)  NP,
integer, intent(in)  ID 
)

Definition at line 1534 of file mod_utils.f90.

1534  IMPLICIT NONE
1535  LOGICAL, INTENT(IN) :: PAR
1536  INTEGER, INTENT(IN) :: NP,ID
1537  character(len=8) :: chr_np,chr_id
1538 
1539 ! CREATED USING: http://www.mudmagic.com/figlet-server/
1540 ! 146) rounded.flf :
1541 
1542 WRITE(ipt,*)'!================================================================!'
1543 WRITE(ipt,*)' _______ _ _ _______ _______ _______ ______ _____ '
1544 WRITE(ipt,*)' (_______)(_) (_)(_______)(_______)(_______)(_____ \ (_____) '
1545 WRITE(ipt,*)' _____ _ _ _ _ _ _ _ _ _____) ) _ __ _ '
1546 WRITE(ipt,*)' | ___) | | | || | | | | || ||_|| |(_____ ( | |/ /| |'
1547 WRITE(ipt,*)' | | \ \ / / | |_____ | |___| || | | | _____) )_| /_| |'
1548 WRITE(ipt,*)' |_| \___/ \______) \_____/ |_| |_|(______/(_)\_____/ '
1549 WRITE(ipt,*)' -- Beta Release'
1550 WRITE(ipt,*)'!================================================================!'
1551 WRITE(ipt,*)'! !'
1552 WRITE(ipt,*)'!========DOMAIN DECOMPOSITION USING: METIS 4.0.1 ================!'
1553 WRITE(ipt,*)'!======Copyright 1998, Regents of University of Minnesota========!'
1554 WRITE(ipt,*)'! !'
1555 IF(par) THEN
1556  WRITE(chr_np,'(I3.3)') np
1557  WRITE(chr_id,'(I3.3)') id
1558 WRITE(ipt,*)'!================================================================!'
1559 WRITE(ipt,*)'! !'
1560 WRITE(ipt,*)'! RUNNING IN PARALLEL: '//trim(chr_np)//' Processors !'
1561 WRITE(ipt,*)'! MYID is '//trim(chr_id)//' !'
1562 WRITE(ipt,*)'!================================================================!'
1563 END IF
1564  RETURN
Here is the caller graph for this function:

Variable Documentation

◆ dbg_io

integer, parameter mod_utils::dbg_io =1

Definition at line 66 of file mod_utils.f90.

66  integer,parameter::dbg_io=1 ! [enm] IO Filenames

◆ dbg_log

integer, parameter mod_utils::dbg_log =0

Definition at line 65 of file mod_utils.f90.

65  integer,parameter::dbg_log=0 ! [enm] Production mode. Debugging is

◆ dbg_lvl

integer mod_utils::dbg_lvl

Definition at line 75 of file mod_utils.f90.

75  integer :: dbg_lvl ! [enm] Debugging level, initialized in input

◆ dbg_mpi

integer, parameter mod_utils::dbg_mpi =3

Definition at line 68 of file mod_utils.f90.

68  integer,parameter::dbg_mpi=3 ! [enm] Mpi Communication

◆ dbg_nbr

integer, parameter mod_utils::dbg_nbr =9

Definition at line 63 of file mod_utils.f90.

63  integer,parameter::dbg_nbr=9 ! [nbr] Number of different debugging levels

◆ dbg_par

logical mod_utils::dbg_par =.false.

Definition at line 76 of file mod_utils.f90.

76  logical :: dbg_par=.false.

◆ dbg_sbr

integer, parameter mod_utils::dbg_sbr =4

Definition at line 69 of file mod_utils.f90.

69  integer,parameter::dbg_sbr=4 ! [enm] Subroutine names on entry and exit

◆ dbg_sbrio

integer, parameter mod_utils::dbg_sbrio =5

Definition at line 70 of file mod_utils.f90.

70  integer,parameter::dbg_sbrio=5 ! [enm] Subroutine I/O

◆ dbg_scl

integer, parameter mod_utils::dbg_scl =2

Definition at line 67 of file mod_utils.f90.

67  integer,parameter::dbg_scl=2 ! [enm] Scalars

◆ dbg_vec

integer, parameter mod_utils::dbg_vec =6

Definition at line 71 of file mod_utils.f90.

71  integer,parameter::dbg_vec=6 ! [enm] Entire vectors

◆ dbg_vrb

integer, parameter mod_utils::dbg_vrb =7

Definition at line 72 of file mod_utils.f90.

72  integer,parameter::dbg_vrb=7 ! [enm] Everything

◆ ext_code

integer, parameter mod_utils::ext_code =3602

Definition at line 50 of file mod_utils.f90.

50  INTEGER, PARAMETER :: EXT_CODE =3602

◆ init_code

integer mod_utils::init_code

Definition at line 55 of file mod_utils.f90.

55  INTEGER :: INIT_CODE

◆ initnest_code

integer mod_utils::initnest_code

Definition at line 57 of file mod_utils.f90.

57  INTEGER :: INITNEST_CODE

◆ nc_code

integer mod_utils::nc_code

Definition at line 53 of file mod_utils.f90.

53  INTEGER :: NC_CODE

◆ ncav_code

integer mod_utils::ncav_code

Definition at line 54 of file mod_utils.f90.

54  INTEGER :: NCAV_CODE

◆ nesting_code

integer mod_utils::nesting_code

Definition at line 56 of file mod_utils.f90.

56  INTEGER :: NESTING_CODE

◆ prg_nm

character(len=80) mod_utils::prg_nm

Definition at line 74 of file mod_utils.f90.

74  CHARACTER(LEN=80) :: prg_nm ! used in mod_sng and mod_input

◆ restart_code

integer mod_utils::restart_code

Definition at line 52 of file mod_utils.f90.

52  INTEGER :: RESTART_CODE

◆ sync_tag

integer, parameter mod_utils::sync_tag =3601

Definition at line 49 of file mod_utils.f90.

49  INTEGER, PARAMETER :: SYNC_TAG =3601

◆ wait_code

integer, parameter mod_utils::wait_code =3603

Definition at line 51 of file mod_utils.f90.

51  INTEGER, PARAMETER :: WAIT_CODE =3603