Chaste  Release::3.4
VtkMeshWriter.cpp
1 /*
2 
3 Copyright (c) 2005-2016, University of Oxford.
4 All rights reserved.
5 
6 University of Oxford means the Chancellor, Masters and Scholars of the
7 University of Oxford, having an administrative office at Wellington
8 Square, Oxford OX1 2JD, UK.
9 
10 This file is part of Chaste.
11 
12 Redistribution and use in source and binary forms, with or without
13 modification, 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 
23 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
24 AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
25 IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
26 ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
27 LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
28 CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
29 GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
30 HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
31 LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
32 OF 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
50 template <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 
63 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
65 {
66  mpVtkUnstructedMesh->Delete(); // Reference counted
67 }
68 
69 template <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=NULL;
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 
172 template <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 
183 template <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  //Uninitialised stuff arises (see #1079), but you can remove
197  //valgrind problems by removing compression:
198  // **** REMOVE WITH CAUTION *****
199  p_writer->SetCompressor(NULL);
200  // **** REMOVE WITH CAUTION *****
201  std::string vtk_file_name = this->mpOutputFileHandler->GetOutputDirectoryFullPath() + this->mBaseName+".vtu";
202  p_writer->SetFileName(vtk_file_name.c_str());
203  //p_writer->PrintSelf(std::cout, vtkIndent());
204  p_writer->Write();
205  p_writer->Delete(); //Reference counted
206  }
207 
208  AddProvenance(this->mBaseName + ".vtu");
209 }
210 
211 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
212 void VtkMeshWriter<ELEMENT_DIM,SPACE_DIM>::AddCellData(std::string dataName, std::vector<double> dataPayload)
213 {
214  vtkDoubleArray* p_scalars = vtkDoubleArray::New();
215  p_scalars->SetName(dataName.c_str());
216  for (unsigned i=0; i<dataPayload.size(); i++)
217  {
218  p_scalars->InsertNextValue(dataPayload[i]);
219  }
220 
221  vtkCellData* p_cell_data = mpVtkUnstructedMesh->GetCellData();
222  p_cell_data->AddArray(p_scalars);
223  p_scalars->Delete(); //Reference counted
224 }
225 
226 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
228 {
229  unsigned num_cell_arrays = mpVtkUnstructedMesh->GetCellData()->GetNumberOfArrays();
230  for (unsigned i = 0; i < num_cell_arrays; i++)
231  {
232  vtkDataArray* array = mpVtkUnstructedMesh->GetCellData()->GetArray(i);
233 
234  //Check data was the correct size before the cables were added
235  unsigned num_cable_pads = this->GetNumCableElements();
236  if (mWriteParallelFiles)
237  {
238  assert((unsigned)array->GetNumberOfTuples() == this->mpDistributedMesh->GetNumLocalElements());
239  num_cable_pads = this->mpMixedMesh->GetNumLocalCableElements();
240  }
241  else
242  {
243  assert((unsigned)array->GetNumberOfTuples() == this->GetNumElements());
244  }
245 
246  //Check that tuples of size 3 will be big enough for padding the rest of the data
247  assert(array->GetNumberOfComponents() <= 3);
248  double null_data[3] = {0.0, 0.0, 0.0};
249 
250  //Pad data
251  for (unsigned new_index = 0; new_index < num_cable_pads; new_index++)
252  {
253  array->InsertNextTuple(null_data);
254  }
255  }
256 }
257 
258 
259 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
260 void VtkMeshWriter<ELEMENT_DIM,SPACE_DIM>::AddCellData(std::string dataName, std::vector<c_vector<double, SPACE_DIM> > dataPayload)
261 {
262  vtkDoubleArray* p_vectors = vtkDoubleArray::New();
263  p_vectors->SetName(dataName.c_str());
264  p_vectors->SetNumberOfComponents(3);
265  for (unsigned i=0; i<dataPayload.size(); i++)
266  {
267  for (unsigned j=0; j<SPACE_DIM; j++)
268  {
269  p_vectors->InsertNextValue(dataPayload[i][j]);
270  }
271  //When SPACE_DIM<3, then pad
272  for (unsigned j=SPACE_DIM; j<3; j++)
273  {
274  p_vectors->InsertNextValue(0.0);
275  }
276  }
277 
278  vtkCellData* p_cell_data = mpVtkUnstructedMesh->GetCellData();
279  p_cell_data->AddArray(p_vectors);
280  p_vectors->Delete(); //Reference counted
281 }
282 
283 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
284 void VtkMeshWriter<ELEMENT_DIM,SPACE_DIM>::AddTensorCellData(std::string dataName, std::vector<c_vector<double,SPACE_DIM*(SPACE_DIM+1)/2> > dataPayload)
285 {
286  assert(SPACE_DIM != 1);
287 
288  vtkDoubleArray* p_vectors = vtkDoubleArray::New();
289  p_vectors->SetName(dataName.c_str());
290  p_vectors->SetNumberOfComponents(SPACE_DIM*SPACE_DIM);
291  for (unsigned i=0; i<dataPayload.size(); i++)
292  {
293  if(SPACE_DIM == 2)
294  {
295  p_vectors->InsertNextValue(dataPayload[i](0)); //a11
296  p_vectors->InsertNextValue(dataPayload[i](1)); //a12
297  p_vectors->InsertNextValue(dataPayload[i](1)); //a21
298  p_vectors->InsertNextValue(dataPayload[i](2)); //a22
299  }
300  else if (SPACE_DIM == 3)
301  {
302  p_vectors->InsertNextValue(dataPayload[i](0)); //a11
303  p_vectors->InsertNextValue(dataPayload[i](1)); //a12
304  p_vectors->InsertNextValue(dataPayload[i](2)); //a13
305  p_vectors->InsertNextValue(dataPayload[i](1)); //a21
306  p_vectors->InsertNextValue(dataPayload[i](3)); //a22
307  p_vectors->InsertNextValue(dataPayload[i](4)); //a23
308  p_vectors->InsertNextValue(dataPayload[i](2)); //a31
309  p_vectors->InsertNextValue(dataPayload[i](4)); //a32
310  p_vectors->InsertNextValue(dataPayload[i](5)); //a33
311  }
312  }
313 
314  vtkCellData* p_cell_data = mpVtkUnstructedMesh->GetCellData();
315  p_cell_data->AddArray(p_vectors);
316  p_vectors->Delete(); //Reference counted
317 }
318 
319 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
320 void VtkMeshWriter<ELEMENT_DIM,SPACE_DIM>::AddTensorCellData(std::string dataName, std::vector<c_matrix<double,SPACE_DIM,SPACE_DIM> > dataPayload)
321 {
322  assert(SPACE_DIM != 1);
323 
324  vtkDoubleArray* p_vectors = vtkDoubleArray::New();
325  p_vectors->SetName(dataName.c_str());
326  p_vectors->SetNumberOfComponents(SPACE_DIM*SPACE_DIM);
327  for (unsigned i=0; i<dataPayload.size(); i++)
328  {
329  if(SPACE_DIM == 2)
330  {
331  p_vectors->InsertNextValue(dataPayload[i](0,0)); //a11
332  p_vectors->InsertNextValue(dataPayload[i](0,1)); //a12
333  p_vectors->InsertNextValue(dataPayload[i](1,0)); //a21
334  p_vectors->InsertNextValue(dataPayload[i](1,1)); //a22
335  }
336  else if (SPACE_DIM == 3)
337  {
338  p_vectors->InsertNextValue(dataPayload[i](0,0)); //a11
339  p_vectors->InsertNextValue(dataPayload[i](0,1)); //a12
340  p_vectors->InsertNextValue(dataPayload[i](0,2)); //a13
341  p_vectors->InsertNextValue(dataPayload[i](1,0)); //a21
342  p_vectors->InsertNextValue(dataPayload[i](1,1)); //a22
343  p_vectors->InsertNextValue(dataPayload[i](1,2)); //a23
344  p_vectors->InsertNextValue(dataPayload[i](2,0)); //a31
345  p_vectors->InsertNextValue(dataPayload[i](2,1)); //a32
346  p_vectors->InsertNextValue(dataPayload[i](2,2)); //a33
347  }
348  }
349 
350  vtkCellData* p_cell_data = mpVtkUnstructedMesh->GetCellData();
351  p_cell_data->AddArray(p_vectors);
352  p_vectors->Delete(); //Reference counted
353 }
354 
355 
356 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
357 void VtkMeshWriter<ELEMENT_DIM,SPACE_DIM>::AddPointData(std::string dataName, std::vector<double> dataPayload)
358 {
359  vtkDoubleArray* p_scalars = vtkDoubleArray::New();
360  p_scalars->SetName(dataName.c_str());
361 
362  if (mWriteParallelFiles && this->mpDistributedMesh != NULL)
363  {
364  // In parallel, the vector we pass will only contain the values from the privately owned nodes.
365  // To get the values from the halo nodes (which will be inserted at the end of the vector we need to
366  // communicate with the equivalent vectors on other processes.
367 
368  // resize the payload data to include halos
369  assert( dataPayload.size() == this->mpDistributedMesh->GetNumLocalNodes() );
370  dataPayload.resize( this->mpDistributedMesh->GetNumLocalNodes() + this->mpDistributedMesh->GetNumHaloNodes() );
371 
372 
373  // then do the communication
374  for ( unsigned rank_offset = 1; rank_offset < PetscTools::GetNumProcs(); rank_offset++ )
375  {
376  unsigned send_to = (PetscTools::GetMyRank() + rank_offset) % (PetscTools::GetNumProcs());
377  unsigned receive_from = (PetscTools::GetMyRank() + PetscTools::GetNumProcs()- rank_offset ) % (PetscTools::GetNumProcs());
378 
379  unsigned number_of_nodes_to_send = mNodesToSendPerProcess[send_to].size();
380  unsigned number_of_nodes_to_receive = mNodesToReceivePerProcess[receive_from].size();
381 
382  boost::scoped_array<double> send_data(new double[number_of_nodes_to_send]);
383  boost::scoped_array<double> receive_data(new double[number_of_nodes_to_receive]);
384  // Pack
385  for (unsigned node = 0; node < number_of_nodes_to_send; node++)
386  {
387  unsigned global_node_index = mNodesToSendPerProcess[send_to][node];
388  unsigned local_node_index = global_node_index
389  - this->mpDistributedMesh->GetDistributedVectorFactory()->GetLow();
390  send_data[node] = dataPayload[local_node_index];
391  }
392  {
393  // Send
394  int ret;
395  MPI_Status status;
396  ret = MPI_Sendrecv(send_data.get(), number_of_nodes_to_send,
397  MPI_DOUBLE,
398  send_to, 0,
399  receive_data.get(), number_of_nodes_to_receive,
400  MPI_DOUBLE,
401  receive_from, 0,
402  PETSC_COMM_WORLD, &status);
403  UNUSED_OPT(ret);
404  assert ( ret == MPI_SUCCESS );
405  }
406 
407  // Unpack
408  for ( unsigned node = 0; node < number_of_nodes_to_receive; node++ )
409  {
410  unsigned global_node_index = mNodesToReceivePerProcess[receive_from][node];
411  unsigned halo_index = mGlobalToNodeIndexMap[global_node_index];
412  assert( halo_index >= this->mpDistributedMesh->GetNumLocalNodes() );
413  dataPayload[halo_index] = receive_data[node];
414  }
415 
416  }
417  }
418 
419  for (unsigned i=0; i<dataPayload.size(); i++)
420  {
421  p_scalars->InsertNextValue(dataPayload[i]);
422  }
423 
424  vtkPointData* p_point_data = mpVtkUnstructedMesh->GetPointData();
425  p_point_data->AddArray(p_scalars);
426  p_scalars->Delete(); //Reference counted
427 }
428 
429 
430 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
431 void VtkMeshWriter<ELEMENT_DIM,SPACE_DIM>::AddPointData(std::string dataName, std::vector<c_vector<double, SPACE_DIM> > dataPayload)
432 {
433  vtkDoubleArray* p_vectors = vtkDoubleArray::New();
434  p_vectors->SetName(dataName.c_str());
435 
436  if (mWriteParallelFiles)
437  {
438  // In parallel, the vector we pass will only contain the values from the privately owned nodes.
439  // To get the values from the halo nodes (which will be inserted at the end of the vector we need to
440  // communicate with the equivalent vectors on other processes.
441 
442  // resize the payload data to include halos
443  assert( dataPayload.size() == this->mpDistributedMesh->GetNumLocalNodes() );
444  dataPayload.resize( this->mpDistributedMesh->GetNumLocalNodes() + this->mpDistributedMesh->GetNumHaloNodes() );
445 
446  // then do the communication
447  for ( unsigned rank_offset = 1; rank_offset < PetscTools::GetNumProcs(); rank_offset++ )
448  {
449  unsigned send_to = (PetscTools::GetMyRank() + rank_offset) % (PetscTools::GetNumProcs());
450  unsigned receive_from = (PetscTools::GetMyRank() + PetscTools::GetNumProcs()- rank_offset ) % (PetscTools::GetNumProcs());
451 
452  unsigned number_of_nodes_to_send = mNodesToSendPerProcess[send_to].size();
453  unsigned number_of_nodes_to_receive = mNodesToReceivePerProcess[receive_from].size();
454 
455  boost::scoped_array<double> send_data(new double[number_of_nodes_to_send * SPACE_DIM]);
456  boost::scoped_array<double> receive_data(new double[number_of_nodes_to_receive * SPACE_DIM]);
457 
458  for (unsigned node = 0; node < number_of_nodes_to_send; node++)
459  {
460  unsigned global_node_index = mNodesToSendPerProcess[send_to][node];
461  unsigned local_node_index = global_node_index
462  - this->mpDistributedMesh->GetDistributedVectorFactory()->GetLow();
463  for (unsigned j=0; j<SPACE_DIM; j++)
464  {
465  send_data[ node*SPACE_DIM + j ] = dataPayload[local_node_index][j];
466  }
467  }
468 
469  int ret;
470  MPI_Status status;
471  ret = MPI_Sendrecv(send_data.get(), number_of_nodes_to_send * SPACE_DIM,
472  MPI_DOUBLE,
473  send_to, 0,
474  receive_data.get(), number_of_nodes_to_receive * SPACE_DIM,
475  MPI_DOUBLE,
476  receive_from, 0,
477  PETSC_COMM_WORLD, &status);
478  UNUSED_OPT(ret);
479  assert ( ret == MPI_SUCCESS );
480 
481  // Unpack
482  for ( unsigned node = 0; node < number_of_nodes_to_receive; node++ )
483  {
484  unsigned global_node_index = mNodesToReceivePerProcess[receive_from][node];
485  unsigned halo_index = mGlobalToNodeIndexMap[global_node_index];
486  assert( halo_index >= this->mpDistributedMesh->GetNumLocalNodes() );
487  for (unsigned j=0; j<SPACE_DIM; j++)
488  {
489  dataPayload[halo_index][j] = receive_data[ node*SPACE_DIM + j ];
490  }
491  }
492  }
493  }
494 
495  p_vectors->SetNumberOfComponents(3);
496  for (unsigned i=0; i<dataPayload.size(); i++)
497  {
498  for (unsigned j=0; j<SPACE_DIM; j++)
499  {
500  p_vectors->InsertNextValue(dataPayload[i][j]);
501  }
502  //When SPACE_DIM<3, then pad
503  for (unsigned j=SPACE_DIM; j<3; j++)
504  {
505  p_vectors->InsertNextValue(0.0);
506  }
507  }
508 
509  vtkPointData* p_point_data = mpVtkUnstructedMesh->GetPointData();
510  p_point_data->AddArray(p_vectors);
511  p_vectors->Delete(); //Reference counted
512 }
513 
514 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
515 void VtkMeshWriter<ELEMENT_DIM,SPACE_DIM>::AddTensorPointData(std::string dataName, std::vector<c_matrix<double,SPACE_DIM,SPACE_DIM> > dataPayload)
516 {
517  assert(SPACE_DIM != 1);
518 
519  vtkDoubleArray* p_vectors = vtkDoubleArray::New();
520  p_vectors->SetName(dataName.c_str());
521  p_vectors->SetNumberOfComponents(SPACE_DIM*SPACE_DIM);
522  for (unsigned i=0; i<dataPayload.size(); i++)
523  {
524  if(SPACE_DIM == 2)
525  {
526  p_vectors->InsertNextValue(dataPayload[i](0,0)); //a11
527  p_vectors->InsertNextValue(dataPayload[i](0,1)); //a12
528  p_vectors->InsertNextValue(dataPayload[i](1,0)); //a21
529  p_vectors->InsertNextValue(dataPayload[i](1,1)); //a22
530  }
531  else if (SPACE_DIM == 3)
532  {
533  p_vectors->InsertNextValue(dataPayload[i](0,0)); //a11
534  p_vectors->InsertNextValue(dataPayload[i](0,1)); //a12
535  p_vectors->InsertNextValue(dataPayload[i](0,2)); //a13
536  p_vectors->InsertNextValue(dataPayload[i](1,0)); //a21
537  p_vectors->InsertNextValue(dataPayload[i](1,1)); //a22
538  p_vectors->InsertNextValue(dataPayload[i](1,2)); //a23
539  p_vectors->InsertNextValue(dataPayload[i](2,0)); //a31
540  p_vectors->InsertNextValue(dataPayload[i](2,1)); //a32
541  p_vectors->InsertNextValue(dataPayload[i](2,2)); //a33
542  }
543  }
544 
545  vtkPointData* p_point_data = mpVtkUnstructedMesh->GetPointData();
546  p_point_data->AddArray(p_vectors);
547  p_vectors->Delete(); //Reference counted
548 }
549 
550 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
552 {
553  //Have we got a distributed mesh?
554  this->mpDistributedMesh = dynamic_cast<DistributedTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* >(&rMesh);
555  mpNodesOnlyMesh = dynamic_cast<NodesOnlyMesh<SPACE_DIM>* >(&rMesh);
556 
557  if (this->mpDistributedMesh == NULL && mpNodesOnlyMesh == NULL)
558  {
559  EXCEPTION("Cannot write parallel files using a sequential mesh");
560  }
561 
563  {
564  return; // mWriteParallelFiles is not set sequentially (so we don't set up data exchange machinery)
565  }
566 
567  mWriteParallelFiles = true;
568 
569  // Populate the global to node index map (as this will be required to add point data)
570 
571  //Node index that we are writing to VTK (index into mNodes and mHaloNodes as if they were concatenated)
572  unsigned index = 0;
573 
574  // Owned nodes
576  node_iter != rMesh.GetNodeIteratorEnd();
577  ++node_iter)
578  {
579  mGlobalToNodeIndexMap[node_iter->GetIndex()] = index;
580  index++;
581  }
582 
583  // Halo nodes
584  if (this->mpDistributedMesh)
585  {
586  for (typename DistributedTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::HaloNodeIterator halo_iter=this->mpDistributedMesh->GetHaloNodeIteratorBegin();
587  halo_iter != this->mpDistributedMesh->GetHaloNodeIteratorEnd();
588  ++halo_iter)
589  {
590  mGlobalToNodeIndexMap[(*halo_iter)->GetIndex()] = index;
591  index++;
592  }
593 
594  //Calculate the halo exchange so that node-wise payloads can be communicated
595  this->mpDistributedMesh->CalculateNodeExchange( mNodesToSendPerProcess, mNodesToReceivePerProcess );
596  }
597 }
598 
600 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
603  bool keepOriginalElementIndexing)
604 {
605  //Have we got a parallel mesh?
606  this->mpDistributedMesh = dynamic_cast<DistributedTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* >(&rMesh);
607  this->mpMixedMesh = dynamic_cast<MixedDimensionMesh<ELEMENT_DIM,SPACE_DIM>* >(&rMesh);
608 
609  if ( PetscTools::IsSequential() || !mWriteParallelFiles || (this->mpDistributedMesh == NULL && mpNodesOnlyMesh == NULL) )
610  {
612  }
613  else
614  {
615  //Make the local mesh into a VtkMesh
616  vtkPoints* p_pts = vtkPoints::New(VTK_DOUBLE);
617  p_pts->GetData()->SetName("Vertex positions");
618 
619  // Owned nodes
621  node_iter != rMesh.GetNodeIteratorEnd();
622  ++node_iter)
623  {
624  c_vector<double, SPACE_DIM> current_item = node_iter->rGetLocation();
625  if (SPACE_DIM == 3)
626  {
627  p_pts->InsertNextPoint(current_item[0], current_item[1], current_item[2]);
628  }
629  else if (SPACE_DIM == 2)
630  {
631  p_pts->InsertNextPoint(current_item[0], current_item[1], 0.0);
632  }
633  else // (SPACE_DIM == 1)
634  {
635  p_pts->InsertNextPoint(current_item[0], 0.0, 0.0);
636  }
637  }
638 
639  // Halo nodes
640  if (this->mpDistributedMesh)
641  {
642  for (typename DistributedTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::HaloNodeIterator halo_iter=this->mpDistributedMesh->GetHaloNodeIteratorBegin();
643  halo_iter != this->mpDistributedMesh->GetHaloNodeIteratorEnd();
644  ++halo_iter)
645  {
646  c_vector<double, SPACE_DIM> current_item = (*halo_iter)->rGetLocation();
647  if (SPACE_DIM == 3)
648  {
649  p_pts->InsertNextPoint(current_item[0], current_item[1], current_item[2]);
650  }
651  else if (SPACE_DIM == 2)
652  {
653  p_pts->InsertNextPoint(current_item[0], current_item[1], 0.0);
654  }
655  else // (SPACE_DIM == 1)
656  {
657  p_pts->InsertNextPoint(current_item[0], 0.0, 0.0);
658  }
659  }
660  }
661 
662  mpVtkUnstructedMesh->SetPoints(p_pts);
663  p_pts->Delete(); //Reference counted
664 
666  elem_iter != rMesh.GetElementIteratorEnd();
667  ++elem_iter)
668  {
669 
670  vtkCell* p_cell=NULL;
672  if (ELEMENT_DIM == 3)
673  {
674  p_cell = vtkTetra::New();
675  }
676  else if (ELEMENT_DIM == 2)
677  {
678  p_cell = vtkTriangle::New();
679  }
680  else //(ELEMENT_DIM == 1)
681  {
682  p_cell = vtkLine::New();
683  }
684  vtkIdList* p_cell_id_list = p_cell->GetPointIds();
685  for (unsigned j = 0; j < ELEMENT_DIM+1; ++j)
686  {
687  unsigned global_node_index = elem_iter->GetNodeGlobalIndex(j);
688  p_cell_id_list->SetId(j, mGlobalToNodeIndexMap[global_node_index]);
689  }
690  mpVtkUnstructedMesh->InsertNextCell(p_cell->GetCellType(), p_cell_id_list);
691  p_cell->Delete(); //Reference counted
692  }
693  //If necessary, construct cables
694  if (this->mpMixedMesh )
695  {
696  AugmentCellData();
697  //Make a blank cell radius data for the regular elements
698  std::vector<double> radii(this->mpMixedMesh->GetNumLocalElements(), 0.0);
699  for (typename MixedDimensionMesh<ELEMENT_DIM,SPACE_DIM>::CableElementIterator elem_iter = this->mpMixedMesh->GetCableElementIteratorBegin();
700  elem_iter != this->mpMixedMesh->GetCableElementIteratorEnd();
701  ++elem_iter)
702  {
703  radii.push_back((*elem_iter)->GetAttribute());
704  vtkCell* p_cell=vtkLine::New();
705  vtkIdList* p_cell_id_list = p_cell->GetPointIds();
706  for (unsigned j = 0; j < 2; ++j)
707  {
708  unsigned global_node_index = (*elem_iter)->GetNodeGlobalIndex(j);
709  p_cell_id_list->SetId(j, mGlobalToNodeIndexMap[global_node_index]);
710  }
711  mpVtkUnstructedMesh->InsertNextCell(p_cell->GetCellType(), p_cell_id_list);
712  p_cell->Delete(); //Reference counted
713  }
714  AddCellData("Cable radius", radii);
715  }
716 
717 
718  //This block is to guard the mesh writers (vtkXMLPUnstructuredGridWriter) so that they
719  //go out of scope, flush buffers and close files
720  {
721  assert(mpVtkUnstructedMesh->CheckAttributes() == 0);
722  vtkXMLPUnstructuredGridWriter* p_writer = vtkXMLPUnstructuredGridWriter::New();
723 
724  p_writer->SetDataModeToBinary();
725 
726  p_writer->SetNumberOfPieces(PetscTools::GetNumProcs());
727  //p_writer->SetGhostLevel(-1);
728  p_writer->SetStartPiece(PetscTools::GetMyRank());
729  p_writer->SetEndPiece(PetscTools::GetMyRank());
730 
731 
732 #if VTK_MAJOR_VERSION >= 6
733  p_writer->SetInputData(mpVtkUnstructedMesh);
734 #else
735  p_writer->SetInput(mpVtkUnstructedMesh);
736 #endif
737  //Uninitialised stuff arises (see #1079), but you can remove
738  //valgrind problems by removing compression:
739  // **** REMOVE WITH CAUTION *****
740  p_writer->SetCompressor(NULL);
741  // **** REMOVE WITH CAUTION *****
742  std::string pvtk_file_name = this->mpOutputFileHandler->GetOutputDirectoryFullPath() + this->mBaseName+ ".pvtu";
743  p_writer->SetFileName(pvtk_file_name.c_str());
744  //p_writer->PrintSelf(std::cout, vtkIndent());
745  p_writer->Write();
746  p_writer->Delete(); //Reference counted
747  }
748 
749  // Add provenance to the individual files
750  std::stringstream filepath;
751  filepath << this->mBaseName << "_" << PetscTools::GetMyRank() << ".vtu";
752  AddProvenance(filepath.str());
754  if (PetscTools::AmMaster())
755  {
756  AddProvenance(this->mBaseName+ ".pvtu");
757  }
758  }
759 }
760 
762 // Explicit instantiation
764 
765 template class VtkMeshWriter<1,1>;
766 template class VtkMeshWriter<1,2>;
767 template class VtkMeshWriter<1,3>;
768 template class VtkMeshWriter<2,2>; // Actually used
769 template class VtkMeshWriter<2,3>;
770 template class VtkMeshWriter<3,3>; // Actually used
771 
772 #endif //CHASTE_VTK
ElementIterator GetElementIteratorBegin(bool skipDeletedElements=true)
#define EXCEPTION(message)
Definition: Exception.hpp:143
vtkUnstructuredGrid * mpVtkUnstructedMesh
std::vector< Node< SPACE_DIM > * >::const_iterator HaloNodeIterator
void SetParallelFiles(AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > &rMesh)
static bool AmMaster()
Definition: PetscTools.cpp:120
void AddCellData(std::string name, std::vector< double > data)
void AddProvenance(std::string fileName)
NodeIterator GetNodeIteratorEnd()
void AugmentCellData()
std::vector< unsigned > NodeIndices
void WriteFilesUsingMesh(AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > &rMesh, bool keepOriginalElementIndexing=true)
static bool IsSequential()
Definition: PetscTools.cpp:91
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)
virtual ~VtkMeshWriter()
#define UNUSED_OPT(var)
Definition: Exception.hpp:216
static std::string GetProvenanceString()
static unsigned GetNumProcs()
Definition: PetscTools.cpp:108
static unsigned GetMyRank()
Definition: PetscTools.cpp:114
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