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