PolyFEM
Loading...
Searching...
No Matches
polyfem::assembler::StokesVelocity Class Reference

#include <Stokes.hpp>

Inheritance diagram for polyfem::assembler::StokesVelocity:
[legend]
Collaboration diagram for polyfem::assembler::StokesVelocity:
[legend]

Public Member Functions

 StokesVelocity ()
 
VectorNd compute_rhs (const AutodiffHessianPt &pt) const override
 
Eigen::Matrix< double, Eigen::Dynamic, 1, 0, 9, 1 > assemble (const LinearAssemblerData &data) const override
 local assembly function that defines the bilinear form (LHS) computes and returns a single local stiffness value
 
void add_multimaterial (const int index, const json &params, const Units &units) override
 
const GenericMatParamviscosity () const
 
virtual std::string name () const override
 
std::map< std::string, ParamFuncparameters () const override
 
bool is_fluid () const override
 
bool is_tensor () const override
 
void assemble (const bool is_volume, const int n_basis, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const AssemblyValsCache &cache, const double t, StiffnessMatrix &stiffness, const bool is_mass=false) const override
 assembles the stiffness matrix for the given basis the bilinear form (local assembler) is encoded by the overloaded assemble (see below) function that the subclass (eg Laplacian) defines sets stiffness and modifies cache if it has not already been computed
 
virtual Eigen::Matrix< double, Eigen::Dynamic, 1, 0, 9, 1 > assemble (const LinearAssemblerData &data) const=0
 local assembly function that defines the bilinear form (LHS) computes and returns a single local stiffness value
 
- Public Member Functions inherited from polyfem::assembler::LinearAssembler
 LinearAssembler ()
 
virtual ~LinearAssembler ()=default
 
virtual bool is_linear () const override
 
- Public Member Functions inherited from polyfem::assembler::Assembler
virtual ~Assembler ()=default
 
int size () const
 
virtual void set_size (const int size)
 
virtual double assemble_energy (const bool is_volume, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const AssemblyValsCache &cache, const double t, const double dt, const Eigen::MatrixXd &displacement, const Eigen::MatrixXd &displacement_prev) const
 
virtual Eigen::VectorXd assemble_energy_per_element (const bool is_volume, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const AssemblyValsCache &cache, const double t, const double dt, const Eigen::MatrixXd &displacement, const Eigen::MatrixXd &displacement_prev) const
 
virtual void assemble_gradient (const bool is_volume, const int n_basis, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const AssemblyValsCache &cache, const double t, const double dt, const Eigen::MatrixXd &displacement, const Eigen::MatrixXd &displacement_prev, Eigen::MatrixXd &rhs) const
 
virtual void assemble_hessian (const bool is_volume, const int n_basis, const bool project_to_psd, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const AssemblyValsCache &cache, const double t, const double dt, const Eigen::MatrixXd &displacement, const Eigen::MatrixXd &displacement_prev, utils::MatrixCache &mat_cache, StiffnessMatrix &grad) const
 
virtual void compute_scalar_value (const OutputData &data, std::vector< NamedMatrix > &result) const
 
virtual void compute_tensor_value (const OutputData &data, std::vector< NamedMatrix > &result) const
 
virtual void compute_stiffness_value (const double t, const assembler::ElementAssemblyValues &vals, const Eigen::MatrixXd &local_pts, const Eigen::MatrixXd &displacement, Eigen::MatrixXd &tensor) const
 
virtual void compute_dstress_dmu_dlambda (const OptAssemblerData &data, Eigen::MatrixXd &dstress_dmu, Eigen::MatrixXd &dstress_dlambda) const
 
virtual void compute_stress_grad_multiply_mat (const OptAssemblerData &data, const Eigen::MatrixXd &mat, Eigen::MatrixXd &stress, Eigen::MatrixXd &result) const
 
virtual void compute_stress_grad_multiply_stress (const OptAssemblerData &data, Eigen::MatrixXd &stress, Eigen::MatrixXd &result) const
 
virtual void compute_stress_grad_multiply_vect (const OptAssemblerData &data, const Eigen::MatrixXd &vect, Eigen::MatrixXd &stress, Eigen::MatrixXd &result) const
 
virtual void compute_stress_grad (const OptAssemblerData &data, const Eigen::MatrixXd &prev_grad_u_i, Eigen::MatrixXd &stress, Eigen::MatrixXd &result) const
 
virtual void compute_stress_prev_grad (const OptAssemblerData &data, const Eigen::MatrixXd &prev_grad_u_i, Eigen::MatrixXd &result) const
 
virtual Eigen::Matrix< AutodiffScalarGrad, Eigen::Dynamic, 1, 0, 3, 1 > kernel (const int dim, const AutodiffGradPt &rvect, const AutodiffScalarGrad &r) const
 
void set_materials (const std::vector< int > &body_ids, const json &body_params, const Units &units)
 
virtual void update_lame_params (const Eigen::MatrixXd &lambdas, const Eigen::MatrixXd &mus)
 
virtual bool is_solution_displacement () const
 

Private Attributes

GenericMatParam viscosity_
 

Additional Inherited Members

- Public Types inherited from polyfem::assembler::Assembler
typedef std::pair< std::string, Eigen::MatrixXd > NamedMatrix
 
typedef std::function< double(const RowVectorNd &, const RowVectorNd &, double, int)> ParamFunc
 
- Protected Attributes inherited from polyfem::assembler::Assembler
int size_ = -1
 

Detailed Description

Definition at line 11 of file Stokes.hpp.

Constructor & Destructor Documentation

◆ StokesVelocity()

polyfem::assembler::StokesVelocity::StokesVelocity ( )

Definition at line 5 of file Stokes.cpp.

Member Function Documentation

◆ add_multimaterial()

void polyfem::assembler::StokesVelocity::add_multimaterial ( const int  index,
const json params,
const Units units 
)
overridevirtual

Reimplemented from polyfem::assembler::Assembler.

Definition at line 9 of file Stokes.cpp.

References polyfem::assembler::GenericMatParam::add_multimaterial(), polyfem::assembler::Assembler::size(), polyfem::Units::viscosity(), and viscosity_.

Here is the call graph for this function:

◆ assemble() [1/3]

void polyfem::assembler::LinearAssembler::assemble ( const bool  is_volume,
const int  n_basis,
const std::vector< basis::ElementBases > &  bases,
const std::vector< basis::ElementBases > &  gbases,
const AssemblyValsCache cache,
const double  t,
StiffnessMatrix stiffness,
const bool  is_mass = false 
) const
overridevirtual

assembles the stiffness matrix for the given basis the bilinear form (local assembler) is encoded by the overloaded assemble (see below) function that the subclass (eg Laplacian) defines sets stiffness and modifies cache if it has not already been computed

Reimplemented from polyfem::assembler::LinearAssembler.

◆ assemble() [2/3]

Eigen::Matrix< double, Eigen::Dynamic, 1, 0, 9, 1 > polyfem::assembler::StokesVelocity::assemble ( const LinearAssemblerData data) const
overridevirtual

local assembly function that defines the bilinear form (LHS) computes and returns a single local stiffness value

Implements polyfem::assembler::LinearAssembler.

Definition at line 17 of file Stokes.cpp.

References polyfem::assembler::ElementAssemblyValues::basis_values, polyfem::assembler::LinearAssemblerData::da, polyfem::assembler::ElementAssemblyValues::element_id, polyfem::assembler::LinearAssemblerData::i, polyfem::assembler::LinearAssemblerData::j, polyfem::assembler::Assembler::size(), polyfem::assembler::LinearAssemblerData::t, polyfem::assembler::ElementAssemblyValues::val, polyfem::assembler::LinearAssemblerData::vals, and viscosity_.

Here is the call graph for this function:

◆ assemble() [3/3]

virtual Eigen::Matrix< double, Eigen::Dynamic, 1, 0, 9, 1 > polyfem::assembler::LinearAssembler::assemble ( const LinearAssemblerData data) const
virtual

local assembly function that defines the bilinear form (LHS) computes and returns a single local stiffness value

Implements polyfem::assembler::LinearAssembler.

◆ compute_rhs()

VectorNd polyfem::assembler::StokesVelocity::compute_rhs ( const AutodiffHessianPt pt) const
overridevirtual

Reimplemented from polyfem::assembler::Assembler.

Definition at line 38 of file Stokes.cpp.

References polyfem::assembler::Assembler::size(), val, and viscosity_.

Here is the call graph for this function:

◆ is_fluid()

bool polyfem::assembler::StokesVelocity::is_fluid ( ) const
inlineoverridevirtual

Reimplemented from polyfem::assembler::Assembler.

Definition at line 30 of file Stokes.hpp.

◆ is_tensor()

bool polyfem::assembler::StokesVelocity::is_tensor ( ) const
inlineoverridevirtual

Reimplemented from polyfem::assembler::Assembler.

Definition at line 31 of file Stokes.hpp.

◆ name()

virtual std::string polyfem::assembler::StokesVelocity::name ( ) const
inlineoverridevirtual

Implements polyfem::assembler::Assembler.

Reimplemented in polyfem::assembler::OperatorSplitting.

Definition at line 27 of file Stokes.hpp.

◆ parameters()

std::map< std::string, Assembler::ParamFunc > polyfem::assembler::StokesVelocity::parameters ( ) const
overridevirtual

Implements polyfem::assembler::Assembler.

Definition at line 56 of file Stokes.cpp.

References viscosity_.

◆ viscosity()

const GenericMatParam & polyfem::assembler::StokesVelocity::viscosity ( ) const
inline

Definition at line 25 of file Stokes.hpp.

References viscosity_.

Member Data Documentation

◆ viscosity_

GenericMatParam polyfem::assembler::StokesVelocity::viscosity_
private

Definition at line 34 of file Stokes.hpp.

Referenced by add_multimaterial(), assemble(), compute_rhs(), parameters(), and viscosity().


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