Metafor

ULiege - Aerospace & Mechanical Engineering

User Tools

Site Tools


doc:user:meshtransfer:datatransferbetweenmeshes

Transfert de données entre deux maillages

Introduction

Le transfert de données est utilisé dans le cadre du remaillage (le nouveau maillage n'a aucun lien avec l'ancien maillage). Pour des raisons extérieures aux algorithmes de transfert, le transfert de données est effectué entre deux objets Metafor par l'intermédiaire de l'objet MetaforTransferOperator. Ainsi, les données (les champs de température, de vitesse, de pression, de contrainte, etc.) sont transférées d'un ancien Metafor dans lequel un calcul a été effectué sur l'ancien maillage vers un nouveau Metafor dans lequel un nouveau maillage a été défini.

Procédure de transfert

La procédure de transfert est donc la suivante.

  • À la fin d'un calcul effectué sur un objet Metafor, on souhaite effectuer un remaillage.
  • Un nouvel objet Metafor est créé avec un nouveau maillage.
  • L'état d'équilibre du système calculé au dernier pas de temps avec l'ancien Metafor est transféré vers le nouveau Metafor.
  • Le calcul peut être poursuivi avec le nouveau Metafor.
Cette procédure de transfert de données est utilisée dans le formalisme X-ALE-FEM.

Cette partie est focalisée sur les commandes permettant de transférer les données calculées avec l'ancien Metafor (l'ancien maillage) vers le nouveau Metafor (le nouveau maillage). On suppose que le calcul avec l'ancien Metafor a été correctement exécuté et que le nouveau Metafor a été créé.

La procédure est identique pour un problème 2D ou 3D.

Données transférées

De quelles données s'agit-il ? Ces données sont fonction de plusieurs facteurs. Elles dépendent:

  • du matériau : il faut par exemple transférer les contraintes si le matériau est mécanique, les déformations plastiques si le matériau est élastoplastique, etc.
  • du type d'élément : on doit transférer la pression séparément des grandeurs déviatoriques lorsqu'on utilise un élément sous-intégré (type SRI).
  • du schéma d'intégration : si le schéma est dynamique, il faut transférer la densité, les vitesses et éventuellement les accélérations si le schéma est implicite.

Elles sont donc de deux types: des données “aux points d'intégration” (par exemple les contraintes) et les données “nodales” (par exemple la température).
Tous ces choix sont automatiques et l'utilisateur ne doit pas s'en préoccuper.

Commandes python nécessaires

On suppose que oldMetafor réfère à l'ancien Metafor et que newMetafor réfère au nouveau Metafor.
On souhaite transférer les données de l'interaction oldInteraction de oldMetafor vers l'interaction newInteraction de newMetafor. Pour utiliser le calcul des transfert de donnée il faut importer wrap.mtDataTransfer.

from wrap.mtDataTransfer import *
transferOperator = MetaforTransferOperator(oldMetafor, newMetafor)
transferOperator.setCouplingInteraction(oldInteraction, newInteraction)
transferOperator.execute()

Par défaut, la méthode de transfert de données utilisée est celle appelée Finite Volume Transfer Method (FVTM) avec la reconstruction linéaire des champs sur l'ancien maillage auxiliaire (LR) et 5 points de Gauss par direction sont utilisés dans l'approximation numérique du calcul du couplage entre les éléments. Cette méthode est recommandée dans la conclusion de l'article suivant pour être utilisée par défaut.

P. Bussetta, R. Boman and J.-P. Ponthot. Efficient 3D data transfer operators based on numerical integration. International Journal for Numerical Methods in Engineering 2014; XX:xxx-xxx, doi: 10.1002/nme.4821

Commandes python optionnelles

Les commandes optionnelles sont à placer avant la commande execute() de l'objet transferOperator

L'objet region est utilisé pour gérer le transfert de donnée entre deux Interaction (de oldInteraction vers newInteraction). Il est créé automatiquement par l'objet MetaforTransferOperator. Pour pouvoir le récupérer, il faut utiliser l'objet oldInteraction.

region = transferOperator.get(oldInteraction)
region.add(transferElementProperties)

Ces commandes permettent de choisir la méthode de transfert de données utilisée par l'intermédiaire de l'objet transferElementProperties (voir la partie suivante: “Methode de transfert de données”). Tous les champs à transférer le seront avec cette méthode de transfert.
Il est également possible de spécifier une méthode de transfert pour un champ particulier:

region.add(IF_DEV_SIG_XX, transElemPropSpec)

Le transfert des contraintes et de toutes les autres grandeurs définies sur les mêmes points d'intégration utilisera les propriétés de transElemPropSpec (objet du même type que transferElementProperties).

region.ignore(IF_FTOTAL)

Cette commande permet de ne pas transférer le tenseur F, ceci aura pour conséquence que toutes les mesures de déformations seront fausses en sortie. Si on n'a pas besoin de ces grandeurs, on économise 9 transferts! Cette commande est dangereuse puisqu'il est ainsi possible de désactiver par exemple le transfert des contraintes (ce qui, dans tous les cas, entrainera une solution erronée).

Méthodes de transfert de données (optionnelle)

