Line data Source code
1 : // The libMesh Finite Element Library.
2 : // Copyright (C) 2002-2025 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 :
4 : // This library is free software; you can redistribute it and/or
5 : // modify it under the terms of the GNU Lesser General Public
6 : // License as published by the Free Software Foundation; either
7 : // version 2.1 of the License, or (at your option) any later version.
8 :
9 : // This library is distributed in the hope that it will be useful,
10 : // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 : // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 : // Lesser General Public License for more details.
13 :
14 : // You should have received a copy of the GNU Lesser General Public
15 : // License along with this library; if not, write to the Free Software
16 : // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 :
18 : #include "libmesh/petsc_auto_fieldsplit.h"
19 :
20 : #ifdef LIBMESH_HAVE_PETSC
21 :
22 : #include <petscksp.h>
23 :
24 : // Local includes
25 : #include "libmesh/dof_map.h"
26 : #include "libmesh/system.h"
27 :
28 : namespace {
29 : using namespace libMesh;
30 :
31 0 : void indices_to_fieldsplit (const Parallel::Communicator & comm,
32 : const std::vector<dof_id_type> & indices,
33 : PC my_pc,
34 : const std::string & field_name)
35 : {
36 0 : const PetscInt * idx = LIBMESH_PETSC_NULLPTR;
37 0 : if (!indices.empty())
38 0 : idx = reinterpret_cast<const PetscInt *>(indices.data());
39 :
40 : IS is;
41 0 : LibmeshPetscCall2(comm, ISCreateGeneral(comm.get(), cast_int<PetscInt>(indices.size()),
42 : idx, PETSC_COPY_VALUES, &is));
43 :
44 0 : LibmeshPetscCall2(comm, PCFieldSplitSetIS(my_pc, field_name.c_str(), is));
45 0 : }
46 :
47 : }
48 :
49 : namespace libMesh
50 : {
51 :
52 219656 : void petsc_auto_fieldsplit (PC my_pc,
53 : const System & sys)
54 : {
55 226252 : std::string sys_prefix = "--solver_group_";
56 :
57 432716 : if (libMesh::on_command_line("--solver-system-names"))
58 : {
59 0 : sys_prefix = sys_prefix + sys.name() + "_";
60 : }
61 :
62 13192 : std::map<std::string, std::vector<dof_id_type>> group_indices;
63 :
64 432716 : if (libMesh::on_command_line("--solver-variable-names"))
65 : {
66 0 : for (auto v : make_range(sys.n_vars()))
67 : {
68 0 : const std::string & var_name = sys.variable_name(v);
69 :
70 0 : std::vector<dof_id_type> var_idx;
71 0 : sys.get_dof_map().local_variable_indices
72 0 : (var_idx, sys.get_mesh(), v);
73 :
74 0 : std::string group_command = sys_prefix + var_name;
75 :
76 0 : const std::string empty_string;
77 :
78 : std::string group_name = libMesh::command_line_value
79 0 : (group_command, empty_string);
80 :
81 0 : if (group_name != empty_string)
82 : {
83 : std::vector<dof_id_type> & indices =
84 0 : group_indices[group_name];
85 0 : const bool prior_indices = !indices.empty();
86 0 : indices.insert(indices.end(), var_idx.begin(),
87 0 : var_idx.end());
88 0 : if (prior_indices)
89 0 : std::sort(indices.begin(), indices.end());
90 : }
91 : else
92 : {
93 0 : indices_to_fieldsplit (sys.comm(), var_idx, my_pc, var_name);
94 : }
95 : }
96 : }
97 :
98 219656 : for (const auto & [field_name, indices] : group_indices)
99 0 : indices_to_fieldsplit(sys.comm(), indices, my_pc, field_name);
100 219656 : }
101 :
102 : } // namespace libMesh
103 :
104 :
105 : #endif // #ifdef LIBMESH_HAVE_PETSC
|