libMesh
Functions
vector_fe_ex1.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 int component, const Real x, const Real y, 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 176 of file vector_fe_ex1.C.

178 {
179 
180  // It is a good idea to make sure we are assembling
181  // the proper system.
182  libmesh_assert_equal_to (system_name, "Poisson");
183 
184  // Get a constant reference to the mesh object.
185  const MeshBase & mesh = es.get_mesh();
186 
187  // The dimension that we are running
188  const unsigned int dim = mesh.mesh_dimension();
189 
190  // Get a reference to the LinearImplicitSystem we are solving
191  LinearImplicitSystem & system = es.get_system<LinearImplicitSystem> ("Poisson");
192 
193  // A reference to the DofMap object for this system. The DofMap
194  // object handles the index translation from node and element numbers
195  // to degree of freedom numbers. We will talk more about the DofMap
196  // in future examples.
197  const DofMap & dof_map = system.get_dof_map();
198 
199  // Get a constant reference to the Finite Element type
200  // for the first (and only) variable in the system.
201  FEType fe_type = dof_map.variable_type(0);
202 
203  // Build a Finite Element object of the specified type.
204  // Note that FEVectorBase is a typedef for the templated FE
205  // class.
206  std::unique_ptr<FEVectorBase> fe (FEVectorBase::build(dim, fe_type));
207 
208  // A 5th order Gauss quadrature rule for numerical integration.
209  QGauss qrule (dim, FIFTH);
210 
211  // Tell the finite element object to use our quadrature rule.
212  fe->attach_quadrature_rule (&qrule);
213 
214  // Declare a special finite element object for
215  // boundary integration.
216  std::unique_ptr<FEVectorBase> fe_face (FEVectorBase::build(dim, fe_type));
217 
218  // Boundary integration requires one quadrature rule,
219  // with dimensionality one less than the dimensionality
220  // of the element.
221  QGauss qface(dim-1, FIFTH);
222 
223  // Tell the finite element object to use our
224  // quadrature rule.
225  fe_face->attach_quadrature_rule (&qface);
226 
227  // Here we define some references to cell-specific data that
228  // will be used to assemble the linear system.
229  //
230  // The element Jacobian * quadrature weight at each integration point.
231  const std::vector<Real> & JxW = fe->get_JxW();
232 
233  // The physical XY locations of the quadrature points on the element.
234  // These might be useful for evaluating spatially varying material
235  // properties at the quadrature points.
236  const std::vector<Point> & q_point = fe->get_xyz();
237 
238  // The element shape functions evaluated at the quadrature points.
239  // Notice the shape functions are a vector rather than a scalar.
240  const std::vector<std::vector<RealGradient>> & phi = fe->get_phi();
241 
242  // The element shape function gradients evaluated at the quadrature
243  // points. Notice that the shape function gradients are a tensor.
244  const std::vector<std::vector<RealTensor>> & dphi = fe->get_dphi();
245 
246  // Define data structures to contain the element matrix
247  // and right-hand-side vector contribution. Following
248  // basic finite element terminology we will denote these
249  // "Ke" and "Fe". These datatypes are templated on
250  // Number, which allows the same code to work for real
251  // or complex numbers.
254 
255  // This vector will hold the degree of freedom indices for
256  // the element. These define where in the global system
257  // the element degrees of freedom get mapped.
258  std::vector<dof_id_type> dof_indices;
259 
260  // Now we will loop over all the elements in the mesh.
261  // We will compute the element matrix and right-hand-side
262  // contribution.
263  //
264  // Element iterators are a nice way to iterate through all the
265  // elements, or all the elements that have some property. The
266  // iterator el will iterate from the first to the last element on
267  // the local processor. The iterator end_el tells us when to stop.
268  // It is smart to make this one const so that we don't accidentally
269  // mess it up! In case users later modify this program to include
270  // refinement, we will be safe and will only consider the active
271  // elements; hence we use a variant of the active_elem_iterator.
272  for (const auto & elem : mesh.active_local_element_ptr_range())
273  {
274  // Get the degree of freedom indices for the
275  // current element. These define where in the global
276  // matrix and right-hand-side this element will
277  // contribute to.
278  dof_map.dof_indices (elem, dof_indices);
279 
280  // Compute the element-specific data for the current
281  // element. This involves computing the location of the
282  // quadrature points (q_point) and the shape functions
283  // (phi, dphi) for the current element.
284  fe->reinit (elem);
285 
286  // Zero the element matrix and right-hand side before
287  // summing them. We use the resize member here because
288  // the number of degrees of freedom might have changed from
289  // the last element. Note that this will be the case if the
290  // element type is different (i.e. the last element was a
291  // triangle, now we are on a quadrilateral).
292 
293  // The DenseMatrix::resize() and the DenseVector::resize()
294  // members will automatically zero out the matrix and vector.
295  Ke.resize (dof_indices.size(),
296  dof_indices.size());
297 
298  Fe.resize (dof_indices.size());
299 
300  // Now loop over the quadrature points. This handles
301  // the numeric integration.
302  for (unsigned int qp=0; qp<qrule.n_points(); qp++)
303  {
304  // Now we will build the element matrix. This involves
305  // a double loop to integrate the test functions (i) against
306  // the trial functions (j).
307  for (std::size_t i=0; i<phi.size(); i++)
308  for (std::size_t j=0; j<phi.size(); j++)
309  Ke(i,j) += JxW[qp] * dphi[i][qp].contract(dphi[j][qp]);
310 
311  // This is the end of the matrix summation loop
312  // Now we build the element right-hand-side contribution.
313  // This involves a single loop in which we integrate the
314  // "forcing function" in the PDE against the test functions.
315  {
316  const Real x = q_point[qp](0);
317  const Real y = q_point[qp](1);
318  const Real eps = 1.e-3;
319 
320  // "f" is the forcing function for the Poisson equation.
321  // In this case we set f to be a finite difference
322  // Laplacian approximation to the (known) exact solution.
323  //
324  // We will use the second-order accurate FD Laplacian
325  // approximation, which in 2D is
326  //
327  // u_xx + u_yy = (u(i,j-1) + u(i,j+1) +
328  // u(i-1,j) + u(i+1,j) +
329  // -4*u(i,j))/h^2
330  //
331  // Since the value of the forcing function depends only
332  // on the location of the quadrature point (q_point[qp])
333  // we will compute it here, outside of the i-loop
334  const Real fx = -(exact_solution(0, x, y-eps) +
335  exact_solution(0, x, y+eps) +
336  exact_solution(0, x-eps, y) +
337  exact_solution(0, x+eps, y) -
338  4.*exact_solution(0, x, y))/eps/eps;
339 
340  const Real fy = -(exact_solution(1, x, y-eps) +
341  exact_solution(1, x, y+eps) +
342  exact_solution(1, x-eps, y) +
343  exact_solution(1, x+eps, y) -
344  4.*exact_solution(1, x, y))/eps/eps;
345 
346  const RealGradient f(fx, fy);
347 
348  for (std::size_t i=0; i<phi.size(); i++)
349  Fe(i) += JxW[qp]*f*phi[i][qp];
350  }
351  }
352 
353  // We have now reached the end of the RHS summation,
354  // and the end of quadrature point loop, so
355  // the interior element integration has
356  // been completed. However, we have not yet addressed
357  // boundary conditions. For this example we will only
358  // consider simple Dirichlet boundary conditions.
359  //
360  // There are several ways Dirichlet boundary conditions
361  // can be imposed. A simple approach, which works for
362  // interpolary bases like the standard Lagrange polynomials,
363  // is to assign function values to the
364  // degrees of freedom living on the domain boundary. This
365  // works well for interpolary bases, but is more difficult
366  // when non-interpolary (e.g Legendre or Hierarchic) bases
367  // are used.
368  //
369  // Dirichlet boundary conditions can also be imposed with a
370  // "penalty" method. In this case essentially the L2 projection
371  // of the boundary values are added to the matrix. The
372  // projection is multiplied by some large factor so that, in
373  // floating point arithmetic, the existing (smaller) entries
374  // in the matrix and right-hand-side are effectively ignored.
375  //
376  // This amounts to adding a term of the form (in latex notation)
377  //
378  // \frac{1}{\epsilon} \int_{\delta \Omega} \phi_i \phi_j = \frac{1}{\epsilon} \int_{\delta \Omega} u \phi_i
379  //
380  // where
381  //
382  // \frac{1}{\epsilon} is the penalty parameter, defined such that \epsilon << 1
383  {
384  // The following loop is over the sides of the element.
385  // If the element has no neighbor on a side then that
386  // side MUST live on a boundary of the domain.
387  for (auto side : elem->side_index_range())
388  if (elem->neighbor_ptr(side) == nullptr)
389  {
390  // The value of the shape functions at the quadrature
391  // points.
392  const std::vector<std::vector<RealGradient>> & phi_face = fe_face->get_phi();
393 
394  // The Jacobian * Quadrature Weight at the quadrature
395  // points on the face.
396  const std::vector<Real> & JxW_face = fe_face->get_JxW();
397 
398  // The XYZ locations (in physical space) of the
399  // quadrature points on the face. This is where
400  // we will interpolate the boundary value function.
401  const std::vector<Point> & qface_point = fe_face->get_xyz();
402 
403  // Compute the shape function values on the element
404  // face.
405  fe_face->reinit(elem, side);
406 
407  // Loop over the face quadrature points for integration.
408  for (unsigned int qp=0; qp<qface.n_points(); qp++)
409  {
410  // The location on the boundary of the current
411  // face quadrature point.
412  const Real xf = qface_point[qp](0);
413  const Real yf = qface_point[qp](1);
414 
415  // The penalty value. \frac{1}{\epsilon}
416  // in the discussion above.
417  const Real penalty = 1.e10;
418 
419  // The boundary values.
420  const RealGradient f(exact_solution(0, xf, yf),
421  exact_solution(1, xf, yf));
422 
423  // Matrix contribution of the L2 projection.
424  for (std::size_t i=0; i<phi_face.size(); i++)
425  for (std::size_t j=0; j<phi_face.size(); j++)
426  Ke(i,j) += JxW_face[qp]*penalty*phi_face[i][qp]*phi_face[j][qp];
427 
428  // Right-hand-side contribution of the L2
429  // projection.
430  for (std::size_t i=0; i<phi_face.size(); i++)
431  Fe(i) += JxW_face[qp]*penalty*f*phi_face[i][qp];
432  }
433  }
434  }
435 
436  // We have now finished the quadrature point loop,
437  // and have therefore applied all the boundary conditions.
438 
439  // If this assembly program were to be used on an adaptive mesh,
440  // we would have to apply any hanging node constraint equations
441  //dof_map.constrain_element_matrix_and_vector (Ke, Fe, dof_indices);
442 
443  // The element matrix and right-hand-side are now built
444  // for this element. Add them to the global matrix and
445  // right-hand-side vector. The SparseMatrix::add_matrix()
446  // and NumericVector::add_vector() members do this for us.
447  system.matrix->add_matrix (Ke, dof_indices);
448  system.rhs->add_vector (Fe, dof_indices);
449  }
450 
451  // All done!
452 }

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

◆ assemble_poisson() [2/2]

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

Referenced by main().

◆ exact_solution()

Real exact_solution ( const int  component,
const Real  x,
const Real  y,
const Real  z = 0. 
)

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 39 of file exact_solution.C.

43 {
44  static const Real pi = acos(-1.);
45 
46  switch (component)
47  {
48  case 0:
49  return cos(.5*pi*x)*sin(.5*pi*y)*cos(.5*pi*z);
50  case 1:
51  return sin(.5*pi*x)*cos(.5*pi*y)*cos(.5*pi*z);
52  case 2:
53  return sin(.5*pi*x)*cos(.5*pi*y)*cos(.5*pi*z)*cos(.5*pi*x*y*z);
54  default:
55  libmesh_error_msg("Invalid component = " << component);
56  }
57 
58  // dummy
59  return 0.0;
60 }

Referenced by assemble_poisson().

◆ main()

int main ( int  argc,
char **  argv 
)

Definition at line 79 of file vector_fe_ex1.C.

80 {
81  // Initialize libraries.
82  LibMeshInit init (argc, argv);
83 
84  // This example requires a linear solver package.
85  libmesh_example_requires(libMesh::default_solver_package() != INVALID_SOLVER_PACKAGE,
86  "--enable-petsc, --enable-trilinos, or --enable-eigen");
87 
88  // Brief message to the user regarding the program name
89  // and command line arguments.
90  libMesh::out << "Running " << argv[0];
91 
92  for (int i=1; i<argc; i++)
93  libMesh::out << " " << argv[i];
94 
95  libMesh::out << std::endl << std::endl;
96 
97  // Skip this 2D example if libMesh was compiled as 1D-only.
98  libmesh_example_requires(2 <= LIBMESH_DIM, "2D support");
99 
100  // Create a mesh, with dimension to be overridden later, on the
101  // default MPI communicator.
102  Mesh mesh(init.comm());
103 
104  // Use the MeshTools::Generation mesh generator to create a uniform
105  // 2D grid on the square [-1,1]^2. We instruct the mesh generator
106  // to build a mesh of 15x15 QUAD9 elements.
108  15, 15,
109  -1., 1.,
110  -1., 1.,
111  QUAD9);
112 
113  // Print information about the mesh to the screen.
114  mesh.print_info();
115 
116  // Create an equation systems object.
117  EquationSystems equation_systems (mesh);
118 
119  // Declare the Poisson system and its variables.
120  // The Poisson system is another example of a steady system.
121  equation_systems.add_system<LinearImplicitSystem> ("Poisson");
122 
123  // Adds the variable "u" to "Poisson". "u"
124  // will be approximated using second-order approximation
125  // using vector Lagrange elements. Since the mesh is 2-D, "u" will
126  // have two components.
127  equation_systems.get_system("Poisson").add_variable("u", SECOND, LAGRANGE_VEC);
128 
129  // Give the system a pointer to the matrix assembly
130  // function. This will be called when needed by the
131  // library.
132  equation_systems.get_system("Poisson").attach_assemble_function (assemble_poisson);
133 
134  // Initialize the data structures for the equation system.
135  equation_systems.init();
136 
137  // Prints information about the system to the screen.
138  equation_systems.print_info();
139 
140  // Solve the system "Poisson". Note that calling this
141  // member will assemble the linear system and invoke
142  // the default numerical solver. With PETSc the solver can be
143  // controlled from the command line. For example,
144  // you can invoke conjugate gradient with:
145  //
146  // ./vector_fe_ex1 -ksp_type cg
147  //
148  // You can also get a nice X-window that monitors the solver
149  // convergence with:
150  //
151  // ./vector_fe_ex1 -ksp_xmonitor
152  //
153  // if you linked against the appropriate X libraries when you
154  // built PETSc.
155  equation_systems.get_system("Poisson").solve();
156 
157 #ifdef LIBMESH_HAVE_EXODUS_API
158  ExodusII_IO(mesh).write_equation_systems("out.e", equation_systems);
159 #endif
160 
161 #ifdef LIBMESH_HAVE_GMV
162  GMVIO(mesh).write_equation_systems("out.gmv", equation_systems);
163 #endif
164 
165  // All done.
166  return 0;
167 }

