This file is indexed.

/usr/include/rheolef/solver.h is in librheolef-dev 6.7-6.

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
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
#ifndef _RHEOLEF_SOLVER_H
#define _RHEOLEF_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
///
/// =========================================================================
// direct & iterative solver interface
//
#include "rheolef/csr.h"
#include "rheolef/solver_option_type.h"

namespace rheolef {

// =======================================================================
// solver_abstract_rep: an abstract interface for solvers
// =======================================================================
// forward declaration:
template <class T, class M> class solver_basic;

template <class T, class M>
class solver_abstract_rep {
public:
  typedef typename csr<T,M>::size_type                 size_type;
  struct determinant_type;
  solver_abstract_rep (const solver_option_type& opt) : _opt(opt) {}
  virtual ~solver_abstract_rep () {}
  const solver_option_type& option() const { return _opt; }
  virtual void update_values (const csr<T,M>& a) = 0;
  virtual vec<T,M> trans_solve (const vec<T,M>& b) const = 0;
  virtual vec<T,M> solve       (const vec<T,M>& b) const = 0;
  virtual void set_preconditionner (const solver_basic<T,M>&);
  virtual determinant_type det() const;
// data:
  mutable solver_option_type _opt;
};
// det = mantissa*2^(exp2)
template <class T, class M>
struct solver_abstract_rep<T,M>::determinant_type {
  T       mantissa, exponant, base;
  determinant_type(): mantissa(0), exponant(0), base(10) {} 
};
// =======================================================================
// solver_rep: a pointer to the abstract representation
// =======================================================================
template <class T, class M>
class solver_rep {
public:
  typedef typename solver_abstract_rep<T,M>::size_type size_type;
  typedef typename solver_abstract_rep<T,M>::determinant_type determinant_type;
  explicit solver_rep   ();
  explicit solver_rep   (const csr<T,M>& a, const solver_option_type& opt = solver_option_type());
  void set_preconditionner (const solver_basic<T,M>&);
  determinant_type det() const;
  void build_eye ();
  void build_lu  (const csr<T,M>& a, const solver_option_type& opt);
  void build_ilu (const csr<T,M>& a, const solver_option_type& opt);
  ~solver_rep ();
  const solver_option_type& option() const;
  void update_values (const csr<T,M>& a);
  vec<T,M> trans_solve (const vec<T,M>& b) const;
  vec<T,M> solve       (const vec<T,M>& b) const;
protected:
  typedef solver_abstract_rep<T,M>  rep;
  rep* _ptr;
}; 
// =======================================================================
// the direct & iterative solver interface
// =======================================================================
/*Class:solver
NAME:   @code{solver} - direct or interative solver interface
@clindex solver
@clindex csr
@clindex vec
DESCRIPTION:
  @noindent
  The class implements a matrix factorization:
  LU factorization for an unsymmetric matrix and
  Choleski fatorisation for a symmetric one.

  Let @var{a} be a square invertible matrix in 
  @code{csr} format (@pxref{csr class}).
  @example
     	csr<Float> a;
  @end example
  @noindent
  We get the factorization by:
  @example
     	solver sa (a);
  @end example
  @noindent
  Each call to the direct solver for a*x = b writes either:
  @example
     	vec<Float> x = sa.solve(b);
  @end example
  When the matrix is modified in a computation loop but
  conserves its sparsity pattern, an efficient re-factorization
  writes:
  @example
     	sa.update_values (new_a);
     	x = sa.solve(b);
  @end example

AUTOMATIC CHOICE AND CUSTOMIZATION:
@cindex direct solver
@cindex iterative solver
  The symmetry of the matrix is tested via the a.is_symmetric() property
  (@pxref{csr class}) while the choice between direct or iterative solver
  is switched by default from the a.pattern_dimension() value. When the pattern
  is 3D, an iterative method is faster and less memory consuming.
  Otherwhise, for 1D or 2D problems, the direct method is prefered.
 
  These default choices can be supersetted by using explicit options:
  @example
	solver_option_type opt;
	opt.iterative = true;
     	solver sa (a, opt);
  @end example
  When using an iterative method, the sa.solve(b) call the conjugate gradient
  when the matrix is symmetric, or the generalized minimum residual
  algorithm when the matrix is unsymmetric.

COMPUTATION OF THE DETERMINANT:
@cindex determinant
  When using a direct method, the determinant of the @code{a} matrix
  can be computed as:
  @example
	solver_option_type opt;
	opt.iterative = false;
     	solver sa (a, opt);
        cout << sa.det().mantissa << "*2^" << sa.det().exponant << endl;
  @end example
  The @code{sa.det()} method returns an object of type @code{solver::determinant_type}
  that contains a mantissa and an exponent in base 2.
  This feature is usefull e.g. when tracking a change of sign in the determinant
  of a matrix.
 
PRECONDITIONNERS FOR ITERATIVE SOLVER:
@cindex preconditionner
  When using an iterative method, the default is to do no preconditionning.
  A suitable preconditionner can be supplied via:
  @example
	solver_option_type opt;
	opt.iterative = true;
     	solver sa (a, opt);
        sa.set_preconditionner (pa);
     	x = sa.solve(b);
  @end example
  The supplied @code{pa} variable is also of type @code{solver}.
  A library of commonly used preconditionners is still in development.

AUTHOR: Pierre.Saramito@imag.fr
DATE:   4 march 2011
End:
*/
//<verbatim:
template <class T, class M = rheo_default_memory_model>
class solver_basic : public smart_pointer<solver_rep<T,M> > {
public:
// typedefs:

  typedef solver_rep<T,M>                 rep;
  typedef smart_pointer<rep>              base;
  typedef typename rep::size_type         size_type;
  typedef typename rep::determinant_type  determinant_type;

// allocator:

  solver_basic ();
  explicit solver_basic (const csr<T,M>& a, const solver_option_type& opt = solver_option_type());
  const solver_option_type& option() const;
  void update_values (const csr<T,M>& a);
  void set_preconditionner (const solver_basic<T,M>&);

// accessors:

  vec<T,M> trans_solve (const vec<T,M>& b) const;
  vec<T,M> solve       (const vec<T,M>& b) const;
  determinant_type det() const;
};
// factorizations:
template <class T, class M>
solver_basic<T,M> ldlt(const csr<T,M>& a, const solver_option_type& opt = solver_option_type());
template <class T, class M>
solver_basic<T,M> lu  (const csr<T,M>& a, const solver_option_type& opt = solver_option_type());

// preconditionners:
template <class T, class M = rheo_default_memory_model>
solver_basic<T,M> eye_basic();
inline solver_basic<Float> eye() { return eye_basic<Float>(); }
template <class T, class M>
solver_basic<T,M> ic0 (const csr<T,M>& a, const solver_option_type& opt = solver_option_type());
template <class T, class M>
solver_basic<T,M> ilu0(const csr<T,M>& a, const solver_option_type& opt = solver_option_type());
template <class T, class M>
solver_basic<T,M> ldlt_seq(const csr<T,M>& a, const solver_option_type& opt = solver_option_type());
template <class T, class M>
solver_basic<T,M> lu_seq  (const csr<T,M>& a, const solver_option_type& opt = solver_option_type());

typedef solver_basic<Float> solver;
//>verbatim:

// =======================================================================
// solver_basic: inlined
// =======================================================================
template <class T, class M>
inline
solver_basic<T,M>::solver_basic ()
 : base (new_macro(rep))
{
}
template <class T, class M>
inline
solver_basic<T,M>::solver_basic (const csr<T,M>& a, const solver_option_type& opt)
 : base (new_macro(rep(a,opt)))
{
}
template <class T, class M>
inline
const solver_option_type&
solver_basic<T,M>::option() const
{
  return base::data().option();
}
template <class T, class M>
inline
solver_basic<T,M>
lu (const csr<T,M>& a, const solver_option_type& opt)
{
  solver_basic<T,M> sa;
  sa.data().build_lu (a, opt);
  return sa;
}
template <class T, class M>
inline
solver_basic<T,M>
ldlt (const csr<T,M>& a, const solver_option_type& opt)
{
  return lu(a,opt); // symmetry flag is indicated in the matrix a
}
template <class T, class M>
inline
solver_basic<T,M>
eye_basic ()
{
  solver_basic<T,M> sa;
  sa.data().build_eye ();
  return sa;
}
template <class T, class M>
inline
solver_basic<T,M>
ilu0 (const csr<T,M>& a, const solver_option_type& opt)
{
  solver_basic<T,M> sa;
  sa.data().build_ilu (a, opt);
  return sa;
}
template <class T, class M>
inline
solver_basic<T,M>
ic0 (const csr<T,M>& a, const solver_option_type& opt)
{
  return ilu0(a,opt); // symmetry flag is indicated in the matrix a
}
template <class T, class M>
inline
solver_basic<T,M>
lu_seq (const csr<T,M>& a, const solver_option_type& opt)
{
  solver_option_type opt1 = opt;
  opt1.force_seq = true;
  return lu (a,opt1);
}
template <class T, class M>
inline
solver_basic<T,M>
ldlt_seq (const csr<T,M>& a, const solver_option_type& opt)
{
  return lu_seq(a,opt); // symmetry flag is indicated in the matrix a
}
template <class T, class M>
inline
void
solver_basic<T,M>::update_values (const csr<T,M>& a)
{
  base::data().update_values (a);
}
template <class T, class M>
inline
void
solver_basic<T,M>::set_preconditionner (const solver_basic<T,M>& sa)
{
  base::data().set_preconditionner (sa);
}
template <class T, class M>
inline
typename solver_basic<T,M>::determinant_type
solver_basic<T,M>::det() const
{
  return base::data().det();
}
template <class T, class M>
inline
vec<T,M>
solver_basic<T,M>::solve (const vec<T,M>& b) const
{
  return base::data().solve (b);
}
template <class T, class M>
inline
vec<T,M>
solver_basic<T,M>::trans_solve (const vec<T,M>& b) const
{
  return base::data().trans_solve (b);
}

} // namespace rheolef
#endif // _RHEOLEF_SOLVER_H