libMesh
Public Member Functions | Public Attributes | Protected Member Functions | List of all members
libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::ProjectEdges Struct Reference

#include <generic_projector.h>

Inheritance diagram for libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::ProjectEdges:
[legend]

Public Member Functions

 ProjectEdges (GenericProjector &p)
 
 ProjectEdges (ProjectEdges &p_e, Threads::split)
 
void operator() (const node_range &range)
 
void insert_id (dof_id_type id, const FValue &val, processor_id_type pid)
 
void insert_ids (const std::vector< dof_id_type > &ids, const std::vector< FValue > &vals, processor_id_type pid)
 
void construct_projection (const std::vector< dof_id_type > &dof_indices_var, const std::vector< unsigned int > &involved_dofs, unsigned int var_component, const Node *node, const FEBase &fe)
 
void insert_id (dof_id_type id, const FValue &val, processor_id_type pid)
 
void insert_ids (const std::vector< dof_id_type > &ids, const std::vector< FValue > &vals, processor_id_type pid)
 
void find_dofs_to_send (const Node &node, const Elem &elem, unsigned short node_num, const var_set &vars)
 
void join (const SubFunctor &other)
 

Public Attributes

ProjectionAction action
 
FFunctor f
 
std::unique_ptr< GFunctor > g
 
const Systemsystem
 
FEMContext context
 
std::vector< FEContinuityconts
 
GenericProjectorprojector
 
std::unordered_map< dof_id_type, std::pair< FValue, processor_id_type > > new_ids_to_push
 
std::unordered_map< dof_id_type, FValue > new_ids_to_save
 

Protected Member Functions

void construct_projection (const std::vector< dof_id_type > &dof_indices_var, const std::vector< unsigned int > &involved_dofs, unsigned int var_component, const Node *node, const FEBase &fe)
 

Detailed Description

template<typename FFunctor, typename GFunctor, typename FValue, typename ProjectionAction>
struct libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::ProjectEdges

Definition at line 287 of file generic_projector.h.

Constructor & Destructor Documentation

◆ ProjectEdges() [1/2]

template<typename FFunctor , typename GFunctor , typename FValue , typename ProjectionAction >
libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::ProjectEdges::ProjectEdges ( GenericProjector p)
inline

Definition at line 288 of file generic_projector.h.

288 : SubProjector(p) {}

◆ ProjectEdges() [2/2]

template<typename FFunctor , typename GFunctor , typename FValue , typename ProjectionAction >
libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::ProjectEdges::ProjectEdges ( ProjectEdges p_e,
Threads::split   
)
inline

Definition at line 290 of file generic_projector.h.

290 : SubProjector(p_e.projector) {}

Member Function Documentation

◆ construct_projection() [1/2]

template<typename FFunctor , typename GFunctor , typename FValue , typename ProjectionAction >
void libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::SubProjector::construct_projection ( const std::vector< dof_id_type > &  dof_indices_var,
const std::vector< unsigned int > &  involved_dofs,
unsigned int  var_component,
const Node node,
const FEBase fe 
)
protectedinherited

Definition at line 2501 of file generic_projector.h.

2506 {
2507  const std::vector<Real> & JxW = fe.get_JxW();
2508  const std::vector<std::vector<Real>> & phi = fe.get_phi();
2509  const std::vector<std::vector<RealGradient>> * dphi = nullptr;
2510  const std::vector<Point> & xyz_values = fe.get_xyz();
2511  const FEContinuity cont = fe.get_continuity();
2512  const std::unordered_map<dof_id_type, FValue> & ids_to_save =
2513  this->projector.ids_to_save;
2514 
2515  if (cont == C_ONE)
2516  dphi = &(fe.get_dphi());
2517 
2518  const unsigned int n_involved_dofs =
2519  cast_int<unsigned int>(involved_dofs.size());
2520 
2521  std::vector<dof_id_type> free_dof_ids;
2522  DenseVector<FValue> Uinvolved(n_involved_dofs);
2523  std::vector<char> dof_is_fixed(n_involved_dofs, false); // bools
2524 
2525  for (auto i : IntRange<unsigned int>(0, n_involved_dofs))
2526  {
2527  const dof_id_type id = dof_indices_var[involved_dofs[i]];
2528  auto iter = ids_to_save.find(id);
2529  if (iter == ids_to_save.end())
2530  free_dof_ids.push_back(id);
2531  else
2532  {
2533  dof_is_fixed[i] = true;
2534  Uinvolved(i) = iter->second;
2535  }
2536  }
2537 
2538  const unsigned int free_dofs = free_dof_ids.size();
2539 
2540  // There may be nothing to project
2541  if (!free_dofs)
2542  return;
2543 
2544  // The element matrix and RHS for projections.
2545  // Note that Ke is always real-valued, whereas
2546  // Fe may be complex valued if complex number
2547  // support is enabled
2548  DenseMatrix<Real> Ke(free_dofs, free_dofs);
2549  DenseVector<FValue> Fe(free_dofs);
2550  // The new degree of freedom coefficients to solve for
2551  DenseVector<FValue> Ufree(free_dofs);
2552 
2553  const unsigned int n_qp =
2554  cast_int<unsigned int>(xyz_values.size());
2555 
2556  // Loop over the quadrature points
2557  for (unsigned int qp=0; qp<n_qp; qp++)
2558  {
2559  // solution at the quadrature point
2560  FValue fineval = f.eval_at_point(context,
2561  var_component,
2562  xyz_values[qp],
2563  system.time);
2564  // solution grad at the quadrature point
2565  VectorValue<FValue> finegrad;
2566  if (cont == C_ONE)
2567  finegrad = g->eval_at_point(context,
2568  var_component,
2569  xyz_values[qp],
2570  system.time);
2571 
2572  // Form edge projection matrix
2573  for (unsigned int sidei=0, freei=0;
2574  sidei != n_involved_dofs; ++sidei)
2575  {
2576  unsigned int i = involved_dofs[sidei];
2577  // fixed DoFs aren't test functions
2578  if (dof_is_fixed[sidei])
2579  continue;
2580  for (unsigned int sidej=0, freej=0;
2581  sidej != n_involved_dofs; ++sidej)
2582  {
2583  unsigned int j = involved_dofs[sidej];
2584  if (dof_is_fixed[sidej])
2585  Fe(freei) -= phi[i][qp] * phi[j][qp] *
2586  JxW[qp] * Uinvolved(sidej);
2587  else
2588  Ke(freei,freej) += phi[i][qp] *
2589  phi[j][qp] * JxW[qp];
2590  if (cont == C_ONE)
2591  {
2592  if (dof_is_fixed[sidej])
2593  Fe(freei) -= ( (*dphi)[i][qp] *
2594  (*dphi)[j][qp] ) *
2595  JxW[qp] * Uinvolved(sidej);
2596  else
2597  Ke(freei,freej) += ( (*dphi)[i][qp] *
2598  (*dphi)[j][qp] )
2599  * JxW[qp];
2600  }
2601  if (!dof_is_fixed[sidej])
2602  freej++;
2603  }
2604  Fe(freei) += phi[i][qp] * fineval * JxW[qp];
2605  if (cont == C_ONE)
2606  Fe(freei) += (finegrad * (*dphi)[i][qp] ) *
2607  JxW[qp];
2608  freei++;
2609  }
2610  }
2611 
2612  Ke.cholesky_solve(Fe, Ufree);
2613 
2614  // Transfer new edge solutions to element
2615  const processor_id_type pid = node ?
2616  node->processor_id() : DofObject::invalid_processor_id;
2617  insert_ids(free_dof_ids, Ufree.get_values(), pid);
2618 }

References libMesh::C_ONE, libMesh::DenseMatrix< T >::cholesky_solve(), libMesh::FEAbstract::get_continuity(), libMesh::FEGenericBase< OutputType >::get_dphi(), libMesh::FEAbstract::get_JxW(), libMesh::FEGenericBase< OutputType >::get_phi(), libMesh::DenseVector< T >::get_values(), libMesh::FEAbstract::get_xyz(), libMesh::DofObject::invalid_processor_id, libMesh::DofObject::processor_id(), and libMesh::DenseVector< T >::size().

◆ construct_projection() [2/2]

template<typename FFunctor , typename GFunctor , typename FValue , typename ProjectionAction >
void libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::SubProjector::construct_projection

Definition at line 2501 of file generic_projector.h.

2506 {
2507  const std::vector<Real> & JxW = fe.get_JxW();
2508  const std::vector<std::vector<Real>> & phi = fe.get_phi();
2509  const std::vector<std::vector<RealGradient>> * dphi = nullptr;
2510  const std::vector<Point> & xyz_values = fe.get_xyz();
2511  const FEContinuity cont = fe.get_continuity();
2512  const std::unordered_map<dof_id_type, FValue> & ids_to_save =
2513  this->projector.ids_to_save;
2514 
2515  if (cont == C_ONE)
2516  dphi = &(fe.get_dphi());
2517 
2518  const unsigned int n_involved_dofs =
2519  cast_int<unsigned int>(involved_dofs.size());
2520 
2521  std::vector<dof_id_type> free_dof_ids;
2522  DenseVector<FValue> Uinvolved(n_involved_dofs);
2523  std::vector<char> dof_is_fixed(n_involved_dofs, false); // bools
2524 
2525  for (auto i : IntRange<unsigned int>(0, n_involved_dofs))
2526  {
2527  const dof_id_type id = dof_indices_var[involved_dofs[i]];
2528  auto iter = ids_to_save.find(id);
2529  if (iter == ids_to_save.end())
2530  free_dof_ids.push_back(id);
2531  else
2532  {
2533  dof_is_fixed[i] = true;
2534  Uinvolved(i) = iter->second;
2535  }
2536  }
2537 
2538  const unsigned int free_dofs = free_dof_ids.size();
2539 
2540  // There may be nothing to project
2541  if (!free_dofs)
2542  return;
2543 
2544  // The element matrix and RHS for projections.
2545  // Note that Ke is always real-valued, whereas
2546  // Fe may be complex valued if complex number
2547  // support is enabled
2548  DenseMatrix<Real> Ke(free_dofs, free_dofs);
2549  DenseVector<FValue> Fe(free_dofs);
2550  // The new degree of freedom coefficients to solve for
2551  DenseVector<FValue> Ufree(free_dofs);
2552 
2553  const unsigned int n_qp =
2554  cast_int<unsigned int>(xyz_values.size());
2555 
2556  // Loop over the quadrature points
2557  for (unsigned int qp=0; qp<n_qp; qp++)
2558  {
2559  // solution at the quadrature point
2560  FValue fineval = f.eval_at_point(context,
2561  var_component,
2562  xyz_values[qp],
2563  system.time);
2564  // solution grad at the quadrature point
2565  VectorValue<FValue> finegrad;
2566  if (cont == C_ONE)
2567  finegrad = g->eval_at_point(context,
2568  var_component,
2569  xyz_values[qp],
2570  system.time);
2571 
2572  // Form edge projection matrix
2573  for (unsigned int sidei=0, freei=0;
2574  sidei != n_involved_dofs; ++sidei)
2575  {
2576  unsigned int i = involved_dofs[sidei];
2577  // fixed DoFs aren't test functions
2578  if (dof_is_fixed[sidei])
2579  continue;
2580  for (unsigned int sidej=0, freej=0;
2581  sidej != n_involved_dofs; ++sidej)
2582  {
2583  unsigned int j = involved_dofs[sidej];
2584  if (dof_is_fixed[sidej])
2585  Fe(freei) -= phi[i][qp] * phi[j][qp] *
2586  JxW[qp] * Uinvolved(sidej);
2587  else
2588  Ke(freei,freej) += phi[i][qp] *
2589  phi[j][qp] * JxW[qp];
2590  if (cont == C_ONE)
2591  {
2592  if (dof_is_fixed[sidej])
2593  Fe(freei) -= ( (*dphi)[i][qp] *
2594  (*dphi)[j][qp] ) *
2595  JxW[qp] * Uinvolved(sidej);
2596  else
2597  Ke(freei,freej) += ( (*dphi)[i][qp] *
2598  (*dphi)[j][qp] )
2599  * JxW[qp];
2600  }
2601  if (!dof_is_fixed[sidej])
2602  freej++;
2603  }
2604  Fe(freei) += phi[i][qp] * fineval * JxW[qp];
2605  if (cont == C_ONE)
2606  Fe(freei) += (finegrad * (*dphi)[i][qp] ) *
2607  JxW[qp];
2608  freei++;
2609  }
2610  }
2611 
2612  Ke.cholesky_solve(Fe, Ufree);
2613 
2614  // Transfer new edge solutions to element
2615  const processor_id_type pid = node ?
2616  node->processor_id() : DofObject::invalid_processor_id;
2617  insert_ids(free_dof_ids, Ufree.get_values(), pid);
2618 }

◆ find_dofs_to_send()

template<typename FFunctor , typename GFunctor , typename FValue , typename ProjectionAction >
void libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::SubFunctor::find_dofs_to_send ( const Node node,
const Elem elem,
unsigned short  node_num,
const var_set vars 
)
inherited

Definition at line 2345 of file generic_projector.h.

2346 {
2347  libmesh_assert (&node == elem.node_ptr(node_num));
2348 
2349  // Any ghosted node in our range might have an owner who needs our
2350  // data
2351  const processor_id_type owner = node.processor_id();
2352  if (owner != system.processor_id())
2353  {
2354  const MeshBase & mesh = system.get_mesh();
2355  const DofMap & dof_map = system.get_dof_map();
2356 
2357  // But let's check and see if we can be certain the owner can
2358  // compute any or all of its own dof coefficients on that node.
2359  std::vector<dof_id_type> node_dof_ids, patch_dof_ids;
2360  for (const auto & var : vars)
2361  {
2362  const Variable & variable = system.variable(var);
2363 
2364  if (!variable.active_on_subdomain(elem.subdomain_id()))
2365  continue;
2366 
2367  dof_map.dof_indices(elem, node_num, node_dof_ids, var);
2368  }
2369  libmesh_assert(std::is_sorted(node_dof_ids.begin(),
2370  node_dof_ids.end()));
2371  const std::vector<dof_id_type> & patch =
2372  (*this->projector.nodes_to_elem)[node.id()];
2373  for (const auto & elem_id : patch)
2374  {
2375  const Elem & patch_elem = mesh.elem_ref(elem_id);
2376  if (!patch_elem.active() || owner != patch_elem.processor_id())
2377  continue;
2378  dof_map.dof_indices(&patch_elem, patch_dof_ids);
2379  std::sort(patch_dof_ids.begin(), patch_dof_ids.end());
2380 
2381  // Remove any node_dof_ids that we see can be calculated on
2382  // this element
2383  std::vector<dof_id_type> diff_ids(node_dof_ids.size());
2384  auto it = std::set_difference(node_dof_ids.begin(), node_dof_ids.end(),
2385  patch_dof_ids.begin(), patch_dof_ids.end(), diff_ids.begin());
2386  diff_ids.resize(it-diff_ids.begin());
2387  node_dof_ids.swap(diff_ids);
2388  if (node_dof_ids.empty())
2389  break;
2390  }
2391 
2392  // Give ids_to_push default invalid pid: not yet computed
2393  for (auto id : node_dof_ids)
2395  }
2396 }

References libMesh::Elem::active(), libMesh::Variable::active_on_subdomain(), libMesh::DofMap::dof_indices(), libMesh::MeshBase::elem_ref(), libMesh::DofObject::id(), libMesh::DofObject::invalid_processor_id, libMesh::Parallel::Utils::is_sorted(), libMesh::libmesh_assert(), mesh, libMesh::Elem::node_ptr(), libMesh::DofObject::processor_id(), and libMesh::Elem::subdomain_id().

◆ insert_id() [1/2]

template<typename FFunctor , typename GFunctor , typename FValue , typename ProjectionAction >
void libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::SubFunctor::insert_id ( dof_id_type  id,
const FValue &  val,
processor_id_type  pid 
)
inherited

Definition at line 1117 of file generic_projector.h.

1118 {
1119  auto iter = new_ids_to_push.find(id);
1120  if (iter == new_ids_to_push.end())
1121  action.insert(id, val);
1122  else
1123  {
1125  iter->second = std::make_pair(val, pid);
1126  }
1127  if (!this->projector.done_saving_ids)
1128  {
1129  libmesh_assert(!new_ids_to_save.count(id));
1130  new_ids_to_save[id] = val;
1131  }
1132 }

References libMesh::DofObject::invalid_processor_id, and libMesh::libmesh_assert().

◆ insert_id() [2/2]

template<typename FFunctor , typename GFunctor , typename FValue , typename ProjectionAction >
void libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::SubFunctor::insert_id

Definition at line 1117 of file generic_projector.h.

1118 {
1119  auto iter = new_ids_to_push.find(id);
1120  if (iter == new_ids_to_push.end())
1121  action.insert(id, val);
1122  else
1123  {
1125  iter->second = std::make_pair(val, pid);
1126  }
1127  if (!this->projector.done_saving_ids)
1128  {
1129  libmesh_assert(!new_ids_to_save.count(id));
1130  new_ids_to_save[id] = val;
1131  }
1132 }

◆ insert_ids() [1/2]

template<typename FFunctor , typename GFunctor , typename FValue , typename ProjectionAction >
void libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::SubFunctor::insert_ids ( const std::vector< dof_id_type > &  ids,
const std::vector< FValue > &  vals,
processor_id_type  pid 
)
inherited

Definition at line 1138 of file generic_projector.h.

1139 {
1140  libmesh_assert_equal_to(ids.size(), vals.size());
1141  for (auto i : index_range(ids))
1142  {
1143  const dof_id_type id = ids[i];
1144  const FValue & val = vals[i];
1145 
1146  auto iter = new_ids_to_push.find(id);
1147  if (iter == new_ids_to_push.end())
1148  action.insert(id, val);
1149  else
1150  {
1152  iter->second = std::make_pair(val, pid);
1153  }
1154  if (!this->projector.done_saving_ids)
1155  {
1156  libmesh_assert(!new_ids_to_save.count(id));
1157  new_ids_to_save[id] = val;
1158  }
1159  }
1160 }

References libMesh::index_range(), libMesh::DofObject::invalid_processor_id, and libMesh::libmesh_assert().

◆ insert_ids() [2/2]

template<typename FFunctor , typename GFunctor , typename FValue , typename ProjectionAction >
void libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::SubFunctor::insert_ids

Definition at line 1138 of file generic_projector.h.

1139 {
1140  libmesh_assert_equal_to(ids.size(), vals.size());
1141  for (auto i : index_range(ids))
1142  {
1143  const dof_id_type id = ids[i];
1144  const FValue & val = vals[i];
1145 
1146  auto iter = new_ids_to_push.find(id);
1147  if (iter == new_ids_to_push.end())
1148  action.insert(id, val);
1149  else
1150  {
1152  iter->second = std::make_pair(val, pid);
1153  }
1154  if (!this->projector.done_saving_ids)
1155  {
1156  libmesh_assert(!new_ids_to_save.count(id));
1157  new_ids_to_save[id] = val;
1158  }
1159  }
1160 }

◆ join()

template<typename FFunctor , typename GFunctor , typename FValue , typename ProjectionAction >
void libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::SubFunctor::join ( const SubFunctor other)
inherited

Definition at line 1166 of file generic_projector.h.

1167 {
1168  new_ids_to_push.insert(other.new_ids_to_push.begin(), other.new_ids_to_push.end());
1169  new_ids_to_save.insert(other.new_ids_to_save.begin(), other.new_ids_to_save.end());
1170 }

References libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::SubFunctor::new_ids_to_push, and libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::SubFunctor::new_ids_to_save.

◆ operator()()

template<typename FFunctor , typename GFunctor , typename FValue , typename ProjectionAction >
void libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::ProjectEdges::operator() ( const node_range range)

Definition at line 1914 of file generic_projector.h.

1915 {
1916  LOG_SCOPE ("project_edges","GenericProjector");
1917 
1918  const unsigned int sys_num = system.number();
1919 
1920  for (const auto & e_pair : range)
1921  {
1922  const Elem & elem = *std::get<0>(e_pair.second);
1923 
1924  // If this is an unchanged element then we already copied all
1925  // its dofs
1926 #ifdef LIBMESH_ENABLE_AMR
1927  if (f.is_grid_projection() &&
1928  (elem.refinement_flag() != Elem::JUST_REFINED &&
1929  elem.refinement_flag() != Elem::JUST_COARSENED &&
1930  elem.p_refinement_flag() != Elem::JUST_REFINED &&
1931  elem.p_refinement_flag() != Elem::JUST_COARSENED))
1932  continue;
1933 #endif // LIBMESH_ENABLE_AMR
1934 
1935  const Node & edge_node = *e_pair.first;
1936  const int dim = elem.dim();
1937  const var_set & edge_vars = std::get<2>(e_pair.second);
1938 
1939  const unsigned short edge_num = std::get<1>(e_pair.second);
1940  const unsigned short node_num = elem.n_vertices() + edge_num;
1941  context.edge = edge_num;
1942  context.pre_fe_reinit(system, &elem);
1943 
1944  this->find_dofs_to_send(edge_node, elem, node_num, edge_vars);
1945 
1946  // Look at all the variables we're supposed to interpolate from
1947  // this element on this edge
1948  for (const auto & var : edge_vars)
1949  {
1950  const Variable & variable = system.variable(var);
1951  const FEType & base_fe_type = variable.type();
1952  const unsigned int var_component =
1954 
1955  if (base_fe_type.family == SCALAR)
1956  continue;
1957 
1958  FEType fe_type = base_fe_type;
1959 
1960  // This may be a p refined element
1961  fe_type.order =
1962  libMesh::Order (fe_type.order + elem.p_level());
1963 
1964  // If this is a Lagrange element with DoFs on edges then by
1965  // convention we interpolate at the node rather than project
1966  // along the edge.
1967  if (fe_type.family == LAGRANGE)
1968  {
1969  if (fe_type.order > 1)
1970  {
1971  const dof_id_type dof_id =
1972  edge_node.dof_number(sys_num, var, 0);
1973  const processor_id_type pid =
1974  edge_node.processor_id();
1975  FValue fval = f.eval_at_point
1976  (context, var_component, edge_node, system.time);
1977  insert_id(dof_id, fval, pid);
1978  }
1979  continue;
1980  }
1981 
1982 #ifdef LIBMESH_ENABLE_AMR
1983  // If this is a low order monomial element which has merely
1984  // been h refined then we already copied all its dofs
1985  if (fe_type.family == MONOMIAL &&
1986  fe_type.order == CONSTANT &&
1987  elem.refinement_flag() != Elem::JUST_COARSENED &&
1988  elem.p_refinement_flag() != Elem::JUST_COARSENED)
1989  continue;
1990 #endif // LIBMESH_ENABLE_AMR
1991 
1992  // FIXME: Need to generalize this to vector-valued elements. [PB]
1993  FEBase * fe = nullptr;
1994  context.get_element_fe( var, fe, dim );
1995  FEBase * edge_fe = nullptr;
1996  context.get_edge_fe( var, edge_fe );
1997 
1998  // If we're JUST_COARSENED we'll need a custom
1999  // evaluation, not just the standard edge FE
2000  const FEBase & proj_fe =
2001 #ifdef LIBMESH_ENABLE_AMR
2002  (elem.refinement_flag() == Elem::JUST_COARSENED) ?
2003  *fe :
2004 #endif
2005  *edge_fe;
2006 
2007 #ifdef LIBMESH_ENABLE_AMR
2008  if (elem.refinement_flag() == Elem::JUST_COARSENED)
2009  {
2010  std::vector<Point> fine_points;
2011 
2012  std::unique_ptr<FEBase> fine_fe
2013  (FEBase::build (dim, base_fe_type));
2014 
2015  std::unique_ptr<QBase> qrule
2016  (base_fe_type.default_quadrature_rule(1));
2017  fine_fe->attach_quadrature_rule(qrule.get());
2018 
2019  const std::vector<Point> & child_xyz =
2020  fine_fe->get_xyz();
2021 
2022  for (unsigned int c = 0, nc = elem.n_children();
2023  c != nc; ++c)
2024  {
2025  if (!elem.is_child_on_edge(c, context.edge))
2026  continue;
2027 
2028  fine_fe->edge_reinit(elem.child_ptr(c), context.edge);
2029  fine_points.insert(fine_points.end(),
2030  child_xyz.begin(),
2031  child_xyz.end());
2032  }
2033 
2034  std::vector<Point> fine_qp;
2035  FEMap::inverse_map (dim, &elem, fine_points, fine_qp);
2036 
2037  context.elem_fe_reinit(&fine_qp);
2038  }
2039  else
2040 #endif // LIBMESH_ENABLE_AMR
2042 
2043  const std::vector<dof_id_type> & dof_indices =
2044  context.get_dof_indices(var);
2045 
2046  std::vector<unsigned int> edge_dofs;
2047  FEInterface::dofs_on_edge(&elem, dim, base_fe_type,
2048  context.edge, edge_dofs);
2049 
2050  this->construct_projection
2051  (dof_indices, edge_dofs, var_component,
2052  &edge_node, proj_fe);
2053  }
2054  }
2055 }

References libMesh::FEGenericBase< OutputType >::build(), libMesh::Elem::child_ptr(), libMesh::CONSTANT, libMesh::FEType::default_quadrature_rule(), dim, libMesh::Elem::dim(), libMesh::DofObject::dof_number(), libMesh::FEInterface::dofs_on_edge(), libMesh::FEType::family, libMesh::FEMap::inverse_map(), libMesh::Elem::is_child_on_edge(), libMesh::Elem::JUST_COARSENED, libMesh::Elem::JUST_REFINED, libMesh::LAGRANGE, libMesh::MONOMIAL, libMesh::Elem::n_children(), libMesh::Elem::n_vertices(), libMesh::FEType::order, libMesh::Elem::p_level(), libMesh::Elem::p_refinement_flag(), libMesh::DofObject::processor_id(), libMesh::Elem::refinement_flag(), libMesh::SCALAR, and libMesh::Variable::type().

Member Data Documentation

◆ action

template<typename FFunctor , typename GFunctor , typename FValue , typename ProjectionAction >
ProjectionAction libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::SubFunctor::action

Definition at line 193 of file generic_projector.h.

◆ context

template<typename FFunctor , typename GFunctor , typename FValue , typename ProjectionAction >
FEMContext libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::SubFunctor::context

Definition at line 198 of file generic_projector.h.

◆ conts

template<typename FFunctor , typename GFunctor , typename FValue , typename ProjectionAction >
std::vector<FEContinuity> libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::SubFunctor::conts

Definition at line 201 of file generic_projector.h.

◆ f

template<typename FFunctor , typename GFunctor , typename FValue , typename ProjectionAction >
FFunctor libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::SubFunctor::f

Definition at line 194 of file generic_projector.h.

◆ g

template<typename FFunctor , typename GFunctor , typename FValue , typename ProjectionAction >
std::unique_ptr<GFunctor> libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::SubProjector::g

Definition at line 227 of file generic_projector.h.

◆ new_ids_to_push

template<typename FFunctor , typename GFunctor , typename FValue , typename ProjectionAction >
std::unordered_map<dof_id_type, std::pair<FValue, processor_id_type> > libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::SubFunctor::new_ids_to_push
inherited

◆ new_ids_to_save

template<typename FFunctor , typename GFunctor , typename FValue , typename ProjectionAction >
std::unordered_map<dof_id_type, FValue> libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::SubFunctor::new_ids_to_save
inherited

◆ projector

template<typename FFunctor , typename GFunctor , typename FValue , typename ProjectionAction >
GenericProjector& libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::SubFunctor::projector
inherited

Definition at line 154 of file generic_projector.h.

◆ system

template<typename FFunctor , typename GFunctor , typename FValue , typename ProjectionAction >
const System& libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::SubFunctor::system

Definition at line 203 of file generic_projector.h.


The documentation for this struct was generated from the following file:
libMesh::C_ONE
Definition: enum_fe_family.h:80
libMesh::dof_id_type
uint8_t dof_id_type
Definition: id_types.h:67
libMesh::GenericProjector::SubFunctor::new_ids_to_save
std::unordered_map< dof_id_type, FValue > new_ids_to_save
Definition: generic_projector.h:169
libMesh::Elem::JUST_REFINED
Definition: elem.h:1172
libMesh::GenericProjector::var_set
std::set< unsigned int > var_set
Definition: generic_projector.h:149
libMesh::DiffContext::get_dof_indices
const std::vector< dof_id_type > & get_dof_indices() const
Accessor for element dof indices.
Definition: diff_context.h:367
libMesh::GenericProjector::SubFunctor::insert_id
void insert_id(dof_id_type id, const FValue &val, processor_id_type pid)
Definition: generic_projector.h:1117
libMesh::index_range
IntRange< std::size_t > index_range(const std::vector< T > &vec)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
Definition: int_range.h:106
libMesh::FEMContext::pre_fe_reinit
virtual void pre_fe_reinit(const System &, const Elem *e)
Reinitializes local data vectors/matrices on the current geometric element.
Definition: fem_context.C:1642
libMesh::Order
Order
Definition: enum_order.h:40
libMesh::FEMContext::get_edge_fe
void get_edge_fe(unsigned int var, FEGenericBase< OutputShape > *&fe) const
Accessor for edge (3D only!) finite element object for variable var.
Definition: fem_context.h:1282
libMesh::System::number
unsigned int number() const
Definition: system.h:2075
libMesh::Variable::type
const FEType & type() const
Definition: variable.h:119
mesh
MeshBase & mesh
Definition: mesh_communication.C:1257
libMesh::GenericProjector::SubFunctor::new_ids_to_push
std::unordered_map< dof_id_type, std::pair< FValue, processor_id_type > > new_ids_to_push
Definition: generic_projector.h:165
libMesh::System::variable
const Variable & variable(unsigned int var) const
Return a constant reference to Variable var.
Definition: system.h:2183
libMesh::GenericProjector::SubProjector::SubProjector
SubProjector(GenericProjector &p)
Definition: generic_projector.h:1096
libMesh::Parallel::Utils::is_sorted
bool is_sorted(const std::vector< KeyType > &v)
Definition: parallel_conversion_utils.h:52
dim
unsigned int dim
Definition: adaptivity_ex3.C:113
libMesh::libmesh_assert
libmesh_assert(ctx)
libMesh::GenericProjector::SubFunctor::insert_ids
void insert_ids(const std::vector< dof_id_type > &ids, const std::vector< FValue > &vals, processor_id_type pid)
Definition: generic_projector.h:1138
libMesh::GenericProjector::SubFunctor::find_dofs_to_send
void find_dofs_to_send(const Node &node, const Elem &elem, unsigned short node_num, const var_set &vars)
Definition: generic_projector.h:2345
libMesh::GenericProjector::done_saving_ids
bool done_saving_ids
Definition: generic_projector.h:93
libMesh::FEMap::inverse_map
static Point inverse_map(const unsigned int dim, const Elem *elem, const Point &p, const Real tolerance=TOLERANCE, const bool secure=true)
Definition: fe_map.C:1622
libMesh::GenericProjector::SubProjector::construct_projection
void construct_projection(const std::vector< dof_id_type > &dof_indices_var, const std::vector< unsigned int > &involved_dofs, unsigned int var_component, const Node *node, const FEBase &fe)
Definition: generic_projector.h:2501
libMesh::GenericProjector::SubFunctor::context
FEMContext context
Definition: generic_projector.h:198
libMesh::System::variable_scalar_number
unsigned int variable_scalar_number(const std::string &var, unsigned int component) const
Definition: system.h:2214
libMesh::FEBase
FEGenericBase< Real > FEBase
Definition: exact_error_estimator.h:39
libMesh::System::get_mesh
const MeshBase & get_mesh() const
Definition: system.h:2083
libMesh::ParallelObject::processor_id
processor_id_type processor_id() const
Definition: parallel_object.h:106
libMesh::FEGenericBase::build
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
libMesh::FEInterface::dofs_on_edge
static void dofs_on_edge(const Elem *const elem, const unsigned int dim, const FEType &fe_t, unsigned int e, std::vector< unsigned int > &di)
Fills the vector di with the local degree of freedom indices associated with edge e of element elem A...
Definition: fe_interface.C:566
libMesh::processor_id_type
uint8_t processor_id_type
Definition: id_types.h:104
libMesh::FEMContext::edge
unsigned char edge
Current edge for edge_* to examine.
Definition: fem_context.h:1015
libMesh::CONSTANT
Definition: enum_order.h:41
libMesh::GenericProjector::nodes_to_elem
std::unordered_map< dof_id_type, std::vector< dof_id_type > > * nodes_to_elem
Definition: generic_projector.h:92
libMesh::Elem::JUST_COARSENED
Definition: elem.h:1173
libMesh::MONOMIAL
Definition: enum_fe_family.h:39
libMesh::System::time
Real time
For time-dependent problems, this is the time t at the beginning of the current timestep.
Definition: system.h:1561
libMesh::FEMContext::elem_fe_reinit
virtual void elem_fe_reinit(const std::vector< Point > *const pts=nullptr)
Reinitializes interior FE objects on the current geometric element.
Definition: fem_context.C:1436
libMesh::FEMContext::edge_fe_reinit
virtual void edge_fe_reinit()
Reinitializes edge FE objects on the current geometric element.
Definition: fem_context.C:1473
libMesh::FEMContext::get_element_fe
void get_element_fe(unsigned int var, FEGenericBase< OutputShape > *&fe) const
Accessor for interior finite element object for variable var for the largest dimension in the mesh.
Definition: fem_context.h:275
libMesh::DofObject::invalid_processor_id
static const processor_id_type invalid_processor_id
An invalid processor_id to distinguish DoFs that have not been assigned to a processor.
Definition: dof_object.h:432
libMesh::GenericProjector::SubFunctor::action
ProjectionAction action
Definition: generic_projector.h:193
libMesh::GenericProjector::SubFunctor::projector
GenericProjector & projector
Definition: generic_projector.h:154
libMesh::GenericProjector::ids_to_save
std::unordered_map< dof_id_type, FValue > ids_to_save
Definition: generic_projector.h:145
libMesh::System::get_dof_map
const DofMap & get_dof_map() const
Definition: system.h:2099
libMesh::GenericProjector::SubFunctor::f
FFunctor f
Definition: generic_projector.h:194
libMesh::LAGRANGE
Definition: enum_fe_family.h:36
libMesh::SCALAR
Definition: enum_fe_family.h:58
libMesh::GenericProjector::SubFunctor::system
const System & system
Definition: generic_projector.h:203
libMesh::FEContinuity
FEContinuity
Definition: enum_fe_family.h:77
libMesh::GenericProjector::SubProjector::g
std::unique_ptr< GFunctor > g
Definition: generic_projector.h:227