Skip to content

Use AutoDiffScalar Concept and fix checkFESByAutoDiff test #329

New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Merged
merged 18 commits into from
Oct 21, 2024
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
13 changes: 7 additions & 6 deletions ikarus/finiteelements/mechanics/enhancedassumedstrains.hh
Original file line number Diff line number Diff line change
Expand Up @@ -18,6 +18,7 @@
#include <ikarus/finiteelements/ferequirements.hh>
#include <ikarus/finiteelements/mechanics/easvariants.hh>
#include <ikarus/finiteelements/mixin.hh>
#include <ikarus/utils/concepts.hh>

namespace Ikarus {

Expand Down Expand Up @@ -165,10 +166,10 @@ protected:
return;

decltype(auto) LMat = [this]() -> decltype(auto) {
if constexpr (std::is_same_v<ScalarType, double>)
return [this]() -> Eigen::MatrixXd& { return L_; }();
else
if constexpr (Concepts::AutodiffScalar<ScalarType>)
return Eigen::MatrixX<ScalarType>{};
else
return [this]() -> Eigen::MatrixXd& { return L_; }();
}();

auto calculateMatrixContribution = [&]<typename EAST>(const EAST& easFunction) {
Expand Down Expand Up @@ -208,10 +209,10 @@ protected:
calculateDAndLMatrix(easFunction, par, D, L_);

decltype(auto) LMat = [this]() -> decltype(auto) {
if constexpr (std::is_same_v<ScalarType, double>)
return [this]() -> Eigen::MatrixXd& { return L_; }();
else
if constexpr (Concepts::AutodiffScalar<ScalarType>)
return Eigen::MatrixX<ScalarType>{};
else
return [this]() -> Eigen::MatrixXd& { return L_; }();
}();
const auto disp = Dune::viewAsFlatEigenVector(uFunction.coefficientsRef());
const auto alpha = (-D.inverse() * L_ * disp).eval();
Expand Down
25 changes: 19 additions & 6 deletions ikarus/finiteelements/mechanics/materials/neohooke.hh
Original file line number Diff line number Diff line change
Expand Up @@ -72,11 +72,13 @@ struct NeoHookeT : public Material<NeoHookeT<ST>>
* \return ScalarType The stored energy.
*/
template <typename Derived>
ScalarType storedEnergyImpl(const Eigen::MatrixBase<Derived>& C) const noexcept {
ScalarType storedEnergyImpl(const Eigen::MatrixBase<Derived>& C) const {
static_assert(Concepts::EigenMatrixOrVoigtNotation3<Derived>);
if constexpr (!Concepts::EigenVector<Derived>) {
const auto traceC = C.trace();
const auto logdetF = log(sqrt(C.determinant()));
const auto traceC = C.trace();
const auto detC = C.determinant();
checkPositiveDetC(detC);
const auto logdetF = log(sqrt(detC));
return materialParameter_.mu / 2.0 * (traceC - 3 - 2 * logdetF) +
materialParameter_.lambda / 2.0 * logdetF * logdetF;
} else
Expand All @@ -96,7 +98,9 @@ struct NeoHookeT : public Material<NeoHookeT<ST>>
static_assert(Concepts::EigenMatrixOrVoigtNotation3<Derived>);
if constexpr (!voigt) {
if constexpr (!Concepts::EigenVector<Derived>) {
const auto logdetF = log(sqrt(C.determinant()));
const auto detC = C.determinant();
checkPositiveDetC(detC);
const auto logdetF = log(sqrt(detC));
const auto invC = C.inverse().eval();
return (materialParameter_.mu * (StrainMatrix::Identity() - invC) + materialParameter_.lambda * logdetF * invC)
.eval();
Expand All @@ -118,8 +122,10 @@ struct NeoHookeT : public Material<NeoHookeT<ST>>
auto tangentModuliImpl(const Eigen::MatrixBase<Derived>& C) const {
static_assert(Concepts::EigenMatrixOrVoigtNotation3<Derived>);
if constexpr (!voigt) {
const auto invC = C.inverse().eval();
const auto logdetF = log(sqrt(C.determinant()));
const auto invC = C.inverse().eval();
const auto detC = C.determinant();
checkPositiveDetC(detC);
const auto logdetF = log(sqrt(detC));
const auto CTinv = tensorView(invC, std::array<Eigen::Index, 2>({3, 3}));
static_assert(Eigen::TensorFixedSize<ScalarType, Eigen::Sizes<3, 3>>::NumIndices == 2);
Eigen::TensorFixedSize<ScalarType, Eigen::Sizes<3, 3, 3, 3>> moduli =
Expand All @@ -144,6 +150,13 @@ struct NeoHookeT : public Material<NeoHookeT<ST>>

private:
MaterialParameters materialParameter_;

void checkPositiveDetC(ScalarType detC) const {
if (Dune::FloatCmp::le(static_cast<double>(detC), 0.0, 1e-10))
DUNE_THROW(Dune::InvalidStateException,
"Determinant of right Cauchy Green tensor C must be greater than zero. detC = " +
std::to_string(static_cast<double>(detC)));
}
};

/**
Expand Down
29 changes: 11 additions & 18 deletions ikarus/finiteelements/mechanics/nonlinearelastic.hh
Original file line number Diff line number Diff line change
Expand Up @@ -152,12 +152,7 @@ public:
*/
template <typename ScalarType, int strainDim, bool voigt = true>
auto materialTangent(const Eigen::Vector<ScalarType, strainDim>& strain) const {
if constexpr (std::is_same_v<ScalarType, double>)
return mat_.template tangentModuli<strainType, voigt>(strain);
else {
decltype(auto) matAD = mat_.template rebind<ScalarType>();
return matAD.template tangentModuli<strainType, voigt>(strain);
}
return material<ScalarType>().template tangentModuli<strainType, voigt>(strain);
}

/**
Expand All @@ -170,12 +165,7 @@ public:
*/
template <typename ScalarType, int strainDim>
auto getInternalEnergy(const Eigen::Vector<ScalarType, strainDim>& strain) const {
if constexpr (std::is_same_v<ScalarType, double>)
return mat_.template storedEnergy<strainType>(strain);
else {
decltype(auto) matAD = mat_.template rebind<ScalarType>();
return matAD.template storedEnergy<strainType>(strain);
}
return material<ScalarType>().template storedEnergy<strainType>(strain);
}

/**
Expand All @@ -189,12 +179,7 @@ public:
*/
template <typename ScalarType, int strainDim, bool voigt = true>
auto getStress(const Eigen::Vector<ScalarType, strainDim>& strain) const {
if constexpr (std::is_same_v<ScalarType, double>)
return mat_.template stresses<strainType, voigt>(strain);
else {
decltype(auto) matAD = mat_.template rebind<ScalarType>();
return matAD.template stresses<strainType, voigt>(strain);
}
return material<ScalarType>().template stresses<strainType, voigt>(strain);
}

const Geometry& geometry() const { return *geo_; }
Expand Down Expand Up @@ -237,6 +222,14 @@ private:
size_t numberOfNodes_{0};
int order_{};

template <typename ScalarType>
decltype(auto) material() const {
if constexpr (Concepts::AutodiffScalar<ScalarType>)
return mat_.template rebind<ScalarType>();
else
return mat_;
}

protected:
/**
* \brief Calculate the matrix associated with the given Requirement.
Expand Down
13 changes: 7 additions & 6 deletions tests/src/checkfebyautodiff.hh
Original file line number Diff line number Diff line change
Expand Up @@ -16,7 +16,8 @@

template <typename GridView, typename BasisHandler, typename Skills, typename AffordanceColl, typename VectorType>
auto checkFESByAutoDiffImpl(const GridView& gridView, const BasisHandler& basis, Skills&& skills,
AffordanceColl affordance, VectorType& d, const std::string& messageIfFailed = "") {
AffordanceColl affordance, VectorType& d, const std::string& messageIfFailed = "",
double tol = 1e-10) {
double lambda = 7.3;
auto fe = Ikarus::makeFE(basis, std::forward<Skills>(skills));
using FE = decltype(fe);
Expand All @@ -28,8 +29,7 @@ auto checkFESByAutoDiffImpl(const GridView& gridView, const BasisHandler& basis,
for (auto element : elements(gridView)) {
auto localView = basis.flat().localView();
localView.bind(element);
auto nDOF = localView.size();
const double tol = 1e-10;
auto nDOF = localView.size();

fe.bind(element);

Expand Down Expand Up @@ -80,13 +80,14 @@ auto checkFESByAutoDiffImpl(const GridView& gridView, const BasisHandler& basis,

template <typename GridView, typename PreBasis, typename Skills, typename AffordanceColl>
auto checkFESByAutoDiff(const GridView& gridView, const PreBasis& pb, Skills&& skills, AffordanceColl affordance,
const std::string& testName = "") {
const std::string& testName = "", double tol = 1e-10) {
Dune::TestSuite t("AutoDiff Test" + testName);
auto basis = Ikarus::makeBasis(gridView, pb);
Eigen::VectorXd d;
d.setZero(basis.flat().dimension());
t.subTest(checkFESByAutoDiffImpl(gridView, basis, skills, affordance, d, " Zero Displacements"));
t.subTest(checkFESByAutoDiffImpl(gridView, basis, skills, affordance, d, " Zero Displacements", tol));
d.setRandom(basis.flat().dimension());
t.subTest(checkFESByAutoDiffImpl(gridView, basis, skills, affordance, d, " Non-zero Displacements"));
d *= 0.2; // to avoid a negative determinant of deformation gradient
t.subTest(checkFESByAutoDiffImpl(gridView, basis, skills, affordance, d, " Non-zero Displacements", tol));
return t;
}
24 changes: 24 additions & 0 deletions tests/src/testmaterial.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -228,6 +228,29 @@ auto testPlaneStrainAgainstPlaneStress(const double tol = 1e-10) {
return t;
}

template <typename MAT>
auto checkThrowNeoHooke(const MAT& matNH) {
static_assert(std::is_same_v<MAT, NeoHookeT<double>>,
"checkThrowNeoHooke is only implemented for Neo-Hooke material law.");
TestSuite t("NeoHooke Test - Checks the throw message for negative determinant of C");
Eigen::Vector3d E;
E << 2.045327969583023, 0.05875570522766141, 0.3423966429644326;
auto reducedMat = planeStress(matNH, 1e-8);

t.checkThrow<Dune::InvalidStateException>(
[&]() { const auto moduli = (reducedMat.template tangentModuli<StrainTags::greenLagrangian, true>(E)); },
"Neo-Hooke test (tangentModuli) should have failed with negative detC for the given E");

t.checkThrow<Dune::InvalidStateException>(
[&]() { const auto stress = (reducedMat.template stresses<StrainTags::greenLagrangian, true>(E)); },
"Neo-Hooke test (stresses) should have failed with negative detC for the given E");

t.checkThrow<Dune::InvalidStateException>(
[&]() { const auto energy = (reducedMat.template storedEnergy<StrainTags::greenLagrangian>(E)); },
"Neo-Hooke test (stresses) should have failed with negative detC for the given E");
return t;
}

int main(int argc, char** argv) {
Ikarus::init(argc, argv);
TestSuite t;
Expand All @@ -239,6 +262,7 @@ int main(int argc, char** argv) {

auto nh = NeoHooke(matPar);
t.subTest(testMaterial(nh));
t.subTest(checkThrowNeoHooke(nh));

auto le = LinearElasticity(matPar);
t.subTest(testMaterial(le));
Expand Down
4 changes: 2 additions & 2 deletions tests/src/testnonlinearelasticity.hh
Original file line number Diff line number Diff line change
Expand Up @@ -290,7 +290,7 @@ auto SingleElementTest(const Material& mat) {
}

template <int gridDim, typename MAT, typename TestSuiteType>
void autoDiffTest(TestSuiteType& t, const MAT& mat, const std::string& testName = "") {
void autoDiffTest(TestSuiteType& t, const MAT& mat, const std::string& testName = "", double tol = 1e-10) {
using namespace Ikarus;
using namespace Dune::Functions::BasisFactory;
auto vL = []<typename VectorType>([[maybe_unused]] const VectorType& globalCoord, auto& lamb) {
Expand All @@ -317,5 +317,5 @@ void autoDiffTest(TestSuiteType& t, const MAT& mat, const std::string& testName
t.subTest(checkFESByAutoDiff(
gridView, power<gridDim>(lagrange<1>()),
skills(Ikarus::nonLinearElastic(mat), volumeLoad<gridDim>(vL), neumannBoundaryLoad(&neumannBoundary, nBL)),
Ikarus::AffordanceCollections::elastoStatics, testName));
Ikarus::AffordanceCollections::elastoStatics, testName, tol));
}
4 changes: 2 additions & 2 deletions tests/src/testnonlinearelasticityneohooke.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -28,8 +28,8 @@ int main(int argc, char** argv) {
t.subTest(NonLinearElasticityLoadControlNRandTR<Grids::Yasp>(matNH1));
t.subTest(NonLinearElasticityLoadControlNRandTR<Grids::IgaSurfaceIn2D>(matNH1));

autoDiffTest<2>(t, planeStressMat1, " nu != 0");
autoDiffTest<2>(t, planeStressMat2, " nu = 0");
autoDiffTest<2>(t, planeStressMat1, " nu != 0", 1e-7);
Copy link
Collaborator

@tarun-mitruka tarun-mitruka Oct 16, 2024

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

This tolerance isn't always necessary but only for certain cases depending on the randomness in the element geometry and the displacement vector (see comment). I am not sure if this is good enough or it could be resolved in a better way? May be fixing the random displacement vector in checkFESByAutoDiff such that this tolerance can be reduced further? We can also fix the tolerance while working on the Issue #211. Do you think increasing minIter to 2 or 3 for NewtonRaphson in VanishingStress would help here?

Copy link
Collaborator

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Increasing minIter might help but I'm not sure here. But using tolerances here is not to bad. But I think there is no nice way to fix this. We might refactor isApproxSame which then uses the objects to compute something reasonable: Maybe something like that https://godbolt.org/z/xEPraWKEf

Copy link
Collaborator

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Okay then I leave the tolerance here as such and we can refactor isApproxSame while perhaps working on #211

autoDiffTest<2>(t, planeStressMat2, " nu = 0", 1e-7);
autoDiffTest<3>(t, matNH1, " nu != 0");
autoDiffTest<3>(t, matNH2, " nu = 0");
return t.exit();
Expand Down
Loading