OGS
detail Namespace Reference

Classes

struct  EigenMatrixType
 
struct  EigenMatrixType< 0, 1 >
 
struct  EigenMatrixType< 0, M >
 
struct  EigenMatrixType< N, 1 >
 

Functions

void rotateToLocal (const MeshLib::RotationMatrix &matR2local, std::vector< MathLib::Point3d > &points)
 rotate points to local coordinates
 
MeshLib::RotationMatrix getRotationMatrixToGlobal (const unsigned element_dimension, const unsigned global_dim, const std::vector< MathLib::Point3d > &points)
 
template<typename Solutions , typename Vector >
void applyKnownSolutions (std::vector< Solutions > const *const known_solutions, Vector &x)
 Applies known solutions to the solution vector x.
 

Function Documentation

◆ applyKnownSolutions()

template<typename Solutions , typename Vector >
void detail::applyKnownSolutions ( std::vector< Solutions > const *const known_solutions,
Vector & x )

Applies known solutions to the solution vector x.

Definition at line 25 of file TimeDiscretizedODESystem.cpp.

27{
28 if (!known_solutions)
29 {
30 return;
31 }
32
33 for (auto const& bc : *known_solutions)
34 {
35 for (std::size_t i = 0; i < bc.ids.size(); ++i)
36 {
37 // TODO that might have bad performance for some Vector types, e.g.,
38 // PETSc.
39 MathLib::setVector(x, bc.ids[i], bc.values[i]);
40 }
41 }
43}
void finalizeAssembly(PETScMatrix &A)
Definition LinAlg.cpp:198
void setVector(PETScVector &v, std::initializer_list< double > values)

References MathLib::LinAlg::finalizeAssembly(), and MathLib::setVector().

Referenced by NumLib::TimeDiscretizedODESystem< ODESystemTag::FirstOrderImplicitQuasilinear, NonlinearSolverTag::Newton >::applyKnownSolutions(), and NumLib::TimeDiscretizedODESystem< ODESystemTag::FirstOrderImplicitQuasilinear, NonlinearSolverTag::Picard >::applyKnownSolutions().

◆ getRotationMatrixToGlobal()

MeshLib::RotationMatrix detail::getRotationMatrixToGlobal ( const unsigned element_dimension,
const unsigned global_dim,
const std::vector< MathLib::Point3d > & points )

get a rotation matrix to the global coordinates it computes R in x=R*x' where x is original coordinates and x' is local coordinates

Definition at line 34 of file ElementCoordinatesMappingLocal.cpp.

38{
39 Eigen::Matrix3d matR;
40 // compute R in x=R*x' where x are original coordinates and x' are local
41 // coordinates
42 if (element_dimension == 1)
43 {
44 Eigen::Vector3d const xx =
45 (points[1].asEigenVector3d() - points[0].asEigenVector3d())
46 .normalized();
47 if (global_dim == 2)
48 {
50 }
51 else
52 {
54 }
55 matR.transposeInPlace();
56 }
57 else if (global_dim == 3 && element_dimension == 2)
58 {
59 // get plane normal
60 auto const [plane_normal, d] = GeoLib::getNewellPlane(points);
61 // compute a rotation matrix to XY
62 matR = GeoLib::computeRotationMatrixToXY(plane_normal);
63 // set a transposed matrix
64 matR.transposeInPlace();
65 }
66 return matR;
67}
Eigen::Matrix3d compute3DRotationMatrixToX(Eigen::Vector3d const &v)
Eigen::Matrix3d computeRotationMatrixToXY(Eigen::Vector3d const &n)
std::pair< Eigen::Vector3d, double > getNewellPlane(InputIterator pnts_begin, InputIterator pnts_end)
Eigen::Matrix3d compute2DRotationMatrixToX(Eigen::Vector3d const &v)

References GeoLib::compute2DRotationMatrixToX(), GeoLib::compute3DRotationMatrixToX(), GeoLib::computeRotationMatrixToXY(), and GeoLib::getNewellPlane().

Referenced by MeshLib::ElementCoordinatesMappingLocal::ElementCoordinatesMappingLocal().

◆ rotateToLocal()

void detail::rotateToLocal ( const MeshLib::RotationMatrix & matR2local,
std::vector< MathLib::Point3d > & points )

rotate points to local coordinates

Definition at line 24 of file ElementCoordinatesMappingLocal.cpp.

26{
27 std::transform(points.begin(), points.end(), points.begin(),
28 [&matR2local](auto const& pnt) { return matR2local * pnt; });
29}

Referenced by MeshLib::ElementCoordinatesMappingLocal::ElementCoordinatesMappingLocal().