libMesh
Functions
transient_ex2.C File Reference

Go to the source code of this file.

Functions

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

Function Documentation

◆ apply_initial()

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

Definition at line 538 of file transient_ex2.C.

References libMesh::EquationSystems::get_system(), libMesh::System::get_vector(), and libMesh::NumericVector< T >::zero().

Referenced by main().

540 {
541  // Get a reference to our system, as before
542  NewmarkSystem & t_system = es.get_system<NewmarkSystem> (system_name);
543 
544  // Numeric vectors for the pressure, velocity and acceleration
545  // values.
546  NumericVector<Number> & pres_vec = t_system.get_vector("displacement");
547  NumericVector<Number> & vel_vec = t_system.get_vector("velocity");
548  NumericVector<Number> & acc_vec = t_system.get_vector("acceleration");
549 
550  // Assume our fluid to be at rest, which would
551  // also be the default conditions in class NewmarkSystem,
552  // but let us do it explicitly here.
553  pres_vec.zero();
554  vel_vec.zero();
555  acc_vec.zero();
556 }
const T_sys & get_system(std::string_view name) const
virtual void zero()=0
Set all entries to zero.
This class contains a specific system class.
const NumericVector< Number > & get_vector(std::string_view vec_name) const
Definition: system.C:918

◆ assemble_wave()

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

Definition at line 318 of file transient_ex2.C.

References libMesh::SparseMatrix< T >::add_matrix(), libMesh::NumericVector< T >::add_vector(), libMesh::FEGenericBase< OutputType >::build(), dim, libMesh::DofMap::dof_indices(), libMesh::Parameters::get(), libMesh::EquationSystems::get_mesh(), libMesh::EquationSystems::get_system(), 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().

