ROMS
Loading...
Searching...
No Matches
rp_t3drelax.F
Go to the documentation of this file.
1#include "cppdefs.h"
3
4#if defined TL_IOMS && defined RPM_RELAXATION && defined SOLVE3D
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 relaxes current representer tangent linear tracer !
14! type variables to previous Picard iteration solution (basic !
15! state) to improve stability and convergence. !
16! !
17!=======================================================================
18!
19 implicit none
20!
21 PRIVATE
22!
23 PUBLIC :: rp_t3drelax
24!
25 CONTAINS
26!
27!***********************************************************************
28 SUBROUTINE rp_t3drelax (ng, tile)
29!***********************************************************************
30!
31 USE mod_param
32 USE mod_grid
33 USE mod_ocean
34 USE mod_stepping
35!
36! Imported variable declarations.
37!
38 integer, intent(in) :: ng, tile
39!
40! Local variable declarations.
41!
42 character (len=*), parameter :: myfile = &
43 & __FILE__
44!
45# include "tile.h"
46!
47# ifdef PROFILE
48 CALL wclock_on (ng, irpm, 24, __line__, myfile)
49# endif
50 CALL rp_t3drelax_tile (ng, tile, &
51 & lbi, ubi, lbj, ubj, &
52 & imins, imaxs, jmins, jmaxs, &
53 & nrhs(ng), nnew(ng), &
54# ifdef MASKING
55 & grid(ng) % umask, &
56 & grid(ng) % vmask, &
57# endif
58 & grid(ng) % Hz, &
59 & grid(ng) % pmon_u, &
60 & grid(ng) % pnom_v, &
61 & grid(ng) % pm, &
62 & grid(ng) % pn, &
63 & ocean(ng) % t, &
64 & ocean(ng) % tl_t)
65# ifdef PROFILE
66 CALL wclock_off (ng, irpm, 24, __line__, myfile)
67# endif
68!
69 RETURN
70 END SUBROUTINE rp_t3drelax
71!
72!***********************************************************************
73 SUBROUTINE rp_t3drelax_tile (ng, tile, &
74 & LBi, UBi, LBj, UBj, &
75 & IminS, ImaxS, JminS, JmaxS, &
76 & nrhs, nnew, &
77# ifdef MASKING
78 & umask, vmask, &
79# endif
80 & Hz, &
81 & pmon_u, pnom_v, pm, pn, &
82 & t, tl_t)
83!***********************************************************************
84!
85 USE mod_param
86 USE mod_scalars
87!
88! Imported variable declarations.
89!
90 integer, intent(in) :: ng, tile
91 integer, intent(in) :: LBi, UBi, LBj, UBj
92 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
93 integer, intent(in) :: nrhs, nnew
94
95# ifdef ASSUMED_SHAPE
96# ifdef MASKING
97 real(r8), intent(in) :: umask(LBi:,LBj:)
98 real(r8), intent(in) :: vmask(LBi:,LBj:)
99# endif
100 real(r8), intent(in) :: Hz(LBi:,LBj:,:)
101 real(r8), intent(in) :: pmon_u(LBi:,LBj:)
102 real(r8), intent(in) :: pnom_v(LBi:,LBj:)
103 real(r8), intent(in) :: pm(LBi:,LBj:)
104 real(r8), intent(in) :: pn(LBi:,LBj:)
105 real(r8), intent(in) :: t(LBi:,LBj:,:,:,:)
106
107 real(r8), intent(inout) :: tl_t(LBi:,LBj:,:,:,:)
108# else
109# ifdef MASKING
110 real(r8), intent(in) :: umask(LBi:UBi,LBj:UBj)
111 real(r8), intent(in) :: vmask(LBi:UBi,LBj:UBj)
112# endif
113 real(r8), intent(in) :: Hz(LBi:UBi,LBj:UBj,N(ng))
114 real(r8), intent(in) :: pmon_u(LBi:UBi,LBj:UBj)
115 real(r8), intent(in) :: pnom_v(LBi:UBi,LBj:UBj)
116 real(r8), intent(in) :: pm(LBi:UBi,LBj:UBj)
117 real(r8), intent(in) :: pn(LBi:UBi,LBj:UBj)
118 real(r8), intent(in) :: t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
119
120 real(r8), intent(inout) :: tl_t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
121# endif
122!
123! Local variable declarations.
124!
125 integer :: i, itrc, j, k
126
127 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: FE
128 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: FX
129
130# include "set_bounds.h"
131!
132!----------------------------------------------------------------------
133! Compute horizontal diffusion relaxation of tracers between current
134! and previous representer tangent linear Picard iteration trajectory
135! (basic state).
136!----------------------------------------------------------------------
137!
138 DO itrc=1,nt(ng)
139 IF (tl_tdiff(itrc,ng).gt.0.0_r8) THEN
140 DO k=1,n(ng)
141!
142! Compute XI- and ETA-components of diffusive tracer flux (T m3/s).
143!
144 DO j=jstr,jend
145 DO i=istr,iend+1
146 fx(i,j)=0.5_r8*tl_tdiff(itrc,ng)*pmon_u(i,j)* &
147 & (hz(i,j,k)+hz(i-1,j,k))* &
148 & (tl_t(i ,j,k,nrhs,itrc)-t(i ,j,k,nrhs,itrc)- &
149 & tl_t(i-1,j,k,nrhs,itrc)+t(i-1,j,k,nrhs,itrc))
150# ifdef MASKING
151 fx(i,j)=fx(i,j)*umask(i,j)
152# endif
153 END DO
154 END DO
155 DO j=jstr,jend+1
156 DO i=istr,iend
157 fe(i,j)=0.5_r8*tl_tdiff(itrc,ng)*pnom_v(i,j)* &
158 & (hz(i,j,k)+hz(i,j-1,k))* &
159 & (tl_t(i,j ,k,nrhs,itrc)-t(i,j ,k,nrhs,itrc)- &
160 & tl_t(i,j-1,k,nrhs,itrc)+t(i,j-1,k,nrhs,itrc))
161# ifdef MASKING
162 fe(i,j)=fe(i,j)*vmask(i,j)
163# endif
164 END DO
165 END DO
166!
167! Time-step horizontal diffusion relaxation term (m Tunits).
168!
169 DO j=jstr,jend
170 DO i=istr,iend
171 tl_t(i,j,k,nnew,itrc)=tl_t(i,j,k,nnew,itrc)+ &
172 & dt(ng)*pm(i,j)*pn(i,j)* &
173 & (fx(i+1,j)-fx(i,j)+ &
174 & fe(i,j+1)-fe(i,j))
175 END DO
176 END DO
177 END DO
178 END IF
179 END DO
180!
181 RETURN
182 END SUBROUTINE rp_t3drelax_tile
183#endif
184 END MODULE rp_t3drelax_mod
type(t_grid), dimension(:), allocatable grid
Definition mod_grid.F:365
type(t_ocean), dimension(:), allocatable ocean
Definition mod_ocean.F:351
integer, parameter irpm
Definition mod_param.F:664
real(r8), dimension(:,:), allocatable tl_tdiff
real(dp), dimension(:), allocatable dt
integer, dimension(:), allocatable nrhs
integer, dimension(:), allocatable nnew
subroutine, public rp_t3drelax(ng, tile)
Definition rp_t3drelax.F:29
subroutine rp_t3drelax_tile(ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, nrhs, nnew, umask, vmask, hz, pmon_u, pnom_v, pm, pn, t, tl_t)
Definition rp_t3drelax.F:83
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