SCOREC core
Parallel unstructured mesh tools
Classes | Typedefs | Enumerations | Functions | Variables
apf Namespace Reference

All APF symbols are contained in this namespace. More...

Classes

class  ReductionOp
 Base class for applying operations to make a Field consistent in parallel. More...
 
class  Integrator
 A virtual base for user-defined integrators. More...
 
struct  Function
 User-defined Analytic Function. More...
 
struct  Up
 statically sized container for upward adjacency queries. More...
 
struct  Copy
 a reference to an object representing the same entity More...
 
class  Mesh
 Interface to a mesh part. More...
 
class  Migration
 Migration plan object: local elements to destinations. More...
 
struct  Sharing
 abstract description of entity copy sharing More...
 
class  Mesh2
 Extended mesh interface for modification. More...
 
class  BuildCallback
 User-defined entity creation callback. More...
 
class  Matrix
 template-generic matrix of M by N doubles More...
 
class  Matrix3x3
 convenience wrapper over apf::Matrix<3,3> More...
 
class  Vector
 template-generic vector of N doubles More...
 
class  Vector3
 convenience wrapper over apf::Vector<3> More...
 
class  DynamicMatrix
 A runtime-sized dense matrix. More...
 
class  DynamicVector
 A runtime-sized linear algebra vector of doubles. More...
 
class  CavityOp
 user-defined mesh cavity operator More...
 
class  EntityShape
 Shape functions over this element. More...
 
class  FieldShape
 Describes field distribution and shape functions. More...
 
struct  Node
 Node identifier. More...
 
class  Splitter
 Splits a mesh part into many. More...
 
class  Balancer
 Load balance over all mesh parts. More...
 
struct  Remap
 a map from old part ids to new part ids More...
 
struct  Divide
 divide the part id More...
 
struct  Multiply
 multiply the part id More...
 
struct  Modulo
 return part id modulo n More...
 
struct  Unmodulo
 inverse of apf::Modulo More...
 
struct  Round
 map to nearest multiple of n More...
 

Typedefs

typedef VectorElement MeshElement
 Mesh Elements represent the mesh coordinate vector field.
 
typedef NumberingOf< int > Numbering
 Numbering is meant to be a 32-bit local numbering.
 
typedef NumberingOf< long > GlobalNumbering
 Global numberings use 64-bit integers.
 
typedef std::map< int, MeshEntity * > Copies
 Remote copy container. More...
 
typedef std::set< int > Parts
 Set of unique part ids.
 
typedef DynamicArray< MeshEntity * > Adjacent
 Set of adjacent mesh entities. More...
 
typedef MeshEntity * Downward[12]
 a static array type downward adjacency queries. More...
 
typedef Copy Match
 matches are just a special case of copies
 
typedef DynamicArray< CopyCopyArray
 a set of copies, possibly multiple copies per part
 
typedef CopyArray Matches
 a set of periodic copies
 
typedef CopyArray DgCopies
 a set of DG copies
 
typedef std::map< Gid, MeshEntity * > GlobalToVert
 a map from global ids to vertex objects
 

Enumerations

enum  ValueType {
  SCALAR , VECTOR , MATRIX , PACKED ,
  VALUE_TYPES
}
 The type of value the field stores. More...
 
enum  ZoltanMethod {
  RCB , RIB , HYPERGRAPH , PARMETIS ,
  GRAPH
}
 Zoltan partitioning method. More...
 
enum  ZoltanApproach {
  PARTITION , REPARTITION , REFINE , PART_KWAY ,
  PART_GEOM , PART_GEOM_KWAY , ADAPT_REPART , REFINE_KWAY
}
 Zoltan partitioning approach. More...
 

Functions

void destroyMesh (Mesh *m)
 Destroys an apf::Mesh. More...
 
MeshElementcreateMeshElement (Mesh *m, MeshEntity *e)
 Creates a Mesh Element over an entity. More...
 
MeshElementcreateMeshElement (Field *c, MeshEntity *e)
 Creates a non-standard Mesh Element over an entity. More...
 
MeshEntity * getMeshEntity (MeshElement *me)
 Retrieve the mesh entity associated with an apf::MeshElement.
 
void destroyMeshElement (MeshElement *e)
 Destroys a Mesh Element. More...
 
Field * createLagrangeField (Mesh *m, const char *name, int valueType, int order)
 Create an apf::Field using a Lagrange distribution. More...
 
Field * createStepField (Mesh *m, const char *name, int valueType)
 Create an apf::Field using a step distribution. More...
 
Field * createIPField (Mesh *m, const char *name, int valueType, int order)
 Create an apf::Field of integration point data. More...
 
Field * createField (Mesh *m, const char *name, int valueType, FieldShape *shape)
 Create a Field from any builtin or user defined FieldShape.
 
Field * createFieldOn (Mesh *m, const char *name, int valueType)
 Create a field using the mesh's coordinate nodal distribution.
 
Field * createPackedField (Mesh *m, const char *name, int components, apf::FieldShape *shape=0)
 Create a field of N components without a tensor type. More...
 
Field * createGeneralField (Mesh *m, const char *name, int valueType, int components, FieldShape *shape)
 Encompasses both packed and typed fields.
 
Field * cloneField (Field *f, Mesh *onto)
 Declare a copy of a field on another apf::Mesh. More...
 
MeshgetMesh (Field *f)
 Retrieve the Mesh over which a Field is defined.
 
bool hasEntity (Field *f, MeshEntity *e)
 Returns true iff an entity has data associate with a field.
 
const char * getName (Field *f)
 Get the name of a Field. More...
 
int getValueType (Field *f)
 Retrieve the type of value a field distributes.
 
void destroyField (Field *f)
 Destroy an apf::Field. More...
 
void setScalar (Field *f, MeshEntity *e, int node, double value)
 Set a nodal value of a scalar field. More...
 
double getScalar (Field *f, MeshEntity *e, int node)
 Get the node value of a scalar field. More...
 
void setVector (Field *f, MeshEntity *e, int node, Vector3 const &value)
 Set the nodal value of a vector field. More...
 
void getVector (Field *f, MeshEntity *e, int node, Vector3 &value)
 Get the nodal value of a vector field.
 
void setMatrix (Field *f, MeshEntity *e, int node, Matrix3x3 const &value)
 Set the nodal value of a matrix field. More...
 
void getMatrix (Field *f, MeshEntity *e, int node, Matrix3x3 &value)
 Get the nodal value of a matrix field.
 
void setComponents (Field *f, MeshEntity *e, int node, double const *components)
 Set the nodal value from an array of component values.
 
void getComponents (Field *f, MeshEntity *e, int node, double *components)
 Copy the nodal value into an array of component values. More...
 
Element * createElement (Field *f, MeshElement *e)
 Create a Field Element from a Mesh Element. More...
 
Element * createElement (Field *f, MeshEntity *e)
 Create a Field Element without a parent Mesh Element. More...
 
void destroyElement (Element *e)
 Destroy a Field Element.
 
MeshElementgetMeshElement (Element *e)
 Get the Mesh Element of a Field Element. More...
 
MeshEntity * getMeshEntity (Element *e)
 Retrieve the mesh entity of an apf::Element.
 
double getScalar (Element *e, Vector3 const &param)
 Evaluate a scalar field at a point. More...
 
void getGrad (Element *e, Vector3 const &param, Vector3 &grad)
 Get the gradient of a scalar field w.r.t. global coordinates. More...
 
void getVector (Element *e, Vector3 const &param, Vector3 &value)
 Evaluate a vector field at a point. More...
 
double getDiv (Element *e, Vector3 const &param)
 Evaluate the divergence of a vector field at a point. More...
 
void getCurl (Element *e, Vector3 const &param, Vector3 &curl)
 Evaluate the curl of a vector field at a point. More...
 
void getVectorGrad (Element *e, Vector3 const &param, Matrix3x3 &deriv)
 Get the gradient of a vector field w.r.t. global coordinates. More...
 
void getMatrix (Element *e, Vector3 const &param, Matrix3x3 &value)
 Evaluate the value of a matrix field. More...
 
void getMatrixGrad (Element *e, Vector3 const &param, Vector< 27 > &value)
 get the gradient of a matrix field More...
 
void getComponents (Element *e, Vector3 const &param, double *components)
 Evaluate a field into an array of component values.
 
int countIntPoints (MeshElement *e, int order)
 Get the number of integration points for an element. More...
 
void getIntPoint (MeshElement *e, int order, int point, Vector3 &param)
 Get an integration point in an element. More...
 
double getIntWeight (MeshElement *e, int order, int point)
 Get the weight of an integration point in an element. More...
 
void mapLocalToGlobal (MeshElement *e, Vector3 const &local, Vector3 &global)
 Map a local coordinate to a global coordinate.
 
double getDV (MeshElement *e, Vector3 const &param)
 Get the differential volume at a point. More...
 
double measure (MeshElement *e)
 Measures the volume, area, or length of a Mesh Element. More...
 
double measure (Mesh *m, MeshEntity *e)
 Measures the volume, area, or length of a Mesh Entity. More...
 
int getOrder (MeshElement *e)
 Returns the polynomial order of the coordinate field.
 
void getJacobian (MeshElement *e, Vector3 const &local, Matrix3x3 &j)
 Returns the Jacobian at a local point.
 
void getJacobianInv (MeshElement *e, Vector3 const &local, Matrix3x3 &jinv)
 Returns the Jacobian inverse at a local point. More...
 
double computeCosAngle (Mesh *m, MeshEntity *pe, MeshEntity *e1, MeshEntity *e2, const Matrix3x3 &Q)
 Returns the cosine of the angle between 2 entities of the parent entity. More...
 
double computeShortestHeightInTet (Mesh *m, MeshEntity *tet, const Matrix3x3 &Q=Matrix3x3(1., 0., 0., 0., 1., 0., 0., 0., 1.))
 Returns the shortest height in a tet. More...
 
double computeLargestHeightInTet (Mesh *m, MeshEntity *tet, const Matrix3x3 &Q=Matrix3x3(1., 0., 0., 0., 1., 0., 0., 0., 1.))
 Returns the largest height in a tet. More...
 
double computeShortestHeightInTri (Mesh *m, MeshEntity *tri, const Matrix3x3 &Q=Matrix3x3(1., 0., 0., 0., 1., 0., 0., 0., 1.))
 Returns the shortest height in a tri. More...
 
double computeLargestHeightInTri (Mesh *m, MeshEntity *tri, const Matrix3x3 &Q=Matrix3x3(1., 0., 0., 0., 1., 0., 0., 0., 1.))
 Returns the largest height in a tri. More...
 
int countNodes (Element *e)
 Returns the number of element nodes. More...
 
void getScalarNodes (Element *e, NewArray< double > &values)
 Returns the element nodal values for a scalar field.
 
void getVectorNodes (Element *e, NewArray< Vector3 > &values)
 Returns the element nodal values for a vector field.
 
void getMatrixNodes (Element *e, NewArray< Matrix3x3 > &values)
 Returns the element nodal values for a matrix field.
 
void getShapeValues (Element *e, Vector3 const &local, NewArray< double > &values)
 Returns the shape function values at a point.
 
void getShapeGrads (Element *e, Vector3 const &local, NewArray< Vector3 > &grads)
 Returns the shape function gradients at a point. More...
 
void getVectorShapeValues (Element *e, Vector3 const &local, NewArray< Vector3 > &values)
 Returns the vector shape function values at a point. More...
 
void getCurlShapeValues (Element *e, Vector3 const &local, NewArray< Vector3 > &values)
 Returns the vector curl shape function values at a point. More...
 
FieldShapegetShape (Field *f)
 Retrieve the apf::FieldShape used by a field.
 
int countComponents (Field *f)
 Count the number of scalar components in the field's value type.
 
void writeVtkFiles (const char *prefix, Mesh *m, int cellDim=-1)
 Write a set of parallel VTK Unstructured Mesh files from an apf::Mesh with binary (base64) encoding and zlib compression (if LION_COMPRESS=ON) More...
 
void writeVtkFiles (const char *prefix, Mesh *m, std::vector< std::string > writeFields, int cellDim=-1)
 Write a set of parallel VTK Unstructured Mesh files from an apf::Mesh with binary (base64) encoding and zlib compression (if LION_COMPRESS=ON) More...
 
void writeOneVtkFile (const char *prefix, Mesh *m)
 Output just the .vtu file with ASCII encoding for this part. More...
 
void writeASCIIVtkFiles (const char *prefix, Mesh *m)
 Write a set of parallel VTK Unstructured Mesh files from an apf::Mesh with ASCII encoding. More...
 
void writeASCIIVtkFiles (const char *prefix, Mesh *m, std::vector< std::string > writeFields)
 Write a set of parallel VTK Unstructured Mesh files from an apf::Mesh with ASCII encoding. More...
 
void writeNedelecVtkFiles (const char *prefix, Mesh *m)
 Output .vtk files with ASCII encoding for this part. More...
 
void getGaussPoint (int type, int order, int point, Vector3 &param)
 Return the location of a gaussian integration point. More...
 
int countGaussPoints (int type, int order)
 Return the number of Gaussian integration points.
 
double getJacobianDeterminant (Matrix3x3 const &J, int dimension)
 Return the Jacobian determinant or dimensional equivalent. More...
 
int getDimension (MeshElement *me)
 Return the dimension of a MeshElement's MeshEntity.
 
void synchronize (Field *f, Sharing *shr=0)
 Synchronize field values along partition boundary. More...
 
void accumulate (Field *f, Sharing *shr=0, bool delete_shr=false)
 Add field values along partition boundary. More...
 
void sharedReduction (Field *f, Sharing *shr, bool delete_shr, const ReductionOp< double > &sum=ReductionSum< double >())
 Apply a reudction operator along partition boundaries. More...
 
bool isPrintable (Field *f)
 Checks whether a Field/Numbering/GlobalNumbering is complete and therefore printable to visualization files. This is a collective operation.
 
void fail (const char *why) __attribute__((noreturn))
 Declare failure of code inside APF. More...
 
void freeze (Field *f)
 Convert a Field from Tag to array storage.
 
void unfreeze (Field *f)
 Convert a Field from array to Tag storage.
 
bool isFrozen (Field *f)
 Returns true iff the Field uses array storage.
 
double * getArrayData (Field *f)
 Return the contiguous array storing this field. More...
 
void zeroField (Field *f)
 Initialize all nodal values with all-zero components.
 
Field * recoverGradientByVolume (Field *f)
 Compute a nodal gradient field from a nodal input field. More...
 
void projectField (Field *to, Field *from)
 Project a field from an existing field.
 
void getBF (FieldShape *s, MeshElement *e, Vector3 const &p, NewArray< double > &BF)
 Get the basis functions over a mesh element. More...
 
void getGradBF (FieldShape *s, MeshElement *e, Vector3 const &p, NewArray< Vector3 > &gradBF)
 Get global gradients of basis functions over a mesh element. More...
 
void verify (Mesh *m, bool abort_on_error=true)
 run consistency checks on an apf::Mesh structure More...
 
int getDimension (Mesh *m, MeshEntity *e)
 get the dimension of a mesh entity
 
void unite (Parts &into, Parts const &from)
 unite two sets of unique part ids More...
 
void removeTagFromDimension (Mesh *m, MeshTag *tag, int d)
 removes a tag from all entities of dimension (d)
 
MeshEntity * findUpward (Mesh *m, int type, MeshEntity **down)
 find an entity from one-level downward adjacencies More...
 
MeshEntity * findElement (Mesh *m, int type, MeshEntity **verts)
 finds an entity from a set of vertices
 
MeshEntity * getEdgeVertOppositeVert (Mesh *m, MeshEntity *edge, MeshEntity *v)
 get the other vertex of an edge
 
void getBridgeAdjacent (Mesh *m, MeshEntity *origin, int bridgeDimension, int targetDimension, Adjacent &result)
 get 2nd-order adjacent entities
 
int countEntitiesOfType (Mesh *m, int type)
 count all on-part entities of one topological type
 
bool isSimplex (int type)
 return true if the topological type is a simplex
 
Vector3 getLinearCentroid (Mesh *m, MeshEntity *e)
 get the average of the entity's vertex coordinates More...
 
SharinggetSharing (Mesh *m)
 create a default sharing object for this mesh More...
 
void getPeers (Mesh *m, int d, Parts &peers)
 scan the part for [vtx|edge|face]-adjacent part ids
 
int findIn (MeshEntity **a, int n, MeshEntity *e)
 find pointer (e) in array (a) of length (n) More...
 
void findTriDown (Mesh *m, MeshEntity **verts, MeshEntity **down)
 given the vertices of a triangle, find its edges More...
 
void changeMeshShape (Mesh *m, FieldShape *newShape, bool project=true)
 deprecated wrapper for apf::Mesh::changeShape
 
void unfreezeFields (Mesh *m)
 unfreeze all associated fields More...
 
void freezeFields (Mesh *m)
 freeze all associated fields More...
 
int countEntitiesOn (Mesh *m, ModelEntity *me, int dim)
 count the number of mesh entities classified on a model entity
 
int countOwned (Mesh *m, int dim, Sharing *shr=NULL)
 count the number of owned entities of dimension (dim) using sharing shr the default sharing is used if none is provided
 
void printStats (Mesh *m)
 print global mesh entity counts per dimension
 
void warnAboutEmptyParts (Mesh *m)
 print to stderr the number of empty parts, if any
 
Copy getOtherCopy (Mesh *m, MeshEntity *s)
 given a mesh face, return its remote copy
 
int getFirstType (Mesh *m, int dim)
 get the type of the first entity in this dimension
 
void getAlignment (Mesh *m, MeshEntity *elem, MeshEntity *boundary, int &which, bool &flip, int &rotate)
 boundary entity alignment to an element More...
 
MeshTag * tagOpposites (GlobalNumbering *gn, const char *name)
 Tag boundary faces with global ids of opposite elements. More...
 
void migrate (Mesh2 *m, Migration *plan)
 APF's migration function, works on apf::Mesh2. More...
 
void setMigrationLimit (size_t maxElements, pcu::PCU *PCUObj)
 set the maximum elements that apf::migrate moves at once More...
 
void displaceMesh (Mesh2 *m, Field *d, double factor=1.0)
 add a field (times a factor) to the mesh coordinates More...
 
MeshEntity * makeOrFind (Mesh2 *m, ModelEntity *c, int type, MeshEntity **down, BuildCallback *cb=0, bool *p_made=0)
 like apf::Mesh2::createEntity, but returns already existing entities
 
MeshEntity * buildElement (Mesh2 *m, ModelEntity *c, int type, MeshEntity **verts, BuildCallback *cb=0)
 build an entity from its vertices More...
 
MeshEntity * buildOneElement (Mesh2 *m, ModelEntity *c, int type, Vector3 const *points)
 build a one-element mesh More...
 
void initResidence (Mesh2 *m, int dim)
 Set entity residence based on remote copies. More...
 
void stitchMesh (Mesh2 *m)
 infer all remote copies from those of vertices More...
 
void clear (Mesh2 *m)
 removes all entities and fields.
 
template<std::size_t M, std::size_t N>
Matrix< N, M > transpose (Matrix< M, N > const &m)
 transpose a matrix
 
template<std::size_t M, std::size_t N>
Matrix< M, N > tensorProduct (Vector< M > const &a, Vector< N > const &b)
 tensor product of two vectors
 
template<std::size_t M, std::size_t N>
Matrix< M-1, N-1 > getMinor (Matrix< M, N > const &A, std::size_t i, std::size_t j)
 get the minor matrix associated with entry (i,j) of matrix A More...
 
template<std::size_t M, std::size_t N>
double getCofactor (Matrix< M, N > const &A, std::size_t i, std::size_t j)
 get the cofactor associated with entry (i,j) of matrix A More...
 
template<std::size_t M, std::size_t N>
double getDeterminant (Matrix< M, N > const &A)
 get the determinant of a matrix A More...
 
Matrix< 3, 3 > cofactor (Matrix< 3, 3 > const &m)
 get the matrix of cofactors for a given matrix
 
Matrix< 2, 2 > invert (Matrix< 2, 2 > const &m)
 invert a 2 by 2 matrix
 
Matrix< 3, 3 > invert (Matrix< 3, 3 > const &m)
 invert a 3 by 3 matrix
 
template<std::size_t M, std::size_t N>
double getInnerProduct (Matrix< M, N > const &a, Matrix< M, N > const &b)
 get the component-wise inner product of two matrices
 
Matrix3x3 cross (Vector3 const &u)
 get the skew-symmetric cross product matrix of a vector
 
Matrix3x3 rotate (Vector3 const &u, double a)
 get the rotation matrix around an axis More...
 
Matrix3x3 getFrame (Vector3 const &v)
 derive an orthogonal frame whose x axis is the given vector More...
 
int eigen (Matrix3x3 const &A, Vector< 3 > *eigenVectors, double *eigenValues)
 get the eigenvectors and eigenvalues of a 3 by 3 matrix
 
template<std::size_t M>
void applyMatrixFunc (Matrix< M, M > const &A, double(*callback)(double), Matrix< M, M > &newMat)
 
Vector< 3 > cross (Vector< 3 > const &a, Vector< 3 > const &b)
 3D vector cross product
 
template<std::size_t N>
Vector< N > project (Vector< N > const &a, Vector< N > const &b)
 Returns vector (a) projected onto vector (b)
 
template<std::size_t N>
Vector< N > reject (Vector< N > const &a, Vector< N > const &b)
 vector rejection
 
