libMesh
mesh_function.C
Go to the documentation of this file.
1 #include <libmesh/equation_systems.h>
2 #include <libmesh/replicated_mesh.h>
3 #include <libmesh/mesh_generation.h>
4 #include <libmesh/dof_map.h>
5 #include <libmesh/system.h>
6 #include <libmesh/mesh_function.h>
7 #include <libmesh/numeric_vector.h>
8 #include <libmesh/elem.h>
9 
10 #include "test_comm.h"
11 #include "libmesh_cppunit.h"
12 
13 
14 using namespace libMesh;
15 
17  const Parameters &,
18  const std::string &,
19  const std::string &)
20 {
21  return
22  cos(.5*libMesh::pi*p(0)) *
23  sin(.5*libMesh::pi*p(1)) *
24  cos(.5*libMesh::pi*p(2));
25 }
26 
28  const Parameters &,
29  const std::string &,
30  const std::string &)
31 {
32  return 8*p(0) + 80*p(1) + 800*p(2);
33 }
34 
35 
36 class MeshFunctionTest : public CppUnit::TestCase
37 {
41 public:
42  LIBMESH_CPPUNIT_TEST_SUITE( MeshFunctionTest );
43 
44 #if LIBMESH_DIM > 1
45  CPPUNIT_TEST( test_subdomain_id_sets );
46 #endif
47 #if LIBMESH_DIM > 2
48 #ifdef LIBMESH_ENABLE_AMR
49  CPPUNIT_TEST( test_p_level );
50 #endif
51 #endif
52 
53  CPPUNIT_TEST_SUITE_END();
54 
55 protected:
56 
57 public:
58  void setUp() {}
59 
60  void tearDown() {}
61 
62  // test that mesh function works correctly with subdomain id sets.
64  {
65  LOG_UNIT_TEST;
66 
68 
70  4, 4,
71  0., 1.,
72  0., 1.,
73  QUAD4);
74 
75  // Set a subdomain id for all elements, based on location.
76  for (auto & elem : mesh.active_element_ptr_range())
77  {
78  Point c = elem->vertex_average();
79  elem->subdomain_id() =
80  subdomain_id_type(c(0)*4) + subdomain_id_type(c(1)*4)*10;
81  }
82 
83  // Add a discontinuous variable so we can easily see what side of
84  // an interface we're querying
86  System & sys = es.add_system<System> ("SimpleSystem");
87  unsigned int u_var = sys.add_variable("u", CONSTANT, MONOMIAL);
88 
89  es.init();
91 
92  MeshFunction mesh_function (sys.get_equation_systems(),
94  sys.get_dof_map(),
95  u_var);
96 
97  // Checkerboard pattern
98  const std::set<subdomain_id_type> sbdids1 {0,2,11,13,20,22,31,33};
99  const std::set<subdomain_id_type> sbdids2 {1,3,10,12,21,23,30,32};
100  mesh_function.init();
101  mesh_function.enable_out_of_mesh_mode(DenseVector<Number>());
102  mesh_function.set_subdomain_ids(&sbdids1);
103 
104  DenseVector<Number> vec_values;
105  const std::string dummy;
106 
107  // Make sure the MeshFunction's values interpolate the projected solution
108  // at the nodes
109  for (auto & elem : mesh.active_local_element_ptr_range())
110  {
111  const Point c = elem->vertex_average();
112  const Real expected_value =
113  libmesh_real(trilinear_function(c, es.parameters, dummy, dummy));
114  const std::vector<Point> offsets
115  {{0,-1/8.}, {1/8.,0}, {0,1/8.}, {-1/8.,0}};
116  for (Point offset : offsets)
117  {
118  const Point p = c + offset;
119  mesh_function(p, 0, vec_values, &sbdids1);
120  const Number retval1 = vec_values.empty() ? -12345 : vec_values(0);
121  mesh_function(p, 0, vec_values, &sbdids2);
122  const Number retval2 = vec_values.empty() ? -12345 : vec_values(0);
123  mesh_function(c, 0, vec_values, nullptr);
124  const Number retval3 = vec_values.empty() ? -12345 : vec_values(0);
125 
126  LIBMESH_ASSERT_NUMBERS_EQUAL
127  (retval3, expected_value, TOLERANCE * TOLERANCE);
128 
129  if (sbdids1.count(elem->subdomain_id()))
130  {
131  CPPUNIT_ASSERT(!sbdids2.count(elem->subdomain_id()));
132  LIBMESH_ASSERT_NUMBERS_EQUAL
133  (retval1, expected_value, TOLERANCE * TOLERANCE);
134 
135  mesh_function(c, 0, vec_values, &sbdids2);
136  CPPUNIT_ASSERT(vec_values.empty());
137  }
138  else
139  {
140  LIBMESH_ASSERT_NUMBERS_EQUAL
141  (retval2, expected_value, TOLERANCE * TOLERANCE);
142 
143  mesh_function(c, 0, vec_values, &sbdids1);
144  CPPUNIT_ASSERT(vec_values.empty());
145  }
146  }
147  }
148  }
149 
150  // test that mesh function works correctly with non-zero
151  // Elem::p_level() values.
152 #ifdef LIBMESH_ENABLE_AMR
154  {
155  LOG_UNIT_TEST;
156 
158 
160  5, 5, 5,
161  0., 1.,
162  0., 1.,
163  0., 1.,
164  HEX20);
165 
166  // Bump the p-level for all elements. We will create a system with
167  // a FIRST, LAGRANGE variable and then use the additional p-level
168  // to solve with quadratic elements. Note: set_p_level(1) actually
169  // _increases_ the existing variable order by 1, it does not set
170  // it to 1.
171  for (auto & elem : mesh.active_element_ptr_range())
172  elem->set_p_level(1);
173 
174  EquationSystems es(mesh);
175  System & sys = es.add_system<System> ("SimpleSystem");
176  unsigned int u_var = sys.add_variable("u", FIRST, LAGRANGE);
177 
178  es.init();
180 
181 
182  std::unique_ptr<NumericVector<Number>> mesh_function_vector
184  mesh_function_vector->init(sys.n_dofs(), sys.n_local_dofs(),
185  sys.get_dof_map().get_send_list(), false,
186  GHOSTED);
187 
188  sys.solution->localize(*mesh_function_vector,
189  sys.get_dof_map().get_send_list());
190 
191 
192  // So the MeshFunction knows which variables to compute values for.
193  // std::make_unique doesn't like if we try to use {u_var} in-place?
194  std::vector<unsigned int> variables {u_var};
195 
196  auto mesh_function =
197  std::make_unique<MeshFunction>(sys.get_equation_systems(),
198  *mesh_function_vector,
199  sys.get_dof_map(),
200  variables);
201 
202  mesh_function->init();
203  mesh_function->set_point_locator_tolerance(0.0001);
204  DenseVector<Number> vec_values;
205  std::string dummy;
206 
207  // Make sure the MeshFunction's values interpolate the projected solution
208  // at the nodes
209  for (const auto & node : mesh.local_node_ptr_range())
210  {
211  (*mesh_function)(*node, /*time=*/ 0., vec_values);
212  Number mesh_function_value =
213  projection_function(*node,
214  es.parameters,
215  dummy,
216  dummy);
217 
218  LIBMESH_ASSERT_NUMBERS_EQUAL
219  (vec_values(0), mesh_function_value, TOLERANCE*TOLERANCE);
220  }
221  }
222 #endif // LIBMESH_ENABLE_AMR
223 };
224 
225 
T libmesh_real(T a)
This is the EquationSystems class.
The ReplicatedMesh class is derived from the MeshBase class, and is used to store identical copies of...
This class provides the ability to map between arbitrary, user-defined strings and several data types...
Definition: parameters.h:67
Number trilinear_function(const Point &p, const Parameters &, const std::string &, const std::string &)
Definition: mesh_function.C:27
libMesh::Parallel::Communicator * TestCommWorld
Definition: driver.C:171
static constexpr Real TOLERANCE
TestClass subdomain_id_type
Based on the 4-byte comment warning above, this probably doesn&#39;t work with exodusII at all...
Definition: id_types.h:43
Number projection_function(const Point &p, const Parameters &, const std::string &, const std::string &)
Definition: mesh_function.C:16
const EquationSystems & get_equation_systems() const
Definition: system.h:721
MeshBase & mesh
const Parallel::Communicator & comm() const
void build_square(UnstructuredMesh &mesh, const unsigned int nx, const unsigned int ny, const Real xmin=0., const Real xmax=1., const Real ymin=0., const Real ymax=1., const ElemType type=INVALID_ELEM, const bool gauss_lobatto_grid=false)
A specialized build_cube() for 2D meshes.
The libMesh namespace provides an interface to certain functionality in the library.
dof_id_type n_local_dofs() const
Definition: system.C:158
dof_id_type n_dofs() const
Definition: system.C:121
void project_solution(FunctionBase< Number > *f, FunctionBase< Gradient > *g=nullptr) const
Projects arbitrary functions onto the current solution.
Manages consistently variables, degrees of freedom, and coefficient vectors.
Definition: system.h:96
virtual bool empty() const override final
Definition: dense_vector.h:109
std::unique_ptr< NumericVector< Number > > solution
Data structure to hold solution values.
Definition: system.h:1593
virtual void init() override
Override the FunctionBase::init() member function.
unsigned int add_variable(std::string_view var, const FEType &type, const std::set< subdomain_id_type > *const active_subdomains=nullptr)
Adds the variable var to the list of variables for this system.
Definition: system.C:1357
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void test_subdomain_id_sets()
Definition: mesh_function.C:63
static std::unique_ptr< NumericVector< T > > build(const Parallel::Communicator &comm, SolverPackage solver_package=libMesh::default_solver_package(), ParallelType parallel_type=AUTOMATIC)
Builds a NumericVector on the processors in communicator comm using the linear solver package specifi...
Parameters parameters
Data structure holding arbitrary parameters.
std::unique_ptr< NumericVector< Number > > current_local_solution
All the values I need to compute my contribution to the simulation at hand.
Definition: system.h:1605
virtual void init()
Initialize all the systems.
This class provides function-like objects for data distributed over a mesh.
Definition: mesh_function.h:54
virtual System & add_system(std::string_view system_type, std::string_view name)
Add the system of type system_type named name to the systems array.
CPPUNIT_TEST_SUITE_REGISTRATION(MeshFunctionTest)
const DofMap & get_dof_map() const
Definition: system.h:2374
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:39
const std::vector< dof_id_type > & get_send_list() const
Definition: dof_map.h:526
void build_cube(UnstructuredMesh &mesh, const unsigned int nx=0, const unsigned int ny=0, const unsigned int nz=0, const Real xmin=0., const Real xmax=1., const Real ymin=0., const Real ymax=1., const Real zmin=0., const Real zmax=1., const ElemType type=INVALID_ELEM, const bool gauss_lobatto_grid=false)
Builds a (elements) cube.
const Real pi
.
Definition: libmesh.h:299