libMesh
Functions
miscellaneous_ex2.C File Reference

Go to the source code of this file.

Functions

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

Function Documentation

◆ add_M_C_K_helmholtz()

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

Definition at line 502 of file miscellaneous_ex2.C.

504 {
505  START_LOG("init phase", "add_M_C_K_helmholtz");
506 
507  // Verify that we are assembling the system we think we are.
508  libmesh_assert_equal_to (system_name, "Helmholtz");
509 
510  // Get a reference to the FrequencySystem.
511  FrequencySystem & f_system =
512  es.get_system<FrequencySystem> (system_name);
513 
514  // Get the frequency, fluid density, and sound speed of the current solve.
515  const Number frequency = es.parameters.get<Number> ("current frequency");
516  const Real rho = es.parameters.get<Real> ("rho");
517  const Real speed = es.parameters.get<Real> ("wave speed");
518 
519  // Compute angular frequency omega and wave number k
520  const Number omega = 2.0*libMesh::pi*frequency;
521  const Number k = omega / speed;
522 
523  // Get writable references to the system matrix and vector, where the
524  // frequency-dependent system is to be collected.
525  SparseMatrix<Number> & matrix = *f_system.matrix;
526  NumericVector<Number> & rhs = *f_system.rhs;
527 
528  // Get writable references to the frequency-independent matrices
529  // and rhs, though we only need to extract values. This write access
530  // is necessary, since solver packages have to close the data structure
531  // before they can extract values for computation.
532  SparseMatrix<Number> & stiffness = f_system.get_matrix("stiffness");
533  SparseMatrix<Number> & damping = f_system.get_matrix("damping");
534  SparseMatrix<Number> & mass = f_system.get_matrix("mass");
535  NumericVector<Number> & freq_indep_rhs = f_system.get_vector("rhs");
536 
537  Number unit_i (0,1);
538 
539  // Compute the scale values for the different operators.
540  const Number scale_stiffness (1., 0.);
541  const Number scale_damping=unit_i*omega; // I is imaginary unit (from complex.h)
542  const Number scale_mass=-k*k;
543  const Number scale_rhs=-unit_i*rho*omega;
544 
545  // Now simply add the matrices together, store the result
546  // in matrix and rhs. Clear them first.
547  matrix.close ();
548  matrix.zero ();
549  rhs.close ();
550  rhs.zero ();
551 
552  // The matrices from which values are added to another matrix
553  // have to be closed. The add() method does take care of
554  // that, but let us do it explicitly.
555  stiffness.close ();
556  damping.close ();
557  mass.close ();
558 
559  STOP_LOG("init phase", "add_M_C_K_helmholtz");
560 
561  START_LOG("global matrix & vector additions", "add_M_C_K_helmholtz");
562 
563  // Add the stiffness and mass with the proper frequency to the
564  // overall system. For this to work properly, the matrix has
565  // to be not only initialized, but filled with the identical
566  // sparsity structure as the matrix added to it, otherwise
567  // solver packages like PETSc crash.
568  //
569  // Note that we have to add the mass and stiffness contributions one
570  // at a time, otherwise the real part of matrix would be fine, but
571  // the imaginary part would be cluttered with unwanted products.
572  matrix.add (scale_stiffness, stiffness);
573  matrix.add (scale_mass, mass);
574  matrix.add (scale_damping, damping);
575  rhs.add (scale_rhs, freq_indep_rhs);
576 
577  STOP_LOG("global matrix & vector additions", "add_M_C_K_helmholtz");
578 
579  // The linear system involving "matrix" and "rhs" is now ready to be solved.
580 }

References libMesh::SparseMatrix< T >::add(), libMesh::NumericVector< T >::add(), libMesh::SparseMatrix< T >::close(), libMesh::NumericVector< T >::close(), libMesh::Parameters::get(), libMesh::ImplicitSystem::get_matrix(), libMesh::EquationSystems::get_system(), libMesh::System::get_vector(), libMesh::ImplicitSystem::matrix, libMesh::EquationSystems::parameters, libMesh::pi, libMesh::Real, libMesh::ExplicitSystem::rhs, libMesh::SparseMatrix< T >::zero(), and libMesh::NumericVector< T >::zero().