void multiply (DynamicMatrix const &a, DynamicVector const &b, DynamicVector &r)
 multiply a DynamicMatrix by a DynamicVector
 
void multiply (DynamicVector const &b, DynamicMatrix const &a, DynamicVector &r)
 multiply a DynamicVector by a DynamicMatrix
 
void multiply (DynamicMatrix const &a, DynamicMatrix const &b, DynamicMatrix &r)
 multiply two DynamicMatrix objects
 
void transpose (DynamicMatrix const &a, DynamicMatrix &r)
 get the transpose of a DynamicMatrix
 
template<std::size_t N, std::size_t M>
DynamicMatrix fromMatrix (Matrix< N, M > other)
 convert an apf::Matrix into an apf::DynamicMatrix
 
template<std::size_t N>
DynamicVector fromVector (Vector< N > other)
 convert an apf::Vector into an apf::DynamicVector
 
FieldShapegetLagrange (int order)
 Get the Lagrangian shape function of some polynomial order. More...
 
FieldShapegetSerendipity ()
 Get the Serendipity shape functions of second order.
 
FieldShapegetConstant (int dimension)
 Get the Constant shape function over some dimension. More...
 
FieldShapegetIPShape (int dimension, int order)
 Get the Integration Point distribution. More...
 
FieldShapegetVoronoiShape (int dimension, int order)
 Get the Voronoi shape function. More...
 
FieldShapegetIPFitShape (int dimension, int order)
 Get the IP Fit shape function. More...
 
FieldShapegetHierarchic (int order)
 Get the quadratic hierarchic shape function. More...
 
FieldShapegetNedelec (int order)
 Get the Nedelec shape function of a given order. More...
 
FieldShapegetL2Shape (int order, int type)
 Get the L2 shapes of a given order and entity type.
 
FieldShapegetH1Shape (int order)
 Get the H1 shapes of a given order. More...
 
void projectHierarchicField (Field *to, Field *from)
 Project a hierarchic field.
 
void projectNedelecField (Field *to, Field *from)
 Project a Nedelec field.
 
void projectL2Field (Field *to, Field *from)
 Project a L2 field.
 
int countElementNodes (FieldShape *s, int type)
 count the number of nodes affecting an element type More...
 
void getElementNodeXis (FieldShape *s, int type, apf::NewArray< apf::Vector3 > &xis)
 gets the xi coordinates for all the nodes More...
 
void getElementNodeXis (FieldShape *s, Mesh *m, MeshEntity *e, apf::NewArray< apf::Vector3 > &xis)
 gets the xi coordinates for all the nodes More...
 
Vector3 boundaryToElementXi (Mesh *m, MeshEntity *boundary, MeshEntity *element, Vector3 const &xi)
 Reparameterize from boundary entity to element. More...
 
NumberingcreateNumbering (Field *f)
 Create a Numbering of degrees of freedom of a Field.
 
NumberingcreateNumbering (Mesh *mesh, const char *name, FieldShape *shape, int components)
 Create a generally-defined Numbering. More...
 
void destroyNumbering (Numbering *n)
 Destroy a Numbering.
 
void fix (Numbering *n, MeshEntity *e, int node, int component, bool fixed)
 Set the fixed/free status of a degree of freedom,. More...
 
bool isFixed (Numbering *n, MeshEntity *e, int node, int component)
 Check whether a degree of freedom is fixed.
 
bool isNumbered (Numbering *n, MeshEntity *e, int node, int component)
 Check whether a degree of freedom is numbered.
 
void number (Numbering *n, MeshEntity *e, int node, int component, int number)
 number a degree of freedom
 
int getNumber (Numbering *n, MeshEntity *e, int node, int component)
 get a degree of freedom number
 
Field * getField (Numbering *n)
 get the field being numbered
 
FieldShapegetShape (Numbering *n)
 get the FieldShape used by a Numbering
 
const char * getName (Numbering *n)
 get the name of a Numbering
 
MeshgetMesh (Numbering *n)
 get the mesh associated with a Numbering
 
int getElementNumbers (Numbering *n, MeshEntity *e, NewArray< int > &numbers)
 returns the node numbers of an element More...
 
int countFixed (Numbering *n)
 return the number of fixed degrees of freedom
 
void synchronize (Numbering *n, Sharing *shr=0, bool delete_shr=false)
 numbers non-owned nodes with the values from their owners More...
 
NumberingnumberOwnedDimension (Mesh *mesh, const char *name, int dim, Sharing *shr=0)
 number the local owned entities of a given dimension
 
NumberingnumberOverlapDimension (Mesh *mesh, const char *name, int dim)
 number all local entities of a given dimension
 
NumberingnumberElements (Mesh *mesh, const char *name)
 number the local elements
 
NumberingnumberOverlapNodes (Mesh *mesh, const char *name, FieldShape *s=0)
 number all local nodes More...
 
NumberingnumberOwnedNodes (Mesh *mesh, const char *name, FieldShape *s=0, Sharing *shr=0)
 number the local owned nodes More...
 
int countNodes (Numbering *n)
 count the number of nodes that have been numbered
 
void getNodes (Numbering *n, DynamicArray< Node > &nodes)
 get an array of numbered nodes More...
 
void getNodesOnClosure (Mesh *m, ModelEntity *me, DynamicArray< Node > &on, FieldShape *sh=0)
 get nodes on the closure of a model entity More...
 
GlobalNumberingcreateGlobalNumbering (Mesh *mesh, const char *name, FieldShape *shape, int components=1)
 create global numbering More...
 
GlobalNumberingcreateGlobalNumbering (Field *f)
 Create a Numbering of degrees of freedom of a Field.
 
MeshgetMesh (GlobalNumbering *n)
 get the mesh associated with a global numbering
 
int countComponents (GlobalNumbering *n)
 get the components associated with a global numbering
 
void number (GlobalNumbering *n, MeshEntity *e, int node, long number)
 assign a global number
 
void number (GlobalNumbering *n, Node node, long number, int component=0)
 assign a global number
 
long getNumber (GlobalNumbering *n, Node node, int component=0)
 get a global number
 
long getNumber (GlobalNumbering *n, MeshEntity *e, int node, int component=0)
 get a global number
 
int getElementNumbers (GlobalNumbering *n, MeshEntity *e, NewArray< long > &numbers)
 get an element's global node numbers
 
Field * getField (GlobalNumbering *n)
 get the field being numbered
 
GlobalNumberingmakeGlobal (Numbering *n, bool destroy=true)
 converts a local numbering into a global numbering. More...
 
void synchronize (GlobalNumbering *n, Sharing *shr=0)
 see the Numbering equivalent and apf::makeGlobal
 
void destroyGlobalNumbering (GlobalNumbering *n)
 destroy a global numbering
 
void getNodes (GlobalNumbering *n, DynamicArray< Node > &nodes)
 see the Numbering equivalent
 
MeshTag * reorder (Mesh *mesh, const char *name)
 Number by adjacency graph traversal. More...
 
int naiveOrder (Numbering *num, Sharing *sharing=NULL)
 number all components by simple iteration
 
int NaiveOrder (Numbering *num)
 todo : mark as deprecated
 
int adjReorder (Numbering *num, Sharing *sharing=NULL)
 like apf::reorder, but numbers all free nodal components
 
int AdjReorder (Numbering *num)
 todo : mark as deprecated
 
void setNumberingOffset (Numbering *num, int off, Sharing *sharing=NULL)
 add an offset to all free nodal component numbers
 
void SetNumberingOffset (Numbering *num, int off)
 todo : mark as deprecated
 
int countDOFs (std::vector< Numbering * > const &n)
 Count the total numbered degrees of freedom. More...
 
int countDOFs (std::vector< GlobalNumbering * > const &n)
 Count the total numbered degrees of freedom. More...
 
void getElementNumbers (std::vector< Numbering * > const &n, MeshEntity *e, std::vector< int > &numbers)
 Get the element numbers for multiple numberings. More...
 
void getElementNumbers (std::vector< GlobalNumbering * > const &n, MeshEntity *e, std::vector< long > &numbers)
 Get the element numbers for multiple global numberings. More...
 
int numberOwned (std::vector< Field * > const &fields, std::vector< Numbering * > &owned)
 Number the owned nodes of multiple fields. More...
 
int numberGhost (std::vector< Field * > const &fields, std::vector< Numbering * > &ghost)
 Number the ghost (overlapped/shared) nodes of multiple fields. More...
 
void makeGlobal (std::vector< Numbering * > &owned, std::vector< GlobalNumbering * > &global, pcu::PCU *PCUObj)
 Globalize a mixed numbering scheme. More...
 
void remapPartition (apf::Mesh2 *m, Remap &remap)
 remap all part ids in the mesh structure More...
 
void convert (Mesh *in, Mesh2 *out, MeshEntity **nodes=NULL, MeshEntity **elems=NULL, bool copy_data=true)
 convert one mesh data structure to another More...
 
NewElements assemble (Mesh2 *m, const apf::Gid *conn, int nelem, int etype, GlobalToVert &globalToVert)
 assemble a mixed-cell-type mesh from just a connectivity array More...
 
void finalise (Mesh2 *m, GlobalToVert &globalToVert)
 finalise construction of a mixed-cell-type mesh from just a connectivity array More...
 
NewElements construct (Mesh2 *m, const apf::Gid *conn, int nelem, int etype, GlobalToVert &globalToVert)
 construct a mesh from just a connectivity array More...
 
void setCoords (Mesh2 *m, const double *coords, int nverts, GlobalToVert &globalToVert)
 Assign coordinates to the mesh. More...
 
void setMatches (Mesh2 *m, const Gid *matches, int nverts, GlobalToVert &globalToVert)
 Assign matching to the mesh. More...
 
void destruct (Mesh2 *m, Gid *&conn, int &nelem, int &etype, int cellDim=-1)
 convert an apf::Mesh2 object into a connectivity array More...
 
void extractCoords (Mesh2 *m, double *&coords, int &nverts)
 get a contiguous set of global vertex coordinates More...
 
Mesh2makeEmptyMdsMesh (gmi_model *model, int dim, bool isMatched, pcu::PCU *PCUObj)
 create an empty MDS part More...
 
Mesh2loadMdsMesh (const char *modelfile, const char *meshfile, pcu::PCU *PCUObj)
 load an MDS mesh and model from file More...
 
Mesh2loadMdsMesh (gmi_model *model, const char *meshfile, pcu::PCU *PCUObj)
 load an MDS mesh from files More...
 
Mesh2createMdsMesh (gmi_model *model, Mesh *from, bool reorder=false, bool copy_data=true)
 create an MDS mesh from an existing mesh More...
 
void reorderMdsMesh (Mesh2 *mesh, MeshTag *t=0)
 apply adjacency-based reordering More...
 
bool alignMdsMatches (Mesh2 *in)
 align the downward adjacencies of matched entities
 
