Ifpack Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Ifpack_LinePartitioner.h
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 
43 #ifndef IFPACK_LINEPARTITIONER_H
44 #define IFPACK_LINEPARTITIONER_H
45 
46 #include "Ifpack_ConfigDefs.h"
47 #include "Ifpack_Partitioner.h"
50 #include "Teuchos_ParameterList.hpp"
51 class Epetra_Comm;
52 class Ifpack_Graph;
53 class Epetra_Map;
54 class Epetra_BlockMap;
55 class Epetra_Import;
56 class Epetra_RowMatrix;
57 
58 /* \brief Ifpack_LinePartitioner: A class to define partitions into a set of lines.
59 
60 These "line" partitions could then be used in to do block Gauss-Seidel smoothing, for instance.
61 
62 The current implementation uses a local (on processor) line detection in one of two forms.
63 Both are inspired by the work of Mavriplis for convection-diffusion (AIAA Journal, Vol 37(10), 1999).
64 
65 
66 Algorithm 1: Matrix entries
67 
68 Here we use the matrix entries, running lines across the "large" matrix entries, if they are sufficiently
69 large relative to the small ones.
70 
71 
72 Algorithms 2: Coordinates
73 
74 Here we use coordinate information to pick "close" points if they are sufficiently far away
75 from the "far" points. We also make sure the line can never double back on itself, so that
76 the associated sub-matrix could (in theory) be handed off to a fast triangular solver. This
77 implementation doesn't actual do that, however.
78 
79 This implementation is deived from the related routine in ML.
80 
81 
82 Supported parameters:
83  \c "partitioner: line type": "matrix entries" or "coordinates" (char *)
84  \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)
85  \c "partitioner: x-coordinates": x coordinates of local nodes (double *)
86  \c "partitioner: y-coordinates": y coordinates of local nodes (double *)
87  \c "partitioner: z-coordinates": z coordinates of local nodes (double *)
88  \c "partitioner: PDE equations": number of equations per node (integer)
89 
90 */
91 
93 
94 public:
95  // Useful typedef
96  typedef enum {COORDINATES=0, MATRIX_ENTRIES,} LINE_MODE;
97 
98 
102  Matrix_(0),
104  NumEqns_(1),
105  xcoord_(0),
106  ycoord_(0),
107  zcoord_(0),
108  threshold_(0.0)
109  {
110 
111  }
112 
115  Matrix_(Matrix),
117  NumEqns_(1),
118  xcoord_(0),
119  ycoord_(0),
120  zcoord_(0),
121  threshold_(0.0)
122  {
124  Graph_ = &*GraphWrapper_;
125 
126  }
127 
128 
131 
134  {
135  std::string mymode;
137  mymode = List.get("partitioner: line mode",mymode);
138  if(mymode=="coordinates") mode_=COORDINATES;
139  else if(mymode=="matrix entries") mode_=MATRIX_ENTRIES;
140 
141  threshold_ = List.get("partitioner: line detection threshold",threshold_);
142  if(threshold_ < 0.0) IFPACK_CHK_ERR(-1);
143  if(mode_==COORDINATES && threshold_ > 1.0) IFPACK_CHK_ERR(-1);
144 
145 
146  NumEqns_ = List.get("partitioner: PDE equations",NumEqns_);
147  if(NumEqns_ < 1 ) IFPACK_CHK_ERR(-2);
148 
149  xcoord_ = List.get("partitioner: x-coordinates",xcoord_);
150  ycoord_ = List.get("partitioner: y-coordinates",ycoord_);
151  zcoord_ = List.get("partitioner: z-coordinates",zcoord_);
152  if(mode_==COORDINATES && !xcoord_ && !ycoord_ && !zcoord_) IFPACK_CHK_ERR(-3);
153 
154  return(0);
155  }
156 
158  int ComputePartitions();
159 
160 private:
161 
162  // Useful functions
163  int Compute_Blocks_AutoLine(int * blockIndices) const;
164  void local_automatic_line_search(int NumEqns, int * blockIndices, int last, int next, int LineID, double tol, int *itemp, double * dtemp) const;
165 
166  // Stuff I allocated
168 
169  // User data
172  int NumEqns_;
173  double * xcoord_;
174  double * ycoord_;
175  double * zcoord_;
176  double threshold_;
177 
178 
179  // State data
180 };
181 
182 #endif // IFPACK_LINEPARTITIONER_H
Ifpack_LinePartitioner(const Epetra_RowMatrix *Matrix)
Teuchos::RCP< const Ifpack_Graph > GraphWrapper_
virtual ~Ifpack_LinePartitioner()
Destructor.
T & get(ParameterList &l, const std::string &name)
int ComputePartitions()
Computes the partitions. Returns 0 if successful.
int Compute_Blocks_AutoLine(int *blockIndices) const
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
const double tol
Ifpack_LinePartitioner(const Ifpack_Graph *Graph)
Constructor.
adjacency_list< vecS, vecS, undirectedS, no_property, property< edge_weight_t, double > > Graph
Ifpack_Graph_Epetra_RowMatrix: a class to define Ifpack_Graph as a light-weight conversion of Epetra_...
const Ifpack_Graph * Graph_
Reference to the graph to be partitioned.
void local_automatic_line_search(int NumEqns, int *blockIndices, int last, int next, int LineID, double tol, int *itemp, double *dtemp) const
const Epetra_RowMatrix * Matrix_
Ifpack_Graph: a pure virtual class that defines graphs for IFPACK.
Definition: Ifpack_Graph.h:61
int SetPartitionParameters(Teuchos::ParameterList &List)
Sets all the parameters for the partitioner.
#define IFPACK_CHK_ERR(ifpack_err)