/usr/include/rheolef/mixed_solver.h is in librheolef-dev 6.5-1+b1.
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 | # ifndef _RHEO_MIXED_SOLVER_H
# define _RHEO_MIXED_SOLVER_H
/// This file is part of Rheolef.
///
/// Copyright (C) 2000-2009 Pierre Saramito <Pierre.Saramito@imag.fr>
///
/// Rheolef is free software; you can redistribute it and/or modify
/// it under the terms of the GNU General Public License as published by
/// the Free Software Foundation; either version 2 of the License, or
/// (at your option) any later version.
///
/// Rheolef 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 General Public License for more details.
///
/// You should have received a copy of the GNU General Public License
/// along with Rheolef; if not, write to the Free Software
/// Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
///
/// =========================================================================
#include "puzawa.h"
#include "pcg.h"
#include "pminres.h"
namespace rheolef {
template <class Matrix, class Vector, class Solver>
struct abtbc_schur_complement { // S = B*inv(A)*B^T + C
Solver iA;
Matrix B, C;
abtbc_schur_complement (const Solver& iA1,const Matrix& B1, const Matrix& C1) : iA(iA1), B(B1), C(C1) {}
Vector operator* (const Vector& p) const {
Vector v = iA.solve(B.trans_mult(p));
return B*v + C*p;
}
};
template <class Matrix, class Vector, class Solver>
struct abtb_schur_complement { // S = B*inv(A)*B^T
Solver iA;
Matrix B;
abtb_schur_complement (const Solver& iA1,const Matrix& B1) : iA(iA1), B(B1) {}
Vector operator* (const Vector& p) const {
Vector v = iA.solve(B.trans_mult(p));
return B*v;
}
};
template <class Matrix, class Vector, class Solver, class Preconditioner, class Size, class Real>
int puzawa_abtbc (const Matrix& A, const Matrix& B, const Matrix& C, Vector& u, Vector& p,
const Vector& Mf, const Vector& Mg, const Preconditioner& S1,
const Solver& inner_solver_A, Size& max_iter, Real& tol, const Float& rho,
odiststream *p_derr = 0, std::string label = "puzawa_abtbc")
{
abtbc_schur_complement<Matrix,Vector,Solver> S (inner_solver_A, B, C);
Vector v = inner_solver_A.solve (Mf);
Vector Mh = B*v - Mg;
int status = puzawa(S, p, Mh, S1, max_iter, tol, rho, p_derr, label);
u = inner_solver_A.solve (Mf - B.trans_mult(p));
return status;
}
template <class Matrix, class Vector, class Solver, class Preconditioner, class Size, class Real>
int puzawa_abtb (const Matrix& A, const Matrix& B, Vector& u, Vector& p,
const Vector& Mf, const Vector& Mg, const Preconditioner& S1,
const Solver& inner_solver_A, Size& max_iter, Real& tol, const Float& rho,
odiststream *p_derr = 0, std::string label = "puzawa_abtb")
{
abtb_schur_complement<Matrix,Vector,Solver> S (inner_solver_A, B);
Vector v = inner_solver_A.solve (Mf);
Vector Mh = B*v - Mg;
int status = puzawa(S, p, Mh, S1, max_iter, tol, rho, p_derr, label);
u = inner_solver_A.solve (Mf - B.trans_mult(p));
return status;
}
/*Class:mixed_solver
NAME: @code{pcg_abtb}, @code{pcg_abtbc}, @code{pminres_abtb}, @code{pminres_abtbc} -- solvers for mixed linear problems
@findex pcg
@findex pminres
@findex pcg\_abtb
@findex pcg\_abtbc
@findex pminres\_abtb
@findex pminres\_abtbc
@cindex mixed linear problem
@cindex conjugate gradien algorithm
@cindex finite element method
@cindex stabilized mixed finite element method
@cindex Stokes problem
@cindex incompresible elasticity
SYNOPSIS:
@example
template <class Matrix, class Vector, class Solver, class Preconditioner, class Size, class Real>
int pcg_abtb (const Matrix& A, const Matrix& B, Vector& u, Vector& p,
const Vector& Mf, const Vector& Mg, const Preconditioner& S1,
const Solver& inner_solver_A, Size& max_iter, Real& tol,
odiststream *p_derr = 0, std::string label = "pcg_abtb");
template <class Matrix, class Vector, class Solver, class Preconditioner, class Size, class Real>
int pcg_abtbc (const Matrix& A, const Matrix& B, const Matrix& C, Vector& u, Vector& p,
const Vector& Mf, const Vector& Mg, const Preconditioner& S1,
const Solver& inner_solver_A, Size& max_iter, Real& tol,
odiststream *p_derr = 0, std::string label = "pcg_abtbc");
@end example
The synopsis is the same with the pminres algorithm.
EXAMPLES:
See the user's manual for practical examples for the nearly incompressible
elasticity, the Stokes and the Navier-Stokes problems.
DESCRIPTION:
@noindent
Preconditioned conjugate gradient algorithm on the pressure p applied to
the stabilized stokes problem:
@example
[ A B^T ] [ u ] [ Mf ]
[ ] [ ] = [ ]
[ B -C ] [ p ] [ Mg ]
@end example
where A is symmetric positive definite and C is symmetric positive
and semi-definite.
Such mixed linear problems appears for instance with the discretization
of Stokes problems with stabilized P1-P1 element, or with nearly
incompressible elasticity.
Formaly u = inv(A)*(Mf - B^T*p) and the reduced system writes for
all non-singular matrix S1:
@example
inv(S1)*(B*inv(A)*B^T)*p = inv(S1)*(B*inv(A)*Mf - Mg)
@end example
Uzawa or conjugate gradient algorithms are considered on the
reduced problem.
Here, S1 is some preconditioner for the Schur complement S=B*inv(A)*B^T.
Both direct or iterative solvers for S1*q = t are supported.
Application of inv(A) is performed via a call to a solver
for systems such as A*v = b.
This last system may be solved either by direct or iterative algorithms,
thus, a general matrix solver class is submitted to the algorithm.
For most applications, such as the Stokes problem,
the mass matrix for the p variable is a good S1 preconditioner
for the Schur complement.
The stoping criteria is expressed using the S1 matrix, i.e. in L2 norm
when this choice is considered.
It is scaled by the L2 norm of the right-hand side of the reduced system,
also in S1 norm.
AUTHOR:
| Pierre.Saramito@imag.fr
LJK-IMAG, 38041 Grenoble cedex 9, France
DATE: 15 april 2009
METHODS: @mixed_solver
End:
*/
template <class Matrix, class Vector, class VectorExpr1, class VectorExpr2,
class Solver, class Preconditioner, class Size, class Real>
int pcg_abtbc (const Matrix& A, const Matrix& B, const Matrix& C, Vector& u, Vector& p,
const VectorExpr1& Mf, const VectorExpr2& Mg, const Preconditioner& S1,
const Solver& inner_solver_A, Size& max_iter, Real& tol,
odiststream *p_derr = 0, std::string label = "pcg_abtbc")
{
abtbc_schur_complement<Matrix,Vector,Solver> S (inner_solver_A, B, C);
Vector v = inner_solver_A.solve (Mf);
Vector Mh = B*v - Mg;
int status = pcg(S, p, Mh, S1, max_iter, tol, p_derr, label);
u = inner_solver_A.solve (Mf - B.trans_mult(p));
return status;
}
template <class Matrix, class Vector, class VectorExpr1, class VectorExpr2,
class Solver, class Preconditioner, class Size, class Real>
int pcg_abtb (const Matrix& A, const Matrix& B, Vector& u, Vector& p,
const VectorExpr1& Mf, const VectorExpr2& Mg, const Preconditioner& S1,
const Solver& inner_solver_A, Size& max_iter, Real& tol,
odiststream *p_derr = 0, std::string label = "pcg_abtb")
{
abtb_schur_complement<Matrix,Vector,Solver> S (inner_solver_A, B);
Vector v = inner_solver_A.solve (Mf);
Vector Mh = B*v - Mg;
int status = pcg(S, p, Mh, S1, max_iter, tol, p_derr, label);
u = inner_solver_A.solve (Mf - B.trans_mult(p));
return status;
}
template <class Matrix, class Vector, class Solver, class Preconditioner, class Size, class Real>
int pminres_abtbc (const Matrix& A, const Matrix& B, const Matrix& C, Vector& u, Vector& p,
const Vector& Mf, const Vector& Mg, const Preconditioner& S1,
const Solver& inner_solver_A, Size& max_iter, Real& tol,
odiststream *p_derr = 0, std::string label = "pminres_abtbc")
{
abtbc_schur_complement<Matrix,Vector,Solver> S (inner_solver_A, B, C);
Vector v = inner_solver_A.solve (Mf);
Vector Mh = B*v - Mg;
int status = pminres(S, p, Mh, S1, max_iter, tol, p_derr, label);
u = inner_solver_A.solve (Mf - B.trans_mult(p));
return status;
}
template <class Matrix, class Vector, class Solver, class Preconditioner, class Size, class Real>
int pminres_abtb (const Matrix& A, const Matrix& B, Vector& u, Vector& p,
const Vector& Mf, const Vector& Mg, const Preconditioner& S1,
const Solver& inner_solver_A, Size& max_iter, Real& tol,
odiststream *p_derr = 0, std::string label = "pminres_abtb")
{
abtb_schur_complement<Matrix,Vector,Solver> S (inner_solver_A, B);
Vector v = inner_solver_A.solve (Mf);
Vector Mh = B*v - Mg;
int status = pminres(S, p, Mh, S1, max_iter, tol, p_derr, label);
u = inner_solver_A.solve (Mf - B.trans_mult(p));
return status;
}
}// namespace rheolef
# endif // _RHEO_MIXED_SOLVER_H
|