1180 lines
72 KiB
C++
1180 lines
72 KiB
C++
|
|
// 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-2022 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/>.
|
|
|
|
/************************************************************************
|
|
* DATE: 23/01/97 *
|
|
* $ *
|
|
* AUTEUR: G RIO (mailto:gerardrio56@free.fr) *
|
|
* $ *
|
|
* PROJET: Herezh++ *
|
|
* $ *
|
|
************************************************************************
|
|
* BUT: Defini l'element generique de mecanique. *
|
|
* $ *
|
|
* '''''''''''''''''''''''''''''''''''''''''''''''''''''''''''' * *
|
|
* VERIFICATION: *
|
|
* *
|
|
* ! date ! auteur ! but ! *
|
|
* ------------------------------------------------------------ *
|
|
* ! ! ! ! *
|
|
* $ *
|
|
* '''''''''''''''''''''''''''''''''''''''''''''''''''''''''''' *
|
|
* MODIFICATIONS: *
|
|
* ! date ! auteur ! but ! *
|
|
* ------------------------------------------------------------ *
|
|
* $ *
|
|
************************************************************************/
|
|
#ifndef ELEMMECA_H
|
|
#define ELEMMECA_H
|
|
|
|
#include "Element.h"
|
|
#include "Tenseur.h"
|
|
#include "NevezTenseur.h"
|
|
#include "Deformation.h"
|
|
#include "Loi_comp_abstraite.h"
|
|
#include "Enum_calcul_masse.h"
|
|
#include "Basiques.h"
|
|
#include "Enum_dure.h"
|
|
#include "CompThermoPhysiqueAbstraite.h"
|
|
#include "CompFrotAbstraite.h"
|
|
#include "LesPtIntegMecaInterne.h"
|
|
#include "Enum_StabHourglass.h"
|
|
#include "LesChargeExtSurElement.h"
|
|
#include "Temps_CPU_HZpp.h"
|
|
#include "Enum_StabMembrane.h"
|
|
|
|
|
|
/// @addtogroup groupe_des_elements_finis
|
|
/// @{
|
|
///
|
|
|
|
|
|
class ElemMeca : public Element
|
|
{
|
|
public :
|
|
// VARIABLES PUBLIQUES :
|
|
|
|
// CONSTRUCTEURS :
|
|
ElemMeca ();
|
|
// Constructeur utile quand le numero de maillage et d'identification de l'element est connu
|
|
ElemMeca (int num_maill,int num_id) ;
|
|
// Constructeur utile quand le numero de maillage et d'identification et le tableau des noeuds
|
|
// de l'element sont connu
|
|
ElemMeca (int num_maill,int num_id,const Tableau<Noeud *>& tab);
|
|
// Constructeur utile quand le numero de maillage et d'identification est connu,
|
|
// ainsi que la geometrie et le type d'interpolation de l'element
|
|
ElemMeca (int num_maill,int num_id,Enum_interpol id_interp_elt,Enum_geom id_geom_elt,string info="");
|
|
// Constructeur utile quand le numero de maillage et d'identification est connu,
|
|
// ainsi que la geometrie et le type d'interpolation de l'element
|
|
ElemMeca (int num_maill,int num_id,char* nom_interpol,char* nom_geom,string info="");
|
|
// Constructeur utile quand toutes les donnees de la classe Element sont connues
|
|
ElemMeca (int num_maill,int num_id,const Tableau<Noeud *>& tab,Enum_interpol id_interp_elt,
|
|
Enum_geom id_geom_elt,string info="");
|
|
// Constructeur utile quand toutes les donnees de la classe Element sont connues
|
|
ElemMeca (int num_maill,int num_id,const Tableau<Noeud *>& tab,char* nom_interpol,
|
|
char* nom_geom,string info="");
|
|
// Constructeur de copie
|
|
ElemMeca (const ElemMeca& elt);
|
|
|
|
// DESTRUCTEUR :
|
|
~ElemMeca ();
|
|
|
|
// METHODES PUBLIQUES :
|
|
// test si l'element est complet
|
|
// = 1 tout est ok, =0 element incomplet
|
|
int TestComplet();
|
|
|
|
// calcul si un point est a l'interieur de l'element ou non
|
|
// il faut que M est la dimension globale
|
|
// les trois fonctions sont pour l'etude a t=0, t et tdt
|
|
// retour : =0 le point est externe, =1 le point est interne ,
|
|
// = 2 le point est sur la frontière à la précision près
|
|
// coor_locales : s'il est différent de NULL, est affecté des coordonnées locales calculées,
|
|
// uniquement précises si le point est interne
|
|
int Interne_0(const Coordonnee& M,Coordonnee* coor_locales=NULL);
|
|
int Interne_t(const Coordonnee& M,Coordonnee* coor_locales=NULL);
|
|
int Interne_tdt(const Coordonnee& M,Coordonnee* coor_locales=NULL);
|
|
|
|
// récupération des énergies intégrées sur l'éléments, résultants d'un précédent calcul
|
|
// explicite, ou implicite
|
|
const EnergieMeca& EnergieTotaleElement() const {return energie_totale;};
|
|
|
|
// test pour savoir si le calcul de contrainte en absolu est possible
|
|
bool ContrainteAbsoluePossible();
|
|
|
|
// METHODES VIRTUELLES:
|
|
|
|
// retourne la liste de tous les types de ddl interne actuellement utilisés
|
|
// par l'élément (actif ou non), sont exclu de cette liste les ddl des noeuds
|
|
// reliés à l'élément (ddl implique grandeur uniquement scalaire !)
|
|
// par contre sont inclus les ddl venant de l'élément, qui sont représenté directement aux noeuds
|
|
// absolue: indique si oui ou non on sort les tenseurs dans la base absolue ou une base particulière
|
|
virtual List_io <Ddl_enum_etendu> Les_type_de_ddl_internes(bool absolue) const ;
|
|
// idem pour les grandeurs évoluées c'est-à-dire directement sous forme de vecteur, tenseurs ....
|
|
// absolue: indique si oui ou non on sort les tenseurs dans la base absolue ou une base particulière
|
|
virtual List_io <TypeQuelconque> Les_type_evolues_internes(bool absolue) const ;
|
|
// idem pour les données particulières
|
|
virtual List_io <TypeQuelconque> Les_types_particuliers_internes(bool absolue) const;
|
|
|
|
// retourne la liste de toutes les grandeurs quelconques relatives aux faces de
|
|
// l'élément (actif ou non),
|
|
// absolue: indique si oui ou non on sort les tenseurs dans la base absolue ou une base particulière
|
|
virtual List_io <TypeQuelconque> Les_type_quelconque_de_face(bool absolue) const;
|
|
// retourne la liste de toutes les grandeurs quelconques relatives aux arêtes de
|
|
// l'élément (actif ou non),
|
|
// absolue: indique si oui ou non on sort les tenseurs dans la base absolue ou une base particulière
|
|
virtual List_io <TypeQuelconque> Les_type_quelconque_de_arete(bool absolue) const;
|
|
|
|
|
|
// --------- calculs utils dans le cadre de la recherche du flambement linéaire
|
|
// dans un premier temps uniquement virtuelles, ensuite se sera virtuelle pure pour éviter
|
|
// les oubli de définition ----> AMODIFIER !!!!!!!!
|
|
// Calcul de la matrice géométrique et initiale
|
|
class MatGeomInit // pour le retour des pointeurs sur des matrices stockées
|
|
// une paire par classe d'éléments
|
|
{ public : MatGeomInit(Mat_pleine * matG,Mat_pleine * matI) :
|
|
matGeom(matG),matInit(matI) {};
|
|
Mat_pleine * matGeom ;Mat_pleine * matInit ;
|
|
};
|
|
virtual MatGeomInit MatricesGeometrique_Et_Initiale (const ParaAlgoControle & pa) ;
|
|
|
|
// --------- calcul d'erreur, calculs du champs de contrainte continu ---------
|
|
// ajout des ddl de contraintes pour les noeuds de l'élément
|
|
virtual void Plus_ddl_Sigma() = 0;
|
|
// inactive les ddl du problème de recherche d'erreur : les contraintes
|
|
virtual void Inactive_ddl_Sigma() = 0;
|
|
// active les ddl du problème de recherche d'erreur : les contraintes
|
|
virtual void Active_ddl_Sigma() = 0 ;
|
|
// active le premier ddl du problème de recherche d'erreur : SIGMA11
|
|
virtual void Active_premier_ddl_Sigma() = 0 ;
|
|
// retourne un tableau de ddl element, correspondant à la
|
|
// composante de sigma -> SIG11, pour chaque noeud qui contiend
|
|
// des ddl de contrainte
|
|
// -> utilisé pour l'assemblage de la raideur d'erreur
|
|
//!!!!!!!!!!!!! pour l'instant en virtuelle il faudra après en
|
|
// virtuelle pure !!!!!!!!!!!!!!!!!!!!!!!!!!!
|
|
virtual DdlElement& Tableau_de_Sig1() const ;
|
|
// retourne un tableau de ddl element, correspondant à la
|
|
// composante d'erreur -> ERREUR, pour chaque noeud
|
|
// -> utilisé pour l'assemblage de la raideur d'erreur
|
|
//dans cette version tous les noeuds sont supposés avoi un ddl erreur
|
|
// dans le cas contraire il faut redéfinir la fonction dans l'élément terminal
|
|
virtual DdlElement Tableau_de_ERREUR() const ;
|
|
// calcul de l'erreur sur l'élément. Ce calcul n'est disponible
|
|
// qu'une fois la remontée aux contraintes effectuées sinon aucune
|
|
// action. En retour la valeur de l'erreur sur l'élément
|
|
// type indique le type de calcul d'erreur :
|
|
// = 1 : erreur = (int (delta sigma):(delta sigma) dv)/(int sigma:sigma dv)
|
|
// le numerateur et le denominateur sont tel que :
|
|
// errElemRelative = numerateur / denominateur , si denominateur different de 0
|
|
// sinon denominateur = numerateur si numerateur est different de 0, sinon
|
|
// tous sont nuls mais on n'effectue pas la division
|
|
//!!!!!!!!!!!!! pour l'instant en virtuelle il faudra après en
|
|
// virtuelle pure !!!!!!!!!!!!!!!!!!!!!!!!!!!
|
|
virtual void ErreurElement(int type,double& errElemRelative
|
|
,double& numerateur, double& denominateur);
|
|
// les 3 routines qui suivent sont virtuelles, car la définition
|
|
//qui est faite dans ElemMeca.cp considère qu'il y a un ddl erreur
|
|
// par noeud de l'élément de manière systématique,
|
|
// avec le status virtuel on peut définir dans la classe dérivée
|
|
// un cas particulier
|
|
// ajout des ddl d'erreur pour les noeuds de l'élément
|
|
virtual void Plus_ddl_Erreur() ;
|
|
// inactive les ddl d'erreur
|
|
virtual void Inactive_ddl_Erreur() ;
|
|
// active les ddl d'erreur
|
|
virtual void Active_ddl_Erreur() ;
|
|
// test pour savoir si l'erreur a été calculée
|
|
bool ErreurDejaCalculee()
|
|
{ if (sigErreur == NULL)
|
|
return false;
|
|
else return true;};
|
|
// sortie de l'erreur à l'élément
|
|
double Erreur( )
|
|
{ return (*sigErreur);};
|
|
|
|
// lecture de données diverses sur le flot d'entrée
|
|
// l'implantation est faite dans les classe dérivées
|
|
virtual void LectureContraintes(UtilLecture * entreePrinc) =0 ;
|
|
|
|
// retour des contraintes en absolu retour true si elle existe sinon false
|
|
virtual bool ContraintesAbsolues(Tableau <Vecteur>& tabSig) = 0;
|
|
|
|
// --------- calcul dynamique ---------
|
|
|
|
// calcul de la longueur d'arrête de l'élément minimal
|
|
// divisé par la célérité dans le matériau
|
|
virtual double Long_arrete_mini_sur_c(Enum_dure temps) = 0;
|
|
// cas du bulk viscosity
|
|
double E_elem_bulk_t,E_elem_bulk_tdt,P_elem_bulk;
|
|
static void ActiveBulkViscosity(int choix) {bulk_viscosity=choix;};
|
|
static void InactiveBulkViscosity() {bulk_viscosity=false;};
|
|
static void ChangeCoefsBulkViscosity(const DeuxDoubles & coef)
|
|
{ c_traceBulk=coef.un;c_carreBulk=coef.deux;};
|
|
|
|
// initialisation pour le calcul de la matrice masse dans le cas de l'algorithme
|
|
// de relaxation dynamique avec optimisation en continu de la matrice masse
|
|
// casMass_relax: permet de choisir entre différentes méthodes de calcul de la masse
|
|
void InitCalculMatriceMassePourRelaxationDynamique(int casMass_relax);
|
|
// phase de calcul de la matrice masse dans le cas de l'algo de relaxation dynamique
|
|
// mi=fonction de (alpha*K+beta*mu+gamma*Isig/3+theta/2*Sig_mises
|
|
// ep: epaisseur, K module de compressibilite, mu: module de cisaillement, Isig trace de sigma,
|
|
// Sig_mises la contrainte de mises
|
|
// casMass_relax: permet de choisir entre différentes méthodes de calcul de la masse
|
|
void CalculMatriceMassePourRelaxationDynamique
|
|
(const double& alph, const double& beta, const double & lambda
|
|
,const double & gamma,const double & theta, int casMass_relax);
|
|
|
|
// ---------- informations annexes ----------------
|
|
|
|
// récupération de la base locales au noeud noe, pour le temps: temps
|
|
const BaseB & Gib_elemeca(Enum_dure temps, const Noeud * noe);
|
|
|
|
// récupération de grandeurs particulières au numéro d'ordre = iteg
|
|
// celles-ci peuvent être quelconques
|
|
// en retour liTQ est modifié et contiend les infos sur les grandeurs particulières
|
|
// absolue: indique si oui ou non on sort les tenseurs dans la base absolue ou une base particulière
|
|
void Grandeur_particuliere (bool absolue,List_io<TypeQuelconque>& liTQ,int iteg);
|
|
|
|
// récupération de grandeurs particulières pour une face au numéro d'ordre = iteg
|
|
// celles-ci peuvent être quelconques
|
|
// en retour liTQ est modifié et contiend les infos sur les grandeurs particulières
|
|
// absolue: indique si oui ou non on sort les tenseurs dans la base absolue ou une base particulière
|
|
void Grandeur_particuliere_face (bool absolue,List_io<TypeQuelconque>& liTQ,int face, int iteg);
|
|
|
|
// récupération de grandeurs particulières pour une arête au numéro d'ordre = iteg
|
|
// celles-ci peuvent être quelconques
|
|
// en retour liTQ est modifié et contiend les infos sur les grandeurs particulières
|
|
// absolue: indique si oui ou non on sort les tenseurs dans la base absolue ou une base particulière
|
|
void Grandeur_particuliere_arete (bool absolue,List_io<TypeQuelconque>& liTQ,int arete, int iteg);
|
|
|
|
// récupération de la loi de frottement, dans le cas où elle n'existe pas
|
|
// retour d'un pointeur nul
|
|
CompFrotAbstraite* LoiDeFrottement() const {return loiFrot;};
|
|
|
|
// récupération de la loi de comportement pour information, dans le cas où elle n'existe pas
|
|
// retour d'un pointeur nul
|
|
const Loi_comp_abstraite* LoiDeComportement() const {return loiComp;};
|
|
|
|
// récup du module de compressibilité moyen de l'élément
|
|
double CompressibiliteMoyenne() const {return lesPtIntegMecaInterne->CompressibiliteMoyenne();};
|
|
|
|
// --- transfert des grandeurs des points d'intégration aux noeuds
|
|
// transfert de ddl des points d'intégrations (de tous) d'un éléments (on ajoute aux noeuds, on ne remplace pas)
|
|
// les ddl doivent déjà exister aux noeuds sinon erreur
|
|
// il doit s'agir du même type de répartition de pt d'integ pour toutes les grandeurs
|
|
// tab_val(i)(j) : valeur associée au i ième pt d'integ et au j ième ddl_enum_etendu
|
|
void TransfertAjoutAuNoeuds(const List_io < Ddl_enum_etendu >& lietendu
|
|
,const Tableau <Tableau <double> > & tab_val,int cas);
|
|
// transfert de type quelconque des points d'intégrations (de tous) aux noeuds d'un éléments (on ajoute aux noeuds,
|
|
// on ne remplace pas). Les types quelconques doivent déjà exister
|
|
// un tableau dans tab_liQ correspondent aux grandeurs quelconque pour tous les pt integ,
|
|
// tab_liQ(i) pour le pt d'integ i
|
|
// liQ_travail: est une liste de travail qui sera utilisée dans le transfert
|
|
void TransfertAjoutAuNoeuds(const Tableau <List_io < TypeQuelconque > >& tab_liQ
|
|
,List_io < TypeQuelconque > & liQ_travail,int cas);
|
|
|
|
// accumulation aux noeuds de grandeurs venant de l'éléments vers ses noeuds (exemple la pression appliquée)
|
|
// autres que celles aux pti classiques, mais directements disponibles
|
|
// le contenu du conteneur stockées dans liQ est utilisé en variable intermédiaire
|
|
void Accumul_aux_noeuds(const List_io < Ddl_enum_etendu >& lietendu
|
|
,List_io < TypeQuelconque > & liQ,int cas);
|
|
|
|
// activation du calcul des invariants de contraintes, qui seront calculé à chaque
|
|
// fois que l'on calcul les contraintes au travers de la loi de comportement
|
|
void ActivCalculInvariantsContraintes();
|
|
// idem pour la déformation
|
|
void ActivCalculInvariantsDeformation();
|
|
// idem pour la vitesse de déformation
|
|
void ActivCalculInvariantsVitesseDeformation();
|
|
|
|
// modification de l'orientation de l'élément en fonction de cas_orientation
|
|
// =0: inversion simple (sans condition) de l'orientation
|
|
// si cas_orientation est diff de 0: on calcul le jacobien aux différents points d'intégration
|
|
// 1. si tous les jacobiens sont négatifs on change d'orientation
|
|
// 2. si tous les jacobiens sont positifs on ne fait rien
|
|
// 3. si certains jacobiens sont positifs et d'autres négatifs message
|
|
// d'erreur et on ne fait rien
|
|
// ramène true: s'il y a eu changement effectif, sinon false
|
|
bool Modif_orient_elem(int cas_orientation);
|
|
|
|
// calcul éventuel de la normale à un noeud
|
|
// ce calcul existe pour les éléments 2D, 1D axi, et aussi pour les éléments 1D
|
|
// qui possède un repère d'orientation
|
|
// en retour coor = la normale si coor.Dimension() est = à la dimension de l'espace
|
|
// si le calcul n'existe pas --> coor.Dimension() = 0
|
|
// ramène un entier :
|
|
// == 1 : calcul normal
|
|
// == 0 : problème de calcul -> coor.Dimension() = 0
|
|
// == 2 : indique que le calcul n'est pas licite pour le noeud passé en paramètre
|
|
// c'est le cas par exemple des noeuds exterieurs pour les éléments SFE
|
|
// mais il n'y a pas d'erreur, c'est seulement que l'élément n'est pas ad hoc pour
|
|
// calculer la normale à ce noeud là
|
|
// temps: indique à quel moment on veut le calcul
|
|
// pour des éléments particulier (ex: SFE) la méthode est surchargée
|
|
virtual int CalculNormale_noeud(Enum_dure temps, const Noeud& noe,Coordonnee& coor);
|
|
|
|
// calcul si un point est a l'interieur de l'element ou non
|
|
// il faut que M est la dimension globale
|
|
// retour : =0 le point est externe, =1 le point est interne ,
|
|
// = 2 le point est sur la frontière à la précision près
|
|
// coor_locales : s'il est différent de NULL, est affecté des coordonnées locales calculées,
|
|
// uniquement précises si le point est interne
|
|
int Interne(Enum_dure temps,const Coordonnee& M,Coordonnee* coor_locales=NULL);
|
|
// -- connaissances particulières sur l'élément
|
|
// ramène l'épaisseur de l'élément
|
|
// =0. si la notion d'épaisseurs ne veut rien dire pour l'élément
|
|
virtual double Epaisseurs(Enum_dure , const Coordonnee& ) {return 0.;};
|
|
// ramène l'épaisseur moyenne de l'élément (indépendante du point)
|
|
// =0. si la notion d'épaisseurs ne veut rien dire pour l'élément
|
|
virtual double EpaisseurMoyenne(Enum_dure ) {return 0.;};
|
|
// ramène la section de l'élément
|
|
// =0. si la notion de section ne veut rien dire pour l'élément
|
|
virtual double Section(Enum_dure , const Coordonnee& ) {return 0.;};
|
|
// ramène la section moyenne de l'élément (indépendante du point)
|
|
// =0. si la notion de section ne veut rien dire pour l'élément
|
|
virtual double SectionMoyenne(Enum_dure ) {return 0.;};
|
|
// // modifie l'épaisseur Moyenne à tdt
|
|
// virtual void Modifie_epaisseur_moyenne_tdt(const double& h_t)
|
|
// {cout << "\n erreur** la methode Modifie_epaisseur_moyenne(.. n'existe pas "; Sortie(1);};
|
|
|
|
// fonction a renseigner par les classes dérivées, concernant les répercutions
|
|
// éventuelles due à la suppression de tous les frontières
|
|
// nums_i : donnent les listes de frontières supprimées
|
|
virtual void Prise_en_compte_des_consequences_suppression_tous_frontieres();
|
|
// idem pour une frontière (avant qu'elle soit supprimée)
|
|
virtual void Prise_en_compte_des_consequences_suppression_une_frontiere(ElFrontiere* elemFront);
|
|
|
|
// -------------------- calcul de frontières -------------------
|
|
|
|
// ramène la frontière point
|
|
// éventuellement création des frontieres points de l'element et stockage dans l'element
|
|
// si c'est la première fois sinon il y a seulement retour de l'elements
|
|
// a moins que le paramètre force est mis a true
|
|
// dans ce dernier cas la frontière effacéee est recréée
|
|
// num indique le numéro du point à créer (numérotation EF)
|
|
virtual ElFrontiere* const Frontiere_points(int num,bool force);
|
|
|
|
// ramène la frontière linéique
|
|
// éventuellement création des frontieres linéique de l'element et stockage dans l'element
|
|
// si c'est la première fois et en 3D sinon il y a seulement retour de l'elements
|
|
// a moins que le paramètre force est mis a true
|
|
// dans ce dernier cas la frontière effacéee est recréée
|
|
// num indique le numéro de l'arête à créer (numérotation EF)
|
|
virtual ElFrontiere* const Frontiere_lineique(int num,bool force);
|
|
|
|
// ramène la frontière surfacique
|
|
// éventuellement création des frontieres surfacique de l'element et stockage dans l'element
|
|
// si c'est la première fois sinon il y a seulement retour de l'elements
|
|
// a moins que le paramètre force est mis a true
|
|
// dans ce dernier cas la frontière effacéee est recréée
|
|
// num indique le numéro de la surface à créer (numérotation EF)
|
|
virtual ElFrontiere* const Frontiere_surfacique(int num,bool force);
|
|
|
|
// Calcul spécifiques des frontieres de l'element et retour d'un tableau tabb des frontières
|
|
// creation des elements frontieres et stockage dans l'element
|
|
// la création n'a lieu qu'au premier appel
|
|
// ou lorsque l'on force le paramètre force a true
|
|
// dans ce dernier cas seul les frontière effacées sont recréée
|
|
// cas :
|
|
// = 0 -> on veut toutes les frontières
|
|
// = 1 -> on veut uniquement les surfaces
|
|
// = 2 -> on veut uniquement les lignes
|
|
// = 3 -> on veut uniquement les points
|
|
// = 4 -> on veut les surfaces + les lignes
|
|
// = 5 -> on veut les surfaces + les points
|
|
// = 6 -> on veut les lignes + les points
|
|
// virtual Tableau <ElFrontiere*> const & Frontiere_specifique(int cas, bool force = false)
|
|
// {return Frontiere_elemeca(cas,force);};
|
|
|
|
// -------------------- init éventuelle avant le chargement -------------------
|
|
// initialisation éventuelle, nécessaire avant d'appliquer l'ensemble des charges
|
|
// par exemple des stockages intermédiaires
|
|
virtual void Initialisation_avant_chargement();
|
|
|
|
// mise à jour éventuel de repère d'anisotropie
|
|
virtual void Mise_a_jour_repere_anisotropie
|
|
(BlocGen & bloc,LesFonctions_nD* lesFonctionsnD);
|
|
|
|
|
|
//=====================================================================================
|
|
protected :
|
|
//=====================================================================================
|
|
// METHODES PROTEGEES utilisables par les classes derivees :
|
|
|
|
// Calcul des frontieres de l'element
|
|
// creation des elements frontieres et retour du tableau de ces elements
|
|
// la création n'a lieu qu'au premier appel
|
|
// ou lorsque l'on force le paramètre force a true
|
|
// dans ce dernier cas seul les frontière effacées sont recréée
|
|
// cas :
|
|
// = 0 -> on veut toutes les frontières
|
|
// = 1 -> on veut uniquement les surfaces
|
|
// = 2 -> on veut uniquement les lignes
|
|
// = 3 -> on veut uniquement les points
|
|
// = 4 -> on veut les surfaces + les lignes
|
|
// = 5 -> on veut les surfaces + les points
|
|
// = 6 -> on veut les lignes + les points
|
|
Tableau <ElFrontiere*> const & Frontiere_elemeca(int cas, bool force = false);
|
|
|
|
// ---------------------------------------------------------------------
|
|
// cas où l'on intègre que selon une liste (un axe, un plan, un volume)
|
|
// ---------------------------------------------------------------------
|
|
// Calcul du residu local et de la raideur locale,
|
|
// pour le schema implicite d'ou a l'instant t + dt
|
|
// ddl represente les degres de liberte specifiques a l'element
|
|
// tabDepsBB = vitesse de déformation, tabDeltaEpsBB = incrément de def entre t et t+dt
|
|
// cald_Dvirtuelle = indique si l'on doit calculer la dérivée de la vitesse de déformation virtuelle
|
|
void Cal_implicit (DdlElement & tab_ddl,Tableau <TenseurBB *> & d_epsBB
|
|
,Tableau < Tableau2 <TenseurBB *> > d2_epsBB,Tableau <TenseurHH *>& d_sigHH,int nbint
|
|
,const Vecteur& poids,const ParaAlgoControle & pa,bool cald_Dvirtuelle);
|
|
|
|
// Calcul du residu local a l'instant t ou tdt
|
|
// atdt = true : calcul à tdt, valeur par défaut
|
|
// = false: calcul à t
|
|
// ddl represente les degres de liberte specifiques a l'element
|
|
// nbint = nb de pt d'integration , poids = poids d'integration
|
|
void Cal_explicit (DdlElement & ddl,Tableau <TenseurBB *>& d_epsBB,int nbint
|
|
,const Vecteur& poids,const ParaAlgoControle & pa,bool atdt=true);
|
|
|
|
// Calcul de la matrice géométrique et de la matrice initiale
|
|
// cette fonction est éventuellement appelée par les classes dérivées
|
|
// ddl represente les degres de liberte specifiques a l'element
|
|
// nbint = nb de pt d'integration , poids = poids d'integration
|
|
// cald_Dvirtuelle = indique si l'on doit calculer la dérivée de la vitesse de déformation virtuelle
|
|
void Cal_matGeom_Init (Mat_pleine & matGeom, Mat_pleine & matInit
|
|
,DdlElement & ddl,Tableau <TenseurBB *>& d_epsBB,Tableau < Tableau2 <TenseurBB *> > d2_epsBB
|
|
,Tableau <TenseurHH *>& d_sigHH,int nbint,const Vecteur& poids
|
|
,const ParaAlgoControle & pa,bool cald_Dvirtuelle);
|
|
|
|
// Calcul de la matrice masse selon différent choix donné par type_matrice_masse,
|
|
// a l'instant initial.
|
|
void Cal_Mat_masse (DdlElement & tab_ddl,Enum_calcul_masse type_matrice_masse,
|
|
int nbint,const Tableau <Vecteur>& taphi,int nbne
|
|
,const Vecteur& poids);
|
|
|
|
// ---------------------------------------------------------------------
|
|
// cas où l'on intègre selon deux listes (ex un axe et un plan, etc.)
|
|
// ---------------------------------------------------------------------
|
|
// Calcul du residu local et de la raideur locale,
|
|
// pour le schema implicite d'ou a l'instant t + dt
|
|
// ddl represente les degres de liberte specifiques a l'element
|
|
void Cal_implicitap (DdlElement & tab_ddl,Tableau <TenseurBB *> & d_epsBB
|
|
,Tableau < Tableau2 <TenseurBB *> > d2_epsBB,Tableau <TenseurHH *>& d_sigHH
|
|
,int nbint1,Vecteur& poids1,int nbint2,const Vecteur& poids2
|
|
,const ParaAlgoControle & pa);
|
|
|
|
// Calcul du residu local a l'instant t ou tdt
|
|
// atdt = true : calcul à tdt, valeur par défaut
|
|
// = false: calcul à t
|
|
// ddl represente les degres de liberte specifiques a l'element
|
|
// d_epsbb = variation des def
|
|
// nbint = nb de pt d'integration , poids = poids d'integration
|
|
void Cal_explicitap (DdlElement & ddl,Tableau <TenseurBB *>& d_epsBB
|
|
,int nbint,const Vecteur& poids,bool atdt=true);
|
|
|
|
// Calcul de la matrice géométrique et de la matrice initiale
|
|
// cette fonction est éventuellement appelée par les classes dérivées
|
|
// ddl represente les degres de liberte specifiques a l'element
|
|
// d_epsbb = variation des def
|
|
// nbint = nb de pt d'integration , poids = poids d'integration
|
|
void Cal_matGeom_Initap (Mat_pleine & matGeom, Mat_pleine & matInit
|
|
,DdlElement & ddl,Tableau <TenseurBB *>& d_epsBB,Tableau < Tableau2 <TenseurBB *> > d2_epsBB
|
|
,Tableau <TenseurHH *>& d_sigHH,int nbint,const Vecteur& poids);
|
|
|
|
// -------------------------- calcul de second membre ---------------------------
|
|
// calcul des seconds membres suivant les chargements
|
|
// cas d'un chargement surfacique, sur les frontières des éléments
|
|
// force indique la force surfacique appliquée
|
|
// retourne le second membre résultant
|
|
// nSurf : le numéro de la surface externe
|
|
// calcul à l'instant tdt ou t en fonction de la variable atdt
|
|
Vecteur& SM_charge_surf_E (DdlElement & ddls,int nSurf
|
|
,const Tableau <Vecteur>& taphi,int nbne
|
|
,const Vecteur& poids,const Coordonnee& force,Fonction_nD* pt_fonct
|
|
,const ParaAlgoControle & pa,bool atdt=true);
|
|
// idem SM_charge_surf_E mais -> implicite,
|
|
// pa : permet de déterminer si l'on fait ou non le calcul de la contribution à la raideur
|
|
// retourne le second membre et la matrice de raideur correspondant
|
|
Element::ResRaid SMR_charge_surf_I (DdlElement & ddls,int nSurf
|
|
,const Tableau <Vecteur>& taphi,int nbne
|
|
,const Vecteur& poids,const Coordonnee& force,Fonction_nD* pt_fonct
|
|
,const ParaAlgoControle & pa);
|
|
// calcul des seconds membres suivant les chargements
|
|
// cas d'un chargement pression, sur les frontières des éléments
|
|
// pression indique la pression appliquée
|
|
// la fonction nD est utilisée que si elle ne dépend pas strictement de grandeurs globales
|
|
// retourne le second membre résultant
|
|
// nSurf : le numéro de la surface externe
|
|
// calcul à l'instant tdt ou t en fonction de la variable atdt
|
|
Vecteur& SM_charge_pres_E (DdlElement & ddls,int nSurf
|
|
,const Tableau <Vecteur>& taphi,int nbne
|
|
,const Vecteur& poids,double pression,Fonction_nD* pt_fonct
|
|
,const ParaAlgoControle & pa,bool atdt=true);
|
|
// idem SM_charge_pres_E mais -> implicite,
|
|
// pa: permet de déterminer si oui ou non on calcul la contribution à la raideur
|
|
// retourne le second membre et la matrice de raideur correspondant
|
|
// la fonction nD est utilisée que si elle ne dépend pas strictement de grandeurs globales
|
|
Element::ResRaid SMR_charge_pres_I (DdlElement & ddls,int nSurf
|
|
,const Tableau <Vecteur>& taphi,int nbne
|
|
,const Vecteur& poids,double pression,Fonction_nD* pt_fonct
|
|
,const ParaAlgoControle & pa);
|
|
// cas d'un chargement lineique, sur les arêtes frontières des éléments
|
|
// force indique la force lineique appliquée
|
|
// retourne le second membre résultant
|
|
// nArete : le numéro de l'arête externe
|
|
// calcul à l'instant tdt ou t en fonction de la variable atdt
|
|
Vecteur& SM_charge_line_E (DdlElement & ddls,int nArete
|
|
,const Tableau <Vecteur>& taphi,int nbne
|
|
,const Vecteur& poids,const Coordonnee& force,Fonction_nD* pt_fonct
|
|
,const ParaAlgoControle & pa,bool atdt=true);
|
|
// idem SM_charge_line_E mais -> implicite,
|
|
// pa: permet de déterminer si oui ou non on calcul la contribution à la raideur
|
|
// retourne le second membre et la matrice de raideur correspondant
|
|
Element::ResRaid SMR_charge_line_I (DdlElement & ddlA,int nArete
|
|
,const Tableau <Vecteur>& taphi,int nbne
|
|
,const Vecteur& poids,const Coordonnee& force,Fonction_nD* pt_fonct
|
|
,const ParaAlgoControle & pa);
|
|
|
|
// cas d'un chargement lineique suiveur, sur les arêtes frontières des éléments
|
|
// pas valable pour des éléments 3D !
|
|
// force indique la force lineique appliquée
|
|
// retourne le second membre résultant
|
|
// nArete : le numéro de l'arête externe
|
|
// calcul à l'instant tdt ou t en fonction de la variable atdt
|
|
Vecteur& SM_charge_line_Suiv_E (DdlElement & ddls,int nArete
|
|
,const Tableau <Vecteur>& taphi,int nbne
|
|
,const Vecteur& poids,const Coordonnee& force,Fonction_nD* pt_fonct
|
|
,const ParaAlgoControle & pa,bool atdt=true);
|
|
// idem SM_charge_line_E mais -> implicite,
|
|
// pa: permet de déterminer si oui ou non on calcul la contribution à la raideur
|
|
// retourne le second membre et la matrice de raideur correspondant
|
|
Element::ResRaid SMR_charge_line_Suiv_I (DdlElement & ddlA,int nArete
|
|
,const Tableau <Vecteur>& taphi,int nbne
|
|
,const Vecteur& poids,const Coordonnee& force,Fonction_nD* pt_fonct
|
|
,const ParaAlgoControle & pa);
|
|
|
|
|
|
// cas d'un chargement surfacique suiveur, sur les surfaces de l'élément
|
|
// la direction varie selon le système suivant: on définit les coordonnées matérielles
|
|
// de la direction, ce qui sert ensuite à calculer les nouvelles directions. L'intensité
|
|
// elle est constante.
|
|
// force indique la force surfacique appliquée
|
|
// retourne le second membre résultant
|
|
// nSurf : le numéro de la surface externe
|
|
// calcul à l'instant tdt ou t en fonction de la variable atdt
|
|
Vecteur& SM_charge_surf_Suiv_E (DdlElement & ddls,int nSurf
|
|
,const Tableau <Vecteur>& taphi,int nbne
|
|
,const Vecteur& poids,const Coordonnee& force,Fonction_nD* pt_fonct
|
|
,const ParaAlgoControle & pa,bool atdt=true);
|
|
// idem SM_charge_surf_Suiv_E mais -> implicite,
|
|
// pa: permet de déterminer si oui ou non on calcul la contribution à la raideur
|
|
// retourne le second membre et la matrice de raideur correspondant
|
|
Element::ResRaid SMR_charge_surf_Suiv_I (DdlElement & ddlA,int nSurf
|
|
,const Tableau <Vecteur>& taphi,int nbne
|
|
,const Vecteur& poids,const Coordonnee& force,Fonction_nD* pt_fonct
|
|
,const ParaAlgoControle & pa);
|
|
|
|
// cas d'un chargement volumique, sur l'élément
|
|
// force indique la force volumique appliquée
|
|
// retourne le second membre résultant
|
|
// calcul à l'instant tdt ou t en fonction de la variable atdt
|
|
Vecteur& SM_charge_vol_E (DdlElement & ddls
|
|
,const Tableau <Vecteur>& taphi,int nbne
|
|
,const Vecteur& poids,const Coordonnee& force,Fonction_nD* pt_fonct
|
|
,const ParaAlgoControle & pa,bool sur_volume_finale_,bool atdt=true);
|
|
|
|
// idem SM_charge_vol_E mais -> implicite,
|
|
// pa: permet de déterminer si oui ou non on calcul la contribution à la raideur
|
|
// retourne le second membre et la matrice de raideur correspondant
|
|
void SMR_charge_vol_I (DdlElement & ddls
|
|
,const Tableau <Vecteur>& taphi,int nbne
|
|
,const Vecteur& poids,const Coordonnee& force,Fonction_nD* pt_fonct
|
|
,const ParaAlgoControle & pa,bool sur_volume_finale_);
|
|
|
|
// cas d'un chargement hydrostatique, sur les surfaces de l'élément
|
|
// la charge dépend de la hauteur à la surface libre du liquide déterminée par un point
|
|
// et une direction normale à la surface libre:
|
|
// nSurf : le numéro de la surface externe
|
|
// poidvol: indique le poids volumique du liquide
|
|
// M_liquide : un point de la surface libre
|
|
// dir_normal_liquide : direction normale à la surface libre
|
|
// sans_limitation: indique si l'on doit limiter aux positions négatives
|
|
// retourne le second membre résultant
|
|
// calcul à l'instant tdt ou t en fonction de la variable atdt
|
|
Vecteur& SM_charge_hydro_E (DdlElement & ddls,int nSurf
|
|
,const Tableau <Vecteur>& taphi,int nbne
|
|
,const Vecteur& poids
|
|
,const Coordonnee& dir_normal_liquide,const double& poidvol
|
|
,const Coordonnee& M_liquide,bool sans_limitation
|
|
,const ParaAlgoControle & pa,bool atdt=true);
|
|
// idem SM_charge_hydro_E mais -> implicite,
|
|
// pa: permet de déterminer si oui ou non on calcul la contribution à la raideur
|
|
// retourne le second membre et la matrice de raideur correspondant
|
|
Element::ResRaid SMR_charge_hydro_I (DdlElement & ddlA,int nSurf
|
|
,const Tableau <Vecteur>& taphi,int nbne
|
|
,const Vecteur& poids
|
|
,const Coordonnee& dir_normal_liquide,const double& poidvol
|
|
,const Coordonnee& M_liquide,bool sans_limitation
|
|
,const ParaAlgoControle & pa);
|
|
|
|
// cas d'un chargement aero-hydrodynamique, sur les frontières de l'élément
|
|
// Il y a trois forces: une suivant la direction de la vitesse: de type traînée aerodynamique
|
|
// Fn = poids_volu * fn(V) * S * (normale*u) * u, u étant le vecteur directeur de V (donc unitaire)
|
|
// une suivant la direction normale à la vitesse de type portance
|
|
// Ft = poids_volu * ft(V) * S * (normale*u) * w, w unitaire, normal à V, et dans le plan n et V
|
|
// une suivant la vitesse tangente de type frottement visqueux
|
|
// T = to(Vt) * S * ut, Vt étant la vitesse tangentielle et ut étant le vecteur directeur de Vt
|
|
// retourne le second membre résultant
|
|
// calcul à l'instant tdt ou t en fonction de la variable atdt
|
|
// coef_mul: est un coefficient multiplicateur global (de tout)
|
|
Vecteur& SM_charge_hydrodyn_E (const double& poidvol,const Tableau <Vecteur>& taphi,int nbne
|
|
,Courbe1D* frot_fluid,const Vecteur& poids
|
|
,Courbe1D* coef_aero_n,int numfront,const double& coef_mul
|
|
,Courbe1D* coef_aero_t,const ParaAlgoControle & ,bool atdt=true);
|
|
// idem SM_charge_hydrodyn_E mais -> implicite,
|
|
// pa: permet de déterminer si oui ou non on calcul la contribution à la raideur
|
|
// retourne le second membre et la matrice de raideur correspondant
|
|
Element::ResRaid SM_charge_hydrodyn_I (const double& poidvol,const Tableau <Vecteur>& taphi,int nbne
|
|
,Courbe1D* frot_fluid,const Vecteur& poids,DdlElement & ddls
|
|
,Courbe1D* coef_aero_n,int numfront,const double& coef_mul
|
|
,Courbe1D* coef_aero_t,const ParaAlgoControle & pa);
|
|
|
|
// -------------------- calcul de frontières en protected -------------------
|
|
|
|
// --- fonction nécessaire pour la construction des Frontières linéiques ou surfaciques particulière à l'élément
|
|
// adressage des frontières linéiques et surfacique
|
|
// définit dans les classes dérivées, et utilisées pour la construction des frontières
|
|
virtual ElFrontiere* new_frontiere_lin(int num,Tableau <Noeud *> & tab, DdlElement& ddelem) = 0;
|
|
virtual ElFrontiere* new_frontiere_surf(int num,Tableau <Noeud *> & tab, DdlElement& ddelem) = 0;
|
|
|
|
// -------------------- stabilisation d'hourglass -------------------
|
|
// calcul d'élément de contrôle d'hourglass associée à un comportement
|
|
void Init_hourglass_comp(const ElemGeomC0& elgeHour, const string & str_precision
|
|
,LoiAbstraiteGeneral * loiHourglass,const BlocGen & bloc);
|
|
|
|
// stabilisation pour un calcul implicit
|
|
void Cal_implicit_hourglass();
|
|
|
|
// stabilisation pour un calcul explicit
|
|
void Cal_explicit_hourglass(bool atdt);
|
|
|
|
// -------------------- stabilisation transversale éventuelle de membrane ou de biel -------------------
|
|
// utilisée et alimentées par les classes dérivées: ElemMeca sert de conteneur ce qui permet
|
|
// d'éviter de redéfinir des conteneurs locaux aux éléments membranes et biel
|
|
class StabMembBiel // conteneur local pour un stockage globale de toutes les grandeurs concernées
|
|
{ public :
|
|
// par défaut
|
|
StabMembBiel();
|
|
// fonction du nombre de noeud ou pti
|
|
StabMembBiel(int nbnoe);
|
|
// fonction d'une valeur numérique et des noms éventuels de fonction
|
|
StabMembBiel(double v,string * sgamma,string* sponder);
|
|
// de copie
|
|
StabMembBiel(const StabMembBiel& a);
|
|
~StabMembBiel();
|
|
StabMembBiel& operator= (const StabMembBiel& a);
|
|
// --- acces aux données en inline
|
|
Enum_StabMembraneBiel& Type_stabMembrane() {return type_stabMembrane;}
|
|
bool& Aa_calculer() {return a_calculer;};
|
|
double& Valgamma() {return gamma;};
|
|
Fonction_nD* Pt_fct_gamma() {return pt_fct_gamma;}
|
|
void Change_pt_fct_gamma(Fonction_nD* fct) {pt_fct_gamma=fct;};
|
|
// affichage
|
|
const int& Affichage_stabilisation() const {return affichage_stabilisation;};
|
|
void Change_affichage_stabilisation(int niveau) {affichage_stabilisation = niveau;};
|
|
// gestion éventuelle des mini maxi
|
|
const bool& Gestion_maxi_mini() const {return gestion_maxi_mini;};
|
|
void Change_gestion_maxi_mini(bool gestion) {gestion_maxi_mini = gestion;};
|
|
//beta
|
|
const double& Beta() const {return beta;};
|
|
void Change_beta(double bet) {beta = bet;};
|
|
//f_mini
|
|
const double& F_mini() const {return f_mini;};
|
|
void Change_f_mini(double f) {f_mini = f;};
|
|
// d_maxi
|
|
const double& D_maxi() const {return d_maxi;};
|
|
void Change_d_maxi(double d) {d_maxi = d;};
|
|
// inactif_en_explicite
|
|
const bool Activite_en_explicite() const {return activite_en_explicite;};
|
|
void change_activite_en_explicite(bool change) {activite_en_explicite=change;};
|
|
// pondération de F_t, utilisé dans le calcul de F(t+dt)
|
|
Fonction_nD* Pt_fct_ponder_Ft() {return pt_fct_ponder_Ft;}
|
|
void Change_pt_fct_ponder_Ft(Fonction_nD* fct) {pt_fct_ponder_Ft=fct;};
|
|
|
|
// l'intensité de la force de stabilisation au noeud i ou au pti i
|
|
double& FF_StabMembBiel(int i) {return F_StabMembBiel(i);};
|
|
double& FF_StabMembBiel_t(int i) {return F_StabMembBiel_t(i);};
|
|
// l'énergie développée par la stabilisation au noeud i ou au pti i
|
|
double& EE_StabMembBiel(int i) {return E_StabMembBiel(i);};
|
|
|
|
// l'énergie développée par la stabilisation au noeud i ou au pti i
|
|
double& EE_StabMembBiel_t(int i) {return E_StabMembBiel_t(i);};
|
|
|
|
// énergie totale développée sur l'élément pour la stabilisation
|
|
double EE_total_StabMembBiel() const ;
|
|
double EE_total_StabMembBiel_t() const ;
|
|
|
|
// valeur de référence pour calculer l'intensité de stabilisation
|
|
// (identique pour tous les pti), typiquement == le max de la raideur par exemple
|
|
double& Stab_ref() {return stab_ref;};
|
|
|
|
string * Nom_fctnD() {return nom_fctnD;};
|
|
// actualisation de la force de stabilisation de t+dt vers t
|
|
void TdtversT();
|
|
// actualisation de la force de stabilisation de t vers tdt
|
|
void TversTdt();
|
|
|
|
// cas donne le niveau de la récupération
|
|
// = 1 : on récupère tout
|
|
// = 2 : on récupère uniquement les données variables (supposées comme telles)
|
|
void Lecture_bas_inf(ifstream& ent,const int cas) ;
|
|
// cas donne le niveau de sauvegarde
|
|
// = 1 : on sauvegarde tout
|
|
// = 2 : on sauvegarde uniquement les données variables (supposées comme telles)
|
|
void Ecriture_bas_inf(ofstream& sort,const int cas) ;
|
|
// complète en attribuant les fct nD si besoin, ceci après une lecture
|
|
// initiale base_info
|
|
void Complete_StabMembBiel(LesFonctions_nD* lesFonctionsnD);
|
|
|
|
// changement du nombre de pti ou de noeuds, suivant le type de stabilisation
|
|
void Change_nb_pti(int nbnoe);
|
|
// le nb de pti ou de noeud du conteneur
|
|
int Taille() const {return F_StabMembBiel.Taille();};
|
|
|
|
// --- données
|
|
protected:
|
|
Enum_StabMembraneBiel type_stabMembrane; // indique le type de stabilisation
|
|
bool a_calculer; // indique si oui ou non il faut calculer la force de stabilisation
|
|
double gamma; // si pt_fct_gamma == NULL -> valeur numérique lue
|
|
// si pt_fct_gamma != NULL -> contient la dernière valeur d'alpha calculée via pt_fct_gamma
|
|
Fonction_nD* pt_fct_gamma; // si non null, c'est cette fonction que l'on utilise sinon
|
|
// c'est la valeur numérique gamma
|
|
double stab_ref; // valeur de référence pour calculer l'intensité de stabilisation
|
|
// (identique pour tous les pti ou noeuds), typiquement == le max de la raideur par exemple
|
|
// peut changer ensuite pendant le calcul, suivant le type de stabilisation
|
|
|
|
bool activite_en_explicite; // par défaut false, indique si oui ou non actif en explicite
|
|
|
|
Fonction_nD* pt_fct_ponder_Ft; // si non nulle, la fonction est utilisée pour pondérer
|
|
// la force de stabilisation à t, utilisée dans le calcul
|
|
// à t+dt
|
|
|
|
int affichage_stabilisation; // niveau d'affichage spécifique
|
|
//--- gestion éventuelle des maxi mini -----
|
|
bool gestion_maxi_mini; // indique la gestion éventuelle (false par défaut)
|
|
// si l'intensité de la stabilisation est supérieure à beta * intensite_Fext
|
|
// on limite (cas où intensite_Fext >= F_mini)
|
|
double beta;
|
|
double f_mini; // force mini pour détection de force externe
|
|
double d_maxi; // limitation éventuelle sur un déplacement maxi
|
|
|
|
// --- stockage des infos -----
|
|
// stockage aux noeuds "ou" aux pti: dépend du type de stabilisation
|
|
Tableau <double> F_StabMembBiel; // la force généralisé de stabilisation actuelle
|
|
Tableau <double> F_StabMembBiel_t; // la force généralisé de stabilisation à t
|
|
Tableau <double> E_StabMembBiel; // l'énergie développée par la stabilisation
|
|
Tableau <double> E_StabMembBiel_t; // l'énergie développée par la stabilisation à t
|
|
|
|
string * nom_fctnD; // sert uniquement lors de la lecture base info, en attendant de définir
|
|
//pt_fct_gamma
|
|
string * nom_fctnD_2; // idem pour pt_fct_ponder_Ft
|
|
};
|
|
|
|
public : // retour de l'énergie totale de stabilisation éventuelle de membrane et ou de biel
|
|
double Energie_stab_membBiel()
|
|
{if (pt_StabMembBiel == NULL)
|
|
{return 0.;}
|
|
else
|
|
{return pt_StabMembBiel->EE_total_StabMembBiel();};
|
|
};
|
|
|
|
protected :
|
|
StabMembBiel* pt_StabMembBiel; // pointeur éventuellement non nul
|
|
double maxi_F_t; // grandeur de travail utiliser par Cal_implicit_StabMembBiel
|
|
Mat_pleine* matD; // raideur éventuelle, définie et supprimée dans la classe dérivée relative à un EF particulier
|
|
Vecteur* resD; // résidu éventuelle, défini et supprimé dans la classe dérivée relative à un EF particulier
|
|
// mise à jour de "a_calculer" en fonction du contexte
|
|
void Mise_a_jour_A_calculer_force_stab();
|
|
// stabilisation pour un calcul implicit,
|
|
// iteg -> donne le numéro du dernier pti sur lequel on a travaillé, auquel met est associé
|
|
// iteg == 0 : un seul calcul global, et on est à la suite d'une boucle
|
|
// correspond au calcul d'alpha: == stab_ref
|
|
// iteg == -1 : fin d'un calcul avec boucle sur tous les pti, et on est à la suite de la boucle
|
|
// iteg entre 1 et nbint: on est dans une boucle de pti
|
|
// nbint : le nombre maxi de pti
|
|
// poid_volume : si iteg > 0 : le poids d'intégration
|
|
// si iteg <= 0 : l'intégrale de l'élément (ici la surface totale)
|
|
// noeud_a_prendre_en_compte: si diff de null indique les noeuds à prendre en compte
|
|
// : noeud_a_prendre_en_compte(i) = 0 --> le noeud i n'est pas à prendre en compte
|
|
// doit avoir la taille du nombre total de noeud de l'élément
|
|
void Cal_implicit_StabMembBiel(int iteg,const Met_abstraite::Impli& met, const int& nbint
|
|
,const double& poid_volume,Tableau <int> * noeud_a_prendre_en_compte);
|
|
// stabilisation pour un calcul explicit
|
|
// iteg -> donne le numéro de pti sur lequel on travaille
|
|
void Cal_explicit_StabMembBiel(int iteg,const Met_abstraite::Expli_t_tdt met, const int& nbint
|
|
,const double& poids_volume,Tableau <int> * noeud_a_prendre_en_compte);
|
|
private :
|
|
// static Tableau< Vecteur> vec_StabMembBiel; // tableau de vecteurs de travail, éventuellement défini
|
|
// -------------------- fin stabilisation transversale éventuelle de membrane ou de biel -------------------
|
|
|
|
public :
|
|
// récupération de l'energie d'hourglass éventuelle
|
|
// uniquement pour un pas de temps
|
|
double Energie_Hourglass()const {return E_Hourglass;};
|
|
// idem pour l'énergie et la puissance de bulk viscosity
|
|
double Energie_Bulk() const {return E_elem_bulk_tdt;};
|
|
double Puissance_Bulk() const {return P_elem_bulk;};
|
|
// idem éventuellement pour l'intensité de la force de stabilisation transversale au pti i
|
|
double Force_StabMembBiel(int i) const
|
|
{if (pt_StabMembBiel != NULL) return pt_StabMembBiel->FF_StabMembBiel(i); else return 0.;};
|
|
// idem éventuellement pour l'énergie de stabilisation transversale
|
|
double Energie_StabMembBiel(int i) const
|
|
{if (pt_StabMembBiel != NULL) return pt_StabMembBiel->EE_StabMembBiel(i); else return 0.;};
|
|
// récupération des temps cpu
|
|
Temps_CPU_HZpp Temps_lois_comportement() const
|
|
{ Temps_CPU_HZpp inter; int taill = lesPtIntegMecaInterne->NbPti()+1;
|
|
for (int i=1;i<taill;i++) inter += (*lesPtIntegMecaInterne)(i).Tps_cpu_loi_comp_const();
|
|
return inter;};
|
|
Temps_CPU_HZpp Temps_metrique_K_SM()const
|
|
{ Temps_CPU_HZpp inter; int taill = lesPtIntegMecaInterne->NbPti()+1;
|
|
for (int i=1;i<taill;i++) inter += (*lesPtIntegMecaInterne)(i).TpsMetrique_const();
|
|
return inter;};
|
|
|
|
// modification de l'erreur relative
|
|
void Change_erreur_relative(double& nevez_erreur_rel)
|
|
{ if (sigErreur_relative == NULL) sigErreur_relative=new double;
|
|
*sigErreur_relative = nevez_erreur_rel;
|
|
};
|
|
// récup de l'erreur et l'erreur relative
|
|
const double Erreur_sigma() const
|
|
{ if (sigErreur == NULL) return 0.;
|
|
else return *sigErreur;
|
|
};
|
|
const double Erreur_sigma_relative() const
|
|
{ if (sigErreur_relative == NULL) return 0.;
|
|
else return *sigErreur_relative;
|
|
};
|
|
|
|
protected :
|
|
// -------------------- calcul d'erreur, remontée des contraintes -------------------
|
|
|
|
// calcul du résidu et de la matrice de raideur pour le calcul d'erreur
|
|
// cas d'une intégration suivant une seule liste
|
|
void SigmaAuNoeud_ResRaid(const int nbne,const Tableau <Vecteur>& taphi
|
|
,const Vecteur& poids,Tableau <Vecteur *>& resErr,Mat_pleine& raidErr
|
|
,const Tableau <Vecteur>& taphiEr,const Vecteur& poidsEr );
|
|
|
|
// calcul de l'erreur sur l'élément. Ce calcul n'est disponible
|
|
// qu'une fois la remontée aux contraintes effectuées sinon aucune
|
|
// action. pour les autres arguments de retour, idem ErreurElement
|
|
// qui est la fonction generique, les autres variables sont spécifiques
|
|
// a l'element.
|
|
void Cal_ErrElem(int type,double& errElemRelative
|
|
,double& numerateur, double& denominateur,const int nbne,const Tableau <Vecteur>& taphi
|
|
,const Vecteur& poids,const Tableau <Vecteur>& taphiEr,const Vecteur& poidsEr);
|
|
|
|
// calcul de l'erreur aux noeuds. Contrairement au cas des contraintes
|
|
// seul le résidu est calculé. Cas d'une intégration suivant une liste
|
|
void Cal_ErrAuxNoeuds(const int nbne, const Tableau <Vecteur>& taphi,
|
|
const Vecteur& poids,Tableau <Vecteur *>& resErr );
|
|
|
|
|
|
// ---------------------------------------------------------------------
|
|
// affichage dans la sortie transmise, des variables duales "nom"
|
|
// aux differents points d'integration
|
|
// dans le cas ou nom est vide, affichage de "toute" les variables
|
|
// cas =1 -> premier passage pour de l'implicit
|
|
// cas = 2 -> premier passage pour de l'explicit
|
|
// cas = 11 -> passage autre que le premier pour de l'implicit
|
|
// cas = 12 -> passage autre que le premier pour de l'explicit
|
|
void VarDualSort(ofstream& sort, Tableau<string>& nom,int nbint,int cas);
|
|
// utilitaires de VarDualSort
|
|
// affiche en fonction d'indic les differentes variables et appel
|
|
// AffDefContiD en fonction de la dimension i
|
|
//
|
|
void AffDefCont( ofstream& sort,Loi_comp_abstraite::SaveResul * saveDon,
|
|
TenseurBB& eps0BB,TenseurBB& epsBB,TenseurBB& DepsBB,TenseurBB& DeltaEpsBB,
|
|
TenseurHH& sigHH,
|
|
TenseurHB& epsHB,TenseurHB& sigHB,
|
|
Coordonnee& valPropreEps,Coordonnee& valPropreDeps,Coordonnee& valPropreSig,
|
|
double Mises,TenseurBB& epsAlmBB,TenseurBB& epslogBB, int indic);
|
|
// cas 1D
|
|
void AffDefCont1D( ofstream& sort,Loi_comp_abstraite::SaveResul * saveDon,
|
|
TenseurBB& eps0BB,TenseurBB& epsBB,TenseurBB& DepsBB,TenseurBB& DeltaEpsBB,
|
|
TenseurHH& sigHH,
|
|
TenseurHB& epsHB,TenseurHB& sigHB,
|
|
Coordonnee& valPropreEps,Coordonnee& valPropreDeps,Coordonnee& valPropreSig,
|
|
double Mises,TenseurBB& epsAlmBB,TenseurBB& epslogBB,int indic);
|
|
// cas 2D
|
|
void AffDefCont2D( ofstream& sort,Loi_comp_abstraite::SaveResul * saveDon,
|
|
TenseurBB& eps0BB,TenseurBB& epsBB,TenseurBB& DepsBB,TenseurBB& DeltaEpsBB,
|
|
TenseurHH& sigHH,
|
|
TenseurHB& epsHB,TenseurHB& sigHB,
|
|
Coordonnee& valPropreEps,Coordonnee& valPropreDeps,Coordonnee& valPropreSig,
|
|
double Mises,TenseurBB& epsAlmBB,TenseurBB& epslogBB,int indic);
|
|
// cas 3D
|
|
void AffDefCont3D( ofstream& sort,Loi_comp_abstraite::SaveResul * saveDon,
|
|
TenseurBB& eps0BB,TenseurBB& epsBB,TenseurBB& DepsBB,TenseurBB& DeltaEpsBB,
|
|
TenseurHH& sigHH,
|
|
TenseurHB& epsHB,TenseurHB& sigHB,
|
|
Coordonnee& valPropreEps,Coordonnee& valPropreDeps,Coordonnee& valPropreSig,
|
|
double Mises,TenseurBB& epsAlmBB,TenseurBB& epslogBB,int indic);
|
|
|
|
|
|
// ---------------------------------------------------------------------
|
|
|
|
// --- méthodes relatives aux calculs d'erreurs -----------
|
|
// ---- utilisées par les classes dérivées ---------
|
|
// ajout des ddl relatif aux contraintes pour les noeuds de l'élément
|
|
void Ad_ddl_Sigma(const DdlElement& tab_ddlErr);
|
|
// inactive les ddl du problème primaire de mécanique
|
|
void Inact_ddl_primaire(DdlElement& tab_ddl);
|
|
// active les ddl du problème primaire de mécanique
|
|
void Act_ddl_primaire(DdlElement& tab_ddl);
|
|
// inactive les ddl du problème de recherche d'erreur : les contraintes
|
|
void Inact_ddl_Sigma(DdlElement& tab_ddlErr);
|
|
// active les ddl du problème de recherche d'erreur : les contraintes
|
|
void Act_ddl_Sigma(DdlElement& tab_ddlErr);
|
|
// active le premier ddl du problème de recherche d'erreur : SIGMA11
|
|
void Act_premier_ddl_Sigma();
|
|
|
|
// ---------------------------------------------------------------------
|
|
// lecture des contraintes sur le flot d'entrée
|
|
void LectureDesContraintes
|
|
(bool cas,UtilLecture * entreePrinc,Tableau <TenseurHH *>& tabSigHH);
|
|
|
|
// retour des contraintes en absolu retour true si elle existe sinon false
|
|
void ContraintesEnAbsolues
|
|
(bool cas,Tableau <TenseurHH *>& tabSigHH,Tableau <Vecteur>& tabSig);
|
|
|
|
|
|
// ---------------- lecture écriture dans base info ----------------
|
|
// programmes utilisés par les classes dérivées
|
|
|
|
// cas donne le niveau de la récupération
|
|
// = 1 : on récupère tout
|
|
// = 2 : on récupère uniquement les données variables (supposées comme telles)
|
|
void Lecture_bas_inf
|
|
(ifstream& ent,const Tableau<Noeud *> * tabMaillageNoeud,const int cas) ;
|
|
// cas donne le niveau de sauvegarde
|
|
// = 1 : on sauvegarde tout
|
|
// = 2 : on sauvegarde uniquement les données variables (supposées comme telles)
|
|
void Ecriture_bas_inf(ofstream& sort,const int cas) ;
|
|
|
|
// --------- utilitaires pour la dynamique
|
|
// calcul de la longueur d'arrête de l'élément minimal
|
|
// divisé par la célérité la plus rapide dans le matériau
|
|
// appelé par les classes dérivées
|
|
// nb_noeud : =0 indique que l'on utilise tous les noeuds du tableau de noeuds
|
|
// = un nombre > 0, indique le nombre de noeuds à utiliser au début du tableau
|
|
double Interne_Long_arrete_mini_sur_c(Enum_dure temps,int nb_noeud=0);
|
|
|
|
// recuperation des coordonnées du point d'intégration numéro = iteg pour
|
|
// la grandeur enu
|
|
// temps: dit si c'est à 0 ou t ou tdt
|
|
// si erreur retourne erreur à true
|
|
Coordonnee CoordPtInt(Enum_dure temps,Enum_ddl enu,int iteg,bool& erreur);
|
|
|
|
// recuperation des coordonnées du point d'intégration numéro = iteg pour
|
|
// la face : face
|
|
// temps: dit si c'est à 0 ou t ou tdt
|
|
// si erreur retourne erreur à true
|
|
Coordonnee CoordPtIntFace(int face, Enum_dure temps,int iteg,bool& erreur);
|
|
// recuperation des coordonnées du point d'intégration numéro = iteg pour
|
|
// la face : face
|
|
// temps: dit si c'est à 0 ou t ou tdt
|
|
// si erreur retourne erreur à true
|
|
Coordonnee CoordPtIntArete(int arete, Enum_dure temps,int iteg,bool& erreur);
|
|
|
|
// procesure permettant de completer l'element, en ce qui concerne les variables gérés
|
|
// par ElemMeca, apres sa creation avec les donnees du bloc transmis
|
|
// peut etre appeler plusieurs fois
|
|
Element* Complete_ElemMeca(BlocGen & bloc,LesFonctions_nD* lesFonctionsnD);
|
|
|
|
// retourne le numero du pt d'ing le plus près ou est exprimé la grandeur enum
|
|
// temps: dit si c'est à 0 ou t ou tdt
|
|
int PtLePlusPres(Enum_dure temps,Enum_ddl enu, const Coordonnee& M);
|
|
|
|
// récupération des valeurs au numéro d'ordre = iteg pour
|
|
// les grandeur enu
|
|
// absolue: indique si oui ou non on sort les tenseurs dans la base absolue ou une base particulière
|
|
// NB Importante: il faut faire attention à ce que ces métriques soient identiques à celles qui ont servit
|
|
// pour le calcul des tenseurs: en particulier si c'est utilisé pour calculer les grandeurs pour le chargement
|
|
// il faut s'assurer que ce sont les "mêmes pti" qui servent pour la charge et pour la raideur !!
|
|
Tableau <double> Valeur_multi(bool absolue,Enum_dure enu_t,const List_io<Ddl_enum_etendu>& enu
|
|
,int iteg,int cas
|
|
) ;
|
|
|
|
// récupération des valeurs Tensorielles (et non scalaire comme avec Valeur_multi)
|
|
// au numéro d'ordre = iteg pour les grandeur enu
|
|
// enu contient les grandeurs de retour
|
|
// absolue: indique si oui ou non on sort les tenseurs dans la base absolue ou une base particulière
|
|
// NB Importante: il faut faire attention à ce que ces métriques soient identiques à celles qui ont servit
|
|
// pour le calcul des tenseurs: en particulier si c'est utilisé pour calculer les grandeurs pour le chargement
|
|
// il faut s'assurer que ce sont les "mêmes pti" qui servent pour la charge et pour la raideur !!
|
|
void Valeurs_Tensorielles(bool absolue, Enum_dure enu_t,List_io<TypeQuelconque>& enu
|
|
,int iteg,int cas
|
|
) ;
|
|
// récupération de valeurs interpolées pour les grandeur enu ou directement calculées
|
|
// absolue: indique si oui ou non on sort les tenseurs dans la base absolue ou une base particulière
|
|
// une seule des 3 métriques doit-être renseigné, les autres doivent être un pointeur nul
|
|
Tableau <double> Valeur_multi_interpoler_ou_calculer
|
|
(bool absolue, Enum_dure temps,const List_io<Ddl_enum_etendu>& enu
|
|
,Deformation & defor
|
|
,const Met_abstraite::Impli* ex_impli
|
|
,const Met_abstraite::Expli_t_tdt* ex_expli_tdt
|
|
,const Met_abstraite::Expli* ex_expli
|
|
);
|
|
// récupération de valeurs interpolées pour les grandeur enu
|
|
// absolue: indique si oui ou non on sort les tenseurs dans la base absolue ou une base particulière
|
|
// une seule des 3 métriques doit-être renseigné, les autres doivent être un pointeur nul
|
|
void Valeurs_Tensorielles_interpoler_ou_calculer
|
|
(bool absolue, Enum_dure temps,List_io<TypeQuelconque>& enu
|
|
,Deformation & defor
|
|
,const Met_abstraite::Impli* ex_impli
|
|
,const Met_abstraite::Expli_t_tdt* ex_expli_tdt
|
|
,const Met_abstraite::Expli* ex_expli
|
|
);
|
|
|
|
// ==== >>>> methodes virtuelles definis dans les classes mères ============
|
|
// ramene l'element geometrique correspondant au ddl passé en paramètre
|
|
virtual ElemGeomC0& ElementGeometrie(Enum_ddl ddl) const ;
|
|
// ramène le nombre de grandeurs génératrices pour un pt d'integ, correspondant à un type enuméré
|
|
// peut-être surchargé pour des éléments particuliers
|
|
virtual int NbGrandeurGene(Enum_ddl enu) const;
|
|
|
|
// ==== >>>> methodes virtuelles definis dans les classes dérivées ============
|
|
// ramene la dimension des tenseurs contraintes et déformations de l'élément
|
|
virtual int Dim_sig_eps() const = 0;
|
|
|
|
// VARIABLES PROTEGEES :
|
|
|
|
//---- variables: contrainte, déformations etc.. aux points d'intégrations -------
|
|
// lesPtIntegMecaInterne pointe vers une entitée entièrement gérés par les classes dérivées
|
|
// contenant les grandeurs mécaniques (contraintes, déformations etc.) stockées aux points d'intégration
|
|
LesPtIntegMecaInterne * lesPtIntegMecaInterne;
|
|
|
|
//----------- les charges externes éventuelles ----------
|
|
LesChargeExtSurElement * lesChargeExtSurEle;
|
|
|
|
|
|
Loi_comp_abstraite * loiComp; // loi de comportement mécanique defini dans les classes derivees
|
|
CompThermoPhysiqueAbstraite * loiTP; // éventuellement une loi de comportement thermo physique
|
|
// defini dans les classes derivees
|
|
CompFrotAbstraite * loiFrot; // éventuellement une loi de comportement au frottement pour ses frontières
|
|
|
|
int dilatation; // indique si oui ou non on tiend compte de la dilatation thermique
|
|
Met_abstraite * met; // definition specifique dans les classes derivees
|
|
Deformation * def; // definition specifique dans les classes derivees
|
|
// mais relative à la déformation mécanique
|
|
Deformation * defEr; // idem que def mais pour la remonte aux contraintes
|
|
Tableau <Deformation *> defSurf; // idem mais pour les déformations des surfaces
|
|
// externes (frontières) si cela est nécessaire
|
|
Tableau <Deformation *> defArete; // idem mais pour les déformations des arretes
|
|
// externes (frontières) si cela est nécessaire
|
|
Deformation * defMas; // idem que def mais pour le calcul de la matrice masse
|
|
double* sigErreur; // erreur sur l'élément dans le calcul des contraintes;
|
|
double* sigErreur_relative; // idem en valeur relative
|
|
Tableau <Loi_comp_abstraite::SaveResul *> tabSaveDon; // donnée particulière à la loi mécanique
|
|
Tableau <CompThermoPhysiqueAbstraite::SaveResul *> tabSaveTP; // donnée parti à la loi thermo physique
|
|
Tableau <Deformation::SaveDefResul *> tabSaveDefDon; // donnée particulière pour la déf
|
|
Tableau <EnergieMeca > tab_energ,tab_energ_t; //les différentes énergies mécaniques mises en jeux
|
|
// dimensionné dans les classes dérivées, a t+dt, et a t
|
|
EnergieMeca energie_totale_t,energie_totale; // le bilan sur l'élément de l'énergie
|
|
bool premier_calcul_meca_impli_expli; // au premier calcul soit en implicite, explicite, mat geom
|
|
// toutes les grandeurs à t=0 de la métrique sont calculées et stockées dans la déformation
|
|
// ensuite elles ne sont plus calculées -> premier_calcul_meca_impli_expli sert pour l'indiquer
|
|
double masse_volumique; // masse volumique de l'élément
|
|
//--- stabilisation d'hourglass éventuelle
|
|
Enum_StabHourglass type_stabHourglass; // méthode de stabilisation
|
|
double E_Hourglass;
|
|
|
|
private:
|
|
//--- stabilisation d'hourglass éventuelle
|
|
Tableau <ElemMeca*> tab_elHourglass; // tableau de pointeurs éventuels sur des éléments servant à la stabilisation
|
|
// les éléments sont définies au travers de l'appel à la méthode Cal_hourglass_comp, par une classe dérivée,
|
|
double coefStabHourglass ; // coef de stabilisation
|
|
Mat_pleine* raid_hourglass_transitoire; // raideur transitoire d'hourglass, éventuellement définie
|
|
static Tableau< Vecteur> vec_hourglass; // tableau de vecteurs de travail, éventuellement défini
|
|
|
|
|
|
// -------- pour faciliter les routines Interne_t _0 et _tdt + paramètres de contrôle ----
|
|
static const Coordonnee & (Met_abstraite::*PointM)
|
|
(const Tableau<Noeud *>& tab_noeud,const Vecteur& phi);
|
|
static void (Met_abstraite::*BaseND)
|
|
(const Tableau<Noeud *>& tab_noeud,
|
|
const Mat_pleine& dphi,const Vecteur& phi,BaseB& bB,BaseH& bH);
|
|
|
|
// ---------- cas du bulk viscosity -------------
|
|
private: // pour éviter les modifications par les classes dérivées
|
|
static int bulk_viscosity; // indique si oui ou non on utilise le bulk viscosity
|
|
// en plus = 1: fonctionnement normal, =2 fonctionnement quelque soit
|
|
// le signe de IDeps
|
|
static double c_traceBulk; // coeff de la trace de D pour le bulk
|
|
static double c_carreBulk; // coeff du carré de la trace de D pour le bulk
|
|
static TenseurHH * sig_bulkHH; // variable de travail pour le bulk
|
|
// ---------- pour le calcul de la masse pour la méthode de relaxation dynamique
|
|
static Ddl_enum_etendu masse_relax_dyn; // simplement pour le définir une fois pour toute
|
|
|
|
protected:
|
|
// pour une utilisation par les classes dérivées
|
|
inline bool Bulk_visco_actif() {return bulk_viscosity;};
|
|
inline double& C_traceBulk() {return c_traceBulk;};
|
|
inline double& C_carreBulk() {return c_carreBulk;};
|
|
// calcul de la contrainte modifiée en fonction du bulk viscosity
|
|
// ramène l'énergie et la puissance par unité de volume, dépensée par le bulk
|
|
DeuxDoubles ModifContrainteBulk(Enum_dure temps,const TenseurHH& gijHH_tdt,TenseurHH & sigHH,const TenseurBB & DepsBB_tdt
|
|
,TenseurHH & sigbulkHH);
|
|
// actualisation des ddl et des grandeurs actives de t+dt vers t
|
|
void TdtversT_();
|
|
// actualisation des ddl et des grandeurs actives de t vers tdt
|
|
void TversTdt_();
|
|
// calcul des intégrales de volume et de volume + temps
|
|
// cas = 1 : pour un appel après un calcul implicit
|
|
// cas = 2 : pour un appel après un calcul explicit
|
|
void Calcul_Integ_vol_et_temps(int cas,const double& poid_jacobien,int iteg);
|
|
// au premier pti init éventuel des intégrales de volume et temps
|
|
void Init_Integ_vol_et_temps();
|
|
|
|
// ---- utilitaires interne
|
|
private:
|
|
// calcul du mini du module d'young équivalent pour tous les points d'intégrations
|
|
double Calcul_mini_E_qui(Enum_dure temps);
|
|
// calcul du maxi du module d'young équivalent pour tous les points d'intégrations
|
|
double Calcul_maxi_E_qui(Enum_dure temps);
|
|
// définition d'un repère d'anisotropie à un point d'intégration
|
|
BaseH DefRepereAnisotropie(int nb_pti,LesFonctions_nD* lesFonctionsnD,const BlocGen & bloc);
|
|
|
|
};
|
|
/// @} // end of group
|
|
|
|
|
|
|
|
#endif
|