1 |
C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_map_ini_gen.F,v 1.6 2009/10/25 02:04:52 gforget Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "CTRL_CPPOPTIONS.h" |
5 |
|
6 |
|
7 |
subroutine ctrl_map_ini_gen3D(xxFileCur, wFileCur, xxDummyCur, |
8 |
& boundsVec, paramFld3d, maskFld3d, paramSmooth, mythid ) |
9 |
|
10 |
c ================================================================== |
11 |
c SUBROUTINE ctrl_map_ini_gen3D |
12 |
c ================================================================== |
13 |
c |
14 |
c started: Gael Forget gforget@mit.edu 8-Feb-2008 |
15 |
c |
16 |
c - Generetic routine for an individual 3D control term |
17 |
c (to be called from ctrl_map_ini in a loop e.g.) |
18 |
c |
19 |
c ================================================================== |
20 |
c SUBROUTINE ctrl_map_ini_gen3D |
21 |
c ================================================================== |
22 |
|
23 |
implicit none |
24 |
|
25 |
c == global variables == |
26 |
|
27 |
#include "EEPARAMS.h" |
28 |
#include "SIZE.h" |
29 |
#include "PARAMS.h" |
30 |
#include "GRID.h" |
31 |
#include "DYNVARS.h" |
32 |
#include "FFIELDS.h" |
33 |
|
34 |
#include "ctrl.h" |
35 |
#include "ctrl_dummy.h" |
36 |
#include "optim.h" |
37 |
#ifdef ALLOW_ECCO |
38 |
#include "ecco_cost.h" |
39 |
#endif |
40 |
|
41 |
#ifdef ALLOW_AUTODIFF_TAMC |
42 |
#include "tamc.h" |
43 |
#include "tamc_keys.h" |
44 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
45 |
|
46 |
c == routine arguments == |
47 |
|
48 |
integer mythid |
49 |
character*(*) wFileCur,xxFileCur |
50 |
_RL boundsVec(5),tmpMax,xxDummyCur |
51 |
|
52 |
_RL wFld3d(1-olx:snx+olx,1-oly:sny+oly,nr,nsx,nsy) |
53 |
_RL xxFld3d(1-olx:snx+olx,1-oly:sny+oly,nr,nsx,nsy) |
54 |
_RL paramFld3d(1-olx:snx+olx,1-oly:sny+oly,nr,nsx,nsy) |
55 |
_RL maskFld3d(1-olx:snx+olx,1-oly:sny+oly,nr,nsx,nsy) |
56 |
integer paramSmooth |
57 |
|
58 |
c == local variables == |
59 |
|
60 |
integer bi,bj |
61 |
integer i,j,k |
62 |
integer itlo,ithi |
63 |
integer jtlo,jthi |
64 |
integer jmin,jmax |
65 |
integer imin,imax |
66 |
integer il |
67 |
|
68 |
logical doglobalread |
69 |
logical ladinit |
70 |
|
71 |
character*( 80) fnamegeneric |
72 |
|
73 |
c == external == |
74 |
|
75 |
integer ilnblnk |
76 |
external ilnblnk |
77 |
|
78 |
c == end of interface == |
79 |
|
80 |
#ifdef ALLOW_AUTODIFF_TAMC |
81 |
act3 = myThid - 1 |
82 |
max3 = nTx*nTy |
83 |
act4 = 0 |
84 |
ikey = (act3 + 1) + act4*max3 |
85 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
86 |
|
87 |
jtlo = mybylo(mythid) |
88 |
jthi = mybyhi(mythid) |
89 |
itlo = mybxlo(mythid) |
90 |
ithi = mybxhi(mythid) |
91 |
c-- only do interior, and exchange outside |
92 |
jmin = 1 |
93 |
jmax = sny |
94 |
imin = 1 |
95 |
imax = snx |
96 |
|
97 |
doglobalread = .false. |
98 |
ladinit = .false. |
99 |
|
100 |
|
101 |
call mdsreadfield(wFileCur,32,'RL',nR,wFld3d,1,mythid) |
102 |
_EXCH_XYZ_RL( wFld3d, mythid ) |
103 |
|
104 |
il=ilnblnk( xxFileCur ) |
105 |
write(fnamegeneric(1:80),'(2a,i10.10)') |
106 |
& xxFileCur(1:il),'.',optimcycle |
107 |
call active_read_xyz( fnamegeneric, xxFld3d, 1, |
108 |
& doglobalread, ladinit, optimcycle, mythid, xxDummyCur ) |
109 |
|
110 |
|
111 |
#ifndef ALLOW_SMOOTH_CORREL3D |
112 |
c avoid xx larger than boundsVec(5) X uncertainty |
113 |
if ( boundsVec(5).GT.0.) then |
114 |
do bj = jtlo,jthi |
115 |
do bi = itlo,ithi |
116 |
do k = 1,nr |
117 |
do j = jmin,jmax |
118 |
do i = imin,imax |
119 |
if ( (maskFld3d(i,j,k,bi,bj).NE.0.).AND. |
120 |
& (wFld3d(i,j,k,bi,bj).GT.0.) ) then |
121 |
tmpMax=boundsVec(5)/sqrt(wFld3d(i,j,k,bi,bj)) |
122 |
if ( abs(xxFld3d(i,j,k,bi,bj)).GT.tmpMax ) then |
123 |
xxFld3d(i,j,k,bi,bj)=sign(tmpMax,xxFld3d(i,j,k,bi,bj)) |
124 |
else |
125 |
xxFld3d(i,j,k,bi,bj)=xxFld3d(i,j,k,bi,bj) |
126 |
endif |
127 |
endif |
128 |
enddo |
129 |
enddo |
130 |
enddo |
131 |
enddo |
132 |
enddo |
133 |
endif |
134 |
#else |
135 |
c apply Weaver And Courtier correlation operator |
136 |
if (paramSmooth.NE.0) then |
137 |
call smooth_correl3D(xxFld3d,paramSmooth,mythid) |
138 |
endif |
139 |
#endif |
140 |
|
141 |
do bj = jtlo,jthi |
142 |
do bi = itlo,ithi |
143 |
do k = 1,nr |
144 |
do j = jmin,jmax |
145 |
do i = imin,imax |
146 |
#ifdef ALLOW_SMOOTH_CORREL3D |
147 |
c scale param adjustment |
148 |
if ( (maskFld3d(i,j,k,bi,bj).NE.0.) |
149 |
& .AND. (wFld3d(i,j,k,bi,bj).GT.0.) ) then |
150 |
xxFld3d(i,j,k,bi,bj)=xxFld3d(i,j,k,bi,bj) |
151 |
& /sqrt( wFld3d(i,j,k,bi,bj) ) |
152 |
else |
153 |
xxFld3d(i,j,k,bi,bj)=0. |
154 |
endif |
155 |
#endif |
156 |
paramFld3d(i,j,k,bi,bj) = paramFld3d(i,j,k,bi,bj) |
157 |
& + xxFld3d(i,j,k,bi,bj) |
158 |
enddo |
159 |
enddo |
160 |
enddo |
161 |
enddo |
162 |
enddo |
163 |
|
164 |
c avoid param out of [boundsVec(1) boundsVec(4)] |
165 |
CALL CTRL_BOUND_3D(paramFld3d,maskFld3d,boundsVec,myThid) |
166 |
|
167 |
#ifdef ALLOW_SMOOTH_CORREL3D |
168 |
write(fnamegeneric(1:80),'(2a,i10.10)') |
169 |
& xxFileCur(1:il),'.effective.',optimcycle |
170 |
call mdswritefield(fnamegeneric,ctrlprec,.FALSE.,'RL', |
171 |
& nr, paramFld3d, 1, optimcycle, mythid) |
172 |
#endif |
173 |
|
174 |
|
175 |
end |
176 |
|
177 |
|
178 |
|
179 |
subroutine ctrl_map_ini_gen2D(xxFileCur, wFileCur, xxDummyCur, |
180 |
& boundsVec, paramFld2d, maskFld3d, paramSmooth, mythid ) |
181 |
|
182 |
c ================================================================== |
183 |
c SUBROUTINE ctrl_map_ini_gen2D |
184 |
c ================================================================== |
185 |
c |
186 |
c started: Gael Forget gforget@mit.edu 8-Feb-2008 |
187 |
c |
188 |
c - Generetic routine for an individual 2D control term |
189 |
c (to be called from ctrl_map_ini in a loop e.g.) |
190 |
c |
191 |
c ================================================================== |
192 |
c SUBROUTINE ctrl_map_ini_gen3D |
193 |
c ================================================================== |
194 |
|
195 |
implicit none |
196 |
|
197 |
c == global variables == |
198 |
|
199 |
#include "EEPARAMS.h" |
200 |
#include "SIZE.h" |
201 |
#include "PARAMS.h" |
202 |
#include "GRID.h" |
203 |
#include "DYNVARS.h" |
204 |
#include "FFIELDS.h" |
205 |
|
206 |
#include "ctrl.h" |
207 |
#include "ctrl_dummy.h" |
208 |
#include "optim.h" |
209 |
#ifdef ALLOW_ECCO |
210 |
#include "ecco_cost.h" |
211 |
#endif |
212 |
|
213 |
#ifdef ALLOW_AUTODIFF_TAMC |
214 |
#include "tamc.h" |
215 |
#include "tamc_keys.h" |
216 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
217 |
|
218 |
c == routine arguments == |
219 |
|
220 |
integer mythid |
221 |
character*(*) wFileCur,xxFileCur |
222 |
_RL boundsVec(5),tmpMax,xxDummyCur |
223 |
|
224 |
_RL wFld2d(1-olx:snx+olx,1-oly:sny+oly,nsx,nsy) |
225 |
_RL xxFld2d(1-olx:snx+olx,1-oly:sny+oly,nsx,nsy) |
226 |
_RL paramFld2d(1-olx:snx+olx,1-oly:sny+oly,nsx,nsy) |
227 |
_RL maskFld3d(1-olx:snx+olx,1-oly:sny+oly,nr,nsx,nsy) |
228 |
integer paramSmooth |
229 |
|
230 |
c == local variables == |
231 |
|
232 |
integer bi,bj |
233 |
integer i,j,k |
234 |
integer itlo,ithi |
235 |
integer jtlo,jthi |
236 |
integer jmin,jmax |
237 |
integer imin,imax |
238 |
integer il |
239 |
|
240 |
logical doglobalread |
241 |
logical ladinit |
242 |
|
243 |
character*( 80) fnamegeneric |
244 |
|
245 |
c == external == |
246 |
|
247 |
integer ilnblnk |
248 |
external ilnblnk |
249 |
|
250 |
c == end of interface == |
251 |
|
252 |
#ifdef ALLOW_AUTODIFF_TAMC |
253 |
act3 = myThid - 1 |
254 |
max3 = nTx*nTy |
255 |
act4 = 0 |
256 |
ikey = (act3 + 1) + act4*max3 |
257 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
258 |
|
259 |
jtlo = mybylo(mythid) |
260 |
jthi = mybyhi(mythid) |
261 |
itlo = mybxlo(mythid) |
262 |
ithi = mybxhi(mythid) |
263 |
c-- only do interior, and exchange outside |
264 |
jmin = 1 |
265 |
jmax = sny |
266 |
imin = 1 |
267 |
imax = snx |
268 |
|
269 |
doglobalread = .false. |
270 |
ladinit = .false. |
271 |
|
272 |
|
273 |
call mdsreadfield(wFileCur,32,'RL',1,wFld2d,1,mythid) |
274 |
_EXCH_XY_RL( wFld2d, mythid ) |
275 |
|
276 |
il=ilnblnk( xxFileCur ) |
277 |
write(fnamegeneric(1:80),'(2a,i10.10)') |
278 |
& xxFileCur(1:il),'.',optimcycle |
279 |
call active_read_xy( fnamegeneric, xxFld2d, 1, |
280 |
& doglobalread, ladinit, optimcycle, mythid, xxDummyCur ) |
281 |
|
282 |
|
283 |
#ifndef ALLOW_SMOOTH_CORREL2D |
284 |
c avoid xx larger than boundsVec(5) X uncertainty |
285 |
if ( boundsVec(5).GT.0.) then |
286 |
do bj = jtlo,jthi |
287 |
do bi = itlo,ithi |
288 |
do j = jmin,jmax |
289 |
do i = imin,imax |
290 |
if ( (maskFld3d(i,j,1,bi,bj).NE.0.).AND. |
291 |
& (wFld2d(i,j,bi,bj).GT.0.) ) then |
292 |
tmpMax=boundsVec(5)/sqrt(wFld2d(i,j,bi,bj)) |
293 |
if ( abs(xxFld2d(i,j,bi,bj)).GT.tmpMax ) then |
294 |
xxFld2d(i,j,bi,bj)=sign(tmpMax,xxFld2d(i,j,bi,bj)) |
295 |
else |
296 |
xxFld2d(i,j,bi,bj)=xxFld2d(i,j,bi,bj) |
297 |
endif |
298 |
endif |
299 |
enddo |
300 |
enddo |
301 |
enddo |
302 |
enddo |
303 |
endif |
304 |
#else |
305 |
c apply Weaver And Courtier correlation operator |
306 |
if (paramSmooth.NE.0) then |
307 |
call smooth_correl2d(xxFld2d,maskFld3d,paramSmooth,mythid) |
308 |
endif |
309 |
#endif |
310 |
|
311 |
do bj = jtlo,jthi |
312 |
do bi = itlo,ithi |
313 |
do j = jmin,jmax |
314 |
do i = imin,imax |
315 |
#ifdef ALLOW_SMOOTH_CORREL2D |
316 |
c scale param adjustment |
317 |
if ( (maskFld3d(i,j,1,bi,bj).NE.0.) |
318 |
& .AND. (wFld2d(i,j,bi,bj).GT.0.) ) then |
319 |
xxFld2d(i,j,bi,bj)=xxFld2d(i,j,bi,bj) |
320 |
& /sqrt( wFld2d(i,j,bi,bj) ) |
321 |
else |
322 |
xxFld2d(i,j,bi,bj)=0. |
323 |
endif |
324 |
#endif |
325 |
paramFld2d(i,j,bi,bj) = paramFld2d(i,j,bi,bj) |
326 |
& + xxFld2d(i,j,bi,bj) |
327 |
enddo |
328 |
enddo |
329 |
enddo |
330 |
enddo |
331 |
|
332 |
CALL CTRL_BOUND_2D(paramFld2d,maskFld3d,boundsVec,myThid) |
333 |
|
334 |
#ifdef ALLOW_SMOOTH_CORREL2D |
335 |
write(fnamegeneric(1:80),'(2a,i10.10)') |
336 |
& xxFileCur(1:il),'.effective.',optimcycle |
337 |
call mdswritefield(fnamegeneric,ctrlprec,.FALSE.,'RL', |
338 |
& 1, paramFld2d, 1, optimcycle, mythid) |
339 |
#endif |
340 |
|
341 |
|
342 |
end |
343 |
|
344 |
|
345 |
|
346 |
|
347 |
|
348 |
|