https://mooseframework.inl.gov
Reconstructions.h
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
10 #pragma once
11 
12 #include "MooseMesh.h"
13 #include "FaceInfo.h"
14 #include "CellCenteredMapFunctor.h"
15 #include "MathFVUtils.h"
16 #include "libmesh/elem.h"
17 
18 #include <unordered_map>
19 #include <utility>
20 
21 namespace Moose
22 {
23 namespace FV
24 {
39 template <typename T, typename Map>
40 void
42  const Moose::FunctorBase<T> & input_functor,
43  const unsigned int num_int_recs,
44  const bool weight_with_sf,
45  const std::vector<const FaceInfo *> & faces,
46  const Moose::StateArg & time)
47 {
48  if (!num_int_recs)
49  return;
50 
51  std::unordered_map<dof_id_type, std::pair<T, Real>> elem_to_num_denom;
52 
53  for (const auto * const face : faces)
54  {
55  mooseAssert(face, "This must be non-null");
56  const Real weight = weight_with_sf ? face->faceArea() * face->faceCoord() : 1;
57  const Moose::FaceArg face_arg{
58  face,
60  true,
61  false,
62  input_functor.hasFaceSide(*face, true)
63  ? (input_functor.hasFaceSide(*face, false) ? nullptr : face->elemPtr())
64  : face->neighborPtr(),
65  nullptr};
66  auto face_value = input_functor(face_arg, time);
67  std::pair<T, Real> * neighbor_pair = nullptr;
68  if (face->neighborPtr() && face->neighborPtr() != libMesh::remote_elem)
69  {
70  neighbor_pair = &elem_to_num_denom[face->neighbor().id()];
71  neighbor_pair->first += face_value * weight;
72  neighbor_pair->second += weight;
73  }
74  auto & elem_pair = elem_to_num_denom[face->elem().id()];
75  elem_pair.first += std::move(face_value) * weight;
76  elem_pair.second += weight;
77  }
78 
79  for (const auto & pair_ : elem_to_num_denom)
80  {
81  const auto & data_pair = pair_.second;
82  output_functor[pair_.first] = data_pair.first / data_pair.second;
83  }
84 
86  output_functor, output_functor, num_int_recs - 1, weight_with_sf, faces, time);
87 }
88 }
89 }
void interpolateReconstruct(CellCenteredMapFunctor< T, Map > &output_functor, const Moose::FunctorBase< T > &input_functor, const unsigned int num_int_recs, const bool weight_with_sf, const std::vector< const FaceInfo *> &faces, const Moose::StateArg &time)
Takes an input functor that can be evaluated at faces, typically by linearly interpolating between ad...
dof_id_type weight(const MeshBase &mesh, const processor_id_type pid)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
A functor whose evaluation relies on querying a map where the keys are element ids and the values cor...
virtual bool hasFaceSide(const FaceInfo &fi, const bool fi_elem_side) const override
const RemoteElem * remote_elem