Line data Source code
1 : // Copyright Amazon.com, Inc. or its affiliates. All Rights Reserved.
2 : // SPDX-License-Identifier: Apache-2.0
3 :
4 : #ifndef PALACE_LINALG_AMS_HPP
5 : #define PALACE_LINALG_AMS_HPP
6 :
7 : #include <memory>
8 : #include <mfem.hpp>
9 : #include "linalg/operator.hpp"
10 : #include "utils/iodata.hpp"
11 :
12 : namespace palace
13 : {
14 :
15 : class FiniteElementSpace;
16 :
17 : //
18 : // A wrapper for Hypre's AMS solver.
19 : //
20 : class HypreAmsSolver : public mfem::HypreSolver
21 : {
22 : private:
23 : // The Hypre solver object.
24 : HYPRE_Solver ams;
25 :
26 : // Parameters used for preconditioner construction.
27 : const int cycle_type, space_dim, ams_it, ams_smooth_it;
28 : const bool ams_singular, agg_coarsen;
29 :
30 : // Control print level for debugging.
31 : const int print;
32 :
33 : // Discrete gradient matrix (not owned).
34 : const mfem::HypreParMatrix *G;
35 :
36 : // Nedelec interpolation matrix and its components, or, for p = 1, the mesh vertex
37 : // coordinates.
38 : std::unique_ptr<mfem::HypreParMatrix> Pi, Pix, Piy, Piz;
39 : std::unique_ptr<mfem::HypreParVector> x, y, z;
40 :
41 : // Helper function to set up the auxiliary objects required by the AMS solver.
42 : void ConstructAuxiliaryMatrices(FiniteElementSpace &nd_fespace,
43 : FiniteElementSpace &h1_fespace);
44 :
45 : // Helper function to construct and configure the AMS solver.
46 : void InitializeSolver();
47 :
48 : public:
49 : // Constructor requires the ND space, but will construct the H1 and (H1)ᵈ spaces
50 : // internally as needed.
51 : HypreAmsSolver(FiniteElementSpace &nd_fespace, FiniteElementSpace &h1_fespace,
52 : int cycle_it, int smooth_it, bool vector_interp, bool singular_op,
53 : bool agg_coarsen, int print);
54 0 : HypreAmsSolver(const IoData &iodata, bool coarse_solver, FiniteElementSpace &nd_fespace,
55 : FiniteElementSpace &h1_fespace, int print)
56 0 : : HypreAmsSolver(
57 : nd_fespace, h1_fespace, coarse_solver ? 1 : iodata.solver.linear.mg_cycle_it,
58 0 : iodata.solver.linear.mg_smooth_it, iodata.solver.linear.ams_vector_interp,
59 0 : iodata.solver.linear.ams_singular_op, iodata.solver.linear.amg_agg_coarsen, print)
60 : {
61 0 : }
62 : ~HypreAmsSolver() override;
63 :
64 : void SetOperator(const Operator &op) override;
65 :
66 0 : operator HYPRE_Solver() const override { return ams; }
67 :
68 0 : HYPRE_PtrToParSolverFcn SetupFcn() const override
69 : {
70 0 : return (HYPRE_PtrToParSolverFcn)HYPRE_AMSSetup;
71 : }
72 :
73 0 : HYPRE_PtrToParSolverFcn SolveFcn() const override
74 : {
75 0 : return (HYPRE_PtrToParSolverFcn)HYPRE_AMSSolve;
76 : }
77 : };
78 :
79 : } // namespace palace
80 :
81 : #endif // PALACE_LINALG_AMS_HPP
|