1 | #include "cppdefs.h"
|
---|
2 | MODULE set_vbc_mod
|
---|
3 | #ifdef NONLINEAR
|
---|
4 | !
|
---|
5 | !svn $Id: set_vbc.F 8 2007-02-06 19:00:29Z arango $
|
---|
6 | !================================================== Hernan G. Arango ===
|
---|
7 | ! Copyright (c) 2002-2007 The ROMS/TOMS Group !
|
---|
8 | ! Licensed under a MIT/X style license !
|
---|
9 | ! See License_ROMS.txt !
|
---|
10 | !=======================================================================
|
---|
11 | ! !
|
---|
12 | ! This module sets vertical boundary conditons for momentum and !
|
---|
13 | ! tracers. !
|
---|
14 | ! !
|
---|
15 | !=======================================================================
|
---|
16 | !
|
---|
17 | implicit none
|
---|
18 |
|
---|
19 | PRIVATE
|
---|
20 | PUBLIC :: set_vbc
|
---|
21 |
|
---|
22 | CONTAINS
|
---|
23 |
|
---|
24 | # ifdef SOLVE3D
|
---|
25 | !
|
---|
26 | !***********************************************************************
|
---|
27 | SUBROUTINE set_vbc (ng, tile)
|
---|
28 | !***********************************************************************
|
---|
29 | !
|
---|
30 | USE mod_param
|
---|
31 | USE mod_grid
|
---|
32 | USE mod_forces
|
---|
33 | USE mod_ocean
|
---|
34 | USE mod_stepping
|
---|
35 | !
|
---|
36 | ! Imported variable declarations.
|
---|
37 | !
|
---|
38 | integer, intent(in) :: ng, tile
|
---|
39 | !
|
---|
40 | ! Local variable declarations.
|
---|
41 | !
|
---|
42 | # include "tile.h"
|
---|
43 | !
|
---|
44 | # ifdef PROFILE
|
---|
45 | CALL wclock_on (ng, iNLM, 6)
|
---|
46 | # endif
|
---|
47 | CALL set_vbc_tile (ng, Istr, Iend, Jstr, Jend, &
|
---|
48 | & LBi, UBi, LBj, UBj, &
|
---|
49 | & nrhs(ng), &
|
---|
50 | & GRID(ng) % Hz, &
|
---|
51 | # if !(defined BBL_MODEL || defined ICESHELF)
|
---|
52 | & GRID(ng) % z_r, &
|
---|
53 | & GRID(ng) % z_w, &
|
---|
54 | # endif
|
---|
55 | # if defined ICESHELF
|
---|
56 | & GRID(ng) % zice, &
|
---|
57 | # endif
|
---|
58 | & OCEAN(ng) % t, &
|
---|
59 | # if !(defined BBL_MODEL || defined ICESHELF)
|
---|
60 | & OCEAN(ng) % u, &
|
---|
61 | & OCEAN(ng) % v, &
|
---|
62 | # endif
|
---|
63 | # ifdef QCORRECTION
|
---|
64 | & FORCES(ng) % dqdt, &
|
---|
65 | & FORCES(ng) % sst, &
|
---|
66 | # endif
|
---|
67 | # if defined SCORRECTION || defined SRELAXATION
|
---|
68 | & FORCES(ng) % sss, &
|
---|
69 | # endif
|
---|
70 | # if defined ICESHELF
|
---|
71 | # ifdef SHORTWAVE
|
---|
72 | & FORCES(ng) % srflx, &
|
---|
73 | # endif
|
---|
74 | & FORCES(ng) % sustr, &
|
---|
75 | & FORCES(ng) % svstr, &
|
---|
76 | # endif
|
---|
77 | # ifndef BBL_MODEL
|
---|
78 | & FORCES(ng) % bustr, &
|
---|
79 | & FORCES(ng) % bvstr, &
|
---|
80 | # endif
|
---|
81 | & FORCES(ng) % stflx, &
|
---|
82 | & FORCES(ng) % btflx)
|
---|
83 | # ifdef PROFILE
|
---|
84 | CALL wclock_off (ng, iNLM, 6)
|
---|
85 | # endif
|
---|
86 | RETURN
|
---|
87 | END SUBROUTINE set_vbc
|
---|
88 | !
|
---|
89 | !***********************************************************************
|
---|
90 | SUBROUTINE set_vbc_tile (ng, Istr, Iend, Jstr, Jend, &
|
---|
91 | & LBi, UBi, LBj, UBj, &
|
---|
92 | & nrhs, &
|
---|
93 | & Hz, &
|
---|
94 | # if !(defined BBL_MODEL || defined ICESHELF)
|
---|
95 | & z_r, z_w, &
|
---|
96 | # endif
|
---|
97 | # if defined ICESHELF
|
---|
98 | & zice, &
|
---|
99 | # endif
|
---|
100 | & t, &
|
---|
101 | # if !(defined BBL_MODEL || defined ICESHELF)
|
---|
102 | & u, v, &
|
---|
103 | # endif
|
---|
104 | # ifdef QCORRECTION
|
---|
105 | & dqdt, sst, &
|
---|
106 | # endif
|
---|
107 | # if defined SCORRECTION || defined SRELAXATION
|
---|
108 | & sss, &
|
---|
109 | # endif
|
---|
110 | # if defined ICESHELF
|
---|
111 | # ifdef SHORTWAVE
|
---|
112 | & srflx, &
|
---|
113 | # endif
|
---|
114 | & sustr, svstr, &
|
---|
115 | # endif
|
---|
116 | # ifndef BBL_MODEL
|
---|
117 | & bustr, bvstr, &
|
---|
118 | # endif
|
---|
119 | & stflx, btflx)
|
---|
120 | !***********************************************************************
|
---|
121 | !
|
---|
122 | USE mod_param
|
---|
123 | USE mod_scalars
|
---|
124 | !
|
---|
125 | USE bc_2d_mod
|
---|
126 | # ifdef DISTRIBUTE
|
---|
127 | USE mp_exchange_mod, ONLY : mp_exchange2d
|
---|
128 | # endif
|
---|
129 | !
|
---|
130 | ! Imported variable declarations.
|
---|
131 | !
|
---|
132 | integer, intent(in) :: ng, Iend, Istr, Jend, Jstr
|
---|
133 | integer, intent(in) :: LBi, UBi, LBj, UBj
|
---|
134 | integer, intent(in) :: nrhs
|
---|
135 | !
|
---|
136 | # ifdef ASSUMED_SHAPE
|
---|
137 | real(r8), intent(in) :: Hz(LBi:,LBj:,:)
|
---|
138 | # if !(defined BBL_MODEL || defined ICESHELF)
|
---|
139 | real(r8), intent(in) :: z_r(LBi:,LBj:,:)
|
---|
140 | real(r8), intent(in) :: z_w(LBi:,LBj:,0:)
|
---|
141 | # endif
|
---|
142 | # if defined ICESHELF
|
---|
143 | real(r8), intent(in) :: zice(LBi:,LBj:)
|
---|
144 | # endif
|
---|
145 | real(r8), intent(in) :: t(LBi:,LBj:,:,:,:)
|
---|
146 | # if !(defined BBL_MODEL || defined ICESHELF)
|
---|
147 | real(r8), intent(in) :: u(LBi:,LBj:,:,:)
|
---|
148 | real(r8), intent(in) :: v(LBi:,LBj:,:,:)
|
---|
149 | # endif
|
---|
150 | # ifdef QCORRECTION
|
---|
151 | real(r8), intent(in) :: dqdt(LBi:,LBj:)
|
---|
152 | real(r8), intent(in) :: sst(LBi:,LBj:)
|
---|
153 | # endif
|
---|
154 | # if defined SCORRECTION || defined SRELAXATION
|
---|
155 | real(r8), intent(in) :: sss(LBi:,LBj:)
|
---|
156 | # endif
|
---|
157 | # if defined ICESHELF
|
---|
158 | # ifdef SHORTWAVE
|
---|
159 | real(r8), intent(inout) :: srflx(LBi:,LBj:)
|
---|
160 | # endif
|
---|
161 | real(r8), intent(inout) :: sustr(LBi:,LBj:)
|
---|
162 | real(r8), intent(inout) :: svstr(LBi:,LBj:)
|
---|
163 | # endif
|
---|
164 | # ifndef BBL_MODEL
|
---|
165 | real(r8), intent(inout) :: bustr(LBi:,LBj:)
|
---|
166 | real(r8), intent(inout) :: bvstr(LBi:,LBj:)
|
---|
167 | # endif
|
---|
168 | real(r8), intent(inout) :: stflx(LBi:,LBj:,:)
|
---|
169 | real(r8), intent(inout) :: btflx(LBi:,LBj:,:)
|
---|
170 | # else
|
---|
171 | real(r8), intent(in) :: Hz(LBi:UBi,LBj:UBj,N(ng))
|
---|
172 | # if !(defined BBL_MODEL || defined ICESHELF)
|
---|
173 | real(r8), intent(in) :: z_r(LBi:UBi,LBj:UBj,N(ng))
|
---|
174 | real(r8), intent(in) :: z_w(LBi:UBi,LBj:UBj,0:N(ng))
|
---|
175 | # endif
|
---|
176 | # if defined ICESHELF
|
---|
177 | real(r8), intent(in) :: zice(LBi:UBi,LBj:UBj)
|
---|
178 | # endif
|
---|
179 | real(r8), intent(in) :: t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
|
---|
180 | # if !(defined BBL_MODEL || defined ICESHELF)
|
---|
181 | real(r8), intent(in) :: u(LBi:UBi,LBj:UBj,N(ng),2)
|
---|
182 | real(r8), intent(in) :: v(LBi:UBi,LBj:UBj,N(ng),2)
|
---|
183 | # endif
|
---|
184 | # ifdef QCORRECTION
|
---|
185 | real(r8), intent(in) :: dqdt(LBi:UBi,LBj:UBj)
|
---|
186 | real(r8), intent(in) :: sst(LBi:UBi,LBj:UBj)
|
---|
187 | # endif
|
---|
188 | # if defined SCORRECTION || defined SRELAXATION
|
---|
189 | real(r8), intent(in) :: sss(LBi:UBi,LBj:UBj)
|
---|
190 | # endif
|
---|
191 | # if defined ICESHELF
|
---|
192 | # ifdef SHORTWAVE
|
---|
193 | real(r8), intent(inout) :: srflx(LBi:UBi,LBj:UBj)
|
---|
194 | # endif
|
---|
195 | real(r8), intent(inout) :: sustr(LBi:UBi,LBj:UBj)
|
---|
196 | real(r8), intent(inout) :: svstr(LBi:UBi,LBj:UBj)
|
---|
197 | # endif
|
---|
198 | # ifndef BBL_MODEL
|
---|
199 | real(r8), intent(inout) :: bustr(LBi:UBi,LBj:UBj)
|
---|
200 | real(r8), intent(inout) :: bvstr(LBi:UBi,LBj:UBj)
|
---|
201 | # endif
|
---|
202 | real(r8), intent(inout) :: stflx(LBi:UBi,LBj:UBj,NT(ng))
|
---|
203 | real(r8), intent(inout) :: btflx(LBi:UBi,LBj:UBj,NT(ng))
|
---|
204 | # endif
|
---|
205 |
|
---|
206 | !
|
---|
207 | ! Local variable declarations.
|
---|
208 | !
|
---|
209 | # ifdef DISTRIBUTE
|
---|
210 | # ifdef EW_PERIODIC
|
---|
211 | logical :: EWperiodic=.TRUE.
|
---|
212 | # else
|
---|
213 | logical :: EWperiodic=.FALSE.
|
---|
214 | # endif
|
---|
215 | # ifdef NS_PERIODIC
|
---|
216 | logical :: NSperiodic=.TRUE.
|
---|
217 | # else
|
---|
218 | logical :: NSperiodic=.FALSE.
|
---|
219 | # endif
|
---|
220 | # endif
|
---|
221 | integer :: IstrR, IendR, JstrR, JendR, IstrU, JstrV
|
---|
222 | integer :: i, j, itrc
|
---|
223 |
|
---|
224 | # if !(defined BBL_MODEL || defined ICESHELF)
|
---|
225 | real(r8) :: cff1, cff2
|
---|
226 | # endif
|
---|
227 |
|
---|
228 | # if !(defined BBL_MODEL || defined ICESHELF) && defined UV_LOGDRAG
|
---|
229 | real(r8), dimension(PRIVATE_2D_SCRATCH_ARRAY) :: wrk
|
---|
230 | # endif
|
---|
231 |
|
---|
232 | # include "set_bounds.h"
|
---|
233 |
|
---|
234 | # ifdef QCORRECTION
|
---|
235 | !
|
---|
236 | !-----------------------------------------------------------------------
|
---|
237 | ! Add in flux correction to surface net heat flux (degC m/s).
|
---|
238 | !-----------------------------------------------------------------------
|
---|
239 | !
|
---|
240 | ! Add in net heat flux correction.
|
---|
241 | !
|
---|
242 | DO j=JstrR,JendR
|
---|
243 | DO i=IstrR,IendR
|
---|
244 | stflx(i,j,itemp)=stflx(i,j,itemp)+ &
|
---|
245 | & dqdt(i,j)*(t(i,j,N(ng),nrhs,itemp)-sst(i,j))
|
---|
246 | END DO
|
---|
247 | END DO
|
---|
248 | # endif
|
---|
249 | # ifdef SALINITY
|
---|
250 | !
|
---|
251 | !-----------------------------------------------------------------------
|
---|
252 | ! Multiply fresh water flux with surface salinity. If appropriate,
|
---|
253 | ! apply correction.
|
---|
254 | !-----------------------------------------------------------------------
|
---|
255 | !
|
---|
256 | DO j=JstrR,JendR
|
---|
257 | DO i=IstrR,IendR
|
---|
258 | # if defined SCORRECTION
|
---|
259 | stflx(i,j,isalt)=stflx(i,j,isalt)*t(i,j,N(ng),nrhs,isalt)- &
|
---|
260 | & Tnudg(isalt,ng)*Hz(i,j,N(ng))* &
|
---|
261 | & (t(i,j,N(ng),nrhs,isalt)-sss(i,j))
|
---|
262 | # elif defined SRELAXATION
|
---|
263 | stflx(i,j,isalt)=-Tnudg(isalt,ng)*Hz(i,j,N(ng))* &
|
---|
264 | & (t(i,j,N(ng),nrhs,isalt)-sss(i,j))
|
---|
265 | # else
|
---|
266 | stflx(i,j,isalt)=stflx(i,j,isalt)*t(i,j,N(ng),nrhs,isalt)
|
---|
267 | # endif
|
---|
268 | btflx(i,j,isalt)=btflx(i,j,isalt)*t(i,j,1,nrhs,isalt)
|
---|
269 | END DO
|
---|
270 | END DO
|
---|
271 | # endif
|
---|
272 | # ifdef ICESHELF
|
---|
273 | !
|
---|
274 | !-----------------------------------------------------------------------
|
---|
275 | ! If ice shelf cavities, zero out for now the surface tracer flux
|
---|
276 | ! over the ice.
|
---|
277 | !-----------------------------------------------------------------------
|
---|
278 | !
|
---|
279 | DO itrc=1,NT(ng)
|
---|
280 | DO j=JstrR,JendR
|
---|
281 | DO i=IstrR,IendR
|
---|
282 | IF (zice(i,j).ne.0.0_r8) THEN
|
---|
283 | stflx(i,j,itrc)=0.0_r8
|
---|
284 | END IF
|
---|
285 | END DO
|
---|
286 | END DO
|
---|
287 | END DO
|
---|
288 | # ifdef SHORTWAVE
|
---|
289 | DO j=JstrR,JendR
|
---|
290 | DO i=IstrR,IendR
|
---|
291 | IF (zice(i,j).ne.0.0_r8) THEN
|
---|
292 | srflx(i,j)=0.0_r8
|
---|
293 | END IF
|
---|
294 | END DO
|
---|
295 | END DO
|
---|
296 | # endif
|
---|
297 | !
|
---|
298 | !-----------------------------------------------------------------------
|
---|
299 | ! If ice shelf cavities, replace surface wind stress with ice shelf
|
---|
300 | ! cavity stress (m2/s2).
|
---|
301 | !-----------------------------------------------------------------------
|
---|
302 |
|
---|
303 | # if defined UV_LOGDRAG
|
---|
304 | !
|
---|
305 | ! Set logarithmic ice shelf cavity stress.
|
---|
306 | !
|
---|
307 | DO j=JstrV-1,Jend
|
---|
308 | DO i=IstrU-1,Iend
|
---|
309 | cff1=1.0_r8/LOG((z_w(i,j,N(ng))-z_r(i,j,N(ng)))/Zob(ng))
|
---|
310 | cff2=vonKar*vonKar*cff1*cff1
|
---|
311 | wrk(i,j)=MIN(Cdb_max,MAX(Cdb_min,cff2))
|
---|
312 | END DO
|
---|
313 | END DO
|
---|
314 | DO j=Jstr,Jend
|
---|
315 | DO i=IstrU,Iend
|
---|
316 | IF (zice(i,j)*zice(i-1,j).ne.0.0_r8) THEN
|
---|
317 | cff1=0.25_r8*(v(i ,j ,N(ng),nrhs)+ &
|
---|
318 | & v(i ,j+1,N(ng),nrhs)+ &
|
---|
319 | & v(i-1,j ,N(ng),nrhs)+ &
|
---|
320 | & v(i-1,j+1,N(ng),nrhs))
|
---|
321 | cff2=SQRT(u(i,j,N(ng),nrhs)*u(i,j,N(ng),nrhs)+cff1*cff1)
|
---|
322 | sustr(i,j)=-0.5_r8*(wrk(i-1,j)+wrk(i,j))* &
|
---|
323 | & u(i,j,N(ng),nrhs)*cff2
|
---|
324 | END IF
|
---|
325 | END DO
|
---|
326 | END DO
|
---|
327 | DO j=JstrV,Jend
|
---|
328 | DO i=Istr,Iend
|
---|
329 | IF (zice(i,j)*zice(i,j-1).ne.0.0_r8) THEN
|
---|
330 | cff1=0.25_r8*(u(i ,j ,N(ng),nrhs)+ &
|
---|
331 | & u(i+1,j ,N(ng),nrhs)+ &
|
---|
332 | & u(i ,j-1,N(ng),nrhs)+ &
|
---|
333 | & u(i+1,j-1,N(ng),nrhs))
|
---|
334 | cff2=SQRT(cff1*cff1+v(i,j,N(ng),nrhs)*v(i,j,N(ng),nrhs))
|
---|
335 | svstr(i,j)=-0.5_r8*(wrk(i,j-1)+wrk(i,j))* &
|
---|
336 | & v(i,j,N(ng),nrhs)*cff2
|
---|
337 | END IF
|
---|
338 | END DO
|
---|
339 | END DO
|
---|
340 | # elif defined UV_QDRAG
|
---|
341 | !
|
---|
342 | ! Set quadratic ice shelf cavity stress.
|
---|
343 | !
|
---|
344 | DO j=Jstr,Jend
|
---|
345 | DO i=IstrU,Iend
|
---|
346 | IF (zice(i,j)*zice(i-1,j).ne.0.0_r8) THEN
|
---|
347 | cff1=0.25_r8*(v(i ,j ,N(ng),nrhs)+ &
|
---|
348 | & v(i ,j+1,N(ng),nrhs)+ &
|
---|
349 | & v(i-1,j ,N(ng),nrhs)+ &
|
---|
350 | & v(i-1,j+1,N(ng),nrhs))
|
---|
351 | cff2=SQRT(u(i,j,N(ng),nrhs)*u(i,j,N(ng),nrhs)+cff2*cff2)
|
---|
352 | sustr(i,j)=-rdrg2(ng)*u(i,j,N(ng),nrhs)*cff2
|
---|
353 | END IF
|
---|
354 | END DO
|
---|
355 | END DO
|
---|
356 | DO j=JstrV,Jend
|
---|
357 | DO i=Istr,Iend
|
---|
358 | IF (zice(i,j)*zice(i,j-1).ne.0.0_r8) THEN
|
---|
359 | cff1=0.25_r8*(u(i ,j ,N(ng),nrhs)+ &
|
---|
360 | & u(i+1,j ,N(ng),nrhs)+ &
|
---|
361 | & u(i ,j-1,N(ng),nrhs)+ &
|
---|
362 | & u(i+1,j-1,N(ng),nrhs))
|
---|
363 | cff2=SQRT(cff1*cff1+v(i,j,N(ng),nrhs)*v(i,j,N(ng),nrhs))
|
---|
364 | svstr(i,j)=-rdrg2(ng)*v(i,j,N(ng),nrhs)*cff2
|
---|
365 | END IF
|
---|
366 | END DO
|
---|
367 | END DO
|
---|
368 | # elif defined UV_LDRAG
|
---|
369 | !
|
---|
370 | ! Set linear ice shelf cavity stress.
|
---|
371 | !
|
---|
372 | DO j=Jstr,Jend
|
---|
373 | DO i=IstrU,Iend
|
---|
374 | IF (zice(i,j)*zice(i-1,j).ne.0.0_r8) THEN
|
---|
375 | sustr(i,j)=-rdrg(ng)*u(i,j,N(ng),nrhs)
|
---|
376 | END IF
|
---|
377 | END DO
|
---|
378 | END DO
|
---|
379 | DO j=JstrV,Jend
|
---|
380 | DO i=Istr,Iend
|
---|
381 | IF (zice(i,j)*zice(i,j-1).ne.0.0_r8) THEN
|
---|
382 | svstr(i,j)=-rdrg(ng)*v(i,j,N(ng),nrhs)
|
---|
383 | END IF
|
---|
384 | END DO
|
---|
385 | END DO
|
---|
386 | # else
|
---|
387 | DO j=Jstr,Jend
|
---|
388 | DO i=IstrU,Iend
|
---|
389 | IF (zice(i,j)*zice(i-1,j).ne.0.0_r8) THEN
|
---|
390 | sustr(i,j)=0.0_r8
|
---|
391 | END IF
|
---|
392 | END DO
|
---|
393 | END DO
|
---|
394 | DO j=JstrV,Jend
|
---|
395 | DO i=Istr,Iend
|
---|
396 | IF (zice(i,j)*zice(i,j-1).ne.0.0_r8) THEN
|
---|
397 | svstr(i,j)=0.0_r8
|
---|
398 | END IF
|
---|
399 | END DO
|
---|
400 | END DO
|
---|
401 | # endif
|
---|
402 | !
|
---|
403 | ! Apply periodic or gradient boundary conditions for output
|
---|
404 | ! purposes only.
|
---|
405 | !
|
---|
406 | CALL bc_u2d_tile (ng, Istr, Iend, Jstr, Jend, &
|
---|
407 | & LBi, UBi, LBj, UBj, &
|
---|
408 | & sustr)
|
---|
409 | CALL bc_v2d_tile (ng, Istr, Iend, Jstr, Jend, &
|
---|
410 | & LBi, UBi, LBj, UBj, &
|
---|
411 | & svstr)
|
---|
412 | # ifdef DISTRIBUTE
|
---|
413 | CALL mp_exchange2d (ng, iNLM, 2, Istr, Iend, Jstr, Jend, &
|
---|
414 | & LBi, UBi, LBj, UBj, &
|
---|
415 | & NghostPoints, EWperiodic, NSperiodic, &
|
---|
416 | & sustr, svstr)
|
---|
417 | # endif
|
---|
418 | # endif
|
---|
419 | # ifndef BBL_MODEL
|
---|
420 | !
|
---|
421 | !-----------------------------------------------------------------------
|
---|
422 | ! Set kinematic bottom momentum flux (m2/s2).
|
---|
423 | !-----------------------------------------------------------------------
|
---|
424 |
|
---|
425 | # if defined UV_LOGDRAG
|
---|
426 | !
|
---|
427 | ! Set logarithmic bottom stress.
|
---|
428 | !
|
---|
429 | DO j=JstrV-1,Jend
|
---|
430 | DO i=IstrU-1,Iend
|
---|
431 | cff1=1.0_r8/LOG((z_r(i,j,1)-z_w(i,j,0))/Zob(ng))
|
---|
432 | cff2=vonKar*vonKar*cff1*cff1
|
---|
433 | wrk(i,j)=MIN(Cdb_max,MAX(Cdb_min,cff2))
|
---|
434 | END DO
|
---|
435 | END DO
|
---|
436 | DO j=Jstr,Jend
|
---|
437 | DO i=IstrU,Iend
|
---|
438 | cff1=0.25_r8*(v(i ,j ,1,nrhs)+ &
|
---|
439 | & v(i ,j+1,1,nrhs)+ &
|
---|
440 | & v(i-1,j ,1,nrhs)+ &
|
---|
441 | & v(i-1,j+1,1,nrhs))
|
---|
442 | cff2=SQRT(u(i,j,1,nrhs)*u(i,j,1,nrhs)+cff1*cff1)
|
---|
443 | bustr(i,j)=0.5_r8*(wrk(i-1,j)+wrk(i,j))* &
|
---|
444 | & u(i,j,1,nrhs)*cff2
|
---|
445 | END DO
|
---|
446 | END DO
|
---|
447 | DO j=JstrV,Jend
|
---|
448 | DO i=Istr,Iend
|
---|
449 | cff1=0.25_r8*(u(i ,j ,1,nrhs)+ &
|
---|
450 | & u(i+1,j ,1,nrhs)+ &
|
---|
451 | & u(i ,j-1,1,nrhs)+ &
|
---|
452 | & u(i+1,j-1,1,nrhs))
|
---|
453 | cff2=SQRT(cff1*cff1+v(i,j,1,nrhs)*v(i,j,1,nrhs))
|
---|
454 | bvstr(i,j)=0.5_r8*(wrk(i,j-1)+wrk(i,j))* &
|
---|
455 | & v(i,j,1,nrhs)*cff2
|
---|
456 | END DO
|
---|
457 | END DO
|
---|
458 | # elif defined UV_QDRAG
|
---|
459 | !
|
---|
460 | ! Set quadratic bottom stress.
|
---|
461 | !
|
---|
462 | DO j=Jstr,Jend
|
---|
463 | DO i=IstrU,Iend
|
---|
464 | cff1=0.25_r8*(v(i ,j ,1,nrhs)+ &
|
---|
465 | & v(i ,j+1,1,nrhs)+ &
|
---|
466 | & v(i-1,j ,1,nrhs)+ &
|
---|
467 | & v(i-1,j+1,1,nrhs))
|
---|
468 | cff2=SQRT(u(i,j,1,nrhs)*u(i,j,1,nrhs)+cff1*cff1)
|
---|
469 | bustr(i,j)=rdrg2(ng)*u(i,j,1,nrhs)*cff2
|
---|
470 | END DO
|
---|
471 | END DO
|
---|
472 | DO j=JstrV,Jend
|
---|
473 | DO i=Istr,Iend
|
---|
474 | cff1=0.25_r8*(u(i ,j ,1,nrhs)+ &
|
---|
475 | & u(i+1,j ,1,nrhs)+ &
|
---|
476 | & u(i ,j-1,1,nrhs)+ &
|
---|
477 | & u(i+1,j-1,1,nrhs))
|
---|
478 | cff2=SQRT(cff1*cff1+v(i,j,1,nrhs)*v(i,j,1,nrhs))
|
---|
479 | bvstr(i,j)=rdrg2(ng)*v(i,j,1,nrhs)*cff2
|
---|
480 | END DO
|
---|
481 | END DO
|
---|
482 | # elif defined UV_LDRAG
|
---|
483 | !
|
---|
484 | ! Set linear bottom stress.
|
---|
485 | !
|
---|
486 | DO j=Jstr,Jend
|
---|
487 | DO i=IstrU,Iend
|
---|
488 | bustr(i,j)=rdrg(ng)*u(i,j,1,nrhs)
|
---|
489 | END DO
|
---|
490 | END DO
|
---|
491 | DO j=JstrV,Jend
|
---|
492 | DO i=Istr,Iend
|
---|
493 | bvstr(i,j)=rdrg(ng)*v(i,j,1,nrhs)
|
---|
494 | END DO
|
---|
495 | END DO
|
---|
496 | # endif
|
---|
497 | !
|
---|
498 | ! Apply boundary conditions.
|
---|
499 | !
|
---|
500 | CALL bc_u2d_tile (ng, Istr, Iend, Jstr, Jend, &
|
---|
501 | & LBi, UBi, LBj, UBj, &
|
---|
502 | & bustr)
|
---|
503 | CALL bc_v2d_tile (ng, Istr, Iend, Jstr, Jend, &
|
---|
504 | & LBi, UBi, LBj, UBj, &
|
---|
505 | & bvstr)
|
---|
506 | # ifdef DISTRIBUTE
|
---|
507 | CALL mp_exchange2d (ng, iNLM, 2, Istr, Iend, Jstr, Jend, &
|
---|
508 | & LBi, UBi, LBj, UBj, &
|
---|
509 | & NghostPoints, EWperiodic, NSperiodic, &
|
---|
510 | & bustr, bvstr)
|
---|
511 | # endif
|
---|
512 | # endif
|
---|
513 | RETURN
|
---|
514 | END SUBROUTINE set_vbc_tile
|
---|
515 |
|
---|
516 | # else
|
---|
517 |
|
---|
518 | !
|
---|
519 | !***********************************************************************
|
---|
520 | SUBROUTINE set_vbc (ng, tile)
|
---|
521 | !***********************************************************************
|
---|
522 | !
|
---|
523 | USE mod_param
|
---|
524 | USE mod_grid
|
---|
525 | USE mod_forces
|
---|
526 | USE mod_ocean
|
---|
527 | USE mod_stepping
|
---|
528 | !
|
---|
529 | ! Imported variable declarations.
|
---|
530 | !
|
---|
531 | integer, intent(in) :: ng, tile
|
---|
532 | !
|
---|
533 | ! Local variable declarations.
|
---|
534 | !
|
---|
535 | # include "tile.h"
|
---|
536 | !
|
---|
537 | # ifdef PROFILE
|
---|
538 | CALL wclock_on (ng, iNLM, 6)
|
---|
539 | # endif
|
---|
540 | CALL set_vbc_tile (ng, Istr, Iend, Jstr, Jend, &
|
---|
541 | & LBi, UBi, LBj, UBj, &
|
---|
542 | & krhs(ng), kstp(ng), knew(ng), &
|
---|
543 | & GRID(ng) % om_u, &
|
---|
544 | & GRID(ng) % om_v, &
|
---|
545 | & GRID(ng) % on_u, &
|
---|
546 | & GRID(ng) % on_v, &
|
---|
547 | & OCEAN(ng) % ubar, &
|
---|
548 | & OCEAN(ng) % vbar, &
|
---|
549 | & FORCES(ng) % bustr, &
|
---|
550 | & FORCES(ng) % bvstr)
|
---|
551 | # ifdef PROFILE
|
---|
552 | CALL wclock_off (ng, iNLM, 6)
|
---|
553 | # endif
|
---|
554 | RETURN
|
---|
555 | END SUBROUTINE set_vbc
|
---|
556 | !
|
---|
557 | !***********************************************************************
|
---|
558 | SUBROUTINE set_vbc_tile (ng, Istr, Iend, Jstr, Jend, &
|
---|
559 | & LBi, UBi, LBj, UBj, &
|
---|
560 | & krhs, kstp, knew, &
|
---|
561 | & om_u, om_v, on_u, on_v, &
|
---|
562 | & ubar, vbar, bustr, bvstr)
|
---|
563 | !***********************************************************************
|
---|
564 | !
|
---|
565 | USE mod_param
|
---|
566 | USE mod_scalars
|
---|
567 | !
|
---|
568 | USE bc_2d_mod
|
---|
569 | # ifdef DISTRIBUTE
|
---|
570 | USE mp_exchange_mod, ONLY : mp_exchange2d
|
---|
571 | # endif
|
---|
572 | !
|
---|
573 | ! Imported variable declarations.
|
---|
574 | !
|
---|
575 | integer, intent(in) :: ng, Iend, Istr, Jend, Jstr
|
---|
576 | integer, intent(in) :: LBi, UBi, LBj, UBj
|
---|
577 | integer, intent(in) :: krhs, kstp, knew
|
---|
578 | !
|
---|
579 | # ifdef ASSUMED_SHAPE
|
---|
580 | real(r8), intent(in) :: om_u(LBi:,LBj:)
|
---|
581 | real(r8), intent(in) :: om_v(LBi:,LBj:)
|
---|
582 | real(r8), intent(in) :: on_u(LBi:,LBj:)
|
---|
583 | real(r8), intent(in) :: on_v(LBi:,LBj:)
|
---|
584 | real(r8), intent(in) :: ubar(LBi:,LBj:,:)
|
---|
585 | real(r8), intent(in) :: vbar(LBi:,LBj:,:)
|
---|
586 | real(r8), intent(inout) :: bustr(LBi:,LBj:)
|
---|
587 | real(r8), intent(inout) :: bvstr(LBi:,LBj:)
|
---|
588 | # else
|
---|
589 | real(r8), intent(in) :: om_u(LBi:UBi,LBj:UBj)
|
---|
590 | real(r8), intent(in) :: om_v(LBi:UBi,LBj:UBj)
|
---|
591 | real(r8), intent(in) :: on_u(LBi:UBi,LBj:UBj)
|
---|
592 | real(r8), intent(in) :: on_v(LBi:UBi,LBj:UBj)
|
---|
593 | real(r8), intent(in) :: ubar(LBi:UBi,LBj:UBj,3)
|
---|
594 | real(r8), intent(in) :: vbar(LBi:UBi,LBj:UBj,3)
|
---|
595 | real(r8), intent(inout) :: bustr(LBi:UBi,LBj:UBj)
|
---|
596 | real(r8), intent(inout) :: bvstr(LBi:UBi,LBj:UBj)
|
---|
597 | # endif
|
---|
598 | !
|
---|
599 | ! Local variable declarations.
|
---|
600 | !
|
---|
601 | # ifdef DISTRIBUTE
|
---|
602 | # ifdef EW_PERIODIC
|
---|
603 | logical :: EWperiodic=.TRUE.
|
---|
604 | # else
|
---|
605 | logical :: EWperiodic=.FALSE.
|
---|
606 | # endif
|
---|
607 | # ifdef NS_PERIODIC
|
---|
608 | logical :: NSperiodic=.TRUE.
|
---|
609 | # else
|
---|
610 | logical :: NSperiodic=.FALSE.
|
---|
611 | # endif
|
---|
612 | # endif
|
---|
613 | integer :: IstrR, IendR, JstrR, JendR, IstrU, JstrV
|
---|
614 | integer :: i, j
|
---|
615 |
|
---|
616 | real(r8) :: cff1, cff2
|
---|
617 |
|
---|
618 | # include "set_bounds.h"
|
---|
619 | !
|
---|
620 | !-----------------------------------------------------------------------
|
---|
621 | ! Set kinematic barotropic bottom momentum stress (m2/s2).
|
---|
622 | !-----------------------------------------------------------------------
|
---|
623 |
|
---|
624 | # if defined UV_LDRAG
|
---|
625 | !
|
---|
626 | ! Set linear bottom stress.
|
---|
627 | !
|
---|
628 | DO j=Jstr,Jend
|
---|
629 | DO i=IstrU,Iend
|
---|
630 | bustr(i,j)=rdrg(ng)*ubar(i,j,krhs)
|
---|
631 | END DO
|
---|
632 | END DO
|
---|
633 | DO j=JstrV,Jend
|
---|
634 | DO i=Istr,Iend
|
---|
635 | bvstr(i,j)=rdrg(ng)*vbar(i,j,krhs)
|
---|
636 | END DO
|
---|
637 | END DO
|
---|
638 | # elif defined UV_QDRAG
|
---|
639 | !
|
---|
640 | ! Set quadratic bottom stress.
|
---|
641 | !
|
---|
642 | DO j=Jstr,Jend
|
---|
643 | DO i=IstrU,Iend
|
---|
644 | cff1=0.25_r8*(vbar(i ,j ,krhs)+ &
|
---|
645 | & vbar(i ,j+1,krhs)+ &
|
---|
646 | & vbar(i-1,j ,krhs)+ &
|
---|
647 | & vbar(i-1,j+1,krhs))
|
---|
648 | cff2=SQRT(ubar(i,j,krhs)*ubar(i,j,krhs)+cff1*cff1)
|
---|
649 | bustr(i,j)=rdrg2(ng)*ubar(i,j,krhs)*cff2
|
---|
650 | END DO
|
---|
651 | END DO
|
---|
652 | DO j=JstrV,Jend
|
---|
653 | DO i=Istr,Iend
|
---|
654 | cff1=0.25_r8*(ubar(i ,j ,krhs)+ &
|
---|
655 | & ubar(i+1,j ,krhs)+ &
|
---|
656 | & ubar(i ,j-1,krhs)+ &
|
---|
657 | & ubar(i+1,j-1,krhs))
|
---|
658 | cff2=SQRT(cff1*cff1+vbar(i,j,krhs)*vbar(i,j,krhs))
|
---|
659 | bvstr(i,j)=rdrg2(ng)*vbar(i,j,krhs)*cff2
|
---|
660 | END DO
|
---|
661 | END DO
|
---|
662 | # endif
|
---|
663 | !
|
---|
664 | ! Apply boundary conditions.
|
---|
665 | !
|
---|
666 | CALL bc_u2d_tile (ng, Istr, Iend, Jstr, Jend, &
|
---|
667 | & LBi, UBi, LBj, UBj, &
|
---|
668 | & bustr)
|
---|
669 | CALL bc_v2d_tile (ng, Istr, Iend, Jstr, Jend, &
|
---|
670 | & LBi, UBi, LBj, UBj, &
|
---|
671 | & bvstr)
|
---|
672 | # ifdef DISTRIBUTE
|
---|
673 | CALL mp_exchange2d (ng, iNLM, 2, Istr, Iend, Jstr, Jend, &
|
---|
674 | & LBi, UBi, LBj, UBj, &
|
---|
675 | & NghostPoints, EWperiodic, NSperiodic, &
|
---|
676 | & bustr, bvstr)
|
---|
677 | # endif
|
---|
678 | RETURN
|
---|
679 | END SUBROUTINE set_vbc_tile
|
---|
680 | # endif
|
---|
681 | #endif
|
---|
682 | END MODULE set_vbc_mod
|
---|