This subroutine allocates space for the static variables used by this module. The metrics may be effectively 0, 1, or 2-D arrays, while fields like the background viscosities are 2-D arrays. ALLOC is a macro defined in MOM_memory.h for allocate or nothing with static memory.
113 type(time_type),
intent(in) :: time
115 type(param_file_type),
intent(in) :: param_file
117 type(tidal_forcing_cs),
pointer :: cs
132 real,
dimension(SZI_(G), SZJ_(G)) :: &
136 real,
dimension(MAX_CONSTITUENTS) :: freq_def, phase0_def, amp_def, love_def
138 logical :: use_m2, use_s2, use_n2, use_k2, use_k1, use_o1, use_p1, use_q1
139 logical :: use_mf, use_mm
141 logical :: fail_if_missing = .true.
143 #include "version_variable.h" 144 character(len=40) :: mdl =
"MOM_tidal_forcing" 145 character(len=128) :: mesg
146 character(len=200) :: tidal_input_files(4*max_constituents)
147 integer :: i, j, c, is, ie, js, je, isd, ied, jsd, jed, nc
148 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
149 isd = g%isd ; ied = g%ied ; jsd = g%jsd; jed = g%jed
151 if (
associated(cs))
then 152 call mom_error(warning,
"tidal_forcing_init called with an associated "// &
153 "control structure.")
158 call log_version(param_file, mdl, version,
"")
159 call get_param(param_file, mdl,
"TIDES", tides, &
160 "If true, apply tidal momentum forcing.", default=.false.)
162 if (.not.tides)
return 168 allocate(cs%sin_struct(isd:ied,jsd:jed,3)) ; cs%sin_struct(:,:,:) = 0.0
169 allocate(cs%cos_struct(isd:ied,jsd:jed,3)) ; cs%cos_struct(:,:,:) = 0.0
170 deg_to_rad = 4.0*atan(1.0)/180.0
171 do j=js-1,je+1 ;
do i=is-1,ie+1
172 lat_rad(i,j) = g%geoLatT(i,j)*deg_to_rad
173 lon_rad(i,j) = g%geoLonT(i,j)*deg_to_rad
175 do j=js-1,je+1 ;
do i=is-1,ie+1
176 cs%sin_struct(i,j,1) = -sin(2.0*lat_rad(i,j)) * sin(lon_rad(i,j))
177 cs%cos_struct(i,j,1) = sin(2.0*lat_rad(i,j)) * cos(lon_rad(i,j))
178 cs%sin_struct(i,j,2) = -cos(lat_rad(i,j))**2 * sin(2.0*lon_rad(i,j))
179 cs%cos_struct(i,j,2) = cos(lat_rad(i,j))**2 * cos(2.0*lon_rad(i,j))
180 cs%sin_struct(i,j,3) = 0.0
181 cs%cos_struct(i,j,3) = (0.5-1.5*sin(lat_rad(i,j))**2)
184 call get_param(param_file, mdl,
"TIDE_M2", use_m2, &
185 "If true, apply tidal momentum forcing at the M2 \n"//&
186 "frequency. This is only used if TIDES is true.", &
188 call get_param(param_file, mdl,
"TIDE_S2", use_s2, &
189 "If true, apply tidal momentum forcing at the S2 \n"//&
190 "frequency. This is only used if TIDES is true.", &
192 call get_param(param_file, mdl,
"TIDE_N2", use_n2, &
193 "If true, apply tidal momentum forcing at the N2 \n"//&
194 "frequency. This is only used if TIDES is true.", &
196 call get_param(param_file, mdl,
"TIDE_K2", use_k2, &
197 "If true, apply tidal momentum forcing at the K2 \n"//&
198 "frequency. This is only used if TIDES is true.", &
200 call get_param(param_file, mdl,
"TIDE_K1", use_k1, &
201 "If true, apply tidal momentum forcing at the K1 \n"//&
202 "frequency. This is only used if TIDES is true.", &
204 call get_param(param_file, mdl,
"TIDE_O1", use_o1, &
205 "If true, apply tidal momentum forcing at the O1 \n"//&
206 "frequency. This is only used if TIDES is true.", &
208 call get_param(param_file, mdl,
"TIDE_P1", use_p1, &
209 "If true, apply tidal momentum forcing at the P1 \n"//&
210 "frequency. This is only used if TIDES is true.", &
212 call get_param(param_file, mdl,
"TIDE_Q1", use_q1, &
213 "If true, apply tidal momentum forcing at the Q1 \n"//&
214 "frequency. This is only used if TIDES is true.", &
216 call get_param(param_file, mdl,
"TIDE_MF", use_mf, &
217 "If true, apply tidal momentum forcing at the MF \n"//&
218 "frequency. This is only used if TIDES is true.", &
220 call get_param(param_file, mdl,
"TIDE_MM", use_mm, &
221 "If true, apply tidal momentum forcing at the MM \n"//&
222 "frequency. This is only used if TIDES is true.", &
227 if (use_m2) nc=nc+1 ;
if (use_s2) nc=nc+1
228 if (use_n2) nc=nc+1 ;
if (use_k2) nc=nc+1
229 if (use_k1) nc=nc+1 ;
if (use_o1) nc=nc+1
230 if (use_p1) nc=nc+1 ;
if (use_q1) nc=nc+1
231 if (use_mf) nc=nc+1 ;
if (use_mm) nc=nc+1
235 call mom_error(fatal,
"tidal_forcing_init: "// &
236 "TIDES are defined, but no tidal constituents are used.")
240 call get_param(param_file, mdl,
"TIDAL_SAL_FROM_FILE", cs%tidal_sal_from_file, &
241 "If true, read the tidal self-attraction and loading \n"//&
242 "from input files, specified by TIDAL_INPUT_FILE. \n"//&
243 "This is only used if TIDES is true.", default=.false.)
244 call get_param(param_file, mdl,
"USE_PREVIOUS_TIDES", cs%use_prev_tides, &
245 "If true, use the SAL from the previous iteration of the \n"//&
246 "tides to facilitate convergent iteration. \n"//&
247 "This is only used if TIDES is true.", default=.false.)
248 call get_param(param_file, mdl,
"TIDE_USE_SAL_SCALAR", cs%use_sal_scalar, &
249 "If true and TIDES is true, use the scalar approximation \n"//&
250 "when calculating self-attraction and loading.", &
251 default=.not.cs%tidal_sal_from_file)
253 if (cs%use_sal_scalar .or. cs%use_prev_tides) &
254 call get_param(param_file, mdl,
"TIDE_SAL_SCALAR_VALUE", cs%sal_scalar, &
255 "The constant of proportionality between sea surface \n"//&
256 "height (really it should be bottom pressure) anomalies \n"//&
257 "and bottom geopotential anomalies. This is only used if \n"//&
258 "TIDES and TIDE_USE_SAL_SCALAR are true.", units=
"m m-1", &
259 fail_if_missing=.true.)
261 if (nc > max_constituents)
then 262 write(mesg,
'("Increase MAX_CONSTITUENTS in MOM_tidal_forcing.F90 to at least",I3, & 263 &"to accomodate all the registered tidal constituents.")') nc
264 call mom_error(fatal,
"MOM_tidal_forcing"//mesg)
267 do c=1,4*max_constituents ; tidal_input_files(c) =
"" ;
enddo 269 if (cs%tidal_sal_from_file .or. cs%use_prev_tides)
then 270 call get_param(param_file, mdl,
"TIDAL_INPUT_FILE", tidal_input_files, &
271 "A list of input files for tidal information.", &
272 default =
"", fail_if_missing=.true.)
278 c=c+1 ; cs%const_name(c) =
"M2" ; cs%freq(c) = 1.4051890e-4 ; cs%struct(c) = 2
279 cs%love_no(c) = 0.693 ; cs%amp(c) = 0.242334 ; cs%phase0(c) = 0.0
280 freq_def(c) = cs%freq(c) ; love_def(c) = cs%love_no(c)
281 amp_def(c) = cs%amp(c) ; phase0_def(c) = cs%phase0(c)
285 c=c+1 ; cs%const_name(c) =
"S2" ; cs%freq(c) = 1.4544410e-4 ; cs%struct(c) = 2
286 cs%love_no(c) = 0.693 ; cs%amp(c) = 0.112743 ; cs%phase0(c) = 0.0
287 freq_def(c) = cs%freq(c) ; love_def(c) = cs%love_no(c)
288 amp_def(c) = cs%amp(c) ; phase0_def(c) = cs%phase0(c)
292 c=c+1 ; cs%const_name(c) =
"N2" ; cs%freq(c) = 1.3787970e-4 ; cs%struct(c) = 2
293 cs%love_no(c) = 0.693 ; cs%amp(c) = 0.046397 ; cs%phase0(c) = 0.0
294 freq_def(c) = cs%freq(c) ; love_def(c) = cs%love_no(c)
295 amp_def(c) = cs%amp(c) ; phase0_def(c) = cs%phase0(c)
299 c=c+1 ; cs%const_name(c) =
"K2" ; cs%freq(c) = 1.4584234e-4 ; cs%struct(c) = 2
300 cs%love_no(c) = 0.693 ; cs%amp(c) = 0.030684 ; cs%phase0(c) = 0.0
301 freq_def(c) = cs%freq(c) ; love_def(c) = cs%love_no(c)
302 amp_def(c) = cs%amp(c) ; phase0_def(c) = cs%phase0(c)
306 c=c+1 ; cs%const_name(c) =
"K1" ; cs%freq(c) = 0.7292117e-4 ; cs%struct(c) = 1
307 cs%love_no(c) = 0.736 ; cs%amp(c) = 0.141565 ; cs%phase0(c) = 0.0
308 freq_def(c) = cs%freq(c) ; love_def(c) = cs%love_no(c)
309 amp_def(c) = cs%amp(c) ; phase0_def(c) = cs%phase0(c)
313 c=c+1 ; cs%const_name(c) =
"O1" ; cs%freq(c) = 0.6759774e-4 ; cs%struct(c) = 1
314 cs%love_no(c) = 0.695 ; cs%amp(c) = 0.100661 ; cs%phase0(c) = 0.0
315 freq_def(c) = cs%freq(c) ; love_def(c) = cs%love_no(c)
316 amp_def(c) = cs%amp(c) ; phase0_def(c) = cs%phase0(c)
320 c=c+1 ; cs%const_name(c) =
"P1" ; cs%freq(c) = 0.7252295e-4 ; cs%struct(c) = 1
321 cs%love_no(c) = 0.706 ; cs%amp(c) = 0.046848 ; cs%phase0(c) = 0.0
322 freq_def(c) = cs%freq(c) ; love_def(c) = cs%love_no(c)
323 amp_def(c) = cs%amp(c) ; phase0_def(c) = cs%phase0(c)
327 c=c+1 ; cs%const_name(c) =
"Q1" ; cs%freq(c) = 0.6495854e-4 ; cs%struct(c) = 1
328 cs%love_no(c) = 0.695 ; cs%amp(c) = 0.019273 ; cs%phase0(c) = 0.0
329 freq_def(c) = cs%freq(c) ; love_def(c) = cs%love_no(c)
330 amp_def(c) = cs%amp(c) ; phase0_def(c) = cs%phase0(c)
334 c=c+1 ; cs%const_name(c) =
"MF" ; cs%freq(c) = 0.053234e-4 ; cs%struct(c) = 3
335 cs%love_no(c) = 0.693 ; cs%amp(c) = 0.042041 ; cs%phase0(c) = 0.0
336 freq_def(c) = cs%freq(c) ; love_def(c) = cs%love_no(c)
337 amp_def(c) = cs%amp(c) ; phase0_def(c) = cs%phase0(c)
341 c=c+1 ; cs%const_name(c) =
"MM" ; cs%freq(c) = 0.026392e-4 ; cs%struct(c) = 3
342 cs%love_no(c) = 0.693 ; cs%amp(c) = 0.022191 ; cs%phase0(c) = 0.0
343 freq_def(c) = cs%freq(c) ; love_def(c) = cs%love_no(c)
344 amp_def(c) = cs%amp(c) ; phase0_def(c) = cs%phase0(c)
351 call get_param(param_file, mdl,
"TIDE_"//trim(cs%const_name(c))//
"_FREQ", cs%freq(c), &
352 "Frequency of the "//trim(cs%const_name(c))//
" tidal constituent. \n"//&
353 "This is only used if TIDES and TIDE_"//trim(cs%const_name(c))// &
354 " are true.", units=
"s-1", default=freq_def(c))
355 call get_param(param_file, mdl,
"TIDE_"//trim(cs%const_name(c))//
"_AMP", cs%amp(c), &
356 "Amplitude of the "//trim(cs%const_name(c))//
" tidal constituent. \n"//&
357 "This is only used if TIDES and TIDE_"//trim(cs%const_name(c))// &
358 " are true.", units=
"m", default=amp_def(c))
359 call get_param(param_file, mdl,
"TIDE_"//trim(cs%const_name(c))//
"_PHASE_T0", cs%phase0(c), &
360 "Phase of the "//trim(cs%const_name(c))//
" tidal constituent at time 0. \n"//&
361 "This is only used if TIDES and TIDE_"//trim(cs%const_name(c))// &
362 " are true.", units=
"radians", default=phase0_def(c))
365 if (cs%tidal_sal_from_file)
then 366 allocate(cs%cosphasesal(isd:ied,jsd:jed,nc))
367 allocate(cs%sinphasesal(isd:ied,jsd:jed,nc))
368 allocate(cs%ampsal(isd:ied,jsd:jed,nc))
371 call find_in_files(tidal_input_files,
"PHASE_SAL_"//trim(cs%const_name(c)),phase,g)
372 call find_in_files(tidal_input_files,
"AMP_SAL_"//trim(cs%const_name(c)),cs%ampsal(:,:,c),g)
373 call pass_var(phase, g%domain,complete=.false.)
374 call pass_var(cs%ampsal(:,:,c),g%domain,complete=.true.)
375 do j=js-1,je+1 ;
do i=is-1,ie+1
376 cs%cosphasesal(i,j,c) = cos(phase(i,j)*deg_to_rad)
377 cs%sinphasesal(i,j,c) = sin(phase(i,j)*deg_to_rad)
382 if (cs%USE_PREV_TIDES)
then 383 allocate(cs%cosphase_prev(isd:ied,jsd:jed,nc))
384 allocate(cs%sinphase_prev(isd:ied,jsd:jed,nc))
385 allocate(cs%amp_prev(isd:ied,jsd:jed,nc))
388 call find_in_files(tidal_input_files,
"PHASE_PREV_"//trim(cs%const_name(c)),phase,g)
389 call find_in_files(tidal_input_files,
"AMP_PREV_"//trim(cs%const_name(c)),cs%amp_prev(:,:,c),g)
390 call pass_var(phase, g%domain,complete=.false.)
391 call pass_var(cs%amp_prev(:,:,c),g%domain,complete=.true.)
392 do j=js-1,je+1 ;
do i=is-1,ie+1
393 cs%cosphase_prev(i,j,c) = cos(phase(i,j)*deg_to_rad)
394 cs%sinphase_prev(i,j,c) = sin(phase(i,j)*deg_to_rad)
399 id_clock_tides = cpu_clock_id(
'(Ocean tides)', grain=clock_module)
Ocean grid type. See mom_grid for details.