25 #include "libmesh/libmesh_config.h" 
   27 #ifdef LIBMESH_HAVE_METAPHYSICL 
   31 #include "libmesh/libmesh_common.h" 
   34 #include "libmesh/ignore_warnings.h" 
   37 #include "metaphysicl/dynamicsparsenumberarray_decl.h" 
   38 #include "libmesh/restore_warnings.h" 
   40 #include "libmesh/compare_types.h" 
   41 #include "libmesh/int_range.h" 
   43 using MetaPhysicL::DynamicSparseNumberArray;
 
   49 template <
typename T, 
typename IndexType>
 
   57 template <
typename T, 
typename IndexType, 
typename T2>
 
   61   MetaPhysicL::DynamicSparseNumberArray
 
   67 template <
typename T, 
typename IndexType>
 
   69   typedef std::vector<std::pair<IndexType,T>> 
type;
 
   72 template <
typename T, 
typename IndexType>
 
   73 const std::vector<std::pair<IndexType,T>>
 
   76   const std::size_t in_size = in.size();
 
   77   std::vector<std::pair<IndexType,T>> returnval(in_size);
 
   79   for (std::size_t i=0; i != in_size; ++i)
 
   81       returnval[i].first = in.raw_index(i);
 
   82       returnval[i].second = in.raw_at(i);
 
   87 template <
typename SendT, 
typename T, 
typename IndexType>
 
   89                            MetaPhysicL::DynamicSparseNumberArray<T,IndexType> & converted)
 
   91   const std::size_t received_size = received.size();
 
   92   converted.resize(received_size);
 
   93   for (std::size_t i=0; i != received_size; ++i)
 
   95       converted.raw_index(i) = received[i].first;
 
   96       converted.raw_at(i) = received[i].second;
 
  105 #include "libmesh/boundary_info.h" 
  106 #include "libmesh/dense_matrix.h" 
  107 #include "libmesh/dense_vector.h" 
  108 #include "libmesh/dof_map.h" 
  109 #include "libmesh/elem.h" 
  110 #include "libmesh/fe_base.h" 
  111 #include "libmesh/fe_interface.h" 
  112 #include "libmesh/generic_projector.h" 
  113 #include "libmesh/libmesh_logging.h" 
  114 #include "libmesh/mesh_base.h" 
  115 #include "libmesh/numeric_vector.h" 
  116 #include "libmesh/quadrature.h" 
  117 #include "libmesh/sparse_matrix.h" 
  118 #include "libmesh/system.h" 
  119 #include "libmesh/threads.h" 
  120 #include "libmesh/wrapped_function.h" 
  121 #include "libmesh/wrapped_functor.h" 
  125 #ifdef LIBMESH_HAVE_METAPHYSICL 
  127 #include "libmesh/ignore_warnings.h" 
  129 #include "metaphysicl/dynamicsparsenumberarray.h" 
  130 #include "libmesh/restore_warnings.h" 
  133 #include "libmesh/dense_matrix_impl.h" 
  136 typedef DynamicSparseNumberArray<Real, dof_id_type> 
DSNAN;
 
  155 #ifdef LIBMESH_ENABLE_AMR 
  179     system(other.system),
 
  189 #endif // LIBMESH_ENABLE_AMR 
  200   const std::set<boundary_id_type> & 
b;
 
  203   std::unique_ptr<FunctionBase<Number>>   
f;
 
  204   std::unique_ptr<FunctionBase<Gradient>> 
g;
 
  210                            const std::vector<unsigned int> & variables_in,
 
  217     variables(variables_in),
 
  221     parameters(parameters_in),
 
  232     variables(in.variables),
 
  236     parameters(in.parameters),
 
  237     new_vector(in.new_vector)
 
  253                              int is_adjoint)
 const 
  257   std::unique_ptr<NumericVector<Number>>
 
  258     old_vector (vector.
clone());
 
  261   this->project_vector (*old_vector, vector, is_adjoint);
 
  272                              int is_adjoint)
 const 
  274   LOG_SCOPE (
"project_vector(old,new)", 
"System");
 
  284 #ifdef LIBMESH_ENABLE_AMR 
  288   std::unique_ptr<NumericVector<Number>> new_vector_built;
 
  290   std::unique_ptr<NumericVector<Number>> local_old_vector_built;
 
  294     (this->get_mesh().active_local_elements_begin(),
 
  295      this->get_mesh().active_local_elements_end());
 
  302       new_vector_ptr = &new_v;
 
  303       old_vector_ptr = &old_v;
 
  318       new_v.
init (this->n_dofs(), this->n_local_dofs(), 
false, 
PARALLEL);
 
  321       new_vector_ptr = new_vector_built.
get();
 
  322       local_old_vector = local_old_vector_built.
get();
 
  323       new_vector_ptr->
init(this->n_dofs(), 
false, 
SERIAL);
 
  326       local_old_vector->
close();
 
  327       old_vector_ptr = local_old_vector;
 
  339       new_v.
init (this->n_dofs(), this->n_local_dofs(),
 
  340                   this->get_dof_map().get_send_list(), 
false, 
GHOSTED);
 
  343       new_vector_ptr = &new_v;
 
  344       local_old_vector = local_old_vector_built.
get();
 
  348       local_old_vector->
close();
 
  349       old_vector_ptr = local_old_vector;
 
  352     libmesh_error_msg(
"ERROR: Unknown old_v.type() == " << old_v.
type());
 
  363   const unsigned int n_variables = this->
n_vars();
 
  367       std::vector<unsigned int> vars(n_variables);
 
  380       FEMProjector projector(*
this, f, &g, setter, vars);
 
  381       projector.project(active_local_elem_range);
 
  386       if (this->processor_id() == (this->n_processors()-1))
 
  388           const DofMap & dof_map = this->get_dof_map();
 
  390             if (this->variable(var).type().family == 
SCALAR)
 
  393                 std::vector<dof_id_type> new_SCALAR_indices, old_SCALAR_indices;
 
  397                   new_vector.
set(new_SCALAR_indices[i], old_vector(old_SCALAR_indices[i]));
 
  412       dist_v->init(this->n_dofs(), this->n_local_dofs(), 
false, 
PARALLEL);
 
  416         if (new_vector(i) != 0.0)
 
  417           dist_v->set(i, new_vector(i));
 
  421       dist_v->localize (new_v, this->get_dof_map().get_send_list());
 
  432         if (new_vector(i) != 0.0)
 
  433           new_v.
set(i, new_vector(i));
 
  437   if (is_adjoint == -1)
 
  438     this->get_dof_map().enforce_constraints_exactly(*
this, &new_v);
 
  439   else if (is_adjoint >= 0)
 
  440     this->get_dof_map().enforce_adjoint_constraints_exactly(new_v,
 
  450 #endif // #ifdef LIBMESH_ENABLE_AMR 
  454 #ifdef LIBMESH_ENABLE_AMR 
  455 #ifdef LIBMESH_HAVE_METAPHYSICL 
  457 template <
typename Output>
 
  461   typedef DynamicSparseNumberArray<Output, dof_id_type> 
type;
 
  464 template <
typename InnerOutput>
 
  480 template <
typename Output,
 
  493     this->old_context.set_algebraic_type(FEMContext::OLD_DOFS_ONLY);
 
  499     this->old_context.set_algebraic_type(FEMContext::OLD_DOFS_ONLY);
 
  504                      unsigned int elem_dim,
 
  506                      bool extra_hanging_dofs,
 
  515                       unsigned int node_num,
 
  516                       unsigned int var_num,
 
  517                       std::vector<dof_id_type> & indices,
 
  518                       std::vector<DSNA> & values)
 
  520     LOG_SCOPE (
"eval_old_dofs(node)", 
"OldSolutionCoefs");
 
  522     this->sys.get_dof_map().dof_indices(elem, node_num, indices, var_num);
 
  524     std::vector<dof_id_type> old_indices;
 
  526     this->sys.get_dof_map().old_dof_indices(elem, node_num, old_indices, var_num);
 
  528     libmesh_assert_equal_to (old_indices.size(), indices.size());
 
  530     values.resize(old_indices.size());
 
  535         values[i].raw_at(0) = 1;
 
  536         values[i].raw_index(0) = old_indices[i];
 
  543                       unsigned int sys_num,
 
  544                       unsigned int var_num,
 
  545                       std::vector<dof_id_type> & indices,
 
  546                       std::vector<DSNA> & values)
 
  548     LOG_SCOPE (
"eval_old_dofs(elem)", 
"OldSolutionCoefs");
 
  552     const Elem & old_elem =
 
  557     const unsigned int nc = FEInterface::n_dofs_per_elem(elem.
dim(),
 
  561     std::vector<dof_id_type> old_dof_indices(nc);
 
  571         const std::pair<unsigned int, unsigned int>
 
  573         const unsigned int vg = vg_and_offset.first;
 
  574         const unsigned int vig = vg_and_offset.second;
 
  576         const unsigned int n_comp = elem.
n_comp_group(sys_num,vg);
 
  577         libmesh_assert_greater(elem.
n_systems(), sys_num);
 
  578         libmesh_assert_greater_equal(n_comp, nc);
 
  580         for (
unsigned int i=0; i<nc; i++)
 
  586             libmesh_assert_not_equal_to (d_old, DofObject::invalid_id);
 
  587             libmesh_assert_not_equal_to (d_new, DofObject::invalid_id);
 
  589             old_dof_indices[i] = d_old;
 
  594     values.resize(old_dof_indices.size());
 
  599         values[i].raw_at(0) = 1;
 
  600         values[i].raw_index(0) = old_dof_indices[i];
 
  609 DynamicSparseNumberArray<Real, dof_id_type>
 
  610 OldSolutionCoefs<Real, &FEMContext::point_value>::
 
  616   LOG_SCOPE (
"eval_at_point()", 
"OldSolutionCoefs");
 
  618   if (!this->check_old_context(c, p))
 
  623   this->old_context.get_element_fe<
Real>
 
  624     (i, fe, this->old_context.get_elem_dim());
 
  628     this->old_context.build_new_fe(fe, p);
 
  631   const std::vector<std::vector<Real> > & phi = fe_new->
get_phi();
 
  632   const std::vector<dof_id_type> & dof_indices =
 
  633     this->old_context.get_dof_indices(i);
 
  635   const std::size_t n_dofs = phi.size();
 
  636   libmesh_assert_equal_to(n_dofs, dof_indices.size());
 
  638   DynamicSparseNumberArray<Real, dof_id_type> returnval;
 
  639   returnval.resize(n_dofs);
 
  643       returnval.raw_at(j) = phi[j][0];
 
  644       returnval.raw_index(j) = dof_indices[j];
 
  661   LOG_SCOPE (
"eval_at_point()", 
"OldSolutionCoefs");
 
  663   if (!this->check_old_context(c, p))
 
  668   this->old_context.get_element_fe<
Real>
 
  669     (i, fe, this->old_context.get_elem_dim());
 
  673     this->old_context.build_new_fe(fe, p);
 
  676   const std::vector<std::vector<RealGradient> > & dphi = fe_new->
get_dphi();
 
  677   const std::vector<dof_id_type> & dof_indices =
 
  678     this->old_context.get_dof_indices(i);
 
  680   const std::size_t n_dofs = dphi.size();
 
  681   libmesh_assert_equal_to(n_dofs, dof_indices.size());
 
  685   for (
unsigned int d = 0; d != LIBMESH_DIM; ++d)
 
  686     returnval(d).resize(n_dofs);
 
  689     for (
int d = 0; d != LIBMESH_DIM; ++d)
 
  691         returnval(d).raw_at(j) = dphi[j][0](d);
 
  692         returnval(d).raw_index(j) = dof_indices[j];
 
  701 DynamicSparseNumberArray<Real, dof_id_type>
 
  707              bool extra_hanging_dofs,
 
  710   LOG_SCOPE (
"Real eval_at_node()", 
"OldSolutionCoefs");
 
  725       (!extra_hanging_dofs ||
 
  726        flag == Elem::JUST_COARSENED ||
 
  727        flag == Elem::DO_NOTHING) &&
 
  731       DynamicSparseNumberArray<Real, dof_id_type> returnval;
 
  735       returnval.raw_at(0) = 1;
 
  736       returnval.raw_index(0) = old_id;
 
  740   return this->eval_at_point(c, i, n, 0);
 
  751              unsigned int elem_dim,
 
  753              bool extra_hanging_dofs,
 
  756   LOG_SCOPE (
"RealGradient eval_at_node()", 
"OldSolutionCoefs");
 
  771       (!extra_hanging_dofs ||
 
  772        flag == Elem::JUST_COARSENED ||
 
  773        flag == Elem::DO_NOTHING) &&
 
  778       for (
unsigned int d = 0; d != elem_dim; ++d)
 
  784           g(d).raw_index(0) = old_id;
 
  789   return this->eval_at_point(c, i, n, 0);
 
  802 template <
typename ValIn, 
typename ValOut>
 
  810     target_matrix(target_mat) {}
 
  813               const DynamicSparseNumberArray<ValIn, dof_id_type> & val)
 
  819       const std::size_t dnsa_size = val.size();
 
  820       for (
unsigned int j = 0; j != dnsa_size; ++j)
 
  823           const ValIn dof_val = val.raw_at(j);
 
  824           target_matrix.
set(
id, dof_j, dof_val);
 
  830   void insert(
const std::vector<dof_id_type> & dof_indices,
 
  831               const std::vector<DynamicSparseNumberArray<ValIn, dof_id_type> > & Ue)
 
  837     unsigned int size = Ue.size();
 
  839     libmesh_assert_equal_to(size, dof_indices.size());
 
  845       for (
unsigned int i = 0; i != size; ++i)
 
  848           if ((dof_i >= begin_dof) && (dof_i < end_dof))
 
  850               const DynamicSparseNumberArray<ValIn,dof_id_type> & dnsa = Ue[i];
 
  851               const std::size_t dnsa_size = dnsa.size();
 
  852               for (
unsigned int j = 0; j != dnsa_size; ++j)
 
  855                   const ValIn dof_val = dnsa.raw_at(j);
 
  856                   target_matrix.
set(dof_i, dof_j, dof_val);
 
  872   LOG_SCOPE (
"projection_matrix()", 
"System");
 
  874   const unsigned int n_variables = this->
n_vars();
 
  879         (this->get_mesh().active_local_elements_begin(),
 
  880          this->get_mesh().active_local_elements_end());
 
  882       std::vector<unsigned int> vars(n_variables);
 
  891                          OldSolutionGradientCoefs,
 
  892                          DynamicSparseNumberArray<Real,dof_id_type>,
 
  895       OldSolutionValueCoefs    f(*
this);
 
  896       OldSolutionGradientCoefs g(*
this);
 
  899       ProjMatFiller mat_filler(*
this, f, &g, setter, vars);
 
  900       mat_filler.project(active_local_elem_range);
 
  905       if (this->processor_id() == (this->n_processors()-1))
 
  907           const DofMap & dof_map = this->get_dof_map();
 
  909             if (this->variable(var).type().family == 
SCALAR)
 
  912                 std::vector<dof_id_type> new_SCALAR_indices, old_SCALAR_indices;
 
  915                 const unsigned int new_n_dofs =
 
  916                   cast_int<unsigned int>(new_SCALAR_indices.size());
 
  918                 for (
unsigned int i=0; i<new_n_dofs; i++)
 
  920                     proj_mat.
set( new_SCALAR_indices[i],
 
  921                                   old_SCALAR_indices[i], 1);
 
  927 #endif // LIBMESH_HAVE_METAPHYSICL 
  928 #endif // LIBMESH_ENABLE_AMR 
  936 void System::project_solution (ValueFunctionPointer 
fptr,
 
  937                                GradientFunctionPointer 
gptr,
 
  942   this->project_solution(&f, &g);
 
  953   this->project_vector(*solution, f, g);
 
  955   solution->localize(*current_local_solution, _dof_map->get_send_list());
 
  966   this->project_vector(*solution, f, g);
 
  968   solution->localize(*current_local_solution, _dof_map->get_send_list());
 
  976 void System::project_vector (ValueFunctionPointer 
fptr,
 
  977                              GradientFunctionPointer 
gptr,
 
  980                              int is_adjoint)
 const 
  984   this->project_vector(new_vector, &f, &g, is_adjoint);
 
  994                              int is_adjoint)
 const 
  996   LOG_SCOPE (
"project_vector(FunctionBase)", 
"System");
 
 1006       this->project_vector(new_vector, &f_fem, &g_fem, is_adjoint);
 
 1009     this->project_vector(new_vector, &f_fem, 
nullptr, is_adjoint);
 
 1020                              int is_adjoint)
 const 
 1022   LOG_SCOPE (
"project_fem_vector()", 
"System");
 
 1027     (this->get_mesh().active_local_elements_begin(),
 
 1028      this->get_mesh().active_local_elements_end() );
 
 1032   const unsigned int n_variables = this->
n_vars();
 
 1034   std::vector<unsigned int> vars(n_variables);
 
 1048       FEMProjector projector(*
this, fw, &gw, setter, vars);
 
 1049       projector.project(active_local_range);
 
 1053       FEMProjector projector(*
this, fw, 
nullptr, setter, vars);
 
 1054       projector.project(active_local_range);
 
 1060   if (this->processor_id() == (this->n_processors()-1))
 
 1065       const DofMap & dof_map = this->get_dof_map();
 
 1067         if (this->variable(var).type().family == 
SCALAR)
 
 1072             context.
pre_fe_reinit(*
this, *(this->get_mesh().active_local_elements_begin()));
 
 1074             std::vector<dof_id_type> SCALAR_indices;
 
 1076             const unsigned int n_SCALAR_dofs =
 
 1077               cast_int<unsigned int>(SCALAR_indices.size());
 
 1079             for (
unsigned int i=0; i<n_SCALAR_dofs; i++)
 
 1081                 const dof_id_type global_index = SCALAR_indices[i];
 
 1082                 const unsigned int component_index =
 
 1083                   this->variable_scalar_number(var,i);
 
 1085                 new_vector.
set(global_index, f->
component(context, component_index, 
Point(), this->time));
 
 1092 #ifdef LIBMESH_ENABLE_CONSTRAINTS 
 1093   if (is_adjoint == -1)
 
 1094     this->get_dof_map().enforce_constraints_exactly(*
this, &new_vector);
 
 1095   else if (is_adjoint >= 0)
 
 1096     this->get_dof_map().enforce_adjoint_constraints_exactly(new_vector,
 
 1109 void System::boundary_project_solution (
const std::set<boundary_id_type> & b,
 
 1110                                         const std::vector<unsigned int> & variables,
 
 1111                                         ValueFunctionPointer 
fptr,
 
 1112                                         GradientFunctionPointer 
gptr,
 
 1117   this->boundary_project_solution(b, variables, &f, &g);
 
 1126 void System::boundary_project_solution (
const std::set<boundary_id_type> & b,
 
 1127                                         const std::vector<unsigned int> & variables,
 
 1131   this->boundary_project_vector(b, variables, *solution, f, g);
 
 1133   solution->localize(*current_local_solution);
 
 1144 void System::boundary_project_vector (
const std::set<boundary_id_type> & b,
 
 1145                                       const std::vector<unsigned int> & variables,
 
 1146                                       ValueFunctionPointer 
fptr,
 
 1147                                       GradientFunctionPointer 
gptr,
 
 1150                                       int is_adjoint)
 const 
 1154   this->boundary_project_vector(b, variables, new_vector, &f, &g,
 
 1162 void System::boundary_project_vector (
const std::set<boundary_id_type> & b,
 
 1163                                       const std::vector<unsigned int> & variables,
 
 1167                                       int is_adjoint)
 const 
 1169   LOG_SCOPE (
"boundary_project_vector()", 
"System");
 
 1173                      this->get_mesh().active_local_elements_end() ),
 
 1175                              this->get_equation_systems().parameters,
 
 1184 #ifdef LIBMESH_ENABLE_CONSTRAINTS 
 1185   if (is_adjoint == -1)
 
 1186     this->get_dof_map().enforce_constraints_exactly(*
this, &new_vector);
 
 1187   else if (is_adjoint >= 0)
 
 1188     this->get_dof_map().enforce_adjoint_constraints_exactly(new_vector,
 
 1197 #ifdef LIBMESH_ENABLE_AMR 
 1198 void BuildProjectionList::unique()
 
 1202   std::sort(this->send_list.begin(),
 
 1203             this->send_list.end());
 
 1206   std::vector<dof_id_type>::iterator new_end =
 
 1207     std::unique (this->send_list.begin(),
 
 1208                  this->send_list.end());
 
 1212   std::vector<dof_id_type>
 
 1213     (this->send_list.begin(), new_end).
swap (this->send_list);
 
 1221   const DofMap & dof_map = system.get_dof_map();
 
 1228   std::vector<dof_id_type> di;
 
 1231   for (
const auto & elem : range)
 
 1246       if (!elem->old_dof_object &&
 
 1247           elem->refinement_flag() != Elem::JUST_REFINED &&
 
 1248           elem->refinement_flag() != Elem::JUST_COARSENED)
 
 1253       if (elem->refinement_flag() == Elem::JUST_REFINED)
 
 1263           for (
auto & node : elem->node_ref_range())
 
 1269                   const unsigned int sysnum = system.number();
 
 1270                   const unsigned int nvg = old_dofs->
n_var_groups(sysnum);
 
 1272                   for (
unsigned int vg=0; vg != nvg; ++vg)
 
 1274                       const unsigned int nvig =
 
 1275                         old_dofs->
n_vars(sysnum, vg);
 
 1276                       for (
unsigned int vig=0; vig != nvig; ++vig)
 
 1278                           const unsigned int n_comp =
 
 1280                           for (
unsigned int c=0; c != n_comp; ++c)
 
 1282                               (old_dofs->
dof_number(sysnum, vg, vig, c, n_comp));
 
 1288           std::sort(di.begin(), di.end());
 
 1289           std::vector<dof_id_type>::iterator new_end =
 
 1290             std::unique(di.begin(), di.end());
 
 1291           std::vector<dof_id_type>(di.begin(), new_end).
swap(di);
 
 1293       else if (elem->refinement_flag() == Elem::JUST_COARSENED)
 
 1295           std::vector<dof_id_type> di_child;
 
 1297           for (
auto & child : elem->child_ref_range())
 
 1300               di.insert(di.end(), di_child.begin(), di_child.end());
 
 1306       for (
auto di_i : di)
 
 1307         if (di_i < first_old_dof || di_i >= end_old_dof)
 
 1308           this->send_list.push_back(di_i);
 
 1317   this->send_list.insert(this->send_list.end(),
 
 1321 #endif // LIBMESH_ENABLE_AMR 
 1338   const unsigned int dim = system.get_mesh().mesh_dimension();
 
 1341   const DofMap & dof_map = system.get_dof_map();
 
 1345     system.get_mesh().get_boundary_info();
 
 1360       const unsigned int var = variables[v];
 
 1369       const unsigned int var_component =
 
 1370         system.variable_scalar_number(var, 0);
 
 1373       std::unique_ptr<FEBase> fe (FEBase::build(
dim, fe_type));
 
 1381       const std::vector<std::vector<Real>> & phi = fe->get_phi();
 
 1385       const std::vector<std::vector<RealGradient>> * dphi = 
nullptr;
 
 1394           const std::vector<std::vector<RealGradient>> &
 
 1395             ref_dphi = fe->get_dphi();
 
 1400       const std::vector<Real> & JxW =
 
 1404       const std::vector<Point> & xyz_values =
 
 1408       std::vector<dof_id_type> dof_indices;
 
 1410       std::vector<unsigned int> side_dofs;
 
 1413       std::vector<boundary_id_type> bc_ids;
 
 1416       for (
const auto & elem : range)
 
 1423           const unsigned short n_nodes = elem->n_nodes();
 
 1424           const unsigned short n_edges = elem->n_edges();
 
 1425           const unsigned short n_sides = elem->n_sides();
 
 1429           std::vector<bool> is_boundary_node(
n_nodes, 
false),
 
 1430             is_boundary_edge(n_edges, 
false),
 
 1431             is_boundary_side(n_sides, 
false);
 
 1434           std::vector<bool> is_boundary_nodeset(
n_nodes, 
false);
 
 1436           for (
unsigned char s=0; s != n_sides; ++s)
 
 1440               bool do_this_side = 
false;
 
 1441               for (
const auto & bc_id : bc_ids)
 
 1444                     do_this_side = 
true;
 
 1450               is_boundary_side[s] = 
true;
 
 1453               for (
unsigned int n=0; n != 
n_nodes; ++n)
 
 1454                 if (elem->is_node_on_side(n,s))
 
 1455                   is_boundary_node[n] = 
true;
 
 1456               for (
unsigned int e=0; e != n_edges; ++e)
 
 1457                 if (elem->is_edge_on_side(e,s))
 
 1458                   is_boundary_edge[e] = 
true;
 
 1463             for (
unsigned int n=0; n != 
n_nodes; ++n)
 
 1465                 boundary_info.
boundary_ids (elem->node_ptr(n), bc_ids);
 
 1467                 for (
const auto & bc_id : bc_ids)
 
 1470                       is_boundary_node[n] = 
true;
 
 1471                       is_boundary_nodeset[n] = 
true;
 
 1477             for (
unsigned short e=0; e != n_edges; ++e)
 
 1481                 for (
const auto & bc_id : bc_ids)
 
 1483                     is_boundary_edge[e] = 
true;
 
 1491           const unsigned int n_dofs =
 
 1492             cast_int<unsigned int>(dof_indices.size());
 
 1495           std::vector<char> dof_is_fixed(n_dofs, 
false); 
 
 1496           std::vector<int> free_dof(n_dofs, 0);
 
 1499           const ElemType elem_type = elem->type();
 
 1511           unsigned int current_dof = 0;
 
 1512           for (
unsigned short n = 0; n != 
n_nodes; ++n)
 
 1516               const unsigned int nc =
 
 1517                 FEInterface::n_dofs_at_node (
dim, fe_type, elem_type,
 
 1519               if ((!elem->is_vertex(n) || !is_boundary_node[n]) &&
 
 1520                   !is_boundary_nodeset[n])
 
 1527                   libmesh_assert_equal_to (nc, 0);
 
 1533                   libmesh_assert_equal_to (nc, 1);
 
 1534                   Ue(current_dof) = f->component(var_component,
 
 1537                   dof_is_fixed[current_dof] = 
true;
 
 1543                   Ue(current_dof) = f->component(var_component,
 
 1546                   dof_is_fixed[current_dof] = 
true;
 
 1548                   Gradient grad = g->component(var_component,
 
 1552                   Ue(current_dof) = grad(0);
 
 1553                   dof_is_fixed[current_dof] = 
true;
 
 1559                       Point nxminus = elem->point(n),
 
 1560                         nxplus = elem->point(n);
 
 1563                       Gradient gxminus = g->component(var_component,
 
 1566                       Gradient gxplus = g->component(var_component,
 
 1570                       Ue(current_dof) = grad(1);
 
 1571                       dof_is_fixed[current_dof] = 
true;
 
 1574                       Ue(current_dof) = (gxplus(1) - gxminus(1))
 
 1576                       dof_is_fixed[current_dof] = 
true;
 
 1583                           Ue(current_dof) = grad(2);
 
 1584                           dof_is_fixed[current_dof] = 
true;
 
 1587                           Ue(current_dof) = (gxplus(2) - gxminus(2))
 
 1589                           dof_is_fixed[current_dof] = 
true;
 
 1592                           Point nyminus = elem->point(n),
 
 1593                             nyplus = elem->point(n);
 
 1596                           Gradient gyminus = g->component(var_component,
 
 1599                           Gradient gyplus = g->component(var_component,
 
 1603                           Ue(current_dof) = (gyplus(2) - gyminus(2))
 
 1605                           dof_is_fixed[current_dof] = 
true;
 
 1608                           Point nxmym = elem->point(n),
 
 1609                             nxmyp = elem->point(n),
 
 1610                             nxpym = elem->point(n),
 
 1611                             nxpyp = elem->point(n);
 
 1620                           Gradient gxmym = g->component(var_component,
 
 1623                           Gradient gxmyp = g->component(var_component,
 
 1626                           Gradient gxpym = g->component(var_component,
 
 1629                           Gradient gxpyp = g->component(var_component,
 
 1632                           Number gxzplus = (gxpyp(2) - gxmyp(2))
 
 1634                           Number gxzminus = (gxpym(2) - gxmym(2))
 
 1637                           Ue(current_dof) = (gxzplus - gxzminus)
 
 1639                           dof_is_fixed[current_dof] = 
true;
 
 1642 #endif // LIBMESH_DIM > 2 
 1644 #endif // LIBMESH_DIM > 1 
 1649               else if (cont == 
C_ONE)
 
 1651                   libmesh_assert_equal_to (nc, 1 + 
dim);
 
 1652                   Ue(current_dof) = f->component(var_component,
 
 1655                   dof_is_fixed[current_dof] = 
true;
 
 1657                   Gradient grad = g->component(var_component,
 
 1660                   for (
unsigned int i=0; i!= 
dim; ++i)
 
 1662                       Ue(current_dof) = grad(i);
 
 1663                       dof_is_fixed[current_dof] = 
true;
 
 1668                 libmesh_error_msg(
"Unknown continuity " << cont);
 
 1673             for (
unsigned short e = 0; e != n_edges; ++e)
 
 1675                 if (!is_boundary_edge[e])
 
 1678                 FEInterface::dofs_on_edge(elem, 
dim, fe_type, e,
 
 1681                 const unsigned int n_side_dofs =
 
 1682                   cast_int<unsigned int>(side_dofs.size());
 
 1686                 unsigned int free_dofs = 0;
 
 1688                   if (!dof_is_fixed[side_dofs[i]])
 
 1689                     free_dof[free_dofs++] = i;
 
 1701                 fe->attach_quadrature_rule (qedgerule.get());
 
 1702                 fe->edge_reinit (elem, e);
 
 1703                 const unsigned int n_qp = qedgerule->n_points();
 
 1706                 for (
unsigned int qp=0; qp<n_qp; qp++)
 
 1709                     Number fineval = f->component(var_component,
 
 1715                       finegrad = g->component(var_component,
 
 1720                     for (
unsigned int sidei=0, freei=0;
 
 1721                          sidei != n_side_dofs; ++sidei)
 
 1723                         unsigned int i = side_dofs[sidei];
 
 1725                         if (dof_is_fixed[i])
 
 1727                         for (
unsigned int sidej=0, freej=0;
 
 1728                              sidej != n_side_dofs; ++sidej)
 
 1730                             unsigned int j = side_dofs[sidej];
 
 1731                             if (dof_is_fixed[j])
 
 1732                               Fe(freei) -= phi[i][qp] * phi[j][qp] *
 
 1735                               Ke(freei,freej) += phi[i][qp] *
 
 1736                                 phi[j][qp] * JxW[qp];
 
 1739                                 if (dof_is_fixed[j])
 
 1740                                   Fe(freei) -= ((*dphi)[i][qp] *
 
 1744                                   Ke(freei,freej) += ((*dphi)[i][qp] *
 
 1748                             if (!dof_is_fixed[j])
 
 1751                         Fe(freei) += phi[i][qp] * fineval * JxW[qp];
 
 1753                           Fe(freei) += (finegrad * (*dphi)[i][qp]) *
 
 1762                 for (
unsigned int i=0; i != free_dofs; ++i)
 
 1764                     Number & ui = Ue(side_dofs[free_dof[i]]);
 
 1768                     dof_is_fixed[side_dofs[free_dof[i]]] = 
true;
 
 1774             for (
unsigned short s = 0; s != n_sides; ++s)
 
 1776                 if (!is_boundary_side[s])
 
 1779                 FEInterface::dofs_on_side(elem, 
dim, fe_type, s,
 
 1784                 unsigned int free_dofs = 0;
 
 1786                   if (!dof_is_fixed[side_dofs[i]])
 
 1787                     free_dof[free_dofs++] = i;
 
 1799                 fe->attach_quadrature_rule (qsiderule.get());
 
 1800                 fe->reinit (elem, s);
 
 1801                 const unsigned int n_qp = qsiderule->n_points();
 
 1803                 const unsigned int n_side_dofs =
 
 1804                   cast_int<unsigned int>(side_dofs.size());
 
 1807                 for (
unsigned int qp=0; qp<n_qp; qp++)
 
 1810                     Number fineval = f->component(var_component,
 
 1816                       finegrad = g->component(var_component,
 
 1821                     for (
unsigned int sidei=0, freei=0;
 
 1822                          sidei != n_side_dofs; ++sidei)
 
 1824                         unsigned int i = side_dofs[sidei];
 
 1826                         if (dof_is_fixed[i])
 
 1828                         for (
unsigned int sidej=0, freej=0;
 
 1829                              sidej != n_side_dofs; ++sidej)
 
 1831                             unsigned int j = side_dofs[sidej];
 
 1832                             if (dof_is_fixed[j])
 
 1833                               Fe(freei) -= phi[i][qp] * phi[j][qp] *
 
 1836                               Ke(freei,freej) += phi[i][qp] *
 
 1837                                 phi[j][qp] * JxW[qp];
 
 1840                                 if (dof_is_fixed[j])
 
 1841                                   Fe(freei) -= ((*dphi)[i][qp] *
 
 1845                                   Ke(freei,freej) += ((*dphi)[i][qp] *
 
 1849                             if (!dof_is_fixed[j])
 
 1852                         Fe(freei) += (fineval * phi[i][qp]) * JxW[qp];
 
 1854                           Fe(freei) += (finegrad * (*dphi)[i][qp]) *
 
 1863                 for (
unsigned int i=0; i != free_dofs; ++i)
 
 1865                     Number & ui = Ue(side_dofs[free_dof[i]]);
 
 1869                     dof_is_fixed[side_dofs[free_dof[i]]] = 
true;
 
 1874             first = new_vector.first_local_index(),
 
 1875             last  = new_vector.last_local_index();
 
 1881             for (
unsigned int i = 0; i < n_dofs; i++)
 
 1882               if (dof_is_fixed[i] &&
 
 1883                   (dof_indices[i] >= first) &&
 
 1884                   (dof_indices[i] <  last))
 
 1886                   new_vector.set(dof_indices[i], Ue(i));