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

Last change on this file since 1243 was 948, checked in by ymipsl, 8 years ago

Bug fix : read unstructured grid with one or more axis.

This is a fix but the method for determine grid type from file must be improved.

YM

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