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

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

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


Revision 1.16 - (show annotations) (download)
Tue Sep 18 20:21:23 2012 UTC (11 years, 7 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint65b, checkpoint65c, checkpoint65a, checkpoint65d, checkpoint65e, checkpoint64, checkpoint65, checkpoint64y, checkpoint64x, checkpoint64z, checkpoint64q, checkpoint64p, checkpoint64s, checkpoint64r, checkpoint64u, checkpoint64t, checkpoint64w, checkpoint64v, checkpoint64i, checkpoint64h, checkpoint64k, checkpoint64j, checkpoint64m, checkpoint64l, checkpoint64o, checkpoint64n, checkpoint64a, checkpoint64c, checkpoint64b, checkpoint64e, checkpoint64d, checkpoint64g, checkpoint64f
Changes since 1.15: +4 -1 lines
use new parameter OB_indexNone for null index value (instead of hard-coded 0)

1 C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_getobcsn.F,v 1.15 2012/08/10 19:38:57 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_getobcsn(
10 I mytime,
11 I myiter,
12 I mythid
13 & )
14
15 c ==================================================================
16 c SUBROUTINE ctrl_getobcsn
17 c ==================================================================
18 c
19 c o Get northern 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_getobcsn
27 c ==================================================================
28
29 implicit none
30
31 c == global variables ==
32 #ifdef ALLOW_OBCSN_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 "optim.h"
44 #endif /* ALLOW_OBCSN_CONTROL */
45
46 c == routine arguments ==
47 _RL mytime
48 integer myiter
49 integer mythid
50
51 #ifdef ALLOW_OBCSN_CONTROL
52 c == local variables ==
53
54 integer bi,bj
55 integer i,j,k
56 integer itlo,ithi
57 integer jtlo,jthi
58 integer jmin,jmax
59 integer imin,imax
60 integer ilobcsn
61 integer iobcs
62 integer jp1
63
64 _RL dummy
65 _RL obcsnfac
66 logical obcsnfirst
67 logical obcsnchanged
68 integer obcsncount0
69 integer obcsncount1
70
71 cgg _RL maskxz (1-olx:snx+olx,nr,nsx,nsy)
72 _RL tmpfldxz (1-olx:snx+olx,nr,nsx,nsy)
73
74 logical doglobalread
75 logical ladinit
76
77 character*(80) fnameobcsn
78
79 #ifdef ALLOW_OBCS_CONTROL_MODES
80 integer nk,nz
81 _RL tmpz (nr,nsx,nsy)
82 _RL stmp
83 #endif
84
85 c == external functions ==
86
87 integer ilnblnk
88 external ilnblnk
89
90
91 c == end of interface ==
92
93 jtlo = mybylo(mythid)
94 jthi = mybyhi(mythid)
95 itlo = mybxlo(mythid)
96 ithi = mybxhi(mythid)
97 cgg jmin = 1-oly
98 cgg jmax = sny+oly
99 cgg imin = 1-olx
100 cgg imax = snx+olx
101
102 jmin = 1
103 jmax = sny
104 imin = 1
105 imax = snx
106 jp1 = 0
107
108 c-- Now, read the control vector.
109 doglobalread = .false.
110 ladinit = .false.
111
112 if (optimcycle .ge. 0) then
113 ilobcsn=ilnblnk( xx_obcsn_file )
114 write(fnameobcsn(1:80),'(2a,i10.10)')
115 & xx_obcsn_file(1:ilobcsn), '.', optimcycle
116 endif
117
118 c-- Get the counters, flags, and the interpolation factor.
119 call ctrl_get_gen_rec(
120 I xx_obcsnstartdate, xx_obcsnperiod,
121 O obcsnfac, obcsnfirst, obcsnchanged,
122 O obcsncount0,obcsncount1,
123 I mytime, myiter, mythid )
124
125 do iobcs = 1,nobcs
126 if ( obcsnfirst ) then
127 call active_read_xz( fnameobcsn, tmpfldxz,
128 & (obcsncount0-1)*nobcs+iobcs,
129 & doglobalread, ladinit, optimcycle,
130 & mythid, xx_obcsn_dummy )
131
132 do bj = jtlo,jthi
133 do bi = itlo,ithi
134 #ifdef ALLOW_OBCS_CONTROL_MODES
135 if (iobcs .gt. 2) then
136 do i = imin,imax
137 j = OB_Jn(i,bi,bj)
138 IF ( j.EQ.OB_indexNone ) j = 1
139 cih Determine number of open vertical layers.
140 nz = 0
141 do k = 1,Nr
142 if (iobcs .eq. 3) then
143 nz = nz + maskS(i,j+jp1,k,bi,bj)
144 else
145 nz = nz + maskW(i,j,k,bi,bj)
146 endif
147 end do
148 cih Compute absolute velocities from the barotropic-baroclinic modes.
149 do k = 1,Nr
150 if (k.le.nz) then
151 stmp = 0.
152 do nk = 1,nz
153 stmp = stmp +
154 & modesv(k,nk,nz)*tmpfldxz(i,nk,bi,bj)
155 end do
156 tmpz(k,bi,bj) = stmp
157 else
158 tmpz(k,bi,bj) = 0.
159 end if
160 end do
161 do k = 1,Nr
162 if (iobcs .eq. 3) then
163 tmpfldxz(i,k,bi,bj) = tmpz(k,bi,bj)
164 & *recip_hFacS(i,j+jp1,k,bi,bj)
165 else
166 tmpfldxz(i,k,bi,bj) = tmpz(k,bi,bj)
167 & *recip_hFacW(i,j,k,bi,bj)
168 endif
169 end do
170 enddo
171 endif
172 #endif
173 do k = 1,nr
174 do i = imin,imax
175 xx_obcsn1(i,k,bi,bj,iobcs) = tmpfldxz (i,k,bi,bj)
176 cgg & * maskxz (i,k,bi,bj)
177 enddo
178 enddo
179 enddo
180 enddo
181 endif
182
183 if ( (obcsnfirst) .or. (obcsnchanged)) then
184
185 do bj = jtlo,jthi
186 do bi = itlo,ithi
187 do k = 1,nr
188 do i = imin,imax
189 xx_obcsn0(i,k,bi,bj,iobcs) = xx_obcsn1(i,k,bi,bj,iobcs)
190 tmpfldxz (i,k,bi,bj) = 0. _d 0
191 enddo
192 enddo
193 enddo
194 enddo
195
196 call active_read_xz( fnameobcsn, tmpfldxz,
197 & (obcsncount1-1)*nobcs+iobcs,
198 & doglobalread, ladinit, optimcycle,
199 & mythid, xx_obcsn_dummy )
200
201 do bj = jtlo,jthi
202 do bi = itlo,ithi
203 #ifdef ALLOW_OBCS_CONTROL_MODES
204 if (iobcs .gt. 2) then
205 do i = imin,imax
206 j = OB_Jn(i,bi,bj)
207 IF ( j.EQ.OB_indexNone ) j = 1
208 cih Determine number of open vertical layers.
209 nz = 0
210 do k = 1,Nr
211 if (iobcs .eq. 3) then
212 nz = nz + maskS(i,j+jp1,k,bi,bj)
213 else
214 nz = nz + maskW(i,j,k,bi,bj)
215 endif
216 end do
217 cih Compute absolute velocities from the barotropic-baroclinic modes.
218 do k = 1,Nr
219 if (k.le.nz) then
220 stmp = 0.
221 do nk = 1,nz
222 stmp = stmp +
223 & modesv(k,nk,nz)*tmpfldxz(i,nk,bi,bj)
224 end do
225 tmpz(k,bi,bj) = stmp
226 else
227 tmpz(k,bi,bj) = 0.
228 end if
229 end do
230 do k = 1,Nr
231 if (iobcs .eq. 3) then
232 tmpfldxz(i,k,bi,bj) = tmpz(k,bi,bj)
233 & *recip_hFacS(i,j+jp1,k,bi,bj)
234 else
235 tmpfldxz(i,k,bi,bj) = tmpz(k,bi,bj)
236 & *recip_hFacW(i,j,k,bi,bj)
237 endif
238 end do
239 enddo
240 endif
241 #endif
242 do k = 1,nr
243 do i = imin,imax
244 xx_obcsn1 (i,k,bi,bj,iobcs) = tmpfldxz (i,k,bi,bj)
245 cgg & * maskxz (i,k,bi,bj)
246 enddo
247 enddo
248 enddo
249 enddo
250
251 endif
252
253 c-- Add control to model variable.
254 do bj = jtlo,jthi
255 do bi = itlo,ithi
256 c-- Calculate mask for tracer cells (0 => land, 1 => water).
257 do k = 1,nr
258 do i = 1,snx
259 j = OB_Jn(I,bi,bj)
260 IF ( j.EQ.OB_indexNone ) j = 1
261 if (iobcs .EQ. 1) then
262 OBNt(i,k,bi,bj) = OBNt (i,k,bi,bj)
263 & + obcsnfac *xx_obcsn0(i,k,bi,bj,iobcs)
264 & + (1. _d 0 - obcsnfac)*xx_obcsn1(i,k,bi,bj,iobcs)
265 OBNt(i,k,bi,bj) = OBNt(i,k,bi,bj)
266 & *maskS(i,j+jp1,k,bi,bj)
267 else if (iobcs .EQ. 2) then
268 OBNs(i,k,bi,bj) = OBNs (i,k,bi,bj)
269 & + obcsnfac *xx_obcsn0(i,k,bi,bj,iobcs)
270 & + (1. _d 0 - obcsnfac)*xx_obcsn1(i,k,bi,bj,iobcs)
271 OBNs(i,k,bi,bj) = OBNs(i,k,bi,bj)
272 & *maskS(i,j+jp1,k,bi,bj)
273 else if (iobcs .EQ. 4) then
274 OBNu(i,k,bi,bj) = OBNu (i,k,bi,bj)
275 & + obcsnfac *xx_obcsn0(i,k,bi,bj,iobcs)
276 & + (1. _d 0 - obcsnfac)*xx_obcsn1(i,k,bi,bj,iobcs)
277 OBNu(i,k,bi,bj) = OBNu(i,k,bi,bj)
278 & *maskW(i,j,k,bi,bj)
279 else if (iobcs .EQ. 3) then
280 OBNv(i,k,bi,bj) = OBNv (i,k,bi,bj)
281 & + obcsnfac *xx_obcsn0(i,k,bi,bj,iobcs)
282 & + (1. _d 0 - obcsnfac)*xx_obcsn1(i,k,bi,bj,iobcs)
283 OBNv(i,k,bi,bj) = OBNv(i,k,bi,bj)
284 & *maskS(i,j+jp1,k,bi,bj)
285 endif
286 enddo
287 enddo
288 enddo
289 enddo
290
291 C-- End over iobcs loop
292 enddo
293
294 #endif /* ALLOW_OBCSN_CONTROL */
295
296 return
297 end

  ViewVC Help
Powered by ViewVC 1.1.22