ROMS
Loading...
Searching...
No Matches
lanc_resid_mod Module Reference

Functions/Subroutines

subroutine, public lanc_resid (ng, tile, linp, lout)
 
subroutine lanc_resid_tile (ng, tile, lbi, ubi, lbj, ubj, lbij, ubij, imins, imaxs, jmins, jmaxs, linp, lout, tl_t_obc, tl_u_obc, tl_v_obc, tl_ubar_obc, tl_vbar_obc, tl_zeta_obc, tl_ustr, tl_vstr, tl_tflux, tl_t, tl_u, tl_v, tl_zeta, ad_t_obc, ad_u_obc, ad_v_obc, ad_ubar_obc, ad_vbar_obc, ad_zeta_obc, ad_ustr, ad_vstr, ad_tflux, ad_t, ad_u, ad_v, ad_zeta)
 

Function/Subroutine Documentation

◆ lanc_resid()

subroutine, public lanc_resid_mod::lanc_resid ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) linp,
integer, intent(in) lout )

Definition at line 26 of file lanc_resid.F.

27!***********************************************************************
28!
29 USE mod_param
30# ifdef ADJUST_BOUNDARY
31 USE mod_boundary
32# endif
33# if defined ADJUST_STFLUX || defined ADJUST_WSTRESS
34 USE mod_forces
35# endif
36 USE mod_ocean
37!
38! Imported variable declarations.
39!
40 integer, intent(in) :: ng, tile, Linp, Lout
41!
42! Local variable declarations.
43!
44# include "tile.h"
45!
46 CALL lanc_resid_tile (ng, tile, &
47 & lbi, ubi, lbj, ubj, lbij, ubij, &
48 & imins, imaxs, jmins, jmaxs, &
49 & linp, lout, &
50# ifdef ADJUST_BOUNDARY
51# ifdef SOLVE3D
52 & boundary(ng) % tl_t_obc, &
53 & boundary(ng) % tl_u_obc, &
54 & boundary(ng) % tl_v_obc, &
55# endif
56 & boundary(ng) % tl_ubar_obc, &
57 & boundary(ng) % tl_vbar_obc, &
58 & boundary(ng) % tl_zeta_obc, &
59# endif
60# ifdef ADJUST_WSTRESS
61 & forces(ng) % tl_ustr, &
62 & forces(ng) % tl_vstr, &
63# endif
64# ifdef SOLVE3D
65# ifdef ADJUST_STFLUX
66 & forces(ng) % tl_tflux, &
67# endif
68 & ocean(ng) % tl_t, &
69 & ocean(ng) % tl_u, &
70 & ocean(ng) % tl_v, &
71# else
72 & ocean(ng) % tl_ubar, &
73 & ocean(ng) % tl_vbar, &
74# endif
75 & ocean(ng) % tl_zeta, &
76# ifdef ADJUST_BOUNDARY
77# ifdef SOLVE3D
78 & boundary(ng) % ad_t_obc, &
79 & boundary(ng) % ad_u_obc, &
80 & boundary(ng) % ad_v_obc, &
81# endif
82 & boundary(ng) % ad_ubar_obc, &
83 & boundary(ng) % ad_vbar_obc, &
84 & boundary(ng) % ad_zeta_obc, &
85# endif
86# ifdef ADJUST_WSTRESS
87 & forces(ng) % ad_ustr, &
88 & forces(ng) % ad_vstr, &
89# endif
90# ifdef SOLVE3D
91# ifdef ADJUST_STFLUX
92 & forces(ng) % ad_tflux, &
93# endif
94 & ocean(ng) % ad_t, &
95 & ocean(ng) % ad_u, &
96 & ocean(ng) % ad_v, &
97# else
98 & ocean(ng) % ad_ubar, &
99 & ocean(ng) % ad_vbar, &
100# endif
101 & ocean(ng) % ad_zeta)
102 RETURN
type(t_boundary), dimension(:), allocatable boundary
type(t_forces), dimension(:), allocatable forces
Definition mod_forces.F:554
type(t_ocean), dimension(:), allocatable ocean
Definition mod_ocean.F:351

