Mmg
Simplicial remeshers (mesh adaptation, isovalue discretization, lagrangian movement)
|
#include <vtkSmartPointer.h>
#include <vtkXMLReader.h>
#include <vtkXMLWriter.h>
#include <vtkXMLUnstructuredGridReader.h>
#include <vtkXMLPolyDataReader.h>
#include <vtkDataSetReader.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 "vtkparser.hpp"
#include "libmmgtypes.h"
#include "mmgcommon_private.h"
Go to the source code of this file.
Functions | |
template<class TReader > | |
vtkDataSet * | MMG5_load_vtkXMLFile (const char *fileName) |
I/O at VTK format. | |
static int | MMG5_count_vtkEntities (vtkDataSet *dataset, MMG5_pMesh mesh, int8_t *ptMeditRef, int8_t *eltMeditRef, int *nsols, int8_t *metricData, int8_t *lsData) |
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) |
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) |
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) |
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) |
|
static |
dataset | pointer toward a vtkDataSet structure |
mesh | pointer toward a MMG5 mesh |
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 |
Count the number of entities of eache type (points, triangles...) in the mesh as well as the number of node data (solutions). If a data field name contains the "medit:ref" string, it will be used as entity reference instead of a solution. If node references are detected, ptMeditRef is incermented. If cell references are detected, eltMeditRef is incremented.
Definition at line 103 of file vtkparser.cpp.
vtkDataSet * MMG5_load_vtkXMLFile | ( | const char * | fileName | ) |
I/O at VTK format.
TReader | one of the VTK reader class. |
filename | name of the input file. |
Open and read a vtk xml file (.vtp, .vtu or .vtk). Store the data in a vtkDataSet object.
Definition at line 63 of file vtkparser.cpp.
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.