Go to the documentation of this file.
   42 #include "libmesh/exodusII_io.h" 
   43 #include "libmesh/libmesh.h" 
   44 #include "libmesh/mesh.h" 
   45 #include "libmesh/mesh_generation.h" 
   46 #include "libmesh/linear_implicit_system.h" 
   47 #include "libmesh/equation_systems.h" 
   48 #include "libmesh/enum_xdr_mode.h" 
   51 #include "libmesh/fe.h" 
   52 #include "libmesh/inf_fe.h" 
   53 #include "libmesh/inf_elem_builder.h" 
   56 #include "libmesh/quadrature_gauss.h" 
   60 #include "libmesh/sparse_matrix.h" 
   61 #include "libmesh/numeric_vector.h" 
   62 #include "libmesh/dense_matrix.h" 
   63 #include "libmesh/dense_vector.h" 
   67 #include "libmesh/dof_map.h" 
   70 #include "libmesh/node.h" 
   73 #include "libmesh/elem.h" 
   81                     const std::string & system_name);
 
   84 int main (
int argc, 
char ** argv)
 
   90 #ifndef LIBMESH_ENABLE_INFINITE_ELEMENTS 
   91   libmesh_example_requires(
false, 
"--enable-ifem");
 
   95   libmesh_example_requires(3 <= LIBMESH_DIM, 
"3D support");
 
   98   libMesh::out << 
"Running ex6 with dim = 3" << std::endl << std::endl;
 
  117 #ifdef LIBMESH_HAVE_EXODUS_API 
  142     if (elem->infinite())
 
  143       elem->subdomain_id() = 1;
 
  150 #ifdef LIBMESH_HAVE_EXODUS_API 
  181   equation_systems.
get_system(
"Wave").add_variable(
"p", fe_type);
 
  194   equation_systems.
init();
 
  212   equation_systems.
write (
"eqn_sys.dat", 
WRITE);
 
  219 #endif // else part of ifndef LIBMESH_ENABLE_INFINITE_ELEMENTS 
  225                    const std::string & libmesh_dbg_var(system_name))
 
  229   libmesh_assert_equal_to (system_name, 
"Wave");
 
  235 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 
  256   const FEType & fe_type = dof_map.variable_type(0);
 
  270   fe->attach_quadrature_rule (&qrule);
 
  277   inf_fe->attach_quadrature_rule (&qrule);
 
  296   std::vector<dof_id_type> dof_indices;
 
  307       dof_map.dof_indices (elem, dof_indices);
 
  309       const unsigned int n_dofs =
 
  310         cast_int<unsigned int>(dof_indices.size());
 
  328       if (elem->infinite())
 
  359       const std::vector<Real> & JxW = cfe->get_JxW();
 
  362       const std::vector<std::vector<Real>> & phi = cfe->get_phi();
 
  366       const std::vector<std::vector<RealGradient>> & dphi = cfe->get_dphi();
 
  376       const std::vector<RealGradient> & dphase  = cfe->get_dphase();
 
  377       const std::vector<Real> &         
weight  = cfe->get_Sobolev_weight();
 
  378       const std::vector<RealGradient> & dweight = cfe->get_Sobolev_dweight();
 
  389       Ke.
resize (n_dofs, n_dofs);
 
  390       Ce.
resize (n_dofs, n_dofs);
 
  391       Me.
resize (n_dofs, n_dofs);
 
  399       unsigned int max_qp = cfe->n_quadrature_points();
 
  402       for (
unsigned int qp=0; qp<max_qp; qp++)
 
  408           const unsigned int n_sf = cfe->n_shape_functions();
 
  425           for (
unsigned int i=0; i<n_sf; i++)
 
  426             for (
unsigned int j=0; j<n_sf; j++)
 
  431                    (dweight[qp] * phi[i][qp] 
 
  440                    (dphase[qp] * dphi[j][qp]) 
 
  443                    (dweight[qp] * dphase[qp]) 
 
  444                    * phi[i][qp] * phi[j][qp]  
 
  446                    (dphi[i][qp] * dphase[qp]) 
 
  453                    (1. - (dphase[qp] * dphase[qp]))       
 
  454                    * phi[i][qp] * phi[j][qp] * 
weight[qp] 
 
  463       Ke.
add(1./speed        , Ce);
 
  464       Ke.
add(1./(speed*speed), Me);
 
  468       dof_map.constrain_element_matrix(Ke, dof_indices);
 
  481         unsigned int dn = node->dof_number(0,0,0);
 
  491 #endif //ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 
  
const Point build_inf_elem(const bool be_verbose=false)
Build infinite elements atop a volume-based mesh, determine origin automatically.
 
virtual void add(const numeric_index_type i, const T value)=0
Adds value to each entry of the vector.
 
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
 
virtual System & add_system(const std::string &system_type, const std::string &name)
Add the system of type system_type named name to the systems array.
 
const MeshBase & get_mesh() const
 
NumericVector< Number > * rhs
The system matrix.
 
This class implements specific orders of Gauss quadrature.
 
virtual SimpleRange< element_iterator > active_local_element_ptr_range()=0
 
The libMesh namespace provides an interface to certain functionality in the library.
 
boostcopy::enable_if_c< ScalarTraits< T2 >::value, void >::type add(const T2 factor, const DenseMatrix< T3 > &mat)
Adds factor times mat to this matrix.
 
This class forms the foundation from which generic finite elements may be derived.
 
const T_sys & get_system(const std::string &name) const
 
static const Real TOLERANCE
 
void assemble_wave(EquationSystems &es, const std::string &system_name)
 
unsigned int mesh_dimension() const
 
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
 
static std::unique_ptr< FEGenericBase > build_InfFE(const unsigned int dim, const FEType &type)
Builds a specific infinite element type.
 
void resize(const unsigned int new_m, const unsigned int new_n)
Resize the matrix.
 
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
 
virtual SimpleRange< element_iterator > element_ptr_range()=0
 
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 MeshBase class.
 
MetaPhysicL::DualNumber< T, D > abs(const MetaPhysicL::DualNumber< T, D > &in)
 
void libmesh_ignore(const Args &...)
 
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
 
virtual SimpleRange< node_iterator > local_node_ptr_range()=0
 
int main(int argc, char **argv)
 
virtual void init()
Initialize all the systems.
 
This class is used to build infinite elements on top of an existing mesh.
 
SparseMatrix< Number > * matrix
The system matrix.
 
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
 
void print_info(std::ostream &os=libMesh::out) const
Prints information about the equation systems, by default to libMesh::out.
 
This is the EquationSystems class.
 
T & set(const std::string &)
 
void resize(const unsigned int n)
Resize the vector.
 
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.
 
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
 
This class handles the numbering of degrees of freedom on a mesh.
 
void print_info(std::ostream &os=libMesh::out) const
Prints relevant information about the mesh.
 
void write(const std::string &name, const XdrMODE, const unsigned int write_flags=(WRITE_DATA), bool partition_agnostic=true) const
Write the systems to disk using the XDR data format.
 
virtual void write(const std::string &fname) override
This method implements writing a mesh to a specified file.
 
const DofMap & get_dof_map() const
 
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
 
Manages consistently variables, degrees of freedom, coefficient vectors, matrices and linear solvers ...
 
const T & get(const std::string &) const
 
virtual void find_neighbors(const bool reset_remote_elements=false, const bool reset_current_list=true)=0
Locate element face (edge in 2D) neighbors.
 
Parameters parameters
Data structure holding arbitrary parameters.