PetscVecTools.cpp

00001 /*
00002 
00003 Copyright (C) University of Oxford, 2005-2011
00004 
00005 University of Oxford means the Chancellor, Masters and Scholars of the
00006 University of Oxford, having an administrative office at Wellington
00007 Square, Oxford OX1 2JD, UK.
00008 
00009 This file is part of Chaste.
00010 
00011 Chaste is free software: you can redistribute it and/or modify it
00012 under the terms of the GNU Lesser General Public License as published
00013 by the Free Software Foundation, either version 2.1 of the License, or
00014 (at your option) any later version.
00015 
00016 Chaste is distributed in the hope that it will be useful, but WITHOUT
00017 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018 FITNESS FOR A PARTICULAR PURPOSE.  See the GNU Lesser General Public
00019 License for more details. The offer of Chaste under the terms of the
00020 License is subject to the License being interpreted in accordance with
00021 English Law and subject to any action against the University of Oxford
00022 being under the jurisdiction of the English Courts.
00023 
00024 You should have received a copy of the GNU Lesser General Public License
00025 along with Chaste. If not, see <http://www.gnu.org/licenses/>.
00026 
00027 */
00028 
00029 #include "PetscVecTools.hpp"
00030 #include "PetscTools.hpp"
00031 #include <petscviewer.h>
00032 #include <cassert>
00033 #include "DistributedVectorFactory.hpp"
00034 #include "DistributedVector.hpp"
00035 #include "PetscException.hpp"
00036 
00038 // Implementation
00040 
00041 void PetscVecTools::Assemble(Vec vector)
00042 {
00043     VecAssemblyBegin(vector);
00044     VecAssemblyEnd(vector);
00045 }
00046 
00047 void PetscVecTools::SetElement(Vec vector, PetscInt row, double value)
00048 {
00049     PetscInt lo, hi;
00050     GetOwnershipRange(vector, lo, hi);
00051 
00052     if (row >= lo && row < hi)
00053     {
00054         VecSetValues(vector, 1, &row, &value, INSERT_VALUES);
00055     }
00056 }
00057 
00058 void PetscVecTools::AddToElement(Vec vector, PetscInt row, double value)
00059 {
00060     PetscInt lo, hi;
00061     GetOwnershipRange(vector, lo, hi);
00062 
00063     if (row >= lo && row < hi)
00064     {
00065         VecSetValues(vector, 1, &row, &value, ADD_VALUES);
00066     }
00067 }
00068 
00069 void PetscVecTools::Display(Vec vector)
00070 {
00071     VecView(vector, PETSC_VIEWER_STDOUT_WORLD);
00072 }
00073 
00074 void PetscVecTools::Zero(Vec vector)
00075 {
00076 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00077     PetscScalar zero = 0.0;
00078     VecSet(&zero, vector);
00079 #else
00080     VecZeroEntries(vector);
00081 #endif
00082 }
00083 
00084 unsigned PetscVecTools::GetSize(Vec vector) 
00085 {
00086     PetscInt size;
00087     VecGetSize(vector, &size);
00088     return (unsigned) size;
00089 }
00090 
00091 void PetscVecTools::GetOwnershipRange(Vec vector, PetscInt& lo, PetscInt& hi)
00092 {
00093     VecGetOwnershipRange(vector, &lo, &hi);
00094 }
00095 
00096 double PetscVecTools::GetElement(Vec vector, PetscInt row)
00097 {
00098     PetscInt lo, hi;
00099     GetOwnershipRange(vector, lo, hi);
00100     assert(lo <= row && row < hi);
00101 
00102     double* p_vector;
00103     PetscInt local_index = row-lo;
00104     VecGetArray(vector, &p_vector);
00105     double answer = p_vector[local_index];
00106     VecRestoreArray(vector, &p_vector);
00107 
00108     return answer;
00109 }
00110 
00111 void PetscVecTools::AddScaledVector(Vec y, Vec x, double scaleFactor)
00112 {
00113 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00114     VecAXPY(&scaleFactor, x, y);
00115 #else
00116     VecAXPY(y, scaleFactor, x);
00117 #endif
00118 }
00119 
00120 
00121 
00122 
00123 void PetscVecTools::Scale(Vec vector, double scaleFactor)
00124 {
00125 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00126     PETSCEXCEPT( VecScale(&scaleFactor, vector) );
00127 #else
00128     PETSCEXCEPT( VecScale(vector, scaleFactor) );
00129 #endif    
00130 }
00131 
00132 void PetscVecTools::WAXPY(Vec w, double a, Vec x, Vec y)
00133 {
00134 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00135     PETSCEXCEPT( VecWAXPY(&a, x, y, w) );
00136 #else
00137     PETSCEXCEPT( VecWAXPY(w, a, x, y) );
00138 #endif
00139 }
00140 
00141 
00142 void PetscVecTools::SetupInterleavedVectorScatterGather(Vec interleavedVec, VecScatter& rFirstVariableScatterContext, VecScatter& rSecondVariableScatterContext)
00143 {
00144     PetscInt num_rows, num_local_rows;
00145 
00146     VecGetSize(interleavedVec, &num_rows);
00147     VecGetLocalSize(interleavedVec, &num_local_rows);
00148 
00149     IS A11_rows, A22_rows;
00150     ISCreateStride(PETSC_COMM_WORLD, num_rows/2, 0, 2, &A11_rows);
00151     ISCreateStride(PETSC_COMM_WORLD, num_rows/2, 1, 2, &A22_rows);
00152 
00153     IS all_vector;
00154     ISCreateStride(PETSC_COMM_WORLD, num_rows/2, 0, 1, &all_vector);
00155 
00156     unsigned subvector_num_rows = num_rows/2;
00157     unsigned subvector_local_rows = num_local_rows/2;
00158     Vec x1_subvector = PetscTools::CreateVec(subvector_num_rows, subvector_local_rows);
00159     Vec x2_subvector = PetscTools::CreateVec(subvector_num_rows, subvector_local_rows);
00160 
00161     VecScatterCreate(interleavedVec, A11_rows, x1_subvector, all_vector, &rFirstVariableScatterContext);
00162     VecScatterCreate(interleavedVec, A22_rows, x2_subvector, all_vector, &rSecondVariableScatterContext);
00163 
00164     VecDestroy(x1_subvector);
00165     VecDestroy(x2_subvector);
00166 
00167     ISDestroy(A11_rows);
00168     ISDestroy(A22_rows);
00169     ISDestroy(all_vector);
00170 }
00171 
00172 void PetscVecTools::DoInterleavedVecScatter(Vec interleavedVec, VecScatter firstVariableScatterContext, Vec firstVariableVec, VecScatter secondVariableScatterContefirstVariableScatterContextt, Vec secondVariableVec)
00173 {
00174     PetscScalar *p_interleaved_vec;
00175     PetscScalar *p_1st_variable_vec;
00176     PetscScalar *p_2nd_variable_vec;
00177 
00178     VecGetArray(interleavedVec, &p_interleaved_vec);
00179     VecGetArray(firstVariableVec, &p_1st_variable_vec);
00180     VecGetArray(secondVariableVec, &p_2nd_variable_vec);
00181 
00182     PetscInt vec_local_size;
00183     VecGetLocalSize(interleavedVec, &vec_local_size);
00184     assert(vec_local_size%2 == 0);
00185 
00186     for (PetscInt local_index=0; local_index<vec_local_size/2; local_index++)
00187     {
00188         p_1st_variable_vec[local_index] = p_interleaved_vec[2*local_index];
00189         p_2nd_variable_vec[local_index] = p_interleaved_vec[2*local_index+1];
00190     }
00191 
00192     VecRestoreArray(interleavedVec, &p_interleaved_vec);
00193     VecRestoreArray(firstVariableVec, &p_1st_variable_vec);
00194     VecRestoreArray(secondVariableVec, &p_2nd_variable_vec);
00195 
00196 
00197 //    DistributedVectorFactory factory(vec_size/2);
00198 //
00199 //    DistributedVector dist_inter_vec = factory.CreateDistributedVector(interleavedVec);
00200 //    DistributedVector::Stripe dist_inter_vec_1st_var(dist_inter_vec, 0);
00201 //    DistributedVector::Stripe dist_inter_vec_2nd_var(dist_inter_vec, 1);
00202 //
00203 //    DistributedVector dist_1st_var_vec = factory.CreateDistributedVector(firstVariableVec);
00204 //    DistributedVector dist_2nd_var_vec = factory.CreateDistributedVector(secondVariableVec);
00205 //
00206 //    for (DistributedVector::Iterator index = dist_1st_var_vec.Begin();
00207 //         index!= dist_1st_var_vec.End();
00208 //         ++index)
00209 //    {
00210 //        dist_1st_var_vec[index] = dist_inter_vec_1st_var[index];
00211 //        dist_2nd_var_vec[index] = dist_inter_vec_2nd_var[index];
00212 //    }
00213 
00214 
00215 //    //PETSc-3.x.x or PETSc-2.3.3
00216 //    #if ( (PETSC_VERSION_MAJOR == 3) || (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 3 && PETSC_VERSION_SUBMINOR == 3)) //2.3.3 or 3.x.x
00217 //        VecScatterBegin(firstVariableScatterContext, interleavedVec, firstVariableVec, INSERT_VALUES, SCATTER_FORWARD);
00218 //        VecScatterEnd(firstVariableScatterContext, interleavedVec, firstVariableVec, INSERT_VALUES, SCATTER_FORWARD);
00219 //    #else
00220 //        VecScatterBegin(interleavedVec, firstVariableVec, INSERT_VALUES, SCATTER_FORWARD, firstVariableScatterContext);
00221 //        VecScatterEnd(interleavedVec, firstVariableVec, INSERT_VALUES, SCATTER_FORWARD, firstVariableScatterContext);
00222 //    #endif
00223 //
00224 //    //PETSc-3.x.x or PETSc-2.3.3
00225 //    #if ( (PETSC_VERSION_MAJOR == 3) || (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 3 && PETSC_VERSION_SUBMINOR == 3)) //2.3.3 or 3.x.x
00226 //        VecScatterBegin(secondVariableScatterContefirstVariableScatterContextt, interleavedVec, secondVariableVec, INSERT_VALUES, SCATTER_FORWARD);
00227 //        VecScatterEnd(secondVariableScatterContefirstVariableScatterContextt, interleavedVec, secondVariableVec, INSERT_VALUES, SCATTER_FORWARD);
00228 //    #else
00229 //        VecScatterBegin(interleavedVec, secondVariableVec, INSERT_VALUES, SCATTER_FORWARD, secondVariableScatterContefirstVariableScatterContextt);
00230 //        VecScatterEnd(interleavedVec, secondVariableVec, INSERT_VALUES, SCATTER_FORWARD, secondVariableScatterContefirstVariableScatterContextt);
00231 //   #endif
00232 }
00233 
00234 void PetscVecTools::DoInterleavedVecGather(Vec interleavedVec, VecScatter firstVariableScatterContext, Vec firstVariableVec, VecScatter secondVariableScatterContext, Vec secondVariableVec)
00235 {
00236     PetscScalar *p_interleaved_vec;
00237     PetscScalar *p_1st_variable_vec;
00238     PetscScalar *p_2nd_variable_vec;
00239 
00240     VecGetArray(interleavedVec, &p_interleaved_vec);
00241     VecGetArray(firstVariableVec, &p_1st_variable_vec);
00242     VecGetArray(secondVariableVec, &p_2nd_variable_vec);
00243 
00244     PetscInt vec_local_size;
00245     VecGetLocalSize(interleavedVec, &vec_local_size);
00246     assert(vec_local_size%2 == 0);
00247 
00248     for (PetscInt local_index=0; local_index<vec_local_size/2; local_index++)
00249     {
00250         p_interleaved_vec[2*local_index] = p_1st_variable_vec[local_index];
00251         p_interleaved_vec[2*local_index+1] = p_2nd_variable_vec[local_index];
00252     }
00253 
00254     VecRestoreArray(interleavedVec, &p_interleaved_vec);
00255     VecRestoreArray(firstVariableVec, &p_1st_variable_vec);
00256     VecRestoreArray(secondVariableVec, &p_2nd_variable_vec);
00257 
00258 //#if ( (PETSC_VERSION_MAJOR == 3) || (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 3 && PETSC_VERSION_SUBMINOR == 3)) //2.3.3 or 3.x.x
00259 //    VecScatterBegin(firstVariableScatterContext, firstVariableVec, interleavedVec, INSERT_VALUES, SCATTER_REVERSE);
00260 //    VecScatterEnd(firstVariableScatterContext, firstVariableVec, interleavedVec, INSERT_VALUES, SCATTER_REVERSE);
00261 //#else
00262 //    VecScatterBegin(firstVariableVec, interleavedVec, INSERT_VALUES, SCATTER_REVERSE, firstVariableScatterContext);
00263 //    VecScatterEnd(firstVariableVec, interleavedVec, INSERT_VALUES, SCATTER_REVERSE, firstVariableScatterContext);
00264 //#endif
00265 //
00267 //#if ( (PETSC_VERSION_MAJOR == 3) || (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 3 && PETSC_VERSION_SUBMINOR == 3)) //2.3.3 or 3.x.x
00268 //    VecScatterBegin(secondVariableScatterContext, secondVariableVec, interleavedVec, INSERT_VALUES, SCATTER_REVERSE);
00269 //    VecScatterEnd(secondVariableScatterContext, secondVariableVec, interleavedVec, INSERT_VALUES, SCATTER_REVERSE);
00270 //#else
00271 //    VecScatterBegin(secondVariableVec, interleavedVec, INSERT_VALUES, SCATTER_REVERSE, secondVariableScatterContext);
00272 //    VecScatterEnd(secondVariableVec, interleavedVec, INSERT_VALUES, SCATTER_REVERSE, secondVariableScatterContext);
00273 //#endif
00274 }

Generated on Mon Apr 18 11:35:28 2011 for Chaste by  doxygen 1.5.5