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

Functions/Subroutines

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)
 

Function/Subroutine Documentation

◆ state_addition()

subroutine, public state_addition_mod::state_addition ( 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) lin1,
integer, intent(in) lin2,
integer, intent(in) lout,
real(r8), intent(in) fac1,
real(r8), intent(in) fac2,
real(r8), dimension(lbi:,lbj:), intent(in) rmask,
real(r8), dimension(lbi:,lbj:), intent(in) umask,
real(r8), dimension(lbi:,lbj:), intent(in) vmask,
real(r8), dimension(lbij:,:,:,:,:,:), intent(inout) s1_t_obc,
real(r8), dimension(lbij:,:,:,:,:,:), intent(in) s2_t_obc,
real(r8), dimension(lbij:,:,:,:,:), intent(inout) s1_u_obc,
real(r8), dimension(lbij:,:,:,:,:), intent(in) s2_u_obc,
real(r8), dimension(lbij:,:,:,:,:), intent(inout) s1_v_obc,
real(r8), dimension(lbij:,:,:,:,:), intent(in) s2_v_obc,
real(r8), dimension(lbij:,:,:,:), intent(inout) s1_ubar_obc,
real(r8), dimension(lbij:,:,:,:), intent(in) s2_ubar_obc,
real(r8), dimension(lbij:,:,:,:), intent(inout) s1_vbar_obc,
real(r8), dimension(lbij:,:,:,:), intent(in) s2_vbar_obc,
real(r8), dimension(lbij:,:,:,:), intent(inout) s1_zeta_obc,
real(r8), dimension(lbij:,:,:,:), intent(in) s2_zeta_obc,
real(r8), dimension(lbi:,lbj:,:,:), intent(inout) s1_sustr,
real(r8), dimension(lbi:,lbj:,:,:), intent(in) s2_sustr,
real(r8), dimension(lbi:,lbj:,:,:), intent(inout) s1_svstr,
real(r8), dimension(lbi:,lbj:,:,:), intent(in) s2_svstr,
real(r8), dimension(lbi:,lbj:,:,:,:), intent(inout) s1_tflux,
real(r8), dimension(lbi:,lbj:,:,:,:), intent(in) s2_tflux,
real(r8), dimension(lbi:,lbj:,:,:,:), intent(inout) s1_t,
real(r8), dimension(lbi:,lbj:,:,:,:), intent(in) s2_t,
real(r8), dimension(lbi:,lbj:,:,:), intent(inout) s1_u,
real(r8), dimension(lbi:,lbj:,:,:), intent(in) s2_u,
real(r8), dimension(lbi:,lbj:,:,:), intent(inout) s1_v,
real(r8), dimension(lbi:,lbj:,:,:), intent(in) s2_v,
real(r8), dimension(lbi:,lbj:,:), intent(inout) s1_ubar,
real(r8), dimension(lbi:,lbj:,:), intent(in) s2_ubar,
real(r8), dimension(lbi:,lbj:,:), intent(inout) s1_vbar,
real(r8), dimension(lbi:,lbj:,:), intent(in) s2_vbar,
real(r8), dimension(lbi:,lbj:,:), intent(inout) s1_zeta,
real(r8), dimension(lbi:,lbj:,:), intent(in) s2_zeta )

Definition at line 27 of file state_addition.F.

