149 libmesh_assert_equal_to (system_name,
"EllipticDG");
160 std::string refinement_type = es.
parameters.
get<std::string> (
"refinement");
175 std::unique_ptr<FEBase> fe (FEBase::build(
dim, fe_type));
176 std::unique_ptr<FEBase> fe_elem_face(FEBase::build(
dim, fe_type));
177 std::unique_ptr<FEBase> fe_neighbor_face(FEBase::build(
dim, fe_type));
185 fe->attach_quadrature_rule (&qrule);
194 fe_elem_face->attach_quadrature_rule(&qface);
195 fe_neighbor_face->attach_quadrature_rule(&qface);
200 const std::vector<Real> & JxW = fe->get_JxW();
201 const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
204 const std::vector<std::vector<Real>> & phi_face = fe_elem_face->get_phi();
205 const std::vector<std::vector<RealGradient>> & dphi_face = fe_elem_face->get_dphi();
206 const std::vector<Real> & JxW_face = fe_elem_face->get_JxW();
207 const std::vector<Point> & qface_normals = fe_elem_face->get_normals();
208 const std::vector<Point> & qface_points = fe_elem_face->get_xyz();
211 const std::vector<std::vector<Real>> & phi_neighbor_face = fe_neighbor_face->get_phi();
212 const std::vector<std::vector<RealGradient>> & dphi_neighbor_face = fe_neighbor_face->get_dphi();
233 std::vector<dof_id_type> dof_indices;
247 for (
const auto & elem :
mesh.active_local_element_ptr_range())
253 dof_map.dof_indices (elem, dof_indices);
254 const unsigned int n_dofs =
255 cast_int<unsigned int>(dof_indices.size());
267 Ke.
resize (n_dofs, n_dofs);
273 for (
unsigned int qp=0; qp<qrule.n_points(); qp++)
274 for (
unsigned int i=0; i<n_dofs; i++)
275 for (
unsigned int j=0; j<n_dofs; j++)
276 Ke(i,j) += JxW[qp]*(dphi[i][qp]*dphi[j][qp]);
283 for (
auto side : elem->side_index_range())
285 if (elem->neighbor_ptr(side) ==
nullptr)
288 fe_elem_face->reinit(elem, side);
291 const auto side_volume = side_builder(*elem, side).volume();
292 const unsigned int elem_b_order =
static_cast<unsigned int> (fe_elem_face->get_order());
293 const Real h_elem = elem->volume()/side_volume * 1./
pow(elem_b_order, 2.);
295 for (
unsigned int qp=0; qp<qface.n_points(); qp++)
298 for (
unsigned int i=0; i<n_dofs; i++)
301 for (
unsigned int j=0; j<n_dofs; j++)
304 Ke(i,j) += JxW_face[qp] * penalty/h_elem * phi_face[i][qp] * phi_face[j][qp];
309 (phi_face[i][qp] * (dphi_face[j][qp]*qface_normals[qp]) +
310 phi_face[j][qp] * (dphi_face[i][qp]*qface_normals[qp]));
316 Fe(i) += JxW_face[qp] * bc_value * penalty/h_elem * phi_face[i][qp];
319 Fe(i) -= JxW_face[qp] * dphi_face[i][qp] * (bc_value*qface_normals[qp]);
334 const unsigned int elem_id = elem->
id();
335 const unsigned int neighbor_id = neighbor->
id();
342 if ((neighbor->
active() &&
343 (neighbor->
level() == elem->level()) &&
344 (elem_id < neighbor_id)) ||
345 (neighbor->
level() < elem->level()))
347 const Elem & elem_side = side_builder(*elem, side);
350 const unsigned int elem_b_order =
static_cast<unsigned int>(fe_elem_face->get_order());
351 const unsigned int neighbor_b_order =
static_cast<unsigned int>(fe_neighbor_face->get_order());
352 const double side_order = (elem_b_order + neighbor_b_order)/2.;
353 const Real h_elem = (elem->volume()/elem_side.
volume()) * 1./
pow(side_order,2.);
356 std::vector<Point> qface_neighbor_point;
359 std::vector<Point > qface_point;
362 fe_elem_face->reinit(elem, side);
365 qface_point = fe_elem_face->get_xyz();
369 if (refinement_type ==
"p")
370 fe_neighbor_face->side_map (neighbor,
374 qface_neighbor_point);
376 FEMap::inverse_map (elem->dim(), neighbor,
378 qface_neighbor_point);
381 fe_neighbor_face->reinit(neighbor, &qface_neighbor_point);
386 std::vector<dof_id_type> neighbor_dof_indices;
387 dof_map.dof_indices (neighbor, neighbor_dof_indices);
388 const unsigned int n_neighbor_dofs =
389 cast_int<unsigned int>(neighbor_dof_indices.size());
397 Kne.
resize (n_neighbor_dofs, n_dofs);
398 Ken.
resize (n_dofs, n_neighbor_dofs);
399 Kee.
resize (n_dofs, n_dofs);
400 Knn.
resize (n_neighbor_dofs, n_neighbor_dofs);
406 for (
unsigned int qp=0; qp<qface.n_points(); qp++)
410 for (
unsigned int i=0; i<n_dofs; i++)
412 for (
unsigned int j=0; j<n_dofs; j++)
417 (phi_face[j][qp]*(qface_normals[qp]*dphi_face[i][qp]) +
418 phi_face[i][qp]*(qface_normals[qp]*dphi_face[j][qp]));
421 Kee(i,j) += JxW_face[qp] * penalty/h_elem * phi_face[j][qp]*phi_face[i][qp];
427 for (
unsigned int i=0; i<n_neighbor_dofs; i++)
429 for (
unsigned int j=0; j<n_neighbor_dofs; j++)
434 (phi_neighbor_face[j][qp]*(qface_normals[qp]*dphi_neighbor_face[i][qp]) +
435 phi_neighbor_face[i][qp]*(qface_normals[qp]*dphi_neighbor_face[j][qp]));
439 JxW_face[qp] * penalty/h_elem * phi_neighbor_face[j][qp]*phi_neighbor_face[i][qp];
445 for (
unsigned int i=0; i<n_neighbor_dofs; i++)
447 for (
unsigned int j=0; j<n_dofs; j++)
452 (phi_neighbor_face[i][qp]*(qface_normals[qp]*dphi_face[j][qp]) -
453 phi_face[j][qp]*(qface_normals[qp]*dphi_neighbor_face[i][qp]));
456 Kne(i,j) -= JxW_face[qp] * penalty/h_elem * phi_face[j][qp]*phi_neighbor_face[i][qp];
462 for (
unsigned int i=0; i<n_dofs; i++)
464 for (
unsigned int j=0; j<n_neighbor_dofs; j++)
469 (phi_neighbor_face[j][qp]*(qface_normals[qp]*dphi_face[i][qp]) -
470 phi_face[i][qp]*(qface_normals[qp]*dphi_neighbor_face[j][qp]));
473 Ken(i,j) -= JxW_face[qp] * penalty/h_elem * phi_face[i][qp]*phi_neighbor_face[j][qp];
481 matrix.
add_matrix(Kne, neighbor_dof_indices, dof_indices);
482 matrix.
add_matrix(Ken, dof_indices, neighbor_dof_indices);
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Manages consistently variables, degrees of freedom, coefficient vectors, matrices and linear solvers ...
void resize(const unsigned int n)
Resize the vector.
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.
NumericVector< Number > * rhs
The system matrix.
Order default_quadrature_order() const
const T_sys & get_system(std::string_view name) const
This is the MeshBase class.
This class handles the numbering of degrees of freedom on a mesh.
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.
const T & get(std::string_view) const
unsigned int which_neighbor_am_i(const Elem *e) const
This function tells you which neighbor e is.
Number exact_solution(const Point &p, const Parameters ¶meters, const std::string &, const std::string &)
Helper for building element sides that minimizes the construction of new elements.
const Elem * neighbor_ptr(unsigned int i) const
BasicOStreamProxy & flush()
Flush the associated stream buffer.
unsigned int level() const
const FEType & variable_type(const unsigned int i) const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
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().
virtual Real volume() const
This class implements specific orders of Gauss quadrature.
unsigned int mesh_dimension() const
Parameters parameters
Data structure holding arbitrary parameters.
const DofMap & get_dof_map() const
const SparseMatrix< Number > & get_system_matrix() const