libMesh
Functions
subdomains_ex2.C File Reference

Go to the source code of this file.

Functions

void assemble_poisson (EquationSystems &es, const std::string &system_name)
 
Real exact_solution (const Real x, const Real y=0., const Real z=0.)
 This is the exact solution that we are trying to obtain. More...
 
int main (int argc, char **argv)
 
void assemble_poisson (EquationSystems &es, const std::string &libmesh_dbg_var(system_name))
 

Function Documentation

◆ assemble_poisson() [1/2]

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

Definition at line 288 of file subdomains_ex2.C.

290 {
291  // It is a good idea to make sure we are assembling
292  // the proper system.
293  libmesh_assert_equal_to (system_name, "Poisson");
294 
295  // Declare a performance log. Give it a descriptive
296  // string to identify what part of the code we are
297  // logging, since there may be many PerfLogs in an
298  // application.
299  PerfLog perf_log ("Matrix Assembly");
300 
301  // Get a constant reference to the mesh object.
302  const MeshBase & mesh = es.get_mesh();
303 
304  // The dimension that we are running
305  const unsigned int dim = mesh.mesh_dimension();
306 
307  // Get a reference to the LinearImplicitSystem we are solving
308  LinearImplicitSystem & system = es.get_system<LinearImplicitSystem>("Poisson");
309 
310  // A reference to the DofMap object for this system. The DofMap
311  // object handles the index translation from node and element numbers
312  // to degree of freedom numbers. We will talk more about the DofMap
313  // in future examples.
314  const DofMap & dof_map = system.get_dof_map();
315 
316  // Get a constant reference to the Finite Element type
317  // for the first (and only) variable in the system.
318  FEType fe_type = dof_map.variable_type(0);
319 
320  // Build a Finite Element object of the specified type. Since the
321  // FEBase::build() member dynamically creates memory we will
322  // store the object as a std::unique_ptr<FEBase>. This can be thought
323  // of as a pointer that will clean up after itself.
324  std::unique_ptr<FEBase> fe (FEBase::build(dim, fe_type));
325 
326  // A 5th order Gauss quadrature rule for numerical integration.
327  QGauss qrule (dim, FIFTH);
328 
329  // Tell the finite element object to use our quadrature rule.
330  fe->attach_quadrature_rule (&qrule);
331 
332  // Declare a special finite element object for
333  // boundary integration.
334  std::unique_ptr<FEBase> fe_face (FEBase::build(dim, fe_type));
335 
336  // Boundary integration requires one quadrature rule,
337  // with dimensionality one less than the dimensionality
338  // of the element.
339  QGauss qface(dim-1, FIFTH);
340 
341  // Tell the finite element object to use our
342  // quadrature rule.
343  fe_face->attach_quadrature_rule (&qface);
344 
345  // Here we define some references to cell-specific data that
346  // will be used to assemble the linear system.
347  // We begin with the element Jacobian * quadrature weight at each
348  // integration point.
349  const std::vector<Real> & JxW = fe->get_JxW();
350 
351  // The physical XY locations of the quadrature points on the element.
352  // These might be useful for evaluating spatially varying material
353  // properties at the quadrature points.
354  const std::vector<Point> & q_point = fe->get_xyz();
355 
356  // The element shape functions evaluated at the quadrature points.
357  const std::vector<std::vector<Real>> & phi = fe->get_phi();
358 
359  // The element shape function gradients evaluated at the quadrature
360  // points.
361  const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
362 
363  // Define data structures to contain the element matrix
364  // and right-hand-side vector contribution. Following
365  // basic finite element terminology we will denote these
366  // "Ke" and "Fe". More detail is in example 3.
369 
370  // This vector will hold the degree of freedom indices for
371  // the element. These define where in the global system
372  // the element degrees of freedom get mapped.
373  std::vector<dof_id_type> dof_indices, dof_indices2;
374 
375  // Now we will loop over all the "local" elements in the mesh. We
376  // will compute the element matrix and right-hand-side contribution.
377  // See example 3 for a discussion of the element iterators. Here we
378  // only want to loop over elements that are owned by the local
379  // processor. This allows each processor to compute its components
380  // of the global matrix.
381  //
382  // "PARALLEL CHANGE"
383  for (const auto & elem : as_range(mesh.local_elements_begin(),
385  {
386  // Start logging the shape function initialization.
387  // This is done through a simple function call with
388  // the name of the event to log.
389  perf_log.push("elem init");
390 
391  // Get the degree of freedom indices for the
392  // current element. These define where in the global
393  // matrix and right-hand-side this element will
394  // contribute to.
395  dof_map.dof_indices (elem, dof_indices, 0);
396  dof_map.dof_indices (elem, dof_indices2, 1);
397 
398  // libMesh::out << "dof_indices.size()="
399  // << dof_indices.size()
400  // << ", dof_indices2.size()="
401  // << dof_indices2.size()
402  // << std::endl;
403 
404  // Compute the element-specific data for the current
405  // element. This involves computing the location of the
406  // quadrature points (q_point) and the shape functions
407  // (phi, dphi) for the current element.
408  fe->reinit (elem);
409 
410  // Zero the element matrix and right-hand side before
411  // summing them. We use the resize member here because
412  // the number of degrees of freedom might have changed from
413  // the last element. Note that this will be the case if the
414  // element type is different (i.e. the last element was a
415  // triangle, now we are on a quadrilateral).
416  Ke.resize (std::max(dof_indices.size(), dof_indices2.size()),
417  std::max(dof_indices.size(), dof_indices2.size()));
418 
419  Fe.resize (std::max(dof_indices.size(), dof_indices2.size()));
420 
421  // Stop logging the shape function initialization.
422  // If you forget to stop logging an event the PerfLog
423  // object will probably catch the error and abort.
424  perf_log.pop("elem init");
425 
426  // Now we will build the element matrix. This involves
427  // a double loop to integrate the test functions (i) against
428  // the trial functions (j).
429  //
430  // We have split the numeric integration into two loops
431  // so that we can log the matrix and right-hand-side
432  // computation separately.
433  //
434  // Now start logging the element matrix computation
435  perf_log.push ("Ke");
436 
437  for (unsigned int qp=0; qp<qrule.n_points(); qp++)
438  for (std::size_t i=0; i<phi.size(); i++)
439  for (std::size_t j=0; j<phi.size(); j++)
440  Ke(i,j) += JxW[qp]*(dphi[i][qp]*dphi[j][qp]);
441 
442 
443  // Stop logging the matrix computation
444  perf_log.pop ("Ke");
445 
446  // Now we build the element right-hand-side contribution.
447  // This involves a single loop in which we integrate the
448  // "forcing function" in the PDE against the test functions.
449  //
450  // Start logging the right-hand-side computation
451  perf_log.push ("Fe");
452 
453  for (unsigned int qp=0; qp<qrule.n_points(); qp++)
454  {
455  // fxy is the forcing function for the Poisson equation.
456  // In this case we set fxy to be a finite difference
457  // Laplacian approximation to the (known) exact solution.
458  //
459  // We will use the second-order accurate FD Laplacian
460  // approximation, which in 2D on a structured grid is
461  //
462  // u_xx + u_yy = (u(i-1,j) + u(i+1,j) +
463  // u(i,j-1) + u(i,j+1) +
464  // -4*u(i,j))/h^2
465  //
466  // Since the value of the forcing function depends only
467  // on the location of the quadrature point (q_point[qp])
468  // we will compute it here, outside of the i-loop
469  const Real x = q_point[qp](0);
470 #if LIBMESH_DIM > 1
471  const Real y = q_point[qp](1);
472 #else
473  const Real y = 0;
474 #endif
475 #if LIBMESH_DIM > 2
476  const Real z = q_point[qp](2);
477 #else
478  const Real z = 0;
479 #endif
480  const Real eps = 1.e-3;
481 
482  const Real uxx = (exact_solution(x-eps, y, z) +
483  exact_solution(x+eps, y, z) +
484  -2.*exact_solution(x, y, z))/eps/eps;
485 
486  const Real uyy = (exact_solution(x, y-eps, z) +
487  exact_solution(x, y+eps, z) +
488  -2.*exact_solution(x, y, z))/eps/eps;
489 
490  const Real uzz = (exact_solution(x, y, z-eps) +
491  exact_solution(x, y, z+eps) +
492  -2.*exact_solution(x, y, z))/eps/eps;
493 
494  Real fxy;
495  if (dim==1)
496  {
497  // In 1D, compute the rhs by differentiating the
498  // exact solution twice.
499  const Real pi = libMesh::pi;
500  fxy = (0.25*pi*pi)*sin(.5*pi*x);
501  }
502  else
503  {
504  fxy = - (uxx + uyy + ((dim==2) ? 0. : uzz));
505  }
506 
507  // Add the RHS contribution
508  for (std::size_t i=0; i<phi.size(); i++)
509  Fe(i) += JxW[qp]*fxy*phi[i][qp];
510  }
511 
512  // Stop logging the right-hand-side computation
513  perf_log.pop ("Fe");
514 
515  // At this point the interior element integration has
516  // been completed. However, we have not yet addressed
517  // boundary conditions. For this example we will only
518  // consider simple Dirichlet boundary conditions imposed
519  // via the penalty method. This is discussed at length in
520  // example 3.
521  {
522  // Start logging the boundary condition computation. We use a
523  // macro to log everything in this scope.
524  LOG_SCOPE_WITH("BCs", "", perf_log);
525 
526  // The following loops over the sides of the element.
527  // If the element has no neighbor on a side then that
528  // side MUST live on a boundary of the domain.
529  for (auto side : elem->side_index_range())
530  if ((elem->neighbor_ptr(side) == nullptr) ||
531  (elem->neighbor_ptr(side)->subdomain_id() != elem->subdomain_id()))
532  {
533 
534  // The penalty value. \frac{1}{\epsilon}
535  // in the discussion above.
536  const Real penalty = 1.e10;
537 
538  // The value of the shape functions at the quadrature
539  // points.
540  const std::vector<std::vector<Real>> & phi_face = fe_face->get_phi();
541 
542  // The Jacobian * Quadrature Weight at the quadrature
543  // points on the face.
544  const std::vector<Real> & JxW_face = fe_face->get_JxW();
545 
546  // The XYZ locations (in physical space) of the
547  // quadrature points on the face. This is where
548  // we will interpolate the boundary value function.
549  const std::vector<Point> & qface_point = fe_face->get_xyz();
550 
551  // Compute the shape function values on the element
552  // face.
553  fe_face->reinit(elem, side);
554 
555  // Loop over the face quadrature points for integration.
556  for (unsigned int qp=0; qp<qface.n_points(); qp++)
557  {
558  // The location on the boundary of the current
559  // face quadrature point.
560  const Real xf = qface_point[qp](0);
561 #if LIBMESH_DIM > 1
562  const Real yf = qface_point[qp](1);
563 #else
564  const Real yf = 0.;
565 #endif
566 #if LIBMESH_DIM > 2
567  const Real zf = qface_point[qp](2);
568 #else
569  const Real zf = 0.;
570 #endif
571 
572 
573  // The boundary value.
574  const Real value = exact_solution(xf, yf, zf);
575 
576  // Matrix contribution of the L2 projection.
577  for (std::size_t i=0; i<phi_face.size(); i++)
578  for (std::size_t j=0; j<phi_face.size(); j++)
579  Ke(i,j) += JxW_face[qp]*penalty*phi_face[i][qp]*phi_face[j][qp];
580 
581  // Right-hand-side contribution of the L2
582  // projection.
583  for (std::size_t i=0; i<phi_face.size(); i++)
584  Fe(i) += JxW_face[qp]*penalty*value*phi_face[i][qp];
585  }
586  }
587  }
588 
589 
590  // The element matrix and right-hand-side are now built
591  // for this element. Add them to the global matrix and
592  // right-hand-side vector. The PetscMatrix::add_matrix()
593  // and PetscVector::add_vector() members do this for us.
594  // Start logging the insertion of the local (element)
595  // matrix and vector into the global matrix and vector
596  LOG_SCOPE_WITH("matrix insertion", "", perf_log);
597 
598  if (dof_indices.size())
599  {
600  system.matrix->add_matrix (Ke, dof_indices);
601  system.rhs->add_vector (Fe, dof_indices);
602  }
603 
604  if (dof_indices2.size())
605  {
606  system.matrix->add_matrix (Ke, dof_indices2);
607  system.rhs->add_vector (Fe, dof_indices2);
608  }
609  }
610 
611  // That's it. We don't need to do anything else to the
612  // PerfLog. When it goes out of scope (at this function return)
613  // it will print its log to the screen. Pretty easy, huh?
614 }

References libMesh::SparseMatrix< T >::add_matrix(), libMesh::NumericVector< T >::add_vector(), libMesh::as_range(), libMesh::FEGenericBase< OutputType >::build(), dim, exact_solution(), libMesh::FIFTH, libMesh::System::get_dof_map(), libMesh::EquationSystems::get_mesh(), libMesh::EquationSystems::get_system(), libMesh::MeshBase::local_elements_begin(), libMesh::MeshBase::local_elements_end(), libMesh::ImplicitSystem::matrix, mesh, libMesh::MeshBase::mesh_dimension(), libMesh::QBase::n_points(), libMesh::pi, libMesh::PerfLog::pop(), libMesh::PerfLog::push(), libMesh::Real, libMesh::DenseVector< T >::resize(), libMesh::DenseMatrix< T >::resize(), libMesh::ExplicitSystem::rhs, and value.

◆ assemble_poisson() [2/2]

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

Referenced by main().

◆ exact_solution()

Real exact_solution ( const Real  x,
const Real  y,
const Real  t 
)

This is the exact solution that we are trying to obtain.

We will solve

  • (u_xx + u_yy) = f

and take a finite difference approximation using this function to get f. This is the well-known "method of manufactured solutions".

Definition at line 43 of file exact_solution.C.

46 {
47  static const Real pi = acos(-1.);
48 
49  return cos(.5*pi*x)*sin(.5*pi*y)*cos(.5*pi*z);
50 }

Referenced by assemble_poisson().

◆ main()

int main ( int  argc,
char **  argv 
)

Definition at line 92 of file subdomains_ex2.C.

93 {
94  // Initialize libMesh and any dependent libraries, like in example 2.
95  LibMeshInit init (argc, argv);
96 
97  // This example requires a linear solver package.
98  libmesh_example_requires(libMesh::default_solver_package() != INVALID_SOLVER_PACKAGE,
99  "--enable-petsc, --enable-trilinos, or --enable-eigen");
100 
101  // Declare a performance log for the main program
102  // PerfLog perf_main("Main Program");
103 
104  // Create a GetPot object to parse the command line
105  GetPot command_line (argc, argv);
106 
107  // Check for proper calling arguments.
108  if (argc < 3)
109  {
110  // This handy function will print the file name, line number,
111  // specified message, and then throw an exception.
112  libmesh_error_msg("Usage:\n" << "\t " << argv[0] << " -d 2(3)" << " -n 15");
113  }
114 
115  // Brief message to the user regarding the program name
116  // and command line arguments.
117  else
118  {
119  libMesh::out << "Running " << argv[0];
120 
121  for (int i=1; i<argc; i++)
122  libMesh::out << " " << argv[i];
123 
124  libMesh::out << std::endl << std::endl;
125  }
126 
127 
128  // Read problem dimension from command line. Use int
129  // instead of unsigned since the GetPot overload is ambiguous
130  // otherwise.
131  int dim = 2;
132  if (command_line.search(1, "-d"))
133  dim = command_line.next(dim);
134 
135  // Skip higher-dimensional examples on a lower-dimensional libMesh build
136  libmesh_example_requires(dim <= LIBMESH_DIM, "2D/3D support");
137 
138  // Create a mesh with user-defined dimension on the default MPI
139  // communicator.
140  Mesh mesh (init.comm(), dim);
141 
142  // Read number of elements from command line
143  int ps = 15;
144  if (command_line.search(1, "-n"))
145  ps = command_line.next(ps);
146 
147  // Read FE order from command line
148  std::string order = "SECOND";
149  if (command_line.search(2, "-Order", "-o"))
150  order = command_line.next(order);
151 
152  // Read FE Family from command line
153  std::string family = "LAGRANGE";
154  if (command_line.search(2, "-FEFamily", "-f"))
155  family = command_line.next(family);
156 
157  // Cannot use discontinuous basis.
158  if ((family == "MONOMIAL") || (family == "XYZ"))
159  {
160  if (mesh.processor_id() == 0)
161  libmesh_error_msg("This example requires a C^0 (or higher) FE basis.");
162  }
163 
164  // Use the MeshTools::Generation mesh generator to create a uniform
165  // grid on the square [-1,1]^D. We instruct the mesh generator
166  // to build a mesh of 8x8 Quad9 elements in 2D, or Hex27
167  // elements in 3D. Building these higher-order elements allows
168  // us to use higher-order approximation, as in example 3.
169 
170  Real halfwidth = dim > 1 ? 1. : 0.;
171  Real halfheight = dim > 2 ? 1. : 0.;
172 
173  if ((family == "LAGRANGE") && (order == "FIRST"))
174  {
175  // No reason to use high-order geometric elements if we are
176  // solving with low-order finite elements.
178  ps,
179  (dim>1) ? ps : 0,
180  (dim>2) ? ps : 0,
181  -1., 1.,
182  -halfwidth, halfwidth,
183  -halfheight, halfheight,
184  (dim==1) ? EDGE2 :
185  ((dim == 2) ? QUAD4 : HEX8));
186  }
187 
188  else
189  {
191  ps,
192  (dim>1) ? ps : 0,
193  (dim>2) ? ps : 0,
194  -1., 1.,
195  -halfwidth, halfwidth,
196  -halfheight, halfheight,
197  (dim==1) ? EDGE3 :
198  ((dim == 2) ? QUAD9 : HEX27));
199  }
200 
201  for (auto & elem : mesh.element_ptr_range())
202  {
203  const Point cent = elem->centroid();
204  if (dim > 1)
205  {
206  if ((cent(0) > 0) == (cent(1) > 0))
207  elem->subdomain_id() = 1;
208  }
209  else if (cent(0) > 0)
210  elem->subdomain_id() = 1;
211  }
212 
213  // Print information about the mesh to the screen.
214  mesh.print_info();
215 
216  // Create an equation systems object.
217  EquationSystems equation_systems (mesh);
218 
219  // Declare the system and its variables.
220  // Create a system named "Poisson"
221  LinearImplicitSystem & system =
222  equation_systems.add_system<LinearImplicitSystem> ("Poisson");
223 
224 
225  std::set<subdomain_id_type> active_subdomains;
226 
227 
228  // Add the variable "u" to "Poisson". "u"
229  // will be approximated using second-order approximation.
230  active_subdomains.clear(); active_subdomains.insert(0);
231  system.add_variable("u",
232  Utility::string_to_enum<Order> (order),
233  Utility::string_to_enum<FEFamily>(family),
234  &active_subdomains);
235 
236  // Add the variable "v" to "Poisson". "v"
237  // will be approximated using second-order approximation.
238  active_subdomains.clear(); active_subdomains.insert(1);
239  system.add_variable("v",
240  Utility::string_to_enum<Order> (order),
241  Utility::string_to_enum<FEFamily>(family),
242  &active_subdomains);
243 
244  // Give the system a pointer to the matrix assembly
245  // function.
247 
248  // Initialize the data structures for the equation system.
249  equation_systems.init();
250 
251  // Print information about the system to the screen.
252  equation_systems.print_info();
253  mesh.print_info();
254 
255  // Solve the system "Poisson", just like example 2.
256  equation_systems.get_system("Poisson").solve();
257 
258  // After solving the system write the solution
259  // to a GMV-formatted plot file.
260  if (dim == 1)
261  {
262  GnuPlotIO plot(mesh, "Subdomains Example 2, 1D", GnuPlotIO::GRID_ON);
263  plot.write_equation_systems("gnuplot_script", equation_systems);
264  }
265  else
266  {
267 #ifdef LIBMESH_HAVE_EXODUS_API
269  "out_3.e" : "out_2.e", equation_systems);
270 #endif // #ifdef LIBMESH_HAVE_EXODUS_API
271  }
272 
273  // All done.
274  return 0;
275 }

References libMesh::EquationSystems::add_system(), libMesh::System::add_variable(), assemble_poisson(), libMesh::System::attach_assemble_function(), libMesh::MeshTools::Generation::build_cube(), libMesh::LinearImplicitSystem::clear(), libMesh::default_solver_package(), dim, libMesh::EDGE2, libMesh::EDGE3, libMesh::MeshBase::element_ptr_range(), libMesh::EquationSystems::get_system(), libMesh::GnuPlotIO::GRID_ON, libMesh::HEX27, libMesh::HEX8, libMesh::TriangleWrapper::init(), libMesh::EquationSystems::init(), libMesh::INVALID_SOLVER_PACKAGE, mesh, libMesh::out, libMesh::EquationSystems::print_info(), libMesh::MeshBase::print_info(), libMesh::ParallelObject::processor_id(), libMesh::QUAD4, libMesh::QUAD9, libMesh::Real, and libMesh::MeshOutput< MT >::write_equation_systems().

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::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::HEX8
Definition: enum_elem_type.h:47
libMesh::LinearImplicitSystem::clear
virtual void clear() override
Clear all the data structures associated with the system.
Definition: linear_implicit_system.C:61
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::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
libMesh::MeshBase::local_elements_begin
virtual element_iterator local_elements_begin()=0
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::MeshBase::element_ptr_range
virtual SimpleRange< element_iterator > element_ptr_range()=0
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::INVALID_SOLVER_PACKAGE
Definition: enum_solver_package.h:43
libMesh::ParallelObject::processor_id
processor_id_type processor_id() const
Definition: parallel_object.h:106
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::QUAD4
Definition: enum_elem_type.h:41
libMesh::PerfLog
The PerfLog class allows monitoring of specific events.
Definition: perf_log.h:125
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
libMesh::as_range
SimpleRange< IndexType > as_range(const std::pair< IndexType, IndexType > &p)
Helper function that allows us to treat a homogenous pair as a range.
Definition: simple_range.h:57
libMesh::LibMeshInit
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
Definition: libmesh.h:83
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
assemble_poisson
void assemble_poisson(EquationSystems &es, const std::string &system_name)
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::MeshBase::local_elements_end
virtual element_iterator local_elements_end()=0
value
static const bool value
Definition: xdr_io.C:56
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
exact_solution
Real exact_solution(const Real x, const Real y=0., const Real z=0.)
This is the exact solution that we are trying to obtain.
Definition: exact_solution.C:43
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