libMesh
Functions
vector_fe_ex10.C File Reference

Go to the source code of this file.

Functions

void assemble_graddiv (EquationSystems &es, const std::string &system_name)
 
int main (int argc, char **argv)
 
void assemble_graddiv (EquationSystems &es, const std::string &libmesh_dbg_var(system_name))
 

Function Documentation

◆ assemble_graddiv() [1/2]

void assemble_graddiv ( EquationSystems es,
const std::string &  system_name 
)

Referenced by main().

◆ assemble_graddiv() [2/2]

void assemble_graddiv ( EquationSystems es,
const std::string &  libmesh_dbg_varsystem_name 
)

Definition at line 231 of file vector_fe_ex10.C.

References libMesh::SparseMatrix< T >::add_matrix(), libMesh::NumericVector< T >::add_vector(), libMesh::FEGenericBase< OutputType >::build(), libMesh::FEType::default_quadrature_order(), dim, GradDivExactSolution::forcing(), libMesh::System::get_dof_map(), libMesh::EquationSystems::get_mesh(), libMesh::EquationSystems::get_system(), libMesh::ImplicitSystem::get_system_matrix(), mesh, libMesh::MeshBase::mesh_dimension(), libMesh::Real, libMesh::DenseVector< T >::resize(), libMesh::DenseMatrix< T >::resize(), libMesh::ExplicitSystem::rhs, and libMesh::System::variable_number().

