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

Detailed Description

Definition at line 39 of file orthoPoly1D.H.

Public Member Functions

 orthoPoly1D (int _order, const std::vector< double > &x)
 
 orthoPoly1D ()=default
 
 orthoPoly1D (const orthoPoly1D &op)
 
orthoPoly1Doperator= (const orthoPoly1D &op)
 
 ~orthoPoly1D ()=default
 
double EvaluateOrthogonal (int power, double xk) const
 
void ComputeAB (const std::vector< double > &x)
 
void EvaluateOrthogonals (const std::vector< double > &x)
 
void ComputePhiTPhiInv (const std::vector< double > &x)
 

Public Attributes

int order
 
std::vector< double > a
 
std::vector< double > b
 
Eigen::MatrixXd phi
 
Eigen::MatrixXd phiTphiInv
 

Constructor & Destructor Documentation

◆ orthoPoly1D() [1/3]

orthoPoly1D::orthoPoly1D ( int  _order,
const std::vector< double > &  x 
)

Definition at line 37 of file orthoPoly1D.C.

References a, b, ComputePhiTPhiInv(), phi, and phiTphiInv.

38  : order(_order)
39 {
40  a.resize(_order + 1);
41  b.resize(_order);
42  phi.resize(x.size(), _order + 1);
43  phiTphiInv.resize(_order + 1, _order + 1);
45 }
std::vector< double > a
Definition: orthoPoly1D.H:77
void ComputePhiTPhiInv(const std::vector< double > &x)
Definition: orthoPoly1D.C:131
std::vector< double > b
Definition: orthoPoly1D.H:79
Eigen::MatrixXd phiTphiInv
Definition: orthoPoly1D.H:83
Eigen::MatrixXd phi
Definition: orthoPoly1D.H:81

◆ orthoPoly1D() [2/3]

orthoPoly1D::orthoPoly1D ( )
default

◆ orthoPoly1D() [3/3]

orthoPoly1D::orthoPoly1D ( const orthoPoly1D op)

Definition at line 47 of file orthoPoly1D.C.

References a, b, order, phi, and phiTphiInv.

48 {
49  order = op.order;
50  a = op.a;
51  b = op.b;
52  phi = op.phi;
54 }
std::vector< double > a
Definition: orthoPoly1D.H:77
std::vector< double > b
Definition: orthoPoly1D.H:79
Eigen::MatrixXd phiTphiInv
Definition: orthoPoly1D.H:83
Eigen::MatrixXd phi
Definition: orthoPoly1D.H:81

◆ ~orthoPoly1D()

orthoPoly1D::~orthoPoly1D ( )
default

Member Function Documentation

◆ ComputeAB()

void orthoPoly1D::ComputeAB ( const std::vector< double > &  x)

Definition at line 91 of file orthoPoly1D.C.

References a, b, EvaluateOrthogonal(), and order.

Referenced by EvaluateOrthogonals().

92 {
93  a[0] = 0.;
94  int n = x.size();
95  for (int i = 0; i < n; ++i)
96  {
97  a[0] += x[i];
98  }
99  a[0] /= n;
100 
101  for (int i = 1; i < order; ++i)
102  {
103  double sum0, sum1, sum2;
104  sum0 = sum1 = sum2 = 0.0;
105  for (int j = 0; j < n; ++j)
106  {
107  double tmp0 = EvaluateOrthogonal(i - 1, x[j]);
108  double tmp1 = EvaluateOrthogonal(i, x[j]);
109  sum0 += tmp0 * tmp0;
110  sum1 += tmp1 * tmp1;
111  sum2 += x[j] * tmp1 * tmp1;
112  }
113  a[i] = sum2 / sum1;
114  b[i - 1] = sum1 / sum0;
115  }
116 }
std::vector< double > a
Definition: orthoPoly1D.H:77
std::vector< double > b
Definition: orthoPoly1D.H:79
double EvaluateOrthogonal(int power, double xk) const
Definition: orthoPoly1D.C:69

◆ ComputePhiTPhiInv()

void orthoPoly1D::ComputePhiTPhiInv ( const std::vector< double > &  x)

Definition at line 131 of file orthoPoly1D.C.

References EvaluateOrthogonals(), phi, and phiTphiInv.

Referenced by orthoPoly1D().

132 {
134  VectorXd tmp = (phi.transpose() * phi).diagonal();
135  for (int i = 0; i < tmp.rows(); ++i)
136  tmp(i) = 1. / tmp(i);
137  phiTphiInv = tmp.asDiagonal();
138 }
void EvaluateOrthogonals(const std::vector< double > &x)
Definition: orthoPoly1D.C:120
Eigen::MatrixXd phiTphiInv
Definition: orthoPoly1D.H:83
Eigen::MatrixXd phi
Definition: orthoPoly1D.H:81

◆ EvaluateOrthogonal()

double orthoPoly1D::EvaluateOrthogonal ( int  power,
double  xk 
) const

Definition at line 69 of file orthoPoly1D.C.

References a, and b.

Referenced by ComputeAB(), and EvaluateOrthogonals().

70 {
71  double p0 = 1.;
72  if (power == 0)
73  return p0;
74  double p1 = xk - a[0];
75  if (power == 1)
76  return p1;
77  double p2 = std::numeric_limits<double>::quiet_NaN();
78  for (int i = 1; i < power; ++i)
79  {
80  p2 = (xk - a[i]) * p1 - b[i - 1] * p0;
81  p0 = p1;
82  p1 = p2;
83  }
84  return p2;
85 }
std::vector< double > a
Definition: orthoPoly1D.H:77
std::vector< double > b
Definition: orthoPoly1D.H:79

◆ EvaluateOrthogonals()

void orthoPoly1D::EvaluateOrthogonals ( const std::vector< double > &  x)

Definition at line 120 of file orthoPoly1D.C.

References ComputeAB(), EvaluateOrthogonal(), order, and phi.

Referenced by ComputePhiTPhiInv().

121 {
122  ComputeAB(x);
123  for (int i = 0; i < order + 1; ++i)
124  {
125  for (int j = 0; j < x.size(); ++j)
126  phi(j, i) = EvaluateOrthogonal(i, x[j]);
127  }
128 }
void ComputeAB(const std::vector< double > &x)
Definition: orthoPoly1D.C:91
double EvaluateOrthogonal(int power, double xk) const
Definition: orthoPoly1D.C:69
Eigen::MatrixXd phi
Definition: orthoPoly1D.H:81

◆ operator=()

orthoPoly1D & orthoPoly1D::operator= ( const orthoPoly1D op)

Definition at line 56 of file orthoPoly1D.C.

References a, b, order, phi, and phiTphiInv.

57 {
58  if (this != &op)
59  {
60  order = op.order;
61  a = op.a;
62  b = op.b;
63  phi = op.phi;
65  }
66  return *this;
67 }
std::vector< double > a
Definition: orthoPoly1D.H:77
std::vector< double > b
Definition: orthoPoly1D.H:79
Eigen::MatrixXd phiTphiInv
Definition: orthoPoly1D.H:83
Eigen::MatrixXd phi
Definition: orthoPoly1D.H:81

Member Data Documentation

◆ a

std::vector<double> orthoPoly1D::a

Definition at line 77 of file orthoPoly1D.H.

Referenced by ComputeAB(), EvaluateOrthogonal(), operator=(), and orthoPoly1D().

◆ b

std::vector<double> orthoPoly1D::b

Definition at line 79 of file orthoPoly1D.H.

Referenced by ComputeAB(), EvaluateOrthogonal(), operator=(), and orthoPoly1D().

◆ order

int orthoPoly1D::order

Definition at line 73 of file orthoPoly1D.H.

Referenced by ComputeAB(), EvaluateOrthogonals(), operator=(), and orthoPoly1D().

◆ phi

Eigen::MatrixXd orthoPoly1D::phi

Definition at line 81 of file orthoPoly1D.H.

Referenced by ComputePhiTPhiInv(), EvaluateOrthogonals(), operator=(), and orthoPoly1D().

◆ phiTphiInv

Eigen::MatrixXd orthoPoly1D::phiTphiInv

Definition at line 83 of file orthoPoly1D.H.

Referenced by ComputePhiTPhiInv(), operator=(), and orthoPoly1D().


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