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

  ViewVC Help
Powered by ViewVC 1.1.22