Go to the source code of this file.
◆ assemble()
Definition at line 233 of file miscellaneous_ex4.C.
  239 #ifdef LIBMESH_ENABLE_AMR 
  242   libmesh_assert_equal_to (system_name, 
"System");
 
  262   std::unique_ptr<FEBase> fe      (FEBase::build(
dim, fe_type));
 
  263   std::unique_ptr<FEBase> fe_face (FEBase::build(
dim, fe_type));
 
  267   QGauss qrule (
dim,   fe_type.default_quadrature_order());
 
  268   QGauss qface (
dim-1, fe_type.default_quadrature_order());
 
  271   fe->attach_quadrature_rule      (&qrule);
 
  272   fe_face->attach_quadrature_rule (&qface);
 
  277   const std::vector<Real> & JxW      = fe->get_JxW();
 
  278   const std::vector<Real> & JxW_face = fe_face->get_JxW();
 
  281   const std::vector<std::vector<Real>> & phi = fe->get_phi();
 
  282   const std::vector<std::vector<Real>> & psi = fe_face->get_phi();
 
  286   const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
 
  312   std::vector<dof_id_type> dof_indices;
 
  339       const unsigned int n_dofs =
 
  340         cast_int<unsigned int>(dof_indices.size());
 
  342       Ke.
resize (n_dofs, n_dofs);
 
  354       for (
unsigned int qp=0; qp<qrule.n_points(); qp++)
 
  357           for (
unsigned int i=0; i<n_dofs; i++)
 
  360               Fe(i) += JxW[qp]*phi[i][qp];
 
  362               for (
unsigned int j=0; j<n_dofs; j++)
 
  367                                       (dphi[i][qp]*dphi[j][qp])
 
  372               Ve(i) += JxW[qp]*phi[i][qp];
 
  373               We(i) += JxW[qp]*phi[i][qp];
 
  388         const Real penalty = 1.e10;
 
  393         for (
auto s : elem->side_index_range())
 
  394           if (elem->neighbor_ptr(s) == 
nullptr)
 
  396               fe_face->reinit(elem, s);
 
  398               for (
unsigned int qp=0; qp<qface.n_points(); qp++)
 
  401                   for (
unsigned int i=0; i<n_dofs; i++)
 
  402                     for (
unsigned int j=0; j<n_dofs; j++)
 
  403                       Ke(i,j) += penalty*JxW_face[qp]*psi[i][qp]*psi[j][qp];
 
  423       std::vector<dof_id_type> dof_indices_backup(dof_indices);
 
  425       dof_indices = dof_indices_backup;
 
  449 #endif // #ifdef LIBMESH_ENABLE_AMR 
 
References libMesh::MeshBase::active_local_element_ptr_range(), libMesh::SparseMatrix< T >::add_matrix(), libMesh::NumericVector< T >::add_vector(), libMesh::FEGenericBase< OutputType >::build(), libMesh::SparseMatrix< T >::close(), libMesh::NumericVector< T >::close(), libMesh::DofMap::constrain_element_dyad_matrix(), libMesh::DofMap::constrain_element_matrix_and_vector(), dim, libMesh::DofMap::dof_indices(), libMesh::System::get_dof_map(), libMesh::ImplicitSystem::get_matrix(), libMesh::EquationSystems::get_mesh(), libMesh::EquationSystems::get_system(), libMesh::System::get_vector(), libMesh::libmesh_ignore(), libMesh::ImplicitSystem::matrix, mesh, libMesh::MeshBase::mesh_dimension(), libMesh::Real, libMesh::DenseVector< T >::resize(), libMesh::DenseMatrix< T >::resize(), libMesh::ExplicitSystem::rhs, and libMesh::System::variable_type().
Referenced by main().
 
 
◆ main()
      
        
          | int main  | 
          ( | 
          int  | 
          argc,  | 
        
        
           | 
           | 
          char **  | 
          argv  | 
        
        
           | 
          ) | 
           |  | 
        
      
 
Definition at line 90 of file miscellaneous_ex4.C.
   95 #if !defined(LIBMESH_ENABLE_AMR) 
   96   libmesh_example_requires(
false, 
"--enable-amr");
 
  105   for (
int i=1; i<argc; i++)
 
  111   libmesh_example_requires(2 <= LIBMESH_DIM, 
"2D support");
 
  148   equation_systems.init ();
 
  151   equation_systems.print_info();
 
  153   equation_systems.parameters.set<
unsigned int>
 
  154     (
"linear solver maximum iterations") = 250;
 
  155   equation_systems.parameters.set<
Real>
 
  159   for (
unsigned int i=0; i<2; i++)
 
  166               if ((elem->id()%20)>8)
 
  167                 elem->set_refinement_flag(Elem::REFINE);
 
  169                 elem->set_refinement_flag(Elem::DO_NOTHING);
 
  172             elem->set_refinement_flag(Elem::INACTIVE);
 
  174       mesh_refinement.refine_elements();
 
  175       equation_systems.reinit();
 
  179   equation_systems.print_info();
 
  192   shellMatrix.matrices.push_back(&shellMatrix0);
 
  194   shellMatrix.matrices.push_back(&shellMatrix1);
 
  213                << 
" iterations, residual norm is " 
  218 #if defined(LIBMESH_HAVE_VTK) && !defined(LIBMESH_ENABLE_PARMESH) 
  221 #endif // #ifdef LIBMESH_HAVE_VTK 
  223 #endif // #ifndef LIBMESH_ENABLE_AMR 
 
References libMesh::ImplicitSystem::add_matrix(), libMesh::System::add_variable(), libMesh::System::add_vector(), assemble(), libMesh::System::attach_assemble_function(), libMesh::LinearImplicitSystem::attach_shell_matrix(), libMesh::MeshTools::Generation::build_square(), libMesh::ParallelObject::comm(), libMesh::default_solver_package(), libMesh::LinearImplicitSystem::detach_shell_matrix(), libMesh::Elem::DO_NOTHING, libMesh::MeshBase::element_ptr_range(), libMesh::LinearImplicitSystem::final_linear_residual(), libMesh::FIRST, libMesh::ImplicitSystem::get_matrix(), libMesh::System::get_vector(), libMesh::Elem::INACTIVE, libMesh::TriangleWrapper::init(), libMesh::NumericVector< T >::init(), libMesh::ImplicitSystem::matrix, mesh, libMesh::System::n_dofs(), libMesh::LinearImplicitSystem::n_linear_iterations(), libMesh::System::n_local_dofs(), libMesh::out, libMesh::PETSC_SOLVERS, libMesh::QUAD4, libMesh::Real, libMesh::Elem::REFINE, libMesh::MeshRefinement::refine_elements(), libMesh::LinearImplicitSystem::solve(), libMesh::TOLERANCE, libMesh::MeshOutput< MT >::write_equation_systems(), and libMesh::SparseMatrix< T >::zero().
 
 
 
virtual void close()=0
Calls the SparseMatrix's internal assembly routines, ensuring that the values are consistent across p...
 
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
 
const MeshBase & get_mesh() const
 
NumericVector< Number > * rhs
The system matrix.
 
This class implements specific orders of Gauss quadrature.
 
void dof_indices(const Elem *const elem, std::vector< dof_id_type > &di) const
Fills the vector di with the global degree of freedom indices for the element.
 
virtual SimpleRange< element_iterator > active_local_element_ptr_range()=0
 
virtual void close()=0
Calls the NumericVector's internal assembly routines, ensuring that the values are consistent across ...
 
SparseMatrix< Number > & add_matrix(const std::string &mat_name)
Adds the additional matrix mat_name to this system.
 
const T_sys & get_system(const std::string &name) const
 
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.
 
static const Real TOLERANCE
 
const Parallel::Communicator & comm() const
 
virtual void solve() override
Assembles & solves the linear system A*x=b.
 
SolverPackage default_solver_package()
 
unsigned int mesh_dimension() const
 
This class implements reading and writing meshes in the VTK format.
 
Shell matrix that is given by a tensor product of two vectors, i.e.
 
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
 
Implements (adaptive) mesh refinement algorithms for a MeshBase.
 
void attach_assemble_function(void fptr(EquationSystems &es, const std::string &name))
Register a user function to use in assembling the system matrix and RHS.
 
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.
 
void constrain_element_dyad_matrix(DenseVector< Number > &v, DenseVector< Number > &w, std::vector< dof_id_type > &row_dofs, bool asymmetric_constraint_rows=true) const
Constrains a dyadic element matrix B = v w'.
 
dof_id_type n_local_dofs() const
 
void libmesh_ignore(const Args &...)
 
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.
 
void assemble(EquationSystems &es, const std::string &system_name)
 
SparseMatrix< Number > * matrix
The system matrix.
 
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
 
This is the EquationSystems class.
 
void constrain_element_matrix_and_vector(DenseMatrix< Number > &matrix, DenseVector< Number > &rhs, std::vector< dof_id_type > &elem_dofs, bool asymmetric_constraint_rows=true) const
Constrains the element matrix and vector.
 
unsigned int n_linear_iterations() const
 
virtual void write_equation_systems(const std::string &, const EquationSystems &, const std::set< std::string > *system_names=nullptr)
This method implements writing a mesh with data to a specified file where the data is taken from the ...
 
const FEType & variable_type(const unsigned int i) const
 
void resize(const unsigned int n)
Resize the vector.
 
void detach_shell_matrix()
Detaches a shell matrix.
 
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.
 
This class allows to use any SparseMatrix object as a shell matrix.
 
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
 
This class handles the numbering of degrees of freedom on a mesh.
 
NumericVector< Number > & add_vector(const std::string &vec_name, const bool projections=true, const ParallelType type=PARALLEL)
Adds the additional vector vec_name to this system.
 
Real final_linear_residual() const
 
const DofMap & get_dof_map() const
 
dof_id_type n_dofs() const
 
This class combines any number of shell matrices to a single shell matrix by summing them together.
 
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
 
Manages consistently variables, degrees of freedom, coefficient vectors, matrices and linear solvers ...
 
const SparseMatrix< Number > & get_matrix(const std::string &mat_name) const
 
void attach_shell_matrix(ShellMatrix< Number > *shell_matrix)
This function enables the user to provide a shell matrix, i.e.
 
const NumericVector< Number > & get_vector(const std::string &vec_name) const
 
virtual void zero()=0
Set all entries to 0.