31 use mom_time_manager, only : time_type,
operator(+),
operator(/),
operator(-)
40 implicit none ;
private 42 #include <MOM_memory.h> 48 logical :: use_temperature
50 logical :: do_integrated
71 real,
pointer,
dimension(:) :: &
73 real,
pointer,
dimension(:,:) :: &
76 real,
pointer,
dimension(:,:,:) :: &
78 precip_cyc => null(), &
79 avg_sst_anom => null(), &
80 avg_sss_anom => null(), &
84 integer :: id_heat_0 = -1
94 real,
dimension(SZI_(G),SZJ_(G)),
intent(in) :: SST_anom
96 real,
dimension(SZI_(G),SZJ_(G)),
intent(in) :: SSS_anom
98 real,
dimension(SZI_(G),SZJ_(G)),
intent(in) :: SSS_mean
100 real,
dimension(SZI_(G),SZJ_(G)),
intent(inout) :: virt_heat
103 real,
dimension(SZI_(G),SZJ_(G)),
intent(inout) :: virt_precip
107 type(time_type),
intent(in) :: day_start
108 real,
intent(in) :: dt
114 real,
dimension(SZIB_(G),SZJ_(G)) :: &
117 real,
dimension(SZI_(G),SZJB_(G)) :: &
120 type(time_type) :: day_end
122 real :: mr_st, mr_end, mr_mid, mr_prev, mr_next
123 real :: dt_wt, dt_heat_rate, dt_prec_rate
124 real :: dt1_heat_rate, dt1_prec_rate, dt2_heat_rate, dt2_prec_rate
125 real :: wt_per1, wt_st, wt_end, wt_mid
126 integer :: m_st, m_end, m_mid, m_u1, m_u2, m_u3
127 integer :: yr, mon, day, hr, min, sec
128 integer :: i, j, is, ie, js, je
130 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
132 if (.not.
associated(cs))
return 133 if ((cs%num_cycle <= 0) .and. (.not.cs%do_integrated))
return 135 day_end = day_start + set_time(floor(dt+0.5))
137 do j=js,je ;
do i=is,ie
138 virt_heat(i,j) = 0.0 ; virt_precip(i,j) = 0.0
141 if (cs%do_integrated)
then 142 dt_heat_rate = dt * cs%heat_int_rate
143 dt_prec_rate = dt * cs%prec_int_rate
144 call pass_var(cs%heat_0, g%Domain, complete=.false.)
145 call pass_var(cs%precip_0, g%Domain)
147 do j=js,je ;
do i=is-1,ie
148 coef = cs%Len2 * (g%dy_Cu(i,j)*g%IdxCu(i,j))
149 flux_heat_x(i,j) = coef * (cs%heat_0(i,j) - cs%heat_0(i+1,j))
150 flux_prec_x(i,j) = coef * (cs%precip_0(i,j) - cs%precip_0(i+1,j))
152 do j=js-1,je ;
do i=is,ie
153 coef = cs%Len2 * (g%dx_Cv(i,j)*g%IdyCv(i,j))
154 flux_heat_y(i,j) = coef * (cs%heat_0(i,j) - cs%heat_0(i,j+1))
155 flux_prec_y(i,j) = coef * (cs%precip_0(i,j) - cs%precip_0(i,j+1))
157 do j=js,je ;
do i=is,ie
158 cs%heat_0(i,j) = cs%heat_0(i,j) + dt_heat_rate * ( &
159 -cs%lam_heat*g%mask2dT(i,j)*sst_anom(i,j) + &
160 (g%IareaT(i,j) * ((flux_heat_x(i-1,j) - flux_heat_x(i,j)) + &
161 (flux_heat_y(i,j-1) - flux_heat_y(i,j))) ) )
163 cs%precip_0(i,j) = cs%precip_0(i,j) + dt_prec_rate * ( &
164 cs%lam_prec * g%mask2dT(i,j)*(sss_anom(i,j) / sss_mean(i,j)) + &
165 (g%IareaT(i,j) * ((flux_prec_x(i-1,j) - flux_prec_x(i,j)) + &
166 (flux_prec_y(i,j-1) - flux_prec_y(i,j))) ) )
168 virt_heat(i,j) = virt_heat(i,j) + cs%heat_0(i,j)
169 virt_precip(i,j) = virt_precip(i,j) + cs%precip_0(i,j)
173 if (cs%num_cycle > 0)
then 175 call get_date(day_start, yr, mon, day, hr, min, sec)
176 mr_st = cs%num_cycle * (time_type_to_real(day_start - set_date(yr, 1, 1)) / &
177 time_type_to_real(set_date(yr+1, 1, 1) - set_date(yr, 1, 1)))
179 call get_date(day_end, yr, mon, day, hr, min, sec)
180 mr_end = cs%num_cycle * (time_type_to_real(day_end - set_date(yr, 1, 1)) / &
181 time_type_to_real(set_date(yr+1, 1, 1) - set_date(yr, 1, 1)))
189 m_end =
periodic_int(
real(ceiling(mr_end)), CS%num_cycle)
190 m_mid =
periodic_int(
real(ceiling(mr_st)), CS%num_cycle)
198 if (m_mid == m_end)
then ; mr_mid = mr_end
199 else ; mr_mid =
periodic_real(
real(m_mid), CS%num_cycle) ; endif
206 if (mr_st < mr_prev) mr_prev = mr_prev - cs%num_cycle
207 if (mr_mid < mr_st) mr_mid = mr_mid + cs%num_cycle
208 if (mr_end < mr_st) mr_end = mr_end + cs%num_cycle
209 if (mr_next < mr_prev) mr_next = mr_next + cs%num_cycle
212 if ((mr_mid < mr_st) .or. (mr_mid > mr_prev + 1.))
call mom_error(fatal, &
213 "apply ctrl_forcing: m_mid interpolation out of bounds; fix the code.")
214 if ((mr_end < mr_st) .or. (mr_end > mr_prev + 2.))
call mom_error(fatal, &
215 "apply ctrl_forcing: m_end interpolation out of bounds; fix the code.")
216 if (mr_end > mr_next)
call mom_error(fatal, &
217 "apply ctrl_forcing: mr_next interpolation out of bounds; fix the code.")
220 if (mr_mid < mr_end) wt_per1 = (mr_mid - mr_st) / (mr_end - mr_st)
223 wt_st = wt_per1 * (1. + (mr_prev - 0.5*(mr_st + mr_mid)))
224 wt_end = (1.0-wt_per1) * (1. + (0.5*(mr_end + mr_mid) - mr_next))
225 wt_mid = 1.0 - (wt_st + wt_end)
226 if ((wt_st < 0.0) .or. (wt_end < 0.0) .or. (wt_mid < 0.0)) &
227 call mom_error(fatal,
"apply_ctrl_forcing: Negative m weights")
228 if ((wt_st > 1.0) .or. (wt_end > 1.0) .or. (wt_mid > 1.0)) &
229 call mom_error(fatal,
"apply_ctrl_forcing: Excessive m weights")
232 do j=js,je ;
do i=is,ie
233 virt_heat(i,j) = virt_heat(i,j) + (wt_st * cs%heat_cyc(i,j,m_st) + &
234 (wt_mid * cs%heat_cyc(i,j,m_mid) + &
235 wt_end * cs%heat_cyc(i,j,m_end)))
236 virt_precip(i,j) = virt_precip(i,j) + (wt_st * cs%precip_cyc(i,j,m_st) + &
237 (wt_mid * cs%precip_cyc(i,j,m_mid) + &
238 wt_end * cs%precip_cyc(i,j,m_end)))
251 if (cs%avg_time(m_end) <= 0.0)
then 252 cs%avg_time(m_end) = 0.0
253 do j=js,je ;
do i=is,ie
254 cs%avg_SST_anom(i,j,m_end) = 0.0
255 cs%avg_SSS_anom(i,j,m_end) = 0.0 ; cs%avg_SSS(i,j,m_end) = 0.0
258 if (cs%avg_time(m_mid) <= 0.0)
then 259 cs%avg_time(m_mid) = 0.0
260 do j=js,je ;
do i=is,ie
261 cs%avg_SST_anom(i,j,m_mid) = 0.0
262 cs%avg_SSS_anom(i,j,m_mid) = 0.0 ; cs%avg_SSS(i,j,m_mid) = 0.0
268 cs%avg_time(m_mid) = cs%avg_time(m_mid) + dt_wt
269 do j=js,je ;
do i=is,ie
270 cs%avg_SST_anom(i,j,m_mid) = cs%avg_SST_anom(i,j,m_mid) + &
271 dt_wt * g%mask2dT(i,j) * sst_anom(i,j)
272 cs%avg_SSS_anom(i,j,m_mid) = cs%avg_SSS_anom(i,j,m_mid) + &
273 dt_wt * g%mask2dT(i,j) * sss_anom(i,j)
274 cs%avg_SSS(i,j,m_mid) = cs%avg_SSS(i,j,m_mid) + dt_wt * sss_mean(i,j)
276 if (wt_per1 < 1.0)
then 277 dt_wt = (1.0-wt_per1) * dt
278 cs%avg_time(m_end) = cs%avg_time(m_end) + dt_wt
279 do j=js,je ;
do i=is,ie
280 cs%avg_SST_anom(i,j,m_end) = cs%avg_SST_anom(i,j,m_end) + &
281 dt_wt * g%mask2dT(i,j) * sst_anom(i,j)
282 cs%avg_SSS_anom(i,j,m_end) = cs%avg_SSS_anom(i,j,m_end) + &
283 dt_wt * g%mask2dT(i,j) * sss_anom(i,j)
284 cs%avg_SSS(i,j,m_end) = cs%avg_SSS(i,j,m_end) + dt_wt * sss_mean(i,j)
293 if (cs%avg_time(m_u1) > 0.0)
then 294 do j=js,je ;
do i=is,ie
295 cs%avg_SST_anom(i,j,m_u1) = cs%avg_SST_anom(i,j,m_u1) / cs%avg_time(m_u1)
296 cs%avg_SSS_anom(i,j,m_u1) = cs%avg_SSS_anom(i,j,m_u1) / cs%avg_time(m_u1)
297 cs%avg_SSS(i,j,m_u1) = cs%avg_SSS(i,j,m_u1) / cs%avg_time(m_u1)
299 cs%avg_time(m_u1) = -1.0
301 if (cs%avg_time(m_u2) > 0.0)
then 302 do j=js,je ;
do i=is,ie
303 cs%avg_SST_anom(i,j,m_u2) = cs%avg_SST_anom(i,j,m_u2) / cs%avg_time(m_u2)
304 cs%avg_SSS_anom(i,j,m_u2) = cs%avg_SSS_anom(i,j,m_u2) / cs%avg_time(m_u2)
305 cs%avg_SSS(i,j,m_u2) = cs%avg_SSS(i,j,m_u2) / cs%avg_time(m_u2)
307 cs%avg_time(m_u2) = -1.0
309 if (cs%avg_time(m_u3) > 0.0)
then 310 do j=js,je ;
do i=is,ie
311 cs%avg_SST_anom(i,j,m_u3) = cs%avg_SST_anom(i,j,m_u3) / cs%avg_time(m_u3)
312 cs%avg_SSS_anom(i,j,m_u3) = cs%avg_SSS_anom(i,j,m_u3) / cs%avg_time(m_u3)
313 cs%avg_SSS(i,j,m_u3) = cs%avg_SSS(i,j,m_u3) / cs%avg_time(m_u3)
315 cs%avg_time(m_u3) = -1.0
318 dt1_heat_rate = wt_per1 * dt * cs%heat_cyc_rate
319 dt1_prec_rate = wt_per1 * dt * cs%prec_cyc_rate
320 dt2_heat_rate = (1.0-wt_per1) * dt * cs%heat_cyc_rate
321 dt2_prec_rate = (1.0-wt_per1) * dt * cs%prec_cyc_rate
323 if (wt_per1 < 1.0)
then 324 call pass_var(cs%heat_cyc(:,:,m_u2), g%Domain, complete=.false.)
325 call pass_var(cs%precip_cyc(:,:,m_u2), g%Domain, complete=.false.)
327 call pass_var(cs%heat_cyc(:,:,m_u1), g%Domain, complete=.false.)
328 call pass_var(cs%precip_cyc(:,:,m_u1), g%Domain)
330 if ((cs%avg_time(m_u1) == -1.0) .and. (cs%avg_time(m_u2) == -1.0))
then 331 do j=js,je ;
do i=is-1,ie
332 coef = cs%Len2 * (g%dy_Cu(i,j)*g%IdxCu(i,j))
333 flux_heat_x(i,j) = coef * (cs%heat_cyc(i,j,m_u1) - cs%heat_cyc(i+1,j,m_u1))
334 flux_prec_x(i,j) = coef * (cs%precip_cyc(i,j,m_u1) - cs%precip_cyc(i+1,j,m_u1))
336 do j=js-1,je ;
do i=is,ie
337 coef = cs%Len2 * (g%dx_Cv(i,j)*g%IdyCv(i,j))
338 flux_heat_y(i,j) = coef * (cs%heat_cyc(i,j,m_u1) - cs%heat_cyc(i,j+1,m_u1))
339 flux_prec_y(i,j) = coef * (cs%precip_cyc(i,j,m_u1) - cs%precip_cyc(i,j+1,m_u1))
341 do j=js,je ;
do i=is,ie
342 cs%heat_cyc(i,j,m_u1) = cs%heat_cyc(i,j,m_u1) + dt1_heat_rate * ( &
343 -cs%lam_cyc_heat*(cs%avg_SST_anom(i,j,m_u2) - cs%avg_SST_anom(i,j,m_u1)) + &
344 (g%IareaT(i,j) * ((flux_heat_x(i-1,j) - flux_heat_x(i,j)) + &
345 (flux_heat_y(i,j-1) - flux_heat_y(i,j))) ) )
347 cs%precip_cyc(i,j,m_u1) = cs%precip_cyc(i,j,m_u1) + dt1_prec_rate * ( &
348 cs%lam_cyc_prec * (cs%avg_SSS_anom(i,j,m_u2) - cs%avg_SSS_anom(i,j,m_u1)) / &
349 (0.5*(cs%avg_SSS(i,j,m_u2) + cs%avg_SSS(i,j,m_u1))) + &
350 (g%IareaT(i,j) * ((flux_prec_x(i-1,j) - flux_prec_x(i,j)) + &
351 (flux_prec_y(i,j-1) - flux_prec_y(i,j))) ) )
355 if ((wt_per1 < 1.0) .and. (cs%avg_time(m_u1) == -1.0) .and. (cs%avg_time(m_u2) == -1.0))
then 356 do j=js,je ;
do i=is-1,ie
357 coef = cs%Len2 * (g%dy_Cu(i,j)*g%IdxCu(i,j))
358 flux_heat_x(i,j) = coef * (cs%heat_cyc(i,j,m_u2) - cs%heat_cyc(i+1,j,m_u2))
359 flux_prec_x(i,j) = coef * (cs%precip_cyc(i,j,m_u2) - cs%precip_cyc(i+1,j,m_u2))
361 do j=js-1,je ;
do i=is,ie
362 coef = cs%Len2 * (g%dx_Cv(i,j)*g%IdyCv(i,j))
363 flux_heat_y(i,j) = coef * (cs%heat_cyc(i,j,m_u2) - cs%heat_cyc(i,j+1,m_u2))
364 flux_prec_y(i,j) = coef * (cs%precip_cyc(i,j,m_u2) - cs%precip_cyc(i,j+1,m_u2))
366 do j=js,je ;
do i=is,ie
367 cs%heat_cyc(i,j,m_u2) = cs%heat_cyc(i,j,m_u2) + dt1_heat_rate * ( &
368 -cs%lam_cyc_heat*(cs%avg_SST_anom(i,j,m_u3) - cs%avg_SST_anom(i,j,m_u2)) + &
369 (g%IareaT(i,j) * ((flux_heat_x(i-1,j) - flux_heat_x(i,j)) + &
370 (flux_heat_y(i,j-1) - flux_heat_y(i,j))) ) )
372 cs%precip_cyc(i,j,m_u2) = cs%precip_cyc(i,j,m_u2) + dt1_prec_rate * ( &
373 cs%lam_cyc_prec * (cs%avg_SSS_anom(i,j,m_u3) - cs%avg_SSS_anom(i,j,m_u2)) / &
374 (0.5*(cs%avg_SSS(i,j,m_u3) + cs%avg_SSS(i,j,m_u2))) + &
375 (g%IareaT(i,j) * ((flux_prec_x(i-1,j) - flux_prec_x(i,j)) + &
376 (flux_prec_y(i,j-1) - flux_prec_y(i,j))) ) )
386 real,
intent(in) :: rval
387 integer,
intent(in) :: num_period
392 m = m + num_period * (1 + (abs(m) / num_period))
393 elseif (m > num_period)
then 394 m = m - num_period * ((m-1) / num_period)
401 real,
intent(in) :: rval
402 integer,
intent(in) :: num_period
406 if (rval < 0)
then ; nshft = floor(abs(rval) / num_period) + 1
407 elseif (rval < num_period)
then ; nshft = 0
408 else ; nshft = -1*floor(rval / num_period) ;
endif 410 val_out = rval + nshft * num_period
425 logical :: controlled, use_temperature
426 character (len=8) :: period_str
428 integer :: isd, ied, jsd, jed, IsdB, IedB, JsdB, JedB
429 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
430 isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
432 if (
associated(cs))
then 433 call mom_error(warning,
"register_ctrl_forcing_restarts called "//&
434 "with an associated control structure.")
439 call read_param(param_file,
"CONTROLLED_FORCING", controlled)
440 if (.not.controlled)
return 442 use_temperature = .true.
443 call read_param(param_file,
"ENABLE_THERMODYNAMICS", use_temperature)
444 if (.not.use_temperature)
call mom_error(fatal, &
445 "register_ctrl_forcing_restarts: CONTROLLED_FORCING only works with "//&
446 "ENABLE_THERMODYNAMICS defined.")
450 cs%do_integrated = .true. ; cs%num_cycle = 0
451 call read_param(param_file,
"CTRL_FORCE_INTEGRATED", cs%do_integrated)
452 call read_param(param_file,
"CTRL_FORCE_NUM_CYCLE", cs%num_cycle)
454 if (cs%do_integrated)
then 455 call safe_alloc_ptr(cs%heat_0,isd,ied,jsd,jed) ; cs%heat_0(:,:) = 0.0
456 call safe_alloc_ptr(cs%precip_0,isd,ied,jsd,jed) ; cs%precip_0(:,:) = 0.0
457 vd =
var_desc(
"Ctrl_heat",
"W m-2",
"Control Integrative Heating",z_grid=
'1')
459 vd =
var_desc(
"Ctrl_precip",
"kg m-2 s-1",
"Control Integrative Precipitation",z_grid=
'1')
463 if (cs%num_cycle > 0)
then 464 write (period_str,
'(i8)') cs%num_cycle
465 period_str = trim(
'p ')//trim(adjustl(period_str))
466 call safe_alloc_ptr(cs%heat_cyc,isd,ied,jsd,jed,cs%num_cycle) ; cs%heat_cyc(:,:,:) = 0.0
467 call safe_alloc_ptr(cs%precip_cyc,isd,ied,jsd,jed,cs%num_cycle) ; cs%precip_cyc(:,:,:) = 0.0
468 vd =
var_desc(
"Ctrl_heat_cycle",
"W m-2",
"Cyclical Control Heating",&
469 z_grid=
'1', t_grid=period_str)
471 vd =
var_desc(
"Ctrl_precip_cycle",
"kg m-2 s-1",
"Cyclical Control Precipitation", &
472 z_grid=
'1', t_grid=period_str)
475 call safe_alloc_ptr(cs%avg_time,cs%num_cycle) ; cs%avg_time(:) = 0.0
476 vd =
var_desc(
"avg_time",
"sec",
"Cyclical accumulated averaging time", &
477 '1',z_grid=
'1',t_grid=period_str)
480 call safe_alloc_ptr(cs%avg_SST_anom,isd,ied,jsd,jed,cs%num_cycle) ; cs%avg_SST_anom(:,:,:) = 0.0
481 call safe_alloc_ptr(cs%avg_SSS_anom,isd,ied,jsd,jed,cs%num_cycle) ; cs%avg_SSS_anom(:,:,:) = 0.0
482 vd =
var_desc(
"avg_SST_anom",
"deg C",
"Cyclical average SST Anomaly", &
483 z_grid=
'1',t_grid=period_str)
485 vd =
var_desc(
"avg_SSS_anom",
"g kg-1",
"Cyclical average SSS Anomaly", &
486 z_grid=
'1',t_grid=period_str)
494 type(time_type),
intent(in) :: Time
499 type(
diag_ctrl),
target,
intent(in) :: diag
504 logical :: do_integrated
507 #include "version_variable.h" 508 character(len=40) :: mdl =
"MOM_controlled_forcing" 514 if (
associated(cs))
then 515 do_integrated = cs%do_integrated ; num_cycle = cs%num_cycle
517 do_integrated = .false. ; num_cycle = 0
522 call log_param(param_file, mdl,
"CTRL_FORCE_INTEGRATED", do_integrated, &
523 "If true, use a PI controller to determine the surface \n"//&
524 "forcing that is consistent with the observed mean properties.", &
526 call log_param(param_file, mdl,
"CTRL_FORCE_NUM_CYCLE", num_cycle, &
527 "The number of cycles per year in the controlled forcing, \n"//&
528 "or 0 for no cyclic forcing.", default=0)
530 if (.not.
associated(cs))
return 534 call get_param(param_file, mdl,
"CTRL_FORCE_HEAT_INT_RATE", cs%heat_int_rate, &
535 "The integrated rate at which heat flux anomalies are \n"//&
536 "accumulated.", units=
"s-1", default=0.0)
537 call get_param(param_file, mdl,
"CTRL_FORCE_PREC_INT_RATE", cs%prec_int_rate, &
538 "The integrated rate at which precipitation anomalies \n"//&
539 "are accumulated.", units=
"s-1", default=0.0)
540 call get_param(param_file, mdl,
"CTRL_FORCE_HEAT_CYC_RATE", cs%heat_cyc_rate, &
541 "The integrated rate at which cyclical heat flux \n"//&
542 "anomalies are accumulated.", units=
"s-1", default=0.0)
543 call get_param(param_file, mdl,
"CTRL_FORCE_PREC_CYC_RATE", cs%prec_cyc_rate, &
544 "The integrated rate at which cyclical precipitation \n"//&
545 "anomalies are accumulated.", units=
"s-1", default=0.0)
546 call get_param(param_file, mdl,
"CTRL_FORCE_SMOOTH_LENGTH", smooth_len, &
547 "The length scales over which controlled forcing \n"//&
548 "anomalies are smoothed.", units=
"m", default=0.0)
549 call get_param(param_file, mdl,
"CTRL_FORCE_LAMDA_HEAT", cs%lam_heat, &
550 "A constant of proportionality between SST anomalies \n"//&
551 "and controlling heat fluxes",
"W m-2 K-1", default=0.0)
552 call get_param(param_file, mdl,
"CTRL_FORCE_LAMDA_PREC", cs%lam_prec, &
553 "A constant of proportionality between SSS anomalies \n"//&
554 "(normalised by mean SSS) and controlling precipitation.", &
555 "kg m-2", default=0.0)
556 call get_param(param_file, mdl,
"CTRL_FORCE_LAMDA_CYC_HEAT", cs%lam_cyc_heat, &
557 "A constant of proportionality between SST anomalies \n"//&
558 "and cyclical controlling heat fluxes",
"W m-2 K-1", default=0.0)
559 call get_param(param_file, mdl,
"CTRL_FORCE_LAMDA_CYC_PREC", cs%lam_cyc_prec, &
560 "A constant of proportionality between SSS anomalies \n"//&
561 "(normalised by mean SSS) and cyclical controlling \n"//&
562 "precipitation.",
"kg m-2", default=0.0)
564 cs%Len2 = smooth_len**2
579 if (
associated(cs))
then 580 if (
associated(cs%heat_0))
deallocate(cs%heat_0)
581 if (
associated(cs%precip_0))
deallocate(cs%precip_0)
582 if (
associated(cs%heat_cyc))
deallocate(cs%heat_cyc)
583 if (
associated(cs%precip_cyc))
deallocate(cs%precip_cyc)
584 if (
associated(cs%avg_SST_anom))
deallocate(cs%avg_SST_anom)
585 if (
associated(cs%avg_SSS_anom))
deallocate(cs%avg_SSS_anom)
586 if (
associated(cs%avg_SSS))
deallocate(cs%avg_SSS)
The following structure contains pointers to various fields which may be used describe the surface st...
type(vardesc) function, public var_desc(name, units, longname, hor_grid, z_grid, t_grid, cmor_field_name, cmor_units, conversion, caller)
Returns a vardesc type whose elements have been filled with the provided fields. The argument name is...
This module implements boundary forcing for MOM6.
integer, parameter, public to_all
Ocean grid type. See mom_grid for details.
subroutine, public controlled_forcing_end(CS)
Clean up this modules control structure.
Provides the ocean grid type.
This module contains I/O framework code.
* By Robert Hallberg, July 2011 *This program contains the subroutines that use control-theory * to adjust the surface heat flux and precipitation, based on the * time-mean or periodically (seasonally) varying anomalies from the * observed state. The techniques behind this are described in * Hallberg and Adcroft (2011, in prep.). *Macros written all in capital letters are defined in MOM_memory.h. *A small fragment of the grid is shown below: *j+1 x ^ x ^ x At x: q * j+1 > o > o > At ^: v, tauy * j x ^ x ^ x At >: u, taux * j > o > o > At o: h, fluxes. * j-1 x ^ x ^ x * i-1 i i+1 At x & ^: * i i+1 At > & o: *The boundaries always run through q grid points (x). *
subroutine, public controlled_forcing_init(Time, G, param_file, diag, CS)
Set up this modules control structure.
logical function, public is_root_pe()
Structure that contains pointers to the boundary forcing used to drive the liquid ocean simulated by ...
subroutine, public apply_ctrl_forcing(SST_anom, SSS_anom, SSS_mean, virt_heat, virt_precip, day_start, dt, G, CS)
This subroutine calls any of the other subroutines in this file that are needed to specify the curren...
subroutine, public mom_mesg(message, verb, all_print)
Type for describing a variable, typically a tracer.
integer function periodic_int(rval, num_period)
This function maps rval into an integer in the range from 1 to num_period.
subroutine, public register_ctrl_forcing_restarts(G, param_file, CS, restart_CS)
This subroutine is used to allocate and register any fields in this module that should be written to ...
real function periodic_real(rval, num_period)
This function shifts rval by an integer multiple of num_period so that 0 <= val_out < num_period...
subroutine, public mom_error(level, message, all_print)