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 330 of file subdomains_ex1.C.
335 libmesh_assert_equal_to (system_name,
"Poisson");
341 PerfLog perf_log (
"Matrix Assembly");
360 FEType fe_type = dof_map.variable_type(0);
366 std::unique_ptr<FEBase> fe (FEBase::build(
dim, fe_type));
372 fe->attach_quadrature_rule (&qrule);
376 std::unique_ptr<FEBase> fe_face (FEBase::build(
dim, fe_type));
385 fe_face->attach_quadrature_rule (&qface);
391 const std::vector<Real> & JxW = fe->get_JxW();
396 const std::vector<Point> & q_point = fe->get_xyz();
399 const std::vector<std::vector<Real>> & phi = fe->get_phi();
403 const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
415 std::vector<dof_id_type> dof_indices;
425 if (elem->subdomain_id()==1)
430 perf_log.push(
"elem init");
436 dof_map.dof_indices (elem, dof_indices);
450 Ke.
resize (dof_indices.size(),
453 Fe.
resize (dof_indices.size());
458 perf_log.pop(
"elem init");
469 perf_log.push (
"Ke");
471 for (
unsigned int qp=0; qp<qrule.n_points(); qp++)
472 for (std::size_t i=0; i<phi.size(); i++)
473 for (std::size_t j=0; j<phi.size(); j++)
474 Ke(i,j) += JxW[qp]*(dphi[i][qp]*dphi[j][qp]);
485 perf_log.push (
"Fe");
487 for (
unsigned int qp=0; qp<qrule.n_points(); qp++)
503 const Real x = q_point[qp](0);
505 const Real y = q_point[qp](1);
510 const Real z = q_point[qp](2);
514 const Real eps = 1.e-3;
534 fxy = (0.25*
pi*
pi)*sin(.5*
pi*x);
538 fxy = - (uxx + uyy + ((
dim==2) ? 0. : uzz));
542 for (std::size_t i=0; i<phi.size(); i++)
543 Fe(i) += JxW[qp]*fxy*phi[i][qp];
558 LOG_SCOPE_WITH(
"BCs",
"", perf_log);
566 for (
auto side : elem->side_index_range())
567 if ((elem->neighbor_ptr(side) ==
nullptr) ||
568 (elem->neighbor_ptr(side)->subdomain_id()!=1))
573 const Real penalty = 1.e10;
577 const std::vector<std::vector<Real>> & phi_face = fe_face->get_phi();
581 const std::vector<Real> & JxW_face = fe_face->get_JxW();
586 const std::vector<Point> & qface_point = fe_face->get_xyz();
590 fe_face->reinit(elem, side);
593 for (
unsigned int qp=0; qp<qface.n_points(); qp++)
597 const Real xf = qface_point[qp](0);
599 const Real yf = qface_point[qp](1);
604 const Real zf = qface_point[qp](2);
614 for (std::size_t i=0; i<phi_face.size(); i++)
615 for (std::size_t j=0; j<phi_face.size(); j++)
616 Ke(i,j) += JxW_face[qp]*penalty*phi_face[i][qp]*phi_face[j][qp];
620 for (std::size_t i=0; i<phi_face.size(); i++)
621 Fe(i) += JxW_face[qp]*penalty*
value*phi_face[i][qp];
628 dof_map.constrain_element_matrix_and_vector (Ke, Fe, dof_indices);
636 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(), libMesh::ExplicitSystem::rhs, and value.
◆ 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().
◆ main()
int main |
( |
int |
argc, |
|
|
char ** |
argv |
|
) |
| |
Definition at line 97 of file subdomains_ex1.C.
107 #ifndef LIBMESH_ENABLE_AMR
108 libmesh_example_requires(
false,
"--enable-amr");
115 GetPot command_line (argc, argv);
122 libmesh_error_msg(
"Usage:\n" <<
"\t " << argv[0] <<
" -d 2(3)" <<
" -n 15");
131 for (
int i=1; i<argc; i++)
142 if (command_line.search(1,
"-d"))
143 dim = command_line.next(
dim);
146 libmesh_example_requires(
dim <= LIBMESH_DIM,
"2D/3D support");
151 if (command_line.search(1,
"-n"))
152 ps = command_line.next(ps);
155 std::string order =
"FIRST";
156 if (command_line.search(2,
"-Order",
"-o"))
157 order = command_line.next(order);
160 std::string family =
"LAGRANGE";
161 if (command_line.search(2,
"-FEFamily",
"-f"))
162 family = command_line.next(family);
165 if ((family ==
"MONOMIAL") || (family ==
"XYZ"))
166 libmesh_error_msg(
"This example requires a C^0 (or higher) FE basis.");
178 Real halfwidth =
dim > 1 ? 1. : 0.;
179 Real halfheight =
dim > 2 ? 1. : 0.;
181 if ((family ==
"LAGRANGE") && (order ==
"FIRST"))
190 -halfwidth, halfwidth,
191 -halfheight, halfheight,
203 -halfwidth, halfwidth,
204 -halfheight, halfheight,
225 bool node_in =
false;
226 bool node_out =
false;
227 for (
auto & n : elem->node_ref_range())
235 if (node_in && node_out)
236 elem->set_refinement_flag(Elem::REFINE);
238 elem->set_refinement_flag(Elem::DO_NOTHING);
241 elem->set_refinement_flag(Elem::INACTIVE);
244 meshRefinement.refine_elements();
254 Real d = elem->centroid().norm();
256 elem->subdomain_id() = 1;
271 Utility::string_to_enum<Order> (order),
272 Utility::string_to_enum<FEFamily>(family));
279 equation_systems.init();
282 equation_systems.print_info();
286 std::set<subdomain_id_type> id_list;
297 equation_systems.get_system(
"Poisson").solve();
303 GnuPlotIO plot(
mesh,
"Subdomains Example 1, 1D", GnuPlotIO::GRID_ON);
304 plot.write_equation_systems(
"gnuplot_script", equation_systems);
309 "out_3.gmv" :
"out_2.gmv", equation_systems);
310 #ifdef LIBMESH_HAVE_EXODUS_API
312 "out_3.e" :
"out_2.e", equation_systems);
313 #endif // #ifdef LIBMESH_HAVE_EXODUS_API
316 #endif // #ifndef LIBMESH_ENABLE_AMR
References 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::Elem::DO_NOTHING, libMesh::EDGE2, libMesh::EDGE3, libMesh::MeshBase::element_ptr_range(), libMesh::EquationSystems::get_system(), libMesh::GnuPlotIO::GRID_ON, libMesh::HEX27, libMesh::HEX8, libMesh::Elem::INACTIVE, libMesh::TriangleWrapper::init(), libMesh::EquationSystems::init(), mesh, libMesh::out, libMesh::PETSC_SOLVERS, libMesh::EquationSystems::print_info(), libMesh::MeshBase::print_info(), libMesh::QUAD4, libMesh::QUAD9, libMesh::Real, libMesh::Elem::REFINE, libMesh::MeshRefinement::refine_elements(), libMesh::LinearImplicitSystem::restrict_solve_to(), libMesh::SUBSET_ZERO, and libMesh::MeshOutput< MT >::write_equation_systems().
This class represents a subset of the dofs of a System, selected by the subdomain_id and possible the...
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.
virtual SimpleRange< element_iterator > active_local_element_ptr_range()=0
Real exact_solution(const Real x, const Real y=0., const Real z=0.)
This is the exact solution that we are trying to obtain.
const T_sys & get_system(const std::string &name) const
SolverPackage default_solver_package()
unsigned int mesh_dimension() const
virtual void restrict_solve_to(const SystemSubset *subset, const SubsetSolveMode subset_solve_mode=SUBSET_ZERO) override
After calling this method, any solve will be limited to the given subset.
This class implements writing meshes in the GMV format.
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
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].
void assemble_poisson(EquationSystems &es, const std::string &system_name)
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.
The PerfLog class allows monitoring of specific events.
SparseMatrix< Number > * matrix
The system matrix.
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
Selection of subdomain ids by a list.
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.
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 ...