next_inactive up previous


T1phvmn

Tutorial 3 :
Calcul thermomécanique
exploitant la symétrie de la pièce &
Utilisation d'un maillage pour le calcul thermique
différent de celui du calcul mécanique
Version 1


Date: 21 février 2007


Table des matières

Table des évolutions

Version Date Modifications Auteur(s)
01 28/01/2007 Création PC
       
       
       

Mots clefs : thermomécanique, projection, maillage tétra, maillage héxa, MECA_STATIQUE.

Avant propos

Dans l'esprit des exercices présents sur CAELINUX, le propos est ici de guider les débutants dans l'utilisation de Code Aster $ ^{\tiny {\textregistered }}$ aussi bien dans la mise en donnée que dans la démarche. Nombre d'explications apparaîtrons triviales rendant le style « lourd », nous espérons cependant qu'elles seront claires.

Ce document se voulant interactif, le signalement de toute erreur est le bienvenu $ \ldots$ il en va de même pour les améliorations.

Le lecteur est invité à lire le tutorial 1 (CAELINUX_post-traitement.pdf) traitant spécifiquement du post-traitement (graphique) des résultats ; il y trouvera quelques façons de visualiser les résultats (déplacements, champs de contraintes, courbes, etc. ...). De cette façon, les différents tutoriaux se concentreront sur le thème de la note.

Introduction

Dans le présent tutorial, nous réaliserons un calcul thermomécanique en exploitant la symétrie cyclique de la pièce.

L'autre particularité du document réside dans l'utilisation de 2 maillages notablement différents :

  1. un premier maillage à base de tétraèdres linéaires pour le calcul thermique,
  2. un second maillage constitué d'hexaèdres quadratiques pour le calcul mécanique1.

Les calculs ont été réalisés sous Code Aster $ ^{\tiny {\textregistered }}$ V8.3.

Préparation du modèle numérique/définition du maillage

La figure 1 présente la maquette 3D de la pièce d'étude. Cette dernière présente une symétrie cyclique à 60$ ^{\circ}$ (Schéma élémentaire : $ \frac{1}{6}$ de la pièce).

Figure 1: Maquette 3D

Maillages

Maillage pour le calcul thermique

Nous avons fait le choix d'utiliser un maillage différent pour le calcul thermique pour les raisons suivantes :
$ \triangleright$
Les éléments thermiques n'ont qu'un seul point d'intégration ; dès lors il est inutile d'utiliser des éléments du second ordre,
$ \triangleright$
La cartographie des températures fait qu'il n'est pas toujours nécessaire d'avoir un maillage fin (si peu ou pas de gradient thermique) alors que cela peut être indispensable pour le calcul mécanique,
$ \triangleright$
Enfin, certaines parties de la pièce présentes dans le calcul thermique sont supprimées dans le calcul mécanique $ \ldots$
$ \triangleright$
$ \ldots$ sur ce dernier point, imaginons que la pièce fait partie d'une canalisation d'eau soumise à un fort gradient thermique dans le plan radial ; il est évident que le champ des températures est conditionné par la présence du fluide $ \Rightarrow$ il est maillé pour le calcul thermique et il sera supprimé lors du calcul mécanique.)

La figure 2 présente le maillage "thermique".

Figure 2: Maillage tétra utilisé pour le calcul thermique
Figure 3: Maillage tétra utilisé pour le calcul thermique
Image maillage_tetra_vue1

Image maillage_hexa_vue1

Maillage pour le calcul mécanique

Classiquement, nous privilégions un maillage structuré du second ordre pour le calcul mécanique (cf. figure 3).

Déroulement des calculs

Je pense qu'il est plus pertinent de lancer les simulations en plusieurs temps :
$ \triangleright$
Le calcul thermique à proprement parlé (THER_LINEAIRE),
$ \triangleright$
Le post-traitement des champs thermiques,
$ \triangleright$
Le calcul mécanique à proprement parlé (MECA_STATIQUE),
$ \triangleright$
Le post-traitement des grandeurs mécaniques,

Il est en effet dommage de relancer tout un calcul (particulièrement si celui-ci est long) pour une problématique de post-traitement (recherche d'un instant particulier, calcul d'une grandeur supplémentaire $ \ldots$ ou tout simplement une erreur de post-traitement).

Etape 1 : calcul thermique

F mess /symetrie_cyclique_TH.mess R 6
F erre /symetrie_cyclique_TH.erre R 9
R base /base_TH R 0
F unv  /symetrie_cyclique_tetra.unv D 19
F unv  /symetrie_cyclique_hexa.unv D 21
F comm /symetrie_cyclique_TH.comm D 1

