NEMoSys  0.63.0
A modular, extensible resource with robust automated mesh generation, mesh quality analysis, adaptive mesh refinement, and data transfer between arbitrary meshes.
polyApprox Class Reference

Detailed Description

Definition at line 39 of file polyApprox.H.

Public Member Functions

 polyApprox (int _order, const std::vector< std::vector< double >> &coords)
 
 ~polyApprox ()=default
 
 polyApprox (const polyApprox &)=delete
 
polyApproxoperator= (const polyApprox &)=delete
 
void computeCoeff (const Eigen::VectorXd &data)
 
void resetCoeff ()
 
double eval (const std::vector< double > &coord) const
 

Static Public Member Functions

static std::unique_ptr< polyApproxCreateUnique (int order, const std::vector< std::vector< double >> &coords)
 

Private Member Functions

Eigen::VectorXd computeBasis (const std::vector< double > &coord) const
 

Private Attributes

int order
 
std::vector< Eigen::VectorXd > basis
 
Eigen::MatrixXd A
 
Eigen::VectorXd b
 
Eigen::VectorXd a
 

Constructor & Destructor Documentation

◆ polyApprox() [1/2]

polyApprox::polyApprox ( int  _order,
const std::vector< std::vector< double >> &  coords 
)

Definition at line 41 of file polyApprox.C.

References A, a, b, basis, computeBasis(), and order.

Referenced by CreateUnique().

43  : order(_order)//, coords(_coords)
44 {
45  switch (order)
46  {
47  case 1:
48  {
49  A = MatrixXd::Zero(4, 4);
50  b = VectorXd::Zero(4);
51  a.resize(4);
52  break;
53  }
54  case 2:
55  {
56  A = MatrixXd::Zero(10, 10);//.resize(10, 10);
57  b = VectorXd::Zero(10);//resize(10);
58  a.resize(10);
59  break;
60  }
61  default:
62  {
63  std::cerr << "Error: order: " << order << " is not supported"
64  << std::endl;
65  exit(1);
66  }
67  }
68  basis.resize(coords.size());
69  for (int i = 0; i < coords.size(); ++i)
70  {
71  basis[i] = computeBasis(coords[i]);
72  A = A + basis[i] * (basis[i].transpose()); // * basis[i];
73  }
74 }
Eigen::VectorXd a
Definition: polyApprox.H:70
std::vector< Eigen::VectorXd > basis
Definition: polyApprox.H:64
Eigen::MatrixXd A
Definition: polyApprox.H:66
Eigen::VectorXd b
Definition: polyApprox.H:68
Eigen::VectorXd computeBasis(const std::vector< double > &coord) const
Definition: polyApprox.C:109

◆ ~polyApprox()

polyApprox::~polyApprox ( )
default

◆ polyApprox() [2/2]

polyApprox::polyApprox ( const polyApprox )
delete

Member Function Documentation

◆ computeBasis()

VectorXd polyApprox::computeBasis ( const std::vector< double > &  coord) const
private

Definition at line 109 of file polyApprox.C.

References order.

Referenced by eval(), and polyApprox().

110 {
111  switch (order)
112  {
113  case 1:
114  {
115  Vector4d basisVec;
116  basisVec(0) = 1;
117  basisVec(1) = coord[0];
118  basisVec(2) = coord[1];
119  basisVec(3) = coord[2];
120  return basisVec;
121  }
122  case 2:
123  {
124  Vector10d basisVec;
125  basisVec(0) = 1;
126  basisVec(1) = coord[0];
127  basisVec(2) = coord[1];
128  basisVec(3) = coord[2];
129  basisVec(4) = coord[0] * coord[0];
130  basisVec(5) = coord[0] * coord[1];
131  basisVec(6) = coord[0] * coord[2];
132  basisVec(7) = coord[1] * coord[1];
133  basisVec(8) = coord[1] * coord[2];
134  basisVec(9) = coord[2] * coord[2];
135  return basisVec;
136  }
137  default:
138  {
139  std::cerr << "Error: order " << order << " is not supported."
140  << std::endl;
141  exit(1);
142  }
143  }
144 }
Eigen::Matrix< double, 10, 1 > Vector10d
Definition: polyApprox.C:38

◆ computeCoeff()

void polyApprox::computeCoeff ( const Eigen::VectorXd &  data)

Definition at line 85 of file polyApprox.C.

References A, a, b, basis, and data.

Referenced by PatchRecovery::computeNodalError(), and PatchRecovery::recoverNodalSolution().

86 {
87  for (int i = 0; i < basis.size(); ++i)
88  {
89  b = b + basis[i] * data(i); // basis[i].transpose() * data(i);
90  }
91  a = A.partialPivLu().solve(b);
92 }
data_type data
Edge/face with sorted point ids (a, b, c, ...) is located at some index i in data[b], with data[b][i].first == [a, c] (for edges, third point id treated as -1).
Eigen::VectorXd a
Definition: polyApprox.H:70
std::vector< Eigen::VectorXd > basis
Definition: polyApprox.H:64
Eigen::MatrixXd A
Definition: polyApprox.H:66
Eigen::VectorXd b
Definition: polyApprox.H:68

◆ CreateUnique()

std::unique_ptr< polyApprox > polyApprox::CreateUnique ( int  order,
const std::vector< std::vector< double >> &  coords 
)
static

Definition at line 78 of file polyApprox.C.

References polyApprox().

Referenced by PatchRecovery::computeNodalError(), and PatchRecovery::recoverNodalSolution().

80 {
81  return std::unique_ptr<polyApprox>(new polyApprox(order, coords));
82 }
polyApprox(int _order, const std::vector< std::vector< double >> &coords)
Definition: polyApprox.C:41

◆ eval()

double polyApprox::eval ( const std::vector< double > &  coord) const

Definition at line 102 of file polyApprox.C.

References a, and computeBasis().

Referenced by PatchRecovery::computeNodalError(), and PatchRecovery::recoverNodalSolution().

103 {
104  VectorXd P = computeBasis(coord);
105  return (P.transpose() * a)(0, 0);
106 }
Eigen::VectorXd a
Definition: polyApprox.H:70
Eigen::VectorXd computeBasis(const std::vector< double > &coord) const
Definition: polyApprox.C:109

◆ operator=()

polyApprox& polyApprox::operator= ( const polyApprox )
delete

◆ resetCoeff()

void polyApprox::resetCoeff ( )

Definition at line 95 of file polyApprox.C.

References a, and b.

Referenced by PatchRecovery::computeNodalError(), and PatchRecovery::recoverNodalSolution().

96 {
97  a.setZero();
98  b.setZero();
99 }
Eigen::VectorXd a
Definition: polyApprox.H:70
Eigen::VectorXd b
Definition: polyApprox.H:68

Member Data Documentation

◆ A

Eigen::MatrixXd polyApprox::A
private

Definition at line 66 of file polyApprox.H.

Referenced by computeCoeff(), and polyApprox().

◆ a

Eigen::VectorXd polyApprox::a
private

Definition at line 70 of file polyApprox.H.

Referenced by computeCoeff(), eval(), polyApprox(), and resetCoeff().

◆ b

Eigen::VectorXd polyApprox::b
private

Definition at line 68 of file polyApprox.H.

Referenced by computeCoeff(), polyApprox(), and resetCoeff().

◆ basis

std::vector<Eigen::VectorXd> polyApprox::basis
private

Definition at line 64 of file polyApprox.H.

Referenced by computeCoeff(), and polyApprox().

◆ order

int polyApprox::order
private

Definition at line 62 of file polyApprox.H.

Referenced by computeBasis(), and polyApprox().


The documentation for this class was generated from the following files: