27 #ifdef LIBMESH_HAVE_METAPHYSICL 33 #include "metaphysicl/dynamicsparsenumberarray_decl.h" 39 using MetaPhysicL::DynamicSparseNumberArray;
45 template <
typename T,
typename I>
53 template <
typename T,
typename I,
typename T2>
54 struct CompareTypes<MetaPhysicL::DynamicSparseNumberArray<T,I>, T2>
57 MetaPhysicL::DynamicSparseNumberArray
85 #ifdef LIBMESH_HAVE_METAPHYSICL 89 #include "metaphysicl/dynamicsparsenumberarray.h" 96 typedef DynamicSparseNumberArray<Real, dof_id_type>
DSNAN;
115 #ifdef LIBMESH_ENABLE_AMR 139 system(other.system),
149 #endif // LIBMESH_ENABLE_AMR 160 const std::set<boundary_id_type> &
b;
163 std::unique_ptr<FunctionBase<Number>>
f;
164 std::unique_ptr<FunctionBase<Gradient>>
g;
170 const std::vector<unsigned int> & variables_in,
177 variables(variables_in),
181 parameters(parameters_in),
184 libmesh_assert(f.get());
192 variables(in.variables),
196 parameters(in.parameters),
197 new_vector(in.new_vector)
199 libmesh_assert(f.get());
213 int is_adjoint)
const 217 std::unique_ptr<NumericVector<Number>>
218 old_vector (vector.
clone());
221 this->project_vector (*old_vector, vector, is_adjoint);
233 #ifdef LIBMESH_ENABLE_AMR
238 LOG_SCOPE (
"project_vector(old,new)",
"System");
248 #ifdef LIBMESH_ENABLE_AMR 252 std::unique_ptr<NumericVector<Number>> new_vector_built;
254 std::unique_ptr<NumericVector<Number>> local_old_vector_built;
258 (this->get_mesh().active_local_elements_begin(),
259 this->get_mesh().active_local_elements_end());
266 new_vector_ptr = &new_v;
267 old_vector_ptr = &old_v;
275 BuildProjectionList projection_list(*
this);
280 projection_list.unique();
282 new_v.
init (this->n_dofs(), this->n_local_dofs(),
false,
PARALLEL);
285 new_vector_ptr = new_vector_built.
get();
286 local_old_vector = local_old_vector_built.
get();
287 new_vector_ptr->
init(this->n_dofs(),
false,
SERIAL);
289 old_v.
localize(*local_old_vector, projection_list.send_list);
290 local_old_vector->
close();
291 old_vector_ptr = local_old_vector;
296 BuildProjectionList projection_list(*
this);
301 projection_list.unique();
303 new_v.
init (this->n_dofs(), this->n_local_dofs(),
304 this->get_dof_map().get_send_list(),
false,
GHOSTED);
307 new_vector_ptr = &new_v;
308 local_old_vector = local_old_vector_built.
get();
310 projection_list.send_list,
false,
GHOSTED);
311 old_v.
localize(*local_old_vector, projection_list.send_list);
312 local_old_vector->
close();
313 old_vector_ptr = local_old_vector;
316 libmesh_error_msg(
"ERROR: Unknown old_v.type() == " << old_v.
type());
321 libmesh_assert(new_vector_ptr);
322 libmesh_assert(old_vector_ptr);
324 NumericVector<Number> & new_vector = *new_vector_ptr;
325 const NumericVector<Number> & old_vector = *old_vector_ptr;
327 const unsigned int n_variables = this->
n_vars();
331 std::vector<unsigned int> vars(n_variables);
336 GenericProjector<OldSolutionValue<Number, &FEMContext::point_value>,
337 OldSolutionValue<Gradient, &FEMContext::point_gradient>,
338 Number, VectorSetAction<Number>> FEMProjector;
340 OldSolutionValue<Number, &FEMContext::point_value> f(*
this, old_vector);
341 OldSolutionValue<Gradient, &FEMContext::point_gradient> g(*
this, old_vector);
342 VectorSetAction<Number> setter(new_vector);
345 FEMProjector(*
this, f, &g, setter, vars));
350 if (this->processor_id() == (this->n_processors()-1))
352 const DofMap & dof_map = this->get_dof_map();
353 for (
unsigned int var=0; var<this->
n_vars(); var++)
354 if (this->variable(var).type().family ==
SCALAR)
357 std::vector<dof_id_type> new_SCALAR_indices, old_SCALAR_indices;
358 dof_map.SCALAR_dof_indices (new_SCALAR_indices, var,
false);
359 dof_map.SCALAR_dof_indices (old_SCALAR_indices, var,
true);
361 new_vector.set(new_SCALAR_indices[i], old_vector(old_SCALAR_indices[i]));
376 dist_v->init(this->n_dofs(), this->n_local_dofs(),
false,
PARALLEL);
380 if (new_vector(i) != 0.0)
381 dist_v->set(i, new_vector(i));
385 dist_v->localize (new_v, this->get_dof_map().get_send_list());
396 if (new_vector(i) != 0.0)
397 new_v.
set(i, new_vector(i));
401 if (is_adjoint == -1)
402 this->get_dof_map().enforce_constraints_exactly(*
this, &new_v);
403 else if (is_adjoint >= 0)
404 this->get_dof_map().enforce_adjoint_constraints_exactly(new_v,
412 #endif // #ifdef LIBMESH_ENABLE_AMR 416 #ifdef LIBMESH_ENABLE_AMR 417 #ifdef LIBMESH_HAVE_METAPHYSICL 419 template <
typename Output>
423 typedef DynamicSparseNumberArray<Output, dof_id_type>
type;
426 template <
typename InnerOutput>
442 template <
typename Output,
443 void (
FEMContext::*point_output) (
unsigned int,
455 this->old_context.set_algebraic_type(FEMContext::OLD_DOFS_ONLY);
461 this->old_context.set_algebraic_type(FEMContext::OLD_DOFS_ONLY);
466 unsigned int elem_dim,
477 std::vector<DSNA> & values)
479 LOG_SCOPE (
"eval_old_dofs()",
"OldSolutionValue");
481 this->check_old_context(c);
483 const std::vector<dof_id_type> & old_dof_indices =
484 this->old_context.get_dof_indices(var);
486 libmesh_assert_equal_to (old_dof_indices.size(), values.size());
491 values[i].raw_at(0) = 1;
492 values[i].raw_index(0) = old_dof_indices[i];
501 DynamicSparseNumberArray<Real, dof_id_type>
502 OldSolutionCoefs<Real, &FEMContext::point_value>::
508 LOG_SCOPE (
"eval_at_point()",
"OldSolutionCoefs");
510 if (!this->check_old_context(c, p))
515 this->old_context.get_element_fe<
Real>
516 (i, fe, this->old_context.get_elem_dim());
520 this->old_context.build_new_fe(fe, p);
523 const std::vector<std::vector<Real> > & phi = fe_new->
get_phi();
524 const std::vector<dof_id_type> & dof_indices =
525 this->old_context.get_dof_indices(i);
527 const std::size_t n_dofs = phi.size();
528 libmesh_assert_equal_to(n_dofs, dof_indices.size());
530 DynamicSparseNumberArray<Real, dof_id_type> returnval;
531 returnval.resize(n_dofs);
535 returnval.raw_at(j) = phi[j][0];
536 returnval.raw_index(j) = dof_indices[j];
553 LOG_SCOPE (
"eval_at_point()",
"OldSolutionCoefs");
555 if (!this->check_old_context(c, p))
560 this->old_context.get_element_fe<
Real>
561 (i, fe, this->old_context.get_elem_dim());
565 this->old_context.build_new_fe(fe, p);
568 const std::vector<std::vector<RealGradient> > & dphi = fe_new->
get_dphi();
569 const std::vector<dof_id_type> & dof_indices =
570 this->old_context.get_dof_indices(i);
572 const std::size_t n_dofs = dphi.size();
573 libmesh_assert_equal_to(n_dofs, dof_indices.size());
577 for (
unsigned int d = 0; d != LIBMESH_DIM; ++d)
578 returnval(d).resize(n_dofs);
581 for (
int d = 0; d != LIBMESH_DIM; ++d)
583 returnval(d).raw_at(j) = dphi[j][0](d);
584 returnval(d).raw_index(j) = dof_indices[j];
593 DynamicSparseNumberArray<Real, dof_id_type>
601 LOG_SCOPE (
"Real eval_at_node()",
"OldSolutionCoefs");
614 DynamicSparseNumberArray<Real, dof_id_type> returnval;
618 returnval.raw_at(0) = 1;
619 returnval.raw_index(0) = old_id;
623 return this->eval_at_point(c, i, n, 0);
634 unsigned int elem_dim,
638 LOG_SCOPE (
"RealGradient eval_at_node()",
"OldSolutionCoefs");
652 for (
unsigned int d = 0; d != elem_dim; ++d)
658 g(d).raw_index(0) = old_id;
663 return this->eval_at_point(c, i, n, 0);
676 template <
typename ValIn,
typename ValOut>
684 target_matrix(target_mat) {}
687 unsigned int var_num,
688 const DenseVector<DynamicSparseNumberArray<ValIn, dof_id_type> > & Ue)
694 const std::vector<dof_id_type> & dof_indices =
697 unsigned int size = Ue.size();
699 libmesh_assert_equal_to(size, dof_indices.size());
705 for (
unsigned int i = 0; i != size; ++i)
708 if ((dof_i >= begin_dof) && (dof_i < end_dof))
710 const DynamicSparseNumberArray<ValIn,dof_id_type> & dnsa = Ue(i);
711 const std::size_t dnsa_size = dnsa.size();
712 for (
unsigned int j = 0; j != dnsa_size; ++j)
715 const ValIn dof_val = dnsa.raw_at(j);
716 target_matrix.
set(dof_i, dof_j, dof_val);
732 LOG_SCOPE (
"projection_matrix()",
"System");
734 const unsigned int n_variables = this->
n_vars();
739 (this->get_mesh().active_local_elements_begin(),
740 this->get_mesh().active_local_elements_end());
742 std::vector<unsigned int> vars(n_variables);
751 OldSolutionGradientCoefs,
752 DynamicSparseNumberArray<Real,dof_id_type>,
755 OldSolutionValueCoefs f(*
this);
756 OldSolutionGradientCoefs g(*
this);
760 ProjMatFiller(*
this, f, &g, setter, vars));
765 if (this->processor_id() == (this->n_processors()-1))
767 const DofMap & dof_map = this->get_dof_map();
768 for (
unsigned int var=0; var<this->
n_vars(); var++)
769 if (this->variable(var).type().family ==
SCALAR)
772 std::vector<dof_id_type> new_SCALAR_indices, old_SCALAR_indices;
775 const unsigned int new_n_dofs =
776 cast_int<unsigned int>(new_SCALAR_indices.size());
778 for (
unsigned int i=0; i<new_n_dofs; i++)
780 proj_mat.
set( new_SCALAR_indices[i],
781 old_SCALAR_indices[i], 1);
787 #endif // LIBMESH_HAVE_METAPHYSICL 788 #endif // LIBMESH_ENABLE_AMR 796 void System::project_solution (ValueFunctionPointer fptr,
797 GradientFunctionPointer gptr,
802 this->project_solution(&f, &g);
813 this->project_vector(*solution, f, g);
815 solution->localize(*current_local_solution, _dof_map->get_send_list());
826 this->project_vector(*solution, f, g);
828 solution->localize(*current_local_solution, _dof_map->get_send_list());
836 void System::project_vector (ValueFunctionPointer fptr,
837 GradientFunctionPointer gptr,
840 int is_adjoint)
const 844 this->project_vector(new_vector, &f, &g, is_adjoint);
854 int is_adjoint)
const 856 LOG_SCOPE (
"project_vector(FunctionBase)",
"System");
866 this->project_vector(new_vector, &f_fem, &g_fem, is_adjoint);
869 this->project_vector(new_vector, &f_fem,
nullptr, is_adjoint);
880 int is_adjoint)
const 882 LOG_SCOPE (
"project_fem_vector()",
"System");
887 (this->get_mesh().active_local_elements_begin(),
888 this->get_mesh().active_local_elements_end() );
892 const unsigned int n_variables = this->
n_vars();
894 std::vector<unsigned int> vars(n_variables);
910 FEMProjector(*
this, fw, &gw, setter, vars));
915 FEMProjector(*
this, fw,
nullptr, setter, vars));
920 if (this->processor_id() == (this->n_processors()-1))
925 const DofMap & dof_map = this->get_dof_map();
926 for (
unsigned int var=0; var<this->
n_vars(); var++)
927 if (this->variable(var).type().family ==
SCALAR)
932 Elem * el =
const_cast<Elem *
>(*(this->get_mesh().active_local_elements_begin()));
935 std::vector<dof_id_type> SCALAR_indices;
937 const unsigned int n_SCALAR_dofs =
938 cast_int<unsigned int>(SCALAR_indices.size());
940 for (
unsigned int i=0; i<n_SCALAR_dofs; i++)
942 const dof_id_type global_index = SCALAR_indices[i];
943 const unsigned int component_index =
944 this->variable_scalar_number(var,i);
946 new_vector.
set(global_index, f->
component(context, component_index,
Point(), this->time));
953 #ifdef LIBMESH_ENABLE_CONSTRAINTS 954 if (is_adjoint == -1)
955 this->get_dof_map().enforce_constraints_exactly(*
this, &new_vector);
956 else if (is_adjoint >= 0)
957 this->get_dof_map().enforce_adjoint_constraints_exactly(new_vector,
968 void System::boundary_project_solution (
const std::set<boundary_id_type> & b,
969 const std::vector<unsigned int> & variables,
970 ValueFunctionPointer fptr,
971 GradientFunctionPointer gptr,
976 this->boundary_project_solution(b, variables, &f, &g);
985 void System::boundary_project_solution (
const std::set<boundary_id_type> & b,
986 const std::vector<unsigned int> & variables,
990 this->boundary_project_vector(b, variables, *solution, f, g);
992 solution->localize(*current_local_solution);
1003 void System::boundary_project_vector (
const std::set<boundary_id_type> & b,
1004 const std::vector<unsigned int> & variables,
1005 ValueFunctionPointer fptr,
1006 GradientFunctionPointer gptr,
1009 int is_adjoint)
const 1013 this->boundary_project_vector(b, variables, new_vector, &f, &g,
1021 void System::boundary_project_vector (
const std::set<boundary_id_type> & b,
1022 const std::vector<unsigned int> & variables,
1026 int is_adjoint)
const 1028 LOG_SCOPE (
"boundary_project_vector()",
"System");
1032 this->get_mesh().active_local_elements_end() ),
1034 this->get_equation_systems().parameters,
1043 #ifdef LIBMESH_ENABLE_CONSTRAINTS 1044 if (is_adjoint == -1)
1045 this->get_dof_map().enforce_constraints_exactly(*
this, &new_vector);
1046 else if (is_adjoint >= 0)
1047 this->get_dof_map().enforce_adjoint_constraints_exactly(new_vector,
1054 #ifdef LIBMESH_ENABLE_AMR 1055 void BuildProjectionList::unique()
1059 std::sort(this->send_list.begin(),
1060 this->send_list.end());
1063 std::vector<dof_id_type>::iterator new_end =
1064 std::unique (this->send_list.begin(),
1065 this->send_list.end());
1069 std::vector<dof_id_type>
1070 (this->send_list.begin(), new_end).
swap (this->send_list);
1078 const DofMap & dof_map = system.get_dof_map();
1085 std::vector<dof_id_type> di;
1088 for (
const auto & elem : range)
1103 if (!elem->old_dof_object &&
1104 elem->refinement_flag() != Elem::JUST_REFINED &&
1105 elem->refinement_flag() != Elem::JUST_COARSENED)
1110 if (elem->refinement_flag() == Elem::JUST_REFINED)
1112 libmesh_assert(parent);
1120 for (
auto & node : elem->node_ref_range())
1126 const unsigned int sysnum = system.number();
1127 const unsigned int nvg = old_dofs->
n_var_groups(sysnum);
1129 for (
unsigned int vg=0; vg != nvg; ++vg)
1131 const unsigned int nvig =
1132 old_dofs->
n_vars(sysnum, vg);
1133 for (
unsigned int vig=0; vig != nvig; ++vig)
1135 const unsigned int n_comp =
1137 for (
unsigned int c=0; c != n_comp; ++c)
1139 (old_dofs->
dof_number(sysnum, vg, vig, c, n_comp));
1145 std::sort(di.begin(), di.end());
1146 std::vector<dof_id_type>::iterator new_end =
1147 std::unique(di.begin(), di.end());
1148 std::vector<dof_id_type>(di.begin(), new_end).
swap(di);
1150 else if (elem->refinement_flag() == Elem::JUST_COARSENED)
1152 std::vector<dof_id_type> di_child;
1154 for (
auto & child : elem->child_ref_range())
1157 di.insert(di.end(), di_child.begin(), di_child.end());
1163 for (std::size_t i=0; i != di.size(); ++i)
1164 if (di[i] < first_old_dof || di[i] >= end_old_dof)
1165 this->send_list.push_back(di[i]);
1174 this->send_list.insert(this->send_list.end(),
1178 #endif // LIBMESH_ENABLE_AMR 1185 libmesh_assert(f.get());
1195 const unsigned int dim = system.get_mesh().mesh_dimension();
1198 const DofMap & dof_map = system.get_dof_map();
1202 system.get_mesh().get_boundary_info();
1215 for (std::size_t v=0; v!=variables.size(); v++)
1217 const unsigned int var = variables[v];
1226 const unsigned int var_component =
1227 system.variable_scalar_number(var, 0);
1230 std::unique_ptr<FEBase> fe (FEBase::build(dim, fe_type));
1238 const std::vector<std::vector<Real>> & phi = fe->get_phi();
1242 const std::vector<std::vector<RealGradient>> * dphi =
nullptr;
1249 libmesh_assert(g.get());
1251 const std::vector<std::vector<RealGradient>> &
1252 ref_dphi = fe->get_dphi();
1257 const std::vector<Real> & JxW =
1261 const std::vector<Point> & xyz_values =
1265 std::vector<dof_id_type> dof_indices;
1267 std::vector<unsigned int> side_dofs;
1270 std::vector<boundary_id_type> bc_ids;
1273 for (
const auto & elem : range)
1280 const unsigned short n_nodes = elem->n_nodes();
1281 const unsigned short n_edges = elem->n_edges();
1282 const unsigned short n_sides = elem->n_sides();
1286 std::vector<bool> is_boundary_node(
n_nodes,
false),
1287 is_boundary_edge(n_edges,
false),
1288 is_boundary_side(n_sides,
false);
1291 std::vector<bool> is_boundary_nodeset(
n_nodes,
false);
1293 for (
unsigned char s=0; s != n_sides; ++s)
1297 bool do_this_side =
false;
1298 for (
const auto &
bc_id : bc_ids)
1301 do_this_side =
true;
1307 is_boundary_side[s] =
true;
1310 for (
unsigned int n=0; n !=
n_nodes; ++n)
1311 if (elem->is_node_on_side(n,s))
1312 is_boundary_node[n] =
true;
1313 for (
unsigned int e=0; e != n_edges; ++e)
1314 if (elem->is_edge_on_side(e,s))
1315 is_boundary_edge[e] =
true;
1320 for (
unsigned int n=0; n !=
n_nodes; ++n)
1322 boundary_info.
boundary_ids (elem->node_ptr(n), bc_ids);
1324 for (
const auto &
bc_id : bc_ids)
1327 is_boundary_node[n] =
true;
1328 is_boundary_nodeset[n] =
true;
1334 for (
unsigned short e=0; e != n_edges; ++e)
1338 for (
const auto &
bc_id : bc_ids)
1340 is_boundary_edge[e] =
true;
1348 const unsigned int n_dofs =
1349 cast_int<unsigned int>(dof_indices.size());
1352 std::vector<char> dof_is_fixed(n_dofs,
false);
1353 std::vector<int> free_dof(n_dofs, 0);
1356 const ElemType elem_type = elem->type();
1368 unsigned int current_dof = 0;
1369 for (
unsigned short n = 0; n !=
n_nodes; ++n)
1373 const unsigned int nc =
1374 FEInterface::n_dofs_at_node (dim, fe_type, elem_type,
1376 if ((!elem->is_vertex(n) || !is_boundary_node[n]) &&
1377 !is_boundary_nodeset[n])
1384 libmesh_assert_equal_to (nc, 0);
1390 libmesh_assert_equal_to (nc, 1);
1391 Ue(current_dof) = f->component(var_component,
1394 dof_is_fixed[current_dof] =
true;
1400 Ue(current_dof) = f->component(var_component,
1403 dof_is_fixed[current_dof] =
true;
1405 Gradient grad = g->component(var_component,
1409 Ue(current_dof) = grad(0);
1410 dof_is_fixed[current_dof] =
true;
1415 Point nxminus = elem->point(n),
1416 nxplus = elem->point(n);
1419 Gradient gxminus = g->component(var_component,
1422 Gradient gxplus = g->component(var_component,
1426 Ue(current_dof) = grad(1);
1427 dof_is_fixed[current_dof] =
true;
1430 Ue(current_dof) = (gxplus(1) - gxminus(1))
1432 dof_is_fixed[current_dof] =
true;
1438 Ue(current_dof) = grad(2);
1439 dof_is_fixed[current_dof] =
true;
1442 Ue(current_dof) = (gxplus(2) - gxminus(2))
1444 dof_is_fixed[current_dof] =
true;
1447 Point nyminus = elem->point(n),
1448 nyplus = elem->point(n);
1451 Gradient gyminus = g->component(var_component,
1454 Gradient gyplus = g->component(var_component,
1458 Ue(current_dof) = (gyplus(2) - gyminus(2))
1460 dof_is_fixed[current_dof] =
true;
1463 Point nxmym = elem->point(n),
1464 nxmyp = elem->point(n),
1465 nxpym = elem->point(n),
1466 nxpyp = elem->point(n);
1475 Gradient gxmym = g->component(var_component,
1478 Gradient gxmyp = g->component(var_component,
1481 Gradient gxpym = g->component(var_component,
1484 Gradient gxpyp = g->component(var_component,
1487 Number gxzplus = (gxpyp(2) - gxmyp(2))
1489 Number gxzminus = (gxpym(2) - gxmym(2))
1492 Ue(current_dof) = (gxzplus - gxzminus)
1494 dof_is_fixed[current_dof] =
true;
1502 else if (cont ==
C_ONE)
1504 libmesh_assert_equal_to (nc, 1 + dim);
1505 Ue(current_dof) = f->component(var_component,
1508 dof_is_fixed[current_dof] =
true;
1510 Gradient grad = g->component(var_component,
1513 for (
unsigned int i=0; i!= dim; ++i)
1515 Ue(current_dof) = grad(i);
1516 dof_is_fixed[current_dof] =
true;
1521 libmesh_error_msg(
"Unknown continuity " << cont);
1526 for (
unsigned short e = 0; e != n_edges; ++e)
1528 if (!is_boundary_edge[e])
1531 FEInterface::dofs_on_edge(elem, dim, fe_type, e,
1534 const unsigned int n_side_dofs =
1535 cast_int<unsigned int>(side_dofs.size());
1539 unsigned int free_dofs = 0;
1541 if (!dof_is_fixed[side_dofs[i]])
1542 free_dof[free_dofs++] = i;
1554 fe->attach_quadrature_rule (qedgerule.get());
1555 fe->edge_reinit (elem, e);
1556 const unsigned int n_qp = qedgerule->n_points();
1559 for (
unsigned int qp=0; qp<n_qp; qp++)
1562 Number fineval = f->component(var_component,
1568 finegrad = g->component(var_component,
1573 for (
unsigned int sidei=0, freei=0;
1574 sidei != n_side_dofs; ++sidei)
1576 unsigned int i = side_dofs[sidei];
1578 if (dof_is_fixed[i])
1580 for (
unsigned int sidej=0, freej=0;
1581 sidej != n_side_dofs; ++sidej)
1583 unsigned int j = side_dofs[sidej];
1584 if (dof_is_fixed[j])
1585 Fe(freei) -= phi[i][qp] * phi[j][qp] *
1588 Ke(freei,freej) += phi[i][qp] *
1589 phi[j][qp] * JxW[qp];
1592 if (dof_is_fixed[j])
1593 Fe(freei) -= ((*dphi)[i][qp] *
1597 Ke(freei,freej) += ((*dphi)[i][qp] *
1601 if (!dof_is_fixed[j])
1604 Fe(freei) += phi[i][qp] * fineval * JxW[qp];
1606 Fe(freei) += (finegrad * (*dphi)[i][qp]) *
1615 for (
unsigned int i=0; i != free_dofs; ++i)
1617 Number & ui = Ue(side_dofs[free_dof[i]]);
1621 dof_is_fixed[side_dofs[free_dof[i]]] =
true;
1627 for (
unsigned short s = 0; s != n_sides; ++s)
1629 if (!is_boundary_side[s])
1632 FEInterface::dofs_on_side(elem, dim, fe_type, s,
1637 unsigned int free_dofs = 0;
1639 if (!dof_is_fixed[side_dofs[i]])
1640 free_dof[free_dofs++] = i;
1652 fe->attach_quadrature_rule (qsiderule.get());
1653 fe->reinit (elem, s);
1654 const unsigned int n_qp = qsiderule->n_points();
1656 const unsigned int n_side_dofs =
1657 cast_int<unsigned int>(side_dofs.size());
1660 for (
unsigned int qp=0; qp<n_qp; qp++)
1663 Number fineval = f->component(var_component,
1669 finegrad = g->component(var_component,
1674 for (
unsigned int sidei=0, freei=0;
1675 sidei != n_side_dofs; ++sidei)
1677 unsigned int i = side_dofs[sidei];
1679 if (dof_is_fixed[i])
1681 for (
unsigned int sidej=0, freej=0;
1682 sidej != n_side_dofs; ++sidej)
1684 unsigned int j = side_dofs[sidej];
1685 if (dof_is_fixed[j])
1686 Fe(freei) -= phi[i][qp] * phi[j][qp] *
1689 Ke(freei,freej) += phi[i][qp] *
1690 phi[j][qp] * JxW[qp];
1693 if (dof_is_fixed[j])
1694 Fe(freei) -= ((*dphi)[i][qp] *
1698 Ke(freei,freej) += ((*dphi)[i][qp] *
1702 if (!dof_is_fixed[j])
1705 Fe(freei) += (fineval * phi[i][qp]) * JxW[qp];
1707 Fe(freei) += (finegrad * (*dphi)[i][qp]) *
1716 for (
unsigned int i=0; i != free_dofs; ++i)
1718 Number & ui = Ue(side_dofs[free_dof[i]]);
1722 dof_is_fixed[side_dofs[free_dof[i]]] =
true;
1727 first = new_vector.first_local_index(),
1728 last = new_vector.last_local_index();
1734 for (
unsigned int i = 0; i < n_dofs; i++)
1735 if (dof_is_fixed[i] &&
1736 (dof_indices[i] >= first) &&
1737 (dof_indices[i] < last))
1739 new_vector.set(dof_indices[i], Ue(i));
Manages the family, order, etc. parameters for a given FE.
void eval_old_dofs(const FEMContext &c, unsigned int var, std::vector< DSNA > &values)
Wrap a libMesh-style function pointer into a FunctionBase object.
dof_id_type dof_number(const unsigned int s, const unsigned int var, const unsigned int comp) const
const Elem * parent() const
A geometric point in (x,y,z) space associated with a DOF.
virtual void pre_fe_reinit(const System &, const Elem *e)
unsigned int n_comp(const unsigned int s, const unsigned int var) const
BuildProjectionList(const System &system_in)
SparseMatrix< ValOut > & target_matrix
unsigned int n_var_groups(const unsigned int s) const
virtual void get(const std::vector< numeric_index_type > &index, T *values) const
void dof_indices(const Elem *const elem, std::vector< dof_id_type > &di) const
void parallel_for(const Range &range, const Body &body)
std::unique_ptr< FunctionBase< Number > > f
BoundaryProjectSolution(const std::set< boundary_id_type > &b_in, const std::vector< unsigned int > &variables_in, const System &system_in, FunctionBase< Number > *f_in, FunctionBase< Gradient > *g_in, const Parameters ¶meters_in, NumericVector< Number > &new_v_in)
virtual numeric_index_type size() const =0
void resize(const unsigned int n)
DynamicSparseNumberArray< Output, dof_id_type > type
The base class for all geometric element types.
IntRange< std::size_t > index_range(const std::vector< T > &vec)
virtual std::unique_ptr< NumericVector< T > > clone() const =0
static const Real TOLERANCE
Utility class for defining generic ranges for threading.
const unsigned int n_vars
const std::vector< unsigned int > & variables
virtual void init(const numeric_index_type n, const numeric_index_type n_local, const bool fast=false, const ParallelType ptype=AUTOMATIC)=0
void old_dof_indices(const Elem *const elem, std::vector< dof_id_type > &di, const unsigned int vn=libMesh::invalid_uint) const
void iota(ForwardIter first, ForwardIter last, T value)
void SCALAR_dof_indices(std::vector< dof_id_type > &di, const unsigned int vn, const bool old_dofs=false) const
virtual numeric_index_type row_stop() const =0
virtual void set(const numeric_index_type i, const numeric_index_type j, const T value)=0
std::vector< boundary_id_type > boundary_ids(const Node *node) const
BuildProjectionList(BuildProjectionList &other, Threads::split)
const Parameters & parameters
Manages the degrees of freedom (DOFs) in a simulation.
std::unique_ptr< QBase > default_quadrature_rule(const unsigned int dim, const int extraorder=0) const
const std::vector< std::vector< OutputGradient > > & get_dphi() const
const dof_id_type n_nodes
A variable which is solved for in a System of equations.
const Variable & variable(const unsigned int c) const
dof_id_type numeric_index_type
OldSolutionCoefs(const libMesh::System &sys_in)
NumericVector< Number > & new_vector
unsigned int n_vars(const unsigned int s, const unsigned int vg) const
VectorValue< DynamicSparseNumberArray< InnerOutput, dof_id_type > > type
Manages consistently variables, degrees of freedom, and coefficient vectors.
void _cholesky_back_substitute(const DenseVector< T2 > &b, DenseVector< T2 > &x) const
Used by the Mesh to keep track of boundary nodes and elements.
dof_id_type end_old_dof(const processor_id_type proc) const
const std::vector< dof_id_type > & get_dof_indices() const
bool active_on_subdomain(subdomain_id_type sid) const
BoundaryProjectSolution(const BoundaryProjectSolution &in)
static std::unique_ptr< NumericVector< Number > > build(const Parallel::Communicator &comm, const SolverPackage solver_package=libMesh::default_solver_package())
OldSolutionCoefs(const OldSolutionCoefs &in)
ParallelType type() const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
DofObject * old_dof_object
DynamicSparseNumberArray< Real, dof_id_type > DSNAN
std::vector< dof_id_type > send_list
virtual numeric_index_type row_start() const =0
virtual void zero() override
void swap(Iterator &lhs, Iterator &rhs)
virtual numeric_index_type local_size() const =0
virtual Output component(const FEMContext &, unsigned int i, const Point &p, Real time=0.)
void resize(const unsigned int new_m, const unsigned int new_n)
void parallel_reduce(const Range &range, Body &body)
MatrixFillAction(SparseMatrix< ValOut > &target_mat)
virtual void set(const numeric_index_type i, const T value)=0
void cholesky_solve(const DenseVector< T2 > &b, DenseVector< T2 > &x)
void insert(const FEMContext &c, unsigned int var_num, const DenseVector< DynamicSparseNumberArray< ValIn, dof_id_type > > &Ue)
const std::set< boundary_id_type > & b
std::unique_ptr< FunctionBase< Gradient > > g
unsigned int n_comp_group(const unsigned int s, const unsigned int vg) const
A matrix object used for finite element assembly and numerics.
A geometric point in (x,y,z) space.
DSNAOutput< Output >::type DSNA
const Elem & get(const ElemType type_in)
dof_id_type first_old_dof(const processor_id_type proc) const
virtual void zero() override
const std::vector< std::vector< OutputShape > > & get_phi() const
const FEType & type() const
std::vector< boundary_id_type > edge_boundary_ids(const Elem *const elem, const unsigned short int edge) const
virtual void localize(std::vector< T > &v_local) const =0