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

#include <generic_projector.h>

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

Public Member Functions

 ProjectInteriors (GenericProjector &p)
 
 ProjectInteriors (ProjectInteriors &p_i, Threads::split)
 
void operator() (const interior_range &range)
 
template<typename InsertInput , typename std::enable_if< std::is_same< typename ProjectionAction::InsertInput, InsertInput >::value, int >::type = 0>
void insert_id (dof_id_type id, const InsertInput &val, processor_id_type pid)
 
template<typename InsertInput , typename std::enable_if< !std::is_same< typename ProjectionAction::InsertInput, InsertInput >::value, int >::type = 0>
void insert_id (dof_id_type id, const InsertInput &val, processor_id_type pid)
 
template<typename InsertInput , typename std::enable_if< std::is_same< typename ProjectionAction::InsertInput, InsertInput >::value, int >::type = 0>
void insert_ids (const std::vector< dof_id_type > &ids, const std::vector< InsertInput > &vals, processor_id_type pid)
 
template<typename InsertInput , typename std::enable_if< !std::is_same< typename ProjectionAction::InsertInput, InsertInput >::value, int >::type = 0>
void insert_ids (const std::vector< dof_id_type > &ids, const std::vector< InsertInput > &vals, processor_id_type pid)
 
template<typename InsertInput , typename std::enable_if< std::is_same< typename ProjectionAction::InsertInput, InsertInput >::value, int >::type = 0>
void insert_id (dof_id_type id, const InsertInput &val, processor_id_type pid)
 
template<typename InsertInput , typename std::enable_if< std::is_same< typename ProjectionAction::InsertInput, InsertInput >::value, int >::type = 0>
void insert_ids (const std::vector< dof_id_type > &ids, const std::vector< InsertInput > &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
 
std::vector< FEFieldTypefield_types
 
GenericProjectorprojector
 
std::unordered_map< dof_id_type, std::pair< typename FFunctor::ValuePushType, processor_id_type > > new_ids_to_push
 
std::unordered_map< dof_id_type, typename FFunctor::ValuePushType > 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 FEGenericBase< typename FFunctor::RealType > &fe)
 

Detailed Description

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

Definition at line 369 of file generic_projector.h.

Constructor & Destructor Documentation

◆ ProjectInteriors() [1/2]

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

Definition at line 370 of file generic_projector.h.

370 : SubProjector(p) {}

◆ ProjectInteriors() [2/2]

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

Definition at line 372 of file generic_projector.h.

372 : SubProjector(p_i.projector) {}

Member Function Documentation

◆ construct_projection()

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 FEGenericBase< typename FFunctor::RealType > &  fe 
)
protectedinherited

Definition at line 2998 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::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::ids_to_save, libMesh::TensorTools::inner_product(), libMesh::DofObject::invalid_processor_id, libMesh::make_range(), libMesh::DofObject::processor_id(), libMesh::DenseVector< T >::size(), libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::system, and libMesh::System::time.

