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

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

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


Revision 1.1 - (show annotations) (download)
Thu Nov 6 22:10:07 2003 UTC (20 years, 7 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint52l_pre, checkpoint52e_pre, hrcube4, checkpoint52n_post, checkpoint52j_post, checkpoint53d_post, checkpoint54a_pre, checkpoint55c_post, checkpoint54e_post, checkpoint52e_post, checkpoint54a_post, checkpoint53c_post, checkpoint55d_pre, hrcube_1, branch-netcdf, checkpoint52d_pre, checkpoint52l_post, checkpoint52k_post, checkpoint52b_pre, checkpoint54b_post, checkpoint53b_pre, checkpoint55b_post, checkpoint54d_post, checkpoint52m_post, checkpoint55, checkpoint53a_post, checkpoint54, checkpoint54f_post, checkpoint53b_post, checkpoint52a_pre, checkpoint53, checkpoint52, checkpoint52d_post, checkpoint52a_post, checkpoint52b_post, checkpoint53g_post, checkpoint52f_post, checkpoint52c_post, ecco_c52_e35, hrcube5, checkpoint52i_post, checkpoint52j_pre, checkpoint53f_post, checkpoint55a_post, checkpoint53d_pre, checkpoint54c_post, checkpoint52i_pre, checkpoint51u_post, checkpoint52h_pre, checkpoint52f_pre, hrcube_2, hrcube_3
Branch point for: netcdf-sm0
o merging from ecco-branch
o pkg/ecco now containes ecco-specific part of cost function
o top level routines the_main_loop, forward_step
  supersede those in model/src/
  previous input data.cost now in data.ecco
  (new namelist ecco_cost_nml)

1
2 #include "COST_CPPOPTIONS.h"
3 #ifdef ALLOW_OBCS
4 # include "OBCS_OPTIONS.h"
5 #endif
6
7 subroutine cost_obcsn(
8 I myiter,
9 I mytime,
10 I startrec,
11 I endrec,
12 I mythid
13 & )
14
15 c ==================================================================
16 c SUBROUTINE cost_obcsn
17 c ==================================================================
18 c
19 c o cost function contribution obc
20 c
21 c o G. Gebbie, gebbie@mit.edu, 18-Mar-2003
22 c ==================================================================
23 c SUBROUTINE cost_obcsn
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 #include "DYNVARS.h"
35 #ifdef ALLOW_OBCS
36 # include "OBCS.h"
37 #endif
38
39 #include "cal.h"
40 #include "ecco_cost.h"
41 #include "ctrl.h"
42 #include "ctrl_dummy.h"
43 #include "optim.h"
44
45 c == routine arguments ==
46
47 integer myiter
48 _RL mytime
49 integer mythid
50 cgg(
51 integer startrec
52 integer endrec
53 cgg)
54
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
64 integer il
65 integer iobcs
66 integer jp1
67 integer nrec
68 integer ilfld
69 integer igg
70 integer ihh
71
72 _RL fctile
73 _RL fcthread
74 _RL dummy
75 _RL gg
76 _RL tmpx
77 _RL tmpfield (1-olx:snx+olx,nr,nsx,nsy)
78 _RL barofield (1-olx:snx+olx,nsx,nsy)
79 _RL maskxz (1-olx:snx+olx,nr,nsx,nsy)
80 _RL vtop
81 _RL area
82 _RL volflux
83
84 character*(80) fnamefld
85
86 logical doglobalread
87 logical ladinit
88
89 #ifdef ECCO_VERBOSE
90 character*(MAX_LEN_MBUF) msgbuf
91 #endif
92
93 c == external functions ==
94
95 integer ilnblnk
96 external ilnblnk
97
98 c == end of interface ==
99
100 jtlo = mybylo(mythid)
101 jthi = mybyhi(mythid)
102 itlo = mybxlo(mythid)
103 ithi = mybxhi(mythid)
104 jmin = 1
105 jmax = sny
106 imin = 1
107 imax = snx
108
109 c-- Read tiled data.
110 doglobalread = .false.
111 ladinit = .false.
112
113 c Number of records to be used.
114 nrec = endrec-startrec+1
115
116 #ifdef ALLOW_OBCSN_COST_CONTRIBUTION
117
118 jp1 = 0
119 fcthread = 0. _d 0
120
121 #ifdef ECCO_VERBOSE
122 _BEGIN_MASTER( mythid )
123 write(msgbuf,'(a)') ' '
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 write(msgbuf,'(a,i9.8)')
130 & ' cost_obcsn: number of records to process: ',nrec
131 call print_message( msgbuf, standardmessageunit,
132 & SQUEEZE_RIGHT , mythid)
133 write(msgbuf,'(a)') ' '
134 call print_message( msgbuf, standardmessageunit,
135 & SQUEEZE_RIGHT , mythid)
136 _END_MASTER( mythid )
137 #endif
138
139 if (optimcycle .ge. 0) then
140 ilfld=ilnblnk( xx_obcsn_file )
141 write(fnamefld(1:80),'(2a,i10.10)')
142 & xx_obcsn_file(1:ilfld), '.', optimcycle
143 endif
144
145 c-- Loop over records.
146 do irec = 1,nrec
147
148 area = 0. _d 0
149 volflux = 0. _d 0
150 do bj = jtlo,jthi
151 do bi = itlo,ithi
152 do i = imin,imax
153 barofield(i,bi,bj) = 0. _d 0
154 enddo
155 enddo
156 enddo
157
158 call active_read_xz( fnamefld, tmpfield, irec, doglobalread,
159 & ladinit, optimcycle, mythid
160 & , xx_obcsn_dummy )
161
162 cgg Need to solve for iobcs would have been.
163 gg = (irec-1)/nobcs
164 igg = int(gg)
165 iobcs = irec - igg*nobcs
166
167 call active_read_xz( 'maskobcsn', maskxz,
168 & iobcs,
169 & doglobalread, ladinit, 0,
170 & mythid, dummy )
171
172 cgg Need to change xx field to baroclinic vel.
173 cgg Level Nr contains barotropic vel, remove it.
174 cgg The deepest wet level velocity should be
175 cgg computed in order to ensure zero integrated flux.
176 if (iobcs .eq. 3 .or. iobcs .eq. 4) then
177 if (iobcs .eq .3) then
178 ihh = igg+1
179 call ctrl_volflux( ihh, area, volflux, mythid)
180 _GLOBAL_SUM_R8( volflux, mythid )
181 _GLOBAL_SUM_R8( area,mythid )
182 print*,'volflux,area',volflux,area
183 endif
184
185 c-- Loop over this thread's tiles.
186 do bj = jtlo,jthi
187 do bi = itlo,ithi
188 do i = imin,imax
189
190 cgg The barotropic velocity is stored in level 1.
191 cgg Penalize for deviation in it.
192 barofield(i,bi,bj) = tmpfield(i,1,bi,bj)
193 & * maskxz(i,1,bi,bj)
194
195 vtop = 0.d0
196 tmpfield(i,1,bi,bj) = 0.d0
197
198 j = OB_Jn(I,bi,bj)
199 if (iobcs .eq. 3) then
200 if (ob_iw(j,bi,bj).eq.(i-1).and.
201 & ob_iw(j,bi,bj).ne. 0) then
202 barofield(i,bi,bj)= volflux/area
203 & * maskxz(i,1,bi,bj)
204 print*,'volflux2,area2',volflux,area
205 print*,'barofield',barofield(i,bi,bj)
206 endif
207 endif
208
209 do k = 1,Nr
210 cgg If cells are not full, this should be modified with hFac.
211 cgg
212 cgg The xx field (tmpfldxz) does not contain the velocity at the
213 cgg lowest wet level. This velocity is not independent; it must
214 cgg exactly balance the volume flux, since we are dealing with
215 cgg the baroclinic velocity structure..
216 vtop = tmpfield(i,k,bi,bj)*
217 & maskxz(i,k,bi,bj) * delZ(k) + vtop
218 enddo
219 tmpfield(i,1,bi,bj) = -vtop / delZ(1)
220 enddo
221 enddo
222 enddo
223
224 endif
225
226 c-- Loop over this thread's tiles.
227 do bj = jtlo,jthi
228 do bi = itlo,ithi
229
230 c-- Determine the weights to be used.
231 fctile = 0. _d 0
232
233 if (iobcs .eq. 3 .or. iobcs .eq. 4) then
234 do i = imin,imax
235 j = OB_Jn(I,bi,bj)
236 tmpx = barofield(i,bi,bj)
237 fctile = fctile
238 & + wbaro*cosphi(i,j,bi,bj)
239 & *tmpx*tmpx
240 & *maskxz(i,1,bi,bj)
241 cgg print*,'S bt fctile',i,j,wbaro,
242 cgg & barofield(i,bi,bj),fctile
243 enddo
244 endif
245
246 do k = 1, Nr
247 do i = imin,imax
248 j = OB_Jn(I,bi,bj)
249 cgg if (maskS(i,j+jp1,k,bi,bj) .ne. 0.) then
250 tmpx = tmpfield(i,k,bi,bj)
251 fctile = fctile
252 & + wobcsnLev(i,k,bi,bj,iobcs)*cosphi(i,j,bi,bj)
253 & *tmpx*tmpx
254 & *maskxz(i,k,bi,bj)
255 cgg endif
256 cgg print*,'S fctile',fctile
257 enddo
258 enddo
259
260 objf_obcsn(bi,bj) = objf_obcsn(bi,bj) + fctile
261 fcthread = fcthread + fctile
262 enddo
263 enddo
264
265 #ifdef ECCO_VERBOSE
266 c-- Print cost function for all tiles.
267 _GLOBAL_SUM_R8( fcthread , myThid )
268 write(msgbuf,'(a)') ' '
269 call print_message( msgbuf, standardmessageunit,
270 & SQUEEZE_RIGHT , mythid)
271 write(msgbuf,'(a,i8.8)')
272 & ' cost_obcsn: irec = ',irec
273 call print_message( msgbuf, standardmessageunit,
274 & SQUEEZE_RIGHT , mythid)
275 write(msgbuf,'(a,a,d22.15)')
276 & ' global cost function value',
277 & ' (obcsn) = ',fcthread
278 call print_message( msgbuf, standardmessageunit,
279 & SQUEEZE_RIGHT , mythid)
280 write(msgbuf,'(a)') ' '
281 call print_message( msgbuf, standardmessageunit,
282 & SQUEEZE_RIGHT , mythid)
283 #endif
284
285 enddo
286 c-- End of loop over records.
287
288 #endif
289
290 return
291 end
292
293
294
295
296
297
298

  ViewVC Help
Powered by ViewVC 1.1.22