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
GitLab community forum
Contribute to GitLab
Provide feedback
Keyboard shortcuts
?
Snippets
Groups
Projects
This is an archived project. Repository and other project resources are read-only.
Show more breadcrumbs
Ansgar Burchardt
dune-elasticity
Commits
a40e365e
Commit
a40e365e
authored
Oct 7, 2017
by
Jonathan Youett
Browse files
Options
Downloads
Patches
Plain Diff
Adjust to upstream changes
parent
c63d18ed
No related branches found
No related tags found
No related merge requests found
Changes
2
Show whitespace changes
Inline
Side-by-side
Showing
2 changed files
linelast.cc
+17
-16
17 additions, 16 deletions
linelast.cc
viscoelast.cc
+15
-16
15 additions, 16 deletions
viscoelast.cc
with
32 additions
and
32 deletions
linelast.cc
+
17
−
16
View file @
a40e365e
...
@@ -27,7 +27,7 @@
...
@@ -27,7 +27,7 @@
#ifdef HAVE_IPOPT
#ifdef HAVE_IPOPT
#include
<dune/solvers/solvers/quadraticipopt.hh>
#include
<dune/solvers/solvers/quadraticipopt.hh>
#endif
#endif
#include
<dune/solvers/iterationsteps/blockgsstep.hh>
#include
<dune/solvers/iterationsteps/blockgsstep
s
.hh>
#include
<dune/solvers/iterationsteps/multigridstep.hh>
#include
<dune/solvers/iterationsteps/multigridstep.hh>
#include
<dune/solvers/transferoperators/compressedmultigridtransfer.hh>
#include
<dune/solvers/transferoperators/compressedmultigridtransfer.hh>
#include
<dune/solvers/solvers/loopsolver.hh>
#include
<dune/solvers/solvers/loopsolver.hh>
...
@@ -182,48 +182,49 @@ int main (int argc, char *argv[]) try
...
@@ -182,48 +182,49 @@ int main (int argc, char *argv[]) try
#if !defined HAVE_IPOPT
#if !defined HAVE_IPOPT
#error You need to have IPOpt installed if you want to use it as the base solver!
#error You need to have IPOpt installed if you want to use it as the base solver!
#endif
#endif
QuadraticIPOptSolver
<
OperatorType
,
VectorType
>
baseSolver
;
QuadraticIPOptSolver
<
OperatorType
,
VectorType
>
baseSolver
(
baseTolerance
,
100
,
Solver
::
QUIET
);
baseSolver
.
tolerance_
=
baseTolerance
;
baseSolver
.
verbosity_
=
Solver
::
QUIET
;
#else // Gauss-Seidel is the base solver
#else // Gauss-Seidel is the base solver
BlockGSStep
<
OperatorType
,
VectorType
>
baseSolverStep
;
auto
baseSolverStep
=
Dune
::
Solvers
::
BlockGSStepFactory
<
OperatorType
,
VectorType
>::
create
(
Dune
::
Solvers
::
BlockGS
::
LocalSolvers
::
direct
(
0.0
));
EnergyNorm
<
OperatorType
,
VectorType
>
baseEnergyNorm
(
baseSolverStep
);
EnergyNorm
<
OperatorType
,
VectorType
>
baseEnergyNorm
(
baseSolverStep
);
LoopSolver
<
VectorType
>
baseSolver
(
&
baseSolverStep
,
LoopSolver
<
VectorType
>
baseSolver
(
baseSolverStep
,
baseIt
,
baseIt
,
baseTolerance
,
baseTolerance
,
&
baseEnergyNorm
,
baseEnergyNorm
,
Solver
::
QUIET
);
Solver
::
QUIET
);
#endif
#endif
// Make pre and postsmoothers
// Make pre and postsmoothers
BlockGSStep
<
OperatorType
,
VectorType
>
presmoother
;
auto
presmoother
=
Dune
::
Solvers
::
BlockGSStepFactory
<
OperatorType
,
VectorType
>::
create
(
BlockGSStep
<
OperatorType
,
VectorType
>
postsmoother
;
Dune
::
Solvers
::
BlockGS
::
LocalSolvers
::
direct
(
0.0
));
auto
postsmoother
=
Dune
::
Solvers
::
BlockGSStepFactory
<
OperatorType
,
VectorType
>::
create
(
Dune
::
Solvers
::
BlockGS
::
LocalSolvers
::
direct
(
0.0
));
MultigridStep
<
OperatorType
,
VectorType
>
multigridStep
(
stiffnessMatrix
,
x
,
rhs
);
MultigridStep
<
OperatorType
,
VectorType
>
multigridStep
(
stiffnessMatrix
,
x
,
rhs
);
multigridStep
.
setMGType
(
mu
,
nu1
,
nu2
);
multigridStep
.
setMGType
(
mu
,
nu1
,
nu2
);
multigridStep
.
i
gnore
Nodes_
=
&
dirichletNodes
;
multigridStep
.
setI
gnore
(
dirichletNodes
)
;
multigridStep
.
b
ase
s
olver
_
=
&
baseSolver
;
multigridStep
.
setB
ase
S
olver
(
baseSolver
)
;
multigridStep
.
setSmoother
(
&
presmoother
,
&
postsmoother
);
multigridStep
.
setSmoother
(
&
presmoother
,
&
postsmoother
);
std
::
vector
<
CompressedMultigridTransfer
<
VectorType
>
*
>
mgTransfers
(
grid
->
maxLevel
());
std
::
vector
<
std
::
shared_ptr
<
CompressedMultigridTransfer
<
VectorType
>
>
>
mgTransfers
(
grid
->
maxLevel
());
for
(
size_t
i
=
0
;
i
<
mgTransfers
.
size
();
i
++
)
{
for
(
size_t
i
=
0
;
i
<
mgTransfers
.
size
();
i
++
)
{
mgTransfers
[
i
]
=
new
CompressedMultigridTransfer
<
VectorType
>
;
mgTransfers
[
i
]
=
std
::
make_shared
<
CompressedMultigridTransfer
<
VectorType
>
>
()
;
mgTransfers
[
i
]
->
setup
(
*
grid
,
i
,
i
+
1
);
mgTransfers
[
i
]
->
setup
(
*
grid
,
i
,
i
+
1
);
}
}
multigridStep
.
setTransferOperators
(
mgTransfers
);
multigridStep
.
setTransferOperators
(
mgTransfers
);
EnergyNorm
<
OperatorType
,
VectorType
>
energyNorm
(
multigridStep
);
EnergyNorm
<
OperatorType
,
VectorType
>
energyNorm
(
multigridStep
);
LoopSolver
<
VectorType
>
solver
(
&
multigridStep
,
LoopSolver
<
VectorType
>
solver
(
multigridStep
,
numIt
,
numIt
,
tolerance
,
tolerance
,
&
energyNorm
,
energyNorm
,
Solver
::
FULL
);
Solver
::
FULL
);
...
...
This diff is collapsed.
Click to expand it.
viscoelast.cc
+
15
−
16
View file @
a40e365e
...
@@ -10,7 +10,7 @@
...
@@ -10,7 +10,7 @@
#include
<dune/istl/io.hh>
#include
<dune/istl/io.hh>
#include
<dune/solvers/iterationsteps/blockgsstep.hh>
#include
<dune/solvers/iterationsteps/blockgsstep
s
.hh>
#include
<dune/solvers/iterationsteps/multigridstep.hh>
#include
<dune/solvers/iterationsteps/multigridstep.hh>
#include
<dune/solvers/transferoperators/compressedmultigridtransfer.hh>
#include
<dune/solvers/transferoperators/compressedmultigridtransfer.hh>
#include
<dune/solvers/norms/energynorm.hh>
#include
<dune/solvers/norms/energynorm.hh>
...
@@ -177,24 +177,23 @@ int main (int argc, char *argv[]) try {
...
@@ -177,24 +177,23 @@ int main (int argc, char *argv[]) try {
// ///////////////////////////
// ///////////////////////////
// First create a base solver
// First create a base solver
BlockGSStep
<
MatrixType
,
VectorType
>
baseSolverStep
;
auto
baseSolverStep
=
Dune
::
Solvers
::
BlockGSStepFactory
<
MatrixType
,
VectorType
>::
create
(
Dune
::
Solvers
::
BlockGS
::
LocalSolvers
::
direct
(
0.0
));
EnergyNorm
<
MatrixType
,
VectorType
>
baseEnergyNorm
(
baseSolverStep
);
EnergyNorm
<
MatrixType
,
VectorType
>
baseEnergyNorm
(
baseSolverStep
);
::
LoopSolver
<
VectorType
>
baseSolver
(
&
baseSolverStep
,
::
LoopSolver
<
VectorType
>
baseSolver
(
baseSolverStep
,
baseIt
,
baseIt
,
baseTolerance
,
baseTolerance
,
&
baseEnergyNorm
,
baseEnergyNorm
,
Solver
::
QUIET
);
Solver
::
QUIET
);
baseSolver
.
verbosity_
=
Solver
::
QUIET
;
baseSolver
.
tolerance_
=
baseTolerance
;
// Make pre and postsmoothers
// Make pre and postsmoothers
BlockGSStep
<
MatrixType
,
VectorType
>
presmoother
;
BlockGSStep
<
MatrixType
,
VectorType
>
postsmoother
;
auto
presmoother
=
Dune
::
Solvers
::
BlockGSStepFactory
<
MatrixType
,
VectorType
>::
create
(
Dune
::
Solvers
::
BlockGS
::
LocalSolvers
::
direct
(
0.0
));
auto
postsmoother
=
Dune
::
Solvers
::
BlockGSStepFactory
<
MatrixType
,
VectorType
>::
create
(
Dune
::
Solvers
::
BlockGS
::
LocalSolvers
::
direct
(
0.0
));
//FEM leads to:
//FEM leads to:
//Quasistatic equation of the form : N*D'(t)+K*D(t)=rhs(t) N viscous_stiffness K elastic_stiffness D displacement
//Quasistatic equation of the form : N*D'(t)+K*D(t)=rhs(t) N viscous_stiffness K elastic_stiffness D displacement
...
@@ -215,23 +214,23 @@ int main (int argc, char *argv[]) try {
...
@@ -215,23 +214,23 @@ int main (int argc, char *argv[]) try {
MultigridStep
<
MatrixType
,
VectorType
>
multigridStep
(
lgs
,
x
,
erhs
);
MultigridStep
<
MatrixType
,
VectorType
>
multigridStep
(
lgs
,
x
,
erhs
);
multigridStep
.
setMGType
(
mu
,
nu1
,
nu2
);
multigridStep
.
setMGType
(
mu
,
nu1
,
nu2
);
multigridStep
.
i
gnore
Nodes_
=
&
dirichletNodes
;
multigridStep
.
setI
gnore
(
dirichletNodes
)
;
multigridStep
.
b
ase
s
olver
_
=
&
baseSolver
;
multigridStep
.
setB
ase
S
olver
(
baseSolver
)
;
multigridStep
.
setSmoother
(
&
presmoother
,
&
postsmoother
);
multigridStep
.
setSmoother
(
&
presmoother
,
&
postsmoother
);
std
::
vector
<
CompressedMultigridTransfer
<
VectorType
>
*
>
transfers
(
grid
->
maxLevel
());
std
::
vector
<
std
::
shared_ptr
<
CompressedMultigridTransfer
<
VectorType
>
>
>
transfers
(
grid
->
maxLevel
());
for
(
size_t
i
=
0
;
i
<
transfers
.
size
();
i
++
)
{
for
(
size_t
i
=
0
;
i
<
transfers
.
size
();
i
++
)
{
transfers
[
i
]
=
new
CompressedMultigridTransfer
<
VectorType
>
;
transfers
[
i
]
=
std
::
make_shared
<
CompressedMultigridTransfer
<
VectorType
>
>
()
;
transfers
[
i
]
->
setup
(
*
grid
,
i
,
i
+
1
);
transfers
[
i
]
->
setup
(
*
grid
,
i
,
i
+
1
);
}
}
multigridStep
.
setTransferOperators
(
transfers
);
multigridStep
.
setTransferOperators
(
transfers
);
EnergyNorm
<
MatrixType
,
VectorType
>
energyNorm
(
multigridStep
);
EnergyNorm
<
MatrixType
,
VectorType
>
energyNorm
(
multigridStep
);
::
LoopSolver
<
VectorType
>
solver
(
&
multigridStep
,
::
LoopSolver
<
VectorType
>
solver
(
multigridStep
,
numIt
,
numIt
,
tolerance
,
tolerance
,
&
energyNorm
,
energyNorm
,
Solver
::
FULL
);
Solver
::
FULL
);
solver
.
preprocess
();
solver
.
preprocess
();
...
...
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