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