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