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

Contents of /MITgcm/pkg/ecco/cost_obcs_ageos.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: checkpoint58l_post, checkpoint57t_post, checkpoint57o_post, checkpoint52l_pre, checkpoint52e_pre, hrcube4, checkpoint58e_post, checkpoint57v_post, checkpoint52n_post, checkpoint52j_post, checkpoint53d_post, checkpoint58u_post, checkpoint58w_post, checkpoint54a_pre, checkpoint57m_post, checkpoint55c_post, checkpoint54e_post, checkpoint52e_post, checkpoint57s_post, checkpoint54a_post, checkpoint53c_post, checkpoint57k_post, checkpoint55d_pre, checkpoint57d_post, checkpoint57g_post, checkpoint57b_post, checkpoint57c_pre, checkpoint58r_post, checkpoint55j_post, checkpoint56b_post, checkpoint57i_post, checkpoint57y_post, hrcube_1, checkpoint57e_post, branch-netcdf, checkpoint52d_pre, checkpoint52l_post, checkpoint55h_post, checkpoint58n_post, checkpoint58x_post, checkpoint52k_post, checkpoint52b_pre, checkpoint57g_pre, checkpoint54b_post, checkpoint53b_pre, checkpoint55b_post, checkpoint58t_post, checkpoint58h_post, checkpoint54d_post, checkpoint56c_post, checkpoint52m_post, checkpoint57y_pre, checkpoint55, checkpoint53a_post, checkpoint57f_pre, checkpoint57a_post, checkpoint54, checkpoint58q_post, checkpoint54f_post, checkpoint53b_post, checkpoint55g_post, checkpoint58j_post, checkpoint52a_pre, checkpoint59e, checkpoint59d, checkpoint59g, checkpoint59f, checkpoint59a, checkpoint55f_post, checkpoint59c, checkpoint59b, checkpoint59h, checkpoint57r_post, checkpoint59, checkpoint58, checkpoint57a_pre, checkpoint55i_post, checkpoint57, checkpoint56, checkpoint53, checkpoint52, checkpoint52d_post, eckpoint57e_pre, checkpoint52a_post, checkpoint57h_done, checkpoint58f_post, checkpoint52b_post, checkpoint53g_post, checkpoint52f_post, checkpoint57x_post, checkpoint57n_post, checkpoint52c_post, checkpoint58d_post, checkpoint58c_post, checkpoint57w_post, checkpoint57p_post, checkpint57u_post, checkpoint57f_post, checkpoint58a_post, checkpoint58i_post, checkpoint57q_post, checkpoint58g_post, ecco_c52_e35, hrcube5, checkpoint58o_post, checkpoint57z_post, checkpoint57c_post, checkpoint58y_post, checkpoint55e_post, checkpoint58k_post, checkpoint52i_post, checkpoint52j_pre, checkpoint58v_post, checkpoint53f_post, checkpoint55a_post, checkpoint53d_pre, checkpoint54c_post, checkpoint58s_post, checkpoint58p_post, checkpoint57j_post, checkpoint58b_post, checkpoint57h_pre, checkpoint58m_post, checkpoint57l_post, checkpoint52i_pre, checkpoint51u_post, checkpoint52h_pre, checkpoint52f_pre, checkpoint57h_post, hrcube_2, hrcube_3, checkpoint56a_post, checkpoint55d_post
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 C $Header: /u/gcmpack/MITgcm/pkg/cost/Attic/cost_obcs_ageos.F,v 1.1.2.2 2003/07/07 16:18:17 heimbach Exp $
2
3 #include "COST_CPPOPTIONS.h"
4 #ifdef ALLOW_OBCS
5 # include "OBCS_OPTIONS.h"
6 #endif
7
8 subroutine cost_obcs_ageos( myiter, mytime, mythid )
9
10 c ==================================================================
11 c SUBROUTINE cost_obcs_ageos
12 c ==================================================================
13 c
14 c o cost function contribution obc -- Ageostrophic boundary flow.
15 c
16 c started: G. Gebbie gebbie@mit.edu 4-Feb-2003
17 c
18 c warning: masks redundantly assume no gradient of topography at
19 c boundary.
20 c
21 c ==================================================================
22 c SUBROUTINE cost_obcs_ageos
23 c ==================================================================
24
25 implicit none
26
27 c == global variables ==
28
29 #include "EEPARAMS.h"
30 #include "SIZE.h"
31 #include "GRID.h"
32 #include "DYNVARS.h"
33 #include "PARAMS.h"
34 #ifdef ALLOW_OBCS
35 # include "OBCS.h"
36 #endif
37
38 #include "cal.h"
39 #include "ecco_cost.h"
40 #include "ctrl.h"
41 #include "ctrl_dummy.h"
42 #include "optim.h"
43
44 c == routine arguments ==
45
46 integer myiter
47 _RL mytime
48 integer mythid
49
50 c == local variables ==
51
52 _RS one_rs
53 parameter( one_rs = 1. )
54
55 integer bi,bj
56 integer i,j,k
57 integer itlo,ithi
58 integer jtlo,jthi
59 integer jmin,jmax
60 integer imin,imax
61 integer irec
62 integer levmon
63 integer levoff
64 integer iltheta
65 integer ilsalt
66 integer iluvel
67 integer ilvvel
68 integer ip1, jp1
69
70 _RL fctile
71 _RL fcthread
72
73 _RL rholoc (1-olx:snx+olx,1-oly:sny+oly,nsx,nsy)
74 _RL xzgrdrho(1-olx:snx+olx,Nr,nsx,nsy)
75 _RL yzgrdrho(1-oly:sny+oly,Nr,nsx,nsy)
76 _RL xzdvel1 (1-olx:snx+olx,nr,nsx,nsy)
77 _RL xzdvel2 (1-olx:snx+olx,nr,nsx,nsy)
78 _RL yzdvel1 (1-oly:sny+oly,nr,nsx,nsy)
79 _RL yzdvel2 (1-oly:sny+oly,nr,nsx,nsy)
80 _RL maskxzageos (1-olx:snx+olx,nr,nsx,nsy)
81 _RL maskyzageos (1-oly:sny+oly,nr,nsx,nsy)
82 _RL dummy
83
84 character*(80) fnametheta
85 character*(80) fnamesalt
86 character*(80) fnameuvel
87 character*(80) fnamevvel
88
89 logical doglobalread
90 logical ladinit
91
92 character*(MAX_LEN_MBUF) msgbuf
93
94 c == external functions ==
95
96 integer ilnblnk
97 external ilnblnk
98
99 c == end of interface ==
100
101 jtlo = mybylo(mythid)
102 jthi = mybyhi(mythid)
103 itlo = mybxlo(mythid)
104 ithi = mybxhi(mythid)
105 jmin = 1
106 jmax = sny
107 imin = 1
108 imax = snx
109
110 c-- Read tiled data.
111 doglobalread = .false.
112 ladinit = .false.
113
114 #ifdef OBCS_AGEOS_COST_CONTRIBUTION
115
116 #ifdef ECCO_VERBOSE
117 _BEGIN_MASTER( mythid )
118 write(msgbuf,'(a)') ' '
119 call print_message( msgbuf, standardmessageunit,
120 & SQUEEZE_RIGHT , mythid)
121 write(msgbuf,'(a,i8.8)')
122 & ' cost_obcs_ageos: number of records to process = ',nmonsrec
123 call print_message( msgbuf, standardmessageunit,
124 & SQUEEZE_RIGHT , mythid)
125 write(msgbuf,'(a)') ' '
126 call print_message( msgbuf, standardmessageunit,
127 & SQUEEZE_RIGHT , mythid)
128 _END_MASTER( mythid )
129 #endif
130
131 if (optimcycle .ge. 0) then
132 iltheta = ilnblnk( tbarfile )
133 write(fnametheta(1:80),'(2a,i10.10)')
134 & tbarfile(1:iltheta),'.',optimcycle
135 ilsalt = ilnblnk( sbarfile )
136 write(fnamesalt(1:80),'(2a,i10.10)')
137 & sbarfile(1:ilsalt),'.',optimcycle
138 iluvel = ilnblnk( ubarfile )
139 write(fnameuvel(1:80),'(2a,i10.10)')
140 & ubarfile(1:iluvel),'.',optimcycle
141 ilvvel = ilnblnk( vbarfile )
142 write(fnamevvel(1:80),'(2a,i10.10)')
143 & vbarfile(1:ilvvel),'.',optimcycle
144 endif
145
146 fcthread = 0. _d 0
147 fctile = 0. _d 0
148
149 cgg Code safe: always initialize to zero.
150 do bj = jtlo,jthi
151 do bi = itlo,ithi
152 do k = 1,nr
153 do i = 1-olx,snx+olx
154 maskxzageos(i,k,bi,bj) = 0. _d 0
155 xzdvel1(i,k,bi,bj) = 0. _d 0
156 xzdvel2(i,k,bi,bj) = 0. _d 0
157 xzgrdrho(i,k,bi,bj) = 0. _d 0
158 enddo
159 do j = 1-oly,sny+oly
160 maskyzageos(j,k,bi,bj) = 0. _d 0
161 yzdvel1(j,k,bi,bj) = 0. _d 0
162 yzdvel2(j,k,bi,bj) = 0. _d 0
163 yzgrdrho(j,k,bi,bj) = 0. _d 0
164 enddo
165 enddo
166 enddo
167 enddo
168
169 do bj = jtlo,jthi
170 do bi = itlo,ithi
171 do j = 1-oly,sny+oly
172 do i = 1-olx,snx+olx
173 rholoc(i,j,bi,bj) = 0. _d 0
174 enddo
175 enddo
176 enddo
177 enddo
178
179 c-- Loop over records.
180 do irec = 1,nmonsrec
181
182 c-- Read time averages and the monthly mean data.
183 call active_read_xyz( fnametheta, tbar, irec,
184 & doglobalread, ladinit,
185 & optimcycle, mythid,
186 & xx_tbar_mean_dummy )
187
188 c-- Read time averages and the monthly mean data.
189 call active_read_xyz( fnamesalt, sbar, irec,
190 & doglobalread, ladinit,
191 & optimcycle, mythid,
192 & xx_sbar_mean_dummy )
193
194 c-- Read time averages and the monthly mean data.
195 call active_read_xyz( fnameuvel, ubar, irec,
196 & doglobalread, ladinit,
197 & optimcycle, mythid,
198 & xx_ubar_mean_dummy )
199
200 c-- Read time averages and the monthly mean data.
201 call active_read_xyz( fnamevvel, vbar, irec,
202 & doglobalread, ladinit,
203 & optimcycle, mythid,
204 & xx_vbar_mean_dummy )
205
206 cgg Minor problem : grad T,S needs overlap values.
207 _BARRIER
208 _EXCH_XYZ_R8(tbar,myThid)
209 _EXCH_XYZ_R8(sbar,myThid)
210
211
212 do bj = jtlo,jthi
213 do bi = itlo,ithi
214
215 #ifdef ALLOW_OBCSN_CONTROL
216 jp1 = 0
217
218 cgg Make a mask for the velocity shear comparison.
219 do k = 1,nr-1
220 do i = imin, imax
221 j = ob_jn(i,bi,bj)
222 cgg All these points need to be wet.
223 if (j .eq. 0) then
224 maskxzageos(i,k,bi,bj) = 0.
225 else
226 maskxzageos(i,k,bi,bj) =
227 & hfacC(i,j+jp1,k,bi,bj)*hfacC(i+1,j+jp1,k,bi,bj) *
228 & hfacC(i-1,j+jp1,k,bi,bj)*hfacC(i,j+jp1,k+1,bi,bj)*
229 & hfacC(i-1,j+jp1,k+1,bi,bj)*hfacC(i+1,j+jp1,k+1,bi,bj)*
230 & hfacS(i,j+jp1,k,bi,bj)*hfacS(i,j+jp1,k+1,bi,bj)
231 endif
232 enddo
233 enddo
234
235 do k = 1,nr
236 call find_rho(bi,bj,imin,imax,jmin,jmax,k,k,
237 & tbar,sbar,rholoc,mythid)
238 _EXCH_XY_R8(rholoc , myThid)
239
240 cgg Compute centered difference horizontal gradient on bdy.
241 do i = imin, imax
242 j = ob_jn(i,bi,bj)
243 xzgrdrho(i,k,bi,bj) =
244 & (rholoc(i-1,j+jp1,bi,bj)-rholoc(i+1,j+jp1,bi,bj))
245 & /(2.*dxc(i,j+jp1,bi,bj))
246 enddo
247 enddo
248
249 cgg Compute vertical shear from geostrophy/thermal wind.
250 cgg Above level 4 needs not be geostrophic.
251 cgg Please get rid of the "4" ASAP. Ridiculous.
252 do k = 4,Nr-1
253 do i = imin,imax
254 j = ob_jn(i,bi,bj)
255 xzdvel1(i,k,bi,bj) = vbar(i,j+jp1,k ,bi,bj)
256 & - vbar(i,j+jp1,k+1,bi,bj)
257 xzdvel2(i,k,bi,bj)=((xzgrdrho(i,k,bi,bj)*delz(k)/2.)+
258 & (xzgrdrho(i,k+1,bi,bj)*delz(k+1)/2.))
259 & * gravity / f0 / rhonil
260
261 fctile = fctile + 100*wcurrent(k,bi,bj) *
262 & maskxzageos(i,k,bi,bj)*
263 & (xzdvel2(i,k,bi,bj) - xzdvel1(i,k,bi,bj))*
264 & (xzdvel2(i,k,bi,bj) - xzdvel1(i,k,bi,bj))
265 if (maskxzageos(i,k,bi,bj) .ne. 0) then
266 cgg print*,'XX shear,geos shear',xzdvel1(i,k,bi,bj),xzdvel2(i,k,bi,bj)
267 cgg print*,'i,j,k,fctile N',i,j,k,fctile
268 endif
269 enddo
270 enddo
271 c-- End of loop over layers.
272 #endif /* ALLOW_OBCSN_CONTROL */
273
274 #ifdef ALLOW_OBCSS_CONTROL
275 jp1 = 1
276
277 cgg Make a mask for the velocity shear comparison.
278 do k = 1,nr-1
279 do i = imin, imax
280 j = ob_js(i,bi,bj)
281 if (j .eq. 0) then
282 maskxzageos(i,k,bi,bj) = 0.
283 else
284 cgg All these points need to be wet.
285 maskxzageos(i,k,bi,bj) =
286 & hfacC(i,j+jp1,k,bi,bj)*hfacC(i+1,j+jp1,k,bi,bj) *
287 & hfacC(i-1,j+jp1,k,bi,bj)*hfacC(i,j+jp1,k+1,bi,bj)*
288 & hfacC(i-1,j+jp1,k+1,bi,bj)*hfacC(i+1,j+jp1,k+1,bi,bj)*
289 & hfacS(i,j+jp1,k,bi,bj)*hfacS(i,j+jp1,k+1,bi,bj)
290 endif
291 enddo
292 enddo
293
294 do k = 1,nr
295
296 call find_rho(bi,bj,imin,imax,jmin,jmax,k,k,
297 & tbar,sbar,rholoc,mythid)
298 _EXCH_XY_R8(rholoc , myThid)
299
300 cgg Compute centered difference horizontal gradient on bdy.
301 do i = imin, imax
302 j = ob_js(i,bi,bj)
303 xzgrdrho(i,k,bi,bj) =
304 & (rholoc(i-1,j+jp1,bi,bj)-rholoc(i+1,j+jp1,bi,bj))
305 & /(2.*dxc(i,j+jp1,bi,bj))
306 enddo
307 enddo
308
309 cgg Compute vertical shear from geostrophy/thermal wind.
310 do k = 4,Nr-1
311 do i = imin,imax
312 j = ob_js(i,bi,bj)
313 cgg Retrieve the model vertical shear.
314 xzdvel1(i,k,bi,bj) = vbar(i,j+jp1,k ,bi,bj)
315 & - vbar(i,j+jp1,k+1,bi,bj)
316
317 cgg Compute vertical shear from geostrophy/thermal wind.
318 xzdvel2(i,k,bi,bj) =((xzgrdrho(i,k,bi,bj)*delz(k)/2.)+
319 & (xzgrdrho(i,k+1,bi,bj)*delz(k+1)/2.))
320 & * gravity / f0 /rhonil
321
322 cgg Make a comparison.
323 fctile = fctile + 100*wcurrent(k,bi,bj) *
324 & maskxzageos(i,k,bi,bj)*
325 & (xzdvel2(i,k,bi,bj) - xzdvel1(i,k,bi,bj))*
326 & (xzdvel2(i,k,bi,bj) - xzdvel1(i,k,bi,bj))
327 cgg print*,'fctile S',fctile
328 enddo
329 enddo
330 c-- End of loop over layers.
331 #endif
332
333 #ifdef ALLOW_OBCSW_CONTROL
334 ip1 = 1
335
336 cgg Make a mask for the velocity shear comparison.
337 do k = 1,nr-1
338 do j = jmin, jmax
339 i = ob_iw(j,bi,bj)
340 cgg All these points need to be wet.
341 if (i.eq. 0) then
342 maskyzageos(j,k,bi,bj) = 0.
343 else
344 maskyzageos(j,k,bi,bj) =
345 & hfacC(i+ip1,j,k,bi,bj)*hfacC(i+ip1,j+1,k,bi,bj) *
346 & hfacC(i+ip1,j-1,k,bi,bj)*hfacC(i+ip1,j,k+1,bi,bj)*
347 & hfacC(i+ip1,j-1,k+1,bi,bj)*hfacC(i+ip1,j+1,k+1,bi,bj)*
348 & hfacW(i+ip1,j,k,bi,bj)*hfacW(i+ip1,j,k+1,bi,bj)
349 endif
350 enddo
351 enddo
352
353 do k = 1,nr
354
355 call find_rho(bi,bj,imin,imax,jmin,jmax,k,k,
356 & tbar,sbar,rholoc,mythid)
357 _EXCH_XY_R8(rholoc , myThid)
358
359 cgg Compute centered difference horizontal gradient on bdy.
360 do j = jmin, jmax
361 i = ob_iw(j,bi,bj)
362 cgg Negative sign due to geostrophy.
363 yzgrdrho(j,k,bi,bj) =
364 & (rholoc(i+ip1,j+1,bi,bj)-rholoc(i+ip1,j-1,bi,bj))
365 & /(2.*dyc(i+ip1,j,bi,bj))
366 enddo
367 enddo
368
369 cgg Compute vertical shear from geostrophy/thermal wind.
370 do k = 4,Nr-1
371 do j = jmin,jmax
372 i = ob_iw(j,bi,bj)
373 cgg Retrieve the model vertical shear.
374 yzdvel1(j,k,bi,bj) = ubar(i+ip1,j,k ,bi,bj)
375 & - ubar(i+ip1,j,k+1,bi,bj)
376
377 cgg Compute vertical shear from geostrophy/thermal wind.
378 yzdvel2(j,k,bi,bj) =((yzgrdrho(j,k ,bi,bj)*delz(k)/2.)+
379 & (yzgrdrho(j,k+1,bi,bj)*delz(k+1)/2.))
380 & * gravity / f0 / rhonil
381
382 cgg Make a comparison.
383 fctile = fctile + 100*wcurrent(k,bi,bj) *
384 & maskyzageos(j,k,bi,bj) *
385 & (yzdvel2(j,k,bi,bj) - yzdvel1(j,k,bi,bj))*
386 & (yzdvel2(j,k,bi,bj) - yzdvel1(j,k,bi,bj))
387 enddo
388 enddo
389 c-- End of loop over layers.
390 #endif
391
392 #ifdef ALLOW_OBCSE_CONTROL
393 ip1 = 0
394
395 cgg Make a mask for the velocity shear comparison.
396 do k = 1,nr-1
397 do j = jmin, jmax
398 i = ob_ie(j,bi,bj)
399 if (i.eq.0) then
400 maskyzageos(j,k,bi,bj) =0.
401 else
402 cgg All these points need to be wet.
403 maskyzageos(j,k,bi,bj) =
404 & hfacC(i+ip1,j,k,bi,bj)*hfacC(i+ip1,j+1,k,bi,bj) *
405 & hfacC(i+ip1,j-1,k,bi,bj)*hfacC(i+ip1,j,k+1,bi,bj)*
406 & hfacC(i+ip1,j-1,k+1,bi,bj)*hfacC(i+ip1,j+1,k+1,bi,bj)*
407 & hfacW(i+ip1,j,k,bi,bj)*hfacW(i+ip1,j,k+1,bi,bj)
408 endif
409 enddo
410 enddo
411
412 do k = 1,nr
413
414 call find_rho(bi,bj,imin,imax,jmin,jmax,k,k,
415 & tbar,sbar,rholoc,mythid)
416 _EXCH_XY_R8(rholoc , myThid)
417
418 cgg Compute centered difference horizontal gradient on bdy.
419 do j = jmin, jmax
420 i = ob_ie(j,bi,bj)
421 cgg Negative sign due to geostrophy.
422 yzgrdrho(j,k,bi,bj) =
423 & (rholoc(i+ip1,,j+1,bi,bj)-rholoc(i+ip1,j-1,bi,bj))
424 & /(2.*dyc(i+ip1,j,bi,bj))
425 enddo
426 enddo
427
428 cgg Compute vertical shear from geostrophy/thermal wind.
429 do k = 4,Nr-1
430 do j = jmin,jmax
431 i = ob_ie(j,bi,bj)
432 cgg Retrieve the model vertical shear.
433 yzdvel1(j,k,bi,bj) = ubar(i+ip1,j,k ,bi,bj)
434 & - ubar(i+ip1,j,k+1,bi,bj)
435
436 cgg Compute vertical shear from geostrophy/thermal wind.
437 yzdvel2(j,k,bi,bj) =((yzgrdrho(j,k ,bi,bj)*delz(k)/2.)+
438 & (yzgrdrho(j,k+1,bi,bj)*delz(k+1)/2.))
439 & * gravity / f0 /rhonil
440
441 cgg Make a comparison.
442 fctile = fctile + 100*wcurrent(k,bi,bj) *
443 & maskyzageos(j,k,bi,bj) *
444 & (yzdvel2(j,k,bi,bj) - yzdvel1(j,k,bi,bj))*
445 & (yzdvel2(j,k,bi,bj) - yzdvel1(j,k,bi,bj))
446
447 enddo
448 enddo
449 c-- End of loop over layers.
450 #endif
451
452 fcthread = fcthread + fctile
453 objf_ageos(bi,bj) = objf_ageos(bi,bj) + fctile
454
455 #ifdef ECCO_VERBOSE
456 c-- Print cost function for each tile in each thread.
457 write(msgbuf,'(a)') ' '
458 call print_message( msgbuf, standardmessageunit,
459 & SQUEEZE_RIGHT , mythid)
460 write(msgbuf,'(a,i8.8,1x,i3.3,1x,i3.3)')
461 & ' cost_Theta: irec,bi,bj = ',irec,bi,bj
462 call print_message( msgbuf, standardmessageunit,
463 & SQUEEZE_RIGHT , mythid)
464 write(msgbuf,'(a,d22.15)')
465 & ' cost function (temperature) = ',
466 & fctile
467 call print_message( msgbuf, standardmessageunit,
468 & SQUEEZE_RIGHT , mythid)
469 write(msgbuf,'(a)') ' '
470 call print_message( msgbuf, standardmessageunit,
471 & SQUEEZE_RIGHT , mythid)
472 #endif
473
474 enddo
475 enddo
476
477 #ifdef ECCO_VERBOSE
478 c-- Print cost function for all tiles.
479 _GLOBAL_SUM_R8( fcthread , myThid )
480 write(msgbuf,'(a)') ' '
481 call print_message( msgbuf, standardmessageunit,
482 & SQUEEZE_RIGHT , mythid)
483 write(msgbuf,'(a,i8.8)')
484 & ' cost_Theta: irec = ',irec
485 call print_message( msgbuf, standardmessageunit,
486 & SQUEEZE_RIGHT , mythid)
487 write(msgbuf,'(a,a,d22.15)')
488 & ' global cost function value',
489 & ' (temperature) = ',fcthread
490 call print_message( msgbuf, standardmessageunit,
491 & SQUEEZE_RIGHT , mythid)
492 write(msgbuf,'(a)') ' '
493 call print_message( msgbuf, standardmessageunit,
494 & SQUEEZE_RIGHT , mythid)
495 #endif
496
497 enddo
498 c-- End of loop over records.
499
500 #endif
501
502 return
503 end
504
505
506

  ViewVC Help
Powered by ViewVC 1.1.22