Skip to content
GitLab
Explore
Sign in
Primary navigation
Search or go to…
Project
D
dune-gfe
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
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
Sander, Oliver
dune-gfe
Commits
408c0310
Commit
408c0310
authored
11 years ago
by
Oliver Sander
Committed by
sander
11 years ago
Browse files
Options
Downloads
Patches
Plain Diff
Start using a dedicated implementation of a symmetric matrix
[[Imported from SVN: r9582]]
parent
26d119b2
No related branches found
No related tags found
No related merge requests found
Changes
3
Hide whitespace changes
Inline
Side-by-side
Showing
3 changed files
dune/gfe/gramschmidtsolver.hh
+6
-16
6 additions, 16 deletions
dune/gfe/gramschmidtsolver.hh
dune/gfe/symmetricmatrix.hh
+68
-0
68 additions, 0 deletions
dune/gfe/symmetricmatrix.hh
dune/gfe/targetspacertrsolver.cc
+7
-1
7 additions, 1 deletion
dune/gfe/targetspacertrsolver.cc
with
81 additions
and
17 deletions
dune/gfe/gramschmidtsolver.hh
+
6
−
16
View file @
408c0310
...
...
@@ -2,8 +2,8 @@
#define DUNE_GFE_GRAMSCHMIDTSOLVER_HH
#include
<dune/common/fvector.hh>
#include
<dune/common/fmatrix.hh>
#include
<dune/gfe/symmetricmatrix.hh>
/** \brief Direct solver for a dense symmetric linear system, using an orthonormal basis
*
...
...
@@ -24,16 +24,10 @@ class GramSchmidtSolver
* \param matrix The matrix inducing the matrix norm
* \param[in,out] v The vector to normalize
*/
static
void
normalize
(
const
Dune
::
Field
Matrix
<
field_type
,
embeddedDim
,
embeddedDim
>&
matrix
,
static
void
normalize
(
const
Dune
::
Symmetric
Matrix
<
field_type
,
embeddedDim
>&
matrix
,
Dune
::
FieldVector
<
field_type
,
embeddedDim
>&
v
)
{
field_type
energyNormSquared
=
0
;
for
(
int
i
=
0
;
i
<
v
.
size
();
i
++
)
for
(
int
j
=
0
;
j
<
v
.
size
();
j
++
)
energyNormSquared
+=
matrix
[
i
][
j
]
*
v
[
i
]
*
v
[
j
];
v
/=
std
::
sqrt
(
energyNormSquared
);
v
/=
std
::
sqrt
(
matrix
.
energyScalarProduct
(
v
,
v
));
}
...
...
@@ -41,16 +35,12 @@ class GramSchmidtSolver
*
* \param matrix The matrix the defines the scalar product
*/
static
void
project
(
const
Dune
::
Field
Matrix
<
field_type
,
embeddedDim
,
embeddedDim
>&
matrix
,
static
void
project
(
const
Dune
::
Symmetric
Matrix
<
field_type
,
embeddedDim
>&
matrix
,
const
Dune
::
FieldVector
<
field_type
,
embeddedDim
>&
vi
,
Dune
::
FieldVector
<
field_type
,
embeddedDim
>&
vj
)
{
field_type
energyScalarProduct
=
0
;
for
(
int
i
=
0
;
i
<
vi
.
size
();
i
++
)
for
(
int
j
=
0
;
j
<
vj
.
size
();
j
++
)
energyScalarProduct
+=
matrix
[
i
][
j
]
*
vi
[
i
]
*
vj
[
j
];
field_type
energyScalarProduct
=
matrix
.
energyScalarProduct
(
vi
,
vj
);
for
(
int
i
=
0
;
i
<
vj
.
size
();
i
++
)
vj
[
i
]
-=
energyScalarProduct
*
vi
[
i
];
...
...
@@ -59,7 +49,7 @@ class GramSchmidtSolver
public
:
/** Solve linear system by constructing an energy-orthonormal basis */
static
void
solve
(
const
Dune
::
Field
Matrix
<
field_type
,
embeddedDim
,
embeddedDim
>&
matrix
,
static
void
solve
(
const
Dune
::
Symmetric
Matrix
<
field_type
,
embeddedDim
>&
matrix
,
Dune
::
FieldVector
<
field_type
,
embeddedDim
>&
x
,
const
Dune
::
FieldVector
<
field_type
,
embeddedDim
>&
rhs
,
const
Dune
::
FieldMatrix
<
field_type
,
dim
,
embeddedDim
>&
basis
)
...
...
This diff is collapsed.
Click to expand it.
dune/gfe/symmetricmatrix.hh
0 → 100644
+
68
−
0
View file @
408c0310
#ifndef DUNE_GFE_SYMMETRICMATRIX_HH
#define DUNE_GFE_SYMMETRICMATRIX_HH
#include
<dune/common/fvector.hh>
namespace
Dune
{
/** \brief A class implementing a symmetric matrix with compile-time size
*
* A \f$ dim\times dim \f$ matrix is stored internally as a <tt>Dune::FieldVector<double, dim*(dim+1)/2></tt>
* The components are assumed to be \f$ [ E(1,1),\ E(2,1),\ E(2,2),\ E(3,1),\ E(3,2),\ E(3,3) ]\f$
* and analogous for other dimensions
* \tparam dim number of lines/columns of the matrix
*/
template
<
class
T
,
int
N
>
class
SymmetricMatrix
{
public:
/** \brief Default constructor
*
* Tensor is initialized containing zeros if no argument is given.
* \param eye if true tensor is initialized as identity
*/
SymmetricMatrix
()
{}
/** \brief Matrix style random read/write access to components
* \param i line index
* \param j column index
* \note You need to know what you are doing: You can only access the lower
* left triangular submatrix using this methods. It requires i<=j.
*/
T
&
operator
()
(
int
i
,
int
j
)
{
assert
(
i
>=
j
);
return
data_
[
i
*
(
i
+
1
)
/
2
+
j
];
}
/** \brief Matrix style random read access to components
* \param i line index
* \param j column index
* \note You need to know what you are doing: You can only access the lower
* left triangular submatrix using this methods. It requires i<=j.
*/
const
T
&
operator
()
(
int
i
,
int
j
)
const
{
assert
(
i
>=
j
);
return
data_
[
i
*
(
i
+
1
)
/
2
+
j
];
}
T
energyScalarProduct
(
const
FieldVector
<
T
,
N
>&
v1
,
const
FieldVector
<
T
,
N
>&
v2
)
const
{
T
result
=
0
;
for
(
size_t
i
=
0
;
i
<
N
;
i
++
)
for
(
size_t
j
=
0
;
j
<=
i
;
j
++
)
result
+=
(
1
+
(
i
!=
j
))
*
operator
()(
i
,
j
)
*
v1
[
i
]
*
v2
[
j
];
return
result
;
}
private
:
Dune
::
FieldVector
<
T
,
N
*
(
N
+
1
)
/
2
>
data_
;
};
}
#endif
This diff is collapsed.
Click to expand it.
dune/gfe/targetspacertrsolver.cc
+
7
−
1
View file @
408c0310
...
...
@@ -107,8 +107,14 @@ void TargetSpaceRiemannianTRSolver<TargetSpace>::solve()
#ifdef USE_GAUSS_SEIDEL_SOLVER
innerSolver_
->
solve
();
#else
Dune
::
SymmetricMatrix
<
field_type
,
embeddedBlocksize
>
symmetricHessian
;
for
(
size_t
j
=
0
;
j
<
embeddedBlocksize
;
j
++
)
for
(
size_t
k
=
0
;
k
<=
j
;
k
++
)
symmetricHessian
(
j
,
k
)
=
hesseMatrix
[
0
][
0
][
j
][
k
];
Dune
::
FieldMatrix
<
field_type
,
blocksize
,
embeddedBlocksize
>
basis
=
x_
.
orthonormalFrame
();
GramSchmidtSolver
<
field_type
,
blocksize
,
embeddedBlocksize
>::
solve
(
hesseMatrix
[
0
][
0
]
,
corr
[
0
],
rhs
[
0
],
basis
);
GramSchmidtSolver
<
field_type
,
blocksize
,
embeddedBlocksize
>::
solve
(
symmetricHessian
,
corr
[
0
],
rhs
[
0
],
basis
);
#endif
#ifdef USE_GAUSS_SEIDEL_SOLVER
...
...
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