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.
230 type(time_type),
intent(in) :: time
231 type(ocean_grid_type),
intent(inout) :: g
232 type(param_file_type),
intent(in) :: param_file
233 type(tidal_forcing_cs),
pointer :: cs
236 real,
dimension(SZI_(G), SZJ_(G)) :: &
240 real,
dimension(MAX_CONSTITUENTS) :: freq_def, phase0_def, amp_def, love_def
241 integer,
dimension(3) :: tide_ref_date
243 logical :: use_m2, use_s2, use_n2, use_k2, use_k1, use_o1, use_p1, use_q1
244 logical :: use_mf, use_mm
246 logical :: fail_if_missing = .true.
248 #include "version_variable.h" 249 character(len=40) :: mdl =
"MOM_tidal_forcing" 250 character(len=128) :: mesg
251 character(len=200) :: tidal_input_files(4*max_constituents)
252 integer :: i, j, c, is, ie, js, je, isd, ied, jsd, jed, nc
253 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
254 isd = g%isd ; ied = g%ied ; jsd = g%jsd; jed = g%jed
256 if (
associated(cs))
then 257 call mom_error(warning,
"tidal_forcing_init called with an associated "// &
258 "control structure.")
263 call log_version(param_file, mdl, version,
"")
264 call get_param(param_file, mdl,
"TIDES", tides, &
265 "If true, apply tidal momentum forcing.", default=.false.)
267 if (.not.tides)
return 273 allocate(cs%sin_struct(isd:ied,jsd:jed,3)) ; cs%sin_struct(:,:,:) = 0.0
274 allocate(cs%cos_struct(isd:ied,jsd:jed,3)) ; cs%cos_struct(:,:,:) = 0.0
275 deg_to_rad = 4.0*atan(1.0)/180.0
276 do j=js-1,je+1 ;
do i=is-1,ie+1
277 lat_rad(i,j) = g%geoLatT(i,j)*deg_to_rad
278 lon_rad(i,j) = g%geoLonT(i,j)*deg_to_rad
280 do j=js-1,je+1 ;
do i=is-1,ie+1
281 cs%sin_struct(i,j,1) = -sin(2.0*lat_rad(i,j)) * sin(lon_rad(i,j))
282 cs%cos_struct(i,j,1) = sin(2.0*lat_rad(i,j)) * cos(lon_rad(i,j))
283 cs%sin_struct(i,j,2) = -cos(lat_rad(i,j))**2 * sin(2.0*lon_rad(i,j))
284 cs%cos_struct(i,j,2) = cos(lat_rad(i,j))**2 * cos(2.0*lon_rad(i,j))
285 cs%sin_struct(i,j,3) = 0.0
286 cs%cos_struct(i,j,3) = (0.5-1.5*sin(lat_rad(i,j))**2)
289 call get_param(param_file, mdl,
"TIDE_M2", use_m2, &
290 "If true, apply tidal momentum forcing at the M2 "//&
291 "frequency. This is only used if TIDES is true.", &
293 call get_param(param_file, mdl,
"TIDE_S2", use_s2, &
294 "If true, apply tidal momentum forcing at the S2 "//&
295 "frequency. This is only used if TIDES is true.", &
297 call get_param(param_file, mdl,
"TIDE_N2", use_n2, &
298 "If true, apply tidal momentum forcing at the N2 "//&
299 "frequency. This is only used if TIDES is true.", &
301 call get_param(param_file, mdl,
"TIDE_K2", use_k2, &
302 "If true, apply tidal momentum forcing at the K2 "//&
303 "frequency. This is only used if TIDES is true.", &
305 call get_param(param_file, mdl,
"TIDE_K1", use_k1, &
306 "If true, apply tidal momentum forcing at the K1 "//&
307 "frequency. This is only used if TIDES is true.", &
309 call get_param(param_file, mdl,
"TIDE_O1", use_o1, &
310 "If true, apply tidal momentum forcing at the O1 "//&
311 "frequency. This is only used if TIDES is true.", &
313 call get_param(param_file, mdl,
"TIDE_P1", use_p1, &
314 "If true, apply tidal momentum forcing at the P1 "//&
315 "frequency. This is only used if TIDES is true.", &
317 call get_param(param_file, mdl,
"TIDE_Q1", use_q1, &
318 "If true, apply tidal momentum forcing at the Q1 "//&
319 "frequency. This is only used if TIDES is true.", &
321 call get_param(param_file, mdl,
"TIDE_MF", use_mf, &
322 "If true, apply tidal momentum forcing at the MF "//&
323 "frequency. This is only used if TIDES is true.", &
325 call get_param(param_file, mdl,
"TIDE_MM", use_mm, &
326 "If true, apply tidal momentum forcing at the MM "//&
327 "frequency. This is only used if TIDES is true.", &
332 if (use_m2) nc=nc+1 ;
if (use_s2) nc=nc+1
333 if (use_n2) nc=nc+1 ;
if (use_k2) nc=nc+1
334 if (use_k1) nc=nc+1 ;
if (use_o1) nc=nc+1
335 if (use_p1) nc=nc+1 ;
if (use_q1) nc=nc+1
336 if (use_mf) nc=nc+1 ;
if (use_mm) nc=nc+1
340 call mom_error(fatal,
"tidal_forcing_init: "// &
341 "TIDES are defined, but no tidal constituents are used.")
345 call get_param(param_file, mdl,
"TIDAL_SAL_FROM_FILE", cs%tidal_sal_from_file, &
346 "If true, read the tidal self-attraction and loading "//&
347 "from input files, specified by TIDAL_INPUT_FILE. "//&
348 "This is only used if TIDES is true.", default=.false.)
349 call get_param(param_file, mdl,
"USE_PREVIOUS_TIDES", cs%use_prev_tides, &
350 "If true, use the SAL from the previous iteration of the "//&
351 "tides to facilitate convergent iteration. "//&
352 "This is only used if TIDES is true.", default=.false.)
353 call get_param(param_file, mdl,
"TIDE_USE_SAL_SCALAR", cs%use_sal_scalar, &
354 "If true and TIDES is true, use the scalar approximation "//&
355 "when calculating self-attraction and loading.", &
356 default=.not.cs%tidal_sal_from_file)
358 if (cs%use_sal_scalar .or. cs%use_prev_tides) &
359 call get_param(param_file, mdl,
"TIDE_SAL_SCALAR_VALUE", cs%sal_scalar, &
360 "The constant of proportionality between sea surface "//&
361 "height (really it should be bottom pressure) anomalies "//&
362 "and bottom geopotential anomalies. This is only used if "//&
363 "TIDES and TIDE_USE_SAL_SCALAR are true.", units=
"m m-1", &
364 fail_if_missing=.true.)
366 if (nc > max_constituents)
then 367 write(mesg,
'("Increase MAX_CONSTITUENTS in MOM_tidal_forcing.F90 to at least",I3, & 368 &"to accommodate all the registered tidal constituents.")') nc
369 call mom_error(fatal,
"MOM_tidal_forcing"//mesg)
372 do c=1,4*max_constituents ; tidal_input_files(c) =
"" ;
enddo 374 if (cs%tidal_sal_from_file .or. cs%use_prev_tides)
then 375 call get_param(param_file, mdl,
"TIDAL_INPUT_FILE", tidal_input_files, &
376 "A list of input files for tidal information.", &
377 default =
"", fail_if_missing=.true.)
380 call get_param(param_file, mdl,
"TIDE_REF_DATE", tide_ref_date, &
381 "Year,month,day to use as reference date for tidal forcing. "//&
382 "If not specified, defaults to 0.", &
385 call get_param(param_file, mdl,
"TIDE_USE_EQ_PHASE", cs%use_eq_phase, &
386 "Correct phases by calculating equilibrium phase arguments for TIDE_REF_DATE. ", &
387 default=.false., fail_if_missing=.false.)
389 if (sum(tide_ref_date) == 0)
then 390 cs%time_ref = set_date(1, 1, 1)
392 if(.not. cs%use_eq_phase)
then 396 call mom_mesg(
'Tidal phases will *not* be corrected with equilibrium arguments.')
398 cs%time_ref = set_date(tide_ref_date(1), tide_ref_date(2), tide_ref_date(3))
403 if (cs%use_eq_phase)
call astro_longitudes_init(cs%time_ref, cs%tidal_longitudes)
406 c=c+1 ; cs%const_name(c) =
"M2" ; cs%struct(c) = 2
407 cs%love_no(c) = 0.693 ; cs%amp(c) = 0.242334
411 c=c+1 ; cs%const_name(c) =
"S2" ; cs%struct(c) = 2
412 cs%love_no(c) = 0.693 ; cs%amp(c) = 0.112743
416 c=c+1 ; cs%const_name(c) =
"N2" ; cs%struct(c) = 2
417 cs%love_no(c) = 0.693 ; cs%amp(c) = 0.046397
421 c=c+1 ; cs%const_name(c) =
"K2" ; cs%struct(c) = 2
422 cs%love_no(c) = 0.693 ; cs%amp(c) = 0.030684
426 c=c+1 ; cs%const_name(c) =
"K1" ; cs%struct(c) = 1
427 cs%love_no(c) = 0.736 ; cs%amp(c) = 0.141565
431 c=c+1 ; cs%const_name(c) =
"O1" ; cs%struct(c) = 1
432 cs%love_no(c) = 0.695 ; cs%amp(c) = 0.100661
436 c=c+1 ; cs%const_name(c) =
"P1" ; cs%struct(c) = 1
437 cs%love_no(c) = 0.706 ; cs%amp(c) = 0.046848
441 c=c+1 ; cs%const_name(c) =
"Q1" ; cs%struct(c) = 1
442 cs%love_no(c) = 0.695 ; cs%amp(c) = 0.019273
446 c=c+1 ; cs%const_name(c) =
"MF" ; cs%struct(c) = 3
447 cs%love_no(c) = 0.693 ; cs%amp(c) = 0.042041
451 c=c+1 ; cs%const_name(c) =
"MM" ; cs%struct(c) = 3
452 cs%love_no(c) = 0.693 ; cs%amp(c) = 0.022191
458 cs%freq(c) = tidal_frequency(cs%const_name(c))
459 freq_def(c) = cs%freq(c)
460 love_def(c) = cs%love_no(c)
461 amp_def(c) = cs%amp(c)
463 if (cs%use_eq_phase)
then 464 phase0_def(c) = eq_phase(cs%const_name(c), cs%tidal_longitudes)
474 call get_param(param_file, mdl,
"TIDE_"//trim(cs%const_name(c))//
"_FREQ", cs%freq(c), &
475 "Frequency of the "//trim(cs%const_name(c))//
" tidal constituent. "//&
476 "This is only used if TIDES and TIDE_"//trim(cs%const_name(c))// &
477 " are true, or if OBC_TIDE_N_CONSTITUENTS > 0 and "//trim(cs%const_name(c))// &
478 " is in OBC_TIDE_CONSTITUENTS.", units=
"s-1", default=freq_def(c))
479 call get_param(param_file, mdl,
"TIDE_"//trim(cs%const_name(c))//
"_AMP", cs%amp(c), &
480 "Amplitude of the "//trim(cs%const_name(c))//
" tidal constituent. "//&
481 "This is only used if TIDES and TIDE_"//trim(cs%const_name(c))// &
482 " are true.", units=
"m", default=amp_def(c))
483 call get_param(param_file, mdl,
"TIDE_"//trim(cs%const_name(c))//
"_PHASE_T0", cs%phase0(c), &
484 "Phase of the "//trim(cs%const_name(c))//
" tidal constituent at time 0. "//&
485 "This is only used if TIDES and TIDE_"//trim(cs%const_name(c))// &
486 " are true.", units=
"radians", default=phase0_def(c))
489 if (cs%tidal_sal_from_file)
then 490 allocate(cs%cosphasesal(isd:ied,jsd:jed,nc))
491 allocate(cs%sinphasesal(isd:ied,jsd:jed,nc))
492 allocate(cs%ampsal(isd:ied,jsd:jed,nc))
495 call find_in_files(tidal_input_files,
"PHASE_SAL_"//trim(cs%const_name(c)),phase,g)
496 call find_in_files(tidal_input_files,
"AMP_SAL_"//trim(cs%const_name(c)),cs%ampsal(:,:,c),g)
497 call pass_var(phase, g%domain,complete=.false.)
498 call pass_var(cs%ampsal(:,:,c),g%domain,complete=.true.)
499 do j=js-1,je+1 ;
do i=is-1,ie+1
500 cs%cosphasesal(i,j,c) = cos(phase(i,j)*deg_to_rad)
501 cs%sinphasesal(i,j,c) = sin(phase(i,j)*deg_to_rad)
506 if (cs%USE_PREV_TIDES)
then 507 allocate(cs%cosphase_prev(isd:ied,jsd:jed,nc))
508 allocate(cs%sinphase_prev(isd:ied,jsd:jed,nc))
509 allocate(cs%amp_prev(isd:ied,jsd:jed,nc))
512 call find_in_files(tidal_input_files,
"PHASE_PREV_"//trim(cs%const_name(c)),phase,g)
513 call find_in_files(tidal_input_files,
"AMP_PREV_"//trim(cs%const_name(c)),cs%amp_prev(:,:,c),g)
514 call pass_var(phase, g%domain,complete=.false.)
515 call pass_var(cs%amp_prev(:,:,c),g%domain,complete=.true.)
516 do j=js-1,je+1 ;
do i=is-1,ie+1
517 cs%cosphase_prev(i,j,c) = cos(phase(i,j)*deg_to_rad)
518 cs%sinphase_prev(i,j,c) = sin(phase(i,j)*deg_to_rad)
523 id_clock_tides = cpu_clock_id(
'(Ocean tides)', grain=clock_module)