31 use coupler_types_mod,
only : coupler_type_set_data, ind_csurf
34 implicit none ;
private
36 #include <MOM_memory.h>
39 public register_isomip_tracer, initialize_isomip_tracer
40 public isomip_tracer_column_physics, isomip_tracer_surface_state, isomip_tracer_end
42 integer,
parameter :: ntr = 1
46 logical :: coupled_tracers = .false.
47 character(len = 200) :: tracer_ic_file
48 type(time_type),
pointer :: time
50 real,
pointer :: tr(:,:,:,:) => null()
51 real :: land_val(ntr) = -1.0
54 integer,
dimension(NTR) :: ind_tr
68 function register_isomip_tracer(HI, GV, param_file, CS, tr_Reg, restart_CS)
78 character(len=80) :: name, longname
80 #include "version_variable.h"
81 character(len=40) :: mdl =
"ISOMIP_tracer"
82 character(len=200) :: inputdir
83 character(len=48) :: flux_units
85 real,
pointer :: tr_ptr(:,:,:) => null()
86 logical :: register_isomip_tracer
87 integer :: isd, ied, jsd, jed, nz, m
88 isd = hi%isd ; ied = hi%ied ; jsd = hi%jsd ; jed = hi%jed ; nz = gv%ke
90 if (
associated(cs))
then
91 call mom_error(warning,
"ISOMIP_register_tracer called with an "// &
92 "associated control structure.")
99 call get_param(param_file, mdl,
"ISOMIP_TRACER_IC_FILE", cs%tracer_IC_file, &
100 "The name of a file from which to read the initial "//&
101 "conditions for the ISOMIP tracers, or blank to initialize "//&
102 "them internally.", default=
" ")
103 if (len_trim(cs%tracer_IC_file) >= 1)
then
104 call get_param(param_file, mdl,
"INPUTDIR", inputdir, default=
".")
105 inputdir = slasher(inputdir)
106 cs%tracer_IC_file = trim(inputdir)//trim(cs%tracer_IC_file)
107 call log_param(param_file, mdl,
"INPUTDIR/ISOMIP_TRACER_IC_FILE", &
110 call get_param(param_file, mdl,
"SPONGE", cs%use_sponge, &
111 "If true, sponges may be applied anywhere in the domain. "//&
112 "The exact location and properties of those sponges are "//&
113 "specified from MOM_initialization.F90.", default=.false.)
115 allocate(cs%tr(isd:ied,jsd:jed,nz,ntr)) ; cs%tr(:,:,:,:) = 0.0
118 if (m < 10)
then ;
write(name,
'("tr_D",I1.1)') m
119 else ;
write(name,
'("tr_D",I2.2)') m ;
endif
120 write(longname,
'("Concentration of ISOMIP Tracer ",I2.2)') m
121 cs%tr_desc(m) = var_desc(name, units=
"kg kg-1", longname=longname, caller=mdl)
122 if (gv%Boussinesq)
then ; flux_units =
"kg kg-1 m3 s-1"
123 else ; flux_units =
"kg s-1" ;
endif
127 tr_ptr => cs%tr(:,:,:,m)
129 call register_tracer(tr_ptr, tr_reg, param_file, hi, gv, &
130 name=name, longname=longname, units=
"kg kg-1", &
131 registry_diags=.true., flux_units=flux_units, &
132 restart_cs=restart_cs)
137 if (cs%coupled_tracers) &
138 cs%ind_tr(m) = aof_set_coupler_flux(trim(name)//
'_flux', &
139 flux_type=
' ', implementation=
' ', caller=
"register_ISOMIP_tracer")
143 register_isomip_tracer = .true.
144 end function register_isomip_tracer
148 subroutine initialize_isomip_tracer(restart, day, G, GV, h, diag, OBC, CS, &
153 logical,
intent(in) :: restart
155 type(time_type),
target,
intent(in) :: day
156 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: h
157 type(
diag_ctrl),
target,
intent(in) :: diag
168 real,
allocatable :: temp(:,:,:)
169 real,
pointer,
dimension(:,:,:) :: &
170 obc_tr1_u => null(), &
174 character(len=16) :: name
175 character(len=72) :: longname
176 character(len=48) :: units
177 character(len=48) :: flux_units
179 real,
pointer :: tr_ptr(:,:,:) => null()
182 real :: e(szk_(g)+1), e_top, e_bot, d_tr
183 integer :: i, j, k, is, ie, js, je, isd, ied, jsd, jed, nz, m
184 integer :: isdb, iedb, jsdb, jedb
186 if (.not.
associated(cs))
return
187 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
188 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
189 isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
190 h_neglect = gv%H_subroundoff
195 if (.not.restart)
then
196 if (len_trim(cs%tracer_IC_file) >= 1)
then
198 if (.not.
file_exists(cs%tracer_IC_file, g%Domain)) &
199 call mom_error(fatal,
"ISOMIP_initialize_tracer: Unable to open "// &
202 call query_vardesc(cs%tr_desc(m), name, caller=
"initialize_ISOMIP_tracer")
203 call mom_read_data(cs%tracer_IC_file, trim(name), cs%tr(:,:,:,m), g%Domain)
207 do k=1,nz ;
do j=js,je ;
do i=is,ie
209 enddo ;
enddo ;
enddo
244 end subroutine initialize_isomip_tracer
248 subroutine isomip_tracer_column_physics(h_old, h_new, ea, eb, fluxes, dt, G, GV, US, CS, &
249 evap_CFL_limit, minimum_forcing_depth)
252 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
254 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
256 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
260 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
264 type(
forcing),
intent(in) :: fluxes
266 real,
intent(in) :: dt
270 real,
optional,
intent(in) :: evap_cfl_limit
272 real,
optional,
intent(in) :: minimum_forcing_depth
280 real :: c1(szi_(g),szk_(g))
281 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)) :: h_work
282 real :: melt(szi_(g),szj_(g))
284 character(len=256) :: mesg
285 integer :: i, j, k, is, ie, js, je, nz, m
286 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
288 if (.not.
associated(cs))
return
290 melt(:,:) = fluxes%iceshelf_melt(:,:)
293 mmax = maxval(melt(is:ie,js:je))
294 call max_across_pes(mmax)
299 do j=js,je ;
do i=is,ie
300 if (melt(i,j) > 0.0)
then
301 cs%tr(i,j,1:2,m) = melt(i,j)/mmax
303 cs%tr(i,j,1:2,m) = 0.0
308 if (
present(evap_cfl_limit) .and.
present(minimum_forcing_depth))
then
310 do k=1,nz ;
do j=js,je ;
do i=is,ie
311 h_work(i,j,k) = h_old(i,j,k)
312 enddo ;
enddo ;
enddo
313 call applytracerboundaryfluxesinout(g, gv, cs%tr(:,:,:,m), dt, fluxes, h_work, &
314 evap_cfl_limit, minimum_forcing_depth)
315 call tracer_vertdiff(h_work, ea, eb, dt, cs%tr(:,:,:,m), g, gv)
319 call tracer_vertdiff(h_old, ea, eb, dt, cs%tr(:,:,:,m), g, gv)
323 end subroutine isomip_tracer_column_physics
328 subroutine isomip_tracer_surface_state(sfc_state, h, G, CS)
330 type(
surface),
intent(inout) :: sfc_state
332 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
340 integer :: m, is, ie, js, je, isd, ied, jsd, jed
341 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
342 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
344 if (.not.
associated(cs))
return
346 if (cs%coupled_tracers)
then
350 call coupler_type_set_data(cs%tr(:,:,1,m), cs%ind_tr(m), ind_csurf, &
351 sfc_state%tr_fields, idim=(/isd, is, ie, ied/), &
352 jdim=(/jsd, js, je, jed/) )
356 end subroutine isomip_tracer_surface_state
359 subroutine isomip_tracer_end(CS)
364 if (
associated(cs))
then
365 if (
associated(cs%tr))
deallocate(cs%tr)
368 end subroutine isomip_tracer_end