libMesh
curl_curl_system.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2024 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 #ifndef CURL_CURL_SYSTEM_H
19 #define CURL_CURL_SYSTEM_H
20 
21 // DiffSystem framework files
22 #include "libmesh/fem_system.h"
23 #include "libmesh/vector_value.h"
24 #include "libmesh/tensor_value.h"
25 #include "libmesh/dirichlet_boundaries.h"
26 
27 #include "curl_curl_exact_solution.h"
28 
29 using namespace libMesh;
30 
35 class CurlCurlSystem : public FEMSystem
36 {
37 public:
38  // Constructor
40  const std::string & name,
41  const unsigned int number);
42 
43  // System initialization
44  virtual void init_data ();
45 
46  // Context initialization
47  virtual void init_context(DiffContext & context);
48 
49  // Element residual and jacobian calculations
50  // Time dependent parts
51  virtual bool element_time_derivative (bool request_jacobian,
52  DiffContext & context);
53 
54  virtual bool side_time_derivative(bool request_jacobian,
55  DiffContext & context);
56 
57 protected:
58  // Indices for each variable;
59  unsigned int u_var;
60 
61  void init_dirichlet_bcs();
62 
63  // Returns the value of a forcing function at point p.
64  RealGradient forcing(const Point & p);
65 
66  // Returns the value of the exact solution for this model problem.
68 
70 };
71 
72 #endif // CURL_CURL_SYSTEM_H
std::string name(const ElemQuality q)
This function returns a string containing some name for q.
Definition: elem_quality.C:42
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 &)
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 u_var
CurlCurlExactSolution soln
FEMSystem, TimeSolver and NewtonSolver will handle most tasks, but we must specify element residuals...
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:39