MarmotMath.h
Go to the documentation of this file.
1 /* ---------------------------------------------------------------------
2  * _
3  * _ __ ___ __ _ _ __ _ __ ___ ___ | |_
4  * | '_ ` _ \ / _` | '__| '_ ` _ \ / _ \| __|
5  * | | | | | | (_| | | | | | | | | (_) | |_
6  * |_| |_| |_|\__,_|_| |_| |_| |_|\___/ \__|
7  *
8  * Unit of Strength of Materials and Structural Analysis
9  * University of Innsbruck,
10  * 2020 - today
11  *
12  * festigkeitslehre@uibk.ac.at
13  *
14  * Matthias Neuner matthias.neuner@uibk.ac.at
15  * Magdalena Schreter magdalena.schreter@uibk.ac.at
16  * Alexander Dummer alexander.dummer@uibk.ac.at
17  *
18  * This file is part of the MAteRialMOdellingToolbox (marmot).
19  *
20  * This library is free software; you can redistribute it and/or
21  * modify it under the terms of the GNU Lesser General Public
22  * License as published by the Free Software Foundation; either
23  * version 2.1 of the License, or (at your option) any later version.
24  *
25  * The full text of the license can be found in the file LICENSE.md at
26  * the top level directory of marmot.
27  * ---------------------------------------------------------------------
28  */
29 
30 #pragma once
31 #include "Marmot/MarmotConstants.h"
32 #include "Marmot/MarmotTypedefs.h"
33 #include "autodiff/forward/dual.hpp"
34 #include "autodiff/forward/real.hpp"
35 #include <algorithm>
36 #include <autodiff/forward/dual/dual.hpp>
37 #include <complex>
38 
39 namespace Marmot {
40  namespace Math {
43  template < typename T >
44  bool isNaN( T x )
45  {
46  return x != x;
47  }
48 
51  double linearInterpolation( double x, double x0, double x1, double y0, double y1 );
52 
55  double exp( double x );
56 
59  int getExponentPowerTen( const double x );
60 
63  constexpr double radToDeg( const double alpha )
64  {
65  return alpha * 180 / Marmot::Constants::Pi;
66  }
67 
70  constexpr double degToRad( const double alpha )
71  {
72  return alpha / 180 * Marmot::Constants::Pi;
73  }
74 
77  constexpr double macauly( double scalar )
78  {
79  return scalar >= 0 ? scalar : 0.0;
80  }
81 
84  constexpr int heaviside( double scalar )
85  {
86  return scalar >= 0 ? 1 : 0;
87  }
88 
89  constexpr int heavisideExclude0( double scalar )
90  {
91  return scalar > 0 ? 1 : 0;
92  }
93 
96  template < typename T >
97  constexpr int sgn( const T& val ) noexcept
98  {
99  return ( T( 0 ) < val ) - ( val < T( 0 ) );
100  }
101 
102  double makeReal( const double& value );
103  double makeReal( const std::complex< double >& value );
104  double makeReal( const autodiff::real& value );
105 
106  template < typename T, typename G >
107  double makeReal( const autodiff::detail::Dual< T, G >& number )
108  {
109  return double( number );
110  }
111 
112  template < typename T, int... Rest >
113  Eigen::Matrix< double, Rest... > makeReal( const Eigen::Matrix< T, Rest... > mat )
114  {
115  Eigen::Matrix< double, Rest... > out;
116  const size_t m = static_cast< size_t >( mat.rows() );
117  const size_t n = static_cast< size_t >( mat.cols() );
118 
119  for ( size_t i = 0; i < m; i++ ) {
120  for ( size_t j = 0; j < n; j++ ) {
121  out( i, j ) = makeReal( mat( i, j ) );
122  }
123  }
124  return out;
125  }
126  template < typename T >
127  Eigen::VectorXd makeReal( Eigen::Vector< T, Eigen::Dynamic > in )
128  {
129 
130  int inSize = in.size();
131  Eigen::VectorXd out( inSize );
132  for ( int i = 0; i < inSize; i++ ) {
133  out( i ) = double( in( i ) );
134  }
135  return out;
136  }
137 
141  template < int nRows, int nCols >
142  Eigen::Matrix< double, nRows, nCols > macaulyMatrix( const Eigen::Matrix< double, nRows, nCols >& mat )
143  {
144  Eigen::Matrix< double, nRows, nCols > positivePart = Eigen::Matrix< double, nRows, nCols >::Zero();
145  for ( int i = 0; i < nRows; i++ ) {
146  for ( int j = 0; j < nCols; j++ ) {
147  positivePart( i, j ) = macauly( mat( i, j ) );
148  }
149  }
150  return positivePart;
151  }
152 
155  template < typename functionType, typename yType, typename... Args >
156  yType explicitEuler( yType yN, const double dt, functionType fRate, Args&&... fRateArgs )
157  {
158  return yN + fRate( yN, fRateArgs... ) * dt;
159  }
160 
166  template < int ySize, typename functionType, typename... Args >
167  Eigen::Matrix< double, ySize, 1 > semiImplicitEuler( Eigen::Matrix< double, ySize, 1 > yN,
168  const double dt,
169  functionType fRate,
170  Args&&... fRateArgs )
171  {
172  Eigen::Matrix< double, ySize, ySize > fS = Eigen::Matrix< double, ySize, ySize >::Zero();
173  Eigen::Matrix< double, ySize, ySize > Iy = Eigen::Matrix< double, ySize, ySize >::Identity();
174 
175  Eigen::VectorXd leftX( ySize );
176  Eigen::VectorXd rightX( ySize );
177 
178  for ( size_t i = 0; i < ySize; i++ ) {
179  double volatile h = std::max( 1.0, std::abs( yN( i ) ) ) * Constants::cubicRootEps();
180  leftX = yN;
181  leftX( i ) -= h;
182  rightX = yN;
183  rightX( i ) += h;
184  fS.col( i ) = 1. / ( 2. * h ) * ( fRate( rightX, fRateArgs... ) - fRate( leftX, fRateArgs... ) );
185  }
186 
187  return yN + ( Iy - dt * fS ).inverse() * dt * fRate( yN, fRateArgs... );
188  }
189 
193  template < int nRows, int nCols, typename functionType >
194  Eigen::Matrix< double, nRows, nCols > centralDiff( functionType f, const Eigen::Matrix< double, nCols, 1 >& X )
195  {
196  Eigen::Matrix< double, nRows, nCols > dXdY = Eigen::Matrix< double, nRows, nCols >::Zero();
197 
198  Eigen::VectorXd leftX( nCols );
199  Eigen::VectorXd rightX( nCols );
200 
201  for ( size_t i = 0; i < nCols; i++ ) {
202  double volatile h = std::max( 1.0, std::abs( X( i ) ) ) * Constants::cubicRootEps();
203  leftX = X;
204  leftX( i ) -= h;
205  rightX = X;
206  rightX( i ) += h;
207  dXdY.col( i ) = 1. / ( 2. * h ) * ( f( rightX ) - f( leftX ) );
208  }
209 
210  return dXdY;
211  }
212 
217  template < typename functionType, typename yType, typename... Args >
218  yType explicitEulerRichardson( yType yN, const double dt, functionType fRate, Args&&... fRateArgs )
219  {
220  yType fN = fRate( yN, fRateArgs... );
221  yType u = yN + fN * dt;
222  yType v = yN + fN * dt / 2.;
223  yType w = v + fRate( v, fRateArgs... ) * dt / 2.;
224 
225  return 2. * w - u;
226  }
227 
232  template < int ySize, typename functionType, typename... Args >
233  std::tuple< Eigen::Matrix< double, ySize, 1 >, double > explicitEulerRichardsonWithErrorEstimator(
234  Eigen::Matrix< double, ySize, 1 > yN,
235  const double dt,
236  const double TOL,
237  functionType fRate,
238  Args&&... fRateArgs )
239  {
240 
241  typedef Eigen::Matrix< double, ySize, 1 > ySized;
242  ySized fN = fRate( yN, fRateArgs... );
243  ySized u = yN + fN * dt;
244  ySized v = yN + fN * dt / 2.;
245  ySized w = v + fRate( v, fRateArgs... ) * dt / 2.;
246  ySized yNew = 2. * w - u;
247 
248  // error estimator
249  const double AERR = 1.0;
250  const double aI = AERR / TOL;
251  const double rI = 1.0;
252  double scaling = 0;
253  ySized ESTVec = ySized::Zero();
254 
255  for ( int i = 0; i < ySize; i++ ) {
256  scaling = aI + rI * std::max( abs( yNew( i ) ), abs( yN( i ) ) );
257  ESTVec( i ) = abs( w( i ) - u( i ) ) / abs( scaling );
258  }
259 
260  const double EST = ESTVec.maxCoeff();
261  const double tauNew = dt * std::min( 2., std::max( 0.2, 0.9 * std::sqrt( TOL / EST ) ) );
262 
263  return std::make_tuple( yNew, tauNew );
264  }
265 
269  Matrix3d directionCosines( const Matrix3d& transformedCoordinateSystem );
270 
275 
279  Matrix3d orthonormalCoordinateSystem( const Vector3d& n1, const Vector3d& n2 );
280  } // namespace Math
281 } // namespace Marmot
Marmot::Math::radToDeg
constexpr double radToDeg(const double alpha)
Definition: MarmotMath.h:63
Marmot::Math::explicitEulerRichardsonWithErrorEstimator
std::tuple< Eigen::Matrix< double, ySize, 1 >, double > explicitEulerRichardsonWithErrorEstimator(Eigen::Matrix< double, ySize, 1 > yN, const double dt, const double TOL, functionType fRate, Args &&... fRateArgs)
Definition: MarmotMath.h:233
Marmot::Math::macaulyMatrix
Eigen::Matrix< double, nRows, nCols > macaulyMatrix(const Eigen::Matrix< double, nRows, nCols > &mat)
Definition: MarmotMath.h:142
MarmotConstants.h
Marmot::Math::exp
double exp(double x)
Definition: MarmotMath.cpp:13
Marmot::Math::degToRad
constexpr double degToRad(const double alpha)
Definition: MarmotMath.h:70
MarmotTypedefs.h
Marmot::Math::isNaN
bool isNaN(T x)
Definition: MarmotMath.h:44
Marmot::Math::heavisideExclude0
constexpr int heavisideExclude0(double scalar)
Definition: MarmotMath.h:89
Marmot::Vector3d
Eigen::Matrix< double, 3, 1 > Vector3d
Definition: MarmotTypedefs.h:42
Marmot::Math::makeReal
double makeReal(const double &value)
Definition: MarmotMath.cpp:37
Marmot::Math::linearInterpolation
double linearInterpolation(double x, double x0, double x1, double y0, double y1)
Definition: MarmotMath.cpp:7
Marmot::Matrix3d
Eigen::Matrix< double, 3, 3 > Matrix3d
Definition: MarmotTypedefs.h:40
Marmot
This file includes functions needed for calculations with stress and strain tensors written in voigt ...
Definition: MarmotTesting.h:32
Marmot::Math::explicitEulerRichardson
yType explicitEulerRichardson(yType yN, const double dt, functionType fRate, Args &&... fRateArgs)
Definition: MarmotMath.h:218
Marmot::Math::explicitEuler
yType explicitEuler(yType yN, const double dt, functionType fRate, Args &&... fRateArgs)
Definition: MarmotMath.h:156
Marmot::Math::heaviside
constexpr int heaviside(double scalar)
Definition: MarmotMath.h:84
Marmot::Math::centralDiff
Eigen::Matrix< double, nRows, nCols > centralDiff(functionType f, const Eigen::Matrix< double, nCols, 1 > &X)
Definition: MarmotMath.h:194
Marmot::Math::orthonormalCoordinateSystem
Matrix3d orthonormalCoordinateSystem(Vector3d &normalVector)
Definition: MarmotMath.cpp:53
Marmot::Math::macauly
constexpr double macauly(double scalar)
Definition: MarmotMath.h:77
Marmot::Math::directionCosines
Matrix3d directionCosines(const Matrix3d &transformedCoordinateSystem)
Definition: MarmotMath.cpp:96
Marmot::Constants::cubicRootEps
double cubicRootEps()
Definition: MarmotConstants.h:36
Marmot::Constants::Pi
constexpr double Pi
Definition: MarmotConstants.h:34
Marmot::Math::semiImplicitEuler
Eigen::Matrix< double, ySize, 1 > semiImplicitEuler(Eigen::Matrix< double, ySize, 1 > yN, const double dt, functionType fRate, Args &&... fRateArgs)
Definition: MarmotMath.h:167
Marmot::Math::getExponentPowerTen
int getExponentPowerTen(const double x)
Definition: MarmotMath.cpp:43
Marmot::Math::sgn
constexpr int sgn(const T &val) noexcept
Definition: MarmotMath.h:97