25 use coupler_types_mod,
only : coupler_type_set_data, ind_csurf
28 implicit none ;
private
30 #include <MOM_memory.h>
32 public register_boundary_impulse_tracer, initialize_boundary_impulse_tracer
33 public boundary_impulse_tracer_column_physics, boundary_impulse_tracer_surface_state
34 public boundary_impulse_stock, boundary_impulse_tracer_end
37 integer,
parameter :: ntr_max = 1
41 integer :: ntr=ntr_max
42 logical :: coupled_tracers = .false.
43 type(time_type),
pointer :: time => null()
45 real,
pointer :: tr(:,:,:,:) => null()
46 logical :: tracers_may_reinit
47 integer,
dimension(NTR_MAX) :: ind_tr
51 real,
dimension(NTR_MAX) :: land_val = -1.0
53 real :: remaining_source_time
66 function register_boundary_impulse_tracer(HI, GV, param_file, CS, tr_Reg, restart_CS)
78 character(len=40) :: mdl =
"boundary_impulse_tracer"
79 character(len=200) :: inputdir
80 character(len=48) :: var_name
81 character(len=3) :: name_tag
82 character(len=48) :: flux_units
85 #include "version_variable.h"
86 real,
pointer :: tr_ptr(:,:,:) => null()
87 real,
pointer :: rem_time_ptr => null()
88 logical :: register_boundary_impulse_tracer
89 integer :: isd, ied, jsd, jed, nz, m, i, j
90 isd = hi%isd ; ied = hi%ied ; jsd = hi%jsd ; jed = hi%jed ; nz = gv%ke
92 if (
associated(cs))
then
93 call mom_error(warning,
"register_boundary_impulse_tracer called with an "// &
94 "associated control structure.")
101 call get_param(param_file, mdl,
"IMPULSE_SOURCE_TIME", cs%remaining_source_time, &
102 "Length of time for the boundary tracer to be injected "//&
103 "into the mixed layer. After this time has elapsed, the "//&
104 "surface becomes a sink for the boundary impulse tracer.", &
106 call get_param(param_file, mdl,
"TRACERS_MAY_REINIT", cs%tracers_may_reinit, &
107 "If true, tracers may go through the initialization code "//&
108 "if they are not found in the restart files. Otherwise "//&
109 "it is a fatal error if the tracers are not found in the "//&
110 "restart files of a restarted run.", default=.false.)
112 allocate(cs%tr(isd:ied,jsd:jed,nz,cs%ntr)) ; cs%tr(:,:,:,:) = 0.0
114 cs%nkml = max(gv%nkml,1)
119 cs%tr_desc(m) = var_desc(trim(
"boundary_impulse"),
"kg kg-1", &
120 "Boundary impulse tracer", caller=mdl)
121 if (gv%Boussinesq)
then ; flux_units =
"kg kg-1 m3 s-1"
122 else ; flux_units =
"kg s-1" ;
endif
124 tr_ptr => cs%tr(:,:,:,m)
125 call query_vardesc(cs%tr_desc(m), name=var_name, caller=
"register_boundary_impulse_tracer")
127 call register_tracer(tr_ptr, tr_reg, param_file, hi, gv, tr_desc=cs%tr_desc(m), &
128 registry_diags=.true., flux_units=flux_units, &
129 restart_cs=restart_cs, mandatory=.not.cs%tracers_may_reinit)
134 if (cs%coupled_tracers) &
135 cs%ind_tr(m) = aof_set_coupler_flux(trim(var_name)//
'_flux', &
136 flux_type=
' ', implementation=
' ', caller=
"register_boundary_impulse_tracer")
139 rem_time_ptr => cs%remaining_source_time
141 .not.cs%tracers_may_reinit, restart_cs, &
142 "Remaining time to apply BIR source",
"s")
145 cs%restart_CSp => restart_cs
146 register_boundary_impulse_tracer = .true.
148 end function register_boundary_impulse_tracer
151 subroutine initialize_boundary_impulse_tracer(restart, day, G, GV, h, diag, OBC, CS, &
153 logical,
intent(in) :: restart
155 type(time_type),
target,
intent(in) :: day
158 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
160 type(
diag_ctrl),
target,
intent(in) :: diag
171 character(len=16) :: name
172 character(len=72) :: longname
173 character(len=48) :: units
174 character(len=48) :: flux_units
177 integer :: i, j, k, is, ie, js, je, isd, ied, jsd, jed, nz, m
178 integer :: isdb, iedb, jsdb, jedb
180 if (.not.
associated(cs))
return
181 if (cs%ntr < 1)
return
182 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
183 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
184 isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
188 name =
"boundary_impulse"
191 call query_vardesc(cs%tr_desc(m), name=name, caller=
"initialize_boundary_impulse_tracer")
192 if ((.not.restart) .or. (.not. &
194 do k=1,cs%nkml ;
do j=jsd,jed ;
do i=isd,ied
196 enddo ;
enddo ;
enddo
200 if (
associated(obc))
then
204 end subroutine initialize_boundary_impulse_tracer
207 subroutine boundary_impulse_tracer_column_physics(h_old, h_new, ea, eb, fluxes, dt, G, GV, US, CS, &
208 tv, debug, evap_CFL_limit, minimum_forcing_depth)
211 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
213 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
215 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
219 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
223 type(
forcing),
intent(in) :: fluxes
225 real,
intent(in) :: dt
231 logical,
intent(in) :: debug
232 real,
optional,
intent(in) :: evap_cfl_limit
234 real,
optional,
intent(in) :: minimum_forcing_depth
245 real :: isecs_per_year = 1.0 / (365.0*86400.0)
246 real :: year, h_total, scale, htot, ih_limit
247 integer :: secs, days
248 integer :: i, j, k, is, ie, js, je, nz, m, k_max
249 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)) :: h_work
251 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
253 if (.not.
associated(cs))
return
254 if (cs%ntr < 1)
return
257 if (
present(evap_cfl_limit) .and.
present(minimum_forcing_depth))
then
258 do k=1,nz ;
do j=js,je ;
do i=is,ie
259 h_work(i,j,k) = h_old(i,j,k)
260 enddo ;
enddo ;
enddo
261 call applytracerboundaryfluxesinout(g, gv, cs%tr(:,:,:,1), dt, fluxes, h_work, &
262 evap_cfl_limit, minimum_forcing_depth)
263 call tracer_vertdiff(h_work, ea, eb, dt, cs%tr(:,:,:,1), g, gv)
265 call tracer_vertdiff(h_old, ea, eb, dt, cs%tr(:,:,:,1), g, gv)
270 if (cs%remaining_source_time>0.0)
then
271 do k=1,cs%nkml ;
do j=js,je ;
do i=is,ie
273 enddo ;
enddo ;
enddo
274 cs%remaining_source_time = cs%remaining_source_time-us%T_to_s*dt
276 do k=1,cs%nkml ;
do j=js,je ;
do i=is,ie
278 enddo ;
enddo ;
enddo
283 end subroutine boundary_impulse_tracer_column_physics
286 function boundary_impulse_stock(h, stocks, G, GV, CS, names, units, stock_index)
289 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in ) :: h
290 real,
dimension(:),
intent( out) :: stocks
294 character(len=*),
dimension(:),
intent( out) :: names
295 character(len=*),
dimension(:),
intent( out) :: units
296 integer,
optional,
intent(in ) :: stock_index
298 integer :: boundary_impulse_stock
305 integer :: i, j, k, is, ie, js, je, nz, m
306 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
308 boundary_impulse_stock = 0
309 if (.not.
associated(cs))
return
310 if (cs%ntr < 1)
return
312 if (
present(stock_index))
then ;
if (stock_index > 0)
then
320 call query_vardesc(cs%tr_desc(m), name=names(m), units=units(m), caller=
"boundary_impulse_stock")
321 units(m) = trim(units(m))//
" kg"
323 do k=1,nz ;
do j=js,je ;
do i=is,ie
324 stocks(m) = stocks(m) + cs%tr(i,j,k,m) * &
325 (g%mask2dT(i,j) * g%US%L_to_m**2*g%areaT(i,j) * h(i,j,k))
326 enddo ;
enddo ;
enddo
327 stocks(m) = gv%H_to_kg_m2 * stocks(m)
330 boundary_impulse_stock = cs%ntr
332 end function boundary_impulse_stock
337 subroutine boundary_impulse_tracer_surface_state(sfc_state, h, G, CS)
339 type(
surface),
intent(inout) :: sfc_state
341 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
349 integer :: m, is, ie, js, je, isd, ied, jsd, jed
350 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
351 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
353 if (.not.
associated(cs))
return
355 if (cs%coupled_tracers)
then
359 call coupler_type_set_data(cs%tr(:,:,1,m), cs%ind_tr(m), ind_csurf, &
360 sfc_state%tr_fields, idim=(/isd, is, ie, ied/), &
361 jdim=(/jsd, js, je, jed/) )
365 end subroutine boundary_impulse_tracer_surface_state
368 subroutine boundary_impulse_tracer_end(CS)
373 if (
associated(cs))
then
374 if (
associated(cs%tr))
deallocate(cs%tr)
377 end subroutine boundary_impulse_tracer_end