]> https://gitweb.dealii.org/ - dealii.git/commitdiff
enable computation of embedding matrices for elements not supporting anisotropic...
authorGuido Kanschat <dr.guido.kanschat@gmail.com>
Wed, 17 Dec 2008 21:17:45 +0000 (21:17 +0000)
committerGuido Kanschat <dr.guido.kanschat@gmail.com>
Wed, 17 Dec 2008 21:17:45 +0000 (21:17 +0000)
git-svn-id: https://svn.dealii.org/trunk@17980 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/deal.II/include/fe/fe_tools.h
deal.II/deal.II/source/fe/fe.cc
deal.II/deal.II/source/fe/fe_tools.cc

index 956b945acc4533602f8beb4f486e607e6bc2e64b..2b96003a8560923b94409a2cbb757ec8821a0143 100644 (file)
@@ -358,10 +358,12 @@ class FETools
                                      * FiniteElement::prolongation
                                      * matrices.
                                      *
-                                     * @param fe The finite element
+                                     * @param <tt>fe</tt> The finite element
                                      * class for which we compute the
-                                     * embedding matrices.  @param
-                                     * matrices A reference to
+                                     * embedding matrices.
+                                     *
+                                     * @param
+                                     * <tt>matrices</tt> A reference to
                                      * <tt>RefinementCase<dim>::isotropic_refinement</tt>
                                      * vectors of FullMatrix
                                      * objects. Each vector
@@ -373,10 +375,17 @@ class FETools
                                      * is the format used in
                                      * FiniteElement, where we want
                                      * to use this function mostly.
+                                     *
+                                     * @param
+                                     * <tt>isotropic_only</tt>: set
+                                     * this <tt>true</tt> if you only
+                                     * want to compute matrices for
+                                     * isotropic refinement.
                                      */
     template <int dim, typename number, int spacedim>
     static void compute_embedding_matrices(const FiniteElement<dim,spacedim> &fe,
-                                          std::vector<std::vector<FullMatrix<number> > >& matrices);
+                                          std::vector<std::vector<FullMatrix<number> > >& matrices,
+                                          const bool isotropic_only = false);
 
                                     /**
                                      * Compute the embedding matrices
index 7fa616e2e5e188cf7cdf5388c01294a40467a287..0e17e536dbc4cd7a124c3dcfa4094591b032ac9f 100644 (file)
@@ -313,7 +313,7 @@ FiniteElement<dim,spacedim>::reinit_restriction_and_prolongation_matrices (
   const bool isotropic_prolongation_only)
 {
   for (unsigned int ref_case=RefinementCase<dim>::cut_x;
-       ref_case<RefinementCase<dim>::isotropic_refinement+1; ++ref_case)
+       ref_case <= RefinementCase<dim>::isotropic_refinement; ++ref_case)
     {
       const unsigned int nc = GeometryInfo<dim>::n_children(RefinementCase<dim>(ref_case));
       
index 214f8d993a8e3d8cdc5e4aebeb1436b9d0bf71ac..e2a6d91b412f35cba254d04280dd50067a49a1b7 100644 (file)
@@ -583,7 +583,8 @@ FETools::compute_embedding_matrices(const FiniteElement<2,3>&,
 template <int dim, typename number, int spacedim>
 void
 FETools::compute_embedding_matrices(const FiniteElement<dim,spacedim>& fe,
-                                   std::vector<std::vector<FullMatrix<number> > >& matrices)
+                                   std::vector<std::vector<FullMatrix<number> > >& matrices,
+                                   const bool isotropic_only)
 {
 
   const unsigned int n  = fe.dofs_per_cell;
@@ -591,8 +592,11 @@ FETools::compute_embedding_matrices(const FiniteElement<dim,spacedim>& fe,
   const unsigned int degree = fe.degree;
 
                                   // loop over all possible refinement cases
-  for (unsigned int ref_case=RefinementCase<dim>::cut_x;
-       ref_case<RefinementCase<dim>::isotropic_refinement+1; ++ref_case)
+  unsigned int ref_case = (isotropic_only)
+                         ? RefinementCase<dim>::isotropic_refinement
+                         : RefinementCase<dim>::cut_x;
+  
+  for (;ref_case <= RefinementCase<dim>::isotropic_refinement; ++ref_case)
     {
       const unsigned int nc = GeometryInfo<dim>::n_children(RefinementCase<dim>(ref_case));
       for (unsigned int i=0;i<nc;++i)
@@ -1967,7 +1971,7 @@ void FETools::get_projection_matrix<deal_II_dimension>
 
 template
 void FETools::compute_embedding_matrices<deal_II_dimension>
-(const FiniteElement<deal_II_dimension> &, std::vector<std::vector<FullMatrix<double> > >&);
+(const FiniteElement<deal_II_dimension> &, std::vector<std::vector<FullMatrix<double> > >&,bool);
 
 template
 void FETools::compute_face_embedding_matrices<deal_II_dimension,double>

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.