MOM6
RGC_initialization.F90
1 module rgc_initialization
2 !***********************************************************************
3 !* GNU General Public License *
4 !* This file is a part of MOM. *
5 !* *
6 !* MOM is free software; you can redistribute it and/or modify it and *
7 !* are expected to follow the terms of the GNU General Public License *
8 !* as published by the Free Software Foundation; either version 2 of *
9 !* the License, or (at your option) any later version. *
10 !* *
11 !* MOM is distributed in the hope that it will be useful, but WITHOUT *
12 !* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY *
13 !* or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public *
14 !* License for more details. *
15 !* *
16 !* For the full text of the GNU General Public License, *
17 !* write to: Free Software Foundation, Inc., *
18 !* 675 Mass Ave, Cambridge, MA 02139, USA. *
19 !* or see: http://www.gnu.org/licenses/gpl.html *
20 !* By Elizabeth Yankovsky, May 2018 *
21 !***********************************************************************
22 
25 use mom_domains, only : pass_var
27 use mom_error_handler, only : mom_mesg, mom_error, fatal, is_root_pe, warning
29 use mom_get_input, only : directories
30 use mom_grid, only : ocean_grid_type
31 use mom_io, only : file_exists, read_data
32 use mom_io, only : slasher
33 use mom_sponge, only : sponge_cs, set_up_sponge_field, initialize_sponge
34 use mom_sponge, only : set_up_sponge_ml_density
39 implicit none ; private
40 
41 #include <MOM_memory.h>
42 
43 character(len=40) :: mod = "RGC_initialization" ! This module's name.
44 public rgc_initialize_sponges
45 
46 contains
47 
48 !> Sets up the the inverse restoration time, and the values towards which the interface heights,
49 !! velocities and tracers should be restored within the sponges for the RGC test case.
50 subroutine rgc_initialize_sponges(G, GV, US, tv, u, v, PF, use_ALE, CSp, ACSp)
51  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure.
52  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure.
53  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
54  type(thermo_var_ptrs), intent(in) :: tv !< A structure containing pointers
55  !! to any available thermodynamic
56  !! fields, potential temperature and
57  !! salinity or mixed layer density.
58  !! Absent fields have NULL ptrs.
59  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
60  target, intent(in) :: u !< Array with the u velocity [L T-1 ~> m s-1]
61  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
62  target, intent(in) :: v !< Array with the v velocity [L T-1 ~> m s-1]
63  type(param_file_type), intent(in) :: pf !< A structure indicating the
64  !! open file to parse for model
65  !! parameter values.
66  logical, intent(in) :: use_ale !< If true, indicates model is in ALE mode
67  type(sponge_cs), pointer :: csp !< Layer-mode sponge structure
68  type(ale_sponge_cs), pointer :: acsp !< ALE-mode sponge structure
69 
70 ! Local variables
71  real :: t(szi_(g),szj_(g),szk_(g)) ! A temporary array for temp
72  real :: s(szi_(g),szj_(g),szk_(g)) ! A temporary array for salt
73  real :: u1(szib_(g),szj_(g),szk_(g)) ! A temporary array for u [L T-1 ~> m s-1]
74  real :: v1(szi_(g),szjb_(g),szk_(g)) ! A temporary array for v [L T-1 ~> m s-1]
75  real :: rho(szi_(g),szj_(g),szk_(g)) ! A temporary array for RHO
76  real :: tmp(szi_(g),szj_(g)) ! A temporary array for tracers.
77  real :: h(szi_(g),szj_(g),szk_(g)) ! A temporary array for thickness at h points
78  real :: idamp(szi_(g),szj_(g)) ! The inverse damping rate at h points [T-1 ~> s-1].
79  real :: tnudg ! Nudging time scale [T ~> s]
80  real :: pres(szi_(g)) ! An array of the reference pressure [R L2 T-2 ~> Pa]
81  real :: e0(szk_(g)+1) ! The resting interface heights, in m, usually !
82  ! negative because it is positive upward. !
83  real :: eta(szi_(g),szj_(g),szk_(g)+1) ! A temporary array for eta.
84  ! positive upward, in m.
85  logical :: sponge_uv ! Nudge velocities (u and v) towards zero
86  real :: min_depth, dummy1, z, delta_h
87  real :: rho_dummy, min_thickness, rho_tmp, xi0
88  real :: lenlat, lenlon, lensponge
89  character(len=40) :: filename, state_file
90  character(len=40) :: temp_var, salt_var, eta_var, inputdir, h_var
91 
92  character(len=40) :: mod = "RGC_initialize_sponges" ! This subroutine's name.
93  integer, dimension(2) :: eosdom ! The i-computational domain for the equation of state
94  integer :: i, j, k, is, ie, js, je, isd, ied, jsd, jed, nz, iscb, iecb, jscb, jecb
95 
96  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
97  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
98  iscb = g%iscB ; iecb = g%iecB; jscb = g%jscB ; jecb = g%jecB
99 
100  call get_param(pf,mod,"MIN_THICKNESS",min_thickness,'Minimum layer thickness',units='m',default=1.e-3)
101 
102  call get_param(pf, mod, "RGC_TNUDG", tnudg, 'Nudging time scale for sponge layers (days)', &
103  default=0.0, scale=86400.0*us%s_to_T)
104 
105  call get_param(pf, mod, "LENLAT", lenlat, &
106  "The latitudinal or y-direction length of the domain", &
107  fail_if_missing=.true., do_not_log=.true.)
108 
109  call get_param(pf, mod, "LENLON", lenlon, &
110  "The longitudinal or x-direction length of the domain", &
111  fail_if_missing=.true., do_not_log=.true.)
112 
113  call get_param(pf, mod, "LENSPONGE", lensponge, &
114  "The length of the sponge layer (km).", &
115  default=10.0)
116 
117  call get_param(pf, mod, "SPONGE_UV", sponge_uv, &
118  "Nudge velocities (u and v) towards zero in the sponge layer.", &
119  default=.false., do_not_log=.true.)
120 
121  t(:,:,:) = 0.0 ; s(:,:,:) = 0.0 ; idamp(:,:) = 0.0; rho(:,:,:) = 0.0
122 
123  call get_param(pf, mod, "MINIMUM_DEPTH", min_depth, &
124  "The minimum depth of the ocean.", units="m", default=0.0)
125 
126  if (associated(csp)) call mom_error(fatal, &
127  "RGC_initialize_sponges called with an associated control structure.")
128  if (associated(acsp)) call mom_error(fatal, &
129  "RGC_initialize_sponges called with an associated ALE-sponge control structure.")
130 
131  ! Here the inverse damping time [T-1 ~> s-1], is set. Set Idamp to 0
132  ! wherever there is no sponge, and the subroutines that are called
133  ! will automatically set up the sponges only where Idamp is positive
134  ! and mask2dT is 1.
135 
136  do i=is,ie ; do j=js,je
137  if ((g%bathyT(i,j) <= min_depth) .or. (g%geoLonT(i,j) <= lensponge)) then
138  idamp(i,j) = 0.0
139  elseif (g%geoLonT(i,j) >= (lenlon - lensponge) .AND. g%geoLonT(i,j) <= lenlon) then
140  dummy1 = (g%geoLonT(i,j)-(lenlon - lensponge))/(lensponge)
141  idamp(i,j) = (1.0/tnudg) * max(0.0,dummy1)
142  else
143  idamp(i,j) = 0.0
144  endif
145  enddo ; enddo
146 
147 
148  ! 1) Read eta, salt and temp from IC file
149  call get_param(pf, mod, "INPUTDIR", inputdir, default=".")
150  inputdir = slasher(inputdir)
151  ! GM: get two different files, one with temp and one with salt values
152  ! this is work around to avoid having wrong values near the surface
153  ! because of the FIT_SALINITY option. To get salt values right in the
154  ! sponge, FIT_SALINITY=False. The oposite is true for temp. One can
155  ! combined the *correct* temp and salt values in one file instead.
156  call get_param(pf, mod, "RGC_SPONGE_FILE", state_file, &
157  "The name of the file with temps., salts. and interfaces to \n"// &
158  " damp toward.", fail_if_missing=.true.)
159  call get_param(pf, mod, "SPONGE_PTEMP_VAR", temp_var, &
160  "The name of the potential temperature variable in \n"//&
161  "SPONGE_STATE_FILE.", default="Temp")
162  call get_param(pf, mod, "SPONGE_SALT_VAR", salt_var, &
163  "The name of the salinity variable in \n"//&
164  "SPONGE_STATE_FILE.", default="Salt")
165  call get_param(pf, mod, "SPONGE_ETA_VAR", eta_var, &
166  "The name of the interface height variable in \n"//&
167  "SPONGE_STATE_FILE.", default="eta")
168  call get_param(pf, mod, "SPONGE_H_VAR", h_var, &
169  "The name of the layer thickness variable in \n"//&
170  "SPONGE_STATE_FILE.", default="h")
171 
172  !read temp and eta
173  filename = trim(inputdir)//trim(state_file)
174  if (.not.file_exists(filename, g%Domain)) &
175  call mom_error(fatal, " RGC_initialize_sponges: Unable to open "//trim(filename))
176  call read_data(filename,temp_var,t(:,:,:), domain=g%Domain%mpp_domain)
177  call read_data(filename,salt_var,s(:,:,:), domain=g%Domain%mpp_domain)
178 
179  if (use_ale) then
180 
181  call read_data(filename,h_var,h(:,:,:), domain=g%Domain%mpp_domain)
182  call pass_var(h, g%domain)
183 
184  !call initialize_ALE_sponge(Idamp, h, nz, G, PF, ACSp)
185  call initialize_ale_sponge(idamp, g, pf, acsp, h, nz)
186 
187  ! The remaining calls to set_up_sponge_field can be in any order. !
188  if ( associated(tv%T) ) then
189  call set_up_ale_sponge_field(t,g,tv%T,acsp)
190  endif
191  if ( associated(tv%S) ) then
192  call set_up_ale_sponge_field(s,g,tv%S,acsp)
193  endif
194 
195  if (sponge_uv) then
196  u1(:,:,:) = 0.0; v1(:,:,:) = 0.0
197  call set_up_ale_sponge_vel_field(u1,v1,g,u,v,acsp)
198  endif
199 
200 
201  else ! layer mode
202 
203  !read eta
204  call read_data(filename,eta_var,eta(:,:,:), domain=g%Domain%mpp_domain)
205 
206  ! Set the inverse damping rates so that the model will know where to
207  ! apply the sponges, along with the interface heights.
208  call initialize_sponge(idamp, eta, g, pf, csp, gv)
209 
210  if ( gv%nkml>0 ) then
211  ! This call to set_up_sponge_ML_density registers the target values of the
212  ! mixed layer density, which is used in determining which layers can be
213  ! inflated without causing static instabilities.
214  do i=is-1,ie ; pres(i) = tv%P_Ref ; enddo
215  eosdom(:) = eos_domain(g%HI)
216  do j=js,je
217  call calculate_density(t(:,j,1), s(:,j,1), pres, tmp(:,j), tv%eqn_of_state, eosdom)
218  enddo
219 
220  call set_up_sponge_ml_density(tmp, g, csp)
221  endif
222 
223  ! Apply sponge in tracer fields
224  call set_up_sponge_field(t, tv%T, g, nz, csp)
225  call set_up_sponge_field(s, tv%S, g, nz, csp)
226 
227  endif
228 
229 end subroutine rgc_initialize_sponges
230 
231 end module rgc_initialization
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:26
Calculates density of sea water from T, S and P.
Definition: MOM_EOS.F90:68
Reads the only Fortran name list needed to boot-strap the model.
A structure that can be parsed to read and document run-time parameters.
Provides the ocean grid type.
Definition: MOM_grid.F90:2
Ddetermine the number of points which are within sponges in this computational domain.
Describes the horizontal ocean grid with only dynamic memory arrays.
Store the reference profile at h points for a variable.
This module contains I/O framework code.
Definition: MOM_io.F90:2
The MOM6 facility to parse input files for runtime parameters.
This module contains the routines used to apply sponge layers when using the ALE mode.
ALE sponge control structure.
Describes various unit conversion factors.
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
Describes the decomposed MOM domain and has routines for communications across PEs.
Definition: MOM_domains.F90:2
Routines for error handling and I/O management.
This control structure holds memory and parameters for the MOM_sponge module.
Definition: MOM_sponge.F90:41
Implements sponge regions in isopycnal mode.
Definition: MOM_sponge.F90:2
Calculate the derivatives of density with temperature and salinity from T, S, and P...
Definition: MOM_EOS.F90:81
Provides subroutines for quantities specific to the equation of state.
Definition: MOM_EOS.F90:2
An overloaded interface to log version information about modules.
Describes the vertical ocean grid, including unit conversion factors.
Indicate whether a file exists, perhaps with domain decomposition.
Definition: MOM_io.F90:68
This subroutine stores the reference profile at u and v points for a vector.
Container for paths and parameter file names.
Pointers to an assortment of thermodynamic fields that may be available, including potential temperat...
Contains a shareable dynamic type for describing horizontal grids and metric data and utilty routines...
Provides a transparent vertical ocean grid type and supporting routines.
Provides transparent structures with groups of MOM6 variables and supporting routines.
Do a halo update on an array.
Definition: MOM_domains.F90:54
An overloaded interface to read and log the values of various types of parameters.
A control structure for the equation of state.
Definition: MOM_EOS.F90:108