libMesh
Functions
amr.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)
 
void assemble (EquationSystems &es, const std::string &libmesh_dbg_var(system_name))
 

Function Documentation

◆ assemble() [1/2]

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().

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:196
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:2343
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:2330
void dof_indices(const Elem *const elem, std::vector< dof_id_type > &di) const
Definition: dof_map.C:2164
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:396
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:75
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:179
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:2508
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:372
const DofMap & get_dof_map() const
Definition: system.h:2374
const SparseMatrix< Number > & get_matrix(std::string_view mat_name) const
Definition: system.C:1124
const SparseMatrix< Number > & get_system_matrix() const
const NumericVector< Number > & get_vector(std::string_view vec_name) const
Definition: system.C:943

◆ assemble() [2/2]

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

Definition at line 131 of file amr.C.

References libMesh::SparseMatrix< T >::add_matrix(), libMesh::NumericVector< T >::add_vector(), libMesh::FEGenericBase< OutputType >::build(), libMesh::DofMap::constrain_element_matrix_and_vector(), dim, libMesh::DofMap::dof_indices(), libMesh::FIFTH, libMesh::FIRST, libMesh::System::get_dof_map(), libMesh::EquationSystems::get_mesh(), libMesh::EquationSystems::get_system(), libMesh::ImplicitSystem::get_system_matrix(), mesh, libMesh::MeshBase::mesh_dimension(), libMesh::QBase::n_points(), libMesh::out, libMesh::Real, libMesh::DenseVector< T >::resize(), libMesh::DenseMatrix< T >::resize(), and libMesh::ExplicitSystem::rhs.

133 {
134  libmesh_assert_equal_to (system_name, "primary");
135 
136  const MeshBase & mesh = es.get_mesh();
137  const unsigned int dim = mesh.mesh_dimension();
138 
139  // Also use a 3x3x3 quadrature rule (3D). Then tell the FE
140  // about the geometry of the problem and the quadrature rule
141  FEType fe_type (FIRST);
142 
143  std::unique_ptr<FEBase> fe(FEBase::build(dim, fe_type));
144  QGauss qrule(dim, FIFTH);
145 
146  fe->attach_quadrature_rule (&qrule);
147 
148  std::unique_ptr<FEBase> fe_face(FEBase::build(dim, fe_type));
149  QGauss qface(dim-1, FIFTH);
150 
151  fe_face->attach_quadrature_rule(&qface);
152 
153  LinearImplicitSystem & system =
154  es.get_system<LinearImplicitSystem>("primary");
155 
156 
157  // These are references to cell-specific data
158  const std::vector<Real> & JxW_face = fe_face->get_JxW();
159  const std::vector<Real> & JxW = fe->get_JxW();
160  const std::vector<Point> & q_point = fe->get_xyz();
161  const std::vector<std::vector<Real>> & phi = fe->get_phi();
162  const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
163 
164  std::vector<dof_id_type> dof_indices_U;
165  std::vector<dof_id_type> dof_indices_V;
166  const DofMap & dof_map = system.get_dof_map();
167 
172 
173  Real vol=0., area=0.;
174 
175  SparseMatrix<Number> & matrix = system.get_system_matrix();
176 
177  for (const auto & elem : mesh.active_local_element_ptr_range())
178  {
179  // recompute the element-specific data for the current element
180  fe->reinit (elem);
181 
182 
183  //fe->print_info();
184 
185  dof_map.dof_indices(elem, dof_indices_U, 0);
186  dof_map.dof_indices(elem, dof_indices_V, 1);
187 
188  const unsigned int n_phi = cast_int<unsigned int>(phi.size());
189 
190  // zero the element matrix and vector
191  Kuu.resize (n_phi, n_phi);
192 
193  Kvv.resize (n_phi, n_phi);
194 
195  Fu.resize (n_phi);
196  Fv.resize (n_phi);
197 
198  // standard stuff... like in code 1.
199  for (unsigned int gp=0; gp<qrule.n_points(); gp++)
200  {
201  for (unsigned int i=0; i<n_phi; ++i)
202  {
203  // this is tricky. ig is the _global_ dof index corresponding
204  // to the _global_ vertex number elem->node_id(i). Note that
205  // in general these numbers will not be the same (except for
206  // the case of one unknown per node on one subdomain) so
207  // we need to go through the dof_map
208 
209  const Real f = q_point[gp]*q_point[gp];
210  // const Real f = (q_point[gp](0) +
211  // q_point[gp](1) +
212  // q_point[gp](2));
213 
214  // add jac*weight*f*phi to the RHS in position ig
215 
216  Fu(i) += JxW[gp]*f*phi[i][gp];
217  Fv(i) += JxW[gp]*f*phi[i][gp];
218 
219  for (unsigned int j=0; j != n_phi; ++j)
220  {
221 
222  Kuu(i,j) += JxW[gp]*((phi[i][gp])*(phi[j][gp]));
223 
224  Kvv(i,j) += JxW[gp]*((phi[i][gp])*(phi[j][gp]) +
225  1.*((dphi[i][gp])*(dphi[j][gp])));
226  };
227  };
228  vol += JxW[gp];
229  };
230 
231 
232  // You can't compute "area" (perimeter) if you are in 2D
233  if (dim == 3)
234  {
235  for (auto side : elem->side_index_range())
236  if (elem->neighbor_ptr(side) == nullptr)
237  {
238  fe_face->reinit (elem, side);
239 
240  for (const auto & val : JxW_face)
241  area += val;
242  }
243  }
244 
245  // Constrain the DOF indices.
246  dof_map.constrain_element_matrix_and_vector(Kuu, Fu, dof_indices_U);
247  dof_map.constrain_element_matrix_and_vector(Kvv, Fv, dof_indices_V);
248 
249 
250  system.rhs->add_vector(Fu, dof_indices_U);
251  system.rhs->add_vector(Fv, dof_indices_V);
252 
253  matrix.add_matrix(Kuu, dof_indices_U);
254  matrix.add_matrix(Kvv, dof_indices_V);
255  }
256 
257  libMesh::out << "Vol=" << vol << std::endl;
258 
259  if (dim == 3)
260  libMesh::out << "Area=" << area << std::endl;
261 }
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:196
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:2330
void dof_indices(const Elem *const elem, std::vector< dof_id_type > &di) const
Definition: dof_map.C:2164
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:396
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:75
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:179
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.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
OStreamProxy out
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:372
const DofMap & get_dof_map() const
Definition: system.h:2374
const SparseMatrix< Number > & get_system_matrix() const

