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 int ier;
44
45 if ( !MMGS_zaldy(mesh) ) {
46 return 0;
47 }
48
49 mesh->ne = mesh->nprism = 0;
50
51 if ( !mesh->nt ) {
52 fprintf(stderr," ** MISSING DATA.\n");
53 fprintf(stderr," Check that your mesh contains triangles.\n");
54 fprintf(stderr," Exit program.\n");
55 return -1;
56 }
57
58 if (mesh->npmax < mesh->np || mesh->ntmax < mesh->nt ) {
59 return -1;
60 }
61
62 ier = MMG5_loadVtkMesh_part2(mesh,sol,dataset,ptMeditRef,eltMeditRef,nsols);
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
72#endif
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;
83 int nsols;
84 int8_t ptMeditRef,eltMeditRef,metricData;
85 vtkDataSet *dataset;
86
87 mesh->dim = 3;
88
89 ier = MMG5_loadVtpMesh_part1(mesh,filename,&dataset,&ptMeditRef,&eltMeditRef,
90 &nsols,&metricData);
91 if ( ier < 1 ) return ier;
92
93 /* Check data fields */
94 if ( nsols > metricData ) {
95 fprintf(stderr,"Error: %d UNEXPECTED DATA FIELD(S)\n",nsols);
96 return -1;
97 }
98
99 // Mesh alloc and transfer of the mesh from dataset toward the MMG5 Mesh Sol
100 ier = MMGS_loadVtkMesh_part2(mesh,&sol,&dataset,ptMeditRef,eltMeditRef,nsols);
101 if ( ier < 1 ) {
102 fprintf(stderr," ** ERROR WHEN PARSING THE INPUT FILE\n");
103 return ier;
104 }
105
106 if ( sol ) {
107 /* Check the metric type */
108 ier = MMG5_chkMetricType(mesh,&sol->type,&sol->entities,NULL);
109 }
110
111 return ier;
112#endif
113}
114
116
117#ifndef USE_VTK
118
119 fprintf(stderr," ** VTK library not founded. Unavailable file format.\n");
120 return -1;
121
122#else
123 int ier;
124 int nsols;
125 int8_t ptMeditRef,eltMeditRef,metricData;
126 vtkDataSet *dataset;
127
128 mesh->dim = 3;
129
130 ier = MMG5_loadVtpMesh_part1(mesh,filename,&dataset,&ptMeditRef,&eltMeditRef,
131 &nsols,&metricData);
132 if ( ier < 1 ) return ier;
133
134 mesh->nsols = nsols;
135 if ( *sol ) MMG5_DEL_MEM(mesh,*sol);
136 MMG5_ADD_MEM(mesh,nsols*sizeof(MMG5_Sol),"solutions array",
137 printf(" Exit program.\n");
138 return -1);
139 MMG5_SAFE_CALLOC(*sol,nsols,MMG5_Sol,return -1);
140
141 // Mesh alloc and transfer of the mesh from dataset toward the MMG5 Mesh Sol
142 ier = MMGS_loadVtkMesh_part2(mesh,sol,&dataset,ptMeditRef,eltMeditRef,nsols);
143
144 return ier;
145#endif
146}
147
149
150#ifndef USE_VTK
151
152 fprintf(stderr," ** VTK library not founded. Unavailable file format.\n");
153 return -1;
154
155#else
156 int ier;
157 int nsols;
158 int8_t ptMeditRef,eltMeditRef,metricData;
159 vtkDataSet *dataset;
160
161 mesh->dim = 3;
162
163 ier = MMG5_loadVtkMesh_part1(mesh,filename,&dataset,&ptMeditRef,&eltMeditRef,
164 &nsols,&metricData);
165 if ( ier < 1 ) return ier;
166
167 /* Check data fields */
168 if ( nsols > metricData ) {
169 fprintf(stderr,"Error: %d UNEXPECTED DATA FIELD(S)\n",nsols);
170 return -1;
171 }
172
173 // Mesh alloc and transfer of the mesh from dataset toward the MMG5 Mesh Sol
174 ier = MMGS_loadVtkMesh_part2(mesh,&sol,&dataset,ptMeditRef,eltMeditRef,nsols);
175 if ( ier < 1 ) {
176 fprintf(stderr," ** ERROR WHEN PARSING THE INPUT FILE\n");
177 return ier;
178 }
179
180 if ( sol ) {
181 /* Check the metric type */
182 ier = MMG5_chkMetricType(mesh,&sol->type,&sol->entities,NULL);
183 }
184
185 return ier;
186#endif
187}
188
190
191#ifndef USE_VTK
192
193 fprintf(stderr," ** VTK library not founded. Unavailable file format.\n");
194 return -1;
195
196#else
197 int ier;
198 int nsols;
199 int8_t ptMeditRef,eltMeditRef,metricData;
200 vtkDataSet *dataset;
201
202 mesh->dim = 3;
203
204 ier = MMG5_loadVtkMesh_part1(mesh,filename,&dataset,&ptMeditRef,&eltMeditRef,
205 &nsols,&metricData);
206 if ( ier < 1 ) return ier;
207
208 mesh->nsols = nsols;
209 if ( *sol ) MMG5_DEL_MEM(mesh,*sol);
210 MMG5_ADD_MEM(mesh,nsols*sizeof(MMG5_Sol),"solutions array",
211 printf(" Exit program.\n");
212 return -1);
213 MMG5_SAFE_CALLOC(*sol,nsols,MMG5_Sol,return -1);
214
215 // Mesh alloc and transfer of the mesh from dataset toward the MMG5 Mesh Sol
216 ier = MMGS_loadVtkMesh_part2(mesh,sol,&dataset,ptMeditRef,eltMeditRef,nsols);
217
218 return ier;
219#endif
220}
221
223
224#ifndef USE_VTK
225
226 fprintf(stderr," ** VTK library not founded. Unavailable file format.\n");
227 return -1;
228
229#else
230 int ier;
231 int nsols;
232 int8_t ptMeditRef,eltMeditRef,metricData;
233 vtkDataSet *dataset;
234
235 mesh->dim = 3;
236
237 ier = MMG5_loadVtuMesh_part1(mesh,filename,&dataset,&ptMeditRef,&eltMeditRef,
238 &nsols,&metricData);
239 if ( ier < 1 ) return ier;
240
241 /* Check data fields */
242 if ( nsols > metricData ) {
243 fprintf(stderr,"Error: %d UNEXPECTED DATA FIELD(S)\n",nsols);
244 return -1;
245 }
246
247 // Mesh alloc and transfer of the mesh from dataset toward the MMG5 Mesh Sol
248 ier = MMGS_loadVtkMesh_part2(mesh,&sol,&dataset,ptMeditRef,eltMeditRef,nsols);
249 if ( ier < 1 ) {
250 fprintf(stderr," ** ERROR WHEN PARSING THE INPUT FILE\n");
251 return ier;
252 }
253
254 if ( sol ) {
255 /* Check the metric type */
256 ier = MMG5_chkMetricType(mesh,&sol->type,&sol->entities,NULL);
257 }
258
259 return ier;
260#endif
261}
262
264
265#ifndef USE_VTK
266
267 fprintf(stderr," ** VTK library not founded. Unavailable file format.\n");
268 return -1;
269
270#else
271 int ier;
272 int nsols;
273 int8_t ptMeditRef,eltMeditRef,metricData;
274 vtkDataSet *dataset;
275
276 mesh->dim = 3;
277
278 ier = MMG5_loadVtuMesh_part1(mesh,filename,&dataset,&ptMeditRef,&eltMeditRef,
279 &nsols,&metricData);
280 if ( ier < 1 ) return ier;
281
282 mesh->nsols = nsols;
283 if ( *sol ) MMG5_DEL_MEM(mesh,*sol);
284 MMG5_ADD_MEM(mesh,nsols*sizeof(MMG5_Sol),"solutions array",
285 printf(" Exit program.\n");
286 return -1);
287 MMG5_SAFE_CALLOC(*sol,nsols,MMG5_Sol,return -1);
288
289 // Mesh alloc and transfer of the mesh from dataset toward the MMG5 Mesh Sol
290 ier = MMGS_loadVtkMesh_part2(mesh,sol,&dataset,ptMeditRef,eltMeditRef,nsols);
291
292 return ier;
293#endif
294}
295
297
298#ifndef USE_VTK
299
300 fprintf(stderr," ** VTK library not found. Unavailable file format.\n");
301 return -1;
302
303#else
304
305 return MMG5_saveVtkMesh<vtkUnstructuredGrid,vtkXMLUnstructuredGridWriter,
306 vtkXMLPUnstructuredGridWriter>(mesh,&sol,filename,1,1);
307
308#endif
309}
310
312
313#ifndef USE_VTK
314
315 fprintf(stderr," ** VTK library not found. Unavailable file format.\n");
316 return -1;
317
318#else
319
320 return MMG5_saveVtkMesh<vtkUnstructuredGrid,vtkXMLUnstructuredGridWriter,
321 vtkXMLPUnstructuredGridWriter>(mesh,sol,filename,0,1);
322
323#endif
324}
325
327
328#ifndef USE_VTK
329
330 fprintf(stderr," ** VTK library not found. Unavailable file format.\n");
331 return -1;
332
333#else
334
335 return MMG5_saveVtkMesh<vtkUnstructuredGrid,vtkDataSetWriter,
336 vtkPDataSetWriter>(mesh,&sol,filename,1,0);
337
338#endif
339}
340
342
343#ifndef USE_VTK
344
345 fprintf(stderr," ** VTK library not found. Unavailable file format.\n");
346 return -1;
347
348#else
349
350 return MMG5_saveVtkMesh<vtkUnstructuredGrid,vtkDataSetWriter,
351 vtkPDataSetWriter>(mesh,sol,filename,0,0);
352
353#endif
354}
355
357
358#ifndef USE_VTK
359
360 fprintf(stderr," ** VTK library not found. Unavailable file format.\n");
361 return -1;
362
363#else
364
365 return MMG5_saveVtkMesh<vtkPolyData,vtkXMLPolyDataWriter,
366 vtkXMLPPolyDataWriter>(mesh,&sol,filename,1,1);
367
368#endif
369}
370
372
373#ifndef USE_VTK
374
375 fprintf(stderr," ** VTK library not found. Unavailable file format.\n");
376 return -1;
377
378#else
379
380 return MMG5_saveVtkMesh<vtkPolyData,vtkXMLPolyDataWriter,
381 vtkXMLPPolyDataWriter>(mesh,sol,filename,0,1);
382
383#endif
384}
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:2661
int MMGS_saveVtuMesh_and_allData(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename)
Definition: inoutcpp_s.cpp:311
int MMGS_loadVtpMesh(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
Definition: inoutcpp_s.cpp:74
int MMGS_saveVtkMesh_and_allData(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename)
Definition: inoutcpp_s.cpp:341
int MMGS_loadVtuMesh(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
Definition: inoutcpp_s.cpp:222
int MMGS_saveVtuMesh(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
Definition: inoutcpp_s.cpp:296
static int MMGS_loadVtkMesh_part2(MMG5_pMesh mesh, MMG5_pSol *sol, vtkDataSet **dataset, int8_t ptMeditRef, int8_t eltMeditRef, MMG5_int nsols)
Input / Output Functions that needs cpp features.
Definition: inoutcpp_s.cpp:40
int MMGS_loadVtkMesh(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
Definition: inoutcpp_s.cpp:148
int MMGS_saveVtkMesh(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
Definition: inoutcpp_s.cpp:326
int MMGS_loadVtpMesh_and_allData(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename)
Definition: inoutcpp_s.cpp:115
int MMGS_loadVtuMesh_and_allData(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename)
Definition: inoutcpp_s.cpp:263
int MMGS_loadVtkMesh_and_allData(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename)
Definition: inoutcpp_s.cpp:189
int MMGS_saveVtpMesh(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
Definition: inoutcpp_s.cpp:356
int MMGS_saveVtpMesh_and_allData(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename)
Definition: inoutcpp_s.cpp:371
API headers 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:605
MMG5_int ntmax
Definition: libmmgtypes.h:612
MMG5_int ne
Definition: libmmgtypes.h:612
MMG5_int npmax
Definition: libmmgtypes.h:612
MMG5_int nt
Definition: libmmgtypes.h:612
MMG5_int np
Definition: libmmgtypes.h:612
MMG5_int nprism
Definition: libmmgtypes.h:613
int MMG5_loadVtkMesh_part2(MMG5_pMesh mesh, MMG5_pSol *sol, vtkDataSet **dataset, int8_t ptMeditRef, int8_t eltMeditRef, int nsols)
Definition: vtkparser.cpp:346
int MMG5_loadVtpMesh_part1(MMG5_pMesh mesh, const char *filename, vtkDataSet **dataset, int8_t *ptMeditRef, int8_t *eltMeditRef, int *nsols, int8_t *metricData)
Definition: vtkparser.cpp:230
int MMG5_loadVtuMesh_part1(MMG5_pMesh mesh, const char *filename, vtkDataSet **dataset, int8_t *ptMeditRef, int8_t *eltMeditRef, int *nsols, int8_t *metricData)
Definition: vtkparser.cpp:308
int MMG5_loadVtkMesh_part1(MMG5_pMesh mesh, const char *filename, vtkDataSet **dataset, int8_t *ptMeditRef, int8_t *eltMeditRef, int *nsols, int8_t *metricData)
Definition: vtkparser.cpp:269
int MMG5_saveVtkMesh(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename, int metricData, int binary)
Definition: vtkparser.hpp:516