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

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

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


Revision 1.14 - (hide annotations) (download)
Tue Nov 16 05:42:12 2004 UTC (19 years, 10 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint57b_post, checkpoint56b_post, checkpoint57d_post, checkpoint57, checkpoint56, checkpoint57a_post, checkpoint57c_post, checkpoint57c_pre, eckpoint57e_pre, checkpoint56a_post, checkpoint56c_post, checkpoint57a_pre
Changes since 1.13: +47 -3 lines
More on dsvd vs. MITgcm interfacing
o handling of g_, ad, via admtlm_vector (mds...vector)
o use ctrl_pack/unpack for admtlm_vector I/O
o use optimcycle for dsvd iteration
o make sure norm is w.r.t. derived quantities

1 heimbach 1.14 C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_map_ini.F,v 1.13 2004/09/17 23:02:01 heimbach Exp $
2 heimbach 1.1
3     #include "CTRL_CPPOPTIONS.h"
4    
5 heimbach 1.5 CBOP
6     C !ROUTINE: ctrl_map_ini
7     C !INTERFACE:
8     subroutine ctrl_map_ini( mythid )
9    
10     C !DESCRIPTION: \bv
11     c *=================================================================
12     c | SUBROUTINE ctrl_map_ini
13     c | Add the temperature, salinity, and diffusivity parts of the
14     c | control vector to the model state and update the tile halos.
15     c | The control vector is defined in the header file "ctrl.h".
16     c *=================================================================
17     C \ev
18 heimbach 1.1
19 heimbach 1.5 C !USES:
20 heimbach 1.1 implicit none
21    
22     c == global variables ==
23 heimbach 1.6 #include "SIZE.h"
24 heimbach 1.1 #include "EEPARAMS.h"
25 heimbach 1.6 #include "PARAMS.h"
26 heimbach 1.1 #include "DYNVARS.h"
27 heimbach 1.6 #include "GRID.h"
28 heimbach 1.1 #include "ctrl.h"
29     #include "ctrl_dummy.h"
30 heimbach 1.2 #include "optim.h"
31 heimbach 1.13 #ifdef ALLOW_PTRACERS
32     # include "PTRACERS_SIZE.h"
33     # include "PTRACERS.h"
34     #endif
35 heimbach 1.12 #ifdef ALLOW_ECCO
36     # include "ecco_cost.h"
37     #endif
38 heimbach 1.1
39 heimbach 1.5 C !INPUT/OUTPUT PARAMETERS:
40 heimbach 1.1 c == routine arguments ==
41     integer mythid
42    
43 heimbach 1.5 C !LOCAL VARIABLES:
44 heimbach 1.1 c == local variables ==
45    
46     integer bi,bj
47     integer i,j,k
48     integer itlo,ithi
49     integer jtlo,jthi
50     integer jmin,jmax
51     integer imin,imax
52     integer il
53    
54     logical equal
55     logical doglobalread
56     logical ladinit
57    
58     character*( 80) fnametheta
59     character*( 80) fnamesalt
60 heimbach 1.2 character*( 80) fnametr1
61 heimbach 1.3 character*( 80) fnamediffkr
62     character*( 80) fnamekapgm
63 heimbach 1.6 character*( 80) fnameefluxy
64     character*( 80) fnameefluxp
65 heimbach 1.7 character*( 80) fnamebottomdrag
66 heimbach 1.14 character*( 80) fnamesss
67     character*( 80) fnamesst
68 heimbach 1.1
69 heimbach 1.5 _RL fac
70 heimbach 1.12 _RL tmptest
71 heimbach 1.5
72 heimbach 1.1 c == external ==
73     integer ilnblnk
74     external ilnblnk
75    
76     c == end of interface ==
77 heimbach 1.5 CEOP
78 heimbach 1.1
79     jtlo = mybylo(mythid)
80     jthi = mybyhi(mythid)
81     itlo = mybxlo(mythid)
82     ithi = mybxhi(mythid)
83 heimbach 1.8 jmin = 1
84     jmax = sny
85     imin = 1
86     imax = snx
87 heimbach 1.1
88     doglobalread = .false.
89     ladinit = .false.
90    
91     equal = .true.
92    
93     if ( equal ) then
94     fac = 1. _d 0
95     else
96     fac = 0. _d 0
97     endif
98    
99     #ifdef ALLOW_THETA0_CONTROL
100     c-- Temperature field.
101     il=ilnblnk( xx_theta_file )
102     write(fnametheta(1:80),'(2a,i10.10)')
103     & xx_theta_file(1:il),'.',optimcycle
104 heimbach 1.11 call active_read_xyz_loc( fnametheta, tmpfld3d, 1,
105 heimbach 1.1 & doglobalread, ladinit, optimcycle,
106     & mythid, xx_theta_dummy )
107    
108     do bj = jtlo,jthi
109     do bi = itlo,ithi
110     do k = 1,nr
111     do j = jmin,jmax
112     do i = imin,imax
113 heimbach 1.12 #ifdef ALLOW_ECCO
114     IF (abs(tmpfld3d(i,j,k,bi,bj)).gt.
115     $ 2.0/sqrt(wtheta(k,bi,bj)))
116     $ tmpfld3d(i,j,k,bi,bj)=
117     $ sign(2.0/sqrt(wtheta(k,bi,bj)),tmpfld3d(i,j,k,bi,bj))
118     #endif
119 heimbach 1.1 theta(i,j,k,bi,bj) = theta(i,j,k,bi,bj) +
120     & fac*tmpfld3d(i,j,k,bi,bj)
121 heimbach 1.9 if(theta(i,j,k,bi,bj).lt.-2.0)
122 heimbach 1.12 & theta(i,j,k,bi,bj)= -2.0
123 heimbach 1.1 enddo
124     enddo
125     enddo
126     enddo
127     enddo
128 heimbach 1.12
129 heimbach 1.1 #endif
130    
131     #ifdef ALLOW_SALT0_CONTROL
132     c-- Temperature field.
133     il=ilnblnk( xx_salt_file )
134     write(fnamesalt(1:80),'(2a,i10.10)')
135     & xx_salt_file(1:il),'.',optimcycle
136 heimbach 1.11 call active_read_xyz_loc( fnamesalt, tmpfld3d, 1,
137 heimbach 1.1 & doglobalread, ladinit, optimcycle,
138     & mythid, xx_salt_dummy )
139    
140     do bj = jtlo,jthi
141     do bi = itlo,ithi
142     do k = 1,nr
143     do j = jmin,jmax
144     do i = imin,imax
145 heimbach 1.12 #ifdef ALLOW_ECCO
146     IF (abs(tmpfld3d(i,j,k,bi,bj)).gt.
147     $ 2.0/sqrt(wsalt(k,bi,bj)))
148     $ tmpfld3d(i,j,k,bi,bj)=
149     $ sign(2.0/sqrt(wsalt(k,bi,bj)),tmpfld3d(i,j,k,bi,bj))
150     #endif
151 heimbach 1.1 salt(i,j,k,bi,bj) = salt(i,j,k,bi,bj) +
152     & fac*tmpfld3d(i,j,k,bi,bj)
153 heimbach 1.12
154 heimbach 1.1 enddo
155     enddo
156     enddo
157     enddo
158 heimbach 1.12 enddO
159 heimbach 1.1 #endif
160    
161 heimbach 1.2 #ifdef ALLOW_TR10_CONTROL
162 heimbach 1.13 #ifdef ALLOW_PTRACERS
163 heimbach 1.2 c-- Temperature field.
164     il=ilnblnk( xx_tr1_file )
165     write(fnametr1(1:80),'(2a,i10.10)')
166     & xx_tr1_file(1:il),'.',optimcycle
167 heimbach 1.11 call active_read_xyz_loc( fnametr1, tmpfld3d, 1,
168 heimbach 1.2 & doglobalread, ladinit, optimcycle,
169     & mythid, xx_tr1_dummy )
170    
171     do bj = jtlo,jthi
172     do bi = itlo,ithi
173     do k = 1,nr
174     do j = jmin,jmax
175     do i = imin,imax
176 heimbach 1.13 ptracer(i,j,k,bi,bj,1) = ptracer(i,j,k,bi,bj,1) +
177 heimbach 1.2 & fac*tmpfld3d(i,j,k,bi,bj)
178     enddo
179     enddo
180     enddo
181     enddo
182     enddo
183     #endif
184 heimbach 1.13 #endif
185 heimbach 1.2
186 heimbach 1.14 #ifdef ALLOW_SST0_CONTROL
187     c-- sst0.
188     il=ilnblnk( xx_sst_file )
189     write(fnamesst(1:80),'(2a,i10.10)')
190     & xx_sst_file(1:il),'.',optimcycle
191     call active_read_xy_loc ( fnamesst, tmpfld2d, 1,
192     & doglobalread, ladinit, optimcycle,
193     & mythid, xx_sst_dummy )
194     do bj = jtlo,jthi
195     do bi = itlo,ithi
196     do j = jmin,jmax
197     do i = imin,imax
198     cph sst(i,j,bi,bj) = sst(i,j,bi,bj) + tmpfld2d(i,j,bi,bj)
199     theta(i,j,1,bi,bj) = theta(i,j,1,bi,bj)
200     & + tmpfld2d(i,j,bi,bj)
201     enddo
202     enddo
203     enddo
204     enddo
205     #endif
206    
207     #ifdef ALLOW_SSS0_CONTROL
208     c-- sss0.
209     il=ilnblnk( xx_sss_file )
210     write(fnamesss(1:80),'(2a,i10.10)')
211     & xx_sss_file(1:il),'.',optimcycle
212     call active_read_xy_loc ( fnamesss, tmpfld2d, 1,
213     & doglobalread, ladinit, optimcycle,
214     & mythid, xx_sss_dummy )
215     do bj = jtlo,jthi
216     do bi = itlo,ithi
217     do j = jmin,jmax
218     do i = imin,imax
219     cph sss(i,j,bi,bj) = sss(i,j,bi,bj) + tmpfld2d(i,j,bi,bj)
220     salt(i,j,1,bi,bj) = salt(i,j,1,bi,bj)
221     & + tmpfld2d(i,j,bi,bj)
222     enddo
223     enddo
224     enddo
225     enddo
226     #endif
227    
228 heimbach 1.3 #ifdef ALLOW_DIFFKR_CONTROL
229     c-- diffkr.
230     il=ilnblnk( xx_diffkr_file )
231     write(fnamediffkr(1:80),'(2a,i10.10)')
232     & xx_diffkr_file(1:il),'.',optimcycle
233 heimbach 1.11 call active_read_xyz_loc( fnamediffkr, tmpfld3d, 1,
234 heimbach 1.3 & doglobalread, ladinit, optimcycle,
235     & mythid, xx_diffkr_dummy )
236     do bj = jtlo,jthi
237     do bi = itlo,ithi
238     do k = 1,nr
239     do j = jmin,jmax
240     do i = imin,imax
241     diffkr(i,j,k,bi,bj) = diffkr(i,j,k,bi,bj) +
242     & tmpfld3d(i,j,k,bi,bj)
243     enddo
244     enddo
245     enddo
246     enddo
247     enddo
248     #endif
249    
250     #ifdef ALLOW_KAPGM_CONTROL
251     c-- kapgm.
252     il=ilnblnk( xx_kapgm_file )
253     write(fnamekapgm(1:80),'(2a,i10.10)')
254     & xx_kapgm_file(1:il),'.',optimcycle
255 heimbach 1.11 call active_read_xyz_loc( fnamekapgm, tmpfld3d, 1,
256 heimbach 1.3 & doglobalread, ladinit, optimcycle,
257     & mythid, xx_kapgm_dummy )
258     do bj = jtlo,jthi
259     do bi = itlo,ithi
260     do k = 1,nr
261     do j = jmin,jmax
262     do i = imin,imax
263     kapgm(i,j,k,bi,bj) = kapgm(i,j,k,bi,bj) +
264     & tmpfld3d(i,j,k,bi,bj)
265     enddo
266     enddo
267     enddo
268     enddo
269     enddo
270     #endif
271    
272 heimbach 1.6 #ifdef ALLOW_EFLUXY0_CONTROL
273     c-- y-component EP-flux field.
274     il=ilnblnk( xx_efluxy_file )
275     write(fnameefluxy(1:80),'(2a,i10.10)')
276     & xx_efluxy_file(1:il),'.',optimcycle
277 heimbach 1.11 call active_read_xyz_loc( fnameefluxy, tmpfld3d, 1,
278 heimbach 1.6 & doglobalread, ladinit, optimcycle,
279     & mythid, xx_efluxy_dummy )
280    
281     do bj = jtlo,jthi
282     do bi = itlo,ithi
283     do k = 1,nr
284     do j = jmin,jmax
285     do i = imin,imax
286     EfluxY(i,j,k,bi,bj) = EfluxY(i,j,k,bi,bj)
287     & - fac*tmpfld3d(i,j,k,bi,bj)
288     & *maskS(i,j,k,bi,bj)
289     cph EfluxY(i,j,k,bi,bj) = EfluxY(i,j,k,bi,bj)
290     cph & - rSphere*cosFacU(J,bi,bj)
291     cph & *fac*tmpfld3d(i,j,k,bi,bj)
292     enddo
293     enddo
294     enddo
295     enddo
296     enddo
297     #endif
298    
299     #ifdef ALLOW_EFLUXP0_CONTROL
300     c-- p-component EP-flux field.
301     il=ilnblnk( xx_efluxp_file )
302     write(fnameefluxp(1:80),'(2a,i10.10)')
303     & xx_efluxp_file(1:il),'.',optimcycle
304 heimbach 1.11 call active_read_xyz_loc( fnameefluxp, tmpfld3d, 1,
305 heimbach 1.6 & doglobalread, ladinit, optimcycle,
306     & mythid, xx_efluxp_dummy )
307    
308     do bj = jtlo,jthi
309     do bi = itlo,ithi
310     do k = 1,nr
311     do j = jmin,jmax
312     do i = imin,imax
313     EfluxP(i,j,k,bi,bj) = EfluxP(i,j,k,bi,bj)
314     & + fCori(i,j,bi,bj)
315     & *fac*tmpfld3d(i,j,k,bi,bj)
316     & *hFacV(i,j,k,bi,bj)
317     cph EfluxP(i,j,k,bi,bj) = EfluxP(i,j,k,bi,bj)
318     cph & + fCori(i,j,bi,bj)
319     cph & *rSphere*cosFacU(J,bi,bj)
320     cph & *fac*tmpfld3d(i,j,k,bi,bj)
321     enddo
322     enddo
323     enddo
324     enddo
325     enddo
326     #endif
327    
328 heimbach 1.7 #ifdef ALLOW_BOTTOMDRAG_CONTROL
329     c-- bottom drag
330     il=ilnblnk( xx_bottomdrag_file )
331     write(fnamebottomdrag(1:80),'(2a,i10.10)')
332     & xx_bottomdrag_file(1:il),'.',optimcycle
333 heimbach 1.11 call active_read_xy_loc ( fnamebottomdrag, tmpfld2d, 1,
334 heimbach 1.7 & doglobalread, ladinit, optimcycle,
335     & mythid, xx_bottomdrag_dummy )
336     do bj = jtlo,jthi
337     do bi = itlo,ithi
338     do j = jmin,jmax
339     do i = imin,imax
340     bottomdragfld(i,j,bi,bj) = bottomdragfld(i,j,bi,bj)
341     & + tmpfld2d(i,j,bi,bj)
342     enddo
343     enddo
344     enddo
345     enddo
346     #endif
347    
348 heimbach 1.1
349     c-- Update the tile edges.
350    
351 heimbach 1.14 #if (defined (ALLOW_THETA0_CONTROL) || defined (ALLOW_SST0_CONTROL))
352 heimbach 1.1 _EXCH_XYZ_R8( theta, mythid )
353     #endif
354 heimbach 1.14 #if (defined (ALLOW_SALT0_CONTROL) || defined (ALLOW_SSS0_CONTROL))
355 heimbach 1.1 _EXCH_XYZ_R8( salt, mythid )
356 heimbach 1.2 #endif
357     #ifdef ALLOW_TR10_CONTROL
358 heimbach 1.13 #ifdef ALLOW_PTRACERS
359     _EXCH_XYZ_R8(pTracer(1-Olx,1-Oly,1,1,1,1),myThid)
360     #endif
361 heimbach 1.1 #endif
362 heimbach 1.3 #ifdef ALLOW_DIFFKR_CONTROL
363     _EXCH_XYZ_R8( diffkr, mythid)
364     #endif
365     #ifdef ALLOW_KAPGM_CONTROL
366     _EXCH_XYZ_R8( kapgm, mythid)
367 heimbach 1.6 #endif
368     #ifdef ALLOW_EFLUXY0_CONTROL
369     _EXCH_XYZ_R8( EfluxY, mythid )
370     #endif
371     #ifdef ALLOW_EFLUXP0_CONTROL
372     _EXCH_XYZ_R8( EfluxP, mythid )
373 heimbach 1.7 #endif
374     #ifdef ALLOW_BOTTOMDRAG_CONTROL
375     _EXCH_XY_R8( bottomdragfld, mythid )
376 heimbach 1.3 #endif
377    
378 heimbach 1.1
379     return
380     end
381    

  ViewVC Help
Powered by ViewVC 1.1.22