1 |
|
2 |
#include "CTRL_CPPOPTIONS.h" |
3 |
|
4 |
|
5 |
subroutine ctrl_set_unpack_xz( |
6 |
& cunit, ivartype, fname, masktype, weighttype, |
7 |
& weightfld, nwetglobal, mythid) |
8 |
|
9 |
c ================================================================== |
10 |
c SUBROUTINE ctrl_set_unpack_xz |
11 |
c ================================================================== |
12 |
c |
13 |
c o Unpack the control vector such that land points are filled in. |
14 |
c |
15 |
c o Open boundary packing added : |
16 |
c gebbie@mit.edu, 18-Mar-2003 |
17 |
c |
18 |
c changed: heimbach@mit.edu 17-Jun-2003 |
19 |
c merged Armin's changes to replace write of |
20 |
c nr * globfld2d by 1 * globfld3d |
21 |
c (ad hoc fix to speed up global I/O) |
22 |
c |
23 |
c ================================================================== |
24 |
|
25 |
implicit none |
26 |
|
27 |
c == global variables == |
28 |
|
29 |
#include "EEPARAMS.h" |
30 |
#include "SIZE.h" |
31 |
#include "PARAMS.h" |
32 |
#include "GRID.h" |
33 |
|
34 |
#include "ctrl.h" |
35 |
#include "cost.h" |
36 |
|
37 |
#ifdef ALLOW_ECCO_OPTIMIZATION |
38 |
#include "optim.h" |
39 |
#endif |
40 |
|
41 |
c == routine arguments == |
42 |
|
43 |
integer cunit |
44 |
integer ivartype |
45 |
character*( 80) fname |
46 |
character* (9) masktype |
47 |
character*( 80) weighttype |
48 |
_RL weightfld( nr,nobcs ) |
49 |
integer nwetglobal(nr,nobcs) |
50 |
integer mythid |
51 |
|
52 |
c == local variables == |
53 |
|
54 |
#ifndef ALLOW_ECCO_OPTIMIZATION |
55 |
integer optimcycle |
56 |
#endif |
57 |
|
58 |
integer bi,bj |
59 |
integer ip,jp |
60 |
integer i,j,k |
61 |
integer ii |
62 |
integer il |
63 |
integer irec,iobcs,nrec_nl |
64 |
integer itlo,ithi |
65 |
integer jtlo,jthi |
66 |
integer jmin,jmax |
67 |
integer imin,imax |
68 |
|
69 |
integer cbuffindex |
70 |
|
71 |
_RL cbuff ( snx*nsx*npx*nsy*npy ) |
72 |
_RL globfldxz( snx,nsx,npx,nsy,npy,nr ) |
73 |
_RL globfld3d( snx,nsx,npx,sny,nsy,npy,nr ) |
74 |
_RL globmskxz( snx,nsx,npx,nsy,npy,nr,nobcs ) |
75 |
#ifdef CTRL_UNPACK_PRECISE |
76 |
_RL weightfldxz( snx,nsx,npx,nsy,npy,nr,nobcs ) |
77 |
#endif |
78 |
|
79 |
integer filenvartype |
80 |
integer filenvarlength |
81 |
character*(10) fileExpId |
82 |
integer fileOptimCycle |
83 |
integer filencbuffindex |
84 |
_RL fileDummy |
85 |
integer fileIg |
86 |
integer fileJg |
87 |
integer fileI |
88 |
integer fileJ |
89 |
integer filensx |
90 |
integer filensy |
91 |
integer filek |
92 |
integer filencvarindex(maxcvars) |
93 |
integer filencvarrecs(maxcvars) |
94 |
integer filencvarxmax(maxcvars) |
95 |
integer filencvarymax(maxcvars) |
96 |
integer filencvarnrmax(maxcvars) |
97 |
character*( 1) filencvargrd(maxcvars) |
98 |
cgg( |
99 |
integer igg |
100 |
_RL gg |
101 |
character*(80) weightname |
102 |
cgg) |
103 |
|
104 |
c == external == |
105 |
|
106 |
integer ilnblnk |
107 |
external ilnblnk |
108 |
|
109 |
cc == end of interface == |
110 |
|
111 |
jtlo = 1 |
112 |
jthi = nsy |
113 |
itlo = 1 |
114 |
ithi = nsx |
115 |
jmin = 1 |
116 |
jmax = sny |
117 |
imin = 1 |
118 |
imax = snx |
119 |
|
120 |
c Initialise temporary file |
121 |
do k = 1,nr |
122 |
do jp = 1,nPy |
123 |
do bj = jtlo,jthi |
124 |
do ip = 1,nPx |
125 |
do bi = itlo,ithi |
126 |
do i = imin,imax |
127 |
globfldxz(i,bi,ip,bj,jp,k) = 0. _d 0 |
128 |
do iobcs=1,nobcs |
129 |
globmskxz(i,bi,ip,bj,jp,k,iobcs) = 0. _d 0 |
130 |
enddo |
131 |
enddo |
132 |
enddo |
133 |
enddo |
134 |
enddo |
135 |
enddo |
136 |
enddo |
137 |
c Initialise temporary file |
138 |
do k = 1,nr |
139 |
do jp = 1,nPy |
140 |
do bj = jtlo,jthi |
141 |
do j = jmin,jmax |
142 |
do ip = 1,nPx |
143 |
do bi = itlo,ithi |
144 |
do i = imin,imax |
145 |
globfld3d(i,bi,ip,j,bj,jp,k) = 0. _d 0 |
146 |
enddo |
147 |
enddo |
148 |
enddo |
149 |
enddo |
150 |
enddo |
151 |
enddo |
152 |
enddo |
153 |
|
154 |
#ifndef ALLOW_ECCO_OPTIMIZATION |
155 |
optimcycle = 0 |
156 |
#endif |
157 |
|
158 |
c-- Only the master thread will do I/O. |
159 |
_BEGIN_MASTER( mythid ) |
160 |
|
161 |
do iobcs=1,nobcs |
162 |
call MDSREADFIELD_XZ_GL( |
163 |
& masktype, ctrlprec, 'RL', |
164 |
& Nr, globmskxz(1,1,1,1,1,1,iobcs), iobcs,mythid) |
165 |
#ifdef CTRL_UNPACK_PRECISE |
166 |
il=ilnblnk( weighttype) |
167 |
write(weightname(1:80),'(80a)') ' ' |
168 |
write(weightname(1:80),'(a)') weighttype(1:il) |
169 |
call MDSREADFIELD_XZ_GL( |
170 |
& weightname, ctrlprec, 'RL', |
171 |
& Nr, weightfldxz(1,1,1,1,1,1,iobcs), iobcs, mythid) |
172 |
CGG One special exception: barotropic velocity should be nondimensionalized |
173 |
cgg differently. Probably introduce new variable. |
174 |
if (iobcs .eq. 3 .or. iobcs .eq. 4) then |
175 |
k = 1 |
176 |
do jp = 1,nPy |
177 |
do bj = jtlo,jthi |
178 |
do ip = 1,nPx |
179 |
do bi = itlo,ithi |
180 |
do i = imin,imax |
181 |
weightfldxz(i,bi,ip,bj,jp,k,iobcs) = wbaro |
182 |
enddo |
183 |
enddo |
184 |
enddo |
185 |
enddo |
186 |
enddo |
187 |
endif |
188 |
#endif /* CTRL_UNPACK_PRECISE */ |
189 |
enddo |
190 |
|
191 |
nrec_nl=int(ncvarrecs(ivartype)/sny) |
192 |
do irec = 1, nrec_nl |
193 |
cgg do iobcs = 1, nobcs |
194 |
cgg And now back-calculate what iobcs should be. |
195 |
do j=1,sny |
196 |
iobcs= mod((irec-1)*sny+j-1,nobcs)+1 |
197 |
|
198 |
read(cunit) filencvarindex(ivartype) |
199 |
if (filencvarindex(ivartype) .NE. ncvarindex(ivartype)) |
200 |
& then |
201 |
print *, 'ctrl-set_unpack:xz:WARNING: wrong ncvarindex ', |
202 |
& filencvarindex(ivartype), ncvarindex(ivartype) |
203 |
STOP 'in S/R ctrl_unpack' |
204 |
endif |
205 |
read(cunit) filej |
206 |
read(cunit) filei |
207 |
do k = 1, Nr |
208 |
cbuffindex = nwetglobal(k,iobcs) |
209 |
if ( cbuffindex .gt. 0 ) then |
210 |
read(cunit) filencbuffindex |
211 |
if (filencbuffindex .NE. cbuffindex) then |
212 |
print *, 'WARNING: wrong cbuffindex ', |
213 |
& filencbuffindex, cbuffindex |
214 |
STOP 'in S/R ctrl_unpack' |
215 |
endif |
216 |
read(cunit) filek |
217 |
if (filek .NE. k) then |
218 |
print *, 'WARNING: wrong k ', |
219 |
& filek, k |
220 |
STOP 'in S/R ctrl_unpack' |
221 |
endif |
222 |
read(cunit) (cbuff(ii), ii=1,cbuffindex) |
223 |
endif |
224 |
cbuffindex = 0 |
225 |
do jp = 1,nPy |
226 |
do bj = jtlo,jthi |
227 |
do ip = 1,nPx |
228 |
do bi = itlo,ithi |
229 |
do i = imin,imax |
230 |
if ( globmskxz(i,bi,ip,bj,jp,k,iobcs) .ne. 0. ) then |
231 |
cbuffindex = cbuffindex + 1 |
232 |
globfld3d(i,bi,ip,j,bj,jp,k) = cbuff(cbuffindex) |
233 |
#ifdef ALLOW_NONDIMENSIONAL_CONTROL_IO |
234 |
globfld3d(i,bi,ip,j,bj,jp,k) = |
235 |
& globfld3d(i,bi,ip,j,bj,jp,k)/ |
236 |
# ifdef CTRL_UNPACK_PRECISE |
237 |
& sqrt(weightfldxz(i,bi,ip,bj,jp,k,iobcs)) |
238 |
# else |
239 |
& sqrt(weightfld(k,iobcs)) |
240 |
# endif |
241 |
#endif /* ALLOW_NONDIMENSIONAL_CONTROL_IO */ |
242 |
else |
243 |
globfld3d(i,bi,ip,j,bj,jp,k) = 0. _d 0 |
244 |
endif |
245 |
enddo |
246 |
enddo |
247 |
enddo |
248 |
enddo |
249 |
enddo |
250 |
c |
251 |
c -- end of k loop -- |
252 |
enddo |
253 |
c -- end of j loop -- |
254 |
enddo |
255 |
|
256 |
call MDSWRITEFIELD_3D_GL( fname, ctrlprec, 'RL', |
257 |
& Nr, globfld3d, irec, |
258 |
& optimcycle, mythid) |
259 |
|
260 |
c -- end of iobcs loop -- This loop removed. 3-28-02. |
261 |
cgg enddo |
262 |
c -- end of irec loop -- |
263 |
enddo |
264 |
|
265 |
do irec = nrec_nl*sny+1, ncvarrecs(ivartype) |
266 |
cgg do iobcs = 1, nobcs |
267 |
cgg And now back-calculate what iobcs should be. |
268 |
iobcs= mod(irec-1,nobcs)+1 |
269 |
|
270 |
read(cunit) filencvarindex(ivartype) |
271 |
if (filencvarindex(ivartype) .NE. ncvarindex(ivartype)) |
272 |
& then |
273 |
print *, 'ctrl-set_unpack:xz:WARNING: wrong ncvarindex ', |
274 |
& filencvarindex(ivartype), ncvarindex(ivartype) |
275 |
STOP 'in S/R ctrl_unpack' |
276 |
endif |
277 |
read(cunit) filej |
278 |
read(cunit) filei |
279 |
do k = 1, Nr |
280 |
cbuffindex = nwetglobal(k,iobcs) |
281 |
if ( cbuffindex .gt. 0 ) then |
282 |
read(cunit) filencbuffindex |
283 |
if (filencbuffindex .NE. cbuffindex) then |
284 |
print *, 'WARNING: wrong cbuffindex ', |
285 |
& filencbuffindex, cbuffindex |
286 |
STOP 'in S/R ctrl_unpack' |
287 |
endif |
288 |
read(cunit) filek |
289 |
if (filek .NE. k) then |
290 |
print *, 'WARNING: wrong k ', |
291 |
& filek, k |
292 |
STOP 'in S/R ctrl_unpack' |
293 |
endif |
294 |
read(cunit) (cbuff(ii), ii=1,cbuffindex) |
295 |
endif |
296 |
cbuffindex = 0 |
297 |
do jp = 1,nPy |
298 |
do bj = jtlo,jthi |
299 |
do ip = 1,nPx |
300 |
do bi = itlo,ithi |
301 |
do i = imin,imax |
302 |
if ( globmskxz(i,bi,ip,bj,jp,k,iobcs) .ne. 0. ) then |
303 |
cbuffindex = cbuffindex + 1 |
304 |
globfldxz(i,bi,ip,bj,jp,k) = cbuff(cbuffindex) |
305 |
#ifdef ALLOW_NONDIMENSIONAL_CONTROL_IO |
306 |
globfldxz(i,bi,ip,bj,jp,k) = |
307 |
& globfldxz(i,bi,ip,bj,jp,k)/ |
308 |
# ifdef CTRL_UNPACK_PRECISE |
309 |
& sqrt(weightfldxz(i,bi,ip,bj,jp,k,iobcs)) |
310 |
# else |
311 |
& sqrt(weightfld(k,iobcs)) |
312 |
# endif |
313 |
#endif /* ALLOW_NONDIMENSIONAL_CONTROL_IO */ |
314 |
else |
315 |
globfldxz(i,bi,ip,bj,jp,k) = 0. _d 0 |
316 |
endif |
317 |
enddo |
318 |
enddo |
319 |
enddo |
320 |
enddo |
321 |
enddo |
322 |
c |
323 |
c -- end of k loop -- |
324 |
enddo |
325 |
|
326 |
call MDSWRITEFIELD_XZ_GL( fname, ctrlprec, 'RL', |
327 |
& Nr, globfldxz, irec, |
328 |
& optimcycle, mythid) |
329 |
|
330 |
c -- end of iobcs loop -- This loop removed. 3-28-02. |
331 |
cgg enddo |
332 |
c -- end of irec loop -- |
333 |
enddo |
334 |
|
335 |
_END_MASTER( mythid ) |
336 |
|
337 |
return |
338 |
end |
339 |
|
340 |
|
341 |
|
342 |
|
343 |
|