1635 lines
60 KiB
C++
Executable file
1635 lines
60 KiB
C++
Executable file
|
|
|
|
// 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.
|
|
//
|
|
// Copyright (C) 1997-2021 Université Bretagne Sud (France)
|
|
// 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 "TenseurQ-2.h"
|
|
#include "ConstMath.h"
|
|
#include "MathUtil.h"
|
|
#include "Tenseur2.h"
|
|
#include "CharUtil.h"
|
|
|
|
#ifndef TenseurQ2_H_deja_inclus
|
|
|
|
// 1) pour le passage repère absolu en repère local
|
|
// T1 et T2 sont des tenseurs, Base représente une base
|
|
// pour des tenseurs T1 et T2 symétriques HHHH ou BBBB
|
|
template <class T1,class T2,class Base>
|
|
T1& produit44_2(T1& A,T2& Nous, const Base & gi)
|
|
{ int dim = 2;
|
|
if ((A.Dimension()!=22)|(Nous.Dimension()!=22))
|
|
{cout << "\n erreur la methode n'est valide que pour la dimension 2"
|
|
<< "\n ici pour des tenseurs Tenseur2HHHH ou Tenseur2BBBB"
|
|
<< "\n T1& produit44_2(T1& A,T2& Nous, const Base & gi)";
|
|
Sortie(1);
|
|
};
|
|
int dimint = gi.Dimension(); //peut être de dim 2 ou 3
|
|
for (int i=1;i<= dim; i++)
|
|
for (int j=1; j<= i; j++)
|
|
for (int k=1;k<= dim; k++)
|
|
for (int l=1; l<= k; l++)
|
|
{int ijkl = (Tenseur2HHHH::cdex2HHHH.odVect(i,j)-1)*3+Tenseur2HHHH::cdex2HHHH.odVect(k,l)-1;
|
|
A.t[ijkl]=0.;
|
|
for (int al=1; al<= dimint; al++)
|
|
for (int be=1; be<= dimint; be++)
|
|
for (int cl=1; cl<= dimint; cl++)
|
|
for (int de=1; de<= dimint; de++)
|
|
{ int abcd=(Tenseur2HHHH::cdex2HHHH.odVect(al,be)-1)*3+Tenseur2HHHH::cdex2HHHH.odVect(cl,de)-1;
|
|
A.t[ijkl] += Nous.t[abcd] * gi.Coordo(i)(al) * gi.Coordo(j)(be)
|
|
* gi.Coordo(k)(cl) * gi.Coordo(l)(de);
|
|
}
|
|
}
|
|
return A;
|
|
};
|
|
|
|
// 2) pour le changement de base
|
|
// T1 et T2 sont des tenseurs, Base représente une base
|
|
// pour des tenseurs T1 et T2 symétriques HHHH ou BBBB
|
|
// A = A?{ijkl) g?i rond g?j rond g?k rond g?l = A'?{efgh) gp?i rond gp?j rond gp?k rond gp?l
|
|
// g?i = beta?i?j gp?j --> A'?{efgh) = A?{ijkl) beta?i?e beta?j?f beta?k?g beta?l?h
|
|
// ici seul g?1 et g?2 est utilisé et les coordonnées: g?al?i sont celles dans gp?j
|
|
// c'est à dire: coordonnées des anciens vecteurs dans la nouvelle base
|
|
|
|
// --- par exemple en contravariant ---
|
|
// A = A^{abcd) g_a rond g_b rond g_c rond g_d = A'^{ijkl) gp_i rond gp_j rond gp_k rond gp_l
|
|
// g_a = beta_a^i gp_i --> A'^{ijkl) = A^{abcd) beta_a^i beta_b^j beta_c^k beta_d^l
|
|
// ici seul g_1 et g_2 sont utilisés et les coordonnées: g_al^i sont celles dans gp_i
|
|
// c'est à dire: coordonnées des anciens vecteurs dans la nouvelle base
|
|
|
|
|
|
template <class T1,class T2,class Base>
|
|
T1& produit55_2(T1& A,T2& Nous, const Base & gi)
|
|
{ //cout << "\n debug: T1& produit55_2(T1& A,T2& Nous, const Base & gi) "
|
|
// << " debut1 " << flush;
|
|
int dim = 2;
|
|
//cout << "\n debug: T1& produit55_2(T1& A,T2& Nous, const Base & gi) "
|
|
// << " debut2 " << flush;
|
|
int dim_gi = gi(1).Dimension(); // dim de l'espace d'arrivé
|
|
//cout << "\n debug: T1& produit55_2(T1& A,T2& Nous, const Base & gi) "
|
|
// << " debut3 " << flush;
|
|
#ifdef MISE_AU_POINT
|
|
// on vérifie que la dimention de retour convient
|
|
// on doit avoir : gal = beta(al,j) g'j
|
|
bool erreur_dim = false;
|
|
switch (A.Dimension())
|
|
{case 22: if (dim_gi != 2) erreur_dim = true; break;
|
|
case 33: if (dim_gi != 3) erreur_dim = true; break;
|
|
default:
|
|
cout << "\n erreur*** cas non pris en compte: A.Dimension()= "
|
|
<< A.Dimension()
|
|
<< "\n T1& produit55_2(T1& A,T2& Nous, const Base & gi) ";
|
|
Sortie(1);
|
|
}
|
|
if (erreur_dim)
|
|
{cout << "\n erreur*** de dimension: A.Dimension()= "
|
|
<< A.Dimension() << ", dim_gi= "<< dim_gi
|
|
<< "\n T1& produit55_2(T1& A,T2& Nous, const Base & gi) ";
|
|
Sortie(1);
|
|
};
|
|
#endif
|
|
// if ((A.Dimension()!=22)|(Nous.Dimension()!=22))
|
|
// {cout << "\n erreur la methode n'est valide que pour la dimension 2"
|
|
// << "\n ici pour des tenseurs Tenseur2HHHH ou Tenseur2BBBB"
|
|
// << "\n T1& produit55_2(T1& A,T2& Nous, const Base & gi)";
|
|
// Sortie(1);
|
|
// };
|
|
|
|
// if (ParaGlob::NiveauImpression() > 8)
|
|
// {cout << "\n debug: T1& produit55_2(T1& A,T2& Nous, const Base & gi) "
|
|
// << " debut " << flush;
|
|
// };
|
|
|
|
switch (A.Dimension())
|
|
{case 22:
|
|
// if (ParaGlob::NiveauImpression() > 8)
|
|
// {cout << "\n debug: T1& produit55_2(T1& A,T2& Nous, const Base & gi) "
|
|
// << " cas 22, dim_gi= " << dim_gi << flush;
|
|
// };
|
|
{int dimint = 2;
|
|
for (int i=1;i<= dim_gi; i++)
|
|
for (int j=1; j<= i; j++)
|
|
for (int k=1;k<= dim_gi; k++)
|
|
for (int l=1; l<= k; l++)
|
|
{int ijkl =
|
|
(Tenseur2HHHH::cdex2HHHH.odVect(i,j)-1)*3+Tenseur2HHHH::cdex2HHHH.odVect(k,l)-1;
|
|
A.t[ijkl]=0.;
|
|
for (int al=1; al<= dimint; al++)
|
|
for (int be=1; be<= dimint; be++)
|
|
for (int cl=1; cl<= dimint; cl++)
|
|
for (int de=1; de<= dimint; de++)
|
|
{ int abcd=(Tenseur2HHHH::cdex2HHHH.odVect(al,be)-1)*3
|
|
+Tenseur2HHHH::cdex2HHHH.odVect(cl,de)-1;
|
|
A.t[ijkl] += Nous.t[abcd] * gi.Coordo(al)(i) * gi.Coordo(be)(j)
|
|
* gi.Coordo(cl)(k) * gi.Coordo(de)(l);
|
|
}
|
|
}
|
|
break;
|
|
}
|
|
case 33:
|
|
{int dimint = 2;
|
|
for (int i=1;i<= dim_gi; i++)
|
|
for (int j=1; j<= i; j++)
|
|
for (int k=1;k<= dim_gi; k++)
|
|
for (int l=1; l<= k; l++)
|
|
{ A.Change(i,j,k,l,0.);
|
|
for (int al=1; al<= dimint; al++)
|
|
for (int be=1; be<= dimint; be++)
|
|
for (int cl=1; cl<= dimint; cl++)
|
|
for (int de=1; de<= dimint; de++)
|
|
{ int abcd=(Tenseur2HHHH::cdex2HHHH.odVect(al,be)-1)*3
|
|
+Tenseur2HHHH::cdex2HHHH.odVect(cl,de)-1;
|
|
double val = Nous.t[abcd] * gi.Coordo(al)(i) * gi.Coordo(be)(j)
|
|
* gi.Coordo(cl)(k) * gi.Coordo(de)(l);
|
|
|
|
A.ChangePlus(i,j,k,l,val);
|
|
}
|
|
}
|
|
break;
|
|
}
|
|
default:
|
|
cout << "\n erreur*** cas non pris en compte: A.Dimension()= "
|
|
<< A.Dimension()
|
|
<< "\n T1& produit55_2(T1& A,T2& Nous, const Base & gi) ";
|
|
Sortie(1);
|
|
}
|
|
|
|
|
|
return A;
|
|
};
|
|
|
|
// variables globales
|
|
// initialisation dans EnteteTenseur.h , utilisé dans le progr principal
|
|
|
|
//------------------------------------------------------------------
|
|
// cas des composantes 4 fois contravariantes 2HHHH
|
|
//------------------------------------------------------------------
|
|
// --- gestion de changement d'index ----
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
Tenseur2HHHH::ChangementIndex::ChangementIndex() :
|
|
idx_i(3),idx_j(3),odVect(2)
|
|
{ idx_i(1)=1;idx_i(2)=2; idx_j(1)=1;idx_j(2)=2;
|
|
idx_i(3)=1; idx_j(3)=2;
|
|
odVect(1,1)=1;odVect(1,2)=3;
|
|
odVect(2,1)=3;odVect(2,2)=2;
|
|
};
|
|
|
|
// Constructeur
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
Tenseur2HHHH::Tenseur2HHHH() :
|
|
ipointe() // par défaut
|
|
{ dimension = 22;
|
|
listdouble9.push_front(Reels9()); // allocation
|
|
ipointe = listdouble9.begin(); // recup de la position de la maille dans la liste
|
|
t = (ipointe)->donnees; // recup de la position des datas dans la maille
|
|
for (int i=0;i<9;i++) t[i]=0.;
|
|
};
|
|
// initialisation de toutes les composantes a une meme valeur
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
Tenseur2HHHH::Tenseur2HHHH( const double val) :
|
|
ipointe()
|
|
{ dimension = 22;
|
|
listdouble9.push_front(Reels9()); // allocation
|
|
ipointe = listdouble9.begin(); // recup de la position de la maille dans la liste
|
|
t = (ipointe)->donnees; // recup de la position des datas dans la maille
|
|
for (int i=0;i<9;i++) t[i]=val;
|
|
};
|
|
// initialisation à partir d'un produit tensoriel avec 2 cas
|
|
// booleen = true : produit tensoriel normal
|
|
// *this=aHH(i,j).bHH(k,l) gBi gBj gBk gBl
|
|
// booleen = false : produit tensoriel barre
|
|
// *this(i,j,k,l)
|
|
// = 1/4 * (a(i,k).b(j,l) + a(j,k).b(i,l) + a(i,l).b(j,k) + a(j,l).b(i,k))
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
Tenseur2HHHH::Tenseur2HHHH(bool normal, const TenseurHH & aHH, const TenseurHH & bHH) :
|
|
ipointe()
|
|
{ dimension = 22;
|
|
listdouble9.push_front(Reels9()); // allocation
|
|
ipointe = listdouble9.begin(); // recup de la position de la maille dans la liste
|
|
t = (ipointe)->donnees; // recup de la position des datas dans la maille
|
|
const Tenseur2HH & a2HH = *((Tenseur2HH*) &aHH); // passage en dim 2
|
|
const Tenseur2HH & b2HH = *((Tenseur2HH*) &bHH); // passage en dim 2
|
|
#ifdef MISE_AU_POINT
|
|
if (Dabs(a2HH.Dimension()) != 2)
|
|
Message(2,string("produit tensoriel a partir d'un premier tenseur non symétriques \n")
|
|
+"Tenseur2HHHH::Tenseur2HHHH(bool normal, const"
|
|
+" TenseurHH & aHH, const TenseurHH & bHH);");
|
|
if (Dabs(b2HH.Dimension()) != 2)
|
|
Message(2,string("produit tensoriel a partir d'un second tenseur non symétriques \n")
|
|
+"Tenseur2HHHH::Tenseur2HHHH(bool normal, const"
|
|
+" TenseurHH & aHH, const TenseurHH & bHH);");
|
|
#endif
|
|
if (normal)
|
|
{ for (int ij=1;ij < 4;ij++)
|
|
for (int kl=1;kl < 4;kl++)
|
|
t[(ij-1)*3+kl-1] = a2HH(cdex2HHHH.idx_i(ij),cdex2HHHH.idx_j(ij))
|
|
* b2HH(cdex2HHHH.idx_i(kl),cdex2HHHH.idx_j(kl));
|
|
}
|
|
else
|
|
{ for (int ij=1;ij < 4;ij++)
|
|
for (int kl=1;kl < 4;kl++)
|
|
{ int i = cdex2HHHH.idx_i(ij); int j = cdex2HHHH.idx_j(ij);
|
|
int k = cdex2HHHH.idx_i(kl); int l = cdex2HHHH.idx_j(kl);
|
|
t[(ij-1)*3+kl-1] = 0.25 * (
|
|
a2HH(i,k) * b2HH(j,l) + a2HH(j,k) * b2HH(i,l)
|
|
+a2HH(i,l) * b2HH(j,k) + a2HH(j,l) * b2HH(i,k)
|
|
);
|
|
};
|
|
}
|
|
};
|
|
// DESTRUCTEUR :
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
Tenseur2HHHH::~Tenseur2HHHH()
|
|
{ listdouble9.erase(ipointe);} ; // suppression de l'élément de la liste
|
|
|
|
// constructeur a partir d'une instance non differenciee
|
|
// dans le cas d'un tenseur quelconque, celui-ci
|
|
// est converti à condition que les symétries existent sinon erreur en debug
|
|
// opération longue dans ce cas !
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
Tenseur2HHHH::Tenseur2HHHH ( const TenseurHHHH & B) :
|
|
ipointe()
|
|
{ dimension = 22;
|
|
// #ifdef MISE_AU_POINT
|
|
// if (Dabs(dimension) != 22)
|
|
// { cout << "\n erreur de dimension, elle devrait etre = 22 ";
|
|
// cout << "\n Tenseur2HHHH::Tenseur2HHHH ( TenseurHHHH &) " << endl;
|
|
// Sortie(1);
|
|
// }
|
|
// #endif
|
|
listdouble9.push_front(Reels9()); // allocation
|
|
ipointe = listdouble9.begin(); // recup de la position de la maille dans la liste
|
|
t = (ipointe)->donnees; // recup de la position des datas dans la maille
|
|
if (Dabs(B.dimension) == 22 ) // cas d'un tenseur du même type
|
|
{ for (int i=0;i< 9;i++)
|
|
t[i] = B.t[i];
|
|
}
|
|
else
|
|
{// cas d'un tenseur quelconque
|
|
double Z=B.MaxiComposante();
|
|
for (int i=1;i < 3;i++)
|
|
for (int j=1;j<=i;j++)
|
|
for (int k=1;k < 3;k++)
|
|
for (int l=1;l<=k;l++)
|
|
{// on teste les symétries et on affecte
|
|
double a = B(i,j,k,l);
|
|
#ifdef MISE_AU_POINT
|
|
if ((!diffpourcent(a,B(j,i,k,l),Z,ConstMath::unpeupetit)
|
|
&& !diffpourcent(a,B(i,j,l,k),Z,ConstMath::unpeupetit))
|
|
|| (Abs(Z) < ConstMath::trespetit) )
|
|
// erreur d'affectation
|
|
if (ParaGlob::NiveauImpression() > 5)
|
|
cout << "\n tenseurHHHH (ijkl= " << i << "," << j << "," << k << "," << l << ")= "
|
|
<< a << " " << B(j,i,k,l) << " " <<B(i,j,l,k) ;
|
|
if (ParaGlob::NiveauImpression() > 1)
|
|
cout << "WARNING ** erreur constructeur, tenseur non symetrique, Tenseur2HHHH::Tenseur2HHHH(const TenseurHHHH & B)";
|
|
#endif
|
|
// si il y a un pb normalement il y a eu un message
|
|
this->Change(i,j,k,l,a);
|
|
};
|
|
};
|
|
};
|
|
// constructeur de copie
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
Tenseur2HHHH::Tenseur2HHHH ( const Tenseur2HHHH & B) :
|
|
ipointe()
|
|
{ this->dimension = B.dimension;
|
|
listdouble9.push_front(Reels9()); // allocation
|
|
ipointe = listdouble9.begin(); // recup de la position de la maille dans la liste
|
|
t = (ipointe)->donnees; // recup de la position des datas dans la maille
|
|
for (int i=0;i< 9;i++)
|
|
this->t[i] = B.t[i];
|
|
};
|
|
// METHODES PUBLIQUES :
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
// initialise toutes les composantes à val
|
|
void Tenseur2HHHH::Inita(double val)
|
|
{ for (int i=0;i< 9;i++)
|
|
t[i] = val;
|
|
};
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
TenseurHHHH & Tenseur2HHHH::operator + ( const TenseurHHHH & B) const
|
|
{ TenseurHHHH * res;
|
|
#ifdef MISE_AU_POINT
|
|
if (B.Dimension() != 22) Message(2,"Tenseur2HHHH::operator + ( etc..");
|
|
#endif
|
|
res = new Tenseur2HHHH;
|
|
LesMaillonsHHHH::NouveauMaillon( res); // ajout d'un tenseur intermediaire
|
|
for (int i = 0; i< 9; i++)
|
|
res->t[i] = this->t[i] + B.t[i]; //somme des données
|
|
return *res ;
|
|
};
|
|
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
void Tenseur2HHHH::operator += ( const TenseurHHHH & B)
|
|
{
|
|
#ifdef MISE_AU_POINT
|
|
if (Dabs(B.Dimension()) != 22) Message(2,"Tenseur2HHHH::operator += ( etc..");
|
|
#endif
|
|
for (int i = 0; i< 9; i++)
|
|
this->t[i] += B.t[i];
|
|
LesMaillonsHHHH::Libere(); // destruction des tenseurs intermediaires
|
|
}; //somme des données
|
|
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
TenseurHHHH & Tenseur2HHHH::operator - () const
|
|
{ TenseurHHHH * res;
|
|
res = new Tenseur2HHHH;
|
|
LesMaillonsHHHH::NouveauMaillon( res); // ajout d'un tenseur intermediaire
|
|
for (int i = 0; i< 9; i++)
|
|
res->t[i] = - this->t[i]; //oppose
|
|
return *res ;
|
|
};
|
|
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
TenseurHHHH & Tenseur2HHHH::operator - ( const TenseurHHHH & B) const
|
|
{ TenseurHHHH * res;
|
|
#ifdef MISE_AU_POINT
|
|
if (Dabs(B.Dimension()) != 22) Message(2,"Tenseur2HHHH::operator - ( etc..");
|
|
#endif
|
|
res = new Tenseur2HHHH;
|
|
LesMaillonsHHHH::NouveauMaillon( res); // ajout d'un tenseur intermediaire
|
|
for (int i = 0; i< 9; i++)
|
|
res->t[i] = this->t[i] - B.t[i]; //soustraction des données
|
|
return *res ;
|
|
};
|
|
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
void Tenseur2HHHH::operator -= ( const TenseurHHHH & B)
|
|
{
|
|
#ifdef MISE_AU_POINT
|
|
if (Dabs(B.Dimension()) != 22) Message(2,"Tenseur2HHHH::operator -= ( etc..");
|
|
#endif
|
|
for (int i = 0; i< 9; i++)
|
|
this->t[i] -= B.t[i];
|
|
LesMaillonsHHHH::Libere(); // destruction des tenseurs intermediaires
|
|
}; //soustraction des données
|
|
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
TenseurHHHH & Tenseur2HHHH::operator = ( const TenseurHHHH & B)
|
|
{
|
|
#ifdef MISE_AU_POINT
|
|
if (Dabs(B.Dimension()) != 22) Message(2,"Tenseur2HHHH::operator = ( etc..");
|
|
#endif
|
|
for (int i = 0; i< 9; i++)
|
|
this->t[i] = B.t[i];
|
|
LesMaillonsHHHH::Libere(); // destruction des tenseurs intermediaires
|
|
return *this;
|
|
}; //affectation des données;
|
|
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
TenseurHHHH & Tenseur2HHHH::operator * ( const double & b) const
|
|
{ TenseurHHHH * res;
|
|
res = new Tenseur2HHHH;
|
|
LesMaillonsHHHH::NouveauMaillon( res); // ajout d'un tenseur intermediaire
|
|
for (int i = 0; i< 9; i++)
|
|
res->t[i] = this->t[i] * b; //multiplication des données
|
|
return *res ;
|
|
};
|
|
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
void Tenseur2HHHH::operator *= ( const double & b)
|
|
{for (int i = 0; i< 9; i++)
|
|
this->t[i] *= b ;}; //multiplication des données
|
|
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
TenseurHHHH & Tenseur2HHHH::operator / ( const double & b) const
|
|
{ TenseurHHHH * res;
|
|
res = new Tenseur2HHHH;
|
|
LesMaillonsHHHH::NouveauMaillon( res); // ajout d'un tenseur intermediaire
|
|
#ifdef MISE_AU_POINT
|
|
if (Dabs(b) < ConstMath::trespetit)
|
|
{ cout << "\n erreur le diviseur est trop petit = " << b;
|
|
cout << "\n Tenseur2HHHH::operator / ( const double & b) " << endl;
|
|
Sortie(1);
|
|
}
|
|
#endif
|
|
for (int i = 0; i< 9; i++)
|
|
res->t[i] = this->t[i] / b; //division des données
|
|
return *res ;
|
|
};
|
|
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
void Tenseur2HHHH::operator /= ( const double & b)
|
|
{
|
|
#ifdef MISE_AU_POINT
|
|
if (Dabs(b) < ConstMath::trespetit)
|
|
{ cout << "\n erreur le diviseur est trop petit = " << b;
|
|
cout << "\n Tenseur2HHHH::operator /= ( const double & b) " << endl;
|
|
Sortie(1);
|
|
}
|
|
#endif
|
|
for (int i = 0; i< 9; i++)
|
|
this->t[i] /= b ;}; //division des données
|
|
|
|
|
|
// produit contracte à droite avec un tenseur du second ordre
|
|
// différent à gauche !!
|
|
// différent à gauche !!
|
|
// A(i,j,k,l)*B(l,k)=A..B
|
|
// on commence par contracter l'indice du milieu puis externe
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
TenseurHH& Tenseur2HHHH::operator && ( const TenseurBB & aBB) const
|
|
{
|
|
#ifdef MISE_AU_POINT
|
|
if (Dabs(aBB.Dimension()) != 2)
|
|
Message(2,"Tenseur2HHHH::operator && ( const TenseurBB & aBB)");
|
|
#endif
|
|
TenseurHH * res;
|
|
res = new Tenseur2HH;
|
|
LesMaillonsHH::NouveauMaillon( res); // ajout d'un tenseur intermediaire
|
|
const Tenseur2BB & a2BB = *((Tenseur2BB*) &aBB); // passage en dim 2
|
|
for (int ij=1;ij < 4;ij++)
|
|
{for (int kl=1;kl < 3;kl++) // partie simple produit : la partie diagonale
|
|
(*res).Coor(cdex2HHHH.idx_i(ij),cdex2HHHH.idx_j(ij)) += t[(ij-1)*3+kl-1]
|
|
* a2BB(cdex2HHHH.idx_i(kl),cdex2HHHH.idx_j(kl));
|
|
int kl = 3; // partie produit doublée : partie extra-diagonale
|
|
// comme le tenseur du second ordre n'est pas forcément symétrique on utilise les
|
|
// 2 coordonnées de chaque coté de la diagonale,
|
|
(*res).Coor(cdex2HHHH.idx_i(ij),cdex2HHHH.idx_j(ij)) += t[(ij-1)*3+kl-1]
|
|
* (a2BB(cdex2HHHH.idx_i(kl),cdex2HHHH.idx_j(kl))
|
|
+a2BB(cdex2HHHH.idx_j(kl),cdex2HHHH.idx_i(kl)));
|
|
};
|
|
return *res ;
|
|
};
|
|
//fonctions définissant le produit tensoriel normal de deux tenseurs
|
|
// *this=aHH(i,j).bHH(k,l) gBi gBj gBk gBl
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
TenseurHHHH & Tenseur2HHHH::Prod_tensoriel(const TenseurHH & aHH, const TenseurHH & bHH)
|
|
{ TenseurHHHH * res;
|
|
res = new Tenseur2HHHH;
|
|
LesMaillonsHHHH::NouveauMaillon( res); // ajout d'un tenseur intermediaire
|
|
const Tenseur2HH & a2HH = *((Tenseur2HH*) &aHH); // passage en dim 3
|
|
const Tenseur2HH & b2HH = *((Tenseur2HH*) &bHH); // passage en dim 3
|
|
#ifdef MISE_AU_POINT
|
|
if (Dabs(a2HH.Dimension()) != 2)
|
|
{ cout << "\n produit tensoriel a partir d'un premier tenseur non symétriques \n"
|
|
<< "Tenseur2HHHH::Prod_tensoriel(const TenseurHH & aHH, const TenseurHH & bHH)";
|
|
Sortie(2);
|
|
}
|
|
if (Dabs(b2HH.Dimension()) != 2)
|
|
{ cout << "\n produit tensoriel a partir d'un second tenseur non symétriques \n"
|
|
<< "Tenseur2HHHH::Prod_tensoriel(const TenseurHH & aHH, const TenseurHH & bHH)";
|
|
Sortie(2);
|
|
}
|
|
#endif
|
|
for (int ij=1;ij < 4;ij++)
|
|
for (int kl=1;kl < 4;kl++)
|
|
res->t[(ij-1)*3+kl-1] = a2HH(cdex2HHHH.idx_i(ij),cdex2HHHH.idx_j(ij))
|
|
* b2HH(cdex2HHHH.idx_i(kl),cdex2HHHH.idx_j(kl));
|
|
return *res;
|
|
};
|
|
//fonctions définissant le produit tensoriel barre de deux tenseurs
|
|
// concervant les symétries !!
|
|
// *this(i,j,k,l)
|
|
// = 1/4 * (a(i,k).b(j,l) + a(j,k).b(i,l) + a(i,l).b(j,k) + a(j,l).b(i,k))
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
TenseurHHHH & Tenseur2HHHH::Prod_tensoriel_barre(const TenseurHH & aHH, const TenseurHH & bHH)
|
|
{ TenseurHHHH * res;
|
|
res = new Tenseur2HHHH;
|
|
LesMaillonsHHHH::NouveauMaillon( res); // ajout d'un tenseur intermediaire
|
|
const Tenseur2HH & a2HH = *((Tenseur2HH*) &aHH); // passage en dim 2
|
|
const Tenseur2HH & b2HH = *((Tenseur2HH*) &bHH); // passage en dim 2
|
|
#ifdef MISE_AU_POINT
|
|
if (Dabs(a2HH.Dimension()) != 2)
|
|
{ cout << "\n produit tensoriel a partir d'un premier tenseur non symétriques \n"
|
|
<< "Tenseur2HHHH::Prod_tensoriel(const TenseurHH & aHH, const TenseurHH & bHH)";
|
|
Sortie(2);
|
|
}
|
|
if (Dabs(b2HH.Dimension()) != 2)
|
|
{ cout << "\n produit tensoriel a partir d'un second tenseur non symétriques \n"
|
|
<< "Tenseur2HHHH::Prod_tensoriel(const TenseurHH & aHH, const TenseurHH & bHH)";
|
|
Sortie(2);
|
|
}
|
|
#endif
|
|
for (int ij=1;ij < 4;ij++)
|
|
for (int kl=1;kl < 4;kl++)
|
|
{ int i = cdex2HHHH.idx_i(ij); int j = cdex2HHHH.idx_j(ij);
|
|
int k = cdex2HHHH.idx_i(kl); int l = cdex2HHHH.idx_j(kl);
|
|
res->t[(ij-1)*3+kl-1] = 0.25 * (
|
|
a2HH(i,k) * b2HH(j,l) + a2HH(j,k) * b2HH(i,l)
|
|
+a2HH(i,l) * b2HH(j,k) + a2HH(j,l) * b2HH(i,k)
|
|
);
|
|
};
|
|
return *res;
|
|
};
|
|
|
|
// ATTENTION creation d'un tenseur transpose qui est supprime par Libere
|
|
// les 2 premiers indices sont échangés avec les deux derniers indices
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
TenseurHHHH & Tenseur2HHHH::Transpose1et2avec3et4() const
|
|
{ TenseurHHHH * res;
|
|
res = new Tenseur2HHHH;
|
|
LesMaillonsHHHH::NouveauMaillon(res); // ajout d'un tenseur intermediaire
|
|
for (int ij=1;ij < 4;ij++)
|
|
for (int kl=1;kl < 4;kl++)
|
|
res->t[(kl-1)*3+ij-1] = t[(ij-1)*3+kl-1] ;
|
|
return *res;
|
|
};
|
|
|
|
|
|
// il s'agit ici de calculer la variation d'un tenseur dans une nouvelle base
|
|
// --> variation par rapport aux composantes covariantes d'un tenseur (ex: composantes eps_ij)
|
|
|
|
// d sigma^ij / d eps_kl = d gamma^i_{.a} / d eps_kl . sigma^ab . gamma^j_{.b}
|
|
// + gamma^i_{.a} . d sigma^ab / d eps_kl. gamma^j_{.b}
|
|
// + gamma^i^_{.a} . sigma^ab . d gamma^j_{.b} / d eps_kl
|
|
|
|
// connaissant sa variation dans la base actuelle
|
|
// var_tensHHHH : en entrée: la variation du tenseur dans la base initiale qu'on appelle g_i
|
|
// ex: var_tensHHHH(ijkl) = d A^ij / d eps_kl
|
|
// : en sortie: la variation du tenseur dans la base finale qu'on appelle gp_i
|
|
// gamma : en entrée gpH(i) = gamma(i,j) * gH(j)
|
|
// var_gamma : en entrée : la variation de gamma
|
|
// ex: var_gamma(i,j,k,l) = d gamma^i_{.j} / d eps_kl
|
|
// tensHH : le tenseur dont on cherche la variation
|
|
/// -- pour mémoire ---
|
|
// changement de base (cf. théorie) : la matrice beta est telle que:
|
|
// gpB(i) = beta(i,j) * gB(j) <==> gp_i = beta_i^j * g_j
|
|
// et la matrice gamma telle que:
|
|
// gamma(i,j) represente les coordonnees de la nouvelle base duale gpH dans l'ancienne gH
|
|
// gpH(i) = gamma(i,j) * gH(j), i indice de ligne, j indice de colonne
|
|
// c-a-d= gp^i = gamma^i_j * g^j
|
|
// rappel des différentes relations entre beta et gamma
|
|
// [beta]^{-1} = [gamma]^T ; [beta]^{-1T} = [gamma]
|
|
// [beta] = [gamma]^{-1T} ; [beta]^{T} = [gamma]^{-1}
|
|
// changement de base pour un tenseur en deux fois covariants:
|
|
// [Ap^kl] = [gamma] * [A^ij] * [gamma]^T
|
|
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
TenseurHHHH & Tenseur2HHHH::Var_tenseur_dans_nouvelle_base
|
|
(const Mat_pleine& gamma,Tenseur2HHHH& var_tensHHHH, const Tableau2 <Tenseur2HH>& var_gamma
|
|
,const Tenseur2HH& tensHH)
|
|
{
|
|
TenseurHHHH * res;
|
|
res = new Tenseur2HHHH;
|
|
LesMaillonsHHHH::NouveauMaillon( res); // ajout d'un tenseur intermediaire
|
|
// d sigma^ij / d eps_kl = d gamma^i_{.a} / d eps_kl . sigma^ab . gamma^j_{.b}
|
|
// + gamma^i_{.a} . d sigma^ab / d eps_kl. gamma^j_{.b}
|
|
// + gamma^i^_{.a} . sigma^ab . d gamma^j_{.b} / d eps_kl
|
|
for (int i=1;i<3;i++)
|
|
for (int j=1;j<3;j++)
|
|
for (int k=1;k<3;k++)
|
|
for (int l=1;l<3;l++)
|
|
{ double d_sig_ij_d_eps_kl = 0.;
|
|
for (int a=1;a<3;a++)
|
|
for (int b=1;b<3;b++)
|
|
d_sig_ij_d_eps_kl += var_gamma(i,a)(k,l) * tensHH(a,b) * gamma(j,b)
|
|
+ gamma(i,a) * var_tensHHHH(a,b,k,l) * gamma(j,b)
|
|
+ gamma(i,a) * tensHH(a,b) * var_gamma(j,b)(k,l);
|
|
res->Change(i,j,k,l,d_sig_ij_d_eps_kl);
|
|
};
|
|
|
|
return *res;
|
|
};
|
|
|
|
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
// affectation de B dans this, plusZero = false: les données manquantes sont inchangées,
|
|
// plusZero = true: les données manquantes sont mises à 0
|
|
// si au contraire la dimension de B est plus grande que *this, il y a uniquement affectation
|
|
// des données possibles
|
|
void Tenseur2HHHH::Affectation_trans_dimension(const TenseurHHHH & aHHHH,bool plusZero)
|
|
{ switch (abs(aHHHH.Dimension()))
|
|
{ case 33 : case 22 :
|
|
for (int ij=1;ij < 4;ij++)
|
|
for (int kl=1;kl < 4;kl++)
|
|
t[(ij-1)*3+kl-1] = aHHHH(cdex2HHHH.idx_i(ij),cdex2HHHH.idx_j(ij)
|
|
,cdex2HHHH.idx_i(kl),cdex2HHHH.idx_j(kl));
|
|
break;
|
|
case 11 :
|
|
if (plusZero)
|
|
this->Inita(0.); // on commence par mettre à 0 si besoin
|
|
// ensuite on affecte
|
|
t[0] = aHHHH(1,1,1,1);
|
|
break;
|
|
default:
|
|
Message(2,string(" *** erreur, la dimension: ")
|
|
+ ChangeEntierSTring(abs(aHHHH.Dimension()))
|
|
+"n'est pas prise en compte \n Tenseur2HHHH::Affectation_trans_dimension(");
|
|
};
|
|
};
|
|
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
// transférer un tenseur général de même dimension accessible en indice, dans un tenseur 9 Tenseur2HHHH
|
|
// il n'y a pas de symétrisation !, seule certaines composantes sont prises en compte
|
|
void Tenseur2HHHH::TransfertDunTenseurGeneral(const TenseurHHHH & aHHHH)
|
|
{ switch (aHHHH.Dimension())
|
|
{ case 11 :
|
|
{int ij=1;
|
|
int kl=1;
|
|
t[(ij-1)*3+kl-1] = aHHHH(1,1,1,1);
|
|
break;
|
|
}
|
|
case 33 : case 22 :
|
|
{for (int ij=1;ij < 4;ij++)
|
|
for (int kl=1;kl < 4;kl++)
|
|
t[(ij-1)*3+kl-1] = aHHHH(cdex2HHHH.idx_i(ij),cdex2HHHH.idx_j(ij)
|
|
,cdex2HHHH.idx_i(kl),cdex2HHHH.idx_j(kl));
|
|
break;
|
|
}
|
|
default:
|
|
Message(2,"\n erreur *** aHHHH.Dimension()= " + ChangeEntierSTring(aHHHH.Dimension())
|
|
+ "\n Tenseur2HHHH::TransfertDunTenseurGeneral(...");
|
|
};
|
|
};
|
|
|
|
// ---- manipulation d'indice ----
|
|
// on baisse les deux premiers indices -> création d'un tenseur TenseurBBHH
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
TenseurBBHH& Tenseur2HHHH::Baisse2premiersIndices()
|
|
{ TenseurBBHH * res;
|
|
res = new Tenseur2BBHH;
|
|
LesMaillonsBBHH::NouveauMaillon(res); // ajout d'un tenseur intermediaire
|
|
for (int ijkl=0;ijkl<9;ijkl++)
|
|
res->t[ijkl] = t[ijkl] ;
|
|
return *res;
|
|
};
|
|
// on baisse les deux derniers indices -> création d'un tenseur TenseurHHBB
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
TenseurHHBB& Tenseur2HHHH::Baisse2derniersIndices()
|
|
{ TenseurHHBB * res;
|
|
res = new Tenseur2HHBB;
|
|
LesMaillonsHHBB::NouveauMaillon(res); // ajout d'un tenseur intermediaire
|
|
for (int ijkl=0;ijkl<9;ijkl++)
|
|
res->t[ijkl] = t[ijkl] ;
|
|
return *res;
|
|
};
|
|
|
|
// calcul des composantes locales du tenseur considéré s'exprimé dans une base absolue
|
|
// en argument : A -> une reference sur le tenseur résultat qui peut avoir une dimension
|
|
// différente du tenseur courant suivant que la dimension absolue et la dimension locale
|
|
// sont égales ou différentes , retour d'une reference sur A
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
TenseurHHHH & Tenseur2HHHH::Baselocale(TenseurHHHH & A,const BaseH & gi) const
|
|
{ return produit44_2(A,*this,gi);
|
|
};
|
|
|
|
// changement des composantes du tenseur, retour donc dans la même variance
|
|
// en argument : A -> une reference sur le tenseur résultat qui a la même dimension
|
|
// retour d'une reference sur A
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
TenseurHHHH & Tenseur2HHHH::ChangeBase(TenseurHHHH & A,const BaseB & gi) const
|
|
{//cout << "\n debug: Tenseur2HHHH::ChangeBase(... "
|
|
// << " debut " << flush;
|
|
return produit55_2(A,*this,gi);
|
|
//cout << "\n debug: Tenseur2HHHH::ChangeBase(... "
|
|
// << " fin " << flush;
|
|
};
|
|
|
|
// test
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
int Tenseur2HHHH::operator == ( const TenseurHHHH & B) const
|
|
{ int res = 1;
|
|
#ifdef MISE_AU_POINT
|
|
if (Dabs(B.Dimension()) != 22) Message(2,"Tenseur2HHHH::operator == ( etc..");
|
|
#endif
|
|
for (int i = 0; i< 9; i++)
|
|
if (this->t[i] != B.t[i]) res = 0 ;
|
|
return res;
|
|
};
|
|
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
// change la composante i,j,k,l du tenseur
|
|
// acces en ecriture,
|
|
void Tenseur2HHHH::Change (int i, int j, int k, int l, const double& val)
|
|
{ t[(cdex2HHHH.odVect(i,j)-1)*3+cdex2HHHH.odVect(k,l)-1] = val;};
|
|
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
// change en cumulant la composante i,j,k,l du tenseur
|
|
// acces en ecriture,
|
|
void Tenseur2HHHH::ChangePlus (int i, int j, int k, int l, const double& val)
|
|
{ t[(cdex2HHHH.odVect(i,j)-1)*3+cdex2HHHH.odVect(k,l)-1] += val;};
|
|
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
// Retourne la composante i,j,k,l du tenseur
|
|
// acces en lecture seule
|
|
double Tenseur2HHHH::operator () (int i, int j, int k, int l) const
|
|
{ return t[(cdex2HHHH.odVect(i,j)-1)*3+cdex2HHHH.odVect(k,l)-1]; };
|
|
|
|
// calcul du maximum en valeur absolu des composantes du tenseur
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
double Tenseur2HHHH::MaxiComposante() const
|
|
{ return DabsMaxiTab(t, 9) ;
|
|
};
|
|
|
|
|
|
// lecture et écriture de données
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
istream & Tenseur2HHHH::Lecture(istream & entree)
|
|
{ // lecture et vérification du type
|
|
string nom_type;
|
|
entree >> nom_type;
|
|
if (nom_type != "Tenseur2HHHH")
|
|
{ Sortie(1);
|
|
return entree;
|
|
};
|
|
// lecture des coordonnées
|
|
for (int i = 0; i< 9; i++)
|
|
entree >> this->t[i];
|
|
return entree;
|
|
};
|
|
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
ostream & Tenseur2HHHH::Ecriture(ostream & sort) const
|
|
{ // écriture du type
|
|
sort << "Tenseur2HHHH ";
|
|
// puis les datas
|
|
for (int i = 0; i< 9; i++)
|
|
sort << setprecision(ParaGlob::NbdigdoCA()) << this->t[i] << " ";
|
|
return sort;
|
|
};
|
|
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
// affichage sous forme de tableau bidim
|
|
void Tenseur2HHHH::Affiche_bidim(ostream & sort) const
|
|
{
|
|
sort << "\n" ;
|
|
for (int kl=1;kl < 4;kl++)
|
|
sort << setw(15) << kl ;
|
|
for (int ij=1;ij < 4;ij++)
|
|
{sort << '\n'<< setw(4) << ij;
|
|
for (int kl=1;kl < 4;kl++)
|
|
{ int i= cdex2HHHH.idx_i(ij); int j= cdex2HHHH.idx_j(ij);
|
|
int k= cdex2HHHH.idx_i(kl); int l= cdex2HHHH.idx_j(kl);
|
|
sort << setw(15) << setprecision(7)
|
|
<< t[(cdex2HHHH.odVect(i,j)-1)*3+cdex2HHHH.odVect(k,l)-1];
|
|
}
|
|
sort << '\n';
|
|
}
|
|
cout << endl ;
|
|
};
|
|
|
|
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
// surcharge de l'operator de lecture
|
|
istream & operator >> (istream & entree, Tenseur2HHHH & A)
|
|
{ int dim = A.Dimension();
|
|
#ifdef MISE_AU_POINT
|
|
if (dim != 22) A.Message(2,"operator >> (istream & entree, Tenseur2HHHH & A)");
|
|
#endif
|
|
// lecture et vérification du type
|
|
string nom_type;
|
|
entree >> nom_type;
|
|
if (nom_type != "Tenseur2HHHH")
|
|
{ Sortie(1);
|
|
return entree;
|
|
}
|
|
// lecture des coordonnées
|
|
for (int i = 0; i< 9; i++)
|
|
entree >> A.t[i];
|
|
return entree;
|
|
};
|
|
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
// surcharge de l'operator d'ecriture
|
|
ostream & operator << (ostream & sort , const Tenseur2HHHH & A)
|
|
{ int dim = A.Dimension();
|
|
// écriture du type
|
|
sort << "Tenseur2HHHH ";
|
|
// puis les datas
|
|
for (int i = 0; i< 9; i++)
|
|
sort << setprecision(ParaGlob::NbdigdoCA()) << A.t[i] << " ";
|
|
return sort;
|
|
};
|
|
|
|
//=========== fonction protected ======================
|
|
// fonction pour le produit contracté à gauche
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
TenseurHH& Tenseur2HHHH::Prod_gauche( const TenseurBB & aBB) const
|
|
{
|
|
#ifdef MISE_AU_POINT
|
|
if (Dabs(aBB.Dimension()) != 2)
|
|
Message(2,"TenseurQ2geneHHHH::Prod_gauche( const TenseurBB & F)");
|
|
#endif
|
|
TenseurHH * res;
|
|
res = new Tenseur2HH;
|
|
LesMaillonsHH::NouveauMaillon( res); // ajout d'un tenseur intermediaire
|
|
const Tenseur2BB & a2BB = *((Tenseur2BB*) &aBB); // passage en dim 2
|
|
|
|
for (int kl=1;kl < 4;kl++)
|
|
{for (int ij=1;ij < 3;ij++) // partie simple produit : la partie diagonale
|
|
(*res).Coor(cdex2HHHH.idx_i(kl),cdex2HHHH.idx_j(kl)) +=
|
|
a2BB(cdex2HHHH.idx_i(ij),cdex2HHHH.idx_j(ij)) * t[(ij-1)*3+kl-1];
|
|
int ij = 3; // partie produit doublée : partie extra-diagonale
|
|
// comme le tenseur du second ordre n'est pas forcément symétrique on utilise les
|
|
// 2 coordonnées de chaque coté de la diagonale,
|
|
(*res).Coor(cdex2HHHH.idx_i(kl),cdex2HHHH.idx_j(kl)) +=
|
|
(a2BB(cdex2HHHH.idx_i(ij),cdex2HHHH.idx_j(ij))
|
|
+a2BB(cdex2HHHH.idx_j(ij),cdex2HHHH.idx_i(ij))
|
|
)
|
|
* t[(ij-1)*3+kl-1] ;
|
|
};
|
|
|
|
// // pour vérif
|
|
// for (int i=1;i<3;i++) for (int j=1;j<3;j++)
|
|
// {res->Coor(i,j)=0;
|
|
// for (int k=1;k<3;k++) for (int l=1;l<3;l++)
|
|
// res->Coor(i,j) += a2BB(i,j) * (*this)(i,j,k,l) ;
|
|
// };
|
|
|
|
return *res ;
|
|
};
|
|
//=========== fin fonction protected ======================
|
|
|
|
//
|
|
//------------------------------------------------------------------
|
|
// cas des composantes 4 fois covariantes
|
|
//------------------------------------------------------------------
|
|
// --- gestion de changement d'index ----
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
Tenseur2BBBB::ChangementIndex::ChangementIndex() :
|
|
idx_i(3),idx_j(3),odVect(2)
|
|
{ idx_i(1)=1;idx_i(2)=2; idx_j(1)=1;idx_j(2)=2;
|
|
idx_i(3)=1; idx_j(3)=2;
|
|
odVect(1,1)=1;odVect(1,2)=3;
|
|
odVect(2,1)=3;odVect(2,2)=2;
|
|
};
|
|
// variables globales
|
|
//Tenseur2BBBB::ChangementIndex Tenseur2BBBB::cdex2BBBB;
|
|
|
|
|
|
// Constructeur
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
Tenseur2BBBB::Tenseur2BBBB() :
|
|
ipointe() // par défaut
|
|
{ dimension = 22;
|
|
listdouble9.push_front(Reels9()); // allocation
|
|
ipointe = listdouble9.begin(); // recup de la position de la maille dans la liste
|
|
t = (ipointe)->donnees; // recup de la position des datas dans la maille
|
|
for (int i=0;i<9;i++) t[i]=0.;
|
|
};
|
|
// initialisation de toutes les composantes a une meme valeur
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
Tenseur2BBBB::Tenseur2BBBB( const double val) :
|
|
ipointe()
|
|
{ dimension = 22;
|
|
listdouble9.push_front(Reels9()); // allocation
|
|
ipointe = listdouble9.begin(); // recup de la position de la maille dans la liste
|
|
t = (ipointe)->donnees; // recup de la position des datas dans la maille
|
|
for (int i=0;i<9;i++) t[i]=val;
|
|
};
|
|
// initialisation à partir d'un produit tensoriel avec 2 cas
|
|
// booleen = true : produit tensoriel normal
|
|
// *this=aBB(i,j).bBB(k,l) gHi gHj gHk gHl
|
|
// booleen = false : produit tensoriel barre
|
|
// *this(i,j,k,l)
|
|
// = 1/4 * (a(i,k).b(j,l) + a(j,k).b(i,l) + a(i,l).b(j,k) + a(j,l).b(i,k))
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
Tenseur2BBBB::Tenseur2BBBB(bool normal, const TenseurBB & aBB, const TenseurBB & bBB) :
|
|
ipointe()
|
|
{ dimension = 22;
|
|
listdouble9.push_front(Reels9()); // allocation
|
|
ipointe = listdouble9.begin(); // recup de la position de la maille dans la liste
|
|
t = (ipointe)->donnees; // recup de la position des datas dans la maille
|
|
const Tenseur2BB & a2BB = *((Tenseur2BB*) &aBB); // passage en dim 2
|
|
const Tenseur2BB & b2BB = *((Tenseur2BB*) &bBB); // passage en dim 2
|
|
#ifdef MISE_AU_POINT
|
|
if (Dabs(a2BB.Dimension()) != 2)
|
|
Message(2,string("produit tensoriel a partir d'un premier tenseur non symetriques \n")+
|
|
"Tenseur2BBBB::Tenseur2BBBB(bool normal, const"
|
|
+ " TenseurBB & aBB, const TenseurBB & bBB);");
|
|
if (Dabs(b2BB.Dimension()) != 2)
|
|
Message(2,string("produit tensoriel a partir d'un second tenseur non symétriques \n")+
|
|
"Tenseur2BBBB::Tenseur2BBBB(bool normal, const"
|
|
+ " TenseurBB & aBB, const TenseurBB & bBB);");
|
|
#endif
|
|
if (normal)
|
|
{ for (int ij=1;ij < 4;ij++)
|
|
for (int kl=1;kl < 4;kl++)
|
|
t[(ij-1)*3+kl-1] = a2BB(cdex2BBBB.idx_i(ij),cdex2BBBB.idx_j(ij))
|
|
* b2BB(cdex2BBBB.idx_i(kl),cdex2BBBB.idx_j(kl));
|
|
}
|
|
else
|
|
{ for (int ij=1;ij < 4;ij++)
|
|
for (int kl=1;kl < 4;kl++)
|
|
{ int i = cdex2BBBB.idx_i(ij); int j = cdex2BBBB.idx_j(ij);
|
|
int k = cdex2BBBB.idx_i(kl); int l = cdex2BBBB.idx_j(kl);
|
|
t[(ij-1)*3+kl-1] = 0.25 * (
|
|
a2BB(i,k) * b2BB(j,l) + a2BB(j,k) * b2BB(i,l)
|
|
+a2BB(i,l) * b2BB(j,k) + a2BB(j,l) * b2BB(i,k)
|
|
);
|
|
};
|
|
}
|
|
};
|
|
// DESTRUCTEUR :
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
Tenseur2BBBB::~Tenseur2BBBB()
|
|
{ listdouble9.erase(ipointe);} ; // suppression de l'élément de la liste
|
|
// constructeur a partir d'une instance non differenciee
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
Tenseur2BBBB::Tenseur2BBBB ( const TenseurBBBB & B) :
|
|
ipointe()
|
|
{ dimension = 22;
|
|
// #ifdef MISE_AU_POINT
|
|
// if (Dabs(dimension) != 22)
|
|
// { cout << "\n erreur de dimension, elle devrait etre = 22 ";
|
|
// cout << "\n Tenseur2BBBB::Tenseur2BBBB ( TenseurBBBB &) " << endl;
|
|
// Sortie(1);
|
|
// }
|
|
// #endif
|
|
listdouble9.push_front(Reels9()); // allocation
|
|
ipointe = listdouble9.begin(); // recup de la position de la maille dans la liste
|
|
t = (ipointe)->donnees; // recup de la position des datas dans la maille
|
|
if (Dabs(B.dimension) == 22 ) // cas d'un tenseur du même type
|
|
{ for (int i=0;i< 9;i++)
|
|
t[i] = B.t[i];
|
|
}
|
|
else
|
|
{// cas d'un tenseur quelconque
|
|
double Z=B.MaxiComposante();
|
|
for (int i=1;i < 3;i++)
|
|
for (int j=1;j<=i;j++)
|
|
for (int k=1;k < 3;k++)
|
|
for (int l=1;l<=k;l++)
|
|
{// on teste les symétries et on affecte
|
|
double a = B(i,j,k,l);
|
|
#ifdef MISE_AU_POINT
|
|
if ((!diffpourcent(a,B(j,i,k,l),Z,ConstMath::unpeupetit)
|
|
&& !diffpourcent(a,B(i,j,l,k),Z,ConstMath::unpeupetit))
|
|
|| (Abs(Z) < ConstMath::trespetit) )
|
|
// erreur d'affectation
|
|
if (ParaGlob::NiveauImpression() > 5)
|
|
cout << "\n tenseurBBBB (ijkl= " << i << "," << j << "," << k << "," << l << ")= "
|
|
<< a << " " << B(j,i,k,l) << " " <<B(i,j,l,k) ;
|
|
if (ParaGlob::NiveauImpression() > 1)
|
|
cout << "WARNING ** erreur constructeur, tenseur non symetrique, Tenseur2BBBB::Tenseur2BBBB(const TenseurBBBB & B)";
|
|
#endif
|
|
// si il y a un pb normalement il y a eu un message
|
|
this->Change(i,j,k,l,a);
|
|
}
|
|
};
|
|
};
|
|
// constructeur de copie
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
Tenseur2BBBB::Tenseur2BBBB ( const Tenseur2BBBB & B) :
|
|
ipointe()
|
|
{ this->dimension = B.dimension;
|
|
listdouble9.push_front(Reels9()); // allocation
|
|
ipointe = listdouble9.begin(); // recup de la position de la maille dans la liste
|
|
t = (ipointe)->donnees; // recup de la position des datas dans la maille
|
|
for (int i=0;i< 9;i++)
|
|
this->t[i] = B.t[i];
|
|
};
|
|
// METHODES PUBLIQUES :
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
// initialise toutes les composantes à val
|
|
void Tenseur2BBBB::Inita(double val)
|
|
{ for (int i=0;i< 9;i++)
|
|
t[i] = val;
|
|
};
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
TenseurBBBB & Tenseur2BBBB::operator + ( const TenseurBBBB & B) const
|
|
{ TenseurBBBB * res;
|
|
#ifdef MISE_AU_POINT
|
|
if (B.Dimension() != 22) Message(2,"Tenseur2BBBB::operator + ( etc..");
|
|
#endif
|
|
res = new Tenseur2BBBB;
|
|
LesMaillonsBBBB::NouveauMaillon( res); // ajout d'un tenseur intermediaire
|
|
for (int i = 0; i< 9; i++)
|
|
res->t[i] = this->t[i] + B.t[i]; //somme des données
|
|
return *res ;};
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
void Tenseur2BBBB::operator += ( const TenseurBBBB & B)
|
|
{
|
|
#ifdef MISE_AU_POINT
|
|
if (Dabs(B.Dimension()) != 22) Message(2,"Tenseur2BBBB::operator += ( etc..");
|
|
#endif
|
|
for (int i = 0; i< 9; i++)
|
|
this->t[i] += B.t[i];
|
|
LesMaillonsBBBB::Libere(); // destruction des tenseurs intermediaires
|
|
}; //somme des données
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
TenseurBBBB & Tenseur2BBBB::operator - () const
|
|
{ TenseurBBBB * res;
|
|
res = new Tenseur2BBBB;
|
|
LesMaillonsBBBB::NouveauMaillon( res); // ajout d'un tenseur intermediaire
|
|
for (int i = 0; i< 9; i++)
|
|
res->t[i] = - this->t[i]; //oppose
|
|
return *res ;};
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
TenseurBBBB & Tenseur2BBBB::operator - ( const TenseurBBBB & B) const
|
|
{ TenseurBBBB * res;
|
|
#ifdef MISE_AU_POINT
|
|
if (Dabs(B.Dimension()) != 22) Message(2,"Tenseur2BBBB::operator - ( etc..");
|
|
#endif
|
|
res = new Tenseur2BBBB;
|
|
LesMaillonsBBBB::NouveauMaillon( res); // ajout d'un tenseur intermediaire
|
|
for (int i = 0; i< 9; i++)
|
|
res->t[i] = this->t[i] - B.t[i]; //soustraction des données
|
|
return *res ;};
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
void Tenseur2BBBB::operator -= ( const TenseurBBBB & B)
|
|
{
|
|
#ifdef MISE_AU_POINT
|
|
if (Dabs(B.Dimension()) != 22) Message(2,"Tenseur2BBBB::operator -= ( etc..");
|
|
#endif
|
|
for (int i = 0; i< 9; i++)
|
|
this->t[i] -= B.t[i];
|
|
LesMaillonsBBBB::Libere(); // destruction des tenseurs intermediaires
|
|
}; //soustraction des données
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
TenseurBBBB & Tenseur2BBBB::operator = ( const TenseurBBBB & B)
|
|
{
|
|
#ifdef MISE_AU_POINT
|
|
if (Dabs(B.Dimension()) != 22) Message(2,"Tenseur2BBBB::operator = ( etc..");
|
|
#endif
|
|
for (int i = 0; i< 9; i++)
|
|
this->t[i] = B.t[i];
|
|
LesMaillonsBBBB::Libere(); // destruction des tenseurs intermediaires
|
|
return *this;
|
|
}; //affectation des données;
|
|
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
TenseurBBBB & Tenseur2BBBB::operator * ( const double & b) const
|
|
{ TenseurBBBB * res;
|
|
res = new Tenseur2BBBB;
|
|
LesMaillonsBBBB::NouveauMaillon( res); // ajout d'un tenseur intermediaire
|
|
for (int i = 0; i< 9; i++)
|
|
res->t[i] = this->t[i] * b; //multiplication des données
|
|
return *res ;};
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
void Tenseur2BBBB::operator *= ( const double & b)
|
|
{for (int i = 0; i< 9; i++)
|
|
this->t[i] *= b ;}; //multiplication des données
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
TenseurBBBB & Tenseur2BBBB::operator / ( const double & b) const
|
|
{ TenseurBBBB * res;
|
|
res = new Tenseur2BBBB;
|
|
LesMaillonsBBBB::NouveauMaillon( res); // ajout d'un tenseur intermediaire
|
|
#ifdef MISE_AU_POINT
|
|
if (Dabs(b) < ConstMath::trespetit)
|
|
{ cout << "\n erreur le diviseur est trop petit = " << b;
|
|
cout << "\n Tenseur2BBBB::operator / ( const double & b) " << endl;
|
|
Sortie(1);
|
|
}
|
|
#endif
|
|
for (int i = 0; i< 9; i++)
|
|
res->t[i] = this->t[i] / b; //division des données
|
|
return *res ;};
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
void Tenseur2BBBB::operator /= ( const double & b)
|
|
{
|
|
#ifdef MISE_AU_POINT
|
|
if (Dabs(b) < ConstMath::trespetit)
|
|
{ cout << "\n erreur le diviseur est trop petit = " << b;
|
|
cout << "\n Tenseur2BBBB::operator /= ( const double & b) " << endl;
|
|
Sortie(1);
|
|
}
|
|
#endif
|
|
for (int i = 0; i< 9; i++)
|
|
this->t[i] /= b ;}; //division des données
|
|
|
|
|
|
// produit contracte à droite avec un tenseur du second ordre
|
|
// différent à gauche !!
|
|
// différent à gauche !!
|
|
// A(i,j,k,l)*B(l,k)=A..B
|
|
// on commence par contracter l'indice du milieu puis externe
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
TenseurBB& Tenseur2BBBB::operator && ( const TenseurHH & aHH) const
|
|
{
|
|
#ifdef MISE_AU_POINT
|
|
if (Dabs(aHH.Dimension()) != 2)
|
|
Message(2,"Tenseur2BBBB::operator && ( const TenseurHH & aHH)");
|
|
#endif
|
|
TenseurBB * res;
|
|
res = new Tenseur2BB;
|
|
LesMaillonsBB::NouveauMaillon( res); // ajout d'un tenseur intermediaire
|
|
const Tenseur2HH & a2HH = *((Tenseur2HH*) &aHH); // passage en dim 2
|
|
for (int ij=1;ij < 4;ij++)
|
|
{for (int kl=1;kl < 3;kl++) // partie simple produit : la partie diagonale
|
|
(*res).Coor(cdex2BBBB.idx_i(ij),cdex2BBBB.idx_j(ij)) += t[(ij-1)*3+kl-1]
|
|
* a2HH(cdex2BBBB.idx_i(kl),cdex2BBBB.idx_j(kl));
|
|
int kl = 3; // partie produit doublée : partie extra-diagonale
|
|
// comme le tenseur du second ordre n'est pas forcément symétrique on utilise les
|
|
// 2 coordonnées de chaque coté de la diagonale,
|
|
(*res).Coor(cdex2BBBB.idx_i(ij),cdex2BBBB.idx_j(ij)) += t[(ij-1)*3+kl-1]
|
|
* (a2HH(cdex2BBBB.idx_i(kl),cdex2BBBB.idx_j(kl))
|
|
+a2HH(cdex2BBBB.idx_j(kl),cdex2BBBB.idx_i(kl)));
|
|
};
|
|
return *res ;
|
|
};
|
|
|
|
|
|
//fonctions définissant le produit tensoriel normal de deux tenseurs
|
|
// *this=aBB(i,j).bBB(k,l) gHi gHj gHk gHl
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
TenseurBBBB & Tenseur2BBBB::Prod_tensoriel(const TenseurBB & aBB, const TenseurBB & bBB)
|
|
{ TenseurBBBB * res;
|
|
res = new Tenseur2BBBB;
|
|
LesMaillonsBBBB::NouveauMaillon( res); // ajout d'un tenseur intermediaire
|
|
const Tenseur2BB & a2BB = *((Tenseur2BB*) &aBB); // passage en dim 2
|
|
const Tenseur2BB & b2BB = *((Tenseur2BB*) &bBB); // passage en dim 2
|
|
#ifdef MISE_AU_POINT
|
|
if (Dabs(a2BB.Dimension()) != 2)
|
|
{ cout << "\n produit tensoriel a partir d'un premier tenseur non symétriques \n"
|
|
<< "Tenseur2BBBB::Prod_tensoriel(const TenseurBB & aBB, const TenseurBB & bBB)";
|
|
Sortie(2);
|
|
}
|
|
if (Dabs(b2BB.Dimension()) != 2)
|
|
{ cout << "\n produit tensoriel a partir d'un second tenseur non symétriques \n"
|
|
<< "Tenseur2BBBB::Prod_tensoriel(const TenseurBB & aBB, const TenseurBB & bBB)";
|
|
Sortie(2);
|
|
}
|
|
#endif
|
|
for (int ij=1;ij < 4;ij++)
|
|
for (int kl=1;kl < 4;kl++)
|
|
res->t[(ij-1)*3+kl-1] = a2BB(cdex2BBBB.idx_i(ij),cdex2BBBB.idx_j(ij))
|
|
* b2BB(cdex2BBBB.idx_i(kl),cdex2BBBB.idx_j(kl));
|
|
return *res;
|
|
};
|
|
//fonctions définissant le produit tensoriel barre de deux tenseurs
|
|
// concervant les symétries !!
|
|
// *this(i,j,k,l)
|
|
// = 1/4 * (a(i,k).b(j,l) + a(j,k).b(i,l) + a(i,l).b(j,k) + a(j,l).b(i,k))
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
TenseurBBBB & Tenseur2BBBB::Prod_tensoriel_barre(const TenseurBB & aBB, const TenseurBB & bBB)
|
|
{ TenseurBBBB * res;
|
|
res = new Tenseur2BBBB;
|
|
LesMaillonsBBBB::NouveauMaillon( res); // ajout d'un tenseur intermediaire
|
|
const Tenseur2BB & a2BB = *((Tenseur2BB*) &aBB); // passage en dim 2
|
|
const Tenseur2BB & b2BB = *((Tenseur2BB*) &bBB); // passage en dim 2
|
|
#ifdef MISE_AU_POINT
|
|
if (Dabs(a2BB.Dimension()) != 2)
|
|
{ cout << "\n produit tensoriel a partir d'un premier tenseur non symétriques \n"
|
|
<< "Tenseur2BBBB::Prod_tensoriel(const TenseurBB & aBB, const TenseurBB & bBB)";
|
|
Sortie(2);
|
|
}
|
|
if (Dabs(b2BB.Dimension()) != 2)
|
|
{ cout << "\n produit tensoriel a partir d'un second tenseur non symétriques \n"
|
|
<< "Tenseur2BBBB::Prod_tensoriel(const TenseurBB & aBB, const TenseurBB & bBB)";
|
|
Sortie(2);
|
|
}
|
|
#endif
|
|
for (int ij=1;ij < 4;ij++)
|
|
for (int kl=1;kl < 4;kl++)
|
|
{ int i = cdex2BBBB.idx_i(ij); int j = cdex2BBBB.idx_j(ij);
|
|
int k = cdex2BBBB.idx_i(kl); int l = cdex2BBBB.idx_j(kl);
|
|
res->t[(ij-1)*3+kl-1] = 0.25 * (
|
|
a2BB(i,k) * b2BB(j,l) + a2BB(j,k) * b2BB(i,l)
|
|
+a2BB(i,l) * b2BB(j,k) + a2BB(j,l) * b2BB(i,k)
|
|
);
|
|
};
|
|
return *res;
|
|
};
|
|
|
|
// ATTENTION creation d'un tenseur transpose qui est supprime par Libere
|
|
// les 2 premiers indices sont échangés avec les deux derniers indices
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
TenseurBBBB & Tenseur2BBBB::Transpose1et2avec3et4() const
|
|
{ TenseurBBBB * res;
|
|
res = new Tenseur2BBBB;
|
|
LesMaillonsBBBB::NouveauMaillon(res); // ajout d'un tenseur intermediaire
|
|
for (int ij=1;ij < 4;ij++)
|
|
for (int kl=1;kl < 4;kl++)
|
|
res->t[(kl-1)*3+ij-1] = t[(ij-1)*3+kl-1] ;
|
|
return *res;
|
|
};
|
|
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
// affectation de B dans this, plusZero = false: les données manquantes sont inchangées,
|
|
// plusZero = true: les données manquantes sont mises à 0
|
|
// si au contraire la dimension de B est plus grande que *this, il y a uniquement affectation
|
|
// des données possibles
|
|
void Tenseur2BBBB::Affectation_trans_dimension(const TenseurBBBB & aBBBB,bool plusZero)
|
|
{ switch (abs(aBBBB.Dimension()))
|
|
{ case 33 : case 22 :
|
|
for (int ij=1;ij < 4;ij++)
|
|
for (int kl=1;kl < 4;kl++)
|
|
t[(ij-1)*3+kl-1] = aBBBB(cdex2BBBB.idx_i(ij),cdex2BBBB.idx_j(ij)
|
|
,cdex2BBBB.idx_i(kl),cdex2BBBB.idx_j(kl));
|
|
break;
|
|
case 11 :
|
|
if (plusZero)
|
|
this->Inita(0.); // on commence par mettre à 0 si besoin
|
|
// ensuite on affecte
|
|
t[0] = aBBBB(1,1,1,1);
|
|
break;
|
|
default:
|
|
Message(2,string(" *** erreur, la dimension: ")
|
|
+ ChangeEntierSTring(abs(aBBBB.Dimension()))
|
|
+"n'est pas prise en compte \n Tenseur2BBBB::Affectation_trans_dimension(");
|
|
};
|
|
};
|
|
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
// transférer un tenseur général de même dimension accessible en indice, dans un tenseur 9 Tenseur2HHHH
|
|
// il n'y a pas de symétrisation !, seule certaines composantes sont prises en compte
|
|
void Tenseur2BBBB::TransfertDunTenseurGeneral(const TenseurBBBB & aBBBB)
|
|
{ switch (aBBBB.Dimension())
|
|
{ case 11 :
|
|
{int ij=1;
|
|
int kl=1;
|
|
t[(ij-1)*3+kl-1] = aBBBB(1,1,1,1);
|
|
break;
|
|
}
|
|
case 33 : case 22 :
|
|
{for (int ij=1;ij < 4;ij++)
|
|
for (int kl=1;kl < 4;kl++)
|
|
t[(ij-1)*3+kl-1] = aBBBB(cdex2BBBB.idx_i(ij),cdex2BBBB.idx_j(ij)
|
|
,cdex2BBBB.idx_i(kl),cdex2BBBB.idx_j(kl));
|
|
break;
|
|
}
|
|
default:
|
|
Message(2,"\n erreur *** aBBBB.Dimension()= " + ChangeEntierSTring(aBBBB.Dimension())
|
|
+ "\n Tenseur2BBBB::TransfertDunTenseurGeneral(...");
|
|
};
|
|
};
|
|
|
|
// ---- manipulation d'indice ----
|
|
// on monte les deux premiers indices -> création d'un tenseur TenseurHHBB
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
TenseurHHBB& Tenseur2BBBB::Monte2premiersIndices()
|
|
{ TenseurHHBB * res;
|
|
res = new Tenseur2HHBB;
|
|
LesMaillonsHHBB::NouveauMaillon(res); // ajout d'un tenseur intermediaire
|
|
for (int ijkl=0;ijkl<9;ijkl++)
|
|
res->t[ijkl] = t[ijkl] ;
|
|
return *res;
|
|
};
|
|
// on monte les deux derniers indices -> création d'un tenseur TenseurBBHH
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
TenseurBBHH& Tenseur2BBBB::Monte2derniersIndices()
|
|
{ TenseurBBHH * res;
|
|
res = new Tenseur2BBHH;
|
|
LesMaillonsBBHH::NouveauMaillon(res); // ajout d'un tenseur intermediaire
|
|
for (int ijkl=0;ijkl<9;ijkl++)
|
|
res->t[ijkl] = t[ijkl] ;
|
|
return *res;
|
|
};
|
|
|
|
// on monte les 4 indices -> création d'un tenseur TenseurHHHH
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
TenseurHHHH& Tenseur2BBBB::Monte4Indices()
|
|
{ TenseurHHHH * res;
|
|
res = new Tenseur2HHHH;
|
|
LesMaillonsHHHH::NouveauMaillon(res); // ajout d'un tenseur intermediaire
|
|
for (int ijkl=0;ijkl<9;ijkl++)
|
|
res->t[ijkl] = t[ijkl] ;
|
|
return *res;
|
|
};
|
|
|
|
|
|
// calcul des composantes locales du tenseur considéré s'exprimé dans une base absolue
|
|
// en argument : A -> une reference sur le tenseur résultat qui peut avoir une dimension
|
|
// différente du tenseur courant suivant que la dimension absolue et la dimension locale
|
|
// sont égales ou différentes , retour d'une reference sur A
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
TenseurBBBB & Tenseur2BBBB::Baselocale(TenseurBBBB & A,const BaseB & gi) const
|
|
{ return produit44_2(A,*this,gi);
|
|
};
|
|
// changement des composantes du tenseur, retour donc dans la même variance
|
|
// en argument : A -> une reference sur le tenseur résultat qui a la même dimension
|
|
// retour d'une reference sur A
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
TenseurBBBB & Tenseur2BBBB::ChangeBase(TenseurBBBB & A,const BaseH & gi) const
|
|
{ return produit55_2(A,*this,gi);
|
|
};
|
|
|
|
// test
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
int Tenseur2BBBB::operator == ( const TenseurBBBB & B) const
|
|
{ int res = 1;
|
|
#ifdef MISE_AU_POINT
|
|
if (Dabs(B.Dimension()) != 22) Message(2,"Tenseur2BBBB::operator == ( etc..");
|
|
#endif
|
|
for (int i = 0; i< 9; i++)
|
|
if (this->t[i] != B.t[i]) res = 0 ;
|
|
return res;
|
|
};
|
|
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
// change la composante i,j,k,l du tenseur
|
|
// acces en ecriture,
|
|
void Tenseur2BBBB::Change (int i, int j, int k, int l,const double& val)
|
|
{ t[(cdex2BBBB.odVect(i,j)-1)*3+cdex2BBBB.odVect(k,l)-1] = val;};
|
|
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
// change en cumulant la composante i,j,k,l du tenseur
|
|
// acces en ecriture,
|
|
void Tenseur2BBBB::ChangePlus (int i, int j, int k, int l,const double& val)
|
|
{ t[(cdex2BBBB.odVect(i,j)-1)*3+cdex2BBBB.odVect(k,l)-1] += val;};
|
|
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
// Retourne la composante i,j,k,l du tenseur
|
|
// acces en lecture seule
|
|
double Tenseur2BBBB::operator () (int i, int j, int k, int l) const
|
|
{ return t[(cdex2BBBB.odVect(i,j)-1)*3+cdex2BBBB.odVect(k,l)-1]; };
|
|
|
|
// calcul du maximum en valeur absolu des composantes du tenseur
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
double Tenseur2BBBB::MaxiComposante() const
|
|
{ return DabsMaxiTab(t, 9) ;
|
|
};
|
|
|
|
|
|
// lecture et écriture de données
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
istream & Tenseur2BBBB::Lecture(istream & entree)
|
|
{ // lecture et vérification du type
|
|
string nom_type;
|
|
entree >> nom_type;
|
|
if (nom_type != "Tenseur2BBBB")
|
|
{ Sortie(1);
|
|
return entree;
|
|
}
|
|
// lecture des coordonnées
|
|
for (int i = 0; i< 9; i++)
|
|
entree >> this->t[i];
|
|
return entree;
|
|
};
|
|
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
ostream & Tenseur2BBBB::Ecriture(ostream & sort) const
|
|
{ // écriture du type
|
|
sort << "Tenseur2BBBB ";
|
|
// puis les datas
|
|
for (int i = 0; i< 9; i++)
|
|
sort << setprecision(ParaGlob::NbdigdoCA()) << this->t[i] << " ";
|
|
return sort;
|
|
};
|
|
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
// affichage sous forme de tableau bidim
|
|
void Tenseur2BBBB::Affiche_bidim(ostream & sort) const
|
|
{
|
|
sort << "\n" ;
|
|
for (int kl=1;kl < 4;kl++)
|
|
sort << setw(15) << kl ;
|
|
for (int ij=1;ij < 4;ij++)
|
|
{sort << '\n'<< setw(4) << ij;
|
|
for (int kl=1;kl < 4;kl++)
|
|
{ int i= cdex2BBBB.idx_i(ij); int j= cdex2BBBB.idx_j(ij);
|
|
int k= cdex2BBBB.idx_i(kl); int l= cdex2BBBB.idx_j(kl);
|
|
sort << setw(15) << setprecision(7)
|
|
<< t[(cdex2BBBB.odVect(i,j)-1)*3+cdex2BBBB.odVect(k,l)-1];
|
|
}
|
|
sort << '\n';
|
|
}
|
|
cout << endl ;
|
|
};
|
|
|
|
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
// surcharge de l'operator de lecture
|
|
istream & operator >> (istream & entree, Tenseur2BBBB & A)
|
|
{ int dim = A.Dimension();
|
|
#ifdef MISE_AU_POINT
|
|
if (dim != 22) A.Message(2,"operator >> (istream & entree, Tenseur2BBBB & A)");
|
|
#endif
|
|
// lecture et vérification du type
|
|
string nom_type;
|
|
entree >> nom_type;
|
|
if (nom_type != "Tenseur2BBBB")
|
|
{ Sortie(1);
|
|
return entree;
|
|
}
|
|
// lecture des coordonnées
|
|
for (int i = 0; i< 9; i++)
|
|
entree >> A.t[i];
|
|
return entree;
|
|
};
|
|
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
// surcharge de l'operator d'ecriture
|
|
ostream & operator << (ostream & sort , const Tenseur2BBBB & A)
|
|
{ int dim = A.Dimension();
|
|
// écriture du type
|
|
sort << "Tenseur2BBBB ";
|
|
// puis les datas
|
|
for (int i = 0; i< 9; i++)
|
|
sort << setprecision(ParaGlob::NbdigdoCA()) << A.t[i] << " ";
|
|
return sort;
|
|
};
|
|
|
|
|
|
//=========== fonction protected ======================
|
|
// fonction pour le produit contracté à gauche
|
|
#ifndef MISE_AU_POINT
|
|
inline
|
|
#endif
|
|
TenseurBB& Tenseur2BBBB::Prod_gauche( const TenseurHH & aHH) const
|
|
{
|
|
#ifdef MISE_AU_POINT
|
|
if (Dabs(aHH.Dimension()) != 2)
|
|
Message(2,"Tenseur2BBBB::Prod_gauche( const TenseurHH & F)");
|
|
#endif
|
|
TenseurBB * res;
|
|
res = new Tenseur2BB;
|
|
LesMaillonsBB::NouveauMaillon( res); // ajout d'un tenseur intermediaire
|
|
const Tenseur2HH & a2HH = *((Tenseur2HH*) &aHH); // passage en dim 2
|
|
|
|
for (int kl=1;kl < 4;kl++)
|
|
{for (int ij=1;ij < 3;ij++) // partie simple produit : la partie diagonale
|
|
(*res).Coor(cdex2BBBB.idx_i(kl),cdex2BBBB.idx_j(kl)) +=
|
|
a2HH(cdex2BBBB.idx_i(ij),cdex2BBBB.idx_j(ij)) * t[(ij-1)*3+kl-1];
|
|
int ij = 3; // partie produit doublée : partie extra-diagonale
|
|
// comme le tenseur du second ordre n'est pas forcément symétrique on utilise les
|
|
// 2 coordonnées de chaque coté de la diagonale,
|
|
(*res).Coor(cdex2BBBB.idx_i(kl),cdex2BBBB.idx_j(kl)) +=
|
|
(a2HH(cdex2BBBB.idx_i(ij),cdex2BBBB.idx_j(ij))
|
|
+a2HH(cdex2BBBB.idx_j(ij),cdex2BBBB.idx_i(ij))
|
|
)
|
|
* t[(ij-1)*3+kl-1] ;
|
|
};
|
|
|
|
// // pour vérif
|
|
// for (int i=1;i<3;i++) for (int j=1;j<3;j++)
|
|
// {res->Coor(i,j)=0;
|
|
// for (int k=1;k<3;k++) for (int l=1;l<3;l++)
|
|
// res->Coor(i,j) += a2BB(i,j) * (*this)(i,j,k,l) ;
|
|
// };
|
|
|
|
|
|
return *res ;
|
|
};
|
|
//=========== fin fonction protected ======================
|
|
|
|
|
|
|
|
#endif
|