36 #include <boost/scoped_array.hpp>
37 #include "VtkMeshWriter.hpp"
38 #include "DistributedTetrahedralMesh.hpp"
39 #include "MixedDimensionMesh.hpp"
40 #include "NodesOnlyMesh.hpp"
43 #include "vtkQuadraticTetra.h"
44 #include "vtkQuadraticTriangle.h"
50 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
52 const std::string& rBaseName,
53 const bool& rCleanDirectory)
55 mWriteParallelFiles(false)
63 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
66 mpVtkUnstructedMesh->Delete();
69 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
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++)
77 std::vector<double> current_item = this->GetNextNode();
79 for (
unsigned dim=SPACE_DIM; dim<3; dim++)
81 current_item.push_back(0.0);
83 assert(current_item.size() == 3);
84 p_pts->InsertPoint(item_num, current_item[0], current_item[1], current_item[2]);
86 mpVtkUnstructedMesh->SetPoints(p_pts);
90 for (
unsigned item_num=0; item_num<this->GetNumElements(); item_num++)
92 std::vector<unsigned> current_element = this->GetNextElement().NodeIndices;
94 assert((current_element.size() == ELEMENT_DIM + 1) || (current_element.size() == (ELEMENT_DIM+1)*(ELEMENT_DIM+2)/2));
96 vtkCell* p_cell=
nullptr;
97 if (ELEMENT_DIM == 3 && current_element.size() == 4)
99 p_cell = vtkTetra::New();
101 else if (ELEMENT_DIM == 3 && current_element.size() == 10)
103 p_cell = vtkQuadraticTetra::New();
105 else if (ELEMENT_DIM == 2 && current_element.size() == 3)
107 p_cell = vtkTriangle::New();
109 else if (ELEMENT_DIM == 2 && current_element.size() == 6)
111 p_cell = vtkQuadraticTriangle::New();
113 else if (ELEMENT_DIM == 1)
115 p_cell = vtkLine::New();
119 vtkIdList* p_cell_id_list = p_cell->GetPointIds();
120 for (
unsigned j = 0; j < current_element.size(); ++j)
122 p_cell_id_list->SetId(j, current_element[j]);
126 if (SPACE_DIM == 2 && current_element.size() == 6)
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]);
133 mpVtkUnstructedMesh->InsertNextCell(p_cell->GetCellType(), p_cell_id_list);
140 for (
unsigned item_num=0; item_num<this->GetNumBoundaryFaces(); item_num++)
142 this->GetNextBoundaryElement();
147 if (this->GetNumCableElements() > 0)
151 std::vector<double> radii(this->GetNumElements(), 0.0);
152 for (
unsigned item_num=0; item_num<this->GetNumCableElements(); item_num++)
154 ElementData cable_element_data = this->GetNextCableElement();
155 std::vector<unsigned> current_element = cable_element_data.
NodeIndices;
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)
162 p_cell_id_list->SetId(j, current_element[j]);
164 mpVtkUnstructedMesh->InsertNextCell(p_cell->GetCellType(), p_cell_id_list);
167 AddCellData(
"Cable radius", radii);
172 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
177 out_stream p_vtu_file = this->mpOutputFileHandler->OpenOutputFile(fileName, std::ios::out | std::ios::app);
179 *p_vtu_file <<
"\n" << comment <<
"\n";
183 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
189 assert(mpVtkUnstructedMesh->CheckAttributes() == 0);
190 vtkXMLUnstructuredGridWriter* p_writer = vtkXMLUnstructuredGridWriter::New();
191 #if VTK_MAJOR_VERSION >= 6
192 p_writer->SetInputData(mpVtkUnstructedMesh);
194 p_writer->SetInput(mpVtkUnstructedMesh);
196 std::string vtk_file_name = this->mpOutputFileHandler->GetOutputDirectoryFullPath() + this->mBaseName+
".vtu";
197 p_writer->SetFileName(vtk_file_name.c_str());
203 AddProvenance(this->mBaseName +
".vtu");
206 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
209 vtkDoubleArray* p_scalars = vtkDoubleArray::New();
210 p_scalars->SetName(dataName.c_str());
211 for (
unsigned i=0; i<dataPayload.size(); i++)
213 p_scalars->InsertNextValue(dataPayload[i]);
216 vtkCellData* p_cell_data = mpVtkUnstructedMesh->GetCellData();
217 p_cell_data->AddArray(p_scalars);
221 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
224 unsigned num_cell_arrays = mpVtkUnstructedMesh->GetCellData()->GetNumberOfArrays();
225 for (
unsigned i = 0; i < num_cell_arrays; i++)
227 vtkDataArray* array = mpVtkUnstructedMesh->GetCellData()->GetArray(i);
230 unsigned num_cable_pads = this->GetNumCableElements();
231 if (mWriteParallelFiles)
233 assert((
unsigned)array->GetNumberOfTuples() == this->mpDistributedMesh->GetNumLocalElements());
234 num_cable_pads = this->mpMixedMesh->GetNumLocalCableElements();
238 assert((
unsigned)array->GetNumberOfTuples() == this->GetNumElements());
242 assert(array->GetNumberOfComponents() <= 3);
243 double null_data[3] = {0.0, 0.0, 0.0};
246 for (
unsigned new_index = 0; new_index < num_cable_pads; new_index++)
248 array->InsertNextTuple(null_data);
254 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
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++)
262 for (
unsigned j=0; j<SPACE_DIM; j++)
264 p_vectors->InsertNextValue(dataPayload[i][j]);
267 for (
unsigned j=SPACE_DIM; j<3; j++)
269 p_vectors->InsertNextValue(0.0);
273 vtkCellData* p_cell_data = mpVtkUnstructedMesh->GetCellData();
274 p_cell_data->AddArray(p_vectors);
278 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
281 assert(SPACE_DIM != 1);
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++)
290 p_vectors->InsertNextValue(dataPayload[i](0));
291 p_vectors->InsertNextValue(dataPayload[i](1));
292 p_vectors->InsertNextValue(dataPayload[i](1));
293 p_vectors->InsertNextValue(dataPayload[i](2));
295 else if (SPACE_DIM == 3)
297 p_vectors->InsertNextValue(dataPayload[i](0));
298 p_vectors->InsertNextValue(dataPayload[i](1));
299 p_vectors->InsertNextValue(dataPayload[i](2));
300 p_vectors->InsertNextValue(dataPayload[i](1));
301 p_vectors->InsertNextValue(dataPayload[i](3));
302 p_vectors->InsertNextValue(dataPayload[i](4));
303 p_vectors->InsertNextValue(dataPayload[i](2));
304 p_vectors->InsertNextValue(dataPayload[i](4));
305 p_vectors->InsertNextValue(dataPayload[i](5));
309 vtkCellData* p_cell_data = mpVtkUnstructedMesh->GetCellData();
310 p_cell_data->AddArray(p_vectors);
314 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
317 assert(SPACE_DIM != 1);
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++)
326 p_vectors->InsertNextValue(dataPayload[i](0,0));
327 p_vectors->InsertNextValue(dataPayload[i](0,1));
328 p_vectors->InsertNextValue(dataPayload[i](1,0));
329 p_vectors->InsertNextValue(dataPayload[i](1,1));
331 else if (SPACE_DIM == 3)
333 p_vectors->InsertNextValue(dataPayload[i](0,0));
334 p_vectors->InsertNextValue(dataPayload[i](0,1));
335 p_vectors->InsertNextValue(dataPayload[i](0,2));
336 p_vectors->InsertNextValue(dataPayload[i](1,0));
337 p_vectors->InsertNextValue(dataPayload[i](1,1));
338 p_vectors->InsertNextValue(dataPayload[i](1,2));
339 p_vectors->InsertNextValue(dataPayload[i](2,0));
340 p_vectors->InsertNextValue(dataPayload[i](2,1));
341 p_vectors->InsertNextValue(dataPayload[i](2,2));
345 vtkCellData* p_cell_data = mpVtkUnstructedMesh->GetCellData();
346 p_cell_data->AddArray(p_vectors);
351 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
354 vtkDoubleArray* p_scalars = vtkDoubleArray::New();
355 p_scalars->SetName(dataName.c_str());
357 if (mWriteParallelFiles && this->mpDistributedMesh !=
nullptr)
364 assert( dataPayload.size() == this->mpDistributedMesh->GetNumLocalNodes() );
365 dataPayload.resize( this->mpDistributedMesh->GetNumLocalNodes() + this->mpDistributedMesh->GetNumHaloNodes() );
374 unsigned number_of_nodes_to_send = mNodesToSendPerProcess[send_to].size();
375 unsigned number_of_nodes_to_receive = mNodesToReceivePerProcess[receive_from].size();
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]);
380 for (
unsigned node = 0; node < number_of_nodes_to_send; node++)
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];
391 ret = MPI_Sendrecv(send_data.get(), number_of_nodes_to_send,
394 receive_data.get(), number_of_nodes_to_receive,
397 PETSC_COMM_WORLD, &status);
399 assert ( ret == MPI_SUCCESS );
403 for (
unsigned node = 0; node < number_of_nodes_to_receive; node++ )
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];
414 for (
unsigned i=0; i<dataPayload.size(); i++)
416 p_scalars->InsertNextValue(dataPayload[i]);
419 vtkPointData* p_point_data = mpVtkUnstructedMesh->GetPointData();
420 p_point_data->AddArray(p_scalars);
425 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
428 vtkDoubleArray* p_vectors = vtkDoubleArray::New();
429 p_vectors->SetName(dataName.c_str());
431 if (mWriteParallelFiles)
438 assert( dataPayload.size() == this->mpDistributedMesh->GetNumLocalNodes() );
439 dataPayload.resize( this->mpDistributedMesh->GetNumLocalNodes() + this->mpDistributedMesh->GetNumHaloNodes() );
447 unsigned number_of_nodes_to_send = mNodesToSendPerProcess[send_to].size();
448 unsigned number_of_nodes_to_receive = mNodesToReceivePerProcess[receive_from].size();
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]);
453 for (
unsigned node = 0; node < number_of_nodes_to_send; node++)
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++)
460 send_data[ node*SPACE_DIM + j ] = dataPayload[local_node_index][j];
466 ret = MPI_Sendrecv(send_data.get(), number_of_nodes_to_send * SPACE_DIM,
469 receive_data.get(), number_of_nodes_to_receive * SPACE_DIM,
472 PETSC_COMM_WORLD, &status);
474 assert ( ret == MPI_SUCCESS );
477 for (
unsigned node = 0; node < number_of_nodes_to_receive; node++ )
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++)
484 dataPayload[halo_index][j] = receive_data[ node*SPACE_DIM + j ];
490 p_vectors->SetNumberOfComponents(3);
491 for (
unsigned i=0; i<dataPayload.size(); i++)
493 for (
unsigned j=0; j<SPACE_DIM; j++)
495 p_vectors->InsertNextValue(dataPayload[i][j]);
498 for (
unsigned j=SPACE_DIM; j<3; j++)
500 p_vectors->InsertNextValue(0.0);
504 vtkPointData* p_point_data = mpVtkUnstructedMesh->GetPointData();
505 p_point_data->AddArray(p_vectors);
509 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
512 assert(SPACE_DIM != 1);
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++)
521 p_vectors->InsertNextValue(dataPayload[i](0,0));
522 p_vectors->InsertNextValue(dataPayload[i](0,1));
523 p_vectors->InsertNextValue(dataPayload[i](1,0));
524 p_vectors->InsertNextValue(dataPayload[i](1,1));
526 else if (SPACE_DIM == 3)
528 p_vectors->InsertNextValue(dataPayload[i](0,0));
529 p_vectors->InsertNextValue(dataPayload[i](0,1));
530 p_vectors->InsertNextValue(dataPayload[i](0,2));
531 p_vectors->InsertNextValue(dataPayload[i](1,0));
532 p_vectors->InsertNextValue(dataPayload[i](1,1));
533 p_vectors->InsertNextValue(dataPayload[i](1,2));
534 p_vectors->InsertNextValue(dataPayload[i](2,0));
535 p_vectors->InsertNextValue(dataPayload[i](2,1));
536 p_vectors->InsertNextValue(dataPayload[i](2,2));
540 vtkPointData* p_point_data = mpVtkUnstructedMesh->GetPointData();
541 p_point_data->AddArray(p_vectors);
545 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
552 if (this->mpDistributedMesh ==
nullptr && mpNodesOnlyMesh ==
nullptr)
554 EXCEPTION(
"Cannot write parallel files using a sequential mesh");
562 mWriteParallelFiles =
true;
574 mGlobalToNodeIndexMap[node_iter->GetIndex()] = index;
579 if (this->mpDistributedMesh)
582 halo_iter != this->mpDistributedMesh->GetHaloNodeIteratorEnd();
585 mGlobalToNodeIndexMap[(*halo_iter)->GetIndex()] = index;
590 this->mpDistributedMesh->CalculateNodeExchange( mNodesToSendPerProcess, mNodesToReceivePerProcess );
595 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
598 bool keepOriginalElementIndexing)
604 if (
PetscTools::IsSequential() || !mWriteParallelFiles || (this->mpDistributedMesh ==
nullptr && mpNodesOnlyMesh ==
nullptr))
611 vtkPoints* p_pts = vtkPoints::New(VTK_DOUBLE);
612 p_pts->GetData()->SetName(
"Vertex positions");
619 c_vector<double, SPACE_DIM> current_item = node_iter->rGetLocation();
622 p_pts->InsertNextPoint(current_item[0], current_item[1], current_item[2]);
624 else if (SPACE_DIM == 2)
626 p_pts->InsertNextPoint(current_item[0], current_item[1], 0.0);
630 p_pts->InsertNextPoint(current_item[0], 0.0, 0.0);
635 if (this->mpDistributedMesh)
638 halo_iter != this->mpDistributedMesh->GetHaloNodeIteratorEnd();
641 c_vector<double, SPACE_DIM> current_item = (*halo_iter)->rGetLocation();
644 p_pts->InsertNextPoint(current_item[0], current_item[1], current_item[2]);
646 else if (SPACE_DIM == 2)
648 p_pts->InsertNextPoint(current_item[0], current_item[1], 0.0);
652 p_pts->InsertNextPoint(current_item[0], 0.0, 0.0);
657 mpVtkUnstructedMesh->SetPoints(p_pts);
665 vtkCell* p_cell=
nullptr;
667 if (ELEMENT_DIM == 3)
669 p_cell = vtkTetra::New();
671 else if (ELEMENT_DIM == 2)
673 p_cell = vtkTriangle::New();
677 p_cell = vtkLine::New();
679 vtkIdList* p_cell_id_list = p_cell->GetPointIds();
680 for (
unsigned j = 0; j < ELEMENT_DIM+1; ++j)
682 unsigned global_node_index = elem_iter->GetNodeGlobalIndex(j);
683 p_cell_id_list->SetId(j, mGlobalToNodeIndexMap[global_node_index]);
685 mpVtkUnstructedMesh->InsertNextCell(p_cell->GetCellType(), p_cell_id_list);
689 if (this->mpMixedMesh )
695 elem_iter != this->mpMixedMesh->GetCableElementIteratorEnd();
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)
703 unsigned global_node_index = (*elem_iter)->GetNodeGlobalIndex(j);
704 p_cell_id_list->SetId(j, mGlobalToNodeIndexMap[global_node_index]);
706 mpVtkUnstructedMesh->InsertNextCell(p_cell->GetCellType(), p_cell_id_list);
709 AddCellData(
"Cable radius", radii);
716 assert(mpVtkUnstructedMesh->CheckAttributes() == 0);
717 vtkXMLPUnstructuredGridWriter* p_writer = vtkXMLPUnstructuredGridWriter::New();
719 p_writer->SetDataModeToBinary();
727 #if VTK_MAJOR_VERSION >= 6
728 p_writer->SetInputData(mpVtkUnstructedMesh);
730 p_writer->SetInput(mpVtkUnstructedMesh);
732 std::string pvtk_file_name = this->mpOutputFileHandler->GetOutputDirectoryFullPath() + this->mBaseName+
".pvtu";
733 p_writer->SetFileName(pvtk_file_name.c_str());
740 std::stringstream filepath;
742 AddProvenance(filepath.str());
746 AddProvenance(this->mBaseName+
".pvtu");
ElementIterator GetElementIteratorBegin(bool skipDeletedElements=true)
ElementIterator GetElementIteratorEnd()
#define EXCEPTION(message)
vtkUnstructuredGrid * mpVtkUnstructedMesh
std::vector< Node< SPACE_DIM > * >::const_iterator HaloNodeIterator
void SetParallelFiles(AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > &rMesh)
void AddCellData(std::string name, std::vector< double > data)
void AddProvenance(std::string fileName)
NodeIterator GetNodeIteratorEnd()
std::vector< unsigned > NodeIndices
void WriteFilesUsingMesh(AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > &rMesh, bool keepOriginalElementIndexing=true)
unsigned GetNumLocalElements() const
void AddTensorCellData(std::string name, std::vector< c_vector< double, SPACE_DIM *(SPACE_DIM+1)/2 > > data)
void AddTensorPointData(std::string name, std::vector< c_matrix< double, SPACE_DIM, SPACE_DIM > > data)
NodeIterator GetNodeIteratorBegin(bool skipDeletedNodes=true)
static std::string GetProvenanceString()
virtual void WriteFilesUsingMesh(AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > &rMesh, bool keepOriginalElementIndexing=true)
VtkMeshWriter(const std::string &rDirectory, const std::string &rBaseName, const bool &rCleanDirectory=true)
void AddPointData(std::string name, std::vector< double > data)
std::vector< Element< 1, SPACE_DIM > * >::const_iterator CableElementIterator