/[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.12 - (show annotations) (download)
Tue Oct 9 00:00:01 2007 UTC (16 years, 7 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint59m, checkpoint59l, checkpoint59i, checkpoint59k, checkpoint59j
Changes since 1.11: +8 -6 lines
add missing cvs $Header:$ or $Name:$

1 C $Header: $
2 C $Name: $
3
4 #include "CTRL_CPPOPTIONS.h"
5
6 subroutine ctrl_set_unpack_yz(
7 & cunit, ivartype, fname, masktype, weighttype,
8 & weightfld, nwetglobal, mythid)
9
10 c ==================================================================
11 c SUBROUTINE ctrl_set_unpack_yz
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 ( nsx*npx*sny*nsy*npy )
66 _RL globfldyz( nsx,npx,sny,nsy,npy,nr )
67 _RL globfld3d( snx,nsx,npx,sny,nsy,npy,nr )
68 _RL globmskyz( nsx,npx,sny,nsy,npy,nr,nobcs )
69 #ifdef CTRL_UNPACK_PRECISE
70 _RL weightfldyz( nsx,npx,sny,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 j = jmin,jmax
100 do ip = 1,nPx
101 do bi = itlo,ithi
102 globfldyz(bi,ip,j,bj,jp,k) = 0. _d 0
103 do iobcs=1,nobcs
104 globmskyz(bi,ip,j,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_YZ_GL(
134 & masktype, ctrlprec, 'RL',
135 & Nr, globmskyz(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_YZ_GL(
141 & weightname, ctrlprec, 'RL',
142 & Nr, weightfldyz(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 j = jmin,jmax
150 do ip = 1,nPx
151 do bi = itlo,ithi
152 cph weightfldyz(bi,ip,j,bj,jp,k,iobcs) = wbaro
153 enddo
154 enddo
155 enddo
156 enddo
157 enddo
158 endif
159 #endif
160 enddo
161
162 nrec_nl=int(ncvarrecs(ivartype)/snx)
163 do irec = 1, nrec_nl
164 cgg do iobcs = 1, nobcs
165 cgg And now back-calculate what iobcs should be.
166 do i=1,snx
167 iobcs= mod((irec-1)*snx+i-1,nobcs)+1
168
169 read(cunit) filencvarindex(ivartype)
170 if (filencvarindex(ivartype) .NE. ncvarindex(ivartype))
171 & then
172 print *, 'ctrl_set_unpack_yz: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 do jp = 1,nPy
197 do bj = jtlo,jthi
198 do j = jmin,jmax
199 do ip = 1,nPx
200 do bi = itlo,ithi
201 ii=mod((i-1)*nr*sny+(k-1)*sny+j-1,snx)+1
202 jj=mod(((i-1)*nr*sny+(k-1)*sny+j-1)/snx,sny)+1
203 kk=int((i-1)*nr*sny+(k-1)*sny+j-1)/(snx*sny)+1
204 if ( globmskyz(bi,ip,j,bj,jp,k,iobcs) .ne. 0. ) then
205 cbuffindex = cbuffindex + 1
206 globfld3d(ii,bi,ip,jj,bj,jp,kk) =
207 & cbuff(cbuffindex)
208 #ifdef ALLOW_NONDIMENSIONAL_CONTROL_IO
209 globfld3d(ii,bi,ip,jj,bj,jp,kk) =
210 & globfld3d(ii,bi,ip,jj,bj,jp,kk)/
211 # ifdef CTRL_UNPACK_PRECISE
212 & sqrt(weightfldyz(bi,ip,j,bj,jp,k,iobcs))
213 # else
214 & sqrt(weightfld(k,iobcs))
215 # endif
216 #endif /* ALLOW_NONDIMENSIONAL_CONTROL_IO */
217 else
218 globfld3d(ii,bi,ip,jj,bj,jp,kk) = 0. _d 0
219 endif
220 enddo
221 enddo
222 enddo
223 enddo
224 enddo
225 c
226 c -- end of k loop --
227 enddo
228 c -- end of i loop --
229 enddo
230
231 call MDSWRITEFIELD_3d_GL( fname, ctrlprec, 'RL',
232 & Nr, globfld3d, irec,
233 & optimcycle, mythid)
234
235 c -- end of iobcs loop -- This loop has been removed.
236 cgg enddo
237 c -- end of irec loop --
238 enddo
239
240 do irec = nrec_nl*snx+1,ncvarrecs(ivartype)
241 iobcs= mod(irec-1,nobcs)+1
242
243 read(cunit) filencvarindex(ivartype)
244 if (filencvarindex(ivartype) .NE. ncvarindex(ivartype))
245 & then
246 print *, 'ctrl_set_unpack_yz:WARNING: wrong ncvarindex ',
247 & filencvarindex(ivartype), ncvarindex(ivartype)
248 STOP 'in S/R ctrl_unpack'
249 endif
250 read(cunit) filej
251 read(cunit) filei
252 do k = 1, Nr
253 cbuffindex = nwetglobal(k,iobcs)
254 if ( cbuffindex .gt. 0 ) then
255 read(cunit) filencbuffindex
256 if (filencbuffindex .NE. cbuffindex) then
257 print *, 'WARNING: wrong cbuffindex ',
258 & filencbuffindex, cbuffindex
259 STOP 'in S/R ctrl_unpack'
260 endif
261 read(cunit) filek
262 if (filek .NE. k) then
263 print *, 'WARNING: wrong k ',
264 & filek, k
265 STOP 'in S/R ctrl_unpack'
266 endif
267 read(cunit) (cbuff(ii), ii=1,cbuffindex)
268 endif
269 cbuffindex = 0
270 do jp = 1,nPy
271 do bj = jtlo,jthi
272 do j = jmin,jmax
273 do ip = 1,nPx
274 do bi = itlo,ithi
275 if ( globmskyz(bi,ip,j,bj,jp,k,iobcs) .ne. 0. ) then
276 cbuffindex = cbuffindex + 1
277 globfldyz(bi,ip,j,bj,jp,k) = cbuff(cbuffindex)
278 #ifdef ALLOW_NONDIMENSIONAL_CONTROL_IO
279 globfldyz(bi,ip,j,bj,jp,k) =
280 & globfldyz(bi,ip,j,bj,jp,k)/
281 # ifdef CTRL_UNPACK_PRECISE
282 & sqrt(weightfldyz(bi,ip,j,bj,jp,k,iobcs))
283 # else
284 & sqrt(weightfld(k,iobcs))
285 # endif
286 #endif /* ALLOW_NONDIMENSIONAL_CONTROL_IO */
287 else
288 globfldyz(bi,ip,j,bj,jp,k) = 0. _d 0
289 endif
290 enddo
291 enddo
292 enddo
293 enddo
294 enddo
295 c
296 c -- end of k loop
297 enddo
298
299 call MDSWRITEFIELD_YZ_GL( fname, ctrlprec, 'RL',
300 & Nr, globfldyz, irec,
301 & optimcycle, mythid)
302
303 c -- end of iobcs loop -- This loop has been removed.
304 cgg enddo
305 c -- end of irec loop --
306 enddo
307
308 _END_MASTER( mythid )
309
310 #endif
311
312 return
313 end
314

  ViewVC Help
Powered by ViewVC 1.1.22