]> https://gitweb.dealii.org/ - dealii.git/commitdiff
changed extract_dofs so that it works in parallel
authorMathias Anselmann <mathias.anselmann@posteo.de>
Tue, 4 Jun 2019 11:47:42 +0000 (13:47 +0200)
committerMathias Anselmann <mathias.anselmann@posteo.de>
Tue, 4 Jun 2019 12:31:09 +0000 (14:31 +0200)
source/dofs/dof_tools.cc

index 0cb7e10617ae151a00badf2872ad358ea73e8b17..c7665b4a25e10da469a3f132ed6233f4f888b458 100644 (file)
@@ -451,40 +451,41 @@ namespace DoFTools
                const ComponentMask &                component_mask,
                std::vector<bool> &                  selected_dofs)
   {
-    const FiniteElement<dim, spacedim> &fe = dof.begin_active()->get_fe();
-    (void)fe;
-
-    Assert(component_mask.represents_n_components(fe.n_components()),
+    Assert(component_mask.represents_n_components(
+             dof.get_fe_collection().n_components()),
            ExcMessage(
              "The given component mask is not sized correctly to represent the "
              "components of the given finite element."));
-    Assert(selected_dofs.size() == dof.n_dofs(),
-           ExcDimensionMismatch(selected_dofs.size(), dof.n_dofs()));
+    Assert(selected_dofs.size() == dof.n_locally_owned_dofs(),
+           ExcDimensionMismatch(selected_dofs.size(),
+                                dof.n_locally_owned_dofs()));
 
     // two special cases: no component is selected, and all components are
     // selected; both rather stupid, but easy to catch
-    if (component_mask.n_selected_components(n_components(dof)) == 0)
+    if (component_mask.n_selected_components(
+          dof.get_fe_collection().n_components()) == 0)
       {
-        std::fill_n(selected_dofs.begin(), dof.n_dofs(), false);
+        std::fill_n(selected_dofs.begin(), dof.n_locally_owned_dofs(), false);
         return;
       }
-    else if (component_mask.n_selected_components(n_components(dof)) ==
-             n_components(dof))
+    else if (component_mask.n_selected_components(
+               dof.get_fe_collection().n_components()) ==
+             dof.get_fe_collection().n_components())
       {
-        std::fill_n(selected_dofs.begin(), dof.n_dofs(), true);
+        std::fill_n(selected_dofs.begin(), dof.n_locally_owned_dofs(), true);
         return;
       }
 
 
     // preset all values by false
-    std::fill_n(selected_dofs.begin(), dof.n_dofs(), false);
+    std::fill_n(selected_dofs.begin(), dof.n_locally_owned_dofs(), false);
 
     // get the component association of each DoF and then select the ones
     // that match the given set of components
-    std::vector<unsigned char> dofs_by_component(dof.n_dofs());
+    std::vector<unsigned char> dofs_by_component(dof.n_locally_owned_dofs());
     internal::get_component_association(dof, component_mask, dofs_by_component);
 
-    for (types::global_dof_index i = 0; i < dof.n_dofs(); ++i)
+    for (types::global_dof_index i = 0; i < dof.n_locally_owned_dofs(); ++i)
       if (component_mask[dofs_by_component[i]] == true)
         selected_dofs[i] = true;
   }

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.