1 |
C $Header: /u/gcmpack/MITgcm/pkg/ecco/cost_kapgm.F,v 1.2 2004/10/13 07:05:51 heimbach Exp $ |
2 |
|
3 |
#include "COST_CPPOPTIONS.h" |
4 |
|
5 |
|
6 |
subroutine cost_kapgm( |
7 |
I myiter, |
8 |
I mytime, |
9 |
I mythid |
10 |
& ) |
11 |
|
12 |
c ================================================================== |
13 |
c SUBROUTINE cost_kapgm |
14 |
c ================================================================== |
15 |
c |
16 |
c o Calculate the Kappa GM contribution to the cost function. |
17 |
c |
18 |
c started:Armin Koehl akoehl@ucsd.edu |
19 |
c - Restructured the code in order to create a package |
20 |
c for the MITgcmUV. |
21 |
c |
22 |
c ================================================================== |
23 |
c SUBROUTINE cost_kapgm |
24 |
c ================================================================== |
25 |
|
26 |
implicit none |
27 |
|
28 |
c == global variables == |
29 |
|
30 |
#include "EEPARAMS.h" |
31 |
#include "SIZE.h" |
32 |
#include "GRID.h" |
33 |
#include "DYNVARS.h" |
34 |
#ifdef ALLOW_GMREDI |
35 |
# include "GMREDI.h" |
36 |
#endif |
37 |
|
38 |
#ifdef ALLOW_ECCO |
39 |
# include "ecco_cost.h" |
40 |
#else |
41 |
# include "ecco_cost.h" |
42 |
#endif |
43 |
#include "ctrl.h" |
44 |
#include "ctrl_dummy.h" |
45 |
#include "optim.h" |
46 |
c == routine arguments == |
47 |
|
48 |
integer myiter |
49 |
_RL mytime |
50 |
integer mythid |
51 |
|
52 |
c == local variables == |
53 |
|
54 |
integer bi,bj |
55 |
integer i,j,k |
56 |
integer itlo,ithi |
57 |
integer jtlo,jthi |
58 |
integer jmin,jmax |
59 |
integer imin,imax |
60 |
integer nrec |
61 |
integer irec |
62 |
integer ilfld |
63 |
|
64 |
_RL fctile |
65 |
_RL fcthread |
66 |
_RL tmpx |
67 |
|
68 |
logical doglobalread |
69 |
logical ladinit |
70 |
|
71 |
character*(80) fnamefld |
72 |
|
73 |
character*(MAX_LEN_MBUF) msgbuf |
74 |
|
75 |
c == external functions == |
76 |
|
77 |
integer ilnblnk |
78 |
external ilnblnk |
79 |
|
80 |
c == end of interface == |
81 |
|
82 |
jtlo = mybylo(mythid) |
83 |
jthi = mybyhi(mythid) |
84 |
itlo = mybxlo(mythid) |
85 |
ithi = mybxhi(mythid) |
86 |
jmin = 1 |
87 |
jmax = sny |
88 |
imin = 1 |
89 |
imax = snx |
90 |
|
91 |
c-- Read state record from global file. |
92 |
doglobalread = .false. |
93 |
ladinit = .false. |
94 |
|
95 |
irec = 1 |
96 |
|
97 |
#ifdef ALLOW_KAPGM_COST_CONTRIBUTION |
98 |
|
99 |
if (optimcycle .ge. 0) then |
100 |
ilfld = ilnblnk( xx_kapgm_file ) |
101 |
write(fnamefld(1:80),'(2a,i10.10)') |
102 |
& xx_kapgm_file(1:ilfld),'.',optimcycle |
103 |
endif |
104 |
|
105 |
fcthread = 0. _d 0 |
106 |
|
107 |
call active_read_xyz( fnamefld, tmpfld3d, irec, doglobalread, |
108 |
& ladinit, optimcycle, mythid |
109 |
& , xx_kapgm_dummy ) |
110 |
|
111 |
c-- Loop over this thread's tiles. |
112 |
do bj = jtlo,jthi |
113 |
do bi = itlo,ithi |
114 |
|
115 |
c-- Determine the weights to be used. |
116 |
|
117 |
fctile = 0. _d 0 |
118 |
do k = 1,nr |
119 |
do j = jmin,jmax |
120 |
do i = imin,imax |
121 |
if (_hFacC(i,j,k,bi,bj) .ne. 0.) then |
122 |
tmpx = (tmpfld3d(i,j,k,bi,bj)-GM_background_K) |
123 |
fctile = fctile |
124 |
& + 1e-8*cosphi(i,j,bi,bj) |
125 |
& *tmpx*tmpx |
126 |
endif |
127 |
enddo |
128 |
enddo |
129 |
enddo |
130 |
|
131 |
objf_kapgm(bi,bj) = objf_kapgm(bi,bj) + fctile |
132 |
fcthread = fcthread + fctile |
133 |
|
134 |
#ifdef ECCO_VERBOSE |
135 |
c-- Print cost function for each tile in each thread. |
136 |
write(msgbuf,'(a)') ' ' |
137 |
call print_message( msgbuf, standardmessageunit, |
138 |
& SQUEEZE_RIGHT , mythid) |
139 |
write(msgbuf,'(a,i8.8,1x,i3.3,1x,i3.3)') |
140 |
& ' cost_kapgm: irec,bi,bj = ',irec,bi,bj |
141 |
call print_message( msgbuf, standardmessageunit, |
142 |
& SQUEEZE_RIGHT , mythid) |
143 |
write(msgbuf,'(a,d22.15)') |
144 |
& ' cost function (dT(0)) = ', |
145 |
& fctile |
146 |
call print_message( msgbuf, standardmessageunit, |
147 |
& SQUEEZE_RIGHT , mythid) |
148 |
#endif |
149 |
enddo |
150 |
enddo |
151 |
|
152 |
#ifdef ECCO_VERBOSE |
153 |
c-- Print cost function for all tiles. |
154 |
_GLOBAL_SUM_R8( fcthread , myThid ) |
155 |
write(msgbuf,'(a)') ' ' |
156 |
call print_message( msgbuf, standardmessageunit, |
157 |
& SQUEEZE_RIGHT , mythid) |
158 |
write(msgbuf,'(a,i8.8)') |
159 |
& ' cost_: irec = ',irec |
160 |
call print_message( msgbuf, standardmessageunit, |
161 |
& SQUEEZE_RIGHT , mythid) |
162 |
write(msgbuf,'(a,d22.15)') |
163 |
& ' global cost function value = ', |
164 |
& fcthread |
165 |
call print_message( msgbuf, standardmessageunit, |
166 |
& SQUEEZE_RIGHT , mythid) |
167 |
write(msgbuf,'(a)') ' ' |
168 |
call print_message( msgbuf, standardmessageunit, |
169 |
& SQUEEZE_RIGHT , mythid) |
170 |
#endif |
171 |
|
172 |
#else |
173 |
c-- Do not enter the calculation of the salinity increment |
174 |
c-- contribution to the final cost function. |
175 |
|
176 |
fctile = 0. _d 0 |
177 |
fcthread = 0. _d 0 |
178 |
|
179 |
#ifdef ECCO_VERBOSE |
180 |
_BEGIN_MASTER( mythid ) |
181 |
write(msgbuf,'(a)') ' ' |
182 |
call print_message( msgbuf, standardmessageunit, |
183 |
& SQUEEZE_RIGHT , mythid) |
184 |
write(msgbuf,'(a,a)') |
185 |
& ' cost_kapgm : no contribution of the I.C. in salin. ', |
186 |
& ' to cost function.' |
187 |
call print_message( msgbuf, standardmessageunit, |
188 |
& SQUEEZE_RIGHT , mythid) |
189 |
write(msgbuf,'(a)') ' ' |
190 |
call print_message( msgbuf, standardmessageunit, |
191 |
& SQUEEZE_RIGHT , mythid) |
192 |
_END_MASTER( mythid ) |
193 |
#endif |
194 |
|
195 |
#endif |
196 |
|
197 |
return |
198 |
end |
199 |
|
200 |
|