MOM6
MOM_controlled_forcing.F90
1 !> Use control-theory to adjust the surface heat flux and precipitation.
2 !!
3 !! Adjustments are based on the time-mean or periodically (seasonally) varying
4 !! anomalies from the observed state.
5 !!
6 !! The techniques behind this are described in Hallberg and Adcroft (2018, in prep.).
8 
9 ! This file is part of MOM6. See LICENSE.md for the license.
10 
11 use mom_diag_mediator, only : post_data, query_averaging_enabled
12 use mom_diag_mediator, only : register_diag_field, diag_ctrl, safe_alloc_ptr
13 use mom_domains, only : pass_var, pass_vector, agrid, to_south, to_west, to_all
14 use mom_error_handler, only : mom_error, fatal, warning, mom_mesg, is_root_pe
16 use mom_forcing_type, only : forcing
17 use mom_grid, only : ocean_grid_type
18 use mom_io, only : vardesc, var_desc
20 use mom_time_manager, only : time_type, operator(+), operator(/), operator(-)
21 use mom_time_manager, only : get_date, set_date
22 use mom_time_manager, only : time_type_to_real, real_to_time
23 use mom_variables, only : surface
24 
25 implicit none ; private
26 
27 #include <MOM_memory.h>
28 
29 public apply_ctrl_forcing, register_ctrl_forcing_restarts
30 public controlled_forcing_init, controlled_forcing_end
31 
32 !> Control structure for MOM_controlled_forcing
33 type, public :: ctrl_forcing_cs ; private
34  logical :: use_temperature !< If true, temperature and salinity are used as
35  !! state variables.
36  logical :: do_integrated !< If true, use time-integrated anomalies to control
37  !! the surface state.
38  integer :: num_cycle !< The number of elements in the forcing cycle.
39  real :: heat_int_rate !< The rate at which heating anomalies accumulate [s-1].
40  real :: prec_int_rate !< The rate at which precipitation anomalies accumulate [s-1].
41  real :: heat_cyc_rate !< The rate at which cyclical heating anomaliess
42  !! accumulate [s-1].
43  real :: prec_cyc_rate !< The rate at which cyclical precipitation anomaliess
44  !! accumulate [s-1].
45  real :: len2 !< The square of the length scale over which the anomalies
46  !! are smoothed via a Laplacian filter [m2].
47  real :: lam_heat !< A constant of proportionality between SST anomalies
48  !! and heat fluxes [W m-2 degC-1].
49  real :: lam_prec !< A constant of proportionality between SSS anomalies
50  !! (normalised by mean SSS) and precipitation [kg m-2].
51  real :: lam_cyc_heat !< A constant of proportionality between cyclical SST
52  !! anomalies and corrective heat fluxes [W m-2 degC-1].
53  real :: lam_cyc_prec !< A constant of proportionality between cyclical SSS
54  !! anomalies (normalised by mean SSS) and corrective
55  !! precipitation [kg m-2].
56 
57  !>@{ Pointers for data.
58  !! \todo Needs more complete documentation.
59  real, pointer, dimension(:) :: &
60  avg_time => null()
61  real, pointer, dimension(:,:) :: &
62  heat_0 => null(), &
63  precip_0 => null()
64  real, pointer, dimension(:,:,:) :: &
65  heat_cyc => null(), &
66  precip_cyc => null(), &
67  avg_sst_anom => null(), &
68  avg_sss_anom => null(), &
69  avg_sss => null()
70  !!@}
71  type(diag_ctrl), pointer :: diag => null() !< A structure that is used to
72  !! regulate the timing of diagnostic output.
73  integer :: id_heat_0 = -1 !< Diagnostic handle
74 end type ctrl_forcing_cs
75 
76 contains
77 
78 !> This subroutine calls any of the other subroutines in this file
79 !! that are needed to specify the current surface forcing fields.
80 subroutine apply_ctrl_forcing(SST_anom, SSS_anom, SSS_mean, virt_heat, virt_precip, &
81  day_start, dt, G, CS)
82  type(ocean_grid_type), intent(inout) :: g !< The ocean's grid structure.
83  real, dimension(SZI_(G),SZJ_(G)), intent(in) :: sst_anom !< The sea surface temperature
84  !! anomalies [degC].
85  real, dimension(SZI_(G),SZJ_(G)), intent(in) :: sss_anom !< The sea surface salinity
86  !! anomlies [ppt].
87  real, dimension(SZI_(G),SZJ_(G)), intent(in) :: sss_mean !< The mean sea surface
88  !! salinity [ppt].
89  real, dimension(SZI_(G),SZJ_(G)), intent(inout) :: virt_heat !< Virtual (corrective) heat
90  !! fluxes that are augmented
91  !! in this subroutine [W m-2].
92  real, dimension(SZI_(G),SZJ_(G)), intent(inout) :: virt_precip !< Virtual (corrective)
93  !! precipitation fluxes that
94  !! are augmented in this
95  !! subroutine [kg m-2 s-1].
96  type(time_type), intent(in) :: day_start !< Start time of the fluxes.
97  real, intent(in) :: dt !< Length of time over which these
98  !! fluxes will be applied [s].
99  type(ctrl_forcing_cs), pointer :: cs !< A pointer to the control structure
100  !! returned by a previous call to
101  !! ctrl_forcing_init.
102 !
103  real, dimension(SZIB_(G),SZJ_(G)) :: &
104  flux_heat_x, &
105  flux_prec_x
106  real, dimension(SZI_(G),SZJB_(G)) :: &
107  flux_heat_y, &
108  flux_prec_y
109  type(time_type) :: day_end
110  real :: coef ! A heat-flux coefficient [m2].
111  real :: mr_st, mr_end, mr_mid, mr_prev, mr_next
112  real :: dt_wt, dt_heat_rate, dt_prec_rate
113  real :: dt1_heat_rate, dt1_prec_rate, dt2_heat_rate, dt2_prec_rate
114  real :: wt_per1, wt_st, wt_end, wt_mid
115  integer :: m_st, m_end, m_mid, m_u1, m_u2, m_u3
116  integer :: yr, mon, day, hr, min, sec
117  integer :: i, j, is, ie, js, je
118 
119  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
120 
121  if (.not.associated(cs)) return
122  if ((cs%num_cycle <= 0) .and. (.not.cs%do_integrated)) return
123 
124  day_end = day_start + real_to_time(dt)
125 
126  do j=js,je ; do i=is,ie
127  virt_heat(i,j) = 0.0 ; virt_precip(i,j) = 0.0
128  enddo ; enddo
129 
130  if (cs%do_integrated) then
131  dt_heat_rate = dt * cs%heat_int_rate
132  dt_prec_rate = dt * cs%prec_int_rate
133  call pass_var(cs%heat_0, g%Domain, complete=.false.)
134  call pass_var(cs%precip_0, g%Domain)
135 
136  do j=js,je ; do i=is-1,ie
137  coef = cs%Len2 * (g%dy_Cu(i,j)*g%IdxCu(i,j))
138  flux_heat_x(i,j) = coef * (cs%heat_0(i,j) - cs%heat_0(i+1,j))
139  flux_prec_x(i,j) = coef * (cs%precip_0(i,j) - cs%precip_0(i+1,j))
140  enddo ; enddo
141  do j=js-1,je ; do i=is,ie
142  coef = cs%Len2 * (g%dx_Cv(i,j)*g%IdyCv(i,j))
143  flux_heat_y(i,j) = coef * (cs%heat_0(i,j) - cs%heat_0(i,j+1))
144  flux_prec_y(i,j) = coef * (cs%precip_0(i,j) - cs%precip_0(i,j+1))
145  enddo ; enddo
146  do j=js,je ; do i=is,ie
147  cs%heat_0(i,j) = cs%heat_0(i,j) + dt_heat_rate * ( &
148  -cs%lam_heat*g%mask2dT(i,j)*sst_anom(i,j) + &
149  (g%IareaT(i,j) * ((flux_heat_x(i-1,j) - flux_heat_x(i,j)) + &
150  (flux_heat_y(i,j-1) - flux_heat_y(i,j))) ) )
151 
152  cs%precip_0(i,j) = cs%precip_0(i,j) + dt_prec_rate * ( &
153  cs%lam_prec * g%mask2dT(i,j)*(sss_anom(i,j) / sss_mean(i,j)) + &
154  (g%IareaT(i,j) * ((flux_prec_x(i-1,j) - flux_prec_x(i,j)) + &
155  (flux_prec_y(i,j-1) - flux_prec_y(i,j))) ) )
156 
157  virt_heat(i,j) = virt_heat(i,j) + cs%heat_0(i,j)
158  virt_precip(i,j) = virt_precip(i,j) + cs%precip_0(i,j)
159  enddo ; enddo
160  endif
161 
162  if (cs%num_cycle > 0) then
163  ! Determine the current period, with values that run from 0 to CS%num_cycle.
164  call get_date(day_start, yr, mon, day, hr, min, sec)
165  mr_st = cs%num_cycle * (time_type_to_real(day_start - set_date(yr, 1, 1)) / &
166  time_type_to_real(set_date(yr+1, 1, 1) - set_date(yr, 1, 1)))
167 
168  call get_date(day_end, yr, mon, day, hr, min, sec)
169  mr_end = cs%num_cycle * (time_type_to_real(day_end - set_date(yr, 1, 1)) / &
170  time_type_to_real(set_date(yr+1, 1, 1) - set_date(yr, 1, 1)))
171 
172  ! The Chapeau functions are centered at whole integer values that are nominally
173  ! the end of the month to enable simple conversion from the fractional-years times
174  ! CS%num_cycle.
175 
176  ! The month-average temperatures have as an index the month number.
177 
178  m_end = periodic_int(real(ceiling(mr_end)), cs%num_cycle)
179  m_mid = periodic_int(real(ceiling(mr_st)), cs%num_cycle)
180  m_st = periodic_int(mr_st, cs%num_cycle)
181 
182  mr_st = periodic_real(mr_st, cs%num_cycle)
183  mr_end = periodic_real(mr_end, cs%num_cycle)
184  ! mr_mid = periodic_real(ceiling(mr_st), CS%num_cycle)
185  mr_prev = periodic_real(real(floor(mr_st)), cs%num_cycle)
186  mr_next = periodic_real(real(m_end), cs%num_cycle)
187  if (m_mid == m_end) then ; mr_mid = mr_end ! There is only one cell.
188  else ; mr_mid = periodic_real(real(m_mid), cs%num_cycle) ; endif
189 
190  ! There may be two cells that run from mr_st to mr_mid and mr_mid to mr_end.
191 
192  ! The values of m for weights are all calculated relative to mr_prev, so
193  ! check whether mr_mid, etc., need to be shifted by CS%num_cycle, so that these
194  ! values satisfiy mr_prev <= mr_st < mr_mid <= mr_end <= mr_next.
195  if (mr_st < mr_prev) mr_prev = mr_prev - cs%num_cycle
196  if (mr_mid < mr_st) mr_mid = mr_mid + cs%num_cycle
197  if (mr_end < mr_st) mr_end = mr_end + cs%num_cycle
198  if (mr_next < mr_prev) mr_next = mr_next + cs%num_cycle
199 
200  !### These might be removed later - they are to check the coding.
201  if ((mr_mid < mr_st) .or. (mr_mid > mr_prev + 1.)) call mom_error(fatal, &
202  "apply ctrl_forcing: m_mid interpolation out of bounds; fix the code.")
203  if ((mr_end < mr_st) .or. (mr_end > mr_prev + 2.)) call mom_error(fatal, &
204  "apply ctrl_forcing: m_end interpolation out of bounds; fix the code.")
205  if (mr_end > mr_next) call mom_error(fatal, &
206  "apply ctrl_forcing: mr_next interpolation out of bounds; fix the code.")
207 
208  wt_per1 = 1.0
209  if (mr_mid < mr_end) wt_per1 = (mr_mid - mr_st) / (mr_end - mr_st)
210 
211  ! Find the 3 Chapeau-function weights, bearing in mind that m_end may be m_mid.
212  wt_st = wt_per1 * (1. + (mr_prev - 0.5*(mr_st + mr_mid)))
213  wt_end = (1.0-wt_per1) * (1. + (0.5*(mr_end + mr_mid) - mr_next))
214  wt_mid = 1.0 - (wt_st + wt_end)
215  if ((wt_st < 0.0) .or. (wt_end < 0.0) .or. (wt_mid < 0.0)) &
216  call mom_error(fatal, "apply_ctrl_forcing: Negative m weights")
217  if ((wt_st > 1.0) .or. (wt_end > 1.0) .or. (wt_mid > 1.0)) &
218  call mom_error(fatal, "apply_ctrl_forcing: Excessive m weights")
219 
220  ! Add to vert_heat and vert_precip.
221  do j=js,je ; do i=is,ie
222  virt_heat(i,j) = virt_heat(i,j) + (wt_st * cs%heat_cyc(i,j,m_st) + &
223  (wt_mid * cs%heat_cyc(i,j,m_mid) + &
224  wt_end * cs%heat_cyc(i,j,m_end)))
225  virt_precip(i,j) = virt_precip(i,j) + (wt_st * cs%precip_cyc(i,j,m_st) + &
226  (wt_mid * cs%precip_cyc(i,j,m_mid) + &
227  wt_end * cs%precip_cyc(i,j,m_end)))
228  enddo ; enddo
229 
230  ! If different from the last period, take the average and determine the
231  ! chapeau weighting
232 
233  ! The Chapeau functions are centered at whole integer values that are nominally
234  ! the end of the month to enable simple conversion from the fractional-years times
235  ! CS%num_cycle.
236 
237  ! The month-average temperatures have as an index the month number, so the averages
238  ! apply to indicies m_end and m_mid.
239 
240  if (cs%avg_time(m_end) <= 0.0) then ! zero out the averages.
241  cs%avg_time(m_end) = 0.0
242  do j=js,je ; do i=is,ie
243  cs%avg_SST_anom(i,j,m_end) = 0.0
244  cs%avg_SSS_anom(i,j,m_end) = 0.0 ; cs%avg_SSS(i,j,m_end) = 0.0
245  enddo ; enddo
246  endif
247  if (cs%avg_time(m_mid) <= 0.0) then ! zero out the averages.
248  cs%avg_time(m_mid) = 0.0
249  do j=js,je ; do i=is,ie
250  cs%avg_SST_anom(i,j,m_mid) = 0.0
251  cs%avg_SSS_anom(i,j,m_mid) = 0.0 ; cs%avg_SSS(i,j,m_mid) = 0.0
252  enddo ; enddo
253  endif
254 
255  ! Accumulate the average anomalies for this period.
256  dt_wt = wt_per1 * dt
257  cs%avg_time(m_mid) = cs%avg_time(m_mid) + dt_wt
258  do j=js,je ; do i=is,ie
259  cs%avg_SST_anom(i,j,m_mid) = cs%avg_SST_anom(i,j,m_mid) + &
260  dt_wt * g%mask2dT(i,j) * sst_anom(i,j)
261  cs%avg_SSS_anom(i,j,m_mid) = cs%avg_SSS_anom(i,j,m_mid) + &
262  dt_wt * g%mask2dT(i,j) * sss_anom(i,j)
263  cs%avg_SSS(i,j,m_mid) = cs%avg_SSS(i,j,m_mid) + dt_wt * sss_mean(i,j)
264  enddo ; enddo
265  if (wt_per1 < 1.0) then
266  dt_wt = (1.0-wt_per1) * dt
267  cs%avg_time(m_end) = cs%avg_time(m_end) + dt_wt
268  do j=js,je ; do i=is,ie
269  cs%avg_SST_anom(i,j,m_end) = cs%avg_SST_anom(i,j,m_end) + &
270  dt_wt * g%mask2dT(i,j) * sst_anom(i,j)
271  cs%avg_SSS_anom(i,j,m_end) = cs%avg_SSS_anom(i,j,m_end) + &
272  dt_wt * g%mask2dT(i,j) * sss_anom(i,j)
273  cs%avg_SSS(i,j,m_end) = cs%avg_SSS(i,j,m_end) + dt_wt * sss_mean(i,j)
274  enddo ; enddo
275  endif
276 
277  ! Update the Chapeau magnitudes for 4 cycles ago.
278  m_u1 = periodic_int(m_st - 4.0, cs%num_cycle)
279  m_u2 = periodic_int(m_st - 3.0, cs%num_cycle)
280  m_u3 = periodic_int(m_st - 2.0, cs%num_cycle)
281 
282  if (cs%avg_time(m_u1) > 0.0) then
283  do j=js,je ; do i=is,ie
284  cs%avg_SST_anom(i,j,m_u1) = cs%avg_SST_anom(i,j,m_u1) / cs%avg_time(m_u1)
285  cs%avg_SSS_anom(i,j,m_u1) = cs%avg_SSS_anom(i,j,m_u1) / cs%avg_time(m_u1)
286  cs%avg_SSS(i,j,m_u1) = cs%avg_SSS(i,j,m_u1) / cs%avg_time(m_u1)
287  enddo ; enddo
288  cs%avg_time(m_u1) = -1.0
289  endif
290  if (cs%avg_time(m_u2) > 0.0) then
291  do j=js,je ; do i=is,ie
292  cs%avg_SST_anom(i,j,m_u2) = cs%avg_SST_anom(i,j,m_u2) / cs%avg_time(m_u2)
293  cs%avg_SSS_anom(i,j,m_u2) = cs%avg_SSS_anom(i,j,m_u2) / cs%avg_time(m_u2)
294  cs%avg_SSS(i,j,m_u2) = cs%avg_SSS(i,j,m_u2) / cs%avg_time(m_u2)
295  enddo ; enddo
296  cs%avg_time(m_u2) = -1.0
297  endif
298  if (cs%avg_time(m_u3) > 0.0) then
299  do j=js,je ; do i=is,ie
300  cs%avg_SST_anom(i,j,m_u3) = cs%avg_SST_anom(i,j,m_u3) / cs%avg_time(m_u3)
301  cs%avg_SSS_anom(i,j,m_u3) = cs%avg_SSS_anom(i,j,m_u3) / cs%avg_time(m_u3)
302  cs%avg_SSS(i,j,m_u3) = cs%avg_SSS(i,j,m_u3) / cs%avg_time(m_u3)
303  enddo ; enddo
304  cs%avg_time(m_u3) = -1.0
305  endif
306 
307  dt1_heat_rate = wt_per1 * dt * cs%heat_cyc_rate
308  dt1_prec_rate = wt_per1 * dt * cs%prec_cyc_rate
309  dt2_heat_rate = (1.0-wt_per1) * dt * cs%heat_cyc_rate
310  dt2_prec_rate = (1.0-wt_per1) * dt * cs%prec_cyc_rate
311 
312  if (wt_per1 < 1.0) then
313  call pass_var(cs%heat_cyc(:,:,m_u2), g%Domain, complete=.false.)
314  call pass_var(cs%precip_cyc(:,:,m_u2), g%Domain, complete=.false.)
315  endif
316  call pass_var(cs%heat_cyc(:,:,m_u1), g%Domain, complete=.false.)
317  call pass_var(cs%precip_cyc(:,:,m_u1), g%Domain)
318 
319  if ((cs%avg_time(m_u1) == -1.0) .and. (cs%avg_time(m_u2) == -1.0)) then
320  do j=js,je ; do i=is-1,ie
321  coef = cs%Len2 * (g%dy_Cu(i,j)*g%IdxCu(i,j))
322  flux_heat_x(i,j) = coef * (cs%heat_cyc(i,j,m_u1) - cs%heat_cyc(i+1,j,m_u1))
323  flux_prec_x(i,j) = coef * (cs%precip_cyc(i,j,m_u1) - cs%precip_cyc(i+1,j,m_u1))
324  enddo ; enddo
325  do j=js-1,je ; do i=is,ie
326  coef = cs%Len2 * (g%dx_Cv(i,j)*g%IdyCv(i,j))
327  flux_heat_y(i,j) = coef * (cs%heat_cyc(i,j,m_u1) - cs%heat_cyc(i,j+1,m_u1))
328  flux_prec_y(i,j) = coef * (cs%precip_cyc(i,j,m_u1) - cs%precip_cyc(i,j+1,m_u1))
329  enddo ; enddo
330  do j=js,je ; do i=is,ie
331  cs%heat_cyc(i,j,m_u1) = cs%heat_cyc(i,j,m_u1) + dt1_heat_rate * ( &
332  -cs%lam_cyc_heat*(cs%avg_SST_anom(i,j,m_u2) - cs%avg_SST_anom(i,j,m_u1)) + &
333  (g%IareaT(i,j) * ((flux_heat_x(i-1,j) - flux_heat_x(i,j)) + &
334  (flux_heat_y(i,j-1) - flux_heat_y(i,j))) ) )
335 
336  cs%precip_cyc(i,j,m_u1) = cs%precip_cyc(i,j,m_u1) + dt1_prec_rate * ( &
337  cs%lam_cyc_prec * (cs%avg_SSS_anom(i,j,m_u2) - cs%avg_SSS_anom(i,j,m_u1)) / &
338  (0.5*(cs%avg_SSS(i,j,m_u2) + cs%avg_SSS(i,j,m_u1))) + &
339  (g%IareaT(i,j) * ((flux_prec_x(i-1,j) - flux_prec_x(i,j)) + &
340  (flux_prec_y(i,j-1) - flux_prec_y(i,j))) ) )
341  enddo ; enddo
342  endif
343 
344  if ((wt_per1 < 1.0) .and. (cs%avg_time(m_u1) == -1.0) .and. (cs%avg_time(m_u2) == -1.0)) then
345  do j=js,je ; do i=is-1,ie
346  coef = cs%Len2 * (g%dy_Cu(i,j)*g%IdxCu(i,j))
347  flux_heat_x(i,j) = coef * (cs%heat_cyc(i,j,m_u2) - cs%heat_cyc(i+1,j,m_u2))
348  flux_prec_x(i,j) = coef * (cs%precip_cyc(i,j,m_u2) - cs%precip_cyc(i+1,j,m_u2))
349  enddo ; enddo
350  do j=js-1,je ; do i=is,ie
351  coef = cs%Len2 * (g%dx_Cv(i,j)*g%IdyCv(i,j))
352  flux_heat_y(i,j) = coef * (cs%heat_cyc(i,j,m_u2) - cs%heat_cyc(i,j+1,m_u2))
353  flux_prec_y(i,j) = coef * (cs%precip_cyc(i,j,m_u2) - cs%precip_cyc(i,j+1,m_u2))
354  enddo ; enddo
355  do j=js,je ; do i=is,ie
356  cs%heat_cyc(i,j,m_u2) = cs%heat_cyc(i,j,m_u2) + dt1_heat_rate * ( &
357  -cs%lam_cyc_heat*(cs%avg_SST_anom(i,j,m_u3) - cs%avg_SST_anom(i,j,m_u2)) + &
358  (g%IareaT(i,j) * ((flux_heat_x(i-1,j) - flux_heat_x(i,j)) + &
359  (flux_heat_y(i,j-1) - flux_heat_y(i,j))) ) )
360 
361  cs%precip_cyc(i,j,m_u2) = cs%precip_cyc(i,j,m_u2) + dt1_prec_rate * ( &
362  cs%lam_cyc_prec * (cs%avg_SSS_anom(i,j,m_u3) - cs%avg_SSS_anom(i,j,m_u2)) / &
363  (0.5*(cs%avg_SSS(i,j,m_u3) + cs%avg_SSS(i,j,m_u2))) + &
364  (g%IareaT(i,j) * ((flux_prec_x(i-1,j) - flux_prec_x(i,j)) + &
365  (flux_prec_y(i,j-1) - flux_prec_y(i,j))) ) )
366  enddo ; enddo
367  endif
368 
369  endif ! (CS%num_cycle > 0)
370 
371 end subroutine apply_ctrl_forcing
372 
373 !> This function maps rval into an integer in the range from 1 to num_period.
374 function periodic_int(rval, num_period) result (m)
375  real, intent(in) :: rval !< Input for mapping.
376  integer, intent(in) :: num_period !< Maximum output.
377  integer :: m !< Return value.
378 
379  m = floor(rval)
380  if (m <= 0) then
381  m = m + num_period * (1 + (abs(m) / num_period))
382  elseif (m > num_period) then
383  m = m - num_period * ((m-1) / num_period)
384  endif
385 end function
386 
387 !> This function shifts rval by an integer multiple of num_period so that
388 !! 0 <= val_out < num_period.
389 function periodic_real(rval, num_period) result(val_out)
390  real, intent(in) :: rval !< Input to be shifted into valid range.
391  integer, intent(in) :: num_period !< Maximum valid value.
392  real :: val_out !< Return value.
393  integer :: nshft
394 
395  if (rval < 0) then ; nshft = floor(abs(rval) / num_period) + 1
396  elseif (rval < num_period) then ; nshft = 0
397  else ; nshft = -1*floor(rval / num_period) ; endif
398 
399  val_out = rval + nshft * num_period
400 end function
401 
402 
403 !> This subroutine is used to allocate and register any fields in this module
404 !! that should be written to or read from the restart file.
405 subroutine register_ctrl_forcing_restarts(G, param_file, CS, restart_CS)
406  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure.
407  type(param_file_type), intent(in) :: param_file !< A structure indicating the
408  !! open file to parse for model
409  !! parameter values.
410  type(ctrl_forcing_cs), pointer :: cs !< A pointer that is set to point to the
411  !! control structure for this module.
412  type(mom_restart_cs), pointer :: restart_cs !< A pointer to the restart control structure.
413 
414  logical :: controlled, use_temperature
415  character (len=8) :: period_str
416  type(vardesc) :: vd
417  integer :: isd, ied, jsd, jed, isdb, iedb, jsdb, jedb
418  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
419  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
420 
421  if (associated(cs)) then
422  call mom_error(warning, "register_ctrl_forcing_restarts called "//&
423  "with an associated control structure.")
424  return
425  endif
426 
427  controlled = .false.
428  call read_param(param_file, "CONTROLLED_FORCING", controlled)
429  if (.not.controlled) return
430 
431  use_temperature = .true.
432  call read_param(param_file, "ENABLE_THERMODYNAMICS", use_temperature)
433  if (.not.use_temperature) call mom_error(fatal, &
434  "register_ctrl_forcing_restarts: CONTROLLED_FORCING only works with "//&
435  "ENABLE_THERMODYNAMICS defined.")
436 
437  allocate(cs)
438 
439  cs%do_integrated = .true. ; cs%num_cycle = 0
440  call read_param(param_file, "CTRL_FORCE_INTEGRATED", cs%do_integrated)
441  call read_param(param_file, "CTRL_FORCE_NUM_CYCLE", cs%num_cycle)
442 
443  if (cs%do_integrated) then
444  call safe_alloc_ptr(cs%heat_0,isd,ied,jsd,jed) ; cs%heat_0(:,:) = 0.0
445  call safe_alloc_ptr(cs%precip_0,isd,ied,jsd,jed) ; cs%precip_0(:,:) = 0.0
446  vd = var_desc("Ctrl_heat","W m-2","Control Integrative Heating",z_grid='1')
447  call register_restart_field(cs%heat_0, vd, .false., restart_cs)
448  vd = var_desc("Ctrl_precip","kg m-2 s-1","Control Integrative Precipitation",z_grid='1')
449  call register_restart_field(cs%precip_0, vd, .false., restart_cs)
450  endif
451 
452  if (cs%num_cycle > 0) then
453  write (period_str, '(i8)') cs%num_cycle
454  period_str = trim('p ')//trim(adjustl(period_str))
455  call safe_alloc_ptr(cs%heat_cyc,isd,ied,jsd,jed,cs%num_cycle) ; cs%heat_cyc(:,:,:) = 0.0
456  call safe_alloc_ptr(cs%precip_cyc,isd,ied,jsd,jed,cs%num_cycle) ; cs%precip_cyc(:,:,:) = 0.0
457  vd = var_desc("Ctrl_heat_cycle", "W m-2","Cyclical Control Heating",&
458  z_grid='1', t_grid=period_str)
459  call register_restart_field(cs%heat_cyc, vd, .false., restart_cs)
460  vd = var_desc("Ctrl_precip_cycle","kg m-2 s-1","Cyclical Control Precipitation", &
461  z_grid='1', t_grid=period_str)
462  call register_restart_field(cs%precip_cyc, vd, .false., restart_cs)
463 
464  call safe_alloc_ptr(cs%avg_time,cs%num_cycle) ; cs%avg_time(:) = 0.0
465  vd = var_desc("avg_time","sec","Cyclical accumulated averaging time", &
466  '1',z_grid='1',t_grid=period_str)
467  call register_restart_field(cs%avg_time, vd, .false., restart_cs)
468 
469  call safe_alloc_ptr(cs%avg_SST_anom,isd,ied,jsd,jed,cs%num_cycle) ; cs%avg_SST_anom(:,:,:) = 0.0
470  call safe_alloc_ptr(cs%avg_SSS_anom,isd,ied,jsd,jed,cs%num_cycle) ; cs%avg_SSS_anom(:,:,:) = 0.0
471  vd = var_desc("avg_SST_anom","deg C","Cyclical average SST Anomaly", &
472  z_grid='1',t_grid=period_str)
473  call register_restart_field(cs%avg_SST_anom, vd, .false., restart_cs)
474  vd = var_desc("avg_SSS_anom","g kg-1","Cyclical average SSS Anomaly", &
475  z_grid='1',t_grid=period_str)
476  call register_restart_field(cs%avg_SSS_anom, vd, .false., restart_cs)
477  endif
478 
479 end subroutine register_ctrl_forcing_restarts
480 
481 !> Set up this modules control structure.
482 subroutine controlled_forcing_init(Time, G, param_file, diag, CS)
483  type(time_type), intent(in) :: time !< The current model time.
484  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure.
485  type(param_file_type), intent(in) :: param_file !< A structure indicating the
486  !! open file to parse for model
487  !! parameter values.
488  type(diag_ctrl), target, intent(in) :: diag !< A structure that is used to regulate
489  !! diagnostic output.
490  type(ctrl_forcing_cs), pointer :: cs !< A pointer that is set to point to the
491  !! control structure for this module.
492  real :: smooth_len
493  logical :: do_integrated
494  integer :: num_cycle
495 ! This include declares and sets the variable "version".
496 #include "version_variable.h"
497  character(len=40) :: mdl = "MOM_controlled_forcing" ! This module's name.
498 
499  ! These should have already been called.
500  ! call read_param(param_file, "CTRL_FORCE_INTEGRATED", CS%do_integrated)
501  ! call read_param(param_file, "CTRL_FORCE_NUM_CYCLE", CS%num_cycle)
502 
503  if (associated(cs)) then
504  do_integrated = cs%do_integrated ; num_cycle = cs%num_cycle
505  else
506  do_integrated = .false. ; num_cycle = 0
507  endif
508 
509  ! Read all relevant parameters and write them to the model log.
510  call log_version(param_file, mdl, version, "")
511  call log_param(param_file, mdl, "CTRL_FORCE_INTEGRATED", do_integrated, &
512  "If true, use a PI controller to determine the surface "//&
513  "forcing that is consistent with the observed mean properties.", &
514  default=.false.)
515  call log_param(param_file, mdl, "CTRL_FORCE_NUM_CYCLE", num_cycle, &
516  "The number of cycles per year in the controlled forcing, "//&
517  "or 0 for no cyclic forcing.", default=0)
518 
519  if (.not.associated(cs)) return
520 
521  cs%diag => diag
522 
523  call get_param(param_file, mdl, "CTRL_FORCE_HEAT_INT_RATE", cs%heat_int_rate, &
524  "The integrated rate at which heat flux anomalies are "//&
525  "accumulated.", units="s-1", default=0.0)
526  call get_param(param_file, mdl, "CTRL_FORCE_PREC_INT_RATE", cs%prec_int_rate, &
527  "The integrated rate at which precipitation anomalies "//&
528  "are accumulated.", units="s-1", default=0.0)
529  call get_param(param_file, mdl, "CTRL_FORCE_HEAT_CYC_RATE", cs%heat_cyc_rate, &
530  "The integrated rate at which cyclical heat flux "//&
531  "anomalies are accumulated.", units="s-1", default=0.0)
532  call get_param(param_file, mdl, "CTRL_FORCE_PREC_CYC_RATE", cs%prec_cyc_rate, &
533  "The integrated rate at which cyclical precipitation "//&
534  "anomalies are accumulated.", units="s-1", default=0.0)
535  call get_param(param_file, mdl, "CTRL_FORCE_SMOOTH_LENGTH", smooth_len, &
536  "The length scales over which controlled forcing "//&
537  "anomalies are smoothed.", units="m", default=0.0)
538  call get_param(param_file, mdl, "CTRL_FORCE_LAMDA_HEAT", cs%lam_heat, &
539  "A constant of proportionality between SST anomalies "//&
540  "and controlling heat fluxes", "W m-2 K-1", default=0.0)
541  call get_param(param_file, mdl, "CTRL_FORCE_LAMDA_PREC", cs%lam_prec, &
542  "A constant of proportionality between SSS anomalies "//&
543  "(normalised by mean SSS) and controlling precipitation.", &
544  "kg m-2", default=0.0)
545  call get_param(param_file, mdl, "CTRL_FORCE_LAMDA_CYC_HEAT", cs%lam_cyc_heat, &
546  "A constant of proportionality between SST anomalies "//&
547  "and cyclical controlling heat fluxes", "W m-2 K-1", default=0.0)
548  call get_param(param_file, mdl, "CTRL_FORCE_LAMDA_CYC_PREC", cs%lam_cyc_prec, &
549  "A constant of proportionality between SSS anomalies "//&
550  "(normalised by mean SSS) and cyclical controlling "//&
551  "precipitation.", "kg m-2", default=0.0)
552 
553  cs%Len2 = smooth_len**2
554 
555 ! ### REPLACE THIS WITH ANY DIAGNOSTICS FROM THIS MODULE.
556 ! CS%id_taux = register_diag_field('ocean_model', 'taux', diag%axesu1, Time, &
557 ! 'Zonal Wind Stress', 'Pascal')
558 
559 end subroutine controlled_forcing_init
560 
561 !> Clean up this modules control structure.
562 subroutine controlled_forcing_end(CS)
563  type(ctrl_forcing_cs), pointer :: cs !< A pointer to the control structure
564  !! returned by a previous call to
565  !! controlled_forcing_init, it will be
566  !! deallocated here.
567 
568  if (associated(cs)) then
569  if (associated(cs%heat_0)) deallocate(cs%heat_0)
570  if (associated(cs%precip_0)) deallocate(cs%precip_0)
571  if (associated(cs%heat_cyc)) deallocate(cs%heat_cyc)
572  if (associated(cs%precip_cyc)) deallocate(cs%precip_cyc)
573  if (associated(cs%avg_SST_anom)) deallocate(cs%avg_SST_anom)
574  if (associated(cs%avg_SSS_anom)) deallocate(cs%avg_SSS_anom)
575  if (associated(cs%avg_SSS)) deallocate(cs%avg_SSS)
576 
577  deallocate(cs)
578  endif
579  cs => null()
580 
581 end subroutine controlled_forcing_end
582 
583 !> \namespace mom_controlled_forcing
584 !! *
585 !! By Robert Hallberg, July 2011 *
586 !! *
587 !! This program contains the subroutines that use control-theory *
588 !! to adjust the surface heat flux and precipitation, based on the *
589 !! time-mean or periodically (seasonally) varying anomalies from the *
590 !! observed state. The techniques behind this are described in *
591 !! Hallberg and Adcroft (2011, in prep.). *
592 !! *
593 !! Macros written all in capital letters are defined in MOM_memory.h. *
594 !! *
595 !! A small fragment of the grid is shown below: *
596 !! *
597 !! j+1 x ^ x ^ x At x: q *
598 !! j+1 > o > o > At ^: v, tauy *
599 !! j x ^ x ^ x At >: u, taux *
600 !! j > o > o > At o: h, fluxes. *
601 !! j-1 x ^ x ^ x *
602 !! i-1 i i+1 At x & ^: *
603 !! i i+1 At > & o: *
604 !! *
605 !! The boundaries always run through q grid points (x). *
606 end module mom_controlled_forcing
mom_time_manager
Wraps the FMS time manager functions.
Definition: MOM_time_manager.F90:2
mom_variables::surface
Pointers to various fields which may be used describe the surface state of MOM, and which will be ret...
Definition: MOM_variables.F90:38
mom_file_parser::log_version
An overloaded interface to log version information about modules.
Definition: MOM_file_parser.F90:109
mom_controlled_forcing
Use control-theory to adjust the surface heat flux and precipitation.
Definition: MOM_controlled_forcing.F90:7
mom_diag_mediator
The subroutines here provide convenient wrappers to the fms diag_manager interfaces with additional d...
Definition: MOM_diag_mediator.F90:3
mom_domains::pass_var
Do a halo update on an array.
Definition: MOM_domains.F90:49
mom_file_parser::param_file_type
A structure that can be parsed to read and document run-time parameters.
Definition: MOM_file_parser.F90:54
mom_file_parser::get_param
An overloaded interface to read and log the values of various types of parameters.
Definition: MOM_file_parser.F90:102
mom_io
This module contains I/O framework code.
Definition: MOM_io.F90:2
mom_restart::mom_restart_cs
A restart registry and the control structure for restarts.
Definition: MOM_restart.F90:72
mom_diag_mediator::post_data
Make a diagnostic available for averaging or output.
Definition: MOM_diag_mediator.F90:70
mom_domains::pass_vector
Do a halo update on a pair of arrays representing the two components of a vector.
Definition: MOM_domains.F90:54
mom_forcing_type
This module implements boundary forcing for MOM6.
Definition: MOM_forcing_type.F90:2
mom_restart
The MOM6 facility for reading and writing restart files, and querying what has been read.
Definition: MOM_restart.F90:2
mom_domains
Describes the decomposed MOM domain and has routines for communications across PEs.
Definition: MOM_domains.F90:2
mom_variables
Provides transparent structures with groups of MOM6 variables and supporting routines.
Definition: MOM_variables.F90:2
mom_file_parser
The MOM6 facility to parse input files for runtime parameters.
Definition: MOM_file_parser.F90:2
mom_grid
Provides the ocean grid type.
Definition: MOM_grid.F90:2
mom_restart::register_restart_field
Register fields for restarts.
Definition: MOM_restart.F90:107
mom_forcing_type::forcing
Structure that contains pointers to the boundary forcing used to drive the liquid ocean simulated by ...
Definition: MOM_forcing_type.F90:49
mom_io::vardesc
Type for describing a variable, typically a tracer.
Definition: MOM_io.F90:53
mom_file_parser::log_param
An overloaded interface to log the values of various types of parameters.
Definition: MOM_file_parser.F90:96
mom_controlled_forcing::ctrl_forcing_cs
Control structure for MOM_controlled_forcing.
Definition: MOM_controlled_forcing.F90:33
mom_error_handler
Routines for error handling and I/O management.
Definition: MOM_error_handler.F90:2
mom_grid::ocean_grid_type
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:25
mom_diag_mediator::diag_ctrl
The following data type a list of diagnostic fields an their variants, as well as variables that cont...
Definition: MOM_diag_mediator.F90:239
mom_file_parser::read_param
An overloaded interface to read various types of parameters.
Definition: MOM_file_parser.F90:90