source: XIOS/trunk/src/io/inetcdf4.cpp @ 686

Last change on this file since 686 was 686, checked in by rlacroix, 9 years ago

Use the NetCDF wrapper for inputs for better error checking.

  • Property copyright set to
    Software name : XIOS (Xml I/O Server)
    http://forge.ipsl.jussieu.fr/ioserver
    Creation date : January 2009
    Licence : CeCCIL version2
    see license file in root directory : Licence_CeCILL_V2-en.txt
    or http://www.cecill.info/licences/Licence_CeCILL_V2-en.html
    Holder : CEA/LSCE (Laboratoire des Sciences du CLimat et de l'Environnement)
    CNRS/IPSL (Institut Pierre Simon Laplace)
    Project Manager : Yann Meurdesoif
    yann.meurdesoif@cea.fr
File size: 24.8 KB
Line 
1#include "inetcdf4.hpp"
2#include "netCdfInterface.hpp"
3
4#include <boost/algorithm/string.hpp>
5
6namespace xios
7{
8  CINetCDF4::CINetCDF4(const StdString& filename, const MPI_Comm* comm /*= NULL*/, bool multifile /*= true*/)
9    : mpi(comm != NULL && !multifile)
10  {
11    if (comm != NULL)
12    {
13      if (!multifile) CNetCdfInterface::openPar(filename, NC_NOWRITE, *comm, MPI_INFO_NULL, this->ncidp);
14      else CNetCdfInterface::open(filename, NC_NOWRITE, this->ncidp);
15    }
16      else CNetCdfInterface::open(filename, NC_NOWRITE, this->ncidp);
17  }
18
19  CINetCDF4::~CINetCDF4(void)
20  { /* Nothing to do */ }
21
22  //---------------------------------------------------------------
23
24  void CINetCDF4::close(void)
25  {
26    CNetCdfInterface::close(this->ncidp);
27  }
28
29  //---------------------------------------------------------------
30
31  int CINetCDF4::getGroup(const CVarPath* const path)
32  {
33    int retvalue = this->ncidp;
34    if (path == NULL) return retvalue;
35    CVarPath::const_iterator it = path->begin(), end = path->end();
36
37    for (; it != end; it++)
38    {
39      const StdString& groupid = *it;
40      CNetCdfInterface::inqNcId(retvalue, groupid, retvalue);
41    }
42
43    return retvalue;
44  }
45
46  int CINetCDF4::getVariable(const StdString& varname,
47                             const CVarPath* const path)
48  {
49    int varid = 0;
50    int grpid = this->getGroup(path);
51    CNetCdfInterface::inqVarId(grpid, varname, varid);
52    return varid;
53  }
54
55  int CINetCDF4::getDimension(const StdString& dimname,
56                              const CVarPath* const path)
57  {
58    int dimid = 0;
59    int grpid = this->getGroup(path);
60    CNetCdfInterface::inqDimId(grpid, dimname, dimid);
61    return dimid;
62  }
63
64  std::pair<nc_type, StdSize> CINetCDF4::getAttribute(const StdString& attname,
65                                                      const StdString* const var,
66                                                      const CVarPath* const path)
67  {
68    std::pair<nc_type, StdSize> retvalue;
69    int grpid = this->getGroup(path);
70    int varid = (var != NULL) ? this->getVariable(*var, path) : NC_GLOBAL;
71    CNetCdfInterface::inqAtt(grpid, varid, attname, retvalue.first, retvalue.second);
72    return retvalue;
73  }
74
75  int CINetCDF4::getUnlimitedDimension(const CVarPath* const path)
76  {
77    int dimid = 0;
78    int grpid = this->getGroup(path);
79    CNetCdfInterface::inqUnLimDim(grpid, dimid);
80    return dimid;
81  }
82
83  StdString CINetCDF4::getUnlimitedDimensionName(const CVarPath* const path)
84  {
85    int grpid = this->getGroup(path);
86    int dimid = this->getUnlimitedDimension(path);
87
88    StdString dimname;
89    if (dimid != -1)
90      CNetCdfInterface::inqDimName(grpid, dimid, dimname);
91    return dimname;
92  }
93
94  //---------------------------------------------------------------
95
96  StdSize CINetCDF4::getNbVertex(const StdString& name,
97                                 const CVarPath* const path)
98  {
99
100    if (this->isRectilinear(name, path) ||
101       this->isCurvilinear(name, path))
102    {
103      if (this->is3Dim(name, path)) return 8;
104      else return 4;
105    }
106    if (this->isUnstructured(name, path))
107    {
108      StdString bound = this->getBoundsId
109            (this->getCoordinatesIdList(name, path).back(), path);
110      StdString dim = this->getDimensionsList(&bound, path).back();
111      return this->getDimensions(&bound, path)[dim];
112    }
113    return size_t(-1);
114  }
115
116  //---------------------------------------------------------------
117
118  std::list<StdString> CINetCDF4::getGroups(const CVarPath* const path)
119  {
120    int nbgroup = 0, *groupid = NULL;
121    int grpid = this->getGroup(path);
122    std::list<StdString> retvalue;
123
124    CNetCdfInterface::inqGrpIds(grpid, nbgroup, NULL);
125    groupid = new int[nbgroup]();
126    CNetCdfInterface::inqGrpIds(grpid, nbgroup, groupid);
127
128    for (int i = 0; i < nbgroup; i++)
129    {
130      StdString fullGrpName;
131      CNetCdfInterface::inqGrpFullName(groupid[i], fullGrpName);
132      retvalue.push_back(fullGrpName);
133    }
134
135    delete [] groupid;
136    return retvalue;
137  }
138
139  std::list<StdString> CINetCDF4::getVariables(const CVarPath* const path)
140  {
141    int nbvar = 0, *varid = NULL;
142    int grpid = this->getGroup(path);
143    std::list<StdString> retvalue;
144
145    CNetCdfInterface::inqVarIds(grpid, nbvar, NULL);
146    varid = new int[nbvar]();
147    CNetCdfInterface::inqVarIds(grpid, nbvar, varid);
148
149    for (int i = 0; i < nbvar; i++)
150    {
151      StdString varName;
152      CNetCdfInterface::inqVarName(grpid, varid[i], varName);
153      retvalue.push_back(varName);
154    }
155
156    delete [] varid;
157    return retvalue;
158  }
159
160  StdSize CINetCDF4::getNbOfTimestep(const CVarPath* const path)
161  {
162    return this->getDimensions(NULL, path)[this->getUnlimitedDimensionName(path)];
163  }
164
165  std::set<StdString> CINetCDF4::getBoundVariables(const CVarPath* const path)
166  {
167    std::set<StdString> retvalue;
168    std::list<StdString> variables = this->getVariables(path);
169    std::list<StdString>::const_iterator it = variables.begin(), end = variables.end();
170    for (; it != end; it++)
171    {
172      const StdString& var = *it;
173      if (this->hasBounds(var, path))
174        retvalue.insert(retvalue.end(), this->getBoundsId(var, path));
175    }
176    return retvalue;
177  }
178
179  std::set<StdString> CINetCDF4::getCoordVariables(const CVarPath* const path)
180  {
181    std::set<StdString> retvalue;
182    std::list<StdString> variables = this->getVariables(path);
183    std::list<StdString>::const_iterator it = variables.begin(), end = variables.end();
184    for (; it != end; it++)
185    {
186      const StdString& var = *it;
187      std::list<StdString> coords = this->getCoordinatesIdList(var, path);
188      std::list<StdString>::const_iterator it = coords.begin(), end = coords.end();
189      for (; it != end; it++)
190      {
191        const StdString& coord = *it;
192        if (this->hasVariable(coord, path))
193          retvalue.insert(retvalue.end(), coord);
194      }
195    }
196    return retvalue;
197  }
198
199  std::list<StdString> CINetCDF4::getDimensionsList(const StdString* const var, const CVarPath* const path)
200  {
201    int nbdim = 0, *dimid = NULL;
202    int grpid = this->getGroup(path);
203    int varid = (var != NULL) ? this->getVariable(*var, path) : NC_GLOBAL;
204    std::list<StdString> retvalue;
205
206    if (var != NULL)
207    {
208      CNetCdfInterface::inqVarNDims(grpid, varid, nbdim);
209      dimid = new int[nbdim]();
210      CNetCdfInterface::inqVarDimId(grpid, varid, dimid);
211    }
212    else
213    {
214      CNetCdfInterface::inqDimIds(grpid, nbdim, NULL, 1);
215      dimid = new int[nbdim]();
216      CNetCdfInterface::inqDimIds(grpid, nbdim, dimid, 1);
217    }
218
219    for (int i = 0; i < nbdim; i++)
220    {
221      std::string dimname;
222      CNetCdfInterface::inqDimName(grpid, dimid[i], dimname);
223      retvalue.push_back(dimname);
224    }
225    delete [] dimid;
226
227    return retvalue;
228  }
229
230  std::map<StdString, StdSize> CINetCDF4::getDimensions(const StdString* const var, const CVarPath* const path)
231  {
232    int nbdim = 0, *dimid = NULL;
233    int grpid = this->getGroup(path);
234    int varid = (var != NULL) ? this->getVariable(*var, path) : NC_GLOBAL;
235    std::map<StdString, StdSize> retvalue;
236
237    if (var != NULL)
238    {
239      CNetCdfInterface::inqVarNDims(grpid, varid, nbdim);
240      dimid = new int[nbdim]();
241      CNetCdfInterface::inqVarDimId(grpid, varid, dimid);
242    }
243    else
244    {
245      CNetCdfInterface::inqDimIds(grpid, nbdim, NULL, 1);
246      dimid = new int[nbdim]();
247      CNetCdfInterface::inqDimIds(grpid, nbdim, dimid, 1);
248    }
249
250    for (int i = 0; i < nbdim; i++)
251    {
252      std::string dimname;
253      CNetCdfInterface::inqDimName(grpid, dimid[i], dimname);
254      StdSize size = 0;
255      CNetCdfInterface::inqDimLen(grpid, dimid[i], size);
256
257      retvalue.insert(retvalue.end(), std::make_pair(dimname, size));
258    }
259    delete [] dimid;
260
261    return retvalue;
262  }
263
264  std::list<StdString> CINetCDF4::getAttributes(const StdString* const var, const CVarPath* const path)
265  {
266    int nbatt = 0;
267    std::list<StdString> retvalue;
268    int grpid = this->getGroup(path);
269    int varid = (var != NULL) ? this->getVariable(*var, path) : NC_GLOBAL;
270
271    if (var != NULL)
272      CNetCdfInterface::inqVarNAtts(grpid, varid, nbatt);
273    else
274      CNetCdfInterface::inqNAtts(grpid, nbatt);
275
276    for (int i = 0; i < nbatt; i++)
277    {
278      StdString attname;
279      CNetCdfInterface::inqAttName(grpid, varid, i, attname);
280      retvalue.push_back(attname);
281    }
282    return retvalue;
283  }
284
285  int CINetCDF4::getAttributeId(const StdString& name,
286                                const StdString* const var,
287                                const CVarPath* const path)
288  {
289    int retvalue = 0;
290    std::list<StdString> atts = this->getAttributes(var, path);
291    std::list<StdString>::const_iterator it = atts.begin(), end = atts.end();
292    for (; it != end; it++)
293    {
294      const StdString& attname = *it;
295      if (attname.compare(name) == 0)
296        return retvalue;
297      retvalue++;
298    }
299    return -1;
300  }
301
302  //---------------------------------------------------------------
303
304  bool CINetCDF4::hasMissingValue(const StdString& name,
305                                  const CVarPath* const path)
306  {
307    return (this->hasAttribute("missing_value", &name, path) || this->hasAttribute("_FillValue", &name, path));
308  }
309
310  bool CINetCDF4::hasAttribute(const StdString& name,
311                               const StdString* const var ,
312                               const CVarPath* const path)
313  {
314    std::list<StdString> atts = this->getAttributes(var, path);
315    std::list<StdString>::const_iterator it = atts.begin(), end = atts.end();
316    for (; it != end; it++)
317    {
318      const StdString& attname = *it;
319      if (attname.compare(name) == 0) return true;
320    }
321    return false;
322  }
323
324  bool CINetCDF4::hasVariable(const StdString& name,
325                              const CVarPath* const path)
326  {
327    std::list<StdString> variables = this->getVariables(path);
328    std::list<StdString>::const_iterator it = variables.begin(), end = variables.end();
329    for (; it != end; it++)
330    {
331      const StdString& varname = *it;
332      if (varname.compare(name) == 0) return true;
333    }
334    return false;
335  }
336
337  bool CINetCDF4::hasCoordinates(const StdString& name,
338                                 const CVarPath* const path)
339  {
340    return this->hasAttribute("coordinates", &name, path);
341  }
342
343  bool CINetCDF4::hasBounds(const StdString& name,
344                            const CVarPath* const path)
345  {
346    return this->hasAttribute("bounds", &name, path);
347  }
348
349  bool CINetCDF4::hasTemporalDim(const CVarPath* const path)
350  {
351    return (this->getUnlimitedDimension(path) != -1);
352  }
353
354  //---------------------------------------------------------------
355
356  template <class T>
357  std::vector<T> CINetCDF4::getAttributeValue(const StdString& name,
358                                              const StdString* const var,
359                                              const CVarPath* const path)
360  {
361    int grpid = this->getGroup(path);
362    int varid = (var != NULL) ? this->getVariable(*var, path) : NC_GLOBAL;
363    std::pair<nc_type , StdSize> attinfos = this->getAttribute(name, var, path);
364    std::vector<T> retvalue(attinfos.second);
365    nc_type type = CNetCdfInterface::getNcType<T>();
366    if (attinfos.first != type)
367      ERROR("CINetCDF4::getAttributeValue<T>(name, var, path)",
368            << "[ name : " << name
369            << ", type requested :" << attinfos.first
370            << ", type stored : " << type << "]"
371            << " Invalid type !");
372    CNetCdfInterface::getAttType(grpid, varid, name.c_str(), &retvalue[0]);
373    return retvalue;
374  }
375
376  template std::vector<double> CINetCDF4::getAttributeValue(const StdString& name,
377                                                            const StdString* const var,
378                                                            const CVarPath* const path);
379  template std::vector<float> CINetCDF4::getAttributeValue(const StdString& name,
380                                                           const StdString* const var,
381                                                           const CVarPath* const path);
382  template std::vector<int> CINetCDF4::getAttributeValue(const StdString& name,
383                                                         const StdString* const var,
384                                                         const CVarPath* const path);
385  template std::vector<char> CINetCDF4::getAttributeValue(const StdString& name,
386                                                          const StdString* const var,
387                                                          const CVarPath* const path);
388
389  StdString CINetCDF4::getAttributeValue(const StdString& name,
390                                         const StdString* const var,
391                                         const CVarPath* const path)
392  {
393    std::vector<char> data = this->getAttributeValue<char>(name, var, path);
394
395    return StdString(data.begin(), data.end());
396  }
397
398  template <class T>
399  T CINetCDF4::getMissingValue(const StdString& name, const CVarPath* const path)
400  {
401    if (this->hasAttribute("missing_value", &name, path))
402      return this->getAttributeValue<T>("missing_value", &name, path)[0];
403    if (this->hasAttribute("_FillValue", &name, path))
404      return this->getAttributeValue<T>("_FillValue", &name, path)[0];
405    return 0;
406  }
407
408  template double CINetCDF4::getMissingValue(const StdString& name, const CVarPath* const path);
409  template float CINetCDF4::getMissingValue(const StdString& name, const CVarPath* const path);
410  template int CINetCDF4::getMissingValue(const StdString& name, const CVarPath* const path);
411  template char CINetCDF4::getMissingValue(const StdString& name, const CVarPath* const path);
412
413  //---------------------------------------------------------------
414
415  std::list<StdString> CINetCDF4::getCoordinatesIdList(const StdString& name, const CVarPath* const path)
416  {
417    std::list<StdString> retvalue;
418    StdString value = this->getCoordinatesId(name, path);
419
420    boost::split(retvalue, value, boost::is_any_of(" "));
421
422    std::list<StdString>::iterator it = retvalue.begin(), end = retvalue.end();
423    for (; it != end; it++)
424    {
425      StdString& coord = *it;
426      coord.assign(coord.data());
427    }
428    return retvalue;
429  }
430
431  StdString CINetCDF4::getCoordinatesId(const StdString& name, const CVarPath* const path)
432  {
433    StdString retvalue;
434    if (this->hasAttribute("coordinates", &name, path))
435    {
436      return this->getAttributeValue("coordinates", &name, path);
437    }
438    else
439    {
440      std::list<StdString> dims = this->getDimensionsList(&name, path);
441      std::list<StdString>::const_iterator it = dims.begin(), end = dims.end();
442      for (; it != end; it++)
443      {
444        const StdString& value = *it;
445        retvalue.append(value).push_back(' ');
446      }
447      retvalue.erase(retvalue.end() - 1) ;
448    }
449
450    return retvalue;
451  }
452
453  StdString CINetCDF4::getBoundsId(const StdString& name,
454                                   const CVarPath* const path)
455  {
456    StdString retvalue;
457    if (this->hasAttribute("bounds", &name, path))
458      retvalue = this->getAttributeValue("bounds", &name, path);
459    return retvalue;
460  }
461
462  //---------------------------------------------------------------
463
464  bool CINetCDF4::isBound(const StdString& name,
465                          const CVarPath* const path)
466  {
467    std::set<StdString> bounds = this->getBoundVariables(path);
468    return (bounds.find(name) != bounds.end());
469  }
470
471  bool CINetCDF4::isCoordinate(const StdString& name,
472                               const CVarPath* const path)
473  {
474    std::set<StdString> coords = this->getCoordVariables(path);
475    return (coords.find(name) != coords.end());
476  }
477
478  bool CINetCDF4::isRectilinear(const StdString& name, const CVarPath* const path)
479  {
480    std::list<StdString> coords = this->getCoordinatesIdList(name, path);
481    std::list<StdString>::const_iterator it = coords.begin(), end = coords.end();
482    for (; it != end; it++)
483    {
484      const StdString& coord = *it;
485      if (this->hasVariable(coord, path) && !this->isTemporal(coord, path))
486      {
487        std::map<StdString, StdSize> dimvar = this->getDimensions(&coord, path);
488        if ((dimvar.size() == 1) && (dimvar.find(coord) != dimvar.end()))
489          continue;
490        else
491          return false;
492      }
493    }
494    return true;
495  }
496
497  bool CINetCDF4::isCurvilinear(const StdString& name, const CVarPath* const path)
498  {
499    if (this->isRectilinear(name, path) || !this->hasCoordinates(name, path))
500      return false;
501
502    std::list<StdString> coords = this->getCoordinatesIdList(name, path);
503    std::list<StdString>::const_iterator it = coords.begin(), end = coords.end();
504    for (; it != end; it++)
505    {
506      const StdString& coord = *it;
507      if (this->hasVariable(coord, path))
508      {
509        std::map<StdString, StdSize> dimvar = this->getDimensions(&coord, path);
510        if (dimvar.size() != 2) return false;
511      }
512      else return false;
513    }
514    return true;
515  }
516
517  bool CINetCDF4::isUnstructured(const StdString& name, const CVarPath* const path)
518  {
519    if (this->isRectilinear(name, path) ||
520        this->isCurvilinear(name, path) ||
521        !this->hasCoordinates(name, path))
522       return false;
523
524    StdString dimname = this->getDimensionsList(&name, path).back();
525
526    std::list<StdString> coords = this->getCoordinatesIdList(name, path);
527    std::list<StdString>::const_iterator it = coords.begin(), end = coords.end();
528    for (; it != end; it++)
529    {
530      const StdString& coord = *it;
531      if (this->hasVariable(coord, path))
532      {
533        std::map<StdString, StdSize> dimvar = this->getDimensions(&coord, path);
534        if ((dimvar.size() == 1) &&
535            (dimvar.find(dimname) != dimvar.end()))
536          continue;
537        else
538          return false;
539      }
540      else return false;
541    }
542
543    return true;
544  }
545
546  bool CINetCDF4::isUnknown(const StdString& name, const CVarPath* const path)
547  {
548    return !(this->isRectilinear(name, path) || this->isCurvilinear(name, path) || this->isUnstructured(name, path));
549  }
550
551  bool CINetCDF4::isTemporal(const StdString& name, const CVarPath* const path)
552  {
553    if (!this->hasTemporalDim(path)) return false;
554    std::map<StdString, StdSize> dims = this->getDimensions(&name, path);
555    return (dims.find(this->getUnlimitedDimensionName(path)) != dims.end());
556  }
557
558  bool CINetCDF4::is3Dim(const StdString& name, const CVarPath* const path)
559  {
560    int i = 0;
561    std::list<StdString> coords = this->getCoordinatesIdList(name, path);
562    std::list<StdString>::const_iterator it = coords.begin(), end = coords.end();
563    for (; it != end; it++)
564    {
565      const StdString& coord = *it;
566      if (this->hasVariable(coord, path))
567      {
568        if (this->isTemporal(coord, path))
569          continue;
570        i++;
571      }
572      else
573      {
574        if (coord.compare(this->getUnlimitedDimensionName()) == 0)
575          continue;
576        i++;
577      }
578    }
579    return (i == 3);
580  }
581
582  bool CINetCDF4::isCellGrid(const StdString& name, const CVarPath* const path)
583  {
584    if (this->isCoordinate(name, path))
585    {
586      return this->hasBounds(name, path);
587    }
588    else
589    {
590      std::list<StdString> coords = this->getCoordinatesIdList(name, path);
591      std::list<StdString>::const_iterator it = coords.begin(), end = coords.end();
592      for (; it != end; it++)
593      {
594        const StdString& coord = *it;
595        if (this->hasVariable(coord, path))
596        {
597          if (this->isTemporal(coord, path))
598            continue;
599          if (this->isCellGrid(coord, path))
600            continue;
601          return false;
602        }
603        else
604        {
605          if (coord.compare(this->getUnlimitedDimensionName()) == 0)
606            continue;
607          return false;
608        }
609      }
610    }
611
612    return true;
613  }
614
615  //---------------------------------------------------------------
616
617  std::list<StdString> CINetCDF4::getDataVariables(bool _is3D,       bool _isRecti,
618                                                   bool _isCurvi,    bool _isUnstr,
619                                                   bool _isCellData, bool _isTemporal,
620                                                   const CVarPath* const path)
621  {
622    std::list<StdString> retvalue;
623    std::list<StdString> allvars  = this->getVariables(path);
624    std::set<StdString> allcoords = this->getCoordVariables(path);
625
626    std::list<StdString>::const_iterator it = allvars.begin(), end = allvars.end();
627    for (; it != end; it++)
628    {
629      const StdString& var = *it;
630      if (this->isCoordinate(var, path)) continue;
631
632      if (!_isRecti && this->isRectilinear(var, path))  continue;
633      if (!_isCurvi && this->isCurvilinear(var, path))  continue;
634      if (!_isUnstr && this->isUnstructured(var, path)) continue;
635
636      if (!_isTemporal && this->isTemporal(var, path)) continue;
637      if (!_is3D       && this->is3Dim(var, path))     continue;
638      if (!_isCellData && this->isCellGrid(var, path)) continue;
639
640      if (this->isUnknown(var, path)) continue;
641
642      retvalue.push_back(var);
643    }
644    return retvalue;
645  }
646
647  //---------------------------------------------------------------
648
649  void CINetCDF4::getDataInfo(const StdString& var, const CVarPath* const path, StdSize record,
650                              std::vector<StdSize>& sstart, std::vector<StdSize>& scount, StdSize& array_size,
651                              const std::vector<StdSize>* start /*= NULL*/, const std::vector<StdSize>* count /*= NULL*/)
652  {
653    std::list<StdString> dimlist = this->getDimensionsList(&var, path);
654    std::map<StdString, StdSize> dimmap = this->getDimensions(&var, path);
655    std::list<StdString>::iterator it = dimlist.begin();
656    if (this->isTemporal(var, path))
657    {
658      if (record != UNLIMITED_DIM)
659        sstart.push_back(record);
660      else
661        sstart.push_back(0);
662      scount.push_back(1);
663      it++;
664    }
665    for (int i = 0; it != dimlist.end(); it++, i++)
666    {
667      if (start && count)
668      {
669        sstart.push_back((*start)[i]);
670        scount.push_back((*count)[i]);
671        array_size *= (*count)[i];
672      }
673      else
674      {
675        sstart.push_back(0);
676        scount.push_back(dimmap[*it]);
677        array_size *= dimmap[*it];
678      }
679    }
680  }
681
682  template <class T>
683  void CINetCDF4::getData(CArray<T, 1>& data, const StdString& var,
684                          const CVarPath* const path, StdSize record)
685  {
686    std::vector<StdSize> start, count;
687    int grpid = this->getGroup(path);
688    int varid = this->getVariable(var, path);
689    StdSize array_size = 1;
690    this->getDataInfo(var, path, record, start, count, array_size);
691    data.resize(array_size);
692    CNetCdfInterface::getVaraType(grpid, varid, &start[0], &count[0], data.dataFirst());
693  }
694
695  template <>
696  void CINetCDF4::getData(CArray<int, 1>& data, const StdString& var,
697                          const CVarPath* const path, StdSize record);
698  template <>
699  void CINetCDF4::getData(CArray<double, 1>& data, const StdString& var,
700                          const CVarPath* const path, StdSize record);
701  template <>
702  void CINetCDF4::getData(CArray<float, 1>& data, const StdString& var,
703                          const CVarPath* const path, StdSize record);
704
705  template <>
706  void CINetCDF4::getData(CArray<double, 1>& data, const StdString& var,
707                          bool collective, StdSize record,
708                          const std::vector<StdSize>* start /*= NULL*/,
709                          const std::vector<StdSize>* count /*= NULL*/)
710  {
711    int varid = this->getVariable(var);
712
713    if (this->mpi && collective)
714      CNetCdfInterface::varParAccess(ncidp, varid, NC_COLLECTIVE);
715    else if (this->mpi && !collective)
716      CNetCdfInterface::varParAccess(ncidp, varid, NC_INDEPENDENT);
717
718    std::vector<StdSize> sstart, scount;
719    StdSize array_size = 1;
720    this->getDataInfo(var, NULL, record, sstart, scount, array_size, start, count);
721
722    if (data.numElements() != array_size)
723    {
724      ERROR("CONetCDF4::getData(...)",
725            << "[ Array size = " << data.numElements()
726            << ", Data size = "  << array_size
727            << " ] Invalid array size");
728    }
729
730    CNetCdfInterface::getVaraType(ncidp, varid, &sstart[0], &scount[0], data.dataFirst());
731  }
732
733  //---------------------------------------------------------------
734
735  StdString CINetCDF4::getLonCoordName(const StdString& varname,
736                                       const CVarPath* const path)
737  {
738    std::list<StdString> clist = this->getCoordinatesIdList(varname, path);
739    if (this->hasCoordinates(varname, path))
740      return *clist.begin();
741    else
742      return *clist.rbegin();
743  }
744
745  StdString CINetCDF4::getLatCoordName(const StdString& varname,
746                                       const CVarPath* const path)
747  {
748    std::list<StdString> clist = this->getCoordinatesIdList(varname, path);
749    if (this->hasCoordinates(varname, path))
750      return *(++clist.begin());
751    else
752      return *(++clist.rbegin());
753  }
754
755  StdString CINetCDF4::getVertCoordName(const StdString& varname,
756                                        const CVarPath* const path)
757  {
758    if (!this->is3Dim(varname, path)) return "";
759    std::list<StdString> clist = this->getCoordinatesIdList(varname, path);
760    if (this->hasCoordinates(varname, path))
761      return *(++(++clist.begin()));
762    else
763      return *(++(++clist.rbegin()));
764  }
765} // namespace xios
Note: See TracBrowser for help on using the repository browser.