]> https://gitweb.dealii.org/ - dealii.git/commitdiff
split off dimension independent stuff
authorGuido Kanschat <dr.guido.kanschat@gmail.com>
Tue, 18 Apr 2000 22:00:13 +0000 (22:00 +0000)
committerGuido Kanschat <dr.guido.kanschat@gmail.com>
Tue, 18 Apr 2000 22:00:13 +0000 (22:00 +0000)
git-svn-id: https://svn.dealii.org/trunk@2746 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/deal.II/include/multigrid/multigrid.templates.h
deal.II/deal.II/source/multigrid/multigrid.all_dimensions.cc [new file with mode: 0644]
deal.II/deal.II/source/multigrid/multigrid.cc

index f93caa8f465ce6e0c6aea137ba313fb52b64f4f8..820d2bdc7d02369eca920ade19e8755e8c674cbc 100644 (file)
@@ -81,6 +81,12 @@ Multigrid<dim>::copy_to_mg (const Vector<number>& src)
       const MGDoFHandler<dim>::active_cell_iterator
        level_end  = mg_dof_handler->end_active(level);
 
+//TODO: Treat hanging nodes properly
+// The single-level vector is not an FE-function, because the values at
+// hanging nodes are set to zero. This should be treated before the restriction.
+
+                                      // Compute coarse level right hand side
+                                      // by restricting from fine level.
       for (; level_cell!=level_end; ++level_cell, ++global_cell)
        {
                                           // get the dof numbers of
@@ -96,8 +102,8 @@ Multigrid<dim>::copy_to_mg (const Vector<number>& src)
          for (unsigned int i=0; i<dofs_per_cell; ++i)
            defect[level](level_dof_indices[i]) = src(global_dof_indices[i]);
 
-//TODO: what happens here?
-                                          // Delete values on refinement edge
+                                          // Delete values on refinement edge,
+                                          // since restriction will add them again.
          for (unsigned int face_n=0; face_n<GeometryInfo<dim>::faces_per_cell; ++face_n)
            {
              const MGDoFHandler<dim>::face_iterator face = level_cell->face(face_n);
diff --git a/deal.II/deal.II/source/multigrid/multigrid.all_dimensions.cc b/deal.II/deal.II/source/multigrid/multigrid.all_dimensions.cc
new file mode 100644 (file)
index 0000000..6a2e53d
--- /dev/null
@@ -0,0 +1,52 @@
+//----------------------------  multigrid.all_dimensions.cc  ---------------------------
+//    $Id$
+//    Version: $Name$
+//
+//    Copyright (C) 1998, 1999, 2000 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.
+//
+//----------------------------  multigrid.all_dimensions.cc  ---------------------------
+
+
+
+#include <grid/tria.h>
+#include <multigrid/mg_dof_handler.h>
+#include <multigrid/mg_dof_accessor.h>
+#include <grid/tria_iterator.h>
+#include <fe/fe.h>
+#include <multigrid/multigrid.h>
+#include <multigrid/multigrid.templates.h>
+#include <multigrid/mg_smoother.h>
+#include <lac/vector.h>
+
+
+MGTransferPrebuilt::~MGTransferPrebuilt () 
+{};
+
+
+void MGTransferPrebuilt::prolongate (const unsigned int   to_level,
+                                    Vector<double>       &dst,
+                                    const Vector<double> &src) const 
+{
+  Assert ((to_level >= 1) && (to_level<=prolongation_matrices.size()),
+         ExcIndexRange (to_level, 1, prolongation_matrices.size()+1));
+
+  prolongation_matrices[to_level-1].vmult (dst, src);
+};
+
+
+void MGTransferPrebuilt::restrict_and_add (const unsigned int   from_level,
+                                          Vector<double>       &dst,
+                                          const Vector<double> &src) const 
+{
+  Assert ((from_level >= 1) && (from_level<=prolongation_matrices.size()),
+         ExcIndexRange (from_level, 1, prolongation_matrices.size()+1));
+
+  prolongation_matrices[from_level-1].Tvmult_add (dst, src);
+};
+
+
index 3dfefa448f3be6f433e80bc4a0b0a6dac8ad486a..4a7c0442af7c84fb6cc9670da43322ece448794a 100644 (file)
@@ -59,10 +59,6 @@ Multigrid<dim>::level_vmult (const unsigned int    level,
 /* ----------------------------- MGTransferPrebuilt ------------------------ */
 
 
-MGTransferPrebuilt::~MGTransferPrebuilt () 
-{};
-
-
 template <int dim>
 void MGTransferPrebuilt::build_matrices (const MGDoFHandler<dim> &mg_dof) 
 {
@@ -174,28 +170,6 @@ void MGTransferPrebuilt::build_matrices (const MGDoFHandler<dim> &mg_dof)
 };
 
 
-void MGTransferPrebuilt::prolongate (const unsigned int   to_level,
-                                    Vector<double>       &dst,
-                                    const Vector<double> &src) const 
-{
-  Assert ((to_level >= 1) && (to_level<=prolongation_matrices.size()),
-         ExcIndexRange (to_level, 1, prolongation_matrices.size()+1));
-
-  prolongation_matrices[to_level-1].vmult (dst, src);
-};
-
-
-void MGTransferPrebuilt::restrict_and_add (const unsigned int   from_level,
-                                          Vector<double>       &dst,
-                                          const Vector<double> &src) const 
-{
-  Assert ((from_level >= 1) && (from_level<=prolongation_matrices.size()),
-         ExcIndexRange (from_level, 1, prolongation_matrices.size()+1));
-
-  prolongation_matrices[from_level-1].Tvmult_add (dst, src);
-};
-
-
 // explicit instatations
 template class Multigrid<deal_II_dimension>;
 template

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.