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