/[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.13 - (hide annotations) (download)
Fri Sep 17 23:02:01 2004 UTC (20 years, 9 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint55, checkpoint55i_post, checkpoint55c_post, checkpoint55g_post, checkpoint55d_post, checkpoint55d_pre, checkpoint55j_post, checkpoint55h_post, checkpoint55b_post, checkpoint55f_post, checkpoint55a_post, checkpoint55e_post
Changes since 1.12: +11 -7 lines
o bringing adjoint up to date for sheduled c55

1 heimbach 1.13 C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_map_ini.F,v 1.12 2004/03/04 19:49:47 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.1
67 heimbach 1.5 _RL fac
68 heimbach 1.12 _RL tmptest
69 heimbach 1.5
70 heimbach 1.1 c == external ==
71     integer ilnblnk
72     external ilnblnk
73    
74     c == end of interface ==
75 heimbach 1.5 CEOP
76 heimbach 1.1
77     jtlo = mybylo(mythid)
78     jthi = mybyhi(mythid)
79     itlo = mybxlo(mythid)
80     ithi = mybxhi(mythid)
81 heimbach 1.8 jmin = 1
82     jmax = sny
83     imin = 1
84     imax = snx
85 heimbach 1.1
86     doglobalread = .false.
87     ladinit = .false.
88    
89     equal = .true.
90    
91     if ( equal ) then
92     fac = 1. _d 0
93     else
94     fac = 0. _d 0
95     endif
96    
97     #ifdef ALLOW_THETA0_CONTROL
98     c-- Temperature field.
99     il=ilnblnk( xx_theta_file )
100     write(fnametheta(1:80),'(2a,i10.10)')
101     & xx_theta_file(1:il),'.',optimcycle
102 heimbach 1.11 call active_read_xyz_loc( fnametheta, tmpfld3d, 1,
103 heimbach 1.1 & doglobalread, ladinit, optimcycle,
104     & mythid, xx_theta_dummy )
105    
106     do bj = jtlo,jthi
107     do bi = itlo,ithi
108     do k = 1,nr
109     do j = jmin,jmax
110     do i = imin,imax
111 heimbach 1.12 #ifdef ALLOW_ECCO
112     IF (abs(tmpfld3d(i,j,k,bi,bj)).gt.
113     $ 2.0/sqrt(wtheta(k,bi,bj)))
114     $ tmpfld3d(i,j,k,bi,bj)=
115     $ sign(2.0/sqrt(wtheta(k,bi,bj)),tmpfld3d(i,j,k,bi,bj))
116     #endif
117 heimbach 1.1 theta(i,j,k,bi,bj) = theta(i,j,k,bi,bj) +
118     & fac*tmpfld3d(i,j,k,bi,bj)
119 heimbach 1.9 if(theta(i,j,k,bi,bj).lt.-2.0)
120 heimbach 1.12 & theta(i,j,k,bi,bj)= -2.0
121 heimbach 1.1 enddo
122     enddo
123     enddo
124     enddo
125     enddo
126 heimbach 1.12
127 heimbach 1.1 #endif
128    
129     #ifdef ALLOW_SALT0_CONTROL
130     c-- Temperature field.
131     il=ilnblnk( xx_salt_file )
132     write(fnamesalt(1:80),'(2a,i10.10)')
133     & xx_salt_file(1:il),'.',optimcycle
134 heimbach 1.11 call active_read_xyz_loc( fnamesalt, tmpfld3d, 1,
135 heimbach 1.1 & doglobalread, ladinit, optimcycle,
136     & mythid, xx_salt_dummy )
137    
138     do bj = jtlo,jthi
139     do bi = itlo,ithi
140     do k = 1,nr
141     do j = jmin,jmax
142     do i = imin,imax
143 heimbach 1.12 #ifdef ALLOW_ECCO
144     IF (abs(tmpfld3d(i,j,k,bi,bj)).gt.
145     $ 2.0/sqrt(wsalt(k,bi,bj)))
146     $ tmpfld3d(i,j,k,bi,bj)=
147     $ sign(2.0/sqrt(wsalt(k,bi,bj)),tmpfld3d(i,j,k,bi,bj))
148     #endif
149 heimbach 1.1 salt(i,j,k,bi,bj) = salt(i,j,k,bi,bj) +
150     & fac*tmpfld3d(i,j,k,bi,bj)
151 heimbach 1.12
152 heimbach 1.1 enddo
153     enddo
154     enddo
155     enddo
156 heimbach 1.12 enddO
157 heimbach 1.1 #endif
158    
159 heimbach 1.2 #ifdef ALLOW_TR10_CONTROL
160 heimbach 1.13 #ifdef ALLOW_PTRACERS
161 heimbach 1.2 c-- Temperature field.
162     il=ilnblnk( xx_tr1_file )
163     write(fnametr1(1:80),'(2a,i10.10)')
164     & xx_tr1_file(1:il),'.',optimcycle
165 heimbach 1.11 call active_read_xyz_loc( fnametr1, tmpfld3d, 1,
166 heimbach 1.2 & doglobalread, ladinit, optimcycle,
167     & mythid, xx_tr1_dummy )
168    
169     do bj = jtlo,jthi
170     do bi = itlo,ithi
171     do k = 1,nr
172     do j = jmin,jmax
173     do i = imin,imax
174 heimbach 1.13 ptracer(i,j,k,bi,bj,1) = ptracer(i,j,k,bi,bj,1) +
175 heimbach 1.2 & fac*tmpfld3d(i,j,k,bi,bj)
176     enddo
177     enddo
178     enddo
179     enddo
180     enddo
181     #endif
182 heimbach 1.13 #endif
183 heimbach 1.2
184 heimbach 1.3 #ifdef ALLOW_DIFFKR_CONTROL
185     c-- diffkr.
186     il=ilnblnk( xx_diffkr_file )
187     write(fnamediffkr(1:80),'(2a,i10.10)')
188     & xx_diffkr_file(1:il),'.',optimcycle
189 heimbach 1.11 call active_read_xyz_loc( fnamediffkr, tmpfld3d, 1,
190 heimbach 1.3 & doglobalread, ladinit, optimcycle,
191     & mythid, xx_diffkr_dummy )
192     do bj = jtlo,jthi
193     do bi = itlo,ithi
194     do k = 1,nr
195     do j = jmin,jmax
196     do i = imin,imax
197     diffkr(i,j,k,bi,bj) = diffkr(i,j,k,bi,bj) +
198     & tmpfld3d(i,j,k,bi,bj)
199     enddo
200     enddo
201     enddo
202     enddo
203     enddo
204     #endif
205    
206     #ifdef ALLOW_KAPGM_CONTROL
207     c-- kapgm.
208     il=ilnblnk( xx_kapgm_file )
209     write(fnamekapgm(1:80),'(2a,i10.10)')
210     & xx_kapgm_file(1:il),'.',optimcycle
211 heimbach 1.11 call active_read_xyz_loc( fnamekapgm, tmpfld3d, 1,
212 heimbach 1.3 & doglobalread, ladinit, optimcycle,
213     & mythid, xx_kapgm_dummy )
214     do bj = jtlo,jthi
215     do bi = itlo,ithi
216     do k = 1,nr
217     do j = jmin,jmax
218     do i = imin,imax
219     kapgm(i,j,k,bi,bj) = kapgm(i,j,k,bi,bj) +
220     & tmpfld3d(i,j,k,bi,bj)
221     enddo
222     enddo
223     enddo
224     enddo
225     enddo
226     #endif
227    
228 heimbach 1.6 #ifdef ALLOW_EFLUXY0_CONTROL
229     c-- y-component EP-flux field.
230     il=ilnblnk( xx_efluxy_file )
231     write(fnameefluxy(1:80),'(2a,i10.10)')
232     & xx_efluxy_file(1:il),'.',optimcycle
233 heimbach 1.11 call active_read_xyz_loc( fnameefluxy, tmpfld3d, 1,
234 heimbach 1.6 & doglobalread, ladinit, optimcycle,
235     & mythid, xx_efluxy_dummy )
236    
237     do bj = jtlo,jthi
238     do bi = itlo,ithi
239     do k = 1,nr
240     do j = jmin,jmax
241     do i = imin,imax
242     EfluxY(i,j,k,bi,bj) = EfluxY(i,j,k,bi,bj)
243     & - fac*tmpfld3d(i,j,k,bi,bj)
244     & *maskS(i,j,k,bi,bj)
245     cph EfluxY(i,j,k,bi,bj) = EfluxY(i,j,k,bi,bj)
246     cph & - rSphere*cosFacU(J,bi,bj)
247     cph & *fac*tmpfld3d(i,j,k,bi,bj)
248     enddo
249     enddo
250     enddo
251     enddo
252     enddo
253     #endif
254    
255     #ifdef ALLOW_EFLUXP0_CONTROL
256     c-- p-component EP-flux field.
257     il=ilnblnk( xx_efluxp_file )
258     write(fnameefluxp(1:80),'(2a,i10.10)')
259     & xx_efluxp_file(1:il),'.',optimcycle
260 heimbach 1.11 call active_read_xyz_loc( fnameefluxp, tmpfld3d, 1,
261 heimbach 1.6 & doglobalread, ladinit, optimcycle,
262     & mythid, xx_efluxp_dummy )
263    
264     do bj = jtlo,jthi
265     do bi = itlo,ithi
266     do k = 1,nr
267     do j = jmin,jmax
268     do i = imin,imax
269     EfluxP(i,j,k,bi,bj) = EfluxP(i,j,k,bi,bj)
270     & + fCori(i,j,bi,bj)
271     & *fac*tmpfld3d(i,j,k,bi,bj)
272     & *hFacV(i,j,k,bi,bj)
273     cph EfluxP(i,j,k,bi,bj) = EfluxP(i,j,k,bi,bj)
274     cph & + fCori(i,j,bi,bj)
275     cph & *rSphere*cosFacU(J,bi,bj)
276     cph & *fac*tmpfld3d(i,j,k,bi,bj)
277     enddo
278     enddo
279     enddo
280     enddo
281     enddo
282     #endif
283    
284 heimbach 1.7 #ifdef ALLOW_BOTTOMDRAG_CONTROL
285     c-- bottom drag
286     il=ilnblnk( xx_bottomdrag_file )
287     write(fnamebottomdrag(1:80),'(2a,i10.10)')
288     & xx_bottomdrag_file(1:il),'.',optimcycle
289 heimbach 1.11 call active_read_xy_loc ( fnamebottomdrag, tmpfld2d, 1,
290 heimbach 1.7 & doglobalread, ladinit, optimcycle,
291     & mythid, xx_bottomdrag_dummy )
292     do bj = jtlo,jthi
293     do bi = itlo,ithi
294     do j = jmin,jmax
295     do i = imin,imax
296     bottomdragfld(i,j,bi,bj) = bottomdragfld(i,j,bi,bj)
297     & + tmpfld2d(i,j,bi,bj)
298     enddo
299     enddo
300     enddo
301     enddo
302     #endif
303    
304 heimbach 1.1
305     c-- Update the tile edges.
306    
307     #ifdef ALLOW_THETA0_CONTROL
308     _EXCH_XYZ_R8( theta, mythid )
309     #endif
310     #ifdef ALLOW_SALT0_CONTROL
311     _EXCH_XYZ_R8( salt, mythid )
312 heimbach 1.2 #endif
313     #ifdef ALLOW_TR10_CONTROL
314 heimbach 1.13 #ifdef ALLOW_PTRACERS
315     _EXCH_XYZ_R8(pTracer(1-Olx,1-Oly,1,1,1,1),myThid)
316     #endif
317 heimbach 1.1 #endif
318 heimbach 1.3 #ifdef ALLOW_DIFFKR_CONTROL
319     _EXCH_XYZ_R8( diffkr, mythid)
320     #endif
321     #ifdef ALLOW_KAPGM_CONTROL
322     _EXCH_XYZ_R8( kapgm, mythid)
323 heimbach 1.6 #endif
324     #ifdef ALLOW_EFLUXY0_CONTROL
325     _EXCH_XYZ_R8( EfluxY, mythid )
326     #endif
327     #ifdef ALLOW_EFLUXP0_CONTROL
328     _EXCH_XYZ_R8( EfluxP, mythid )
329 heimbach 1.7 #endif
330     #ifdef ALLOW_BOTTOMDRAG_CONTROL
331     _EXCH_XY_R8( bottomdragfld, mythid )
332 heimbach 1.3 #endif
333    
334 heimbach 1.1
335     return
336     end
337    

  ViewVC Help
Powered by ViewVC 1.1.22