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 282 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)

Definition at line 283 of file generic_projector.h.

283 : 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   
)

Definition at line 285 of file generic_projector.h.

285 : 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 2478 of file generic_projector.h.

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().

2483 {
2484  const std::vector<Real> & JxW = fe.get_JxW();
2485  const std::vector<std::vector<Real>> & phi = fe.get_phi();
2486  const std::vector<std::vector<RealGradient>> * dphi = nullptr;
2487  const std::vector<Point> & xyz_values = fe.get_xyz();
2488  const FEContinuity cont = fe.get_continuity();
2489  const std::unordered_map<dof_id_type, FValue> & ids_to_save =
2490  this->projector.ids_to_save;
2491 
2492  if (cont == C_ONE)
2493  dphi = &(fe.get_dphi());
2494 
2495  const unsigned int n_involved_dofs =
2496  cast_int<unsigned int>(involved_dofs.size());
2497 
2498  std::vector<dof_id_type> free_dof_ids;
2499  DenseVector<FValue> Uinvolved(n_involved_dofs);
2500  std::vector<char> dof_is_fixed(n_involved_dofs, false); // bools
2501 
2502  for (auto i : IntRange<unsigned int>(0, n_involved_dofs))
2503  {
2504  const dof_id_type id = dof_indices_var[involved_dofs[i]];
2505  auto iter = ids_to_save.find(id);
2506  if (iter == ids_to_save.end())
2507  free_dof_ids.push_back(id);
2508  else
2509  {
2510  dof_is_fixed[i] = true;
2511  Uinvolved(i) = iter->second;
2512  }
2513  }
2514 
2515  const unsigned int free_dofs = free_dof_ids.size();
2516 
2517  // There may be nothing to project
2518  if (!free_dofs)
2519  return;
2520 
2521  // The element matrix and RHS for projections.
2522  // Note that Ke is always real-valued, whereas
2523  // Fe may be complex valued if complex number
2524  // support is enabled
2525  DenseMatrix<Real> Ke(free_dofs, free_dofs);
2526  DenseVector<FValue> Fe(free_dofs);
2527  // The new degree of freedom coefficients to solve for
2528  DenseVector<FValue> Ufree(free_dofs);
2529 
2530  const unsigned int n_qp =
2531  cast_int<unsigned int>(xyz_values.size());
2532 
2533  // Loop over the quadrature points
2534  for (unsigned int qp=0; qp<n_qp; qp++)
2535  {
2536  // solution at the quadrature point
2537  FValue fineval = f.eval_at_point(context,
2538  var_component,
2539  xyz_values[qp],
2540  system.time);
2541  // solution grad at the quadrature point
2542  VectorValue<FValue> finegrad;
2543  if (cont == C_ONE)
2544  finegrad = g->eval_at_point(context,
2545  var_component,
2546  xyz_values[qp],
2547  system.time);
2548 
2549  // Form edge projection matrix
2550  for (unsigned int sidei=0, freei=0;
2551  sidei != n_involved_dofs; ++sidei)
2552  {
2553  unsigned int i = involved_dofs[sidei];
2554  // fixed DoFs aren't test functions
2555  if (dof_is_fixed[sidei])
2556  continue;
2557  for (unsigned int sidej=0, freej=0;
2558  sidej != n_involved_dofs; ++sidej)
2559  {
2560  unsigned int j = involved_dofs[sidej];
2561  if (dof_is_fixed[sidej])
2562  Fe(freei) -= phi[i][qp] * phi[j][qp] *
2563  JxW[qp] * Uinvolved(sidej);
2564  else
2565  Ke(freei,freej) += phi[i][qp] *
2566  phi[j][qp] * JxW[qp];
2567  if (cont == C_ONE)
2568  {
2569  if (dof_is_fixed[sidej])
2570  Fe(freei) -= ( (*dphi)[i][qp] *
2571  (*dphi)[j][qp] ) *
2572  JxW[qp] * Uinvolved(sidej);
2573  else
2574  Ke(freei,freej) += ( (*dphi)[i][qp] *
2575  (*dphi)[j][qp] )
2576  * JxW[qp];
2577  }
2578  if (!dof_is_fixed[sidej])
2579  freej++;
2580  }
2581  Fe(freei) += phi[i][qp] * fineval * JxW[qp];
2582  if (cont == C_ONE)
2583  Fe(freei) += (finegrad * (*dphi)[i][qp] ) *
2584  JxW[qp];
2585  freei++;
2586  }
2587  }
2588 
2589  Ke.cholesky_solve(Fe, Ufree);
2590 
2591  // Transfer new edge solutions to element
2592  const processor_id_type pid = node ?
2593  node->processor_id() : DofObject::invalid_processor_id;
2594  insert_ids(free_dof_ids, Ufree.get_values(), pid);
2595 }
Real time
For time-dependent problems, this is the time t at the beginning of the current timestep.
Definition: system.h:1545
uint8_t processor_id_type
Definition: id_types.h:104
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:404
std::unordered_map< dof_id_type, FValue > ids_to_save
FEContinuity
defines an enum for finite element types to libmesh_assert a certain level (or type? Hcurl?) of continuity.
void insert_ids(const std::vector< dof_id_type > &ids, const std::vector< FValue > &vals, processor_id_type pid)
uint8_t dof_id_type
Definition: id_types.h:67

◆ 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 2478 of file generic_projector.h.

2483 {
2484  const std::vector<Real> & JxW = fe.get_JxW();
2485  const std::vector<std::vector<Real>> & phi = fe.get_phi();
2486  const std::vector<std::vector<RealGradient>> * dphi = nullptr;
2487  const std::vector<Point> & xyz_values = fe.get_xyz();
2488  const FEContinuity cont = fe.get_continuity();
2489  const std::unordered_map<dof_id_type, FValue> & ids_to_save =
2490  this->projector.ids_to_save;
2491 
2492  if (cont == C_ONE)
2493  dphi = &(fe.get_dphi());
2494 
2495  const unsigned int n_involved_dofs =
2496  cast_int<unsigned int>(involved_dofs.size());
2497 
2498  std::vector<dof_id_type> free_dof_ids;
2499  DenseVector<FValue> Uinvolved(n_involved_dofs);
2500  std::vector<char> dof_is_fixed(n_involved_dofs, false); // bools
2501 
2502  for (auto i : IntRange<unsigned int>(0, n_involved_dofs))
2503  {
2504  const dof_id_type id = dof_indices_var[involved_dofs[i]];
2505  auto iter = ids_to_save.find(id);
2506  if (iter == ids_to_save.end())
2507  free_dof_ids.push_back(id);
2508  else
2509  {
2510  dof_is_fixed[i] = true;
2511  Uinvolved(i) = iter->second;
2512  }
2513  }
2514 
2515  const unsigned int free_dofs = free_dof_ids.size();
2516 
2517  // There may be nothing to project
2518  if (!free_dofs)
2519  return;
2520 
2521  // The element matrix and RHS for projections.
2522  // Note that Ke is always real-valued, whereas
2523  // Fe may be complex valued if complex number
2524  // support is enabled
2525  DenseMatrix<Real> Ke(free_dofs, free_dofs);
2526  DenseVector<FValue> Fe(free_dofs);
2527  // The new degree of freedom coefficients to solve for
2528  DenseVector<FValue> Ufree(free_dofs);
2529 
2530  const unsigned int n_qp =
2531  cast_int<unsigned int>(xyz_values.size());
2532 
2533  // Loop over the quadrature points
2534  for (unsigned int qp=0; qp<n_qp; qp++)
2535  {
2536  // solution at the quadrature point
2537  FValue fineval = f.eval_at_point(context,
2538  var_component,
2539  xyz_values[qp],
2540  system.time);
2541  // solution grad at the quadrature point
2542  VectorValue<FValue> finegrad;
2543  if (cont == C_ONE)
2544  finegrad = g->eval_at_point(context,
2545  var_component,
2546  xyz_values[qp],
2547  system.time);
2548 
2549  // Form edge projection matrix
2550  for (unsigned int sidei=0, freei=0;
2551  sidei != n_involved_dofs; ++sidei)
2552  {
2553  unsigned int i = involved_dofs[sidei];
2554  // fixed DoFs aren't test functions
2555  if (dof_is_fixed[sidei])
2556  continue;
2557  for (unsigned int sidej=0, freej=0;
2558  sidej != n_involved_dofs; ++sidej)
2559  {
2560  unsigned int j = involved_dofs[sidej];
2561  if (dof_is_fixed[sidej])
2562  Fe(freei) -= phi[i][qp] * phi[j][qp] *
2563  JxW[qp] * Uinvolved(sidej);
2564  else
2565  Ke(freei,freej) += phi[i][qp] *
2566  phi[j][qp] * JxW[qp];
2567  if (cont == C_ONE)
2568  {
2569  if (dof_is_fixed[sidej])
2570  Fe(freei) -= ( (*dphi)[i][qp] *
2571  (*dphi)[j][qp] ) *
2572  JxW[qp] * Uinvolved(sidej);
2573  else
2574  Ke(freei,freej) += ( (*dphi)[i][qp] *
2575  (*dphi)[j][qp] )
2576  * JxW[qp];
2577  }
2578  if (!dof_is_fixed[sidej])
2579  freej++;
2580  }
2581  Fe(freei) += phi[i][qp] * fineval * JxW[qp];
2582  if (cont == C_ONE)
2583  Fe(freei) += (finegrad * (*dphi)[i][qp] ) *
2584  JxW[qp];
2585  freei++;
2586  }
2587  }
2588 
2589  Ke.cholesky_solve(Fe, Ufree);
2590 
2591  // Transfer new edge solutions to element
2592  const processor_id_type pid = node ?
2593  node->processor_id() : DofObject::invalid_processor_id;
2594  insert_ids(free_dof_ids, Ufree.get_values(), pid);
2595 }
Real time
For time-dependent problems, this is the time t at the beginning of the current timestep.
Definition: system.h:1545
uint8_t processor_id_type
Definition: id_types.h:104
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:404
std::unordered_map< dof_id_type, FValue > ids_to_save
FEContinuity
defines an enum for finite element types to libmesh_assert a certain level (or type? Hcurl?) of continuity.
void insert_ids(const std::vector< dof_id_type > &ids, const std::vector< FValue > &vals, processor_id_type pid)
uint8_t dof_id_type
Definition: id_types.h:67

