IFPACK  Development
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 #include "Ifpack_ConfigDefs.h"
47 #include "Ifpack_Partitioner.h"
48 #include "Ifpack_OverlappingPartitioner.h"
49 #include "Teuchos_ParameterList.hpp"
50 class Epetra_Comm;
51 class Ifpack_Graph;
52 class Epetra_Map;
53 class Epetra_BlockMap;
54 class Epetra_Import;
55 
56 /* \brief Ifpack_LinePartitioner: A class to define partitions into a set of lines.
57 
58 These "line" partitions could then be used in to do block Gauss-Seidel smoothing, for instance.
59 
60 The current implementation uses a local line detection inspired by the work of Mavriplis
61 for convection-diffusion (AIAA Journal, Vol 37(10), 1999).
62 
63 Here we use coordinate information to pick "close" points if they are sufficiently far away
64 from the "far" points. We also make sure the line can never double back on itself, so that
65 the associated sub-matrix could (in theory) be handed off to a fast triangular solver. This
66 implementation doesn't actual do that, however.
67 
68 This implementation is deived from the related routine in ML.
69 
70 Supported parameters:
71  \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)
72  \c "partitioner: x-coordinates": x coordinates of local nodes (double *)
73  \c "partitioner: y-coordinates": y coordinates of local nodes (double *)
74  \c "partitioner: z-coordinates": z coordinates of local nodes (double *)
75  \c "partitioner: PDE equations": number of equations per node (integer)
76 
77 */
78 
80 
81 public:
82 
86  NumEqns_(1),
87  xcoord_(0),
88  ycoord_(0),
89  zcoord_(0),
90  threshold_(0.0)
91  {}
92 
95 
97  int SetPartitionParameters(Teuchos::ParameterList& List)
98  {
99  threshold_ = List.get("partitioner: line detection threshold",threshold_);
100  if(threshold_ < 0.0 || threshold_ > 1.0) IFPACK_CHK_ERR(-1);
101 
102  NumEqns_ = List.get("partitioner: PDE equations",NumEqns_);
103  if(NumEqns_ < 1 ) IFPACK_CHK_ERR(-2);
104 
105  xcoord_ = List.get("partitioner: x-coordinates",xcoord_);
106  ycoord_ = List.get("partitioner: y-coordinates",ycoord_);
107  zcoord_ = List.get("partitioner: z-coordinates",zcoord_);
108  if(!xcoord_ && !ycoord_ && !zcoord_) IFPACK_CHK_ERR(-3);
109 
110  return(0);
111  }
112 
114  int ComputePartitions();
115 
116 private:
117  // Useful functions
118  int Compute_Blocks_AutoLine(int * blockIndices) const;
119  void local_automatic_line_search(int NumEqns, int * blockIndices, int last, int next, int LineID, double tol, int *itemp, double * dtemp) const;
120 
121 
122  // User data
123  int NumEqns_;
124  double * xcoord_;
125  double * ycoord_;
126  double * zcoord_;
127  double threshold_;
128 
129  // State data
130 };
131 
132 #endif // IFPACK_LINEPARTITIONER_H
virtual ~Ifpack_LinePartitioner()
Destructor.
int ComputePartitions()
Computes the partitions. Returns 0 if successful.
Ifpack_LinePartitioner(const Ifpack_Graph *Graph)
Constructor.
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.