Skip to content

Commit 59d3e3e

Browse files
Mariona Claretjkrasting
authored andcommitted
Merged mlres tracer from user/mc/bling13
- Taken from internal repository at commit ed2f242
1 parent 8dda39f commit 59d3e3e

File tree

2 files changed

+421
-7
lines changed

2 files changed

+421
-7
lines changed

generic_tracers/generic_mlres.F90

Lines changed: 395 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,395 @@
1+
!----------------------------------------------------------------
2+
! <CONTACT EMAIL="[email protected]"> Niki Zadeh
3+
! </CONTACT>
4+
!
5+
! <REVIEWER EMAIL="[email protected]"> William Cooke
6+
! </REVIEWER>
7+
!<OVERVIEW>
8+
!</OVERVIEW>
9+
!
10+
!<DESCRIPTION>
11+
! Implementation of routines to solve the amount of time that a water parcel
12+
! at any give location has last resided in the mixed layer
13+
!</DESCRIPTION>
14+
!
15+
! </DEVELOPER_NOTES>
16+
! </INFO>
17+
!
18+
!----------------------------------------------------------------
19+
20+
module generic_mlres
21+
22+
use coupler_types_mod, only: coupler_2d_bc_type
23+
use field_manager_mod, only: fm_string_len
24+
!use mpp_mod, only : mpp_error, NOTE, WARNING, FATAL, stdout
25+
use mpp_mod, only: input_nml_file, mpp_error, stdlog, NOTE, WARNING, FATAL, stdout, mpp_chksum
26+
use fms_mod, only: write_version_number, open_namelist_file, check_nml_error, close_file
27+
use time_manager_mod, only : time_type
28+
use fm_util_mod, only: fm_util_start_namelist, fm_util_end_namelist
29+
30+
use g_tracer_utils, only : g_tracer_type,g_tracer_start_param_list,g_tracer_end_param_list
31+
use g_tracer_utils, only : g_tracer_add,g_tracer_add_param, g_tracer_set_files
32+
use g_tracer_utils, only : g_tracer_set_values,g_tracer_get_pointer,g_tracer_get_common
33+
use g_tracer_utils, only : g_tracer_coupler_set,g_tracer_coupler_get
34+
use g_tracer_utils, only : g_tracer_send_diag, g_tracer_get_values
35+
36+
37+
implicit none ; private
38+
39+
character(len=fm_string_len), parameter :: mod_name = 'generic_mlres'
40+
character(len=fm_string_len), parameter :: package_name = 'generic_mlres'
41+
42+
public do_generic_mlres
43+
public generic_mlres_register
44+
public generic_mlres_init
45+
public generic_mlres_update_from_coupler
46+
public generic_mlres_update_from_source
47+
public generic_mlres_set_boundary_values
48+
public generic_mlres_end
49+
50+
!The following logical for using this module is overwritten
51+
! by generic_tracer_nml namelist
52+
logical, save :: do_generic_mlres = .false.
53+
54+
real, parameter :: epsln=1.0e-30
55+
56+
real :: reset_time=1.0
57+
58+
namelist /generic_mlres_nml/ reset_time
59+
!
60+
!This type contains all the parameters and arrays used in this module.
61+
!
62+
!Note that there is no programatic reason for treating
63+
!the following as a type. These are the parameters used only in this module.
64+
!It suffices for varables to be a declared at the top of the module.
65+
!nnz: Find out about the timing overhead for using type%x rather than x
66+
67+
type generic_mlres_params
68+
real :: Rho_0
69+
character(len=fm_string_len) :: ice_restart_file
70+
character(len=fm_string_len) :: ocean_restart_file,IC_file
71+
end type generic_mlres_params
72+
73+
74+
type(generic_mlres_params) :: param
75+
76+
contains
77+
78+
subroutine generic_mlres_register(tracer_list)
79+
type(g_tracer_type), pointer :: tracer_list
80+
81+
integer :: ioun
82+
integer :: ierr
83+
integer :: io_status
84+
character(len=fm_string_len) :: name
85+
integer :: stdoutunit,stdlogunit
86+
87+
character(len=fm_string_len), parameter :: sub_name = 'generic_mlres_register'
88+
89+
! provide for namelist over-ride
90+
! This needs to go before the add_tracers in order to allow the namelist
91+
! settings to switch tracers on and off.
92+
!
93+
stdoutunit=stdout();stdlogunit=stdlog()
94+
95+
#ifdef INTERNAL_FILE_NML
96+
read (input_nml_file, nml=generic_mlres_nml, iostat=io_status)
97+
ierr = check_nml_error(io_status,'generic_mlres_nml')
98+
#else
99+
ioun = open_namelist_file()
100+
read (ioun, generic_bling_nml,iostat=io_status)
101+
ierr = check_nml_error(io_status,'generic_mlres_nml')
102+
call close_file (ioun)
103+
#endif
104+
105+
write (stdoutunit,'(/)')
106+
write (stdoutunit, generic_mlres_nml)
107+
write (stdlogunit, generic_mlres_nml)
108+
109+
!Specify all prognostic and diagnostic tracers of this modules.
110+
call user_add_tracers(tracer_list)
111+
112+
113+
end subroutine generic_mlres_register
114+
115+
! <SUBROUTINE NAME="generic_mlres_init">
116+
! <OVERVIEW>
117+
! Initialize the generic mixed layer residence tracer
118+
! </OVERVIEW>
119+
! <DESCRIPTION>
120+
! This subroutine:
121+
! Adds all the residence tracers to the list of generic Tracers passed to it via utility subroutine g_tracer_add().
122+
! Adds all the parameters used by this module via utility subroutine g_tracer_add_param().
123+
! Allocates all work arrays used in the module.
124+
! </DESCRIPTION>
125+
! <TEMPLATE>
126+
! call generic_mlres_init(tracer_list)
127+
! </TEMPLATE>
128+
! <IN NAME="tracer_list" TYPE="type(g_tracer_type), pointer">
129+
! Pointer to the head of generic tracer list.
130+
! </IN>
131+
! </SUBROUTINE>
132+
133+
subroutine generic_mlres_init(tracer_list)
134+
type(g_tracer_type), pointer :: tracer_list
135+
136+
character(len=fm_string_len), parameter :: sub_name = 'generic_mlres_init'
137+
138+
!Specify and initialize all parameters used by this package
139+
call user_add_params
140+
141+
!Allocate and initiate all the private work arrays used by this module.
142+
!call user_allocate_arrays
143+
144+
end subroutine generic_mlres_init
145+
146+
subroutine user_allocate_arrays
147+
148+
! Nothing to do for mixed layer tracer
149+
! call g_tracer_get_common(isc,iec,jsc,jec,isd,ied,jsd,jed,nk,ntau)
150+
151+
end subroutine user_allocate_arrays
152+
153+
subroutine user_deallocate_arrays
154+
! This is an internal sub, not a public interface.
155+
! Deallocate all the work arrays allocated by user_allocate_arrays.
156+
!
157+
end subroutine user_deallocate_arrays
158+
!
159+
! This is an internal sub, not a public interface.
160+
! Add all the parameters to be used in this module.
161+
!
162+
subroutine user_add_params
163+
164+
!Specify all parameters used in this modules.
165+
!==============================================================
166+
!User adds one call for each parameter below!
167+
!User also adds the definition of each parameter in generic_mlres_params type
168+
!==============================================================
169+
170+
!=============
171+
!Block Starts: g_tracer_add_param
172+
!=============
173+
!Add the known experimental parameters used for calculations
174+
!in this module.
175+
!All the g_tracer_add_param calls must happen between
176+
!g_tracer_start_param_list and g_tracer_end_param_list calls.
177+
!This implementation enables runtime overwrite via field_table.
178+
179+
call g_tracer_start_param_list(package_name)
180+
! Rho_0 is used in the Boussinesq
181+
! approximation to calculations of pressure and
182+
! pressure gradients, in units of kg m-3.
183+
call g_tracer_add_param('RHO_0', param%Rho_0, 1035.0)
184+
185+
call g_tracer_end_param_list(package_name)
186+
!===========
187+
!Block Ends: g_tracer_add_param
188+
!===========
189+
190+
end subroutine user_add_params
191+
192+
!
193+
! This is an internal sub, not a public interface.
194+
! Add all the tracers to be used in this module.
195+
!
196+
subroutine user_add_tracers(tracer_list)
197+
type(g_tracer_type), pointer :: tracer_list
198+
199+
character(len=fm_string_len), parameter :: sub_name = 'user_add_tracers'
200+
201+
call g_tracer_start_param_list(package_name)!nnz: Does this append?
202+
call g_tracer_add_param('ice_restart_file' , param%ice_restart_file , 'ice_ocmip_mlres.res.nc')
203+
call g_tracer_add_param('ocean_restart_file' , param%ocean_restart_file , 'ocmip_mlres.res.nc' )
204+
call g_tracer_add_param('IC_file' , param%IC_file , '')
205+
call g_tracer_end_param_list(package_name)
206+
207+
! Set Restart files
208+
call g_tracer_set_files(ice_restart_file=param%ice_restart_file, ocean_restart_file=param%ocean_restart_file )
209+
210+
!=====================================================
211+
!Specify all prognostic tracers of this modules.
212+
!=====================================================
213+
!User adds one call for each prognostic tracer below!
214+
!User should specify if fluxes must be extracted from boundary
215+
!by passing one or more of the following methods as .true.
216+
!and provide the corresponding parameters array
217+
!methods: flux_gas,flux_runoff,flux_wetdep,flux_drydep
218+
!
219+
!prog_tracers: mlres
220+
!diag_tracers: none
221+
!
222+
!age
223+
call g_tracer_add(tracer_list,package_name, &
224+
name = 'mlres', &
225+
longname = 'residence time inside mixed layer', &
226+
units = 'years', &
227+
init_value = 0.0, &
228+
prog = .true.)
229+
230+
call g_tracer_add(tracer_list,package_name, &
231+
name = 'mlres_inv', &
232+
longname = 'residence time outside mixed layer', &
233+
units = 'years', &
234+
init_value = 0.0, &
235+
prog = .true.)
236+
237+
238+
end subroutine user_add_tracers
239+
240+
! <SUBROUTINE NAME="generic_mlres_update_from_coupler">
241+
! <OVERVIEW>
242+
! Modify the values obtained from the coupler if necessary.
243+
! </OVERVIEW>
244+
! <DESCRIPTION>
245+
! Currently an empty stub for ages.
246+
! Some tracer fields need to be modified after values are obtained from the coupler.
247+
! This subroutine is the place for specific tracer manipulations.
248+
! </DESCRIPTION>
249+
! <TEMPLATE>
250+
! call generic_mlres_update_from_coupler(tracer_list)
251+
! </TEMPLATE>
252+
! <IN NAME="tracer_list" TYPE="type(g_tracer_type), pointer">
253+
! Pointer to the head of generic tracer list.
254+
! </IN>
255+
! </SUBROUTINE>
256+
subroutine generic_mlres_update_from_coupler(tracer_list)
257+
type(g_tracer_type), pointer :: tracer_list
258+
character(len=fm_string_len), parameter :: sub_name = 'generic_mlres_update_from_coupler'
259+
!
260+
!Nothing specific to be done for mixed layer residence tracers
261+
!
262+
return
263+
end subroutine generic_mlres_update_from_coupler
264+
265+
! <SUBROUTINE NAME="generic_mlres_update_from_source">
266+
! <OVERVIEW>
267+
! Update tracer concentration fields due to the source/sink contributions.
268+
! </OVERVIEW>
269+
! <DESCRIPTION>
270+
! Sets age to zero in uppermost level and increments it by the time step elsewhere
271+
! </DESCRIPTION>
272+
! </SUBROUTINE>
273+
subroutine generic_mlres_update_from_source( tracer_list,tau,dt, &
274+
hblt_depth,dzt,ilb,jlb )
275+
276+
type(g_tracer_type), pointer :: tracer_list
277+
integer, intent(in) :: tau, ilb, jlb
278+
real, intent(in) :: dt
279+
real, dimension(ilb:,jlb:), intent(in) :: hblt_depth
280+
real, dimension(ilb:,jlb:,:), intent(in) :: dzt
281+
282+
character(len=fm_string_len), parameter :: sub_name = 'generic_mlres_update_from_source'
283+
integer :: isc,iec,jsc,jec,isd,ied,jsd,jed,nk,ntau,i,j,k
284+
real, parameter :: secs_in_year_r = 1.0 / (86400.0 * 365.25)
285+
real, dimension(:,:,:) , pointer :: grid_tmask
286+
real, dimension(:,:,:,:), pointer :: p_mlres_field, p_mlres_inv_field
287+
288+
real, dimension(:,:,:), allocatable :: zt
289+
290+
call g_tracer_get_common(isc,iec,jsc,jec,isd,ied,jsd,jed,nk,ntau,grid_tmask=grid_tmask)
291+
292+
allocate( zt(isd:ied,jsd:jed,nk) ); zt=0.0
293+
294+
call g_tracer_get_pointer(tracer_list, 'mlres' , 'field', p_mlres_field )
295+
call g_tracer_get_pointer(tracer_list, 'mlres_inv', 'field', p_mlres_inv_field)
296+
297+
! Compute depth of the bottom of the cell
298+
zt = 0.0
299+
do j = jsc, jec ; do i = isc, iec !{
300+
zt(i,j,1) = dzt(i,j,1)
301+
enddo; enddo !} i,j
302+
do k = 2, nk ; do j = jsc, jec ; do i = isc, iec !{
303+
zt(i,j,k) = zt(i,j,k-1) + dzt(i,j,k)
304+
enddo; enddo; enddo !} i,j
305+
306+
do k = 1, nk; do j = jsc, jec ; do i = isc, iec
307+
308+
if (zt(i,j,k) .le. hblt_depth(i,j) ) then
309+
!if ( p_mlres_inv_field(i,j,k,tau) .gt. 1 ) then
310+
if ( p_mlres_inv_field(i,j,k,tau) .gt. reset_time ) then
311+
! Reset tracers
312+
p_mlres_field (i,j,k,tau) = 0.0
313+
p_mlres_inv_field(i,j,k,tau) = 0.0
314+
else
315+
p_mlres_field(i,j,k,tau) = p_mlres_field(i,j,k,tau) + dt*secs_in_year_r*grid_tmask(i,j,k)
316+
endif
317+
else
318+
p_mlres_inv_field(i,j,k,tau) = p_mlres_inv_field(i,j,k,tau) + dt*secs_in_year_r*grid_tmask(i,j,k)
319+
endif
320+
321+
enddo; enddo; enddo !} i,j,k
322+
323+
deallocate (zt)
324+
325+
return
326+
end subroutine generic_mlres_update_from_source
327+
328+
! <SUBROUTINE NAME="generic_mlres_set_boundary_values">
329+
! <OVERVIEW>
330+
! Calculate and set coupler values at the surface / bottom
331+
! </OVERVIEW>
332+
! <DESCRIPTION>
333+
!
334+
! </DESCRIPTION>
335+
! <TEMPLATE>
336+
! call generic_mlres_set_boundary_values(tracer_list,SST,SSS,rho,ilb,jlb,tau)
337+
! </TEMPLATE>
338+
! <IN NAME="tracer_list" TYPE="type(g_tracer_type), pointer">
339+
! Pointer to the head of generic tracer list.
340+
! </IN>
341+
! <IN NAME="ilb,jlb" TYPE="integer">
342+
! Lower bounds of x and y extents of input arrays on data domain
343+
! </IN>
344+
! <IN NAME="SST" TYPE="real, dimension(ilb:,jlb:)">
345+
! Sea Surface Temperature
346+
! </IN>
347+
! <IN NAME="SSS" TYPE="real, dimension(ilb:,jlb:)">
348+
! Sea Surface Salinity
349+
! </IN>
350+
! <IN NAME="rho" TYPE="real, dimension(ilb:,jlb:,:,:)">
351+
! Ocean density
352+
! </IN>
353+
! <IN NAME="tau" TYPE="integer">
354+
! Time step index of %field
355+
! </IN>
356+
! </SUBROUTINE>
357+
358+
!User must provide the calculations for these boundary values.
359+
subroutine generic_mlres_set_boundary_values(tracer_list,ST,SSS,rho,ilb,jlb,taum1)
360+
type(g_tracer_type), pointer :: tracer_list
361+
real, dimension(ilb:,jlb:), intent(in) :: ST, SSS
362+
real, dimension(ilb:,jlb:,:,:), intent(in) :: rho
363+
integer, intent(in) :: ilb,jlb,taum1
364+
365+
integer :: isc,iec, jsc,jec,isd,ied,jsd,jed,nk,ntau , i, j
366+
real :: conv_fac,sal,ta,SST,alpha,sc
367+
real, dimension(:,:,:) ,pointer :: grid_tmask
368+
real, dimension(:,:,:,:), pointer :: g_mlres_field
369+
real, dimension(:,:), ALLOCATABLE :: g_mlres_alpha,g_mlres_csurf
370+
real, dimension(:,:), ALLOCATABLE :: sc_no
371+
372+
character(len=fm_string_len), parameter :: sub_name = 'generic_mlres_set_boundary_values'
373+
374+
return
375+
end subroutine generic_mlres_set_boundary_values
376+
377+
! <SUBROUTINE NAME="generic_mlres_end">
378+
! <OVERVIEW>
379+
! End the module.
380+
! </OVERVIEW>
381+
! <DESCRIPTION>
382+
! Deallocate all work arrays
383+
! </DESCRIPTION>
384+
! <TEMPLATE>
385+
! call generic_mlres_end
386+
! </TEMPLATE>
387+
! </SUBROUTINE>
388+
389+
subroutine generic_mlres_end
390+
character(len=fm_string_len), parameter :: sub_name = 'generic_mlres_end'
391+
!call user_deallocate_arrays
392+
393+
end subroutine generic_mlres_end
394+
395+
end module generic_mlres

0 commit comments

Comments
 (0)