Go to the source code of this file.
◆ assemble_poisson() [1/2]
      
        
          | void assemble_poisson  | 
          ( | 
          EquationSystems &  | 
          es,  | 
        
        
           | 
           | 
          const std::string &  | 
          libmesh_dbg_varsystem_name  | 
        
        
           | 
          ) | 
           |  | 
        
      
 
Definition at line 332 of file introduction_ex4.C.
  337   libmesh_assert_equal_to (system_name, 
"Poisson");
 
  343   PerfLog perf_log (
"Matrix Assembly");
 
  362   FEType fe_type = dof_map.variable_type(0);
 
  368   std::unique_ptr<FEBase> fe (FEBase::build(
dim, fe_type));
 
  374   fe->attach_quadrature_rule (&qrule);
 
  378   std::unique_ptr<FEBase> fe_face (FEBase::build(
dim, fe_type));
 
  387   fe_face->attach_quadrature_rule (&qface);
 
  393   const std::vector<Real> & JxW = fe->get_JxW();
 
  398   const std::vector<Point> & q_point = fe->get_xyz();
 
  401   const std::vector<std::vector<Real>> & phi = fe->get_phi();
 
  405   const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
 
  417   std::vector<dof_id_type> dof_indices;
 
  428       perf_log.push(
"elem init");
 
  434       dof_map.dof_indices (elem, dof_indices);
 
  441       const unsigned int n_dofs =
 
  442         cast_int<unsigned int>(dof_indices.size());
 
  452       libmesh_assert_equal_to (n_dofs, phi.size());
 
  460       Ke.
resize (n_dofs, n_dofs);
 
  467       perf_log.pop(
"elem init");
 
  478       perf_log.push (
"Ke");
 
  480       for (
unsigned int qp=0; qp<qrule.n_points(); qp++)
 
  481         for (
unsigned int i=0; i != n_dofs; i++)
 
  482           for (
unsigned int j=0; j != n_dofs; j++)
 
  483             Ke(i,j) += JxW[qp]*(dphi[i][qp]*dphi[j][qp]);
 
  494       perf_log.push (
"Fe");
 
  496       for (
unsigned int qp=0; qp<qrule.n_points(); qp++)
 
  512           const Real x = q_point[qp](0);
 
  514           const Real y = q_point[qp](1);
 
  519           const Real z = q_point[qp](2);
 
  523           const Real eps = 1.e-3;
 
  543               fxy = (0.25*
pi*
pi)*sin(.5*
pi*x);
 
  547               fxy = - (uxx + uyy + ((
dim==2) ? 0. : uzz));
 
  551           for (
unsigned int i=0; i != n_dofs; i++)
 
  552             Fe(i) += JxW[qp]*fxy*phi[i][qp];
 
  562       dof_map.heterogenously_constrain_element_matrix_and_vector (Ke, Fe, dof_indices);
 
  570       LOG_SCOPE_WITH(
"matrix insertion", 
"", perf_log);
 
 
References libMesh::MeshBase::active_local_element_ptr_range(), libMesh::SparseMatrix< T >::add_matrix(), libMesh::NumericVector< T >::add_vector(), libMesh::FEGenericBase< OutputType >::build(), dim, exact_solution(), libMesh::FIFTH, libMesh::System::get_dof_map(), libMesh::EquationSystems::get_mesh(), libMesh::EquationSystems::get_system(), libMesh::ImplicitSystem::matrix, mesh, libMesh::MeshBase::mesh_dimension(), libMesh::QBase::n_points(), libMesh::pi, libMesh::PerfLog::pop(), libMesh::PerfLog::push(), libMesh::Real, libMesh::DenseVector< T >::resize(), libMesh::DenseMatrix< T >::resize(), and libMesh::ExplicitSystem::rhs.
 
 
◆ assemble_poisson() [2/2]
      
        
          | void assemble_poisson  | 
          ( | 
          EquationSystems &  | 
          es,  | 
        
        
           | 
           | 
          const std::string &  | 
          system_name  | 
        
        
           | 
          ) | 
           |  | 
        
      
 
 
◆ exact_solution()
      
        
          | Real exact_solution  | 
          ( | 
          const Real  | 
          x,  | 
        
        
           | 
           | 
          const Real  | 
          y,  | 
        
        
           | 
           | 
          const Real  | 
          t  | 
        
        
           | 
          ) | 
           |  | 
        
      
 
This is the exact solution that we are trying to obtain. 
We will solve
and take a finite difference approximation using this function to get f. This is the well-known "method of
manufactured solutions". 
Definition at line 43 of file exact_solution.C.
   47   static const Real pi = acos(-1.);
 
   49   return cos(.5*
pi*x)*sin(.5*
pi*y)*cos(.5*
pi*z);
 
 
Referenced by assemble_poisson(), and exact_solution_wrapper().
 
 
◆ exact_solution_wrapper()
      
        
          | void exact_solution_wrapper  | 
          ( | 
          DenseVector< Number > &  | 
          output,  | 
        
        
           | 
           | 
          const Point &  | 
          p,  | 
        
        
           | 
           | 
          const  | 
          Real  | 
        
        
           | 
          ) | 
           |  | 
        
      
 
 
◆ main()
      
        
          | int main  | 
          ( | 
          int  | 
          argc,  | 
        
        
           | 
           | 
          char **  | 
          argv  | 
        
        
           | 
          ) | 
           |  | 
        
      
 
