source: XIOS/dev/branch_openmp/src/io/inetcdf4.cpp @ 1491

Last change on this file since 1491 was 1491, checked in by yushan, 6 years ago

Branch EP merged with Dev_cmip6 @r1490

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