Ifpack Package Browser (Single Doxygen Collection)
Development
Toggle main menu visibility
Loading...
Searching...
No Matches
example
Ifpack_ex_Reordering.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_Reordering.h
"
51
#include "
Ifpack_RCMReordering.h
"
52
#include "
Ifpack_ReorderFilter.h
"
53
#include "
Ifpack_Utils.h
"
54
#include "Teuchos_RefCountPtr.hpp"
55
56
//==============================================================================
57
int
main
(
int
argc,
char
*argv[])
58
{
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
// only one process
70
if
(Comm.
NumProc
() != 1) {
71
#ifdef HAVE_MPI
72
MPI_Finalize();
73
#endif
74
if
(Comm.
MyPID
() == 0)
75
cout <<
"Please run this test with one process only"
<< endl;
76
// return success not to break the tests
77
exit(EXIT_SUCCESS);
78
}
79
80
// ======================================================== //
81
// now create the famous "upper arrow" matrix, which //
82
// should be reordered as a "lower arrow". Sparsity pattern //
83
// will be printed on screen. //
84
// ======================================================== //
85
86
int
NumPoints = 16;
87
88
#if !defined(EPETRA_NO_32BIT_GLOBAL_INDICES) || !defined(EPETRA_NO_64BIT_GLOBAL_INDICES)
89
Epetra_Map
Map(-1,NumPoints,0,Comm);
90
#else
91
Epetra_Map
Map;
92
#endif
93
94
std::vector<int> Indices(NumPoints);
95
std::vector<double> Values(NumPoints);
96
97
Teuchos::RefCountPtr<Epetra_CrsMatrix> A = Teuchos::rcp(
new
Epetra_CrsMatrix
(
Copy
,Map,0) );
98
for
(
int
i = 0 ; i < NumPoints ; ++i) {
99
100
int
NumEntries;
101
if
(i == 0) {
102
NumEntries = NumPoints;
103
for
(
int
j = 0 ; j < NumPoints ; ++j) {
104
Indices[j] = j;
105
Values[j] = 1.0;
106
}
107
}
108
else
{
109
NumEntries = 2;
110
Indices[0] = 0;
111
Indices[1] = i;
112
Values[0] = 1.0;
113
Values[1] = 1.0;
114
}
115
116
#if !defined(EPETRA_NO_32BIT_GLOBAL_INDICES) || !defined(EPETRA_NO_64BIT_GLOBAL_INDICES)
117
A->InsertGlobalValues(i, NumEntries, &Values[0], &Indices[0]);
118
#endif
119
}
120
121
A->FillComplete();
122
123
// print the sparsity to file, postscript format
125
126
// create the reordering...
127
Teuchos::RefCountPtr<Ifpack_RCMReordering> Reorder = Teuchos::rcp(
new
Ifpack_RCMReordering
() );
128
// and compute is on A
129
IFPACK_CHK_ERR
(Reorder->Compute(*A));
130
131
// cout information
132
cout << *Reorder;
133
134
// create a reordered matrix
135
Ifpack_ReorderFilter
ReordA(A, Reorder);
136
137
// print the sparsity to file, postscript format
139
140
#ifdef HAVE_MPI
141
MPI_Finalize();
142
#endif
143
return
(EXIT_SUCCESS);
144
145
}
Epetra_CrsMatrix.h
Copy
Copy
Epetra_Map.h
Epetra_MpiComm.h
Epetra_SerialComm.h
Ifpack_ConfigDefs.h
IFPACK_CHK_ERR
#define IFPACK_CHK_ERR(ifpack_err)
Definition
Ifpack_ConfigDefs.h:125
Ifpack_RCMReordering.h
Ifpack_ReorderFilter.h
Ifpack_Reordering.h
Ifpack_Utils.h
main
int main(int argc, char *argv[])
Definition
Ifpack_ex_Reordering.cpp:57
Epetra_CrsMatrix
Epetra_Map
Epetra_MpiComm
Epetra_MpiComm::NumProc
int NumProc() const
Epetra_MpiComm::MyPID
int MyPID() const
Epetra_SerialComm
Ifpack_RCMReordering
Ifpack_RCMReordering: reverse Cuthill-McKee reordering.
Definition
Ifpack_RCMReordering.h:58
Ifpack_ReorderFilter
Ifpack_ReorderFilter: a class for light-weight reorder of local rows and columns of an Epetra_RowMatr...
Definition
Ifpack_ReorderFilter.h:81
Generated by
1.17.0