My Project
Functions/Subroutines
advection_edge_gcy.f90 File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine advection_edge_gcy (XFLUX, YFLUX)
 

Function/Subroutine Documentation

◆ advection_edge_gcy()

subroutine advection_edge_gcy ( real(sp), dimension(0:nt,kb), intent(out)  XFLUX,
real(sp), dimension(0:nt,kb), intent(out)  YFLUX 
)

Definition at line 42 of file advection_edge_gcy.f90.

42 !==============================================================================|
43 ! Calculate the Advection and Diffusion Terms of 3D Velocity Field |
44 ! These Terms will be vertically integrated to form the Mean Terms in |
45 ! the Gx and Gy Terms of the External Mode Equation |
46 ! Ghost cell boundary conditions are used here |
47 !==============================================================================|
48 
49  USE mod_utils
50  USE all_vars
51  USE bcs
52  USE mod_spherical
53  USE mod_northpole
54  USE mod_wd
55 
56  IMPLICIT NONE
57  REAL(SP), INTENT(OUT), DIMENSION(0:NT,KB) :: XFLUX,YFLUX
58  REAL(SP) :: DIJ
59  REAL(SP) :: COFA1,COFA2,COFA3,COFA4,COFA5,COFA6,COFA7,COFA8
60  REAL(SP) :: XADV,YADV,TXXIJ,TYYIJ,TXYIJ,UN
61  REAL(SP) :: VISCOF,VISCOF1,VISCOF2,TEMP,TPA,TPB
62  REAL(SP) :: XIJA,YIJA,XIJB,YIJB,UIJ,VIJ
63  REAL(SP) :: FACT,FM1
64  INTEGER :: I,IA,IB,J1,J2,K1,K2,K3,K4,K5,K6,K,II,J,I1,I2
65  REAL(SP) :: ISWETTMP
66 
67  REAL(SP) :: UIJ1,VIJ1,UIJ2,VIJ2,FXX,FYY
68  REAL(SP) :: UK1(KB),UK2(KB),UK3(KB),UK4(KB),UK5(KB),UK6(KB), &
69  VK1(KB),VK2(KB),VK3(KB),VK4(KB),VK5(KB),VK6(KB)
70 !------------------------------------------------------------------------------|
71 
72  if(dbg_set(dbg_sbr)) write(ipt,*) "Start: advection_edge_gcy.F"
73 
74 
75  SELECT CASE(horizontal_mixing_type)
76  CASE ('closure')
77  fact = 1.0_sp
78  fm1 = 0.0_sp
79  CASE('constant')
80  fact = 0.0_sp
81  fm1 = 1.0_sp
82  CASE DEFAULT
83  CALL fatal_error("UNKNOW HORIZONTAL MIXING TYPE:",&
84  & trim(horizontal_mixing_type) )
85  END SELECT
86 
87 !
88 !--Initialize Variables--------------------------------------------------------|
89 !
90  xflux = 0.0_sp
91  yflux = 0.0_sp
92 
93 !
94 !--Loop Over Edges and Accumulate Fluxes-For Each Element----------------------|
95 !
96 
97  DO i=1,ne
98  ia=iec(i,1)
99  ib=iec(i,2)
100  IF(iswetct(ia)*iswetc(ia) == 1 .OR. iswetct(ib)*iswetc(ib) == 1)THEN
101  j1=ienode(i,1)
102  j2=ienode(i,2)
103 
104  k1=nbe(ia,1)
105  k2=nbe(ia,2)
106  k3=nbe(ia,3)
107  k4=nbe(ib,1)
108  k5=nbe(ib,2)
109  k6=nbe(ib,3)
110  xija=xijc(i)-xc(ia)
111  yija=yijc(i)-yc(ia)
112  xijb=xijc(i)-xc(ib)
113  yijb=yijc(i)-yc(ib)
114 
115  uk1 = u(k1,:)
116  uk2 = u(k2,:)
117  uk3 = u(k3,:)
118  uk4 = u(k4,:)
119  uk5 = u(k5,:)
120  uk6 = u(k6,:)
121  vk1 = v(k1,:)
122  vk2 = v(k2,:)
123  vk3 = v(k3,:)
124  vk4 = v(k4,:)
125  vk5 = v(k5,:)
126  vk6 = v(k6,:)
127 
128  IF(k1 == 0) CALL ghostuv3(ia,1,uk1,vk1)
129  IF(k2 == 0) CALL ghostuv3(ia,2,uk2,vk2)
130  IF(k3 == 0) CALL ghostuv3(ia,3,uk3,vk3)
131  IF(k4 == 0) CALL ghostuv3(ib,1,uk4,vk4)
132  IF(k5 == 0) CALL ghostuv3(ib,2,uk5,vk5)
133  IF(k6 == 0) CALL ghostuv3(ib,3,uk6,vk6)
134 
135  DO k=1,kbm1
136 
137  dij= 0.5_sp*(dt(j1)*dz(j1,k)+dt(j2)*dz(j2,k))
138 
139  !!FORM THE LEFT FLUX
140  cofa1=a1u(ia,1)*u(ia,k)+a1u(ia,2)*uk1(k)+a1u(ia,3)*uk2(k)+a1u(ia,4)*uk3(k)
141  cofa2=a2u(ia,1)*u(ia,k)+a2u(ia,2)*uk1(k)+a2u(ia,3)*uk2(k)+a2u(ia,4)*uk3(k)
142  cofa5=a1u(ia,1)*v(ia,k)+a1u(ia,2)*vk1(k)+a1u(ia,3)*vk2(k)+a1u(ia,4)*vk3(k)
143  cofa6=a2u(ia,1)*v(ia,k)+a2u(ia,2)*vk1(k)+a2u(ia,3)*vk2(k)+a2u(ia,4)*vk3(k)
144  uij1=u(ia,k)+cofa1*xija+cofa2*yija
145  vij1=v(ia,k)+cofa5*xija+cofa6*yija
146 
147  !!FORM THE RIGHT FLUX
148  cofa3=a1u(ib,1)*u(ib,k)+a1u(ib,2)*uk4(k)+a1u(ib,3)*uk5(k)+a1u(ib,4)*uk6(k)
149  cofa4=a2u(ib,1)*u(ib,k)+a2u(ib,2)*uk4(k)+a2u(ib,3)*uk5(k)+a2u(ib,4)*uk6(k)
150  cofa7=a1u(ib,1)*v(ib,k)+a1u(ib,2)*vk4(k)+a1u(ib,3)*vk5(k)+a1u(ib,4)*vk6(k)
151  cofa8=a2u(ib,1)*v(ib,k)+a2u(ib,2)*vk4(k)+a2u(ib,3)*vk5(k)+a2u(ib,4)*vk6(k)
152  uij2=u(ib,k)+cofa3*xijb+cofa4*yijb
153  vij2=v(ib,k)+cofa7*xijb+cofa8*yijb
154 
155  !!COMPUTE THE NORMAL VELOCITY ACROSS THE EDGE
156  uij=0.5_sp*(uij1+uij2)
157  vij=0.5_sp*(vij1+vij2)
158  un=vij*dltxc(i) - uij*dltyc(i)
159 
160  viscof1=art(ia)*sqrt(cofa1**2+cofa6**2+0.5_sp*(cofa2+cofa5)**2)
161  viscof2=art(ib)*sqrt(cofa3**2+cofa8**2+0.5_sp*(cofa4+cofa7)**2)
162 
163 ! VISCOF = HORCON*(FACT*0.5_SP*(VISCOF1+VISCOF2)/HPRNU + FM1)
164 ! VISCOF = HORCON*(FACT*0.5_SP*(VISCOF1+VISCOF2) + FM1)
165  ! David moved HPRNU and added HVC
166  viscof=(fact*0.5_sp*(viscof1*cc_hvc(ia)+viscof2*cc_hvc(ib)) + fm1*0.5_sp*(cc_hvc(ia)+cc_hvc(ib)))/hprnu
167 
168  txxij=(cofa1+cofa3)*viscof
169  tyyij=(cofa6+cofa8)*viscof
170  txyij=0.5_sp*(cofa2+cofa4+cofa5+cofa7)*viscof
171  fxx=dij*(txxij*dltyc(i)-txyij*dltxc(i))
172  fyy=dij*(txyij*dltyc(i)-tyyij*dltxc(i))
173 
174 
175  !!UPWIND THE ADVECTIVE FLUX
176  xadv=dij*un*((1.0_sp-sign(1.0_sp,un))*uij2+(1.0_sp+sign(1.0_sp,un))*uij1)*0.5_sp
177  yadv=dij*un*((1.0_sp-sign(1.0_sp,un))*vij2+(1.0_sp+sign(1.0_sp,un))*vij1)*0.5_sp
178 
179 
180  !!COMPUTE BOUNDARY FLUX AUGMENTERS
181  tpa = float(1-isbc(i))*epor(ia)
182  tpb = float(1-isbc(i))*epor(ib)
183 
184 
185  !!ACCUMULATE THE FLUX
186  xflux(ia,k)=xflux(ia,k)+xadv*tpa+fxx*tpa
187  yflux(ia,k)=yflux(ia,k)+yadv*tpa+fyy*tpa
188  xflux(ib,k)=xflux(ib,k)-xadv*tpb-fxx*tpb
189  yflux(ib,k)=yflux(ib,k)-yadv*tpb-fyy*tpb
190 
191 
192  END DO
193  END IF
194  END DO
195 
196 
197  DO i=1,n
198  iswettmp = iswetct(i)*iswetc(i)
199  DO k=1,kbm1
200  xflux(i,k) = xflux(i,k)*iswettmp
201  yflux(i,k) = yflux(i,k)*iswettmp
202  END DO
203  END DO
204 
205 
206 !
207 !--Boundary Conditions on Flux-------------------------------------------------|
208 !
209  DO i=1,n
210  IF(isbce(i) == 2)THEN
211  DO k=1,kbm1
212  xflux(i,k)=0.0_sp
213  yflux(i,k)=0.0_sp
214  END DO
215  END IF
216  END DO
217 
218 !
219 !--Adjust Flux for Fresh Water Inflow------------------------------------------|
220 !
221 
222  IF(numqbc > 0) THEN
223  IF(river_inflow_location == 'node') THEN
224  DO ii=1,numqbc
225  j=inodeq(ii)
226  i1=nbve(j,1)
227  i2=nbve(j,ntve(j))
228  DO k=1,kbm1
229  vlctyq(ii)=qdis(ii)/qarea(ii)
230 ! TEMP=0.5_SP*QDIS(II)*VQDIST(II,K)*VLCTYQ(II)
231  temp=0.5_sp*qdis(ii)*vqdist(ii,k)*vqdist(ii,k)*vlctyq(ii)/dz(j,k)
232 ! XFLUX(I1,K)=XFLUX(I1,K)-TEMP/DZ(J,K)*COS(ANGLEQ(II))
233 ! XFLUX(I2,K)=XFLUX(I2,K)-TEMP/DZ(J,K)*COS(ANGLEQ(II))
234 ! YFLUX(I1,K)=YFLUX(I1,K)-TEMP/DZ(J,K)*SIN(ANGLEQ(II))
235 ! YFLUX(I2,K)=YFLUX(I2,K)-TEMP/DZ(J,K)*SIN(ANGLEQ(II))
236  xflux(i1,k)=xflux(i1,k)-temp*cos(angleq(ii))
237  xflux(i2,k)=xflux(i2,k)-temp*cos(angleq(ii))
238  yflux(i1,k)=yflux(i1,k)-temp*sin(angleq(ii))
239  yflux(i2,k)=yflux(i2,k)-temp*sin(angleq(ii))
240  END DO
241  END DO
242  ELSE IF(river_inflow_location == 'edge') THEN
243  DO ii=1,numqbc
244  i1=icellq(ii)
245  DO k=1,kbm1
246  vlctyq(ii)=qdis(ii)/qarea(ii)
247 ! TEMP=QDIS(II)*VQDIST(II,K)*VLCTYQ(II)
248  temp=qdis(ii)*vqdist(ii,k)*vqdist(ii,k)*vlctyq(ii)/dz1(i1,k)
249 ! XFLUX(I1,K)=XFLUX(I1,K)-TEMP/DZ1(I1,K)*COS(ANGLEQ(II))
250 ! YFLUX(I1,K)=YFLUX(I1,K)-TEMP/DZ1(I1,K)*SIN(ANGLEQ(II))
251  xflux(i1,k)=xflux(i1,k)-temp*cos(angleq(ii))
252  yflux(i1,k)=yflux(i1,k)-temp*sin(angleq(ii))
253  END DO
254  END DO
255  END IF
256  END IF
257 
258  if(dbg_set(dbg_sbr)) write(ipt,*) "End: advection_edge_gcy.F"
259 
260  RETURN
integer, dimension(:,:), allocatable, target ienode
Definition: mod_main.f90:1029
real(sp), dimension(:), allocatable, target epor
Definition: mod_main.f90:1056
real(sp), dimension(:), allocatable, target art
Definition: mod_main.f90:1009
real(sp), dimension(:,:), allocatable, target v
Definition: mod_main.f90:1269
real(sp), dimension(:,:), allocatable, target vqdist
Definition: mod_main.f90:1217
logical function dbg_set(vrb)
Definition: mod_utils.f90:182
real(sp), dimension(:), allocatable, target qdis
Definition: mod_main.f90:1220
real(sp), dimension(:), allocatable, target yc
Definition: mod_main.f90:1004
real(sp), dimension(:), allocatable cc_hvc
Definition: mod_main.f90:1302
real(sp), dimension(:), allocatable, target dltxc
Definition: mod_main.f90:1040
subroutine ghostuv3(I, JJ, UAKK, VAKK)
Definition: ghostuv.f90:86
real(sp), dimension(:,:), allocatable, target a1u
Definition: mod_main.f90:1331
integer, dimension(:), allocatable iswetct
Definition: mod_wd.f90:54
real(sp), dimension(:,:), allocatable, target u
Definition: mod_main.f90:1268
integer, dimension(:), allocatable, target isbc
Definition: mod_main.f90:1026
real(sp), dimension(:), allocatable, target angleq
Definition: mod_main.f90:1228
integer, dimension(:,:), allocatable, target iec
Definition: mod_main.f90:1028
integer, dimension(:), allocatable, target ntve
Definition: mod_main.f90:1022
integer, dimension(:,:), allocatable, target nbe
Definition: mod_main.f90:1020
real(sp), dimension(:), allocatable, target xijc
Definition: mod_main.f90:1044
real(sp), dimension(:), allocatable, target qarea
Definition: mod_main.f90:1226
real(sp), dimension(:), allocatable, target dltyc
Definition: mod_main.f90:1041
real(sp), dimension(:,:), allocatable, target dz
Definition: mod_main.f90:1092
integer, dimension(:), allocatable, target icellq
Definition: mod_main.f90:1215
real(sp), dimension(:), allocatable, target yijc
Definition: mod_main.f90:1045
real(sp), dimension(:,:), allocatable, target a2u
Definition: mod_main.f90:1332
integer, dimension(:), allocatable iswetc
Definition: mod_wd.f90:52
integer, dimension(:,:), allocatable, target nbve
Definition: mod_main.f90:1034
subroutine fatal_error(ER1, ER2, ER3, ER4)
Definition: mod_utils.f90:230
real(sp), dimension(:), allocatable, target xc
Definition: mod_main.f90:1003
real(sp), dimension(:,:), allocatable, target dz1
Definition: mod_main.f90:1096
integer, dimension(:), allocatable, target isbce
Definition: mod_main.f90:1027
real(sp), dimension(:), allocatable, target vlctyq
Definition: mod_main.f90:1229
integer, parameter dbg_sbr
Definition: mod_utils.f90:69
integer, dimension(:), allocatable, target inodeq
Definition: mod_main.f90:1214
real(sp), dimension(:), allocatable, target dt
Definition: mod_main.f90:1133
Here is the call graph for this function: