MarmotMaterialFiniteStrain.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  * Alexander Dummer alexander.dummer@uibk.ac.at
16  *
17  * This file is part of the MAteRialMOdellingToolbox (marmot).
18  *
19  * This library is free software; you can redistribute it and/or
20  * modify it under the terms of the GNU Lesser General Public
21  * License as published by the Free Software Foundation; either
22  * version 2.1 of the License, or (at your option) any later version.
23  *
24  * The full text of the license can be found in the file LICENSE.md at
25  * the top level directory of marmot.
26  * ---------------------------------------------------------------------
27  */
28 #pragma once
29 #include "Fastor/Fastor.h"
30 #include "Marmot/MarmotMaterial.h"
31 #include <Fastor/tensor/Tensor.h>
32 
34 
35  /*
36  Abstract basic class for mechanical materials in the finite strain regime
37  */
38 public:
39  template < int nDim >
41  Fastor::Tensor< double, nDim, nDim > tau; // kirchhoff stress
42  double rho; // density
43  double elasticEnergyDensity; // elastic energy per unit volume
44  };
45 
46  template < int nDim >
48  Fastor::Tensor< double, nDim, nDim, nDim, nDim > dTau_dF; // tangent operator w.r.t. deformation gradient
49  };
50 
51  template < int nDim >
52  struct Deformation {
53  Fastor::Tensor< double, nDim, nDim > F;
54  };
55 
56  struct TimeIncrement {
57  const double time;
58  const double dT;
59  };
60 
62 
63  virtual void computeStress( ConstitutiveResponse< 3 >& response,
64  AlgorithmicModuli< 3 >& tangents,
65  const Deformation< 3 >&,
66  const TimeIncrement& ) = 0;
67 
71  virtual void computeStress( ConstitutiveResponse< 3 >& response,
72  AlgorithmicModuli< 3 >& tangents,
73  const Deformation< 3 >& deformation,
74  const TimeIncrement& timeIncrement,
75  const std::tuple< double, double, double >& eigenDeformation );
76 
77  virtual void computePlaneStrain( ConstitutiveResponse< 3 >& response,
78  AlgorithmicModuli< 3 >& algorithmicModuli,
79  const Deformation< 3 >& deformation,
80  const TimeIncrement& timeIncrement );
81  /***/
82  /* * Compute stress under plane strain conditions, but account for eigen deformations (e.g, geostatic stress
83  * states).*/
84  /* * Modifies algorithmic tangent accordingly.*/
85  virtual void computePlaneStrain( ConstitutiveResponse< 3 >& response,
86  AlgorithmicModuli< 3 >& algorithmicModuli,
87  const Deformation< 3 >& deformation,
88  const TimeIncrement& timeIncrement,
89  const std::tuple< double, double, double >& eigenDeformation );
90 
91  virtual void computePlaneStress( ConstitutiveResponse< 2 >& response,
92  AlgorithmicModuli< 2 >& algorithmicModuli,
93  const Deformation< 2 >& deformation,
94  const TimeIncrement& timeIncrement );
97  std::tuple< double, double, double > findEigenDeformationForEigenStress(
98  const std::tuple< double, double, double >& initialGuess,
99  const std::tuple< double, double, double >& eigenStress );
100 };
MarmotMaterialFiniteStrain
Definition: MarmotMaterialFiniteStrain.h:33
MarmotMaterialFiniteStrain::Deformation::F
Fastor::Tensor< double, nDim, nDim > F
Definition: MarmotMaterialFiniteStrain.h:53
MarmotMaterialFiniteStrain::AlgorithmicModuli::dTau_dF
Fastor::Tensor< double, nDim, nDim, nDim, nDim > dTau_dF
Definition: MarmotMaterialFiniteStrain.h:48
MarmotMaterial::MarmotMaterial
MarmotMaterial(const double *materialProperties, int nMaterialProperties, int materialNumber)
MarmotMaterialFiniteStrain::ConstitutiveResponse
Definition: MarmotMaterialFiniteStrain.h:40
MarmotMaterialFiniteStrain::Deformation
Definition: MarmotMaterialFiniteStrain.h:52
MarmotMaterialFiniteStrain::findEigenDeformationForEigenStress
std::tuple< double, double, double > findEigenDeformationForEigenStress(const std::tuple< double, double, double > &initialGuess, const std::tuple< double, double, double > &eigenStress)
Definition: MarmotMaterialFiniteStrain.cpp:78
MarmotMaterialFiniteStrain::ConstitutiveResponse::rho
double rho
Definition: MarmotMaterialFiniteStrain.h:42
MarmotMaterialFiniteStrain::TimeIncrement
Definition: MarmotMaterialFiniteStrain.h:56
MarmotMaterialFiniteStrain::AlgorithmicModuli
Definition: MarmotMaterialFiniteStrain.h:47
MarmotMaterialFiniteStrain::computePlaneStrain
virtual void computePlaneStrain(ConstitutiveResponse< 3 > &response, AlgorithmicModuli< 3 > &algorithmicModuli, const Deformation< 3 > &deformation, const TimeIncrement &timeIncrement)
Definition: MarmotMaterialFiniteStrain.cpp:51
MarmotMaterial
Definition: MarmotMaterial.h:32
MarmotMaterialFiniteStrain::TimeIncrement::time
const double time
Definition: MarmotMaterialFiniteStrain.h:57
MarmotMaterialFiniteStrain::TimeIncrement::dT
const double dT
Definition: MarmotMaterialFiniteStrain.h:58
MarmotMaterialFiniteStrain::ConstitutiveResponse::elasticEnergyDensity
double elasticEnergyDensity
Definition: MarmotMaterialFiniteStrain.h:43
MarmotMaterial.h
MarmotMaterialFiniteStrain::computePlaneStress
virtual void computePlaneStress(ConstitutiveResponse< 2 > &response, AlgorithmicModuli< 2 > &algorithmicModuli, const Deformation< 2 > &deformation, const TimeIncrement &timeIncrement)
Definition: MarmotMaterialFiniteStrain.cpp:68
MarmotMaterialFiniteStrain::computeStress
virtual void computeStress(ConstitutiveResponse< 3 > &response, AlgorithmicModuli< 3 > &tangents, const Deformation< 3 > &, const TimeIncrement &)=0
MarmotMaterialFiniteStrain::ConstitutiveResponse::tau
Fastor::Tensor< double, nDim, nDim > tau
Definition: MarmotMaterialFiniteStrain.h:41