Ifpack Package Browser (Single Doxygen Collection)
Development
Toggle main menu visibility
Loading...
Searching...
No Matches
example
Ifpack_ex_Filtering_LL.cpp
Go to the documentation of this file.
1
// @HEADER
2
// ***********************************************************************
3
//
4
// IFPACK
5
// Copyright (2004) Sandia Corporation
6
//
7
// Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8
// license for use of this work by or on behalf of the U.S. Government.
9
//
10
// This library is free software; you can redistribute it and/or modify
11
// it under the terms of the GNU Lesser General Public License as
12
// published by the Free Software Foundation; either version 2.1 of the
13
// License, or (at your option) any later version.
14
//
15
// This library is distributed in the hope that it will be useful, but
16
// WITHOUT ANY WARRANTY; without even the implied warranty of
17
// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18
// Lesser General Public License for more details.
19
//
20
// You should have received a copy of the GNU Lesser General Public
21
// License along with this library; if not, write to the Free Software
22
// Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301
23
// USA
24
// Questions? Contact Michael A. Heroux (maherou@sandia.gov)
25
//
26
// ***********************************************************************
27
// @HEADER
28
29
#include "
Ifpack_ConfigDefs.h
"
30
#ifdef HAVE_MPI
31
#include "
Epetra_MpiComm.h
"
32
#else
33
#include "
Epetra_SerialComm.h
"
34
#endif
35
#include "
Epetra_Map.h
"
36
#include "
Epetra_CrsMatrix.h
"
37
#include "
Ifpack_DropFilter.h
"
38
#include "
Ifpack_SparsityFilter.h
"
39
#include "
Ifpack_SingletonFilter.h
"
40
#include "
Ifpack_Utils.h
"
41
#include "Teuchos_RefCountPtr.hpp"
42
43
int
main
(
int
argc,
char
*argv[])
44
{
45
using
std::cerr;
46
using
std::cout;
47
using
std::endl;
48
49
#ifdef HAVE_MPI
50
MPI_Init(&argc,&argv);
51
Epetra_MpiComm
Comm( MPI_COMM_WORLD );
52
#else
53
Epetra_SerialComm
Comm;
54
#endif
55
56
if
(Comm.
NumProc
() != 1) {
57
cerr <<
"This example must be run with one process only."
<< endl;
58
// exit with success not to break the test harness
59
#ifdef HAVE_MPI
60
MPI_Finalize();
61
#endif
62
exit(EXIT_SUCCESS);
63
}
64
65
long
long
NumPoints = 5;
66
#if !defined(EPETRA_NO_32BIT_GLOBAL_INDICES) || !defined(EPETRA_NO_64BIT_GLOBAL_INDICES)
67
Epetra_Map
Map(NumPoints,0LL,Comm);
68
#else
69
Epetra_Map
Map;
70
#endif
71
72
Teuchos::RefCountPtr<Epetra_CrsMatrix> Matrix = Teuchos::rcp(
new
Epetra_CrsMatrix
(
Copy
,Map,0) );
73
74
std::vector<long long> Indices(NumPoints);
75
std::vector<double> Values(NumPoints);
76
double
Diag = 0.0;
77
78
for
(
long
long
i = 0 ; i < NumPoints ; ++i) {
79
// add a diagonal
80
Matrix->InsertGlobalValues(i,1,&Diag,&i);
81
82
// add off-diagonals
83
int
NumEntries = 0;
84
for
(
long
long
j = i + 1 ; j < NumPoints ; ++j) {
85
Indices[NumEntries] = j;
86
Values[NumEntries] = 1.0 * (j - i);
87
++NumEntries;
88
}
89
#if !defined(EPETRA_NO_32BIT_GLOBAL_INDICES) || !defined(EPETRA_NO_64BIT_GLOBAL_INDICES)
90
Matrix->InsertGlobalValues(i,NumEntries,&Values[0],&Indices[0]);
91
#endif
92
}
93
Matrix->FillComplete();
94
95
// ================================= //
96
// print sparsity of original matrix //
97
// ================================= //
98
99
cout <<
"Sparsity, non-dropped matrix"
<< endl;
100
Ifpack_PrintSparsity_Simple
(*Matrix);
101
102
// ====================================== //
103
// create a new matrix, dropping by value //
104
// ====================================== //
105
//
106
// drop all elements below 4.0. Only the upper-right element
107
// is maintained, plus all the diagonals that are not
108
// considering in dropping.
109
Ifpack_DropFilter
DropA(Matrix,4.0);
110
assert (DropA.
MaxNumEntries
() == 2);
111
112
cout <<
"Sparsity, dropping by value"
<< endl;
113
Ifpack_PrintSparsity_Simple
(DropA);
114
115
// ========================================= //
116
// create a new matrix, dropping by sparsity //
117
// ========================================= //
118
//
119
// Mantain 2 off-diagonal elements.
120
Ifpack_SparsityFilter
SparsityA(Matrix,2);
121
122
cout <<
"Sparsity, dropping by sparsity"
<< endl;
123
Ifpack_PrintSparsity_Simple
(SparsityA);
124
assert (SparsityA.
MaxNumEntries
() == 3);
125
126
// ======================================== //
127
// create new matrices, dropping singletons //
128
// ======================================== //
129
//
130
// If we apply this filter NumPoints - 1 times,
131
// we end up with a one-row matrix
132
Ifpack_SingletonFilter
Filter1(Matrix);
133
Ifpack_SingletonFilter
Filter2(Teuchos::rcp(&Filter1,
false
));
134
Ifpack_SingletonFilter
Filter3(Teuchos::rcp(&Filter2,
false
));
135
Ifpack_SingletonFilter
Filter4(Teuchos::rcp(&Filter3,
false
));
136
137
cout <<
"Sparsity, dropping singletons 4 times"
<< endl;
138
Ifpack_PrintSparsity_Simple
(Filter4);
139
assert (Filter4.
NumMyRows
() == 1);
140
141
#ifdef HAVE_MPI
142
MPI_Finalize();
143
#endif
144
return
(EXIT_SUCCESS);
145
}
146
Epetra_CrsMatrix.h
Copy
Copy
Epetra_Map.h
Epetra_MpiComm.h
Epetra_SerialComm.h
Ifpack_ConfigDefs.h
Ifpack_DropFilter.h
Ifpack_SingletonFilter.h
Ifpack_SparsityFilter.h
Ifpack_PrintSparsity_Simple
void Ifpack_PrintSparsity_Simple(const Epetra_RowMatrix &A)
Definition
Ifpack_Utils.cpp:299
Ifpack_Utils.h
main
int main(int argc, char *argv[])
Definition
Ifpack_ex_Filtering_LL.cpp:43
Epetra_CrsMatrix
Epetra_Map
Epetra_MpiComm
Epetra_MpiComm::NumProc
int NumProc() const
Epetra_SerialComm
Ifpack_DropFilter
Ifpack_DropFilter: Filter based on matrix entries.
Definition
Ifpack_DropFilter.h:81
Ifpack_DropFilter::MaxNumEntries
virtual int MaxNumEntries() const
Returns the maximum number of entries.
Definition
Ifpack_DropFilter.h:99
Ifpack_SingletonFilter
Ifpack_SingletonFilter: Filter based on matrix entries.
Definition
Ifpack_SingletonFilter.h:58
Ifpack_SingletonFilter::NumMyRows
virtual int NumMyRows() const
Definition
Ifpack_SingletonFilter.h:178
Ifpack_SparsityFilter
Ifpack_SparsityFilter: a class to drop based on sparsity.
Definition
Ifpack_SparsityFilter.h:58
Ifpack_SparsityFilter::MaxNumEntries
virtual int MaxNumEntries() const
Definition
Ifpack_SparsityFilter.h:73
Generated by
1.17.0