3003 {
3004  const auto & JxW = fe.get_JxW();
3005  const auto & phi = fe.get_phi();
3006  const std::vector<std::vector<typename FEGenericBase<typename FFunctor::RealType>::OutputGradient>> * dphi = nullptr;
3007  const std::vector<Point> & xyz_values = fe.get_xyz();
3008  const FEContinuity cont = fe.get_continuity();
3009  const std::unordered_map<dof_id_type, typename FFunctor::ValuePushType> & ids_to_save =
3010  this->projector.ids_to_save;
3011 
3012  if (cont == C_ONE)
3013  dphi = &(fe.get_dphi());
3014 
3015  const unsigned int n_involved_dofs =
3016  cast_int<unsigned int>(involved_dofs.size());
3017 
3018  std::vector<dof_id_type> free_dof_ids;
3019  DenseVector<typename FFunctor::ValuePushType> Uinvolved(n_involved_dofs);
3020  std::vector<char> dof_is_fixed(n_involved_dofs, false); // bools
3021 
3022  for (auto i : make_range(n_involved_dofs))
3023  {
3024  const dof_id_type id = dof_indices_var[involved_dofs[i]];
3025  auto iter = ids_to_save.find(id);
3026  if (iter == ids_to_save.end())
3027  free_dof_ids.push_back(id);
3028  else
3029  {
3030  dof_is_fixed[i] = true;
3031  Uinvolved(i) = iter->second;
3032  }
3033  }
3034 
3035  const unsigned int free_dofs = free_dof_ids.size();
3036 
3037  // There may be nothing to project
3038  if (!free_dofs)
3039  return;
3040 
3041  // The element matrix and RHS for projections.
3042  // Note that Ke is always real-valued, whereas
3043  // Fe may be complex valued if complex number
3044  // support is enabled
3045  DenseMatrix<Real> Ke(free_dofs, free_dofs);
3046  DenseVector<typename FFunctor::ValuePushType> Fe(free_dofs);
3047  // The new degree of freedom coefficients to solve for
3048  DenseVector<typename FFunctor::ValuePushType> Ufree(free_dofs);
3049 
3050  const unsigned int n_qp =
3051  cast_int<unsigned int>(xyz_values.size());
3052 
3053  // Loop over the quadrature points
3054  for (unsigned int qp=0; qp<n_qp; qp++)
3055  {
3056  // solution at the quadrature point
3057  FValue fineval = f.eval_at_point(context,
3058  var_component,
3059  xyz_values[qp],
3060  system.time,
3061  false);
3062  // solution grad at the quadrature point
3063  typename GFunctor::FunctorValue finegrad;
3064  if (cont == C_ONE)
3065  finegrad = g->eval_at_point(context,
3066  var_component,
3067  xyz_values[qp],
3068  system.time,
3069  false);
3070 
3071  // Form edge projection matrix
3072  for (unsigned int sidei=0, freei=0;
3073  sidei != n_involved_dofs; ++sidei)
3074  {
3075  unsigned int i = involved_dofs[sidei];
3076  // fixed DoFs aren't test functions
3077  if (dof_is_fixed[sidei])
3078  continue;
3079  for (unsigned int sidej=0, freej=0;
3080  sidej != n_involved_dofs; ++sidej)
3081  {
3082  unsigned int j = involved_dofs[sidej];
3083  if (dof_is_fixed[sidej])
3084  Fe(freei) -= phi[i][qp] * phi[j][qp] *
3085  JxW[qp] * Uinvolved(sidej);
3086  else
3087  Ke(freei,freej) += phi[i][qp] *
3088  phi[j][qp] * JxW[qp];
3089  if (cont == C_ONE)
3090  {
3091  if (dof_is_fixed[sidej])
3092  Fe(freei) -= ( TensorTools::inner_product((*dphi)[i][qp],
3093  (*dphi)[j][qp]) ) *
3094  JxW[qp] * Uinvolved(sidej);
3095  else
3096  Ke(freei,freej) += ( TensorTools::inner_product((*dphi)[i][qp],
3097  (*dphi)[j][qp]) )
3098  * JxW[qp];
3099  }
3100  if (!dof_is_fixed[sidej])
3101  freej++;
3102  }
3103  Fe(freei) += phi[i][qp] * fineval * JxW[qp];
3104  if (cont == C_ONE)
3105  Fe(freei) += (TensorTools::inner_product(finegrad,
3106  (*dphi)[i][qp]) ) *
3107  JxW[qp];
3108  freei++;
3109  }
3110  }
3111 
3112  Ke.cholesky_solve(Fe, Ufree);
3113 
3114  // Transfer new edge solutions to element
3115  const processor_id_type pid = node ?
3116  node->processor_id() : DofObject::invalid_processor_id;
3117  insert_ids(free_dof_ids, Ufree.get_values(), pid);
3118 }
Real time
For time-dependent problems, this is the time t at the beginning of the current timestep.
Definition: system.h:1615
uint8_t processor_id_type
std::unordered_map< dof_id_type, typename FFunctor::ValuePushType > ids_to_save
void insert_ids(const std::vector< dof_id_type > &ids, const std::vector< InsertInput > &vals, processor_id_type pid)
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:493
template class LIBMESH_EXPORT DenseMatrix< Real >
Definition: dense_matrix.C:35
FEContinuity
defines an enum for finite element types to libmesh_assert a certain level (or type? Hcurl?) of continuity.
IntRange< T > make_range(T beg, T end)
The 2-parameter make_range() helper function returns an IntRange<T> when both input parameters are of...
Definition: int_range.h:140
boostcopy::enable_if_c< ScalarTraits< T >::value &&ScalarTraits< T2 >::value, typename CompareTypes< T, T2 >::supertype >::type inner_product(const T &a, const T2 &b)
Definition: tensor_tools.h:51
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 2839 of file generic_projector.h.

References libMesh::Elem::active(), libMesh::Variable::active_on_subdomain(), libMesh::DofMap::dof_indices(), libMesh::System::get_dof_map(), libMesh::System::get_mesh(), libMesh::DofObject::id(), libMesh::DofObject::invalid_processor_id, libMesh::Parallel::Utils::is_sorted(), libMesh::libmesh_assert(), mesh, libMesh::Elem::node_ptr(), libMesh::ParallelObject::processor_id(), libMesh::DofObject::processor_id(), libMesh::Elem::subdomain_id(), libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::system, and libMesh::System::variable().

2840 {
2841  libmesh_assert (&node == elem.node_ptr(node_num));
2842 
2843  // Any ghosted node in our range might have an owner who needs our
2844  // data
2845  const processor_id_type owner = node.processor_id();
2846  if (owner != system.processor_id())
2847  {
2848  const MeshBase & mesh = system.get_mesh();
2849  const DofMap & dof_map = system.get_dof_map();
2850 
2851  // But let's check and see if we can be certain the owner can
2852  // compute any or all of its own dof coefficients on that node.
2853  std::vector<dof_id_type> node_dof_ids, patch_dof_ids;
2854  for (const auto & var : vars)
2855  {
2856  const Variable & variable = system.variable(var);
2857 
2858  if (!variable.active_on_subdomain(elem.subdomain_id()))
2859  continue;
2860 
2861  dof_map.dof_indices(elem, node_num, node_dof_ids, var);
2862  }
2863  libmesh_assert(std::is_sorted(node_dof_ids.begin(),
2864  node_dof_ids.end()));
2865  const std::vector<dof_id_type> & patch =
2866  (*this->projector.nodes_to_elem)[node.id()];
2867  for (const auto & elem_id : patch)
2868  {
2869  const Elem & patch_elem = mesh.elem_ref(elem_id);
2870  if (!patch_elem.active() || owner != patch_elem.processor_id())
2871  continue;
2872  dof_map.dof_indices(&patch_elem, patch_dof_ids);
2873  std::sort(patch_dof_ids.begin(), patch_dof_ids.end());
2874 
2875  // Remove any node_dof_ids that we see can be calculated on
2876  // this element
2877  std::vector<dof_id_type> diff_ids(node_dof_ids.size());
2878  auto it = std::set_difference(node_dof_ids.begin(), node_dof_ids.end(),
2879  patch_dof_ids.begin(), patch_dof_ids.end(), diff_ids.begin());
2880  diff_ids.resize(it-diff_ids.begin());
2881  node_dof_ids.swap(diff_ids);
2882  if (node_dof_ids.empty())
2883  break;
2884  }
2885 
2886  // Give ids_to_push default invalid pid: not yet computed
2887  for (auto id : node_dof_ids)
2889  }
2890 }
std::unordered_map< dof_id_type, std::pair< typename FFunctor::ValuePushType, processor_id_type > > new_ids_to_push
const Variable & variable(unsigned int var) const
Return a constant reference to Variable var.
Definition: system.h:2458
NodesToElemMap * nodes_to_elem
MeshBase & mesh
const MeshBase & get_mesh() const
Definition: system.h:2358
uint8_t processor_id_type
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:493
libmesh_assert(ctx)
processor_id_type processor_id() const
const DofMap & get_dof_map() const
Definition: system.h:2374

