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

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

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


Revision 1.12 - (hide annotations) (download)
Wed Apr 20 19:14:07 2011 UTC (14 years, 2 months ago) by mmazloff
Branch: MAIN
CVS Tags: checkpoint62w, checkpoint62x
Changes since 1.11: +82 -10 lines
Ability to decompose obcs controls into modes using ALLOW_OBCS_CONTROL_MODES

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

  ViewVC Help
Powered by ViewVC 1.1.22