MOM6
mom_dyn_horgrid Module Reference

Data Types

type  dyn_horgrid_type
 

Functions/Subroutines

subroutine, public create_dyn_horgrid (G, HI, bathymetry_at_vel)
 Allocate memory used by the dyn_horgrid_type and related structures. More...
 
subroutine, public set_derived_dyn_horgrid (G)
 set_derived_dyn_horgrid calculates metric terms that are derived from other metrics. More...
 
real function adcroft_reciprocal (val)
 Adcroft_reciprocal(x) = 1/x for |x|>0 or 0 for x=0. More...
 
subroutine, public destroy_dyn_horgrid (G)
 Release memory used by the dyn_horgrid_type and related structures. More...
 

Function/Subroutine Documentation

◆ adcroft_reciprocal()

real function mom_dyn_horgrid::adcroft_reciprocal ( real, intent(in)  val)
private

Adcroft_reciprocal(x) = 1/x for |x|>0 or 0 for x=0.

Parameters
[in]valThe value being inverted.
Returns
The Adcroft reciprocal of val.

Definition at line 305 of file MOM_dyn_horgrid.F90.

Referenced by set_derived_dyn_horgrid().

305  real, intent(in) :: val !< The value being inverted.
306  real :: i_val !< The Adcroft reciprocal of val.
307 
308  i_val = 0.0 ; if (val /= 0.0) i_val = 1.0/val
Here is the caller graph for this function:

◆ create_dyn_horgrid()

subroutine, public mom_dyn_horgrid::create_dyn_horgrid ( type(dyn_horgrid_type), pointer  G,
type(hor_index_type), intent(in)  HI,
logical, intent(in), optional  bathymetry_at_vel 
)

Allocate memory used by the dyn_horgrid_type and related structures.

Parameters
gA pointer to the dynamic horizontal grid type
[in]hiA hor_index_type for array extents
[in]bathymetry_at_velIf true, there are separate values for the basin depths at velocity points. Otherwise the effects of topography are entirely determined from thickness points.

Definition at line 149 of file MOM_dyn_horgrid.F90.

References mom_error_handler::mom_error().

Referenced by mom::initialize_mom().

