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

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

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


Revision 1.8 - (show annotations) (download)
Fri May 28 16:04:42 2004 UTC (20 years ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint54d_post, checkpoint54e_post, checkpoint55, checkpoint54, checkpoint54f_post, checkpoint55i_post, checkpoint55c_post, checkpoint53d_post, checkpoint54b_post, checkpoint55g_post, checkpoint55d_post, checkpoint54a_pre, checkpoint55d_pre, checkpoint55j_post, checkpoint54a_post, checkpoint55h_post, checkpoint55b_post, checkpoint55f_post, checkpoint53g_post, checkpoint53f_post, checkpoint55a_post, checkpoint55e_post, checkpoint54c_post
Changes since 1.7: +0 -19 lines
Use ctrl_pack/unpack as standalone to map back and forth
between xx_/adxx_ and vector
(useful when analysing wetpoint gradient- and control-VECTOR)
Needs modified the_model_main.F

1
2 #include "CTRL_CPPOPTIONS.h"
3
4 subroutine ctrl_set_unpack_xz(
5 & cunit, ivartype, fname, masktype, weighttype,
6 & weightfld, nwetglobal, mythid)
7
8 c ==================================================================
9 c SUBROUTINE ctrl_set_unpack_xz
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 ( snx*nsx*npx*nsy*npy )
70 _RL globfldxz( snx,nsx,npx,nsy,npy,nr )
71 _RL globfld3d( snx,nsx,npx,sny,nsy,npy,nr )
72 _RL globmskxz( snx,nsx,npx,nsy,npy,nr,nobcs )
73 #ifdef CTRL_UNPACK_PRECISE
74 _RL weightfldxz( snx,nsx,npx,nsy,npy,nr,nobcs )
75 #endif
76
77 cgg(
78 integer igg
79 _RL gg
80 character*(80) weightname
81 cgg)
82
83 c == external ==
84
85 integer ilnblnk
86 external ilnblnk
87
88 cc == end of interface ==
89
90 jtlo = 1
91 jthi = nsy
92 itlo = 1
93 ithi = nsx
94 jmin = 1
95 jmax = sny
96 imin = 1
97 imax = snx
98
99 c Initialise temporary file
100 do k = 1,nr
101 do jp = 1,nPy
102 do bj = jtlo,jthi
103 do ip = 1,nPx
104 do bi = itlo,ithi
105 do i = imin,imax
106 globfldxz(i,bi,ip,bj,jp,k) = 0. _d 0
107 do iobcs=1,nobcs
108 globmskxz(i,bi,ip,bj,jp,k,iobcs) = 0. _d 0
109 enddo
110 enddo
111 enddo
112 enddo
113 enddo
114 enddo
115 enddo
116 c Initialise temporary file
117 do k = 1,nr
118 do jp = 1,nPy
119 do bj = jtlo,jthi
120 do j = jmin,jmax
121 do ip = 1,nPx
122 do bi = itlo,ithi
123 do i = imin,imax
124 globfld3d(i,bi,ip,j,bj,jp,k) = 0. _d 0
125 enddo
126 enddo
127 enddo
128 enddo
129 enddo
130 enddo
131 enddo
132
133 #ifndef ALLOW_ECCO_OPTIMIZATION
134 optimcycle = 0
135 #endif
136
137 c-- Only the master thread will do I/O.
138 _BEGIN_MASTER( mythid )
139
140 do iobcs=1,nobcs
141 call MDSREADFIELD_XZ_GL(
142 & masktype, ctrlprec, 'RL',
143 & Nr, globmskxz(1,1,1,1,1,1,iobcs), iobcs,mythid)
144 #ifdef CTRL_UNPACK_PRECISE
145 il=ilnblnk( weighttype)
146 write(weightname(1:80),'(80a)') ' '
147 write(weightname(1:80),'(a)') weighttype(1:il)
148 call MDSREADFIELD_XZ_GL(
149 & weightname, ctrlprec, 'RL',
150 & Nr, weightfldxz(1,1,1,1,1,1,iobcs), iobcs, mythid)
151 CGG One special exception: barotropic velocity should be nondimensionalized
152 cgg differently. Probably introduce new variable.
153 if (iobcs .eq. 3 .or. iobcs .eq. 4) then
154 k = 1
155 do jp = 1,nPy
156 do bj = jtlo,jthi
157 do ip = 1,nPx
158 do bi = itlo,ithi
159 do i = imin,imax
160 weightfldxz(i,bi,ip,bj,jp,k,iobcs) = wbaro
161 enddo
162 enddo
163 enddo
164 enddo
165 enddo
166 endif
167 #endif /* CTRL_UNPACK_PRECISE */
168 enddo
169
170 nrec_nl=int(ncvarrecs(ivartype)/sny)
171 do irec = 1, nrec_nl
172 cgg do iobcs = 1, nobcs
173 cgg And now back-calculate what iobcs should be.
174 do j=1,sny
175 iobcs= mod((irec-1)*sny+j-1,nobcs)+1
176
177 read(cunit) filencvarindex(ivartype)
178 if (filencvarindex(ivartype) .NE. ncvarindex(ivartype))
179 & then
180 print *, 'ctrl-set_unpack:xz:WARNING: wrong ncvarindex ',
181 & filencvarindex(ivartype), ncvarindex(ivartype)
182 STOP 'in S/R ctrl_unpack'
183 endif
184 read(cunit) filej
185 read(cunit) filei
186 do k = 1, Nr
187 cbuffindex = nwetglobal(k,iobcs)
188 if ( cbuffindex .gt. 0 ) then
189 read(cunit) filencbuffindex
190 if (filencbuffindex .NE. cbuffindex) then
191 print *, 'WARNING: wrong cbuffindex ',
192 & filencbuffindex, cbuffindex
193 STOP 'in S/R ctrl_unpack'
194 endif
195 read(cunit) filek
196 if (filek .NE. k) then
197 print *, 'WARNING: wrong k ',
198 & filek, k
199 STOP 'in S/R ctrl_unpack'
200 endif
201 read(cunit) (cbuff(ii), ii=1,cbuffindex)
202 endif
203 cbuffindex = 0
204 jj=mod((j-1)*nr+k-1,sny)+1
205 kk=int((j-1)*nr+k-1)/sny+1
206 do jp = 1,nPy
207 do bj = jtlo,jthi
208 do ip = 1,nPx
209 do bi = itlo,ithi
210 do i = imin,imax
211 if ( globmskxz(i,bi,ip,bj,jp,k,iobcs) .ne. 0. ) then
212 cbuffindex = cbuffindex + 1
213 globfld3d(i,bi,ip,jj,bj,jp,kk) =
214 & cbuff(cbuffindex)
215 #ifdef ALLOW_NONDIMENSIONAL_CONTROL_IO
216 globfld3d(i,bi,ip,jj,bj,jp,kk) =
217 & globfld3d(i,bi,ip,jj,bj,jp,kk)/
218 # ifdef CTRL_UNPACK_PRECISE
219 & sqrt(weightfldxz(i,bi,ip,bj,jp,k,iobcs))
220 # else
221 & sqrt(weightfld(k,iobcs))
222 # endif
223 #endif /* ALLOW_NONDIMENSIONAL_CONTROL_IO */
224 else
225 globfld3d(i,bi,ip,jj,bj,jp,kk) = 0. _d 0
226 endif
227 enddo
228 enddo
229 enddo
230 enddo
231 enddo
232 c
233 c -- end of k loop --
234 enddo
235 c -- end of j loop --
236 enddo
237
238 call MDSWRITEFIELD_3D_GL( fname, ctrlprec, 'RL',
239 & Nr, globfld3d, irec,
240 & optimcycle, mythid)
241
242 c -- end of iobcs loop -- This loop removed. 3-28-02.
243 cgg enddo
244 c -- end of irec loop --
245 enddo
246
247 do irec = nrec_nl*sny+1, ncvarrecs(ivartype)
248 cgg do iobcs = 1, nobcs
249 cgg And now back-calculate what iobcs should be.
250 iobcs= mod(irec-1,nobcs)+1
251
252 read(cunit) filencvarindex(ivartype)
253 if (filencvarindex(ivartype) .NE. ncvarindex(ivartype))
254 & then
255 print *, 'ctrl-set_unpack:xz:WARNING: wrong ncvarindex ',
256 & filencvarindex(ivartype), ncvarindex(ivartype)
257 STOP 'in S/R ctrl_unpack'
258 endif
259 read(cunit) filej
260 read(cunit) filei
261 do k = 1, Nr
262 cbuffindex = nwetglobal(k,iobcs)
263 if ( cbuffindex .gt. 0 ) then
264 read(cunit) filencbuffindex
265 if (filencbuffindex .NE. cbuffindex) then
266 print *, 'WARNING: wrong cbuffindex ',
267 & filencbuffindex, cbuffindex
268 STOP 'in S/R ctrl_unpack'
269 endif
270 read(cunit) filek
271 if (filek .NE. k) then
272 print *, 'WARNING: wrong k ',
273 & filek, k
274 STOP 'in S/R ctrl_unpack'
275 endif
276 read(cunit) (cbuff(ii), ii=1,cbuffindex)
277 endif
278 cbuffindex = 0
279 do jp = 1,nPy
280 do bj = jtlo,jthi
281 do ip = 1,nPx
282 do bi = itlo,ithi
283 do i = imin,imax
284 if ( globmskxz(i,bi,ip,bj,jp,k,iobcs) .ne. 0. ) then
285 cbuffindex = cbuffindex + 1
286 globfldxz(i,bi,ip,bj,jp,k) = cbuff(cbuffindex)
287 #ifdef ALLOW_NONDIMENSIONAL_CONTROL_IO
288 globfldxz(i,bi,ip,bj,jp,k) =
289 & globfldxz(i,bi,ip,bj,jp,k)/
290 # ifdef CTRL_UNPACK_PRECISE
291 & sqrt(weightfldxz(i,bi,ip,bj,jp,k,iobcs))
292 # else
293 & sqrt(weightfld(k,iobcs))
294 # endif
295 #endif /* ALLOW_NONDIMENSIONAL_CONTROL_IO */
296 else
297 globfldxz(i,bi,ip,bj,jp,k) = 0. _d 0
298 endif
299 enddo
300 enddo
301 enddo
302 enddo
303 enddo
304 c
305 c -- end of k loop --
306 enddo
307
308 call MDSWRITEFIELD_XZ_GL( fname, ctrlprec, 'RL',
309 & Nr, globfldxz, irec,
310 & optimcycle, mythid)
311
312 c -- end of iobcs loop -- This loop removed. 3-28-02.
313 cgg enddo
314 c -- end of irec loop --
315 enddo
316
317 _END_MASTER( mythid )
318
319 return
320 end
321
322
323
324
325

  ViewVC Help
Powered by ViewVC 1.1.22