Skip to content

Commit 7e41fc5

Browse files
committed
Add test for ShellMatrix.
1 parent f2ddeeb commit 7e41fc5

File tree

1 file changed

+167
-0
lines changed

1 file changed

+167
-0
lines changed

tests/numerics/shell_matrix.C

Lines changed: 167 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,167 @@
1+
#include "libmesh/shell_matrix.h"
2+
#include "libmesh/enum_solver_type.h"
3+
#include "libmesh/sparse_shell_matrix.h"
4+
#ifdef LIBMESH_HAVE_PETSC
5+
#include "libmesh/petsc_linear_solver.h"
6+
#include "libmesh/petsc_vector.h"
7+
#include "libmesh/petsc_matrix.h"
8+
#include "libmesh/petsc_shell_matrix.h"
9+
10+
#include "libmesh_cppunit.h"
11+
#include "test_comm.h"
12+
13+
using namespace libMesh;
14+
15+
class ShellMatrixTest : public CppUnit::TestCase
16+
{
17+
public:
18+
LIBMESH_CPPUNIT_TEST_SUITE(ShellMatrixTest);
19+
20+
CPPUNIT_TEST(testPetscShell);
21+
CPPUNIT_TEST(testSparseShell_petsc);
22+
23+
CPPUNIT_TEST_SUITE_END();
24+
25+
public:
26+
void setUp()
27+
{
28+
_comm = TestCommWorld;
29+
}
30+
31+
void tearDown() {}
32+
33+
void testPetscShell()
34+
{
35+
LOG_UNIT_TEST;
36+
37+
// Depending on the config: Use different vector-types here!
38+
PetscLinearSolver<Number> linear_solver(*_comm);
39+
linear_solver.init();
40+
auto ierr = KSPSetType (linear_solver.ksp(), const_cast<KSPType>(KSPGMRES));
41+
CHKERRABORT(linear_solver.comm().get(), ierr);
42+
ierr = PCSetType (linear_solver.pc(), const_cast<PCType>(PCNONE));
43+
CHKERRABORT(linear_solver.comm().get(), ierr);
44+
KSPSetInitialGuessNonzero(linear_solver.ksp(), PETSC_TRUE);
45+
46+
unsigned int maxits=20;
47+
Real tol = 1e-3;
48+
49+
UnityShellMat mat = (*_comm);
50+
51+
// Depending on the config: Use different vector-types here!
52+
PetscVector<Number> rhs(*_comm, 20);
53+
rhs.set(5,2.);
54+
rhs.close();
55+
auto solution = rhs.clone();
56+
solution->set(7,2.);
57+
solution->close();
58+
//auto rval = linear_solver.solve(mat, *solution, rhs, tol, maxits);
59+
linear_solver.solve(mat, *solution, rhs, tol, maxits);
60+
61+
for (numeric_index_type m=0; m < rhs.size(); ++m)
62+
CPPUNIT_ASSERT_EQUAL(solution->el(m), rhs.el(m));
63+
}
64+
65+
void testSparseShell_petsc()
66+
{
67+
LOG_UNIT_TEST;
68+
69+
// Depending on the config: Use different vector-types here!
70+
PetscLinearSolver<Number> linear_solver(*_comm);
71+
linear_solver.init();
72+
auto ierr = KSPSetType (linear_solver.ksp(), const_cast<KSPType>(KSPGMRES));
73+
CHKERRABORT(linear_solver.comm().get(), ierr);
74+
ierr = PCSetType (linear_solver.pc(), const_cast<PCType>(PCNONE));
75+
CHKERRABORT(linear_solver.comm().get(), ierr);
76+
KSPSetInitialGuessNonzero(linear_solver.ksp(), PETSC_TRUE);
77+
78+
PetscMatrix<Number> petsc_mat(*_comm);
79+
petsc_mat.init(20,20,20,20,1,0,1); // hard-coded dimension.
80+
set_unity(petsc_mat);
81+
petsc_mat.close();
82+
SparseShellMatrix<Number> mat(petsc_mat);
83+
84+
unsigned int maxits=20;
85+
Real tol = 1e-3;
86+
87+
// Depending on the config: Use different vector-types here!
88+
PetscVector<Number> rhs(*_comm, 20);
89+
rhs.set(5,2.);
90+
rhs.close();
91+
auto solution = rhs.clone();
92+
solution->set(7,2.);
93+
solution->close();
94+
//auto rval = linear_solver.solve(mat, *solution, rhs, tol, maxits);
95+
linear_solver.solve(mat, *solution, rhs, tol, maxits);
96+
97+
for (numeric_index_type m=0; m < rhs.size(); ++m)
98+
CPPUNIT_ASSERT_EQUAL(solution->el(m), rhs.el(m));
99+
100+
// Make the matrix singular:
101+
// a) set an irrelevant element to '0':
102+
solution->set(7,2.);
103+
petsc_mat.set(2,2,0.);
104+
auto rval = linear_solver.solve(mat, *solution, rhs, tol, maxits);
105+
// In this case, we can solve the equation.
106+
CPPUNIT_ASSERT_EQUAL(rval.second, 0.00);
107+
// b) set the relevant element to '0':
108+
petsc_mat.set(5,5,0.);
109+
rval = linear_solver.solve(mat, *solution, rhs, tol, maxits);
110+
// In this case, we can not solve the equation.
111+
// The best solution is the 0-vector and we have an error of 2.0
112+
CPPUNIT_ASSERT_EQUAL(rval.second, 2.00);
113+
114+
}
115+
116+
private:
117+
Parallel::Communicator * _comm;
118+
119+
void set_unity(SparseMatrix<Number> & M)
120+
{
121+
for (libMesh::numeric_index_type n=0; n < M.m(); ++n)
122+
M.set(n,n,1.);
123+
}
124+
125+
class UnityShellMat : public libMesh::ShellMatrix<libMesh::Number>
126+
{
127+
public:
128+
UnityShellMat(Parallel::Communicator & comm):
129+
ShellMatrix(comm){}
130+
131+
virtual ~UnityShellMat() = default;
132+
133+
virtual libMesh::numeric_index_type m () const override
134+
{
135+
return 20;
136+
}
137+
virtual libMesh::numeric_index_type n () const override
138+
{
139+
return 20;
140+
}
141+
142+
virtual void vector_mult(libMesh::NumericVector<libMesh::Number> &dest,
143+
const libMesh::NumericVector< libMesh::Number > &arg) const override
144+
{
145+
dest = arg;
146+
}
147+
148+
virtual void vector_mult_add (libMesh::NumericVector<libMesh::Number> &dest,
149+
const libMesh::NumericVector< libMesh::Number > &arg) const override
150+
{
151+
dest += arg;
152+
}
153+
154+
virtual void get_diagonal (libMesh::NumericVector<libMesh::Number> &dest) const override
155+
{
156+
157+
for (numeric_index_type m=0; m < dest.local_size(); ++m)
158+
dest.set(m, 1.);
159+
}
160+
161+
};
162+
163+
};
164+
165+
CPPUNIT_TEST_SUITE_REGISTRATION(ShellMatrixTest);
166+
167+
#endif // #ifdef LIBMESH_HAVE_PETSC

0 commit comments

Comments
 (0)