/[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.18 - (hide annotations) (download)
Mon Feb 28 17:29:38 2005 UTC (19 years, 3 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint57f_post, checkpoint57e_post, checkpoint57f_pre
Changes since 1.17: +22 -1 lines
Adding eddy stress controls a la Ferreira et al.

1 heimbach 1.18 C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_pack.F,v 1.17 2005/01/04 22:02:31 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.18 call ctrl_set_fname(xx_edtaux_file, fname_edtaux, mythid)
131     call ctrl_set_fname(xx_edtauy_file, fname_edtauy, mythid)
132 heimbach 1.5
133 heimbach 1.1 c
134 heimbach 1.5 c-- Only the master thread will do I/O.
135 heimbach 1.1 _BEGIN_MASTER( mythid )
136    
137 heimbach 1.13 if ( first ) then
138 heimbach 1.11 c >>> Initialise control vector for optimcycle=0 <<<
139     lxxadxx = .TRUE.
140     ictrlgrad = 1
141     fcloc = fmin
142     write(cfile(1:128),'(4a,i4.4)')
143 heimbach 1.13 & ctrlname(1:9),'_',yctrlid(1:10),
144     & yctrlpospack, optimcycle
145 heimbach 1.17 print *, 'ph-pack: packing ', ctrlname(1:9)
146 heimbach 1.11 else
147 heimbach 1.1 c >>> Write gradient vector <<<
148 heimbach 1.11 lxxadxx = .FALSE.
149     ictrlgrad = 2
150     fcloc = fc
151 heimbach 1.5 write(cfile(1:128),'(4a,i4.4)')
152 heimbach 1.13 & costname(1:9),'_',yctrlid(1:10),
153     & yctrlpospack, optimcycle
154 heimbach 1.17 print *, 'ph-pack: packing ', costname(1:9)
155 heimbach 1.11 endif
156 heimbach 1.1
157 heimbach 1.11 call mdsfindunit( cunit, mythid )
158     open( cunit, file = cfile,
159     & status = 'unknown',
160     & form = 'unformatted',
161     & access = 'sequential' )
162 heimbach 1.1
163     c-- Header information.
164 mlosch 1.15 write(cunit) nvartype
165     write(cunit) nvarlength
166     write(cunit) yctrlid
167     write(cunit) optimCycle
168     write(cunit) fc
169     C place holder of obsolete variable iG
170     write(cunit) 1
171     C place holder of obsolete variable jG
172     write(cunit) 1
173     write(cunit) nsx
174     write(cunit) nsy
175     write(cunit) (nWetcGlobal(k), k=1,nr)
176     write(cunit) (nWetsGlobal(k), k=1,nr)
177     write(cunit) (nWetwGlobal(k), k=1,nr)
178 heimbach 1.7 #ifdef ALLOW_CTRL_WETV
179 mlosch 1.15 write(cunit) (nWetvGlobal(k), k=1,nr)
180 heimbach 1.7 #endif
181 heimbach 1.11
182 heimbach 1.5 #ifdef ALLOW_OBCSN_CONTROL
183     write(cunit) ((nWetobcsnGlo(k,iobcs), k=1,nr),iobcs= 1,nobcs)
184     #endif
185     #ifdef ALLOW_OBCSS_CONTROL
186     write(cunit) ((nWetobcssGlo(k,iobcs), k=1,nr),iobcs= 1,nobcs)
187     #endif
188     #ifdef ALLOW_OBCSW_CONTROL
189     write(cunit) ((nWetobcswGlo(k,iobcs), k=1,nr),iobcs= 1,nobcs)
190     #endif
191     #ifdef ALLOW_OBCSE_CONTROL
192     write(cunit) ((nWetobcseGlo(k,iobcs), k=1,nr),iobcs= 1,nobcs)
193     #endif
194 mlosch 1.15 write(cunit) (ncvarindex(i), i=1,maxcvars)
195     write(cunit) (ncvarrecs(i), i=1,maxcvars)
196     write(cunit) (ncvarxmax(i), i=1,maxcvars)
197     write(cunit) (ncvarymax(i), i=1,maxcvars)
198     write(cunit) (ncvarnrmax(i), i=1,maxcvars)
199     write(cunit) (ncvargrd(i), i=1,maxcvars)
200 heimbach 1.1 write(cunit)
201    
202     #ifdef ALLOW_THETA0_CONTROL
203 heimbach 1.5 ivartype = 1
204 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
205     write(weighttype(1:80),'(a)') "wtheta"
206 heimbach 1.5 call ctrl_set_pack_xyz(
207 heimbach 1.11 & cunit, ivartype, fname_theta(ictrlgrad), "hFacC",
208 heimbach 1.8 & weighttype, wtheta, lxxadxx, mythid)
209 heimbach 1.1 #endif
210    
211     #ifdef ALLOW_SALT0_CONTROL
212 heimbach 1.5 ivartype = 2
213 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
214     write(weighttype(1:80),'(a)') "wsalt"
215 heimbach 1.5 call ctrl_set_pack_xyz(
216 heimbach 1.11 & cunit, ivartype, fname_salt(ictrlgrad), "hFacC",
217 heimbach 1.8 & weighttype, wsalt, lxxadxx, mythid)
218 heimbach 1.5 #endif
219    
220     #if (defined (ALLOW_HFLUX_CONTROL) || \
221     defined (ALLOW_HFLUX0_CONTROL))
222     ivartype = 3
223     write(weighttype(1:80),'(80a)') ' '
224     write(weighttype(1:80),'(a)') "whflux"
225     call ctrl_set_pack_xy(
226 heimbach 1.11 & cunit, ivartype, fname_hflux(ictrlgrad), "hFacC",
227     & weighttype, lxxadxx, mythid)
228 heimbach 1.5 #endif
229    
230     #if (defined (ALLOW_SFLUX_CONTROL) || \
231     defined (ALLOW_SFLUX0_CONTROL))
232     ivartype = 4
233     write(weighttype(1:80),'(80a)') ' '
234     write(weighttype(1:80),'(a)') "wsflux"
235     call ctrl_set_pack_xy(
236 heimbach 1.11 & cunit, ivartype, fname_sflux(ictrlgrad), "hFacC",
237     & weighttype, lxxadxx, mythid)
238 heimbach 1.5 #endif
239    
240     #if (defined (ALLOW_USTRESS_CONTROL) || \
241     defined (ALLOW_TAUU0_CONTROL))
242     ivartype = 5
243     write(weighttype(1:80),'(80a)') ' '
244     write(weighttype(1:80),'(a)') "wtauu"
245     call ctrl_set_pack_xy(
246 heimbach 1.11 & cunit, ivartype, fname_tauu(ictrlgrad), "maskW",
247     & weighttype, lxxadxx, mythid)
248 heimbach 1.5 #endif
249    
250     #if (defined (ALLOW_VSTRESS_CONTROL) || \
251     defined (ALLOW_TAUV0_CONTROL))
252     ivartype = 6
253     write(weighttype(1:80),'(80a)') ' '
254     write(weighttype(1:80),'(a)') "wtauv"
255     call ctrl_set_pack_xy(
256 heimbach 1.11 & cunit, ivartype, fname_tauv(ictrlgrad), "maskS",
257     & weighttype, lxxadxx, mythid)
258 heimbach 1.5 #endif
259    
260     #ifdef ALLOW_ATEMP_CONTROL
261     ivartype = 7
262     write(weighttype(1:80),'(80a)') ' '
263     write(weighttype(1:80),'(a)') "watemp"
264     call ctrl_set_pack_xy(
265 heimbach 1.11 & cunit, ivartype, fname_atemp(ictrlgrad), "hFacC",
266     & weighttype, lxxadxx, mythid)
267 heimbach 1.5 #endif
268    
269     #ifdef ALLOW_AQH_CONTROL
270     ivartype = 8
271     write(weighttype(1:80),'(80a)') ' '
272     write(weighttype(1:80),'(a)') "waqh"
273     call ctrl_set_pack_xy(
274 heimbach 1.11 & cunit, ivartype, fname_aqh(ictrlgrad), "hFacC",
275     & weighttype, lxxadxx, mythid)
276 heimbach 1.5 #endif
277    
278     #ifdef ALLOW_UWIND_CONTROL
279     ivartype = 9
280     write(weighttype(1:80),'(80a)') ' '
281     write(weighttype(1:80),'(a)') "wuwind"
282     call ctrl_set_pack_xy(
283 heimbach 1.11 & cunit, ivartype, fname_uwind(ictrlgrad), "maskW",
284     & weighttype, lxxadxx, mythid)
285 heimbach 1.5 #endif
286    
287     #ifdef ALLOW_VWIND_CONTROL
288     ivartype = 10
289     write(weighttype(1:80),'(80a)') ' '
290     write(weighttype(1:80),'(a)') "wvwind"
291     call ctrl_set_pack_xy(
292 heimbach 1.11 & cunit, ivartype, fname_vwind(ictrlgrad), "maskS",
293     & weighttype, lxxadxx, mythid)
294 heimbach 1.5 #endif
295    
296     #ifdef ALLOW_OBCSN_CONTROL
297     ivartype = 11
298 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
299     write(weighttype(1:80),'(a)') "wobcsn"
300 heimbach 1.5 call ctrl_set_pack_xz(
301 heimbach 1.11 & cunit, ivartype, fname_obcsn(ictrlgrad), "maskobcsn",
302 heimbach 1.8 & weighttype, wobcsn, lxxadxx, mythid)
303 heimbach 1.5 #endif
304    
305     #ifdef ALLOW_OBCSS_CONTROL
306     ivartype = 12
307 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
308     write(weighttype(1:80),'(a)') "wobcss"
309 heimbach 1.5 call ctrl_set_pack_xz(
310 heimbach 1.11 & cunit, ivartype, fname_obcss(ictrlgrad), "maskobcss",
311 heimbach 1.8 & weighttype, wobcss, lxxadxx, mythid)
312 heimbach 1.5 #endif
313    
314     #ifdef ALLOW_OBCSW_CONTROL
315     ivartype = 13
316 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
317     write(weighttype(1:80),'(a)') "wobcsw"
318 heimbach 1.5 call ctrl_set_pack_yz(
319 heimbach 1.11 & cunit, ivartype, fname_obcsw(ictrlgrad), "maskobcsw",
320 heimbach 1.8 & weighttype, wobcsw, lxxadxx, mythid)
321 heimbach 1.5 #endif
322    
323     #ifdef ALLOW_OBCSE_CONTROL
324     ivartype = 14
325 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
326     write(weighttype(1:80),'(a)') "wobcse"
327 heimbach 1.5 call ctrl_set_pack_yz(
328 heimbach 1.11 & cunit, ivartype, fname_obcse(ictrlgrad), "maskobcse",
329 heimbach 1.8 & weighttype, wobcse, lxxadxx, mythid)
330 heimbach 1.1 #endif
331 heimbach 1.3
332     #ifdef ALLOW_DIFFKR_CONTROL
333 heimbach 1.5 ivartype = 15
334 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
335     write(weighttype(1:80),'(a)') "wdiffkr"
336 heimbach 1.5 call ctrl_set_pack_xyz(
337 heimbach 1.11 & cunit, ivartype, fname_diffkr(ictrlgrad), "hFacC",
338 heimbach 1.8 & weighttype, wunit, lxxadxx, mythid)
339 heimbach 1.3 #endif
340    
341     #ifdef ALLOW_KAPGM_CONTROL
342 heimbach 1.5 ivartype = 16
343 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
344     write(weighttype(1:80),'(a)') "wkapgm"
345 heimbach 1.5 call ctrl_set_pack_xyz(
346 heimbach 1.11 & cunit, ivartype, fname_kapgm(ictrlgrad), "hFacC",
347 heimbach 1.8 & weighttype, wunit, lxxadxx, mythid)
348 heimbach 1.3 #endif
349    
350 heimbach 1.5 #ifdef ALLOW_TR10_CONTROL
351     ivartype = 17
352 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
353     write(weighttype(1:80),'(a)') "wtr1"
354 heimbach 1.5 call ctrl_set_pack_xyz(
355 heimbach 1.11 & cunit, ivartype, fname_tr1(ictrlgrad), "hFacC",
356 heimbach 1.8 & weighttype, wunit, lxxadxx, mythid)
357 heimbach 1.5 #endif
358    
359 heimbach 1.6 #ifdef ALLOW_SST0_CONTROL
360     ivartype = 18
361     write(weighttype(1:80),'(80a)') ' '
362     write(weighttype(1:80),'(a)') "wsst0"
363     call ctrl_set_pack_xy(
364 heimbach 1.11 & cunit, ivartype, fname_sst(ictrlgrad), "hFacC",
365     & weighttype, lxxadxx, mythid)
366 heimbach 1.6 #endif
367    
368     #ifdef ALLOW_SSS0_CONTROL
369     ivartype = 19
370     write(weighttype(1:80),'(80a)') ' '
371     write(weighttype(1:80),'(a)') "wsss0"
372     call ctrl_set_pack_xy(
373 heimbach 1.11 & cunit, ivartype, fname_sss(ictrlgrad), "hFacC",
374     & weighttype, lxxadxx, mythid)
375 heimbach 1.6 #endif
376    
377     #ifdef ALLOW_HFACC_CONTROL
378     ivartype = 20
379 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
380     write(weighttype(1:80),'(a)') "whfacc"
381     # ifdef ALLOW_HFACC3D_CONTROL
382 heimbach 1.6 call ctrl_set_pack_xyz(
383 heimbach 1.11 & cunit, ivartype, fname_hfacc(ictrlgrad), "hFacC",
384 heimbach 1.8 & weighttype, wunit, lxxadxx, mythid)
385     # else
386 heimbach 1.6 call ctrl_set_pack_xy(
387 heimbach 1.11 & cunit, ivartype, fname_hfacc(ictrlgrad), "hFacC",
388     & weighttype, lxxadxx, mythid)
389 heimbach 1.8 # endif
390 heimbach 1.6 #endif
391    
392 heimbach 1.5 #ifdef ALLOW_EFLUXY0_CONTROL
393     ivartype = 21
394 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
395     write(weighttype(1:80),'(a)') "wefluxy0"
396 heimbach 1.5 call ctrl_set_pack_xyz(
397 heimbach 1.11 & cunit, ivartype, fname_efluxy(ictrlgrad), "hFacS",
398 heimbach 1.8 & weighttype, wunit, lxxadxx, mythid)
399 heimbach 1.5 #endif
400    
401     #ifdef ALLOW_EFLUXP0_CONTROL
402     ivartype = 22
403 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
404     write(weighttype(1:80),'(a)') "wefluxp0"
405 heimbach 1.5 call ctrl_set_pack_xyz(
406 heimbach 1.11 & cunit, ivartype, fname_efluxp(ictrlgrad), "hFacV",
407 heimbach 1.8 & weighttype, wunit, lxxadxx, mythid)
408 heimbach 1.6 #endif
409    
410     #ifdef ALLOW_BOTTOMDRAG_CONTROL
411     ivartype = 23
412     write(weighttype(1:80),'(80a)') ' '
413     write(weighttype(1:80),'(a)') "wbottomdrag"
414     call ctrl_set_pack_xy(
415 heimbach 1.11 & cunit, ivartype, fname_bottomdrag(ictrlgrad), "hFacC",
416     & weighttype, lxxadxx, mythid)
417 heimbach 1.5 #endif
418    
419 heimbach 1.18 #ifdef ALLOW_EDTAUX_CONTROL
420     ivartype = 25
421     write(weighttype(1:80),'(80a)') ' '
422     write(weighttype(1:80),'(a)') "wedtaux"
423     call ctrl_set_pack_xyz(
424     & cunit, ivartype, fname_edtaux(ictrlgrad), "maskW",
425     & weighttype, wunit, lxxadxx, mythid)
426     #endif
427    
428     #ifdef ALLOW_EDTAUY_CONTROL
429     ivartype = 26
430     write(weighttype(1:80),'(80a)') ' '
431     write(weighttype(1:80),'(a)') "wedtauy"
432     call ctrl_set_pack_xyz(
433     & cunit, ivartype, fname_edtauy(ictrlgrad), "maskS",
434     & weighttype, wunit, lxxadxx, mythid)
435     #endif
436    
437    
438 heimbach 1.5 close ( cunit )
439 heimbach 1.1
440 heimbach 1.5 _END_MASTER( mythid )
441 heimbach 1.11
442     #endif /* EXCLUDE_CTRL_PACK */
443 heimbach 1.1
444     return
445     end
446    

  ViewVC Help
Powered by ViewVC 1.1.22