/[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.25 - (show annotations) (download)
Tue Oct 9 00:00:00 2007 UTC (17 years, 8 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint59i
Changes since 1.24: +13 -12 lines
add missing cvs $Header:$ or $Name:$

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

  ViewVC Help
Powered by ViewVC 1.1.22