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

Functions/Subroutines

subroutine, public rp_frc_adjust (ng, tile, linp)
 
subroutine rp_frc_adjust_tile (ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, tl_ustr, tl_vstr, tl_sustr, tl_svstr, tl_tflux, tl_stflx, linp)
 

Function/Subroutine Documentation

◆ rp_frc_adjust()

subroutine, public rp_frc_adjust_mod::rp_frc_adjust ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) linp )

Definition at line 33 of file rp_frc_adjust.F.

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, irpm, 7, __line__, myfile)
52# endif
53 CALL rp_frc_adjust_tile (ng, tile, &
54 & lbi, ubi, lbj, ubj, &
55 & imins, imaxs, jmins, jmaxs, &
56# ifdef ADJUST_WSTRESS
57 & forces(ng) % tl_ustr, &
58 & forces(ng) % tl_vstr, &
59 & forces(ng) % tl_sustr, &
60 & forces(ng) % tl_svstr, &
61# endif
62# if defined ADJUST_STFLUX && defined SOLVE3D
63 & forces(ng) % tl_tflux, &
64 & forces(ng) % tl_stflx, &
65# endif
66 & linp)
67# ifdef PROFILE
68 CALL wclock_off (ng, irpm, 7, __line__, myfile)
69# endif
70!
71 RETURN
type(t_forces), dimension(:), allocatable forces
Definition mod_forces.F:554
integer, parameter irpm
Definition mod_param.F:664
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

References mod_forces::forces, mod_param::irpm, rp_frc_adjust_tile(), wclock_off(), and wclock_on().

Referenced by rp_main3d().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ rp_frc_adjust_tile()

subroutine rp_frc_adjust_mod::rp_frc_adjust_tile ( 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) imins,
integer, intent(in) imaxs,
integer, intent(in) jmins,
integer, intent(in) jmaxs,
real(r8), dimension(lbi:,lbj:,:,:), intent(in) tl_ustr,
real(r8), dimension(lbi:,lbj:,:,:), intent(in) tl_vstr,
real(r8), dimension(lbi:,lbj:), intent(inout) tl_sustr,
real(r8), dimension(lbi:,lbj:), intent(inout) tl_svstr,
real(r8), dimension(lbi:,lbj:,:,:,:), intent(inout) tl_tflux,
real(r8), dimension(lbi:,lbj:,:), intent(inout) tl_stflx,
integer, intent(in) linp )
private

Definition at line 75 of file rp_frc_adjust.F.

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(in) :: tl_ustr(LBi:,LBj:,:,:)
101 real(r8), intent(in) :: tl_vstr(LBi:,LBj:,:,:)
102# endif
103# if defined ADJUST_STFLUX && defined SOLVE3D
104 real(r8), intent(inout) :: tl_tflux(LBi:,LBj:,:,:,:)
105# endif
106# ifdef ADJUST_WSTRESS
107 real(r8), intent(inout) :: tl_sustr(LBi:,LBj:)
108 real(r8), intent(inout) :: tl_svstr(LBi:,LBj:)
109# endif
110# if defined ADJUST_STFLUX && defined SOLVE3D
111 real(r8), intent(inout) :: tl_stflx(LBi:,LBj:,:)
112# endif
113# else
114# ifdef ADJUST_WSTRESS
115 real(r8), intent(in) :: tl_ustr(LBi:UBi,LBj:UBj,Nfrec(ng),2)
116 real(r8), intent(in) :: tl_vstr(LBi:UBi,LBj:UBj,Nfrec(ng),2)
117# endif
118# if defined ADJUST_STFLUX && defined SOLVE3D
119 real(r8), intent(in) :: tl_tflux(LBi:UBi,LBj:UBj, &
120 & Nfrec(ng),2,NT(ng))
121# endif
122# ifdef ADJUST_WSTRESS
123 real(r8), intent(inout) :: tl_sustr(LBi:UBi,LBj:UBj)
124 real(r8), intent(inout) :: tl_svstr(LBi:UBi,LBj:UBj)
125# endif
126# if defined ADJUST_STFLUX && defined SOLVE3D
127 real(r8), intent(inout) :: tl_stflx(LBi:UBi,LBj:UBj,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 representer model surface forcing fields with 4DVar
143! increments.
144!-----------------------------------------------------------------------
145!
146! Set time records and interpolation factor, if any.
147!
148 IF (nfrec(ng).eq.1) THEN
149 it1=1
150 it2=1
151 fac1=1.0_r8
152 fac2=0.0_r8
153 ELSE
154# ifdef GENERIC_DSTART
155 it1=max(0,(iic(ng)-ntstart(ng))/nsff(ng))+1
156# else
157 it1=max(0,(iic(ng)-1)/nsff(ng))+1
158# endif
159 it2=min(it1+1,nfrec(ng))
160 fac1=sf_time(it2,ng)-(time(ng)+dt(ng))
161 fac2=(time(ng)+dt(ng))-sf_time(it1,ng)
162 fac=1.0_r8/(fac1+fac2)
163 fac1=fac*fac1
164 fac2=fac*fac2
165 END IF
166
167# ifdef ADJUST_WSTRESS
168!
169! Adjust surface wind stress. Interpolate between surface forcing
170! increments, if appropriate.
171!
172 DO j=jstrr,jendr
173 DO i=istr,iendr
174 tl_sustr(i,j)=tl_sustr(i,j)+ &
175 & fac1*tl_ustr(i,j,it1,linp)+ &
176 & fac2*tl_ustr(i,j,it2,linp)
177 END DO
178 END DO
179 DO j=jstr,jendr
180 DO i=istrr,iendr
181 tl_svstr(i,j)=tl_svstr(i,j)+ &
182 & fac1*tl_vstr(i,j,it1,linp)+ &
183 & fac2*tl_vstr(i,j,it2,linp)
184 END DO
185 END DO
186# endif
187# if defined ADJUST_STFLUX && defined SOLVE3D
188!
189! Adjust surface tracer fluxes. Interpolate between surface forcing
190! increments, if appropriate.
191!
192 DO itrc=1,nt(ng)
193 IF (lstflux(itrc,ng)) THEN
194 DO j=jstrr,jendr
195 DO i=istrr,iendr
196 tl_stflx(i,j,itrc)=tl_stflx(i,j,itrc)+ &
197 & fac1*tl_tflux(i,j,it1,linp,itrc)+ &
198 & fac2*tl_tflux(i,j,it2,linp,itrc)
199 END DO
200 END DO
201 END IF
202 END DO
203# endif
204!
205 RETURN
integer, dimension(:), allocatable nt
Definition mod_param.F:489
integer, dimension(:), allocatable iic
real(dp), dimension(:), allocatable dt
logical, dimension(:,:), allocatable lstflux
integer, dimension(:), allocatable nfrec
real(dp), dimension(:), allocatable time
integer, dimension(:), allocatable nsff
integer, dimension(:), allocatable ntstart
real(dp), dimension(:,:), allocatable sf_time

References mod_scalars::dt, mod_scalars::iic, mod_scalars::lstflux, mod_scalars::nsff, mod_scalars::ntstart, mod_scalars::sf_time, and mod_scalars::time.

Referenced by rp_frc_adjust().

Here is the caller graph for this function: