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

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

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


Revision 1.16 - (hide annotations) (download)
Thu Oct 9 00:49:26 2014 UTC (9 years, 7 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.15: +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 gforget 1.16 C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_getobcse.F,v 1.15 2012/09/18 20:21:23 jmc Exp $
2 jmc 1.7 C $Name: $
3 heimbach 1.2
4 jmc 1.14 #include "CTRL_OPTIONS.h"
5 heimbach 1.2 #ifdef ALLOW_OBCS
6     # include "OBCS_OPTIONS.h"
7     #endif
8    
9     subroutine ctrl_getobcse(
10     I mytime,
11     I myiter,
12     I mythid
13     & )
14    
15     c ==================================================================
16     c SUBROUTINE ctrl_getobcse
17     c ==================================================================
18     c
19     c o Get eastern 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 ==================================================================
25     c SUBROUTINE ctrl_getobcse
26     c ==================================================================
27    
28     implicit none
29    
30 jmc 1.12 c == global variables ==
31 heimbach 1.2 #ifdef ALLOW_OBCSE_CONTROL
32     #include "EEPARAMS.h"
33     #include "SIZE.h"
34     #include "PARAMS.h"
35     #include "GRID.h"
36 jmc 1.12 c#include "OBCS_PARAMS.h"
37     #include "OBCS_GRID.h"
38     #include "OBCS_FIELDS.h"
39 heimbach 1.13 #include "CTRL_SIZE.h"
40 heimbach 1.2 #include "ctrl.h"
41     #include "ctrl_dummy.h"
42     #include "optim.h"
43 gforget 1.16 #include "CTRL_OBCS.h"
44 jmc 1.12 #endif /* ALLOW_OBCSE_CONTROL */
45 heimbach 1.2
46     c == routine arguments ==
47     _RL mytime
48     integer myiter
49     integer mythid
50    
51 jmc 1.12 #ifdef ALLOW_OBCSE_CONTROL
52 heimbach 1.2 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 ilobcse
61     integer iobcs
62    
63     _RL dummy
64     _RL obcsefac
65     logical obcsefirst
66     logical obcsechanged
67     integer obcsecount0
68     integer obcsecount1
69     integer ip1
70    
71     cgg _RL maskyz (1-oly:sny+oly,nr,nsx,nsy)
72 mlosch 1.9 _RL tmpfldyz (1-oly:sny+oly,nr,nsx,nsy)
73 heimbach 1.2
74     logical doglobalread
75     logical ladinit
76    
77     character*(80) fnameobcse
78    
79 mmazloff 1.11 #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 = 0
102    
103     c-- Now, read the control vector.
104     doglobalread = .false.
105     ladinit = .false.
106    
107     if (optimcycle .ge. 0) then
108 mlosch 1.10 ilobcse=ilnblnk( xx_obcse_file )
109     write(fnameobcse(1:80),'(2a,i10.10)')
110     & xx_obcse_file(1:ilobcse), '.', 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_obcsestartdate, xx_obcseperiod,
116     O obcsefac, obcsefirst, obcsechanged,
117     O obcsecount0,obcsecount1,
118     I mytime, myiter, mythid )
119    
120     do iobcs = 1,nobcs
121    
122 mlosch 1.10 if ( obcsefirst ) then
123     call active_read_yz( fnameobcse, tmpfldyz,
124     & (obcsecount0-1)*nobcs+iobcs,
125     & doglobalread, ladinit, optimcycle,
126     & mythid, xx_obcse_dummy )
127 heimbach 1.2
128 mlosch 1.10 do bj = jtlo,jthi
129     do bi = itlo,ithi
130 mmazloff 1.11 #ifdef ALLOW_OBCS_CONTROL_MODES
131     if (iobcs .gt. 2) then
132     do j = jmin,jmax
133     i = OB_Ie(j,bi,bj)
134 jmc 1.15 IF ( i.EQ.OB_indexNone ) i = 1
135 mmazloff 1.11 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 mlosch 1.10 do k = 1,nr
170     do j = jmin,jmax
171     xx_obcse1(j,k,bi,bj,iobcs) = tmpfldyz (j,k,bi,bj)
172     cgg & * maskyz (j,k,bi,bj)
173 mlosch 1.8 enddo
174 heimbach 1.2 enddo
175 mlosch 1.10 enddo
176     enddo
177     endif
178 jmc 1.12
179 mlosch 1.10 if ( (obcsefirst) .or. (obcsechanged)) then
180 jmc 1.12
181 mlosch 1.10 do bj = jtlo,jthi
182     do bi = itlo,ithi
183     do j = jmin,jmax
184     do k = 1,nr
185     xx_obcse0(j,k,bi,bj,iobcs) = xx_obcse1(j,k,bi,bj,iobcs)
186     tmpfldyz (j,k,bi,bj) = 0. _d 0
187     enddo
188 heimbach 1.2 enddo
189 mlosch 1.10 enddo
190     enddo
191 jmc 1.12
192 mlosch 1.10 call active_read_yz( fnameobcse, tmpfldyz,
193     & (obcsecount1-1)*nobcs+iobcs,
194     & doglobalread, ladinit, optimcycle,
195     & mythid, xx_obcse_dummy )
196 heimbach 1.2
197     do bj = jtlo,jthi
198 mlosch 1.9 do bi = itlo,ithi
199 mmazloff 1.11 #ifdef ALLOW_OBCS_CONTROL_MODES
200     if (iobcs .gt. 2) then
201     do j = jmin,jmax
202     i = OB_Ie(j,bi,bj)
203 jmc 1.15 IF ( i.EQ.OB_indexNone ) i = 1
204 mmazloff 1.11 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     endif
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 mlosch 1.9 do k = 1,nr
239 mlosch 1.10 do j = jmin,jmax
240     xx_obcse1 (j,k,bi,bj,iobcs) = tmpfldyz (j,k,bi,bj)
241     cgg & * maskyz (j,k,bi,bj)
242 heimbach 1.2 enddo
243 mlosch 1.9 enddo
244     enddo
245 heimbach 1.2 enddo
246 mlosch 1.10 endif
247 jmc 1.12
248 mlosch 1.10 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_Ie(j,bi,bj)
255 jmc 1.15 IF ( i.EQ.OB_indexNone ) i = 1
256 mlosch 1.10 if (iobcs .EQ. 1) then
257     OBEt(j,k,bi,bj) = OBEt(j,k,bi,bj)
258     & + obcsefac *xx_obcse0(j,k,bi,bj,iobcs)
259     & + (1. _d 0 - obcsefac)*xx_obcse1(j,k,bi,bj,iobcs)
260     OBEt(j,k,bi,bj) = OBEt(j,k,bi,bj)
261     & *maskW(i+ip1,j,k,bi,bj)
262     else if (iobcs .EQ. 2) then
263     OBEs(j,k,bi,bj) = OBEs(j,k,bi,bj)
264     & + obcsefac *xx_obcse0(j,k,bi,bj,iobcs)
265     & + (1. _d 0 - obcsefac)*xx_obcse1(j,k,bi,bj,iobcs)
266     OBEs(j,k,bi,bj) = OBEs(j,k,bi,bj)
267     & *maskW(i+ip1,j,k,bi,bj)
268     else if (iobcs .EQ. 3) then
269     OBEu(j,k,bi,bj) = OBEu(j,k,bi,bj)
270     & + obcsefac *xx_obcse0(j,k,bi,bj,iobcs)
271     & + (1. _d 0 - obcsefac)*xx_obcse1(j,k,bi,bj,iobcs)
272     OBEu(j,k,bi,bj) = OBEu(j,k,bi,bj)
273     & *maskW(i+ip1,j,k,bi,bj)
274     else if (iobcs .EQ. 4) then
275     OBEv(j,k,bi,bj) = OBEv(j,k,bi,bj)
276     & + obcsefac *xx_obcse0(j,k,bi,bj,iobcs)
277     & + (1. _d 0 - obcsefac)*xx_obcse1(j,k,bi,bj,iobcs)
278     OBEv(j,k,bi,bj) = OBEv(j,k,bi,bj)
279     & *maskS(i,j,k,bi,bj)
280     endif
281     enddo
282     enddo
283     enddo
284     enddo
285 jmc 1.12
286 heimbach 1.2 C-- End over iobcs loop
287     enddo
288    
289     #endif /* ALLOW_OBCSE_CONTROL */
290    
291 jmc 1.12 return
292 heimbach 1.2 end

  ViewVC Help
Powered by ViewVC 1.1.22