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

Contents of /MITgcm/pkg/ecco/cost_obcss.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_obcss(
8 I myiter,
9 I mytime,
10 I startrec,
11 I endrec,
12 I mythid
13 & )
14
15 c ==================================================================
16 c SUBROUTINE cost_obcss
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_obcss
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
71 _RL fctile
72 _RL fcthread
73 _RL dummy
74 _RL gg
75 _RL tmpx
76 _RL tmpfield (1-olx:snx+olx,nr,nsx,nsy)
77 _RL barofield (1-olx:snx+olx,nsx,nsy)
78 _RL maskxz (1-olx:snx+olx,nr,nsx,nsy)
79 _RL vtop
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_OBCSS_COST_CONTRIBUTION
114
115 jp1 = 1
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_obcss: 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_obcss_file )
138 write(fnamefld(1:80),'(2a,i10.10)')
139 & xx_obcss_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 i = imin,imax
148 barofield(i,bi,bj) = 0. _d 0
149 enddo
150 enddo
151 enddo
152
153 call active_read_xz( fnamefld, tmpfield, irec, doglobalread,
154 & ladinit, optimcycle, mythid
155 & , xx_obcss_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 cgg print*,'S IREC, IOBCS',irec, iobcs
162
163 call active_read_xz( 'maskobcss', maskxz,
164 & iobcs,
165 & doglobalread, ladinit, 0,
166 & mythid, dummy )
167
168 cgg Need to change xx field to baroclinic vel.
169 cgg Level Nr contains barotropic vel, remove it.
170 cgg The deepest wet level velocity should be
171 cgg computed in order to ensure zero integrated flux.
172 if (iobcs .eq. 3 .or. iobcs .eq. 4) then
173 c-- Loop over this thread's tiles.
174 do bj = jtlo,jthi
175 do bi = itlo,ithi
176 do i = imin,imax
177
178 cgg The barotropic velocity is stored in level 1.
179 cgg Penalize for deviation in it.
180 barofield(i,bi,bj) = tmpfield(i,1,bi,bj)
181 & * maskxz(i,1,bi,bj)
182
183 vtop = 0.d0
184 tmpfield(i,1,bi,bj) = 0.d0
185
186 do k = 1,Nr
187 cgg If cells are not full, this should be modified with hFac.
188 cgg
189 cgg The xx field (tmpfldxz) does not contain the velocity at the
190 cgg lowest wet level. This velocity is not independent; it must
191 cgg exactly balance the volume flux, since we are dealing with
192 cgg the baroclinic velocity structure..
193 vtop = tmpfield(i,k,bi,bj)*
194 & maskxz(i,k,bi,bj) * delZ(k) + vtop
195 enddo
196 tmpfield(i,1,bi,bj) = -vtop / delZ(1)
197 enddo
198 enddo
199 enddo
200
201 endif
202
203 c-- Loop over this thread's tiles.
204 do bj = jtlo,jthi
205 do bi = itlo,ithi
206
207 c-- Determine the weights to be used.
208 fctile = 0. _d 0
209
210 if (iobcs .eq. 3 .or. iobcs .eq. 4) then
211 do i = imin,imax
212 j = OB_Js(I,bi,bj)
213 tmpx = barofield(i,bi,bj)
214 fctile = fctile
215 & + wbaro*cosphi(i,j,bi,bj)
216 & *tmpx*tmpx
217 & *maskxz(i,1,bi,bj)
218 cgg print*,'S bt fctile',fctile
219 enddo
220 endif
221
222 do k = 1, Nr
223 do i = imin,imax
224 j = OB_Js(I,bi,bj)
225 cgg if (maskS(i,j+jp1,k,bi,bj) .ne. 0.) then
226 tmpx = tmpfield(i,k,bi,bj)
227 fctile = fctile
228 & + wobcssLev(i,k,bi,bj,iobcs)*cosphi(i,j,bi,bj)
229 & *tmpx*tmpx
230 & *maskxz(i,k,bi,bj)
231 cgg endif
232 cgg print*,'S fctile',fctile
233 enddo
234 enddo
235
236 objf_obcss(bi,bj) = objf_obcss(bi,bj) + fctile
237 fcthread = fcthread + fctile
238 enddo
239 enddo
240
241 #ifdef ECCO_VERBOSE
242 c-- Print cost function for all tiles.
243 _GLOBAL_SUM_R8( fcthread , myThid )
244 write(msgbuf,'(a)') ' '
245 call print_message( msgbuf, standardmessageunit,
246 & SQUEEZE_RIGHT , mythid)
247 write(msgbuf,'(a,i8.8)')
248 & ' cost_obcss: irec = ',irec
249 call print_message( msgbuf, standardmessageunit,
250 & SQUEEZE_RIGHT , mythid)
251 write(msgbuf,'(a,a,d22.15)')
252 & ' global cost function value',
253 & ' (obcss) = ',fcthread
254 call print_message( msgbuf, standardmessageunit,
255 & SQUEEZE_RIGHT , mythid)
256 write(msgbuf,'(a)') ' '
257 call print_message( msgbuf, standardmessageunit,
258 & SQUEEZE_RIGHT , mythid)
259 #endif
260
261 enddo
262 c-- End of loop over records.
263
264 #endif
265
266 return
267 end
268
269
270
271
272
273
274

  ViewVC Help
Powered by ViewVC 1.1.22