/[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.9 - (hide annotations) (download)
Tue Jun 24 16:07:06 2003 UTC (20 years, 11 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint51, checkpoint51b_pre, checkpoint51a_post
Changes since 1.8: +3 -7 lines
Merging for c51 vs. e34

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

  ViewVC Help
Powered by ViewVC 1.1.22