Skip to content
GitLab
Explore
Sign in
Primary navigation
Search or go to…
Project
D
dune-istl
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
Model registry
Operate
Environments
Monitor
Incidents
Service Desk
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
Core Modules
dune-istl
Commits
3abefe26
Commit
3abefe26
authored
18 years ago
by
Markus Blatt
Browse files
Options
Downloads
Patches
Plain Diff
Check whether diagonals are present.
[[Imported from SVN: r695]]
parent
68895878
Branches
Branches containing commit
Tags
Tags containing commit
No related merge requests found
Changes
2
Hide whitespace changes
Inline
Side-by-side
Showing
2 changed files
istl/matrixutils.hh
+36
-0
36 additions, 0 deletions
istl/matrixutils.hh
istl/preconditioners.hh
+13
-4
13 additions, 4 deletions
istl/preconditioners.hh
with
49 additions
and
4 deletions
istl/matrixutils.hh
+
36
−
0
View file @
3abefe26
...
...
@@ -5,6 +5,7 @@
#include
<dune/common/typetraits.hh>
#include
<dune/common/helpertemplates.hh>
#include
"istlexception.hh"
namespace
Dune
{
...
...
@@ -54,6 +55,41 @@ namespace Dune
}
/**
* @brief Check whether the a matrix has diagonal values
* on blocklevel recursion levels.
*/
template
<
std
::
size_t
blocklevel
,
std
::
size_t
l
=
blocklevel
>
struct
CheckIfDiagonalPresent
{
/**
* @brief Check whether the a matrix has diagonal values
* on blocklevel recursion levels.
*/
template
<
class
Matrix
>
static
void
check
(
const
Matrix
&
mat
)
{
#ifdef DUNE_ISTL_WITH_CHECKING
CheckIfDiagonalPresent
<
blocklevel
-
1
,
l
>::
check
(
mat
);
#endif
}
};
template
<
std
::
size_t
l
>
struct
CheckIfDiagonalPresent
<
0
,
l
>
{
template
<
class
Matrix
>
static
void
check
(
const
Matrix
&
mat
)
{
typedef
typename
Matrix
::
ConstRowIterator
Row
;
for
(
Row
row
=
mat
.
begin
();
row
!=
mat
.
end
();
++
row
)
{
if
(
row
->
find
(
row
.
index
())
==
row
->
end
())
DUNE_THROW
(
ISTLException
,
"Missing diagonal value in row "
<<
row
.
index
()
<<
" at block recursion level "
<<
l
);
}
}
};
/**
* @brief Get the number of nonzero fields in the matrix.
*
...
...
This diff is collapsed.
Click to expand it.
istl/preconditioners.hh
+
13
−
4
View file @
3abefe26
...
...
@@ -11,6 +11,7 @@
#include
"solvercategory.hh"
#include
"istlexception.hh"
#include
"matrixutils.hh"
#include
"io.hh"
#include
"gsetc.hh"
#include
"ilu.hh"
...
...
@@ -134,7 +135,9 @@ namespace Dune {
*/
SeqSSOR
(
const
M
&
A
,
int
n
,
field_type
w
)
:
_A_
(
A
),
_n
(
n
),
_w
(
w
)
{
}
{
CheckIfDiagonalPresent
<
l
>::
check
(
_A_
);
}
/*!
\brief Prepare the preconditioner.
...
...
@@ -207,7 +210,9 @@ namespace Dune {
*/
SeqSOR
(
const
M
&
A
,
int
n
,
field_type
w
)
:
_A_
(
A
),
_n
(
n
),
_w
(
w
)
{
}
{
CheckIfDiagonalPresent
<
l
>::
check
(
_A_
);
}
/*!
\brief Prepare the preconditioner.
...
...
@@ -277,7 +282,9 @@ namespace Dune {
*/
SeqGS
(
const
M
&
A
,
int
n
,
field_type
w
)
:
_A_
(
A
),
_n
(
n
),
_w
(
w
)
{
}
{
CheckIfDiagonalPresent
<
l
>::
check
(
_A_
);
}
/*!
\brief Prepare the preconditioner.
...
...
@@ -347,7 +354,9 @@ namespace Dune {
*/
SeqJac
(
const
M
&
A
,
int
n
,
field_type
w
)
:
_A_
(
A
),
_n
(
n
),
_w
(
w
)
{
}
{
CheckIfDiagonalPresent
<
l
>::
check
(
_A_
);
}
/*!
\brief Prepare the preconditioner.
...
...
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