320 {
321  // It is a good idea to make sure we are assembling
322  // the proper system.
323  libmesh_assert_equal_to (system_name, "Wave");
324 
325  // Get a constant reference to the mesh object.
326  const MeshBase & mesh = es.get_mesh();
327 
328  // The dimension that we are running.
329  const unsigned int dim = mesh.mesh_dimension();
330 
331  // Copy the speed of sound to a local variable.
332  const Real speed = es.parameters.get<Real>("speed");
333 
334  // If we added Neumann conditions we would need density too
335  // const Real rho = es.parameters.get<Real>("fluid density");
336 
337  // Get a reference to our system, as before.
338  NewmarkSystem & t_system = es.get_system<NewmarkSystem> (system_name);
339 
340  // Get a constant reference to the Finite Element type
341  // for the first (and only) variable in the system.
342  FEType fe_type = t_system.get_dof_map().variable_type(0);
343 
344  // In here, we will add the element matrices to the
345  // @e additional matrices "stiffness_mass" and "damping"
346  // and the additional vector "force", not to the members
347  // "matrix" and "rhs". Therefore, get writable
348  // references to them.
349  SparseMatrix<Number> & stiffness = t_system.get_matrix("stiffness");
350  SparseMatrix<Number> & damping = t_system.get_matrix("damping");
351  SparseMatrix<Number> & mass = t_system.get_matrix("mass");
352  NumericVector<Number> & force = t_system.get_vector("force");
353 
354  // Some solver packages (PETSc) are especially picky about
355  // allocating sparsity structure and truly assigning values
356  // to this structure. Namely, matrix additions, as performed
357  // later, exhibit acceptable performance only for identical
358  // sparsity structures. Therefore, explicitly zero the
359  // values in the collective matrix, so that matrix additions
360  // encounter identical sparsity structures.
361  SparseMatrix<Number> & matrix = *t_system.matrix;
362  DenseMatrix<Number> zero_matrix;
363 
364  // Build a Finite Element object of the specified type. Since the
365  // FEBase::build() member dynamically creates memory we will
366  // store the object as a std::unique_ptr<FEBase>. This can be thought
367  // of as a pointer that will clean up after itself.
368  std::unique_ptr<FEBase> fe (FEBase::build(dim, fe_type));
369 
370  // A 2nd order Gauss quadrature rule for numerical integration.
371  QGauss qrule (dim, SECOND);
372 
373  // Tell the finite element object to use our quadrature rule.
374  fe->attach_quadrature_rule (&qrule);
375 
376  // The element Jacobian * quadrature weight at each integration point.
377  const std::vector<Real> & JxW = fe->get_JxW();
378 
379  // The element shape functions evaluated at the quadrature points.
380  const std::vector<std::vector<Real>> & phi = fe->get_phi();
381 
382  // The element shape function gradients evaluated at the quadrature
383  // points.
384  const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
385 
386  // A reference to the DofMap object for this system. The DofMap
387  // object handles the index translation from node and element numbers
388  // to degree of freedom numbers.
389  const DofMap & dof_map = t_system.get_dof_map();
390 
391  // The element mass, damping and stiffness matrices
392  // and the element contribution to the rhs.
393  DenseMatrix<Number> Ke, Ce, Me;
395 
396  // This vector will hold the degree of freedom indices for
397  // the element. These define where in the global system
398  // the element degrees of freedom get mapped.
399  std::vector<dof_id_type> dof_indices;
400 
401  // Now we will loop over all the elements in the mesh.
402  // We will compute the element matrix and right-hand-side
403  // contribution.
404  for (const auto & elem : mesh.active_local_element_ptr_range())
405  {
406  // Get the degree of freedom indices for the
407  // current element. These define where in the global
408  // matrix and right-hand-side this element will
409  // contribute to.
410  dof_map.dof_indices (elem, dof_indices);
411 
412  // Compute the element-specific data for the current
413  // element. This involves computing the location of the
414  // quadrature points (q_point) and the shape functions
415  // (phi, dphi) for the current element.
416  fe->reinit (elem);
417 
418  // Zero the element matrices and rhs before
419  // summing them. We use the resize member here because
420  // the number of degrees of freedom might have changed from
421  // the last element. Note that this will be the case if the
422  // element type is different (i.e. the last element was HEX8
423  // and now have a PRISM6).
424  {
425  const unsigned int n_dof_indices = dof_indices.size();
426 
427  Ke.resize (n_dof_indices, n_dof_indices);
428  Ce.resize (n_dof_indices, n_dof_indices);
429  Me.resize (n_dof_indices, n_dof_indices);
430  zero_matrix.resize (n_dof_indices, n_dof_indices);
431  Fe.resize (n_dof_indices);
432  }
433 
434  // Now loop over the quadrature points. This handles
435  // the numeric integration.
436  for (unsigned int qp=0; qp<qrule.n_points(); qp++)
437  {
438  // Now we will build the element matrix. This involves
439  // a double loop to integrate the test functions (i) against
440  // the trial functions (j).
441  for (std::size_t i=0; i<phi.size(); i++)
442  for (std::size_t j=0; j<phi.size(); j++)
443  {
444  Ke(i,j) += JxW[qp]*(dphi[i][qp]*dphi[j][qp]);
445  Me(i,j) += JxW[qp]*phi[i][qp]*phi[j][qp]
446  *1./(speed*speed);
447  } // end of the matrix summation loop
448  } // end of quadrature point loop
449 
450  // Now compute the contribution to the element matrix and the
451  // right-hand-side vector if the current element lies on the
452  // boundary.
453  {
454  // In this example no natural boundary conditions will
455  // be considered. The code is left here so it can easily
456  // be extended.
457  //
458  // don't do this for any side
459 #if 0
460  for (auto side : elem->side_index_range())
461  if (elem->neighbor_ptr(side) == nullptr)
462  {
463  // Declare a special finite element object for
464  // boundary integration.
465  std::unique_ptr<FEBase> fe_face (FEBase::build(dim, fe_type));
466 
467  // Boundary integration requires one quadrature rule,
468  // with dimensionality one less than the dimensionality
469  // of the element.
470  QGauss qface(dim-1, SECOND);
471 
472  // Tell the finite element object to use our
473  // quadrature rule.
474  fe_face->attach_quadrature_rule (&qface);
475 
476  // The value of the shape functions at the quadrature
477  // points.
478  const std::vector<std::vector<Real>> & phi_face = fe_face->get_phi();
479 
480  // The Jacobian * Quadrature Weight at the quadrature
481  // points on the face.
482  const std::vector<Real> & JxW_face = fe_face->get_JxW();
483 
484  // Compute the shape function values on the element
485  // face.
486  fe_face->reinit(elem, side);
487 
488  // Here we consider a normal acceleration acc_n=1 applied to
489  // the whole boundary of our mesh.
490  const Real acc_n_value = 1.0;
491 
492  // Loop over the face quadrature points for integration.
493  for (unsigned int qp=0; qp<qface.n_points(); qp++)
494  {
495  // Right-hand-side contribution due to prescribed
496  // normal acceleration.
497  for (std::size_t i=0; i<phi_face.size(); i++)
498  {
499  Fe(i) += acc_n_value*rho
500  *phi_face[i][qp]*JxW_face[qp];
501  }
502  } // end face quadrature point loop
503  } // end if (elem->neighbor_ptr(side) == nullptr)
504 #endif // 0
505 
506  // In this example the Dirichlet boundary conditions will be
507  // imposed via penalty method after the
508  // system is assembled.
509 
510  } // end boundary condition section
511 
512  // If this assembly program were to be used on an adaptive mesh,
513  // we would have to apply any hanging node constraint equations
514  // by uncommenting the following lines:
515  // std::vector<unsigned int> dof_indicesC = dof_indices;
516  // std::vector<unsigned int> dof_indicesM = dof_indices;
517  // dof_map.constrain_element_matrix_and_vector (Ke, Fe, dof_indices);
518  // dof_map.constrain_element_matrix (Ce, dof_indicesC);
519  // dof_map.constrain_element_matrix (Me, dof_indicesM);
520 
521  // Finally, simply add the contributions to the additional
522  // matrices and vector.
523  stiffness.add_matrix (Ke, dof_indices);
524  damping.add_matrix (Ce, dof_indices);
525  mass.add_matrix (Me, dof_indices);
526 
527  force.add_vector (Fe, dof_indices);
528 
529  // For the overall matrix, explicitly zero the entries where
530  // we added values in the other ones, so that we have
531  // identical sparsity footprints.
532  matrix.add_matrix(zero_matrix, dof_indices);
533 
534  } // end of element loop
535 }
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:182
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
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
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
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.
const T & get(std::string_view) const
Definition: parameters.h:426
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 contains a specific system class.
This class implements specific orders of Gauss quadrature.
unsigned int mesh_dimension() const
Definition: mesh_base.C:324
Parameters parameters
Data structure holding arbitrary parameters.

◆ fill_dirichlet_bc()

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

Definition at line 559 of file transient_ex2.C.

References std::abs(), libMesh::NumericVector< T >::add(), libMesh::DofObject::dof_number(), libMesh::Parameters::get(), libMesh::EquationSystems::get_mesh(), libMesh::EquationSystems::get_system(), libMesh::ImplicitSystem::matrix, mesh, n_nodes, libMesh::MeshBase::n_nodes(), libMesh::MeshBase::node_ref(), libMesh::EquationSystems::parameters, libMesh::pi, libMesh::Real, libMesh::ExplicitSystem::rhs, libMesh::System::time, and libMesh::TOLERANCE.

Referenced by main().

561 {
562  // It is a good idea to make sure we are assembling
563  // the proper system.
564  libmesh_assert_equal_to (system_name, "Wave");
565 
566  // Get a reference to our system, as before.
567  NewmarkSystem & t_system = es.get_system<NewmarkSystem> (system_name);
568 
569  // Get writable references to the overall matrix and vector.
570  SparseMatrix<Number> & matrix = *t_system.matrix;
571  NumericVector<Number> & rhs = *t_system.rhs;
572 
573  // Get a constant reference to the mesh object.
574  const MeshBase & mesh = es.get_mesh();
575 
576  // Get libMesh's pi
577  const Real pi = libMesh::pi;
578 
579  // Ask the EquationSystems flag whether
580  // we should do this also for the matrix
581  const bool do_for_matrix =
582  es.parameters.get<bool>("Newmark set BC for Matrix");
583 
584  // Number of nodes in the mesh.
585  unsigned int n_nodes = mesh.n_nodes();
586 
587  for (unsigned int n_cnt=0; n_cnt<n_nodes; n_cnt++)
588  {
589  // Get a reference to the current node.
590  const Node & curr_node = mesh.node_ref(n_cnt);
591 
592  // Check if Dirichlet BCs should be applied to this node.
593  // Use the TOLERANCE from mesh_common.h as tolerance.
594  // Here a pressure value is applied if the z-coord.
595  // is equal to 4, which corresponds to one end of the
596  // pipe-mesh in this directory.
597  const Real z_coo = 4.;
598 
599  if (std::abs(curr_node(2)-z_coo) < TOLERANCE)
600  {
601  // The global number of the respective degree of freedom.
602  unsigned int dn = curr_node.dof_number(0, 0, 0);
603 
604  // The penalty parameter.
605  const Real penalty = 1.e10;
606 
607  // Here we apply sinusoidal pressure values for 0<t<0.002
608  // at one end of the pipe-mesh.
609  Real p_value;
610  if (t_system.time < .002)
611  p_value = sin(2*pi*t_system.time/.002);
612  else
613  p_value = .0;
614 
615  // Now add the contributions to the matrix and the rhs.
616  rhs.add(dn, p_value*penalty);
617 
618  // Add the penalty parameter to the global matrix
619  // if desired.
620  if (do_for_matrix)
621  matrix.add(dn, dn, penalty);
622  }
623  } // loop n_cnt
624 }
Real time
For time-dependent problems, this is the time t at the beginning of the current timestep.
Definition: system.h:1595
dof_id_type dof_number(const unsigned int s, const unsigned int var, const unsigned int comp) const
Definition: dof_object.h:1025
A Node is like a Point, but with more information.
Definition: node.h:52
static constexpr Real TOLERANCE
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
ADRealEigenVector< T, D, asd > abs(const ADRealEigenVector< T, D, asd > &)
Definition: type_vector.h:57
const dof_id_type n_nodes
Definition: tecplot_io.C:67
const T & get(std::string_view) const
Definition: parameters.h:426
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
SparseMatrix< Number > * matrix
The system matrix.
const MeshBase & get_mesh() const
This class contains a specific system class.
Parameters parameters
Data structure holding arbitrary parameters.
virtual const Node & node_ref(const dof_id_type i) const
Definition: mesh_base.h:575
virtual void add(const numeric_index_type i, const T value)=0
Adds value to the vector entry specified by i.
virtual dof_id_type n_nodes() const =0
const Real pi
.
Definition: libmesh.h:274