References libMesh::EquationSystems::add_system(), assemble_poisson(), libMesh::MeshTools::Generation::build_square(), libMesh::default_solver_package(), libMesh::EquationSystems::get_system(), libMesh::TriangleWrapper::init(), libMesh::EquationSystems::init(), libMesh::INVALID_SOLVER_PACKAGE, libMesh::LAGRANGE_VEC, mesh, libMesh::out, libMesh::EquationSystems::print_info(), libMesh::MeshBase::print_info(), libMesh::QUAD9, libMesh::SECOND, 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::MeshBase::active_local_element_ptr_range
virtual SimpleRange< element_iterator > active_local_element_ptr_range()=0
libMesh::EquationSystems::get_system
const T_sys & get_system(const std::string &name) const
Definition: equation_systems.h:757
assemble_poisson
void assemble_poisson(EquationSystems &es, const std::string &system_name)
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::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
libMesh::SECOND
Definition: enum_order.h: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::VectorValue< Real >
libMesh::TriangleWrapper::init
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
libMesh::MeshTools::Generation::build_square
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.
Definition: mesh_generation.C:1501
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::INVALID_SOLVER_PACKAGE
Definition: enum_solver_package.h:43
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::ImplicitSystem::matrix
SparseMatrix< Number > * matrix
The system matrix.
Definition: implicit_system.h:393
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
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::QUAD9
Definition: enum_elem_type.h:43
libMesh::System::get_dof_map
const DofMap & get_dof_map() const
Definition: system.h:2099
libMesh::LAGRANGE_VEC
Definition: enum_fe_family.h:60
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
exact_solution
Real exact_solution(const int component, 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:39
libMesh::DenseVector< Number >