PhysicsBasedAnimationToolkit 0.0.10
Cross-platform C++20 library of algorithms and data structures commonly used in computer graphics research on physically-based simulation.
|
FEM gradient operator. More...
#include "Concepts.h"
#include "ShapeFunctions.h"
#include "pbat/Aliases.h"
#include "pbat/common/Eigen.h"
#include "pbat/profiling/Profiling.h"
#include <exception>
#include <fmt/core.h>
#include <tbb/parallel_for.h>
Go to the source code of this file.
Namespaces | |
namespace | pbat |
The main namespace of the library. | |
namespace | pbat::fem |
Finite Element Method (FEM) | |
Functions | |
template<CElement TElement, int Dims, class TDerivedE, class TDerivedeg, class TDerivedGNeg, class TDerivedIn, class TDerivedOut> | |
void | pbat::fem::GemmGradient (Eigen::DenseBase< TDerivedE > const &E, Eigen::Index nNodes, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::MatrixBase< TDerivedGNeg > const &GNeg, Eigen::MatrixBase< TDerivedIn > const &X, Eigen::DenseBase< TDerivedOut > &Y) |
Compute gradient matrix-vector multiply \( Y += \mathbf{G} X \). | |
template<CMesh TMesh, class TDerivedeg, class TDerivedGNeg, class TDerivedIn, class TDerivedOut> | |
void | pbat::fem::GemmGradient (TMesh const &mesh, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::MatrixBase< TDerivedGNeg > const &GNeg, Eigen::MatrixBase< TDerivedIn > const &X, Eigen::DenseBase< TDerivedOut > &Y) |
Compute gradient matrix-vector multiply \( Y += \mathbf{G} X \) using mesh. | |
template<CElement TElement, int Dims, Eigen::StorageOptions Options, class TDerivedE, class TDerivedeg, class TDerivedGNeg> | |
auto | pbat::fem::GradientMatrix (Eigen::DenseBase< TDerivedE > const &E, Eigen::Index nNodes, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::MatrixBase< TDerivedGNeg > const &GNeg) -> Eigen::SparseMatrix< typename TDerivedGNeg::Scalar, Options, typename TDerivedE::Scalar > |
Construct the gradient operator's sparse matrix representation. | |
template<Eigen::StorageOptions Options, CMesh TMesh, class TDerivedeg, class TDerivedGNeg> | |
auto | pbat::fem::GradientMatrix (TMesh const &mesh, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::MatrixBase< TDerivedGNeg > const &GNeg) -> Eigen::SparseMatrix< typename TDerivedGNeg::Scalar, Options, typename TMesh::IndexType > |
Construct the gradient operator's sparse matrix representation using mesh. | |
FEM gradient operator.
Given a function \( u(X) \) discretized at mesh nodes \( X_i \),
\[\nabla u(X) = \sum_j u_j \nabla \phi_j(X) \; \forall \; e \in E \]
where \( \nabla \phi_j(X) \) is the gradient of the shape function \( \phi_j(X) \) at element \( e \). The gradient operator \( \mathbf{G} \in \mathbb{R}^{d|Q| \times n} \) maps the nodal values to the gradient at quadrature points \( Q \) in dimensions \( d \).
This file provides functions to compute gradient operator related quantities.