Mmg
Simplicial remeshers (mesh adaptation, isovalue discretization, lagrangian movement)
inoutcpp_s.cpp
Go to the documentation of this file.
1/* =============================================================================
2** This file is part of the mmg software package for the tetrahedral
3** mesh modification.
4** Copyright (c) Bx INP/CNRS/Inria/UBordeaux/UPMC, 2004-
5**
6** mmg is free software: you can redistribute it and/or modify it
7** under the terms of the GNU Lesser General Public License as published
8** by the Free Software Foundation, either version 3 of the License, or
9** (at your option) any later version.
10**
11** mmg is distributed in the hope that it will be useful, but WITHOUT
12** ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
13** FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
14** License for more details.
15**
16** You should have received a copy of the GNU Lesser General Public
17** License and of the GNU General Public License along with mmg (in
18** files COPYING.LESSER and COPYING). If not, see
19** <http://www.gnu.org/licenses/>. Please read their terms carefully and
20** use this copy of the mmg distribution only if you accept them.
21** =============================================================================
22*/
23
34#include "vtkparser.hpp"
35#include "libmmgs.h"
36#include "libmmgs_private.h"
37
38#ifdef USE_VTK
39
41 vtkDataSet **dataset, int8_t ptMeditRef,
42 int8_t eltMeditRef,MMG5_int nsols,
43 int8_t metricData, int8_t lsData ) {
44 int ier;
45
46 if ( !MMGS_zaldy(mesh) ) {
47 return 0;
48 }
49
50 mesh->ne = mesh->nprism = 0;
51
52 if ( !mesh->nt ) {
53 fprintf(stderr," ** MISSING DATA.\n");
54 fprintf(stderr," Check that your mesh contains triangles.\n");
55 fprintf(stderr," Exit program.\n");
56 return -1;
57 }
58
59 if (mesh->npmax < mesh->np || mesh->ntmax < mesh->nt ) {
60 return -1;
61 }
62
63 ier = MMG5_loadVtkMesh_part2(mesh,sol,dataset,ptMeditRef,eltMeditRef,nsols,metricData,lsData);
64
65 if ( ier < 1 ) {
66 fprintf(stderr," ** ERROR WHEN PARSING THE INPUT FILE\n");
67 return ier;
68 }
69
70 return ier;
71}
72
73#endif
74
76
77#ifndef USE_VTK
78
79 fprintf(stderr," ** VTK library not founded. Unavailable file format.\n");
80 return -1;
81
82#else
83 int ier;
84 int nsols;
85 int8_t ptMeditRef,eltMeditRef,metricData,lsData;
86 vtkDataSet *dataset;
87 MMG5_pSol allSol[2];
88
89 mesh->dim = 3;
90
91 ier = MMG5_loadVtpMesh_part1(mesh,filename,&dataset,&ptMeditRef,&eltMeditRef,
92 &nsols,&metricData,&lsData);
93 if ( ier < 1 ) return ier;
94
95 /* Check data fields */
96 if ( nsols > (metricData+lsData) ) {
97 fprintf(stderr,"Error: %d UNEXPECTED DATA FIELD(S)\n",nsols-metricData-lsData);
98 return -1;
99 }
100
101 // Mesh alloc and transfer of the mesh from dataset toward the MMG5 Mesh Sol
102 allSol[0] = met;
103 allSol[1] = sol;
104 ier = MMGS_loadVtkMesh_part2(mesh,allSol,&dataset,ptMeditRef,eltMeditRef,nsols,metricData,lsData);
105 if ( ier < 1 ) {
106 fprintf(stderr," ** ERROR WHEN PARSING THE INPUT FILE\n");
107 return ier;
108 }
109
110 if ( sol ) {
111 /* Check the metric type */
112 ier = MMG5_chkMetricType(mesh,&sol->type,&sol->entities,NULL);
113 }
114
115 return ier;
116#endif
117}
118
120
121#ifndef USE_VTK
122
123 fprintf(stderr," ** VTK library not founded. Unavailable file format.\n");
124 return -1;
125
126#else
127 int ier;
128 int nsols;
129 int8_t ptMeditRef,eltMeditRef,metricData,lsData;
130 vtkDataSet *dataset;
131 MMG5_pSol allSol[2];
132
133 mesh->dim = 3;
134
135 ier = MMG5_loadVtpMesh_part1(mesh,filename,&dataset,&ptMeditRef,&eltMeditRef,
136 &nsols,&metricData,&lsData);
137 if ( ier < 1 ) return ier;
138
139 mesh->nsols = nsols;
140 if ( *sol ) MMG5_DEL_MEM(mesh,*sol);
141 MMG5_ADD_MEM(mesh,nsols*sizeof(MMG5_Sol),"solutions array",
142 printf(" Exit program.\n");
143 return -1);
144 MMG5_SAFE_CALLOC(*sol,nsols,MMG5_Sol,return -1);
145
146 // Mesh alloc and transfer of the mesh from dataset toward the MMG5 Mesh Sol
147 allSol[0] = NULL;
148 allSol[1] = *sol;
149 ier = MMGS_loadVtkMesh_part2(mesh,allSol,&dataset,ptMeditRef,eltMeditRef,nsols,metricData,0);
150
151 return ier;
152#endif
153}
154
156
157#ifndef USE_VTK
158
159 fprintf(stderr," ** VTK library not founded. Unavailable file format.\n");
160 return -1;
161
162#else
163 int ier;
164 int nsols;
165 int8_t ptMeditRef,eltMeditRef,metricData,lsData;
166 vtkDataSet *dataset;
167 MMG5_pSol allSol[2];
168
169 mesh->dim = 3;
170
171 ier = MMG5_loadVtkMesh_part1(mesh,filename,&dataset,&ptMeditRef,&eltMeditRef,
172 &nsols,&metricData,&lsData);
173 if ( ier < 1 ) return ier;
174
175 /* Check data fields */
176 if ( nsols > (metricData+lsData) ) {
177 fprintf(stderr,"Error: %d UNEXPECTED DATA FIELD(S)\n",nsols-metricData-lsData);
178 return -1;
179 }
180
181 // Mesh alloc and transfer of the mesh from dataset toward the MMG5 Mesh Sol
182 allSol[0] = met;
183 allSol[1] = sol;
184 ier = MMGS_loadVtkMesh_part2(mesh,allSol,&dataset,ptMeditRef,eltMeditRef,nsols,metricData,lsData);
185 if ( ier < 1 ) {
186 fprintf(stderr," ** ERROR WHEN PARSING THE INPUT FILE\n");
187 return ier;
188 }
189
190 if ( sol ) {
191 /* Check the metric type */
192 ier = MMG5_chkMetricType(mesh,&sol->type,&sol->entities,NULL);
193 }
194
195 return ier;
196#endif
197}
198
200
201#ifndef USE_VTK
202
203 fprintf(stderr," ** VTK library not founded. Unavailable file format.\n");
204 return -1;
205
206#else
207 int ier;
208 int nsols;
209 int8_t ptMeditRef,eltMeditRef,metricData,lsData;
210 vtkDataSet *dataset;
211 MMG5_pSol allSol[2];
212
213 mesh->dim = 3;
214
215 ier = MMG5_loadVtkMesh_part1(mesh,filename,&dataset,&ptMeditRef,&eltMeditRef,
216 &nsols,&metricData,&lsData);
217 if ( ier < 1 ) return ier;
218
219 mesh->nsols = nsols;
220 if ( *sol ) MMG5_DEL_MEM(mesh,*sol);
221 MMG5_ADD_MEM(mesh,nsols*sizeof(MMG5_Sol),"solutions array",
222 printf(" Exit program.\n");
223 return -1);
224 MMG5_SAFE_CALLOC(*sol,nsols,MMG5_Sol,return -1);
225
226 // Mesh alloc and transfer of the mesh from dataset toward the MMG5 Mesh Sol
227 allSol[0] = NULL;
228 allSol[1] = *sol;
229 ier = MMGS_loadVtkMesh_part2(mesh,allSol,&dataset,ptMeditRef,eltMeditRef,nsols,metricData,0);
230
231 return ier;
232#endif
233}
234
236
237#ifndef USE_VTK
238
239 fprintf(stderr," ** VTK library not founded. Unavailable file format.\n");
240 return -1;
241
242#else
243 int ier;
244 int nsols;
245 int8_t ptMeditRef,eltMeditRef,metricData,lsData;
246 vtkDataSet *dataset;
247 MMG5_pSol allSol[2];
248
249 mesh->dim = 3;
250
251 ier = MMG5_loadVtuMesh_part1(mesh,filename,&dataset,&ptMeditRef,&eltMeditRef,
252 &nsols,&metricData,&lsData);
253 if ( ier < 1 ) return ier;
254
255 /* Check data fields */
256 if ( nsols > (metricData+lsData) ) {
257 fprintf(stderr,"Error: %d UNEXPECTED DATA FIELD(S)\n",nsols-metricData-lsData);
258 return -1;
259 }
260
261 // Mesh alloc and transfer of the mesh from dataset toward the MMG5 Mesh Sol
262 allSol[0] = met;
263 allSol[1] = sol;
264 ier = MMGS_loadVtkMesh_part2(mesh,allSol,&dataset,ptMeditRef,eltMeditRef,nsols,metricData,lsData);
265 if ( ier < 1 ) {
266 fprintf(stderr," ** ERROR WHEN PARSING THE INPUT FILE\n");
267 return ier;
268 }
269
270 if ( sol ) {
271 /* Check the metric type */
272 ier = MMG5_chkMetricType(mesh,&sol->type,&sol->entities,NULL);
273 }
274
275 return ier;
276#endif
277}
278
280
281#ifndef USE_VTK
282
283 fprintf(stderr," ** VTK library not founded. Unavailable file format.\n");
284 return -1;
285
286#else
287 int ier;
288 int nsols;
289 int8_t ptMeditRef,eltMeditRef,metricData,lsData;
290 vtkDataSet *dataset;
291 MMG5_pSol allSol[2];
292
293 mesh->dim = 3;
294
295 ier = MMG5_loadVtuMesh_part1(mesh,filename,&dataset,&ptMeditRef,&eltMeditRef,
296 &nsols,&metricData,&lsData);
297 if ( ier < 1 ) return ier;
298
299 mesh->nsols = nsols;
300 if ( *sol ) MMG5_DEL_MEM(mesh,*sol);
301 MMG5_ADD_MEM(mesh,nsols*sizeof(MMG5_Sol),"solutions array",
302 printf(" Exit program.\n");
303 return -1);
304 MMG5_SAFE_CALLOC(*sol,nsols,MMG5_Sol,return -1);
305
306 // Mesh alloc and transfer of the mesh from dataset toward the MMG5 Mesh Sol
307 allSol[0] = NULL;
308 allSol[1] = *sol;
309 ier = MMGS_loadVtkMesh_part2(mesh,allSol,&dataset,ptMeditRef,eltMeditRef,nsols,metricData,0);
310
311 return ier;
312#endif
313}
314
316
317#ifndef USE_VTK
318
319 fprintf(stderr," ** VTK library not found. Unavailable file format.\n");
320 return -1;
321
322#else
323
324 return MMG5_saveVtkMesh<vtkUnstructuredGrid,vtkXMLUnstructuredGridWriter,
325 vtkXMLPUnstructuredGridWriter>(mesh,&sol,filename,1,1);
326
327#endif
328}
329
331
332 MMG5_pSol allSol[2];
333
334#ifndef USE_VTK
335
336 fprintf(stderr," ** VTK library not found. Unavailable file format.\n");
337 return -1;
338
339#else
340
341 allSol[0] = NULL;
342 allSol[1] = *sol;
343
344 return MMG5_saveVtkMesh<vtkUnstructuredGrid,vtkXMLUnstructuredGridWriter,
345 vtkXMLPUnstructuredGridWriter>(mesh,allSol,filename,0,1);
346
347#endif
348}
349
351
352#ifndef USE_VTK
353
354 fprintf(stderr," ** VTK library not found. Unavailable file format.\n");
355 return -1;
356
357#else
358
359 return MMG5_saveVtkMesh<vtkUnstructuredGrid,vtkDataSetWriter,
360 vtkPDataSetWriter>(mesh,&sol,filename,1,0);
361
362#endif
363}
364
366
367 MMG5_pSol allSol[2];
368
369#ifndef USE_VTK
370
371 fprintf(stderr," ** VTK library not found. Unavailable file format.\n");
372 return -1;
373
374#else
375
376 allSol[0] = NULL;
377 allSol[1] = *sol;
378
379 return MMG5_saveVtkMesh<vtkUnstructuredGrid,vtkDataSetWriter,
380 vtkPDataSetWriter>(mesh,allSol,filename,0,0);
381
382#endif
383}
384
386
387#ifndef USE_VTK
388
389 fprintf(stderr," ** VTK library not found. Unavailable file format.\n");
390 return -1;
391
392#else
393
394 return MMG5_saveVtkMesh<vtkPolyData,vtkXMLPolyDataWriter,
395 vtkXMLPPolyDataWriter>(mesh,&sol,filename,1,1);
396
397#endif
398}
399
401
402 MMG5_pSol allSol[2];
403
404#ifndef USE_VTK
405
406 fprintf(stderr," ** VTK library not found. Unavailable file format.\n");
407 return -1;
408
409#else
410
411 allSol[0] = NULL;
412 allSol[1] = *sol;
413
414 return MMG5_saveVtkMesh<vtkPolyData,vtkXMLPolyDataWriter,
415 vtkXMLPPolyDataWriter>(mesh,allSol,filename,0,1);
416
417#endif
418}
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 MMGS_saveVtuMesh_and_allData(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename)
Write a mesh and multiple data fields in vtu Vtk file format (.vtu extension).
Definition: inoutcpp_s.cpp:330
int MMGS_loadVtuMesh(MMG5_pMesh mesh, MMG5_pSol met, MMG5_pSol sol, const char *filename)
Load a mesh and possibly data in VTU (VTK) format from file.
Definition: inoutcpp_s.cpp:235
int MMGS_saveVtkMesh_and_allData(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename)
Save a mesh and multiple data fields in VTK format.
Definition: inoutcpp_s.cpp:365
int MMGS_saveVtuMesh(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
Write mesh and optionally one data field vtu Vtk file format (.vtu extension).
Definition: inoutcpp_s.cpp:315
int MMGS_loadVtpMesh(MMG5_pMesh mesh, MMG5_pSol met, MMG5_pSol sol, const char *filename)
Load a mesh and optionally a solution in VTP (VTK) format from file.
Definition: inoutcpp_s.cpp:75
int MMGS_saveVtkMesh(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
Write mesh and optionally one data field in Vtk file format (.vtk extension).
Definition: inoutcpp_s.cpp:350
int MMGS_loadVtpMesh_and_allData(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename)
Load a mesh and multiple solutions in VTP (VTK) format from file.
Definition: inoutcpp_s.cpp:119
static int MMGS_loadVtkMesh_part2(MMG5_pMesh mesh, MMG5_pSol *sol, vtkDataSet **dataset, int8_t ptMeditRef, int8_t eltMeditRef, MMG5_int nsols, int8_t metricData, int8_t lsData)
Input / Output Functions that needs cpp features.
Definition: inoutcpp_s.cpp:40
int MMGS_loadVtuMesh_and_allData(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename)
Load a mesh and multiple solutions in VTU (VTK) format from file.
Definition: inoutcpp_s.cpp:279
int MMGS_loadVtkMesh_and_allData(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename)
Load a mesh and multiple solutions in VTK format from file.
Definition: inoutcpp_s.cpp:199
int MMGS_loadVtkMesh(MMG5_pMesh mesh, MMG5_pSol met, MMG5_pSol sol, const char *filename)
Load a mesh and possibly data in VTK format from file.
Definition: inoutcpp_s.cpp:155
int MMGS_saveVtpMesh(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
Save a mesh and optionally one data field in VTP format.
Definition: inoutcpp_s.cpp:385
int MMGS_saveVtpMesh_and_allData(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename)
Save a mesh and multiple data fields in VTP format.
Definition: inoutcpp_s.cpp:400
API headers and documentation for the mmgs library.
int MMGS_zaldy(MMG5_pMesh mesh)
Definition: zaldy_s.c:263
#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 npmax
Definition: libmmgtypes.h:620
MMG5_int nt
Definition: libmmgtypes.h:620
MMG5_int np
Definition: libmmgtypes.h:620
MMG5_int nprism
Definition: libmmgtypes.h:621
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_loadVtpMesh_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:245
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