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.
orthoPoly3D Class Reference

Detailed Description

Definition at line 77 of file orthoPoly3D.H.

Public Member Functions

 orthoPoly3D ()
 
 orthoPoly3D (int _order, const Eigen::VectorXd &_sigma, const std::vector< double > &_x, const std::vector< double > &_y, const std::vector< double > &_z)
 
 orthoPoly3D (int _order, const std::vector< std::vector< double >> &coords)
 
 orthoPoly3D (orthoPoly3D &&op)
 
orthoPoly3Doperator= (orthoPoly3D &&op)
 
 orthoPoly3D (const orthoPoly3D &)=delete
 
orthoPoly3Doperator= (const orthoPoly3D &)=delete
 
 ~orthoPoly3D ()=default
 
void computeA (const Eigen::VectorXd &sigma)
 
double operator() (const std::vector< double > &coord)
 
double eval (const std::vector< double > &coord)
 
bool status () const
 
void resetA ()
 
void Reset ()
 
Eigen::VectorXd getCoeffs () const
 

Static Public Member Functions

static orthoPoly3DCreate (int _order, const std::vector< std::vector< double >> &coords)
 
static std::unique_ptr< orthoPoly3DCreateUnique (int _order, const std::vector< std::vector< double >> &coords)
 

Private Member Functions

void initCheck ()
 

Private Attributes

int order
 
std::unique_ptr< orthoPoly1Dopx
 
std::unique_ptr< orthoPoly1Dopy
 
std::unique_ptr< orthoPoly1Dopz
 
Eigen::VectorXd a
 
std::vector< unsigned int > toRemove
 
bool finished
 

Constructor & Destructor Documentation

◆ orthoPoly3D() [1/5]

orthoPoly3D::orthoPoly3D ( )

Definition at line 110 of file orthoPoly3D.C.

References a.

Referenced by Create(), and initCheck().

110  : opx(new orthoPoly1D()),
111  opy(new orthoPoly1D()),
112  opz(new orthoPoly1D()),
113  finished(false)
114 {
115  a.resize(0);
116 }
Eigen::VectorXd a
Definition: orthoPoly3D.H:146
std::unique_ptr< orthoPoly1D > opy
Definition: orthoPoly3D.H:142
std::unique_ptr< orthoPoly1D > opx
Definition: orthoPoly3D.H:140
std::unique_ptr< orthoPoly1D > opz
Definition: orthoPoly3D.H:144

◆ orthoPoly3D() [2/5]

orthoPoly3D::orthoPoly3D ( int  _order,
const Eigen::VectorXd &  _sigma,
const std::vector< double > &  _x,
const std::vector< double > &  _y,
const std::vector< double > &  _z 
)

◆ orthoPoly3D() [3/5]

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

Definition at line 153 of file orthoPoly3D.C.

References initCheck(), opx, opy, opz, order, and nemAux::printVec().

155  : order(_order), finished(false)
156 {
157  //int root(round(cbrt(coords.size())));
158  //if (coords.size() != root * root * root &&
159  // coords.size() != (root + 1) * (root + 1) * (root + 1))
160  //{
161  // std::cerr << "coords are definitely not in a cartesian product set"
162  // << std::endl;
163  // exit(1);
164  //}
165 
166  initCheck();
167 
168  std::vector<double> x(coords.size());
169  std::vector<double> y(coords.size());
170  std::vector<double> z(coords.size());
171 
172  for (int i = 0; i < coords.size(); ++i)
173  {
174  x[i] = coords[i][0];
175  y[i] = coords[i][1];
176  z[i] = coords[i][2];
177  }
178  std::cout << "non-unique coords" << std::endl;
179 #ifdef NEMOSYS_DEBUG
180  nemAux::printVec(x); //nemAux::printVec(y); nemAux::printVec(z);
181 #endif
182  //std::vector<double> x_unique;
183  //for (int i = 0; i < x.size(); ++i)
184  //{
185  // x_unique.push_back(x[i]);
186  // for (int j = i + 1; j < x.size(); ++j)
187  // {
188  // if (x_unique[i] != x[j])
189  // }
190  //}
191 
192  //std::sort(x.begin(), x.end());
193  //x.erase(std::unique(x.begin(), x.end(), Pred), x.end());
194  //std::cout << "unique coords" << std::endl;
195  //nemAux::printVec(x); //nemAux::printVec(y); nemAux::printVec(z);
196 
197  opx = std::unique_ptr<orthoPoly1D>(new orthoPoly1D(order, x));
198  opy = std::unique_ptr<orthoPoly1D>(new orthoPoly1D(order, y));
199  opz = std::unique_ptr<orthoPoly1D>(new orthoPoly1D(order, z));
200 }
std::unique_ptr< orthoPoly1D > opy
Definition: orthoPoly3D.H:142
std::unique_ptr< orthoPoly1D > opx
Definition: orthoPoly3D.H:140
void initCheck()
Definition: orthoPoly3D.C:118
void printVec(const std::vector< T > &v)
std::unique_ptr< orthoPoly1D > opz
Definition: orthoPoly3D.H:144

◆ orthoPoly3D() [4/5]

orthoPoly3D::orthoPoly3D ( orthoPoly3D &&  op)

Definition at line 203 of file orthoPoly3D.C.

204  : order(op.order),
205  opx(std::move(op.opx)),
206  opy(std::move(op.opy)),
207  opz(std::move(op.opz)),
208  a(op.a),
209  toRemove(op.toRemove),
210  finished(op.finished)
211 {
212  op.Reset();
213 }
Eigen::VectorXd a
Definition: orthoPoly3D.H:146
std::unique_ptr< orthoPoly1D > opy
Definition: orthoPoly3D.H:142
std::unique_ptr< orthoPoly1D > opx
Definition: orthoPoly3D.H:140
void Reset()
Definition: orthoPoly3D.C:376
std::vector< unsigned int > toRemove
Definition: orthoPoly3D.H:148
std::unique_ptr< orthoPoly1D > opz
Definition: orthoPoly3D.H:144

◆ orthoPoly3D() [5/5]

orthoPoly3D::orthoPoly3D ( const orthoPoly3D )
delete

◆ ~orthoPoly3D()

orthoPoly3D::~orthoPoly3D ( )
default

Member Function Documentation

◆ computeA()

void orthoPoly3D::computeA ( const Eigen::VectorXd &  sigma)

Definition at line 254 of file orthoPoly3D.C.

References a, nemAux::Timer::elapsed(), finished, opx, opy, opz, nemAux::printVec(), removeRowT(), nemAux::Timer::start(), nemAux::Timer::stop(), and toRemove.

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

255 {
256  if (!finished)
257  {
258 #ifdef NEMOSYS_DEBUG
259  nemAux::Timer T;
260 #endif
261 
262  MatrixXd tmp1 = opx->phi * opx->phiTphiInv;
263  MatrixXd tmp2 = opy->phi * opy->phiTphiInv;
264  MatrixXd tmp3 = opz->phi * opz->phiTphiInv;
265 #ifdef NEMOSYS_DEBUG
266  nemAux::Timer T1;
267  T1.start();
268 #endif
269  MatrixXd kronProd_full = kroneckerProduct(tmp1, tmp2);
270  kronProd_full = kroneckerProduct(kronProd_full, tmp3).eval();
271 #ifdef NEMOSYS_DEBUG
272  T1.stop();
273  std::cout << "Time for building kronprod in constructor: " << T1.elapsed()
274  << "\n";
275  std::cout << "kronProd_full: " << kronProd_full.rows() << " "
276  << kronProd_full.cols() << "\n";
277  std::cout << "Removing rows: ";
279 #endif
280  MatrixXdRM kronProd_red(kronProd_full.cols() - toRemove.size(),
281  kronProd_full.rows());
282 #ifdef NEMOSYS_DEBUG
283  T.start();
284 #endif
285  removeRowT(kronProd_full, kronProd_red, toRemove);
286  std::cout << "kronProd_red: " << kronProd_red.rows() << " "
287  << kronProd_red.cols() << "\n";
288  std::cout << "sigma size: " << sigma.size() << "\n";
289  a = kronProd_red * sigma;
290 #ifdef NEMOSYS_DEBUG
291  T.stop();
292  std::cout << "Time: " << T.elapsed() << std::endl;
293  std::cout << "shape of a: " << a.rows() << " " << a.cols() << std::endl;
294 #endif
295  finished = true;
296  }
297 }
Eigen::VectorXd a
Definition: orthoPoly3D.H:146
std::unique_ptr< orthoPoly1D > opy
Definition: orthoPoly3D.H:142
std::unique_ptr< orthoPoly1D > opx
Definition: orthoPoly3D.H:140
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor > MatrixXdRM
Definition: orthoPoly3D.H:61
void removeRowT(const MatrixXd &matrix, MatrixXdRM &matrix_red, const std::vector< unsigned int > &toRemove)
Definition: orthoPoly3D.C:85
std::vector< unsigned int > toRemove
Definition: orthoPoly3D.H:148
void printVec(const std::vector< T > &v)
std::unique_ptr< orthoPoly1D > opz
Definition: orthoPoly3D.H:144

◆ Create()

orthoPoly3D * orthoPoly3D::Create ( int  _order,
const std::vector< std::vector< double >> &  coords 
)
static

Definition at line 239 of file orthoPoly3D.C.

References orthoPoly3D().

Referenced by CreateUnique().

240 {
241  return new orthoPoly3D(_order, coords);
242 }

◆ CreateUnique()

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

Definition at line 245 of file orthoPoly3D.C.

References Create().

Referenced by PatchRecovery::recoverNodalSolution().

247 {
248  return std::unique_ptr<orthoPoly3D>(
249  orthoPoly3D::Create(_order, coords));
250 }
static orthoPoly3D * Create(int _order, const std::vector< std::vector< double >> &coords)
Definition: orthoPoly3D.C:239

◆ eval()

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

Definition at line 309 of file orthoPoly3D.C.

References a, finished, opx, opy, opz, order, removeColumn(), and toRemove.

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

310 {
311  if (!finished)
312  {
313  std::cerr << "Coefficients to basis have not been generated" << std::endl;
314  exit(1);
315  }
316 
317  MatrixXd phix(1, order + 1);
318  MatrixXd phiy(1, order + 1);
319  MatrixXd phiz(1, order + 1);
320  for (int i = 0; i < order + 1; ++i)
321  {
322  phix(0, i) = opx->EvaluateOrthogonal(i, coord[0]);
323  phiy(0, i) = opy->EvaluateOrthogonal(i, coord[1]);
324  phiz(0, i) = opz->EvaluateOrthogonal(i, coord[2]);
325  }
326 
327  MatrixXd Phi = kroneckerProduct(phix, phiy);
328  Phi = kroneckerProduct(Phi, phiz).eval();
329 #ifdef NEMOSYS_DEBUG
330  std::cout << Phi.rows() << " " << Phi.cols() << std::endl;
331 #endif
332  MatrixXd Phi_red(Phi.rows(), Phi.cols() - toRemove.size());
333  removeColumn(Phi, Phi_red, toRemove);
334 
335 #ifdef NEMOSYS_DEBUG
336  std::cout << Phi_red.rows() << " " << Phi_red.cols() << std::endl;
337  std::cout << a.rows() << " " << a.cols() << std::endl;
338 #endif
339 
340  return (Phi_red * a)(0, 0);
341 }
Eigen::VectorXd a
Definition: orthoPoly3D.H:146
std::unique_ptr< orthoPoly1D > opy
Definition: orthoPoly3D.H:142
std::unique_ptr< orthoPoly1D > opx
Definition: orthoPoly3D.H:140
std::vector< unsigned int > toRemove
Definition: orthoPoly3D.H:148
void removeColumn(MatrixXd &matrix, unsigned int colToRemove)
Definition: orthoPoly3D.C:42
std::unique_ptr< orthoPoly1D > opz
Definition: orthoPoly3D.H:144

◆ getCoeffs()

Eigen::VectorXd orthoPoly3D::getCoeffs ( ) const
inline

Definition at line 134 of file orthoPoly3D.H.

134 { return a; }
Eigen::VectorXd a
Definition: orthoPoly3D.H:146

◆ initCheck()

void orthoPoly3D::initCheck ( )
private

Definition at line 118 of file orthoPoly3D.C.

References computeA(), finished, opx, opy, opz, order, orthoPoly3D(), and toRemove.

Referenced by orthoPoly3D().

119 {
120  if (order == 1)
121  toRemove = {3, 5, 6, 7};
122  else if (order == 2)
123  toRemove = {5, 7, 8, 11, 13, 14, 15, 16, 17, 19, 20, 21, 22, 23, 24, 25,
124  26};
125  else if (order == 3)
126  toRemove = {7, 10, 11, 13, 14, 15, 19, 22, 23, 25, 26,
127  27, 28, 29, 30, 31, 34, 35, 37, 38, 39, 40,
128  41, 42, 43, 44, 45, 46, 47, 49, 50, 51, 52,
129  53, 54, 55, 56, 57, 58, 59, 60, 61, 62, 63};
130  else
131  {
132  std::cerr << "Polynomial order greater than 3 is not supported"
133  << std::endl;
134  exit(1);
135  }
136 }
std::vector< unsigned int > toRemove
Definition: orthoPoly3D.H:148

◆ operator()()

double orthoPoly3D::operator() ( const std::vector< double > &  coord)

Definition at line 344 of file orthoPoly3D.C.

References a, nemAux::Timer::elapsed(), eval(), opx, opy, opz, removeColumn(), nemAux::Timer::start(), nemAux::Timer::stop(), and toRemove.

345 {
346  return eval(coord);
347 }
double eval(const std::vector< double > &coord)
Definition: orthoPoly3D.C:309

◆ operator=() [1/2]

orthoPoly3D & orthoPoly3D::operator= ( orthoPoly3D &&  op)

Definition at line 217 of file orthoPoly3D.C.

References a, finished, opx, opy, opz, order, Reset(), and toRemove.

218 {
219  // check against self assignment
220  if (this != &op)
221  {
222  // release held resources
223  this->Reset();
224  // move resources from op to this (i.e. pilfer)
225  opx = std::move(op.opx);
226  opy = std::move(op.opy);
227  opz = std::move(op.opz);
228  a = op.a;
229  toRemove = op.toRemove;
230  finished = op.finished;
231  order = op.order;
232  // set op to default state
233  op.Reset();
234  }
235  return *this;
236 }
Eigen::VectorXd a
Definition: orthoPoly3D.H:146
std::unique_ptr< orthoPoly1D > opy
Definition: orthoPoly3D.H:142
std::unique_ptr< orthoPoly1D > opx
Definition: orthoPoly3D.H:140
void Reset()
Definition: orthoPoly3D.C:376
std::vector< unsigned int > toRemove
Definition: orthoPoly3D.H:148
std::unique_ptr< orthoPoly1D > opz
Definition: orthoPoly3D.H:144

◆ operator=() [2/2]

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

◆ Reset()

void orthoPoly3D::Reset ( )

Definition at line 376 of file orthoPoly3D.C.

References a, finished, opx, opy, opz, and toRemove.

Referenced by operator=().

377 {
378  if (opx) opx.reset();
379  if (opy) opy.reset();
380  if (opz) opz.reset();
381  a.resize(0);
382  finished = false;
383  toRemove.resize(0);
384 }
Eigen::VectorXd a
Definition: orthoPoly3D.H:146
std::unique_ptr< orthoPoly1D > opy
Definition: orthoPoly3D.H:142
std::unique_ptr< orthoPoly1D > opx
Definition: orthoPoly3D.H:140
std::vector< unsigned int > toRemove
Definition: orthoPoly3D.H:148
std::unique_ptr< orthoPoly1D > opz
Definition: orthoPoly3D.H:144

◆ resetA()

void orthoPoly3D::resetA ( )

Definition at line 299 of file orthoPoly3D.C.

References a, and finished.

Referenced by PatchRecovery::recoverNodalSolution().

300 {
301  if (finished)
302  {
303  a.resize(0);
304  finished = false;
305  }
306 }
Eigen::VectorXd a
Definition: orthoPoly3D.H:146

◆ status()

bool orthoPoly3D::status ( ) const
inline

Definition at line 127 of file orthoPoly3D.H.

127 { return finished; }

Member Data Documentation

◆ a

Eigen::VectorXd orthoPoly3D::a
private

Definition at line 146 of file orthoPoly3D.H.

Referenced by computeA(), eval(), operator()(), operator=(), orthoPoly3D(), Reset(), and resetA().

◆ finished

bool orthoPoly3D::finished
private

Definition at line 150 of file orthoPoly3D.H.

Referenced by computeA(), eval(), initCheck(), operator=(), Reset(), and resetA().

◆ opx

std::unique_ptr<orthoPoly1D> orthoPoly3D::opx
private

Definition at line 140 of file orthoPoly3D.H.

Referenced by computeA(), eval(), initCheck(), operator()(), operator=(), orthoPoly3D(), and Reset().

◆ opy

std::unique_ptr<orthoPoly1D> orthoPoly3D::opy
private

Definition at line 142 of file orthoPoly3D.H.

Referenced by computeA(), eval(), initCheck(), operator()(), operator=(), orthoPoly3D(), and Reset().

◆ opz

std::unique_ptr<orthoPoly1D> orthoPoly3D::opz
private

Definition at line 144 of file orthoPoly3D.H.

Referenced by computeA(), eval(), initCheck(), operator()(), operator=(), orthoPoly3D(), and Reset().

◆ order

int orthoPoly3D::order
private

Definition at line 138 of file orthoPoly3D.H.

Referenced by eval(), initCheck(), operator=(), and orthoPoly3D().

◆ toRemove

std::vector<unsigned int> orthoPoly3D::toRemove
private

Definition at line 148 of file orthoPoly3D.H.

Referenced by computeA(), eval(), initCheck(), operator()(), operator=(), and Reset().


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