Chaste Commit::675f9facbe008c5eacb9006feaeb6423206579ea
VtkMeshWriter.cpp
1/*
2
3Copyright (c) 2005-2025, University of Oxford.
4All rights reserved.
5
6University of Oxford means the Chancellor, Masters and Scholars of the
7University of Oxford, having an administrative office at Wellington
8Square, Oxford OX1 2JD, UK.
9
10This file is part of Chaste.
11
12Redistribution and use in source and binary forms, with or without
13modification, are permitted provided that the following conditions are met:
14 * Redistributions of source code must retain the above copyright notice,
15 this list of conditions and the following disclaimer.
16 * Redistributions in binary form must reproduce the above copyright notice,
17 this list of conditions and the following disclaimer in the documentation
18 and/or other materials provided with the distribution.
19 * Neither the name of the University of Oxford nor the names of its
20 contributors may be used to endorse or promote products derived from this
21 software without specific prior written permission.
22
23THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
24AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
25IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
26ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
27LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
28CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
29GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
30HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
31LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
32OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
33
34*/
35
36#include <boost/scoped_array.hpp>
37#include "VtkMeshWriter.hpp"
38#include "DistributedTetrahedralMesh.hpp"
39#include "MixedDimensionMesh.hpp"
40#include "NodesOnlyMesh.hpp"
41
42#ifdef CHASTE_VTK
43#include "vtkQuadraticTetra.h"
44#include "vtkQuadraticTriangle.h"
45
46
48// Implementation
50template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
52 const std::string& rBaseName,
53 const bool& rCleanDirectory)
54 : AbstractTetrahedralMeshWriter<ELEMENT_DIM, SPACE_DIM>(rDirectory, rBaseName, rCleanDirectory),
55 mWriteParallelFiles(false),
56 mWriteMeshCells(true)
57{
58 this->mIndexFromZero = true;
59
60 // Dubious, since we shouldn't yet know what any details of the mesh are.
61 mpVtkUnstructedMesh = vtkUnstructuredGrid::New();
62
63
64}
65
66template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
68{
69 mpVtkUnstructedMesh->Delete(); // Reference counted
70}
71
72template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
74{
75 //Construct nodes aka as Points
76 vtkPoints* p_pts = vtkPoints::New(VTK_DOUBLE);
77 p_pts->GetData()->SetName("Vertex positions");
78 for (unsigned item_num=0; item_num<this->GetNumNodes(); item_num++)
79 {
80 std::vector<double> current_item = this->GetNextNode(); //this->mNodeData[item_num];
81 // Add zeroes if the dimension is below 3
82 for (unsigned dim=SPACE_DIM; dim<3; dim++)
83 {
84 current_item.push_back(0.0);//For y and z-coordinates if necessary
85 }
86 assert(current_item.size() == 3);
87 p_pts->InsertPoint(item_num, current_item[0], current_item[1], current_item[2]);
88 }
89 mpVtkUnstructedMesh->SetPoints(p_pts);
90 p_pts->Delete(); //Reference counted
91
92 if (mWriteMeshCells == true)
93 {
94 //Construct elements aka Cells
95 for (unsigned item_num=0; item_num<this->GetNumElements(); item_num++)
96 {
97 std::vector<unsigned> current_element = this->GetNextElement().NodeIndices; // this->mElementData[item_num];
98
99 assert((current_element.size() == ELEMENT_DIM + 1) || (current_element.size() == (ELEMENT_DIM+1)*(ELEMENT_DIM+2)/2));
100
101 vtkCell* p_cell=nullptr;
102 if (ELEMENT_DIM == 3 && current_element.size() == 4)
103 {
104 p_cell = vtkTetra::New();
105 }
106 else if (ELEMENT_DIM == 3 && current_element.size() == 10)
107 {
108 p_cell = vtkQuadraticTetra::New();
109 }
110 else if (ELEMENT_DIM == 2 && current_element.size() == 3)
111 {
112 p_cell = vtkTriangle::New();
113 }
114 else if (ELEMENT_DIM == 2 && current_element.size() == 6)
116 p_cell = vtkQuadraticTriangle::New();
117 }
118 else if (ELEMENT_DIM == 1)
119 {
120 p_cell = vtkLine::New();
121 }
122
123 //Set the linear nodes
124 vtkIdList* p_cell_id_list = p_cell->GetPointIds();
125 for (unsigned j = 0; j < current_element.size(); ++j)
127 p_cell_id_list->SetId(j, current_element[j]);
128 }
129
130 //VTK defines the node ordering in quadratic triangles differently to Chaste, so they must be treated as a special case
131 if (SPACE_DIM == 2 && current_element.size() == 6)
132 {
133 p_cell_id_list->SetId(3, current_element[5]);
134 p_cell_id_list->SetId(4, current_element[3]);
135 p_cell_id_list->SetId(5, current_element[4]);
136 }
137
138 mpVtkUnstructedMesh->InsertNextCell(p_cell->GetCellType(), p_cell_id_list);
139 p_cell->Delete(); //Reference counted
140 }
142
143 if (SPACE_DIM > 1)
144 {
146 for (unsigned item_num=0; item_num<this->GetNumBoundaryFaces(); item_num++)
147 {
148 this->GetNextBoundaryElement();
149 }
151
152 //If necessary, construct cables
153 if (this->GetNumCableElements() > 0)
154 {
155 AugmentCellData();
156 //Make a blank cell radius data for the regular elements
157 std::vector<double> radii(this->GetNumElements(), 0.0);
158 for (unsigned item_num=0; item_num<this->GetNumCableElements(); item_num++)
159 {
160 ElementData cable_element_data = this->GetNextCableElement();
161 std::vector<unsigned> current_element = cable_element_data.NodeIndices;
162 radii.push_back(cable_element_data.AttributeValue);
163 assert(current_element.size() == 2);
164 vtkCell* p_cell=vtkLine::New();
165 vtkIdList* p_cell_id_list = p_cell->GetPointIds();
166 for (unsigned j = 0; j < 2; ++j)
167 {
168 p_cell_id_list->SetId(j, current_element[j]);
170 mpVtkUnstructedMesh->InsertNextCell(p_cell->GetCellType(), p_cell_id_list);
171 p_cell->Delete(); //Reference counted
172 }
173 AddCellData("Cable radius", radii);
174
175 }
176
177}
178
179template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
181{
182 mWriteMeshCells = writeMeshCells;
183}
184
185template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
187{
188 std::string comment = "<!-- " + ChasteBuildInfo::GetProvenanceString() + "-->";
189
190 out_stream p_vtu_file = this->mpOutputFileHandler->OpenOutputFile(fileName, std::ios::out | std::ios::app);
191
192 *p_vtu_file << "\n" << comment << "\n";
193 p_vtu_file->close();
194}
195
196template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
198{
199 // Using separate scope here to make sure file is properly closed before re-opening it to add provenance info.
200 {
201 MakeVtkMesh();
202 assert(mpVtkUnstructedMesh->CheckAttributes() == 0);
203 vtkXMLUnstructuredGridWriter* p_writer = vtkXMLUnstructuredGridWriter::New();
204#if VTK_MAJOR_VERSION >= 6
205 p_writer->SetInputData(mpVtkUnstructedMesh);
206#else
207 p_writer->SetInput(mpVtkUnstructedMesh);
208#endif
209 std::string vtk_file_name = this->mpOutputFileHandler->GetOutputDirectoryFullPath() + this->mBaseName+".vtu";
210 p_writer->SetFileName(vtk_file_name.c_str());
211 //p_writer->PrintSelf(std::cout, vtkIndent());
212 p_writer->Write();
213 p_writer->Delete(); //Reference counted
214 }
215
216 AddProvenance(this->mBaseName + ".vtu");
217}
218
219template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
220void VtkMeshWriter<ELEMENT_DIM,SPACE_DIM>::AddCellData(std::string dataName, std::vector<double> dataPayload)
221{
222 vtkDoubleArray* p_scalars = vtkDoubleArray::New();
223 p_scalars->SetName(dataName.c_str());
224 for (unsigned i=0; i<dataPayload.size(); i++)
225 {
226 p_scalars->InsertNextValue(dataPayload[i]);
227 }
228
229 vtkCellData* p_cell_data = mpVtkUnstructedMesh->GetCellData();
230 p_cell_data->AddArray(p_scalars);
231 p_scalars->Delete(); //Reference counted
233
234template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
236{
237 unsigned num_cell_arrays = mpVtkUnstructedMesh->GetCellData()->GetNumberOfArrays();
238 for (unsigned i = 0; i < num_cell_arrays; i++)
239 {
240 vtkDataArray* array = mpVtkUnstructedMesh->GetCellData()->GetArray(i);
241
242 //Check data was the correct size before the cables were added
243 unsigned num_cable_pads = this->GetNumCableElements();
244 if (mWriteParallelFiles)
246 assert((unsigned)array->GetNumberOfTuples() == this->mpDistributedMesh->GetNumLocalElements());
247 num_cable_pads = this->mpMixedMesh->GetNumLocalCableElements();
248 }
249 else
250 {
251 assert((unsigned)array->GetNumberOfTuples() == this->GetNumElements());
252 }
253
254 //Check that tuples of size 3 will be big enough for padding the rest of the data
255 assert(array->GetNumberOfComponents() <= 3);
256 double null_data[3] = {0.0, 0.0, 0.0};
257
258 //Pad data
259 for (unsigned new_index = 0; new_index < num_cable_pads; new_index++)
260 {
261 array->InsertNextTuple(null_data);
262 }
263 }
264}
265
266
267template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
268void VtkMeshWriter<ELEMENT_DIM,SPACE_DIM>::AddCellData(std::string dataName, std::vector<c_vector<double, SPACE_DIM> > dataPayload)
269{
270 vtkDoubleArray* p_vectors = vtkDoubleArray::New();
271 p_vectors->SetName(dataName.c_str());
272 p_vectors->SetNumberOfComponents(3);
273 for (unsigned i=0; i<dataPayload.size(); i++)
274 {
275 for (unsigned j=0; j<SPACE_DIM; j++)
276 {
277 p_vectors->InsertNextValue(dataPayload[i][j]);
278 }
279 //When SPACE_DIM<3, then pad
280 for (unsigned j=SPACE_DIM; j<3; j++)
281 {
282 p_vectors->InsertNextValue(0.0);
283 }
284 }
285
286 vtkCellData* p_cell_data = mpVtkUnstructedMesh->GetCellData();
287 p_cell_data->AddArray(p_vectors);
288 p_vectors->Delete(); //Reference counted
289}
290
291template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
292void VtkMeshWriter<ELEMENT_DIM,SPACE_DIM>::AddTensorCellData(std::string dataName, std::vector<c_vector<double,SPACE_DIM*(SPACE_DIM+1)/2> > dataPayload)
293{
294 assert(SPACE_DIM != 1); // LCOV_EXCL_LINE
295
296 vtkDoubleArray* p_vectors = vtkDoubleArray::New();
297 p_vectors->SetName(dataName.c_str());
298 p_vectors->SetNumberOfComponents(SPACE_DIM*SPACE_DIM);
299 for (unsigned i=0; i<dataPayload.size(); i++)
300 {
301 if (SPACE_DIM == 2)
302 {
303 p_vectors->InsertNextValue(dataPayload[i](0)); //a11
304 p_vectors->InsertNextValue(dataPayload[i](1)); //a12
305 p_vectors->InsertNextValue(dataPayload[i](1)); //a21
306 p_vectors->InsertNextValue(dataPayload[i](2)); //a22
307 }
308 else if (SPACE_DIM == 3)
309 {
310 p_vectors->InsertNextValue(dataPayload[i](0)); //a11
311 p_vectors->InsertNextValue(dataPayload[i](1)); //a12
312 p_vectors->InsertNextValue(dataPayload[i](2)); //a13
313 p_vectors->InsertNextValue(dataPayload[i](1)); //a21
314 p_vectors->InsertNextValue(dataPayload[i](3)); //a22
315 p_vectors->InsertNextValue(dataPayload[i](4)); //a23
316 p_vectors->InsertNextValue(dataPayload[i](2)); //a31
317 p_vectors->InsertNextValue(dataPayload[i](4)); //a32
318 p_vectors->InsertNextValue(dataPayload[i](5)); //a33
319 }
320 }
321
322 vtkCellData* p_cell_data = mpVtkUnstructedMesh->GetCellData();
323 p_cell_data->AddArray(p_vectors);
324 p_vectors->Delete(); //Reference counted
325}
326
327template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
328void VtkMeshWriter<ELEMENT_DIM,SPACE_DIM>::AddTensorCellData(std::string dataName, std::vector<c_matrix<double,SPACE_DIM,SPACE_DIM> > dataPayload)
329{
330 assert(SPACE_DIM != 1); // LCOV_EXCL_LINE
331
332 vtkDoubleArray* p_vectors = vtkDoubleArray::New();
333 p_vectors->SetName(dataName.c_str());
334 p_vectors->SetNumberOfComponents(SPACE_DIM*SPACE_DIM);
335 for (unsigned i=0; i<dataPayload.size(); i++)
336 {
337 if (SPACE_DIM == 2)
338 {
339 p_vectors->InsertNextValue(dataPayload[i](0,0)); //a11
340 p_vectors->InsertNextValue(dataPayload[i](0,1)); //a12
341 p_vectors->InsertNextValue(dataPayload[i](1,0)); //a21
342 p_vectors->InsertNextValue(dataPayload[i](1,1)); //a22
343 }
344 else if (SPACE_DIM == 3)
345 {
346 p_vectors->InsertNextValue(dataPayload[i](0,0)); //a11
347 p_vectors->InsertNextValue(dataPayload[i](0,1)); //a12
348 p_vectors->InsertNextValue(dataPayload[i](0,2)); //a13
349 p_vectors->InsertNextValue(dataPayload[i](1,0)); //a21
350 p_vectors->InsertNextValue(dataPayload[i](1,1)); //a22
351 p_vectors->InsertNextValue(dataPayload[i](1,2)); //a23
352 p_vectors->InsertNextValue(dataPayload[i](2,0)); //a31
353 p_vectors->InsertNextValue(dataPayload[i](2,1)); //a32
354 p_vectors->InsertNextValue(dataPayload[i](2,2)); //a33
355 }
356 }
357
358 vtkCellData* p_cell_data = mpVtkUnstructedMesh->GetCellData();
359 p_cell_data->AddArray(p_vectors);
360 p_vectors->Delete(); //Reference counted
361}
362
363
364template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
365void VtkMeshWriter<ELEMENT_DIM,SPACE_DIM>::AddPointData(std::string dataName, std::vector<double> dataPayload)
366{
367 vtkDoubleArray* p_scalars = vtkDoubleArray::New();
368 p_scalars->SetName(dataName.c_str());
369
370 if (mWriteParallelFiles && this->mpDistributedMesh != nullptr)
371 {
372 // In parallel, the vector we pass will only contain the values from the privately owned nodes.
373 // To get the values from the halo nodes (which will be inserted at the end of the vector we need to
374 // communicate with the equivalent vectors on other processes.
375
376 // resize the payload data to include halos
377 assert( dataPayload.size() == this->mpDistributedMesh->GetNumLocalNodes() );
378 dataPayload.resize( this->mpDistributedMesh->GetNumLocalNodes() + this->mpDistributedMesh->GetNumHaloNodes() );
379
380
381 // then do the communication
382 for ( unsigned rank_offset = 1; rank_offset < PetscTools::GetNumProcs(); rank_offset++ )
383 {
384 unsigned send_to = (PetscTools::GetMyRank() + rank_offset) % (PetscTools::GetNumProcs());
385 unsigned receive_from = (PetscTools::GetMyRank() + PetscTools::GetNumProcs()- rank_offset ) % (PetscTools::GetNumProcs());
386
387 unsigned number_of_nodes_to_send = mNodesToSendPerProcess[send_to].size();
388 unsigned number_of_nodes_to_receive = mNodesToReceivePerProcess[receive_from].size();
389
390 boost::scoped_array<double> send_data(new double[number_of_nodes_to_send]);
391 boost::scoped_array<double> receive_data(new double[number_of_nodes_to_receive]);
392 // Pack
393 for (unsigned node = 0; node < number_of_nodes_to_send; node++)
394 {
395 unsigned global_node_index = mNodesToSendPerProcess[send_to][node];
396 unsigned local_node_index = global_node_index
397 - this->mpDistributedMesh->GetDistributedVectorFactory()->GetLow();
398 send_data[node] = dataPayload[local_node_index];
399 }
400 {
401 // Send
402 int ret;
403 MPI_Status status;
404 ret = MPI_Sendrecv(send_data.get(), number_of_nodes_to_send,
405 MPI_DOUBLE,
406 send_to, 0,
407 receive_data.get(), number_of_nodes_to_receive,
408 MPI_DOUBLE,
409 receive_from, 0,
410 PETSC_COMM_WORLD, &status);
411 UNUSED_OPT(ret);
412 assert ( ret == MPI_SUCCESS );
413 }
414
415 // Unpack
416 for ( unsigned node = 0; node < number_of_nodes_to_receive; node++ )
417 {
418 unsigned global_node_index = mNodesToReceivePerProcess[receive_from][node];
419 unsigned halo_index = mGlobalToNodeIndexMap[global_node_index];
420 assert( halo_index >= this->mpDistributedMesh->GetNumLocalNodes() );
421 dataPayload[halo_index] = receive_data[node];
422 }
423
424 }
425 }
426
427 for (unsigned i=0; i<dataPayload.size(); i++)
428 {
429 p_scalars->InsertNextValue(dataPayload[i]);
430 }
431
432 vtkPointData* p_point_data = mpVtkUnstructedMesh->GetPointData();
433 p_point_data->AddArray(p_scalars);
434 p_scalars->Delete(); //Reference counted
435}
436
437
438template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
439void VtkMeshWriter<ELEMENT_DIM,SPACE_DIM>::AddPointData(std::string dataName, std::vector<c_vector<double, SPACE_DIM> > dataPayload)
440{
441 vtkDoubleArray* p_vectors = vtkDoubleArray::New();
442 p_vectors->SetName(dataName.c_str());
443
444 if (mWriteParallelFiles)
445 {
446 // In parallel, the vector we pass will only contain the values from the privately owned nodes.
447 // To get the values from the halo nodes (which will be inserted at the end of the vector we need to
448 // communicate with the equivalent vectors on other processes.
449
450 // resize the payload data to include halos
451 assert( dataPayload.size() == this->mpDistributedMesh->GetNumLocalNodes() );
452 dataPayload.resize( this->mpDistributedMesh->GetNumLocalNodes() + this->mpDistributedMesh->GetNumHaloNodes() );
453
454 // then do the communication
455 for ( unsigned rank_offset = 1; rank_offset < PetscTools::GetNumProcs(); rank_offset++ )
456 {
457 unsigned send_to = (PetscTools::GetMyRank() + rank_offset) % (PetscTools::GetNumProcs());
458 unsigned receive_from = (PetscTools::GetMyRank() + PetscTools::GetNumProcs()- rank_offset ) % (PetscTools::GetNumProcs());
459
460 unsigned number_of_nodes_to_send = mNodesToSendPerProcess[send_to].size();
461 unsigned number_of_nodes_to_receive = mNodesToReceivePerProcess[receive_from].size();
462
463 boost::scoped_array<double> send_data(new double[number_of_nodes_to_send * SPACE_DIM]);
464 boost::scoped_array<double> receive_data(new double[number_of_nodes_to_receive * SPACE_DIM]);
465
466 for (unsigned node = 0; node < number_of_nodes_to_send; node++)
467 {
468 unsigned global_node_index = mNodesToSendPerProcess[send_to][node];
469 unsigned local_node_index = global_node_index
470 - this->mpDistributedMesh->GetDistributedVectorFactory()->GetLow();
471 for (unsigned j=0; j<SPACE_DIM; j++)
472 {
473 send_data[ node*SPACE_DIM + j ] = dataPayload[local_node_index][j];
474 }
475 }
476
477 int ret;
478 MPI_Status status;
479 ret = MPI_Sendrecv(send_data.get(), number_of_nodes_to_send * SPACE_DIM,
480 MPI_DOUBLE,
481 send_to, 0,
482 receive_data.get(), number_of_nodes_to_receive * SPACE_DIM,
483 MPI_DOUBLE,
484 receive_from, 0,
485 PETSC_COMM_WORLD, &status);
486 UNUSED_OPT(ret);
487 assert ( ret == MPI_SUCCESS );
488
489 // Unpack
490 for ( unsigned node = 0; node < number_of_nodes_to_receive; node++ )
491 {
492 unsigned global_node_index = mNodesToReceivePerProcess[receive_from][node];
493 unsigned halo_index = mGlobalToNodeIndexMap[global_node_index];
494 assert( halo_index >= this->mpDistributedMesh->GetNumLocalNodes() );
495 for (unsigned j=0; j<SPACE_DIM; j++)
496 {
497 dataPayload[halo_index][j] = receive_data[ node*SPACE_DIM + j ];
498 }
499 }
500 }
501 }
502
503 p_vectors->SetNumberOfComponents(3);
504 for (unsigned i=0; i<dataPayload.size(); i++)
505 {
506 for (unsigned j=0; j<SPACE_DIM; j++)
507 {
508 p_vectors->InsertNextValue(dataPayload[i][j]);
509 }
510 //When SPACE_DIM<3, then pad
511 for (unsigned j=SPACE_DIM; j<3; j++)
512 {
513 p_vectors->InsertNextValue(0.0);
514 }
515 }
516
517 vtkPointData* p_point_data = mpVtkUnstructedMesh->GetPointData();
518 p_point_data->AddArray(p_vectors);
519 p_vectors->Delete(); //Reference counted
520}
521
522template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
523void VtkMeshWriter<ELEMENT_DIM,SPACE_DIM>::AddTensorPointData(std::string dataName, std::vector<c_matrix<double,SPACE_DIM,SPACE_DIM> > dataPayload)
524{
525 assert(SPACE_DIM != 1); // LCOV_EXCL_LINE
526
527 vtkDoubleArray* p_vectors = vtkDoubleArray::New();
528 p_vectors->SetName(dataName.c_str());
529 p_vectors->SetNumberOfComponents(SPACE_DIM*SPACE_DIM);
530 for (unsigned i=0; i<dataPayload.size(); i++)
531 {
532 if (SPACE_DIM == 2)
533 {
534 p_vectors->InsertNextValue(dataPayload[i](0,0)); //a11
535 p_vectors->InsertNextValue(dataPayload[i](0,1)); //a12
536 p_vectors->InsertNextValue(dataPayload[i](1,0)); //a21
537 p_vectors->InsertNextValue(dataPayload[i](1,1)); //a22
538 }
539 else if (SPACE_DIM == 3)
540 {
541 p_vectors->InsertNextValue(dataPayload[i](0,0)); //a11
542 p_vectors->InsertNextValue(dataPayload[i](0,1)); //a12
543 p_vectors->InsertNextValue(dataPayload[i](0,2)); //a13
544 p_vectors->InsertNextValue(dataPayload[i](1,0)); //a21
545 p_vectors->InsertNextValue(dataPayload[i](1,1)); //a22
546 p_vectors->InsertNextValue(dataPayload[i](1,2)); //a23
547 p_vectors->InsertNextValue(dataPayload[i](2,0)); //a31
548 p_vectors->InsertNextValue(dataPayload[i](2,1)); //a32
549 p_vectors->InsertNextValue(dataPayload[i](2,2)); //a33
550 }
551 }
552
553 vtkPointData* p_point_data = mpVtkUnstructedMesh->GetPointData();
554 p_point_data->AddArray(p_vectors);
555 p_vectors->Delete(); //Reference counted
556}
557
558template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
560{
561 //Have we got a distributed mesh?
562 this->mpDistributedMesh = dynamic_cast<DistributedTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* >(&rMesh);
563 mpNodesOnlyMesh = dynamic_cast<NodesOnlyMesh<SPACE_DIM>* >(&rMesh);
564
565 if (this->mpDistributedMesh == nullptr && mpNodesOnlyMesh == nullptr)
566 {
567 EXCEPTION("Cannot write parallel files using a sequential mesh");
568 }
569
571 {
572 return; // mWriteParallelFiles is not set sequentially (so we don't set up data exchange machinery)
573 }
574
575 mWriteParallelFiles = true;
576
577 // Populate the global to node index map (as this will be required to add point data)
578
579 //Node index that we are writing to VTK (index into mNodes and mHaloNodes as if they were concatenated)
580 unsigned index = 0;
581
582 // Owned nodes
584 node_iter != rMesh.GetNodeIteratorEnd();
585 ++node_iter)
586 {
587 mGlobalToNodeIndexMap[node_iter->GetIndex()] = index;
588 index++;
589 }
590
591 // Halo nodes
592 if (this->mpDistributedMesh)
593 {
594 for (typename DistributedTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::HaloNodeIterator halo_iter=this->mpDistributedMesh->GetHaloNodeIteratorBegin();
595 halo_iter != this->mpDistributedMesh->GetHaloNodeIteratorEnd();
596 ++halo_iter)
597 {
598 mGlobalToNodeIndexMap[(*halo_iter)->GetIndex()] = index;
599 index++;
600 }
601
602 //Calculate the halo exchange so that node-wise payloads can be communicated
603 this->mpDistributedMesh->CalculateNodeExchange( mNodesToSendPerProcess, mNodesToReceivePerProcess );
604 }
605}
606
608template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
611 bool keepOriginalElementIndexing)
612{
613 // Have we got a parallel mesh?
614 this->mpDistributedMesh = dynamic_cast<DistributedTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* >(&rMesh);
615 this->mpMixedMesh = dynamic_cast<MixedDimensionMesh<ELEMENT_DIM,SPACE_DIM>* >(&rMesh);
616
617 if (PetscTools::IsSequential() || !mWriteParallelFiles || (this->mpDistributedMesh == nullptr && mpNodesOnlyMesh == nullptr))
618 {
620 }
621 else
622 {
623 //Make the local mesh into a VtkMesh
624 vtkPoints* p_pts = vtkPoints::New(VTK_DOUBLE);
625 p_pts->GetData()->SetName("Vertex positions");
626
627 // Owned nodes
629 node_iter != rMesh.GetNodeIteratorEnd();
630 ++node_iter)
631 {
632 c_vector<double, SPACE_DIM> current_item = node_iter->rGetLocation();
633 if (SPACE_DIM == 3)
634 {
635 p_pts->InsertNextPoint(current_item[0], current_item[1], current_item[2]);
636 }
637 else if (SPACE_DIM == 2)
638 {
639 p_pts->InsertNextPoint(current_item[0], current_item[1], 0.0);
640 }
641 else // (SPACE_DIM == 1)
642 {
643 p_pts->InsertNextPoint(current_item[0], 0.0, 0.0);
644 }
645 }
646
647 // Halo nodes
648 if (this->mpDistributedMesh)
649 {
650 for (typename DistributedTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::HaloNodeIterator halo_iter=this->mpDistributedMesh->GetHaloNodeIteratorBegin();
651 halo_iter != this->mpDistributedMesh->GetHaloNodeIteratorEnd();
652 ++halo_iter)
653 {
654 c_vector<double, SPACE_DIM> current_item = (*halo_iter)->rGetLocation();
655 if (SPACE_DIM == 3)
656 {
657 p_pts->InsertNextPoint(current_item[0], current_item[1], current_item[2]);
658 }
659 else if (SPACE_DIM == 2)
660 {
661 p_pts->InsertNextPoint(current_item[0], current_item[1], 0.0);
662 }
663 else // (SPACE_DIM == 1)
664 {
665 p_pts->InsertNextPoint(current_item[0], 0.0, 0.0);
666 }
667 }
668 }
669
670 mpVtkUnstructedMesh->SetPoints(p_pts);
671 p_pts->Delete(); //Reference counted
672
674 elem_iter != rMesh.GetElementIteratorEnd();
675 ++elem_iter)
676 {
677
678 vtkCell* p_cell=nullptr;
680 if (ELEMENT_DIM == 3)
681 {
682 p_cell = vtkTetra::New();
683 }
684 else if (ELEMENT_DIM == 2)
685 {
686 p_cell = vtkTriangle::New();
687 }
688 else //(ELEMENT_DIM == 1)
689 {
690 p_cell = vtkLine::New();
691 }
692 vtkIdList* p_cell_id_list = p_cell->GetPointIds();
693 for (unsigned j = 0; j < ELEMENT_DIM+1; ++j)
694 {
695 unsigned global_node_index = elem_iter->GetNodeGlobalIndex(j);
696 p_cell_id_list->SetId(j, mGlobalToNodeIndexMap[global_node_index]);
697 }
698 mpVtkUnstructedMesh->InsertNextCell(p_cell->GetCellType(), p_cell_id_list);
699 p_cell->Delete(); //Reference counted
700 }
701 //If necessary, construct cables
702 if (this->mpMixedMesh )
703 {
704 AugmentCellData();
705 //Make a blank cell radius data for the regular elements
706 std::vector<double> radii(this->mpMixedMesh->GetNumLocalElements(), 0.0);
707 for (typename MixedDimensionMesh<ELEMENT_DIM,SPACE_DIM>::CableElementIterator elem_iter = this->mpMixedMesh->GetCableElementIteratorBegin();
708 elem_iter != this->mpMixedMesh->GetCableElementIteratorEnd();
709 ++elem_iter)
710 {
711 radii.push_back((*elem_iter)->GetAttribute());
712 vtkCell* p_cell=vtkLine::New();
713 vtkIdList* p_cell_id_list = p_cell->GetPointIds();
714 for (unsigned j = 0; j < 2; ++j)
715 {
716 unsigned global_node_index = (*elem_iter)->GetNodeGlobalIndex(j);
717 p_cell_id_list->SetId(j, mGlobalToNodeIndexMap[global_node_index]);
718 }
719 mpVtkUnstructedMesh->InsertNextCell(p_cell->GetCellType(), p_cell_id_list);
720 p_cell->Delete(); //Reference counted
721 }
722 AddCellData("Cable radius", radii);
723 }
724
725
726 //This block is to guard the mesh writers (vtkXMLPUnstructuredGridWriter) so that they
727 //go out of scope, flush buffers and close files
728 {
729 assert(mpVtkUnstructedMesh->CheckAttributes() == 0);
730 vtkXMLPUnstructuredGridWriter* p_writer = vtkXMLPUnstructuredGridWriter::New();
731
732 p_writer->SetDataModeToBinary();
733
734 p_writer->SetNumberOfPieces(PetscTools::GetNumProcs());
735 //p_writer->SetGhostLevel(-1);
736 p_writer->SetStartPiece(PetscTools::GetMyRank());
737 p_writer->SetEndPiece(PetscTools::GetMyRank());
738
739
740#if VTK_MAJOR_VERSION >= 6
741 p_writer->SetInputData(mpVtkUnstructedMesh);
742#else
743 p_writer->SetInput(mpVtkUnstructedMesh);
744#endif
745 std::string pvtk_file_name = this->mpOutputFileHandler->GetOutputDirectoryFullPath() + this->mBaseName+ ".pvtu";
746 p_writer->SetFileName(pvtk_file_name.c_str());
747 //p_writer->PrintSelf(std::cout, vtkIndent());
748 p_writer->Write();
749 p_writer->Delete(); //Reference counted
750 }
751
752 // Add provenance to the individual files
753 std::stringstream filepath;
754 filepath << this->mBaseName << "_" << PetscTools::GetMyRank() << ".vtu";
755 AddProvenance(filepath.str());
758 {
759 AddProvenance(this->mBaseName+ ".pvtu");
760 }
761 }
762}
763
764// Explicit instantiation
765template class VtkMeshWriter<1,1>;
766template class VtkMeshWriter<1,2>;
767template class VtkMeshWriter<1,3>;
768template class VtkMeshWriter<2,2>; // Actually used
769template class VtkMeshWriter<2,3>;
770template class VtkMeshWriter<3,3>; // Actually used
771
772#endif //CHASTE_VTK
#define EXCEPTION(message)
#define UNUSED_OPT(var)
NodeIterator GetNodeIteratorEnd()
NodeIterator GetNodeIteratorBegin(bool skipDeletedNodes=true)
virtual void WriteFilesUsingMesh(AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > &rMesh, bool keepOriginalElementIndexing=true)
ElementIterator GetElementIteratorBegin(bool skipDeletedElements=true)
static std::string GetProvenanceString()
std::vector< Node< SPACE_DIM > * >::const_iterator HaloNodeIterator
std::vector< Element< 1, SPACE_DIM > * >::const_iterator CableElementIterator
static bool AmMaster()
static bool IsSequential()
static unsigned GetMyRank()
static unsigned GetNumProcs()
void AddPointData(std::string name, std::vector< double > data)
void SetParallelFiles(AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > &rMesh)
void AddTensorPointData(std::string name, std::vector< c_matrix< double, SPACE_DIM, SPACE_DIM > > data)
void SetWriteMeshCells(bool writeMeshCells)
void AddProvenance(std::string fileName)
void AddTensorCellData(std::string name, std::vector< c_vector< double, SPACE_DIM *(SPACE_DIM+1)/2 > > data)
void AddCellData(std::string name, std::vector< double > data)
VtkMeshWriter(const std::string &rDirectory, const std::string &rBaseName, const bool &rCleanDirectory=true)
void WriteFilesUsingMesh(AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > &rMesh, bool keepOriginalElementIndexing=true)
vtkUnstructuredGrid * mpVtkUnstructedMesh
virtual ~VtkMeshWriter()
std::vector< unsigned > NodeIndices