149  type(dyn_horgrid_type), pointer :: g !< A pointer to the dynamic horizontal grid type
150  type(hor_index_type), intent(in) :: hi !< A hor_index_type for array extents
151  logical, optional, intent(in) :: bathymetry_at_vel !< If true, there are
152  !! separate values for the basin depths at velocity
153  !! points. Otherwise the effects of topography are
154  !! entirely determined from thickness points.
155  integer :: isd, ied, jsd, jed, isdb, iedb, jsdb, jedb, isg, ieg, jsg, jeg
156 
157  ! This subroutine allocates the lateral elements of the dyn_horgrid_type that
158  ! are always used and zeros them out.
159 
160  if (associated(g)) then
161  call mom_error(warning, "create_dyn_horgrid called with an associated horgrid_type.")
162  else
163  allocate(g)
164  endif
165 
166  g%HI = hi
167 
168  g%isc = hi%isc ; g%iec = hi%iec ; g%jsc = hi%jsc ; g%jec = hi%jec
169  g%isd = hi%isd ; g%ied = hi%ied ; g%jsd = hi%jsd ; g%jed = hi%jed
170  g%isg = hi%isg ; g%ieg = hi%ieg ; g%jsg = hi%jsg ; g%jeg = hi%jeg
171 
172  g%IscB = hi%IscB ; g%IecB = hi%IecB ; g%JscB = hi%JscB ; g%JecB = hi%JecB
173  g%IsdB = hi%IsdB ; g%IedB = hi%IedB ; g%JsdB = hi%JsdB ; g%JedB = hi%JedB
174  g%IsgB = hi%IsgB ; g%IegB = hi%IegB ; g%JsgB = hi%JsgB ; g%JegB = hi%JegB
175 
176  g%idg_offset = hi%idg_offset ; g%jdg_offset = hi%jdg_offset
177  g%isd_global = g%isd + hi%idg_offset ; g%jsd_global = g%jsd + hi%jdg_offset
178  g%symmetric = hi%symmetric
179 
180  g%bathymetry_at_vel = .false.
181  if (present(bathymetry_at_vel)) g%bathymetry_at_vel = bathymetry_at_vel
182 
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
185  isg = g%isg ; ieg = g%ieg ; jsg = g%jsg ; jeg = g%jeg
186 
187  allocate(g%dxT(isd:ied,jsd:jed)) ; g%dxT(:,:) = 0.0
188  allocate(g%dxCu(isdb:iedb,jsd:jed)) ; g%dxCu(:,:) = 0.0
189  allocate(g%dxCv(isd:ied,jsdb:jedb)) ; g%dxCv(:,:) = 0.0
190  allocate(g%dxBu(isdb:iedb,jsdb:jedb)) ; g%dxBu(:,:) = 0.0
191  allocate(g%IdxT(isd:ied,jsd:jed)) ; g%IdxT(:,:) = 0.0
192  allocate(g%IdxCu(isdb:iedb,jsd:jed)) ; g%IdxCu(:,:) = 0.0
193  allocate(g%IdxCv(isd:ied,jsdb:jedb)) ; g%IdxCv(:,:) = 0.0
194  allocate(g%IdxBu(isdb:iedb,jsdb:jedb)) ; g%IdxBu(:,:) = 0.0
195 
196  allocate(g%dyT(isd:ied,jsd:jed)) ; g%dyT(:,:) = 0.0
197  allocate(g%dyCu(isdb:iedb,jsd:jed)) ; g%dyCu(:,:) = 0.0
198  allocate(g%dyCv(isd:ied,jsdb:jedb)) ; g%dyCv(:,:) = 0.0
199  allocate(g%dyBu(isdb:iedb,jsdb:jedb)) ; g%dyBu(:,:) = 0.0
200  allocate(g%IdyT(isd:ied,jsd:jed)) ; g%IdyT(:,:) = 0.0
201  allocate(g%IdyCu(isdb:iedb,jsd:jed)) ; g%IdyCu(:,:) = 0.0
202  allocate(g%IdyCv(isd:ied,jsdb:jedb)) ; g%IdyCv(:,:) = 0.0
203  allocate(g%IdyBu(isdb:iedb,jsdb:jedb)) ; g%IdyBu(:,:) = 0.0
204 
205  allocate(g%areaT(isd:ied,jsd:jed)) ; g%areaT(:,:) = 0.0
206  allocate(g%IareaT(isd:ied,jsd:jed)) ; g%IareaT(:,:) = 0.0
207  allocate(g%areaBu(isdb:iedb,jsdb:jedb)) ; g%areaBu(:,:) = 0.0
208  allocate(g%IareaBu(isdb:iedb,jsdb:jedb)) ; g%IareaBu(:,:) = 0.0
209 
210  allocate(g%mask2dT(isd:ied,jsd:jed)) ; g%mask2dT(:,:) = 0.0
211  allocate(g%mask2dCu(isdb:iedb,jsd:jed)) ; g%mask2dCu(:,:) = 0.0
212  allocate(g%mask2dCv(isd:ied,jsdb:jedb)) ; g%mask2dCv(:,:) = 0.0
213  allocate(g%mask2dBu(isdb:iedb,jsdb:jedb)) ; g%mask2dBu(:,:) = 0.0
214  allocate(g%geoLatT(isd:ied,jsd:jed)) ; g%geoLatT(:,:) = 0.0
215  allocate(g%geoLatCu(isdb:iedb,jsd:jed)) ; g%geoLatCu(:,:) = 0.0
216  allocate(g%geoLatCv(isd:ied,jsdb:jedb)) ; g%geoLatCv(:,:) = 0.0
217  allocate(g%geoLatBu(isdb:iedb,jsdb:jedb)) ; g%geoLatBu(:,:) = 0.0
218  allocate(g%geoLonT(isd:ied,jsd:jed)) ; g%geoLonT(:,:) = 0.0
219  allocate(g%geoLonCu(isdb:iedb,jsd:jed)) ; g%geoLonCu(:,:) = 0.0
220  allocate(g%geoLonCv(isd:ied,jsdb:jedb)) ; g%geoLonCv(:,:) = 0.0
221  allocate(g%geoLonBu(isdb:iedb,jsdb:jedb)) ; g%geoLonBu(:,:) = 0.0
222 
223  allocate(g%dx_Cv(isd:ied,jsdb:jedb)) ; g%dx_Cv(:,:) = 0.0
224  allocate(g%dy_Cu(isdb:iedb,jsd:jed)) ; g%dy_Cu(:,:) = 0.0
225  allocate(g%dx_Cv_obc(isd:ied,jsdb:jedb)) ; g%dx_Cv_obc(:,:) = 0.0
226  allocate(g%dy_Cu_obc(isdb:iedb,jsd:jed)) ; g%dy_Cu_obc(:,:) = 0.0
227 
228  allocate(g%areaCu(isdb:iedb,jsd:jed)) ; g%areaCu(:,:) = 0.0
229  allocate(g%areaCv(isd:ied,jsdb:jedb)) ; g%areaCv(:,:) = 0.0
230  allocate(g%IareaCu(isdb:iedb,jsd:jed)) ; g%IareaCu(:,:) = 0.0
231  allocate(g%IareaCv(isd:ied,jsdb:jedb)) ; g%IareaCv(:,:) = 0.0
232 
233  allocate(g%bathyT(isd:ied, jsd:jed)) ; g%bathyT(:,:) = 0.0
234  allocate(g%CoriolisBu(isdb:iedb, jsdb:jedb)) ; g%CoriolisBu(:,:) = 0.0
235  allocate(g%dF_dx(isd:ied, jsd:jed)) ; g%dF_dx(:,:) = 0.0
236  allocate(g%dF_dy(isd:ied, jsd:jed)) ; g%dF_dy(:,:) = 0.0
237 
238  allocate(g%sin_rot(isd:ied,jsd:jed)) ; g%sin_rot(:,:) = 0.0
239  allocate(g%cos_rot(isd:ied,jsd:jed)) ; g%cos_rot(:,:) = 1.0
240 
241  if (g%bathymetry_at_vel) then
242  allocate(g%Dblock_u(isdb:iedb, jsd:jed)) ; g%Dblock_u(:,:) = 0.0
243  allocate(g%Dopen_u(isdb:iedb, jsd:jed)) ; g%Dopen_u(:,:) = 0.0
244  allocate(g%Dblock_v(isd:ied, jsdb:jedb)) ; g%Dblock_v(:,:) = 0.0
245  allocate(g%Dopen_v(isd:ied, jsdb:jedb)) ; g%Dopen_v(:,:) = 0.0
246  endif
247 
248  ! gridLonB and gridLatB are used as edge values in some cases, so they
249  ! always need to use symmetric memory allcoations.
250  allocate(g%gridLonT(isg:ieg)) ; g%gridLonT(:) = 0.0
251  allocate(g%gridLonB(isg-1:ieg)) ; g%gridLonB(:) = 0.0
252  allocate(g%gridLatT(jsg:jeg)) ; g%gridLatT(:) = 0.0
253  allocate(g%gridLatB(jsg-1:jeg)) ; g%gridLatB(:) = 0.0
254 
Here is the call graph for this function:
Here is the caller graph for this function:

◆ destroy_dyn_horgrid()

subroutine, public mom_dyn_horgrid::destroy_dyn_horgrid ( type(dyn_horgrid_type), pointer  G)

Release memory used by the dyn_horgrid_type and related structures.

Parameters
gThe dynamic horizontal grid type

Definition at line 314 of file MOM_dyn_horgrid.F90.

References mom_error_handler::mom_error().

Referenced by mom::initialize_mom().

314  type(dyn_horgrid_type), pointer :: g !< The dynamic horizontal grid type
315 
316  if (.not.associated(g)) then
317  call mom_error(fatal, "destroy_dyn_horgrid called with an unassociated horgrid_type.")
318  endif
319 
320  deallocate(g%dxT) ; deallocate(g%dxCu) ; deallocate(g%dxCv) ; deallocate(g%dxBu)
321  deallocate(g%IdxT) ; deallocate(g%IdxCu) ; deallocate(g%IdxCv) ; deallocate(g%IdxBu)
322 
323  deallocate(g%dyT) ; deallocate(g%dyCu) ; deallocate(g%dyCv) ; deallocate(g%dyBu)
324  deallocate(g%IdyT) ; deallocate(g%IdyCu) ; deallocate(g%IdyCv) ; deallocate(g%IdyBu)
325 
326  deallocate(g%areaT) ; deallocate(g%IareaT)
327  deallocate(g%areaBu) ; deallocate(g%IareaBu)
328  deallocate(g%areaCu) ; deallocate(g%IareaCu)
329  deallocate(g%areaCv) ; deallocate(g%IareaCv)
330 
331  deallocate(g%mask2dT) ; deallocate(g%mask2dCu)
332  deallocate(g%mask2dCv) ; deallocate(g%mask2dBu)
333 
334  deallocate(g%geoLatT) ; deallocate(g%geoLatCu)
335  deallocate(g%geoLatCv) ; deallocate(g%geoLatBu)
336  deallocate(g%geoLonT) ; deallocate(g%geoLonCu)
337  deallocate(g%geoLonCv) ; deallocate(g%geoLonBu)
338 
339  deallocate(g%dx_Cv) ; deallocate(g%dy_Cu)
340  deallocate(g%dx_Cv_obc) ; deallocate(g%dy_Cu_obc)
341 
342  deallocate(g%bathyT) ; deallocate(g%CoriolisBu)
343  deallocate(g%dF_dx) ; deallocate(g%dF_dy)
344  deallocate(g%sin_rot) ; deallocate(g%cos_rot)
345 
346  if (allocated(g%Dblock_u)) deallocate(g%Dblock_u)
347  if (allocated(g%Dopen_u)) deallocate(g%Dopen_u)
348  if (allocated(g%Dblock_v)) deallocate(g%Dblock_v)
349  if (allocated(g%Dopen_v)) deallocate(g%Dopen_v)
350 
351  deallocate(g%gridLonT) ; deallocate(g%gridLatT)
352  deallocate(g%gridLonB) ; deallocate(g%gridLatB)
353 
354  deallocate(g%Domain%mpp_domain)
355  deallocate(g%Domain)
356 
357  deallocate(g)
358 
Here is the call graph for this function:
Here is the caller graph for this function:

