/[MITgcm]/MITgcm/pkg/cost/cost_ssh.F
ViewVC logotype

Contents of /MITgcm/pkg/cost/cost_ssh.F

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


Revision 1.1.2.3 - (show annotations) (download)
Thu Apr 4 10:58:59 2002 UTC (22 years, 2 months ago) by heimbach
Branch: ecco-branch
CVS Tags: ecco_c51_e34d, ecco_c51_e34e, ecco_c51_e34f, ecco_c51_e34g, ecco_c51_e34a, ecco_c51_e34b, ecco_c51_e34c, ecco_c51_e34, ecco_c50_e29, ecco_c50_e28, ecco_c44_e22, ecco_c44_e23, ecco_c44_e26, ecco_c44_e27, ecco_c44_e24, ecco_c44_e25, icebear5, icebear4, icebear3, icebear2, ecco_c50_e33a, ecco_c50_e32, ecco_c50_e33, ecco_c50_e30, ecco_c50_e31, ecco_ice2, ecco_ice1
Branch point for: icebear, c24_e25_ice
Changes since 1.1.2.2: +14 -1 lines
Merged new ECCO cost function terms into repository.
Various additions, modifications, quality checks.

1 C $Header: /u/gcmpack/MITgcm/pkg/cost/Attic/cost_ssh.F,v 1.1.2.2 2002/02/15 22:27:47 heimbach Exp $
2
3 #include "COST_CPPOPTIONS.h"
4
5
6 subroutine cost_ssh(
7 I myiter,
8 I mytime,
9 I mythid
10 & )
11
12 c ==================================================================
13 c SUBROUTINE cost_ssh
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
19 c started: Detlef Stammer, Ralf Giering Jul-1996
20 c
21 c changed: Christian Eckert eckert@mit.edu 25-Feb-2000
22 c
23 c - Restructured the code in order to create a package
24 c for the MITgcmUV.
25 c
26 c changed: Ralf Giering Ralf.Giering@FastOpt.de 12-Jun-2001
27 c
28 c - totally rewrite for parallel processing
29 c
30 c ==================================================================
31 c SUBROUTINE cost_ssh
32 c ==================================================================
33
34 implicit none
35
36 c == global variables ==
37
38 #include "EEPARAMS.h"
39 #include "SIZE.h"
40 #include "PARAMS.h"
41
42 #include "cost.h"
43 #include "ctrl.h"
44 #include "ctrl_dummy.h"
45 #include "optim.h"
46
47 c == routine arguments ==
48
49 integer myiter
50 _RL mytime
51 integer mythid
52
53 #ifdef ALLOW_SSH_COST_CONTRIBUTION
54 c == local variables ==
55
56 integer bi,bj
57 integer i,j
58 integer itlo,ithi
59 integer jtlo,jthi
60 integer jmin,jmax
61 integer imin,imax
62 integer irec
63 integer ilps
64
65 logical doglobalread
66 logical ladinit
67
68 _RL offset
69 _RL erscost
70 _RL tpcost
71 _RL costmean
72 _RL offset_sum
73 _RL psmean ( 1-olx:snx+olx, 1-oly:sny+oly, nsx, nsy )
74 _RL wwwtp ( 1-olx:snx+olx, 1-oly:sny+oly )
75 _RL wwwers ( 1-olx:snx+olx, 1-oly:sny+oly )
76 _RL junk
77
78 character*(80) fname
79 #ifdef ECCO_VERBOSE
80 character*(MAX_LEN_MBUF) msgbuf
81 #endif
82
83 c == external functions ==
84
85 integer ilnblnk
86 external ilnblnk
87
88 c == end of interface ==
89
90 jtlo = mybylo(mythid)
91 jthi = mybyhi(mythid)
92 itlo = mybxlo(mythid)
93 ithi = mybxhi(mythid)
94 jmin = 1
95 jmax = sny
96 imin = 1
97 imax = snx
98
99 c-- Initialise local variables.
100 costmean = 0. _d 0
101
102 do j = jmin, jmax
103 do i = imin, imax
104 wwwtp(i,j) = 0. _d 0
105 wwwers(i,j) = 0. _d 0
106 enddo
107 enddo
108
109 c-- First, read tiled data.
110 doglobalread = .false.
111 ladinit = .false.
112
113 write(fname(1:80),'(80a)') ' '
114 ilps=ilnblnk( psbarfile )
115 write(fname(1:80),'(2a,i10.10)')
116 & psbarfile(1:ilps),'.',optimcycle
117
118 #ifdef ECCO_VERBOSE
119 write(msgbuf,'(a)') ' '
120 call print_message( msgbuf, standardmessageunit,
121 & SQUEEZE_RIGHT , mythid)
122 write(msgbuf,'(a,i8.8)')
123 & ' cost_SSH: number of records to process = ',ndaysrec
124 call print_message( msgbuf, standardmessageunit,
125 & SQUEEZE_RIGHT , mythid)
126 write(msgbuf,'(a)') ' '
127 call print_message( msgbuf, standardmessageunit,
128 & SQUEEZE_RIGHT , mythid)
129 #endif
130
131 c-- ============
132 c-- Mean values.
133 c-- ============
134
135 do bj = jtlo,jthi
136 do bi = itlo,ithi
137 do j = jmin,jmax
138 do i = imin,imax
139 psmean(i,j,bi,bj) = 0. _d 0
140 enddo
141 enddo
142 enddo
143 enddo
144
145 c-- Loop over records for the first time.
146 do irec = 1, ndaysrec
147
148 c-- Compute the mean over all psbar records.
149 call active_read_xy( fname, psbar, irec, doglobalread,
150 & ladinit, optimcycle, mythid,
151 & xx_psbar_mean_dummy )
152
153 do bj = jtlo,jthi
154 do bi = itlo,ithi
155 do j = jmin,jmax
156 do i = imin,imax
157 psmean(i,j,bi,bj) = psmean(i,j,bi,bj) +
158 & psbar(i,j,bi,bj)*frame(i,j)/float(ndaysrec)
159 enddo
160 enddo
161 enddo
162 enddo
163
164 enddo
165
166 call cost_ReadTopexMean( mythid )
167
168 c-- Compute and remove offset for current tile and sum over all
169 c-- tiles of this instance.
170 offset = 0. _d 0
171 offset_sum = 0. _d 0
172
173 c-- Sum over this thread's tiles.
174 do bj = jtlo,jthi
175 do bi = itlo,ithi
176 do j = 1,sny
177 do i = 1,snx
178 offset = offset +
179 & tpmeanmask(i,j,bi,bj)*cosphi(i,j,bi,bj)*
180 & (tpmean(i,j,bi,bj) - psmean(i,j,bi,bj))
181 offset_sum = offset_sum +
182 & tpmeanmask(i,j,bi,bj)*cosphi(i,j,bi,bj)
183 enddo
184 enddo
185 enddo
186 enddo
187
188 c-- Do a global summation.
189 _GLOBAL_SUM_R8( offset , mythid )
190 _GLOBAL_SUM_R8( offset_sum , mythid )
191
192 if (offset_sum .eq. 0.0) then
193 _BEGIN_MASTER( mythid )
194 print*
195 print*,' cost_ssh: offset_sum = zero!'
196 print*
197 _END_MASTER( mythid )
198 stop ' ... stopped in cost_ssh.'
199 else
200 _BEGIN_MASTER( mythid )
201 print*
202 print*,' cost_ssh: offset_sum = ',offset_sum
203 print*
204 _END_MASTER( mythid )
205 endif
206
207 offset = offset / offset_sum
208
209 #ifdef ALLOW_SSH_MEAN_COST_CONTRIBUTION
210
211 c-- ==========
212 c-- Mean
213 c-- ==========
214 c-- compute mean ssh difference cost contribution
215 call cost_ssh_mean(
216 I psmean, offset
217 O , costmean
218 I , mythid
219 & )
220
221
222 objf_hmean = costmean
223
224 #endif /* ALLOW_SSH_MEAN_COST_CONTRIBUTION */
225
226 c-- ==========
227 c-- Anomalies.
228 c-- ==========
229
230 erscost = 0. _d 0
231 tpcost = 0. _d 0
232
233 c-- Loop over records for the second time.
234 do irec = 1, ndaysrec
235
236 call active_read_xy( fname, psbar, irec, doglobalread,
237 & ladinit, optimcycle, mythid,
238 & xx_psbar_mean_dummy )
239
240 #ifdef ALLOW_SSH_TPANOM_COST_CONTRIBUTION
241 call cost_ReadTopex( irec, mythid )
242 #endif
243
244 #ifdef ALLOW_SSH_ERSANOM_COST_CONTRIBUTION
245 call cost_ReadERS( irec, mythid )
246 #endif
247
248 do bj = jtlo,jthi
249 do bi = itlo,ithi
250
251 erscost = 0. _d 0
252 tpcost = 0. _d 0
253
254 #ifdef ALLOW_SSH_TPANOM_COST_CONTRIBUTION
255 do j = jmin,jmax
256 do i = imin,imax
257 c-- The array psobs contains SSH anomalies.
258 wwwtp(i,j) = wtp(i,j,bi,bj) *cosphi(i,j,bi,bj)
259 junk = ((psbar(i,j,bi,bj) - psmean(i,j,bi,bj)) -
260 & tpobs(i,j,bi,bj))*tpmask(i,j,bi,bj)
261 tpcost = tpcost + junk*junk*wwwtp(i,j)
262 enddo
263 enddo
264 #endif
265
266 #ifdef ALLOW_SSH_ERSANOM_COST_CONTRIBUTION
267 do j = jmin,jmax
268 do i = imin,imax
269 c-- The array ersobs contains SSH anomalies.
270 wwwers(i,j) = wers(i,j,bi,bj)*cosphi(i,j,bi,bj)
271 junk = ((psbar(i,j,bi,bj) - psmean(i,j,bi,bj)) -
272 & ersobs(i,j,bi,bj))*ersmask(i,j,bi,bj)
273 erscost = erscost + junk*junk*wwwers(i,j)
274 enddo
275
276 enddo
277 #endif
278
279 objf_h(bi,bj) = objf_h(bi,bj) + tpcost + erscost
280
281 #ifdef ECCO_VERBOSE
282 write(msgbuf,'(a)') ' '
283 call print_message( msgbuf, standardmessageunit,
284 & SQUEEZE_RIGHT , mythid)
285 write(msgbuf,'(a,i8.8,1x,i3.3,1x,i3.3)')
286 & ' COST_SSH: irec,bi,bj = ',irec,bi,bj
287 call print_message( msgbuf, standardmessageunit,
288 & SQUEEZE_RIGHT , mythid)
289 write(msgbuf,'(a,d22.15)')
290 & ' COST_SSH: cost function (TOPEX) = ',tpcost
291 call print_message( msgbuf, standardmessageunit,
292 & SQUEEZE_RIGHT , mythid)
293 write(msgbuf,'(a,d22.15)')
294 & ' COST_SSH: cost function (ERS) = ',erscost
295 call print_message( msgbuf, standardmessageunit,
296 & SQUEEZE_RIGHT , mythid)
297 write(msgbuf,'(a,d22.15)')
298 & ' COST_SSH: cost function (Total) = ',objf_h(bi,bj)
299 call print_message( msgbuf, standardmessageunit,
300 & SQUEEZE_RIGHT , mythid)
301 write(msgbuf,'(a)') ' '
302 call print_message( msgbuf, standardmessageunit,
303 & SQUEEZE_RIGHT , mythid)
304 #endif
305
306 enddo
307 enddo
308
309 enddo
310 c-- End of second loop over records.
311
312 #endif /* ifdef ALLOW_SSH_COST_CONTRIBUTION */
313
314 end

  ViewVC Help
Powered by ViewVC 1.1.22