◆ main()

int main ( int  argc,
char **  argv 
)

Definition at line 100 of file transient_ex2.C.

References libMesh::EquationSystems::add_system(), libMesh::System::add_variable(), apply_initial(), assemble_wave(), libMesh::default_solver_package(), libMesh::DofObject::dof_number(), libMesh::EIGEN_SOLVERS, fill_dirichlet_bc(), libMesh::FIRST, libMesh::EquationSystems::get_system(), libMesh::TriangleWrapper::init(), libMesh::EquationSystems::init(), libMesh::INVALID_SOLVER_PACKAGE, libMesh::NumericVector< T >::localize(), mesh, libMesh::MeshBase::node_ref(), libMesh::out, libMesh::EquationSystems::parameters, libMesh::PETSC_SOLVERS, libMesh::EquationSystems::print_info(), libMesh::MeshBase::print_info(), libMesh::MeshBase::read(), libMesh::Real, libMesh::Parameters::set(), libMesh::NumericVector< T >::size(), and libMesh::MeshOutput< MT >::write_equation_systems().

101 {
102  // Initialize libraries, like in example 2.
103  LibMeshInit init (argc, argv);
104 
105  // This example requires a linear solver package.
106  libmesh_example_requires(libMesh::default_solver_package() != INVALID_SOLVER_PACKAGE,
107  "--enable-petsc, --enable-trilinos, or --enable-eigen");
108 
109  // Check for proper usage.
110  libmesh_error_msg_if(argc < 2, "Usage: " << argv[0] << " [meshfile]");
111 
112  // Tell the user what we are doing.
113  libMesh::out << "Running " << argv[0];
114 
115  for (int i=1; i<argc; i++)
116  libMesh::out << " " << argv[i];
117 
118  libMesh::out << std::endl << std::endl;
119 
120  // LasPack solvers don't work so well for this example, Trilinos doesn't work at all.
121  // PETSc and Eigen both work...
122  libmesh_example_requires(libMesh::default_solver_package() == PETSC_SOLVERS || \
123  libMesh::default_solver_package() == EIGEN_SOLVERS, "--enable-petsc");
124 
125  // Get the name of the mesh file
126  // from the command line.
127  std::string mesh_file = argv[1];
128  libMesh::out << "Mesh file is: " << mesh_file << std::endl;
129 
130  // Skip this 3D example if libMesh was compiled as 1D or 2D-only.
131  libmesh_example_requires(3 <= LIBMESH_DIM, "3D support");
132 
133  // Create a mesh.
134  // This example directly references all mesh nodes and is
135  // incompatible with DistributedMesh use.
136  //
137  // Create a ReplicatedMesh object, with dimension to be overridden
138  // later, distributed across the default MPI communicator.
139  ReplicatedMesh mesh(init.comm());
140 
141  // Read the meshfile specified on the command line.
142  mesh.read(mesh_file);
143 
144  // Print information about the mesh to the screen.
145  mesh.print_info();
146 
147  // The node that should be monitored.
148  const unsigned int result_node = 274;
149 
150 
151  // Time stepping issues
152  //
153  // Note that the total current time is stored as a parameter
154  // in the \pEquationSystems object.
155  //
156  // the time step size
157  const Real delta_t = .0000625;
158 
159  // The number of time steps.
160  unsigned int n_time_steps = 300;
161 
162  // Create an equation systems object.
163  EquationSystems equation_systems (mesh);
164 
165  // Declare the system and its variables.
166  // Create a NewmarkSystem named "Wave"
167  equation_systems.add_system<NewmarkSystem> ("Wave");
168 
169  // Use a handy reference to this system
170  NewmarkSystem & t_system = equation_systems.get_system<NewmarkSystem> ("Wave");
171 
172  // Add the variable "p" to "Wave". "p"
173  // will be approximated using first-order approximation.
174  t_system.add_variable("p", FIRST);
175 
176  // Give the system a pointer to the matrix assembly
177  // function and the initial condition function defined
178  // below.
179  t_system.attach_assemble_function (assemble_wave);
180  t_system.attach_init_function (apply_initial);
181 
182  // Set the time step size, and optionally the
183  // Newmark parameters, so that NewmarkSystem can
184  // compute integration constants. Here we simply use
185  // pass only the time step and use default values
186  // for alpha=.25 and delta=.5.
187  t_system.set_newmark_parameters(delta_t);
188 
189  // Set the speed of sound and fluid density
190  // as EquationSystems parameter,
191  // so that assemble_wave() can access it.
192  equation_systems.parameters.set<Real>("speed") = 1000.;
193  equation_systems.parameters.set<Real>("fluid density") = 1000.;
194 
195  // Start time integration from t=0
196  t_system.time = 0.;
197 
198  // Initialize the data structures for the equation system.
199  equation_systems.init();
200 
201  // Prints information about the system to the screen.
202  equation_systems.print_info();
203 
204  // A file to store the results at certain nodes.
205  std::ofstream res_out("pressure_node.res");
206 
207  // get the dof_numbers for the nodes that
208  // should be monitored.
209  const unsigned int res_node_no = result_node;
210  const Node & res_node = mesh.node_ref(res_node_no-1);
211  unsigned int dof_no = res_node.dof_number(0, 0, 0);
212 
213  // Assemble the time independent system matrices and rhs.
214  // This function will also compute the effective system matrix
215  // K~=K+a_0*M+a_1*C and apply user specified initial
216  // conditions.
217  t_system.assemble();
218 
219  // Now solve for each time step.
220  // For convenience, use a local buffer of the
221  // current time. But once this time is updated,
222  // also update the EquationSystems parameter
223  // Start with t_time = 0 and write a short header
224  // to the nodal result file
225  res_out << "# pressure at node " << res_node_no << "\n"
226  << "# time\tpressure\n"
227  << t_system.time << "\t" << 0 << std::endl;
228 
229 
230  for (unsigned int time_step=0; time_step<n_time_steps; time_step++)
231  {
232  // Update the time. Both here and in the
233  // EquationSystems object
234  t_system.time += delta_t;
235 
236  // Update the rhs.
237  t_system.update_rhs();
238 
239  // Impose essential boundary conditions.
240  // Not that since the matrix is only assembled once,
241  // the penalty parameter should be added to the matrix
242  // only in the first time step. The applied
243  // boundary conditions may be time-dependent and hence
244  // the rhs vector is considered in each time step.
245  if (time_step == 0)
246  {
247  // The local function fill_dirichlet_bc()
248  // may also set Dirichlet boundary conditions for the
249  // matrix. When you set the flag as shown below,
250  // the flag will return true. If you want it to return
251  // false, simply do not set it.
252  equation_systems.parameters.set<bool>("Newmark set BC for Matrix") = true;
253 
254  fill_dirichlet_bc(equation_systems, "Wave");
255 
256  // unset the flag, so that it returns false
257  equation_systems.parameters.set<bool>("Newmark set BC for Matrix") = false;
258  }
259  else
260  fill_dirichlet_bc(equation_systems, "Wave");
261 
262  // Solve the system "Wave".
263  t_system.solve();
264 
265  // After solving the system, write the solution
266  // to a GMV-formatted plot file.
267  // Do only for a few time steps.
268  if (time_step == 30 || time_step == 60 ||
269  time_step == 90 || time_step == 120)
270  {
271  std::ostringstream file_name;
272 
273 #ifdef LIBMESH_HAVE_VTK
274  file_name << "out_"
275  << std::setw(3)
276  << std::setfill('0')
277  << std::right
278  << time_step
279  << ".pvtu";
280 
281  VTKIO(mesh).write_equation_systems (file_name.str(), equation_systems);
282 #else
283 
284  file_name << "out."
285  << std::setw(3)
286  << std::setfill('0')
287  << std::right
288  << time_step
289  << ".gmv";
290 
291  GMVIO(mesh).write_equation_systems (file_name.str(), equation_systems);
292 #endif
293  }
294 
295  // Update the p, v and a.
296  t_system.update_u_v_a();
297 
298  // dof_no may not be local in parallel runs, so we may need a
299  // global displacement vector
300  NumericVector<Number> & displacement = t_system.get_vector("displacement");
301  std::vector<Number> global_displacement(displacement.size());
302  displacement.localize(global_displacement);
303 
304  // Write nodal results to file. The results can then
305  // be viewed with e.g. gnuplot (run gnuplot and type
306  // 'plot "pressure_node.res" with lines' in the command line)
307  res_out << t_system.time << "\t"
308  << global_displacement[dof_no]
309  << std::endl;
310  }
311 
312  // All done.
313  return 0;
314 }
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.
dof_id_type dof_number(const unsigned int s, const unsigned int var, const unsigned int comp) const
Definition: dof_object.h:1025
The ReplicatedMesh class is derived from the MeshBase class, and is used to store identical copies of...
A Node is like a Point, but with more information.
Definition: node.h:52
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
virtual numeric_index_type size() const =0
MeshBase & mesh
void apply_initial(EquationSystems &es, const std::string &system_name)
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
This class implements reading and writing meshes in the VTK format.
Definition: vtk_io.h:60
SolverPackage default_solver_package()
Definition: libmesh.C:1050
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:1489
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 fill_dirichlet_bc(EquationSystems &es, const std::string &system_name)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
OStreamProxy out
This class contains a specific system class.
virtual const Node & node_ref(const dof_id_type i) const
Definition: mesh_base.h:575
void assemble_wave(EquationSystems &es, const std::string &system_name)
virtual void localize(std::vector< T > &v_local) const =0
Creates a copy of the global vector in the local vector v_local.