source: trunk/SRC/Grid/ncdf_meshread.pro @ 302

Last change on this file since 302 was 302, checked in by smasson, 17 years ago

compatibility with new meshmask using gdept_O, ...

  • Property svn:keywords set to Id
File size: 23.6 KB
Line 
1;+
2;
3; @file_comments
4; read NetCDF meshmask file created by OPA
5;
6; @categories
7; Grid
8;
9; @examples
10; IDL> ncdf_meshread [,' filename']
11;
12; @param filename {in}{optional}{default='meshmask.nc'}{type=scalar string}
13;    Name of the meshmask file to read. If this name does not contain any "/"
14;    and if iodirectory keyword is not specified, then the common variable
15;    iodir will be used to define the mesh file path.
16;
17; @keyword GLAMBOUNDARY {default=those defined in the file}{type=2 elements vector}
18;    Longitude boundaries that should be used to visualize the data.
19;      lon2 > lon1
20;      lon2 - lon1 le 360
21;    By default, the common (cm_4mesh) variable key_shift will be automatically
22;    defined according to GLAMBOUNDARY.
23;
24; @keyword CHECKDAT
25; Suppressed. Use <pro>micromeshmask</pro> to create an appropriate meshmask.
26;
27; @keyword ONEARTH {default=1}{type=scalar: 0 or 1}
28;    Force the manual definition of data localization on the earth or not
29;       0) if the data are not on the earth
30;       1) if the data are on earth (in that case we can for example use
31;          the labels 'longitude', 'latitude' in plots).
32;    The resulting value will be stored in the common (cm_4mesh) variable key_onearth
33;    ONEARTH = 0 forces PERIODIC = 0, SHIFT = 0 and is cancelling GLAMBOUNDARY
34;
35; @keyword GETDIMENSIONS {default=0}{type=scalar: 0 or 1}
36;    Activate this keywords if you only want to know the dimension
37;    of the domain stored in the mesh file. This dimension will be
38;    defined in jpiglo, jpjglo, jpkglo (cm_4mesh common variables)
39;
40; @keyword PERIODIC {default=computed by using the first line of glamt}{type=scalar: 0 or 1}
41;    Force the manual definition of the grid zonal periodicity.
42;    The resulting value will be stored in the common (cm_4mesh) variable key_periodic
43;    PERIODIC = 0 forces SHIFT = 0
44;
45; @keyword SHIFT {default=computed according to glamboundary}{type=scalar}
46;    Force the manual definition of the zonal shift that must be apply to the data.
47;    The resulting value will be stored in the common (cm_4mesh) variable key_shift
48;    Note that if key_periodic=0 then in any case key_shift = 0.
49;
50; @keyword STRCALLING {type=scalar string}
51;    the calling command used to call <pro>computegrid</pro> (this is used by <pro>xxx</pro>)
52;
53; @keyword STRIDE {default=[1, 1, 1]}{type=3 elements vector}
54;    Specify the stride in x, y and z direction. The resulting
55;    value will be stored in the common (cm_4mesh) variable key_stride
56;
57; @keyword _EXTRA
58; Used to pass keywords to <pro>isafile</pro>
59;
60; @uses
61; cm_4mesh
62; cm_4data
63; cm_4cal
64;
65; @restrictions
66; ixminmesh, ixmaxmesh, iyminmesh, iymaxmesh, izminmesh, izmaxmesh must
67; be defined before calling ncdf_meshread. If some of those values
68; are equal to -1 they will be automatically defined
69;
70; @history
71; Sebastien Masson (smasson\@lodyc.jussieu.fr)
72;                      12/1999
73; July 2004, Sebastien Masson: Several modifications (micromeshmask,
74; clean partial steps, clean use of key_stride, automatic definition
75; of key_shift, ...)
76; Oct. 2004, Sebastien Masson: add PERIODIC and SHIFT
77; Aug. 2005, Sebastien Masson: some cleaning + english
78;
79; @version
80; $Id$
81;
82;-
83;
84PRO ncdf_meshread, filename, GLAMBOUNDARY = glamboundary, CHECKDAT = checkdat $
85                  , ONEARTH = onearth, GETDIMENSIONS = getdimensions $
86                  , PERIODIC = periodic, SHIFT = shift, STRIDE = stride $
87                  , STRCALLING = strcalling, _EXTRA = ex
88;
89  compile_opt idl2, strictarrsubs
90;
91@cm_4mesh
92@cm_4data
93@cm_4cal
94  IF NOT keyword_set(key_forgetold) THEN BEGIN
95@updatenew
96@updatekwd
97  ENDIF
98;---------------------------------------------------------
99;
100  tempsun = systime(1)          ; for key_performance
101  IF keyword_set(CHECKDAT) THEN BEGIN
102    ras = report([' The keyword CHECKDAT has been suppressed (it could create bugs).', $
103    ' Remove it from the call of ncdf_meshread', $
104    ' Please use smallmeshmask.pro or micromeshmask.pro to create a', $
105    ' meshmask that has manageable size'])
106    return
107  ENDIF
108;-------------------------------------------------------
109; find meshfile name and open it!
110;-------------------------------------------------------
111; def of filename by default
112  IF n_params() EQ 0 then filename = 'meshmask.nc'
113  meshname = isafile(file = filename, iodirectory = iodir, _EXTRA = ex)
114  meshname = meshname[0]
115;
116  noticebase = xnotice('Reading file !C '+meshname+'!C ...')
117; if the meshmask is on tape archive ... get it back
118  IF !version.OS_FAMILY EQ 'unix' THEN spawn, '\file '+meshname+' > /dev/null'
119  cdfid = ncdf_open(meshname)
120  inq = ncdf_inquire(cdfid)
121;------------------------------------------------------------
122; dimensions
123;------------------------------------------------------------
124  ncdf_diminq, cdfid, 'x', name, jpiglo
125  ncdf_diminq, cdfid, 'y', name, jpjglo
126  listdims = strlowcase(ncdf_listdims(cdfid))
127  IF (where(listdims EQ 'z'))[0] NE -1 THEN ncdf_diminq, cdfid, 'z', name, jpkglo ELSE BEGIN
128    dimid = (where(strmid(listdims, 0, 5) EQ 'depth'))[0]
129    IF dimid NE -1 THEN ncdf_diminq, cdfid, dimid, name, jpkglo ELSE BEGIN
130      dummy = report('We could not find the vertical dimension..., its name must be z or start with depth')
131      stop
132    ENDELSE
133  ENDELSE
134;
135  if keyword_set(getdimensions) then begin
136    widget_control, noticebase, bad_id = nothing, /destroy
137    ncdf_close,  cdfid
138    return
139  endif
140;-------------------------------------------------------
141; check that all i[xyz]min[ax]mesh are well defined
142;-------------------------------------------------------
143  if n_elements(ixminmesh) EQ 0 THEN ixminmesh = 0
144  if n_elements(ixmaxmesh) EQ 0 then ixmaxmesh = jpiglo-1
145  if ixminmesh EQ -1 THEN ixminmesh = 0
146  IF ixmaxmesh EQ -1 then ixmaxmesh = jpiglo-1
147  if n_elements(iyminmesh) EQ 0 THEN iyminmesh = 0
148  IF n_elements(iymaxmesh) EQ 0 then iymaxmesh = jpjglo-1
149  if iyminmesh EQ -1 THEN iyminmesh = 0
150  IF iymaxmesh EQ -1 then iymaxmesh = jpjglo-1
151  if n_elements(izminmesh) EQ 0 THEN izminmesh = 0
152  IF n_elements(izmaxmesh) EQ 0 then izmaxmesh = jpkglo-1
153  if izminmesh EQ -1 THEN izminmesh = 0
154  IF izmaxmesh EQ -1 then izmaxmesh = jpkglo-1
155; definition of jpi,jpj,jpj
156  jpi = long(ixmaxmesh-ixminmesh+1)
157  jpj = long(iymaxmesh-iyminmesh+1)
158  jpk = long(izmaxmesh-izminmesh+1)
159;-------------------------------------------------------
160; check onearth and its consequences
161;-------------------------------------------------------
162  IF n_elements(onearth) EQ 0 THEN key_onearth = 1 $
163  ELSE key_onearth = keyword_set(onearth)
164  IF NOT key_onearth THEN BEGIN
165    periodic = 0
166    shift = 0
167  ENDIF
168;-------------------------------------------------------
169; automatic definition of key_periodic
170;-------------------------------------------------------
171  IF n_elements(periodic) EQ 0 THEN BEGIN
172    IF jpi GT 1 THEN BEGIN
173      varinq = ncdf_varinq(cdfid, 'glamt')
174      CASE varinq.ndims OF
175        2:ncdf_varget, cdfid, 'glamt', xaxis $
176                       , offset = [ixminmesh, iyminmesh], count = [jpi, 1]
177        3:ncdf_varget, cdfid, 'glamt', xaxis $
178                       , offset = [ixminmesh, iyminmesh, 0], count = [jpi, 1, 1]
179        4:ncdf_varget, cdfid, 'glamt', xaxis $
180                       , offset = [ixminmesh, iyminmesh, 0, 0], count = [jpi, 1, 1, 1]
181      ENDCASE
182      xaxis = (xaxis+720) MOD 360
183      xaxis = xaxis[sort(xaxis)]
184      key_periodic = (xaxis[jpi-1]+2*(xaxis[jpi-1]-xaxis[jpi-2])) $
185                     GE (xaxis[0]+360)
186    ENDIF ELSE key_periodic = 0
187  ENDIF ELSE key_periodic = keyword_set(periodic)
188;-------------------------------------------------------
189; automatic definition of key_shift
190;-------------------------------------------------------
191  IF n_elements(shift) EQ 0 THEN BEGIN
192    key_shift = long(testvar(var = key_shift))
193;  key_shift will be defined according to the first line of glamt.
194    if keyword_set(glamboundary) AND jpi GT 1 AND key_periodic EQ 1 $
195    THEN BEGIN
196      varinq = ncdf_varinq(cdfid, 'glamt')
197      CASE varinq.ndims OF
198        2:ncdf_varget, cdfid, 'glamt', xaxis $
199                       , offset = [ixminmesh, iyminmesh], count = [jpi, 1]
200        3:ncdf_varget, cdfid, 'glamt', xaxis $
201                       , offset = [ixminmesh, iyminmesh, 0], count = [jpi, 1, 1]
202        4:ncdf_varget, cdfid, 'glamt', xaxis $
203                       , offset = [ixminmesh, iyminmesh, 0, 0], count = [jpi, 1, 1, 1]
204      ENDCASE
205; xaxis between glamboundary[0] and glamboundary[1]
206      xaxis = xaxis MOD 360
207      smaller = where(xaxis LT glamboundary[0])
208      if smaller[0] NE -1 then xaxis[smaller] = xaxis[smaller]+360
209      bigger = where(xaxis GE glamboundary[1])
210      if bigger[0] NE -1 then xaxis[bigger] = xaxis[bigger]-360
211;
212      key_shift = (where(xaxis EQ min(xaxis)))[0]
213      IF key_shift NE 0 THEN BEGIN
214        key_shift = jpi-key_shift
215        xaxis = shift(xaxis, key_shift)
216      ENDIF
217;
218      IF array_equal(sort(xaxis), lindgen(jpi)) NE 1 THEN BEGIN
219        ras = report(['the x axis (1st line of glamt) is not sorted in the increasing order after the automatic definition of key_shift', $
220        'please use the keyword shift (and periodic) to suppress the automatic definition of key_shift (and key_periodic) and define by hand a more suitable value...'])
221        widget_control, noticebase, bad_id = nothing, /destroy
222        return
223      ENDIF
224;
225    ENDIF ELSE key_shift = 0
226  ENDIF ELSE key_shift = long(shift)*(key_periodic EQ 1)
227;-------------------------------------------------------
228; check key_stride and related things
229;-------------------------------------------------------
230  if n_elements(stride) eq 3 then key_stride = stride
231  if n_elements(key_stride) LE 2 then key_stride = [1, 1, 1]
232  key_stride = 1l > long(key_stride)
233  IF total(key_stride) NE 3  THEN BEGIN
234    IF key_shift NE 0 THEN BEGIN
235; for explanation, see header of read_ncdf_varget.pro
236      jpiright = key_shift
237      jpileft = jpi - key_shift - ( (key_stride[0]-1)-((key_shift-1) MOD key_stride[0]) )
238      jpi = ((jpiright-1)/key_stride[0]+1) + ((jpileft-1)/key_stride[0]+1)
239    ENDIF ELSE jpi = (jpi-1)/key_stride[0]+1
240    jpj = (jpj-1)/key_stride[1]+1
241    jpk = (jpk-1)/key_stride[2]+1
242  ENDIF
243;-------------------------------------------------------
244; default definitions to be able to use read_ncdf_varget
245;-------------------------------------------------------
246; default definitions to be able to use read_ncdf_varget
247  ixmindtasauve = testvar(var = ixmindta)
248  iymindtasauve = testvar(var = iymindta)
249  izmindtasauve = testvar(var = izmindta)
250;
251  ixmindta = 0l
252  iymindta = 0l
253  izmindta = 0l
254;
255  jpt = 1
256  time = 1
257  firsttps = 0
258;
259  firstx = 0
260  lastx = jpi-1
261  firsty = 0
262  lasty = jpj-1
263  firstz = 0
264  lastz = jpk-1
265  nx = jpi
266  ny = jpj
267  nz = 1
268  izminmeshsauve = izminmesh
269  izminmesh = 0
270;
271  key_yreverse = 0
272  key_zreverse = 0
273  key_gridtype = 'c'
274;-------------------------------------------------------
275; 2d arrays:
276;-------------------------------------------------------
277; list the 2d variables that must be read
278  namevar = ['glamt', 'glamu', 'glamv', 'glamf' $
279             , 'gphit', 'gphiu', 'gphiv', 'gphif' $
280             , 'e1t', 'e1u', 'e1v', 'e1f' $
281             , 'e2t', 'e2u', 'e2v', 'e2f']
282; for the variables related to the partial steps
283  allvarname =  ncdf_listvars(cdfid)
284  IF (where(allvarname EQ 'hdept'))[0] NE -1 THEN BEGIN
285    key_partialstep = 1
286    namevar = [namevar, 'hdept', 'hdepw']
287  ENDIF ELSE BEGIN
288    key_partialstep = 0
289    hdept = -1
290    hdepw = -1
291  ENDELSE
292; for compatibility with old versions of meshmask/partial steps
293  IF (where(allvarname EQ 'e3tp'))[0] NE -1 THEN $
294    namevar = [namevar, 'e3tp', 'e3wp'] ELSE BEGIN
295    e3t_ps = -1
296    e3w_ps = -1
297  ENDELSE
298  IF (where(allvarname EQ 'e3t_ps'))[0] NE -1 $
299  THEN namevar = [namevar, 'e3t_ps', 'e3w_ps' ]ELSE BEGIN
300    e3t_ps = -1
301    e3w_ps = -1
302  ENDELSE
303  IF (where(allvarname EQ 'e3u_ps'))[0] NE -1 $
304  THEN namevar = [namevar, 'e3u_ps', 'e3v_ps'] ELSE BEGIN
305    e3u_ps = -1
306    e3v_ps = -1
307  ENDELSE
308;
309; read all the 2d variables
310;
311  for i = 0, n_elements(namevar)-1 do begin
312    varinq = ncdf_varinq(cdfid, namevar[i])
313    name = varinq.name
314@read_ncdf_varget
315    command = namevar[i]+'=float(temporary(res))'
316    nothing = execute(command)
317  ENDFOR
318; for compatibility with old versions of meshmask/partial steps
319; change e3[tw]p to e3[tw]_ps
320  IF n_elements(e3tp) NE 0 THEN e3t_ps = temporary(e3tp)
321  IF n_elements(e3wp) NE 0 THEN e3w_ps = temporary(e3wp)
322; in the case of key_stride ne [1, 1, 1] redefine f points
323; coordinates: they must be in the middle of 3 T points
324  if key_stride[0] NE 1 OR key_stride[1] NE 1 then BEGIN
325; we must recompute glamf and gphif...
326    IF jpi GT 1 THEN BEGIN
327      if (keyword_set(key_onearth) AND keyword_set(xnotsorted)) $
328        OR (keyword_set(key_periodic) AND key_irregular) then BEGIN
329        stepxf = (glamt + 720) MOD 360
330        stepxf = shift(stepxf, -1, -1) - stepxf
331        stepxf = [ [[stepxf]], [[stepxf + 360]], [[stepxf - 360]] ]
332        stepxf = min(abs(stepxf), dimension = 3)
333        IF NOT keyword_set(key_periodic) THEN $
334          stepxf[jpi-1, *] = stepxf[jpi-2, *]
335      ENDIF ELSE BEGIN
336        stepxf = shift(glamt, -1, -1) - glamt
337        IF keyword_set(key_periodic) THEN $
338          stepxf[jpi-1, *] = 360 + stepxf[jpi-1, *] $
339        ELSE stepxf[jpi-1, *] = stepxf[jpi-2, *]
340      ENDELSE
341      IF jpj GT 1 THEN BEGIN
342        stepxf[*, jpj-1] = stepxf[*, jpj-2]
343        stepxf[jpi-1, jpj-1] = stepxf[jpi-2, jpj-2]
344      ENDIF
345      glamf = glamt + 0.5 * stepxf
346    ENDIF ELSE glamf = glamt + 0.5
347    IF jpj GT 1 THEN BEGIN
348; we must compute stepyf: y distance between T(i,j) T(i+1,j+1)
349      stepyf = shift(gphit, -1, -1) - gphit
350      stepyf[*, jpj-1] = stepyf[*, jpj-2]
351      IF jpi GT 1 THEN BEGIN
352        if NOT keyword_set(key_periodic) THEN $
353          stepyf[jpi-1, *] = stepyf[jpi-2, *]
354        stepyf[jpi-1, jpj-1] = stepyf[jpi-2, jpj-2]
355      ENDIF
356      gphif = gphit + 0.5 * stepyf
357    ENDIF ELSE gphif = gphit + 0.5
358  ENDIF
359;-------------------------------------------------------
360; 3d arrays:
361;-------------------------------------------------------
362  nz = jpk
363  izminmesh = izminmeshsauve
364;
365  listdims = ncdf_listdims(cdfid)
366  micromask = (where(listdims EQ 'y_m'))[0]
367;
368  varinq = ncdf_varinq(cdfid, 'tmask')
369  name = varinq.name
370  IF micromask NE -1 THEN BEGIN
371; keep original values
372    iyminmeshtrue = iyminmesh
373    key_stridetrue = key_stride
374    yyy1 = firsty*key_stridetrue[1]+iyminmeshtrue
375    yyy2 = lasty*key_stridetrue[1]+iyminmeshtrue
376; the mask is stored as the bit values of the byte array (along the y
377; dimension, see micromeshmask.pro)...
378; we must modify several parameters...
379    iyminmesh = 0L
380    firsty = yyy1/8
381    lasty = yyy2/8
382    ny = lasty-firsty+1
383    key_stride = [key_stride[0], 1, key_stride[2]]
384@read_ncdf_varget
385    tmask = bytarr(jpi, jpj, jpk)
386; now we must get back the mask
387; loop on the level to save memory (the loop is short and, thus,
388; should be fast enough)
389    FOR k = 0, jpk-1 DO BEGIN
390      zzz = transpose(res[*, *, k])
391      zzz = reform(binary(zzz), 8*ny, nx, /over)
392      zzz = transpose(temporary(zzz))
393      zzz = zzz[*, yyy1 MOD 8: 8*ny - 8 + yyy2 MOD 8]
394      IF key_stridetrue[1] NE 1 THEN BEGIN
395;        IF float(strmid(!version.release,0,3)) LT 5.6 THEN BEGIN
396        nnny = (size(zzz))[2]
397        yind = key_stridetrue[1]*lindgen((nnny-1)/key_stridetrue[1]+1)
398        tmask[*, *, k] = temporary(zzz[*, yind])
399;        ENDIF ELSE tmask[*, *, k] = temporary(zzz[*, 0:*:key_stridetrue[1]])
400      ENDIF ELSE tmask[*, *, k] = temporary(zzz)
401    ENDFOR
402  ENDIF ELSE BEGIN
403@read_ncdf_varget
404    tmask = byte(res)
405  ENDELSE
406; boundary conditions used to compute umask.
407  varinq = ncdf_varinq(cdfid, 'umask')
408  name = varinq.name
409  nx = 1L
410  firstx = jpi-1
411  lastx = jpi-1
412  IF micromask NE -1 THEN BEGIN
413@read_ncdf_varget
414    umaskred = reform(binary(res), 8*ny, jpk, /over)
415    umaskred = umaskred[yyy1 MOD 8: 8*ny - 8 + yyy2 MOD 8, *]
416    IF key_stridetrue[1] NE 1 THEN umaskred = temporary(umaskred[yind, *])
417  ENDIF ELSE BEGIN
418@read_ncdf_varget
419    umaskred = reform(byte(res), /over)
420  ENDELSE
421; boundary conditions used to compute fmask (1).
422  varinq = ncdf_varinq(cdfid, 'fmask')
423  name = varinq.name
424  IF micromask NE -1 THEN BEGIN
425@read_ncdf_varget
426    fmaskredy = reform(binary(res), 8*ny, jpk, /over)
427    fmaskredy = fmaskredy[yyy1 MOD 8: 8*ny - 8 + yyy2 MOD 8, *]
428    IF key_stridetrue[1] NE 1 THEN fmaskredy = temporary(fmaskredy[yind, *])
429  ENDIF ELSE BEGIN
430@read_ncdf_varget
431    fmaskredy = reform(byte(res), /over)
432    fmaskredy = temporary(fmaskredy) MOD 2
433  ENDELSE
434; boundary conditions used to compute vmask
435  varinq = ncdf_varinq(cdfid, 'vmask')
436  name = varinq.name
437  nx = jpi
438  firstx = 0L
439  lastx = jpi-1L
440  ny = 1L
441  firsty = jpj-1
442  lasty = jpj-1
443  IF micromask NE -1 THEN BEGIN
444    yyy1 = firsty*key_stridetrue[1]+iyminmeshtrue
445    yyy2 = lasty*key_stridetrue[1]+iyminmeshtrue
446    iyminmesh = 0L
447    firsty = yyy1/8
448    lasty = yyy2/8
449    ny = lasty-firsty+1
450@read_ncdf_varget
451    IF jpk EQ 1 THEN res = reform(res, jpi, 1, jpk, /over)
452    vmaskred = transpose(temporary(res), [1, 0, 2])
453    vmaskred = reform(binary(vmaskred), 8*ny, nx, nz, /over)
454    vmaskred = transpose(temporary(vmaskred), [1, 0, 2])
455    vmaskred = reform(vmaskred[*, yyy1 MOD 8: 8*ny - 8 + yyy2 MOD 8, *])
456  ENDIF ELSE BEGIN
457@read_ncdf_varget
458    vmaskred = reform(byte(res), /over)
459  ENDELSE
460; boundary conditions used to compute fmask (2).
461  varinq = ncdf_varinq(cdfid, 'fmask')
462  name = varinq.name
463  IF micromask NE -1 THEN BEGIN
464@read_ncdf_varget
465    IF jpk EQ 1 THEN res = reform(res, jpi, 1, jpk, /over)
466    fmaskredx = transpose(temporary(res), [1, 0, 2])
467    fmaskredx = reform(binary(fmaskredx), 8*ny, nx, nz, /over)
468    fmaskredx = transpose(temporary(fmaskredx), [1, 0, 2])
469    fmaskredx = reform(fmaskredx[*, yyy1 MOD 8: 8*ny - 8 + yyy2 MOD 8, *])
470;
471    iyminmesh = iyminmeshtrue
472    key_stride = key_stridetrue
473  ENDIF ELSE BEGIN
474@read_ncdf_varget
475    fmaskredx = reform(byte(res), /over)
476    fmaskredx = fmaskredx MOD 2
477  ENDELSE
478;-------------------------------------------------------
479; 1d arrays
480;-------------------------------------------------------
481  IF (where(allvarname EQ 'e3t_O'))[0] NE -1 THEN fnamevar = ['e3t_O', 'e3w_O', 'gdept_O', 'gdepw_O'] $
482  ELSE                                            fnamevar = ['e3t', 'e3w', 'gdept', 'gdepw']
483  namevar = ['e3t', 'e3w', 'gdept', 'gdepw']
484  for i = 0, n_elements(namevar)-1 do begin
485    varinq = ncdf_varinq(cdfid, namevar[i])
486    CASE n_elements(varinq.dim) OF
487      4:BEGIN
488        command = 'ncdf_varget,cdfid,fnamevar[i],'+namevar[i] $
489                   +',offset = [0,0,izminmesh,0], count = [1,1,jpk,1]'
490        if key_stride[2] NE 1 then command = command+', stride=[1,1,key_stride[2],1]'
491      END
492      2:BEGIN
493        command = 'ncdf_varget,cdfid,fnamevar[i],'+namevar[i] $
494                   +',offset = [izminmesh,0], count = [jpk,1]'
495        if key_stride[2] NE 1 then command = command+', stride=key_stride[2]'
496      END
497      1:BEGIN
498        command = 'ncdf_varget,cdfid,fnamevar[i],'+namevar[i] $
499                   +',offset = [izminmesh], count = [jpk]'
500        if key_stride[2] NE 1 then command = command+', stride=key_stride[2]'
501      END
502    ENDCASE
503    nothing = execute(command)
504    command = namevar[i]+'=float('+namevar[i]+')'
505    nothing = execute(command)
506    command = 'if size('+namevar[i]+', /n_dimension) gt 0 then '+namevar[i]+' = reform('+namevar[i]+', /over)'
507    nothing = execute(command)
508  ENDFOR
509;-------------------------------------------------------
510  ncdf_close,  cdfid
511;-------------------------------------------------------
512; Apply Glamboundary
513;-------------------------------------------------------
514  if keyword_set(glamboundary) AND key_onearth then BEGIN
515    if glamboundary[0] NE glamboundary[1] then BEGIN
516      glamt = temporary(glamt) MOD 360
517      smaller = where(glamt LT glamboundary[0])
518      if smaller[0] NE -1 then glamt[smaller] = glamt[smaller]+360
519      bigger = where(glamt GE glamboundary[1])
520      if bigger[0] NE -1 then glamt[bigger] = glamt[bigger]-360
521      glamu = temporary(glamu) MOD 360
522      smaller = where(glamu LT glamboundary[0])
523      if smaller[0] NE -1 then glamu[smaller] = glamu[smaller]+360
524      bigger = where(glamu GE glamboundary[1])
525      if bigger[0] NE -1 then glamu[bigger] = glamu[bigger]-360
526      glamv = temporary(glamv) MOD 360
527      smaller = where(glamv LT glamboundary[0])
528      if smaller[0] NE -1 then glamv[smaller] = glamv[smaller]+360
529      bigger = where(glamv GE glamboundary[1])
530      if bigger[0] NE -1 then glamv[bigger] = glamv[bigger]-360
531      glamf = temporary(glamf) MOD 360
532      smaller = where(glamf LT glamboundary[0])
533      if smaller[0] NE -1 then glamf[smaller] = glamf[smaller]+360
534      bigger = where(glamf GE glamboundary[1])
535      if bigger[0] NE -1 then glamf[bigger] = glamf[bigger]-360
536      toosmall = where(glamu EQ glamboundary[0])
537      IF toosmall[0] NE -1 THEN glamu[toosmall] = glamu[toosmall] + 360
538      toosmall = where(glamf EQ glamboundary[0])
539      IF toosmall[0] NE -1 THEN glamf[toosmall] = glamf[toosmall] + 360
540    endif
541  endif
542;-------------------------------------------------------
543; make sure we do have 2d arrays when jpj eq 1
544;-------------------------------------------------------
545  IF jpj EQ 1 THEN BEGIN
546    glamt = reform(glamt, jpi, jpj, /over)
547    gphit = reform(gphit, jpi, jpj, /over)
548    e1t = reform(e1t, jpi, jpj, /over)
549    e2t = reform(e2t, jpi, jpj, /over)
550    glamu = reform(glamu, jpi, jpj, /over)
551    gphiu = reform(gphiu, jpi, jpj, /over)
552    e1u = reform(e1u, jpi, jpj, /over)
553    e2u = reform(e2u, jpi, jpj, /over)
554    glamv = reform(glamv, jpi, jpj, /over)
555    gphiv = reform(gphiv, jpi, jpj, /over)
556    e1v = reform(e1v, jpi, jpj, /over)
557    e2v = reform(e2v, jpi, jpj, /over)
558    glamf = reform(glamf, jpi, jpj, /over)
559    gphif = reform(gphif, jpi, jpj, /over)
560    e1f = reform(e1f, jpi, jpj, /over)
561    e2f = reform(e2f, jpi, jpj, /over)
562    IF keyword_set(key_partialstep) THEN BEGIN
563      hdept = reform(hdept, jpi, jpj, /over)
564      hdepw = reform(hdepw, jpi, jpj, /over)
565      e3t_ps = reform(e3t_ps, jpi, jpj, /over)
566      e3w_ps = reform(e3w_ps, jpi, jpj, /over)
567    ENDIF
568  ENDIF
569;-------------------------------------------------------
570  ixmindta = ixmindtasauve
571  iymindta = iymindtasauve
572  izmindta = izmindtasauve
573;-------------------------------------------------------
574  widget_control, noticebase, bad_id = nothing, /destroy
575;
576;====================================================
577; grid parameters used by xxx
578;====================================================
579;
580  IF NOT keyword_set(strcalling) THEN BEGIN
581    IF n_elements(ccmeshparameters) EQ 0 THEN strcalling = 'ncdf_meshread' $
582    ELSE strcalling = ccmeshparameters.filename
583  ENDIF
584  IF n_elements(glamt) GE 2 THEN BEGIN
585    glaminfo = moment(glamt)
586    IF finite(glaminfo[2]) EQ 0 THEN glaminfo = glaminfo[0:1]
587    gphiinfo = moment(gphit)
588    IF finite(gphiinfo[2]) EQ 0 THEN gphiinfo = gphiinfo[0:1]
589  ENDIF ELSE BEGIN
590    glaminfo = glamt
591    gphiinfo = gphit
592  ENDELSE
593  romszinfos = {h:-1, zeta:-1, theta_s:-1, theta_b:-1, hc:-1}
594  ccmeshparameters = {filename:strcalling  $
595          , glaminfo:float(string(glaminfo, format = '(E11.4)')) $
596          , gphiinfo:float(string(gphiinfo, format = '(E11.4)')) $
597          , jpiglo:jpiglo, jpjglo:jpjglo, jpkglo:jpkglo $
598          , jpi:jpi, jpj:jpj, jpk:jpk $
599          , ixminmesh:ixminmesh, ixmaxmesh:ixmaxmesh $
600          , iyminmesh:iyminmesh, iymaxmesh:iymaxmesh $
601          , izminmesh:izminmesh, izmaxmesh:izmaxmesh $
602          , key_shift:key_shift, key_periodic:key_periodic $
603          , key_stride:key_stride, key_gridtype:key_gridtype $
604          , key_yreverse:key_yreverse, key_zreverse:key_zreverse $
605          , key_partialstep:key_partialstep, key_onearth:key_onearth}
606;
607  if keyword_set(key_performance) THEN $
608    print, 'time ncdf_meshread', systime(1)-tempsun
609
610;-------------------------------------------------------
611   @updateold
612;-------------------------------------------------------
613   return
614 end
Note: See TracBrowser for help on using the repository browser.