libMesh::OldSolutionCoefs< Output, point_output > Class Template Reference
Inheritance diagram for libMesh::OldSolutionCoefs< Output, point_output >:

Public Types

typedef DSNAOutput< Output >::type DSNA
 

Public Member Functions

 OldSolutionCoefs (const libMesh::System &sys_in)
 
 OldSolutionCoefs (const OldSolutionCoefs &in)
 
DSNA eval_at_node (const FEMContext &c, unsigned int i, unsigned int elem_dim, const Node &n, Real=0.)
 
DSNA eval_at_point (const FEMContext &c, unsigned int i, const Point &p, Real=0.)
 
void eval_old_dofs (const FEMContext &c, unsigned int var, std::vector< DSNA > &values)
 
template<>
DynamicSparseNumberArray< Real, dof_id_typeeval_at_point (const FEMContext &c, unsigned int i, const Point &p, Real)
 
template<>
VectorValue< DynamicSparseNumberArray< Real, dof_id_type > > eval_at_point (const FEMContext &c, unsigned int i, const Point &p, Real)
 
template<>
DynamicSparseNumberArray< Real, dof_id_typeeval_at_node (const FEMContext &c, unsigned int i, unsigned int, const Node &n, Real)
 
template<>
VectorValue< DynamicSparseNumberArray< Real, dof_id_type > > eval_at_node (const FEMContext &c, unsigned int i, unsigned int elem_dim, const Node &n, Real)
 
template<>
void get_shape_outputs (FEBase &fe)
 
template<>
void get_shape_outputs (FEBase &fe)
 
template<>
void get_shape_outputs (FEBase &fe)
 
template<>
void get_shape_outputs (FEBase &fe)
 
void init_context (FEMContext &c)
 
bool is_grid_projection ()
 

Static Public Member Functions

static void get_shape_outputs (FEBase &fe)
 

Protected Member Functions

template<>
const Real out_of_elem_tol
 
template<>
const Real out_of_elem_tol
 
template<>
const Real out_of_elem_tol
 
template<>
const Real out_of_elem_tol
 
void check_old_context (const FEMContext &c)
 
bool check_old_context (const FEMContext &c, const Point &p)
 

Protected Attributes

const Elemlast_elem
 
const Systemsys
 
FEMContext old_context
 
std::vector< unsigned int > component_to_var
 

Static Protected Attributes

static const Real out_of_elem_tol
 

Detailed Description

template<typename Output, void(FEMContext::*)(unsigned int, const Point &, Output &, const Real) const point_output>
class libMesh::OldSolutionCoefs< Output, point_output >

The OldSolutionCoefs input functor class can be used with GenericProjector to read solution transfer coefficients on a just-refined-and-coarsened mesh.

Author
Roy H. Stogner
Date
2017

Definition at line 447 of file system_projection.C.

Member Typedef Documentation

◆ DSNA

template<typename Output , void(FEMContext::*)(unsigned int, const Point &, Output &, const Real) const point_output>
typedef DSNAOutput<Output>::type libMesh::OldSolutionCoefs< Output, point_output >::DSNA

Definition at line 450 of file system_projection.C.

Constructor & Destructor Documentation

◆ OldSolutionCoefs() [1/2]

template<typename Output , void(FEMContext::*)(unsigned int, const Point &, Output &, const Real) const point_output>
libMesh::OldSolutionCoefs< Output, point_output >::OldSolutionCoefs ( const libMesh::System sys_in)
inline

Definition at line 452 of file system_projection.C.

452  :
453  OldSolutionBase<Output, point_output>(sys_in)
454  {
456  }
void set_algebraic_type(const AlgebraicType atype)
Definition: fem_context.h:948

◆ OldSolutionCoefs() [2/2]

template<typename Output , void(FEMContext::*)(unsigned int, const Point &, Output &, const Real) const point_output>
libMesh::OldSolutionCoefs< Output, point_output >::OldSolutionCoefs ( const OldSolutionCoefs< Output, point_output > &  in)
inline

Definition at line 458 of file system_projection.C.

