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_density_integrals, only : int_specific_vol_dp
11 use mom_domains, only : pass_var, pass_vector, sum_across_pes, broadcast
12 use mom_domains, only : root_pe, to_all, scalar_pair, cgrid_ne, agrid
13 use mom_error_handler, only : mom_mesg, mom_error, fatal, warning, is_root_pe
14 use mom_error_handler, only : calltree_enter, calltree_leave, calltree_waypoint
16 use mom_file_parser, only : log_version
17 use mom_get_input, only : directories
18 use mom_grid, only : ocean_grid_type, ispointincell
20 use mom_regridding, only : regridding_cs
21 use mom_remapping, only : remapping_cs, initialize_remapping
22 use mom_remapping, only : remapping_core_h
23 use mom_string_functions, only : uppercase
26 use mom_verticalgrid, only : verticalgrid_type, setverticalgridaxes
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 :: answers_2018, default_2018_answers, hor_regrid_answers_2018
94  logical :: reentrant_x, tripolar_n
95 
96  id_clock_routine = cpu_clock_id('(Initialize tracer from Z)', grain=clock_routine)
97  id_clock_ale = cpu_clock_id('(Initialize tracer from Z) ALE', grain=clock_loop)
98 
99  call cpu_clock_begin(id_clock_routine)
100 
101  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
102  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
103 
104  call calltree_enter(trim(mdl)//"(), MOM_state_initialization.F90")
105 
106  call get_param(pf, mdl, "Z_INIT_HOMOGENIZE", homog, &
107  "If True, then horizontally homogenize the interpolated "//&
108  "initial conditions.", default=.false.)
109  call get_param(pf, mdl, "Z_INIT_ALE_REMAPPING", useale, &
110  "If True, then remap straight to model coordinate from file.",&
111  default=.true.)
112  call get_param(pf, mdl, "Z_INIT_REMAPPING_SCHEME", remapscheme, &
113  "The remapping scheme to use if using Z_INIT_ALE_REMAPPING is True.", &
114  default="PLM")
115  call get_param(pf, mdl, "DEFAULT_2018_ANSWERS", default_2018_answers, &
116  "This sets the default value for the various _2018_ANSWERS parameters.", &
117  default=.false.)
118  if (useale) then
119  call get_param(pf, mdl, "REMAPPING_2018_ANSWERS", answers_2018, &
120  "If true, use the order of arithmetic and expressions that recover the "//&
121  "answers from the end of 2018. Otherwise, use updated and more robust "//&
122  "forms of the same expressions.", default=default_2018_answers)
123  endif
124  call get_param(pf, mdl, "HOR_REGRID_2018_ANSWERS", hor_regrid_answers_2018, &
125  "If true, use the order of arithmetic for horizonal regridding that recovers "//&
126  "the answers from the end of 2018. Otherwise, use rotationally symmetric "//&
127  "forms of the same expressions.", default=default_2018_answers)
128 
129  ! These are model grid properties, but being applied to the data grid for now.
130  ! need to revisit this (mjh)
131  reentrant_x = .false. ; call get_param(pf, mdl, "REENTRANT_X", reentrant_x,default=.true.)
132  tripolar_n = .false. ; call get_param(pf, mdl, "TRIPOLAR_N", tripolar_n, default=.false.)
133 
134  if (PRESENT(homogenize)) homog=homogenize
135  if (PRESENT(usealeremapping)) useale=usealeremapping
136  if (PRESENT(remappingscheme)) remapscheme=remappingscheme
137  recnum=1
138  if (PRESENT(src_var_record)) recnum = src_var_record
139  convert=1.0
140  if (PRESENT(src_var_unit_conversion)) convert = src_var_unit_conversion
141 
142  call horiz_interp_and_extrap_tracer(src_file, src_var_nam, convert, recnum, &
143  g, tr_z, mask_z, z_in, z_edges_in, missing_value, reentrant_x, tripolar_n, &
144  homog, m_to_z=us%m_to_Z, answers_2018=hor_regrid_answers_2018)
145 
146  kd = size(z_edges_in,1)-1
147  call pass_var(tr_z,g%Domain)
148  call pass_var(mask_z,g%Domain)
149 
150 ! Done with horizontal interpolation.
151 ! Now remap to model coordinates
152  if (useale) then
153  call cpu_clock_begin(id_clock_ale)
154  ! First we reserve a work space for reconstructions of the source data
155  allocate( h1(kd) )
156  allocate( hsrc(isd:ied,jsd:jed,kd) )
157  ! Set parameters for reconstructions
158  call initialize_remapping( remapcs, remapscheme, boundary_extrapolation=.false., answers_2018=answers_2018 )
159  ! Next we initialize the regridding package so that it knows about the target grid
160 
161  do j = js, je ; do i = is, ie
162  if (g%mask2dT(i,j)>0.) then
163  ! Build the source grid
164  ztopofcell = 0. ; zbottomofcell = 0. ; npoints = 0
165  z_bathy = g%bathyT(i,j)
166  do k = 1, kd
167  if (mask_z(i,j,k) > 0.) then
168  zbottomofcell = -min( z_edges_in(k+1), z_bathy )
169  elseif (k>1) then
170  zbottomofcell = -z_bathy
171  endif
172  h1(k) = ztopofcell - zbottomofcell
173  if (h1(k)>0.) npoints = npoints + 1
174  ztopofcell = zbottomofcell ! Bottom becomes top for next value of k
175  enddo
176  h1(kd) = h1(kd) + ( ztopofcell + z_bathy ) ! In case data is deeper than model
177  else
178  tr(i,j,:) = 0.
179  endif ! mask2dT
180  hsrc(i,j,:) = gv%Z_to_H * h1(:)
181  enddo ; enddo
182 
183  call ale_remap_scalar(remapcs, g, gv, kd, hsrc, tr_z, h, tr, all_cells=.false., answers_2018=answers_2018 )
184 
185  deallocate( hsrc )
186  deallocate( h1 )
187 
188  do k=1,nz
189  call mystats(tr(:,:,k), missing_value, is, ie, js, je, k, 'Tracer from ALE()')
190  enddo
191  call cpu_clock_end(id_clock_ale)
192  endif ! useALEremapping
193 
194 ! Fill land values
195  do k=1,nz ; do j=js,je ; do i=is,ie
196  if (tr(i,j,k) == missing_value) then
197  tr(i,j,k)=land_fill
198  endif
199  enddo ; enddo ; enddo
200 
201  call calltree_leave(trim(mdl)//'()')
202  call cpu_clock_end(id_clock_routine)
203 
204 end subroutine mom_initialize_tracer_from_z
205 
This module contains the main regridding routines.
Definition: MOM_ALE.F90:9
Horizontal interpolation.
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:26
Calculates density of sea water from T, S and P.
Definition: MOM_EOS.F90:68
Reads the only Fortran name list needed to boot-strap the model.
A structure that can be parsed to read and document run-time parameters.
Provides the ocean grid type.
Definition: MOM_grid.F90:2
Generates vertical grids as part of the ALE algorithm.
Provides column-wise vertical remapping functions.
Wraps the MPP cpu clock functions.
The MOM6 facility to parse input files for runtime parameters.
An overloaded interface to log the values of various types of parameters.
Do a halo update on a pair of arrays representing the two components of a vector. ...
Definition: MOM_domains.F90:59
Container for remapping parameters.
Interfaces to non-domain-oriented communication subroutines, including the MOM6 reproducing sums faci...
Definition: MOM_coms.F90:3
Describes various unit conversion factors.
Regridding control structure.
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
Describes the decomposed MOM domain and has routines for communications across PEs.
Definition: MOM_domains.F90:2
Routines for error handling and I/O management.
Calculate the derivatives of density with temperature and salinity from T, S, and P...
Definition: MOM_EOS.F90:81
Provides subroutines for quantities specific to the equation of state.
Definition: MOM_EOS.F90:2
An overloaded interface to log version information about modules.
Describes the vertical ocean grid, including unit conversion factors.
Provides integrals of density.
An overloaded interface to read various types of parameters.
Container for paths and parameter file names.
Pointers to an assortment of thermodynamic fields that may be available, including potential temperat...
Handy functions for manipulating strings.
Provides a transparent vertical ocean grid type and supporting routines.
Provides transparent structures with groups of MOM6 variables and supporting routines.
Initializes hydrography from z-coordinate climatology files.
Do a halo update on an array.
Definition: MOM_domains.F90:54
An overloaded interface to read and log the values of various types of parameters.
Provides checksumming functions for debugging.
A control structure for the equation of state.
Definition: MOM_EOS.F90:108