MOM6
mom_controlled_forcing Module Reference

Detailed Description

Use control-theory to adjust the surface heat flux and precipitation.

Adjustments are based on the time-mean or periodically (seasonally) varying anomalies from the observed state.

The techniques behind this are described in Hallberg and Adcroft (2018, in prep.).

  • 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). *

Data Types

type  ctrl_forcing_cs
 Control structure for MOM_controlled_forcing. More...
 

Functions/Subroutines

subroutine, public apply_ctrl_forcing (SST_anom, SSS_anom, SSS_mean, virt_heat, virt_precip, day_start, dt, G, US, CS)
 This subroutine calls any of the other subroutines in this file that are needed to specify the current surface forcing fields. More...
 
integer function periodic_int (rval, num_period)
 This function maps rval into an integer in the range from 1 to num_period. More...
 
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. More...
 
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 or read from the restart file. More...
 
subroutine, public controlled_forcing_init (Time, G, param_file, diag, CS)
 Set up this modules control structure. More...
 
subroutine, public controlled_forcing_end (CS)
 Clean up this modules control structure. More...
 

Function/Subroutine Documentation

◆ apply_ctrl_forcing()

subroutine, public mom_controlled_forcing::apply_ctrl_forcing ( real, dimension(szi_(g),szj_(g)), intent(in)  SST_anom,
real, dimension(szi_(g),szj_(g)), intent(in)  SSS_anom,
real, dimension(szi_(g),szj_(g)), intent(in)  SSS_mean,
real, dimension(szi_(g),szj_(g)), intent(inout)  virt_heat,
real, dimension(szi_(g),szj_(g)), intent(inout)  virt_precip,
type(time_type), intent(in)  day_start,
real, intent(in)  dt,
type(ocean_grid_type), intent(inout)  G,
type(unit_scale_type), intent(in)  US,
type(ctrl_forcing_cs), pointer  CS 
)

This subroutine calls any of the other subroutines in this file that are needed to specify the current surface forcing fields.

Parameters
[in,out]gThe ocean's grid structure.
[in]sst_anomThe sea surface temperature anomalies [degC].
[in]sss_anomThe sea surface salinity anomlies [ppt].
[in]sss_meanThe mean sea surface salinity [ppt].
[in,out]virt_heatVirtual (corrective) heat fluxes that are augmented in this subroutine [W m-2].
[in,out]virt_precipVirtual (corrective) precipitation fluxes that are augmented in this subroutine [kg m-2 s-1].
[in]day_startStart time of the fluxes.
[in]dtLength of time over which these fluxes will be applied [s].
[in]usA dimensional unit scaling type
csA pointer to the control structure returned by a previous call to ctrl_forcing_init.

Definition at line 83 of file MOM_controlled_forcing.F90.