233 {
234 
235  // It is a good idea to make sure we are assembling
236  // the proper system.
237  libmesh_assert_equal_to (system_name, "GradDiv");
238 
239  // Get a constant reference to the mesh object.
240  const MeshBase & mesh = es.get_mesh();
241 
242  // The dimension that we are running.
243  const unsigned int dim = mesh.mesh_dimension();
244 
245  // Get a reference to the LinearImplicitSystem we are solving.
246  LinearImplicitSystem & system = es.get_system<LinearImplicitSystem>("GradDiv");
247 
248  // A reference to the DofMap object for this system. The DofMap
249  // object handles the index translation from node and element numbers
250  // to degree of freedom numbers.
251  const DofMap & dof_map = system.get_dof_map();
252 
253  // Get a constant reference to the Finite Element type
254  // for the variable in the system.
255  FEType vector_fe_type = dof_map.variable_type(system.variable_number("u"));
256 
257  // Build the Finite Element object. Since the
258  // FEBase::build() member dynamically creates memory we will
259  // store the object as a std::unique_ptr<FEBase>. This can be thought
260  // of as a pointer that will clean up after itself. Introduction Example 4
261  // describes some advantages of std::unique_ptr's in the context of
262  // quadrature rules.
263  std::unique_ptr<FEVectorBase> vector_fe (FEVectorBase::build(dim, vector_fe_type));
264 
265  // A just-high-enough Gauss quadrature rule for numerical integration.
266  QGauss qrule (dim, vector_fe_type.default_quadrature_order());
267 
268  // Tell the finite element object to use our quadrature rule.
269  vector_fe->attach_quadrature_rule (&qrule);
270 
271  // Declare a special finite element object for boundary integration.
272  std::unique_ptr<FEVectorBase> vector_fe_face (FEVectorBase::build(dim, vector_fe_type));
273 
274  // Boundary integration requires one quadrature rule with dimensionality one
275  // less than the dimensionality of the element.
276  QGauss qface(dim-1, vector_fe_type.default_quadrature_order());
277 
278  // Tell the finite element object to use our quadrature rule.
279  vector_fe_face->attach_quadrature_rule (&qface);
280 
281  // Here we define some references to cell-specific data that
282  // will be used to assemble the linear system.
283  //
284  // The element Jacobian * quadrature weight at each integration point.
285  const std::vector<Real> & JxW = vector_fe->get_JxW();
286 
287  // The physical XY locations of the quadrature points on the element.
288  // These might be useful for evaluating spatially varying material
289  // properties at the quadrature points.
290  const std::vector<Point> & q_point = vector_fe->get_xyz();
291 
292  // The element shape functions evaluated at the quadrature points.
293  const std::vector<std::vector<RealGradient>> & vector_phi = vector_fe->get_phi();
294 
295  // The divergence of the element vector shape functions evaluated at the
296  // quadrature points.
297  const std::vector<std::vector<Real>> & div_vector_phi = vector_fe->get_div_phi();
298 
299  // Define data structures to contain the element matrix
300  // and right-hand-side vector contribution. Following
301  // basic finite element terminology we will denote these
302  // "Ke" and "Fe". These datatypes are templated on
303  // Number, which allows the same code to work for real
304  // or complex numbers.
307 
308  // These vectors will hold the degree of freedom indices for
309  // the element. These define where in the global system
310  // the element degrees of freedom get mapped.
311  std::vector<dof_id_type> dof_indices;
312 
313  // The global system matrix
314  SparseMatrix<Number> & matrix = system.get_system_matrix();
315 
316  // Now we will loop over all the elements in the mesh.
317  // We will compute the element matrix and right-hand-side
318  // contribution.
319  //
320  // Element ranges are a nice way to iterate through all the
321  // elements, or all the elements that have some property. The
322  // range will iterate from the first to the last element on
323  // the local processor.
324  // It is smart to make this one const so that we don't accidentally
325  // mess it up! In case users later modify this program to include
326  // refinement, we will be safe and will only consider the active
327  // elements; hence we use a variant of the
328  // active_local_element_ptr_range.
329  for (const auto & elem : mesh.active_local_element_ptr_range())
330  {
331  // Get the degree of freedom indices for the
332  // current element. These define where in the global
333  // matrix and right-hand-side this element will
334  // contribute to.
335  dof_map.dof_indices (elem, dof_indices);
336 
337  // Cache the total number of degrees of freedom on this element,
338  // for use as array and loop bounds later.
339  // We use cast_int to explicitly convert from size() (which may be
340  // 64-bit) to unsigned int (which may be 32-bit but which is definitely
341  // enough to count *local* degrees of freedom.
342  const unsigned int n_dofs =
343  cast_int<unsigned int>(dof_indices.size());
344 
345  // Compute the element-specific data for the current
346  // element. This involves computing the location of the
347  // quadrature points (q_point) and the shape functions
348  // and their divergences for the current element.
349  vector_fe->reinit (elem);
350 
351  // We should also have the same number of degrees of freedom as
352  // shape functions for our variable.
353  libmesh_assert_equal_to (n_dofs, vector_phi.size());
354 
355  // Zero the element matrix and right-hand side before
356  // summing them. We use the resize member here because
357  // the number of degrees of freedom might have changed from
358  // the last element. Note that this will be the case if the
359  // element type is different (i.e. the last element was a
360  // triangle, now we are on a quadrilateral).
361 
362  // The DenseMatrix::resize() and the DenseVector::resize()
363  // members will automatically zero out the matrix and vector.
364  Ke.resize (n_dofs, n_dofs);
365  Fe.resize (n_dofs);
366 
367  // Now loop over the quadrature points. This handles
368  // the numeric integration.
369  for (unsigned int qp=0; qp<qrule.n_points(); qp++)
370  {
371 
372  // Now we will build the element matrix.
373  // This a double loop to integrate the vector test functions (i)
374  // against the vector trial functions (j) and their divergences
375  for (unsigned int i = 0; i != n_dofs; i++)
376  for (unsigned int j = 0; j != n_dofs; j++)
377  {
378  Ke(i, j) += JxW[qp]*(div_vector_phi[i][qp]*div_vector_phi[j][qp]+
379  vector_phi[i][qp]*vector_phi[j][qp]);
380  }
381 
382  // This is the end of the matrix summation loop
383  // Now we build the element right-hand-side contribution.
384  // This involves a single loop in which we integrate the "forcing
385  // function" in the PDE against the vector test functions (k).
386  {
387  // "f" is the forcing function, given by the well-known
388  // "method of manufactured solutions".
389  RealGradient f = GradDivExactSolution().forcing(q_point[qp]);
390 
391  // Loop to integrate the vector test functions (k) against the
392  // forcing function.
393  for (unsigned int k = 0; k != n_dofs; k++)
394  {
395  Fe(k) += JxW[qp]*f*vector_phi[k][qp];
396  }
397  }
398  }
399 
400  // We have now reached the end of the quadrature point loop, so
401  // the interior element integration has been completed. However, we have
402  // not yet addressed boundary conditions.
403  {
404 
405  // The following loop is over the sides of the element.
406  // If the element has no neighbor on a side then that
407  // side MUST live on a boundary of the domain.
408  for (auto side : elem->side_index_range())
409  if (elem->neighbor_ptr(side) == nullptr)
410  {
411  // The value of the shape functions at the quadrature points.
412  const std::vector<std::vector<RealGradient>> & vector_phi_face = vector_fe_face->get_phi();
413 
414  // The Jacobian * Quadrature Weight at the quadrature
415  // points on the face.
416  const std::vector<Real> & JxW_face = vector_fe_face->get_JxW();
417 
418  // The XYZ locations (in physical space) of, and the normals at,
419  // the quadrature points on the face. This is where
420  // we will interpolate the boundary value function.
421  const std::vector<Point> & qface_point = vector_fe_face->get_xyz();
422  const std::vector<Point> & normals = vector_fe_face->get_normals();
423 
424  // Compute the vector shape function values on the element face.
425  vector_fe_face->reinit(elem, side);
426 
427  // Some shape functions will be 0 on the face, but for ease of
428  // indexing and generality of code we loop over them anyway.
429  libmesh_assert_equal_to (n_dofs, vector_phi_face.size());
430 
431  // Loop over the face quadrature points for integration.
432  for (unsigned int qp=0; qp<qface.n_points(); qp++)
433  {
434  // The boundary value for the vector variable.
435  RealGradient vector_value = GradDivExactSolution()(qface_point[qp]);
436 
437  // We use the penalty method to set the flux of the vector
438  // variable at the boundary, i.e. the RT vector boundary dof.
439  const Real penalty = 1.e10;
440 
441  // A double loop to integrate the normal component of the
442  // vector test functions (i) against the normal component of
443  // the vector trial functions (j).
444  for (unsigned int i = 0; i != n_dofs; i++)
445  for (unsigned int j = 0; j != n_dofs; j++)
446  {
447  Ke(i, j) += JxW_face[qp]*penalty*vector_phi_face[i][qp]*
448  normals[qp]*vector_phi_face[j][qp]*normals[qp];
449  }
450 
451  // Loop to integrate the normal component of the vector test
452  // functions (i) against the normal component of the
453  // exact solution for the vector variable.
454  for (unsigned int i = 0; i != n_dofs; i++)
455  {
456  Fe(i) += JxW_face[qp]*penalty*vector_phi_face[i][qp]*normals[qp]*
457  vector_value*normals[qp];
458  }
459  }
460  }
461  }
462 
463  // We have now finished the quadrature point loop,
464  // and have therefore applied all the boundary conditions.
465 
466  // If this assembly program were to be used on an adaptive mesh,
467  // we would have to apply any hanging node constraint equations.
468  dof_map.constrain_element_matrix_and_vector (Ke, Fe, dof_indices);
469 
470  // The element matrix and right-hand-side are now built
471  // for this element. Add them to the global matrix and
472  // right-hand-side vector. The SparseMatrix::add_matrix()
473  // and NumericVector::add_vector() members do this for us.
474  matrix.add_matrix (Ke, dof_indices);
475  system.rhs->add_vector (Fe, dof_indices);
476  }
477 
478  // All done!
479 }
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:196
unsigned int dim
Manages consistently variables, degrees of freedom, coefficient vectors, matrices and linear solvers ...
void resize(const unsigned int n)
Resize the vector.
Definition: dense_vector.h:396
virtual void add_vector(const T *v, const std::vector< numeric_index_type > &dof_indices)
Computes , where v is a pointer and each dof_indices[i] specifies where to add value v[i]...
MeshBase & mesh
NumericVector< Number > * rhs
The system matrix.
Order default_quadrature_order() const
Definition: fe_type.h:371
const T_sys & get_system(std::string_view name) const
This is the MeshBase class.
Definition: mesh_base.h:75
unsigned int variable_number(std::string_view var) const
Definition: system.C:1609
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:179
virtual void add_matrix(const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &rows, const std::vector< numeric_index_type > &cols)=0
Add the full matrix dm to the SparseMatrix.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const MeshBase & get_mesh() const
void resize(const unsigned int new_m, const unsigned int new_n)
Resizes the matrix to the specified size and calls zero().
Definition: dense_matrix.h:895
This class implements specific orders of Gauss quadrature.
unsigned int mesh_dimension() const
Definition: mesh_base.C:372
RealGradient forcing(Point p)
const DofMap & get_dof_map() const
Definition: system.h:2374
const SparseMatrix< Number > & get_system_matrix() const

◆ main()

int main ( int  argc,
char **  argv 
)

Definition at line 83 of file vector_fe_ex10.C.

References libMesh::EquationSystems::add_system(), libMesh::System::add_variable(), assemble_graddiv(), libMesh::System::attach_assemble_function(), libMesh::ExactSolution::attach_exact_derivs(), libMesh::ExactSolution::attach_exact_values(), libMesh::MeshTools::Generation::build_cube(), libMesh::MeshTools::Generation::build_square(), libMesh::ExactSolution::compute_error(), libMesh::default_solver_package(), libMesh::ExactSolution::error_norm(), libMesh::ExactSolution::extra_quadrature_order(), libMesh::FIFTH, libMesh::FIRST, libMesh::ExactSolution::hdiv_error(), libMesh::HDIV_SEMINORM, libMesh::TriangleWrapper::init(), libMesh::EquationSystems::init(), libMesh::INVALID_SOLVER_PACKAGE, libMesh::ExactSolution::l2_error(), mesh, libMesh::out, libMesh::MeshTools::Modification::permute_elements(), libMesh::EquationSystems::print_info(), libMesh::MeshBase::print_info(), libMesh::RAVIART_THOMAS, libMesh::Real, GradDivExactSolution::RM(), libMesh::MeshTools::Modification::rotate(), libMesh::LinearImplicitSystem::solve(), and libMesh::ExodusII_IO::write_equation_systems().

84 {
85  // Initialize libMesh.
86  LibMeshInit init (argc, argv);
87 
88  // This example requires a linear solver package.
89  libmesh_example_requires(libMesh::default_solver_package() != INVALID_SOLVER_PACKAGE,
90  "--enable-petsc, --enable-trilinos, or --enable-eigen");
91 
92  // Parse the input file.
93  GetPot infile("vector_fe_ex10.in");
94 
95  // But allow the command line to override it.
96  infile.parse_command_line(argc, argv);
97 
98  // Read in parameters from the command line and the input file.
99  const unsigned int dimension = infile("dim", 2);
100  const unsigned int grid_size = infile("grid_size", 15);
101 
102  // Skip higher-dimensional examples on a lower-dimensional libMesh build.
103  libmesh_example_requires(dimension <= LIBMESH_DIM, dimension << "D support");
104 
105  // Create a mesh, with dimension to be overridden later, distributed
106  // across the default MPI communicator.
107  Mesh mesh(init.comm());
108 
109  // Use the MeshTools::Generation mesh generator to create a uniform
110  // grid on the cube [-1,1]^D. To accomodate Raviart-Thomas elements, we must
111  // use TRI6/7 or QUAD8/9 elements in 2d, or TET14 or HEX27 in 3d.
112  const std::string elem_str = infile("element_type", std::string("TRI6"));
113 
114  libmesh_error_msg_if((dimension == 2 && elem_str != "TRI6" && elem_str != "TRI7" && elem_str != "QUAD8" && elem_str != "QUAD9") ||
115  (dimension == 3 && elem_str != "TET14" && elem_str != "HEX27"),
116  "You selected " << elem_str <<
117  " but this example must be run with TRI6, TRI7, QUAD8, or QUAD9 in 2d" <<
118  " or with TET14, or HEX27 in 3d.");
119 
120  if (dimension == 2)
122  grid_size,
123  grid_size,
124  -1., 1.,
125  -1., 1.,
126  Utility::string_to_enum<ElemType>(elem_str));
127  else if (dimension == 3)
129  grid_size,
130  grid_size,
131  grid_size,
132  -1., 1.,
133  -1., 1.,
134  -1., 1.,
135  Utility::string_to_enum<ElemType>(elem_str));
136 
137  // Make sure the code is robust against nodal reorderings.
139 
140  // Make sure the code is robust against solves on 2d meshes rotated out of
141  // the xy plane. By default, all Euler angles are zero, the rotation matrix
142  // is the identity, and the mesh stays in place.
143  const Real phi = infile("phi", 0.), theta = infile("theta", 0.), psi = infile("psi", 0.);
145 
146  // Print information about the mesh to the screen.
147  mesh.print_info();
148 
149  // Create an equation systems object.
150  EquationSystems equation_systems (mesh);
151 
152  // Declare the system "GradDiv" and its variable.
153  LinearImplicitSystem & system = equation_systems.add_system<LinearImplicitSystem>("GradDiv");
154 
155  // Set the FE approximation order for the vector field variable.
156  const Order vector_order = static_cast<Order>(infile("order", 1u));
157 
158  libmesh_error_msg_if(vector_order < FIRST || vector_order > ((dimension == 3) ? FIRST : FIFTH),
159  "You selected: " << vector_order <<
160  " but this example must be run with either 1 <= order <= 5 in 2d"
161  " or with order 1 in 3d.");
162 
163  // Adds the variable "u" to "GradDiv". "u" will be our vector field.
164  system.add_variable("u", vector_order, RAVIART_THOMAS);
165 
166  // Give the system a pointer to the matrix assembly
167  // function. This will be called when needed by the library.
169 
170  // Initialize the data structures for the equation system.
171  equation_systems.init();
172 
173  // Prints information about the system to the screen.
174  equation_systems.print_info();
175 
176  // Solve the system "GradDiv". Note that calling this
177  // member will assemble the linear system and invoke
178  // the default numerical solver.
179  system.solve();
180 
181  ExactSolution exact_sol(equation_systems);
182 
183  SolutionFunction soln_func;
184  SolutionGradient soln_grad;
185 
186  // Build FunctionBase* containers to attach to the ExactSolution object.
187  std::vector<FunctionBase<Number> *> sols(1, &soln_func);
188  std::vector<FunctionBase<Gradient> *> grads(1, &soln_grad);
189 
190  exact_sol.attach_exact_values(sols);
191  exact_sol.attach_exact_derivs(grads);
192 
193  // Use higher quadrature order for more accurate error results.
194  int extra_error_quadrature = infile("extra_error_quadrature", 2);
195  exact_sol.extra_quadrature_order(extra_error_quadrature);
196 
197  // Compute the error.
198  exact_sol.compute_error("GradDiv", "u");
199 
200  // Print out the error values.
201  libMesh::out << "~~ Vector field (u) ~~"
202  << std::endl;
203  libMesh::out << "L2 error is: "
204  << exact_sol.l2_error("GradDiv", "u")
205  << std::endl;
206  libMesh::out << "HDiv semi-norm error is: "
207  << exact_sol.error_norm("GradDiv", "u", HDIV_SEMINORM)
208  << std::endl;
209  libMesh::out << "HDiv error is: "
210  << exact_sol.hdiv_error("GradDiv", "u")
211  << std::endl;
212 
213 #ifdef LIBMESH_HAVE_EXODUS_API
214 
215  // We write the file in the ExodusII format.
216  ExodusII_IO(mesh).write_equation_systems("out.e", equation_systems);
217 
218 #endif // #ifdef LIBMESH_HAVE_EXODUS_API
219 
220  // All done.
221  return 0;
222 }
This class handles the computation of the L2 and/or H1 error for the Systems in the EquationSystems o...
This is the EquationSystems class.
Order
defines an enum for polynomial orders.
Definition: enum_order.h:40
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
Definition: exodusII_io.h:52
void assemble_graddiv(EquationSystems &es, const std::string &system_name)
Manages consistently variables, degrees of freedom, coefficient vectors, matrices and linear solvers ...
MeshBase & mesh
void build_square(UnstructuredMesh &mesh, const unsigned int nx, const unsigned int ny, const Real xmin=0., const Real xmax=1., const Real ymin=0., const Real ymax=1., const ElemType type=INVALID_ELEM, const bool gauss_lobatto_grid=false)
A specialized build_cube() for 2D meshes.
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
Definition: libmesh.h:90
virtual void solve() override
Assembles & solves the linear system A*x=b.
void permute_elements(MeshBase &mesh)
Randomly permute the nodal ordering of each element (without twisting the element mapping)...
SolverPackage default_solver_package()
Definition: libmesh.C:1117
virtual void write_equation_systems(const std::string &fname, const EquationSystems &es, const std::set< std::string > *system_names=nullptr) override
Writes out the solution for no specific time or timestep.
Definition: exodusII_io.C:2033
void print_info(std::ostream &os=libMesh::out, const unsigned int verbosity=0, const bool global=true) const
Prints relevant information about the mesh.
Definition: mesh_base.C:1562
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
unsigned int add_variable(std::string_view 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:1357
static void RM(RealTensor T)
void attach_assemble_function(void fptr(EquationSystems &es, const std::string &name))
Register a user function to use in assembling the system matrix and RHS.
Definition: system.C:2161
RealTensorValue rotate(MeshBase &mesh, const Real phi, const Real theta=0., const Real psi=0.)
Rotates the mesh in the xy plane.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
OStreamProxy out
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
Definition: mesh.h:50
void build_cube(UnstructuredMesh &mesh, const unsigned int nx=0, const unsigned int ny=0, const unsigned int nz=0, const Real xmin=0., const Real xmax=1., const Real ymin=0., const Real ymax=1., const Real zmin=0., const Real zmax=1., const ElemType type=INVALID_ELEM, const bool gauss_lobatto_grid=false)
Builds a (elements) cube.