Mmg
Simplicial remeshers (mesh adaptation, isovalue discretization, lagrangian movement)
inoutcpp_2d.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
36#include "libmmg2d.h"
37#include "libmmg2d_private.h"
38
39#ifdef USE_VTK
40
42 vtkDataSet **dataset, int8_t ptMeditRef,
43 int8_t eltMeditRef,int nsols ) {
44 int ier;
45
46 if ( !MMG2D_zaldy(mesh) ) {
47 return 0;
48 }
49
50 if ( mesh->ne || mesh->nprism ) {
51 fprintf(stderr,"\n ## Error: %s: Input mesh must be a two-dimensional mesh.\n",
52 __func__);
53 return -1;
54 }
55
56 if ( !mesh->nt ) {
57 fprintf(stderr," ** MISSING DATA.\n");
58 fprintf(stderr," Check that your mesh contains triangles.\n");
59 fprintf(stderr," Exit program.\n");
60 return -1;
61 }
62
63 if (mesh->npmax < mesh->np || mesh->ntmax < mesh->nt ) {
64 return -1;
65 }
66
67 ier = MMG5_loadVtkMesh_part2(mesh,sol,dataset,ptMeditRef,eltMeditRef,nsols);
68
69 if ( ier < 1 ) {
70 fprintf(stderr," ** ERROR WHEN PARSING THE INPUT FILE\n");
71 return ier;
72 }
73
74 /* Mark all points as used in case of mesh generation and check the
75 * z-componant */
76 if ( !MMG2D_2dMeshCheck(mesh) ) return -1;
77
78 return ier;
79}
80
81#endif
82
84#ifndef USE_VTK
85
86 fprintf(stderr," ** VTK library not founded. Unavailable file format.\n");
87 return -1;
88
89#else
90 int ier,nsols;
91 int8_t ptMeditRef,eltMeditRef,metricData;
92 vtkDataSet *dataset;
93
94 mesh->dim = 2;
95
96 ier = MMG5_loadVtpMesh_part1(mesh,filename,&dataset,&ptMeditRef,&eltMeditRef,
97 &nsols,&metricData);
98 if ( ier < 1 ) return ier;
99
100 /* Check data fields */
101 if ( nsols > metricData ) {
102 fprintf(stderr,"Error: %d UNEXPECTED DATA FIELD(S)\n",nsols);
103 return -1;
104 }
105
106 // Mesh alloc and transfer of the mesh from dataset toward the MMG5 Mesh Sol
107 ier = MMG2D_loadVtkMesh_part2(mesh,&sol,&dataset,ptMeditRef,eltMeditRef,nsols);
108 if ( ier < 1 ) {
109 fprintf(stderr," ** ERROR WHEN PARSING THE INPUT FILE\n");
110 return ier;
111 }
112
113 if ( sol ) {
114 /* Check the metric type */
115 ier = MMG5_chkMetricType(mesh,&sol->type,&sol->entities,NULL);
116 if ( ier <1 ) {
117 fprintf(stderr," ** UNEXPECTED METRIC TYPE\n");
118 return ier;
119 }
120 }
121
122 return ier;
123#endif
124}
125
127#ifndef USE_VTK
128
129 fprintf(stderr," ** VTK library not founded. Unavailable file format.\n");
130 return -1;
131
132#else
133 int ier,nsols;
134 int8_t ptMeditRef,eltMeditRef,metricData;
135 vtkDataSet *dataset;
136
137 mesh->dim = 2;
138
139 ier = MMG5_loadVtpMesh_part1(mesh,filename,&dataset,&ptMeditRef,&eltMeditRef,
140 &nsols,&metricData);
141 if ( ier < 1 ) return ier;
142
143 mesh->nsols = nsols;
144 if ( *sol ) MMG5_DEL_MEM(mesh,*sol);
145
146 MMG5_ADD_MEM(mesh,nsols*sizeof(MMG5_Sol),"solutions array",
147 printf(" Exit program.\n");
148 return -1);
149 MMG5_SAFE_CALLOC(*sol,nsols,MMG5_Sol,return -1);
150
151 // Mesh alloc and transfer of the mesh from dataset toward the MMG5 Mesh Sol
152 ier = MMG2D_loadVtkMesh_part2(mesh,sol,&dataset,ptMeditRef,eltMeditRef,nsols);
153
154 return ier;
155#endif
156}
157
159#ifndef USE_VTK
160
161 fprintf(stderr," ** VTK library not founded. Unavailable file format.\n");
162 return -1;
163
164#else
165 int ier,nsols;
166 int8_t ptMeditRef,eltMeditRef,metricData;
167 vtkDataSet *dataset;
168
169 mesh->dim = 2;
170
171 ier = MMG5_loadVtkMesh_part1(mesh,filename,&dataset,&ptMeditRef,&eltMeditRef,
172 &nsols,&metricData);
173 if ( ier < 1 ) return ier;
174
175 /* Check data fields */
176 if ( nsols > metricData ) {
177 fprintf(stderr,"Error: %d UNEXPECTED DATA FIELD(S)\n",nsols);
178 return -1;
179 }
180
181 // Mesh alloc and transfer of the mesh from dataset toward the MMG5 Mesh Sol
182 ier = MMG2D_loadVtkMesh_part2(mesh,&sol,&dataset,ptMeditRef,eltMeditRef,nsols);
183 if ( ier < 1 ) {
184 fprintf(stderr," ** ERROR WHEN PARSING THE INPUT FILE\n");
185 return ier;
186 }
187
188 if ( sol ) {
189 /* Check the metric type */
190 ier = MMG5_chkMetricType(mesh,&sol->type,&sol->entities,NULL);
191 if ( ier <1 ) {
192 fprintf(stderr," ** UNEXPECTED METRIC TYPE\n");
193 return ier;
194 }
195 }
196
197 return ier;
198#endif
199}
200
202#ifndef USE_VTK
203
204 fprintf(stderr," ** VTK library not founded. Unavailable file format.\n");
205 return -1;
206
207#else
208 int ier,nsols;
209 int8_t ptMeditRef,eltMeditRef,metricData;
210 vtkDataSet *dataset;
211
212 mesh->dim = 2;
213
214 ier = MMG5_loadVtkMesh_part1(mesh,filename,&dataset,&ptMeditRef,&eltMeditRef,
215 &nsols,&metricData);
216 if ( ier < 1 ) return ier;
217
218 mesh->nsols = nsols;
219 if ( *sol ) MMG5_DEL_MEM(mesh,*sol);
220
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 ier = MMG2D_loadVtkMesh_part2(mesh,sol,&dataset,ptMeditRef,eltMeditRef,nsols);
228
229 return ier;
230#endif
231}
232
234#ifndef USE_VTK
235
236 fprintf(stderr," ** VTK library not founded. Unavailable file format.\n");
237 return -1;
238
239#else
240 int ier,nsols;
241 int8_t ptMeditRef,eltMeditRef,metricData;
242 vtkDataSet *dataset;
243
244 mesh->dim = 2;
245
246 ier = MMG5_loadVtuMesh_part1(mesh,filename,&dataset,&ptMeditRef,&eltMeditRef,
247 &nsols,&metricData);
248 if ( ier < 1 ) return ier;
249
250 /* Check data fields */
251 if ( nsols > metricData ) {
252 fprintf(stderr,"Error: %d UNEXPECTED DATA FIELD(S)\n",nsols);
253 return -1;
254 }
255
256 // Mesh alloc and transfer of the mesh from dataset toward the MMG5 Mesh Sol
257 ier = MMG2D_loadVtkMesh_part2(mesh,&sol,&dataset,ptMeditRef,eltMeditRef,nsols);
258 if ( ier < 1 ) {
259 fprintf(stderr," ** ERROR WHEN PARSING THE INPUT FILE\n");
260 return ier;
261 }
262
263 if ( sol ) {
264 /* Check the metric type */
265 ier = MMG5_chkMetricType(mesh,&sol->type,&sol->entities,NULL);
266 if ( ier <1 ) {
267 fprintf(stderr," ** UNEXPECTED METRIC TYPE\n");
268 return ier;
269 }
270 }
271
272 return ier;
273#endif
274}
275
277
278#ifndef USE_VTK
279
280 fprintf(stderr," ** VTK library not founded. Unavailable file format.\n");
281 return -1;
282
283#else
284 int ier,nsols;
285 int8_t ptMeditRef,eltMeditRef,metricData;
286 vtkDataSet *dataset;
287
288 mesh->dim = 2;
289
290 ier = MMG5_loadVtuMesh_part1(mesh,filename,&dataset,&ptMeditRef,&eltMeditRef,
291 &nsols,&metricData);
292 if ( ier < 1 ) return ier;
293
294 /* Check element count */
295 if ( nsols>1 ) {
296 fprintf(stderr,"Error: SEVERAL SOLUTIONS FOUND (%d)\n",nsols);
297 return -1;
298 }
299
300 mesh->nsols = nsols;
301 if ( *sol ) MMG5_DEL_MEM(mesh,*sol);
302
303 MMG5_ADD_MEM(mesh,nsols*sizeof(MMG5_Sol),"solutions array",
304 printf(" Exit program.\n");
305 return -1);
306 MMG5_SAFE_CALLOC(*sol,nsols,MMG5_Sol,return -1);
307
308 // Mesh alloc and transfer of the mesh from dataset toward the MMG5 Mesh Sol
309 ier = MMG2D_loadVtkMesh_part2(mesh,sol,&dataset,ptMeditRef,eltMeditRef,nsols);
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#ifndef USE_VTK
333
334 fprintf(stderr," ** VTK library not found. Unavailable file format.\n");
335 return -1;
336
337#else
338
339 return MMG5_saveVtkMesh<vtkUnstructuredGrid,vtkXMLUnstructuredGridWriter,
340 vtkXMLPUnstructuredGridWriter>(mesh,sol,filename,0,1);
341
342#endif
343}
344
346
347#ifndef USE_VTK
348
349 fprintf(stderr," ** VTK library not found. Unavailable file format.\n");
350 return -1;
351
352#else
353
354 return MMG5_saveVtkMesh<vtkUnstructuredGrid,vtkDataSetWriter,
355 vtkPDataSetWriter>(mesh,&sol,filename,1,0);
356
357#endif
358}
359
361
362#ifndef USE_VTK
363
364 fprintf(stderr," ** VTK library not found. Unavailable file format.\n");
365 return -1;
366
367#else
368
369 return MMG5_saveVtkMesh<vtkUnstructuredGrid,vtkDataSetWriter,
370 vtkPDataSetWriter>(mesh,sol,filename,0,0);
371
372#endif
373}
374
376
377#ifndef USE_VTK
378
379 fprintf(stderr," ** VTK library not found. Unavailable file format.\n");
380 return -1;
381
382#else
383
384 return MMG5_saveVtkMesh<vtkPolyData,vtkXMLPolyDataWriter,
385 vtkXMLPPolyDataWriter>(mesh,&sol,filename,1,1);
386
387#endif
388
389}
390
392
393#ifndef USE_VTK
394
395 fprintf(stderr," ** VTK library not found. Unavailable file format.\n");
396 return -1;
397
398#else
399
400 return MMG5_saveVtkMesh<vtkPolyData,vtkXMLPolyDataWriter,
401 vtkXMLPPolyDataWriter>(mesh,sol,filename,0,1);
402
403#endif
404
405}
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 MMG2D_2dMeshCheck(MMG5_pMesh mesh)
Definition: inout_2d.c:674
static int MMG2D_loadVtkMesh_part2(MMG5_pMesh mesh, MMG5_pSol *sol, vtkDataSet **dataset, int8_t ptMeditRef, int8_t eltMeditRef, int nsols)
Input / Output Functions that needs cpp features.
Definition: inoutcpp_2d.cpp:41
int MMG2D_loadVtkMesh_and_allData(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename)
int MMG2D_loadVtkMesh(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
int MMG2D_loadVtpMesh_and_allData(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename)
int MMG2D_saveVtpMesh_and_allData(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename)
int MMG2D_saveVtkMesh_and_allData(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename)
int MMG2D_saveVtkMesh(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
int MMG2D_saveVtuMesh(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
int MMG2D_loadVtuMesh(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
int MMG2D_saveVtpMesh(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
int MMG2D_saveVtuMesh_and_allData(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename)
int MMG2D_loadVtpMesh(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
Definition: inoutcpp_2d.cpp:83
int MMG2D_loadVtuMesh_and_allData(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename)
API headers for the mmg2d library.
int MMG2D_zaldy(MMG5_pMesh mesh)
Definition: zaldy_2d.c:303
#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