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

Functions/Subroutines

subroutine, public tl_variability (ng, tile, linp, lweak)
 
subroutine tl_variability_tile (ng, tile, lbi, ubi, lbj, ubj, lbij, ubij, imins, imaxs, jmins, jmaxs, linp, lweak, t_obc_std, u_obc_std, v_obc_std, ubar_obc_std, vbar_obc_std, zeta_obc_std, sustr_std, svstr_std, stflx_std, t_std, u_std, v_std, ubar_std, vbar_std, zeta_std, 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_ubar, tl_vbar, tl_zeta)
 

Function/Subroutine Documentation

◆ tl_variability()

subroutine, public tl_variability_mod::tl_variability ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) linp,
logical, intent(in) lweak )

Definition at line 56 of file tl_variability.F.

57!***********************************************************************
58!
59 USE mod_param
60# ifdef ADJUST_BOUNDARY
61 USE mod_boundary
62# endif
63# if defined ADJUST_STFLUX || defined ADJUST_WSTRESS
64 USE mod_forces
65# endif
66 USE mod_ocean
67!
68! Imported variable declarations.
69!
70 logical, intent(in) :: Lweak
71
72 integer, intent(in) :: ng, tile, Linp
73!
74! Local variable declarations.
75!
76# include "tile.h"
77!
78 CALL tl_variability_tile (ng, tile, &
79 & lbi, ubi, lbj, ubj, lbij, ubij, &
80 & imins, imaxs, jmins, jmaxs, &
81 & linp, lweak, &
82# ifdef ADJUST_BOUNDARY
83# ifdef SOLVE3D
84 & boundary(ng) % e_t_obc, &
85 & boundary(ng) % e_u_obc, &
86 & boundary(ng) % e_v_obc, &
87# endif
88 & boundary(ng) % e_ubar_obc, &
89 & boundary(ng) % e_vbar_obc, &
90 & boundary(ng) % e_zeta_obc, &
91# endif
92# ifdef ADJUST_WSTRESS
93 & forces(ng) % e_sustr, &
94 & forces(ng) % e_svstr, &
95# endif
96# if defined ADJUST_STFLUX && defined SOLVE3D
97 & forces(ng) % e_stflx, &
98# endif
99# ifdef SOLVE3D
100 & ocean(ng) % e_t, &
101 & ocean(ng) % e_u, &
102 & ocean(ng) % e_v, &
103# if defined WEAK_CONSTRAINT && defined TIME_CONV
104 & ocean(ng) % e_ubar, &
105 & ocean(ng) % e_vbar, &
106# endif
107# else
108 & ocean(ng) % e_ubar, &
109 & ocean(ng) % e_vbar, &
110# endif
111 & ocean(ng) % e_zeta, &
112# ifdef ADJUST_BOUNDARY
113# ifdef SOLVE3D
114 & boundary(ng) % tl_t_obc, &
115 & boundary(ng) % tl_u_obc, &
116 & boundary(ng) % tl_v_obc, &
117# endif
118 & boundary(ng) % tl_ubar_obc, &
119 & boundary(ng) % tl_vbar_obc, &
120 & boundary(ng) % tl_zeta_obc, &
121# endif
122# ifdef ADJUST_WSTRESS
123 & forces(ng) % tl_ustr, &
124 & forces(ng) % tl_vstr, &
125# endif
126# if defined ADJUST_STFLUX && defined SOLVE3D
127 & forces(ng) % tl_tflux, &
128# endif
129# ifdef SOLVE3D
130 & ocean(ng) % tl_t, &
131 & ocean(ng) % tl_u, &
132 & ocean(ng) % tl_v, &
133# if defined WEAK_CONSTRAINT && defined TIME_CONV
134 & ocean(ng) % tl_ubar, &
135 & ocean(ng) % tl_vbar, &
136# endif
137# else
138 & ocean(ng) % tl_ubar, &
139 & ocean(ng) % tl_vbar, &
140# endif
141 & ocean(ng) % tl_zeta)
142
143 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, mod_ocean::ocean, and tl_variability_tile().

Referenced by i4dvar_mod::analysis(), convolve_mod::convolve(), convolve_mod::error_covariance(), i4dvar_mod::increment(), and convolve_mod::saddlec().

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

◆ tl_variability_tile()

subroutine tl_variability_mod::tl_variability_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,
logical, intent(in) lweak,
real(r8), dimension(lbij:,:,:,:), intent(in) t_obc_std,
real(r8), dimension(lbij:,:,:), intent(in) u_obc_std,
real(r8), dimension(lbij:,:,:), intent(in) v_obc_std,
real(r8), dimension(lbij:,:), intent(in) ubar_obc_std,
real(r8), dimension(lbij:,:), intent(in) vbar_obc_std,
real(r8), dimension(lbij:,:), intent(in) zeta_obc_std,
real(r8), dimension(lbi:,lbj:), intent(in) sustr_std,
real(r8), dimension(lbi:,lbj:), intent(in) svstr_std,
real(r8), dimension(lbi:,lbj:,:), intent(in) stflx_std,
real(r8), dimension(lbi:,lbj:,:,:,:), intent(in) t_std,
real(r8), dimension(lbi:,lbj:,:,:), intent(in) u_std,
real(r8), dimension(lbi:,lbj:,:,:), intent(in) v_std,
real(r8), dimension(lbi:,lbj:,:), intent(in) ubar_std,
real(r8), dimension(lbi:,lbj:,:), intent(in) vbar_std,
real(r8), dimension(lbi:,lbj:,:), intent(in) zeta_std,
real(r8), dimension(lbij:,:,:,:,:,:), intent(inout) tl_t_obc,
real(r8), dimension(lbij:,:,:,:,:), intent(inout) tl_u_obc,
real(r8), dimension(lbij:,:,:,:,:), intent(inout) tl_v_obc,
real(r8), dimension(lbij:,:,:,:), intent(inout) tl_ubar_obc,
real(r8), dimension(lbij:,:,:,:), intent(inout) tl_vbar_obc,
real(r8), dimension(lbij:,:,:,:), intent(inout) tl_zeta_obc,
real(r8), dimension(lbi:,lbj:,:,:), intent(inout) tl_ustr,
real(r8), dimension(lbi:,lbj:,:,:), intent(inout) tl_vstr,
real(r8), dimension(lbi:,lbj:,:,:,:), intent(inout) tl_tflux,
real(r8), dimension(lbi:,lbj:,:,:,:), intent(inout) tl_t,
real(r8), dimension(lbi:,lbj:,:,:), intent(inout) tl_u,
real(r8), dimension(lbi:,lbj:,:,:), intent(inout) tl_v,
real(r8), dimension(lbi:,lbj:,:), intent(inout) tl_ubar,
real(r8), dimension(lbi:,lbj:,:), intent(inout) tl_vbar,
real(r8), dimension(lbi:,lbj:,:), intent(inout) tl_zeta )
private

Definition at line 147 of file tl_variability.F.

195!***********************************************************************
196!
197 USE mod_param
198 USE mod_ncparam
199 USE mod_scalars
200
201# ifdef DISTRIBUTE
202!
204# endif
205!
206! Imported variable declarations.
207!
208 logical, intent(in) :: Lweak
209
210 integer, intent(in) :: ng, tile
211 integer, intent(in) :: LBi, UBi, LBj, UBj, LBij, UBij
212 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
213 integer, intent(in) :: Linp
214!
215# ifdef ASSUMED_SHAPE
216# ifdef ADJUST_BOUNDARY
217# ifdef SOLVE3D
218 real(r8), intent(in) :: t_obc_std(LBij:,:,:,:)
219 real(r8), intent(in) :: u_obc_std(LBij:,:,:)
220 real(r8), intent(in) :: v_obc_std(LBij:,:,:)
221# endif
222 real(r8), intent(in) :: ubar_obc_std(LBij:,:)
223 real(r8), intent(in) :: vbar_obc_std(LBij:,:)
224 real(r8), intent(in) :: zeta_obc_std(LBij:,:)
225# endif
226# ifdef ADJUST_WSTRESS
227 real(r8), intent(in) :: sustr_std(LBi:,LBj:)
228 real(r8), intent(in) :: svstr_std(LBi:,LBj:)
229# endif
230# if defined ADJUST_STFLUX && defined SOLVE3D
231 real(r8), intent(in) :: stflx_std(LBi:,LBj:,:)
232# endif
233# ifdef SOLVE3D
234 real(r8), intent(in) :: t_std(LBi:,LBj:,:,:,:)
235 real(r8), intent(in) :: u_std(LBi:,LBj:,:,:)
236 real(r8), intent(in) :: v_std(LBi:,LBj:,:,:)
237# if defined WEAK_CONSTRAINT && defined TIME_CONV
238 real(r8), intent(in) :: ubar_std(LBi:,LBj:,:)
239 real(r8), intent(in) :: vbar_std(LBi:,LBj:,:)
240# endif
241# else
242 real(r8), intent(in) :: ubar_std(LBi:,LBj:,:)
243 real(r8), intent(in) :: vbar_std(LBi:,LBj:,:)
244# endif
245 real(r8), intent(in) :: zeta_std(LBi:,LBj:,:)
246# ifdef ADJUST_BOUNDARY
247# ifdef SOLVE3D
248 real(r8), intent(inout) :: tl_t_obc(LBij:,:,:,:,:,:)
249 real(r8), intent(inout) :: tl_u_obc(LBij:,:,:,:,:)
250 real(r8), intent(inout) :: tl_v_obc(LBij:,:,:,:,:)
251# endif
252 real(r8), intent(inout) :: tl_ubar_obc(LBij:,:,:,:)
253 real(r8), intent(inout) :: tl_vbar_obc(LBij:,:,:,:)
254 real(r8), intent(inout) :: tl_zeta_obc(LBij:,:,:,:)
255# endif
256# ifdef ADJUST_WSTRESS
257 real(r8), intent(inout) :: tl_ustr(LBi:,LBj:,:,:)
258 real(r8), intent(inout) :: tl_vstr(LBi:,LBj:,:,:)
259# endif
260# if defined ADJUST_STFLUX && defined SOLVE3D
261 real(r8), intent(inout) :: tl_tflux(LBi:,LBj:,:,:,:)
262# endif
263# ifdef SOLVE3D
264 real(r8), intent(inout) :: tl_t(LBi:,LBj:,:,:,:)
265 real(r8), intent(inout) :: tl_u(LBi:,LBj:,:,:)
266 real(r8), intent(inout) :: tl_v(LBi:,LBj:,:,:)
267# if defined WEAK_CONSTRAINT && defined TIME_CONV
268 real(r8), intent(inout) :: tl_ubar(LBi:,LBj:,:)
269 real(r8), intent(inout) :: tl_vbar(LBi:,LBj:,:)
270# endif
271# else
272 real(r8), intent(inout) :: tl_ubar(LBi:,LBj:,:)
273 real(r8), intent(inout) :: tl_vbar(LBi:,LBj:,:)
274# endif
275 real(r8), intent(inout) :: tl_zeta(LBi:,LBj:,:)
276
277# else
278
279# ifdef ADJUST_BOUNDARY
280# ifdef SOLVE3D
281 real(r8), intent(in) :: t_obc_std(LBij:UBij,N(ng),4,NT(ng))
282 real(r8), intent(in) :: u_obc_std(LBij:UBij,N(ng),4)
283 real(r8), intent(in) :: v_obc_std(LBij:UBij,N(ng),4)
284# endif
285 real(r8), intent(in) :: ubar_obc_std(LBij:UBij,4)
286 real(r8), intent(in) :: vbar_obc_std(LBij:UBij,4)
287 real(r8), intent(in) :: zeta_obc_std(LBij:UBij,4)
288# endif
289# ifdef ADJUST_WSTRESS
290 real(r8), intent(in) :: sustr_std(LBi:,LBj:)
291 real(r8), intent(in) :: svstr_std(LBi:,LBj:)
292# endif
293# if defined ADJUST_STFLUX && defined SOLVE3D
294 real(r8), intent(in) :: stflx_std(LBi:UBi,LBj:UBj,NT(ng))
295# endif
296# ifdef SOLVE3D
297 real(r8), intent(in) :: t_std(LBi:UBi,LBj:UBj,N(ng),NSA,NT(ng))
298 real(r8), intent(in) :: u_std(LBi:UBi,LBj:UBj,N(ng),NSA)
299 real(r8), intent(in) :: v_std(LBi:UBi,LBj:UBj,N(ng),NSA)
300# if defined WEAK_CONSTRAINT && defined TIME_CONV
301 real(r8), intent(in) :: ubar_std(LBi:UBi,LBj:UBj,NSA)
302 real(r8), intent(in) :: vbar_std(LBi:UBi,LBj:UBj,NSA)
303# endif
304# else
305 real(r8), intent(in) :: ubar_std(LBi:UBi,LBj:UBj,NSA)
306 real(r8), intent(in) :: vbar_std(LBi:UBi,LBj:UBj,NSA)
307# endif
308 real(r8), intent(in) :: zeta_std(LBi:UBi,LBj:UBj,NSA)
309# ifdef ADJUST_BOUNDARY
310# ifdef SOLVE3D
311 real(r8), intent(inout) :: tl_t_obc(LBij:UBij,N(ng),4, &
312 & Nbrec(ng),2,NT(ng))
313 real(r8), intent(inout) :: tl_u_obc(LBij:UBij,N(ng),4,Nbrec(ng),2)
314 real(r8), intent(inout) :: tl_v_obc(LBij:UBij,N(ng),4,Nbrec(ng),2)
315# endif
316 real(r8), intent(inout) :: tl_ubar_obc(LBij:UBij,4,Nbrec(ng),2)
317 real(r8), intent(inout) :: tl_vbar_obc(LBij:UBij,4,Nbrec(ng),2)
318 real(r8), intent(inout) :: tl_zeta_obc(LBij:UBij,4,Nbrec(ng),2)
319# endif
320# ifdef ADJUST_WSTRESS
321 real(r8), intent(inout) :: tl_ustr(LBi:UBi,LBj:UBj,Nfrec(ng),2)
322 real(r8), intent(inout) :: tl_vstr(LBi:UBi,LBj:UBj,Nfrec(ng),2)
323# endif
324# if defined ADJUST_STFLUX && defined SOLVE3D
325 real(r8), intent(inout) :: tl_tflux(LBi:UBi,LBj:UBj, &
326 & Nfrec(ng),2,NT(ng))
327# endif
328# ifdef SOLVE3D
329 real(r8), intent(inout) :: tl_t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
330 real(r8), intent(inout) :: tl_u(LBi:UBi,LBj:UBj,N(ng),2)
331 real(r8), intent(inout) :: tl_v(LBi:UBi,LBj:UBj,N(ng),2)
332# if defined WEAK_CONSTRAINT && defined TIME_CONV
333 real(r8), intent(inout) :: tl_ubar(LBi:UBi,LBj:UBj,:)
334 real(r8), intent(inout) :: tl_vbar(LBi:UBi,LBj:UBj,:)
335# endif
336# else
337 real(r8), intent(inout) :: tl_ubar(LBi:UBi,LBj:UBj,:)
338 real(r8), intent(inout) :: tl_vbar(LBi:UBi,LBj:UBj,:)
339# endif
340 real(r8), intent(inout) :: tl_zeta(LBi:UBi,LBj:UBj,:)
341# endif
342!
343! Local variable declarations.
344!
345 integer :: i, ib, ir, j, rec
346# ifdef SOLVE3D
347 integer :: it, k
348# endif
349
350# include "set_bounds.h"
351!
352!-----------------------------------------------------------------------
353! Determine error covariance standard deviation factors to use.
354!-----------------------------------------------------------------------
355!
356 IF (lweak) THEN
357 rec=2 ! weak constraint
358 ELSE
359 rec=1 ! strong constraint
360 END IF
361!
362!-----------------------------------------------------------------------
363! Initial conditions and model error covariance: Multiply tangent
364! linear state by its corresponding standard deviation.
365!-----------------------------------------------------------------------
366!
367! Free-surface.
368!
369 DO j=jstrt,jendt
370 DO i=istrt,iendt
371 tl_zeta(i,j,linp)=tl_zeta(i,j,linp)*zeta_std(i,j,rec)
372 END DO
373 END DO
374# ifdef DISTRIBUTE
375 CALL mp_exchange2d (ng, tile, itlm, 1, &
376 & lbi, ubi, lbj, ubj, &
377 & nghostpoints, &
378 & ewperiodic(ng), nsperiodic(ng), &
379 & tl_zeta(:,:,linp))
380# endif
381# if !defined SOLVE3D || (defined WEAK_CONSTRAINT && defined TIME_CONV)
382!
383! 2D U-momentum.
384!
385 DO j=jstrt,jendt
386 DO i=istrp,iendt
387 tl_ubar(i,j,linp)=tl_ubar(i,j,linp)*ubar_std(i,j,rec)
388 END DO
389 END DO
390 DO j=jstrp,jendt
391 DO i=istrt,iendt
392 tl_vbar(i,j,linp)=tl_vbar(i,j,linp)*vbar_std(i,j,rec)
393 END DO
394 END DO
395# ifdef DISTRIBUTE
396 CALL mp_exchange2d (ng, tile, itlm, 2, &
397 & lbi, ubi, lbj, ubj, &
398 & nghostpoints, &
399 & ewperiodic(ng), nsperiodic(ng), &
400 & tl_ubar(:,:,linp), &
401 & tl_vbar(:,:,linp))
402# endif
403# endif
404# ifdef SOLVE3D
405!
406! 3D U-momentum.
407!
408 DO k=1,n(ng)
409 DO j=jstrt,jendt
410 DO i=istrp,iendt
411 tl_u(i,j,k,linp)=tl_u(i,j,k,linp)*u_std(i,j,k,rec)
412 END DO
413 END DO
414 DO j=jstrp,jendt
415 DO i=istrt,iendt
416 tl_v(i,j,k,linp)=tl_v(i,j,k,linp)*v_std(i,j,k,rec)
417 END DO
418 END DO
419 END DO
420# ifdef DISTRIBUTE
421 CALL mp_exchange3d (ng, tile, itlm, 2, &
422 & lbi, ubi, lbj, ubj, 1, n(ng), &
423 & nghostpoints, &
424 & ewperiodic(ng), nsperiodic(ng), &
425 & tl_u(:,:,:,linp), &
426 & tl_v(:,:,:,linp))
427# endif
428!
429! Tangent linear tracers.
430!
431 DO it=1,nt(ng)
432 DO k=1,n(ng)
433 DO j=jstrt,jendt
434 DO i=istrt,iendt
435 tl_t(i,j,k,linp,it)=tl_t(i,j,k,linp,it)* &
436 & t_std(i,j,k,rec,it)
437 END DO
438 END DO
439 END DO
440 END DO
441# ifdef DISTRIBUTE
442 CALL mp_exchange4d (ng, tile, itlm, 1, &
443 & lbi, ubi, lbj, ubj, 1, n(ng), 1, nt(ng), &
444 & nghostpoints, &
445 & ewperiodic(ng), nsperiodic(ng), &
446 & tl_t(:,:,:,linp,:))
447# endif
448# endif
449
450# ifdef ADJUST_BOUNDARY
451!
452!-----------------------------------------------------------------------
453! Boundary conditions and model error covariance: Multiply tangent
454! linear state by its corresponding standard deviation.
455!-----------------------------------------------------------------------
456!
457! Free-surface open boundaries.
458!
459 IF (.not.lweak.and.(any(lobc(:,isfsur,ng)))) THEN
460 DO ir=1,nbrec(ng)
461 IF ((lobc(iwest,isfsur,ng)).and. &
462 & domain(ng)%Western_Edge(tile)) THEN
463 ib=iwest
464 DO j=jstr,jend
465 tl_zeta_obc(j,ib,ir,linp)=tl_zeta_obc(j,ib,ir,linp)* &
466 & zeta_obc_std(j,ib)
467 END DO
468 END IF
469 IF ((lobc(ieast,isfsur,ng)).and. &
470 & domain(ng)%Eastern_Edge(tile)) THEN
471 ib=ieast
472 DO j=jstr,jend
473 tl_zeta_obc(j,ib,ir,linp)=tl_zeta_obc(j,ib,ir,linp)* &
474 & zeta_obc_std(j,ib)
475 END DO
476 END IF
477 IF ((lobc(isouth,isfsur,ng)).and. &
478 & domain(ng)%Southern_Edge(tile)) THEN
479 ib=isouth
480 DO i=istr,iend
481 tl_zeta_obc(i,ib,ir,linp)=tl_zeta_obc(i,ib,ir,linp)* &
482 & zeta_obc_std(i,ib)
483 END DO
484 END IF
485 IF ((lobc(inorth,isfsur,ng)).and. &
486 & domain(ng)%Northern_Edge(tile)) THEN
487 ib=inorth
488 DO i=istr,iend
489 tl_zeta_obc(i,ib,ir,linp)=tl_zeta_obc(i,ib,ir,linp)* &
490 & zeta_obc_std(i,ib)
491 END DO
492 END IF
493# ifdef DISTRIBUTE
494 DO ib=1,4
495 IF (lobc(ib,isfsur,ng)) THEN
496 CALL mp_exchange2d_bry (ng, tile, itlm, 1, ib, &
497 & lbij, ubij, &
498 & nghostpoints, &
499 & ewperiodic(ng), nsperiodic(ng), &
500 & tl_zeta_obc(:,ib,ir,linp))
501 END IF
502 END DO
503# endif
504 END DO
505 END IF
506!
507! 2D U-momentum open boundaries.
508!
509 IF (.not.lweak.and.(any(lobc(:,isubar,ng)))) THEN
510 DO ir=1,nbrec(ng)
511 IF ((lobc(iwest,isubar,ng)).and. &
512 & domain(ng)%Western_Edge(tile)) THEN
513 ib=iwest
514 DO j=jstr,jend
515 tl_ubar_obc(j,ib,ir,linp)=tl_ubar_obc(j,ib,ir,linp)* &
516 & ubar_obc_std(j,ib)
517 END DO
518 END IF
519 IF ((lobc(ieast,isubar,ng)).and. &
520 & domain(ng)%Eastern_Edge(tile)) THEN
521 ib=ieast
522 DO j=jstr,jend
523 tl_ubar_obc(j,ib,ir,linp)=tl_ubar_obc(j,ib,ir,linp)* &
524 & ubar_obc_std(j,ib)
525 END DO
526 END IF
527 IF ((lobc(isouth,isubar,ng)).and. &
528 & domain(ng)%Southern_Edge(tile)) THEN
529 ib=isouth
530 DO i=istru,iend
531 tl_ubar_obc(i,ib,ir,linp)=tl_ubar_obc(i,ib,ir,linp)* &
532 & ubar_obc_std(i,ib)
533 END DO
534 END IF
535 IF ((lobc(inorth,isubar,ng)).and. &
536 & domain(ng)%Northern_Edge(tile)) THEN
537 ib=inorth
538 DO i=istru,iend
539 tl_ubar_obc(i,ib,ir,linp)=tl_ubar_obc(i,ib,ir,linp)* &
540 & ubar_obc_std(i,ib)
541 END DO
542 END IF
543# ifdef DISTRIBUTE
544 DO ib=1,4
545 IF (lobc(ib,isubar,ng)) THEN
546 CALL mp_exchange2d_bry (ng, tile, itlm, 1, ib, &
547 & lbij, ubij, &
548 & nghostpoints, &
549 & ewperiodic(ng), nsperiodic(ng), &
550 & tl_ubar_obc(:,ib,ir,linp))
551 END IF
552 END DO
553# endif
554 END DO
555 END IF
556!
557! 2D V-momentum open boundaries.
558!
559 IF (.not.lweak.and.(any(lobc(:,isvbar,ng)))) THEN
560 DO ir=1,nbrec(ng)
561 IF ((lobc(iwest,isvbar,ng)).and. &
562 & domain(ng)%Western_Edge(tile)) THEN
563 ib=iwest
564 DO j=jstrv,jend
565 tl_vbar_obc(j,ib,ir,linp)=tl_vbar_obc(j,ib,ir,linp)* &
566 & vbar_obc_std(j,ib)
567 END DO
568 END IF
569 IF ((lobc(ieast,isvbar,ng)).and. &
570 & domain(ng)%Eastern_Edge(tile)) THEN
571 ib=ieast
572 DO j=jstrv,jend
573 tl_vbar_obc(j,ib,ir,linp)=tl_vbar_obc(j,ib,ir,linp)* &
574 & vbar_obc_std(j,ib)
575 END DO
576 END IF
577 IF ((lobc(isouth,isvbar,ng)).and. &
578 & domain(ng)%Southern_Edge(tile)) THEN
579 ib=isouth
580 DO i=istr,iend
581 tl_vbar_obc(i,ib,ir,linp)=tl_vbar_obc(i,ib,ir,linp)* &
582 & vbar_obc_std(i,ib)
583 END DO
584 END IF
585 IF ((lobc(inorth,isvbar,ng)).and. &
586 & domain(ng)%Northern_Edge(tile)) THEN
587 ib=inorth
588 DO i=istr,iend
589 tl_vbar_obc(i,ib,ir,linp)=tl_vbar_obc(i,ib,ir,linp)* &
590 & vbar_obc_std(i,ib)
591 END DO
592 END IF
593# ifdef DISTRIBUTE
594 DO ib=1,4
595 IF (lobc(ib,isvbar,ng)) THEN
596 CALL mp_exchange2d_bry (ng, tile, itlm, 1, ib, &
597 & lbij, ubij, &
598 & nghostpoints, &
599 & ewperiodic(ng), nsperiodic(ng), &
600 & tl_vbar_obc(:,ib,ir,linp))
601 END IF
602 END DO
603# endif
604 END DO
605 END IF
606
607# ifdef SOLVE3D
608!
609! 3D U-momentum open boundaries.
610!
611 IF (.not.lweak.and.(any(lobc(:,isuvel,ng)))) THEN
612 DO ir=1,nbrec(ng)
613 IF ((lobc(iwest,isuvel,ng)).and. &
614 & domain(ng)%Western_Edge(tile)) THEN
615 ib=iwest
616 DO k=1,n(ng)
617 DO j=jstr,jend
618 tl_u_obc(j,k,ib,ir,linp)=tl_u_obc(j,k,ib,ir,linp)* &
619 & u_obc_std(j,k,ib)
620 END DO
621 END DO
622 END IF
623 IF ((lobc(ieast,isuvel,ng)).and. &
624 & domain(ng)%Eastern_Edge(tile)) THEN
625 ib=ieast
626 DO k=1,n(ng)
627 DO j=jstr,jend
628 tl_u_obc(j,k,ib,ir,linp)=tl_u_obc(j,k,ib,ir,linp)* &
629 & u_obc_std(j,k,ib)
630 END DO
631 END DO
632 END IF
633 IF ((lobc(isouth,isuvel,ng)).and. &
634 & domain(ng)%Southern_Edge(tile)) THEN
635 ib=isouth
636 DO k=1,n(ng)
637 DO i=istru,iend
638 tl_u_obc(i,k,ib,ir,linp)=tl_u_obc(i,k,ib,ir,linp)* &
639 & u_obc_std(i,k,ib)
640 END DO
641 END DO
642 END IF
643 IF ((lobc(inorth,isuvel,ng)).and. &
644 & domain(ng)%Northern_Edge(tile)) THEN
645 ib=inorth
646 DO k=1,n(ng)
647 DO i=istru,iend
648 tl_u_obc(i,k,ib,ir,linp)=tl_u_obc(i,k,ib,ir,linp)* &
649 & u_obc_std(i,k,ib)
650 END DO
651 END DO
652 END IF
653# ifdef DISTRIBUTE
654 DO ib=1,4
655 IF (lobc(ib,isuvel,ng)) THEN
656 CALL mp_exchange3d_bry (ng, tile, itlm, 1, ib, &
657 & lbij, ubij, 1, n(ng), &
658 & nghostpoints, &
659 & ewperiodic(ng), nsperiodic(ng), &
660 & tl_u_obc(:,:,ib,ir,linp))
661 END IF
662 END DO
663# endif
664 END DO
665 END IF
666!
667! 3D V-momentum open boundaries.
668!
669 IF (.not.lweak.and.(any(lobc(:,isvvel,ng)))) THEN
670 DO ir=1,nbrec(ng)
671 IF ((lobc(iwest,isvvel,ng)).and. &
672 & domain(ng)%Western_Edge(tile)) THEN
673 ib=iwest
674 DO k=1,n(ng)
675 DO j=jstrv,jend
676 tl_v_obc(j,k,ib,ir,linp)=tl_v_obc(j,k,ib,ir,linp)* &
677 & v_obc_std(j,k,ib)
678 END DO
679 END DO
680 END IF
681 IF ((lobc(ieast,isvvel,ng)).and. &
682 & domain(ng)%Eastern_Edge(tile)) THEN
683 ib=ieast
684 DO k=1,n(ng)
685 DO j=jstrv,jend
686 tl_v_obc(j,k,ib,ir,linp)=tl_v_obc(j,k,ib,ir,linp)* &
687 & v_obc_std(j,k,ib)
688 END DO
689 END DO
690 END IF
691 IF ((lobc(isouth,isvvel,ng)).and. &
692 & domain(ng)%Southern_Edge(tile)) THEN
693 ib=isouth
694 DO k=1,n(ng)
695 DO i=istr,iend
696 tl_v_obc(i,k,ib,ir,linp)=tl_v_obc(i,k,ib,ir,linp)* &
697 & v_obc_std(i,k,ib)
698 END DO
699 END DO
700 END IF
701 IF ((lobc(inorth,isvvel,ng)).and. &
702 & domain(ng)%Northern_Edge(tile)) THEN
703 ib=inorth
704 DO k=1,n(ng)
705 DO i=istr,iend
706 tl_v_obc(i,k,ib,ir,linp)=tl_v_obc(i,k,ib,ir,linp)* &
707 & v_obc_std(i,k,ib)
708 END DO
709 END DO
710 END IF
711# ifdef DISTRIBUTE
712 DO ib=1,4
713 IF (lobc(ib,isvvel,ng)) THEN
714 CALL mp_exchange3d_bry (ng, tile, itlm, 1, ib, &
715 & lbij, ubij, 1, n(ng), &
716 & nghostpoints, &
717 & ewperiodic(ng), nsperiodic(ng), &
718 & tl_v_obc(:,:,ib,ir,linp))
719 END IF
720 END DO
721# endif
722 END DO
723 END IF
724!
725! Tracers open boundaries.
726!
727 DO it=1,nt(ng)
728 IF (.not.lweak.and.(any(lobc(:,istvar(it),ng)))) THEN
729 DO ir=1,nbrec(ng)
730 IF ((lobc(iwest,istvar(it),ng)).and. &
731 & domain(ng)%Western_Edge(tile)) THEN
732 ib=iwest
733 DO k=1,n(ng)
734 DO j=jstr,jend
735 tl_t_obc(j,k,ib,ir,linp,it)= &
736 & tl_t_obc(j,k,ib,ir,linp,it)* &
737 & t_obc_std(j,k,ib,it)
738 END DO
739 END DO
740 END IF
741 IF ((lobc(ieast,istvar(it),ng)).and. &
742 & domain(ng)%Eastern_Edge(tile)) THEN
743 ib=ieast
744 DO k=1,n(ng)
745 DO j=jstr,jend
746 tl_t_obc(j,k,ib,ir,linp,it)= &
747 & tl_t_obc(j,k,ib,ir,linp,it)* &
748 & t_obc_std(j,k,ib,it)
749 END DO
750 END DO
751 END IF
752 IF ((lobc(isouth,istvar(it),ng)).and. &
753 & domain(ng)%Southern_Edge(tile)) THEN
754 ib=isouth
755 DO k=1,n(ng)
756 DO i=istr,iend
757 tl_t_obc(i,k,ib,ir,linp,it)= &
758 & tl_t_obc(i,k,ib,ir,linp,it)* &
759 & t_obc_std(i,k,ib,it)
760 END DO
761 END DO
762 END IF
763 IF ((lobc(inorth,istvar(it),ng)).and. &
764 & domain(ng)%Northern_Edge(tile)) THEN
765 ib=inorth
766 DO k=1,n(ng)
767 DO i=istr,iend
768 tl_t_obc(i,k,ib,ir,linp,it)= &
769 & tl_t_obc(i,k,ib,ir,linp,it)* &
770 & t_obc_std(i,k,ib,it)
771 END DO
772 END DO
773 END IF
774# ifdef DISTRIBUTE
775 DO ib=1,4
776 IF (lobc(ib,istvar(it),ng)) THEN
777 CALL mp_exchange3d_bry (ng, tile, itlm, 1, ib, &
778 & lbij, ubij, 1, n(ng), &
779 & nghostpoints, &
780 & ewperiodic(ng), nsperiodic(ng), &
781 & tl_t_obc(:,:,ib,ir,linp,it))
782 END IF
783 END DO
784# endif
785 END DO
786 END IF
787 END DO
788# endif
789# endif
790
791# if defined ADJUST_WSTRESS || defined ADJUST_STFLUX
792!
793!-----------------------------------------------------------------------
794! Initial conditions and model error covariance: Multiply tangent
795! linear state by its corresponding standard deviation.
796!-----------------------------------------------------------------------
797
798# ifdef ADJUST_WSTRESS
799!
800! Surface momentum stress.
801!
802 IF (.not.lweak) THEN
803 DO ir=1,nfrec(ng)
804 DO j=jstrt,jendt
805 DO i=istrp,iendt
806 tl_ustr(i,j,ir,linp)=tl_ustr(i,j,ir,linp)*sustr_std(i,j)
807 END DO
808 END DO
809 DO j=jstrp,jendt
810 DO i=istrt,iendt
811 tl_vstr(i,j,ir,linp)=tl_vstr(i,j,ir,linp)*svstr_std(i,j)
812 END DO
813 END DO
814 END DO
815# ifdef DISTRIBUTE
816 CALL mp_exchange3d (ng, tile, itlm, 2, &
817 & lbi, ubi, lbj, ubj, 1, nfrec(ng), &
818 & nghostpoints, &
819 & ewperiodic(ng), nsperiodic(ng), &
820 & tl_ustr(:,:,:,linp), &
821 & tl_vstr(:,:,:,linp))
822# endif
823 END IF
824# endif
825# if defined ADJUST_STFLUX && defined SOLVE3D
826!
827! Tangent linear surface tracers flux.
828!
829 IF (.not.lweak) THEN
830 DO it=1,nt(ng)
831 IF (lstflux(it,ng)) THEN
832 DO ir=1,nfrec(ng)
833 DO j=jstrt,jendt
834 DO i=istrt,iendt
835 tl_tflux(i,j,ir,linp,it)=tl_tflux(i,j,ir,linp,it)* &
836 & stflx_std(i,j,it)
837 END DO
838 END DO
839 END DO
840# ifdef DISTRIBUTE
841 CALL mp_exchange3d (ng, tile, itlm, 1, &
842 & lbi, ubi, lbj, ubj, 1, nfrec(ng), &
843 & nghostpoints, &
844 & ewperiodic(ng), nsperiodic(ng), &
845 & tl_tflux(:,:,:,linp,it))
846# endif
847 END IF
848 END DO
849 END IF
850# endif
851# endif
852
853 RETURN
integer isvvel
integer isvbar
integer, dimension(:), allocatable istvar
integer isuvel
integer isfsur
integer isubar
integer, dimension(:), allocatable n
Definition mod_param.F:479
integer nghostpoints
Definition mod_param.F:710
type(t_domain), dimension(:), allocatable domain
Definition mod_param.F:329
integer, parameter itlm
Definition mod_param.F:663
integer, dimension(:), allocatable nt
Definition mod_param.F:489
logical, dimension(:,:,:), allocatable lobc
logical, dimension(:), allocatable ewperiodic
integer, parameter iwest
logical, dimension(:), allocatable nsperiodic
logical, dimension(:,:), allocatable lstflux
integer, dimension(:), allocatable nfrec
integer, parameter isouth
integer, parameter ieast
integer, parameter inorth
integer, dimension(:), allocatable nbrec
subroutine mp_exchange2d_bry(ng, tile, model, nvar, boundary, lbij, ubij, nghost, ew_periodic, ns_periodic, a, b, c, d)
subroutine mp_exchange4d(ng, tile, model, nvar, lbi, ubi, lbj, ubj, lbk, ubk, lbt, ubt, nghost, ew_periodic, ns_periodic, a, b, c)
subroutine mp_exchange2d(ng, tile, model, nvar, lbi, ubi, lbj, ubj, nghost, ew_periodic, ns_periodic, a, b, c, d)
subroutine mp_exchange3d(ng, tile, model, nvar, lbi, ubi, lbj, ubj, lbk, ubk, nghost, ew_periodic, ns_periodic, a, b, c, d)
subroutine mp_exchange3d_bry(ng, tile, model, nvar, boundary, lbij, ubij, lbk, ubk, nghost, ew_periodic, ns_periodic, a, b, c, d)

References mod_param::domain, mod_scalars::ewperiodic, 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_param::itlm, mod_scalars::iwest, mod_scalars::lobc, mod_scalars::lstflux, mp_exchange_mod::mp_exchange2d(), mp_exchange_mod::mp_exchange2d_bry(), mp_exchange_mod::mp_exchange3d(), mp_exchange_mod::mp_exchange3d_bry(), mp_exchange_mod::mp_exchange4d(), mod_param::nghostpoints, and mod_scalars::nsperiodic.

Referenced by tl_variability().

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