]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Augment documentation.
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Wed, 20 Mar 2013 23:27:52 +0000 (23:27 +0000)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Wed, 20 Mar 2013 23:27:52 +0000 (23:27 +0000)
git-svn-id: https://svn.dealii.org/trunk@28964 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/include/deal.II/lac/arpack_solver.h

index 79cf542f8a68b4b8ca9319b6ade8d8c6341c1436..75d6888072534840c436ed8f0f584b89acf9a347 100644 (file)
 
 DEAL_II_NAMESPACE_OPEN
 
+
+extern "C" void dnaupd_(int *ido, char *bmat, const unsigned int *n, char *which,
+                        const unsigned int *nev, const double *tol, double *resid, int *ncv,
+                        double *v, int *ldv, int *iparam, int *ipntr,
+                        double *workd, double *workl, int *lworkl,
+                        int *info);
+
+extern "C" void dneupd_(int *rvec, char *howmany, int *select, double *d,
+                        double *di, double *z, int *ldz, double *sigmar,
+                        double *sigmai, double *workev, char *bmat,const unsigned int *n, char *which,
+                        const unsigned int *nev, const double *tol, double *resid, int *ncv,
+                        double *v, int *ldv, int *iparam, int *ipntr,
+                        double *workd, double *workl, int *lworkl, int *info);
+
 /**
- * Interface for using ARPACK. ARPACK is a collection of Fortran77
- * subroutines designed to solve large scale eigenvalue problems.
- * Here we interface to the routines dnaupd and dneupd of ARPACK.
- * The package is designed to compute a few eigenvalues and
- * corresponding eigenvectors of a general n by n matrix A. It is
- * most appropriate for large sparse matrices A.
+ * Interface for using ARPACK. ARPACK is a collection of Fortran77 subroutines
+ * designed to solve large scale eigenvalue problems.  Here we interface to
+ * the routines <code>dneupd</code> and <code>dnaupd</code> of ARPACK.  The
+ * package is designed to compute a few eigenvalues and corresponding
+ * eigenvectors of a general n by n matrix A. It is most appropriate for large
+ * sparse matrices A.
  *
  * In this class we make use of the method applied to the
  * generalized eigenspectrum problem $(A-\lambda B)x=0$, for
@@ -39,40 +53,25 @@ DEAL_II_NAMESPACE_OPEN
  @code
   SolverControl solver_control (1000, 1e-9);
   ArpackSolver (solver_control);
-  system.solve (A, B, lambda, x, size_of_spectrum);
+  system.solve (A, B, P, lambda, x, size_of_spectrum);
  @endcode
  * for the generalized eigenvalue problem $Ax=B\lambda x$, where
- * the variable <code>const unsigned int size_of_spectrum</code>
+ * the variable <code>size_of_spectrum</code>
  * tells ARPACK the number of eigenvector/eigenvalue pairs to
- * solve for.
+ * solve for. Here, <code>lambda</code> is a vector that will contain
+ * the eigenvalues computed, <code>x</code> a vector that will
+ * contain the eigenvectors computed, and <code>P</code> is
+ * a preconditioner for the matrix <code>A</code>.
  *
  * Through the AdditionalData the user can specify some of the
  * parameters to be set.
  *
- * For further information on how the ARPACK routines dnaupd and
- * dneupd work and also how to set the parameters appropriately
+ * For further information on how the ARPACK routines <code>dneupd</code> and
+ * <code>dnaupd</code> work and also how to set the parameters appropriately
  * please take a look into the ARPACK manual.
  *
- * @author Bรคrbel Janssen, Agnieszka Miedlar, 2010.
+ * @author Baerbel Janssen, Agnieszka Miedlar, 2010.
  */
-
-extern "C" void dnaupd_(int *ido, char *bmat, const unsigned int *n, char *which,
-                        const unsigned int *nev, const double *tol, double *resid, int *ncv,
-                        double *v, int *ldv, int *iparam, int *ipntr,
-                        double *workd, double *workl, int *lworkl,
-                        int *info);
-
-extern "C" void dneupd_(int *rvec, char *howmany, int *select, double *d,
-                        double *di, double *z, int *ldz, double *sigmar,
-                        double *sigmai, double *workev, char *bmat,const unsigned int *n, char *which,
-                        const unsigned int *nev, const double *tol, double *resid, int *ncv,
-                        double *v, int *ldv, int *iparam, int *ipntr,
-                        double *workd, double *workl, int *lworkl, int *info);
-
-/**
- * Class to interface with the ARPACK routines.
- */
-
 class ArpackSolver : public Subscriptor
 {
 public:
@@ -123,7 +122,8 @@ public:
   /**
    * Solve the generalized eigensprectrum
    * problem $A x=\lambda B x$ by calling
-   * dneupd and dnaupd of ARPACK.
+   * the <code>dneupd</code> and <code>dnaupd</code>
+   * functions of ARPACK.
    */
   template <typename VECTOR, typename MATRIX1,
            typename MATRIX2, typename INVERSE>

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.