/[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.10 - (show annotations) (download)
Tue Apr 1 08:03:42 2014 UTC (10 years, 2 months ago) by atn
Branch: MAIN
CVS Tags: checkpoint64y, checkpoint64x, checkpoint64z, checkpoint64w, checkpoint64v, checkpoint65, checkpoint65b, checkpoint65c, checkpoint65a, checkpoint65d, checkpoint65e
Changes since 1.9: +4 -2 lines
o move smooth2Ddiffnbt param from pkg/smooth to pkg/ecco
o minor bug fix in set gencost default values in ecco_readparms.F

1 C $Header: /u/gcmpack/MITgcm/pkg/ecco/cost_gencost_sstv4.F,v 1.8 2012/10/02 01:40:03 gforget Exp $
2 C $Name: $
3
4 #include "ECCO_OPTIONS.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_SIZE.h"
37 #include "ctrl.h"
38 #include "ctrl_dummy.h"
39 #include "optim.h"
40 #include "DYNVARS.h"
41 #include "cal.h"
42 #ifdef ALLOW_SMOOTH
43 #include "SMOOTH.h"
44 #endif
45
46 c == routine arguments ==
47
48 integer myiter
49 _RL mytime
50 integer mythid
51
52 #ifdef ALLOW_DAILYSST_COST_CONTRIBUTION
53 #ifdef ALLOW_SMOOTH
54 #ifdef ALLOW_GENCOST_CONTRIBUTION
55 c == local variables ==
56
57 integer bi,bj
58 integer i,j,k
59 integer itlo,ithi
60 integer jtlo,jthi
61 integer jmin,jmax
62 integer imin,imax
63 integer irec,jrec,krec
64 integer ilps
65 integer gwunit
66
67 logical doglobalread
68 logical ladinit
69
70 _RL anom_sst(1-olx:snx+olx,1-oly:sny+oly,nsx,nsy)
71 _RL obs_sst (1-olx:snx+olx,1-oly:sny+oly,nsx,nsy)
72 _RL nb_sst (1-olx:snx+olx,1-oly:sny+oly,nsx,nsy)
73 _RL msk_sst (1-olx:snx+olx,1-oly:sny+oly,nsx,nsy)
74 _RL tmp_sst (1-olx:snx+olx,1-oly:sny+oly,nsx,nsy)
75 _RL spval
76
77 integer num_var
78
79 _RL junk,junkweight
80
81 integer ndaysave
82 _RL ndaysaveRL
83
84 character*(80) fname
85 character*(80) fname2
86
87 #ifdef ALLOW_ECCO_DEBUG
88 character*(MAX_LEN_MBUF) msgBuf
89 INTEGER ioUnit
90 #endif
91 logical exst
92
93 _RL daytime
94 _RL diffsecs
95 integer il, localrec
96 integer dayiter
97 integer daydate(4)
98 integer difftime(4)
99 integer tempDate_1
100 integer middate(4)
101 integer locstartdate(4)
102 integer yday, ymod
103 integer md, dd, sd, ld, wd
104 integer mody, modm
105
106 integer igen_amsre, igen_amsre_lsc
107
108 c == external functions ==
109
110 integer ilnblnk
111 external ilnblnk
112
113 c == end of interface ==
114
115 jtlo = mybylo(mythid)
116 jthi = mybyhi(mythid)
117 itlo = mybxlo(mythid)
118 ithi = mybxhi(mythid)
119 jmin = 1
120 jmax = sny
121 imin = 1
122 imax = snx
123
124 #ifdef ALLOW_ECCO_DEBUG
125 ioUnit=standardMessageUnit
126 #endif
127
128 call cal_FullDate(19920101,0,locstartdate,mythid)
129
130 c-- detect the relevant gencost indices
131 igen_amsre=0
132 igen_amsre_lsc=0
133 do k=1,NGENCOST
134 if (gencost_name(k).EQ.'sstv4-amsre') igen_amsre=k
135 if (gencost_name(k).EQ.'sstv4-amsre-lsc') igen_amsre_lsc=k
136 enddo
137
138 c-- First, read tiled data.
139 doglobalread = .false.
140 ladinit = .false.
141
142 write(fname(1:80),'(80a)') ' '
143 ilps=ilnblnk( sstbarfile )
144 write(fname(1:80),'(2a,i10.10)')
145 & sstbarfile(1:ilps),'.',optimcycle
146
147 spval = -9990.
148
149 cgf =======================================================
150 cgf PART 1: compute smooth SST cost term
151 cgf =======================================================
152
153
154 ndaysave=7
155 ndaysaveRL=ndaysave
156
157 do irec = 1, ndaysrec-ndaysave+1, 7
158
159 do bj = jtlo,jthi
160 do bi = itlo,ithi
161 do j = jmin,jmax
162 do i = imin,imax
163 anom_sst(i,j,bi,bj) = 0. _d 0
164 obs_sst(i,j,bi,bj) = 0. _d 0
165 nb_sst(i,j,bi,bj) = 0. _d 0
166 msk_sst(i,j,bi,bj) = 0. _d 0
167 enddo
168 enddo
169 enddo
170 enddo
171
172 c PART 1.1: compute running sample average over ndaysave
173 c ------------------------------------------------------
174
175 do jrec=1,ndaysave
176
177 krec=irec+jrec-1
178
179 c get modeled sst:
180 call active_read_xy( fname, sstbar, krec, doglobalread,
181 & ladinit, optimcycle, mythid,
182 & xx_sstbar_mean_dummy )
183
184 c get observed sst:
185 daytime = FLOAT(secondsperday*(krec-1)) + modelstart
186 dayiter = hoursperday*(krec-1)+modeliter0
187 call cal_getdate( dayiter, daytime, daydate, mythid )
188 call cal_convdate( daydate,yday,md,dd,sd,ld,wd,mythid )
189 ymod = locstartdate(1)/10000
190 if ( ymod .GE. yday ) then
191 middate(1)=1
192 call cal_FullDate(locstartdate(1),0,middate,mythid)
193 else
194 middate(1)=1
195 tempDate_1 = yday*10000+100+1
196 call cal_FullDate( tempDate_1, 0, middate, mythid)
197 endif
198 call cal_TimePassed( middate, daydate, difftime, mythid )
199 call cal_ToSeconds( difftime, diffsecs, mythid )
200 c localrec = floor(diffsecs/86400.) + 1
201 localrec = int(diffsecs/86400.) + 1
202
203 il=ilnblnk(gencost_datafile(igen_amsre))
204 write(fname2(1:80),'(2a,i4)')
205 & gencost_datafile(igen_amsre)(1:il), '_', yday
206 inquire( file=fname2, exist=exst )
207
208 #ifdef ALLOW_ECCO_DEBUG
209 WRITE(msgBuf,'(A,I4,A,I4,A,I10,A,1PE15.2)') 'sstv4 reading ',
210 & yday,' ',ymod,' ',localrec,' ',diffsecs
211 CALL PRINT_MESSAGE( msgBuf, ioUnit, SQUEEZE_RIGHT, myThid )
212 c
213 CALL CAL_PRINTDATE(middate,mythid)
214 CALL CAL_PRINTDATE(daydate,mythid)
215 CALL CAL_PRINTDATE(difftime,mythid)
216 #endif
217
218 if ( ( localrec .GT. 0 ).AND.(diffsecs .GT. 0.d0) ) then
219 call mdsreadfield( fname2, cost_iprec, cost_yftype, 1,
220 & tmp_sst, localrec, mythid )
221 else
222 do bj = jtlo,jthi
223 do bi = itlo,ithi
224 do j = jmin,jmax
225 do i = imin,imax
226 tmp_sst(i,j,bi,bj) = spval
227 enddo
228 enddo
229 enddo
230 enddo
231 endif
232
233 c accumulate obs and misfit:
234 do bj = jtlo,jthi
235 do bi = itlo,ithi
236 do j = jmin,jmax
237 do i = imin,imax
238 if ( (tmp_sst(i,j,bi,bj).GT.spval).AND.
239 & (maskc(i,j,1,bi,bj).EQ.1.) ) then
240 anom_sst(i,j,bi,bj)= anom_sst(i,j,bi,bj)+
241 & sstbar(i,j,bi,bj)-tmp_sst(i,j,bi,bj)
242 obs_sst(i,j,bi,bj)= obs_sst(i,j,bi,bj)+
243 & tmp_sst(i,j,bi,bj)
244 nb_sst(i,j,bi,bj)=nb_sst(i,j,bi,bj)+1. _d 0
245 endif
246 enddo
247 enddo
248 enddo
249 enddo
250
251 enddo !do jrec=1,ndaysave
252
253 c average obs and misfit:
254 do bj = jtlo,jthi
255 do bi = itlo,ithi
256 do j = jmin,jmax
257 do i = imin,imax
258 if ( nb_sst(i,j,bi,bj) .NE. 0. ) then
259 obs_sst(i,j,bi,bj) =
260 & obs_sst(i,j,bi,bj)/nb_sst(i,j,bi,bj)
261 anom_sst(i,j,bi,bj) =
262 & anom_sst(i,j,bi,bj)/nb_sst(i,j,bi,bj)
263 msk_sst(i,j,bi,bj) = 1. _d 0
264 endif
265 enddo
266 enddo
267 enddo
268 enddo
269
270
271 c PART 1.2: smooth anom_sst in space
272 c ----------------------------------------
273
274 #ifdef ALLOW_GENCOST_SSTV4_OUTPUT
275 write(fname2(1:80),'(1a)') 'sstdiff_raw'
276 call mdswritefield(fname2,32,.false.,'RL',
277 & 1,anom_sst,irec,1,mythid)
278
279 write(fname2(1:80),'(1a)') 'sstobs_raw'
280 call mdswritefield(fname2,32,.false.,'RL',
281 & 1,obs_sst,irec,1,mythid)
282 #endif
283
284 call smooth_hetero2d(anom_sst,maskc,
285 & gencost_scalefile(igen_amsre_lsc),
286 & gencost_smooth2Ddiffnbt(igen_amsre_lsc),mythid)
287
288 #ifdef ALLOW_GENCOST_SSTV4_OUTPUT
289 call smooth_hetero2d(obs_sst,maskc,
290 & gencost_scalefile(igen_amsre_lsc),
291 & gencost_smooth2Ddiffnbt(igen_amsre_lsc),mythid)
292
293 write(fname2(1:80),'(1a)') 'sstdiff_smooth'
294 call mdswritefield(fname2,32,.false.,'RL',
295 & 1,anom_sst,irec,1,mythid)
296
297 write(fname2(1:80),'(1a)') 'sstobs_smooth'
298 call mdswritefield(fname2,32,.false.,'RL',
299 & 1,obs_sst,irec,1,mythid)
300 #endif
301
302 c PART 1.3: compute cost function term
303 c ------------------------------------
304
305 do bj = jtlo,jthi
306 do bi = itlo,ithi
307 do j = jmin,jmax
308 do i = imin,imax
309 junk = anom_sst(i,j,bi,bj)
310 junkweight = gencost_weight(i,j,bi,bj,igen_amsre_lsc)*
311 & maskc(i,j,1,bi,bj)
312 objf_gencost(igen_amsre_lsc,bi,bj) =
313 & objf_gencost(igen_amsre_lsc,bi,bj)
314 & +junk*junk*junkweight/ndaysaveRL
315 if ( (junkweight.GT.0.).AND.(nb_sst(i,j,bi,bj).GT.0.) )
316 & num_gencost(igen_amsre_lsc,bi,bj) =
317 & num_gencost(igen_amsre_lsc,bi,bj) + 1. _d 0 /ndaysaveRL
318 enddo
319 enddo
320 enddo
321 enddo
322
323 enddo
324
325 cgf =======================================================
326 cgf PART 2: compute raw SST cost term
327 cgf =======================================================
328
329 do irec = 1, ndaysrec
330
331 c get modeled sst:
332 call active_read_xy( fname, sstbar, irec, doglobalread,
333 & ladinit, optimcycle, mythid,
334 & xx_sstbar_mean_dummy )
335
336 c get observed sst:
337 daytime = FLOAT(secondsperday*(irec-1)) + modelstart
338 dayiter = hoursperday*(irec-1)+modeliter0
339 call cal_getdate( dayiter, daytime, daydate, mythid )
340 call cal_convdate( daydate,yday,md,dd,sd,ld,wd,mythid )
341 ymod = locstartdate(1)/10000
342 if ( ymod .GE. yday ) then
343 middate(1)=1
344 call cal_FullDate(locstartdate(1),0,middate,mythid)
345 else
346 middate(1)=1
347 tempDate_1 = yday*10000+100+1
348 call cal_FullDate( tempDate_1, 0, middate, mythid)
349 endif
350 call cal_TimePassed( middate, daydate, difftime, mythid )
351 call cal_ToSeconds( difftime, diffsecs, mythid )
352 c localrec = floor(diffsecs/86400.) + 1
353 localrec = int(diffsecs/86400.) + 1
354
355 il=ilnblnk(gencost_datafile(igen_amsre))
356 write(fname2(1:80),'(2a,i4)')
357 & gencost_datafile(igen_amsre)(1:il), '_', yday
358 inquire( file=fname2, exist=exst )
359
360 #ifdef ALLOW_ECCO_DEBUG
361 WRITE(msgBuf,'(A,I4,A,I4,A,I10,A,1PE15.2)') 'sstv4 reading ',
362 & yday,' ',ymod,' ',localrec,' ',diffsecs
363 CALL PRINT_MESSAGE( msgBuf, ioUnit, SQUEEZE_RIGHT, myThid )
364 c
365 CALL CAL_PRINTDATE(middate,mythid)
366 CALL CAL_PRINTDATE(daydate,mythid)
367 CALL CAL_PRINTDATE(difftime,mythid)
368 #endif
369
370 if ( ( localrec .GT. 0 ).AND.(diffsecs .GT. 0.d0) ) then
371 call mdsreadfield( fname2, cost_iprec, cost_yftype, 1,
372 & tmp_sst, localrec, mythid )
373 else
374 do bj = jtlo,jthi
375 do bi = itlo,ithi
376 do j = jmin,jmax
377 do i = imin,imax
378 tmp_sst(i,j,bi,bj) = spval
379 enddo
380 enddo
381 enddo
382 enddo
383 endif
384
385 c compute misfit:
386 do bj = jtlo,jthi
387 do bi = itlo,ithi
388 do j = jmin,jmax
389 do i = imin,imax
390 if ( (tmp_sst(i,j,bi,bj).GT.spval).AND.
391 & (maskc(i,j,1,bi,bj).EQ.1.) ) then
392 anom_sst(i,j,bi,bj) =
393 & sstbar(i,j,bi,bj)-tmp_sst(i,j,bi,bj)
394 msk_sst(i,j,bi,bj) = 1. _d 0
395 else
396 anom_sst(i,j,bi,bj) = 0. _d 0
397 msk_sst(i,j,bi,bj) = 0. _d 0
398 endif
399 enddo
400 enddo
401 enddo
402 enddo
403
404 #ifdef ALLOW_GENCOST_SSTV4_OUTPUT
405 write(fname2(1:80),'(1a)') 'sstdiff_point'
406 call mdswritefield(fname2,32,.false.,'RL',
407 & 1,anom_sst,irec,1,mythid)
408 #endif
409
410 c compute cost:
411
412 do bj = jtlo,jthi
413 do bi = itlo,ithi
414 do j = jmin,jmax
415 do i = imin,imax
416 junk = anom_sst(i,j,bi,bj)
417 junkweight = gencost_weight(i,j,bi,bj,igen_amsre)*
418 & maskc(i,j,1,bi,bj)*msk_sst(i,j,bi,bj)
419 objf_gencost(igen_amsre,bi,bj) =
420 & objf_gencost(igen_amsre,bi,bj)+junk*junk*junkweight
421 if (junkweight.GT.0.)
422 & num_gencost(igen_amsre,bi,bj) =
423 & num_gencost(igen_amsre,bi,bj) + 1. _d 0
424 enddo
425 enddo
426 enddo
427 enddo
428
429 enddo
430
431 #endif /* ifdef ALLOW_GENCOST_CONTRIBUTION */
432 #endif /* ifdef ALLOW_SMOOTH */
433 #endif /* ifdef ALLOW_DAILYSST_COST_CONTRIBUTION */
434
435 end

  ViewVC Help
Powered by ViewVC 1.1.22