1 module rgc_initialization
34 use mom_sponge,
only : set_up_sponge_ml_density
39 implicit none ;
private
41 #include <MOM_memory.h>
43 character(len=40) :: mod =
"RGC_initialization"
44 public rgc_initialize_sponges
50 subroutine rgc_initialize_sponges(G, GV, US, tv, u, v, PF, use_ALE, CSp, ACSp)
51 type(ocean_grid_type),
intent(in) :: G
52 type(verticalGrid_type),
intent(in) :: GV
53 type(unit_scale_type),
intent(in) :: US
54 type(thermo_var_ptrs),
intent(in) :: tv
59 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
60 target,
intent(in) :: u
61 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)), &
62 target,
intent(in) :: v
63 type(param_file_type),
intent(in) :: PF
66 logical,
intent(in) :: use_ALE
67 type(sponge_CS),
pointer :: CSp
68 type(ALE_sponge_CS),
pointer :: ACSp
71 real :: T(SZI_(G),SZJ_(G),SZK_(G))
72 real :: S(SZI_(G),SZJ_(G),SZK_(G))
73 real :: U1(SZIB_(G),SZJ_(G),SZK_(G))
74 real :: V1(SZI_(G),SZJB_(G),SZK_(G))
75 real :: RHO(SZI_(G),SZJ_(G),SZK_(G))
76 real :: tmp(SZI_(G),SZJ_(G))
77 real :: h(SZI_(G),SZJ_(G),SZK_(G))
78 real :: Idamp(SZI_(G),SZJ_(G))
83 real :: eta(SZI_(G),SZJ_(G),SZK_(G)+1)
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
92 character(len=40) :: mod =
"RGC_initialize_sponges"
93 integer,
dimension(2) :: EOSdom
94 integer :: i, j, k, is, ie, js, je, isd, ied, jsd, jed, nz, iscB, iecB, jscB, jecB
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
100 call get_param(pf,mod,
"MIN_THICKNESS",min_thickness,
'Minimum layer thickness',units=
'm',default=1.e-3)
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)
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.)
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.)
113 call get_param(pf, mod,
"LENSPONGE", lensponge, &
114 "The length of the sponge layer (km).", &
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.)
121 t(:,:,:) = 0.0 ; s(:,:,:) = 0.0 ; idamp(:,:) = 0.0; rho(:,:,:) = 0.0
123 call get_param(pf, mod,
"MINIMUM_DEPTH", min_depth, &
124 "The minimum depth of the ocean.", units=
"m", default=0.0)
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.")
136 do i=is,ie ;
do j=js,je
137 if ((g%bathyT(i,j) <= min_depth) .or. (g%geoLonT(i,j) <= lensponge))
then
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)
149 call get_param(pf, mod,
"INPUTDIR", inputdir, default=
".")
150 inputdir = slasher(inputdir)
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")
173 filename = trim(inputdir)//trim(state_file)
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)
181 call read_data(filename,h_var,h(:,:,:), domain=g%Domain%mpp_domain)
188 if (
associated(tv%T) )
then
191 if (
associated(tv%S) )
then
196 u1(:,:,:) = 0.0; v1(:,:,:) = 0.0
204 call read_data(filename,eta_var,eta(:,:,:), domain=g%Domain%mpp_domain)
208 call initialize_sponge(idamp, eta, g, pf, csp, gv)
210 if ( gv%nkml>0 )
then
214 do i=is-1,ie ; pres(i) = tv%P_Ref ;
enddo
215 eosdom(:) = eos_domain(g%HI)
217 call calculate_density(t(:,j,1), s(:,j,1), pres, tmp(:,j), tv%eqn_of_state, eosdom)
220 call set_up_sponge_ml_density(tmp, g, csp)
224 call set_up_sponge_field(t, tv%T, g, nz, csp)
225 call set_up_sponge_field(s, tv%S, g, nz, csp)
229 end subroutine rgc_initialize_sponges
231 end module rgc_initialization