83  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure.
84  real, dimension(SZI_(G),SZJ_(G)), intent(in) :: SST_anom !< The sea surface temperature
85  !! anomalies [degC].
86  real, dimension(SZI_(G),SZJ_(G)), intent(in) :: SSS_anom !< The sea surface salinity
87  !! anomlies [ppt].
88  real, dimension(SZI_(G),SZJ_(G)), intent(in) :: SSS_mean !< The mean sea surface
89  !! salinity [ppt].
90  real, dimension(SZI_(G),SZJ_(G)), intent(inout) :: virt_heat !< Virtual (corrective) heat
91  !! fluxes that are augmented
92  !! in this subroutine [W m-2].
93  real, dimension(SZI_(G),SZJ_(G)), intent(inout) :: virt_precip !< Virtual (corrective)
94  !! precipitation fluxes that
95  !! are augmented in this
96  !! subroutine [kg m-2 s-1].
97  type(time_type), intent(in) :: day_start !< Start time of the fluxes.
98  real, intent(in) :: dt !< Length of time over which these
99  !! fluxes will be applied [s].
100  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
101  type(ctrl_forcing_CS), pointer :: CS !< A pointer to the control structure
102  !! returned by a previous call to
103  !! ctrl_forcing_init.
104 !
105  real, dimension(SZIB_(G),SZJ_(G)) :: &
106  flux_heat_x, &
107  flux_prec_x
108  real, dimension(SZI_(G),SZJB_(G)) :: &
109  flux_heat_y, &
110  flux_prec_y
111  type(time_type) :: day_end
112  real :: coef ! A heat-flux coefficient [m2].
113  real :: mr_st, mr_end, mr_mid, mr_prev, mr_next
114  real :: dt_wt, dt_heat_rate, dt_prec_rate
115  real :: dt1_heat_rate, dt1_prec_rate, dt2_heat_rate, dt2_prec_rate
116  real :: wt_per1, wt_st, wt_end, wt_mid
117  integer :: m_st, m_end, m_mid, m_u1, m_u2, m_u3
118  integer :: yr, mon, day, hr, min, sec
119  integer :: i, j, is, ie, js, je
120 
121  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
122 
123  if (.not.associated(cs)) return
124  if ((cs%num_cycle <= 0) .and. (.not.cs%do_integrated)) return
125 
126  day_end = day_start + real_to_time(dt)
127 
128  do j=js,je ; do i=is,ie
129  virt_heat(i,j) = 0.0 ; virt_precip(i,j) = 0.0
130  enddo ; enddo
131 
132  if (cs%do_integrated) then
133  dt_heat_rate = dt * cs%heat_int_rate
134  dt_prec_rate = dt * cs%prec_int_rate
135  call pass_var(cs%heat_0, g%Domain, complete=.false.)
136  call pass_var(cs%precip_0, g%Domain)
137 
138  do j=js,je ; do i=is-1,ie
139  coef = cs%Len2 * (g%dy_Cu(i,j)*g%IdxCu(i,j))
140  flux_heat_x(i,j) = coef * (cs%heat_0(i,j) - cs%heat_0(i+1,j))
141  flux_prec_x(i,j) = coef * (cs%precip_0(i,j) - cs%precip_0(i+1,j))
142  enddo ; enddo
143  do j=js-1,je ; do i=is,ie
144  coef = cs%Len2 * (g%dx_Cv(i,j)*g%IdyCv(i,j))
145  flux_heat_y(i,j) = coef * (cs%heat_0(i,j) - cs%heat_0(i,j+1))
146  flux_prec_y(i,j) = coef * (cs%precip_0(i,j) - cs%precip_0(i,j+1))
147  enddo ; enddo
148  do j=js,je ; do i=is,ie
149  cs%heat_0(i,j) = cs%heat_0(i,j) + dt_heat_rate * ( &
150  -cs%lam_heat*g%mask2dT(i,j)*sst_anom(i,j) + &
151  (us%m_to_L**2*g%IareaT(i,j) * ((flux_heat_x(i-1,j) - flux_heat_x(i,j)) + &
152  (flux_heat_y(i,j-1) - flux_heat_y(i,j))) ) )
153 
154  cs%precip_0(i,j) = cs%precip_0(i,j) + dt_prec_rate * ( &
155  cs%lam_prec * g%mask2dT(i,j)*(sss_anom(i,j) / sss_mean(i,j)) + &
156  (us%m_to_L**2*g%IareaT(i,j) * ((flux_prec_x(i-1,j) - flux_prec_x(i,j)) + &
157  (flux_prec_y(i,j-1) - flux_prec_y(i,j))) ) )
158 
159  virt_heat(i,j) = virt_heat(i,j) + cs%heat_0(i,j)
160  virt_precip(i,j) = virt_precip(i,j) + cs%precip_0(i,j)
161  enddo ; enddo
162  endif
163 
164  if (cs%num_cycle > 0) then
165  ! Determine the current period, with values that run from 0 to CS%num_cycle.
166  call get_date(day_start, yr, mon, day, hr, min, sec)
167  mr_st = cs%num_cycle * (time_type_to_real(day_start - set_date(yr, 1, 1)) / &
168  time_type_to_real(set_date(yr+1, 1, 1) - set_date(yr, 1, 1)))
169 
170  call get_date(day_end, yr, mon, day, hr, min, sec)
171  mr_end = cs%num_cycle * (time_type_to_real(day_end - set_date(yr, 1, 1)) / &
172  time_type_to_real(set_date(yr+1, 1, 1) - set_date(yr, 1, 1)))
173 
174  ! The Chapeau functions are centered at whole integer values that are nominally
175  ! the end of the month to enable simple conversion from the fractional-years times
176  ! CS%num_cycle.
177 
178  ! The month-average temperatures have as an index the month number.
179 
180  m_end = periodic_int(real(ceiling(mr_end)), cs%num_cycle)
181  m_mid = periodic_int(real(ceiling(mr_st)), cs%num_cycle)
182  m_st = periodic_int(mr_st, cs%num_cycle)
183 
184  mr_st = periodic_real(mr_st, cs%num_cycle)
185  mr_end = periodic_real(mr_end, cs%num_cycle)
186  ! mr_mid = periodic_real(ceiling(mr_st), CS%num_cycle)
187  mr_prev = periodic_real(real(floor(mr_st)), cs%num_cycle)
188  mr_next = periodic_real(real(m_end), cs%num_cycle)
189  if (m_mid == m_end) then ; mr_mid = mr_end ! There is only one cell.
190  else ; mr_mid = periodic_real(real(m_mid), cs%num_cycle) ; endif
191 
192  ! There may be two cells that run from mr_st to mr_mid and mr_mid to mr_end.
193 
194  ! The values of m for weights are all calculated relative to mr_prev, so
195  ! check whether mr_mid, etc., need to be shifted by CS%num_cycle, so that these
196  ! values satisfiy mr_prev <= mr_st < mr_mid <= mr_end <= mr_next.
197  if (mr_st < mr_prev) mr_prev = mr_prev - cs%num_cycle
198  if (mr_mid < mr_st) mr_mid = mr_mid + cs%num_cycle
199  if (mr_end < mr_st) mr_end = mr_end + cs%num_cycle
200  if (mr_next < mr_prev) mr_next = mr_next + cs%num_cycle
201 
202  !### These might be removed later - they are to check the coding.
203  if ((mr_mid < mr_st) .or. (mr_mid > mr_prev + 1.)) call mom_error(fatal, &
204  "apply ctrl_forcing: m_mid interpolation out of bounds; fix the code.")
205  if ((mr_end < mr_st) .or. (mr_end > mr_prev + 2.)) call mom_error(fatal, &
206  "apply ctrl_forcing: m_end interpolation out of bounds; fix the code.")
207  if (mr_end > mr_next) call mom_error(fatal, &
208  "apply ctrl_forcing: mr_next interpolation out of bounds; fix the code.")
209 
210  wt_per1 = 1.0
211  if (mr_mid < mr_end) wt_per1 = (mr_mid - mr_st) / (mr_end - mr_st)
212 
213  ! Find the 3 Chapeau-function weights, bearing in mind that m_end may be m_mid.
214  wt_st = wt_per1 * (1. + (mr_prev - 0.5*(mr_st + mr_mid)))
215  wt_end = (1.0-wt_per1) * (1. + (0.5*(mr_end + mr_mid) - mr_next))
216  wt_mid = 1.0 - (wt_st + wt_end)
217  if ((wt_st < 0.0) .or. (wt_end < 0.0) .or. (wt_mid < 0.0)) &
218  call mom_error(fatal, "apply_ctrl_forcing: Negative m weights")
219  if ((wt_st > 1.0) .or. (wt_end > 1.0) .or. (wt_mid > 1.0)) &
220  call mom_error(fatal, "apply_ctrl_forcing: Excessive m weights")
221 
222  ! Add to vert_heat and vert_precip.
223  do j=js,je ; do i=is,ie
224  virt_heat(i,j) = virt_heat(i,j) + (wt_st * cs%heat_cyc(i,j,m_st) + &
225  (wt_mid * cs%heat_cyc(i,j,m_mid) + &
226  wt_end * cs%heat_cyc(i,j,m_end)))
227  virt_precip(i,j) = virt_precip(i,j) + (wt_st * cs%precip_cyc(i,j,m_st) + &
228  (wt_mid * cs%precip_cyc(i,j,m_mid) + &
229  wt_end * cs%precip_cyc(i,j,m_end)))
230  enddo ; enddo
231 
232  ! If different from the last period, take the average and determine the
233  ! chapeau weighting
234 
235  ! The Chapeau functions are centered at whole integer values that are nominally
236  ! the end of the month to enable simple conversion from the fractional-years times
237  ! CS%num_cycle.
238 
239  ! The month-average temperatures have as an index the month number, so the averages
240  ! apply to indicies m_end and m_mid.
241 
242  if (cs%avg_time(m_end) <= 0.0) then ! zero out the averages.
243  cs%avg_time(m_end) = 0.0
244  do j=js,je ; do i=is,ie
245  cs%avg_SST_anom(i,j,m_end) = 0.0
246  cs%avg_SSS_anom(i,j,m_end) = 0.0 ; cs%avg_SSS(i,j,m_end) = 0.0
247  enddo ; enddo
248  endif
249  if (cs%avg_time(m_mid) <= 0.0) then ! zero out the averages.
250  cs%avg_time(m_mid) = 0.0
251  do j=js,je ; do i=is,ie
252  cs%avg_SST_anom(i,j,m_mid) = 0.0
253  cs%avg_SSS_anom(i,j,m_mid) = 0.0 ; cs%avg_SSS(i,j,m_mid) = 0.0
254  enddo ; enddo
255  endif
256 
257  ! Accumulate the average anomalies for this period.
258  dt_wt = wt_per1 * dt
259  cs%avg_time(m_mid) = cs%avg_time(m_mid) + dt_wt
260  do j=js,je ; do i=is,ie
261  cs%avg_SST_anom(i,j,m_mid) = cs%avg_SST_anom(i,j,m_mid) + &
262  dt_wt * g%mask2dT(i,j) * sst_anom(i,j)
263  cs%avg_SSS_anom(i,j,m_mid) = cs%avg_SSS_anom(i,j,m_mid) + &
264  dt_wt * g%mask2dT(i,j) * sss_anom(i,j)
265  cs%avg_SSS(i,j,m_mid) = cs%avg_SSS(i,j,m_mid) + dt_wt * sss_mean(i,j)
266  enddo ; enddo
267  if (wt_per1 < 1.0) then
268  dt_wt = (1.0-wt_per1) * dt
269  cs%avg_time(m_end) = cs%avg_time(m_end) + dt_wt
270  do j=js,je ; do i=is,ie
271  cs%avg_SST_anom(i,j,m_end) = cs%avg_SST_anom(i,j,m_end) + &
272  dt_wt * g%mask2dT(i,j) * sst_anom(i,j)
273  cs%avg_SSS_anom(i,j,m_end) = cs%avg_SSS_anom(i,j,m_end) + &
274  dt_wt * g%mask2dT(i,j) * sss_anom(i,j)
275  cs%avg_SSS(i,j,m_end) = cs%avg_SSS(i,j,m_end) + dt_wt * sss_mean(i,j)
276  enddo ; enddo
277  endif
278 
279  ! Update the Chapeau magnitudes for 4 cycles ago.
280  m_u1 = periodic_int(m_st - 4.0, cs%num_cycle)
281  m_u2 = periodic_int(m_st - 3.0, cs%num_cycle)
282  m_u3 = periodic_int(m_st - 2.0, cs%num_cycle)
283 
284  if (cs%avg_time(m_u1) > 0.0) then
285  do j=js,je ; do i=is,ie
286  cs%avg_SST_anom(i,j,m_u1) = cs%avg_SST_anom(i,j,m_u1) / cs%avg_time(m_u1)
287  cs%avg_SSS_anom(i,j,m_u1) = cs%avg_SSS_anom(i,j,m_u1) / cs%avg_time(m_u1)
288  cs%avg_SSS(i,j,m_u1) = cs%avg_SSS(i,j,m_u1) / cs%avg_time(m_u1)
289  enddo ; enddo
290  cs%avg_time(m_u1) = -1.0
291  endif
292  if (cs%avg_time(m_u2) > 0.0) then
293  do j=js,je ; do i=is,ie
294  cs%avg_SST_anom(i,j,m_u2) = cs%avg_SST_anom(i,j,m_u2) / cs%avg_time(m_u2)
295  cs%avg_SSS_anom(i,j,m_u2) = cs%avg_SSS_anom(i,j,m_u2) / cs%avg_time(m_u2)
296  cs%avg_SSS(i,j,m_u2) = cs%avg_SSS(i,j,m_u2) / cs%avg_time(m_u2)
297  enddo ; enddo
298  cs%avg_time(m_u2) = -1.0
299  endif
300  if (cs%avg_time(m_u3) > 0.0) then
301  do j=js,je ; do i=is,ie
302  cs%avg_SST_anom(i,j,m_u3) = cs%avg_SST_anom(i,j,m_u3) / cs%avg_time(m_u3)
303  cs%avg_SSS_anom(i,j,m_u3) = cs%avg_SSS_anom(i,j,m_u3) / cs%avg_time(m_u3)
304  cs%avg_SSS(i,j,m_u3) = cs%avg_SSS(i,j,m_u3) / cs%avg_time(m_u3)
305  enddo ; enddo
306  cs%avg_time(m_u3) = -1.0
307  endif
308 
309  dt1_heat_rate = wt_per1 * dt * cs%heat_cyc_rate
310  dt1_prec_rate = wt_per1 * dt * cs%prec_cyc_rate
311  dt2_heat_rate = (1.0-wt_per1) * dt * cs%heat_cyc_rate
312  dt2_prec_rate = (1.0-wt_per1) * dt * cs%prec_cyc_rate
313 
314  if (wt_per1 < 1.0) then
315  call pass_var(cs%heat_cyc(:,:,m_u2), g%Domain, complete=.false.)
316  call pass_var(cs%precip_cyc(:,:,m_u2), g%Domain, complete=.false.)
317  endif
318  call pass_var(cs%heat_cyc(:,:,m_u1), g%Domain, complete=.false.)
319  call pass_var(cs%precip_cyc(:,:,m_u1), g%Domain)
320 
321  if ((cs%avg_time(m_u1) == -1.0) .and. (cs%avg_time(m_u2) == -1.0)) then
322  do j=js,je ; do i=is-1,ie
323  coef = cs%Len2 * (g%dy_Cu(i,j)*g%IdxCu(i,j))
324  flux_heat_x(i,j) = coef * (cs%heat_cyc(i,j,m_u1) - cs%heat_cyc(i+1,j,m_u1))
325  flux_prec_x(i,j) = coef * (cs%precip_cyc(i,j,m_u1) - cs%precip_cyc(i+1,j,m_u1))
326  enddo ; enddo
327  do j=js-1,je ; do i=is,ie
328  coef = cs%Len2 * (g%dx_Cv(i,j)*g%IdyCv(i,j))
329  flux_heat_y(i,j) = coef * (cs%heat_cyc(i,j,m_u1) - cs%heat_cyc(i,j+1,m_u1))
330  flux_prec_y(i,j) = coef * (cs%precip_cyc(i,j,m_u1) - cs%precip_cyc(i,j+1,m_u1))
331  enddo ; enddo
332  do j=js,je ; do i=is,ie
333  cs%heat_cyc(i,j,m_u1) = cs%heat_cyc(i,j,m_u1) + dt1_heat_rate * ( &
334  -cs%lam_cyc_heat*(cs%avg_SST_anom(i,j,m_u2) - cs%avg_SST_anom(i,j,m_u1)) + &
335  (us%m_to_L**2*g%IareaT(i,j) * ((flux_heat_x(i-1,j) - flux_heat_x(i,j)) + &
336  (flux_heat_y(i,j-1) - flux_heat_y(i,j))) ) )
337 
338  cs%precip_cyc(i,j,m_u1) = cs%precip_cyc(i,j,m_u1) + dt1_prec_rate * ( &
339  cs%lam_cyc_prec * (cs%avg_SSS_anom(i,j,m_u2) - cs%avg_SSS_anom(i,j,m_u1)) / &
340  (0.5*(cs%avg_SSS(i,j,m_u2) + cs%avg_SSS(i,j,m_u1))) + &
341  (us%m_to_L**2*g%IareaT(i,j) * ((flux_prec_x(i-1,j) - flux_prec_x(i,j)) + &
342  (flux_prec_y(i,j-1) - flux_prec_y(i,j))) ) )
343  enddo ; enddo
344  endif
345 
346  if ((wt_per1 < 1.0) .and. (cs%avg_time(m_u1) == -1.0) .and. (cs%avg_time(m_u2) == -1.0)) then
347  do j=js,je ; do i=is-1,ie
348  coef = cs%Len2 * (g%dy_Cu(i,j)*g%IdxCu(i,j))
349  flux_heat_x(i,j) = coef * (cs%heat_cyc(i,j,m_u2) - cs%heat_cyc(i+1,j,m_u2))
350  flux_prec_x(i,j) = coef * (cs%precip_cyc(i,j,m_u2) - cs%precip_cyc(i+1,j,m_u2))
351  enddo ; enddo
352  do j=js-1,je ; do i=is,ie
353  coef = cs%Len2 * (g%dx_Cv(i,j)*g%IdyCv(i,j))
354  flux_heat_y(i,j) = coef * (cs%heat_cyc(i,j,m_u2) - cs%heat_cyc(i,j+1,m_u2))
355  flux_prec_y(i,j) = coef * (cs%precip_cyc(i,j,m_u2) - cs%precip_cyc(i,j+1,m_u2))
356  enddo ; enddo
357  do j=js,je ; do i=is,ie
358  cs%heat_cyc(i,j,m_u2) = cs%heat_cyc(i,j,m_u2) + dt1_heat_rate * ( &
359  -cs%lam_cyc_heat*(cs%avg_SST_anom(i,j,m_u3) - cs%avg_SST_anom(i,j,m_u2)) + &
360  (us%m_to_L**2*g%IareaT(i,j) * ((flux_heat_x(i-1,j) - flux_heat_x(i,j)) + &
361  (flux_heat_y(i,j-1) - flux_heat_y(i,j))) ) )
362 
363  cs%precip_cyc(i,j,m_u2) = cs%precip_cyc(i,j,m_u2) + dt1_prec_rate * ( &
364  cs%lam_cyc_prec * (cs%avg_SSS_anom(i,j,m_u3) - cs%avg_SSS_anom(i,j,m_u2)) / &
365  (0.5*(cs%avg_SSS(i,j,m_u3) + cs%avg_SSS(i,j,m_u2))) + &
366  (us%m_to_L**2*g%IareaT(i,j) * ((flux_prec_x(i-1,j) - flux_prec_x(i,j)) + &
367  (flux_prec_y(i,j-1) - flux_prec_y(i,j))) ) )
368  enddo ; enddo
369  endif
370 
371  endif ! (CS%num_cycle > 0)
372 

