libMesh
vector_fe_ex3.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>Vector Finite Elements Example 3 - Nedelec Elements</h1>
21 // \author Paul Bauman
22 // \date 2012
23 //
24 // This example shows an example of using the Nedelec elements of the
25 // first type to solve a model problem in H(curl).
26 
27 // Basic include files
28 #include "libmesh/equation_systems.h"
29 #include "libmesh/getpot.h"
30 #include "libmesh/exodusII_io.h"
31 #include "libmesh/mesh.h"
32 #include "libmesh/mesh_generation.h"
33 #include "libmesh/mesh_modification.h"
34 #include "libmesh/mesh_refinement.h"
35 #include "libmesh/exact_solution.h"
36 #include "libmesh/string_to_enum.h"
37 #include "libmesh/enum_solver_package.h"
38 #include "libmesh/enum_norm_type.h"
39 
40 // The systems and solvers we may use
41 #include "curl_curl_system.h"
42 #include "libmesh/diff_solver.h"
43 #include "libmesh/steady_solver.h"
44 #include "solution_function.h"
45 
46 // Bring in everything from the libMesh namespace
47 using namespace libMesh;
48 
49 // Define static data member holding (optional) rotation matrix
51 
52 // The main program.
53 int main (int argc, char ** argv)
54 {
55  // Initialize libMesh.
56  LibMeshInit init (argc, argv);
57 
58  // This example requires a linear solver package.
59  libmesh_example_requires(libMesh::default_solver_package() != INVALID_SOLVER_PACKAGE,
60  "--enable-petsc, --enable-trilinos, or --enable-eigen");
61 
62  // Parse the input file
63  GetPot infile("vector_fe_ex3.in");
64 
65  // But allow the command line to override it.
66  infile.parse_command_line(argc, argv);
67 
68  // Read in parameters from the input file
69  const unsigned int grid_size = infile("grid_size", 2);
70 
71  // Skip higher-dimensional examples on a lower-dimensional libMesh build
72  libmesh_example_requires(2 <= LIBMESH_DIM, "2D support");
73 
74  // Create a mesh, with dimension to be overridden later, on the
75  // default MPI communicator.
76  Mesh mesh(init.comm());
77 
78  // Use the MeshTools::Generation mesh generator to create a uniform
79  // grid on the square [-1,1]^D. We must use at least TRI6 or QUAD8 elements
80  // for the Nedelec triangle or quadrilateral elements, respectively.
81  const std::string elem_str = infile("element_type", std::string("TRI6"));
82 
83  libmesh_error_msg_if(elem_str != "TRI6" && elem_str != "TRI7" && elem_str != "QUAD8" && elem_str != "QUAD9",
84  "You selected: " << elem_str <<
85  " but this example must be run with TRI6, TRI7, QUAD8, or QUAD9.");
86 
88  grid_size,
89  grid_size,
90  -1., 1.,
91  -1., 1.,
92  Utility::string_to_enum<ElemType>(elem_str));
93 
94  // Make sure the code is robust against nodal reorderings.
96 
97 #ifdef LIBMESH_ENABLE_AMR
98  // Make sure the code is robust against mesh refinements.
99  MeshRefinement mesh_refinement(mesh);
100  mesh_refinement.uniformly_refine(infile("refine", 0));
101 #endif
102 
103  // Make sure the code is robust against solves on 2d meshes rotated out of
104  // the xy plane. By default, all Euler angles are zero, the rotation matrix
105  // is the identity, and the mesh stays in place.
106  const Real phi = infile("phi", 0.), theta = infile("theta", 0.), psi = infile("psi", 0.);
108 
109  // Print information about the mesh to the screen.
110  mesh.print_info();
111 
112  // Create an equation systems object.
113  EquationSystems equation_systems (mesh);
114 
115  // Declare the system "CurlCurl" and its variables.
116  CurlCurlSystem & system =
117  equation_systems.add_system<CurlCurlSystem> ("CurlCurl");
118 
119  // Set the FE approximation order.
120  const Order order = static_cast<Order>(infile("order", 1u));
121 
122  libmesh_error_msg_if(order < FIRST || order > FIFTH,
123  "You selected: " << order <<
124  " but this example must be run with 1 <= order <= 5.");
125 
126  system.order(order);
127 
128  // This example only implements the steady-state problem
129  system.time_solver = std::make_unique<SteadySolver>(system);
130 
131  // Initialize the system
132  equation_systems.init();
133 
134  // And the nonlinear solver options
135  DiffSolver & solver = *(system.time_solver->diff_solver().get());
136  solver.quiet = infile("solver_quiet", true);
137  solver.verbose = !solver.quiet;
138  solver.max_nonlinear_iterations = infile("max_nonlinear_iterations", 15);
139  solver.relative_step_tolerance = infile("relative_step_tolerance", 1.e-3);
140  solver.relative_residual_tolerance = infile("relative_residual_tolerance", 1.0e-13);
141  solver.absolute_residual_tolerance = infile("absolute_residual_tolerance", 0.0);
142 
143  // And the linear solver options
144  solver.max_linear_iterations = infile("max_linear_iterations", 50000);
145  solver.initial_linear_tolerance = infile("initial_linear_tolerance", 1.e-10);
146 
147  // Print information about the system to the screen.
148  equation_systems.print_info();
149 
150  // Solve the system.
151  system.solve();
152 
153  ExactSolution exact_sol(equation_systems);
154 
155  SolutionFunction soln_func(system.variable_number("u"));
156  SolutionGradient soln_grad(system.variable_number("u"));
157 
158  // Build FunctionBase* containers to attach to the ExactSolution object.
159  std::vector<FunctionBase<Number> *> sols(1, &soln_func);
160  std::vector<FunctionBase<Gradient> *> grads(1, &soln_grad);
161 
162  exact_sol.attach_exact_values(sols);
163  exact_sol.attach_exact_derivs(grads);
164 
165  // Use higher quadrature order for more accurate error results
166  int extra_error_quadrature = infile("extra_error_quadrature", 2);
167  exact_sol.extra_quadrature_order(extra_error_quadrature);
168 
169  // Compute the error.
170  exact_sol.compute_error("CurlCurl", "u");
171 
172  // Print out the error values
173  libMesh::out << "L2-Error is: "
174  << exact_sol.l2_error("CurlCurl", "u")
175  << std::endl;
176  libMesh::out << "HCurl semi-norm error is: "
177  << exact_sol.error_norm("CurlCurl", "u", HCURL_SEMINORM)
178  << std::endl;
179  libMesh::out << "HCurl-Error is: "
180  << exact_sol.hcurl_error("CurlCurl", "u")
181  << std::endl;
182 
183 #ifdef LIBMESH_HAVE_EXODUS_API
184 
185  // We write the file in the ExodusII format.
186  ExodusII_IO(mesh).write_equation_systems("out.e", equation_systems);
187 
188 #endif // #ifdef LIBMESH_HAVE_EXODUS_API
189 
190  // All done.
191  return 0;
192 }
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
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 ...
int main(int argc, char **argv)
Definition: vector_fe_ex3.C:53
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
Real hcurl_error(std::string_view sys_name, std::string_view unknown_name)
void extra_quadrature_order(const int extraorder)
Increases or decreases the order of the quadrature rule used for numerical integration.
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
static void RM(RealTensor T)
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
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
This is a generic class that defines a solver to handle ImplicitSystem classes, including NonlinearIm...
Definition: diff_solver.h:68
Implements (adaptive) mesh refinement algorithms for a MeshBase.
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)
RealTensorValue rotate(MeshBase &mesh, const Real phi, const Real theta=0., const Real psi=0.)
Rotates the mesh in the xy plane.
FEMSystem, TimeSolver and NewtonSolver will handle most tasks, but we must specify element residuals...
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
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
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
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
Real relative_residual_tolerance
Definition: diff_solver.h:190
This class defines a tensor in LIBMESH_DIM dimensional Real or Complex space.
void order(const Order &order)
void uniformly_refine(unsigned int n=1)
Uniformly refines the mesh n times.
Real error_norm(std::string_view sys_name, std::string_view unknown_name, const FEMNormType &norm)