ROMS
Loading...
Searching...
No Matches
set_grid.F
Go to the documentation of this file.
1#include "cppdefs.h"
2 SUBROUTINE set_grid (ng, model)
3!
4!git $Id$
5!=======================================================================
6! Copyright (c) 2002-2025 The ROMS Group !
7! Licensed under a MIT/X style license !
8! See License_ROMS.md !
9!================================================== Hernan G. Arango ===
10! !
11! This routine sets application grid and associated variables and !
12! parameters. It called only once during the initialization stage. !
13! !
14!=======================================================================
15!
16 USE mod_param
17 USE mod_parallel
18#ifdef NESTING
19 USE mod_nesting
20#endif
21 USE mod_scalars
22!
24#ifdef DISTRIBUTE
25 USE distribute_mod, ONLY : mp_bcasti
26#endif
27#ifdef TIDE_GENERATING_FORCES
28 USE equilibrium_tide_mod, ONLY : harmonic_constituents
29#endif
30#ifdef GRID_EXTRACT
31 USE get_extract_mod, ONLY : get_extract
32#endif
33#ifndef ANA_GRID
34 USE get_grid_mod, ONLY : get_grid
35#endif
36#ifndef ANA_NUDGCOEF
38#endif
39 USE metrics_mod, ONLY : metrics
40#ifdef NESTING
41 USE nesting_mod, ONLY : nesting
42#endif
43 USE strings_mod, ONLY : founderror
44!
45 implicit none
46!
47! Imported variable declarations.
48!
49 integer, intent(in) :: ng, model
50!
51! Local variable declarations.
52!
53 integer :: tile
54!
55 character (len=*), parameter :: MyFile = &
56 & __FILE__
57!
58!-----------------------------------------------------------------------
59! Set horizontal grid, bathymetry, and Land/Sea masking (if any).
60! Use analytical functions or read in from a grid NetCDF.
61!-----------------------------------------------------------------------
62!
63#ifdef ANA_GRID
64 DO tile=first_tile(ng),last_tile(ng),+1
65 CALL ana_grid (ng, tile, model)
66# ifdef MASKING
67 CALL ana_mask (ng, tile, model)
68# endif
69# if defined AD_SENSITIVITY || defined I4DVAR_ANA_SENSITIVITY || \
70 defined opt_observations || defined so_semi || \
71 defined sensitivity_4dvar
72 CALL ana_scope (ng, tile, model)
73# endif
74 END DO
75!$OMP BARRIER
76#else
77!$OMP MASTER
78# ifdef DISTRIBUTE
79 CALL get_grid (ng, myrank, model)
80# else
81 CALL get_grid (ng, -1, model)
82# endif
83!$OMP END MASTER
84# ifdef DISTRIBUTE
85 CALL mp_bcasti (ng, model, exit_flag)
86# endif
87!$OMP BARRIER
88 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
89#endif
90
91#ifdef GRID_EXTRACT
92!
93!-----------------------------------------------------------------------
94! Get grid geometry for field extraction through interpolation or
95! decimation.
96!-----------------------------------------------------------------------
97!
98# ifdef DISTRIBUTE
99 CALL get_extract (ng, myrank, model)
100# else
101 CALL get_extract (ng, -1, model)
102# endif
103#endif
104
105#ifdef SOLVE3D
106!
107!-----------------------------------------------------------------------
108! Set vertical terrain-following coordinate transformation function.
109!-----------------------------------------------------------------------
110!
111!$OMP MASTER
112 CALL set_scoord (ng)
113!$OMP END MASTER
114!$OMP BARRIER
115#endif
116
117#ifdef SOLVE3D
118!
119!-----------------------------------------------------------------------
120! Set barotropic time-steps average weighting function.
121!-----------------------------------------------------------------------
122!
123!$OMP MASTER
124 CALL set_weights (ng)
125!$OMP END MASTER
126!$OMP BARRIER
127#endif
128
129!
130!-----------------------------------------------------------------------
131! Compute various metric term combinations.
132!-----------------------------------------------------------------------
133!
134 DO tile=first_tile(ng),last_tile(ng),+1
135 CALL metrics (ng, tile, model)
136 END DO
137!$OMP BARRIER
138
139#ifdef NESTING
140!
141!-----------------------------------------------------------------------
142! If nesting, initialize grid spacing (on_u and om_v) in REFINED(:)
143! structure. They are used to impose mass flux at the finer grid
144! physical boundaries.
145!-----------------------------------------------------------------------
146!
147 CALL nesting (ng, model, ndxdy)
148#endif
149
150#if defined WTYPE_GRID && defined ANA_WTYPE && \
151 (defined lmd_skpp || defined solar_source)
152!
153!-----------------------------------------------------------------------
154! Set spatially varying Jerlov water type.
155!-----------------------------------------------------------------------
156!
157 DO tile=first_tile(ng),last_tile(ng),+1
158 CALL ana_wtype (ng, tile, inlm)
159 END DO
160!$OMP BARRIER
161#endif
162
163#ifndef CORRELATION
164!
165!-----------------------------------------------------------------------
166! If appropriate, set spatially varying nudging coefficients time
167! scales.
168!-----------------------------------------------------------------------
169!
170# ifdef ANA_NUDGCOEF
171 IF (lnudging(ng)) THEN
172 DO tile=first_tile(ng),last_tile(ng),+1
173 CALL ana_nudgcoef (ng, tile, model)
174 END DO
175!$OMP BARRIER
176 END IF
177# else
178 IF (lnudging(ng)) THEN
179!$OMP MASTER
180# ifdef DISTRIBUTE
181 CALL get_nudgcoef (ng, myrank, model)
182# else
183 CALL get_nudgcoef (ng, -1, model)
184# endif
185!$OMP END MASTER
186# ifdef DISTRIBUTE
187 CALL mp_bcasti (ng, model, exit_flag)
188# endif
189!$OMP BARRIER
190 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
191 END IF
192# endif
193
194# ifdef TIDE_GENERATING_FORCES
195!
196!-----------------------------------------------------------------------
197! If applying tide generating forces, compute harmonic constituent
198! parameters of the equilibrium tide on Greenwich meridian for the
199! reference tide date number.
200!-----------------------------------------------------------------------
201!
202 CALL harmonic_constituents (lnodal)
203# endif
204#endif
205!
206 RETURN
207 END SUBROUTINE set_grid
subroutine ana_grid(ng, tile, model)
Definition ana_grid.h:3
subroutine lmd_skpp(ng, tile)
Definition lmd_skpp.F:45
subroutine ana_mask(ng, tile, model)
Definition ana_mask.h:3
subroutine ana_wtype(ng, tile, model)
Definition ana_wtype.h:3
subroutine ana_scope(ng, tile, model)
Definition ana_scope.h:3
subroutine ana_nudgcoef(ng, tile, model)
Definition ana_nudgcoef.h:3
subroutine, public get_grid(ng, tile, model)
Definition get_grid.F:55
subroutine, public get_nudgcoef(ng, tile, model)
subroutine, public metrics(ng, tile, model)
Definition metrics.F:24
integer, parameter ndxdy
Definition mod_nesting.F:76
integer, dimension(:), allocatable first_tile
integer, dimension(:), allocatable last_tile
integer, parameter inlm
Definition mod_param.F:662
logical lnodal
logical, dimension(:), allocatable lnudging
integer exit_flag
integer noerror
subroutine, public nesting(ng, model, isection)
Definition nesting.F:140
logical function, public founderror(flag, noerr, line, routine)
Definition strings.F:52
subroutine set_grid(ng, model)
Definition set_grid.F:3
subroutine set_scoord(ng)
Definition set_scoord.F:4
subroutine set_weights(ng)
Definition set_weights.F:4