My Project
Functions/Subroutines | Variables
mod_setup Module Reference

Functions/Subroutines

subroutine setup_center_coords
 
subroutine setup_horizontal_mixing_coefficient
 
subroutine setup_bottom_roughness
 
subroutine setup_depth
 
subroutine setup_coriolis
 
subroutine setup_gravity
 
subroutine setup_sponge
 
subroutine coordinate_units (XL, YL)
 
subroutine setup_sigma
 
subroutine setup_sigma_derivatives
 

Variables

real(sp), dimension(:), allocatable, target x_gbl
 
real(sp), dimension(:), allocatable, target y_gbl
 
real(sp), dimension(:), allocatable, target x_lcl
 
real(sp), dimension(:), allocatable, target y_lcl
 
real(sp), dimension(:), allocatable, target c_lcl
 
real(sp), dimension(:), allocatable, target h_lcl
 
integer, dimension(:), allocatable, target n_spg
 
real(sp), dimension(:), allocatable, target r_spg
 
real(sp), dimension(:), allocatable, target c_spg
 
real(sp), dimension(:), allocatable, target x_spg
 
real(sp), dimension(:), allocatable, target y_spg
 
integer nsponge
 

Function/Subroutine Documentation

◆ coordinate_units()

subroutine mod_setup::coordinate_units ( real(sp), dimension(:), allocatable  XL,
real(sp), dimension(:), allocatable  YL 
)

Definition at line 300 of file mod_setup.f90.

300  IMPLICIT NONE
301  REAL(SP), ALLOCATABLE :: XL(:),YL(:)
302  integer status, ierr
303 
304 
305  REAL(SP), allocatable :: x_buf(:), lon_buf(:)
306  REAL(SP), allocatable :: y_buf(:), lat_buf(:)
307 
308  IF(dbg_set(dbg_sbr)) write(ipt,*) "COODINATE_UNITS: START"
309 
310 !! IF NOT SPHERICAL CASE
311 
312  IF (grid_file_units == 'degrees') THEN
313 
314  IF (use_proj) THEN
315 
316  lon = xl
317  lat = yl
318 
319  ! USE PROJECTION TOOL BOX TO CONVERT TO METERS
320  IF (serial) THEN
321  CALL degrees2meters(xl(1:mt),yl(1:mt), &
322  & projection_reference,vx(1:mt),vy(1:mt),mt)
323  END IF
324  xm = vx
325  ym = vy
326 
327  ELSE
328 
329  CALL fatal_error('You must specify a valid projection reference'&
330  &,'and compile with PROJ to use files with latitude and longitue in cartesian mode')
331  END IF
332 
333 
334  ELSE IF(grid_file_units == 'meters') THEN
335 
336  vx = xl
337  vy = yl
338 
339  xm = xl
340  ym = yl
341 
342  ! USE PROJECTION TOOL BOX TO CONVERT TO DEGREES
343  IF (use_proj ) THEN
344  IF (serial) THEN
345  CALL meters2degrees(xl(1:mt),yl(1:mt), &
347  END IF
348 
349  END IF
350 
351  ELSE
352  CALL fatal_error('UNRECOGNIZED GRID_FILE_UNITS: '//trim(grid_file_units))
353 
354  END IF
355 
356 
357  IF(dbg_set(dbg_sbr)) write(ipt,*) "COODINATE_UNITS: END"
358 
logical use_proj
Definition: mod_main.f90:633
logical serial
Definition: mod_main.f90:100
integer mt
Definition: mod_main.f90:78
logical function dbg_set(vrb)
Definition: mod_utils.f90:182
character(len=80) grid_file_units
Definition: mod_main.f90:626
character(len=200) projection_reference
Definition: mod_main.f90:627
real(sp), dimension(:), allocatable, target vx
Definition: mod_main.f90:1001
real(sp), dimension(:), allocatable, target vy
Definition: mod_main.f90:1002
real(sp), dimension(:), allocatable, target xm
Definition: mod_main.f90:991
real(sp), dimension(:), allocatable, target lat
Definition: mod_main.f90:996
subroutine fatal_error(ER1, ER2, ER3, ER4)
Definition: mod_utils.f90:230
real(sp), dimension(:), allocatable, target lon
Definition: mod_main.f90:995
integer ipt
Definition: mod_main.f90:922
integer, parameter dbg_sbr
Definition: mod_utils.f90:69
real(sp), dimension(:), allocatable, target ym
Definition: mod_main.f90:992
Here is the call graph for this function:
Here is the caller graph for this function:

◆ setup_bottom_roughness()

subroutine mod_setup::setup_bottom_roughness ( )

Definition at line 166 of file mod_setup.f90.

166  IMPLICIT NONE
167 
168  if (bottom_roughness_kind .eq. sttc) THEN
169 
170  if(dbg_set(dbg_log)) then
171  write(ipt,*) "! "
172  write(ipt,*) "! Setting Staticly Variable Bottom Roughness"
173  end if
174 
176 
177  else if(bottom_roughness_kind .eq. cnstnt) THEN
178 
180  cc_z0b(0) = 0.0_sp
181 
182  else
183  CALL fatal_error&
184  &("HORIZONTAL_MIXING_KIND ERROR",&
185  & "This should not happen")
186 
187  end if
188 
189 
character(len=80), parameter sttc
Definition: mod_main.f90:489
logical function dbg_set(vrb)
Definition: mod_utils.f90:182
real(sp), dimension(:), allocatable, target cc_z0b
Definition: mod_main.f90:1171
character(len=80), parameter cnstnt
Definition: mod_main.f90:488
subroutine load_bottom_roughness(Z0)
Definition: mod_input.f90:1747
subroutine fatal_error(ER1, ER2, ER3, ER4)
Definition: mod_utils.f90:230
integer ipt
Definition: mod_main.f90:922
real(sp) bottom_roughness_lengthscale
Definition: mod_main.f90:372
character(len=80) bottom_roughness_kind
Definition: mod_main.f90:368
integer, parameter dbg_log
Definition: mod_utils.f90:65
Here is the call graph for this function:
Here is the caller graph for this function:

◆ setup_center_coords()

subroutine mod_setup::setup_center_coords ( )

Definition at line 74 of file mod_setup.f90.

74  USE mod_spherical
75  IMPLICIT NONE
76  INTEGER I,J,IERR,STATUS,SENDID
77  REAL(SP) SBUF
78 
79  INTEGER K,ITMP
80  REAL(DP) VX1,VY1,VX2,VY2,VX3,VY3
81  REAL(DP) EVX12,EVX13,EVX23,EVY12,EVY13,EVY23,EVXY
82  REAL(DP) VX12,VY12,VX23,VY23,VX31,VY31
83  INTEGER :: SENDER
84 
85  ! David Added:
86  REAL(SP), allocatable :: xc_buf(:), lonc_buf(:)
87  REAL(SP), allocatable :: yc_buf(:), latc_buf(:)
88 
89  IF(dbg_set(dbg_sbr)) write(ipt,*) "SETUP_CENTER_COORDS: START"
90 
91 
92 !==============================================================================|
93 ! SET UP LOCAL MESH (HORIZONTAL COORDINATES) |
94 !==============================================================================|
95 !--------------CALCULATE GLOBAL MINIMUMS AND MAXIMUMS--------------------------!
96 
97 ! IF THIS IS NOT A SPHERICAL MODEL REMOVE MIN COORDINATE VALUE
98 
99  vxmin = minval(vx(1:mt)) ; vxmax = maxval(vx(1:mt))
100  vymin = minval(vy(1:mt)) ; vymax = maxval(vy(1:mt))
101 
102  !--------------SHIFT GRID TO UPPER RIGHT CARTESIAN-----------------------------!
103 
104  vx = vx - vxmin
105  vy = vy - vymin
106  vx(0) = 0.0_sp ; vy(0) = 0.0_sp
107 
108  !--------------CALCULATE GLOBAL ELEMENT CENTER GRID COORDINATES----------------!
109  CALL n2e2d(vx,xc)
110  CALL n2e2d(vy,yc)
111 
112  xmc = xc + vxmin
113  xmc(0)= 0.0_sp
114  ymc = yc + vymin
115  ymc(0)= 0.0_sp
116 
117  IF (use_proj ) THEN
118  IF (serial) THEN
119  CALL meters2degrees(xmc(1:nt),ymc(1:nt)&
121  END IF
122 
123  END IF
124 
125 
126 
127  IF(dbg_set(dbg_sbr)) write(ipt,*) "SETUP_CENTER_COORDS: END"
128 
logical use_proj
Definition: mod_main.f90:633
logical serial
Definition: mod_main.f90:100
real(sp) vymax
Definition: mod_main.f90:989
integer mt
Definition: mod_main.f90:78
logical function dbg_set(vrb)
Definition: mod_utils.f90:182
real(sp), dimension(:), allocatable, target yc
Definition: mod_main.f90:1004
real(sp), dimension(:), allocatable, target ymc
Definition: mod_main.f90:994
real(sp), dimension(:), allocatable, target latc
Definition: mod_main.f90:998
real(sp) vymin
Definition: mod_main.f90:989
character(len=200) projection_reference
Definition: mod_main.f90:627
real(sp), dimension(:), allocatable, target vx
Definition: mod_main.f90:1001
real(sp), dimension(:), allocatable, target vy
Definition: mod_main.f90:1002
real(sp), dimension(:), allocatable, target xmc
Definition: mod_main.f90:993
real(sp), dimension(:), allocatable, target lonc
Definition: mod_main.f90:997
subroutine n2e2d(NVAR, EVAR)
Definition: mod_main.f90:1390
real(sp) vxmax
Definition: mod_main.f90:989
real(sp), dimension(:), allocatable, target xc
Definition: mod_main.f90:1003
integer ipt
Definition: mod_main.f90:922
integer, parameter dbg_sbr
Definition: mod_utils.f90:69
integer nt
Definition: mod_main.f90:77
real(sp) vxmin
Definition: mod_main.f90:989
Here is the call graph for this function:
Here is the caller graph for this function:

◆ setup_coriolis()

subroutine mod_setup::setup_coriolis ( )

Definition at line 214 of file mod_setup.f90.

214  IMPLICIT NONE
215  integer:: I
216  !--------------READ IN CORIOLIS PARAMETER--------------------------------------!
217 
218  IF(dbg_set(dbg_sbr)) write(ipt,*) "SETUP_CORIOLIS: START"
219 
220  CALL n2e2d(c_lcl,cor)
221 
222 
223  cor = 2.*7.292e-5_sp * sin(cor * deg2rad)
224 
225  !! ggao for equatoral min (4deg)
226  IF(.NOT. equator_beta_plane)THEN
227  WHERE(cor < 1.e-5_sp .AND. cor > 0.0_sp) cor = 1.e-5_sp
228  WHERE(cor > -1.e-5_sp .AND. cor < 0.0_sp) cor = -1.e-5_sp
229  END IF
230 
231  IF(dbg_set(dbg_sbr)) write(ipt,*) "SETUP_CORIOLIS: END"
232 
logical equator_beta_plane
Definition: mod_main.f90:404
real(sp), dimension(:), allocatable, target cor
Definition: mod_main.f90:1113
logical function dbg_set(vrb)
Definition: mod_utils.f90:182
real(sp), dimension(:), allocatable, target c_lcl
Definition: mod_setup.f90:58
subroutine n2e2d(NVAR, EVAR)
Definition: mod_main.f90:1390
real(dp), parameter deg2rad
Definition: mod_main.f90:885
integer ipt
Definition: mod_main.f90:922
integer, parameter dbg_sbr
Definition: mod_utils.f90:69
Here is the call graph for this function:
Here is the caller graph for this function:

◆ setup_depth()

subroutine mod_setup::setup_depth ( )

Definition at line 195 of file mod_setup.f90.

195  IMPLICIT NONE
196  INTEGER :: IERR
197  REAL(SP):: SBUF
198 
199  IF(dbg_set(dbg_sbr)) write(ipt,*) "SETUP_DEPTH: START"
200  hmax = maxval(h_lcl(1:mt))
201  hmin = minval(h_lcl(1:mt))
202 
203 
204  h = h_lcl
205 
206  IF(dbg_set(dbg_sbr)) write(ipt,*) "SETUP_DEPTH: END"
207 
real(sp), dimension(:), allocatable, target h
Definition: mod_main.f90:1131
integer mt
Definition: mod_main.f90:78
real(sp) hmin
Definition: mod_main.f90:918
real(sp) hmax
Definition: mod_main.f90:917
logical function dbg_set(vrb)
Definition: mod_utils.f90:182
real(sp), dimension(:), allocatable, target h_lcl
Definition: mod_setup.f90:61
integer ipt
Definition: mod_main.f90:922
integer, parameter dbg_sbr
Definition: mod_utils.f90:69
Here is the call graph for this function:
Here is the caller graph for this function:

◆ setup_gravity()

subroutine mod_setup::setup_gravity ( )

Definition at line 241 of file mod_setup.f90.

241  IMPLICIT NONE
242  INTEGER I
243 
244  IF(dbg_set(dbg_sbr)) write(ipt,*) "SETUP_GRAVITY: START"
245 
246  grav_n = grav
247  grav_e = grav
248 
249  IF(dbg_set(dbg_sbr)) write(ipt,*) "SETUP_GRAVITY: END"
logical function dbg_set(vrb)
Definition: mod_utils.f90:182
real(dp), parameter grav
Definition: mod_main.f90:879
real(sp), dimension(:), allocatable, target grav_e
Definition: mod_main.f90:1013
real(sp), dimension(:), allocatable, target grav_n
Definition: mod_main.f90:1013
integer ipt
Definition: mod_main.f90:922
integer, parameter dbg_sbr
Definition: mod_utils.f90:69
Here is the call graph for this function:
Here is the caller graph for this function:

◆ setup_horizontal_mixing_coefficient()

subroutine mod_setup::setup_horizontal_mixing_coefficient ( )

Definition at line 134 of file mod_setup.f90.

134  IMPLICIT NONE
135 
136  if (horizontal_mixing_kind .eq. sttc) THEN
137 
138  if(dbg_set(dbg_log)) then
139  write(ipt,*) "! "
140  write(ipt,*) "! Setting Staticly Variable viscosity"
141  end if
142 
144 
145 
146  else if(horizontal_mixing_kind .eq. cnstnt) THEN
147 
149  cc_hvc(0)=0.0_sp
151  nn_hvc(0)=0.0_sp
152 
153  else
154  CALL fatal_error&
155  &("HORIZONTAL_MIXING_KIND ERROR",&
156  & "This should not happen")
157 
158  end if
159 
160 
character(len=80), parameter sttc
Definition: mod_main.f90:489
logical function dbg_set(vrb)
Definition: mod_utils.f90:182
subroutine load_horizontal_mixing_coefficient(NN, CC)
Definition: mod_input.f90:1690
real(sp), dimension(:), allocatable cc_hvc
Definition: mod_main.f90:1302
character(len=80) horizontal_mixing_kind
Definition: mod_main.f90:353
real(sp), dimension(:), allocatable nn_hvc
Definition: mod_main.f90:1303
character(len=80), parameter cnstnt
Definition: mod_main.f90:488
subroutine fatal_error(ER1, ER2, ER3, ER4)
Definition: mod_utils.f90:230
real(sp) horizontal_mixing_coefficient
Definition: mod_main.f90:354
integer ipt
Definition: mod_main.f90:922
integer, parameter dbg_log
Definition: mod_utils.f90:65
Here is the call graph for this function:
Here is the caller graph for this function:

◆ setup_sigma()

subroutine mod_setup::setup_sigma ( )

Definition at line 403 of file mod_setup.f90.

403  !==============================================================================|
404  IMPLICIT NONE
405  INTEGER :: K,KK
406  INTEGER :: I
407  REAL(SP):: ZTMP(KB)
408  REAL(SP):: X1,X2,X3
409  REAL(SP):: DR,RCL,RCU
410  !==============================================================================|
411 
412  IF(dbg_set(dbg_sbr)) WRITE(ipt,*)"SETUP_SIGMA: START"
413 
414  IF(dbg_set(dbg_sbrio)) THEN
415  WRITE(ipt,*)"==================="
416  WRITE(ipt,*)" SET_SIGMA IO"
417  WRITE(ipt,*)" STYPE = "//trim(stype)
418  WRITE(ipt,*)" P_SIGMA = ", p_sigma
419  WRITE(ipt,*)" KB = ", kb
420  WRITE(ipt,*)"==================="
421  END IF
422 
423  !---------------------------------------------
424  !---------------------------------------------
425  SELECT CASE(stype)
426  !---------------------------------------------
427  !SIGMA_COORDINATE_TYPE = UNIFORM (DEGENERATE CASE OF GEOMETRIC)
428  CASE(stype_uniform)
429  !---------------------------------------------
430  IF(p_sigma > 1 .AND. mod(kb,2) == 0) &
431  CALL fatal_error('SETUP_SIGMA: COORDINATE TYPE:'//trim(stype)&
432  &,'kb shoude be an odd number for this type of sigma coordinates....' )
433  CALL sigma_geometric
434  !---------------------------------------------
435  !SIGMA_COORDINATE_TYPE = GEOMETRIC
436  CASE(stype_geometric)
437  !---------------------------------------------
438 
439  IF(p_sigma > 1 .AND. mod(kb,2) == 0) &
440  CALL fatal_error('SETUP_SIGMA: COORDINATE TYPE:'//trim(stype)&
441  &,'kb shoude be an odd number for this type of sigma coordinates....' )
442  CALL sigma_geometric
443  !---------------------------------------------
444  ! SIGMA_COORDINATE_TYPE = TANH
445  CASE(stype_tanh)
446  !---------------------------------------------
447  CALL sigma_tanh
448  !---------------------------------------------
449  !SIGMA_COORDINATE_TYPE = GENERALIZED
450  CASE(stype_generalized)
451  !---------------------------------------------
452  CALL sigma_generalized
453 
454  ! THIS IS A CURRENTLY UNUSED METHOD
455 !!$!---------------------------------------------
456 !!$!SIGMA_COORDINATE_TYPE = WHAT THE HELL IS THIS?
457 !!$ CASE("UNKNOWN")
458 !!$!---------------------------------------------
459 !!$
460 !!$ CALL SIGMA_UNKNOWN
461  CASE DEFAULT
462  CALL fatal_error("SET_SIGMA: REACHED DEFAULT CASE FOR SIGMA COOR&
463  &DINATE TYPE")
464  END SELECT
465 
466  !---------COMPUTE SIGMA DERIVATIVES AND INTRA SIGMA LEVELS---------------------!
467 
468  IF(dbg_set(dbg_sbr)) WRITE(ipt,*)"SETUP_SIGMA: END"
469 
470  RETURN
471 
character(len=80), parameter stype_uniform
Definition: mod_main.f90:894
logical function dbg_set(vrb)
Definition: mod_utils.f90:182
integer kb
Definition: mod_main.f90:64
integer, parameter dbg_sbrio
Definition: mod_utils.f90:70
subroutine fatal_error(ER1, ER2, ER3, ER4)
Definition: mod_utils.f90:230
real(sp) p_sigma
Definition: mod_main.f90:900
character(len=80), parameter stype_geometric
Definition: mod_main.f90:895
character(len=80), parameter stype_tanh
Definition: mod_main.f90:896
integer ipt
Definition: mod_main.f90:922
integer, parameter dbg_sbr
Definition: mod_utils.f90:69
character(len=80) stype
Definition: mod_main.f90:893
character(len=80), parameter stype_generalized
Definition: mod_main.f90:897
Here is the call graph for this function:
Here is the caller graph for this function:

◆ setup_sigma_derivatives()

subroutine mod_setup::setup_sigma_derivatives ( )

Definition at line 636 of file mod_setup.f90.

636  USE all_vars
637  IMPLICIT NONE
638  INTEGER :: K, I
639 
640 
641  IF(dbg_set(dbg_sbr)) WRITE(ipt,*)"SETUP_SIGMA_DERIVATIVES: START"
642 
643  DO k=1,kb-1
644  DO i=1,mt
645  dz(i,k) = z(i,k)-z(i,k+1)
646  zz(i,k) = .5_sp*(z(i,k)+z(i,k+1))
647  END DO
648  DO i=1,nt
649  dz1(i,k) = z1(i,k)-z1(i,k+1)
650  zz1(i,k) = .5_sp*(z1(i,k)+z1(i,k+1))
651  END DO
652  END DO
653 
654  DO i=1,mt
655  zz(i,kb) = 2.0_sp*zz(i,kb-1)-zz(i,kb-2)
656  END DO
657  DO i=1,nt
658  zz1(i,kb) = 2.0_sp*zz1(i,kb-1)-zz1(i,kb-2)
659  END DO
660 
661  DO k=1,kbm2
662  DO i=1,mt
663  dzz(i,k) = zz(i,k)-zz(i,k+1)
664  END DO
665  DO i=1,nt
666  dzz1(i,k) = zz1(i,k)-zz1(i,k+1)
667  END DO
668  END DO
669 
670  dzz(:,kbm1) = 0.0_sp
671  dz(:,kb) = 0.0_sp
672  dzz1(:,kbm1) = 0.0_sp
673  dz1(:,kb) = 0.0_sp
674 
675 
676  !----------OUTPUT VALUES-TO INFOFILE-------------------------------------------!
677 
678  IF(dbg_set(dbg_log)) THEN
679  WRITE(ipt,* )'!'
680  WRITE(ipt,* )'!'
681  WRITE(ipt,*)'! SIGMA LAYER INFO '
682  WRITE(ipt,*) "SIGMA TYPE:",trim(stype)
683  WRITE(ipt,70)
684  SELECT CASE(stype)
685  CASE(stype_uniform)
686  DO k=1,kb
687  WRITE(ipt,80) k,z(1,k),zz(1,k),dz(1,k),dzz(1,k)
688  END DO
689  CASE(stype_restart)
690  DO k=1,kb
691  WRITE(ipt,80) k,z(1,k),zz(1,k),dz(1,k),dzz(1,k)
692  END DO
693  CASE(stype_geometric)
694  DO k=1,kb
695  WRITE(ipt,80) k,z(1,k),zz(1,k),dz(1,k),dzz(1,k)
696  END DO
697  CASE(stype_tanh)
698  DO k=1,kb
699  WRITE(ipt,80) k,z(1,k),zz(1,k),dz(1,k),dzz(1,k)
700  END DO
701  CASE(stype_generalized) ! THIS IS CASE SPECIFIC
702  WRITE(ipt,*)"SET CASE SPECIFIC GENERALIZED SIGMA LAYER OUTPUT TO SCREEN &
703  &IN mod_setup.F"
704  END SELECT
705  WRITE(ipt,* )'!'
706  END IF
707 
708  IF(dbg_set(dbg_sbr)) WRITE(ipt,*)"END SETUP_SIGMA_DERIVATIVES"
709 
710  !----------FORMAT STATEMENTS---------------------------------------------------!
711 
712 70 FORMAT(2x,'k',13x,'z',11x,'zz',11x,'dz',11x,'dzz')
713 80 FORMAT(' ',i5,4f13.8)
714 
715 
716  IF(dbg_set(dbg_sbr)) WRITE(ipt,*)"SETUP_SIGMA_DERIVATIVES: END"
integer mt
Definition: mod_main.f90:78
character(len=80), parameter stype_uniform
Definition: mod_main.f90:894
logical function dbg_set(vrb)
Definition: mod_utils.f90:182
real(sp), dimension(:,:), allocatable, target dzz1
Definition: mod_main.f90:1097
integer kb
Definition: mod_main.f90:64
integer kbm2
Definition: mod_main.f90:66
real(sp), dimension(:,:), allocatable, target zz1
Definition: mod_main.f90:1095
real(sp), dimension(:,:), allocatable, target dzz
Definition: mod_main.f90:1093
real(sp), dimension(:,:), allocatable, target dz
Definition: mod_main.f90:1092
real(sp), dimension(:,:), allocatable, target z
Definition: mod_main.f90:1090
real(sp), dimension(:,:), allocatable, target dz1
Definition: mod_main.f90:1096
real(sp), dimension(:,:), allocatable, target z1
Definition: mod_main.f90:1094
character(len=80), parameter stype_restart
Definition: mod_main.f90:898
character(len=80), parameter stype_geometric
Definition: mod_main.f90:895
character(len=80), parameter stype_tanh
Definition: mod_main.f90:896
integer ipt
Definition: mod_main.f90:922
integer, parameter dbg_sbr
Definition: mod_utils.f90:69
integer kbm1
Definition: mod_main.f90:65
real(sp), dimension(:,:), allocatable, target zz
Definition: mod_main.f90:1091
integer nt
Definition: mod_main.f90:77
character(len=80) stype
Definition: mod_main.f90:893
integer, parameter dbg_log
Definition: mod_utils.f90:65
character(len=80), parameter stype_generalized
Definition: mod_main.f90:897
Here is the call graph for this function:
Here is the caller graph for this function:

◆ setup_sponge()

subroutine mod_setup::setup_sponge ( )

Definition at line 257 of file mod_setup.f90.

257  USE mod_spherical
258  IMPLICIT NONE
259  REAL(SP) TEMP,DTMP,C_SPONGE
260  INTEGER :: I1, I, SENDER, IERR
261  REAL(DP) X1_DP,Y1_DP,X2_DP,Y2_DP,DTMP_DP
262 
263  !--SET SPONGE PARAMETERS-------------------------------------------------------|
264 
265  IF(dbg_set(dbg_sbr)) write(ipt,*) "SETUP_SPONGE: START"
266 
267  IF (nsponge ==0) RETURN
268 
269 ! NOTE: X_SPG/Y_SPG COORDINATES MUST BE AJUSTED FOR VXMIN/VYMIN
270 
271  x_spg = x_spg - vxmin
272  y_spg = y_spg - vymin
273 
274  DO i=1,nt
275  DO i1=1,nsponge
276  dtmp=(xc(i)-x_spg(i1))**2+(yc(i)-y_spg(i1))**2
277  dtmp=sqrt(dtmp)/r_spg(i1)
278 
279  IF(dtmp <= 1.0_sp) THEN
280  c_sponge=c_spg(i1)*(1.0_sp-dtmp)
281  cc_sponge(i)=max(c_sponge,cc_sponge(i))
282  END IF
283  END DO
284  END DO
285 
286 
287  DEALLOCATE(n_spg,r_spg,c_spg,x_spg,y_spg)
288 
289  IF(dbg_set(dbg_sbr)) write(ipt,*) "SETUP_SPONGE: END"
290 
291  RETURN
logical function dbg_set(vrb)
Definition: mod_utils.f90:182
integer nsponge
Definition: mod_setup.f90:66
real(sp), dimension(:), allocatable, target yc
Definition: mod_main.f90:1004
real(sp), dimension(:), allocatable, target c_spg
Definition: mod_setup.f90:65
real(sp), dimension(:), allocatable, target r_spg
Definition: mod_setup.f90:65
real(sp), dimension(:), allocatable, target y_spg
Definition: mod_setup.f90:65
real(sp) vymin
Definition: mod_main.f90:989
real(sp), dimension(:), allocatable, target x_spg
Definition: mod_setup.f90:65
integer, dimension(:), allocatable, target n_spg
Definition: mod_setup.f90:64
real(sp), dimension(:), allocatable, target xc
Definition: mod_main.f90:1003
real(sp), dimension(:), allocatable, target cc_sponge
Definition: mod_main.f90:1127
integer ipt
Definition: mod_main.f90:922
integer, parameter dbg_sbr
Definition: mod_utils.f90:69
integer nt
Definition: mod_main.f90:77
real(sp) vxmin
Definition: mod_main.f90:989
Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ c_lcl

real(sp), dimension(:), allocatable, target mod_setup::c_lcl

Definition at line 58 of file mod_setup.f90.

58  REAL(SP), ALLOCATABLE, TARGET :: C_LCL(:) !!CORIOLIS PARAMETER AT NODES

◆ c_spg

real(sp), dimension(:), allocatable, target mod_setup::c_spg

Definition at line 65 of file mod_setup.f90.

◆ h_lcl

real(sp), dimension(:), allocatable, target mod_setup::h_lcl

Definition at line 61 of file mod_setup.f90.

61  REAL(SP), ALLOCATABLE, TARGET :: H_LCL(:) !! DEPTH PARAMETER AT NODES

◆ n_spg

integer, dimension(:), allocatable, target mod_setup::n_spg

Definition at line 64 of file mod_setup.f90.

64  INTEGER, ALLOCATABLE, TARGET :: N_SPG(:)

◆ nsponge

integer mod_setup::nsponge

Definition at line 66 of file mod_setup.f90.

66  INTEGER :: NSPONGE

◆ r_spg

real(sp), dimension(:), allocatable, target mod_setup::r_spg

Definition at line 65 of file mod_setup.f90.

65  REAL(SP), ALLOCATABLE, TARGET :: R_SPG(:),C_SPG(:),X_SPG(:),Y_SPG(:)

◆ x_gbl

real(sp), dimension(:), allocatable, target mod_setup::x_gbl

Definition at line 54 of file mod_setup.f90.

54  REAL(SP), ALLOCATABLE, TARGET :: X_GBL(:),Y_GBL(:)

◆ x_lcl

real(sp), dimension(:), allocatable, target mod_setup::x_lcl

Definition at line 55 of file mod_setup.f90.

55  REAL(SP), ALLOCATABLE, TARGET :: X_LCL(:),Y_LCL(:)

◆ x_spg

real(sp), dimension(:), allocatable, target mod_setup::x_spg

Definition at line 65 of file mod_setup.f90.

◆ y_gbl

real(sp), dimension(:), allocatable, target mod_setup::y_gbl

Definition at line 54 of file mod_setup.f90.

◆ y_lcl

real(sp), dimension(:), allocatable, target mod_setup::y_lcl

Definition at line 55 of file mod_setup.f90.

◆ y_spg

real(sp), dimension(:), allocatable, target mod_setup::y_spg

Definition at line 65 of file mod_setup.f90.