VtkMeshWriter.cpp

00001 /*
00002 
00003 Copyright (c) 2005-2015, University of Oxford.
00004 All rights reserved.
00005 
00006 University of Oxford means the Chancellor, Masters and Scholars of the
00007 University of Oxford, having an administrative office at Wellington
00008 Square, Oxford OX1 2JD, UK.
00009 
00010 This file is part of Chaste.
00011 
00012 Redistribution and use in source and binary forms, with or without
00013 modification, are permitted provided that the following conditions are met:
00014  * Redistributions of source code must retain the above copyright notice,
00015    this list of conditions and the following disclaimer.
00016  * Redistributions in binary form must reproduce the above copyright notice,
00017    this list of conditions and the following disclaimer in the documentation
00018    and/or other materials provided with the distribution.
00019  * Neither the name of the University of Oxford nor the names of its
00020    contributors may be used to endorse or promote products derived from this
00021    software without specific prior written permission.
00022 
00023 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
00024 AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
00025 IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
00026 ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
00027 LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
00028 CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
00029 GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
00030 HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
00031 LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
00032 OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
00033 
00034 */
00035 
00036 #include <boost/scoped_array.hpp>
00037 #include "VtkMeshWriter.hpp"
00038 #include "DistributedTetrahedralMesh.hpp"
00039 #include "MixedDimensionMesh.hpp"
00040 #include "NodesOnlyMesh.hpp"
00041 
00042 #ifdef CHASTE_VTK
00043 #include "vtkQuadraticTetra.h"
00044 #include "vtkQuadraticTriangle.h"
00045 
00046 
00048 // Implementation
00050 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00051 VtkMeshWriter<ELEMENT_DIM, SPACE_DIM>::VtkMeshWriter(const std::string& rDirectory,
00052                      const std::string& rBaseName,
00053                      const bool& rCleanDirectory)
00054     : AbstractTetrahedralMeshWriter<ELEMENT_DIM, SPACE_DIM>(rDirectory, rBaseName, rCleanDirectory),
00055       mWriteParallelFiles(false)
00056 {
00057     this->mIndexFromZero = true;
00058 
00059     // Dubious, since we shouldn't yet know what any details of the mesh are.
00060     mpVtkUnstructedMesh = vtkUnstructuredGrid::New();
00061 }
00062 
00063 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00064 VtkMeshWriter<ELEMENT_DIM,SPACE_DIM>::~VtkMeshWriter()
00065 {
00066     mpVtkUnstructedMesh->Delete(); // Reference counted
00067 }
00068 
00069 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00070 void VtkMeshWriter<ELEMENT_DIM,SPACE_DIM>::MakeVtkMesh()
00071 {
00072     //Construct nodes aka as Points
00073     vtkPoints* p_pts = vtkPoints::New(VTK_DOUBLE);
00074     p_pts->GetData()->SetName("Vertex positions");
00075     for (unsigned item_num=0; item_num<this->GetNumNodes(); item_num++)
00076     {
00077         std::vector<double> current_item = this->GetNextNode(); //this->mNodeData[item_num];
00078         // Add zeroes if the dimension is below 3
00079         for (unsigned dim=SPACE_DIM; dim<3; dim++)
00080         {
00081             current_item.push_back(0.0);//For y and z-coordinates if necessary
00082         }
00083         assert(current_item.size() == 3);
00084         p_pts->InsertPoint(item_num, current_item[0], current_item[1], current_item[2]);
00085     }
00086     mpVtkUnstructedMesh->SetPoints(p_pts);
00087     p_pts->Delete(); //Reference counted
00088 
00089     //Construct elements aka Cells
00090     for (unsigned item_num=0; item_num<this->GetNumElements(); item_num++)
00091     {
00092         std::vector<unsigned> current_element = this->GetNextElement().NodeIndices; // this->mElementData[item_num];
00093 
00094         assert((current_element.size() == ELEMENT_DIM + 1) || (current_element.size() == (ELEMENT_DIM+1)*(ELEMENT_DIM+2)/2));
00095 
00096         vtkCell* p_cell=NULL;
00097         if (ELEMENT_DIM == 3 && current_element.size() == 4)
00098         {
00099             p_cell = vtkTetra::New();
00100         }
00101         else if(ELEMENT_DIM == 3 && current_element.size() == 10)
00102         {
00103             p_cell = vtkQuadraticTetra::New();
00104         }
00105         else if (ELEMENT_DIM == 2 && current_element.size() == 3)
00106         {
00107             p_cell = vtkTriangle::New();
00108         }
00109         else if (ELEMENT_DIM == 2 && current_element.size() == 6)
00110         {
00111             p_cell = vtkQuadraticTriangle::New();
00112         }
00113         else if (ELEMENT_DIM == 1)
00114         {
00115             p_cell = vtkLine::New();
00116         }
00117 
00118         //Set the linear nodes
00119         vtkIdList* p_cell_id_list = p_cell->GetPointIds();
00120         for (unsigned j = 0; j < current_element.size(); ++j)
00121         {
00122             p_cell_id_list->SetId(j, current_element[j]);
00123         }
00124 
00125         //VTK defines the node ordering in quadratic triangles differently to Chaste, so they must be treated as a special case
00126         if (SPACE_DIM == 2 && current_element.size() == 6)
00127         {
00128             p_cell_id_list->SetId(3, current_element[5]);
00129             p_cell_id_list->SetId(4, current_element[3]);
00130             p_cell_id_list->SetId(5, current_element[4]);
00131         }
00132 
00133         mpVtkUnstructedMesh->InsertNextCell(p_cell->GetCellType(), p_cell_id_list);
00134         p_cell->Delete(); //Reference counted
00135     }
00136 
00137     if (SPACE_DIM > 1)
00138     {
00140         for (unsigned item_num=0; item_num<this->GetNumBoundaryFaces(); item_num++)
00141         {
00142             this->GetNextBoundaryElement();
00143         }
00144     }
00145 
00146     //If necessary, construct cables
00147     if (this->GetNumCableElements() > 0)
00148     {
00149         AugmentCellData();
00150         //Make a blank cell radius data for the regular elements
00151         std::vector<double> radii(this->GetNumElements(), 0.0);
00152         for (unsigned item_num=0; item_num<this->GetNumCableElements(); item_num++)
00153         {
00154             ElementData cable_element_data = this->GetNextCableElement();
00155             std::vector<unsigned> current_element = cable_element_data.NodeIndices;
00156             radii.push_back(cable_element_data.AttributeValue);
00157             assert(current_element.size() == 2);
00158             vtkCell* p_cell=vtkLine::New();
00159             vtkIdList* p_cell_id_list = p_cell->GetPointIds();
00160             for (unsigned j = 0; j < 2; ++j)
00161             {
00162                 p_cell_id_list->SetId(j, current_element[j]);
00163             }
00164             mpVtkUnstructedMesh->InsertNextCell(p_cell->GetCellType(), p_cell_id_list);
00165             p_cell->Delete(); //Reference counted
00166         }
00167         AddCellData("Cable radius", radii);
00168 
00169     }
00170 }
00171 
00172 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00173 void VtkMeshWriter<ELEMENT_DIM,SPACE_DIM>::AddProvenance(std::string fileName)
00174 {
00175     std::string comment = "<!-- " + ChasteBuildInfo::GetProvenanceString() + "-->";
00176 
00177     out_stream p_vtu_file = this->mpOutputFileHandler->OpenOutputFile(fileName, std::ios::out | std::ios::app);
00178 
00179     *p_vtu_file << "\n" << comment << "\n";
00180     p_vtu_file->close();
00181 }
00182 
00183 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00184 void VtkMeshWriter<ELEMENT_DIM,SPACE_DIM>::WriteFiles()
00185 {
00186     // Using separate scope here to make sure file is properly closed before re-opening it to add provenance info.
00187     {
00188         MakeVtkMesh();
00189         assert(mpVtkUnstructedMesh->CheckAttributes() == 0);
00190         vtkXMLUnstructuredGridWriter* p_writer = vtkXMLUnstructuredGridWriter::New();
00191 #if VTK_MAJOR_VERSION >= 6
00192         p_writer->SetInputData(mpVtkUnstructedMesh);
00193 #else
00194         p_writer->SetInput(mpVtkUnstructedMesh);
00195 #endif
00196         //Uninitialised stuff arises (see #1079), but you can remove
00197         //valgrind problems by removing compression:
00198         // **** REMOVE WITH CAUTION *****
00199         p_writer->SetCompressor(NULL);
00200         // **** REMOVE WITH CAUTION *****
00201         std::string vtk_file_name = this->mpOutputFileHandler->GetOutputDirectoryFullPath() + this->mBaseName+".vtu";
00202         p_writer->SetFileName(vtk_file_name.c_str());
00203         //p_writer->PrintSelf(std::cout, vtkIndent());
00204         p_writer->Write();
00205         p_writer->Delete(); //Reference counted
00206     }
00207 
00208     AddProvenance(this->mBaseName + ".vtu");
00209 }
00210 
00211 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00212 void VtkMeshWriter<ELEMENT_DIM,SPACE_DIM>::AddCellData(std::string dataName, std::vector<double> dataPayload)
00213 {
00214     vtkDoubleArray* p_scalars = vtkDoubleArray::New();
00215     p_scalars->SetName(dataName.c_str());
00216     for (unsigned i=0; i<dataPayload.size(); i++)
00217     {
00218         p_scalars->InsertNextValue(dataPayload[i]);
00219     }
00220 
00221     vtkCellData* p_cell_data = mpVtkUnstructedMesh->GetCellData();
00222     p_cell_data->AddArray(p_scalars);
00223     p_scalars->Delete(); //Reference counted
00224 }
00225 
00226 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00227 void VtkMeshWriter<ELEMENT_DIM,SPACE_DIM>::AugmentCellData()
00228 {
00229     unsigned num_cell_arrays = mpVtkUnstructedMesh->GetCellData()->GetNumberOfArrays();
00230     for (unsigned i = 0; i < num_cell_arrays; i++)
00231     {
00232         vtkDataArray* array = mpVtkUnstructedMesh->GetCellData()->GetArray(i);
00233 
00234         //Check data was the correct size before the cables were added
00235         unsigned num_cable_pads = this->GetNumCableElements();
00236         if (mWriteParallelFiles)
00237         {
00238             assert((unsigned)array->GetNumberOfTuples() == this->mpDistributedMesh->GetNumLocalElements());
00239             num_cable_pads =  this->mpMixedMesh->GetNumLocalCableElements();
00240         }
00241         else
00242         {
00243             assert((unsigned)array->GetNumberOfTuples() == this->GetNumElements());
00244         }
00245 
00246         //Check that tuples of size 3 will be big enough for padding the rest of the data
00247         assert(array->GetNumberOfComponents() <= 3);
00248         double null_data[3] = {0.0, 0.0, 0.0};
00249 
00250         //Pad data
00251         for (unsigned new_index = 0; new_index <  num_cable_pads; new_index++)
00252         {
00253             array->InsertNextTuple(null_data);
00254         }
00255     }
00256 }
00257 
00258 
00259 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00260 void VtkMeshWriter<ELEMENT_DIM,SPACE_DIM>::AddCellData(std::string dataName, std::vector<c_vector<double, SPACE_DIM> > dataPayload)
00261 {
00262     vtkDoubleArray* p_vectors = vtkDoubleArray::New();
00263     p_vectors->SetName(dataName.c_str());
00264     p_vectors->SetNumberOfComponents(3);
00265     for (unsigned i=0; i<dataPayload.size(); i++)
00266     {
00267         for (unsigned j=0; j<SPACE_DIM; j++)
00268         {
00269             p_vectors->InsertNextValue(dataPayload[i][j]);
00270         }
00271         //When SPACE_DIM<3, then pad
00272         for (unsigned j=SPACE_DIM; j<3; j++)
00273         {
00274             p_vectors->InsertNextValue(0.0);
00275         }
00276     }
00277 
00278     vtkCellData* p_cell_data = mpVtkUnstructedMesh->GetCellData();
00279     p_cell_data->AddArray(p_vectors);
00280     p_vectors->Delete(); //Reference counted
00281 }
00282 
00283 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00284 void VtkMeshWriter<ELEMENT_DIM,SPACE_DIM>::AddTensorCellData(std::string dataName, std::vector<c_vector<double,SPACE_DIM*(SPACE_DIM+1)/2> > dataPayload)
00285 {
00286     assert(SPACE_DIM != 1);
00287 
00288     vtkDoubleArray* p_vectors = vtkDoubleArray::New();
00289     p_vectors->SetName(dataName.c_str());
00290     p_vectors->SetNumberOfComponents(SPACE_DIM*SPACE_DIM);
00291     for (unsigned i=0; i<dataPayload.size(); i++)
00292     {
00293         if(SPACE_DIM == 2)
00294         {
00295             p_vectors->InsertNextValue(dataPayload[i](0)); //a11
00296             p_vectors->InsertNextValue(dataPayload[i](1)); //a12
00297             p_vectors->InsertNextValue(dataPayload[i](1)); //a21
00298             p_vectors->InsertNextValue(dataPayload[i](2)); //a22
00299         }
00300         else if (SPACE_DIM == 3)
00301         {
00302             p_vectors->InsertNextValue(dataPayload[i](0)); //a11
00303             p_vectors->InsertNextValue(dataPayload[i](1)); //a12
00304             p_vectors->InsertNextValue(dataPayload[i](2)); //a13
00305             p_vectors->InsertNextValue(dataPayload[i](1)); //a21
00306             p_vectors->InsertNextValue(dataPayload[i](3)); //a22
00307             p_vectors->InsertNextValue(dataPayload[i](4)); //a23
00308             p_vectors->InsertNextValue(dataPayload[i](2)); //a31
00309             p_vectors->InsertNextValue(dataPayload[i](4)); //a32
00310             p_vectors->InsertNextValue(dataPayload[i](5)); //a33
00311         }
00312     }
00313 
00314     vtkCellData* p_cell_data = mpVtkUnstructedMesh->GetCellData();
00315     p_cell_data->AddArray(p_vectors);
00316     p_vectors->Delete(); //Reference counted
00317 }
00318 
00319 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00320 void VtkMeshWriter<ELEMENT_DIM,SPACE_DIM>::AddTensorCellData(std::string dataName, std::vector<c_matrix<double,SPACE_DIM,SPACE_DIM> > dataPayload)
00321 {
00322     assert(SPACE_DIM != 1);
00323 
00324     vtkDoubleArray* p_vectors = vtkDoubleArray::New();
00325     p_vectors->SetName(dataName.c_str());
00326     p_vectors->SetNumberOfComponents(SPACE_DIM*SPACE_DIM);
00327     for (unsigned i=0; i<dataPayload.size(); i++)
00328     {
00329         if(SPACE_DIM == 2)
00330         {
00331             p_vectors->InsertNextValue(dataPayload[i](0,0)); //a11
00332             p_vectors->InsertNextValue(dataPayload[i](0,1)); //a12
00333             p_vectors->InsertNextValue(dataPayload[i](1,0)); //a21
00334             p_vectors->InsertNextValue(dataPayload[i](1,1)); //a22
00335         }
00336         else if (SPACE_DIM == 3)
00337         {
00338             p_vectors->InsertNextValue(dataPayload[i](0,0)); //a11
00339             p_vectors->InsertNextValue(dataPayload[i](0,1)); //a12
00340             p_vectors->InsertNextValue(dataPayload[i](0,2)); //a13
00341             p_vectors->InsertNextValue(dataPayload[i](1,0)); //a21
00342             p_vectors->InsertNextValue(dataPayload[i](1,1)); //a22
00343             p_vectors->InsertNextValue(dataPayload[i](1,2)); //a23
00344             p_vectors->InsertNextValue(dataPayload[i](2,0)); //a31
00345             p_vectors->InsertNextValue(dataPayload[i](2,1)); //a32
00346             p_vectors->InsertNextValue(dataPayload[i](2,2)); //a33
00347         }
00348     }
00349 
00350     vtkCellData* p_cell_data = mpVtkUnstructedMesh->GetCellData();
00351     p_cell_data->AddArray(p_vectors);
00352     p_vectors->Delete(); //Reference counted
00353 }
00354 
00355 
00356 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00357 void VtkMeshWriter<ELEMENT_DIM,SPACE_DIM>::AddPointData(std::string dataName, std::vector<double> dataPayload)
00358 {
00359     vtkDoubleArray* p_scalars = vtkDoubleArray::New();
00360     p_scalars->SetName(dataName.c_str());
00361 
00362     if (mWriteParallelFiles && this->mpDistributedMesh != NULL)
00363     {
00364         // In parallel, the vector we pass will only contain the values from the privately owned nodes.
00365         // To get the values from the halo nodes (which will be inserted at the end of the vector we need to
00366         // communicate with the equivalent vectors on other processes.
00367 
00368         // resize the payload data to include halos
00369         assert( dataPayload.size() == this->mpDistributedMesh->GetNumLocalNodes() );
00370         dataPayload.resize( this->mpDistributedMesh->GetNumLocalNodes() + this->mpDistributedMesh->GetNumHaloNodes() );
00371 
00372 
00373         // then do the communication
00374         for ( unsigned rank_offset = 1; rank_offset < PetscTools::GetNumProcs(); rank_offset++ )
00375         {
00376             unsigned send_to      = (PetscTools::GetMyRank() + rank_offset) % (PetscTools::GetNumProcs());
00377             unsigned receive_from = (PetscTools::GetMyRank() + PetscTools::GetNumProcs()- rank_offset ) % (PetscTools::GetNumProcs());
00378 
00379             unsigned number_of_nodes_to_send    = mNodesToSendPerProcess[send_to].size();
00380             unsigned number_of_nodes_to_receive = mNodesToReceivePerProcess[receive_from].size();
00381 
00382             boost::scoped_array<double> send_data(new double[number_of_nodes_to_send]);
00383             boost::scoped_array<double> receive_data(new double[number_of_nodes_to_receive]);
00384             // Pack
00385             for (unsigned node = 0; node < number_of_nodes_to_send; node++)
00386             {
00387                 unsigned global_node_index = mNodesToSendPerProcess[send_to][node];
00388                 unsigned local_node_index = global_node_index
00389                             - this->mpDistributedMesh->GetDistributedVectorFactory()->GetLow();
00390                 send_data[node] = dataPayload[local_node_index];
00391             }
00392             {
00393                 // Send
00394                 int ret;
00395                 MPI_Status status;
00396                 ret = MPI_Sendrecv(send_data.get(), number_of_nodes_to_send,
00397                                    MPI_DOUBLE,
00398                                    send_to, 0,
00399                                    receive_data.get(),  number_of_nodes_to_receive,
00400                                    MPI_DOUBLE,
00401                                    receive_from, 0,
00402                                    PETSC_COMM_WORLD, &status);
00403                 UNUSED_OPT(ret);
00404                 assert ( ret == MPI_SUCCESS );
00405             }
00406 
00407             // Unpack
00408             for ( unsigned node = 0; node < number_of_nodes_to_receive; node++ )
00409             {
00410                 unsigned global_node_index = mNodesToReceivePerProcess[receive_from][node];
00411                 unsigned halo_index = mGlobalToNodeIndexMap[global_node_index];
00412                 assert( halo_index >= this->mpDistributedMesh->GetNumLocalNodes() );
00413                 dataPayload[halo_index] = receive_data[node];
00414             }
00415 
00416         }
00417     }
00418 
00419     for (unsigned i=0; i<dataPayload.size(); i++)
00420     {
00421         p_scalars->InsertNextValue(dataPayload[i]);
00422     }
00423 
00424     vtkPointData* p_point_data = mpVtkUnstructedMesh->GetPointData();
00425     p_point_data->AddArray(p_scalars);
00426     p_scalars->Delete(); //Reference counted
00427 }
00428 
00429 
00430 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00431 void VtkMeshWriter<ELEMENT_DIM,SPACE_DIM>::AddPointData(std::string dataName, std::vector<c_vector<double, SPACE_DIM> > dataPayload)
00432 {
00433     vtkDoubleArray* p_vectors = vtkDoubleArray::New();
00434     p_vectors->SetName(dataName.c_str());
00435 
00436     if (mWriteParallelFiles)
00437     {
00438         // In parallel, the vector we pass will only contain the values from the privately owned nodes.
00439         // To get the values from the halo nodes (which will be inserted at the end of the vector we need to
00440         // communicate with the equivalent vectors on other processes.
00441 
00442         // resize the payload data to include halos
00443         assert( dataPayload.size() == this->mpDistributedMesh->GetNumLocalNodes() );
00444         dataPayload.resize( this->mpDistributedMesh->GetNumLocalNodes() + this->mpDistributedMesh->GetNumHaloNodes() );
00445 
00446         // then do the communication
00447         for ( unsigned rank_offset = 1; rank_offset < PetscTools::GetNumProcs(); rank_offset++ )
00448         {
00449             unsigned send_to      = (PetscTools::GetMyRank() + rank_offset) % (PetscTools::GetNumProcs());
00450             unsigned receive_from = (PetscTools::GetMyRank() + PetscTools::GetNumProcs()- rank_offset ) % (PetscTools::GetNumProcs());
00451 
00452             unsigned number_of_nodes_to_send    = mNodesToSendPerProcess[send_to].size();
00453             unsigned number_of_nodes_to_receive = mNodesToReceivePerProcess[receive_from].size();
00454 
00455             boost::scoped_array<double> send_data(new double[number_of_nodes_to_send * SPACE_DIM]);
00456             boost::scoped_array<double> receive_data(new double[number_of_nodes_to_receive * SPACE_DIM]);
00457 
00458             for (unsigned node = 0; node < number_of_nodes_to_send; node++)
00459             {
00460                 unsigned global_node_index = mNodesToSendPerProcess[send_to][node];
00461                 unsigned local_node_index = global_node_index
00462                             - this->mpDistributedMesh->GetDistributedVectorFactory()->GetLow();
00463                 for (unsigned j=0; j<SPACE_DIM; j++)
00464                 {
00465                     send_data[ node*SPACE_DIM + j ] = dataPayload[local_node_index][j];
00466                 }
00467             }
00468 
00469                 int ret;
00470                 MPI_Status status;
00471                 ret = MPI_Sendrecv(send_data.get(), number_of_nodes_to_send * SPACE_DIM,
00472                                    MPI_DOUBLE,
00473                                    send_to, 0,
00474                                    receive_data.get(),  number_of_nodes_to_receive * SPACE_DIM,
00475                                    MPI_DOUBLE,
00476                                    receive_from, 0,
00477                                    PETSC_COMM_WORLD, &status);
00478                 UNUSED_OPT(ret);
00479                 assert ( ret == MPI_SUCCESS );
00480 
00481             // Unpack
00482             for ( unsigned node = 0; node < number_of_nodes_to_receive; node++ )
00483             {
00484                 unsigned global_node_index = mNodesToReceivePerProcess[receive_from][node];
00485                 unsigned halo_index = mGlobalToNodeIndexMap[global_node_index];
00486                 assert( halo_index >= this->mpDistributedMesh->GetNumLocalNodes() );
00487                 for (unsigned j=0; j<SPACE_DIM; j++)
00488                 {
00489                     dataPayload[halo_index][j] = receive_data[ node*SPACE_DIM + j ];
00490                 }
00491             }
00492         }
00493     }
00494 
00495     p_vectors->SetNumberOfComponents(3);
00496     for (unsigned i=0; i<dataPayload.size(); i++)
00497     {
00498         for (unsigned j=0; j<SPACE_DIM; j++)
00499         {
00500             p_vectors->InsertNextValue(dataPayload[i][j]);
00501         }
00502         //When SPACE_DIM<3, then pad
00503         for (unsigned j=SPACE_DIM; j<3; j++)
00504         {
00505             p_vectors->InsertNextValue(0.0);
00506         }
00507     }
00508 
00509     vtkPointData* p_point_data = mpVtkUnstructedMesh->GetPointData();
00510     p_point_data->AddArray(p_vectors);
00511     p_vectors->Delete(); //Reference counted
00512 }
00513 
00514 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00515 void VtkMeshWriter<ELEMENT_DIM,SPACE_DIM>::AddTensorPointData(std::string dataName, std::vector<c_matrix<double,SPACE_DIM,SPACE_DIM> > dataPayload)
00516 {
00517     assert(SPACE_DIM != 1);
00518 
00519     vtkDoubleArray* p_vectors = vtkDoubleArray::New();
00520     p_vectors->SetName(dataName.c_str());
00521     p_vectors->SetNumberOfComponents(SPACE_DIM*SPACE_DIM);
00522     for (unsigned i=0; i<dataPayload.size(); i++)
00523     {
00524         if(SPACE_DIM == 2)
00525         {
00526             p_vectors->InsertNextValue(dataPayload[i](0,0)); //a11
00527             p_vectors->InsertNextValue(dataPayload[i](0,1)); //a12
00528             p_vectors->InsertNextValue(dataPayload[i](1,0)); //a21
00529             p_vectors->InsertNextValue(dataPayload[i](1,1)); //a22
00530         }
00531         else if (SPACE_DIM == 3)
00532         {
00533             p_vectors->InsertNextValue(dataPayload[i](0,0)); //a11
00534             p_vectors->InsertNextValue(dataPayload[i](0,1)); //a12
00535             p_vectors->InsertNextValue(dataPayload[i](0,2)); //a13
00536             p_vectors->InsertNextValue(dataPayload[i](1,0)); //a21
00537             p_vectors->InsertNextValue(dataPayload[i](1,1)); //a22
00538             p_vectors->InsertNextValue(dataPayload[i](1,2)); //a23
00539             p_vectors->InsertNextValue(dataPayload[i](2,0)); //a31
00540             p_vectors->InsertNextValue(dataPayload[i](2,1)); //a32
00541             p_vectors->InsertNextValue(dataPayload[i](2,2)); //a33
00542         }
00543     }
00544 
00545     vtkPointData* p_point_data = mpVtkUnstructedMesh->GetPointData();
00546     p_point_data->AddArray(p_vectors);
00547     p_vectors->Delete(); //Reference counted
00548 }
00549 
00550 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00551 void VtkMeshWriter<ELEMENT_DIM,SPACE_DIM>::SetParallelFiles( AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>& rMesh )
00552 {
00553     //Have we got a distributed mesh?
00554     this->mpDistributedMesh = dynamic_cast<DistributedTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* >(&rMesh);
00555     mpNodesOnlyMesh = dynamic_cast<NodesOnlyMesh<SPACE_DIM>* >(&rMesh);
00556 
00557     if (this->mpDistributedMesh == NULL && mpNodesOnlyMesh == NULL)
00558     {
00559         EXCEPTION("Cannot write parallel files using a sequential mesh");
00560     }
00561 
00562     if (PetscTools::IsSequential())
00563     {
00564         return;     // mWriteParallelFiles is not set sequentially (so we don't set up data exchange machinery)
00565     }
00566 
00567     mWriteParallelFiles = true;
00568 
00569     // Populate the global to node index map (as this will be required to add point data)
00570 
00571     //Node index that we are writing to VTK (index into mNodes and mHaloNodes as if they were concatenated)
00572     unsigned index = 0;
00573 
00574     // Owned nodes
00575     for (typename AbstractMesh<ELEMENT_DIM,SPACE_DIM>::NodeIterator node_iter = rMesh.GetNodeIteratorBegin();
00576          node_iter != rMesh.GetNodeIteratorEnd();
00577          ++node_iter)
00578     {
00579         mGlobalToNodeIndexMap[node_iter->GetIndex()] = index;
00580         index++;
00581     }
00582 
00583     // Halo nodes
00584     if (this->mpDistributedMesh)
00585     {
00586         for (typename DistributedTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::HaloNodeIterator halo_iter=this->mpDistributedMesh->GetHaloNodeIteratorBegin();
00587                 halo_iter != this->mpDistributedMesh->GetHaloNodeIteratorEnd();
00588                 ++halo_iter)
00589         {
00590             mGlobalToNodeIndexMap[(*halo_iter)->GetIndex()] = index;
00591             index++;
00592         }
00593 
00594         //Calculate the halo exchange so that node-wise payloads can be communicated
00595         this->mpDistributedMesh->CalculateNodeExchange( mNodesToSendPerProcess, mNodesToReceivePerProcess );
00596     }
00597 }
00598 
00600 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00601 void VtkMeshWriter<ELEMENT_DIM, SPACE_DIM>::WriteFilesUsingMesh(
00602       AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>& rMesh,
00603       bool keepOriginalElementIndexing)
00604 {
00605     //Have we got a parallel mesh?
00606     this->mpDistributedMesh = dynamic_cast<DistributedTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* >(&rMesh);
00607     this->mpMixedMesh = dynamic_cast<MixedDimensionMesh<ELEMENT_DIM,SPACE_DIM>* >(&rMesh);
00608 
00609     if ( PetscTools::IsSequential() || !mWriteParallelFiles || (this->mpDistributedMesh == NULL && mpNodesOnlyMesh == NULL) )
00610     {
00611         AbstractTetrahedralMeshWriter<ELEMENT_DIM,SPACE_DIM>::WriteFilesUsingMesh( rMesh,keepOriginalElementIndexing );
00612     }
00613     else
00614     {
00615         //Make the local mesh into a VtkMesh
00616         vtkPoints* p_pts = vtkPoints::New(VTK_DOUBLE);
00617         p_pts->GetData()->SetName("Vertex positions");
00618 
00619         // Owned nodes
00620         for (typename AbstractMesh<ELEMENT_DIM,SPACE_DIM>::NodeIterator node_iter = rMesh.GetNodeIteratorBegin();
00621              node_iter != rMesh.GetNodeIteratorEnd();
00622              ++node_iter)
00623         {
00624             c_vector<double, SPACE_DIM> current_item = node_iter->rGetLocation();
00625             if (SPACE_DIM == 3)
00626             {
00627                 p_pts->InsertNextPoint(current_item[0], current_item[1], current_item[2]);
00628             }
00629             else if (SPACE_DIM == 2)
00630             {
00631                 p_pts->InsertNextPoint(current_item[0], current_item[1], 0.0);
00632             }
00633             else // (SPACE_DIM == 1)
00634             {
00635                 p_pts->InsertNextPoint(current_item[0], 0.0, 0.0);
00636             }
00637         }
00638 
00639         // Halo nodes
00640         if (this->mpDistributedMesh)
00641         {
00642             for (typename DistributedTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::HaloNodeIterator halo_iter=this->mpDistributedMesh->GetHaloNodeIteratorBegin();
00643                     halo_iter != this->mpDistributedMesh->GetHaloNodeIteratorEnd();
00644                     ++halo_iter)
00645             {
00646                 c_vector<double, SPACE_DIM> current_item = (*halo_iter)->rGetLocation();
00647                 if (SPACE_DIM == 3)
00648                 {
00649                     p_pts->InsertNextPoint(current_item[0], current_item[1], current_item[2]);
00650                 }
00651                 else if (SPACE_DIM == 2)
00652                 {
00653                     p_pts->InsertNextPoint(current_item[0], current_item[1], 0.0);
00654                 }
00655                 else // (SPACE_DIM == 1)
00656                 {
00657                     p_pts->InsertNextPoint(current_item[0], 0.0, 0.0);
00658                 }
00659             }
00660         }
00661 
00662         mpVtkUnstructedMesh->SetPoints(p_pts);
00663         p_pts->Delete(); //Reference counted
00664 
00665         for (typename AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>::ElementIterator elem_iter = rMesh.GetElementIteratorBegin();
00666              elem_iter != rMesh.GetElementIteratorEnd();
00667              ++elem_iter)
00668         {
00669 
00670             vtkCell* p_cell=NULL;
00672             if (ELEMENT_DIM == 3)
00673             {
00674                 p_cell = vtkTetra::New();
00675             }
00676             else if (ELEMENT_DIM == 2)
00677             {
00678                 p_cell = vtkTriangle::New();
00679             }
00680             else //(ELEMENT_DIM == 1)
00681             {
00682                 p_cell = vtkLine::New();
00683             }
00684             vtkIdList* p_cell_id_list = p_cell->GetPointIds();
00685             for (unsigned j = 0; j < ELEMENT_DIM+1; ++j)
00686             {
00687                 unsigned global_node_index = elem_iter->GetNodeGlobalIndex(j);
00688                 p_cell_id_list->SetId(j, mGlobalToNodeIndexMap[global_node_index]);
00689             }
00690             mpVtkUnstructedMesh->InsertNextCell(p_cell->GetCellType(), p_cell_id_list);
00691             p_cell->Delete(); //Reference counted
00692         }
00693         //If necessary, construct cables
00694         if (this->mpMixedMesh )
00695         {
00696             AugmentCellData();
00697             //Make a blank cell radius data for the regular elements
00698             std::vector<double> radii(this->mpMixedMesh->GetNumLocalElements(), 0.0);
00699             for (typename MixedDimensionMesh<ELEMENT_DIM,SPACE_DIM>::CableElementIterator elem_iter = this->mpMixedMesh->GetCableElementIteratorBegin();
00700                  elem_iter != this->mpMixedMesh->GetCableElementIteratorEnd();
00701                  ++elem_iter)
00702             {
00703                 radii.push_back((*elem_iter)->GetAttribute());
00704                 vtkCell* p_cell=vtkLine::New();
00705                 vtkIdList* p_cell_id_list = p_cell->GetPointIds();
00706                 for (unsigned j = 0; j < 2; ++j)
00707                 {
00708                     unsigned global_node_index = (*elem_iter)->GetNodeGlobalIndex(j);
00709                     p_cell_id_list->SetId(j, mGlobalToNodeIndexMap[global_node_index]);
00710                 }
00711                 mpVtkUnstructedMesh->InsertNextCell(p_cell->GetCellType(), p_cell_id_list);
00712                 p_cell->Delete(); //Reference counted
00713             }
00714             AddCellData("Cable radius", radii);
00715         }
00716 
00717 
00718         //This block is to guard the mesh writers (vtkXMLPUnstructuredGridWriter) so that they
00719         //go out of scope, flush buffers and close files
00720         {
00721             assert(mpVtkUnstructedMesh->CheckAttributes() == 0);
00722             vtkXMLPUnstructuredGridWriter* p_writer = vtkXMLPUnstructuredGridWriter::New();
00723 
00724             p_writer->SetDataModeToBinary();
00725 
00726             p_writer->SetNumberOfPieces(PetscTools::GetNumProcs());
00727             //p_writer->SetGhostLevel(-1);
00728             p_writer->SetStartPiece(PetscTools::GetMyRank());
00729             p_writer->SetEndPiece(PetscTools::GetMyRank());
00730 
00731 
00732 #if VTK_MAJOR_VERSION >= 6
00733             p_writer->SetInputData(mpVtkUnstructedMesh);
00734 #else
00735             p_writer->SetInput(mpVtkUnstructedMesh);
00736 #endif
00737             //Uninitialised stuff arises (see #1079), but you can remove
00738             //valgrind problems by removing compression:
00739             // **** REMOVE WITH CAUTION *****
00740             p_writer->SetCompressor(NULL);
00741             // **** REMOVE WITH CAUTION *****
00742             std::string pvtk_file_name = this->mpOutputFileHandler->GetOutputDirectoryFullPath() + this->mBaseName+ ".pvtu";
00743             p_writer->SetFileName(pvtk_file_name.c_str());
00744             //p_writer->PrintSelf(std::cout, vtkIndent());
00745             p_writer->Write();
00746             p_writer->Delete(); //Reference counted
00747         }
00748 
00749         // Add provenance to the individual files
00750         std::stringstream filepath;
00751         filepath << this->mBaseName << "_" << PetscTools::GetMyRank() << ".vtu";
00752         AddProvenance(filepath.str());
00754         if (PetscTools::AmMaster())
00755         {
00756             AddProvenance(this->mBaseName+ ".pvtu");
00757         }
00758     }
00759 }
00760 
00762 // Explicit instantiation
00764 
00765 template class VtkMeshWriter<1,1>;
00766 template class VtkMeshWriter<1,2>;
00767 template class VtkMeshWriter<1,3>;
00768 template class VtkMeshWriter<2,2>; // Actually used
00769 template class VtkMeshWriter<2,3>;
00770 template class VtkMeshWriter<3,3>; // Actually used
00771 
00772 #endif //CHASTE_VTK

Generated by  doxygen 1.6.2