![]()  | 
  
    Mmg
    
   Simplicial remeshers (mesh adaptation, isovalue discretization, lagrangian movement) 
   | 
 
#include "mmgcommon_private.h"#include <vtkSmartPointer.h>#include <vtkXMLReader.h>#include <vtkXMLWriter.h>#include <vtkXMLUnstructuredGridReader.h>#include <vtkXMLUnstructuredGridWriter.h>#include <vtkXMLPUnstructuredGridWriter.h>#include <vtkXMLPolyDataReader.h>#include <vtkXMLPolyDataWriter.h>#include <vtkXMLPPolyDataWriter.h>#include <vtkDataSetReader.h>#include <vtkDataSetWriter.h>#include <vtkPDataSetWriter.h>#include <vtkDataSet.h>#include <vtkUnstructuredGrid.h>#include <vtkPolyData.h>#include <vtkStructuredGrid.h>#include <vtkPointData.h>#include <vtkCellData.h>#include <vtkFieldData.h>#include <vtkCellTypes.h>#include <vtkDataArray.h>#include <vtkFloatArray.h>#include <vtkDoubleArray.h>#include <vtkLine.h>#include <vtkTriangle.h>#include <vtkQuad.h>#include <vtkTetra.h>#include <vtkWedge.h>#include <vtkCellArray.h>#include <typeinfo>

Go to the source code of this file.
Functions | |
| int | MMG5_loadVtpMesh_part1 (MMG5_pMesh, const char *, vtkDataSet **, int8_t *, int8_t *, int *, int8_t *, int8_t *) | 
| int | MMG5_loadVtuMesh_part1 (MMG5_pMesh, const char *, vtkDataSet **, int8_t *, int8_t *, int *, int8_t *, int8_t *) | 
| int | MMG5_loadVtkMesh_part1 (MMG5_pMesh, const char *, vtkDataSet **, int8_t *, int8_t *, int *, int8_t *, int8_t *) | 
| int | MMG5_loadVtkMesh_part2 (MMG5_pMesh, MMG5_pSol *, vtkDataSet **, int8_t, int8_t, int, int8_t, int8_t) | 
| static void | MMG5_internal_VTKSetLine (vtkSmartPointer< vtkPolyData > d, vtkSmartPointer< vtkCellArray > *ca) | 
| static void | MMG5_internal_VTKSetLine (vtkSmartPointer< vtkUnstructuredGrid > d, vtkSmartPointer< vtkCellArray > *ca) | 
| static void | MMG5_internal_VTKSetCells (int *t, vtkSmartPointer< vtkPolyData > d, vtkSmartPointer< vtkCellArray > ca) | 
| static void | MMG5_internal_VTKSetCells (int *t, vtkSmartPointer< vtkUnstructuredGrid > d, vtkSmartPointer< vtkCellArray > ca) | 
| static void | MMG5_internal_VTKbinary (vtkXMLUnstructuredGridWriter *w, int binary) | 
| static void | MMG5_internal_VTKbinary (vtkXMLPolyDataWriter *w, int binary) | 
| static void | MMG5_internal_VTKbinary (vtkDataSetWriter *w, int binary) | 
| template<class T , class TWriter , class PWriter > | |
| int | MMG5_saveVtkMesh_i (MMG5_pMesh mesh, MMG5_pSol *sol, const char *mfilename, int metricData, int binary, int npart, int myid, int master) | 
| template<class T , class TWriter , class PWriter > | |
| int | MMG5_saveVtkMesh (MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename, int metricData, int binary) | 
      
  | 
  static | 
| w | vtk writer | 
| binary | 1 if we want to save in binary format | 
Try to set the suitable file format to the vtk writer
Definition at line 148 of file vtkparser.hpp.
      
  | 
  static | 
| w | vtk writer | 
| binary | 1 if we want to save in binary format | 
Try to set the suitable file format to the vtk writer
Definition at line 134 of file vtkparser.hpp.
      
  | 
  static | 
| w | vtk writer | 
| binary | 1 if we want to save in binary format | 
Try to set the suitable file format to the vtk writer
Definition at line 120 of file vtkparser.hpp.
      
  | 
  static | 
| t | array of integer storing the cell types | 
| d | vtk data type in which we want to store the array ca | 
| ca | vtk cell array containing the cells connectivity | 
Store a list of vtk cells into a vtkPolyData
Definition at line 101 of file vtkparser.hpp.

      
  | 
  static | 
| t | array of integer storing the cell types | 
| d | vtk data type in which we want to store the array ca | 
| ca | vtk cell array containing the cells connectivity | 
Store a list of vtk cells into a vtkUnstructuredGrid
Definition at line 111 of file vtkparser.hpp.
      
  | 
  static | 
| d | vtk data type in which we want to store the array ca | 
| ca | vtk cell array containing the lines connectivity | 
Store a list of vtk cells containing only lines into a vtkPolyData and reset the list of cells. Otherwise, lines are casted into polygons by the SetPoly method.
Definition at line 76 of file vtkparser.hpp.

      
  | 
  static | 
| d | vtk data type in which we want to store the array ca | 
| ca | vtk cell array containing the lines connectivity | 
Nothing to do here, the SetCell method allow to handle both lines and elements.
Definition at line 90 of file vtkparser.hpp.
| int MMG5_loadVtkMesh_part1 | ( | MMG5_pMesh | mesh, | 
| const char * | filename, | ||
| vtkDataSet ** | dataset, | ||
| int8_t * | ptMeditRef, | ||
| int8_t * | eltMeditRef, | ||
| int * | nsols, | ||
| int8_t * | metricData, | ||
| int8_t * | lsData | ||
| ) | 
| mesh | pointer toward a MMG5 mesh | 
| filename | name of the input file. | 
| dataset | vtkdataset structure | 
| ptMeditRef | index of point data field that contains references (field named medit:ref), -1 if no references | 
| eltMeditRef | index of a cell data field that contains references (field named medit:ref), -1 if no references | 
| nsols | number of point data (except the medit:ref ones) | 
| metricData | 1 if file contains a metric data highlighted by the :metric name | 
| lsData | 1 if file contains a metric data highlighted by the :ls name | 
I/O at Vtk VTK file format.
Definition at line 285 of file vtkparser.cpp.


| int MMG5_loadVtkMesh_part2 | ( | MMG5_pMesh | mesh, | 
| MMG5_pSol * | sol, | ||
| vtkDataSet ** | dataset, | ||
| int8_t | ptMeditRef, | ||
| int8_t | eltMeditRef, | ||
| int | nsols, | ||
| int8_t | metricData, | ||
| int8_t | lsData | ||
| ) | 
| mesh | pointer toward a MMG5 mesh | 
| sol | pointer toward a list of solution structures | 
| ptMeditRef | 1 if a point data field contains references (field named medit:ref) | 
| eltMeditRef | 1 if a cell data field contains references (field named medit:ref) | 
| nsols | number of point data (except the medit:ref ones) | 
| metricData | 1 if file contains a metric data highlighted by the :metric name | 
| lsData | 1 if file contains a metric data highlighted by the :ls name | 
I/O at Vtu VTK file format, part 2: mesh and solution storing
Definition at line 366 of file vtkparser.cpp.


| int MMG5_loadVtpMesh_part1 | ( | MMG5_pMesh | mesh, | 
| const char * | filename, | ||
| vtkDataSet ** | dataset, | ||
| int8_t * | ptMeditRef, | ||
| int8_t * | eltMeditRef, | ||
| int * | nsols, | ||
| int8_t * | metricData, | ||
| int8_t * | lsData | ||
| ) | 
| mesh | pointer toward a MMG5 mesh | 
| filename | name of the input file. | 
| dataset | vtkdataset structure | 
| ptMeditRef | index of point data field that contains references (field named medit:ref), -1 if no references | 
| eltMeditRef | index of a cell data field that contains references (field named medit:ref), -1 if no references | 
| nsols | number of point data (except the medit:ref ones) | 
| metricData | 1 if file contains a metric data highlighted by the :metric name | 
| lsData | 1 if file contains a metric data highlighted by the :ls name | 
I/O at Vtp VTK file format.
Definition at line 245 of file vtkparser.cpp.


| int MMG5_loadVtuMesh_part1 | ( | MMG5_pMesh | mesh, | 
| const char * | filename, | ||
| vtkDataSet ** | dataset, | ||
| int8_t * | ptMeditRef, | ||
| int8_t * | eltMeditRef, | ||
| int * | nsols, | ||
| int8_t * | metricData, | ||
| int8_t * | lsData | ||
| ) | 
| mesh | pointer toward a MMG5 mesh | 
| filename | pointer toward the filename | 
| dataset | vtkdataset structure | 
| ptMeditRef | index of point data field that contains references (field named medit:ref), -1 if no references | 
| eltMeditRef | index of a cell data field that contains references (field named medit:ref), -1 if no references | 
| nsols | number of point data (except the medit:ref ones) | 
| metricData | 1 if file contains a metric data highlighted by the :metric name | 
| lsData | 1 if file contains a metric data highlighted by the :ls name | 
I/O at Vtu VTK file format, part 1: file reading + count of the number of entities.
Definition at line 325 of file vtkparser.cpp.


| int MMG5_saveVtkMesh | ( | MMG5_pMesh | mesh, | 
| MMG5_pSol * | sol, | ||
| const char * | filename, | ||
| int | metricData, | ||
| int | binary | ||
| ) | 
| T | one of the VTK data class. | 
| TWriter | one of the VTK writer class. | 
| mesh | pointer toward a MMG5 mesh | 
| sol | pointer toward a MMG5 solution array | 
| filename | name of the input file. | 
| metricData | 1 if sol contains a metric array | 
| binary | 1 to save file at binary format (if available in TWriter class) | 
Save a vtk file at .vtp, .vtu or .vtk format.
Definition at line 539 of file vtkparser.hpp.

| int MMG5_saveVtkMesh_i | ( | MMG5_pMesh | mesh, | 
| MMG5_pSol * | sol, | ||
| const char * | mfilename, | ||
| int | metricData, | ||
| int | binary, | ||
| int | npart, | ||
| int | myid, | ||
| int | master | ||
| ) | 
| T | one of the VTK data class. | 
| TWriter | one of the VTK writer class. | 
| PWriter | one of the parallel VTK writer class. | 
| mesh | pointer toward a MMG5 mesh | 
| sol | pointer toward a MMG5 solution array | 
| mfilename | name of the master file to save (if call by master). | 
| metricData | 1 if sol contains a metric array | 
| binary | 1 to save file at binary format (if available in TWriter class) | 
| npart | number of parts of the saving | 
| myid | id of the process (save the file part number myid) | 
| master | id of the master process (save its part file + the master file) | 
Save a vtk file at .(p)vtp, .(p)vtu or .(p)vtk format.
Definition at line 170 of file vtkparser.hpp.
