]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
shifted matrices in a separate file
authorguido <guido@0785d39b-7218-0410-832d-ea1e28bc413d>
Wed, 27 Jun 2001 17:10:10 +0000 (17:10 +0000)
committerguido <guido@0785d39b-7218-0410-832d-ea1e28bc413d>
Wed, 27 Jun 2001 17:10:10 +0000 (17:10 +0000)
git-svn-id: https://svn.dealii.org/trunk@4766 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/lac/include/lac/eigen.h
deal.II/lac/include/lac/shifted_matrix.h [new file with mode: 0644]

index c2886da3ba73b2f2e37003ff25120d85d157744f..fb43f77a546ddb7c7808bc7c51c346c1581e03cc 100644 (file)
 #define __deal2__eigen_h
 
 
+#include <lac/shifted_matrix.h>
 #include <lac/solver.h>
 #include <lac/solver_control.h>
 #include <lac/solver_cg.h>
 #include <lac/vector_memory.h>
 #include <lac/precondition.h>
 
-/**
- * Matrix with shifted diagonal values.
- *
- * Given a matrix $A$, this class implements a matrix-vector product with
- * $A+\sigma I$, where sigma is a provided shift parameter.
- *
- * @author Guido Kanschat, 2000
- */
-template<class MATRIX>
-class ShiftedMatrix
-{
-  public:
-                                    /**
-                                     * Constructor.
-                                     * Provide the base matrix and a shift parameter.
-                                     */
-    ShiftedMatrix (const MATRIX& A, const double sigma);
-
-                                    /**
-                                     * Set the shift parameter.
-                                     */
-    void shift (const double sigma);
-
-                                    /**
-                                     * Access to the shift parameter.
-                                     */
-    double shift () const;
-
-                                    /**
-                                     * Matrix-vector-product.
-                                     */
-    template <class VECTOR>
-    void vmult (VECTOR& dst, const VECTOR& src) const;
-
-                                    /**
-                                     * Residual.
-                                     */
-    template <class VECTOR>
-    double residual (VECTOR& dst, const VECTOR& src, const VECTOR& rhs) const;
-    
-  private:
-                                    /**
-                                     * Storage for base matrix.
-                                     */
-    const MATRIX& A;
-
-                                    /**
-                                     * Shift parameter.
-                                     */
-    double sigma;
-};
-    
-
 /**
  * Power method (von Mises).
  *
@@ -230,58 +178,6 @@ class EigenInverse : private Solver<VECTOR>
     AdditionalData additional_data;
 };
 
-//----------------------------------------------------------------------//
-
-template <class MATRIX>
-inline
-ShiftedMatrix<MATRIX>::ShiftedMatrix (const MATRIX& A, const double sigma)
-               :
-               A(A), sigma(sigma)
-{}
-
-
-
-template <class MATRIX>
-inline void
-ShiftedMatrix<MATRIX>::shift (const double s)
-{
-  sigma = s;
-}
-
-
-template <class MATRIX>
-inline double
-ShiftedMatrix<MATRIX>::shift () const
-{
-  return sigma;
-}
-
-
-
-template <class MATRIX>
-template <class VECTOR>
-inline void
-ShiftedMatrix<MATRIX>::vmult (VECTOR& dst, const VECTOR& src) const
-{
-  A.vmult(dst, src);
-  dst.add(sigma, src);
-}
-
-
-template <class MATRIX>
-template <class VECTOR>
-inline double
-ShiftedMatrix<MATRIX>::residual (VECTOR& dst,
-                                const VECTOR& src,
-                                const VECTOR& rhs) const
-{
-  A.vmult(dst, src);
-  dst.add(sigma, src);
-  dst.sadd(-1.,1.,rhs);
-  return dst.l2_norm ();
-}
-
-
 //----------------------------------------------------------------------
 
 
diff --git a/deal.II/lac/include/lac/shifted_matrix.h b/deal.II/lac/include/lac/shifted_matrix.h
new file mode 100644 (file)
index 0000000..5784926
--- /dev/null
@@ -0,0 +1,242 @@
+//----------------------------  eigen.h  ---------------------------
+//    $Id$
+//    Version: $Name$
+//
+//    Copyright (C) 2001 by the deal.II authors
+//
+//    This file is subject to QPL and may not be  distributed
+//    without copyright and license information. Please refer
+//    to the file deal.II/doc/license.html for the  text  and
+//    further information on this license.
+//
+//----------------------------  eigen.h  ---------------------------
+#ifndef __deal2__shifted_matrix_h
+#define __deal2__shifted_matrix_h
+
+
+#include <base/smartpointer.h>
+
+/**
+ * Matrix with shifted diagonal values.
+ *
+ * Given a matrix p{A}, this class implements a matrix-vector product with
+ * @p{A+\sigma I}, where sigma is a provided shift parameter.
+ *
+ * @author Guido Kanschat, 2000, 2001
+ */
+template<class MATRIX>
+class ShiftedMatrix
+{
+  public:
+                                    /**
+                                     * Constructor.
+                                     * Provide the base matrix and a shift parameter.
+                                     */
+    ShiftedMatrix (const MATRIX& A, const double sigma);
+
+                                    /**
+                                     * Set the shift parameter.
+                                     */
+    void shift (const double sigma);
+
+                                    /**
+                                     * Access to the shift parameter.
+                                     */
+    double shift () const;
+
+                                    /**
+                                     * Matrix-vector-product.
+                                     */
+    template <class VECTOR>
+    void vmult (VECTOR& dst, const VECTOR& src) const;
+
+                                    /**
+                                     * Residual.
+                                     */
+    template <class VECTOR>
+    double residual (VECTOR& dst, const VECTOR& src, const VECTOR& rhs) const;
+    
+  private:
+                                    /**
+                                     * Storage for base matrix.
+                                     */
+    SmartPointer<const MATRIX> A;
+
+                                    /**
+                                     * Shift parameter.
+                                     */
+    double sigma;
+};
+    
+
+
+/**
+ * Matrix with shifted diagonal values with respect to a certain scalar product.
+ *
+ * Given a matrix @p{A}, this class implements a matrix-vector product
+ * with @p{A+\sigma M}, where sigma is a provided shift parameter and
+ * @p{M} is the matrix representing the identity
+ *
+ * @author Guido Kanschat, 2001
+ */
+template<class MATRIX, class MASSMATRIX>
+class ShiftedMatrixGeneralized
+{
+  public:
+                                    /**
+                                     * Constructor.
+                                     * Provide the base matrix and a shift parameter.
+                                     */
+    ShiftedMatrixGeneralized (const MATRIX& A,
+                             const MASSMATRIX& M,
+                             const double sigma);
+
+                                    /**
+                                     * Set the shift parameter.
+                                     */
+    void shift (const double sigma);
+
+                                    /**
+                                     * Access to the shift parameter.
+                                     */
+    double shift () const;
+
+                                    /**
+                                     * Matrix-vector-product.
+                                     */
+    template <class VECTOR>
+    void vmult (VECTOR& dst, const VECTOR& src) const;
+
+                                    /**
+                                     * Residual.
+                                     */
+    template <class VECTOR>
+    double residual (VECTOR& dst, const VECTOR& src, const VECTOR& rhs) const;
+    
+  private:
+                                    /**
+                                     * Storage for base matrix.
+                                     */
+    SmartPointer<const MATRIX> A;
+                                    /**
+                                     * Storage for mass matrix.
+                                     */
+    SmartPointer<const MASSMATRIX> M;
+
+                                    /**
+                                     * Shift parameter.
+                                     */
+    double sigma;
+};
+    
+
+
+//----------------------------------------------------------------------//
+
+template <class MATRIX>
+inline
+ShiftedMatrix<MATRIX>::ShiftedMatrix (const MATRIX& A, const double sigma)
+               :
+               A(&A), sigma(sigma)
+{}
+
+
+
+template <class MATRIX>
+inline void
+ShiftedMatrix<MATRIX>::shift (const double s)
+{
+  sigma = s;
+}
+
+
+template <class MATRIX>
+inline double
+ShiftedMatrix<MATRIX>::shift () const
+{
+  return sigma;
+}
+
+
+
+template <class MATRIX>
+template <class VECTOR>
+inline void
+ShiftedMatrix<MATRIX>::vmult (VECTOR& dst, const VECTOR& src) const
+{
+  A.vmult(dst, src);
+  dst.add(sigma, src);
+}
+
+
+template <class MATRIX>
+template <class VECTOR>
+inline double
+ShiftedMatrix<MATRIX>::residual (VECTOR& dst,
+                                const VECTOR& src,
+                                const VECTOR& rhs) const
+{
+  A.vmult(dst, src);
+  dst.add(sigma, src);
+  dst.sadd(-1.,1.,rhs);
+  return dst.l2_norm ();
+}
+
+
+//----------------------------------------------------------------------//
+
+template <class MATRIX, class MASSMATRIX>
+inline
+ShiftedMatrixGeneralized<MATRIX, MASSMATRIX>
+::ShiftedMatrixGeneralized (const MATRIX& A,
+                           const MASSMATRIX& M,
+                           const double sigma)
+               :
+               A(&A), M(&M), sigma(sigma)
+{}
+
+
+
+template <class MATRIX, class MASSMATRIX>
+inline void
+ShiftedMatrixGeneralized<MATRIX, MASSMATRIX>::shift (const double s)
+{
+  sigma = s;
+}
+
+
+template <class MATRIX, class MASSMATRIX>
+inline double
+ShiftedMatrixGeneralized<MATRIX, MASSMATRIX>::shift () const
+{
+  return sigma;
+}
+
+
+
+template <class MATRIX, class MASSMATRIX>
+template <class VECTOR>
+inline void
+ShiftedMatrixGeneralized<MATRIX, MASSMATRIX>::vmult (VECTOR& dst,
+                                                    const VECTOR& src) const
+{
+  A.vmult(dst, src);
+  dst.add(sigma, src);
+}
+
+
+template <class MATRIX, class MASSMATRIX>
+template <class VECTOR>
+inline double
+ShiftedMatrixGeneralized<MATRIX, MASSMATRIX>::residual (VECTOR& dst,
+                                                       const VECTOR& src,
+                                                       const VECTOR& rhs) const
+{
+  A.vmult(dst, src);
+  dst.add(sigma, src);
+  dst.sadd(-1.,1.,rhs);
+  return dst.l2_norm ();
+}
+
+
+#endif

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.