ParallelColumnDataWriter.cpp
00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029 #include "ParallelColumnDataWriter.hpp"
00030 #include "Exception.hpp"
00031 #include "DistributedVectorFactory.hpp"
00032
00033 ParallelColumnDataWriter::ParallelColumnDataWriter(const std::string& rDirectory,
00034 const std::string& rBaseName,
00035 bool cleanDirectory)
00036 : ColumnDataWriter::ColumnDataWriter(rDirectory, rBaseName, cleanDirectory),
00037 mConcentrated(NULL)
00038 {
00039 int num_procs;
00040 MPI_Comm_size(PETSC_COMM_WORLD, &num_procs);
00041 if (num_procs==1)
00042 {
00043 mIsParallel = false;
00044 }
00045 else
00046 {
00047 mIsParallel = true;
00048 }
00049 }
00050
00051 void ParallelColumnDataWriter::PutVector(int variableID, Vec petscVector)
00052 {
00053 int size;
00054 VecGetSize(petscVector,&size);
00055
00056 if (size != mFixedDimensionSize)
00057 {
00058 EXCEPTION("Size of vector does not match FixedDimensionSize.");
00059 }
00060
00061
00062 if (mConcentrated==NULL)
00063 {
00064 VecScatterCreateToZero(petscVector, &mToMaster, &mConcentrated);
00065 }
00066
00067
00068
00069
00070
00071
00072 #if ( (PETSC_VERSION_MAJOR == 3) || (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 3 && PETSC_VERSION_SUBMINOR == 3)) //2.3.3 or 3.x.x
00073 VecScatterBegin(mToMaster, petscVector, mConcentrated, INSERT_VALUES, SCATTER_FORWARD);
00074 VecScatterEnd(mToMaster, petscVector, mConcentrated, INSERT_VALUES, SCATTER_FORWARD);
00075 #else
00076 VecScatterBegin(petscVector, mConcentrated, INSERT_VALUES, SCATTER_FORWARD, mToMaster);
00077 VecScatterEnd(petscVector, mConcentrated, INSERT_VALUES, SCATTER_FORWARD, mToMaster);
00078 #endif
00079
00080
00081
00082 if (PetscTools::AmMaster())
00083 {
00084 double *concentrated_vector;
00085 VecGetArray(mConcentrated, &concentrated_vector);
00086 for (int i=0; i<size; i++)
00087 {
00088 ColumnDataWriter::PutVariable(variableID, concentrated_vector[i], i);
00089 }
00090 VecRestoreArray(mConcentrated, &concentrated_vector);
00091 }
00092 }
00093
00094 void ParallelColumnDataWriter::PutVectorStripe(int variableId, DistributedVector::Stripe& rStripe)
00095 {
00096
00097 DistributedVectorFactory* p_factory = rStripe.GetFactory();
00098 Vec unstriped_petsc = p_factory->CreateVec();
00099 DistributedVector unstriped = p_factory->CreateDistributedVector(unstriped_petsc);
00100 for (DistributedVector::Iterator index = unstriped.Begin();
00101 index!= unstriped.End();
00102 ++index)
00103 {
00104 unstriped[index] = rStripe[index];
00105 }
00106
00107
00108 ParallelColumnDataWriter::PutVector(variableId, unstriped_petsc);
00109 VecDestroy(unstriped_petsc);
00110 }
00111
00112 void ParallelColumnDataWriter::EndDefineMode()
00113 {
00114 if (PetscTools::AmMaster())
00115 {
00116 ColumnDataWriter::EndDefineMode();
00117 }
00118 else
00119 {
00120 mIsInDefineMode = false;
00121 }
00122 }
00123
00132 void ParallelColumnDataWriter::PutVariable(int variableID, double variableValue, long dimensionPosition)
00133 {
00134 if (PetscTools::AmMaster())
00135 {
00136
00137 ColumnDataWriter::PutVariable(variableID, variableValue, dimensionPosition);
00138 }
00139 }
00140
00141 ParallelColumnDataWriter::~ParallelColumnDataWriter()
00142 {
00143 if (mConcentrated != NULL)
00144 {
00145 VecScatterDestroy(mToMaster);
00146 VecDestroy(mConcentrated);
00147 }
00148 Close();
00149 }
00150
00151 void ParallelColumnDataWriter::AdvanceAlongUnlimitedDimension()
00152 {
00153
00154 PetscTools::Barrier("ParallelColumnDataWriter::AdvanceAlongUnlimitedDimension");
00155
00156 if (PetscTools::AmMaster())
00157 {
00158 ColumnDataWriter::DoAdvanceAlongUnlimitedDimension();
00159 }
00160 }
00161
00162 void ParallelColumnDataWriter::Close()
00163 {
00164
00165 PetscTools::Barrier("ParallelColumnDataWriter::Close");
00166
00167 if (PetscTools::AmMaster())
00168 {
00169 ColumnDataWriter::Close();
00170 }
00171 }