References mod_boundary::boundary, mod_forces::forces, lanc_resid_tile(), and mod_ocean::ocean.

Here is the call graph for this function:

◆ lanc_resid_tile()

subroutine lanc_resid_mod::lanc_resid_tile ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj,
integer, intent(in) lbij,
integer, intent(in) ubij,
integer, intent(in) imins,
integer, intent(in) imaxs,
integer, intent(in) jmins,
integer, intent(in) jmaxs,
integer, intent(in) linp,
integer, intent(in) lout,
real(r8), dimension(lbij:,:,:,:,:,:), intent(in) tl_t_obc,
real(r8), dimension(lbij:,:,:,:,:), intent(in) tl_u_obc,
real(r8), dimension(lbij:,:,:,:,:), intent(in) tl_v_obc,
real(r8), dimension(lbij:,:,:,:), intent(in) tl_ubar_obc,
real(r8), dimension(lbij:,:,:,:), intent(in) tl_vbar_obc,
real(r8), dimension(lbij:,:,:,:), intent(in) tl_zeta_obc,
real(r8), dimension(lbi:,lbj:,:,:), intent(in) tl_ustr,
real(r8), dimension(lbi:,lbj:,:,:), intent(in) tl_vstr,
real(r8), dimension(lbi:,lbj:,:,:,:), intent(in) tl_tflux,
real(r8), dimension(lbi:,lbj:,:,:,:), intent(in) tl_t,
real(r8), dimension(lbi:,lbj:,:,:), intent(in) tl_u,
real(r8), dimension(lbi:,lbj:,:,:), intent(in) tl_v,
real(r8), dimension(lbi:,lbj:,:), intent(in) tl_zeta,
real(r8), dimension(lbij:,:,:,:,:,:), intent(inout) ad_t_obc,
real(r8), dimension(lbij:,:,:,:,:), intent(inout) ad_u_obc,
real(r8), dimension(lbij:,:,:,:,:), intent(inout) ad_v_obc,
real(r8), dimension(lbij:,:,:,:), intent(inout) ad_ubar_obc,
real(r8), dimension(lbij:,:,:,:), intent(inout) ad_vbar_obc,
real(r8), dimension(lbij:,:,:,:), intent(inout) ad_zeta_obc,
real(r8), dimension(lbi:,lbj:,:,:), intent(inout) ad_ustr,
real(r8), dimension(lbi:,lbj:,:,:), intent(inout) ad_vstr,
real(r8), dimension(lbi:,lbj:,:,:,:), intent(inout) ad_tflux,
real(r8), dimension(lbi:,lbj:,:,:,:), intent(inout) ad_t,
real(r8), dimension(lbi:,lbj:,:,:), intent(inout) ad_u,
real(r8), dimension(lbi:,lbj:,:,:), intent(inout) ad_v,
real(r8), dimension(lbi:,lbj:,:), intent(inout) ad_zeta )
private

Definition at line 106 of file lanc_resid.F.

