MOM6
MOM_tracer_initialization_from_Z.F90
1 !> Initializes hydrography from z-coordinate climatology files
3 
4 ! This file is part of MOM6. See LICENSE.md for the license.
5 
6 use mom_debugging, only : hchksum
7 use mom_coms, only : max_across_pes, min_across_pes
8 use mom_cpu_clock, only : cpu_clock_id, cpu_clock_begin, cpu_clock_end
9 use mom_cpu_clock, only : clock_routine, clock_loop
10 use mom_domains, only : pass_var, pass_vector, sum_across_pes, broadcast
11 use mom_domains, only : root_pe, to_all, scalar_pair, cgrid_ne, agrid
12 use mom_error_handler, only : mom_mesg, mom_error, fatal, warning, is_root_pe
13 use mom_error_handler, only : calltree_enter, calltree_leave, calltree_waypoint
15 use mom_file_parser, only : log_version
16 use mom_get_input, only : directories
17 use mom_grid, only : ocean_grid_type, ispointincell
19 use mom_regridding, only : regridding_cs
20 use mom_remapping, only : remapping_cs, initialize_remapping
21 use mom_remapping, only : remapping_core_h
22 use mom_string_functions, only : uppercase
25 use mom_verticalgrid, only : verticalgrid_type, setverticalgridaxes
27 use mom_eos, only : int_specific_vol_dp
28 use mom_ale, only : ale_remap_scalar
29 
30 implicit none ; private
31 
32 #include <MOM_memory.h>
33 
34 public :: mom_initialize_tracer_from_z
35 
36 ! A note on unit descriptions in comments: MOM6 uses units that can be rescaled for dimensional
37 ! consistency testing. These are noted in comments with units like Z, H, L, and T, along with
38 ! their mks counterparts with notation like "a velocity [Z T-1 ~> m s-1]". If the units
39 ! vary with the Boussinesq approximation, the Boussinesq variant is given first.
40 
41 character(len=40) :: mdl = "MOM_tracer_initialization_from_Z" !< This module's name.
42 
43 contains
44 
45 !> Initializes a tracer from a z-space data file.
46 subroutine mom_initialize_tracer_from_z(h, tr, G, GV, US, PF, src_file, src_var_nam, &
47  src_var_unit_conversion, src_var_record, homogenize, &
48  useALEremapping, remappingScheme, src_var_gridspec )
49  type(ocean_grid_type), intent(inout) :: g !< Ocean grid structure.
50  type(verticalgrid_type), intent(in) :: gv !< Ocean vertical grid structure.
51  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
52  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
53  intent(in) :: h !< Layer thickness [H ~> m or kg m-2].
54  real, dimension(:,:,:), pointer :: tr !< Pointer to array to be initialized
55  type(param_file_type), intent(in) :: pf !< parameter file
56  character(len=*), intent(in) :: src_file !< source filename
57  character(len=*), intent(in) :: src_var_nam !< variable name in file
58  real, optional, intent(in) :: src_var_unit_conversion !< optional multiplicative unit conversion
59  integer, optional, intent(in) :: src_var_record !< record to read for multiple time-level files
60  logical, optional, intent(in) :: homogenize !< optionally homogenize to mean value
61  logical, optional, intent(in) :: usealeremapping !< to remap or not (optional)
62  character(len=*), optional, intent(in) :: remappingscheme !< remapping scheme to use.
63  character(len=*), optional, intent(in) :: src_var_gridspec !< Source variable name in a gridspec file.
64  !! This is not implemented yet.
65  ! Local variables
66  real :: land_fill = 0.0
67  character(len=200) :: inputdir ! The directory where NetCDF input files are.
68  character(len=200) :: mesg
69  real :: convert
70  integer :: recnum
71  character(len=10) :: remapscheme
72  logical :: homog,useale
73 
74  ! This include declares and sets the variable "version".
75 # include "version_variable.h"
76  character(len=40) :: mdl = "MOM_initialize_tracers_from_Z" ! This module's name.
77 
78  integer :: is, ie, js, je, nz ! compute domain indices
79  integer :: isd, ied, jsd, jed ! data domain indices
80  integer :: i, j, k, kd
81  real, allocatable, dimension(:,:,:), target :: tr_z, mask_z
82  real, allocatable, dimension(:), target :: z_edges_in, z_in
83 
84  ! Local variables for ALE remapping
85  real, dimension(:,:,:), allocatable :: hsrc ! Source thicknesses [H ~> m or kg m-2].
86  real, dimension(:), allocatable :: h1 ! A 1-d column of source thicknesses [Z ~> m].
87  real :: ztopofcell, zbottomofcell, z_bathy ! Heights [Z ~> m].
88  type(remapping_cs) :: remapcs ! Remapping parameters and work arrays
89 
90  real :: missing_value
91  integer :: npoints
92  integer :: id_clock_routine, id_clock_ale
93  logical :: reentrant_x, tripolar_n
94 
95  id_clock_routine = cpu_clock_id('(Initialize tracer from Z)', grain=clock_routine)
96  id_clock_ale = cpu_clock_id('(Initialize tracer from Z) ALE', grain=clock_loop)
97 
98  call cpu_clock_begin(id_clock_routine)
99 
100  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
101  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
102 
103  call calltree_enter(trim(mdl)//"(), MOM_state_initialization.F90")
104 
105  call get_param(pf, mdl, "Z_INIT_HOMOGENIZE", homog, &
106  "If True, then horizontally homogenize the interpolated "//&
107  "initial conditions.", default=.false.)
108  call get_param(pf, mdl, "Z_INIT_ALE_REMAPPING", useale, &
109  "If True, then remap straight to model coordinate from file.",&
110  default=.true.)
111  call get_param(pf, mdl, "Z_INIT_REMAPPING_SCHEME", remapscheme, &
112  "The remapping scheme to use if using Z_INIT_ALE_REMAPPING is True.", &
113  default="PLM")
114 
115  ! These are model grid properties, but being applied to the data grid for now.
116  ! need to revisit this (mjh)
117  reentrant_x = .false. ; call get_param(pf, mdl, "REENTRANT_X", reentrant_x,default=.true.)
118  tripolar_n = .false. ; call get_param(pf, mdl, "TRIPOLAR_N", tripolar_n, default=.false.)
119 
120  if (PRESENT(homogenize)) homog=homogenize
121  if (PRESENT(usealeremapping)) useale=usealeremapping
122  if (PRESENT(remappingscheme)) remapscheme=remappingscheme
123  recnum=1
124  if (PRESENT(src_var_record)) recnum = src_var_record
125  convert=1.0
126  if (PRESENT(src_var_unit_conversion)) convert = src_var_unit_conversion
127 
128  call horiz_interp_and_extrap_tracer(src_file, src_var_nam, convert, recnum, &
129  g, tr_z, mask_z, z_in, z_edges_in, missing_value, reentrant_x, tripolar_n, &
130  homog, m_to_z=us%m_to_Z)
131 
132  kd = size(z_edges_in,1)-1
133  call pass_var(tr_z,g%Domain)
134  call pass_var(mask_z,g%Domain)
135 
136 ! Done with horizontal interpolation.
137 ! Now remap to model coordinates
138  if (useale) then
139  call cpu_clock_begin(id_clock_ale)
140  ! First we reserve a work space for reconstructions of the source data
141  allocate( h1(kd) )
142  allocate( hsrc(isd:ied,jsd:jed,kd) )
143  call initialize_remapping( remapcs, remapscheme, boundary_extrapolation=.false. ) ! Data for reconstructions
144  ! Next we initialize the regridding package so that it knows about the target grid
145 
146  do j = js, je ; do i = is, ie
147  if (g%mask2dT(i,j)>0.) then
148  ! Build the source grid
149  ztopofcell = 0. ; zbottomofcell = 0. ; npoints = 0
150  z_bathy = g%bathyT(i,j)
151  do k = 1, kd
152  if (mask_z(i,j,k) > 0.) then
153  zbottomofcell = -min( z_edges_in(k+1), z_bathy )
154  elseif (k>1) then
155  zbottomofcell = -z_bathy
156  endif
157  h1(k) = ztopofcell - zbottomofcell
158  if (h1(k)>0.) npoints = npoints + 1
159  ztopofcell = zbottomofcell ! Bottom becomes top for next value of k
160  enddo
161  h1(kd) = h1(kd) + ( ztopofcell + z_bathy ) ! In case data is deeper than model
162  else
163  tr(i,j,:) = 0.
164  endif ! mask2dT
165  hsrc(i,j,:) = gv%Z_to_H * h1(:)
166  enddo ; enddo
167 
168  call ale_remap_scalar(remapcs, g, gv, kd, hsrc, tr_z, h, tr, all_cells=.false. )
169 
170  deallocate( hsrc )
171  deallocate( h1 )
172 
173  do k=1,nz
174  call mystats(tr(:,:,k), missing_value, is, ie, js, je, k, 'Tracer from ALE()')
175  enddo
176  call cpu_clock_end(id_clock_ale)
177  endif ! useALEremapping
178 
179 ! Fill land values
180  do k=1,nz ; do j=js,je ; do i=is,ie
181  if (tr(i,j,k) == missing_value) then
182  tr(i,j,k)=land_fill
183  endif
184  enddo ; enddo ; enddo
185 
186  call calltree_leave(trim(mdl)//'()')
187  call cpu_clock_end(id_clock_routine)
188 
189 end subroutine mom_initialize_tracer_from_z
190 
mom_regridding::regridding_cs
Regridding control structure.
Definition: MOM_regridding.F90:45
mom_horizontal_regridding::horiz_interp_and_extrap_tracer
Extrapolate and interpolate data.
Definition: MOM_horizontal_regridding.F90:52
mom_verticalgrid
Provides a transparent vertical ocean grid type and supporting routines.
Definition: MOM_verticalGrid.F90:2
mom_file_parser::log_version
An overloaded interface to log version information about modules.
Definition: MOM_file_parser.F90:109
mom_variables::thermo_var_ptrs
Pointers to an assortment of thermodynamic fields that may be available, including potential temperat...
Definition: MOM_variables.F90:82
mom_get_input::directories
Container for paths and parameter file names.
Definition: MOM_get_input.F90:20
mom_string_functions
Handy functions for manipulating strings.
Definition: MOM_string_functions.F90:2
mom_eos
Provides subroutines for quantities specific to the equation of state.
Definition: MOM_EOS.F90:2
mom_domains::pass_var
Do a halo update on an array.
Definition: MOM_domains.F90:49
mom_remapping::remapping_cs
Container for remapping parameters.
Definition: MOM_remapping.F90:22
mom_file_parser::param_file_type
A structure that can be parsed to read and document run-time parameters.
Definition: MOM_file_parser.F90:54
mom_file_parser::get_param
An overloaded interface to read and log the values of various types of parameters.
Definition: MOM_file_parser.F90:102
mom_tracer_initialization_from_z
Initializes hydrography from z-coordinate climatology files.
Definition: MOM_tracer_initialization_from_Z.F90:2
mom_horizontal_regridding
Horizontal interpolation.
Definition: MOM_horizontal_regridding.F90:2
mom_get_input
Reads the only Fortran name list needed to boot-strap the model.
Definition: MOM_get_input.F90:6
mom_unit_scaling::unit_scale_type
Describes various unit conversion factors.
Definition: MOM_unit_scaling.F90:14
mom_coms
Interfaces to non-domain-oriented communication subroutines, including the MOM6 reproducing sums faci...
Definition: MOM_coms.F90:3
mom_remapping
Provides column-wise vertical remapping functions.
Definition: MOM_remapping.F90:2
mom_domains::pass_vector
Do a halo update on a pair of arrays representing the two components of a vector.
Definition: MOM_domains.F90:54
mom_eos::eos_type
A control structure for the equation of state.
Definition: MOM_EOS.F90:86
mom_verticalgrid::verticalgrid_type
Describes the vertical ocean grid, including unit conversion factors.
Definition: MOM_verticalGrid.F90:24
mom_domains
Describes the decomposed MOM domain and has routines for communications across PEs.
Definition: MOM_domains.F90:2
mom_ale
This module contains the main regridding routines.
Definition: MOM_ALE.F90:9
mom_variables
Provides transparent structures with groups of MOM6 variables and supporting routines.
Definition: MOM_variables.F90:2
mom_cpu_clock
Wraps the MPP cpu clock functions.
Definition: MOM_cpu_clock.F90:2
mom_file_parser
The MOM6 facility to parse input files for runtime parameters.
Definition: MOM_file_parser.F90:2
mom_eos::calculate_density_derivs
Calculate the derivatives of density with temperature and salinity from T, S, and P.
Definition: MOM_EOS.F90:70
mom_grid
Provides the ocean grid type.
Definition: MOM_grid.F90:2
mom_unit_scaling
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
Definition: MOM_unit_scaling.F90:2
mom_debugging
Provides checksumming functions for debugging.
Definition: MOM_debugging.F90:7
mom_regridding
Generates vertical grids as part of the ALE algorithm.
Definition: MOM_regridding.F90:2
mom_file_parser::log_param
An overloaded interface to log the values of various types of parameters.
Definition: MOM_file_parser.F90:96
mom_error_handler
Routines for error handling and I/O management.
Definition: MOM_error_handler.F90:2
mom_grid::ocean_grid_type
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:25
mom_eos::calculate_density
Calculates density of sea water from T, S and P.
Definition: MOM_EOS.F90:60
mom_file_parser::read_param
An overloaded interface to read various types of parameters.
Definition: MOM_file_parser.F90:90