oil_tracer.F90
1! This file is part of MOM6, the Modular Ocean Model version 6.
2! See the LICENSE file for licensing information.
3! SPDX-License-Identifier: Apache-2.0
4
5!> A tracer package to mimic dissolved oil.
6module oil_tracer
7
9use mom_coupler_types, only : set_coupler_type_data, atmos_ocn_coupler_flux
10use mom_diag_mediator, only : diag_ctrl
11use mom_error_handler, only : mom_error, fatal, warning
12use mom_file_parser, only : get_param, log_param, log_version, param_file_type
13use mom_forcing_type, only : forcing
14use mom_grid, only : ocean_grid_type
15use mom_hor_index, only : hor_index_type
16use mom_io, only : file_exists, mom_read_data, slasher
17use mom_io, only : vardesc, var_desc, query_vardesc
18use mom_open_boundary, only : ocean_obc_type
19use mom_restart, only : query_initialized, set_initialized, mom_restart_cs
20use mom_spatial_means, only : global_mass_int_efp
21use mom_sponge, only : set_up_sponge_field, sponge_cs
22use mom_time_manager, only : time_type, time_type_to_real
23use mom_tracer_registry, only : register_tracer, tracer_registry_type
24use mom_tracer_diabatic, only : tracer_vertdiff, applytracerboundaryfluxesinout
25use mom_tracer_z_init, only : tracer_z_init
26use mom_unit_scaling, only : unit_scale_type
27use mom_variables, only : surface, thermo_var_ptrs
28use mom_verticalgrid, only : verticalgrid_type
29
30implicit none ; private
31
32#include <MOM_memory.h>
33
34public register_oil_tracer, initialize_oil_tracer
36public oil_stock, oil_tracer_end
37
38integer, parameter :: ntr_max = 20 !< the maximum number of tracers in this module.
39
40!> The control structure for the oil tracer package
41type, public :: oil_tracer_cs ; private
42 integer :: ntr !< The number of tracers that are actually used.
43 logical :: coupled_tracers = .false. !< These tracers are not offered to the coupler.
44 character(len=200) :: ic_file !< The file in which the age-tracer initial values
45 !! can be found, or an empty string for internal initialization.
46 logical :: z_ic_file !< If true, the IC_file is in Z-space. The default is false.
47 real :: oil_source_longitude !< Latitude of source location (geographic) [degrees_N]
48 real :: oil_source_latitude !< Longitude of source location (geographic) [degrees_E]
49 integer :: oil_source_i=-999 !< Local i of source location (computational index location)
50 integer :: oil_source_j=-999 !< Local j of source location (computational index location)
51 real :: oil_source_rate !< Rate of oil injection [kg T-1 ~> kg s-1]
52 real :: oil_start_year !< The time at which the oil source starts [years]
53 real :: oil_end_year !< The time at which the oil source ends [years]
54 type(time_type), pointer :: time => null() !< A pointer to the ocean model's clock.
55 type(tracer_registry_type), pointer :: tr_reg => null() !< A pointer to the MOM tracer registry
56 real, pointer :: tr(:,:,:,:) => null() !< The array of tracers used in this subroutine, [kg m-3]
57 real, dimension(NTR_MAX) :: ic_val = 0.0 !< The (uniform) initial condition value [kg m-3]
58 real, dimension(NTR_MAX) :: land_val = -1.0 !< The value of tr used where land is masked out [kg m-3]
59 real, dimension(NTR_MAX) :: oil_decay_rate !< Decay rate of oil [T-1 ~> s-1] calculated from oil_decay_days
60 integer, dimension(NTR_MAX) :: oil_source_k !< Layer of source
61 logical :: oil_may_reinit !< If true, oil tracers may be reset by the initialization code
62 !! if they are not found in the restart files.
63 integer, dimension(NTR_MAX) :: ind_tr !< Indices returned by atmos_ocn_coupler_flux if it is used and the
64 !! surface tracer concentrations are to be provided to the coupler.
65 type(vardesc) :: tr_desc(ntr_max) !< Descriptions and metadata for the tracers
66
67 type(diag_ctrl), pointer :: diag => null() !< A structure that is used to
68 !! regulate the timing of diagnostic output.
69 type(mom_restart_cs), pointer :: restart_csp => null() !< A pointer to the restart control structure
70end type oil_tracer_cs
71
72contains
73
74!> Register oil tracer fields and subroutines to be used with MOM.
75function register_oil_tracer(HI, GV, US, param_file, CS, tr_Reg, restart_CS)
76 type(hor_index_type), intent(in) :: hi !< A horizontal index type structure
77 type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure
78 type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
79 type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
80 type(oil_tracer_cs), pointer :: cs !< A pointer that is set to point to the control
81 !! structure for this module
82 type(tracer_registry_type), pointer :: tr_reg !< A pointer that is set to point to the control
83 !! structure for the tracer advection and
84 !! diffusion module
85 type(mom_restart_cs), target, intent(inout) :: restart_cs !< MOM restart control struct
86
87 ! Local variables
88 character(len=40) :: mdl = "oil_tracer" ! This module's name.
89 ! This include declares and sets the variable "version".
90# include "version_variable.h"
91 real, dimension(NTR_MAX) :: oil_decay_days !< Decay time scale of oil [days]
92 character(len=200) :: inputdir ! The directory where the input files are.
93 character(len=48) :: var_name ! The variable's name.
94 character(len=3) :: name_tag ! String for creating identifying oils
95 character(len=48) :: flux_units ! The units for tracer fluxes, here
96 ! kg(oil) s-1 or kg(oil) m-3 kg(water) s-1.
97 real, pointer :: tr_ptr(:,:,:) => null() ! The tracer concentration [kg m-3]
98 logical :: register_oil_tracer
99 integer :: isd, ied, jsd, jed, nz, m
100 isd = hi%isd ; ied = hi%ied ; jsd = hi%jsd ; jed = hi%jed ; nz = gv%ke
101
102 if (associated(cs)) then
103 call mom_error(fatal, "register_oil_tracer called with an "// &
104 "associated control structure.")
105 endif
106 allocate(cs)
107
108 ! Read all relevant parameters and write them to the model log.
109 call log_version(param_file, mdl, version, "")
110 call get_param(param_file, mdl, "OIL_IC_FILE", cs%IC_file, &
111 "The file in which the oil tracer initial values can be "//&
112 "found, or an empty string for internal initialization.", &
113 default=" ")
114 if ((len_trim(cs%IC_file) > 0) .and. (scan(cs%IC_file,'/') == 0)) then
115 ! Add the directory if CS%IC_file is not already a complete path.
116 call get_param(param_file, mdl, "INPUTDIR", inputdir, default=".")
117 cs%IC_file = trim(slasher(inputdir))//trim(cs%IC_file)
118 call log_param(param_file, mdl, "INPUTDIR/OIL_IC_FILE", cs%IC_file)
119 endif
120 call get_param(param_file, mdl, "OIL_IC_FILE_IS_Z", cs%Z_IC_file, &
121 "If true, OIL_IC_FILE is in depth space, not layer space", &
122 default=.false.)
123
124 call get_param(param_file, mdl, "OIL_MAY_REINIT", cs%oil_may_reinit, &
125 "If true, oil tracers may go through the initialization "//&
126 "code if they are not found in the restart files. "//&
127 "Otherwise it is a fatal error if the oil tracers are not "//&
128 "found in the restart files of a restarted run.", &
129 default=.false.)
130 call get_param(param_file, mdl, "OIL_SOURCE_LONGITUDE", cs%oil_source_longitude, &
131 "The geographic longitude of the oil source.", units="degrees_E", &
132 fail_if_missing=.true.)
133 call get_param(param_file, mdl, "OIL_SOURCE_LATITUDE", cs%oil_source_latitude, &
134 "The geographic latitude of the oil source.", units="degrees_N", &
135 fail_if_missing=.true.)
136 call get_param(param_file, mdl, "OIL_SOURCE_LAYER", cs%oil_source_k, &
137 "The layer into which the oil is introduced, or a "//&
138 "negative number for a vertically uniform source, "//&
139 "or 0 not to use this tracer.", units="Layer", default=0)
140 call get_param(param_file, mdl, "OIL_SOURCE_RATE", cs%oil_source_rate, &
141 "The rate of oil injection.", &
142 units="kg s-1", scale=us%T_to_s, default=1.0)
143 call get_param(param_file, mdl, "OIL_DECAY_DAYS", oil_decay_days, &
144 "The decay timescale in days (if positive), or no decay "//&
145 "if 0, or use the temperature dependent decay rate of "//&
146 "Adcroft et al. (GRL, 2010) if negative.", units="days", &
147 default=0.0)
148 call get_param(param_file, mdl, "OIL_DATED_START_YEAR", cs%oil_start_year, &
149 "The time at which the oil source starts", units="years", &
150 default=0.0)
151 call get_param(param_file, mdl, "OIL_DATED_END_YEAR", cs%oil_end_year, &
152 "The time at which the oil source ends", units="years", &
153 default=1.0e99)
154
155 cs%ntr = 0
156 cs%oil_decay_rate(:) = 0.
157 do m=1,ntr_max
158 if (cs%oil_source_k(m)/=0) then
159 write(name_tag(1:3),'("_",I2.2)') m
160 cs%ntr = cs%ntr + 1
161 cs%tr_desc(m) = var_desc("oil"//trim(name_tag), "kg m-3", "Oil Tracer", caller=mdl)
162 cs%IC_val(m) = 0.0
163 if (oil_decay_days(m) > 0.) then
164 cs%oil_decay_rate(m) = 1. / (86400.0*us%s_to_T * oil_decay_days(m))
165 elseif (oil_decay_days(m) < 0.) then
166 cs%oil_decay_rate(m) = -1.
167 endif
168 endif
169 enddo
170 call log_param(param_file, mdl, "OIL_DECAY_RATE", cs%oil_decay_rate(1:cs%ntr), &
171 units="s-1", unscale=us%s_to_T)
172
173 ! This needs to be changed if the units of tracer are changed above.
174 if (gv%Boussinesq) then ; flux_units = "kg s-1"
175 else ; flux_units = "kg m-3 kg s-1" ; endif
176
177 allocate(cs%tr(isd:ied,jsd:jed,nz,cs%ntr), source=0.0)
178
179 do m=1,cs%ntr
180 ! This is needed to force the compiler not to do a copy in the registration
181 ! calls. Curses on the designers and implementers of Fortran90.
182 tr_ptr => cs%tr(:,:,:,m)
183 call query_vardesc(cs%tr_desc(m), name=var_name, caller="register_oil_tracer")
184 ! Register the tracer for horizontal advection, diffusion, and restarts.
185 call register_tracer(tr_ptr, tr_reg, param_file, hi, gv, tr_desc=cs%tr_desc(m), &
186 registry_diags=.true., flux_units=flux_units, restart_cs=restart_cs, &
187 mandatory=.not.cs%oil_may_reinit)
188
189 ! Set coupled_tracers to be true (hard-coded above) to provide the surface
190 ! values to the coupler (if any). This is meta-code and its arguments will
191 ! currently (deliberately) give fatal errors if it is used.
192 if (cs%coupled_tracers) &
193 cs%ind_tr(m) = atmos_ocn_coupler_flux(trim(var_name)//'_flux', &
194 flux_type=' ', implementation=' ', caller="register_oil_tracer")
195 enddo
196
197 cs%tr_Reg => tr_reg
198 cs%restart_CSp => restart_cs
199 register_oil_tracer = .true.
200
201end function register_oil_tracer
202
203!> Initialize the oil tracers and set up tracer output
204subroutine initialize_oil_tracer(restart, day, G, GV, US, h, diag, OBC, CS, &
205 sponge_CSp)
206 logical, intent(in) :: restart !< .true. if the fields have already
207 !! been read from a restart file.
208 type(time_type), target, intent(in) :: day !< Time of the start of the run.
209 type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
210 type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure
211 type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
212 real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
213 intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
214 type(diag_ctrl), target, intent(in) :: diag !< A structure that is used to regulate
215 !! diagnostic output.
216 type(ocean_obc_type), pointer :: obc !< This open boundary condition type specifies
217 !! whether, where, and what open boundary
218 !! conditions are used.
219 type(oil_tracer_cs), pointer :: cs !< The control structure returned by a previous
220 !! call to register_oil_tracer.
221 type(sponge_cs), pointer :: sponge_csp !< Pointer to the control structure for the sponges.
222
223 ! Local variables
224 character(len=16) :: name ! A variable's name in a NetCDF file.
225 logical :: ok
226 integer :: i, j, k, is, ie, js, je, isd, ied, jsd, jed, nz, m
227 integer :: isdb, iedb, jsdb, jedb
228
229 if (.not.associated(cs)) return
230 if (cs%ntr < 1) return
231 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
232 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
233 isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
234
235 ! Establish location of source
236 do j=g%jsdB+1,g%jed ; do i=g%isdB+1,g%ied
237 ! This test for i,j index is specific to a lat/lon (non-rotated grid).
238 ! and needs to be generalized to work properly on the tri-polar grid.
239 if (cs%oil_source_longitude<g%geoLonBu(i,j) .and. &
240 cs%oil_source_longitude>=g%geoLonBu(i-1,j) .and. &
241 cs%oil_source_latitude<g%geoLatBu(i,j) .and. &
242 cs%oil_source_latitude>=g%geoLatBu(i,j-1) ) then
243 cs%oil_source_i=i
244 cs%oil_source_j=j
245 endif
246 enddo ; enddo
247
248 cs%Time => day
249 cs%diag => diag
250
251 do m=1,cs%ntr
252 call query_vardesc(cs%tr_desc(m), name=name, caller="initialize_oil_tracer")
253 if ((.not.restart) .or. (cs%oil_may_reinit .and. .not. &
254 query_initialized(cs%tr(:,:,:,m), name, cs%restart_CSp))) then
255
256 if (len_trim(cs%IC_file) > 0) then
257 ! Read the tracer concentrations from a netcdf file.
258 if (.not.file_exists(cs%IC_file, g%Domain)) &
259 call mom_error(fatal, "initialize_oil_tracer: Unable to open "//cs%IC_file)
260
261 if (cs%Z_IC_file) then
262 ok = tracer_z_init(cs%tr(:,:,:,m), h, cs%IC_file, name, &
263 g, gv, us, -1e34, 0.0) ! CS%land_val(m))
264 if (.not.ok) then
265 ok = tracer_z_init(cs%tr(:,:,:,m), h, cs%IC_file, &
266 trim(name), g, gv, us, -1e34, 0.0) ! CS%land_val(m))
267 if (.not.ok) call mom_error(fatal,"initialize_oil_tracer: "//&
268 "Unable to read "//trim(name)//" from "//&
269 trim(cs%IC_file)//".")
270 endif
271 else
272 call mom_read_data(cs%IC_file, trim(name), cs%tr(:,:,:,m), g%Domain)
273 endif
274 else
275 do k=1,nz ; do j=js,je ; do i=is,ie
276 if (g%mask2dT(i,j) < 0.5) then
277 cs%tr(i,j,k,m) = cs%land_val(m)
278 else
279 cs%tr(i,j,k,m) = cs%IC_val(m)
280 endif
281 enddo ; enddo ; enddo
282 endif
283 call set_initialized(cs%tr(:,:,:,m), name, cs%restart_CSp)
284 endif ! restart
285 enddo ! Tracer loop
286
287 if (associated(obc)) then
288 ! Put something here...
289 endif
290
291end subroutine initialize_oil_tracer
292
293!> Apply sources, sinks, diapycnal mixing and rising motions to the oil tracers
294subroutine oil_tracer_column_physics(h_old, h_new, ea, eb, fluxes, dt, G, GV, US, CS, tv, &
295 evap_CFL_limit, minimum_forcing_depth)
296 type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
297 type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure
298 real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
299 intent(in) :: h_old !< Layer thickness before entrainment [H ~> m or kg m-2].
300 real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
301 intent(in) :: h_new !< Layer thickness after entrainment [H ~> m or kg m-2].
302 real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
303 intent(in) :: ea !< an array to which the amount of fluid entrained
304 !! from the layer above during this call will be
305 !! added [H ~> m or kg m-2].
306 real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
307 intent(in) :: eb !< an array to which the amount of fluid entrained
308 !! from the layer below during this call will be
309 !! added [H ~> m or kg m-2].
310 type(forcing), intent(in) :: fluxes !< A structure containing pointers to thermodynamic
311 !! and tracer forcing fields. Unused fields have NULL ptrs.
312 real, intent(in) :: dt !< The amount of time covered by this call [T ~> s]
313 type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
314 type(oil_tracer_cs), pointer :: cs !< The control structure returned by a previous
315 !! call to register_oil_tracer.
316 type(thermo_var_ptrs), intent(in) :: tv !< A structure pointing to various thermodynamic variables
317 real, optional, intent(in) :: evap_cfl_limit !< Limit on the fraction of the water that can
318 !! be fluxed out of the top layer in a timestep [nondim]
319 real, optional, intent(in) :: minimum_forcing_depth !< The smallest depth over which
320 !! fluxes can be applied [H ~> m or kg m-2]
321! This subroutine applies diapycnal diffusion and any other column
322! tracer physics or chemistry to the tracers from this file.
323! This is a simple example of a set of advected passive tracers.
324
325! The arguments to this subroutine are redundant in that
326! h_new(k) = h_old(k) + ea(k) - eb(k-1) + eb(k) - ea(k+1)
327
328 ! Local variables
329 real, dimension(SZI_(G),SZJ_(G),SZK_(GV)) :: h_work ! Used so that h can be modified [H ~> m or kg m-2]
330 real :: isecs_per_year = 1.0 / (365.0*86400.0) ! Conversion factor from seconds to year [year s-1]
331 real :: vol_scale ! A conversion factor for volumes into m3 [m3 H-1 L-2 ~> 1 or m3 kg-1]
332 real :: year ! Time in fractional years [years]
333 real :: h_total ! A running sum of thicknesses [H ~> m or kg m-2]
334 real :: decay_timescale ! Chemical decay timescale for oil [T ~> s]
335 real :: ldecay ! Chemical decay rate of oil [T-1 ~> s-1]
336 integer :: i, j, k, is, ie, js, je, nz, m, k_max
337 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
338
339 if (.not.associated(cs)) return
340 if (cs%ntr < 1) return
341
342 if (present(evap_cfl_limit) .and. present(minimum_forcing_depth)) then
343 do m=1,cs%ntr
344 do k=1,nz ;do j=js,je ; do i=is,ie
345 h_work(i,j,k) = h_old(i,j,k)
346 enddo ; enddo ; enddo
347 call applytracerboundaryfluxesinout(g, gv, cs%tr(:,:,:,m), dt, fluxes, h_work, &
348 evap_cfl_limit, minimum_forcing_depth)
349 call tracer_vertdiff(h_work, ea, eb, dt, cs%tr(:,:,:,m), g, gv)
350 enddo
351 else
352 do m=1,cs%ntr
353 call tracer_vertdiff(h_old, ea, eb, dt, cs%tr(:,:,:,m), g, gv)
354 enddo
355 endif
356
357 year = time_type_to_real(cs%Time) * isecs_per_year
358
359 ! Decay tracer (limit decay rate to 1/dt - just in case)
360 do m=2,cs%ntr
361 do k=1,nz ; do j=js,je ; do i=is,ie
362 !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
363 !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
364 if (cs%oil_decay_rate(m)>0.) then
365 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
366 elseif (cs%oil_decay_rate(m)<0.) then
367 decay_timescale = (12.0 * (3.0**(-(tv%T(i,j,k)-20.0*us%degC_to_C)/10.0*us%degC_to_C))) * &
368 (86400.0*us%s_to_T) ! Timescale [T ~> s]
369 ldecay = 1. / decay_timescale ! Rate [T-1 ~> s-1]
370 cs%tr(i,j,k,m) = g%mask2dT(i,j)*max(1. - dt*ldecay,0.)*cs%tr(i,j,k,m)
371 endif
372 enddo ; enddo ; enddo
373 enddo
374
375 ! Add oil at the source location
376 if (year>=cs%oil_start_year .and. year<=cs%oil_end_year .and. &
377 cs%oil_source_i>-999 .and. cs%oil_source_j>-999) then
378 i = cs%oil_source_i ; j = cs%oil_source_j
379 k_max = nz ; h_total = 0.
380 vol_scale = gv%H_to_m * us%L_to_m**2
381 do k=nz, 2, -1
382 h_total = h_total + h_new(i,j,k)
383 if (h_total < 10.*gv%m_to_H) k_max=k-1 ! Find bottom most interface that is 10 m above bottom
384 enddo
385 do m=1,cs%ntr
386 k = cs%oil_source_k(m)
387 if (k>0) then
388 k = min(k,k_max) ! Only insert k or first layer with interface 10 m above bottom
389 cs%tr(i,j,k,m) = cs%tr(i,j,k,m) + cs%oil_source_rate*dt / &
390 (vol_scale * (h_new(i,j,k)+gv%H_subroundoff) * g%areaT(i,j) )
391 elseif (k<0) then
392 h_total = gv%H_subroundoff
393 do k=1, nz
394 h_total = h_total + h_new(i,j,k)
395 enddo
396 do k=1, nz
397 cs%tr(i,j,k,m) = cs%tr(i,j,k,m) + cs%oil_source_rate*dt / (vol_scale * h_total * g%areaT(i,j) )
398 enddo
399 endif
400 enddo
401 endif
402
403end subroutine oil_tracer_column_physics
404
405!> Calculate the mass-weighted integral of the oil tracer stocks, returning the number of stocks it
406!! has calculated. If the stock_index is present, only the stock corresponding to that coded index is returned.
407function oil_stock(h, stocks, G, GV, CS, names, units, stock_index)
408 type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
409 type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure
410 real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
411 type(efp_type), dimension(:), intent(out) :: stocks !< The mass-weighted integrated amount of each
412 !! tracer, in kg times concentration units [kg conc]
413 type(oil_tracer_cs), pointer :: cs !< The control structure returned by a previous
414 !! call to register_oil_tracer.
415 character(len=*), dimension(:), intent(out) :: names !< the names of the stocks calculated.
416 character(len=*), dimension(:), intent(out) :: units !< the units of the stocks calculated.
417 integer, optional, intent(in) :: stock_index !< the coded index of a specific stock
418 !! being sought.
419 integer :: oil_stock !< The number of stocks calculated here.
420
421 ! Local variables
422 integer :: m
423
424 oil_stock = 0
425 if (.not.associated(cs)) return
426 if (cs%ntr < 1) return
427
428 if (present(stock_index)) then ; if (stock_index > 0) then
429 ! Check whether this stock is available from this routine.
430
431 ! No stocks from this routine are being checked yet. Return 0.
432 return
433 endif ; endif
434
435 do m=1,cs%ntr
436 call query_vardesc(cs%tr_desc(m), name=names(m), units=units(m), caller="oil_stock")
437 units(m) = trim(units(m))//" kg"
438 stocks(m) = global_mass_int_efp(h, g, gv, cs%tr(:,:,:,m), on_pe_only=.true.)
439 enddo
440 oil_stock = cs%ntr
441
442end function oil_stock
443
444!> This subroutine extracts the surface fields from this tracer package that
445!! are to be shared with the atmosphere in coupled configurations.
446!! This particular tracer package does not report anything back to the coupler.
447subroutine oil_tracer_surface_state(sfc_state, h, G, GV, CS)
448 type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure.
449 type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure
450 type(surface), intent(inout) :: sfc_state !< A structure containing fields that
451 !! describe the surface state of the ocean.
452 real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
453 intent(in) :: h !< Layer thickness [H ~> m or kg m-2].
454 type(oil_tracer_cs), pointer :: cs !< The control structure returned by a previous
455 !! call to register_oil_tracer.
456
457 ! This particular tracer package does not report anything back to the coupler.
458 ! The code that is here is just a rough guide for packages that would.
459
460 integer :: m, is, ie, js, je, isd, ied, jsd, jed
461 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
462 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
463
464 if (.not.associated(cs)) return
465
466 if (cs%coupled_tracers) then
467 do m=1,cs%ntr
468 ! This call loads the surface values into the appropriate array in the
469 ! coupler-type structure.
470 call set_coupler_type_data(cs%tr(:,:,1,m), cs%ind_tr(m), sfc_state%tr_fields, &
471 idim=(/isd, is, ie, ied/), jdim=(/jsd, js, je, jed/), turns=g%HI%turns)
472 enddo
473 endif
474
475end subroutine oil_tracer_surface_state
476
477!> Deallocate memory associated with this tracer package
478subroutine oil_tracer_end(CS)
479 type(oil_tracer_cs), pointer :: cs !< The control structure returned by a previous
480 !! call to register_oil_tracer.
481
482 if (associated(cs)) then
483 if (associated(cs%tr)) deallocate(cs%tr)
484 deallocate(cs)
485 endif
486end subroutine oil_tracer_end
487
488!> \namespace oil_tracer
489!!
490!! By Alistair Adcroft and Robert Hallberg, 2010 *
491!!
492!! In the midst of the Deepwater Horizon oil spill, it became evident that
493!! models were needed to predict the long-term fate of dissolved oil in the
494!! open ocean. This tracer packages mimics the transport, dilution and decay
495!! of dissolved oil plumes in the ocean.
496!!
497!! This tracer package was central to the simulations used by Adcroft et al.,
498!! GRL 2010, to prove that the Deepwater Horizon spill was an important regional
499!! event, with implications for dissolved oxygen levels in certains regions,
500!! see above reference for details.
501
502end module oil_tracer