moab
moab::ScdNCHelper Class Reference

Child helper class for scd mesh, e.g. CAM_EL or CAM_FV. More...

#include <NCHelper.hpp>

Inheritance diagram for moab::ScdNCHelper:
moab::NCHelper moab::NCHelperEuler moab::NCHelperFV

List of all members.

Public Member Functions

 ScdNCHelper (ReadNC *readNC, int fileId, const FileOptions &opts, EntityHandle fileSet)
virtual ~ScdNCHelper ()

Protected Attributes

int gDims [6]
 Dimensions of global grid in file.
int lDims [6]
 Dimensions of my local part of grid.
int gCDims [6]
 Center dimensions of global grid in file.
int lCDims [6]
 Center dimensions of my local part of grid.
std::vector< double > ilVals
 Values for i/j.
std::vector< double > jlVals
std::vector< double > ilCVals
 Center values for i/j.
std::vector< double > jlCVals
int iDim
 Dimension numbers for i/j.
int jDim
int iCDim
 Center dimension numbers for i/j.
int jCDim
int locallyPeriodic [3]
 Whether mesh is locally periodic in i or j or k.
int globallyPeriodic [3]
 Whether mesh is globally periodic in i or j or k.

Private Member Functions

virtual ErrorCode check_existing_mesh ()
 Implementation of NCHelper::check_existing_mesh()
virtual ErrorCode create_mesh (Range &faces)
 Implementation of NCHelper::create_mesh()
virtual ErrorCode read_variables (std::vector< std::string > &var_names, std::vector< int > &tstep_nums)
 Implementation of NCHelper::read_variables()
ErrorCode read_scd_variable_to_nonset_allocate (std::vector< ReadNC::VarData > &vdatas, std::vector< int > &tstep_nums)
 Read non-set variables for scd mesh.
ErrorCode read_scd_variable_to_nonset (std::vector< ReadNC::VarData > &vdatas, std::vector< int > &tstep_nums)
ErrorCode create_quad_coordinate_tag ()
 Create COORDS tag for quads coordinate.
template<typename T >
ErrorCode kji_to_jik (size_t ni, size_t nj, size_t nk, void *dest, T *source)

Detailed Description

Child helper class for scd mesh, e.g. CAM_EL or CAM_FV.

Definition at line 101 of file NCHelper.hpp.


Constructor & Destructor Documentation

moab::ScdNCHelper::ScdNCHelper ( ReadNC readNC,
int  fileId,
const FileOptions opts,
EntityHandle  fileSet 
) [inline]

Definition at line 104 of file NCHelper.hpp.

: NCHelper(readNC, fileId, opts, fileSet),
  iDim(-1), jDim(-1), iCDim(-1), jCDim(-1)
  {
    for (unsigned int i = 0; i < 6; i++) {
      gDims[i] = -1;
      lDims[i] = -1;
      gCDims[i] = -1;
      lCDims[i] = -1;
    }

    locallyPeriodic[0] = locallyPeriodic[1] = locallyPeriodic[2] = 0;
    globallyPeriodic[0] = globallyPeriodic[1] = globallyPeriodic[2] = 0;
  }
virtual moab::ScdNCHelper::~ScdNCHelper ( ) [inline, virtual]

Definition at line 118 of file NCHelper.hpp.

{}

Member Function Documentation

Implementation of NCHelper::check_existing_mesh()

Implements moab::NCHelper.

Definition at line 872 of file NCHelper.cpp.

                                           {
  Interface*& mbImpl = _readNC->mbImpl;

  // Get the number of vertices
  int num_verts;
  ErrorCode rval = mbImpl->get_number_entities_by_dimension(_fileSet, 0, num_verts);
  ERRORR(rval, "Trouble getting number of vertices.");

  /*
  // Check against parameters
  // When ghosting is used, this check might fail (to be updated later)
  if (num_verts > 0) {
    int expected_verts = (lDims[3] - lDims[0] + 1) * (lDims[4] - lDims[1] + 1) * (-1 == lDims[2] ? 1 : lDims[5] - lDims[2] + 1);
    if (num_verts != expected_verts) {
      ERRORR(MB_FAILURE, "Number of vertices doesn't match.");
    }
  }
  */

  // Check the number of elements too
  int num_elems;
  rval = mbImpl->get_number_entities_by_dimension(_fileSet, (-1 == lCDims[2] ? 2 : 3), num_elems);
  ERRORR(rval, "Trouble getting number of elements.");

  /*
  // Check against parameters
  // When ghosting is used, this check might fail (to be updated later)
  if (num_elems > 0) {
    int expected_elems = (lCDims[3] - lCDims[0] + 1) * (lCDims[4] - lCDims[1] + 1) * (-1 == lCDims[2] ? 1 : (lCDims[5] - lCDims[2] + 1));
    if (num_elems != expected_elems) {
      ERRORR(MB_FAILURE, "Number of elements doesn't match.");
    }
  }
  */

  return MB_SUCCESS;
}
ErrorCode moab::ScdNCHelper::create_mesh ( Range faces) [private, virtual]

Implementation of NCHelper::create_mesh()

Implements moab::NCHelper.

Definition at line 910 of file NCHelper.cpp.

{
  Interface*& mbImpl = _readNC->mbImpl;
  Tag& mGlobalIdTag = _readNC->mGlobalIdTag;
  const Tag*& mpFileIdTag = _readNC->mpFileIdTag;
  DebugOutput& dbgOut = _readNC->dbgOut;
  ScdInterface* scdi = _readNC->scdi;
  ScdParData& parData = _readNC->parData;

  Range tmp_range;
  ScdBox* scd_box;

  ErrorCode rval = scdi->construct_box(HomCoord(lDims[0], lDims[1], lDims[2], 1), HomCoord(lDims[3], lDims[4], lDims[5], 1), 
                                       NULL, 0, scd_box, locallyPeriodic, &parData, true);
  ERRORR(rval, "Trouble creating scd vertex sequence.");

  // Add verts to tmp_range first, so we can duplicate global ids in vertex ids
  tmp_range.insert(scd_box->start_vertex(), scd_box->start_vertex() + scd_box->num_vertices() - 1);

  if (mpFileIdTag) {
    int count;
    void* data;
    rval = mbImpl->tag_iterate(*mpFileIdTag, tmp_range.begin(), tmp_range.end(), count, data);
    ERRORR(rval, "Failed to get tag iterator on file id tag.");
    assert(count == scd_box->num_vertices());
    int* fid_data = (int*) data;
    rval = mbImpl->tag_iterate(mGlobalIdTag, tmp_range.begin(), tmp_range.end(), count, data);
    ERRORR(rval, "Failed to get tag iterator on global id tag.");
    assert(count == scd_box->num_vertices());
    int* gid_data = (int*) data;
    for (int i = 0; i < count; i++)
      fid_data[i] = gid_data[i];
  }

  // Then add box set and elements to the range, then to the file set
  tmp_range.insert(scd_box->start_element(), scd_box->start_element() + scd_box->num_elements() - 1);
  tmp_range.insert(scd_box->box_set());
  rval = mbImpl->add_entities(_fileSet, tmp_range);
  ERRORR(rval, "Couldn't add new vertices to file set.");

  dbgOut.tprintf(1, "scdbox %d quads, %d vertices\n", scd_box->num_elements(), scd_box->num_vertices());

  // Set the vertex coordinates
  double *xc, *yc, *zc;
  rval = scd_box->get_coordinate_arrays(xc, yc, zc);
  ERRORR(rval, "Couldn't get vertex coordinate arrays.");

  int i, j, k, il, jl, kl;
  int dil = lDims[3] - lDims[0] + 1;
  int djl = lDims[4] - lDims[1] + 1;
  assert(dil == (int)ilVals.size() && djl == (int)jlVals.size() &&
      (-1 == lDims[2] || lDims[5] - lDims[2] + 1 == (int)levVals.size()));
//#define INDEX(i, j, k) ()
  for (kl = lDims[2]; kl <= lDims[5]; kl++) {
    k = kl - lDims[2];
    for (jl = lDims[1]; jl <= lDims[4]; jl++) {
      j = jl - lDims[1];
      for (il = lDims[0]; il <= lDims[3]; il++) {
        i = il - lDims[0];
        unsigned int pos = i + j * dil + k * dil * djl;
        xc[pos] = ilVals[i];
        yc[pos] = jlVals[j];
        zc[pos] = (-1 == lDims[2] ? 0.0 : levVals[k]);
      }
    }
  }
//#undef INDEX

#ifndef NDEBUG
  int num_verts = (lDims[3] - lDims[0] + 1) * (lDims[4] - lDims[1] + 1) * (-1 == lDims[2] ? 1 : lDims[5] - lDims[2] + 1);
  std::vector<int> gids(num_verts);
  Range verts(scd_box->start_vertex(), scd_box->start_vertex() + scd_box->num_vertices() - 1);
  rval = mbImpl->tag_get_data(mGlobalIdTag, verts, &gids[0]);
  ERRORR(rval, "Trouble getting gid values.");
  int vmin = *(std::min_element(gids.begin(), gids.end())), vmax = *(std::max_element(gids.begin(), gids.end()));
  dbgOut.tprintf(1, "Vertex gids %d-%d\n", vmin, vmax);
#endif

  // Add elements to the range passed in
  faces.insert(scd_box->start_element(), scd_box->start_element() + scd_box->num_elements() - 1);

  if (2 <= dbgOut.get_verbosity()) {
    assert(scd_box->boundary_complete());
    EntityHandle dum_ent = scd_box->start_element();
    rval = mbImpl->list_entities(&dum_ent, 1);
    ERRORR(rval, "Trouble listing first hex.");

    std::vector<EntityHandle> connect;
    rval = mbImpl->get_connectivity(&dum_ent, 1, connect);
    ERRORR(rval, "Trouble getting connectivity.");

    rval = mbImpl->list_entities(&connect[0], connect.size());
    ERRORR(rval, "Trouble listing element connectivity.");
  }

  Range edges;
  mbImpl->get_adjacencies(faces, 1, true, edges, Interface::UNION);

  return MB_SUCCESS;
}

