libMesh
laplace_system.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2019 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 #include "libmesh/getpot.h"
19 
20 #include "laplace_system.h"
21 
22 #include "libmesh/dof_map.h"
23 #include "libmesh/fe_base.h"
24 #include "libmesh/fe_interface.h"
25 #include "libmesh/fem_context.h"
26 #include "libmesh/mesh.h"
27 #include "libmesh/quadrature.h"
28 #include "libmesh/string_to_enum.h"
29 
30 
31 // Bring in everything from the libMesh namespace
32 using namespace libMesh;
33 
35  const std::string & name_in,
36  const unsigned int number_in) :
37  FEMSystem(es, name_in, number_in)
38 {}
39 
41 {
42  // Check the input file for Reynolds number, application type,
43  // approximation type
44  GetPot infile("laplace.in");
45 
46  // Add the solution variable
47  u_var = this->add_variable ("u", FIRST, LAGRANGE_VEC);
48 
49  // The solution is evolving, with a first order time derivative
50  this->time_evolving(u_var, 1);
51 
52  // Useful debugging options
53  // Set verify_analytic_jacobians to 1e-6 to use
54  this->verify_analytic_jacobians = infile("verify_analytic_jacobians", 0.);
55  this->print_jacobians = infile("print_jacobians", false);
56  this->print_element_jacobians = infile("print_element_jacobians", false);
57 
58  this->init_dirichlet_bcs();
59 
60  // Do the parent's initialization after variables and boundary constraints are defined
61  FEMSystem::init_data();
62 }
63 
65 {
66 #ifdef LIBMESH_ENABLE_DIRICHLET
67  const boundary_id_type all_ids[6] = {0, 1, 2, 3, 4, 5};
68  std::set<boundary_id_type> boundary_ids(all_ids, all_ids+6);
69 
70  std::vector<unsigned int> vars;
71  vars.push_back(u_var);
72 
73  // Note that for vector-valued variables, it is assumed each
74  // component is stored contiguously. For 2-D elements in 3-D space,
75  // only two components should be returned.
76  SolutionFunction func(u_var);
77 
78  // In general, when reusing a system-indexed exact solution, we want
79  // to use the default system-ordering constructor for
80  // DirichletBoundary, so we demonstrate that here. In this case,
81  // though, we have only one (vector-valued!) variable, so system-
82  // and local- orderings are the same.
84  (libMesh::DirichletBoundary(boundary_ids, vars, func));
85 #endif // LIBMESH_ENABLE_DIRICHLET
86 }
87 
89 {
90  FEMContext & c = cast_ref<FEMContext &>(context);
91 
92  // Get finite element object
95 
96  // We should prerequest all the data
97  // we will need to build the linear system.
98  fe->get_JxW();
99  fe->get_phi();
100  fe->get_dphi();
101  fe->get_xyz();
102 
103  FEGenericBase<RealGradient> * side_fe;
104  c.get_side_fe<RealGradient>(u_var, side_fe);
105 
106  side_fe->get_JxW();
107  side_fe->get_phi();
108 }
109 
110 
111 bool LaplaceSystem::element_time_derivative (bool request_jacobian,
112  DiffContext & context)
113 {
114  FEMContext & c = cast_ref<FEMContext &>(context);
115 
116  // Get finite element object
117  FEGenericBase<RealGradient> * fe = nullptr;
119 
120  // First we get some references to cell-specific data that
121  // will be used to assemble the linear system.
122 
123  // Element Jacobian * quadrature weights for interior integration
124  const std::vector<Real> & JxW = fe->get_JxW();
125 
126  // The velocity shape functions at interior quadrature points.
127  const std::vector<std::vector<RealGradient>> & phi = fe->get_phi();
128 
129  // The velocity shape function gradients at interior
130  // quadrature points.
131  const std::vector<std::vector<RealTensor>> & grad_phi = fe->get_dphi();
132 
133  const std::vector<Point> & qpoint = fe->get_xyz();
134 
135  // The number of local degrees of freedom in each variable
136  const unsigned int n_u_dofs = c.n_dof_indices(u_var);
137 
139 
141 
142  // Now we will build the element Jacobian and residual.
143  // Constructing the residual requires the solution and its
144  // gradient from the previous timestep. This must be
145  // calculated at each quadrature point by summing the
146  // solution degree-of-freedom values by the appropriate
147  // weight functions.
148  const unsigned int n_qpoints = c.get_element_qrule().n_points();
149 
150  for (unsigned int qp=0; qp != n_qpoints; qp++)
151  {
152  Tensor grad_u;
153 
154  c.interior_gradient(u_var, qp, grad_u);
155 
156  // Value of the forcing function at this quadrature point
157  RealGradient f = this->forcing(qpoint[qp]);
158 
159  // First, an i-loop over the velocity degrees of freedom.
160  // We know that n_u_dofs == n_v_dofs so we can compute contributions
161  // for both at the same time.
162  for (unsigned int i=0; i != n_u_dofs; i++)
163  {
164  Fu(i) += (grad_u.contract(grad_phi[i][qp]) - f*phi[i][qp])*JxW[qp];
165 
166  // Matrix contributions for the uu and vv couplings.
167  if (request_jacobian)
168  for (unsigned int j=0; j != n_u_dofs; j++)
169  Kuu(i,j) += grad_phi[j][qp].contract(grad_phi[i][qp])*JxW[qp];
170  }
171  } // end of the quadrature point qp-loop
172 
173  return request_jacobian;
174 }
175 
176 
177 
178 // bool LaplaceSystem::side_constraint (bool request_jacobian,
179 // DiffContext & context)
180 // {
181 // FEMContext & c = cast_ref<FEMContext &>(context);
182 //
183 // // Get finite element object
184 // FEGenericBase<RealGradient> * side_fe = nullptr;
185 // c.get_side_fe<RealGradient>(u_var, side_fe);
186 //
187 // // First we get some references to cell-specific data that
188 // // will be used to assemble the linear system.
189 //
190 // // Element Jacobian * quadrature weights for interior integration
191 // const std::vector<Real> & JxW = side_fe->get_JxW();
192 //
193 // // The velocity shape functions at interior quadrature points.
194 // const std::vector<std::vector<RealGradient>> & phi = side_fe->get_phi();
195 //
196 // // The number of local degrees of freedom in each variable
197 // const unsigned int n_u_dofs = c.dof_indices_var[u_var].size();
198 //
199 // const std::vector<Point> & qpoint = side_fe->get_xyz();
200 //
201 // // The penalty value. \frac{1}{\epsilon}
202 // // in the discussion above.
203 // const Real penalty = 1.e10;
204 //
205 // DenseSubMatrix<Number> & Kuu = *c.elem_subjacobians[u_var][u_var];
206 // DenseSubVector<Number> & Fu = *c.elem_subresiduals[u_var];
207 //
208 // const unsigned int n_qpoints = (c.get_side_qrule())->n_points();
209 //
210 // for (unsigned int qp=0; qp != n_qpoints; qp++)
211 // {
212 // Gradient u;
213 // c.side_value(u_var, qp, u);
214 //
215 // Gradient u_exact(this->exact_solution(0, qpoint[qp](0), qpoint[qp](1)),
216 // this->exact_solution(1, qpoint[qp](0), qpoint[qp](1)));
217 //
218 // for (unsigned int i=0; i != n_u_dofs; i++)
219 // {
220 // Fu(i) += penalty*(u - u_exact)*phi[i][qp]*JxW[qp];
221 //
222 // if (request_jacobian)
223 // for (unsigned int j=0; j != n_u_dofs; j++)
224 // Kuu(i,j) += penalty*phi[j][qp]*phi[i][qp]*JxW[qp];
225 // }
226 // }
227 //
228 // return request_jacobian;
229 // }
230 
231 
232 
234 {
235  Real x = p(0);
236  Real y = p(1);
237  Real z = p(2);
238 
239  const Real eps = 1.e-3;
240 
241  const Real fx = -(exact_solution(0, x, y, z-eps) +
242  exact_solution(0, x, y, z+eps) +
243  exact_solution(0, x, y-eps, z) +
244  exact_solution(0, x, y+eps, z) +
245  exact_solution(0, x-eps, y, z) +
246  exact_solution(0, x+eps, y, z) -
247  6.*exact_solution(0, x, y, z))/eps/eps;
248 
249  const Real fy = -(exact_solution(1, x, y, z-eps) +
250  exact_solution(1, x, y, z+eps) +
251  exact_solution(1, x, y-eps, z) +
252  exact_solution(1, x, y+eps, z) +
253  exact_solution(1, x-eps, y, z) +
254  exact_solution(1, x+eps, y, z) -
255  6.*exact_solution(1, x, y, z))/eps/eps;
256 
257  const Real fz = -(exact_solution(2, x, y, z-eps) +
258  exact_solution(2, x, y, z+eps) +
259  exact_solution(2, x, y-eps, z) +
260  exact_solution(2, x, y+eps, z) +
261  exact_solution(2, x-eps, y, z) +
262  exact_solution(2, x+eps, y, z) -
263  6.*exact_solution(2, x, y, z))/eps/eps;
264 
265  return RealGradient(fx, fy, fz);
266 }
LaplaceSystem::u_var
unsigned int u_var
Definition: laplace_system.h:58
libMesh::DiffContext::get_elem_residual
const DenseVector< Number > & get_elem_residual() const
Const accessor for element residual.
Definition: diff_context.h:249
libMesh::FEMContext::get_element_qrule
const QBase & get_element_qrule() const
Accessor for element interior quadrature rule for the dimension of the current _elem.
Definition: fem_context.h:790
libMesh::DifferentiableSystem::print_element_jacobians
bool print_element_jacobians
Set print_element_jacobians to true to print each J_elem contribution.
Definition: diff_system.h:361
libMesh::DifferentiablePhysics::time_evolving
virtual void time_evolving(unsigned int var)
Tells the DiffSystem that variable var is evolving with respect to time.
Definition: diff_physics.h:250
libMesh::FEAbstract::get_xyz
const std::vector< Point > & get_xyz() const
Definition: fe_abstract.h:237
libMesh::DofMap::add_dirichlet_boundary
void add_dirichlet_boundary(const DirichletBoundary &dirichlet_boundary)
Adds a copy of the specified Dirichlet boundary to the system.
Definition: dof_map_constraints.C:4390
libMesh::DenseSubMatrix
Defines a dense submatrix for use in Finite Element-type computations.
Definition: dense_submatrix.h:45
libMesh::QBase::n_points
unsigned int n_points() const
Definition: quadrature.h:126
libMesh::RealGradient
RealVectorValue RealGradient
Definition: hp_coarsentest.h:49
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::FEGenericBase
This class forms the foundation from which generic finite elements may be derived.
Definition: exact_error_estimator.h:39
libMesh::FEAbstract::get_JxW
const std::vector< Real > & get_JxW() const
Definition: fe_abstract.h:244
libMesh::FEGenericBase::get_dphi
const std::vector< std::vector< OutputGradient > > & get_dphi() const
Definition: fe_base.h:214
LaplaceSystem::LaplaceSystem
LaplaceSystem(EquationSystems &es, const std::string &name_in, const unsigned int number_in)
Definition: L-shaped.h:15
LaplaceSystem::exact_solution
LaplaceExactSolution exact_solution
Definition: laplace_system.h:65
libMesh::boundary_id_type
int8_t boundary_id_type
Definition: id_types.h:51
LaplaceSystem::element_time_derivative
virtual bool element_time_derivative(bool request_jacobian, DiffContext &context)
Adds the time derivative contribution on elem to elem_residual.
Definition: L-shaped.C:58
libMesh::FEMSystem
This class provides a specific system class.
Definition: fem_system.h:53
libMesh::VectorValue< Real >
LaplaceSystem::forcing
RealGradient forcing(const Point &p)
Definition: laplace_system.C:233
SolutionFunction
Definition: solution_function.h:30
libMesh::System::add_variable
unsigned int add_variable(const std::string &var, const FEType &type, const std::set< subdomain_id_type > *const active_subdomains=nullptr)
Adds the variable var to the list of variables for this system.
Definition: system.C:1069
LaplaceSystem::init_dirichlet_bcs
void init_dirichlet_bcs()
Definition: laplace_system.C:64
libMesh::DifferentiableSystem::print_jacobians
bool print_jacobians
Set print_jacobians to true to print J whenever it is assembled.
Definition: diff_system.h:346
libMesh::TensorValue
This class defines a tensor in LIBMESH_DIM dimensional Real or Complex space.
Definition: exact_solution.h:53
libMesh::Point
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:38
libMesh::TypeTensor::contract
CompareTypes< T, T2 >::supertype contract(const TypeTensor< T2 > &) const
Multiply 2 tensors together to return a scalar, i.e.
Definition: type_tensor.h:1286
libMesh::DiffContext::n_dof_indices
unsigned int n_dof_indices() const
Total number of dof indices on the element.
Definition: diff_context.h:399
libMesh::DiffContext
This class provides all data required for a physics package (e.g.
Definition: diff_context.h:55
libMesh::EquationSystems
This is the EquationSystems class.
Definition: equation_systems.h:74
libMesh::DenseSubVector< Number >
laplace_system.h
libMesh::FEMContext::get_element_fe
void get_element_fe(unsigned int var, FEGenericBase< OutputShape > *&fe) const
Accessor for interior finite element object for variable var for the largest dimension in the mesh.
Definition: fem_context.h:275
libMesh::FEMContext::interior_gradient
Gradient interior_gradient(unsigned int var, unsigned int qp) const
Definition: fem_context.C:426
libMesh::FEMSystem::verify_analytic_jacobians
Real verify_analytic_jacobians
If verify_analytic_jacobian is equal to zero (as it is by default), no numeric jacobians will be calc...
Definition: fem_system.h:209
LaplaceSystem::init_context
virtual void init_context(DiffContext &context)
Definition: L-shaped.C:34
libMesh::DiffContext::get_elem_jacobian
const DenseMatrix< Number > & get_elem_jacobian() const
Const accessor for element Jacobian.
Definition: diff_context.h:283
libMesh::DirichletBoundary
This class allows one to associate Dirichlet boundary values with a given set of mesh boundary ids an...
Definition: dirichlet_boundaries.h:88
libMesh::FEMContext::get_side_fe
void get_side_fe(unsigned int var, FEGenericBase< OutputShape > *&fe) const
Accessor for edge/face (2D/3D) finite element object for variable var for the largest dimension in th...
Definition: fem_context.h:312
libMesh::System::get_dof_map
const DofMap & get_dof_map() const
Definition: system.h:2099
libMesh::LAGRANGE_VEC
Definition: enum_fe_family.h:60
LaplaceSystem::init_data
virtual void init_data()
Initializes the member data fields associated with the system, so that, e.g., assemble() may be used.
Definition: L-shaped.C:17
libMesh::FEGenericBase::get_phi
const std::vector< std::vector< OutputShape > > & get_phi() const
Definition: fe_base.h:206
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::FIRST
Definition: enum_order.h:42
libMesh::FEMContext
This class provides all data required for a physics package (e.g.
Definition: fem_context.h:62