source: trunk/SRC/Grid/ncdf_meshroms.pro @ 327

Last change on this file since 327 was 327, checked in by pinsard, 17 years ago

modification of headers : mainly blanks around = sign for keywords in declaration of function and pro

  • Property svn:keywords set to Id
File size: 19.9 KB
Line 
1;+
2;
3; @file_comments
4; read NetCDF grid file created by ROMS
5;
6; @categories
7; Grid
8;
9; @examples
10; IDL> ncdf_meshroms [,' filename']
11;
12; @param filename {in}{optional}{default='roms_grd.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 specify, 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 ONEARTH {default=1}{type=scalar: 0 or 1}
25;    Force the manual definition of data localization on the earth or not
26;       0) if the data are not on the earth
27;       1) if the data are on earth (in that case we can for example use
28;          the labels 'longitude', 'latitude' in plots).
29;    The resulting value will be stored in the common (cm_4mesh) variable key_onearth
30;    ONEARTH=0 forces PERIODIC=0, SHIFT=0 and is cancelling GLAMBOUNDARY
31;
32; @keyword GETDIMENSIONS {default=0}{type=scalar: 0 or 1}
33;    Activate this keywords if you only want to know the dimension
34;    of the domain stored in the mesh file. This dimension will be
35;    defined in jpiglo, jpjglo, jpkglo (cm_4mesh common variables)
36;
37; @keyword PERIODIC {default=computed by using the first line of glamt}{type=scalar: 0 or 1}
38;    Force the manual definition of the grid zonal periodicity.
39;    The resulting value will be stored in the common (cm_4mesh) variable key_periodic
40;    PERIODIC=0 forces SHIFT=0
41;
42; @keyword NRHO {default=1}{type=scalar}
43;    Specify the number of rho level that contain the data we want to explore.
44;    This is mainly useful when using <pro>xxx</pro> to get access to the deeper levers and vertical sections.
45;
46; @keyword SHIFT {default=computed according to glamboundary}{type=scalar}
47;    Force the manual definition of the zonal shift that must be apply to the data.
48;    The resulting value will be stored in the common (cm_4mesh) variable key_shift
49;    Note that if key_periodic=0 then in any case key_shift=0.
50;
51; @keyword STRCALLING {type=scalar string}
52;    the calling command used to call <pro>computegrid</pro> (this is used by <pro>xxx</pro>)
53;
54; @keyword STRIDE {default=[1, 1, 1]}{type=3 elements vector}
55;    Specify the stride in x, y and z direction. The resulting
56;    value will be stored in the common (cm_4mesh) variable key_stride
57;
58; @keyword _EXTRA
59; Used to pass keywords to <pro>isafile</pro>
60;
61; @uses
62; cm_4mesh
63; cm_4data
64; cm_4cal
65;
66; @restrictions
67; ixminmesh, ixmaxmesh, iyminmesh, iymaxmesh, izminmesh, izmaxmesh must
68; be defined before calling <pro>ncdf_meshread</pro>. If some of those values
69; are equal to -1 they will be automatically defined
70;
71; In the original ROMS grid, if F grid has (jpi,jpj) points then T
72; grid will have (jpi+1,jpj+1) points, U grid will have (jpi,jpj+1)
73; points and V grid will have (jpi+1,jpj) points.
74; By default C-grid used in this package needs the same number of
75; points for T,U,V and F grid, with a T point at the bottom left
76; corner of the grid. We therefore ignore the last column of T and
77; V points and the last line of T and U points.
78;
79; Scale factors are computed using the distance between the points
80; (which is not the exact definition for irregular grid).
81;
82; @history
83; Sebastien Masson (smasson\@lodyc.jussieu.fr) September 2006
84;
85; @version
86; $Id$
87;
88;-
89PRO ncdf_meshroms, filename, NRHO=nrho, GLAMBOUNDARY=glamboundary $
90                  , ONEARTH=onearth, GETDIMENSIONS=getdimensions $
91                  , PERIODIC=periodic, SHIFT=shift, STRIDE=stride $
92                  , STRCALLING=strcalling, _EXTRA=ex
93;
94  compile_opt idl2, strictarrsubs
95;
96@cm_4mesh
97@cm_4data
98@cm_4cal
99  IF NOT keyword_set(key_forgetold) THEN BEGIN
100@updatenew
101@updatekwd
102  ENDIF
103;---------------------------------------------------------
104;
105  tempsun = systime(1)          ; for key_performance
106;-------------------------------------------------------
107; find meshfile name and open it!
108;-------------------------------------------------------
109; def of filename by default
110  IF n_params() EQ 0 then filename = 'roms_grd.nc'
111  meshname = isafile(file = filename, iodirectory = iodir, _EXTRA = ex)
112  meshname = meshname[0]
113;
114  noticebase = xnotice('Reading file !C '+meshname+'!C ...')
115; if the meshmask is on tape archive ... get it back
116  IF !version.OS_FAMILY EQ 'unix' THEN spawn, '\file '+meshname+' > /dev/null'
117  cdfid = ncdf_open(meshname)
118  inq = ncdf_inquire(cdfid)
119;------------------------------------------------------------
120; dimensions
121;------------------------------------------------------------
122  ncdf_diminq, cdfid, 'xi_rho', name, jpiglo
123  ncdf_diminq, cdfid, 'eta_rho', name, jpjglo
124  IF n_elements(nrho) NE 0 THEN jpkglo = long(nrho[0]) $
125    ELSE jpkglo = 1L
126;
127  if keyword_set(getdimensions) then begin
128    widget_control, noticebase, bad_id = nothing, /destroy
129    ncdf_close,  cdfid
130    return
131  endif
132;-------------------------------------------------------
133; check that all i[xyz]min[ax]mesh are well defined
134;-------------------------------------------------------
135  if n_elements(ixminmesh) EQ 0 THEN ixminmesh = 0
136  if n_elements(ixmaxmesh) EQ 0 then ixmaxmesh = jpiglo-1
137  if ixminmesh EQ -1 THEN ixminmesh = 0
138  IF ixmaxmesh EQ -1 then ixmaxmesh = jpiglo-1
139  if n_elements(iyminmesh) EQ 0 THEN iyminmesh = 0
140  IF n_elements(iymaxmesh) EQ 0 then iymaxmesh = jpjglo-1
141  if iyminmesh EQ -1 THEN iyminmesh = 0
142  IF iymaxmesh EQ -1 then iymaxmesh = jpjglo-1
143  if n_elements(izminmesh) EQ 0 THEN izminmesh = 0
144  IF n_elements(izmaxmesh) EQ 0 then izmaxmesh = jpkglo-1
145  if izminmesh EQ -1 THEN izminmesh = 0
146  IF izmaxmesh EQ -1 then izmaxmesh = jpkglo-1
147; definition of jpi,jpj,jpj
148  jpi = long(ixmaxmesh-ixminmesh+1)
149  jpj = long(iymaxmesh-iyminmesh+1)
150  jpk = long(izmaxmesh-izminmesh+1)
151;-------------------------------------------------------
152; check onearth and its consequences
153;-------------------------------------------------------
154  IF n_elements(onearth) EQ 0 THEN key_onearth = 1 $
155  ELSE key_onearth = keyword_set(onearth)
156  IF NOT key_onearth THEN BEGIN
157    periodic = 0
158    shift = 0
159  ENDIF
160;-------------------------------------------------------
161; automatic definition of key_periodic
162;-------------------------------------------------------
163  IF n_elements(periodic) EQ 0 THEN BEGIN
164    IF jpi GT 1 THEN BEGIN
165      ncdf_varget, cdfid, 'lon_rho', xaxis $
166                   , offset = [ixminmesh, iyminmesh], count = [jpi, 1]
167      xaxis = (xaxis+720) MOD 360
168      xaxis = xaxis[sort(xaxis)]
169      key_periodic = (xaxis[jpi-1]+2*(xaxis[jpi-1]-xaxis[jpi-2])) $
170                     GE (xaxis[0]+360)
171    ENDIF ELSE key_periodic = 0
172  ENDIF ELSE key_periodic = keyword_set(periodic)
173;-------------------------------------------------------
174; automatic definition of key_shift
175;-------------------------------------------------------
176  IF n_elements(shift) EQ 0 THEN BEGIN
177    key_shift = long(testvar(var = key_shift))
178;  key_shift will be defined according to the first line of glamt.
179    if keyword_set(glamboundary) AND jpi GT 1 AND key_periodic EQ 1 $
180    THEN BEGIN
181      ncdf_varget, cdfid, 'lon_rho', xaxis $
182                   , offset = [ixminmesh, iyminmesh], count = [jpi, 1]
183; xaxis between glamboundary[0] and glamboundary[1]
184      xaxis = xaxis MOD 360
185      smaller = where(xaxis LT glamboundary[0])
186      if smaller[0] NE -1 then xaxis[smaller] = xaxis[smaller]+360
187      bigger = where(xaxis GE glamboundary[1])
188      if bigger[0] NE -1 then xaxis[bigger] = xaxis[bigger]-360
189;
190      key_shift = (where(xaxis EQ min(xaxis)))[0]
191      IF key_shift NE 0 THEN BEGIN
192        key_shift = jpi-key_shift
193        xaxis = shift(xaxis, key_shift)
194      ENDIF
195;
196      IF array_equal(sort(xaxis), lindgen(jpi)) NE 1 THEN BEGIN
197        ras = report (['the x axis (1st line of glamt) is not sorted in the increasing order after the automatic definition of key_shift', $
198        '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...'])
199        widget_control, noticebase, bad_id = nothing, /destroy
200        return
201      ENDIF
202;
203    ENDIF ELSE key_shift = 0
204  ENDIF ELSE key_shift = long(shift)*(key_periodic EQ 1)
205;-------------------------------------------------------
206; check key_stride and related things
207;-------------------------------------------------------
208  if n_elements(stride) eq 3 then key_stride = stride
209  if n_elements(key_stride) LE 2 then key_stride = [1, 1, 1]
210  key_stride = 1l > long(key_stride)
211  IF total(key_stride) NE 3  THEN BEGIN
212    IF key_shift NE 0 THEN BEGIN
213; for explanation, see header of read_ncdf_varget.pro
214      jpiright = key_shift
215      jpileft = jpi - key_shift - ( (key_stride[0]-1)-((key_shift-1) MOD key_stride[0]) )
216      jpi = ((jpiright-1)/key_stride[0]+1) + ((jpileft-1)/key_stride[0]+1)
217    ENDIF ELSE jpi = (jpi-1)/key_stride[0]+1
218    jpj = (jpj-1)/key_stride[1]+1
219    jpk = (jpk-1)/key_stride[2]+1
220  ENDIF
221;-------------------------------------------------------
222; for the variables related to the partial steps
223;-------------------------------------------------------
224  key_partialstep = 0
225  hdept = -1
226  hdepw = -1
227;-------------------------------------------------------
228; default definitions to be able to use read_ncdf_varget
229;-------------------------------------------------------
230; default definitions to be able to use read_ncdf_varget
231  ixmindtasauve = testvar(var = ixmindta)
232  iymindtasauve = testvar(var = iymindta)
233  izmindtasauve = testvar(var = izmindta)
234;
235  ixmindta = 0l
236  iymindta = 0l
237  izmindta = 0l
238;
239  jpt = 1
240  time = 1
241  firsttps = 0
242;
243  firstx = 0
244  lastx = jpi-1
245  firsty = 0
246  lasty = jpj-1
247  firstz = 0
248  lastz = jpk-1
249  nx = jpi
250  ny = jpj
251  nz = 1
252  izminmeshsauve = izminmesh
253  izminmesh = 0
254;
255  key_yreverse = 0
256  key_zreverse = 1
257  key_gridtype = 'c'
258;-------------------------------------------------------
259; 2d arrays:
260;-------------------------------------------------------
261; list the 2d variables that must be read
262  namebase  = ['lon_', 'lat_', 'mask_', 'x_', 'y_']
263  namebase2 = ['glam', 'gphi', 'mask' , 'd1', 'd2']
264;
265; read all grid T variables
266;
267  for i = 0, n_elements(namebase)-1 do begin
268    varinq = ncdf_varinq(cdfid, namebase[i]+'rho')
269    name = varinq.name
270@read_ncdf_varget
271    command = namebase2[i]+'t = float(temporary(res))'
272    nothing = execute(command)
273  ENDFOR
274  d1t = 1.e3*(shift(d1t, -1, 0) - d1t)
275  d2t = 1.e3*(shift(d2t, 0, -1) - d2t)
276  for i = 0, n_elements(namebase2)-1 do begin
277    command = namebase2[i]+'t = '+namebase2[i]+'t[0:jpi-2, 0:jpj-2]'
278    nothing = execute(command)
279  ENDFOR
280  tmask = byte(temporary(maskt))
281  IF jpk GT 1 THEN tmask = reform(tmask[*]#replicate(1b, jpk), jpi-1, jpj-1, jpk, /overwrite)
282  e1u = temporary(d1t)
283  e2v = temporary(d2t)
284; h: Final bathymetry at RHO-points
285    varinq = ncdf_varinq(cdfid, 'h')
286    name = varinq.name
287@read_ncdf_varget
288    hroms = float(temporary(res))
289    hroms = hroms[0:jpi-2, 0:jpj-2]
290;
291; read all grid U variables
292;
293  jpiglo = jpiglo - 1
294  jpi = jpi - 1
295  ixmaxmesh = ixmaxmesh - 1
296  firstx = 0
297  lastx = jpi-1
298  nx = jpi
299;
300  for i = 0, n_elements(namebase)-1 do begin
301    varinq = ncdf_varinq(cdfid, namebase[i]+'u')
302    name = varinq.name
303@read_ncdf_varget
304    command = namebase2[i]+'u = float(temporary(res))'
305    nothing = execute(command)
306  ENDFOR
307  tmpsave = 2. * 1.e3 * d1u[0, 0:jpj-2]
308  d1u = 1.e3*(shift(d1u, -1, 0) - d1u)
309  d2u = 1.e3*(shift(d2u, 0, -1) - d2u)
310  for i = 0, n_elements(namebase2)-1 do begin
311    command = namebase2[i]+'u = '+namebase2[i]+'u[*, 0:jpj-2]'
312    nothing = execute(command)
313  ENDFOR
314  umaskred = byte((temporary(masku))[jpi-1, *])
315  IF jpk GT 1 THEN umaskred = reform(umaskred[*]#replicate(1b, jpk), 1, jpj-1, jpk, /overwrite)
316  e1t = temporary(d1u)
317  e1t = shift(temporary(e1t), 1, 0)
318  e1t[0, *] = temporary(tmpsave)
319  e2f = temporary(d2u)
320;
321; read all grid V variables
322;
323  jpiglo = jpiglo + 1
324  jpi = jpi + 1
325  ixmaxmesh = ixmaxmesh + 1
326  firstx = 0
327  lastx = jpi-1
328  nx = jpi
329  jpjglo = jpjglo - 1
330  jpj = jpj - 1
331  iymaxmesh = iymaxmesh - 1
332  firsty = 0
333  lasty = jpj-1
334  ny = jpj
335;
336  for i = 0, n_elements(namebase)-1 do begin
337    varinq = ncdf_varinq(cdfid, namebase[i]+'v')
338    name = varinq.name
339@read_ncdf_varget
340    command = namebase2[i]+'v = float(temporary(res))'
341    nothing = execute(command)
342  ENDFOR
343  d1v = 1.e3*(shift(d1v, -1, 0) - d1v)
344  tmpsave = 2. * 1.e3 * d2v[0:jpi-2, 0]
345  d2v = 1.e3*(shift(d2v, 0, -1) - d2v)
346  for i = 0, n_elements(namebase2)-1 do begin
347    command = namebase2[i]+'v = '+namebase2[i]+'v[0:jpi-2, *]'
348    nothing = execute(command)
349  ENDFOR
350  vmaskred = byte((temporary(maskv))[*, jpj-1])
351  IF jpk GT 1 THEN vmaskred = reform(vmaskred[*]#replicate(1b, jpk), jpi-1, 1, jpk, /overwrite)
352  e1f = temporary(d1v)
353  e2t = temporary(d2v)
354  e2t = shift(temporary(e2t), 0, 1)
355  e2t[*, 0] = temporary(tmpsave)
356;
357; read all grid F variables
358;
359  jpiglo = jpiglo - 1
360  jpi = jpi - 1
361  ixmaxmesh = ixmaxmesh - 1
362  firstx = 0
363  lastx = jpi-1
364  nx = jpi
365;
366  for i = 0, n_elements(namebase)-1 do begin
367    varinq = ncdf_varinq(cdfid, namebase[i]+'psi')
368    name = varinq.name
369@read_ncdf_varget
370    command = namebase2[i]+'f = float(temporary(res))'
371    nothing = execute(command)
372  ENDFOR
373  tmpsave1 = 2. * 1.e3 * d1f[0, *]
374  d1f = 1.e3*(shift(d1f, -1, 0) - d1f)
375  tmpsave2 = 2. * 1.e3 * d2f[*, 0]
376  d2f = 1.e3*(shift(d2f, 0, -1) - d2f)
377  fmaskredy = byte(maskf[jpi-1, *])
378  IF jpk GT 1 THEN fmaskredy = reform(fmaskredy[*]#replicate(1b, jpk), 1, jpj, jpk, /overwrite)
379  fmaskredx = byte((temporary(maskf))[*, jpj-1])
380  IF jpk GT 1 THEN fmaskredx = reform(fmaskredx[*]#replicate(1b, jpk), jpi, 1, jpk, /overwrite)
381  e1v = temporary(d1f)
382  e1v = shift(temporary(e1v), 1, 0)
383  e1v[0, *] = temporary(tmpsave1)
384  e2u = temporary(d2f)
385  e2u = shift(temporary(e2u), 0, 1)
386  e2u[*, 0] = temporary(tmpsave2)
387;-------------------------------------------------------
388; in the case of key_stride ne [1, 1, 1] redefine f points
389; coordinates: they must be in the middle of 3 T points
390;-------------------------------------------------------
391  if key_stride[0] NE 1 OR key_stride[1] NE 1 then BEGIN
392; we must recompute glamf and gphif...
393    IF jpi GT 1 THEN BEGIN
394      if (keyword_set(key_onearth) AND keyword_set(xnotsorted)) $
395        OR (keyword_set(key_periodic) AND key_irregular) then BEGIN
396        stepxf = (glamt + 720) MOD 360
397        stepxf = shift(stepxf, -1, -1) - stepxf
398        stepxf = [ [[stepxf]], [[stepxf + 360]], [[stepxf - 360]] ]
399        stepxf = min(abs(stepxf), dimension = 3)
400        IF NOT keyword_set(key_periodic) THEN $
401          stepxf[jpi-1, *] = stepxf[jpi-2, *]
402      ENDIF ELSE BEGIN
403        stepxf = shift(glamt, -1, -1) - glamt
404        IF keyword_set(key_periodic) THEN $
405          stepxf[jpi-1, *] = 360 + stepxf[jpi-1, *] $
406        ELSE stepxf[jpi-1, *] = stepxf[jpi-2, *]
407      ENDELSE
408      IF jpj GT 1 THEN BEGIN
409        stepxf[*, jpj-1] = stepxf[*, jpj-2]
410        stepxf[jpi-1, jpj-1] = stepxf[jpi-2, jpj-2]
411      ENDIF
412      glamf = glamt + 0.5 * stepxf
413    ENDIF ELSE glamf = glamt + 0.5
414    IF jpj GT 1 THEN BEGIN
415; we must compute stepyf: y distance between T(i,j) T(i+1,j+1)
416      stepyf = shift(gphit, -1, -1) - gphit
417      stepyf[*, jpj-1] = stepyf[*, jpj-2]
418      IF jpi GT 1 THEN BEGIN
419        if NOT keyword_set(key_periodic) THEN $
420          stepyf[jpi-1, *] = stepyf[jpi-2, *]
421        stepyf[jpi-1, jpj-1] = stepyf[jpi-2, jpj-2]
422      ENDIF
423      gphif = gphit + 0.5 * stepyf
424    ENDIF ELSE gphif = gphit + 0.5
425  ENDIF
426;-------------------------------------------------------
427; 1d arrays
428;-------------------------------------------------------
429  gdept = findgen(jpk)
430  gdepw = findgen(jpk)
431  e3t = replicate(1., jpk)
432  e3w = replicate(1., jpk)
433;-------------------------------------------------------
434  ncdf_close, cdfid
435;-------------------------------------------------------
436; Apply Glamboudary
437;-------------------------------------------------------
438  if keyword_set(glamboundary) AND key_onearth then BEGIN
439    if glamboundary[0] NE glamboundary[1] then BEGIN
440      glamt = temporary(glamt) MOD 360
441      smaller = where(glamt LT glamboundary[0])
442      if smaller[0] NE -1 then glamt[smaller] = glamt[smaller]+360
443      bigger = where(glamt GE glamboundary[1])
444      if bigger[0] NE -1 then glamt[bigger] = glamt[bigger]-360
445      glamu = temporary(glamu) MOD 360
446      smaller = where(glamu LT glamboundary[0])
447      if smaller[0] NE -1 then glamu[smaller] = glamu[smaller]+360
448      bigger = where(glamu GE glamboundary[1])
449      if bigger[0] NE -1 then glamu[bigger] = glamu[bigger]-360
450      glamv = temporary(glamv) MOD 360
451      smaller = where(glamv LT glamboundary[0])
452      if smaller[0] NE -1 then glamv[smaller] = glamv[smaller]+360
453      bigger = where(glamv GE glamboundary[1])
454      if bigger[0] NE -1 then glamv[bigger] = glamv[bigger]-360
455      glamf = temporary(glamf) MOD 360
456      smaller = where(glamf LT glamboundary[0])
457      if smaller[0] NE -1 then glamf[smaller] = glamf[smaller]+360
458      bigger = where(glamf GE glamboundary[1])
459      if bigger[0] NE -1 then glamf[bigger] = glamf[bigger]-360
460      toosmall = where(glamu EQ glamboundary[0])
461      IF toosmall[0] NE -1 THEN glamu[toosmall] = glamu[toosmall] + 360
462      toosmall = where(glamf EQ glamboundary[0])
463      IF toosmall[0] NE -1 THEN glamf[toosmall] = glamf[toosmall] + 360
464    endif
465  endif
466;-------------------------------------------------------
467; make sure we do have 2d arrays when jpj eq 1
468;-------------------------------------------------------
469  IF jpj EQ 1 THEN BEGIN
470    glamt = reform(glamt, jpi, jpj, /over)
471    gphit = reform(gphit, jpi, jpj, /over)
472    e1t = reform(e1t, jpi, jpj, /over)
473    e2t = reform(e2t, jpi, jpj, /over)
474    glamu = reform(glamu, jpi, jpj, /over)
475    gphiu = reform(gphiu, jpi, jpj, /over)
476    e1u = reform(e1u, jpi, jpj, /over)
477    e2u = reform(e2u, jpi, jpj, /over)
478    glamv = reform(glamv, jpi, jpj, /over)
479    gphiv = reform(gphiv, jpi, jpj, /over)
480    e1v = reform(e1v, jpi, jpj, /over)
481    e2v = reform(e2v, jpi, jpj, /over)
482    glamf = reform(glamf, jpi, jpj, /over)
483    gphif = reform(gphif, jpi, jpj, /over)
484    e1f = reform(e1f, jpi, jpj, /over)
485    e2f = reform(e2f, jpi, jpj, /over)
486  ENDIF
487;-------------------------------------------------------
488  ixmindta = ixmindtasauve
489  iymindta = iymindtasauve
490  izmindta = izmindtasauve
491;-------------------------------------------------------
492  widget_control, noticebase, bad_id = nothing, /destroy
493;
494;====================================================
495; grid parameters used by xxx
496;====================================================
497;
498  IF NOT keyword_set(strcalling) THEN BEGIN
499    IF n_elements(ccmeshparameters) EQ 0 THEN strcalling = 'ncdf_meshroms' $
500    ELSE strcalling = ccmeshparameters.filename
501  ENDIF
502  IF n_elements(glamt) GE 2 THEN BEGIN
503    glaminfo = moment(glamt)
504    IF finite(glaminfo[2]) EQ 0 THEN glaminfo = glaminfo[0:1]
505    gphiinfo = moment(gphit)
506    IF finite(gphiinfo[2]) EQ 0 THEN gphiinfo = gphiinfo[0:1]
507  ENDIF ELSE BEGIN
508    glaminfo = glamt
509    gphiinfo = gphit
510  ENDELSE
511  romszinfos = {h:hroms, zeta:replicate(0., jpi, jpj), theta_s:-1, theta_b:-1, hc:-1}
512  ccmeshparameters = {filename:strcalling  $
513          , glaminfo:float(string(glaminfo, format = '(E11.4)')) $
514          , gphiinfo:float(string(gphiinfo, format = '(E11.4)')) $
515          , jpiglo:jpiglo, jpjglo:jpjglo, jpkglo:jpkglo $
516          , jpi:jpi, jpj:jpj, jpk:jpk $
517          , ixminmesh:ixminmesh, ixmaxmesh:ixmaxmesh $
518          , iyminmesh:iyminmesh, iymaxmesh:iymaxmesh $
519          , izminmesh:izminmesh, izmaxmesh:izmaxmesh $
520          , key_shift:key_shift, key_periodic:key_periodic $
521          , key_stride:key_stride, key_gridtype:key_gridtype $
522          , key_yreverse:key_yreverse, key_zreverse:key_zreverse $
523          , key_partialstep:key_partialstep, key_onearth:key_onearth}
524;
525  if keyword_set(key_performance) THEN $
526    print, 'time ncdf_meshread', systime(1)-tempsun
527
528;-------------------------------------------------------
529   @updateold
530;-------------------------------------------------------
531   return
532 end
Note: See TracBrowser for help on using the repository browser.