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

Contents of /MITgcm/pkg/ecco/cost_obcse.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, 6 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_obcse(
8 I myiter,
9 I mytime,
10 I startrec,
11 I endrec,
12 I mythid
13 & )
14
15 c ==================================================================
16 c SUBROUTINE cost_obcse
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_obcse
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 ip1
67 integer nrec
68 integer ilfld
69 integer igg
70
71 _RL fctile
72 _RL fcthread
73 _RL dummy
74 _RL gg
75 _RL tmpx
76 _RL tmpfield (1-oly:sny+oly,nr,nsx,nsy)
77 _RL barofield (1-oly:sny+oly,nsx,nsy)
78 _RL maskyz (1-oly:sny+oly,nr,nsx,nsy)
79 _RL ubottom
80
81 character*(80) fnamefld
82
83 logical doglobalread
84 logical ladinit
85
86 #ifdef ECCO_VERBOSE
87 character*(MAX_LEN_MBUF) msgbuf
88 #endif
89
90 c == external functions ==
91
92 integer ilnblnk
93 external ilnblnk
94
95 c == end of interface ==
96
97 jtlo = mybylo(mythid)
98 jthi = mybyhi(mythid)
99 itlo = mybxlo(mythid)
100 ithi = mybxhi(mythid)
101 jmin = 1
102 jmax = sny
103 imin = 1
104 imax = snx
105
106 c-- Read tiled data.
107 doglobalread = .false.
108 ladinit = .false.
109
110 c Number of records to be used.
111 nrec = endrec-startrec+1
112
113 #ifdef ALLOW_OBCSE_COST_CONTRIBUTION
114
115 ip1 = 0
116 fcthread = 0. _d 0
117
118 #ifdef ECCO_VERBOSE
119 _BEGIN_MASTER( mythid )
120 write(msgbuf,'(a)') ' '
121 call print_message( msgbuf, standardmessageunit,
122 & SQUEEZE_RIGHT , mythid)
123 write(msgbuf,'(a)') ' '
124 call print_message( msgbuf, standardmessageunit,
125 & SQUEEZE_RIGHT , mythid)
126 write(msgbuf,'(a,i9.8)')
127 & ' cost_obcse: number of records to process: ',nrec
128 call print_message( msgbuf, standardmessageunit,
129 & SQUEEZE_RIGHT , mythid)
130 write(msgbuf,'(a)') ' '
131 call print_message( msgbuf, standardmessageunit,
132 & SQUEEZE_RIGHT , mythid)
133 _END_MASTER( mythid )
134 #endif
135
136 if (optimcycle .ge. 0) then
137 ilfld=ilnblnk( xx_obcse_file )
138 write(fnamefld(1:80),'(2a,i10.10)')
139 & xx_obcse_file(1:ilfld), '.', optimcycle
140 endif
141
142 c-- Loop over records.
143 do irec = 1,nrec
144
145 do bj = jtlo,jthi
146 do bi = itlo,ithi
147 do j = jmin,jmax
148 barofield(j,bi,bj) = 0. _d 0
149 enddo
150 enddo
151 enddo
152
153 call active_read_yz( fnamefld, tmpfield, irec, doglobalread,
154 & ladinit, optimcycle, mythid
155 & , xx_obcse_dummy )
156
157 cgg Need to solve for iobcs would have been.
158 gg = (irec-1)/nobcs
159 igg = int(gg)
160 iobcs = irec - igg*nobcs
161
162 c-- Loop over this thread's tiles.
163 do bj = jtlo,jthi
164 do bi = itlo,ithi
165
166 call active_read_yz( 'maskobcse', maskyz,
167 & iobcs,
168 & doglobalread, ladinit, 0,
169 & mythid, dummy )
170
171 cgg Need to change xx field to baroclinic vel.
172 cgg Level Nr contains barotropic vel, remove it.
173 cgg The deepest wet level velocity should be
174 cgg computed in order to ensure zero integrated flux.
175 if (iobcs .eq. 3 .or. iobcs .eq. 4) then
176 c-- Loop over this thread's tiles.
177 do bj = jtlo,jthi
178 do bi = itlo,ithi
179 do j = jmin,jmax
180
181 cgg The barotropic velocity is stored in level 1.
182 cgg Penalize for deviation in it.
183 barofield(j,bi,bj) = tmpfield(j,1,bi,bj)
184 & * maskyz(j,1,bi,bj)
185
186 tmpfield(j,1,bi,bj) = 0.d0
187 utop = 0.d0
188
189 do k = 1,Nr
190 cgg If cells are not full, this should be modified with hFac.
191 cgg
192 cgg The xx field (tmpfldxz) does not contain the velocity at the
193 cgg lowest wet level. This velocity is not independent; it must
194 cgg exactly balance the volume flux, since we are dealing with
195 cgg the baroclinic velocity structure..
196 utop = tmpfield(j,k,bi,bj)*
197 & maskyz(j,k,bi,bj) * delZ(k) + utop
198 enddo
199 tmpfield(j,1,bi,bj) = -utop / delZ(1)
200 enddo
201 enddo
202 enddo
203 endif
204
205 c-- Determine the weights to be used.
206 fctile = 0. _d 0
207
208 if (iobcs .eq. 3 .or. iobcs .eq. 4) then
209 do j = jmin,jmax
210 i = OB_Ie(j,bi,bj)
211 tmpx = barofield(j,bi,bj)
212 fctile = fctile
213 & + wbaro*cosphi(i,j,bi,bj)
214 & *tmpx*tmpx
215 & *maskyz(j,1,bi,bj)
216 print*,'fctile barotropic E',fctile
217 enddo
218 endif
219
220 do k = 1, Nr
221 do j = jmin,jmax
222 i = OB_Ie(j,bi,bj)
223 cgg if (maskW(i+ip1,j,k,bi,bj) .ne. 0.) then
224 tmpx = tmpfield(j,k,bi,bj)
225 fctile = fctile
226 & + wobcseLev(j,k,bi,bj,iobcs)*cosphi(i,j,bi,bj)
227 & *tmpx*tmpx
228 & *maskyz(j,k,bi,bj)
229 cgg endif
230 enddo
231 enddo
232
233 objf_obcse(bi,bj) = objf_obcse(bi,bj) + fctile
234 fcthread = fcthread + fctile
235 enddo
236 enddo
237
238 #ifdef ECCO_VERBOSE
239 c-- Print cost function for all tiles.
240 _GLOBAL_SUM_R8( fcthread , myThid )
241 write(msgbuf,'(a)') ' '
242 call print_message( msgbuf, standardmessageunit,
243 & SQUEEZE_RIGHT , mythid)
244 write(msgbuf,'(a,i8.8)')
245 & ' cost_obcse: irec = ',irec
246 call print_message( msgbuf, standardmessageunit,
247 & SQUEEZE_RIGHT , mythid)
248 write(msgbuf,'(a,a,d22.15)')
249 & ' global cost function value',
250 & ' (obcse) = ',fcthread
251 call print_message( msgbuf, standardmessageunit,
252 & SQUEEZE_RIGHT , mythid)
253 write(msgbuf,'(a)') ' '
254 call print_message( msgbuf, standardmessageunit,
255 & SQUEEZE_RIGHT , mythid)
256 #endif
257
258 enddo
259 c-- End of loop over records.
260
261 #endif
262
263 return
264 end
265
266
267
268
269
270
271

  ViewVC Help
Powered by ViewVC 1.1.22