Herezh_dev/Resolin/Matrices/Mat_pleine.cc

2277 lines
78 KiB
C++
Raw Permalink Normal View History

2021-09-28 17:28:23 +02:00
// FICHIER : Mat_pleine.cc
// CLASSE : Mat_pleine
// This file is part of the Herezh++ application.
//
// The finite element software Herezh++ is dedicated to the field
// of mechanics for large transformations of solid structures.
// It is developed by Gérard Rio (APP: IDDN.FR.010.0106078.000.R.P.2006.035.20600)
// INSTITUT DE RECHERCHE DUPUY DE LÔME (IRDL) <https://www.irdl.fr/>.
//
// Herezh++ is distributed under GPL 3 license ou ultérieure.
//
2023-05-03 17:23:49 +02:00
// Copyright (C) 1997-2022 Université Bretagne Sud (France)
2021-09-28 17:28:23 +02:00
// AUTHOR : Gérard Rio
// E-MAIL : gerardrio56@free.fr
//
// This program 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 3 of the License,
// or (at your option) any later version.
//
// This program 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 this program. If not, see <https://www.gnu.org/licenses/>.
//
// For more information, please consult: <https://herezh.irdl.fr/>.
//#include "Debug.h"
# include <iostream>
using namespace std; //introduces namespace std
#include <math.h>
#include <stdlib.h>
#include "Sortie.h"
#include "ConstMath.h"
#include "MathUtil.h"
#include "Mat_pleine.h"
#include "MatDiag.h"
2023-05-03 17:23:49 +02:00
#include "ParaGlob.h"
2021-09-28 17:28:23 +02:00
#ifndef Mat_pleine_H_deja_inclus
// Constructeur par defaut
#ifndef MISE_AU_POINT
inline
#endif
Mat_pleine::Mat_pleine ():
2023-05-03 17:23:49 +02:00
Mat_abstraite(RECTANGLE,CHOLESKY,RIEN_PRECONDITIONNEMENT)
,val_de_base(),val()
2021-09-28 17:28:23 +02:00
{ };
// Constructeur utile quand les dimensions de la matrice sont connues
#ifndef MISE_AU_POINT
inline
#endif
Mat_pleine::Mat_pleine (int nb_ligne,int nb_colonne,double val_init):
Mat_abstraite(RECTANGLE,CHOLESKY,RIEN_PRECONDITIONNEMENT)
2023-05-03 17:23:49 +02:00
,val_de_base(nb_ligne*nb_colonne)
2021-09-28 17:28:23 +02:00
{ if (nb_ligne == nb_colonne)
{type_matrice = CARREE;
if (nb_ligne < 4) // par défaut pour de petite matrice on utilise cramer
type_resolution = CRAMER;
};
2023-05-03 17:23:49 +02:00
Creation_tab_val(nb_ligne,nb_colonne);
Initialise(val_init);
// Initialise(nb_ligne,nb_colonne,val_init);
2021-09-28 17:28:23 +02:00
};
// Constructeur de copie
#ifndef MISE_AU_POINT
inline
#endif
Mat_pleine::Mat_pleine (const Mat_pleine& mat_pl) :
2023-05-03 17:23:49 +02:00
Mat_abstraite(mat_pl),val_de_base(mat_pl.val_de_base)
{Creation_tab_val(mat_pl.Nb_ligne(),mat_pl.Nb_colonne());
};
/// --------- cas très particulier d'un vecteur val_de_base lié
/// --------- utilisé uniquement par des classes friend
// c-a-d qui ne gère pas sa mémoire -> memoire = false
/// il s'agit d'un tuillage, val_de_base est positionné a "vec"
/// aucun test n'est fait concernant la taille disponible à partir de "vec"
/// routine dangereuse !!!!!
/// avec_recopie == true : les infos de mat_pl sont recopiées, sinon non
#ifndef MISE_AU_POINT
inline
#endif
Mat_pleine::Mat_pleine(bool avec_recopie,double* vec, const Mat_pleine& mat_pl) :
Mat_abstraite(mat_pl),val_de_base(mat_pl.val_de_base.Taille(),vec,false)
{Creation_tab_val(mat_pl.Nb_ligne(),mat_pl.Nb_colonne());
if (avec_recopie)
val_de_base = mat_pl.val_de_base;
};
2021-09-28 17:28:23 +02:00
#ifndef MISE_AU_POINT
inline
#endif
Mat_pleine::~Mat_pleine ()
// Destructeur
2023-05-03 17:23:49 +02:00
{ // val_de_base est automatiquement détruit
int nb_ligne = val.Taille();
for (int i =1;i<=nb_ligne;i++)
delete val(i); // en fait ne fait pas grand chose !! car il s'agit
// de vecteurs liés, on pourrait zapper !
};
2021-09-28 17:28:23 +02:00
// fonction permettant de creer une nouvelle instance d'element
#ifndef MISE_AU_POINT
inline
#endif
Mat_abstraite * Mat_pleine::NouvelElement() const
{ Mat_abstraite * a;
a = new Mat_pleine(*this);
return a;
};
// surcharge de l'opérateur d'affectation: cas des matrices abstraites
// il y a ajustement automatique de taille de la matrice en fonction de mat_pl
#ifndef MISE_AU_POINT
inline
#endif
Mat_abstraite & Mat_pleine::operator = ( const Mat_abstraite & b)
{
// #ifdef MISE_AU_POINT
// dans le cas où les matrices ne sont pas de caractéristiques identiques
// il y a un message d'erreur
if (type_matrice != b.Type_matrice())
{cout << "erreur dans l'operation d'affectation";
cout << " les matrices sont de types différents "
<< Nom_matrice(type_matrice) << " " << Nom_matrice(b.Type_matrice()) << '\n';
cout << "Mat_pleine::operator = ( const Mat_abstraite & b)" << endl;
Sortie(1);
2023-05-03 17:23:49 +02:00
};
2021-09-28 17:28:23 +02:00
// #endif
const Mat_pleine & mat_pl = *((Mat_pleine*) & b);
2023-05-03 17:23:49 +02:00
// et on affecte entre matrices Mat_pleine
*this = mat_pl;
// if ( val.Taille() != 0 )
// // cas ou la matrice se trouvant a gauche du signe = n'est pas vide
// {
// Libere(); // desallocation de ce vecteur
// };
// // traitement----
// // cas ou la matrice se trouvant a droite du signe = est vide
// // on a rien a faire car de toute manière la matrice this viens d'être
// // supprimé dans le cas où elle était non nul
// if ( mat_pl.Nb_ligne()!=0 )
// // autres cas
// {
// int taill_ligne = mat_pl.Nb_ligne();
// val.Change_taille(taill_ligne);
// for (int i=1;i<=taill_ligne;i++)
// // affectation de la ieme ligne de la matrice mat_pl a la ieme ligne
// // de la matrice courante
// val(i)=mat_pl(i);
// };
return (*this);
2021-09-28 17:28:23 +02:00
};
// transfert des informations de *this dans la matrice passée en paramètre
// la matrice paramètre est au préalable, mise à 0.
#ifndef MISE_AU_POINT
inline
#endif
void Mat_pleine::Transfert_vers_mat( Mat_abstraite & b )
{b.Initialise(0.); // init de la matrice
// puis on transfert
int nb_ligne = val.Taille();
int nb_col = Nb_colonne();
for (int i=1;i<=nb_ligne;i++)
2023-05-03 17:23:49 +02:00
{ Vecteur& vee = *(val(i));
2021-09-28 17:28:23 +02:00
for (int j=1;j<= nb_col; j++)
b(i,j)=vee(j);
};
};
#ifndef MISE_AU_POINT
inline
#endif
Mat_pleine& Mat_pleine::operator= ( const Mat_pleine& mat_pl)
// Surcharge de l'operateur = : affectation d'une matrice a une autre
// il y a ajustement automatique de taille de la matrice en fonction de mat_pl
2023-05-03 17:23:49 +02:00
{ if (val_de_base.Taille() == mat_pl.val_de_base.Taille())
// l'ensemble à la même taille, on regarde le nombre de colonne
{if (mat_pl.Nb_ligne() != val.Taille())
// on a la même taille globale mais pas le même nombre de ligne
// il faut changer le tableau des lignes
{Change_tailles_tab_val(mat_pl.Nb_ligne(),mat_pl.Nb_colonne());
};
// sinon on ne fait rien car tout est ok
}
else // cas ou val_de_base n'a pas la bonne taille
// on réaffecte
{ val_de_base.Change_taille(mat_pl.val_de_base.Taille());
// on change le tableau val des vecteurs liés
Change_tailles_tab_val(mat_pl.Nb_ligne(),mat_pl.Nb_colonne());
};
// toutes les tailles sont ok on peut affecter globalement
val_de_base = mat_pl.val_de_base;
//if ( val.Taille() != 0 )
// // cas ou la matrice se trouvant a gauche du signe = n'est pas vide
// {
// Libere(); // desallocation de ce vecteur
// };
// // traitement----
// // cas ou la matrice se trouvant a droite du signe = est vide
// // on a rien a faire car de toute manière la matrice this viens d'être
// // supprimé dans le cas où elle était non nul
// if ( mat_pl.Nb_ligne()!=0 )
// // autres cas
// {
// int taill_ligne = mat_pl.Nb_ligne();
// val.Change_taille(taill_ligne);
// for (int i=1;i<=taill_ligne;i++)
// // affectation de la ieme ligne de la matrice mat_pl a la ieme ligne
// // de la matrice courante
// val(i)=mat_pl(i);
// };
return (*this);
2021-09-28 17:28:23 +02:00
};
// Surcharge de l'operateur += : addition d'une matrice a la matrice courante
#ifndef MISE_AU_POINT
inline
#endif
void Mat_pleine::operator+= (const Mat_abstraite& b)
{
#ifdef MISE_AU_POINT
// dans le cas où les matrices ne sont pas de caractéristiques identiques
// il y a un message d'erreur
if ((type_matrice != b.Type_matrice())&&(b.Type_matrice() != DIAGONALE))
{cout << " erreur les matrices sont de types differents "
<< Nom_matrice(type_matrice) << " " << Nom_matrice(b.Type_matrice()) << '\n';
cout << "Mat_pleine::operator+= (const Mat_abstraite& b)" << endl;
Sortie(1);
}
// si les tailles ne sont pas identique pb
if (this->Nb_ligne() != b.Nb_ligne())
{cout << " les matrices ont un nombre de ligne different "
<< this->Nb_ligne() << " " << b.Nb_ligne() << '\n';
cout << "Mat_pleine::operator+= (const Mat_abstraite& b)" << endl;
Sortie(1);
}
// on ne vérifie pas le nombre de colonne pour les matrices diagonales
if ((this->Nb_colonne() != b.Nb_colonne())&&(b.Type_matrice() != DIAGONALE))
{cout << " les matrices ont un nombre de ligne different "
<< this->Nb_colonne() << " " << b.Nb_colonne() << '\n';
cout << "Mat_pleine::operator+= (const Mat_abstraite& b)" << endl;
Sortie(1);
}
#endif
if (b.Type_matrice() != DIAGONALE)
{ const Mat_pleine & mat_pl = *((Mat_pleine*) & b);
this->operator+= (mat_pl);
}
else // cas d'une matrice argument diagonale
{ int nb_lign = val.Taille();
const MatDiag & mat_d = *((MatDiag*) & b);
for (int i=1;i<=nb_lign;i++)
(*this)(i,i) += mat_d(i,i);
}
};
// Surcharge de l'operateur -= : soustraction d'une matrice a la matrice courante
#ifndef MISE_AU_POINT
inline
#endif
void Mat_pleine::operator-= (const Mat_abstraite& b)
{
#ifdef MISE_AU_POINT
// dans le cas où les matrices ne sont pas de caractéristiques identiques
// il y a un message d'erreur
if ((type_matrice != b.Type_matrice())&&(b.Type_matrice() != DIAGONALE))
{cout << " erreur les matrices sont de types différents "
<< Nom_matrice(type_matrice) << " " << Nom_matrice(b.Type_matrice()) << '\n';
cout << "Mat_pleine::operator-= (const Mat_abstraite& b)" << endl;
Sortie(1);
}
// si les tailles ne sont pas identique pb
if (this->Nb_ligne() != b.Nb_ligne())
{cout << " les matrices ont un nombre de ligne différent "
<< this->Nb_ligne() << " " << b.Nb_ligne() << '\n';
cout << "Mat_pleine::operator-= (const Mat_abstraite& b)" << endl;
Sortie(1);
}
// on ne vérifie pas le nombre de colonne pour les matrices diagonales
if ((this->Nb_colonne() != b.Nb_colonne())&&(b.Type_matrice() != DIAGONALE))
{cout << " les matrices ont un nombre de ligne différent "
<< this->Nb_colonne() << " " << b.Nb_colonne() << '\n';
cout << "Mat_pleine::operator-= (const Mat_abstraite& b)" << endl;
Sortie(1);
}
#endif
if (b.Type_matrice() != DIAGONALE)
{ const Mat_pleine & mat_pl = *((Mat_pleine*) & b);
this->operator-= (mat_pl);
}
else // cas d'une matrice argument diagonale
{ int nb_lign = val.Taille();
const MatDiag & mat_d = *((MatDiag*) & b);
for (int i=1;i<=nb_lign;i++)
(*this)(i,i) -= mat_d(i,i);
}
};
// Surcharge de l'operateur == : test d'egalite entre deux matrices
#ifndef MISE_AU_POINT
inline
#endif
int Mat_pleine::operator== (const Mat_abstraite& b) const
{
#ifdef MISE_AU_POINT
// dans le cas où les matrices ne sont pas de caractéristiques identiques
// il y a un message d'erreur
if (type_matrice != b.Type_matrice())
{cout << " erreur les matrices sont de types différents "
<< Nom_matrice(type_matrice) << " " << Nom_matrice(b.Type_matrice()) << '\n';
cout << "Mat_pleine::operator+= (const Mat_abstraite& b)" << endl;
Sortie(1);
}
#endif
const Mat_pleine & mat_pl = *((Mat_pleine*) & b);
#ifdef MISE_AU_POINT
// si les tailles ne sont pas identique pb
if (this->Nb_ligne() != b.Nb_ligne())
{cout << " les matrices ont un nombre de ligne différent "
<< this->Nb_ligne() << " " << b.Nb_ligne() << '\n';
cout << "Mat_pleine::operator+= (const Mat_abstraite& b)" << endl;
Sortie(1);
}
if (this->Nb_colonne() != b.Nb_colonne())
{cout << " les matrices ont un nombre de ligne différent "
<< this->Nb_colonne() << " " << b.Nb_colonne() << '\n';
cout << "Mat_pleine::operator+= (const Mat_abstraite& b)" << endl;
Sortie(1);
}
#endif
return this->operator==(mat_pl);
};
#ifndef MISE_AU_POINT
inline
#endif
void Mat_pleine::Affiche () const
// Affiche les donnees de la matrice : dimensions, composantes
{int nbligne = Nb_ligne();
int nbcol = Nb_colonne();
cout << "\tNombre de lignes : " << nbligne ; //<< "\n";
cout << "\tNombre de colonnes : " << nbcol << "\n";
// for (int i=1;i<=nbligne;i++)
// {
// cout << "Ligne " << i << " :\t{\t";
// for (int j=1;j<=nbcol;j++)
// cout << (*this)(i,j) << "\t";
// if (i!= nbligne) {cout << "}\n";}
// else {cout << "}";}
// };
int i,j;
cout << " ";
for (j=1;j<=nbcol;j++) cout << "col " << setw(4) << j << " " ;
cout << '\n';
for (i=1;i<=nbligne;i++)
{ cout << "lig " << setw(4) << i << " * ";
for (j=1;j<=nbcol;j++)
cout << setw(14) << setprecision(7) << (*this)(i,j);
cout << '\n';
}
cout << endl;
};
2023-05-03 17:23:49 +02:00
//#ifndef MISE_AU_POINT
// inline
//#endif
//void Mat_pleine::Initialise (double val_init)
//// Initialise les composantes de la matrice a val_init
//{ int nb_ligne = val.Taille();
// int nb_col = Nb_colonne();
// for (int i=1;i<=nb_ligne;i++)
// // initialisation des composantes de la matrice a val_init
// for (int j=1;j<=nb_col;j++)
// (*this)(i,j)=val_init;
//};
2021-09-28 17:28:23 +02:00
#ifndef MISE_AU_POINT
inline
#endif
void Mat_pleine::Initialise (int nb_ligne,int nb_colonne,double val_init)
// Initialise les composantes de la matrice de nombre de ligne nb_ligne
// et de nombre de colonne nb_colonne a val_init
{
#ifdef MISE_AU_POINT
if ( (nb_ligne<0) || (nb_colonne<0) )
{ cout << "Erreur sur les dimensions de la matrice !\n";
cout << "MAT_PLEINE::MAT_PLEINE (int ,int ,double )\n";
Sortie(1);
};
#endif
2023-05-03 17:23:49 +02:00
int taille_glob = nb_ligne * nb_colonne;
// changement de taille éventuelle
val_de_base.Change_taille(taille_glob);
Change_tailles_tab_val(nb_ligne,nb_colonne);
// initialisation globalle
val_de_base.Inita(val_init);
//
//
// val.Change_taille(nb_ligne);
// for (int i=1;i<=nb_ligne;i++)
// val(i).Change_taille(nb_colonne);
// for (int i=1;i<=nb_ligne;i++)
// // initialisation des composantes de la matrice a val_init
// for (int j=1;j<=nb_colonne;j++)
// (*this)(i,j)=val_init;
2021-09-28 17:28:23 +02:00
};
#ifndef MISE_AU_POINT
inline
#endif
int Mat_pleine::Symetrie () const
// Test sur la symetrie de la matrice
// Retourne 1 si la matrice est symetrique
// 0 sinon
2023-05-03 17:23:49 +02:00
{if ( Nb_ligne()!=Nb_colonne() )
return 0;
2021-09-28 17:28:23 +02:00
int nb_lign = val.Taille();
for (int i=1;i<=nb_lign;i++)
for (int j=1;j<i;j++)
// if ( (*this)(i,j)!=(*this)(j,i) ) // test d'egalite de la ieme jieme composante
if ( diffpetit((*this)(i,j),(*this)(j,i))) // test d'egalite de la ieme jieme composante
{
// affichage -- débug
// cout << " ("<<i<<","<<j<<") et ("<<j<<","<<i<<") = "<<(*this)(i,j)<<" "<< (*this)(j,i)<<", " << endl;
// fin affichage -- débug
return 0;
};
return 1;
};
#ifndef MISE_AU_POINT
inline
#endif
// affiche les termes non symétriques de la matrice s'il y en a
void Mat_pleine::AfficheNonSymetries() const
{ int nb_lign = val.Taille();
cout << "\n affichage des termes non symetriques de la matrice: ";
for (int i=1;i<=nb_lign;i++)
for (int j=1;j<i;j++)
if ( diffpetit((*this)(i,j),(*this)(j,i))) // test d'egalite de la ieme jieme composante
cout << " ("<<i<<","<<j<<") et ("<<j<<","<<i<<") = "<<(*this)(i,j)<<" "<< (*this)(j,i)<<", ";
};
#ifndef MISE_AU_POINT
inline
#endif
Vecteur Mat_pleine::Resol_syst
(const Vecteur& b,const double &tole, const int maxi, const int rest)
// Resolution du systeme Ax=b
// N.B. : dans le cas de la methode de Cholesky, aucun test n'est realise
// pour verifier que la matrice A est definie positive ( condition necessaire a
// l'utilisation de la methode de Cholesky )
// La matrice est decomposee de la facon suivante : A=B.Bt
// ou :
// B : matrice triangulaire inferieure reguliere
// Bt: transposee de la matrice B ( triangulaire superieure )
{
#ifdef MISE_AU_POINT
switch (type_resolution)
{case CHOLESKY:
// si c'est cholesky normal, il faut une matrice symétrique
if ( !(*this).Symetrie() )
{ if (ParaGlob::NiveauImpression() > 2)
{ cout << "\n warning : attention ne fonctionne qu'avec une matrice symetrique !"
<< " ce qui n'est pas le cas ici !! "
<< "\n on utilise ici le triangle inferieure pour la resolution ! \n";
cout << "MAT_PLEINE::RESOL_SYST (const Vecteur& )\n";
AfficheNonSymetries();
};
};
// de plus si c'est cramer ou cholesky après symétrisation, on doit avoir une matrice carrée
case CRAMER : case SYMETRISATION_PUIS_CHOLESKY :
if ( Nb_ligne()!=Nb_colonne() )
{ cout << "Erreur : la matrice doit etre carree !\n";
cout << "MAT_PLEINE::RESOL_SYST (const Vecteur& )\n";
Sortie(1);
};
break;
default: break;// rien
}
#endif
// dans le cas ou la matrice est uniquement un nombre !! on passe en CRAMER
// car cholesky est d'une part compliqué et pose pb !!
Enum_type_resolution_matri type_resolution_en_pratique = type_resolution; // init
if (Nb_ligne()==1)
type_resolution_en_pratique = CRAMER;
switch (type_resolution_en_pratique)
{case SYMETRISATION_PUIS_CHOLESKY:
Symetrisation(); // on commence par symétriser
case CHOLESKY:
{int N=Nb_ligne();
Vecteur resul(N); // Vecteur solution
Mat_pleine B(N,N); // B est une matrice triangulaire
// inferieure reguliere
Triangulation (N,B);
Resolution (N,B,b,resul);
return resul;
break;
}
case CRAMER:
{ Vecteur resul(b);Cramer(b,resul);
return resul;
break;
}
default:
{ Vecteur resul(b);
Mat_abstraite::Resolution_syst(b,resul,tole,maxi,rest);
return resul;
};
};
};
//2) avec en sortie le vecteur d'entree
#ifndef MISE_AU_POINT
inline
#endif
Vecteur& Mat_pleine::Resol_systID
(Vecteur& b,const double &tole, const int maxi, const int rest)
{
#ifdef MISE_AU_POINT
switch (type_resolution)
{case CHOLESKY:
// si c'est cholesky normal, il faut une matrice symétrique
if ( !(*this).Symetrie() )
{ if (ParaGlob::NiveauImpression() > 2)
{ cout << "\n warning : attention ne fonction qu'avec une matrice symetrique !"
<< " ce qui n'est pas le cas ici !!"
<< "\n on utilise ici le triangle inferieure pour la resolution ! \n";
cout << "MAT_PLEINE::RESOL_SYST ( Vecteur& )\n";
AfficheNonSymetries();
};
};
// de plus si c'est cramer ou cholesky après symétrisation, on doit avoir une matrice carrée
case CRAMER : case SYMETRISATION_PUIS_CHOLESKY :
if ( Nb_ligne()!=Nb_colonne() )
{ cout << "Erreur : la matrice doit etre carree !\n";
cout << "MAT_PLEINE::RESOL_SYST ( Vecteur& )\n";
Sortie(1);
};
break;
default: break;// rien
};
#endif
// dans le cas ou la matrice est uniquement un nombre !! on passe en CRAMER
// car cholesky est d'une part compliqué et pose pb !!
Enum_type_resolution_matri type_resolution_en_pratique = type_resolution; // init
if (Nb_ligne()==1)
type_resolution_en_pratique = CRAMER;
switch (type_resolution_en_pratique)
{case SYMETRISATION_PUIS_CHOLESKY:
Symetrisation(); // on commence par symétriser
case CHOLESKY:
{int N=Nb_ligne();
Mat_pleine B(N,N); // B est une matrice triangulaire
// inferieure reguliere
Triangulation (N,B);
Resolution (N,B,b,b);
break;
}
case CRAMER:
{ Cramer(b,b);
break;
}
default:
{ Vecteur bb(b);
Mat_abstraite::Resolution_syst(b,bb,tole,maxi,rest);
b = bb;
}
};
// retour
return b;
};
//3) avec en entrée un tableau de vecteur second membre et
// en sortie un nouveau tableau de vecteurs
#ifndef MISE_AU_POINT
inline
#endif
Tableau <Vecteur> Mat_pleine::Resol_syst
(const Tableau <Vecteur>& b,const double &tole, const int maxi, const int rest)
{
#ifdef MISE_AU_POINT
switch (type_resolution)
{case CHOLESKY:
// si c'est cholesky normal, il faut une matrice symétrique
if ( !(*this).Symetrie() )
{ if (ParaGlob::NiveauImpression() > 2)
{ cout << "\n warning : attention ne fonction qu'avec une matrice symetrique !"
<< " ce qui n'est pas le cas ici !!"
<< "\n on utilise ici le triangle inferieure pour la resolution ! \n";
cout << "MAT_PLEINE::RESOL_SYST ( const Tableau <Vecteur>..)\n";
AfficheNonSymetries();
};
};
// de plus si c'est cramer ou cholesky après symétrisation, on doit avoir une matrice carrée
case CRAMER : case SYMETRISATION_PUIS_CHOLESKY :
if ( Nb_ligne()!=Nb_colonne() )
{ cout << "Erreur : la matrice doit etre carree !\n";
cout << "MAT_PLEINE::RESOL_SYST ( const Tableau <Vecteur>..)\n";
Sortie(1);
};
break;
default: break;// rien
}
#endif
// dans le cas ou la matrice est uniquement un nombre !! on passe en CRAMER
// car cholesky est d'une part compliqué et pose pb !!
Enum_type_resolution_matri type_resolution_en_pratique = type_resolution; // init
if (Nb_ligne()==1)
type_resolution_en_pratique = CRAMER;
switch (type_resolution_en_pratique)
{case SYMETRISATION_PUIS_CHOLESKY:
Symetrisation(); // on commence par symétriser
case CHOLESKY:
{int N=Nb_ligne();
Mat_pleine B(N,N); // B est une matrice triangulaire
// inferieure reguliere
Triangulation (N,B);
int dimtab = b.Taille();
Tableau <Vecteur> c(dimtab);
for (int i=1;i<= dimtab;i++)
{ c(i).Change_taille(b(i).Taille());
Resolution (N,B,b(i),c(i));
};
return c;
break;
}
case CRAMER:
{int dimtab = b.Taille();
Tableau <Vecteur> c(dimtab);
Cramer(b,c);
return c;
break;
}
default:
{Tableau <Vecteur> bb(b);
Mat_abstraite::Resolution_syst(b,bb,tole,maxi,rest);
return bb;
}
}
};
//4) avec en entrée un tableau de vecteur second membre et
// en sortie le tableau de vecteurs d'entree
#ifndef MISE_AU_POINT
inline
#endif
Tableau <Vecteur>& Mat_pleine::Resol_systID
(Tableau <Vecteur>& b,const double &tole, const int maxi, const int rest)
{
#ifdef MISE_AU_POINT
switch (type_resolution)
{case CHOLESKY:
// si c'est cholesky normal, il faut une matrice symétrique
if ( !(*this).Symetrie() )
{ if (ParaGlob::NiveauImpression() > 2)
{ cout << "\n warning : attention ne fonction qu'avec une matrice symetrique !"
<< " ce qui n'est pas le cas ici !!"
<< "\n on utilise ici le triangle inferieure pour la resolution ! \n";
cout << "MAT_PLEINE::RESOL_SYST ( Tableau <Vecteur>..)\n";
AfficheNonSymetries();
};
};
// de plus si c'est cramer ou cholesky après symétrisation, on doit avoir une matrice carrée
case CRAMER : case SYMETRISATION_PUIS_CHOLESKY :
if ( Nb_ligne()!=Nb_colonne() )
{ cout << "Erreur : la matrice doit etre carree !\n";
cout << "MAT_PLEINE::RESOL_SYST ( Tableau <Vecteur>..)\n";
Sortie(1);
};
default: break; // sinon pas de pb ici
};
#endif
// dans le cas ou la matrice est uniquement un nombre !! on passe en CRAMER
// car cholesky est d'une part compliqué et pose pb !!
Enum_type_resolution_matri type_resolution_en_pratique = type_resolution; // init
if (Nb_ligne()==1)
type_resolution_en_pratique = CRAMER;
switch (type_resolution_en_pratique)
{case SYMETRISATION_PUIS_CHOLESKY:
Symetrisation(); // on commence par symétriser
case CHOLESKY:
{int N=Nb_ligne();
Mat_pleine B(N,N); // B est une matrice triangulaire
// inferieure reguliere
Triangulation (N,B);
int dimtab = b.Taille();
for (int i=1;i<= dimtab;i++)
Resolution (N,B,b(i),b(i));
break;
}
case CRAMER:
{Cramer (b,b);
break;
}
default:
{Tableau <Vecteur> bb(b);
Mat_abstraite::Resolution_syst(b,bb,tole,maxi,rest);
b = bb;
}
};
return b;
};
#ifndef MISE_AU_POINT
inline
#endif
Vecteur Mat_pleine::Prod_vec_mat (const Vecteur& vec) const
// Realise le produit du vecteur vec par la matrice
{ Vecteur result(Nb_colonne());
#ifdef MISE_AU_POINT
if (vec.Taille() != Colonne(1).Taille())
{cout << "erreur de taille, la première colonne de la matrice n'a pas la même dimension";
cout << " que le vecteur à multiplier "
<< "\n Mat_pleine::Prod_vec_mat (const Vecteur& vec)" << endl;
Sortie(1);
}
#endif
int nb_col = Nb_colonne();
for (int j=1;j<=nb_col;j++)
result(j)=Colonne(j)*vec;
return result;
};
#ifndef MISE_AU_POINT
inline
#endif
Vecteur& Mat_pleine::Prod_vec_mat (const Vecteur& vec,Vecteur& result) const
// Realise le produit du vecteur vec par la matrice
// ici on se sert du second vecteur pour le résultat
{ result.Change_taille(Nb_colonne());
#ifdef MISE_AU_POINT
if (vec.Taille() != Colonne(1).Taille())
{cout << "erreur de taille, la première colonne de la matrice n'a pas la même dimension";
cout << " que le vecteur à multiplier "
<< "\n Mat_pleine::Prod_vec_mat (const Vecteur& vec,Vecteur& res)" << endl;
Sortie(1);
}
#endif
int nb_col = Nb_colonne();
for (int j=1;j<=nb_col;j++)
result(j)=Colonne(j)*vec;
return result;
};
//5) avec en sortie le dernier vecteur d'entree, le premier étant le second membre
// et restant inchangé, en sortie c'est donc soit le retour ou soit vortie, les
// deux étant identiques
#ifndef MISE_AU_POINT
inline
#endif
Vecteur& Mat_pleine::Resol_systID_2 (const Vecteur& b,Vecteur& vortie
, const double &tole,const int maxi,const int rest)
{
#ifdef MISE_AU_POINT
switch (type_resolution)
{case CHOLESKY:
// si c'est cholesky normal, il faut une matrice symétrique
if ( !(*this).Symetrie() )
{ if (ParaGlob::NiveauImpression() > 2)
{ cout << "\n warning : attention ne fonction qu'avec une matrice symetrique !"
<< " ce qui n'est pas le cas ici !!"
<< "\n on utilise ici le triangle inferieure pour la resolution ! \n";
cout << "Mat_pleine::Resol_systID_2 (...\n";
AfficheNonSymetries();
};
};
// de plus si c'est cramer ou cholesky après symétrisation, on doit avoir une matrice carrée
case CRAMER : case SYMETRISATION_PUIS_CHOLESKY :
if ( Nb_ligne()!=Nb_colonne() )
{ cout << "Erreur : la matrice doit etre carree !\n";
cout << "Mat_pleine::Resol_systID_2 (...\n";
Sortie(1);
};
default: break; // sinon pas de pb ici
}
#endif
// dans le cas ou la matrice est uniquement un nombre !! on passe en CRAMER
// car cholesky est d'une part compliqué et pose pb !!
Enum_type_resolution_matri type_resolution_en_pratique = type_resolution; // init
if (Nb_ligne()==1)
type_resolution_en_pratique = CRAMER;
switch (type_resolution_en_pratique)
{case SYMETRISATION_PUIS_CHOLESKY:
Symetrisation(); // on commence par symétriser
case CHOLESKY:
{int N=Nb_ligne();
Mat_pleine B(*this); // B est une matrice triangulaire
// inferieure reguliere
Triangulation (N,B);
Resolution (N,B,b,vortie);
return vortie;
break;
}
case CRAMER:
{Cramer(b,vortie); return vortie; break;
}
default:
{Mat_abstraite::Resolution_syst(b,vortie,tole,maxi,rest);
return vortie;
}
};
};
// ===== RÉSOLUTION EN DEUX TEMPS ================ :
// 1) préparation de la matrice donc modification de la matrice éventuellement
// par exemple pour les matrices bandes avec cholesky : triangulation
#ifndef MISE_AU_POINT
inline
#endif
void Mat_pleine::Preparation_resol()
{
// dans le cas ou la matrice est uniquement un nombre !! on passe en CRAMER
// car cholesky est d'une part compliqué et pose pb !!
Enum_type_resolution_matri type_resolution_en_pratique = type_resolution; // init
if (Nb_ligne()==1)
type_resolution_en_pratique = CRAMER;
#ifdef MISE_AU_POINT
switch (type_resolution_en_pratique)
{case CHOLESKY:
2023-05-03 17:23:49 +02:00
{// si c'est cholesky normal, il faut une matrice symétrique
if ( !(*this).Symetrie() )
2021-09-28 17:28:23 +02:00
{ if (ParaGlob::NiveauImpression() > 2)
{ cout << "\n warning : attention ne fonction qu'avec une matrice symetrique !"
<< " ce qui n'est pas le cas ici !!"
<< "\n on utilise ici le triangle inferieure pour la resolution ! \n";
cout << "Mat_pleine::Preparation_resol() \n";
AfficheNonSymetries();
};
};
2023-05-03 17:23:49 +02:00
}
2021-09-28 17:28:23 +02:00
// de plus si c'est cramer ou cholesky après symétrisation, on doit avoir une matrice carrée
case CRAMER : case SYMETRISATION_PUIS_CHOLESKY :
2023-05-03 17:23:49 +02:00
{if ( Nb_ligne()!=Nb_colonne() )
2021-09-28 17:28:23 +02:00
{ cout << "Erreur : la matrice doit etre carree !\n";
cout << "Mat_pleine::Preparation_resol() \n";
Sortie(1);
};
2023-05-03 17:23:49 +02:00
}
2021-09-28 17:28:23 +02:00
default: break; // sinon ok ici
};
#endif
// dans le cas de cholesky on triangule
switch (type_resolution_en_pratique)
{case SYMETRISATION_PUIS_CHOLESKY:
Symetrisation(); // on commence par symétriser
case CHOLESKY:
{int N=Nb_ligne();
Triangulation (N,*this);
}
2023-05-03 17:23:49 +02:00
default: break; // sinon ok ici
};
2021-09-28 17:28:23 +02:00
};
// 2) *** résolution sans modification de la matrice DOIT ÊTRE PRÉCÉDÉ DE L'APPEL DE
// Preparation_resol
// a) avec en sortie un new vecteur
#ifndef MISE_AU_POINT
inline
#endif
Vecteur Mat_pleine::Simple_Resol_syst (const Vecteur& b,const double &tol
,const int maxi,const int rest) const
{
// dans le cas ou la matrice est uniquement un nombre !! on passe en CRAMER
// car cholesky est d'une part compliqué et pose pb !!
Enum_type_resolution_matri type_resolution_en_pratique = type_resolution; // init
if (Nb_ligne()==1)
type_resolution_en_pratique = CRAMER;
switch (type_resolution_en_pratique)
{case CHOLESKY: case SYMETRISATION_PUIS_CHOLESKY:
{int N=Nb_ligne();
Vecteur resul(N); // Vecteur solution
Resolution (N,*this,b,resul);
return resul; break;
}
case CRAMER:
{ Vecteur resul(b); Cramer(b,resul); return resul; break; }
default:
{ Vecteur resul(b);
Mat_abstraite::Resolution_syst(b,resul,tol,maxi,rest);
return resul;
}
};
};
// b) avec en sortie le vecteur d'entree
#ifndef MISE_AU_POINT
inline
#endif
Vecteur& Mat_pleine::Simple_Resol_systID (Vecteur& b,const double &tol
,const int maxi,const int restart ) const
{
// dans le cas ou la matrice est uniquement un nombre !! on passe en CRAMER
// car cholesky est d'une part compliqué et pose pb !!
Enum_type_resolution_matri type_resolution_en_pratique = type_resolution; // init
if (Nb_ligne()==1)
type_resolution_en_pratique = CRAMER;
switch (type_resolution_en_pratique)
{case CHOLESKY: case SYMETRISATION_PUIS_CHOLESKY:
{int N=Nb_ligne();
Resolution (N,*this,b,b);
}
case CRAMER:
{ Cramer(b,b);
break;
}
default:
{ Vecteur bb(b);
Mat_abstraite::Resolution_syst(b,bb,tol,maxi,restart);
b = bb;
}
};
return b;
};
// c) avec en entrée un tableau de vecteur second membre et
// en sortie un nouveau tableau de vecteurs
#ifndef MISE_AU_POINT
inline
#endif
Tableau <Vecteur> Mat_pleine::Simple_Resol_syst
(const Tableau <Vecteur>& b,const double &tol
,const int maxi,const int restart ) const
{
// dans le cas ou la matrice est uniquement un nombre !! on passe en CRAMER
// car cholesky est d'une part compliqué et pose pb !!
Enum_type_resolution_matri type_resolution_en_pratique = type_resolution; // init
if (Nb_ligne()==1)
type_resolution_en_pratique = CRAMER;
switch (type_resolution_en_pratique)
{case CHOLESKY: case SYMETRISATION_PUIS_CHOLESKY:
{int N=Nb_ligne();
int dimtab = b.Taille();
Tableau <Vecteur> c(dimtab);
for (int i=1;i<= dimtab;i++)
{ c(i).Change_taille(b(i).Taille());
Resolution (N,*this,b(i),c(i));
};
return c;
break;
}
case CRAMER:
{int dimtab = b.Taille();
Tableau <Vecteur> c(dimtab);
Cramer(b,c);
return c;
break;
}
default:
{ Tableau <Vecteur> bb(b);
Mat_abstraite::Resolution_syst(b,bb,tol,maxi,restart);
return bb;
}
};
};
// d) avec en entrée un tableau de vecteur second membre et
// en sortie le tableau de vecteurs d'entree
#ifndef MISE_AU_POINT
inline
#endif
Tableau <Vecteur>& Mat_pleine::Simple_Resol_systID
(Tableau <Vecteur>& b,const double &tol
,const int maxi,const int restart ) const
{
// dans le cas ou la matrice est uniquement un nombre !! on passe en CRAMER
// car cholesky est d'une part compliqué et pose pb !!
Enum_type_resolution_matri type_resolution_en_pratique = type_resolution; // init
if (Nb_ligne()==1)
type_resolution_en_pratique = CRAMER;
switch (type_resolution_en_pratique)
{case CHOLESKY: case SYMETRISATION_PUIS_CHOLESKY:
{int N=Nb_ligne();
int dimtab = b.Taille();
for (int i=1;i<= dimtab;i++)
Resolution (N,*this,b(i),b(i));
break;
}
case CRAMER:
{ Cramer (b,b);
break;
}
default:
{ Tableau <Vecteur> bb(b);
Mat_abstraite::Resolution_syst(b,bb,tol,maxi,restart);
b = bb;
}
}
return b;
};
// e) avec en sortie le dernier vecteur d'entree, le premier étant le second membre
// et restant inchangé, en sortie c'est donc soit le retour ou soit vortie, les
// deux étant identiques
#ifndef MISE_AU_POINT
inline
#endif
Vecteur& Mat_pleine::Simple_Resol_systID_2 (const Vecteur& b,Vecteur& vortie
, const double &tol,const int maxi,const int restart) const
{
// dans le cas ou la matrice est uniquement un nombre !! on passe en CRAMER
// car cholesky est d'une part compliqué et pose pb !!
Enum_type_resolution_matri type_resolution_en_pratique = type_resolution; // init
if (Nb_ligne()==1)
type_resolution_en_pratique = CRAMER;
switch (type_resolution_en_pratique)
{case CHOLESKY: case SYMETRISATION_PUIS_CHOLESKY:
{int N=Nb_ligne();
Resolution (N,*this,b,vortie);
return vortie;
break;
}
case CRAMER:
{ Cramer(b,vortie); return vortie; break;
}
default:
{ Mat_abstraite::Resolution_syst(b,vortie,tol,maxi,restart);
return vortie;
}
};
};
// ===== FIN RÉSOLUTION EN DEUX TEMPS ================ :
#ifndef MISE_AU_POINT
inline
#endif
Vecteur Mat_pleine::Prod_mat_vec (const Vecteur& vec) const
// Realise le produit de la matrice par le vecteur vec
{ Vecteur resul(Nb_ligne());
#ifdef MISE_AU_POINT
// il faut que les lignes et vec aient la même dimension
if (vec.Taille() != Ligne(1).Taille())
{cout << "erreur de taille, la première ligne de la matrice n'a pas la même dimension";
cout << " que le vecteur à multiplier "
<< "\n Mat_pleine::Prod_mat_vec (const Vecteur& vec)" << endl;
Sortie(1);
};
#endif
int nb_lign = val.Taille();
for (int i=1;i<=nb_lign;i++)
resul(i)=Ligne(i)*vec;
return resul;
};
#ifndef MISE_AU_POINT
inline
#endif
Vecteur& Mat_pleine::Prod_mat_vec (const Vecteur& vec,Vecteur& result) const
// Realise le produit de la matrice par le vecteur vec
// ici on se sert du second vecteur pour le résultat
{ result.Change_taille(Nb_ligne());
#ifdef MISE_AU_POINT
// il faut que les lignes et vec aient la même dimension
if (vec.Taille() != Ligne(1).Taille())
{cout << "erreur de taille, la première ligne de la matrice n'a pas la même dimension";
cout << " que le vecteur à multiplier "
<< "\n Mat_pleine::Prod_mat_vec (const Vecteur& vec,Vecteur& result)" << endl;
Sortie(1);
};
#endif
int nb_lign = val.Taille();
for (int i=1;i<=nb_lign;i++)
result(i)=Ligne(i)*vec;
return result;
};
#ifndef MISE_AU_POINT
inline
#endif
Mat_pleine Mat_pleine::Transpose () const
// Retourne la transposee de la matrice mat
{
2023-05-03 17:23:49 +02:00
int nb_col = Nb_colonne(); int nb_lign = val.Taille();
#ifdef MISE_AU_POINT
if ( ( nb_lign ==0 ) || ( nb_col ==0 ) )
2021-09-28 17:28:23 +02:00
{
cout << "Erreur : dimensions de la matrice non valables\n";
cout << "MAT_PLEINE::TRANSPOSE ( )\n";
Sortie(1);
};
#endif
2023-05-03 17:23:49 +02:00
2021-09-28 17:28:23 +02:00
Mat_pleine result(nb_col,nb_lign);
for (int i=1;i<=nb_col;i++)
{
for (int j=1;j<=nb_lign;j++)
{
result(i,j)=(*this)(j,i);
};
};
return result;
};
// determine l'inverse d'une matrice carre
// actuellement uniquement implemente pour dim = 1,2,3
#ifndef MISE_AU_POINT
inline
#endif
Mat_pleine Mat_pleine::Inverse() const
{ int iL = this->Nb_ligne(); int cC = this->Nb_colonne();
#ifdef MISE_AU_POINT
if (iL != cC)
{ cout << "\nErreur : le nombre de ligne est different du nb de colonne !\n";
cout << "Mat_pleine::Inverse() \n" << endl;
Sortie(1);
}
else if ((iL != 1) && (iL != 2) && (iL!=3))
{ cout << "\nErreur : desole, actuellement seule les cas dim = 1 ou 2 ou 3"
<< " sont implante !\n";
cout << "Mat_pleine::Inverse() \n" << endl;
Sortie(1);
}
#endif
const Mat_pleine & mat = *this;
Mat_pleine res (cC,cC);
Mat_pleine inter (cC,cC); // pour contenir la matrice régularisée
double det;
if (cC == 1)
{if ( Dabs(mat(1,1)) <= ConstMath::trespetit)
{ cout << "\nErreur1 : le determinant est nulle !\n";
if (ParaGlob::NiveauImpression() > 3)
mat.Affiche();
cout << "Mat_pleine::Inverse() \n" << endl;
throw ErrResolve_system_lineaire(1);
Sortie(1);
}
else
res(1,1) = 1./mat(1,1);
}
else if (cC == 2)
{ // on va regulariser la matrice pour limiter les erreurs de dépassement de capacité
// en espérant que cela fonctionne
// 1-- on récupère le maxi des valeurs de la matrice
double maxival = mat(1,1); // init
maxival = DabsMaX(maxival, mat(1,2));
maxival = DabsMaX(maxival, mat(2,1));
maxival = DabsMaX(maxival, mat(2,2));
if ( maxival <= ConstMath::trespetit)
{ cout << "\nErreur11 : la matrice est nulle !\n";
if (ParaGlob::NiveauImpression() > 3)
mat.Affiche();
cout << "Mat_pleine::Inverse() \n" << endl;
throw ErrResolve_system_lineaire(1);
Sortie(1);
};
// 2 --- on régularise
double unSurMaxival = 1./maxival;
inter (1,1) = mat(1,1) * unSurMaxival;
inter (1,2) = mat(1,2) * unSurMaxival;
inter (2,1) = mat(2,1) * unSurMaxival;
inter (2,2) = mat(2,2) * unSurMaxival;
// suite calcul normal
det = inter(1,1) * inter(2,2) - inter(1,2) * inter (2,1);
if ( Dabs(det) <= ConstMath::trespetit)
{ cout << "\nErreur2 : le determinant est nulle !\n";
if (ParaGlob::NiveauImpression() > 3)
mat.Affiche();
cout << "Mat_pleine::Inverse() \n" << endl;
throw ErrResolve_system_lineaire(1);
Sortie(1);
};
// 3-- on tient compte de la régularisation
res (1,1) = unSurMaxival * inter(2,2) / det;
res(1,2) = - unSurMaxival * inter(1,2) / det;
res(2,1) = - unSurMaxival * inter(2,1) / det;
res(2,2) = unSurMaxival * inter(1,1) / det;
}
else if (cC == 3)
{ // on va regulariser la matrice pour limiter les erreurs de dépassement de capacité
// en espérant que cela fonctionne
// 1-- on récupère le maxi des valeurs de la matrice
double maxival = mat(1,1); // init
maxival = DabsMaX(maxival, mat(1,2));
maxival = DabsMaX(maxival, mat(1,3));
maxival = DabsMaX(maxival, mat(2,1));
maxival = DabsMaX(maxival, mat(2,2));
maxival = DabsMaX(maxival, mat(2,3));
maxival = DabsMaX(maxival, mat(3,1));
maxival = DabsMaX(maxival, mat(3,2));
maxival = DabsMaX(maxival, mat(3,3));
if ( maxival <= ConstMath::trespetit)
{ cout << "\nErreur12 : la matrice est nulle !\n";
if (ParaGlob::NiveauImpression() > 3)
mat.Affiche();
cout << "Mat_pleine::Inverse() \n" << endl;
throw ErrResolve_system_lineaire(1);
Sortie(1);
};
// 2 --- on régularise
double unSurMaxival = 1./maxival;
inter (1,1) = mat(1,1) * unSurMaxival;
inter (1,2) = mat(1,2) * unSurMaxival;
inter (1,3) = mat(1,3) * unSurMaxival;
inter (2,1) = mat(2,1) * unSurMaxival;
inter (2,2) = mat(2,2) * unSurMaxival;
inter (2,3) = mat(2,3) * unSurMaxival;
inter (3,1) = mat(3,1) * unSurMaxival;
inter (3,2) = mat(3,2) * unSurMaxival;
inter (3,3) = mat(3,3) * unSurMaxival;
//det = ((a11 a22 - a12 a21) a33 + (a13 a21 - a11 a23) a32 + (a12 a23 - a13 a22) a31);
det = inter(1,1) *( inter(2,2) * inter(3,3) - inter(2,3) * inter (3,2));
det -= inter(1,2) *( inter(2,1) * inter(3,3) - inter(2,3) * inter (3,1));
det += inter(1,3) *( inter(2,1) * inter(3,2) - inter(2,2) * inter (3,1));
if ( Dabs(det) <= ConstMath::trespetit)
{ cout << "\nErreur3 : le determinant est nulle !\n";
if (ParaGlob::NiveauImpression() > 3)
mat.Affiche();
cout << "Mat_pleine::Inverse() \n" << endl;
throw ErrResolve_system_lineaire(1);
Sortie(1);
};
// 3-- on tient compte de la régularisation
res(1,1) = unSurMaxival * (inter(2,2) * inter(3,3) - inter(2,3) * inter (3,2)) / det;
res(2,2) = unSurMaxival * (inter(1,1) * inter(3,3) - inter(1,3) * inter (3,1)) / det;
res(3,3) = unSurMaxival * (inter(1,1) * inter(2,2) - inter(1,2) * inter (2,1)) / det;
res(1,2) = -unSurMaxival * (inter(1,2) * inter(3,3) - inter(1,3) * inter (3,2)) / det;
res(1,3) = unSurMaxival * (inter(1,2) * inter(2,3) - inter(1,3) * inter (2,2)) / det;
res(2,1) = -unSurMaxival * (inter(2,1) * inter(3,3) - inter(2,3) * inter (3,1)) / det;
res(3,1) = unSurMaxival * (inter(2,1) * inter(3,2) - inter(2,2) * inter (3,1)) / det;
res(3,2) = -unSurMaxival * (inter(1,1) * inter(3,2) - inter(1,2) * inter (3,1)) / det;
res(2,3) = -unSurMaxival * (inter(1,1) * inter(2,3) - inter(1,3) * inter (2,1)) / det;
};
return res;
};
#ifndef MISE_AU_POINT
inline
#endif
// determine l'inverse d'une matrice carre
// actuellement uniquement implemente pour dim = 1,2,3
// ici retourne l'inverse dans la matrice passée en paramètre
// qui doit être de même type et de même dimension que this
Mat_pleine& Mat_pleine::Inverse(Mat_pleine& res) const
{ int iL = this->Nb_ligne(); int cC = this->Nb_colonne();
#ifdef MISE_AU_POINT
if (iL != cC)
{ cout << "\nErreur : le nombre de ligne est different du nb de colonne !\n";
cout << "Mat_pleine::Inverse() \n" << endl;
Sortie(1);
}
else if ((iL != 1) && (iL != 2) && (iL!=3))
{ cout << "\nErreur : desole, actuellement seule les cas dim = 1 ou 2 ou 3"
<< " sont implante !\n";
cout << "Mat_pleine::Inverse() \n" << endl;
Sortie(1);
};
if (iL != res.Nb_ligne())
{ cout << "\nErreur : le nombre de ligne est different entre this " << iL
<< " et res " << res.Nb_ligne() << "!\n";
cout << "Mat_pleine::Inverse() \n" << endl;
Sortie(1);
}
else if (cC != res.Nb_colonne())
{ cout << "\nErreur : le nombre de colonne est different entre this " << cC
<< " et res " << res.Nb_colonne() << "!\n";
cout << "Mat_pleine::Inverse() \n" << endl;
Sortie(1);
};
#endif
const Mat_pleine & mat = *this;
Mat_pleine inter (cC,cC); // pour contenir la matrice régularisée
double det;
if (cC == 1)
{if ( Dabs(mat(1,1)) <= ConstMath::trespetit)
{ cout << "\nErreur : le determinant est nulle !\n";
cout << "Mat_pleine::Inverse() \n" << endl;
throw ErrResolve_system_lineaire(1);
Sortie(1);
}
else
res(1,1) = 1./mat(1,1);
}
else if (cC == 2)
{ // on va regulariser la matrice pour limiter les erreurs de dépassement de capacité
// en espérant que cela fonctionne
// 1-- on récupère le maxi des valeurs de la matrice
double maxival = mat(1,1); // init
maxival = DabsMaX(maxival, mat(1,2));
maxival = DabsMaX(maxival, mat(2,1));
maxival = DabsMaX(maxival, mat(2,2));
if ( maxival <= ConstMath::trespetit)
{ cout << "\nErreur110 : la matrice est nulle !\n";
if (ParaGlob::NiveauImpression() > 3)
mat.Affiche();
cout << "Mat_pleine::Inverse() \n" << endl;
throw ErrResolve_system_lineaire(1);
Sortie(1);
};
// 2 --- on régularise
double unSurMaxival = 1./maxival;
inter (1,1) = mat(1,1) * unSurMaxival;
inter (1,2) = mat(1,2) * unSurMaxival;
inter (2,1) = mat(2,1) * unSurMaxival;
inter (2,2) = mat(2,2) * unSurMaxival;
// suite calcul normal
det = inter(1,1) * inter(2,2) - inter(1,2) * inter (2,1);
if ( Dabs(det) <= ConstMath::trespetit)
{ cout << "\nErreur : le determinant est nulle !\n";
cout << "Mat_pleine::Inverse() \n" << endl;
throw ErrResolve_system_lineaire(1);
Sortie(1);
};
// 3-- on tient compte de la régularisation
res (1,1) = unSurMaxival * inter(2,2) / det;
res(1,2) = - unSurMaxival * inter(1,2) / det;
res(2,1) = - unSurMaxival * inter(2,1) / det;
res(2,2) = unSurMaxival * inter(1,1) / det;
// res (1,1) = mat(2,2) / det;
// res(1,2) = -mat(1,2) / det;
// res(2,1) = -mat(2,1) / det;
// res(2,2) = mat(1,1) / det;
}
else if (cC == 3)
{ // on va regulariser la matrice pour limiter les erreurs de dépassement de capacité
// en espérant que cela fonctionne
// 1-- on récupère le maxi des valeurs de la matrice
double maxival = mat(1,1); // init
maxival = DabsMaX(maxival, mat(1,2));
maxival = DabsMaX(maxival, mat(1,3));
maxival = DabsMaX(maxival, mat(2,1));
maxival = DabsMaX(maxival, mat(2,2));
maxival = DabsMaX(maxival, mat(2,3));
maxival = DabsMaX(maxival, mat(3,1));
maxival = DabsMaX(maxival, mat(3,2));
maxival = DabsMaX(maxival, mat(3,3));
if ( maxival <= ConstMath::trespetit)
{ cout << "\nErreur12 : la matrice est nulle !\n";
if (ParaGlob::NiveauImpression() > 3)
mat.Affiche();
cout << "Mat_pleine::Inverse() \n" << endl;
throw ErrResolve_system_lineaire(1);
Sortie(1);
};
// 2 --- on régularise
double unSurMaxival = 1./maxival;
inter (1,1) = mat(1,1) * unSurMaxival;
inter (1,2) = mat(1,2) * unSurMaxival;
inter (1,3) = mat(1,3) * unSurMaxival;
inter (2,1) = mat(2,1) * unSurMaxival;
inter (2,2) = mat(2,2) * unSurMaxival;
inter (2,3) = mat(2,3) * unSurMaxival;
inter (3,1) = mat(3,1) * unSurMaxival;
inter (3,2) = mat(3,2) * unSurMaxival;
inter (3,3) = mat(3,3) * unSurMaxival;
//det = ((a11 a22 - a12 a21) a33 + (a13 a21 - a11 a23) a32 + (a12 a23 - a13 a22) a31);
det = inter(1,1) *( inter(2,2) * inter(3,3) - inter(2,3) * inter (3,2));
det -= inter(1,2) *( inter(2,1) * inter(3,3) - inter(2,3) * inter (3,1));
det += inter(1,3) *( inter(2,1) * inter(3,2) - inter(2,2) * inter (3,1));
if ( Dabs(det) <= ConstMath::trespetit)
{ cout << "\nErreur : le determinant est nulle !\n";
cout << "Mat_pleine::Inverse() \n" << endl;
throw ErrResolve_system_lineaire(1);
Sortie(1);
};
// 3-- on tient compte de la régularisation
res(1,1) = unSurMaxival * (inter(2,2) * inter(3,3) - inter(2,3) * inter (3,2)) / det;
res(2,2) = unSurMaxival * (inter(1,1) * inter(3,3) - inter(1,3) * inter (3,1)) / det;
res(3,3) = unSurMaxival * (inter(1,1) * inter(2,2) - inter(1,2) * inter (2,1)) / det;
res(1,2) = -unSurMaxival * (inter(1,2) * inter(3,3) - inter(1,3) * inter (3,2)) / det;
res(1,3) = unSurMaxival * (inter(1,2) * inter(2,3) - inter(1,3) * inter (2,2)) / det;
res(2,1) = -unSurMaxival * (inter(2,1) * inter(3,3) - inter(2,3) * inter (3,1)) / det;
res(3,1) = unSurMaxival * (inter(2,1) * inter(3,2) - inter(2,2) * inter (3,1)) / det;
res(3,2) = -unSurMaxival * (inter(1,1) * inter(3,2) - inter(1,2) * inter (3,1)) / det;
res(2,3) = -unSurMaxival * (inter(1,1) * inter(2,3) - inter(1,3) * inter (2,1)) / det;
// res(1,1) = (mat(2,2) * mat(3,3) - mat(2,3) * mat (3,2)) / det;
// res(2,2) = (mat(1,1) * mat(3,3) - mat(1,3) * mat (3,1)) / det;
// res(3,3) = (mat(1,1) * mat(2,2) - mat(1,2) * mat (2,1)) / det;
// res(1,2) = -(mat(1,2) * mat(3,3) - mat(1,3) * mat (3,2)) / det;
// res(1,3) = (mat(1,2) * mat(2,3) - mat(1,3) * mat (2,2)) / det;
// res(2,1) = -(mat(2,1) * mat(3,3) - mat(2,3) * mat (3,1)) / det;
// res(3,1) = (mat(2,1) * mat(3,2) - mat(2,2) * mat (3,1)) / det;
// res(3,2) = -(mat(1,1) * mat(3,2) - mat(1,2) * mat (3,1)) / det;
// res(2,3) = -(mat(1,1) * mat(2,3) - mat(1,3) * mat (2,1)) / det;
};
return res;
};
#ifndef MISE_AU_POINT
inline
#endif
// détermine le déterminant d'une matrice
// actuellement uniquement implemente pour dim = 1,2,3
double Mat_pleine::Determinant() const
{ int iL = this->Nb_ligne(); int cC = this->Nb_colonne();
#ifdef MISE_AU_POINT
if (iL != cC)
{ cout << "\nErreur : le nombre de ligne est different du nb de colonne !\n";
cout << "Mat_pleine::Determinant() \n" << endl;
Sortie(1);
}
else if ((iL != 1) && (iL != 2) && (iL!=3))
{ cout << "\nErreur : desole, actuellement seule les cas dim = 1 ou 2 ou 3"
<< " sont implante !\n";
cout << "Mat_pleine::Determinant() \n" << endl;
Sortie(1);
}
#endif
const Mat_pleine & mat = *this;
double det;
if (cC == 1)
{det = mat(1,1);
}
else if (cC == 2)
{det = mat(1,1) * mat(2,2) - mat(1,2) * mat (2,1);
}
else if (cC == 3)
{ //det = ((a11 a22 - a12 a21) a33 + (a13 a21 - a11 a23) a32 + (a12 a23 - a13 a22) a31);
det = mat(1,1) *( mat(2,2) * mat(3,3) - mat(2,3) * mat (3,2));
det -= mat(1,2) *( mat(2,1) * mat(3,3) - mat(2,3) * mat (3,1));
det += mat(1,3) *( mat(2,1) * mat(3,2) - mat(2,2) * mat (3,1));
};
return det;
};
#ifndef MISE_AU_POINT
inline
#endif
Mat_pleine Mat_pleine::operator+ (const Mat_pleine& mat_pl) const
// Surcharge de l'operateur + : addition de deux matrices
2023-05-03 17:23:49 +02:00
{ int nb_lign = val.Taille();
int nb_col = Nb_colonne();
#ifdef MISE_AU_POINT
if ( ( nb_lign != mat_pl.Nb_ligne() )
|| ( nb_col != mat_pl.Nb_colonne() ) )
2021-09-28 17:28:23 +02:00
{ cout << "Erreur : dimensions des matrices non egales\n";
cout << "MAT_PLEINE::OPERATOR+ (const Mat_pleine& )\n";
Sortie(1);
};
#endif
2023-05-03 17:23:49 +02:00
Mat_pleine result(nb_lign,nb_col);
result.val_de_base = val_de_base + mat_pl.val_de_base;
// result.val.Change_taille(nb_lign);
// for (int i=1;i<=nb_lign;i++)
// // addition des composantes a l'aide de la surcharge de l'operateur
// // + pour les vecteurs
// result.val(i)=val(i)+mat_pl(i);
2021-09-28 17:28:23 +02:00
return result;
};
#ifndef MISE_AU_POINT
inline
#endif
Mat_pleine Mat_pleine::operator- (const Mat_pleine& mat_pl) const
// Surcharge de l'operateur - : soustraction entre deux matrices
{ int nb_lign = val.Taille();
2023-05-03 17:23:49 +02:00
int nb_col = Nb_colonne();
#ifdef MISE_AU_POINT
if ( ( nb_lign != mat_pl.Nb_ligne() )
|| ( nb_col != mat_pl.Nb_colonne() ) )
2021-09-28 17:28:23 +02:00
{
cout << "Erreur : dimensions des matrices non egales\n";
cout << "MAT_PLEINE::OPERATOR- (const Mat_pleine& )\n";
Sortie(1);
};
#endif
2023-05-03 17:23:49 +02:00
Mat_pleine result(nb_lign,nb_col);
result.val_de_base = val_de_base - mat_pl.val_de_base;
//
// result.val.Change_taille(nb_lign);
// for (int i=1;i<=nb_lign;i++)
// // soustraction des composantes a l'aide de la surcharge de l'operateur
// // - pour les vecteurs
// result.val(i)=val(i)-mat_pl(i);
2021-09-28 17:28:23 +02:00
return result;
};
#ifndef MISE_AU_POINT
inline
#endif
Mat_pleine Mat_pleine::operator- () const
// Surcharge de l'operateur - : oppose d'une matrice
2023-05-03 17:23:49 +02:00
{ Mat_pleine result(this->Nb_ligne(),this->Nb_colonne());
result.val_de_base = -val_de_base;
// int nb_lign = val.Taille();
// result.val.Change_taille(nb_lign);
// for (int i=1;i<=nb_lign;i++)
// // determination de l'oppose des composantes de la matrice courante a l'aide de la surcharge
// // de l'operateur -= pour les vecteurs
// result.val(i)=-val(i);
2021-09-28 17:28:23 +02:00
return result;
};
#ifndef MISE_AU_POINT
inline
#endif
void Mat_pleine::operator+= ( const Mat_pleine& mat_pl)
// Addition a la matrice courante de la matrice mat_pl
{ int nb_lign = val.Taille();
#ifdef MISE_AU_POINT
if ( ( nb_lign!=mat_pl.Nb_ligne() )
|| ( Nb_colonne()!=mat_pl.Nb_colonne() ) )
{
cout << "Erreur : dimensions des matrices non egales\n";
cout << "MAT_PLEINE::OPERATOR+= (const Mat_pleine& )\n";
Sortie(1);
};
#endif
2023-05-03 17:23:49 +02:00
val_de_base += mat_pl.val_de_base;
//
// for (int i=1;i<=nb_lign;i++)
// // l'operateur += est l'operateur += des vecteurs
// val(i)+=mat_pl(i);
2021-09-28 17:28:23 +02:00
};
#ifndef MISE_AU_POINT
inline
#endif
void Mat_pleine::operator-= ( const Mat_pleine& mat_pl)
{ // Soustrait a la matrice courante la matrice mat_pl
int nb_lign = val.Taille();
#ifdef MISE_AU_POINT
if ( ( nb_lign!=mat_pl.Nb_ligne() )
|| ( Nb_colonne()!=mat_pl.Nb_colonne() ) )
{
cout << "Erreur : dimensions des matrices non egales\n";
cout << "MAT_PLEINE::OPERATOR-= (const Mat_pleine& )\n";
Sortie(1);
};
#endif
2023-05-03 17:23:49 +02:00
val_de_base -= mat_pl.val_de_base;
// for (int i=1;i<=nb_lign;i++)
// // l'operateur -= est l'operateur -= des vecteurs
// val(i)-=mat_pl(i);
2021-09-28 17:28:23 +02:00
};
#ifndef MISE_AU_POINT
inline
#endif
Mat_pleine Mat_pleine::operator* (const Mat_pleine& mat_pl) const
// Surcharge de l'operateur * : multiplication de deux matrices
{ int mat_nb_col = mat_pl.Nb_colonne(); int nb_lign = val.Taille();
#ifdef MISE_AU_POINT
if ( Nb_colonne()!=mat_pl.Nb_ligne() )
{
cout << "Erreur : dimensions des matrices non valables\n";
cout << "MAT_PLEINE::OPERATOR* (const Mat_pleine& )\n";
Sortie(1);
};
#endif
2023-05-03 17:23:49 +02:00
Mat_pleine result(nb_lign,mat_nb_col);
// result.Initialise(nb_lign,mat_nb_col);
2021-09-28 17:28:23 +02:00
for (int i=1;i<=nb_lign;i++)
{
for (int j=1;j<=mat_nb_col;j++)
// la ieme jieme composante de la matrice result est obtenue par le produit de la
// ieme ligne de la matrice courante avec la jieme colonne de la matrice mat_pl
result(i,j)=Ligne(i)*mat_pl.Colonne(j);
};
return result;
};
#ifndef MISE_AU_POINT
inline
#endif
void Mat_pleine::operator*= (double coeff)
//multiplication de la matrice courante par un scalaire
2023-05-03 17:23:49 +02:00
{ val_de_base *= coeff;
//
// int nb_lign = val.Taille();
// for (int i=1;i<=nb_lign;i++)
// // l'operateur *= est l'operateur *= des vecteurs
// val(i) *= coeff;
2021-09-28 17:28:23 +02:00
};
#ifndef MISE_AU_POINT
inline
#endif
Mat_pleine Mat_pleine::operator* (double coeff) const
// Surcharge de l'operateur * : multiplication d'une matrice par un scalaire
2023-05-03 17:23:49 +02:00
{ Mat_pleine result(*this);
result *= coeff;
return result;
//
//
//
//Mat_pleine result;
// int nb_lign = val.Taille();int nb_col = Nb_colonne();
// result.Initialise(nb_lign,nb_col);
// for (int i=1;i<=nb_lign;i++)
// {
// for (int j=1;j<=nb_col;j++)
// result(i,j)=coeff*(*this)(i,j);
// };
// return result;
2021-09-28 17:28:23 +02:00
};
#ifndef MISE_AU_POINT
inline
#endif
int Mat_pleine::operator== (const Mat_pleine& mat_pl) const
// Surcharge de l'operateur == : test d'egalite entre deux matrices
// Retourne 1 si les deux matrices sont egales
// 0 sinon
{
int nb_lign = val.Taille();int nb_col = Nb_colonne();
if ( ( nb_lign != mat_pl.Nb_ligne() ) || ( nb_col != mat_pl.Nb_colonne() ) )
return 0;
2023-05-03 17:23:49 +02:00
int taille_glob = val_de_base.Taille();
for (int i=1;i<=taille_glob;i++)
{if (val_de_base(i) != mat_pl.val_de_base(i))
return 0;
}
//
// for (int i=1;i <= nb_lign;i++)
// {for (int j=1;j <= nb_col;j++)
// {if ( (*this)(i,j) != mat_pl(i,j) )
// return 0;
// };
// };
2021-09-28 17:28:23 +02:00
return 1;
};
// remplace la ligne de la matrice par la ligne fournie
#ifndef MISE_AU_POINT
inline
#endif
void Mat_pleine::RemplaceLigneSpe(int i,const Vecteur & v)
{
#ifdef MISE_AU_POINT
2023-05-03 17:23:49 +02:00
if (val(i)->Taille() != v.Taille())
2021-09-28 17:28:23 +02:00
{ cout << " \nla taille de la ligne i= " << i <<" a remplacer n est pas bonne " ;
2023-05-03 17:23:49 +02:00
cout << " taille = " << v.Taille() << " tailleLigneMatrice = " << val(i)->Taille();
2021-09-28 17:28:23 +02:00
cout << " void Mat_pleine::RemplaceLigne(int i,Vecteur & v)" << endl;
Sortie (1);
}
#endif
2023-05-03 17:23:49 +02:00
*(val(i)) = v;
2021-09-28 17:28:23 +02:00
};
//met une valeur identique sur toute la ligne
#ifndef MISE_AU_POINT
inline
#endif
void Mat_pleine::MetValLigne(int i,double x)
2023-05-03 17:23:49 +02:00
{ val(i)->Inita(x);
// Vecteur& v=*(val(i));
// int vali_taille = val(i)->Taille();
// for (int j=1;j <= vali_taille;j++)
// v(j) = x;
2021-09-28 17:28:23 +02:00
};
// remplace la Colonne de la matrice par la colonne fournie
#ifndef MISE_AU_POINT
inline
#endif
void Mat_pleine::RemplaceColonneSpe(int j, const Vecteur & v)
{
#ifdef MISE_AU_POINT
if (val.Taille() != v.Taille())
{ cout << " \nla taille de la colonne j= " << j <<" a remplacer n est pas bonne " ;
cout << " taille = " << v.Taille() << " tailleColMatrice = " << val.Taille();
cout << " void Mat_pleine::RemplaceColonne(int j,Vecteur & v) " << endl;
Sortie (1);
}
#endif
int nb_lign = val.Taille();
for (int i=1;i <= nb_lign;i++)
(*this)(i,j)=v(i);
};
//met une valeur identique sur toute la colonne
#ifndef MISE_AU_POINT
inline
#endif
void Mat_pleine::MetValColonne(int j,double y)
{ int nb_lign = val.Taille();
for (int i=1;i <= nb_lign;i++)
(*this)(i,j) = y;
};
// Multiplication d'une ligne iligne de la matrice avec un vecteur de
// dimension = le nombre de colonne de la matrice
#ifndef MISE_AU_POINT
inline
#endif
double Mat_pleine::Prod_Ligne_vec ( int iligne, const Vecteur& vec) const
{
#ifdef MISE_AU_POINT
2023-05-03 17:23:49 +02:00
if (val(iligne)->Taille() != vec.Taille())
2021-09-28 17:28:23 +02:00
{ cout << " \n erreur la taille de la ligne i= " << iligne <<" n est pas correcte ";
2023-05-03 17:23:49 +02:00
cout << " tailleligneMatrice = " <<val(iligne)->Taille()
2021-09-28 17:28:23 +02:00
<< " taillevecteur = " << vec.Taille();
cout << " double Mat_pleine::Prod_Ligne_vec ( int iligne,Vecteur& vec) " << endl;
Sortie (1);
}
#endif
double res = 0;
2023-05-03 17:23:49 +02:00
int val_iligne_taille = val(iligne)->Taille();
2021-09-28 17:28:23 +02:00
for (int j=1;j <= val_iligne_taille;j++)
res += (*this)(iligne,j) * vec(j);
return res;
};
// Multiplication d'un vecteur avec une colonne icol de la matrice
// dimension = le nombre de ligne de la matrice
#ifndef MISE_AU_POINT
inline
#endif
double Mat_pleine::Prod_vec_col( int jcol, const Vecteur& vec) const
{
#ifdef MISE_AU_POINT
if (val.Taille() != vec.Taille())
{ cout << " \n erreur la taille de la colonne j= " << jcol <<" n est pas correcte ";
cout << " tailleColMatrice = " <<val.Taille()
<< " taillevecteur = " << vec.Taille();
cout << " double Mat_pleine::Prod_vec_col( int jcol,Vecteur& vec) " << endl;
Sortie (1);
}
#endif
double res = 0;
int val_taille = val.Taille();
for (int i=1;i <= val_taille;i++)
res += vec(i) * (*this)(i,jcol) ;
return res;
};
// calcul du produit : (vec_1)^T * A * (vect_2)
#ifndef MISE_AU_POINT
inline
#endif
double Mat_pleine::vectT_mat_vec(const Vecteur& vec1, const Vecteur& vec2) const
{ double resu=0.; // valeur de retour
#ifdef MISE_AU_POINT
2023-05-03 17:23:49 +02:00
if (vec1.Taille() != val(1)->Taille())
2021-09-28 17:28:23 +02:00
{cout << "erreur de taille, la dimension de (vec1)^T= " << vec1.Taille() << " n'a pas la même dimension";
2023-05-03 17:23:49 +02:00
cout << " que le le nombre de ligne de la matrice= " << val(1)->Taille()
2021-09-28 17:28:23 +02:00
<< "\n Mat_pleine::vectT_mat_vec(const Vecteur& vec1, const Vecteur& vec2)" << endl;
Sortie(1);
}
if (vec2.Taille() != val.Taille())
{cout << "erreur de taille, la dimension de (vec2)= " << vec2.Taille() << " n'a pas la même dimension";
cout << " que le le nombre de colonne de la matrice= " << val.Taille()
<< "\n Mat_pleine::vectT_mat_vec(const Vecteur& vec1, const Vecteur& vec2)" << endl;
Sortie(1);
}
#endif
int nb_col = Nb_colonne();
int nb_lig = Nb_ligne();
for (int i=1;i<=nb_col ;i++)
for (int j=1;j<=nb_lig ;j++)
resu += vec1(i) * (*this)(i,j) * vec2(j);
return resu;
};
// ramène un tableau de coordonnées (avec 3 variances possibles) correspondant à la matrice
// tab(i) = la colonne i de la matrice
#ifndef MISE_AU_POINT
inline
#endif
Tableau <Coordonnee > Mat_pleine::Coordonnee_Base_associee() const
{ int nb_col = Nb_colonne();
Tableau <Coordonnee > tab(nb_col);
int nb_lig = Nb_ligne();
for (int i=1;i<=nb_col ;i++)
{tab(i).Change_dim(nb_lig);
for (int j=1;j<=nb_lig ;j++)
tab(i)(j) = (*this)(j,i);
};
return tab;
};
#ifndef MISE_AU_POINT
inline
#endif
Tableau <CoordonneeH > Mat_pleine::CoordonneeH_Base_associee() const
{ int nb_col = Nb_colonne();
Tableau <CoordonneeH > tab(nb_col);
int nb_lig = Nb_ligne();
for (int i=1;i<=nb_col ;i++)
{tab(i).Change_dim(nb_lig);
for (int j=1;j<=nb_lig ;j++)
tab(i)(j) = (*this)(j,i);
};
return tab;
};
#ifndef MISE_AU_POINT
inline
#endif
Tableau <CoordonneeB > Mat_pleine::CoordonneeB_Base_associee() const
{ int nb_col = Nb_colonne();
Tableau <CoordonneeB > tab(nb_col);
int nb_lig = Nb_ligne();
for (int i=1;i<=nb_col ;i++)
{tab(i).Change_dim(nb_lig);
for (int j=1;j<=nb_lig ;j++)
tab(i)(j) = (*this)(j,i);
};
return tab;
};
#ifndef MISE_AU_POINT
inline
#endif
// idem mais une seule colonne i
Coordonnee Mat_pleine::Coordonnee_Base_associee(int i) const
{ Coordonnee coo;
int nb_lig = Nb_ligne();
coo.Change_dim(nb_lig);
for (int j=1;j<=nb_lig ;j++)
coo(j) = (*this)(j,i);
return coo;
};
#ifndef MISE_AU_POINT
inline
#endif
CoordonneeH Mat_pleine::CoordonneeH_Base_associee(int i) const
{ CoordonneeH coo;
int nb_lig = Nb_ligne();
coo.Change_dim(nb_lig);
for (int j=1;j<=nb_lig ;j++)
coo(j) = (*this)(j,i);
return coo;
};
#ifndef MISE_AU_POINT
inline
#endif
CoordonneeB Mat_pleine::CoordonneeB_Base_associee(int i) const
{ CoordonneeB coo;
int nb_lig = Nb_ligne();
coo.Change_dim(nb_lig);
for (int j=1;j<=nb_lig ;j++)
coo(j) = (*this)(j,i);
return coo;
};
// ramène le maxi en valeur absolue des valeurs de la matrice, et les indices associées
#ifndef MISE_AU_POINT
inline
#endif
double Mat_pleine::MaxiValAbs(int & imax, int & jmax) const
{ double lemaxi=0.;
imax = jmax = 0;
int nb_lig = val.Taille();
int nb_col = Nb_colonne();
for (int i=1;i<= nb_lig;i++)
for (int j=1;j<= nb_col;j++)
2023-05-03 17:23:49 +02:00
{double valeur = Dabs((*this)(i,j));
if (valeur > lemaxi)
{ lemaxi = valeur; imax = i; jmax = j; };
2021-09-28 17:28:23 +02:00
};
return lemaxi;
};
// ramène le maxi en valeur absolue de la somme des valeurs absolues de chaque ligne,
// et l'indice de ligne associé
// permet d'avoir une approximation de la valeur propre maximale de la matrice via le
// théorème de Gerschgorin : |lambda| < Max_i (|lambda_i|)
// avec: |lambda_i| < Max_j somme_j |k_ij|)
#ifndef MISE_AU_POINT
inline
#endif
double Mat_pleine::Maxi_ligne_ValAbs(int & imax) const
{ double lemaxi=0.;
imax = 0;
int nb_lig = val.Taille();
for (int i=1;i<= nb_lig;i++)
2023-05-03 17:23:49 +02:00
{double somme_ligne = val(i)->SommeAbs();
2021-09-28 17:28:23 +02:00
if (somme_ligne > lemaxi)
{ lemaxi = somme_ligne; imax = i;};
};
return lemaxi;
};
// =================== méthodes protégées =======================
// symétrisation de la matrice (il faut qu'elle soit carrée
#ifndef MISE_AU_POINT
inline
#endif
void Mat_pleine::Symetrisation()
{// arrivée ici on considère que l'on a une matrice carrée
int N=Nb_ligne();
for (int i=1;i<= N;i++)
for (int j=i+1;j<=N;j++)
(*this)(i,j)=(*this)(j,i)=0.5*((*this)(i,j)+(*this)(j,i));
};
// calcul de la matrice triangulée dans le cadre de la méthode de cholesky
// utilisation particulière : la matrice résultat B est défini dans le programme
// appelant.
// On considère que la matrice B est symétrique, et que l'on utilise uniquement la partie supérieure
#ifndef MISE_AU_POINT
inline
#endif
void Mat_pleine::Triangulation (int N,Mat_pleine& B)
{ double somme;
// Determination des composantes de la matrice
// triangulaire inferieure B :
for (int k=1;k<=N;k++)
{
somme=0;
for (int j=1;j<=k-1;j++)
{
// somme=somme+ B(k,j)*B(k,j); // pow(B(k,j),2.0);
somme += B(k,j)*B(k,j); // pow(B(k,j),2.0);
// somme += B(j,k)*B(j,k); // pow(B(k,j),2.0); on utilise que la sup
};
B(k,k)=sqrt((*this)(k,k)-somme);
for (int i=k+1;i<=N;i++)
{
somme=0;
for (int j=1;j<=k-1;j++)
{
// somme=somme+(B(i,j)*B(k,j));
somme += (B(i,j)*B(k,j));
// somme += (B(j,i)*B(j,k)); // on utilise que le triangle sup
};
B(i,k)=((*this)(i,k)-somme)/B(k,k);
};
};
};
// résolution du problème triangularisé
#ifndef MISE_AU_POINT
inline
#endif
void Mat_pleine::Resolution
(int N, const Mat_pleine& B,const Vecteur& b,Vecteur& resul) const
{ // Resolution du systeme B.y=b ( avec y=Bt.x )
// Les valeurs de y sont stockees dans le vecteur result
for (int i=1;i<=N;i++)
{
resul(i)=b(i);
// for (int j=1;j<=i-1;j++)
// resul(i)=resul(i)-B(i,j)*resul(j);
// resul(i)=resul(i)/B(i,i);
for (int j=1;j<=i-1;j++)
resul(i) -= B(i,j)*resul(j);
resul(i) /= B(i,i);
};
// Resolution du systeme Bt.x=y
// Les valeurs de x sont stockees dans le vecteur result
// B=B.Transpose();
// for (int i=N;i>=1;i--)
// {
// for (int j=i+1;j<=N;j++)
// result(i)=result(i)-B(j,i)*result(j);
// result(i)=result(i)/B(i,i);
// };
for (int i=N;i>=1;i--)
{
for (int j=i+1;j<=N;j++)
resul(i) -= B(j,i)*resul(j);
// resul(i)=resul(i)-B(j,i)*resul(j);
// resul(i)=resul(i)/B(i,i);
resul(i) /= B(i,i);
};
};
// résolution directe par la méthode de cramer (pour les petits systèmes dim 1,2,3)
#ifndef MISE_AU_POINT
inline
#endif
void Mat_pleine::Cramer( const Vecteur& b,Vecteur& res) const
{ int iL = this->Nb_ligne(); int cC = this->Nb_colonne();
#ifdef MISE_AU_POINT
if (iL != cC)
{ cout << "\nErreur : le nombre de ligne est different du nb de colonne !\n";
cout << "void Mat_pleine::Cramer(.. \n" << endl;
Sortie(1);
}
else if ((iL != 1) && (iL != 2) && (iL!=3))
{ cout << "\nErreur : desole, actuellement seule les cas dim = 1 ou 2 ou 3"
<< " sont implante !\n";
cout << "void Mat_pleine::Cramer(.. \n" << endl;
Sortie(1);
};
#endif
const Mat_pleine & mat = *this;
double det;
switch (cC)
{ case 1:
{if ( Dabs(mat(1,1)) <= ConstMath::trespetit)
{ cout << "\nErreur : le determinant est nulle !\n";
cout << "Mat_pleine::Cramer(.. \n" << endl;
throw ErrResolve_system_lineaire(1);
Sortie(1);
}
else
{ double res11 = 1./mat(1,1); res(1)=res11 * b(1);}
break;
}
case 2:
{det = mat(1,1) * mat(2,2) - mat(1,2) * mat (2,1);
if ( Dabs(det) <= ConstMath::trespetit)
{ cout << "\nErreur : le determinant est nulle !\n";
cout << "Mat_pleine::Cramer(.. \n" << endl;
throw ErrResolve_system_lineaire(1);
Sortie(1);
};
double res11 = mat(2,2) / det; double res12 = mat(1,2) / det;
double res21 = mat(2,1) / det; double res22 = mat(1,1) / det;
double res1 = res11 * b(1) - res12 * b(2);
double res2 = -res21 * b(1) + res22 * b(2);
res(1) = res1; res(2) = res2; // car le vecteur res peut-être égale à b
break;
}
case 3:
{det = mat(1,1) *( mat(2,2) * mat(3,3) - mat(2,3) * mat (3,2));
det -= mat(1,2) *( mat(2,1) * mat(3,3) - mat(2,3) * mat (3,1));
det += mat(1,3) *( mat(2,1) * mat(3,2) - mat(2,2) * mat (3,1));
if ( Dabs(det) <= ConstMath::trespetit)
{ cout << "\nErreur : le determinant est nulle !\n";
cout << "Mat_pleine::Cramer(.. \n" << endl;
throw ErrResolve_system_lineaire(1);
Sortie(1);
};
double res11 = (mat(2,2) * mat(3,3) - mat(2,3) * mat (3,2)) / det;
double res22 = (mat(1,1) * mat(3,3) - mat(1,3) * mat (3,1)) / det;
double res33 = (mat(1,1) * mat(2,2) - mat(1,2) * mat (2,1)) / det;
double res12 = -(mat(1,2) * mat(3,3) + mat(1,3) * mat (3,2)) / det;
double res13 = (mat(1,2) * mat(2,3) - mat(1,3) * mat (2,2)) / det;
double res21 = -(mat(2,1) * mat(3,3) + mat(2,3) * mat (3,1)) / det;
double res31 = (mat(2,1) * mat(3,2) - mat(2,2) * mat (3,1)) / det;
double res32 = -(mat(1,1) * mat(3,2) + mat(1,2) * mat (3,1)) / det;
double res23 = -(mat(1,1) * mat(2,3) - mat(1,3) * mat (2,1)) / det;
double res1 = res11 * b(1) + res12 * b(2) + res13 * b(3);
double res2 = res21 * b(1) + res22 * b(2) + res23 * b(3);
double res3 = res31 * b(1) + res32 * b(2) + res33 * b(3);
res(1) = res1; res(2) = res2; res(3) = res3; // car res peut-être égale à b
break;
}
};
};
// résolution directe par la méthode de cramer (pour les petits systèmes dim 1,2,3)
// idem précédent mais pour un tableau de vecteurs
#ifndef MISE_AU_POINT
inline
#endif
void Mat_pleine::Cramer( const Tableau <Vecteur>& b,Tableau <Vecteur>& res) const
{ int iL = this->Nb_ligne(); int cC = this->Nb_colonne();
#ifdef MISE_AU_POINT
if (iL != cC)
{ cout << "\nErreur : le nombre de ligne est different du nb de colonne !\n";
cout << "void Mat_pleine::Cramer(.. \n" << endl;
Sortie(1);
}
else if ((iL != 1) && (iL != 2) && (iL!=3))
{ cout << "\nErreur : desole, actuellement seule les cas dim = 1 ou 2 ou 3"
<< " sont implante !\n";
cout << "void Mat_pleine::Cramer(.. \n" << endl;
Sortie(1);
};
#endif
const Mat_pleine & mat = *this;
double det;
switch (cC)
{ case 1:
{if ( Dabs(mat(1,1)) <= ConstMath::trespetit)
{ cout << "\nErreur : le determinant est nulle !\n";
cout << "Mat_pleine::Cramer(.. \n" << endl;
throw ErrResolve_system_lineaire(1);
Sortie(1);
};
double res11 = 1./mat(1,1);
int dimtab = b.Taille();
res.Change_taille(dimtab);
for (int i=1;i<= dimtab;i++)
{ res(i).Change_taille(1); res(i)(1)= b(i)(1) * res11;}
break;
}
case 2:
{det = mat(1,1) * mat(2,2) - mat(1,2) * mat (2,1);
if ( Dabs(det) <= ConstMath::trespetit)
{ cout << "\nErreur : le determinant est nulle !\n";
cout << "Mat_pleine::Cramer(.. \n" << endl;
throw ErrResolve_system_lineaire(1);
Sortie(1);
};
double res11 = mat(2,2) / det; double res12 = mat(2,1) / det;
double res21 = mat(1,2) / det; double res22 = mat(1,1) / det;
int dimtab = b.Taille();
res.Change_taille(dimtab);
for (int i=1;i<= dimtab;i++)
{ res(i).Change_taille(2); Vecteur& bi=b(i);
double res1 = res11 * bi(1) + res12 * bi(2);
double res2 = res21 * bi(1) + res22 * bi(2);
res(i)(1) = res1; res(i)(2) = res2; // car le vecteur res peut-être égale à b
};
break;
}
case 3:
{det = mat(1,1) *( mat(2,2) * mat(3,3) - mat(2,3) * mat (3,2));
det -= mat(1,2) *( mat(2,1) * mat(3,3) - mat(2,3) * mat (3,1));
det += mat(1,3) *( mat(2,1) * mat(3,2) - mat(2,2) * mat (3,1));
if ( Dabs(det) <= ConstMath::trespetit)
{ cout << "\nErreur : le determinant est nulle !\n";
cout << "Mat_pleine::Cramer(.. \n" << endl;
throw ErrResolve_system_lineaire(1);
Sortie(1);
};
double res11 = (mat(2,2) * mat(3,3) - mat(2,3) * mat (3,2)) / det;
double res22 = (mat(1,1) * mat(3,3) - mat(1,3) * mat (3,1)) / det;
double res33 = (mat(1,1) * mat(2,2) - mat(1,2) * mat (2,1)) / det;
double res12 = -(mat(1,2) * mat(3,3) + mat(1,3) * mat (3,2)) / det;
double res13 = (mat(1,2) * mat(2,3) - mat(1,3) * mat (2,2)) / det;
double res21 = -(mat(2,1) * mat(3,3) + mat(2,3) * mat (3,1)) / det;
double res31 = (mat(2,1) * mat(3,2) - mat(2,2) * mat (3,1)) / det;
double res32 = -(mat(1,1) * mat(3,2) + mat(1,2) * mat (3,1)) / det;
double res23 = -(mat(1,1) * mat(2,3) - mat(1,3) * mat (2,1)) / det;
int dimtab = b.Taille();
res.Change_taille(dimtab);
for (int i=1;i<= dimtab;i++)
{ res(i).Change_taille(3); Vecteur& bi=b(i); Vecteur& resi=res(i);
double res1 = res11 * bi(1) + res12 * bi(2) + res13 * bi(3);
double res2 = res21 * bi(1) + res22 * bi(2) + res23 * bi(3);
double res3 = res31 * bi(1) + res32 * bi(2) + res33 * bi(3);
resi(1) = res1; resi(2) = res2; resi(3) = res3; // car res peut-être égale à b
};
break;
}
};
};
2023-05-03 17:23:49 +02:00
// méthode interne pour modifier éventuellement les tailles
// en liaison entre les deux stockages
#ifndef MISE_AU_POINT
inline
#endif
void Mat_pleine::Change_tailles_tab_val(int nb_lig,int nb_col)
{ // on considère que val_de_base est déjà ok
#ifdef MISE_AU_POINT
if (nb_lig*nb_col != val_de_base.Taille())
{ cout << "\nErreur : la taille de val_de_base: " << val_de_base.Taille()
<< " est differente de nb_lig*nb_col: " << nb_lig*nb_col;
cout << "void Change_tailles_tab_val(.. \n" << endl;
Sortie(1);
}
#endif
// l'objectif est la modification éventuelle de val
val.Change_taille(nb_lig);
double* pt = val_de_base.Pointeur_vect(); // init
for (int i=1;i<= nb_lig;i++)
{bool creation = false;
if (val(i) != NULL)
{ double* pt_val = val(i)->Pointeur_vect();
if ((pt_val != pt ) // il y a discordance
|| (val(i)->Taille() != nb_col)
)
{delete val(i); // on doit reconstruire le vecteur
creation = true;
};
// sinon c'est ok, on ne fait rien
}
else
{creation = true;};
if (creation) // on crée un vecteur lié
{// il n'y a pas d'initialisation de composantes des vecteur
// elles gardent les valeurs de val_de_base
bool memoire = false;
val(i)=new Vecteur(nb_col,pt,memoire);
};
// on déplace le pointeur de base sur le vecteur suivant
for (int j = 1; j<= nb_col; j++)
pt++;
};
};
// méthode interne pour créer val en fonction des dimensions
// ce qui crée une liaison entre les deux stockages
#ifndef MISE_AU_POINT
inline
#endif
void Mat_pleine::Creation_tab_val(int nb_lig,int nb_col)
{ // on considère que val_de_base est déjà ok
#ifdef MISE_AU_POINT
if (nb_lig*nb_col != val_de_base.Taille())
{ cout << "\nErreur : la taille de val_de_base: " << val_de_base.Taille()
<< " est differente de nb_lig*nb_col: " << nb_lig*nb_col;
cout << "void Mat_pleine::Creation_tab_val(.. \n" << endl;
Sortie(1);
}
#endif
// l'objectif est la construction de val
if (val.Taille() != 0)
// si val est déjà construit on voit s'il faut changer quelque chose
{Change_tailles_tab_val(nb_lig,nb_col);}
else
// sinon on construit
{bool memoire = false;
val.Change_taille(nb_lig);
double* pt = val_de_base.Pointeur_vect(); // init
for (int i=1;i<= nb_lig;i++)
{// il n'y a pas d'initialisation de composantes des vecteur
// elles gardent les valeurs de val_de_base
val(i)=new Vecteur(nb_col,pt,memoire);
for (int j = 1; j<= nb_col; j++)
pt++;
};
};
};
2021-09-28 17:28:23 +02:00
// surcharge de l'operateur de lecture typée
#ifndef MISE_AU_POINT
inline
#endif
istream & operator >> (istream & entree, Mat_pleine & mat)
{ // vérification du type
string type;
entree >> type;
if (type != "Mat_pleine")
2023-05-03 17:23:49 +02:00
{cout << "\n*** erreur en lecture d'une matrice, on attendait Mat_pleine "
<< " et on a lue: "<< type;
Sortie (1);
2021-09-28 17:28:23 +02:00
return entree;
};
2023-05-03 17:23:49 +02:00
// les tailles
int nb_ligne, nb_colonne;
string toto;
entree >> nb_ligne >> toto >> nb_colonne;
// redimentionnement éventuelle de la matrice
mat.Initialise (nb_ligne,nb_colonne,0.0);
// on lit les composantes
entree >> mat.val_de_base ;
2021-09-28 17:28:23 +02:00
return entree;
};
// surcharge de l'operateur d'ecriture typée
#ifndef MISE_AU_POINT
inline
#endif
ostream & operator << ( ostream & sort,const Mat_pleine & mat)
{ // un indicateur donnant le type puis les datas
2023-05-03 17:23:49 +02:00
sort << "Mat_pleine " << mat.Nb_ligne() << " X "<< mat.Nb_colonne() << " "
<< mat.val_de_base ;
2021-09-28 17:28:23 +02:00
return sort;
};
#endif