Definition at line 42 of file mesh_function_dfem.C.
◆ build_mesh()
void MeshfunctionDFEM::build_mesh |
( |
| ) |
|
|
inlineprotected |
◆ CPPUNIT_TEST() [1/3]
◆ CPPUNIT_TEST() [2/3]
◆ CPPUNIT_TEST() [3/3]
◆ CPPUNIT_TEST_SUITE()
The goal of this test is to verify proper operation of mesh function for discontinuous shape functions and points close to the boundary.
◆ CPPUNIT_TEST_SUITE_END()
MeshfunctionDFEM::CPPUNIT_TEST_SUITE_END |
( |
| ) |
|
◆ setUp()
void MeshfunctionDFEM::setUp |
( |
| ) |
|
|
inline |
◆ tearDown()
void MeshfunctionDFEM::tearDown |
( |
| ) |
|
|
inline |
◆ test_mesh_function_dfem()
void MeshfunctionDFEM::test_mesh_function_dfem |
( |
| ) |
|
|
inline |
Definition at line 164 of file mesh_function_dfem.C.
167 Point top(0.5, 0.5, 0.0);
170 Point bottom(0.5, -0.5, 0.0);
173 Point face(0.5, 0.0, 0.0);
183 std::vector<unsigned int> variables;
185 std::unique_ptr<NumericVector<Number>> mesh_function_vector =
188 sys.
solution->localize(*mesh_function_vector);
192 mesh_function.
init();
195 std::map<const Elem *, Number> top_val = mesh_function.discontinuous_value(top);
198 CPPUNIT_ASSERT (top_val.size() == 1);
201 for (std::map<const Elem *, Number>::const_iterator it = top_val.begin(); it != top_val.end(); ++it)
202 CPPUNIT_ASSERT (it->first->id() == 0 &&
std::abs(it->second) < 1.0e-10);
205 std::map<const Elem *, Number> bottom_val = mesh_function.discontinuous_value(bottom);
208 CPPUNIT_ASSERT (bottom_val.size() == 1);
211 for (std::map<const Elem *, Number>::const_iterator it = bottom_val.begin(); it != bottom_val.end(); ++it)
212 CPPUNIT_ASSERT (it->first->id() == 1 &&
std::abs(it->second - 1.0) < 1.0e-10);
215 std::map<const Elem *, Number> face_val = mesh_function.discontinuous_value(face);
218 CPPUNIT_ASSERT (face_val.size() == 2);
221 for (std::map<const Elem *, Number>::const_iterator it = face_val.begin(); it != face_val.end(); ++it)
222 if (it->first->id() == 0)
223 CPPUNIT_ASSERT (
std::abs(it->second) < 1.0e-10);
224 else if (it->first->id() == 1)
225 CPPUNIT_ASSERT (
std::abs(it->second - 1.0) < 1.0e-10);
227 CPPUNIT_ASSERT (
false);
References std::abs(), libMesh::EquationSystems::add_system(), libMesh::System::add_variable(), libMesh::NumericVector< T >::build(), libMesh::ParallelObject::comm(), libMesh::CONSTANT, libMesh::System::get_all_variable_numbers(), libMesh::System::get_dof_map(), libMesh::MeshFunction::init(), libMesh::EquationSystems::init(), libMesh::MONOMIAL, libMesh::System::n_dofs(), libMesh::EquationSystems::parameters, position_function(), libMesh::System::project_solution(), libMesh::SERIAL, and libMesh::System::solution.
◆ test_mesh_function_dfem_grad()
void MeshfunctionDFEM::test_mesh_function_dfem_grad |
( |
| ) |
|
|
inline |
Definition at line 232 of file mesh_function_dfem.C.
235 Point top(0.5, 0.5, 0.0);
238 Point bottom(0.5, -0.5, 0.0);
241 Point face(0.5, 0.0, 0.0);
251 std::vector<unsigned int> variables;
253 std::unique_ptr<NumericVector<Number>> mesh_function_vector =
256 sys.
solution->localize(*mesh_function_vector);
260 mesh_function.
init();
263 std::map<const Elem *, Gradient> top_val = mesh_function.discontinuous_gradient(top);
266 CPPUNIT_ASSERT (top_val.size() == 1);
269 for (std::map<const Elem *, Gradient>::const_iterator it = top_val.begin(); it != top_val.end(); ++it)
270 CPPUNIT_ASSERT (it->first->id() == 0 &&
std::abs(it->second(1) - 2.0) < 1.0e-10);
273 std::map<const Elem *, Gradient> bottom_val = mesh_function.discontinuous_gradient(bottom);
276 CPPUNIT_ASSERT (bottom_val.size() == 1);
279 for (std::map<const Elem *, Gradient>::const_iterator it = bottom_val.begin(); it != bottom_val.end(); ++it)
280 CPPUNIT_ASSERT (it->first->id() == 1 &&
std::abs(it->second(1) - 1.0) < 1.0e-10);
283 std::map<const Elem *, Gradient> face_val = mesh_function.discontinuous_gradient(face);
286 CPPUNIT_ASSERT (face_val.size() == 2);
289 for (std::map<const Elem *, Gradient>::const_iterator it = face_val.begin(); it != face_val.end(); ++it)
290 if (it->first->id() == 0)
291 CPPUNIT_ASSERT (
std::abs(it->second(1) - 2.0) < 1.0e-10);
292 else if (it->first->id() == 1)
293 CPPUNIT_ASSERT (
std::abs(it->second(1) - 1.0) < 1.0e-10);
295 CPPUNIT_ASSERT (
false);
References std::abs(), libMesh::EquationSystems::add_system(), libMesh::System::add_variable(), libMesh::NumericVector< T >::build(), libMesh::ParallelObject::comm(), libMesh::FIRST, libMesh::System::get_all_variable_numbers(), libMesh::System::get_dof_map(), libMesh::MeshFunction::init(), libMesh::EquationSystems::init(), libMesh::LAGRANGE, libMesh::System::n_dofs(), libMesh::EquationSystems::parameters, position_function2(), libMesh::System::project_solution(), libMesh::SERIAL, and libMesh::System::solution.
◆ test_point_locator_dfem()
void MeshfunctionDFEM::test_point_locator_dfem |
( |
| ) |
|
|
inline |
Definition at line 133 of file mesh_function_dfem.C.
136 Point interior(0.5, -0.5, 0.0);
139 Point face(0.5, 0.0, 0.0);
142 std::set<const Elem *> int_cand;
143 (*_point_locator)(interior, int_cand,
nullptr);
144 CPPUNIT_ASSERT (int_cand.size() == 1);
145 for (std::set<const Elem *>::iterator it = int_cand.begin(); it != int_cand.end(); ++it)
146 CPPUNIT_ASSERT ((*it)->id() == 1);
149 std::set<const Elem *> face_cand;
150 (*_point_locator)(face, face_cand,
nullptr);
151 CPPUNIT_ASSERT (face_cand.size() == 2);
152 int array[2] = {0, 0};
153 for (std::set<const Elem *>::iterator it = face_cand.begin(); it != face_cand.end(); ++it)
156 CPPUNIT_ASSERT (array[(*it)->id()] == 0);
157 array[(*it)->id()] = 1;
159 CPPUNIT_ASSERT (array[0] == 1);
160 CPPUNIT_ASSERT (array[1] == 1);
◆ _mesh
◆ _point_locator
The documentation for this class was generated from the following file:
Manages consistently variables, degrees of freedom, and coefficient vectors.
virtual void init() override
Override the FunctionBase::init() member function by calling our own and specifying the Trees::NODES ...
virtual const Node * node_ptr(const dof_id_type i) const override
This class provides function-like objects for data distributed over a mesh.
virtual void init(const numeric_index_type n, const numeric_index_type n_local, const bool fast=false, const ParallelType ptype=AUTOMATIC)=0
Change the dimension of the vector to n.
virtual Elem * add_elem(Elem *e) override
Add elem e to the end of the element array.
Number position_function(const Point &p, const Parameters &, const std::string &, const std::string &)
The ReplicatedMesh class is derived from the MeshBase class, and is used to store identical copies of...
Provides a uniform interface to vector storage schemes for different linear algebra libraries.
MetaPhysicL::DualNumber< T, D > abs(const MetaPhysicL::DualNumber< T, D > &in)
unsigned int add_variable(const std::string &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.
libMesh::Parallel::Communicator * TestCommWorld
A Point defines a location in LIBMESH_DIM dimensional Real space.
Number position_function2(const Point &p, const Parameters &, const std::string &, const std::string &)
The QUAD4 is an element in 2D composed of 4 nodes.
This is the EquationSystems class.
virtual Node *& set_node(const unsigned int i)
std::unique_ptr< PointLocatorBase > sub_point_locator() const
std::unique_ptr< NumericVector< Number > > solution
Data structure to hold solution values.
This is the base class from which all geometric element types are derived.
const DofMap & get_dof_map() const
dof_id_type n_dofs() const
void prepare_for_use(const bool skip_renumber_nodes_and_elements=false, const bool skip_find_neighbors=false)
Prepare a newly ecreated (or read) mesh for use.
std::unique_ptr< PointLocatorBase > _point_locator
void get_all_variable_numbers(std::vector< unsigned int > &all_variable_numbers) const
Fills all_variable_numbers with all the variable numbers for the variables that have been added to th...
void project_solution(FunctionBase< Number > *f, FunctionBase< Gradient > *g=nullptr) const
Projects arbitrary functions onto the current solution.
void set_mesh_dimension(unsigned char d)
Resets the logical dimension of the mesh.
virtual Node * add_point(const Point &p, const dof_id_type id=DofObject::invalid_id, const processor_id_type proc_id=DofObject::invalid_processor_id) override
functions for adding /deleting nodes elements.