]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
Constraint matrices are now also computed on the fly for 3d, so that we don't need...
authorwolf <wolf@0785d39b-7218-0410-832d-ea1e28bc413d>
Wed, 17 Nov 2004 03:43:31 +0000 (03:43 +0000)
committerwolf <wolf@0785d39b-7218-0410-832d-ea1e28bc413d>
Wed, 17 Nov 2004 03:43:31 +0000 (03:43 +0000)
git-svn-id: https://svn.dealii.org/trunk@9782 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/deal.II/include/fe/fe_q.h
deal.II/deal.II/source/fe/fe_q_1d.cc [deleted file]
deal.II/deal.II/source/fe/fe_q_2d.cc [deleted file]
deal.II/deal.II/source/fe/fe_q_3d.cc [deleted file]

index 351aa3a406f8bbb5d36f47cb2ff15f4efd881c5e..fe4899ee2b984d8ebd09ca2f4d13dabbbc39faae 100644 (file)
@@ -305,43 +305,6 @@ class FE_Q : public FE_Poly<TensorProductPolynomials<dim>,dim>
                                      */
     virtual unsigned int memory_consumption () const;
 
-                                    /**
-                                     * Declare a nested class which
-                                     * will hold static definitions
-                                     * of various matrices such as
-                                     * constraint. All other
-                                     * information about this
-                                     * particular finite element,
-                                     * such as embedding matrices,
-                                     * will be computed on-the-fly
-                                     * during construction. The
-                                     * definition of the various
-                                     * static fields are in the files
-                                     * <tt>fe_q_[123]d.cc</tt> in the
-                                     * source directory.
-                                     */
-    struct Matrices
-    {
-                                        /**
-                                         * As the
-                                         * @p embedding_matrices
-                                         * field, but for the
-                                         * interface constraints. One
-                                         * for each element for which
-                                         * it has been computed.
-                                         */
-       static const double * const constraint_matrices[];
-
-                                        /**
-                                         * Like
-                                         * @p n_embedding_matrices,
-                                         * but for the number of
-                                         * interface constraint
-                                         * matrices.
-                                         */
-       static const unsigned int n_constraint_matrices;
-    };
-
   protected:    
                                     /**
                                      * @p clone function instead of
@@ -533,27 +496,4 @@ template <>
 void FE_Q<3>::initialize_constraints ();
 
 
-// declaration of explicit specializations of member variables, if the
-// compiler allows us to do that (the standard says we must)
-#ifndef DEAL_II_MEMBER_VAR_SPECIALIZATION_BUG
-template <>
-const double * const FE_Q<1>::Matrices::constraint_matrices[];
-
-template <>
-const unsigned int FE_Q<1>::Matrices::n_constraint_matrices;
-
-template <>
-const double * const FE_Q<2>::Matrices::constraint_matrices[];
-
-template <>
-const unsigned int FE_Q<2>::Matrices::n_constraint_matrices;
-
-template <>
-const double * const FE_Q<3>::Matrices::constraint_matrices[];
-
-template <>
-const unsigned int FE_Q<3>::Matrices::n_constraint_matrices;
-
-#endif
-
 #endif
diff --git a/deal.II/deal.II/source/fe/fe_q_1d.cc b/deal.II/deal.II/source/fe/fe_q_1d.cc
deleted file mode 100644 (file)
index fac6527..0000000
+++ /dev/null
@@ -1,38 +0,0 @@
-//----------------------------------------------------------------
-//    $Id$
-//    Version: $Name$
-//
-//    Copyright (C) 2001, 2002, 2003 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.
-//
-//----------------------------------------------------------------
-
-
-
-// only compile this file if in 1d
-#if deal_II_dimension == 1
-
-
-#include <fe/fe_q.h>
-
-
-// No constraints in 1d
-template <>
-const unsigned int 
-FE_Q<1>::Matrices::n_constraint_matrices = 0;
-
-
-template <>
-const double * const
-FE_Q<1>::Matrices::constraint_matrices[] = { 0 };
-
-
-#else // #if deal_II_dimension
-// On gcc2.95 on Alpha OSF1, the native assembler does not like empty
-// files, so provide some dummy code
-namespace { void dummy () {} }
-#endif // #if deal_II_dimension == 1
diff --git a/deal.II/deal.II/source/fe/fe_q_2d.cc b/deal.II/deal.II/source/fe/fe_q_2d.cc
deleted file mode 100644 (file)
index ec34cef..0000000
+++ /dev/null
@@ -1,42 +0,0 @@
-//----------------------------------------------------------------
-//    $Id$
-//    Version: $Name$
-//
-//    Copyright (C) 2001, 2002, 2003, 2004 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.
-//
-//----------------------------------------------------------------
-
-
-// only compile this file if in 2d
-#if deal_II_dimension == 2
-
-
-#include <fe/fe_q.h>
-
-// constraint matrices in 2d are now implemented by computing them on the fly
-// for all polynomial degrees. the array is thus empty. unfortunately, some
-// compilers dislike empty initializers for arrays of unknown size
-// (particularly the hp compiler), so we simply initialize a single element
-// with a null pointer. in fact, the compiler behavior is supported by clause
-// 8.5.1/4
-template <>
-const double * const 
-FE_Q<2>::Matrices::constraint_matrices[] = { 0 };
-
-
-template <>
-const unsigned int 
-FE_Q<2>::Matrices::n_constraint_matrices = 0;
-
-
-
-#else // #if deal_II_dimension
-// On gcc2.95 on Alpha OSF1, the native assembler does not like empty
-// files, so provide some dummy code
-namespace { void dummy () {} }
-#endif // #if deal_II_dimension == 2
diff --git a/deal.II/deal.II/source/fe/fe_q_3d.cc b/deal.II/deal.II/source/fe/fe_q_3d.cc
deleted file mode 100644 (file)
index 6b0de5a..0000000
+++ /dev/null
@@ -1,85 +0,0 @@
-//----------------------------------------------------------------
-//    $Id$
-//    Version: $Name$
-//
-//    Copyright (C) 2001, 2002, 2003 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.
-//
-//----------------------------------------------------------------
-
-// Transfer matrices for finite elements
-
-
-// only compile this file if in 3d
-#if deal_II_dimension == 3
-
-#include <fe/fe_q.h>
-
-// Constraint matrices taken from Wolfgangs old version
-
-namespace FE_Q_3d 
-{
-  static const double constraint_q1[] =
-  {
-       .25,.25,.25,.25,
-       .5,.5,0.,0.,
-       0.,.5,.5,0.,
-       0.,0.,.5,.5,
-       .5,0.,0.,.5
-  };
-
-  static const double constraint_q2[] =
-  {
-       0,0,0,0,0,0,0,0,1,
-       0,0,0,0,1,0,0,0,0,
-       0,0,0,0,0,1,0,0,0,
-       0,0,0,0,0,0,1,0,0,
-       0,0,0,0,0,0,0,1,0,
-       0,0,0,0,.375,0,-.125,0,.75,
-       0,0,0,0,0,.375,0,-.125,.75,
-       0,0,0,0,-.125,0,.375,0,.75,
-       0,0,0,0,0,-.125,0,.375,.75,
-       .375,-.125,0,0,.75,0,0,0,0,
-       -.125,.375,0,0,.75,0,0,0,0,
-       0,.375,-.125,0,0,.75,0,0,0,
-       0,-.125,.375,0,0,.75,0,0,0,
-       0,0,-.125,.375,0,0,.75,0,0,
-       0,0,.375,-.125,0,0,.75,0,0,
-       .375,0,0,-.125,0,0,0,.75,0,
-       -.125,0,0,.375,0,0,0,.75,0,
-       .140625,-.046875,.015625,-.046875,.28125,-.09375,-.09375,.28125,.5625,
-       -.046875,.140625,-.046875,.015625,.28125,.28125,-.09375,-.09375,.5625,
-       .015625,-.046875,.140625,-.046875,-.09375,.28125,.28125,-.09375,.5625,
-       -.046875,.015625,-.046875,.140625,-.09375,-.09375,.28125,.28125,.5625
-  };
-}
-
-
-
-template <>
-const double * const 
-FE_Q<3>::Matrices::constraint_matrices[] =
-{
-  FE_Q_3d::constraint_q1,
-  FE_Q_3d::constraint_q2
-};
-
-
-
-template <>
-const unsigned int 
-FE_Q<3>::Matrices::n_constraint_matrices
-  = sizeof(FE_Q<3>::Matrices::constraint_matrices) /
-    sizeof(FE_Q<3>::Matrices::constraint_matrices[0]);
-
-
-
-#else // #if deal_II_dimension
-// On gcc2.95 on Alpha OSF1, the native assembler does not like empty
-// files, so provide some dummy code
-namespace { void dummy () {} }
-#endif // #if deal_II_dimension == 3

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.