◆ 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 2322 of file generic_projector.h.

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::Elem::node_ptr(), libMesh::DofObject::processor_id(), and libMesh::Elem::subdomain_id().

2323 {
2324  libmesh_assert (&node == elem.node_ptr(node_num));
2325 
2326  // Any ghosted node in our range might have an owner who needs our
2327  // data
2328  const processor_id_type owner = node.processor_id();
2329  if (owner != system.processor_id())
2330  {
2331  const MeshBase & mesh = system.get_mesh();
2332  const DofMap & dof_map = system.get_dof_map();
2333 
2334  // But let's check and see if we can be certain the owner can
2335  // compute any or all of its own dof coefficients on that node.
2336  std::vector<dof_id_type> node_dof_ids, patch_dof_ids;
2337  for (const auto & var : vars)
2338  {
2339  const Variable & variable = system.variable(var);
2340 
2341  if (!variable.active_on_subdomain(elem.subdomain_id()))
2342  continue;
2343 
2344  dof_map.dof_indices(elem, node_num, node_dof_ids, var);
2345  }
2346  libmesh_assert(std::is_sorted(node_dof_ids.begin(),
2347  node_dof_ids.end()));
2348  const std::vector<dof_id_type> & patch =
2349  (*this->projector.nodes_to_elem)[node.id()];
2350  for (const auto & elem_id : patch)
2351  {
2352  const Elem & patch_elem = mesh.elem_ref(elem_id);
2353  if (!patch_elem.active() || owner != patch_elem.processor_id())
2354  continue;
2355  dof_map.dof_indices(&patch_elem, patch_dof_ids);
2356  std::sort(patch_dof_ids.begin(), patch_dof_ids.end());
2357 
2358  // Remove any node_dof_ids that we see can be calculated on
2359  // this element
2360  std::vector<dof_id_type> diff_ids(node_dof_ids.size());
2361  auto it = std::set_difference(node_dof_ids.begin(), node_dof_ids.end(),
2362  patch_dof_ids.begin(), patch_dof_ids.end(), diff_ids.begin());
2363  diff_ids.resize(it-diff_ids.begin());
2364  node_dof_ids.swap(diff_ids);
2365  if (node_dof_ids.empty())
2366  break;
2367  }
2368 
2369  // Give ids_to_push default invalid pid: not yet computed
2370  for (auto id : node_dof_ids)
2372  }
2373 }
const Variable & variable(unsigned int var) const
Return a constant reference to Variable var.
Definition: system.h:2167
const MeshBase & get_mesh() const
Definition: system.h:2067
uint8_t processor_id_type
Definition: id_types.h:104
bool is_sorted(const std::vector< KeyType > &v)
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:404
std::unordered_map< dof_id_type, std::vector< dof_id_type > > * nodes_to_elem
std::unordered_map< dof_id_type, std::pair< FValue, processor_id_type > > new_ids_to_push
processor_id_type processor_id() const
const DofMap & get_dof_map() const
Definition: system.h:2083

