https://mooseframework.inl.gov
TestFaceCenteredMapFunctor.C
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 #include "gtest/gtest.h"
11 #include "MeshGeneratorMesh.h"
12 #include "GeneratedMeshGenerator.h"
13 #include "FaceCenteredMapFunctor.h"
14 #include "AppFactory.h"
15 #include "libmesh/quadrature_gauss.h"
16 #include "MooseMain.h"
17 
18 using namespace libMesh;
19 using namespace Moose;
20 using namespace FV;
21 
22 template <typename T, typename Map>
24 {
25 public:
27  const std::set<SubdomainID> & sub_ids,
28  const std::string & name)
29  : FaceCenteredMapFunctor<T, Map>(mesh, sub_ids, name)
30  {
31  }
32 
33  bool hasBlocks(SubdomainID) const override { return true; }
34 };
35 
36 TEST(FaceCenteredMapFunctorTest, testArgs)
37 {
38  const char * argv[2] = {"foo", "\0"};
39 
40  // First we create a simple mesh
41  auto app = Moose::createMooseApp("NavierStokesUnitApp", 1, (char **)argv);
42  auto * factory = &app->getFactory();
43  std::string mesh_type = "MeshGeneratorMesh";
44 
45  std::shared_ptr<MeshGeneratorMesh> mesh;
46  {
47  InputParameters params = factory->getValidParams(mesh_type);
48  mesh = factory->create<MeshGeneratorMesh>(mesh_type, "moose_mesh", params);
49  }
50 
51  app->actionWarehouse().mesh() = mesh;
52 
53  {
54  std::unique_ptr<MeshBase> lm_mesh;
55  InputParameters params = factory->getValidParams("GeneratedMeshGenerator");
56  params.set<unsigned int>("nx") = 2;
57  params.set<unsigned int>("ny") = 2;
58  params.set<MooseEnum>("dim") = "2";
59  auto mesh_gen =
60  factory->create<GeneratedMeshGenerator>("GeneratedMeshGenerator", "mesh_gen", params);
61  lm_mesh = mesh_gen->generate();
62  mesh->setMeshBase(std::move(lm_mesh));
63  }
64 
65  mesh->prepare(nullptr);
66  MultiMooseEnum coord_type_enum("XYZ RZ RSPHERICAL", "XYZ");
67  mesh->setCoordSystem({}, coord_type_enum);
68  mesh->buildFiniteVolumeInfo();
69  mesh->computeFiniteVolumeCoords();
70  const auto & all_fi = mesh->allFaceInfo();
71 
72  // We create a face-centered functor
74  "u");
75 
76  // We fill up the functor with known values
77  for (auto & fi : all_fi)
78  {
79  const auto & face_center = fi.faceCentroid();
80  u[fi.id()] = RealVectorValue(
81  -sin(face_center(0)) * cos(face_center(1)), cos(face_center(0)) * sin(face_center(1)), 0);
82  }
83 
84  // We check if the functor has the right face values
85  for (auto & fi : all_fi)
86  {
87  const auto & face_center = fi.faceCentroid();
88  const auto face_arg = Moose::FaceArg{
89  &fi, Moose::FV::LimiterType::CentralDifference, true, false, nullptr, nullptr};
90 
91  const auto result = u(face_arg, Moose::currentState());
92 
93  EXPECT_NEAR(result(0), -sin(face_center(0)) * cos(face_center(1)), 1e-14);
94  EXPECT_NEAR(result(1), cos(face_center(0)) * sin(face_center(1)), 1e-14);
95  EXPECT_EQ(result(2), 0);
96  }
97 
98  // Next, we test for the error messages of not-implemented methods
99  // This is a lambda for testing errors when requesting a gradient evaluation
100  auto test_gradient = [&u](const auto & arg)
101  {
102  try
103  {
104  u.gradient(arg, Moose::currentState());
105  EXPECT_TRUE(false);
106  }
107  catch (std::runtime_error & e)
108  {
109  EXPECT_TRUE(std::string(e.what()).find("not implemented") != std::string::npos);
110  }
111  };
112 
113  // This is a lambda for testing errors when requesting a regular evaluation with an unsupported
114  // argument
115  auto test_evaluate = [&u](const auto & arg)
116  {
117  try
118  {
119  u(arg, Moose::currentState());
120  EXPECT_TRUE(false);
121  }
122  catch (std::runtime_error & e)
123  {
124  EXPECT_TRUE(std::string(e.what()).find("not implemented") != std::string::npos);
125  }
126  };
127 
128  // Arguments for the simple error checks, we use the first face and the corresponding
129  // owner element
130  QGauss qrule(1, CONSTANT);
131  const auto face_arg = Moose::FaceArg{
132  &all_fi[0], Moose::FV::LimiterType::CentralDifference, true, false, nullptr, nullptr};
133  const auto elem_arg = ElemArg{all_fi[0].elemPtr(), false};
134  const auto elem_qp_arg = ElemQpArg({all_fi[0].elemPtr(), 0, &qrule, Point(0)});
135  const auto elem_side_qp_arg = ElemSideQpArg({all_fi[0].elemPtr(), 0, 0, &qrule, Point(0)});
136  const auto elem_point_arg = ElemPointArg({all_fi[0].elemPtr(), Point(0), false});
137 
138  test_gradient(elem_arg);
139  test_gradient(face_arg);
140 
141  test_evaluate(elem_qp_arg);
142  test_evaluate(elem_side_qp_arg);
143  test_evaluate(elem_point_arg);
144 
145  // Lastly, we check for errors when encountering faces with incorrect subdomains
147  unrestricted_error_test(*mesh, "not_restricted");
148  try
149  {
150  unrestricted_error_test(
151  FaceArg{&all_fi[2], LimiterType::CentralDifference, true, false, nullptr, nullptr},
153  EXPECT_TRUE(false);
154  }
155  catch (std::runtime_error & e)
156  {
157  EXPECT_TRUE(std::string(e.what()).find("not_restricted") != std::string::npos);
158  EXPECT_TRUE(std::string(e.what()).find("Make sure to fill") != std::string::npos);
159  }
160 
162  restricted_error_test(*mesh, {1}, "is_restricted");
163  try
164  {
165  restricted_error_test(
166  FaceArg{&all_fi[2], LimiterType::CentralDifference, true, false, nullptr, nullptr},
168  EXPECT_TRUE(false);
169  }
170  catch (std::runtime_error & e)
171  {
172  EXPECT_TRUE(std::string(e.what()).find("is_restricted") != std::string::npos);
173  EXPECT_TRUE(std::string(e.what()).find("0") != std::string::npos);
174  EXPECT_TRUE(
175  std::string(e.what()).find(
176  "that subdomain id is not one of the subdomain ids the functor is restricted to") !=
177  std::string::npos);
178  }
179 }
FictionalFaceCenteredMapFunctor(const MooseMesh &mesh, const std::set< SubdomainID > &sub_ids, const std::string &name)
std::shared_ptr< MooseApp > createMooseApp(const std::string &default_app_type, int argc, char *argv[])
CTSub CT_OPERATOR_BINARY CTMul CTCompareLess CTCompareGreater CTCompareEqual _arg template * sin(_arg) *_arg.template D< dtag >()) CT_SIMPLE_UNARY_FUNCTION(tan
A functor whose evaluation relies on querying a map where the keys are face info ids and the values c...
VectorValue< Real > RealVectorValue
T & set(const std::string &name, bool quiet_mode=false)
MeshBase & mesh
The following methods are specializations for using the Parallel::packed_range_* routines for a vecto...
CTSub CT_OPERATOR_BINARY CTMul CTCompareLess CTCompareGreater CTCompareEqual _arg template cos(_arg) *_arg.template D< dtag >()) CT_SIMPLE_UNARY_FUNCTION(cos
const std::string name
Definition: Setup.h:20
GradientType gradient(const ElemArg &elem, const StateArg &state) const
TEST(FaceCenteredMapFunctorTest, testArgs)
bool hasBlocks(SubdomainID) const override
StateArg currentState()