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

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

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


Revision 1.3 - (show annotations) (download)
Fri Mar 29 14:53:09 2013 UTC (11 years, 2 months ago) by gforget
Branch: MAIN
CVS Tags: checkpoint64o, checkpoint64q, checkpoint64p, checkpoint64s, checkpoint64r, checkpoint64n, checkpoint64g, checkpoint64f, checkpoint64i, checkpoint64h, checkpoint64k, checkpoint64j, checkpoint64m, checkpoint64l
Changes since 1.2: +5 -4 lines
- bug fix.

1 C $Header: /u/gcmpack/MITgcm/pkg/ecco/cost_gencost_seaicev4.F,v 1.2 2013/03/28 13:40:09 jmc Exp $
2 C $Name: $
3
4 #include "ECCO_OPTIONS.h"
5
6 subroutine cost_gencost_seaicev4(mythid)
7
8 c ==================================================================
9 c SUBROUTINE cost_gencost_seaicev4
10 c ==================================================================
11 c
12 c o Evaluate cost function contributions of ice concentration.
13 c
14 c ==================================================================
15 c SUBROUTINE cost_gencost_seaicev4
16 c ==================================================================
17
18 implicit none
19
20 c == global variables ==
21
22 #include "EEPARAMS.h"
23 #include "SIZE.h"
24 #include "PARAMS.h"
25 #include "GRID.h"
26 #ifdef ALLOW_CAL
27 # include "cal.h"
28 #endif
29 #ifdef ALLOW_COST
30 # include "ecco_cost.h"
31 # include "optim.h"
32 # ifdef ALLOW_SEAICE
33 # include "SEAICE_COST.h"
34 # include "SEAICE_PARAMS.h"
35 # endif
36 #endif
37
38 c == routine arguments ==
39 integer mythid
40
41 #ifdef ALLOW_SEAICE_COST_CONTRIBUTION
42 #ifdef ALLOW_GENCOST_CONTRIBUTION
43
44 c == local variables ==
45
46 integer nnzobs
47 parameter (nnzobs = 1 )
48 integer nrecloc
49 integer localstartdate(4)
50
51 _RL areabbbar (1-olx:snx+olx,1-oly:sny+oly,1,nsx,nsy)
52 _RL heffbbbar (1-olx:snx+olx,1-oly:sny+oly,1,nsx,nsy)
53 _RL sstbbbar (1-olx:snx+olx,1-oly:sny+oly,1,nsx,nsy)
54 _RL localweight(1-olx:snx+olx,1-oly:sny+oly,1,nsx,nsy)
55 _RL xx_areabbbar_mean_dummy
56 _RL xx_heffbbbar_mean_dummy
57 _RL xx_sstbbbar_mean_dummy
58 _RL mult_local
59 _RL localperiod
60 _RL spminloc
61 _RL spmaxloc
62 _RL spzeroloc
63 _RL objf_local(nsx,nsy)
64 _RL num_local(nsx,nsy)
65
66 character*(MAX_LEN_FNAM) areabbbarfile
67 character*(MAX_LEN_FNAM) heffbbbarfile
68 character*(MAX_LEN_FNAM) sstbbbarfile
69 character*(MAX_LEN_FNAM) localobsfile
70
71 integer igen_conc, igen_sst, igen_vol
72
73 integer bi,bj
74 integer i,j,k
75 integer itlo,ithi
76 integer jtlo,jthi
77 integer jmin,jmax
78 integer imin,imax
79 integer irec
80 integer il
81 integer localrec
82 integer obsrec
83
84 logical doglobalread
85 logical ladinit
86
87 _RL spval
88 parameter (spval = -9999. )
89 _RL localwww
90 _RL junk
91
92 _RL localmask (1-olx:snx+olx,1-oly:sny+oly,nr,nsx,nsy)
93 _RL localobs (1-olx:snx+olx,1-oly:sny+oly,nnzobs,nsx,nsy)
94 _RL cmask (1-olx:snx+olx,1-oly:sny+oly,nnzobs)
95
96 character*(128) fname0, fname1, fname2, fname3
97 character*(MAX_LEN_MBUF) msgbuf
98
99 _RL daytime
100 _RL diffsecs
101 integer dayiter
102 integer daydate(4)
103 integer difftime(4)
104 integer middate(4)
105 integer yday, ymod
106 integer md, dd, sd, ld, wd
107 integer mody, modm
108 integer beginmodel, beginlocal
109 logical exst
110
111 c == external functions ==
112
113 integer ilnblnk
114 external ilnblnk
115
116 c == end of interface ==
117
118 jtlo = mybylo(mythid)
119 jthi = mybyhi(mythid)
120 itlo = mybxlo(mythid)
121 ithi = mybxhi(mythid)
122 jmin = 1
123 jmax = sny
124 imin = 1
125 imax = snx
126
127 c-- detect the relevant gencost indices
128 igen_conc=0
129 igen_sst=0
130 igen_vol=0
131 do k=1,NGENCOST
132 if (gencost_name(k).EQ.'siv4-conc') igen_conc=k
133 if (gencost_name(k).EQ.'siv4-sst') igen_sst=k
134 if (gencost_name(k).EQ.'siv4-vol') igen_vol=k
135 enddo
136
137 if ((igen_conc.NE.0).AND.(igen_sst.NE.0).AND.(igen_vol.NE.0))
138 & then
139
140 c-- Initialise local variables.
141
142 localwww = 0. _d 0
143
144 do bj = jtlo,jthi
145 do bi = itlo,ithi
146 objf_gencost(bi,bj,igen_conc) = 0. _d 0
147 objf_gencost(bi,bj,igen_vol) = 0. _d 0
148 objf_gencost(bi,bj,igen_sst) = 0. _d 0
149 num_gencost(bi,bj,igen_conc) = 0. _d 0
150 num_gencost(bi,bj,igen_vol) = 0. _d 0
151 num_gencost(bi,bj,igen_sst) = 0. _d 0
152 do k = 1,nnzobs
153 do j = jmin,jmax
154 do i = imin,imax
155 localobs(i,j,k,bi,bj) = 0. _d 0
156 enddo
157 enddo
158 enddo
159 enddo
160 enddo
161
162 c-- Assign mask
163 do bj = jtlo,jthi
164 do bi = itlo,ithi
165 do k = 1,Nr
166 do j = 1-oly,sny+oly
167 do i = 1-olx,snx+olx
168 localmask(i,j,k,bi,bj) = maskC(i,j,k,bi,bj)
169 enddo
170 enddo
171 enddo
172 enddo
173 enddo
174
175 nrecloc=nmonsrec
176 localperiod=0.
177
178 areabbbarfile=gencost_barfile(igen_conc)
179 heffbbbarfile=gencost_barfile(igen_vol)
180 sstbbbarfile=gencost_barfile(igen_sst)
181 localobsfile=gencost_datafile(igen_conc)
182 xx_areabbbar_mean_dummy=xx_genbar_dummy(igen_conc)
183 xx_heffbbbar_mean_dummy=xx_genbar_dummy(igen_vol)
184 xx_sstbbbar_mean_dummy=xx_genbar_dummy(igen_sst)
185 localstartdate(1)=modelstartdate(1)
186 localstartdate(2)=modelstartdate(2)
187 localstartdate(3)=modelstartdate(3)
188 localstartdate(4)=modelstartdate(4)
189 spminloc=gencost_spmin(igen_conc)
190 spmaxloc=gencost_spmax(igen_conc)
191 spzeroloc=gencost_spzero(igen_conc)
192
193 c-- First, read tiled data.
194 doglobalread = .false.
195 ladinit = .false.
196
197 write(fname1(1:128),'(80a)') ' '
198 il=ilnblnk( areabbbarfile )
199 write(fname1(1:128),'(2a,i10.10)')
200 & areabbbarfile(1:il),'.',optimcycle
201
202 il=ilnblnk( heffbbbarfile )
203 write(fname2(1:128),'(2a,i10.10)')
204 & heffbbbarfile(1:il),'.',optimcycle
205
206 il=ilnblnk( sstbbbarfile )
207 write(fname3(1:128),'(2a,i10.10)')
208 & sstbbbarfile(1:il),'.',optimcycle
209
210 if ( .NOT. ( localobsfile.EQ.' ' ) ) then
211
212 c-- Loop over records for the second time.
213 do irec = 1, nrecloc
214
215 call active_read_xy( fname1, areabbbar, irec, doglobalread,
216 & ladinit, optimcycle, mythid,
217 & xx_areabbbar_mean_dummy )
218
219 call active_read_xy( fname2, heffbbbar, irec, doglobalread,
220 & ladinit, optimcycle, mythid,
221 & xx_heffbbbar_mean_dummy )
222
223 call active_read_xy( fname3, sstbbbar, irec, doglobalread,
224 & ladinit, optimcycle, mythid,
225 & xx_sstbbbar_mean_dummy )
226
227 if ( localperiod .EQ. 86400. ) then
228 c-- assume daily fields
229 obsrec = irec
230 daytime = FLOAT(secondsperday*(irec-1))
231 dayiter = hoursperday*(irec-1)
232 call cal_getdate( dayiter, daytime, daydate, mythid )
233 call cal_convdate( daydate,yday,md,dd,sd,ld,wd,mythid )
234 ymod = localstartdate(1)/10000
235 if ( ymod .EQ. yday ) then
236 middate(1) = modelstartdate(1)
237 else
238 middate(1) = yday*10000+100+1
239 endif
240 middate(2) = 0
241 middate(3) = modelstartdate(3)
242 middate(4) = modelstartdate(4)
243 call cal_TimePassed( middate, daydate, difftime, mythid )
244 call cal_ToSeconds( difftime, diffsecs, mythid )
245 localrec = int(diffsecs/localperiod) + 1
246 else
247 c-- assume monthly fields
248 beginlocal = localstartdate(1)/10000
249 beginmodel = modelstartdate(1)/10000
250 obsrec =
251 & ( beginmodel - beginlocal )*nmonthyear
252 & + ( mod(modelstartdate(1)/100,100)
253 & -mod(localstartdate(1)/100,100) )
254 & + irec
255 mody = modelstartdate(1)/10000
256 modm = modelstartdate(1)/100 - mody*100
257 yday = mody + INT((modm-1+irec-1)/12)
258 localrec = 1 + MOD(modm-1+irec-1,12)
259 endif
260
261 il=ilnblnk(localobsfile)
262 write(fname0(1:128),'(2a,i4)')
263 & localobsfile(1:il), '_', yday
264 inquire( file=fname0, exist=exst )
265 if (.NOT. exst) then
266 write(fname0(1:128),'(a)') localobsfile(1:il)
267 c to use the data in a repreated cycle, comment next line?
268 localrec = obsrec
269 endif
270
271 if ( localrec .GT. 0 ) then
272 call mdsreadfield( fname0, cost_iprec, cost_yftype, nnzobs,
273 & localobs, localrec, mythid )
274 else
275 do bj = jtlo,jthi
276 do bi = itlo,ithi
277 do k = 1,nnzobs
278 do j = jmin,jmax
279 do i = imin,imax
280 localobs(i,j,k,bi,bj) = spval
281 c not sure why this is not spzeroloc
282 enddo
283 enddo
284 enddo
285 enddo
286 enddo
287 endif
288
289 do bj = jtlo,jthi
290 do bi = itlo,ithi
291
292 c-- Determine the mask on weights
293 do k = 1,nnzobs
294 do j = jmin,jmax
295 do i = imin,imax
296 cmask(i,j,k) = cosphi(i,j,bi,bj)*localmask(i,j,k,bi,bj)
297 if ( localobs(i,j,k,bi,bj) .lt. spminloc .or.
298 & localobs(i,j,k,bi,bj) .gt. spmaxloc .or.
299 & localobs(i,j,k,bi,bj) .eq. spzeroloc ) then
300 cmask(i,j,k) = 0. _d 0
301 endif
302 localweight(i,j,k,bi,bj)=
303 & gencost_weight(i,j,bi,bj,igen_conc)
304 enddo
305 enddo
306 enddo
307 c--
308 do k = 1,nnzobs
309 do j = jmin,jmax
310 do i = imin,imax
311
312 c area term
313 localwww = localweight(i,j,k,bi,bj)*cmask(i,j,k)
314 junk = ( areabbbar(i,j,k,bi,bj) -
315 & localobs(i,j,k,bi,bj) )
316 objf_gencost(bi,bj,igen_conc) =
317 & objf_gencost(bi,bj,igen_conc) + junk*junk*localwww
318
319 if ( localwww .ne. 0. )
320 & num_gencost(bi,bj,igen_conc) =
321 & num_gencost(bi,bj,igen_conc) + 1. _d 0
322
323 c heff term
324 if ( (localobs(i,j,k,bi,bj) .EQ. 0.).AND.
325 & (heffbbbar(i,j,k,bi,bj) .GT. 0.) ) then
326 junk=10. _d 0 *cmask(i,j,k)*heffbbbar(i,j,k,bi,bj)
327 num_gencost(bi,bj,igen_vol) =
328 & num_gencost(bi,bj,igen_vol) + 1. _d 0
329 else
330 junk = 0. _d 0
331 num_gencost(bi,bj,igen_vol) =
332 & num_gencost(bi,bj,igen_vol) + 0. _d 0
333 endif
334
335 objf_gencost(bi,bj,igen_vol) =
336 & objf_gencost(bi,bj,igen_vol) + junk
337
338 c sst term
339 if ( (areabbbar(i,j,1,bi,bj) .LE. 0.).AND.
340 & (localobs(i,j,1,bi,bj) .GT. 0.) ) then
341 junk=1. _d 0 *cmask(i,j,k)*
342 & ( 3. _d 0 + sstbbbar(i,j,k,bi,bj) )
343 num_gencost(bi,bj,igen_sst) =
344 & num_gencost(bi,bj,igen_sst) + 1. _d 0
345 else
346 junk = 0. _d 0
347 num_gencost(bi,bj,igen_sst) =
348 & num_gencost(bi,bj,igen_sst) + 0. _d 0
349 endif
350
351 objf_gencost(bi,bj,igen_sst) =
352 & objf_gencost(bi,bj,igen_sst) + junk
353
354 enddo
355 enddo
356 enddo
357
358 enddo
359 enddo
360
361 enddo
362
363 endif !if ( .NOT. ( localobsfile.EQ.' ' ) ) then
364 endif !if ((igen_conc.NE.0).AND.(igen_sst.NE.0).AND.(igen_vol.NE.0)) then
365
366 #endif /* ALLOW_GENCOST_CONTRIBUTION */
367 #endif /* ALLOW_SEAICE_COST_CONTRIBUTION */
368
369 end

  ViewVC Help
Powered by ViewVC 1.1.22