00001 /*@HEADER 00002 // *********************************************************************** 00003 // 00004 // Ifpack: Object-Oriented Algebraic Preconditioner Package 00005 // Copyright (2002) Sandia Corporation 00006 // 00007 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive 00008 // license for use of this work by or on behalf of the U.S. Government. 00009 // 00010 // Redistribution and use in source and binary forms, with or without 00011 // modification, are permitted provided that the following conditions are 00012 // met: 00013 // 00014 // 1. Redistributions of source code must retain the above copyright 00015 // notice, this list of conditions and the following disclaimer. 00016 // 00017 // 2. Redistributions in binary form must reproduce the above copyright 00018 // notice, this list of conditions and the following disclaimer in the 00019 // documentation and/or other materials provided with the distribution. 00020 // 00021 // 3. Neither the name of the Corporation nor the names of the 00022 // contributors may be used to endorse or promote products derived from 00023 // this software without specific prior written permission. 00024 // 00025 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY 00026 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE 00027 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR 00028 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE 00029 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, 00030 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, 00031 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR 00032 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF 00033 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING 00034 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS 00035 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE. 00036 // 00037 // Questions? Contact Michael A. Heroux (maherou@sandia.gov) 00038 // 00039 // *********************************************************************** 00040 //@HEADER 00041 */ 00042 00043 #ifndef IFPACK_LINEPARTITIONER_H 00044 #define IFPACK_LINEPARTITIONER_H 00045 00046 #include "Ifpack_ConfigDefs.h" 00047 #include "Ifpack_Partitioner.h" 00048 #include "Ifpack_OverlappingPartitioner.h" 00049 #include "Teuchos_ParameterList.hpp" 00050 class Epetra_Comm; 00051 class Ifpack_Graph; 00052 class Epetra_Map; 00053 class Epetra_BlockMap; 00054 class Epetra_Import; 00055 00056 /* \brief Ifpack_LinePartitioner: A class to define partitions into a set of lines. 00057 00058 These "line" partitions could then be used in to do block Gauss-Seidel smoothing, for instance. 00059 00060 The current implementation uses a local line detection inspired by the work of Mavriplis 00061 for convection-diffusion (AIAA Journal, Vol 37(10), 1999). 00062 00063 Here we use coordinate information to pick "close" points if they are sufficiently far away 00064 from the "far" points. We also make sure the line can never double back on itself, so that 00065 the associated sub-matrix could (in theory) be handed off to a fast triangular solver. This 00066 implementation doesn't actual do that, however. 00067 00068 This implementation is deived from the related routine in ML. 00069 00070 Supported parameters: 00071 \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) 00072 \c "partitioner: x-coordinates": x coordinates of local nodes (double *) 00073 \c "partitioner: y-coordinates": y coordinates of local nodes (double *) 00074 \c "partitioner: z-coordinates": z coordinates of local nodes (double *) 00075 \c "partitioner: PDE equations": number of equations per node (integer) 00076 00077 */ 00078 00079 class Ifpack_LinePartitioner : public Ifpack_OverlappingPartitioner { 00080 00081 public: 00082 00084 Ifpack_LinePartitioner(const Ifpack_Graph* Graph) : 00085 Ifpack_OverlappingPartitioner(Graph), 00086 NumEqns_(1), 00087 xcoord_(0), 00088 ycoord_(0), 00089 zcoord_(0), 00090 threshold_(0.0) 00091 {} 00092 00094 virtual ~Ifpack_LinePartitioner() {}; 00095 00097 int SetPartitionParameters(Teuchos::ParameterList& List) 00098 { 00099 threshold_ = List.get("partitioner: line detection threshold",threshold_); 00100 if(threshold_ < 0.0 || threshold_ > 1.0) IFPACK_CHK_ERR(-1); 00101 00102 NumEqns_ = List.get("partitioner: PDE equations",NumEqns_); 00103 if(NumEqns_ < 1 ) IFPACK_CHK_ERR(-2); 00104 00105 xcoord_ = List.get("partitioner: x-coordinates",xcoord_); 00106 ycoord_ = List.get("partitioner: y-coordinates",ycoord_); 00107 zcoord_ = List.get("partitioner: z-coordinates",zcoord_); 00108 if(!xcoord_ && !ycoord_ && !zcoord_) IFPACK_CHK_ERR(-3); 00109 00110 return(0); 00111 } 00112 00114 int ComputePartitions(); 00115 00116 private: 00117 // Useful functions 00118 int Compute_Blocks_AutoLine(int * blockIndices) const; 00119 void local_automatic_line_search(int NumEqns, int * blockIndices, int last, int next, int LineID, double tol, int *itemp, double * dtemp) const; 00120 00121 00122 // User data 00123 int NumEqns_; 00124 double * xcoord_; 00125 double * ycoord_; 00126 double * zcoord_; 00127 double threshold_; 00128 00129 // State data 00130 }; 00131 00132 #endif // IFPACK_LINEPARTITIONER_H
1.7.6.1