version 0.4.1
linearalgebrahelper.hh
Go to the documentation of this file.
1// SPDX-FileCopyrightText: 2021-2024 The Ikarus Developers mueller@ibb.uni-stuttgart.de
2// SPDX-License-Identifier: LGPL-3.0-or-later
3
9#pragma once
10#include "concepts.hh"
11#include "traits.hh"
12
13#include <iosfwd>
14#include <random>
15#include <ranges>
16
17#include <dune/common/tuplevector.hh>
18#include <dune/functions/functionspacebases/lagrangebasis.hh>
19#include <dune/istl/bvector.hh>
20
21#include <Eigen/Core>
22
23#include <autodiff/forward/dual/dual.hpp>
24
25namespace Ikarus {
26
34template <typename Derived>
35auto orthonormalizeMatrixColumns(const Eigen::MatrixBase<Derived>& A) {
36 // Gram Schmidt Ortho
37 auto Q = A.eval();
38
39 Q.col(0).normalize();
40
41 for (int colIndex = 1; colIndex < Q.cols(); colIndex++) {
42 Q.col(colIndex) -= Q.leftCols(colIndex) * (Q.leftCols(colIndex).transpose() * A.col(colIndex));
43 Q.col(colIndex).normalize();
44 }
45
46 return Q;
47}
48
56template <typename ValueType>
57auto viewAsFlatEigenVector(Dune::BlockVector<ValueType>& blockedVector) {
58 Eigen::Map<Eigen::VectorX<typename ValueType::field_type>> vec(&blockedVector.begin()->begin().operator*(),
59 blockedVector.size() * blockedVector[0].size());
60
61 return vec;
62}
63
71template <typename ValueType>
72auto viewAsFlatEigenVector(const Dune::BlockVector<ValueType>& blockedVector) {
73 Eigen::Map<const Eigen::VectorX<typename ValueType::field_type>> vec(&blockedVector.begin()->begin().operator*(),
74 blockedVector.size() * blockedVector[0].size());
75
76 return vec;
77}
78
87template <typename ValueType>
88auto viewAsEigenMatrixAsDynFixed(Dune::BlockVector<ValueType>& blockedVector) {
89 Eigen::Map<Eigen::Matrix<typename ValueType::field_type, Eigen::Dynamic, ValueType::valueSize, Eigen::RowMajor>> vec(
90 &blockedVector.begin()->begin().operator*(), blockedVector.size(), blockedVector[0].size());
91
92 return vec;
93}
94
103template <typename ValueType>
104auto viewAsEigenMatrixAsDynFixed(const Dune::BlockVector<ValueType>& blockedVector) {
105 Eigen::Map<const Eigen::Matrix<typename ValueType::field_type, Eigen::Dynamic, ValueType::valueSize, Eigen::RowMajor>>
106 vec(&blockedVector.begin()->begin().operator*(), blockedVector.size(), blockedVector[0].size());
107
108 return vec;
109}
110
119template <typename ValueType>
120auto viewAsEigenMatrixFixedDyn(Dune::BlockVector<ValueType>& blockedVector) {
121 Eigen::Map<Eigen::Matrix<typename ValueType::field_type, ValueType::valueSize, Eigen::Dynamic>> vec(
122 &blockedVector.begin()->begin().operator*(), blockedVector[0].size(), blockedVector.size());
123
124 return vec;
125}
126
134template <typename ValueType>
135auto viewAsEigenMatrixFixedDyn(const Dune::BlockVector<ValueType>& blockedVector) {
136 Eigen::Map<const Eigen::Matrix<typename ValueType::field_type, ValueType::valueSize, Eigen::Dynamic>> vec(
137 &blockedVector.begin()->begin().operator*(), blockedVector[0].size(), blockedVector.size());
138
139 return vec;
140}
141
149template <typename Type>
150size_t correctionSize(const Dune::BlockVector<Type>& a)
151requires requires { Type::correctionSize; }
152{
153 return a.size() * Type::correctionSize;
154}
155
163template <typename Type>
164size_t valueSize(const Dune::BlockVector<Type>& a)
165requires requires { Type::valueSize; }
166{
167 return a.size() * Type::valueSize;
168}
169
179template <typename Type, typename Derived>
180Dune::BlockVector<Type>& operator+=(Dune::BlockVector<Type>& a, const Eigen::MatrixBase<Derived>& b)
181requires(Ikarus::Concepts::AddAssignAble<Type, decltype(b.template segment<Type::correctionSize>(0))> and
182 requires() { Type::correctionSize; })
183{
184 assert(correctionSize(a) == static_cast<size_t>(b.size()) && " The passed vector has wrong size");
185 for (auto i = 0U; i < a.size(); ++i)
186 a[i] += b.template segment<Type::correctionSize>(i * Type::correctionSize);
187 return a;
188}
189
199template <typename Type, typename Derived>
200Dune::BlockVector<Type>& operator-=(Dune::BlockVector<Type>& a, const Eigen::MatrixBase<Derived>& b)
201requires(Ikarus::Concepts::AddAssignAble<Type, decltype(b.template segment<Type::correctionSize>(0))> and
202 requires() { Type::correctionSize; })
203{
204 return a += (-b);
205}
206
216template <typename... Types, typename Derived>
217Dune::TupleVector<Types...>& operator+=(Dune::TupleVector<Types...>& a, const Eigen::MatrixBase<Derived>& b) {
218 using namespace Dune::Indices;
219 size_t posStart = 0;
220 Dune::Hybrid::forEach(Dune::Hybrid::integralRange(Dune::index_constant<sizeof...(Types)>()), [&](const auto i) {
221 const size_t size = correctionSize(a[i]);
222 a[i] += b(Eigen::seqN(posStart, size));
223 posStart += size;
224 });
225
226 return a;
227}
228
238template <typename ManifoldPoint, typename Derived>
239Dune::BlockVector<ManifoldPoint>& addInEmbedding(Dune::BlockVector<ManifoldPoint>& a,
240 const Eigen::MatrixBase<Derived>& b)
241requires(Ikarus::Concepts::AddAssignAble<ManifoldPoint, decltype(b.template segment<ManifoldPoint::valueSize>(0))> and
242 requires() { ManifoldPoint::valueSize; })
243{
244 assert(valueSize(a) == static_cast<size_t>(b.size()) && " The passed vector has wrong size");
245 for (auto i = 0U; i < a.size(); ++i)
246 a[i].addInEmbedding(b.template segment<ManifoldPoint::valueSize>(i * ManifoldPoint::valueSize));
247 return a;
248}
249
257template <typename Derived>
258requires(!std::floating_point<Derived>)
259auto norm(const Eigen::MatrixBase<Derived>& v) {
260 return v.norm();
261}
262
269auto norm(const std::floating_point auto& v) { return std::abs(v); }
270
280template <typename Scalar, int size>
281auto operator*(const Eigen::DiagonalMatrix<Scalar, size>& a, const Eigen::DiagonalMatrix<Scalar, size>& b) {
282 return (a.diagonal().cwiseProduct(b.diagonal())).asDiagonal();
283}
284
294template <typename Scalar, int size>
295auto operator+=(Eigen::DiagonalMatrix<Scalar, size>& a, const Eigen::DiagonalMatrix<Scalar, size>& b) {
296 a.diagonal() += b.diagonal();
297 return a;
298}
299
310template <typename Derived, typename Scalar, int size>
311auto operator+(const Eigen::MatrixBase<Derived>& a, const Eigen::DiagonalMatrix<Scalar, size>& b) {
312 auto c = a.derived().eval();
313 c.diagonal() += b.diagonal();
314 return c;
315}
316
327template <typename Derived, typename Scalar, int size>
328auto operator+(const Eigen::DiagonalMatrix<Scalar, size>& a, const Eigen::MatrixBase<Derived>& b) {
329 return b + a;
330}
331
340template <typename Scalar, int size>
341auto operator-(const Eigen::DiagonalMatrix<Scalar, size>& a) {
342 return (-a.diagonal()).asDiagonal();
343}
344
354template <typename Derived, typename Derived2>
355auto operator+(const Eigen::MatrixBase<Derived>& a, const Eigen::DiagonalWrapper<Derived2>& b) {
356 auto c = a.derived().eval();
357 c.diagonal() += b.diagonal();
358 return c;
359}
360
370template <typename Derived, typename Derived2>
371auto operator+(const Eigen::DiagonalWrapper<Derived>& a, const Eigen::MatrixBase<Derived2>& b) {
372 return b + a;
373}
374
384template <typename Scalar, int size>
385std::ostream& operator<<(std::ostream& os, const Eigen::DiagonalMatrix<Scalar, size>& a) {
386 os << Eigen::Matrix<Scalar, size, size>(a);
387 return os;
388}
389
397template <typename Derived>
398Derived sym(const Eigen::MatrixBase<Derived>& A) {
399 return 0.5 * (A + A.transpose());
400}
401
409template <typename Derived>
410Derived skew(const Eigen::MatrixBase<Derived>& A) {
411 return 0.5 * (A - A.transpose());
412}
413
420template <typename Derived>
422 Eigen::IOFormat mapleFmt(Eigen::FullPrecision, 0, ", ", "|\n", "<", ">", "<", ">");
423 if constexpr (std::convertible_to<Derived, const Eigen::MatrixBase<Derived>&>) {
424 std::cout << "\n" << A.derived().format(mapleFmt) << std::endl;
425 } else { // branch for Eigen::DiagonalMatrix
426 using Scalar = typename Derived::Scalar;
427 using namespace Eigen;
428 constexpr int diag_size = EIGEN_SIZE_MIN_PREFER_DYNAMIC(Derived::RowsAtCompileTime, Derived::ColsAtCompileTime);
429 std::cout << "\n"
430 << Eigen::Matrix<Scalar, diag_size, diag_size>(A.derived().diagonal().asDiagonal()).format(mapleFmt)
431 << std::endl;
432 }
433}
434
443template <typename FieldVectorT>
444auto createRandomVector(typename FieldVectorT::value_type lower = -1, typename FieldVectorT::value_type upper = 1) {
445 std::random_device rd;
446 std::mt19937 mt(rd());
447 std::uniform_real_distribution<typename FieldVectorT::value_type> dist(lower, upper);
448 auto rand = [&dist, &mt]() { return dist(mt); };
449 FieldVectorT vec;
450 std::generate(vec.begin(), vec.end(), rand);
451 return vec;
452}
453
461template <typename ScalarType>
462Eigen::Matrix<ScalarType, 3, 3> skew(const Eigen::Vector<ScalarType, 3>& a) {
463 Eigen::Matrix<ScalarType, 3, 3> A;
464 A << 0, -a(2), a(1), a(2), 0, -a(0), -a(1), a(0), 0;
465 return A;
466}
467
468namespace Impl {
469 constexpr std::tuple<std::array<std::array<int, 2>, 1>, std::array<std::array<int, 2>, 3>,
470 std::array<std::array<int, 2>, 6>>
471 voigtIndices = {{{{0, 0}}}, {{{0, 0}, {1, 1}, {0, 1}}}, {{{0, 0}, {1, 1}, {2, 2}, {1, 2}, {0, 2}, {0, 1}}}};
472}
473
482template <int dim>
483constexpr auto voigtNotationContainer = std::get<dim - 1>(Impl::voigtIndices);
484
497template <typename Derived, size_t sizeOfCondensedIndices>
498auto staticCondensation(const Eigen::MatrixBase<Derived>& E,
499 const std::array<size_t, sizeOfCondensedIndices>& indices) {
500 constexpr size_t colsFull = std::remove_cvref_t<Derived>::ColsAtCompileTime;
501 constexpr size_t rowsFull = std::remove_cvref_t<Derived>::RowsAtCompileTime;
502 static_assert(colsFull == rowsFull, "staticCondensation only allowed for square matrices");
503 std::array<size_t, rowsFull - sizeOfCondensedIndices> remainingIndices{};
504 std::ranges::set_difference(std::ranges::iota_view(size_t(0), size_t(colsFull)), indices, remainingIndices.begin());
505
506 auto K11 = E(remainingIndices, remainingIndices);
507 auto K12 = E(indices, remainingIndices);
508 auto K22 = E(indices, indices);
509
510 return (K11 - K12.transpose() * K22.inverse() * K12).eval();
511}
512
525template <typename Derived, size_t sizeOfRemovedCols>
526auto removeCol(const Eigen::MatrixBase<Derived>& E, const std::array<size_t, sizeOfRemovedCols>& indices) {
527 constexpr size_t colsFull = std::remove_cvref_t<Derived>::ColsAtCompileTime;
528 constexpr size_t rowsFull = std::remove_cvref_t<Derived>::RowsAtCompileTime;
529 static_assert(colsFull == 1);
530
531 std::array<size_t, rowsFull - sizeOfRemovedCols> remainingIndices{};
532 std::ranges::set_difference(std::ranges::iota_view(size_t(0), size_t(rowsFull)), indices, remainingIndices.begin());
533
534 return (E(remainingIndices)).eval();
535}
536
551template <typename ST, typename MaterialImpl>
552auto toVoigtAndMaybeReduce(const Eigen::Matrix<ST, 3, 3>& E, [[maybe_unused]] const MaterialImpl& material,
553 bool isStrain = true) {
554 if constexpr (!MaterialImpl::isReduced)
555 return toVoigt(E, isStrain);
556 else {
557 auto ev = toVoigt(E, isStrain);
558 static_assert(decltype(ev)::RowsAtCompileTime == 6 and decltype(ev)::ColsAtCompileTime == 1);
559
560 auto evRed = removeCol(ev, MaterialImpl::fixedVoigtIndices);
561 static_assert(decltype(evRed)::RowsAtCompileTime == 6 - MaterialImpl::fixedVoigtIndices.size() and
562 decltype(evRed)::ColsAtCompileTime == 1);
563 return evRed;
564 }
565}
566
579template <typename M, typename Derived>
580decltype(auto) enlargeIfReduced(const Eigen::MatrixBase<Derived>& E) {
581 if constexpr (!Concepts::EigenVector6<Derived> and Concepts::EigenVector<Derived>) {
582 static_assert(M::isReduced, "You should only end up here, if your material is reduced");
583
584 auto freeindices = M::MaterialImpl::freeVoigtIndices;
585 auto p_E = Eigen::Vector<typename M::MaterialImpl::ScalarType, 6>::Zero().eval();
586 for (int ri = 0; auto i : freeindices)
587 p_E(i) = E(ri++);
588 return p_E;
589
590 } else if constexpr (Concepts::EigenMatrix66<Derived> or Concepts::EigenMatrix33<Derived> or
591 Concepts::EigenVector6<Derived>) {
592 return E.derived();
593 } else {
594 static_assert(M::isReduced, "You should only end up here, if your material is reduced");
595
596 auto freeindices = M::MaterialImpl::freeVoigtIndices;
597 auto p_E = Eigen::Matrix<typename M::MaterialImpl::ScalarType, 6, 6>::Zero().eval();
598 for (int ri = 0; auto i : freeindices) {
599 for (int rj = 0; auto j : freeindices)
600 p_E(i, j) = E(ri, rj++);
601 ++ri;
602 }
603 return p_E;
604 }
605}
606
607} // namespace Ikarus
Several concepts.
Contains stl-like type traits.
auto operator+(const Eigen::MatrixBase< Derived > &a, const Eigen::DiagonalMatrix< Scalar, size > &b)
Eigen::Matrix + Eigen::DiagonalMatrix addition missing in Eigen.
Definition: linearalgebrahelper.hh:311
Dune::BlockVector< Type > & operator+=(Dune::BlockVector< Type > &a, const Eigen::MatrixBase< Derived > &b)
Enables the += operator for Dune::BlockVector += Eigen::Vector.
Definition: linearalgebrahelper.hh:180
auto norm(const Eigen::MatrixBase< Derived > &v)
Adding free norm function to Eigen types.
Definition: linearalgebrahelper.hh:259
Derived skew(const Eigen::MatrixBase< Derived > &A)
Returns the skew part of a matrix.
Definition: linearalgebrahelper.hh:410
auto staticCondensation(const Eigen::MatrixBase< Derived > &E, const std::array< size_t, sizeOfCondensedIndices > &indices)
Performs static condensation on a square matrix.
Definition: linearalgebrahelper.hh:498
auto removeCol(const Eigen::MatrixBase< Derived > &E, const std::array< size_t, sizeOfRemovedCols > &indices)
Removes specified columns from a matrix.
Definition: linearalgebrahelper.hh:526
decltype(auto) enlargeIfReduced(const Eigen::MatrixBase< Derived > &E)
Enlarges a matrix if it reduced in the context of material laws, i.e., VanishingStress If the materia...
Definition: linearalgebrahelper.hh:580
auto toVoigtAndMaybeReduce(const Eigen::Matrix< ST, 3, 3 > &E, const MaterialImpl &material, bool isStrain=true)
Converts a 3x3 matrix to Voigt notation, possibly reducing it based on material properties.
Definition: linearalgebrahelper.hh:552
void printForMaple(const Eigen::EigenBase< Derived > &A)
Method to print the matrix in a format that can directly be copied to Maple.
Definition: linearalgebrahelper.hh:421
auto viewAsFlatEigenVector(Dune::BlockVector< ValueType > &blockedVector)
View Dune::BlockVector as an Eigen::Vector.
Definition: linearalgebrahelper.hh:57
size_t correctionSize(const Dune::BlockVector< Type > &a)
Returns the total correction size of a block vector with a Manifold as the underlying type.
Definition: linearalgebrahelper.hh:150
auto createRandomVector(typename FieldVectorT::value_type lower=-1, typename FieldVectorT::value_type upper=1)
Creates a random vector of the specified type within a given range.
Definition: linearalgebrahelper.hh:444
constexpr auto voigtNotationContainer
Container for Voigt notation indices based on dimension.1D: 0,0 2D: 0,0; 1,1; 0,1 3D: 0,...
Definition: linearalgebrahelper.hh:483
Dune::BlockVector< Type > & operator-=(Dune::BlockVector< Type > &a, const Eigen::MatrixBase< Derived > &b)
Enables the -= operator for Dune::BlockVector += Eigen::Vector.
Definition: linearalgebrahelper.hh:200
auto operator-(const Eigen::DiagonalMatrix< Scalar, size > &a)
Unary minus for Eigen::DiagonalMatrix.
Definition: linearalgebrahelper.hh:341
auto orthonormalizeMatrixColumns(const Eigen::MatrixBase< Derived > &A)
Orthonormalizes all Matrix columns using Gram-Schmidt Orthogonalization.
Definition: linearalgebrahelper.hh:35
Derived sym(const Eigen::MatrixBase< Derived > &A)
Returns the symmetric part of a matrix.
Definition: linearalgebrahelper.hh:398
size_t valueSize(const Dune::BlockVector< Type > &a)
Returns the total value size of a block vector with a Manifold as the underlying type.
Definition: linearalgebrahelper.hh:164
auto viewAsEigenMatrixFixedDyn(Dune::BlockVector< ValueType > &blockedVector)
View Dune::BlockVector as an Eigen::Matrix with fixed rows depending on the size of the ValueType and...
Definition: linearalgebrahelper.hh:120
std::ostream & operator<<(std::ostream &os, const Eigen::DiagonalMatrix< Scalar, size > &a)
Output stream operator for Eigen::DiagonalMatrix.
Definition: linearalgebrahelper.hh:385
auto operator*(const Eigen::DiagonalMatrix< Scalar, size > &a, const Eigen::DiagonalMatrix< Scalar, size > &b)
Eigen::DiagonalMatrix Product Missing in Eigen.
Definition: linearalgebrahelper.hh:281
Dune::BlockVector< ManifoldPoint > & addInEmbedding(Dune::BlockVector< ManifoldPoint > &a, const Eigen::MatrixBase< Derived > &b)
Enables the addition in the embedding space of a vector in the space M^n, where M is a manifold with ...
Definition: linearalgebrahelper.hh:239
auto viewAsEigenMatrixAsDynFixed(Dune::BlockVector< ValueType > &blockedVector)
View Dune::BlockVector as an Eigen::Matrix with dynamic rows and fixed columns depending on the size ...
Definition: linearalgebrahelper.hh:88
constexpr Eigen::Index toVoigt(Eigen::Index i, Eigen::Index j) noexcept
Converts 2D indices to Voigt notation index.
Definition: tensorutils.hh:166
Definition: simpleassemblers.hh:22
Definition: truncatedconjugategradient.hh:24
Definition: concepts.hh:25
Concept defining the requirements for types that support in-place addition.
Definition: concepts.hh:294
Concept defining the requirements for Eigen vectors.
Definition: concepts.hh:353