Skip to content
GitLab
Explore
Sign in
Primary navigation
Search or go to…
Project
D
dune-tectonic
Manage
Activity
Members
Labels
Plan
Issues
Issue boards
Milestones
Iterations
Wiki
Requirements
Code
Merge requests
Repository
Branches
Commits
Tags
Repository graph
Compare revisions
Locked files
Deploy
Releases
Model registry
Monitor
Incidents
Analyze
Value stream analytics
Contributor 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
agnumpde
dune-tectonic
Commits
3360ff17
Commit
3360ff17
authored
11 years ago
by
Elias Pipping
Browse files
Options
Downloads
Patches
Plain Diff
[Cleanup] Outsource more assembly; Turn header into class
parent
279be7cb
No related branches found
No related tags found
No related merge requests found
Changes
4
Hide whitespace changes
Inline
Side-by-side
Showing
4 changed files
src/assemblers.cc
+91
-34
91 additions, 34 deletions
src/assemblers.cc
src/assemblers.hh
+49
-19
49 additions, 19 deletions
src/assemblers.hh
src/assemblers_tmpl.cc
+1
-22
1 addition, 22 deletions
src/assemblers_tmpl.cc
src/one-body-sample.cc
+46
-86
46 additions, 86 deletions
src/one-body-sample.cc
with
187 additions
and
161 deletions
src/assemblers.cc
+
91
−
34
View file @
3360ff17
...
@@ -2,56 +2,113 @@
...
@@ -2,56 +2,113 @@
#include
"config.h"
#include
"config.h"
#endif
#endif
#include
<dune/istl/scaledidmatrix.hh>
#include
<dune/fufem/boundarypatch.hh>
#include
<dune/fufem/boundarypatch.hh>
#include
<dune/fufem/functions/constantfunction.hh>
#include
<dune/fufem/functions/constantfunction.hh>
#include
<dune/fufem/assemblers/localassemblers/boundarymassassembler.hh>
#include
<dune/fufem/assemblers/localassemblers/l2functionalassembler.hh>
#include
<dune/fufem/assemblers/localassemblers/massassembler.hh>
#include
<dune/fufem/assemblers/localassemblers/neumannboundaryassembler.hh>
#include
<dune/fufem/assemblers/localassemblers/neumannboundaryassembler.hh>
#include
<dune/fufem/assemblers/localassemblers/stvenantkirchhoffassembler.hh>
#include
<dune/fufem/assemblers/localassemblers/viscosityassembler.hh>
#include
<dune/tectonic/globalruinanonlinearity.hh>
#include
<dune/tectonic/globalruinanonlinearity.hh>
#include
"assemblers.hh"
#include
"assemblers.hh"
// Assembles Neumann boundary term in f
template
<
class
GridView
,
int
dimension
>
template
<
class
GridView
,
class
LocalVector
,
class
Assembler
>
MyAssembler
<
GridView
,
dimension
>::
MyAssembler
(
GridView
const
&
_gridView
)
void
assembleNeumann
(
GridView
const
&
gridView
,
Assembler
const
&
assembler
,
:
vertexBasis
(
_gridView
),
Dune
::
BitSetVector
<
1
>
const
&
neumannNodes
,
gridView
(
_gridView
),
Dune
::
BlockVector
<
LocalVector
>
&
f
,
assembler
(
vertexBasis
,
vertexBasis
)
{}
Dune
::
VirtualFunction
<
double
,
double
>
const
&
neumann
,
double
relativeTime
)
{
// constant sample function on
template
<
class
GridView
,
int
dimension
>
// neumann boundary
void
MyAssembler
<
GridView
,
dimension
>::
assembleFrictionalBoundaryMass
(
BoundaryPatch
<
GridView
>
const
neumannBoundary
(
gridView
,
neumannNodes
);
BoundaryPatch
<
GridView
>
const
&
frictionalBoundary
,
ScalarMatrix
&
frictionalBoundaryMass
)
{
BoundaryMassAssembler
<
Grid
,
BoundaryPatch
<
GridView
>
,
LocalVertexBasis
,
LocalVertexBasis
,
Dune
::
FieldMatrix
<
double
,
1
,
1
>>
const
frictionalBoundaryMassAssembler
(
frictionalBoundary
);
assembler
.
assembleOperator
(
frictionalBoundaryMassAssembler
,
frictionalBoundaryMass
);
}
template
<
class
GridView
,
int
dimension
>
void
MyAssembler
<
GridView
,
dimension
>::
assembleMass
(
double
density
,
Matrix
&
M
)
{
MassAssembler
<
Grid
,
LocalVertexBasis
,
LocalVertexBasis
,
Dune
::
ScaledIdentityMatrix
<
double
,
dimension
>>
const
localMass
;
assembler
.
assembleOperator
(
localMass
,
M
);
M
*=
density
;
}
template
<
class
GridView
,
int
dimension
>
void
MyAssembler
<
GridView
,
dimension
>::
assembleElasticity
(
double
E
,
double
nu
,
Matrix
&
A
)
{
StVenantKirchhoffAssembler
<
Grid
,
LocalVertexBasis
,
LocalVertexBasis
>
const
localStiffness
(
E
,
nu
);
assembler
.
assembleOperator
(
localStiffness
,
A
);
}
template
<
class
GridView
,
int
dimension
>
void
MyAssembler
<
GridView
,
dimension
>::
assembleViscosity
(
double
shearViscosity
,
double
bulkViscosity
,
Matrix
&
C
)
{
ViscosityAssembler
<
Grid
,
LocalVertexBasis
,
LocalVertexBasis
>
const
localViscosity
(
shearViscosity
,
bulkViscosity
);
assembler
.
assembleOperator
(
localViscosity
,
C
);
}
template
<
class
GridView
,
int
dimension
>
void
MyAssembler
<
GridView
,
dimension
>::
assembleBodyForce
(
double
gravity
,
double
density
,
LocalVector
zenith
,
Vector
&
f
)
{
LocalVector
weightedGravitationalDirection
=
zenith
;
weightedGravitationalDirection
*=
density
*
gravity
;
weightedGravitationalDirection
*=
-
1
;
ConstantFunction
<
LocalVector
,
LocalVector
>
const
gravityFunction
(
weightedGravitationalDirection
);
L2FunctionalAssembler
<
Grid
,
LocalVertexBasis
,
LocalVector
>
gravityFunctionalAssembler
(
gravityFunction
);
assembler
.
assembleFunctional
(
gravityFunctionalAssembler
,
f
);
}
template
<
class
GridView
,
int
dimension
>
void
MyAssembler
<
GridView
,
dimension
>::
assembleNeumann
(
BoundaryPatch
<
GridView
>
const
&
neumannBoundary
,
Vector
&
f
,
Dune
::
VirtualFunction
<
double
,
double
>
const
&
neumann
,
double
relativeTime
)
{
LocalVector
localNeumann
(
0
);
LocalVector
localNeumann
(
0
);
neumann
.
evaluate
(
relativeTime
,
localNeumann
[
0
]);
neumann
.
evaluate
(
relativeTime
,
localNeumann
[
0
]);
ConstantFunction
<
LocalVector
,
LocalVector
>
const
fNeumann
(
localNeumann
);
ConstantFunction
<
LocalVector
,
LocalVector
>
const
fNeumann
(
localNeumann
);
NeumannBoundaryAssembler
<
typename
GridView
::
Grid
,
LocalVector
>
NeumannBoundaryAssembler
<
Grid
,
LocalVector
>
neumannBoundaryAssembler
(
neumannBoundaryAssembler
(
fNeumann
);
fNeumann
);
assembler
.
assembleBoundaryFunctional
(
neumannBoundaryAssembler
,
f
,
assembler
.
assembleBoundaryFunctional
(
neumannBoundaryAssembler
,
f
,
neumannBoundary
);
neumannBoundary
);
}
}
// Assembles constant 1-function on frictional boundary in nodalIntegrals
template
<
class
GridView
,
int
dimension
>
template
<
class
GridView
,
class
LocalVector
,
class
Assembler
>
auto
MyAssembler
<
GridView
,
dimension
>::
assembleFrictionNonlinearity
(
Dune
::
shared_ptr
<
Dune
::
BlockVector
<
Dune
::
FieldVector
<
double
,
1
>>>
BoundaryPatch
<
GridView
>
const
&
frictionalBoundary
,
FrictionData
const
&
fd
)
assembleFrictionWeightsal
(
GridView
const
&
gridView
,
Assembler
const
&
assembler
,
->
std
::
shared_ptr
<
Dune
::
GlobalNonlinearity
<
Matrix
,
Vector
>>
{
Dune
::
BitSetVector
<
1
>
const
&
frictionalNodes
)
{
// Lump negative normal stress (kludge)
using
Scalar
=
Dune
::
FieldVector
<
double
,
1
>
;
ScalarVector
weights
;
BoundaryPatch
<
GridView
>
const
frictionalBoundary
(
gridView
,
frictionalNodes
);
{
ConstantFunction
<
LocalVector
,
Scalar
>
const
constantOneFunction
(
1
);
ConstantFunction
<
LocalVector
,
typename
ScalarVector
::
block_type
>
const
NeumannBoundaryAssembler
<
typename
GridView
::
Grid
,
Scalar
>
constantOneFunction
(
1
);
frictional
BoundaryAssembler
(
constantOneFunction
);
Neumann
BoundaryAssembler
<
Grid
,
typename
ScalarVector
::
block_type
>
frictionalBoundaryAssembler
(
constantOneFunction
);
auto
const
nodalIntegrals
=
Dune
::
make_shared
<
Dune
::
BlockVector
<
Scalar
>>
();
assembler
.
assembleBoundaryFunctional
(
frictionalBoundaryAssembler
,
weights
,
assembler
.
assembleBoundaryFunctional
(
frictionalBoundary
Assembler
,
frictionalBoundary
);
*
nodalIntegrals
,
frictionalBoundary
);
}
return
nodalIntegrals
;
for
(
size_t
i
=
0
;
i
<
weights
.
size
();
++
i
)
}
assert
(
weights
[
i
]
>=
0.0
);
template
<
class
Matrix
,
class
Vector
>
Dune
::
BitSetVector
<
1
>
frictionalNodes
;
Dune
::
shared_ptr
<
Dune
::
GlobalNonlinearity
<
Matrix
,
Vector
>>
assembleNonlinearity
(
frictionalBoundary
.
getVertices
(
frictionalNodes
);
Dune
::
BitSetVector
<
1
>
const
&
frictionalNodes
,
return
std
::
make_shared
<
Dune
::
GlobalRuinaNonlinearity
<
Matrix
,
Vector
>>
(
Dune
::
BlockVector
<
Dune
::
FieldVector
<
double
,
1
>>
const
&
nodalIntegrals
,
frictionalNodes
,
weights
,
fd
);
FrictionData
const
&
fd
)
{
return
Dune
::
make_shared
<
Dune
::
GlobalRuinaNonlinearity
<
Matrix
,
Vector
>>
(
frictionalNodes
,
nodalIntegrals
,
fd
);
}
}
#include
"assemblers_tmpl.cc"
#include
"assemblers_tmpl.cc"
This diff is collapsed.
Click to expand it.
src/assemblers.hh
+
49
−
19
View file @
3360ff17
...
@@ -3,29 +3,59 @@
...
@@ -3,29 +3,59 @@
#include
<dune/common/bitsetvector.hh>
#include
<dune/common/bitsetvector.hh>
#include
<dune/common/function.hh>
#include
<dune/common/function.hh>
#include
<dune/common/fvector.hh>
#include
<dune/istl/bcrsmatrix.hh>
#include
<dune/common/shared_ptr.hh>
#include
<dune/istl/bvector.hh>
#include
<dune/istl/bvector.hh>
#include
<dune/fufem/assemblers/assembler.hh>
#include
<dune/fufem/assemblers/assembler.hh>
#include
<dune/fufem/functionspacebases/p1nodalbasis.hh>
#include
<dune/tectonic/globalnonlinearity.hh>
#include
<dune/tectonic/globalnonlinearity.hh>
template
<
class
GridView
,
class
LocalVector
,
class
Assembler
>
template
<
class
GridView
,
int
dimension
>
class
MyAssembler
{
void
assembleNeumann
(
GridView
const
&
gridView
,
Assembler
const
&
assembler
,
public:
Dune
::
BitSetVector
<
1
>
const
&
neumannNodes
,
using
ScalarMatrix
=
Dune
::
BCRSMatrix
<
Dune
::
FieldMatrix
<
double
,
1
,
1
>>
;
Dune
::
BlockVector
<
LocalVector
>
&
f
,
using
Matrix
=
Dune
::
VirtualFunction
<
double
,
double
>
const
&
neumann
,
Dune
::
BCRSMatrix
<
Dune
::
FieldMatrix
<
double
,
dimension
,
dimension
>>
;
double
relativeTime
);
using
ScalarVector
=
Dune
::
BlockVector
<
Dune
::
FieldVector
<
double
,
1
>>
;
using
Vector
=
Dune
::
BlockVector
<
Dune
::
FieldVector
<
double
,
dimension
>>
;
template
<
class
GridView
,
class
LocalVector
,
class
Assembler
>
Dune
::
shared_ptr
<
Dune
::
BlockVector
<
Dune
::
FieldVector
<
double
,
1
>>>
using
VertexBasis
=
P1NodalBasis
<
GridView
,
double
>
;
assembleFrictionWeightsal
(
GridView
const
&
gridView
,
Assembler
const
&
assembler
,
VertexBasis
vertexBasis
;
Dune
::
BitSetVector
<
1
>
const
&
frictionalNodes
);
private:
template
<
class
Matrix
,
class
Vector
>
using
Grid
=
typename
GridView
::
Grid
;
Dune
::
shared_ptr
<
Dune
::
GlobalNonlinearity
<
Matrix
,
Vector
>>
assembleNonlinearity
(
using
LocalVector
=
typename
Vector
::
block_type
;
Dune
::
BitSetVector
<
1
>
const
&
frictionalNodes
,
using
LocalVertexBasis
=
typename
VertexBasis
::
LocalFiniteElement
;
Dune
::
BlockVector
<
Dune
::
FieldVector
<
double
,
1
>>
const
&
nodalIntegrals
,
FrictionData
const
&
fd
);
GridView
const
&
gridView
;
Assembler
<
VertexBasis
,
VertexBasis
>
assembler
;
public:
MyAssembler
(
GridView
const
&
gridView
);
void
assembleFrictionalBoundaryMass
(
BoundaryPatch
<
GridView
>
const
&
frictionalBoundary
,
ScalarMatrix
&
frictionalBoundaryMass
);
void
assembleMass
(
double
density
,
Matrix
&
M
);
void
assembleElasticity
(
double
E
,
double
nu
,
Matrix
&
A
);
void
assembleViscosity
(
double
shearViscosity
,
double
bulkViscosity
,
Matrix
&
C
);
void
assembleBodyForce
(
double
gravity
,
double
density
,
LocalVector
zenith
,
Vector
&
f
);
void
assembleNeumann
(
BoundaryPatch
<
GridView
>
const
&
neumannBoundary
,
Vector
&
f
,
Dune
::
VirtualFunction
<
double
,
double
>
const
&
neumann
,
double
relativeTime
);
std
::
shared_ptr
<
Dune
::
GlobalNonlinearity
<
Matrix
,
Vector
>>
assembleFrictionNonlinearity
(
BoundaryPatch
<
GridView
>
const
&
frictionalBoundary
,
FrictionData
const
&
fd
);
};
#endif
#endif
This diff is collapsed.
Click to expand it.
src/assemblers_tmpl.cc
+
1
−
22
View file @
3360ff17
...
@@ -2,27 +2,6 @@
...
@@ -2,27 +2,6 @@
#error DIM unset
#error DIM unset
#endif
#endif
#include
<dune/fufem/functionspacebases/p1nodalbasis.hh>
#include
"explicitgrid.hh"
#include
"explicitgrid.hh"
#include
"explicitvectors.hh"
using
P1Basis
=
P1NodalBasis
<
GridView
,
double
>
;
using
MyAssembler
=
Assembler
<
P1Basis
,
P1Basis
>
;
template
void
assembleNeumann
<
GridView
,
SmallVector
,
MyAssembler
>(
GridView
const
&
gridView
,
MyAssembler
const
&
assembler
,
Dune
::
BitSetVector
<
1
>
const
&
neumannNodes
,
Dune
::
BlockVector
<
SmallVector
>
&
f
,
Dune
::
VirtualFunction
<
double
,
double
>
const
&
neumann
,
double
relativeTime
);
template
Dune
::
shared_ptr
<
Dune
::
BlockVector
<
Dune
::
FieldVector
<
double
,
1
>
>>
assembleFrictionWeightsal
<
GridView
,
SmallVector
,
MyAssembler
>
(
GridView
const
&
gridView
,
MyAssembler
const
&
assembler
,
Dune
::
BitSetVector
<
1
>
const
&
frictionalNodes
);
template
Dune
::
shared_ptr
<
Dune
::
GlobalNonlinearity
<
Matrix
,
Vector
>
>
template
class
MyAssembler
<
GridView
,
DIM
>;
assembleNonlinearity
<
Matrix
,
Vector
>
(
Dune
::
BitSetVector
<
1
>
const
&
frictionalNodes
,
Dune
::
BlockVector
<
Dune
::
FieldVector
<
double
,
1
>>
const
&
nodalIntegrals
,
FrictionData
const
&
fd
);
This diff is collapsed.
Click to expand it.
src/one-body-sample.cc
+
46
−
86
View file @
3360ff17
...
@@ -49,19 +49,12 @@
...
@@ -49,19 +49,12 @@
#include
<dune/grid/utility/structuredgridfactory.hh>
#include
<dune/grid/utility/structuredgridfactory.hh>
#include
<dune/istl/bcrsmatrix.hh>
#include
<dune/istl/bcrsmatrix.hh>
#include
<dune/istl/bvector.hh>
#include
<dune/istl/bvector.hh>
#include
<dune/istl/scaledidmatrix.hh>
#include
<dune/fufem/assemblers/assembler.hh>
#include
<dune/fufem/assemblers/assembler.hh>
#include
<dune/fufem/assemblers/localassemblers/boundarymassassembler.hh>
#include
<dune/fufem/assemblers/localassemblers/l2functionalassembler.hh>
#include
<dune/fufem/assemblers/localassemblers/massassembler.hh>
#include
<dune/fufem/assemblers/localassemblers/stvenantkirchhoffassembler.hh>
#include
<dune/fufem/assemblers/localassemblers/viscosityassembler.hh>
#include
<dune/fufem/assemblers/localassemblers/vonmisesstressassembler.hh>
#include
<dune/fufem/assemblers/localassemblers/vonmisesstressassembler.hh>
#include
<dune/fufem/boundarypatch.hh>
#include
<dune/fufem/boundarypatch.hh>
#include
<dune/fufem/dunepython.hh>
#include
<dune/fufem/dunepython.hh>
#include
<dune/fufem/functions/basisgridfunction.hh>
#include
<dune/fufem/functions/basisgridfunction.hh>
#include
<dune/fufem/functions/constantfunction.hh>
#pragma clang diagnostic push
#pragma clang diagnostic push
#pragma clang diagnostic ignored "-Wsign-compare"
#pragma clang diagnostic ignored "-Wsign-compare"
...
@@ -163,6 +156,8 @@ int main(int argc, char *argv[]) {
...
@@ -163,6 +156,8 @@ int main(int argc, char *argv[]) {
auto
const
E
=
parset
.
get
<
double
>
(
"body.E"
);
auto
const
E
=
parset
.
get
<
double
>
(
"body.E"
);
auto
const
nu
=
parset
.
get
<
double
>
(
"body.nu"
);
auto
const
nu
=
parset
.
get
<
double
>
(
"body.nu"
);
auto
const
density
=
parset
.
get
<
double
>
(
"body.density"
);
double
const
gravity
=
9.81
;
// {{{ Set up grid
// {{{ Set up grid
using
Grid
=
Dune
::
ALUGrid
<
dims
,
dims
,
Dune
::
simplex
,
Dune
::
nonconforming
>
;
using
Grid
=
Dune
::
ALUGrid
<
dims
,
dims
,
Dune
::
simplex
,
Dune
::
nonconforming
>
;
...
@@ -182,6 +177,8 @@ int main(int argc, char *argv[]) {
...
@@ -182,6 +177,8 @@ int main(int argc, char *argv[]) {
grid
=
Dune
::
StructuredGridFactory
<
Grid
>::
createSimplexGrid
(
grid
=
Dune
::
StructuredGridFactory
<
Grid
>::
createSimplexGrid
(
lowerLeft
,
upperRight
,
elements
);
lowerLeft
,
upperRight
,
elements
);
}
}
SmallVector
zenith
(
0
);
zenith
[
1
]
=
1
;
auto
const
refinements
=
parset
.
get
<
size_t
>
(
"grid.refinements"
);
auto
const
refinements
=
parset
.
get
<
size_t
>
(
"grid.refinements"
);
grid
->
globalRefine
(
refinements
);
grid
->
globalRefine
(
refinements
);
...
@@ -193,11 +190,8 @@ int main(int argc, char *argv[]) {
...
@@ -193,11 +190,8 @@ int main(int argc, char *argv[]) {
// Set up bases
// Set up bases
using
P0Basis
=
P0Basis
<
GridView
,
double
>
;
using
P0Basis
=
P0Basis
<
GridView
,
double
>
;
using
P1Basis
=
P1NodalBasis
<
GridView
,
double
>
;
P0Basis
const
p0Basis
(
leafView
);
P0Basis
const
p0Basis
(
leafView
);
P1Basis
const
p1Basis
(
leafView
);
Assembler
<
P0Basis
,
P0Basis
>
const
p0Assembler
(
p0Basis
,
p0Basis
);
Assembler
<
P0Basis
,
P0Basis
>
const
p0Assembler
(
p0Basis
,
p0Basis
);
Assembler
<
P1Basis
,
P1Basis
>
const
p1Assembler
(
p1Basis
,
p1Basis
);
// Set up the boundary
// Set up the boundary
Dune
::
BitSetVector
<
dims
>
velocityDirichletNodes
(
fineVertexCount
,
false
);
Dune
::
BitSetVector
<
dims
>
velocityDirichletNodes
(
fineVertexCount
,
false
);
...
@@ -205,9 +199,10 @@ int main(int argc, char *argv[]) {
...
@@ -205,9 +199,10 @@ int main(int argc, char *argv[]) {
velocityDirichletNodes
;
velocityDirichletNodes
;
Dune
::
BitSetVector
<
dims
>
accelerationDirichletNodes
(
fineVertexCount
,
false
);
Dune
::
BitSetVector
<
dims
>
accelerationDirichletNodes
(
fineVertexCount
,
false
);
Dune
::
BitSetVector
<
1
>
neumannNodes
(
fineVertexCount
,
false
);
Dune
::
BitSetVector
<
1
>
neumannNodes
(
fineVertexCount
,
false
);
Dune
::
BitSetVector
<
1
>
frictionalNodes
(
fineVertexCount
,
false
);
BoundaryPatch
<
GridView
>
const
neumannBoundary
(
leafView
,
neumannNodes
);
Vector
vertexCoordinates
(
fineVertexCount
);
Vector
vertexCoordinates
(
fineVertexCount
);
Dune
::
BitSetVector
<
1
>
frictionalNodes
(
fineVertexCount
,
false
);
{
{
Dune
::
MultipleCodimMultipleGeomTypeMapper
<
Dune
::
MultipleCodimMultipleGeomTypeMapper
<
GridView
,
Dune
::
MCMGVertexLayout
>
const
vertexMapper
(
leafView
);
GridView
,
Dune
::
MCMGVertexLayout
>
const
vertexMapper
(
leafView
);
...
@@ -254,89 +249,54 @@ int main(int argc, char *argv[]) {
...
@@ -254,89 +249,54 @@ int main(int argc, char *argv[]) {
// Set up normal stress, mass matrix, and gravity functional
// Set up normal stress, mass matrix, and gravity functional
double
normalStress
;
double
normalStress
;
Matrix
M
;
EnergyNorm
<
Matrix
,
Vector
>
const
MNorm
(
M
);
Vector
gravityFunctional
;
{
{
double
const
gravity
=
9.81
;
double
volume
=
1.0
;
double
const
density
=
parset
.
get
<
double
>
(
"body.density"
);
for
(
size_t
i
=
0
;
i
<
dims
;
++
i
)
{
volume
*=
(
upperRight
[
i
]
-
lowerLeft
[
i
]);
MassAssembler
<
Grid
,
P1Basis
::
LocalFiniteElement
,
P1Basis
::
LocalFiniteElement
,
double
area
=
1.0
;
Dune
::
ScaledIdentityMatrix
<
double
,
dims
>>
const
localMass
;
for
(
size_t
i
=
0
;
i
<
dims
;
++
i
)
p1Assembler
.
assembleOperator
(
localMass
,
M
);
if
(
i
!=
1
)
M
*=
density
;
area
*=
(
upperRight
[
i
]
-
lowerLeft
[
i
]);
}
{
// volume * gravity * density / area = normal stress
double
volume
=
1.0
;
// V * g * rho / A = sigma_n
for
(
size_t
i
=
0
;
i
<
dims
;
++
i
)
// m^d * N/kg * kg/m^d / m^(d-1) = N/m^(d-1)
volume
*=
(
upperRight
[
i
]
-
lowerLeft
[
i
]);
normalStress
=
volume
*
gravity
*
density
/
area
;
double
area
=
1.0
;
for
(
size_t
i
=
0
;
i
<
dims
;
++
i
)
if
(
i
!=
1
)
area
*=
(
upperRight
[
i
]
-
lowerLeft
[
i
]);
// volume * gravity * density / area = normal stress
// V * g * rho / A = sigma_n
// m^d * N/kg * kg/m^d / m^(d-1) = N/m^(d-1)
normalStress
=
volume
*
gravity
*
density
/
area
;
}
{
SmallVector
weightedGravitationalDirection
(
0
);
weightedGravitationalDirection
[
1
]
=
-
density
*
gravity
;
ConstantFunction
<
SmallVector
,
SmallVector
>
const
gravityFunction
(
weightedGravitationalDirection
);
L2FunctionalAssembler
<
Grid
,
P1Basis
::
LocalFiniteElement
,
SmallVector
>
gravityFunctionalAssembler
(
gravityFunction
);
p1Assembler
.
assembleFunctional
(
gravityFunctionalAssembler
,
gravityFunctional
);
}
}
}
FrictionData
const
frictionData
(
parset
.
sub
(
"boundary.friction"
),
FrictionData
const
frictionData
(
parset
.
sub
(
"boundary.friction"
),
normalStress
);
normalStress
);
Matrix
A
;
using
MyAssembler
=
MyAssembler
<
GridView
,
dims
>
;
EnergyNorm
<
Matrix
,
Vector
>
const
ANorm
(
A
);
{
StVenantKirchhoffAssembler
<
Grid
,
P1Basis
::
LocalFiniteElement
,
P1Basis
::
LocalFiniteElement
>
const
localStiffness
(
E
,
nu
);
p1Assembler
.
assembleOperator
(
localStiffness
,
A
);
}
Matrix
C
;
{
ViscosityAssembler
<
Grid
,
P1Basis
::
LocalFiniteElement
,
P1Basis
::
LocalFiniteElement
>
const
localViscosity
(
parset
.
get
<
double
>
(
"body.shearViscosity"
),
parset
.
get
<
double
>
(
"body.bulkViscosity"
));
p1Assembler
.
assembleOperator
(
localViscosity
,
C
);
}
ScalarMatrix
frictionalBoundaryMassMatrix
;
MyAssembler
myAssembler
(
leafView
);
EnergyNorm
<
ScalarMatrix
,
ScalarVector
>
const
stateEnergyNorm
(
frictionalBoundaryMassMatrix
);
Matrix
A
,
C
,
M
;
{
myAssembler
.
assembleElasticity
(
E
,
nu
,
A
);
BoundaryMassAssembler
<
myAssembler
.
assembleViscosity
(
parset
.
get
<
double
>
(
"body.shearViscosity"
),
Grid
,
BoundaryPatch
<
GridView
>
,
P1Basis
::
LocalFiniteElement
,
parset
.
get
<
double
>
(
"body.bulkViscosity"
),
C
);
P1Basis
::
LocalFiniteElement
,
SmallScalarMatrix
>
const
myAssembler
.
assembleMass
(
density
,
M
);
frictionalBoundaryMassAssembler
(
frictionalBoundary
);
EnergyNorm
<
Matrix
,
Vector
>
const
ANorm
(
A
),
MNorm
(
M
);
p1Assembler
.
assembleOperator
(
frictionalBoundaryMassAssembler
,
frictionalBoundaryMassMatrix
);
}
// Q: Does it make sense to weigh them in this manner?
// Q: Does it make sense to weigh them in this manner?
SumNorm
<
Vector
>
const
AMNorm
(
1.0
,
ANorm
,
1.0
,
MNorm
);
SumNorm
<
Vector
>
const
AMNorm
(
1.0
,
ANorm
,
1.0
,
MNorm
);
auto
const
nodalIntegrals
=
ScalarMatrix
frictionalBoundaryMass
;
assembleFrictionWeightsal
<
GridView
,
SmallVector
>
(
leafView
,
p1Assembler
,
myAssembler
.
assembleFrictionalBoundaryMass
(
frictionalBoundary
,
frictionalNodes
);
frictionalBoundaryMass
);
auto
myGlobalNonlinearity
=
assembleNonlinearity
<
Matrix
,
Vector
>
(
EnergyNorm
<
ScalarMatrix
,
ScalarVector
>
const
stateEnergyNorm
(
frictionalNodes
,
*
nodalIntegrals
,
frictionData
);
frictionalBoundaryMass
);
// Assemble forces
Vector
gravityFunctional
;
myAssembler
.
assembleBodyForce
(
gravity
,
density
,
zenith
,
gravityFunctional
);
auto
myGlobalNonlinearity
=
myAssembler
.
assembleFrictionNonlinearity
(
frictionalBoundary
,
frictionData
);
// Problem formulation: right-hand side
// Problem formulation: right-hand side
auto
const
createRHS
=
[
&
](
double
_relativeTime
,
Vector
&
_ell
)
{
auto
const
createRHS
=
[
&
](
double
_relativeTime
,
Vector
&
_ell
)
{
assembleNeumann
(
leafView
,
p1Assembler
,
neumannNodes
,
_ell
,
myAssembler
.
assembleNeumann
(
neumannBoundary
,
_ell
,
neumannFunction
,
neumannFunction
,
_relativeTime
);
_relativeTime
);
_ell
+=
gravityFunctional
;
_ell
+=
gravityFunctional
;
};
};
Vector
ell
(
fineVertexCount
);
Vector
ell
(
fineVertexCount
);
...
@@ -569,15 +529,15 @@ int main(int argc, char *argv[]) {
...
@@ -569,15 +529,15 @@ int main(int argc, char *argv[]) {
relaxationWriter
<<
std
::
endl
;
relaxationWriter
<<
std
::
endl
;
if
(
parset
.
get
<
bool
>
(
"io.writeVTK"
))
{
if
(
parset
.
get
<
bool
>
(
"io.writeVTK"
))
{
auto
const
gridDisplacement
=
auto
const
gridDisplacement
=
Dune
::
make_shared
<
Dune
::
make_shared
<
BasisGridFunction
<
P1
Basis
,
Vector
>
const
>
(
p1Basis
,
BasisGridFunction
<
typename
MyAssembler
::
Vertex
Basis
,
Vector
>
const
>
(
u
);
myAssembler
.
vertexBasis
,
u
);
ScalarVector
vonMisesStress
;
ScalarVector
vonMisesStress
;
VonMisesStressAssembler
<
Grid
,
P0Basis
::
LocalFiniteElement
>
VonMisesStressAssembler
<
Grid
,
P0Basis
::
LocalFiniteElement
>
localStressAssembler
(
E
,
nu
,
gridDisplacement
);
localStressAssembler
(
E
,
nu
,
gridDisplacement
);
p0Assembler
.
assembleFunctional
(
localStressAssembler
,
vonMisesStress
);
p0Assembler
.
assembleFunctional
(
localStressAssembler
,
vonMisesStress
);
writeVtk
(
p1
Basis
,
u
,
alpha
,
p0Basis
,
vonMisesStress
,
writeVtk
(
myAssembler
.
vertex
Basis
,
u
,
alpha
,
p0Basis
,
vonMisesStress
,
(
boost
::
format
(
"obs%d"
)
%
run
).
str
());
(
boost
::
format
(
"obs%d"
)
%
run
).
str
());
}
}
}
}
...
...
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