Definition at line 117 of file introduction_ex4.C.
  124                            "--enable-petsc, --enable-trilinos, or --enable-eigen");
 
  130   GetPot command_line (argc, argv);
 
  137       libmesh_error_msg(
"Usage:\n" << 
"\t " << argv[0] << 
" -d 2(3)" << 
" -n 15");
 
  146       for (
int i=1; i<argc; i++)
 
  156   if (command_line.search(1, 
"-d"))
 
  157     dim = command_line.next(
dim);
 
  160   libmesh_example_requires(
dim <= LIBMESH_DIM, 
"2D/3D support");
 
  163 #ifndef LIBMESH_ENABLE_DIRICHLET 
  164   libmesh_example_requires(
false, 
"--enable-dirichlet");
 
  170   if (command_line.search(1, 
"-n"))
 
  171     ps = command_line.next(ps);
 
  174   std::string order = 
"SECOND";
 
  175   if (command_line.search(2, 
"-Order", 
"-o"))
 
  176     order = command_line.next(order);
 
  179   std::string family = 
"LAGRANGE";
 
  180   if (command_line.search(2, 
"-FEFamily", 
"-f"))
 
  181     family = command_line.next(family);
 
  184   if ((family == 
"MONOMIAL") || (family == 
"XYZ"))
 
  185     libmesh_error_msg(
"ex4 currently requires a C^0 (or higher) FE basis.");
 
  197   Real halfwidth = 
dim > 1 ? 1. : 0.;
 
  198   Real halfheight = 
dim > 2 ? 1. : 0.;
 
  200   if ((family == 
"LAGRANGE") && (order == 
"FIRST"))
 
  209                                          -halfwidth, halfwidth,
 
  210                                          -halfheight, halfheight,
 
  222                                          -halfwidth, halfwidth,
 
  223                                          -halfheight, halfheight,
 
  245                                            Utility::string_to_enum<Order>   (order),
 
  246                                            Utility::string_to_enum<FEFamily>(family));
 
  257   std::set<boundary_id_type> boundary_ids;
 
  259   boundary_ids.insert(0);
 
  260   boundary_ids.insert(1);
 
  264       boundary_ids.insert(2);
 
  265       boundary_ids.insert(3);
 
  270       boundary_ids.insert(4);
 
  271       boundary_ids.insert(5);
 
  275   std::vector<unsigned int> variables(1);
 
  276   variables[0] = u_var;
 
  282 #ifdef LIBMESH_ENABLE_DIRICHLET 
  289     (boundary_ids, variables, exact_solution_object);
 
  297   equation_systems.init();
 
  300   equation_systems.print_info();
 
  310       GnuPlotIO plot(
mesh, 
"Introduction Example 4, 1D", GnuPlotIO::GRID_ON);
 
  311       plot.write_equation_systems(
"gnuplot_script", equation_systems);
 
  313 #ifdef LIBMESH_HAVE_EXODUS_API 
  317                                                  "out_3.e" : 
"out_2.e", equation_systems);
 
  319 #endif // #ifdef LIBMESH_HAVE_EXODUS_API 
 
References libMesh::DofMap::add_dirichlet_boundary(), libMesh::EquationSystems::add_system(), libMesh::System::add_variable(), assemble_poisson(), libMesh::System::attach_assemble_function(), libMesh::MeshTools::Generation::build_cube(), libMesh::default_solver_package(), dim, libMesh::EDGE2, libMesh::EDGE3, exact_solution_wrapper(), libMesh::System::get_dof_map(), libMesh::GnuPlotIO::GRID_ON, libMesh::HEX27, libMesh::HEX8, libMesh::TriangleWrapper::init(), libMesh::EquationSystems::init(), libMesh::INVALID_SOLVER_PACKAGE, mesh, libMesh::out, libMesh::EquationSystems::print_info(), libMesh::MeshBase::print_info(), libMesh::QUAD4, libMesh::QUAD9, libMesh::Real, libMesh::LinearImplicitSystem::solve(), and libMesh::MeshOutput< MT >::write_equation_systems().
 
 
 
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 add_dirichlet_boundary(const DirichletBoundary &dirichlet_boundary)
Adds a copy of the specified Dirichlet boundary to the system.
 
virtual SimpleRange< element_iterator > active_local_element_ptr_range()=0
 
const T_sys & get_system(const std::string &name) const
 
virtual void solve() override
Assembles & solves the linear system A*x=b.
 
SolverPackage default_solver_package()
 
unsigned int mesh_dimension() const
 
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
 
Real exact_solution(const Real x, const Real y, const Real z=0.)
This is the exact solution that we are trying to obtain.
 
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.
 
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.
 
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.
 
Wraps a function pointer into a FunctionBase object.
 
The PerfLog class allows monitoring of specific events.
 
SparseMatrix< Number > * matrix
The system matrix.
 
void assemble_poisson(EquationSystems &es, const std::string &system_name)
 
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
 
void exact_solution_wrapper(DenseVector< Number > &output, const Point &p, const Real)
 
This is the EquationSystems class.
 
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 ...
 
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.
 
This class allows one to associate Dirichlet boundary values with a given set of mesh boundary ids an...
 
const DofMap & get_dof_map() const
 
This class implements writing meshes using GNUplot, designed for use only with 1D meshes.
 
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
 
Manages consistently variables, degrees of freedom, coefficient vectors, matrices and linear solvers ...