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 #include "ReplicatableVector.hpp"
00031
00032 #include <cassert>
00033
00034
00035
00036 void ReplicatableVector::RemovePetscContext()
00037 {
00038 if (mToAll != NULL)
00039 {
00040 VecScatterDestroy(mToAll);
00041 mToAll = NULL;
00042 }
00043
00044 if (mReplicated != NULL)
00045 {
00046 VecDestroy(mReplicated);
00047 mReplicated = NULL;
00048 }
00049
00050 if (mDistributed != NULL)
00051 {
00052 VecDestroy(mDistributed);
00053 mDistributed = NULL;
00054 }
00055 }
00056
00057
00058
00059 ReplicatableVector::ReplicatableVector()
00060 : mToAll(NULL),
00061 mReplicated(NULL),
00062 mDistributed(NULL)
00063 {
00064 }
00065
00066 ReplicatableVector::ReplicatableVector(Vec vec)
00067 : mToAll(NULL),
00068 mReplicated(NULL),
00069 mDistributed(NULL)
00070 {
00071 ReplicatePetscVector(vec);
00072 }
00073
00074 ReplicatableVector::ReplicatableVector(unsigned size)
00075 : mToAll(NULL),
00076 mReplicated(NULL),
00077 mDistributed(NULL)
00078 {
00079 resize(size);
00080 }
00081
00082 ReplicatableVector::~ReplicatableVector()
00083 {
00084 RemovePetscContext();
00085 }
00086
00087
00088
00089
00090 unsigned ReplicatableVector::size()
00091 {
00092 return mData.size();
00093 }
00094
00095 void ReplicatableVector::resize(unsigned size)
00096 {
00097
00098 RemovePetscContext();
00099 mData.resize(size);
00100 }
00101
00102 double& ReplicatableVector::operator[](unsigned index)
00103 {
00104 assert(index < mData.size());
00105 return mData[index];
00106 }
00107
00108
00109
00110
00111 void ReplicatableVector::Replicate(unsigned lo, unsigned hi)
00112 {
00113
00114 if (mDistributed==NULL)
00115 {
00116 VecCreateMPI(PETSC_COMM_WORLD, hi-lo, this->size(), &mDistributed);
00117 }
00118
00119 double *p_distributed;
00120 VecGetArray(mDistributed, &p_distributed);
00121 for (unsigned global_index=lo; global_index<hi; global_index++)
00122 {
00123 p_distributed[ (global_index-lo) ] = mData[global_index];
00124 }
00125 VecAssemblyBegin(mDistributed);
00126 VecAssemblyEnd(mDistributed);
00127
00128
00129 ReplicatePetscVector(mDistributed);
00130 }
00131
00132 void ReplicatableVector::ReplicatePetscVector(Vec vec)
00133 {
00134
00135 PetscInt isize;
00136 VecGetSize(vec, &isize);
00137 unsigned size=isize;
00138
00139 if (this->size() != size)
00140 {
00141 resize(size);
00142 }
00143 if (mReplicated == NULL)
00144 {
00145
00146 VecScatterCreateToAll(vec, &mToAll, &mReplicated);
00147 }
00148
00149
00150 #if (PETSC_VERSION_MINOR == 3 && PETSC_VERSION_SUBMINOR == 3)
00151 VecScatterBegin(mToAll, vec, mReplicated, INSERT_VALUES, SCATTER_FORWARD);
00152 VecScatterEnd (mToAll, vec, mReplicated, INSERT_VALUES, SCATTER_FORWARD);
00153 #else
00154 VecScatterBegin(vec, mReplicated, INSERT_VALUES, SCATTER_FORWARD, mToAll);
00155 VecScatterEnd (vec, mReplicated, INSERT_VALUES, SCATTER_FORWARD, mToAll);
00156 #endif
00157
00158
00159
00160 double *p_replicated;
00161 VecGetArray(mReplicated, &p_replicated);
00162 for (unsigned i=0; i<size; i++)
00163 {
00164 mData[i] = p_replicated[i];
00165 }
00166 }