Etape 2 : Post-traitement des résultats thermiques - sous Gibi ici

F mess /symetrie_cyclique_TH_post.mess R 6
F erre /symetrie_cyclique_TH_post.erre R 9
R base /base_TH D 0
F comm /symetrie_cyclique_TH_post.comm D 1
F cast /symetrie_cyclique_TH.cast R 50

Etape 3 : calcul mécanique

F mess /symetrie_cyclique_MEC.mess R 6
F erre /symetrie_cyclique_MEC.erre R 9
R base /base_TH D 0
R base /base_MECA R 0
F comm /symetrie_cyclique_MEC.comm D 1

Etape 4 : Post-traitement des résultats mécaniques - sous Gibi ici

F comm /symetrie_cyclique_MEC_post.comm D 1
R base /base_MECA D 0
F mess /symetrie_cyclique_MEC_post.mess R 6
F resu /symetrie_cyclique_MEC_post.resu R 8
F erre /symetrie_cyclique_MEC_post.erre R 9
F cast /symetrie_cyclique_MEC_post.cast R 50

définition et affectation du modèle

Nous lisons ici des maillages au format IDEAS $ ^{\tiny {\textregistered }}$(.unv). Les éléments choisis sont des éléments de volume 3D pour la thermique et 3D sous-intégrés (3D_SI) pour la mécanique2.

PRE_IDEAS(UNITE_IDEAS=19,
          UNITE_MAILLAGE=20,);

MAIL_TH=LIRE_MAILLAGE(UNITE=20,);
#Lecture du maillage mécanique

PRE_IDEAS(UNITE_IDEAS=21,
          UNITE_MAILLAGE=22,);

MAIL_MEC=LIRE_MAILLAGE(UNITE=22,);


#Definition & affectation des modèles

MOD_TH=AFFE_MODELE(MAILLAGE=MAIL_TH,
                   AFFE=_F(TOUT='OUI',
                           PHENOMENE='THERMIQUE',
                           MODELISATION='3D',),);

MOD_MEC=AFFE_MODELE(MAILLAGE=MAIL_MEC,
                    AFFE=_F(TOUT='OUI',
                            PHENOMENE='MECANIQUE',
                            MODELISATION='3D_SI',),);

Calcul thermique

Définition du matériau

Nous définissons les caractéristiques thermiques du matériau, à savoir la conductivité thermique $ \lambda$=14.9 mW/mm-$ ^{\circ}$C et la chaleur volumique $ \rho C_{p}$=3.76 mJ/mm$ ^{3}$.

BIDON_TH=DEFI_MATERIAU(THER=_F(LAMBDA=14.9,
                               RHO_CP=3.76,),);

MAT_TH=AFFE_MATERIAU(MAILLAGE=MAIL_TH,
                     AFFE=_F(TOUT='OUI',
                             MATER=BIDON_TH,),);

Définition des conditions aux limites

Nous avons choisi d'imposer une température de 200$ ^\circ$C sur la face intérieure et des conditions de convection sur la face extérieure (coefficient d'échange h=0.5 mW/mm$ ^{2}$ et température extérieure de 20$ ^\circ$C).

A l'instar du calcul mécanique, il est possible d'exploiter la symétrie cyclique de la pièce en imposant des conditions adiabatiques aux 2 plans de symétrie ; cela se traduit par un flux imposé nul ($ \Phi=0$).

Rappelons que les températures s'appliquent aux noeuds (le groupe de mailles T_IMPO est transformé en groupe de noeuds avec DEFI_GROUP), les flux et les échanges convectifs sur les faces (cf. figure 4).

Figure 4: Groupes de mailles utilisés pour les CL thermiques
Image GROUP_MA_CL_TH

MAIL_TH=DEFI_GROUP(reuse =MAIL_TH,
                   MAILLAGE=MAIL_TH,
                   CREA_GROUP_NO=_F(GROUP_MA='T_IMPO',),);

CL_TH=AFFE_CHAR_THER(MODELE=MOD_TH,
                     TEMP_IMPO=_F(GROUP_NO='T_IMPO',
                                  TEMP=200,),
                     FLUX_REP=_F(GROUP_MA='ADIAB',
                                 FLUN=0.0,),
                     ECHANGE=_F(GROUP_MA='CONV',
                                COEF_H=0.5,
                                TEMP_EXT=20,),);

Solveur thermique

Nous réalisons un calcul thermique linéaire en condition stationnaire.
RESOL_TH=THER_LINEAIRE(MODELE=MOD_TH,
                       CHAM_MATER=MAT_TH,
                       EXCIT=_F(CHARGE=CL_TH,),
                       TEMP_INIT=_F(STATIONNAIRE='OUI',),);

Post-traitement (sous GIBI $ ^{\tiny {\textregistered }}$)

Le fichier de post-traitement sous Code Aster $ ^{\tiny {\textregistered }}$ est de la forme :
POURSUITE();
RESOL_TH=CALC_ELEM(reuse =RESOL_TH,
                   RESULTAT=RESOL_TH,
                   OPTION='FLUX_ELNO_TEMP',);

IMPR_RESU(FORMAT='CASTEM',
          UNITE=50,
          RESU=_F(MAILLAGE=MAIL_TH,
                  RESULTAT=RESOL_TH,
                  NOM_CHAM=('TEMP','FLUX_ELNO_TEMP',),),);
FIN();

A noter que le maillage DOIT ABSOLUMENT être imprimé dans le fichier résultat (cf. MAILLAGE=MAIL_TH,).

Les figures 5 et 6 présentent le résultat graphique.

Figure: Cartographie des températures visualisée sous GIBI $ ^{\tiny {\textregistered }}$ - vue 1
Figure: Cartographie des températures visualisée sous GIBI $ ^{\tiny {\textregistered }}$ - vue 2
Image champ_thermique_gibi_vue1

Image champ_thermique_gibi_vue2

Calcul mécanique

Chainage thermique $ \rightarrow $ mécanique

Les 2 maillages étant différents, il est nécessaire d'indiquer au maillage "mécanique" les températures nodales ; pour ce faire l'instruction PROJ_CHAMP projette le champ thermique (valeurs nodales) de l'un vers l'autre.
PROJ_TH=PROJ_CHAMP(RESULTAT=RESOL_TH,
                   MODELE_1=MOD_TH,
                   MODELE_2=MOD_MEC,
                   VIS_A_VIS=_F(TOUT_1='OUI',
                                TOUT_2='OUI',
                                CAS_FIGURE='3D',),);

Le reste de la mise en données est assez classique.

Définition du matériau

Nous avons choisi comme matériau un acier inoxydable X2 CrNi 18-10 (AISI 304L) dont les propriétés varient avec la température. Pour ce faire, des fonctions sont créées :

Température Module d'YOUNG E Coefficient de Coefficient de POISSON
    dilatation thermique  
[$ ^{\circ}$C] [MPa] [/$ ^{\circ}$C] [adimensionel]
20 193000 ref. 0.3
100 186000 $ 16.4\cdot10^{-6}$ 0.3
200 178000 $ 17.0\cdot10^{-6}$ 0.3
300 170000 $ 17.5\cdot10^{-6}$ 0.3

Ces propriétés sont affectées à l'ensemble du maillage.

ALPHA_T=DEFI_FONCTION(NOM_PARA='TEMP',
                      VALE=
                      (100,16.4E-6,200,17E-6,300,17.5E-6,),
                      PROL_DROITE='CONSTANT',
                      PROL_GAUCHE='CONSTANT',);

YOUNG_T=DEFI_FONCTION(NOM_PARA='TEMP',
                      VALE=
                      (20,193000,100,186000,200,178000,300,170000,),
                      PROL_DROITE='CONSTANT',
                      PROL_GAUCHE='CONSTANT',);

POISSON=DEFI_CONSTANTE(VALE=0.3,);

SS304L=DEFI_MATERIAU(ELAS_FO=_F(E=YOUNG_T,
                                NU=POISSON,
                                RHO=8.47e-009,
                                TEMP_DEF_ALPHA=20.0,
                                ALPHA=ALPHA_T,),);
                                
MATERIAU=AFFE_MATERIAU(MAILLAGE=MAIL_MEC,
                       AFFE=_F(TOUT='OUI',
                               MATER=SS304L,
                               TEMP_REF=20.0,),);

Conditions aux limites et chargement

Les figures suivantes montrent les groupes de mailles sur lesquelles les conditions aux limites et le chargement en pression seront appliqués (figure 7), de même que les repères global (X$ _{1}$) et local (X$ _{2}$) (figure 8).

La définition du repère local se fait au travers de LIAISON_OBLIQUE dans AFFE_CHAR_MECA (nous retrouvons les 60$ ^{\circ}$ dans la définition de l'angle nautique3).

Nous avons choisi d'utiliser ici des groupes de mailles (que nous transformons en groupe de noeuds avec CREA_GROUP_NO dans DEFI_GROUP).

Figure 7: Vue des groupes de mailles servant aux CL
Figure 8: Conditions aux limites - Repères global et local
Image groupes_ma_meca

Image reperes

Contrairement aux éléments de coque ou de poutre (6 DDL), les éléments volumiques sont définis par 3 DDL (les 3 translations - comme pour les membranes et les barres) qu'il convient de bloquer sous peine d'avoir des pivots nuls (signe que le système est hypostatique et donc insoluble).

Le système reste isostatique en bloquant les noeuds des groupes (cf. 7) :

$ \triangleright$
PLAN_YOZ selon DX (repère global),
$ \triangleright$
ENCASTRG selon DY (repère global),
$ \triangleright$
PLAN_12 selon DX (repère local),
$ \triangleright$
ENCAS_LO selon DY (repère local),

NOTA :

  1. les groupes de mailles ENCASTRE et PLAN_12 possèdent des noeuds communs alors qu'ils sont définis l'un dans le repère local, l'autre dans le repère global ; bloquer des translations aux 2 sans précaution conduit à une surabondance des CL4,
  2. Bloquer DX dans le repère local équivaut à bloquer DX & DZ dans le repère global $ \ldots$ sachant que DY$ _{local}$ équivaut à DY$ _{global}$, les 3 translations sont biens bloquées.

Suite aux remarques précédentes, des opérations booléennes sont opérées de façon à retirer les noeuds communs aux 2 groupes de ENCASTRE selon le schéma :

$ \triangleright$
ENCAS_LO = PLAN_12 intersection ENCASTRE,
$ \triangleright$
ENCASTRG = ENCASTRE moins ENCAS_LO,
$ \triangleright$
il est nécessaire de créer de nouveaux groupes $ \ldots$ ici ENCAS_LO pour ENCAS_LO(cal) et ENCASTRG pour ENCASTRG(lobal) respectivement définis dans le repère local et global.

Ce qui donne :

MAIL=DEFI_GROUP(reuse =MAIL,
                MAILLAGE=MAIL,
                CREA_GROUP_NO=(_F(GROUP_MA='ENCASTRE',),
                               _F(GROUP_MA='PLAN_12',),
                               _F(GROUP_MA='PLAN_YOZ',),),);

MAIL=DEFI_GROUP(reuse =MAIL,
                MAILLAGE=MAIL,
                CREA_GROUP_NO=(_F(INTERSEC=('ENCASTRE','PLAN_12',),
                                  NOM='ENCAS_LO',),
                               _F(DIFFE=('ENCASTRE','ENCAS_LO',),
                                  NOM='ENCASTRG',),),);

Le blocage des CL devient alors :

CL=AFFE_CHAR_MECA(MODELE=MOD_MEC,
                  DDL_IMPO=(_F(GROUP_NO='ENCASTRG',
                               DY=0,),
                            _F(GROUP_NO='PLAN_YOZ',
                               DX=0,),),
                  LIAISON_OBLIQUE=(_F(GROUP_NO='PLAN_12',
                                      ANGL_NAUT=(0,-60.0,0,),
                                      DX=0,),
                                   _F(GROUP_NO='ENCAS_LO',
                                      ANGL_NAUT=(0,-60,0,),
                                      DY=0,),),);

CHARGMNT=AFFE_CHAR_MECA(MODELE=MOD_MEC,
                        TEMP_CALCULEE=PROJ_TH,
                        PESANTEUR=(0,-1,0,981000,),
                        PRES_REP=_F(GROUP_MA='PRESSION',
                                    PRES=172.81,),);

NOTA : Pour calculer les dilatations de la pièce, il est nécessaire d'indiquer à Code Aster $ ^{\tiny {\textregistered }}$ l'existence d'un champ thermique grâce à TEMP_CALCULEE dans AFFE_CHAR_MECA.

Pour le chargement, nous cherchons à modéliser un pseudo-serrage5 sous la forme d'une pression (correspondant à une vis ChC de M6 serrée à 10 NM). Nous supposerons que la pièce se dilate librement sous la tête de vis. Pour le fun, nous supposons une gravité de 100G.

Rappelons que toute pression s'exerce en sens inverse de la normale (sortante) des faces (cette dernière est définie par la topologie des éléments) ; il est préférable de demander à Code Aster $ ^{\tiny {\textregistered }}$ de réorienter les normales des groupes de mailles utilisés (au pire il ne réoriente aucune maille ; au mieux il évite une erreur) :

MAIL=MODI_MAILLAGE(reuse =MAIL,
                   MAILLAGE=MAIL,
                   ORIE_PEAU_3D=_F(GROUP_MA='PRESSION',),
                   MODELE=MODELE,);


CHARGMNT=AFFE_CHAR_MECA(MODELE=MODELE,
                        PESANTEUR=(0,-1,0,981000,),
                        PRES_REP=_F(GROUP_MA='PRESSION',
                                    PRES=172.81,),);

Type de calcul

Nous faisons ici un classique calcul statique linéaire et nous pointons sur le modèle (MODELE), les caractéristiques du matériau (CHAM_MATER) et les chargements/CL (EXIT).

RESOL_ME=MECA_STATIQUE(MODELE=MOD_MEC,
                       CHAM_MATER=MATERIAU,
                       EXCIT=(_F(CHARGE=CHARGMNT,),
                              _F(CHARGE=CL,),),);

Fichier de post-traitement de Code Aster $ ^{\tiny {\textregistered }}$

La mise en données sous Code ASTER $ ^{\tiny {\textregistered }}$ est du type :
POURSUITE();

RESOL=CALC_ELEM(reuse =RESOL,
                RESULTAT=RESOL,
                OPTION=('SIGM_ELNO_DEPL','EQUI_ELNO_SIGM',),);

IMPR_RESU(FORMAT='CASTEM',
          UNITE=51,
          RESU=_F(MAILLAGE=MAIL,
                  RESULTAT=RESOL,
                  NOM_CHAM=('SIGM_ELNO_DEPL','EQUI_ELNO_SIGM','DEPL',),),);

FIN();

A noter que le maillage DOIT ABSOLUMENT être imprimé dans le fichier résultat (cf. MAILLAGE=MAIL_MEC,).

Les figures 9 et 10 présentent les isovaleurs de déplacements ; les 2 suivantes (fig 12 & 13) sont intéressantes parce qu'elles visualisent les vecteurs déplacements.

Figure: Isovaleurs des déplacements sous GIBI $ ^{\tiny {\textregistered }}$ - vue 1
Figure: Isovaleurs des déplacements sous GIBI $ ^{\tiny {\textregistered }}$ - vue 2
Image champ_depl_gibi_vue1

Image champ_depl_gibi_vue2

Une fonctionnalité intéressante : il est possible de superposer le maillage déformé sur le maillage initial (cf. figure 11).

Figure: Maillage déformé VS maillage initial sous GIBI $ ^{\tiny {\textregistered }}$ - vue 1
Image deformee_gibi

Figure: Vecteurs déplacements sous GIBI $ ^{\tiny {\textregistered }}$ - vue 1
Figure: Vecteurs déplacements sous GIBI $ ^{\tiny {\textregistered }}$ - vue 2
Image champ_depl_gibi_fleches_vue1

Image champ_depl_gibi_fleches_vue2

Les cartographies des contraintes sont visibles aux figures 14 et 15.

Figure 14: Visualisation des $ \sigma _{eq}$ de Von Mises sous GIBI $ ^{\tiny {\textregistered }}$ - vue 1
Figure 15: Visualisation des $ \sigma _{eq}$ de Von Mises sous GIBI $ ^{\tiny {\textregistered }}$ - vue 2
Image vmis_gibi_vue1

Image vmis_gibi_vue2

De même, il est possible de visualiser les champs de contraintes sur le maillage déformée (cf. figure 16.

Figure 16: Visualisation des $ \sigma _{eq}$ de Von Mises sur le maillage déformé (GIBI $ ^{\tiny {\textregistered }}$)
Image vmis_gibi_vue_deformee

Fichiers Code Aster $ ^{\tiny {\textregistered }}$ complets

Calcul thermique

DEBUT();
#Lecture des maillages : 
#Lecture maillage thermique

PRE_IDEAS(UNITE_IDEAS=19,
          UNITE_MAILLAGE=20,);

MAIL_TH=LIRE_MAILLAGE(UNITE=20,);
#Lecture du maillage mécanique

PRE_IDEAS(UNITE_IDEAS=21,
          UNITE_MAILLAGE=22,);

MAIL_MEC=LIRE_MAILLAGE(UNITE=22,);


#Definition & affectation des modèles

MOD_TH=AFFE_MODELE(MAILLAGE=MAIL_TH,
                   AFFE=_F(TOUT='OUI',
                           PHENOMENE='THERMIQUE',
                           MODELISATION='3D',),);

MOD_MEC=AFFE_MODELE(MAILLAGE=MAIL_MEC,
                    AFFE=_F(TOUT='OUI',
                            PHENOMENE='MECANIQUE',
                            MODELISATION='3D_SI',),);

BIDON_TH=DEFI_MATERIAU(THER=_F(LAMBDA=14.9,
                               RHO_CP=3.76,),);

MAT_TH=AFFE_MATERIAU(MAILLAGE=MAIL_TH,
                     AFFE=_F(TOUT='OUI',
                             MATER=BIDON_TH,),);
#Défintion du chargement thermique :
#- GROUP_MA : ADIAB => flux nul (symétrie)
#- GROUP_MA : CONV => convection forcée
#- GROUP_MA : T_IMPO => T imposée

MAIL_TH=DEFI_GROUP(reuse =MAIL_TH,
                   MAILLAGE=MAIL_TH,
                   CREA_GROUP_NO=_F(GROUP_MA='T_IMPO',),);

CL_TH=AFFE_CHAR_THER(MODELE=MOD_TH,
                     TEMP_IMPO=_F(GROUP_NO='T_IMPO',
                                  TEMP=200,),
                     FLUX_REP=_F(GROUP_MA='ADIAB',
                                 FLUN=0.0,),
                     ECHANGE=_F(GROUP_MA='CONV',
                                COEF_H=0.5,
                                TEMP_EXT=20,),);
#Solveur thermique

RESOL_TH=THER_LINEAIRE(MODELE=MOD_TH,
                       CHAM_MATER=MAT_TH,
                       EXCIT=_F(CHARGE=CL_TH,),
                       TEMP_INIT=_F(STATIONNAIRE='OUI',),);

FIN();

post-traitement thermique

POURSUITE();

RESOL_TH=CALC_ELEM(reuse =RESOL_TH,
                   RESULTAT=RESOL_TH,
                   OPTION='FLUX_ELNO_TEMP',);

IMPR_RESU(FORMAT='CASTEM',
          UNITE=50,
          RESU=_F(MAILLAGE=MAIL_TH,
                  RESULTAT=RESOL_TH,
                  NOM_CHAM=('TEMP','FLUX_ELNO_TEMP',),),);


RESOL_TH=CALC_NO(reuse =RESOL_TH,
                 RESULTAT=RESOL_TH,
                 OPTION='FLUX_NOEU_TEMP',);

IMPR_RESU(FORMAT='GMSH',
          UNITE=51,
          RESU=_F(RESULTAT=RESOL_TH,
                  NOM_CHAM='TEMP',),);

FIN();

Calcul thermique

POURSUITE();
#Projection du champ thermique du maillage thermique (tétra linéaire) 
# sur le maillage mécanique (hexa quadratique)

PROJ_TH=PROJ_CHAMP(RESULTAT=RESOL_TH,
                   MODELE_1=MOD_TH,
                   MODELE_2=MOD_MEC,
                   VIS_A_VIS=_F(TOUT_1='OUI',
                                TOUT_2='OUI',
                                CAS_FIGURE='3D',),);
#Préparation du modèle mécanique :
#- réorientation des normales pour l'application d'une pression

MAIL_MEC=MODI_MAILLAGE(reuse =MAIL_MEC,
                       MAILLAGE=MAIL_MEC,
                       ORIE_PEAU_3D=_F(GROUP_MA='PRESSION',),
                       MODELE=MOD_MEC,);
#- Création des groupes de noeuds pour les CL

MAIL_MEC=DEFI_GROUP(reuse =MAIL_MEC,
                    MAILLAGE=MAIL_MEC,
                    CREA_GROUP_NO=(_F(GROUP_MA='ENCASTRE',),
                                   _F(GROUP_MA='PLAN_12',),
                                   _F(GROUP_MA='PLAN_YOZ',),),);

MAIL_MEC=DEFI_GROUP(reuse =MAIL_MEC,
                    MAILLAGE=MAIL_MEC,
                    CREA_GROUP_NO=(_F(INTERSEC=('ENCASTRE','PLAN_12',),
                                      NOM='ENCAS_LO',),
                                   _F(DIFFE=('ENCASTRE','ENCAS_LO',),
                                      NOM='ENCASTRG',),),);
#Propriétés des caractéristiques mécaniques du matériaux

ALPHA_T=DEFI_FONCTION(NOM_PARA='TEMP',
                      VALE=
                      (100,16.4E-6,200,17E-6,300,17.5E-6,),
                      PROL_DROITE='CONSTANT',
                      PROL_GAUCHE='CONSTANT',);

YOUNG_T=DEFI_FONCTION(NOM_PARA='TEMP',
                      VALE=
                      (20,193000,100,186000,200,178000,300,170000,),
                      PROL_DROITE='CONSTANT',
                      PROL_GAUCHE='CONSTANT',);

POISSON=DEFI_CONSTANTE(VALE=0.3,);

SS304L=DEFI_MATERIAU(ELAS_FO=_F(E=YOUNG_T,
                                NU=POISSON,
                                RHO=8.47e-009,
                                TEMP_DEF_ALPHA=20.0,
                                ALPHA=ALPHA_T,),);

MATERIAU=AFFE_MATERIAU(MAILLAGE=MAIL_MEC,
                       AFFE=_F(TOUT='OUI',
                               MATER=SS304L,
                               TEMP_REF=20.0,),);

#Définition des conditions aux limites :
#- LIAISON_OBLIQUE => Chargement de repère


CL=AFFE_CHAR_MECA(MODELE=MOD_MEC,
                  TEMP_CALCULEE=PROJ_TH,
                  DDL_IMPO=(_F(GROUP_NO='ENCASTRG',
                               DY=0,),
                            _F(GROUP_NO='PLAN_YOZ',
                               DX=0,),),
                  LIAISON_OBLIQUE=(_F(GROUP_NO='PLAN_12',
                                      ANGL_NAUT=(0,-60.0,0,),
                                      DX=0,),
                                   _F(GROUP_NO='ENCAS_LO',
                                      ANGL_NAUT=(0,-60,0,),
                                      DY=0,),),);

CHARGMNT=AFFE_CHAR_MECA(MODELE=MOD_MEC,
                        PESANTEUR=(0,-1,0,981000,),
                        PRES_REP=_F(GROUP_MA='PRESSION',
                                    PRES=172.81,),);

RESOL_ME=MECA_STATIQUE(MODELE=MOD_MEC,
                       CHAM_MATER=MATERIAU,
                       EXCIT=(_F(CHARGE=CHARGMNT,),
                              _F(CHARGE=CL,),),);

FIN();

post-traitement mécanique

POURSUITE();

RESOL_ME=CALC_ELEM(reuse =RESOL_ME,
                   RESULTAT=RESOL_ME,
                   OPTION=('SIGM_ELNO_DEPL','EQUI_ELNO_SIGM',),);

RESOL_ME=CALC_NO(reuse =RESOL_ME,
                 RESULTAT=RESOL_ME,
                 OPTION=('SIGM_NOEU_DEPL','EQUI_NOEU_SIGM',),);

IMPR_RESU(FORMAT='CASTEM',
          UNITE=50,
          RESU=_F(MAILLAGE=MAIL_MEC,
                  RESULTAT=RESOL_ME,
                  NOM_CHAM=('DEPL','SIGM_ELNO_DEPL','EQUI_ELNO_SIGM',),),);

IMPR_RESU(FORMAT='GMSH',
          UNITE=51,
          RESU=_F(RESULTAT=RESOL_ME,
                  NOM_CHAM='EQUI_NOEU_SIGM',
                  NOM_CMP='VMIS',),);

IMPR_RESU(FORMAT='GMSH',
          UNITE=52,
          RESU=_F(RESULTAT=RESOL_ME,
                  NOM_CHAM='DEPL',),);

Conclusion, remerciements, auteur(s)

Droit d'auteur(s) : L'utilisation de ce document sous quelque forme que ce soit est absolument libre au sens que la licence GPL donne à ce terme. Nous souhaitons simplement, si de larges extraits de cette publication sont utilisés dans d'autres documents, qu'il soit fait mention du nom du(des) auteur(s) du document initial.

Remerciements : L'auteur(s) souhaite(nt) remercier les personnes de la communauté de Code Aster $ ^{\tiny {\textregistered }}$ et de Cast3M $ ^{\tiny {\textregistered }}$ qui ont fourni aides et informations précieuses : qu'ils en soient chaleureusement remercier. L'ensemble des notes présentes sur le site ne sont que la continuation de cet esprit "libre" qui vise entre-autres choses partager son savoir et son travail $ \ldots$

L'auteur : Paul CARRICO6.

Les commentaires sont à adresser à :


paul.carrico_at_free.fr

Correspondance dans les unités

Unités METRIQUE METRIQUE ANGLO-SAXONNE ANGLO-SAXONNE
  MKS mmNS FPS IPS
Longueur m mm ft in
Temps sec sec sec sec
Masse Kg tonne slug lbf-sec$ ^{2}$
Force N N lbf lbf
Température $ ^{\circ}$C $ ^{\circ}$C $ ^{\circ}$F $ ^{\circ}$F
Aire m$ ^{2}$ mm$ ^{2}$ ft$ ^{2}$ in$ ^{2}$
Volume m$ ^{3}$ mm$ ^{3}$ ft$ ^{3}$ (cu-ft) in$ ^{3}$ (cu-in)
Vitesse m/sec mm/sec ft/sec in/sec
Accélération m/sec$ ^{2}$ mm/sec$ ^{2}$ ft/sec$ ^{2}$ in/sec$ ^{2}$
Angle, rotation rad rad rad rad
Vitesse angulaire rad/sec$ ^{2}$ rad/sec$ ^{2}$ rad/sec$ ^{2}$ rad/sec$ ^{2}$
Masse volumique Kg/m$ ^{3}$ Tonne/mm$ ^{3}$ slug/ft$ ^{3}$ lbf-sec$ ^{2}$/in$ ^{4}$
Moment, couple N-m N-mm ft-lbf in-lbf
Force linéïque N/m N/mm lbf/ft lbf/in
Force répartie sur une surface N/m$ ^{2}$ N/mm$ ^{2}$ lbf/ft$ ^{2}$ lbf/in$ ^{2}$
(Contrainte, pression,Module d'Young) (Pa) (MPa)   (Psi)
Coefficient de Dilatation thermique /$ ^{\circ}$C /$ ^{\circ}$C /$ ^{\circ}$F /$ ^{\circ}$F
  (/K) (/K) (/K) (/K)
Moment Quadratique d'une poutre $ IG_{z}$ m$ ^{4}$ mm$ ^{4}$ ft$ ^{4}$ in$ ^{4}$
Moment d'inertie transverse d'une poutre Kg-m$ ^{2}$ tonne-mm$ ^{2}$ slug-ft$ ^{2}$ lbf-in-sec$ ^{2}$
Energie, Travail, Chaleur J mJ ft-lbf in-lbf
Puissance, taux de transfert thermique W mW ft-lbf/sec in-lbf/sec
Gradient de température $ ^{\circ}$C/m $ ^{\circ}$C/mm $ ^{\circ}$F/ft $ ^{\circ}$F/in
Flux thermique W/m$ ^{2}$ mW/mm$ ^{2}$ lbf/ft-sec lbf/in-sec
Conductivité thermique W/m-$ ^{\circ}$C mW/mm-$ ^{\circ}$C lbf/sec-$ ^{\circ}$F lbf/sec-$ ^{\circ}$F
Chaleur spécifique C$ _{p}$ J/Kg-$ ^{\circ}$C mJ/tonne-$ ^{\circ}$C ft-lbf/slug-$ ^{\circ}$F in$ ^{2}$/sec$ ^{2}$-$ ^{\circ}$F

Rappels :

$ \triangleright$
1 W $ \equiv$ 1 N-m/sec,
$ \triangleright$
1 mJ $ \equiv$ 1 N-mm,
$ \triangleright$
1 mW $ \equiv$ 1 N-mm/sec,
$ \triangleright$
1 N/m$ ^{2}$ $ \equiv$ 1 Pa,
$ \triangleright$
1 N/mm$ ^{2}$ $ \equiv$ 1 MPa,

À propos de ce document...

Tutorial 3 :
Calcul thermomécanique
exploitant la symétrie de la pièce &
Utilisation d'un maillage pour le calcul thermique
différent de celui du calcul mécanique
Version 1

This document was generated using the LaTeX2HTML translator Version 2002-2-1 (1.70)

Copyright © 1993, 1994, 1995, 1996, Nikos Drakos, Computer Based Learning Unit, University of Leeds.
Copyright © 1997, 1998, 1999, Ross Moore, Mathematics Department, Macquarie University, Sydney.

The command line arguments were:
latex2html -split 0 CAELINUX_thermomecanique.tex

The translation was initiated by on 2007-02-21


Notes

... mécanique1
Il s'agit du maillage utilisé au cours du tutorial 2
... mécanique2
Ces éléments donnent les mêmes résultats que les éléments 3D tout en présentant des gains de temps calculs appréciables $ \ldots$ particulièrement en dynamique ; ils ne sont cependant pas utilisables avec des éléments tétra ou penta.
... nautique3
cf. U4.42.01.pdf - mot clef ORIENTATION pour la définition et le signe de l'angle
... CL4
ce que Code Aster $ ^{\tiny {\textregistered }}$ annoncera avec une erreur.
... pseudo-serrage5
Nous ne cherchons pas à réaliser un serrage dans les règles de l'art et nous supposons ici que la surface à bloquer sous la pièce est la même que celle sous la tête de vis $ \ldots$ c'est-à-dire un diamètre $ \phi$=10 mm.
... CARRICO6
Avant de travailler sur des outils propriétaires, j'ai fait mes premières armes sur Code Aster $ ^{\textregistered}$ que j'utilise encore ponctuellement. Ces tutoriaux sont pour moi l'occasion de coucher sur le papier une démarche (pas nécessairement exhaustive) qui est la mienne à la fois pour la partager et l'améliorer/corriger au grès des retours, et l'occasion de garder la main sur un code que considère comme rigoureux (par opposition aux boîtes noires qui fleurissent depuis un certain nombre d'années).

next_inactive up previous
2007-02-21