]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
Patch by Juan Carlos Araujo Cabarcas: Explain how to use ARPACK as an alternative...
authorbangerth <bangerth@0785d39b-7218-0410-832d-ea1e28bc413d>
Mon, 25 Mar 2013 14:53:34 +0000 (14:53 +0000)
committerbangerth <bangerth@0785d39b-7218-0410-832d-ea1e28bc413d>
Mon, 25 Mar 2013 14:53:34 +0000 (14:53 +0000)
git-svn-id: https://svn.dealii.org/trunk@29027 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/doc/news/changes.h
deal.II/examples/step-36/doc/results.dox

index f1a687f08905189a1a61ba93f26b625d90d546f5..c81a48e3400b93930b72d3a05886a625ebba52c8 100644 (file)
@@ -42,6 +42,11 @@ this function.
 
 
 <ol>
+  <li> New: The results section of step-36 now explains how to use ARPACK
+  as an alternative to SLEPc as eigenvalue solver.
+  <br>
+  (Juan Carlos Araujo Cabarcas, 2013/03/25)
+
   <li> New: deal.II now uses <a href="http://www.cmake.org/">CMake</a>
   as its configuration and build tool. Please read through the
   readme and other installation files for information about how the
@@ -58,7 +63,6 @@ this function.
   </ul>
   <br>
   (Matthias Maier, 2013/03/07)
-
 </ol>
 
 
@@ -78,10 +82,10 @@ sense.
 (Guido Kanschat, 2013/03/21)
 </li>
 
-<li> Added GridOut::write_svg to allow for the output of two-dimensional 
-triangulations in two space dimensions in the SVG format (Scalable Vector 
-Graphics, an XML-based vector image format recommended by the World 
-Wide Web Consortium W3C). This function also provides cell coloring 
+<li> Added GridOut::write_svg to allow for the output of two-dimensional
+triangulations in two space dimensions in the SVG format (Scalable Vector
+Graphics, an XML-based vector image format recommended by the World
+Wide Web Consortium W3C). This function also provides cell coloring
 and cell labeling for the visualization of basic cell properties.
 <br>
 (Christian Wülker, 2013/03/21)
index fbbb313f63ad4ad89ab05273fa46c7c85da08273..9684b72a4e973a60dc410e68e716b19c96cbb69a 100644 (file)
@@ -190,11 +190,94 @@ PETScWrappers and SLEPcWrappers and is suitable for running on serial
 machine architecture. However, for larger grids and with a larger
 number of degrees-of-freedom, we may want to run our application on
 parallel architectures. A parallel implementation of the above code
-can be particuarily useful here since the generalized eigenspectrum
+can be particularily useful here since the generalized eigenspectrum
 problem is somewhat more expensive to solve than the standard problems
-considered most of the earlier tutorials. Fortunately, modifying the above
-program to be MPI complient is a relatively straightforward
+considered in most of the earlier tutorials. Fortunately, modifying the above
+program to be MPI compliant is a relatively straightforward
 procedure. A sketch of how this can be done can be found in @ref
 step_17 "step-17".
 
+<li> Finally, there are alternatives to using the SLEPc eigenvalue
+solvers. deal.II has interfaces to one of them, ARPACK (see
+http://www.dealii.org/developer/external-libs/arpack.html), implemented in the
+ArpackSolver class. Here is a short and quick overview of what one would need
+to change to use it, provided you have a working installation of ARPACK and
+deal.II has been configured properly for it (see the deal.II ReadMe file
+at http://www.dealii.org/readme.html):
+
+First, in order to use the ARPACK interfaces, we can go back to using standard
+deal.II matrices and vectors, so we start by replacing the PETSc and SLEPc
+headers
+@code
+#include <deal.II/lac/petsc_sparse_matrix.h>
+#include <deal.II/lac/petsc_vector.h>
+#include <deal.II/lac/slepc_solver.h>
+@endcode
+with these:
+@code
+#include <deal.II/lac/arpack_solver.h>
+#include <deal.II/lac/sparse_direct.h>
+#include <deal.II/lac/sparse_matrix.h>
+#include <deal.II/lac/compressed_sparsity_pattern.h>
+@endcode
+ARPACK allows complex eigenvalues, so we will also need
+@code
+#include <complex>
+@endcode
+
+Secondly, we switch back to the deal.II matrix and vector definitions in the
+main class:
+@code
+    SparsityPattern                     sparsity_pattern;
+    SparseMatrix<double>                stiffness_matrix, mass_matrix;
+    std::vector<Vector<double> >        eigenfunctions;
+    std::vector<complex<double>>      eigenvalues;
+@endcode
+and initialize them as usual in <code>make_grid_and_dofs()</code>:
+@code
+    sparsity_pattern.reinit (dof_handler.n_dofs(),
+                             dof_handler.n_dofs(),
+                             dof_handler.max_couplings_between_dofs());
+
+    DoFTools::make_sparsity_pattern (dof_handler, sparsity_pattern);
+    constraints.condense (sparsity_pattern);
+    sparsity_pattern.compress();
+
+    stiffness_matrix.reinit (sparsity_pattern);
+    mass_matrix.reinit (sparsity_pattern);
+@endcode
+
+For solving the eigenvalue problem with ARPACK, we finally need to modify
+<code>solve()</code>:
+@code
+  template <int dim>
+  unsigned int EigenvalueProblem<dim>::solve ()
+  {
+    const unsigned int num_eigs =
+      parameters.get_integer ("Number of eigenvalues/eigenfunctions");
+
+    SolverControl solver_control (dof_handler.n_dofs(), 1e-9);
+
+    SparseDirectUMFPACK inverse;
+    inverse.initialize (stiffness_matrix);
+
+    const unsigned int num_arnoldi_vectors = 2*num_eigs + 2
+    ArpackSolver::AdditionalData additional_data(num_arnoldi_vectors);
+
+    ArpackSolver eigensolver (solver_control, additional_data);
+    eigensolver.solve (stiffness_matrix,
+                       mass_matrix,
+                       preconditioner,
+                       eigenvalues,
+                       eigenfunctions,
+                       num_eigs);
+
+    for (unsigned int i=0; i<eigenfunctions.size(); ++i)
+      eigenfunctions[i] /= eigenfunctions[i].linfty_norm ();
+
+    return solver_control.last_step ();
+  }
+@endcode
+Note how we have used an exact decomposition (using SparseDirectUMFPACK) as a
+preconditioner to ARPACK.
 </ul>

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.