64!***********************************************************************
65!
66 USE mod_param
67#if defined ADJUST_BOUNDARY || defined ADJUST_STFLUX || \
68 defined adjust_wstress
69 USE mod_ncparam
70 USE mod_scalars
71#endif
72!
73! Imported variable declarations.
74!
75 integer, intent(in) :: ng, tile
76 integer, intent(in) :: LBi, UBi, LBj, UBj, LBij, UBij
77 integer, intent(in) :: Lin1, Lin2, Lout
78!
79 real(r8), intent(in) :: fac1, fac2
80!
81#ifdef ASSUMED_SHAPE
82# ifdef MASKING
83 real(r8), intent(in) :: rmask(LBi:,LBj:)
84 real(r8), intent(in) :: umask(LBi:,LBj:)
85 real(r8), intent(in) :: vmask(LBi:,LBj:)
86# endif
87# ifdef ADJUST_BOUNDARY
88# ifdef SOLVE3D
89 real(r8), intent(in) :: s2_t_obc(LBij:,:,:,:,:,:)
90 real(r8), intent(in) :: s2_u_obc(LBij:,:,:,:,:)
91 real(r8), intent(in) :: s2_v_obc(LBij:,:,:,:,:)
92# endif
93 real(r8), intent(in) :: s2_ubar_obc(LBij:,:,:,:)
94 real(r8), intent(in) :: s2_vbar_obc(LBij:,:,:,:)
95 real(r8), intent(in) :: s2_zeta_obc(LBij:,:,:,:)
96# endif
97# ifdef ADJUST_WSTRESS
98 real(r8), intent(in) :: s2_sustr(LBi:,LBj:,:,:)
99 real(r8), intent(in) :: s2_svstr(LBi:,LBj:,:,:)
100# endif
101# ifdef SOLVE3D
102# ifdef ADJUST_STFLUX
103 real(r8), intent(in) :: s2_tflux(LBi:,LBj:,:,:,:)
104# endif
105 real(r8), intent(in) :: s2_t(LBi:,LBj:,:,:,:)
106 real(r8), intent(in) :: s2_u(LBi:,LBj:,:,:)
107 real(r8), intent(in) :: s2_v(LBi:,LBj:,:,:)
108# if defined WEAK_CONSTRAINT && defined TIME_CONV
109 real(r8), intent(in) :: s2_ubar(LBi:,LBj:,:)
110 real(r8), intent(in) :: s2_vbar(LBi:,LBj:,:)
111# endif
112# else
113 real(r8), intent(in) :: s2_ubar(LBi:,LBj:,:)
114 real(r8), intent(in) :: s2_vbar(LBi:,LBj:,:)
115# endif
116 real(r8), intent(in) :: s2_zeta(LBi:,LBj:,:)
117
118# ifdef ADJUST_BOUNDARY
119# ifdef SOLVE3D
120 real(r8), intent(inout) :: s1_t_obc(LBij:,:,:,:,:,:)
121 real(r8), intent(inout) :: s1_u_obc(LBij:,:,:,:,:)
122 real(r8), intent(inout) :: s1_v_obc(LBij:,:,:,:,:)
123# endif
124 real(r8), intent(inout) :: s1_ubar_obc(LBij:,:,:,:)
125 real(r8), intent(inout) :: s1_vbar_obc(LBij:,:,:,:)
126 real(r8), intent(inout) :: s1_zeta_obc(LBij:,:,:,:)
127# endif
128# ifdef ADJUST_WSTRESS
129 real(r8), intent(inout) :: s1_sustr(LBi:,LBj:,:,:)
130 real(r8), intent(inout) :: s1_svstr(LBi:,LBj:,:,:)
131# endif
132# ifdef SOLVE3D
133# ifdef ADJUST_STFLUX
134 real(r8), intent(inout) :: s1_tflux(LBi:,LBj:,:,:,:)
135# endif
136 real(r8), intent(inout) :: s1_t(LBi:,LBj:,:,:,:)
137 real(r8), intent(inout) :: s1_u(LBi:,LBj:,:,:)
138 real(r8), intent(inout) :: s1_v(LBi:,LBj:,:,:)
139# if defined WEAK_CONSTRAINT && defined TIME_CONV
140 real(r8), intent(inout) :: s1_ubar(LBi:,LBj:,:)
141 real(r8), intent(inout) :: s1_vbar(LBi:,LBj:,:)
142# endif
143# else
144 real(r8), intent(inout) :: s1_ubar(LBi:,LBj:,:)
145 real(r8), intent(inout) :: s1_vbar(LBi:,LBj:,:)
146# endif
147 real(r8), intent(inout) :: s1_zeta(LBi:,LBj:,:)
148
149#else
150
151# ifdef MASKING
152 real(r8), intent(in) :: rmask(LBi:UBi,LBj:UBj)
153 real(r8), intent(in) :: umask(LBi:UBi,LBj:UBj)
154 real(r8), intent(in) :: vmask(LBi:UBi,LBj:UBj)
155# endif
156# ifdef ADJUST_BOUNDARY
157# ifdef SOLVE3D
158 real(r8), intent(in) :: s2_t_obc(LBij:UBij,N(ng),4, &
159 & Nbrec(ng),2,NT(ng))
160 real(r8), intent(in) :: s2_u_obc(LBij:UBij,N(ng),4,Nbrec(ng),2)
161 real(r8), intent(in) :: s2_v_obc(LBij:UBij,N(ng),4,Nbrec(ng),2)
162# endif
163 real(r8), intent(in) :: s2_ubar_obc(LBij:UBij,4,Nbrec(ng),2)
164 real(r8), intent(in) :: s2_vbar_obc(LBij:UBij,4,Nbrec(ng),2)
165 real(r8), intent(in) :: s2_zeta_obc(LBij:UBij,4,Nbrec(ng),2)
166# endif
167# ifdef ADJUST_WSTRESS
168 real(r8), intent(in) :: s2_sustr(LBi:UBi,LBj:UBj,Nfrec(ng),2)
169 real(r8), intent(in) :: s2_svstr(LBi:UBi,LBj:UBj,Nfrec(ng),2)
170# endif
171# ifdef SOLVE3D
172# ifdef ADJUST_STFLUX
173 real(r8), intent(in) :: s2_tflux(LBi:UBi,LBj:UBj, &
174 & Nfrec(ng),2,NT(ng))
175# endif
176 real(r8), intent(in) :: s2_t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
177 real(r8), intent(in) :: s2_u(LBi:UBi,LBj:UBj,N(ng),2)
178 real(r8), intent(in) :: s2_v(LBi:UBi,LBj:UBj,N(ng),2)
179# if defined WEAK_CONSTRAINT && defined TIME_CONV
180 real(r8), intent(in) :: s2_ubar(LBi:UBi,LBj:UBj,:)
181 real(r8), intent(in) :: s2_vbar(LBi:UBi,LBj:UBj,:)
182# endif
183# else
184 real(r8), intent(in) :: s2_ubar(LBi:UBi,LBj:UBj,:)
185 real(r8), intent(in) :: s2_vbar(LBi:UBi,LBj:UBj,:)
186# endif
187 real(r8), intent(in) :: s2_zeta(LBi:UBi,LBj:UBj,:)
188
189# ifdef ADJUST_BOUNDARY
190# ifdef SOLVE3D
191 real(r8), intent(inout) :: s1_t_obc(LBij:UBij,N(ng),4, &
192 & Nbrec(ng),2,NT(ng))
193 real(r8), intent(inout) :: s1_u_obc(LBij:UBij,N(ng),4,Nbrec(ng),2)
194 real(r8), intent(inout) :: s1_v_obc(LBij:UBij,N(ng),4,Nbrec(ng),2)
195# endif
196 real(r8), intent(inout) :: s1_ubar_obc(LBij:UBij,4,Nbrec(ng),2)
197 real(r8), intent(inout) :: s1_vbar_obc(LBij:UBij,4,Nbrec(ng),2)
198 real(r8), intent(inout) :: s1_zeta_obc(LBij:UBij,4,Nbrec(ng),2)
199# endif
200# ifdef ADJUST_WSTRESS
201 real(r8), intent(inout) :: s1_sustr(LBi:UBi,LBj:UBj,Nfrec(ng),2)
202 real(r8), intent(inout) :: s1_svstr(LBi:UBi,LBj:UBj,Nfrec(ng),2)
203# endif
204# ifdef SOLVE3D
205# ifdef ADJUST_STFLUX
206 real(r8), intent(inout) :: s1_tflux(LBi:UBi,LBj:UBj, &
207 & Nfrec(ng),2,NT(ng))
208# endif
209 real(r8), intent(inout) :: s1_t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
210 real(r8), intent(inout) :: s1_u(LBi:UBi,LBj:UBj,N(ng),2)
211 real(r8), intent(inout) :: s1_v(LBi:UBi,LBj:UBj,N(ng),2)
212# if defined WEAK_CONSTRAINT && defined TIME_CONV
213 real(r8), intent(inout) :: s1_ubar(LBi:UBi,LBj:UBj,:)
214 real(r8), intent(inout) :: s1_vbar(LBi:UBi,LBj:UBj,:)
215# endif
216# else
217 real(r8), intent(inout) :: s1_ubar(LBi:UBi,LBj:UBj,:)
218 real(r8), intent(inout) :: s1_vbar(LBi:UBi,LBj:UBj,:)
219# endif
220 real(r8), intent(inout) :: s1_zeta(LBi:UBi,LBj:UBj,:)
221#endif
222!
223! Local variable declarations.
224!
225 integer :: i, j, k
226 integer :: ib, ir, it
227
228#include "set_bounds.h"
229!
230!-----------------------------------------------------------------------
231! Compute the following operation between S1 and S2 model state
232! trajectories:
233! S1(Lout) = fac1 * S1(Lin1) + fac2 * S2(Lin2)
234!-----------------------------------------------------------------------
235!
236! Free-surface.
237!
238 DO j=jstrt,jendt
239 DO i=istrt,iendt
240 s1_zeta(i,j,lout)=fac1*s1_zeta(i,j,lin1)+ &
241 & fac2*s2_zeta(i,j,lin2)
242#ifdef MASKING
243 s1_zeta(i,j,lout)=s1_zeta(i,j,lout)*rmask(i,j)
244#endif
245 END DO
246 END DO
247
248#ifdef ADJUST_BOUNDARY
249!
250! Free-surface open boundaries.
251!
252 IF (any(lobc(:,isfsur,ng))) THEN
253 DO ir=1,nbrec(ng)
254 IF ((lobc(iwest,isfsur,ng)).and. &
255 & domain(ng)%Western_Edge(tile)) THEN
256 ib=iwest
257 DO j=jstr,jend
258 s1_zeta_obc(j,ib,ir,lout)=fac1*s1_zeta_obc(j,ib,ir,lin1)+ &
259 & fac2*s2_zeta_obc(j,ib,ir,lin2)
260# ifdef MASKING
261 s1_zeta_obc(j,ib,ir,lout)=s1_zeta_obc(j,ib,ir,lout)* &
262 & rmask(istr-1,j)
263# endif
264 END DO
265 END IF
266 IF ((lobc(ieast,isfsur,ng)).and. &
267 & domain(ng)%Eastern_Edge(tile)) THEN
268 ib=ieast
269 DO j=jstr,jend
270 s1_zeta_obc(j,ib,ir,lout)=fac1*s1_zeta_obc(j,ib,ir,lin1)+ &
271 & fac2*s2_zeta_obc(j,ib,ir,lin2)
272# ifdef MASKING
273 s1_zeta_obc(j,ib,ir,lout)=s1_zeta_obc(j,ib,ir,lout)* &
274 & rmask(iend+1,j)
275# endif
276 END DO
277 END IF
278 IF ((lobc(isouth,isfsur,ng)).and. &
279 & domain(ng)%Southern_Edge(tile)) THEN
280 ib=isouth
281 DO i=istr,iend
282 s1_zeta_obc(i,ib,ir,lout)=fac1*s1_zeta_obc(i,ib,ir,lin1)+ &
283 & fac2*s2_zeta_obc(i,ib,ir,lin2)
284# ifdef MASKING
285 s1_zeta_obc(i,ib,ir,lout)=s1_zeta_obc(i,ib,ir,lout)* &
286 & rmask(i,jstr-1)
287# endif
288 END DO
289 END IF
290 IF ((lobc(inorth,isfsur,ng)).and. &
291 & domain(ng)%Northern_Edge(tile)) THEN
292 ib=inorth
293 DO i=istr,iend
294 s1_zeta_obc(i,ib,ir,lout)=fac1*s1_zeta_obc(i,ib,ir,lin1)+ &
295 & fac2*s2_zeta_obc(i,ib,ir,lin2)
296# ifdef MASKING
297 s1_zeta_obc(i,ib,ir,lout)=s1_zeta_obc(i,ib,ir,lout)* &
298 & rmask(i,jend+1)
299# endif
300 END DO
301 END IF
302 END DO
303 END IF
304#endif
305
306#if !defined SOLVE3D || (defined WEAK_CONSTRAINT && defined TIME_CONV )
307!
308! 2D U-momentum component.
309!
310 DO j=jstrt,jendt
311 DO i=istrp,iendt
312 s1_ubar(i,j,lout)=fac1*s1_ubar(i,j,lin1)+ &
313 & fac2*s2_ubar(i,j,lin2)
314# ifdef MASKING
315 s1_ubar(i,j,lout)=s1_ubar(i,j,lout)*umask(i,j)
316# endif
317 END DO
318 END DO
319#endif
320
321#ifdef ADJUST_BOUNDARY
322!
323! 2D U-momentum open boundaries.
324!
325 IF (any(lobc(:,isubar,ng))) THEN
326 DO ir=1,nbrec(ng)
327 IF ((lobc(iwest,isubar,ng)).and. &
328 & domain(ng)%Western_Edge(tile)) THEN
329 ib=iwest
330 DO j=jstr,jend
331 s1_ubar_obc(j,ib,ir,lout)=fac1*s1_ubar_obc(j,ib,ir,lin1)+ &
332 & fac2*s2_ubar_obc(j,ib,ir,lin2)
333# ifdef MASKING
334 s1_ubar_obc(j,ib,ir,lout)=s1_ubar_obc(j,ib,ir,lout)* &
335 & umask(istr,j)
336# endif
337 END DO
338 END IF
339 IF ((lobc(ieast,isubar,ng)).and. &
340 & domain(ng)%Eastern_Edge(tile)) THEN
341 ib=ieast
342 DO j=jstr,jend
343 s1_ubar_obc(j,ib,ir,lout)=fac1*s1_ubar_obc(j,ib,ir,lin1)+ &
344 & fac2*s2_ubar_obc(j,ib,ir,lin2)
345# ifdef MASKING
346 s1_ubar_obc(j,ib,ir,lout)=s1_ubar_obc(j,ib,ir,lout)* &
347 & umask(iend+1,j)
348# endif
349 END DO
350 END IF
351 IF ((lobc(isouth,isubar,ng)).and. &
352 & domain(ng)%Southern_Edge(tile)) THEN
353 ib=isouth
354 DO i=istru,iend
355 s1_ubar_obc(i,ib,ir,lout)=fac1*s1_ubar_obc(i,ib,ir,lin1)+ &
356 & fac2*s2_ubar_obc(i,ib,ir,lin2)
357# ifdef MASKING
358 s1_ubar_obc(i,ib,ir,lout)=s1_ubar_obc(i,ib,ir,lout)* &
359 & umask(i,jstr-1)
360# endif
361 END DO
362 END IF
363 IF ((lobc(inorth,isubar,ng)).and. &
364 & domain(ng)%Northern_Edge(tile)) THEN
365 ib=inorth
366 DO i=istru,iend
367 s1_ubar_obc(i,ib,ir,lout)=fac1*s1_ubar_obc(i,ib,ir,lin1)+ &
368 & fac2*s2_ubar_obc(i,ib,ir,lin2)
369# ifdef MASKING
370 s1_ubar_obc(i,ib,ir,lout)=s1_ubar_obc(i,ib,ir,lout)* &
371 & umask(i,jend+1)
372# endif
373 END DO
374 END IF
375 END DO
376 END IF
377#endif
378
379#if !defined SOLVE3D || (defined WEAK_CONSTRAINT && defined TIME_CONV )
380!
381! 2D V-momentum component.
382!
383 DO j=jstrp,jendt
384 DO i=istrt,iendt
385 s1_vbar(i,j,lout)=fac1*s1_vbar(i,j,lin1)+ &
386 & fac2*s2_vbar(i,j,lin2)
387# ifdef MASKING
388 s1_vbar(i,j,lout)=s1_vbar(i,j,lout)*vmask(i,j)
389# endif
390 END DO
391 END DO
392#endif
393
394#ifdef ADJUST_BOUNDARY
395!
396! 2D V-momentum open boundaries.
397!
398 IF (any(lobc(:,isvbar,ng))) THEN
399 DO ir=1,nbrec(ng)
400 IF ((lobc(iwest,isvbar,ng)).and. &
401 & domain(ng)%Western_Edge(tile)) THEN
402 ib=iwest
403 DO j=jstrv,jend
404 s1_vbar_obc(j,ib,ir,lout)=fac1*s1_vbar_obc(j,ib,ir,lin1)+ &
405 & fac2*s2_vbar_obc(j,ib,ir,lin2)
406# ifdef MASKING
407 s1_vbar_obc(j,ib,ir,lout)=s1_vbar_obc(j,ib,ir,lout)* &
408 & vmask(istr-1,j)
409# endif
410 END DO
411 END IF
412 IF ((lobc(ieast,isvbar,ng)).and. &
413 & domain(ng)%Eastern_Edge(tile)) THEN
414 ib=ieast
415 DO j=jstrv,jend
416 s1_vbar_obc(j,ib,ir,lout)=fac1*s1_vbar_obc(j,ib,ir,lin1)+ &
417 & fac2*s2_vbar_obc(j,ib,ir,lin2)
418# ifdef MASKING
419 s1_vbar_obc(j,ib,ir,lout)=s1_vbar_obc(j,ib,ir,lout)* &
420 & vmask(iend+1,j)
421# endif
422 END DO
423 END IF
424 IF ((lobc(isouth,isvbar,ng)).and. &
425 & domain(ng)%Southern_Edge(tile)) THEN
426 ib=isouth
427 DO i=istr,iend
428 s1_vbar_obc(i,ib,ir,lout)=fac1*s1_vbar_obc(i,ib,ir,lin1)+ &
429 & fac2*s2_vbar_obc(i,ib,ir,lin2)
430# ifdef MASKING
431 s1_vbar_obc(i,ib,ir,lout)=s1_vbar_obc(i,ib,ir,lout)* &
432 & vmask(i,jstr)
433# endif
434 END DO
435 END IF
436 IF ((lobc(inorth,isvbar,ng)).and. &
437 & domain(ng)%Northern_Edge(tile)) THEN
438 ib=inorth
439 DO i=istr,iend
440 s1_vbar_obc(i,ib,ir,lout)=fac1*s1_vbar_obc(i,ib,ir,lin1)+ &
441 & fac2*s2_vbar_obc(i,ib,ir,lin2)
442# ifdef MASKING
443 s1_vbar_obc(i,ib,ir,lout)=s1_vbar_obc(i,ib,ir,lout)* &
444 & vmask(i,jend+1)
445# endif
446 END DO
447 END IF
448 END DO
449 END IF
450#endif
451
452#ifdef ADJUST_WSTRESS
453!
454! Surface momentum stress.
455!
456 DO ir=1,nfrec(ng)
457 DO j=jstrt,jendt
458 DO i=istrp,iendt
459 s1_sustr(i,j,ir,lout)=fac1*s1_sustr(i,j,ir,lin1)+ &
460 & fac2*s2_sustr(i,j,ir,lin2)
461# ifdef MASKING
462 s1_sustr(i,j,ir,lout)=s1_sustr(i,j,ir,lout)*umask(i,j)
463# endif
464 END DO
465 END DO
466 DO j=jstrp,jendt
467 DO i=istrt,iendt
468 s1_svstr(i,j,ir,lout)=fac1*s1_svstr(i,j,ir,lin1)+ &
469 & fac2*s2_svstr(i,j,ir,lin2)
470# ifdef MASKING
471 s1_svstr(i,j,ir,lout)=s1_svstr(i,j,ir,lout)*vmask(i,j)
472# endif
473 END DO
474 END DO
475 END DO
476#endif
477
478#ifdef SOLVE3D
479!
480! 3D U-momentum component.
481!
482 DO k=1,n(ng)
483 DO j=jstrt,jendt
484 DO i=istrp,iendt
485 s1_u(i,j,k,lout)=fac1*s1_u(i,j,k,lin1)+ &
486 & fac2*s2_u(i,j,k,lin2)
487# ifdef MASKING
488 s1_u(i,j,k,lout)=s1_u(i,j,k,lout)*umask(i,j)
489# endif
490 END DO
491 END DO
492 END DO
493
494# ifdef ADJUST_BOUNDARY
495!
496! 3D U-momentum open boundaries.
497!
498 IF (any(lobc(:,isuvel,ng))) THEN
499 DO ir=1,nbrec(ng)
500 IF ((lobc(iwest,isuvel,ng)).and. &
501 & domain(ng)%Western_Edge(tile)) THEN
502 ib=iwest
503 DO k=1,n(ng)
504 DO j=jstr,jend
505 s1_u_obc(j,k,ib,ir,lout)=fac1*s1_u_obc(j,k,ib,ir,lin1)+ &
506 & fac2*s2_u_obc(j,k,ib,ir,lin2)
507# ifdef MASKING
508 s1_u_obc(j,k,ib,ir,lout)=s1_u_obc(j,k,ib,ir,lout)* &
509 & umask(istr,j)
510# endif
511 END DO
512 END DO
513 END IF
514 IF ((lobc(ieast,isuvel,ng)).and. &
515 & domain(ng)%Eastern_Edge(tile)) THEN
516 ib=ieast
517 DO k=1,n(ng)
518 DO j=jstr,jend
519 s1_u_obc(j,k,ib,ir,lout)=fac1*s1_u_obc(j,k,ib,ir,lin1)+ &
520 & fac2*s2_u_obc(j,k,ib,ir,lin2)
521# ifdef MASKING
522 s1_u_obc(j,k,ib,ir,lout)=s1_u_obc(j,k,ib,ir,lout)* &
523 & umask(iend+1,j)
524# endif
525 END DO
526 END DO
527 END IF
528 IF ((lobc(isouth,isuvel,ng)).and. &
529 & domain(ng)%Southern_Edge(tile)) THEN
530 ib=isouth
531 DO k=1,n(ng)
532 DO i=istru,iend
533 s1_u_obc(i,k,ib,ir,lout)=fac1*s1_u_obc(i,k,ib,ir,lin1)+ &
534 & fac2*s2_u_obc(i,k,ib,ir,lin2)
535# ifdef MASKING
536 s1_u_obc(i,k,ib,ir,lout)=s1_u_obc(i,k,ib,ir,lout)* &
537 & umask(i,jstr-1)
538# endif
539 END DO
540 END DO
541 END IF
542 IF ((lobc(inorth,isuvel,ng)).and. &
543 & domain(ng)%Northern_Edge(tile)) THEN
544 ib=inorth
545 DO k=1,n(ng)
546 DO i=istru,iend
547 s1_u_obc(i,k,ib,ir,lout)=fac1*s1_u_obc(i,k,ib,ir,lin1)+ &
548 & fac2*s2_u_obc(i,k,ib,ir,lin2)
549# ifdef MASKING
550 s1_u_obc(i,k,ib,ir,lout)=s1_u_obc(i,k,ib,ir,lout)* &
551 & umask(i,jend+1)
552# endif
553 END DO
554 END DO
555 END IF
556 END DO
557 END IF
558# endif
559!
560! 3D V-momentum component.
561!
562 DO k=1,n(ng)
563 DO j=jstrp,jendt
564 DO i=istrt,iendt
565 s1_v(i,j,k,lout)=fac1*s1_v(i,j,k,lin1)+ &
566 & fac2*s2_v(i,j,k,lin2)
567# ifdef MASKING
568 s1_v(i,j,k,lout)=s1_v(i,j,k,lout)*vmask(i,j)
569# endif
570 END DO
571 END DO
572 END DO
573
574# ifdef ADJUST_BOUNDARY
575!
576! 3D V-momentum open boundaries.
577!
578 IF (any(lobc(:,isvvel,ng))) THEN
579 DO ir=1,nbrec(ng)
580 IF ((lobc(iwest,isvvel,ng)).and. &
581 & domain(ng)%Western_Edge(tile)) THEN
582 ib=iwest
583 DO k=1,n(ng)
584 DO j=jstrv,jend
585 s1_v_obc(j,k,ib,ir,lout)=fac1*s1_v_obc(j,k,ib,ir,lin1)+ &
586 & fac2*s2_v_obc(j,k,ib,ir,lin2)
587# ifdef MASKING
588 s1_v_obc(j,k,ib,ir,lout)=s1_v_obc(j,k,ib,ir,lout)* &
589 & vmask(istr-1,j)
590# endif
591 END DO
592 END DO
593 END IF
594 IF ((lobc(ieast,isvvel,ng)).and. &
595 & domain(ng)%Eastern_Edge(tile)) THEN
596 ib=ieast
597 DO k=1,n(ng)
598 DO j=jstrv,jend
599 s1_v_obc(j,k,ib,ir,lout)=fac1*s1_v_obc(j,k,ib,ir,lin1)+ &
600 & fac2*s2_v_obc(j,k,ib,ir,lin2)
601# ifdef MASKING
602 s1_v_obc(j,k,ib,ir,lout)=s1_v_obc(j,k,ib,ir,lout)* &
603 & vmask(iend+1,j)
604# endif
605 END DO
606 END DO
607 END IF
608 IF ((lobc(isouth,isvvel,ng)).and. &
609 & domain(ng)%Southern_Edge(tile)) THEN
610 ib=isouth
611 DO k=1,n(ng)
612 DO i=istr,iend
613 s1_v_obc(i,k,ib,ir,lout)=fac1*s1_v_obc(i,k,ib,ir,lin1)+ &
614 & fac2*s2_v_obc(i,k,ib,ir,lin2)
615# ifdef MASKING
616 s1_v_obc(i,k,ib,ir,lout)=s1_v_obc(i,k,ib,ir,lout)* &
617 & vmask(i,jstr)
618# endif
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 s1_v_obc(i,k,ib,ir,lout)=fac1*s1_v_obc(i,k,ib,ir,lin1)+ &
628 & fac2*s2_v_obc(i,k,ib,ir,lin2)
629# ifdef MASKING
630 s1_v_obc(i,k,ib,ir,lout)=s1_v_obc(i,k,ib,ir,lout)* &
631 & vmask(i,jend+1)
632# endif
633 END DO
634 END DO
635 END IF
636 END DO
637 END IF
638# endif
639!
640! Tracers.
641!
642 DO it=1,nt(ng)
643 DO k=1,n(ng)
644 DO j=jstrt,jendt
645 DO i=istrt,iendt
646 s1_t(i,j,k,lout,it)=fac1*s1_t(i,j,k,lin1,it)+ &
647 & fac2*s2_t(i,j,k,lin2,it)
648# ifdef MASKING
649 s1_t(i,j,k,lout,it)=s1_t(i,j,k,lout,it)*rmask(i,j)
650# endif
651 END DO
652 END DO
653 END DO
654 END DO
655
656# ifdef ADJUST_BOUNDARY
657!
658! Tracers open boundaries.
659!
660 DO it=1,nt(ng)
661 IF (any(lobc(:,istvar(it),ng))) THEN
662 DO ir=1,nbrec(ng)
663 IF ((lobc(iwest,istvar(it),ng)).and. &
664 & domain(ng)%Western_Edge(tile)) THEN
665 ib=iwest
666 DO k=1,n(ng)
667 DO j=jstr,jend
668 s1_t_obc(j,k,ib,ir,lout,it)= &
669 & fac1*s1_t_obc(j,k,ib,ir,lin1,it)+ &
670 & fac2*s2_t_obc(j,k,ib,ir,lin2,it)
671# ifdef MASKING
672 s1_t_obc(j,k,ib,ir,lout,it)= &
673 & s1_t_obc(j,k,ib,ir,lout,it)*rmask(istr-1,j)
674# endif
675 END DO
676 END DO
677 END IF
678 IF ((lobc(ieast,istvar(it),ng)).and. &
679 & domain(ng)%Eastern_Edge(tile)) THEN
680 ib=ieast
681 DO k=1,n(ng)
682 DO j=jstr,jend
683 s1_t_obc(j,k,ib,ir,lout,it)= &
684 & fac1*s1_t_obc(j,k,ib,ir,lin1,it)+ &
685 & fac2*s2_t_obc(j,k,ib,ir,lin2,it)
686# ifdef MASKING
687 s1_t_obc(j,k,ib,ir,lout,it)= &
688 & s1_t_obc(j,k,ib,ir,lout,it)*rmask(iend+1,j)
689# endif
690 END DO
691 END DO
692 END IF
693 IF ((lobc(isouth,istvar(it),ng)).and. &
694 & domain(ng)%Southern_Edge(tile)) THEN
695 ib=isouth
696 DO k=1,n(ng)
697 DO i=istr,iend
698 s1_t_obc(i,k,ib,ir,lout,it)= &
699 & fac1*s1_t_obc(i,k,ib,ir,lin1,it)+ &
700 & fac2*s2_t_obc(i,k,ib,ir,lin2,it)
701# ifdef MASKING
702 s1_t_obc(i,k,ib,ir,lout,it)= &
703 & s1_t_obc(i,k,ib,ir,lout,it)*rmask(i,jstr-1)
704# endif
705 END DO
706 END DO
707 END IF
708 IF ((lobc(inorth,istvar(it),ng)).and. &
709 & domain(ng)%Northern_Edge(tile)) THEN
710 ib=inorth
711 DO k=1,n(ng)
712 DO i=istr,iend
713 s1_t_obc(i,k,ib,ir,lout,it)= &
714 & fac1*s1_t_obc(i,k,ib,ir,lin1,it)+ &
715 & fac2*s2_t_obc(i,k,ib,ir,lin2,it)
716# ifdef MASKING
717 s1_t_obc(i,k,ib,ir,lout,it)= &
718 & s1_t_obc(i,k,ib,ir,lout,it)*rmask(i,jend+1)
719# endif
720 END DO
721 END DO
722 END IF
723 END DO
724 END IF
725 END DO
726# endif
727
728# ifdef ADJUST_STFLUX
729!
730! Surface tracers flux.
731!
732 DO it=1,nt(ng)
733 IF (lstflux(it,ng)) THEN
734 DO ir=1,nfrec(ng)
735 DO j=jstrt,jendt
736 DO i=istrt,iendt
737 s1_tflux(i,j,ir,lout,it)=fac1*s1_tflux(i,j,ir,lin1,it)+ &
738 & fac2*s2_tflux(i,j,ir,lin2,it)
739# ifdef MASKING
740 s1_tflux(i,j,ir,lout,it)=s1_tflux(i,j,ir,lout,it)* &
741 & rmask(i,j)
742# endif
743 END DO
744 END DO
745 END DO
746 END IF
747 END DO
748# endif
749
750#endif
751!
752 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, mod_scalars::lstflux, mod_param::n, mod_scalars::nbrec, mod_scalars::nfrec, and mod_param::nt.

Referenced by posterior_mod::analysis_error(), cgradient_mod::hessian_evecs(), ini_lanczos_mod::ini_lanczos_tile(), cgradient_mod::lanczos(), posterior_mod::lanczos(), ini_adjust_mod::load_adtotl_tile(), ini_adjust_mod::load_tltoad_tile(), cgradient_mod::new_cost(), cgradient_mod::new_gradient(), posterior_mod::posterior_eofs(), posterior_var_mod::posterior_var_tile(), cgradient_mod::precond(), inner2state_mod::tl_inner2state_tile(), and cgradient_mod::tl_new_state().

Here is the caller graph for this function: