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

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

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


Revision 1.8 - (hide 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 heimbach 1.2
2     #include "CTRL_CPPOPTIONS.h"
3    
4     subroutine ctrl_set_unpack_yz(
5 heimbach 1.3 & cunit, ivartype, fname, masktype, weighttype,
6 heimbach 1.2 & weightfld, nwetglobal, mythid)
7    
8     c ==================================================================
9     c SUBROUTINE ctrl_set_unpack_yz
10     c ==================================================================
11     c
12 heimbach 1.3 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 heimbach 1.2 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 heimbach 1.3 character*( 80) weighttype
46 heimbach 1.2 _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 heimbach 1.7 integer ii,jj,kk
60 heimbach 1.2 integer il
61 heimbach 1.3 integer irec,iobcs,nrec_nl
62 heimbach 1.2 integer itlo,ithi
63     integer jtlo,jthi
64     integer jmin,jmax
65     integer imin,imax
66    
67     integer cbuffindex
68    
69 heimbach 1.7 real*4 cbuff ( nsx*npx*sny*nsy*npy )
70 heimbach 1.2 _RL globfldyz( nsx,npx,sny,nsy,npy,nr )
71 heimbach 1.3 _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 heimbach 1.2
77     cgg(
78     integer igg
79     _RL gg
80 heimbach 1.3 character*(80) weightname
81 heimbach 1.2 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 j = jmin,jmax
104     do ip = 1,nPx
105     do bi = itlo,ithi
106     globfldyz(bi,ip,j,bj,jp,k) = 0. _d 0
107 heimbach 1.3 do iobcs=1,nobcs
108     globmskyz(bi,ip,j,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 heimbach 1.2 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 heimbach 1.3 do iobcs=1,nobcs
141     call MDSREADFIELD_YZ_GL(
142     & masktype, ctrlprec, 'RL',
143     & Nr, globmskyz(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_YZ_GL(
149     & weightname, ctrlprec, 'RL',
150     & Nr, weightfldyz(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 j = jmin,jmax
158     do ip = 1,nPx
159     do bi = itlo,ithi
160     weightfldyz(bi,ip,j,bj,jp,k,iobcs) = wbaro
161     enddo
162     enddo
163     enddo
164     enddo
165     enddo
166     endif
167     #endif
168     enddo
169    
170     nrec_nl=int(ncvarrecs(ivartype)/snx)
171     do irec = 1, nrec_nl
172 heimbach 1.2 cgg do iobcs = 1, nobcs
173     cgg And now back-calculate what iobcs should be.
174 heimbach 1.3 do i=1,snx
175     iobcs= mod((irec-1)*snx+i-1,nobcs)+1
176    
177     read(cunit) filencvarindex(ivartype)
178     if (filencvarindex(ivartype) .NE. ncvarindex(ivartype))
179     & then
180     print *, 'ctrl_set_unpack_yz: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     do jp = 1,nPy
205     do bj = jtlo,jthi
206     do j = jmin,jmax
207     do ip = 1,nPx
208     do bi = itlo,ithi
209 heimbach 1.7 ii=mod((i-1)*nr*sny+(k-1)*sny+j-1,snx)+1
210     jj=mod(((i-1)*nr*sny+(k-1)*sny+j-1)/snx,sny)+1
211     kk=int((i-1)*nr*sny+(k-1)*sny+j-1)/(snx*sny)+1
212 heimbach 1.3 if ( globmskyz(bi,ip,j,bj,jp,k,iobcs) .ne. 0. ) then
213     cbuffindex = cbuffindex + 1
214 heimbach 1.7 globfld3d(ii,bi,ip,jj,bj,jp,kk) =
215 heimbach 1.4 & cbuff(cbuffindex)
216 heimbach 1.3 #ifdef ALLOW_NONDIMENSIONAL_CONTROL_IO
217 heimbach 1.7 globfld3d(ii,bi,ip,jj,bj,jp,kk) =
218     & globfld3d(ii,bi,ip,jj,bj,jp,kk)/
219 heimbach 1.3 # ifdef CTRL_UNPACK_PRECISE
220     & sqrt(weightfldyz(bi,ip,j,bj,jp,k,iobcs))
221     # else
222     & sqrt(weightfld(k,iobcs))
223     # endif
224     #endif /* ALLOW_NONDIMENSIONAL_CONTROL_IO */
225     else
226 heimbach 1.7 globfld3d(ii,bi,ip,jj,bj,jp,kk) = 0. _d 0
227 heimbach 1.3 endif
228     enddo
229     enddo
230     enddo
231     enddo
232     enddo
233     c
234     c -- end of k loop --
235     enddo
236     c -- end of i loop --
237     enddo
238 heimbach 1.2
239 heimbach 1.3 call MDSWRITEFIELD_3d_GL( fname, ctrlprec, 'RL',
240     & Nr, globfld3d, irec,
241     & optimcycle, mythid)
242    
243     c -- end of iobcs loop -- This loop has been removed.
244     cgg enddo
245     c -- end of irec loop --
246     enddo
247    
248     do irec = nrec_nl*snx+1,ncvarrecs(ivartype)
249     iobcs= mod(irec-1,nobcs)+1
250 heimbach 1.2
251     read(cunit) filencvarindex(ivartype)
252     if (filencvarindex(ivartype) .NE. ncvarindex(ivartype))
253     & then
254     print *, 'ctrl_set_unpack_yz:WARNING: wrong ncvarindex ',
255     & filencvarindex(ivartype), ncvarindex(ivartype)
256     STOP 'in S/R ctrl_unpack'
257     endif
258     read(cunit) filej
259     read(cunit) filei
260     do k = 1, Nr
261 heimbach 1.3 cbuffindex = nwetglobal(k,iobcs)
262 heimbach 1.2 if ( cbuffindex .gt. 0 ) then
263     read(cunit) filencbuffindex
264     if (filencbuffindex .NE. cbuffindex) then
265     print *, 'WARNING: wrong cbuffindex ',
266     & filencbuffindex, cbuffindex
267     STOP 'in S/R ctrl_unpack'
268     endif
269     read(cunit) filek
270     if (filek .NE. k) then
271     print *, 'WARNING: wrong k ',
272     & filek, k
273     STOP 'in S/R ctrl_unpack'
274     endif
275     read(cunit) (cbuff(ii), ii=1,cbuffindex)
276     endif
277     cbuffindex = 0
278     do jp = 1,nPy
279     do bj = jtlo,jthi
280     do j = jmin,jmax
281     do ip = 1,nPx
282     do bi = itlo,ithi
283 heimbach 1.3 if ( globmskyz(bi,ip,j,bj,jp,k,iobcs) .ne. 0. ) then
284 heimbach 1.2 cbuffindex = cbuffindex + 1
285     globfldyz(bi,ip,j,bj,jp,k) = cbuff(cbuffindex)
286     #ifdef ALLOW_NONDIMENSIONAL_CONTROL_IO
287     globfldyz(bi,ip,j,bj,jp,k) =
288     & globfldyz(bi,ip,j,bj,jp,k)/
289 heimbach 1.3 # ifdef CTRL_UNPACK_PRECISE
290     & sqrt(weightfldyz(bi,ip,j,bj,jp,k,iobcs))
291     # else
292 heimbach 1.2 & sqrt(weightfld(k,iobcs))
293 heimbach 1.3 # endif
294     #endif /* ALLOW_NONDIMENSIONAL_CONTROL_IO */
295 heimbach 1.2 else
296     globfldyz(bi,ip,j,bj,jp,k) = 0. _d 0
297     endif
298     enddo
299     enddo
300     enddo
301     enddo
302     enddo
303     c
304 heimbach 1.3 c -- end of k loop
305 heimbach 1.2 enddo
306    
307     call MDSWRITEFIELD_YZ_GL( fname, ctrlprec, 'RL',
308     & Nr, globfldyz, irec,
309     & optimcycle, mythid)
310    
311     c -- end of iobcs loop -- This loop has been removed.
312     cgg enddo
313     c -- end of irec loop --
314     enddo
315    
316     _END_MASTER( mythid )
317    
318     return
319     end
320    

  ViewVC Help
Powered by ViewVC 1.1.22