/[MITgcm]/MITgcm/pkg/ecco/cost_ssh_mean.F
ViewVC logotype

Contents of /MITgcm/pkg/ecco/cost_ssh_mean.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.7 - (show annotations) (download)
Tue Apr 28 17:55:53 2009 UTC (15 years ago) by jmc
Branch: MAIN
CVS Tags: checkpoint62, checkpoint62c, checkpoint62b, checkpoint62a, checkpoint62e, checkpoint62d, checkpoint61n, checkpoint61o, checkpoint61m, checkpoint61v, checkpoint61w, checkpoint61t, checkpoint61u, checkpoint61r, checkpoint61s, checkpoint61p, checkpoint61q, checkpoint61z, checkpoint61x, checkpoint61y
Changes since 1.6: +4 -4 lines
this is the only place I found where a _R8 macro is applied to a real*8
 variable; replace with direct call to GLOBAL_SUM_R8.

1 C $Header: /u/gcmpack/MITgcm/pkg/ecco/cost_ssh_mean.F,v 1.6 2007/10/09 00:02:51 jmc Exp $
2 C $Name: $
3
4 #include "COST_CPPOPTIONS.h"
5
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 #include "PARAMS.h"
43
44 #include "ecco_cost.h"
45 #include "optim.h"
46
47 #ifdef ALLOW_EGM96_ERROR_COV
48 # if (defined (ALLOW_USE_MPI) || defined (ALWAYS_USE_MPI))
49 # include "EESUPPORT.h"
50 # else
51 INTEGER nProcs
52 PARAMETER (nProcs=1)
53 # endif
54 # include "sphere.h"
55 #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 cph(
74 _RL diagnosfld(1-olx:snx+olx,1-oly:sny+oly,nsx,nsy)
75 cph)
76
77 #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 Real*8 sumc
96 character*(max_len_mbuf) msgbuf
97
98 c == end of interface ==
99
100 c-- Initialise local variables.
101 jtlo = mybylo(mythid)
102 jthi = mybyhi(mythid)
103 itlo = mybxlo(mythid)
104 ithi = mybxhi(mythid)
105
106 jmin = 1
107 jmax = sny
108 imin = 1
109 imax = snx
110
111 costmean = 0.
112
113 #ifdef ALLOW_EGM96_ERROR_COV
114
115 c-- compute misfit on local domain
116 do bj = jtlo,jthi
117 do bi = itlo,ithi
118 do j = jmin,jmax
119 do i = imin,imax
120 if ( tpmeanmask(i,j,bi,bj) .gt. 0. ) then
121 misfit(i,j,bi,bj) =
122 & tpmean(i,j,bi,bj) - offset - psmean(i,j,bi,bj)
123 else
124 misfit(i,j,bi,bj) = 0.
125 endif
126 enddo
127 enddo
128 enddo
129 enddo
130
131 c-- communicate to get misfit on full 2d model surface grid
132 #if (defined (ALLOW_USE_MPI) || defined (ALWAYS_USE_MPI))
133 call exch_all_2d_rl( misfit, misfitgl, mythid )
134 #else
135 do bj = jtlo,jthi
136 do bi = itlo,ithi
137 do j = jmin,jmax
138 do i = imin,imax
139 misfitgl(i,j,bi,bj,1,1) = misfit(i,j,bi,bj)
140 enddo
141 enddo
142 enddo
143 enddo
144 #endif
145
146 c-- set meridional index offset,
147 c-- it is non zero if the grid does no reach the poles
148 joffset = (latshc - ny)/2
149
150 write(msgbuf,'(a,I10)')
151 & 'cost_ssh_mean: grid starts at point ', joffset
152 call print_message( msgbuf, standardmessageunit,
153 & SQUEEZE_RIGHT , mythid)
154
155 c-- preset field
156 c-- necessary if the grid does no reach the poles
157 do j = 1, latshc
158 do i = 1, lonshc
159 diffearth(i,j) = 0.
160 enddo
161 enddo
162
163 c-- -interpolate- from model grid to regular grid
164 c-- so far we assume that the model grid is already regular
165 _BEGIN_MASTER( mythid )
166
167 do iproc = 1, nprocs
168
169 #if (defined (ALLOW_USE_MPI) || defined (ALWAYS_USE_MPI))
170 c-- get coordinates of processor (iporc-1)
171 call MPI_Cart_coords(
172 I MPI_COMM_MODEL, iproc-1, 2, mpicrd
173 O , mpirc
174 & )
175 ipx = 1 + mpicrd(1)
176 ipy = 1 + mpicrd(2)
177 #else
178 ipx = 1
179 ipy = 1
180 #endif
181
182 do bj = jtlo,jthi
183 do bi = itlo,ithi
184
185 do j = jmin,jmax
186 do i = imin,imax
187 jglobe = joffset+ j + sNy*(bj-1) + sNy*nSy*(ipy-1)
188 iglobe = i + sNx*(bi-1) + sNx*nSx*(ipx-1)
189
190 diffearth(iglobe,jglobe) = misfitgl(i,j,bi,bj,ipx,ipy)
191 enddo
192 enddo
193
194 enddo
195 enddo
196
197 enddo
198
199 c-- Project regular grid misfit onto sperical harmonics
200 call shc4grid(
201 I lshcmax
202 O , shc
203 I , latshc, lonshc
204 I , diffearth
205 W , pwork
206 & )
207
208 c-- Remove the C(0,0) component, i.e. the global mean.
209 shc(1) = 0.
210
211 C-- Compute the cost function for the mean SSH.
212 call cost_geoid(
213 O costmean
214 I , shc
215 I , mythid
216 & )
217
218 _END_MASTER( mythid )
219
220 _BARRIER
221
222 #else
223
224 c-- Compute cost function for SSH by using the diagonal
225 c-- of the error covariance only.
226 c--
227 c-- Note: wp is assumed to include latitude dependence
228 c-- of error due to meridian convergence;
229 c-- --> no weighting by cosphi.
230
231 sumc = 0.
232 do bj = jtlo,jthi
233 do bi = itlo,ithi
234 do j = jmin,jmax
235 do i = imin,imax
236 diff = psmean(i,j,bi,bj) - tpmean(i,j,bi,bj) + offset
237 sumc = sumc + diff*diff
238 & * wp(i,j,bi,bj)*tpmeanmask(i,j,bi,bj)
239 if ( wp(i,j,bi,bj)*tpmeanmask(i,j,bi,bj) .ne. 0. )
240 & num_hmean = num_hmean + 1. _d 0
241 diagnosfld(i,j,bi,bj) = diff*diff
242 & * wp(i,j,bi,bj)*tpmeanmask(i,j,bi,bj)
243 enddo
244 enddo
245 enddo
246 enddo
247
248 cph(
249 CALL WRITE_FLD_XY_RL( 'DiagnosSSHmean', ' ', diagnosfld,
250 & optimcycle, mythid )
251 cph)
252
253 c-- Do the global summation.
254 CALL GLOBAL_SUM_R8( sumc, mythid )
255 _GLOBAL_SUM_RL( num_hmean, mythid )
256 costmean = sumc
257
258 #endif
259
260 end

  ViewVC Help
Powered by ViewVC 1.1.22