1 | SUBROUTINE ana_initial (ng, tile, model)
|
---|
2 | !
|
---|
3 | !! svn $Id: ana_initial.h 42 2007-05-04 20:58:11Z arango $
|
---|
4 | !!======================================================================
|
---|
5 | !! Copyright (c) 2002-2007 The ROMS/TOMS Group !
|
---|
6 | !! Licensed under a MIT/X style license !
|
---|
7 | !! See License_ROMS.txt !
|
---|
8 | !! !
|
---|
9 | !=======================================================================
|
---|
10 | ! !
|
---|
11 | ! This subroutine sets initial conditions for momentum and tracer !
|
---|
12 | ! type variables using analytical expressions. !
|
---|
13 | ! !
|
---|
14 | !=======================================================================
|
---|
15 | !
|
---|
16 | USE mod_param
|
---|
17 | USE mod_grid
|
---|
18 | USE mod_ncparam
|
---|
19 | USE mod_ocean
|
---|
20 | USE mod_stepping
|
---|
21 | !
|
---|
22 | ! Imported variable declarations.
|
---|
23 | !
|
---|
24 | integer, intent(in) :: ng, tile, model
|
---|
25 |
|
---|
26 | #include "tile.h"
|
---|
27 | !
|
---|
28 | IF (model.eq.iNLM) THEN
|
---|
29 | CALL ana_NLMinitial_tile (ng, model, Istr, Iend, Jstr, Jend, &
|
---|
30 | & LBi, UBi, LBj, UBj, &
|
---|
31 | & GRID(ng) % h, &
|
---|
32 | #ifdef SPHERICAL
|
---|
33 | & GRID(ng) % lonr, &
|
---|
34 | & GRID(ng) % latr, &
|
---|
35 | #else
|
---|
36 | & GRID(ng) % xr, &
|
---|
37 | & GRID(ng) % yr, &
|
---|
38 | #endif
|
---|
39 | #ifdef SOLVE3D
|
---|
40 | & GRID(ng) % z_r, &
|
---|
41 | & OCEAN(ng) % u, &
|
---|
42 | & OCEAN(ng) % v, &
|
---|
43 | & OCEAN(ng) % t, &
|
---|
44 | #endif
|
---|
45 | & OCEAN(ng) % ubar, &
|
---|
46 | & OCEAN(ng) % vbar, &
|
---|
47 | & OCEAN(ng) % zeta)
|
---|
48 | #ifdef TANGENT
|
---|
49 | ELSE IF ((model.eq.iTLM).or.(model.eq.iRPM)) THEN
|
---|
50 | CALL ana_TLMinitial_tile (ng, model, Istr, Iend, Jstr, Jend, &
|
---|
51 | & LBi, UBi, LBj, UBj, &
|
---|
52 | & kstp(ng), &
|
---|
53 | # ifdef SOLVE3D
|
---|
54 | & nstp(ng), &
|
---|
55 | & OCEAN(ng) % tl_u, &
|
---|
56 | & OCEAN(ng) % tl_v, &
|
---|
57 | & OCEAN(ng) % tl_t, &
|
---|
58 | # endif
|
---|
59 | & OCEAN(ng) % tl_ubar, &
|
---|
60 | & OCEAN(ng) % tl_vbar, &
|
---|
61 | & OCEAN(ng) % tl_zeta)
|
---|
62 | #endif
|
---|
63 | #ifdef ADJOINT
|
---|
64 | ELSE IF (model.eq.iADM) THEN
|
---|
65 | CALL ana_ADMinitial_tile (ng, model, Istr, Iend, Jstr, Jend, &
|
---|
66 | & LBi, UBi, LBj, UBj, &
|
---|
67 | & knew(ng), &
|
---|
68 | # ifdef SOLVE3D
|
---|
69 | & nstp(ng), &
|
---|
70 | & OCEAN(ng) % ad_u, &
|
---|
71 | & OCEAN(ng) % ad_v, &
|
---|
72 | & OCEAN(ng) % ad_t, &
|
---|
73 | # endif
|
---|
74 | & OCEAN(ng) % ad_ubar, &
|
---|
75 | & OCEAN(ng) % ad_vbar, &
|
---|
76 | & OCEAN(ng) % ad_zeta)
|
---|
77 | #endif
|
---|
78 | END IF
|
---|
79 | !
|
---|
80 | ! Set analytical header file name used.
|
---|
81 | !
|
---|
82 | IF (Lanafile) THEN
|
---|
83 | ANANAME(10)='ROMS/Functionals/ana_initial.h'
|
---|
84 | END IF
|
---|
85 |
|
---|
86 | RETURN
|
---|
87 | END SUBROUTINE ana_initial
|
---|
88 | !
|
---|
89 | !***********************************************************************
|
---|
90 | SUBROUTINE ana_NLMinitial_tile (ng, model, Istr, Iend, Jstr, Jend,&
|
---|
91 | & LBi, UBi, LBj, UBj, &
|
---|
92 | & h, &
|
---|
93 | #ifdef SPHERICAL
|
---|
94 | & lonr, latr, &
|
---|
95 | #else
|
---|
96 | & xr, yr, &
|
---|
97 | #endif
|
---|
98 | #ifdef SOLVE3D
|
---|
99 | & z_r, &
|
---|
100 | & u, v, t, &
|
---|
101 | #endif
|
---|
102 | & ubar, vbar, zeta)
|
---|
103 | !***********************************************************************
|
---|
104 | !
|
---|
105 | USE mod_param
|
---|
106 | USE mod_scalars
|
---|
107 | !
|
---|
108 | ! Imported variable declarations.
|
---|
109 | !
|
---|
110 | integer, intent(in) :: ng, model, Iend, Istr, Jend, Jstr
|
---|
111 | integer, intent(in) :: LBi, UBi, LBj, UBj
|
---|
112 | !
|
---|
113 | #ifdef ASSUMED_SHAPE
|
---|
114 | real(r8), intent(in) :: h(LBi:,LBj:)
|
---|
115 | # ifdef SPHERICAL
|
---|
116 | real(r8), intent(in) :: lonr(LBi:,LBj:)
|
---|
117 | real(r8), intent(in) :: latr(LBi:,LBj:)
|
---|
118 | # else
|
---|
119 | real(r8), intent(in) :: xr(LBi:,LBj:)
|
---|
120 | real(r8), intent(in) :: yr(LBi:,LBj:)
|
---|
121 | # endif
|
---|
122 | # ifdef SOLVE3D
|
---|
123 | real(r8), intent(in) :: z_r(LBi:,LBj:,:)
|
---|
124 |
|
---|
125 | real(r8), intent(out) :: u(LBi:,LBj:,:,:)
|
---|
126 | real(r8), intent(out) :: v(LBi:,LBj:,:,:)
|
---|
127 | real(r8), intent(out) :: t(LBi:,LBj:,:,:,:)
|
---|
128 | # endif
|
---|
129 | real(r8), intent(out) :: ubar(LBi:,LBj:,:)
|
---|
130 | real(r8), intent(out) :: vbar(LBi:,LBj:,:)
|
---|
131 | real(r8), intent(out) :: zeta(LBi:,LBj:,:)
|
---|
132 | #else
|
---|
133 | # ifdef SPHERICAL
|
---|
134 | real(r8), intent(in) :: lonr(LBi:UBi,LBj:UBj)
|
---|
135 | real(r8), intent(in) :: latr(LBi:UBi,LBj:UBj)
|
---|
136 | # else
|
---|
137 | real(r8), intent(in) :: xr(LBi:UBi,LBj:UBj)
|
---|
138 | real(r8), intent(in) :: yr(LBi:UBi,LBj:UBj)
|
---|
139 | # endif
|
---|
140 | real(r8), intent(in) :: h(LBi:UBi,LBj:UBj)
|
---|
141 | # ifdef SOLVE3D
|
---|
142 | real(r8), intent(in) :: z_r(LBi:UBi,LBj:UBj,N(ng))
|
---|
143 |
|
---|
144 | real(r8), intent(out) :: u(LBi:UBi,LBj:UBj,N(ng),2)
|
---|
145 | real(r8), intent(out) :: v(LBi:UBi,LBj:UBj,N(ng),2)
|
---|
146 | real(r8), intent(out) :: t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
|
---|
147 | # endif
|
---|
148 | real(r8), intent(out) :: ubar(LBi:UBi,LBj:UBj,3)
|
---|
149 | real(r8), intent(out) :: vbar(LBi:UBi,LBj:UBj,3)
|
---|
150 | real(r8), intent(out) :: zeta(LBi:UBi,LBj:UBj,3)
|
---|
151 | #endif
|
---|
152 | !
|
---|
153 | ! Local variable declarations.
|
---|
154 | !
|
---|
155 | integer :: IstrR, IendR, JstrR, JendR, IstrU, JstrV
|
---|
156 | integer :: Iless, Iplus, i, itrc, j, k
|
---|
157 | real(r8) :: depth, dx, val1, val2, val3, val4, x, x0, y, y0
|
---|
158 |
|
---|
159 | #include "set_bounds.h"
|
---|
160 | !
|
---|
161 | !-----------------------------------------------------------------------
|
---|
162 | ! Initial conditions for 2D momentum (m/s) components.
|
---|
163 | !-----------------------------------------------------------------------
|
---|
164 | !
|
---|
165 | #if defined SOLITON
|
---|
166 | x0=2.0_r8*xl(ng)/3.0_r8
|
---|
167 | y0=0.5_r8*el(ng)
|
---|
168 | val1=0.395_r8
|
---|
169 | val2=0.771_r8*(val1*val1)
|
---|
170 | DO j=JstrR,JendR
|
---|
171 | DO i=Istr,IendR
|
---|
172 | x=0.5_r8*(xr(i-1,j)+xr(i,j))-x0
|
---|
173 | y=0.5_r8*(yr(i-1,j)+yr(i,j))-y0
|
---|
174 | val3=EXP(-val1*x)
|
---|
175 | val4=val2*((2.0_r8*val3/(1.0_r8+(val3*val3)))**2)
|
---|
176 | ubar(i,j,1)=0.25_r8*val4*(6.0_r8*y*y-9.0_r8)* &
|
---|
177 | & EXP(-0.5_r8*y*y)
|
---|
178 | END DO
|
---|
179 | END DO
|
---|
180 | DO j=Jstr,JendR
|
---|
181 | DO i=IstrR,IendR
|
---|
182 | x=0.5_r8*(xr(i,j-1)+xr(i,j))-x0
|
---|
183 | y=0.5_r8*(yr(i,j-1)+yr(i,j))-y0
|
---|
184 | val3=EXP(-val1*x)
|
---|
185 | val4=val2*((2.0_r8*val3/(1.0_r8+(val3*val3)))**2)
|
---|
186 | vbar(i,j,1)=2.0_r8*val4*y*(-2.0_r8*val1*TANH(val1*x))* &
|
---|
187 | & EXP(-0.5_r8*y*y)
|
---|
188 | END DO
|
---|
189 | END DO
|
---|
190 | #elif defined RIVERPLUME2
|
---|
191 | DO j=JstrR,JendR
|
---|
192 | DO i=Istr,IendR
|
---|
193 | ubar(i,j,1)=0.0_r8
|
---|
194 | END DO
|
---|
195 | END DO
|
---|
196 | DO j=Jstr,JendR
|
---|
197 | DO i=IstrR,IendR
|
---|
198 | vbar(i,j,1)=-0.05_r8
|
---|
199 | END DO
|
---|
200 | END DO
|
---|
201 | #elif defined SED_TEST1
|
---|
202 | val1=100.0_r8
|
---|
203 | DO j=JstrR,JendR
|
---|
204 | DO i=Istr,IendR
|
---|
205 | ubar(i,j,1)=-10.0_r8/(10.0_r8+9.0E-06_r8*REAL(i,r8)*val1)
|
---|
206 | END DO
|
---|
207 | END DO
|
---|
208 | DO j=Jstr,JendR
|
---|
209 | DO i=IstrR,IendR
|
---|
210 | vbar(i,j,1)=0.0_r8
|
---|
211 | END DO
|
---|
212 | END DO
|
---|
213 | #elif defined SED_TOY
|
---|
214 | DO j=JstrR,JendR
|
---|
215 | DO i=Istr,IendR
|
---|
216 | ubar(i,j,1)=0.0_r8
|
---|
217 | END DO
|
---|
218 | END DO
|
---|
219 | DO j=Jstr,JendR
|
---|
220 | DO i=IstrR,IendR
|
---|
221 | vbar(i,j,1)=0.0_r8
|
---|
222 | END DO
|
---|
223 | END DO
|
---|
224 | #elif defined TEST_CHAN
|
---|
225 | val1=100.0_r8
|
---|
226 | DO j=JstrR,JendR
|
---|
227 | DO i=Istr,IendR
|
---|
228 | ubar(i,j,1)=0.0_r8
|
---|
229 | END DO
|
---|
230 | END DO
|
---|
231 | DO j=Jstr,JendR
|
---|
232 | DO i=IstrR,IendR
|
---|
233 | vbar(i,j,1)=0.0_r8
|
---|
234 | END DO
|
---|
235 | END DO
|
---|
236 | #else
|
---|
237 | DO j=JstrR,JendR
|
---|
238 | DO i=Istr,IendR
|
---|
239 | ubar(i,j,1)=0.0_r8
|
---|
240 | END DO
|
---|
241 | END DO
|
---|
242 | DO j=Jstr,JendR
|
---|
243 | DO i=IstrR,IendR
|
---|
244 | vbar(i,j,1)=0.0_r8
|
---|
245 | END DO
|
---|
246 | END DO
|
---|
247 | #endif
|
---|
248 | !
|
---|
249 | !-----------------------------------------------------------------------
|
---|
250 | ! Initial conditions for free-surface (m).
|
---|
251 | !-----------------------------------------------------------------------
|
---|
252 | !
|
---|
253 | #if defined KELVIN
|
---|
254 | !! val1=1.0_r8 ! zeta0
|
---|
255 | !! val2=2.0_r8*pi/(12.42_r8*3600.0_r8) ! M2 Tide period
|
---|
256 | DO j=JstrR,JendR
|
---|
257 | DO i=IstrR,IendR
|
---|
258 | !! zeta(i,j,1)=val1* &
|
---|
259 | !! & EXP(-GRID(ng)%f(i,j)*GRID(ng)%yp(i,j)/ &
|
---|
260 | !! & SQRT(g*GRID(ng)%h(i,j)))* &
|
---|
261 | !! & COS(val2*GRID(ng)%xp(i,j)/ &
|
---|
262 | !! & SQRT(g*GRID(ng)%h(i,j)))
|
---|
263 | zeta(i,j,1)=0.0_r8
|
---|
264 | END DO
|
---|
265 | END DO
|
---|
266 | #elif defined SOLITON
|
---|
267 | x0=2.0_r8*xl(ng)/3.0_r8
|
---|
268 | y0=0.5_r8*el(ng)
|
---|
269 | val1=0.395_r8
|
---|
270 | val2=0.771_r8*(val1*val1)
|
---|
271 | DO j=JstrR,JendR
|
---|
272 | DO i=IstrR,IendR
|
---|
273 | x=xr(i,j)-x0
|
---|
274 | y=yr(i,j)-y0
|
---|
275 | val3=EXP(-val1*x)
|
---|
276 | val4=val2*((2.0_r8*val3/(1.0_r8+(val3*val3)))**2)
|
---|
277 | zeta(i,j,1)=0.25_r8*val4*(6.0_r8*y*y+3.0_r8)* &
|
---|
278 | & EXP(-0.5_r8*y*y)
|
---|
279 | END DO
|
---|
280 | END DO
|
---|
281 | #elif defined SED_TEST1
|
---|
282 | val1=100.0_r8
|
---|
283 | DO j=JstrR,JendR
|
---|
284 | DO i=IstrR,IendR
|
---|
285 | zeta(i,j,1)=9.0E-06_r8*REAL(i,r8)*val1
|
---|
286 | END DO
|
---|
287 | END DO
|
---|
288 | #elif defined TEST_CHAN
|
---|
289 | DO j=JstrR,JendR
|
---|
290 | DO i=IstrR,IendR
|
---|
291 | zeta(i,j,1)=-0.4040_r8*REAL(i,r8)/REAL(Lm(ng)+1,r8)
|
---|
292 | END DO
|
---|
293 | END DO
|
---|
294 | #else
|
---|
295 | DO j=JstrR,JendR
|
---|
296 | DO i=IstrR,IendR
|
---|
297 | zeta(i,j,1)=0.0_r8
|
---|
298 | END DO
|
---|
299 | END DO
|
---|
300 | #endif
|
---|
301 | #ifdef SOLVE3D
|
---|
302 | !
|
---|
303 | !-----------------------------------------------------------------------
|
---|
304 | ! Initial conditions for 3D momentum components (m/s).
|
---|
305 | !-----------------------------------------------------------------------
|
---|
306 | !
|
---|
307 | # if defined RIVERPLUME2
|
---|
308 | DO k=1,N(ng)
|
---|
309 | DO j=JstrR,JendR
|
---|
310 | DO i=Istr,IendR
|
---|
311 | u(i,j,k,1)=0.0_r8
|
---|
312 | END DO
|
---|
313 | END DO
|
---|
314 | DO j=Jstr,JendR
|
---|
315 | DO i=IstrR,IendR
|
---|
316 | v(i,j,k,1)=-0.05_r8*LOG((h(i,j)+z_r(i,j,k))/Zob(ng))/ &
|
---|
317 | & (LOG(h(i,j)/Zob(ng))-1.0_r8+Zob(ng)/h(i,j))
|
---|
318 | END DO
|
---|
319 | END DO
|
---|
320 | END DO
|
---|
321 | # elif defined SED_TEST1
|
---|
322 | DO k=1,N(ng)
|
---|
323 | DO j=JstrR,JendR
|
---|
324 | DO i=Istr,IendR
|
---|
325 | u(i,j,k,1)=-1.0_r8*LOG((h(i,j)+z_r(i,j,k))/Zob(ng))/ &
|
---|
326 | & (LOG(h(i,j)/Zob(ng))-1.0_r8+Zob(ng)/h(i,j))
|
---|
327 | END DO
|
---|
328 | END DO
|
---|
329 | DO j=Jstr,JendR
|
---|
330 | DO i=IstrR,IendR
|
---|
331 | v(i,j,k,1)=0.0_r8
|
---|
332 | END DO
|
---|
333 | END DO
|
---|
334 | END DO
|
---|
335 | # elif defined SED_TOY
|
---|
336 | DO k=1,N(ng)
|
---|
337 | DO j=JstrR,JendR
|
---|
338 | DO i=Istr,IendR
|
---|
339 | u(i,j,k,1)=1.0_r8
|
---|
340 | !! u(i,j,k,1)=-1.0_r8
|
---|
341 | !! u(i,j,k,1)=0.0_r8
|
---|
342 | END DO
|
---|
343 | END DO
|
---|
344 | DO j=Jstr,JendR
|
---|
345 | DO i=IstrR,IendR
|
---|
346 | v(i,j,k,1)=0.0_r8
|
---|
347 | END DO
|
---|
348 | END DO
|
---|
349 | END DO
|
---|
350 | # else
|
---|
351 | DO k=1,N(ng)
|
---|
352 | DO j=JstrR,JendR
|
---|
353 | DO i=Istr,IendR
|
---|
354 | u(i,j,k,1)=0.0_r8
|
---|
355 | END DO
|
---|
356 | END DO
|
---|
357 | DO j=Jstr,JendR
|
---|
358 | DO i=IstrR,IendR
|
---|
359 | v(i,j,k,1)=0.0_r8
|
---|
360 | END DO
|
---|
361 | END DO
|
---|
362 | END DO
|
---|
363 | # endif
|
---|
364 | !
|
---|
365 | !-----------------------------------------------------------------------
|
---|
366 | ! Initial conditions for tracer type variables.
|
---|
367 | !-----------------------------------------------------------------------
|
---|
368 | !
|
---|
369 | ! Set initial conditions for potential temperature (Celsius) and
|
---|
370 | ! salinity (PSU).
|
---|
371 | !
|
---|
372 | # if defined BENCHMARK
|
---|
373 | val1=(44.69_r8/39.382_r8)**2
|
---|
374 | val2=val1*(rho0*800.0_r8/g)*(5.0E-05_r8/((42.689_r8/44.69_r8)**2))
|
---|
375 | DO k=1,N(ng)
|
---|
376 | DO j=JstrR,JendR
|
---|
377 | DO i=IstrR,IendR
|
---|
378 | t(i,j,k,1,itemp)=val2*EXP(z_r(i,j,k)/800.0_r8)* &
|
---|
379 | & (0.6_r8-0.4_r8*TANH(z_r(i,j,k)/800.0_r8))
|
---|
380 | t(i,j,k,1,isalt)=35.0_r8
|
---|
381 | END DO
|
---|
382 | END DO
|
---|
383 | END DO
|
---|
384 | # elif defined BASIN
|
---|
385 | val1=(44.69_r8/39.382_r8)**2
|
---|
386 | val2=val1*(rho0*800.0_r8/g)*(5.0E-05_r8/((42.689_r8/44.69_r8)**2))
|
---|
387 | DO k=1,N(ng)
|
---|
388 | DO j=JstrR,JendR
|
---|
389 | DO i=IstrR,IendR
|
---|
390 | t(i,j,k,1,itemp)=val2*EXP(z_r(i,j,k)/800.0_r8)* &
|
---|
391 | & (0.6_r8-0.4_r8*TANH(z_r(i,j,k)/800.0_r8))
|
---|
392 | END DO
|
---|
393 | END DO
|
---|
394 | END DO
|
---|
395 | # elif defined BL_TEST
|
---|
396 | DO k=1,N(ng)
|
---|
397 | DO j=JstrR,JendR
|
---|
398 | DO i=IstrR,IendR
|
---|
399 | val1=TANH(1.1_r8*z_r(i,j,k)+11.0_r8)
|
---|
400 | t(i,j,k,1,itemp)=T0(ng)+6.25_r8*val1
|
---|
401 | t(i,j,k,1,isalt)=S0(ng)-0.75_r8*val1
|
---|
402 | END DO
|
---|
403 | END DO
|
---|
404 | END DO
|
---|
405 | # elif defined CANYON
|
---|
406 | DO k=1,N(ng)
|
---|
407 | DO j=JstrR,JendR
|
---|
408 | DO i=IstrR,IendR
|
---|
409 | t(i,j,k,1,itemp)=3.488_r8*EXP(z_r(i,j,k)/800.0_r8)* &
|
---|
410 | & (1.0_r8-(2.0_r8/3.0_r8)* &
|
---|
411 | & TANH(z_r(i,j,k)/800.0_r8))
|
---|
412 | END DO
|
---|
413 | END DO
|
---|
414 | END DO
|
---|
415 | # elif defined CHANNEL_NECK
|
---|
416 | DO k=1,N(ng)
|
---|
417 | DO j=JstrR,JendR
|
---|
418 | DO i=IstrR,IendR
|
---|
419 | t(i,j,k,1,itemp)=20.0_r8
|
---|
420 | !! t(i,j,k,1,itemp)=14.0_r8+8.0_r8*EXP(z_r(i,j,k)/50.0_r8)
|
---|
421 | END DO
|
---|
422 | END DO
|
---|
423 | END DO
|
---|
424 | # elif defined COUPLING_TEST
|
---|
425 | val1=40.0_r8
|
---|
426 | DO k=1,N(ng)
|
---|
427 | DO j=JstrR,JendR
|
---|
428 | DO i=IstrR,IendR
|
---|
429 | t(i,j,k,1,itemp)=val1*EXP(z_r(i,j,k)/800.0_r8)* &
|
---|
430 | & (0.6_r8-0.4_r8*TANH(z_r(i,j,k)/800.0_r8))+ &
|
---|
431 | & 1.5_r8
|
---|
432 | t(i,j,k,1,isalt)=35.0_r8
|
---|
433 | END DO
|
---|
434 | END DO
|
---|
435 | END DO
|
---|
436 | # elif defined DOUBLE_GYRE
|
---|
437 | val1=(44.69_r8/39.382_r8)**2
|
---|
438 | val2=val1*(rho0*100.0_r8/g)*(5.0E-5_r8/((42.689_r8/44.69_r8)**2))
|
---|
439 | DO k=1,N(ng)
|
---|
440 | DO j=JstrR,JendR
|
---|
441 | DO i=IstrR,IendR
|
---|
442 | val3=T0(ng)+val2*EXP(z_r(i,j,k)/100.0_r8)* &
|
---|
443 | & (10.0_r8-0.4_r8*TANH(z_r(i,j,k)/100.0_r8))
|
---|
444 | val4=yr(i,j)/el(ng)
|
---|
445 | t(i,j,k,1,itemp)=val3-3.0_r8*val4
|
---|
446 | # ifdef SALINITY
|
---|
447 | t(i,j,k,1,isalt)=34.5_r8-0.001_r8*z_r(i,j,k)-val4
|
---|
448 | # endif
|
---|
449 | END DO
|
---|
450 | END DO
|
---|
451 | END DO
|
---|
452 | # elif defined ESTUARY_TEST
|
---|
453 | DO k=1,N(ng)
|
---|
454 | DO j=JstrR,JendR
|
---|
455 | DO i=IstrR,IendR
|
---|
456 | t(i,j,k,1,itemp)=10.0_r8
|
---|
457 | IF (xr(i,j).le.30000.0_r8) then
|
---|
458 | t(i,j,k,1,isalt)=30.0_r8
|
---|
459 | ELSEIF (xr(i,j).le.80000.0_r8) then
|
---|
460 | t(i,j,k,1,isalt)=(80000.0_r8-xr(i,j))/50000.0_r8*30.0_r8
|
---|
461 | ELSE
|
---|
462 | t(i,j,k,1,isalt)=0.0_r8
|
---|
463 | END IF
|
---|
464 | END DO
|
---|
465 | END DO
|
---|
466 | END DO
|
---|
467 | # elif defined FLT_TEST
|
---|
468 | DO k=1,N(ng)
|
---|
469 | DO j=JstrR,JendR
|
---|
470 | DO i=IstrR,IendR
|
---|
471 | t(i,j,k,1,itemp)=T0(ng)
|
---|
472 | END DO
|
---|
473 | END DO
|
---|
474 | END DO
|
---|
475 | # elif defined GRAV_ADJ
|
---|
476 | DO k=1,N(ng)
|
---|
477 | DO j=JstrR,JendR
|
---|
478 | DO i=IstrR,MIN((Lm(ng)+1)/2,IendR)
|
---|
479 | t(i,j,k,1,itemp)=T0(ng)+5.0_r8
|
---|
480 | t(i,j,k,1,isalt)=0.0_r8
|
---|
481 | END DO
|
---|
482 | DO i=MAX((Lm(ng)+1)/2+1,IstrR),IendR
|
---|
483 | t(i,j,k,1,itemp)=T0(ng)
|
---|
484 | t(i,j,k,1,isalt)=S0(ng)
|
---|
485 | END DO
|
---|
486 | !! DO i=IstrR,IendR
|
---|
487 | !! IF (i.lt.Lm(ng)/2) THEN
|
---|
488 | !! t(i,j,k,1,itemp)=T0(ng)+5.0_r8
|
---|
489 | !! ELSE IF (i.eq.Lm(ng)/2) THEN
|
---|
490 | !! t(i,j,k,1,itemp)=T0(ng)+4.0_r8
|
---|
491 | !! ELSE IF (i.eq.Lm(ng)/2+1) THEN
|
---|
492 | !! t(i,j,k,1,itemp)=T0(ng)+1.0_r8
|
---|
493 | !! ELSE
|
---|
494 | !! t(i,j,k,1,itemp)=T0(ng)
|
---|
495 | !! END IF
|
---|
496 | !! END DO
|
---|
497 | END DO
|
---|
498 | END DO
|
---|
499 | # elif defined LAB_CANYON
|
---|
500 | DO k=1,N(ng)
|
---|
501 | DO j=JstrR,JendR
|
---|
502 | DO i=IstrR,IendR
|
---|
503 | t(i,j,k,1,itemp)=-659.34183_r8*z_r(i,j,k)
|
---|
504 | END DO
|
---|
505 | END DO
|
---|
506 | END DO
|
---|
507 | # elif defined LAKE_SIGNELL
|
---|
508 | DO k=1,N(ng)
|
---|
509 | DO j=JstrR,JendR
|
---|
510 | DO i=IstrR,IendR
|
---|
511 | t(i,j,k,1,itemp)=10.0_r8
|
---|
512 | t(i,j,k,1,isalt)=30.0_r8
|
---|
513 | END DO
|
---|
514 | END DO
|
---|
515 | END DO
|
---|
516 | # elif defined LMD_TEST
|
---|
517 | DO k=1,N(ng)
|
---|
518 | DO j=JstrR,JendR
|
---|
519 | DO i=IstrR,IendR
|
---|
520 | t(i,j,k,1,itemp)=MIN(13.0_r8, &
|
---|
521 | & 7.0_r8+0.2_r8*(z_r(i,j,k)+50.0_r8))
|
---|
522 | t(i,j,k,1,isalt)=35.0_r8
|
---|
523 | END DO
|
---|
524 | END DO
|
---|
525 | END DO
|
---|
526 | # elif defined MIXED_LAYER
|
---|
527 | DO k=1,N(ng)
|
---|
528 | DO j=JstrR,JendR
|
---|
529 | DO i=IstrR,IendR
|
---|
530 | t(i,j,k,1,itemp)=10.0_r8+3.0_r8*(z_r(i,j,k)+h(i,j))/ &
|
---|
531 | & h(i,j)
|
---|
532 | t(i,j,k,1,isalt)=S0(ng)
|
---|
533 | END DO
|
---|
534 | END DO
|
---|
535 | END DO
|
---|
536 | # elif defined NJ_BIGHT
|
---|
537 | DO k=1,N(ng)
|
---|
538 | DO j=JstrR,JendR
|
---|
539 | DO i=IstrR,IendR
|
---|
540 | depth=z_r(i,j,k)
|
---|
541 | IF (depth.ge.-15.0_r8) THEN
|
---|
542 | t(i,j,k,1,itemp)= 2.049264257728403E+01_r8-depth* &
|
---|
543 | & (2.640850848793918E-01_r8+depth* &
|
---|
544 | & (2.751125328535212E-01_r8+depth* &
|
---|
545 | & (9.207489761648872E-02_r8+depth* &
|
---|
546 | & (1.449075725742839E-02_r8+depth* &
|
---|
547 | & (1.078215685912076E-03_r8+depth* &
|
---|
548 | & (3.240318053903974E-05_r8+ &
|
---|
549 | & 1.262826857690271E-07_r8*depth))))))
|
---|
550 | t(i,j,k,1,isalt)= 3.066489149193135E+01_r8-depth* &
|
---|
551 | & (1.476725262946735E-01_r8+depth* &
|
---|
552 | & (1.126455760313399E-01_r8+depth* &
|
---|
553 | & (3.900923281871022E-02_r8+depth* &
|
---|
554 | & (6.939014937447098E-03_r8+depth* &
|
---|
555 | & (6.604436696792939E-04_r8+depth* &
|
---|
556 | & (3.191792361954220E-05_r8+ &
|
---|
557 | & 6.177352634409320E-07_r8*depth))))))
|
---|
558 | ELSE
|
---|
559 | t(i,j,k,1,itemp)=14.6_r8+ &
|
---|
560 | & 6.70_r8*TANH(1.1_r8*depth+15.9_r8)
|
---|
561 | t(i,j,k,1,isalt)=31.3_r8- &
|
---|
562 | & 0.55_r8*TANH(1.1_r8*depth+15.9_r8)
|
---|
563 | END IF
|
---|
564 | END DO
|
---|
565 | END DO
|
---|
566 | END DO
|
---|
567 | # elif defined OVERFLOW
|
---|
568 | DO k=1,N(ng)
|
---|
569 | DO j=JstrR,JendR
|
---|
570 | DO i=IstrR,IendR
|
---|
571 | t(i,j,k,1,itemp)=T0(ng)-0.5_r8*T0(ng)*(1.0_r8+ &
|
---|
572 | & TANH((yr(i,j)-60000.0_r8)/2000.0_r8))
|
---|
573 | END DO
|
---|
574 | END DO
|
---|
575 | END DO
|
---|
576 | # elif defined RIVERPLUME1
|
---|
577 | DO k=1,N(ng)
|
---|
578 | DO j=JstrR,JendR
|
---|
579 | DO i=IstrR,IendR
|
---|
580 | t(i,j,k,1,itemp)=T0(ng)+0.01_r8*REAL(k,r8)
|
---|
581 | t(i,j,k,1,isalt)=S0(ng)
|
---|
582 | END DO
|
---|
583 | END DO
|
---|
584 | END DO
|
---|
585 | # elif defined RIVERPLUME2
|
---|
586 | DO k=1,N(ng)
|
---|
587 | DO j=JstrR,JendR
|
---|
588 | DO i=IstrR,IendR
|
---|
589 | t(i,j,k,1,itemp)=T0(ng)
|
---|
590 | t(i,j,k,1,isalt)=S0(ng)
|
---|
591 | END DO
|
---|
592 | END DO
|
---|
593 | END DO
|
---|
594 | # elif defined SEAMOUNT
|
---|
595 | DO k=1,N(ng)
|
---|
596 | DO j=JstrR,JendR
|
---|
597 | DO i=IstrR,IendR
|
---|
598 | t(i,j,k,1,itemp)=T0(ng)+7.5_r8*EXP(z_r(i,j,k)/1000.0_r8)
|
---|
599 | END DO
|
---|
600 | END DO
|
---|
601 | END DO
|
---|
602 | # elif defined SED_TEST1
|
---|
603 | DO k=1,N(ng)
|
---|
604 | DO j=JstrR,JendR
|
---|
605 | DO i=IstrR,IendR
|
---|
606 | t(i,j,k,1,itemp)=20.0_r8
|
---|
607 | t(i,j,k,1,isalt)=0.0_r8
|
---|
608 | END DO
|
---|
609 | END DO
|
---|
610 | END DO
|
---|
611 | # elif defined UPWELLING
|
---|
612 | DO k=1,N(ng)
|
---|
613 | DO j=JstrR,JendR
|
---|
614 | DO i=IstrR,IendR
|
---|
615 | t(i,j,k,1,itemp)=T0(ng)+8.0_r8*EXP(z_r(i,j,k)/50.0_r8)
|
---|
616 | !! t(i,j,k,1,itemp)=T0(ng)+(z_r(i,j,k)+75.0_r8)/150.0_r8+
|
---|
617 | !! & 4.0_r8*(1.0_r8+TANH((z_r(i,j,k)+35.0_r8)/
|
---|
618 | !! & 6.5_r8))
|
---|
619 | !! t(i,j,k,1,isalt)=1.0E-04_r8*yr(i,j)-S0(ng)
|
---|
620 | t(i,j,k,1,isalt)=S0(ng)
|
---|
621 | !! IF (j.lt.Mm(ng)/2) THEN
|
---|
622 | !! t(i,j,k,1,isalt)=0.0_r8
|
---|
623 | !! ELSE IF (j.eq.Mm(ng)/2) THEN
|
---|
624 | !! t(i,j,k,1,isalt)=0.5_r8
|
---|
625 | !! ELSE IF (j.gt.Mm(ng)/2) THEN
|
---|
626 | !! t(i,j,k,1,isalt)=1.0_r8
|
---|
627 | !! END IF
|
---|
628 | END DO
|
---|
629 | END DO
|
---|
630 | END DO
|
---|
631 | # elif defined WINDBASIN
|
---|
632 | DO k=1,N(ng)
|
---|
633 | DO j=JstrR,JendR
|
---|
634 | DO i=IstrR,IendR
|
---|
635 | t(i,j,k,1,itemp)=20.0_r8 ! homogeneous
|
---|
636 | !! t(i,j,k,1,itemp)=14.0_r8+8.0_r8*EXP(z_r(i,j,k)/50.0_r8)- &
|
---|
637 | !! & T0(ng) ! stratified
|
---|
638 | END DO
|
---|
639 | END DO
|
---|
640 | END DO
|
---|
641 | # else
|
---|
642 | DO k=1,N(ng)
|
---|
643 | DO j=JstrR,JendR
|
---|
644 | DO i=IstrR,IendR
|
---|
645 | t(i,j,k,1,itemp)=T0(ng)
|
---|
646 | # ifdef SALINITY
|
---|
647 | t(i,j,k,1,isalt)=S0(ng)
|
---|
648 | # endif
|
---|
649 | END DO
|
---|
650 | END DO
|
---|
651 | END DO
|
---|
652 | # endif
|
---|
653 | #endif
|
---|
654 | RETURN
|
---|
655 | END SUBROUTINE ana_NLMinitial_tile
|
---|
656 |
|
---|
657 | #ifdef TANGENT
|
---|
658 | !
|
---|
659 | !***********************************************************************
|
---|
660 | SUBROUTINE ana_TLMinitial_tile (ng, model, Istr, Iend, Jstr, Jend,&
|
---|
661 | & LBi, UBi, LBj, UBj, &
|
---|
662 | & kstp, &
|
---|
663 | # ifdef SOLVE3D
|
---|
664 | & nstp, &
|
---|
665 | & tl_u, tl_v, tl_t, &
|
---|
666 | # endif
|
---|
667 | & tl_ubar, tl_vbar, tl_zeta)
|
---|
668 | !***********************************************************************
|
---|
669 | !
|
---|
670 | USE mod_param
|
---|
671 | USE mod_scalars
|
---|
672 | !
|
---|
673 | ! Imported variable declarations.
|
---|
674 | !
|
---|
675 | integer, intent(in) :: ng, model, Iend, Istr, Jend, Jstr
|
---|
676 | integer, intent(in) :: LBi, UBi, LBj, UBj
|
---|
677 | integer, intent(in) :: kstp
|
---|
678 | # ifdef SOLVE3D
|
---|
679 | integer, intent(in) :: nstp
|
---|
680 | # endif
|
---|
681 | !
|
---|
682 | # ifdef ASSUMED_SHAPE
|
---|
683 | # ifdef SOLVE3D
|
---|
684 | real(r8), intent(out) :: tl_u(LBi:,LBj:,:,:)
|
---|
685 | real(r8), intent(out) :: tl_v(LBi:,LBj:,:,:)
|
---|
686 | real(r8), intent(out) :: tl_t(LBi:,LBj:,:,:,:)
|
---|
687 | # endif
|
---|
688 | real(r8), intent(out) :: tl_ubar(LBi:,LBj:,:)
|
---|
689 | real(r8), intent(out) :: tl_vbar(LBi:,LBj:,:)
|
---|
690 | real(r8), intent(out) :: tl_zeta(LBi:,LBj:,:)
|
---|
691 | # else
|
---|
692 | # ifdef SOLVE3D
|
---|
693 | real(r8), intent(out) :: tl_u(LBi:UBi,LBj:UBj,N(ng),2)
|
---|
694 | real(r8), intent(out) :: tl_v(LBi:UBi,LBj:UBj,N(ng),2)
|
---|
695 | real(r8), intent(out) :: tl_t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
|
---|
696 | # endif
|
---|
697 | real(r8), intent(out) :: tl_ubar(LBi:UBi,LBj:UBj,3)
|
---|
698 | real(r8), intent(out) :: tl_vbar(LBi:UBi,LBj:UBj,3)
|
---|
699 | real(r8), intent(out) :: tl_zeta(LBi:UBi,LBj:UBj,3)
|
---|
700 | # endif
|
---|
701 | !
|
---|
702 | ! Local variable declarations.
|
---|
703 | !
|
---|
704 | integer :: IstrR, IendR, JstrR, JendR, IstrU, JstrV
|
---|
705 | integer :: i, itrc, j, k
|
---|
706 |
|
---|
707 | # include "set_bounds.h"
|
---|
708 | !
|
---|
709 | !-----------------------------------------------------------------------
|
---|
710 | ! Initial conditions for tangent linear 2D momentum (s/m) components.
|
---|
711 | !-----------------------------------------------------------------------
|
---|
712 | !
|
---|
713 | # if defined MY_APPLICATION
|
---|
714 | # else
|
---|
715 | DO j=JstrR,JendR
|
---|
716 | DO i=Istr,IendR
|
---|
717 | tl_ubar(i,j,kstp)=0.0_r8
|
---|
718 | END DO
|
---|
719 | END DO
|
---|
720 | DO j=Jstr,JendR
|
---|
721 | DO i=IstrR,IendR
|
---|
722 | tl_vbar(i,j,kstp)=0.0_r8
|
---|
723 | END DO
|
---|
724 | END DO
|
---|
725 | # endif
|
---|
726 | !
|
---|
727 | !-----------------------------------------------------------------------
|
---|
728 | ! Initial conditions for tangent linear free-surface (1/m).
|
---|
729 | !-----------------------------------------------------------------------
|
---|
730 | !
|
---|
731 | # if defined MY_APPLICATION
|
---|
732 | # else
|
---|
733 | DO j=JstrR,JendR
|
---|
734 | DO i=IstrR,IendR
|
---|
735 | tl_zeta(i,j,kstp)=0.0_r8
|
---|
736 | END DO
|
---|
737 | END DO
|
---|
738 | # endif
|
---|
739 | # ifdef SOLVE3D
|
---|
740 | !
|
---|
741 | !-----------------------------------------------------------------------
|
---|
742 | ! Initial conditions for tangent linear 3D momentum components (s/m).
|
---|
743 | !-----------------------------------------------------------------------
|
---|
744 | !
|
---|
745 | # if defined MY_APPLICATION
|
---|
746 | # else
|
---|
747 | DO k=1,N(ng)
|
---|
748 | DO j=JstrR,JendR
|
---|
749 | DO i=Istr,IendR
|
---|
750 | tl_u(i,j,k,nstp)=0.0_r8
|
---|
751 | END DO
|
---|
752 | END DO
|
---|
753 | DO j=Jstr,JendR
|
---|
754 | DO i=IstrR,IendR
|
---|
755 | tl_v(i,j,k,nstp)=0.0_r8
|
---|
756 | END DO
|
---|
757 | END DO
|
---|
758 | END DO
|
---|
759 | # endif
|
---|
760 | !
|
---|
761 | !-----------------------------------------------------------------------
|
---|
762 | ! Initial conditions for tangent linear active tracers (1/Tunits).
|
---|
763 | !-----------------------------------------------------------------------
|
---|
764 | !
|
---|
765 | # if defined MY_APLICATION
|
---|
766 | # else
|
---|
767 | DO itrc=1,NAT
|
---|
768 | DO k=1,N(ng)
|
---|
769 | DO j=JstrR,JendR
|
---|
770 | DO i=IstrR,IendR
|
---|
771 | tl_t(i,j,k,nstp,itrc)=0.0_r8
|
---|
772 | END DO
|
---|
773 | END DO
|
---|
774 | END DO
|
---|
775 | END DO
|
---|
776 | # endif
|
---|
777 | # endif
|
---|
778 | RETURN
|
---|
779 | END SUBROUTINE ana_TLMinitial_tile
|
---|
780 | #endif
|
---|
781 |
|
---|
782 | #ifdef ADJOINT
|
---|
783 | !
|
---|
784 | !***********************************************************************
|
---|
785 | SUBROUTINE ana_ADMinitial_tile (ng, model, Istr, Iend, Jstr, Jend,&
|
---|
786 | & LBi, UBi, LBj, UBj, &
|
---|
787 | & knew, &
|
---|
788 | # ifdef SOLVE3D
|
---|
789 | & nstp, &
|
---|
790 | & ad_u, ad_v, ad_t, &
|
---|
791 | # endif
|
---|
792 | & ad_ubar, ad_vbar, ad_zeta)
|
---|
793 | !***********************************************************************
|
---|
794 | !
|
---|
795 | USE mod_param
|
---|
796 | USE mod_scalars
|
---|
797 | !
|
---|
798 | ! Imported variable declarations.
|
---|
799 | !
|
---|
800 | integer, intent(in) :: ng, model, Iend, Istr, Jend, Jstr
|
---|
801 | integer, intent(in) :: LBi, UBi, LBj, UBj
|
---|
802 | integer, intent(in) :: knew
|
---|
803 | # ifdef SOLVE3D
|
---|
804 | integer, intent(in) :: nstp
|
---|
805 | # endif
|
---|
806 | !
|
---|
807 | # ifdef ASSUMED_SHAPE
|
---|
808 | # ifdef SOLVE3D
|
---|
809 | real(r8), intent(out) :: ad_u(LBi:,LBj:,:,:)
|
---|
810 | real(r8), intent(out) :: ad_v(LBi:,LBj:,:,:)
|
---|
811 | real(r8), intent(out) :: ad_t(LBi:,LBj:,:,:,:)
|
---|
812 | # endif
|
---|
813 | real(r8), intent(out) :: ad_ubar(LBi:,LBj:,:)
|
---|
814 | real(r8), intent(out) :: ad_vbar(LBi:,LBj:,:)
|
---|
815 | real(r8), intent(out) :: ad_zeta(LBi:,LBj:,:)
|
---|
816 | # else
|
---|
817 | # ifdef SOLVE3D
|
---|
818 | real(r8), intent(out) :: ad_u(LBi:UBi,LBj:UBj,N(ng),2)
|
---|
819 | real(r8), intent(out) :: ad_v(LBi:UBi,LBj:UBj,N(ng),2)
|
---|
820 | real(r8), intent(out) :: ad_t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
|
---|
821 | # endif
|
---|
822 | real(r8), intent(out) :: ad_ubar(LBi:UBi,LBj:UBj,3)
|
---|
823 | real(r8), intent(out) :: ad_vbar(LBi:UBi,LBj:UBj,3)
|
---|
824 | real(r8), intent(out) :: ad_zeta(LBi:UBi,LBj:UBj,3)
|
---|
825 | # endif
|
---|
826 | !
|
---|
827 | ! Local variable declarations.
|
---|
828 | !
|
---|
829 | integer :: IstrR, IendR, JstrR, JendR, IstrU, JstrV
|
---|
830 | integer :: i, itrc, j, k
|
---|
831 |
|
---|
832 | # include "set_bounds.h"
|
---|
833 | !
|
---|
834 | !-----------------------------------------------------------------------
|
---|
835 | ! Initial conditions for adjoint 2D momentum (s/m) components.
|
---|
836 | !-----------------------------------------------------------------------
|
---|
837 | !
|
---|
838 | # if defined MY_APPLICATION
|
---|
839 | # else
|
---|
840 | DO j=JstrR,JendR
|
---|
841 | DO i=Istr,IendR
|
---|
842 | ad_ubar(i,j,knew)=0.0_r8
|
---|
843 | END DO
|
---|
844 | END DO
|
---|
845 | DO j=Jstr,JendR
|
---|
846 | DO i=IstrR,IendR
|
---|
847 | ad_vbar(i,j,knew)=0.0_r8
|
---|
848 | END DO
|
---|
849 | END DO
|
---|
850 | # endif
|
---|
851 | !
|
---|
852 | !-----------------------------------------------------------------------
|
---|
853 | ! Initial conditions for adjoint free-surface (1/m).
|
---|
854 | !-----------------------------------------------------------------------
|
---|
855 | !
|
---|
856 | # if defined MY_APPLICATION
|
---|
857 | # else
|
---|
858 | DO j=JstrR,JendR
|
---|
859 | DO i=IstrR,IendR
|
---|
860 | ad_zeta(i,j,knew)=0.0_r8
|
---|
861 | END DO
|
---|
862 | END DO
|
---|
863 | # endif
|
---|
864 | # ifdef SOLVE3D
|
---|
865 | !
|
---|
866 | !-----------------------------------------------------------------------
|
---|
867 | ! Initial conditions for adjoint 3D momentum components (s/m).
|
---|
868 | !-----------------------------------------------------------------------
|
---|
869 | !
|
---|
870 | # if defined MY_APPLICATION
|
---|
871 | # else
|
---|
872 | DO k=1,N(ng)
|
---|
873 | DO j=JstrR,JendR
|
---|
874 | DO i=Istr,IendR
|
---|
875 | ad_u(i,j,k,nstp)=0.0_r8
|
---|
876 | END DO
|
---|
877 | END DO
|
---|
878 | DO j=Jstr,JendR
|
---|
879 | DO i=IstrR,IendR
|
---|
880 | ad_v(i,j,k,nstp)=0.0_r8
|
---|
881 | END DO
|
---|
882 | END DO
|
---|
883 | END DO
|
---|
884 | # endif
|
---|
885 | !
|
---|
886 | !-----------------------------------------------------------------------
|
---|
887 | ! Initial conditions for adjoint active tracers (1/Tunits).
|
---|
888 | !-----------------------------------------------------------------------
|
---|
889 | !
|
---|
890 | # if defined MY_APLICATION
|
---|
891 | # else
|
---|
892 | DO itrc=1,NAT
|
---|
893 | DO k=1,N(ng)
|
---|
894 | DO j=JstrR,JendR
|
---|
895 | DO i=IstrR,IendR
|
---|
896 | ad_t(i,j,k,nstp,itrc)=0.0_r8
|
---|
897 | END DO
|
---|
898 | END DO
|
---|
899 | END DO
|
---|
900 | END DO
|
---|
901 | # endif
|
---|
902 | # endif
|
---|
903 | RETURN
|
---|
904 | END SUBROUTINE ana_ADMinitial_tile
|
---|
905 | #endif
|
---|