Skip to content
GitLab
Explore
Sign in
Primary navigation
Search or go to…
Project
D
dune-elasticity
Manage
Activity
Members
Labels
Plan
Issues
Issue boards
Milestones
Wiki
Requirements
Code
Merge requests
Repository
Branches
Commits
Tags
Repository graph
Compare revisions
Snippets
Locked files
Build
Pipelines
Jobs
Pipeline schedules
Test cases
Artifacts
Deploy
Releases
Package registry
Container registry
Model registry
Operate
Environments
Terraform modules
Monitor
Incidents
Analyze
Value stream analytics
Contributor analytics
CI/CD analytics
Repository analytics
Code review analytics
Issue analytics
Insights
Model experiments
Help
Help
Support
GitLab documentation
Compare GitLab plans
Community forum
Contribute to GitLab
Provide feedback
Keyboard shortcuts
?
Snippets
Groups
Projects
Show more breadcrumbs
Patrick Jaap
dune-elasticity
Commits
1190401c
Commit
1190401c
authored
11 years ago
by
akbib
Committed by
akbib
11 years ago
Browse files
Options
Downloads
Patches
Plain Diff
Inherit from new abstract base class + cleanup and renaming.
[[Imported from SVN: r12666]]
parent
76ee7b69
Branches
Branches containing commit
Tags
Tags containing commit
No related merge requests found
Changes
1
Hide whitespace changes
Inline
Side-by-side
Showing
1 changed file
dune/elasticity/materials/geomexactstvenantkirchhoffmaterial.hh
+56
-95
56 additions, 95 deletions
...lasticity/materials/geomexactstvenantkirchhoffmaterial.hh
with
56 additions
and
95 deletions
dune/elasticity/materials/geomexactstvenantkirchhoffmaterial.hh
+
56
−
95
View file @
1190401c
#ifndef GEOMETRICALLY_EXACT_ST_VENANT_KIRCHHOFF_MATERIAL
#define GEOMETRICALLY_EXACT_ST_VENANT_KIRCHHOFF_MATERIAL
// -*- tab-width: 8; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set ts=8 sw=4 et sts=4:
#ifndef DUNE_ELASTICITY_MATERIALS_GEOMETRICALLY_EXACT_ST_VENANT_KIRCHHOFF_MATERIAL_HH
#define DUNE_ELASTICITY_MATERIALS_GEOMETRICALLY_EXACT_ST_VENANT_KIRCHHOFF_MATERIAL_HH
#include
<dune/fufem/quadraturerules/quadraturerulecache.hh>
#include
<dune/fufem/symmetrictensor.hh>
#include
<dune/fufem/functions/virtualgridfunction.hh>
#include
<dune/fufem/functions/basisgridfunction.hh>
#include
<dune/fufem/assemblers/localassemblers/geononlinstvenantfunctionalassembler.hh>
#include
<dune/fufem/assemblers/localassemblers/geononlinlinearizedstvenantassembler.hh>
#include
<dune/elasticity/materials/material.hh>
#include
<dune/elasticity/common/elasticityhelpers.hh>
/* \brief Class representing a hyperelastic geometrically exact St.Venant Kirchhoff material.
*
* \tparam Basis Global basis that is used for the spatial discretization.
* (This is not nice but I need a LocalFiniteElement for the local Hessian assembler :-( )
*/
template
<
class
Basis
>
class
GeomExactStVenantMaterial
class
GeomExactStVenantMaterial
:
public
Material
<
Basis
>
{
typedef
Material
<
Basis
>
Base
;
public:
typedef
typename
Basis
::
GridView
::
Grid
GridType
;
typedef
typename
Basis
::
LocalFiniteElement
Lfe
;
typedef
typename
Base
::
GlobalBasis
GlobalBasis
;
typedef
typename
Base
::
Lfe
Lfe
;
typedef
typename
Base
::
LocalLinearization
LocalLinearization
;
typedef
typename
Base
::
LocalHessian
LocalHessian
;
typedef
typename
Base
::
VectorType
VectorType
;
typedef
typename
Base
::
GridFunction
GridFunction
;
typedef
typename
Base
::
ReturnType
ReturnType
;
typedef
GeomNonlinElasticityFunctionalAssembler
<
GridType
,
Lfe
>
LocalLinearization
;
typedef
GeomNonlinLinearizedStVenantAssembler
<
GridType
,
Lfe
,
Lfe
>
LocalHessian
;
typedef
Basis
GlobalBasis
;
private:
using
Base
::
dim
;
typedef
typename
GridType
::
ctype
ctype
;
static
const
int
dim
=
GridType
::
dimension
;
static
const
int
dimworld
=
GridType
::
dimensionworld
;
typedef
GeomNonlinElasticityFunctionalAssembler
<
GridType
,
Lfe
>
GeomLinearization
;
typedef
GeomNonlinLinearizedStVenantAssembler
<
GridType
,
Lfe
,
Lfe
>
GeomHessian
;
typedef
typename
GridType
::
template
Codim
<
0
>
::
Geometry
::
LocalCoordinate
LocalCoordinate
;
typedef
typename
GridType
::
template
Codim
<
0
>
::
LeafIterator
ElementIterator
;
...
...
@@ -41,51 +46,31 @@ public:
nu_
(
0.3
)
{}
GeomExactStVenantMaterial
(
const
Basis
&
basis
,
ctype
E
,
ctype
nu
)
:
localLinearization_
(
E
,
nu
),
localHessian_
(
E
,
nu
),
basis_
(
basis
),
GeomExactStVenantMaterial
(
const
Basis
&
basis
,
ReturnType
E
,
ReturnType
nu
)
:
Base
(
basis
),
E_
(
E
),
nu_
(
nu
)
{}
void
getMaterialParameters
(
ctype
&
E
,
ctype
&
nu
)
{
E
=
E_
;
nu
=
nu_
;
localLinearization_
=
Dune
::
shared_ptr
<
GeomLinearization
>
(
new
GeomLinearization
(
E
,
nu
))
;
localHessian_
=
Dune
::
shared_ptr
<
GeomHessian
>
(
new
GeomHessian
(
E
,
nu
))
;
}
void
setMaterialParameters
(
ctype
E
,
ctype
nu
)
{
E_
=
E
;
nu_
=
nu
;
}
void
setup
(
ctype
E
,
ctype
nu
,
const
Basis
&
basis
)
void
setup
(
ReturnType
E
,
ReturnType
nu
,
const
Basis
&
basis
)
{
E_
=
E
;
nu_
=
nu
;
if
(
localLinearization_
)
localLinearization_
.
reset
();
if
(
localHessian_
)
localHessian_
.
reset
();
localLinearization_
=
Dune
::
shared_ptr
<
LocalLinearization
>
(
new
LocalLinearization
(
E
,
nu
));
localHessian_
=
Dune
::
shared_ptr
<
LocalHessian
>
(
new
LocalHessian
(
E
,
nu
));
basis_
=
&
basis
;
localLinearization_
=
Dune
::
shared_ptr
<
GeomLinearization
>
(
new
GeomLinearization
(
E
,
nu
));
localHessian_
=
Dune
::
shared_ptr
<
GeomHessian
>
(
new
GeomHessian
(
E
,
nu
));
this
->
basis_
=
&
basis
;
}
//! Evaluate the strain energy
template
<
class
CoeffType
>
ctype
energy
(
const
CoeffType
&
coeff
)
ReturnType
energy
(
const
GridFunction
&
displace
)
{
// make grid function
BasisGridFunction
<
Basis
,
CoeffType
>
configuration
(
*
basis_
,
coeff
);
ctype
energy
=
0
;
const
GridType
&
grid
=
configuration
.
grid
();
const
GridType
&
grid
=
this
->
basis_
->
getGridView
()
.
grid
();
ElementIterator
eIt
=
grid
.
template
leafbegin
<
0
>();
ElementIterator
eItEnd
=
grid
.
template
leafend
<
0
>();
...
...
@@ -93,7 +78,7 @@ public:
for
(;
eIt
!=
eItEnd
;
++
eIt
)
{
// get quadrature rule
QuadratureRuleKey
quad1
(
basis_
->
getLocalFiniteElement
(
*
eIt
));
QuadratureRuleKey
quad1
(
this
->
basis_
->
getLocalFiniteElement
(
*
eIt
));
QuadratureRuleKey
quadKey
=
quad1
.
derivative
().
square
().
square
();
const
Dune
::
template
QuadratureRule
<
ctype
,
dim
>
&
quad
=
QuadratureRuleCache
<
ctype
,
dim
>::
rule
(
quadKey
);
...
...
@@ -108,23 +93,22 @@ public:
const
ctype
integrationElement
=
eIt
->
geometry
().
integrationElement
(
quadPos
);
// evaluate displacement gradient at the quadrature point
typename
Basis
GridFunction
<
Basis
,
CoeffType
>
::
DerivativeType
local
Conf
Grad
;
typename
GridFunction
::
DerivativeType
local
Disp
Grad
;
if
(
configuration
.
isDefinedOn
(
*
eIt
))
configuration
.
evaluateDerivativeLocal
(
*
eIt
,
quadPos
,
local
Conf
Grad
);
if
(
displace
->
isDefinedOn
(
*
eIt
))
displace
->
evaluateDerivativeLocal
(
*
eIt
,
quadPos
,
local
Disp
Grad
);
else
configuration
.
evaluateDerivative
(
eIt
->
geometry
().
global
(
quadPos
),
local
Conf
Grad
);
displace
->
evaluateDerivative
(
eIt
->
geometry
().
global
(
quadPos
),
local
Disp
Grad
);
// Compute the nonlinear strain tensor from the deformation gradient
SymmetricTensor
<
dim
>
strain
;
compute
Strain
(
local
Conf
Grad
,
strain
);
SymmetricTensor
<
dim
,
ReturnType
>
strain
;
Dune
::
Elasticity
::
computeNonlinear
Strain
(
local
Disp
Grad
,
strain
);
// and the stress
SymmetricTensor
<
dim
>
stress
=
hookeTimesStrain
(
strain
);
SymmetricTensor
<
dim
,
ReturnType
>
stress
=
hookeTimesStrain
(
strain
);
ctype
z
=
quad
[
pt
].
weight
()
*
integrationElement
;
energy
+=
(
stress
*
strain
)
*
z
;
}
}
...
...
@@ -132,67 +116,44 @@ public:
}
//! Return the local assembler of the first derivative of the strain energy
LocalLinearization
&
firstDerivative
()
{
return
*
localLinearization_
;}
LocalLinearization
&
firstDerivative
(
const
GridFunction
&
displace
)
{
localLinearization_
->
setConfiguration
(
displace
);
return
*
localLinearization_
;
}
//! Return the local assembler of the second derivative of the strain energy
LocalHessian
&
secondDerivative
()
{
return
*
localHessian_
;}
//! Return the global basis
const
Basis
&
basis
()
{
return
*
basis_
;}
LocalHessian
&
secondDerivative
(
const
GridFunction
&
displace
)
{
localHessian_
->
setConfiguration
(
displace
);
return
*
localHessian_
;
}
private
:
//! First derivative of the strain energy
Dune
::
shared_ptr
<
Local
Linearization
>
localLinearization_
;
Dune
::
shared_ptr
<
Geom
Linearization
>
localLinearization_
;
//! Second derivative of the strain energy
Dune
::
shared_ptr
<
LocalHessian
>
localHessian_
;
//! Global basis used for the spatial discretization
const
Basis
*
basis_
;
Dune
::
shared_ptr
<
GeomHessian
>
localHessian_
;
//! Elasticity modulus
ct
ype
E_
;
ReturnT
ype
E_
;
//! Poisson ratio
ct
ype
nu_
;
ReturnT
ype
nu_
;
//! Compute nonlinear strain tensor from the deformation gradient.
void
computeStrain
(
const
Dune
::
FieldMatrix
<
ctype
,
dim
,
dim
>&
grad
,
SymmetricTensor
<
dim
>&
strain
)
const
{
strain
=
0
;
for
(
int
i
=
0
;
i
<
dim
;
++
i
)
{
strain
(
i
,
i
)
+=
grad
[
i
][
i
];
for
(
int
k
=
0
;
k
<
dim
;
++
k
)
strain
(
i
,
i
)
+=
0.5
*
grad
[
k
][
i
]
*
grad
[
k
][
i
];
for
(
int
j
=
i
+
1
;
j
<
dim
;
++
j
)
{
strain
(
i
,
j
)
+=
0.5
*
(
grad
[
i
][
j
]
+
grad
[
j
][
i
]);
for
(
int
k
=
0
;
k
<
dim
;
++
k
)
strain
(
i
,
j
)
+=
0.5
*
grad
[
k
][
i
]
*
grad
[
k
][
j
];
}
}
}
//! Compute linear elastic stress from the strain
SymmetricTensor
<
dim
>
hookeTimesStrain
(
const
SymmetricTensor
<
dim
>&
strain
)
const
{
SymmetricTensor
<
dim
,
ReturnType
>
hookeTimesStrain
(
const
SymmetricTensor
<
dim
,
ReturnType
>&
strain
)
const
{
SymmetricTensor
<
dim
>
stress
=
strain
;
SymmetricTensor
<
dim
,
ReturnType
>
stress
=
strain
;
stress
*=
E_
/
(
1
+
nu_
);
// Compute the trace of the strain
ctype
trace
=
0
;
for
(
int
i
=
0
;
i
<
dim
;
i
++
)
trace
+=
strain
(
i
,
i
);
ctype
f
=
E_
*
nu_
/
((
1
+
nu_
)
*
(
1
-
2
*
nu_
))
*
trace
;
ReturnType
f
=
E_
*
nu_
/
((
1
+
nu_
)
*
(
1
-
2
*
nu_
))
*
strain
.
trace
();
for
(
int
i
=
0
;
i
<
dim
;
i
++
)
stress
(
i
,
i
)
+=
f
;
return
stress
;
}
};
#endif
This diff is collapsed.
Click to expand it.
Preview
0%
Loading
Try again
or
attach a new file
.
Cancel
You are about to add
0
people
to the discussion. Proceed with caution.
Finish editing this message first!
Save comment
Cancel
Please
register
or
sign in
to comment