/[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.13 - (show annotations) (download)
Tue May 24 20:48:28 2011 UTC (12 years, 11 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint62z, checkpoint62y, checkpoint63g, checkpoint63, checkpoint63p, checkpoint63q, checkpoint63l, checkpoint63m, checkpoint63n, checkpoint63o, checkpoint63h, checkpoint63i, checkpoint63j, checkpoint63k, checkpoint63d, checkpoint63e, checkpoint63f, checkpoint63a, checkpoint63b, checkpoint63c
Changes since 1.12: +12 -26 lines
split "OBCS.h" into 4 separated header files (OBCS_PARAMS,GRID,FIELDS,SEAICE)

1 C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_getobcsn.F,v 1.12 2011/04/20 19:14:07 mmazloff Exp $
2 C $Name: $
3
4 #include "CTRL_CPPOPTIONS.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
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 cih Determine number of open vertical layers.
139 nz = 0
140 do k = 1,Nr
141 if (iobcs .eq. 3) then
142 nz = nz + maskS(i,j+jp1,k,bi,bj)
143 else
144 nz = nz + maskW(i,j,k,bi,bj)
145 endif
146 end do
147 cih Compute absolute velocities from the barotropic-baroclinic modes.
148 do k = 1,Nr
149 if (k.le.nz) then
150 stmp = 0.
151 do nk = 1,nz
152 stmp = stmp +
153 & modesv(k,nk,nz)*tmpfldxz(i,nk,bi,bj)
154 end do
155 tmpz(k,bi,bj) = stmp
156 else
157 tmpz(k,bi,bj) = 0.
158 end if
159 end do
160 do k = 1,Nr
161 if (iobcs .eq. 3) then
162 tmpfldxz(i,k,bi,bj) = tmpz(k,bi,bj)
163 & *recip_hFacS(i,j+jp1,k,bi,bj)
164 else
165 tmpfldxz(i,k,bi,bj) = tmpz(k,bi,bj)
166 & *recip_hFacW(i,j,k,bi,bj)
167 endif
168 end do
169 enddo
170 endif
171 #endif
172 do k = 1,nr
173 do i = imin,imax
174 xx_obcsn1(i,k,bi,bj,iobcs) = tmpfldxz (i,k,bi,bj)
175 cgg & * maskxz (i,k,bi,bj)
176 enddo
177 enddo
178 enddo
179 enddo
180 endif
181
182 if ( (obcsnfirst) .or. (obcsnchanged)) then
183
184 do bj = jtlo,jthi
185 do bi = itlo,ithi
186 do k = 1,nr
187 do i = imin,imax
188 xx_obcsn0(i,k,bi,bj,iobcs) = xx_obcsn1(i,k,bi,bj,iobcs)
189 tmpfldxz (i,k,bi,bj) = 0. _d 0
190 enddo
191 enddo
192 enddo
193 enddo
194
195 call active_read_xz( fnameobcsn, tmpfldxz,
196 & (obcsncount1-1)*nobcs+iobcs,
197 & doglobalread, ladinit, optimcycle,
198 & mythid, xx_obcsn_dummy )
199
200 do bj = jtlo,jthi
201 do bi = itlo,ithi
202 #ifdef ALLOW_OBCS_CONTROL_MODES
203 if (iobcs .gt. 2) then
204 do i = imin,imax
205 j = OB_Jn(i,bi,bj)
206 cih Determine number of open vertical layers.
207 nz = 0
208 do k = 1,Nr
209 if (iobcs .eq. 3) then
210 nz = nz + maskS(i,j+jp1,k,bi,bj)
211 else
212 nz = nz + maskW(i,j,k,bi,bj)
213 endif
214 end do
215 cih Compute absolute velocities from the barotropic-baroclinic modes.
216 do k = 1,Nr
217 if (k.le.nz) then
218 stmp = 0.
219 do nk = 1,nz
220 stmp = stmp +
221 & modesv(k,nk,nz)*tmpfldxz(i,nk,bi,bj)
222 end do
223 tmpz(k,bi,bj) = stmp
224 else
225 tmpz(k,bi,bj) = 0.
226 end if
227 end do
228 do k = 1,Nr
229 if (iobcs .eq. 3) then
230 tmpfldxz(i,k,bi,bj) = tmpz(k,bi,bj)
231 & *recip_hFacS(i,j+jp1,k,bi,bj)
232 else
233 tmpfldxz(i,k,bi,bj) = tmpz(k,bi,bj)
234 & *recip_hFacW(i,j,k,bi,bj)
235 endif
236 end do
237 enddo
238 endif
239 #endif
240 do k = 1,nr
241 do i = imin,imax
242 xx_obcsn1 (i,k,bi,bj,iobcs) = tmpfldxz (i,k,bi,bj)
243 cgg & * maskxz (i,k,bi,bj)
244 enddo
245 enddo
246 enddo
247 enddo
248
249 endif
250
251 c-- Add control to model variable.
252 do bj = jtlo,jthi
253 do bi = itlo,ithi
254 c-- Calculate mask for tracer cells (0 => land, 1 => water).
255 do k = 1,nr
256 do i = 1,snx
257 j = OB_Jn(I,bi,bj)
258 if (iobcs .EQ. 1) then
259 OBNt(i,k,bi,bj) = OBNt (i,k,bi,bj)
260 & + obcsnfac *xx_obcsn0(i,k,bi,bj,iobcs)
261 & + (1. _d 0 - obcsnfac)*xx_obcsn1(i,k,bi,bj,iobcs)
262 OBNt(i,k,bi,bj) = OBNt(i,k,bi,bj)
263 & *maskS(i,j+jp1,k,bi,bj)
264 else if (iobcs .EQ. 2) then
265 OBNs(i,k,bi,bj) = OBNs (i,k,bi,bj)
266 & + obcsnfac *xx_obcsn0(i,k,bi,bj,iobcs)
267 & + (1. _d 0 - obcsnfac)*xx_obcsn1(i,k,bi,bj,iobcs)
268 OBNs(i,k,bi,bj) = OBNs(i,k,bi,bj)
269 & *maskS(i,j+jp1,k,bi,bj)
270 else if (iobcs .EQ. 4) then
271 OBNu(i,k,bi,bj) = OBNu (i,k,bi,bj)
272 & + obcsnfac *xx_obcsn0(i,k,bi,bj,iobcs)
273 & + (1. _d 0 - obcsnfac)*xx_obcsn1(i,k,bi,bj,iobcs)
274 OBNu(i,k,bi,bj) = OBNu(i,k,bi,bj)
275 & *maskW(i,j,k,bi,bj)
276 else if (iobcs .EQ. 3) then
277 OBNv(i,k,bi,bj) = OBNv (i,k,bi,bj)
278 & + obcsnfac *xx_obcsn0(i,k,bi,bj,iobcs)
279 & + (1. _d 0 - obcsnfac)*xx_obcsn1(i,k,bi,bj,iobcs)
280 OBNv(i,k,bi,bj) = OBNv(i,k,bi,bj)
281 & *maskS(i,j+jp1,k,bi,bj)
282 endif
283 enddo
284 enddo
285 enddo
286 enddo
287
288 C-- End over iobcs loop
289 enddo
290
291 #endif /* ALLOW_OBCSN_CONTROL */
292
293 return
294 end

  ViewVC Help
Powered by ViewVC 1.1.22