libMesh
L-shaped.h
Go to the documentation of this file.
1 #include "libmesh/enum_fe_family.h"
2 #include "libmesh/fem_system.h"
3 #include "libmesh/qoi_set.h"
4 #include "libmesh/system.h"
5 
6 // Bring in everything from the libMesh namespace
7 using namespace libMesh;
8 
9 // FEMSystem, TimeSolver and NewtonSolver will handle most tasks,
10 // but we must specify element residuals
11 class LaplaceSystem : public FEMSystem
12 {
13 public:
14  // Constructor
16  const std::string & name_in,
17  const unsigned int number_in) :
18  FEMSystem(es, name_in, number_in),
19  _fe_family("LAGRANGE"),
20  _fe_order(1),
21  _analytic_jacobians(true)
22  { this->init_qois(2); }
23 
24  std::string & fe_family() { return _fe_family; }
25  unsigned int & fe_order() { return _fe_order; }
26  bool & analytic_jacobians() { return _analytic_jacobians; }
27 
28  // Postprocessing function which we are going to override for this application
29 
30  virtual void postprocess();
31 
32  Number & get_QoI_value(std::string type,
33  unsigned int QoI_index)
34  {
35  if (type == "exact")
36  return exact_QoI[QoI_index];
37  else
38  return computed_QoI[QoI_index];
39  }
40 
41 protected:
42  // System initialization
43  virtual void init_data ();
44 
45  // Context initialization
46  virtual void init_context (DiffContext & context);
47 
48  // Element residual and jacobian calculations
49  // Time dependent parts
50  virtual bool element_time_derivative (bool request_jacobian,
51  DiffContext & context);
52 
53  // Constraint parts
54  virtual bool side_constraint (bool request_jacobian,
55  DiffContext & context);
56 
57  // Overloading the postprocess function
58 
59  virtual void element_postprocess(DiffContext & context);
60 
61  virtual void side_postprocess(DiffContext & context);
62 
63  // Overloading the qoi function on elements
64 
65  virtual void element_qoi_derivative (DiffContext & context,
66  const QoISet & qois);
67 
68  // Overloading the qoi function on sides
69 
70  virtual void side_qoi_derivative (DiffContext & context,
71  const QoISet & qois);
72 
73  Number exact_solution (const Point &);
74 
75  // Variables to hold the computed QoIs
76 
77  Number computed_QoI[2];
78 
79  // Variables to read in the exact QoIs from l-shaped.in
80 
81  Number exact_QoI[2];
82 
83  // The FE type to use
84  std::string _fe_family;
85  unsigned int _fe_order;
86 
87  // Calculate Jacobians analytically or not?
89 };
This is the EquationSystems class.
This class provides all data required for a physics package (e.g.
Definition: diff_context.h:55
Number(* exact_solution)(const Point &p, const Parameters &, const std::string &, const std::string &)
Data structure for specifying which Quantities of Interest should be calculated in an adjoint or a pa...
Definition: qoi_set.h:45
Number & get_QoI_value(std::string type, unsigned int QoI_index)
Definition: L-shaped.h:32
std::string & fe_family()
Definition: L-shaped.h:24
LaplaceSystem(EquationSystems &es, const std::string &name_in, const unsigned int number_in)
Definition: L-shaped.h:15
bool _analytic_jacobians
Definition: L-shaped.h:88
unsigned int _fe_order
Definition: L-shaped.h:85
bool & analytic_jacobians()
Definition: L-shaped.h:26
The libMesh namespace provides an interface to certain functionality in the library.
This class provides a specific system class.
Definition: fem_system.h:53
std::string _fe_family
Definition: L-shaped.h:84
unsigned int & fe_order()
Definition: L-shaped.h:25
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:39