petsc_vector.C
Go to the documentation of this file.00001 // The libMesh Finite Element Library. 00002 // Copyright (C) 2002-2012 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner 00003 00004 // This library is free software; you can redistribute it and/or 00005 // modify it under the terms of the GNU Lesser General Public 00006 // License as published by the Free Software Foundation; either 00007 // version 2.1 of the License, or (at your option) any later version. 00008 00009 // This library is distributed in the hope that it will be useful, 00010 // but WITHOUT ANY WARRANTY; without even the implied warranty of 00011 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU 00012 // Lesser General Public License for more details. 00013 00014 // You should have received a copy of the GNU Lesser General Public 00015 // License along with this library; if not, write to the Free Software 00016 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA 00017 00018 00019 00020 // C++ includes 00021 00022 // Local Includes 00023 #include "libmesh/petsc_vector.h" 00024 #include "libmesh/petsc_matrix.h" 00025 00026 #ifdef LIBMESH_HAVE_PETSC 00027 00028 #include "libmesh/dense_subvector.h" 00029 #include "libmesh/dense_vector.h" 00030 #include "libmesh/parallel.h" 00031 #include "libmesh/petsc_macro.h" 00032 #include "libmesh/utility.h" 00033 00034 namespace libMesh 00035 { 00036 00037 00038 00039 //----------------------------------------------------------------------- 00040 // PetscVector members 00041 00042 // void PetscVector<T>::init (const NumericVector<T>& v, const bool fast) 00043 // { 00044 // libmesh_error(); 00045 00046 // init (v.local_size(), v.size(), fast); 00047 00048 // vec = libmesh_cast_ref<const PetscVector<T>&>(v).vec; 00049 // } 00050 00051 template <typename T> 00052 T PetscVector<T>::sum () const 00053 { 00054 this->_restore_array(); 00055 libmesh_assert(this->closed()); 00056 00057 PetscErrorCode ierr=0; 00058 PetscScalar value=0.; 00059 00060 ierr = VecSum (_vec, &value); 00061 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00062 00063 return static_cast<T>(value); 00064 } 00065 00066 00067 template <typename T> 00068 Real PetscVector<T>::l1_norm () const 00069 { 00070 this->_restore_array(); 00071 libmesh_assert(this->closed()); 00072 00073 PetscErrorCode ierr=0; 00074 PetscReal value=0.; 00075 00076 ierr = VecNorm (_vec, NORM_1, &value); 00077 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00078 00079 return static_cast<Real>(value); 00080 } 00081 00082 00083 00084 template <typename T> 00085 Real PetscVector<T>::l2_norm () const 00086 { 00087 this->_restore_array(); 00088 libmesh_assert(this->closed()); 00089 00090 PetscErrorCode ierr=0; 00091 PetscReal value=0.; 00092 00093 ierr = VecNorm (_vec, NORM_2, &value); 00094 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00095 00096 return static_cast<Real>(value); 00097 } 00098 00099 00100 00101 00102 template <typename T> 00103 Real PetscVector<T>::linfty_norm () const 00104 { 00105 this->_restore_array(); 00106 libmesh_assert(this->closed()); 00107 00108 PetscErrorCode ierr=0; 00109 PetscReal value=0.; 00110 00111 ierr = VecNorm (_vec, NORM_INFINITY, &value); 00112 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00113 00114 return static_cast<Real>(value); 00115 } 00116 00117 00118 00119 00120 template <typename T> 00121 NumericVector<T>& 00122 PetscVector<T>::operator += (const NumericVector<T>& v) 00123 { 00124 this->_restore_array(); 00125 libmesh_assert(this->closed()); 00126 00127 this->add(1., v); 00128 00129 return *this; 00130 } 00131 00132 00133 00134 template <typename T> 00135 NumericVector<T>& 00136 PetscVector<T>::operator -= (const NumericVector<T>& v) 00137 { 00138 this->_restore_array(); 00139 libmesh_assert(this->closed()); 00140 00141 this->add(-1., v); 00142 00143 return *this; 00144 } 00145 00146 00147 00148 template <typename T> 00149 void PetscVector<T>::set (const numeric_index_type i, const T value) 00150 { 00151 this->_restore_array(); 00152 libmesh_assert_less (i, size()); 00153 00154 PetscErrorCode ierr=0; 00155 PetscInt i_val = static_cast<PetscInt>(i); 00156 PetscScalar petsc_value = static_cast<PetscScalar>(value); 00157 00158 ierr = VecSetValues (_vec, 1, &i_val, &petsc_value, INSERT_VALUES); 00159 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00160 00161 this->_is_closed = false; 00162 } 00163 00164 00165 00166 template <typename T> 00167 void PetscVector<T>::reciprocal() 00168 { 00169 PetscErrorCode ierr = 0; 00170 00171 // VecReciprocal has been in PETSc since at least 2.3.3 days 00172 ierr = VecReciprocal(_vec); 00173 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00174 } 00175 00176 00177 00178 template <typename T> 00179 void PetscVector<T>::add (const numeric_index_type i, const T value) 00180 { 00181 this->_restore_array(); 00182 libmesh_assert_less (i, size()); 00183 00184 PetscErrorCode ierr=0; 00185 PetscInt i_val = static_cast<PetscInt>(i); 00186 PetscScalar petsc_value = static_cast<PetscScalar>(value); 00187 00188 ierr = VecSetValues (_vec, 1, &i_val, &petsc_value, ADD_VALUES); 00189 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00190 00191 this->_is_closed = false; 00192 } 00193 00194 00195 00196 template <typename T> 00197 void PetscVector<T>::add_vector (const std::vector<T>& v, 00198 const std::vector<numeric_index_type>& dof_indices) 00199 { 00200 this->_restore_array(); 00201 libmesh_assert_equal_to (v.size(), dof_indices.size()); 00202 00203 for (unsigned int i=0; i<v.size(); i++) 00204 this->add (dof_indices[i], v[i]); 00205 } 00206 00207 00208 00209 template <typename T> 00210 void PetscVector<T>::add_vector (const NumericVector<T>& V, 00211 const std::vector<numeric_index_type>& dof_indices) 00212 { 00213 libmesh_assert_equal_to (V.size(), dof_indices.size()); 00214 00215 for (unsigned int i=0; i<V.size(); i++) 00216 this->add (dof_indices[i], V(i)); 00217 } 00218 00219 00220 00221 template <typename T> 00222 void PetscVector<T>::add_vector (const NumericVector<T>& V_in, 00223 const SparseMatrix<T>& A_in) 00224 { 00225 this->_restore_array(); 00226 // Make sure the data passed in are really of Petsc types 00227 const PetscVector<T>* V = libmesh_cast_ptr<const PetscVector<T>*>(&V_in); 00228 const PetscMatrix<T>* A = libmesh_cast_ptr<const PetscMatrix<T>*>(&A_in); 00229 00230 PetscErrorCode ierr=0; 00231 00232 A->close(); 00233 00234 // The const_cast<> is not elegant, but it is required since PETSc 00235 // is not const-correct. 00236 ierr = MatMultAdd(const_cast<PetscMatrix<T>*>(A)->mat(), V->_vec, _vec, _vec); 00237 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00238 } 00239 00240 00241 00242 template <typename T> 00243 void PetscVector<T>::add_vector (const DenseVector<T>& V, 00244 const std::vector<numeric_index_type>& dof_indices) 00245 { 00246 libmesh_assert_equal_to (V.size(), dof_indices.size()); 00247 00248 for (unsigned int i=0; i<V.size(); i++) 00249 this->add (dof_indices[i], V(i)); 00250 } 00251 00252 00253 template <typename T> 00254 void PetscVector<T>::add_vector_transpose (const NumericVector<T>& V_in, 00255 const SparseMatrix<T>& A_in) 00256 { 00257 this->_restore_array(); 00258 // Make sure the data passed in are really of Petsc types 00259 const PetscVector<T>* V = libmesh_cast_ptr<const PetscVector<T>*>(&V_in); 00260 const PetscMatrix<T>* A = libmesh_cast_ptr<const PetscMatrix<T>*>(&A_in); 00261 00262 PetscErrorCode ierr=0; 00263 00264 A->close(); 00265 00266 // The const_cast<> is not elegant, but it is required since PETSc 00267 // is not const-correct. 00268 ierr = MatMultTransposeAdd(const_cast<PetscMatrix<T>*>(A)->mat(), V->_vec, _vec, _vec); 00269 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00270 } 00271 00272 00273 #if PETSC_VERSION_LESS_THAN(3,1,0) 00274 template <typename T> 00275 void PetscVector<T>::add_vector_conjugate_transpose (const NumericVector<T>&, 00276 const SparseMatrix<T>&) 00277 { 00278 00279 libMesh::out << "MatMultHermitianTranspose was introduced in PETSc 3.1.0," 00280 << "No one has made it backwards compatible with older " 00281 << "versions of PETSc so far." << std::endl; 00282 libmesh_error(); 00283 } 00284 00285 #else 00286 00287 template <typename T> 00288 void PetscVector<T>::add_vector_conjugate_transpose (const NumericVector<T>& V_in, 00289 const SparseMatrix<T>& A_in) 00290 { 00291 this->_restore_array(); 00292 // Make sure the data passed in are really of Petsc types 00293 const PetscVector<T>* V = libmesh_cast_ptr<const PetscVector<T>*>(&V_in); 00294 const PetscMatrix<T>* A = libmesh_cast_ptr<const PetscMatrix<T>*>(&A_in); 00295 00296 A->close(); 00297 00298 // Store a temporary copy since MatMultHermitianTransposeAdd doesn't seem to work 00299 // TODO: Find out why MatMultHermitianTransposeAdd doesn't work, might be a PETSc bug? 00300 AutoPtr< NumericVector<Number> > this_clone = this->clone(); 00301 00302 // The const_cast<> is not elegant, but it is required since PETSc 00303 // is not const-correct. 00304 PetscErrorCode ierr = MatMultHermitianTranspose(const_cast<PetscMatrix<T>*>(A)->mat(), V->_vec, _vec); 00305 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00306 00307 // Add the temporary copy to the matvec result 00308 this->add(1., *this_clone); 00309 } 00310 #endif 00311 00312 00313 template <typename T> 00314 void PetscVector<T>::add (const T v_in) 00315 { 00316 this->_restore_array(); 00317 00318 PetscErrorCode ierr=0; 00319 PetscScalar* values; 00320 const PetscScalar v = static_cast<PetscScalar>(v_in); 00321 00322 if(this->type() != GHOSTED) 00323 { 00324 const PetscInt n = static_cast<PetscInt>(this->local_size()); 00325 const PetscInt fli = static_cast<PetscInt>(this->first_local_index()); 00326 00327 for (PetscInt i=0; i<n; i++) 00328 { 00329 ierr = VecGetArray (_vec, &values); 00330 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00331 00332 PetscInt ig = fli + i; 00333 00334 PetscScalar value = (values[i] + v); 00335 00336 ierr = VecRestoreArray (_vec, &values); 00337 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00338 00339 ierr = VecSetValues (_vec, 1, &ig, &value, INSERT_VALUES); 00340 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00341 } 00342 } 00343 else 00344 { 00345 /* Vectors that include ghost values require a special 00346 handling. */ 00347 Vec loc_vec; 00348 ierr = VecGhostGetLocalForm (_vec,&loc_vec); 00349 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00350 00351 PetscInt n=0; 00352 ierr = VecGetSize(loc_vec, &n); 00353 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00354 00355 for (PetscInt i=0; i<n; i++) 00356 { 00357 ierr = VecGetArray (loc_vec, &values); 00358 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00359 00360 PetscScalar value = (values[i] + v); 00361 00362 ierr = VecRestoreArray (loc_vec, &values); 00363 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00364 00365 ierr = VecSetValues (loc_vec, 1, &i, &value, INSERT_VALUES); 00366 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00367 } 00368 00369 ierr = VecGhostRestoreLocalForm (_vec,&loc_vec); 00370 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00371 } 00372 00373 this->_is_closed = false; 00374 } 00375 00376 00377 00378 template <typename T> 00379 void PetscVector<T>::add (const NumericVector<T>& v) 00380 { 00381 this->add (1., v); 00382 } 00383 00384 00385 00386 template <typename T> 00387 void PetscVector<T>::add (const T a_in, const NumericVector<T>& v_in) 00388 { 00389 this->_restore_array(); 00390 00391 PetscErrorCode ierr = 0; 00392 PetscScalar a = static_cast<PetscScalar>(a_in); 00393 00394 // Make sure the NumericVector passed in is really a PetscVector 00395 const PetscVector<T>* v = libmesh_cast_ptr<const PetscVector<T>*>(&v_in); 00396 v->_restore_array(); 00397 00398 libmesh_assert_equal_to (this->size(), v->size()); 00399 00400 if(this->type() != GHOSTED) 00401 { 00402 #if PETSC_VERSION_LESS_THAN(2,3,0) 00403 // 2.2.x & earlier style 00404 ierr = VecAXPY(&a, v->_vec, _vec); 00405 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00406 #else 00407 // 2.3.x & later style 00408 ierr = VecAXPY(_vec, a, v->_vec); 00409 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00410 #endif 00411 } 00412 else 00413 { 00414 Vec loc_vec; 00415 Vec v_loc_vec; 00416 ierr = VecGhostGetLocalForm (_vec,&loc_vec); 00417 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00418 ierr = VecGhostGetLocalForm (v->_vec,&v_loc_vec); 00419 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00420 00421 #if PETSC_VERSION_LESS_THAN(2,3,0) 00422 // 2.2.x & earlier style 00423 ierr = VecAXPY(&a, v_loc_vec, loc_vec); 00424 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00425 #else 00426 // 2.3.x & later style 00427 ierr = VecAXPY(loc_vec, a, v_loc_vec); 00428 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00429 #endif 00430 00431 ierr = VecGhostRestoreLocalForm (v->_vec,&v_loc_vec); 00432 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00433 ierr = VecGhostRestoreLocalForm (_vec,&loc_vec); 00434 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00435 } 00436 } 00437 00438 00439 00440 template <typename T> 00441 void PetscVector<T>::insert (const std::vector<T>& v, 00442 const std::vector<numeric_index_type>& dof_indices) 00443 { 00444 libmesh_assert_equal_to (v.size(), dof_indices.size()); 00445 00446 for (unsigned int i=0; i<v.size(); i++) 00447 this->set (dof_indices[i], v[i]); 00448 } 00449 00450 00451 00452 template <typename T> 00453 void PetscVector<T>::insert (const NumericVector<T>& V, 00454 const std::vector<numeric_index_type>& dof_indices) 00455 { 00456 libmesh_assert_equal_to (V.size(), dof_indices.size()); 00457 00458 for (unsigned int i=0; i<V.size(); i++) 00459 this->set (dof_indices[i], V(i)); 00460 } 00461 00462 00463 00464 template <typename T> 00465 void PetscVector<T>::insert (const DenseVector<T>& V, 00466 const std::vector<numeric_index_type>& dof_indices) 00467 { 00468 libmesh_assert_equal_to (V.size(), dof_indices.size()); 00469 00470 for (unsigned int i=0; i<V.size(); i++) 00471 this->set (dof_indices[i], V(i)); 00472 } 00473 00474 00475 00476 template <typename T> 00477 void PetscVector<T>::insert (const DenseSubVector<T>& V, 00478 const std::vector<numeric_index_type>& dof_indices) 00479 { 00480 libmesh_assert_equal_to (V.size(), dof_indices.size()); 00481 00482 for (unsigned int i=0; i<V.size(); i++) 00483 this->set (dof_indices[i], V(i)); 00484 } 00485 00486 00487 00488 template <typename T> 00489 void PetscVector<T>::scale (const T factor_in) 00490 { 00491 this->_restore_array(); 00492 00493 PetscErrorCode ierr = 0; 00494 PetscScalar factor = static_cast<PetscScalar>(factor_in); 00495 00496 if(this->type() != GHOSTED) 00497 { 00498 #if PETSC_VERSION_LESS_THAN(2,3,0) 00499 // 2.2.x & earlier style 00500 ierr = VecScale(&factor, _vec); 00501 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00502 #else 00503 // 2.3.x & later style 00504 ierr = VecScale(_vec, factor); 00505 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00506 #endif 00507 } 00508 else 00509 { 00510 Vec loc_vec; 00511 ierr = VecGhostGetLocalForm (_vec,&loc_vec); 00512 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00513 00514 #if PETSC_VERSION_LESS_THAN(2,3,0) 00515 // 2.2.x & earlier style 00516 ierr = VecScale(&factor, loc_vec); 00517 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00518 #else 00519 // 2.3.x & later style 00520 ierr = VecScale(loc_vec, factor); 00521 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00522 #endif 00523 00524 ierr = VecGhostRestoreLocalForm (_vec,&loc_vec); 00525 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00526 } 00527 } 00528 00529 template <typename T> 00530 void PetscVector<T>::abs() 00531 { 00532 this->_restore_array(); 00533 00534 PetscErrorCode ierr = 0; 00535 00536 if(this->type() != GHOSTED) 00537 { 00538 ierr = VecAbs(_vec); 00539 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00540 } 00541 else 00542 { 00543 Vec loc_vec; 00544 ierr = VecGhostGetLocalForm (_vec,&loc_vec); 00545 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00546 00547 ierr = VecAbs(loc_vec); 00548 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00549 00550 ierr = VecGhostRestoreLocalForm (_vec,&loc_vec); 00551 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00552 } 00553 } 00554 00555 template <typename T> 00556 T PetscVector<T>::dot (const NumericVector<T>& V) const 00557 { 00558 this->_restore_array(); 00559 00560 // Error flag 00561 PetscErrorCode ierr = 0; 00562 00563 // Return value 00564 PetscScalar value=0.; 00565 00566 // Make sure the NumericVector passed in is really a PetscVector 00567 const PetscVector<T>* v = libmesh_cast_ptr<const PetscVector<T>*>(&V); 00568 00569 // 2.3.x (at least) style. Untested for previous versions. 00570 ierr = VecDot(this->_vec, v->_vec, &value); 00571 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00572 00573 return static_cast<T>(value); 00574 } 00575 00576 template <typename T> 00577 T PetscVector<T>::indefinite_dot (const NumericVector<T>& V) const 00578 { 00579 this->_restore_array(); 00580 00581 // Error flag 00582 PetscErrorCode ierr = 0; 00583 00584 // Return value 00585 PetscScalar value=0.; 00586 00587 // Make sure the NumericVector passed in is really a PetscVector 00588 const PetscVector<T>* v = libmesh_cast_ptr<const PetscVector<T>*>(&V); 00589 00590 // 2.3.x (at least) style. Untested for previous versions. 00591 ierr = VecTDot(this->_vec, v->_vec, &value); 00592 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00593 00594 return static_cast<T>(value); 00595 } 00596 00597 00598 template <typename T> 00599 NumericVector<T>& 00600 PetscVector<T>::operator = (const T s_in) 00601 { 00602 this->_restore_array(); 00603 libmesh_assert(this->closed()); 00604 00605 PetscErrorCode ierr = 0; 00606 PetscScalar s = static_cast<PetscScalar>(s_in); 00607 00608 if (this->size() != 0) 00609 { 00610 if(this->type() != GHOSTED) 00611 { 00612 #if PETSC_VERSION_LESS_THAN(2,3,0) 00613 // 2.2.x & earlier style 00614 ierr = VecSet(&s, _vec); 00615 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00616 #else 00617 // 2.3.x & later style 00618 ierr = VecSet(_vec, s); 00619 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00620 #endif 00621 } 00622 else 00623 { 00624 Vec loc_vec; 00625 ierr = VecGhostGetLocalForm (_vec,&loc_vec); 00626 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00627 00628 #if PETSC_VERSION_LESS_THAN(2,3,0) 00629 // 2.2.x & earlier style 00630 ierr = VecSet(&s, loc_vec); 00631 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00632 #else 00633 // 2.3.x & later style 00634 ierr = VecSet(loc_vec, s); 00635 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00636 #endif 00637 00638 ierr = VecGhostRestoreLocalForm (_vec,&loc_vec); 00639 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00640 } 00641 } 00642 00643 return *this; 00644 } 00645 00646 00647 00648 template <typename T> 00649 NumericVector<T>& 00650 PetscVector<T>::operator = (const NumericVector<T>& v_in) 00651 { 00652 // Make sure the NumericVector passed in is really a PetscVector 00653 const PetscVector<T>* v = libmesh_cast_ptr<const PetscVector<T>*>(&v_in); 00654 00655 *this = *v; 00656 00657 return *this; 00658 } 00659 00660 00661 00662 template <typename T> 00663 PetscVector<T>& 00664 PetscVector<T>::operator = (const PetscVector<T>& v) 00665 { 00666 this->_restore_array(); 00667 v._restore_array(); 00668 00669 libmesh_assert_equal_to (this->size(), v.size()); 00670 libmesh_assert_equal_to (this->local_size(), v.local_size()); 00671 libmesh_assert (v.closed()); 00672 00673 PetscErrorCode ierr = 0; 00674 00675 if (((this->type()==PARALLEL) && (v.type()==GHOSTED)) || 00676 ((this->type()==GHOSTED) && (v.type()==PARALLEL)) || 00677 ((this->type()==GHOSTED) && (v.type()==SERIAL)) || 00678 ((this->type()==SERIAL) && (v.type()==GHOSTED))) 00679 { 00680 /* Allow assignment of a ghosted to a parallel vector since this 00681 causes no difficulty. See discussion in libmesh-devel of 00682 June 24, 2010. */ 00683 ierr = VecCopy (v._vec, this->_vec); 00684 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00685 } 00686 else 00687 { 00688 /* In all other cases, we assert that both vectors are of equal 00689 type. */ 00690 libmesh_assert_equal_to (this->_type, v._type); 00691 libmesh_assert (this->_global_to_local_map == v._global_to_local_map); 00692 00693 if (v.size() != 0) 00694 { 00695 if(this->type() != GHOSTED) 00696 { 00697 ierr = VecCopy (v._vec, this->_vec); 00698 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00699 } 00700 else 00701 { 00702 Vec loc_vec; 00703 Vec v_loc_vec; 00704 ierr = VecGhostGetLocalForm (_vec,&loc_vec); 00705 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00706 ierr = VecGhostGetLocalForm (v._vec,&v_loc_vec); 00707 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00708 00709 ierr = VecCopy (v_loc_vec, loc_vec); 00710 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00711 00712 ierr = VecGhostRestoreLocalForm (v._vec,&v_loc_vec); 00713 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00714 ierr = VecGhostRestoreLocalForm (_vec,&loc_vec); 00715 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00716 } 00717 } 00718 } 00719 00720 close(); 00721 00722 return *this; 00723 } 00724 00725 00726 00727 template <typename T> 00728 NumericVector<T>& 00729 PetscVector<T>::operator = (const std::vector<T>& v) 00730 { 00731 this->_restore_array(); 00732 00733 const numeric_index_type nl = this->local_size(); 00734 const numeric_index_type ioff = this->first_local_index(); 00735 PetscErrorCode ierr=0; 00736 PetscScalar* values; 00737 00742 if (this->size() == v.size()) 00743 { 00744 ierr = VecGetArray (_vec, &values); 00745 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00746 00747 for (numeric_index_type i=0; i<nl; i++) 00748 values[i] = static_cast<PetscScalar>(v[i+ioff]); 00749 00750 ierr = VecRestoreArray (_vec, &values); 00751 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00752 } 00753 00758 else 00759 { 00760 libmesh_assert_equal_to (this->local_size(), v.size()); 00761 00762 ierr = VecGetArray (_vec, &values); 00763 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00764 00765 for (numeric_index_type i=0; i<nl; i++) 00766 values[i] = static_cast<PetscScalar>(v[i]); 00767 00768 ierr = VecRestoreArray (_vec, &values); 00769 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00770 } 00771 00772 // Make sure ghost dofs are up to date 00773 if (this->type() == GHOSTED) 00774 this->close(); 00775 00776 return *this; 00777 } 00778 00779 00780 00781 template <typename T> 00782 void PetscVector<T>::localize (NumericVector<T>& v_local_in) const 00783 { 00784 this->_restore_array(); 00785 00786 // Make sure the NumericVector passed in is really a PetscVector 00787 PetscVector<T>* v_local = libmesh_cast_ptr<PetscVector<T>*>(&v_local_in); 00788 00789 libmesh_assert(v_local); 00790 libmesh_assert_equal_to (v_local->size(), this->size()); 00791 00792 PetscErrorCode ierr = 0; 00793 const PetscInt n = this->size(); 00794 00795 IS is; 00796 VecScatter scatter; 00797 00798 // Create idx, idx[i] = i; 00799 std::vector<PetscInt> idx(n); Utility::iota (idx.begin(), idx.end(), 0); 00800 00801 // Create the index set & scatter object 00802 ierr = ISCreateLibMesh(libMesh::COMM_WORLD, n, &idx[0], PETSC_USE_POINTER, &is); 00803 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00804 00805 ierr = VecScatterCreate(_vec, is, 00806 v_local->_vec, is, 00807 &scatter); 00808 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00809 00810 // Perform the scatter 00811 #if PETSC_VERSION_LESS_THAN(2,3,3) 00812 00813 ierr = VecScatterBegin(_vec, v_local->_vec, INSERT_VALUES, 00814 SCATTER_FORWARD, scatter); 00815 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00816 00817 ierr = VecScatterEnd (_vec, v_local->_vec, INSERT_VALUES, 00818 SCATTER_FORWARD, scatter); 00819 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00820 #else 00821 // API argument order change in PETSc 2.3.3 00822 ierr = VecScatterBegin(scatter, _vec, v_local->_vec, 00823 INSERT_VALUES, SCATTER_FORWARD); 00824 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00825 00826 ierr = VecScatterEnd (scatter, _vec, v_local->_vec, 00827 INSERT_VALUES, SCATTER_FORWARD); 00828 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00829 #endif 00830 00831 // Clean up 00832 ierr = LibMeshISDestroy (&is); 00833 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00834 00835 ierr = LibMeshVecScatterDestroy(&scatter); 00836 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00837 00838 // Make sure ghost dofs are up to date 00839 if (v_local->type() == GHOSTED) 00840 v_local->close(); 00841 } 00842 00843 00844 00845 template <typename T> 00846 void PetscVector<T>::localize (NumericVector<T>& v_local_in, 00847 const std::vector<numeric_index_type>& send_list) const 00848 { 00849 // FIXME: Workaround for a strange bug at large-scale. 00850 // If we have ghosting, PETSc lets us just copy the solution, and 00851 // doing so avoids a segfault? 00852 if (v_local_in.type() == GHOSTED && 00853 this->type() == PARALLEL) 00854 { 00855 v_local_in = *this; 00856 return; 00857 } 00858 00859 // Normal code path begins here 00860 00861 this->_restore_array(); 00862 00863 // Make sure the NumericVector passed in is really a PetscVector 00864 PetscVector<T>* v_local = libmesh_cast_ptr<PetscVector<T>*>(&v_local_in); 00865 00866 libmesh_assert(v_local); 00867 libmesh_assert_equal_to (v_local->size(), this->size()); 00868 libmesh_assert_less_equal (send_list.size(), v_local->size()); 00869 00870 PetscErrorCode ierr=0; 00871 const numeric_index_type n_sl = send_list.size(); 00872 00873 IS is; 00874 VecScatter scatter; 00875 00876 std::vector<PetscInt> idx(n_sl + this->local_size()); 00877 00878 for (numeric_index_type i=0; i<n_sl; i++) 00879 idx[i] = static_cast<PetscInt>(send_list[i]); 00880 for (numeric_index_type i = 0; i != this->local_size(); ++i) 00881 idx[n_sl+i] = i + this->first_local_index(); 00882 00883 // Create the index set & scatter object 00884 if (idx.empty()) 00885 ierr = ISCreateLibMesh(libMesh::COMM_WORLD, 00886 n_sl+this->local_size(), PETSC_NULL, PETSC_USE_POINTER, &is); 00887 else 00888 ierr = ISCreateLibMesh(libMesh::COMM_WORLD, 00889 n_sl+this->local_size(), &idx[0], PETSC_USE_POINTER, &is); 00890 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00891 00892 ierr = VecScatterCreate(_vec, is, 00893 v_local->_vec, is, 00894 &scatter); 00895 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00896 00897 00898 // Perform the scatter 00899 #if PETSC_VERSION_LESS_THAN(2,3,3) 00900 00901 ierr = VecScatterBegin(_vec, v_local->_vec, INSERT_VALUES, 00902 SCATTER_FORWARD, scatter); 00903 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00904 00905 ierr = VecScatterEnd (_vec, v_local->_vec, INSERT_VALUES, 00906 SCATTER_FORWARD, scatter); 00907 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00908 00909 #else 00910 00911 // API argument order change in PETSc 2.3.3 00912 ierr = VecScatterBegin(scatter, _vec, v_local->_vec, 00913 INSERT_VALUES, SCATTER_FORWARD); 00914 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00915 00916 ierr = VecScatterEnd (scatter, _vec, v_local->_vec, 00917 INSERT_VALUES, SCATTER_FORWARD); 00918 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00919 00920 #endif 00921 00922 00923 // Clean up 00924 ierr = LibMeshISDestroy (&is); 00925 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00926 00927 ierr = LibMeshVecScatterDestroy(&scatter); 00928 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00929 00930 // Make sure ghost dofs are up to date 00931 if (v_local->type() == GHOSTED) 00932 v_local->close(); 00933 } 00934 00935 00936 template <typename T> 00937 void PetscVector<T>::localize (const numeric_index_type first_local_idx, 00938 const numeric_index_type last_local_idx, 00939 const std::vector<numeric_index_type>& send_list) 00940 { 00941 this->_restore_array(); 00942 00943 libmesh_assert_less_equal (send_list.size(), this->size()); 00944 libmesh_assert_less_equal (last_local_idx+1, this->size()); 00945 00946 const numeric_index_type my_size = this->size(); 00947 const numeric_index_type my_local_size = (last_local_idx - first_local_idx + 1); 00948 PetscErrorCode ierr=0; 00949 00950 // Don't bother for serial cases 00951 // if ((first_local_idx == 0) && 00952 // (my_local_size == my_size)) 00953 // But we do need to stay in sync for degenerate cases 00954 if (libMesh::n_processors() == 1) 00955 return; 00956 00957 00958 // Build a parallel vector, initialize it with the local 00959 // parts of (*this) 00960 PetscVector<T> parallel_vec; 00961 00962 parallel_vec.init (my_size, my_local_size, true, PARALLEL); 00963 00964 00965 // Copy part of *this into the parallel_vec 00966 { 00967 IS is; 00968 VecScatter scatter; 00969 00970 // Create idx, idx[i] = i+first_local_idx; 00971 std::vector<PetscInt> idx(my_local_size); 00972 Utility::iota (idx.begin(), idx.end(), first_local_idx); 00973 00974 // Create the index set & scatter object 00975 ierr = ISCreateLibMesh(libMesh::COMM_WORLD, my_local_size, 00976 my_local_size ? &idx[0] : NULL, PETSC_USE_POINTER, &is); 00977 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00978 00979 ierr = VecScatterCreate(_vec, is, 00980 parallel_vec._vec, is, 00981 &scatter); 00982 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00983 00984 // Perform the scatter 00985 #if PETSC_VERSION_LESS_THAN(2,3,3) 00986 00987 ierr = VecScatterBegin(_vec, parallel_vec._vec, INSERT_VALUES, 00988 SCATTER_FORWARD, scatter); 00989 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00990 00991 ierr = VecScatterEnd (_vec, parallel_vec._vec, INSERT_VALUES, 00992 SCATTER_FORWARD, scatter); 00993 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00994 00995 #else 00996 00997 // API argument order change in PETSc 2.3.3 00998 ierr = VecScatterBegin(scatter, _vec, parallel_vec._vec, 00999 INSERT_VALUES, SCATTER_FORWARD); 01000 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01001 01002 ierr = VecScatterEnd (scatter, _vec, parallel_vec._vec, 01003 INSERT_VALUES, SCATTER_FORWARD); 01004 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01005 01006 #endif 01007 01008 // Clean up 01009 ierr = LibMeshISDestroy (&is); 01010 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01011 01012 ierr = LibMeshVecScatterDestroy(&scatter); 01013 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01014 } 01015 01016 // localize like normal 01017 parallel_vec.close(); 01018 parallel_vec.localize (*this, send_list); 01019 this->close(); 01020 } 01021 01022 01023 01024 template <typename T> 01025 void PetscVector<T>::localize (std::vector<T>& v_local) const 01026 { 01027 this->_restore_array(); 01028 01029 // This function must be run on all processors at once 01030 parallel_only(); 01031 01032 PetscErrorCode ierr=0; 01033 const PetscInt n = this->size(); 01034 const PetscInt nl = this->local_size(); 01035 PetscScalar *values; 01036 01037 v_local.clear(); 01038 v_local.resize(n, 0.); 01039 01040 ierr = VecGetArray (_vec, &values); 01041 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01042 01043 numeric_index_type ioff = first_local_index(); 01044 01045 for (PetscInt i=0; i<nl; i++) 01046 v_local[i+ioff] = static_cast<T>(values[i]); 01047 01048 ierr = VecRestoreArray (_vec, &values); 01049 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01050 01051 CommWorld.sum(v_local); 01052 } 01053 01054 01055 01056 // Full specialization for Real datatypes 01057 #ifdef LIBMESH_USE_REAL_NUMBERS 01058 01059 template <> 01060 void PetscVector<Real>::localize_to_one (std::vector<Real>& v_local, 01061 const processor_id_type pid) const 01062 { 01063 this->_restore_array(); 01064 01065 PetscErrorCode ierr=0; 01066 const PetscInt n = size(); 01067 const PetscInt nl = local_size(); 01068 PetscScalar *values; 01069 01070 01071 v_local.resize(n); 01072 01073 01074 // only one processor 01075 if (n == nl) 01076 { 01077 ierr = VecGetArray (_vec, &values); 01078 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01079 01080 for (PetscInt i=0; i<n; i++) 01081 v_local[i] = static_cast<Real>(values[i]); 01082 01083 ierr = VecRestoreArray (_vec, &values); 01084 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01085 } 01086 01087 // otherwise multiple processors 01088 else 01089 { 01090 numeric_index_type ioff = this->first_local_index(); 01091 std::vector<Real> local_values (n, 0.); 01092 01093 { 01094 ierr = VecGetArray (_vec, &values); 01095 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01096 01097 for (PetscInt i=0; i<nl; i++) 01098 local_values[i+ioff] = static_cast<Real>(values[i]); 01099 01100 ierr = VecRestoreArray (_vec, &values); 01101 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01102 } 01103 01104 01105 MPI_Reduce (&local_values[0], &v_local[0], n, MPI_REAL, MPI_SUM, 01106 pid, libMesh::COMM_WORLD); 01107 } 01108 } 01109 01110 #endif 01111 01112 01113 // Full specialization for Complex datatypes 01114 #ifdef LIBMESH_USE_COMPLEX_NUMBERS 01115 01116 template <> 01117 void PetscVector<Complex>::localize_to_one (std::vector<Complex>& v_local, 01118 const processor_id_type pid) const 01119 { 01120 this->_restore_array(); 01121 01122 PetscErrorCode ierr=0; 01123 const PetscInt n = size(); 01124 const PetscInt nl = local_size(); 01125 PetscScalar *values; 01126 01127 01128 v_local.resize(n); 01129 01130 01131 for (PetscInt i=0; i<n; i++) 01132 v_local[i] = 0.; 01133 01134 // only one processor 01135 if (n == nl) 01136 { 01137 ierr = VecGetArray (_vec, &values); 01138 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01139 01140 for (PetscInt i=0; i<n; i++) 01141 v_local[i] = static_cast<Complex>(values[i]); 01142 01143 ierr = VecRestoreArray (_vec, &values); 01144 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01145 } 01146 01147 // otherwise multiple processors 01148 else 01149 { 01150 numeric_index_type ioff = this->first_local_index(); 01151 01152 /* in here the local values are stored, acting as send buffer for MPI 01153 * initialize to zero, since we collect using MPI_SUM 01154 */ 01155 std::vector<Real> real_local_values(n, 0.); 01156 std::vector<Real> imag_local_values(n, 0.); 01157 01158 { 01159 ierr = VecGetArray (_vec, &values); 01160 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01161 01162 // provide my local share to the real and imag buffers 01163 for (PetscInt i=0; i<nl; i++) 01164 { 01165 real_local_values[i+ioff] = static_cast<Complex>(values[i]).real(); 01166 imag_local_values[i+ioff] = static_cast<Complex>(values[i]).imag(); 01167 } 01168 01169 ierr = VecRestoreArray (_vec, &values); 01170 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01171 } 01172 01173 /* have buffers of the real and imaginary part of v_local. 01174 * Once MPI_Reduce() collected all the real and imaginary 01175 * parts in these std::vector<Real>, the values can be 01176 * copied to v_local 01177 */ 01178 std::vector<Real> real_v_local(n); 01179 std::vector<Real> imag_v_local(n); 01180 01181 // collect entries from other proc's in real_v_local, imag_v_local 01182 MPI_Reduce (&real_local_values[0], &real_v_local[0], n, 01183 MPI_REAL, MPI_SUM, 01184 pid, libMesh::COMM_WORLD); 01185 01186 MPI_Reduce (&imag_local_values[0], &imag_v_local[0], n, 01187 MPI_REAL, MPI_SUM, 01188 pid, libMesh::COMM_WORLD); 01189 01190 // copy real_v_local and imag_v_local to v_local 01191 for (PetscInt i=0; i<n; i++) 01192 v_local[i] = Complex(real_v_local[i], imag_v_local[i]); 01193 } 01194 } 01195 01196 #endif 01197 01198 01199 01200 template <typename T> 01201 void PetscVector<T>::pointwise_mult (const NumericVector<T>& vec1, 01202 const NumericVector<T>& vec2) 01203 { 01204 this->_restore_array(); 01205 01206 PetscErrorCode ierr = 0; 01207 01208 // Convert arguments to PetscVector*. 01209 const PetscVector<T>* vec1_petsc = libmesh_cast_ptr<const PetscVector<T>*>(&vec1); 01210 const PetscVector<T>* vec2_petsc = libmesh_cast_ptr<const PetscVector<T>*>(&vec2); 01211 01212 // Call PETSc function. 01213 01214 #if PETSC_VERSION_LESS_THAN(2,3,1) 01215 01216 libMesh::out << "This method has been developed with PETSc 2.3.1. " 01217 << "No one has made it backwards compatible with older " 01218 << "versions of PETSc so far; however, it might work " 01219 << "without any change with some older version." << std::endl; 01220 libmesh_error(); 01221 01222 #else 01223 01224 if(this->type() != GHOSTED) 01225 { 01226 ierr = VecPointwiseMult(this->vec(), 01227 const_cast<PetscVector<T>*>(vec1_petsc)->vec(), 01228 const_cast<PetscVector<T>*>(vec2_petsc)->vec()); 01229 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01230 } 01231 else 01232 { 01233 Vec loc_vec; 01234 Vec v1_loc_vec; 01235 Vec v2_loc_vec; 01236 ierr = VecGhostGetLocalForm (_vec,&loc_vec); 01237 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01238 ierr = VecGhostGetLocalForm (const_cast<PetscVector<T>*>(vec1_petsc)->vec(),&v1_loc_vec); 01239 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01240 ierr = VecGhostGetLocalForm (const_cast<PetscVector<T>*>(vec2_petsc)->vec(),&v2_loc_vec); 01241 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01242 01243 ierr = VecPointwiseMult(loc_vec,v1_loc_vec,v2_loc_vec); 01244 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01245 01246 ierr = VecGhostRestoreLocalForm (const_cast<PetscVector<T>*>(vec1_petsc)->vec(),&v1_loc_vec); 01247 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01248 ierr = VecGhostRestoreLocalForm (const_cast<PetscVector<T>*>(vec2_petsc)->vec(),&v2_loc_vec); 01249 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01250 ierr = VecGhostRestoreLocalForm (_vec,&loc_vec); 01251 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01252 } 01253 01254 #endif 01255 01256 } 01257 01258 01259 01260 template <typename T> 01261 void PetscVector<T>::print_matlab (const std::string name) const 01262 { 01263 this->_restore_array(); 01264 libmesh_assert (this->closed()); 01265 01266 PetscErrorCode ierr=0; 01267 PetscViewer petsc_viewer; 01268 01269 01270 ierr = PetscViewerCreate (libMesh::COMM_WORLD, 01271 &petsc_viewer); 01272 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01273 01278 if (name != "NULL") 01279 { 01280 ierr = PetscViewerASCIIOpen( libMesh::COMM_WORLD, 01281 name.c_str(), 01282 &petsc_viewer); 01283 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01284 01285 ierr = PetscViewerSetFormat (petsc_viewer, 01286 PETSC_VIEWER_ASCII_MATLAB); 01287 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01288 01289 ierr = VecView (_vec, petsc_viewer); 01290 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01291 } 01292 01296 else 01297 { 01298 ierr = PetscViewerSetFormat (PETSC_VIEWER_STDOUT_WORLD, 01299 PETSC_VIEWER_ASCII_MATLAB); 01300 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01301 01302 ierr = VecView (_vec, PETSC_VIEWER_STDOUT_WORLD); 01303 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01304 } 01305 01306 01310 ierr = LibMeshPetscViewerDestroy (&petsc_viewer); 01311 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01312 } 01313 01314 01315 01316 01317 01318 template <typename T> 01319 void PetscVector<T>::create_subvector(NumericVector<T>& subvector, 01320 const std::vector<numeric_index_type>& rows) const 01321 { 01322 this->_restore_array(); 01323 01324 // PETSc data structures 01325 IS parent_is, subvector_is; 01326 VecScatter scatter; 01327 PetscErrorCode ierr = 0; 01328 01329 // Make sure the passed in subvector is really a PetscVector 01330 PetscVector<T>* petsc_subvector = libmesh_cast_ptr<PetscVector<T>*>(&subvector); 01331 01332 // If the petsc_subvector is already initialized, we assume that the 01333 // user has already allocated the *correct* amount of space for it. 01334 // If not, we use the appropriate PETSc routines to initialize it. 01335 if (!petsc_subvector->initialized()) 01336 { 01337 // Initialize the petsc_subvector to have enough space to hold 01338 // the entries which will be scattered into it. Note: such an 01339 // init() function (where we let PETSc decide the number of local 01340 // entries) is not currently offered by the PetscVector 01341 // class. Should we differentiate here between sequential and 01342 // parallel vector creation based on libMesh::n_processors() ? 01343 ierr = VecCreateMPI(libMesh::COMM_WORLD, 01344 PETSC_DECIDE, // n_local 01345 rows.size(), // n_global 01346 &(petsc_subvector->_vec)); CHKERRABORT(libMesh::COMM_WORLD,ierr); 01347 01348 ierr = VecSetFromOptions (petsc_subvector->_vec); CHKERRABORT(libMesh::COMM_WORLD,ierr); 01349 01350 // Mark the subvector as initialized 01351 petsc_subvector->_is_initialized = true; 01352 } 01353 else 01354 { 01355 petsc_subvector->_restore_array(); 01356 } 01357 01358 // Use iota to fill an array with entries [0,1,2,3,4,...rows.size()] 01359 std::vector<PetscInt> idx(rows.size()); 01360 Utility::iota (idx.begin(), idx.end(), 0); 01361 01362 // Construct index sets 01363 ierr = ISCreateLibMesh(libMesh::COMM_WORLD, 01364 rows.size(), 01365 (PetscInt*) &rows[0], 01366 PETSC_USE_POINTER, 01367 &parent_is); CHKERRABORT(libMesh::COMM_WORLD,ierr); 01368 01369 ierr = ISCreateLibMesh(libMesh::COMM_WORLD, 01370 rows.size(), 01371 (PetscInt*) &idx[0], 01372 PETSC_USE_POINTER, 01373 &subvector_is); CHKERRABORT(libMesh::COMM_WORLD,ierr); 01374 01375 // Construct the scatter object 01376 ierr = VecScatterCreate(this->_vec, 01377 parent_is, 01378 petsc_subvector->_vec, 01379 subvector_is, 01380 &scatter); CHKERRABORT(libMesh::COMM_WORLD,ierr); 01381 01382 // Actually perform the scatter 01383 #if PETSC_VERSION_LESS_THAN(2,3,3) 01384 ierr = VecScatterBegin(this->_vec, 01385 petsc_subvector->_vec, 01386 INSERT_VALUES, 01387 SCATTER_FORWARD, 01388 scatter); CHKERRABORT(libMesh::COMM_WORLD,ierr); 01389 01390 ierr = VecScatterEnd(this->_vec, 01391 petsc_subvector->_vec, 01392 INSERT_VALUES, 01393 SCATTER_FORWARD, 01394 scatter); CHKERRABORT(libMesh::COMM_WORLD,ierr); 01395 #else 01396 // API argument order change in PETSc 2.3.3 01397 ierr = VecScatterBegin(scatter, 01398 this->_vec, 01399 petsc_subvector->_vec, 01400 INSERT_VALUES, 01401 SCATTER_FORWARD); CHKERRABORT(libMesh::COMM_WORLD,ierr); 01402 01403 ierr = VecScatterEnd(scatter, 01404 this->_vec, 01405 petsc_subvector->_vec, 01406 INSERT_VALUES, 01407 SCATTER_FORWARD); CHKERRABORT(libMesh::COMM_WORLD,ierr); 01408 #endif 01409 01410 // Clean up 01411 ierr = LibMeshISDestroy(&parent_is); CHKERRABORT(libMesh::COMM_WORLD,ierr); 01412 ierr = LibMeshISDestroy(&subvector_is); CHKERRABORT(libMesh::COMM_WORLD,ierr); 01413 ierr = LibMeshVecScatterDestroy(&scatter); CHKERRABORT(libMesh::COMM_WORLD,ierr); 01414 01415 } 01416 01417 01418 01419 01420 //------------------------------------------------------------------ 01421 // Explicit instantiations 01422 template class PetscVector<Number>; 01423 01424 } // namespace libMesh 01425 01426 01427 01428 #endif // #ifdef LIBMESH_HAVE_PETSC
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:48 UTC
Hosted By: