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_MODELS_LAPLACE_OPERATOR_HPP
5 : #define PALACE_MODELS_LAPLACE_OPERATOR_HPP
6 :
7 : #include <map>
8 : #include <memory>
9 : #include <vector>
10 : #include <mfem.hpp>
11 : #include "fem/fespace.hpp"
12 : #include "linalg/operator.hpp"
13 : #include "linalg/vector.hpp"
14 : #include "models/materialoperator.hpp"
15 :
16 : namespace palace
17 : {
18 :
19 : class IoData;
20 : class Mesh;
21 :
22 : //
23 : // A class handling discretization of Laplace problems for electrostatics.
24 : //
25 : class LaplaceOperator
26 : {
27 : private:
28 : // Helper variable for log file printing.
29 : bool print_hdr;
30 :
31 : // Essential boundary condition markers.
32 : mfem::Array<int> dbc_attr;
33 : std::vector<mfem::Array<int>> dbc_tdof_lists;
34 :
35 : // Objects defining the finite element spaces for the electrostatic potential (H1) and
36 : // electric field (Nedelec) on the given mesh. The RT spaces are used for error
37 : // estimation.
38 : std::vector<std::unique_ptr<mfem::H1_FECollection>> h1_fecs;
39 : std::unique_ptr<mfem::ND_FECollection> nd_fec;
40 : std::vector<std::unique_ptr<mfem::RT_FECollection>> rt_fecs;
41 : FiniteElementSpaceHierarchy h1_fespaces;
42 : FiniteElementSpace nd_fespace;
43 : FiniteElementSpaceHierarchy rt_fespaces;
44 :
45 : // Operator for domain material properties.
46 : MaterialOperator mat_op;
47 :
48 : // Boundary attributes for each terminal index.
49 : std::map<int, mfem::Array<int>> source_attr_lists;
50 :
51 : mfem::Array<int> SetUpBoundaryProperties(const IoData &iodata, const mfem::ParMesh &mesh);
52 : std::map<int, mfem::Array<int>> ConstructSources(const IoData &iodata);
53 :
54 : public:
55 : LaplaceOperator(const IoData &iodata, const std::vector<std::unique_ptr<Mesh>> &mesh);
56 :
57 : // Return material operator for postprocessing.
58 12 : const MaterialOperator &GetMaterialOp() const { return mat_op; }
59 :
60 : // Access source attribute lists.
61 : const auto &GetSources() const { return source_attr_lists; }
62 :
63 : // Return the parallel finite element space objects.
64 0 : auto &GetH1Spaces() { return h1_fespaces; }
65 : const auto &GetH1Spaces() const { return h1_fespaces; }
66 : auto &GetH1Space() { return h1_fespaces.GetFinestFESpace(); }
67 : const auto &GetH1Space() const { return h1_fespaces.GetFinestFESpace(); }
68 6 : auto &GetNDSpace() { return nd_fespace; }
69 : const auto &GetNDSpace() const { return nd_fespace; }
70 : auto &GetRTSpaces() { return rt_fespaces; }
71 : const auto &GetRTSpaces() const { return rt_fespaces; }
72 : auto &GetRTSpace() { return rt_fespaces.GetFinestFESpace(); }
73 : const auto &GetRTSpace() const { return rt_fespaces.GetFinestFESpace(); }
74 :
75 : // Access the underlying mesh object.
76 : const auto &GetMesh() const { return GetH1Space().GetMesh(); }
77 :
78 : // Return the number of true (conforming) dofs on the finest H1 space.
79 : auto GlobalTrueVSize() const { return GetH1Space().GlobalTrueVSize(); }
80 :
81 : // Construct and return system matrix representing discretized Laplace operator for
82 : // Gauss's law.
83 : std::unique_ptr<Operator> GetStiffnessMatrix();
84 :
85 : // Construct and return the discrete gradient matrix.
86 : const Operator &GetGradMatrix() const
87 : {
88 3 : return GetNDSpace().GetDiscreteInterpolator(GetH1Space());
89 : }
90 :
91 : // Assemble the solution boundary conditions and right-hand side vector for a nonzero
92 : // prescribed voltage on the specified surface index.
93 : void GetExcitationVector(int idx, const Operator &K, Vector &X, Vector &RHS);
94 :
95 : // Get the associated MPI communicator.
96 : MPI_Comm GetComm() const { return GetH1Space().GetComm(); }
97 : };
98 :
99 : } // namespace palace
100 :
101 : #endif // PALACE_MODELS_LAPLACE_OPERATOR_HPP
|