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

#include <generic_projector.h>

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

Public Member Functions

 SubFunctor (GenericProjector &p)
 
void find_dofs_to_send (const Node &node, const Elem &elem, unsigned short node_num, const var_set &vars)
 
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 join (const SubFunctor &other)
 

Public Attributes

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 Attributes

ProjectionAction action
 
FFunctor f
 
FEMContext context
 
std::vector< FEContinuityconts
 
const Systemsystem
 

Detailed Description

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

Definition at line 148 of file generic_projector.h.

Constructor & Destructor Documentation

◆ SubFunctor()

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

Definition at line 1026 of file generic_projector.h.

References libMesh::C_ONE, libMesh::FEType::family, libMesh::FEAbstract::get_continuity(), libMesh::FEGenericBase< OutputType >::get_dphi(), libMesh::FEAbstract::get_fe_type(), libMesh::FEAbstract::get_JxW(), libMesh::FEGenericBase< OutputType >::get_phi(), libMesh::FEAbstract::get_xyz(), and libMesh::SCALAR.

1026  :
1027  projector(p), action(p.master_action), f(p.master_f),
1028  context(p.system), conts(p.system.n_vars()), system(p.system)
1029 {
1030  // Loop over all the variables we've been requested to project, to
1031  // pre-request
1032  for (const auto & var : this->projector.variables)
1033  {
1034  // FIXME: Need to generalize this to vector-valued elements. [PB]
1035  FEBase * fe = nullptr;
1036  FEBase * side_fe = nullptr;
1037  FEBase * edge_fe = nullptr;
1038 
1039  const std::set<unsigned char> & elem_dims =
1041 
1042  for (const auto & dim : elem_dims)
1043  {
1044  context.get_element_fe( var, fe, dim );
1045  if (fe->get_fe_type().family == SCALAR)
1046  continue;
1047  if (dim > 1)
1048  context.get_side_fe( var, side_fe, dim );
1049  if (dim > 2)
1050  context.get_edge_fe( var, edge_fe );
1051 
1052  fe->get_JxW();
1053  fe->get_xyz();
1054  fe->get_JxW();
1055 
1056  fe->get_phi();
1057  if (dim > 1)
1058  {
1059  side_fe->get_JxW();
1060  side_fe->get_xyz();
1061  side_fe->get_phi();
1062  }
1063  if (dim > 2)
1064  {
1065  edge_fe->get_JxW();
1066  edge_fe->get_xyz();
1067  edge_fe->get_phi();
1068  }
1069 
1070  const FEContinuity cont = fe->get_continuity();
1071  this->conts[var] = cont;
1072  if (cont == C_ONE)
1073  {
1074  fe->get_dphi();
1075  if (dim > 1)
1076  side_fe->get_dphi();
1077  if (dim > 2)
1078  edge_fe->get_dphi();
1079  }
1080  }
1081  }
1082 
1083  // Now initialize any user requested shape functions, xyz vals, etc
1084  f.init_context(context);
1085 }
void get_side_fe(unsigned int var, FEGenericBase< OutputShape > *&fe) const
Accessor for edge/face (2D/3D) finite element object for variable var for the largest dimension in th...
Definition: fem_context.h:299
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
FEGenericBase< Real > FEBase
const std::vector< unsigned int > & variables
FEContinuity
defines an enum for finite element types to libmesh_assert a certain level (or type? Hcurl?) of continuity.
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
std::vector< FEContinuity > conts
const std::set< unsigned char > & elem_dimensions() const
Definition: fem_context.h:925

Member Function Documentation

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

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

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 
)

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

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 
)

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

◆ join()

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

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

Member Data Documentation

◆ action

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

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
protected

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
protected

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
protected

Definition at line 189 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

◆ 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

◆ projector

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

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
protected

Definition at line 198 of file generic_projector.h.


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