source: Sophya/trunk/SophyaLib/TArray/sopemtx.cc@ 2615

Last change on this file since 2615 was 2615, checked in by cmv, 21 years ago

using namespace sophya enleve de machdefs.h, nouveau sopnamsp.h cmv 10/09/2004

File size: 25.7 KB
RevLine 
[772]1// C.Magneville, R.Ansari 03/2000
2
[2615]3#include "sopnamsp.h"
[772]4#include "machdefs.h"
5#include <stdio.h>
[2322]6#include <iostream>
[772]7#include <complex>
8#include <math.h>
9#include "sopemtx.h"
10#include "smathconst.h"
11
12////////////////////////////////////////////////////////////////
[939]13// ---------------------------------------------------------- //
14// La classe de gestion des lignes et colonnes d'une matrice //
15// ---------------------------------------------------------- //
[772]16////////////////////////////////////////////////////////////////
[926]17
[939]18/*!
19 \class SOPHYA::TMatrixRC
20 \ingroup TArray
21 Class for representing rows, columns and diagonal of a matrix.
22 \sa TMatrix TArray
23 */
24
[926]25//! Class of line, column or diagonal of a TMatrix
26/*!
27 A TMatrixRC represents a line, a column or the diagonal of a TMatrix
28 */
[804]29template <class T>
30class TMatrixRC {
31public:
[926]32 //! Define type of TMatrixRC
33 enum TRCKind {
34 TmatrixRow=0, //!< TMatrixRC ligne
35 TmatrixCol=1, //!< TMatrixRC column
36 TmatrixDiag=2 //!< TMatrixRC diagonal
37 };
[804]38 TMatrixRC();
39 TMatrixRC(TMatrix<T>& m, TRCKind kind, uint_4 index=0);
40 virtual ~TMatrixRC() {}
[772]41
[804]42 // Acces aux rangees et colonnes de matrices
43 static TMatrixRC<T> Row(TMatrix<T> & m, uint_4 r);
44 static TMatrixRC<T> Col(TMatrix<T> & m, uint_4 c);
45 static TMatrixRC<T> Diag(TMatrix<T> & m);
46
47 // ---- A virer $CHECK$ Reza 03/2000
48 // int_4 Next();
49 // int_4 Prev();
50 int_4 SetCol(int_4 c);
51 int_4 SetRow(int_4 r);
52 int_4 SetDiag();
53
54 static uint_4 Step(const TMatrix<T>& m, TRCKind rckind);
55 static T* Org(const TMatrix<T>&, TRCKind rckind, uint_4 ind=0);
56
[926]57 //! Return the kind of TMatrix (line,column,diagonal)
[804]58 TRCKind Kind() const { return kind; }
59 uint_4 NElts() const;
60 T& operator()(uint_4 i);
61 T operator()(uint_4 i) const;
62
63
64 TMatrixRC<T>& operator = (const TMatrixRC<T>& rc);
65
66// ---- A virer $CHECK$ Reza 03/2000
67// TVector<T> GetVect() const;
68// TMatrixRC<T>& operator += (const TMatrixRC<T>& rc);
69// TMatrixRC<T>& operator -= (const TMatrixRC<T>& rc);
70
71 TMatrixRC<T>& operator *= (T x);
72 TMatrixRC<T>& operator /= (T x);
73
74// TMatrixRC<T>& operator -= (T x);
75// TMatrixRC<T>& operator += (T x);
76
77
78 TMatrixRC<T>& LinComb(T a, T b, const TMatrixRC& rc, uint_4 first=0);
79 TMatrixRC<T>& LinComb(T b, const TMatrixRC<T>& rc, uint_4 first=0);
80
81 uint_4 IMaxAbs(uint_4 first=0);
[813]82 void Print(ostream & os) const ;
[804]83
84 static void Swap(TMatrixRC<T>& rc1, TMatrixRC<T>& rc2);
85
[926]86 //! Define Absolute value for uint_1
[804]87 inline static double Abs_Value(uint_1 v) {return (double) v;}
[926]88 //! Define Absolute value for uint_2
[804]89 inline static double Abs_Value(uint_2 v) {return (double) v;}
[926]90 //! Define Absolute value for int_2
[804]91 inline static double Abs_Value(int_2 v) {return (v>0)? (double) v: (double) -v;}
[926]92 //! Define Absolute value for int_4
[804]93 inline static double Abs_Value(int_4 v) {return (v>0)? (double) v: (double) -v;}
[926]94 //! Define Absolute value for int_8
[804]95 inline static double Abs_Value(int_8 v) {return (v>0)? (double) v: (double) -v;}
[926]96 //! Define Absolute value for uint_4
[804]97 inline static double Abs_Value(uint_4 v) {return (double) v;}
[926]98 //! Define Absolute value for uint_8
[804]99 inline static double Abs_Value(uint_8 v) {return (double) v;}
[926]100 //! Define Absolute value for r_4
[2147]101 inline static double Abs_Value(r_4 v) {return (double) fabs((double)v);}
[926]102 //! Define Absolute value for r_8
[804]103 inline static double Abs_Value(r_8 v) {return fabs(v);}
[926]104 //! Define Absolute value for complex r_4
105 inline static double Abs_Value(complex<r_4> v)
[804]106 {return sqrt(v.real()*v.real()+v.imag()*v.imag());}
[926]107 //! Define Absolute value for complex r_8
108 inline static double Abs_Value(complex<r_8> v)
[804]109 {return sqrt(v.real()*v.real()+v.imag()*v.imag());}
110
111protected:
[926]112 TMatrix<T>* matrix; //!< pointer to the TMatrix
113 T* data; //!< pointer to the beginnig of interesting datas
114 int_4 index; //!< index of the line/column
115 uint_4 step; //!< step of the line/column
116 TRCKind kind; //!< type: line, column or diagonal
[804]117};
118
[935]119//! Scalar product of two TMatrixRC
120/*!
121 \return sum[ a(i) * b(i) ]
122 */
[804]123template <class T>
124inline T operator * (const TMatrixRC<T>& a, const TMatrixRC<T>& b)
125 {
126 if ( a.NElts() != b.NElts() )
127 throw(SzMismatchError("TMatrixRC::operator * size mismatch\n"));
128 if ( a.Kind() != b.Kind() )
129 throw(SzMismatchError("TMatrixRC::operator * type mismatch\n"));
130 T sum = 0;
131 for(uint_4 i=0; i<a.NElts(); i++) sum += a(i)*b(i);
132 return sum;
133 }
134
[935]135//! Get the step in datas for a TMatrix for type rckind
136/*!
137 \param rckind : line, column or diagonal
138 \return step in TMatrix along TMatrixRC
139 */
[804]140template <class T>
141inline uint_4 TMatrixRC<T>::Step(const TMatrix<T>& m, TRCKind rckind)
[813]142 { switch (rckind) { case TmatrixRow : return m.Step(m.ColsKA());
143 case TmatrixCol : return m.Step(m.RowsKA());
[804]144 case TmatrixDiag : return (m.Step(m.RowsKA())+m.Step(m.ColsKA())); }
145 return 0; }
146
[935]147//! Get the origin of datas.
148/*!
149 Get the origin of datas in the TMatrix for a TMatrixRC for type
150 \b rckind and index \b index .
151 \param rckind : line, column or diagonal
152 \param index : index of the line or column.
153 \return adress of the first element in datas.
154 */
[804]155template <class T>
156inline T* TMatrixRC<T>::Org(const TMatrix<T>& m, TRCKind rckind, uint_4 index)
157 { switch (rckind) { case TmatrixRow : return const_cast<T *>(&(m(index,0)));
158 case TmatrixCol : return const_cast<T *>(&(m(0,index)));
159 case TmatrixDiag : return const_cast<T *>(m.Data()); }
160 return NULL; }
161
[926]162//! return number of elements for a TMatrixRC
[804]163template <class T> inline uint_4 TMatrixRC<T>::NElts() const
164 { if (!matrix) return 0;
165 switch (kind) { case TmatrixRow : return matrix->NCols();
166 case TmatrixCol : return matrix->NRows();
167 case TmatrixDiag : return matrix->NCols(); }
168 return 0; }
169
[926]170//! access of element \b i
[804]171template <class T>
172inline T& TMatrixRC<T>::operator()(uint_4 i) {return data[i*step];}
[926]173//! access of element \b i
[804]174template <class T>
175inline T TMatrixRC<T>::operator()(uint_4 i) const {return data[i*step];}
176
177////////////////////////////////////////////////////////////////
[935]178//! Typedef to simplify TMatrixRC<r_8> writing
[804]179typedef TMatrixRC<r_8> MatrixRC;
180
181
[926]182//! Default constructor
[772]183template <class T> TMatrixRC<T>::TMatrixRC()
184: matrix(NULL), data(NULL), index(0), step(0)
185{}
186
[926]187//! Constructor
188/*!
189 \param m : matrix
190 \param rckind : select line, column or diagonal
191 \param ind : number of the line or column
192*/
[772]193template <class T> TMatrixRC<T>::TMatrixRC(TMatrix<T>& m,TRCKind rckind,uint_4 ind)
194: matrix(&m), data(Org(m,rckind,ind)),
195 index(ind), step(Step(m,rckind)), kind(rckind)
196{
197if (kind == TmatrixDiag && m.NCols() != m.NRows())
198 throw(SzMismatchError("TMatrixRC::TMatrixRC(...,TmatrixDiag,...) size mismatch\n"));
199}
200
201////////////////////////////////////////////////////////////////
202// Acces aux rangees et colonnes de matrices
203
[926]204//! Return TMatrixRC for line \b r of matrix \b m
[772]205template <class T>
206TMatrixRC<T> TMatrixRC<T>::Row(TMatrix<T> & m, uint_4 r)
207{
208TMatrixRC<T> rc(m, TmatrixRow, r);
209return rc;
210}
211
[926]212//! Return TMatrixRC for column \b r of matrix \b m
[772]213template <class T>
214TMatrixRC<T> TMatrixRC<T>::Col(TMatrix<T> & m, uint_4 c)
215{
216TMatrixRC<T> rc(m, TmatrixCol, c);
217return rc;
218}
219
[926]220//! Return TMatrixRC for diagonal of matrix \b m
[772]221template <class T>
222TMatrixRC<T> TMatrixRC<T>::Diag(TMatrix<T> & m)
223{
224TMatrixRC<T> rc(m, TmatrixDiag);
225return rc;
226}
227
228
[804]229// ---- A virer $CHECK$ Reza 03/2000
230// template <class T> int_4 TMatrixRC<T>::Next()
231// {
232// if (!matrix || kind==TmatrixDiag) return -1;
233// index++;
234// if(kind == TmatrixRow) {
235// if(index > (int_4)matrix->NRows()) {index = (int_4)matrix->NRows(); return -1;}
236// data += matrix->NCols();
237// } else {
238// if (index > (int_4)matrix->NCols()) {index = (int_4)matrix->NCols(); return -1;}
239// data++;
240// }
241// return index;
242// }
[772]243
[804]244// template <class T> int_4 TMatrixRC<T>::Prev()
245// {
246// if (!matrix || kind == TmatrixDiag) return -1;
247// index--;
248// if(index < 0) {index = 0; return -1;}
249// if(kind == TmatrixRow) data -= matrix->NCols();
250// else data--;
251// return index;
252// }
[772]253
[926]254//! Set column \b c for this TMatrixRC
[772]255template <class T> int_4 TMatrixRC<T>::SetCol(int_4 c)
256{
257if(!matrix) return -1;
258if(c<0 || c>(int_4)matrix->NCols()) return -1;
259kind = TmatrixCol;
260index = c;
261step = Step(*matrix, TmatrixCol);
262data = Org(*matrix, TmatrixCol, c);
263return c;
264}
265
[926]266//! Set line \b r for this TMatrixRC
[772]267template <class T> int_4 TMatrixRC<T>::SetRow(int_4 r)
268{
269if(!matrix) return -1;
270if(r<0 && r>(int_4)matrix->NRows()) return -1;
271kind = TmatrixRow;
272index = r;
273step = Step(*matrix, TmatrixRow);
274data = Org(*matrix, TmatrixRow, r);
275return r;
276}
277
[926]278//! Set line diaginal for this TMatrixRC
[772]279template <class T> int_4 TMatrixRC<T>::SetDiag()
280{
281if (!matrix) return -1;
282if (matrix->NCols() != matrix->NRows())
283 throw(SzMismatchError("TMatrixRC::SetDiag size mismatch\n"));
284kind = TmatrixDiag;
285index = 0;
286step = Step(*matrix, TmatrixDiag);
287data = Org(*matrix, TmatrixDiag);
288return 0;
289}
290
[926]291//! Operator =
[772]292template <class T> TMatrixRC<T>& TMatrixRC<T>::operator = (const TMatrixRC<T>& rc)
293{
294matrix = rc.matrix;
295data = rc.data;
296index = rc.index;
297step = rc.step;
298kind = rc.kind;
299return *this;
300}
301
[804]302// ---- A virer $CHECK$ Reza 03/2000
303// template <class T> TVector<T> TMatrixRC<T>::GetVect() const
304// {
305// TVector<T> v(NElts());
306// for (uint_4 i=0; i<NElts(); i++) v(i) = (*this)(i);
307// return v;
308// }
[772]309
[804]310// template <class T> TMatrixRC<T>& TMatrixRC<T>::operator += (const TMatrixRC<T>& rc)
311// {
312// if ( NElts() != rc.NElts() )
313// throw(SzMismatchError("TMatrixRC::operator+= size mismatch\n"));
314// if ( kind != rc.kind )
315// throw(SzMismatchError("TMatrixRC::operator+= type mismatch\n"));
316// for (uint_4 i=0; i<NElts(); i++) (*this)(i) += rc(i);
317// return *this;
318// }
[772]319
[804]320// template <class T> TMatrixRC<T>& TMatrixRC<T>::operator -= (const TMatrixRC<T>& rc)
321// {
322// if( NElts() != rc.NElts() )
323// throw(SzMismatchError("TMatrixRC::operator-= size mismatch\n"));
324// if( kind != rc.kind )
325// throw(SzMismatchError("TMatrixRC::operator-= type mismatch\n"));
326// for(uint_4 i=0; i<NElts(); i++) (*this)(i) -= rc(i);
327// return *this;
328// }
[772]329
[926]330//! Operator to multiply by constant \b x
[772]331template <class T> TMatrixRC<T>& TMatrixRC<T>::operator *= (T x)
332{
333for(uint_4 i=0; i<NElts(); i++) (*this)(i) *= x;
334return *this;
335}
336
[926]337//! Operator to divide by constant \b x
[772]338template <class T> TMatrixRC<T>& TMatrixRC<T>::operator /= (T x)
339{
340for(uint_4 i=0; i<NElts(); i++) (*this)(i) /= x;
341return *this;
342}
343
344
[804]345// ---- A virer $CHECK$ Reza 03/2000
346// template <class T> TMatrixRC<T>& TMatrixRC<T>::operator -= (T x)
347// {
348// for(uint_4 i=0; i<NElts(); i++) (*this)(i) -= x;
349// return *this;
350// }
[772]351
[804]352// template <class T> TMatrixRC<T>& TMatrixRC<T>::operator += (T x)
353// {
354// for(uint_4 i=0; i<NElts(); i++) (*this)(i) += x;
355// return *this;
356// }
[772]357
[926]358//! Linear combination
359/*!
360 Do : \f$ MRC(i) = MRC(i)*a + rc(i)*b \f$
361 \return *this
362 */
[772]363template <class T>
364TMatrixRC<T>& TMatrixRC<T>::LinComb(T a, T b, const TMatrixRC<T>& rc, uint_4 first)
365{
366if ( NElts() != rc.NElts() )
367 throw(SzMismatchError("TMatrixRC::LinComb size mismatch\n"));
368if ( kind != rc.kind )
369 throw(SzMismatchError("TMatrixRC::LinComb type mismatch\n"));
370for(uint_4 i=first; i<NElts(); i++) (*this)(i) = (*this)(i)*a + rc(i)*b;
371return *this;
372}
373
[926]374//! Linear combination
375/*!
376 Do : \f$ MRC(i) = MRC(i) + rc(i)*b \f$
377 */
[772]378template <class T>
379TMatrixRC<T>& TMatrixRC<T>::LinComb(T b, const TMatrixRC<T>& rc, uint_4 first)
380{
381if ( NElts() != rc.NElts() )
382 throw(SzMismatchError("TMatrixRC::LinComb size mismatch\n"));
383if ( kind != rc.kind )
384 throw(SzMismatchError("TMatrixRC::LinComb type mismatch\n"));
385for(uint_4 i=first; i<NElts(); i++) (*this)(i) += rc(i)*b;
386return *this;
387}
388
[926]389//! Find maximum absolute value in TMatrixRC, search begin at \b first
[772]390template <class T> uint_4 TMatrixRC<T>::IMaxAbs(uint_4 first)
391{
392if (first>NElts())
393 throw(SzMismatchError("TMatrixRC::IMaxAbs size mismatch\n"));
394uint_4 imax = first;
395double vmax = Abs_Value((*this)(first));
396for(uint_4 i=first+1; i<NElts(); i++) {
397 double v = Abs_Value((*this)(i));
398 if(v > vmax) {vmax = v; imax = i;}
399}
400return imax;
401}
402
[926]403//! Print on stream \b os
[772]404template <class T>
[813]405void TMatrixRC<T>::Print(ostream & os) const
406{
407 os << " TMatrixRC<T>::Print(ostream & os) " << NElts() << " Kind="
408 << kind << " Index=" << index << " Step= " << step << endl;
409 for(uint_4 i=0; i<NElts(); i++) {
410 os << (*this)(i);
411 if (kind == TmatrixCol) os << endl;
412 else os << ", ";
413 }
414 os << endl;
415}
416
[926]417//! Swap two TMatrixRC of the same kind
[813]418template <class T>
[772]419void TMatrixRC<T>::Swap(TMatrixRC<T>& rc1, TMatrixRC<T>& rc2)
420{
421if(rc1.NElts() != rc2.NElts())
422 throw(SzMismatchError("TMatrixRC::Swap size mismatch\n"));
423if(rc1.kind != rc2.kind)
424 throw(SzMismatchError("TMatrixRC::Swap type mismatch\n"));
425if(rc1.data == rc2.data) return;
426for(uint_4 i=0; i<rc1.NElts(); i++)
427 {T tmp = rc1(i); rc1(i) = rc2(i); rc2(i) = tmp;}
428}
429
430
[926]431////////////////////////////////////////////////////////////////
[939]432// ---------------------------------------------------------- //
433// La classe de calcul simple sur les TMatrix //
434// ---------------------------------------------------------- //
[926]435////////////////////////////////////////////////////////////////
[772]436
437//**** Pour inversion
438#ifndef M_LN2
439#define M_LN2 0.69314718055994530942
440#endif
[1007]441//// CMV BUG BUG : sur OSF 5.0 DMINEXP est deconnant (~1.e+19 !!!)
[772]442#ifndef LN_MINDOUBLE
443#define LN_MINDOUBLE (M_LN2 * (DMINEXP - 1))
444#endif
445#ifndef LN_MAXDOUBLE
446#define LN_MAXDOUBLE (M_LN2 * DMAXEXP)
447#endif
448
[1007]449template <class T>
450int SimpleMatrixOperation<T>::GausPivScaling = PIV_GLOB_SCALE;
451
[926]452//! Gaussian pivoting
453/*!
[999]454 Do Gauss pivoting of \b a, doing the same operations on matrix \b b
[1007]455 \param computedet = true : return determimant of \b a (beware of over/underfloat)
456 \param computedet = false : return 1 if OK, 0 if not.
[999]457 \verbatim
458 Solve linear system A(n,n) * X(n,m) = B(n,m)
459 and put solution X in B for return.
460 \endverbatim
461 \warning If \b b is identity matrix, return inverse of \b a
462 \warning matrix \b a and \b b are modified.
[926]463 */
[772]464template <class T>
[1007]465T SimpleMatrixOperation<T>::GausPiv(TMatrix<T>& a, TMatrix<T>& b,bool computedet)
[772]466// Pivot de Gauss
467// * Attention: egcs impose que cette fonction soit mise dans le .cc
468// avant ::Inverse() (car Inverse() l'utilise)
469// {TMatrix A(a); TMatrix B(b); return (T) TMatrix::GausPiv(A,B);}
470{
471uint_4 n = a.NRows();
472if(n!=b.NRows())
473 throw(SzMismatchError("TMatrix::GausPiv size mismatch\n"));
474
[926]475T det = 1;
[1007]476
477//////////////////
478// Data scaling //
479//////////////////
480
481// Pas de normalisation des donnees
482if(GausPivScaling==PIV_NO_SCALE) {
483 if(computedet) det = (T) 1;
484// normalisation des donnees ligne par ligne
485} else if(GausPivScaling==PIV_ROW_SCALE) {
486 double nrm = 0.;
487 for(uint_4 iii=0; iii<a.NRows(); iii++) {
488 uint_4 jjj; double vmax = -1.;
489 for(jjj=0; jjj<a.NCols(); jjj++) {
490 double v = TMatrixRC<T>::Abs_Value(a(iii,jjj));
491 if(v>vmax) vmax = v;
492 }
493 if(vmax>0.) {
494 for(jjj=0; jjj<a.NCols(); jjj++) a(iii,jjj) /= (T) vmax;
495 for(jjj=0; jjj<b.NCols(); jjj++) b(iii,jjj) /= (T) vmax;
496 nrm += log(vmax);
497 } else return (T) 0;
[772]498 }
[1007]499 if(computedet) {
500 if(nrm <= LN_MINDOUBLE || nrm >= LN_MAXDOUBLE) {
501 cerr<<"GausPiv_Row: normalisation failure, "
502 <<"determinant has to be multiplied by exp("<<nrm<<")"<<endl;
503 } else det = (T) exp(nrm);
504 }
505// On fait une normalisation un peu brutale globale
506} else {
507 double vmin=MAXDOUBLE, vmax=0;
508 for(uint_4 iii=0; iii<a.NRows(); iii++)
509 for(uint_4 jjj=0; jjj<a.NCols(); jjj++) {
510 double v = TMatrixRC<T>::Abs_Value(a(iii,jjj));
511 if(v>vmax) vmax = v;
512 if(v<vmin && v>0.) vmin = v;
513 }
514 double nrm = sqrt(vmin*vmax);
515 if(nrm>0.) { a /= (T) nrm; b /= (T) nrm;} else return (T) 0;
516 if(computedet) {
517 nrm = a.NRows() * log(nrm);
518 if (nrm <= LN_MINDOUBLE || nrm >= LN_MAXDOUBLE) {
519 cerr<<"GausPiv_Glob: normalisation failure, "
520 <<"determinant has to be multiplied by exp("<<nrm<<")"<<endl;
521 } else det = (T) exp(nrm);
522 }
[772]523}
524
[1007]525////////////////////////////////////////
526// Gaussian elimination with pivoting //
527////////////////////////////////////////
528
[926]529TMatrixRC<T> pivRowa(a,TMatrixRC<T>::TmatrixRow);
530TMatrixRC<T> pivRowb(b,TMatrixRC<T>::TmatrixRow);
[772]531
532for(uint_4 k=0; k<n-1; k++) {
533 uint_4 iPiv = TMatrixRC<T>::Col(a, k).IMaxAbs(k);
534 if(iPiv != k) {
535 TMatrixRC<T> aIPiv(TMatrixRC<T>::Row(a,iPiv));
536 TMatrixRC<T> aK(TMatrixRC<T>::Row(a, k));
537 TMatrixRC<T>::Swap(aIPiv,aK);
538 TMatrixRC<T> bIPiv(TMatrixRC<T>::Row(b, iPiv));
539 TMatrixRC<T> bK(TMatrixRC<T>::Row(b, k));
540 TMatrixRC<T>::Swap(bIPiv,bK);
541 }
[926]542 T pivot = a(k,k);
543 if( TMatrixRC<T>::Abs_Value(pivot) < 1.e-50 ) return (T) 0;
[1007]544 if(computedet) det *= pivot;
[772]545 pivRowa.SetRow(k); // to avoid constructors
546 pivRowb.SetRow(k);
547 for (uint_4 i=k+1; i<n; i++) {
[926]548 T r = -a(i,k)/pivot;
[772]549 TMatrixRC<T>::Row(a, i).LinComb(r, pivRowa); // + rapide que -= r * pivRowa
550 TMatrixRC<T>::Row(b, i).LinComb(r, pivRowb);
551 }
552}
[1007]553if(computedet) det *= a(n-1, n-1);
[772]554
555// on remonte
556for(uint_4 kk=n-1; kk>0; kk--) {
[926]557 T pivot = a(kk,kk);
558 if( TMatrixRC<T>::Abs_Value(pivot) <= 1.e-50 ) return (T) 0;
[772]559 pivRowa.SetRow(kk); // to avoid constructors
560 pivRowb.SetRow(kk);
561 for(uint_4 jj=0; jj<kk; jj++) {
[926]562 T r = -a(jj,kk)/pivot;
[772]563 TMatrixRC<T>::Row(a, jj).LinComb(r, pivRowa);
564 TMatrixRC<T>::Row(b, jj).LinComb(r, pivRowb);
565 }
566}
567
568for(uint_4 l=0; l<n; l++) {
[926]569 if( TMatrixRC<T>::Abs_Value(a(l,l)) <= 1.e-50 ) return (T) 0;
[772]570 TMatrixRC<T>::Row(b, l) /= a(l,l);
571}
572
573return det;
574}
575
[926]576//! Return the inverse matrix of \b A
[772]577template <class T>
578TMatrix<T> SimpleMatrixOperation<T>::Inverse(TMatrix<T> const & A)
579{
[976]580TMatrix<T> a(A,false);
[813]581TMatrix<T> b(a.NCols(),a.NRows()); b = IdentityMatrix(1.);
[1007]582if(GausPiv(a,b)==(T) 0)
[926]583 throw(MathExc("TMatrix Inverse() Singular Matrix"));
[772]584return b;
585}
586
587
[926]588////////////////////////////////////////////////////////////////
[939]589// ---------------------------------------------------------- //
590// La classe fit lineaire //
591// ---------------------------------------------------------- //
[926]592////////////////////////////////////////////////////////////////
593
[939]594//! Creator
595template <class T>
596LinFitter<T>::LinFitter()
[804]597{
598}
[772]599
[939]600//! Destructor
601template <class T>
602LinFitter<T>::~LinFitter()
[804]603{
604}
605
[939]606// fit lineaire des y(k) en tant que somme de c(i)f(i,x(k)), i=0..nf-1;
607//! Linear fitting
608/*!
609 Linear fit of y(k) as the sum of \f$ c(i)f(i,x(k)), i=0..nf-1 \f$
610 \param x : vector of X values
611 \param y : vector of datas
612 \param nf: number of functions
613 \param f : f(i,x(k)), i=0..nf-1
614 \return c : vector of solutions
615 \return return chisquare
616 */
617template <class T>
618r_8 LinFitter<T>::LinFit(const TVector<T>& x, const TVector<T>& y,
619 uint_4 nf, T (*f)(uint_4,T), TVector<T>& c)
[804]620{
[939]621uint_4 n = x.NElts();
[2421]622if (n != y.NElts())
623 throw SzMismatchError("LinFitter<T>::LinFit(x,y,nf,f,c)/Error x.NElts() <> y.Nelts() ");
[939]624
625TMatrix<T> fx(nf,n);
[804]626
[939]627for(uint_4 i=0; i<nf; i++)
628 for(uint_4 j=0; j<n; j++) fx(i,j) = f(i,x(j));
629
630return LinFit(fx,y,c);
[804]631}
632
[939]633// fit lineaire des y(k) en tant que somme de c(i)f(i,x(k)), i=0..nf-1,
634// la matrice fx contient les valeurs des f: fx(i,j) = f(i, x(j)).
635//! Linear fitting
636/*!
637 Linear fit of y(k) as the sum of \f$ c(i)f(i,x(k)), i=0..nf-1 \f$.
638 \param fx : matrix which contains \f$ fx(i,j) = f(i, x(j)) \f$.
639 \param y : vector of datas
640 \return c : vector of solutions
641 \return return chisquare
642 */
643template <class T>
644r_8 LinFitter<T>::LinFit(const TMatrix<T>& fx, const TVector<T>& y, TVector<T>& c)
[804]645{
[939]646uint_4 n = y.NElts();
[2421]647if (n != fx.NCol())
648 throw SzMismatchError("LinFitter<T>::LinFit(fx, y, c)/Error y.NElts() <> fx.Nelts() ");
[804]649
[939]650uint_4 nf = fx.NRows();
[804]651
[939]652TMatrix<T> a(nf,nf);
[804]653
[939]654for(uint_4 j=0; j<nf; j++)
655 for(uint_4 k=j; k<nf; k++)
656 a(j,k) = a(k,j) = TMatrixRC<T>::Row(const_cast<TMatrix<T> &>(fx), j)
657 * TMatrixRC<T>::Row(const_cast<TMatrix<T> &>(fx), k);
[804]658
[939]659c = fx * y;
[804]660
[2421]661if(SimpleMatrixOperation<T>::GausPiv(a,c)==(T) 0)
662 throw SingMatrixExc("LinFitter<T>::LinFit(fx, y, c) - Non invertible matrix (by GausPiv())");
[804]663
[939]664r_8 xi2 = 0., ax;
665T x;
666for(uint_4 k=0; k<n; k++) {
667 x = (T) 0;
668 for(uint_4 i=0; i<nf; i++) x += c(i) * fx(i,k);
669 x -= y(k);
670 ax = TMatrixRC<T>::Abs_Value(x);
671 xi2 += ax*ax;
[804]672}
[939]673return xi2;
674}
[804]675
[939]676// fit lineaire des y(k) en tant que somme de c(i)f(i,x(k)), i=0..nf-1,
677// errY2 contient les carres des erreurs sur les Y.
678// au retour, errC contient les erreurs sur les coefs.
679//! Linear fitting with errors
680/*!
681 Linear fit with errors of y(k) as the sum of \f$ c(i)f(i,x(k)), i=0..nf-1 \f$.
682 \param x : vector of X values
683 \param y : vector of datas
684 \param errY2 : vector of errors square on Y
685 \param nf: number of functions
686 \param f : f(i,x(k)), i=0..nf-1
687 \return c : vector of solutions
688 \return errC : vector of errors on solutions C
689 \return return chisquare
690 */
691template <class T>
692r_8 LinFitter<T>::LinFit(const TVector<T>& x, const TVector<T>& y,
693 const TVector<T>& errY2, uint_4 nf, T (*f)(uint_4,T),
694 TVector<T>& c, TVector<T>& errC)
695{
696uint_4 n = x.NElts();
[2421]697if (n != y.NElts())
698 throw SzMismatchError("LinFitter<T>::LinFit(x,y,errY2,nf,f,c,errC)/Error x.NElts() <> y.Nelts() ");
[804]699
[939]700TMatrix<T> fx(nf, n);
701for(uint_4 i=0; i<nf; i++)
702 for(uint_4 j=0; j<n; j++)
703 fx(i,j) = f(i,x(j));
[804]704
[939]705return LinFit(fx,y,errY2,c,errC);
[804]706}
707
[939]708// fit lineaire des y(k) en tant que somme de c(i)f(i,x(k)), i=0..nf-1,
709// la matrice fx contient les valeurs des f:
710// fx(i,j) = f(i, x(j)).
711// errY2 contient les carres des erreurs sur les Y.
712// au retour, errC contient les erreurs sur les coefs.
713//! Linear fitting with errors
714/*!
715 \param fx : matrix which contains \f$ fx(i,j) = f(i, x(j)) \f$.
716 \param y : vector of datas
717 \param errY2 : vector of errors square on Y
718 \return c : vector of solutions
719 \return errC : vector of errors on solutions on C
720 \return return chisquare
721 */
722template <class T>
723r_8 LinFitter<T>::LinFit(const TMatrix<T>& fx, const TVector<T>& y,
724 const TVector<T>& errY2,TVector<T>& c, TVector<T>& errC)
[804]725{
[939]726uint_4 n = y.NElts();
[2421]727if( n != errY2.NElts())
728 throw SzMismatchError("LinFitter<T>::LinFit(fx,y,errY2,c,errC)/Error y.NElts() <> errY2.Nelts() ");
729if( n != fx.NCol())
730 throw SzMismatchError("LinFitter<T>::LinFit(fx,y,errY2,c,errC)/Error y.NElts() <> fx.NCols() ");
[939]731uint_4 nf = fx.NRows();
[804]732
[939]733TMatrix<T> a(nf,nf);
[804]734
[939]735c.Realloc(nf);
736errC.Realloc(nf);
[804]737
[939]738for(uint_4 j=0; j<nf; j++)
739 for(uint_4 k=j; k<nf; k++) {
740 T x=0;
741 // Matrice a inverser
742 for(uint_4 l=0; l<n; l++) x += fx(j,l)*fx(k,l)/errY2(l);
743 a(j,k) = a(k,j) = x;
744 }
[804]745
[939]746TMatrix<T> d(nf,nf+1);
747for(uint_4 k=0; k<nf; k++) {
748 T x = (T) 0;
749 // Second membre 1ere colonne
750 for(uint_4 l=0; l<n; l++) x += fx(k,l)*y(l)/errY2(l);
751 d(k,0) = x;
752 // Reste second membre = Id.
753 for(uint_4 m=1; m<=nf; m++) d(k,m) = (k==m)?1:0;
754}
[804]755
[2421]756if(SimpleMatrixOperation<T>::GausPiv(a,d)==(T) 0)
757 throw SingMatrixExc("LinFitter<T>::LinFit(...ErrY2...) - Non invertible matrix (by GausPiv())");
[804]758
[2421]759
[939]760for(uint_4 l=0; l<nf; l++) {
761 c(l) = d(l,0); // Parametres = 1ere colonne
762 errC(l) = d(l,l+1); // Erreurs = diag inverse.
763}
[804]764
[939]765r_8 xi2 = 0., ax;
766T x;
767for(uint_4 jj=0; jj<n; jj++) {
768 x = (T) 0;
769 for(uint_4 ii=0; ii<nf; ii++) x += c(ii) * fx(ii,jj);
[804]770 x -= y(jj);
[939]771 ax = TMatrixRC<T>::Abs_Value(x);
772 xi2 += ax*ax/TMatrixRC<T>::Abs_Value(errY2(jj));
[804]773 }
774 return xi2;
775}
776
[939]777///////////////////////////////////////////////////////////////////////////
778///////////////////////////////////////////////////////////////////////////
779///////////////////////////////////////////////////////////////////////////
780///////////////////////////////////////////////////////////////////////////
[813]781void _ZZ_TestTMatrixRC_YY_(TMatrix<r_8> & m)
782{
783 cout << " ::::: _ZZ_TestTMatrixRC_YY_ :::: M= \n" << m << endl;
784 TMatrixRC<r_8> l0 = TMatrixRC<r_8>::Row(m,0);
785 cout << "TMatrixRC<r_8>::Row(m,0) = \n" ;
786 l0.Print(cout);
787 TMatrixRC<r_8> l1 = TMatrixRC<r_8>::Row(m,1);
788 cout << "TMatrixRC<r_8>::Row(m,1) = \n" ;
789 l1.Print(cout);
790 l0.SetRow(2);
791 cout << "TMatrixRC<r_8>::l0.SetRow(2 = \n" ;
792 l0.Print(cout);
[804]793
[813]794 TMatrixRC<r_8> c0 = TMatrixRC<r_8>::Col(m,0);
795 cout << "TMatrixRC<r_8>::Col(m,0) = \n" ;
796 c0.Print(cout);
797 TMatrixRC<r_8> c1 = TMatrixRC<r_8>::Col(m,1);
798 cout << "TMatrixRC<r_8>::Col(m,1) = \n" ;
799 c1.Print(cout);
800 c0.SetCol(2);
801 cout << "TMatrixRC<r_8>::c0.SetCol(2) = \n" ;
802 c0.Print(cout);
803 TMatrixRC<r_8> c00 = TMatrixRC<r_8>::Col(m,0);
804 TMatrixRC<r_8>::Swap(c0, c00);
805 cout << " ::::: M Apres Swap = \n" << m << endl;
806
807}
[804]808
[772]809///////////////////////////////////////////////////////////////
810#ifdef __CXX_PRAGMA_TEMPLATES__
811// Instances gestion lignes/colonnes
812#pragma define_template TMatrixRC<int_4>
813#pragma define_template TMatrixRC<r_4>
814#pragma define_template TMatrixRC<r_8>
[926]815#pragma define_template TMatrixRC< complex<r_4> >
816#pragma define_template TMatrixRC< complex<r_8> >
[804]817#pragma define_template SimpleMatrixOperation<int_4>
818#pragma define_template SimpleMatrixOperation<r_4>
[772]819#pragma define_template SimpleMatrixOperation<r_8>
[926]820#pragma define_template SimpleMatrixOperation< complex<r_4> >
821#pragma define_template SimpleMatrixOperation< complex<r_8> >
[939]822#pragma define_template LinFitter<r_4>
823#pragma define_template LinFitter<r_8>
824#pragma define_template LinFitter< complex<r_4> >
825#pragma define_template LinFitter< complex<r_8> >
[772]826#endif
827
828#if defined(ANSI_TEMPLATES) || defined(GNU_TEMPLATES)
829// Instances gestion lignes/colonnes
830template class TMatrixRC<int_4>;
831template class TMatrixRC<r_4>;
832template class TMatrixRC<r_8>;
[926]833template class TMatrixRC< complex<r_4> >;
834template class TMatrixRC< complex<r_8> >;
[804]835template class SimpleMatrixOperation<int_4>;
[772]836template class SimpleMatrixOperation<r_4>;
837template class SimpleMatrixOperation<r_8>;
[926]838template class SimpleMatrixOperation< complex<r_4> >;
839template class SimpleMatrixOperation< complex<r_8> >;
[939]840template class LinFitter<r_4>;
841template class LinFitter<r_8>;
842template class LinFitter< complex<r_4> >;
843template class LinFitter< complex<r_8> >;
[772]844#endif
[935]845
Note: See TracBrowser for help on using the repository browser.