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

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

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


Revision 1.17 - (hide annotations) (download)
Tue Jan 4 22:02:31 2005 UTC (19 years, 5 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint57d_post, checkpoint57c_post, checkpoint57c_pre, eckpoint57e_pre
Changes since 1.16: +3 -3 lines
o Add ctrlvec diagnostics in pack/unpack for nondimensional I/O
o May be enabled via doPackDiag

1 heimbach 1.17 C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_pack.F,v 1.16 2004/12/04 17:57:37 heimbach Exp $
2 heimbach 1.16 C $Name: $
3 heimbach 1.1
4 heimbach 1.12 #include "PACKAGES_CONFIG.h"
5 heimbach 1.1 #include "CTRL_CPPOPTIONS.h"
6    
7 heimbach 1.11 subroutine ctrl_pack( first, mythid )
8 heimbach 1.8
9     c ==================================================================
10     c SUBROUTINE ctrl_pack
11     c ==================================================================
12     c
13     c o Compress the control vector such that only ocean points are
14     c written to file.
15     c
16     c started: Christian Eckert eckert@mit.edu 10-Mar=2000
17     c
18     c changed: Patrick Heimbach heimbach@mit.edu 06-Jun-2000
19     c - Transferred some filename declarations
20     c from here to namelist in ctrl_init
21     c
22     c Patrick Heimbach heimbach@mit.edu 16-Jun-2000
23     c - single file name convention with or without
24     c ALLOW_ECCO_OPTIMIZATION
25     c
26     c G. Gebbie, added open boundary control packing,
27     c gebbie@mit.edu 18 -Mar- 2003
28     c
29 heimbach 1.12 c heimbach@mit.edu totally restructured 28-Oct-2003
30 heimbach 1.11 c
31 heimbach 1.8 c ==================================================================
32     c SUBROUTINE ctrl_pack
33     c ==================================================================
34    
35 heimbach 1.1 implicit none
36    
37     c == global variables ==
38 heimbach 1.5
39 heimbach 1.1 #include "EEPARAMS.h"
40     #include "SIZE.h"
41     #include "PARAMS.h"
42     #include "GRID.h"
43 heimbach 1.5
44 heimbach 1.1 #include "ctrl.h"
45 heimbach 1.14 #include "optim.h"
46 heimbach 1.5
47 heimbach 1.16 #ifdef ALLOW_COST
48     # include "cost.h"
49     #endif
50 heimbach 1.12 #ifdef ALLOW_ECCO
51     # include "ecco_cost.h"
52     #else
53     # include "ctrl_weights.h"
54     #endif
55    
56 heimbach 1.1 c == routine arguments ==
57 heimbach 1.5
58 heimbach 1.11 logical first
59 heimbach 1.1 integer mythid
60    
61 heimbach 1.11 #ifndef EXCLUDE_CTRL_PACK
62 heimbach 1.1 c == local variables ==
63    
64 heimbach 1.11 _RL fcloc
65    
66 heimbach 1.5 integer i, j, k
67 heimbach 1.1 integer ii
68     integer il
69     integer irec
70 heimbach 1.5 integer ig,jg
71     integer ivartype
72     integer iobcs
73 heimbach 1.1
74     logical doglobalread
75     logical ladinit
76 heimbach 1.5 integer cbuffindex
77 heimbach 1.11 logical lxxadxx
78    
79 heimbach 1.5 integer cunit
80 heimbach 1.11 integer ictrlgrad
81 heimbach 1.1
82     character*(128) cfile
83 heimbach 1.5 character*( 80) weighttype
84    
85 heimbach 1.1 c == external ==
86 heimbach 1.5
87 heimbach 1.1 integer ilnblnk
88     external ilnblnk
89    
90     c == end of interface ==
91    
92 heimbach 1.5 #ifndef ALLOW_ECCO_OPTIMIZATION
93 heimbach 1.11 fmin = 0. _d 0
94 heimbach 1.5 #endif
95 heimbach 1.1
96     c-- Tiled files are used.
97     doglobalread = .false.
98    
99     c-- Initialise adjoint variables on active files.
100     ladinit = .false.
101    
102 heimbach 1.14 c-- Initialise global buffer index
103     nbuffglobal = 0
104    
105 heimbach 1.5 c-- Assign file names.
106    
107 heimbach 1.11 call ctrl_set_fname(xx_theta_file, fname_theta, mythid)
108     call ctrl_set_fname(xx_salt_file, fname_salt, mythid)
109     call ctrl_set_fname(xx_hflux_file, fname_hflux, mythid)
110     call ctrl_set_fname(xx_sflux_file, fname_sflux, mythid)
111     call ctrl_set_fname(xx_tauu_file, fname_tauu, mythid)
112     call ctrl_set_fname(xx_tauv_file, fname_tauv, mythid)
113     call ctrl_set_fname(xx_atemp_file, fname_atemp, mythid)
114     call ctrl_set_fname(xx_aqh_file, fname_aqh, mythid)
115     call ctrl_set_fname(xx_uwind_file, fname_uwind, mythid)
116     call ctrl_set_fname(xx_vwind_file, fname_vwind, mythid)
117     call ctrl_set_fname(xx_obcsn_file, fname_obcsn, mythid)
118     call ctrl_set_fname(xx_obcss_file, fname_obcss, mythid)
119     call ctrl_set_fname(xx_obcsw_file, fname_obcsw, mythid)
120     call ctrl_set_fname(xx_obcse_file, fname_obcse, mythid)
121     call ctrl_set_fname(xx_diffkr_file, fname_diffkr, mythid)
122     call ctrl_set_fname(xx_kapgm_file, fname_kapgm, mythid)
123     call ctrl_set_fname(xx_tr1_file, fname_tr1, mythid)
124     call ctrl_set_fname(xx_sst_file, fname_sst, mythid)
125     call ctrl_set_fname(xx_sss_file, fname_sss, mythid)
126     call ctrl_set_fname(xx_hfacc_file, fname_hfacc, mythid)
127     call ctrl_set_fname(xx_efluxy_file, fname_efluxy, mythid)
128     call ctrl_set_fname(xx_efluxp_file, fname_efluxp, mythid)
129     call ctrl_set_fname(xx_bottomdrag_file, fname_bottomdrag, mythid)
130 heimbach 1.5
131 heimbach 1.1 c
132 heimbach 1.5 c-- Only the master thread will do I/O.
133 heimbach 1.1 _BEGIN_MASTER( mythid )
134    
135 heimbach 1.13 if ( first ) then
136 heimbach 1.11 c >>> Initialise control vector for optimcycle=0 <<<
137     lxxadxx = .TRUE.
138     ictrlgrad = 1
139     fcloc = fmin
140     write(cfile(1:128),'(4a,i4.4)')
141 heimbach 1.13 & ctrlname(1:9),'_',yctrlid(1:10),
142     & yctrlpospack, optimcycle
143 heimbach 1.17 print *, 'ph-pack: packing ', ctrlname(1:9)
144 heimbach 1.11 else
145 heimbach 1.1 c >>> Write gradient vector <<<
146 heimbach 1.11 lxxadxx = .FALSE.
147     ictrlgrad = 2
148     fcloc = fc
149 heimbach 1.5 write(cfile(1:128),'(4a,i4.4)')
150 heimbach 1.13 & costname(1:9),'_',yctrlid(1:10),
151     & yctrlpospack, optimcycle
152 heimbach 1.17 print *, 'ph-pack: packing ', costname(1:9)
153 heimbach 1.11 endif
154 heimbach 1.1
155 heimbach 1.11 call mdsfindunit( cunit, mythid )
156     open( cunit, file = cfile,
157     & status = 'unknown',
158     & form = 'unformatted',
159     & access = 'sequential' )
160 heimbach 1.1
161     c-- Header information.
162 mlosch 1.15 write(cunit) nvartype
163     write(cunit) nvarlength
164     write(cunit) yctrlid
165     write(cunit) optimCycle
166     write(cunit) fc
167     C place holder of obsolete variable iG
168     write(cunit) 1
169     C place holder of obsolete variable jG
170     write(cunit) 1
171     write(cunit) nsx
172     write(cunit) nsy
173     write(cunit) (nWetcGlobal(k), k=1,nr)
174     write(cunit) (nWetsGlobal(k), k=1,nr)
175     write(cunit) (nWetwGlobal(k), k=1,nr)
176 heimbach 1.7 #ifdef ALLOW_CTRL_WETV
177 mlosch 1.15 write(cunit) (nWetvGlobal(k), k=1,nr)
178 heimbach 1.7 #endif
179 heimbach 1.11
180 heimbach 1.5 #ifdef ALLOW_OBCSN_CONTROL
181     write(cunit) ((nWetobcsnGlo(k,iobcs), k=1,nr),iobcs= 1,nobcs)
182     #endif
183     #ifdef ALLOW_OBCSS_CONTROL
184     write(cunit) ((nWetobcssGlo(k,iobcs), k=1,nr),iobcs= 1,nobcs)
185     #endif
186     #ifdef ALLOW_OBCSW_CONTROL
187     write(cunit) ((nWetobcswGlo(k,iobcs), k=1,nr),iobcs= 1,nobcs)
188     #endif
189     #ifdef ALLOW_OBCSE_CONTROL
190     write(cunit) ((nWetobcseGlo(k,iobcs), k=1,nr),iobcs= 1,nobcs)
191     #endif
192 mlosch 1.15 write(cunit) (ncvarindex(i), i=1,maxcvars)
193     write(cunit) (ncvarrecs(i), i=1,maxcvars)
194     write(cunit) (ncvarxmax(i), i=1,maxcvars)
195     write(cunit) (ncvarymax(i), i=1,maxcvars)
196     write(cunit) (ncvarnrmax(i), i=1,maxcvars)
197     write(cunit) (ncvargrd(i), i=1,maxcvars)
198 heimbach 1.1 write(cunit)
199    
200     #ifdef ALLOW_THETA0_CONTROL
201 heimbach 1.5 ivartype = 1
202 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
203     write(weighttype(1:80),'(a)') "wtheta"
204 heimbach 1.5 call ctrl_set_pack_xyz(
205 heimbach 1.11 & cunit, ivartype, fname_theta(ictrlgrad), "hFacC",
206 heimbach 1.8 & weighttype, wtheta, lxxadxx, mythid)
207 heimbach 1.1 #endif
208    
209     #ifdef ALLOW_SALT0_CONTROL
210 heimbach 1.5 ivartype = 2
211 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
212     write(weighttype(1:80),'(a)') "wsalt"
213 heimbach 1.5 call ctrl_set_pack_xyz(
214 heimbach 1.11 & cunit, ivartype, fname_salt(ictrlgrad), "hFacC",
215 heimbach 1.8 & weighttype, wsalt, lxxadxx, mythid)
216 heimbach 1.5 #endif
217    
218     #if (defined (ALLOW_HFLUX_CONTROL) || \
219     defined (ALLOW_HFLUX0_CONTROL))
220     ivartype = 3
221     write(weighttype(1:80),'(80a)') ' '
222     write(weighttype(1:80),'(a)') "whflux"
223     call ctrl_set_pack_xy(
224 heimbach 1.11 & cunit, ivartype, fname_hflux(ictrlgrad), "hFacC",
225     & weighttype, lxxadxx, mythid)
226 heimbach 1.5 #endif
227    
228     #if (defined (ALLOW_SFLUX_CONTROL) || \
229     defined (ALLOW_SFLUX0_CONTROL))
230     ivartype = 4
231     write(weighttype(1:80),'(80a)') ' '
232     write(weighttype(1:80),'(a)') "wsflux"
233     call ctrl_set_pack_xy(
234 heimbach 1.11 & cunit, ivartype, fname_sflux(ictrlgrad), "hFacC",
235     & weighttype, lxxadxx, mythid)
236 heimbach 1.5 #endif
237    
238     #if (defined (ALLOW_USTRESS_CONTROL) || \
239     defined (ALLOW_TAUU0_CONTROL))
240     ivartype = 5
241     write(weighttype(1:80),'(80a)') ' '
242     write(weighttype(1:80),'(a)') "wtauu"
243     call ctrl_set_pack_xy(
244 heimbach 1.11 & cunit, ivartype, fname_tauu(ictrlgrad), "maskW",
245     & weighttype, lxxadxx, mythid)
246 heimbach 1.5 #endif
247    
248     #if (defined (ALLOW_VSTRESS_CONTROL) || \
249     defined (ALLOW_TAUV0_CONTROL))
250     ivartype = 6
251     write(weighttype(1:80),'(80a)') ' '
252     write(weighttype(1:80),'(a)') "wtauv"
253     call ctrl_set_pack_xy(
254 heimbach 1.11 & cunit, ivartype, fname_tauv(ictrlgrad), "maskS",
255     & weighttype, lxxadxx, mythid)
256 heimbach 1.5 #endif
257    
258     #ifdef ALLOW_ATEMP_CONTROL
259     ivartype = 7
260     write(weighttype(1:80),'(80a)') ' '
261     write(weighttype(1:80),'(a)') "watemp"
262     call ctrl_set_pack_xy(
263 heimbach 1.11 & cunit, ivartype, fname_atemp(ictrlgrad), "hFacC",
264     & weighttype, lxxadxx, mythid)
265 heimbach 1.5 #endif
266    
267     #ifdef ALLOW_AQH_CONTROL
268     ivartype = 8
269     write(weighttype(1:80),'(80a)') ' '
270     write(weighttype(1:80),'(a)') "waqh"
271     call ctrl_set_pack_xy(
272 heimbach 1.11 & cunit, ivartype, fname_aqh(ictrlgrad), "hFacC",
273     & weighttype, lxxadxx, mythid)
274 heimbach 1.5 #endif
275    
276     #ifdef ALLOW_UWIND_CONTROL
277     ivartype = 9
278     write(weighttype(1:80),'(80a)') ' '
279     write(weighttype(1:80),'(a)') "wuwind"
280     call ctrl_set_pack_xy(
281 heimbach 1.11 & cunit, ivartype, fname_uwind(ictrlgrad), "maskW",
282     & weighttype, lxxadxx, mythid)
283 heimbach 1.5 #endif
284    
285     #ifdef ALLOW_VWIND_CONTROL
286     ivartype = 10
287     write(weighttype(1:80),'(80a)') ' '
288     write(weighttype(1:80),'(a)') "wvwind"
289     call ctrl_set_pack_xy(
290 heimbach 1.11 & cunit, ivartype, fname_vwind(ictrlgrad), "maskS",
291     & weighttype, lxxadxx, mythid)
292 heimbach 1.5 #endif
293    
294     #ifdef ALLOW_OBCSN_CONTROL
295     ivartype = 11
296 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
297     write(weighttype(1:80),'(a)') "wobcsn"
298 heimbach 1.5 call ctrl_set_pack_xz(
299 heimbach 1.11 & cunit, ivartype, fname_obcsn(ictrlgrad), "maskobcsn",
300 heimbach 1.8 & weighttype, wobcsn, lxxadxx, mythid)
301 heimbach 1.5 #endif
302    
303     #ifdef ALLOW_OBCSS_CONTROL
304     ivartype = 12
305 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
306     write(weighttype(1:80),'(a)') "wobcss"
307 heimbach 1.5 call ctrl_set_pack_xz(
308 heimbach 1.11 & cunit, ivartype, fname_obcss(ictrlgrad), "maskobcss",
309 heimbach 1.8 & weighttype, wobcss, lxxadxx, mythid)
310 heimbach 1.5 #endif
311    
312     #ifdef ALLOW_OBCSW_CONTROL
313     ivartype = 13
314 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
315     write(weighttype(1:80),'(a)') "wobcsw"
316 heimbach 1.5 call ctrl_set_pack_yz(
317 heimbach 1.11 & cunit, ivartype, fname_obcsw(ictrlgrad), "maskobcsw",
318 heimbach 1.8 & weighttype, wobcsw, lxxadxx, mythid)
319 heimbach 1.5 #endif
320    
321     #ifdef ALLOW_OBCSE_CONTROL
322     ivartype = 14
323 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
324     write(weighttype(1:80),'(a)') "wobcse"
325 heimbach 1.5 call ctrl_set_pack_yz(
326 heimbach 1.11 & cunit, ivartype, fname_obcse(ictrlgrad), "maskobcse",
327 heimbach 1.8 & weighttype, wobcse, lxxadxx, mythid)
328 heimbach 1.1 #endif
329 heimbach 1.3
330     #ifdef ALLOW_DIFFKR_CONTROL
331 heimbach 1.5 ivartype = 15
332 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
333     write(weighttype(1:80),'(a)') "wdiffkr"
334 heimbach 1.5 call ctrl_set_pack_xyz(
335 heimbach 1.11 & cunit, ivartype, fname_diffkr(ictrlgrad), "hFacC",
336 heimbach 1.8 & weighttype, wunit, lxxadxx, mythid)
337 heimbach 1.3 #endif
338    
339     #ifdef ALLOW_KAPGM_CONTROL
340 heimbach 1.5 ivartype = 16
341 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
342     write(weighttype(1:80),'(a)') "wkapgm"
343 heimbach 1.5 call ctrl_set_pack_xyz(
344 heimbach 1.11 & cunit, ivartype, fname_kapgm(ictrlgrad), "hFacC",
345 heimbach 1.8 & weighttype, wunit, lxxadxx, mythid)
346 heimbach 1.3 #endif
347    
348 heimbach 1.5 #ifdef ALLOW_TR10_CONTROL
349     ivartype = 17
350 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
351     write(weighttype(1:80),'(a)') "wtr1"
352 heimbach 1.5 call ctrl_set_pack_xyz(
353 heimbach 1.11 & cunit, ivartype, fname_tr1(ictrlgrad), "hFacC",
354 heimbach 1.8 & weighttype, wunit, lxxadxx, mythid)
355 heimbach 1.5 #endif
356    
357 heimbach 1.6 #ifdef ALLOW_SST0_CONTROL
358     ivartype = 18
359     write(weighttype(1:80),'(80a)') ' '
360     write(weighttype(1:80),'(a)') "wsst0"
361     call ctrl_set_pack_xy(
362 heimbach 1.11 & cunit, ivartype, fname_sst(ictrlgrad), "hFacC",
363     & weighttype, lxxadxx, mythid)
364 heimbach 1.6 #endif
365    
366     #ifdef ALLOW_SSS0_CONTROL
367     ivartype = 19
368     write(weighttype(1:80),'(80a)') ' '
369     write(weighttype(1:80),'(a)') "wsss0"
370     call ctrl_set_pack_xy(
371 heimbach 1.11 & cunit, ivartype, fname_sss(ictrlgrad), "hFacC",
372     & weighttype, lxxadxx, mythid)
373 heimbach 1.6 #endif
374    
375     #ifdef ALLOW_HFACC_CONTROL
376     ivartype = 20
377 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
378     write(weighttype(1:80),'(a)') "whfacc"
379     # ifdef ALLOW_HFACC3D_CONTROL
380 heimbach 1.6 call ctrl_set_pack_xyz(
381 heimbach 1.11 & cunit, ivartype, fname_hfacc(ictrlgrad), "hFacC",
382 heimbach 1.8 & weighttype, wunit, lxxadxx, mythid)
383     # else
384 heimbach 1.6 call ctrl_set_pack_xy(
385 heimbach 1.11 & cunit, ivartype, fname_hfacc(ictrlgrad), "hFacC",
386     & weighttype, lxxadxx, mythid)
387 heimbach 1.8 # endif
388 heimbach 1.6 #endif
389    
390 heimbach 1.5 #ifdef ALLOW_EFLUXY0_CONTROL
391     ivartype = 21
392 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
393     write(weighttype(1:80),'(a)') "wefluxy0"
394 heimbach 1.5 call ctrl_set_pack_xyz(
395 heimbach 1.11 & cunit, ivartype, fname_efluxy(ictrlgrad), "hFacS",
396 heimbach 1.8 & weighttype, wunit, lxxadxx, mythid)
397 heimbach 1.5 #endif
398    
399     #ifdef ALLOW_EFLUXP0_CONTROL
400     ivartype = 22
401 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
402     write(weighttype(1:80),'(a)') "wefluxp0"
403 heimbach 1.5 call ctrl_set_pack_xyz(
404 heimbach 1.11 & cunit, ivartype, fname_efluxp(ictrlgrad), "hFacV",
405 heimbach 1.8 & weighttype, wunit, lxxadxx, mythid)
406 heimbach 1.6 #endif
407    
408     #ifdef ALLOW_BOTTOMDRAG_CONTROL
409     ivartype = 23
410     write(weighttype(1:80),'(80a)') ' '
411     write(weighttype(1:80),'(a)') "wbottomdrag"
412     call ctrl_set_pack_xy(
413 heimbach 1.11 & cunit, ivartype, fname_bottomdrag(ictrlgrad), "hFacC",
414     & weighttype, lxxadxx, mythid)
415 heimbach 1.5 #endif
416    
417     close ( cunit )
418 heimbach 1.1
419 heimbach 1.5 _END_MASTER( mythid )
420 heimbach 1.11
421     #endif /* EXCLUDE_CTRL_PACK */
422 heimbach 1.1
423     return
424     end
425    

  ViewVC Help
Powered by ViewVC 1.1.22