bool alignMdsRemotes (Mesh2 *in)
 align the downward adjacencies of remote copies
 
void deriveMdsModel (Mesh2 *in)
 build a null model such that apf::verify accepts the mesh. More...
 
void deriveMdlFromManifold (Mesh2 *mesh, bool *isModelVert, int nBFaces, int(*bFaces)[5], GlobalToVert &globalToVert, std::map< int, apf::MeshEntity * > &globalToRegion)
 Given the mesh vertices that are also model vertices, and the classification on boundary mesh faces, constructs the classification on the rest of the boundary entities. More...
 
void derive2DMdlFromManifold (Mesh2 *mesh, bool *isModelVert, int nBEdges, int(*bEdges)[4], GlobalToVert &globalToVert, std::map< int, apf::MeshEntity * > &globalToFace)
 Given the mesh vertices that are also model vertices, and the classification on boundary mesh edges, constructs the classification on the rest of the boundary entities for a 2-D mesh. More...
 
void changeMdsDimension (Mesh2 *in, int d)
 change the dimension of an MDS mesh More...
 
int getMdsIndex (Mesh2 *in, MeshEntity *e)
 returns the dimension-unique index for this entity More...
 
MeshEntity * getMdsEntity (Mesh2 *in, int dimension, int index)
 retrieve an entity by dimension and index More...
 
Mesh2loadMdsFromANSYS (const char *nodefile, const char *elemfile, pcu::PCU *PCUObj)
 load an MDS mesh from ANSYS .node and .elem files More...
 
Mesh2makeMdsBox (int nx, int ny, int nz, double wx, double wy, double wz, bool is, pcu::PCU *PCUObj)
 create a box from an MDS mesh More...
 
gmi_modelmakeMdsBoxModel (int nx, int ny, int nz, double wx, double wy, double wz, bool is, pcu::PCU *PCUObj)
 see makeMdsBox - only creates geometric model
 
SplittermakeZoltanSplitter (Mesh *mesh, int method, int approach, bool debug=false, bool sync=true)
 Make a Zoltan Splitter object. More...
 
SplittermakeZoltanGlobalSplitter (Mesh *mesh, int method, int approach, bool debug=false)
 Make a Zoltan Splitter object. More...
 
BalancermakeZoltanBalancer (Mesh *mesh, int method, int approach, bool debug=false)
 Make a Zoltan Balancer object. More...
 
int * getElementToElement (apf::Mesh *m)
 Get an element-to-element connectivity array. More...
 
SplittermakeMETISsplitter (Mesh *mesh)
 Make an apf::Splitter that calls METIS for the underlying algorithm. More...
 
BalancermakeMETISbalancer (Mesh *mesh)
 Make an apf::Balancer that calls METIS for the underlying algorithm. More...
 
bool hasMETIS ()
 Query whether METIS is supported. More...
 

Variables

int const tri_edge_verts [3][2]
 map from triangle edge order to triangle vertex order
 
int const quad_edge_verts [4][2]
 map from quad edge order to quad vertex order
 
int const tet_edge_verts [6][2]
 map from tet edge order to tet vertex order
 
int const prism_edge_verts [9][2]
 map from prism edge order to prism vertex order
 
int const pyramid_edge_verts [8][2]
 map from pyramid edge order to pyramid vertex order
 
int const tet_tri_verts [4][3]
 map from tet triangle order to tet vertex order
 
int const hex_quad_verts [6][4]
 map from hex quad order to hex vertex order
 
int const prism_tri_verts [2][3]
 map from prism triangle order to prism vertex order
 
int const prism_quad_verts [3][4]
 map from prism quad order to prism vertex order
 
int const pyramid_tri_verts [4][3]
 map from pyramid triangle order to pyramid vertex order
 
double const pi
 The mathematical constant pi. More...
 

Detailed Description

All APF symbols are contained in this namespace.

Wrapping a namespace over everything gives reasonable insurance against future symbol conflicts with other packages while very common names are being used, like Mesh or VECTOR. Users will be able to use the simple names directly with a using statement, likewise all APF code is written without apf::

Typedef Documentation

◆ Adjacent

typedef DynamicArray<MeshEntity*> apf::Adjacent

Set of adjacent mesh entities.

see also apf::Downward and apf::Up

Definition at line 47 of file apfMesh.h.

◆ Copies

typedef std::map<int,MeshEntity*> apf::Copies

Remote copy container.

the key is the part id, the value is the on-part pointer to the remote copy

Definition at line 42 of file apfMesh.h.

◆ Downward

typedef MeshEntity* apf::Downward[12]

a static array type downward adjacency queries.

using statically sized arrays saves time by avoiding dynamic allocation, and downward adjacencies have a guaranteed bound.

Definition at line 53 of file apfMesh.h.

Enumeration Type Documentation

◆ ValueType

The type of value the field stores.

The near future may bring more complex tensors.

Enumerator
SCALAR 

a single scalar value.

VECTOR 

a 3D vector value

MATRIX 

a 3x3 matrix

PACKED 

a user-defined set of components

VALUE_TYPES 

placeholder used to set array sizes

Definition at line 133 of file apf.h.

133  {
135  SCALAR,
137  VECTOR,
139  MATRIX,
141  PACKED,
144 };
@ PACKED
a user-defined set of components
Definition: apf.h:141
@ VECTOR
a 3D vector value
Definition: apf.h:137
@ SCALAR
a single scalar value.
Definition: apf.h:135
@ VALUE_TYPES
placeholder used to set array sizes
Definition: apf.h:143
@ MATRIX
a 3x3 matrix
Definition: apf.h:139

◆ ZoltanApproach

Zoltan partitioning approach.

Enumerator
PARTITION 

(Hyper)Graph - does not consider the initial distribution

REPARTITION 

(Hyper)Graph - considers the initial distribution

REFINE 

(HYPER)Graph - targets partitions needing only small changes

PART_KWAY 

Graph - multilevel.

PART_GEOM 

Graph - space filling curves.

PART_GEOM_KWAY 

Graph - hybrid method combining PART_KWAY and PART_GEOM.

ADAPT_REPART 

Graph - targets graphs generated from adaptively refined meshes.

REFINE_KWAY 

Graph - targets partitions needing only small changes.

Definition at line 50 of file apfZoltan.h.

50  {
52  PARTITION,
56  REFINE,
58  PART_KWAY,
60  PART_GEOM,
67 };
@ PART_GEOM_KWAY
Graph - hybrid method combining PART_KWAY and PART_GEOM.
Definition: apfZoltan.h:62
@ PART_GEOM
Graph - space filling curves.
Definition: apfZoltan.h:60
@ REPARTITION
(Hyper)Graph - considers the initial distribution
Definition: apfZoltan.h:54
@ REFINE_KWAY
Graph - targets partitions needing only small changes.
Definition: apfZoltan.h:66
@ PARTITION
(Hyper)Graph - does not consider the initial distribution
Definition: apfZoltan.h:52
@ REFINE
(HYPER)Graph - targets partitions needing only small changes
Definition: apfZoltan.h:56
@ ADAPT_REPART
Graph - targets graphs generated from adaptively refined meshes.
Definition: apfZoltan.h:64
@ PART_KWAY
Graph - multilevel.
Definition: apfZoltan.h:58

◆ ZoltanMethod

Zoltan partitioning method.

Enumerator
RCB 

Recursive Coordinate Bisection.

RIB 

Recursive Inertial Bisection.

HYPERGRAPH 

Hyper-graph partitioning.

PARMETIS 

Use ParMetis.

GRAPH 

General graph partitionig.

Definition at line 36 of file apfZoltan.h.

36  {
38  RCB,
40  RIB,
42  HYPERGRAPH,
44  PARMETIS,
46  GRAPH
47 };
@ RCB
Recursive Coordinate Bisection.
Definition: apfZoltan.h:38
@ RIB
Recursive Inertial Bisection.
Definition: apfZoltan.h:40
@ GRAPH
General graph partitionig.
Definition: apfZoltan.h:46
@ PARMETIS
Use ParMetis.
Definition: apfZoltan.h:44
@ HYPERGRAPH
Hyper-graph partitioning.
Definition: apfZoltan.h:42

Function Documentation

◆ accumulate()

void apf::accumulate ( Field *  f,
Sharing shr = 0,
bool  delete_shr = false 
)

Add field values along partition boundary.

Using the copies described by an apf::Sharing object, add up the field values of all copies of an entity and assign the sum as the value for all copies.

◆ applyMatrixFunc()

template<std::size_t M>
void apf::applyMatrixFunc ( Matrix< M, M > const &  A,
double(*)(double)  callback,
Matrix< M, M > &  newMat 
)

This function applies a function to a symmetric tensor. Using the polar decomposition theorem. One example is taking the sqrt of a tensor

Parameters
Aa square matrix
callbacka function pointer that takes a double and returns a double. This is the function that will be applied to the matrix.
newMata matrix that contains the results of the operation

◆ assemble()

NewElements apf::assemble ( Mesh2 m,
const apf::Gid *  conn,
int  nelem,
int  etype,
GlobalToVert globalToVert 
)

assemble a mixed-cell-type mesh from just a connectivity array

construct is now split into two functions, assemble and finalise. The premise of assemble being that it is called multiple times for a given cell type, across several different cell types in the input mesh.

◆ boundaryToElementXi()

Vector3 apf::boundaryToElementXi ( Mesh m,
MeshEntity *  boundary,
MeshEntity *  element,
Vector3 const &  xi 
)

Reparameterize from boundary entity to element.

This function converts a point in the local parametric space of a boundary mesh entity into the equivalent point in the local parametric space of an adjacent element.

◆ buildElement()

MeshEntity* apf::buildElement ( Mesh2 m,
ModelEntity *  c,
int  type,
MeshEntity **  verts,
BuildCallback cb = 0 
)

build an entity from its vertices

any missing intermediate entities will also be built, and all will be classified to (c). If a non-zero BuildCallback pointer is given, it will be called for each entity created, including intermediate ones.

◆ buildOneElement()

MeshEntity* apf::buildOneElement ( Mesh2 m,
ModelEntity *  c,
int  type,
Vector3 const *  points 
)

build a one-element mesh

this is mostly useful for debugging

Todo:
this doesn't get used much, maybe remove it

◆ changeMdsDimension()

void apf::changeMdsDimension ( Mesh2 in,
int  d 
)

change the dimension of an MDS mesh

this should be called before adding entities of dimension higher than the previous mesh dimension (when building a higher dimensional mesh from a lower one), or after removing all entities of higher dimension (when reducing a high dimensional mesh to a lower one)

◆ cloneField()

Field* apf::cloneField ( Field *  f,
Mesh onto 
)

Declare a copy of a field on another apf::Mesh.

This will just make a Field object with the same properties, but not fill in any data.

◆ computeCosAngle()

double apf::computeCosAngle ( Mesh m,
MeshEntity *  pe,
MeshEntity *  e1,
MeshEntity *  e2,
const Matrix3x3 Q 
)

Returns the cosine of the angle between 2 entities of the parent entity.

Parameters
mmesh
peparent entity
e1first entity
e2second entity
Qmetric

All the angles are computed based on the Jacobian, and therefore, this function will work for curved elements as well. For anisotropic meshes the user can provide Q (he metric value of the parent) to get the angles in the metric space. Otherwise and Identity matrix must be given for Q.

◆ computeLargestHeightInTet()

double apf::computeLargestHeightInTet ( Mesh m,
MeshEntity *  tet,
const Matrix3x3 Q = Matrix3x3(1., 0., 0., 0., 1., 0., 0., 0., 1.) 
)

Returns the largest height in a tet.

Parameters
mmesh
tettet
Qmetric (default Identity)

◆ computeLargestHeightInTri()

double apf::computeLargestHeightInTri ( Mesh m,
MeshEntity *  tri,
const Matrix3x3 Q = Matrix3x3(1., 0., 0., 0., 1., 0., 0., 0., 1.) 
)

Returns the largest height in a tri.

Parameters
mmesh
tritri
Qmetric (default Identity)

◆ computeShortestHeightInTet()

double apf::computeShortestHeightInTet ( Mesh m,
MeshEntity *  tet,
const Matrix3x3 Q = Matrix3x3(1., 0., 0., 0., 1., 0., 0., 0., 1.) 
)

Returns the shortest height in a tet.

Parameters
mmesh
tettet
Qmetric (default Identity)

◆ computeShortestHeightInTri()

double apf::computeShortestHeightInTri ( Mesh m,
MeshEntity *  tri,
const Matrix3x3 Q = Matrix3x3(1., 0., 0., 0., 1., 0., 0., 0., 1.) 
)

Returns the shortest height in a tri.

Parameters
mmesh
tritri
Qmetric (default Identity)

◆ construct()

NewElements apf::construct ( Mesh2 m,
const apf::Gid *  conn,
int  nelem,
int  etype,
GlobalToVert globalToVert 
)

construct a mesh from just a connectivity array

this function is here to interface with very simple mesh formats. Given a set of elements described only in terms of the ordered global ids of their vertices, this function builds a reasonable apf::Mesh2 structure and as a side effect returns a map from global ids to local vertices. This functions assumes a uniform cell type. Use a combination of assemble and finalise for meshes loaded with mixed cell types.

This is a fully scalable parallel mesh construction algorithm, no processor incurs memory or runtime costs proportional to the global mesh size.

Note that all vertices will have zero coordinates, so it is often good to use apf::setCoords after this.

◆ convert()

void apf::convert ( Mesh in,
Mesh2 out,
MeshEntity **  nodes = NULL,
MeshEntity **  elems = NULL,
bool  copy_data = true 
)

convert one mesh data structure to another

this function will fill in a structure that fully implements apf::Mesh2 by using information from an implementation of apf::Mesh. This is a fully scalable parallel mesh conversion tool.

◆ countDOFs() [1/2]

int apf::countDOFs ( std::vector< GlobalNumbering * > const &  n)

Count the total numbered degrees of freedom.

Parameters
nThe input local numberings.
Returns
The number of on-part numbered degrees of freedom

◆ countDOFs() [2/2]

int apf::countDOFs ( std::vector< Numbering * > const &  n)

Count the total numbered degrees of freedom.

Parameters
nThe input local numberings.
Returns
The number of on-part numbered degrees of freedom

◆ countElementNodes()

int apf::countElementNodes ( FieldShape s,
int  type 
)

count the number of nodes affecting an element type

Parameters
typeselect from apf::Mesh::Type

◆ countIntPoints()

int apf::countIntPoints ( MeshElement e,
int  order 
)

Get the number of integration points for an element.

Parameters
orderthe polynomial order of accuracy desired for the integration (not to be confused with the polynomial order of shape functions).

◆ countNodes()

int apf::countNodes ( Element *  e)

Returns the number of element nodes.

This is the number of nodes affecting an element, as opposed to the nodes tagged to an entity.

◆ createElement() [1/2]

Element* apf::createElement ( Field *  f,
MeshElement e 
)

Create a Field Element from a Mesh Element.

A Field Element object caches elemental data for use in evaluation, mapping, and integration. use destroyElement to free this data.

Parameters
fThe field which the Element will represent
eAn existing MeshElement for the desired entity
Returns
The new field Element

◆ createElement() [2/2]

Element* apf::createElement ( Field *  f,
MeshEntity *  e 
)

Create a Field Element without a parent Mesh Element.

Warning: most users should call the version which takes a MeshElement as input. Only call this function if you know the other one isn't right for you.

◆ createGlobalNumbering()

GlobalNumbering* apf::createGlobalNumbering ( Mesh mesh,
const char *  name,
FieldShape shape,
int  components = 1 
)

create global numbering

see apf::createNumbering. so far global numberings have one component

◆ createIPField()

Field* apf::createIPField ( Mesh m,
const char *  name,
int  valueType,
int  order 
)

Create an apf::Field of integration point data.

Parameters
mthe mesh over which the field is defined
namea unique name for this field
valueTypethe type of field data
orderpolynomial order of accuracy

◆ createLagrangeField()

Field* apf::createLagrangeField ( Mesh m,
const char *  name,
int  valueType,
int  order 
)

Create an apf::Field using a Lagrange distribution.

Parameters
mthe mesh over which the field is defined
namea unique name for this field
valueTypethe type of field data
orderthe polynomial order of the shape functions (so far 1 or 2)

◆ createMdsMesh()

Mesh2* apf::createMdsMesh ( gmi_model model,
Mesh from,
bool  reorder = false,
bool  copy_data = true 
)

create an MDS mesh from an existing mesh

Parameters
fromthe mesh to copy
reorderif true reorder mesh vertices and elements (start from a vertex with minimum Y)
copy_dataif true (default), copy Fields/Numberings/Tags

this function uses apf::convert to copy any apf::Mesh

◆ createMeshElement() [1/2]

MeshElement* apf::createMeshElement ( Field *  c,
MeshEntity *  e 
)

Creates a non-standard Mesh Element over an entity.

Warning: most users should use the above function unless you know this is the right one for you. This function is useful for creating a mesh element over a coordinate field c that may be distinct from the mesh's underlying coordinate field. This may be useful for updated-Lagrangian simulations, where quantities must be computed on reference and updated configurations.

◆ createMeshElement() [2/2]

MeshElement* apf::createMeshElement ( Mesh m,
MeshEntity *  e 
)

Creates a Mesh Element over an entity.

A Mesh Element allows queries to the coordinate field, including mapping, differential and total volume, as well as gauss integration point data. A Mesh Element is also required to build a Field Element.

◆ createNumbering()

Numbering* apf::createNumbering ( Mesh mesh,
const char *  name,
FieldShape shape,
int  components 
)

Create a generally-defined Numbering.

This numbering will be available via mesh->findNumbering, etc. The shape determines where the nodes are, and the component count determines how many integers there are per node.

◆ createPackedField()

Field* apf::createPackedField ( Mesh m,
const char *  name,
int  components,
apf::FieldShape shape = 0 
)

Create a field of N components without a tensor type.

Packed fields are used to interface with applications whose fields are not easily categorized as tensors of order 0,1,2. They contain enough information to interpolate values in an element and such, but some higher-level functionality is left out.

◆ createStepField()

Field* apf::createStepField ( Mesh m,
const char *  name,
int  valueType 
)

Create an apf::Field using a step distribution.

A step-wise distribution is a C-1 continuous field defined by one node at each element, with the field value being constant over the element and discontinuous at element boundaries

Parameters
mthe mesh over which the field is defined
namea unique name for this field
valueTypethe type of field data

◆ derive2DMdlFromManifold()

void apf::derive2DMdlFromManifold ( Mesh2 mesh,
bool *  isModelVert,
int  nBEdges,
int(*)  bEdges[4],
GlobalToVert globalToVert,
std::map< int, apf::MeshEntity * > &  globalToFace 
)

Given the mesh vertices that are also model vertices, and the classification on boundary mesh edges, constructs the classification on the rest of the boundary entities for a 2-D mesh.

Only for triangular mesh with single model region. The tags provided for edge classification are treated as reserved, and all newly generated tags are distinct regardless of dimension. It is assumed that both mesh vertices are indexed from 0 to (n_verts - 1) and mesh faces from 0 to (n_faces -1).

Parameters
meshThe mesh in consideration
isModelVertArray of bools, one per mesh vertex, telling if that vertex is also a model vertex
nBEdgesnumber of boundary faces
bEdges2D Array of size (nBEdges x 4). For each face, the row is [model_edge_tag, adj_face_tag, global_vtx_id_1, global_vtx_id_2]
globalToVertMaps mesh vertex ID to the mesh vertex. Typically output from apf::construct
globalToFaceMaps mesh face ID to the mesh face

◆ deriveMdlFromManifold()

void apf::deriveMdlFromManifold ( Mesh2 mesh,
bool *  isModelVert,
int  nBFaces,
int(*)  bFaces[5],
GlobalToVert globalToVert,
std::map< int, apf::MeshEntity * > &  globalToRegion 
)

Given the mesh vertices that are also model vertices, and the classification on boundary mesh faces, constructs the classification on the rest of the boundary entities.

Only for tetrahedral mesh with single model region. The tags provided for face classification are treated as reserved, and all newly generated tags are distinct regardless of dimension. It is assumed that the mesh was created using apf::construct, which by default assigns a tag 0 to the model region. Due to this, it is advised to provide face tags starting from 1 if uniqueness is desired. It is assumed that both mesh vertices are indexed from 0 to (n_verts - 1) and mesh regions from 0 to (n_regions -1).

Parameters
meshThe mesh in consideration
isModelVertArray of bools, one per mesh vertex, telling if that vertex is also a model vertex
nBFacesnumber of boundary faces
bFaces2D Array of size (n_bfaces x 5). For each face, the row is [model_face_tag, adj_region_tag, global_vtx_id_1, global_vtx_id_2, global_vtx_id_3]
globalToVertMaps mesh vertex ID to the mesh vertex. Typically output from apf::construct
globalToRegionMaps mesh region ID to the mesh region

◆ deriveMdsModel()

void apf::deriveMdsModel ( Mesh2 in)

build a null model such that apf::verify accepts the mesh.

given an MDS mesh that is (wrongly) classified on a null model, this algorithm will classify all interior entities onto a model region and all boundary entities onto a boundary model entity, as defined by mesh upward adjacencies.

◆ destroyField()

void apf::destroyField ( Field *  f)

Destroy an apf::Field.

This function will also remove any field data that this field attached to its Mesh domain.

◆ destroyMesh()

void apf::destroyMesh ( Mesh m)

Destroys an apf::Mesh.

This only destroys the apf::Mesh object, the underlying mesh database is unaffected. Mesh objects are wrappers over mesh databases, and are created by functions provided outside the APF core.

◆ destroyMeshElement()

void apf::destroyMeshElement ( MeshElement e)

Destroys a Mesh Element.

This only destroys the apf::MeshElement object, the underlying mesh entity and field data are unaffected.

◆ destruct()

void apf::destruct ( Mesh2 m,
Gid *&  conn,
int &  nelem,
int &  etype,
int  cellDim = -1 
)

convert an apf::Mesh2 object into a connectivity array

this is useful for debugging the apf::convert function

Parameters
meshthe apf mesh
nelemnumber of elements
etypeapf::Mesh::Type
cellDimdimension of elements (if embedded in a higher dimension manifold)

◆ displaceMesh()

void apf::displaceMesh ( Mesh2 m,
Field *  d,
double  factor = 1.0 
)

add a field (times a factor) to the mesh coordinates

this is useful in mechanical deformation problems to transform the mesh from reference space to deformed space. Setting the factor to -1 will undo the deformation

◆ extractCoords()

void apf::extractCoords ( Mesh2 m,
double *&  coords,
int &  nverts 
)

get a contiguous set of global vertex coordinates

this is used for debugging apf::setCoords

◆ fail()

void apf::fail ( const char *  why)

Declare failure of code inside APF.

This function prints the string as an APF failure to stderr and then calls abort. It can be called from code that is part of the apf namespace, but not outside of that.

◆ finalise()

void apf::finalise ( Mesh2 m,
GlobalToVert globalToVert 
)

finalise construction of a mixed-cell-type mesh from just a connectivity array

construct is now split into two functions, assemble and finalise. Once the mixed cell type mesh is assembled finalise should be called. Doing it this way provides non-breaking changes for current users of construct, which now just calls assemble and finalise.

◆ findIn()

int apf::findIn ( MeshEntity **  a,
int  n,
MeshEntity *  e 
)

find pointer (e) in array (a) of length (n)

Returns
-1 if not found, otherwise i such that a[i] = e

◆ findTriDown()

void apf::findTriDown ( Mesh m,
MeshEntity **  verts,
MeshEntity **  down 
)

given the vertices of a triangle, find its edges

Parameters
downthe resulting array of edges

◆ findUpward()

MeshEntity* apf::findUpward ( Mesh m,
int  type,
MeshEntity **  down 
)

find an entity from one-level downward adjacencies

this function ignores the ordering of adjacent entities

◆ fix()

void apf::fix ( Numbering n,
MeshEntity *  e,
int  node,
int  component,
bool  fixed 
)

Set the fixed/free status of a degree of freedom,.

must be called prior to making any isFixed() calls on the same node/component.

Parameters
nthe numbering object
ethe mesh entity with which the node is associated
nodethe node number withing the mesh entity
componentthe component number within the nodal tensor

◆ freezeFields()

void apf::freezeFields ( Mesh m)

freeze all associated fields

see apf::freezeField

◆ getAlignment()

void apf::getAlignment ( Mesh m,
MeshEntity *  elem,
MeshEntity *  boundary,
int &  which,
bool &  flip,
int &  rotate 
)

boundary entity alignment to an element

Parameters
mthe mesh
elemthe element
boundaryan entity on the boundary of elem
whichindex of (boundary) in getDownward((elem)...)
fliptrue iff orientation of (boundary) is opposite canonical
rotateposition of canonical vertex 0 in boundary vertices, or in boundary vertices reversed if (flip)==true

◆ getArrayData()

double* apf::getArrayData ( Field *  f)

Return the contiguous array storing this field.

This function is only defined for fields which are using array storage, for which apf::isFrozen returns true.

◆ getBF()

void apf::getBF ( FieldShape s,
MeshElement e,
Vector3 const &  p,
NewArray< double > &  BF 
)

Get the basis functions over a mesh element.

Parameters
sthe field shape that defines the basis functions
ethe mesh element over which the basis functions are defined
pthe local coordinates at which the basis functions are evaluated
BFnodal array of basis functions evaluated at p

◆ getCofactor()

template<std::size_t M, std::size_t N>
double apf::getCofactor ( Matrix< M, N > const &  A,
std::size_t  i,
std::size_t  j 
)

get the cofactor associated with entry (i,j) of matrix A

this is only instantiated for square matrices up to 4 by 4

◆ getComponents()

void apf::getComponents ( Field *  f,
MeshEntity *  e,
int  node,
double *  components 
)

Copy the nodal value into an array of component values.

the output array must already be allocated, apf::countComponents can help with this.

◆ getConstant()

FieldShape* apf::getConstant ( int  dimension)

Get the Constant shape function over some dimension.

this pseudo-shape function places a node on every element of the given dimension. Dimensions up to 3 are available

◆ getCurl()

void apf::getCurl ( Element *  e,
Vector3 const &  param,
Vector3 curl 
)

Evaluate the curl of a vector field at a point.

Parameters
paramThe local coordinates in the element.
curlThe curl vector at that point.

◆ getCurlShapeValues()

void apf::getCurlShapeValues ( Element *  e,
Vector3 const &  local,
NewArray< Vector3 > &  values 
)

Returns the vector curl shape function values at a point.

used only for Nedelec shapes (Piola transformation used to map from parent to physical coordinates)

◆ getDeterminant()

template<std::size_t M, std::size_t N>
double apf::getDeterminant ( Matrix< M, N > const &  A)

get the determinant of a matrix A

this is only instantiated for square matrices up to 4 by 4

◆ getDiv()

double apf::getDiv ( Element *  e,
Vector3 const &  param 
)

Evaluate the divergence of a vector field at a point.

Parameters
paramThe local coordinates in the element.
Returns
The divergence at that point.

◆ getDV()

double apf::getDV ( MeshElement e,
Vector3 const &  param 
)

Get the differential volume at a point.

This function is meant to provide the differential measure of an element at a point, based on the Jacobian determinant in the case of regions, and equivalent terms for lower dimensions.

Returns
The differential volume

◆ getElementNodeXis() [1/2]

void apf::getElementNodeXis ( FieldShape s,
int  type,
apf::NewArray< apf::Vector3 > &  xis 
)

gets the xi coordinates for all the nodes

order follows canonical notation. See tables apf::Mesh::tri_edge_verts, apf::Mesh::tet_edge_verts, and apf::Mesh::tet_tri_verts

Parameters
typeselect from apf::Mesh::Type

◆ getElementNodeXis() [2/2]

void apf::getElementNodeXis ( FieldShape s,
Mesh m,
MeshEntity *  e,
apf::NewArray< apf::Vector3 > &  xis 
)

gets the xi coordinates for all the nodes

order follows downward adjacency and global directions for the bounding entities. xi coordinates will be with respect to the entity e

Parameters
typeselect from apf::Mesh::Type

◆ getElementNumbers() [1/3]

int apf::getElementNumbers ( Numbering n,
MeshEntity *  e,
NewArray< int > &  numbers 
)

returns the node numbers of an element

numbers are returned in the standard element node ordering for its shape functions

Returns
the number of element nodes

◆ getElementNumbers() [2/3]

void apf::getElementNumbers ( std::vector< GlobalNumbering * > const &  n,
MeshEntity *  e,
std::vector< long > &  numbers 
)

Get the element numbers for multiple global numberings.

Parameters
nThe input global numberings.
eThe mesh entity for which to get element numbers.
numbersThe output element numbers.

◆ getElementNumbers() [3/3]

void apf::getElementNumbers ( std::vector< Numbering * > const &  n,
MeshEntity *  e,
std::vector< int > &  numbers 
)

Get the element numbers for multiple numberings.

Parameters
nThe input numberings.
eThe mesh entity for which to get element numbers.
numbersthe output element numbers.

◆ getElementToElement()

int* apf::getElementToElement ( apf::Mesh m)

Get an element-to-element connectivity array.

this function assumes the mesh has one element type. the resulting array is created with new int[nelements * nsides]. nsides is the number of faces of an element. entry [i * nsides + j] is the global id of the j'th adjacent element to local element i, which can be -1 for a geometric boundary.

◆ getFrame()

Matrix3x3 apf::getFrame ( Vector3 const &  v)

derive an orthogonal frame whose x axis is the given vector

this is a robust algorithm for choosing an arbitrary frame that is aligned with a given vector while avoiding the numerical issues that arise when the vector is near global axes.

Note, however, that the resulting frame is not normalized. some uses of this function require that the x basis vector be exactly the input, so that is the default behavior. It is the user's responsibility to normalize the result if desired

◆ getGaussPoint()

void apf::getGaussPoint ( int  type,
int  order,
int  point,
Vector3 param 
)

Return the location of a gaussian integration point.

Parameters
typethe element type, from apf::Mesh::getType
orderthe order of the integration rule
pointwhich point of the integration rule
paramthe resulting parent element coordinates

◆ getGrad()

void apf::getGrad ( Element *  e,
Vector3 const &  param,
Vector3 grad 
)

Get the gradient of a scalar field w.r.t. global coordinates.

Parameters
paramThe local coordinates in the element.
gradThe gradient vector at that point.

◆ getGradBF()

void apf::getGradBF ( FieldShape s,
MeshElement e,
Vector3 const &  p,
NewArray< Vector3 > &  gradBF 
)

Get global gradients of basis functions over a mesh element.

Parameters
gradBFnodal array of gradients of basis functions evaluated at p

◆ getH1Shape()

FieldShape* apf::getH1Shape ( int  order)

Get the H1 shapes of a given order.

These are hierarchic shapes that are compatible with MFEM's impl.

◆ getHierarchic()

FieldShape* apf::getHierarchic ( int  order)

Get the quadratic hierarchic shape function.

only first and second order so far

◆ getIntPoint()

void apf::getIntPoint ( MeshElement e,
int  order,
int  point,
Vector3 param 
)

Get an integration point in an element.

Parameters
orderThe polynomial order of accuracy.
pointThe integration point number.
paramThe resulting local coordinate of the integration point.

◆ getIntWeight()

double apf::getIntWeight ( MeshElement e,
int  order,
int  point 
)

Get the weight of an integration point in an element.

All integration point tables in APF are scaled such that the sum of the weights equals the area of of the parent element.

◆ getIPFitShape()

FieldShape* apf::getIPFitShape ( int  dimension,
int  order 
)

Get the IP Fit shape function.

the IP Fit FieldShape is equivalent to the IPShape except that it is capable of evaluating as a shape function whose value at any point in the element is a polynomial fit to the integration point data in that element.

◆ getIPShape()

FieldShape* apf::getIPShape ( int  dimension,
int  order 
)

Get the Integration Point distribution.

Parameters
dimensionThe dimensionality of the elements
orderThe order of accuracy, determines the integration points

This allows users to create a field which has values at the integration points of elements. Orders 1 to 3 for dimension 2 or 3 are available

◆ getJacobianDeterminant()

double apf::getJacobianDeterminant ( Matrix3x3 const &  J,
int  dimension 
)

Return the Jacobian determinant or dimensional equivalent.

this is a special function taking into account the various formulae for differential volume, area, lenght, etc.

Parameters
JJacobian matrix, vector gradient of coordinate field with respect to parent element coordinates
dimensionspacial dimension of the entity being measured

◆ getJacobianInv()

void apf::getJacobianInv ( MeshElement e,
Vector3 const &  local,
Matrix3x3 jinv 
)

Returns the Jacobian inverse at a local point.

computes the Moore-Penrose psuedo-inverse of J. This is needed to handle non-square Jacobians that arise from edges embedded in 2D or higher and faces embeded in 3D. If J is invertible, the Moore-Penrose inverse equals the inverse.

◆ getLagrange()

FieldShape* apf::getLagrange ( int  order)

Get the Lagrangian shape function of some polynomial order.

we have only first and second order so far

◆ getLinearCentroid()

Vector3 apf::getLinearCentroid ( Mesh m,
MeshEntity *  e 
)

get the average of the entity's vertex coordinates

this also works if given just a vertex, so its a convenient way to get the centroid of any entity

◆ getMatrix()

void apf::getMatrix ( Element *  e,
Vector3 const &  param,
Matrix3x3 value 
)

Evaluate the value of a matrix field.

Parameters
paramThe local coordinates in the element.
valueThe field value at that point.

◆ getMatrixGrad()

void apf::getMatrixGrad ( Element *  e,
Vector3 const &  param,
Vector< 27 > &  value 
)

get the gradient of a matrix field

Parameters
paramThe local coordinate in the element.
valuethe gradient of the field at a point. If the matrix is defined by A_{ij}, the gradient $\frac{\partial A_{ij}}{\partial X_k}=value[i*3+j+9*k]$ .

◆ getMdsEntity()

MeshEntity* apf::getMdsEntity ( Mesh2 in,
int  dimension,
int  index 
)

retrieve an entity by dimension and index

indices follow iteration order, so this function is equivalent to iterating (index) times, but is actually much faster than that. this function only works when the arrays have no gaps, so call apf::reorderMdsMesh after any mesh modification.

◆ getMdsIndex()

int apf::getMdsIndex ( Mesh2 in,
MeshEntity *  e 
)

returns the dimension-unique index for this entity

this function only works when the arrays have no gaps, so call apf::reorderMdsMesh after any mesh modification.

◆ getMeshElement()

MeshElement* apf::getMeshElement ( Element *  e)

Get the Mesh Element of a Field Element.

Each apf::Element operates over an apf::MeshElement, which must be maintained as long as the apf::Element exists. Multiple apf::Elements may share an apf::MeshElement.

◆ getMinor()

template<std::size_t M, std::size_t N>
Matrix<M-1,N-1> apf::getMinor ( Matrix< M, N > const &  A,
std::size_t  i,
std::size_t  j 
)

get the minor matrix associated with entry (i,j) of matrix A

this is only instantiated for square matrices up to 4 by 4

◆ getName()

const char* apf::getName ( Field *  f)

Get the name of a Field.

Both for use convenience and for technical reasons related to tagging, each Field should have a unique name.

◆ getNedelec()

FieldShape* apf::getNedelec ( int  order)

Get the Nedelec shape function of a given order.

TODO: complete later

◆ getNodes()

void apf::getNodes ( Numbering n,
DynamicArray< Node > &  nodes 
)

get an array of numbered nodes

the array is sorted by dimension first and then by mesh iterator order

◆ getNodesOnClosure()

void apf::getNodesOnClosure ( Mesh m,
ModelEntity *  me,
DynamicArray< Node > &  on,
FieldShape sh = 0 
)

get nodes on the closure of a model entity

all local nodes associated with mesh entities classified on the closure of the model entity are returned. this is useful for applying boundary conditions, and correctly handles cases when a part has, for example, vertices classified on the edge of a model face, but none of the triangles classified on the model face itself.

◆ getScalar() [1/2]

double apf::getScalar ( Element *  e,
Vector3 const &  param 
)

Evaluate a scalar field at a point.

Parameters
paramThe local coordinates in the element.
Returns
The field value at that point.

◆ getScalar() [2/2]

double apf::getScalar ( Field *  f,
MeshEntity *  e,
int  node 
)

Get the node value of a scalar field.

Returns
The field value at that node.

◆ getShapeGrads()

void apf::getShapeGrads ( Element *  e,
Vector3 const &  local,
NewArray< Vector3 > &  grads 
)

Returns the shape function gradients at a point.

these are gradients with respect to global coordinates.

◆ getSharing()

Sharing* apf::getSharing ( Mesh m)

create a default sharing object for this mesh

for normal meshes, the sharing object just describes remote copies. For matched meshes, the sharing object describes matches for matched entities and remote copies for other entities

◆ getVector()

void apf::getVector ( Element *  e,
Vector3 const &  param,
Vector3 value 
)

Evaluate a vector field at a point.

Parameters
paramThe local coordinates in the element.
valueThe field value at that point.

◆ getVectorGrad()

void apf::getVectorGrad ( Element *  e,
Vector3 const &  param,
Matrix3x3 deriv 
)

Get the gradient of a vector field w.r.t. global coordinates.

Parameters
paramThe local coordinates in the element.
derivThe gradient matrix at that point.

Note: The return parameter component deriv[j][i] stores the value $(grad u)_{i,j} = \frac{\partial u_i} / \frac{\partial x_j}$, where $u$ is the vector field of interest, i is the vector index, and j is the derivative index.

◆ getVectorShapeValues()

void apf::getVectorShapeValues ( Element *  e,
Vector3 const &  local,
NewArray< Vector3 > &  values 
)

Returns the vector shape function values at a point.

used only for Nedelec shapes (Piola transformation used to map from parent to physical coordinates)

◆ getVoronoiShape()

FieldShape* apf::getVoronoiShape ( int  dimension,
int  order 
)

Get the Voronoi shape function.

the Voronoi FieldShape is equivalent to the IPShape except that it is capable of evaluating as a shape function whose value at any point in the element is the value of the closest integration point in that element.

◆ initResidence()

void apf::initResidence ( Mesh2 m,
int  dim 
)

Set entity residence based on remote copies.

this function acts on all entities of one dimension

◆ loadMdsFromANSYS()

Mesh2* apf::loadMdsFromANSYS ( const char *  nodefile,
const char *  elemfile,
pcu::PCU PCUObj 
)

load an MDS mesh from ANSYS .node and .elem files

this call takes two filenames, one for a .node and another for a .elem file.

the resulting MDS mesh will be constructed with a null geometric model via gmi_load(".null"), so be sure to call gmi_register_null before this function.

currently, ANSYS element types SOLID72 and SOLID92 are supported, which become linear and quadratic tetrahedra, respectively.

◆ loadMdsMesh() [1/2]

Mesh2* apf::loadMdsMesh ( const char *  modelfile,
const char *  meshfile,
pcu::PCU PCUObj 
)

load an MDS mesh and model from file

Parameters
modelfilewill be passed to gmi_load to get the model
Note
gmi_register_mesh and gmi_register_null need to be called before this function. Also, gmi_register_sim may also be called to enable loading of GeomSim, Parasolid, and ACIS models

◆ loadMdsMesh() [2/2]

Mesh2* apf::loadMdsMesh ( gmi_model model,
const char *  meshfile,
pcu::PCU PCUObj 
)

load an MDS mesh from files

Parameters
modelthe geometric model interface
meshfileThe path to an SMB format mesh. If the path is "something"".smb", then the file "somethingN.smb" will be loaded where N is the part number. If the path is "something/", then the file "something/N.smb" will be loaded. For both of these cases, if the path is prepended with "bz2:", then it will be uncompressed using PCU file IO functions. Calling apf::Mesh::writeNative on the resulting object will do the same in reverse.

◆ makeEmptyMdsMesh()

Mesh2* apf::makeEmptyMdsMesh ( gmi_model model,
int  dim,
bool  isMatched,
pcu::PCU PCUObj 
)

create an empty MDS part

Parameters
modelthe geometric model interface
dimthe eventual mesh dimension. MDS needs to allocate arrays based on this before users add entities.
isMatchedwhether or not there will be matched entities

◆ makeGlobal() [1/2]

GlobalNumbering* apf::makeGlobal ( Numbering n,
bool  destroy = true 
)

converts a local numbering into a global numbering.

Parameters
destroyShould the input Numbering* be destroyed?

the original local numbering is destroyed. All local numbers are increased by an offset; the offset on part P is the sum of the numbered nodes on parts [0,P-1].

this is done in O(log N) time in parallel, where N is the part count.

this function has no intrinsic knowledge of ownership, it operates simply on nodes which have been explicitly numbered. the input to this function is usually a numbering produced by a numberOwned* function, and the result is a global numbering of all the owned nodes. subsequently calling synchronize on the global numbering completes the typical process which leaves all nodes (owned and not) with a global number attached.

◆ makeGlobal() [2/2]

void apf::makeGlobal ( std::vector< Numbering * > &  owned,
std::vector< GlobalNumbering * > &  global,
pcu::PCU PCUObj 
)

Globalize a mixed numbering scheme.

Parameters
ownedThe local owned on-part numberings.
globalThe output global numberings.

◆ makeMdsBox()

Mesh2* apf::makeMdsBox ( int  nx,
int  ny,
int  nz,
double  wx,
double  wy,
double  wz,
bool  is,
pcu::PCU PCUObj 
)

create a box from an MDS mesh

Parameters
nxnumber of x elements
nynumber of y elements
nznumber of z elements
wxx dimension width
wyy dimension width
wzz dimension width
istrue = simplical mesh, false = quad/hex

set ny,nz=0 for a 1D mesh, set nz=0 for a 2D mesh

◆ makeZoltanBalancer()

Balancer* apf::makeZoltanBalancer ( Mesh mesh,
int  method,
int  approach,
bool  debug = false 
)

Make a Zoltan Balancer object.

this Balancer will apply Zoltan to the global mesh to improve load balance.

Also note that this Balancer will create a Zoltan edge between two elements that share matched faces.

Parameters
methodselect from apf::ZoltanMethod
approachselect from apf::ZoltanApproach
debugprint the full Zoltan configuration

◆ makeZoltanGlobalSplitter()

Splitter* apf::makeZoltanGlobalSplitter ( Mesh mesh,
int  method,
int  approach,
bool  debug = false 
)

Make a Zoltan Splitter object.

the resulting splitter will apply Zoltan to the global mesh part to break it into several new parts.

Parameters
methodselect from apf::ZoltanMethod
approachselect from apf::ZoltanApproach
debugprint the full Zoltan configuration when splitting

◆ makeZoltanSplitter()

Splitter* apf::makeZoltanSplitter ( Mesh mesh,
int  method,
int  approach,
bool  debug = false,
bool  sync = true 
)

Make a Zoltan Splitter object.

the resulting splitter will apply Zoltan to the local mesh part to break it into several new parts.

Parameters
methodselect from apf::ZoltanMethod
approachselect from apf::ZoltanApproach
debugprint the full Zoltan configuration when splitting
syncall parts are splitting by the same factor, multiply the part ids in the resulting apf::Migration accordingly

◆ measure() [1/2]

double apf::measure ( Mesh m,
MeshEntity *  e 
)

Measures the volume, area, or length of a Mesh Entity.

Returns
The measure of the element

◆ measure() [2/2]

double apf::measure ( MeshElement e)

Measures the volume, area, or length of a Mesh Element.

By integrating the differential volume over the element, a general measure is obtained. This correctly measures curved meshes.

Returns
The measure of the element

◆ migrate()

void apf::migrate ( Mesh2 m,
Migration plan 
)

APF's migration function, works on apf::Mesh2.

if your database implements apf::Mesh2 (and residence is separate from remote copies) then you may use this to implement most of apf::Mesh::migrate. Users of APF are encouraged to call apf::Mesh::migrate instead of calling this directly

◆ numberGhost()

int apf::numberGhost ( std::vector< Field * > const &  fields,
std::vector< Numbering * > &  ghost 
)

Number the ghost (overlapped/shared) nodes of multiple fields.

Parameters
fieldsThe input fields to be numbered.
ghostThe output local ghost numberings.
Returns
The number of on-part ghost nodes across all fields.

◆ numberOverlapNodes()

Numbering* apf::numberOverlapNodes ( Mesh mesh,
const char *  name,
FieldShape s = 0 
)

number all local nodes

Parameters
sif non-zero, use nodes from this FieldShape, otherwise use the mesh's coordinate nodes

◆ numberOwned()

int apf::numberOwned ( std::vector< Field * > const &  fields,
std::vector< Numbering * > &  owned 
)

Number the owned nodes of multiple fields.

Parameters
fieldsThe input fields to be numbered.
ownedThe output local owned numberings.
Returns
The number of on-part owned nodes across all fields.

◆ numberOwnedNodes()

Numbering* apf::numberOwnedNodes ( Mesh mesh,
const char *  name,
FieldShape s = 0,
Sharing shr = 0 
)

number the local owned nodes

Parameters
sif non-zero, use nodes from this FieldShape, otherwise use the mesh's coordinate nodes
shrif non-zero, use this Sharing to determine ownership, otherwise call apf::getSharing

◆ recoverGradientByVolume()

Field* apf::recoverGradientByVolume ( Field *  f)

Compute a nodal gradient field from a nodal input field.

given a nodal field, compute approximate nodal gradient values by giving each node a volume-weighted average of the gradients computed at each element around it.

◆ remapPartition()

void apf::remapPartition ( apf::Mesh2 m,
Remap remap 
)

remap all part ids in the mesh structure

when using sub-group partitioning schemes or splitting meshes (see Parma_SplitPartition or apf::Splitter), it is useful to be able to update all partition model structures in a mesh to reflect a transition from one partitioning scheme to the next.

this function applies the given map to all part ids in the remote copies, resident sets, and matching using the apf::Mesh2 interface

◆ reorder()

MeshTag* apf::reorder ( Mesh mesh,
const char *  name 
)

Number by adjacency graph traversal.

a plain single-integer tag is used to number the vertices and elements of a mesh

◆ reorderMdsMesh()

void apf::reorderMdsMesh ( Mesh2 mesh,
MeshTag *  t = 0 
)

apply adjacency-based reordering

Parameters
tOptional user-defined ordering of the vertices. Set this to NULL to use the internal ordering system. Otherwise, attach a unique integer to each vertex in the range [0, #vertices). this will indicate the order in which they appear after reordering.

similar to the algorithm for apf::reorder, this function will traverse adjacencies to reorder each topological type. Then all MDS arrays are re-formed in this new order. An important side effect of this function is that there are no gaps in the MDS arrays after this

◆ rotate()

Matrix3x3 apf::rotate ( Vector3 const &  u,
double  a 
)

get the rotation matrix around an axis

Parameters
uthe axis to rotate around
athe amount of rotation in radians

◆ setCoords()

void apf::setCoords ( Mesh2 m,
const double *  coords,
int  nverts,
GlobalToVert globalToVert 
)

Assign coordinates to the mesh.

Each peer provides a set of the coordinates. The coords most be ordered according to the global ids of the vertices. Peer 0 provides the coords for vertices 0 to m-1, peer to for m to n-1, ... After this call, all vertices in the apf::Mesh2 object have correct coordinates assigned.

◆ setMatches()

void apf::setMatches ( Mesh2 m,
const Gid *  matches,
int  nverts,
GlobalToVert globalToVert 
)

Assign matching to the mesh.

Each peer provides a set of the matched entity global ids. An id set to -1 indicates that the vertex is not matched. The ids most be ordered according to the global ids of the vertices. Peer 0 provides the ids for vertices 0 to m-1, peer to for m to n-1, ... After this call, all vertices in the apf::Mesh2 object have correct coordinates assigned.

◆ setMatrix()

void apf::setMatrix ( Field *  f,
MeshEntity *  e,
int  node,
Matrix3x3 const &  value 
)

Set the nodal value of a matrix field.

Parameters
valueThe vector value.

◆ setMigrationLimit()

void apf::setMigrationLimit ( size_t  maxElements,
pcu::PCU PCUObj 
)

set the maximum elements that apf::migrate moves at once

apf::migrate implements gradual limited migration in an effort to help applications keep memory use to a minimum. This function globally sets the limit on migration, which causes any migration requests greater than the limit to be performed as several consecutive migrations.

◆ setScalar()

void apf::setScalar ( Field *  f,
MeshEntity *  e,
int  node,
double  value 
)

Set a nodal value of a scalar field.

Parameters
nodeThe node number within the entity. So far, it is just 0 for vertices and for edges in 2nd order. higher order will bring more nodes per edge and onto faces and such.

◆ setVector()

void apf::setVector ( Field *  f,
MeshEntity *  e,
int  node,
Vector3 const &  value 
)

Set the nodal value of a vector field.

Parameters
valueThe vector value.

◆ sharedReduction()

void apf::sharedReduction ( Field *  f,
Sharing shr,
bool  delete_shr,
const ReductionOp< double > &  sum = ReductionSum< double >() 
)

Apply a reudction operator along partition boundaries.

Using the copies described by an apf::Sharing object, applied the specified operation pairwise to the values of the field on each partition. No guarantee is made about hte order of the pairwise application

◆ stitchMesh()

void apf::stitchMesh ( Mesh2 m)

infer all remote copies from those of vertices

given that the remote copies of the vertices are set up correctly, this function will synchronize the remote copies and resident part sets for all other entities correctly.

◆ synchronize() [1/2]

void apf::synchronize ( Field *  f,
Sharing shr = 0 
)

Synchronize field values along partition boundary.

Using the ownership and copies described by an apf::Sharing object, copy values from the owned nodes to their copies, possibly assigning them values for the first time.

◆ synchronize() [2/2]

void apf::synchronize ( Numbering n,
Sharing shr = 0,
bool  delete_shr = false 
)

numbers non-owned nodes with the values from their owners

Works even if the non-owned nodes have no number currently assigned, after this they are all numbered

Parameters
shrif non-zero, use this Sharing model to determine ownership and copies, otherwise call apf::getSharing

◆ tagOpposites()

MeshTag* apf::tagOpposites ( GlobalNumbering gn,
const char *  name 
)

Tag boundary faces with global ids of opposite elements.

This function creates a LONG tag of one value and attaches to all mesh faces (edges) classified on a partition model face (edge) (i.e., a face in 3d or edge in 2d that has a remote copy in another part) the global id of the mesh region (face) on the other side.

Parameters
gnglobal element numbering
namethe name of the resulting tag
Returns
a new MeshTag called name with global IDs of opposite elements

◆ unfreezeFields()

void apf::unfreezeFields ( Mesh m)

unfreeze all associated fields

see apf::unfreezeField

◆ unite()

void apf::unite ( Parts into,
Parts const &  from 
)

unite two sets of unique part ids

Parameters
intobecomes the union

◆ verify()

void apf::verify ( Mesh m,
bool  abort_on_error = true 
)

run consistency checks on an apf::Mesh structure

this can be used to implement apf::Mesh::verify. Other implementations may define their own.

◆ writeASCIIVtkFiles() [1/2]

void apf::writeASCIIVtkFiles ( const char *  prefix,
Mesh m 
)

Write a set of parallel VTK Unstructured Mesh files from an apf::Mesh with ASCII encoding.

Nodal fields whose shape differs from the mesh shape will not be output. Fields with incomplete data will not be output.

◆ writeASCIIVtkFiles() [2/2]

void apf::writeASCIIVtkFiles ( const char *  prefix,
Mesh m,
std::vector< std::string >  writeFields 
)

Write a set of parallel VTK Unstructured Mesh files from an apf::Mesh with ASCII encoding.

Only fields whose name appears in the vector writeFields will be output. Nodal fields whose shape differs from the mesh shape will not be output. Fields with incomplete data will not be output.

◆ writeNedelecVtkFiles()

void apf::writeNedelecVtkFiles ( const char *  prefix,
Mesh m 
)

Output .vtk files with ASCII encoding for this part.

this function is useful for debugging meshes with Nedelec fields on them.

◆ writeOneVtkFile()

void apf::writeOneVtkFile ( const char *  prefix,
Mesh m 
)

Output just the .vtu file with ASCII encoding for this part.

this function is useful for debugging large parallel meshes.

◆ writeVtkFiles() [1/2]

void apf::writeVtkFiles ( const char *  prefix,
Mesh m,
int  cellDim = -1 
)

Write a set of parallel VTK Unstructured Mesh files from an apf::Mesh with binary (base64) encoding and zlib compression (if LION_COMPRESS=ON)

Nodal fields whose shape differs from the mesh shape will not be output. Fields with incomplete data will not be output.

◆ writeVtkFiles() [2/2]

void apf::writeVtkFiles ( const char *  prefix,
Mesh m,
std::vector< std::string >  writeFields,
int  cellDim = -1 
)

Write a set of parallel VTK Unstructured Mesh files from an apf::Mesh with binary (base64) encoding and zlib compression (if LION_COMPRESS=ON)

Only fields whose name appears in the vector writeFields will be output. Nodal fields whose shape differs from the mesh shape will not be output. Fields with incomplete data will not be output.

Variable Documentation

◆ pi

double const apf::pi
extern

The mathematical constant pi.

although it doesn't fit perfectly in this header, there is no more appropriate place for this in APF.