Referenced by main().

◆ assemble_helmholtz()

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

Definition at line 293 of file miscellaneous_ex2.C.

295 {
296 
297  // It is a good idea to make sure we are assembling
298  // the proper system.
299  libmesh_assert_equal_to (system_name, "Helmholtz");
300 
301  // Get a constant reference to the mesh object.
302  const MeshBase & mesh = es.get_mesh();
303 
304  // The maximum dimension of the elements stored in the mesh.
305  const unsigned int dim = mesh.mesh_dimension();
306 
307  // Get a reference to our system, as before
308  FrequencySystem & f_system =
309  es.get_system<FrequencySystem> (system_name);
310 
311  // A const reference to the DofMap object for this system. The DofMap
312  // object handles the index translation from node and element numbers
313  // to degree of freedom numbers.
314  const DofMap & dof_map = f_system.get_dof_map();
315 
316  // Get a constant reference to the finite element type
317  // for the first (and only) variable in the system.
318  const FEType & fe_type = dof_map.variable_type(0);
319 
320  // The fluid density is used by the admittance boundary condition.
321  const Real rho = es.parameters.get<Real>("rho");
322 
323  // In this example, we assemble the element matrices into the
324  // additional matrices "stiffness_mass", "damping", and the element
325  // vectors into "rhs". We get writable references to these objects
326  // now.
327  SparseMatrix<Number> & stiffness = f_system.get_matrix("stiffness");
328  SparseMatrix<Number> & damping = f_system.get_matrix("damping");
329  SparseMatrix<Number> & mass = f_system.get_matrix("mass");
330 
331  // Build a finite element object of the specified type. Since the
332  // FEBase::build() member dynamically creates memory we will
333  // store the object as a std::unique_ptr<FEBase>. This can be thought
334  // of as a pointer that will clean up after itself.
335  std::unique_ptr<FEBase> fe (FEBase::build(dim, fe_type));
336 
337  // A 5th-order Gauss quadrature rule for numerical integration.
338  QGauss qrule (dim, FIFTH);
339 
340  // Tell the finite element object to use our quadrature rule.
341  fe->attach_quadrature_rule (&qrule);
342 
343  // The element Jacobian times the quadrature weight at each integration point.
344  const std::vector<Real> & JxW = fe->get_JxW();
345 
346  // The element shape functions evaluated at the quadrature points.
347  const std::vector<std::vector<Real>> & phi = fe->get_phi();
348 
349  // The element shape function gradients evaluated at the quadrature
350  // points.
351  const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
352 
353  // Here we do not assemble directly in the System matrix, but to the
354  // additional matrices "stiffness_mass" and "damping". The same
355  // approach is followed for the right-hand-side vector Fe, which we
356  // will later on store in the additional vector "rhs".
357  // The zero_matrix is used to explicitly induce the same sparsity
358  // structure in the overall matrix. The mass and stiffness matrices
359  // are real-valued. Therefore, it would be possible to store them
360  // as a single complex-valued matrix to save on memory, but we do
361  // not follow this approach here.
362  DenseMatrix<Number> Ke, Ce, Me, zero_matrix;
364 
365  // This vector will hold the degree of freedom indices for
366  // the element. These define where in the global system
367  // the element degrees of freedom get mapped.
368  std::vector<dof_id_type> dof_indices;
369 
370  // Now we will loop over all the elements in the mesh, and compute
371  // the element matrix and right-hand-side contributions.
372  for (const auto & elem : mesh.active_local_element_ptr_range())
373  {
374  // Start logging the element initialization.
375  START_LOG("elem init", "assemble_helmholtz");
376 
377  // Get the degree of freedom indices for the
378  // current element. These define where in the global
379  // matrix and right-hand-side this element will
380  // contribute to.
381  dof_map.dof_indices (elem, dof_indices);
382 
383  // Compute the element-specific data for the current
384  // element. This involves computing the location of the
385  // quadrature points (q_point) and the shape functions
386  // (phi, dphi) for the current element.
387  fe->reinit (elem);
388 
389  // Zero & resize the element matrix and right-hand side before
390  // summing them, with different element types in the mesh this
391  // is quite necessary.
392  {
393  const unsigned int n_dof_indices = dof_indices.size();
394 
395  Ke.resize (n_dof_indices, n_dof_indices);
396  Ce.resize (n_dof_indices, n_dof_indices);
397  Me.resize (n_dof_indices, n_dof_indices);
398  zero_matrix.resize (n_dof_indices, n_dof_indices);
399  Fe.resize (n_dof_indices);
400  }
401 
402  // Stop logging the element initialization.
403  STOP_LOG("elem init", "assemble_helmholtz");
404 
405  // Now loop over the quadrature points. This handles
406  // the numeric integration.
407  START_LOG("stiffness & mass", "assemble_helmholtz");
408 
409  for (unsigned int qp=0; qp<qrule.n_points(); qp++)
410  {
411  // Now we will build the element matrix. This involves
412  // a double loop to integrate the test functions (i) against
413  // the trial functions (j).
414  for (std::size_t i=0; i<phi.size(); i++)
415  for (std::size_t j=0; j<phi.size(); j++)
416  {
417  Ke(i,j) += JxW[qp]*(dphi[i][qp]*dphi[j][qp]);
418  Me(i,j) += JxW[qp]*(phi[i][qp]*phi[j][qp]);
419  }
420  }
421 
422  STOP_LOG("stiffness & mass", "assemble_helmholtz");
423 
424  // Now compute the contribution to the element matrix
425  // (due to mixed boundary conditions) if the current
426  // element lies on the boundary.
427  //
428  // The following loops over the sides of the element.
429  // If the element has no neighbor on a side then that
430  // side MUST live on a boundary of the domain.
431  for (auto side : elem->side_index_range())
432  if (elem->neighbor_ptr(side) == nullptr)
433  {
434  LOG_SCOPE("damping", "assemble_helmholtz");
435 
436  // Declare a special finite element object for
437  // boundary integration.
438  std::unique_ptr<FEBase> fe_face (FEBase::build(dim, fe_type));
439 
440  // Boundary integration requires one quadrature rule,
441  // with dimensionality one less than the dimensionality
442  // of the element.
443  QGauss qface(dim-1, SECOND);
444 
445  // Tell the finite element object to use our
446  // quadrature rule.
447  fe_face->attach_quadrature_rule (&qface);
448 
449  // The value of the shape functions at the quadrature
450  // points.
451  const std::vector<std::vector<Real>> & phi_face =
452  fe_face->get_phi();
453 
454  // The Jacobian times the quadrature weight at the quadrature
455  // points on the face.
456  const std::vector<Real> & JxW_face = fe_face->get_JxW();
457 
458  // Compute the shape function values on the element
459  // face.
460  fe_face->reinit(elem, side);
461 
462  // For the Robin BCs, consider a normal admittance an=1
463  // at some parts of the boundary
464  const Real an_value = 1.;
465 
466  // Loop over the face quadrature points for integration.
467  for (unsigned int qp=0; qp<qface.n_points(); qp++)
468  {
469  // Element matrix contribution due to prescribed
470  // admittance boundary conditions.
471  for (std::size_t i=0; i<phi_face.size(); i++)
472  for (std::size_t j=0; j<phi_face.size(); j++)
473  Ce(i,j) += rho * an_value * JxW_face[qp] * phi_face[i][qp] * phi_face[j][qp];
474  }
475  }
476 
477  // If this assembly program were to be used on an adaptive mesh,
478  // we would have to apply any hanging node constraint equations
479  // by uncommenting the following lines:
480  // std::vector<unsigned int> dof_indicesC = dof_indices;
481  // std::vector<unsigned int> dof_indicesM = dof_indices;
482  // dof_map.constrain_element_matrix (Ke, dof_indices);
483  // dof_map.constrain_element_matrix (Ce, dof_indicesC);
484  // dof_map.constrain_element_matrix (Me, dof_indicesM);
485 
486  // Finally, simply add the contributions to the additional
487  // matrices and vector.
488  stiffness.add_matrix (Ke, dof_indices);
489  damping.add_matrix (Ce, dof_indices);
490  mass.add_matrix (Me, dof_indices);
491 
492  // For the overall matrix, explicitly zero the entries where
493  // we added values in the other ones, so that we have
494  // identical sparsity footprints.
495  f_system.matrix->add_matrix(zero_matrix, dof_indices);
496  } // end loop over elements
497 }

References libMesh::MeshBase::active_local_element_ptr_range(), libMesh::SparseMatrix< T >::add_matrix(), libMesh::FEGenericBase< OutputType >::build(), dim, libMesh::FIFTH, libMesh::Parameters::get(), libMesh::System::get_dof_map(), libMesh::ImplicitSystem::get_matrix(), libMesh::EquationSystems::get_mesh(), libMesh::EquationSystems::get_system(), libMesh::ImplicitSystem::matrix, mesh, libMesh::MeshBase::mesh_dimension(), libMesh::QBase::n_points(), libMesh::EquationSystems::parameters, libMesh::Real, libMesh::DenseVector< T >::resize(), libMesh::DenseMatrix< T >::resize(), and libMesh::SECOND.

Referenced by main().

◆ main()

int main ( int  argc,
char **  argv 
)

Definition at line 103 of file miscellaneous_ex2.C.

104 {
105  // The libMeshInit object initializes MPI, PETSc, etc, and must be
106  // constructed before all other objects.
107  LibMeshInit init (argc, argv);
108 
109  // This example is designed for complex numbers.
110 #ifndef LIBMESH_USE_COMPLEX_NUMBERS
111  libmesh_example_requires(false, "--enable-complex");
112 #else
113 
114  // This example requires at least 2D support.
115  libmesh_example_requires(2 <= LIBMESH_DIM, "2D support");
116 
117  // Check for proper usage. frequency has two terms:
118  if ( argc <4 )
119  libmesh_error_msg("Usage: " << argv[0] << " -f [real_frequency imag_frequency]");
120 
121  if (init.comm().size() > 1)
122  {
123  if (init.comm().rank() == 0)
124  {
125  libMesh::err << "TODO: This example should be able to run in parallel."
126  << std::endl;
127  }
128  return 0;
129  }
130 
131  // Tell the user what we are doing.
132  else
133  {
134  libMesh::out << "Running " << argv[0];
135 
136  for (int i=1; i<argc; i++)
137  libMesh::out << " " << argv[i];
138 
139  libMesh::out << std::endl << std::endl;
140  }
141 
142  // Get the frequency from argv[2] as a float, solve for 1/3rd, 2/3rd
143  // and 1/1th of the given frequency.
144 
145  const Number frequency_in (atof(argv[2]), atof(argv[3]));
146 
147  // this must work as well, but is deprecated.
148  // const Real frequency_in = atof(argv[2]);
149 
150  const unsigned int n_frequencies = 3;
151 
152  // Create a mesh, with dimension to be overridden later, distributed
153  // across the default MPI communicator.
154  Mesh mesh(init.comm());
155 
156  // Read the mesh file. Here the file lshape.unv contains
157  // an L-shaped domain in .unv format.
158  UNVIO unvio(mesh);
159  unvio.read("lshape.unv");
160 
161  // Manually prepare the mesh for use, this is not done automatically
162  // by the UNVIO reader, so we have to do it here. Note that calling
163  // this function renumbers the nodes and elements of the Mesh, but
164  // the original numbering is still stored in the UNVIO object for
165  // correctly mapping dataset values (see below) to the correct nodes.
167 
168  // Read the dataset accompanying this problem. The load on the
169  // boundary of the domain is stored in the .unv data file
170  // lshape_data.unv. The data are given as complex-valued normal
171  // velocities.
172  unvio.read_dataset("lshape_data.unv");
173 
174  // Print information about the mesh to the screen.
175  mesh.print_info();
176 
177  // Create an EquationSystems object.
178  EquationSystems equation_systems (mesh);
179 
180  // Create a FrequencySystem named "Helmholtz" and store a
181  // reference to it.
182  FrequencySystem & f_system =
183  equation_systems.add_system<FrequencySystem> ("Helmholtz");
184 
185  // Add the variable "p" to the "Helmholtz" system. "p"
186  // will be approximated using second-order approximation.
187  f_system.add_variable("p", SECOND);
188 
189  // The FrequencySystem requires two user-provided functions: one for
190  // assembling the different operators and another that specifies how
191  // to combine them before the solve.
194 
195  // To enable the fast solution scheme, additional sparse matrices
196  // and one vector have to be added. The FrequencySystem object
197  // takes care of sizing the additional objects. The user should
198  // still set the sparsity structure of the f_system.matrix, so that
199  // the fast matrix addition method can be used. The procedure for
200  // this is shown in detail in the assembly function.
201  f_system.add_matrix ("stiffness");
202  f_system.add_matrix ("damping");
203  f_system.add_matrix ("mass");
204  f_system.add_vector ("rhs");
205 
206  // Communicate the frequencies to the system. Note that the
207  // frequency system stores the frequencies as parameters in the
208  // EquationSystems object, so that our assemble and solve functions
209  // may directly access them. Must be called before the
210  // EquationSystems object is initialized. Will solve for 1/3,
211  // 2/3, and 1 times the given frequency.
212  f_system.set_frequencies_by_steps (frequency_in/static_cast<Number>(n_frequencies),
213  frequency_in,
214  n_frequencies);
215 
216  // Set the wave velocity and fluid density parameter values,
217  // otherwise default values will be used.
218  equation_systems.parameters.set<Real> ("wave speed") = 1.;
219  equation_systems.parameters.set<Real> ("rho") = 1.;
220 
221  // Initialize the EquationSystems object.
222  equation_systems.init ();
223 
224  // Set values in the "rhs" vector based on the entries stored in the
225  // UNVIO object from the dataset we read in. These values only need
226  // to be set once, as they are the same for every frequency. We can
227  // only do this once equation_systems.init() has been called...
228  {
229  NumericVector<Number> & freq_indep_rhs = f_system.get_vector("rhs");
230 
231  for (const auto & node : mesh.node_ptr_range())
232  {
233  // Get the data read in from the dataset for the current Node, if any.
234  const std::vector<Number> * nodal_data = unvio.get_data(node);
235 
236  // Set the rhs value based on values read in from the dataset.
237  if (nodal_data)
238  {
239  unsigned int dn = node->dof_number(/*system=*/0,
240  /*variable=*/0,
241  /*component=*/0);
242  freq_indep_rhs.set(dn, (*nodal_data)[0]);
243  }
244  }
245  }
246 
247  // Print information about the system to the screen.
248  equation_systems.print_info ();
249 
250  for (unsigned int n=0; n < n_frequencies; n++)
251  {
252  // Solve the "Helmholtz" system for the n-th frequency.
253  // Since we attached an assemble() function to the system,
254  // the mass, damping, and stiffness contributions will only
255  // be assembled once. Then, the system is solved for the
256  // given frequencies. Note that solve() may also solve
257  // the system only for specific frequencies.
258  f_system.solve (n, n);
259 
260  // After solving the system, write the solution to an ExodusII
261  // file for every frequency.
262 #ifdef LIBMESH_HAVE_EXODUS_API
263  std::ostringstream file_name;
264 
265  file_name << "out_"
266  << std::setw(4)
267  << std::setfill('0')
268  << std::right
269  << n
270  << ".e";
271 
272  ExodusII_IO(mesh).write_equation_systems (file_name.str(),
273  equation_systems);
274 #endif
275  }
276 
277  // Alternatively, the whole EquationSystems object can be
278  // written to disk. By default, the additional vectors are also
279  // saved.
280  equation_systems.write ("eqn_sys.dat", WRITE);
281 
282  // All done.
283  return 0;
284 
285 #endif
286 }

References add_M_C_K_helmholtz(), libMesh::ImplicitSystem::add_matrix(), libMesh::EquationSystems::add_system(), libMesh::System::add_variable(), libMesh::System::add_vector(), assemble_helmholtz(), libMesh::System::attach_assemble_function(), libMesh::FrequencySystem::attach_solve_function(), libMesh::err, libMesh::System::get_vector(), libMesh::TriangleWrapper::init(), libMesh::EquationSystems::init(), mesh, libMesh::MeshBase::node_ptr_range(), libMesh::out, libMesh::EquationSystems::parameters, libMesh::MeshBase::prepare_for_use(), libMesh::EquationSystems::print_info(), libMesh::MeshBase::print_info(), libMesh::Real, libMesh::SECOND, libMesh::Parameters::set(), libMesh::NumericVector< T >::set(), libMesh::FrequencySystem::set_frequencies_by_steps(), libMesh::FrequencySystem::solve(), libMesh::WRITE, libMesh::EquationSystems::write(), and libMesh::MeshOutput< MT >::write_equation_systems().

libMesh::SparseMatrix::close
virtual void close()=0
Calls the SparseMatrix's internal assembly routines, ensuring that the values are consistent across p...
libMesh::NumericVector::zero
virtual void zero()=0
Set all entries to zero.
libMesh::Number
Real Number
Definition: libmesh_common.h:195
libMesh::pi
const Real pi
.
Definition: libmesh.h:237
libMesh::NumericVector::add
virtual void add(const numeric_index_type i, const T value)=0
Adds value to each entry of the vector.
libMesh::Mesh
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
Definition: mesh.h:50
libMesh::EquationSystems::get_mesh
const MeshBase & get_mesh() const
Definition: equation_systems.h:637
libMesh::ExplicitSystem::rhs
NumericVector< Number > * rhs
The system matrix.
Definition: explicit_system.h:114
libMesh::QGauss
This class implements specific orders of Gauss quadrature.
Definition: quadrature_gauss.h:39
libMesh::FrequencySystem::attach_solve_function
void attach_solve_function(void fptr(EquationSystems &es, const std::string &name))
Register a required user function to use in assembling/solving the system.
Definition: frequency_system.C:415
libMesh::MeshBase::active_local_element_ptr_range
virtual SimpleRange< element_iterator > active_local_element_ptr_range()=0
libMesh::NumericVector::close
virtual void close()=0
Calls the NumericVector's internal assembly routines, ensuring that the values are consistent across ...
libMesh::ImplicitSystem::add_matrix
SparseMatrix< Number > & add_matrix(const std::string &mat_name)
Adds the additional matrix mat_name to this system.
Definition: implicit_system.C:202
libMesh::EquationSystems::get_system
const T_sys & get_system(const std::string &name) const
Definition: equation_systems.h:757
libMesh::DenseMatrix< Number >
libMesh::WRITE
Definition: enum_xdr_mode.h:40
mesh
MeshBase & mesh
Definition: mesh_communication.C:1257
libMesh::MeshBase::mesh_dimension
unsigned int mesh_dimension() const
Definition: mesh_base.C:135
libMesh::FIFTH
Definition: enum_order.h:46
libMesh::ExodusII_IO
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
Definition: exodusII_io.h:51
libMesh::SECOND
Definition: enum_order.h:43
libMesh::SparseMatrix< Number >
dim
unsigned int dim
Definition: adaptivity_ex3.C:113
libMesh::DenseMatrix::resize
void resize(const unsigned int new_m, const unsigned int new_n)
Resize the matrix.
Definition: dense_matrix.h:822
libMesh::NumericVector< Number >
libMesh::TriangleWrapper::init
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
libMesh::System::attach_assemble_function
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:1755
libMesh::MeshBase
This is the MeshBase class.
Definition: mesh_base.h:78
libMesh::UNVIO
The UNVIO class implements the Ideas UNV universal file format.
Definition: unv_io.h:52
libMesh::MeshBase::node_ptr_range
virtual SimpleRange< node_iterator > node_ptr_range()=0
libMesh::System::add_variable
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.
Definition: system.C:1069
libMesh::ImplicitSystem::matrix
SparseMatrix< Number > * matrix
The system matrix.
Definition: implicit_system.h:393
libMesh::FrequencySystem
FrequencySystem provides a specific system class for frequency-dependent (linear) systems.
Definition: frequency_system.h:63
libMesh::LibMeshInit
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
Definition: libmesh.h:83
libMesh::SparseMatrix::add
virtual void add(const numeric_index_type i, const numeric_index_type j, const T value)=0
Add value to the element (i,j).
libMesh::EquationSystems
This is the EquationSystems class.
Definition: equation_systems.h:74
libMesh::MeshOutput::write_equation_systems
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
libMesh::DenseVector::resize
void resize(const unsigned int n)
Resize the vector.
Definition: dense_vector.h:355
libMesh::SparseMatrix::add_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.
libMesh::FEType
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:178
libMesh::FrequencySystem::set_frequencies_by_steps
void set_frequencies_by_steps(const Number base_freq, const Number freq_step=0., const unsigned int n_freq=1, const bool allocate_solution_duplicates=true)
Set the frequency range for which the system should be solved.
Definition: frequency_system.C:180
libMesh::NumericVector::set
virtual void set(const numeric_index_type i, const T value)=0
Sets v(i) = value.
libMesh::DofMap
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:176
libMesh::MeshBase::print_info
void print_info(std::ostream &os=libMesh::out) const
Prints relevant information about the mesh.
Definition: mesh_base.C:585
assemble_helmholtz
void assemble_helmholtz(EquationSystems &es, const std::string &system_name)
Definition: miscellaneous_ex2.C:293
libMesh::FrequencySystem::solve
virtual void solve() override
Solves the system for all frequencies.
Definition: frequency_system.C:338
add_M_C_K_helmholtz
void add_M_C_K_helmholtz(EquationSystems &es, const std::string &system_name)
Definition: miscellaneous_ex2.C:502
libMesh::System::add_vector
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.
Definition: system.C:661
libMesh::System::get_dof_map
const DofMap & get_dof_map() const
Definition: system.h:2099
libMesh::err
OStreamProxy err
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::MeshBase::prepare_for_use
void prepare_for_use(const bool skip_renumber_nodes_and_elements=false, const bool skip_find_neighbors=false)
Prepare a newly ecreated (or read) mesh for use.
Definition: mesh_base.C:318
libMesh::Parameters::get
const T & get(const std::string &) const
Definition: parameters.h:421
libMesh::out
OStreamProxy out
libMesh::ImplicitSystem::get_matrix
const SparseMatrix< Number > & get_matrix(const std::string &mat_name) const
Definition: implicit_system.C:262
libMesh::System::get_vector
const NumericVector< Number > & get_vector(const std::string &vec_name) const
Definition: system.C:774
libMesh::SparseMatrix::zero
virtual void zero()=0
Set all entries to 0.
libMesh::DenseVector< Number >
libMesh::EquationSystems::parameters
Parameters parameters
Data structure holding arbitrary parameters.
Definition: equation_systems.h:557