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

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

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


Revision 1.9 - (hide annotations) (download)
Tue Jan 4 22:02:31 2005 UTC (19 years, 4 months ago) by heimbach
Branch: MAIN
Changes since 1.8: +90 -8 lines
o Add ctrlvec diagnostics in pack/unpack for nondimensional I/O
o May be enabled via doPackDiag

1 heimbach 1.2
2     #include "CTRL_CPPOPTIONS.h"
3    
4     subroutine ctrl_set_unpack_xyz(
5 heimbach 1.9 & lxxadxx, cunit, ivartype, fname, masktype, weighttype,
6 heimbach 1.2 & weightfld, nwetglobal, mythid)
7    
8     c ==================================================================
9     c SUBROUTINE ctrl_set_unpack_xyz
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 Use a more precise nondimensionalization that depends on (x,y)
15     c Added weighttype to the argument list so that I can geographically
16     c vary the nondimensionalization.
17     c gebbie@mit.edu, 18-Mar-2003
18 heimbach 1.2 c
19     c ==================================================================
20    
21     implicit none
22    
23     c == global variables ==
24    
25     #include "EEPARAMS.h"
26     #include "SIZE.h"
27     #include "PARAMS.h"
28     #include "GRID.h"
29    
30     #include "ctrl.h"
31     #include "optim.h"
32    
33     c == routine arguments ==
34    
35 heimbach 1.9 logical lxxadxx
36 heimbach 1.2 integer cunit
37     integer ivartype
38     character*( 80) fname
39     character* (5) masktype
40 heimbach 1.3 character*( 80) weighttype
41 heimbach 1.2 _RL weightfld( nr,nsx,nsy )
42     integer nwetglobal(nr)
43     integer mythid
44    
45     c == local variables ==
46    
47     integer bi,bj
48     integer ip,jp
49     integer i,j,k
50     integer ii
51     integer il
52     integer irec
53     integer itlo,ithi
54     integer jtlo,jthi
55     integer jmin,jmax
56     integer imin,imax
57    
58     integer cbuffindex
59    
60     _RL globmsk ( snx,nsx,npx,sny,nsy,npy,nr )
61     _RL globfld3d( snx,nsx,npx,sny,nsy,npy,nr )
62 heimbach 1.3 #ifdef CTRL_UNPACK_PRECISE
63     _RL weightfld3d( snx,nsx,npx,sny,nsy,npy,nr )
64     #endif
65 heimbach 1.9 real*4 cbuff ( snx*nsx*npx*sny*nsy*npy )
66     real*4 globfldtmp2( snx,nsx,npx,sny,nsy,npy )
67     real*4 globfldtmp3( snx,nsx,npx,sny,nsy,npy )
68 heimbach 1.2
69     character*(128) cfile
70 heimbach 1.3 character*(80) weightname
71 heimbach 1.2
72 heimbach 1.9 #ifdef CTRL_DELZNORM
73     _RL delZnorm
74     #endif
75     integer reclen, irectrue
76     integer cunit2, cunit3
77     character*(80) cfile2, cfile3
78    
79 heimbach 1.2 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 heimbach 1.9 #ifdef CTRL_DELZNORM
96     delZnorm = 0.
97     do k = 1, Nr
98     delZnorm = delZnorm + delR(k)/FLOAT(Nr)
99     enddo
100     #endif
101    
102 heimbach 1.2 c Initialise temporary file
103     do k = 1,nr
104     do jp = 1,nPy
105     do bj = jtlo,jthi
106     do j = jmin,jmax
107     do ip = 1,nPx
108     do bi = itlo,ithi
109     do i = imin,imax
110 heimbach 1.9 globfld3d (i,bi,ip,j,bj,jp,k) = 0. _d 0
111     globmsk (i,bi,ip,j,bj,jp,k) = 0. _d 0
112     globfldtmp2(i,bi,ip,j,bj,jp) = 0.
113     globfldtmp3(i,bi,ip,j,bj,jp) = 0.
114 heimbach 1.2 enddo
115     enddo
116     enddo
117     enddo
118     enddo
119     enddo
120     enddo
121    
122     c-- Only the master thread will do I/O.
123     _BEGIN_MASTER( mythid )
124    
125 heimbach 1.9 #ifdef CTRL_DELZNORM
126     do k = 1, nr
127     print *, 'ph-delznorm ', k, delZnorm, delR(k)
128     print *, 'ph-weight ', weightfld(k,1,1)
129     enddo
130     #endif
131    
132     if ( doPackDiag ) then
133     write(cfile2(1:80),'(80a)') ' '
134     write(cfile3(1:80),'(80a)') ' '
135     if ( lxxadxx ) then
136     write(cfile2(1:80),'(a,I2.2,a,I4.4,a)')
137     & 'diag_unpack_nondim_ctrl_',
138     & ivartype, '_', optimcycle, '.bin'
139     write(cfile3(1:80),'(a,I2.2,a,I4.4,a)')
140     & 'diag_unpack_dimens_ctrl_',
141     & ivartype, '_', optimcycle, '.bin'
142     else
143     write(cfile2(1:80),'(a,I2.2,a,I4.4,a)')
144     & 'diag_unpack_nondim_grad_',
145     & ivartype, '_', optimcycle, '.bin'
146     write(cfile3(1:80),'(a,I2.2,a,I4.4,a)')
147     & 'diag_unpack_dimens_grad_',
148     & ivartype, '_', optimcycle, '.bin'
149     endif
150    
151     reclen = FLOAT(snx*nsx*npx*sny*nsy*npy*4)
152     call mdsfindunit( cunit2, mythid )
153     open( cunit2, file=cfile2, status='new',
154     & access='direct', recl=reclen )
155     call mdsfindunit( cunit3, mythid )
156     open( cunit3, file=cfile3, status='new',
157     & access='direct', recl=reclen )
158     endif
159    
160 heimbach 1.3 #ifdef CTRL_UNPACK_PRECISE
161     il=ilnblnk( weighttype)
162     write(weightname(1:80),'(80a)') ' '
163     write(weightname(1:80),'(a)') weighttype(1:il)
164    
165     call MDSREADFIELD_3D_GL(
166     & weightname, ctrlprec, 'RL',
167     & Nr, weightfld3d, 1, mythid)
168     #endif
169    
170 heimbach 1.2 call MDSREADFIELD_3D_GL(
171     & masktype, ctrlprec, 'RL',
172     & Nr, globmsk, 1, mythid)
173    
174     do irec = 1, ncvarrecs(ivartype)
175     read(cunit) filencvarindex(ivartype)
176     if (filencvarindex(ivartype) .NE. ncvarindex(ivartype))
177     & then
178     print *, 'ctrl_set_unpack_xyz:WARNING: wrong ncvarindex ',
179     & filencvarindex(ivartype), ncvarindex(ivartype)
180     STOP 'in S/R ctrl_unpack'
181     endif
182     read(cunit) filej
183     read(cunit) filei
184     do k = 1, Nr
185 heimbach 1.9 irectrue = (irec-1)*nr + k
186 heimbach 1.2 cbuffindex = nwetglobal(k)
187     if ( cbuffindex .gt. 0 ) then
188     read(cunit) filencbuffindex
189     if (filencbuffindex .NE. cbuffindex) then
190     print *, 'WARNING: wrong cbuffindex ',
191     & filencbuffindex, cbuffindex
192     STOP 'in S/R ctrl_unpack'
193     endif
194     read(cunit) filek
195     if (filek .NE. k) then
196     print *, 'WARNING: wrong k ',
197     & filek, k
198     STOP 'in S/R ctrl_unpack'
199     endif
200     read(cunit) (cbuff(ii), ii=1,cbuffindex)
201     endif
202     cbuffindex = 0
203     do jp = 1,nPy
204     do bj = jtlo,jthi
205     do j = jmin,jmax
206     do ip = 1,nPx
207     do bi = itlo,ithi
208     do i = imin,imax
209     if ( globmsk(i,bi,ip,j,bj,jp,k) .ne. 0. ) then
210     cbuffindex = cbuffindex + 1
211     globfld3d(i,bi,ip,j,bj,jp,k) = cbuff(cbuffindex)
212 heimbach 1.9 cph(
213     globfldtmp2(i,bi,ip,j,bj,jp) = cbuff(cbuffindex)
214     cph)
215 heimbach 1.2 #ifdef ALLOW_NONDIMENSIONAL_CONTROL_IO
216 heimbach 1.9 if ( lxxadxx ) then
217     globfld3d(i,bi,ip,j,bj,jp,k) =
218     & globfld3d(i,bi,ip,j,bj,jp,k)
219     # ifdef CTRL_UNPACK_PRECISE
220     & / sqrt(weightfld3d(i,bi,ip,j,bj,jp,k))
221     # else
222     & / sqrt(weightfld(k,bi,bj))
223     # endif
224     else
225     globfld3d(i,bi,ip,j,bj,jp,k) =
226     & globfld3d(i,bi,ip,j,bj,jp,k)
227 heimbach 1.3 # ifdef CTRL_UNPACK_PRECISE
228 heimbach 1.9 & * sqrt(weightfld3d(i,bi,ip,j,bj,jp,k))
229 heimbach 1.3 # else
230 heimbach 1.9 & * sqrt(weightfld(k,bi,bj))
231 heimbach 1.3 # endif
232 heimbach 1.9 endif
233 heimbach 1.3 #endif /* ALLOW_NONDIMENSIONAL_CONTROL_IO */
234 heimbach 1.2 else
235     globfld3d(i,bi,ip,j,bj,jp,k) = 0. _d 0
236     endif
237 heimbach 1.9 cph(
238     globfldtmp3(i,bi,ip,j,bj,jp) =
239     & globfld3d(i,bi,ip,j,bj,jp,k)
240     cph)
241 heimbach 1.2 enddo
242     enddo
243     enddo
244     enddo
245     enddo
246     enddo
247     c
248 heimbach 1.9 if ( doPackDiag ) then
249     write(cunit2,rec=irectrue) globfldtmp2
250     write(cunit3,rec=irectrue) globfldtmp3
251     endif
252     c
253 heimbach 1.2 enddo
254    
255     call MDSWRITEFIELD_3D_GL( fname, ctrlprec, 'RL',
256     & Nr, globfld3d,
257     & irec, optimcycle, mythid)
258    
259     enddo
260    
261 heimbach 1.9 if ( doPackDiag ) then
262     close ( cunit2 )
263     close ( cunit3 )
264     endif
265    
266 heimbach 1.2 _END_MASTER( mythid )
267    
268     return
269     end
270    

  ViewVC Help
Powered by ViewVC 1.1.22