Skip to content

Commit

Permalink
fixes + format
Browse files Browse the repository at this point in the history
  • Loading branch information
henrij22 committed Dec 2, 2024
1 parent b0421b0 commit f6fbfb9
Show file tree
Hide file tree
Showing 9 changed files with 34 additions and 18 deletions.
2 changes: 1 addition & 1 deletion ikarus/finiteelements/mechanics/kirchhoffloveshell.hh
Original file line number Diff line number Diff line change
Expand Up @@ -289,7 +289,7 @@ protected:
template <typename ScalarType>
void calculateMassImpl(const Requirement& par, const MatrixAffordance& affordance,
typename Traits::template MatrixType<> M) const {
const auto geo = underlying().localView().element().geometry();
const auto geo = underlying().localView().element().geometry();
const auto rhoT = thickness_ * density_;

for (const auto& [gpIndex, gp] : localBasis_.viewOverIntegrationPoints()) {
Expand Down
24 changes: 17 additions & 7 deletions ikarus/solver/eigenvaluesolver/generalizedeigensolver.hh
Original file line number Diff line number Diff line change
Expand Up @@ -24,6 +24,13 @@

namespace Ikarus {

namespace Impl {
void assertNev(std::optional<Eigen::Index> nev, Eigen::Index nevMax) {
if (nev.has_value() and nev.value() > nevMax)
DUNE_THROW(Dune::InvalidStateException, "You can not ask for more eigenvalues or eigenvectors then calculated.");
}
} // namespace Impl

/**
* \brief A strongly typed enum class representing the type of solver to use for the eigenvalue problem
*/
Expand Down Expand Up @@ -101,25 +108,27 @@ struct GeneralizedSymEigenSolver<EigenValueSolverType::Spectra, MT>
* \brief Starts the computation of the eigenvalue solver
*
* \param tolerance given tolerance for iterative eigenvalue solving (default: 1e-10)
* \param maxit givenn maximum iterations for eigenvalue solving (default: 1000)
* \param maxit given maximum iterations for eigenvalue solving (default: 1000)
* \return true solving was successful
* \return false solving was not successful
*/
bool compute(Eigen::Index maxit = 1000, ScalarType tolerance = 1e-10) {
solverSmallest_.init();
solverSmallest_.compute(Spectra::SortRule::SmallestAlge, maxit, tolerance, Spectra::SortRule::SmallestAlge);

eigenvalues_.head(nevsPartition_.first) = solverSmallest_.eigenvalues();
eigenvectors_.leftCols(nevsPartition_.first) = solverSmallest_.eigenvectors();

solverGreatest_.init();
solverGreatest_.compute(Spectra::SortRule::LargestAlge, maxit, tolerance, Spectra::SortRule::SmallestAlge);

eigenvalues_.tail(nevsPartition_.second) = solverGreatest_.eigenvalues();
eigenvectors_.rightCols(nevsPartition_.second) = solverGreatest_.eigenvectors();

computed_ = solverSmallest_.info() == Spectra::CompInfo::Successful and
solverGreatest_.info() == Spectra::CompInfo::Successful;
if (not computed_)
return false;

eigenvalues_.head(nevsPartition_.first) = solverSmallest_.eigenvalues();
eigenvectors_.leftCols(nevsPartition_.first) = solverSmallest_.eigenvectors();

eigenvalues_.tail(nevsPartition_.second) = solverGreatest_.eigenvalues();
eigenvectors_.rightCols(nevsPartition_.second) = solverGreatest_.eigenvectors();

return computed_;
}
Expand Down Expand Up @@ -375,6 +384,7 @@ struct PartialGeneralizedSymEigenSolver
* \return auto matrix with the eigevectors as columns
*/
auto eigenvectors(std::optional<Eigen::Index> _nev = std::nullopt) const {
Impl::assertNev(_nev, nev_);
assertCompute();
return solver_.eigenvectors(_nev.value_or(nev_));
}
Expand Down
1 change: 1 addition & 0 deletions ikarus/utils/concepts.hh
Original file line number Diff line number Diff line change
Expand Up @@ -636,6 +636,7 @@ namespace Concepts {
{ es.compute() } -> std::same_as<bool>;
{ es.eigenvalues() } -> std::convertible_to<Eigen::VectorX<typename ES::ScalarType>>;
{ es.eigenvectors() } -> std::convertible_to<Eigen::MatrixX<typename ES::ScalarType>>;
{ es.nev() } -> std::convertible_to<int>;
};

} // namespace Concepts
Expand Down
3 changes: 2 additions & 1 deletion ikarus/utils/modalanalysis/lumpingschemes.hh
Original file line number Diff line number Diff line change
Expand Up @@ -15,7 +15,8 @@
namespace Ikarus::Dynamics::LumpingSchemes {

/**
* \brief Implements the operator() for a row-sum lumping scheme that can be passed to the assemblermanipulator to lump the already assembled global mass matrix.
* \brief Implements the operator() for a row-sum lumping scheme that can be passed to the assemblermanipulator to lump
* the already assembled global mass matrix.
*/
struct RowSumLumping
{
Expand Down
9 changes: 5 additions & 4 deletions ikarus/utils/modalanalysis/modalanalysis.hh
Original file line number Diff line number Diff line change
Expand Up @@ -98,9 +98,9 @@ struct ModalAnalysis
* \return true solving was successful
* \return false solving was not successful
*/
bool compute(ScalarType tolerance = 1e-10, Eigen::Index maxit = 1000) {
bool compute(Eigen::Index maxit = 1000, ScalarType tolerance = 1e-10) {
solver_.emplace(stiffAssembler_, lumpedMassAssembler_);
return solver_->compute(tolerance, maxit);
return solver_->compute(maxit, tolerance);
}

/**
Expand Down Expand Up @@ -174,9 +174,10 @@ struct ModalAnalysis
* \param filename filename of the output pvd file.
* \param nev_ optionally specify how many eigenmodes should be written out, defaults to all.
*/
void writeEigenModes(const std::string& filename, std::optional<Eigen::Index> nev_ = std::nullopt) const {
void writeEigenModes(const std::string& filename, std::optional<Eigen::Index> _nev = std::nullopt) const {
assertCompute();
writeEigenmodesAsTimeSeries(solver_.value(), stiffAssembler_, filename, nev_);
Impl::assertNev(_nev, nev());
writeEigenmodesAsTimeSeries(solver_.value(), stiffAssembler_, filename, _nev);
}

/** \brief Returns the number of eigenvalues of the problem */
Expand Down
5 changes: 4 additions & 1 deletion ikarus/utils/modalanalysis/modalanalysishelper.hh
Original file line number Diff line number Diff line change
Expand Up @@ -14,6 +14,7 @@

#include <ikarus/assembler/assemblermanipulatorfuser.hh>
#include <ikarus/io/vtkwriter.hh>
#include <ikarus/solver/eigenvaluesolver/generalizedeigensolver.hh>
#include <ikarus/utils/concepts.hh>
#include <ikarus/utils/makeenum.hh>
#include <ikarus/utils/modalanalysis/lumpingschemes.hh>
Expand Down Expand Up @@ -51,6 +52,7 @@ auto makeLumpedFlatAssembler(const std::shared_ptr<AS>& assembler) {
template <Concepts::EigenValueSolver Eigensolver, Concepts::FlatAssembler Assembler>
void writeEigenmodesToVTK(const Eigensolver& solver, std::shared_ptr<Assembler> assembler, const std::string& filename,
std::optional<Eigen::Index> nev_ = std::nullopt) {
Impl::assertNev(nev_, solver.nev());
auto nev = nev_.value_or(solver.nev());
auto eigenvectors = solver.eigenvectors();
auto basis = assembler->basis();
Expand Down Expand Up @@ -107,8 +109,9 @@ void writeEigenmodesAsTimeSeries(const Eigen::MatrixBase<Derived>& eigenvectors,
template <Concepts::EigenValueSolver Eigensolver, Concepts::FlatAssembler Assembler>
void writeEigenmodesAsTimeSeries(const Eigensolver& solver, std::shared_ptr<Assembler> assembler,
const std::string& filename, std::optional<Eigen::Index> nev_ = std::nullopt) {
Impl::assertNev(nev_, solver.nev());
auto nev = nev_.value_or(solver.nev());
auto eigenvectors = solver.eigenvectors();
auto eigenvectors = solver.eigenvectors().leftCols(nev).eval();
writeEigenmodesAsTimeSeries(eigenvectors, assembler, filename);
}

Expand Down
3 changes: 2 additions & 1 deletion tests/src/dummyproblem.hh
Original file line number Diff line number Diff line change
Expand Up @@ -43,7 +43,8 @@ struct DummyProblem

// YASPGrid needs an int, structuresgridfactory an unsigned int
explicit DummyProblem(
const std::array<std::conditional_t<useYASP, int, unsigned int>, 2>& elementsPerDirection = {10, 10}, double Lx = 4.0, double Ly = 4.0)
const std::array<std::conditional_t<useYASP, int, unsigned int>, 2>& elementsPerDirection = {10, 10},
double Lx = 4.0, double Ly = 4.0)
: grid_([&]() {
const Dune::FieldVector<double, 2> bbox = {Lx, Ly};

Expand Down
1 change: 0 additions & 1 deletion tests/src/testeigenvaluesolver.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -121,7 +121,6 @@ auto testEigenSolvers() {
return t;
}


auto checkThrows() {
TestSuite t;
using Ikarus::EigenValueSolverType;
Expand Down
4 changes: 2 additions & 2 deletions tests/src/testenhancedassumedstrains.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -46,12 +46,12 @@ int main(int argc, char** argv) {
t.subTest(testFEElement(firstOrderLagrangePrePower2Basis, "EAS", randomlyDistorted,
Dune::ReferenceElements<double, 2>::cube(), linearElasticFuncPlaneStress,
Ikarus::skills(Ikarus::eas()), Ikarus::AffordanceCollections::elastoStatics,
checkJacobianFunctor));
checkJacobianFunctor, singleElementTestFunctor));

t.subTest(testFEElement(firstOrderLagrangePrePower3Basis, "EAS", randomlyDistorted,
Dune::ReferenceElements<double, 3>::cube(), linearElasticFunc3D,
Ikarus::skills(Ikarus::eas()), Ikarus::AffordanceCollections::elastoStatics,
checkJacobianFunctor));
checkJacobianFunctor, singleElementTestFunctor));

// Plane Stress
t.subTest(testFEElement(
Expand Down

0 comments on commit f6fbfb9

Please sign in to comment.