◆ set_derived_dyn_horgrid()

subroutine, public mom_dyn_horgrid::set_derived_dyn_horgrid ( type(dyn_horgrid_type), intent(inout)  G)

set_derived_dyn_horgrid calculates metric terms that are derived from other metrics.

Parameters
[in,out]gThe dynamic horizontal grid type

Definition at line 259 of file MOM_dyn_horgrid.F90.

References adcroft_reciprocal().

Referenced by mom_transcribe_grid::copy_mom_grid_to_dyngrid(), and mom_grid_initialize::set_grid_metrics().

259  type(dyn_horgrid_type), intent(inout) :: g !< The dynamic horizontal grid type
260 ! Various inverse grid spacings and derived areas are calculated within this
261 ! subroutine.
262  integer :: i, j, isd, ied, jsd, jed
263  integer :: isdb, iedb, jsdb, jedb
264 
265  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
266  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
267 
268  do j=jsd,jed ; do i=isd,ied
269  if (g%dxT(i,j) < 0.0) g%dxT(i,j) = 0.0
270  if (g%dyT(i,j) < 0.0) g%dyT(i,j) = 0.0
271  g%IdxT(i,j) = adcroft_reciprocal(g%dxT(i,j))
272  g%IdyT(i,j) = adcroft_reciprocal(g%dyT(i,j))
273  g%IareaT(i,j) = adcroft_reciprocal(g%areaT(i,j))
274  enddo ; enddo
275 
276  do j=jsd,jed ; do i=isdb,iedb
277  if (g%dxCu(i,j) < 0.0) g%dxCu(i,j) = 0.0
278  if (g%dyCu(i,j) < 0.0) g%dyCu(i,j) = 0.0
279  g%IdxCu(i,j) = adcroft_reciprocal(g%dxCu(i,j))
280  g%IdyCu(i,j) = adcroft_reciprocal(g%dyCu(i,j))
281  enddo ; enddo
282 
283  do j=jsdb,jedb ; do i=isd,ied
284  if (g%dxCv(i,j) < 0.0) g%dxCv(i,j) = 0.0
285  if (g%dyCv(i,j) < 0.0) g%dyCv(i,j) = 0.0
286  g%IdxCv(i,j) = adcroft_reciprocal(g%dxCv(i,j))
287  g%IdyCv(i,j) = adcroft_reciprocal(g%dyCv(i,j))
288  enddo ; enddo
289 
290  do j=jsdb,jedb ; do i=isdb,iedb
291  if (g%dxBu(i,j) < 0.0) g%dxBu(i,j) = 0.0
292  if (g%dyBu(i,j) < 0.0) g%dyBu(i,j) = 0.0
293 
294  g%IdxBu(i,j) = adcroft_reciprocal(g%dxBu(i,j))
295  g%IdyBu(i,j) = adcroft_reciprocal(g%dyBu(i,j))
296  ! areaBu has usually been set to a positive area elsewhere.
297  if (g%areaBu(i,j) <= 0.0) g%areaBu(i,j) = g%dxBu(i,j) * g%dyBu(i,j)
298  g%IareaBu(i,j) = adcroft_reciprocal(g%areaBu(i,j))
299  enddo ; enddo
300 
Here is the call graph for this function:
Here is the caller graph for this function: