]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Improve conditioning of project_bv_curl_conforming.
authorMarkus Buerg <buerg@math.tamu.edu>
Thu, 5 Sep 2013 21:26:45 +0000 (21:26 +0000)
committerMarkus Buerg <buerg@math.tamu.edu>
Thu, 5 Sep 2013 21:26:45 +0000 (21:26 +0000)
git-svn-id: https://svn.dealii.org/trunk@30613 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/include/deal.II/numerics/vector_tools.templates.h

index 6044c92994d0f997180e7d1e758be084f422ff28..8acb1e4eacf08d299cb93eed29a9bddb93fec855 100644 (file)
@@ -2787,7 +2787,7 @@ namespace VectorTools
                              std::vector<double> &dof_values,
                              std::vector<bool> &dofs_processed)
     {
-      const double tol = 0.5 * cell->get_fe ().degree * 1e-13 / cell->face (face)->line (line)->diameter ();
+      const double tol = 0.5 * cell->face (face)->line (line)->diameter () / cell->get_fe ().degree;
       const unsigned int dim = 3;
       const unsigned int spacedim = 3;
 
@@ -2991,6 +2991,7 @@ namespace VectorTools
         {
         case 2:
         {
+          const double tol = 0.5 * cell->face (face)->diameter () / cell->get_fe ().degree;
           std::vector<Point<dim> >
           tangentials (fe_values.n_quadrature_points);
 
@@ -3024,18 +3025,18 @@ namespace VectorTools
                 = reference_quadrature_points[q_point];
 
               shifted_reference_point_1 (face_coordinate_direction[face])
-              += 1e-13;
+              += tol;
               shifted_reference_point_2 (face_coordinate_direction[face])
-              -= 1e-13;
+              -= tol;
               tangentials[q_point]
-                = 2e13
-                  * (fe_values.get_mapping ()
+                = (fe_values.get_mapping ()
                      .transform_unit_to_real_cell (cell,
                                                    shifted_reference_point_1)
                      -
                      fe_values.get_mapping ()
                      .transform_unit_to_real_cell (cell,
-                                                   shifted_reference_point_2));
+                                                   shifted_reference_point_2))
+                  / tol;
               tangentials[q_point]
               /= std::sqrt (tangentials[q_point].square ());
               // Compute the degrees
@@ -3412,15 +3413,13 @@ namespace VectorTools
                   // if the degree of
                   // freedom is not
                   // already constrained.
-                  const double tol = 0.5 * cell->get_fe ().degree * 1e-13  / cell->face (face)->diameter ();
-
                   for (unsigned int dof = 0; dof < dofs_per_face; ++dof)
                     if (dofs_processed[dof] && constraints.can_store_line (face_dof_indices[dof])
                         && !(constraints.is_constrained (face_dof_indices[dof])))
                       {
                         constraints.add_line (face_dof_indices[dof]);
 
-                        if (std::abs (dof_values[dof]) > tol)
+                        if (std::abs (dof_values[dof]) > 1e-13)
                           constraints.set_inhomogeneity (face_dof_indices[dof], dof_values[dof]);
                       }
                 }
@@ -3507,19 +3506,16 @@ namespace VectorTools
                   // Store the computed
                   // values in the global
                   // vector.
-
                   cell->face (face)->get_dof_indices (face_dof_indices,
                                                       cell->active_fe_index ());
 
-                  const double tol = 0.5 * superdegree * 1e-13 / cell->face (face)->diameter ();
-
                   for (unsigned int dof = 0; dof < dofs_per_face; ++dof)
                     if (dofs_processed[dof] && constraints.can_store_line (face_dof_indices[dof])
                         && !(constraints.is_constrained (face_dof_indices[dof])))
                       {
                         constraints.add_line (face_dof_indices[dof]);
 
-                        if (std::abs (dof_values[dof]) > tol)
+                        if (std::abs (dof_values[dof]) > 1e-13)
                           constraints.set_inhomogeneity (face_dof_indices[dof], dof_values[dof]);
                       }
                 }
@@ -3611,15 +3607,13 @@ namespace VectorTools
                   cell->face (face)->get_dof_indices (face_dof_indices,
                                                       cell->active_fe_index ());
 
-                  const double tol = 0.5 * cell->get_fe ().degree * 1e-13  / cell->face (face)->diameter ();
-
                   for (unsigned int dof = 0; dof < dofs_per_face; ++dof)
                     if (dofs_processed[dof] && constraints.can_store_line (face_dof_indices[dof])
                         && !(constraints.is_constrained (face_dof_indices[dof])))
                       {
                         constraints.add_line (face_dof_indices[dof]);
 
-                        if (std::abs (dof_values[dof]) > tol)
+                        if (std::abs (dof_values[dof]) > 1e-13)
                           constraints.set_inhomogeneity (face_dof_indices[dof], dof_values[dof]);
                       }
                 }
@@ -3705,15 +3699,13 @@ namespace VectorTools
                   cell->face (face)->get_dof_indices (face_dof_indices,
                                                       cell->active_fe_index ());
 
-                  const double tol = 0.5 * superdegree * 1e-13  / cell->face (face)->diameter ();
-
                   for (unsigned int dof = 0; dof < dofs_per_face; ++dof)
                     if (dofs_processed[dof] && constraints.can_store_line (face_dof_indices[dof])
                         && !(constraints.is_constrained (face_dof_indices[dof])))
                       {
                         constraints.add_line (face_dof_indices[dof]);
 
-                        if (std::abs (dof_values[dof]) > tol)
+                        if (std::abs (dof_values[dof]) > 1e-13)
                           constraints.set_inhomogeneity (face_dof_indices[dof], dof_values[dof]);
                       }
                 }

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.