libMesh
Functions
miscellaneous_ex4.C File Reference

Go to the source code of this file.

Functions

void assemble (EquationSystems &es, const std::string &system_name)
 
int main (int argc, char **argv)
 

Function Documentation

◆ assemble()

void assemble ( EquationSystems es,
const std::string &  system_name 
)

Definition at line 239 of file miscellaneous_ex4.C.

References libMesh::SparseMatrix< T >::add_matrix(), libMesh::NumericVector< T >::add_vector(), libMesh::FEGenericBase< OutputType >::build(), libMesh::NumericVector< T >::close(), libMesh::SparseMatrix< 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::System::get_matrix(), libMesh::EquationSystems::get_mesh(), libMesh::EquationSystems::get_system(), libMesh::ImplicitSystem::get_system_matrix(), libMesh::System::get_vector(), libMesh::libmesh_ignore(), 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(), HDGProblem::precheck(), and HDGProblem::residual_and_jacobian().

241 {
242  // Ignore unused parameter warnings when !LIBMESH_ENABLE_AMR.
243  libmesh_ignore(es, system_name);
244 
245 #ifdef LIBMESH_ENABLE_AMR
246  // It is a good idea to make sure we are assembling
247  // the proper system.
248  libmesh_assert_equal_to (system_name, "System");
249 
250  // Get a constant reference to the mesh object.
251  const MeshBase & mesh = es.get_mesh();
252 
253  // The dimension that we are running
254  const unsigned int dim = mesh.mesh_dimension();
255 
256  // Get a reference to the Convection-Diffusion system object.
257  LinearImplicitSystem & system =
258  es.get_system<LinearImplicitSystem> ("System");
259 
260  // Get the Finite Element type for the first (and only)
261  // variable in the system.
262  FEType fe_type = system.variable_type(0);
263 
264  // Build a Finite Element object of the specified type. Since the
265  // FEBase::build() member dynamically creates memory we will
266  // store the object as a std::unique_ptr<FEBase>. This can be thought
267  // of as a pointer that will clean up after itself.
268  std::unique_ptr<FEBase> fe (FEBase::build(dim, fe_type));
269  std::unique_ptr<FEBase> fe_face (FEBase::build(dim, fe_type));
270 
271  // A Gauss quadrature rule for numerical integration.
272  // Let the FEType object decide what order rule is appropriate.
273  QGauss qrule (dim, fe_type.default_quadrature_order());
274  QGauss qface (dim-1, fe_type.default_quadrature_order());
275 
276  // Tell the finite element object to use our quadrature rule.
277  fe->attach_quadrature_rule (&qrule);
278  fe_face->attach_quadrature_rule (&qface);
279 
280  // Here we define some references to cell-specific data that
281  // will be used to assemble the linear system. We will start
282  // with the element Jacobian * quadrature weight at each integration point.
283  const std::vector<Real> & JxW = fe->get_JxW();
284  const std::vector<Real> & JxW_face = fe_face->get_JxW();
285 
286  // The element shape functions evaluated at the quadrature points.
287  const std::vector<std::vector<Real>> & phi = fe->get_phi();
288  const std::vector<std::vector<Real>> & psi = fe_face->get_phi();
289 
290  // The element shape function gradients evaluated at the quadrature
291  // points.
292  const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
293 
294  // The XY locations of the quadrature points used for face integration
295  //const std::vector<Point>& qface_points = fe_face->get_xyz();
296 
297  // A reference to the DofMap object for this system. The DofMap
298  // object handles the index translation from node and element numbers
299  // to degree of freedom numbers. We will talk more about the DofMap
300  // in future examples.
301  const DofMap & dof_map = system.get_dof_map();
302 
303  // Define data structures to contain the element matrix
304  // and right-hand-side vector contribution. Following
305  // basic finite element terminology we will denote these
306  // "Ke" and "Fe".
309 
310  // Analogous data structures for thw two vectors v and w that form
311  // the tensor shell matrix.
314 
315  // This vector will hold the degree of freedom indices for
316  // the element. These define where in the global system
317  // the element degrees of freedom get mapped.
318  std::vector<dof_id_type> dof_indices;
319 
320  // The global system matrix
321  SparseMatrix<Number> & matrix = system.get_system_matrix();
322 
323  // Now we will loop over all the elements in the mesh that
324  // live on the local processor. We will compute the element
325  // matrix and right-hand-side contribution. Since the mesh
326  // will be refined we want to only consider the ACTIVE elements,
327  // hence we use a variant of the active_elem_iterator.
328  for (const auto & elem : mesh.active_local_element_ptr_range())
329  {
330  // Get the degree of freedom indices for the
331  // current element. These define where in the global
332  // matrix and right-hand-side this element will
333  // contribute to.
334  dof_map.dof_indices (elem, dof_indices);
335 
336  // Compute the element-specific data for the current
337  // element. This involves computing the location of the
338  // quadrature points (q_point) and the shape functions
339  // (phi, dphi) for the current element.
340  fe->reinit (elem);
341 
342  // Zero the element matrix and right-hand side before
343  // summing them. We use the resize member here because
344  // the number of degrees of freedom might have changed from
345  // the last element. Note that this will be the case if the
346  // element type is different (i.e. the last element was a
347  // triangle, now we are on a quadrilateral).
348  const unsigned int n_dofs =
349  cast_int<unsigned int>(dof_indices.size());
350 
351  Ke.resize (n_dofs, n_dofs);
352 
353  Fe.resize (n_dofs);
354  Ve.resize (n_dofs);
355  We.resize (n_dofs);
356 
357  // Now we will build the element matrix and right-hand-side.
358  // Constructing the RHS requires the solution and its
359  // gradient from the previous timestep. This myst be
360  // calculated at each quadrature point by summing the
361  // solution degree-of-freedom values by the appropriate
362  // weight functions.
363  for (unsigned int qp=0; qp<qrule.n_points(); qp++)
364  {
365  // Now compute the element matrix and RHS contributions.
366  for (unsigned int i=0; i<n_dofs; i++)
367  {
368  // The RHS contribution
369  Fe(i) += JxW[qp]*phi[i][qp];
370 
371  for (unsigned int j=0; j<n_dofs; j++)
372  {
373  // The matrix contribution
374  Ke(i,j) += JxW[qp]*(
375  // Stiffness matrix
376  (dphi[i][qp]*dphi[j][qp])
377  );
378  }
379 
380  // V and W are the same for this example.
381  Ve(i) += JxW[qp]*phi[i][qp];
382  We(i) += JxW[qp]*phi[i][qp];
383  }
384  }
385 
386  // At this point the interior element integration has
387  // been completed. However, we have not yet addressed
388  // boundary conditions. For this example we will only
389  // consider simple Dirichlet boundary conditions imposed
390  // via the penalty method.
391  //
392  // The following loops over the sides of the element.
393  // If the element has no neighbor on a side then that
394  // side MUST live on a boundary of the domain.
395  {
396  // The penalty value.
397  const Real penalty = 1.e10;
398 
399  // The following loops over the sides of the element.
400  // If the element has no neighbor on a side then that
401  // side MUST live on a boundary of the domain.
402  for (auto s : elem->side_index_range())
403  if (elem->neighbor_ptr(s) == nullptr)
404  {
405  fe_face->reinit(elem, s);
406 
407  for (unsigned int qp=0; qp<qface.n_points(); qp++)
408  {
409  // Matrix contribution
410  for (unsigned int i=0; i<n_dofs; i++)
411  for (unsigned int j=0; j<n_dofs; j++)
412  Ke(i,j) += penalty*JxW_face[qp]*psi[i][qp]*psi[j][qp];
413  }
414  }
415  }
416 
417 
418  // We have now built the element matrix and RHS vector in terms
419  // of the element degrees of freedom. However, it is possible
420  // that some of the element DOFs are constrained to enforce
421  // solution continuity, i.e. they are not really "free". We need
422  // to constrain those DOFs in terms of non-constrained DOFs to
423  // ensure a continuous solution. The
424  // DofMap::constrain_element_matrix_and_vector() method does
425  // just that.
426 
427  // However, constraining both the sparse matrix (and right hand
428  // side) plus the rank 1 matrix is tricky. The dof_indices
429  // vector has to be backed up for that because the constraining
430  // functions modify it.
431 
432  std::vector<dof_id_type> dof_indices_backup(dof_indices);
433  dof_map.constrain_element_matrix_and_vector (Ke, Fe, dof_indices);
434  dof_indices = dof_indices_backup;
435  dof_map.constrain_element_dyad_matrix(Ve, We, dof_indices);
436 
437  // The element matrix and right-hand-side are now built
438  // for this element. Add them to the global matrix and
439  // right-hand-side vector. The SparseMatrix::add_matrix()
440  // and NumericVector::add_vector() members do this for us.
441  matrix.add_matrix (Ke, dof_indices);
442  system.get_matrix("Preconditioner").add_matrix (Ke, dof_indices);
443  system.rhs->add_vector (Fe, dof_indices);
444  system.get_vector("v").add_vector(Ve, dof_indices);
445  system.get_vector("w").add_vector(We, dof_indices);
446  }
447  // Finished computing the system matrix and right-hand side.
448 
449  // Matrices and vectors must be closed manually. This is necessary
450  // because the matrix is not directly used as the system matrix (in
451  // which case the solver closes it) but as a part of a shell matrix.
452  matrix.close();
453  system.get_matrix("Preconditioner").close();
454  system.rhs->close();
455  system.get_vector("v").close();
456  system.get_vector("w").close();
457 
458 #endif // #ifdef LIBMESH_ENABLE_AMR
459 }
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:182
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&#39;.
Definition: dof_map.h:2267
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.
Definition: dof_map.h:2254
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.
Definition: dof_map.C:1992
unsigned int dim
Manages consistently variables, degrees of freedom, coefficient vectors, matrices and linear solvers ...
void resize(const unsigned int n)
Resize the vector.
Definition: dense_vector.h:374
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]...
MeshBase & mesh
NumericVector< Number > * rhs
The system matrix.
const T_sys & get_system(std::string_view name) const
This is the MeshBase class.
Definition: mesh_base.h:74
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:169
void libmesh_ignore(const Args &...)
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.
virtual void close()=0
Calls the NumericVector&#39;s internal assembly routines, ensuring that the values are consistent across ...
virtual void close()=0
Calls the SparseMatrix&#39;s internal assembly routines, ensuring that the values are consistent across p...
const FEType & variable_type(const unsigned int i) const
Definition: system.h:2427
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().
Definition: dense_matrix.h:895
This class implements specific orders of Gauss quadrature.
unsigned int mesh_dimension() const
Definition: mesh_base.C:324
const DofMap & get_dof_map() const
Definition: system.h:2293
const SparseMatrix< Number > & get_matrix(std::string_view mat_name) const
Definition: system.C:1073
const SparseMatrix< Number > & get_system_matrix() const
const NumericVector< Number > & get_vector(std::string_view vec_name) const
Definition: system.C:918

