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

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

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


Revision 1.7 - (show annotations) (download)
Thu Nov 6 22:05:08 2003 UTC (20 years, 7 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint52l_pre, hrcube4, hrcube5, checkpoint52d_pre, checkpoint52j_pre, checkpoint52l_post, checkpoint52k_post, checkpoint53, checkpoint52, checkpoint52f_post, checkpoint52i_pre, hrcube_1, hrcube_2, hrcube_3, checkpoint52e_pre, checkpoint52e_post, checkpoint52b_pre, checkpoint52m_post, checkpoint52b_post, checkpoint52c_post, checkpoint52f_pre, checkpoint53c_post, checkpoint53a_post, checkpoint52d_post, checkpoint52a_pre, checkpoint52i_post, checkpoint52h_pre, checkpoint52j_post, branch-netcdf, checkpoint52n_post, checkpoint53b_pre, checkpoint53b_post, checkpoint52a_post, ecco_c52_e35, checkpoint53d_pre, checkpoint51u_post
Branch point for: netcdf-sm0
Changes since 1.6: +9 -13 lines
o merging from ecco-branch
o cleaned some cross-dependencies and updated CPP options

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

  ViewVC Help
Powered by ViewVC 1.1.22