Mmg
Simplicial remeshers (mesh adaptation, isovalue discretization, lagrangian movement)
Functions
vtkparser.cpp File Reference
#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 dependency graph for vtkparser.cpp:

Go to the source code of this file.

Functions

template<class TReader >
vtkDataSet * MMG5_load_vtkXMLFile (const char *fileName)
 I/O at VTK format. More...
 
static int MMG5_count_vtkEntities (vtkDataSet *dataset, MMG5_pMesh mesh, int8_t *ptMeditRef, int8_t *eltMeditRef, int *nsols, int8_t *metricData)
 
int MMG5_loadVtpMesh_part1 (MMG5_pMesh mesh, const char *filename, vtkDataSet **dataset, int8_t *ptMeditRef, int8_t *eltMeditRef, int *nsols, int8_t *metricData)
 
int MMG5_loadVtkMesh_part1 (MMG5_pMesh mesh, const char *filename, vtkDataSet **dataset, int8_t *ptMeditRef, int8_t *eltMeditRef, int *nsols, int8_t *metricData)
 
int MMG5_loadVtuMesh_part1 (MMG5_pMesh mesh, const char *filename, vtkDataSet **dataset, int8_t *ptMeditRef, int8_t *eltMeditRef, int *nsols, int8_t *metricData)
 
int MMG5_loadVtkMesh_part2 (MMG5_pMesh mesh, MMG5_pSol *sol, vtkDataSet **dataset, int8_t ptMeditRef, int8_t eltMeditRef, int nsols)
 

Function Documentation

◆ MMG5_count_vtkEntities()

static int MMG5_count_vtkEntities ( vtkDataSet *  dataset,
MMG5_pMesh  mesh,
int8_t *  ptMeditRef,
int8_t *  eltMeditRef,
int *  nsols,
int8_t *  metricData 
)
static
Parameters
datasetpointer toward a vtkDataSet structure
meshpointer toward a MMG5 mesh
ptMeditRefindex of point data field that contains references (field named medit:ref), -1 if no references
eltMeditRefindex of a cell data field that contains references (field named medit:ref), -1 if no references
nsolsnumber of point data (except the medit:ref ones)
metricData1 if file contains a metric data highlighted by the :metric name
Returns
1 if success, -1 otherwise

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.

Remarks
Note that Mmg supports only 1 reference by point and 1 reference per cell.

Definition at line 100 of file vtkparser.cpp.

Here is the caller graph for this function:

◆ MMG5_load_vtkXMLFile()

template<class TReader >
vtkDataSet * MMG5_load_vtkXMLFile ( const char *  fileName)

I/O at VTK format.

Author
Charles Dapogny (UPMC)
Pascal Frey (UPMC)
Algiane Froehly (Inria)
Version
5
Template Parameters
TReaderone of the VTK reader class.
Parameters
filenamename of the input file.
Returns
pointer toward a vtkDataSet object that contains the mesh and the associated data.

Open and read a vtk xml file (.vtp, .vtu or .vtk). Store the data in a vtkDataSet object.

Definition at line 61 of file vtkparser.cpp.

◆ MMG5_loadVtkMesh_part1()

int MMG5_loadVtkMesh_part1 ( MMG5_pMesh  mesh,
const char *  filename,
vtkDataSet **  dataset,
int8_t *  ptMeditRef,
int8_t *  eltMeditRef,
int *  nsols,
int8_t *  metricData 
)
Parameters
meshpointer toward a MMG5 mesh
filenamename of the input file.
datasetvtkdataset structure
ptMeditRefindex of point data field that contains references (field named medit:ref), -1 if no references
eltMeditRefindex of a cell data field that contains references (field named medit:ref), -1 if no references
nsolsnumber of point data (except the medit:ref ones)
metricData1 if file contains a metric data highlighted by the :metric name
Returns
1 if success, 0 if fail to open/load the file, -1 otherwise;

I/O at Vtk VTK file format.

Definition at line 269 of file vtkparser.cpp.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ MMG5_loadVtkMesh_part2()

int MMG5_loadVtkMesh_part2 ( MMG5_pMesh  mesh,
MMG5_pSol sol,
vtkDataSet **  dataset,
int8_t  ptMeditRef,
int8_t  eltMeditRef,
int  nsols 
)
Parameters
meshpointer toward a MMG5 mesh
ptMeditRef1 if a point data field contains references (field named medit:ref)
eltMeditRef1 if a cell data field contains references (field named medit:ref)
nsolsnumber of point data (except the medit:ref ones)
Returns
1 if success, -1 if fail.

I/O at Vtu VTK file format, part 2: mesh and solution storing

Definition at line 346 of file vtkparser.cpp.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ MMG5_loadVtpMesh_part1()

int MMG5_loadVtpMesh_part1 ( MMG5_pMesh  mesh,
const char *  filename,
vtkDataSet **  dataset,
int8_t *  ptMeditRef,
int8_t *  eltMeditRef,
int *  nsols,
int8_t *  metricData 
)
Parameters
meshpointer toward a MMG5 mesh
filenamename of the input file.
datasetvtkdataset structure
ptMeditRefindex of point data field that contains references (field named medit:ref), -1 if no references
eltMeditRefindex of a cell data field that contains references (field named medit:ref), -1 if no references
nsolsnumber of point data (except the medit:ref ones)
metricData1 if file contains a metric data highlighted by the :metric name
Returns
1 if success, 0 if fail to open/load the file, -1 otherwise;

I/O at Vtp VTK file format.

Definition at line 230 of file vtkparser.cpp.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ MMG5_loadVtuMesh_part1()

int MMG5_loadVtuMesh_part1 ( MMG5_pMesh  mesh,
const char *  filename,
vtkDataSet **  dataset,
int8_t *  ptMeditRef,
int8_t *  eltMeditRef,
int *  nsols,
int8_t *  metricData 
)
Parameters
meshpointer toward a MMG5 mesh
filenamepointer toward the filename
datasetvtkdataset structure
ptMeditRefindex of point data field that contains references (field named medit:ref), -1 if no references
eltMeditRefindex of a cell data field that contains references (field named medit:ref), -1 if no references
nsolsnumber of point data (except the medit:ref ones)
metricData1 if file contains a metric data highlighted by the :metric name
Returns
1 if success, 0 if fail to open/load the file, -1 other errors;

I/O at Vtu VTK file format, part 1: file reading + count of the number of entities.

Definition at line 308 of file vtkparser.cpp.

Here is the call graph for this function:
Here is the caller graph for this function: