1 |
C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_getobcsw.F,v 1.10 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_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 |
_RL tmpfldyz (1-oly:sny+oly,nr,nsx,nsy) |
73 |
|
74 |
logical doglobalread |
75 |
logical ladinit |
76 |
|
77 |
character*(80) fnameobcsw |
78 |
|
79 |
cgg( Variables for splitting barotropic/baroclinic vels. |
80 |
_RL ubaro |
81 |
_RL utop |
82 |
cgg) |
83 |
|
84 |
c == external functions == |
85 |
|
86 |
integer ilnblnk |
87 |
external ilnblnk |
88 |
|
89 |
|
90 |
c == end of interface == |
91 |
|
92 |
jtlo = mybylo(mythid) |
93 |
jthi = mybyhi(mythid) |
94 |
itlo = mybxlo(mythid) |
95 |
ithi = mybxhi(mythid) |
96 |
jmin = 1-oly |
97 |
jmax = sny+oly |
98 |
imin = 1-olx |
99 |
imax = snx+olx |
100 |
ip1 = 1 |
101 |
|
102 |
cgg( Initialize variables for balancing volume flux. |
103 |
ubaro = 0.d0 |
104 |
utop = 0.d0 |
105 |
cgg) |
106 |
|
107 |
c-- Now, read the control vector. |
108 |
doglobalread = .false. |
109 |
ladinit = .false. |
110 |
|
111 |
if (optimcycle .ge. 0) then |
112 |
ilobcsw=ilnblnk( xx_obcsw_file ) |
113 |
write(fnameobcsw(1:80),'(2a,i10.10)') |
114 |
& xx_obcsw_file(1:ilobcsw), '.', optimcycle |
115 |
endif |
116 |
|
117 |
c-- Get the counters, flags, and the interpolation factor. |
118 |
call ctrl_get_gen_rec( |
119 |
I xx_obcswstartdate, xx_obcswperiod, |
120 |
O obcswfac, obcswfirst, obcswchanged, |
121 |
O obcswcount0,obcswcount1, |
122 |
I mytime, myiter, mythid ) |
123 |
|
124 |
do iobcs = 1,nobcs |
125 |
if ( obcswfirst ) then |
126 |
call active_read_yz( fnameobcsw, tmpfldyz, |
127 |
& (obcswcount0-1)*nobcs+iobcs, |
128 |
& doglobalread, ladinit, optimcycle, |
129 |
& mythid, xx_obcsw_dummy ) |
130 |
|
131 |
do bj = jtlo,jthi |
132 |
do bi = itlo,ithi |
133 |
do k = 1,nr |
134 |
do j = jmin,jmax |
135 |
xx_obcsw1(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 ( (obcswfirst) .or. (obcswchanged)) then |
144 |
|
145 |
do bj = jtlo,jthi |
146 |
do bi = itlo,ithi |
147 |
do k = 1,nr |
148 |
do j = jmin,jmax |
149 |
xx_obcsw0(j,k,bi,bj,iobcs) = xx_obcsw1(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( fnameobcsw, tmpfldyz, |
157 |
& (obcswcount1-1)*nobcs+iobcs, |
158 |
& doglobalread, ladinit, optimcycle, |
159 |
& mythid, xx_obcsw_dummy ) |
160 |
|
161 |
do bj = jtlo,jthi |
162 |
do bi = itlo,ithi |
163 |
do k = 1,nr |
164 |
do j = jmin,jmax |
165 |
xx_obcsw1 (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_Iw(j,bi,bj) |
180 |
if (iobcs .EQ. 1) then |
181 |
OBWt(j,k,bi,bj) = OBWt (j,k,bi,bj) |
182 |
& + obcswfac *xx_obcsw0(j,k,bi,bj,iobcs) |
183 |
& + (1. _d 0 - obcswfac)*xx_obcsw1(j,k,bi,bj,iobcs) |
184 |
OBWt(j,k,bi,bj) = OBWt(j,k,bi,bj) |
185 |
& *maskW(i+ip1,j,k,bi,bj) |
186 |
else if (iobcs .EQ. 2) then |
187 |
OBWs(j,k,bi,bj) = OBWs (j,k,bi,bj) |
188 |
& + obcswfac *xx_obcsw0(j,k,bi,bj,iobcs) |
189 |
& + (1. _d 0 - obcswfac)*xx_obcsw1(j,k,bi,bj,iobcs) |
190 |
OBWs(j,k,bi,bj) = OBWs(j,k,bi,bj) |
191 |
& *maskW(i+ip1,j,k,bi,bj) |
192 |
else if (iobcs .EQ. 3) then |
193 |
OBWu(j,k,bi,bj) = OBWu (j,k,bi,bj) |
194 |
& + obcswfac *xx_obcsw0(j,k,bi,bj,iobcs) |
195 |
& + (1. _d 0 - obcswfac)*xx_obcsw1(j,k,bi,bj,iobcs) |
196 |
OBWu(j,k,bi,bj) = OBWu(j,k,bi,bj) |
197 |
& *maskW(i+ip1,j,k,bi,bj) |
198 |
else if (iobcs .EQ. 4) then |
199 |
OBWv(j,k,bi,bj) = OBWv (j,k,bi,bj) |
200 |
& + obcswfac *xx_obcsw0(j,k,bi,bj,iobcs) |
201 |
& + (1. _d 0 - obcswfac)*xx_obcsw1(j,k,bi,bj,iobcs) |
202 |
OBWv(j,k,bi,bj) = OBWv(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_OBCSW_CONTROL undefined */ |
214 |
|
215 |
c == routine arguments == |
216 |
|
217 |
_RL mytime |
218 |
integer myiter |
219 |
integer mythid |
220 |
|
221 |
c-- CPP flag ALLOW_OBCSW_CONTROL undefined. |
222 |
|
223 |
#endif /* ALLOW_OBCSW_CONTROL */ |
224 |
|
225 |
end |
226 |
|