◆ main()

int main ( int  argc,
char **  argv 
)

Definition at line 90 of file miscellaneous_ex4.C.

References libMesh::System::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::LinearImplicitSystem::final_linear_residual(), libMesh::FIRST, libMesh::System::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().

91 {
92  // Initialize libMesh.
93  LibMeshInit init (argc, argv);
94 
95 #if !defined(LIBMESH_ENABLE_AMR)
96  libmesh_example_requires(false, "--enable-amr");
97 #else
98  libmesh_example_requires(libMesh::default_solver_package() == PETSC_SOLVERS, "--enable-petsc");
99 
100  // Brief message to the user regarding the program name
101  // and command line arguments.
102 
103  libMesh::out << "Running: " << argv[0];
104 
105  for (int i=1; i<argc; i++)
106  libMesh::out << " " << argv[i];
107 
108  libMesh::out << std::endl << std::endl;
109 
110  // Skip this 2D example if libMesh was compiled as 1D-only.
111  libmesh_example_requires(2 <= LIBMESH_DIM, "2D support");
112 
113  // Create a mesh, with dimension to be overridden later, distributed
114  // across the default MPI communicator.
115  Mesh mesh(init.comm());
116 
117  // Create an equation systems object.
118  EquationSystems equation_systems (mesh);
119 
120  // Get command line arguments for mesh size
121  GetPot input(argc, argv);
122 
123  const unsigned int nx = input("nx", 16);
124  const unsigned int ny = input("ny", 16);
125 
127  nx,
128  ny,
129  -1., 1.,
130  -1., 1.,
131  QUAD4);
132 
133  LinearImplicitSystem & system =
134  equation_systems.add_system<LinearImplicitSystem>
135  ("System");
136 
137  // Adds the variable "u" to "System". "u"
138  // will be approximated using first-order approximation.
139  system.add_variable ("u", FIRST);
140 
141  // Also, we need to add two vectors. The tensor matrix v*w^T of
142  // these two vectors will be part of the system matrix.
143  system.add_vector("v", false);
144  system.add_vector("w", false);
145 
146  // We need an additional matrix to be used for preconditioning since
147  // a shell matrix is not suitable for that.
148  system.add_matrix("Preconditioner");
149 
150  // Give the system a pointer to the matrix assembly function.
152 
153  // Initialize the data structures for the equation system.
154  equation_systems.init ();
155 
156  // Prints information about the system to the screen.
157  equation_systems.print_info();
158 
159  equation_systems.parameters.set<unsigned int>
160  ("linear solver maximum iterations") = 250;
161  equation_systems.parameters.set<Real>
162  ("linear solver tolerance") = TOLERANCE;
163 
164  // Refine arbitrarily some elements.
165  for (unsigned int i=0; i<2; i++)
166  {
167  MeshRefinement mesh_refinement(mesh);
168  for (auto & elem : mesh.element_ptr_range())
169  {
170  if (elem->active())
171  {
172  if ((elem->id()%20)>8)
173  elem->set_refinement_flag(Elem::REFINE);
174  else
175  elem->set_refinement_flag(Elem::DO_NOTHING);
176  }
177  else
178  elem->set_refinement_flag(Elem::INACTIVE);
179  }
180  mesh_refinement.refine_elements();
181  equation_systems.reinit();
182  }
183 
184  // Prints information about the system to the screen.
185  equation_systems.print_info();
186 
187  // Before assembling the matrix, we have to clear the two
188  // vectors that form the tensor matrix (since this is not performed
189  // automatically).
190  system.get_vector("v").init(system.n_dofs(), system.n_local_dofs());
191  system.get_vector("w").init(system.n_dofs(), system.n_local_dofs());
192 
193  // We need a shell matrix to solve. There is currently no way to
194  // store the shell matrix in the system. We just create it locally
195  // here (a shell matrix does not occupy much memory).
196  SumShellMatrix<Number> shellMatrix(system.comm());
197  TensorShellMatrix<Number> shellMatrix0(system.get_vector("v"), system.get_vector("w"));
198  shellMatrix.matrices.push_back(&shellMatrix0);
199  SparseShellMatrix<Number> shellMatrix1(*system.matrix);
200  shellMatrix.matrices.push_back(&shellMatrix1);
201 
202  // Attach that to the system.
203  system.attach_shell_matrix(&shellMatrix);
204 
205  // Reset the preconditioning matrix to zero (for the system matrix,
206  // the same thing is done automatically).
207  system.get_matrix("Preconditioner").zero();
208 
209  // Assemble & solve the linear system
210  system.solve();
211 
212  // Detach the shell matrix from the system since it will go out of
213  // scope. Nobody should solve the system outside this function.
214  system.detach_shell_matrix();
215 
216  // Print a nice message.
217  libMesh::out << "Solved linear system in "
218  << system.n_linear_iterations()
219  << " iterations, residual norm is "
220  << system.final_linear_residual()
221  << "."
222  << std::endl;
223 
224 #if defined(LIBMESH_HAVE_VTK) && !defined(LIBMESH_ENABLE_PARMESH)
225  // Write result to file.
226  VTKIO(mesh).write_equation_systems ("out.pvtu", equation_systems);
227 #endif // #ifdef LIBMESH_HAVE_VTK
228 
229 #endif // #ifndef LIBMESH_ENABLE_AMR
230 
231  return 0;
232 }
void attach_shell_matrix(ShellMatrix< Number > *shell_matrix)
This function enables the user to provide a shell matrix, i.e.
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 ...
Definition: mesh_output.C:31
static constexpr Real TOLERANCE
Manages consistently variables, degrees of freedom, coefficient vectors, matrices and linear solvers ...
MeshBase & mesh
const Parallel::Communicator & comm() const
NumericVector< Number > & add_vector(std::string_view vec_name, const bool projections=true, const ParallelType type=PARALLEL)
Adds the additional vector vec_name to this system.
Definition: system.C:751
This class combines any number of shell matrices to a single shell matrix by summing them together...
void build_square(UnstructuredMesh &mesh, const unsigned int nx, const unsigned int ny, const Real xmin=0., const Real xmax=1., const Real ymin=0., const Real ymax=1., const ElemType type=INVALID_ELEM, const bool gauss_lobatto_grid=false)
A specialized build_cube() for 2D meshes.
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
Definition: libmesh.h:90
dof_id_type n_local_dofs() const
Definition: system.C:150
This class implements reading and writing meshes in the VTK format.
Definition: vtk_io.h:60
virtual void solve() override
Assembles & solves the linear system A*x=b.
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.
void assemble(EquationSystems &es, const std::string &system_name)
dof_id_type n_dofs() const
Definition: system.C:113
SolverPackage default_solver_package()
Definition: libmesh.C:1050
Shell matrix that is given by a tensor product of two vectors, i.e.
Implements (adaptive) mesh refinement algorithms for a MeshBase.
virtual void zero()=0
Set all entries to 0.
unsigned int n_linear_iterations() const
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
unsigned int add_variable(std::string_view 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.
Definition: system.C:1305
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.
Definition: system.C:2109
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
OStreamProxy out
SparseMatrix< Number > * matrix
The system matrix.
This class allows to use any SparseMatrix object as a shell matrix.
void detach_shell_matrix()
Detaches a shell matrix.
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
Definition: mesh.h:50
SparseMatrix< Number > & add_matrix(std::string_view mat_name, ParallelType type=PARALLEL, MatrixBuildType mat_build_type=MatrixBuildType::AUTOMATIC)
Adds the additional matrix mat_name to this system.
Definition: system.C:985
const SparseMatrix< Number > & get_matrix(std::string_view mat_name) const
Definition: system.C:1073
const NumericVector< Number > & get_vector(std::string_view vec_name) const
Definition: system.C:918