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

Last change on this file since 1583 was 1007, checked in by ansari, 25 years ago

gestion det dans GausPiv et new norm cmv 16/5/00

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