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 int8_t metricData, int8_t lsData) {
45 int ier;
46
47 if ( !MMG2D_zaldy(mesh) ) {
48 return 0;
49 }
50
51 if ( mesh->ne || mesh->nprism ) {
52 fprintf(stderr,"\n ## Error: %s: Input mesh must be a two-dimensional mesh.\n",
53 __func__);
54 return -1;
55 }
56
57 if ( !mesh->nt ) {
58 fprintf(stderr," ** MISSING DATA.\n");
59 fprintf(stderr," Check that your mesh contains triangles.\n");
60 fprintf(stderr," Exit program.\n");
61 return -1;
62 }
63
64 if (mesh->npmax < mesh->np || mesh->ntmax < mesh->nt ) {
65 return -1;
66 }
67
68 ier = MMG5_loadVtkMesh_part2(mesh,sol,dataset,ptMeditRef,eltMeditRef,nsols,metricData,lsData);
69
70 if ( ier < 1 ) {
71 fprintf(stderr," ** ERROR WHEN PARSING THE INPUT FILE\n");
72 return ier;
73 }
74
75 /* Mark all points as used in case of mesh generation and check the
76 * z-componant */
77 if ( !MMG2D_2dMeshCheck(mesh) ) return -1;
78
79 return ier;
80}
81
82#endif
83
85#ifndef USE_VTK
86
87 fprintf(stderr," ** VTK library not founded. Unavailable file format.\n");
88 return -1;
89
90#else
91 int ier,nsols;
92 int8_t ptMeditRef,eltMeditRef,metricData,lsData;
93 vtkDataSet *dataset;
94 MMG5_pSol allSol[2];
95
96 mesh->dim = 2;
97
98 ier = MMG5_loadVtpMesh_part1(mesh,filename,&dataset,&ptMeditRef,&eltMeditRef,
99 &nsols,&metricData,&lsData);
100 if ( ier < 1 ) return ier;
101
102 /* Check data fields */
103 if ( nsols > (metricData+lsData) ) {
104 fprintf(stderr,"Error: %d UNEXPECTED DATA FIELD(S)\n",nsols-metricData-lsData);
105 return -1;
106 }
107
108 // Mesh alloc and transfer of the mesh from dataset toward the MMG5 Mesh Sol
109 allSol[0] = met;
110 allSol[1] = sol;
111 ier = MMG2D_loadVtkMesh_part2(mesh,allSol,&dataset,ptMeditRef,eltMeditRef,nsols,metricData,lsData);
112 if ( ier < 1 ) {
113 fprintf(stderr," ** ERROR WHEN PARSING THE INPUT FILE\n");
114 return ier;
115 }
116
117 if ( sol ) {
118 /* Check the metric type */
119 ier = MMG5_chkMetricType(mesh,&sol->type,&sol->entities,NULL);
120 if ( ier <1 ) {
121 fprintf(stderr," ** UNEXPECTED METRIC TYPE\n");
122 return ier;
123 }
124 }
125
126 return ier;
127#endif
128}
129
131#ifndef USE_VTK
132
133 fprintf(stderr," ** VTK library not founded. Unavailable file format.\n");
134 return -1;
135
136#else
137 int ier,nsols;
138 int8_t ptMeditRef,eltMeditRef,metricData,lsData;
139 vtkDataSet *dataset;
140 MMG5_pSol allSol[2];
141
142 mesh->dim = 2;
143
144 ier = MMG5_loadVtpMesh_part1(mesh,filename,&dataset,&ptMeditRef,&eltMeditRef,
145 &nsols,&metricData,&lsData);
146 if ( ier < 1 ) return ier;
147
148 mesh->nsols = nsols;
149 if ( *sol ) MMG5_DEL_MEM(mesh,*sol);
150
151 MMG5_ADD_MEM(mesh,nsols*sizeof(MMG5_Sol),"solutions array",
152 printf(" Exit program.\n");
153 return -1);
154 MMG5_SAFE_CALLOC(*sol,nsols,MMG5_Sol,return -1);
155
156 // Mesh alloc and transfer of the mesh from dataset toward the MMG5 Mesh Sol
157 allSol[0] = NULL;
158 allSol[1] = *sol;
159 ier = MMG2D_loadVtkMesh_part2(mesh,allSol,&dataset,ptMeditRef,eltMeditRef,nsols,metricData,0);
160
161 return ier;
162#endif
163}
164
166#ifndef USE_VTK
167
168 fprintf(stderr," ** VTK library not founded. Unavailable file format.\n");
169 return -1;
170
171#else
172 int ier,nsols;
173 int8_t ptMeditRef,eltMeditRef,metricData,lsData;
174 vtkDataSet *dataset;
175 MMG5_pSol allSol[2];
176
177 mesh->dim = 2;
178
179 ier = MMG5_loadVtkMesh_part1(mesh,filename,&dataset,&ptMeditRef,&eltMeditRef,
180 &nsols,&metricData,&lsData);
181 if ( ier < 1 ) return ier;
182
183 /* Check data fields */
184 if ( nsols > (metricData+lsData) ) {
185 fprintf(stderr,"Error: %d UNEXPECTED DATA FIELD(S)\n",nsols-metricData-lsData);
186 return -1;
187 }
188
189 // Mesh alloc and transfer of the mesh from dataset toward the MMG5 Mesh Sol
190 allSol[0] = met;
191 allSol[1] = sol;
192 ier = MMG2D_loadVtkMesh_part2(mesh,allSol,&dataset,ptMeditRef,eltMeditRef,nsols,metricData,lsData);
193 if ( ier < 1 ) {
194 fprintf(stderr," ** ERROR WHEN PARSING THE INPUT FILE\n");
195 return ier;
196 }
197
198 if ( sol ) {
199 /* Check the metric type */
200 ier = MMG5_chkMetricType(mesh,&sol->type,&sol->entities,NULL);
201 if ( ier <1 ) {
202 fprintf(stderr," ** UNEXPECTED METRIC TYPE\n");
203 return ier;
204 }
205 }
206
207 return ier;
208#endif
209}
210
212#ifndef USE_VTK
213
214 fprintf(stderr," ** VTK library not founded. Unavailable file format.\n");
215 return -1;
216
217#else
218 int ier,nsols;
219 int8_t ptMeditRef,eltMeditRef,metricData,lsData;
220 vtkDataSet *dataset;
221 MMG5_pSol allSol[2];
222
223 mesh->dim = 2;
224
225 ier = MMG5_loadVtkMesh_part1(mesh,filename,&dataset,&ptMeditRef,&eltMeditRef,
226 &nsols,&metricData,&lsData);
227 if ( ier < 1 ) return ier;
228
229 mesh->nsols = nsols;
230 if ( *sol ) MMG5_DEL_MEM(mesh,*sol);
231
232 MMG5_ADD_MEM(mesh,nsols*sizeof(MMG5_Sol),"solutions array",
233 printf(" Exit program.\n");
234 return -1);
235 MMG5_SAFE_CALLOC(*sol,nsols,MMG5_Sol,return -1);
236
237 // Mesh alloc and transfer of the mesh from dataset toward the MMG5 Mesh Sol
238 allSol[0] = NULL;
239 allSol[1] = *sol;
240 ier = MMG2D_loadVtkMesh_part2(mesh,allSol,&dataset,ptMeditRef,eltMeditRef,nsols,metricData,0);
241
242 return ier;
243#endif
244}
245
247#ifndef USE_VTK
248
249 fprintf(stderr," ** VTK library not founded. Unavailable file format.\n");
250 return -1;
251
252#else
253 int ier,nsols;
254 int8_t ptMeditRef,eltMeditRef,metricData,lsData;
255 vtkDataSet *dataset;
256 MMG5_pSol allSol[2];
257
258 mesh->dim = 2;
259
260 ier = MMG5_loadVtuMesh_part1(mesh,filename,&dataset,&ptMeditRef,&eltMeditRef,
261 &nsols,&metricData,&lsData);
262 if ( ier < 1 ) return ier;
263
264 /* Check data fields */
265 if ( nsols > (metricData+lsData) ) {
266 fprintf(stderr,"Error: %d UNEXPECTED DATA FIELD(S)\n",nsols-metricData-lsData);
267 return -1;
268 }
269
270 // Mesh alloc and transfer of the mesh from dataset toward the MMG5 Mesh Sol
271 allSol[0] = met;
272 allSol[1] = sol;
273 ier = MMG2D_loadVtkMesh_part2(mesh,allSol,&dataset,ptMeditRef,eltMeditRef,nsols,metricData,lsData);
274 if ( ier < 1 ) {
275 fprintf(stderr," ** ERROR WHEN PARSING THE INPUT FILE\n");
276 return ier;
277 }
278
279 if ( sol ) {
280 /* Check the metric type */
281 ier = MMG5_chkMetricType(mesh,&sol->type,&sol->entities,NULL);
282 if ( ier <1 ) {
283 fprintf(stderr," ** UNEXPECTED METRIC TYPE\n");
284 return ier;
285 }
286 }
287
288 return ier;
289#endif
290}
291
293
294#ifndef USE_VTK
295
296 fprintf(stderr," ** VTK library not founded. Unavailable file format.\n");
297 return -1;
298
299#else
300 int ier,nsols;
301 int8_t ptMeditRef,eltMeditRef,metricData,lsData;
302 vtkDataSet *dataset;
303 MMG5_pSol allSol[2];
304
305 mesh->dim = 2;
306
307 ier = MMG5_loadVtuMesh_part1(mesh,filename,&dataset,&ptMeditRef,&eltMeditRef,
308 &nsols,&metricData,&lsData);
309 if ( ier < 1 ) return ier;
310
311 /* Check element count */
312 if ( nsols>1 ) {
313 fprintf(stderr,"Error: SEVERAL SOLUTIONS FOUND (%d)\n",nsols);
314 return -1;
315 }
316
317 mesh->nsols = nsols;
318 if ( *sol ) MMG5_DEL_MEM(mesh,*sol);
319
320 MMG5_ADD_MEM(mesh,nsols*sizeof(MMG5_Sol),"solutions array",
321 printf(" Exit program.\n");
322 return -1);
323 MMG5_SAFE_CALLOC(*sol,nsols,MMG5_Sol,return -1);
324
325 // Mesh alloc and transfer of the mesh from dataset toward the MMG5 Mesh Sol
326 allSol[0] = NULL;
327 allSol[1] = *sol;
328 ier = MMG2D_loadVtkMesh_part2(mesh,allSol,&dataset,ptMeditRef,eltMeditRef,nsols,metricData,0);
329
330 return ier;
331#endif
332}
333
335
336#ifndef USE_VTK
337
338 fprintf(stderr," ** VTK library not found. Unavailable file format.\n");
339 return -1;
340
341#else
342
343 return MMG5_saveVtkMesh<vtkUnstructuredGrid,vtkXMLUnstructuredGridWriter,
344 vtkXMLPUnstructuredGridWriter>(mesh,&sol,filename,1,1);
345
346#endif
347}
348
350
351 MMG5_pSol allSol[2];
352
353#ifndef USE_VTK
354
355 fprintf(stderr," ** VTK library not found. Unavailable file format.\n");
356 return -1;
357
358#else
359
360 allSol[0] = NULL;
361 allSol[1] = *sol;
362
363 return MMG5_saveVtkMesh<vtkUnstructuredGrid,vtkXMLUnstructuredGridWriter,
364 vtkXMLPUnstructuredGridWriter>(mesh,allSol,filename,0,1);
365
366#endif
367}
368
370
371#ifndef USE_VTK
372
373 fprintf(stderr," ** VTK library not found. Unavailable file format.\n");
374 return -1;
375
376#else
377
378 return MMG5_saveVtkMesh<vtkUnstructuredGrid,vtkDataSetWriter,
379 vtkPDataSetWriter>(mesh,&sol,filename,1,0);
380
381#endif
382}
383
385
386 MMG5_pSol allSol[2];
387
388#ifndef USE_VTK
389
390 fprintf(stderr," ** VTK library not found. Unavailable file format.\n");
391 return -1;
392
393#else
394
395 allSol[0] = NULL;
396 allSol[1] = *sol;
397
398 return MMG5_saveVtkMesh<vtkUnstructuredGrid,vtkDataSetWriter,
399 vtkPDataSetWriter>(mesh,allSol,filename,0,0);
400
401#endif
402}
403
405
406#ifndef USE_VTK
407
408 fprintf(stderr," ** VTK library not found. Unavailable file format.\n");
409 return -1;
410
411#else
412
413 return MMG5_saveVtkMesh<vtkPolyData,vtkXMLPolyDataWriter,
414 vtkXMLPPolyDataWriter>(mesh,&sol,filename,1,1);
415
416#endif
417
418}
419
421
422 MMG5_pSol allSol[2];
423
424#ifndef USE_VTK
425
426 fprintf(stderr," ** VTK library not found. Unavailable file format.\n");
427 return -1;
428
429#else
430
431 allSol[0] = NULL;
432 allSol[1] = *sol;
433
434 return MMG5_saveVtkMesh<vtkPolyData,vtkXMLPolyDataWriter,
435 vtkXMLPPolyDataWriter>(mesh,allSol,filename,0,1);
436
437#endif
438
439}
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 MMG2D_2dMeshCheck(MMG5_pMesh mesh)
Definition: inout_2d.c:679
int MMG2D_loadVtpMesh(MMG5_pMesh mesh, MMG5_pSol met, MMG5_pSol sol, const char *filename)
Load a mesh and possibly a solution in VTP (VTK) format from file.
Definition: inoutcpp_2d.cpp:84
static int MMG2D_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_2d.cpp:41
int MMG2D_loadVtkMesh_and_allData(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename)
Load a mesh and multiple solutions in VTK format from file.
int MMG2D_loadVtkMesh(MMG5_pMesh mesh, MMG5_pSol met, MMG5_pSol sol, const char *filename)
Load a mesh and possibly data in VTK format from file.
int MMG2D_loadVtpMesh_and_allData(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename)
Load a mesh and multiple solutions in VTP (VTK) format from file.
int MMG2D_saveVtpMesh_and_allData(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename)
Save a mesh and multiple data fields in VTP format.
int MMG2D_saveVtkMesh_and_allData(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename)
Save a mesh and multiple data fields in VTK format.
int MMG2D_saveVtkMesh(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
Save a mesh and optionally one solution in VTK format.
int MMG2D_saveVtuMesh(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
Save a mesh and optionally one data field in VTU format.
int MMG2D_saveVtpMesh(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
Save a mesh and optionally one data field in VTP format.
int MMG2D_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.
int MMG2D_saveVtuMesh_and_allData(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename)
Save a mesh and multiple data fields in VTU format.
int MMG2D_loadVtuMesh_and_allData(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename)
Load a mesh and multiple solutions in VTU (VTK) format from file.
API headers and documentation for the mmg2d library.
int MMG2D_zaldy(MMG5_pMesh mesh)
Definition: zaldy_2d.c:307
#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