◆ main()

int main ( int  argc,
char **  argv 
)

Definition at line 46 of file amr.C.

References libMesh::DofMap::_dof_coupling, libMesh::EquationSystems::add_system(), libMesh::System::add_variable(), assemble(), libMesh::System::attach_assemble_function(), dim, libMesh::MeshBase::elem_ref(), libMesh::FIRST, libMesh::System::get_dof_map(), libMesh::TriangleWrapper::init(), libMesh::EquationSystems::init(), mesh, libMesh::DofMap::print_dof_constraints(), libMesh::EquationSystems::print_info(), libMesh::MeshBase::print_info(), libMesh::MeshBase::read(), libMesh::Elem::REFINE, libMesh::MeshRefinement::refine_and_coarsen_elements(), libMesh::EquationSystems::reinit(), libMesh::CouplingMatrix::resize(), libMesh::Elem::set_refinement_flag(), libMesh::LinearImplicitSystem::solve(), libMesh::MeshRefinement::uniformly_refine(), libMesh::GMVIO::write(), and libMesh::MeshOutput< MT >::write_equation_systems().

47 {
48  LibMeshInit init(argc, argv);
49 
50  libmesh_error_msg_if(argc < 4, "Usage: ./prog -d DIM filename");
51 
52  // Variables to get us started
53  const unsigned char dim = cast_int<unsigned char>(atoi(argv[2]));
54 
55  std::string meshname (argv[3]);
56 
57  // declare a mesh...
58  Mesh mesh(init.comm(), dim);
59 
60  // Read a mesh
61  mesh.read(meshname);
62 
63  GMVIO(mesh).write ("out_0.gmv");
64 
65  mesh.elem_ref(0).set_refinement_flag (Elem::REFINE);
66 
67  MeshRefinement mesh_refinement (mesh);
68 
69  mesh_refinement.refine_and_coarsen_elements ();
70  mesh_refinement.uniformly_refine (2);
71 
72  mesh.print_info();
73 
74 
75  // Set up the equation system(s)
76  EquationSystems es (mesh);
77 
78  LinearImplicitSystem & primary =
79  es.add_system<LinearImplicitSystem>("primary");
80 
81  primary.add_variable ("U", FIRST);
82  primary.add_variable ("V", FIRST);
83 
84  primary.get_dof_map()._dof_coupling->resize(2);
85  (*primary.get_dof_map()._dof_coupling)(0,0) = 1;
86  (*primary.get_dof_map()._dof_coupling)(1,1) = 1;
87 
89 
90  es.init ();
91 
92  es.print_info ();
94 
95  // call the solver.
96  primary.solve ();
97 
98  GMVIO(mesh).write_equation_systems ("out_1.gmv",
99  es);
100 
101 
102 
103  // Refine uniformly
104  mesh_refinement.uniformly_refine (1);
105  es.reinit ();
106 
107  // Write out the projected solution
108  GMVIO(mesh).write_equation_systems ("out_2.gmv",
109  es);
110 
111  // Solve again. Output the refined solution
112  primary.solve ();
113  GMVIO(mesh).write_equation_systems ("out_3.gmv",
114  es);
115 
116  return 0;
117 }
This is the EquationSystems class.
virtual void read(const std::string &name, void *mesh_data=nullptr, bool skip_renumber_nodes_and_elements=false, bool skip_find_neighbors=false)=0
Interfaces for reading/writing a mesh to/from a file.
virtual void write(const std::string &) override
This method implements writing a mesh to a specified file.
Definition: gmv_io.C:271
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
unsigned int dim
Manages consistently variables, degrees of freedom, coefficient vectors, matrices and linear solvers ...
MeshBase & mesh
void set_refinement_flag(const RefinementState rflag)
Sets the value of the refinement flag for the element.
Definition: elem.h:3218
void resize(const std::size_t n)
Resizes the matrix and initializes all entries to be 0.
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
Definition: libmesh.h:90
This class implements writing meshes in the GMV format.
Definition: gmv_io.h:46
virtual void solve() override
Assembles & solves the linear system A*x=b.
void assemble(EquationSystems &es, const std::string &system_name)
Implements (adaptive) mesh refinement algorithms for a MeshBase.
void print_info(std::ostream &os=libMesh::out, const unsigned int verbosity=0, const bool global=true) const
Prints relevant information about the mesh.
Definition: mesh_base.C:1562
void print_dof_constraints(std::ostream &os=libMesh::out, bool print_nonlocal=false) const
Prints (from processor 0) all DoF and Node constraints.
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:1357
CouplingMatrix * _dof_coupling
Degree of freedom coupling.
Definition: dof_map.h:1613
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:2161
virtual const Elem & elem_ref(const dof_id_type i) const
Definition: mesh_base.h:639
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
Definition: mesh.h:50
const DofMap & get_dof_map() const
Definition: system.h:2374