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
bb9f42a1
Commit
bb9f42a1
authored
6 years ago
by
Oliver Sander
Browse files
Options
Downloads
Patches
Plain Diff
Implement MatrixDimension for dense Matrices
parent
9e485731
No related branches found
No related tags found
1 merge request
!268
Implement MatrixDimension for dense Matrices
Pipeline
#15596
passed
6 years ago
Stage: test
Changes
2
Pipelines
1
Hide whitespace changes
Inline
Side-by-side
Showing
2 changed files
dune/istl/matrixutils.hh
+34
-0
34 additions, 0 deletions
dune/istl/matrixutils.hh
dune/istl/test/iotest.cc
+22
-3
22 additions, 3 deletions
dune/istl/test/iotest.cc
with
56 additions
and
3 deletions
dune/istl/matrixutils.hh
+
34
−
0
View file @
bb9f42a1
...
...
@@ -219,6 +219,40 @@ namespace Dune
}
};
// Default implementation for scalar types
template
<
typename
B
,
typename
TA
>
struct
MatrixDimension
<
Matrix
<
B
,
TA
>
>
{
using
block_type
=
typename
Matrix
<
B
,
TA
>::
block_type
;
using
size_type
=
typename
Matrix
<
B
,
TA
>::
size_type
;
static
size_type
rowdim
(
const
Matrix
<
B
,
TA
>&
A
,
size_type
i
)
{
return
MatrixDimension
<
block_type
>::
rowdim
(
A
[
i
][
0
]);
}
static
size_type
coldim
(
const
Matrix
<
B
,
TA
>&
A
,
size_type
c
)
{
return
MatrixDimension
<
block_type
>::
coldim
(
A
[
0
][
c
]);
}
static
size_type
rowdim
(
const
Matrix
<
B
,
TA
>&
A
)
{
size_type
nn
=
0
;
for
(
size_type
i
=
0
;
i
<
A
.
N
();
i
++
)
nn
+=
rowdim
(
A
,
i
);
return
nn
;
}
static
size_type
coldim
(
const
Matrix
<
B
,
TA
>&
A
)
{
size_type
nn
=
0
;
for
(
size_type
i
=
0
;
i
<
A
.
M
();
i
++
)
nn
+=
coldim
(
A
,
i
);
return
nn
;
}
};
template
<
typename
B
,
typename
TA
>
struct
MatrixDimension
<
BCRSMatrix
<
B
,
TA
>
>
...
...
This diff is collapsed.
Click to expand it.
dune/istl/test/iotest.cc
+
22
−
3
View file @
bb9f42a1
...
...
@@ -5,6 +5,7 @@
#include
<dune/common/diagonalmatrix.hh>
#include
<dune/istl/scaledidmatrix.hh>
#include
<dune/istl/bcrsmatrix.hh>
#include
<dune/istl/matrix.hh>
#include
<dune/istl/io.hh>
#include
"laplacian.hh"
...
...
@@ -85,19 +86,37 @@ int main(int argc, char** argv)
Dune
::
writeVectorToMatlabHelper
(
v3
,
std
::
cout
);
// Test the printmatrix method
// BCRSMatrix
{
Dune
::
BCRSMatrix
<
double
>
matrix
;
setupLaplacian
(
matrix
,
3
);
Dune
::
printmatrix
(
std
::
cout
,
matrix
,
"
matrix
"
,
"--"
);
Dune
::
printmatrix
(
std
::
cout
,
matrix
,
"
BCRSMatrix<double>
"
,
"--"
);
}
{
Dune
::
BCRSMatrix
<
Dune
::
FieldMatrix
<
double
,
1
,
1
>
>
matrix
;
setupLaplacian
(
matrix
,
3
);
Dune
::
printmatrix
(
std
::
cout
,
matrix
,
"
matrix
"
,
"--"
);
Dune
::
printmatrix
(
std
::
cout
,
matrix
,
"
BCRSMatrix<FieldMatrix<double,1,1> >
"
,
"--"
);
}
{
Dune
::
BCRSMatrix
<
Dune
::
FieldMatrix
<
double
,
2
,
3
>
>
matrix
;
setupLaplacian
(
matrix
,
3
);
Dune
::
printmatrix
(
std
::
cout
,
matrix
,
"matrix"
,
"--"
);
Dune
::
printmatrix
(
std
::
cout
,
matrix
,
"BCRSMatrix<FieldMatrix<double,2,3> >"
,
"--"
);
}
// Matrix
{
Dune
::
Matrix
<
double
>
matrix
(
3
,
3
);
matrix
=
0
;
Dune
::
printmatrix
(
std
::
cout
,
matrix
,
"Matrix<double>"
,
"--"
);
}
{
Dune
::
Matrix
<
Dune
::
FieldMatrix
<
double
,
1
,
1
>
>
matrix
(
3
,
3
);
matrix
=
0
;
Dune
::
printmatrix
(
std
::
cout
,
matrix
,
"Matrix<FieldMatrix<double,1,1> >"
,
"--"
);
}
{
Dune
::
Matrix
<
Dune
::
FieldMatrix
<
double
,
2
,
3
>
>
matrix
(
3
,
3
);
matrix
=
0
;
Dune
::
printmatrix
(
std
::
cout
,
matrix
,
"Matrix<FieldMatrix<double,2,3> >"
,
"--"
);
}
}
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