libMesh
L2system.h
Go to the documentation of this file.
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 // libMesh includes
19 #include "libmesh/enum_fe_family.h"
20 #include "libmesh/fdm_gradient.h"
21 #include "libmesh/fem_function_base.h"
22 #include "libmesh/fem_system.h"
23 #include "libmesh/libmesh_common.h"
24 
25 // C++ includes
26 #include <map>
27 #include <memory>
28 
29 
30 // FEMSystem, TimeSolver and NewtonSolver will handle most tasks,
31 // but we must specify element residuals
33 {
34 public:
35  // Constructor
37  const std::string & name,
38  const unsigned int number)
39  : libMesh::FEMSystem(es, name, number),
40  input_system(nullptr),
41  _fe_family("LAGRANGE"),
42  _fe_order(1),
43  _hilbert_order(0),
44  _subdomains_list() {}
45 
46  // Default destructor
48 
49  std::string & fe_family() { return _fe_family; }
50  unsigned int & fe_order() { return _fe_order; }
51  std::set<libMesh::subdomain_id_type> & subdomains_list() { return _subdomains_list; }
52 
53  unsigned int & hilbert_order() { return _hilbert_order; }
55  _fdm_eps = eps;
56  if (_goal_func.get())
58  }
59 
61  {
62  _goal_func = goal.clone();
63  _goal_grad = std::make_unique<libMesh::FDMGradient<libMesh::Gradient>>(*_goal_func, _fdm_eps);
64  }
65 
66  // We want to be able to project functions based on *other* systems'
67  // values. For that we need not only a FEMFunction but also a
68  // reference to the system where it applies and a separate context
69  // object (or multiple separate context objects, in the threaded
70  // case) for that system.
72 
73 protected:
74  std::unique_ptr<libMesh::FEMFunctionBase<libMesh::Number> > _goal_func;
75 
76  std::map<libMesh::FEMContext *, std::unique_ptr<libMesh::FEMContext>>
78 
79  // System initialization
80  virtual void init_data ();
81 
82  // Context initialization
83  virtual void init_context (libMesh::DiffContext & context);
84 
85  // Element residual and jacobian calculations
86  // Time dependent parts
87  virtual bool element_time_derivative (bool request_jacobian,
88  libMesh::DiffContext & context);
89 
90  // The FE type to use
91  std::string _fe_family;
92  unsigned int _fe_order;
93 
94  // The Hilbert order our subclass will project with
95  unsigned int _hilbert_order;
96 
97  // The function we will call to finite difference our goal
98  // function
99  std::unique_ptr<libMesh::FDMGradient<libMesh::Gradient>> _goal_grad;
100 
101  // The perturbation we will use when finite differencing our goal
102  // function
104 
105  // Which subdomains to integrate on (all subdomains, if empty())
106  std::set<libMesh::subdomain_id_type> _subdomains_list;
107 };
std::set< libMesh::subdomain_id_type > & subdomains_list()
Definition: L2system.h:51
This is the EquationSystems class.
std::set< libMesh::subdomain_id_type > _subdomains_list
Definition: L2system.h:106
unsigned int _fe_order
Definition: L2system.h:92
This class provides all data required for a physics package (e.g.
Definition: diff_context.h:55
HilbertSystem(libMesh::EquationSystems &es, const std::string &name, const unsigned int number)
Definition: L2system.h:36
std::map< libMesh::FEMContext *, std::unique_ptr< libMesh::FEMContext > > input_contexts
Definition: L2system.h:77
virtual bool element_time_derivative(bool request_jacobian, libMesh::DiffContext &context)
Adds the time derivative contribution on elem to elem_residual.
Definition: L2system.C:88
unsigned int _hilbert_order
Definition: L2system.h:95
The libMesh namespace provides an interface to certain functionality in the library.
This class provides a specific system class.
Definition: fem_system.h:53
unsigned int number() const
Definition: system.h:2350
void set_goal_func(libMesh::FEMFunctionBase< libMesh::Number > &goal)
Definition: L2system.h:60
libMesh::Real _fdm_eps
Definition: L2system.h:103
Manages consistently variables, degrees of freedom, and coefficient vectors.
Definition: system.h:96
FEMSystem(EquationSystems &es, const std::string &name, const unsigned int number)
Constructor.
Definition: fem_system.C:858
std::string & fe_family()
Definition: L2system.h:49
unsigned int & fe_order()
Definition: L2system.h:50
unsigned int & hilbert_order()
Definition: L2system.h:53
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
libMesh::System * input_system
Definition: L2system.h:71
const std::string & name() const
Definition: system.h:2342
std::string _fe_family
Definition: L2system.h:91
virtual std::unique_ptr< FEMFunctionBase< Output > > clone() const =0
FEMFunctionBase is a base class from which users can derive in order to define "function-like" object...
std::unique_ptr< libMesh::FDMGradient< libMesh::Gradient > > _goal_grad
Definition: L2system.h:99
std::unique_ptr< libMesh::FEMFunctionBase< libMesh::Number > > _goal_func
Definition: L2system.h:74
virtual void init_data()
Initializes the member data fields associated with the system, so that, e.g., assemble() may be used...
Definition: L2system.C:34
void set_fdm_eps(libMesh::Real eps)
Definition: L2system.h:54
virtual void init_context(libMesh::DiffContext &context)
Definition: L2system.C:45