ROMS
Loading...
Searching...
No Matches
ana_nudgcoef.h
Go to the documentation of this file.
1!!
2 SUBROUTINE ana_nudgcoef (ng, tile, model)
3!
4!! git $Id$
5!!================================================= Hernan G. Arango ===
6!! Copyright (c) 2002-2025 The ROMS Group !
7!! Licensed under a MIT/X style license !
8!! See License_ROMS.md !
9!=======================================================================
10! !
11! This routine sets spatially varying nudging coefficients time- !
12! scales (1/s). They are used for nuding to climatology in the !
13! governing equations. !
14! !
15! It is HIGHLY recommended to write these nudging coefficients into !
16! input NetCDF NUDNAME instead of using analytical expressions !
17! below. It is very easy to introduce parallel bugs. Also, Users !
18! can plot their spatial distribution and fine tune their values !
19! during the pre-proccessing stage for a particular application. !
20! !
21! REMARK: Nudging of free-surface in the vertically integrated !
22! ====== continuity equation is NOT allowed because it VIOLATES !
23! mass/volume conservation. If such nudging effects are required, !
24! it needs to be specified on the momentum equations for (u,v) !
25! and/or (ubar,vbar). If done on (u,v) only, its effects enter !
26! the 2D momentum equations via the residual vertically integrated !
27! forcing term. !
28! !
29!=======================================================================
30!
31 USE mod_param
32 USE mod_ncparam
33!
34! Imported variable declarations.
35!
36 integer, intent(in) :: ng, tile, model
37!
38! Local variable declarations.
39!
40 character (len=*), parameter :: MyFile = &
41 & __FILE__
42!
43#include "tile.h"
44!
45 CALL ana_nudgcoef_tile (ng, tile, model, &
46 & lbi, ubi, lbj, ubj, &
47 & imins, imaxs, jmins, jmaxs)
48!
49! Set analytical header file name used.
50!
51#ifdef DISTRIBUTE
52 IF (lanafile) THEN
53#else
54 IF (lanafile.and.(tile.eq.0)) THEN
55#endif
56 ananame(16)=myfile
57 END IF
58!
59 RETURN
60 END SUBROUTINE ana_nudgcoef
61!
62!***********************************************************************
63 SUBROUTINE ana_nudgcoef_tile (ng, tile, model, &
64 & LBi, UBi, LBj, UBj, &
65 & IminS, ImaxS, JminS, JmaxS)
66!***********************************************************************
67!
68 USE mod_param
69 USE mod_parallel
70 USE mod_clima
71 USE mod_grid
72 USE mod_ncparam
73 USE mod_scalars
74#ifdef DISTRIBUTE
75!
76 USE distribute_mod, ONLY : mp_collect
78# ifdef SOLVE3D
81# endif
82#endif
83!
84 implicit none
85!
86! Imported variable declarations.
87!
88 integer, intent(in) :: ng, tile, model
89 integer, intent(in) :: LBi, UBi, LBj, UBj
90 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
91!
92! Local variable declarations.
93!
94 integer :: Iwrk, i, itrc, j, k
95!
96 real(r8) :: cff1, cff2, cff3
97!
98 real(r8), parameter :: IniVal = 0.0_r8
99!
100 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: wrk
101
102#include "set_bounds.h"
103!
104!-----------------------------------------------------------------------
105! Set up nudging towards data time-scale coefficients (1/s).
106!-----------------------------------------------------------------------
107!
108! Initialize.
109!
110 DO j=jstrt,jendt
111 DO i=istrt,iendt
112 wrk(i,j)=0.0_r8
113 END DO
114 END DO
115
116#if defined ADRIA02
117!
118! Set tracer nudging coefficients in the southern edges from a three
119! days time scale at the boundary point to decrease linearly to 30 days
120! six grids points away from the boundary.
121!
122 cff1=1.0_r8/(3.0_r8*86400.0_r8)
123 cff2=1.0_r8/(30.0_r8*86400.0_r8)
124 DO j=jstrt,min(6,jendt)
125 DO i=istrt,iendt
126 wrk(i,j)=cff2+real(6-j,r8)*(cff1-cff2)/6.0_r8
127 END DO
128 END DO
129
130 IF (any(lnudgetclm(:,ng))) THEN
131 DO itrc=1,ntclm(ng)
132 DO k=1,n(ng)
133 DO j=jstrt,jendt
134 DO i=istrt,iendt
135 clima(ng)%Tnudgcof(i,j,k,itrc)=wrk(i,j)
136 END DO
137 END DO
138 END DO
139 END DO
140 END IF
141
142#elif defined DAMEE_4
143!
144! Set tracer nudging coefficients in the southern and northern edges
145! from a five days time scale at the boundary point to decrease
146! linearly to 60 days seven grids points away from the boundary.
147!
148 cff1=1.0_r8/(5.0_r8*86400.0_r8)
149 cff2=1.0_r8/(60.0_r8*86400.0_r8)
150 cff3=(7.0_r8*cff1-cff2)/6.0_r8
151 DO j=jstrt,min(8,jendt)
152 DO i=istrt,iendt
153 wrk(i,j)=cff2+real(8-j,r8)*(cff1-cff2)/7.0_r8
154 END DO
155 END DO
156
157 DO j=max(jstrt,mm(ng)-7),jendt
158 DO i=istrt,iendt
159 wrk(i,j)=cff1+real(mm(ng)-j,r8)*(cff2-cff1)/7.0_r8
160 END DO
161 END DO
162
163 DO j=max(jstrt,74),min(80,jendt)
164 DO i=max(istrt,102),min(108,iendt)
165 cff1=sqrt(real((i-109)*(i-109)+(j-77)*(j-77),r8))
166 wrk(i,j)=max(0.0_r8,(cff3+cff1*(cff2-cff3)/6.0_r8))
167 END DO
168 END DO
169
170 IF (any(lnudgetclm(:,ng))) THEN
171 DO itrc=1,ntclm(ng)
172 DO k=1,n(ng)
173 DO j=jstrt,jendt
174 DO i=istrt,iendt
175 clima(ng)%Tnudgcof(i,j,k,itrc)=wrk(i,j)
176 END DO
177 END DO
178 END DO
179 END DO
180 END IF
181
182#else
183!
184! Default nudging coefficients. Set nudging coefficients uniformly to
185! the values specified in the standard input file.
186!
187 IF (lnudgem2clm(ng)) THEN
188 DO j=jstrt,jendt
189 DO i=istrt,iendt
190 clima(ng)%M2nudgcof(i,j)=m2nudg(ng)
191 END DO
192 END DO
193 END IF
194
195# ifdef SOLVE3D
196!
197 IF (lnudgem3clm(ng)) THEN
198 DO k=1,n(ng)
199 DO j=jstrt,jendt
200 DO i=istrt,iendt
201 clima(ng)%M3nudgcof(i,j,k)=m3nudg(ng)
202 END DO
203 END DO
204 END DO
205 END IF
206!
207 IF (any(lnudgetclm(:,ng))) THEN
208 DO itrc=1,ntclm(ng)
209 DO k=1,n(ng)
210 DO j=jstrt,jendt
211 DO i=istrt,iendt
212 clima(ng)%Tnudgcof(i,j,k,itrc)=tnudg(itrc,ng)
213 END DO
214 END DO
215 END DO
216 END DO
217 END IF
218# endif
219#endif
220#ifdef DISTRIBUTE
221!
222!-----------------------------------------------------------------------
223! Exchage nudging coefficients information.
224!-----------------------------------------------------------------------
225!
226 IF (lnudgem2clm(ng)) THEN
227 CALL mp_exchange2d (ng, tile, model, 1, &
228 & lbi, ubi, lbj, ubj, &
229 & nghostpoints, .false., .false., &
230 & clima(ng)%M2nudgcof)
231 END IF
232
233# ifdef SOLVE3D
234!
235 IF (lnudgem3clm(ng)) THEN
236 CALL mp_exchange3d (ng, tile, model, 1, &
237 & lbi, ubi, lbj, ubj, 1, n(ng), &
238 & nghostpoints, .false., .false., &
239 & clima(ng)%M3nudgcof)
240 END IF
241!
242 IF (any(lnudgetclm(:,ng))) THEN
243 CALL mp_exchange4d (ng, tile, model, 1, &
244 & lbi, ubi, lbj, ubj, 1, n(ng), 1, ntclm(ng), &
245 & nghostpoints, .false., .false., &
246 & clima(ng)%Tnudgcof)
247 END IF
248# endif
249#endif
250!
251 RETURN
252 END SUBROUTINE ana_nudgcoef_tile
subroutine ana_nudgcoef_tile(ng, tile, model, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs)
subroutine ana_nudgcoef(ng, tile, model)
Definition ana_nudgcoef.h:3
type(t_clima), dimension(:), allocatable clima
Definition mod_clima.F:153
logical lanafile
character(len=256), dimension(39) ananame
integer, dimension(:), allocatable ntclm
Definition mod_param.F:494
integer, dimension(:), allocatable n
Definition mod_param.F:479
integer nghostpoints
Definition mod_param.F:710
integer, dimension(:), allocatable mm
Definition mod_param.F:456
logical, dimension(:), allocatable lnudgem2clm
real(dp), dimension(:), allocatable m2nudg
real(dp), dimension(:,:), allocatable tnudg
logical, dimension(:), allocatable lnudgem3clm
real(dp), dimension(:), allocatable m3nudg
logical, dimension(:,:), allocatable lnudgetclm
subroutine mp_exchange4d(ng, tile, model, nvar, lbi, ubi, lbj, ubj, lbk, ubk, lbt, ubt, nghost, ew_periodic, ns_periodic, a, b, c)
subroutine mp_exchange2d(ng, tile, model, nvar, lbi, ubi, lbj, ubj, nghost, ew_periodic, ns_periodic, a, b, c, d)
subroutine mp_exchange3d(ng, tile, model, nvar, lbi, ubi, lbj, ubj, lbk, ubk, nghost, ew_periodic, ns_periodic, a, b, c, d)