]> https://gitweb.dealii.org/ - dealii.git/commitdiff
parameters exchanged, flux_pattern in MGDoFTools
authorGuido Kanschat <dr.guido.kanschat@gmail.com>
Tue, 22 Feb 2000 16:37:38 +0000 (16:37 +0000)
committerGuido Kanschat <dr.guido.kanschat@gmail.com>
Tue, 22 Feb 2000 16:37:38 +0000 (16:37 +0000)
git-svn-id: https://svn.dealii.org/trunk@2463 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/deal.II/include/multigrid/mg_dof_tools.h
deal.II/deal.II/source/dofs/dof_renumbering.cc
deal.II/deal.II/source/multigrid/mg_dof_tools.cc
deal.II/deal.II/source/numerics/data_out.cc

index ee586c0aa59f07f46c5ac1dfa0e6c116e652d151..8926c91bb980c3e844370126fb2d30aef955a019 100644 (file)
@@ -8,6 +8,7 @@
 #include <lac/forward_declarations.h>
 #include <grid/forward_declarations.h>
 
+//TODO: Consider incorporating these functions in DoFTools
 
 /**
  * This is a collection of functions operating on, and manipulating
@@ -18,7 +19,7 @@
  * for more information.
  *
  * All member functions are static, so there is no need to create an
- * object of class #DoFTools#.
+ * object of class #MGDoFTools#.
  *
  * @author Wolfgang Bangerth and others, 1999
  */
@@ -28,29 +29,36 @@ class MGDoFTools
                                     /**
                                      * Write the sparsity structure
                                      * of the matrix belonging to the
-                                     * specified #level# including
-                                     * constrained degrees of freedom
-                                     * into the matrix structure. The
-                                     * sparsity pattern does not
-                                     * include entries introduced by
-                                     * the elimination of constrained
-                                     * nodes.  The sparsity pattern
-                                     * is not compressed, since if
-                                     * you want to call
-                                     * #ConstraintMatrix::condense(1)#
-                                     * afterwards, new entries have
-                                     * to be added. However, if you
-                                     * don't want to call
-                                     * #ConstraintMatrix::condense(1)#,
+                                     * specified #level#. The sparsity pattern
+                                     * is not compressed, so before 
+                                     * creating the actual matrix
                                      * you have to compress the
                                      * matrix yourself, using
                                      * #SparseMatrixStruct::compress()#.
+                                     *
+                                     * There is no need to consider
+                                     * hanging nodes here, since only
+                                     * one level is considered.
                                      */
     template <int dim>
     static void
     make_sparsity_pattern (const MGDoFHandler<dim> &dof_handler,
-                          const unsigned int       level,
-                          SparsityPattern         &sparsity);
+                          SparsityPattern         &sparsity,
+                          const unsigned int       level);
+
+                                    /**
+                                     * Make a sparsity pattern including fluxes
+                                     * of discontinuous Galerkin methods.
+                                     * @see make_sparsity_pattern
+                                     * $see DoFTools
+                                     */
+    template <int dim>
+    static void
+    make_flux_sparsity_pattern (const MGDoFHandler<dim> &dof_handler,
+                               SparsityPattern         &sparsity,
+                               const unsigned int       level);
+
+
 };
 
 
index 4cb578d7960b8bcd2c52f6dedf6d7e620aeaff08..6394aa15e2cb76359b471a947b6fa6572dde2c0e 100644 (file)
@@ -217,7 +217,7 @@ void DoFRenumbering::Cuthill_McKee (MGDoFHandler<dim>          &dof_handler,
                                   // make the connection graph
   SparsityPattern sparsity (dof_handler.n_dofs(level),
                            dof_handler.max_couplings_between_dofs());
-  MGDoFTools::make_sparsity_pattern (dof_handler, level, sparsity);
+  MGDoFTools::make_sparsity_pattern (dof_handler, sparsity, level);
     
   const unsigned int n_dofs = sparsity.n_rows();
                                   // store the new dof numbers; invalid_dof_index means
index 3e16b5da0400cd580ba486d4ac88bafef1df22af..7cc6d0b60142186236e04d0263689a8abef2ab7c 100644 (file)
 
 
 template <int dim>
-void MGDoFTools::make_sparsity_pattern (const MGDoFHandler<dim> &mg_dof_handler,
-                                       const unsigned int       level,
-                                       SparsityPattern         &sparsity)
+void MGDoFTools::make_sparsity_pattern (const MGDoFHandler<dim> &dof,
+                                       SparsityPattern         &sparsity,
+                                       const unsigned int       level)
 {
-  Assert (sparsity.n_rows() == mg_dof_handler.n_dofs(level),
-         ExcDimensionMismatch (sparsity.n_rows(), mg_dof_handler.n_dofs(level)));
-  Assert (sparsity.n_cols() == mg_dof_handler.n_dofs(level),
-         ExcDimensionMismatch (sparsity.n_cols(), mg_dof_handler.n_dofs(level)));
+  const unsigned int n_dofs = dof.n_dofs(level);
 
-  const unsigned int dofs_per_cell = mg_dof_handler.get_fe().dofs_per_cell;
+  Assert (sparsity.n_rows() == n_dofs,
+         ExcDimensionMismatch (sparsity.n_rows(), n_dofs));
+  Assert (sparsity.n_cols() == n_dofs,
+         ExcDimensionMismatch (sparsity.n_cols(), n_dofs));
+
+  const unsigned int dofs_per_cell = dof.get_fe().dofs_per_cell;
   vector<unsigned int> dofs_on_this_cell(dofs_per_cell);
-  MGDoFHandler<dim>::cell_iterator cell = mg_dof_handler.begin(level),
-                                  endc = mg_dof_handler.end(level);
+  MGDoFHandler<dim>::cell_iterator cell = dof.begin(level),
+                                  endc = dof.end(level);
   for (; cell!=endc; ++cell) 
     {
       cell->get_mg_dof_indices (dofs_on_this_cell);
@@ -31,13 +33,67 @@ void MGDoFTools::make_sparsity_pattern (const MGDoFHandler<dim> &mg_dof_handler,
        for (unsigned int j=0; j<dofs_per_cell; ++j)
          sparsity.add (dofs_on_this_cell[i],
                        dofs_on_this_cell[j]);
-    };
-};
+    }
+}
+
+template<int dim>
+void
+MGDoFTools::make_flux_sparsity_pattern (const MGDoFHandler<dim> &dof,
+                                       SparsityPattern       &sparsity,
+                                       const unsigned int level)
+{
+  const unsigned int n_dofs = dof.n_dofs(level);
+  
+  Assert (sparsity.n_rows() == n_dofs,
+         ExcDimensionMismatch (sparsity.n_rows(), n_dofs));
+  Assert (sparsity.n_cols() == n_dofs,
+         ExcDimensionMismatch (sparsity.n_cols(), n_dofs));
+
+  const unsigned int dofs_per_cell = dof.get_fe().dofs_per_cell;
+  vector<unsigned int> dofs_on_this_cell(dofs_per_cell);
+  vector<unsigned int> dofs_on_other_cell(dofs_per_cell);
+  DoFHandler<dim>::cell_iterator cell = dof.begin(level),
+                                       endc = dof.end(level);
+  for (; cell!=endc; ++cell)
+    {
+      cell->get_dof_indices (dofs_on_this_cell);
+                                      // make sparsity pattern for this cell
+      for (unsigned int i=0; i<dofs_per_cell; ++i)
+       for (unsigned int j=0; j<dofs_per_cell; ++j)
+         sparsity.add (dofs_on_this_cell[i],
+                       dofs_on_this_cell[j]);
 
+                                      // Loop over all interior neighbors
+      for (unsigned int face = 0;
+          face < GeometryInfo<dim>::faces_per_cell;
+          ++face)
+       {
+         if ( (! cell->at_boundary(face)) && (cell->neighbor_level(face) == level) )
+           {
+             DoFHandler<dim>::cell_iterator neighbor = cell->neighbor(face);
+             neighbor->get_dof_indices (dofs_on_other_cell);
+                                              // only add one direction
+                                              // The other is taken care of
+                                              // by neighbor.
+             for (unsigned int i=0; i<dofs_per_cell; ++i)
+               {
+                 for (unsigned int j=0; j<dofs_per_cell; ++j)
+                   {
+                     sparsity.add (dofs_on_this_cell[i],
+                                   dofs_on_other_cell[j]);
+                   }
+               }
+           }
+       } 
+    }
+}
 
 
 // explicit instantiations
 template void MGDoFTools::make_sparsity_pattern (const MGDoFHandler<deal_II_dimension> &,
-                                                const unsigned int,
-                                                SparsityPattern &);
+                                                SparsityPattern &,
+                                                const unsigned int);
+template void MGDoFTools::make_flux_sparsity_pattern (const MGDoFHandler<deal_II_dimension> &,
+                                                     SparsityPattern &,
+                                                     const unsigned int);
 
index f521a9e7929ea25abe5b1ecb6b19c0171a6e5cf2..c8d83e7f11f036d378c8bf21ce91109b7782e16e 100644 (file)
@@ -85,8 +85,18 @@ template <int dim>
 void DataOut_DoFData<dim>::add_data_vector (const Vector<double> &vec,
                                            const string         &name)
 {
+//  unsigned int n = dofs->get_fe().n_components ();
+//    if (n > 1)
+//      {
+//        vector<string> names (dofs->get_fe().n_components ());
+//        for (unsigned int i=0;i<n;++i)
+//     {
+//       names[i] = name + string("_");
+//     }
+//    }
+//TODO: Murks hier
   add_data_vector (vec, vector<string>(1,name));
-};
+}
 
 
 

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.