The source file exact_solution.C with comments:
#include <math.h>
Mesh library includes
#include "libmesh/libmesh_common.h"
Bring in everything from the libMesh namespace
using namespace libMesh;
/**
* 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".
*/
Real exact_solution (const Real x,
const Real y,
const Real z = 0.)
{
static const Real pi = acos(-1.);
return cos(.5*pi*x)*sin(.5*pi*y)*cos(.5*pi*z);
}
The source file subdomains_ex2.C with comments:
Subdomains Example 2 - Subdomain-Restricted Variables
This example builds on the fourth example program by showing how to restrict solution fields to a subdomain (or union of subdomains).
C++ include files that we need
#include <iostream>
#include <algorithm>
#include <math.h>
Basic include file needed for the mesh functionality.
#include "libmesh/libmesh.h"
#include "libmesh/mesh.h"
#include "libmesh/mesh_generation.h"
#include "libmesh/exodusII_io.h"
#include "libmesh/gnuplot_io.h"
#include "libmesh/linear_implicit_system.h"
#include "libmesh/equation_systems.h"
Define the Finite Element object.
#include "libmesh/fe.h"
Define Gauss quadrature rules.
#include "libmesh/quadrature_gauss.h"
Define the DofMap, which handles degree of freedom
indexing.
#include "libmesh/dof_map.h"
Define useful datatypes for finite element
matrix and vector components.
#include "libmesh/sparse_matrix.h"
#include "libmesh/numeric_vector.h"
#include "libmesh/dense_matrix.h"
#include "libmesh/dense_vector.h"
Define the PerfLog, a performance logging utility.
It is useful for timing events in a code and giving
you an idea where bottlenecks lie.
#include "libmesh/perf_log.h"
The definition of a geometric element
#include "libmesh/elem.h"
#include "libmesh/string_to_enum.h"
#include "libmesh/getpot.h"
Bring in everything from the libMesh namespace
using namespace libMesh;
Function prototype. This is the function that will assemble
the linear system for our Poisson problem. Note that the
function will take the \p EquationSystems object and the
name of the system we are assembling as input. From the
\p EquationSystems object we have acess to the \p Mesh and
other objects we might need.
void assemble_poisson(EquationSystems& es,
const std::string& system_name);
Exact solution function prototype.
Real exact_solution (const Real x,
const Real y = 0.,
const Real z = 0.);
Begin the main program.
int main (int argc, char** argv)
{
Initialize libMesh and any dependent libaries, like in example 2.
LibMeshInit init (argc, argv);
Declare a performance log for the main program
PerfLog perf_main("Main Program");
Create a GetPot object to parse the command line
Create a GetPot object to parse the command line
GetPot command_line (argc, argv);
Check for proper calling arguments.
if (argc < 3)
{
if (libMesh::processor_id() == 0)
std::cerr << "Usage:\n"
<<"\t " << argv[0] << " -d 2(3)" << " -n 15"
<< std::endl;
This handy function will print the file name, line number,
and then abort. Currrently the library does not use C++
exception handling.
libmesh_error();
}
Brief message to the user regarding the program name
and command line arguments.
else
{
std::cout << "Running " << argv[0];
for (int i=1; i<argc; i++)
std::cout << " " << argv[i];
std::cout << std::endl << std::endl;
}
Read problem dimension from command line. Use int
instead of unsigned since the GetPot overload is ambiguous
otherwise.
int dim = 2;
if ( command_line.search(1, "-d") )
dim = command_line.next(dim);
Skip higher-dimensional examples on a lower-dimensional libMesh build
libmesh_example_assert(dim <= LIBMESH_DIM, "2D/3D support");
Create a mesh with user-defined dimension.
Mesh mesh (dim);
Read number of elements from command line
int ps = 15;
if ( command_line.search(1, "-n") )
ps = command_line.next(ps);
Read FE order from command line
std::string order = "SECOND";
if ( command_line.search(2, "-Order", "-o") )
order = command_line.next(order);
Read FE Family from command line
std::string family = "LAGRANGE";
if ( command_line.search(2, "-FEFamily", "-f") )
family = command_line.next(family);
Cannot use discontinuous basis.
if ((family == "MONOMIAL") || (family == "XYZ"))
{
if (libMesh::processor_id() == 0)
std::cerr << "ex28 currently requires a C^0 (or higher) FE basis." << std::endl;
libmesh_error();
}
Use the MeshTools::Generation mesh generator to create a uniform
grid on the square [-1,1]^D. We instruct the mesh generator
to build a mesh of 8x8 \p Quad9 elements in 2D, or \p Hex27
elements in 3D. Building these higher-order elements allows
us to use higher-order approximation, as in example 3.
Real halfwidth = dim > 1 ? 1. : 0.;
Real halfheight = dim > 2 ? 1. : 0.;
if ((family == "LAGRANGE") && (order == "FIRST"))
{
No reason to use high-order geometric elements if we are
solving with low-order finite elements.
MeshTools::Generation::build_cube (mesh,
ps,
(dim>1) ? ps : 0,
(dim>2) ? ps : 0,
-1., 1.,
-halfwidth, halfwidth,
-halfheight, halfheight,
(dim==1) ? EDGE2 :
((dim == 2) ? QUAD4 : HEX8));
}
else
{
MeshTools::Generation::build_cube (mesh,
ps,
(dim>1) ? ps : 0,
(dim>2) ? ps : 0,
-1., 1.,
-halfwidth, halfwidth,
-halfheight, halfheight,
(dim==1) ? EDGE3 :
((dim == 2) ? QUAD9 : HEX27));
}
{
MeshBase::element_iterator el = mesh.elements_begin();
const MeshBase::element_iterator end_el = mesh.elements_end();
for ( ; el != end_el; ++el)
{
Elem* elem = *el;
const Point cent = elem->centroid();
if (dim > 1)
{
if ((cent(0) > 0) == (cent(1) > 0))
elem->subdomain_id() = 1;
}
else
{
if (cent(0) > 0)
elem->subdomain_id() = 1;
}
}
}
Print information about the mesh to the screen.
mesh.print_info();
Create an equation systems object.
EquationSystems equation_systems (mesh);
Declare the system and its variables.
Create a system named "Poisson"
LinearImplicitSystem& system =
equation_systems.add_system<LinearImplicitSystem> ("Poisson");
std::set<subdomain_id_type> active_subdomains;
Add the variable "u" to "Poisson". "u"
will be approximated using second-order approximation.
active_subdomains.clear(); active_subdomains.insert(0);
system.add_variable("u",
Utility::string_to_enum<Order> (order),
Utility::string_to_enum<FEFamily>(family),
&active_subdomains);
Add the variable "v" to "Poisson". "v"
will be approximated using second-order approximation.
active_subdomains.clear(); active_subdomains.insert(1);
system.add_variable("v",
Utility::string_to_enum<Order> (order),
Utility::string_to_enum<FEFamily>(family),
&active_subdomains);
Give the system a pointer to the matrix assembly
function.
system.attach_assemble_function (assemble_poisson);
Initialize the data structures for the equation system.
equation_systems.init();
Print information about the system to the screen.
equation_systems.print_info();
mesh.print_info();
Solve the system "Poisson", just like example 2.
equation_systems.get_system("Poisson").solve();
After solving the system write the solution
to a GMV-formatted plot file.
if(dim == 1)
{
GnuPlotIO plot(mesh,"Subdomains Example 2, 1D",GnuPlotIO::GRID_ON);
plot.write_equation_systems("gnuplot_script",equation_systems);
}
else
{
#ifdef LIBMESH_HAVE_EXODUS_API
ExodusII_IO (mesh).write_equation_systems ((dim == 3) ?
"out_3.e" : "out_2.e",equation_systems);
#endif // #ifdef LIBMESH_HAVE_EXODUS_API
}
All done.
return 0;
}
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.
libmesh_assert_equal_to (system_name, "Poisson");
Declare a performance log. Give it a descriptive
string to identify what part of the code we are
logging, since there may be many PerfLogs in an
application.
PerfLog perf_log ("Matrix Assembly");
Get a constant reference to the mesh object.
const MeshBase& 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 \p DofMap object for this system. The \p DofMap
object handles the index translation from node and element numbers
to degree of freedom numbers. We will talk more about the \p 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
\p FEBase::build() member dynamically creates memory we will
store the object as an \p AutoPtr. This can be thought
of as a pointer that will clean up after itself.
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 finte 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.
We begin with 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". More detail is in example 3.
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<dof_id_type> dof_indices, dof_indices2;
Now we will loop over all the elements in the mesh.
We will compute the element matrix and right-hand-side
contribution. See example 3 for a discussion of the
element iterators. Here we use the \p const_local_elem_iterator
to indicate we only want to loop over elements that are assigned
to the local processor. This allows each processor to compute
its components of the global matrix.
"PARALLEL CHANGE"
"PARALLEL CHANGE"
MeshBase::const_element_iterator el = mesh.local_elements_begin();
const MeshBase::const_element_iterator end_el = mesh.local_elements_end();
for ( ; el != end_el; ++el)
{
Start logging the shape function initialization.
This is done through a simple function call with
the name of the event to log.
perf_log.push("elem init");
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,0);
dof_map.dof_indices (elem, dof_indices2,1);
std::cout << "dof_indices.size()="
<< dof_indices.size()
<< ", dof_indices2.size()="
<< dof_indices2.size()
<< std::endl;
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.
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).
Ke.resize (std::max(dof_indices.size(), dof_indices2.size()),
std::max(dof_indices.size(), dof_indices2.size()));
Fe.resize (std::max(dof_indices.size(), dof_indices2.size()));
Stop logging the shape function initialization.
If you forget to stop logging an event the PerfLog
object will probably catch the error and abort.
perf_log.pop("elem init");
Now we will build the element matrix. This involves
a double loop to integrate the test funcions (i) against
the trial functions (j).
We have split the numeric integration into two loops so that we can log the matrix and right-hand-side computation seperately.
Now start logging the element matrix computation
We have split the numeric integration into two loops so that we can log the matrix and right-hand-side computation seperately.
Now start logging the element matrix computation
perf_log.push ("Ke");
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]);
Stop logging the matrix computation
perf_log.pop ("Ke");
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.
Start logging the right-hand-side computation
Start logging the right-hand-side computation
perf_log.push ("Fe");
for (unsigned int qp=0; qp<qrule.n_points(); qp++)
{
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 on a structured grid is
u_xx + u_yy = (u(i-1,j) + u(i+1,j) + u(i,j-1) + u(i,j+1) + -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 on a structured grid is
u_xx + u_yy = (u(i-1,j) + u(i+1,j) + u(i,j-1) + u(i,j+1) + -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 x = q_point[qp](0);
#if LIBMESH_DIM > 1
const Real y = q_point[qp](1);
#else
const Real y = 0;
#endif
#if LIBMESH_DIM > 2
const Real z = q_point[qp](2);
#else
const Real z = 0;
#endif
const Real eps = 1.e-3;
const Real uxx = (exact_solution(x-eps,y,z) +
exact_solution(x+eps,y,z) +
-2.*exact_solution(x,y,z))/eps/eps;
const Real uyy = (exact_solution(x,y-eps,z) +
exact_solution(x,y+eps,z) +
-2.*exact_solution(x,y,z))/eps/eps;
const Real uzz = (exact_solution(x,y,z-eps) +
exact_solution(x,y,z+eps) +
-2.*exact_solution(x,y,z))/eps/eps;
Real fxy;
if(dim==1)
{
In 1D, compute the rhs by differentiating the
exact solution twice.
const Real pi = libMesh::pi;
fxy = (0.25*pi*pi)*sin(.5*pi*x);
}
else
{
fxy = - (uxx + uyy + ((dim==2) ? 0. : uzz));
}
Add the RHS contribution
for (unsigned int i=0; i<phi.size(); i++)
Fe(i) += JxW[qp]*fxy*phi[i][qp];
}
Stop logging the right-hand-side computation
perf_log.pop ("Fe");
At this point 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 imposed
via the penalty method. This is discussed at length in
example 3.
{
Start logging the boundary condition computation
perf_log.push ("BCs");
The following loops 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) ||
(elem->neighbor(side)->subdomain_id() != elem->subdomain_id()))
{
The penalty value. \frac{1}{\epsilon}
in the discussion above.
const Real penalty = 1.e10;
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);
#if LIBMESH_DIM > 1
const Real yf = qface_point[qp](1);
#else
const Real yf = 0.;
#endif
#if LIBMESH_DIM > 2
const Real zf = qface_point[qp](2);
#else
const Real zf = 0.;
#endif
The boundary value.
const Real value = exact_solution(xf, yf, zf);
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];
}
}
Stop logging the boundary condition computation
perf_log.pop ("BCs");
}
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 \p PetscMatrix::add_matrix()
and \p PetscVector::add_vector() members do this for us.
Start logging the insertion of the local (element)
matrix and vector into the global matrix and vector
perf_log.push ("matrix insertion");
if (dof_indices.size())
{
system.matrix->add_matrix (Ke, dof_indices);
system.rhs->add_vector (Fe, dof_indices);
}
if (dof_indices2.size())
{
system.matrix->add_matrix (Ke, dof_indices2);
system.rhs->add_vector (Fe, dof_indices2);
}
Start logging the insertion of the local (element)
matrix and vector into the global matrix and vector
perf_log.pop ("matrix insertion");
}
That's it. We don't need to do anything else to the
PerfLog. When it goes out of scope (at this function return)
it will print its log to the screen. Pretty easy, huh?
}
The source file exact_solution.C without comments:
#include <math.h>
#include "libmesh/libmesh_common.h"
using namespace libMesh;
/**
* 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".
*/
Real exact_solution (const Real x,
const Real y,
const Real z = 0.)
{
static const Real pi = acos(-1.);
return cos(.5*pi*x)*sin(.5*pi*y)*cos(.5*pi*z);
}
The source file subdomains_ex2.C without comments:
#include <iostream>
#include <algorithm>
#include <math.h>
#include "libmesh/libmesh.h"
#include "libmesh/mesh.h"
#include "libmesh/mesh_generation.h"
#include "libmesh/exodusII_io.h"
#include "libmesh/gnuplot_io.h"
#include "libmesh/linear_implicit_system.h"
#include "libmesh/equation_systems.h"
#include "libmesh/fe.h"
#include "libmesh/quadrature_gauss.h"
#include "libmesh/dof_map.h"
#include "libmesh/sparse_matrix.h"
#include "libmesh/numeric_vector.h"
#include "libmesh/dense_matrix.h"
#include "libmesh/dense_vector.h"
#include "libmesh/perf_log.h"
#include "libmesh/elem.h"
#include "libmesh/string_to_enum.h"
#include "libmesh/getpot.h"
using namespace libMesh;
void assemble_poisson(EquationSystems& es,
const std::string& system_name);
Real exact_solution (const Real x,
const Real y = 0.,
const Real z = 0.);
int main (int argc, char** argv)
{
LibMeshInit init (argc, argv);
GetPot command_line (argc, argv);
if (argc < 3)
{
if (libMesh::processor_id() == 0)
std::cerr << "Usage:\n"
<<"\t " << argv[0] << " -d 2(3)" << " -n 15"
<< std::endl;
libmesh_error();
}
else
{
std::cout << "Running " << argv[0];
for (int i=1; i<argc; i++)
std::cout << " " << argv[i];
std::cout << std::endl << std::endl;
}
int dim = 2;
if ( command_line.search(1, "-d") )
dim = command_line.next(dim);
libmesh_example_assert(dim <= LIBMESH_DIM, "2D/3D support");
Mesh mesh (dim);
int ps = 15;
if ( command_line.search(1, "-n") )
ps = command_line.next(ps);
std::string order = "SECOND";
if ( command_line.search(2, "-Order", "-o") )
order = command_line.next(order);
std::string family = "LAGRANGE";
if ( command_line.search(2, "-FEFamily", "-f") )
family = command_line.next(family);
if ((family == "MONOMIAL") || (family == "XYZ"))
{
if (libMesh::processor_id() == 0)
std::cerr << "ex28 currently requires a C^0 (or higher) FE basis." << std::endl;
libmesh_error();
}
Real halfwidth = dim > 1 ? 1. : 0.;
Real halfheight = dim > 2 ? 1. : 0.;
if ((family == "LAGRANGE") && (order == "FIRST"))
{
MeshTools::Generation::build_cube (mesh,
ps,
(dim>1) ? ps : 0,
(dim>2) ? ps : 0,
-1., 1.,
-halfwidth, halfwidth,
-halfheight, halfheight,
(dim==1) ? EDGE2 :
((dim == 2) ? QUAD4 : HEX8));
}
else
{
MeshTools::Generation::build_cube (mesh,
ps,
(dim>1) ? ps : 0,
(dim>2) ? ps : 0,
-1., 1.,
-halfwidth, halfwidth,
-halfheight, halfheight,
(dim==1) ? EDGE3 :
((dim == 2) ? QUAD9 : HEX27));
}
{
MeshBase::element_iterator el = mesh.elements_begin();
const MeshBase::element_iterator end_el = mesh.elements_end();
for ( ; el != end_el; ++el)
{
Elem* elem = *el;
const Point cent = elem->centroid();
if (dim > 1)
{
if ((cent(0) > 0) == (cent(1) > 0))
elem->subdomain_id() = 1;
}
else
{
if (cent(0) > 0)
elem->subdomain_id() = 1;
}
}
}
mesh.print_info();
EquationSystems equation_systems (mesh);
LinearImplicitSystem& system =
equation_systems.add_system<LinearImplicitSystem> ("Poisson");
std::set<subdomain_id_type> active_subdomains;
active_subdomains.clear(); active_subdomains.insert(0);
system.add_variable("u",
Utility::string_to_enum<Order> (order),
Utility::string_to_enum<FEFamily>(family),
&active_subdomains);
active_subdomains.clear(); active_subdomains.insert(1);
system.add_variable("v",
Utility::string_to_enum<Order> (order),
Utility::string_to_enum<FEFamily>(family),
&active_subdomains);
system.attach_assemble_function (assemble_poisson);
equation_systems.init();
equation_systems.print_info();
mesh.print_info();
equation_systems.get_system("Poisson").solve();
if(dim == 1)
{
GnuPlotIO plot(mesh,"Subdomains Example 2, 1D",GnuPlotIO::GRID_ON);
plot.write_equation_systems("gnuplot_script",equation_systems);
}
else
{
#ifdef LIBMESH_HAVE_EXODUS_API
ExodusII_IO (mesh).write_equation_systems ((dim == 3) ?
"out_3.e" : "out_2.e",equation_systems);
#endif // #ifdef LIBMESH_HAVE_EXODUS_API
}
return 0;
}
void assemble_poisson(EquationSystems& es,
const std::string& system_name)
{
libmesh_assert_equal_to (system_name, "Poisson");
PerfLog perf_log ("Matrix Assembly");
const MeshBase& 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<dof_id_type> dof_indices, dof_indices2;
MeshBase::const_element_iterator el = mesh.local_elements_begin();
const MeshBase::const_element_iterator end_el = mesh.local_elements_end();
for ( ; el != end_el; ++el)
{
perf_log.push("elem init");
const Elem* elem = *el;
dof_map.dof_indices (elem, dof_indices,0);
dof_map.dof_indices (elem, dof_indices2,1);
fe->reinit (elem);
Ke.resize (std::max(dof_indices.size(), dof_indices2.size()),
std::max(dof_indices.size(), dof_indices2.size()));
Fe.resize (std::max(dof_indices.size(), dof_indices2.size()));
perf_log.pop("elem init");
perf_log.push ("Ke");
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]);
perf_log.pop ("Ke");
perf_log.push ("Fe");
for (unsigned int qp=0; qp<qrule.n_points(); qp++)
{
const Real x = q_point[qp](0);
#if LIBMESH_DIM > 1
const Real y = q_point[qp](1);
#else
const Real y = 0;
#endif
#if LIBMESH_DIM > 2
const Real z = q_point[qp](2);
#else
const Real z = 0;
#endif
const Real eps = 1.e-3;
const Real uxx = (exact_solution(x-eps,y,z) +
exact_solution(x+eps,y,z) +
-2.*exact_solution(x,y,z))/eps/eps;
const Real uyy = (exact_solution(x,y-eps,z) +
exact_solution(x,y+eps,z) +
-2.*exact_solution(x,y,z))/eps/eps;
const Real uzz = (exact_solution(x,y,z-eps) +
exact_solution(x,y,z+eps) +
-2.*exact_solution(x,y,z))/eps/eps;
Real fxy;
if(dim==1)
{
const Real pi = libMesh::pi;
fxy = (0.25*pi*pi)*sin(.5*pi*x);
}
else
{
fxy = - (uxx + uyy + ((dim==2) ? 0. : uzz));
}
for (unsigned int i=0; i<phi.size(); i++)
Fe(i) += JxW[qp]*fxy*phi[i][qp];
}
perf_log.pop ("Fe");
{
perf_log.push ("BCs");
for (unsigned int side=0; side<elem->n_sides(); side++)
if ((elem->neighbor(side) == NULL) ||
(elem->neighbor(side)->subdomain_id() != elem->subdomain_id()))
{
const Real penalty = 1.e10;
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);
#if LIBMESH_DIM > 1
const Real yf = qface_point[qp](1);
#else
const Real yf = 0.;
#endif
#if LIBMESH_DIM > 2
const Real zf = qface_point[qp](2);
#else
const Real zf = 0.;
#endif
const Real value = exact_solution(xf, yf, zf);
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];
}
}
perf_log.pop ("BCs");
}
perf_log.push ("matrix insertion");
if (dof_indices.size())
{
system.matrix->add_matrix (Ke, dof_indices);
system.rhs->add_vector (Fe, dof_indices);
}
if (dof_indices2.size())
{
system.matrix->add_matrix (Ke, dof_indices2);
system.rhs->add_vector (Fe, dof_indices2);
}
perf_log.pop ("matrix insertion");
}
}
The console output of the program:
***************************************************************
* Running Example subdomains_ex2:
* mpirun -np 2 example-devel -d 1 -n 20 -pc_type bjacobi -sub_pc_type ilu -sub_pc_factor_levels 4 -sub_pc_factor_zeropivot 0 -ksp_right_pc -log_summary
***************************************************************
Running /workspace/libmesh/examples/subdomains/subdomains_ex2/.libs/lt-example-devel -d 1 -n 20 -pc_type bjacobi -sub_pc_type ilu -sub_pc_factor_levels 4 -sub_pc_factor_zeropivot 0 -ksp_right_pc -log_summary
Mesh Information:
mesh_dimension()=1
spatial_dimension()=3
n_nodes()=41
n_local_nodes()=21
n_elem()=20
n_local_elem()=10
n_active_elem()=20
n_subdomains()=2
n_partitions()=2
n_processors()=2
n_threads()=1
processor_id()=0
EquationSystems
n_systems()=1
System #0, "Poisson"
Type "LinearImplicit"
Variables="u" "v"
Finite Element Types="LAGRANGE", "JACOBI_20_00" "LAGRANGE", "JACOBI_20_00"
Infinite Element Mapping="CARTESIAN" "CARTESIAN"
Approximation Orders="SECOND", "THIRD" "SECOND", "THIRD"
n_dofs()=42
n_local_dofs()=22
n_constrained_dofs()=0
n_local_constrained_dofs()=0
n_vectors()=1
n_matrices()=1
DofMap Sparsity
Average On-Processor Bandwidth <= 3.7619
Average Off-Processor Bandwidth <= 0.0952381
Maximum On-Processor Bandwidth <= 5
Maximum Off-Processor Bandwidth <= 2
DofMap Constraints
Number of DoF Constraints = 0
Number of Node Constraints = 0
Mesh Information:
mesh_dimension()=1
spatial_dimension()=3
n_nodes()=41
n_local_nodes()=21
n_elem()=20
n_local_elem()=10
n_active_elem()=20
n_subdomains()=2
n_partitions()=2
n_processors()=2
n_threads()=1
processor_id()=0
----------------------------------------------------------------------------------------------------------------------
| Processor id: 0 |
| Num Processors: 2 |
| Time: Tue Feb 5 13:42:09 2013 |
| OS: Linux |
| HostName: hbar.ices.utexas.edu |
| OS Release: 2.6.32-279.1.1.el6.x86_64 |
| OS Version: #1 SMP Tue Jul 10 11:24:23 CDT 2012 |
| Machine: x86_64 |
| Username: benkirk |
| Configuration: ./configure '--enable-everything' |
| '--prefix=/workspace/libmesh/install' |
| 'CXX=mpicxx' |
| 'CC=mpicc' |
| 'F77=mpif77' |
| 'FC=mpif90' |
| 'PETSC_DIR=/opt/apps/ossw/libraries/petsc/petsc-3.3-p2' |
| 'PETSC_ARCH=intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt' |
| 'SLEPC_DIR=/opt/apps/ossw/libraries/slepc/slepc-3.3-p2-petsc-3.3-p2-cxx-opt' |
| 'TRILINOS_DIR=/opt/apps/ossw/libraries/trilinos/trilinos-10.12.2/sl6/intel-12.1/mpich2-1.4.1p1/mkl-intel-10.3.12.361'|
| 'VTK_DIR=/opt/apps/ossw/libraries/vtk/vtk-5.10.0/sl6/intel-12.1' |
----------------------------------------------------------------------------------------------------------------------
-----------------------------------------------------------------------------------------------------------
| Matrix Assembly Performance: Alive time=0.001499, Active time=0.000965 |
-----------------------------------------------------------------------------------------------------------
| Event nCalls Total Time Avg Time Total Time Avg Time % of Active Time |
| w/o Sub w/o Sub With Sub With Sub w/o S With S |
|-----------------------------------------------------------------------------------------------------------|
| |
| BCs 10 0.0001 0.000014 0.0001 0.000014 14.92 14.92 |
| Fe 10 0.0001 0.000007 0.0001 0.000007 7.15 7.15 |
| Ke 10 0.0000 0.000003 0.0000 0.000003 2.59 2.59 |
| elem init 10 0.0007 0.000066 0.0007 0.000066 68.81 68.81 |
| matrix insertion 10 0.0001 0.000006 0.0001 0.000006 6.53 6.53 |
-----------------------------------------------------------------------------------------------------------
| Totals: 50 0.0010 100.00 |
-----------------------------------------------------------------------------------------------------------
************************************************************************************************************************
*** WIDEN YOUR WINDOW TO 120 CHARACTERS. Use 'enscript -r -fCourier9' to print this document ***
************************************************************************************************************************
---------------------------------------------- PETSc Performance Summary: ----------------------------------------------
/workspace/libmesh/examples/subdomains/subdomains_ex2/.libs/lt-example-devel on a intel-12. named hbar.ices.utexas.edu with 2 processors, by benkirk Tue Feb 5 13:42:09 2013
Using Petsc Release Version 3.3.0, Patch 2, Fri Jul 13 15:42:00 CDT 2012
Max Max/Min Avg Total
Time (sec): 2.702e-02 1.00004 2.702e-02
Objects: 4.300e+01 1.00000 4.300e+01
Flops: 3.018e+03 1.07021 2.919e+03 5.838e+03
Flops/sec: 1.117e+05 1.07017 1.080e+05 2.161e+05
MPI Messages: 1.150e+01 1.00000 1.150e+01 2.300e+01
MPI Message Lengths: 1.300e+02 1.00000 1.130e+01 2.600e+02
MPI Reductions: 5.300e+01 1.00000
Flop counting convention: 1 flop = 1 real number operation of type (multiply/divide/add/subtract)
e.g., VecAXPY() for real vectors of length N --> 2N flops
and VecAXPY() for complex vectors of length N --> 8N flops
Summary of Stages: ----- Time ------ ----- Flops ----- --- Messages --- -- Message Lengths -- -- Reductions --
Avg %Total Avg %Total counts %Total Avg %Total counts %Total
0: Main Stage: 2.6988e-02 99.9% 5.8380e+03 100.0% 2.300e+01 100.0% 1.130e+01 100.0% 5.200e+01 98.1%
------------------------------------------------------------------------------------------------------------------------
See the 'Profiling' chapter of the users' manual for details on interpreting output.
Phase summary info:
Count: number of times phase was executed
Time and Flops: Max - maximum over all processors
Ratio - ratio of maximum to minimum over all processors
Mess: number of messages sent
Avg. len: average message length
Reduct: number of global reductions
Global: entire computation
Stage: stages of a computation. Set stages with PetscLogStagePush() and PetscLogStagePop().
%T - percent time in this phase %f - percent flops in this phase
%M - percent messages in this phase %L - percent message lengths in this phase
%R - percent reductions in this phase
Total Mflop/s: 10e-6 * (sum of flops over all processors)/(max time over all processors)
------------------------------------------------------------------------------------------------------------------------
Event Count Time (sec) Flops --- Global --- --- Stage --- Total
Max Ratio Max Ratio Max Ratio Mess Avg len Reduct %T %f %M %L %R %T %f %M %L %R Mflop/s
------------------------------------------------------------------------------------------------------------------------
--- Event Stage 0: Main Stage
VecMDot 3 1.0 2.4796e-05 1.0 2.58e+02 1.1 0.0e+00 0.0e+00 3.0e+00 0 8 0 0 6 0 8 0 0 6 20
VecNorm 5 1.0 3.2902e-05 1.0 2.20e+02 1.1 0.0e+00 0.0e+00 5.0e+00 0 7 0 0 9 0 7 0 0 10 13
VecScale 4 1.0 1.4782e-05 1.0 8.80e+01 1.1 0.0e+00 0.0e+00 0.0e+00 0 3 0 0 0 0 3 0 0 0 11
VecCopy 2 1.0 3.0994e-06 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00 0 0 0 0 0 0 0 0 0 0 0
VecSet 11 1.0 5.0068e-06 1.2 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00 0 0 0 0 0 0 0 0 0 0 0
VecAXPY 2 1.0 6.1152e-03 1.0 8.80e+01 1.1 0.0e+00 0.0e+00 0.0e+00 23 3 0 0 0 23 3 0 0 0 0
VecMAXPY 4 1.0 1.9073e-06 2.0 3.96e+02 1.1 0.0e+00 0.0e+00 0.0e+00 0 13 0 0 0 0 13 0 0 0 396
VecAssemblyBegin 3 1.0 5.7697e-05 1.0 0.00e+00 0.0 2.0e+00 6.0e+00 9.0e+00 0 0 9 5 17 0 0 9 5 17 0
VecAssemblyEnd 3 1.0 1.8358e-05 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00 0 0 0 0 0 0 0 0 0 0 0
VecScatterBegin 5 1.0 4.0531e-05 1.0 0.00e+00 0.0 1.0e+01 1.4e+01 0.0e+00 0 0 43 55 0 0 0 43 55 0 0
VecScatterEnd 5 1.0 1.0490e-05 1.2 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00 0 0 0 0 0 0 0 0 0 0 0
VecNormalize 4 1.0 7.8917e-05 1.0 2.64e+02 1.1 0.0e+00 0.0e+00 4.0e+00 0 9 0 0 8 0 9 0 0 8 6
MatMult 4 1.0 6.9618e-05 1.0 5.84e+02 1.1 8.0e+00 1.2e+01 0.0e+00 0 19 35 37 0 0 19 35 37 0 16
MatSolve 5 1.0 1.0014e-05 1.1 1.07e+03 1.0 0.0e+00 0.0e+00 0.0e+00 0 36 0 0 0 0 36 0 0 0 209
MatLUFactorNum 1 1.0 2.0027e-05 1.0 3.14e+02 1.0 0.0e+00 0.0e+00 0.0e+00 0 11 0 0 0 0 11 0 0 0 31
MatILUFactorSym 1 1.0 6.6996e-05 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 3.0e+00 0 0 0 0 6 0 0 0 0 6 0
MatAssemblyBegin 2 1.0 3.1996e-04 2.8 0.00e+00 0.0 3.0e+00 1.7e+01 4.0e+00 1 0 13 20 8 1 0 13 20 8 0
MatAssemblyEnd 2 1.0 2.7895e-04 1.0 0.00e+00 0.0 4.0e+00 5.0e+00 8.0e+00 1 0 17 8 15 1 0 17 8 15 0
MatGetRowIJ 1 1.0 5.0068e-06 1.6 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00 0 0 0 0 0 0 0 0 0 0 0
MatGetOrdering 1 1.0 4.9114e-05 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 2.0e+00 0 0 0 0 4 0 0 0 0 4 0
MatZeroEntries 3 1.0 1.8120e-05 2.6 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00 0 0 0 0 0 0 0 0 0 0 0
KSPGMRESOrthog 3 1.0 4.6730e-05 1.0 5.22e+02 1.1 0.0e+00 0.0e+00 3.0e+00 0 17 0 0 6 0 17 0 0 6 21
KSPSetUp 2 1.0 6.9857e-05 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00 0 0 0 0 0 0 0 0 0 0 0
KSPSolve 1 1.0 7.2110e-03 1.0 3.02e+03 1.1 8.0e+00 1.2e+01 1.5e+01 27100 35 37 28 27100 35 37 29 1
PCSetUp 2 1.0 5.0306e-04 1.0 3.14e+02 1.0 0.0e+00 0.0e+00 7.0e+00 2 11 0 0 13 2 11 0 0 13 1
PCSetUpOnBlocks 1 1.0 2.3198e-04 1.0 3.14e+02 1.0 0.0e+00 0.0e+00 5.0e+00 1 11 0 0 9 1 11 0 0 10 3
PCApply 5 1.0 9.1791e-05 1.0 1.07e+03 1.0 0.0e+00 0.0e+00 0.0e+00 0 36 0 0 0 0 36 0 0 0 23
------------------------------------------------------------------------------------------------------------------------
Memory usage is given in bytes:
Object Type Creations Destructions Memory Descendants' Mem.
Reports information only for process 0.
--- Event Stage 0: Main Stage
Vector 24 24 39024 0
Vector Scatter 2 2 2072 0
Index Set 7 7 5296 0
IS L to G Mapping 1 1 564 0
Matrix 4 4 13268 0
Krylov Solver 2 2 19360 0
Preconditioner 2 2 1784 0
Viewer 1 0 0 0
========================================================================================================================
Average time to get PetscTime(): 9.53674e-08
Average time for MPI_Barrier(): 8.10623e-07
Average time for zero size MPI_Send(): 1.29938e-05
#PETSc Option Table entries:
-d 1
-ksp_right_pc
-log_summary
-n 20
-pc_type bjacobi
-sub_pc_factor_levels 4
-sub_pc_factor_zeropivot 0
-sub_pc_type ilu
#End of PETSc Option Table entries
Compiled without FORTRAN kernels
Compiled with full precision matrices (default)
sizeof(short) 2 sizeof(int) 4 sizeof(long) 8 sizeof(void*) 8 sizeof(PetscScalar) 8 sizeof(PetscInt) 4
Configure run at: Thu Nov 8 11:21:02 2012
Configure options: --with-debugging=false --COPTFLAGS=-O3 --CXXOPTFLAGS=-O3 --FOPTFLAGS=-O3 --with-clanguage=C++ --with-shared-libraries=1 --with-mpi-dir=/opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1 --with-mumps=true --download-mumps=1 --with-metis=true --download-metis=1 --with-parmetis=true --download-parmetis=1 --with-superlu=true --download-superlu=1 --with-superludir=true --download-superlu_dist=1 --with-blacs=true --download-blacs=1 --with-scalapack=true --download-scalapack=1 --with-hypre=true --download-hypre=1 --with-blas-lib="[/opt/apps/sysnet/intel/12.1/mkl/10.3.12.361/lib/intel64/libmkl_intel_lp64.so,/opt/apps/sysnet/intel/12.1/mkl/10.3.12.361/lib/intel64/libmkl_sequential.so,/opt/apps/sysnet/intel/12.1/mkl/10.3.12.361/lib/intel64/libmkl_core.so]" --with-lapack-lib="[/opt/apps/sysnet/intel/12.1/mkl/10.3.12.361/lib/intel64/libmkl_lapack95_lp64.a]"
-----------------------------------------
Libraries compiled on Thu Nov 8 11:21:02 2012 on daedalus.ices.utexas.edu
Machine characteristics: Linux-2.6.32-279.1.1.el6.x86_64-x86_64-with-redhat-6.3-Carbon
Using PETSc directory: /opt/apps/ossw/libraries/petsc/petsc-3.3-p2
Using PETSc arch: intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt
-----------------------------------------
Using C compiler: /opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1/bin/mpicxx -wd1572 -O3 -fPIC ${COPTFLAGS} ${CFLAGS}
Using Fortran compiler: /opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1/bin/mpif90 -fPIC -O3 ${FOPTFLAGS} ${FFLAGS}
-----------------------------------------
Using include paths: -I/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt/include -I/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/include -I/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/include -I/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt/include -I/opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1/include
-----------------------------------------
Using C linker: /opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1/bin/mpicxx
Using Fortran linker: /opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1/bin/mpif90
Using libraries: -Wl,-rpath,/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt/lib -L/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt/lib -lpetsc -lX11 -Wl,-rpath,/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt/lib -L/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt/lib -lcmumps -ldmumps -lsmumps -lzmumps -lmumps_common -lpord -lHYPRE -lpthread -lsuperlu_dist_3.0 -lparmetis -lmetis -lscalapack -lblacs -lsuperlu_4.3 -Wl,-rpath,/opt/apps/sysnet/intel/12.1/mkl/10.3.12.361/lib/intel64 -L/opt/apps/sysnet/intel/12.1/mkl/10.3.12.361/lib/intel64 -lmkl_lapack95_lp64 -lmkl_intel_lp64 -lmkl_sequential -lmkl_core -Wl,-rpath,/opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1/lib -L/opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1/lib -Wl,-rpath,/opt/apps/sysnet/intel/12.1/composer_xe_2011_sp1.7.256/compiler/lib/intel64 -L/opt/apps/sysnet/intel/12.1/composer_xe_2011_sp1.7.256/compiler/lib/intel64 -Wl,-rpath,/usr/lib/gcc/x86_64-redhat-linux/4.4.6 -L/usr/lib/gcc/x86_64-redhat-linux/4.4.6 -lmpichf90 -lifport -lifcore -lm -lm -lmpichcxx -ldl -lmpich -lopa -lmpl -lrt -lpthread -limf -lsvml -lipgo -ldecimal -lcilkrts -lstdc++ -lgcc_s -lirc -lirc_s -ldl
-----------------------------------------
----------------------------------------------------------------------------------------------------------------
| libMesh Performance: Alive time=0.092489, Active time=0.02021 |
----------------------------------------------------------------------------------------------------------------
| Event nCalls Total Time Avg Time Total Time Avg Time % of Active Time |
| w/o Sub w/o Sub With Sub With Sub w/o S With S |
|----------------------------------------------------------------------------------------------------------------|
| |
| |
| DofMap |
| add_neighbors_to_send_list() 1 0.0004 0.000383 0.0005 0.000456 1.90 2.26 |
| build_sparsity() 1 0.0006 0.000585 0.0012 0.001195 2.89 5.91 |
| create_dof_constraints() 1 0.0000 0.000012 0.0000 0.000012 0.06 0.06 |
| distribute_dofs() 1 0.0007 0.000668 0.0015 0.001497 3.31 7.41 |
| dof_indices() 42 0.0009 0.000021 0.0009 0.000021 4.43 4.43 |
| prepare_send_list() 1 0.0000 0.000008 0.0000 0.000008 0.04 0.04 |
| reinit() 1 0.0007 0.000666 0.0007 0.000666 3.30 3.30 |
| |
| EquationSystems |
| build_solution_vector() 1 0.0002 0.000200 0.0010 0.001028 0.99 5.09 |
| |
| FE |
| compute_shape_functions() 12 0.0001 0.000006 0.0001 0.000006 0.33 0.33 |
| init_shape_functions() 3 0.0001 0.000034 0.0001 0.000034 0.50 0.50 |
| |
| FEMap |
| compute_affine_map() 12 0.0001 0.000006 0.0001 0.000006 0.34 0.34 |
| compute_face_map() 2 0.0000 0.000005 0.0000 0.000005 0.05 0.05 |
| init_face_shape_functions() 1 0.0000 0.000008 0.0000 0.000008 0.04 0.04 |
| init_reference_to_physical_map() 3 0.0001 0.000022 0.0001 0.000022 0.32 0.32 |
| |
| GnuPlotIO |
| write_nodal_data() 1 0.0007 0.000664 0.0007 0.000664 3.29 3.29 |
| |
| Mesh |
| find_neighbors() 1 0.0007 0.000684 0.0007 0.000735 3.38 3.64 |
| renumber_nodes_and_elem() 2 0.0001 0.000040 0.0001 0.000040 0.40 0.40 |
| |
| MeshCommunication |
| compute_hilbert_indices() 2 0.0003 0.000134 0.0003 0.000134 1.33 1.33 |
| find_global_indices() 2 0.0004 0.000182 0.0015 0.000767 1.81 7.60 |
| parallel_sort() 2 0.0005 0.000267 0.0006 0.000308 2.64 3.05 |
| |
| MeshOutput |
| write_equation_systems() 1 0.0001 0.000088 0.0018 0.001831 0.44 9.06 |
| |
| MeshTools::Generation |
| build_cube() 1 0.0003 0.000334 0.0003 0.000334 1.65 1.65 |
| |
| MetisPartitioner |
| partition() 1 0.0007 0.000686 0.0011 0.001137 3.39 5.63 |
| |
| Parallel |
| allgather() 11 0.0002 0.000016 0.0002 0.000018 0.89 0.98 |
| max(bool) 1 0.0000 0.000007 0.0000 0.000007 0.03 0.03 |
| max(scalar) 115 0.0003 0.000003 0.0003 0.000003 1.46 1.46 |
| max(vector) 26 0.0002 0.000007 0.0004 0.000014 0.84 1.82 |
| min(bool) 133 0.0003 0.000002 0.0003 0.000002 1.62 1.62 |
| min(scalar) 109 0.0011 0.000010 0.0011 0.000010 5.24 5.24 |
| min(vector) 26 0.0002 0.000008 0.0004 0.000015 1.00 1.98 |
| probe() 12 0.0001 0.000009 0.0001 0.000009 0.54 0.54 |
| receive() 12 0.0001 0.000007 0.0002 0.000016 0.43 0.97 |
| send() 12 0.0001 0.000004 0.0001 0.000004 0.25 0.25 |
| send_receive() 16 0.0002 0.000011 0.0005 0.000030 0.91 2.37 |
| sum() 22 0.0001 0.000006 0.0006 0.000026 0.65 2.87 |
| |
| Parallel::Request |
| wait() 12 0.0000 0.000003 0.0000 0.000003 0.21 0.21 |
| |
| Partitioner |
| set_node_processor_ids() 1 0.0003 0.000271 0.0003 0.000347 1.34 1.72 |
| set_parent_processor_ids() 1 0.0001 0.000058 0.0001 0.000058 0.29 0.29 |
| |
| PetscLinearSolver |
| solve() 1 0.0085 0.008509 0.0085 0.008509 42.10 42.10 |
| |
| System |
| assemble() 1 0.0011 0.001089 0.0017 0.001707 5.39 8.45 |
----------------------------------------------------------------------------------------------------------------
| Totals: 606 0.0202 100.00 |
----------------------------------------------------------------------------------------------------------------
***************************************************************
* Done Running Example subdomains_ex2:
* mpirun -np 2 example-devel -d 1 -n 20 -pc_type bjacobi -sub_pc_type ilu -sub_pc_factor_levels 4 -sub_pc_factor_zeropivot 0 -ksp_right_pc -log_summary
***************************************************************
***************************************************************
* Running Example subdomains_ex2:
* mpirun -np 2 example-devel -d 2 -n 15 -pc_type bjacobi -sub_pc_type ilu -sub_pc_factor_levels 4 -sub_pc_factor_zeropivot 0 -ksp_right_pc -log_summary
***************************************************************
Running /workspace/libmesh/examples/subdomains/subdomains_ex2/.libs/lt-example-devel -d 2 -n 15 -pc_type bjacobi -sub_pc_type ilu -sub_pc_factor_levels 4 -sub_pc_factor_zeropivot 0 -ksp_right_pc -log_summary
Mesh Information:
mesh_dimension()=2
spatial_dimension()=3
n_nodes()=961
n_local_nodes()=495
n_elem()=225
n_local_elem()=112
n_active_elem()=225
n_subdomains()=2
n_partitions()=2
n_processors()=2
n_threads()=1
processor_id()=0
EquationSystems
n_systems()=1
System #0, "Poisson"
Type "LinearImplicit"
Variables="u" "v"
Finite Element Types="LAGRANGE", "JACOBI_20_00" "LAGRANGE", "JACOBI_20_00"
Infinite Element Mapping="CARTESIAN" "CARTESIAN"
Approximation Orders="SECOND", "THIRD" "SECOND", "THIRD"
n_dofs()=1022
n_local_dofs()=540
n_constrained_dofs()=0
n_local_constrained_dofs()=0
n_vectors()=1
n_matrices()=1
DofMap Sparsity
Average On-Processor Bandwidth <= 14.1341
Average Off-Processor Bandwidth <= 0.493151
Maximum On-Processor Bandwidth <= 25
Maximum Off-Processor Bandwidth <= 14
DofMap Constraints
Number of DoF Constraints = 0
Number of Node Constraints = 0
Mesh Information:
mesh_dimension()=2
spatial_dimension()=3
n_nodes()=961
n_local_nodes()=495
n_elem()=225
n_local_elem()=112
n_active_elem()=225
n_subdomains()=2
n_partitions()=2
n_processors()=2
n_threads()=1
processor_id()=0
----------------------------------------------------------------------------------------------------------------------
| Processor id: 0 |
| Num Processors: 2 |
| Time: Tue Feb 5 13:42:09 2013 |
| OS: Linux |
| HostName: hbar.ices.utexas.edu |
| OS Release: 2.6.32-279.1.1.el6.x86_64 |
| OS Version: #1 SMP Tue Jul 10 11:24:23 CDT 2012 |
| Machine: x86_64 |
| Username: benkirk |
| Configuration: ./configure '--enable-everything' |
| '--prefix=/workspace/libmesh/install' |
| 'CXX=mpicxx' |
| 'CC=mpicc' |
| 'F77=mpif77' |
| 'FC=mpif90' |
| 'PETSC_DIR=/opt/apps/ossw/libraries/petsc/petsc-3.3-p2' |
| 'PETSC_ARCH=intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt' |
| 'SLEPC_DIR=/opt/apps/ossw/libraries/slepc/slepc-3.3-p2-petsc-3.3-p2-cxx-opt' |
| 'TRILINOS_DIR=/opt/apps/ossw/libraries/trilinos/trilinos-10.12.2/sl6/intel-12.1/mpich2-1.4.1p1/mkl-intel-10.3.12.361'|
| 'VTK_DIR=/opt/apps/ossw/libraries/vtk/vtk-5.10.0/sl6/intel-12.1' |
----------------------------------------------------------------------------------------------------------------------
-----------------------------------------------------------------------------------------------------------
| Matrix Assembly Performance: Alive time=0.027621, Active time=0.026716 |
-----------------------------------------------------------------------------------------------------------
| Event nCalls Total Time Avg Time Total Time Avg Time % of Active Time |
| w/o Sub w/o Sub With Sub With Sub w/o S With S |
|-----------------------------------------------------------------------------------------------------------|
| |
| BCs 112 0.0073 0.000065 0.0073 0.000065 27.20 27.20 |
| Fe 112 0.0008 0.000007 0.0008 0.000007 3.02 3.02 |
| Ke 112 0.0039 0.000035 0.0039 0.000035 14.61 14.61 |
| elem init 112 0.0142 0.000127 0.0142 0.000127 53.26 53.26 |
| matrix insertion 112 0.0005 0.000005 0.0005 0.000005 1.91 1.91 |
-----------------------------------------------------------------------------------------------------------
| Totals: 560 0.0267 100.00 |
-----------------------------------------------------------------------------------------------------------
************************************************************************************************************************
*** WIDEN YOUR WINDOW TO 120 CHARACTERS. Use 'enscript -r -fCourier9' to print this document ***
************************************************************************************************************************
---------------------------------------------- PETSc Performance Summary: ----------------------------------------------
/workspace/libmesh/examples/subdomains/subdomains_ex2/.libs/lt-example-devel on a intel-12. named hbar.ices.utexas.edu with 2 processors, by benkirk Tue Feb 5 13:42:09 2013
Using Petsc Release Version 3.3.0, Patch 2, Fri Jul 13 15:42:00 CDT 2012
Max Max/Min Avg Total
Time (sec): 1.452e-01 1.00119 1.451e-01
Objects: 5.300e+01 1.00000 5.300e+01
Flops: 3.206e+06 1.06355 3.111e+06 6.221e+06
Flops/sec: 2.211e+07 1.06481 2.143e+07 4.287e+07
MPI Messages: 3.050e+01 1.00000 3.050e+01 6.100e+01
MPI Message Lengths: 1.447e+04 1.00000 4.744e+02 2.894e+04
MPI Reductions: 9.100e+01 1.00000
Flop counting convention: 1 flop = 1 real number operation of type (multiply/divide/add/subtract)
e.g., VecAXPY() for real vectors of length N --> 2N flops
and VecAXPY() for complex vectors of length N --> 8N flops
Summary of Stages: ----- Time ------ ----- Flops ----- --- Messages --- -- Message Lengths -- -- Reductions --
Avg %Total Avg %Total counts %Total Avg %Total counts %Total
0: Main Stage: 1.4510e-01 100.0% 6.2212e+06 100.0% 6.100e+01 100.0% 4.744e+02 100.0% 9.000e+01 98.9%
------------------------------------------------------------------------------------------------------------------------
See the 'Profiling' chapter of the users' manual for details on interpreting output.
Phase summary info:
Count: number of times phase was executed
Time and Flops: Max - maximum over all processors
Ratio - ratio of maximum to minimum over all processors
Mess: number of messages sent
Avg. len: average message length
Reduct: number of global reductions
Global: entire computation
Stage: stages of a computation. Set stages with PetscLogStagePush() and PetscLogStagePop().
%T - percent time in this phase %f - percent flops in this phase
%M - percent messages in this phase %L - percent message lengths in this phase
%R - percent reductions in this phase
Total Mflop/s: 10e-6 * (sum of flops over all processors)/(max time over all processors)
------------------------------------------------------------------------------------------------------------------------
Event Count Time (sec) Flops --- Global --- --- Stage --- Total
Max Ratio Max Ratio Max Ratio Mess Avg len Reduct %T %f %M %L %R %T %f %M %L %R Mflop/s
------------------------------------------------------------------------------------------------------------------------
--- Event Stage 0: Main Stage
VecMDot 22 1.0 1.7629e-0310.4 2.73e+05 1.1 0.0e+00 0.0e+00 2.2e+01 1 8 0 0 24 1 8 0 0 24 293
VecNorm 24 1.0 5.2240e-03 4.0 2.59e+04 1.1 0.0e+00 0.0e+00 2.4e+01 2 1 0 0 26 2 1 0 0 27 9
VecScale 23 1.0 2.4796e-05 1.1 1.24e+04 1.1 0.0e+00 0.0e+00 0.0e+00 0 0 0 0 0 0 0 0 0 0 948
VecCopy 2 1.0 3.0994e-06 1.1 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00 0 0 0 0 0 0 0 0 0 0 0
VecSet 30 1.0 1.4782e-05 1.1 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00 0 0 0 0 0 0 0 0 0 0 0
VecAXPY 2 1.0 1.2022e-02 1.6 2.16e+03 1.1 0.0e+00 0.0e+00 0.0e+00 7 0 0 0 0 7 0 0 0 0 0
VecMAXPY 23 1.0 7.8201e-05 1.0 2.97e+05 1.1 0.0e+00 0.0e+00 0.0e+00 0 9 0 0 0 0 9 0 0 0 7188
VecAssemblyBegin 3 1.0 2.0480e-04 1.2 0.00e+00 0.0 2.0e+00 2.9e+02 9.0e+00 0 0 3 2 10 0 0 3 2 10 0
VecAssemblyEnd 3 1.0 1.9073e-05 1.1 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00 0 0 0 0 0 0 0 0 0 0 0
VecScatterBegin 24 1.0 6.3419e-05 1.1 0.00e+00 0.0 4.8e+01 4.2e+02 0.0e+00 0 0 79 70 0 0 0 79 70 0 0
VecScatterEnd 24 1.0 1.9503e-0320.5 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00 1 0 0 0 0 1 0 0 0 0 0
VecNormalize 23 1.0 5.2006e-03 3.9 3.73e+04 1.1 0.0e+00 0.0e+00 2.3e+01 2 1 0 0 25 2 1 0 0 26 14
MatMult 23 1.0 2.0742e-03 5.8 3.46e+05 1.1 4.6e+01 4.1e+02 0.0e+00 1 11 75 65 0 1 11 75 65 0 319
MatSolve 24 1.0 5.1165e-04 1.0 1.37e+06 1.0 0.0e+00 0.0e+00 0.0e+00 0 43 0 0 0 0 43 0 0 0 5249
MatLUFactorNum 1 1.0 7.1287e-04 1.1 8.76e+05 1.0 0.0e+00 0.0e+00 0.0e+00 0 28 0 0 0 0 28 0 0 0 2411
MatILUFactorSym 1 1.0 1.8108e-03 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 3.0e+00 1 0 0 0 3 1 0 0 0 3 0
MatAssemblyBegin 2 1.0 8.7190e-04 7.5 0.00e+00 0.0 3.0e+00 2.3e+03 4.0e+00 0 0 5 24 4 0 0 5 24 4 0
MatAssemblyEnd 2 1.0 2.9802e-04 1.0 0.00e+00 0.0 4.0e+00 1.0e+02 8.0e+00 0 0 7 1 9 0 0 7 1 9 0
MatGetRowIJ 1 1.0 2.8610e-06 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00 0 0 0 0 0 0 0 0 0 0 0
MatGetOrdering 1 1.0 5.2929e-05 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 2.0e+00 0 0 0 0 2 0 0 0 0 2 0
MatZeroEntries 3 1.0 3.7193e-05 1.6 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00 0 0 0 0 0 0 0 0 0 0 0
KSPGMRESOrthog 22 1.0 1.8744e-03 6.8 5.46e+05 1.1 0.0e+00 0.0e+00 2.2e+01 1 17 0 0 24 1 17 0 0 24 552
KSPSetUp 2 1.0 7.4863e-05 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00 0 0 0 0 0 0 0 0 0 0 0
KSPSolve 1 1.0 1.9518e-02 1.0 3.21e+06 1.1 4.6e+01 4.1e+02 5.3e+01 13100 75 65 58 13100 75 65 59 319
PCSetUp 2 1.0 2.9299e-03 1.0 8.76e+05 1.0 0.0e+00 0.0e+00 7.0e+00 2 28 0 0 8 2 28 0 0 8 587
PCSetUpOnBlocks 1 1.0 2.6648e-03 1.0 8.76e+05 1.0 0.0e+00 0.0e+00 5.0e+00 2 28 0 0 5 2 28 0 0 6 645
PCApply 24 1.0 7.3051e-04 1.1 1.37e+06 1.0 0.0e+00 0.0e+00 0.0e+00 0 43 0 0 0 0 43 0 0 0 3676
------------------------------------------------------------------------------------------------------------------------
Memory usage is given in bytes:
Object Type Creations Destructions Memory Descendants' Mem.
Reports information only for process 0.
--- Event Stage 0: Main Stage
Vector 34 34 176976 0
Vector Scatter 2 2 2072 0
Index Set 7 7 7896 0
IS L to G Mapping 1 1 564 0
Matrix 4 4 468600 0
Krylov Solver 2 2 19360 0
Preconditioner 2 2 1784 0
Viewer 1 0 0 0
========================================================================================================================
Average time to get PetscTime(): 9.53674e-08
Average time for MPI_Barrier(): 2.61784e-05
Average time for zero size MPI_Send(): 0.00020647
#PETSc Option Table entries:
-d 2
-ksp_right_pc
-log_summary
-n 15
-pc_type bjacobi
-sub_pc_factor_levels 4
-sub_pc_factor_zeropivot 0
-sub_pc_type ilu
#End of PETSc Option Table entries
Compiled without FORTRAN kernels
Compiled with full precision matrices (default)
sizeof(short) 2 sizeof(int) 4 sizeof(long) 8 sizeof(void*) 8 sizeof(PetscScalar) 8 sizeof(PetscInt) 4
Configure run at: Thu Nov 8 11:21:02 2012
Configure options: --with-debugging=false --COPTFLAGS=-O3 --CXXOPTFLAGS=-O3 --FOPTFLAGS=-O3 --with-clanguage=C++ --with-shared-libraries=1 --with-mpi-dir=/opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1 --with-mumps=true --download-mumps=1 --with-metis=true --download-metis=1 --with-parmetis=true --download-parmetis=1 --with-superlu=true --download-superlu=1 --with-superludir=true --download-superlu_dist=1 --with-blacs=true --download-blacs=1 --with-scalapack=true --download-scalapack=1 --with-hypre=true --download-hypre=1 --with-blas-lib="[/opt/apps/sysnet/intel/12.1/mkl/10.3.12.361/lib/intel64/libmkl_intel_lp64.so,/opt/apps/sysnet/intel/12.1/mkl/10.3.12.361/lib/intel64/libmkl_sequential.so,/opt/apps/sysnet/intel/12.1/mkl/10.3.12.361/lib/intel64/libmkl_core.so]" --with-lapack-lib="[/opt/apps/sysnet/intel/12.1/mkl/10.3.12.361/lib/intel64/libmkl_lapack95_lp64.a]"
-----------------------------------------
Libraries compiled on Thu Nov 8 11:21:02 2012 on daedalus.ices.utexas.edu
Machine characteristics: Linux-2.6.32-279.1.1.el6.x86_64-x86_64-with-redhat-6.3-Carbon
Using PETSc directory: /opt/apps/ossw/libraries/petsc/petsc-3.3-p2
Using PETSc arch: intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt
-----------------------------------------
Using C compiler: /opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1/bin/mpicxx -wd1572 -O3 -fPIC ${COPTFLAGS} ${CFLAGS}
Using Fortran compiler: /opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1/bin/mpif90 -fPIC -O3 ${FOPTFLAGS} ${FFLAGS}
-----------------------------------------
Using include paths: -I/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt/include -I/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/include -I/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/include -I/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt/include -I/opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1/include
-----------------------------------------
Using C linker: /opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1/bin/mpicxx
Using Fortran linker: /opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1/bin/mpif90
Using libraries: -Wl,-rpath,/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt/lib -L/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt/lib -lpetsc -lX11 -Wl,-rpath,/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt/lib -L/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt/lib -lcmumps -ldmumps -lsmumps -lzmumps -lmumps_common -lpord -lHYPRE -lpthread -lsuperlu_dist_3.0 -lparmetis -lmetis -lscalapack -lblacs -lsuperlu_4.3 -Wl,-rpath,/opt/apps/sysnet/intel/12.1/mkl/10.3.12.361/lib/intel64 -L/opt/apps/sysnet/intel/12.1/mkl/10.3.12.361/lib/intel64 -lmkl_lapack95_lp64 -lmkl_intel_lp64 -lmkl_sequential -lmkl_core -Wl,-rpath,/opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1/lib -L/opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1/lib -Wl,-rpath,/opt/apps/sysnet/intel/12.1/composer_xe_2011_sp1.7.256/compiler/lib/intel64 -L/opt/apps/sysnet/intel/12.1/composer_xe_2011_sp1.7.256/compiler/lib/intel64 -Wl,-rpath,/usr/lib/gcc/x86_64-redhat-linux/4.4.6 -L/usr/lib/gcc/x86_64-redhat-linux/4.4.6 -lmpichf90 -lifport -lifcore -lm -lm -lmpichcxx -ldl -lmpich -lopa -lmpl -lrt -lpthread -limf -lsvml -lipgo -ldecimal -lcilkrts -lstdc++ -lgcc_s -lirc -lirc_s -ldl
-----------------------------------------
----------------------------------------------------------------------------------------------------------------
| libMesh Performance: Alive time=0.204726, Active time=0.134462 |
----------------------------------------------------------------------------------------------------------------
| Event nCalls Total Time Avg Time Total Time Avg Time % of Active Time |
| w/o Sub w/o Sub With Sub With Sub w/o S With S |
|----------------------------------------------------------------------------------------------------------------|
| |
| |
| DofMap |
| add_neighbors_to_send_list() 1 0.0099 0.009924 0.0124 0.012447 7.38 9.26 |
| build_sparsity() 1 0.0039 0.003915 0.0135 0.013463 2.91 10.01 |
| create_dof_constraints() 1 0.0003 0.000332 0.0003 0.000332 0.25 0.25 |
| distribute_dofs() 1 0.0088 0.008770 0.0195 0.019543 6.52 14.53 |
| dof_indices() 480 0.0292 0.000061 0.0292 0.000061 21.70 21.70 |
| prepare_send_list() 1 0.0001 0.000053 0.0001 0.000053 0.04 0.04 |
| reinit() 1 0.0105 0.010516 0.0105 0.010516 7.82 7.82 |
| |
| EquationSystems |
| build_solution_vector() 1 0.0008 0.000818 0.0210 0.020992 0.61 15.61 |
| |
| ExodusII_IO |
| write_nodal_data() 1 0.0028 0.002817 0.0028 0.002817 2.10 2.10 |
| |
| FE |
| compute_shape_functions() 178 0.0032 0.000018 0.0032 0.000018 2.38 2.38 |
| init_shape_functions() 67 0.0003 0.000004 0.0003 0.000004 0.22 0.22 |
| inverse_map() 198 0.0016 0.000008 0.0016 0.000008 1.18 1.18 |
| |
| FEMap |
| compute_affine_map() 178 0.0015 0.000009 0.0015 0.000009 1.15 1.15 |
| compute_face_map() 66 0.0011 0.000017 0.0027 0.000041 0.83 2.02 |
| init_face_shape_functions() 1 0.0000 0.000012 0.0000 0.000012 0.01 0.01 |
| init_reference_to_physical_map() 67 0.0014 0.000020 0.0014 0.000020 1.01 1.01 |
| |
| Mesh |
| find_neighbors() 1 0.0030 0.002959 0.0030 0.003037 2.20 2.26 |
| renumber_nodes_and_elem() 2 0.0003 0.000156 0.0003 0.000156 0.23 0.23 |
| |
| MeshCommunication |
| compute_hilbert_indices() 2 0.0021 0.001053 0.0021 0.001053 1.57 1.57 |
| find_global_indices() 2 0.0009 0.000470 0.0041 0.002034 0.70 3.03 |
| parallel_sort() 2 0.0006 0.000320 0.0007 0.000360 0.48 0.54 |
| |
| MeshOutput |
| write_equation_systems() 1 0.0001 0.000097 0.0244 0.024393 0.07 18.14 |
| |
| MeshTools::Generation |
| build_cube() 1 0.0015 0.001516 0.0015 0.001516 1.13 1.13 |
| |
| MetisPartitioner |
| partition() 1 0.0030 0.002998 0.0047 0.004736 2.23 3.52 |
| |
| Parallel |
| allgather() 11 0.0002 0.000018 0.0002 0.000020 0.15 0.16 |
| max(bool) 1 0.0000 0.000003 0.0000 0.000003 0.00 0.00 |
| max(scalar) 115 0.0009 0.000008 0.0009 0.000008 0.69 0.69 |
| max(vector) 26 0.0003 0.000011 0.0008 0.000032 0.21 0.63 |
| min(bool) 133 0.0011 0.000008 0.0011 0.000008 0.79 0.79 |
| min(scalar) 109 0.0130 0.000119 0.0130 0.000119 9.65 9.65 |
| min(vector) 26 0.0003 0.000013 0.0010 0.000038 0.26 0.73 |
| probe() 12 0.0002 0.000015 0.0002 0.000015 0.13 0.13 |
| receive() 12 0.0001 0.000008 0.0003 0.000023 0.07 0.21 |
| send() 12 0.0001 0.000004 0.0001 0.000004 0.04 0.04 |
| send_receive() 16 0.0002 0.000013 0.0006 0.000036 0.16 0.43 |
| sum() 22 0.0004 0.000017 0.0105 0.000479 0.28 7.84 |
| |
| Parallel::Request |
| wait() 12 0.0000 0.000003 0.0000 0.000003 0.03 0.03 |
| |
| Partitioner |
| set_node_processor_ids() 1 0.0007 0.000653 0.0007 0.000731 0.49 0.54 |
| set_parent_processor_ids() 1 0.0002 0.000194 0.0002 0.000194 0.14 0.14 |
| |
| PetscLinearSolver |
| solve() 1 0.0204 0.020432 0.0204 0.020432 15.20 15.20 |
| |
| System |
| assemble() 1 0.0094 0.009434 0.0278 0.027798 7.02 20.67 |
----------------------------------------------------------------------------------------------------------------
| Totals: 1766 0.1345 100.00 |
----------------------------------------------------------------------------------------------------------------
***************************************************************
* Done Running Example subdomains_ex2:
* mpirun -np 2 example-devel -d 2 -n 15 -pc_type bjacobi -sub_pc_type ilu -sub_pc_factor_levels 4 -sub_pc_factor_zeropivot 0 -ksp_right_pc -log_summary
***************************************************************
***************************************************************
* Running Example subdomains_ex2:
* mpirun -np 2 example-devel -d 3 -n 6 -pc_type bjacobi -sub_pc_type ilu -sub_pc_factor_levels 4 -sub_pc_factor_zeropivot 0 -ksp_right_pc -log_summary
***************************************************************
Running /workspace/libmesh/examples/subdomains/subdomains_ex2/.libs/lt-example-devel -d 3 -n 6 -pc_type bjacobi -sub_pc_type ilu -sub_pc_factor_levels 4 -sub_pc_factor_zeropivot 0 -ksp_right_pc -log_summary
Mesh Information:
mesh_dimension()=3
spatial_dimension()=3
n_nodes()=2197
n_local_nodes()=1209
n_elem()=216
n_local_elem()=108
n_active_elem()=216
n_subdomains()=2
n_partitions()=2
n_processors()=2
n_threads()=1
processor_id()=0
EquationSystems
n_systems()=1
System #0, "Poisson"
Type "LinearImplicit"
Variables="u" "v"
Finite Element Types="LAGRANGE", "JACOBI_20_00" "LAGRANGE", "JACOBI_20_00"
Infinite Element Mapping="CARTESIAN" "CARTESIAN"
Approximation Orders="SECOND", "THIRD" "SECOND", "THIRD"
n_dofs()=2522
n_local_dofs()=1404
n_constrained_dofs()=0
n_local_constrained_dofs()=0
n_vectors()=1
n_matrices()=1
DofMap Sparsity
Average On-Processor Bandwidth <= 44.6479
Average Off-Processor Bandwidth <= 4.89611
Maximum On-Processor Bandwidth <= 125
Maximum Off-Processor Bandwidth <= 60
DofMap Constraints
Number of DoF Constraints = 0
Number of Node Constraints = 0
Mesh Information:
mesh_dimension()=3
spatial_dimension()=3
n_nodes()=2197
n_local_nodes()=1209
n_elem()=216
n_local_elem()=108
n_active_elem()=216
n_subdomains()=2
n_partitions()=2
n_processors()=2
n_threads()=1
processor_id()=0
----------------------------------------------------------------------------------------------------------------------
| Processor id: 0 |
| Num Processors: 2 |
| Time: Tue Feb 5 13:42:10 2013 |
| OS: Linux |
| HostName: hbar.ices.utexas.edu |
| OS Release: 2.6.32-279.1.1.el6.x86_64 |
| OS Version: #1 SMP Tue Jul 10 11:24:23 CDT 2012 |
| Machine: x86_64 |
| Username: benkirk |
| Configuration: ./configure '--enable-everything' |
| '--prefix=/workspace/libmesh/install' |
| 'CXX=mpicxx' |
| 'CC=mpicc' |
| 'F77=mpif77' |
| 'FC=mpif90' |
| 'PETSC_DIR=/opt/apps/ossw/libraries/petsc/petsc-3.3-p2' |
| 'PETSC_ARCH=intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt' |
| 'SLEPC_DIR=/opt/apps/ossw/libraries/slepc/slepc-3.3-p2-petsc-3.3-p2-cxx-opt' |
| 'TRILINOS_DIR=/opt/apps/ossw/libraries/trilinos/trilinos-10.12.2/sl6/intel-12.1/mpich2-1.4.1p1/mkl-intel-10.3.12.361'|
| 'VTK_DIR=/opt/apps/ossw/libraries/vtk/vtk-5.10.0/sl6/intel-12.1' |
----------------------------------------------------------------------------------------------------------------------
-----------------------------------------------------------------------------------------------------------
| Matrix Assembly Performance: Alive time=0.288688, Active time=0.287576 |
-----------------------------------------------------------------------------------------------------------
| Event nCalls Total Time Avg Time Total Time Avg Time % of Active Time |
| w/o Sub w/o Sub With Sub With Sub w/o S With S |
|-----------------------------------------------------------------------------------------------------------|
| |
| BCs 108 0.1267 0.001173 0.1267 0.001173 44.06 44.06 |
| Fe 108 0.0035 0.000032 0.0035 0.000032 1.21 1.21 |
| Ke 108 0.0950 0.000879 0.0950 0.000879 33.02 33.02 |
| elem init 108 0.0598 0.000553 0.0598 0.000553 20.78 20.78 |
| matrix insertion 108 0.0027 0.000025 0.0027 0.000025 0.93 0.93 |
-----------------------------------------------------------------------------------------------------------
| Totals: 540 0.2876 100.00 |
-----------------------------------------------------------------------------------------------------------
************************************************************************************************************************
*** WIDEN YOUR WINDOW TO 120 CHARACTERS. Use 'enscript -r -fCourier9' to print this document ***
************************************************************************************************************************
---------------------------------------------- PETSc Performance Summary: ----------------------------------------------
/workspace/libmesh/examples/subdomains/subdomains_ex2/.libs/lt-example-devel on a intel-12. named hbar.ices.utexas.edu with 2 processors, by benkirk Tue Feb 5 13:42:10 2013
Using Petsc Release Version 3.3.0, Patch 2, Fri Jul 13 15:42:00 CDT 2012
Max Max/Min Avg Total
Time (sec): 6.141e-01 1.00001 6.141e-01
Objects: 5.300e+01 1.00000 5.300e+01
Flops: 6.176e+07 2.00499 4.628e+07 9.256e+07
Flops/sec: 1.006e+08 2.00500 7.536e+07 1.507e+08
MPI Messages: 2.650e+01 1.00000 2.650e+01 5.300e+01
MPI Message Lengths: 1.614e+05 1.00000 6.089e+03 3.227e+05
MPI Reductions: 8.300e+01 1.00000
Flop counting convention: 1 flop = 1 real number operation of type (multiply/divide/add/subtract)
e.g., VecAXPY() for real vectors of length N --> 2N flops
and VecAXPY() for complex vectors of length N --> 8N flops
Summary of Stages: ----- Time ------ ----- Flops ----- --- Messages --- -- Message Lengths -- -- Reductions --
Avg %Total Avg %Total counts %Total Avg %Total counts %Total
0: Main Stage: 6.1411e-01 100.0% 9.2558e+07 100.0% 5.300e+01 100.0% 6.089e+03 100.0% 8.200e+01 98.8%
------------------------------------------------------------------------------------------------------------------------
See the 'Profiling' chapter of the users' manual for details on interpreting output.
Phase summary info:
Count: number of times phase was executed
Time and Flops: Max - maximum over all processors
Ratio - ratio of maximum to minimum over all processors
Mess: number of messages sent
Avg. len: average message length
Reduct: number of global reductions
Global: entire computation
Stage: stages of a computation. Set stages with PetscLogStagePush() and PetscLogStagePop().
%T - percent time in this phase %f - percent flops in this phase
%M - percent messages in this phase %L - percent message lengths in this phase
%R - percent reductions in this phase
Total Mflop/s: 10e-6 * (sum of flops over all processors)/(max time over all processors)
------------------------------------------------------------------------------------------------------------------------
Event Count Time (sec) Flops --- Global --- --- Stage --- Total
Max Ratio Max Ratio Max Ratio Mess Avg len Reduct %T %f %M %L %R %T %f %M %L %R Mflop/s
------------------------------------------------------------------------------------------------------------------------
--- Event Stage 0: Main Stage
VecMDot 18 1.0 2.3043e-0313.5 4.80e+05 1.3 0.0e+00 0.0e+00 1.8e+01 0 1 0 0 22 0 1 0 0 22 374
VecNorm 20 1.0 3.0279e-04 4.5 5.62e+04 1.3 0.0e+00 0.0e+00 2.0e+01 0 0 0 0 24 0 0 0 0 24 333
VecScale 19 1.0 2.9802e-05 1.1 2.67e+04 1.3 0.0e+00 0.0e+00 0.0e+00 0 0 0 0 0 0 0 0 0 0 1608
VecCopy 2 1.0 5.0068e-06 1.1 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00 0 0 0 0 0 0 0 0 0 0 0
VecSet 26 1.0 2.5988e-05 1.1 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00 0 0 0 0 0 0 0 0 0 0 0
VecAXPY 2 1.0 7.0388e-03 1.0 5.62e+03 1.3 0.0e+00 0.0e+00 0.0e+00 1 0 0 0 0 1 0 0 0 0 1
VecMAXPY 19 1.0 1.4496e-04 1.2 5.31e+05 1.3 0.0e+00 0.0e+00 0.0e+00 0 1 0 0 0 0 1 0 0 0 6576
VecAssemblyBegin 3 1.0 6.3896e-05 1.1 0.00e+00 0.0 2.0e+00 2.6e+03 9.0e+00 0 0 4 2 11 0 0 4 2 11 0
VecAssemblyEnd 3 1.0 2.0981e-05 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00 0 0 0 0 0 0 0 0 0 0 0
VecScatterBegin 20 1.0 8.0347e-05 1.2 0.00e+00 0.0 4.0e+01 3.1e+03 0.0e+00 0 0 75 38 0 0 0 75 38 0 0
VecScatterEnd 20 1.0 4.7434e-021344.3 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00 4 0 0 0 0 4 0 0 0 0 0
VecNormalize 19 1.0 2.7919e-04 2.0 8.00e+04 1.3 0.0e+00 0.0e+00 1.9e+01 0 0 0 0 23 0 0 0 0 23 515
MatMult 19 1.0 4.8417e-0238.6 2.62e+06 1.3 3.8e+01 3.0e+03 0.0e+00 4 5 72 35 0 4 5 72 35 0 95
MatSolve 20 1.0 4.7274e-03 1.8 1.34e+07 1.7 0.0e+00 0.0e+00 0.0e+00 1 23 0 0 0 1 23 0 0 0 4476
MatLUFactorNum 1 1.0 2.5796e-02 2.2 4.46e+07 2.2 0.0e+00 0.0e+00 0.0e+00 3 70 0 0 0 3 70 0 0 0 2513
MatILUFactorSym 1 1.0 7.0440e-02 1.9 0.00e+00 0.0 0.0e+00 0.0e+00 3.0e+00 9 0 0 0 4 9 0 0 0 4 0
MatAssemblyBegin 2 1.0 4.8208e-04 1.8 0.00e+00 0.0 3.0e+00 6.2e+04 4.0e+00 0 0 6 58 5 0 0 6 58 5 0
MatAssemblyEnd 2 1.0 1.2081e-03 1.1 0.00e+00 0.0 4.0e+00 7.4e+02 8.0e+00 0 0 8 1 10 0 0 8 1 10 0
MatGetRowIJ 1 1.0 9.5367e-07 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00 0 0 0 0 0 0 0 0 0 0 0
MatGetOrdering 1 1.0 5.6028e-05 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 2.0e+00 0 0 0 0 2 0 0 0 0 2 0
MatZeroEntries 3 1.0 2.3699e-04 1.3 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00 0 0 0 0 0 0 0 0 0 0 0
KSPGMRESOrthog 18 1.0 2.4524e-03 7.2 9.60e+05 1.3 0.0e+00 0.0e+00 1.8e+01 0 2 0 0 22 0 2 0 0 22 703
KSPSetUp 2 1.0 8.0824e-05 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00 0 0 0 0 0 0 0 0 0 0 0
KSPSolve 1 1.0 1.1080e-01 1.0 6.18e+07 2.0 3.8e+01 3.0e+03 4.5e+01 18100 72 35 54 18100 72 35 55 835
PCSetUp 2 1.0 9.6663e-02 2.0 4.46e+07 2.2 0.0e+00 0.0e+00 7.0e+00 12 70 0 0 8 12 70 0 0 9 671
PCSetUpOnBlocks 1 1.0 9.6388e-02 2.0 4.46e+07 2.2 0.0e+00 0.0e+00 5.0e+00 12 70 0 0 6 12 70 0 0 6 672
PCApply 20 1.0 4.9284e-03 1.7 1.34e+07 1.7 0.0e+00 0.0e+00 0.0e+00 1 23 0 0 0 1 23 0 0 0 4294
------------------------------------------------------------------------------------------------------------------------
Memory usage is given in bytes:
Object Type Creations Destructions Memory Descendants' Mem.
Reports information only for process 0.
--- Event Stage 0: Main Stage
Vector 34 34 384240 0
Vector Scatter 2 2 2072 0
Index Set 7 7 14760 0
IS L to G Mapping 1 1 564 0
Matrix 4 4 4965804 0
Krylov Solver 2 2 19360 0
Preconditioner 2 2 1784 0
Viewer 1 0 0 0
========================================================================================================================
Average time to get PetscTime(): 9.53674e-08
Average time for MPI_Barrier(): 2.95639e-06
Average time for zero size MPI_Send(): 1.60933e-05
#PETSc Option Table entries:
-d 3
-ksp_right_pc
-log_summary
-n 6
-pc_type bjacobi
-sub_pc_factor_levels 4
-sub_pc_factor_zeropivot 0
-sub_pc_type ilu
#End of PETSc Option Table entries
Compiled without FORTRAN kernels
Compiled with full precision matrices (default)
sizeof(short) 2 sizeof(int) 4 sizeof(long) 8 sizeof(void*) 8 sizeof(PetscScalar) 8 sizeof(PetscInt) 4
Configure run at: Thu Nov 8 11:21:02 2012
Configure options: --with-debugging=false --COPTFLAGS=-O3 --CXXOPTFLAGS=-O3 --FOPTFLAGS=-O3 --with-clanguage=C++ --with-shared-libraries=1 --with-mpi-dir=/opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1 --with-mumps=true --download-mumps=1 --with-metis=true --download-metis=1 --with-parmetis=true --download-parmetis=1 --with-superlu=true --download-superlu=1 --with-superludir=true --download-superlu_dist=1 --with-blacs=true --download-blacs=1 --with-scalapack=true --download-scalapack=1 --with-hypre=true --download-hypre=1 --with-blas-lib="[/opt/apps/sysnet/intel/12.1/mkl/10.3.12.361/lib/intel64/libmkl_intel_lp64.so,/opt/apps/sysnet/intel/12.1/mkl/10.3.12.361/lib/intel64/libmkl_sequential.so,/opt/apps/sysnet/intel/12.1/mkl/10.3.12.361/lib/intel64/libmkl_core.so]" --with-lapack-lib="[/opt/apps/sysnet/intel/12.1/mkl/10.3.12.361/lib/intel64/libmkl_lapack95_lp64.a]"
-----------------------------------------
Libraries compiled on Thu Nov 8 11:21:02 2012 on daedalus.ices.utexas.edu
Machine characteristics: Linux-2.6.32-279.1.1.el6.x86_64-x86_64-with-redhat-6.3-Carbon
Using PETSc directory: /opt/apps/ossw/libraries/petsc/petsc-3.3-p2
Using PETSc arch: intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt
-----------------------------------------
Using C compiler: /opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1/bin/mpicxx -wd1572 -O3 -fPIC ${COPTFLAGS} ${CFLAGS}
Using Fortran compiler: /opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1/bin/mpif90 -fPIC -O3 ${FOPTFLAGS} ${FFLAGS}
-----------------------------------------
Using include paths: -I/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt/include -I/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/include -I/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/include -I/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt/include -I/opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1/include
-----------------------------------------
Using C linker: /opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1/bin/mpicxx
Using Fortran linker: /opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1/bin/mpif90
Using libraries: -Wl,-rpath,/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt/lib -L/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt/lib -lpetsc -lX11 -Wl,-rpath,/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt/lib -L/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt/lib -lcmumps -ldmumps -lsmumps -lzmumps -lmumps_common -lpord -lHYPRE -lpthread -lsuperlu_dist_3.0 -lparmetis -lmetis -lscalapack -lblacs -lsuperlu_4.3 -Wl,-rpath,/opt/apps/sysnet/intel/12.1/mkl/10.3.12.361/lib/intel64 -L/opt/apps/sysnet/intel/12.1/mkl/10.3.12.361/lib/intel64 -lmkl_lapack95_lp64 -lmkl_intel_lp64 -lmkl_sequential -lmkl_core -Wl,-rpath,/opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1/lib -L/opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1/lib -Wl,-rpath,/opt/apps/sysnet/intel/12.1/composer_xe_2011_sp1.7.256/compiler/lib/intel64 -L/opt/apps/sysnet/intel/12.1/composer_xe_2011_sp1.7.256/compiler/lib/intel64 -Wl,-rpath,/usr/lib/gcc/x86_64-redhat-linux/4.4.6 -L/usr/lib/gcc/x86_64-redhat-linux/4.4.6 -lmpichf90 -lifport -lifcore -lm -lm -lmpichcxx -ldl -lmpich -lopa -lmpl -lrt -lpthread -limf -lsvml -lipgo -ldecimal -lcilkrts -lstdc++ -lgcc_s -lirc -lirc_s -ldl
-----------------------------------------
----------------------------------------------------------------------------------------------------------------
| libMesh Performance: Alive time=0.624066, Active time=0.600518 |
----------------------------------------------------------------------------------------------------------------
| Event nCalls Total Time Avg Time Total Time Avg Time % of Active Time |
| w/o Sub w/o Sub With Sub With Sub w/o S With S |
|----------------------------------------------------------------------------------------------------------------|
| |
| |
| DofMap |
| add_neighbors_to_send_list() 1 0.0299 0.029933 0.0504 0.050410 4.98 8.39 |
| build_sparsity() 1 0.0250 0.024975 0.0538 0.053844 4.16 8.97 |
| create_dof_constraints() 1 0.0004 0.000354 0.0004 0.000354 0.06 0.06 |
| distribute_dofs() 1 0.0202 0.020246 0.0437 0.043726 3.37 7.28 |
| dof_indices() 528 0.0954 0.000181 0.0954 0.000181 15.88 15.88 |
| prepare_send_list() 1 0.0005 0.000452 0.0005 0.000452 0.08 0.08 |
| reinit() 1 0.0228 0.022791 0.0228 0.022791 3.80 3.80 |
| |
| EquationSystems |
| build_solution_vector() 1 0.0011 0.001122 0.0263 0.026269 0.19 4.37 |
| |
| ExodusII_IO |
| write_nodal_data() 1 0.0041 0.004140 0.0041 0.004140 0.69 0.69 |
| |
| FE |
| compute_shape_functions() 288 0.0339 0.000118 0.0339 0.000118 5.64 5.64 |
| init_shape_functions() 181 0.0011 0.000006 0.0011 0.000006 0.19 0.19 |
| |
| FEMap |
| compute_affine_map() 288 0.0095 0.000033 0.0095 0.000033 1.59 1.59 |
| compute_face_map() 180 0.0043 0.000024 0.0043 0.000024 0.71 0.71 |
| init_face_shape_functions() 1 0.0001 0.000081 0.0001 0.000081 0.01 0.01 |
| init_reference_to_physical_map() 181 0.0652 0.000360 0.0652 0.000360 10.85 10.85 |
| |
| Mesh |
| find_neighbors() 1 0.0041 0.004093 0.0042 0.004202 0.68 0.70 |
| renumber_nodes_and_elem() 2 0.0007 0.000371 0.0007 0.000371 0.12 0.12 |
| |
| MeshCommunication |
| compute_hilbert_indices() 2 0.0020 0.001021 0.0020 0.001021 0.34 0.34 |
| find_global_indices() 2 0.0009 0.000451 0.0040 0.001997 0.15 0.66 |
| parallel_sort() 2 0.0007 0.000325 0.0007 0.000369 0.11 0.12 |
| |
| MeshOutput |
| write_equation_systems() 1 0.0001 0.000113 0.0306 0.030578 0.02 5.09 |
| |
| MeshTools::Generation |
| build_cube() 1 0.0033 0.003323 0.0033 0.003323 0.55 0.55 |
| |
| MetisPartitioner |
| partition() 1 0.0045 0.004545 0.0062 0.006243 0.76 1.04 |
| |
| Parallel |
| allgather() 11 0.0003 0.000023 0.0003 0.000025 0.04 0.05 |
| max(bool) 1 0.0000 0.000003 0.0000 0.000003 0.00 0.00 |
| max(scalar) 115 0.0003 0.000003 0.0003 0.000003 0.05 0.05 |
| max(vector) 26 0.0002 0.000006 0.0004 0.000014 0.03 0.06 |
| min(bool) 133 0.0003 0.000002 0.0003 0.000002 0.05 0.05 |
| min(scalar) 109 0.0040 0.000036 0.0040 0.000036 0.66 0.66 |
| min(vector) 26 0.0002 0.000008 0.0004 0.000017 0.03 0.07 |
| probe() 12 0.0006 0.000053 0.0006 0.000053 0.11 0.11 |
| receive() 12 0.0001 0.000011 0.0008 0.000064 0.02 0.13 |
| send() 12 0.0001 0.000005 0.0001 0.000005 0.01 0.01 |
| send_receive() 16 0.0004 0.000022 0.0012 0.000077 0.06 0.20 |
| sum() 22 0.0002 0.000010 0.0004 0.000019 0.04 0.07 |
| |
| Parallel::Request |
| wait() 12 0.0000 0.000003 0.0000 0.000003 0.01 0.01 |
| |
| Partitioner |
| set_node_processor_ids() 1 0.0011 0.001065 0.0012 0.001151 0.18 0.19 |
| set_parent_processor_ids() 1 0.0002 0.000199 0.0002 0.000199 0.03 0.03 |
| |
| PetscLinearSolver |
| solve() 1 0.1133 0.113285 0.1133 0.113285 18.86 18.86 |
| |
| System |
| assemble() 1 0.1495 0.149480 0.2889 0.288905 24.89 48.11 |
----------------------------------------------------------------------------------------------------------------
| Totals: 2178 0.6005 100.00 |
----------------------------------------------------------------------------------------------------------------
***************************************************************
* Done Running Example subdomains_ex2:
* mpirun -np 2 example-devel -d 3 -n 6 -pc_type bjacobi -sub_pc_type ilu -sub_pc_factor_levels 4 -sub_pc_factor_zeropivot 0 -ksp_right_pc -log_summary
***************************************************************
This library is free software; you can redistribute it and/or modify it under the terms of the GNU Lesser General Public License as published by the Free Software Foundation; either version 2.1 of the License, or (at your option) any later version.
This library is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more details.
You should have received a copy of the GNU Lesser General Public License along with this library; if not, write to the Free Software Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
C++ Includes