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 #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
00068 EXCEPTION("Size of vector does not match FixedDimensionSize.");
00069 }
00070
00071
00072 if (mConcentrated==NULL)
00073 {
00074 VecScatterCreateToZero(petscVector, &mToMaster, &mConcentrated);
00075 }
00076
00077
00078
00079
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
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
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
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
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
00166 MPI_Barrier(PETSC_COMM_WORLD);
00167
00168
00169
00170
00171 if (mAmMaster)
00172 {
00176 ColumnDataWriter::DoAdvanceAlongUnlimitedDimension();
00177 }
00178 }
00179
00180 void ParallelColumnDataWriter::Close()
00181 {
00182
00183 MPI_Barrier(PETSC_COMM_WORLD);
00184
00186 if (mAmMaster)
00187 {
00188 ColumnDataWriter::Close();
00189 }
00190 }
00191
00192