libMesh
vector_fe_ex2.C
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 
19 
20 // <h1>FEMSystem Example 1 - Unsteady Navier-Stokes Equations with
21 // FEMSystem</h1>
22 // \author Paul Bauman
23 // \date 2012
24 //
25 // This example shows how the transient nonlinear problem from
26 // systems_of_equations_ex2 can be solved using the
27 // DifferentiableSystem class framework
28 
29 // Basic include files
30 #include "libmesh/equation_systems.h"
31 #include "libmesh/getpot.h"
32 #include "libmesh/exodusII_io.h"
33 #include "libmesh/mesh.h"
34 #include "libmesh/mesh_generation.h"
35 #include "libmesh/exact_solution.h"
36 #include "libmesh/ucd_io.h"
37 #include "libmesh/enum_solver_package.h"
38 
39 // The systems and solvers we may use
40 #include "laplace_system.h"
41 #include "libmesh/diff_solver.h"
42 #include "libmesh/steady_solver.h"
43 
44 // C++ includes
45 #include <memory>
46 
47 // Bring in everything from the libMesh namespace
48 using namespace libMesh;
49 
50 // The main program.
51 int main (int argc, char** argv)
52 {
53  // Initialize libMesh.
54  LibMeshInit init (argc, argv);
55 
56  // This example requires a linear solver package.
57  libmesh_example_requires(libMesh::default_solver_package() != INVALID_SOLVER_PACKAGE,
58  "--enable-petsc, --enable-trilinos, or --enable-eigen");
59 
60  // Parse the input file
61  GetPot infile("vector_fe_ex2.in");
62 
63  // But allow the command line to override it.
64  infile.parse_command_line(argc, argv);
65 
66  // Read in parameters from the input file
67  const unsigned int grid_size = infile("grid_size", 2);
68 
69  // Skip higher-dimensional examples on a lower-dimensional libMesh build
70  libmesh_example_requires(3 <= LIBMESH_DIM, "2D/3D support");
71 
72  // We use Dirichlet boundary conditions here
73 #ifndef LIBMESH_ENABLE_DIRICHLET
74  libmesh_example_requires(false, "--enable-dirichlet");
75 #endif
76 
77  // Create a mesh, with dimension to be overridden later, on the
78  // default MPI communicator.
79  Mesh mesh(init.comm());
80 
81  // Use the MeshTools::Generation mesh generator to create a uniform
82  // grid on the square [-1,1]^D. We instruct the mesh generator
83  // to build a mesh of 8x8 Quad9 elements in 2D, or Hex27
84  // elements in 3D. Building these higher-order elements allows
85  // us to use higher-order approximation, as in example 3.
87  grid_size,
88  grid_size,
89  grid_size,
90  -1., 1.,
91  -1., 1.,
92  -1., 1.,
93  HEX8);
94 
95  // Print information about the mesh to the screen.
96  mesh.print_info();
97 
98  // Create an equation systems object.
99  EquationSystems equation_systems (mesh);
100 
101  // Declare the system "Laplace" and its variables.
102  LaplaceSystem & system =
103  equation_systems.add_system<LaplaceSystem> ("Laplace");
104 
105  // This example only implements the steady-state problem
106  system.time_solver = std::make_unique<SteadySolver>(system);
107 
108  // Initialize the system
109  equation_systems.init();
110 
111  // And the nonlinear solver options
112  DiffSolver & solver = *(system.time_solver->diff_solver().get());
113  solver.quiet = infile("solver_quiet", true);
114  solver.verbose = !solver.quiet;
115  solver.max_nonlinear_iterations = infile("max_nonlinear_iterations", 15);
116  solver.relative_step_tolerance = infile("relative_step_tolerance", 1.e-3);
117  solver.relative_residual_tolerance = infile("relative_residual_tolerance", 0.0);
118  solver.absolute_residual_tolerance = infile("absolute_residual_tolerance", 0.0);
119 
120  // And the linear solver options
121  solver.max_linear_iterations = infile("max_linear_iterations", 50000);
122  solver.initial_linear_tolerance = infile("initial_linear_tolerance", 1.e-3);
123 
124  // Print information about the system to the screen.
125  equation_systems.print_info();
126 
127  system.solve();
128 
129  ExactSolution exact_sol(equation_systems);
130 
131  SolutionFunction soln_func(system.variable_number("u"));
132  SolutionGradient soln_grad(system.variable_number("u"));
133 
134  // Build FunctionBase* containers to attach to the ExactSolution object.
135  std::vector<FunctionBase<Number> *> sols(1, &soln_func);
136  std::vector<FunctionBase<Gradient> *> grads(1, &soln_grad);
137 
138  exact_sol.attach_exact_values(sols);
139  exact_sol.attach_exact_derivs(grads);
140 
141  // Use higher quadrature order for more accurate error results
142  int extra_error_quadrature = infile("extra_error_quadrature", 2);
143  exact_sol.extra_quadrature_order(extra_error_quadrature);
144 
145  // Compute the error.
146  exact_sol.compute_error("Laplace", "u");
147 
148  // Print out the error values
149  libMesh::out << "L2-Error is: "
150  << exact_sol.l2_error("Laplace", "u")
151  << std::endl;
152  libMesh::out << "H1-Error is: "
153  << exact_sol.h1_error("Laplace", "u")
154  << std::endl;
155 
156 #ifdef LIBMESH_HAVE_EXODUS_API
157 
158  // We write the file in the ExodusII format.
159  ExodusII_IO(mesh).write_equation_systems("out.e", equation_systems);
160 
161 #endif // #ifdef LIBMESH_HAVE_EXODUS_API
162 
163  UCDIO(mesh).write_equation_systems("out.inp", equation_systems);
164 
165  // All done.
166  return 0;
167 }
This class handles the computation of the L2 and/or H1 error for the Systems in the EquationSystems o...
This is the EquationSystems class.
virtual void write_equation_systems(const std::string &, const EquationSystems &, const std::set< std::string > *system_names=nullptr)
This method implements writing a mesh with data to a specified file where the data is taken from the ...
Definition: mesh_output.C:31
bool quiet
The DiffSolver should not print anything to libMesh::out unless quiet is set to false; default is tru...
Definition: diff_solver.h:160
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
Definition: exodusII_io.h:52
void attach_exact_values(const std::vector< FunctionBase< Number > *> &f)
Clone and attach arbitrary functors which compute the exact values of the EquationSystems&#39; solutions ...
unsigned int max_nonlinear_iterations
The DiffSolver should exit in failure if max_nonlinear_iterations is exceeded and continue_after_max_...
Definition: diff_solver.h:154
Real absolute_residual_tolerance
The DiffSolver should exit after the residual is reduced to either less than absolute_residual_tolera...
Definition: diff_solver.h:189
std::unique_ptr< TimeSolver > time_solver
A pointer to the solver object we&#39;re going to use.
Definition: diff_system.h:253
void print_info(std::ostream &os=libMesh::out) const
Prints information about the equation systems, by default to libMesh::out.
MeshBase & mesh
void extra_quadrature_order(const int extraorder)
Increases or decreases the order of the quadrature rule used for numerical integration.
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
Definition: libmesh.h:90
The libMesh namespace provides an interface to certain functionality in the library.
void attach_exact_derivs(const std::vector< FunctionBase< Gradient > *> &g)
Clone and attach arbitrary functors which compute the exact gradients of the EquationSystems&#39; solutio...
unsigned int variable_number(std::string_view var) const
Definition: system.C:1609
SolverPackage default_solver_package()
Definition: libmesh.C:1117
This is a generic class that defines a solver to handle ImplicitSystem classes, including NonlinearIm...
Definition: diff_solver.h:68
This class implements reading & writing meshes in the AVS&#39;s UCD format.
Definition: ucd_io.h:44
unsigned int max_linear_iterations
Each linear solver step should exit after max_linear_iterations is exceeded.
Definition: diff_solver.h:146
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 compute_error(std::string_view sys_name, std::string_view unknown_name)
Computes and stores the error in the solution value e = u-u_h, the gradient grad(e) = grad(u) - grad(...
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
Real l2_error(std::string_view sys_name, std::string_view unknown_name)
bool verbose
The DiffSolver may print a lot more to libMesh::out if verbose is set to true; default is false...
Definition: diff_solver.h:166
OStreamProxy out
virtual void solve() override
Invokes the solver associated with the system.
Definition: fem_system.C:1076
virtual void init()
Initialize all the systems.
virtual System & add_system(std::string_view system_type, std::string_view name)
Add the system of type system_type named name to the systems array.
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
Definition: mesh.h:50
Real h1_error(std::string_view sys_name, std::string_view unknown_name)
double initial_linear_tolerance
Any required linear solves will at first be done with this tolerance; the DiffSolver may tighten the ...
Definition: diff_solver.h:208
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.
Real relative_residual_tolerance
Definition: diff_solver.h:190
int main(int argc, char **argv)
Definition: vector_fe_ex2.C:51