/[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.7 - (hide annotations) (download)
Fri Nov 29 13:38:37 2002 UTC (22 years, 7 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint47e_post, checkpoint47c_post, checkpoint47d_pre, checkpoint47i_post, checkpoint47d_post, checkpoint47g_post, branch-exfmods-tag, checkpoint47b_post, checkpoint47f_post, checkpoint47h_post
Branch point for: branch-exfmods-curt
Changes since 1.6: +25 -1 lines
Controls of sst, sss, hfacc, bottomdrag.
(no ice climbing).

1 heimbach 1.7 C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_map_ini.F,v 1.6 2002/07/13 02:47:32 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     jmin = 1-oly
75     jmax = sny+oly
76     imin = 1-olx
77     imax = snx+olx
78    
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.4 cph gtNm1(i,j,k,bi,bj) = gtNm1(i,j,k,bi,bj) +
107     cph & fac*tmpfld3d(i,j,k,bi,bj)
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 heimbach 1.4 cph gsNm1(i,j,k,bi,bj) = gsNm1(i,j,k,bi,bj) +
132     cph & fac*tmpfld3d(i,j,k,bi,bj)
133 heimbach 1.1 enddo
134     enddo
135     enddo
136     enddo
137     enddo
138     #endif
139    
140 heimbach 1.2 #ifdef ALLOW_TR10_CONTROL
141     c-- Temperature field.
142     il=ilnblnk( xx_tr1_file )
143     write(fnametr1(1:80),'(2a,i10.10)')
144     & xx_tr1_file(1:il),'.',optimcycle
145     call active_read_xyz( fnametr1, tmpfld3d, 1,
146     & doglobalread, ladinit, optimcycle,
147     & mythid, xx_tr1_dummy )
148    
149     do bj = jtlo,jthi
150     do bi = itlo,ithi
151     do k = 1,nr
152     do j = jmin,jmax
153     do i = imin,imax
154     tr1(i,j,k,bi,bj) = tr1(i,j,k,bi,bj) +
155     & fac*tmpfld3d(i,j,k,bi,bj)
156 heimbach 1.4 cph gtr1Nm1(i,j,k,bi,bj) = gtr1Nm1(i,j,k,bi,bj) +
157     cph & fac*tmpfld3d(i,j,k,bi,bj)
158 heimbach 1.2 enddo
159     enddo
160     enddo
161     enddo
162     enddo
163     #endif
164    
165 heimbach 1.3 #ifdef ALLOW_DIFFKR_CONTROL
166     c-- diffkr.
167     il=ilnblnk( xx_diffkr_file )
168     write(fnamediffkr(1:80),'(2a,i10.10)')
169     & xx_diffkr_file(1:il),'.',optimcycle
170     call active_read_xyz( fnamediffkr, tmpfld3d, 1,
171     & doglobalread, ladinit, optimcycle,
172     & mythid, xx_diffkr_dummy )
173     do bj = jtlo,jthi
174     do bi = itlo,ithi
175     do k = 1,nr
176     do j = jmin,jmax
177     do i = imin,imax
178     diffkr(i,j,k,bi,bj) = diffkr(i,j,k,bi,bj) +
179     & tmpfld3d(i,j,k,bi,bj)
180     enddo
181     enddo
182     enddo
183     enddo
184     enddo
185     #endif
186    
187     #ifdef ALLOW_KAPGM_CONTROL
188     c-- kapgm.
189     il=ilnblnk( xx_kapgm_file )
190     write(fnamekapgm(1:80),'(2a,i10.10)')
191     & xx_kapgm_file(1:il),'.',optimcycle
192     call active_read_xyz( fnamekapgm, tmpfld3d, 1,
193     & doglobalread, ladinit, optimcycle,
194     & mythid, xx_kapgm_dummy )
195     do bj = jtlo,jthi
196     do bi = itlo,ithi
197     do k = 1,nr
198     do j = jmin,jmax
199     do i = imin,imax
200     kapgm(i,j,k,bi,bj) = kapgm(i,j,k,bi,bj) +
201     & tmpfld3d(i,j,k,bi,bj)
202     enddo
203     enddo
204     enddo
205     enddo
206     enddo
207     #endif
208    
209 heimbach 1.6 #ifdef ALLOW_EFLUXY0_CONTROL
210     c-- y-component EP-flux field.
211     il=ilnblnk( xx_efluxy_file )
212     write(fnameefluxy(1:80),'(2a,i10.10)')
213     & xx_efluxy_file(1:il),'.',optimcycle
214     call active_read_xyz( fnameefluxy, tmpfld3d, 1,
215     & doglobalread, ladinit, optimcycle,
216     & mythid, xx_efluxy_dummy )
217    
218     do bj = jtlo,jthi
219     do bi = itlo,ithi
220     do k = 1,nr
221     do j = jmin,jmax
222     do i = imin,imax
223     EfluxY(i,j,k,bi,bj) = EfluxY(i,j,k,bi,bj)
224     & - fac*tmpfld3d(i,j,k,bi,bj)
225     & *maskS(i,j,k,bi,bj)
226     cph EfluxY(i,j,k,bi,bj) = EfluxY(i,j,k,bi,bj)
227     cph & - rSphere*cosFacU(J,bi,bj)
228     cph & *fac*tmpfld3d(i,j,k,bi,bj)
229     enddo
230     enddo
231     enddo
232     enddo
233     enddo
234     #endif
235    
236     #ifdef ALLOW_EFLUXP0_CONTROL
237     c-- p-component EP-flux field.
238     il=ilnblnk( xx_efluxp_file )
239     write(fnameefluxp(1:80),'(2a,i10.10)')
240     & xx_efluxp_file(1:il),'.',optimcycle
241     call active_read_xyz( fnameefluxp, tmpfld3d, 1,
242     & doglobalread, ladinit, optimcycle,
243     & mythid, xx_efluxp_dummy )
244    
245     do bj = jtlo,jthi
246     do bi = itlo,ithi
247     do k = 1,nr
248     do j = jmin,jmax
249     do i = imin,imax
250     EfluxP(i,j,k,bi,bj) = EfluxP(i,j,k,bi,bj)
251     & + fCori(i,j,bi,bj)
252     & *fac*tmpfld3d(i,j,k,bi,bj)
253     & *hFacV(i,j,k,bi,bj)
254     cph EfluxP(i,j,k,bi,bj) = EfluxP(i,j,k,bi,bj)
255     cph & + fCori(i,j,bi,bj)
256     cph & *rSphere*cosFacU(J,bi,bj)
257     cph & *fac*tmpfld3d(i,j,k,bi,bj)
258     enddo
259     enddo
260     enddo
261     enddo
262     enddo
263     #endif
264    
265 heimbach 1.7 #ifdef ALLOW_BOTTOMDRAG_CONTROL
266     c-- bottom drag
267     il=ilnblnk( xx_bottomdrag_file )
268     write(fnamebottomdrag(1:80),'(2a,i10.10)')
269     & xx_bottomdrag_file(1:il),'.',optimcycle
270     call active_read_xy ( fnamebottomdrag, tmpfld2d, 1,
271     & doglobalread, ladinit, optimcycle,
272     & mythid, xx_bottomdrag_dummy )
273     do bj = jtlo,jthi
274     do bi = itlo,ithi
275     do j = jmin,jmax
276     do i = imin,imax
277     bottomdragfld(i,j,bi,bj) = bottomdragfld(i,j,bi,bj)
278     & + tmpfld2d(i,j,bi,bj)
279     enddo
280     enddo
281     enddo
282     enddo
283     #endif
284    
285 heimbach 1.1
286     c-- Update the tile edges.
287    
288     #ifdef ALLOW_THETA0_CONTROL
289     _EXCH_XYZ_R8( theta, mythid )
290 heimbach 1.6 cph _EXCH_XYZ_R8( gtNm1, mythid )
291 heimbach 1.1 #endif
292     #ifdef ALLOW_SALT0_CONTROL
293     _EXCH_XYZ_R8( salt, mythid )
294 heimbach 1.6 cph _EXCH_XYZ_R8( gsNm1, mythid )
295 heimbach 1.2 #endif
296     #ifdef ALLOW_TR10_CONTROL
297 heimbach 1.3 _EXCH_XYZ_R8( tr1, mythid )
298 heimbach 1.6 cph _EXCH_XYZ_R8( gTr1Nm1, mythid )
299 heimbach 1.1 #endif
300 heimbach 1.3 #ifdef ALLOW_DIFFKR_CONTROL
301     _EXCH_XYZ_R8( diffkr, mythid)
302     #endif
303     #ifdef ALLOW_KAPGM_CONTROL
304     _EXCH_XYZ_R8( kapgm, mythid)
305 heimbach 1.6 #endif
306     #ifdef ALLOW_EFLUXY0_CONTROL
307     _EXCH_XYZ_R8( EfluxY, mythid )
308     #endif
309     #ifdef ALLOW_EFLUXP0_CONTROL
310     _EXCH_XYZ_R8( EfluxP, mythid )
311 heimbach 1.7 #endif
312     #ifdef ALLOW_BOTTOMDRAG_CONTROL
313     _EXCH_XY_R8( bottomdragfld, mythid )
314 heimbach 1.3 #endif
315    
316 heimbach 1.1
317     return
318     end
319    

  ViewVC Help
Powered by ViewVC 1.1.22