/usr/include/libmesh/linear_solver.h is in libmesh-dev 0.7.1-2ubuntu1.
This file is owned by root:root, with mode 0o644.
The actual contents of the file can be viewed below.
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 | // $Id: linear_solver.h 4278 2011-03-21 15:23:30Z roystgnr $
// The libMesh Finite Element Library.
// Copyright (C) 2002-2008 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
// This library is free software; you can redistribute it and/or
// modify it under the terms of the GNU Lesser General Public
// License as published by the Free Software Foundation; either
// version 2.1 of the License, or (at your option) any later version.
// This library is distributed in the hope that it will be useful,
// but WITHOUT ANY WARRANTY; without even the implied warranty of
// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
// Lesser General Public License for more details.
// You should have received a copy of the GNU Lesser General Public
// License along with this library; if not, write to the Free Software
// Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
#ifndef __linear_solver_h__
#define __linear_solver_h__
// C++ includes
#include <vector>
// Local includes
#include "libmesh_common.h"
#include "enum_solver_package.h"
#include "enum_solver_type.h"
#include "enum_preconditioner_type.h"
#include "enum_subset_solve_mode.h"
#include "reference_counted_object.h"
#include "libmesh.h"
namespace libMesh
{
// forward declarations
template <typename T> class AutoPtr;
template <typename T> class SparseMatrix;
template <typename T> class NumericVector;
template <typename T> class ShellMatrix;
template <typename T> class Preconditioner;
/**
* This class provides a uniform interface for linear solvers. This base
* class is overloaded to provide linear solvers from different packages
* like PETSC or LASPACK.
*
* @author Benjamin Kirk, 2003
*/
template <typename T>
class LinearSolver : public ReferenceCountedObject<LinearSolver<T> >
{
public:
/**
* Constructor. Initializes Solver data structures
*/
LinearSolver ();
/**
* Destructor.
*/
virtual ~LinearSolver ();
/**
* Builds a \p LinearSolver using the linear solver package specified by
* \p solver_package
*/
static AutoPtr<LinearSolver<T> > build(const SolverPackage solver_package =
libMesh::default_solver_package());
/**
* @returns true if the data structures are
* initialized, false otherwise.
*/
bool initialized () const { return _is_initialized; }
/**
* Release all memory and clear data structures.
*/
virtual void clear () {}
/**
* Initialize data structures if not done so already.
*/
virtual void init () = 0;
/**
* Returns the type of solver to use.
*/
SolverType solver_type () const { return _solver_type; }
/**
* Sets the type of solver to use.
*/
void set_solver_type (const SolverType st)
{ _solver_type = st; }
/**
* Returns the type of preconditioner to use.
*/
PreconditionerType preconditioner_type () const;
/**
* Sets the type of preconditioner to use.
*/
void set_preconditioner_type (const PreconditionerType pct);
/**
* Attaches a Preconditioner object to be used
*/
void attach_preconditioner(Preconditioner<T> * preconditioner);
/**
* After calling this method, all successive solves will be
* restricted to the given set of dofs, which must contain local
* dofs on each processor only and not contain any duplicates. This
* mode can be disabled by calling this method with \p dofs being a
* \p NULL pointer.
*/
virtual void restrict_solve_to (const std::vector<unsigned int>* const dofs,
const SubsetSolveMode subset_solve_mode=SUBSET_ZERO);
/**
* This function calls the solver
* "_solver_type" preconditioned with the
* "_preconditioner_type" preconditioner. Note that this method
* will compute the preconditioner from the system matrix.
*/
virtual std::pair<unsigned int, Real> solve (SparseMatrix<T>&, // System Matrix
NumericVector<T>&, // Solution vector
NumericVector<T>&, // RHS vector
const double, // Stopping tolerance
const unsigned int) = 0; // N. Iterations
/**
* This function calls the solver
* "_solver_type" preconditioned with the
* "_preconditioner_type" preconditioner.
*/
virtual std::pair<unsigned int, Real> solve (SparseMatrix<T>&, // System Matrix
SparseMatrix<T>&, // Preconditioning Matrix
NumericVector<T>&, // Solution vector
NumericVector<T>&, // RHS vector
const double, // Stopping tolerance
const unsigned int) = 0; // N. Iterations
/**
* This function calls the solver "_solver_type" preconditioned with
* the "_preconditioner_type" preconditioner. The preconditioning
* matrix is used if it is provided, or the system matrix is used if
* \p precond_matrix is null
*/
std::pair<unsigned int, Real> solve (SparseMatrix<T>& matrix,
SparseMatrix<T>* precond_matrix,
NumericVector<T>&, // Solution vector
NumericVector<T>&, // RHS vector
const double, // Stopping tolerance
const unsigned int); // N. Iterations
/**
* This function solves a system whose matrix is a shell matrix.
*/
virtual std::pair<unsigned int, Real> solve (const ShellMatrix<T>& shell_matrix,
NumericVector<T>&, // Solution vector
NumericVector<T>&, // RHS vector
const double, // Stopping tolerance
const unsigned int) = 0; // N. Iterations
/**
* This function solves a system whose matrix is a shell matrix, but
* a sparse matrix is used as preconditioning matrix, this allowing
* other preconditioners than JACOBI.
*/
virtual std::pair<unsigned int, Real> solve (const ShellMatrix<T>& shell_matrix,
const SparseMatrix<T>& precond_matrix,
NumericVector<T>&, // Solution vector
NumericVector<T>&, // RHS vector
const double, // Stopping tolerance
const unsigned int) = 0; // N. Iterations
/**
* This function solves a system whose matrix is a shell matrix, but
* an optional sparse matrix may be used as preconditioning matrix.
*/
std::pair<unsigned int, Real> solve (const ShellMatrix<T>& matrix,
const SparseMatrix<T>* precond_matrix,
NumericVector<T>&, // Solution vector
NumericVector<T>&, // RHS vector
const double, // Stopping tolerance
const unsigned int); // N. Iterations
/**
* Prints a useful message about why the latest linear solve
* con(di)verged.
*/
virtual void print_converged_reason() = 0;
/**
* Boolean flag to indicate whether we want to use an identical
* preconditioner to the previous solve. This can save
* substantial work in the cases where the system matrix is
* the same for successive solves.
*/
bool same_preconditioner;
protected:
/**
* Enum stating which type of iterative solver to use.
*/
SolverType _solver_type;
/**
* Enum statitng with type of preconditioner to use.
*/
PreconditionerType _preconditioner_type;
/**
* Flag indicating if the data structures have been initialized.
*/
bool _is_initialized;
/**
* Holds the Preconditioner object to be used for the linear solves.
*/
Preconditioner<T> * _preconditioner;
};
/*----------------------- inline functions ----------------------------------*/
template <typename T>
inline
LinearSolver<T>::LinearSolver () :
same_preconditioner (false),
_solver_type (GMRES),
_preconditioner_type (ILU_PRECOND),
_is_initialized (false),
_preconditioner (NULL)
{
}
template <typename T>
inline
LinearSolver<T>::~LinearSolver ()
{
this->clear ();
}
template <typename T>
inline
std::pair<unsigned int, Real>
LinearSolver<T>::solve (SparseMatrix<T>& mat,
SparseMatrix<T>* pc_mat,
NumericVector<T>& sol,
NumericVector<T>& rhs,
const double tol,
const unsigned int n_iter)
{
if (pc_mat)
return this->solve(mat, *pc_mat, sol, rhs, tol, n_iter);
else
return this->solve(mat, sol, rhs, tol, n_iter);
}
template <typename T>
inline
std::pair<unsigned int, Real>
LinearSolver<T>::solve (const ShellMatrix<T>& mat,
const SparseMatrix<T>* pc_mat,
NumericVector<T>& sol,
NumericVector<T>& rhs,
const double tol,
const unsigned int n_iter)
{
if (pc_mat)
return this->solve(mat, *pc_mat, sol, rhs, tol, n_iter);
else
return this->solve(mat, sol, rhs, tol, n_iter);
}
} // namespace libMesh
#endif // #ifdef __solver_h__
|