1 |
jmc |
1.12 |
C $Header: /u/gcmpack/MITgcm/pkg/ecco/cost_ssh_mean.F,v 1.11 2012/08/10 19:45:27 jmc Exp $ |
2 |
jmc |
1.6 |
C $Name: $ |
3 |
heimbach |
1.1 |
|
4 |
jmc |
1.11 |
#include "ECCO_OPTIONS.h" |
5 |
jmc |
1.10 |
#define ORIGINAL_COST_SSH_MEAN |
6 |
heimbach |
1.1 |
|
7 |
|
|
subroutine cost_ssh_mean( |
8 |
|
|
I psmean, offset |
9 |
|
|
O , costmean |
10 |
|
|
I , mythid |
11 |
|
|
& ) |
12 |
|
|
|
13 |
|
|
c ================================================================== |
14 |
|
|
c SUBROUTINE cost_ssh_mean |
15 |
|
|
c ================================================================== |
16 |
|
|
c |
17 |
|
|
c o Evaluate cost function contribution of sea surface height. |
18 |
|
|
c using of geoid error covariances requires regular model grid |
19 |
|
|
c o TODO: interpolate model grid to regular grid |
20 |
|
|
c check mask |
21 |
|
|
c |
22 |
|
|
c started: Detlef Stammer, Ralf Giering Jul-1996 |
23 |
|
|
c changed: Ralf Giering Ralf.Giering@FastOpt.de 12-Jun-2001 |
24 |
|
|
c - totally rewrite for parallel processing |
25 |
|
|
c heimbach@mit.edu 13-Mar-2002 |
26 |
|
|
c - several wrong i-loop boundaries (spotted by G. Gebbie) |
27 |
|
|
c - nprocs need to be replaced by nPx*nPy |
28 |
|
|
c - geoid error does not work as of now |
29 |
heimbach |
1.4 |
c heimbach@mit.edu 05-May-2005 |
30 |
|
|
c - debugged and restructuted |
31 |
heimbach |
1.1 |
c |
32 |
|
|
c ================================================================== |
33 |
|
|
c SUBROUTINE cost_ssh_mean |
34 |
|
|
c ================================================================== |
35 |
|
|
|
36 |
|
|
implicit none |
37 |
|
|
|
38 |
|
|
c == global variables == |
39 |
|
|
|
40 |
|
|
#include "EEPARAMS.h" |
41 |
|
|
#include "SIZE.h" |
42 |
|
|
#include "PARAMS.h" |
43 |
|
|
|
44 |
|
|
#include "ecco_cost.h" |
45 |
heimbach |
1.5 |
#include "optim.h" |
46 |
heimbach |
1.1 |
|
47 |
|
|
#ifdef ALLOW_EGM96_ERROR_COV |
48 |
jmc |
1.12 |
# ifdef ALLOW_USE_MPI |
49 |
heimbach |
1.4 |
# include "EESUPPORT.h" |
50 |
|
|
# else |
51 |
|
|
INTEGER nProcs |
52 |
|
|
PARAMETER (nProcs=1) |
53 |
|
|
# endif |
54 |
|
|
# include "sphere.h" |
55 |
heimbach |
1.1 |
#endif |
56 |
|
|
|
57 |
|
|
c == routine arguments == |
58 |
|
|
|
59 |
|
|
_RL psmean ( 1-olx:snx+olx, 1-oly:sny+oly, nsx, nsy ) |
60 |
|
|
_RL offset |
61 |
|
|
_RL costmean |
62 |
|
|
integer mythid |
63 |
|
|
|
64 |
|
|
c == local variables == |
65 |
|
|
|
66 |
|
|
integer bi,bj |
67 |
|
|
integer i,j |
68 |
|
|
integer itlo,ithi |
69 |
|
|
integer jtlo,jthi |
70 |
|
|
integer jmin,jmax |
71 |
|
|
integer imin,imax |
72 |
|
|
|
73 |
heimbach |
1.5 |
cph( |
74 |
|
|
_RL diagnosfld(1-olx:snx+olx,1-oly:sny+oly,nsx,nsy) |
75 |
|
|
cph) |
76 |
|
|
|
77 |
heimbach |
1.1 |
#ifdef ALLOW_EGM96_ERROR_COV |
78 |
|
|
_RL cphi |
79 |
|
|
_RL xmean |
80 |
|
|
_RL shc( ncshc ) |
81 |
|
|
_RL misfit ( 1-olx:snx+olx, 1-oly:sny+oly, nsx,nsy ) |
82 |
|
|
_RL misfitgl( 1-olx:snx+olx, 1-oly:sny+oly, nsx,nsy, npx,npy ) |
83 |
|
|
|
84 |
|
|
_RL pwork ( (lshcmax+1)*(lshcmax+2)/2 ) |
85 |
|
|
_RL diffearth ( lonshc, latshc ) |
86 |
|
|
integer joffset |
87 |
|
|
integer ipx, ipy |
88 |
|
|
integer iglobe, jglobe |
89 |
|
|
integer iproc |
90 |
|
|
integer mpirc |
91 |
|
|
integer mpicrd(2) |
92 |
|
|
#endif |
93 |
|
|
|
94 |
|
|
_RL diff |
95 |
jmc |
1.8 |
_RL sumc |
96 |
heimbach |
1.4 |
character*(max_len_mbuf) msgbuf |
97 |
jmc |
1.8 |
_RL sumC0, sumC1, sumC2 |
98 |
heimbach |
1.1 |
|
99 |
|
|
c == end of interface == |
100 |
|
|
|
101 |
|
|
c-- Initialise local variables. |
102 |
|
|
jtlo = mybylo(mythid) |
103 |
|
|
jthi = mybyhi(mythid) |
104 |
|
|
itlo = mybxlo(mythid) |
105 |
|
|
ithi = mybxhi(mythid) |
106 |
|
|
|
107 |
|
|
jmin = 1 |
108 |
|
|
jmax = sny |
109 |
|
|
imin = 1 |
110 |
|
|
imax = snx |
111 |
|
|
|
112 |
|
|
costmean = 0. |
113 |
|
|
|
114 |
|
|
#ifdef ALLOW_EGM96_ERROR_COV |
115 |
|
|
|
116 |
|
|
c-- compute misfit on local domain |
117 |
|
|
do bj = jtlo,jthi |
118 |
|
|
do bi = itlo,ithi |
119 |
|
|
do j = jmin,jmax |
120 |
|
|
do i = imin,imax |
121 |
|
|
if ( tpmeanmask(i,j,bi,bj) .gt. 0. ) then |
122 |
|
|
misfit(i,j,bi,bj) = |
123 |
|
|
& tpmean(i,j,bi,bj) - offset - psmean(i,j,bi,bj) |
124 |
|
|
else |
125 |
|
|
misfit(i,j,bi,bj) = 0. |
126 |
|
|
endif |
127 |
|
|
enddo |
128 |
|
|
enddo |
129 |
|
|
enddo |
130 |
|
|
enddo |
131 |
|
|
|
132 |
|
|
c-- communicate to get misfit on full 2d model surface grid |
133 |
jmc |
1.12 |
#ifdef ALLOW_USE_MPI |
134 |
heimbach |
1.1 |
call exch_all_2d_rl( misfit, misfitgl, mythid ) |
135 |
heimbach |
1.4 |
#else |
136 |
|
|
do bj = jtlo,jthi |
137 |
|
|
do bi = itlo,ithi |
138 |
|
|
do j = jmin,jmax |
139 |
|
|
do i = imin,imax |
140 |
|
|
misfitgl(i,j,bi,bj,1,1) = misfit(i,j,bi,bj) |
141 |
|
|
enddo |
142 |
|
|
enddo |
143 |
|
|
enddo |
144 |
|
|
enddo |
145 |
|
|
#endif |
146 |
heimbach |
1.1 |
|
147 |
|
|
c-- set meridional index offset, |
148 |
|
|
c-- it is non zero if the grid does no reach the poles |
149 |
|
|
joffset = (latshc - ny)/2 |
150 |
|
|
|
151 |
heimbach |
1.4 |
write(msgbuf,'(a,I10)') |
152 |
|
|
& 'cost_ssh_mean: grid starts at point ', joffset |
153 |
|
|
call print_message( msgbuf, standardmessageunit, |
154 |
|
|
& SQUEEZE_RIGHT , mythid) |
155 |
heimbach |
1.1 |
|
156 |
|
|
c-- preset field |
157 |
|
|
c-- necessary if the grid does no reach the poles |
158 |
|
|
do j = 1, latshc |
159 |
|
|
do i = 1, lonshc |
160 |
|
|
diffearth(i,j) = 0. |
161 |
|
|
enddo |
162 |
|
|
enddo |
163 |
|
|
|
164 |
|
|
c-- -interpolate- from model grid to regular grid |
165 |
|
|
c-- so far we assume that the model grid is already regular |
166 |
|
|
_BEGIN_MASTER( mythid ) |
167 |
|
|
|
168 |
|
|
do iproc = 1, nprocs |
169 |
|
|
|
170 |
jmc |
1.12 |
#ifdef ALLOW_USE_MPI |
171 |
heimbach |
1.1 |
c-- get coordinates of processor (iporc-1) |
172 |
|
|
call MPI_Cart_coords( |
173 |
|
|
I MPI_COMM_MODEL, iproc-1, 2, mpicrd |
174 |
|
|
O , mpirc |
175 |
|
|
& ) |
176 |
|
|
ipx = 1 + mpicrd(1) |
177 |
|
|
ipy = 1 + mpicrd(2) |
178 |
heimbach |
1.4 |
#else |
179 |
|
|
ipx = 1 |
180 |
|
|
ipy = 1 |
181 |
|
|
#endif |
182 |
heimbach |
1.1 |
|
183 |
|
|
do bj = jtlo,jthi |
184 |
|
|
do bi = itlo,ithi |
185 |
|
|
|
186 |
|
|
do j = jmin,jmax |
187 |
|
|
do i = imin,imax |
188 |
|
|
jglobe = joffset+ j + sNy*(bj-1) + sNy*nSy*(ipy-1) |
189 |
|
|
iglobe = i + sNx*(bi-1) + sNx*nSx*(ipx-1) |
190 |
|
|
|
191 |
|
|
diffearth(iglobe,jglobe) = misfitgl(i,j,bi,bj,ipx,ipy) |
192 |
|
|
enddo |
193 |
|
|
enddo |
194 |
|
|
|
195 |
|
|
enddo |
196 |
|
|
enddo |
197 |
|
|
|
198 |
|
|
enddo |
199 |
|
|
|
200 |
|
|
c-- Project regular grid misfit onto sperical harmonics |
201 |
|
|
call shc4grid( |
202 |
|
|
I lshcmax |
203 |
|
|
O , shc |
204 |
|
|
I , latshc, lonshc |
205 |
|
|
I , diffearth |
206 |
|
|
W , pwork |
207 |
|
|
& ) |
208 |
|
|
|
209 |
|
|
c-- Remove the C(0,0) component, i.e. the global mean. |
210 |
|
|
shc(1) = 0. |
211 |
|
|
|
212 |
|
|
C-- Compute the cost function for the mean SSH. |
213 |
|
|
call cost_geoid( |
214 |
|
|
O costmean |
215 |
|
|
I , shc |
216 |
|
|
I , mythid |
217 |
|
|
& ) |
218 |
|
|
|
219 |
|
|
_END_MASTER( mythid ) |
220 |
|
|
|
221 |
|
|
_BARRIER |
222 |
|
|
|
223 |
jmc |
1.8 |
#else /* else ALLOW_EGM96_ERROR_COV */ |
224 |
heimbach |
1.1 |
|
225 |
|
|
c-- Compute cost function for SSH by using the diagonal |
226 |
|
|
c-- of the error covariance only. |
227 |
|
|
c-- |
228 |
|
|
c-- Note: wp is assumed to include latitude dependence |
229 |
|
|
c-- of error due to meridian convergence; |
230 |
|
|
c-- --> no weighting by cosphi. |
231 |
|
|
|
232 |
jmc |
1.10 |
#ifdef ORIGINAL_COST_SSH_MEAN |
233 |
heimbach |
1.1 |
sumc = 0. |
234 |
|
|
do bj = jtlo,jthi |
235 |
|
|
do bi = itlo,ithi |
236 |
|
|
do j = jmin,jmax |
237 |
|
|
do i = imin,imax |
238 |
|
|
diff = psmean(i,j,bi,bj) - tpmean(i,j,bi,bj) + offset |
239 |
|
|
sumc = sumc + diff*diff |
240 |
|
|
& * wp(i,j,bi,bj)*tpmeanmask(i,j,bi,bj) |
241 |
heimbach |
1.3 |
if ( wp(i,j,bi,bj)*tpmeanmask(i,j,bi,bj) .ne. 0. ) |
242 |
|
|
& num_hmean = num_hmean + 1. _d 0 |
243 |
jmc |
1.7 |
diagnosfld(i,j,bi,bj) = diff*diff |
244 |
heimbach |
1.5 |
& * wp(i,j,bi,bj)*tpmeanmask(i,j,bi,bj) |
245 |
heimbach |
1.1 |
enddo |
246 |
|
|
enddo |
247 |
|
|
enddo |
248 |
|
|
enddo |
249 |
|
|
|
250 |
jmc |
1.8 |
c-- Do the global summation. |
251 |
|
|
_GLOBAL_SUM_RL( sumc, mythid ) |
252 |
|
|
_GLOBAL_SUM_RL( num_hmean, mythid ) |
253 |
|
|
costmean = sumc |
254 |
|
|
|
255 |
jmc |
1.10 |
#else /* ORIGINAL_COST_SSH_MEAN */ |
256 |
|
|
C- Expand explicitly cost expression as quadratic function of "offset" |
257 |
|
|
C to avoid problems with adjoint of global_sum function. |
258 |
|
|
C Because of different truncation, expect also slight difference on 1 proc. |
259 |
|
|
sumC0 = 0. |
260 |
|
|
sumC1 = 0. |
261 |
|
|
sumC2 = 0. |
262 |
|
|
do bj = jtlo,jthi |
263 |
|
|
do bi = itlo,ithi |
264 |
|
|
do j = jmin,jmax |
265 |
|
|
do i = imin,imax |
266 |
|
|
diff = psmean(i,j,bi,bj) - tpmean(i,j,bi,bj) |
267 |
|
|
sumC0 = sumC0 + wp(i,j,bi,bj)*tpmeanmask(i,j,bi,bj) |
268 |
|
|
sumC1 = sumC1 + diff |
269 |
|
|
& * wp(i,j,bi,bj)*tpmeanmask(i,j,bi,bj) |
270 |
|
|
sumC2 = sumC2 + diff*diff |
271 |
|
|
& * wp(i,j,bi,bj)*tpmeanmask(i,j,bi,bj) |
272 |
|
|
if ( wp(i,j,bi,bj)*tpmeanmask(i,j,bi,bj) .ne. 0. ) |
273 |
|
|
& num_hmean = num_hmean + 1. _d 0 |
274 |
|
|
diagnosfld(i,j,bi,bj) = |
275 |
|
|
& ( diff*diff + 2. _d 0*offset*diff + offset*offset ) |
276 |
|
|
& * wp(i,j,bi,bj)*tpmeanmask(i,j,bi,bj) |
277 |
|
|
enddo |
278 |
|
|
enddo |
279 |
|
|
enddo |
280 |
|
|
enddo |
281 |
|
|
|
282 |
|
|
C-- Do the global summation. |
283 |
|
|
_GLOBAL_SUM_RL( sumC0, mythid ) |
284 |
|
|
_GLOBAL_SUM_RL( sumC1, mythid ) |
285 |
|
|
_GLOBAL_SUM_RL( sumC2, mythid ) |
286 |
|
|
_GLOBAL_SUM_RL( num_hmean, mythid ) |
287 |
|
|
|
288 |
|
|
costmean = sumC2 + 2. _d 0*offset*sumC1 + offset*offset*sumC0 |
289 |
|
|
|
290 |
|
|
#endif /* else ORIGINAL_COST_SSH_MEAN */ |
291 |
|
|
|
292 |
heimbach |
1.5 |
cph( |
293 |
|
|
CALL WRITE_FLD_XY_RL( 'DiagnosSSHmean', ' ', diagnosfld, |
294 |
|
|
& optimcycle, mythid ) |
295 |
|
|
cph) |
296 |
|
|
|
297 |
jmc |
1.8 |
#endif /* else ALLOW_EGM96_ERROR_COV */ |
298 |
heimbach |
1.1 |
|
299 |
jmc |
1.8 |
return |
300 |
heimbach |
1.1 |
end |