Ces méthodes de transferts de données sont utilisées et testées dans le répertoire : apps/remeshing
Ces méthodes de transferts peuvent aussi être utilisées dans le cadre du formalisme Arbitraire Eulérien Lagrangien (ALE) : voir dans le répertoire apps/ale

Les trois premières méthodes de transfert sont décrites en détail dans l'article suivant:

P. Bussetta, R. Boman and J.-P. Ponthot. Efficient 3D data transfer operators based on numerical integration. International Journal for Numerical Methods in Engineering 2014; XX:xxx-xxx, doi: 10.1002/nme.4821

Méthode de transfert par élément (interpolation/extrapolation)

La méthode de transfert par éléments est la plus utilisée, mais aussi la plus intuitive. Le champ sur le nouveau maillage est calculé directement aux points caractéristiques (noeuds ou points d’intégration) par interpolation /extrapolation des valeurs de l’élément de l’ancien maillage dont il fait partie. Ainsi, pour un champ défini par des valeurs nodales, la valeur de chacun des noeuds du nouveau maillage est calculée par interpolation des valeurs nodales de l’élément de l’ancien maillage le contenant. De même, pour un champ défini par des valeurs aux points d'intégration, la valeur de chaque point d’intégration du nouveau maillage est calculée par interpolation/extrapolation des valeurs des points d’intégration de l’élément le contenant.

La définition de l' ElementProperties permettant d'utiliser cette méthode est la suivante:

transferElementProperties = ElementProperties(ETMCell)

Méthode de transfert utilisant les éléments joints (Mortar Elements)

La méthode de transfert utilisant les éléments joints est une méthode de transfert sous forme faible. Le couplage entre les deux maillages est effectué grâce à ces éléments joints. Ainsi, les valeurs d’un champ sur le nouveau maillage sont calculées directement aux points caractéristiques à partir des valeurs des éléments joints et des valeurs aux points caractéristiques de l’ancien maillage. Cette méthode rend l'opération de projection implicite. Ainsi, avec cette méthode, deux possibilités sont proposées pour le calcul du transfert :

  • un calcul local (explicite) ou
  • un calcul global (implicite).

Pour chacune de ces possibilités, deux méthodes de calcul des éléments joints ont été proposées :

  • un calcul approché (plus rapide) et
  • un calcul exact (plus couteux en temps de calcul, car nécessitant la construction d’un super-maillage).

La définition de l' ElementProperties permettant d'utiliser cette méthode avec le calcul approché est :

transferElementProperties = ElementProperties(MTMCell)
transferElementProperties.put (MTMSOLVERTYPE, solver  )
transferElementProperties.put (INTPT_NB     , intePtNb)
transferElementProperties.put (QUADRATURETYPE, QuadType)
  • solver = type de solveur:
    • GLOBALSOLVER : calcul global (implicite)
    • LOCALSOLVER : calcul local (explicite)
  • intePtNb = Nombre de points d'intégration pour le calcul du couplage
    • nombre de points par direction pour les quadrangles, hexaèdres, etc.
    • nombre de points d'intégration total pour les triangles, tétraèdres, etc.
  • QuadType = schéma d'intégration utilisé dans le calcul approché du couplage:
    • GAUSSQUADRATURE : méthodes de quadrature de Gauss (valeur par défaut)
    • TRAPEZOIDQUADRATURE : méthode des trapèzes
    • LOBATTOQUADRATURE : méthodes de quadrature de Lobatto

Pour utiliser le calcul exact (nécessite la construction d'un supermaillage) il faut remplacer MTMCell par ExactMTMCell et importer wrap.mtExactDataTransfer_CGAL. L'ElementProperties est donc défini par :

from wrap.mtExactDataTransfer_CGAL import *
transferElementProperties = ElementProperties(ExactMTMCell)
transferElementProperties.put (MTMSOLVERTYPE, solver  )

Méthode de transfert utilisant des volumes finis

La méthode de transfert utilisant des volumes finis est une méthode de transfert sous forme faible. Cette méthode de transfert utilise une reconstruction des champs sur des volumes finis. Avec cette méthode, les champs sont conservés de manière faible, mais le transfert est calculé entre des champs construits sur des volumes finis localisés sur les points caractéristiques (noeuds ou points d’intégration). Le champ reconstruit sur ces volumes finis peut être constant ou défini par une fonction linéaire sur chaque volume fini. Avec cette méthode, le calcul du transfert est explicite. Deux méthodes de calcul du couplage sont proposées :

  • un calcul approché (plus rapide) et
  • un calcul exact (plus couteux en temps de calcul, car nécessitant la construction d’un super-maillage).

La définition de l' ElementProperties permettant d'utiliser cette méthode avec le calcul approché est la suivante:

transferElementProperties = ElementProperties(FVTMCell)
transferElementProperties.put (INTPT_NB     , intePtNb)
transferElementProperties.put (QUADRATURETYPE, QuadType)
  • intePtNb = Nombre de points d'intégration pour le calcul du couplage
    • nombre de points par direction pour les quadrangles, hexaèdres, etc.
    • nombre de points d'intégration total pour un triangle, tétraèdre, etc.
  • QuadType = schéma d'intégration utilisé dans le calcul approché du couplage:
    • GAUSSQUADRATURE : méthodes de quadrature de Gauss (valeur par défaut)
    • TRAPEZOIDQUADRATURE : méthode des trapèzes
    • LOBATTOQUADRATURE : méthodes de quadrature de Lobatto

Pour utiliser le calcul exact (nécessite la construction d'un supermaillage) il faut remplacer FVTMCell par ExactFVTMCell et importer wrap.mtExactDataTransfer_CGAL. L'ElementProperties est donc défini par :

from wrap.mtExactDataTransfer_CGAL import *
transferElementProperties = ElementProperties(ExactFVTMCell)

Option pour la reconstruction constante

transferElementProperties.put (FVCELLTYPE  , GODUNOVCELL)

Option pour la reconstruction linéaire

transferElementProperties.put (FVCELLTYPE  , LINEARRECCELL)
transferElementProperties.put (STENCILTYPE, stencil )
transferElementProperties.put (LIMITERTYPE, limiter )

Cette méthode plus précise permet de capter de manière plus correcte d'éventuelles discontinuités (gradients élevés) et montre beaucoup moins de diffusion. Elle n'est cependant pas TVD si on n'utilise pas un limiteur de flux.

  • stencil = type de stencil:
    • DUMMY_STENCIL : le gradient est nul
    • LEASTSQUARE_STENCIL : calcul du gradient par moindres carrés
    • GREENGAUSS_STENCIL : calcul du gradient par la formule de Green-Gauss
  • limiter = type de limiteur:
    • NO_LIMITER : le gradient n'est pas limité
    • SIMPLE_LIMITER : limite le gradient (la valeur sur les frontières est inférieure à la valeur principale des voisins)

Méthode de transfert par voisinage

Cette méthode de transfert de données est utilisée uniquement dans le cadre d'un remaillage partiel pour transférer les données vers le nouveau Metafor dans les zones non remaillées. Avec cette méthode la valeur en point caractéristique (noeud ou point d’intégration) du nouveau maillage est égale à celle du point caractéristique le plus proche de l'ancien maillage.

La définition de l' ElementProperties permettant d'utiliser cette méthode est la suivante:

transferElementProperties = ElementProperties(NeighbourTMCell)
transferElementProperties.put (FVCELLTYPE, GODUNOVCELL)

Exemple

On souhaite transférer les données de oldInteraction1 vers newInteraction1 et de oldInteraction2 vers newInteraction2

transferElementProperties = ElementProperties(FVTMCell)
transferElementProperties.put (FVCELLTYPE , LINEARRECCELL)
transferElementProperties.put (STENCILTYPE, LEASTSQUARE_STENCIL)
transferElementProperties.put (LIMITERTYPE, SIMPLE_LIMITER )
transferElementProperties.put (INTPT_NB     , 5)
transferOperator = MetaforTransferOperator(oldMetafor, newMetafor)
transferOperator.setCouplingInteraction(oldInteraction1, newInteraction1)
transferOperator.setCouplingInteraction(oldInteraction2, newInteraction2)
region1 = transferOperator.get(oldInteraction1) # optionnelle
region1.add(transferElementProperties) # optionnelle
region2 = transferOperator.get(oldInteraction2) # optionnelle
region2.add(transferElementProperties) # optionnelle
transferOperator.execute()

Extracteur de valeurs

TransferValueExtractor

Cet extracteur permet de calculer une valeur de l'erreur commise durant le transfert de données

region = transferOperator.get(oldInteraction)
TransferValueExtractor(region, Field, type)

where

field Field1D or internal field
oldInteraction est une Interaction de l'ancien Metafor.
  • type
    • 'BEFORE' : calcul de l'intégrale sur la region de la valeur du champ avant remaillage
    • 'AFTER' : calcul de l'intégrale sur la region de la valeur du champ après remaillage
    • 'DIFF' : calcul de la différence entre l'intégrale sur la region de la valeur du champ avant et après remaillage en %
    • 'TRANSFER_ERROR' : calcul de la somme des erreurs de transfert de chacun des éléments de la region (erreur = (valeur avant - valeur après)/valeur avant)

Example:

TransferValueExtractor(region, IF_EPL, 'TRANSFER_ERROR')

TransferDiffValueExtractor

Cet extracteur permet de calculer la différence après transfert entre la valeur d'un champ et une fonction référence fctMultiChamps.

region = transferOperator.get(oldInteraction)
TransferDiffValueExtractor(region, Field, fctMultiChamps , fieldList)
field Field1D or internal field
oldInteraction est une Interaction de l'ancien Metafor.
fctMultiChamps Functions y=f(t) (spatiotemporal dependency)
fieldList list of dependency fields

Example:

TransferDiffValueExtractor(region, IF_EPL, fct, FieldList(Field1D(TX,AB), Field1D(TY,AB)))
doc/user/meshtransfer/datatransferbetweenmeshes.txt · Last modified: 2016/03/30 15:23 by 127.0.0.1

Donate Powered by PHP Valid HTML5 Valid CSS Driven by DokuWiki