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

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

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


Revision 1.10 - (show annotations) (download)
Mon Mar 14 17:08:00 2011 UTC (13 years, 2 months ago) by mlosch
Branch: MAIN
CVS Tags: checkpoint62v, checkpoint62u
Changes since 1.9: +76 -230 lines
remove obsolete and partially broken code, step 2:
remove code within ALLOW_CTRL_OBCS_BALANCE

1 C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_getobcse.F,v 1.9 2011/03/07 09:24:10 mlosch Exp $
2 C $Name: $
3
4 #include "CTRL_CPPOPTIONS.h"
5 #ifdef ALLOW_OBCS
6 # include "OBCS_OPTIONS.h"
7 #endif
8
9
10 subroutine ctrl_getobcse(
11 I mytime,
12 I myiter,
13 I mythid
14 & )
15
16 c ==================================================================
17 c SUBROUTINE ctrl_getobcse
18 c ==================================================================
19 c
20 c o Get eastern 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 ==================================================================
26 c SUBROUTINE ctrl_getobcse
27 c ==================================================================
28
29 implicit none
30
31 #ifdef ALLOW_OBCSE_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 ilobcse
60 integer iobcs
61
62 _RL dummy
63 _RL obcsefac
64 logical obcsefirst
65 logical obcsechanged
66 integer obcsecount0
67 integer obcsecount1
68 integer ip1
69
70 cgg _RL maskyz (1-oly:sny+oly,nr,nsx,nsy)
71 _RL tmpfldyz (1-oly:sny+oly,nr,nsx,nsy)
72
73 logical doglobalread
74 logical ladinit
75
76 character*(80) fnameobcse
77
78 cgg( Variables for splitting barotropic/baroclinic vels.
79 _RL ubaro
80 _RL utop
81 cgg)
82
83 c == external functions ==
84
85 integer ilnblnk
86 external ilnblnk
87
88
89 c == end of interface ==
90
91 jtlo = mybylo(mythid)
92 jthi = mybyhi(mythid)
93 itlo = mybxlo(mythid)
94 ithi = mybxhi(mythid)
95 jmin = 1-oly
96 jmax = sny+oly
97 imin = 1-olx
98 imax = snx+olx
99 ip1 = 0
100
101 cgg( Initialize variables for balancing volume flux.
102 ubaro = 0.d0
103 utop = 0.d0
104 cgg)
105
106 c-- Now, read the control vector.
107 doglobalread = .false.
108 ladinit = .false.
109
110 if (optimcycle .ge. 0) then
111 ilobcse=ilnblnk( xx_obcse_file )
112 write(fnameobcse(1:80),'(2a,i10.10)')
113 & xx_obcse_file(1:ilobcse), '.', optimcycle
114 endif
115
116 c-- Get the counters, flags, and the interpolation factor.
117 call ctrl_get_gen_rec(
118 I xx_obcsestartdate, xx_obcseperiod,
119 O obcsefac, obcsefirst, obcsechanged,
120 O obcsecount0,obcsecount1,
121 I mytime, myiter, mythid )
122
123 do iobcs = 1,nobcs
124
125 if ( obcsefirst ) then
126 call active_read_yz( fnameobcse, tmpfldyz,
127 & (obcsecount0-1)*nobcs+iobcs,
128 & doglobalread, ladinit, optimcycle,
129 & mythid, xx_obcse_dummy )
130
131 do bj = jtlo,jthi
132 do bi = itlo,ithi
133 do k = 1,nr
134 do j = jmin,jmax
135 xx_obcse1(j,k,bi,bj,iobcs) = tmpfldyz (j,k,bi,bj)
136 cgg & * maskyz (j,k,bi,bj)
137 enddo
138 enddo
139 enddo
140 enddo
141 endif
142
143 if ( (obcsefirst) .or. (obcsechanged)) then
144
145 do bj = jtlo,jthi
146 do bi = itlo,ithi
147 do j = jmin,jmax
148 do k = 1,nr
149 xx_obcse0(j,k,bi,bj,iobcs) = xx_obcse1(j,k,bi,bj,iobcs)
150 tmpfldyz (j,k,bi,bj) = 0. _d 0
151 enddo
152 enddo
153 enddo
154 enddo
155
156 call active_read_yz( fnameobcse, tmpfldyz,
157 & (obcsecount1-1)*nobcs+iobcs,
158 & doglobalread, ladinit, optimcycle,
159 & mythid, xx_obcse_dummy )
160
161 do bj = jtlo,jthi
162 do bi = itlo,ithi
163 do k = 1,nr
164 do j = jmin,jmax
165 xx_obcse1 (j,k,bi,bj,iobcs) = tmpfldyz (j,k,bi,bj)
166 cgg & * maskyz (j,k,bi,bj)
167 enddo
168 enddo
169 enddo
170 enddo
171 endif
172
173 c-- Add control to model variable.
174 do bj = jtlo,jthi
175 do bi = itlo,ithi
176 c-- Calculate mask for tracer cells (0 => land, 1 => water).
177 do k = 1,nr
178 do j = 1,sny
179 i = OB_Ie(j,bi,bj)
180 if (iobcs .EQ. 1) then
181 OBEt(j,k,bi,bj) = OBEt(j,k,bi,bj)
182 & + obcsefac *xx_obcse0(j,k,bi,bj,iobcs)
183 & + (1. _d 0 - obcsefac)*xx_obcse1(j,k,bi,bj,iobcs)
184 OBEt(j,k,bi,bj) = OBEt(j,k,bi,bj)
185 & *maskW(i+ip1,j,k,bi,bj)
186 else if (iobcs .EQ. 2) then
187 OBEs(j,k,bi,bj) = OBEs(j,k,bi,bj)
188 & + obcsefac *xx_obcse0(j,k,bi,bj,iobcs)
189 & + (1. _d 0 - obcsefac)*xx_obcse1(j,k,bi,bj,iobcs)
190 OBEs(j,k,bi,bj) = OBEs(j,k,bi,bj)
191 & *maskW(i+ip1,j,k,bi,bj)
192 else if (iobcs .EQ. 3) then
193 OBEu(j,k,bi,bj) = OBEu(j,k,bi,bj)
194 & + obcsefac *xx_obcse0(j,k,bi,bj,iobcs)
195 & + (1. _d 0 - obcsefac)*xx_obcse1(j,k,bi,bj,iobcs)
196 OBEu(j,k,bi,bj) = OBEu(j,k,bi,bj)
197 & *maskW(i+ip1,j,k,bi,bj)
198 else if (iobcs .EQ. 4) then
199 OBEv(j,k,bi,bj) = OBEv(j,k,bi,bj)
200 & + obcsefac *xx_obcse0(j,k,bi,bj,iobcs)
201 & + (1. _d 0 - obcsefac)*xx_obcse1(j,k,bi,bj,iobcs)
202 OBEv(j,k,bi,bj) = OBEv(j,k,bi,bj)
203 & *maskS(i,j,k,bi,bj)
204 endif
205 enddo
206 enddo
207 enddo
208 enddo
209
210 C-- End over iobcs loop
211 enddo
212
213 #else /* ALLOW_OBCSE_CONTROL undefined */
214
215 c == routine arguments ==
216
217 _RL mytime
218 integer myiter
219 integer mythid
220
221 c-- CPP flag ALLOW_OBCSE_CONTROL undefined.
222
223 #endif /* ALLOW_OBCSE_CONTROL */
224
225 end
226

  ViewVC Help
Powered by ViewVC 1.1.22