Ifpack Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Ifpack_ex_Filtering.cpp
Go to the documentation of this file.
1 //@HEADER
2 // ***********************************************************************
3 //
4 // Ifpack: Object-Oriented Algebraic Preconditioner Package
5 // Copyright (2002) 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 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38 //
39 // ***********************************************************************
40 //@HEADER
41 
42 #include "Ifpack_ConfigDefs.h"
43 #ifdef HAVE_MPI
44 #include "Epetra_MpiComm.h"
45 #else
46 #include "Epetra_SerialComm.h"
47 #endif
48 #include "Epetra_Map.h"
49 #include "Epetra_CrsMatrix.h"
50 #include "Ifpack_DropFilter.h"
51 #include "Ifpack_SparsityFilter.h"
52 #include "Ifpack_SingletonFilter.h"
53 #include "Ifpack_Utils.h"
54 #include "Teuchos_RefCountPtr.hpp"
55 
56 int main(int argc, char *argv[])
57 {
58  using std::cerr;
59  using std::cout;
60  using std::endl;
61 
62 #ifdef HAVE_MPI
63  MPI_Init(&argc,&argv);
64  Epetra_MpiComm Comm( MPI_COMM_WORLD );
65 #else
66  Epetra_SerialComm Comm;
67 #endif
68 
69  if (Comm.NumProc() != 1) {
70  cerr << "This example must be run with one process only." << endl;
71  // exit with success not to break the test harness
72 #ifdef HAVE_MPI
73  MPI_Finalize();
74 #endif
75  exit(EXIT_SUCCESS);
76  }
77 
78  int NumPoints = 5;
79 #if !defined(EPETRA_NO_32BIT_GLOBAL_INDICES) || !defined(EPETRA_NO_64BIT_GLOBAL_INDICES)
80  Epetra_Map Map(NumPoints,0,Comm);
81 #else
82  Epetra_Map Map;
83 #endif
84 
85  Teuchos::RefCountPtr<Epetra_CrsMatrix> Matrix = Teuchos::rcp( new Epetra_CrsMatrix(Copy,Map,0) );
86 
87  std::vector<int> Indices(NumPoints);
88  std::vector<double> Values(NumPoints);
89  double Diag = 0.0;
90 
91  for (int i = 0 ; i < NumPoints ; ++i) {
92  // add a diagonal
93  Matrix->InsertGlobalValues(i,1,&Diag,&i);
94 
95  // add off-diagonals
96  int NumEntries = 0;
97  for (int j = i + 1 ; j < NumPoints ; ++j) {
98  Indices[NumEntries] = j;
99  Values[NumEntries] = 1.0 * (j - i);
100  ++NumEntries;
101  }
102 #if !defined(EPETRA_NO_32BIT_GLOBAL_INDICES) || !defined(EPETRA_NO_64BIT_GLOBAL_INDICES)
103  Matrix->InsertGlobalValues(i,NumEntries,&Values[0],&Indices[0]);
104 #endif
105  }
106  Matrix->FillComplete();
107 
108  // ================================= //
109  // print sparsity of original matrix //
110  // ================================= //
111 
112  cout << "Sparsity, non-dropped matrix" << endl;
114 
115  // ====================================== //
116  // create a new matrix, dropping by value //
117  // ====================================== //
118  //
119  // drop all elements below 4.0. Only the upper-right element
120  // is maintained, plus all the diagonals that are not
121  // considering in dropping.
122  Ifpack_DropFilter DropA(Matrix,4.0);
123  assert (DropA.MaxNumEntries() == 2);
124 
125  cout << "Sparsity, dropping by value" << endl;
127 
128  // ========================================= //
129  // create a new matrix, dropping by sparsity //
130  // ========================================= //
131  //
132  // Mantain 2 off-diagonal elements.
133  Ifpack_SparsityFilter SparsityA(Matrix,2);
134 
135  cout << "Sparsity, dropping by sparsity" << endl;
136  Ifpack_PrintSparsity_Simple(SparsityA);
137  assert (SparsityA.MaxNumEntries() == 3);
138 
139  // ======================================== //
140  // create new matrices, dropping singletons //
141  // ======================================== //
142  //
143  // If we apply this filter NumPoints - 1 times,
144  // we end up with a one-row matrix
145  Ifpack_SingletonFilter Filter1(Matrix);
146  Ifpack_SingletonFilter Filter2(Teuchos::rcp(&Filter1, false));
147  Ifpack_SingletonFilter Filter3(Teuchos::rcp(&Filter2, false));
148  Ifpack_SingletonFilter Filter4(Teuchos::rcp(&Filter3, false));
149 
150  cout << "Sparsity, dropping singletons 4 times" << endl;
152  assert (Filter4.NumMyRows() == 1);
153 
154 #ifdef HAVE_MPI
155  MPI_Finalize();
156 #endif
157  return(EXIT_SUCCESS);
158 }
159 
virtual int NumMyRows() const
virtual int MaxNumEntries() const
Returns the maximum number of entries.
Ifpack_DropFilter: Filter based on matrix entries.
virtual int MaxNumEntries() const
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
void Ifpack_PrintSparsity_Simple(const Epetra_RowMatrix &A)
Ifpack_SingletonFilter: Filter based on matrix entries.
int main(int argc, char *argv[])
int NumProc() const
Ifpack_SparsityFilter: a class to drop based on sparsity.