IFPACK Development
Loading...
Searching...
No Matches
Ifpack_LinePartitioner.h
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
43#ifndef IFPACK_LINEPARTITIONER_H
44#define IFPACK_LINEPARTITIONER_H
45
46#if defined(Ifpack_SHOW_DEPRECATED_WARNINGS)
47#ifdef __GNUC__
48#warning "The Ifpack package is deprecated"
49#endif
50#endif
51
52#include "Ifpack_ConfigDefs.h"
53#include "Ifpack_Partitioner.h"
54#include "Ifpack_OverlappingPartitioner.h"
55#include "Ifpack_Graph_Epetra_RowMatrix.h"
56#include "Teuchos_ParameterList.hpp"
57class Epetra_Comm;
58class Ifpack_Graph;
59class Epetra_Map;
60class Epetra_BlockMap;
61class Epetra_Import;
63
64/* \brief Ifpack_LinePartitioner: A class to define partitions into a set of lines.
65
66These "line" partitions could then be used in to do block Gauss-Seidel smoothing, for instance.
67
68The current implementation uses a local (on processor) line detection in one of two forms.
69Both are inspired by the work of Mavriplis for convection-diffusion (AIAA Journal, Vol 37(10), 1999).
70
71
72Algorithm 1: Matrix entries
73
74Here we use the matrix entries, running lines across the "large" matrix entries, if they are sufficiently
75large relative to the small ones.
76
77
78Algorithms 2: Coordinates
79
80Here we use coordinate information to pick "close" points if they are sufficiently far away
81from the "far" points. We also make sure the line can never double back on itself, so that
82the associated sub-matrix could (in theory) be handed off to a fast triangular solver. This
83implementation doesn't actual do that, however.
84
85This implementation is deived from the related routine in ML.
86
87
88Supported parameters:
89 \c "partitioner: line type": "matrix entries" or "coordinates" (char *)
90 \c "partitioner: line detection threshold": if ||x_j - x_i||^2 < thresh * max_k||x_k - x_i||^2, then the points are close enough to line smooth (double)
91 \c "partitioner: x-coordinates": x coordinates of local nodes (double *)
92 \c "partitioner: y-coordinates": y coordinates of local nodes (double *)
93 \c "partitioner: z-coordinates": z coordinates of local nodes (double *)
94 \c "partitioner: PDE equations": number of equations per node (integer)
95
96*/
97
99
100public:
101 // Useful typedef
102 typedef enum {COORDINATES=0, MATRIX_ENTRIES,} LINE_MODE;
103
104
108 Matrix_(0),
109 mode_(COORDINATES),
110 NumEqns_(1),
111 xcoord_(0),
112 ycoord_(0),
113 zcoord_(0),
114 threshold_(0.0)
115 {
116
117 }
118
121 Matrix_(Matrix),
122 mode_(COORDINATES),
123 NumEqns_(1),
124 xcoord_(0),
125 ycoord_(0),
126 zcoord_(0),
127 threshold_(0.0)
128 {
129 GraphWrapper_ = Teuchos::rcp(new Ifpack_Graph_Epetra_RowMatrix(Teuchos::rcp(Matrix,false)));
130 Graph_ = &*GraphWrapper_;
131
132 }
133
134
137
139 int SetPartitionParameters(Teuchos::ParameterList& List)
140 {
141 std::string mymode;
142 mode_=COORDINATES;
143 mymode = List.get("partitioner: line mode",mymode);
144 if(mymode=="coordinates") mode_=COORDINATES;
145 else if(mymode=="matrix entries") mode_=MATRIX_ENTRIES;
146
147 threshold_ = List.get("partitioner: line detection threshold",threshold_);
148 if(threshold_ < 0.0) IFPACK_CHK_ERR(-1);
149 if(mode_==COORDINATES && threshold_ > 1.0) IFPACK_CHK_ERR(-1);
150
151
152 NumEqns_ = List.get("partitioner: PDE equations",NumEqns_);
153 if(NumEqns_ < 1 ) IFPACK_CHK_ERR(-2);
154
155 xcoord_ = List.get("partitioner: x-coordinates",xcoord_);
156 ycoord_ = List.get("partitioner: y-coordinates",ycoord_);
157 zcoord_ = List.get("partitioner: z-coordinates",zcoord_);
158 if(mode_==COORDINATES && !xcoord_ && !ycoord_ && !zcoord_) IFPACK_CHK_ERR(-3);
159
160 return(0);
161 }
162
164 int ComputePartitions();
165
166private:
167
168 // Useful functions
169 int Compute_Blocks_AutoLine(int * blockIndices) const;
170 void local_automatic_line_search(int NumEqns, int * blockIndices, int last, int next, int LineID, double tol, int *itemp, double * dtemp) const;
171
172 // Stuff I allocated
173 Teuchos::RCP<const Ifpack_Graph> GraphWrapper_;
174
175 // User data
176 const Epetra_RowMatrix* Matrix_;
177 LINE_MODE mode_;
178 int NumEqns_;
179 double * xcoord_;
180 double * ycoord_;
181 double * zcoord_;
182 double threshold_;
183
184
185 // State data
186};
187
188#endif // IFPACK_LINEPARTITIONER_H
Ifpack_Graph_Epetra_RowMatrix: a class to define Ifpack_Graph as a light-weight conversion of Epetra_...
Ifpack_Graph: a pure virtual class that defines graphs for IFPACK.
int SetPartitionParameters(Teuchos::ParameterList &List)
Sets all the parameters for the partitioner.
int ComputePartitions()
Computes the partitions. Returns 0 if successful.
Ifpack_LinePartitioner(const Ifpack_Graph *Graph)
Constructor.
virtual ~Ifpack_LinePartitioner()
Destructor.
const Ifpack_Graph * Graph_
Reference to the graph to be partitioned.