MOM6
oil_tracer.F90
Go to the documentation of this file.
1 module oil_tracer
2 !***********************************************************************
3 !* GNU General Public License *
4 !* This file is a part of MOM. *
5 !* *
6 !* MOM is free software; you can redistribute it and/or modify it and *
7 !* are expected to follow the terms of the GNU General Public License *
8 !* as published by the Free Software Foundation; either version 2 of *
9 !* the License, or (at your option) any later version. *
10 !* *
11 !* MOM is distributed in the hope that it will be useful, but WITHOUT *
12 !* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY *
13 !* or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public *
14 !* License for more details. *
15 !* *
16 !* For the full text of the GNU General Public License, *
17 !* write to: Free Software Foundation, Inc., *
18 !* 675 Mass Ave, Cambridge, MA 02139, USA. *
19 !* or see: http://www.gnu.org/licenses/gpl.html *
20 !***********************************************************************
21 
22 !********+*********+*********+*********+*********+*********+*********+**
23 !* *
24 !* By Robert Hallberg, 2002 *
25 !* *
26 !* This file contains an example of the code that is needed to set *
27 !* up and use a set (in this case two) of dynamically passive tracers *
28 !* for diagnostic purposes. The tracers here are an ideal age tracer *
29 !* that ages at a rate of 1/year once it is isolated from the surface,*
30 !* and a vintage tracer, whose surface concentration grows exponen- *
31 !* with time with a 30-year timescale (similar to CFCs). *
32 !* *
33 !* A single subroutine is called from within each file to register *
34 !* each of the tracers for reinitialization and advection and to *
35 !* register the subroutine that initializes the tracers and set up *
36 !* their output and the subroutine that does any tracer physics or *
37 !* chemistry along with diapycnal mixing (included here because some *
38 !* tracers may float or swim vertically or dye diapycnal processes). *
39 !* *
40 !* *
41 !* Macros written all in capital letters are defined in MOM_memory.h. *
42 !* *
43 !* A small fragment of the grid is shown below: *
44 !* *
45 !* j+1 x ^ x ^ x At x: q *
46 !* j+1 > o > o > At ^: v *
47 !* j x ^ x ^ x At >: u *
48 !* j > o > o > At o: h, tr *
49 !* j-1 x ^ x ^ x *
50 !* i-1 i i+1 At x & ^: *
51 !* i i+1 At > & o: *
52 !* *
53 !* The boundaries always run through q grid points (x). *
54 !* *
55 !********+*********+*********+*********+*********+*********+*********+**
56 
57 use mom_diag_mediator, only : post_data, register_diag_field, safe_alloc_ptr
58 use mom_diag_mediator, only : diag_ctrl
60 use mom_error_handler, only : mom_error, fatal, warning
62 use mom_forcing_type, only : forcing
63 use mom_grid, only : ocean_grid_type
64 use mom_hor_index, only : hor_index_type
65 use mom_io, only : file_exists, read_data, slasher, vardesc, var_desc, query_vardesc
69 use mom_time_manager, only : time_type, get_time
74 use mom_variables, only : surface
77 use coupler_util, only : set_coupler_values, ind_csurf
79 
80 implicit none ; private
81 
82 #include <MOM_memory.h>
83 
87 
88 ! NTR_MAX is the maximum number of tracers in this module.
89 integer, parameter :: ntr_max = 20
90 
91 type p3d
92  real, dimension(:,:,:), pointer :: p => null()
93 end type p3d
94 
95 type, public :: oil_tracer_cs ; private
96  integer :: ntr ! The number of tracers that are actually used.
97  logical :: coupled_tracers = .false. ! These tracers are not offered to the
98  ! coupler.
99  character(len=200) :: ic_file ! The file in which the age-tracer initial values
100  ! can be found, or an empty string for internal initialization.
101  logical :: z_ic_file ! If true, the IC_file is in Z-space. The default is false.
102  real :: oil_source_longitude, oil_source_latitude ! Lat,lon of source location (geographic)
103  integer :: oil_source_i=-999, oil_source_j=-999 ! Local i,j of source location (computational)
104  real :: oil_source_rate ! Rate of oil injection (kg/s)
105  real :: oil_start_year ! The year in which tracers start aging, or at which the
106  ! surface value equals young_val, in years.
107  real :: oil_end_year ! The year in which tracers start aging, or at which the
108  ! surface value equals young_val, in years.
109  type(time_type), pointer :: time ! A pointer to the ocean model's clock.
110  type(tracer_registry_type), pointer :: tr_reg => null()
111  real, pointer :: tr(:,:,:,:) => null() ! The array of tracers used in this
112  ! subroutine, in g m-3?
113  type(p3d), dimension(NTR_MAX) :: &
114  tr_adx, &! Tracer zonal advective fluxes in g m-3 m3 s-1.
115  tr_ady, &! Tracer meridional advective fluxes in g m-3 m3 s-1.
116  tr_dfx, &! Tracer zonal diffusive fluxes in g m-3 m3 s-1.
117  tr_dfy ! Tracer meridional diffusive fluxes in g m-3 m3 s-1.
118  real, dimension(NTR_MAX) :: &
119  ic_val = 0.0, & ! The (uniform) initial condition value.
120  young_val = 0.0, & ! The value assigned to tr at the surface.
121  land_val = -1.0, & ! The value of tr used where land is masked out.
122  sfc_growth_rate ! The exponential growth rate for the surface value,
123  ! in units of year-1.
124  real, dimension(NTR_MAX) :: oil_decay_days, & ! Decay time scale of oil (in days)
125  oil_decay_rate ! Decay rate of oil (in s^-1) calculated from oil_decay_days
126  integer, dimension(NTR_MAX) :: oil_source_k ! Layer of source
127  logical :: mask_tracers ! If true, oils are masked out in massless layers.
128  logical :: oil_may_reinit ! If true, oil may go through the
129  ! initialization code if they are not found in the
130  ! restart files.
131  integer, dimension(NTR_MAX) :: &
132  ind_tr, & ! Indices returned by aof_set_coupler_flux if it is used and the
133  ! surface tracer concentrations are to be provided to the coupler.
134  id_tracer = -1, id_tr_adx = -1, id_tr_ady = -1, &
135  id_tr_dfx = -1, id_tr_dfy = -1
136 
137  type(diag_ctrl), pointer :: diag ! A structure that is used to regulate the
138  ! timing of diagnostic output.
139  type(mom_restart_cs), pointer :: restart_csp => null()
140 
141  type(vardesc) :: tr_desc(ntr_max)
142 end type oil_tracer_cs
143 
144 contains
145 
146 function register_oil_tracer(HI, GV, param_file, CS, tr_Reg, restart_CS)
147  type(hor_index_type), intent(in) :: HI
148  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
149  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
150  type(oil_tracer_cs), pointer :: CS
151  type(tracer_registry_type), pointer :: tr_Reg
152  type(mom_restart_cs), pointer :: restart_CS
153 ! This subroutine is used to register tracer fields and subroutines
154 ! to be used with MOM.
155 ! Arguments: HI - A horizontal index type structure.
156 ! (in) GV - The ocean's vertical grid structure.
157 ! (in) param_file - A structure indicating the open file to parse for
158 ! model parameter values.
159 ! (in/out) CS - A pointer that is set to point to the control structure
160 ! for this module
161 ! (in/out) tr_Reg - A pointer that is set to point to the control structure
162 ! for the tracer advection and diffusion module.
163 ! (in) restart_CS - A pointer to the restart control structure.
164 
165 ! This include declares and sets the variable "version".
166 #include "version_variable.h"
167  character(len=40) :: mdl = "oil_tracer" ! This module's name.
168  character(len=200) :: inputdir ! The directory where the input files are.
169  character(len=48) :: var_name ! The variable's name.
170  character(len=3) :: name_tag ! String for creating identifying oils
171  real, pointer :: tr_ptr(:,:,:) => null()
172  logical :: register_oil_tracer
173  integer :: isd, ied, jsd, jed, nz, m, i, j
174  isd = hi%isd ; ied = hi%ied ; jsd = hi%jsd ; jed = hi%jed ; nz = gv%ke
175 
176  if (associated(cs)) then
177  call mom_error(warning, "register_oil_tracer called with an "// &
178  "associated control structure.")
179  return
180  endif
181  allocate(cs)
182 
183  ! Read all relevant parameters and write them to the model log.
184  call log_version(param_file, mdl, version, "")
185  call get_param(param_file, mdl, "OIL_IC_FILE", cs%IC_file, &
186  "The file in which the oil tracer initial values can be \n"//&
187  "found, or an empty string for internal initialization.", &
188  default=" ")
189  if ((len_trim(cs%IC_file) > 0) .and. (scan(cs%IC_file,'/') == 0)) then
190  ! Add the directory if CS%IC_file is not already a complete path.
191  call get_param(param_file, mdl, "INPUTDIR", inputdir, default=".")
192  cs%IC_file = trim(slasher(inputdir))//trim(cs%IC_file)
193  call log_param(param_file, mdl, "INPUTDIR/CFC_IC_FILE", cs%IC_file)
194  endif
195  call get_param(param_file, mdl, "OIL_IC_FILE_IS_Z", cs%Z_IC_file, &
196  "If true, OIL_IC_FILE is in depth space, not layer space", &
197  default=.false.)
198 
199  call get_param(param_file, mdl, "MASK_MASSLESS_TRACERS", cs%mask_tracers, &
200  "If true, the tracers are masked out in massless layer. \n"//&
201  "This can be a problem with time-averages.", default=.false.)
202  call get_param(param_file, mdl, "OIL_MAY_REINIT", cs%oil_may_reinit, &
203  "If true, oil tracers may go through the initialization \n"//&
204  "code if they are not found in the restart files. \n"//&
205  "Otherwise it is a fatal error if the oil tracers are not \n"//&
206  "found in the restart files of a restarted run.", &
207  default=.false.)
208  call get_param(param_file, mdl, "OIL_SOURCE_LONGITUDE", cs%oil_source_longitude, &
209  "The geographic longitude of the oil source.", units="degrees E", &
210  fail_if_missing=.true.)
211  call get_param(param_file, mdl, "OIL_SOURCE_LATITUDE", cs%oil_source_latitude, &
212  "The geographic latitude of the oil source.", units="degrees N", &
213  fail_if_missing=.true.)
214  call get_param(param_file, mdl, "OIL_SOURCE_LAYER", cs%oil_source_k, &
215  "The layer into which the oil is introduced, or a \n"//&
216  "negative number for a vertically uniform source, \n"//&
217  "or 0 not to use this tracer.", units="Layer", default=0)
218  call get_param(param_file, mdl, "OIL_SOURCE_RATE", cs%oil_source_rate, &
219  "The rate of oil injection.", units="kg s-1", default=1.0)
220  call get_param(param_file, mdl, "OIL_DECAY_DAYS", cs%oil_decay_days, &
221  "The decay timescale in days (if positive), or no decay \n"//&
222  "if 0, or use the temperature dependent decay rate of \n"//&
223  "Adcroft et al. (GRL, 2010) if negative.", units="days", &
224  default=0.0)
225  call get_param(param_file, mdl, "OIL_DATED_START_YEAR", cs%oil_start_year, &
226  "The time at which the oil source starts", units="years", &
227  default=0.0)
228  call get_param(param_file, mdl, "OIL_DATED_END_YEAR", cs%oil_end_year, &
229  "The time at which the oil source ends", units="years", &
230  default=1.0e99)
231 
232  cs%ntr = 0
233  cs%oil_decay_rate(:) = 0.
234  do m=1,ntr_max
235  if (cs%oil_source_k(m)/=0) then
236  write(name_tag(1:3),'("_",I2.2)') m
237  cs%ntr = cs%ntr + 1
238  cs%tr_desc(m) = var_desc("oil"//trim(name_tag), "kg/m3", "Oil Tracer", caller=mdl)
239  cs%IC_val(m) = 0.0
240  if (cs%oil_decay_days(m)>0.) then
241  cs%oil_decay_rate(m)=1./(86400.0*cs%oil_decay_days(m))
242  elseif (cs%oil_decay_days(m)<0.) then
243  cs%oil_decay_rate(m)=-1.
244  endif
245  endif
246  enddo
247  call log_param(param_file, mdl, "OIL_DECAY_RATE", cs%oil_decay_rate(1:cs%ntr))
248 
249  allocate(cs%tr(isd:ied,jsd:jed,nz,cs%ntr)) ; cs%tr(:,:,:,:) = 0.0
250 
251  do m=1,cs%ntr
252  ! This is needed to force the compiler not to do a copy in the registration
253  ! calls. Curses on the designers and implementers of Fortran90.
254  tr_ptr => cs%tr(:,:,:,m)
255  call query_vardesc(cs%tr_desc(m), name=var_name, caller="register_oil_tracer")
256  ! Register the tracer for the restart file.
257  call register_restart_field(tr_ptr, cs%tr_desc(m), &
258  .not.cs%oil_may_reinit, restart_cs)
259  ! Register the tracer for horizontal advection & diffusion.
260  call register_tracer(tr_ptr, cs%tr_desc(m), param_file, hi, gv, tr_reg, &
261  tr_desc_ptr=cs%tr_desc(m))
262 
263  ! Set coupled_tracers to be true (hard-coded above) to provide the surface
264  ! values to the coupler (if any). This is meta-code and its arguments will
265  ! currently (deliberately) give fatal errors if it is used.
266  if (cs%coupled_tracers) &
267  cs%ind_tr(m) = aof_set_coupler_flux(trim(var_name)//'_flux', &
268  flux_type=' ', implementation=' ', caller="register_oil_tracer")
269  enddo
270 
271  cs%tr_Reg => tr_reg
272  cs%restart_CSp => restart_cs
273  register_oil_tracer = .true.
274 
275 end function register_oil_tracer
276 
277 subroutine initialize_oil_tracer(restart, day, G, GV, h, diag, OBC, CS, &
278  sponge_CSp, diag_to_Z_CSp)
279  logical, intent(in) :: restart
280  type(time_type), target, intent(in) :: day
281  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
282  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
283  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer thicknesses, in H (usually m or kg m-2)
284  type(diag_ctrl), target, intent(in) :: diag
285  type(ocean_obc_type), pointer :: OBC
286  type(oil_tracer_cs), pointer :: CS
287  type(sponge_cs), pointer :: sponge_CSp
288  type(diag_to_z_cs), pointer :: diag_to_Z_CSp
289 ! This subroutine initializes the CS%ntr tracer fields in tr(:,:,:,:)
290 ! and it sets up the tracer output.
291 
292 ! Arguments: restart - .true. if the fields have already been read from
293 ! a restart file.
294 ! (in) day - Time of the start of the run.
295 ! (in) G - The ocean's grid structure.
296 ! (in) GV - The ocean's vertical grid structure.
297 ! (in) h - Layer thickness, in m or kg m-2.
298 ! (in) diag - A structure that is used to regulate diagnostic output.
299 ! (in) OBC - This open boundary condition type specifies whether, where,
300 ! and what open boundary conditions are used.
301 ! (in/out) CS - The control structure returned by a previous call to
302 ! register_oil_tracer.
303 ! (in/out) sponge_CSp - A pointer to the control structure for the sponges, if
304 ! they are in use. Otherwise this may be unassociated.
305 ! (in/out) diag_to_Z_Csp - A pointer to the control structure for diagnostics
306 ! in depth space.
307  character(len=16) :: name ! A variable's name in a NetCDF file.
308  character(len=72) :: longname ! The long name of that variable.
309  character(len=48) :: units ! The dimensions of the variable.
310  character(len=48) :: flux_units ! The units for age tracer fluxes, either
311  ! years m3 s-1 or years kg s-1.
312  logical :: OK
313  integer :: i, j, k, is, ie, js, je, isd, ied, jsd, jed, nz, m
314  integer :: IsdB, IedB, JsdB, JedB
315 
316  if (.not.associated(cs)) return
317  if (cs%ntr < 1) return
318  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
319  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
320  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
321 
322  ! Establish location of source
323  do j=g%jsdB+1,g%jed ; do i=g%isdB+1,g%ied
324  ! This test for i,j index is specific to a lat/lon (non-rotated grid).
325  ! and needs to be generalized to work properly on the tri-polar grid.
326  if (cs%oil_source_longitude<g%geoLonBu(i,j) .and. &
327  cs%oil_source_longitude>=g%geoLonBu(i-1,j) .and. &
328  cs%oil_source_latitude<g%geoLatBu(i,j) .and. &
329  cs%oil_source_latitude>=g%geoLatBu(i,j-1) ) then
330  cs%oil_source_i=i
331  cs%oil_source_j=j
332  endif
333  enddo; enddo
334 
335  cs%Time => day
336  cs%diag => diag
337 
338  do m=1,cs%ntr
339  call query_vardesc(cs%tr_desc(m), name=name, caller="initialize_oil_tracer")
340  if ((.not.restart) .or. (cs%oil_may_reinit .and. .not. &
341  query_initialized(cs%tr(:,:,:,m), name, cs%restart_CSp))) then
342 
343  if (len_trim(cs%IC_file) > 0) then
344  ! Read the tracer concentrations from a netcdf file.
345  if (.not.file_exists(cs%IC_file, g%Domain)) &
346  call mom_error(fatal, "initialize_oil_tracer: "// &
347  "Unable to open "//cs%IC_file)
348 
349  if (cs%Z_IC_file) then
350  ok = tracer_z_init(cs%tr(:,:,:,m), h, cs%IC_file, name, &
351  g, -1e34, 0.0) ! CS%land_val(m))
352  if (.not.ok) then
353  ok = tracer_z_init(cs%tr(:,:,:,m), h, cs%IC_file, &
354  trim(name), g, -1e34, 0.0) ! CS%land_val(m))
355  if (.not.ok) call mom_error(fatal,"initialize_oil_tracer: "//&
356  "Unable to read "//trim(name)//" from "//&
357  trim(cs%IC_file)//".")
358  endif
359  else
360  call read_data(cs%IC_file, trim(name), cs%tr(:,:,:,m), &
361  domain=g%Domain%mpp_domain)
362  endif
363  else
364  do k=1,nz ; do j=js,je ; do i=is,ie
365  if (g%mask2dT(i,j) < 0.5) then
366  cs%tr(i,j,k,m) = cs%land_val(m)
367  else
368  cs%tr(i,j,k,m) = cs%IC_val(m)
369  endif
370  enddo ; enddo ; enddo
371  endif
372 
373  endif ! restart
374  enddo ! Tracer loop
375 
376  if (associated(obc)) then
377  ! All tracers but the first have 0 concentration in their inflows. As this
378  ! is the default value, the following calls are unnecessary.
379  ! do m=1,CS%ntr
380  ! call add_tracer_OBC_values(trim(CS%tr_desc(m)%name), CS%tr_Reg, 0.0)
381  ! enddo
382  endif
383 
384  ! This needs to be changed if the units of tracer are changed above.
385  if (gv%Boussinesq) then ; flux_units = "years m3 s-1"
386  else ; flux_units = "years kg s-1" ; endif
387 
388  do m=1,cs%ntr
389  ! Register the tracer for the restart file.
390  call query_vardesc(cs%tr_desc(m), name, units=units, longname=longname, &
391  caller="initialize_oil_tracer")
392  cs%id_tracer(m) = register_diag_field("ocean_model", trim(name), cs%diag%axesTL, &
393  day, trim(longname) , trim(units))
394  cs%id_tr_adx(m) = register_diag_field("ocean_model", trim(name)//"_adx", &
395  cs%diag%axesCuL, day, trim(longname)//" advective zonal flux" , &
396  trim(flux_units))
397  cs%id_tr_ady(m) = register_diag_field("ocean_model", trim(name)//"_ady", &
398  cs%diag%axesCvL, day, trim(longname)//" advective meridional flux" , &
399  trim(flux_units))
400  cs%id_tr_dfx(m) = register_diag_field("ocean_model", trim(name)//"_dfx", &
401  cs%diag%axesCuL, day, trim(longname)//" diffusive zonal flux" , &
402  trim(flux_units))
403  cs%id_tr_dfy(m) = register_diag_field("ocean_model", trim(name)//"_dfy", &
404  cs%diag%axesCvL, day, trim(longname)//" diffusive zonal flux" , &
405  trim(flux_units))
406  if (cs%id_tr_adx(m) > 0) call safe_alloc_ptr(cs%tr_adx(m)%p,isdb,iedb,jsd,jed,nz)
407  if (cs%id_tr_ady(m) > 0) call safe_alloc_ptr(cs%tr_ady(m)%p,isd,ied,jsdb,jedb,nz)
408  if (cs%id_tr_dfx(m) > 0) call safe_alloc_ptr(cs%tr_dfx(m)%p,isdb,iedb,jsd,jed,nz)
409  if (cs%id_tr_dfy(m) > 0) call safe_alloc_ptr(cs%tr_dfy(m)%p,isd,ied,jsdb,jedb,nz)
410 
411 ! Register the tracer for horizontal advection & diffusion.
412  if ((cs%id_tr_adx(m) > 0) .or. (cs%id_tr_ady(m) > 0) .or. &
413  (cs%id_tr_dfx(m) > 0) .or. (cs%id_tr_dfy(m) > 0)) &
414  call add_tracer_diagnostics(name, cs%tr_Reg, cs%tr_adx(m)%p, &
415  cs%tr_ady(m)%p,cs%tr_dfx(m)%p,cs%tr_dfy(m)%p)
416 
417  call register_z_tracer(cs%tr(:,:,:,m), trim(name), longname, units, &
418  day, g, diag_to_z_csp)
419  enddo
420 
421 end subroutine initialize_oil_tracer
422 
423 subroutine oil_tracer_column_physics(h_old, h_new, ea, eb, fluxes, dt, G, GV, CS, tv, &
424  evap_CFL_limit, minimum_forcing_depth)
425  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
426  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
427  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h_old, h_new, ea, eb
428  type(forcing), intent(in) :: fluxes
429  real, intent(in) :: dt !< The amount of time covered by this call, in s
430  type(oil_tracer_cs), pointer :: CS
431  type(thermo_var_ptrs), intent(in) :: tv !< A structure pointing to various thermodynamic variables
432  real, optional,intent(in) :: evap_CFL_limit
433  real, optional,intent(in) :: minimum_forcing_depth
434 ! This subroutine applies diapycnal diffusion and any other column
435 ! tracer physics or chemistry to the tracers from this file.
436 ! This is a simple example of a set of advected passive tracers.
437 
438 ! Arguments: h_old - Layer thickness before entrainment, in m or kg m-2.
439 ! (in) h_new - Layer thickness after entrainment, in m or kg m-2.
440 ! (in) ea - an array to which the amount of fluid entrained
441 ! from the layer above during this call will be
442 ! added, in m or kg m-2.
443 ! (in) eb - an array to which the amount of fluid entrained
444 ! from the layer below during this call will be
445 ! added, in m or kg m-2.
446 ! (in) fluxes - A structure containing pointers to any possible
447 ! forcing fields. Unused fields have NULL ptrs.
448 ! (in) dt - The amount of time covered by this call, in s.
449 ! (in) G - The ocean's grid structure.
450 ! (in) GV - The ocean's vertical grid structure.
451 ! (in) CS - The control structure returned by a previous call to
452 ! register_oil_tracer.
453 !
454 ! The arguments to this subroutine are redundant in that
455 ! h_new[k] = h_old[k] + ea[k] - eb[k-1] + eb[k] - ea[k+1]
456 
457  real, dimension(SZI_(G),SZJ_(G),SZK_(G)) :: h_work ! Used so that h can be modified
458  real :: Isecs_per_year = 1.0 / (365.0*86400.0)
459  real :: year, h_total, ldecay
460  integer :: secs, days
461  integer :: i, j, k, is, ie, js, je, nz, m, k_max
462  real, allocatable :: local_tr(:,:,:)
463  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
464 
465  if (.not.associated(cs)) return
466  if (cs%ntr < 1) return
467 
468  if (present(evap_cfl_limit) .and. present(minimum_forcing_depth)) then
469  do m=1,cs%ntr
470  do k=1,nz ;do j=js,je ; do i=is,ie
471  h_work(i,j,k) = h_old(i,j,k)
472  enddo ; enddo ; enddo;
473  call applytracerboundaryfluxesinout(g, gv, cs%tr(:,:,:,m) , dt, fluxes, h_work, &
474  evap_cfl_limit, minimum_forcing_depth)
475  call tracer_vertdiff(h_work, ea, eb, dt, cs%tr(:,:,:,m), g, gv)
476  enddo
477  else
478  do m=1,cs%ntr
479  call tracer_vertdiff(h_old, ea, eb, dt, cs%tr(:,:,:,m), g, gv)
480  enddo
481  endif
482 
483  ! Set the surface value of tracer 1 to increase exponentially
484  ! with a 30 year time scale.
485  call get_time(cs%Time, secs, days)
486  year = (86400.0*days + real(secs)) * isecs_per_year
487 
488  ! Decay tracer (limit decay rate to 1/dt - just in case)
489  do m=2,cs%ntr
490  do k=1,nz ; do j=js,je ; do i=is,ie
491  !CS%tr(i,j,k,m) = CS%tr(i,j,k,m) - dt*CS%oil_decay_rate(m)*CS%tr(i,j,k,m) ! Simple
492  !CS%tr(i,j,k,m) = CS%tr(i,j,k,m) - min(dt*CS%oil_decay_rate(m),1.)*CS%tr(i,j,k,m) ! Safer
493  if (cs%oil_decay_rate(m)>0.) then
494  cs%tr(i,j,k,m) = g%mask2dT(i,j)*max(1.-dt*cs%oil_decay_rate(m),0.)*cs%tr(i,j,k,m) ! Safest
495  elseif (cs%oil_decay_rate(m)<0.) then
496  ldecay = 12.*(3.0**(-(tv%T(i,j,k)-20.)/10.)) ! Timescale in days
497  ldecay = 1./(86400.*ldecay) ! Rate in s^-1
498  cs%tr(i,j,k,m) = g%mask2dT(i,j)*max(1.-dt*ldecay,0.)*cs%tr(i,j,k,m)
499  endif
500  enddo ; enddo ; enddo
501  enddo
502 
503  ! Add oil at the source location
504  if (year>=cs%oil_start_year .and. year<=cs%oil_end_year .and. &
505  cs%oil_source_i>-999 .and. cs%oil_source_j>-999) then
506  i=cs%oil_source_i ; j=cs%oil_source_j
507  k_max=nz ; h_total=0.
508  do k=nz, 2, -1
509  h_total = h_total + h_new(i,j,k)
510  if (h_total<10.) k_max=k-1 ! Find bottom most interface that is 10 m above bottom
511  enddo
512  do m=1,cs%ntr
513  k=cs%oil_source_k(m)
514  if (k>0) then
515  k=min(k,k_max) ! Only insert k or first layer with interface 10 m above bottom
516  cs%tr(i,j,k,m) = cs%tr(i,j,k,m) + cs%oil_source_rate*dt / &
517  ((h_new(i,j,k)+gv%H_subroundoff) * g%areaT(i,j) )
518  elseif (k<0) then
519  h_total=gv%H_subroundoff
520  do k=1, nz
521  h_total = h_total + h_new(i,j,k)
522  enddo
523  do k=1, nz
524  cs%tr(i,j,k,m) = cs%tr(i,j,k,m) + cs%oil_source_rate*dt/(h_total &
525  * g%areaT(i,j) )
526  enddo
527  endif
528  enddo
529  endif
530 
531  allocate(local_tr(g%isd:g%ied,g%jsd:g%jed,nz))
532  do m=1,cs%ntr
533  if (cs%id_tracer(m)>0) then
534  if (cs%mask_tracers) then
535  do k=1,nz ; do j=js,je ; do i=is,ie
536  if (h_new(i,j,k) < 1.1*gv%Angstrom) then
537  local_tr(i,j,k) = cs%land_val(m)
538  else
539  local_tr(i,j,k) = cs%tr(i,j,k,m)
540  endif
541  enddo ; enddo ; enddo
542  else
543  do k=1,nz ; do j=js,je ; do i=is,ie
544  local_tr(i,j,k) = cs%tr(i,j,k,m)
545  enddo ; enddo ; enddo
546  endif ! CS%mask_tracers
547  call post_data(cs%id_tracer(m),local_tr,cs%diag)
548  endif ! CS%id_tracer(m)>0
549  if (cs%id_tr_adx(m)>0) &
550  call post_data(cs%id_tr_adx(m),cs%tr_adx(m)%p(:,:,:),cs%diag)
551  if (cs%id_tr_ady(m)>0) &
552  call post_data(cs%id_tr_ady(m),cs%tr_ady(m)%p(:,:,:),cs%diag)
553  if (cs%id_tr_dfx(m)>0) &
554  call post_data(cs%id_tr_dfx(m),cs%tr_dfx(m)%p(:,:,:),cs%diag)
555  if (cs%id_tr_dfy(m)>0) &
556  call post_data(cs%id_tr_dfy(m),cs%tr_dfy(m)%p(:,:,:),cs%diag)
557  enddo
558  deallocate(local_tr)
559 
560 end subroutine oil_tracer_column_physics
561 
562 function oil_stock(h, stocks, G, GV, CS, names, units, stock_index)
563  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
564  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
565  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer thicknesses, in H (usually m or kg m-2)
566  real, dimension(:), intent(out) :: stocks
567  type(oil_tracer_cs), pointer :: CS
568  character(len=*), dimension(:), intent(out) :: names
569  character(len=*), dimension(:), intent(out) :: units
570  integer, optional, intent(in) :: stock_index
571  integer :: oil_stock
572 ! This function calculates the mass-weighted integral of all tracer stocks,
573 ! returning the number of stocks it has calculated. If the stock_index
574 ! is present, only the stock corresponding to that coded index is returned.
575 
576 ! Arguments: h - Layer thickness, in m or kg m-2.
577 ! (out) stocks - the mass-weighted integrated amount of each tracer,
578 ! in kg times concentration units.
579 ! (in) G - The ocean's grid structure.
580 ! (in) GV - The ocean's vertical grid structure.
581 ! (in) CS - The control structure returned by a previous call to
582 ! register_oil_tracer.
583 ! (out) names - the names of the stocks calculated.
584 ! (out) units - the units of the stocks calculated.
585 ! (in,opt) stock_index - the coded index of a specific stock being sought.
586 ! Return value: the number of stocks calculated here.
587 
588  integer :: i, j, k, is, ie, js, je, nz, m
589  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
590 
591  oil_stock = 0
592  if (.not.associated(cs)) return
593  if (cs%ntr < 1) return
594 
595  if (present(stock_index)) then ; if (stock_index > 0) then
596  ! Check whether this stock is available from this routine.
597 
598  ! No stocks from this routine are being checked yet. Return 0.
599  return
600  endif ; endif
601 
602  do m=1,cs%ntr
603  call query_vardesc(cs%tr_desc(m), name=names(m), units=units(m), caller="oil_stock")
604  units(m) = trim(units(m))//" kg"
605  stocks(m) = 0.0
606  do k=1,nz ; do j=js,je ; do i=is,ie
607  stocks(m) = stocks(m) + cs%tr(i,j,k,m) * &
608  (g%mask2dT(i,j) * g%areaT(i,j) * h(i,j,k))
609  enddo ; enddo ; enddo
610  stocks(m) = gv%H_to_kg_m2 * stocks(m)
611  enddo
612  oil_stock = cs%ntr
613 
614 end function oil_stock
615 
616 subroutine oil_tracer_surface_state(state, h, G, CS)
617  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
618  type(surface), intent(inout) :: state
619  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer thicknesses, in H (usually m or kg m-2)
620  type(oil_tracer_cs), pointer :: CS
621 ! This particular tracer package does not report anything back to the coupler.
622 ! The code that is here is just a rough guide for packages that would.
623 ! Arguments: state - A structure containing fields that describe the
624 ! surface state of the ocean.
625 ! (in) h - Layer thickness, in m or kg m-2.
626 ! (in) G - The ocean's grid structure.
627 ! (in) CS - The control structure returned by a previous call to
628 ! register_oil_tracer.
629  integer :: m, is, ie, js, je
630  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
631 
632  if (.not.associated(cs)) return
633 
634  if (cs%coupled_tracers) then
635  do m=1,cs%ntr
636  ! This call loads the surface vlues into the appropriate array in the
637  ! coupler-type structure.
638  call set_coupler_values(cs%tr(:,:,1,m), state%tr_fields, cs%ind_tr(m), &
639  ind_csurf, is, ie, js, je)
640  enddo
641  endif
642 
643 end subroutine oil_tracer_surface_state
644 
645 subroutine oil_tracer_end(CS)
646  type(oil_tracer_cs), pointer :: CS
647  integer :: m
648 
649  if (associated(cs)) then
650  if (associated(cs%tr)) deallocate(cs%tr)
651  do m=1,cs%ntr
652  if (associated(cs%tr_adx(m)%p)) deallocate(cs%tr_adx(m)%p)
653  if (associated(cs%tr_ady(m)%p)) deallocate(cs%tr_ady(m)%p)
654  if (associated(cs%tr_dfx(m)%p)) deallocate(cs%tr_dfx(m)%p)
655  if (associated(cs%tr_dfy(m)%p)) deallocate(cs%tr_dfy(m)%p)
656  enddo
657 
658  deallocate(cs)
659  endif
660 end subroutine oil_tracer_end
661 
662 end module oil_tracer
The following structure contains pointers to various fields which may be used describe the surface st...
type(vardesc) function, public var_desc(name, units, longname, hor_grid, z_grid, t_grid, cmor_field_name, cmor_units, conversion, caller)
Returns a vardesc type whose elements have been filled with the provided fields. The argument name is...
Definition: MOM_io.F90:585
subroutine, public oil_tracer_column_physics(h_old, h_new, ea, eb, fluxes, dt, G, GV, CS, tv, evap_CFL_limit, minimum_forcing_depth)
Definition: oil_tracer.F90:425
This module implements boundary forcing for MOM6.
integer function, public aof_set_coupler_flux(name, flux_type, implementation, atm_tr_index, param, flag, ice_restart_file, ocean_restart_file, units, caller)
logical function, public tracer_z_init(tr, h, filename, tr_name, G, missing_val, land_val)
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:19
The following data type a list of diagnostic fields an their variants, as well as variables that cont...
Provides the ocean grid type.
Definition: MOM_grid.F90:2
This module contains I/O framework code.
Definition: MOM_io.F90:2
Defines the horizontal index type (hor_index_type) used for providing index ranges.
subroutine, public initialize_oil_tracer(restart, day, G, GV, h, diag, OBC, CS, sponge_CSp, diag_to_Z_CSp)
Definition: oil_tracer.F90:279
Container for horizontal index ranges for data, computational and global domains. ...
subroutine, public register_tracer(tr1, tr_desc, param_file, HI, GV, Reg, tr_desc_ptr, ad_x, ad_y, df_x, df_y, OBC_inflow, OBC_in_u, OBC_in_v, ad_2d_x, ad_2d_y, df_2d_x, df_2d_y, advection_xy)
This subroutine registers a tracer to be advected and laterally diffused.
subroutine, public oil_tracer_surface_state(state, h, G, CS)
Definition: oil_tracer.F90:617
This module contains the tracer_registry_type and the subroutines that handle registration of tracers...
subroutine, public oil_tracer_end(CS)
Definition: oil_tracer.F90:646
subroutine, public applytracerboundaryfluxesinout(G, GV, Tr, dt, fluxes, h, evap_CFL_limit, minimum_forcing_depth, in_flux_optional, out_flux_optional, update_h_opt)
This routine is modeled after applyBoundaryFluxesInOut in MOM_diabatic_aux.F90 NOTE: Please note that...
subroutine, public set_up_sponge_field(sp_val, f_ptr, G, nlay, CS, sp_val_i_mean)
Definition: MOM_sponge.F90:271
Type to carry basic tracer information.
This module contains routines that implement physical fluxes of tracers (e.g. due to surface fluxes o...
subroutine, public register_z_tracer(tr_ptr, name, long_name, units, Time, G, CS, standard_name, cmor_field_name, cmor_long_name, cmor_units, cmor_standard_name)
This subroutine registers a tracer to be output in depth space.
integer function, public oil_stock(h, stocks, G, GV, CS, names, units, stock_index)
Definition: oil_tracer.F90:563
subroutine, public add_tracer_obc_values(name, Reg, OBC_inflow, OBC_in_u, OBC_in_v)
This subroutine adds open boundary condition concentrations for a tracer that has previously been reg...
subroutine, public add_tracer_diagnostics(name, Reg, ad_x, ad_y, df_x, df_y, ad_2d_x, ad_2d_y, df_2d_x, df_2d_y, advection_xy)
This subroutine adds diagnostic arrays for a tracer that has previously been registered by a call to ...
Structure that contains pointers to the boundary forcing used to drive the liquid ocean simulated by ...
subroutine, public set_coupler_values(array_in, BC_struc, BC_index, BC_element, is, ie, js, je, conversion)
Type for describing a variable, typically a tracer.
Definition: MOM_io.F90:51
subroutine, public tracer_vertdiff(h_old, ea, eb, dt, tr, G, GV, sfc_flux, btm_flux, btm_reservoir, sink_rate, convert_flux_in)
This subroutine solves a tridiagonal equation for the final tracer concentrations after the dual-entr...
integer, parameter ntr_max
Definition: oil_tracer.F90:89
logical function, public register_oil_tracer(HI, GV, param_file, CS, tr_Reg, restart_CS)
Definition: oil_tracer.F90:147
The thermo_var_ptrs structure contains pointers to an assortment of thermodynamic fields that may be ...
subroutine, public query_vardesc(vd, name, units, longname, hor_grid, z_grid, t_grid, cmor_field_name, cmor_units, conversion, caller)
This routine queries vardesc.
Definition: MOM_io.F90:664
Controls where open boundary conditions are applied.
integer function, public register_diag_field(module_name, field_name, axes, init_time, long_name, units, missing_value, range, mask_variant, standard_name, verbose, do_not_log, err_msg, interp_method, tile_count, cmor_field_name, cmor_long_name, cmor_units, cmor_standard_name, cell_methods, x_cell_method, y_cell_method, v_cell_method, conversion, v_extensive)
Returns the "diag_mediator" handle for a group (native, CMOR, z-coord, ...) of diagnostics derived fr...
subroutine, public mom_error(level, message, all_print)