MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
pweakform.hpp
Go to the documentation of this file.
1// Copyright (c) 2010-2024, 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_PDPGWEAKFORM
13#define MFEM_PDPGWEAKFORM
14
15#include "mfem.hpp"
16#include "weakform.hpp"
17
18#ifdef MFEM_USE_MPI
19
20namespace mfem
21{
22
23/** @brief Class representing the parallel weak formulation.
24 * (Convenient for DPG Equations) */
26{
27
28protected:
29 // Trial FE spaces
31
32 // ess_tdof list for each space
34
35 // split ess_tdof_list given in global tdof (for all spaces)
36 // to individual lists for each space
37 void FillEssTdofLists(const Array<int> & ess_tdof_list);
38
39 // Block operator of HypreParMatrix
40 BlockOperator * P = nullptr; // Block Prolongation
41 BlockMatrix * R = nullptr; // Block Restriction
42
43 // Block operator of HypreParMatrix
44 BlockOperator * p_mat = nullptr;
46
47 void BuildProlongation();
48
49private:
50
51public:
52
54
55 /// Creates bilinear form associated with FE spaces @a trial_pfes_.
58 : DPGWeakForm()
59 {
60 SetParSpaces(trial_pfes_,fecol_);
61 }
62
65 {
66 trial_pfes = trial_pfes_;
67 ess_tdofs.SetSize(trial_pfes.Size());
68
70 for (int i = 0; i<trial_sfes.Size(); i++)
71 {
72 trial_sfes[i] = (FiniteElementSpace *)trial_pfes[i];
73 ess_tdofs[i] = new Array<int>();
74 }
75
76 SetSpaces(trial_sfes,fecol_);
77 }
78
79
80 /// Assemble the local matrix
81 void Assemble(int skip_zeros = 1);
82
83 /// Returns the matrix assembled on the true dofs, i.e. P^t A P.
84 /** The returned matrix has to be deleted by the caller. */
86
87 void FormLinearSystem(const Array<int> &ess_tdof_list, Vector &x,
89 Vector &B, int copy_interior = 0);
90
91 void FormSystemMatrix(const Array<int> &ess_tdof_list,
93
94 /** Call this method after solving a linear system constructed using the
95 FormLinearSystem method to recover the solution as a ParGridFunction-size
96 vector in x. Use the same arguments as in the FormLinearSystem call. */
97 virtual void RecoverFEMSolution(const Vector &X, Vector &x);
98
99 virtual void Update();
100
101 /// Destroys bilinear form.
102 virtual ~ParDPGWeakForm();
103
104};
105
106} // namespace mfem
107
108
109#endif // MFEM_USE_MPI
110
111
112#endif
int Size() const
Return the logical size of the array.
Definition: array.hpp:144
A class to handle Block systems in a matrix-free implementation.
Class representing the DPG weak formulation. Given the variational formulation a(u,...
Definition: weakform.hpp:34
BlockMatrix * mat
Block matrix to be associated with the Block bilinear form. Owned.
Definition: weakform.hpp:49
void SetSpaces(Array< FiniteElementSpace * > &fes_, Array< FiniteElementCollection * > &fecol_)
Definition: weakform.hpp:125
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:220
Pointer to an Operator of a specified type.
Definition: handle.hpp:34
Class representing the parallel weak formulation. (Convenient for DPG Equations)
Definition: pweakform.hpp:26
BlockOperator * p_mat
Definition: pweakform.hpp:44
Array< ParFiniteElementSpace * > trial_pfes
Definition: pweakform.hpp:30
BlockOperator * P
Definition: pweakform.hpp:40
virtual ~ParDPGWeakForm()
Destroys bilinear form.
Definition: pweakform.cpp:206
void Assemble(int skip_zeros=1)
Assemble the local matrix.
Definition: pweakform.cpp:33
void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
Form the linear system A X = B, corresponding to this DPG weak form.
Definition: pweakform.cpp:100
BlockMatrix * R
Definition: pweakform.hpp:41
ParDPGWeakForm(Array< ParFiniteElementSpace * > &trial_pfes_, Array< FiniteElementCollection * > &fecol_)
Creates bilinear form associated with FE spaces trial_pfes_.
Definition: pweakform.hpp:56
void SetParSpaces(Array< ParFiniteElementSpace * > &trial_pfes_, Array< FiniteElementCollection * > &fecol_)
Definition: pweakform.hpp:63
virtual void RecoverFEMSolution(const Vector &X, Vector &x)
Definition: pweakform.cpp:173
void FillEssTdofLists(const Array< int > &ess_tdof_list)
Definition: pweakform.cpp:19
Array< Array< int > * > ess_tdofs
Definition: pweakform.hpp:33
BlockOperator * p_mat_e
Definition: pweakform.hpp:45
void ParallelAssemble(BlockMatrix *mat)
Returns the matrix assembled on the true dofs, i.e. P^t A P.
Definition: pweakform.cpp:38
void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
Form the linear system matrix A, see FormLinearSystem() for details.
Definition: pweakform.cpp:140
virtual void Update()
Update the DPGWeakForm after mesh modifications (AMR)
Definition: pweakform.cpp:188
Vector data type.
Definition: vector.hpp:80