Create COORDS tag for quads coordinate.

Definition at line 1288 of file NCHelper.cpp.

                                                  {
  Interface*& mbImpl = _readNC->mbImpl;

  Range ents;
  ErrorCode rval = mbImpl->get_entities_by_type(_fileSet, moab::MBQUAD, ents);
  ERRORR(rval, "Trouble getting QUAD entity.");

  std::size_t numOwnedEnts = 0;
#ifdef USE_MPI
  Range ents_owned;
  bool& isParallel = _readNC->isParallel;
  if (isParallel) {
    ParallelComm*& myPcomm = _readNC->myPcomm;
    rval = myPcomm->filter_pstatus(ents, PSTATUS_NOT_OWNED, PSTATUS_NOT, -1, &ents_owned);
    ERRORR(rval, "Trouble getting owned QUAD entity.");
    numOwnedEnts = ents_owned.size();
  }
  else
  {
    numOwnedEnts = ents.size();
    ents_owned = ents;
  }
#else
  numOwnedEnts = ents.size();
#endif

  if (numOwnedEnts == 0)
    return MB_SUCCESS;

  assert(numOwnedEnts == ilCVals.size() * jlCVals.size());
  std::vector<double> coords(numOwnedEnts * 3);
  std::size_t pos = 0;
  for (std::size_t j = 0; j != jlCVals.size(); ++j) {
    for (std::size_t i = 0; i != ilCVals.size(); ++i) {
      pos = j * ilCVals.size() * 3 + i * 3;
      coords[pos] = ilCVals[i];
      coords[pos + 1] = jlCVals[j];
      coords[pos + 2] = 0.0;
    }
  }
  std::string tag_name = "COORDS";
  Tag tagh = 0;
  rval = mbImpl->tag_get_handle(tag_name.c_str(), 3, MB_TYPE_DOUBLE, tagh, MB_TAG_DENSE | MB_TAG_CREAT);
  ERRORR(rval, "Trouble creating COORDS tag.");

  void *data;
  int count;
#ifdef USE_MPI
  rval = mbImpl->tag_iterate(tagh, ents_owned.begin(), ents_owned.end(), count, data);
#else
  rval = mbImpl->tag_iterate(tagh, ents.begin(), ents.end(), count, data);
#endif
  ERRORR(rval, "Failed to get COORDS tag iterator.");
  assert(count == (int)numOwnedEnts);
  double* quad_data = (double*) data;
  std::copy(coords.begin(), coords.end(), quad_data);

  return MB_SUCCESS;
}
template<typename T >
ErrorCode moab::ScdNCHelper::kji_to_jik ( size_t  ni,
size_t  nj,
size_t  nk,
void *  dest,
T source 
) [inline, private]

