1%%%%%%%%%%%%%%%%%%%
2% XLiFE++ is an extended library of finite elements written in C++
3%     Copyright (C) 2014  Lunéville, Eric; Kielbasiewicz, Nicolas; Lafranche, Yvon; Nguyen, Manh-Ha; Chambeyron, Colin
4%
5%     This program is free software: you can redistribute it and/or modify
6%     it under the terms of the GNU General Public License as published by
7%     the Free Software Foundation, either version 3 of the License, or
8%     (at your option) any later version.
9%     This program is distributed in the hope that it will be useful,
10%     but WITHOUT ANY WARRANTY; without even the implied warranty of
11%     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
12%     GNU General Public License for more details.
13%     You should have received a copy of the GNU General Public License
14%     along with this program.  If not, see <http://www.gnu.org/licenses/>.
15%%%%%%%%%%%%%%%%%%%
16
17\section{The {\classtitle TermMatrix} class}
18
19The {\class TermMatrix} class carries numerical representation of a bilinear form  and more generally any matrix attached to a pair of unknowns. It follows a similar organization to {\class TermVector} class but manages row and column information:
20\vspace{.1cm}
21\begin{lstlisting}[deletekeywords={[3] map}]
22class TermMatrix : public Term
23{protected :
24BilinearForm bilinForm_;                   //bilinear form if exists
25std::map<uvPair, SuTermMatrix*> suTerms_;  //list of SuTermMatrix
26MatrixEntry* scalar_entries_p;             //scalar entries representation
27std::vector<DofComponent> cdofs_c;         //column scalar dofs
28std::vector<DofComponent> cdofs_r;         //row scalar dofs
29SetOfConstraints* constraints_u_p;         //constraints to apply to row
30SetOfConstraints* constraints_v_p;         //constraints to apply to column
31MatrixEntry* rhs_matrix_p;      //matrix used by right hand side correction
32}
33\end{lstlisting}
34\vspace{.3cm}
35The bilinear form is a multiple unknown pairs form, may be reduced to a single unknown pair form (see {\lib form} library). This bilinear form is void if {\class TermMatrix} is not explicitly construct from a bilinear form (linear combination, nodal construction, ...). In particular, this bilinear form is not updated when {\class TermMatrix} objects are combined. \\
36The map \verb? suTerms_? carries each matrix block related to a pair of unknowns. When there is only one element in map, it means that the {\class TermMatrix} is actually a single unknown pair term matrix. \\
37As block representation may be too restrictive in certain cases, it is possible to create its non block representation in \emph{scalar\_entries\_p} ({\class MatrixEntry} pointer). The {\class MatrixEntry} object created has to be of real scalar or complex scalar type. In that case, non block row and non block column numbering, consisting in vectors of {\class DofComponent} (Unknown pointer, dof number, component number), are also defined (\verb?cdofs_r? and \verb?cdofs_c?).\\
38
39Essential conditions (conditions on space functions) induce special treatments in computation of {\class TermMatrix}. It is the reason why they are handled with
40{\class SetOfConstraints} objects which are the algebraic representations of essential conditions (see {\lib essentialConditions} for details). {\class SetOfConstraints} on row may differ from  {\class SetOfConstraints} on column (it is often the same). In case of non homogeneous essential condition, the part of matrix that is eliminated contributes to the right hand side; this part is stored into \verb?rhs_matrix_p? {\class MatrixEntry} pointer. Obviously, the pointers \verb?constraints_u_p, constraints_v_p, rhs_matrix_p? are null pointer whenever no essential conditions are imposed.\\
41
42The {\class TermMatrix} class provides useful constructors and assign operators designed for end users :
43\vspace{.1cm}
44\begin{lstlisting}[]{}
45TermMatrix(const string_t & na="?");
46
47//constructors with no essential boundary conditions and options
48TermMatrix(const BilinearForm&, const string_t& na="");
49TermMatrix(const BilinearForm&, TermOption, const string_t& na="");
50TermMatrix(const BilinearForm&, TermOption, TermOption, const string_t& na="");
51...
52//constructor with one essential boundary condition (on u and v)
53TermMatrix(const BilinearForm&, const EssentialConditions&, const string_t& na="");
54TermMatrix(const BilinearForm&, const EssentialConditions&, TermOption,
55           const string_t& na="");
56...
57//constructor with one essential boundary condition and explicit reduction method
58TermMatrix(const BilinearForm&, const EssentialConditions&, const ReductionMethod&,
59           const string_t& na="");
60TermMatrix(const BilinearForm&, const EssentialConditions&, const ReductionMethod&,
61           TermOption, const string_t& na="");
62...
63//constructor with  two essential boundary conditions (on u and v)
64TermMatrix(const BilinearForm&, const EssentialConditions&,
65           const EssentialConditions&, const string_t& na="");
66TermMatrix(const BilinearForm&, const EssentialConditions&,
67           const EssentialConditions&, TermOption, const string_t& na="");
68...
69//constructor with two essential boundary conditions and explicit reduction method
70TermMatrix(const BilinearForm&, const EssentialConditions&, const EssentialConditions&,
71           const ReductionMethod&, const string_t& na="");
72TermMatrix(const BilinearForm&, const EssentialConditions&, const EssentialConditions&,
73           const ReductionMethod&, TermOption, const string_t& na="");
74...
75//effective constructor
76void build(const BilinearForm& blf, const EssentialConditions* ecu,
77           const EssentialConditions* ecv, const std::vector<TermOption>&,
78           const string_t& na);
79
80//copy constructors
81TermMatrix(const TermMatrix &, const string_t& na="");
82// diagonal TermMatrix from a TermVector or a vector
83TermMatrix(const TermVector &, const string_t& na="");
84template<typename T>
85TermMatrix(const Unknown&, const GeomDomain&, const Vector<T> &,
86           const string_t& na="");
87//TermMatrix from a SuTermMatrix
88TermMatrix(const SuTermMatrix &, const string_t& na=""); //full copy
89TermMatrix(SuTermMatrix *, const string_t& na="");       //pointer copy
90//TermMatrix of the form opker(xi,yj)
91TermMatrix(const Unknown&, const GeomDomain&, const Unknown&,
92           const GeomDomain&, const OperatorOnKernel&,
93           const string_t& na="");
94// TermMatrix from a linear combination
95TermMatrix(const LcTerm&,const string_t& na="");
96
97//assign operator
98TermMatrix& operator=(const TermMatrix&);
99TermMatrix& operator=(const LcTerm&);
100\end{lstlisting}
101\vspace{.2cm}
102The {\class ReductionMethod} class allows the user to specify the reduction method he wants, one of the enumeration {\class ReductionMethodType}: {\it \_pseudoReduction, \_realReduction, \_penalizationReduction, \_dualReduction}. It is also possible to provide an additional real/complex parameter (say alpha) : the diagonal coefficient of the pseudo-eliminated part or the penalization coefficient.
103\vspace{.1cm}
104\begin{lstlisting}[]{}
105ReductionMethodType method;  //type of reduction method
106complex_t alpha;            //diagonal or penalization coefficient
107ReductionMethod(ReductionMethodType= _noReduction, const complex_t& = 1.);
108void print(std::ostream& out) const;
109friend :ostream& operator << (ostream&, const ReductionMethod&);
110\end{lstlisting}
111\vspace{.2cm}
112\begin{remark}
113When reduction method is not specified, the pseudo reduction with 1 as diagonal coefficient is used. {\bf Only the pseudo reduction method is available!}
114\end{remark}
115
116\begin{remark}
117When TermMatrix is constructed, it is automatically computed except if it has set to be not computed using a the {\class TermOption} \textit{\_notCompute} in its construction.
118\end{remark}
119
120Developers may used other construction/destruction stuff:
121\vspace{.1cm}
122\begin{lstlisting}[]{}
123void initFromBlf(const BilinearForm &, const String& na ="",
124                 bool noass = false);
125void initTermVector(TermVector&, ValueType, bool col=true);
126template<typename T>
127void initTermVectorT(TermVector&,const T&, bool col=true);
128void insert(const SuTermMatrix &);
129void insert(SuTermMatrix *);
130void copy(const TermMatrix&);
131void clear();
132virtual ~TermMatrix();
133
134\end{lstlisting}
135\vspace{.2cm}
136This class provides a lot of accessors and properties (shorcuts to other class accessors):
137\vspace{.1cm}
138\begin{lstlisting}[deletekeywords={[3] set}]
139bool isSingleUnknown() const;
140std::set<const Unknown*> rowUnknowns() const;
141std::set<const Unknown*> colUnknowns() const;
142Number nbTerms() const;
143ValueType valueType() const;
144AccessType storageAccess() const;
145void setStorage(StorageType, AccessType);
146StorageType storageType() const
147bool hasConstraints() const;
148FactorizationType factorization() const;
149SymType globalSymmetry() const;
150TermMatrix operator()(const Unknown&,const Unknown&) const;
151
152//stuff reserved to developers
153cit_mustm begin() const;
154it_mustm begin();
155cit_mustm end() const;
156it_mustm end();
157std::map<uvPair, SuTermMatrix*>& suTerms();
158MatrixEntry*& scalar_entries();
159const MatrixEntry* scalar_entries() const;
160MatrixEntry*& rhs_matrix();
161const MatrixEntry* rhs_matrix() const;
162const std::vector<DofComponent>& cdofsr() const;
163const std::vector<DofComponent>& cdofsc() const;
164std::vector<DofComponent>& cdofsr();
165std::vector<DofComponent>& cdofsc();
166SuTermMatrix& subMatrix(const Unknown*,const Unknown*);
167const SuTermMatrix& subMatrix(const Unknown*,const Unknown*) const;
168SuTermMatrix* subMatrix_p(const Unknown*,const Unknown*);
169const SuTermMatrix* subMatrix_p(const Unknown*,const Unknown*) const;
170
171TermVector& getRowCol(number_t, AccessType, TermVector&) const;
172TermVector row(number_t) const;
173TermVector column(number_t) const;
174template <typename T>
175LargeMatrix<T>& getLargeMatrix(StrucType =_undefStrucType) const;
176template <typename T>
177HMatrix<T,FeDof>& getHMatrix(StrucType =_undefStrucType) const;
178
179\end{lstlisting}
180\begin{remark}
181The access {\cmd operator()} returns a {\class SuTermMatrix} as a {\class TermMatrix} (obviously a single unknown one). Note that it is a full copy. To avoid copy, use {\cmd subMatrix} member functions.
182\end{remark}
183
184\vspace{.2cm}
185\verb?it_mustm? and \verb?cit_mustm? are iterator aliases:
186\vspace{.1cm}
187\begin{lstlisting}[deletekeywords={[3] map}]
188typedef std::map<uvPair, SuTermMatrix*>::iterator it_mustm;
189typedef std::map<uvPair, SuTermMatrix*>::const_iterator cit_mustm;
190\end{lstlisting}
191\vspace{.2cm}
192It is also possible to get row or column (index starts from 1) of the TermMatrix if it is a single unknown one and to access to the LargeMatrix or the HMatrix object that really stores the matrix values:
193\vspace{.1cm}
194\begin{lstlisting}[deletekeywords={[3] set}]
195TermVector& getRowCol(number_t, AccessType, TermVector&) const;
196TermVector row(number_t) const;
197TermVector column(number_t) const;
198template <typename T>
199LargeMatrix<T>& getLargeMatrix(StrucType =_undefStrucType) const;
200template <typename T>
201HMatrix<T,FeDof>& getHMatrix(StrucType =_undefStrucType) const;
202\end{lstlisting}
203\textdbend \ \ The {\cmd setStorage} function allows the user to impose the matrix storage.\\
204
205The class proposes to users interfaces to computation algorithms (FE computation and algebraic operations):
206\vspace{.1cm}
207\begin{lstlisting}[]{}
208void compute();
209void compute(const LcTerm&,const String& na="");
210
211template<typename T> TermMatrix& operator*=(const T&);
212template<typename T> TermMatrix& operator/=(const T&);
213TermMatrix& operator+=(const TermMatrix&);
214TermMatrix& operator-=(const TermMatrix&);
215friend TermVector& multMatrixVector(const TermMatrix&,const TermVector&, TermVector&);
216friend TermVector& multVectorMatrix(const TermMatrix&, const TermVector&, TermVector&);
217friend TermVector operator*(const TermMatrix&,const TermVector&);
218friend TermVector operator*(const TermVector&, const TermMatrix&);
219friend TermMatrix operator*(const TermMatrix&, const TermMatrix&);
220\end{lstlisting}
221\vspace{.2cm}
222
223The computation of  {\class TermMatrix} is a quite intricate process in case of essential conditions to apply to bilinear form. The implementation of \verb?compute()? looks like (some details are omitted):
224\vspace{.1cm}
225\begin{lstlisting}[]{}
226void TermMatrix::compute()
227{
228  //compute suterms
229  for(it_mustm it = suTerms_.begin(); it != suTerms_.end(); it++)
230      it->second->compute();
231  if(constraints_u_p==0 && constraints_v_p==0) {computed() = true; return;}
232  //deal with essential conditions
233  bool global=false;
234  if(constraints_u_p!=0) global=constraints_u_p->isGlobal();
235  if(constraints_v_p!=0 && !global) global=constraints_u_p->isGlobal();
236  if(global) toGlobal(_noStorage, _noAccess, _noSymmetry, false);
237  else toScalar();
238  switch(computingInfo_.elimination)
239    {
240      case _noElimination             :  break;
241      case _pseudoElimination         :  pseudoElimination(); break;
242    }
243  computed() = true;
244}
245\end{lstlisting}
246\vspace{.2cm}
247The principle of pseudo-elimination consists in combining some rows and columns, and "eliminating" some rows and columns (in fact set to 0 except diagonal coefficient). Two cases are to be considered:
248\begin{itemize}
249\item case of essential conditions which do not couple unknowns
250\item case of essential conditions which couples unknowns (global constraints)
251\end{itemize}
252In the first case, the pseudo-elimination process may be performed on each  matrix of {\class SuTermMatrix} involving unknowns concerned by essential conditions whereas in the second case, the process must be performed and the whole matrix of {\class TermMatrix}. If no such representation exists, we have to create it (\verb?toGlobal? function). See the \verb?pseudoElimination()? function.\\
253\\
254
255This class offers many interfaces to factorization tools and direct solvers :
256\vspace{.1cm}
257\begin{lstlisting}[]{}
258//developers stuff
259void initTermVector(TermVector&, ValueType, bool col=true);
260template<typename T>
261void initTermVectorT(TermVector&, const T&, bool col=true);
262
263void toScalar(bool keepEntries=false);
264void toGlobal(StorageType, AccessType, SymType=_noSymmetry, bool keepSuTerms=false);
265void toSkyline();
266void mergeNumbering(std::map<const Unknown*, std::list<SuTermMatrix* > >&, std::map<SuTermMatrix*, std::vector<Number > >&, std::vector<DofComponent>&, AccessType);
267void addIndices(std::vector<std::set<Number> >&, MatrixStorage*, const std::vector<Number>&, const std::vector<Number>&);
268void mergeBlocks();
269friend TermVector prepareLinearSystem(TermMatrix&, TermVector&, MatrixEntry*&, VectorEntry*&, bool toScal=false);
270void pseudoElimination();
271
272//direct solver member functions
273void ldltSolve(const TermVector&, TermVector&);
274void ldlstarSolve(const TermVector&,TermVector&);
275void luSolve(const TermVector&,TermVector&);
276void umfpackSolve(const TermVector&,TermVector&);
277void sorSolve(const TermVector& V, TermVector& R, const Real w, SorSolverType sType );
278void sorUpperSolve(const TermVector& V, TermVector& R, const Real w );
279void sorDiagonalMatrixVector(const TermVector& V, TermVector& R, const Real w );
280void sorLowerSolve(const TermVector& V, TermVector& R, const Real w );
281void sorDiagonalSolve(const TermVector& V, TermVector& R, const Real w );
282
283\end{lstlisting}
284\vspace{.2cm}
285End users have to use the following external functions:
286\vspace{.1cm}
287\begin{lstlisting}[]{}
288//factorization tool
289void ldltFactorize(TermMatrix&, TermMatrix&);
290void ldlstarFactorize(TermMatrix&, TermMatrix&);
291void luFactorize(TermMatrix&, TermMatrix&);
292void umfpackFactorize(TermMatrix&, TermMatrix&);
293void factorize(TermMatrix&, TermMatrix&, FactorizationType ft=_noFactorization);
294
295//direct solver (one right hand side)
296TermVector factSolve(TermMatrix&, const TermVector&);
297TermVector ldltSolve(TermMatrix&, const TermVector&);
298TermVector ldlstarSolve(TermMatrix&, const TermVector&);
299TermVector luSolve(TermMatrix&, const TermVector&);
300TermVector gaussSolve(TermMatrix&, const TermVector&, bool keepA = false);
301TermVector umfpackSolve(TermMatrix&, const TermVector&, bool keepA = false);
302TermVector magmaSolve(TermMatrix& A, const TermVector& B, bool useGPU=false,
303                      bool keepA=false);
304TermVector lapackSolve(TermMatrix& A, const TermVector& B, bool keepA=false);
305TermVector directSolve(TermMatrix&, const TermVector&, bool keepA = false);
306TermVector schurSolve(TermMatrix&, const TermVector&, const Unknown&, const Unknown&,
307                      bool keepA = false);
308
309//direct solver (few right hand sides)
310TermVectors factSolve(TermMatrix&, const vector<TermVector>&);
311TermVectors ldltSolve(TermMatrix&, const vector<TermVector>&, TermMatrix&);
312TermVectors ldlstarSolve(TermMatrix&, const vector<TermVector>&, TermMatrix&);
313TermVectors luSolve(TermMatrix&, const vector<TermVector>&, TermMatrix&);
314TermVectors gaussSolve(TermMatrix&, const vector<TermVector>&,bool keepA = false);
315TermVectors umfpackSolve(TermMatrix&, const vector<TermVector>&,bool keepA = false);
316TermVectors directSolve(TermMatrix&, const vector<TermVector>&,bool keepA = false);
317
318// direct solver (TermMatrix right hand side)
319TermMatrix factSolve(TermMatrix&, TermMatrix&);
320TermMatrix directSolve(TermMatrix&, TermMatrix&, KeepStatus=_keep);
321TermMatrix inverse(TermMatrix&);
322\end{lstlisting}
323\vspace{.2cm}
324{\cmd factorize()}, {\cmd factSolve()} and {\cmd directSolve} functions are general functions that determine automatically the adapted methods to apply.\\
325
326\begin{remark}
327	to use umfpackSolve, magmaSolve or lapackSolve, libraries umfpack, magma, lapack-blas have to be installed and activated. directSolve choose umfPack solver when matrix is sparse and umfpack available, and Lapack solver when matrix is dense and lapack available.	umfpack solver is always faster than \xlifepp solvers but lapack solver may be slower than the \xlifepp gauss solver if non optimized lapack-blas libraries are used
328\end{remark}\\
329
330\vspace{.2cm}
331Note that the {\cmd factSolve} and {\cmd directSolve} functions can use a {\class TermMatrix} as right hand side. In that case, they compute $\mathbb{A}^{-1}\mathbb{B}$ using a factorization of $\mathbb{A}$. The  {\cmd inverse()} function computes the inverse of a {\class TermMatrix} using the {\cmd directSolve} function with the identity matrix as right hand side. These functions provide some {\class TermMatrix} that are stored as column dense storage. So be care with memory usage!\\
332
333\begin{remark}
334Up to now, the usage of {\cmd factSolve}, {\cmd directSolve} functions with a {\class TermMatrix} as right hand side and {\cmd inverse()} is restricted to single unknown {\class TermMatrix}.
335\end{remark}
336\vspace{.2cm}
337
338In a same way, interfaces to various iterative solvers are provided :
339\vspace{.1cm}
340\begin{lstlisting}
341//! general iterative solver
342TermVector iterativeSolveGen(IterativeSolverType isType, TermMatrix& A,
343       TermVector& B, const TermVector& X0, const PreconditionerTerm& P,
344       const real_t tol, const number_t iterMax, const real_t omega,
345       const number_t krylovDim, const number_t verboseLevel);
346//! general front end for iterativeSolve
347TermVector iterativeSolve(TermMatrix& A, TermVector& B, const TermVector& X0,
348      const PreconditionerTerm& P, const std::vector<Parameter>& ps);
349
350//! user front-end for iterative solver
351inline TermVector iterativeSolve(TermMatrix& A, TermVector& B, const TermVector& X0, const PreconditionerTerm& P, 0 to 5 Parameter);
352inline TermVector iterativeSolve(TermMatrix& A, TermVector& B, const PreconditionerTerm& P, 0 to 6 Parameter);
353inline TermVector iterativeSolve(TermMatrix& A, TermVector& B, const TermVector& X0, 0 to 6 Parameter);
354inline TermVector iterativeSolve(TermMatrix& A, TermVector& B, 0 to 6 Parameter);
355
356//! user front-end for BiCG solver
357inline TermVector bicgSolve(TermMatrix& A, TermVector& B, const TermVector& X0, const PreconditionerTerm& P, 0 to 4 Parameter);
358inline TermVector bicgSolve(TermMatrix& A, TermVector& B, const PreconditionerTerm& P, 0 to 4 Parameter);
359inline TermVector bicgSolve(TermMatrix& A, TermVector& B, const TermVector& X0, 0 to 4 Parameter);
360inline TermVector bicgSolve(TermMatrix& A, TermVector& B, 0 to 4 Parameter);
361
362//! user front-end for BiCGStab solver
363inline TermVector bicgStabSolve(TermMatrix& A, TermVector& B, const TermVector& X0, const PreconditionerTerm& P, 0 to 4 Parameter);
364inline TermVector bicgStabSolve(TermMatrix& A, TermVector& B, const PreconditionerTerm& P, 0 to 4 Parameter);
365inline TermVector bicgStabSolve(TermMatrix& A, TermVector& B, const TermVector& X0, 0 to 4 Parameter);
366inline TermVector bicgStabSolve(TermMatrix& A, TermVector& B, 0 to 4 Parameter);
367
368//! user front-end for CG solver
369inline TermVector cgSolve(TermMatrix& A, TermVector& B, const TermVector& X0, const PreconditionerTerm& P, 0 to 4 Parameter);
370inline TermVector cgSolve(TermMatrix& A, TermVector& B, const PreconditionerTerm& P, 0 to 4 Parameter);
371inline TermVector cgSolve(TermMatrix& A, TermVector& B, const TermVector& X0, 0 to 4 Parameter);
372inline TermVector cgSolve(TermMatrix& A, TermVector& B, 0 to 4 Parameter);
373
374//! user front-end for CGS solver
375inline TermVector cgsSolve(TermMatrix& A, TermVector& B, const TermVector& X0, const PreconditionerTerm& P, 0 to 4 Parameter);
376inline TermVector cgsSolve(TermMatrix& A, TermVector& B, const PreconditionerTerm& P, 0 to 4 Parameter);
377inline TermVector cgsSolve(TermMatrix& A, TermVector& B, const TermVector& X0, 0 to 4 Parameter);
378inline TermVector cgsSolve(TermMatrix& A, TermVector& B, 0 to 4 Parameter);
379
380//! user front-end for GMRes solver
381inline TermVector gmresSolve(TermMatrix& A, TermVector& B, const TermVector& X0, const PreconditionerTerm& P, 0 to 4 Parameter);
382inline TermVector gmresSolve(TermMatrix& A, TermVector& B, const PreconditionerTerm& P, 0 to 4 Parameter);
383inline TermVector gmresSolve(TermMatrix& A, TermVector& B, const TermVector& X0, 0 to 4 Parameter);
384inline TermVector gmresSolve(TermMatrix& A, TermVector& B, 0 to 4 Parameter);
385
386//! user front-end for QMR solver
387inline TermVector qmrSolve(TermMatrix& A, TermVector& B, const TermVector& X0, const PreconditionerTerm& P, 0 to 4 Parameter);
388inline TermVector qmrSolve(TermMatrix& A, TermVector& B, const PreconditionerTerm& P, 0 to 4 Parameter);
389inline TermVector qmrSolve(TermMatrix& A, TermVector& B, const TermVector& X0, 0 to 4 Parameter);
390inline TermVector qmrSolve(TermMatrix& A, TermVector& B, 0 to 4 Parameter);
391\end{lstlisting}
392\vspace{.2cm}
393{\cmd iterativeSolve(\ldots)} are the front end of all other iterative methods.\\
394{\cmd iterativeSolveGen(\ldots)} are the front end of all iterative methods.\\
395
396To compute eigenvalues, internal solvers or interface to arpack solvers are available:
397
398\vspace{.1cm}
399\begin{lstlisting}[]{}
400//! general eigen solver
401EigenElements eigenSolveGen(TermMatrix* pA, TermMatrix* pB, number_t nev, string_t which, real_t tol, EigenComputationalMode eigCompMode,
402                            complex_t sigma, const char mode, EigenSolverType solver);
403//! general eigen solver
404EigenElements eigenInternGen(TermMatrix* pA, TermMatrix* pB, number_t nev, string_t which, real_t tol, EigenComputationalMode eigCompMode,
405                            complex_t sigma, bool isShift);
406
407//! general front-end for generalized eigen solver
408EigenElements eigenSolve(TermMatrix& A, TermMatrix& B, std::vector<Parameter> ps);
409
410//! user front-end for eigen solver
411EigenElements eigenSolve(TermMatrix& A, TermMatrix& B, 0 to 10 Parameter);
412EigenElements eigenSolve(TermMatrix& A, 0 to 10 Parameter);
413
414//! user front-end for internal eigen solver
415inline EigenElements eigenInternSolve(TermMatrix& A, TermMatrix& B, 0 to 9 Parameter);
416inline EigenElements eigenInternSolve(TermMatrix& A, 0 to 9 Parameter);
417
418//! user front-end for external Arpack eigen solver
419inline EigenElements arpackSolve(TermMatrix& A, TermMatrix& B, 0 to 9 Parameter);
420inline EigenElements arpackSolve(TermMatrix& A, 0 to 9 Parameter);
421\end{lstlisting}
422\vspace{.2cm}
423The {\class EigenElements} class handles the list of eigenvalues and the list of eigenvectors. See the {\class TermVector} section.\\
424
425\begin{remark}
426arpack solver are available only if the arpack lib is installed.
427\end{remark}
428
429\vspace{.3cm}
430When addressing integral representations, say
431$$r(x)=\int_{\Sigma}op(K)(x,y)*op(u)(y)\,dy$$
432one can be interested in the matrix with coefficients
433$$r_{ij}=\int_{\Sigma}op(K)(x_i,y)*op(w_j)(y)\,dy$$
434where $(x_i)_i$ are some points and $(w_j)_j$ some shape functions related to the unknown $u$.
435Special functions, named {\var integralFunction} will produce such matrices embedded in a {\class TermMatrix} object:
436\begin{lstlisting}[]{}
437TermMatrix integralRepresentation(const Unknown&, const GeomDomain&,
438                                  const LinearForm&);
439TermMatrix integralRepresentation(const GeomDomain&, const LinearForm&);
440TermMatrix integralRepresentation(const std::vector<Point>&, const LinearForm&);
441\end{lstlisting}
442When no domain is explicitly passed, one shadow {\class PointsDomain} object will be created from the list of points given. When no unknown is explicitly  passed, one (minimal) space and related shadow unknown will be created to represent the row unknown.\\
443These functions call the fundamental computation function:
444\begin{lstlisting}[]{}
445template<typename T, typename K>
446void SuTermMatrix::computeIR(const SuLinearForm&f, T*, K&, const std::vector<Point>& xs);
447\end{lstlisting}
448
449Finally some print and export functions are provided:
450\vspace{.1cm}
451\begin{lstlisting}[]{}
452//as member function
453void print(std::ostream&) const;
454void viewStorage(std::ostream&) const;
455void saveToFile(const String&, StorageType, bool=false) const;
456//as non member function
457void saveToFile(const TermMatrix&, const String&, StorageType,
458                bool enc = false);
459\end{lstlisting}
460\vspace{.2cm}
461
462\displayInfos{
463library=term,
464header=TermMatrix.hpp,
465implementation=TermMatrix.cpp,
466test=test\_TermMatrix.cpp,
467header dep={SuTermMatrix.hpp, termUtils.hpp, Term.hpp, form.h, config.h, utils.h}
468}
469