Loading [MathJax]/extensions/tex2jax.js
Plugins
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
hooke.h
Go to the documentation of this file.
1/*
2 * Copyright (c) 2014-2015: G-CSC, Goethe University Frankfurt
3 * Author: Raphael Prohl
4 *
5 * This file is part of UG4.
6 *
7 * UG4 is free software: you can redistribute it and/or modify it under the
8 * terms of the GNU Lesser General Public License version 3 (as published by the
9 * Free Software Foundation) with the following additional attribution
10 * requirements (according to LGPL/GPL v3 §7):
11 *
12 * (1) The following notice must be displayed in the Appropriate Legal Notices
13 * of covered and combined works: "Based on UG4 (www.ug4.org/license)".
14 *
15 * (2) The following notice must be displayed at a prominent place in the
16 * terminal output of covered works: "Based on UG4 (www.ug4.org/license)".
17 *
18 * (3) The following bibliography is recommended for citation and must be
19 * preserved in all covered files:
20 * "Reiter, S., Vogel, A., Heppner, I., Rupp, M., and Wittum, G. A massively
21 * parallel geometric multigrid solver on hierarchically distributed grids.
22 * Computing and visualization in science 16, 4 (2013), 151-164"
23 * "Vogel, A., Reiter, S., Rupp, M., Nägel, A., and Wittum, G. UG4 -- a novel
24 * flexible software system for simulating pde based models on high performance
25 * computers. Computing and visualization in science 16, 4 (2013), 165-179"
26 *
27 * This program is distributed in the hope that it will be useful,
28 * but WITHOUT ANY WARRANTY; without even the implied warranty of
29 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
30 * GNU Lesser General Public License for more details.
31 */
32
33#ifndef HOOKE_H_
34#define HOOKE_H_
35
36namespace ug{
37namespace SmallStrainMechanics{
38
42
76template <typename TDomain>
78 : public IMaterialLaw<TDomain>
79{
80 private:
83
84 public:
86 static const int dim = base_type::dim;
87
88 public:
90 HookeLaw(): IMaterialLaw<TDomain>()
91 {
93 };
94
97
98 public:
100 // INTERFACE-METHODS
102 virtual void init();
103
105 virtual void stressTensor(MathMatrix<dim,dim>& stressTens, const size_t ip,
106 const MathMatrix<dim, dim>& GradU);
107
111
113 elasticityTensor(const size_t ip, const MathMatrix<dim, dim>& GradU){
114 return elasticityTensor();};
115
116
117 public:
120 void set_hooke_elasticity_tensor(const number lambda, const number mu);
121 void set_hooke_elasticity_tensor_E_nu(const number E, const number nu);
122
126
130 const number C11, const number C12, const number C13,
131 const number C22, const number C23,
132 const number C33,
133 const number C44,
134 const number C55,
135 const number C66 );
136
138 const number E1, const number E2, const number E3,
139 const number G12, const number G13, const number G23,
140 const number v12, const number v13, const number v23);
141
143 const number E1, const number E2,
144 const number G12,
145 const number v12);
146
148 const number E1, const number E2, const number E3,
149 const number G12, const number G13, const number G23,
150 const number v12, const number v13, const number v23
151 );
153
154 public:
156
157 public:
158 void strainTensor(MathMatrix<dim,dim>& strainTens, const MathMatrix<dim, dim>& GradU);
159
160 protected:
163};
164
165
166}// end of namespace SmallStrainMechanics
167}// end of namespace ug
168
169#include "hooke_impl.h"
170
171#endif /* HOOKE_H_ */
Definition mat_law_interface.h:44
static const int dim
World dimension.
Definition mat_law_interface.h:47
bool m_bConstElastTens
flag indicating, if elasticity tensor is constant
Definition mat_law_interface.h:131
std::string m_materialConfiguration
Definition mat_law_interface.h:124
void set_elasticity_tensor_orthotropic_E_G_nu(const number E1, const number E2, const number E3, const number G12, const number G13, const number G23, const number v12, const number v13, const number v23)
Definition hooke_impl.h:240
virtual SmartPtr< MathTensor4< TDomain::dim, TDomain::dim, TDomain::dim, TDomain::dim > > elasticityTensor()
computes the elasticity tensor; commonly denoted by C
Definition hooke_impl.h:99
static const int dim
World dimension.
Definition hooke.h:86
void set_hooke_elasticity_tensor_plain_strain_E_nu(const number E, const number nu)
Definition hooke_impl.h:396
SmartPtr< MathTensor4< dim, dim, dim, dim > > m_spElastTensorFunct
elasticity tensor
Definition hooke.h:162
IMaterialLaw< TDomain > base_type
Base class type.
Definition hooke.h:82
void strainTensor(MathMatrix< dim, dim > &strainTens, const MathMatrix< dim, dim > &GradU)
Definition hooke_impl.h:109
void set_hooke_elasticity_tensor_E_nu(const number E, const number nu)
Definition hooke_impl.h:407
void set_hooke_elasticity_tensor(const number lambda, const number mu)
Definition hooke_impl.h:417
void set_elasticity_tensor_orthotropic(const number C11, const number C12, const number C13, const number C22, const number C23, const number C33, const number C44, const number C55, const number C66)
Definition hooke_impl.h:294
virtual SmartPtr< MathTensor4< TDomain::dim, TDomain::dim, TDomain::dim, TDomain::dim > > elasticityTensor(const size_t ip, const MathMatrix< dim, dim > &GradU)
Definition hooke.h:113
void set_elasticity_tensor_orthotropic_plain_strain_E_G_nu(const number E1, const number E2, const number E3, const number G12, const number G13, const number G23, const number v12, const number v13, const number v23)
Definition hooke_impl.h:180
virtual void stressTensor(MathMatrix< dim, dim > &stressTens, const size_t ip, const MathMatrix< dim, dim > &GradU)
computes the cauchy stress tensor sigma at an integration point 'ip'
Definition hooke_impl.h:71
void set_hooke_elasticity_tensor_plain_stress_E_nu(const number E, const number nu)
Definition hooke_impl.h:355
void set_elasticity_tensor_orthotropic_plain_stress_E_G_nu(const number E1, const number E2, const number G12, const number v12)
Definition hooke_impl.h:121
~HookeLaw()
Destructor.
Definition hooke.h:96
virtual void init()
Definition hooke_impl.h:59
HookeLaw()
constructor
Definition hooke.h:90
double number