PolyFEM
Loading...
Searching...
No Matches
SaintVenantElasticity.hpp
Go to the documentation of this file.
1#pragma once
2
7
8namespace polyfem::assembler
9{
10 // Similar to HookeLinear but with non-linear stress strain: C:½(F+Fᵀ+FᵀF)
12 {
13 public:
15
19
20 double compute_energy(const NonLinearAssemblerData &data) const override;
21 Eigen::VectorXd assemble_gradient(const NonLinearAssemblerData &data) const override;
22 Eigen::MatrixXd assemble_hessian(const NonLinearAssemblerData &data) const override;
23
24 VectorNd compute_rhs(const AutodiffHessianPt &pt) const override;
25
26 void set_size(const int size) override;
27
28 void set_stiffness_tensor(int i, int j, const double val);
29 double stifness_tensor(int i, int j) const;
30
31 void add_multimaterial(const int index, const json &params, const Units &units) override;
32
33 std::string name() const override { return "SaintVenant"; }
34 std::map<std::string, ParamFunc> parameters() const override;
35
36 void assign_stress_tensor(const OutputData &data,
37 const int all_size,
38 const ElasticityTensorType &type,
39 Eigen::MatrixXd &all,
40 const std::function<Eigen::MatrixXd(const Eigen::MatrixXd &)> &fun) const override;
41
42 private:
44
45 template <typename T, unsigned long N>
46 T stress(const std::array<T, N> &strain, const int j) const;
47
48 template <typename T>
49 T compute_energy_aux(const NonLinearAssemblerData &data) const;
50 };
51} // namespace polyfem::assembler
double val
Definition Assembler.cpp:86
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 override
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 override
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 override
void add_multimaterial(const int index, const json &params, const Units &units) override
T stress(const std::array< T, N > &strain, const int j) const
void set_stiffness_tensor(int i, int j, const double val)
double compute_energy(const NonLinearAssemblerData &data) const override
T compute_energy_aux(const NonLinearAssemblerData &data) const
double stifness_tensor(int i, int j) const
Eigen::MatrixXd assemble_hessian(const NonLinearAssemblerData &data) const override
std::map< std::string, ParamFunc > parameters() const override
VectorNd compute_rhs(const AutodiffHessianPt &pt) const override
void assign_stress_tensor(const OutputData &data, const int all_size, const ElasticityTensorType &type, Eigen::MatrixXd &all, const std::function< Eigen::MatrixXd(const Eigen::MatrixXd &)> &fun) const override
Eigen::VectorXd assemble_gradient(const NonLinearAssemblerData &data) const override
Used for test only.
Eigen::Matrix< AutodiffScalarHessian, Eigen::Dynamic, 1, 0, 3, 1 > AutodiffHessianPt
Eigen::Matrix< double, Eigen::Dynamic, 1, 0, 3, 1 > VectorNd
Definition Types.hpp:11
nlohmann::json json
Definition Common.hpp:9