Definition at line 137 of file NCHelper.hpp.

  {
    size_t nik = ni * nk, nij = ni * nj;
    T* tmp_data = reinterpret_cast<T*>(dest);
    for (std::size_t j = 0; j != nj; j++)
      for (std::size_t i = 0; i != ni; i++)
        for (std::size_t k = 0; k != nk; k++)
          tmp_data[j*nik + i*nk + k] = source[k*nij + j*ni + i];
    return MB_SUCCESS;
  }
ErrorCode moab::ScdNCHelper::read_scd_variable_to_nonset ( std::vector< ReadNC::VarData > &  vdatas,
std::vector< int > &  tstep_nums 
) [private]

Definition at line 1162 of file NCHelper.cpp.

{
  DebugOutput& dbgOut = _readNC->dbgOut;

  ErrorCode rval = read_scd_variable_to_nonset_allocate(vdatas, tstep_nums);
  ERRORR(rval, "Trouble allocating read variables.");

  // Finally, read into that space
  int success;
  for (unsigned int i = 0; i < vdatas.size(); i++) {
    std::size_t sz = vdatas[i].sz;

    // A typical supported variable: float T(time, lev, lat, lon)
    // For tag values, need transpose (lev, lat, lon) to (lat, lon, lev)
    size_t ni = vdatas[i].readCounts[3]; // lon or slon
    size_t nj = vdatas[i].readCounts[2]; // lat or slat
    size_t nk = vdatas[i].readCounts[1]; // lev

    for (unsigned int t = 0; t < tstep_nums.size(); t++) {
      // Tag data for this timestep
      void* data = vdatas[i].varDatas[t];

      // Set readStart for each timestep along time dimension
      vdatas[i].readStarts[0] = tstep_nums[t];

      switch (vdatas[i].varDataType) {
        case NC_BYTE:
        case NC_CHAR: {
          std::vector<char> tmpchardata(sz);
          success = NCFUNCAG(_vara_text)(_fileId, vdatas[i].varId, &vdatas[i].readStarts[0], &vdatas[i].readCounts[0],
                                        &tmpchardata[0]);
          if (vdatas[i].numLev != 1)
            // Transpose (lev, lat, lon) to (lat, lon, lev)
            success = kji_to_jik(ni, nj, nk, data, &tmpchardata[0]);
          else {
            for (std::size_t idx = 0; idx != tmpchardata.size(); idx++)
              ((char*) data)[idx] = tmpchardata[idx];
          }
          ERRORS(success, "Failed to read char data.");
          break;
        }
        case NC_DOUBLE: {
          std::vector<double> tmpdoubledata(sz);
          success = NCFUNCAG(_vara_double)(_fileId, vdatas[i].varId, &vdatas[i].readStarts[0], &vdatas[i].readCounts[0],
                                          &tmpdoubledata[0]);
          if (vdatas[i].numLev != 1)
            // Transpose (lev, lat, lon) to (lat, lon, lev)
            success = kji_to_jik(ni, nj, nk, data, &tmpdoubledata[0]);
          else {
            for (std::size_t idx = 0; idx != tmpdoubledata.size(); idx++)
              ((double*) data)[idx] = tmpdoubledata[idx];
          }
          ERRORS(success, "Failed to read double data.");
          break;
        }
        case NC_FLOAT: {
          std::vector<float> tmpfloatdata(sz);
          success = NCFUNCAG(_vara_float)(_fileId, vdatas[i].varId, &vdatas[i].readStarts[0], &vdatas[i].readCounts[0],
                                          &tmpfloatdata[0]);
          if (vdatas[i].numLev != 1)
            // Transpose (lev, lat, lon) to (lat, lon, lev)
            success = kji_to_jik(ni, nj, nk, data, &tmpfloatdata[0]);
          else {
            for (std::size_t idx = 0; idx != tmpfloatdata.size(); idx++)
              ((float*) data)[idx] = tmpfloatdata[idx];
          }
          ERRORS(success, "Failed to read float data.");
          break;
        }
        case NC_INT: {
          std::vector<int> tmpintdata(sz);
          success = NCFUNCAG(_vara_int)(_fileId, vdatas[i].varId, &vdatas[i].readStarts[0], &vdatas[i].readCounts[0],
                                        &tmpintdata[0]);
          if (vdatas[i].numLev != 1)
            // Transpose (lev, lat, lon) to (lat, lon, lev)
            success = kji_to_jik(ni, nj, nk, data, &tmpintdata[0]);
          else {
            for (std::size_t idx = 0; idx != tmpintdata.size(); idx++)
              ((int*) data)[idx] = tmpintdata[idx];
          }
          ERRORS(success, "Failed to read int data.");
          break;
        }
        case NC_SHORT: {
          std::vector<short> tmpshortdata(sz);
          success = NCFUNCAG(_vara_short)(_fileId, vdatas[i].varId, &vdatas[i].readStarts[0], &vdatas[i].readCounts[0],
                                          &tmpshortdata[0]);
          if (vdatas[i].numLev != 1)
            // Transpose (lev, lat, lon) to (lat, lon, lev)
            success = kji_to_jik(ni, nj, nk, data, &tmpshortdata[0]);
          else {
            for (std::size_t idx = 0; idx != tmpshortdata.size(); idx++)
              ((short*) data)[idx] = tmpshortdata[idx];
          }
          ERRORS(success, "Failed to read short data.");
          break;
        }
        default:
          success = 1;
      }

      if (success)
        ERRORR(MB_FAILURE, "Trouble reading variable.");
    }
  }

  for (unsigned int i = 0; i < vdatas.size(); i++) {
    for (unsigned int t = 0; t < tstep_nums.size(); t++) {
      dbgOut.tprintf(2, "Converting variable %s, time step %d\n", vdatas[i].varName.c_str(), tstep_nums[t]);
      ErrorCode tmp_rval = convert_variable(vdatas[i], t);
      if (MB_SUCCESS != tmp_rval)
        rval = tmp_rval;
    }
  }

  // Debug output, if requested
  if (1 == dbgOut.get_verbosity()) {
    dbgOut.printf(1, "Read variables: %s", vdatas.begin()->varName.c_str());
    for (unsigned int i = 1; i < vdatas.size(); i++)
      dbgOut.printf(1, ", %s ", vdatas[i].varName.c_str());
    dbgOut.tprintf(1, "\n");
  }

  return rval;
}
ErrorCode moab::ScdNCHelper::read_scd_variable_to_nonset_allocate ( std::vector< ReadNC::VarData > &  vdatas,
std::vector< int > &  tstep_nums 
) [private]

Read non-set variables for scd mesh.

Definition at line 1036 of file NCHelper.cpp.

{
  Interface*& mbImpl = _readNC->mbImpl;
  std::vector<int>& dimLens = _readNC->dimLens;
  DebugOutput& dbgOut = _readNC->dbgOut;

  ErrorCode rval = MB_SUCCESS;

  Range* range = NULL;

  // Get vertices in set
  Range verts;
  rval = mbImpl->get_entities_by_dimension(_fileSet, 0, verts);
  ERRORR(rval, "Trouble getting vertices in set.");
  assert("Should only have a single vertex subrange, since they were read in one shot" &&
      verts.psize() == 1);

  Range edges;
  rval = mbImpl->get_entities_by_dimension(_fileSet, 1, edges);
  ERRORR(rval, "Trouble getting edges in set.");

  // Get faces in set
  Range faces;
  rval = mbImpl->get_entities_by_dimension(_fileSet, 2, faces);
  ERRORR(rval, "Trouble getting faces in set.");
  assert("Should only have a single face subrange, since they were read in one shot" &&
      faces.psize() == 1);

#ifdef USE_MPI
  moab::Range faces_owned;
  bool& isParallel = _readNC->isParallel;
  if (isParallel) {
    ParallelComm*& myPcomm = _readNC->myPcomm;
    rval = myPcomm->filter_pstatus(faces, PSTATUS_NOT_OWNED, PSTATUS_NOT, -1, &faces_owned);
    ERRORR(rval, "Trouble getting owned faces in set.");
  }
  else
    faces_owned = faces; // Not running in parallel, but still with MPI
#endif

  for (unsigned int i = 0; i < vdatas.size(); i++) {
    // Support non-set variables with 4 dimensions like (time, lev, lat, lon)
    assert(4 == vdatas[i].varDims.size());

    // For a non-set variable, time should be the first dimension
    assert(tDim == vdatas[i].varDims[0]);

    // Set up readStarts and readCounts
    vdatas[i].readStarts.resize(4);
    vdatas[i].readCounts.resize(4);

    // First: time
    vdatas[i].readStarts[0] = 0; // This value is timestep dependent, will be set later
    vdatas[i].readCounts[0] = 1;

    // Next: lev
    vdatas[i].readStarts[1] = 0;
    vdatas[i].readCounts[1] = vdatas[i].numLev;

    // Finally: lat (or slat) and lon (or slon)
    switch (vdatas[i].entLoc) {
      case ReadNC::ENTLOCVERT:
        // Vertices
        vdatas[i].readStarts[2] = lDims[1];
        vdatas[i].readCounts[2] = lDims[4] - lDims[1] + 1;
        vdatas[i].readStarts[3] = lDims[0];
        vdatas[i].readCounts[3] = lDims[3] - lDims[0] + 1;
        range = &verts;
        break;
      case ReadNC::ENTLOCNSEDGE:
        ERRORR(MB_FAILURE, "Reading edge data not implemented yet.");
        break;
      case ReadNC::ENTLOCEWEDGE:
        ERRORR(MB_FAILURE, "Reading edge data not implemented yet.");
        break;
      case ReadNC::ENTLOCFACE:
        // Faces
        vdatas[i].readStarts[2] = lCDims[1];
        vdatas[i].readCounts[2] = lCDims[4] - lCDims[1] + 1;
        vdatas[i].readStarts[3] = lCDims[0];
        vdatas[i].readCounts[3] = lCDims[3] - lCDims[0] + 1;
#ifdef USE_MPI
        range = &faces_owned;
#else
        range = &faces;
#endif
        break;
      case ReadNC::ENTLOCSET:
        // Set
        break;
      default:
        ERRORR(MB_FAILURE, "Unrecognized entity location type.");
        break;
    }

    for (unsigned int t = 0; t < tstep_nums.size(); t++) {
      dbgOut.tprintf(2, "Reading variable %s, time step %d\n", vdatas[i].varName.c_str(), tstep_nums[t]);

      if (tstep_nums[t] >= dimLens[tDim]) {
        ERRORR(MB_INDEX_OUT_OF_RANGE, "Wrong value for a timestep number.");
      }

      // Get the tag to read into
      if (!vdatas[i].varTags[t]) {
        rval = get_tag_to_nonset(vdatas[i], tstep_nums[t], vdatas[i].varTags[t], vdatas[i].numLev);
        ERRORR(rval, "Trouble getting tag.");
      }

      // Get ptr to tag space
      void* data;
      int count;
      rval = mbImpl->tag_iterate(vdatas[i].varTags[t], range->begin(), range->end(), count, data);
      ERRORR(rval, "Failed to get tag iterator.");
      assert((unsigned)count == range->size());
      vdatas[i].varDatas[t] = data;
    }

    // Get variable size
    vdatas[i].sz = 1;
    for (std::size_t idx = 0; idx != vdatas[i].readCounts.size(); idx++)
      vdatas[i].sz *= vdatas[i].readCounts[idx];
  }

  return rval;
}
ErrorCode moab::ScdNCHelper::read_variables ( std::vector< std::string > &  var_names,
std::vector< int > &  tstep_nums 
) [private, virtual]

Implementation of NCHelper::read_variables()

Implements moab::NCHelper.

Definition at line 1011 of file NCHelper.cpp.

{
  std::vector<ReadNC::VarData> vdatas;
  std::vector<ReadNC::VarData> vsetdatas;

  ErrorCode rval = read_variable_setup(var_names, tstep_nums, vdatas, vsetdatas);
  ERRORR(rval, "Trouble setting up read variable.");

  // Create COORDS tag for quads
  rval = create_quad_coordinate_tag();
  ERRORR(rval, "Trouble creating coordinate tags to entities quads");

  if (!vsetdatas.empty()) {
    rval = read_variable_to_set(vsetdatas, tstep_nums);
    ERRORR(rval, "Trouble read variables to set.");
  }

  if (!vdatas.empty()) {
    rval = read_scd_variable_to_nonset(vdatas, tstep_nums);
    ERRORR(rval, "Trouble read variables to entities verts/edges/faces.");
  }

  return MB_SUCCESS;
}

Member Data Documentation

int moab::ScdNCHelper::gCDims[6] [protected]

Center dimensions of global grid in file.

Definition at line 156 of file NCHelper.hpp.

int moab::ScdNCHelper::gDims[6] [protected]

Dimensions of global grid in file.

Definition at line 150 of file NCHelper.hpp.

