/[MITgcm]/MITgcm/pkg/ctrl/ctrl_map_ini.F
ViewVC logotype

Contents of /MITgcm/pkg/ctrl/ctrl_map_ini.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.23 - (show annotations) (download)
Wed Jun 6 16:59:35 2007 UTC (17 years ago) by utke
Branch: MAIN
CVS Tags: checkpoint59d
Changes since 1.22: +6 -6 lines
renamed macro

1 C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_map_ini.F,v 1.22 2007/05/14 22:02:34 heimbach Exp $
2
3 #include "CTRL_CPPOPTIONS.h"
4
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
19 C !USES:
20 implicit none
21
22 c == global variables ==
23 #include "SIZE.h"
24 #include "EEPARAMS.h"
25 #include "PARAMS.h"
26 #include "GRID.h"
27 #include "DYNVARS.h"
28 #include "FFIELDS.h"
29 #include "ctrl.h"
30 #include "ctrl_dummy.h"
31 #include "optim.h"
32 #ifdef ALLOW_PTRACERS
33 # include "PTRACERS_SIZE.h"
34 # include "PTRACERS.h"
35 #endif
36 #ifdef ALLOW_ECCO
37 # include "ecco_cost.h"
38 #endif
39
40 C !INPUT/OUTPUT PARAMETERS:
41 c == routine arguments ==
42 integer mythid
43
44 C !LOCAL VARIABLES:
45 c == local variables ==
46
47 integer bi,bj
48 integer i,j,k
49 integer itlo,ithi
50 integer jtlo,jthi
51 integer jmin,jmax
52 integer imin,imax
53 integer il
54
55 logical equal
56 logical doglobalread
57 logical ladinit
58
59 character*( 80) fnamegeneric
60
61 _RL fac
62 _RL tmptest
63
64 c == external ==
65 integer ilnblnk
66 external ilnblnk
67
68 c == end of interface ==
69 CEOP
70
71 jtlo = mybylo(mythid)
72 jthi = mybyhi(mythid)
73 itlo = mybxlo(mythid)
74 ithi = mybxhi(mythid)
75 jmin = 1
76 jmax = sny
77 imin = 1
78 imax = snx
79
80 doglobalread = .false.
81 ladinit = .false.
82
83 equal = .true.
84
85 if ( equal ) then
86 fac = 1. _d 0
87 else
88 fac = 0. _d 0
89 endif
90
91 #ifdef ALLOW_THETA0_CONTROL
92 c-- Temperature field.
93 il=ilnblnk( xx_theta_file )
94 write(fnamegeneric(1:80),'(2a,i10.10)')
95 & xx_theta_file(1:il),'.',optimcycle
96 call active_read_xyz( fnamegeneric, tmpfld3d, 1,
97 & doglobalread, ladinit, optimcycle,
98 & mythid, xx_theta_dummy )
99
100 do bj = jtlo,jthi
101 do bi = itlo,ithi
102 do k = 1,nr
103 do j = jmin,jmax
104 do i = imin,imax
105 #ifdef ALLOW_ECCO
106 IF (abs(tmpfld3d(i,j,k,bi,bj)).gt.
107 $ 2.0/sqrt(wtheta(k,bi,bj)))
108 $ tmpfld3d(i,j,k,bi,bj)=
109 $ sign(2.0/sqrt(wtheta(k,bi,bj)),tmpfld3d(i,j,k,bi,bj))
110 #endif
111 #ifdef ALLOW_AUTODIFF_OPENAD
112 theta(i,j,k,bi,bj) = theta(i,j,k,bi,bj) +
113 & fac*xx_theta(i,j,k,bi,bj) +
114 & fac*tmpfld3d(i,j,k,bi,bj)
115 #else
116 theta(i,j,k,bi,bj) = theta(i,j,k,bi,bj) +
117 & fac*tmpfld3d(i,j,k,bi,bj)
118 #endif
119 if(theta(i,j,k,bi,bj).lt.-2.0)
120 & theta(i,j,k,bi,bj)= -2.0
121 enddo
122 enddo
123 enddo
124 enddo
125 enddo
126
127 #endif
128
129 #ifdef ALLOW_SALT0_CONTROL
130 c-- Temperature field.
131 il=ilnblnk( xx_salt_file )
132 write(fnamegeneric(1:80),'(2a,i10.10)')
133 & xx_salt_file(1:il),'.',optimcycle
134 call active_read_xyz( fnamegeneric, tmpfld3d, 1,
135 & 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 #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 #ifdef ALLOW_AUTODIFF_OPENAD
150 salt(i,j,k,bi,bj) = salt(i,j,k,bi,bj) +
151 & fac*xx_salt(i,j,k,bi,bj) +
152 & fac*tmpfld3d(i,j,k,bi,bj)
153 #else
154 salt(i,j,k,bi,bj) = salt(i,j,k,bi,bj) +
155 & fac*tmpfld3d(i,j,k,bi,bj)
156 #endif
157
158 enddo
159 enddo
160 enddo
161 enddo
162 enddO
163 #endif
164
165 #ifdef ALLOW_TR10_CONTROL
166 #ifdef ALLOW_PTRACERS
167 c-- Temperature field.
168 il=ilnblnk( xx_tr1_file )
169 write(fnamegeneric(1:80),'(2a,i10.10)')
170 & xx_tr1_file(1:il),'.',optimcycle
171 call active_read_xyz( fnamegeneric, tmpfld3d, 1,
172 & doglobalread, ladinit, optimcycle,
173 & mythid, xx_tr1_dummy )
174
175 do bj = jtlo,jthi
176 do bi = itlo,ithi
177 do k = 1,nr
178 do j = jmin,jmax
179 do i = imin,imax
180 ptracer(i,j,k,bi,bj,1) = ptracer(i,j,k,bi,bj,1) +
181 & fac*tmpfld3d(i,j,k,bi,bj)
182 enddo
183 enddo
184 enddo
185 enddo
186 enddo
187 #endif
188 #endif
189
190 #ifdef ALLOW_SST0_CONTROL
191 c-- sst0.
192 il=ilnblnk( xx_sst_file )
193 write(fnamegeneric(1:80),'(2a,i10.10)')
194 & xx_sst_file(1:il),'.',optimcycle
195 call active_read_xy ( fnamegeneric, tmpfld2d, 1,
196 & doglobalread, ladinit, optimcycle,
197 & mythid, xx_sst_dummy )
198 do bj = jtlo,jthi
199 do bi = itlo,ithi
200 do j = jmin,jmax
201 do i = imin,imax
202 cph sst(i,j,bi,bj) = sst(i,j,bi,bj) + tmpfld2d(i,j,bi,bj)
203 theta(i,j,1,bi,bj) = theta(i,j,1,bi,bj)
204 & + tmpfld2d(i,j,bi,bj)
205 enddo
206 enddo
207 enddo
208 enddo
209 #endif
210
211 #ifdef ALLOW_SSS0_CONTROL
212 c-- sss0.
213 il=ilnblnk( xx_sss_file )
214 write(fnamegeneric(1:80),'(2a,i10.10)')
215 & xx_sss_file(1:il),'.',optimcycle
216 call active_read_xy ( fnamegeneric, tmpfld2d, 1,
217 & doglobalread, ladinit, optimcycle,
218 & mythid, xx_sss_dummy )
219 do bj = jtlo,jthi
220 do bi = itlo,ithi
221 do j = jmin,jmax
222 do i = imin,imax
223 cph sss(i,j,bi,bj) = sss(i,j,bi,bj) + tmpfld2d(i,j,bi,bj)
224 salt(i,j,1,bi,bj) = salt(i,j,1,bi,bj)
225 & + tmpfld2d(i,j,bi,bj)
226 enddo
227 enddo
228 enddo
229 enddo
230 #endif
231
232 #ifdef ALLOW_DIFFKR_CONTROL
233 c-- diffkr.
234 il=ilnblnk( xx_diffkr_file )
235 write(fnamegeneric(1:80),'(2a,i10.10)')
236 & xx_diffkr_file(1:il),'.',optimcycle
237 call active_read_xyz( fnamegeneric, tmpfld3d, 1,
238 & doglobalread, ladinit, optimcycle,
239 & mythid, xx_diffkr_dummy )
240 do bj = jtlo,jthi
241 do bi = itlo,ithi
242 do k = 1,nr
243 do j = jmin,jmax
244 do i = imin,imax
245 diffkr(i,j,k,bi,bj) = diffkr(i,j,k,bi,bj) +
246 & tmpfld3d(i,j,k,bi,bj)
247 enddo
248 enddo
249 enddo
250 enddo
251 enddo
252 #endif
253
254 #ifdef ALLOW_KAPGM_CONTROL
255 c-- kapgm.
256 il=ilnblnk( xx_kapgm_file )
257 write(fnamegeneric(1:80),'(2a,i10.10)')
258 & xx_kapgm_file(1:il),'.',optimcycle
259 call active_read_xyz( fnamegeneric, tmpfld3d, 1,
260 & doglobalread, ladinit, optimcycle,
261 & mythid, xx_kapgm_dummy )
262 do bj = jtlo,jthi
263 do bi = itlo,ithi
264 do k = 1,nr
265 do j = jmin,jmax
266 do i = imin,imax
267 kapgm(i,j,k,bi,bj) = kapgm(i,j,k,bi,bj) +
268 & tmpfld3d(i,j,k,bi,bj)
269 enddo
270 enddo
271 enddo
272 enddo
273 enddo
274 #endif
275
276 #ifdef ALLOW_EFLUXY0_CONTROL
277 c-- y-component EP-flux field.
278 il=ilnblnk( xx_efluxy_file )
279 write(fnamegeneric(1:80),'(2a,i10.10)')
280 & xx_efluxy_file(1:il),'.',optimcycle
281 call active_read_xyz( fnamegeneric, tmpfld3d, 1,
282 & doglobalread, ladinit, optimcycle,
283 & mythid, xx_efluxy_dummy )
284
285 do bj = jtlo,jthi
286 do bi = itlo,ithi
287 do k = 1,nr
288 do j = jmin,jmax
289 do i = imin,imax
290 EfluxY(i,j,k,bi,bj) = EfluxY(i,j,k,bi,bj)
291 & - fac*tmpfld3d(i,j,k,bi,bj)
292 & *maskS(i,j,k,bi,bj)
293 cph EfluxY(i,j,k,bi,bj) = EfluxY(i,j,k,bi,bj)
294 cph & - rSphere*cosFacU(J,bi,bj)
295 cph & *fac*tmpfld3d(i,j,k,bi,bj)
296 enddo
297 enddo
298 enddo
299 enddo
300 enddo
301 #endif
302
303 #ifdef ALLOW_EFLUXP0_CONTROL
304 c-- p-component EP-flux field.
305 il=ilnblnk( xx_efluxp_file )
306 write(fnamegeneric(1:80),'(2a,i10.10)')
307 & xx_efluxp_file(1:il),'.',optimcycle
308 call active_read_xyz( fnamegeneric, tmpfld3d, 1,
309 & doglobalread, ladinit, optimcycle,
310 & mythid, xx_efluxp_dummy )
311
312 do bj = jtlo,jthi
313 do bi = itlo,ithi
314 do k = 1,nr
315 do j = jmin,jmax
316 do i = imin,imax
317 EfluxP(i,j,k,bi,bj) = EfluxP(i,j,k,bi,bj)
318 & + fCori(i,j,bi,bj)
319 & *fac*tmpfld3d(i,j,k,bi,bj)
320 & *hFacV(i,j,k,bi,bj)
321 cph EfluxP(i,j,k,bi,bj) = EfluxP(i,j,k,bi,bj)
322 cph & + fCori(i,j,bi,bj)
323 cph & *rSphere*cosFacU(J,bi,bj)
324 cph & *fac*tmpfld3d(i,j,k,bi,bj)
325 enddo
326 enddo
327 enddo
328 enddo
329 enddo
330 #endif
331
332 #ifdef ALLOW_BOTTOMDRAG_CONTROL
333 c-- bottom drag
334 il=ilnblnk( xx_bottomdrag_file )
335 write(fnamegeneric(1:80),'(2a,i10.10)')
336 & xx_bottomdrag_file(1:il),'.',optimcycle
337 call active_read_xy ( fnamegeneric, tmpfld2d, 1,
338 & doglobalread, ladinit, optimcycle,
339 & mythid, xx_bottomdrag_dummy )
340 do bj = jtlo,jthi
341 do bi = itlo,ithi
342 do j = jmin,jmax
343 do i = imin,imax
344 bottomdragfld(i,j,bi,bj) = bottomdragfld(i,j,bi,bj)
345 & + tmpfld2d(i,j,bi,bj)
346 enddo
347 enddo
348 enddo
349 enddo
350 #endif
351
352 #ifdef ALLOW_EDTAUX_CONTROL
353 c-- zonal eddy stress : edtaux
354 il=ilnblnk( xx_edtaux_file )
355 write(fnamegeneric(1:80),'(2a,i10.10)')
356 & xx_edtaux_file(1:il),'.',optimcycle
357 call active_read_xyz( fnamegeneric, tmpfld3d, 1,
358 & doglobalread, ladinit, optimcycle,
359 & mythid, xx_edtaux_dummy )
360 do bj = jtlo,jthi
361 do bi = itlo,ithi
362 do k = 1,nr
363 do j = jmin,jmax
364 do i = imin,imax
365 eddyTauX(i,j,k,bi,bj) = eddyTauX(i,j,k,bi,bj) +
366 & fCori(i,j,bi,bj)*tmpfld3d(i,j,k,bi,bj)
367 enddo
368 enddo
369 enddo
370 enddo
371 enddo
372 #endif
373
374 #ifdef ALLOW_EDTAUY_CONTROL
375 c-- meridional eddy stress : edtauy
376 il=ilnblnk( xx_edtauy_file )
377 write(fnamegeneric(1:80),'(2a,i10.10)')
378 & xx_edtauy_file(1:il),'.',optimcycle
379 call active_read_xyz( fnamegeneric, tmpfld3d, 1,
380 & doglobalread, ladinit, optimcycle,
381 & mythid, xx_edtauy_dummy )
382 do bj = jtlo,jthi
383 do bi = itlo,ithi
384 do k = 1,nr
385 do j = jmin,jmax
386 do i = imin,imax
387 eddyTauY(i,j,k,bi,bj) = eddyTauY(i,j,k,bi,bj) +
388 & fCoriG(i,j,bi,bj)*tmpfld3d(i,j,k,bi,bj)
389 enddo
390 enddo
391 enddo
392 enddo
393 enddo
394 #endif
395
396 #ifdef ALLOW_UVEL0_CONTROL
397 c-- initial zonal velocity
398 il=ilnblnk( xx_uvel_file )
399 write(fnamegeneric(1:80),'(2a,i10.10)')
400 & xx_uvel_file(1:il),'.',optimcycle
401 call active_read_xyz( fnamegeneric, tmpfld3d, 1,
402 & doglobalread, ladinit, optimcycle,
403 & mythid, xx_uvel_dummy )
404 do bj = jtlo,jthi
405 do bi = itlo,ithi
406 do k = 1,nr
407 do j = jmin,jmax
408 do i = imin,imax
409 #ifdef ALLOW_AUTODIFF_OPENAD
410 uVel(i,j,k,bi,bj) = uVel(i,j,k,bi,bj) +
411 & fac*xx_uvel(i,j,k,bi,bj)
412 #else
413 uVel(i,j,k,bi,bj) = uVel(i,j,k,bi,bj) +
414 & fac*tmpfld3d(i,j,k,bi,bj)
415 #endif
416 enddo
417 enddo
418 enddo
419 enddo
420 enddo
421 #endif
422
423 #ifdef ALLOW_VVEL0_CONTROL
424 c-- initial merid. velocity
425 il=ilnblnk( xx_vvel_file )
426 write(fnamegeneric(1:80),'(2a,i10.10)')
427 & xx_vvel_file(1:il),'.',optimcycle
428 call active_read_xyz( fnamegeneric, tmpfld3d, 1,
429 & doglobalread, ladinit, optimcycle,
430 & mythid, xx_vvel_dummy )
431 do bj = jtlo,jthi
432 do bi = itlo,ithi
433 do k = 1,nr
434 do j = jmin,jmax
435 do i = imin,imax
436 #ifdef ALLOW_AUTODIFF_OPENAD
437 vVel(i,j,k,bi,bj) = vVel(i,j,k,bi,bj) +
438 & fac*xx_vvel(i,j,k,bi,bj)
439 #else
440 vVel(i,j,k,bi,bj) = vVel(i,j,k,bi,bj) +
441 & fac*tmpfld3d(i,j,k,bi,bj)
442 #endif
443 enddo
444 enddo
445 enddo
446 enddo
447 enddo
448 #endif
449
450 #ifdef ALLOW_ETAN0_CONTROL
451 c-- initial Eta.
452 il=ilnblnk( xx_etan_file )
453 write(fnamegeneric(1:80),'(2a,i10.10)')
454 & xx_etan_file(1:il),'.',optimcycle
455 call active_read_xy ( fnamegeneric, tmpfld2d, 1,
456 & doglobalread, ladinit, optimcycle,
457 & mythid, xx_etan_dummy )
458 do bj = jtlo,jthi
459 do bi = itlo,ithi
460 do j = jmin,jmax
461 do i = imin,imax
462 #ifdef ALLOW_AUTODIFF_OPENAD
463 etaN(i,j,bi,bj) = etaN(i,j,bi,bj) +
464 & fac*xx_etan(i,j,bi,bj)
465 #else
466 etaN(i,j,bi,bj) = etaN(i,j,bi,bj) +
467 & fac*tmpfld2d(i,j,bi,bj)
468 #endif
469 enddo
470 enddo
471 enddo
472 enddo
473 #endif
474
475 #ifdef ALLOW_RELAXSST_CONTROL
476 c-- SST relaxation coefficient.
477 il=ilnblnk( xx_relaxsst_file )
478 write(fnamegeneric(1:80),'(2a,i10.10)')
479 & xx_relaxsst_file(1:il),'.',optimcycle
480 call active_read_xy ( fnamegeneric, tmpfld2d, 1,
481 & doglobalread, ladinit, optimcycle,
482 & mythid, xx_relaxsst_dummy )
483 do bj = jtlo,jthi
484 do bi = itlo,ithi
485 do j = jmin,jmax
486 do i = imin,imax
487 lambdaThetaClimRelax(i,j,bi,bj) =
488 & lambdaThetaClimRelax(i,j,bi,bj)
489 & + tmpfld2d(i,j,bi,bj)
490 enddo
491 enddo
492 enddo
493 enddo
494 #endif
495
496 #ifdef ALLOW_RELAXSSS_CONTROL
497 c-- SSS relaxation coefficient.
498 il=ilnblnk( xx_relaxsss_file )
499 write(fnamegeneric(1:80),'(2a,i10.10)')
500 & xx_relaxsss_file(1:il),'.',optimcycle
501 call active_read_xy ( fnamegeneric, tmpfld2d, 1,
502 & doglobalread, ladinit, optimcycle,
503 & mythid, xx_relaxsss_dummy )
504 do bj = jtlo,jthi
505 do bi = itlo,ithi
506 do j = jmin,jmax
507 do i = imin,imax
508 lambdaSaltClimRelax(i,j,bi,bj) =
509 & lambdaSaltClimRelax(i,j,bi,bj)
510 & + tmpfld2d(i,j,bi,bj)
511 enddo
512 enddo
513 enddo
514 enddo
515 #endif
516
517 c-- Update the tile edges.
518
519 #if (defined (ALLOW_THETA0_CONTROL) || defined (ALLOW_SST0_CONTROL))
520 _EXCH_XYZ_R8( theta, mythid )
521 #endif
522 #if (defined (ALLOW_SALT0_CONTROL) || defined (ALLOW_SSS0_CONTROL))
523 _EXCH_XYZ_R8( salt, mythid )
524 #endif
525 #ifdef ALLOW_TR10_CONTROL
526 #ifdef ALLOW_PTRACERS
527 _EXCH_XYZ_R8(pTracer(1-Olx,1-Oly,1,1,1,1),myThid)
528 #endif
529 #endif
530 #ifdef ALLOW_DIFFKR_CONTROL
531 _EXCH_XYZ_R8( diffkr, mythid)
532 #endif
533 #ifdef ALLOW_KAPGM_CONTROL
534 _EXCH_XYZ_R8( kapgm, mythid)
535 #endif
536 #ifdef ALLOW_EFLUXY0_CONTROL
537 _EXCH_XYZ_R8( EfluxY, mythid )
538 #endif
539 #ifdef ALLOW_EFLUXP0_CONTROL
540 _EXCH_XYZ_R8( EfluxP, mythid )
541 #endif
542 #ifdef ALLOW_BOTTOMDRAG_CONTROL
543 _EXCH_XY_R8( bottomdragfld, mythid )
544 #endif
545
546 #if (defined (ALLOW_EDTAUX_CONTROL) && defined (ALLOW_EDTAUY_CONTROL))
547 CALL EXCH_UV_XYZ_RS(eddyTauX,eddyTauY,.TRUE.,myThid)
548 #elif (defined (ALLOW_EDTAUX_CONTROL) || defined (ALLOW_EDTAUY_CONTROL))
549 STOP 'ctrl_map_forcing: need BOTH ALLOW_EDTAU[X,Y]_CONTROL'
550 #endif
551
552 #ifdef ALLOW_UVEL0_CONTROL
553 _EXCH_XYZ_R8( uVel, mythid)
554 #endif
555
556 #ifdef ALLOW_VVEL0_CONTROL
557 _EXCH_XYZ_R8( vVel, mythid)
558 #endif
559
560 #ifdef ALLOW_ETAN0_CONTROL
561 _EXCH_XY_R8( etaN, mythid )
562 #endif
563
564 #ifdef ALLOW_RELAXSST_CONTROL
565 _EXCH_XY_R4( lambdaThetaClimRelax, mythid )
566 #endif
567
568 #ifdef ALLOW_RELAXSSS_CONTROL
569 _EXCH_XY_R4( lambdaThetaClimRelax, mythid )
570 #endif
571
572 return
573 end
574

  ViewVC Help
Powered by ViewVC 1.1.22