/[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.12 - (hide annotations) (download)
Thu Mar 4 19:49:47 2004 UTC (20 years, 3 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint52l_pre, hrcube5, checkpoint54d_post, checkpoint54e_post, checkpoint52l_post, checkpoint54, checkpoint53, checkpoint54f_post, checkpoint53d_post, checkpoint54b_post, checkpoint52m_post, checkpoint54a_pre, checkpoint53c_post, checkpoint54a_post, checkpoint53a_post, checkpoint53g_post, checkpoint53f_post, checkpoint52n_post, checkpoint53b_pre, checkpoint53b_post, checkpoint53d_pre, checkpoint54c_post
Changes since 1.11: +21 -3 lines
Some tricks...

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

  ViewVC Help
Powered by ViewVC 1.1.22