ParallelColumnDataWriter.cpp

00001 /*
00002 
00003 Copyright (C) University of Oxford, 2005-2009
00004 
00005 University of Oxford means the Chancellor, Masters and Scholars of the
00006 University of Oxford, having an administrative office at Wellington
00007 Square, Oxford OX1 2JD, UK.
00008 
00009 This file is part of Chaste.
00010 
00011 Chaste is free software: you can redistribute it and/or modify it
00012 under the terms of the GNU Lesser General Public License as published
00013 by the Free Software Foundation, either version 2.1 of the License, or
00014 (at your option) any later version.
00015 
00016 Chaste is distributed in the hope that it will be useful, but WITHOUT
00017 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018 FITNESS FOR A PARTICULAR PURPOSE.  See the GNU Lesser General Public
00019 License for more details. The offer of Chaste under the terms of the
00020 License is subject to the License being interpreted in accordance with
00021 English Law and subject to any action against the University of Oxford
00022 being under the jurisdiction of the English Courts.
00023 
00024 You should have received a copy of the GNU Lesser General Public License
00025 along with Chaste. If not, see <http://www.gnu.org/licenses/>.
00026 
00027 */
00028 #include "ParallelColumnDataWriter.hpp"
00029 #include "Exception.hpp"
00030 #include <iostream>
00031 
00032 ParallelColumnDataWriter::ParallelColumnDataWriter(std::string directory, std::string baseName, bool cleanDirectory)
00033         : ColumnDataWriter::ColumnDataWriter(directory, baseName, cleanDirectory)
00034 {
00035     mConcentrated=NULL;
00036 
00037     int num_procs, my_rank;
00038     MPI_Comm_size(PETSC_COMM_WORLD, &num_procs);
00039     if (num_procs==1)
00040     {
00041         mIsParallel=false;
00042     }
00043     else
00044     {
00045         mIsParallel=true;
00046     }
00047 
00048     MPI_Comm_rank(PETSC_COMM_WORLD, &my_rank);
00049     if (my_rank==0)
00050     {
00051         mAmMaster=true;
00052     }
00053     else
00054     {
00055         mAmMaster=false;
00056     }
00057 }
00058 
00059 void
00060 ParallelColumnDataWriter::PutVector(int variableID, Vec petscVector)
00061 {
00062     int size;
00063     VecGetSize(petscVector,&size);
00064 
00065     if (size != mFixedDimensionSize)
00066     {
00067         //std::cout << "fixed_dim: " << mFixedDimensionSize << ", vec_size " << size << "\n";
00068         EXCEPTION("Size of vector does not match FixedDimensionSize.");
00069     }
00070 
00071     //Construct the appropriate "scatter" object to concentrate the vector on the master
00072     if (mConcentrated==NULL)
00073     {
00074         VecScatterCreateToZero(petscVector, &mToMaster, &mConcentrated);
00075     }
00076 
00077 //    int size2;
00078 //    VecGetSize(mConcentrated, &size2);
00079 //    std::cout << "Vector size=" << size << "," << size2 << std::endl << std::flush;
00080 
00081 #if (PETSC_VERSION_MINOR == 3 && PETSC_VERSION_SUBMINOR == 3)
00082     VecScatterBegin(mToMaster, petscVector, mConcentrated, INSERT_VALUES, SCATTER_FORWARD);
00083     VecScatterEnd(mToMaster, petscVector, mConcentrated, INSERT_VALUES, SCATTER_FORWARD);
00084 #else
00085     VecScatterBegin(petscVector, mConcentrated, INSERT_VALUES, SCATTER_FORWARD, mToMaster);
00086     VecScatterEnd(petscVector, mConcentrated, INSERT_VALUES, SCATTER_FORWARD, mToMaster);
00087 #endif
00088 
00089 //    std::cout << "Done scatter" << std::endl << std::flush;
00090 
00091     if (mAmMaster)
00092     {
00093         double *concentrated_vector;
00094         VecGetArray(mConcentrated, &concentrated_vector);
00095         for (int i=0; i<size; i++)
00096         {
00097             ColumnDataWriter::PutVariable(variableID, concentrated_vector[i], i);
00098         }
00099         VecRestoreArray(mConcentrated, &concentrated_vector);
00100     }
00101 
00102 }
00103 
00104 void ParallelColumnDataWriter::PutVectorStripe(int variableId, DistributedVector::Stripe stripe)
00105 {
00106     // put the stripe into its own 'unstriped' vector
00107     Vec unstriped_petsc = DistributedVector::CreateVec();
00108     DistributedVector unstriped(unstriped_petsc);
00109     for (DistributedVector::Iterator index = DistributedVector::Begin();
00110          index!= DistributedVector::End();
00111          ++index)
00112         {
00113             unstriped[index] =  stripe[index];
00114         }
00115 
00116     // put the unstriped vector
00117     ParallelColumnDataWriter::PutVector(variableId, unstriped_petsc);
00118     VecDestroy(unstriped_petsc);
00119 }
00120 
00121 void ParallelColumnDataWriter::EndDefineMode()
00122 {
00123     if (mAmMaster)
00124     {
00125         ColumnDataWriter::EndDefineMode();
00126     }
00127     else
00128     {
00129         mIsInDefineMode = false;
00130     }
00131 }
00132 
00141 void ParallelColumnDataWriter::PutVariable(int variableID, double variableValue,long dimensionPosition)
00142 {
00143     if (mAmMaster)
00144     {
00145         //Master process is allowed to write
00146         ColumnDataWriter::PutVariable(variableID, variableValue, dimensionPosition);
00147     }
00148 }
00149 
00150 
00151 ParallelColumnDataWriter::~ParallelColumnDataWriter()
00152 {
00153     if (mConcentrated != NULL)
00154     {
00155         VecScatterDestroy(mToMaster);
00156         VecDestroy(mConcentrated);
00157     }
00158     Close();
00159 }
00160 
00161 
00162 
00163 void ParallelColumnDataWriter::AdvanceAlongUnlimitedDimension()
00164 {
00165     //Make sure that everyone has queued their messages
00166     MPI_Barrier(PETSC_COMM_WORLD);
00167 
00168 //    std::cout<<"In AdvanceAlongUnlimitedDimension mAmMaster="<< mAmMaster<<
00169 //     " mpCurrentOutputFile="<<mpCurrentOutputFile<<"\n"<<std::flush;
00170 
00171     if (mAmMaster)
00172     {
00176         ColumnDataWriter::DoAdvanceAlongUnlimitedDimension();
00177     }
00178 }
00179 
00180 void ParallelColumnDataWriter::Close()
00181 {
00182     //std::cout<<"In Close mMyRank="<< mMyRank<<"\n";
00183     MPI_Barrier(PETSC_COMM_WORLD);
00184 
00186     if (mAmMaster)
00187     {
00188         ColumnDataWriter::Close();
00189     }
00190 }
00191 
00192 

Generated on Wed Mar 18 12:51:52 2009 for Chaste by  doxygen 1.5.5