ROMS
Loading...
Searching...
No Matches
forcing.F
Go to the documentation of this file.
1#include "cppdefs.h"
3
4#if defined NLM_OUTER || defined RBL4DVAR || \
5 defined sensitivity_4dvar || defined sp4dvar || \
6 defined tl_rbl4dvar
7
8!
9!git $Id$
10!================================================== Hernan G. Arango ===
11! Copyright (c) 2002-2025 The ROMS Group !
12! Licensed under a MIT/X style license !
13! See License_ROMS.md !
14!=======================================================================
15! !
16! This routine is used to force the nonlinear state equations in !
17! weak constraint variational data assimilation. !
18! !
19!=======================================================================
20!
21 implicit none
22!
23 PRIVATE
24 PUBLIC :: forcing
25!
26 CONTAINS
27!
28!***********************************************************************
29 SUBROUTINE forcing (ng, tile, Kfrc, Nfrc)
30!***********************************************************************
31!
32 USE mod_param
33 USE mod_ocean
34# ifdef SOLVE3D
35 USE mod_coupling
36# endif
37!
38! Imported variable declarations.
39!
40 integer, intent(in) :: ng, tile, kfrc, nfrc
41!
42! Local variable declarations.
43!
44# include "tile.h"
45!
46 CALL forcing_tile (ng, tile, &
47 & lbi, ubi, lbj, ubj, &
48 & imins, imaxs, jmins, jmaxs, &
49 & kfrc, nfrc, &
50# ifdef SOLVE3D
51 & ocean(ng) % f_t, &
52 & ocean(ng) % f_u, &
53 & ocean(ng) % f_v, &
54# endif
55 & ocean(ng) % f_ubar, &
56 & ocean(ng) % f_vbar, &
57 & ocean(ng) % f_zeta, &
58# ifdef SOLVE3D
59 & ocean(ng) % t, &
60 & ocean(ng) % u, &
61 & ocean(ng) % v, &
62# endif
63 & ocean(ng) % ubar, &
64 & ocean(ng) % vbar, &
65# ifdef SOLVE3D
66 & coupling(ng) % Zt_avg1, &
67# endif
68 & ocean(ng) % zeta)
69!
70 RETURN
71 END SUBROUTINE forcing
72!
73!***********************************************************************
74 SUBROUTINE forcing_tile (ng, tile, &
75 & LBi, UBi, LBj, UBj, &
76 & IminS, ImaxS, JminS, JmaxS, &
77 & Kfrc, Nfrc, &
78# ifdef SOLVE3D
79 & f_t, f_u, f_v, &
80# endif
81 & f_ubar, f_vbar, &
82 & f_zeta, &
83# ifdef SOLVE3D
84 & t, u, v, &
85# endif
86 & ubar, vbar, &
87# ifdef SOLVE3D
88 & Zt_avg1, &
89# endif
90 & zeta)
91!***********************************************************************
92!
93 USE mod_param
94 USE mod_parallel
95 USE mod_iounits
96 USE mod_scalars
97!
98! Imported variable declarations.
99!
100 integer, intent(in) :: ng, tile
101 integer, intent(in) :: LBi, UBi, LBj, UBj
102 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
103 integer, intent(in) :: Kfrc
104 integer, intent(in) :: Nfrc
105!
106# ifdef ASSUMED_SHAPE
107# ifdef SOLVE3D
108 real(r8), intent(in) :: f_t(LBi:,LBj:,:,:)
109 real(r8), intent(in) :: f_u(LBi:,LBj:,:)
110 real(r8), intent(in) :: f_v(LBi:,LBj:,:)
111# endif
112 real(r8), intent(in) :: f_ubar(LBi:,LBj:)
113 real(r8), intent(in) :: f_vbar(LBi:,LBj:)
114 real(r8), intent(in) :: f_zeta(LBi:,LBj:)
115# ifdef SOLVE3D
116 real(r8), intent(inout) :: t(LBi:,LBj:,:,:,:)
117 real(r8), intent(inout) :: u(LBi:,LBj:,:,:)
118 real(r8), intent(inout) :: v(LBi:,LBj:,:,:)
119# endif
120 real(r8), intent(inout) :: ubar(LBi:,LBj:,:)
121 real(r8), intent(inout) :: vbar(LBi:,LBj:,:)
122# ifdef SOLVE3D
123 real(r8), intent(inout) :: Zt_avg1(LBi:,LBj:)
124# endif
125 real(r8), intent(inout) :: zeta(LBi:,LBj:,:)
126# else
127# ifdef SOLVE3D
128 real(r8), intent(in) :: f_t(LBi:UBi,LBj:UBj,N(ng),NT(ng))
129 real(r8), intent(in) :: f_u(LBi:UBi,LBj:UBj,N(ng))
130 real(r8), intent(in) :: f_v(LBi:UBi,LBj:UBj,N(ng))
131# endif
132 real(r8), intent(in) :: f_ubar(LBi:UBi,LBj:UBj)
133 real(r8), intent(in) :: f_vbar(LBi:UBi,LBj:UBj)
134 real(r8), intent(in) :: f_zeta(LBi:UBi,LBj:UBj)
135# ifdef SOLVE3D
136 real(r8), intent(inout) :: t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
137 real(r8), intent(inout) :: u(LBi:UBi,LBj:UBj,N(ng),2)
138 real(r8), intent(inout) :: v(LBi:UBi,LBj:UBj,N(ng),2)
139# endif
140 real(r8), intent(inout) :: ubar(LBi:UBi,LBj:UBj,:)
141 real(r8), intent(inout) :: vbar(LBi:UBi,LBj:UBj,:)
142# ifdef SOLVE3D
143 real(r8), intent(inout) :: Zt_avg1(LBi:UBi,LBj:UBj)
144# endif
145 real(r8), intent(inout) :: zeta(LBi:UBi,LBj:UBj,:)
146# endif
147!
148! Local variable declarations.
149!
150 integer :: i, j
151# ifdef SOLVE3D
152 integer :: itrc, k
153# endif
154
155# include "set_bounds.h"
156!
157!-----------------------------------------------------------------------
158! Add adjoint impulse forcing to nonlinear linear state.
159!-----------------------------------------------------------------------
160!
161 IF (domain(ng)%SouthWest_Corner(tile)) THEN
162 IF (frequentimpulse(ng)) THEN
163 IF (master) WRITE (stdout,10) time_code(ng)
164 END IF
165# ifdef RPCG
166 IF (iauswitch(ng)) THEN
167 IF (master) WRITE (stdout,20) time_code(ng)
168 END IF
169# endif
170 END IF
171!
172! Free-surface.
173!
174# ifdef SOLVE3D
175 IF (iic(ng).eq.ntstart(ng)) THEN
176 DO j=jstrr,jendr
177 DO i=istrr,iendr
178 zeta(i,j,kfrc)=zeta(i,j,kfrc)+f_zeta(i,j)
179 END DO
180 END DO
181 ELSE
182 DO j=jstrr,jendr
183 DO i=istrr,iendr
184 zt_avg1(i,j)=zt_avg1(i,j)+f_zeta(i,j)
185 END DO
186 END DO
187 END IF
188# else
189 DO j=jstrr,jendr
190 DO i=istrr,iendr
191 zeta(i,j,kfrc)=zeta(i,j,kfrc)+f_zeta(i,j)
192 END DO
193 END DO
194# endif
195
196# ifndef SOLVE3D
197!
198! 2D momentum.
199!
200 DO j=jstrr,jendr
201 DO i=istr,iendr
202 ubar(i,j,kfrc)=ubar(i,j,kfrc)+f_ubar(i,j)
203 END DO
204 END DO
205!
206 DO j=jstr,jendr
207 DO i=istrr,iendr
208 vbar(i,j,kfrc)=vbar(i,j,kfrc)+f_vbar(i,j)
209 END DO
210 END DO
211
212# else
213!
214! 3D momentum.
215!
216 DO k=1,n(ng)
217 DO j=jstrr,jendr
218 DO i=istr,iendr
219 u(i,j,k,nfrc)=u(i,j,k,nfrc)+f_u(i,j,k)
220 END DO
221 END DO
222 DO j=jstr,jendr
223 DO i=istrr,iendr
224 v(i,j,k,nfrc)=v(i,j,k,nfrc)+f_v(i,j,k)
225 END DO
226 END DO
227 END DO
228!
229! Tracers.
230!
231 DO itrc=1,nt(ng)
232 DO k=1,n(ng)
233 DO j=jstrr,jendr
234 DO i=istrr,iendr
235 t(i,j,k,nfrc,itrc)=t(i,j,k,nfrc,itrc)+ &
236 & f_t(i,j,k,itrc)
237 END DO
238 END DO
239 END DO
240 END DO
241# endif
242!
243 10 FORMAT (2x,'NL_FORCING - added convolved adjoint impulse,', &
244 & t71,'t = ', a)
245# ifdef RPCG
246 20 FORMAT (2x,'NL_FORCING - incremental analysis update, ', &
247 & t71,'t = ', a)
248# endif
249!
250 RETURN
251 END SUBROUTINE forcing_tile
252#endif
253 END MODULE forcing_mod
subroutine, public forcing(ng, tile, kfrc, nfrc)
Definition forcing.F:30
subroutine forcing_tile(ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, kfrc, nfrc, f_t, f_u, f_v, f_ubar, f_vbar, f_zeta, t, u, v, ubar, vbar, zt_avg1, zeta)
Definition forcing.F:91
type(t_coupling), dimension(:), allocatable coupling
integer stdout
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, dimension(:), allocatable iic
logical, dimension(:), allocatable frequentimpulse
character(len=22), dimension(:), allocatable time_code
logical, dimension(:), allocatable iauswitch
integer, dimension(:), allocatable ntstart