Zoltan2
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Zoltan2_MeshAdapter.hpp
Go to the documentation of this file.
1 // @HEADER
2 //
3 // ***********************************************************************
4 //
5 // Zoltan2: A package of combinatorial algorithms for scientific computing
6 // Copyright 2012 Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Karen Devine (kddevin@sandia.gov)
39 // Erik Boman (egboman@sandia.gov)
40 // Siva Rajamanickam (srajama@sandia.gov)
41 //
42 // ***********************************************************************
43 //
44 // @HEADER
45 
51 #ifndef _ZOLTAN2_MESHADAPTER_HPP_
52 #define _ZOLTAN2_MESHADAPTER_HPP_
53 
54 #include <Zoltan2_Adapter.hpp>
55 #include "TpetraExt_MatrixMatrix.hpp"
56 
57 namespace Zoltan2 {
58 
68 };
69 
76  POINT, // a 0D entity (e.g. a vertex)
77  LINE_SEGMENT, // a 1D entity (e.g. an edge)
78  POLYGON, // a general 2D entity
79  TRIANGLE, // a specific 2D entity bounded by 3 edge entities
80  QUADRILATERAL, // a specific 2D entity bounded by 4 edge entities
81  POLYHEDRON, // a general 3D entity
82  TETRAHEDRON, // a specific 3D entity bounded by 4 triangle entities
83  HEXAHEDRON, // a specific 3D entity bounded by 6 quadrilateral
84  // entities
85  PRISM, // a specific 3D entity bounded by a combination of 3
86  //quadrilateral entities and 2 triangle entities
87  PYRAMID // a specific 3D entity bounded by a combination of 1
88  // quadrilateral entity and 4 triangle entities
89 };
90 
123 template <typename User>
124 class MeshAdapter : public BaseAdapter<User> {
125 public:
126 
127 #ifndef DOXYGEN_SHOULD_SKIP_THIS
128  typedef typename InputTraits<User>::scalar_t scalar_t;
129  typedef typename InputTraits<User>::offset_t offset_t;
130  typedef typename InputTraits<User>::lno_t lno_t;
131  typedef typename InputTraits<User>::gno_t gno_t;
132  typedef typename InputTraits<User>::part_t part_t;
133  typedef typename InputTraits<User>::node_t node_t;
134  typedef User user_t;
135  typedef User userCoord_t;
137 #endif
138 
140 
143  virtual ~MeshAdapter() {};
144 
145  // Default MeshEntityType is MESH_REGION with MESH_FACE-based adjacencies and
146  // second adjacencies and coordinates
147  MeshAdapter() : primaryEntityType(MESH_REGION),
148  adjacencyEntityType(MESH_FACE),
149  secondAdjacencyEntityType(MESH_FACE) {};
150 
152  // Methods to be defined in derived classes.
153 
158  virtual bool areEntityIDsUnique(MeshEntityType etype) const
159  {
160  return etype==this->getPrimaryEntityType();
161  }
162 
165  //virtual size_t getGlobalNumOf(MeshEntityType etype) const = 0;
166 
169  virtual size_t getLocalNumOf(MeshEntityType etype) const = 0;
170 
171 
176  virtual void getIDsViewOf(MeshEntityType etype,
177  gno_t const *&Ids) const = 0;
178 
179 
184  virtual void getTopologyViewOf(MeshEntityType etype,
185  enum EntityTopologyType const *&Types) const
186  {
187  Types = NULL;
189  }
190 
196  virtual int getNumWeightsPerOf(MeshEntityType etype) const { return 0; }
197 
211  virtual void getWeightsViewOf(MeshEntityType etype,
212  const scalar_t *&weights, int &stride, int idx = 0) const
213  {
214  weights = NULL;
215  stride = 0;
217  }
218 
219 
228  virtual int getDimension() const { return 0; }
229 
241  const scalar_t *&coords, int &stride, int coordDim) const
242  {
243  coords = NULL;
244  stride = 0;
246  }
247 
248 
251  virtual bool availAdjs(MeshEntityType source, MeshEntityType target) const {
252  return false;
253  }
254 
255 
258  virtual size_t getLocalNumAdjs(MeshEntityType source,
259  MeshEntityType target) const { return 0;}
260 
261 
272  virtual void getAdjsView(MeshEntityType source, MeshEntityType target,
273  const offset_t *&offsets, const gno_t *& adjacencyIds) const
274  {
275  offsets = NULL;
276  adjacencyIds = NULL;
278  }
279 
280 
285  virtual bool avail2ndAdjs(MeshEntityType sourcetarget,
286  MeshEntityType through) const
287  {
288  return false;
289  }
290 
294  virtual size_t getLocalNum2ndAdjs(MeshEntityType sourcetarget,
295  MeshEntityType through) const
296  {
297  return 0;
298  }
299 
310  virtual void get2ndAdjsView(MeshEntityType sourcetarget,
311  MeshEntityType through,
312  const offset_t *&offsets,
313  const gno_t *&adjacencyIds) const
314  {
315  offsets = NULL;
316  adjacencyIds = NULL;
318  }
319 
323  virtual int getNumWeightsPer2ndAdj(MeshEntityType sourcetarget,
324  MeshEntityType through) const { return 0;}
325 
326 
336  virtual void get2ndAdjWeightsView(MeshEntityType sourcetarget,
337  MeshEntityType through,
338  const scalar_t *&weights,
339  int &stride,
340  int idx) const
341  {
342  weights = NULL;
343  stride = 0;
345  }
346 
348  // Implementations of base-class methods
349 
352  inline enum MeshEntityType getPrimaryEntityType() const {
353  return this->primaryEntityType;
354  }
355 
362  return this->adjacencyEntityType;
363  }
364 
371  return this->secondAdjacencyEntityType;
372  }
373 
380  void setEntityTypes(std::string ptypestr, std::string atypestr,
381  std::string satypestr) {
382 
383  if (ptypestr != atypestr && ptypestr != satypestr) {
384  if (ptypestr == "region")
385  this->primaryEntityType = MESH_REGION;
386  else if (ptypestr == "face")
387  this->primaryEntityType = MESH_FACE;
388  else if (ptypestr == "edge")
389  this->primaryEntityType = MESH_EDGE;
390  else if (ptypestr == "vertex")
391  this->primaryEntityType = MESH_VERTEX;
392  else {
393  std::ostringstream emsg;
394  emsg << __FILE__ << "," << __LINE__
395  << " error: Invalid MeshEntityType " << ptypestr << std::endl;
396  emsg << "Valid values: region face edge vertex" << std::endl;
397  throw std::runtime_error(emsg.str());
398  }
399 
400  if (atypestr == "region")
401  this->adjacencyEntityType = MESH_REGION;
402  else if (atypestr == "face")
403  this->adjacencyEntityType = MESH_FACE;
404  else if (atypestr == "edge")
405  this->adjacencyEntityType = MESH_EDGE;
406  else if (atypestr == "vertex")
407  this->adjacencyEntityType = MESH_VERTEX;
408  else {
409  std::ostringstream emsg;
410  emsg << __FILE__ << "," << __LINE__
411  << " error: Invalid MeshEntityType " << atypestr << std::endl;
412  emsg << "Valid values: region face edge vertex" << std::endl;
413  throw std::runtime_error(emsg.str());
414  }
415 
416  if (satypestr == "region")
417  this->secondAdjacencyEntityType = MESH_REGION;
418  else if (satypestr == "face")
419  this->secondAdjacencyEntityType = MESH_FACE;
420  else if (satypestr == "edge")
421  this->secondAdjacencyEntityType = MESH_EDGE;
422  else if (satypestr == "vertex")
423  this->secondAdjacencyEntityType = MESH_VERTEX;
424  else {
425  std::ostringstream emsg;
426  emsg << __FILE__ << "," << __LINE__
427  << " error: Invalid MeshEntityType " << satypestr << std::endl;
428  emsg << "Valid values: region face edge vertex" << std::endl;
429  throw std::runtime_error(emsg.str());
430  }
431  }
432  else {
433  std::ostringstream emsg;
434  emsg << __FILE__ << "," << __LINE__
435  << " error: PrimaryEntityType " << ptypestr
436  << " matches AdjacencyEntityType " << atypestr
437  << " or SecondAdjacencyEntityType " << satypestr << std::endl;
438  throw std::runtime_error(emsg.str());
439  }
440  }
441 
446  virtual bool useDegreeAsWeightOf(MeshEntityType etype, int idx) const
447  {
448  return false;
449  }
450 
452  // Functions from the BaseAdapter interface
453  size_t getLocalNumIDs() const {
455  }
456 
457  void getIDsView(const gno_t *&Ids) const {
459  }
460 
461  int getNumWeightsPerID() const {
463  }
464 
465  void getWeightsView(const scalar_t *&wgt, int &stride, int idx = 0) const {
466  getWeightsViewOf(getPrimaryEntityType(), wgt, stride, idx);
467  }
468 
469  void getCoordinatesView(const scalar_t *&coords, int &stride,
470  int coordDim) const
471  {
472  getCoordinatesViewOf(getPrimaryEntityType(), coords, stride, coordDim);
473  }
474 
475  bool useDegreeAsWeight(int idx) const
476  {
478  }
479 
480 private:
481  enum MeshEntityType primaryEntityType; // Entity type
482  // to be partitioned, ordered,
483  // colored, matched, etc.
484  enum MeshEntityType adjacencyEntityType; // Entity type defining first-order
485  // adjacencies; adjacencies are of
486  // this type.
487  enum MeshEntityType secondAdjacencyEntityType; // Bridge entity type
488  // defining second-order
489  // adjacencies.
490 };
491 
492 } //namespace Zoltan2
493 
494 #endif
virtual void getAdjsView(MeshEntityType source, MeshEntityType target, const offset_t *&offsets, const gno_t *&adjacencyIds) const
Sets pointers to this process&#39; mesh first adjacencies.
InputTraits< User >::scalar_t scalar_t
Zoltan2::BaseAdapter< userTypes_t > base_adapter_t
virtual bool areEntityIDsUnique(MeshEntityType etype) const
Provide a pointer to the entity topology types.
virtual bool availAdjs(MeshEntityType source, MeshEntityType target) const
Returns whether a first adjacency combination is available.
void setEntityTypes(std::string ptypestr, std::string atypestr, std::string satypestr)
Sets the primary, adjacency, and second adjacency entity types. Called by algorithm based on paramete...
InputTraits< User >::gno_t gno_t
int getNumWeightsPerID() const
Returns the number of weights per object. Number of weights per object should be zero or greater...
MeshAdapter defines the interface for mesh input.
virtual size_t getLocalNumOf(MeshEntityType etype) const =0
Returns the global number of mesh entities of MeshEntityType.
default_part_t part_t
The data type to represent part numbers.
default_offset_t offset_t
The data type to represent offsets.
virtual ~MeshAdapter()
Destructor.
virtual void getWeightsViewOf(MeshEntityType etype, const scalar_t *&weights, int &stride, int idx=0) const
Provide a pointer to one of the number of this process&#39; optional entity weights.
InputTraits< User >::lno_t lno_t
virtual bool avail2ndAdjs(MeshEntityType sourcetarget, MeshEntityType through) const
Returns whether a second adjacency combination is available. If combination is not available in the M...
enum MeshEntityType getAdjacencyEntityType() const
Returns the entity that describes adjacencies between the entities to be partitioned, ordered, colored, etc. That is, a primaryEntityType that contains an adjacencyEntityType are adjacent.
#define Z2_THROW_NOT_IMPLEMENTED
static ArrayRCP< ArrayRCP< zscalar_t > > weights
virtual void get2ndAdjsView(MeshEntityType sourcetarget, MeshEntityType through, const offset_t *&offsets, const gno_t *&adjacencyIds) const
if avail2ndAdjs(), set pointers to this process&#39; second adjacencies
virtual int getDimension() const
Return dimension of the entity coordinates, if any.
void getCoordinatesView(const scalar_t *&coords, int &stride, int coordDim) const
virtual void get2ndAdjWeightsView(MeshEntityType sourcetarget, MeshEntityType through, const scalar_t *&weights, int &stride, int idx) const
Provide a pointer to the second adjacency weights, if any. Note: second-adjacency weights may be used...
virtual size_t getLocalNum2ndAdjs(MeshEntityType sourcetarget, MeshEntityType through) const
if avail2ndAdjs(), returns the number of second adjacencies on this process.
virtual void getIDsViewOf(MeshEntityType etype, gno_t const *&Ids) const =0
Provide a pointer to this process&#39; identifiers.
virtual void getCoordinatesViewOf(MeshEntityType etype, const scalar_t *&coords, int &stride, int coordDim) const
Provide a pointer to one dimension of entity coordinates.
enum MeshEntityType getSecondAdjacencyEntityType() const
Returns the entity that describes second adjacencies between the entities to be partitioned, ordered, colored, etc. That is, two primaryEntityType that share a secondAdjacencyEntityType are adjacent.
default_lno_t lno_t
The ordinal type (e.g., int, long, int64_t) that represents local counts and local indices...
BaseAdapterType
An enum to identify general types of adapters.
virtual size_t getLocalNumAdjs(MeshEntityType source, MeshEntityType target) const
Returns the number of first adjacencies on this process.
void getWeightsView(const scalar_t *&wgt, int &stride, int idx=0) const
Provide pointer to a weight array with stride.
virtual int getNumWeightsPer2ndAdj(MeshEntityType sourcetarget, MeshEntityType through) const
Returns the number (0 or greater) of weights per second adjacency. Note: second-adjacency weights may...
InputTraits< User >::part_t part_t
EntityTopologyType
Enumerate entity topology types for meshes: points,lines,polygons,triangles,quadrilaterals, polyhedrons, tetrahedrons, hexhedrons, prisms, or pyramids.
default_gno_t gno_t
The ordinal type (e.g., int, long, int64_t) that can represent global counts and identifiers.
default_node_t node_t
The Kokkos node type. This is only meaningful for users of Tpetra objects.
MeshEntityType
Enumerate entity types for meshes: Regions, Faces, Edges, or Vertices.
bool useDegreeAsWeight(int idx) const
virtual int getNumWeightsPerOf(MeshEntityType etype) const
Return the number of weights per entity.
virtual void getTopologyViewOf(MeshEntityType etype, enum EntityTopologyType const *&Types) const
Provide a pointer to the entity topology types.
size_t getLocalNumIDs() const
Returns the number of objects on this process.
virtual bool useDegreeAsWeightOf(MeshEntityType etype, int idx) const
Optional method allowing the idx-th weight of entity type etype to be set as the number of neighbors ...
enum MeshEntityType getPrimaryEntityType() const
Returns the entity to be partitioned, ordered, colored, etc.
InputTraits< User >::offset_t offset_t
default_scalar_t scalar_t
The data type for weights and coordinates.
enum BaseAdapterType adapterType() const
Returns the type of adapter.
void getIDsView(const gno_t *&Ids) const
Provide a pointer to this process&#39; identifiers.
Zoltan2::BasicUserTypes< zscalar_t, zlno_t, zgno_t > user_t
Definition: Metric.cpp:74