Compadre  1.2.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Classes | Namespaces | Macros | Enumerations | Variables
Compadre_Operators.hpp File Reference
#include "Compadre_Typedefs.hpp"

Go to the source code of this file.

Classes

struct  Compadre::SamplingFunctional
 

Namespaces

 Compadre
 

Macros

#define make_sampling_functional(input, output, targets, nontrivial, transform)   SamplingFunctional(input, output, targets, nontrivial, transform, __COUNTER__)
 

Enumerations

enum  Compadre::TargetOperation {
  Compadre::ScalarPointEvaluation, Compadre::VectorPointEvaluation, Compadre::LaplacianOfScalarPointEvaluation, Compadre::VectorLaplacianPointEvaluation,
  Compadre::GradientOfScalarPointEvaluation, Compadre::GradientOfVectorPointEvaluation, Compadre::DivergenceOfVectorPointEvaluation, Compadre::CurlOfVectorPointEvaluation,
  Compadre::CurlCurlOfVectorPointEvaluation, Compadre::PartialXOfScalarPointEvaluation, Compadre::PartialYOfScalarPointEvaluation, Compadre::PartialZOfScalarPointEvaluation,
  Compadre::ChainedStaggeredLaplacianOfScalarPointEvaluation, Compadre::GaussianCurvaturePointEvaluation, Compadre::ScalarFaceAverageEvaluation, Compadre::COUNT =15
}
 Available target functionals. More...
 
enum  Compadre::ReconstructionSpace { Compadre::ScalarTaylorPolynomial, Compadre::VectorTaylorPolynomial, Compadre::VectorOfScalarClonesTaylorPolynomial, Compadre::DivergenceFreeVectorTaylorPolynomial }
 Space in which to reconstruct polynomial. More...
 
enum  Compadre::SamplingTransformType { Compadre::Identity, Compadre::SameForAll, Compadre::DifferentEachTarget, Compadre::DifferentEachNeighbor }
 Describes the SamplingFunction relationship to targets, neighbors. More...
 
enum  Compadre::DenseSolverType { Compadre::QR, Compadre::SVD, Compadre::LU }
 Dense solver type. More...
 
enum  Compadre::ProblemType { Compadre::STANDARD, Compadre::MANIFOLD }
 Problem type, that optionally can handle manifolds. More...
 
enum  Compadre::ConstraintType { Compadre::NO_CONSTRAINT, Compadre::NEUMANN_GRAD_SCALAR }
 Constraint type. More...
 
enum  Compadre::WeightingFunctionType { Compadre::Power, Compadre::Gaussian, Compadre::CubicSpline }
 Available weighting kernel function types. More...
 
enum  Compadre::CoordinatesType { Compadre::Ambient, Compadre::Local }
 Coordinate type for input and output format of vector data on manifold problems. More...
 

Variables

constexpr int Compadre::TargetOutputTensorRank []
 Rank of target functional output for each TargetOperation Rank of target functional input for each TargetOperation is based on the output rank of the SamplingFunctional used on the polynomial basis. More...
 
constexpr int Compadre::ActualReconstructionSpaceRank []
 Number of actual components in the ReconstructionSpace. More...
 
constexpr SamplingFunctional Compadre::PointSample = make_sampling_functional(0,0,false,false,(int)Identity)
 Available sampling functionals. More...
 
constexpr SamplingFunctional Compadre::VectorPointSample = make_sampling_functional(1,1,false,false,(int)Identity)
 Point evaluations of the entire vector source function. More...
 
constexpr SamplingFunctional Compadre::ManifoldVectorPointSample = make_sampling_functional(1,1,false,false,(int)DifferentEachTarget)
 Point evaluations of the entire vector source function (but on a manifold, so it includes a transform into local coordinates) More...
 
constexpr SamplingFunctional Compadre::StaggeredEdgeAnalyticGradientIntegralSample = make_sampling_functional(0,0,true,true,(int)SameForAll)
 Analytical integral of a gradient source vector is just a difference of the scalar source at neighbor and target. More...
 
constexpr SamplingFunctional Compadre::StaggeredEdgeIntegralSample = make_sampling_functional(1,0,true,true,(int)DifferentEachNeighbor)
 Samples consist of the result of integrals of a vector dotted with the tangent along edges between neighbor and target. More...
 
constexpr SamplingFunctional Compadre::VaryingManifoldVectorPointSample = make_sampling_functional(1,1,false,false,(int)DifferentEachNeighbor)
 For integrating polynomial dotted with normal over an edge. More...
 
constexpr SamplingFunctional Compadre::FaceNormalIntegralSample = make_sampling_functional(1,0,false,false,(int)Identity)
 For integrating polynomial dotted with normal over an edge. More...
 
constexpr SamplingFunctional Compadre::FaceNormalPointSample = make_sampling_functional(1,0,false,false,(int)Identity)
 For polynomial dotted with normal on edge. More...
 
constexpr SamplingFunctional Compadre::FaceTangentIntegralSample = make_sampling_functional(1,0,false,false,(int)Identity)
 For integrating polynomial dotted with tangent over an edge. More...
 
constexpr SamplingFunctional Compadre::FaceTangentPointSample = make_sampling_functional(1,0,false,false,(int)Identity)
 For polynomial dotted with tangent. More...
 
constexpr SamplingFunctional Compadre::ScalarFaceAverageSample = make_sampling_functional(0,0,false,false,(int)DifferentEachNeighbor)
 For polynomial integrated on faces. More...
 

Macro Definition Documentation

#define make_sampling_functional (   input,
  output,
  targets,
  nontrivial,
  transform 
)    SamplingFunctional(input, output, targets, nontrivial, transform, __COUNTER__)

Definition at line 6 of file Compadre_Operators.hpp.