References mom_error_handler::mom_error(), periodic_int(), periodic_real(), and mom_time_manager::real_to_time().

Here is the call graph for this function:

◆ controlled_forcing_end()

subroutine, public mom_controlled_forcing::controlled_forcing_end ( type(ctrl_forcing_cs), pointer  CS)

Clean up this modules control structure.

Parameters
csA pointer to the control structure returned by a previous call to controlled_forcing_init, it will be deallocated here.

Definition at line 565 of file MOM_controlled_forcing.F90.

565  type(ctrl_forcing_CS), pointer :: CS !< A pointer to the control structure
566  !! returned by a previous call to
567  !! controlled_forcing_init, it will be
568  !! deallocated here.
569 
570  if (associated(cs)) then
571  if (associated(cs%heat_0)) deallocate(cs%heat_0)
572  if (associated(cs%precip_0)) deallocate(cs%precip_0)
573  if (associated(cs%heat_cyc)) deallocate(cs%heat_cyc)
574  if (associated(cs%precip_cyc)) deallocate(cs%precip_cyc)
575  if (associated(cs%avg_SST_anom)) deallocate(cs%avg_SST_anom)
576  if (associated(cs%avg_SSS_anom)) deallocate(cs%avg_SSS_anom)
577  if (associated(cs%avg_SSS)) deallocate(cs%avg_SSS)
578 
579  deallocate(cs)
580  endif
581  cs => null()
582 

