#include <iostream>
#include <algorithm>
#include <math.h>
Basic include files needed for the mesh functionality.
#include "libmesh.h"
#include "mesh.h"
#include "mesh_generation.h"
#include "gmv_io.h"
#include "linear_implicit_system.h"
#include "equation_systems.h"
Define the Finite Element object.
#include "fe.h"
Define Gauss quadrature rules.
#include "quadrature_gauss.h"
Define useful datatypes for finite element
matrix and vector components.
#include "sparse_matrix.h"
#include "numeric_vector.h"
#include "dense_matrix.h"
#include "dense_vector.h"
#include "elem.h"
Define the DofMap, which handles degree of freedom
indexing.
#include "dof_map.h"
#include "boundary_mesh.h"
#include "boundary_info.h"
Function prototype. This is the function that will assemble
the linear system for our Poisson problem. Note that the
function will take the EquationSystems object and the
name of the system we are assembling as input. From the
EquationSystems object we have access to the Mesh and
other objects we might need.
void assemble_poisson(EquationSystems& es,
const std::string& system_name);
Function prototype for the exact solution.
Real exact_solution (const Real x,
const Real y,
const Real z = 0.);
int main (int argc, char** argv)
{
Initialize Petsc, like in example 2.
libMesh::init (argc, argv);
Braces are used to force object scope, like in example 2
{
Brief message to the user regarding the program name
and command line arguments.
std::cout << "Running " << argv[0];
for (int i=1; i<argc; i++)
std::cout << " " << argv[i];
std::cout << std::endl << std::endl;
Create a 2D mesh.
Mesh mesh (2);
Use the MeshTools::Generation mesh generator to create a uniform
grid on the square [-1,1]^2. We instruct the mesh generator
to build a mesh of 15x15 QUAD9 elements. Building QUAD9
elements instead of the default QUAD4's we used in example 2
allow us to use higher-order approximation.
MeshTools::Generation::build_square (mesh,
15, 15,
-1., 1.,
-1., 1.,
QUAD9);
Print information about the mesh to the screen.
Note that 5x5 QUAD9 elements actually has 11x11 nodes,
so this mesh is significantly larger than the one in example 2.
mesh.print_info();
Create an equation systems object.
EquationSystems equation_systems (mesh);
Declare the Poisson system and its variables.
The Poisson system is another example of a steady system.
equation_systems.add_system<LinearImplicitSystem> ("Poisson");
Adds the variable "u" to "Poisson". "u"
will be approximated using second-order approximation.
equation_systems.get_system("Poisson").add_variable("u", SECOND);
Give the system a pointer to the matrix assembly
function. This will be called when needed by the
library.
equation_systems.get_system("Poisson").attach_assemble_function (assemble_poisson);
Initialize the data structures for the equation system.
equation_systems.init();
Prints information about the system to the screen.
equation_systems.print_info();
Solve the system "Poisson". Note that calling this
member will assemble the linear system and invoke
the default Petsc solver, however the solver can be
controlled from the command line. For example,
you can invoke conjugate gradient with:
./ex3 -ksp_type cg
and you can get a nice X-window that monitors the solver convergence with:
./ex3 -ksp_xmonitor
if you linked against the appropriate X libraries when you built PETSc.
./ex3 -ksp_type cg
and you can get a nice X-window that monitors the solver convergence with:
./ex3 -ksp_xmonitor
if you linked against the appropriate X libraries when you built PETSc.
equation_systems.get_system("Poisson").solve();
After solving the system write the solution
to a GMV-formatted plot file.
GMVIO (mesh).write_equation_systems ("out.gmv", equation_systems);
for (unsigned int i=0; i<3; i++)
{
here();
equation_systems.get_system("Poisson").solve();
BoundaryMesh boundary_mesh (mesh.mesh_dimension()-1);
mesh.boundary_info->sync(boundary_mesh, false);
GMVIO(boundary_mesh).write("boundary.gmv");
}
}
All done.
return libMesh::close();
}
We now define the matrix assembly function for the
Poisson system. We need to first compute element
matrices and right-hand sides, and then take into
account the boundary conditions, which will be handled
via a penalty method.
void assemble_poisson(EquationSystems& es,
const std::string& system_name)
{
It is a good idea to make sure we are assembling
the proper system.
assert (system_name == "Poisson");
Get a constant reference to the mesh object.
const Mesh& mesh = es.get_mesh();
The dimension that we are running
const unsigned int dim = mesh.mesh_dimension();
Get a reference to the LinearImplicitSystem we are solving
LinearImplicitSystem& system = es.get_system<LinearImplicitSystem> ("Poisson");
A reference to the DofMap object for this system. The DofMap
object handles the index translation from node and element numbers
to degree of freedom numbers. We will talk more about the DofMap
in future examples.
const DofMap& dof_map = system.get_dof_map();
Get a constant reference to the Finite Element type
for the first (and only) variable in the system.
FEType fe_type = dof_map.variable_type(0);
Build a Finite Element object of the specified type. Since the
FEBase::build() member dynamically creates memory we will
store the object as an AutoPtr. This can be thought
of as a pointer that will clean up after itself. Example 4
describes some advantages of AutoPtr's in the context of
quadrature rules.
AutoPtr<FEBase> fe (FEBase::build(dim, fe_type));
A 5th order Gauss quadrature rule for numerical integration.
QGauss qrule (dim, FIFTH);
Tell the finite element object to use our quadrature rule.
fe->attach_quadrature_rule (&qrule);
Declare a special finite element object for
boundary integration.
AutoPtr<FEBase> fe_face (FEBase::build(dim, fe_type));
Boundary integration requires one quadraure rule,
with dimensionality one less than the dimensionality
of the element.
QGauss qface(dim-1, FIFTH);
Tell the finite element object to use our
quadrature rule.
fe_face->attach_quadrature_rule (&qface);
Here we define some references to cell-specific data that
will be used to assemble the linear system.
The element Jacobian * quadrature weight at each integration point.
The element Jacobian * quadrature weight at each integration point.
const std::vector<Real>& JxW = fe->get_JxW();
The physical XY locations of the quadrature points on the element.
These might be useful for evaluating spatially varying material
properties at the quadrature points.
const std::vector<Point>& q_point = fe->get_xyz();
The element shape functions evaluated at the quadrature points.
const std::vector<std::vector<Real> >& phi = fe->get_phi();
The element shape function gradients evaluated at the quadrature
points.
const std::vector<std::vector<RealGradient> >& dphi = fe->get_dphi();
Define data structures to contain the element matrix
and right-hand-side vector contribution. Following
basic finite element terminology we will denote these
"Ke" and "Fe". These datatypes are templated on
Number, which allows the same code to work for real
or complex numbers.
DenseMatrix<Number> Ke;
DenseVector<Number> Fe;
This vector will hold the degree of freedom indices for
the element. These define where in the global system
the element degrees of freedom get mapped.
std::vector<unsigned int> dof_indices;
Now we will loop over all the elements in the mesh.
We will compute the element matrix and right-hand-side
contribution.
Element iterators are a nice way to iterate through all the elements, or all the elements that have some property. There are many types of element iterators, but here we will use the most basic type, the const_elem_iterator. The iterator el will iterate from the first to the last element. The iterator end_el tells us when to stop. It is smart to make this one const so that we don't accidentally mess it up! const_elem_iterator el (mesh.elements_begin()); const const_elem_iterator end_el (mesh.elements_end());
Element iterators are a nice way to iterate through all the elements, or all the elements that have some property. There are many types of element iterators, but here we will use the most basic type, the const_elem_iterator. The iterator el will iterate from the first to the last element. The iterator end_el tells us when to stop. It is smart to make this one const so that we don't accidentally mess it up! const_elem_iterator el (mesh.elements_begin()); const const_elem_iterator end_el (mesh.elements_end());
MeshBase::const_element_iterator el = mesh.elements_begin();
const MeshBase::const_element_iterator end_el = mesh.elements_end();
Loop over the elements. Note that ++el is preferred to
el++ since the latter requires an unnecessary temporary
object.
for ( ; el != end_el ; ++el)
{
Store a pointer to the element we are currently
working on. This allows for nicer syntax later.
const Elem* elem = *el;
Get the degree of freedom indices for the
current element. These define where in the global
matrix and right-hand-side this element will
contribute to.
dof_map.dof_indices (elem, dof_indices);
Compute the element-specific data for the current
element. This involves computing the location of the
quadrature points (q_point) and the shape functions
(phi, dphi) for the current element.
fe->reinit (elem);
Zero the element matrix and right-hand side before
summing them. We use the resize member here because
the number of degrees of freedom might have changed from
the last element. Note that this will be the case if the
element type is different (i.e. the last element was a
triangle, now we are on a quadrilateral).
The DenseMatrix::resize() and the DenseVector::resize() members will automatically zero out the matrix and vector.
The DenseMatrix::resize() and the DenseVector::resize() members will automatically zero out the matrix and vector.
Ke.resize (dof_indices.size(),
dof_indices.size());
Fe.resize (dof_indices.size());
Now loop over the quadrature points. This handles
the numeric integration.
for (unsigned int qp=0; qp<qrule.n_points(); qp++)
{
Now we will build the element matrix. This involves
a double loop to integrate the test funcions (i) against
the trial functions (j).
for (unsigned int i=0; i<phi.size(); i++)
for (unsigned int j=0; j<phi.size(); j++)
{
Ke(i,j) += JxW[qp]*(dphi[i][qp]*dphi[j][qp]);
}
This is the end of the matrix summation loop
Now we build the element right-hand-side contribution.
This involves a single loop in which we integrate the
"forcing function" in the PDE against the test functions.
{
const Real x = q_point[qp](0);
const Real y = q_point[qp](1);
const Real eps = 1.e-3;
"fxy" is the forcing function for the Poisson equation.
In this case we set fxy to be a finite difference
Laplacian approximation to the (known) exact solution.
We will use the second-order accurate FD Laplacian approximation, which in 2D is
u_xx + u_yy = (u(i,j-1) + u(i,j+1) + u(i-1,j) + u(i+1,j) + -4*u(i,j))/h^2
Since the value of the forcing function depends only on the location of the quadrature point (q_point[qp]) we will compute it here, outside of the i-loop
We will use the second-order accurate FD Laplacian approximation, which in 2D is
u_xx + u_yy = (u(i,j-1) + u(i,j+1) + u(i-1,j) + u(i+1,j) + -4*u(i,j))/h^2
Since the value of the forcing function depends only on the location of the quadrature point (q_point[qp]) we will compute it here, outside of the i-loop
const Real fxy = -(exact_solution(x,y-eps) +
exact_solution(x,y+eps) +
exact_solution(x-eps,y) +
exact_solution(x+eps,y) -
4.*exact_solution(x,y))/eps/eps;
for (unsigned int i=0; i<phi.size(); i++)
Fe(i) += JxW[qp]*fxy*phi[i][qp];
}
}
We have now reached the end of the RHS summation,
and the end of quadrature point loop, so
the interior element integration has
been completed. However, we have not yet addressed
boundary conditions. For this example we will only
consider simple Dirichlet boundary conditions.
There are several ways Dirichlet boundary conditions can be imposed. A simple approach, which works for interpolary bases like the standard Lagrange polynomials, is to assign function values to the degrees of freedom living on the domain boundary. This works well for interpolary bases, but is more difficult when non-interpolary (e.g Legendre or Hierarchic) bases are used.
Dirichlet boundary conditions can also be imposed with a "penalty" method. In this case essentially the L2 projection of the boundary values are added to the matrix. The projection is multiplied by some large factor so that, in floating point arithmetic, the existing (smaller) entries in the matrix and right-hand-side are effectively ignored.
This amounts to adding a term of the form (in latex notation)
\frac{1}{\epsilon} \int_{\delta \Omega} \phi_i \phi_j = \frac{1}{\epsilon} \int_{\delta \Omega} u \phi_i
where
\frac{1}{\epsilon} is the penalty parameter, defined such that \epsilon << 1
There are several ways Dirichlet boundary conditions can be imposed. A simple approach, which works for interpolary bases like the standard Lagrange polynomials, is to assign function values to the degrees of freedom living on the domain boundary. This works well for interpolary bases, but is more difficult when non-interpolary (e.g Legendre or Hierarchic) bases are used.
Dirichlet boundary conditions can also be imposed with a "penalty" method. In this case essentially the L2 projection of the boundary values are added to the matrix. The projection is multiplied by some large factor so that, in floating point arithmetic, the existing (smaller) entries in the matrix and right-hand-side are effectively ignored.
This amounts to adding a term of the form (in latex notation)
\frac{1}{\epsilon} \int_{\delta \Omega} \phi_i \phi_j = \frac{1}{\epsilon} \int_{\delta \Omega} u \phi_i
where
\frac{1}{\epsilon} is the penalty parameter, defined such that \epsilon << 1
{
The following loop is over the sides of the element.
If the element has no neighbor on a side then that
side MUST live on a boundary of the domain.
for (unsigned int side=0; side<elem->n_sides(); side++)
if (elem->neighbor(side) == NULL)
{
The value of the shape functions at the quadrature
points.
const std::vector<std::vector<Real> >& phi_face = fe_face->get_phi();
The Jacobian * Quadrature Weight at the quadrature
points on the face.
const std::vector<Real>& JxW_face = fe_face->get_JxW();
The XYZ locations (in physical space) of the
quadrature points on the face. This is where
we will interpolate the boundary value function.
const std::vector<Point >& qface_point = fe_face->get_xyz();
Compute the shape function values on the element
face.
fe_face->reinit(elem, side);
Loop over the face quadrature points for integration.
for (unsigned int qp=0; qp<qface.n_points(); qp++)
{
The location on the boundary of the current
face quadrature point.
const Real xf = qface_point[qp](0);
const Real yf = qface_point[qp](1);
The penalty value. \frac{1}{\epsilon}
in the discussion above.
const Real penalty = 1.e10;
The boundary value.
const Real value = exact_solution(xf, yf);
Matrix contribution of the L2 projection.
for (unsigned int i=0; i<phi_face.size(); i++)
for (unsigned int j=0; j<phi_face.size(); j++)
Ke(i,j) += JxW_face[qp]*penalty*phi_face[i][qp]*phi_face[j][qp];
Right-hand-side contribution of the L2
projection.
for (unsigned int i=0; i<phi_face.size(); i++)
Fe(i) += JxW_face[qp]*penalty*value*phi_face[i][qp];
}
}
}
We have now finished the quadrature point loop,
and have therefore applied all the boundary conditions.
The element matrix and right-hand-side are now built for this element. Add them to the global matrix and right-hand-side vector. The SparseMatrix::add_matrix() and NumericVector::add_vector() members do this for us.
The element matrix and right-hand-side are now built for this element. Add them to the global matrix and right-hand-side vector. The SparseMatrix::add_matrix() and NumericVector::add_vector() members do this for us.
system.matrix->add_matrix (Ke, dof_indices);
system.rhs->add_vector (Fe, dof_indices);
}
All done!
}
The program without comments:
#include <iostream>
#include <algorithm>
#include <math.h>
#include "libmesh.h"
#include "mesh.h"
#include "mesh_generation.h"
#include "gmv_io.h"
#include "linear_implicit_system.h"
#include "equation_systems.h"
#include "fe.h"
#include "quadrature_gauss.h"
#include "sparse_matrix.h"
#include "numeric_vector.h"
#include "dense_matrix.h"
#include "dense_vector.h"
#include "elem.h"
#include "dof_map.h"
#include "boundary_mesh.h"
#include "boundary_info.h"
void assemble_poisson(EquationSystems& es,
const std::string& system_name);
Real exact_solution (const Real x,
const Real y,
const Real z = 0.);
int main (int argc, char** argv)
{
libMesh::init (argc, argv);
{
std::cout << "Running " << argv[0];
for (int i=1; i<argc; i++)
std::cout << " " << argv[i];
std::cout << std::endl << std::endl;
Mesh mesh (2);
MeshTools::Generation::build_square (mesh,
15, 15,
-1., 1.,
-1., 1.,
QUAD9);
mesh.print_info();
EquationSystems equation_systems (mesh);
equation_systems.add_system<LinearImplicitSystem> ("Poisson");
equation_systems.get_system("Poisson").add_variable("u", SECOND);
equation_systems.get_system("Poisson").attach_assemble_function (assemble_poisson);
equation_systems.init();
equation_systems.print_info();
equation_systems.get_system("Poisson").solve();
GMVIO (mesh).write_equation_systems ("out.gmv", equation_systems);
for (unsigned int i=0; i<3; i++)
{
here();
equation_systems.get_system("Poisson").solve();
BoundaryMesh boundary_mesh (mesh.mesh_dimension()-1);
mesh.boundary_info->sync(boundary_mesh, false);
}
}
return libMesh::close();
}
void assemble_poisson(EquationSystems& es,
const std::string& system_name)
{
assert (system_name == "Poisson");
const Mesh& mesh = es.get_mesh();
const unsigned int dim = mesh.mesh_dimension();
LinearImplicitSystem& system = es.get_system<LinearImplicitSystem> ("Poisson");
const DofMap& dof_map = system.get_dof_map();
FEType fe_type = dof_map.variable_type(0);
AutoPtr<FEBase> fe (FEBase::build(dim, fe_type));
QGauss qrule (dim, FIFTH);
fe->attach_quadrature_rule (&qrule);
AutoPtr<FEBase> fe_face (FEBase::build(dim, fe_type));
QGauss qface(dim-1, FIFTH);
fe_face->attach_quadrature_rule (&qface);
const std::vector<Real>& JxW = fe->get_JxW();
const std::vector<Point>& q_point = fe->get_xyz();
const std::vector<std::vector<Real> >& phi = fe->get_phi();
const std::vector<std::vector<RealGradient> >& dphi = fe->get_dphi();
DenseMatrix<Number> Ke;
DenseVector<Number> Fe;
std::vector<unsigned int> dof_indices;
MeshBase::const_element_iterator el = mesh.elements_begin();
const MeshBase::const_element_iterator end_el = mesh.elements_end();
for ( ; el != end_el ; ++el)
{
const Elem* elem = *el;
dof_map.dof_indices (elem, dof_indices);
fe->reinit (elem);
Ke.resize (dof_indices.size(),
dof_indices.size());
Fe.resize (dof_indices.size());
for (unsigned int qp=0; qp<qrule.n_points(); qp++)
{
for (unsigned int i=0; i<phi.size(); i++)
for (unsigned int j=0; j<phi.size(); j++)
{
Ke(i,j) += JxW[qp]*(dphi[i][qp]*dphi[j][qp]);
}
{
const Real x = q_point[qp](0);
const Real y = q_point[qp](1);
const Real eps = 1.e-3;
const Real fxy = -(exact_solution(x,y-eps) +
exact_solution(x,y+eps) +
exact_solution(x-eps,y) +
exact_solution(x+eps,y) -
4.*exact_solution(x,y))/eps/eps;
for (unsigned int i=0; i<phi.size(); i++)
Fe(i) += JxW[qp]*fxy*phi[i][qp];
}
}
{
for (unsigned int side=0; side<elem->n_sides(); side++)
if (elem->neighbor(side) == NULL)
{
const std::vector<std::vector<Real> >& phi_face = fe_face->get_phi();
const std::vector<Real>& JxW_face = fe_face->get_JxW();
const std::vector<Point >& qface_point = fe_face->get_xyz();
fe_face->reinit(elem, side);
for (unsigned int qp=0; qp<qface.n_points(); qp++)
{
const Real xf = qface_point[qp](0);
const Real yf = qface_point[qp](1);
const Real penalty = 1.e10;
const Real value = exact_solution(xf, yf);
for (unsigned int i=0; i<phi_face.size(); i++)
for (unsigned int j=0; j<phi_face.size(); j++)
Ke(i,j) += JxW_face[qp]*penalty*phi_face[i][qp]*phi_face[j][qp];
for (unsigned int i=0; i<phi_face.size(); i++)
Fe(i) += JxW_face[qp]*penalty*value*phi_face[i][qp];
}
}
}
system.matrix->add_matrix (Ke, dof_indices);
system.rhs->add_vector (Fe, dof_indices);
}
}
The console output of the program:
***************************************************************
* Running Example ./ex3-devel
***************************************************************
Running ./ex3-devel
Mesh Information:
mesh_dimension()=2
spatial_dimension()=3
n_nodes()=961
n_elem()=225
n_local_elem()=225
n_active_elem()=225
n_subdomains()=1
n_processors()=1
processor_id()=0
EquationSystems
n_systems()=1
System "Poisson"
Type "LinearImplicit"
Variables="u"
Finite Element Types="LAGRANGE"
Approximation Orders="SECOND"
n_dofs()=961
n_local_dofs()=961
n_constrained_dofs()=0
n_vectors()=1
[0] ex3.C, line 162, compiled Jun 6 2007 at 11:54:25
[0] ex3.C, line 162, compiled Jun 6 2007 at 11:54:25
[0] ex3.C, line 162, compiled Jun 6 2007 at 11:54:25
***************************************************************
* Done Running Example ./ex3-devel
***************************************************************
Example 3 - Solving a Poisson Problem
This is the third example program. It builds on the second example program by showing how to solve a simple Poisson system. This example also introduces the notion of customized matrix assembly functions, working with an exact solution, and using element iterators. We will not comment on things that were already explained in the second example.
C++ include files that we need