Skip to content
GitLab
Explore
Sign in
Register
Primary navigation
Search or go to…
Project
dumux
Manage
Activity
Members
Labels
Plan
Issues
Issue boards
Milestones
Wiki
Code
Merge requests
Repository
Branches
Commits
Tags
Repository graph
Compare revisions
Build
Pipelines
Jobs
Pipeline schedules
Artifacts
Deploy
Releases
Container Registry
Model registry
Operate
Environments
Monitor
Incidents
Analyze
Value stream analytics
Contributor analytics
CI/CD analytics
Repository analytics
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
dumux-repositories
dumux
Commits
f0e9dcbc
"materialsystem.html" did not exist on "ade8b17b3d9398d9bf2e65dec44a1d5b61c0d26d"
Commit
f0e9dcbc
authored
4 years ago
by
Kilian Weishaupt
Browse files
Options
Downloads
Patches
Plain Diff
[linear][uzawa] Improve docu
parent
47e1b168
No related branches found
No related tags found
1 merge request
!1921
Feature/uzawa 3.2
Changes
2
Hide whitespace changes
Inline
Side-by-side
Showing
2 changed files
doc/handbook/dumux-handbook.bib
+25
-0
25 additions, 0 deletions
doc/handbook/dumux-handbook.bib
dumux/linear/preconditioners.hh
+53
-0
53 additions, 0 deletions
dumux/linear/preconditioners.hh
with
78 additions
and
0 deletions
doc/handbook/dumux-handbook.bib
+
25
−
0
View file @
f0e9dcbc
...
@@ -1882,3 +1882,28 @@ year={1999}
...
@@ -1882,3 +1882,28 @@ year={1999}
pages
=
{37}
,
pages
=
{37}
,
year
=
{1964}
year
=
{1964}
}
}
@article
{
benzi2005
,
doi
=
{10.1017/s0962492904000212}
,
url
=
{https://doi.org/10.1017/s0962492904000212}
,
year
=
{2005}
,
publisher
=
{Cambridge University Press ({CUP})}
,
volume
=
{14}
,
pages
=
{1--137}
,
author
=
{Michele Benzi and Gene H. Golub and J\"{o}rg Liesen}
,
title
=
{Numerical solution of saddle point problems}
,
journal
=
{Acta Numerica}
}
@article
{
ho2017
,
doi
=
{10.1137/16m1076770}
,
url
=
{https://doi.org/10.1137/16m1076770}
,
year
=
{2017}
,
publisher
=
{Society for Industrial {\&} Applied Mathematics ({SIAM})}
,
volume
=
{39}
,
number
=
{5}
,
pages
=
{S461--S476}
,
author
=
{Nguyenho Ho and Sarah D. Olson and Homer F. Walker}
,
title
=
{Accelerating the Uzawa Algorithm}
,
journal
=
{{SIAM} Journal on Scientific Computing}
}
This diff is collapsed.
Click to expand it.
dumux/linear/preconditioners.hh
+
53
−
0
View file @
f0e9dcbc
...
@@ -44,6 +44,40 @@
...
@@ -44,6 +44,40 @@
namespace
Dumux
{
namespace
Dumux
{
#if DUNE_VERSION_GTE(DUNE_ISTL,2,7)
#if DUNE_VERSION_GTE(DUNE_ISTL,2,7)
/*!
* \ingroup Linear
* \brief A preconditioner based on the Uzawa algorithm for saddle-point problems of the form
* \f$
\begin{pmatrix}
A & B \\
C & D
\end{pmatrix}
\begin{pmatrix}
u\\
p
\end{pmatrix}
=
\begin{pmatrix}
f\\
g
\end{pmatrix}
* \f$
*
* This preconditioner is especially suited for solving the incompressible (Navier-)Stokes equations.
* Here, \f$D = 0\f$ and \f$B = C^T\f$ if \f$\rho = 1\f$.
* We do not expect good convergence if energy or mass transport is considered.
*
* See: Benzi, M., Golub, G. H., & Liesen, J. (2005). Numerical solution of saddle point problems. Acta numerica, 14, 1-137 \cite benzi2005 and <BR>
* Ho, N., Olson, S. D., & Walker, H. F. (2017). Accelerating the Uzawa algorithm. SIAM Journal on Scientific Computing, 39(5), S461-S476 \cite ho2017
*
* \tparam M Type of the matrix.
* \tparam X Type of the update.
* \tparam Y Type of the defect.
* \tparam l Preconditioner block level (for compatibility reasons, unused).
*/
template
<
class
M
,
class
X
,
class
Y
,
int
l
=
1
>
template
<
class
M
,
class
X
,
class
Y
,
int
l
=
1
>
class
SeqUzawa
:
public
Dune
::
Preconditioner
<
X
,
Y
>
class
SeqUzawa
:
public
Dune
::
Preconditioner
<
X
,
Y
>
{
{
...
@@ -189,6 +223,25 @@ private:
...
@@ -189,6 +223,25 @@ private:
#endif
#endif
}
}
/*!
* \brief Estimate the relaxation factor omega
*
* The optimal relaxation factor is omega = 2/(lambdaMin + lambdaMax), where lambdaMin and lambdaMax are the smallest and largest
* eigenvalues of the Schur complement -C*Ainv*B (assuming D = 0).
* lambdaMax can be easily determined using the power iteration algorithm (https://en.wikipedia.org/wiki/Power_iteration) and lambdaMin
* could be estimated in a similar manner. We do not consider lambdaMin because for certain cases, e.g., when C contains some rows of zeroes only,
* this estimate will fail.
*
* Instead we assume that lambdaMin is sufficiently close to lambdaMax such that omega = 1/lambdaMax.
* This seems to work rather well for various applications.
* We will underestimate omega by a factor of 2 in the worst case (i.e, lambdaMin = 0).
*
* When facing convergence issues, you may set LinearSolver.Preconditioner.Verbosity = 1 to see the estimate of lambdaMax.
* In a new simulation run, you can then set LinearSolver.Preconditioner.DetermineRelaxationFactor = false and set some other value
* for LinearSolver.Preconditioner.Relaxation based on the estimate of lambdaMax.
*
* See: Benzi, M., Golub, G. H., & Liesen, J. (2005). Numerical solution of saddle point problems. Acta numerica, 14, 1-137.
*/
scalar_field_type
estimateOmega_
()
const
scalar_field_type
estimateOmega_
()
const
{
{
using
namespace
Dune
::
Indices
;
using
namespace
Dune
::
Indices
;
...
...
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