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

Annotation of /MITgcm_contrib/SOSE/BoxAdj/code_ad/ctrl_getobcss.F

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


Revision 1.2 - (hide annotations) (download)
Wed Jan 19 12:39:13 2011 UTC (14 years, 11 months ago) by mmazloff
Branch: MAIN
Changes since 1.1: +18 -5 lines
account for partial cells

1 mmazloff 1.1
2     #include "CTRL_CPPOPTIONS.h"
3     #ifdef ALLOW_OBCS
4     # include "OBCS_OPTIONS.h"
5     #endif
6    
7    
8     subroutine ctrl_getobcss(
9     I mytime,
10     I myiter,
11     I mythid
12     & )
13    
14     c ==================================================================
15     c SUBROUTINE ctrl_getobcss
16     c ==================================================================
17     c
18     c o Get southern obc of the control vector and add it
19     c to dyn. fields
20     c
21     c started: heimbach@mit.edu, 29-Aug-2001
22     c
23     c new flags: gebbie@mit.edu, 25 Jan 2003.
24     c
25     c ==================================================================
26     c SUBROUTINE ctrl_getobcss
27     c ==================================================================
28    
29     implicit none
30    
31     #ifdef ALLOW_OBCSS_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 ilobcss
60     integer iobcs
61    
62     _RL dummy
63     _RL obcssfac
64     logical obcssfirst
65     logical obcsschanged
66     integer obcsscount0
67     integer obcsscount1
68     integer jp1
69     cih
70     Integer nk,nz
71     _RL tmpz (nr,nsx,nsy)
72     _RL stmp
73     character*(80) fnamein
74     cih
75     logical doglobalread
76     logical ladinit
77    
78     character*(80) fnameobcss
79    
80    
81     c == external functions ==
82    
83     integer ilnblnk
84     external ilnblnk
85    
86    
87     c == end of interface ==
88    
89     jtlo = mybylo(mythid)
90     jthi = mybyhi(mythid)
91     itlo = mybxlo(mythid)
92     ithi = mybxhi(mythid)
93     jmin = 1-oly
94     jmax = sny+oly
95     imin = 1-olx
96     imax = snx+olx
97     jp1 = 1
98    
99    
100     c-- Now, read the control vector.
101     doglobalread = .false.
102     ladinit = .false.
103    
104     if (optimcycle .ge. 0) then
105     ilobcss=ilnblnk( xx_obcss_file )
106     write(fnameobcss(1:80),'(2a,i10.10)')
107     & xx_obcss_file(1:ilobcss), '.', optimcycle
108     endif
109    
110     c-- Get the counters, flags, and the interpolation factor.
111     call ctrl_get_gen_rec(
112     I xx_obcssstartdate, xx_obcssperiod,
113     O obcssfac, obcssfirst, obcsschanged,
114     O obcsscount0,obcsscount1,
115     I mytime, myiter, mythid )
116     cih
117     do iobcs = 1,nobcs
118     if ( obcssfirst ) then
119     call active_read_xz( fnameobcss, tmpfldxz,
120     & (obcsscount0-1)*nobcs+iobcs,
121     & doglobalread, ladinit, optimcycle,
122     & mythid, xx_obcss_dummy )
123 mmazloff 1.2 CNMA
124     #ifdef CTRL_OBCS_MODES
125     cih
126 mmazloff 1.1 if ( optimcycle .ge. 0) then
127     cih If normal velocity
128     if (iobcs .eq. 3) then
129     cih Begin loop over x-points.
130     do bj = jtlo,jthi
131     do bi = itlo, ithi
132     do i = imin,imax
133     cih If open boundary.
134     if ( OB_Js(i,bi,bj) .ne. 0. ) then
135     j = OB_Js(i,bi,bj)
136     cih Determine number of open vertical layers.
137     nz = 0
138     do k = 1,Nr
139     nz = nz + maskS(i,j+jp1,k,bi,bj)
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)*tmpfldxz(i,nk,bi,bj)
148     end do
149     tmpz(k,bi,bj) = stmp
150     else
151     tmpz(k,bi,bj) = 0.
152     end if
153     end do
154     do k = 1,Nr
155 mmazloff 1.2 tmpfldxz(i,k,bi,bj) = tmpz(k,bi,bj)
156     & *recip_hFacS(i,j+jp1,k,bi,bj)
157 mmazloff 1.1 end do
158     cih End if open boundary.
159     end if
160     cih End loop over x-points
161     end do
162     end do
163     end do
164     cih End if iobcs = 3.
165     end if
166     cih
167     cih
168     cih If tangential velocity
169     if (iobcs .eq. 4) then
170     cih Begin loop over x-points.
171     do bj = jtlo,jthi
172     do bi = itlo, ithi
173     do i = imin,imax
174     cih If open boundary.
175     if ( OB_Js(i,bi,bj) .ne. 0. ) then
176     j = OB_Js(i,bi,bj)
177     cih Determine number of open vertical layers.
178     nz = 0
179     do k = 1,Nr
180     nz = nz + maskW(i,j,k,bi,bj)
181     end do
182     cih Compute absolute velocities from the barotropic-baroclinic modes.
183     do k = 1,Nr
184     if (k.le.nz) then
185     stmp = 0.
186     do nk = 1,nz
187     stmp = stmp +
188     & modesv(k,nk,nz)*tmpfldxz(i,nk,bi,bj)
189     end do
190     tmpz(k,bi,bj) = stmp
191     else
192     tmpz(k,bi,bj) = 0.
193     end if
194     end do
195     do k = 1,Nr
196 mmazloff 1.2 tmpfldxz(i,k,bi,bj) = tmpz(k,bi,bj)
197     & *recip_hFacW(i,j,k,bi,bj)
198 mmazloff 1.1 end do
199     cih End if open boundary.
200     end if
201     cih End loop over x-points
202     end do
203     end do
204     end do
205     cih End if iobcs = 4.
206     end if
207     cih End if optimcycle > 0 .
208     end if
209 mmazloff 1.2 CNMA
210     #endif
211 mmazloff 1.1 cih
212     do bj = jtlo,jthi
213     do bi = itlo,ithi
214     do k = 1,nr
215     do i = imin,imax
216     xx_obcss1(i,k,bi,bj,iobcs) = tmpfldxz (i,k,bi,bj)
217     enddo
218     enddo
219     enddo
220     enddo
221     cih End if obcssfirst
222     endif
223     cih
224     if ( (obcssfirst) .or. (obcsschanged) ) then
225     cih
226     do bj = jtlo,jthi
227     do bi = itlo,ithi
228     do k = 1,nr
229     do i = imin,imax
230     tmpfldxz(i,k,bi,bj) = xx_obcss1(i,k,bi,bj,iobcs)
231     enddo
232     enddo
233     enddo
234     enddo
235     cih
236     call exf_swapffields_xz( tmpfldxz2, tmpfldxz, mythid)
237     cih
238     do bj = jtlo,jthi
239     do bi = itlo,ithi
240     do k = 1,nr
241     do i = imin,imax
242     xx_obcss0(i,k,bi,bj,iobcs) = tmpfldxz2(i,k,bi,bj)
243     enddo
244     enddo
245     enddo
246     enddo
247     cih
248     call active_read_xz( fnameobcss, tmpfldxz,
249     & (obcsscount1-1)*nobcs+iobcs,
250     & doglobalread, ladinit, optimcycle,
251     & mythid, xx_obcss_dummy )
252 mmazloff 1.2 CNMA
253     #ifdef CTRL_OBCS_MODES
254 mmazloff 1.1 cih
255     if ( optimcycle .ge. 0) then
256     cih If normal velocity
257     if (iobcs .eq. 3) then
258     cih Begin loop over x-points.
259     do bj = jtlo,jthi
260     do bi = itlo, ithi
261     do i = imin,imax
262     cih If open boundary.
263     if ( OB_Js(i,bi,bj) .ne. 0. ) then
264     j = OB_Js(i,bi,bj)
265     cih Determine number of open vertical layers.
266     nz = 0
267     do k = 1,Nr
268     nz = nz + maskS(i,j+jp1,k,bi,bj)
269     end do
270     cih Compute absolute velocities from the barotropic-baroclinic modes.
271     do k = 1,Nr
272     if (k.le.nz) then
273     stmp = 0.
274     do nk = 1,nz
275     stmp = stmp +
276     & modesv(k,nk,nz)*tmpfldxz(i,nk,bi,bj)
277     end do
278     tmpz(k,bi,bj) = stmp
279     else
280     tmpz(k,bi,bj) = 0.
281     end if
282     end do
283     do k = 1,Nr
284 mmazloff 1.2 tmpfldxz(i,k,bi,bj) = tmpz(k,bi,bj)
285     & *recip_hFacS(i,j+jp1,k,bi,bj)
286 mmazloff 1.1 end do
287     cih End if open boundary.
288     end if
289     cih End loop over x-points
290     end do
291     end do
292     end do
293     cih End if iobcs = 3.
294     end if
295     cih
296     cih
297     cih If tangential velocity
298     if (iobcs .eq. 4) then
299     cih Begin loop over x-points.
300     do bj = jtlo,jthi
301     do bi = itlo, ithi
302     do i = imin,imax
303     cih If open boundary.
304     if ( OB_Js(i,bi,bj) .ne. 0. ) then
305     j = OB_Js(i,bi,bj)
306     cih Determine number of open vertical layers.
307     nz = 0
308     do k = 1,Nr
309     nz = nz + maskW(i,j,k,bi,bj)
310     end do
311     cih Compute absolute velocities from the barotropic-baroclinic modes.
312     do k = 1,Nr
313     if (k.le.nz) then
314     stmp = 0.
315     do nk = 1,nz
316     stmp = stmp +
317     & modesv(k,nk,nz)*tmpfldxz(i,nk,bi,bj)
318     end do
319     tmpz(k,bi,bj) = stmp
320     else
321     tmpz(k,bi,bj) = 0.
322     end if
323     end do
324     do k = 1,Nr
325 mmazloff 1.2 c tmpfldxz(i,k,bi,bj) = tmpz(k,bi,bj)
326     tmpfldxz(i,k,bi,bj) = tmpz(k,bi,bj)
327     & *recip_hFacW(i,j,k,bi,bj)
328 mmazloff 1.1 end do
329     cih End if open boundary.
330     end if
331     cih End loop over x-points
332     end do
333     end do
334     end do
335     cih End if iobcs = 4.
336     end if
337     cih End if optimcycle > 0 .
338     end if
339 mmazloff 1.2 CNMA
340     #endif
341 mmazloff 1.1 cih
342     do bj = jtlo,jthi
343     do bi = itlo,ithi
344     do k = 1,nr
345     do i = imin,imax
346     xx_obcss1 (i,k,bi,bj,iobcs) = tmpfldxz (i,k,bi,bj)
347     enddo
348     enddo
349     enddo
350     enddo
351     cih End if obcssfirst .or. obcschanged
352     endif
353    
354     c-- Add control to model variable.
355     do bj = jtlo,jthi
356     do bi = itlo,ithi
357     c-- Calculate mask for tracer cells (0 => land, 1 => water).
358     do k = 1,nr
359     do i = 1,snx
360     j = OB_Js(i,bi,bj)
361     if (iobcs .EQ. 1) then
362     OBSt(i,k,bi,bj) = OBSt (i,k,bi,bj)
363     & + obcssfac *xx_obcss0(i,k,bi,bj,iobcs)
364     & + (1. _d 0 - obcssfac)*xx_obcss1(i,k,bi,bj,iobcs)
365     OBSt(i,k,bi,bj) = OBSt(i,k,bi,bj)
366     & *maskS(i,j+jp1,k,bi,bj)
367     else if (iobcs .EQ. 2) then
368     OBSs(i,k,bi,bj) = OBSs (i,k,bi,bj)
369     & + obcssfac *xx_obcss0(i,k,bi,bj,iobcs)
370     & + (1. _d 0 - obcssfac)*xx_obcss1(i,k,bi,bj,iobcs)
371     OBSs(i,k,bi,bj) = OBSs(i,k,bi,bj)
372     & *maskS(i,j+jp1,k,bi,bj)
373     else if (iobcs .EQ. 4) then
374     OBSu(i,k,bi,bj) = OBSu (i,k,bi,bj)
375     & + obcssfac *xx_obcss0(i,k,bi,bj,iobcs)
376     & + (1. _d 0 - obcssfac)*xx_obcss1(i,k,bi,bj,iobcs)
377     OBSu(i,k,bi,bj) = OBSu(i,k,bi,bj)
378     & *maskW(i,j,k,bi,bj)
379     else if (iobcs .EQ. 3) then
380     OBSv(i,k,bi,bj) = OBSv (i,k,bi,bj)
381     & + obcssfac *xx_obcss0(i,k,bi,bj,iobcs)
382     & + (1. _d 0 - obcssfac)*xx_obcss1(i,k,bi,bj,iobcs)
383     OBSv(i,k,bi,bj) = OBSv(i,k,bi,bj)
384     & *maskS(i,j+jp1,k,bi,bj)
385     endif
386     enddo
387     enddo
388     enddo
389     enddo
390    
391     C-- End over iobcs loop
392     enddo
393    
394     #else /* ALLOW_OBCSS_CONTROL undefined */
395    
396     c == routine arguments ==
397    
398     _RL mytime
399     integer myiter
400     integer mythid
401    
402     c-- CPP flag ALLOW_OBCSS_CONTROL undefined.
403    
404     #endif /* ALLOW_OBCSS_CONTROL */
405    
406     end
407    

  ViewVC Help
Powered by ViewVC 1.1.22