MFEM v4.8.0
Finite element discretization library
Loading...
Searching...
No Matches
hybridization_ext.hpp
Go to the documentation of this file.
1// Copyright (c) 2010-2025, Lawrence Livermore National Security, LLC. Produced
2// at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3// LICENSE and NOTICE for details. LLNL-CODE-806117.
4//
5// This file is part of the MFEM library. For more information and source code
6// availability visit https://mfem.org.
7//
8// MFEM is free software; you can redistribute it and/or modify it under the
9// terms of the BSD-3 license. We welcome feedback and contributions, see file
10// CONTRIBUTING.md for details.
11
12#ifndef MFEM_HYBRIDIZATION_EXT
13#define MFEM_HYBRIDIZATION_EXT
14
15#include "../config/config.hpp"
16#include "../general/array.hpp"
17#include "../linalg/vector.hpp"
18
19namespace mfem
20{
21
22/// @brief Extension class supporting Hybridization on device (GPU).
23///
24/// Similar to BilinearFormExtension and LinearFormExtension, this extension
25/// class provides device execution capabilities for the Hybridization class.
26///
27/// As with the other extension classes, a limitation of this class is that it
28/// requires meshes consisting only of tensor-product elements, and finite
29/// element spaces without variable polynomial degrees.
31{
32 friend class Hybridization;
33public:
34 enum DofType : char
35 {
39 };
40protected:
41 class Hybridization &h; ///< The associated Hybridization object.=
42 int num_hat_dofs; ///< Number of Lagrange multipliers.
43 mutable Vector tmp1, tmp2; ///< Temporary vectors.
44
47
50 Vector Ct_mat; ///< Constraint matrix (transposed) stored element-wise.
51
53
56
57public:
58 /// Construct the constraint matrix.
59 void ConstructC();
60
61 template <int MID, int MBD>
62 void FactorElementMatrices(Vector &AhatInvCt_mat);
63
64 /// Form the Schur complement matrix $H$.
65 void ConstructH();
66
67 /// Compute the action of C^t x.
68 void MultCt(const Vector &x, Vector &y) const;
69
70 /// Compute the action of C x.
71 void MultC(const Vector &x, Vector &y) const;
72
73 /// @brief Assemble the element matrix A into the hybridized system matrix.
74 ///
75 /// @warning Using the interface will be very slow. AssembleElementMatrices()
76 /// should be used instead.
77 void AssembleMatrix(int el, const class DenseMatrix &elmat);
78
79 /// @brief Assemble the boundary element matrix A into the hybridized system
80 /// matrix.
81 ///
82 /// @warning Using the interface will be very slow. AssembleElementMatrices()
83 /// should be used instead.
84 void AssembleBdrMatrix(int bdr_el, const class DenseMatrix &elmat);
85
86 /// Invert and store the element matrices Ahat.
87 void AssembleElementMatrices(const class DenseTensor &el_mats);
88
89 /// Apply the action of R mapping from "hat DOFs" to T-vector
90 void MultR(const Vector &b, Vector &b_hat) const;
91
92 /// Apply the action of R^t mapping into the "hat DOF" space.
93 void MultRt(const Vector &b, Vector &b_hat) const;
94
95 /// Apply the elementwise A_hat^{-1}.
96 void MultAhatInv(Vector &x) const;
97
98 /// Constructor.
99 HybridizationExtension(class Hybridization &hybridization_);
100
101 /// Prepare for assembly; form the constraint matrix.
102 void Init(const Array<int> &ess_tdof_list);
103
104 /// @brief Given a right-hand side on the original space, compute the
105 /// corresponding right-hand side for the Lagrange multipliers.
106 void ReduceRHS(const Vector &b, Vector &b_r) const;
107
108 /// @brief Given Lagrange multipliers @a sol_r and the original right-hand
109 /// side @a b, recover the solution @a sol on the original finite element
110 /// space.
111 void ComputeSolution(const Vector &b, const Vector &sol_r, Vector &sol) const;
112
113 /// Destroys the stored element matrices.
114 void Reset() { Ahat = 0.0; }
115};
116
117}
118
119#endif
Data type dense matrix using column-major storage.
Definition densemat.hpp:24
Rank 3 tensor (array of matrices)
Extension class supporting Hybridization on device (GPU).
Vector Ct_mat
Constraint matrix (transposed) stored element-wise.
void ConstructC()
Construct the constraint matrix.
void Init(const Array< int > &ess_tdof_list)
Prepare for assembly; form the constraint matrix.
void FactorElementMatrices(Vector &AhatInvCt_mat)
void ReduceRHS(const Vector &b, Vector &b_r) const
Given a right-hand side on the original space, compute the corresponding right-hand side for the Lagr...
void ComputeSolution(const Vector &b, const Vector &sol_r, Vector &sol) const
Given Lagrange multipliers sol_r and the original right-hand side b, recover the solution sol on the ...
void Reset()
Destroys the stored element matrices.
HybridizationExtension(class Hybridization &hybridization_)
Constructor.
Vector tmp2
Temporary vectors.
void ConstructH()
Form the Schur complement matrix .
void MultR(const Vector &b, Vector &b_hat) const
Apply the action of R mapping from "hat DOFs" to T-vector.
void MultC(const Vector &x, Vector &y) const
Compute the action of C x.
void MultAhatInv(Vector &x) const
Apply the elementwise A_hat^{-1}.
int num_hat_dofs
Number of Lagrange multipliers.
void AssembleMatrix(int el, const class DenseMatrix &elmat)
Assemble the element matrix A into the hybridized system matrix.
void MultCt(const Vector &x, Vector &y) const
Compute the action of C^t x.
class Hybridization & h
The associated Hybridization object.=.
void MultRt(const Vector &b, Vector &b_hat) const
Apply the action of R^t mapping into the "hat DOF" space.
void AssembleBdrMatrix(int bdr_el, const class DenseMatrix &elmat)
Assemble the boundary element matrix A into the hybridized system matrix.
void AssembleElementMatrices(const class DenseTensor &el_mats)
Invert and store the element matrices Ahat.
Auxiliary class Hybridization, used to implement BilinearForm hybridization.
Vector data type.
Definition vector.hpp:82
real_t b
Definition lissajous.cpp:42