◆ controlled_forcing_init()

subroutine, public mom_controlled_forcing::controlled_forcing_init ( type(time_type), intent(in)  Time,
type(ocean_grid_type), intent(in)  G,
type(param_file_type), intent(in)  param_file,
type(diag_ctrl), intent(in), target  diag,
type(ctrl_forcing_cs), pointer  CS 
)

Set up this modules control structure.

Parameters
[in]timeThe current model time.
[in]gThe ocean's grid structure.
[in]param_fileA structure indicating the open file to parse for model parameter values.
[in]diagA structure that is used to regulate diagnostic output.
csA pointer that is set to point to the control structure for this module.

Definition at line 485 of file MOM_controlled_forcing.F90.

485  type(time_type), intent(in) :: Time !< The current model time.
486  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
487  type(param_file_type), intent(in) :: param_file !< A structure indicating the
488  !! open file to parse for model
489  !! parameter values.
490  type(diag_ctrl), target, intent(in) :: diag !< A structure that is used to regulate
491  !! diagnostic output.
492  type(ctrl_forcing_CS), pointer :: CS !< A pointer that is set to point to the
493  !! control structure for this module.
494  real :: smooth_len
495  logical :: do_integrated
496  integer :: num_cycle
497 ! This include declares and sets the variable "version".
498 #include "version_variable.h"
499  character(len=40) :: mdl = "MOM_controlled_forcing" ! This module's name.
500 
501  ! These should have already been called.
502  ! call read_param(param_file, "CTRL_FORCE_INTEGRATED", CS%do_integrated)
503  ! call read_param(param_file, "CTRL_FORCE_NUM_CYCLE", CS%num_cycle)
504 
505  if (associated(cs)) then
506  do_integrated = cs%do_integrated ; num_cycle = cs%num_cycle
507  else
508  do_integrated = .false. ; num_cycle = 0
509  endif
510 
511  ! Read all relevant parameters and write them to the model log.
512  call log_version(param_file, mdl, version, "")
513  call log_param(param_file, mdl, "CTRL_FORCE_INTEGRATED", do_integrated, &
514  "If true, use a PI controller to determine the surface "//&
515  "forcing that is consistent with the observed mean properties.", &
516  default=.false.)
517  call log_param(param_file, mdl, "CTRL_FORCE_NUM_CYCLE", num_cycle, &
518  "The number of cycles per year in the controlled forcing, "//&
519  "or 0 for no cyclic forcing.", default=0)
520 
521  if (.not.associated(cs)) return
522 
523  cs%diag => diag
524 
525  call get_param(param_file, mdl, "CTRL_FORCE_HEAT_INT_RATE", cs%heat_int_rate, &
526  "The integrated rate at which heat flux anomalies are "//&
527  "accumulated.", units="s-1", default=0.0)
528  call get_param(param_file, mdl, "CTRL_FORCE_PREC_INT_RATE", cs%prec_int_rate, &
529  "The integrated rate at which precipitation anomalies "//&
530  "are accumulated.", units="s-1", default=0.0)
531  call get_param(param_file, mdl, "CTRL_FORCE_HEAT_CYC_RATE", cs%heat_cyc_rate, &
532  "The integrated rate at which cyclical heat flux "//&
533  "anomalies are accumulated.", units="s-1", default=0.0)
534  call get_param(param_file, mdl, "CTRL_FORCE_PREC_CYC_RATE", cs%prec_cyc_rate, &
535  "The integrated rate at which cyclical precipitation "//&
536  "anomalies are accumulated.", units="s-1", default=0.0)
537  call get_param(param_file, mdl, "CTRL_FORCE_SMOOTH_LENGTH", smooth_len, &
538  "The length scales over which controlled forcing "//&
539  "anomalies are smoothed.", units="m", default=0.0)
540  call get_param(param_file, mdl, "CTRL_FORCE_LAMDA_HEAT", cs%lam_heat, &
541  "A constant of proportionality between SST anomalies "//&
542  "and controlling heat fluxes", "W m-2 K-1", default=0.0)
543  call get_param(param_file, mdl, "CTRL_FORCE_LAMDA_PREC", cs%lam_prec, &
544  "A constant of proportionality between SSS anomalies "//&
545  "(normalised by mean SSS) and controlling precipitation.", &
546  "kg m-2", default=0.0)
547  call get_param(param_file, mdl, "CTRL_FORCE_LAMDA_CYC_HEAT", cs%lam_cyc_heat, &
548  "A constant of proportionality between SST anomalies "//&
549  "and cyclical controlling heat fluxes", "W m-2 K-1", default=0.0)
550  call get_param(param_file, mdl, "CTRL_FORCE_LAMDA_CYC_PREC", cs%lam_cyc_prec, &
551  "A constant of proportionality between SSS anomalies "//&
552  "(normalised by mean SSS) and cyclical controlling "//&
553  "precipitation.", "kg m-2", default=0.0)
554 
555  cs%Len2 = smooth_len**2
556 
557 ! ### REPLACE THIS WITH ANY DIAGNOSTICS FROM THIS MODULE.
558 ! CS%id_taux = register_diag_field('ocean_model', 'taux', diag%axesu1, Time, &
559 ! 'Zonal Wind Stress', 'Pascal')
560 

◆ periodic_int()

integer function mom_controlled_forcing::periodic_int ( real, intent(in)  rval,
integer, intent(in)  num_period 
)
private

This function maps rval into an integer in the range from 1 to num_period.

Parameters
[in]rvalInput for mapping.
[in]num_periodMaximum output.
Returns
Return value.

Definition at line 377 of file MOM_controlled_forcing.F90.

377  real, intent(in) :: rval !< Input for mapping.
378  integer, intent(in) :: num_period !< Maximum output.
379  integer :: m !< Return value.
380 
381  m = floor(rval)
382  if (m <= 0) then
383  m = m + num_period * (1 + (abs(m) / num_period))
384  elseif (m > num_period) then
385  m = m - num_period * ((m-1) / num_period)
386  endif

Referenced by apply_ctrl_forcing().

Here is the caller graph for this function:

◆ periodic_real()

real function mom_controlled_forcing::periodic_real ( real, intent(in)  rval,
integer, intent(in)  num_period 
)
private

This function shifts rval by an integer multiple of num_period so that 0 <= val_out < num_period.

Parameters
[in]rvalInput to be shifted into valid range.
[in]num_periodMaximum valid value.
Returns
Return value.

Definition at line 392 of file MOM_controlled_forcing.F90.

392  real, intent(in) :: rval !< Input to be shifted into valid range.
393  integer, intent(in) :: num_period !< Maximum valid value.
394  real :: val_out !< Return value.
395  integer :: nshft
396 
397  if (rval < 0) then ; nshft = floor(abs(rval) / num_period) + 1
398  elseif (rval < num_period) then ; nshft = 0
399  else ; nshft = -1*floor(rval / num_period) ; endif
400 
401  val_out = rval + nshft * num_period

