petsc_linear_solver.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 #include "libmesh/libmesh_common.h" 00019 00020 #ifdef LIBMESH_HAVE_PETSC 00021 00022 00023 // C++ includes 00024 #include <string.h> 00025 00026 // Local Includes 00027 #include "libmesh/libmesh_logging.h" 00028 #include "libmesh/petsc_linear_solver.h" 00029 #include "libmesh/shell_matrix.h" 00030 #include "libmesh/petsc_matrix.h" 00031 #include "libmesh/petsc_preconditioner.h" 00032 #include "libmesh/petsc_vector.h" 00033 00034 namespace libMesh 00035 { 00036 00037 extern "C" 00038 { 00039 #if PETSC_VERSION_LESS_THAN(2,2,1) 00040 typedef int PetscErrorCode; 00041 typedef int PetscInt; 00042 #endif 00043 00044 00045 #if PETSC_VERSION_LESS_THAN(3,0,1) && PETSC_VERSION_RELEASE 00046 PetscErrorCode __libmesh_petsc_preconditioner_setup (void * ctx) 00047 { 00048 Preconditioner<Number> * preconditioner = static_cast<Preconditioner<Number>*>(ctx); 00049 00050 if(!preconditioner->initialized()) 00051 { 00052 err<<"Preconditioner not initialized! Make sure you call init() before solve!"<<std::endl; 00053 libmesh_error(); 00054 } 00055 00056 preconditioner->setup(); 00057 00058 return 0; 00059 } 00060 00061 00062 PetscErrorCode __libmesh_petsc_preconditioner_apply(void *ctx, Vec x, Vec y) 00063 { 00064 Preconditioner<Number> * preconditioner = static_cast<Preconditioner<Number>*>(ctx); 00065 00066 PetscVector<Number> x_vec(x); 00067 PetscVector<Number> y_vec(y); 00068 00069 preconditioner->apply(x_vec,y_vec); 00070 00071 return 0; 00072 } 00073 #else 00074 PetscErrorCode __libmesh_petsc_preconditioner_setup (PC pc) 00075 { 00076 void *ctx; 00077 PetscErrorCode ierr = PCShellGetContext(pc,&ctx);CHKERRQ(ierr); 00078 Preconditioner<Number> * preconditioner = static_cast<Preconditioner<Number>*>(ctx); 00079 00080 if(!preconditioner->initialized()) 00081 { 00082 err<<"Preconditioner not initialized! Make sure you call init() before solve!"<<std::endl; 00083 libmesh_error(); 00084 } 00085 00086 preconditioner->setup(); 00087 00088 return 0; 00089 } 00090 00091 PetscErrorCode __libmesh_petsc_preconditioner_apply(PC pc, Vec x, Vec y) 00092 { 00093 void *ctx; 00094 PetscErrorCode ierr = PCShellGetContext(pc,&ctx);CHKERRQ(ierr); 00095 Preconditioner<Number> * preconditioner = static_cast<Preconditioner<Number>*>(ctx); 00096 00097 PetscVector<Number> x_vec(x); 00098 PetscVector<Number> y_vec(y); 00099 00100 preconditioner->apply(x_vec,y_vec); 00101 00102 return 0; 00103 } 00104 #endif 00105 } // end extern "C" 00106 00107 /*----------------------- functions ----------------------------------*/ 00108 template <typename T> 00109 void PetscLinearSolver<T>::clear () 00110 { 00111 if (this->initialized()) 00112 { 00113 /* If we were restricted to some subset, this restriction must 00114 be removed and the subset index set destroyed. */ 00115 if(_restrict_solve_to_is!=NULL) 00116 { 00117 int ierr = LibMeshISDestroy(&_restrict_solve_to_is); 00118 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00119 _restrict_solve_to_is = NULL; 00120 } 00121 00122 if(_restrict_solve_to_is_complement!=NULL) 00123 { 00124 int ierr = LibMeshISDestroy(&_restrict_solve_to_is_complement); 00125 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00126 _restrict_solve_to_is_complement = NULL; 00127 } 00128 00129 this->_is_initialized = false; 00130 00131 int ierr=0; 00132 00133 #if PETSC_VERSION_LESS_THAN(2,2,0) 00134 00135 // 2.1.x & earlier style 00136 ierr = SLESDestroy(_sles); 00137 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00138 00139 #else 00140 00141 // 2.2.0 & newer style 00142 ierr = LibMeshKSPDestroy(&_ksp); 00143 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00144 00145 #endif 00146 00147 00148 // Mimic PETSc default solver and preconditioner 00149 this->_solver_type = GMRES; 00150 00151 if(!this->_preconditioner) 00152 { 00153 if (libMesh::n_processors() == 1) 00154 this->_preconditioner_type = ILU_PRECOND; 00155 else 00156 this->_preconditioner_type = BLOCK_JACOBI_PRECOND; 00157 } 00158 } 00159 } 00160 00161 00162 00163 template <typename T> 00164 void PetscLinearSolver<T>::init () 00165 { 00166 // Initialize the data structures if not done so already. 00167 if (!this->initialized()) 00168 { 00169 this->_is_initialized = true; 00170 00171 int ierr=0; 00172 00173 // 2.1.x & earlier style 00174 #if PETSC_VERSION_LESS_THAN(2,2,0) 00175 00176 // Create the linear solver context 00177 ierr = SLESCreate (libMesh::COMM_WORLD, &_sles); 00178 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00179 00180 // Create the Krylov subspace & preconditioner contexts 00181 ierr = SLESGetKSP (_sles, &_ksp); 00182 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00183 ierr = SLESGetPC (_sles, &_pc); 00184 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00185 00186 // Set user-specified solver and preconditioner types 00187 this->set_petsc_solver_type(); 00188 00189 // Set the options from user-input 00190 // Set runtime options, e.g., 00191 // -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol> 00192 // These options will override those specified above as long as 00193 // SLESSetFromOptions() is called _after_ any other customization 00194 // routines. 00195 00196 ierr = SLESSetFromOptions (_sles); 00197 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00198 00199 // 2.2.0 & newer style 00200 #else 00201 00202 // Create the linear solver context 00203 ierr = KSPCreate (libMesh::COMM_WORLD, &_ksp); 00204 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00205 00206 // Create the preconditioner context 00207 ierr = KSPGetPC (_ksp, &_pc); 00208 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00209 00210 // Set user-specified solver and preconditioner types 00211 this->set_petsc_solver_type(); 00212 00213 // Set the options from user-input 00214 // Set runtime options, e.g., 00215 // -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol> 00216 // These options will override those specified above as long as 00217 // KSPSetFromOptions() is called _after_ any other customization 00218 // routines. 00219 00220 ierr = KSPSetFromOptions (_ksp); 00221 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00222 00223 // Not sure if this is necessary, or if it is already handled by KSPSetFromOptions? 00224 // NOT NECESSARY!!!! 00225 //ierr = PCSetFromOptions (_pc); 00226 //CHKERRABORT(libMesh::COMM_WORLD,ierr); 00227 00228 00229 #endif 00230 00231 // Have the Krylov subspace method use our good initial guess 00232 // rather than 0, unless the user requested a KSPType of 00233 // preonly, which complains if asked to use initial guesses. 00234 #if PETSC_VERSION_LESS_THAN(3,0,0) || !PETSC_VERSION_RELEASE 00235 // Pre-3.0 and petsc-dev (as of October 2012) use non-const versions 00236 KSPType ksp_type; 00237 #else 00238 const KSPType ksp_type; 00239 #endif 00240 00241 ierr = KSPGetType (_ksp, &ksp_type); 00242 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00243 00244 if (strcmp(ksp_type, "preonly")) 00245 { 00246 ierr = KSPSetInitialGuessNonzero (_ksp, PETSC_TRUE); 00247 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00248 } 00249 00250 // Notify PETSc of location to store residual history. 00251 // This needs to be called before any solves, since 00252 // it sets the residual history length to zero. The default 00253 // behavior is for PETSc to allocate (internally) an array 00254 // of size 1000 to hold the residual norm history. 00255 ierr = KSPSetResidualHistory(_ksp, 00256 PETSC_NULL, // pointer to the array which holds the history 00257 PETSC_DECIDE, // size of the array holding the history 00258 PETSC_TRUE); // Whether or not to reset the history for each solve. 00259 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00260 00261 PetscPreconditioner<T>::set_petsc_preconditioner_type(this->_preconditioner_type,_pc); 00262 00263 //If there is a preconditioner object we need to set the internal setup and apply routines 00264 if(this->_preconditioner) 00265 { 00266 this->_preconditioner->init(); 00267 PCShellSetContext(_pc,(void*)this->_preconditioner); 00268 PCShellSetSetUp(_pc,__libmesh_petsc_preconditioner_setup); 00269 PCShellSetApply(_pc,__libmesh_petsc_preconditioner_apply); 00270 } 00271 } 00272 } 00273 00274 00275 template <typename T> 00276 void PetscLinearSolver<T>::init ( PetscMatrix<T>* matrix ) 00277 { 00278 // Initialize the data structures if not done so already. 00279 if (!this->initialized()) 00280 { 00281 this->_is_initialized = true; 00282 00283 int ierr=0; 00284 00285 // 2.1.x & earlier style 00286 #if PETSC_VERSION_LESS_THAN(2,2,0) 00287 00288 // Create the linear solver context 00289 ierr = SLESCreate (libMesh::COMM_WORLD, &_sles); 00290 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00291 00292 // Create the Krylov subspace & preconditioner contexts 00293 ierr = SLESGetKSP (_sles, &_ksp); 00294 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00295 ierr = SLESGetPC (_sles, &_pc); 00296 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00297 00298 // Set user-specified solver and preconditioner types 00299 this->set_petsc_solver_type(); 00300 00301 // Set the options from user-input 00302 // Set runtime options, e.g., 00303 // -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol> 00304 // These options will override those specified above as long as 00305 // SLESSetFromOptions() is called _after_ any other customization 00306 // routines. 00307 00308 ierr = SLESSetFromOptions (_sles); 00309 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00310 00311 // 2.2.0 & newer style 00312 #else 00313 00314 // Create the linear solver context 00315 ierr = KSPCreate (libMesh::COMM_WORLD, &_ksp); 00316 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00317 00318 00319 //ierr = PCCreate (libMesh::COMM_WORLD, &_pc); 00320 // CHKERRABORT(libMesh::COMM_WORLD,ierr); 00321 00322 // Create the preconditioner context 00323 ierr = KSPGetPC (_ksp, &_pc); 00324 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00325 00326 // Set operators. The input matrix works as the preconditioning matrix 00327 ierr = KSPSetOperators(_ksp, matrix->mat(), matrix->mat(),DIFFERENT_NONZERO_PATTERN); 00328 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00329 00330 // Set user-specified solver and preconditioner types 00331 this->set_petsc_solver_type(); 00332 00333 // Set the options from user-input 00334 // Set runtime options, e.g., 00335 // -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol> 00336 // These options will override those specified above as long as 00337 // KSPSetFromOptions() is called _after_ any other customization 00338 // routines. 00339 00340 ierr = KSPSetFromOptions (_ksp); 00341 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00342 00343 // Not sure if this is necessary, or if it is already handled by KSPSetFromOptions? 00344 // NOT NECESSARY!!!! 00345 //ierr = PCSetFromOptions (_pc); 00346 //CHKERRABORT(libMesh::COMM_WORLD,ierr); 00347 00348 00349 #endif 00350 00351 // Have the Krylov subspace method use our good initial guess 00352 // rather than 0, unless the user requested a KSPType of 00353 // preonly, which complains if asked to use initial guesses. 00354 #if PETSC_VERSION_LESS_THAN(3,0,0) 00355 KSPType ksp_type; 00356 #else 00357 const KSPType ksp_type; 00358 #endif 00359 00360 ierr = KSPGetType (_ksp, &ksp_type); 00361 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00362 00363 if (strcmp(ksp_type, "preonly")) 00364 { 00365 ierr = KSPSetInitialGuessNonzero (_ksp, PETSC_TRUE); 00366 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00367 } 00368 00369 // Notify PETSc of location to store residual history. 00370 // This needs to be called before any solves, since 00371 // it sets the residual history length to zero. The default 00372 // behavior is for PETSc to allocate (internally) an array 00373 // of size 1000 to hold the residual norm history. 00374 ierr = KSPSetResidualHistory(_ksp, 00375 PETSC_NULL, // pointer to the array which holds the history 00376 PETSC_DECIDE, // size of the array holding the history 00377 PETSC_TRUE); // Whether or not to reset the history for each solve. 00378 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00379 00380 PetscPreconditioner<T>::set_petsc_preconditioner_type(this->_preconditioner_type,_pc); 00381 if(this->_preconditioner) 00382 { 00383 this->_preconditioner->set_matrix(*matrix); 00384 this->_preconditioner->init(); 00385 PCShellSetContext(_pc,(void*)this->_preconditioner); 00386 PCShellSetSetUp(_pc,__libmesh_petsc_preconditioner_setup); 00387 PCShellSetApply(_pc,__libmesh_petsc_preconditioner_apply); 00388 } 00389 } 00390 } 00391 00392 00393 00394 template <typename T> 00395 void 00396 PetscLinearSolver<T>::restrict_solve_to (const std::vector<unsigned int>* const dofs, 00397 const SubsetSolveMode subset_solve_mode) 00398 { 00399 int ierr=0; 00400 00401 /* The preconditioner (in particular if a default preconditioner) 00402 will have to be reset. We call this->clear() to do that. This 00403 call will also remove and free any previous subset that may have 00404 been set before. */ 00405 this->clear(); 00406 00407 _subset_solve_mode = subset_solve_mode; 00408 00409 if(dofs!=NULL) 00410 { 00411 PetscInt* petsc_dofs = NULL; 00412 ierr = PetscMalloc(dofs->size()*sizeof(PetscInt), &petsc_dofs); 00413 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00414 00415 for(size_t i=0; i<dofs->size(); i++) 00416 { 00417 petsc_dofs[i] = (*dofs)[i]; 00418 } 00419 00420 ierr = ISCreateLibMesh(libMesh::COMM_WORLD,dofs->size(),petsc_dofs,PETSC_OWN_POINTER,&_restrict_solve_to_is); 00421 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00422 } 00423 } 00424 00425 00426 00427 template <typename T> 00428 std::pair<unsigned int, Real> 00429 PetscLinearSolver<T>::solve (SparseMatrix<T>& matrix_in, 00430 SparseMatrix<T>& precond_in, 00431 NumericVector<T>& solution_in, 00432 NumericVector<T>& rhs_in, 00433 const double tol, 00434 const unsigned int m_its) 00435 { 00436 START_LOG("solve()", "PetscLinearSolver"); 00437 00438 // Make sure the data passed in are really of Petsc types 00439 PetscMatrix<T>* matrix = libmesh_cast_ptr<PetscMatrix<T>*>(&matrix_in); 00440 PetscMatrix<T>* precond = libmesh_cast_ptr<PetscMatrix<T>*>(&precond_in); 00441 PetscVector<T>* solution = libmesh_cast_ptr<PetscVector<T>*>(&solution_in); 00442 PetscVector<T>* rhs = libmesh_cast_ptr<PetscVector<T>*>(&rhs_in); 00443 00444 this->init (matrix); 00445 00446 int ierr=0; 00447 int its=0, max_its = static_cast<int>(m_its); 00448 PetscReal final_resid=0.; 00449 00450 // Close the matrices and vectors in case this wasn't already done. 00451 matrix->close (); 00452 precond->close (); 00453 solution->close (); 00454 rhs->close (); 00455 00456 // // If matrix != precond, then this means we have specified a 00457 // // special preconditioner, so reset preconditioner type to PCMAT. 00458 // if (matrix != precond) 00459 // { 00460 // this->_preconditioner_type = USER_PRECOND; 00461 // this->set_petsc_preconditioner_type (); 00462 // } 00463 00464 // 2.1.x & earlier style 00465 #if PETSC_VERSION_LESS_THAN(2,2,0) 00466 00467 if(_restrict_solve_to_is!=NULL) 00468 { 00469 libmesh_not_implemented(); 00470 } 00471 00472 // Set operators. The input matrix works as the preconditioning matrix 00473 ierr = SLESSetOperators(_sles, matrix->mat(), precond->mat(), 00474 DIFFERENT_NONZERO_PATTERN); 00475 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00476 00477 // Set the tolerances for the iterative solver. Use the user-supplied 00478 // tolerance for the relative residual & leave the others at default values. 00479 ierr = KSPSetTolerances (_ksp, tol, PETSC_DEFAULT, 00480 PETSC_DEFAULT, max_its); 00481 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00482 00483 00484 // Solve the linear system 00485 ierr = SLESSolve (_sles, rhs->vec(), solution->vec(), &its); 00486 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00487 00488 00489 // Get the norm of the final residual to return to the user. 00490 ierr = KSPGetResidualNorm (_ksp, &final_resid); 00491 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00492 00493 // 2.2.0 00494 #elif PETSC_VERSION_LESS_THAN(2,2,1) 00495 00496 if(_restrict_solve_to_is!=NULL) 00497 { 00498 libmesh_not_implemented(); 00499 } 00500 00501 // Set operators. The input matrix works as the preconditioning matrix 00502 //ierr = KSPSetOperators(_ksp, matrix->mat(), precond->mat(), 00503 // SAME_NONZERO_PATTERN); 00504 // CHKERRABORT(libMesh::COMM_WORLD,ierr); 00505 00506 00507 // Set the tolerances for the iterative solver. Use the user-supplied 00508 // tolerance for the relative residual & leave the others at default values. 00509 // Convergence is detected at iteration k if 00510 // ||r_k||_2 < max(rtol*||b||_2 , abstol) 00511 // where r_k is the residual vector and b is the right-hand side. Note that 00512 // it is the *maximum* of the two values, the larger of which will almost 00513 // always be rtol*||b||_2. 00514 ierr = KSPSetTolerances (_ksp, 00515 tol, // rtol = relative decrease in residual (1.e-5) 00516 PETSC_DEFAULT, // abstol = absolute convergence tolerance (1.e-50) 00517 PETSC_DEFAULT, // dtol = divergence tolerance (1.e+5) 00518 max_its); 00519 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00520 00521 00522 // Set the solution vector to use 00523 ierr = KSPSetSolution (_ksp, solution->vec()); 00524 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00525 00526 // Set the RHS vector to use 00527 ierr = KSPSetRhs (_ksp, rhs->vec()); 00528 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00529 00530 // Solve the linear system 00531 ierr = KSPSolve (_ksp); 00532 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00533 00534 // Get the number of iterations required for convergence 00535 ierr = KSPGetIterationNumber (_ksp, &its); 00536 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00537 00538 // Get the norm of the final residual to return to the user. 00539 ierr = KSPGetResidualNorm (_ksp, &final_resid); 00540 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00541 00542 // 2.2.1 & newer style 00543 #else 00544 00545 Mat submat = NULL; 00546 Mat subprecond = NULL; 00547 Vec subrhs = NULL; 00548 Vec subsolution = NULL; 00549 VecScatter scatter = NULL; 00550 PetscMatrix<Number>* subprecond_matrix = NULL; 00551 00552 // Set operators. Also restrict rhs and solution vector to 00553 // subdomain if neccessary. 00554 if(_restrict_solve_to_is!=NULL) 00555 { 00556 size_t is_local_size = this->_restrict_solve_to_is_local_size(); 00557 00558 ierr = VecCreate(libMesh::COMM_WORLD,&subrhs); 00559 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00560 ierr = VecSetSizes(subrhs,is_local_size,PETSC_DECIDE); 00561 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00562 ierr = VecSetFromOptions(subrhs); 00563 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00564 00565 ierr = VecCreate(libMesh::COMM_WORLD,&subsolution); 00566 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00567 ierr = VecSetSizes(subsolution,is_local_size,PETSC_DECIDE); 00568 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00569 ierr = VecSetFromOptions(subsolution); 00570 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00571 00572 ierr = VecScatterCreate(rhs->vec(),_restrict_solve_to_is, subrhs,NULL, &scatter); 00573 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00574 00575 ierr = VecScatterBegin(scatter,rhs->vec(),subrhs,INSERT_VALUES,SCATTER_FORWARD); 00576 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00577 ierr = VecScatterEnd(scatter,rhs->vec(),subrhs,INSERT_VALUES,SCATTER_FORWARD); 00578 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00579 00580 ierr = VecScatterBegin(scatter,solution->vec(),subsolution,INSERT_VALUES,SCATTER_FORWARD); 00581 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00582 ierr = VecScatterEnd(scatter,solution->vec(),subsolution,INSERT_VALUES,SCATTER_FORWARD); 00583 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00584 00585 #if PETSC_VERSION_LESS_THAN(3,1,0) 00586 ierr = MatGetSubMatrix(matrix->mat(), 00587 _restrict_solve_to_is,_restrict_solve_to_is, 00588 PETSC_DECIDE,MAT_INITIAL_MATRIX,&submat); 00589 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00590 ierr = MatGetSubMatrix(precond->mat(), 00591 _restrict_solve_to_is,_restrict_solve_to_is, 00592 PETSC_DECIDE,MAT_INITIAL_MATRIX,&subprecond); 00593 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00594 #else 00595 ierr = MatGetSubMatrix(matrix->mat(), 00596 _restrict_solve_to_is,_restrict_solve_to_is, 00597 MAT_INITIAL_MATRIX,&submat); 00598 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00599 ierr = MatGetSubMatrix(precond->mat(), 00600 _restrict_solve_to_is,_restrict_solve_to_is, 00601 MAT_INITIAL_MATRIX,&subprecond); 00602 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00603 #endif 00604 00605 /* Since removing columns of the matrix changes the equation 00606 system, we will now change the right hand side to compensate 00607 for this. Note that this is not necessary if \p SUBSET_ZERO 00608 has been selected. */ 00609 if(_subset_solve_mode!=SUBSET_ZERO) 00610 { 00611 _create_complement_is(rhs_in); 00612 size_t is_complement_local_size = rhs_in.local_size()-is_local_size; 00613 00614 Vec subvec1 = NULL; 00615 Mat submat1 = NULL; 00616 VecScatter scatter1 = NULL; 00617 00618 ierr = VecCreate(libMesh::COMM_WORLD,&subvec1); 00619 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00620 ierr = VecSetSizes(subvec1,is_complement_local_size,PETSC_DECIDE); 00621 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00622 ierr = VecSetFromOptions(subvec1); 00623 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00624 00625 ierr = VecScatterCreate(rhs->vec(),_restrict_solve_to_is_complement, subvec1,NULL, &scatter1); 00626 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00627 00628 ierr = VecScatterBegin(scatter1,_subset_solve_mode==SUBSET_COPY_RHS ? rhs->vec() : solution->vec(),subvec1,INSERT_VALUES,SCATTER_FORWARD); 00629 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00630 ierr = VecScatterEnd(scatter1,_subset_solve_mode==SUBSET_COPY_RHS ? rhs->vec() : solution->vec(),subvec1,INSERT_VALUES,SCATTER_FORWARD); 00631 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00632 00633 ierr = VecScale(subvec1,-1.0); 00634 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00635 00636 #if PETSC_VERSION_LESS_THAN(3,1,0) 00637 ierr = MatGetSubMatrix(matrix->mat(), 00638 _restrict_solve_to_is,_restrict_solve_to_is_complement, 00639 PETSC_DECIDE,MAT_INITIAL_MATRIX,&submat1); 00640 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00641 #else 00642 ierr = MatGetSubMatrix(matrix->mat(), 00643 _restrict_solve_to_is,_restrict_solve_to_is_complement, 00644 MAT_INITIAL_MATRIX,&submat1); 00645 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00646 #endif 00647 00648 ierr = MatMultAdd(submat1,subvec1,subrhs,subrhs); 00649 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00650 00651 ierr = LibMeshVecScatterDestroy(&scatter1); 00652 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00653 ierr = LibMeshVecDestroy(&subvec1); 00654 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00655 ierr = LibMeshMatDestroy(&submat1); 00656 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00657 } 00658 00659 ierr = KSPSetOperators(_ksp, submat, subprecond, 00660 this->same_preconditioner ? SAME_PRECONDITIONER : DIFFERENT_NONZERO_PATTERN); 00661 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00662 00663 if(this->_preconditioner) 00664 { 00665 subprecond_matrix = new PetscMatrix<Number>(subprecond); 00666 this->_preconditioner->set_matrix(*subprecond_matrix); 00667 this->_preconditioner->init(); 00668 } 00669 } 00670 else 00671 { 00672 ierr = KSPSetOperators(_ksp, matrix->mat(), precond->mat(), 00673 this->same_preconditioner ? SAME_PRECONDITIONER : DIFFERENT_NONZERO_PATTERN); 00674 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00675 00676 if(this->_preconditioner) 00677 { 00678 this->_preconditioner->set_matrix(matrix_in); 00679 this->_preconditioner->init(); 00680 } 00681 } 00682 00683 // Set the tolerances for the iterative solver. Use the user-supplied 00684 // tolerance for the relative residual & leave the others at default values. 00685 ierr = KSPSetTolerances (_ksp, tol, PETSC_DEFAULT, 00686 PETSC_DEFAULT, max_its); 00687 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00688 00689 // Solve the linear system 00690 if(_restrict_solve_to_is!=NULL) 00691 { 00692 ierr = KSPSolve (_ksp, subrhs, subsolution); 00693 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00694 } 00695 else 00696 { 00697 ierr = KSPSolve (_ksp, rhs->vec(), solution->vec()); 00698 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00699 } 00700 00701 // Get the number of iterations required for convergence 00702 ierr = KSPGetIterationNumber (_ksp, &its); 00703 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00704 00705 // Get the norm of the final residual to return to the user. 00706 ierr = KSPGetResidualNorm (_ksp, &final_resid); 00707 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00708 00709 if(_restrict_solve_to_is!=NULL) 00710 { 00711 switch(_subset_solve_mode) 00712 { 00713 case SUBSET_ZERO: 00714 ierr = VecZeroEntries(solution->vec()); 00715 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00716 break; 00717 00718 case SUBSET_COPY_RHS: 00719 ierr = VecCopy(rhs->vec(),solution->vec()); 00720 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00721 break; 00722 00723 case SUBSET_DONT_TOUCH: 00724 /* Nothing to do here. */ 00725 break; 00726 00727 } 00728 ierr = VecScatterBegin(scatter,subsolution,solution->vec(),INSERT_VALUES,SCATTER_REVERSE); 00729 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00730 ierr = VecScatterEnd(scatter,subsolution,solution->vec(),INSERT_VALUES,SCATTER_REVERSE); 00731 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00732 00733 ierr = LibMeshVecScatterDestroy(&scatter); 00734 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00735 00736 if(this->_preconditioner) 00737 { 00738 /* Before we delete subprecond_matrix, we should give the 00739 _preconditioner a different matrix. */ 00740 this->_preconditioner->set_matrix(matrix_in); 00741 this->_preconditioner->init(); 00742 delete subprecond_matrix; 00743 subprecond_matrix = NULL; 00744 } 00745 00746 ierr = LibMeshVecDestroy(&subsolution); 00747 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00748 ierr = LibMeshVecDestroy(&subrhs); 00749 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00750 ierr = LibMeshMatDestroy(&submat); 00751 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00752 ierr = LibMeshMatDestroy(&subprecond); 00753 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00754 } 00755 00756 #endif 00757 00758 STOP_LOG("solve()", "PetscLinearSolver"); 00759 // return the # of its. and the final residual norm. 00760 return std::make_pair(its, final_resid); 00761 } 00762 00763 template <typename T> 00764 std::pair<unsigned int, Real> 00765 PetscLinearSolver<T>::adjoint_solve (SparseMatrix<T>& matrix_in, 00766 NumericVector<T>& solution_in, 00767 NumericVector<T>& rhs_in, 00768 const double tol, 00769 const unsigned int m_its) 00770 { 00771 START_LOG("solve()", "PetscLinearSolver"); 00772 00773 // Make sure the data passed in are really of Petsc types 00774 PetscMatrix<T>* matrix = libmesh_cast_ptr<PetscMatrix<T>*>(&matrix_in); 00775 // Note that the matrix and precond matrix are the same 00776 PetscMatrix<T>* precond = libmesh_cast_ptr<PetscMatrix<T>*>(&matrix_in); 00777 PetscVector<T>* solution = libmesh_cast_ptr<PetscVector<T>*>(&solution_in); 00778 PetscVector<T>* rhs = libmesh_cast_ptr<PetscVector<T>*>(&rhs_in); 00779 00780 this->init (matrix); 00781 00782 int ierr=0; 00783 int its=0, max_its = static_cast<int>(m_its); 00784 PetscReal final_resid=0.; 00785 00786 // Close the matrices and vectors in case this wasn't already done. 00787 matrix->close (); 00788 precond->close (); 00789 solution->close (); 00790 rhs->close (); 00791 00792 // // If matrix != precond, then this means we have specified a 00793 // // special preconditioner, so reset preconditioner type to PCMAT. 00794 // if (matrix != precond) 00795 // { 00796 // this->_preconditioner_type = USER_PRECOND; 00797 // this->set_petsc_preconditioner_type (); 00798 // } 00799 00800 // 2.1.x & earlier style 00801 #if PETSC_VERSION_LESS_THAN(2,2,0) 00802 00803 if(_restrict_solve_to_is!=NULL) 00804 { 00805 libmesh_not_implemented(); 00806 } 00807 00808 // Based on http://wolfgang.math.tamu.edu/svn/public/deal.II/branches/MATH676/2008/deal.II/lac/source/petsc_solver.cc, http://tccc.iesl.forth.gr/AMS_EPEAEK/Elements/doc/in_html/petsc/SLES/index.html 00809 00810 SLES sles; 00811 ierr = SLESCreate (libMesh::COMM_WORLD, &sles); 00812 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00813 00814 ierr = SLESSetOperators (sles, matrix->mat(), precond->mat(), this->same_preconditioner ? SAME_PRECONDITIONER : DIFFERENT_NONZERO_PATTERN); 00815 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00816 00817 KSP ksp; 00818 ierr = SLESGetKSP (sles, &ksp); 00819 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00820 00821 ierr = SLESSetUp (sles, rhs->vec(), solution->vec()); 00822 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00823 00824 // See http://tccc.iesl.forth.gr/AMS_EPEAEK/Elements/doc/in_html/petsc/KSP/KSPSolveTrans.html#KSPSolveTrans 00825 ierr = SLESSolveTrans (ksp, &its); 00826 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00827 00828 // 2.2.0 00829 #elif PETSC_VERSION_LESS_THAN(2,2,1) 00830 00831 if(_restrict_solve_to_is!=NULL) 00832 { 00833 libmesh_not_implemented(); 00834 } 00835 00836 // Set operators. The input matrix works as the preconditioning matrix 00837 // This was commented earlier but it looks like KSPSetOperators is supported 00838 // after PETSc 2.2.0 00839 ierr = KSPSetOperators(_ksp, matrix->mat(), precond->mat(), 00840 this->same_preconditioner ? SAME_PRECONDITIONER : DIFFERENT_NONZERO_PATTERN); 00841 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00842 00843 00844 // Set the tolerances for the iterative solver. Use the user-supplied 00845 // tolerance for the relative residual & leave the others at default values. 00846 // Convergence is detected at iteration k if 00847 // ||r_k||_2 < max(rtol*||b||_2 , abstol) 00848 // where r_k is the residual vector and b is the right-hand side. Note that 00849 // it is the *maximum* of the two values, the larger of which will almost 00850 // always be rtol*||b||_2. 00851 ierr = KSPSetTolerances (_ksp, 00852 tol, // rtol = relative decrease in residual (1.e-5) 00853 PETSC_DEFAULT, // abstol = absolute convergence tolerance (1.e-50) 00854 PETSC_DEFAULT, // dtol = divergence tolerance (1.e+5) 00855 max_its); 00856 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00857 00858 00859 // Set the solution vector to use 00860 ierr = KSPSetSolution (_ksp, solution->vec()); 00861 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00862 00863 // Set the RHS vector to use 00864 ierr = KSPSetRhs (_ksp, rhs->vec()); 00865 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00866 00867 // Solve the linear system 00868 ierr = KSPSolveTranspose (_ksp); 00869 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00870 00871 // Get the number of iterations required for convergence 00872 ierr = KSPGetIterationNumber (_ksp, &its); 00873 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00874 00875 // Get the norm of the final residual to return to the user. 00876 ierr = KSPGetResidualNorm (_ksp, &final_resid); 00877 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00878 00879 // 2.2.1 & newer style 00880 #else 00881 00882 Mat submat = NULL; 00883 Mat subprecond = NULL; 00884 Vec subrhs = NULL; 00885 Vec subsolution = NULL; 00886 VecScatter scatter = NULL; 00887 PetscMatrix<Number>* subprecond_matrix = NULL; 00888 00889 // Set operators. Also restrict rhs and solution vector to 00890 // subdomain if neccessary. 00891 if(_restrict_solve_to_is!=NULL) 00892 { 00893 size_t is_local_size = this->_restrict_solve_to_is_local_size(); 00894 00895 ierr = VecCreate(libMesh::COMM_WORLD,&subrhs); 00896 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00897 ierr = VecSetSizes(subrhs,is_local_size,PETSC_DECIDE); 00898 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00899 ierr = VecSetFromOptions(subrhs); 00900 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00901 00902 ierr = VecCreate(libMesh::COMM_WORLD,&subsolution); 00903 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00904 ierr = VecSetSizes(subsolution,is_local_size,PETSC_DECIDE); 00905 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00906 ierr = VecSetFromOptions(subsolution); 00907 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00908 00909 ierr = VecScatterCreate(rhs->vec(),_restrict_solve_to_is, subrhs,NULL, &scatter); 00910 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00911 00912 ierr = VecScatterBegin(scatter,rhs->vec(),subrhs,INSERT_VALUES,SCATTER_FORWARD); 00913 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00914 ierr = VecScatterEnd(scatter,rhs->vec(),subrhs,INSERT_VALUES,SCATTER_FORWARD); 00915 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00916 00917 ierr = VecScatterBegin(scatter,solution->vec(),subsolution,INSERT_VALUES,SCATTER_FORWARD); 00918 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00919 ierr = VecScatterEnd(scatter,solution->vec(),subsolution,INSERT_VALUES,SCATTER_FORWARD); 00920 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00921 00922 #if PETSC_VERSION_LESS_THAN(3,1,0) 00923 ierr = MatGetSubMatrix(matrix->mat(), 00924 _restrict_solve_to_is,_restrict_solve_to_is, 00925 PETSC_DECIDE,MAT_INITIAL_MATRIX,&submat); 00926 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00927 ierr = MatGetSubMatrix(precond->mat(), 00928 _restrict_solve_to_is,_restrict_solve_to_is, 00929 PETSC_DECIDE,MAT_INITIAL_MATRIX,&subprecond); 00930 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00931 #else 00932 ierr = MatGetSubMatrix(matrix->mat(), 00933 _restrict_solve_to_is,_restrict_solve_to_is, 00934 MAT_INITIAL_MATRIX,&submat); 00935 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00936 ierr = MatGetSubMatrix(precond->mat(), 00937 _restrict_solve_to_is,_restrict_solve_to_is, 00938 MAT_INITIAL_MATRIX,&subprecond); 00939 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00940 #endif 00941 00942 /* Since removing columns of the matrix changes the equation 00943 system, we will now change the right hand side to compensate 00944 for this. Note that this is not necessary if \p SUBSET_ZERO 00945 has been selected. */ 00946 if(_subset_solve_mode!=SUBSET_ZERO) 00947 { 00948 _create_complement_is(rhs_in); 00949 size_t is_complement_local_size = rhs_in.local_size()-is_local_size; 00950 00951 Vec subvec1 = NULL; 00952 Mat submat1 = NULL; 00953 VecScatter scatter1 = NULL; 00954 00955 ierr = VecCreate(libMesh::COMM_WORLD,&subvec1); 00956 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00957 ierr = VecSetSizes(subvec1,is_complement_local_size,PETSC_DECIDE); 00958 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00959 ierr = VecSetFromOptions(subvec1); 00960 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00961 00962 ierr = VecScatterCreate(rhs->vec(),_restrict_solve_to_is_complement, subvec1,NULL, &scatter1); 00963 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00964 00965 ierr = VecScatterBegin(scatter1,_subset_solve_mode==SUBSET_COPY_RHS ? rhs->vec() : solution->vec(),subvec1,INSERT_VALUES,SCATTER_FORWARD); 00966 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00967 ierr = VecScatterEnd(scatter1,_subset_solve_mode==SUBSET_COPY_RHS ? rhs->vec() : solution->vec(),subvec1,INSERT_VALUES,SCATTER_FORWARD); 00968 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00969 00970 ierr = VecScale(subvec1,-1.0); 00971 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00972 00973 #if PETSC_VERSION_LESS_THAN(3,1,0) 00974 ierr = MatGetSubMatrix(matrix->mat(), 00975 _restrict_solve_to_is,_restrict_solve_to_is_complement, 00976 PETSC_DECIDE,MAT_INITIAL_MATRIX,&submat1); 00977 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00978 #else 00979 ierr = MatGetSubMatrix(matrix->mat(), 00980 _restrict_solve_to_is,_restrict_solve_to_is_complement, 00981 MAT_INITIAL_MATRIX,&submat1); 00982 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00983 #endif 00984 00985 ierr = MatMultAdd(submat1,subvec1,subrhs,subrhs); 00986 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00987 00988 ierr = LibMeshVecScatterDestroy(&scatter1); 00989 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00990 ierr = LibMeshVecDestroy(&subvec1); 00991 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00992 ierr = LibMeshMatDestroy(&submat1); 00993 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00994 } 00995 00996 ierr = KSPSetOperators(_ksp, submat, subprecond, 00997 this->same_preconditioner ? SAME_PRECONDITIONER : DIFFERENT_NONZERO_PATTERN); 00998 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00999 01000 if(this->_preconditioner) 01001 { 01002 subprecond_matrix = new PetscMatrix<Number>(subprecond); 01003 this->_preconditioner->set_matrix(*subprecond_matrix); 01004 this->_preconditioner->init(); 01005 } 01006 } 01007 else 01008 { 01009 ierr = KSPSetOperators(_ksp, matrix->mat(), precond->mat(), 01010 this->same_preconditioner ? SAME_PRECONDITIONER : DIFFERENT_NONZERO_PATTERN); 01011 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01012 01013 if(this->_preconditioner) 01014 { 01015 this->_preconditioner->set_matrix(matrix_in); 01016 this->_preconditioner->init(); 01017 } 01018 } 01019 01020 // Set the tolerances for the iterative solver. Use the user-supplied 01021 // tolerance for the relative residual & leave the others at default values. 01022 ierr = KSPSetTolerances (_ksp, tol, PETSC_DEFAULT, 01023 PETSC_DEFAULT, max_its); 01024 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01025 01026 // Solve the linear system 01027 if(_restrict_solve_to_is!=NULL) 01028 { 01029 ierr = KSPSolveTranspose (_ksp, subrhs, subsolution); 01030 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01031 } 01032 else 01033 { 01034 ierr = KSPSolveTranspose (_ksp, rhs->vec(), solution->vec()); 01035 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01036 } 01037 01038 // Get the number of iterations required for convergence 01039 ierr = KSPGetIterationNumber (_ksp, &its); 01040 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01041 01042 // Get the norm of the final residual to return to the user. 01043 ierr = KSPGetResidualNorm (_ksp, &final_resid); 01044 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01045 01046 if(_restrict_solve_to_is!=NULL) 01047 { 01048 switch(_subset_solve_mode) 01049 { 01050 case SUBSET_ZERO: 01051 ierr = VecZeroEntries(solution->vec()); 01052 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01053 break; 01054 01055 case SUBSET_COPY_RHS: 01056 ierr = VecCopy(rhs->vec(),solution->vec()); 01057 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01058 break; 01059 01060 case SUBSET_DONT_TOUCH: 01061 /* Nothing to do here. */ 01062 break; 01063 01064 } 01065 ierr = VecScatterBegin(scatter,subsolution,solution->vec(),INSERT_VALUES,SCATTER_REVERSE); 01066 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01067 ierr = VecScatterEnd(scatter,subsolution,solution->vec(),INSERT_VALUES,SCATTER_REVERSE); 01068 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01069 01070 ierr = LibMeshVecScatterDestroy(&scatter); 01071 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01072 01073 if(this->_preconditioner) 01074 { 01075 /* Before we delete subprecond_matrix, we should give the 01076 _preconditioner a different matrix. */ 01077 this->_preconditioner->set_matrix(matrix_in); 01078 this->_preconditioner->init(); 01079 delete subprecond_matrix; 01080 subprecond_matrix = NULL; 01081 } 01082 01083 ierr = LibMeshVecDestroy(&subsolution); 01084 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01085 ierr = LibMeshVecDestroy(&subrhs); 01086 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01087 ierr = LibMeshMatDestroy(&submat); 01088 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01089 ierr = LibMeshMatDestroy(&subprecond); 01090 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01091 } 01092 01093 #endif 01094 01095 STOP_LOG("solve()", "PetscLinearSolver"); 01096 // return the # of its. and the final residual norm. 01097 return std::make_pair(its, final_resid); 01098 } 01099 01100 01101 template <typename T> 01102 std::pair<unsigned int, Real> 01103 PetscLinearSolver<T>::solve (const ShellMatrix<T>& shell_matrix, 01104 NumericVector<T>& solution_in, 01105 NumericVector<T>& rhs_in, 01106 const double tol, 01107 const unsigned int m_its) 01108 { 01109 01110 #if PETSC_VERSION_LESS_THAN(2,3,1) 01111 // FIXME[JWP]: There will be a bunch of unused variable warnings 01112 // for older PETScs here. 01113 libMesh::out << "This method has been developed with PETSc 2.3.1. " 01114 << "No one has made it backwards compatible with older " 01115 << "versions of PETSc so far; however, it might work " 01116 << "without any change with some older version." << std::endl; 01117 libmesh_error(); 01118 return std::make_pair(0,0.0); 01119 01120 #else 01121 01122 #if PETSC_VERSION_LESS_THAN(3,1,0) 01123 if(_restrict_solve_to_is!=NULL) 01124 { 01125 libMesh::out << "The current implementation of subset solves with " 01126 << "shell matrices requires PETSc version 3.1 or above. " 01127 << "Older PETSc version do not support automatic " 01128 << "submatrix generation of shell matrices." 01129 << std::endl; 01130 libmesh_error(); 01131 } 01132 #endif 01133 01134 START_LOG("solve()", "PetscLinearSolver"); 01135 01136 // Make sure the data passed in are really of Petsc types 01137 PetscVector<T>* solution = libmesh_cast_ptr<PetscVector<T>*>(&solution_in); 01138 PetscVector<T>* rhs = libmesh_cast_ptr<PetscVector<T>*>(&rhs_in); 01139 01140 this->init (); 01141 01142 int ierr=0; 01143 int its=0, max_its = static_cast<int>(m_its); 01144 PetscReal final_resid=0.; 01145 01146 Mat submat = NULL; 01147 Vec subrhs = NULL; 01148 Vec subsolution = NULL; 01149 VecScatter scatter = NULL; 01150 01151 // Close the matrices and vectors in case this wasn't already done. 01152 solution->close (); 01153 rhs->close (); 01154 01155 // Prepare the matrix. 01156 Mat mat; 01157 ierr = MatCreateShell(libMesh::COMM_WORLD, 01158 rhs_in.local_size(), 01159 solution_in.local_size(), 01160 rhs_in.size(), 01161 solution_in.size(), 01162 const_cast<void*>(static_cast<const void*>(&shell_matrix)), 01163 &mat); 01164 /* Note that the const_cast above is only necessary because PETSc 01165 does not accept a const void*. Inside the member function 01166 _petsc_shell_matrix() below, the pointer is casted back to a 01167 const ShellMatrix<T>*. */ 01168 01169 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01170 ierr = MatShellSetOperation(mat,MATOP_MULT,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_mult)); 01171 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01172 ierr = MatShellSetOperation(mat,MATOP_MULT_ADD,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_mult_add)); 01173 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01174 ierr = MatShellSetOperation(mat,MATOP_GET_DIAGONAL,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_get_diagonal)); 01175 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01176 01177 // Restrict rhs and solution vectors and set operators. The input 01178 // matrix works as the preconditioning matrix. 01179 if(_restrict_solve_to_is!=NULL) 01180 { 01181 size_t is_local_size = this->_restrict_solve_to_is_local_size(); 01182 01183 ierr = VecCreate(libMesh::COMM_WORLD,&subrhs); 01184 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01185 ierr = VecSetSizes(subrhs,is_local_size,PETSC_DECIDE); 01186 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01187 ierr = VecSetFromOptions(subrhs); 01188 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01189 01190 ierr = VecCreate(libMesh::COMM_WORLD,&subsolution); 01191 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01192 ierr = VecSetSizes(subsolution,is_local_size,PETSC_DECIDE); 01193 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01194 ierr = VecSetFromOptions(subsolution); 01195 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01196 01197 ierr = VecScatterCreate(rhs->vec(),_restrict_solve_to_is, subrhs,NULL, &scatter); 01198 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01199 01200 ierr = VecScatterBegin(scatter,rhs->vec(),subrhs,INSERT_VALUES,SCATTER_FORWARD); 01201 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01202 ierr = VecScatterEnd(scatter,rhs->vec(),subrhs,INSERT_VALUES,SCATTER_FORWARD); 01203 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01204 01205 ierr = VecScatterBegin(scatter,solution->vec(),subsolution,INSERT_VALUES,SCATTER_FORWARD); 01206 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01207 ierr = VecScatterEnd(scatter,solution->vec(),subsolution,INSERT_VALUES,SCATTER_FORWARD); 01208 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01209 01210 #if PETSC_VERSION_LESS_THAN(3,1,0) 01211 /* This point can't be reached, see above. */ 01212 libmesh_assert(false); 01213 #else 01214 ierr = MatGetSubMatrix(mat, 01215 _restrict_solve_to_is,_restrict_solve_to_is, 01216 MAT_INITIAL_MATRIX,&submat); 01217 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01218 #endif 01219 01220 /* Since removing columns of the matrix changes the equation 01221 system, we will now change the right hand side to compensate 01222 for this. Note that this is not necessary if \p SUBSET_ZERO 01223 has been selected. */ 01224 if(_subset_solve_mode!=SUBSET_ZERO) 01225 { 01226 _create_complement_is(rhs_in); 01227 size_t is_complement_local_size = rhs_in.local_size()-is_local_size; 01228 01229 Vec subvec1 = NULL; 01230 Mat submat1 = NULL; 01231 VecScatter scatter1 = NULL; 01232 01233 ierr = VecCreate(libMesh::COMM_WORLD,&subvec1); 01234 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01235 ierr = VecSetSizes(subvec1,is_complement_local_size,PETSC_DECIDE); 01236 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01237 ierr = VecSetFromOptions(subvec1); 01238 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01239 01240 ierr = VecScatterCreate(rhs->vec(),_restrict_solve_to_is_complement, subvec1,NULL, &scatter1); 01241 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01242 01243 ierr = VecScatterBegin(scatter1,_subset_solve_mode==SUBSET_COPY_RHS ? rhs->vec() : solution->vec(),subvec1,INSERT_VALUES,SCATTER_FORWARD); 01244 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01245 ierr = VecScatterEnd(scatter1,_subset_solve_mode==SUBSET_COPY_RHS ? rhs->vec() : solution->vec(),subvec1,INSERT_VALUES,SCATTER_FORWARD); 01246 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01247 01248 ierr = VecScale(subvec1,-1.0); 01249 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01250 01251 #if PETSC_VERSION_LESS_THAN(3,1,0) 01252 /* This point can't be reached, see above. */ 01253 libmesh_assert(false); 01254 #else 01255 ierr = MatGetSubMatrix(mat, 01256 _restrict_solve_to_is,_restrict_solve_to_is_complement, 01257 MAT_INITIAL_MATRIX,&submat1); 01258 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01259 #endif 01260 01261 // The following lines would be correct, but don't work 01262 // correctly in PETSc up to 3.1.0-p5. See discussion in 01263 // petsc-users of Nov 9, 2010. 01264 // 01265 // ierr = MatMultAdd(submat1,subvec1,subrhs,subrhs); 01266 // CHKERRABORT(libMesh::COMM_WORLD,ierr); 01267 // 01268 // We workaround by using a temporary vector. Note that the 01269 // fix in PETsc 3.1.0-p6 uses a temporary vector internally, 01270 // so this is no effective performance loss. 01271 Vec subvec2 = NULL; 01272 ierr = VecCreate(libMesh::COMM_WORLD,&subvec2); 01273 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01274 ierr = VecSetSizes(subvec2,is_local_size,PETSC_DECIDE); 01275 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01276 ierr = VecSetFromOptions(subvec2); 01277 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01278 ierr = MatMult(submat1,subvec1,subvec2); 01279 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01280 ierr = VecAXPY(subrhs,1.0,subvec2); 01281 01282 ierr = LibMeshVecScatterDestroy(&scatter1); 01283 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01284 ierr = LibMeshVecDestroy(&subvec1); 01285 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01286 ierr = LibMeshMatDestroy(&submat1); 01287 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01288 } 01289 01290 ierr = KSPSetOperators(_ksp, submat, submat, 01291 DIFFERENT_NONZERO_PATTERN); 01292 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01293 } 01294 else 01295 { 01296 ierr = KSPSetOperators(_ksp, mat, mat, 01297 DIFFERENT_NONZERO_PATTERN); 01298 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01299 } 01300 01301 // Set the tolerances for the iterative solver. Use the user-supplied 01302 // tolerance for the relative residual & leave the others at default values. 01303 ierr = KSPSetTolerances (_ksp, tol, PETSC_DEFAULT, 01304 PETSC_DEFAULT, max_its); 01305 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01306 01307 // Solve the linear system 01308 if(_restrict_solve_to_is!=NULL) 01309 { 01310 ierr = KSPSolve (_ksp, subrhs, subsolution); 01311 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01312 } 01313 else 01314 { 01315 ierr = KSPSolve (_ksp, rhs->vec(), solution->vec()); 01316 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01317 } 01318 01319 // Get the number of iterations required for convergence 01320 ierr = KSPGetIterationNumber (_ksp, &its); 01321 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01322 01323 // Get the norm of the final residual to return to the user. 01324 ierr = KSPGetResidualNorm (_ksp, &final_resid); 01325 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01326 01327 if(_restrict_solve_to_is!=NULL) 01328 { 01329 switch(_subset_solve_mode) 01330 { 01331 case SUBSET_ZERO: 01332 ierr = VecZeroEntries(solution->vec()); 01333 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01334 break; 01335 01336 case SUBSET_COPY_RHS: 01337 ierr = VecCopy(rhs->vec(),solution->vec()); 01338 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01339 break; 01340 01341 case SUBSET_DONT_TOUCH: 01342 /* Nothing to do here. */ 01343 break; 01344 01345 } 01346 ierr = VecScatterBegin(scatter,subsolution,solution->vec(),INSERT_VALUES,SCATTER_REVERSE); 01347 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01348 ierr = VecScatterEnd(scatter,subsolution,solution->vec(),INSERT_VALUES,SCATTER_REVERSE); 01349 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01350 01351 ierr = LibMeshVecScatterDestroy(&scatter); 01352 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01353 01354 ierr = LibMeshVecDestroy(&subsolution); 01355 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01356 ierr = LibMeshVecDestroy(&subrhs); 01357 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01358 ierr = LibMeshMatDestroy(&submat); 01359 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01360 } 01361 01362 // Destroy the matrix. 01363 ierr = LibMeshMatDestroy(&mat); 01364 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01365 01366 STOP_LOG("solve()", "PetscLinearSolver"); 01367 // return the # of its. and the final residual norm. 01368 return std::make_pair(its, final_resid); 01369 01370 #endif 01371 01372 } 01373 01374 01375 01376 template <typename T> 01377 std::pair<unsigned int, Real> 01378 PetscLinearSolver<T>::solve (const ShellMatrix<T>& shell_matrix, 01379 const SparseMatrix<T>& precond_matrix, 01380 NumericVector<T> &solution_in, 01381 NumericVector<T> &rhs_in, 01382 const double tol, 01383 const unsigned int m_its) 01384 { 01385 01386 #if PETSC_VERSION_LESS_THAN(2,3,1) 01387 // FIXME[JWP]: There will be a bunch of unused variable warnings 01388 // for older PETScs here. 01389 libMesh::out << "This method has been developed with PETSc 2.3.1. " 01390 << "No one has made it backwards compatible with older " 01391 << "versions of PETSc so far; however, it might work " 01392 << "without any change with some older version." << std::endl; 01393 libmesh_error(); 01394 return std::make_pair(0,0.0); 01395 01396 #else 01397 01398 #if PETSC_VERSION_LESS_THAN(3,1,0) 01399 if(_restrict_solve_to_is!=NULL) 01400 { 01401 libMesh::out << "The current implementation of subset solves with " 01402 << "shell matrices requires PETSc version 3.1 or above. " 01403 << "Older PETSc version do not support automatic " 01404 << "submatrix generation of shell matrices." 01405 << std::endl; 01406 libmesh_error(); 01407 } 01408 #endif 01409 01410 START_LOG("solve()", "PetscLinearSolver"); 01411 01412 // Make sure the data passed in are really of Petsc types 01413 const PetscMatrix<T>* precond = libmesh_cast_ptr<const PetscMatrix<T>*>(&precond_matrix); 01414 PetscVector<T>* solution = libmesh_cast_ptr<PetscVector<T>*>(&solution_in); 01415 PetscVector<T>* rhs = libmesh_cast_ptr<PetscVector<T>*>(&rhs_in); 01416 01417 this->init (); 01418 01419 int ierr=0; 01420 int its=0, max_its = static_cast<int>(m_its); 01421 PetscReal final_resid=0.; 01422 01423 Mat submat = NULL; 01424 Mat subprecond = NULL; 01425 Vec subrhs = NULL; 01426 Vec subsolution = NULL; 01427 VecScatter scatter = NULL; 01428 PetscMatrix<Number>* subprecond_matrix = NULL; 01429 01430 // Close the matrices and vectors in case this wasn't already done. 01431 solution->close (); 01432 rhs->close (); 01433 01434 // Prepare the matrix. 01435 Mat mat; 01436 ierr = MatCreateShell(libMesh::COMM_WORLD, 01437 rhs_in.local_size(), 01438 solution_in.local_size(), 01439 rhs_in.size(), 01440 solution_in.size(), 01441 const_cast<void*>(static_cast<const void*>(&shell_matrix)), 01442 &mat); 01443 /* Note that the const_cast above is only necessary because PETSc 01444 does not accept a const void*. Inside the member function 01445 _petsc_shell_matrix() below, the pointer is casted back to a 01446 const ShellMatrix<T>*. */ 01447 01448 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01449 ierr = MatShellSetOperation(mat,MATOP_MULT,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_mult)); 01450 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01451 ierr = MatShellSetOperation(mat,MATOP_MULT_ADD,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_mult_add)); 01452 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01453 ierr = MatShellSetOperation(mat,MATOP_GET_DIAGONAL,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_get_diagonal)); 01454 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01455 01456 // Restrict rhs and solution vectors and set operators. The input 01457 // matrix works as the preconditioning matrix. 01458 if(_restrict_solve_to_is!=NULL) 01459 { 01460 size_t is_local_size = this->_restrict_solve_to_is_local_size(); 01461 01462 ierr = VecCreate(libMesh::COMM_WORLD,&subrhs); 01463 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01464 ierr = VecSetSizes(subrhs,is_local_size,PETSC_DECIDE); 01465 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01466 ierr = VecSetFromOptions(subrhs); 01467 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01468 01469 ierr = VecCreate(libMesh::COMM_WORLD,&subsolution); 01470 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01471 ierr = VecSetSizes(subsolution,is_local_size,PETSC_DECIDE); 01472 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01473 ierr = VecSetFromOptions(subsolution); 01474 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01475 01476 ierr = VecScatterCreate(rhs->vec(),_restrict_solve_to_is, subrhs,NULL, &scatter); 01477 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01478 01479 ierr = VecScatterBegin(scatter,rhs->vec(),subrhs,INSERT_VALUES,SCATTER_FORWARD); 01480 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01481 ierr = VecScatterEnd(scatter,rhs->vec(),subrhs,INSERT_VALUES,SCATTER_FORWARD); 01482 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01483 01484 ierr = VecScatterBegin(scatter,solution->vec(),subsolution,INSERT_VALUES,SCATTER_FORWARD); 01485 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01486 ierr = VecScatterEnd(scatter,solution->vec(),subsolution,INSERT_VALUES,SCATTER_FORWARD); 01487 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01488 01489 #if PETSC_VERSION_LESS_THAN(3,1,0) 01490 /* This point can't be reached, see above. */ 01491 libmesh_assert(false); 01492 #else 01493 ierr = MatGetSubMatrix(mat, 01494 _restrict_solve_to_is,_restrict_solve_to_is, 01495 MAT_INITIAL_MATRIX,&submat); 01496 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01497 ierr = MatGetSubMatrix(const_cast<PetscMatrix<T>*>(precond)->mat(), 01498 _restrict_solve_to_is,_restrict_solve_to_is, 01499 MAT_INITIAL_MATRIX,&subprecond); 01500 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01501 #endif 01502 01503 /* Since removing columns of the matrix changes the equation 01504 system, we will now change the right hand side to compensate 01505 for this. Note that this is not necessary if \p SUBSET_ZERO 01506 has been selected. */ 01507 if(_subset_solve_mode!=SUBSET_ZERO) 01508 { 01509 _create_complement_is(rhs_in); 01510 size_t is_complement_local_size = rhs_in.local_size()-is_local_size; 01511 01512 Vec subvec1 = NULL; 01513 Mat submat1 = NULL; 01514 VecScatter scatter1 = NULL; 01515 01516 ierr = VecCreate(libMesh::COMM_WORLD,&subvec1); 01517 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01518 ierr = VecSetSizes(subvec1,is_complement_local_size,PETSC_DECIDE); 01519 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01520 ierr = VecSetFromOptions(subvec1); 01521 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01522 01523 ierr = VecScatterCreate(rhs->vec(),_restrict_solve_to_is_complement, subvec1,NULL, &scatter1); 01524 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01525 01526 ierr = VecScatterBegin(scatter1,_subset_solve_mode==SUBSET_COPY_RHS ? rhs->vec() : solution->vec(),subvec1,INSERT_VALUES,SCATTER_FORWARD); 01527 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01528 ierr = VecScatterEnd(scatter1,_subset_solve_mode==SUBSET_COPY_RHS ? rhs->vec() : solution->vec(),subvec1,INSERT_VALUES,SCATTER_FORWARD); 01529 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01530 01531 ierr = VecScale(subvec1,-1.0); 01532 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01533 01534 #if PETSC_VERSION_LESS_THAN(3,1,0) 01535 /* This point can't be reached, see above. */ 01536 libmesh_assert(false); 01537 #else 01538 ierr = MatGetSubMatrix(mat, 01539 _restrict_solve_to_is,_restrict_solve_to_is_complement, 01540 MAT_INITIAL_MATRIX,&submat1); 01541 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01542 #endif 01543 01544 // The following lines would be correct, but don't work 01545 // correctly in PETSc up to 3.1.0-p5. See discussion in 01546 // petsc-users of Nov 9, 2010. 01547 // 01548 // ierr = MatMultAdd(submat1,subvec1,subrhs,subrhs); 01549 // CHKERRABORT(libMesh::COMM_WORLD,ierr); 01550 // 01551 // We workaround by using a temporary vector. Note that the 01552 // fix in PETsc 3.1.0-p6 uses a temporary vector internally, 01553 // so this is no effective performance loss. 01554 Vec subvec2 = NULL; 01555 ierr = VecCreate(libMesh::COMM_WORLD,&subvec2); 01556 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01557 ierr = VecSetSizes(subvec2,is_local_size,PETSC_DECIDE); 01558 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01559 ierr = VecSetFromOptions(subvec2); 01560 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01561 ierr = MatMult(submat1,subvec1,subvec2); 01562 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01563 ierr = VecAXPY(subrhs,1.0,subvec2); 01564 01565 ierr = LibMeshVecScatterDestroy(&scatter1); 01566 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01567 ierr = LibMeshVecDestroy(&subvec1); 01568 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01569 ierr = LibMeshMatDestroy(&submat1); 01570 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01571 } 01572 01573 ierr = KSPSetOperators(_ksp, submat, subprecond, 01574 DIFFERENT_NONZERO_PATTERN); 01575 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01576 01577 if(this->_preconditioner) 01578 { 01579 subprecond_matrix = new PetscMatrix<Number>(subprecond); 01580 this->_preconditioner->set_matrix(*subprecond_matrix); 01581 this->_preconditioner->init(); 01582 } 01583 } 01584 else 01585 { 01586 ierr = KSPSetOperators(_ksp, mat, const_cast<PetscMatrix<T>*>(precond)->mat(), 01587 DIFFERENT_NONZERO_PATTERN); 01588 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01589 01590 if(this->_preconditioner) 01591 { 01592 this->_preconditioner->set_matrix(const_cast<SparseMatrix<Number>&>(precond_matrix)); 01593 this->_preconditioner->init(); 01594 } 01595 } 01596 01597 // Set the tolerances for the iterative solver. Use the user-supplied 01598 // tolerance for the relative residual & leave the others at default values. 01599 ierr = KSPSetTolerances (_ksp, tol, PETSC_DEFAULT, 01600 PETSC_DEFAULT, max_its); 01601 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01602 01603 // Solve the linear system 01604 if(_restrict_solve_to_is!=NULL) 01605 { 01606 ierr = KSPSolve (_ksp, subrhs, subsolution); 01607 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01608 } 01609 else 01610 { 01611 ierr = KSPSolve (_ksp, rhs->vec(), solution->vec()); 01612 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01613 } 01614 01615 // Get the number of iterations required for convergence 01616 ierr = KSPGetIterationNumber (_ksp, &its); 01617 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01618 01619 // Get the norm of the final residual to return to the user. 01620 ierr = KSPGetResidualNorm (_ksp, &final_resid); 01621 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01622 01623 if(_restrict_solve_to_is!=NULL) 01624 { 01625 switch(_subset_solve_mode) 01626 { 01627 case SUBSET_ZERO: 01628 ierr = VecZeroEntries(solution->vec()); 01629 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01630 break; 01631 01632 case SUBSET_COPY_RHS: 01633 ierr = VecCopy(rhs->vec(),solution->vec()); 01634 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01635 break; 01636 01637 case SUBSET_DONT_TOUCH: 01638 /* Nothing to do here. */ 01639 break; 01640 01641 } 01642 ierr = VecScatterBegin(scatter,subsolution,solution->vec(),INSERT_VALUES,SCATTER_REVERSE); 01643 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01644 ierr = VecScatterEnd(scatter,subsolution,solution->vec(),INSERT_VALUES,SCATTER_REVERSE); 01645 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01646 01647 ierr = LibMeshVecScatterDestroy(&scatter); 01648 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01649 01650 if(this->_preconditioner) 01651 { 01652 /* Before we delete subprecond_matrix, we should give the 01653 _preconditioner a different matrix. */ 01654 this->_preconditioner->set_matrix(const_cast<SparseMatrix<Number>&>(precond_matrix)); 01655 this->_preconditioner->init(); 01656 delete subprecond_matrix; 01657 subprecond_matrix = NULL; 01658 } 01659 01660 ierr = LibMeshVecDestroy(&subsolution); 01661 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01662 ierr = LibMeshVecDestroy(&subrhs); 01663 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01664 ierr = LibMeshMatDestroy(&submat); 01665 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01666 ierr = LibMeshMatDestroy(&subprecond); 01667 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01668 } 01669 01670 // Destroy the matrix. 01671 ierr = LibMeshMatDestroy(&mat); 01672 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01673 01674 STOP_LOG("solve()", "PetscLinearSolver"); 01675 // return the # of its. and the final residual norm. 01676 return std::make_pair(its, final_resid); 01677 01678 #endif 01679 01680 } 01681 01682 01683 01684 template <typename T> 01685 void PetscLinearSolver<T>::get_residual_history(std::vector<double>& hist) 01686 { 01687 int ierr = 0; 01688 int its = 0; 01689 01690 // Fill the residual history vector with the residual norms 01691 // Note that GetResidualHistory() does not copy any values, it 01692 // simply sets the pointer p. Note that for some Krylov subspace 01693 // methods, the number of residuals returned in the history 01694 // vector may be different from what you are expecting. For 01695 // example, TFQMR returns two residual values per iteration step. 01696 PetscReal* p; 01697 ierr = KSPGetResidualHistory(_ksp, &p, &its); 01698 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01699 01700 // Check for early return 01701 if (its == 0) return; 01702 01703 // Create space to store the result 01704 hist.resize(its); 01705 01706 // Copy history into the vector provided by the user. 01707 for (int i=0; i<its; ++i) 01708 { 01709 hist[i] = *p; 01710 p++; 01711 } 01712 } 01713 01714 01715 01716 01717 template <typename T> 01718 Real PetscLinearSolver<T>::get_initial_residual() 01719 { 01720 int ierr = 0; 01721 int its = 0; 01722 01723 // Fill the residual history vector with the residual norms 01724 // Note that GetResidualHistory() does not copy any values, it 01725 // simply sets the pointer p. Note that for some Krylov subspace 01726 // methods, the number of residuals returned in the history 01727 // vector may be different from what you are expecting. For 01728 // example, TFQMR returns two residual values per iteration step. 01729 PetscReal* p; 01730 ierr = KSPGetResidualHistory(_ksp, &p, &its); 01731 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01732 01733 // Check no residual history 01734 if (its == 0) 01735 { 01736 libMesh::err << "No iterations have been performed, returning 0." << std::endl; 01737 return 0.; 01738 } 01739 01740 // Otherwise, return the value pointed to by p. 01741 return *p; 01742 } 01743 01744 01745 01746 01747 template <typename T> 01748 void PetscLinearSolver<T>::set_petsc_solver_type() 01749 { 01750 int ierr = 0; 01751 01752 switch (this->_solver_type) 01753 { 01754 01755 case CG: 01756 ierr = KSPSetType (_ksp, (char*) KSPCG); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 01757 01758 case CR: 01759 ierr = KSPSetType (_ksp, (char*) KSPCR); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 01760 01761 case CGS: 01762 ierr = KSPSetType (_ksp, (char*) KSPCGS); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 01763 01764 case BICG: 01765 ierr = KSPSetType (_ksp, (char*) KSPBICG); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 01766 01767 case TCQMR: 01768 ierr = KSPSetType (_ksp, (char*) KSPTCQMR); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 01769 01770 case TFQMR: 01771 ierr = KSPSetType (_ksp, (char*) KSPTFQMR); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 01772 01773 case LSQR: 01774 ierr = KSPSetType (_ksp, (char*) KSPLSQR); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 01775 01776 case BICGSTAB: 01777 ierr = KSPSetType (_ksp, (char*) KSPBCGS); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 01778 01779 case MINRES: 01780 ierr = KSPSetType (_ksp, (char*) KSPMINRES); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 01781 01782 case GMRES: 01783 ierr = KSPSetType (_ksp, (char*) KSPGMRES); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 01784 01785 case RICHARDSON: 01786 ierr = KSPSetType (_ksp, (char*) KSPRICHARDSON); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 01787 01788 case CHEBYSHEV: 01789 #if defined(LIBMESH_HAVE_PETSC) && PETSC_VERSION_LESS_THAN(3,3,0) 01790 ierr = KSPSetType (_ksp, (char*) KSPCHEBYCHEV); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 01791 #else 01792 ierr = KSPSetType (_ksp, (char*) KSPCHEBYSHEV); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 01793 #endif 01794 01795 01796 default: 01797 libMesh::err << "ERROR: Unsupported PETSC Solver: " 01798 << this->_solver_type << std::endl 01799 << "Continuing with PETSC defaults" << std::endl; 01800 } 01801 } 01802 01803 template <typename T> 01804 void PetscLinearSolver<T>::print_converged_reason() 01805 { 01806 KSPConvergedReason reason; 01807 KSPGetConvergedReason(_ksp, &reason); 01808 libMesh::out << "Linear solver convergence/divergence reason: " << KSPConvergedReasons[reason] << std::endl; 01809 } 01810 01811 01812 01813 template <typename T> 01814 PetscErrorCode PetscLinearSolver<T>::_petsc_shell_matrix_mult(Mat mat, Vec arg, Vec dest) 01815 { 01816 /* Get the matrix context. */ 01817 int ierr=0; 01818 void* ctx; 01819 ierr = MatShellGetContext(mat,&ctx); 01820 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01821 01822 /* Get user shell matrix object. */ 01823 const ShellMatrix<T>& shell_matrix = *static_cast<const ShellMatrix<T>*>(ctx); 01824 01825 /* Make \p NumericVector instances around the vectors. */ 01826 PetscVector<T> arg_global(arg); 01827 PetscVector<T> dest_global(dest); 01828 01829 /* Call the user function. */ 01830 shell_matrix.vector_mult(dest_global,arg_global); 01831 01832 return ierr; 01833 } 01834 01835 01836 01837 template <typename T> 01838 PetscErrorCode PetscLinearSolver<T>::_petsc_shell_matrix_mult_add(Mat mat, Vec arg, Vec add, Vec dest) 01839 { 01840 /* Get the matrix context. */ 01841 int ierr=0; 01842 void* ctx; 01843 ierr = MatShellGetContext(mat,&ctx); 01844 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01845 01846 /* Get user shell matrix object. */ 01847 const ShellMatrix<T>& shell_matrix = *static_cast<const ShellMatrix<T>*>(ctx); 01848 01849 /* Make \p NumericVector instances around the vectors. */ 01850 PetscVector<T> arg_global(arg); 01851 PetscVector<T> dest_global(dest); 01852 PetscVector<T> add_global(add); 01853 01854 if(add!=arg) 01855 { 01856 arg_global = add_global; 01857 } 01858 01859 /* Call the user function. */ 01860 shell_matrix.vector_mult_add(dest_global,arg_global); 01861 01862 return ierr; 01863 } 01864 01865 01866 01867 template <typename T> 01868 PetscErrorCode PetscLinearSolver<T>::_petsc_shell_matrix_get_diagonal(Mat mat, Vec dest) 01869 { 01870 /* Get the matrix context. */ 01871 int ierr=0; 01872 void* ctx; 01873 ierr = MatShellGetContext(mat,&ctx); 01874 CHKERRABORT(libMesh::COMM_WORLD,ierr); 01875 01876 /* Get user shell matrix object. */ 01877 const ShellMatrix<T>& shell_matrix = *static_cast<const ShellMatrix<T>*>(ctx); 01878 01879 /* Make \p NumericVector instances around the vector. */ 01880 PetscVector<T> dest_global(dest); 01881 01882 /* Call the user function. */ 01883 shell_matrix.get_diagonal(dest_global); 01884 01885 return ierr; 01886 } 01887 01888 01889 01890 //------------------------------------------------------------------ 01891 // Explicit instantiations 01892 template class PetscLinearSolver<Number>; 01893 01894 } // namespace libMesh 01895 01896 01897 01898 #endif // #ifdef LIBMESH_HAVE_PETSC
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:48 UTC
Hosted By: