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

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

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


Revision 1.2 - (show annotations) (download)
Thu Oct 28 18:20:24 2010 UTC (13 years, 6 months ago) by gforget
Branch: MAIN
CVS Tags: checkpoint63h, checkpoint63i, checkpoint63j, checkpoint63d, checkpoint63e, checkpoint63f, checkpoint63g, checkpoint63a, checkpoint63b, checkpoint63c, checkpoint63, checkpoint62o, checkpoint62n, checkpoint62s, checkpoint62r, checkpoint62q, checkpoint62p, checkpoint62w, checkpoint62v, checkpoint62u, checkpoint62t, checkpoint62z, checkpoint62y, checkpoint62x
Changes since 1.1: +3 -3 lines
bug fix.

1 C $Header: /u/gcmpack/MITgcm/pkg/ecco/cost_gencost_sstv4.F,v 1.1 2010/09/07 21:20:39 gforget Exp $
2 C $Name: $
3
4 #include "COST_CPPOPTIONS.h"
5
6
7 subroutine cost_gencost_sstv4(
8 I myiter,
9 I mytime,
10 I mythid
11 & )
12
13 c ==================================================================
14 c SUBROUTINE cost_gencost_sstv4
15 c ==================================================================
16 c
17 c o Evaluate cost function contributions of sea surface temperature.
18 c (Daily Pointwise and then Large Scale)
19 c
20 c started: Gael Forget, Oct-2009
21 c
22 c ==================================================================
23 c SUBROUTINE cost_gencost_sstv4
24 c ==================================================================
25
26 implicit none
27
28 c == global variables ==
29
30 #include "EEPARAMS.h"
31 #include "SIZE.h"
32 #include "PARAMS.h"
33 #include "GRID.h"
34
35 #include "ecco_cost.h"
36 #include "ctrl.h"
37 #include "ctrl_dummy.h"
38 #include "optim.h"
39 #include "DYNVARS.h"
40 #include "cal.h"
41
42 c == routine arguments ==
43
44 integer myiter
45 _RL mytime
46 integer mythid
47
48 #ifdef ALLOW_DAILYSST_COST_CONTRIBUTION
49 #ifdef ALLOW_SMOOTH
50 #ifdef ALLOW_GENCOST_CONTRIBUTION
51 #ifdef ALLOW_GENCOST_FREEFORM
52 c == local variables ==
53
54 integer bi,bj
55 integer i,j,k
56 integer itlo,ithi
57 integer jtlo,jthi
58 integer jmin,jmax
59 integer imin,imax
60 integer irec,jrec,krec
61 integer ilps
62 integer gwunit
63
64 logical doglobalread
65 logical ladinit
66
67 _RL anom_sst(1-olx:snx+olx,1-oly:sny+oly,nsx,nsy)
68 _RL obs_sst (1-olx:snx+olx,1-oly:sny+oly,nsx,nsy)
69 _RL nb_sst (1-olx:snx+olx,1-oly:sny+oly,nsx,nsy)
70 _RL msk_sst (1-olx:snx+olx,1-oly:sny+oly,nsx,nsy)
71 _RL tmp_sst (1-olx:snx+olx,1-oly:sny+oly,nsx,nsy)
72 _RL spval
73
74 integer num_var
75
76 _RL junk,junkweight
77
78 integer ndaysave
79 _RL ndaysaveRL
80
81 character*(80) fname
82 character*(80) fname2
83 character*(MAX_LEN_MBUF) msgbuf
84
85 _RL daytime
86 _RL diffsecs
87 integer il, localrec
88 integer dayiter
89 integer daydate(4)
90 integer difftime(4)
91 integer middate(4)
92 integer yday, ymod
93 integer md, dd, sd, ld, wd
94 integer mody, modm
95
96 integer igen_amsre, igen_amsre_lsc
97
98 c == external functions ==
99
100 integer ilnblnk
101 external ilnblnk
102
103 c == end of interface ==
104
105 jtlo = mybylo(mythid)
106 jthi = mybyhi(mythid)
107 itlo = mybxlo(mythid)
108 ithi = mybxhi(mythid)
109 jmin = 1
110 jmax = sny
111 imin = 1
112 imax = snx
113
114 c-- detect the relevant gencost indices
115 igen_amsre=0
116 igen_amsre_lsc=0
117 do k=1,NGENCOST
118 if (gencost_name(k).EQ.'sstv4-amsre') igen_amsre=k
119 if (gencost_name(k).EQ.'sstv4-amsre-lsc') igen_amsre_lsc=k
120 enddo
121
122 c-- First, read tiled data.
123 doglobalread = .false.
124 ladinit = .false.
125
126 write(fname(1:80),'(80a)') ' '
127 ilps=ilnblnk( sstbarfile )
128 write(fname(1:80),'(2a,i10.10)')
129 & sstbarfile(1:ilps),'.',optimcycle
130
131 spval = -9990.
132
133 cgf =======================================================
134 cgf PART 1: compute smooth SST cost term
135 cgf =======================================================
136
137
138 ndaysave=7
139 ndaysaveRL=ndaysave
140
141 do irec = 1, ndaysrec-ndaysave+1, 7
142
143 do bj = jtlo,jthi
144 do bi = itlo,ithi
145 do j = jmin,jmax
146 do i = imin,imax
147 anom_sst(i,j,bi,bj) = 0. _d 0
148 obs_sst(i,j,bi,bj) = 0. _d 0
149 nb_sst(i,j,bi,bj) = 0. _d 0
150 msk_sst(i,j,bi,bj) = 0. _d 0
151 enddo
152 enddo
153 enddo
154 enddo
155
156 c PART 1.1: compute running sample average over ndaysave
157 c ------------------------------------------------------
158
159 do jrec=1,ndaysave
160
161 krec=irec+jrec-1
162
163 c get modeled sst:
164 call active_read_xy( fname, sstbar, krec, doglobalread,
165 & ladinit, optimcycle, mythid,
166 & xx_sstbar_mean_dummy )
167
168 c get observed sst:
169 daytime = FLOAT(secondsperday*(krec-1))
170 dayiter = hoursperday*(krec-1)
171 call cal_getdate( dayiter, daytime, daydate, mythid )
172 call cal_convdate( daydate,yday,md,dd,sd,ld,wd,mythid )
173 ymod = modelstartdate(1)/10000
174 if ( ymod .EQ. yday ) then
175 middate(1) = modelstartdate(1)
176 else
177 middate(1) = yday*10000+100+1
178 endif
179 middate(2) = 0
180 middate(3) = modelstartdate(3)
181 middate(4) = modelstartdate(4)
182 call cal_TimePassed( middate, daydate, difftime, mythid )
183 call cal_ToSeconds( difftime, diffsecs, mythid )
184 localrec = int(diffsecs/86400.) + 1
185
186 il=ilnblnk(gencost_datafile(igen_amsre))
187 write(fname2(1:80),'(2a,i4)')
188 & gencost_datafile(igen_amsre)(1:il), '_', yday
189
190 if ( localrec .GT. 0 ) then
191 call mdsreadfield( fname2, cost_iprec, cost_yftype, 1,
192 & tmp_sst, localrec, mythid )
193 else
194 do bj = jtlo,jthi
195 do bi = itlo,ithi
196 do j = jmin,jmax
197 do i = imin,imax
198 tmp_sst(i,j,bi,bj) = spval
199 enddo
200 enddo
201 enddo
202 enddo
203 endif
204
205 c accumulate obs and misfit:
206 do bj = jtlo,jthi
207 do bi = itlo,ithi
208 do j = jmin,jmax
209 do i = imin,imax
210 if ( (tmp_sst(i,j,bi,bj).GT.spval).AND.
211 & (maskc(i,j,1,bi,bj).EQ.1.) ) then
212 anom_sst(i,j,bi,bj)= anom_sst(i,j,bi,bj)+
213 & sstbar(i,j,bi,bj)-tmp_sst(i,j,bi,bj)
214 obs_sst(i,j,bi,bj)= obs_sst(i,j,bi,bj)+
215 & tmp_sst(i,j,bi,bj)
216 nb_sst(i,j,bi,bj)=nb_sst(i,j,bi,bj)+1. _d 0
217 endif
218 enddo
219 enddo
220 enddo
221 enddo
222
223 enddo !do jrec=1,ndaysave
224
225 c average obs and misfit:
226 do bj = jtlo,jthi
227 do bi = itlo,ithi
228 do j = jmin,jmax
229 do i = imin,imax
230 if ( nb_sst(i,j,bi,bj) .NE. 0. ) then
231 obs_sst(i,j,bi,bj) =
232 & obs_sst(i,j,bi,bj)/nb_sst(i,j,bi,bj)
233 anom_sst(i,j,bi,bj) =
234 & anom_sst(i,j,bi,bj)/nb_sst(i,j,bi,bj)
235 msk_sst(i,j,bi,bj) = 1. _d 0
236 endif
237 enddo
238 enddo
239 enddo
240 enddo
241
242
243 c PART 1.2: smooth anom_sst in space
244 c ----------------------------------------
245
246 #ifdef ALLOW_GENCOST_SSTV4_OUTPUT
247 write(fname2(1:80),'(1a)') 'sstdiff_raw'
248 call mdswritefield(fname2,32,.false.,'RL',
249 & 1,anom_sst,irec,1,mythid)
250
251 write(fname2(1:80),'(1a)') 'sstobs_raw'
252 call mdswritefield(fname2,32,.false.,'RL',
253 & 1,obs_sst,irec,1,mythid)
254 #endif
255
256 call smooth_hetero2d(anom_sst,maskc,
257 & gencost_scalefile(igen_amsre_lsc),300,mythid)
258
259 #ifdef ALLOW_GENCOST_SSTV4_OUTPUT
260 call smooth_hetero2d(obs_sst,maskc,
261 & gencost_scalefile(igen_amsre_lsc),300,mythid)
262
263 write(fname2(1:80),'(1a)') 'sstdiff_smooth'
264 call mdswritefield(fname2,32,.false.,'RL',
265 & 1,anom_sst,irec,1,mythid)
266
267 write(fname2(1:80),'(1a)') 'sstobs_smooth'
268 call mdswritefield(fname2,32,.false.,'RL',
269 & 1,obs_sst,irec,1,mythid)
270 #endif
271
272 c PART 1.3: compute cost function term
273 c ------------------------------------
274
275 do bj = jtlo,jthi
276 do bi = itlo,ithi
277 do j = jmin,jmax
278 do i = imin,imax
279 junk = anom_sst(i,j,bi,bj)
280 junkweight = gencost_weight(i,j,bi,bj,igen_amsre_lsc)*
281 & maskc(i,j,1,bi,bj)
282 objf_gencost(igen_amsre_lsc,bi,bj) =
283 & objf_gencost(igen_amsre_lsc,bi,bj)
284 & +junk*junk*junkweight/ndaysaveRL
285 if ( (junkweight.GT.0.).AND.(nb_sst(i,j,bi,bj).GT.0.) )
286 & num_gencost(igen_amsre_lsc,bi,bj) =
287 & num_gencost(igen_amsre_lsc,bi,bj) + 1. _d 0 /ndaysaveRL
288 enddo
289 enddo
290 enddo
291 enddo
292
293 enddo
294
295 cgf =======================================================
296 cgf PART 2: compute raw SST cost term
297 cgf =======================================================
298
299 do irec = 1, ndaysrec
300
301 c get modeled sst:
302 call active_read_xy( fname, sstbar, irec, doglobalread,
303 & ladinit, optimcycle, mythid,
304 & xx_sstbar_mean_dummy )
305
306 c get observed sst:
307 daytime = FLOAT(secondsperday*(irec-1))
308 dayiter = hoursperday*(irec-1)
309 call cal_getdate( dayiter, daytime, daydate, mythid )
310 call cal_convdate( daydate,yday,md,dd,sd,ld,wd,mythid )
311 ymod = modelstartdate(1)/10000
312 if ( ymod .EQ. yday ) then
313 middate(1) = modelstartdate(1)
314 else
315 middate(1) = yday*10000+100+1
316 endif
317 middate(2) = 0
318 middate(3) = modelstartdate(3)
319 middate(4) = modelstartdate(4)
320 call cal_TimePassed( middate, daydate, difftime, mythid )
321 call cal_ToSeconds( difftime, diffsecs, mythid )
322 localrec = int(diffsecs/86400.) + 1
323
324 il=ilnblnk(gencost_datafile(igen_amsre))
325 write(fname2(1:80),'(2a,i4)')
326 & gencost_datafile(igen_amsre)(1:il), '_', yday
327
328 if ( localrec .GT. 0 ) then
329 call mdsreadfield( fname2, cost_iprec, cost_yftype, 1,
330 & tmp_sst, localrec, mythid )
331 else
332 do bj = jtlo,jthi
333 do bi = itlo,ithi
334 do j = jmin,jmax
335 do i = imin,imax
336 tmp_sst(i,j,bi,bj) = spval
337 enddo
338 enddo
339 enddo
340 enddo
341 endif
342
343 c compute misfit:
344 do bj = jtlo,jthi
345 do bi = itlo,ithi
346 do j = jmin,jmax
347 do i = imin,imax
348 if ( (tmp_sst(i,j,bi,bj).GT.spval).AND.
349 & (maskc(i,j,1,bi,bj).EQ.1.) ) then
350 anom_sst(i,j,bi,bj) =
351 & sstbar(i,j,bi,bj)-tmp_sst(i,j,bi,bj)
352 msk_sst(i,j,bi,bj) = 1. _d 0
353 else
354 anom_sst(i,j,bi,bj) = 0. _d 0
355 msk_sst(i,j,bi,bj) = 0. _d 0
356 endif
357 enddo
358 enddo
359 enddo
360 enddo
361
362 #ifdef ALLOW_GENCOST_SSTV4_OUTPUT
363 write(fname2(1:80),'(1a)') 'sstdiff_point'
364 call mdswritefield(fname2,32,.false.,'RL',
365 & 1,anom_sst,irec,1,mythid)
366 #endif
367
368 c compute cost:
369
370 do bj = jtlo,jthi
371 do bi = itlo,ithi
372 do j = jmin,jmax
373 do i = imin,imax
374 junk = anom_sst(i,j,bi,bj)
375 junkweight = gencost_weight(i,j,bi,bj,igen_amsre)*
376 & maskc(i,j,1,bi,bj)*msk_sst(i,j,bi,bj)
377 objf_gencost(igen_amsre,bi,bj) =
378 & objf_gencost(igen_amsre,bi,bj)+junk*junk*junkweight
379 if (junkweight.GT.0.)
380 & num_gencost(igen_amsre,bi,bj) =
381 & num_gencost(igen_amsre,bi,bj) + 1. _d 0
382 enddo
383 enddo
384 enddo
385 enddo
386
387 enddo
388
389 #endif /* ifdef ALLOW_GENCOST_FREEFORM */
390 #endif /* ifdef ALLOW_GENCOST_CONTRIBUTION */
391 #endif /* ifdef ALLOW_SMOOTH */
392 #endif /* ifdef ALLOW_DAILYSST_COST_CONTRIBUTION */
393
394 end

  ViewVC Help
Powered by ViewVC 1.1.22