Whether mesh is globally periodic in i or j or k.

Definition at line 177 of file NCHelper.hpp.

int moab::ScdNCHelper::iCDim [protected]

Center dimension numbers for i/j.

Definition at line 171 of file NCHelper.hpp.

int moab::ScdNCHelper::iDim [protected]

Dimension numbers for i/j.

Definition at line 168 of file NCHelper.hpp.

std::vector<double> moab::ScdNCHelper::ilCVals [protected]

Center values for i/j.

Definition at line 165 of file NCHelper.hpp.

std::vector<double> moab::ScdNCHelper::ilVals [protected]

Values for i/j.

Definition at line 162 of file NCHelper.hpp.

int moab::ScdNCHelper::jCDim [protected]

Definition at line 171 of file NCHelper.hpp.

int moab::ScdNCHelper::jDim [protected]

Definition at line 168 of file NCHelper.hpp.

std::vector<double> moab::ScdNCHelper::jlCVals [protected]

Definition at line 165 of file NCHelper.hpp.

std::vector<double> moab::ScdNCHelper::jlVals [protected]

Definition at line 162 of file NCHelper.hpp.

int moab::ScdNCHelper::lCDims[6] [protected]

Center dimensions of my local part of grid.

Definition at line 159 of file NCHelper.hpp.

int moab::ScdNCHelper::lDims[6] [protected]

Dimensions of my local part of grid.

Definition at line 153 of file NCHelper.hpp.

Whether mesh is locally periodic in i or j or k.

Definition at line 174 of file NCHelper.hpp.


The documentation for this class was generated from the following files:
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines