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

Functions/Subroutines

subroutine, public state_scale (ng, tile, lbi, ubi, lbj, ubj, lbij, ubij, linp, 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)
 

Function/Subroutine Documentation

◆ state_scale()

subroutine, public state_scale_mod::state_scale ( 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) linp,
integer, intent(in) lout,
real(r8), intent(in) fac,
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) s_t_obc,
real(r8), dimension(lbij:,:,:,:,:), intent(inout) s_u_obc,
real(r8), dimension(lbij:,:,:,:,:), intent(inout) s_v_obc,
real(r8), dimension(lbij:,:,:,:), intent(inout) s_ubar_obc,
real(r8), dimension(lbij:,:,:,:), intent(inout) s_vbar_obc,
real(r8), dimension(lbij:,:,:,:), intent(inout) s_zeta_obc,
real(r8), dimension(lbi:,lbj:,:,:), intent(inout) s_sustr,
real(r8), dimension(lbi:,lbj:,:,:), intent(inout) s_svstr,
real(r8), dimension(lbi:,lbj:,:,:,:), intent(inout) s_tflux,
real(r8), dimension(lbi:,lbj:,:,:,:), intent(inout) s_t,
real(r8), dimension(lbi:,lbj:,:,:), intent(inout) s_u,
real(r8), dimension(lbi:,lbj:,:,:), intent(inout) s_v,
real(r8), dimension(lbi:,lbj:,:), intent(inout) s_zeta )

Definition at line 24 of file state_scale.F.

50!***********************************************************************
51!
52 USE mod_param
53#if defined ADJUST_BOUNDARY || defined ADJUST_STFLUX || \
54 defined adjust_wstress
55 USE mod_ncparam
56 USE mod_scalars
57#endif
58!
59! Imported variable declarations.
60!
61 integer, intent(in) :: ng, tile
62 integer, intent(in) :: LBi, UBi, LBj, UBj, LBij, UBij
63 integer, intent(in) :: Linp, Lout
64!
65 real(r8), intent(in) :: fac
66!
67#ifdef ASSUMED_SHAPE
68# ifdef MASKING
69 real(r8), intent(in) :: rmask(LBi:,LBj:)
70 real(r8), intent(in) :: umask(LBi:,LBj:)
71 real(r8), intent(in) :: vmask(LBi:,LBj:)
72# endif
73# ifdef ADJUST_BOUNDARY
74# ifdef SOLVE3D
75 real(r8), intent(inout) :: s_t_obc(LBij:,:,:,:,:,:)
76 real(r8), intent(inout) :: s_u_obc(LBij:,:,:,:,:)
77 real(r8), intent(inout) :: s_v_obc(LBij:,:,:,:,:)
78# endif
79 real(r8), intent(inout) :: s_ubar_obc(LBij:,:,:,:)
80 real(r8), intent(inout) :: s_vbar_obc(LBij:,:,:,:)
81 real(r8), intent(inout) :: s_zeta_obc(LBij:,:,:,:)
82# endif
83# ifdef ADJUST_WSTRESS
84 real(r8), intent(inout) :: s_sustr(LBi:,LBj:,:,:)
85 real(r8), intent(inout) :: s_svstr(LBi:,LBj:,:,:)
86# endif
87# ifdef SOLVE3D
88# ifdef ADJUST_STFLUX
89 real(r8), intent(inout) :: s_tflux(LBi:,LBj:,:,:,:)
90# endif
91 real(r8), intent(inout) :: s_t(LBi:,LBj:,:,:,:)
92 real(r8), intent(inout) :: s_u(LBi:,LBj:,:,:)
93 real(r8), intent(inout) :: s_v(LBi:,LBj:,:,:)
94# else
95 real(r8), intent(inout) :: s_ubar(LBi:,LBj:,:)
96 real(r8), intent(inout) :: s_vbar(LBi:,LBj:,:)
97# endif
98 real(r8), intent(inout) :: s_zeta(LBi:,LBj:,:)
99
100#else
101
102# ifdef MASKING
103 real(r8), intent(in) :: rmask(LBi:UBi,LBj:UBj)
104 real(r8), intent(in) :: umask(LBi:UBi,LBj:UBj)
105 real(r8), intent(in) :: vmask(LBi:UBi,LBj:UBj)
106# endif
107# ifdef ADJUST_BOUNDARY
108# ifdef SOLVE3D
109 real(r8), intent(inout) :: s_t_obc(LBij:UBij,N(ng),4, &
110 & Nbrec(ng),2,NT(ng))
111 real(r8), intent(inout) :: s_u_obc(LBij:UBij,N(ng),4,Nbrec(ng),2)
112 real(r8), intent(inout) :: s_v_obc(LBij:UBij,N(ng),4,Nbrec(ng),2)
113# endif
114 real(r8), intent(inout) :: s_ubar_obc(LBij:UBij,4,Nbrec(ng),2)
115 real(r8), intent(inout) :: s_vbar_obc(LBij:UBij,4,Nbrec(ng),2)
116 real(r8), intent(inout) :: s_zeta_obc(LBij:UBij,4,Nbrec(ng),2)
117# endif
118# ifdef ADJUST_WSTRESS
119 real(r8), intent(inout) :: s_sustr(LBi:UBi,LBj:UBj,Nfrec(ng),2)
120 real(r8), intent(inout) :: s_svstr(LBi:UBi,LBj:UBj,Nfrec(ng),2)
121# endif
122# ifdef SOLVE3D
123# ifdef ADJUST_STFLUX
124 real(r8), intent(inout) :: s_tflux(LBi:UBi,LBj:UBj, &
125 & Nfrec(ng),2,NT(ng))
126# endif
127 real(r8), intent(inout) :: s_t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
128 real(r8), intent(inout) :: s_u(LBi:UBi,LBj:UBj,N(ng),2)
129 real(r8), intent(inout) :: s_v(LBi:UBi,LBj:UBj,N(ng),2)
130# else
131 real(r8), intent(inout) :: s_ubar(LBi:UBi,LBj:UBj,:)
132 real(r8), intent(inout) :: s_vbar(LBi:UBi,LBj:UBj,:)
133# endif
134 real(r8), intent(inout) :: s_zeta(LBi:UBi,LBj:UBj,:)
135#endif
136!
137! Local variable declarations.
138!
139 integer :: i, j, k
140 integer :: ib, ir, it
141
142#include "set_bounds.h"
143!
144!-----------------------------------------------------------------------
145! Scale model state variable by a constant.
146!-----------------------------------------------------------------------
147!
148! Free-surface.
149!
150 DO j=jstrt,jendt
151 DO i=istrt,iendt
152 s_zeta(i,j,lout)=fac*s_zeta(i,j,linp)
153#ifdef MASKING
154 s_zeta(i,j,lout)=s_zeta(i,j,lout)*rmask(i,j)
155#endif
156 END DO
157 END DO
158
159#ifdef ADJUST_BOUNDARY
160!
161! Free-surface open boundaries.
162!
163 IF (any(lobc(:,isfsur,ng))) THEN
164 DO ir=1,nbrec(ng)
165 IF ((lobc(iwest,isfsur,ng)).and. &
166 & domain(ng)%Western_Edge(tile)) THEN
167 ib=iwest
168 DO j=jstr,jend
169 s_zeta_obc(j,ib,ir,lout)=fac*s_zeta_obc(j,ib,ir,linp)
170# ifdef MASKING
171 s_zeta_obc(j,ib,ir,lout)=s_zeta_obc(j,ib,ir,lout)* &
172 & rmask(istr-1,j)
173# endif
174 END DO
175 END IF
176 IF ((lobc(ieast,isfsur,ng)).and. &
177 & domain(ng)%Eastern_Edge(tile)) THEN
178 ib=ieast
179 DO j=jstr,jend
180 s_zeta_obc(j,ib,ir,lout)=fac*s_zeta_obc(j,ib,ir,linp)
181# ifdef MASKING
182 s_zeta_obc(j,ib,ir,lout)=s_zeta_obc(j,ib,ir,lout)* &
183 & rmask(iend+1,j)
184# endif
185 END DO
186 END IF
187 IF ((lobc(isouth,isfsur,ng)).and. &
188 & domain(ng)%Southern_Edge(tile)) THEN
189 ib=isouth
190 DO i=istr,iend
191 s_zeta_obc(i,ib,ir,lout)=fac*s_zeta_obc(i,ib,ir,linp)
192# ifdef MASKING
193 s_zeta_obc(i,ib,ir,lout)=s_zeta_obc(i,ib,ir,lout)* &
194 & rmask(i,jstr-1)
195# endif
196 END DO
197 END IF
198 IF ((lobc(inorth,isfsur,ng)).and. &
199 & domain(ng)%Northern_Edge(tile)) THEN
200 ib=inorth
201 DO i=istr,iend
202 s_zeta_obc(i,ib,ir,lout)=fac*s_zeta_obc(i,ib,ir,linp)
203# ifdef MASKING
204 s_zeta_obc(i,ib,ir,lout)=s_zeta_obc(i,ib,ir,lout)* &
205 & rmask(i,jend+1)
206# endif
207 END DO
208 END IF
209 END DO
210 END IF
211#endif
212
213#ifndef SOLVE3D
214!
215! 2D U-momentum component.
216!
217 DO j=jstrt,jendt
218 DO i=istrp,iendt
219 s_ubar(i,j,lout)=fac*s_ubar(i,j,linp)
220# ifdef MASKING
221 s_ubar(i,j,lout)=s_ubar(i,j,lout)*umask(i,j)
222# endif
223 END DO
224 END DO
225#endif
226
227#ifdef ADJUST_BOUNDARY
228!
229! 2D U-momentum open boundaries.
230!
231 IF (any(lobc(:,isubar,ng))) THEN
232 DO ir=1,nbrec(ng)
233 IF ((lobc(iwest,isubar,ng)).and. &
234 & domain(ng)%Western_Edge(tile)) THEN
235 ib=iwest
236 DO j=jstr,jend
237 s_ubar_obc(j,ib,ir,lout)=fac*s_ubar_obc(j,ib,ir,linp)
238# ifdef MASKING
239 s_ubar_obc(j,ib,ir,lout)=s_ubar_obc(j,ib,ir,lout)* &
240 & umask(istr,j)
241# endif
242 END DO
243 END IF
244 IF ((lobc(ieast,isubar,ng)).and. &
245 & domain(ng)%Eastern_Edge(tile)) THEN
246 ib=ieast
247 DO j=jstr,jend
248 s_ubar_obc(j,ib,ir,lout)=fac*s_ubar_obc(j,ib,ir,linp)
249# ifdef MASKING
250 s_ubar_obc(j,ib,ir,lout)=s_ubar_obc(j,ib,ir,lout)* &
251 & umask(iend+1,j)
252# endif
253 END DO
254 END IF
255 IF ((lobc(isouth,isubar,ng)).and. &
256 & domain(ng)%Southern_Edge(tile)) THEN
257 ib=isouth
258 DO i=istru,iend
259 s_ubar_obc(i,ib,ir,lout)=fac*s_ubar_obc(i,ib,ir,linp)
260# ifdef MASKING
261 s_ubar_obc(i,ib,ir,lout)=s_ubar_obc(i,ib,ir,lout)* &
262 & umask(i,jstr-1)
263# endif
264 END DO
265 END IF
266 IF ((lobc(inorth,isubar,ng)).and. &
267 & domain(ng)%Northern_Edge(tile)) THEN
268 ib=inorth
269 DO i=istru,iend
270 s_ubar_obc(i,ib,ir,lout)=fac*s_ubar_obc(i,ib,ir,linp)
271# ifdef MASKING
272 s_ubar_obc(i,ib,ir,lout)=s_ubar_obc(i,ib,ir,lout)* &
273 & umask(i,jend+1)
274# endif
275 END DO
276 END IF
277 END DO
278 END IF
279#endif
280
281#ifndef SOLVE3D
282!
283! 2D V-momentum component.
284!
285 DO j=jstrp,jendt
286 DO i=istrt,iendt
287 s_vbar(i,j,lout)=fac*s_vbar(i,j,linp)
288# ifdef MASKING
289 s_vbar(i,j,lout)=s_vbar(i,j,lout)*vmask(i,j)
290# endif
291 END DO
292 END DO
293#endif
294
295#ifdef ADJUST_BOUNDARY
296!
297! 2D V-momentum open boundaries.
298!
299 IF (any(lobc(:,isvbar,ng))) THEN
300 DO ir=1,nbrec(ng)
301 IF ((lobc(iwest,isvbar,ng)).and. &
302 & domain(ng)%Western_Edge(tile)) THEN
303 ib=iwest
304 DO j=jstrv,jend
305 s_vbar_obc(j,ib,ir,lout)=fac*s_vbar_obc(j,ib,ir,linp)
306# ifdef MASKING
307 s_vbar_obc(j,ib,ir,lout)=s_vbar_obc(j,ib,ir,lout)* &
308 & vmask(istr-1,j)
309# endif
310 END DO
311 END IF
312 IF ((lobc(ieast,isvbar,ng)).and. &
313 & domain(ng)%Eastern_Edge(tile)) THEN
314 ib=ieast
315 DO j=jstrv,jend
316 s_vbar_obc(j,ib,ir,lout)=fac*s_vbar_obc(j,ib,ir,linp)
317# ifdef MASKING
318 s_vbar_obc(j,ib,ir,lout)=s_vbar_obc(j,ib,ir,lout)* &
319 & vmask(iend+1,j)
320# endif
321 END DO
322 END IF
323 IF ((lobc(isouth,isvbar,ng)).and. &
324 & domain(ng)%Southern_Edge(tile)) THEN
325 ib=isouth
326 DO i=istr,iend
327 s_vbar_obc(i,ib,ir,lout)=fac*s_vbar_obc(i,ib,ir,linp)
328# ifdef MASKING
329 s_vbar_obc(i,ib,ir,lout)=s_vbar_obc(i,ib,ir,lout)* &
330 & vmask(i,jstr)
331# endif
332 END DO
333 END IF
334 IF ((lobc(inorth,isvbar,ng)).and. &
335 & domain(ng)%Northern_Edge(tile)) THEN
336 ib=inorth
337 DO i=istr,iend
338 s_vbar_obc(i,ib,ir,lout)=fac*s_vbar_obc(i,ib,ir,linp)
339# ifdef MASKING
340 s_vbar_obc(i,ib,ir,lout)=s_vbar_obc(i,ib,ir,lout)* &
341 & vmask(i,jend+1)
342# endif
343 END DO
344 END IF
345 END DO
346 END IF
347#endif
348
349#ifdef ADJUST_WSTRESS
350!
351! Surface momentum stress.
352!
353 DO ir=1,nfrec(ng)
354 DO j=jstrt,jendt
355 DO i=istrp,iendt
356 s_sustr(i,j,ir,lout)=fac*s_sustr(i,j,ir,linp)
357# ifdef MASKING
358 s_sustr(i,j,ir,lout)=s_sustr(i,j,ir,lout)*umask(i,j)
359# endif
360 END DO
361 END DO
362 DO j=jstrp,jendt
363 DO i=istrt,iendt
364 s_svstr(i,j,ir,lout)=fac*s_svstr(i,j,ir,linp)
365# ifdef MASKING
366 s_svstr(i,j,ir,lout)=s_svstr(i,j,ir,lout)*vmask(i,j)
367# endif
368 END DO
369 END DO
370 END DO
371#endif
372
373#ifdef SOLVE3D
374!
375! 3D U-momentum component.
376!
377 DO k=1,n(ng)
378 DO j=jstrt,jendt
379 DO i=istrp,iendt
380 s_u(i,j,k,lout)=fac*s_u(i,j,k,linp)
381# ifdef MASKING
382 s_u(i,j,k,lout)=s_u(i,j,k,lout)*umask(i,j)
383# endif
384 END DO
385 END DO
386 END DO
387
388# ifdef ADJUST_BOUNDARY
389!
390! 3D U-momentum open boundaries.
391!
392 IF (any(lobc(:,isuvel,ng))) THEN
393 DO ir=1,nbrec(ng)
394 IF ((lobc(iwest,isuvel,ng)).and. &
395 & domain(ng)%Western_Edge(tile)) THEN
396 ib=iwest
397 DO k=1,n(ng)
398 DO j=jstr,jend
399 s_u_obc(j,k,ib,ir,lout)=fac*s_u_obc(j,k,ib,ir,linp)
400# ifdef MASKING
401 s_u_obc(j,k,ib,ir,lout)=s_u_obc(j,k,ib,ir,lout)* &
402 & umask(istr,j)
403# endif
404 END DO
405 END DO
406 END IF
407 IF ((lobc(ieast,isuvel,ng)).and. &
408 & domain(ng)%Eastern_Edge(tile)) THEN
409 ib=ieast
410 DO k=1,n(ng)
411 DO j=jstr,jend
412 s_u_obc(j,k,ib,ir,lout)=fac*s_u_obc(j,k,ib,ir,linp)
413# ifdef MASKING
414 s_u_obc(j,k,ib,ir,lout)=s_u_obc(j,k,ib,ir,lout)* &
415 & umask(iend+1,j)
416# endif
417 END DO
418 END DO
419 END IF
420 IF ((lobc(isouth,isuvel,ng)).and. &
421 & domain(ng)%Southern_Edge(tile)) THEN
422 ib=isouth
423 DO k=1,n(ng)
424 DO i=istru,iend
425 s_u_obc(i,k,ib,ir,lout)=fac*s_u_obc(i,k,ib,ir,linp)
426# ifdef MASKING
427 s_u_obc(i,k,ib,ir,lout)=s_u_obc(i,k,ib,ir,lout)* &
428 & umask(i,jstr-1)
429# endif
430 END DO
431 END DO
432 END IF
433 IF ((lobc(inorth,isuvel,ng)).and. &
434 & domain(ng)%Northern_Edge(tile)) THEN
435 ib=inorth
436 DO k=1,n(ng)
437 DO i=istru,iend
438 s_u_obc(i,k,ib,ir,lout)=fac*s_u_obc(i,k,ib,ir,linp)
439# ifdef MASKING
440 s_u_obc(i,k,ib,ir,lout)=s_u_obc(i,k,ib,ir,lout)* &
441 & umask(i,jend+1)
442# endif
443 END DO
444 END DO
445 END IF
446 END DO
447 END IF
448# endif
449!
450! 3D V-momentum component.
451!
452 DO k=1,n(ng)
453 DO j=jstrp,jendt
454 DO i=istrt,iendt
455 s_v(i,j,k,lout)=fac*s_v(i,j,k,linp)
456# ifdef MASKING
457 s_v(i,j,k,lout)=s_v(i,j,k,lout)*vmask(i,j)
458# endif
459 END DO
460 END DO
461 END DO
462
463# ifdef ADJUST_BOUNDARY
464!
465! 3D V-momentum open boundaries.
466!
467 IF (any(lobc(:,isvvel,ng))) THEN
468 DO ir=1,nbrec(ng)
469 IF ((lobc(iwest,isvvel,ng)).and. &
470 & domain(ng)%Western_Edge(tile)) THEN
471 ib=iwest
472 DO k=1,n(ng)
473 DO j=jstrv,jend
474 s_v_obc(j,k,ib,ir,lout)=fac*s_v_obc(j,k,ib,ir,linp)
475# ifdef MASKING
476 s_v_obc(j,k,ib,ir,lout)=s_v_obc(j,k,ib,ir,lout)* &
477 & vmask(istr-1,j)
478# endif
479 END DO
480 END DO
481 END IF
482 IF ((lobc(ieast,isvvel,ng)).and. &
483 & domain(ng)%Eastern_Edge(tile)) THEN
484 ib=ieast
485 DO k=1,n(ng)
486 DO j=jstrv,jend
487 s_v_obc(j,k,ib,ir,lout)=fac*s_v_obc(j,k,ib,ir,linp)
488# ifdef MASKING
489 s_v_obc(j,k,ib,ir,lout)=s_v_obc(j,k,ib,ir,lout)* &
490 & vmask(iend+1,j)
491# endif
492 END DO
493 END DO
494 END IF
495 IF ((lobc(isouth,isvvel,ng)).and. &
496 & domain(ng)%Southern_Edge(tile)) THEN
497 ib=isouth
498 DO k=1,n(ng)
499 DO i=istr,iend
500 s_v_obc(i,k,ib,ir,lout)=fac*s_v_obc(i,k,ib,ir,linp)
501# ifdef MASKING
502 s_v_obc(i,k,ib,ir,lout)=s_v_obc(i,k,ib,ir,lout)* &
503 & vmask(i,jstr)
504# endif
505 END DO
506 END DO
507 END IF
508 IF ((lobc(inorth,isvvel,ng)).and. &
509 & domain(ng)%Northern_Edge(tile)) THEN
510 ib=inorth
511 DO k=1,n(ng)
512 DO i=istr,iend
513 s_v_obc(i,k,ib,ir,lout)=fac*s_v_obc(i,k,ib,ir,linp)
514# ifdef MASKING
515 s_v_obc(i,k,ib,ir,lout)=s_v_obc(i,k,ib,ir,lout)* &
516 & vmask(i,jend+1)
517# endif
518 END DO
519 END DO
520 END IF
521 END DO
522 END IF
523# endif
524!
525! Tracers.
526!
527 DO it=1,nt(ng)
528 DO k=1,n(ng)
529 DO j=jstrt,jendt
530 DO i=istrt,iendt
531 s_t(i,j,k,lout,it)=fac*s_t(i,j,k,linp,it)
532# ifdef MASKING
533 s_t(i,j,k,lout,it)=s_t(i,j,k,lout,it)*rmask(i,j)
534# endif
535 END DO
536 END DO
537 END DO
538 END DO
539
540# ifdef ADJUST_BOUNDARY
541!
542! Tracers open boundaries.
543!
544 DO it=1,nt(ng)
545 IF (any(lobc(:,istvar(it),ng))) THEN
546 DO ir=1,nbrec(ng)
547 IF ((lobc(iwest,istvar(it),ng)).and. &
548 & domain(ng)%Western_Edge(tile)) THEN
549 ib=iwest
550 DO k=1,n(ng)
551 DO j=jstr,jend
552 s_t_obc(j,k,ib,ir,lout,it)=fac* &
553 & s_t_obc(j,k,ib,ir,linp,it)
554# ifdef MASKING
555 s_t_obc(j,k,ib,ir,lout,it)=s_t_obc(j,k,ib,ir,lout,it)*&
556 & rmask(istr-1,j)
557# endif
558 END DO
559 END DO
560 END IF
561 IF ((lobc(ieast,istvar(it),ng)).and. &
562 & domain(ng)%Eastern_Edge(tile)) THEN
563 ib=ieast
564 DO k=1,n(ng)
565 DO j=jstr,jend
566 s_t_obc(j,k,ib,ir,lout,it)=fac* &
567 & s_t_obc(j,k,ib,ir,linp,it)
568# ifdef MASKING
569 s_t_obc(j,k,ib,ir,lout,it)=s_t_obc(j,k,ib,ir,lout,it)*&
570 & rmask(iend+1,j)
571# endif
572 END DO
573 END DO
574 END IF
575 IF ((lobc(isouth,istvar(it),ng)).and. &
576 & domain(ng)%Southern_Edge(tile)) THEN
577 ib=isouth
578 DO k=1,n(ng)
579 DO i=istr,iend
580 s_t_obc(i,k,ib,ir,lout,it)=fac* &
581 & s_t_obc(i,k,ib,ir,linp,it)
582# ifdef MASKING
583 s_t_obc(i,k,ib,ir,lout,it)=s_t_obc(i,k,ib,ir,lout,it)*&
584 & rmask(i,jstr-1)
585# endif
586 END DO
587 END DO
588 END IF
589 IF ((lobc(inorth,istvar(it),ng)).and. &
590 & domain(ng)%Northern_Edge(tile)) THEN
591 ib=inorth
592 DO k=1,n(ng)
593 DO i=istr,iend
594 s_t_obc(i,k,ib,ir,lout,it)=fac* &
595 & s_t_obc(i,k,ib,ir,linp,it)
596# ifdef MASKING
597 s_t_obc(i,k,ib,ir,lout,it)=s_t_obc(i,k,ib,ir,lout,it)*&
598 & rmask(i,jend+1)
599# endif
600 END DO
601 END DO
602 END IF
603 END DO
604 END IF
605 END DO
606# endif
607
608# ifdef ADJUST_STFLUX
609!
610! Surface tracers flux.
611!
612 DO it=1,nt(ng)
613 IF (lstflux(it,ng)) THEN
614 DO ir=1,nfrec(ng)
615 DO j=jstrt,jendt
616 DO i=istrt,iendt
617 s_tflux(i,j,ir,lout,it)=fac*s_tflux(i,j,ir,linp,it)
618# ifdef MASKING
619 s_tflux(i,j,ir,lout,it)=s_tflux(i,j,ir,lout,it)* &
620 & rmask(i,j)
621# endif
622 END DO
623 END DO
624 END DO
625 END IF
626 END DO
627# endif
628#endif
629!
630 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 cgradient_mod::hessian_evecs(), cgradient_mod::lanczos(), posterior_mod::lanczos(), posterior_mod::posterior_eofs(), and inner2state_mod::tl_inner2state_tile().

Here is the caller graph for this function: