/[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.4 - (show annotations) (download)
Thu May 5 23:47:57 2005 UTC (19 years, 1 month ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint57h_done, checkpoint57h_pre, checkpoint57h_post
Changes since 1.3: +32 -9 lines
o debugged and taffed geoid/sphere code for single CPU

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

  ViewVC Help
Powered by ViewVC 1.1.22