]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Remove the VectorBase class, which was only needed for the old way of handling multigrid.
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Fri, 4 Sep 1998 11:23:38 +0000 (11:23 +0000)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Fri, 4 Sep 1998 11:23:38 +0000 (11:23 +0000)
git-svn-id: https://svn.dealii.org/trunk@557 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/lac/include/lac/dvector.h
deal.II/lac/include/lac/vectorbase.h [deleted file]
deal.II/lac/source/dvector.cc

index 7c72d92f19707f022d2ec334a5d293f5dde2ccb3..2ed3a1452a1472a733edacf55cdef5e6d3314c96 100644 (file)
@@ -362,63 +362,6 @@ class dVector {
                                     //@}
     
     
-                                    /**
-                                     * @name 4: Modification of components
-                                     */
-                                    //@{
-                                    /**
-                                     *  U(i)=0             . ONE Component only
-                                     */
-    void czero (const unsigned int i);
-    
-                                    /**
-                                     *  U(i)=a*V(j)        . Replacing
-                                     */
-    void cequ(unsigned int i, const dVector& V,
-             double a, unsigned int j);
-    
-                                    /**
-                                     *  U(i)=a*V(j)+b*V(k) . Replacing by sum
-                                     */
-    void cequ (const unsigned int i, const dVector& V,
-              const double a, const unsigned int j,
-              const double b, const unsigned int k);
-    
-                                    /**
-                                     * U(i)=a*V(j)+b*V(k)+c*V(l)+d*V(m).
-                                     * Replacing by sum
-                                     */
-    void cequ (const unsigned int i, const dVector& V,
-              const double a, const unsigned int j,
-              const double b, const unsigned int k,
-              const double c, const unsigned int l,
-              const double d, const unsigned int m);
-    
-                                    /**
-                                     *  U(i)+=a*V(j)       . Simple addition
-                                     */
-    void cadd (const unsigned int i, const dVector& V,
-              const double a, const unsigned int j);
-    
-                                    /**
-                                     *  U(i)+=a*V(j)+b*V(k). Multiple addition
-                                     */ 
-    void cadd (const unsigned int i, const dVector& V,
-              const double a, const unsigned int j,
-              const double b, const unsigned int k);
-    
-                                    /**
-                                     * U(i)+=a*V(j)+b*V(k)+c*V(l)+d*V(m).
-                                     * Multiple addition
-                                     */
-    void cadd (const unsigned int i, const dVector& V,
-              const double a, const unsigned int j,
-              const double b, const unsigned int k,
-              const double c, const unsigned int l,
-              const double d, const unsigned int m);
-                                    //@}
-    
-    
                                     /**
                                      * @name 5: Mixed stuff
                                      */
diff --git a/deal.II/lac/include/lac/vectorbase.h b/deal.II/lac/include/lac/vectorbase.h
deleted file mode 100644 (file)
index c69d9bb..0000000
+++ /dev/null
@@ -1,95 +0,0 @@
-// $Id$
-
-// This file is part of the DEAL Library
-// DEAL is Copyright(1995) by
-// Roland Becker, Guido Kanschat, Franz-Theo Suttmeier
-
-
-#ifndef __lac_vectorbase_h
-#define __lac_vectorbase_h
-
-
-
-/**
- * Vector Baseclass (abstract).
- *  CONVENTIONS for used `equations` : <p>
- *  - THIS vector is always named `U` <p>
- *  - vectors are always uppercase , scalars are lowercase
- */
-class VectorBase
-{
-  public:
-    
-                                    /**@name 1: Basic Object-handling
-                                     */
-                                    //@{
-                                    /**
-                                     *           Destructor. Should clear memory
-                                     */
-    virtual ~VectorBase() {}
-                                    //@}
-    
-    
-                                    /**@name 2: Modification of components
-                                     */
-                                    //@{
-                                    /**
-                                     *  U(i)=0             . ONE Component only
-                                     */
-    virtual void czero (const unsigned int i) = 0;
-    
-                                    /**
-                                     *  U(i)=a*V(j)        . Replacing
-                                     */
-    virtual void cequ (const unsigned int i, const VectorBase& V,
-                      const double a, const unsigned int j) = 0;
-    
-                                    /**
-                                     *  U(i)=a*V(j)+b*V(k) . Replacing by sum
-                                     */
-    virtual void cequ (const unsigned int i, const VectorBase& V,
-                      const double a, const unsigned int j,
-                      const double b, const unsigned int k) = 0;
-    
-                                    /**
-                                     *  U(i)=a*V(j)+b*V(k)+c*V(l)+d*V(m) . Replacing by sum
-                                     */
-    virtual void cequ (const unsigned int i, const VectorBase& V,
-                      const double a, const unsigned int j,
-                      const double b, const unsigned int k,
-                      const double c, const unsigned int l,
-                      const double d, const unsigned int m) = 0;
-    
-                                    /**
-                                     *  U(i)+=a*V(j)       . Simple addition
-                                     */
-    virtual void cadd (const unsigned int i, const VectorBase& V,
-                      const double a, const unsigned int j) = 0;
-    
-                                    /**
-                                     *  U(i)+=a*V(j)+b*V(k). Multiple addition
-                                     */
-    virtual void cadd (const unsigned int i, const VectorBase& V,
-                      const double a, const unsigned int j,
-                      const double b, const unsigned int k) = 0;
-    
-                                    /**
-                                     *  U(i)+=a*V(j)+b*V(k)+c*V(l)+d*V(m) . Multiple addition
-                                     */
-    virtual void cadd (const unsigned int i, const VectorBase& V,
-                      const double a, const unsigned int j,
-                      const double b, const unsigned int k,
-                      const double c, const unsigned int l,
-                      const double d, const unsigned int m) = 0;
-                                    //@}
-    
-    
-                                    /**@name 3: Mixed stuff
-                                     */
-                                    //@{
-                                    ///
-    virtual const char* name() const = 0;
-                                    //@}
-};
-
-#endif
index 05febceda76b77cfbd747ce4453bc37952bfaacb..6001a6b4a10d3284492f2e2321bf1ca6a9b23bc7 100644 (file)
@@ -485,104 +485,6 @@ dVector& dVector::operator = (const dVector& v)
 
 
 
-void dVector::cadd (const unsigned int i, const dVector& v,
-                   const double s, const unsigned int j)
-{
-  Assert (dim!=0, ExcEmptyVector());
-  Assert (i<dim, ExcInvalidIndex(i,dim));
-  Assert (j<v.dim, ExcInvalidIndex(j,v.dim));
-
-  val[i] += s*v.val[j];
-}
-
-
-
-void dVector::cadd(unsigned int i, const dVector &v,
-                  double s, unsigned int j,
-                  double t, unsigned int k)
-{
-  Assert (dim!=0, ExcEmptyVector());
-  Assert (i<dim, ExcInvalidIndex(i,dim));
-  Assert (j<v.dim, ExcInvalidIndex(j,v.dim));
-  Assert (k<v.dim, ExcInvalidIndex(k,v.dim));
-
-  val[i] += s*v.val[j] + t*v.val[k];
-}
-
-
-
-void dVector::cadd (const unsigned int i, const dVector &v,
-                   const double s, const unsigned int j,
-                   const double t, const unsigned int k,
-                   const double q, const unsigned int l,
-                   const double r, const unsigned int m)
-{
-  Assert (dim!=0, ExcEmptyVector());
-  Assert (i<dim, ExcInvalidIndex(i,dim));
-  Assert (j<v.dim, ExcInvalidIndex(j,v.dim));
-  Assert (k<v.dim, ExcInvalidIndex(k,v.dim));
-  Assert (l<v.dim, ExcInvalidIndex(l,v.dim));
-  Assert (m<v.dim, ExcInvalidIndex(m,v.dim));
-
-  val[i] += s*v.val[j] + t*v.val[k] + q*v.val[l] + r*v.val[m];
-}
-
-
-
-void dVector::czero (const unsigned int i)
-{
-  Assert (dim!=0, ExcEmptyVector());
-  Assert (i<dim, ExcInvalidIndex(i,dim));
-  val[i] = 0.;
-}
-
-
-
-void dVector::cequ (const unsigned int i, const dVector &v,
-                   const double s, const unsigned int j)
-{
-  Assert (dim!=0, ExcEmptyVector());
-  Assert (i<dim, ExcInvalidIndex(i,dim));
-  Assert (j<v.dim, ExcInvalidIndex(j,v.dim));
-
-  val[i] = s*v.val[j];
-}
-
-
-
-void dVector::cequ (const unsigned int i, const dVector &v,
-                   const double s, const unsigned int j,
-                   const double t, const unsigned int k)
-{
-  Assert (dim!=0, ExcEmptyVector());
-  Assert (i<dim, ExcInvalidIndex(i,dim));
-  Assert (j<v.dim, ExcInvalidIndex(j,v.dim));
-  Assert (k<v.dim, ExcInvalidIndex(k,v.dim));
-
-  val[i] = s*v.val[j] + t*v.val[k];
-}
-
-
-
-void dVector::cequ (const unsigned int i, const dVector &v,
-                   const double s, const unsigned int j,
-                   const double t, const unsigned int k,
-                   const double q, const unsigned int l,
-                   const double r, const unsigned int m)
-{
-  Assert (dim!=0, ExcEmptyVector());
-  Assert (i<dim, ExcInvalidIndex(i,dim));
-  Assert (j<v.dim, ExcInvalidIndex(j,v.dim));
-  Assert (k<v.dim, ExcInvalidIndex(k,v.dim));
-  Assert (l<v.dim, ExcInvalidIndex(l,v.dim));
-  Assert (m<v.dim, ExcInvalidIndex(m,v.dim));
-
-  val[i] = s*v.val[j] + t*v.val[k] + q*v.val[l] + r*v.val[m];
-}
-
-
-
-
 void dVector::print (FILE* f, const char* format) const
 {
   Assert (dim!=0, ExcEmptyVector());

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.