libMesh
subdomains_ex1.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>Subdomains Example 1 - Solving on a Subdomain</h1>
21 // \author Tim Kroger
22 // \date 2010
23 //
24 // This example builds on the example 4 by showing what to do in
25 // order to solve an equation only on a subdomain.
26 
27 
28 // C++ include files that we need
29 #include <iostream>
30 #include <algorithm>
31 #include <math.h>
32 
33 // Basic include file needed for the mesh functionality.
34 #include "libmesh/libmesh.h"
35 #include "libmesh/mesh.h"
36 #include "libmesh/mesh_generation.h"
37 #include "libmesh/exodusII_io.h"
38 #include "libmesh/gmv_io.h"
39 #include "libmesh/gnuplot_io.h"
40 #include "libmesh/linear_implicit_system.h"
41 #include "libmesh/equation_systems.h"
42 
43 // Define the Finite Element object.
44 #include "libmesh/fe.h"
45 
46 // Define Gauss quadrature rules.
47 #include "libmesh/quadrature_gauss.h"
48 
49 // Define the DofMap, which handles degree of freedom
50 // indexing.
51 #include "libmesh/dof_map.h"
52 
53 // Define useful datatypes for finite element
54 // matrix and vector components.
55 #include "libmesh/sparse_matrix.h"
56 #include "libmesh/numeric_vector.h"
57 #include "libmesh/dense_matrix.h"
58 #include "libmesh/dense_vector.h"
59 
60 // Define the PerfLog, a performance logging utility.
61 // It is useful for timing events in a code and giving
62 // you an idea where bottlenecks lie.
63 #include "libmesh/perf_log.h"
64 
65 // The definition of a geometric element
66 #include "libmesh/elem.h"
67 
68 #include "libmesh/mesh_refinement.h"
69 
70 // Classes needed for subdomain computation.
71 #include "libmesh/system_subset_by_subdomain.h"
72 
73 #include "libmesh/string_to_enum.h"
74 #include "libmesh/getpot.h"
75 #include "libmesh/enum_solver_package.h"
76 
77 // Bring in everything from the libMesh namespace
78 using namespace libMesh;
79 
80 
81 
82 // Function prototype. This is the function that will assemble
83 // the linear system for our Poisson problem. Note that the
84 // function will take the EquationSystems object and the
85 // name of the system we are assembling as input. From the
86 // EquationSystems object we have access to the Mesh and
87 // other objects we might need.
89  const std::string & system_name);
90 
91 // Exact solution function prototype.
92 Real exact_solution (const Real x,
93  const Real y = 0.,
94  const Real z = 0.);
95 
96 // Begin the main program.
97 int main (int argc, char ** argv)
98 {
99  // Initialize libMesh and any dependent libraries, like in example 2.
100  LibMeshInit init (argc, argv);
101 
102  // Only our PETSc interface currently supports solves restricted to
103  // subdomains
104  libmesh_example_requires(libMesh::default_solver_package() == PETSC_SOLVERS, "--enable-petsc");
105 
106  // Skip adaptive examples on a non-adaptive libMesh build
107 #ifndef LIBMESH_ENABLE_AMR
108  libmesh_example_requires(false, "--enable-amr");
109 #else
110 
111  // Declare a performance log for the main program
112  // PerfLog perf_main("Main Program");
113 
114  // Create a GetPot object to parse the command line
115  GetPot command_line (argc, argv);
116 
117  // Check for proper calling arguments.
118  if (argc < 3)
119  {
120  // This handy function will print the file name, line number,
121  // specified message, and then throw an exception.
122  libmesh_error_msg("Usage:\n" << "\t " << argv[0] << " -d 2(3)" << " -n 15");
123  }
124 
125  // Brief message to the user regarding the program name
126  // and command line arguments.
127  else
128  {
129  libMesh::out << "Running " << argv[0];
130 
131  for (int i=1; i<argc; i++)
132  libMesh::out << " " << argv[i];
133 
134  libMesh::out << std::endl << std::endl;
135  }
136 
137 
138  // Read problem dimension from command line. Use int
139  // instead of unsigned since the GetPot overload is ambiguous
140  // otherwise.
141  int dim = 2;
142  if (command_line.search(1, "-d"))
143  dim = command_line.next(dim);
144 
145  // Skip higher-dimensional examples on a lower-dimensional libMesh build
146  libmesh_example_requires(dim <= LIBMESH_DIM, "2D/3D support");
147 
148  // Create a mesh with user-defined dimension.
149  // Read number of elements from command line
150  int ps = 15;
151  if (command_line.search(1, "-n"))
152  ps = command_line.next(ps);
153 
154  // Read FE order from command line
155  std::string order = "FIRST";
156  if (command_line.search(2, "-Order", "-o"))
157  order = command_line.next(order);
158 
159  // Read FE Family from command line
160  std::string family = "LAGRANGE";
161  if (command_line.search(2, "-FEFamily", "-f"))
162  family = command_line.next(family);
163 
164  // Cannot use discontinuous basis.
165  if ((family == "MONOMIAL") || (family == "XYZ"))
166  libmesh_error_msg("This example requires a C^0 (or higher) FE basis.");
167 
168  // Create a mesh, with dimension to be overridden later, on the
169  // default MPI communicator.
170  Mesh mesh(init.comm());
171 
172  // Use the MeshTools::Generation mesh generator to create a uniform
173  // grid on the square [-1,1]^D. We instruct the mesh generator
174  // to build a mesh of 8x8 Quad9 elements in 2D, or Hex27
175  // elements in 3D. Building these higher-order elements allows
176  // us to use higher-order approximation, as in example 3.
177 
178  Real halfwidth = dim > 1 ? 1. : 0.;
179  Real halfheight = dim > 2 ? 1. : 0.;
180 
181  if ((family == "LAGRANGE") && (order == "FIRST"))
182  {
183  // No reason to use high-order geometric elements if we are
184  // solving with low-order finite elements.
186  ps,
187  (dim>1) ? ps : 0,
188  (dim>2) ? ps : 0,
189  -1., 1.,
190  -halfwidth, halfwidth,
191  -halfheight, halfheight,
192  (dim==1) ? EDGE2 :
193  ((dim == 2) ? QUAD4 : HEX8));
194  }
195 
196  else
197  {
199  ps,
200  (dim>1) ? ps : 0,
201  (dim>2) ? ps : 0,
202  -1., 1.,
203  -halfwidth, halfwidth,
204  -halfheight, halfheight,
205  (dim==1) ? EDGE3 :
206  ((dim == 2) ? QUAD9 : HEX27));
207  }
208 
209 
210  // To demonstrate solving on a subdomain, we will solve only on the
211  // interior of a circle (ball in 3d) with radius 0.8. So show that
212  // this also works well on locally refined meshes, we refine once
213  // all elements that are located on the boundary of this circle (or
214  // ball).
215  {
216  // A MeshRefinement object is needed to refine meshes.
217  MeshRefinement meshRefinement(mesh);
218 
219  // Loop over all elements.
220  for (auto & elem : mesh.element_ptr_range())
221  if (elem->active())
222  {
223  // Just check whether the current element has at least one
224  // node inside and one node outside the circle.
225  bool node_in = false;
226  bool node_out = false;
227  for (auto & n : elem->node_ref_range())
228  {
229  Real d = n.norm();
230  if (d<0.8)
231  node_in = true;
232  else
233  node_out = true;
234  }
235  if (node_in && node_out)
236  elem->set_refinement_flag(Elem::REFINE);
237  else
238  elem->set_refinement_flag(Elem::DO_NOTHING);
239  }
240  else
241  elem->set_refinement_flag(Elem::INACTIVE);
242 
243  // Now actually refine.
244  meshRefinement.refine_elements();
245  }
246 
247  // Print information about the mesh to the screen.
248  mesh.print_info();
249 
250  // Now set the subdomain_id of all elements whose centroid is inside
251  // the circle to 1.
252  for (auto elem : mesh.element_ptr_range())
253  {
254  Real d = elem->centroid().norm();
255  if (d < 0.8)
256  elem->subdomain_id() = 1;
257  }
258 
259  // Create an equation systems object.
260  EquationSystems equation_systems (mesh);
261 
262  // Declare the system and its variables.
263  // Create a system named "Poisson"
264  LinearImplicitSystem & system =
265  equation_systems.add_system<LinearImplicitSystem> ("Poisson");
266 
267 
268  // Add the variable "u" to "Poisson". "u"
269  // will be approximated using second-order approximation.
270  system.add_variable("u",
271  Utility::string_to_enum<Order> (order),
272  Utility::string_to_enum<FEFamily>(family));
273 
274  // Give the system a pointer to the matrix assembly
275  // function.
277 
278  // Initialize the data structures for the equation system.
279  equation_systems.init();
280 
281  // Print information about the system to the screen.
282  equation_systems.print_info();
283  mesh.print_info();
284 
285  // Restrict solves to those elements that have subdomain_id set to 1.
286  std::set<subdomain_id_type> id_list;
287  id_list.insert(1);
289  SystemSubsetBySubdomain subset(system, selection);
290  system.restrict_solve_to(&subset, SUBSET_ZERO);
291 
292  // Note that using SUBSET_ZERO will cause all dofs outside the
293  // subdomain to be cleared. This will, however, cause some hanging
294  // nodes outside the subdomain to have inconsistent values.
295 
296  // Solve the system "Poisson", just like example 2.
297  equation_systems.get_system("Poisson").solve();
298 
299  // After solving the system write the solution
300  // to a GMV-formatted plot file.
301  if (dim == 1)
302  {
303  GnuPlotIO plot(mesh, "Subdomains Example 1, 1D", GnuPlotIO::GRID_ON);
304  plot.write_equation_systems("gnuplot_script", equation_systems);
305  }
306  else
307  {
309  "out_3.gmv" : "out_2.gmv", equation_systems);
310 #ifdef LIBMESH_HAVE_EXODUS_API
312  "out_3.e" : "out_2.e", equation_systems);
313 #endif // #ifdef LIBMESH_HAVE_EXODUS_API
314  }
315 
316 #endif // #ifndef LIBMESH_ENABLE_AMR
317 
318  // All done.
319  return 0;
320 }
321 
322 
323 
324 
325 // We now define the matrix assembly function for the
326 // Poisson system. We need to first compute element
327 // matrices and right-hand sides, and then take into
328 // account the boundary conditions, which will be handled
329 // via a penalty method.
331  const std::string & libmesh_dbg_var(system_name))
332 {
333  // It is a good idea to make sure we are assembling
334  // the proper system.
335  libmesh_assert_equal_to (system_name, "Poisson");
336 
337  // Declare a performance log. Give it a descriptive
338  // string to identify what part of the code we are
339  // logging, since there may be many PerfLogs in an
340  // application.
341  PerfLog perf_log ("Matrix Assembly");
342 
343  // Get a constant reference to the mesh object.
344  const MeshBase & mesh = es.get_mesh();
345 
346  // The dimension that we are running
347  const unsigned int dim = mesh.mesh_dimension();
348 
349  // Get a reference to the LinearImplicitSystem we are solving
350  LinearImplicitSystem & system = es.get_system<LinearImplicitSystem>("Poisson");
351 
352  // A reference to the DofMap object for this system. The DofMap
353  // object handles the index translation from node and element numbers
354  // to degree of freedom numbers. We will talk more about the DofMap
355  // in future examples.
356  const DofMap & dof_map = system.get_dof_map();
357 
358  // Get a constant reference to the Finite Element type
359  // for the first (and only) variable in the system.
360  FEType fe_type = dof_map.variable_type(0);
361 
362  // Build a Finite Element object of the specified type. Since the
363  // FEBase::build() member dynamically creates memory we will
364  // store the object as a std::unique_ptr<FEBase>. This can be thought
365  // of as a pointer that will clean up after itself.
366  std::unique_ptr<FEBase> fe (FEBase::build(dim, fe_type));
367 
368  // A 5th order Gauss quadrature rule for numerical integration.
369  QGauss qrule (dim, FIFTH);
370 
371  // Tell the finite element object to use our quadrature rule.
372  fe->attach_quadrature_rule (&qrule);
373 
374  // Declare a special finite element object for
375  // boundary integration.
376  std::unique_ptr<FEBase> fe_face (FEBase::build(dim, fe_type));
377 
378  // Boundary integration requires one quadrature rule,
379  // with dimensionality one less than the dimensionality
380  // of the element.
381  QGauss qface(dim-1, FIFTH);
382 
383  // Tell the finite element object to use our
384  // quadrature rule.
385  fe_face->attach_quadrature_rule (&qface);
386 
387  // Here we define some references to cell-specific data that
388  // will be used to assemble the linear system.
389  // We begin with the element Jacobian * quadrature weight at each
390  // integration point.
391  const std::vector<Real> & JxW = fe->get_JxW();
392 
393  // The physical XY locations of the quadrature points on the element.
394  // These might be useful for evaluating spatially varying material
395  // properties at the quadrature points.
396  const std::vector<Point> & q_point = fe->get_xyz();
397 
398  // The element shape functions evaluated at the quadrature points.
399  const std::vector<std::vector<Real>> & phi = fe->get_phi();
400 
401  // The element shape function gradients evaluated at the quadrature
402  // points.
403  const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
404 
405  // Define data structures to contain the element matrix
406  // and right-hand-side vector contribution. Following
407  // basic finite element terminology we will denote these
408  // "Ke" and "Fe". More detail is in example 3.
411 
412  // This vector will hold the degree of freedom indices for
413  // the element. These define where in the global system
414  // the element degrees of freedom get mapped.
415  std::vector<dof_id_type> dof_indices;
416 
417  // Now we will loop over all the elements in the mesh.
418  // We will compute the element matrix and right-hand-side
419  // contribution. See example 3 for a discussion of the
420  // element iterators.
421  for (const auto & elem : mesh.active_local_element_ptr_range())
422  {
423  // Elements with subdomain_id other than 1 are not in the active
424  // subdomain. We don't assemble anything for them.
425  if (elem->subdomain_id()==1)
426  {
427  // Start logging the shape function initialization.
428  // This is done through a simple function call with
429  // the name of the event to log.
430  perf_log.push("elem init");
431 
432  // Get the degree of freedom indices for the
433  // current element. These define where in the global
434  // matrix and right-hand-side this element will
435  // contribute to.
436  dof_map.dof_indices (elem, dof_indices);
437 
438  // Compute the element-specific data for the current
439  // element. This involves computing the location of the
440  // quadrature points (q_point) and the shape functions
441  // (phi, dphi) for the current element.
442  fe->reinit (elem);
443 
444  // Zero the element matrix and right-hand side before
445  // summing them. We use the resize member here because
446  // the number of degrees of freedom might have changed from
447  // the last element. Note that this will be the case if the
448  // element type is different (i.e. the last element was a
449  // triangle, now we are on a quadrilateral).
450  Ke.resize (dof_indices.size(),
451  dof_indices.size());
452 
453  Fe.resize (dof_indices.size());
454 
455  // Stop logging the shape function initialization.
456  // If you forget to stop logging an event the PerfLog
457  // object will probably catch the error and abort.
458  perf_log.pop("elem init");
459 
460  // Now we will build the element matrix. This involves
461  // a double loop to integrate the test functions (i) against
462  // the trial functions (j).
463  //
464  // We have split the numeric integration into two loops
465  // so that we can log the matrix and right-hand-side
466  // computation separately.
467  //
468  // Now start logging the element matrix computation
469  perf_log.push ("Ke");
470 
471  for (unsigned int qp=0; qp<qrule.n_points(); qp++)
472  for (std::size_t i=0; i<phi.size(); i++)
473  for (std::size_t j=0; j<phi.size(); j++)
474  Ke(i,j) += JxW[qp]*(dphi[i][qp]*dphi[j][qp]);
475 
476 
477  // Stop logging the matrix computation
478  perf_log.pop ("Ke");
479 
480  // Now we build the element right-hand-side contribution.
481  // This involves a single loop in which we integrate the
482  // "forcing function" in the PDE against the test functions.
483  //
484  // Start logging the right-hand-side computation
485  perf_log.push ("Fe");
486 
487  for (unsigned int qp=0; qp<qrule.n_points(); qp++)
488  {
489  // fxy is the forcing function for the Poisson equation.
490  // In this case we set fxy to be a finite difference
491  // Laplacian approximation to the (known) exact solution.
492  //
493  // We will use the second-order accurate FD Laplacian
494  // approximation, which in 2D on a structured grid is
495  //
496  // u_xx + u_yy = (u(i-1,j) + u(i+1,j) +
497  // u(i,j-1) + u(i,j+1) +
498  // -4*u(i,j))/h^2
499  //
500  // Since the value of the forcing function depends only
501  // on the location of the quadrature point (q_point[qp])
502  // we will compute it here, outside of the i-loop
503  const Real x = q_point[qp](0);
504 #if LIBMESH_DIM > 1
505  const Real y = q_point[qp](1);
506 #else
507  const Real y = 0.;
508 #endif
509 #if LIBMESH_DIM > 2
510  const Real z = q_point[qp](2);
511 #else
512  const Real z = 0.;
513 #endif
514  const Real eps = 1.e-3;
515 
516  const Real uxx = (exact_solution(x-eps, y, z) +
517  exact_solution(x+eps, y, z) +
518  -2.*exact_solution(x, y, z))/eps/eps;
519 
520  const Real uyy = (exact_solution(x, y-eps, z) +
521  exact_solution(x, y+eps, z) +
522  -2.*exact_solution(x, y, z))/eps/eps;
523 
524  const Real uzz = (exact_solution(x, y, z-eps) +
525  exact_solution(x, y, z+eps) +
526  -2.*exact_solution(x, y, z))/eps/eps;
527 
528  Real fxy;
529  if (dim==1)
530  {
531  // In 1D, compute the rhs by differentiating the
532  // exact solution twice.
533  const Real pi = libMesh::pi;
534  fxy = (0.25*pi*pi)*sin(.5*pi*x);
535  }
536  else
537  {
538  fxy = - (uxx + uyy + ((dim==2) ? 0. : uzz));
539  }
540 
541  // Add the RHS contribution
542  for (std::size_t i=0; i<phi.size(); i++)
543  Fe(i) += JxW[qp]*fxy*phi[i][qp];
544  }
545 
546  // Stop logging the right-hand-side computation
547  perf_log.pop ("Fe");
548 
549  // At this point the interior element integration has
550  // been completed. However, we have not yet addressed
551  // boundary conditions. For this example we will only
552  // consider simple Dirichlet boundary conditions imposed
553  // via the penalty method. This is discussed at length in
554  // example 3.
555  {
556  // Start logging the boundary condition computation. We use a
557  // macro to log everything in this scope.
558  LOG_SCOPE_WITH("BCs", "", perf_log);
559 
560  // The following loops over the sides of the element. If
561  // the element has no neighbor on a side then that side
562  // MUST live on a boundary of the domain. If there is a
563  // neighbor, check that neighbor's subdomain_id; if that
564  // is different from 1, the side is also located on the
565  // boundary.
566  for (auto side : elem->side_index_range())
567  if ((elem->neighbor_ptr(side) == nullptr) ||
568  (elem->neighbor_ptr(side)->subdomain_id()!=1))
569  {
570 
571  // The penalty value. \frac{1}{\epsilon}
572  // in the discussion above.
573  const Real penalty = 1.e10;
574 
575  // The value of the shape functions at the quadrature
576  // points.
577  const std::vector<std::vector<Real>> & phi_face = fe_face->get_phi();
578 
579  // The Jacobian * Quadrature Weight at the quadrature
580  // points on the face.
581  const std::vector<Real> & JxW_face = fe_face->get_JxW();
582 
583  // The XYZ locations (in physical space) of the
584  // quadrature points on the face. This is where
585  // we will interpolate the boundary value function.
586  const std::vector<Point> & qface_point = fe_face->get_xyz();
587 
588  // Compute the shape function values on the element
589  // face.
590  fe_face->reinit(elem, side);
591 
592  // Loop over the face quadrature points for integration.
593  for (unsigned int qp=0; qp<qface.n_points(); qp++)
594  {
595  // The location on the boundary of the current
596  // face quadrature point.
597  const Real xf = qface_point[qp](0);
598 #if LIBMESH_DIM > 1
599  const Real yf = qface_point[qp](1);
600 #else
601  const Real yf = 0.;
602 #endif
603 #if LIBMESH_DIM > 2
604  const Real zf = qface_point[qp](2);
605 #else
606  const Real zf = 0.;
607 #endif
608 
609 
610  // The boundary value.
611  const Real value = exact_solution(xf, yf, zf);
612 
613  // Matrix contribution of the L2 projection.
614  for (std::size_t i=0; i<phi_face.size(); i++)
615  for (std::size_t j=0; j<phi_face.size(); j++)
616  Ke(i,j) += JxW_face[qp]*penalty*phi_face[i][qp]*phi_face[j][qp];
617 
618  // Right-hand-side contribution of the L2
619  // projection.
620  for (std::size_t i=0; i<phi_face.size(); i++)
621  Fe(i) += JxW_face[qp]*penalty*value*phi_face[i][qp];
622  }
623  }
624  }
625 
626  // If this assembly program were to be used on an adaptive mesh,
627  // we would have to apply any hanging node constraint equations
628  dof_map.constrain_element_matrix_and_vector (Ke, Fe, dof_indices);
629 
630  // The element matrix and right-hand-side are now built
631  // for this element. Add them to the global matrix and
632  // right-hand-side vector. The SparseMatrix::add_matrix()
633  // and NumericVector::add_vector() members do this for us.
634  // Start logging the insertion of the local (element)
635  // matrix and vector into the global matrix and vector
636  LOG_SCOPE_WITH("matrix insertion", "", perf_log);
637 
638  system.matrix->add_matrix (Ke, dof_indices);
639  system.rhs->add_vector (Fe, dof_indices);
640  }
641  }
642 
643  // That's it. We don't need to do anything else to the
644  // PerfLog. When it goes out of scope (at this function return)
645  // it will print its log to the screen. Pretty easy, huh?
646 }
libMesh::pi
const Real pi
.
Definition: libmesh.h:237
libMesh::SystemSubsetBySubdomain
This class represents a subset of the dofs of a System, selected by the subdomain_id and possible the...
Definition: system_subset_by_subdomain.h:45
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::MeshRefinement::refine_elements
bool refine_elements()
Only refines the user-requested elements.
Definition: mesh_refinement.C:683
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::PETSC_SOLVERS
Definition: enum_solver_package.h:36
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
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::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::LinearImplicitSystem::restrict_solve_to
virtual void restrict_solve_to(const SystemSubset *subset, const SubsetSolveMode subset_solve_mode=SUBSET_ZERO) override
After calling this method, any solve will be limited to the given subset.
Definition: linear_implicit_system.C:96
libMesh::GMVIO
This class implements writing meshes in the GMV format.
Definition: gmv_io.h:54
libMesh::ExodusII_IO
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
Definition: exodusII_io.h:51
dim
unsigned int dim
Definition: adaptivity_ex3.C:113
libMesh::SUBSET_ZERO
Definition: enum_subset_solve_mode.h:37
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::MeshRefinement
Implements (adaptive) mesh refinement algorithms for a MeshBase.
Definition: mesh_refinement.h:61
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
assemble_poisson
void assemble_poisson(EquationSystems &es, const std::string &system_name)
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::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::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::ImplicitSystem::matrix
SparseMatrix< Number > * matrix
The system matrix.
Definition: implicit_system.h:393
libMesh::Elem::DO_NOTHING
Definition: elem.h:1170
libMesh::LibMeshInit
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
Definition: libmesh.h:83
libMesh::Elem::REFINE
Definition: elem.h:1171
libMesh::SystemSubsetBySubdomain::SubdomainSelectionByList
Selection of subdomain ids by a list.
Definition: system_subset_by_subdomain.h:97
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
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
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::Elem::INACTIVE
Definition: elem.h:1174
libMesh::out
OStreamProxy out
libMesh::DenseVector< Number >
main
int main(int argc, char **argv)
Definition: subdomains_ex1.C:97
libMesh::EDGE2
Definition: enum_elem_type.h:35
libMesh::GnuPlotIO::GRID_ON
Definition: gnuplot_io.h:50