Mmg
Simplicial remeshers (mesh adaptation, isovalue discretization, lagrangian movement)
inoutcpp_3d.cpp
Go to the documentation of this file.
1
2/* =============================================================================
3** This file is part of the mmg software package for the tetrahedral
4** mesh modification.
5** Copyright (c) Bx INP/CNRS/Inria/UBordeaux/UPMC, 2004-
6**
7** mmg is free software: you can redistribute it and/or modify it
8** under the terms of the GNU Lesser General Public License as published
9** by the Free Software Foundation, either version 3 of the License, or
10** (at your option) any later version.
11**
12** mmg is distributed in the hope that it will be useful, but WITHOUT
13** ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
14** FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
15** License for more details.
16**
17** You should have received a copy of the GNU Lesser General Public
18** License and of the GNU General Public License along with mmg (in
19** files COPYING.LESSER and COPYING). If not, see
20** <http://www.gnu.org/licenses/>. Please read their terms carefully and
21** use this copy of the mmg distribution only if you accept them.
22** =============================================================================
23*/
24
35#include "vtkparser.hpp"
36
37#include "libmmg3d.h"
38#include "libmmg3d_private.h"
39
40#ifdef USE_VTK
41
43 vtkDataSet **dataset, int8_t ptMeditRef,
44 int8_t eltMeditRef,int nsols,
45 int8_t metricData, int8_t lsData) {
46
47 int ier;
48
49 if ( !MMG3D_zaldy(mesh) ) {
50 return -1;
51 }
52 if (mesh->npmax < mesh->np || mesh->ntmax < mesh->nt || mesh->nemax < mesh->ne) {
53 return -1;
54 }
55 if ( !mesh->ne ) {
56 fprintf(stderr," ** MISSING DATA.\n");
57 fprintf(stderr," Check that your mesh contains tetrahedra.\n");
58 fprintf(stderr," Exit program.\n");
59 return -1;
60 }
61
62 ier = MMG5_loadVtkMesh_part2(mesh,sol,dataset,ptMeditRef,eltMeditRef,nsols,metricData,lsData);
63
64 if ( ier < 1 ) {
65 fprintf(stderr," ** ERROR WHEN PARSING THE INPUT FILE\n");
66 return ier;
67 }
68
69 return ier;
70}
71#endif
72
73
75
76#ifndef USE_VTK
77
78 fprintf(stderr," ** VTK library not founded. Unavailable file format.\n");
79 return -1;
80
81#else
82 int ier,nsols;
83 int8_t ptMeditRef,eltMeditRef,metricData,lsData;
84 vtkDataSet *dataset;
85 MMG5_pSol allSol[2];
86
87 mesh->dim = 3;
88
89 ier = MMG5_loadVtuMesh_part1(mesh,filename,&dataset,&ptMeditRef,&eltMeditRef,
90 &nsols,&metricData,&lsData);
91 if ( ier < 1 ) return ier;
92
93 /* Check data fields */
94 if ( nsols > (metricData+lsData) ) {
95 fprintf(stderr,"Error: %d UNEXPECTED DATA FIELD(S)\n",nsols-metricData-lsData);
96 return -1;
97 }
98
99 // Mesh alloc and transfer of the mesh from dataset toward the MMG5 Mesh Sol
100 allSol[0] = met;
101 allSol[1] = sol;
102 ier = MMG3D_loadVtkMesh_part2(mesh,allSol,&dataset,ptMeditRef,eltMeditRef,nsols,metricData,lsData);
103 if ( ier < 1 ) {
104 fprintf(stderr," ** ERROR WHEN PARSING THE INPUT FILE\n");
105 return ier;
106 }
107
108 if ( sol ) {
109 /* Check the metric type */
110 ier = MMG5_chkMetricType(mesh,&sol->type,&sol->entities,NULL);
111 }
112
113 return ier;
114#endif
115}
116
118
119#ifndef USE_VTK
120
121 fprintf(stderr," ** VTK library not founded. Unavailable file format.\n");
122 return -1;
123
124#else
125 int ier,nsols;
126 int8_t ptMeditRef,eltMeditRef,metricData,lsData;
127 vtkDataSet *dataset;
128 MMG5_pSol allSol[2];
129
130 mesh->dim = 3;
131
132 ier = MMG5_loadVtuMesh_part1(mesh,filename,&dataset,&ptMeditRef,&eltMeditRef,
133 &nsols,&metricData,&lsData);
134 if ( ier < 1 ) return ier;
135
136 mesh->nsols = nsols;
137 if ( *sol ) MMG5_DEL_MEM(mesh,*sol);
138
139 MMG5_ADD_MEM(mesh,nsols*sizeof(MMG5_Sol),"solutions array",
140 printf(" Exit program.\n");
141 return -1);
142 MMG5_SAFE_CALLOC(*sol,nsols,MMG5_Sol,return -1);
143
144 // Mesh alloc and transfer of the mesh from dataset toward the MMG5 Mesh Sol
145 allSol[0] = NULL;
146 allSol[1] = *sol;
147 ier = MMG3D_loadVtkMesh_part2(mesh,allSol,&dataset,ptMeditRef,eltMeditRef,nsols,metricData,0);
148
149 return ier;
150#endif
151}
152
154
155#ifndef USE_VTK
156
157 fprintf(stderr," ** VTK library not founded. Unavailable file format.\n");
158 return -1;
159
160#else
161 int ier,nsols;
162 int8_t ptMeditRef,eltMeditRef,metricData,lsData;
163 vtkDataSet *dataset;
164 MMG5_pSol allSol[2];
165
166 mesh->dim = 3;
167
168 ier = MMG5_loadVtkMesh_part1(mesh,filename,&dataset,&ptMeditRef,&eltMeditRef,
169 &nsols,&metricData,&lsData);
170 if ( ier < 1 ) return ier;
171
172 /* Check data fields */
173 if ( nsols > (metricData+lsData) ) {
174 fprintf(stderr,"Error: %d UNEXPECTED DATA FIELD(S)\n",nsols-metricData-lsData);
175 return -1;
176 }
177
178 // Mesh alloc and transfer of the mesh from dataset toward the MMG5 Mesh Sol
179 allSol[0] = met;
180 allSol[1] = sol;
181 ier = MMG3D_loadVtkMesh_part2(mesh,allSol,&dataset,ptMeditRef,eltMeditRef,nsols,metricData,lsData);
182 if ( ier < 1 ) {
183 fprintf(stderr," ** ERROR WHEN PARSING THE INPUT FILE\n");
184 return ier;
185 }
186
187 if ( sol ) {
188 /* Check the metric type */
189 ier = MMG5_chkMetricType(mesh,&sol->type,&sol->entities,NULL);
190 }
191
192 return ier;
193#endif
194}
195
197
198#ifndef USE_VTK
199
200 fprintf(stderr," ** VTK library not founded. Unavailable file format.\n");
201 return -1;
202
203#else
204 int ier,nsols;
205 int8_t ptMeditRef,eltMeditRef,metricData, lsData;
206 vtkDataSet *dataset;
207 MMG5_pSol allSol[2];
208
209 mesh->dim = 3;
210
211 ier = MMG5_loadVtkMesh_part1(mesh,filename,&dataset,&ptMeditRef,&eltMeditRef,
212 &nsols,&metricData,&lsData);
213 if ( ier < 1 ) return ier;
214
215 mesh->nsols = nsols;
216 if ( *sol ) MMG5_DEL_MEM(mesh,*sol);
217
218 MMG5_ADD_MEM(mesh,nsols*sizeof(MMG5_Sol),"solutions array",
219 printf(" Exit program.\n");
220 return -1);
221 MMG5_SAFE_CALLOC(*sol,nsols,MMG5_Sol,return -1);
222
223 // Mesh alloc and transfer of the mesh from dataset toward the MMG5 Mesh Sol
224 allSol[0] = NULL;
225 allSol[1] = *sol;
226 ier = MMG3D_loadVtkMesh_part2(mesh,allSol,&dataset,ptMeditRef,eltMeditRef,nsols,metricData,0);
227
228 return ier;
229#endif
230}
231
233
234#ifndef USE_VTK
235
236 fprintf(stderr," ** VTK library not founded. Unavailable file format.\n");
237 return -1;
238
239#else
240 return MMG5_saveVtkMesh<vtkUnstructuredGrid,vtkXMLUnstructuredGridWriter,
241 vtkXMLPUnstructuredGridWriter>(mesh,&sol,filename,1,1);
242
243#endif
244}
245
247
248 MMG5_pSol allSol[2];
249
250#ifndef USE_VTK
251
252 fprintf(stderr," ** VTK library not founded. Unavailable file format.\n");
253 return -1;
254
255#else
256
257 allSol[0] = NULL;
258 allSol[1] = *sol;
259
260 return MMG5_saveVtkMesh<vtkUnstructuredGrid,vtkXMLUnstructuredGridWriter,
261 vtkXMLPUnstructuredGridWriter>(mesh,allSol,filename,0,1);
262
263#endif
264}
265
267
268#ifndef USE_VTK
269
270 fprintf(stderr," ** VTK library not founded. Unavailable file format.\n");
271 return -1;
272
273#else
274
275 return MMG5_saveVtkMesh<vtkUnstructuredGrid,vtkDataSetWriter,
276 vtkPDataSetWriter>(mesh,&sol,filename,1,0);
277
278#endif
279}
280
282
283 MMG5_pSol allSol[2];
284
285#ifndef USE_VTK
286
287 fprintf(stderr," ** VTK library not founded. Unavailable file format.\n");
288 return -1;
289
290#else
291
292 allSol[0] = NULL;
293 allSol[1] = *sol;
294
295 return MMG5_saveVtkMesh<vtkUnstructuredGrid,vtkDataSetWriter,
296 vtkPDataSetWriter>(mesh,allSol,filename,0,0);
297
298#endif
299}
int ier
MMG5_pMesh MMG5_pSol * sol
MMG5_pMesh * mesh
MMG5_pMesh char * filename
int MMG5_chkMetricType(MMG5_pMesh mesh, int *type, int *entities, FILE *inm)
Definition: inout.c:2663
int MMG3D_loadVtuMesh_and_allData(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename)
Load a mesh and multiple solutions in VTU (VTK) format from file.
int MMG3D_saveVtuMesh_and_allData(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename)
Save a mesh and multiple data fields in VTU format.
int MMG3D_loadVtkMesh_and_allData(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename)
Load a mesh and multiple solutions from a file in VTK format.
int MMG3D_loadVtkMesh(MMG5_pMesh mesh, MMG5_pSol met, MMG5_pSol sol, const char *filename)
Load a mesh and possibly a solution from a file in VTK format.
static int MMG3D_loadVtkMesh_part2(MMG5_pMesh mesh, MMG5_pSol *sol, vtkDataSet **dataset, int8_t ptMeditRef, int8_t eltMeditRef, int nsols, int8_t metricData, int8_t lsData)
Input / Output Functions that needs cpp features.
Definition: inoutcpp_3d.cpp:42
int MMG3D_saveVtkMesh(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
Save a mesh and optionally one solution in VTK format.
int MMG3D_loadVtuMesh(MMG5_pMesh mesh, MMG5_pSol met, MMG5_pSol sol, const char *filename)
Load a mesh and possibly a solution in VTU (VTK) format from file.
Definition: inoutcpp_3d.cpp:74
int MMG3D_saveVtkMesh_and_allData(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename)
Save a mesh and multiple data fields in VTK format.
int MMG3D_saveVtuMesh(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
Save a mesh and optionally one data field in VTU format.
API headers and documentation for the mmg3d library, for volumetric meshes in 3D.
int MMG3D_zaldy(MMG5_pMesh mesh)
Definition: zaldy_3d.c:346
#define MMG5_SAFE_CALLOC(ptr, size, type, law)
#define MMG5_ADD_MEM(mesh, size, message, law)
#define MMG5_DEL_MEM(mesh, ptr)
MMG mesh structure.
Definition: libmmgtypes.h:613
MMG5_int ntmax
Definition: libmmgtypes.h:620
MMG5_int ne
Definition: libmmgtypes.h:620
MMG5_int nemax
Definition: libmmgtypes.h:620
MMG5_int npmax
Definition: libmmgtypes.h:620
MMG5_int nt
Definition: libmmgtypes.h:620
MMG5_int np
Definition: libmmgtypes.h:620
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)
Definition: vtkparser.cpp:285
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)
Definition: vtkparser.cpp:325
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)
Definition: vtkparser.cpp:366
int MMG5_saveVtkMesh(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename, int metricData, int binary)
Definition: vtkparser.hpp:539