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

Functions/Subroutines

subroutine, public cost_grad (ng, tile, linp1, linp2, lout)
 
subroutine cost_grad_tile (ng, tile, lbi, ubi, lbj, ubj, lbij, ubij, imins, imaxs, jmins, jmaxs, linp1, linp2, 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

◆ cost_grad()

subroutine, public cost_grad_mod::cost_grad ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) linp1,
integer, intent(in) linp2,
integer, intent(in) lout )

Definition at line 42 of file cost_grad.F.

43!***********************************************************************
44!
45 USE mod_param
46# ifdef ADJUST_BOUNDARY
47 USE mod_boundary
48# endif
49# if defined ADJUST_STFLUX || defined ADJUST_WSTRESS
50 USE mod_forces
51# endif
52 USE mod_ocean
53!
54! Imported variable declarations.
55!
56 integer, intent(in) :: ng, tile, Linp1, Linp2, Lout
57!
58! Local variable declarations.
59!
60# include "tile.h"
61!
62 CALL cost_grad_tile (ng, tile, &
63 & lbi, ubi, lbj, ubj, lbij, ubij, &
64 & imins, imaxs, jmins, jmaxs, &
65 & linp1, linp2, lout, &
66# ifdef ADJUST_BOUNDARY
67# ifdef SOLVE3D
68 & boundary(ng) % tl_t_obc, &
69 & boundary(ng) % tl_u_obc, &
70 & boundary(ng) % tl_v_obc, &
71# endif
72 & boundary(ng) % tl_ubar_obc, &
73 & boundary(ng) % tl_vbar_obc, &
74 & boundary(ng) % tl_zeta_obc, &
75# endif
76# ifdef ADJUST_WSTRESS
77 & forces(ng) % tl_ustr, &
78 & forces(ng) % tl_vstr, &
79# endif
80# ifdef SOLVE3D
81# ifdef ADJUST_STFLUX
82 & forces(ng) % tl_tflux, &
83# endif
84 & ocean(ng) % tl_t, &
85 & ocean(ng) % tl_u, &
86 & ocean(ng) % tl_v, &
87# else
88 & ocean(ng) % tl_ubar, &
89 & ocean(ng) % tl_vbar, &
90# endif
91 & ocean(ng) % tl_zeta, &
92# ifdef ADJUST_BOUNDARY
93# ifdef SOLVE3D
94 & boundary(ng) % ad_t_obc, &
95 & boundary(ng) % ad_u_obc, &
96 & boundary(ng) % ad_v_obc, &
97# endif
98 & boundary(ng) % ad_ubar_obc, &
99 & boundary(ng) % ad_vbar_obc, &
100 & boundary(ng) % ad_zeta_obc, &
101# endif
102# ifdef ADJUST_WSTRESS
103 & forces(ng) % ad_ustr, &
104 & forces(ng) % ad_vstr, &
105# endif
106# ifdef SOLVE3D
107# ifdef ADJUST_STFLUX
108 & forces(ng) % ad_tflux, &
109# endif
110 & ocean(ng) % ad_t, &
111 & ocean(ng) % ad_u, &
112 & ocean(ng) % ad_v, &
113# else
114 & ocean(ng) % ad_ubar, &
115 & ocean(ng) % ad_vbar, &
116# endif
117 & ocean(ng) % ad_zeta)
118 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, cost_grad_tile(), mod_forces::forces, and mod_ocean::ocean.

Referenced by i4dvar_mod::increment().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ cost_grad_tile()

subroutine cost_grad_mod::cost_grad_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) linp1,
integer, intent(in) linp2,
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 122 of file cost_grad.F.

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

Here is the caller graph for this function: