MAST3
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
MAST::FEBasis::Evaluation Namespace Reference

Classes

class  FEShapeDerivative
 

Functions

template<typename NodalScalarType , uint_t ElemDim, uint_t SpatialDim>
void compute_detJ (const Eigen::Matrix< NodalScalarType, SpatialDim *ElemDim, Eigen::Dynamic > &dx_dxi, Eigen::Matrix< NodalScalarType, Eigen::Dynamic, 1 > &detJ)
 
template<typename NodalScalarType , uint_t ElemDim, uint_t SpatialDim, typename ContextType >
std::enable_if< ElemDim==SpatialDim &&ElemDim==1, void >::type compute_detJ_side (const ContextType &c, const uint_t s, const Eigen::Matrix< NodalScalarType, ElemDim *SpatialDim, Eigen::Dynamic > &dx_dxi, Eigen::Matrix< NodalScalarType, Eigen::Dynamic, 1 > &detJ)
 
template<typename NodalScalarType , uint_t ElemDim, uint_t SpatialDim, typename ContextType >
std::enable_if< ElemDim==SpatialDim &&ElemDim==2, void >::type compute_detJ_side (const ContextType &c, const uint_t s, const Eigen::Matrix< NodalScalarType, ElemDim *SpatialDim, Eigen::Dynamic > &dx_dxi, Eigen::Matrix< NodalScalarType, Eigen::Dynamic, 1 > &detJ)
 
template<typename NodalScalarType , uint_t ElemDim, uint_t SpatialDim, typename ContextType >
std::enable_if< ElemDim==SpatialDim &&ElemDim==3, void >::type compute_detJ_side (const ContextType &c, const uint_t s, const Eigen::Matrix< NodalScalarType, ElemDim *SpatialDim, Eigen::Dynamic > &dx_dxi, Eigen::Matrix< NodalScalarType, Eigen::Dynamic, 1 > &detJ)
 
template<typename NodalScalarType , uint_t ElemDim, uint_t SpatialDim, typename ContextType >
std::enable_if< ElemDim==SpatialDim &&ElemDim==3, void >::type compute_detJ_side_hex (const ContextType &c, const uint_t s, const Eigen::Matrix< NodalScalarType, ElemDim *SpatialDim, Eigen::Dynamic > &dx_dxi, Eigen::Matrix< NodalScalarType, Eigen::Dynamic, 1 > &detJ)
 
template<typename NodalScalarType , uint_t ElemDim, uint_t SpatialDim, typename ContextType >
std::enable_if< ElemDim==SpatialDim &&ElemDim==2, void >::type compute_detJ_side_quad (const ContextType &c, const uint_t s, const Eigen::Matrix< NodalScalarType, ElemDim *SpatialDim, Eigen::Dynamic > &dx_dxi, Eigen::Matrix< NodalScalarType, Eigen::Dynamic, 1 > &detJ)
 
template<typename NodalScalarType , uint_t ElemDim, uint_t SpatialDim, typename FEBasisType , typename ContextType >
void compute_detJxW (const FEBasisType &fe_basis, const Eigen::Matrix< NodalScalarType, Eigen::Dynamic, 1 > &detJ, Eigen::Matrix< NodalScalarType, Eigen::Dynamic, 1 > &detJxW)
 
template<typename NodalScalarType , uint_t ElemDim, uint_t SpatialDim, typename FEBasisType >
void compute_dphi_dx (const FEBasisType &fe_basis, const Eigen::Matrix< NodalScalarType, ElemDim *SpatialDim, Eigen::Dynamic > &dxi_dx, Eigen::Matrix< NodalScalarType, Eigen::Dynamic, Eigen::Dynamic > &dphi_dx)
 
template<typename NodalScalarType , uint_t ElemDim, uint_t SpatialDim, typename ContextType >
std::enable_if< ElemDim==SpatialDim &&ElemDim==3, void >::type compute_hex_side_tangent_and_normal (const ContextType &c, const uint_t s, const Eigen::Matrix< NodalScalarType, ElemDim *SpatialDim, Eigen::Dynamic > &dx_dxi, Eigen::Matrix< NodalScalarType, SpatialDim, Eigen::Dynamic > &tangent, Eigen::Matrix< NodalScalarType, SpatialDim, Eigen::Dynamic > &normal)
 
template<typename NodalScalarType , uint_t ElemDim, uint_t SpatialDim, typename FEBasisType , typename ContextType >
void compute_Jac (const ContextType &c, const FEBasisType &fe_basis, const Eigen::Matrix< NodalScalarType, SpatialDim, Eigen::Dynamic > &node_coord, Eigen::Matrix< NodalScalarType, SpatialDim *ElemDim, Eigen::Dynamic > &dx_dxi)
 
template<typename NodalScalarType , uint_t ElemDim, uint_t SpatialDim>
std::enable_if< ElemDim==SpatialDim, void >::type compute_Jac_inv (const Eigen::Matrix< NodalScalarType, ElemDim *ElemDim, Eigen::Dynamic > &dx_dxi, Eigen::Matrix< NodalScalarType, ElemDim *ElemDim, Eigen::Dynamic > &dxi_dx)
 
template<typename NodalScalarType , uint_t ElemDim, uint_t SpatialDim, typename ContextType >
std::enable_if< ElemDim==SpatialDim &&ElemDim==2, void >::type compute_quad_side_tangent_and_normal (const ContextType &c, const uint_t s, const Eigen::Matrix< NodalScalarType, ElemDim *SpatialDim, Eigen::Dynamic > &dx_dxi, Eigen::Matrix< NodalScalarType, SpatialDim, Eigen::Dynamic > &tangent, Eigen::Matrix< NodalScalarType, SpatialDim, Eigen::Dynamic > &normal)
 
template<typename NodalScalarType , uint_t ElemDim, uint_t SpatialDim, typename ContextType >
std::enable_if< ElemDim==SpatialDim &&ElemDim==1, void >::type compute_side_tangent_and_normal (const ContextType &c, const uint_t s, const Eigen::Matrix< NodalScalarType, SpatialDim *ElemDim, Eigen::Dynamic > &dx_dxi, Eigen::Matrix< NodalScalarType, SpatialDim, Eigen::Dynamic > &tangent, Eigen::Matrix< NodalScalarType, SpatialDim, Eigen::Dynamic > &normal)
 
template<typename NodalScalarType , uint_t ElemDim, uint_t SpatialDim, typename ContextType >
std::enable_if< ElemDim==SpatialDim &&ElemDim==2, void >::type compute_side_tangent_and_normal (const ContextType &c, const uint_t s, const Eigen::Matrix< NodalScalarType, SpatialDim *ElemDim, Eigen::Dynamic > &dx_dxi, Eigen::Matrix< NodalScalarType, SpatialDim, Eigen::Dynamic > &tangent, Eigen::Matrix< NodalScalarType, SpatialDim, Eigen::Dynamic > &normal)
 
template<typename NodalScalarType , uint_t ElemDim, uint_t SpatialDim, typename ContextType >
std::enable_if< ElemDim==SpatialDim &&ElemDim==3, void >::type compute_side_tangent_and_normal (const ContextType &c, const uint_t s, const Eigen::Matrix< NodalScalarType, SpatialDim *ElemDim, Eigen::Dynamic > &dx_dxi, Eigen::Matrix< NodalScalarType, SpatialDim, Eigen::Dynamic > &tangent, Eigen::Matrix< NodalScalarType, SpatialDim, Eigen::Dynamic > &normal)
 
template<typename NodalScalarType , uint_t ElemDim, uint_t SpatialDim, typename FEBasisType , typename ContextType >
void compute_xyz (const ContextType &c, const FEBasisType &fe_basis, Eigen::Matrix< NodalScalarType, SpatialDim, Eigen::Dynamic > &node_coord, Eigen::Matrix< NodalScalarType, SpatialDim, Eigen::Dynamic > &xyz)
 
void hex_side_Jac_cols (const uint_t s, uint_t &c1, uint_t &c2)
 
uint_t quad_side_Jac_col (uint_t s)
 

Function Documentation

◆ compute_detJ()

template<typename NodalScalarType , uint_t ElemDim, uint_t SpatialDim>
void MAST::FEBasis::Evaluation::compute_detJ ( const Eigen::Matrix< NodalScalarType, SpatialDim *ElemDim, Eigen::Dynamic > &  dx_dxi,
Eigen::Matrix< NodalScalarType, Eigen::Dynamic, 1 > &  detJ 
)
inline

Definition at line 100 of file fe_derivative_evaluation.hpp.

101  {
102 
103  uint_t
104  nq = dx_dxi.cols();
105 
106  detJ = Eigen::Matrix<NodalScalarType, Eigen::Dynamic, 1>::Zero(nq);
107 
108  for (uint_t i=0; i<nq; i++) {
109 
110  Eigen::Map<const typename Eigen::Matrix<NodalScalarType, SpatialDim, ElemDim>>
111  dxdxi(dx_dxi.col(i).data(), SpatialDim, ElemDim);
112 
113  // determinant of dx_dxi
114  detJ(i) = dxdxi.determinant();
115  }
116 }
unsigned int uint_t

◆ compute_detJ_side() [1/3]

template<typename NodalScalarType , uint_t ElemDim, uint_t SpatialDim, typename ContextType >
std::enable_if<ElemDim == SpatialDim && ElemDim == 1, void>::type MAST::FEBasis::Evaluation::compute_detJ_side ( const ContextType &  c,
const uint_t  s,
const Eigen::Matrix< NodalScalarType, ElemDim *SpatialDim, Eigen::Dynamic > &  dx_dxi,
Eigen::Matrix< NodalScalarType, Eigen::Dynamic, 1 > &  detJ 
)
inline

Definition at line 127 of file fe_derivative_evaluation.hpp.

130  {
131 
132  Assert1(dx_dxi.cols() == 1, dx_dxi.cols(), "Only one quadrature point on side of 1D element");
133  detJ = Eigen::Matrix<NodalScalarType, Eigen::Dynamic, 1>::Ones(1);
134 }
#define Assert1(cond, v1, msg)
Definition: exceptions.hpp:143

◆ compute_detJ_side() [2/3]

template<typename NodalScalarType , uint_t ElemDim, uint_t SpatialDim, typename ContextType >
std::enable_if<ElemDim == SpatialDim && ElemDim == 2, void>::type MAST::FEBasis::Evaluation::compute_detJ_side ( const ContextType &  c,
const uint_t  s,
const Eigen::Matrix< NodalScalarType, ElemDim *SpatialDim, Eigen::Dynamic > &  dx_dxi,
Eigen::Matrix< NodalScalarType, Eigen::Dynamic, 1 > &  detJ 
)
inline

Definition at line 311 of file fe_derivative_evaluation.hpp.

314  {
315 
316  if (c.elem_is_quad())
317  compute_detJ_side_quad<NodalScalarType, ElemDim, SpatialDim, ContextType>(c, s, dx_dxi, detJ);
318  else
319  Error(false, "Not implemented for element type.");
320  }
#define Error(cond, msg)
Definition: exceptions.hpp:166

◆ compute_detJ_side() [3/3]

template<typename NodalScalarType , uint_t ElemDim, uint_t SpatialDim, typename ContextType >
std::enable_if<ElemDim == SpatialDim && ElemDim == 3, void>::type MAST::FEBasis::Evaluation::compute_detJ_side ( const ContextType &  c,
const uint_t  s,
const Eigen::Matrix< NodalScalarType, ElemDim *SpatialDim, Eigen::Dynamic > &  dx_dxi,
Eigen::Matrix< NodalScalarType, Eigen::Dynamic, 1 > &  detJ 
)
inline

Definition at line 331 of file fe_derivative_evaluation.hpp.

334  {
335 
336  if (c.elem_is_hex())
337  compute_detJ_side_hex<NodalScalarType, ElemDim, SpatialDim, ContextType>(c, s, dx_dxi, detJ);
338  else
339  Error(false, "Not implemented for element type.");
340  }
#define Error(cond, msg)
Definition: exceptions.hpp:166

◆ compute_detJ_side_hex()

template<typename NodalScalarType , uint_t ElemDim, uint_t SpatialDim, typename ContextType >
std::enable_if<ElemDim == SpatialDim && ElemDim == 3, void>::type MAST::FEBasis::Evaluation::compute_detJ_side_hex ( const ContextType &  c,
const uint_t  s,
const Eigen::Matrix< NodalScalarType, ElemDim *SpatialDim, Eigen::Dynamic > &  dx_dxi,
Eigen::Matrix< NodalScalarType, Eigen::Dynamic, 1 > &  detJ 
)
inline

Definition at line 276 of file fe_derivative_evaluation.hpp.

279  {
280 
281  Assert0(c.elem_is_hex(), "Element must be a hex");
282 
283  uint_t
284  nq = dx_dxi.cols(),
285  c1 = -1,
286  c2 = -1;
287 
289 
290  detJ = Eigen::Matrix<NodalScalarType, Eigen::Dynamic, 1>::Zero(nq);
291 
292  for (uint_t i=0; i<nq; i++) {
293 
294  Eigen::Map<const typename Eigen::Matrix<NodalScalarType, 3, 3>>
295  dxdxi(dx_dxi.col(i).data(), SpatialDim, ElemDim);
296 
297  // determinant of dx_dxi
298  detJ(i) = dxdxi.col(c1).cross(dxdxi.col(c2)).norm();
299  }
300 }
void hex_side_Jac_cols(const uint_t s, uint_t &c1, uint_t &c2)
#define Assert0(cond, msg)
Definition: exceptions.hpp:134
unsigned int uint_t