148!***********************************************************************
149!
150 USE mod_param
151 USE mod_ncparam
152# if defined ADJUST_STFLUX || defined ADJUST_WSTRESS || \
153 defined adjust_boundary
154 USE mod_scalars
155# endif
156!
157! Imported variable declarations.
158!
159 integer, intent(in) :: ng, tile
160 integer, intent(in) :: LBi, UBi, LBj, UBj, LBij, UBij
161 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
162 integer, intent(in) :: Linp, Lout
163!
164# ifdef ASSUMED_SHAPE
165# ifdef ADJUST_BOUNDARY
166# ifdef SOLVE3D
167 real(r8), intent(in) :: tl_t_obc(LBij:,:,:,:,:,:)
168 real(r8), intent(in) :: tl_u_obc(LBij:,:,:,:,:)
169 real(r8), intent(in) :: tl_v_obc(LBij:,:,:,:,:)
170# endif
171 real(r8), intent(in) :: tl_ubar_obc(LBij:,:,:,:)
172 real(r8), intent(in) :: tl_vbar_obc(LBij:,:,:,:)
173 real(r8), intent(in) :: tl_zeta_obc(LBij:,:,:,:)
174# endif
175# ifdef ADJUST_WSTRESS
176 real(r8), intent(in) :: tl_ustr(LBi:,LBj:,:,:)
177 real(r8), intent(in) :: tl_vstr(LBi:,LBj:,:,:)
178# endif
179# ifdef SOLVE3D
180# ifdef ADJUST_STFLUX
181 real(r8), intent(in) :: tl_tflux(LBi:,LBj:,:,:,:)
182# endif
183 real(r8), intent(in) :: tl_t(LBi:,LBj:,:,:,:)
184 real(r8), intent(in) :: tl_u(LBi:,LBj:,:,:)
185 real(r8), intent(in) :: tl_v(LBi:,LBj:,:,:)
186# else
187 real(r8), intent(in) :: tl_ubar(LBi:,LBj:,:)
188 real(r8), intent(in) :: tl_vbar(LBi:,LBj:,:)
189# endif
190 real(r8), intent(in) :: tl_zeta(LBi:,LBj:,:)
191# ifdef ADJUST_BOUNDARY
192# ifdef SOLVE3D
193 real(r8), intent(inout) :: ad_t_obc(LBij:,:,:,:,:,:)
194 real(r8), intent(inout) :: ad_u_obc(LBij:,:,:,:,:)
195 real(r8), intent(inout) :: ad_v_obc(LBij:,:,:,:,:)
196# endif
197 real(r8), intent(inout) :: ad_ubar_obc(LBij:,:,:,:)
198 real(r8), intent(inout) :: ad_vbar_obc(LBij:,:,:,:)
199 real(r8), intent(inout) :: ad_zeta_obc(LBij:,:,:,:)
200# endif
201# ifdef ADJUST_WSTRESS
202 real(r8), intent(inout) :: ad_ustr(LBi:,LBj:,:,:)
203 real(r8), intent(inout) :: ad_vstr(LBi:,LBj:,:,:)
204# endif
205# ifdef SOLVE3D
206# ifdef ADJUST_STFLUX
207 real(r8), intent(inout) :: ad_tflux(LBi:,LBj:,:,:,:)
208# endif
209 real(r8), intent(inout) :: ad_t(LBi:,LBj:,:,:,:)
210 real(r8), intent(inout) :: ad_u(LBi:,LBj:,:,:)
211 real(r8), intent(inout) :: ad_v(LBi:,LBj:,:,:)
212# else
213 real(r8), intent(inout) :: ad_ubar(LBi:,LBj:,:)
214 real(r8), intent(inout) :: ad_vbar(LBi:,LBj:,:)
215# endif
216 real(r8), intent(inout) :: ad_zeta(LBi:,LBj:,:)
217# else
218# ifdef ADJUST_WSTRESS
219 real(r8), intent(in) :: tl_ustr(LBi:UBi,LBj:UBj,Nfrec(ng),2)
220 real(r8), intent(in) :: tl_vstr(LBi:UBi,LBj:UBj,Nfrec(ng),2)
221# endif
222# ifdef ADJUST_BOUNDARY
223# ifdef SOLVE3D
224 real(r8), intent(in) :: tl_t_obc(LBij:UBij,N(ng),4, &
225 & Nbrec(ng),2,NT(ng))
226 real(r8), intent(in) :: tl_u_obc(LBij:UBij,N(ng),4,Nbrec(ng),2)
227 real(r8), intent(in) :: tl_v_obc(LBij:UBij,N(ng),4,Nbrec(ng),2)
228# endif
229 real(r8), intent(in) :: tl_ubar_obc(LBij:UBij,4,Nbrec(ng),2)
230 real(r8), intent(in) :: tl_vbar_obc(LBij:UBij,4,Nbrec(ng),2)
231 real(r8), intent(in) :: tl_zeta_obc(LBij:UBij,4,Nbrec(ng),2)
232# endif
233# ifdef SOLVE3D
234# ifdef ADJUST_STFLUX
235 real(r8), intent(in) :: tl_tflux(LBi:UBi,LBj:UBj, &
236 & Nfrec(ng),2,NT(ng))
237# endif
238 real(r8), intent(in) :: tl_t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
239 real(r8), intent(in) :: tl_u(LBi:UBi,LBj:UBj,N(ng),2)
240 real(r8), intent(in) :: tl_v(LBi:UBi,LBj:UBj,N(ng),2)
241# else
242 real(r8), intent(in) :: tl_ubar(LBi:UBi,LBj:UBj,:)
243 real(r8), intent(in) :: tl_vbar(LBi:UBi,LBj:UBj,:)
244# endif
245 real(r8), intent(in) :: tl_zeta(LBi:UBi,LBj:UBj,:)
246# ifdef ADJUST_BOUNDARY
247# ifdef SOLVE3D
248 real(r8), intent(inout) :: ad_t_obc(LBij:UBij,N(ng),4, &
249 & Nbrec(ng),2,NT(ng))
250 real(r8), intent(inout) :: ad_u_obc(LBij:UBij,N(ng),4,Nbrec(ng),2)
251 real(r8), intent(inout) :: ad_v_obc(LBij:UBij,N(ng),4,Nbrec(ng),2)
252# endif
253 real(r8), intent(inout) :: ad_ubar_obc(LBij:UBij,4,Nbrec(ng),2)
254 real(r8), intent(inout) :: ad_vbar_obc(LBij:UBij,4,Nbrec(ng),2)
255 real(r8), intent(inout) :: ad_zeta_obc(LBij:UBij,4,Nbrec(ng),2)
256# endif
257# ifdef ADJUST_WSTRESS
258 real(r8), intent(inout) :: ad_ustr(LBi:UBi,LBj:UBj,Nfrec(ng),2)
259 real(r8), intent(inout) :: ad_vstr(LBi:UBi,LBj:UBj,Nfrec(ng),2)
260# endif
261# ifdef SOLVE3D
262# ifdef ADJUST_STFLUX
263 real(r8), intent(inout) :: ad_tflux(LBi:UBi,LBj:UBj, &
264 & Nfrec(ng),2,NT(ng))
265# endif
266 real(r8), intent(inout) :: ad_t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
267 real(r8), intent(inout) :: ad_u(LBi:UBi,LBj:UBj,N(ng),2)
268 real(r8), intent(inout) :: ad_v(LBi:UBi,LBj:UBj,N(ng),2)
269# else
270 real(r8), intent(inout) :: ad_ubar(LBi:UBi,LBj:UBj,:)
271 real(r8), intent(inout) :: ad_vbar(LBi:UBi,LBj:UBj,:)
272# endif
273 real(r8), intent(inout) :: ad_zeta(LBi:UBi,LBj:UBj,:)
274# endif
275!
276! Local variable declarations.
277!
278 integer :: i, ib, ir, j, k
279# ifdef SOLVE3D
280 integer :: itrc
281# endif
282
283# include "set_bounds.h"
284!
285!-----------------------------------------------------------------------
286! Compute total Lanczos residual at the start of each outer-loop.
287!-----------------------------------------------------------------------
288!
289! Free-surface.
290!
291 DO j=jstrr,jendr
292 DO i=istrr,iendr
293 ad_zeta(i,j,lout)=-tl_zeta(i,j,linp)+ &
294 & ad_zeta(i,j,lout )
295 END DO
296 END DO
297
298# ifdef ADJUST_BOUNDARY
299!
300! Free-surface open boundaries.
301!
302 IF (any(lobc(:,isfsur,ng))) THEN
303 DO ir=1,nbrec(ng)
304 IF ((lobc(iwest,isfsur,ng)).and. &
305 & domain(ng)%Western_Edge(tile)) THEN
306 ib=iwest
307 DO j=jstr,jend
308 ad_zeta_obc(j,ib,ir,lout)=-tl_zeta_obc(j,ib,ir,linp)+ &
309 & ad_zeta_obc(j,ib,ir,lout )
310 END DO
311 END IF
312 IF ((lobc(ieast,isfsur,ng)).and. &
313 & domain(ng)%Eastern_Edge(tile)) THEN
314 ib=ieast
315 DO j=jstr,jend
316 ad_zeta_obc(j,ib,ir,lout)=-tl_zeta_obc(j,ib,ir,linp)+ &
317 & ad_zeta_obc(j,ib,ir,lout )
318 END DO
319 END IF
320 IF ((lobc(isouth,isfsur,ng)).and. &
321 & domain(ng)%Southern_Edge(tile)) THEN
322 ib=isouth
323 DO i=istr,iend
324 ad_zeta_obc(i,ib,ir,lout)=-tl_zeta_obc(i,ib,ir,linp)+ &
325 & ad_zeta_obc(i,ib,ir,lout )
326 END DO
327 END IF
328 IF ((lobc(inorth,isfsur,ng)).and. &
329 & domain(ng)%Northern_Edge(tile)) THEN
330 ib=inorth
331 DO i=istr,iend
332 ad_zeta_obc(i,ib,ir,lout)=-tl_zeta_obc(i,ib,ir,linp)+ &
333 & ad_zeta_obc(i,ib,ir,lout )
334 END DO
335 END IF
336 END DO
337 END IF
338# endif
339
340# ifndef SOLVE3D
341!
342! 2D U-momentum component.
343!
344 DO j=jstrr,jendr
345 DO i=istr,iendr
346 ad_ubar(i,j,lout)=-tl_ubar(i,j,linp)+ &
347 & ad_ubar(i,j,lout )
348 END DO
349 END DO
350# endif
351
352# ifdef ADJUST_BOUNDARY
353!
354! 2D U-momentum open boundaries.
355!
356 IF (any(lobc(:,isubar,ng))) THEN
357 DO ir=1,nbrec(ng)
358 IF ((lobc(iwest,isubar,ng)).and. &
359 & domain(ng)%Western_Edge(tile)) THEN
360 ib=iwest
361 DO j=jstr,jend
362 ad_ubar_obc(j,ib,ir,lout)=-tl_ubar_obc(j,ib,ir,linp)+ &
363 & ad_ubar_obc(j,ib,ir,lout )
364 END DO
365 END IF
366 IF ((lobc(ieast,isubar,ng)).and. &
367 & domain(ng)%Eastern_Edge(tile)) THEN
368 ib=ieast
369 DO j=jstr,jend
370 ad_ubar_obc(j,ib,ir,lout)=-tl_ubar_obc(j,ib,ir,linp)+ &
371 & ad_ubar_obc(j,ib,ir,lout )
372 END DO
373 END IF
374 IF ((lobc(isouth,isubar,ng)).and. &
375 & domain(ng)%Southern_Edge(tile)) THEN
376 ib=isouth
377 DO i=istru,iend
378 ad_ubar_obc(i,ib,ir,lout)=-tl_ubar_obc(i,ib,ir,linp)+ &
379 & ad_ubar_obc(i,ib,ir,lout )
380 END DO
381 END IF
382 IF ((lobc(inorth,isubar,ng)).and. &
383 & domain(ng)%Northern_Edge(tile)) THEN
384 ib=inorth
385 DO i=istru,iend
386 ad_ubar_obc(i,ib,ir,lout)=-tl_ubar_obc(i,ib,ir,linp)+ &
387 & ad_ubar_obc(i,ib,ir,lout )
388 END DO
389 END IF
390 END DO
391 END IF
392# endif
393
394# ifndef SOLVE3D
395!
396! 2D V-momentum component.
397!
398 DO j=jstr,jendr
399 DO i=istrr,iendr
400 ad_vbar(i,j,lout)=-tl_vbar(i,j,linp)+ &
401 & ad_vbar(i,j,lout )
402 END DO
403 END DO
404# endif
405
406# ifdef ADJUST_BOUNDARY
407!
408! 2D V-momentum open boundaries.
409!
410 IF (any(lobc(:,isvbar,ng))) THEN
411 DO ir=1,nbrec(ng)
412 IF ((lobc(iwest,isvbar,ng)).and. &
413 & domain(ng)%Western_Edge(tile)) THEN
414 ib=iwest
415 DO j=jstrv,jend
416 ad_vbar_obc(j,ib,ir,lout)=-tl_vbar_obc(j,ib,ir,linp)+ &
417 & ad_vbar_obc(j,ib,ir,lout )
418 END DO
419 END IF
420 IF ((lobc(ieast,isvbar,ng)).and. &
421 & domain(ng)%Eastern_Edge(tile)) THEN
422 ib=ieast
423 DO j=jstrv,jend
424 ad_vbar_obc(j,ib,ir,lout)=-tl_vbar_obc(j,ib,ir,linp)+ &
425 & ad_vbar_obc(j,ib,ir,lout )
426 END DO
427 END IF
428 IF ((lobc(isouth,isvbar,ng)).and. &
429 & domain(ng)%Southern_Edge(tile)) THEN
430 ib=isouth
431 DO i=istr,iend
432 ad_vbar_obc(i,ib,ir,lout)=-tl_vbar_obc(i,ib,ir,linp)+ &
433 & ad_vbar_obc(i,ib,ir,lout )
434 END DO
435 END IF
436 IF ((lobc(inorth,isvbar,ng)).and. &
437 & domain(ng)%Northern_Edge(tile)) THEN
438 ib=inorth
439 DO i=istr,iend
440 ad_vbar_obc(i,ib,ir,lout)=-tl_vbar_obc(i,ib,ir,linp)+ &
441 & ad_vbar_obc(i,ib,ir,lout )
442 END DO
443 END IF
444 END DO
445 END IF
446# endif
447
448# ifdef ADJUST_WSTRESS
449!
450! Surface momentum stress.
451!
452 DO k=1,nfrec(ng)
453 DO j=jstrr,jendr
454 DO i=istr,iendr
455 ad_ustr(i,j,k,lout)=-tl_ustr(i,j,k,linp)+ &
456 & ad_ustr(i,j,k,lout )
457 END DO
458 END DO
459 DO j=jstr,jendr
460 DO i=istrr,iendr
461 ad_vstr(i,j,k,lout)=-tl_vstr(i,j,k,linp)+ &
462 & ad_vstr(i,j,k,lout )
463 END DO
464 END DO
465 END DO
466# endif
467
468# ifdef SOLVE3D
469!
470! 3D U-momentum component.
471!
472 DO k=1,n(ng)
473 DO j=jstrr,jendr
474 DO i=istr,iendr
475 ad_u(i,j,k,lout)=-tl_u(i,j,k,linp)+ &
476 & ad_u(i,j,k,lout )
477 END DO
478 END DO
479 END DO
480
481# ifdef ADJUST_BOUNDARY
482!
483! 3D U-momentum open boundaries.
484!
485 IF (any(lobc(:,isuvel,ng))) THEN
486 DO ir=1,nbrec(ng)
487 IF ((lobc(iwest,isuvel,ng)).and. &
488 & domain(ng)%Western_Edge(tile)) THEN
489 ib=iwest
490 DO k=1,n(ng)
491 DO j=jstr,jend
492 ad_u_obc(j,k,ib,ir,lout)=-tl_u_obc(j,k,ib,ir,linp)+ &
493 & ad_u_obc(j,k,ib,ir,lout )
494 END DO
495 END DO
496 END IF
497 IF ((lobc(ieast,isuvel,ng)).and. &
498 & domain(ng)%Eastern_Edge(tile)) THEN
499 ib=ieast
500 DO k=1,n(ng)
501 DO j=jstr,jend
502 ad_u_obc(j,k,ib,ir,lout)=-tl_u_obc(j,k,ib,ir,linp)+ &
503 & ad_u_obc(j,k,ib,ir,lout )
504 END DO
505 END DO
506 END IF
507 IF ((lobc(isouth,isuvel,ng)).and. &
508 & domain(ng)%Southern_Edge(tile)) THEN
509 ib=isouth
510 DO k=1,n(ng)
511 DO i=istru,iend
512 ad_u_obc(i,k,ib,ir,lout)=-tl_u_obc(i,k,ib,ir,linp)+ &
513 & ad_u_obc(i,k,ib,ir,lout )
514 END DO
515 END DO
516 END IF
517 IF ((lobc(inorth,isuvel,ng)).and. &
518 & domain(ng)%Northern_Edge(tile)) THEN
519 ib=inorth
520 DO k=1,n(ng)
521 DO i=istru,iend
522 ad_u_obc(i,k,ib,ir,lout)=-tl_u_obc(i,k,ib,ir,linp)+ &
523 & ad_u_obc(i,k,ib,ir,lout )
524 END DO
525 END DO
526 END IF
527 END DO
528 END IF
529# endif
530!
531! 3D V-momentum component.
532!
533 DO k=1,n(ng)
534 DO j=jstr,jendr
535 DO i=istrr,iendr
536 ad_v(i,j,k,lout)=-tl_v(i,j,k,linp)+ &
537 & ad_v(i,j,k,lout )
538 END DO
539 END DO
540 END DO
541
542# ifdef ADJUST_BOUNDARY
543!
544! 3D V-momentum open boundaries.
545!
546 IF (any(lobc(:,isvvel,ng))) THEN
547 DO ir=1,nbrec(ng)
548 IF ((lobc(iwest,isvvel,ng)).and. &
549 & domain(ng)%Western_Edge(tile)) THEN
550 ib=iwest
551 DO k=1,n(ng)
552 DO j=jstrv,jend
553 ad_v_obc(j,k,ib,ir,lout)=-tl_v_obc(j,k,ib,ir,linp)+ &
554 & ad_v_obc(j,k,ib,ir,lout )
555 END DO
556 END DO
557 END IF
558 IF ((lobc(ieast,isvvel,ng)).and. &
559 & domain(ng)%Eastern_Edge(tile)) THEN
560 ib=ieast
561 DO k=1,n(ng)
562 DO j=jstrv,jend
563 ad_v_obc(j,k,ib,ir,lout)=-tl_v_obc(j,k,ib,ir,linp)+ &
564 & ad_v_obc(j,k,ib,ir,lout )
565 END DO
566 END DO
567 END IF
568 IF ((lobc(isouth,isvvel,ng)).and. &
569 & domain(ng)%Southern_Edge(tile)) THEN
570 ib=isouth
571 DO k=1,n(ng)
572 DO i=istr,iend
573 ad_v_obc(i,k,ib,ir,lout)=-tl_v_obc(i,k,ib,ir,linp)+ &
574 & ad_v_obc(i,k,ib,ir,lout )
575 END DO
576 END DO
577 END IF
578 IF ((lobc(inorth,isvvel,ng)).and. &
579 & domain(ng)%Northern_Edge(tile)) THEN
580 ib=inorth
581 DO k=1,n(ng)
582 DO i=istr,iend
583 ad_v_obc(i,k,ib,ir,lout)=-tl_v_obc(i,k,ib,ir,linp)+ &
584 & ad_v_obc(i,k,ib,ir,lout )
585 END DO
586 END DO
587 END IF
588 END DO
589 END IF
590# endif
591!
592! Tracers.
593!
594 DO itrc=1,nt(ng)
595 DO k=1,n(ng)
596 DO j=jstrr,jendr
597 DO i=istrr,iendr
598 ad_t(i,j,k,lout,itrc)=-tl_t(i,j,k,linp,itrc)+ &
599 & ad_t(i,j,k,lout ,itrc)
600 END DO
601 END DO
602 END DO
603 END DO
604
605# ifdef ADJUST_BOUNDARY
606!
607! Tracers open boundaries.
608!
609 DO itrc=1,nt(ng)
610 IF (any(lobc(:,istvar(itrc),ng))) THEN
611 DO ir=1,nbrec(ng)
612 IF ((lobc(iwest,istvar(itrc),ng)).and. &
613 & domain(ng)%Western_Edge(tile)) THEN
614 ib=iwest
615 DO k=1,n(ng)
616 DO j=jstr,jend
617 ad_t_obc(j,k,ib,ir,lout,itrc)= &
618 & -tl_t_obc(j,k,ib,ir,linp,itrc)+ &
619 & ad_t_obc(j,k,ib,ir,lout ,itrc)
620 END DO
621 END DO
622 END IF
623 IF ((lobc(ieast,istvar(itrc),ng)).and. &
624 & domain(ng)%Eastern_Edge(tile)) THEN
625 ib=ieast
626 DO k=1,n(ng)
627 DO j=jstr,jend
628 ad_t_obc(j,k,ib,ir,lout,itrc)= &
629 & -tl_t_obc(j,k,ib,ir,linp,itrc)+ &
630 & ad_t_obc(j,k,ib,ir,lout ,itrc)
631 END DO
632 END DO
633 END IF
634 IF ((lobc(isouth,istvar(itrc),ng)).and. &
635 & domain(ng)%Southern_Edge(tile)) THEN
636 ib=isouth
637 DO k=1,n(ng)
638 DO i=istr,iend
639 ad_t_obc(i,k,ib,ir,lout,itrc)= &
640 & -tl_t_obc(i,k,ib,ir,linp,itrc)+ &
641 & ad_t_obc(i,k,ib,ir,lout ,itrc)
642 END DO
643 END DO
644 END IF
645 IF ((lobc(inorth,istvar(itrc),ng)).and. &
646 & domain(ng)%Northern_Edge(tile)) THEN
647 ib=inorth
648 DO k=1,n(ng)
649 DO i=istr,iend
650 ad_t_obc(i,k,ib,ir,lout,itrc)= &
651 & -tl_t_obc(i,k,ib,ir,linp,itrc)+ &
652 & ad_t_obc(i,k,ib,ir,lout ,itrc)
653 END DO
654 END DO
655 END IF
656 END DO
657 END IF
658 END DO
659# endif
660# ifdef ADJUST_STFLUX
661!
662! Surface tracers flux.
663!
664 DO itrc=1,nt(ng)
665 IF (lstflux(itrc,ng)) THEN
666 DO k=1,nfrec(ng)
667 DO j=jstrr,jendr
668 DO i=istrr,iendr
669 ad_tflux(i,j,k,lout,itrc)=-tl_tflux(i,j,k,linp,itrc)+ &
670 & ad_tflux(i,j,k,lout ,itrc)
671 END DO
672 END DO
673 END DO
674 END IF
675 END DO
676# endif
677# endif
678
679 RETURN
integer isvvel
integer isvbar
integer, dimension(:), allocatable istvar
integer isuvel
integer isfsur
integer isubar
integer, dimension(:), allocatable n
Definition mod_param.F:479
type(t_domain), dimension(:), allocatable domain
Definition mod_param.F:329
integer, dimension(:), allocatable nt
Definition mod_param.F:489
logical, dimension(:,:,:), allocatable lobc
integer, parameter iwest
logical, dimension(:,:), allocatable lstflux
integer, dimension(:), allocatable nfrec
integer, parameter isouth
integer, parameter ieast
integer, parameter inorth
integer, dimension(:), allocatable nbrec

References mod_param::domain, mod_scalars::ieast, mod_scalars::inorth, mod_ncparam::isfsur, mod_scalars::isouth, mod_ncparam::istvar, mod_ncparam::isubar, mod_ncparam::isuvel, mod_ncparam::isvbar, mod_ncparam::isvvel, mod_scalars::iwest, mod_scalars::lobc, and mod_scalars::lstflux.

Referenced by lanc_resid().

Here is the caller graph for this function: