libMesh
Classes | Functions | Variables
periodic_bc_test.C File Reference

Go to the source code of this file.

Classes

struct  PeriodicQuadFunction
 
class  PeriodicBCTest
 

Functions

Number quadratic_solution (const Point &p, const Parameters &, const std::string &, const std::string &)
 
void periodic_bc_test_poisson (EquationSystems &es, const std::string &)
 
 CPPUNIT_TEST_SUITE_REGISTRATION (PeriodicBCTest)
 

Variables

const boundary_id_type top_id = 50
 
const boundary_id_type bottom_id = 60
 
const boundary_id_type side_id = 70
 

Function Documentation

◆ CPPUNIT_TEST_SUITE_REGISTRATION()

CPPUNIT_TEST_SUITE_REGISTRATION ( PeriodicBCTest  )

◆ periodic_bc_test_poisson()

void periodic_bc_test_poisson ( EquationSystems es,
const std::string &   
)

Definition at line 74 of file periodic_bc_test.C.

References libMesh::SparseMatrix< T >::add_matrix(), libMesh::NumericVector< T >::add_vector(), libMesh::FEGenericBase< OutputType >::build(), libMesh::NumericVector< T >::close(), libMesh::SparseMatrix< T >::close(), libMesh::FEType::default_quadrature_order(), libMesh::System::get_dof_map(), libMesh::EquationSystems::get_mesh(), libMesh::EquationSystems::get_system(), libMesh::ImplicitSystem::matrix, mesh, libMesh::DenseVector< T >::resize(), libMesh::DenseMatrix< T >::resize(), libMesh::ExplicitSystem::rhs, libMesh::DenseVector< T >::size(), and libMesh::TOLERANCE.

Referenced by PeriodicBCTest::testPeriodicBC().

76 {
77  const MeshBase& mesh = es.get_mesh();
79  const DofMap& dof_map = system.get_dof_map();
80 
81  // Interior FE to integrate smooth part of poisson problem
82  FEType fe_type = dof_map.variable_type(0);
83  std::unique_ptr<FEBase> fe (FEBase::build(2, fe_type));
84  QGauss qrule (2, fe_type.default_quadrature_order());
85  fe->attach_quadrature_rule(&qrule);
86 
87  const std::vector<Real> & JxW = fe->get_JxW();
88  const std::vector<std::vector<Real>> & phi = fe->get_phi();
89  const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
90 
91  // Face FE to integrate delta function part of forcing
92  std::unique_ptr<FEBase> fe_face (FEBase::build(2, fe_type));
93  QGauss qface (1, fe_type.default_quadrature_order());
94  fe_face->attach_quadrature_rule(&qface);
95 
96  const std::vector<Real> & JxW_face = fe_face->get_JxW();
97  const std::vector<std::vector<Real>> & phi_face = fe_face->get_phi();
98  std::unique_ptr<const Elem> elem_side;
99 
102 
103  std::vector<dof_id_type> dof_indices;
104 
105  for (const Elem * elem : mesh.active_local_element_ptr_range())
106  {
107  dof_map.dof_indices (elem, dof_indices);
108  const unsigned int n_dofs = dof_indices.size();
109 
110  Ke.resize (n_dofs, n_dofs);
111  Fe.resize (n_dofs);
112 
113  fe->reinit (elem);
114 
115  // Integrate the poisson problem over the element interior, where we have
116  // a nice f=2 forcing function aside from the discontinuity
117  for (unsigned int qp=0; qp<qrule.n_points(); qp++)
118  {
119  for (unsigned int i=0; i != n_dofs; i++)
120  {
121  for (unsigned int j=0; j != n_dofs; j++)
122  Ke(i,j) += -JxW[qp]*(dphi[i][qp]*dphi[j][qp]);
123 
124  Fe(i) += JxW[qp]*phi[i][qp]*2;
125  }
126  }
127 
128  for(unsigned int s=0; s != elem->n_sides(); ++s)
129  {
130  elem->build_side_ptr(elem_side, s);
131  Point centroid = elem_side->vertex_average();
132  if (std::abs(centroid(1) - 1) < TOLERANCE)
133  {
134  fe_face->reinit(elem, s);
135 
136  for (unsigned int qp=0; qp<qface.n_points(); qp++)
137  {
138  // delta function divided by 2 since we hit it from
139  // both sides
140  for (unsigned int i=0; i != n_dofs; i++)
141  Fe(i) += JxW_face[qp]*phi_face[i][qp]*(-4);
142  }
143  }
144  }
145 
146  dof_map.heterogenously_constrain_element_matrix_and_vector(Ke, Fe, dof_indices);
147  system.matrix->add_matrix (Ke, dof_indices);
148  system.rhs->add_vector (Fe, dof_indices);
149  }
150 
151  system.rhs->close();
152  system.matrix->close();
153 }
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:196
static constexpr Real TOLERANCE
Manages consistently variables, degrees of freedom, coefficient vectors, matrices and linear solvers ...
void resize(const unsigned int n)
Resize the vector.
Definition: dense_vector.h:396
virtual void add_vector(const T *v, const std::vector< numeric_index_type > &dof_indices)
Computes , where v is a pointer and each dof_indices[i] specifies where to add value v[i]...
This is the base class from which all geometric element types are derived.
Definition: elem.h:94
MeshBase & mesh
NumericVector< Number > * rhs
The system matrix.
Order default_quadrature_order() const
Definition: fe_type.h:371
const T_sys & get_system(std::string_view name) const
This is the MeshBase class.
Definition: mesh_base.h:75
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:179
virtual void add_matrix(const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &rows, const std::vector< numeric_index_type > &cols)=0
Add the full matrix dm to the SparseMatrix.
virtual void close()=0
Calls the NumericVector&#39;s internal assembly routines, ensuring that the values are consistent across ...
virtual void close()=0
Calls the SparseMatrix&#39;s internal assembly routines, ensuring that the values are consistent across p...
SparseMatrix< Number > * matrix
The system matrix.
const MeshBase & get_mesh() const
void resize(const unsigned int new_m, const unsigned int new_n)
Resizes the matrix to the specified size and calls zero().
Definition: dense_matrix.h:895
This class implements specific orders of Gauss quadrature.
virtual unsigned int size() const override final
Definition: dense_vector.h:104
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

◆ quadratic_solution()

Number quadratic_solution ( const Point p,
const Parameters ,
const std::string &  ,
const std::string &   
)

Definition at line 23 of file periodic_bc_test.C.

References libMesh::Real.

Referenced by PeriodicQuadFunction::component(), PeriodicQuadFunction::operator()(), and PeriodicBCTest::testPeriodicBC().

27 {
28  const Real & y = LIBMESH_DIM > 1 ? p(1) : 0;
29 
30  // Discontinuity at y=1 from the forcing function there.
31  // Periodic on -3 < y < 2
32  return (y > 1) ? (y-3)*(y-3) : (y+1)*(y+1);
33 }
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real

Variable Documentation

◆ bottom_id

const boundary_id_type bottom_id = 60

Definition at line 157 of file periodic_bc_test.C.

Referenced by PeriodicBCTest::testPeriodicBC().

◆ side_id

const boundary_id_type side_id = 70

◆ top_id

const boundary_id_type top_id = 50