void SlepcEigenSolver<T>:: set_slepc_problem_type() { PetscErrorCode ierr = 0; switch (this->_eigen_problem_type) { case NHEP: ierr = EPSSetProblemType (_eps, EPS_NHEP); LIBMESH_CHKERR(ierr); return; case GNHEP: ierr = EPSSetProblemType (_eps, EPS_GNHEP); LIBMESH_CHKERR(ierr); return; case HEP: ierr = EPSSetProblemType (_eps, EPS_HEP); LIBMESH_CHKERR(ierr); return; case GHEP: ierr = EPSSetProblemType (_eps, EPS_GHEP); LIBMESH_CHKERR(ierr); return; #if !SLEPC_VERSION_LESS_THAN(3,3,0) // EPS_GHIEP added in 3.3.0 case GHIEP: ierr = EPSSetProblemType (_eps, EPS_GHIEP); LIBMESH_CHKERR(ierr); return; #endif default: libMesh::err << "ERROR: Unsupported SLEPc Eigen Problem: " << this->_eigen_problem_type << std::endl << "Continuing with SLEPc defaults" << std::endl; } }
std::pair<unsigned int, unsigned int> SlepcEigenSolver<T>::solve_standard (ShellMatrix<T> & shell_matrix, int nev, // number of requested eigenpairs int ncv, // number of basis vectors const double tol, // solver tolerance const unsigned int m_its) // maximum number of iterations { this->init (); PetscErrorCode ierr=0; // Prepare the matrix. Note that the const_cast is only necessary // because PETSc does not accept a const void *. Inside the member // function _petsc_shell_matrix() below, the pointer is casted back // to a const ShellMatrix<T> *. Mat mat; ierr = MatCreateShell(this->comm().get(), shell_matrix.m(), // Specify the number of local rows shell_matrix.n(), // Specify the number of local columns PETSC_DETERMINE, PETSC_DETERMINE, const_cast<void *>(static_cast<const void *>(&shell_matrix)), &mat); LIBMESH_CHKERR(ierr); ierr = MatShellSetOperation(mat,MATOP_MULT,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_mult)); LIBMESH_CHKERR(ierr); ierr = MatShellSetOperation(mat,MATOP_GET_DIAGONAL,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_get_diagonal)); LIBMESH_CHKERR(ierr); return _solve_standard_helper(mat, nev, ncv, tol, m_its); }
std::pair<Real, Real> SlepcEigenSolver<T>::get_eigenpair(dof_id_type i, NumericVector<T> & solution_in) { PetscErrorCode ierr=0; PetscReal re, im; // Make sure the NumericVector passed in is really a PetscVector PetscVector<T> * solution = dynamic_cast<PetscVector<T> *>(&solution_in); if (!solution) libmesh_error_msg("Error getting eigenvector: input vector must be a PetscVector."); // real and imaginary part of the ith eigenvalue. PetscScalar kr, ki; solution->close(); ierr = EPSGetEigenpair(_eps, i, &kr, &ki, solution->vec(), PETSC_NULL); LIBMESH_CHKERR(ierr); #ifdef LIBMESH_USE_COMPLEX_NUMBERS re = PetscRealPart(kr); im = PetscImaginaryPart(kr); #else re = kr; im = ki; #endif return std::make_pair(re, im); }
//! Function to give PETSc that sets the Restriction Matrix between two DMs PetscErrorCode libmesh_petsc_DMCreateRestriction (DM dmc /*coarse*/, DM dmf/*fine*/, Mat * mat) { libmesh_assert(dmc); libmesh_assert(dmf); libmesh_assert(mat); PetscErrorCode ierr; // get a communicator from incomming DM MPI_Comm comm; PetscObjectGetComm((PetscObject)dmc, &comm); // extract our fine context from the incoming DM void * ctx_f = nullptr; ierr = DMShellGetContext(dmf, &ctx_f);LIBMESH_CHKERR(ierr); libmesh_assert(ctx_f); PetscDMContext * p_ctx_f = static_cast<PetscDMContext*>(ctx_f); // check / give PETSc its matrix libmesh_assert(p_ctx_f->K_restrict_ptr); *(mat) = p_ctx_f->K_restrict_ptr->mat(); return 0; }
std::pair<unsigned int, unsigned int> SlepcEigenSolver<T>::solve_generalized (SparseMatrix<T> & matrix_A_in, ShellMatrix<T> & shell_matrix_B, int nev, // number of requested eigenpairs int ncv, // number of basis vectors const double tol, // solver tolerance const unsigned int m_its) // maximum number of iterations { this->init (); PetscErrorCode ierr=0; PetscMatrix<T> * matrix_A = dynamic_cast<PetscMatrix<T> *>(&matrix_A_in); if (!matrix_A) libmesh_error_msg("Error: inputs to solve_generalized() must be of type PetscMatrix."); // Close the matrix and vectors in case this wasn't already done. if (this->_close_matrix_before_solve) matrix_A->close (); // Prepare the matrix. Note that the const_cast is only necessary // because PETSc does not accept a const void *. Inside the member // function _petsc_shell_matrix() below, the pointer is casted back // to a const ShellMatrix<T> *. Mat mat_B; ierr = MatCreateShell(this->comm().get(), shell_matrix_B.m(), // Specify the number of local rows shell_matrix_B.n(), // Specify the number of local columns PETSC_DETERMINE, PETSC_DETERMINE, const_cast<void *>(static_cast<const void *>(&shell_matrix_B)), &mat_B); LIBMESH_CHKERR(ierr); ierr = MatShellSetOperation(mat_B,MATOP_MULT,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_mult)); LIBMESH_CHKERR(ierr); ierr = MatShellSetOperation(mat_B,MATOP_GET_DIAGONAL,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_get_diagonal)); LIBMESH_CHKERR(ierr); return _solve_generalized_helper (matrix_A->mat(), mat_B, nev, ncv, tol, m_its); }
int TaoOptimizationSolver<T>::get_converged_reason() { PetscErrorCode ierr=0; if (this->initialized()) { ierr = TaoGetConvergedReason(_tao, &_reason); LIBMESH_CHKERR(ierr); } return static_cast<int>(_reason); }
void TaoOptimizationSolver<T>::clear () { if (this->initialized()) { this->_is_initialized = false; PetscErrorCode ierr=0; ierr = TaoDestroy(&_tao); LIBMESH_CHKERR(ierr); } }
void SlepcEigenSolver<T>::set_slepc_solver_type() { PetscErrorCode ierr = 0; switch (this->_eigen_solver_type) { case POWER: ierr = EPSSetType (_eps, EPSPOWER); LIBMESH_CHKERR(ierr); return; case SUBSPACE: ierr = EPSSetType (_eps, EPSSUBSPACE); LIBMESH_CHKERR(ierr); return; case LAPACK: ierr = EPSSetType (_eps, EPSLAPACK); LIBMESH_CHKERR(ierr); return; case ARNOLDI: ierr = EPSSetType (_eps, EPSARNOLDI); LIBMESH_CHKERR(ierr); return; case LANCZOS: ierr = EPSSetType (_eps, EPSLANCZOS); LIBMESH_CHKERR(ierr); return; case KRYLOVSCHUR: ierr = EPSSetType (_eps, EPSKRYLOVSCHUR); LIBMESH_CHKERR(ierr); return; // case ARPACK: // ierr = EPSSetType (_eps, (char *) EPSARPACK); LIBMESH_CHKERR(ierr); return; default: libMesh::err << "ERROR: Unsupported SLEPc Eigen Solver: " << Utility::enum_to_string(this->_eigen_solver_type) << std::endl << "Continuing with SLEPc defaults" << std::endl; } }
void TaoOptimizationSolver<T>::init () { // Initialize the data structures if not done so already. if (!this->initialized()) { this->_is_initialized = true; PetscErrorCode ierr=0; ierr = TaoCreate(this->comm().get(),&_tao); LIBMESH_CHKERR(ierr); } }
Real SlepcEigenSolver<T>::get_relative_error(unsigned int i) { PetscErrorCode ierr=0; PetscReal error; #if SLEPC_VERSION_LESS_THAN(3,6,0) ierr = EPSComputeRelativeError(_eps, i, &error); #else ierr = EPSComputeError(_eps, i, EPS_ERROR_RELATIVE, &error); #endif LIBMESH_CHKERR(ierr); return error; }
void SlepcEigenSolver<T>::clear () { if (this->initialized()) { this->_is_initialized = false; PetscErrorCode ierr=0; ierr = LibMeshEPSDestroy(&_eps); LIBMESH_CHKERR(ierr); // SLEPc default eigenproblem solver this->_eigen_solver_type = KRYLOVSCHUR; } }
void TaoOptimizationSolver<T>::get_dual_variables() { LOG_SCOPE("get_dual_variables()", "TaoOptimizationSolver"); PetscVector<T> * lambda_eq_petsc = cast_ptr<PetscVector<T> *>(this->system().lambda_eq.get()); PetscVector<T> * lambda_ineq_petsc = cast_ptr<PetscVector<T> *>(this->system().lambda_ineq.get()); Vec lambda_eq_petsc_vec = lambda_eq_petsc->vec(); Vec lambda_ineq_petsc_vec = lambda_ineq_petsc->vec(); PetscErrorCode ierr = 0; ierr = TaoGetDualVariables(_tao, &lambda_eq_petsc_vec, &lambda_ineq_petsc_vec); LIBMESH_CHKERR(ierr); }
void SlepcEigenSolver<T>::init () { PetscErrorCode ierr=0; // Initialize the data structures if not done so already. if (!this->initialized()) { this->_is_initialized = true; // Create the eigenproblem solver context ierr = EPSCreate (this->comm().get(), &_eps); LIBMESH_CHKERR(ierr); // Set user-specified solver set_slepc_solver_type(); } }
//! Help PETSc identify the finer DM given a dmc PetscErrorCode libmesh_petsc_DMRefine(DM dmc, MPI_Comm /*comm*/, DM * dmf) { libmesh_assert(dmc); libmesh_assert(dmf); PetscErrorCode ierr; // extract our context from the incoming dmc void * ctx_c = nullptr; ierr = DMShellGetContext(dmc, & ctx_c);LIBMESH_CHKERR(ierr); libmesh_assert(ctx_c); PetscDMContext * p_ctx = static_cast<PetscDMContext * >(ctx_c); // check / set the finer DM libmesh_assert(p_ctx->finer_dm); libmesh_assert(*(p_ctx->finer_dm)); *(dmf) = *(p_ctx->finer_dm); return 0; }
std::pair<Real, Real> SlepcEigenSolver<T>::get_eigenvalue(dof_id_type i) { PetscErrorCode ierr=0; PetscReal re, im; // real and imaginary part of the ith eigenvalue. PetscScalar kr, ki; ierr = EPSGetEigenvalue(_eps, i, &kr, &ki); LIBMESH_CHKERR(ierr); #ifdef LIBMESH_USE_COMPLEX_NUMBERS re = PetscRealPart(kr); im = PetscImaginaryPart(kr); #else re = kr; im = ki; #endif return std::make_pair(re, im); }
void SlepcEigenSolver<T>::attach_deflation_space(NumericVector<T> & deflation_vector_in) { this->init(); PetscErrorCode ierr = 0; // Make sure the input vector is actually a PetscVector PetscVector<T> * deflation_vector_petsc_vec = dynamic_cast<PetscVector<T> *>(&deflation_vector_in); if (!deflation_vector_petsc_vec) libmesh_error_msg("Error attaching deflation space: input vector must be a PetscVector."); // Get a handle for the underlying Vec. Vec deflation_vector = deflation_vector_petsc_vec->vec(); #if SLEPC_VERSION_LESS_THAN(3,1,0) ierr = EPSAttachDeflationSpace(_eps, 1, &deflation_vector, PETSC_FALSE); #else ierr = EPSSetDeflationSpace(_eps, 1, &deflation_vector); #endif LIBMESH_CHKERR(ierr); }
void SlepcEigenSolver<T>::set_initial_space(NumericVector<T> & initial_space_in) { #if SLEPC_VERSION_LESS_THAN(3,1,0) libmesh_error_msg("SLEPc 3.1 is required to call EigenSolver::set_initial_space()"); #else this->init(); PetscErrorCode ierr = 0; // Make sure the input vector is actually a PetscVector PetscVector<T> * initial_space_petsc_vec = dynamic_cast<PetscVector<T> *>(&initial_space_in); if (!initial_space_petsc_vec) libmesh_error_msg("Error attaching initial space: input vector must be a PetscVector."); // Get a handle for the underlying Vec. Vec initial_vector = initial_space_petsc_vec->vec(); ierr = EPSSetInitialSpace(_eps, 1, &initial_vector); LIBMESH_CHKERR(ierr); #endif }
void TaoOptimizationSolver<T>::solve () { LOG_SCOPE("solve()", "TaoOptimizationSolver"); this->init (); this->system().solution->zero(); PetscMatrix<T> * hessian = cast_ptr<PetscMatrix<T> *>(this->system().matrix); // PetscVector<T> * gradient = cast_ptr<PetscVector<T> *>(this->system().rhs); PetscVector<T> * x = cast_ptr<PetscVector<T> *>(this->system().solution.get()); PetscVector<T> * ceq = cast_ptr<PetscVector<T> *>(this->system().C_eq.get()); PetscMatrix<T> * ceq_jac = cast_ptr<PetscMatrix<T> *>(this->system().C_eq_jac.get()); PetscVector<T> * cineq = cast_ptr<PetscVector<T> *>(this->system().C_ineq.get()); PetscMatrix<T> * cineq_jac = cast_ptr<PetscMatrix<T> *>(this->system().C_ineq_jac.get()); PetscVector<T> * lb = cast_ptr<PetscVector<T> *>(&this->system().get_vector("lower_bounds")); PetscVector<T> * ub = cast_ptr<PetscVector<T> *>(&this->system().get_vector("upper_bounds")); // Set the starting guess to zero. x->zero(); PetscErrorCode ierr = 0; // Workaround for bug where TaoSetFromOptions *reset* // programmatically set tolerance and max. function evaluation // values when "-tao_type ipm" was specified on the command line: we // call TaoSetFromOptions twice (both before and after setting // custom options programatically) ierr = TaoSetFromOptions(_tao); LIBMESH_CHKERR(ierr); // Set convergence tolerances // f(X) - f(X*) (estimated) <= fatol // |f(X) - f(X*)| (estimated) / |f(X)| <= frtol // ||g(X)|| <= gatol // ||g(X)|| / |f(X)| <= grtol // ||g(X)|| / ||g(X0)|| <= gttol // Command line equivalents: -tao_fatol, -tao_frtol, -tao_gatol, -tao_grtol, -tao_gttol ierr = TaoSetTolerances(_tao, #if PETSC_RELEASE_LESS_THAN(3,7,0) // Releases up to 3.X.Y had fatol and frtol, after that they were removed. // Hopefully we'll be able to know X and Y soon. Guessing at 3.7.0. /*fatol=*/PETSC_DEFAULT, /*frtol=*/PETSC_DEFAULT, #endif /*gatol=*/PETSC_DEFAULT, /*grtol=*/this->objective_function_relative_tolerance, /*gttol=*/PETSC_DEFAULT); LIBMESH_CHKERR(ierr); // Set the max-allowed number of objective function evaluations // Command line equivalent: -tao_max_funcs ierr = TaoSetMaximumFunctionEvaluations(_tao, this->max_objective_function_evaluations); LIBMESH_CHKERR(ierr); // Set the max-allowed number of optimization iterations. // Command line equivalent: -tao_max_it // Not implemented for now as it seems fairly similar to // ierr = TaoSetMaximumIterations(_tao, 4); // LIBMESH_CHKERR(ierr); // Set solution vec and an initial guess ierr = TaoSetInitialVector(_tao, x->vec()); LIBMESH_CHKERR(ierr); // We have to have an objective function libmesh_assert( this->objective_object ); // Set routines for objective, gradient, hessian evaluation ierr = TaoSetObjectiveRoutine(_tao, __libmesh_tao_objective, this); LIBMESH_CHKERR(ierr); if ( this->gradient_object ) { ierr = TaoSetGradientRoutine(_tao, __libmesh_tao_gradient, this); LIBMESH_CHKERR(ierr); } if ( this->hessian_object ) { ierr = TaoSetHessianRoutine(_tao, hessian->mat(), hessian->mat(), __libmesh_tao_hessian, this); LIBMESH_CHKERR(ierr); } if ( this->lower_and_upper_bounds_object ) { // Need to actually compute the bounds vectors first this->lower_and_upper_bounds_object->lower_and_upper_bounds(this->system()); ierr = TaoSetVariableBounds(_tao, lb->vec(), ub->vec()); LIBMESH_CHKERR(ierr); } if ( this->equality_constraints_object ) { ierr = TaoSetEqualityConstraintsRoutine(_tao, ceq->vec(), __libmesh_tao_equality_constraints, this); LIBMESH_CHKERR(ierr); } if ( this->equality_constraints_jacobian_object ) { ierr = TaoSetJacobianEqualityRoutine(_tao, ceq_jac->mat(), ceq_jac->mat(), __libmesh_tao_equality_constraints_jacobian, this); LIBMESH_CHKERR(ierr); } // Optionally set inequality constraints if ( this->inequality_constraints_object ) { ierr = TaoSetInequalityConstraintsRoutine(_tao, cineq->vec(), __libmesh_tao_inequality_constraints, this); LIBMESH_CHKERR(ierr); } // Optionally set inequality constraints Jacobian if ( this->inequality_constraints_jacobian_object ) { ierr = TaoSetJacobianInequalityRoutine(_tao, cineq_jac->mat(), cineq_jac->mat(), __libmesh_tao_inequality_constraints_jacobian, this); LIBMESH_CHKERR(ierr); } // Check for Tao command line options ierr = TaoSetFromOptions(_tao); LIBMESH_CHKERR(ierr); // Perform the optimization ierr = TaoSolve(_tao); LIBMESH_CHKERR(ierr); // Store the convergence/divergence reason ierr = TaoGetConvergedReason(_tao, &_reason); LIBMESH_CHKERR(ierr); }
std::pair<unsigned int, unsigned int> SlepcEigenSolver<T>::_solve_generalized_helper (Mat mat_A, Mat mat_B, int nev, // number of requested eigenpairs int ncv, // number of basis vectors const double tol, // solver tolerance const unsigned int m_its) // maximum number of iterations { LOG_SCOPE("solve_generalized()", "SlepcEigenSolver"); PetscErrorCode ierr=0; // converged eigen pairs and number of iterations PetscInt nconv=0; PetscInt its=0; #ifdef DEBUG // The relative error. PetscReal error, re, im; // Pointer to vectors of the real parts, imaginary parts. PetscScalar kr, ki; #endif // Set operators. ierr = EPSSetOperators (_eps, mat_A, mat_B); LIBMESH_CHKERR(ierr); //set the problem type and the position of the spectrum set_slepc_problem_type(); set_slepc_position_of_spectrum(); // Set eigenvalues to be computed. #if SLEPC_VERSION_LESS_THAN(3,0,0) ierr = EPSSetDimensions (_eps, nev, ncv); #else ierr = EPSSetDimensions (_eps, nev, ncv, PETSC_DECIDE); #endif LIBMESH_CHKERR(ierr); // Set the tolerance and maximum iterations. ierr = EPSSetTolerances (_eps, tol, m_its); LIBMESH_CHKERR(ierr); // Set runtime options, e.g., // -eps_type <type>, -eps_nev <nev>, -eps_ncv <ncv> // Similar to PETSc, these options will override those specified // above as long as EPSSetFromOptions() is called _after_ any // other customization routines. ierr = EPSSetFromOptions (_eps); LIBMESH_CHKERR(ierr); // If the SolverConfiguration object is provided, use it to override // solver options. if (this->_solver_configuration) { this->_solver_configuration->configure_solver(); } // Solve the eigenproblem. ierr = EPSSolve (_eps); LIBMESH_CHKERR(ierr); // Get the number of iterations. ierr = EPSGetIterationNumber (_eps, &its); LIBMESH_CHKERR(ierr); // Get number of converged eigenpairs. ierr = EPSGetConverged(_eps,&nconv); LIBMESH_CHKERR(ierr); #ifdef DEBUG // ierr = PetscPrintf(this->comm().get(), // "\n Number of iterations: %d\n" // " Number of converged eigenpairs: %d\n\n", its, nconv); // Display eigenvalues and relative errors. ierr = PetscPrintf(this->comm().get(), " k ||Ax-kx||/|kx|\n" " ----------------- -----------------\n" ); LIBMESH_CHKERR(ierr); for (PetscInt i=0; i<nconv; i++ ) { ierr = EPSGetEigenpair(_eps, i, &kr, &ki, PETSC_NULL, PETSC_NULL); LIBMESH_CHKERR(ierr); #if SLEPC_VERSION_LESS_THAN(3,6,0) ierr = EPSComputeRelativeError(_eps, i, &error); #else ierr = EPSComputeError(_eps, i, EPS_ERROR_RELATIVE, &error); #endif LIBMESH_CHKERR(ierr); #ifdef LIBMESH_USE_COMPLEX_NUMBERS re = PetscRealPart(kr); im = PetscImaginaryPart(kr); #else re = kr; im = ki; #endif if (im != .0) { ierr = PetscPrintf(this->comm().get()," %9f%+9f i %12f\n", re, im, error); LIBMESH_CHKERR(ierr); } else { ierr = PetscPrintf(this->comm().get()," %12f %12f\n", re, error); LIBMESH_CHKERR(ierr); } } ierr = PetscPrintf(this->comm().get(),"\n" ); LIBMESH_CHKERR(ierr); #endif // DEBUG // return the number of converged eigenpairs // and the number of iterations return std::make_pair(nconv, its); }
//! Function to give PETSc that sets the Interpolation Matrix between two DMs PetscErrorCode libmesh_petsc_DMCreateInterpolation (DM dmc /*coarse*/, DM dmf /*fine*/, Mat * mat ,Vec * vec) { libmesh_assert(dmc); libmesh_assert(dmf); libmesh_assert(mat); libmesh_assert(vec); // Optional scaling (not needed for mg) // Get a communicator from incoming DM PetscErrorCode ierr; MPI_Comm comm; PetscObjectGetComm((PetscObject)dmc, &comm); // Extract our coarse context from the incoming DM void * ctx_c = nullptr; ierr = DMShellGetContext(dmc, &ctx_c); LIBMESH_CHKERR(ierr); libmesh_assert(ctx_c); PetscDMContext * p_ctx_c = static_cast<PetscDMContext*>(ctx_c); // Extract our fine context from the incoming DM void * ctx_f = nullptr; ierr = DMShellGetContext(dmf, &ctx_f);LIBMESH_CHKERR(ierr); libmesh_assert(ctx_f); PetscDMContext * p_ctx_f = static_cast<PetscDMContext*>(ctx_f); // Check for existing projection matrix libmesh_assert(p_ctx_c->K_interp_ptr); libmesh_assert(p_ctx_c->K_sub_interp_ptr); // If were doing fieldsplit we need to construct sub projection // matrices. We compare the passed in number of DMs fields to a // global DM in order to determine if a subprojection is needed. PetscInt nfieldsc,nfieldsf, nfieldsg; libmesh_assert(p_ctx_c->global_dm); DM * globaldm = p_ctx_c->global_dm; ierr = DMCreateFieldIS(dmc, &nfieldsc, nullptr, nullptr); LIBMESH_CHKERR(ierr); ierr = DMCreateFieldIS(dmf, &nfieldsf, nullptr, nullptr); LIBMESH_CHKERR(ierr); ierr = DMCreateFieldIS(*globaldm, &nfieldsg, nullptr, nullptr); LIBMESH_CHKERR(ierr); // If subfields are identified, were doing FS so we need to create the subProjectionMatrix if (nfieldsc < nfieldsg) { // Loop over the fields and merge their index sets. std::vector<std::vector<numeric_index_type>> allrows,allcols; std::vector<numeric_index_type> rows,cols; allrows = p_ctx_f->dof_vec; allcols = p_ctx_c->dof_vec; // For internal libmesh submat extraction need to merge all // field dofs and then sort the vectors so that they match // the Projection Matrix ordering const int n_subfields = nfieldsc; if ( n_subfields > 1 ) { for (int i = 0 ; i < n_subfields ; i++) { rows.insert(rows.end(), allrows[i].begin(), allrows[i].end()); cols.insert(cols.end(), allcols[i].begin(), allcols[i].end()); } std::sort(rows.begin(),rows.end()); std::sort(cols.begin(),cols.end()); } // Now that we have merged the fine and coarse index sets // were ready to make the submatrix and pass it off to PETSc p_ctx_c->K_interp_ptr->create_submatrix (*p_ctx_c->K_sub_interp_ptr, rows, cols); *(mat) = p_ctx_c->K_sub_interp_ptr->mat(); } else // We are not doing fieldsplit, so return entire projection *(mat) = p_ctx_c->K_interp_ptr->mat(); // Vec scaling isnt needed so were done. *(vec) = PETSC_NULL; return 0; } // end libmesh_petsc_DMCreateInterpolation
//! Help PETSc identify the coarser DM dmc given the fine DM dmf PetscErrorCode libmesh_petsc_DMCoarsen(DM dmf, MPI_Comm /*comm*/, DM * dmc) { libmesh_assert(dmc); libmesh_assert(dmf); PetscErrorCode ierr; // Extract our context from the incoming dmf void * ctx_f = nullptr; ierr = DMShellGetContext(dmf, &ctx_f);LIBMESH_CHKERR(ierr); libmesh_assert(ctx_f); PetscDMContext * p_ctx = static_cast<PetscDMContext*>(ctx_f); // First, ensure that there exists a coarse DM that we want to // set. There ought to be as we created it while walking the // hierarchy. libmesh_assert(p_ctx->coarser_dm); libmesh_assert(*(p_ctx->coarser_dm)); // In situations using fieldsplit we need to (potentially) // provide a coarser DM which only has the relevant subfields in // it. Since we create global DMs for each mesh level, we need // to extract the section from the DM, and check the number of // fields. When less than all the fields are used, we need to // create the proper subsections. // Get the number of fields and their names from the incomming // fine DM and the global reference DM PetscInt nfieldsf, nfieldsg; char ** fieldnamesf; char ** fieldnamesg; libmesh_assert(p_ctx->global_dm); DM * globaldm = p_ctx->global_dm; ierr = DMCreateFieldIS(dmf, &nfieldsf, &fieldnamesf, nullptr); LIBMESH_CHKERR(ierr); ierr = DMCreateFieldIS(*globaldm, &nfieldsg, &fieldnamesg, nullptr); LIBMESH_CHKERR(ierr); // If the probed number of fields is less than the number of // global fields, this amounts to PETSc 'indicating' to us we // are doing FS. So, we must create subsections for the coarser // DMs. if ( nfieldsf < nfieldsg ) { PetscSection section; PetscSection subsection; std::vector<PetscInt> subfields(nfieldsf); // extracted fields // First, get the section from the coarse DM #if PETSC_VERSION_LESS_THAN(3,10,0) ierr = DMGetDefaultSection(*(p_ctx->coarser_dm), §ion); #else ierr = DMGetSection(*(p_ctx->coarser_dm), §ion); #endif LIBMESH_CHKERR(ierr); // Now, match fine grid DM field names to their global DM // counterparts. Since PETSc can internally reassign field // numbering under a fieldsplit, we must extract // subsections via the field names. This is admittedly // gross, but c'est la vie. for (int i = 0; i < nfieldsf ; i++) { for (int j = 0; j < nfieldsg ;j++) if ( strcmp( fieldnamesg[j], fieldnamesf[i] ) == 0 ) subfields[i] = j; } // Next, for the found fields we now make a subsection and set it for the coarser DM ierr = PetscSectionCreateSubsection(section, nfieldsf, subfields.data(), &subsection); LIBMESH_CHKERR(ierr); #if PETSC_VERSION_LESS_THAN(3,10,0) ierr = DMSetDefaultSection(*(p_ctx->coarser_dm) , subsection); #else ierr = DMSetSection(*(p_ctx->coarser_dm) , subsection); #endif LIBMESH_CHKERR(ierr); ierr = PetscSectionDestroy(&subsection); LIBMESH_CHKERR(ierr); } // Finally, set the coarser DM *(dmc) = *(p_ctx->coarser_dm); return 0; }
void SlepcEigenSolver<T>:: set_slepc_position_of_spectrum() { PetscErrorCode ierr = 0; switch (this->_position_of_spectrum) { case LARGEST_MAGNITUDE: { ierr = EPSSetWhichEigenpairs (_eps, EPS_LARGEST_MAGNITUDE); LIBMESH_CHKERR(ierr); return; } case SMALLEST_MAGNITUDE: { ierr = EPSSetWhichEigenpairs (_eps, EPS_SMALLEST_MAGNITUDE); LIBMESH_CHKERR(ierr); return; } case LARGEST_REAL: { ierr = EPSSetWhichEigenpairs (_eps, EPS_LARGEST_REAL); LIBMESH_CHKERR(ierr); return; } case SMALLEST_REAL: { ierr = EPSSetWhichEigenpairs (_eps, EPS_SMALLEST_REAL); LIBMESH_CHKERR(ierr); return; } case LARGEST_IMAGINARY: { ierr = EPSSetWhichEigenpairs (_eps, EPS_LARGEST_IMAGINARY); LIBMESH_CHKERR(ierr); return; } case SMALLEST_IMAGINARY: { ierr = EPSSetWhichEigenpairs (_eps, EPS_SMALLEST_IMAGINARY); LIBMESH_CHKERR(ierr); return; } // The EPS_TARGET_XXX enums were added in SLEPc 3.1 #if !SLEPC_VERSION_LESS_THAN(3,1,0) case TARGET_MAGNITUDE: { ierr = EPSSetTarget(_eps, this->_target_val); LIBMESH_CHKERR(ierr); ierr = EPSSetWhichEigenpairs (_eps, EPS_TARGET_MAGNITUDE); LIBMESH_CHKERR(ierr); return; } case TARGET_REAL: { ierr = EPSSetTarget(_eps, this->_target_val); LIBMESH_CHKERR(ierr); ierr = EPSSetWhichEigenpairs (_eps, EPS_TARGET_REAL); LIBMESH_CHKERR(ierr); return; } case TARGET_IMAGINARY: { ierr = EPSSetTarget(_eps, this->_target_val); LIBMESH_CHKERR(ierr); ierr = EPSSetWhichEigenpairs (_eps, EPS_TARGET_IMAGINARY); LIBMESH_CHKERR(ierr); return; } #endif default: libmesh_error_msg("ERROR: Unsupported SLEPc position of spectrum: " << this->_position_of_spectrum); } }