◆ 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 1112 of file generic_projector.h.

References libMesh::DofObject::invalid_processor_id.

1113 {
1114  auto iter = new_ids_to_push.find(id);
1115  if (iter == new_ids_to_push.end())
1116  action.insert(id, val);
1117  else
1118  {
1119  libmesh_assert(pid != DofObject::invalid_processor_id);
1120  iter->second = std::make_pair(val, pid);
1121  }
1122  if (!this->projector.done_saving_ids)
1123  {
1124  libmesh_assert(!new_ids_to_save.count(id));
1125  new_ids_to_save[id] = val;
1126  }
1127 }
std::unordered_map< dof_id_type, FValue > new_ids_to_save
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:404
std::unordered_map< dof_id_type, std::pair< FValue, processor_id_type > > new_ids_to_push

◆ 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 1112 of file generic_projector.h.

1113 {
1114  auto iter = new_ids_to_push.find(id);
1115  if (iter == new_ids_to_push.end())
1116  action.insert(id, val);
1117  else
1118  {
1119  libmesh_assert(pid != DofObject::invalid_processor_id);
1120  iter->second = std::make_pair(val, pid);
1121  }
1122  if (!this->projector.done_saving_ids)
1123  {
1124  libmesh_assert(!new_ids_to_save.count(id));
1125  new_ids_to_save[id] = val;
1126  }
1127 }
std::unordered_map< dof_id_type, FValue > new_ids_to_save
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:404
std::unordered_map< dof_id_type, std::pair< FValue, processor_id_type > > new_ids_to_push

◆ 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 1133 of file generic_projector.h.

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

1134 {
1135  libmesh_assert_equal_to(ids.size(), vals.size());
1136  for (auto i : index_range(ids))
1137  {
1138  const dof_id_type id = ids[i];
1139  const FValue & val = vals[i];
1140 
1141  auto iter = new_ids_to_push.find(id);
1142  if (iter == new_ids_to_push.end())
1143  action.insert(id, val);
1144  else
1145  {
1146  libmesh_assert(pid != DofObject::invalid_processor_id);
1147  iter->second = std::make_pair(val, pid);
1148  }
1149  if (!this->projector.done_saving_ids)
1150  {
1151  libmesh_assert(!new_ids_to_save.count(id));
1152  new_ids_to_save[id] = val;
1153  }
1154  }
1155 }
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:104
std::unordered_map< dof_id_type, FValue > new_ids_to_save
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:404
std::unordered_map< dof_id_type, std::pair< FValue, processor_id_type > > new_ids_to_push
uint8_t dof_id_type
Definition: id_types.h:67

◆ 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 1133 of file generic_projector.h.

1134 {
1135  libmesh_assert_equal_to(ids.size(), vals.size());
1136  for (auto i : index_range(ids))
1137  {
1138  const dof_id_type id = ids[i];
1139  const FValue & val = vals[i];
1140 
1141  auto iter = new_ids_to_push.find(id);
1142  if (iter == new_ids_to_push.end())
1143  action.insert(id, val);
1144  else
1145  {
1146  libmesh_assert(pid != DofObject::invalid_processor_id);
1147  iter->second = std::make_pair(val, pid);
1148  }
1149  if (!this->projector.done_saving_ids)
1150  {
1151  libmesh_assert(!new_ids_to_save.count(id));
1152  new_ids_to_save[id] = val;
1153  }
1154  }
1155 }
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:104
std::unordered_map< dof_id_type, FValue > new_ids_to_save
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:404
std::unordered_map< dof_id_type, std::pair< FValue, processor_id_type > > new_ids_to_push
uint8_t dof_id_type
Definition: id_types.h:67

