Line data Source code
1 : //* This file is part of the MOOSE framework 2 : //* https://mooseframework.inl.gov 3 : //* 4 : //* All rights reserved, see COPYRIGHT for full restrictions 5 : //* https://github.com/idaholab/moose/blob/master/COPYRIGHT 6 : //* 7 : //* Licensed under LGPL 2.1, please see LICENSE for details 8 : //* https://www.gnu.org/licenses/lgpl-2.1.html 9 : 10 : #include "libmesh/petsc_macro.h" 11 : #include "FieldSplitPreconditioner.h" 12 : 13 : // MOOSE includes 14 : #include "FEProblem.h" 15 : #include "MooseEnum.h" 16 : #include "MooseVariableFE.h" 17 : #include "NonlinearSystem.h" 18 : #include "PetscSupport.h" 19 : 20 : #include "libmesh/libmesh_common.h" 21 : #include "libmesh/petsc_nonlinear_solver.h" 22 : #include "libmesh/coupling_matrix.h" 23 : 24 : registerMooseObjectAliased("MooseApp", FieldSplitPreconditioner, "FSP"); 25 : 26 : InputParameters 27 14473 : FieldSplitPreconditioner::validParams() 28 : { 29 14473 : InputParameters params = MoosePreconditioner::validParams(); 30 14473 : params.addClassDescription("Preconditioner designed to map onto PETSc's PCFieldSplit."); 31 : 32 14473 : params.addRequiredParam<std::vector<std::string>>( 33 : "topsplit", "Entrance to splits, the top split will specify how splits will go."); 34 : // We should use full coupling Jacobian matrix by default 35 43419 : params.addParam<bool>("full", 36 28946 : true, 37 : "Set to true if you want the full set of couplings between variables " 38 : "simply for convenience so you don't have to set every off_diag_row " 39 : "and off_diag_column combination."); 40 14473 : return params; 41 0 : } 42 : 43 104 : FieldSplitPreconditioner::FieldSplitPreconditioner(const InputParameters & parameters) 44 : : MoosePreconditioner(parameters), 45 104 : _top_split(getParam<std::vector<std::string>>("topsplit")), 46 208 : _nl(_fe_problem.getNonlinearSystemBase(_nl_sys_num)) 47 : { 48 : // number of variables 49 104 : unsigned int n_vars = _nl.nVariables(); 50 : // if we want to construct a full Jacobian? 51 : // it is recommended to have a full Jacobian for using 52 : // the fieldSplit preconditioner 53 104 : bool full = getParam<bool>("full"); 54 : 55 : // how variables couple 56 104 : std::unique_ptr<CouplingMatrix> cm = std::make_unique<CouplingMatrix>(n_vars); 57 104 : if (!full) 58 : { 59 0 : if (isParamValid("off_diag_row") && isParamValid("off_diag_column")) 60 : { 61 : 62 0 : const auto off_diag_rows = getParam<std::vector<NonlinearVariableName>>("off_diag_row"); 63 0 : const auto off_diag_columns = getParam<std::vector<NonlinearVariableName>>("off_diag_column"); 64 : 65 : // put 1s on diagonal 66 0 : for (unsigned int i = 0; i < n_vars; i++) 67 0 : (*cm)(i, i) = 1; 68 : 69 : // off-diagonal entries 70 0 : std::vector<std::vector<unsigned int>> off_diag(n_vars); 71 0 : if (off_diag_rows.size() * off_diag_columns.size() != 0 && 72 0 : off_diag_rows.size() == off_diag_columns.size()) 73 0 : for (const auto i : index_range(off_diag_rows)) 74 : { 75 0 : unsigned int row = _nl.getVariable(0, off_diag_rows[i]).number(); 76 0 : unsigned int column = _nl.getVariable(0, off_diag_columns[i]).number(); 77 0 : (*cm)(row, column) = 1; 78 : } 79 0 : } 80 : } 81 : else 82 : { 83 328 : for (unsigned int i = 0; i < n_vars; i++) 84 720 : for (unsigned int j = 0; j < n_vars; j++) 85 496 : (*cm)(i, j) = 1; // full coupling 86 : } 87 104 : setCouplingMatrix(std::move(cm)); 88 : 89 : // turn on a flag 90 104 : _nl.useFieldSplitPreconditioner(true); 91 : 92 : // set a top splitting 93 104 : _nl.setDecomposition(_top_split); 94 : 95 : // apply prefix and store PETSc options 96 104 : _nl.setupFieldDecomposition(); 97 104 : }