◆ insert_id() [1/3]

template<typename FFunctor , typename GFunctor , typename FValue , typename ProjectionAction >
template<typename InsertInput , typename std::enable_if< std::is_same< typename ProjectionAction::InsertInput, InsertInput >::value, int >::type >
void libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::SubFunctor::insert_id ( dof_id_type  id,
const InsertInput &  val,
processor_id_type  pid 
)
inherited

Definition at line 1378 of file generic_projector.h.

1380 {
1381  libmesh_error_msg("ID insertion should only occur when the provided input matches that "
1382  "expected by the ProjectionAction");
1383 }

◆ insert_id() [2/3]

template<typename FFunctor , typename GFunctor , typename FValue , typename ProjectionAction >
template<typename InsertInput , typename std::enable_if< !std::is_same< typename ProjectionAction::InsertInput, InsertInput >::value, int >::type = 0>
void libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::SubFunctor::insert_id ( typename InsertInput  ,
typename std::enable_if< !std::is_same< typename ProjectionAction::InsertInput, InsertInput >::value, int >::type  = 0 
)

◆ insert_id() [3/3]

template<typename FFunctor , typename GFunctor , typename FValue , typename ProjectionAction >
template<typename InsertInput , typename std::enable_if< std::is_same< typename ProjectionAction::InsertInput, InsertInput >::value, int >::type = 0>
void libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::SubFunctor::insert_id ( typename InsertInput  ,
typename std::enable_if< std::is_same< typename ProjectionAction::InsertInput, InsertInput >::value, int >::type  = 0 
)

Definition at line 1378 of file generic_projector.h.

1380 {
1381  libmesh_error_msg("ID insertion should only occur when the provided input matches that "
1382  "expected by the ProjectionAction");
1383 }

◆ insert_ids() [1/3]

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

Definition at line 1420 of file generic_projector.h.

1424 {
1425  libmesh_error_msg("ID insertion should only occur when the provided input matches that "
1426  "expected by the ProjectionAction");
1427 }

◆ insert_ids() [2/3]

template<typename FFunctor , typename GFunctor , typename FValue , typename ProjectionAction >
template<typename InsertInput , typename std::enable_if< !std::is_same< typename ProjectionAction::InsertInput, InsertInput >::value, int >::type = 0>
void libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::SubFunctor::insert_ids ( typename InsertInput  ,
typename std::enable_if< !std::is_same< typename ProjectionAction::InsertInput, InsertInput >::value, int >::type  = 0 
)

◆ insert_ids() [3/3]

template<typename FFunctor , typename GFunctor , typename FValue , typename ProjectionAction >
template<typename InsertInput , typename std::enable_if< std::is_same< typename ProjectionAction::InsertInput, InsertInput >::value, int >::type = 0>
void libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::SubFunctor::insert_ids ( typename InsertInput  ,
typename std::enable_if< std::is_same< typename ProjectionAction::InsertInput, InsertInput >::value, int >::type  = 0 
)

Definition at line 1420 of file generic_projector.h.

1424 {
1425  libmesh_error_msg("ID insertion should only occur when the provided input matches that "
1426  "expected by the ProjectionAction");
1427 }

◆ 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 1471 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.

Referenced by libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::SortAndCopy::join().

1472 {
1473  new_ids_to_push.insert(other.new_ids_to_push.begin(), other.new_ids_to_push.end());
1474  new_ids_to_save.insert(other.new_ids_to_save.begin(), other.new_ids_to_save.end());
1475 }
std::unordered_map< dof_id_type, std::pair< typename FFunctor::ValuePushType, processor_id_type > > new_ids_to_push
std::unordered_map< dof_id_type, typename FFunctor::ValuePushType > new_ids_to_save

◆ operator()()

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

Definition at line 2693 of file generic_projector.h.

References libMesh::Variable::active_on_subdomain(), libMesh::FEType::default_quadrature_rule(), dim, libMesh::DofObject::dof_number(), libMesh::FEType::family, libMesh::System::get_dof_map(), libMesh::FEMap::inverse_map(), libMesh::Utility::iota(), libMesh::Elem::JUST_COARSENED, libMesh::LAGRANGE, libMesh::LAGRANGE_VEC, n_nodes, libMesh::System::number(), libMesh::FEType::order, libMesh::DofObject::processor_id(), libMesh::SCALAR, libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::system, libMesh::System::time, libMesh::Variable::type(), libMesh::DofMap::var_group_from_var_number(), libMesh::System::variable(), and libMesh::System::variable_scalar_number().

2694 {
2695  LOG_SCOPE ("project_interiors","GenericProjector");
2696 
2697  const unsigned int sys_num = system.number();
2698 
2699  // Iterate over all dof-bearing element interiors in the range
2700  for (const auto & elem : range)
2701  {
2702  unsigned char dim = cast_int<unsigned char>(elem->dim());
2703 
2704  context.pre_fe_reinit(system, elem);
2705 
2706  // Loop over all the variables we've been requested to project, to
2707  // do the projection
2708  for (const auto & var : this->projector.variables)
2709  {
2710  const Variable & variable = system.variable(var);
2711 
2712  if (!variable.active_on_subdomain(elem->subdomain_id()))
2713  continue;
2714 
2715  const FEType & base_fe_type = variable.type();
2716 
2717  if (base_fe_type.family == SCALAR)
2718  continue;
2719 
2720  FEGenericBase<typename FFunctor::RealType> * fe = nullptr;
2721  context.get_element_fe( var, fe, dim );
2722 
2723  FEType fe_type = base_fe_type;
2724  const auto & dof_map = system.get_dof_map();
2725  const auto vg = dof_map.var_group_from_var_number(var);
2726  const bool add_p_level = dof_map.should_p_refine(vg);
2727 
2728  // This may be a p refined element
2729  fe_type.order = fe_type.order + add_p_level * elem->p_level();
2730 
2731  const unsigned int var_component =
2733 
2734  // If this is a Lagrange element with interior DoFs then by
2735  // convention we interpolate at the interior node rather
2736  // than project along the interior.
2737  if (fe_type.family == LAGRANGE || fe_type.family == LAGRANGE_VEC)
2738  {
2739  if (fe_type.order > 1)
2740  {
2741  const unsigned int first_interior_node =
2742  (elem->n_vertices() +
2743  ((elem->dim() > 2) * elem->n_edges()) +
2744  ((elem->dim() > 1) * elem->n_sides()));
2745  const unsigned int n_nodes = elem->n_nodes();
2746 
2747  // < vs != is important here for HEX20, QUAD8!
2748  for (unsigned int n = first_interior_node; n < n_nodes; ++n)
2749  {
2750  const Node & interior_node = elem->node_ref(n);
2751 
2752  FValue fval = f.eval_at_point
2753  (context, var_component, interior_node,
2754  system.time, false);
2755  const processor_id_type pid =
2756  interior_node.processor_id();
2757 
2758  if (fe_type.family == LAGRANGE_VEC)
2759  {
2760  // We will have a number of nodal value DoFs equal to the elem dim
2761  for (unsigned int i = 0; i < elem->dim(); ++i)
2762  {
2763  const dof_id_type dof_id = interior_node.dof_number(sys_num, var, i);
2764 
2765  // Need this conversion so that this method
2766  // will compile for TYPE_SCALAR instantiations
2767  const auto insert_val =
2768  raw_value<typename ProjectionAction::InsertInput>(fval, i);
2769 
2770  insert_id(dof_id, insert_val, pid);
2771  }
2772  }
2773  else // We are LAGRANGE
2774  {
2775  const dof_id_type dof_id =
2776  interior_node.dof_number(sys_num, var, 0);
2777  insert_id(dof_id, fval, pid);
2778  }
2779  }
2780  }
2781  continue;
2782  }
2783 
2784 #ifdef LIBMESH_ENABLE_AMR
2785  if (elem->refinement_flag() == Elem::JUST_COARSENED)
2786  {
2787  std::vector<Point> fine_points;
2788 
2789  std::unique_ptr<FEGenericBase<typename FFunctor::RealType>> fine_fe
2791 
2792  std::unique_ptr<QBase> qrule
2793  (base_fe_type.default_quadrature_rule(dim));
2794  fine_fe->attach_quadrature_rule(qrule.get());
2795 
2796  const std::vector<Point> & child_xyz =
2797  fine_fe->get_xyz();
2798 
2799  for (auto & child : elem->child_ref_range())
2800  {
2801  fine_fe->reinit(&child);
2802  fine_points.insert(fine_points.end(),
2803  child_xyz.begin(),
2804  child_xyz.end());
2805  }
2806 
2807  std::vector<Point> fine_qp;
2808  FEMap::inverse_map (dim, elem, fine_points, fine_qp);
2809 
2810  context.elem_fe_reinit(&fine_qp);
2811  }
2812  else
2813 #endif // LIBMESH_ENABLE_AMR
2815 
2816  const std::vector<dof_id_type> & dof_indices =
2817  context.get_dof_indices(var);
2818 
2819  const unsigned int n_dofs =
2820  cast_int<unsigned int>(dof_indices.size());
2821 
2822  std::vector<unsigned int> all_dofs(n_dofs);
2823  std::iota(all_dofs.begin(), all_dofs.end(), 0);
2824 
2825  this->construct_projection
2826  (dof_indices, all_dofs, var_component,
2827  nullptr, *fe);
2828  } // end variables loop
2829  } // end elements loop
2830 }
Real time
For time-dependent problems, this is the time t at the beginning of the current timestep.
Definition: system.h:1615
unsigned int variable_scalar_number(std::string_view var, unsigned int component) const
Definition: system.h:2489
const Variable & variable(unsigned int var) const
Return a constant reference to Variable var.
Definition: system.h:2458
virtual void pre_fe_reinit(const System &, const Elem *e)
Reinitializes local data vectors/matrices on the current geometric element.
Definition: fem_context.C:1683
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 FEGenericBase< typename FFunctor::RealType > &fe)
unsigned int dim
static Point inverse_map(const unsigned int dim, const Elem *elem, const Point &p, const Real tolerance=TOLERANCE, const bool secure=true, const bool extra_checks=true)
Definition: fe_map.C:1628
void insert_id(dof_id_type id, const InsertInput &val, processor_id_type pid)
void iota(ForwardIter first, ForwardIter last, T value)
Utility::iota was created back when std::iota was just an SGI STL extension.
Definition: utility.h:229
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:1477
unsigned int var_group_from_var_number(unsigned int var_num) const
Definition: dof_map.h:2430
uint8_t processor_id_type
const dof_id_type n_nodes
Definition: tecplot_io.C:67
unsigned int number() const
Definition: system.h:2350
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:363
const std::vector< unsigned int > & variables
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:277
const DofMap & get_dof_map() const
Definition: system.h:2374
uint8_t dof_id_type
Definition: id_types.h:67
const FEType & type() const
Definition: variable.h:144

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 231 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 236 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 239 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 232 of file generic_projector.h.

◆ field_types

template<typename FFunctor , typename GFunctor , typename FValue , typename ProjectionAction >
std::vector<FEFieldType> libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::SubFunctor::field_types

Definition at line 240 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 267 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<typename FFunctor::ValuePushType, 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, typename FFunctor::ValuePushType> 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 170 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 242 of file generic_projector.h.


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