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