MarmotKelvinChain.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  * Alexander Dummer alexander.dummer@uibk.ac.at
15  *
16  * This file is part of the MAteRialMOdellingToolbox (marmot).
17  *
18  * This library is free software; you can redistribute it and/or
19  * modify it under the terms of the GNU Lesser General Public
20  * License as published by the Free Software Foundation; either
21  * version 2.1 of the License, or (at your option) any later version.
22  *
23  * The full text of the license can be found in the file LICENSE.md at
24  * the top level directory of marmot.
25  * ---------------------------------------------------------------------
26  */
27 
28 #pragma once
29 #include "Marmot/MarmotTypedefs.h"
30 #include "autodiff/forward/real.hpp"
31 #include <functional>
32 
33 namespace Marmot::Materials {
34 
35  namespace KelvinChain {
36 
37  typedef Eigen::VectorXd Properties;
38  typedef Eigen::Map< Properties > mapProperties;
39 
40  typedef Eigen::Matrix< double, 6, Eigen::Dynamic > StateVarMatrix;
41  typedef Eigen::Map< StateVarMatrix > mapStateVarMatrix;
42 
43  template < int N >
44  struct Factorial {
45  enum { value = N * Factorial< N - 1 >::value };
46  };
47 
48  template <>
49  struct Factorial< 0 > {
50  enum { value = 1 };
51  };
52 
53  template < int k >
54  double evaluatePostWidderFormula( std::function< autodiff::Real< k, double >( autodiff::Real< k, double > ) > phi,
55  double tau )
56  {
57 
58  autodiff::Real< k, double > tau_( tau * k );
59 
60  double val = -pow( -tau * k, k ) / double( Factorial< k - 1 >::value );
61  val *= autodiff::derivatives( phi, autodiff::along( 1. ), autodiff::at( tau_ ) )[k];
62  return val;
63  }
64 
65  template < int k >
66  Properties computeElasticModuli( std::function< autodiff::Real< k, double >( autodiff::Real< k, double > ) > phi,
67  Properties retardationTimes )
68  {
69  Properties elasticModuli( retardationTimes.size() );
70  double spacing = log( retardationTimes( 1 ) / retardationTimes( 0 ) );
71 
72  for ( int i = 0; i < retardationTimes.size(); i++ ) {
73  double tau = retardationTimes( i );
74  elasticModuli( i ) = 1. / ( spacing * evaluatePostWidderFormula< k >( phi, tau ) );
75  }
76 
77  return elasticModuli;
78  }
79 
80  Properties generateRetardationTimes( int n, double min, double spacing );
81 
82  void updateStateVarMatrix( const double dT,
83  Properties elasticModuli,
84  Properties retardationTimes,
85  Eigen::Ref< StateVarMatrix > stateVars,
86  const Marmot::Vector6d& dStress,
87  const Marmot::Matrix6d& unitComplianceMatrix );
88 
89  void evaluateKelvinChain( const double dT,
90  Properties elasticModuli,
91  Properties retardationTimes,
92  StateVarMatrix stateVars,
93  double& uniaxialCompliance,
94  Marmot::Vector6d& dStrain,
95  const double factor );
96 
97  void computeLambdaAndBeta( double dT, double tau, double& lambda, double& beta );
98 
99  } // namespace KelvinChain
100 } // namespace Marmot::Materials
Marmot::Materials::KelvinChain::evaluatePostWidderFormula
double evaluatePostWidderFormula(std::function< autodiff::Real< k, double >(autodiff::Real< k, double >) > phi, double tau)
Definition: MarmotKelvinChain.h:54
Marmot::Materials::KelvinChain::evaluateKelvinChain
void evaluateKelvinChain(const double dT, Properties elasticModuli, Properties retardationTimes, StateVarMatrix stateVars, double &uniaxialCompliance, Marmot::Vector6d &dStrain, const double factor)
Definition: MarmotKelvinChain.cpp:18
Marmot::Materials::KelvinChain::computeLambdaAndBeta
void computeLambdaAndBeta(double dT, double tau, double &lambda, double &beta)
Definition: MarmotKelvinChain.cpp:57
MarmotTypedefs.h
Marmot::Materials::KelvinChain::mapStateVarMatrix
Eigen::Map< StateVarMatrix > mapStateVarMatrix
Definition: MarmotKelvinChain.h:41
Marmot::Matrix6d
Eigen::Matrix< double, 6, 6 > Matrix6d
Definition: MarmotTypedefs.h:35
Marmot::Materials
Definition: MarmotKelvinChain.h:33
Marmot::Materials::KelvinChain::mapProperties
Eigen::Map< Properties > mapProperties
Definition: MarmotKelvinChain.h:38
Marmot::Materials::KelvinChain::updateStateVarMatrix
void updateStateVarMatrix(const double dT, Properties elasticModuli, Properties retardationTimes, Eigen::Ref< StateVarMatrix > stateVars, const Marmot::Vector6d &dStress, const Marmot::Matrix6d &unitComplianceMatrix)
Marmot::Materials::KelvinChain::generateRetardationTimes
Properties generateRetardationTimes(int n, double min, double spacing)
Definition: MarmotKelvinChain.cpp:10
Marmot::Materials::KelvinChain::StateVarMatrix
Eigen::Matrix< double, 6, Eigen::Dynamic > StateVarMatrix
Definition: MarmotKelvinChain.h:40
Marmot::Materials::KelvinChain::Properties
Eigen::VectorXd Properties
Definition: MarmotKelvinChain.h:37
Marmot::Vector6d
Eigen::Matrix< double, 6, 1 > Vector6d
Definition: MarmotTypedefs.h:43
Marmot::Materials::KelvinChain::Factorial::value
@ value
Definition: MarmotKelvinChain.h:45
Marmot::Materials::KelvinChain::Factorial
Definition: MarmotKelvinChain.h:44
Marmot::Materials::KelvinChain::computeElasticModuli
Properties computeElasticModuli(std::function< autodiff::Real< k, double >(autodiff::Real< k, double >) > phi, Properties retardationTimes)
Definition: MarmotKelvinChain.h:66
Marmot::FiniteElement::Spatial1D::Bar2::N
NSized N(double xi)
Definition: MarmotFiniteElement1D.cpp:15