ROMS
Loading...
Searching...
No Matches
ad_frc_adjust.F
Go to the documentation of this file.
1#include "cppdefs.h"
3
4#if defined ADJUST_STFLUX || defined ADJUST_WSTRESS
5!
6!git $Id$
7!================================================== Hernan G. Arango ===
8! Copyright (c) 2002-2025 The ROMS Group !
9! Licensed under a MIT/X style license !
10! See License_ROMS.md !
11!=======================================================================
12! !
13! This routine is the adjoint of the time-interpolation of 4DVar !
14! surface forcing increments. The increments can be constant !
15! (Nfrec=1) or time interpolated between snapshots (Nfrec>1). !
16! !
17! On Input: !
18! !
19! ng Nested grid number. !
20! tile Domain partition. !
21! Linp 4DVar increment time index to process. !
22! !
23!=======================================================================
24!
25 implicit none
26!
27 PRIVATE
28 PUBLIC :: ad_frc_adjust
29!
30 CONTAINS
31!
32!***********************************************************************
33 SUBROUTINE ad_frc_adjust (ng, tile, Linp)
34!***********************************************************************
35!
36 USE mod_param
37 USE mod_forces
38!
39! Imported variable declarations.
40!
41 integer, intent(in) :: ng, tile, linp
42!
43! Local variable declarations.
44!
45 character (len=*), parameter :: myfile = &
46 & __FILE__
47!
48# include "tile.h"
49!
50# ifdef PROFILE
51 CALL wclock_on (ng, iadm, 7, __line__, myfile)
52# endif
53 CALL ad_frc_adjust_tile (ng, tile, &
54 & lbi, ubi, lbj, ubj, &
55 & imins, imaxs, jmins, jmaxs, &
56# ifdef ADJUST_WSTRESS
57 & forces(ng) % ad_ustr, &
58 & forces(ng) % ad_vstr, &
59 & forces(ng) % ad_sustr, &
60 & forces(ng) % ad_svstr, &
61# endif
62# if defined ADJUST_STFLUX && defined SOLVE3D
63 & forces(ng) % ad_tflux, &
64 & forces(ng) % ad_stflx, &
65# endif
66 & linp)
67# ifdef PROFILE
68 CALL wclock_off (ng, iadm, 7, __line__, myfile)
69# endif
70!
71 RETURN
72 END SUBROUTINE ad_frc_adjust
73!
74!***********************************************************************
75 SUBROUTINE ad_frc_adjust_tile (ng, tile, &
76 & LBi, UBi, LBj, UBj, &
77 & IminS, ImaxS, JminS, JmaxS, &
78# ifdef ADJUST_WSTRESS
79 & ad_ustr, ad_vstr, &
80 & ad_sustr, ad_svstr, &
81# endif
82# if defined ADJUST_STFLUX && defined SOLVE3D
83 & ad_tflux, ad_stflx, &
84# endif
85 & Linp)
86!***********************************************************************
87!
88 USE mod_param
89 USE mod_scalars
90!
91! Imported variable declarations.
92!
93 integer, intent(in) :: ng, tile
94 integer, intent(in) :: LBi, UBi, LBj, UBj
95 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
96 integer, intent(in) :: Linp
97!
98# ifdef ASSUMED_SHAPE
99# ifdef ADJUST_WSTRESS
100 real(r8), intent(inout) :: ad_sustr(LBi:,LBj:)
101 real(r8), intent(inout) :: ad_svstr(LBi:,LBj:)
102# endif
103# if defined ADJUST_STFLUX && defined SOLVE3D
104 real(r8), intent(inout) :: ad_stflx(LBi:,LBj:,:)
105# endif
106# ifdef ADJUST_WSTRESS
107 real(r8), intent(inout) :: ad_ustr(LBi:,LBj:,:,:)
108 real(r8), intent(inout) :: ad_vstr(LBi:,LBj:,:,:)
109# endif
110# if defined ADJUST_STFLUX && defined SOLVE3D
111 real(r8), intent(inout) :: ad_tflux(LBi:,LBj:,:,:,:)
112# endif
113# else
114# ifdef ADJUST_WSTRESS
115 real(r8), intent(inout) :: ad_sustr(LBi:UBi,LBj:UBj)
116 real(r8), intent(inout) :: ad_svstr(LBi:UBi,LBj:UBj)
117# endif
118# if defined ADJUST_STFLUX && defined SOLVE3D
119 real(r8), intent(inout) :: ad_stflx(LBi:UBi,LBj:UBj,NT(ng))
120# endif
121# ifdef ADJUST_WSTRESS
122 real(r8), intent(inout) :: ad_ustr(LBi:UBi,LBj:UBj,Nfrec(ng),2)
123 real(r8), intent(inout) :: ad_vstr(LBi:UBi,LBj:UBj,Nfrec(ng),2)
124# endif
125# if defined ADJUST_STFLUX && defined SOLVE3D
126 real(r8), intent(inout) :: ad_tflux(LBi:UBi,LBj:UBj, &
127 & Nfrec(ng),2,NT(ng))
128# endif
129# endif
130!
131! Local variable declarations.
132!
133 integer :: i, it1, it2, j
134# ifdef SOLVE3D
135 integer :: itrc
136# endif
137 real(r8) :: fac, fac1, fac2
138
139# include "set_bounds.h"
140!
141!-----------------------------------------------------------------------
142! Adjust nonlinear surface forcing fields with 4DVAR increments.
143!-----------------------------------------------------------------------
144!
145! Set time records and interpolation factor, if any.
146!
147 IF (nfrec(ng).eq.1) THEN
148 it1=1
149 it2=1
150 fac1=1.0_r8
151 fac2=0.0_r8
152 ELSE
153# ifdef GENERIC_DSTART
154 it1=max(0,(iic(ng)-ntstart(ng))/nsff(ng))+1
155# else
156 it1=max(0,(iic(ng)-1)/nsff(ng))+1
157# endif
158 it2=min(it1+1,nfrec(ng))
159 fac1=sf_time(it2,ng)-(time(ng)+dt(ng))
160 fac2=(time(ng)+dt(ng))-sf_time(it1,ng)
161 fac=1.0_r8/(fac1+fac2)
162 fac1=fac*fac1
163 fac2=fac*fac2
164 END IF
165
166# ifdef ADJUST_WSTRESS
167!
168! Adjoint of surface wind stress adjustment.
169!
170 DO j=jstrr,jendr
171 DO i=istr,iendr
172!^ tl_sustr(i,j)=fac1*tl_ustr(i,j,it1,Linp)+ &
173!^ & fac2*tl_ustr(i,j,it2,Linp)
174!^
175 ad_ustr(i,j,it1,linp)=ad_ustr(i,j,it1,linp)+ &
176 & fac1*ad_sustr(i,j)
177 ad_ustr(i,j,it2,linp)=ad_ustr(i,j,it2,linp)+ &
178 & fac2*ad_sustr(i,j)
179 ad_sustr(i,j)=0.0_r8
180 END DO
181 END DO
182 DO j=jstr,jendr
183 DO i=istrr,iendr
184!^ tl_svstr(i,j)=fac1*tl_vstr(i,j,it1,Linp)+ &
185!^ & fac2*tl_vstr(i,j,it2,Linp)
186!^
187 ad_vstr(i,j,it1,linp)=ad_vstr(i,j,it1,linp)+ &
188 & fac1*ad_svstr(i,j)
189 ad_vstr(i,j,it2,linp)=ad_vstr(i,j,it2,linp)+ &
190 & fac2*ad_svstr(i,j)
191 ad_svstr(i,j)=0.0_r8
192 END DO
193 END DO
194# endif
195# if defined ADJUST_STFLUX && defined SOLVE3D
196!
197! Adjoint of surface tracer fluxes adjustmenst.
198!
199 IF (lstflux(itemp,ng)) THEN
200# ifdef QCORRECTION
201 DO j=jstrr,jendr
202 DO i=istrr,iendr
203!^ tl_stflx(i,j,itemp)=tl_stflx(i,j,itemp)+ &
204!^ & fac1*tl_tflux(i,j,it1,Linp,itemp)+ &
205!^ & fac2*tl_tflux(i,j,it2,Linp,itemp)
206!^
207 ad_tflux(i,j,it1,linp,itemp)=ad_tflux(i,j,it1,linp,itemp)+ &
208 & fac1*ad_stflx(i,j,itemp)
209 ad_tflux(i,j,it2,linp,itemp)=ad_tflux(i,j,it2,linp,itemp)+ &
210 & fac2*ad_stflx(i,j,itemp)
211 END DO
212 END DO
213# else
214 DO j=jstrr,jendr
215 DO i=istrr,iendr
216!^ tl_stflx(i,j,itemp)=fac1*tl_tflux(i,j,it1,Linp,itemp)+ &
217!^ & fac2*tl_tflux(i,j,it2,Linp,itemp)
218!^
219 ad_tflux(i,j,it1,linp,itemp)=ad_tflux(i,j,it1,linp,itemp)+ &
220 & fac1*ad_stflx(i,j,itemp)
221 ad_tflux(i,j,it2,linp,itemp)=ad_tflux(i,j,it2,linp,itemp)+ &
222 & fac2*ad_stflx(i,j,itemp)
223 ad_stflx(i,j,itemp)=0.0_r8
224 END DO
225 END DO
226# endif
227 END IF
228
229# ifdef SALINITY
230 IF (lstflux(isalt,ng)) THEN
231 DO j=jstrr,jendr
232 DO i=istrr,iendr
233!^ tl_stflx(i,j,isalt)=tl_stflx(i,j,isalt)+ &
234!^ & fac1*tl_tflux(i,j,it1,Linp,isalt)+ &
235!^ & fac2*tl_tflux(i,j,it2,Linp,isalt)
236!^
237 ad_tflux(i,j,it1,linp,isalt)=ad_tflux(i,j,it1,linp,isalt)+ &
238 & fac1*ad_stflx(i,j,isalt)
239 ad_tflux(i,j,it2,linp,isalt)=ad_tflux(i,j,it2,linp,isalt)+ &
240 & fac2*ad_stflx(i,j,isalt)
241 END DO
242 END DO
243 END IF
244# endif
245
246 DO itrc=nat+1,nt(ng)
247 IF (lstflux(itrc,ng)) THEN
248 DO j=jstrr,jendr
249 DO i=istrr,iendr
250!^ tl_stflx(i,j,itrc)=fac1*tl_tflux(i,j,it1,Linp,itrc)+ &
251!^ & fac2*tl_tflux(i,j,it2,Linp,itrc)
252!^
253 ad_tflux(i,j,it1,linp,itrc)=ad_tflux(i,j,it1,linp,itrc)+ &
254 & fac1*ad_stflx(i,j,itrc)
255 ad_tflux(i,j,it2,linp,itrc)=ad_tflux(i,j,it2,linp,itrc)+ &
256 & fac2*ad_stflx(i,j,itrc)
257 ad_stflx(i,j,itrc)=0.0_r8
258 END DO
259 END DO
260 END IF
261 END DO
262# endif
263!
264 RETURN
265 END SUBROUTINE ad_frc_adjust_tile
266#endif
267 END MODULE ad_frc_adjust_mod
subroutine ad_frc_adjust_tile(ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, ad_ustr, ad_vstr, ad_sustr, ad_svstr, ad_tflux, ad_stflx, linp)
subroutine, public ad_frc_adjust(ng, tile, linp)
type(t_forces), dimension(:), allocatable forces
Definition mod_forces.F:554
integer nat
Definition mod_param.F:499
integer, parameter iadm
Definition mod_param.F:665
integer, dimension(:), allocatable iic
real(dp), dimension(:), allocatable dt
logical, dimension(:,:), allocatable lstflux
integer isalt
integer itemp
real(dp), dimension(:), allocatable time
integer, dimension(:), allocatable nsff
integer, dimension(:), allocatable ntstart
real(dp), dimension(:,:), allocatable sf_time
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