ROMS
Loading...
Searching...
No Matches
posterior_var.F
Go to the documentation of this file.
1#include "cppdefs.h"
2
4
5#if defined WEAK_CONSTRAINT && \
6 (defined posterior_error_i || defined posterior_error_f)
7
8!
9!git $Id$
10!================================================== Hernan G. Arango ===
11! Copyright (c) 2002-2025 The ROMS Group Andrew M. Moore !
12! Licensed under a MIT/X style license !
13! See License_ROMS.md !
14!=======================================================================
15! !
16! This routine computes the full (diagonal) posterior analysis error !
17! covariance matrix. !
18! !
19!=======================================================================
20!
21 USE mod_param
22 USE mod_parallel
23# ifdef ADJUST_BOUNDARY
24 USE mod_boundary
25# endif
26# ifdef SOLVE3D
27 USE mod_coupling
28# endif
29# if defined ADJUST_STFLUX || defined ADJUST_WSTRESS
30 USE mod_forces
31# endif
32 USE mod_fourdvar
33 USE mod_grid
34 USE mod_iounits
35 USE mod_ncparam
36 USE mod_ocean
37 USE mod_scalars
38 USE mod_stepping
39!
40# ifdef DISTRIBUTE
42 USE distribute_mod, ONLY : mp_reduce
43# endif
45 USE state_copy_mod, ONLY : state_copy
48 USE state_read_mod, ONLY : state_read
49!
50 implicit none
51!
52 PRIVATE
53 PUBLIC :: posterior_var
54!
55 CONTAINS
56!
57!***********************************************************************
58 SUBROUTINE posterior_var (ng, tile, model, outLoop)
59!***********************************************************************
60!
61! Imported variable declarations.
62!
63 integer, intent(in) :: ng, tile, model, outloop
64!
65! Local variable declarations.
66!
67 character (len=*), parameter :: myfile = &
68 & __FILE__
69!
70# include "tile.h"
71!
72# ifdef PROFILE
73 CALL wclock_on (ng, model, 83, __line__, myfile)
74# endif
75!
76 CALL posterior_var_tile (ng, tile, model, &
77 & lbi, ubi, lbj, ubj, lbij, ubij, &
78 & imins, imaxs, jmins, jmaxs, &
79 & lold(ng), lnew(ng), outloop, &
80# ifdef MASKING
81 & grid(ng) % rmask, &
82 & grid(ng) % umask, &
83 & grid(ng) % vmask, &
84# endif
85# ifdef ADJUST_BOUNDARY
86# ifdef SOLVE3D
87 & boundary(ng) % e_t_obc, &
88 & boundary(ng) % e_u_obc, &
89 & boundary(ng) % e_v_obc, &
90# endif
91 & boundary(ng) % e_ubar_obc, &
92 & boundary(ng) % e_vbar_obc, &
93 & boundary(ng) % e_zeta_obc, &
94# endif
95# ifdef ADJUST_WSTRESS
96 & forces(ng) % e_sustr, &
97 & forces(ng) % e_svstr, &
98# endif
99# if defined ADJUST_STFLUX && defined SOLVE3D
100 & forces(ng) % e_stflx, &
101# endif
102# ifdef SOLVE3D
103 & ocean(ng) % e_t, &
104 & ocean(ng) % e_u, &
105 & ocean(ng) % e_v, &
106# if defined WEAK_CONSTRAINT && defined TIME_CONV
107 & ocean(ng) % e_ubar, &
108 & ocean(ng) % e_vbar, &
109# endif
110# else
111 & ocean(ng) % e_ubar, &
112 & ocean(ng) % e_vbar, &
113# endif
114 & ocean(ng) % e_zeta, &
115# ifdef ADJUST_BOUNDARY
116# ifdef SOLVE3D
117 & boundary(ng) % t_obc, &
118 & boundary(ng) % u_obc, &
119 & boundary(ng) % v_obc, &
120# endif
121 & boundary(ng) % ubar_obc, &
122 & boundary(ng) % vbar_obc, &
123 & boundary(ng) % zeta_obc, &
124# endif
125# ifdef ADJUST_WSTRESS
126 & forces(ng) % ustr, &
127 & forces(ng) % vstr, &
128# endif
129# ifdef SOLVE3D
130# ifdef ADJUST_STFLUX
131 & forces(ng) % tflux, &
132# endif
133 & ocean(ng) % t, &
134 & ocean(ng) % u, &
135 & ocean(ng) % v, &
136# if defined WEAK_CONSTRAINT && defined TIME_CONV
137 & ocean(ng) % ubar, &
138 & ocean(ng) % vbar, &
139# endif
140# else
141 & ocean(ng) % ubar, &
142 & ocean(ng) % vbar, &
143# endif
144 & ocean(ng) % zeta, &
145# ifdef ADJUST_BOUNDARY
146# ifdef SOLVE3D
147 & boundary(ng) % tl_t_obc, &
148 & boundary(ng) % tl_u_obc, &
149 & boundary(ng) % tl_v_obc, &
150# endif
151 & boundary(ng) % tl_ubar_obc, &
152 & boundary(ng) % tl_vbar_obc, &
153 & boundary(ng) % tl_zeta_obc, &
154# endif
155# ifdef ADJUST_WSTRESS
156 & forces(ng) % tl_ustr, &
157 & forces(ng) % tl_vstr, &
158# endif
159# ifdef SOLVE3D
160# ifdef ADJUST_STFLUX
161 & forces(ng) % tl_tflux, &
162# endif
163 & ocean(ng) % tl_t, &
164 & ocean(ng) % tl_u, &
165 & ocean(ng) % tl_v, &
166# if defined WEAK_CONSTRAINT && defined TIME_CONV
167 & ocean(ng) % tl_ubar, &
168 & ocean(ng) % tl_vbar, &
169# endif
170# else
171 & ocean(ng) % tl_ubar, &
172 & ocean(ng) % tl_vbar, &
173# endif
174 & ocean(ng) % tl_zeta, &
175# ifdef ADJUST_BOUNDARY
176# ifdef SOLVE3D
177 & boundary(ng) % ad_t_obc, &
178 & boundary(ng) % ad_u_obc, &
179 & boundary(ng) % ad_v_obc, &
180# endif
181 & boundary(ng) % ad_ubar_obc, &
182 & boundary(ng) % ad_vbar_obc, &
183 & boundary(ng) % ad_zeta_obc, &
184# endif
185# ifdef ADJUST_WSTRESS
186 & forces(ng) % ad_ustr, &
187 & forces(ng) % ad_vstr, &
188# endif
189# ifdef SOLVE3D
190# ifdef ADJUST_STFLUX
191 & forces(ng) % ad_tflux, &
192# endif
193 & ocean(ng) % ad_t, &
194 & ocean(ng) % ad_u, &
195 & ocean(ng) % ad_v, &
196# if defined WEAK_CONSTRAINT && defined TIME_CONV
197 & ocean(ng) % ad_ubar, &
198 & ocean(ng) % ad_vbar, &
199# endif
200# else
201 & ocean(ng) % ad_ubar, &
202 & ocean(ng) % ad_vbar, &
203# endif
204 & ocean(ng) % ad_zeta)
205# ifdef PROFILE
206 CALL wclock_off (ng, model, 83, __line__, myfile)
207# endif
208!
209 RETURN
210 END SUBROUTINE posterior_var
211!
212!***********************************************************************
213 SUBROUTINE posterior_var_tile (ng, tile, model, &
214 & LBi, UBi, LBj, UBj, LBij, UBij, &
215 & IminS, ImaxS, JminS, JmaxS, &
216 & Lold, Lnew, outLoop, &
217# ifdef MASKING
218 & rmask, umask, vmask, &
219# endif
220# ifdef ADJUST_BOUNDARY
221# ifdef SOLVE3D
222 & t_obc_std, u_obc_std, v_obc_std, &
223# endif
224 & ubar_obc_std, vbar_obc_std, &
225 & zeta_obc_std, &
226# endif
227# ifdef ADJUST_WSTRESS
228 & sustr_std, svstr_std, &
229# endif
230# if defined ADJUST_STFLUX && defined SOLVE3D
231 & stflx_std, &
232# endif
233# ifdef SOLVE3D
234 & t_std, u_std, v_std, &
235# if defined WEAK_CONSTRAINT && defined TIME_CONV
236 & ubar_std, vbar_std, &
237# endif
238# else
239 & ubar_std, vbar_std, &
240# endif
241 & zeta_std, &
242# ifdef ADJUST_BOUNDARY
243# ifdef SOLVE3D
244 & nl_t_obc, nl_u_obc, nl_v_obc, &
245# endif
246 & nl_ubar_obc, nl_vbar_obc, &
247 & nl_zeta_obc, &
248# endif
249# ifdef ADJUST_WSTRESS
250 & nl_ustr, nl_vstr, &
251# endif
252# ifdef SOLVE3D
253# ifdef ADJUST_STFLUX
254 & nl_tflux, &
255# endif
256 & nl_t, nl_u, nl_v, &
257# if defined WEAK_CONSTRAINT && defined TIME_CONV
258 & nl_ubar, nl_vbar, &
259# endif
260# else
261 & nl_ubar, nl_vbar, &
262# endif
263 & nl_zeta, &
264# ifdef ADJUST_BOUNDARY
265# ifdef SOLVE3D
266 & tl_t_obc, tl_u_obc, tl_v_obc, &
267# endif
268 & tl_ubar_obc, tl_vbar_obc, &
269 & tl_zeta_obc, &
270# endif
271# ifdef ADJUST_WSTRESS
272 & tl_ustr, tl_vstr, &
273# endif
274# ifdef SOLVE3D
275# ifdef ADJUST_STFLUX
276 & tl_tflux, &
277# endif
278 & tl_t, tl_u, tl_v, &
279# if defined WEAK_CONSTRAINT && defined TIME_CONV
280 & tl_ubar, tl_vbar, &
281# endif
282# else
283 & tl_ubar, tl_vbar, &
284# endif
285 & tl_zeta, &
286# ifdef ADJUST_BOUNDARY
287# ifdef SOLVE3D
288 & ad_t_obc, ad_u_obc, ad_v_obc, &
289# endif
290 & ad_ubar_obc, ad_vbar_obc, &
291 & ad_zeta_obc, &
292# endif
293# ifdef ADJUST_WSTRESS
294 & ad_ustr, ad_vstr, &
295# endif
296# ifdef SOLVE3D
297# ifdef ADJUST_STFLUX
298 & ad_tflux, &
299# endif
300 & ad_t, ad_u, ad_v, &
301# if defined WEAK_CONSTRAINT && defined TIME_CONV
302 & ad_ubar, ad_vbar, &
303# endif
304# else
305 & ad_ubar, ad_vbar, &
306# endif
307 & ad_zeta)
308!***********************************************************************
309!
310! Imported variable declarations.
311!
312 integer, intent(in) :: ng, tile, model
313 integer, intent(in) :: LBi, UBi, LBj, UBj, LBij, UBij
314 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
315 integer, intent(in) :: Lold, Lnew, outLoop
316!
317# ifdef ASSUMED_SHAPE
318# ifdef MASKING
319 real(r8), intent(in) :: rmask(LBi:,LBj:)
320 real(r8), intent(in) :: umask(LBi:,LBj:)
321 real(r8), intent(in) :: vmask(LBi:,LBj:)
322# endif
323# ifdef ADJUST_BOUNDARY
324# ifdef SOLVE3D
325 real(r8), intent(in) :: t_obc_std(LBij:,:,:,:)
326 real(r8), intent(in) :: u_obc_std(LBij:,:,:)
327 real(r8), intent(in) :: v_obc_std(LBij:,:,:)
328# endif
329 real(r8), intent(in) :: ubar_obc_std(LBij:,:)
330 real(r8), intent(in) :: vbar_obc_std(LBij:,:)
331 real(r8), intent(in) :: zeta_obc_std(LBij:,:)
332# endif
333# ifdef ADJUST_WSTRESS
334 real(r8), intent(in) :: sustr_std(LBi:,LBj:)
335 real(r8), intent(in) :: svstr_std(LBi:,LBj:)
336# endif
337# if defined ADJUST_STFLUX && defined SOLVE3D
338 real(r8), intent(in) :: stflx_std(LBi:,LBj:,:)
339# endif
340# ifdef SOLVE3D
341 real(r8), intent(in) :: t_std(LBi:,LBj:,:,:,:)
342 real(r8), intent(in) :: u_std(LBi:,LBj:,:,:)
343 real(r8), intent(in) :: v_std(LBi:,LBj:,:,:)
344# if defined WEAK_CONSTRAINT && defined TIME_CONV
345 real(r8), intent(in) :: ubar_std(LBi:,LBj:,:)
346 real(r8), intent(in) :: vbar_std(LBi:,LBj:,:)
347# endif
348# else
349 real(r8), intent(in) :: ubar_std(LBi:,LBj:,:)
350 real(r8), intent(in) :: vbar_std(LBi:,LBj:,:)
351# endif
352 real(r8), intent(in) :: zeta_std(LBi:,LBj:,:)
353# ifdef ADJUST_BOUNDARY
354# ifdef SOLVE3D
355 real(r8), intent(inout) :: ad_t_obc(LBij:,:,:,:,:,:)
356 real(r8), intent(inout) :: ad_u_obc(LBij:,:,:,:,:)
357 real(r8), intent(inout) :: ad_v_obc(LBij:,:,:,:,:)
358# endif
359 real(r8), intent(inout) :: ad_ubar_obc(LBij:,:,:,:)
360 real(r8), intent(inout) :: ad_vbar_obc(LBij:,:,:,:)
361 real(r8), intent(inout) :: ad_zeta_obc(LBij:,:,:,:)
362# endif
363# ifdef ADJUST_WSTRESS
364 real(r8), intent(inout) :: ad_ustr(LBi:,LBj:,:,:)
365 real(r8), intent(inout) :: ad_vstr(LBi:,LBj:,:,:)
366# endif
367# ifdef SOLVE3D
368# ifdef ADJUST_STFLUX
369 real(r8), intent(inout) :: ad_tflux(LBi:,LBj:,:,:,:)
370# endif
371 real(r8), intent(inout) :: ad_t(LBi:,LBj:,:,:,:)
372 real(r8), intent(inout) :: ad_u(LBi:,LBj:,:,:)
373 real(r8), intent(inout) :: ad_v(LBi:,LBj:,:,:)
374# if defined WEAK_CONSTRAINT && defined TIME_CONV
375 real(r8), intent(inout) :: ad_ubar(LBi:,LBj:,:)
376 real(r8), intent(inout) :: ad_vbar(LBi:,LBj:,:)
377# endif
378# else
379 real(r8), intent(inout) :: ad_ubar(LBi:,LBj:,:)
380 real(r8), intent(inout) :: ad_vbar(LBi:,LBj:,:)
381# endif
382 real(r8), intent(inout) :: ad_zeta(LBi:,LBj:,:)
383# ifdef ADJUST_BOUNDARY
384# ifdef SOLVE3D
385 real(r8), intent(inout) :: nl_t_obc(LBij:,:,:,:,:,:)
386 real(r8), intent(inout) :: nl_u_obc(LBij:,:,:,:,:)
387 real(r8), intent(inout) :: nl_v_obc(LBij:,:,:,:,:)
388# endif
389 real(r8), intent(inout) :: nl_ubar_obc(LBij:,:,:,:)
390 real(r8), intent(inout) :: nl_vbar_obc(LBij:,:,:,:)
391 real(r8), intent(inout) :: nl_zeta_obc(LBij:,:,:,:)
392# endif
393# ifdef ADJUST_WSTRESS
394 real(r8), intent(inout) :: nl_ustr(LBi:,LBj:,:,:)
395 real(r8), intent(inout) :: nl_vstr(LBi:,LBj:,:,:)
396# endif
397# ifdef SOLVE3D
398# ifdef ADJUST_STFLUX
399 real(r8), intent(inout) :: nl_tflux(LBi:,LBj:,:,:,:)
400# endif
401 real(r8), intent(inout) :: nl_t(LBi:,LBj:,:,:,:)
402 real(r8), intent(inout) :: nl_u(LBi:,LBj:,:,:)
403 real(r8), intent(inout) :: nl_v(LBi:,LBj:,:,:)
404# if defined WEAK_CONSTRAINT && defined TIME_CONV
405 real(r8), intent(inout) :: nl_ubar(LBi:,LBj:,:)
406 real(r8), intent(inout) :: nl_vbar(LBi:,LBj:,:)
407# endif
408# else
409 real(r8), intent(inout) :: nl_ubar(LBi:,LBj:,:)
410 real(r8), intent(inout) :: nl_vbar(LBi:,LBj:,:)
411# endif
412 real(r8), intent(inout) :: nl_zeta(LBi:,LBj:,:)
413# ifdef ADJUST_BOUNDARY
414# ifdef SOLVE3D
415 real(r8), intent(inout) :: tl_t_obc(LBij:,:,:,:,:,:)
416 real(r8), intent(inout) :: tl_u_obc(LBij:,:,:,:,:)
417 real(r8), intent(inout) :: tl_v_obc(LBij:,:,:,:,:)
418# endif
419 real(r8), intent(inout) :: tl_ubar_obc(LBij:,:,:,:)
420 real(r8), intent(inout) :: tl_vbar_obc(LBij:,:,:,:)
421 real(r8), intent(inout) :: tl_zeta_obc(LBij:,:,:,:)
422# endif
423# ifdef ADJUST_WSTRESS
424 real(r8), intent(inout) :: tl_ustr(LBi:,LBj:,:,:)
425 real(r8), intent(inout) :: tl_vstr(LBi:,LBj:,:,:)
426# endif
427# ifdef SOLVE3D
428# ifdef ADJUST_STFLUX
429 real(r8), intent(inout) :: tl_tflux(LBi:,LBj:,:,:,:)
430# endif
431 real(r8), intent(inout) :: tl_t(LBi:,LBj:,:,:,:)
432 real(r8), intent(inout) :: tl_u(LBi:,LBj:,:,:)
433 real(r8), intent(inout) :: tl_v(LBi:,LBj:,:,:)
434# if defined WEAK_CONSTRAINT && defined TIME_CONV
435 real(r8), intent(inout) :: tl_ubar(LBi:,LBj:,:)
436 real(r8), intent(inout) :: tl_vbar(LBi:,LBj:,:)
437# endif
438# else
439 real(r8), intent(inout) :: tl_ubar(LBi:,LBj:,:)
440 real(r8), intent(inout) :: tl_vbar(LBi:,LBj:,:)
441# endif
442 real(r8), intent(inout) :: tl_zeta(LBi:,LBj:,:)
443
444# else
445
446# ifdef MASKING
447 real(r8), intent(in) :: rmask(LBi:UBi,LBj:UBj)
448 real(r8), intent(in) :: umask(LBi:UBi,LBj:UBj)
449 real(r8), intent(in) :: vmask(LBi:UBi,LBj:UBj)
450# endif
451# ifdef ADJUST_BOUNDARY
452# ifdef SOLVE3D
453 real(r8), intent(in) :: t_obc_std(LBij:UBij,N(ng),4,NT(ng))
454 real(r8), intent(in) :: u_obc_std(LBij:UBij,N(ng),4)
455 real(r8), intent(in) :: v_obc_std(LBij:UBij,N(ng),4)
456# endif
457 real(r8), intent(in) :: ubar_obc_std(LBij:UBij,4)
458 real(r8), intent(in) :: vbar_obc_std(LBij:UBij,4)
459 real(r8), intent(in) :: zeta_obc_std(LBij:UBij,4)
460# endif
461# ifdef ADJUST_WSTRESS
462 real(r8), intent(in) :: sustr_std(LBi:,LBj:)
463 real(r8), intent(in) :: svstr_std(LBi:,LBj:)
464# endif
465# if defined ADJUST_STFLUX && defined SOLVE3D
466 real(r8), intent(in) :: stflx_std(LBi:UBi,LBj:UBj,NT(ng))
467# endif
468# ifdef SOLVE3D
469 real(r8), intent(in) :: t_std(LBi:UBi,LBj:UBj,N(ng),NSA,NT(ng))
470 real(r8), intent(in) :: u_std(LBi:UBi,LBj:UBj,N(ng),NSA)
471 real(r8), intent(in) :: v_std(LBi:UBi,LBj:UBj,N(ng),NSA)
472# if defined WEAK_CONSTRAINT && defined TIME_CONV
473 real(r8), intent(in) :: ubar_std(LBi:UBi,LBj:UBj,NSA)
474 real(r8), intent(in) :: vbar_std(LBi:UBi,LBj:UBj,NSA)
475# endif
476# else
477 real(r8), intent(in) :: ubar_std(LBi:UBi,LBj:UBj,NSA)
478 real(r8), intent(in) :: vbar_std(LBi:UBi,LBj:UBj,NSA)
479# endif
480 real(r8), intent(in) :: zeta_std(LBi:UBi,LBj:UBj,NSA)
481# ifdef ADJUST_BOUNDARY
482# ifdef SOLVE3D
483 real(r8), intent(inout) :: ad_t_obc(LBij:UBij,N(ng),4, &
484 & Nbrec(ng),2,NT(ng))
485 real(r8), intent(inout) :: ad_u_obc(LBij:UBij,N(ng),4,Nbrec(ng),2)
486 real(r8), intent(inout) :: ad_v_obc(LBij:UBij,N(ng),4,Nbrec(ng),2)
487# endif
488 real(r8), intent(inout) :: ad_ubar_obc(LBij:UBij,4,Nbrec(ng),2)
489 real(r8), intent(inout) :: ad_vbar_obc(LBij:UBij,4,Nbrec(ng),2)
490 real(r8), intent(inout) :: ad_zeta_obc(LBij:UBij,4,Nbrec(ng),2)
491# endif
492# ifdef ADJUST_WSTRESS
493 real(r8), intent(inout) :: ad_ustr(LBi:UBi,LBj:UBj,Nfrec(ng),2)
494 real(r8), intent(inout) :: ad_vstr(LBi:UBi,LBj:UBj,Nfrec(ng),2)
495# endif
496# ifdef SOLVE3D
497# ifdef ADJUST_STFLUX
498 real(r8), intent(inout) :: ad_tflux(LBi:UBi,LBj:UBj, &
499 & Nfrec(ng),2,NT(ng))
500# endif
501 real(r8), intent(inout) :: ad_t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
502 real(r8), intent(inout) :: ad_u(LBi:UBi,LBj:UBj,N(ng),2)
503 real(r8), intent(inout) :: ad_v(LBi:UBi,LBj:UBj,N(ng),2)
504# if defined WEAK_CONSTRAINT && defined TIME_CONV
505 real(r8), intent(inout) :: ad_ubar(LBi:UBi,LBj:UBj,:)
506 real(r8), intent(inout) :: ad_vbar(LBi:UBi,LBj:UBj,:)
507# endif
508# else
509 real(r8), intent(inout) :: ad_ubar(LBi:UBi,LBj:UBj,:)
510 real(r8), intent(inout) :: ad_vbar(LBi:UBi,LBj:UBj,:)
511# endif
512 real(r8), intent(inout) :: ad_zeta(LBi:UBi,LBj:UBj,:)
513# ifdef ADJUST_BOUNDARY
514# ifdef SOLVE3D
515 real(r8), intent(inout) :: nl_t_obc(LBij:UBij,N(ng),4, &
516 & Nbrec(ng),2,NT(ng))
517 real(r8), intent(inout) :: nl_u_obc(LBij:UBij,N(ng),4,Nbrec(ng),2)
518 real(r8), intent(inout) :: nl_v_obc(LBij:UBij,N(ng),4,Nbrec(ng),2)
519# endif
520 real(r8), intent(inout) :: nl_ubar_obc(LBij:UBij,4,Nbrec(ng),2)
521 real(r8), intent(inout) :: nl_vbar_obc(LBij:UBij,4,Nbrec(ng),2)
522 real(r8), intent(inout) :: nl_zeta_obc(LBij:UBij,4,Nbrec(ng),2)
523# endif
524# ifdef ADJUST_WSTRESS
525 real(r8), intent(inout) :: nl_ustr(LBi:UBi,LBj:UBj,Nfrec(ng),2)
526 real(r8), intent(inout) :: nl_vstr(LBi:UBi,LBj:UBj,Nfrec(ng),2)
527# endif
528# ifdef SOLVE3D
529# ifdef ADJUST_STFLUX
530 real(r8), intent(inout) :: nl_tflux(LBi:UBi,LBj:UBj, &
531 & Nfrec(ng),2,NT(ng))
532# endif
533 real(r8), intent(inout) :: nl_t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
534 real(r8), intent(inout) :: nl_u(LBi:UBi,LBj:UBj,N(ng),2)
535 real(r8), intent(inout) :: nl_v(LBi:UBi,LBj:UBj,N(ng),2)
536# if defined WEAK_CONSTRAINT && defined TIME_CONV
537 real(r8), intent(inout) :: nl_ubar(LBi:UBi,LBj:UBj,:)
538 real(r8), intent(inout) :: nl_vbar(LBi:UBi,LBj:UBj,:)
539# endif
540# else
541 real(r8), intent(inout) :: nl_ubar(LBi:UBi,LBj:UBj,:)
542 real(r8), intent(inout) :: nl_vbar(LBi:UBi,LBj:UBj,:)
543# endif
544 real(r8), intent(inout) :: nl_zeta(LBi:UBi,LBj:UBj,:)
545# ifdef ADJUST_BOUNDARY
546# ifdef SOLVE3D
547 real(r8), intent(inout) :: tl_t_obc(LBij:UBij,N(ng),4, &
548 & Nbrec(ng),2,NT(ng))
549 real(r8), intent(inout) :: tl_u_obc(LBij:UBij,N(ng),4,Nbrec(ng),2)
550 real(r8), intent(inout) :: tl_v_obc(LBij:UBij,N(ng),4,Nbrec(ng),2)
551# endif
552 real(r8), intent(inout) :: tl_ubar_obc(LBij:UBij,4,Nbrec(ng),2)
553 real(r8), intent(inout) :: tl_vbar_obc(LBij:UBij,4,Nbrec(ng),2)
554 real(r8), intent(inout) :: tl_zeta_obc(LBij:UBij,4,Nbrec(ng),2)
555# endif
556# ifdef ADJUST_WSTRESS
557 real(r8), intent(inout) :: tl_ustr(LBi:UBi,LBj:UBj,Nfrec(ng),2)
558 real(r8), intent(inout) :: tl_vstr(LBi:UBi,LBj:UBj,Nfrec(ng),2)
559# endif
560# ifdef SOLVE3D
561# ifdef ADJUST_STFLUX
562 real(r8), intent(inout) :: tl_tflux(LBi:UBi,LBj:UBj, &
563 & Nfrec(ng),2,NT(ng))
564# endif
565 real(r8), intent(inout) :: tl_t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
566 real(r8), intent(inout) :: tl_u(LBi:UBi,LBj:UBj,N(ng),2)
567 real(r8), intent(inout) :: tl_v(LBi:UBi,LBj:UBj,N(ng),2)
568# if defined WEAK_CONSTRAINT && defined TIME_CONV
569 real(r8), intent(inout) :: tl_ubar(LBi:UBi,LBj:UBj,:)
570 real(r8), intent(inout) :: tl_vbar(LBi:UBi,LBj:UBj,:)
571# endif
572# else
573 real(r8), intent(inout) :: tl_ubar(LBi:UBi,LBj:UBj,:)
574 real(r8), intent(inout) :: tl_vbar(LBi:UBi,LBj:UBj,:)
575# endif
576 real(r8), intent(inout) :: tl_zeta(LBi:UBi,LBj:UBj,:)
577# endif
578!
579! Local variable declarations.
580!
581 integer :: L1 = 1
582 integer :: L2 = 2
583
584 integer :: Linp, Lout, Lscale, Lwrk, Lwrk1, i, j, ic
585 integer :: info, ivec, jvec, rec
586
587 real(r8) :: fac, fac1, fac2
588
589 logical :: Lweak
590
591 character (len=256) :: ncname
592
593# include "set_bounds.h"
594!
595 calledfrom=myfile
596 sourcefile=myfile
597!
598!-----------------------------------------------------------------------
599! Computes the full (diagonal) posterior analysis error covariance
600# if defined POSTERIOR_ERROR_I
601! matrix at the initial time of the assimilation window.
602# elif defined POSTERIOR_ERROR_F
603! matrix at the final time of the assimilation window.
604# endif
605!-----------------------------------------------------------------------
606!
607! NOTE: In the case of weak constraint ("WEAK_CONSTRAINT") and
608! and time convolutions ("TIME_CONV"), the state arrays
609! tl_ubar, tl_vbar, ad_ubar, and ad_vbar are only passed
610! as required by the "state" operators routines but they
611! are not used in subsequent calculations.
612!
613 IF (master) WRITE (stdout,10)
614 10 FORMAT (/,' <<<< Full Posterior Error Covariance Matrix >>>>',/)
615!
616 lweak=.false.
617!
618! Invert the tridiagonal matrix of inner-loop Lanczos vector
619! coefficients.
620!
621 zlanczos_coef=0.0_r8
622 DO i=1,ninner
623 zlanczos_coef(i,i)=cg_delta(i,outloop)
624 END DO
625 DO i=1,ninner-1
626 zlanczos_coef(i,i+1)=cg_beta(i+1,outloop)
627 END DO
628 DO i=2,ninner
629 zlanczos_coef(i,i-1)=cg_beta(i,outloop)
630 END DO
632!
633! Invert "zLanczos_coef" using LAPACK routines DPOTRF and DPOTRI.
634! On exit "zLanczos_inv" contains the inverse of "zLanczos_coef".
635! The "zLanczos_coef" matrix is saved to check the inversion below.
636!
637 IF (master) THEN
638 CALL dpotrf ('U', ninner, zlanczos_inv, ninner, info)
639 END IF
640# ifdef DISTRIBUTE
641 CALL mp_bcasti (ng, model, info)
642# endif
643 IF (info.ne.0) THEN
644 IF (master) WRITE (stdout,*) ' Error in DPOTRF: info = ', info
645 exit_flag=8
646 RETURN
647 END IF
648!
649 IF (master) THEN
650 CALL dpotri ('U', ninner, zlanczos_inv, ninner, info)
651 END IF
652# ifdef DISTRIBUTE
653 CALL mp_bcasti (ng, model, info)
654 CALL mp_bcastf (ng, model, zlanczos_inv)
655# endif
656 IF (info.ne.0) THEN
657 IF (master) WRITE (stdout,*) ' Error in DPOTRI: info = ', info
658 exit_flag=8
659 RETURN
660 END IF
661!
662! Copy the upper triangle of the inverse into the lower triangle
663! as required by the solution from DPOTRI.
664!
665 DO j=1,ninner
666 DO i=j+1,ninner
667 zlanczos_inv(i,j)=zlanczos_inv(j,i)
668 END DO
669 END DO
670!
671! Test the inverse, we need to get the identity matrix within roundoff.
672!
673 DO j=1,ninner
674 DO i=1,ninner
675 zlanczos_err(i,j)=0.0_r8
676 DO ic=1,ninner
677 zlanczos_err(i,j)=zlanczos_err(i,j)+ &
678 & zlanczos_inv(i,ic)*zlanczos_coef(ic,j)
679 END DO
680 END DO
681 END DO
682!
683! Now compute the diagonal of the posterior analysis error covariance
684! matrix. This will be stored in record L1 of the tl_var arrays.
685!
686 fac=0.0_r8
687 CALL state_initialize (ng, tile, &
688 & lbi, ubi, lbj, ubj, lbij, ubij, &
689 & l1, fac, &
690# ifdef MASKING
691 & rmask, umask, vmask, &
692# endif
693# ifdef ADJUST_BOUNDARY
694# ifdef SOLVE3D
695 & tl_t_obc, tl_u_obc, tl_v_obc, &
696# endif
697 & tl_ubar_obc, tl_vbar_obc, &
698 & tl_zeta_obc, &
699# endif
700# ifdef ADJUST_WSTRESS
701 & tl_ustr, tl_vstr, &
702# endif
703# ifdef SOLVE3D
704# ifdef ADJUST_STFLUX
705 & tl_tflux, &
706# endif
707 & tl_t, tl_u, tl_v, &
708# else
709 & tl_ubar, tl_vbar, &
710# endif
711 & tl_zeta)
712!
713! Read in turn each evolved and convolved Lanczos vector of the
714! stabilized representer matrix and compute the product with all
715! the other vectors. Use the ad_var arrays as temporary storage.
716!
717 ncname=hss(ng)%name
718 DO ivec=1,ninner
719 rec=ivec
720!
721! Read vectors stored in hessian netcdf file using ad_var(L1) as
722! temporary storage.
723!
724 CALL state_read (ng, tile, model, hss(ng)%IOtype, &
725 & lbi, ubi, lbj, ubj, lbij, ubij, &
726 & l1, rec, &
727 & 0, hss(ng)%ncid, ncname, &
728# ifdef MASKING
729 & rmask, umask, vmask, &
730# endif
731# ifdef ADJUST_BOUNDARY
732# ifdef SOLVE3D
733 & ad_t_obc, ad_u_obc, ad_v_obc, &
734# endif
735 & ad_ubar_obc, ad_vbar_obc, &
736 & ad_zeta_obc, &
737# endif
738# ifdef ADJUST_WSTRESS
739 & ad_ustr, ad_vstr, &
740# endif
741# ifdef SOLVE3D
742# ifdef ADJUST_STFLUX
743 & ad_tflux, &
744# endif
745 & ad_t, ad_u, ad_v, &
746# else
747 & ad_ubar, ad_vbar, &
748# endif
749 & ad_zeta)
750!
751! Compute state product of ad_var(L1) with itself using nl_var(L1)
752! as temporary storage.
753!
754 CALL state_product (ng, tile, model, &
755 & lbi, ubi, lbj, ubj, lbij, ubij, &
756# ifdef MASKING
757 & rmask, umask, vmask, &
758# endif
759# ifdef ADJUST_BOUNDARY
760# ifdef SOLVE3D
761 & ad_t_obc(:,:,:,:,l1,:), &
762 & ad_t_obc(:,:,:,:,l1,:), &
763 & nl_t_obc(:,:,:,:,l1,:), &
764 & ad_u_obc(:,:,:,:,l1), &
765 & ad_u_obc(:,:,:,:,l1), &
766 & nl_u_obc(:,:,:,:,l1), &
767 & ad_v_obc(:,:,:,:,l1), &
768 & ad_v_obc(:,:,:,:,l1), &
769 & nl_v_obc(:,:,:,:,l1), &
770# endif
771 & ad_ubar_obc(:,:,:,l1), &
772 & ad_ubar_obc(:,:,:,l1), &
773 & nl_ubar_obc(:,:,:,l1), &
774 & ad_vbar_obc(:,:,:,l1), &
775 & ad_vbar_obc(:,:,:,l1), &
776 & nl_vbar_obc(:,:,:,l1), &
777 & ad_zeta_obc(:,:,:,l1), &
778 & ad_zeta_obc(:,:,:,l1), &
779 & nl_zeta_obc(:,:,:,l1), &
780# endif
781# ifdef ADJUST_WSTRESS
782 & ad_ustr(:,:,:,l1), &
783 & ad_ustr(:,:,:,l1), &
784 & nl_ustr(:,:,:,l1), &
785 & ad_vstr(:,:,:,l1), &
786 & ad_vstr(:,:,:,l1), &
787 & nl_vstr(:,:,:,l1), &
788# endif
789# ifdef SOLVE3D
790# ifdef ADJUST_STFLUX
791 & ad_tflux(:,:,:,l1,:), &
792 & ad_tflux(:,:,:,l1,:), &
793 & nl_tflux(:,:,:,l1,:), &
794# endif
795 & ad_t(:,:,:,l1,:), &
796 & ad_t(:,:,:,l1,:), &
797 & nl_t(:,:,:,l1,:), &
798 & ad_u(:,:,:,l1), &
799 & ad_u(:,:,:,l1), &
800 & nl_u(:,:,:,l1), &
801 & ad_v(:,:,:,l1), &
802 & ad_v(:,:,:,l1), &
803 & nl_v(:,:,:,l1), &
804# else
805 & ad_ubar(:,:,l1), &
806 & ad_ubar(:,:,l1), &
807 & nl_ubar(:,:,l1), &
808 & ad_vbar(:,:,l1), &
809 & ad_vbar(:,:,l1), &
810 & nl_vbar(:,:,l1), &
811# endif
812 & ad_zeta(:,:,l1), &
813 & ad_zeta(:,:,l1), &
814 & nl_zeta(:,:,l1))
815!
816! tl_var(L1) = fac1 * tl_var(L1) + fac2 * nl_var(L1). See NOTE above.
817!
818 fac1=1.0_r8
819 fac2=zlanczos_inv(ivec,ivec)
820
821 CALL state_addition (ng, tile, &
822 & lbi, ubi, lbj, ubj, lbij, ubij, &
823 & l1, l1, l1, fac1, fac2, &
824# ifdef MASKING
825 & rmask, umask, vmask, &
826# endif
827# ifdef ADJUST_BOUNDARY
828# ifdef SOLVE3D
829 & tl_t_obc, nl_t_obc, &
830 & tl_u_obc, nl_u_obc, &
831 & tl_v_obc, nl_v_obc, &
832# endif
833 & tl_ubar_obc, nl_ubar_obc, &
834 & tl_vbar_obc, nl_vbar_obc, &
835 & tl_zeta_obc, nl_zeta_obc, &
836# endif
837# ifdef ADJUST_WSTRESS
838 & tl_ustr, nl_ustr, &
839 & tl_vstr, nl_vstr, &
840# endif
841# ifdef SOLVE3D
842# ifdef ADJUST_STFLUX
843 & tl_tflux, nl_tflux, &
844# endif
845 & tl_t, nl_t, &
846 & tl_u, nl_u, &
847 & tl_v, nl_v, &
848# if defined WEAK_CONSTRAINT && defined TIME_CONV
849 & tl_ubar, nl_ubar, &
850 & tl_vbar, nl_vbar, &
851# endif
852# else
853 & tl_ubar, nl_ubar, &
854 & tl_vbar, nl_vbar, &
855# endif
856 & tl_zeta, nl_zeta)
857!
858 DO jvec=ivec+1,ninner
859 rec=jvec
860!
861! Read vectors stored in hessian netcdf file using ad_var(L2) as
862! temporary storage.
863!
864 CALL state_read (ng, tile, model, hss(ng)%IOtype, &
865 & lbi, ubi, lbj, ubj, lbij, ubij, &
866 & l2, rec, &
867 & 0, hss(ng)%ncid, ncname, &
868# ifdef MASKING
869 & rmask, umask, vmask, &
870# endif
871# ifdef ADJUST_BOUNDARY
872# ifdef SOLVE3D
873 & ad_t_obc, ad_u_obc, ad_v_obc, &
874# endif
875 & ad_ubar_obc, ad_vbar_obc, &
876 & ad_zeta_obc, &
877# endif
878# ifdef ADJUST_WSTRESS
879 & ad_ustr, ad_vstr, &
880# endif
881# ifdef SOLVE3D
882# ifdef ADJUST_STFLUX
883 & ad_tflux, &
884# endif
885 & ad_t, ad_u, ad_v, &
886# else
887 & ad_ubar, ad_vbar, &
888# endif
889 & ad_zeta)
890!
891! Compute state product of ad_var(L2) with ad_var(L1) using nl_var(L1)
892! as temporary storage.
893!
894 CALL state_product (ng, tile, model, &
895 & lbi, ubi, lbj, ubj, lbij, ubij, &
896# ifdef MASKING
897 & rmask, umask, vmask, &
898# endif
899# ifdef ADJUST_BOUNDARY
900# ifdef SOLVE3D
901 & ad_t_obc(:,:,:,:,l1,:), &
902 & ad_t_obc(:,:,:,:,l2,:), &
903 & nl_t_obc(:,:,:,:,l1,:), &
904 & ad_u_obc(:,:,:,:,l1), &
905 & ad_u_obc(:,:,:,:,l2), &
906 & nl_u_obc(:,:,:,:,l1), &
907 & ad_v_obc(:,:,:,:,l1), &
908 & ad_v_obc(:,:,:,:,l2), &
909 & nl_v_obc(:,:,:,:,l1), &
910# endif
911 & ad_ubar_obc(:,:,:,l1), &
912 & ad_ubar_obc(:,:,:,l2), &
913 & nl_ubar_obc(:,:,:,l1), &
914 & ad_vbar_obc(:,:,:,l1), &
915 & ad_vbar_obc(:,:,:,l2), &
916 & nl_vbar_obc(:,:,:,l1), &
917 & ad_zeta_obc(:,:,:,l1), &
918 & ad_zeta_obc(:,:,:,l2), &
919 & nl_zeta_obc(:,:,:,l1), &
920# endif
921# ifdef ADJUST_WSTRESS
922 & ad_ustr(:,:,:,l1), &
923 & ad_ustr(:,:,:,l2), &
924 & nl_ustr(:,:,:,l1), &
925 & ad_vstr(:,:,:,l1), &
926 & ad_vstr(:,:,:,l2), &
927 & nl_vstr(:,:,:,l1), &
928# endif
929# ifdef SOLVE3D
930# ifdef ADJUST_STFLUX
931 & ad_tflux(:,:,:,l1,:), &
932 & ad_tflux(:,:,:,l2,:), &
933 & nl_tflux(:,:,:,l1,:), &
934# endif
935 & ad_t(:,:,:,l1,:), &
936 & ad_t(:,:,:,l2,:), &
937 & nl_t(:,:,:,l1,:), &
938 & ad_u(:,:,:,l1), &
939 & ad_u(:,:,:,l2), &
940 & nl_u(:,:,:,l1), &
941 & ad_v(:,:,:,l1), &
942 & ad_v(:,:,:,l2), &
943 & nl_v(:,:,:,l1), &
944# else
945 & ad_ubar(:,:,l1), &
946 & ad_ubar(:,:,l2), &
947 & nl_ubar(:,:,l1), &
948 & ad_vbar(:,:,l1), &
949 & ad_vbar(:,:,l2), &
950 & nl_vbar(:,:,l1), &
951# endif
952 & ad_zeta(:,:,l1), &
953 & ad_zeta(:,:,l2), &
954 & nl_zeta(:,:,l1))
955!
956! tl_var(L1) = fac1 * tl_var(L1) + fac2 * nl_var(L1). See NOTE above.
957!
958 fac1=1.0_r8
959 fac2=2.0_r8*zlanczos_inv(ivec,jvec)
960
961 CALL state_addition (ng, tile, &
962 & lbi, ubi, lbj, ubj, lbij, ubij, &
963 & l1, l1, l1, fac1, fac2, &
964# ifdef MASKING
965 & rmask, umask, vmask, &
966# endif
967# ifdef ADJUST_BOUNDARY
968# ifdef SOLVE3D
969 & tl_t_obc, nl_t_obc, &
970 & tl_u_obc, nl_u_obc, &
971 & tl_v_obc, nl_v_obc, &
972# endif
973 & tl_ubar_obc, nl_ubar_obc, &
974 & tl_vbar_obc, nl_vbar_obc, &
975 & tl_zeta_obc, nl_zeta_obc, &
976# endif
977# ifdef ADJUST_WSTRESS
978 & tl_ustr, nl_ustr, &
979 & tl_vstr, nl_vstr, &
980# endif
981# ifdef SOLVE3D
982# ifdef ADJUST_STFLUX
983 & tl_tflux, nl_tflux, &
984# endif
985 & tl_t, nl_t, &
986 & tl_u, nl_u, &
987 & tl_v, nl_v, &
988# if defined WEAK_CONSTRAINT && defined TIME_CONV
989 & tl_ubar, nl_ubar, &
990 & tl_vbar, nl_vbar, &
991# endif
992# else
993 & tl_ubar, nl_ubar, &
994 & tl_vbar, nl_vbar, &
995# endif
996 & tl_zeta, nl_zeta)
997 END DO
998 END DO
999!
1000! Finally subtract the result in tl_var(L1) from the background
1001! error variances to get the posterior/analysis error variance.
1002!
1003 CALL analysis_var (ng, tile, model, &
1004 & lbi, ubi, lbj, ubj, lbij, ubij, &
1005 & l1, lweak, &
1006# ifdef MASKING
1007 & rmask, umask, vmask, &
1008# endif
1009# ifdef ADJUST_BOUNDARY
1010# ifdef SOLVE3D
1011 & t_obc_std, u_obc_std, v_obc_std, &
1012# endif
1013 & ubar_obc_std, vbar_obc_std, zeta_obc_std, &
1014# endif
1015# ifdef ADJUST_WSTRESS
1016 & sustr_std, svstr_std, &
1017# endif
1018# ifdef SOLVE3D
1019# ifdef ADJUST_STFLUX
1020 & stflx_std, &
1021# endif
1022 & t_std, u_std, v_std, &
1023# else
1024 & ubar_std, vbar_std, &
1025# endif
1026 & zeta_std, &
1027# ifdef ADJUST_BOUNDARY
1028# ifdef SOLVE3D
1029 & tl_t_obc, tl_u_obc, tl_v_obc, &
1030# endif
1031 & tl_ubar_obc, tl_vbar_obc, tl_zeta_obc, &
1032# endif
1033# ifdef ADJUST_WSTRESS
1034 & tl_ustr, tl_vstr, &
1035# endif
1036# ifdef SOLVE3D
1037# ifdef ADJUST_STFLUX
1038 & tl_tflux, &
1039# endif
1040 & tl_t, tl_u, tl_v, &
1041# else
1042 & tl_ubar, tl_vbar, &
1043# endif
1044 & tl_zeta)
1045!
1046 RETURN
1047 END SUBROUTINE posterior_var_tile
1048!
1049!***********************************************************************
1050 SUBROUTINE analysis_var (ng, tile, model, &
1051 & LBi, UBi, LBj, UBj, LBij, UBij, &
1052 & Lout, Lweak, &
1053# ifdef MASKING
1054 & rmask, umask, vmask, &
1055# endif
1056# ifdef ADJUST_BOUNDARY
1057# ifdef SOLVE3D
1058 & t_obc_std, u_obc_std, v_obc_std, &
1059# endif
1060 & ubar_obc_std, vbar_obc_std, &
1061 & zeta_obc_std, &
1062# endif
1063# ifdef ADJUST_WSTRESS
1064 & sustr_std, svstr_std, &
1065# endif
1066# ifdef SOLVE3D
1067# ifdef ADJUST_STFLUX
1068 & stflx_std, &
1069# endif
1070 & t_std, u_std, v_std, &
1071# else
1072 & ubar_std, vbar_std, &
1073# endif
1074 & zeta_std, &
1075# ifdef ADJUST_BOUNDARY
1076# ifdef SOLVE3D
1077 & s1_t_obc, s1_u_obc, s1_v_obc, &
1078# endif
1079 & s1_ubar_obc, s1_vbar_obc, &
1080 & s1_zeta_obc, &
1081# endif
1082# ifdef ADJUST_WSTRESS
1083 & s1_sustr, s1_svstr, &
1084# endif
1085# ifdef SOLVE3D
1086# ifdef ADJUST_STFLUX
1087 & s1_tflux, &
1088# endif
1089 & s1_t, s1_u, s1_v, &
1090# else
1091 & s1_ubar, s1_vbar, &
1092# endif
1093 & s1_zeta)
1094!***********************************************************************
1095!
1096! Imported variable declarations.
1097!
1098 integer, intent(in) :: ng, tile, model, Lout
1099 integer, intent(in) :: LBi, UBi, LBj, UBj, LBij, UBij
1100 logical, intent(in) :: Lweak
1101!
1102# ifdef ASSUMED_SHAPE
1103# ifdef MASKING
1104 real(r8), intent(in) :: rmask(LBi:,LBj:)
1105 real(r8), intent(in) :: umask(LBi:,LBj:)
1106 real(r8), intent(in) :: vmask(LBi:,LBj:)
1107# endif
1108# ifdef ADJUST_BOUNDARY
1109# ifdef SOLVE3D
1110 real(r8), intent(inout) :: s1_t_obc(LBij:,:,:,:,:,:)
1111 real(r8), intent(inout) :: s1_u_obc(LBij:,:,:,:,:)
1112 real(r8), intent(inout) :: s1_v_obc(LBij:,:,:,:,:)
1113# endif
1114 real(r8), intent(inout) :: s1_ubar_obc(LBij:,:,:,:)
1115 real(r8), intent(inout) :: s1_vbar_obc(LBij:,:,:,:)
1116 real(r8), intent(inout) :: s1_zeta_obc(LBij:,:,:,:)
1117# endif
1118# ifdef ADJUST_WSTRESS
1119 real(r8), intent(inout) :: s1_sustr(LBi:,LBj:,:,:)
1120 real(r8), intent(inout) :: s1_svstr(LBi:,LBj:,:,:)
1121# endif
1122# ifdef SOLVE3D
1123# ifdef ADJUST_STFLUX
1124 real(r8), intent(inout) :: s1_tflux(LBi:,LBj:,:,:,:)
1125# endif
1126 real(r8), intent(inout) :: s1_t(LBi:,LBj:,:,:,:)
1127 real(r8), intent(inout) :: s1_u(LBi:,LBj:,:,:)
1128 real(r8), intent(inout) :: s1_v(LBi:,LBj:,:,:)
1129# else
1130 real(r8), intent(inout) :: s1_ubar(LBi:,LBj:,:)
1131 real(r8), intent(inout) :: s1_vbar(LBi:,LBj:,:)
1132# endif
1133 real(r8), intent(inout) :: s1_zeta(LBi:,LBj:,:)
1134# ifdef ADJUST_BOUNDARY
1135# ifdef SOLVE3D
1136 real(r8), intent(in) :: t_obc_std(LBij:,:,:,:)
1137 real(r8), intent(in) :: u_obc_std(LBij:,:,:)
1138 real(r8), intent(in) :: v_obc_std(LBij:,:,:)
1139# endif
1140 real(r8), intent(in) :: ubar_obc_std(LBij:,:)
1141 real(r8), intent(in) :: vbar_obc_std(LBij:,:)
1142 real(r8), intent(in) :: zeta_obc_std(LBij:,:)
1143# endif
1144# ifdef ADJUST_WSTRESS
1145 real(r8), intent(in) :: sustr_std(LBi:,LBj:)
1146 real(r8), intent(in) :: svstr_std(LBi:,LBj:)
1147# endif
1148# if defined ADJUST_STFLUX && defined SOLVE3D
1149 real(r8), intent(in) :: stflx_std(LBi:,LBj:,:)
1150# endif
1151# ifdef SOLVE3D
1152 real(r8), intent(in) :: t_std(LBi:,LBj:,:,:,:)
1153 real(r8), intent(in) :: u_std(LBi:,LBj:,:,:)
1154 real(r8), intent(in) :: v_std(LBi:,LBj:,:,:)
1155# else
1156 real(r8), intent(in) :: ubar_std(LBi:,LBj:,:)
1157 real(r8), intent(in) :: vbar_std(LBi:,LBj:,:)
1158# endif
1159 real(r8), intent(in) :: zeta_std(LBi:,LBj:,:)
1160
1161# else
1162
1163# ifdef MASKING
1164 real(r8), intent(in) :: rmask(LBi:UBi,LBj:UBj)
1165 real(r8), intent(in) :: umask(LBi:UBi,LBj:UBj)
1166 real(r8), intent(in) :: vmask(LBi:UBi,LBj:UBj)
1167# endif
1168
1169# ifdef ADJUST_BOUNDARY
1170# ifdef SOLVE3D
1171 real(r8), intent(inout) :: s1_t_obc(LBij:UBij,N(ng),4, &
1172 & Nbrec(ng),2,NT(ng))
1173 real(r8), intent(inout) :: s1_u_obc(LBij:UBij,N(ng),4,Nbrec(ng),2)
1174 real(r8), intent(inout) :: s1_v_obc(LBij:UBij,N(ng),4,Nbrec(ng),2)
1175# endif
1176 real(r8), intent(inout) :: s1_ubar_obc(LBij:UBij,4,Nbrec(ng),2)
1177 real(r8), intent(inout) :: s1_vbar_obc(LBij:UBij,4,Nbrec(ng),2)
1178 real(r8), intent(inout) :: s1_zeta_obc(LBij:UBij,4,Nbrec(ng),2)
1179# endif
1180# ifdef ADJUST_WSTRESS
1181 real(r8), intent(inout) :: s1_sustr(LBi:UBi,LBj:UBj,Nfrec(ng),2)
1182 real(r8), intent(inout) :: s1_svstr(LBi:UBi,LBj:UBj,Nfrec(ng),2)
1183# endif
1184# ifdef SOLVE3D
1185# ifdef ADJUST_STFLUX
1186 real(r8), intent(inout) :: s1_tflux(LBi:UBi,LBj:UBj, &
1187 & Nfrec(ng),2,NT(ng))
1188# endif
1189 real(r8), intent(inout) :: s1_t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
1190 real(r8), intent(inout) :: s1_u(LBi:UBi,LBj:UBj,N(ng),2)
1191 real(r8), intent(inout) :: s1_v(LBi:UBi,LBj:UBj,N(ng),2)
1192# else
1193 real(r8), intent(inout) :: s1_ubar(LBi:UBi,LBj:UBj,:)
1194 real(r8), intent(inout) :: s1_vbar(LBi:UBi,LBj:UBj,:)
1195# endif
1196 real(r8), intent(inout) :: s1_zeta(LBi:UBi,LBj:UBj,:)
1197# ifdef ADJUST_BOUNDARY
1198# ifdef SOLVE3D
1199 real(r8), intent(in) :: t_obc_std(LBij:UBij,N(ng),4,NT(ng))
1200 real(r8), intent(in) :: u_obc_std(LBij:UBij,N(ng),4)
1201 real(r8), intent(in) :: v_obc_std(LBij:UBij,N(ng),4)
1202# endif
1203 real(r8), intent(in) :: ubar_obc_std(LBij:UBij,4)
1204 real(r8), intent(in) :: vbar_obc_std(LBij:UBij,4)
1205 real(r8), intent(in) :: zeta_obc_std(LBij:UBij,4)
1206# endif
1207# ifdef ADJUST_WSTRESS
1208 real(r8), intent(in) :: sustr_std(LBi:,LBj:)
1209 real(r8), intent(in) :: svstr_std(LBi:,LBj:)
1210# endif
1211# if defined ADJUST_STFLUX && defined SOLVE3D
1212 real(r8), intent(in) :: stflx_std(LBi:UBi,LBj:UBj,NT(ng))
1213# endif
1214# ifdef SOLVE3D
1215 real(r8), intent(in) :: t_std(LBi:UBi,LBj:UBj,N(ng),NSA,NT(ng))
1216 real(r8), intent(in) :: u_std(LBi:UBi,LBj:UBj,N(ng),NSA)
1217 real(r8), intent(in) :: v_std(LBi:UBi,LBj:UBj,N(ng),NSA)
1218# else
1219 real(r8), intent(in) :: ubar_std(LBi:UBi,LBj:UBj,NSA)
1220 real(r8), intent(in) :: vbar_std(LBi:UBi,LBj:UBj,NSA)
1221# endif
1222 real(r8), intent(in) :: zeta_std(LBi:UBi,LBj:UBj,NSA)
1223# endif
1224!
1225! Local variable declarations.
1226!
1227 integer :: NSUB, i, j, k
1228 integer :: ir, it, rec
1229
1230 real(r8) :: cff
1231
1232# include "set_bounds.h"
1233!
1234!-----------------------------------------------------------------------
1235! Determine error covariance standard deviation factors to use.
1236!-----------------------------------------------------------------------
1237!
1238 IF (lweak) THEN
1239 rec=2 ! weak constraint
1240 ELSE
1241 rec=1 ! strong constraint
1242 END IF
1243!
1244!-----------------------------------------------------------------------
1245! Compute S1=STD*STD-S1.
1246!-----------------------------------------------------------------------
1247!
1248! Free-surface.
1249!
1250 DO j=jstrt,jendt
1251 DO i=istrt,iendt
1252 cff=zeta_std(i,j,rec)*zeta_std(i,j,rec)- &
1253 & s1_zeta(i,j,lout)
1254# ifdef MASKING
1255 cff=cff*rmask(i,j)
1256# endif
1257 s1_zeta(i,j,lout)=cff
1258 END DO
1259 END DO
1260
1261# ifdef ADJUST_BOUNDARY
1262!
1263! Free-surface open boundaries.
1264!
1265 IF (any(lobc(:,isfsur,ng))) THEN
1266 DO ir=1,nbrec(ng)
1267 IF ((lobc(iwest,isfsur,ng)).and. &
1268 & domain(ng)%Western_Edge(tile)) THEN
1269 DO j=jstr,jend
1270 cff=zeta_obc_std(j,iwest)*zeta_obc_std(j,iwest)- &
1271 & s1_zeta_obc(j,iwest,ir,lout)
1272# ifdef MASKING
1273 cff=cff*rmask(istr-1,j)
1274# endif
1275 s1_zeta_obc(j,iwest,ir,lout)=cff
1276 END DO
1277 END IF
1278 IF ((lobc(ieast,isfsur,ng)).and. &
1279 & domain(ng)%Eastern_Edge(tile)) THEN
1280 DO j=jstr,jend
1281 cff=zeta_obc_std(j,ieast)*zeta_obc_std(j,ieast)- &
1282 & s1_zeta_obc(j,ieast,ir,lout)
1283# ifdef MASKING
1284 cff=cff*rmask(iend+1,j)
1285# endif
1286 s1_zeta_obc(j,ieast,ir,lout)=cff
1287 END DO
1288 END IF
1289 IF ((lobc(isouth,isfsur,ng)).and. &
1290 & domain(ng)%Southern_Edge(tile)) THEN
1291 DO i=istr,iend
1292 cff=zeta_obc_std(i,isouth)*zeta_obc_std(i,isouth)- &
1293 & s1_zeta_obc(i,isouth,ir,lout)
1294# ifdef MASKING
1295 cff=cff*rmask(i,jstr-1)
1296# endif
1297 s1_zeta_obc(i,isouth,ir,lout)=cff
1298 END DO
1299 END IF
1300 IF ((lobc(inorth,isfsur,ng)).and. &
1301 & domain(ng)%Northern_Edge(tile)) THEN
1302 DO i=istr,iend
1303 cff=zeta_obc_std(i,inorth)*zeta_obc_std(i,inorth)- &
1304 & s1_zeta_obc(i,inorth,ir,lout)
1305# ifdef MASKING
1306 cff=cff*rmask(i,jend+1)
1307# endif
1308 s1_zeta_obc(i,inorth,ir,lout)=cff
1309 END DO
1310 END IF
1311 END DO
1312 END IF
1313# endif
1314
1315# ifndef SOLVE3D
1316!
1317! 2D U-momentum component.
1318!
1319 DO j=jstrt,jendt
1320 DO i=istrp,iendt
1321 cff=ubar_std(i,j,rec)*ubar_std(i,j,rec)- &
1322 & s1_ubar(i,j,lout)
1323# ifdef MASKING
1324 cff=cff*umask(i,j)
1325# endif
1326 s1_ubar(i,j,lout)=cff
1327 END DO
1328 END DO
1329# endif
1330
1331# ifdef ADJUST_BOUNDARY
1332!
1333! 2D U-momentum open boundaries.
1334!
1335 IF (any(lobc(:,isubar,ng))) THEN
1336 DO ir=1,nbrec(ng)
1337 IF ((lobc(iwest,isubar,ng)).and. &
1338 & domain(ng)%Western_Edge(tile)) THEN
1339 DO j=jstr,jend
1340 cff=ubar_obc_std(j,iwest)*ubar_obc_std(j,iwest)- &
1341 & s1_ubar_obc(j,iwest,ir,lout)
1342# ifdef MASKING
1343 cff=cff*umask(istr,j)
1344# endif
1345 s1_ubar_obc(j,iwest,ir,lout)=cff
1346 END DO
1347 END IF
1348 IF ((lobc(ieast,isubar,ng)).and. &
1349 & domain(ng)%Eastern_Edge(tile)) THEN
1350 DO j=jstr,jend
1351 cff=ubar_obc_std(j,ieast)*ubar_obc_std(j,ieast)- &
1352 & s1_ubar_obc(j,ieast,ir,lout)
1353# ifdef MASKING
1354 cff=cff*umask(iend+1,j)
1355# endif
1356 s1_ubar_obc(j,ieast,ir,lout)=cff
1357 END DO
1358 END IF
1359 IF ((lobc(isouth,isubar,ng)).and. &
1360 & domain(ng)%Southern_Edge(tile)) THEN
1361 DO i=istru,iend
1362 cff=ubar_obc_std(i,isouth)*ubar_obc_std(i,isouth)- &
1363 & s1_ubar_obc(i,isouth,ir,lout)
1364# ifdef MASKING
1365 cff=cff*umask(i,jstr-1)
1366# endif
1367 s1_ubar_obc(i,isouth,ir,lout)=cff
1368 END DO
1369 END IF
1370 IF ((lobc(inorth,isubar,ng)).and. &
1371 & domain(ng)%Northern_Edge(tile)) THEN
1372 DO i=istru,iend
1373 cff=ubar_obc_std(i,inorth)*ubar_obc_std(i,inorth)- &
1374 & s1_ubar_obc(i,inorth,ir,lout)
1375# ifdef MASKING
1376 cff=cff*umask(i,jend+1)
1377# endif
1378 s1_ubar_obc(i,inorth,ir,lout)=cff
1379 END DO
1380 END IF
1381 END DO
1382 END IF
1383# endif
1384
1385# ifndef SOLVE3D
1386!
1387! 2D V-momentum component.
1388!
1389 DO j=jstrp,jendt
1390 DO i=istrt,iendt
1391 cff=vbar_std(i,j,rec)*vbar_std(i,j,rec)- &
1392 & s1_vbar(i,j,lout)
1393# ifdef MASKING
1394 cff=cff*vmask(i,j)
1395# endif
1396 s1_vbar(i,j,lout)=cff
1397 END DO
1398 END DO
1399# endif
1400
1401# ifdef ADJUST_BOUNDARY
1402!
1403! 2D V-momentum open boundaries.
1404!
1405 IF (any(lobc(:,isvbar,ng))) THEN
1406 DO ir=1,nbrec(ng)
1407 IF ((lobc(iwest,isvbar,ng)).and. &
1408 & domain(ng)%Western_Edge(tile)) THEN
1409 DO j=jstrv,jend
1410 cff=vbar_obc_std(j,iwest)*vbar_obc_std(j,iwest)- &
1411 & s1_vbar_obc(j,iwest,ir,lout)
1412# ifdef MASKING
1413 cff=cff*vmask(istr-1,j)
1414# endif
1415 s1_vbar_obc(j,iwest,ir,lout)=cff
1416 END DO
1417 END IF
1418 IF ((lobc(ieast,isvbar,ng)).and. &
1419 & domain(ng)%Eastern_Edge(tile)) THEN
1420 DO j=jstrv,jend
1421 cff=vbar_obc_std(j,ieast)*vbar_obc_std(j,ieast)- &
1422 & s1_vbar_obc(j,ieast,ir,lout)
1423# ifdef MASKING
1424 cff=cff*vmask(iend+1,j)
1425# endif
1426 s1_vbar_obc(j,ieast,ir,lout)=cff
1427 END DO
1428 END IF
1429 IF ((lobc(isouth,isvbar,ng)).and. &
1430 & domain(ng)%Southern_Edge(tile)) THEN
1431 DO i=istr,iend
1432 cff=vbar_obc_std(i,isouth)*vbar_obc_std(i,isouth)- &
1433 & s1_vbar_obc(i,isouth,ir,lout)
1434# ifdef MASKING
1435 cff=cff*vmask(i,jstr)
1436# endif
1437 s1_vbar_obc(i,isouth,ir,lout)=cff
1438 END DO
1439 END IF
1440 IF ((lobc(inorth,isvbar,ng)).and. &
1441 & domain(ng)%Northern_Edge(tile)) THEN
1442 DO i=istr,iend
1443 cff=vbar_obc_std(i,inorth)*vbar_obc_std(i,inorth)- &
1444 & s1_vbar_obc(i,inorth,ir,lout)
1445# ifdef MASKING
1446 cff=cff*vmask(i,jend+1)
1447# endif
1448 s1_vbar_obc(i,inorth,ir,lout)=cff
1449 END DO
1450 END IF
1451 END DO
1452 END IF
1453# endif
1454
1455# ifdef ADJUST_WSTRESS
1456!
1457! Surface momentum stress.
1458!
1459 DO ir=1,nfrec(ng)
1460 DO j=jstrt,jendt
1461 DO i=istrp,iendt
1462 cff=sustr_std(i,j)*sustr_std(i,j)- &
1463 & s1_sustr(i,j,ir,lout)
1464# ifdef MASKING
1465 cff=cff*umask(i,j)
1466# endif
1467 s1_sustr(i,j,ir,lout)=cff
1468 END DO
1469 END DO
1470 DO j=jstrp,jendt
1471 DO i=istrt,iendt
1472 cff=svstr_std(i,j)*svstr_std(i,j)- &
1473 & s1_svstr(i,j,ir,lout)
1474# ifdef MASKING
1475 cff=cff*vmask(i,j)
1476# endif
1477 s1_svstr(i,j,ir,lout)=cff
1478 END DO
1479 END DO
1480 END DO
1481# endif
1482
1483# ifdef SOLVE3D
1484!
1485! 3D U-momentum component.
1486!
1487 DO k=1,n(ng)
1488 DO j=jstrt,jendt
1489 DO i=istrp,iendt
1490 cff=u_std(i,j,k,rec)*u_std(i,j,k,rec)- &
1491 & s1_u(i,j,k,lout)
1492# ifdef MASKING
1493 cff=cff*umask(i,j)
1494# endif
1495 s1_u(i,j,k,lout)=cff
1496 END DO
1497 END DO
1498 END DO
1499
1500# ifdef ADJUST_BOUNDARY
1501!
1502! 3D U-momentum open boundaries.
1503!
1504 IF (any(lobc(:,isuvel,ng))) THEN
1505 DO ir=1,nbrec(ng)
1506 IF ((lobc(iwest,isuvel,ng)).and. &
1507 & domain(ng)%Western_Edge(tile)) THEN
1508 DO k=1,n(ng)
1509 DO j=jstr,jend
1510 cff=u_obc_std(j,k,iwest)*u_obc_std(j,k,iwest)- &
1511 & s1_u_obc(j,k,iwest,ir,lout)
1512# ifdef MASKING
1513 cff=cff*umask(istr,j)
1514# endif
1515 s1_u_obc(j,k,iwest,ir,lout)=cff
1516 END DO
1517 END DO
1518 END IF
1519 IF ((lobc(ieast,isuvel,ng)).and. &
1520 & domain(ng)%Eastern_Edge(tile)) THEN
1521 DO k=1,n(ng)
1522 DO j=jstr,jend
1523 cff=u_obc_std(j,k,ieast)*u_obc_std(j,k,ieast)- &
1524 & s1_u_obc(j,k,ieast,ir,lout)
1525# ifdef MASKING
1526 cff=cff*umask(iend+1,j)
1527# endif
1528 s1_u_obc(j,k,ieast,ir,lout)=cff
1529 END DO
1530 END DO
1531 END IF
1532 IF ((lobc(isouth,isuvel,ng)).and. &
1533 & domain(ng)%Southern_Edge(tile)) THEN
1534 DO k=1,n(ng)
1535 DO i=istru,iend
1536 cff=u_obc_std(i,k,isouth)*u_obc_std(i,k,isouth)- &
1537 & s1_u_obc(i,k,isouth,ir,lout)
1538# ifdef MASKING
1539 cff=cff*umask(i,jstr-1)
1540# endif
1541 s1_u_obc(i,k,isouth,ir,lout)=cff
1542 END DO
1543 END DO
1544 END IF
1545 IF ((lobc(inorth,isuvel,ng)).and. &
1546 & domain(ng)%Northern_Edge(tile)) THEN
1547 DO k=1,n(ng)
1548 DO i=istru,iend
1549 cff=u_obc_std(i,k,inorth)*u_obc_std(i,k,inorth)- &
1550 & s1_u_obc(i,k,inorth,ir,lout)
1551# ifdef MASKING
1552 cff=cff*umask(i,jend+1)
1553# endif
1554 s1_u_obc(i,k,inorth,ir,lout)=cff
1555 END DO
1556 END DO
1557 END IF
1558 END DO
1559 END IF
1560# endif
1561!
1562! 3D V-momentum component.
1563!
1564 DO k=1,n(ng)
1565 DO j=jstrp,jendt
1566 DO i=istrt,iendt
1567 cff=v_std(i,j,k,rec)*v_std(i,j,k,rec)- &
1568 & s1_v(i,j,k,lout)
1569# ifdef MASKING
1570 cff=cff*vmask(i,j)
1571# endif
1572 s1_v(i,j,k,lout)=cff
1573 END DO
1574 END DO
1575 END DO
1576
1577# ifdef ADJUST_BOUNDARY
1578!
1579! 3D V-momentum open boundaries.
1580!
1581 IF (any(lobc(:,isvvel,ng))) THEN
1582 DO ir=1,nbrec(ng)
1583 IF ((lobc(iwest,isvvel,ng)).and. &
1584 & domain(ng)%Western_Edge(tile)) THEN
1585 DO k=1,n(ng)
1586 DO j=jstrv,jend
1587 cff=v_obc_std(j,k,iwest)*v_obc_std(j,k,iwest)- &
1588 & s1_v_obc(j,k,iwest,ir,lout)
1589# ifdef MASKING
1590 cff=cff*vmask(istr-1,j)
1591# endif
1592 s1_v_obc(j,k,iwest,ir,lout)=cff
1593 END DO
1594 END DO
1595 END IF
1596 IF ((lobc(ieast,isvvel,ng)).and. &
1597 & domain(ng)%Eastern_Edge(tile)) THEN
1598 DO k=1,n(ng)
1599 DO j=jstrv,jend
1600 cff=v_obc_std(j,k,ieast)*v_obc_std(j,k,ieast)- &
1601 & s1_v_obc(j,k,ieast,ir,lout)
1602# ifdef MASKING
1603 cff=cff*vmask(iend+1,j)
1604# endif
1605 s1_v_obc(j,k,ieast,ir,lout)=cff
1606 END DO
1607 END DO
1608 END IF
1609 IF ((lobc(isouth,isvvel,ng)).and. &
1610 & domain(ng)%Southern_Edge(tile)) THEN
1611 DO k=1,n(ng)
1612 DO i=istr,iend
1613 cff=v_obc_std(i,k,isouth)*v_obc_std(i,k,isouth)- &
1614 & s1_v_obc(i,k,isouth,ir,lout)
1615# ifdef MASKING
1616 cff=cff*vmask(i,jstr)
1617# endif
1618 s1_v_obc(i,k,isouth,ir,lout)=cff
1619 END DO
1620 END DO
1621 END IF
1622 IF ((lobc(inorth,isvvel,ng)).and. &
1623 & domain(ng)%Northern_Edge(tile)) THEN
1624 DO k=1,n(ng)
1625 DO i=istr,iend
1626 cff=v_obc_std(i,k,inorth)*v_obc_std(i,k,inorth)- &
1627 & s1_v_obc(i,k,inorth,ir,lout)
1628# ifdef MASKING
1629 cff=cff*vmask(i,jend+1)
1630# endif
1631 s1_v_obc(i,k,inorth,ir,lout)=cff
1632 END DO
1633 END DO
1634 END IF
1635 END DO
1636 END IF
1637# endif
1638!
1639! Tracers.
1640!
1641 DO it=1,nt(ng)
1642 DO k=1,n(ng)
1643 DO j=jstrt,jendt
1644 DO i=istrt,iendt
1645 cff=t_std(i,j,k,rec,it)*t_std(i,j,k,rec,it)- &
1646 & s1_t(i,j,k,lout,it)
1647# ifdef MASKING
1648 cff=cff*rmask(i,j)
1649# endif
1650 s1_t(i,j,k,lout,it)=cff
1651 END DO
1652 END DO
1653 END DO
1654 END DO
1655
1656# ifdef ADJUST_BOUNDARY
1657!
1658! Tracers open boundaries.
1659!
1660 DO it=1,nt(ng)
1661 IF (any(lobc(:,istvar(it),ng))) THEN
1662 DO ir=1,nbrec(ng)
1663 IF ((lobc(iwest,istvar(it),ng)).and. &
1664 & domain(ng)%Western_Edge(tile)) THEN
1665 DO k=1,n(ng)
1666 DO j=jstr,jend
1667 cff=t_obc_std(j,k,iwest,it)*t_obc_std(j,k,iwest,it)- &
1668 & s1_t_obc(j,k,iwest,ir,lout,it)
1669# ifdef MASKING
1670 cff=cff*rmask(istr-1,j)
1671# endif
1672 s1_t_obc(j,k,iwest,ir,lout,it)=cff
1673 END DO
1674 END DO
1675 END IF
1676 IF ((lobc(ieast,istvar(it),ng)).and. &
1677 & domain(ng)%Eastern_Edge(tile)) THEN
1678 DO k=1,n(ng)
1679 DO j=jstr,jend
1680 cff=t_obc_std(j,k,ieast,it)*t_obc_std(j,k,ieast,it)- &
1681 & s1_t_obc(j,k,ieast,ir,lout,it)
1682# ifdef MASKING
1683 cff=cff*rmask(iend+1,j)
1684# endif
1685 s1_t_obc(j,k,ieast,ir,lout,it)=cff
1686 END DO
1687 END DO
1688 END IF
1689 IF ((lobc(isouth,istvar(it),ng)).and. &
1690 & domain(ng)%Southern_Edge(tile)) THEN
1691 DO k=1,n(ng)
1692 DO i=istr,iend
1693 cff=t_obc_std(i,k,isouth,it)*t_obc_std(i,k,isouth,it)-&
1694 & s1_t_obc(i,k,isouth,ir,lout,it)
1695# ifdef MASKING
1696 cff=cff*rmask(i,jstr-1)
1697# endif
1698 s1_t_obc(i,k,isouth,ir,lout,it)=cff
1699 END DO
1700 END DO
1701 END IF
1702 IF ((lobc(inorth,istvar(it),ng)).and. &
1703 & domain(ng)%Northern_Edge(tile)) THEN
1704 DO k=1,n(ng)
1705 DO i=istr,iend
1706 cff=t_obc_std(i,k,inorth,it)*t_obc_std(i,k,inorth,it)-&
1707 & s1_t_obc(i,k,inorth,ir,lout,it)
1708# ifdef MASKING
1709 cff=cff*rmask(i,jend+1)
1710# endif
1711 s1_t_obc(i,k,inorth,ir,lout,it)=cff
1712 END DO
1713 END DO
1714 END IF
1715 END DO
1716 END IF
1717 END DO
1718# endif
1719
1720# ifdef ADJUST_STFLUX
1721!
1722! Surface tracers flux.
1723!
1724 DO it=1,nt(ng)
1725 IF (lstflux(it,ng)) THEN
1726 DO ir=1,nfrec(ng)
1727 DO j=jstrt,jendt
1728 DO i=istrt,iendt
1729 cff=stflx_std(i,j,it)*stflx_std(i,j,it)- &
1730 & s1_tflux(i,j,ir,lout,it)
1731# ifdef MASKING
1732 cff=cff*rmask(i,j)
1733# endif
1734 s1_tflux(i,j,ir,lout,it)=cff
1735 END DO
1736 END DO
1737 END DO
1738 END IF
1739 END DO
1740# endif
1741
1742# endif
1743!
1744 RETURN
1745 END SUBROUTINE analysis_var
1746#endif
1747 END MODULE posterior_var_mod
type(t_boundary), dimension(:), allocatable boundary
type(t_forces), dimension(:), allocatable forces
Definition mod_forces.F:554
real(dp), dimension(:,:), allocatable cg_beta
real(dp), dimension(:,:), allocatable cg_delta
real(r8), dimension(:,:), allocatable zlanczos_inv
real(r8), dimension(:,:), allocatable zlanczos_err
real(r8), dimension(:,:), allocatable zlanczos_coef
type(t_grid), dimension(:), allocatable grid
Definition mod_grid.F:365
character(len=256) calledfrom
type(t_io), dimension(:), allocatable hss
integer stdout
character(len=256) sourcefile
integer isvvel
integer isvbar
integer, dimension(:), allocatable istvar
integer isuvel
integer isfsur
integer isubar
type(t_ocean), dimension(:), allocatable ocean
Definition mod_ocean.F:351
logical master
type(t_domain), dimension(:), allocatable domain
Definition mod_param.F:329
integer ninner
logical, dimension(:,:,:), allocatable lobc
integer, parameter iwest
logical, dimension(:,:), allocatable lstflux
integer exit_flag
integer, parameter isouth
integer, parameter ieast
integer, parameter inorth
integer, dimension(:), allocatable lold
integer, dimension(:), allocatable lnew
subroutine analysis_var(ng, tile, model, lbi, ubi, lbj, ubj, lbij, ubij, lout, lweak, rmask, umask, vmask, 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, zeta_std, s1_t_obc, s1_u_obc, s1_v_obc, s1_ubar_obc, s1_vbar_obc, s1_zeta_obc, s1_sustr, s1_svstr, s1_tflux, s1_t, s1_u, s1_v, s1_zeta)
subroutine, public posterior_var(ng, tile, model, outloop)
subroutine posterior_var_tile(ng, tile, model, lbi, ubi, lbj, ubj, lbij, ubij, imins, imaxs, jmins, jmaxs, lold, lnew, outloop, rmask, umask, vmask, 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, nl_t_obc, nl_u_obc, nl_v_obc, nl_ubar_obc, nl_vbar_obc, nl_zeta_obc, nl_ustr, nl_vstr, nl_tflux, nl_t, nl_u, nl_v, nl_ubar, nl_vbar, nl_zeta, 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, 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_ubar, ad_vbar, ad_zeta)
subroutine, public state_addition(ng, tile, lbi, ubi, lbj, ubj, lbij, ubij, lin1, lin2, lout, fac1, fac2, rmask, umask, vmask, s1_t_obc, s2_t_obc, s1_u_obc, s2_u_obc, s1_v_obc, s2_v_obc, s1_ubar_obc, s2_ubar_obc, s1_vbar_obc, s2_vbar_obc, s1_zeta_obc, s2_zeta_obc, s1_sustr, s2_sustr, s1_svstr, s2_svstr, s1_tflux, s2_tflux, s1_t, s2_t, s1_u, s2_u, s1_v, s2_v, s1_ubar, s2_ubar, s1_vbar, s2_vbar, s1_zeta, s2_zeta)
subroutine, public state_copy(ng, tile, lbi, ubi, lbj, ubj, lbij, ubij, linp, lout, s1_t_obc, s2_t_obc, s1_u_obc, s2_u_obc, s1_v_obc, s2_v_obc, s1_ubar_obc, s2_ubar_obc, s1_vbar_obc, s2_vbar_obc, s1_zeta_obc, s2_zeta_obc, s1_sustr, s2_sustr, s1_svstr, s2_svstr, s1_tflux, s2_tflux, s1_t, s2_t, s1_u, s2_u, s1_v, s2_v, s1_ubar, s2_ubar, s1_vbar, s2_vbar, s1_zeta, s2_zeta)
Definition state_copy.F:57
subroutine, public state_initialize(ng, tile, lbi, ubi, lbj, ubj, lbij, ubij, lout, fac, rmask, umask, vmask, s_t_obc, s_u_obc, s_v_obc, s_ubar_obc, s_vbar_obc, s_zeta_obc, s_sustr, s_svstr, s_tflux, s_t, s_u, s_v, s_zeta)
subroutine, public state_product(ng, tile, model, lbi, ubi, lbj, ubj, lbij, ubij, rmask, umask, vmask, s1_t_obc, s2_t_obc, s3_t_obc, s1_u_obc, s2_u_obc, s3_u_obc, s1_v_obc, s2_v_obc, s3_v_obc, s1_ubar_obc, s2_ubar_obc, s3_ubar_obc, s1_vbar_obc, s2_vbar_obc, s3_vbar_obc, s1_zeta_obc, s2_zeta_obc, s3_zeta_obc, s1_sustr, s2_sustr, s3_sustr, s1_svstr, s2_svstr, s3_svstr, s1_tflux, s2_tflux, s3_tflux, s1_t, s2_t, s3_t, s1_u, s2_u, s3_u, s1_v, s2_v, s3_v, s1_zeta, s2_zeta, s3_zeta)
subroutine, public state_read(ng, tile, model, iotype, lbi, ubi, lbj, ubj, lbij, ubij, lout, rec, nopen, ncid, piofile, ncname, rmask, umask, vmask, s_t_obc, s_u_obc, s_v_obc, s_ubar_obc, s_vbar_obc, s_zeta_obc, s_ustr, s_vstr, s_tflux, s_t, s_u, s_v, s_zeta)
Definition state_read.F:137
recursive subroutine wclock_off(ng, model, region, line, routine)
Definition timers.F:148
recursive subroutine wclock_on(ng, model, region, line, routine)
Definition timers.F:3