Skip to content
GitLab
Explore
Sign in
Primary navigation
Search or go to…
Project
D
dune-solvers
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
Build
Pipelines
Jobs
Pipeline schedules
Test cases
Artifacts
Deploy
Releases
Model registry
Operate
Environments
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
agnumpde
dune-solvers
Commits
81a90c95
Commit
81a90c95
authored
16 years ago
by
Oliver Sander
Committed by
sander@PCPOOL.MI.FU-BERLIN.DE
16 years ago
Browse files
Options
Downloads
Patches
Plain Diff
moved to new module dune-solvers
[[Imported from SVN: r2345]]
parent
c408784a
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
dune-solvers/iterationsteps/projectedblockgsstep.cc
+87
-0
87 additions, 0 deletions
dune-solvers/iterationsteps/projectedblockgsstep.cc
dune-solvers/iterationsteps/projectedblockgsstep.hh
+41
-0
41 additions, 0 deletions
dune-solvers/iterationsteps/projectedblockgsstep.hh
with
128 additions
and
0 deletions
dune-solvers/iterationsteps/projectedblockgsstep.cc
0 → 100644
+
87
−
0
View file @
81a90c95
template
<
class
OperatorType
,
class
DiscFuncType
>
inline
DiscFuncType
ProjectedBlockGSStep
<
OperatorType
,
DiscFuncType
>::
getSol
()
{
return
*
(
this
->
x_
);
}
template
<
class
OperatorType
,
class
DiscFuncType
>
inline
void
ProjectedBlockGSStep
<
OperatorType
,
DiscFuncType
>::
iterate
()
{
if
(
hasObstacle_
->
size
()
!=
(
unsigned
int
)
this
->
x_
->
size
())
DUNE_THROW
(
SolverError
,
"Size of hasObstacle ("
<<
hasObstacle_
->
size
()
<<
") doesn't match solution vector ("
<<
this
->
x_
->
size
()
<<
")"
);
const
OperatorType
&
mat
=
*
this
->
mat_
;
for
(
int
i
=
0
;
i
<
this
->
x_
->
size
();
i
++
)
{
/** \todo Handle more general boundary conditions */
if
((
*
this
->
ignoreNodes_
)[
i
][
0
])
continue
;
bool
zeroDiagonal
=
false
;
for
(
int
j
=
0
;
j
<
BlockSize
;
j
++
)
{
// When using this solver as part of a truncated multigrid solver,
// the diagonal entries of the matrix may get completely truncated
// away. In this case we just do nothing here.
zeroDiagonal
|=
(
mat
[
i
][
i
][
j
][
j
]
<
1e-10
);
}
if
(
zeroDiagonal
)
continue
;
VectorBlock
r
;
residual
(
i
,
r
);
// Compute x_i += A_{i,i}^{-1} r[i]
VectorBlock
v
;
VectorBlock
&
x
=
(
*
this
->
x_
)[
i
];
if
((
*
hasObstacle_
)[
i
]
==
false
)
{
// No obstacle --> solve linear problem
mat
[
i
][
i
].
solve
(
v
,
r
);
}
else
{
// Solve the local constraint minimization problem
// We use a projected Gauss-Seidel, for lack of anything better
BoxConstraint
<
field_type
,
BlockSize
>
defectObstacle
=
(
*
obstacles_
)[
i
];
defectObstacle
-=
x
;
// Initial correction
v
=
0
;
for
(
int
j
=
0
;
j
<
20
;
j
++
)
{
for
(
int
k
=
0
;
k
<
BlockSize
;
k
++
)
{
// Compute residual
field_type
sr
=
0
;
for
(
int
l
=
0
;
l
<
BlockSize
;
l
++
)
sr
+=
mat
[
i
][
i
][
k
][
l
]
*
v
[
l
];
sr
=
r
[
k
]
-
sr
;
v
[
k
]
+=
sr
/
mat
[
i
][
i
][
k
][
k
];
// Project
if
(
v
[
k
]
<
defectObstacle
.
lower
(
k
))
v
[
k
]
=
defectObstacle
.
lower
(
k
);
else
if
(
v
[
k
]
>
defectObstacle
.
upper
(
k
))
v
[
k
]
=
defectObstacle
.
upper
(
k
);
}
}
}
// Add correction
x
+=
v
;
}
}
This diff is collapsed.
Click to expand it.
dune-solvers/iterationsteps/projectedblockgsstep.hh
0 → 100644
+
41
−
0
View file @
81a90c95
#ifndef DUNE_PROJECTED_BLOCK_GAUSS_SEIDEL_STEP_HH
#define DUNE_PROJECTED_BLOCK_GAUSS_SEIDEL_STEP_HH
#include
<vector>
#include
"blockgsstep.hh"
#include
"boxconstraint.hh"
template
<
class
OperatorType
,
class
DiscFuncType
>
class
ProjectedBlockGSStep
:
public
BlockGSStep
<
OperatorType
,
DiscFuncType
>
{
typedef
typename
DiscFuncType
::
block_type
VectorBlock
;
typedef
typename
DiscFuncType
::
field_type
field_type
;
enum
{
BlockSize
=
VectorBlock
::
dimension
};
public
:
//! Default constructor. Doesn't init anything
ProjectedBlockGSStep
()
{}
//! Constructor with a linear problem
ProjectedBlockGSStep
(
OperatorType
&
mat
,
DiscFuncType
&
x
,
DiscFuncType
&
rhs
)
:
LinearIterationStep
<
OperatorType
,
DiscFuncType
>
(
mat
,
x
,
rhs
)
{}
//! Perform one iteration
virtual
void
iterate
();
virtual
DiscFuncType
getSol
();
Dune
::
BitSetVector
<
1
>*
hasObstacle_
;
const
std
::
vector
<
BoxConstraint
<
field_type
,
BlockSize
>
>*
obstacles_
;
};
#include
"projectedblockgsstep.cc"
#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