/[MITgcm]/MITgcm/pkg/ctrl/ctrl_getobcsw.F
ViewVC logotype

Contents of /MITgcm/pkg/ctrl/ctrl_getobcsw.F

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


Revision 1.17 - (show annotations) (download)
Thu Oct 9 00:49:26 2014 UTC (10 years, 9 months ago) by gforget
Branch: MAIN
CVS Tags: checkpoint65z, checkpoint65x, checkpoint65y, checkpoint65r, checkpoint65s, checkpoint65p, checkpoint65q, checkpoint65v, checkpoint65w, checkpoint65t, checkpoint65u, checkpoint65j, checkpoint65k, checkpoint65h, checkpoint65i, checkpoint65n, checkpoint65l, checkpoint65m, checkpoint65f, checkpoint65g, checkpoint66g, checkpoint66f, checkpoint66e, checkpoint66d, checkpoint66c, checkpoint66b, checkpoint66a, checkpoint66o, checkpoint66n, checkpoint66m, checkpoint66l, checkpoint66k, checkpoint66j, checkpoint66i, checkpoint66h, checkpoint65o, HEAD
Changes since 1.16: +2 -1 lines
- pkg/ctrl/CTRL_OBCS.h (new) : regroup all obcs ctrl variables.
- pkg/ctrl/ctrl.h, ctrl_dummy.h, ctrl_weights.h : rm obcs
  ctrl variables (now all in CTRL_OBCS.h).

- pkg/ctrl/ctrl_getobcse.F, ctrl_getobcsn.F, ctrl_getobcss.F,
  ctrl_getobcsw.F, ctrl_getrec.F, ctrl_init.F, ctrl_init_obcs_variables.F,
  ctrl_init_wet.F, ctrl_mask_set_xz.F, ctrl_mask_set_yz.F,
  ctrl_pack.F, ctrl_unpack.F, ctrl_readparms.F,
  ctrl_set_pack_xz.F, ctrl_set_pack_yz.F, ctrl_set_unpack_xz.F,
  ctrl_set_unpack_yz.F : add CPP brackets and CTRL_OBCS.h

- pkg/ctrl/ctrl_pack.F, ctrl_unpack.F : add CPP brackets

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

  ViewVC Help
Powered by ViewVC 1.1.22