The source file eigenproblems_ex1.C with comments:
#include "libmesh/libmesh.h"
#include "libmesh/mesh.h"
#include "libmesh/mesh_generation.h"
#include "libmesh/exodusII_io.h"
#include "libmesh/eigen_system.h"
#include "libmesh/equation_systems.h"
#include "libmesh/fe.h"
#include "libmesh/quadrature_gauss.h"
#include "libmesh/dense_matrix.h"
#include "libmesh/sparse_matrix.h"
#include "libmesh/numeric_vector.h"
#include "libmesh/dof_map.h"
Bring in everything from the libMesh namespace
using namespace libMesh;
Function prototype. This is the function that will assemble
the eigen system. Here, we will simply assemble a mass matrix.
void assemble_mass(EquationSystems& es,
const std::string& system_name);
int main (int argc, char** argv)
{
Initialize libMesh and the dependent libraries.
LibMeshInit init (argc, argv);
Skip SLEPc examples on a non-SLEPc libMesh build
#ifndef LIBMESH_HAVE_SLEPC
libmesh_example_assert(false, "--enable-slepc");
}
#else
#ifdef LIBMESH_DEFAULT_SINGLE_PRECISION
SLEPc currently gives us a nasty crash with Real==float
libmesh_example_assert(false, "--disable-singleprecision");
#endif
Check for proper usage.
if (argc < 3)
{
if (libMesh::processor_id() == 0)
std::cerr << "\nUsage: " << argv[0]
<< " -n <number of eigen values>"
<< std::endl;
libmesh_error();
}
Tell the user what we are doing.
else
{
std::cout << "Running " << argv[0];
for (int i=1; i<argc; i++)
std::cout << " " << argv[i];
std::cout << std::endl << std::endl;
}
Get the number of eigen values to be computed from argv[2]
const unsigned int nev = std::atoi(argv[2]);
Skip this 2D example if libMesh was compiled as 1D-only.
libmesh_example_assert(2 <= LIBMESH_DIM, "2D support");
Create a mesh.
Mesh mesh;
Use the internal mesh generator to create a uniform
2D grid on a square.
MeshTools::Generation::build_square (mesh,
20, 20,
-1., 1.,
-1., 1.,
QUAD4);
Print information about the mesh to the screen.
mesh.print_info();
Create an equation systems object.
EquationSystems equation_systems (mesh);
Create a EigenSystem named "Eigensystem" and (for convenience)
use a reference to the system we create.
EigenSystem & eigen_system =
equation_systems.add_system<EigenSystem> ("Eigensystem");
Declare the system variables.
Adds the variable "p" to "Eigensystem". "p"
will be approximated using second-order approximation.
eigen_system.add_variable("p", FIRST);
Give the system a pointer to the matrix assembly
function defined below.
eigen_system.attach_assemble_function (assemble_mass);
Set necessary parametrs used in EigenSystem::solve(),
i.e. the number of requested eigenpairs \p nev and the number
of basis vectors \p ncv used in the solution algorithm. Note that
ncv >= nev must hold and ncv >= 2*nev is recommended.
equation_systems.parameters.set<unsigned int>("eigenpairs") = nev;
equation_systems.parameters.set<unsigned int>("basis vectors") = nev*3;
You may optionally change the default eigensolver used by SLEPc.
The Krylov-Schur method is mathematically equivalent to implicitly
restarted Arnoldi, the method of Arpack, so there is currently no
point in using SLEPc with Arpack.
ARNOLDI = default in SLEPc 2.3.1 and earlier
KRYLOVSCHUR default in SLEPc 2.3.2 and later
eigen_system.eigen_solver->set_eigensolver_type(KRYLOVSCHUR);
Set the solver tolerance and the maximum number of iterations.
Set the solver tolerance and the maximum number of iterations.
equation_systems.parameters.set<Real>
("linear solver tolerance") = pow(TOLERANCE, 5./3.);
equation_systems.parameters.set<unsigned int>
("linear solver maximum iterations") = 1000;
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 "Eigensystem".
eigen_system.solve();
Get the number of converged eigen pairs.
unsigned int nconv = eigen_system.get_n_converged();
std::cout << "Number of converged eigenpairs: " << nconv
<< "\n" << std::endl;
Get the last converged eigenpair
if (nconv != 0)
{
eigen_system.get_eigenpair(nconv-1);
#ifdef LIBMESH_HAVE_EXODUS_API
Write the eigen vector to file.
ExodusII_IO (mesh).write_equation_systems ("out.e", equation_systems);
#endif // #ifdef LIBMESH_HAVE_EXODUS_API
}
else
{
std::cout << "WARNING: Solver did not converge!\n" << nconv << std::endl;
}
All done.
return 0;
}
#endif // LIBMESH_HAVE_SLEPC
void assemble_mass(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, "Eigensystem");
#ifdef LIBMESH_HAVE_SLEPC
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 our system.
EigenSystem & eigen_system = es.get_system<EigenSystem> (system_name);
Get a constant reference to the Finite Element type
for the first (and only) variable in the system.
FEType fe_type = eigen_system.get_dof_map().variable_type(0);
A reference to the system matrix
SparseMatrix<Number>& matrix_A = *eigen_system.matrix_A;
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 Gauss quadrature rule for numerical integration.
Use the default quadrature order.
QGauss qrule (dim, fe_type.default_quadrature_order());
Tell the finite element object to use our quadrature rule.
fe->attach_quadrature_rule (&qrule);
The element Jacobian * quadrature weight at each integration point.
const std::vector<Real>& JxW = fe->get_JxW();
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 >& dphi = fe->get_dphi();
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.
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.
const DofMap& dof_map = eigen_system.get_dof_map();
The element mass matrix.
DenseMatrix<Number> Me;
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;
Now we will loop over all the elements in the mesh that
live on the local processor. We will compute the element
matrix and right-hand-side contribution. In case users
later modify this program to include refinement, we will
be safe and will only consider the active elements;
hence we use a variant of the \p active_elem_iterator.
MeshBase::const_element_iterator el = mesh.active_local_elements_begin();
const MeshBase::const_element_iterator end_el = mesh.active_local_elements_end();
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 matrices and rhs 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).
Me.resize (dof_indices.size(), dof_indices.size());
Now loop over the quadrature points. This handles
the numeric integration.
We will build the element matrix. This involves a double loop to integrate the test funcions (i) against the trial functions (j).
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 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++)
Me(i,j) += JxW[qp]*phi[i][qp]*phi[j][qp];
On an unrefined mesh, constrain_element_matrix does
nothing. If this assembly function is ever repurposed to
run on a refined mesh, getting the hanging node constraints
right will be important. Note that, even with
asymmetric_constraint_rows = false, the constrained dof
diagonals still exist in the matrix, with diagonal entries
that are there to ensure non-singular matrices for linear
solves but which would generate positive non-physical
eigenvalues for eigensolves.
dof_map.constrain_element_matrix(Me, dof_indices, false);
Finally, simply add the element contribution to the
overall matrix.
matrix_A.add_matrix (Me, dof_indices);
} // end of element loop
#endif // LIBMESH_HAVE_SLEPC
/**
* All done!
*/
return;
}
The source file eigenproblems_ex1.C without comments:
#include "libmesh/libmesh.h"
#include "libmesh/mesh.h"
#include "libmesh/mesh_generation.h"
#include "libmesh/exodusII_io.h"
#include "libmesh/eigen_system.h"
#include "libmesh/equation_systems.h"
#include "libmesh/fe.h"
#include "libmesh/quadrature_gauss.h"
#include "libmesh/dense_matrix.h"
#include "libmesh/sparse_matrix.h"
#include "libmesh/numeric_vector.h"
#include "libmesh/dof_map.h"
using namespace libMesh;
void assemble_mass(EquationSystems& es,
const std::string& system_name);
int main (int argc, char** argv)
{
LibMeshInit init (argc, argv);
#ifndef LIBMESH_HAVE_SLEPC
libmesh_example_assert(false, "--enable-slepc");
}
#else
#ifdef LIBMESH_DEFAULT_SINGLE_PRECISION
libmesh_example_assert(false, "--disable-singleprecision");
#endif
if (argc < 3)
{
if (libMesh::processor_id() == 0)
std::cerr << "\nUsage: " << argv[0]
<< " -n <number of eigen values>"
<< 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;
}
const unsigned int nev = std::atoi(argv[2]);
libmesh_example_assert(2 <= LIBMESH_DIM, "2D support");
Mesh mesh;
MeshTools::Generation::build_square (mesh,
20, 20,
-1., 1.,
-1., 1.,
QUAD4);
mesh.print_info();
EquationSystems equation_systems (mesh);
EigenSystem & eigen_system =
equation_systems.add_system<EigenSystem> ("Eigensystem");
eigen_system.add_variable("p", FIRST);
eigen_system.attach_assemble_function (assemble_mass);
equation_systems.parameters.set<unsigned int>("eigenpairs") = nev;
equation_systems.parameters.set<unsigned int>("basis vectors") = nev*3;
equation_systems.parameters.set<Real>
("linear solver tolerance") = pow(TOLERANCE, 5./3.);
equation_systems.parameters.set<unsigned int>
("linear solver maximum iterations") = 1000;
equation_systems.init();
equation_systems.print_info();
eigen_system.solve();
unsigned int nconv = eigen_system.get_n_converged();
std::cout << "Number of converged eigenpairs: " << nconv
<< "\n" << std::endl;
if (nconv != 0)
{
eigen_system.get_eigenpair(nconv-1);
#ifdef LIBMESH_HAVE_EXODUS_API
ExodusII_IO (mesh).write_equation_systems ("out.e", equation_systems);
#endif // #ifdef LIBMESH_HAVE_EXODUS_API
}
else
{
std::cout << "WARNING: Solver did not converge!\n" << nconv << std::endl;
}
return 0;
}
#endif // LIBMESH_HAVE_SLEPC
void assemble_mass(EquationSystems& es,
const std::string& system_name)
{
libmesh_assert_equal_to (system_name, "Eigensystem");
#ifdef LIBMESH_HAVE_SLEPC
const MeshBase& mesh = es.get_mesh();
const unsigned int dim = mesh.mesh_dimension();
EigenSystem & eigen_system = es.get_system<EigenSystem> (system_name);
FEType fe_type = eigen_system.get_dof_map().variable_type(0);
SparseMatrix<Number>& matrix_A = *eigen_system.matrix_A;
AutoPtr<FEBase> fe (FEBase::build(dim, fe_type));
QGauss qrule (dim, fe_type.default_quadrature_order());
fe->attach_quadrature_rule (&qrule);
const std::vector<Real>& JxW = fe->get_JxW();
const std::vector<std::vector<Real> >& phi = fe->get_phi();
const DofMap& dof_map = eigen_system.get_dof_map();
DenseMatrix<Number> Me;
std::vector<dof_id_type> dof_indices;
MeshBase::const_element_iterator el = mesh.active_local_elements_begin();
const MeshBase::const_element_iterator end_el = mesh.active_local_elements_end();
for ( ; el != end_el; ++el)
{
const Elem* elem = *el;
dof_map.dof_indices (elem, dof_indices);
fe->reinit (elem);
Me.resize (dof_indices.size(), 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++)
Me(i,j) += JxW[qp]*phi[i][qp]*phi[j][qp];
dof_map.constrain_element_matrix(Me, dof_indices, false);
matrix_A.add_matrix (Me, dof_indices);
} // end of element loop
#endif // LIBMESH_HAVE_SLEPC
/**
* All done!
*/
return;
}
The console output of the program:
***************************************************************
* Running Example eigenproblems_ex1:
* mpirun -np 2 example-devel -n 5 -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/eigenproblems/eigenproblems_ex1/.libs/lt-example-devel -n 5 -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()=441
n_local_nodes()=232
n_elem()=400
n_local_elem()=200
n_active_elem()=400
n_subdomains()=1
n_partitions()=2
n_processors()=2
n_threads()=1
processor_id()=0
EquationSystems
n_systems()=1
System #0, "Eigensystem"
Type "Eigen"
Variables="p"
Finite Element Types="LAGRANGE", "JACOBI_20_00"
Infinite Element Mapping="CARTESIAN"
Approximation Orders="FIRST", "THIRD"
n_dofs()=441
n_local_dofs()=232
n_constrained_dofs()=0
n_local_constrained_dofs()=0
n_vectors()=0
n_matrices()=1
DofMap Sparsity
Average On-Processor Bandwidth <= 8.30612
Average Off-Processor Bandwidth <= 0.290249
Maximum On-Processor Bandwidth <= 11
Maximum Off-Processor Bandwidth <= 4
DofMap Constraints
Number of DoF Constraints = 0
Number of Node Constraints = 0
Number of converged eigenpairs: 5
************************************************************************************************************************
*** WIDEN YOUR WINDOW TO 120 CHARACTERS. Use 'enscript -r -fCourier9' to print this document ***
************************************************************************************************************************
---------------------------------------------- PETSc Performance Summary: ----------------------------------------------
/workspace/libmesh/examples/eigenproblems/eigenproblems_ex1/.libs/lt-example-devel on a intel-12. named hbar.ices.utexas.edu with 2 processors, by benkirk Tue Feb 5 13:40:44 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.108e-01 1.00000 1.108e-01
Objects: 4.400e+01 1.00000 4.400e+01
Flops: 1.335e+07 1.11152 1.268e+07 2.536e+07
Flops/sec: 1.204e+08 1.11152 1.144e+08 2.288e+08
MPI Messages: 6.055e+02 1.00000 6.055e+02 1.211e+03
MPI Message Lengths: 1.169e+05 1.00000 1.931e+02 2.338e+05
MPI Reductions: 1.276e+03 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.1080e-01 100.0% 2.5356e+07 100.0% 1.211e+03 100.0% 1.931e+02 100.0% 1.275e+03 99.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
STSetUp 1 1.0 1.0967e-05 1.2 0.00e+00 0.0 0.0e+00 0.0e+00 2.0e+00 0 0 0 0 0 0 0 0 0 0 0
STApply 600 1.0 4.6027e-03 1.0 2.22e+06 1.1 1.2e+03 1.9e+02 0.0e+00 4 17 99 99 0 4 17 99 99 0 913
EPSSetUp 1 1.0 4.4394e-04 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 3.6e+01 0 0 0 0 3 0 0 0 0 3 0
EPSSolve 1 1.0 2.8505e-02 1.0 1.33e+07 1.1 1.2e+03 1.9e+02 1.2e+03 26100 99 99 95 26100 99 99 95 889
IPOrthogonalize 601 1.0 1.9284e-02 1.0 1.06e+07 1.1 0.0e+00 0.0e+00 1.2e+03 17 80 0 0 95 17 80 0 0 95 1048
IPInnerProduct 4822 1.0 1.6248e-02 1.0 5.59e+06 1.1 0.0e+00 0.0e+00 1.2e+03 15 42 0 0 95 15 42 0 0 95 654
DSSolve 47 1.0 2.3146e-03 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00 2 0 0 0 0 2 0 0 0 0 0
DSVectors 53 1.0 3.1209e-04 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
DSOther 94 1.0 6.2871e-04 1.0 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
UpdateVectors 48 1.0 5.0783e-04 1.1 3.74e+05 1.1 0.0e+00 0.0e+00 0.0e+00 0 3 0 0 0 0 3 0 0 0 1398
VecMAXPBY 1200 1.0 1.8556e-03 1.0 5.04e+06 1.1 0.0e+00 0.0e+00 0.0e+00 2 38 0 0 0 2 38 0 0 0 5162
VecScale 554 1.0 2.4223e-04 1.1 1.29e+05 1.1 0.0e+00 0.0e+00 0.0e+00 0 1 0 0 0 0 1 0 0 0 1009
VecCopy 2 1.0 4.0531e-06 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
VecSet 3 1.0 2.1458e-06 2.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
VecAssemblyBegin 1 1.0 6.8903e-05 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 3.0e+00 0 0 0 0 0 0 0 0 0 0 0
VecAssemblyEnd 1 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 600 1.0 8.7953e-04 1.0 0.00e+00 0.0 1.2e+03 1.9e+02 0.0e+00 1 0 99 99 0 1 0 99 99 0 0
VecScatterEnd 600 1.0 6.7735e-04 1.4 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
VecReduceArith 2411 1.0 1.1530e-02 1.0 5.59e+06 1.1 0.0e+00 0.0e+00 0.0e+00 10 42 0 0 0 10 42 0 0 0 921
VecReduceComm 1211 1.0 3.3510e-03 1.1 0.00e+00 0.0 0.0e+00 0.0e+00 1.2e+03 3 0 0 0 95 3 0 0 0 95 0
MatMult 600 1.0 4.3874e-03 1.0 2.22e+06 1.1 1.2e+03 1.9e+02 0.0e+00 4 17 99 99 0 4 17 99 99 0 957
MatAssemblyBegin 1 1.0 5.5790e-05 1.0 0.00e+00 0.0 3.0e+00 9.8e+02 2.0e+00 0 0 0 1 0 0 0 0 1 0 0
MatAssemblyEnd 1 1.0 2.6989e-04 1.0 0.00e+00 0.0 4.0e+00 5.0e+01 8.0e+00 0 0 0 0 1 0 0 0 0 1 0
MatZeroEntries 2 1.0 1.0014e-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
------------------------------------------------------------------------------------------------------------------------
Memory usage is given in bytes:
Object Type Creations Destructions Memory Descendants' Mem.
Reports information only for process 0.
--- Event Stage 0: Main Stage
Container 2 2 1096 0
Spectral Transform 1 1 736 0
Eigenproblem Solver 1 1 1680 0
Inner product 1 1 624 0
Direct solver 1 1 9408 0
Vector 24 24 45440 0
Vector Scatter 2 2 2072 0
Index Set 4 4 3152 0
IS L to G Mapping 1 1 564 0
Matrix 3 3 39476 0
PetscRandom 1 1 608 0
Krylov Solver 1 1 1072 0
Preconditioner 1 1 856 0
Viewer 1 0 0 0
========================================================================================================================
Average time to get PetscTime(): 0
Average time for MPI_Barrier(): 7.62939e-07
Average time for zero size MPI_Send(): 1.50204e-05
#PETSc Option Table entries:
-ksp_right_pc
-log_summary
-n 5
-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
-----------------------------------------
----------------------------------------------------------------------------------------------------------------------
| Processor id: 0 |
| Num Processors: 2 |
| Time: Tue Feb 5 13:40:44 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' |
----------------------------------------------------------------------------------------------------------------------
----------------------------------------------------------------------------------------------------------------
| libMesh Performance: Alive time=0.122503, Active time=0.10149 |
----------------------------------------------------------------------------------------------------------------
| 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.0046 0.004633 0.0059 0.005875 4.56 5.79 |
| build_sparsity() 1 0.0031 0.003144 0.0091 0.009116 3.10 8.98 |
| create_dof_constraints() 1 0.0004 0.000429 0.0004 0.000429 0.42 0.42 |
| distribute_dofs() 1 0.0044 0.004437 0.0117 0.011673 4.37 11.50 |
| dof_indices() 646 0.0214 0.000033 0.0214 0.000033 21.05 21.05 |
| prepare_send_list() 1 0.0000 0.000030 0.0000 0.000030 0.03 0.03 |
| reinit() 1 0.0070 0.007003 0.0070 0.007003 6.90 6.90 |
| |
| EquationSystems |
| build_solution_vector() 1 0.0014 0.001404 0.0109 0.010892 1.38 10.73 |
| |
| ExodusII_IO |
| write_nodal_data() 1 0.0047 0.004676 0.0047 0.004676 4.61 4.61 |
| |
| FE |
| compute_shape_functions() 200 0.0008 0.000004 0.0008 0.000004 0.77 0.77 |
| init_shape_functions() 1 0.0000 0.000037 0.0000 0.000037 0.04 0.04 |
| |
| FEMap |
| compute_affine_map() 200 0.0010 0.000005 0.0010 0.000005 0.94 0.94 |
| init_reference_to_physical_map() 1 0.0001 0.000069 0.0001 0.000069 0.07 0.07 |
| |
| Mesh |
| find_neighbors() 1 0.0051 0.005143 0.0052 0.005201 5.07 5.12 |
| renumber_nodes_and_elem() 2 0.0003 0.000125 0.0003 0.000125 0.25 0.25 |
| |
| MeshCommunication |
| compute_hilbert_indices() 2 0.0037 0.001831 0.0037 0.001831 3.61 3.61 |
| find_global_indices() 2 0.0015 0.000771 0.0063 0.003146 1.52 6.20 |
| parallel_sort() 2 0.0007 0.000363 0.0008 0.000404 0.72 0.80 |
| |
| MeshOutput |
| write_equation_systems() 1 0.0001 0.000094 0.0157 0.015730 0.09 15.50 |
| |
| MeshTools::Generation |
| build_cube() 1 0.0012 0.001163 0.0012 0.001163 1.15 1.15 |
| |
| MetisPartitioner |
| partition() 1 0.0046 0.004621 0.0074 0.007428 4.55 7.32 |
| |
| Parallel |
| allgather() 9 0.0002 0.000020 0.0002 0.000022 0.18 0.20 |
| max(bool) 1 0.0000 0.000003 0.0000 0.000003 0.00 0.00 |
| max(scalar) 105 0.0003 0.000003 0.0003 0.000003 0.31 0.31 |
| max(vector) 24 0.0002 0.000007 0.0004 0.000015 0.16 0.37 |
| min(bool) 121 0.0003 0.000003 0.0003 0.000003 0.31 0.31 |
| min(scalar) 99 0.0008 0.000009 0.0008 0.000009 0.84 0.84 |
| min(vector) 24 0.0002 0.000009 0.0004 0.000018 0.20 0.42 |
| probe() 12 0.0002 0.000015 0.0002 0.000015 0.18 0.18 |
| receive() 12 0.0001 0.000007 0.0003 0.000022 0.08 0.27 |
| send() 12 0.0000 0.000004 0.0000 0.000004 0.04 0.04 |
| send_receive() 16 0.0002 0.000013 0.0006 0.000035 0.20 0.55 |
| sum() 20 0.0001 0.000007 0.0002 0.000011 0.14 0.22 |
| |
| Parallel::Request |
| wait() 12 0.0000 0.000003 0.0000 0.000003 0.04 0.04 |
| |
| Partitioner |
| set_node_processor_ids() 1 0.0006 0.000608 0.0007 0.000691 0.60 0.68 |
| set_parent_processor_ids() 1 0.0003 0.000296 0.0003 0.000296 0.29 0.29 |
| |
| SlepcEigenSolver |
| solve_standard() 1 0.0299 0.029904 0.0299 0.029904 29.46 29.46 |
| |
| System |
| assemble() 1 0.0018 0.001790 0.0092 0.009173 1.76 9.04 |
----------------------------------------------------------------------------------------------------------------
| Totals: 1539 0.1015 100.00 |
----------------------------------------------------------------------------------------------------------------
***************************************************************
* Done Running Example eigenproblems_ex1:
* mpirun -np 2 example-devel -n 5 -pc_type bjacobi -sub_pc_type ilu -sub_pc_factor_levels 4 -sub_pc_factor_zeropivot 0 -ksp_right_pc -log_summary
***************************************************************