◆ compute_detJ_side_quad()

template<typename NodalScalarType , uint_t ElemDim, uint_t SpatialDim, typename ContextType >
std::enable_if<ElemDim == SpatialDim && ElemDim == 2, void>::type MAST::FEBasis::Evaluation::compute_detJ_side_quad ( const ContextType &  c,
const uint_t  s,
const Eigen::Matrix< NodalScalarType, ElemDim *SpatialDim, Eigen::Dynamic > &  dx_dxi,
Eigen::Matrix< NodalScalarType, Eigen::Dynamic, 1 > &  detJ 
)
inline

Definition at line 244 of file fe_derivative_evaluation.hpp.

247  {
248 
249  Assert0(c.elem_is_quad(), "Element must be a quadrilateral");
250 
251  uint_t
252  nq = dx_dxi.cols(),
254 
255  detJ = Eigen::Matrix<NodalScalarType, Eigen::Dynamic, 1>::Zero(nq);
256 
257  for (uint_t i=0; i<nq; i++) {
258 
259  Eigen::Map<const typename Eigen::Matrix<NodalScalarType, 2, 2>>
260  dxdxi(dx_dxi.col(i).data(), SpatialDim, ElemDim);
261 
262  // determinant of dx_dxi
263  detJ(i) = dxdxi.col(col).norm();
264  }
265 }
#define Assert0(cond, msg)
Definition: exceptions.hpp:134
unsigned int uint_t

◆ compute_detJxW()

template<typename NodalScalarType , uint_t ElemDim, uint_t SpatialDim, typename FEBasisType , typename ContextType >
void MAST::FEBasis::Evaluation::compute_detJxW ( const FEBasisType &  fe_basis,
const Eigen::Matrix< NodalScalarType, Eigen::Dynamic, 1 > &  detJ,
Eigen::Matrix< NodalScalarType, Eigen::Dynamic, 1 > &  detJxW 
)
inline

Definition at line 519 of file fe_derivative_evaluation.hpp.

521  {
522 
523  Assert2(fe_basis.n_q_points() == detJ.rows(),
524  fe_basis.n_q_points(), detJ.rows(),
525  "Incompatible number of quadrature points of detJ and FEBasis.");
526 
527  uint_t
528  nq = detJ.rows();
529 
530  detJxW = Eigen::Matrix<NodalScalarType, Eigen::Dynamic, 1>::Zero(nq);
531 
532  for (uint_t i=0; i<nq; i++) {
533 
534  // determinant times weight
535  detJxW(i) = detJ(i) * fe_basis.qp_weight(i);
536  }
537 }
unsigned int uint_t
#define Assert2(cond, v1, v2, msg)
Definition: exceptions.hpp:152

◆ compute_dphi_dx()

template<typename NodalScalarType , uint_t ElemDim, uint_t SpatialDim, typename FEBasisType >
void MAST::FEBasis::Evaluation::compute_dphi_dx ( const FEBasisType &  fe_basis,
const Eigen::Matrix< NodalScalarType, ElemDim *SpatialDim, Eigen::Dynamic > &  dxi_dx,
Eigen::Matrix< NodalScalarType, Eigen::Dynamic, Eigen::Dynamic > &  dphi_dx 
)
inline

Definition at line 578 of file fe_derivative_evaluation.hpp.

580  {
581 
582  uint_t
583  nq = fe_basis.n_q_points(),
584  n_basis = fe_basis.n_basis();
585 
586  Assert2(dxi_dx.cols() == nq,
587  dxi_dx.cols(), nq,
588  "Incompatible quadrature points in FEBasis and dxi_dx.");
589  Assert2(dxi_dx.rows() == ElemDim*SpatialDim,
590  dxi_dx.rows(), ElemDim*SpatialDim,
591  "Incompatible rows in dxi_dx.");
592 
593  dphi_dx = Eigen::Matrix<NodalScalarType, Eigen::Dynamic, Eigen::Dynamic>::Zero(SpatialDim*n_basis, nq);
594 
595  for (uint_t i=0; i<nq; i++) {
596 
597  // quadrature point spatial coordinate derivatives dx/dxi
598  Eigen::Map<const typename Eigen::Matrix<NodalScalarType, ElemDim, SpatialDim>>
599  dxidx (dxi_dx.col(i).data(), ElemDim, SpatialDim);
600  Eigen::Map<typename Eigen::Matrix<NodalScalarType, Eigen::Dynamic, SpatialDim>>
601  dphidx(dphi_dx.col(i).data(), n_basis, SpatialDim);
602 
603  for (uint_t l=0; l<n_basis; l++)
604  for (uint_t j=0; j<SpatialDim; j++)
605  for (uint_t k=0; k<ElemDim; k++)
606  dphidx(l, j) += fe_basis.dphi_dxi(i, l, k) * dxidx(k, j);
607  }
608 }
unsigned int uint_t
#define Assert2(cond, v1, v2, msg)
Definition: exceptions.hpp:152

◆ compute_hex_side_tangent_and_normal()

template<typename NodalScalarType , uint_t ElemDim, uint_t SpatialDim, typename ContextType >
std::enable_if<ElemDim == SpatialDim && ElemDim == 3, void>::type MAST::FEBasis::Evaluation::compute_hex_side_tangent_and_normal ( const ContextType &  c,
const uint_t  s,
const Eigen::Matrix< NodalScalarType, ElemDim *SpatialDim, Eigen::Dynamic > &  dx_dxi,
Eigen::Matrix< NodalScalarType, SpatialDim, Eigen::Dynamic > &  tangent,
Eigen::Matrix< NodalScalarType, SpatialDim, Eigen::Dynamic > &  normal 
)
inline

Definition at line 429 of file fe_derivative_evaluation.hpp.

433  {
434 
435  Assert0(c.elem_is_hex(), "Element must be a hexagon");
436 
437  uint_t
438  nq = dx_dxi.cols(),
439  c1 = -1,
440  c2 = -1;
441 
443 
444  tangent = Eigen::Matrix<NodalScalarType, SpatialDim, Eigen::Dynamic>::Zero(SpatialDim, nq);
445  normal = Eigen::Matrix<NodalScalarType, SpatialDim, Eigen::Dynamic>::Zero(SpatialDim, nq);
446 
447  Eigen::Matrix<NodalScalarType, 3, 1>
448  dx = Eigen::Matrix<NodalScalarType, 3, 1>::Zero(3);
449 
450  for (uint_t i=0; i<nq; i++) {
451 
452  Eigen::Map<const typename Eigen::Matrix<NodalScalarType, 3, 3>>
453  dxdxi(dx_dxi.col(i).data(), 3, 3);
454 
455  // normal is dx1 x dx2
456  // dn = | i j k |
457  // | dx1 dy1 dz1 |
458  // | dx2 dy2 dz2 |
459  normal.col(i) = dxdxi.col(c1).cross(dxdxi.col(c2));
460  normal.col(i) /= normal.col(i).norm();
461  }
462 }
void hex_side_Jac_cols(const uint_t s, uint_t &c1, uint_t &c2)
#define Assert0(cond, msg)
Definition: exceptions.hpp:134
unsigned int uint_t

◆ compute_Jac()

template<typename NodalScalarType , uint_t ElemDim, uint_t SpatialDim, typename FEBasisType , typename ContextType >
void MAST::FEBasis::Evaluation::compute_Jac ( const ContextType &  c,
const FEBasisType &  fe_basis,
const Eigen::Matrix< NodalScalarType, SpatialDim, Eigen::Dynamic > &  node_coord,
Eigen::Matrix< NodalScalarType, SpatialDim *ElemDim, Eigen::Dynamic > &  dx_dxi 
)
inline

Definition at line 70 of file fe_derivative_evaluation.hpp.

73  {
74 
75  uint_t
76  nq = fe_basis.n_q_points(),
77  n_nodes = c.n_nodes();
78 
79  dx_dxi = Eigen::Matrix<NodalScalarType, SpatialDim*ElemDim, Eigen::Dynamic>::Zero(SpatialDim*ElemDim, nq);
80 
81  // quadrature point spatial coordinate derivatives dx/dxi
82  for (uint_t i=0; i<nq; i++) {
83 
84  Eigen::Map<typename Eigen::Matrix<NodalScalarType, SpatialDim, ElemDim>>
85  dxdxi(dx_dxi.col(i).data(), SpatialDim, ElemDim);
86 
87  for (uint_t l=0; l<n_nodes; l++)
88  for (uint_t j=0; j<SpatialDim; j++)
89  for (uint_t k=0; k<ElemDim; k++)
90  dxdxi(j, k) += fe_basis.dphi_dxi(i, l, k) * node_coord(j, l);
91  }
92 
93 }
unsigned int uint_t

◆ compute_Jac_inv()

template<typename NodalScalarType , uint_t ElemDim, uint_t SpatialDim>
std::enable_if<ElemDim == SpatialDim, void>::type MAST::FEBasis::Evaluation::compute_Jac_inv ( const Eigen::Matrix< NodalScalarType, ElemDim *ElemDim, Eigen::Dynamic > &  dx_dxi,
Eigen::Matrix< NodalScalarType, ElemDim *ElemDim, Eigen::Dynamic > &  dxi_dx 
)
inline

Definition at line 548 of file fe_derivative_evaluation.hpp.

549  {
550 
551  uint_t
552  nq = dx_dxi.cols();
553 
554  dxi_dx = Eigen::Matrix<NodalScalarType, ElemDim*ElemDim, Eigen::Dynamic>::Zero(ElemDim*ElemDim, nq);
555 
556  for (uint_t i=0; i<nq; i++) {
557 
558  // quadrature point spatial coordinate derivatives dx/dxi
559  Eigen::Map<const typename Eigen::Matrix<NodalScalarType, ElemDim, ElemDim>>
560  dxdxi(dx_dxi.col(i).data(), ElemDim, ElemDim);
561  Eigen::Map<typename Eigen::Matrix<NodalScalarType, ElemDim, ElemDim>>
562  dxidx(dxi_dx.col(i).data(), ElemDim, ElemDim);
563 
564  // compute dx/dxi
565  dxidx = dxdxi.inverse();
566  }
567 }
unsigned int uint_t

◆ compute_quad_side_tangent_and_normal()

template<typename NodalScalarType , uint_t ElemDim, uint_t SpatialDim, typename ContextType >
std::enable_if<ElemDim == SpatialDim && ElemDim == 2, void>::type MAST::FEBasis::Evaluation::compute_quad_side_tangent_and_normal ( const ContextType &  c,
const uint_t  s,
const Eigen::Matrix< NodalScalarType, ElemDim *SpatialDim, Eigen::Dynamic > &  dx_dxi,
Eigen::Matrix< NodalScalarType, SpatialDim, Eigen::Dynamic > &  tangent,
Eigen::Matrix< NodalScalarType, SpatialDim, Eigen::Dynamic > &  normal 
)
inline

Definition at line 377 of file fe_derivative_evaluation.hpp.

381  {
382 
383  Assert0(c.elem_is_quad(), "Element must be a quadrilateral");
384 
385  uint_t
386  nq = dx_dxi.cols(),
388 
389  tangent = Eigen::Matrix<NodalScalarType, SpatialDim, Eigen::Dynamic>::Zero(SpatialDim, nq);
390  normal = Eigen::Matrix<NodalScalarType, SpatialDim, Eigen::Dynamic>::Zero(SpatialDim, nq);
391 
392  // for bottom and right edges, the tangent is d{x, y}/dxi and d{x, y}/deta.
393  // for top and left edges, the tangent is -d{x, y}/dxi and -d{x, y}/deta.
394  real_t
395  v = s>1?-1.:1;
396 
397  Eigen::Matrix<NodalScalarType, 2, 1>
398  dx = Eigen::Matrix<NodalScalarType, 2, 1>::Zero(2);
399 
400  for (uint_t i=0; i<nq; i++) {
401 
402  Eigen::Map<const typename Eigen::Matrix<NodalScalarType, 2, 2>>
403  dxdxi(dx_dxi.col(i).data(), 2, 2);
404 
405  // tangent
406  dx = v * dxdxi.col(col);
407  dx /= dx.norm();
408  tangent.col(i) = dx;
409 
410  // normal is tangent cross k_hat
411  // dn = | i j k | = i ty - j tx
412  // | tx ty 0 |
413  // | 0 0 1 |
414  normal(0, i) = dx(1);
415  normal(1, i) = -dx(0);
416  }
417 }
#define Assert0(cond, msg)
Definition: exceptions.hpp:134
unsigned int uint_t
double real_t

◆ compute_side_tangent_and_normal() [1/3]

template<typename NodalScalarType , uint_t ElemDim, uint_t SpatialDim, typename ContextType >
std::enable_if<ElemDim == SpatialDim && ElemDim == 1, void>::type MAST::FEBasis::Evaluation::compute_side_tangent_and_normal ( const ContextType &  c,
const uint_t  s,
const Eigen::Matrix< NodalScalarType, SpatialDim *ElemDim, Eigen::Dynamic > &  dx_dxi,
Eigen::Matrix< NodalScalarType, SpatialDim, Eigen::Dynamic > &  tangent,
Eigen::Matrix< NodalScalarType, SpatialDim, Eigen::Dynamic > &  normal 
)
inline

Definition at line 352 of file fe_derivative_evaluation.hpp.

356  {
357 
358  // side of 1D is a 0-d element, so shoudl have a single point
359  Assert1(dx_dxi.cols() == 1, dx_dxi.cols(), "Only one quadrature point on side of 1D element");
360 
361  normal = Eigen::Matrix<NodalScalarType, 1, Eigen::Dynamic>::Zero(1, 1);
362  tangent = Eigen::Matrix<NodalScalarType, 1, Eigen::Dynamic>::Zero(1, 1);
363 
364  // left side normal is -1 and right side normal is +1
365  normal(0, 0) = s==0?-1.:1.;
366 }
#define Assert1(cond, v1, msg)
Definition: exceptions.hpp:143

◆ compute_side_tangent_and_normal() [2/3]

template<typename NodalScalarType , uint_t ElemDim, uint_t SpatialDim, typename ContextType >
std::enable_if<ElemDim == SpatialDim && ElemDim == 2, void>::type MAST::FEBasis::Evaluation::compute_side_tangent_and_normal ( const ContextType &  c,
const uint_t  s,
const Eigen::Matrix< NodalScalarType, SpatialDim *ElemDim, Eigen::Dynamic > &  dx_dxi,
Eigen::Matrix< NodalScalarType, SpatialDim, Eigen::Dynamic > &  tangent,
Eigen::Matrix< NodalScalarType, SpatialDim, Eigen::Dynamic > &  normal 
)
inline

Definition at line 473 of file fe_derivative_evaluation.hpp.

477  {
478 
479  if (c.elem_is_quad())
481  <NodalScalarType, ElemDim, SpatialDim, ContextType>
482  (c, s, dx_dxi, tangent, normal);
483  else
484  Error(false, "Not implemented for element type.");
485 }
#define Error(cond, msg)
Definition: exceptions.hpp:166
std::enable_if< ElemDim==SpatialDim &&ElemDim==2, void >::type compute_quad_side_tangent_and_normal(const ContextType &c, const uint_t s, const Eigen::Matrix< NodalScalarType, ElemDim *SpatialDim, Eigen::Dynamic > &dx_dxi, Eigen::Matrix< NodalScalarType, SpatialDim, Eigen::Dynamic > &tangent, Eigen::Matrix< NodalScalarType, SpatialDim, Eigen::Dynamic > &normal)

◆ compute_side_tangent_and_normal() [3/3]

template<typename NodalScalarType , uint_t ElemDim, uint_t SpatialDim, typename ContextType >
std::enable_if<ElemDim == SpatialDim && ElemDim == 3, void>::type MAST::FEBasis::Evaluation::compute_side_tangent_and_normal ( const ContextType &  c,
const uint_t  s,
const Eigen::Matrix< NodalScalarType, SpatialDim *ElemDim, Eigen::Dynamic > &  dx_dxi,
Eigen::Matrix< NodalScalarType, SpatialDim, Eigen::Dynamic > &  tangent,
Eigen::Matrix< NodalScalarType, SpatialDim, Eigen::Dynamic > &  normal 
)
inline

Definition at line 497 of file fe_derivative_evaluation.hpp.

501  {
502 
503  if (c.elem_is_hex())
505  <NodalScalarType, ElemDim, SpatialDim, ContextType>
506  (c, s, dx_dxi, tangent, normal);
507  else
508  Error(false, "Not implemented for element type.");
509 }
#define Error(cond, msg)
Definition: exceptions.hpp:166
std::enable_if< ElemDim==SpatialDim &&ElemDim==3, void >::type compute_hex_side_tangent_and_normal(const ContextType &c, const uint_t s, const Eigen::Matrix< NodalScalarType, ElemDim *SpatialDim, Eigen::Dynamic > &dx_dxi, Eigen::Matrix< NodalScalarType, SpatialDim, Eigen::Dynamic > &tangent, Eigen::Matrix< NodalScalarType, SpatialDim, Eigen::Dynamic > &normal)

◆ compute_xyz()

template<typename NodalScalarType , uint_t ElemDim, uint_t SpatialDim, typename FEBasisType , typename ContextType >
void MAST::FEBasis::Evaluation::compute_xyz ( const ContextType &  c,
const FEBasisType &  fe_basis,
Eigen::Matrix< NodalScalarType, SpatialDim, Eigen::Dynamic > &  node_coord,
Eigen::Matrix< NodalScalarType, SpatialDim, Eigen::Dynamic > &  xyz 
)
inline

Definition at line 38 of file fe_derivative_evaluation.hpp.

41  {
42 
43  uint_t
44  nq = fe_basis.n_q_points(),
45  n_nodes = c.n_nodes();
46 
47  node_coord = Eigen::Matrix<NodalScalarType, SpatialDim, Eigen::Dynamic>::Zero(SpatialDim, n_nodes);
48  xyz = Eigen::Matrix<NodalScalarType, SpatialDim, Eigen::Dynamic>::Zero(SpatialDim, nq);
49 
50  // get the nodal locations
51  for (uint_t i=0; i<n_nodes; i++)
52  for (uint_t j=0; j<SpatialDim; j++)
53  node_coord(j, i) = c.nodal_coord(i, j);
54 
55  // quadrature point coordinates
56  for (uint_t i=0; i<nq; i++)
57  for (uint_t k=0; k<n_nodes; k++)
58  for (uint_t j=0; j<SpatialDim; j++)
59  xyz(j, i) += fe_basis.phi(i, k) * node_coord(j, k);
60 }
unsigned int uint_t

◆ hex_side_Jac_cols()

void MAST::FEBasis::Evaluation::hex_side_Jac_cols ( const uint_t  s,
uint_t c1,
uint_t c2 
)
inline

Definition at line 163 of file fe_derivative_evaluation.hpp.

165  {
166 
167  /*
168  * HEX8: 7 6
169  * o--------z
170  * /: /| zeta
171  * / : / | ^ eta (into page)
172  * 4 / : 5 / | | /
173  * o--------o | |/
174  * | o....|...o 2 o---> xi
175  * | .3 | /
176  * | . | /
177  * |. |/
178  * o--------o
179  * 0 1
180  *
181  * libMesh side numbering:
182  * {0, 3, 2, 1}, // Side 0
183  * {0, 1, 5, 4}, // Side 1
184  * {1, 2, 6, 5}, // Side 2
185  * {2, 3, 7, 6}, // Side 3
186  * {3, 0, 4, 7}, // Side 4
187  * {4, 5, 6, 7} // Side 5
188  */
189 
190  // identify row of the Jacobian matrix that will be used to compute
191  // the size
192  switch (s) {
193 
194  case 0: {
195  c1 = 1; // { dx/deta, dy/deta, dz/deta}
196  c2 = 0; // { dx/xi, dy/xi, dz/xi}
197  }
198  break;
199 
200  case 1: {
201  c1 = 0; // { dx/xi, dy/xi, dz/xi}
202  c2 = 2; // {dx/dzeta, dy/dzeta, dz/dzeta}
203  }
204  break;
205 
206  case 2: {
207  c1 = 1; // { dx/deta, dy/deta, dz/deta}
208  c2 = 2; // {dx/dzeta, dy/dzeta, dz/dzeta}
209  }
210  break;
211 
212  case 3: {
213  c1 = 2; // {dx/dzeta, dy/dzeta, dz/dzeta}
214  c2 = 0; // { dx/xi, dy/xi, dz/xi}
215  }
216  break;
217 
218  case 4: {
219  c1 = 2; // {dx/dzeta, dy/dzeta, dz/dzeta}
220  c2 = 1; // { dx/deta, dy/deta, dz/deta}
221  }
222  break;
223 
224  case 5: {
225  c1 = 0; // { dx/deta, dy/deta, dz/deta}
226  c2 = 1; // { dx/xi, dy/xi, dz/xi}
227  }
228  break;
229 
230  default:
231  Error(false, "Invalid side number for hex");
232  }
233 }
#define Error(cond, msg)
Definition: exceptions.hpp:166

◆ quad_side_Jac_col()

uint_t MAST::FEBasis::Evaluation::quad_side_Jac_col ( uint_t  s)
inline

Definition at line 138 of file fe_derivative_evaluation.hpp.

138  {
139 
140  // identify row of the Jacobian matrix that will be used to compute
141  // the size
142  switch (s) {
143  case 0:
144  case 2:
145  return 0; // (dx/dxi, dy/dxi)
146  break;
147 
148  case 1:
149  case 3:
150  return 1; // (dx/deta, dy/deta)
151  break;
152 
153  default:
154  Error(false, "Invalid side number for quad");
155  }
156 
157  // should not get here
158  return -1;
159 }
#define Error(cond, msg)
Definition: exceptions.hpp:166