◆ 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 1161 of file generic_projector.h.

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

1162 {
1163  new_ids_to_push.insert(other.new_ids_to_push.begin(), other.new_ids_to_push.end());
1164  new_ids_to_save.insert(other.new_ids_to_save.begin(), other.new_ids_to_save.end());
1165 }
std::unordered_map< dof_id_type, FValue > new_ids_to_save
std::unordered_map< dof_id_type, std::pair< FValue, processor_id_type > > new_ids_to_push

◆ 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 1892 of file generic_projector.h.

References libMesh::FEGenericBase< OutputType >::build(), libMesh::Elem::child_ptr(), libMesh::CONSTANT, libMesh::FEType::default_quadrature_rule(), libMesh::Elem::dim(), libMesh::DofObject::dof_number(), libMesh::FEInterface::dofs_on_edge(), libMesh::FEType::family, libMesh::FEInterface::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().

1893 {
1894  LOG_SCOPE ("project_edges","GenericProjector");
1895 
1896  const unsigned int sys_num = system.number();
1897 
1898  for (const auto & e_pair : range)
1899  {
1900  const Elem & elem = *std::get<0>(e_pair.second);
1901 
1902  // If this is an unchanged element then we already copied all
1903  // its dofs
1904 #ifdef LIBMESH_ENABLE_AMR
1905  if (f.is_grid_projection() &&
1906  (elem.refinement_flag() != Elem::JUST_REFINED &&
1907  elem.refinement_flag() != Elem::JUST_COARSENED &&
1908  elem.p_refinement_flag() != Elem::JUST_REFINED &&
1909  elem.p_refinement_flag() != Elem::JUST_COARSENED))
1910  continue;
1911 #endif // LIBMESH_ENABLE_AMR
1912 
1913  const Node & edge_node = *e_pair.first;
1914  const int dim = elem.dim();
1915  const var_set & edge_vars = std::get<2>(e_pair.second);
1916 
1917  const unsigned short edge_num = std::get<1>(e_pair.second);
1918  const unsigned short node_num = elem.n_vertices() + edge_num;
1919  context.edge = edge_num;
1920  context.pre_fe_reinit(system, &elem);
1921 
1922  this->find_dofs_to_send(edge_node, elem, node_num, edge_vars);
1923 
1924  // Look at all the variables we're supposed to interpolate from
1925  // this element on this edge
1926  for (const auto & var : edge_vars)
1927  {
1928  const Variable & variable = system.variable(var);
1929  const FEType & base_fe_type = variable.type();
1930  const unsigned int var_component =
1932 
1933  if (base_fe_type.family == SCALAR)
1934  continue;
1935 
1936  FEType fe_type = base_fe_type;
1937 
1938  // This may be a p refined element
1939  fe_type.order =
1940  libMesh::Order (fe_type.order + elem.p_level());
1941 
1942  // If this is a Lagrange element with DoFs on edges then by
1943  // convention we interpolate at the node rather than project
1944  // along the edge.
1945  if (fe_type.family == LAGRANGE)
1946  {
1947  if (fe_type.order > 1)
1948  {
1949  const dof_id_type dof_id =
1950  edge_node.dof_number(sys_num, var, 0);
1951  const processor_id_type pid =
1952  edge_node.processor_id();
1953  FValue fval = f.eval_at_point
1954  (context, var_component, edge_node, system.time);
1955  insert_id(dof_id, fval, pid);
1956  }
1957  continue;
1958  }
1959 
1960  // If this is a low order monomial element which has merely
1961  // been h refined then we already copied all its dofs
1962  if (fe_type.family == MONOMIAL &&
1963  fe_type.order == CONSTANT &&
1964  elem.refinement_flag() != Elem::JUST_COARSENED &&
1965  elem.p_refinement_flag() != Elem::JUST_COARSENED)
1966  continue;
1967 
1968  // FIXME: Need to generalize this to vector-valued elements. [PB]
1969  FEBase * fe = nullptr;
1970  context.get_element_fe( var, fe, dim );
1971  FEBase * edge_fe = nullptr;
1972  context.get_edge_fe( var, edge_fe );
1973 
1974  // If we're JUST_COARSENED we'll need a custom
1975  // evaluation, not just the standard edge FE
1976  const FEBase & proj_fe =
1977 #ifdef LIBMESH_ENABLE_AMR
1978  (elem.refinement_flag() == Elem::JUST_COARSENED) ?
1979  *fe :
1980 #endif
1981  *edge_fe;
1982 
1983 #ifdef LIBMESH_ENABLE_AMR
1984  if (elem.refinement_flag() == Elem::JUST_COARSENED)
1985  {
1986  std::vector<Point> fine_points;
1987 
1988  std::unique_ptr<FEBase> fine_fe
1989  (FEBase::build (dim, base_fe_type));
1990 
1991  std::unique_ptr<QBase> qrule
1992  (base_fe_type.default_quadrature_rule(1));
1993  fine_fe->attach_quadrature_rule(qrule.get());
1994 
1995  const std::vector<Point> & child_xyz =
1996  fine_fe->get_xyz();
1997 
1998  for (unsigned int c = 0, nc = elem.n_children();
1999  c != nc; ++c)
2000  {
2001  if (!elem.is_child_on_edge(c, context.edge))
2002  continue;
2003 
2004  fine_fe->edge_reinit(elem.child_ptr(c), context.edge);
2005  fine_points.insert(fine_points.end(),
2006  child_xyz.begin(),
2007  child_xyz.end());
2008  }
2009 
2010  std::vector<Point> fine_qp;
2011  FEInterface::inverse_map (dim, base_fe_type, &elem,
2012  fine_points, fine_qp);
2013 
2014  context.elem_fe_reinit(&fine_qp);
2015  }
2016  else
2017 #endif // LIBMESH_ENABLE_AMR
2019 
2020  const std::vector<dof_id_type> & dof_indices =
2021  context.get_dof_indices(var);
2022 
2023  std::vector<unsigned int> edge_dofs;
2024  FEInterface::dofs_on_edge(&elem, dim, base_fe_type,
2025  context.edge, edge_dofs);
2026 
2027  this->construct_projection
2028  (dof_indices, edge_dofs, var_component,
2029  &edge_node, proj_fe);
2030  }
2031  }
2032 }
Real time
For time-dependent problems, this is the time t at the beginning of the current timestep.
Definition: system.h:1545
unsigned int variable_scalar_number(const std::string &var, unsigned int component) const
Definition: system.h:2198
Order
defines an enum for polynomial orders.
Definition: enum_order.h:40
const Variable & variable(unsigned int var) const
Return a constant reference to Variable var.
Definition: system.h:2167
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:1253
void insert_id(dof_id_type id, const FValue &val, processor_id_type pid)
uint8_t processor_id_type
Definition: id_types.h:104
unsigned char edge
Current edge for edge_* to examine.
Definition: fem_context.h:996
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...
unsigned int number() const
Definition: system.h:2059
virtual void elem_fe_reinit(const std::vector< Point > *const pts=nullptr)
Reinitializes interior FE objects on the current geometric element.
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
const std::vector< dof_id_type > & get_dof_indices() const
Accessor for element dof indices.
Definition: diff_context.h:367
static Point inverse_map(const unsigned int dim, const FEType &fe_t, const Elem *elem, const Point &p, const Real tolerance=TOLERANCE, const bool secure=true)
FEGenericBase< Real > FEBase
void find_dofs_to_send(const Node &node, const Elem &elem, unsigned short node_num, const var_set &vars)
std::set< unsigned int > var_set
virtual void edge_fe_reinit()
Reinitializes edge FE objects on the current geometric element.
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:262
virtual void pre_fe_reinit(const System &, const Elem *e)
Reinitializes local data vectors/matrices on the current geometric element.
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)
uint8_t dof_id_type
Definition: id_types.h:67
const FEType & type() const
Definition: variable.h:119

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 188 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 193 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 196 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 189 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 222 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 149 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 198 of file generic_projector.h.


The documentation for this struct was generated from the following file: