/[MITgcm]/MITgcm_contrib/SOSE/BoxAdj/code_ad/ctrl_getobcse.F
ViewVC logotype

Contents of /MITgcm_contrib/SOSE/BoxAdj/code_ad/ctrl_getobcse.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.6 - (show annotations) (download)
Fri Mar 18 18:43:05 2011 UTC (14 years, 9 months ago) by mmazloff
Branch: MAIN
Changes since 1.5: +158 -382 lines
bringing ALLOW_OBCS_CONTROL_MODES up to date

1 C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_getobcse.F,v 1.10 2011/03/14 17:08:00 mlosch Exp $
2 C $Name: $
3
4 #include "CTRL_CPPOPTIONS.h"
5 #ifdef ALLOW_OBCS
6 # include "OBCS_OPTIONS.h"
7 #endif
8
9
10 subroutine ctrl_getobcse(
11 I mytime,
12 I myiter,
13 I mythid
14 & )
15
16 c ==================================================================
17 c SUBROUTINE ctrl_getobcse
18 c ==================================================================
19 c
20 c o Get eastern obc of the control vector and add it
21 c to dyn. fields
22 c
23 c started: heimbach@mit.edu, 29-Aug-2001
24 c
25 c ==================================================================
26 c SUBROUTINE ctrl_getobcse
27 c ==================================================================
28
29 implicit none
30
31 #ifdef ALLOW_OBCSE_CONTROL
32
33 c == global variables ==
34
35 #include "EEPARAMS.h"
36 #include "SIZE.h"
37 #include "PARAMS.h"
38 #include "GRID.h"
39 #include "OBCS.h"
40
41 #include "ctrl.h"
42 #include "ctrl_dummy.h"
43 #include "optim.h"
44
45 c == routine arguments ==
46
47 _RL mytime
48 integer myiter
49 integer mythid
50
51 c == local variables ==
52
53 integer bi,bj
54 integer i,j,k
55 integer itlo,ithi
56 integer jtlo,jthi
57 integer jmin,jmax
58 integer imin,imax
59 integer ilobcse
60 integer iobcs
61
62 _RL dummy
63 _RL obcsefac
64 logical obcsefirst
65 logical obcsechanged
66 integer obcsecount0
67 integer obcsecount1
68 integer ip1
69
70 cgg _RL maskyz (1-oly:sny+oly,nr,nsx,nsy)
71 _RL tmpfldyz (1-oly:sny+oly,nr,nsx,nsy)
72
73 logical doglobalread
74 logical ladinit
75
76 character*(80) fnameobcse
77
78 cmm( modes:
79 integer nk,nz
80 _RL tmpz (nr,nsx,nsy)
81 _RL stmp
82
83 c == external functions ==
84
85 integer ilnblnk
86 external ilnblnk
87
88
89 c == end of interface ==
90
91 jtlo = mybylo(mythid)
92 jthi = mybyhi(mythid)
93 itlo = mybxlo(mythid)
94 ithi = mybxhi(mythid)
95 jmin = 1-oly
96 jmax = sny+oly
97 imin = 1-olx
98 imax = snx+olx
99 ip1 = 0
100
101 c-- Now, read the control vector.
102 doglobalread = .false.
103 ladinit = .false.
104
105 if (optimcycle .ge. 0) then
106 ilobcse=ilnblnk( xx_obcse_file )
107 write(fnameobcse(1:80),'(2a,i10.10)')
108 & xx_obcse_file(1:ilobcse), '.', optimcycle
109 endif
110
111 c-- Get the counters, flags, and the interpolation factor.
112 call ctrl_get_gen_rec(
113 I xx_obcsestartdate, xx_obcseperiod,
114 O obcsefac, obcsefirst, obcsechanged,
115 O obcsecount0,obcsecount1,
116 I mytime, myiter, mythid )
117
118 do iobcs = 1,nobcs
119
120 if ( obcsefirst ) then
121 call active_read_yz( fnameobcse, tmpfldyz,
122 & (obcsecount0-1)*nobcs+iobcs,
123 & doglobalread, ladinit, optimcycle,
124 & mythid, xx_obcse_dummy )
125
126 do bj = jtlo,jthi
127 do bi = itlo,ithi
128 #ifdef ALLOW_OBCS_CONTROL_MODES
129 if (iobcs .gt. 2) then
130 do j = jmin,jmax
131 i = OB_Ie(j,bi,bj)
132 cih Determine number of open vertical layers.
133 nz = 0
134 do k = 1,Nr
135 if (iobcs .eq. 3) then
136 nz = nz + maskW(i+ip1,j,k,bi,bj)
137 else
138 nz = nz + maskS(i,j,k,bi,bj)
139 endif
140 end do
141 cih Compute absolute velocities from the barotropic-baroclinic modes.
142 do k = 1,Nr
143 if (k.le.nz) then
144 stmp = 0.
145 do nk = 1,nz
146 stmp = stmp +
147 & modesv(k,nk,nz)*tmpfldyz(j,nk,bi,bj)
148 end do
149 tmpz(k,bi,bj) = stmp
150 else
151 tmpz(k,bi,bj) = 0.
152 end if
153 enddo
154 do k = 1,Nr
155 if (iobcs .eq. 3) then
156 tmpfldyz(j,k,bi,bj) = tmpz(k,bi,bj)
157 & *recip_hFacW(i+ip1,j,k,bi,bj)
158 else
159 tmpfldyz(j,k,bi,bj) = tmpz(k,bi,bj)
160 & *recip_hFacS(i,j,k,bi,bj)
161 endif
162 end do
163 enddo
164 endif
165 #endif
166 do k = 1,nr
167 do j = jmin,jmax
168 xx_obcse1(j,k,bi,bj,iobcs) = tmpfldyz (j,k,bi,bj)
169 cgg & * maskyz (j,k,bi,bj)
170 enddo
171 enddo
172 enddo
173 enddo
174 endif
175
176 if ( (obcsefirst) .or. (obcsechanged)) then
177
178 do bj = jtlo,jthi
179 do bi = itlo,ithi
180 do j = jmin,jmax
181 do k = 1,nr
182 xx_obcse0(j,k,bi,bj,iobcs) = xx_obcse1(j,k,bi,bj,iobcs)
183 tmpfldyz (j,k,bi,bj) = 0. _d 0
184 enddo
185 enddo
186 enddo
187 enddo
188
189 call active_read_yz( fnameobcse, tmpfldyz,
190 & (obcsecount1-1)*nobcs+iobcs,
191 & doglobalread, ladinit, optimcycle,
192 & mythid, xx_obcse_dummy )
193
194 do bj = jtlo,jthi
195 do bi = itlo,ithi
196 #ifdef ALLOW_OBCS_CONTROL_MODES
197 if (iobcs .gt. 2) then
198 do j = jmin,jmax
199 i = OB_Ie(j,bi,bj)
200 cih Determine number of open vertical layers.
201 nz = 0
202 do k = 1,Nr
203 if (iobcs .eq. 3) then
204 nz = nz + maskW(i+ip1,j,k,bi,bj)
205 else
206 nz = nz + maskS(i,j,k,bi,bj)
207 endif
208 end do
209 cih Compute absolute velocities from the barotropic-baroclinic modes.
210 do k = 1,Nr
211 if (k.le.nz) then
212 stmp = 0.
213 do nk = 1,nz
214 stmp = stmp +
215 & modesv(k,nk,nz)*tmpfldyz(j,nk,bi,bj)
216 end do
217 tmpz(k,bi,bj) = stmp
218 else
219 tmpz(k,bi,bj) = 0.
220 endif
221 enddo
222 do k = 1,Nr
223 if (iobcs .eq. 3) then
224 tmpfldyz(j,k,bi,bj) = tmpz(k,bi,bj)
225 & *recip_hFacW(i+ip1,j,k,bi,bj)
226 else
227 tmpfldyz(j,k,bi,bj) = tmpz(k,bi,bj)
228 & *recip_hFacS(i,j,k,bi,bj)
229 endif
230 end do
231 enddo
232 endif
233 #endif
234 do k = 1,nr
235 do j = jmin,jmax
236 xx_obcse1 (j,k,bi,bj,iobcs) = tmpfldyz (j,k,bi,bj)
237 cgg & * maskyz (j,k,bi,bj)
238 enddo
239 enddo
240 enddo
241 enddo
242 endif
243
244 c-- Add control to model variable.
245 do bj = jtlo,jthi
246 do bi = itlo,ithi
247 c-- Calculate mask for tracer cells (0 => land, 1 => water).
248 do k = 1,nr
249 do j = 1,sny
250 i = OB_Ie(j,bi,bj)
251 if (iobcs .EQ. 1) then
252 OBEt(j,k,bi,bj) = OBEt(j,k,bi,bj)
253 & + obcsefac *xx_obcse0(j,k,bi,bj,iobcs)
254 & + (1. _d 0 - obcsefac)*xx_obcse1(j,k,bi,bj,iobcs)
255 OBEt(j,k,bi,bj) = OBEt(j,k,bi,bj)
256 & *maskW(i+ip1,j,k,bi,bj)
257 else if (iobcs .EQ. 2) then
258 OBEs(j,k,bi,bj) = OBEs(j,k,bi,bj)
259 & + obcsefac *xx_obcse0(j,k,bi,bj,iobcs)
260 & + (1. _d 0 - obcsefac)*xx_obcse1(j,k,bi,bj,iobcs)
261 OBEs(j,k,bi,bj) = OBEs(j,k,bi,bj)
262 & *maskW(i+ip1,j,k,bi,bj)
263 else if (iobcs .EQ. 3) then
264 OBEu(j,k,bi,bj) = OBEu(j,k,bi,bj)
265 & + obcsefac *xx_obcse0(j,k,bi,bj,iobcs)
266 & + (1. _d 0 - obcsefac)*xx_obcse1(j,k,bi,bj,iobcs)
267 OBEu(j,k,bi,bj) = OBEu(j,k,bi,bj)
268 & *maskW(i+ip1,j,k,bi,bj)
269 else if (iobcs .EQ. 4) then
270 OBEv(j,k,bi,bj) = OBEv(j,k,bi,bj)
271 & + obcsefac *xx_obcse0(j,k,bi,bj,iobcs)
272 & + (1. _d 0 - obcsefac)*xx_obcse1(j,k,bi,bj,iobcs)
273 OBEv(j,k,bi,bj) = OBEv(j,k,bi,bj)
274 & *maskS(i,j,k,bi,bj)
275 endif
276 enddo
277 enddo
278 enddo
279 enddo
280
281 C-- End over iobcs loop
282 enddo
283
284 #else /* ALLOW_OBCSE_CONTROL undefined */
285
286 c == routine arguments ==
287
288 _RL mytime
289 integer myiter
290 integer mythid
291
292 c-- CPP flag ALLOW_OBCSE_CONTROL undefined.
293
294 #endif /* ALLOW_OBCSE_CONTROL */
295
296 end
297

  ViewVC Help
Powered by ViewVC 1.1.22