/[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.17 - (show annotations) (download)
Thu Apr 7 23:38:43 2005 UTC (20 years, 2 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint57m_post, checkpoint57g_pre, checkpoint57s_post, checkpoint57g_post, checkpoint57y_post, checkpoint57r_post, checkpoint57i_post, checkpoint57n_post, checkpoint57z_post, checkpoint57l_post, checkpoint57t_post, checkpoint57v_post, checkpoint57h_pre, checkpoint57h_post, checkpoint57y_pre, checkpoint57p_post, checkpint57u_post, checkpoint57q_post, checkpoint57h_done, checkpoint57j_post, checkpoint57o_post, checkpoint57k_post, checkpoint57w_post, checkpoint57x_post
Changes since 1.16: +151 -37 lines
o separate masks used for ctrl_pack/unpack 'from write_grid' output
  (suggested by G. Forget)
o added new control variables
  * init. uVel, vVel, etanN
  * lambda[Theta,Salt]ClimRelax

1 C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_map_ini.F,v 1.16 2005/02/28 19:05:12 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_loc( 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 theta(i,j,k,bi,bj) = theta(i,j,k,bi,bj) +
112 & fac*tmpfld3d(i,j,k,bi,bj)
113 if(theta(i,j,k,bi,bj).lt.-2.0)
114 & theta(i,j,k,bi,bj)= -2.0
115 enddo
116 enddo
117 enddo
118 enddo
119 enddo
120
121 #endif
122
123 #ifdef ALLOW_SALT0_CONTROL
124 c-- Temperature field.
125 il=ilnblnk( xx_salt_file )
126 write(fnamegeneric(1:80),'(2a,i10.10)')
127 & xx_salt_file(1:il),'.',optimcycle
128 call active_read_xyz_loc( fnamegeneric, tmpfld3d, 1,
129 & doglobalread, ladinit, optimcycle,
130 & mythid, xx_salt_dummy )
131
132 do bj = jtlo,jthi
133 do bi = itlo,ithi
134 do k = 1,nr
135 do j = jmin,jmax
136 do i = imin,imax
137 #ifdef ALLOW_ECCO
138 IF (abs(tmpfld3d(i,j,k,bi,bj)).gt.
139 $ 2.0/sqrt(wsalt(k,bi,bj)))
140 $ tmpfld3d(i,j,k,bi,bj)=
141 $ sign(2.0/sqrt(wsalt(k,bi,bj)),tmpfld3d(i,j,k,bi,bj))
142 #endif
143 salt(i,j,k,bi,bj) = salt(i,j,k,bi,bj) +
144 & fac*tmpfld3d(i,j,k,bi,bj)
145
146 enddo
147 enddo
148 enddo
149 enddo
150 enddO
151 #endif
152
153 #ifdef ALLOW_TR10_CONTROL
154 #ifdef ALLOW_PTRACERS
155 c-- Temperature field.
156 il=ilnblnk( xx_tr1_file )
157 write(fnamegeneric(1:80),'(2a,i10.10)')
158 & xx_tr1_file(1:il),'.',optimcycle
159 call active_read_xyz_loc( fnamegeneric, tmpfld3d, 1,
160 & doglobalread, ladinit, optimcycle,
161 & mythid, xx_tr1_dummy )
162
163 do bj = jtlo,jthi
164 do bi = itlo,ithi
165 do k = 1,nr
166 do j = jmin,jmax
167 do i = imin,imax
168 ptracer(i,j,k,bi,bj,1) = ptracer(i,j,k,bi,bj,1) +
169 & fac*tmpfld3d(i,j,k,bi,bj)
170 enddo
171 enddo
172 enddo
173 enddo
174 enddo
175 #endif
176 #endif
177
178 #ifdef ALLOW_SST0_CONTROL
179 c-- sst0.
180 il=ilnblnk( xx_sst_file )
181 write(fnamegeneric(1:80),'(2a,i10.10)')
182 & xx_sst_file(1:il),'.',optimcycle
183 call active_read_xy_loc ( fnamegeneric, tmpfld2d, 1,
184 & doglobalread, ladinit, optimcycle,
185 & mythid, xx_sst_dummy )
186 do bj = jtlo,jthi
187 do bi = itlo,ithi
188 do j = jmin,jmax
189 do i = imin,imax
190 cph sst(i,j,bi,bj) = sst(i,j,bi,bj) + tmpfld2d(i,j,bi,bj)
191 theta(i,j,1,bi,bj) = theta(i,j,1,bi,bj)
192 & + tmpfld2d(i,j,bi,bj)
193 enddo
194 enddo
195 enddo
196 enddo
197 #endif
198
199 #ifdef ALLOW_SSS0_CONTROL
200 c-- sss0.
201 il=ilnblnk( xx_sss_file )
202 write(fnamegeneric(1:80),'(2a,i10.10)')
203 & xx_sss_file(1:il),'.',optimcycle
204 call active_read_xy_loc ( fnamegeneric, tmpfld2d, 1,
205 & doglobalread, ladinit, optimcycle,
206 & mythid, xx_sss_dummy )
207 do bj = jtlo,jthi
208 do bi = itlo,ithi
209 do j = jmin,jmax
210 do i = imin,imax
211 cph sss(i,j,bi,bj) = sss(i,j,bi,bj) + tmpfld2d(i,j,bi,bj)
212 salt(i,j,1,bi,bj) = salt(i,j,1,bi,bj)
213 & + tmpfld2d(i,j,bi,bj)
214 enddo
215 enddo
216 enddo
217 enddo
218 #endif
219
220 #ifdef ALLOW_DIFFKR_CONTROL
221 c-- diffkr.
222 il=ilnblnk( xx_diffkr_file )
223 write(fnamegeneric(1:80),'(2a,i10.10)')
224 & xx_diffkr_file(1:il),'.',optimcycle
225 call active_read_xyz_loc( fnamegeneric, tmpfld3d, 1,
226 & doglobalread, ladinit, optimcycle,
227 & mythid, xx_diffkr_dummy )
228 do bj = jtlo,jthi
229 do bi = itlo,ithi
230 do k = 1,nr
231 do j = jmin,jmax
232 do i = imin,imax
233 diffkr(i,j,k,bi,bj) = diffkr(i,j,k,bi,bj) +
234 & tmpfld3d(i,j,k,bi,bj)
235 enddo
236 enddo
237 enddo
238 enddo
239 enddo
240 #endif
241
242 #ifdef ALLOW_KAPGM_CONTROL
243 c-- kapgm.
244 il=ilnblnk( xx_kapgm_file )
245 write(fnamegeneric(1:80),'(2a,i10.10)')
246 & xx_kapgm_file(1:il),'.',optimcycle
247 call active_read_xyz_loc( fnamegeneric, tmpfld3d, 1,
248 & doglobalread, ladinit, optimcycle,
249 & mythid, xx_kapgm_dummy )
250 do bj = jtlo,jthi
251 do bi = itlo,ithi
252 do k = 1,nr
253 do j = jmin,jmax
254 do i = imin,imax
255 kapgm(i,j,k,bi,bj) = kapgm(i,j,k,bi,bj) +
256 & tmpfld3d(i,j,k,bi,bj)
257 enddo
258 enddo
259 enddo
260 enddo
261 enddo
262 #endif
263
264 #ifdef ALLOW_EFLUXY0_CONTROL
265 c-- y-component EP-flux field.
266 il=ilnblnk( xx_efluxy_file )
267 write(fnamegeneric(1:80),'(2a,i10.10)')
268 & xx_efluxy_file(1:il),'.',optimcycle
269 call active_read_xyz_loc( fnamegeneric, tmpfld3d, 1,
270 & doglobalread, ladinit, optimcycle,
271 & mythid, xx_efluxy_dummy )
272
273 do bj = jtlo,jthi
274 do bi = itlo,ithi
275 do k = 1,nr
276 do j = jmin,jmax
277 do i = imin,imax
278 EfluxY(i,j,k,bi,bj) = EfluxY(i,j,k,bi,bj)
279 & - fac*tmpfld3d(i,j,k,bi,bj)
280 & *maskS(i,j,k,bi,bj)
281 cph EfluxY(i,j,k,bi,bj) = EfluxY(i,j,k,bi,bj)
282 cph & - rSphere*cosFacU(J,bi,bj)
283 cph & *fac*tmpfld3d(i,j,k,bi,bj)
284 enddo
285 enddo
286 enddo
287 enddo
288 enddo
289 #endif
290
291 #ifdef ALLOW_EFLUXP0_CONTROL
292 c-- p-component EP-flux field.
293 il=ilnblnk( xx_efluxp_file )
294 write(fnamegeneric(1:80),'(2a,i10.10)')
295 & xx_efluxp_file(1:il),'.',optimcycle
296 call active_read_xyz_loc( fnamegeneric, tmpfld3d, 1,
297 & doglobalread, ladinit, optimcycle,
298 & mythid, xx_efluxp_dummy )
299
300 do bj = jtlo,jthi
301 do bi = itlo,ithi
302 do k = 1,nr
303 do j = jmin,jmax
304 do i = imin,imax
305 EfluxP(i,j,k,bi,bj) = EfluxP(i,j,k,bi,bj)
306 & + fCori(i,j,bi,bj)
307 & *fac*tmpfld3d(i,j,k,bi,bj)
308 & *hFacV(i,j,k,bi,bj)
309 cph EfluxP(i,j,k,bi,bj) = EfluxP(i,j,k,bi,bj)
310 cph & + fCori(i,j,bi,bj)
311 cph & *rSphere*cosFacU(J,bi,bj)
312 cph & *fac*tmpfld3d(i,j,k,bi,bj)
313 enddo
314 enddo
315 enddo
316 enddo
317 enddo
318 #endif
319
320 #ifdef ALLOW_BOTTOMDRAG_CONTROL
321 c-- bottom drag
322 il=ilnblnk( xx_bottomdrag_file )
323 write(fnamegeneric(1:80),'(2a,i10.10)')
324 & xx_bottomdrag_file(1:il),'.',optimcycle
325 call active_read_xy_loc ( fnamegeneric, tmpfld2d, 1,
326 & doglobalread, ladinit, optimcycle,
327 & mythid, xx_bottomdrag_dummy )
328 do bj = jtlo,jthi
329 do bi = itlo,ithi
330 do j = jmin,jmax
331 do i = imin,imax
332 bottomdragfld(i,j,bi,bj) = bottomdragfld(i,j,bi,bj)
333 & + tmpfld2d(i,j,bi,bj)
334 enddo
335 enddo
336 enddo
337 enddo
338 #endif
339
340 #ifdef ALLOW_EDTAUX_CONTROL
341 c-- zonal eddy stress : edtaux
342 il=ilnblnk( xx_edtaux_file )
343 write(fnamegeneric(1:80),'(2a,i10.10)')
344 & xx_edtaux_file(1:il),'.',optimcycle
345 call active_read_xyz( fnamegeneric, tmpfld3d, 1,
346 & doglobalread, ladinit, optimcycle,
347 & mythid, xx_edtaux_dummy )
348 do bj = jtlo,jthi
349 do bi = itlo,ithi
350 do k = 1,nr
351 do j = jmin,jmax
352 do i = imin,imax
353 Eddytaux(i,j,k,bi,bj) = Eddytaux(i,j,k,bi,bj) +
354 & tmpfld3d(i,j,k,bi,bj)
355 enddo
356 enddo
357 enddo
358 enddo
359 enddo
360 #endif
361
362 #ifdef ALLOW_EDTAUY_CONTROL
363 c-- meridional eddy stress : edtauy
364 il=ilnblnk( xx_edtauy_file )
365 write(fnamegeneric(1:80),'(2a,i10.10)')
366 & xx_edtauy_file(1:il),'.',optimcycle
367 call active_read_xyz( fnamegeneric, tmpfld3d, 1,
368 & doglobalread, ladinit, optimcycle,
369 & mythid, xx_edtauy_dummy )
370 do bj = jtlo,jthi
371 do bi = itlo,ithi
372 do k = 1,nr
373 do j = jmin,jmax
374 do i = imin,imax
375 Eddytauy(i,j,k,bi,bj) = Eddytauy(i,j,k,bi,bj) +
376 & tmpfld3d(i,j,k,bi,bj)
377 enddo
378 enddo
379 enddo
380 enddo
381 enddo
382 #endif
383
384 #ifdef ALLOW_UVEL0_CONTROL
385 c-- initial zonal velocity
386 il=ilnblnk( xx_uvel_file )
387 write(fnamegeneric(1:80),'(2a,i10.10)')
388 & xx_uvel_file(1:il),'.',optimcycle
389 call active_read_xyz( fnamegeneric, tmpfld3d, 1,
390 & doglobalread, ladinit, optimcycle,
391 & mythid, xx_uvel_dummy )
392 do bj = jtlo,jthi
393 do bi = itlo,ithi
394 do k = 1,nr
395 do j = jmin,jmax
396 do i = imin,imax
397 uVel(i,j,k,bi,bj) = uVel(i,j,k,bi,bj) +
398 & tmpfld3d(i,j,k,bi,bj)
399 enddo
400 enddo
401 enddo
402 enddo
403 enddo
404 #endif
405
406 #ifdef ALLOW_VVEL0_CONTROL
407 c-- initial merid. velocity
408 il=ilnblnk( xx_vvel_file )
409 write(fnamegeneric(1:80),'(2a,i10.10)')
410 & xx_vvel_file(1:il),'.',optimcycle
411 call active_read_xyz( fnamegeneric, tmpfld3d, 1,
412 & doglobalread, ladinit, optimcycle,
413 & mythid, xx_vvel_dummy )
414 do bj = jtlo,jthi
415 do bi = itlo,ithi
416 do k = 1,nr
417 do j = jmin,jmax
418 do i = imin,imax
419 vVel(i,j,k,bi,bj) = vVel(i,j,k,bi,bj) +
420 & tmpfld3d(i,j,k,bi,bj)
421 enddo
422 enddo
423 enddo
424 enddo
425 enddo
426 #endif
427
428 #ifdef ALLOW_ETAN0_CONTROL
429 c-- initial Eta.
430 il=ilnblnk( xx_etan_file )
431 write(fnamegeneric(1:80),'(2a,i10.10)')
432 & xx_etan_file(1:il),'.',optimcycle
433 call active_read_xy_loc ( fnamegeneric, tmpfld2d, 1,
434 & doglobalread, ladinit, optimcycle,
435 & mythid, xx_etan_dummy )
436 do bj = jtlo,jthi
437 do bi = itlo,ithi
438 do j = jmin,jmax
439 do i = imin,imax
440 etaN(i,j,bi,bj) = etaN(i,j,bi,bj) + tmpfld2d(i,j,bi,bj)
441 enddo
442 enddo
443 enddo
444 enddo
445 #endif
446
447 #ifdef ALLOW_RELAXSST_CONTROL
448 c-- SST relaxation coefficient.
449 il=ilnblnk( xx_relaxsst_file )
450 write(fnamegeneric(1:80),'(2a,i10.10)')
451 & xx_relaxsst_file(1:il),'.',optimcycle
452 call active_read_xy_loc ( fnamegeneric, tmpfld2d, 1,
453 & doglobalread, ladinit, optimcycle,
454 & mythid, xx_relaxsst_dummy )
455 do bj = jtlo,jthi
456 do bi = itlo,ithi
457 do j = jmin,jmax
458 do i = imin,imax
459 lambdaThetaClimRelax(i,j,bi,bj) =
460 & lambdaThetaClimRelax(i,j,bi,bj)
461 & + tmpfld2d(i,j,bi,bj)
462 enddo
463 enddo
464 enddo
465 enddo
466 #endif
467
468 #ifdef ALLOW_RELAXSSS_CONTROL
469 c-- SSS relaxation coefficient.
470 il=ilnblnk( xx_relaxsss_file )
471 write(fnamegeneric(1:80),'(2a,i10.10)')
472 & xx_relaxsss_file(1:il),'.',optimcycle
473 call active_read_xy_loc ( fnamegeneric, tmpfld2d, 1,
474 & doglobalread, ladinit, optimcycle,
475 & mythid, xx_relaxsss_dummy )
476 do bj = jtlo,jthi
477 do bi = itlo,ithi
478 do j = jmin,jmax
479 do i = imin,imax
480 lambdaSaltClimRelax(i,j,bi,bj) =
481 & lambdaSaltClimRelax(i,j,bi,bj)
482 & + tmpfld2d(i,j,bi,bj)
483 enddo
484 enddo
485 enddo
486 enddo
487 #endif
488
489 c-- Update the tile edges.
490
491 #if (defined (ALLOW_THETA0_CONTROL) || defined (ALLOW_SST0_CONTROL))
492 _EXCH_XYZ_R8( theta, mythid )
493 #endif
494 #if (defined (ALLOW_SALT0_CONTROL) || defined (ALLOW_SSS0_CONTROL))
495 _EXCH_XYZ_R8( salt, mythid )
496 #endif
497 #ifdef ALLOW_TR10_CONTROL
498 #ifdef ALLOW_PTRACERS
499 _EXCH_XYZ_R8(pTracer(1-Olx,1-Oly,1,1,1,1),myThid)
500 #endif
501 #endif
502 #ifdef ALLOW_DIFFKR_CONTROL
503 _EXCH_XYZ_R8( diffkr, mythid)
504 #endif
505 #ifdef ALLOW_KAPGM_CONTROL
506 _EXCH_XYZ_R8( kapgm, mythid)
507 #endif
508 #ifdef ALLOW_EFLUXY0_CONTROL
509 _EXCH_XYZ_R8( EfluxY, mythid )
510 #endif
511 #ifdef ALLOW_EFLUXP0_CONTROL
512 _EXCH_XYZ_R8( EfluxP, mythid )
513 #endif
514 #ifdef ALLOW_BOTTOMDRAG_CONTROL
515 _EXCH_XY_R8( bottomdragfld, mythid )
516 #endif
517
518 #if (defined (ALLOW_EDTAUX_CONTROL) && defined (ALLOW_EDTAUY_CONTROL))
519 CALL EXCH_UV_XYZ_RS(Eddytaux,Eddytauy,.TRUE.,myThid)
520 #elif (defined (ALLOW_EDTAUX_CONTROL) || defined (ALLOW_EDTAUY_CONTROL))
521 STOP 'ctrl_map_forcing: need BOTH ALLOW_EDTAU[X,Y]_CONTROL'
522 #endif
523
524 #ifdef ALLOW_UVEL0_CONTROL
525 _EXCH_XYZ_R8( uVel, mythid)
526 #endif
527
528 #ifdef ALLOW_VVEL0_CONTROL
529 _EXCH_XYZ_R8( vVel, mythid)
530 #endif
531
532 #ifdef ALLOW_ETAN0_CONTROL
533 _EXCH_XY_R8( etaN, mythid )
534 #endif
535
536 #ifdef ALLOW_RELAXSST_CONTROL
537 _EXCH_XY_R4( lambdaThetaClimRelax, mythid )
538 #endif
539
540 #ifdef ALLOW_RELAXSSS_CONTROL
541 _EXCH_XY_R4( lambdaThetaClimRelax, mythid )
542 #endif
543
544 return
545 end
546

  ViewVC Help
Powered by ViewVC 1.1.22