458  :
459  OldSolutionBase<Output, point_output>(in.sys)
460  {
462  }
void set_algebraic_type(const AlgebraicType atype)
Definition: fem_context.h:948

Member Function Documentation

◆ check_old_context() [1/2]

template<typename Output , void(FEMContext::*)(unsigned int, const Point &, Output &, const Real) const point_output>
void libMesh::OldSolutionBase< Output, point_output >::check_old_context ( const FEMContext c)
inlineprotectedinherited

Definition at line 259 of file generic_projector.h.

References libMesh::FEMContext::get_elem(), libMesh::Elem::JUST_COARSENED, libMesh::Elem::JUST_REFINED, libMesh::DofObject::old_dof_object, libMesh::Elem::parent(), and libMesh::Elem::refinement_flag().

260  {
261  LOG_SCOPE ("check_old_context(c)", "OldSolutionBase");
262  const Elem & elem = c.get_elem();
263  if (last_elem != &elem)
264  {
265  if (elem.refinement_flag() == Elem::JUST_REFINED)
266  {
267  old_context.pre_fe_reinit(sys, elem.parent());
268  }
269  else if (elem.refinement_flag() == Elem::JUST_COARSENED)
270  {
271  libmesh_error();
272  }
273  else
274  {
275  if (!elem.old_dof_object)
276  {
277  libmesh_error();
278  }
279 
280  old_context.pre_fe_reinit(sys, &elem);
281  }
282 
283  last_elem = &elem;
284  }
285  else
286  {
287  libmesh_assert(old_context.has_elem());
288  }
289  }
virtual void pre_fe_reinit(const System &, const Elem *e)
Definition: fem_context.C:1552
bool has_elem() const
Definition: fem_context.h:865

◆ check_old_context() [2/2]

template<typename Output , void(FEMContext::*)(unsigned int, const Point &, Output &, const Real) const point_output>
bool libMesh::OldSolutionBase< Output, point_output >::check_old_context ( const FEMContext c,
const Point p 
)
inlineprotectedinherited

Definition at line 292 of file generic_projector.h.

References libMesh::Elem::child_ref_range(), libMesh::FEMContext::get_elem(), libMesh::Elem::hmax(), libMesh::Elem::JUST_COARSENED, libMesh::Elem::JUST_REFINED, libMesh::DofObject::old_dof_object, libMesh::Elem::parent(), libMesh::Real, and libMesh::Elem::refinement_flag().

293  {
294  LOG_SCOPE ("check_old_context(c,p)", "OldSolutionBase");
295  const Elem & elem = c.get_elem();
296  if (last_elem != &elem)
297  {
298  if (elem.refinement_flag() == Elem::JUST_REFINED)
299  {
300  old_context.pre_fe_reinit(sys, elem.parent());
301  }
302  else if (elem.refinement_flag() == Elem::JUST_COARSENED)
303  {
304  // Find the child with this point. Use out_of_elem_tol
305  // (in physical space, which may correspond to a large
306  // tolerance in master space!) to allow for out-of-element
307  // finite differencing of mixed gradient terms. Pray we
308  // have no quadrature locations which are within 1e-5 of
309  // the element subdivision boundary but are not exactly on
310  // that boundary.
311  const Real master_tol = out_of_elem_tol / elem.hmax() * 2;
312 
313  for (auto & child : elem.child_ref_range())
314  if (child.close_to_point(p, master_tol))
315  {
316  old_context.pre_fe_reinit(sys, &child);
317  break;
318  }
319 
320  libmesh_assert
321  (old_context.get_elem().close_to_point(p, master_tol));
322  }
323  else
324  {
325  if (!elem.old_dof_object)
326  return false;
327 
328  old_context.pre_fe_reinit(sys, &elem);
329  }
330 
331  last_elem = &elem;
332  }
333  else
334  {
335  libmesh_assert(old_context.has_elem());
336 
337  const Real master_tol = out_of_elem_tol / elem.hmax() * 2;
338 
339  if (!old_context.get_elem().close_to_point(p, master_tol))
340  {
341  libmesh_assert_equal_to
342  (elem.refinement_flag(), Elem::JUST_COARSENED);
343 
344  for (auto & child : elem.child_ref_range())
345  if (child.close_to_point(p, master_tol))
346  {
347  old_context.pre_fe_reinit(sys, &child);
348  break;
349  }
350 
351  libmesh_assert
352  (old_context.get_elem().close_to_point(p, master_tol));
353  }
354  }
355 
356  return true;
357  }
virtual void pre_fe_reinit(const System &, const Elem *e)
Definition: fem_context.C:1552
const Elem & get_elem() const
Definition: fem_context.h:871
virtual bool close_to_point(const Point &p, Real tol) const
Definition: elem.C:2326
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
bool has_elem() const
Definition: fem_context.h:865

◆ eval_at_node() [1/3]

template<typename Output , void(FEMContext::*)(unsigned int, const Point &, Output &, const Real) const point_output>
DSNA libMesh::OldSolutionCoefs< Output, point_output >::eval_at_node ( const FEMContext c,
unsigned int  i,
unsigned int  elem_dim,
const Node n,
Real  = 0. 
)

◆ eval_at_node() [2/3]

template<>
DynamicSparseNumberArray< Real, dof_id_type > libMesh::OldSolutionCoefs< Real, &FEMContext::point_value >::eval_at_node ( const FEMContext c,
unsigned int  i,
unsigned  int,
const Node n,
Real   
)
inline

Definition at line 595 of file system_projection.C.

References libMesh::DofObject::dof_number(), libMesh::DofObject::n_comp(), libMesh::DofObject::n_vars(), and libMesh::DofObject::old_dof_object.

600 {
601  LOG_SCOPE ("Real eval_at_node()", "OldSolutionCoefs");
602 
603  // Optimize for the common case, where this node was part of the
604  // old solution.
605  //
606  // Be sure to handle cases where the variable wasn't defined on
607  // this node (due to changing subdomain support) or where the
608  // variable has no components on this node (due to Elem order
609  // exceeding FE order)
610  if (n.old_dof_object &&
611  n.old_dof_object->n_vars(sys.number()) &&
612  n.old_dof_object->n_comp(sys.number(), i))
613  {
614  DynamicSparseNumberArray<Real, dof_id_type> returnval;
615  const dof_id_type old_id =
616  n.old_dof_object->dof_number(sys.number(), i, 0);
617  returnval.resize(1);
618  returnval.raw_at(0) = 1;
619  returnval.raw_index(0) = old_id;
620  return returnval;
621  }
622 
623  return this->eval_at_point(c, i, n, 0);
624 }
DSNA eval_at_point(const FEMContext &c, unsigned int i, const Point &p, Real=0.)
unsigned int number() const
Definition: system.h:2025
uint8_t dof_id_type
Definition: id_types.h:64

◆ eval_at_node() [3/3]

template<>
VectorValue< DynamicSparseNumberArray< Real, dof_id_type > > libMesh::OldSolutionCoefs< RealGradient, &FEMContext::point_gradient >::eval_at_node ( const FEMContext c,
unsigned int  i,
unsigned int  elem_dim,
const Node n,
Real   
)
inline

Definition at line 632 of file system_projection.C.

References libMesh::DofObject::dof_number(), libMesh::DofObject::n_comp(), libMesh::DofObject::n_vars(), and libMesh::DofObject::old_dof_object.

637 {
638  LOG_SCOPE ("RealGradient eval_at_node()", "OldSolutionCoefs");
639 
640  // Optimize for the common case, where this node was part of the
641  // old solution.
642  //
643  // Be sure to handle cases where the variable wasn't defined on
644  // this node (due to changing subdomain support) or where the
645  // variable has no components on this node (due to Elem order
646  // exceeding FE order)
647  if (n.old_dof_object &&
648  n.old_dof_object->n_vars(sys.number()) &&
649  n.old_dof_object->n_comp(sys.number(), i))
650  {
651  VectorValue<DynamicSparseNumberArray<Real, dof_id_type> > g;
652  for (unsigned int d = 0; d != elem_dim; ++d)
653  {
654  const dof_id_type old_id =
655  n.old_dof_object->dof_number(sys.number(), i, d+1);
656  g(d).resize(1);
657  g(d).raw_at(0) = 1;
658  g(d).raw_index(0) = old_id;
659  }
660  return g;
661  }
662 
663  return this->eval_at_point(c, i, n, 0);
664 }
DSNA eval_at_point(const FEMContext &c, unsigned int i, const Point &p, Real=0.)
unsigned int number() const
Definition: system.h:2025
uint8_t dof_id_type
Definition: id_types.h:64

◆ eval_at_point() [1/3]

template<typename Output , void(FEMContext::*)(unsigned int, const Point &, Output &, const Real) const point_output>
DSNA libMesh::OldSolutionCoefs< Output, point_output >::eval_at_point ( const FEMContext c,
unsigned int  i,
const Point p,
Real  = 0. 
)

◆ eval_at_point() [2/3]

template<>
DynamicSparseNumberArray< Real, dof_id_type > libMesh::OldSolutionCoefs< Real, &FEMContext::point_value >::eval_at_point ( const FEMContext c,
unsigned int  i,
const Point p,
Real   
)
inline

Definition at line 503 of file system_projection.C.

References libMesh::FEGenericBase< OutputType >::get_phi(), libMesh::index_range(), and libMesh::Real.

507 {
508  LOG_SCOPE ("eval_at_point()", "OldSolutionCoefs");
509 
510  if (!this->check_old_context(c, p))
511  return 0;
512 
513  // Get finite element object
514  FEGenericBase<Real> * fe = nullptr;
516  (i, fe, this->old_context.get_elem_dim());
517 
518  // Build a FE for calculating phi(p)
519  FEGenericBase<Real> * fe_new =
520  this->old_context.build_new_fe(fe, p);
521 
522  // Get the values and global indices of the shape functions
523  const std::vector<std::vector<Real> > & phi = fe_new->get_phi();
524  const std::vector<dof_id_type> & dof_indices =
525  this->old_context.get_dof_indices(i);
526 
527  const std::size_t n_dofs = phi.size();
528  libmesh_assert_equal_to(n_dofs, dof_indices.size());
529 
530  DynamicSparseNumberArray<Real, dof_id_type> returnval;
531  returnval.resize(n_dofs);
532 
533  for (auto j : index_range(phi))
534  {
535  returnval.raw_at(j) = phi[j][0];
536  returnval.raw_index(j) = dof_indices[j];
537  }
538 
539  return returnval;
540 }
void check_old_context(const FEMContext &c)
unsigned char get_elem_dim() const
Definition: fem_context.h:906
IntRange< std::size_t > index_range(const std::vector< T > &vec)
Definition: int_range.h:104
const std::vector< dof_id_type > & get_dof_indices() const
Definition: diff_context.h:367
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void get_element_fe(unsigned int var, FEGenericBase< OutputShape > *&fe) const
Definition: fem_context.h:262
FEGenericBase< OutputShape > * build_new_fe(const FEGenericBase< OutputShape > *fe, const Point &p, const Real tolerance=TOLERANCE) const
Definition: fem_context.C:1862

◆ eval_at_point() [3/3]

template<>
VectorValue< DynamicSparseNumberArray< Real, dof_id_type > > libMesh::OldSolutionCoefs< RealGradient, &FEMContext::point_gradient >::eval_at_point ( const FEMContext c,
unsigned int  i,
const Point p,
Real   
)
inline

Definition at line 548 of file system_projection.C.

References libMesh::FEGenericBase< OutputType >::get_dphi(), libMesh::index_range(), and libMesh::Real.

552 {
553  LOG_SCOPE ("eval_at_point()", "OldSolutionCoefs");
554 
555  if (!this->check_old_context(c, p))
556  return 0;
557 
558  // Get finite element object
559  FEGenericBase<Real> * fe = nullptr;
561  (i, fe, this->old_context.get_elem_dim());
562 
563  // Build a FE for calculating phi(p)
564  FEGenericBase<Real> * fe_new =
565  this->old_context.build_new_fe(fe, p);
566 
567  // Get the values and global indices of the shape functions
568  const std::vector<std::vector<RealGradient> > & dphi = fe_new->get_dphi();
569  const std::vector<dof_id_type> & dof_indices =
570  this->old_context.get_dof_indices(i);
571 
572  const std::size_t n_dofs = dphi.size();
573  libmesh_assert_equal_to(n_dofs, dof_indices.size());
574 
575  VectorValue<DynamicSparseNumberArray<Real, dof_id_type> > returnval;
576 
577  for (unsigned int d = 0; d != LIBMESH_DIM; ++d)
578  returnval(d).resize(n_dofs);
579 
580  for (auto j : index_range(dphi))
581  for (int d = 0; d != LIBMESH_DIM; ++d)
582  {
583  returnval(d).raw_at(j) = dphi[j][0](d);
584  returnval(d).raw_index(j) = dof_indices[j];
585  }
586 
587  return returnval;
588 }
void check_old_context(const FEMContext &c)
unsigned char get_elem_dim() const
Definition: fem_context.h:906
IntRange< std::size_t > index_range(const std::vector< T > &vec)
Definition: int_range.h:104
const std::vector< dof_id_type > & get_dof_indices() const
Definition: diff_context.h:367
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void get_element_fe(unsigned int var, FEGenericBase< OutputShape > *&fe) const
Definition: fem_context.h:262
FEGenericBase< OutputShape > * build_new_fe(const FEGenericBase< OutputShape > *fe, const Point &p, const Real tolerance=TOLERANCE) const
Definition: fem_context.C:1862

◆ eval_old_dofs()

template<typename Output , void(FEMContext::*)(unsigned int, const Point &, Output &, const Real) const point_output>
void libMesh::OldSolutionCoefs< Output, point_output >::eval_old_dofs ( const FEMContext c,
unsigned int  var,
std::vector< DSNA > &  values 
)
inline

Definition at line 475 of file system_projection.C.

References libMesh::index_range().

478  {
479  LOG_SCOPE ("eval_old_dofs()", "OldSolutionValue");
480 
481  this->check_old_context(c);
482 
483  const std::vector<dof_id_type> & old_dof_indices =
484  this->old_context.get_dof_indices(var);
485 
486  libmesh_assert_equal_to (old_dof_indices.size(), values.size());
487 
488  for (auto i : index_range(values))
489  {
490  values[i].resize(1);
491  values[i].raw_at(0) = 1;
492  values[i].raw_index(0) = old_dof_indices[i];
493  }
494  }
void check_old_context(const FEMContext &c)
IntRange< std::size_t > index_range(const std::vector< T > &vec)
Definition: int_range.h:104
const std::vector< dof_id_type > & get_dof_indices() const
Definition: diff_context.h:367

◆ get_shape_outputs() [1/5]

template<typename Output , void(FEMContext::*)(unsigned int, const Point &, Output &, const Real) const point_output>
static void libMesh::OldSolutionBase< Output, point_output >::get_shape_outputs ( FEBase fe)
staticinherited

◆ get_shape_outputs() [2/5]

template<>
void libMesh::OldSolutionBase< Number, &FEMContext::point_value >::get_shape_outputs ( FEBase fe)
inlineinherited

Definition at line 451 of file generic_projector.h.

References libMesh::FEGenericBase< OutputType >::get_phi().

452 {
453  fe.get_phi();
454 }

◆ get_shape_outputs() [3/5]

template<>
void libMesh::OldSolutionBase< Gradient, &FEMContext::point_gradient >::get_shape_outputs ( FEBase fe)
inlineinherited

Definition at line 459 of file generic_projector.h.

References libMesh::FEGenericBase< OutputType >::get_dphi().

460 {
461  fe.get_dphi();
462 }

◆ get_shape_outputs() [4/5]

template<>
void libMesh::OldSolutionBase< Real, &FEMContext::point_value >::get_shape_outputs ( FEBase fe)
inlineinherited

Definition at line 468 of file generic_projector.h.

References libMesh::FEGenericBase< OutputType >::get_phi().

469 {
470  fe.get_phi();
471 }

◆ get_shape_outputs() [5/5]

template<>
void libMesh::OldSolutionBase< RealGradient, &FEMContext::point_gradient >::get_shape_outputs ( FEBase fe)
inlineinherited

Definition at line 476 of file generic_projector.h.

References libMesh::FEGenericBase< OutputType >::get_dphi().

477 {
478  fe.get_dphi();
479 }

◆ init_context()

template<typename Output , void(FEMContext::*)(unsigned int, const Point &, Output &, const Real) const point_output>
void libMesh::OldSolutionBase< Output, point_output >::init_context ( FEMContext c)
inlineinherited

Definition at line 237 of file generic_projector.h.

References libMesh::FEMContext::DOFS_ONLY, and libMesh::FEMContext::set_algebraic_type().

238  {
239  c.set_algebraic_type(FEMContext::DOFS_ONLY);
240 
241  // Loop over variables, to prerequest
242  for (unsigned int var=0; var!=sys.n_vars(); ++var)
243  {
244  FEBase * fe = nullptr;
245  const std::set<unsigned char> & elem_dims =
247 
248  for (const auto & dim : elem_dims)
249  {
250  old_context.get_element_fe(var, fe, dim);
251  get_shape_outputs(*fe);
252  }
253  }
254  }
static void get_shape_outputs(FEBase &fe)
FEGenericBase< Real > FEBase
void get_element_fe(unsigned int var, FEGenericBase< OutputShape > *&fe) const
Definition: fem_context.h:262
unsigned int n_vars() const
Definition: system.h:2105
const std::set< unsigned char > & elem_dimensions() const
Definition: fem_context.h:913

◆ is_grid_projection()

template<typename Output , void(FEMContext::*)(unsigned int, const Point &, Output &, const Real) const point_output>
bool libMesh::OldSolutionBase< Output, point_output >::is_grid_projection ( )
inlineinherited

Definition at line 256 of file generic_projector.h.

256 { return true; }

◆ out_of_elem_tol() [1/4]

template<>
const Real libMesh::OldSolutionBase< Number, &FEMContext::point_value >::out_of_elem_tol ( )
protectedinherited

Definition at line 565 of file generic_projector.h.

◆ out_of_elem_tol() [2/4]

template<>
const Real libMesh::OldSolutionBase< Gradient, &FEMContext::point_gradient >::out_of_elem_tol ( )
protectedinherited

Definition at line 568 of file generic_projector.h.

◆ out_of_elem_tol() [3/4]

template<>
const Real libMesh::OldSolutionBase< Real, &FEMContext::point_value >::out_of_elem_tol ( )
protectedinherited

Definition at line 572 of file generic_projector.h.

◆ out_of_elem_tol() [4/4]

template<>
const Real libMesh::OldSolutionBase< RealGradient, &FEMContext::point_gradient >::out_of_elem_tol ( )
protectedinherited

Definition at line 575 of file generic_projector.h.

Member Data Documentation

◆ component_to_var

template<typename Output , void(FEMContext::*)(unsigned int, const Point &, Output &, const Real) const point_output>
std::vector<unsigned int> libMesh::OldSolutionBase< Output, point_output >::component_to_var
protectedinherited

Definition at line 363 of file generic_projector.h.

◆ last_elem

template<typename Output , void(FEMContext::*)(unsigned int, const Point &, Output &, const Real) const point_output>
const Elem* libMesh::OldSolutionBase< Output, point_output >::last_elem
protectedinherited

Definition at line 360 of file generic_projector.h.

◆ old_context

template<typename Output , void(FEMContext::*)(unsigned int, const Point &, Output &, const Real) const point_output>
FEMContext libMesh::OldSolutionBase< Output, point_output >::old_context
protectedinherited

Definition at line 362 of file generic_projector.h.

◆ out_of_elem_tol

template<typename Output , void(FEMContext::*)(unsigned int, const Point &, Output &, const Real) const point_output>
const Real libMesh::OldSolutionBase< Output, point_output >::out_of_elem_tol
staticprotectedinherited

Definition at line 365 of file generic_projector.h.

◆ sys

template<typename Output , void(FEMContext::*)(unsigned int, const Point &, Output &, const Real) const point_output>
const System& libMesh::OldSolutionBase< Output, point_output >::sys
protectedinherited

Definition at line 361 of file generic_projector.h.


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