Line data Source code
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 : { 25 : /** 26 : * Takes an input functor that can be evaluated at faces, typically by linearly interpolating 27 : * between adjacent cell center values, and then creates an output functor whose cell-center 28 : * evaluations will correspond to weighted averages of the input functor's surrounding face 29 : * evaluations 30 : * @param output_functor the output functor 31 : * @param input_functor the input functor 32 : * @param num_int_recs the total number of interpolation and reconstruction operations to perform. 33 : * If this number is greater than 1, then this function will recurse 34 : * @param weight_with_sf when reconstructing the cell center value, decides whether the face values 35 : * (and maybe gradients) are weighted with the surface vector. If this is false, then the weights 36 : * are simply unity 37 : * @param faces the mesh faces we will be looping over for the interpolations and reconstructions 38 : */ 39 : template <typename T, typename Map> 40 : void 41 2040 : interpolateReconstruct(CellCenteredMapFunctor<T, Map> & output_functor, 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 2040 : if (!num_int_recs) 49 600 : return; 50 : 51 : std::unordered_map<dof_id_type, std::pair<T, Real>> elem_to_num_denom; 52 : 53 9098976 : for (const auto * const face : faces) 54 : { 55 : mooseAssert(face, "This must be non-null"); 56 9097536 : const Real weight = weight_with_sf ? face->faceArea() * face->faceCoord() : 1; 57 18195072 : const Moose::FaceArg face_arg{ 58 : face, 59 : Moose::FV::LimiterType::CentralDifference, 60 : true, 61 : false, 62 9097536 : input_functor.hasFaceSide(*face, true) 63 9097536 : ? (input_functor.hasFaceSide(*face, false) ? nullptr : face->elemPtr()) 64 : : face->neighborPtr(), 65 : nullptr}; 66 9097536 : auto face_value = input_functor(face_arg, time); 67 : std::pair<T, Real> * neighbor_pair = nullptr; 68 8880752 : if (face->neighborPtr() && face->neighborPtr() != libMesh::remote_elem) 69 : { 70 8880752 : neighbor_pair = &elem_to_num_denom[face->neighbor().id()]; 71 7511664 : neighbor_pair->first += face_value * weight; 72 8880752 : neighbor_pair->second += weight; 73 : } 74 9097536 : auto & elem_pair = elem_to_num_denom[face->elem().id()]; 75 7714112 : elem_pair.first += std::move(face_value) * weight; 76 9097536 : elem_pair.second += weight; 77 : } 78 : 79 4496672 : for (const auto & pair_ : elem_to_num_denom) 80 : { 81 : const auto & data_pair = pair_.second; 82 8302336 : output_functor[pair_.first] = data_pair.first / data_pair.second; 83 : } 84 : 85 1440 : interpolateReconstruct( 86 : output_functor, output_functor, num_int_recs - 1, weight_with_sf, faces, time); 87 : } 88 : } 89 : }