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