Referenced by apply_ctrl_forcing().

Here is the caller graph for this function:

◆ register_ctrl_forcing_restarts()

subroutine, public mom_controlled_forcing::register_ctrl_forcing_restarts ( type(ocean_grid_type), intent(in)  G,
type(param_file_type), intent(in)  param_file,
type(ctrl_forcing_cs), pointer  CS,
type(mom_restart_cs), pointer  restart_CS 
)

This subroutine is used to allocate and register any fields in this module that should be written to or read from the restart file.

Parameters
[in]gThe ocean's grid structure.
[in]param_fileA structure indicating the open file to parse for model parameter values.
csA pointer that is set to point to the control structure for this module.
restart_csA pointer to the restart control structure.

Definition at line 408 of file MOM_controlled_forcing.F90.

408  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
409  type(param_file_type), intent(in) :: param_file !< A structure indicating the
410  !! open file to parse for model
411  !! parameter values.
412  type(ctrl_forcing_CS), pointer :: CS !< A pointer that is set to point to the
413  !! control structure for this module.
414  type(MOM_restart_CS), pointer :: restart_CS !< A pointer to the restart control structure.
415 
416  logical :: controlled, use_temperature
417  character (len=8) :: period_str
418  type(vardesc) :: vd
419  integer :: isd, ied, jsd, jed, IsdB, IedB, JsdB, JedB
420  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
421  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
422 
423  if (associated(cs)) then
424  call mom_error(warning, "register_ctrl_forcing_restarts called "//&
425  "with an associated control structure.")
426  return
427  endif
428 
429  controlled = .false.
430  call read_param(param_file, "CONTROLLED_FORCING", controlled)
431  if (.not.controlled) return
432 
433  use_temperature = .true.
434  call read_param(param_file, "ENABLE_THERMODYNAMICS", use_temperature)
435  if (.not.use_temperature) call mom_error(fatal, &
436  "register_ctrl_forcing_restarts: CONTROLLED_FORCING only works with "//&
437  "ENABLE_THERMODYNAMICS defined.")
438 
439  allocate(cs)
440 
441  cs%do_integrated = .true. ; cs%num_cycle = 0
442  call read_param(param_file, "CTRL_FORCE_INTEGRATED", cs%do_integrated)
443  call read_param(param_file, "CTRL_FORCE_NUM_CYCLE", cs%num_cycle)
444 
445  if (cs%do_integrated) then
446  call safe_alloc_ptr(cs%heat_0,isd,ied,jsd,jed) ; cs%heat_0(:,:) = 0.0
447  call safe_alloc_ptr(cs%precip_0,isd,ied,jsd,jed) ; cs%precip_0(:,:) = 0.0
448  vd = var_desc("Ctrl_heat","W m-2","Control Integrative Heating",z_grid='1')
449  call register_restart_field(cs%heat_0, vd, .false., restart_cs)
450  vd = var_desc("Ctrl_precip","kg m-2 s-1","Control Integrative Precipitation",z_grid='1')
451  call register_restart_field(cs%precip_0, vd, .false., restart_cs)
452  endif
453 
454  if (cs%num_cycle > 0) then
455  write (period_str, '(i8)') cs%num_cycle
456  period_str = trim('p ')//trim(adjustl(period_str))
457  call safe_alloc_ptr(cs%heat_cyc,isd,ied,jsd,jed,cs%num_cycle) ; cs%heat_cyc(:,:,:) = 0.0
458  call safe_alloc_ptr(cs%precip_cyc,isd,ied,jsd,jed,cs%num_cycle) ; cs%precip_cyc(:,:,:) = 0.0
459  vd = var_desc("Ctrl_heat_cycle", "W m-2","Cyclical Control Heating",&
460  z_grid='1', t_grid=period_str)
461  call register_restart_field(cs%heat_cyc, vd, .false., restart_cs)
462  vd = var_desc("Ctrl_precip_cycle","kg m-2 s-1","Cyclical Control Precipitation", &
463  z_grid='1', t_grid=period_str)
464  call register_restart_field(cs%precip_cyc, vd, .false., restart_cs)
465 
466  call safe_alloc_ptr(cs%avg_time,cs%num_cycle) ; cs%avg_time(:) = 0.0
467  vd = var_desc("avg_time","sec","Cyclical accumulated averaging time", &
468  '1',z_grid='1',t_grid=period_str)
469  call register_restart_field(cs%avg_time, vd, .false., restart_cs)
470 
471  call safe_alloc_ptr(cs%avg_SST_anom,isd,ied,jsd,jed,cs%num_cycle) ; cs%avg_SST_anom(:,:,:) = 0.0
472  call safe_alloc_ptr(cs%avg_SSS_anom,isd,ied,jsd,jed,cs%num_cycle) ; cs%avg_SSS_anom(:,:,:) = 0.0
473  vd = var_desc("avg_SST_anom","deg C","Cyclical average SST Anomaly", &
474  z_grid='1',t_grid=period_str)
475  call register_restart_field(cs%avg_SST_anom, vd, .false., restart_cs)
476  vd = var_desc("avg_SSS_anom","g kg-1","Cyclical average SSS Anomaly", &
477  z_grid='1',t_grid=period_str)
478  call register_restart_field(cs%avg_SSS_anom, vd, .false., restart_cs)
479  endif
480 

References mom_error_handler::mom_error(), and mom_io::var_desc().

Here is the call graph for this function: