My Project
external_step.f90
Go to the documentation of this file.
1 
2 
3 
4 
5 
6 
7 
8 
9 
10 
11 
12 !/===========================================================================/
13 ! Copyright (c) 2007, The University of Massachusetts Dartmouth
14 ! Produced at the School of Marine Science & Technology
15 ! Marine Ecosystem Dynamics Modeling group
16 ! All rights reserved.
17 !
18 ! FVCOM has been developed by the joint UMASSD-WHOI research team. For
19 ! details of authorship and attribution of credit please see the FVCOM
20 ! technical manual or contact the MEDM group.
21 !
22 !
23 ! This file is part of FVCOM. For details, see http://fvcom.smast.umassd.edu
24 ! The full copyright notice is contained in the file COPYRIGHT located in the
25 ! root directory of the FVCOM code. This original header must be maintained
26 ! in all distributed versions.
27 !
28 ! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
29 ! AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO,
30 ! THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
31 ! PURPOSE ARE DISCLAIMED.
32 !
33 !/---------------------------------------------------------------------------/
34 ! CVS VERSION INFORMATION
35 ! $Id$
36 ! $Name$
37 ! $Revision$
38 !/===========================================================================/
39 
40 SUBROUTINE external_step
41 
42  USE mod_nesting
43  USE mod_utils
44  USE all_vars
45  USE mod_time
46  USE mod_obcs
47  USE mod_wd
48 
49 
50  USE mod_ice
51  USE mod_ice2d
52 
53 
54 
55 
56  IMPLICIT NONE
57  REAL(SP) :: TMP
58  INTEGER :: K, I, J, JN, J1,i1,i2
59 
60 !------------------------------------------------------------------------------|
61 
62  if(dbg_set(dbg_sbr)) write(ipt,*) "Start: external_step"
63 
64  !! David for VISIT
65  !!=======================================================
66  !!=======================================================
67 
68 
69 !----SET RAMP FACTOR TO EASE SPINUP--------------------------------------------!
70  IF(iramp /= 0) THEN
71  tmp = real(iint-1,sp)+real(iext,sp)/real(isplit,sp)
72  ramp=tanh(tmp/real(iramp,sp))
73  ELSE
74  ramp = 1.0_sp
75  END IF
76 
77 !
78 !------SURFACE BOUNDARY CONDITIONS FOR EXTERNAL MODEL--------------------------!
79 !
80  CALL bcond_gcn(9,0)
81 
82 
83 
84 !
85 !------SAVE VALUES FROM CURRENT TIME STEP--------------------------------------!
86 !
87  elrk1 = el1
88  elrk = el
89  uark = ua
90  vark = va
91 
92 
93 ! New Open Boundary Condition ----5
94 
95 
96 !
97 !------BEGIN MAIN LOOP OVER EXTERNAL MODEL 4 STAGE RUNGE-KUTTA INTEGRATION-----!
98 !
99 
100  DO k=1,4
101 
102  rktime = exttime + imdte * (alpha_rk(k) - 1.0_dp)
103 
104 ! CALL PRINT_REAL_TIME(RKTIME,IPT,"RUNGE-KUTTA")
105 
106 
107 ! New Open Boundary Condition ----6
108 
109 
110 !FREE SURFACE AMPLITUDE UPDATE --> ELF
111  CALL extel_edge(k)
112 
113 
114 
115 
116 
117  ! New Open Boundary Condition ----7
118 
119  ! VALUES FOR THE OPEN BOUNDARY ARE ONLY UPDATED IN THE LOCAL DOMAIN
120  ! THE HALO IS NOT SET HERE
121  CALL bcond_gcn(1,k)
122 
123  IF(nesting_on )THEN
124  CALL set_var(exttime,el=elf)
125  END IF
126 
127  DO i=1,ibcn(1)
128  jn = obc_lst(1,i)
129  j=i_obc_n(jn)
130  elf(j)=elrk(j)+alpha_rk(k)*(elf(j)-elrk(j))
131  END DO
132 
133 
134  ! DAVID ADDED THIS EXCHANGE CALL:
135  ! IT SEEMS LIKELY THAT THE HALO VALUES OF ELF WILL BE USED
136  ! BEFORE THEY ARE SET CORRECTLY OTHERWISE
137 
138 !---------------For Dam Model-----------------------
139 ! Jadon
140 
141 
142  CALL n2e2d(elf,elf1)
143 
144  IF(wetting_drying_on)CALL wet_judge
145 
146  CALL flux_obn(k)
147 
148  !CALCULATE ADVECTIVE, DIFFUSIVE, AND BAROCLINIC MODES --> UAF ,VAF
149  CALL advave_edge_gcn(advua,advva) !Compute Ext Mode Adv/Diff
150 
151 
152  CALL extuv_edge(k)
153 
154 
155  CALL bcond_gcn(2,k)
156 
157 
158  IF(nesting_on )THEN
159  CALL set_var(exttime,ua=uaf)
160  CALL set_var(exttime,va=vaf)
161  END IF
162 
163 
164 
165 
166  !UPDATE WATER SURFACE ELEVATION
168 
169  el = elf
170  el1 = elf1
171 
172 
173 
174  !!INTERPOLATE DEPTH FROM NODE-BASED TO ELEMENT-BASED VALUES
175  CALL n2e2d(el,el1)
176 
177  !UPDATE DEPTH AND VERTICALLY AVERAGED VELOCITY FIELD
178  d = h + el
179  d1 = h1 + el1
180  ua = uaf
181  va = vaf
182  dtfa = d
183 
184 
185  ! New Open Boundary Condition ----8
186 
187 
188  !!ENSURE ALL CELLS ARE WET IN NO FLOOD/DRY CASE
189 
190  !EXCHANGE ELEMENT-BASED VALUES ACROSS THE INTERFACE
191 
192  !SAVE VALUES FOR 3D MOMENTUM CORRECTION AND UPDATE
193  IF(k == 3)THEN
194  uard = uard + ua*d1
195  vard = vard + va*d1
196  egf = egf + el/isplit
197 
198 
199 
200 !!# if defined (1)
201 !! UARDS = UARDS + UAS*D1
202 !! VARDS = VARDS + VAS*D1
203 !!# endif
204  END IF
205 
206  !CALCULATE VALUES USED FOR SALINITY/TEMP BOUNDARY CONDITIONS
207  IF(k == 4.AND.iobcn > 0) THEN
208  DO i=1,iobcn
209  j=i_obc_n(i)
210  tmp=-(elf(j)-elrk(j))*art1(j)/dte-xflux_obcn(i)
211  uard_obcn(i)=uard_obcn(i)+tmp/float(isplit)
212  END DO
213  END IF
214 !end !defined (TWO_D_MODEL)
215 
216 
217  !UPDATE WET/DRY FACTORS
218  IF(wetting_drying_on)CALL wd_update(1)
219 
220  END DO !! END RUNGE-KUTTA LOOP
221 
222 
223  if(dbg_set(dbg_sbr)) write(ipt,*) "End: external_step"
224 
225 END SUBROUTINE external_step
real(sp), dimension(:), allocatable, target elrk1
Definition: mod_main.f90:1121
real(sp), dimension(:), allocatable, target va
Definition: mod_main.f90:1104
real(sp), dimension(:), allocatable, target d
Definition: mod_main.f90:1132
real(sp), dimension(:), allocatable, target d1
Definition: mod_main.f90:1116
subroutine flux_obn(K)
Definition: mod_obcs.f90:829
real(sp), dimension(:), allocatable, target elrk
Definition: mod_main.f90:1138
real(sp), dimension(:), allocatable, target h
Definition: mod_main.f90:1131
real(sp), dimension(:), allocatable, target dtfa
Definition: mod_main.f90:1124
real(sp), dimension(:), allocatable, target uark
Definition: mod_main.f90:1108
real(sp), dimension(:), allocatable, target el
Definition: mod_main.f90:1134
real(sp), dimension(:), allocatable, target advua
Definition: mod_main.f90:1256
logical function dbg_set(vrb)
Definition: mod_utils.f90:182
real(sp), dimension(:), allocatable, target art1
Definition: mod_main.f90:1010
integer, dimension(:,:), allocatable obc_lst
Definition: mod_obcs.f90:84
real(sp), dimension(:), allocatable, target egf
Definition: mod_main.f90:1136
subroutine bcond_gcn(IDX, K_RK)
Definition: bcond_gcn.f90:58
subroutine set_var(NOW, UA, VA, EL, U, V, S1, T1, HYW)
logical nesting_on
real(sp), dimension(:), allocatable, target vard
Definition: mod_main.f90:1111
real(sp), dimension(:), allocatable, target el1
Definition: mod_main.f90:1118
real(sp), dimension(:), allocatable, target uard
Definition: mod_main.f90:1110
real(sp), dimension(:), allocatable uard_obcn
Definition: mod_obcs.f90:112
real(sp), dimension(:), allocatable, target vaf
Definition: mod_main.f90:1106
real(sp), dimension(:), allocatable, target elf
Definition: mod_main.f90:1140
real(sp), dimension(:), allocatable xflux_obcn
Definition: mod_obcs.f90:111
integer, dimension(5) ibcn
Definition: mod_obcs.f90:82
real(sp), dimension(:), allocatable, target ua
Definition: mod_main.f90:1103
subroutine n2e2d(NVAR, EVAR)
Definition: mod_main.f90:1390
subroutine wet_judge
Definition: mod_wd.f90:180
subroutine assign_elm1_to_elm2
Definition: mod_obcs.f90:156
subroutine extel_edge(K)
Definition: extel_edge.f90:44
real(sp), dimension(:), allocatable, target h1
Definition: mod_main.f90:1115
subroutine external_step
real(sp), dimension(:), allocatable, target uaf
Definition: mod_main.f90:1105
real(sp), dimension(:), allocatable, target elf1
Definition: mod_main.f90:1123
subroutine advave_edge_gcn(XFLUX, YFLUX)
subroutine extuv_edge(K)
Definition: extuv_edge.f90:45
integer, parameter dbg_sbr
Definition: mod_utils.f90:69
subroutine wd_update(INCASE)
Definition: mod_wd.f90:253
real(sp), dimension(:), allocatable, target vark
Definition: mod_main.f90:1109
real(sp), dimension(:), allocatable, target advva
Definition: mod_main.f90:1257