Mmg
Simplicial remeshers (mesh adaptation, isovalue discretization, lagrangian movement)
vtkparser.hpp
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/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
24#ifndef VTKPARSER_HPP
25#define VTKPARSER_HPP
26
27#ifdef USE_VTK
28
29#include "mmgcommon_private.h"
30
31#include <vtkSmartPointer.h>
32#include <vtkXMLReader.h>
33#include <vtkXMLWriter.h>
34#include <vtkXMLUnstructuredGridReader.h>
35#include <vtkXMLUnstructuredGridWriter.h>
36#include <vtkXMLPUnstructuredGridWriter.h>
37#include <vtkXMLPolyDataReader.h>
38#include <vtkXMLPolyDataWriter.h>
39#include <vtkXMLPPolyDataWriter.h>
40#include <vtkDataSetReader.h>
41#include <vtkDataSetWriter.h>
42#include <vtkPDataSetWriter.h>
43#include <vtkDataSet.h>
44#include <vtkUnstructuredGrid.h>
45#include <vtkPolyData.h>
46#include <vtkStructuredGrid.h>
47#include <vtkPointData.h>
48#include <vtkCellData.h>
49#include <vtkFieldData.h>
50#include <vtkCellTypes.h>
51#include <vtkDataArray.h>
52#include <vtkFloatArray.h>
53#include <vtkDoubleArray.h>
54
55#include <vtkLine.h>
56#include <vtkTriangle.h>
57#include <vtkQuad.h>
58#include <vtkTetra.h>
59#include <vtkWedge.h>
60#include <vtkCellArray.h>
61#include <typeinfo>
62
63int MMG5_loadVtpMesh_part1(MMG5_pMesh,const char*,vtkDataSet**,int8_t*,int8_t*,int*,int8_t*);
64int MMG5_loadVtuMesh_part1(MMG5_pMesh,const char*,vtkDataSet**,int8_t*,int8_t*,int*,int8_t*);
65int MMG5_loadVtkMesh_part1(MMG5_pMesh,const char*,vtkDataSet**,int8_t*,int8_t*,int*,int8_t*);
66
67int MMG5_loadVtkMesh_part2(MMG5_pMesh,MMG5_pSol*,vtkDataSet**,int8_t,int8_t,int);
68
76static void MMG5_internal_VTKSetLine(vtkSmartPointer<vtkPolyData> d,vtkSmartPointer<vtkCellArray> *ca) {
77
78 d->SetLines(*ca);
79
80 *ca = NULL;
81 *ca = vtkSmartPointer< vtkCellArray >::New();
82
83}
84
90static void MMG5_internal_VTKSetLine(vtkSmartPointer<vtkUnstructuredGrid> d,vtkSmartPointer<vtkCellArray> *ca) {
91
92}
93
94
101static void MMG5_internal_VTKSetCells(int * t,vtkSmartPointer<vtkPolyData> d,vtkSmartPointer<vtkCellArray> ca) {
102 d->SetPolys(ca);
103}
104
111static void MMG5_internal_VTKSetCells(int * t,vtkSmartPointer<vtkUnstructuredGrid> d,vtkSmartPointer<vtkCellArray> ca) {
112 d->SetCells(t,ca);
113}
114
120static void MMG5_internal_VTKbinary(vtkXMLUnstructuredGridWriter *w, int binary) {
121 if ( binary ) {
122 w->SetDataModeToBinary();
123 }
124 else {
125 w->SetDataModeToAscii();
126 }
127}
128
134static void MMG5_internal_VTKbinary(vtkXMLPolyDataWriter *w, int binary) {
135 if ( binary ) {
136 w->SetDataModeToBinary();
137 }
138 else {
139 w->SetDataModeToAscii();
140 }
141}
142
148static void MMG5_internal_VTKbinary(vtkDataSetWriter *w, int binary) {
149 return;
150}
151
169template <class T, class TWriter, class PWriter>
171 const char *mfilename,
172 int metricData,int binary,
173 int npart, int myid,int master) {
174 int hasPointRef = 0, hasCellRef = 0;
175
176 // Transfer points from Mmg to VTK
177 vtkSmartPointer<vtkPoints> points = vtkSmartPointer<vtkPoints>::New();
178
179 MMG5_int np = 0;
180 for ( MMG5_int k=1; k<=mesh->np; ++k ) {
181 MMG5_pPoint ppt = &mesh->point[k];
182 if ( !MG_VOK(ppt) ) continue;
183
184 ppt->tmp = np++;
185
186 if ( mesh->dim == 2 ) ppt->c[2] = 0.;
187 if ( ppt->ref ) hasPointRef = 1;
188
189 points->InsertNextPoint(ppt->c[0], ppt->c[1], ppt->c[2]);
190 }
191
192
193 auto dataset = vtkSmartPointer<T> :: New();
194 dataset->SetPoints(points);
195
196
197 // VTK types
198 vtkSmartPointer<vtkCellArray> cellArray = vtkSmartPointer<vtkCellArray> ::New();
199 vtkSmartPointer<vtkWedge> wedge = vtkSmartPointer<vtkWedge> ::New();
200 vtkSmartPointer<vtkTetra> tetra = vtkSmartPointer<vtkTetra> ::New();
201 vtkSmartPointer<vtkTriangle> tria = vtkSmartPointer<vtkTriangle> ::New();
202 vtkSmartPointer<vtkQuad> quadra = vtkSmartPointer<vtkQuad> ::New();
203 vtkSmartPointer<vtkLine> edge = vtkSmartPointer<vtkLine> ::New();
204
205 MMG5_int nc = mesh->na+mesh->nt+mesh->nquad+mesh->ne+mesh->nprism;
206 MMG5_int ic = 0;
207
208 int* types = NULL;
209 MMG5_SAFE_MALLOC ( types, nc, int, return 0 );
210
211 // transfer edges from Mmg to VTK
212 for ( MMG5_int k=1; k<=mesh->na; ++k ) {
213 MMG5_pEdge pa = &mesh->edge[k];
214 if ( !pa || !pa->a ) continue;
215
216 if ( pa->ref ) hasCellRef = 1;
217
218 edge->GetPointIds()->SetId(0, mesh->point[pa->a].tmp);
219 edge->GetPointIds()->SetId(1, mesh->point[pa->b].tmp);
220
221 cellArray->InsertNextCell(edge);
222 types[ic++] = VTK_LINE;
223 }
224
225 MMG5_internal_VTKSetLine(dataset, &cellArray);
226
227 // Transfer trias from Mmg to VTK
228 for ( MMG5_int k=1; k<=mesh->nt; ++k ) {
229 MMG5_pTria ptt = &mesh->tria[k];
230 if ( !MG_EOK(ptt) ) continue;
231
232 if ( ptt->ref ) hasCellRef = 1;
233
234 for ( int i=0; i<3; ++i ) {
235 tria->GetPointIds()->SetId(i, mesh->point[ptt->v[i]].tmp );
236 }
237 cellArray->InsertNextCell(tria);
238 types[ic++] = VTK_TRIANGLE;
239 }
240
241 // Transfer quads from Mmg to VTK
242 for ( MMG5_int k=1; k<=mesh->nquad; ++k ) {
243 MMG5_pQuad pq = &mesh->quadra[k];
244 if ( !MG_EOK(pq) ) continue;
245
246 if ( pq->ref ) hasCellRef = 1;
247
248 for ( int i=0; i<4; ++i ) {
249 quadra->GetPointIds()->SetId(i, mesh->point[pq->v[i]].tmp );
250 }
251 cellArray->InsertNextCell(quadra);
252 types[ic++] = VTK_QUAD;
253 }
254
255 // Transfer tets from Mmg to VTK
256 for ( MMG5_int k=1; k<=mesh->ne; ++k ) {
257 MMG5_pTetra pt = &mesh->tetra[k];
258 if ( !MG_EOK(pt) ) continue;
259
260 if ( pt->ref ) hasCellRef = 1;
261
262 for ( int i=0; i<4; ++i ) {
263 tetra->GetPointIds()->SetId(i, mesh->point[pt->v[i]].tmp );
264 }
265 cellArray->InsertNextCell(tetra);
266 types[ic++] = VTK_TETRA;
267 }
268
269 // Transfer prisms from Mmg to VTK
270 for ( MMG5_int k=1; k<=mesh->nprism; ++k ) {
271 MMG5_pPrism ppr = &mesh->prism[k];
272 if ( !MG_EOK(ppr) ) continue;
273
274 if ( ppr->ref ) hasCellRef = 1;
275
276 for ( int i=0; i<6; ++i ) {
277 wedge->GetPointIds()->SetId(i, mesh->point[ppr->v[i]].tmp );
278 }
279 cellArray->InsertNextCell(wedge);
280 types[ic++] = VTK_WEDGE;
281 }
282
283 assert ( ic == nc );
284
285 // Transfer cell array into data set
286 MMG5_internal_VTKSetCells(types,dataset,cellArray);
287
288 // Transfer references if needed (i.e. the mesh contains non 0 refs)
289 if ( hasPointRef ) {
290 auto *ar = vtkFloatArray::New();
291
292 ar->SetNumberOfComponents(1);
293 ar->SetNumberOfTuples(mesh->np);
294 ar->SetName("medit:ref");
295
296 for ( MMG5_int k = 0; k < mesh->np; k++ ) {
297 MMG5_pPoint ppt = &mesh->point[k+1];
298 if ( !MG_VOK(ppt) ) continue;
299
300 ar->SetTuple1(ppt->tmp,ppt->ref);
301 }
302
303 dataset->GetPointData()->AddArray(ar);
304 }
305 if ( hasCellRef ) {
306 auto *ar = vtkFloatArray::New();
307
308 ar->SetNumberOfComponents(1);
309 ar->SetNumberOfTuples(nc);
310 ar->SetName("medit:ref");
311
312 MMG5_int idx = 0;
313 for ( MMG5_int k = 1; k <= mesh->na; k++ ) {
314 if ( !mesh->edge[k].a ) continue;
315 ar->SetTuple1(idx++,mesh->edge[k].ref);
316 }
317 for ( MMG5_int k = 1; k <= mesh->nt; k++ ) {
318 if ( !MG_EOK(&mesh->tria[k]) ) continue;
319 ar->SetTuple1(idx++,mesh->tria[k].ref);
320 }
321 for ( MMG5_int k = 1; k <= mesh->nquad; k++ ) {
322 if ( !MG_EOK(&mesh->quadra[k]) ) continue;
323 ar->SetTuple1(idx++,mesh->quadra[k].ref);
324 }
325 for ( MMG5_int k = 1; k <= mesh->ne; k++ ) {
326 if ( !MG_EOK(&mesh->tetra[k]) ) continue;
327 ar->SetTuple1(idx++,mesh->tetra[k].ref);
328 }
329 for ( MMG5_int k = 1; k <= mesh->nprism; k++ ) {
330 if ( !MG_EOK(&mesh->prism[k]) ) continue;
331 ar->SetTuple1(idx++,mesh->prism[k].ref);
332 }
333
334 dataset->GetCellData()->AddArray(ar);
335 }
336
337 // Transfer point solutions into data set
338 MMG5_pSol psl = NULL;
339 int nsols;
340
341 if ( metricData==1 ) {
342 if ( sol && *sol && sol[0]->np ) {
343 nsols = 1;
344 }
345 else {
346 /* In analysis mode (-noinsert -noswap -nomove), metric is not allocated */
347 nsols = 0;
348 }
349 }
350 else {
351 nsols = mesh->nsols;
352 }
353
354 static int mmgWarn = 0;
355 for ( int isol=0; isol<nsols; ++isol) {
356 psl = *sol + isol;
357
358 if ( !psl->m ) {
359 if ( !mmgWarn ) {
360 mmgWarn = 1;
361 fprintf(stderr, " ## Warning: %s: missing data for at least 1 solution."
362 " Skipped.\n",__func__);
363 }
364 continue;
365 }
366
367 int ncp;
368 if ( psl->size == 1 ) {
369 ncp = 1;
370 }
371 else if ( psl->size == psl->dim ) {
372 ncp = psl->dim;
373 }
374 else {
375 ncp = psl->dim*psl->dim;
376 }
377
378 auto *ar = vtkDoubleArray::New();
379
380 ar->SetNumberOfComponents(ncp);
381 ar->SetNumberOfTuples(mesh->np);
382
383 if ( psl->namein ) {
384 char *tmp = MMG5_Get_basename(psl->namein);
385 char *data;
386
387 MMG5_SAFE_CALLOC(data,strlen(tmp)+8,char,
388 MMG5_SAFE_FREE ( types ); return 0);
389
390 strcpy(data,tmp);
391 free(tmp); tmp = 0;
392
393 if ( metricData ) {
394 strcat ( data , ":metric");
395 }
396 ar->SetName(data);
397
398 MMG5_DEL_MEM(mesh,data);
399 }
400 else {
401 ar->SetName("no_name");
402 }
403
404 double* dfmt = NULL;
405 MMG5_SAFE_MALLOC ( dfmt, ncp, double, MMG5_SAFE_FREE ( types );return 0 );
406
407 if ( psl->size!= (psl->dim*(psl->dim+1))/2 ) {
408 /* scalar or vector field: data order and size isn't modified */
409 for ( MMG5_int k=1; k<=mesh->np; k++) {
410 MMG5_pPoint ppt = &mesh->point[k];
411 if ( !MG_VOK(ppt) ) continue;
412
413 MMG5_int iadr = k*psl->size;
414 for ( int i=0; i<psl->size; ++i ) {
415 dfmt[i] = psl->m[iadr+i];
416 }
417
418 if ( psl->dim==2 && ncp==3 ) {
419 dfmt[2] = 0; // z-component for a vector field
420 }
421 ar->SetTuple(ppt->tmp,dfmt);
422 }
423 }
424 else {
425 /* Tensor field: we need to reconstruct the entire symetric matrix */
426 for ( MMG5_int k=1; k<=mesh->np; k++) {
427 MMG5_pPoint ppt = &mesh->point[k];
428 if ( !MG_VOK(ppt) ) continue;
429
430 MMG5_int iadr = k*psl->size;
431 double *d = &psl->m[iadr];
432 if ( psl->dim==2 ) {
433 dfmt[0] = d[0];
434 dfmt[1] = d[1];
435 dfmt[2] = d[1];
436 dfmt[3] = d[2];
437 }
438 else {
439 double dbuf[6] = {0,0,0,0,0,0};
440 if ( metricData ) {
441 assert(!mesh->nsols);
442 MMG5_build3DMetric(mesh,psl,k,dbuf);
443 }
444 else {
445 for ( int i=0; i<psl->size; i++) dbuf[i] = psl->m[psl->size*k+i];
446 }
447 dfmt[0] = dbuf[0];
448 dfmt[1] = dbuf[1];
449 dfmt[2] = dbuf[2];
450 dfmt[3] = dbuf[1];
451 dfmt[4] = dbuf[3];
452 dfmt[5] = dbuf[4];
453 dfmt[6] = dbuf[2];
454 dfmt[7] = dbuf[4];
455 dfmt[8] = dbuf[5];
456 }
457 ar->SetTuple(ppt->tmp,dfmt);
458 }
459 }
460 dataset->GetPointData()->AddArray(ar);
461
462 MMG5_SAFE_FREE(dfmt);
463 }
464
465 if ( npart ) {
466 // distributed IO
467 vtkSmartPointer<PWriter> writer = vtkSmartPointer<PWriter>::New();
468
469#if VTK_MAJOR_VERSION <= 5
470 writer->SetInput(dataset);
471#else
472 writer->SetInputData(dataset);
473#endif
474
475 writer->SetFileName(mfilename);
476
477 writer->SetNumberOfPieces(npart);
478 writer->SetStartPiece(myid);
479 writer->SetEndPiece(myid);
480 writer->Write();
481 }
482 else {
483 // centralized IO
484 vtkSmartPointer<TWriter> writer = vtkSmartPointer<TWriter>::New();
485
486 writer->SetFileName(mfilename);
487
488#if VTK_MAJOR_VERSION <= 5
489 writer->SetInput(dataset);
490#else
491 writer->SetInputData(dataset);
492#endif
493
494 //MMG5_internal_VTKbinary(writer,binary);
495 writer->Write();
496 }
497
498 MMG5_SAFE_FREE(types);
499
500 return 1;
501}
502
515template <class T, class TWriter,class PWriter>
517 int metricData,int binary) {
518
519 return MMG5_saveVtkMesh_i<T,TWriter,PWriter>( mesh,sol,filename,
520 metricData,binary,0,0,0);
521
522}
523
524#endif
525
526#endif
char * MMG5_Get_basename(char *path)
tmp[*strlen0]
MMG5_pMesh MMG5_pSol * sol
MMG5_pMesh * mesh
MMG5_pMesh char * filename
void MMG5_build3DMetric(MMG5_pMesh mesh, MMG5_pSol sol, MMG5_int ip, double dbuf[6])
Definition: inout.c:1528
#define MMG5_SAFE_CALLOC(ptr, size, type, law)
#define MG_EOK(pt)
#define MG_VOK(ppt)
#define MMG5_SAFE_MALLOC(ptr, size, type, law)
#define MMG5_SAFE_FREE(ptr)
#define MMG5_DEL_MEM(mesh, ptr)
Structure to store edges of a MMG mesh.
Definition: libmmgtypes.h:305
MMG5_int b
Definition: libmmgtypes.h:306
MMG5_int ref
Definition: libmmgtypes.h:307
MMG5_int a
Definition: libmmgtypes.h:306
MMG mesh structure.
Definition: libmmgtypes.h:605
MMG5_pQuad quadra
Definition: libmmgtypes.h:648
MMG5_int ne
Definition: libmmgtypes.h:612
MMG5_pPoint point
Definition: libmmgtypes.h:641
MMG5_pPrism prism
Definition: libmmgtypes.h:645
MMG5_int nquad
Definition: libmmgtypes.h:613
MMG5_pTetra tetra
Definition: libmmgtypes.h:643
MMG5_int nt
Definition: libmmgtypes.h:612
MMG5_pTria tria
Definition: libmmgtypes.h:647
MMG5_int np
Definition: libmmgtypes.h:612
MMG5_pEdge edge
Definition: libmmgtypes.h:649
MMG5_int nprism
Definition: libmmgtypes.h:613
MMG5_int na
Definition: libmmgtypes.h:612
Structure to store points of a MMG mesh.
Definition: libmmgtypes.h:270
MMG5_int tmp
Definition: libmmgtypes.h:280
double c[3]
Definition: libmmgtypes.h:271
MMG5_int ref
Definition: libmmgtypes.h:278
MMG5_int v[6]
Definition: libmmgtypes.h:463
MMG5_int ref
Definition: libmmgtypes.h:464
MMG5_int ref
Definition: libmmgtypes.h:368
MMG5_int v[4]
Definition: libmmgtypes.h:367
char * namein
Definition: libmmgtypes.h:673
double * m
Definition: libmmgtypes.h:671
MMG5_int v[4]
Definition: libmmgtypes.h:403
MMG5_int ref
Definition: libmmgtypes.h:404
MMG5_int ref
Definition: libmmgtypes.h:335
MMG5_int v[3]
Definition: libmmgtypes.h:334
int MMG5_loadVtkMesh_part1(MMG5_pMesh, const char *, vtkDataSet **, int8_t *, int8_t *, int *, int8_t *)
Definition: vtkparser.cpp:269
int MMG5_loadVtuMesh_part1(MMG5_pMesh, const char *, vtkDataSet **, int8_t *, int8_t *, int *, int8_t *)
Definition: vtkparser.cpp:308
static void MMG5_internal_VTKSetLine(vtkSmartPointer< vtkPolyData > d, vtkSmartPointer< vtkCellArray > *ca)
Definition: vtkparser.hpp:76
int MMG5_saveVtkMesh(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename, int metricData, int binary)
Definition: vtkparser.hpp:516
static void MMG5_internal_VTKbinary(vtkXMLUnstructuredGridWriter *w, int binary)
Definition: vtkparser.hpp:120
int MMG5_saveVtkMesh_i(MMG5_pMesh mesh, MMG5_pSol *sol, const char *mfilename, int metricData, int binary, int npart, int myid, int master)
Definition: vtkparser.hpp:170
int MMG5_loadVtkMesh_part2(MMG5_pMesh, MMG5_pSol *, vtkDataSet **, int8_t, int8_t, int)
Definition: vtkparser.cpp:346
static void MMG5_internal_VTKSetCells(int *t, vtkSmartPointer< vtkPolyData > d, vtkSmartPointer< vtkCellArray > ca)
Definition: vtkparser.hpp:101
int MMG5_loadVtpMesh_part1(MMG5_pMesh, const char *, vtkDataSet **, int8_t *, int8_t *, int *, int8_t *)
Definition: vtkparser.cpp:230