libMesh
introduction_ex4.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 
19 
20 // <h1>Introduction Example 4 - Solving a 1D, 2D or 3D Poisson Problem in Parallel</h1>
21 // \author Benjamin S. Kirk
22 // \date 2003
23 //
24 // This is the fourth example program. It builds on
25 // the third example program by showing how to formulate
26 // the code in a dimension-independent way. Very minor
27 // changes to the example will allow the problem to be
28 // solved in one, two or three dimensions.
29 //
30 // This example will also introduce the PerfLog class
31 // as a way to monitor your code's performance. We will
32 // use it to instrument the matrix assembly code and look
33 // for bottlenecks where we should focus optimization efforts.
34 //
35 // This example also shows how to extend example 3 to run in
36 // parallel. Notice how little has changed!
37 
38 
39 // C++ include files that we need
40 #include <iostream>
41 #include <algorithm>
42 #include <math.h>
43 #include <set>
44 
45 // Basic include file needed for the mesh functionality.
46 #include "libmesh/libmesh.h"
47 #include "libmesh/mesh.h"
48 #include "libmesh/mesh_generation.h"
49 #include "libmesh/exodusII_io.h"
50 #include "libmesh/gnuplot_io.h"
51 #include "libmesh/linear_implicit_system.h"
52 #include "libmesh/equation_systems.h"
53 
54 // Define the Finite Element object.
55 #include "libmesh/fe.h"
56 
57 // Define Gauss quadrature rules.
58 #include "libmesh/quadrature_gauss.h"
59 
60 // Define the DofMap, which handles degree of freedom
61 // indexing.
62 #include "libmesh/dof_map.h"
63 
64 // Define useful datatypes for finite element
65 // matrix and vector components.
66 #include "libmesh/sparse_matrix.h"
67 #include "libmesh/numeric_vector.h"
68 #include "libmesh/dense_matrix.h"
69 #include "libmesh/dense_vector.h"
70 
71 // Define the PerfLog, a performance logging utility.
72 // It is useful for timing events in a code and giving
73 // you an idea where bottlenecks lie.
74 #include "libmesh/perf_log.h"
75 
76 // The definition of a geometric element
77 #include "libmesh/elem.h"
78 
79 // To impose Dirichlet boundary conditions
80 #include "libmesh/dirichlet_boundaries.h"
81 #include "libmesh/analytic_function.h"
82 
83 #include "libmesh/string_to_enum.h"
84 #include "libmesh/getpot.h"
85 #include "libmesh/enum_solver_package.h"
86 
87 // Bring in everything from the libMesh namespace
88 using namespace libMesh;
89 
90 
91 
92 // Function prototype. This is the function that will assemble
93 // the linear system for our Poisson problem. Note that the
94 // function will take the EquationSystems object and the
95 // name of the system we are assembling as input. From the
96 // EquationSystems object we have access to the Mesh and
97 // other objects we might need.
99  const std::string & system_name);
100 
101 // Exact solution function prototype.
102 Real exact_solution (const Real x,
103  const Real y,
104  const Real z = 0.);
105 
106 // Define a wrapper for exact_solution that will be needed below
108  const Point & p,
109  const Real)
110 {
111  output(0) = exact_solution(p(0),
112  (LIBMESH_DIM>1)?p(1):0,
113  (LIBMESH_DIM>2)?p(2):0);
114 }
115 
116 // Begin the main program.
117 int main (int argc, char ** argv)
118 {
119  // Initialize libMesh and any dependent libraries, like in example 2.
120  LibMeshInit init (argc, argv);
121 
122  // This example requires a linear solver package.
123  libmesh_example_requires(libMesh::default_solver_package() != INVALID_SOLVER_PACKAGE,
124  "--enable-petsc, --enable-trilinos, or --enable-eigen");
125 
126  // Declare a performance log for the main program
127  // PerfLog perf_main("Main Program");
128 
129  // Create a GetPot object to parse the command line
130  GetPot command_line (argc, argv);
131 
132  // Check for proper calling arguments.
133  if (argc < 3)
134  {
135  // This handy function will print the file name, line number,
136  // specified message, and then throw an exception.
137  libmesh_error_msg("Usage:\n" << "\t " << argv[0] << " -d 2(3)" << " -n 15");
138  }
139 
140  // Brief message to the user regarding the program name
141  // and command line arguments.
142  else
143  {
144  libMesh::out << "Running " << argv[0];
145 
146  for (int i=1; i<argc; i++)
147  libMesh::out << " " << argv[i];
148 
149  libMesh::out << std::endl << std::endl;
150  }
151 
152  // Read problem dimension from command line. Use int
153  // instead of unsigned since the GetPot overload is ambiguous
154  // otherwise.
155  int dim = 2;
156  if (command_line.search(1, "-d"))
157  dim = command_line.next(dim);
158 
159  // Skip higher-dimensional examples on a lower-dimensional libMesh build
160  libmesh_example_requires(dim <= LIBMESH_DIM, "2D/3D support");
161 
162  // We use Dirichlet boundary conditions here
163 #ifndef LIBMESH_ENABLE_DIRICHLET
164  libmesh_example_requires(false, "--enable-dirichlet");
165 #endif
166 
167  // Create a mesh with user-defined dimension.
168  // Read number of elements from command line
169  int ps = 15;
170  if (command_line.search(1, "-n"))
171  ps = command_line.next(ps);
172 
173  // Read FE order from command line
174  std::string order = "SECOND";
175  if (command_line.search(2, "-Order", "-o"))
176  order = command_line.next(order);
177 
178  // Read FE Family from command line
179  std::string family = "LAGRANGE";
180  if (command_line.search(2, "-FEFamily", "-f"))
181  family = command_line.next(family);
182 
183  // Cannot use discontinuous basis.
184  if ((family == "MONOMIAL") || (family == "XYZ"))
185  libmesh_error_msg("ex4 currently requires a C^0 (or higher) FE basis.");
186 
187  // Create a mesh, with dimension to be overridden later, distributed
188  // across the default MPI communicator.
189  Mesh mesh(init.comm());
190 
191  // Use the MeshTools::Generation mesh generator to create a uniform
192  // grid on the square [-1,1]^D. We instruct the mesh generator
193  // to build a mesh of 8x8 Quad9 elements in 2D, or Hex27
194  // elements in 3D. Building these higher-order elements allows
195  // us to use higher-order approximation, as in example 3.
196 
197  Real halfwidth = dim > 1 ? 1. : 0.;
198  Real halfheight = dim > 2 ? 1. : 0.;
199 
200  if ((family == "LAGRANGE") && (order == "FIRST"))
201  {
202  // No reason to use high-order geometric elements if we are
203  // solving with low-order finite elements.
205  ps,
206  (dim>1) ? ps : 0,
207  (dim>2) ? ps : 0,
208  -1., 1.,
209  -halfwidth, halfwidth,
210  -halfheight, halfheight,
211  (dim==1) ? EDGE2 :
212  ((dim == 2) ? QUAD4 : HEX8));
213  }
214 
215  else
216  {
218  ps,
219  (dim>1) ? ps : 0,
220  (dim>2) ? ps : 0,
221  -1., 1.,
222  -halfwidth, halfwidth,
223  -halfheight, halfheight,
224  (dim==1) ? EDGE3 :
225  ((dim == 2) ? QUAD9 : HEX27));
226  }
227 
228 
229  // Print information about the mesh to the screen.
230  mesh.print_info();
231 
232 
233  // Create an equation systems object.
234  EquationSystems equation_systems (mesh);
235 
236  // Declare the system and its variables.
237  // Create a system named "Poisson"
238  LinearImplicitSystem & system =
239  equation_systems.add_system<LinearImplicitSystem> ("Poisson");
240 
241 
242  // Add the variable "u" to "Poisson". "u"
243  // will be approximated using second-order approximation.
244  unsigned int u_var = system.add_variable("u",
245  Utility::string_to_enum<Order> (order),
246  Utility::string_to_enum<FEFamily>(family));
247 
248  // Give the system a pointer to the matrix assembly
249  // function.
251 
252  // Construct a Dirichlet boundary condition object
253 
254  // Indicate which boundary IDs we impose the BC on
255  // We either build a line, a square or a cube, and
256  // here we indicate the boundaries IDs in each case
257  std::set<boundary_id_type> boundary_ids;
258  // the dim==1 mesh has two boundaries with IDs 0 and 1
259  boundary_ids.insert(0);
260  boundary_ids.insert(1);
261  // the dim==2 mesh has four boundaries with IDs 0, 1, 2 and 3
262  if (dim>=2)
263  {
264  boundary_ids.insert(2);
265  boundary_ids.insert(3);
266  }
267  // the dim==3 mesh has four boundaries with IDs 0, 1, 2, 3, 4 and 5
268  if (dim==3)
269  {
270  boundary_ids.insert(4);
271  boundary_ids.insert(5);
272  }
273 
274  // Create a vector storing the variable numbers which the BC applies to
275  std::vector<unsigned int> variables(1);
276  variables[0] = u_var;
277 
278  // Create an AnalyticFunction object that we use to project the BC
279  // This function just calls the function exact_solution via exact_solution_wrapper
280  AnalyticFunction<> exact_solution_object(exact_solution_wrapper);
281 
282 #ifdef LIBMESH_ENABLE_DIRICHLET
283  // In general, when reusing a system-indexed exact solution, we want
284  // to use the default system-ordering constructor for
285  // DirichletBoundary, so we demonstrate that here. In this case,
286  // though, we have only one variable, so system- and local-
287  // orderings are the same.
288  DirichletBoundary dirichlet_bc
289  (boundary_ids, variables, exact_solution_object);
290 
291  // We must add the Dirichlet boundary condition _before_
292  // we call equation_systems.init()
293  system.get_dof_map().add_dirichlet_boundary(dirichlet_bc);
294 #endif
295 
296  // Initialize the data structures for the equation system.
297  equation_systems.init();
298 
299  // Print information about the system to the screen.
300  equation_systems.print_info();
301  mesh.print_info();
302 
303  // Solve the system "Poisson", just like example 2.
304  system.solve();
305 
306  // After solving the system write the solution
307  // to a GMV-formatted plot file.
308  if (dim == 1)
309  {
310  GnuPlotIO plot(mesh, "Introduction Example 4, 1D", GnuPlotIO::GRID_ON);
311  plot.write_equation_systems("gnuplot_script", equation_systems);
312  }
313 #ifdef LIBMESH_HAVE_EXODUS_API
314  else
315  {
317  "out_3.e" : "out_2.e", equation_systems);
318  }
319 #endif // #ifdef LIBMESH_HAVE_EXODUS_API
320 
321  // All done.
322  return 0;
323 }
324 
325 
326 
327 
328 // We now define the matrix assembly function for the
329 // Poisson system. We need to first compute element
330 // matrices and right-hand sides, and then take into
331 // account the boundary conditions.
333  const std::string & libmesh_dbg_var(system_name))
334 {
335  // It is a good idea to make sure we are assembling
336  // the proper system.
337  libmesh_assert_equal_to (system_name, "Poisson");
338 
339  // Declare a performance log. Give it a descriptive
340  // string to identify what part of the code we are
341  // logging, since there may be many PerfLogs in an
342  // application.
343  PerfLog perf_log ("Matrix Assembly");
344 
345  // Get a constant reference to the mesh object.
346  const MeshBase & mesh = es.get_mesh();
347 
348  // The dimension that we are running
349  const unsigned int dim = mesh.mesh_dimension();
350 
351  // Get a reference to the LinearImplicitSystem we are solving
352  LinearImplicitSystem & system = es.get_system<LinearImplicitSystem>("Poisson");
353 
354  // A reference to the DofMap object for this system. The DofMap
355  // object handles the index translation from node and element numbers
356  // to degree of freedom numbers. We will talk more about the DofMap
357  // in future examples.
358  const DofMap & dof_map = system.get_dof_map();
359 
360  // Get a constant reference to the Finite Element type
361  // for the first (and only) variable in the system.
362  FEType fe_type = dof_map.variable_type(0);
363 
364  // Build a Finite Element object of the specified type. Since the
365  // FEBase::build() member dynamically creates memory we will
366  // store the object as a std::unique_ptr<FEBase>. This can be thought
367  // of as a pointer that will clean up after itself.
368  std::unique_ptr<FEBase> fe (FEBase::build(dim, fe_type));
369 
370  // A 5th order Gauss quadrature rule for numerical integration.
371  QGauss qrule (dim, FIFTH);
372 
373  // Tell the finite element object to use our quadrature rule.
374  fe->attach_quadrature_rule (&qrule);
375 
376  // Declare a special finite element object for
377  // boundary integration.
378  std::unique_ptr<FEBase> fe_face (FEBase::build(dim, fe_type));
379 
380  // Boundary integration requires one quadrature rule,
381  // with dimensionality one less than the dimensionality
382  // of the element.
383  QGauss qface(dim-1, FIFTH);
384 
385  // Tell the finite element object to use our
386  // quadrature rule.
387  fe_face->attach_quadrature_rule (&qface);
388 
389  // Here we define some references to cell-specific data that
390  // will be used to assemble the linear system.
391  // We begin with the element Jacobian * quadrature weight at each
392  // integration point.
393  const std::vector<Real> & JxW = fe->get_JxW();
394 
395  // The physical XY locations of the quadrature points on the element.
396  // These might be useful for evaluating spatially varying material
397  // properties at the quadrature points.
398  const std::vector<Point> & q_point = fe->get_xyz();
399 
400  // The element shape functions evaluated at the quadrature points.
401  const std::vector<std::vector<Real>> & phi = fe->get_phi();
402 
403  // The element shape function gradients evaluated at the quadrature
404  // points.
405  const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
406 
407  // Define data structures to contain the element matrix
408  // and right-hand-side vector contribution. Following
409  // basic finite element terminology we will denote these
410  // "Ke" and "Fe". More detail is in example 3.
413 
414  // This vector will hold the degree of freedom indices for
415  // the element. These define where in the global system
416  // the element degrees of freedom get mapped.
417  std::vector<dof_id_type> dof_indices;
418 
419  // Now we will loop over all the elements in the mesh.
420  // We will compute the element matrix and right-hand-side
421  // contribution. See example 3 for a discussion of the
422  // element iterators.
423  for (const auto & elem : mesh.active_local_element_ptr_range())
424  {
425  // Start logging the shape function initialization.
426  // This is done through a simple function call with
427  // the name of the event to log.
428  perf_log.push("elem init");
429 
430  // Get the degree of freedom indices for the
431  // current element. These define where in the global
432  // matrix and right-hand-side this element will
433  // contribute to.
434  dof_map.dof_indices (elem, dof_indices);
435 
436  // Cache the number of degrees of freedom on this element, for
437  // use as a loop bound later. We use cast_int to explicitly
438  // convert from size() (which may be 64-bit) to unsigned int
439  // (which may be 32-bit but which is definitely enough to count
440  // *local* degrees of freedom.
441  const unsigned int n_dofs =
442  cast_int<unsigned int>(dof_indices.size());
443 
444  // Compute the element-specific data for the current
445  // element. This involves computing the location of the
446  // quadrature points (q_point) and the shape functions
447  // (phi, dphi) for the current element.
448  fe->reinit (elem);
449 
450  // With one variable, we should have the same number of degrees
451  // of freedom as shape functions.
452  libmesh_assert_equal_to (n_dofs, phi.size());
453 
454  // Zero the element matrix and right-hand side before
455  // summing them. We use the resize member here because
456  // the number of degrees of freedom might have changed from
457  // the last element. Note that this will be the case if the
458  // element type is different (i.e. the last element was a
459  // triangle, now we are on a quadrilateral).
460  Ke.resize (n_dofs, n_dofs);
461 
462  Fe.resize (n_dofs);
463 
464  // Stop logging the shape function initialization.
465  // If you forget to stop logging an event the PerfLog
466  // object will probably catch the error and abort.
467  perf_log.pop("elem init");
468 
469  // Now we will build the element matrix. This involves
470  // a double loop to integrate the test functions (i) against
471  // the trial functions (j).
472  //
473  // We have split the numeric integration into two loops
474  // so that we can log the matrix and right-hand-side
475  // computation separately.
476  //
477  // Now start logging the element matrix computation
478  perf_log.push ("Ke");
479 
480  for (unsigned int qp=0; qp<qrule.n_points(); qp++)
481  for (unsigned int i=0; i != n_dofs; i++)
482  for (unsigned int j=0; j != n_dofs; j++)
483  Ke(i,j) += JxW[qp]*(dphi[i][qp]*dphi[j][qp]);
484 
485 
486  // Stop logging the matrix computation
487  perf_log.pop ("Ke");
488 
489  // Now we build the element right-hand-side contribution.
490  // This involves a single loop in which we integrate the
491  // "forcing function" in the PDE against the test functions.
492  //
493  // Start logging the right-hand-side computation
494  perf_log.push ("Fe");
495 
496  for (unsigned int qp=0; qp<qrule.n_points(); qp++)
497  {
498  // fxy is the forcing function for the Poisson equation.
499  // In this case we set fxy to be a finite difference
500  // Laplacian approximation to the (known) exact solution.
501  //
502  // We will use the second-order accurate FD Laplacian
503  // approximation, which in 2D on a structured grid is
504  //
505  // u_xx + u_yy = (u(i-1,j) + u(i+1,j) +
506  // u(i,j-1) + u(i,j+1) +
507  // -4*u(i,j))/h^2
508  //
509  // Since the value of the forcing function depends only
510  // on the location of the quadrature point (q_point[qp])
511  // we will compute it here, outside of the i-loop
512  const Real x = q_point[qp](0);
513 #if LIBMESH_DIM > 1
514  const Real y = q_point[qp](1);
515 #else
516  const Real y = 0.;
517 #endif
518 #if LIBMESH_DIM > 2
519  const Real z = q_point[qp](2);
520 #else
521  const Real z = 0.;
522 #endif
523  const Real eps = 1.e-3;
524 
525  const Real uxx = (exact_solution(x-eps, y, z) +
526  exact_solution(x+eps, y, z) +
527  -2.*exact_solution(x, y, z))/eps/eps;
528 
529  const Real uyy = (exact_solution(x, y-eps, z) +
530  exact_solution(x, y+eps, z) +
531  -2.*exact_solution(x, y, z))/eps/eps;
532 
533  const Real uzz = (exact_solution(x, y, z-eps) +
534  exact_solution(x, y, z+eps) +
535  -2.*exact_solution(x, y, z))/eps/eps;
536 
537  Real fxy;
538  if (dim==1)
539  {
540  // In 1D, compute the rhs by differentiating the
541  // exact solution twice.
542  const Real pi = libMesh::pi;
543  fxy = (0.25*pi*pi)*sin(.5*pi*x);
544  }
545  else
546  {
547  fxy = - (uxx + uyy + ((dim==2) ? 0. : uzz));
548  }
549 
550  // Add the RHS contribution
551  for (unsigned int i=0; i != n_dofs; i++)
552  Fe(i) += JxW[qp]*fxy*phi[i][qp];
553  }
554 
555  // Stop logging the right-hand-side computation
556  perf_log.pop ("Fe");
557 
558  // If this assembly program were to be used on an adaptive mesh,
559  // we would have to apply any hanging node constraint equations
560  // Also, note that here we call heterogenously_constrain_element_matrix_and_vector
561  // to impose a inhomogeneous Dirichlet boundary conditions.
562  dof_map.heterogenously_constrain_element_matrix_and_vector (Ke, Fe, dof_indices);
563 
564  // The element matrix and right-hand-side are now built
565  // for this element. Add them to the global matrix and
566  // right-hand-side vector. The SparseMatrix::add_matrix()
567  // and NumericVector::add_vector() members do this for us.
568  // Start logging the insertion of the local (element)
569  // matrix and vector into the global matrix and vector
570  LOG_SCOPE_WITH("matrix insertion", "", perf_log);
571 
572  system.matrix->add_matrix (Ke, dof_indices);
573  system.rhs->add_vector (Fe, dof_indices);
574  }
575 
576  // That's it. We don't need to do anything else to the
577  // PerfLog. When it goes out of scope (at this function return)
578  // it will print its log to the screen. Pretty easy, huh?
579 }
libMesh::pi
const Real pi
.
Definition: libmesh.h:237
libMesh::Mesh
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
Definition: mesh.h:50
libMesh::EquationSystems::add_system
virtual System & add_system(const std::string &system_type, const std::string &name)
Add the system of type system_type named name to the systems array.
Definition: equation_systems.C:345
libMesh::EquationSystems::get_mesh
const MeshBase & get_mesh() const
Definition: equation_systems.h:637
libMesh::ExplicitSystem::rhs
NumericVector< Number > * rhs
The system matrix.
Definition: explicit_system.h:114
libMesh::QGauss
This class implements specific orders of Gauss quadrature.
Definition: quadrature_gauss.h:39
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::HEX8
Definition: enum_elem_type.h:47
libMesh::MeshBase::active_local_element_ptr_range
virtual SimpleRange< element_iterator > active_local_element_ptr_range()=0
libMesh::QBase::n_points
unsigned int n_points() const
Definition: quadrature.h:126
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::EquationSystems::get_system
const T_sys & get_system(const std::string &name) const
Definition: equation_systems.h:757
libMesh::MeshTools::Generation::build_cube
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.
Definition: mesh_generation.C:298
libMesh::DenseMatrix< Number >
libMesh::LinearImplicitSystem::solve
virtual void solve() override
Assembles & solves the linear system A*x=b.
Definition: linear_implicit_system.C:108
libMesh::default_solver_package
SolverPackage default_solver_package()
Definition: libmesh.C:993
mesh
MeshBase & mesh
Definition: mesh_communication.C:1257
libMesh::MeshBase::mesh_dimension
unsigned int mesh_dimension() const
Definition: mesh_base.C:135
libMesh::FIFTH
Definition: enum_order.h:46
libMesh::ExodusII_IO
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
Definition: exodusII_io.h:51
exact_solution
Real exact_solution(const Real x, const Real y, const Real z=0.)
This is the exact solution that we are trying to obtain.
Definition: exact_solution.C:43
dim
unsigned int dim
Definition: adaptivity_ex3.C:113
libMesh::DenseMatrix::resize
void resize(const unsigned int new_m, const unsigned int new_n)
Resize the matrix.
Definition: dense_matrix.h:822
libMesh::TriangleWrapper::init
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
libMesh::System::attach_assemble_function
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:1755
libMesh::NumericVector::add_vector
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].
Definition: numeric_vector.C:363
libMesh::MeshBase
This is the MeshBase class.
Definition: mesh_base.h:78
libMesh::HEX27
Definition: enum_elem_type.h:49
libMesh::PerfLog::pop
void pop(const char *label, const char *header="")
Pop the event label off the stack, resuming any lower event.
Definition: perf_log.C:163
libMesh::INVALID_SOLVER_PACKAGE
Definition: enum_solver_package.h:43
main
int main(int argc, char **argv)
Definition: introduction_ex4.C:117
libMesh::FEGenericBase::build
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
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
libMesh::AnalyticFunction
Wraps a function pointer into a FunctionBase object.
Definition: analytic_function.h:48
libMesh::QUAD4
Definition: enum_elem_type.h:41
libMesh::PerfLog
The PerfLog class allows monitoring of specific events.
Definition: perf_log.h:125
libMesh::EquationSystems::init
virtual void init()
Initialize all the systems.
Definition: equation_systems.C:96
libMesh::Point
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:38
libMesh::ImplicitSystem::matrix
SparseMatrix< Number > * matrix
The system matrix.
Definition: implicit_system.h:393
assemble_poisson
void assemble_poisson(EquationSystems &es, const std::string &system_name)
libMesh::LibMeshInit
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
Definition: libmesh.h:83
exact_solution_wrapper
void exact_solution_wrapper(DenseVector< Number > &output, const Point &p, const Real)
Definition: introduction_ex4.C:107
libMesh::EquationSystems::print_info
void print_info(std::ostream &os=libMesh::out) const
Prints information about the equation systems, by default to libMesh::out.
Definition: equation_systems.C:1314
libMesh::EquationSystems
This is the EquationSystems class.
Definition: equation_systems.h:74
libMesh::MeshOutput::write_equation_systems
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
libMesh::DenseVector::resize
void resize(const unsigned int n)
Resize the vector.
Definition: dense_vector.h:355
libMesh::SparseMatrix::add_matrix
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.
libMesh::FEType
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:178
libMesh::DofMap
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:176
libMesh::MeshBase::print_info
void print_info(std::ostream &os=libMesh::out) const
Prints relevant information about the mesh.
Definition: mesh_base.C:585
libMesh::EDGE3
Definition: enum_elem_type.h:36
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::PerfLog::push
void push(const char *label, const char *header="")
Push the event label onto the stack, pausing any active event.
Definition: perf_log.C:133
libMesh::QUAD9
Definition: enum_elem_type.h:43
libMesh::System::get_dof_map
const DofMap & get_dof_map() const
Definition: system.h:2099
libMesh::GnuPlotIO
This class implements writing meshes using GNUplot, designed for use only with 1D meshes.
Definition: gnuplot_io.h:43
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::LinearImplicitSystem
Manages consistently variables, degrees of freedom, coefficient vectors, matrices and linear solvers ...
Definition: linear_implicit_system.h:55
libMesh::out
OStreamProxy out
libMesh::DenseVector< Number >
libMesh::EDGE2
Definition: enum_elem_type.h:35
libMesh::GnuPlotIO::GRID_ON
Definition: gnuplot_io.h:50