MOM6
MOM_PointAccel.F90
1 !> Debug accelerations at a given point
2 !!
3 !! The two subroutines in this file write out all of the terms
4 !! in the u- or v-momentum balance at a given point. Usually
5 !! these subroutines are called after the velocities exceed some
6 !! threshold, in order to determine which term is culpable.
7 !! often this is done for debugging purposes.
9 
10 ! This file is part of MOM6. See LICENSE.md for the license.
11 
12 use mom_diag_mediator, only : diag_ctrl
13 use mom_domains, only : pe_here
14 use mom_error_handler, only : mom_error, note
16 use mom_get_input, only : directories
17 use mom_grid, only : ocean_grid_type
18 use mom_io, only : open_file
19 use mom_io, only : append_file, ascii_file, multiple, single_file
20 use mom_time_manager, only : time_type, get_time, get_date, set_date, operator(-)
24 
25 implicit none ; private
26 
27 #include <MOM_memory.h>
28 
29 public write_u_accel, write_v_accel, pointaccel_init
30 
31 !> The control structure for the MOM_PointAccel module
32 type, public :: pointaccel_cs ; private
33  character(len=200) :: u_trunc_file !< The complete path to the file in which a column's worth of
34  !! u-accelerations are written if u-velocity truncations occur.
35  character(len=200) :: v_trunc_file !< The complete path to the file in which a column's worth of
36  !! v-accelerations are written if v-velocity truncations occur.
37  integer :: u_file !< The unit number for an opened u-truncation files, or -1 if it has not yet been opened.
38  integer :: v_file !< The unit number for an opened v-truncation files, or -1 if it has not yet been opened.
39  integer :: cols_written !< The number of columns whose output has been
40  !! written by this PE during the current run.
41  integer :: max_writes !< The maximum number of times any PE can write out
42  !! a column's worth of accelerations during a run.
43  type(time_type), pointer :: time => null() !< A pointer to the ocean model's clock.
44  type(diag_ctrl), pointer :: diag => null() !< A structure that is used to
45  !! regulate the timing of diagnostic output.
46 ! The following are pointers to many of the state variables and accelerations
47 ! that are used to step the physical model forward. They all use the same
48 ! names as the variables they point to in MOM.F90
49  real, pointer, dimension(:,:,:) :: &
50  u_av => null(), & !< Time average u-velocity [m s-1].
51  v_av => null(), & !< Time average velocity [m s-1].
52  u_prev => null(), & !< Previous u-velocity [m s-1].
53  v_prev => null(), & !< Previous v-velocity [m s-1].
54  t => null(), & !< Temperature [degC].
55  s => null(), & !< Salinity [ppt].
56  u_accel_bt => null(), & !< Barotropic u-acclerations [m s-2]
57  v_accel_bt => null() !< Barotropic v-acclerations [m s-2]
58  real, pointer, dimension(:,:,:) :: pbce => null() !< pbce times eta gives the baroclinic
59  !! pressure anomaly in each layer due to free surface height anomalies
60  !! [m2 s-2 H-1 ~> m s-2 or m4 kg-1 s-2].
61 
62 end type pointaccel_cs
63 
64 contains
65 
66 !> This subroutine writes to an output file all of the accelerations
67 !! that have been applied to a column of zonal velocities over the
68 !! previous timestep. This subroutine is called from vertvisc.
69 subroutine write_u_accel(I, j, um, hin, ADp, CDp, dt, G, GV, US, CS, vel_rpt, str, a, hv)
70  integer, intent(in) :: i !< The zonal index of the column to be documented.
71  integer, intent(in) :: j !< The meridional index of the column to be documented.
72  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure.
73  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure.
74  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
75  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
76  intent(in) :: um !< The new zonal velocity [m s-1].
77  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
78  intent(in) :: hin !< The layer thickness [H ~> m or kg m-2].
79  type(accel_diag_ptrs), intent(in) :: adp !< A structure pointing to the various
80  !! accelerations in the momentum equations.
81  type(cont_diag_ptrs), intent(in) :: cdp !< A structure with pointers to various terms
82  !! in the continuity equations.
83  real, intent(in) :: dt !< The ocean dynamics time step [s].
84  type(pointaccel_cs), pointer :: cs !< The control structure returned by a previous
85  !! call to PointAccel_init.
86  real, intent(in) :: vel_rpt !< The velocity magnitude that triggers a report [m s-1].
87  real, optional, intent(in) :: str !< The surface wind stress integrated over a time
88  !! step divided by the Boussinesq density [m2 s-1].
89  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
90  optional, intent(in) :: a !< The layer coupling coefficients from vertvisc [Z s-1 ~> m s-1].
91  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
92  optional, intent(in) :: hv !< The layer thicknesses at velocity grid points,
93  !! from vertvisc [H ~> m or kg m-2].
94  ! Local variables
95  real :: f_eff, cfl
96  real :: angstrom
97  real :: truncvel, du
98  real :: inorm(szk_(g))
99  real :: e(szk_(g)+1)
100  real :: h_scale, uh_scale
101  integer :: yr, mo, day, hr, minute, sec, yearday
102  integer :: k, ks, ke
103  integer :: nz
104  logical :: do_k(szk_(g)+1)
105  logical :: prev_avail
106  integer :: file
107 
108  angstrom = gv%Angstrom_H + gv%H_subroundoff
109  h_scale = gv%H_to_m ; uh_scale = gv%H_to_m
110 
111 ! if (.not.associated(CS)) return
112  nz = g%ke
113  if (cs%cols_written < cs%max_writes) then
114  cs%cols_written = cs%cols_written + 1
115 
116  ks = 1 ; ke = nz
117  do_k(:) = .false.
118 
119  ! Open up the file for output if this is the first call.
120  if (cs%u_file < 0) then
121  if (len_trim(cs%u_trunc_file) < 1) return
122  call open_file(cs%u_file, trim(cs%u_trunc_file), action=append_file, &
123  form=ascii_file, threading=multiple, fileset=single_file)
124  if (cs%u_file < 0) then
125  call mom_error(note, 'Unable to open file '//trim(cs%u_trunc_file)//'.')
126  return
127  endif
128  endif
129  file = cs%u_file
130 
131  prev_avail = (associated(cs%u_prev) .and. associated(cs%v_prev))
132 
133  ! Determine which layers to write out accelerations for.
134  do k=1,nz
135  if (((max(cs%u_av(i,j,k),um(i,j,k)) >= vel_rpt) .or. &
136  (min(cs%u_av(i,j,k),um(i,j,k)) <= -vel_rpt)) .and. &
137  ((hin(i,j,k) + hin(i+1,j,k)) > 3.0*angstrom)) exit
138  enddo
139  ks = k
140  do k=nz,1,-1
141  if (((max(cs%u_av(i,j,k), um(i,j,k)) >= vel_rpt) .or. &
142  (min(cs%u_av(i,j,k), um(i,j,k)) <= -vel_rpt)) .and. &
143  ((hin(i,j,k) + hin(i+1,j,k)) > 3.0*angstrom)) exit
144  enddo
145  ke = k
146  if (ke < ks) then
147  ks = 1; ke = nz; write(file,'("U: Unable to set ks & ke.")')
148  endif
149 
150  call get_date(cs%Time, yr, mo, day, hr, minute, sec)
151  call get_time((cs%Time - set_date(yr, 1, 1, 0, 0, 0)), sec, yearday)
152  write (file,'(/,"--------------------------")')
153  write (file,'(/,"Time ",i5,i4,F6.2," U-velocity violation at ",I4,": ",2(I3), &
154  & " (",F7.2," E "F7.2," N) Layers ",I3," to ",I3,". dt = ",1PG10.4)') &
155  yr, yearday, (real(sec)/3600.0), pe_here(), i, j, &
156  g%geoLonCu(i,j), g%geoLatCu(i,j), ks, ke, dt
157 
158  if (ks <= gv%nk_rho_varies) ks = 1
159  do k=ks,ke
160  if ((hin(i,j,k) + hin(i+1,j,k)) > 3.0*angstrom) do_k(k) = .true.
161  enddo
162 
163  write(file,'(/,"Layers:",$)')
164  do k=ks,ke ; if (do_k(k)) write(file,'(I10," ",$)') (k); enddo
165  write(file,'(/,"u(m): ",$)')
166  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') (um(i,j,k)); enddo
167  if (prev_avail) then
168  write(file,'(/,"u(mp): ",$)')
169  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') (cs%u_prev(i,j,k)); enddo
170  endif
171  write(file,'(/,"u(3): ",$)')
172  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') (cs%u_av(i,j,k)); enddo
173 
174  write(file,'(/,"CFL u: ",$)')
175  do k=ks,ke ; if (do_k(k)) then
176  cfl = abs(um(i,j,k)) * dt * g%dy_Cu(i,j)
177  if (um(i,j,k) < 0.0) then ; cfl = cfl * g%IareaT(i+1,j)
178  else ; cfl = cfl * g%IareaT(i,j) ; endif
179  write(file,'(ES10.3," ",$)') cfl
180  endif ; enddo
181  write(file,'(/,"CFL0 u:",$)')
182  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') &
183  abs(um(i,j,k)) * dt * g%IdxCu(i,j) ; enddo
184 
185  if (prev_avail) then
186  write(file,'(/,"du: ",$)')
187  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') &
188  ((um(i,j,k)-cs%u_prev(i,j,k))); enddo
189  endif
190  write(file,'(/,"CAu: ",$)')
191  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') (dt*adp%CAu(i,j,k)); enddo
192  write(file,'(/,"PFu: ",$)')
193  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') (dt*adp%PFu(i,j,k)); enddo
194  write(file,'(/,"diffu: ",$)')
195  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') (dt*us%s_to_T*adp%diffu(i,j,k)); enddo
196 
197  if (associated(adp%gradKEu)) then
198  write(file,'(/,"KEu: ",$)')
199  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') &
200  (dt*adp%gradKEu(i,j,k)); enddo
201  endif
202  if (associated(adp%rv_x_v)) then
203  write(file,'(/,"Coru: ",$)')
204  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') &
205  dt*(adp%CAu(i,j,k)-adp%rv_x_v(i,j,k)); enddo
206  endif
207  if (associated(adp%du_dt_visc)) then
208  write(file,'(/,"ubv: ",$)')
209  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') &
210  (um(i,j,k)-dt*adp%du_dt_visc(i,j,k)); enddo
211  write(file,'(/,"duv: ",$)')
212  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') &
213  (dt*adp%du_dt_visc(i,j,k)); enddo
214  endif
215  if (associated(adp%du_other)) then
216  write(file,'(/,"du_other: ",$)')
217  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') &
218  (adp%du_other(i,j,k)); enddo
219  endif
220  if (present(a)) then
221  write(file,'(/,"a: ",$)')
222  do k=ks,ke+1 ; if (do_k(k)) write(file,'(ES10.3," ",$)') a(i,j,k)*us%Z_to_m*dt; enddo
223  endif
224  if (present(hv)) then
225  write(file,'(/,"hvel: ",$)')
226  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') hv(i,j,k); enddo
227  endif
228  write(file,'(/,"Stress: ",ES10.3)') str
229 
230  if (associated(cs%u_accel_bt)) then
231  write(file,'("dubt: ",$)')
232  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') &
233  (dt*cs%u_accel_bt(i,j,k)) ; enddo
234  write(file,'(/)')
235  endif
236 
237  write(file,'(/,"h--: ",$)')
238  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') (h_scale*hin(i,j-1,k)); enddo
239  write(file,'(/,"h+-: ",$)')
240  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') (h_scale*hin(i+1,j-1,k)); enddo
241  write(file,'(/,"h-0: ",$)')
242  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') (h_scale*hin(i,j,k)); enddo
243  write(file,'(/,"h+0: ",$)')
244  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') (h_scale*hin(i+1,j,k)); enddo
245  write(file,'(/,"h-+: ",$)')
246  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') (h_scale*hin(i,j+1,k)); enddo
247  write(file,'(/,"h++: ",$)')
248  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') (h_scale*hin(i+1,j+1,k)); enddo
249 
250 
251  e(nz+1) = -us%Z_to_m*g%bathyT(i,j)
252  do k=nz,1,-1 ; e(k) = e(k+1) + h_scale*hin(i,j,k) ; enddo
253  write(file,'(/,"e-: ",$)')
254  write(file,'(ES10.3," ",$)') e(ks)
255  do k=ks+1,ke+1 ; if (do_k(k-1)) write(file,'(ES10.3," ",$)') e(k); enddo
256 
257  e(nz+1) = -us%Z_to_m*g%bathyT(i+1,j)
258  do k=nz,1,-1 ; e(k) = e(k+1) + h_scale*hin(i+1,j,k) ; enddo
259  write(file,'(/,"e+: ",$)')
260  write(file,'(ES10.3," ",$)') e(ks)
261  do k=ks+1,ke+1 ; if (do_k(k-1)) write(file,'(ES10.3," ",$)') e(k) ; enddo
262  if (associated(cs%T)) then
263  write(file,'(/,"T-: ",$)')
264  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') cs%T(i,j,k); enddo
265  write(file,'(/,"T+: ",$)')
266  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') cs%T(i+1,j,k); enddo
267  endif
268  if (associated(cs%S)) then
269  write(file,'(/,"S-: ",$)')
270  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') cs%S(i,j,k); enddo
271  write(file,'(/,"S+: ",$)')
272  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') cs%S(i+1,j,k); enddo
273  endif
274 
275  if (prev_avail) then
276  write(file,'(/,"v--: ",$)')
277  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') (cs%v_prev(i,j-1,k)); enddo
278  write(file,'(/,"v-+: ",$)')
279  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') (cs%v_prev(i,j,k)); enddo
280  write(file,'(/,"v+-: ",$)')
281  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') (cs%v_prev(i+1,j-1,k)); enddo
282  write(file,'(/,"v++: ",$)')
283  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') (cs%v_prev(i+1,j,k)); enddo
284  endif
285 
286  write(file,'(/,"vh--: ",$)')
287  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') &
288  (uh_scale*cdp%vh(i,j-1,k)*g%IdxCv(i,j-1)); enddo
289  write(file,'(/," vhC--:",$)')
290  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') &
291  (0.5*cs%v_av(i,j-1,k)*h_scale*(hin(i,j-1,k) + hin(i,j,k))); enddo
292  if (prev_avail) then
293  write(file,'(/," vhCp--:",$)')
294  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') &
295  (0.5*cs%v_prev(i,j-1,k)*h_scale*(hin(i,j-1,k) + hin(i,j,k))); enddo
296  endif
297 
298  write(file,'(/,"vh-+: ",$)')
299  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') &
300  (uh_scale*cdp%vh(i,j,k)*g%IdxCv(i,j)); enddo
301  write(file,'(/," vhC-+:",$)')
302  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') &
303  (0.5*cs%v_av(i,j,k)*h_scale*(hin(i,j,k) + hin(i,j+1,k))); enddo
304  if (prev_avail) then
305  write(file,'(/," vhCp-+:",$)')
306  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') &
307  (0.5*cs%v_prev(i,j,k)*h_scale*(hin(i,j,k) + hin(i,j+1,k))); enddo
308  endif
309 
310  write(file,'(/,"vh+-: ",$)')
311  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') &
312  (uh_scale*cdp%vh(i+1,j-1,k)*g%IdxCv(i+1,j-1)); enddo
313  write(file,'(/," vhC+-:",$)')
314  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') &
315  (0.5*cs%v_av(i+1,j-1,k)*h_scale*(hin(i+1,j-1,k) + hin(i+1,j,k))); enddo
316  if (prev_avail) then
317  write(file,'(/," vhCp+-:",$)')
318  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') &
319  (0.5*cs%v_prev(i+1,j-1,k)*h_scale*(hin(i+1,j-1,k) + hin(i+1,j,k))); enddo
320  endif
321 
322  write(file,'(/,"vh++: ",$)')
323  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') &
324  (uh_scale*cdp%vh(i+1,j,k)*g%IdxCv(i+1,j)); enddo
325  write(file,'(/," vhC++:",$)')
326  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') &
327  (0.5*cs%v_av(i+1,j,k)*h_scale*(hin(i+1,j,k) + hin(i+1,j+1,k))); enddo
328  if (prev_avail) then
329  write(file,'(/," vhCp++:",$)')
330  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') &
331  (0.5*cs%v_av(i+1,j,k)*h_scale*(hin(i+1,j,k) + hin(i+1,j+1,k))); enddo
332  endif
333 
334  write(file,'(/,"D: ",2(ES10.3))') us%Z_to_m*g%bathyT(i,j),us%Z_to_m*g%bathyT(i+1,j)
335 
336  ! From here on, the normalized accelerations are written.
337  if (prev_avail) then
338  do k=ks,ke
339  du = um(i,j,k)-cs%u_prev(i,j,k)
340  if (abs(du) < 1.0e-6) du = 1.0e-6
341  inorm(k) = 1.0 / du
342  enddo
343 
344  write(file,'(2/,"Norm: ",$)')
345  do k=ks,ke ; if (do_k(k)) write(file,'(F10.6," ",$)') (1.0/inorm(k)); enddo
346 
347  write(file,'(/,"du: ",$)')
348  do k=ks,ke ; if (do_k(k)) write(file,'(F10.6," ",$)') &
349  ((um(i,j,k)-cs%u_prev(i,j,k))*inorm(k)); enddo
350 
351  write(file,'(/,"CAu: ",$)')
352  do k=ks,ke ; if (do_k(k)) write(file,'(F10.6," ",$)') &
353  (dt*adp%CAu(i,j,k)*inorm(k)); enddo
354 
355  write(file,'(/,"PFu: ",$)')
356  do k=ks,ke ; if (do_k(k)) write(file,'(F10.6," ",$)') &
357  (dt*adp%PFu(i,j,k)*inorm(k)); enddo
358 
359  write(file,'(/,"diffu: ",$)')
360  do k=ks,ke ; if (do_k(k)) write(file,'(F10.6," ",$)') &
361  (dt*us%s_to_T*adp%diffu(i,j,k)*inorm(k)); enddo
362 
363  if (associated(adp%gradKEu)) then
364  write(file,'(/,"KEu: ",$)')
365  do k=ks,ke ; if (do_k(k)) write(file,'(F10.6," ",$)') &
366  (dt*adp%gradKEu(i,j,k)*inorm(k)); enddo
367  endif
368  if (associated(adp%rv_x_v)) then
369  write(file,'(/,"Coru: ",$)')
370  do k=ks,ke ; if (do_k(k)) write(file,'(F10.6," ",$)') &
371  dt*(adp%CAu(i,j,k)-adp%rv_x_v(i,j,k))*inorm(k); enddo
372  endif
373  if (associated(adp%du_dt_visc)) then
374  write(file,'(/,"duv: ",$)')
375  do k=ks,ke ; if (do_k(k)) write(file,'(F10.6," ",$)') &
376  (dt*adp%du_dt_visc(i,j,k))*inorm(k); enddo
377  endif
378  if (associated(adp%du_other)) then
379  write(file,'(/,"du_other: ",$)')
380  do k=ks,ke ; if (do_k(k)) write(file,'(F10.6," ",$)') &
381  (adp%du_other(i,j,k))*inorm(k); enddo
382  endif
383  if (associated(cs%u_accel_bt)) then
384  write(file,'(/,"dubt: ",$)')
385  do k=ks,ke ; if (do_k(k)) write(file,'(F10.6," ",$)') &
386  (dt*cs%u_accel_bt(i,j,k)*inorm(k)) ; enddo
387  endif
388  endif
389 
390  write(file,'(2/)')
391 
392  call flush(file)
393  endif
394 
395 end subroutine write_u_accel
396 
397 !> This subroutine writes to an output file all of the accelerations
398 !! that have been applied to a column of meridional velocities over
399 !! the previous timestep. This subroutine is called from vertvisc.
400 subroutine write_v_accel(i, J, vm, hin, ADp, CDp, dt, G, GV, US, CS, vel_rpt, str, a, hv)
401  integer, intent(in) :: i !< The zonal index of the column to be documented.
402  integer, intent(in) :: j !< The meridional index of the column to be documented.
403  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure.
404  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure.
405  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
406  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
407  intent(in) :: vm !< The new meridional velocity [m s-1].
408  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
409  intent(in) :: hin !< The layer thickness [H ~> m or kg m-2].
410  type(accel_diag_ptrs), intent(in) :: adp !< A structure pointing to the various
411  !! accelerations in the momentum equations.
412  type(cont_diag_ptrs), intent(in) :: cdp !< A structure with pointers to various terms in
413  !! the continuity equations.
414  real, intent(in) :: dt !< The ocean dynamics time step [s].
415  type(pointaccel_cs), pointer :: cs !< The control structure returned by a previous
416  !! call to PointAccel_init.
417  real, intent(in) :: vel_rpt !< The velocity magnitude that triggers a report [m s-1].
418  real, optional, intent(in) :: str !< The surface wind stress integrated over a time
419  !! step divided by the Boussinesq density [m2 s-1].
420  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
421  optional, intent(in) :: a !< The layer coupling coefficients from vertvisc [Z s-1 ~> m s-1].
422  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
423  optional, intent(in) :: hv !< The layer thicknesses at velocity grid points,
424  !! from vertvisc [H ~> m or kg m-2].
425  ! Local variables
426  real :: f_eff, cfl
427  real :: angstrom
428  real :: truncvel, dv
429  real :: inorm(szk_(g))
430  real :: e(szk_(g)+1)
431  real :: h_scale, uh_scale
432  integer :: yr, mo, day, hr, minute, sec, yearday
433  integer :: k, ks, ke
434  integer :: nz
435  logical :: do_k(szk_(g)+1)
436  logical :: prev_avail
437  integer :: file
438 
439  angstrom = gv%Angstrom_H + gv%H_subroundoff
440  h_scale = gv%H_to_m ; uh_scale = gv%H_to_m
441 
442 ! if (.not.associated(CS)) return
443  nz = g%ke
444  if (cs%cols_written < cs%max_writes) then
445  cs%cols_written = cs%cols_written + 1
446 
447  ks = 1 ; ke = nz
448  do_k(:) = .false.
449 
450  ! Open up the file for output if this is the first call.
451  if (cs%v_file < 0) then
452  if (len_trim(cs%v_trunc_file) < 1) return
453  call open_file(cs%v_file, trim(cs%v_trunc_file), action=append_file, &
454  form=ascii_file, threading=multiple, fileset=single_file)
455  if (cs%v_file < 0) then
456  call mom_error(note, 'Unable to open file '//trim(cs%v_trunc_file)//'.')
457  return
458  endif
459  endif
460  file = cs%v_file
461 
462  prev_avail = (associated(cs%u_prev) .and. associated(cs%v_prev))
463 
464  do k=1,nz
465  if (((max(cs%v_av(i,j,k), vm(i,j,k)) >= vel_rpt) .or. &
466  (min(cs%v_av(i,j,k), vm(i,j,k)) <= -vel_rpt)) .and. &
467  ((hin(i,j,k) + hin(i,j+1,k)) > 3.0*angstrom)) exit
468  enddo
469  ks = k
470  do k=nz,1,-1
471  if (((max(cs%v_av(i,j,k), vm(i,j,k)) >= vel_rpt) .or. &
472  (min(cs%v_av(i,j,k), vm(i,j,k)) <= -vel_rpt)) .and. &
473  ((hin(i,j,k) + hin(i,j+1,k)) > 3.0*angstrom)) exit
474  enddo
475  ke = k
476  if (ke < ks) then
477  ks = 1; ke = nz; write(file,'("V: Unable to set ks & ke.")')
478  endif
479 
480  call get_date(cs%Time, yr, mo, day, hr, minute, sec)
481  call get_time((cs%Time - set_date(yr, 1, 1, 0, 0, 0)), sec, yearday)
482  write (file,'(/,"--------------------------")')
483  write (file,'(/,"Time ",i5,i4,F6.2," V-velocity violation at ",I4,": ",2(I3), &
484  & " (",F7.2," E ",F7.2," N) Layers ",I3," to ",I3,". dt = ",1PG10.4)') &
485  yr, yearday, (real(sec)/3600.0), pe_here(), i, j, &
486  g%geoLonCv(i,j), g%geoLatCv(i,j), ks, ke, dt
487 
488  if (ks <= gv%nk_rho_varies) ks = 1
489  do k=ks,ke
490  if ((hin(i,j,k) + hin(i,j+1,k)) > 3.0*angstrom) do_k(k) = .true.
491  enddo
492 
493  write(file,'(/,"Layers:",$)')
494  do k=ks,ke ; if (do_k(k)) write(file,'(I10," ",$)') (k); enddo
495  write(file,'(/,"v(m): ",$)')
496  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') (vm(i,j,k)); enddo
497 
498  if (prev_avail) then
499  write(file,'(/,"v(mp): ",$)')
500  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') (cs%v_prev(i,j,k)); enddo
501  endif
502 
503  write(file,'(/,"v(3): ",$)')
504  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') (cs%v_av(i,j,k)); enddo
505  write(file,'(/,"CFL v: ",$)')
506  do k=ks,ke ; if (do_k(k)) then
507  cfl = abs(vm(i,j,k)) * dt * g%dx_Cv(i,j)
508  if (vm(i,j,k) < 0.0) then ; cfl = cfl * g%IareaT(i,j+1)
509  else ; cfl = cfl * g%IareaT(i,j) ; endif
510  write(file,'(ES10.3," ",$)') cfl
511  endif ; enddo
512  write(file,'(/,"CFL0 v:",$)')
513  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') &
514  abs(vm(i,j,k)) * dt * g%IdyCv(i,j) ; enddo
515 
516  if (prev_avail) then
517  write(file,'(/,"dv: ",$)')
518  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') &
519  ((vm(i,j,k)-cs%v_prev(i,j,k))); enddo
520  endif
521 
522  write(file,'(/,"CAv: ",$)')
523  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') (dt*adp%CAv(i,j,k)); enddo
524 
525  write(file,'(/,"PFv: ",$)')
526  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') (dt*adp%PFv(i,j,k)); enddo
527 
528  write(file,'(/,"diffv: ",$)')
529  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') (dt*us%s_to_T*adp%diffv(i,j,k)); enddo
530 
531  if (associated(adp%gradKEv)) then
532  write(file,'(/,"KEv: ",$)')
533  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') &
534  (dt*adp%gradKEv(i,j,k)); enddo
535  endif
536  if (associated(adp%rv_x_u)) then
537  write(file,'(/,"Corv: ",$)')
538  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') &
539  dt*(adp%CAv(i,j,k)-adp%rv_x_u(i,j,k)); enddo
540  endif
541  if (associated(adp%dv_dt_visc)) then
542  write(file,'(/,"vbv: ",$)')
543  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') &
544  (vm(i,j,k)-dt*adp%dv_dt_visc(i,j,k)); enddo
545 
546  write(file,'(/,"dvv: ",$)')
547  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') &
548  (dt*adp%dv_dt_visc(i,j,k)); enddo
549  endif
550  if (associated(adp%dv_other)) then
551  write(file,'(/,"dv_other: ",$)')
552  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') &
553  (adp%dv_other(i,j,k)); enddo
554  endif
555  if (present(a)) then
556  write(file,'(/,"a: ",$)')
557  do k=ks,ke+1 ; if (do_k(k)) write(file,'(ES10.3," ",$)') a(i,j,k)*us%Z_to_m*dt; enddo
558  endif
559  if (present(hv)) then
560  write(file,'(/,"hvel: ",$)')
561  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') hv(i,j,k); enddo
562  endif
563  write(file,'(/,"Stress: ",ES10.3)') str
564 
565  if (associated(cs%v_accel_bt)) then
566  write(file,'("dvbt: ",$)')
567  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') &
568  (dt*cs%v_accel_bt(i,j,k)) ; enddo
569  write(file,'(/)')
570  endif
571 
572  write(file,'("h--: ",$)')
573  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') h_scale*hin(i-1,j,k); enddo
574  write(file,'(/,"h0-: ",$)')
575  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') h_scale*hin(i,j,k); enddo
576  write(file,'(/,"h+-: ",$)')
577  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') h_scale*hin(i+1,j,k); enddo
578  write(file,'(/,"h-+: ",$)')
579  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') h_scale*hin(i-1,j+1,k); enddo
580  write(file,'(/,"h0+: ",$)')
581  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') h_scale*hin(i,j+1,k); enddo
582  write(file,'(/,"h++: ",$)')
583  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') h_scale*hin(i+1,j+1,k); enddo
584 
585  e(nz+1) = -us%Z_to_m*g%bathyT(i,j)
586  do k=nz,1,-1 ; e(k) = e(k+1) + h_scale*hin(i,j,k); enddo
587  write(file,'(/,"e-: ",$)')
588  write(file,'(ES10.3," ",$)') e(ks)
589  do k=ks+1,ke+1 ; if (do_k(k-1)) write(file,'(ES10.3," ",$)') e(k); enddo
590 
591  e(nz+1) = -us%Z_to_m*g%bathyT(i,j+1)
592  do k=nz,1,-1 ; e(k) = e(k+1) + h_scale*hin(i,j+1,k) ; enddo
593  write(file,'(/,"e+: ",$)')
594  write(file,'(ES10.3," ",$)') e(ks)
595  do k=ks+1,ke+1 ; if (do_k(k-1)) write(file,'(ES10.3," ",$)') e(k); enddo
596  if (associated(cs%T)) then
597  write(file,'(/,"T-: ",$)')
598  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') cs%T(i,j,k); enddo
599  write(file,'(/,"T+: ",$)')
600  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') cs%T(i,j+1,k); enddo
601  endif
602  if (associated(cs%S)) then
603  write(file,'(/,"S-: ",$)')
604  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') cs%S(i,j,k); enddo
605  write(file,'(/,"S+: ",$)')
606  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') cs%S(i,j+1,k); enddo
607  endif
608 
609  if (prev_avail) then
610  write(file,'(/,"u--: ",$)')
611  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') cs%u_prev(i-1,j,k); enddo
612  write(file,'(/,"u-+: ",$)')
613  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') cs%u_prev(i-1,j+1,k); enddo
614  write(file,'(/,"u+-: ",$)')
615  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') cs%u_prev(i,j,k); enddo
616  write(file,'(/,"u++: ",$)')
617  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') cs%u_prev(i,j+1,k); enddo
618  endif
619 
620  write(file,'(/,"uh--: ",$)')
621  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') &
622  (uh_scale*cdp%uh(i-1,j,k)*g%IdyCu(i-1,j)); enddo
623  write(file,'(/," uhC--: ",$)')
624  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') &
625  (cs%u_av(i-1,j,k) * h_scale*0.5*(hin(i-1,j,k) + hin(i,j,k))); enddo
626  if (prev_avail) then
627  write(file,'(/," uhCp--:",$)')
628  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') &
629  (cs%u_prev(i-1,j,k) * h_scale*0.5*(hin(i-1,j,k) + hin(i,j,k))); enddo
630  endif
631 
632  write(file,'(/,"uh-+: ",$)')
633  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') &
634  (uh_scale*cdp%uh(i-1,j+1,k)*g%IdyCu(i-1,j+1)); enddo
635  write(file,'(/," uhC-+: ",$)')
636  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') &
637  (cs%u_av(i-1,j+1,k) * h_scale*0.5*(hin(i-1,j+1,k) + hin(i,j+1,k))); enddo
638  if (prev_avail) then
639  write(file,'(/," uhCp-+:",$)')
640  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') &
641  (cs%u_prev(i-1,j+1,k) * h_scale*0.5*(hin(i-1,j+1,k) + hin(i,j+1,k))); enddo
642  endif
643 
644  write(file,'(/,"uh+-: ",$)')
645  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') &
646  (uh_scale*cdp%uh(i,j,k)*g%IdyCu(i,j)); enddo
647  write(file,'(/," uhC+-: ",$)')
648  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') &
649  (cs%u_av(i,j,k) * h_scale*0.5*(hin(i,j,k) + hin(i+1,j,k))); enddo
650  if (prev_avail) then
651  write(file,'(/," uhCp+-:",$)')
652  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') &
653  (cs%u_prev(i,j,k) * h_scale*0.5*(hin(i,j,k) + hin(i+1,j,k))); enddo
654  endif
655 
656  write(file,'(/,"uh++: ",$)')
657  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') &
658  (uh_scale*cdp%uh(i,j+1,k)*g%IdyCu(i,j+1)); enddo
659  write(file,'(/," uhC++: ",$)')
660  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') &
661  (cs%u_av(i,j+1,k) * 0.5*h_scale*(hin(i,j+1,k) + hin(i+1,j+1,k))); enddo
662  if (prev_avail) then
663  write(file,'(/," uhCp++:",$)')
664  do k=ks,ke ; if (do_k(k)) write(file,'(ES10.3," ",$)') &
665  (cs%u_prev(i,j+1,k) * h_scale*0.5*(hin(i,j+1,k) + hin(i+1,j+1,k))); enddo
666  endif
667 
668  write(file,'(/,"D: ",2(ES10.3))') us%Z_to_m*g%bathyT(i,j),us%Z_to_m*g%bathyT(i,j+1)
669 
670  ! From here on, the normalized accelerations are written.
671  if (prev_avail) then
672  do k=ks,ke
673  dv = vm(i,j,k)-cs%v_prev(i,j,k)
674  if (abs(dv) < 1.0e-6) dv = 1.0e-6
675  inorm(k) = 1.0 / dv
676  enddo
677 
678  write(file,'(2/,"Norm: ",$)')
679  do k=ks,ke ; if (do_k(k)) write(file,'(F10.6," ",$)') (1.0/inorm(k)); enddo
680  write(file,'(/,"dv: ",$)')
681  do k=ks,ke ; if (do_k(k)) write(file,'(F10.6," ",$)') &
682  ((vm(i,j,k)-cs%v_prev(i,j,k))*inorm(k)); enddo
683  write(file,'(/,"CAv: ",$)')
684  do k=ks,ke ; if (do_k(k)) write(file,'(F10.6," ",$)') &
685  (dt*adp%CAv(i,j,k)*inorm(k)); enddo
686  write(file,'(/,"PFv: ",$)')
687  do k=ks,ke ; if (do_k(k)) write(file,'(F10.6," ",$)') &
688  (dt*adp%PFv(i,j,k)*inorm(k)); enddo
689  write(file,'(/,"diffv: ",$)')
690  do k=ks,ke ; if (do_k(k)) write(file,'(F10.6," ",$)') &
691  (dt*us%s_to_T*adp%diffv(i,j,k)*inorm(k)); enddo
692 
693  if (associated(adp%gradKEu)) then
694  write(file,'(/,"KEv: ",$)')
695  do k=ks,ke ; if (do_k(k)) write(file,'(F10.6," ",$)') &
696  (dt*adp%gradKEv(i,j,k)*inorm(k)); enddo
697  endif
698  if (associated(adp%rv_x_u)) then
699  write(file,'(/,"Corv: ",$)')
700  do k=ks,ke ; if (do_k(k)) write(file,'(F10.6," ",$)') &
701  dt*(adp%CAv(i,j,k)-adp%rv_x_u(i,j,k))*inorm(k); enddo
702  endif
703  if (associated(adp%dv_dt_visc)) then
704  write(file,'(/,"dvv: ",$)')
705  do k=ks,ke ; if (do_k(k)) write(file,'(F10.6," ",$)') &
706  (dt*adp%dv_dt_visc(i,j,k)*inorm(k)); enddo
707  endif
708  if (associated(adp%dv_other)) then
709  write(file,'(/,"dv_other: ",$)')
710  do k=ks,ke ; if (do_k(k)) write(file,'(F10.6," ",$)') &
711  (adp%dv_other(i,j,k)*inorm(k)); enddo
712  endif
713  if (associated(cs%v_accel_bt)) then
714  write(file,'(/,"dvbt: ",$)')
715  do k=ks,ke ; if (do_k(k)) write(file,'(F10.6," ",$)') &
716  (dt*cs%v_accel_bt(i,j,k)*inorm(k)) ; enddo
717  endif
718  endif
719 
720  write(file,'(2/)')
721 
722  call flush(file)
723  endif
724 
725 end subroutine write_v_accel
726 
727 !> This subroutine initializes the parameters regulating how truncations are logged.
728 subroutine pointaccel_init(MIS, Time, G, param_file, diag, dirs, CS)
730  target, intent(in) :: mis !< For "MOM Internal State" a set of pointers
731  !! to the fields and accelerations that make
732  !! up the ocean's physical state.
733  type(time_type), target, intent(in) :: time !< The current model time.
734  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure.
735  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time
736  !! parameters.
737  type(diag_ctrl), target, intent(inout) :: diag !< A structure that is used to regulate
738  !! diagnostic output.
739  type(directories), intent(in) :: dirs !< A structure containing several relevant
740  !! directory paths.
741  type(pointaccel_cs), pointer :: cs !< A pointer that is set to point to the
742  !! control structure for this module.
743 ! This include declares and sets the variable "version".
744 #include "version_variable.h"
745  character(len=40) :: mdl = "MOM_PointAccel" ! This module's name.
746 
747  if (associated(cs)) return
748  allocate(cs)
749 
750  cs%diag => diag ; cs%Time => time
751 
752  cs%T => mis%T ; cs%S => mis%S ; cs%pbce => mis%pbce
753  cs%u_accel_bt => mis%u_accel_bt ; cs%v_accel_bt => mis%v_accel_bt
754  cs%u_prev => mis%u_prev ; cs%v_prev => mis%v_prev
755  cs%u_av => mis%u_av; if (.not.associated(mis%u_av)) cs%u_av => mis%u(:,:,:)
756  cs%v_av => mis%v_av; if (.not.associated(mis%v_av)) cs%v_av => mis%v(:,:,:)
757 
758  ! Read all relevant parameters and write them to the model log.
759  call log_version(param_file, mdl, version, "")
760  call get_param(param_file, mdl, "U_TRUNC_FILE", cs%u_trunc_file, &
761  "The absolute path to the file where the accelerations "//&
762  "leading to zonal velocity truncations are written. \n"//&
763  "Leave this empty for efficiency if this diagnostic is "//&
764  "not needed.", default="", debuggingparam=.true.)
765  call get_param(param_file, mdl, "V_TRUNC_FILE", cs%v_trunc_file, &
766  "The absolute path to the file where the accelerations "//&
767  "leading to meridional velocity truncations are written. \n"//&
768  "Leave this empty for efficiency if this diagnostic is "//&
769  "not needed.", default="", debuggingparam=.true.)
770  call get_param(param_file, mdl, "MAX_TRUNC_FILE_SIZE_PER_PE", cs%max_writes, &
771  "The maximum number of colums of truncations that any PE "//&
772  "will write out during a run.", default=50, debuggingparam=.true.)
773 
774  if (len_trim(dirs%output_directory) > 0) then
775  if (len_trim(cs%u_trunc_file) > 0) &
776  cs%u_trunc_file = trim(dirs%output_directory)//trim(cs%u_trunc_file)
777  if (len_trim(cs%v_trunc_file) > 0) &
778  cs%v_trunc_file = trim(dirs%output_directory)//trim(cs%v_trunc_file)
779  call log_param(param_file, mdl, "output_dir/U_TRUNC_FILE", cs%u_trunc_file)
780  call log_param(param_file, mdl, "output_dir/V_TRUNC_FILE", cs%v_trunc_file)
781  endif
782  cs%u_file = -1 ; cs%v_file = -1 ; cs%cols_written = 0
783 
784 end subroutine pointaccel_init
785 
786 end module mom_pointaccel
mom_time_manager
Wraps the FMS time manager functions.
Definition: MOM_time_manager.F90:2
mom_verticalgrid
Provides a transparent vertical ocean grid type and supporting routines.
Definition: MOM_verticalGrid.F90:2
mom_file_parser::log_version
An overloaded interface to log version information about modules.
Definition: MOM_file_parser.F90:109
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_get_input::directories
Container for paths and parameter file names.
Definition: MOM_get_input.F90:20
mom_pointaccel::pointaccel_cs
The control structure for the MOM_PointAccel module.
Definition: MOM_PointAccel.F90:32
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_get_input
Reads the only Fortran name list needed to boot-strap the model.
Definition: MOM_get_input.F90:6
mom_unit_scaling::unit_scale_type
Describes various unit conversion factors.
Definition: MOM_unit_scaling.F90:14
mom_verticalgrid::verticalgrid_type
Describes the vertical ocean grid, including unit conversion factors.
Definition: MOM_verticalGrid.F90:24
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_variables::cont_diag_ptrs
Pointers to arrays with transports, which can later be used for derived diagnostics,...
Definition: MOM_variables.F90:185
mom_variables::accel_diag_ptrs
Pointers to arrays with accelerations, which can later be used for derived diagnostics,...
Definition: MOM_variables.F90:155
mom_grid
Provides the ocean grid type.
Definition: MOM_grid.F90:2
mom_unit_scaling
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
Definition: MOM_unit_scaling.F90:2
mom_file_parser::log_param
An overloaded interface to log the values of various types of parameters.
Definition: MOM_file_parser.F90:96
mom_variables::ocean_internal_state
Pointers to all of the prognostic variables allocated in MOM_variables.F90 and MOM....
Definition: MOM_variables.F90:126
mom_pointaccel
Debug accelerations at a given point.
Definition: MOM_PointAccel.F90:8
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