/usr/include/rheolef/solver_abtb.h is in librheolef-dev 6.5-1build1.
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 | #ifndef _SKIT_SOLVER_ABTB_H
#define _SKIT_SOLVER_ABTB_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 "rheolef/solver.h"
#include "rheolef/mixed_solver.h"
#include "rheolef/csr.h"
namespace rheolef {
/*Class:solver_abtb
NAME: @code{solver_abtb} -- direct or iterative solver iterface for mixed linear systems
@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
solver_abtb stokes (a,b,mp);
solver_abtb elasticity (a,b,c,mp);
@end example
DESCRIPTION:
@noindent
The @code{solver_abtb} class provides direct or iterative algorithms for some mixed problem:
@example
[ A B^T ] [ u ] [ Mf ]
[ ] [ ] = [ ]
[ B -C ] [ p ] [ Mg ]
@end example
where A is symmetric positive definite and C is symmetric positive.
By default, iterative algorithms are considered for tridimensional problems
and direct methods otherwise.
Such mixed linear problems appears for instance with the discretization
of Stokes problems.
The C matrix can be zero and then the corresponding argument can be omitted
when invoking the constructor.
Non-zero C matrix appears for of Stokes problems with stabilized P1-P1 element,
or for nearly incompressible elasticity problems.
DIRECT ALGORITHM:
When the kernel of @code{B^T} is not reduced to zero, then the pressure p is defined up to a constant
and the system is singular. In the case of iterative methods, this is not a problem.
But when using direct method, the system is then completed to impose a constraint on the pressure term
and the whole matrix is factored one time for all.
ITERATIVE ALGORITHM:
The preconditionned conjugate gradient algorithm is used, where the @code{mp} matrix
is used as preconditionner. See @pxref{mixed_solver algorithm}.
EXAMPLES:
See the user's manual for practical examples for the nearly incompressible
elasticity, the Stokes and the Navier-Stokes problems.
AUTHOR:
| Pierre.Saramito@imag.fr
LJK-IMAG, 38041 Grenoble cedex 9, France
DATE: 19 january 2012
METHODS: @mixed_solver
End:
*/
//<solver_abtb:
template <class T, class M = rheo_default_memory_model>
class solver_abtb_basic {
public:
// typedefs:
typedef typename csr<T,M>::size_type size_type;
// allocators:
solver_abtb_basic ();
solver_abtb_basic (const csr<T,M>& a, const csr<T,M>& b, const csr<T,M>& mp,
const solver_option_type& opt = solver_option_type());
solver_abtb_basic (const csr<T,M>& a, const csr<T,M>& b, const csr<T,M>& c, const csr<T,M>& mp,
const solver_option_type& opt = solver_option_type());
// accessors:
void solve (const vec<T,M>& f, const vec<T,M>& g, vec<T,M>& u, vec<T,M>& p) const;
protected:
// internal
void init();
// data:
mutable solver_option_type _opt;
csr<T,M> _a;
csr<T,M> _b;
csr<T,M> _c;
csr<T,M> _mp;
solver_basic<T,M> _sA;
solver_basic<T,M> _sa;
solver_basic<T,M> _smp;
bool _need_constraint;
};
typedef solver_abtb_basic<Float,rheo_default_memory_model> solver_abtb;
//>solver_abtb:
} // namespace rheolef
#endif // _SKIT_SOLVER_ABTB_H
|