libMesh
transient_ex2.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2019 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 // <h1>Transient Example 2 - The Newmark System and the Wave Equation</h1>
21 // \author Steffen Petersen
22 // \date 2003
23 //
24 // This is the eighth example program. It builds on
25 // the previous example programs. It introduces the
26 // NewmarkSystem class. In this example the wave equation
27 // is solved using the time integration scheme provided
28 // by the NewmarkSystem class.
29 //
30 // This example comes with a cylindrical mesh given in the
31 // universal file pipe-mesh.unv.
32 // The mesh contains HEX8 and PRISM6 elements.
33 
34 // C++ include files that we need
35 #include <iostream>
36 #include <fstream>
37 #include <sstream>
38 #include <algorithm>
39 #include <stdio.h>
40 #include <math.h>
41 
42 // Basic include file needed for the mesh functionality.
43 #include "libmesh/libmesh.h"
44 #include "libmesh/replicated_mesh.h"
45 #include "libmesh/gmv_io.h"
46 #include "libmesh/vtk_io.h"
47 #include "libmesh/newmark_system.h"
48 #include "libmesh/equation_systems.h"
49 #include "libmesh/enum_solver_package.h"
50 
51 // Define the Finite Element object.
52 #include "libmesh/fe.h"
53 
54 // Define Gauss quadrature rules.
55 #include "libmesh/quadrature_gauss.h"
56 
57 // Define useful datatypes for finite element
58 #include "libmesh/dense_matrix.h"
59 #include "libmesh/dense_vector.h"
60 
61 // Define matrix and vector data types for the global
62 // equation system. These are base classes,
63 // from which specific implementations, like
64 // the PETSc and Eigen implementations, are derived.
65 #include "libmesh/sparse_matrix.h"
66 #include "libmesh/numeric_vector.h"
67 
68 // Define the DofMap, which handles degree of freedom
69 // indexing.
70 #include "libmesh/dof_map.h"
71 
72 // The definition of a vertex associated with a Mesh.
73 #include "libmesh/node.h"
74 
75 // The definition of a geometric element
76 #include "libmesh/elem.h"
77 
78 // Bring in everything from the libMesh namespace
79 using namespace libMesh;
80 
81 // Function prototype. This is the function that will assemble
82 // the linear system for our problem, governed by the linear
83 // wave equation.
85  const std::string & system_name);
86 
87 
88 // Function Prototype. This function will be used to apply the
89 // initial conditions.
91  const std::string & system_name);
92 
93 // Function Prototype. This function imposes
94 // Dirichlet Boundary conditions via the penalty
95 // method after the system is assembled.
97  const std::string & system_name);
98 
99 // The main program
100 int main (int argc, char** argv)
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  if (argc < 2)
111  libmesh_error_msg("Usage: " << argv[0] << " [meshfile]");
112 
113  // Tell the user what we are doing.
114  else
115  {
116  libMesh::out << "Running " << argv[0];
117 
118  for (int i=1; i<argc; i++)
119  libMesh::out << " " << argv[i];
120 
121  libMesh::out << std::endl << std::endl;
122 
123  }
124 
125  // LasPack solvers don't work so well for this example, Trilinos doesn't work at all.
126  // PETSc and Eigen both work...
127  libmesh_example_requires(libMesh::default_solver_package() == PETSC_SOLVERS || \
128  libMesh::default_solver_package() == EIGEN_SOLVERS, "--enable-petsc");
129 
130  // Get the name of the mesh file
131  // from the command line.
132  std::string mesh_file = argv[1];
133  libMesh::out << "Mesh file is: " << mesh_file << std::endl;
134 
135  // Skip this 3D example if libMesh was compiled as 1D or 2D-only.
136  libmesh_example_requires(3 <= LIBMESH_DIM, "3D support");
137 
138  // Create a mesh.
139  // This example directly references all mesh nodes and is
140  // incompatible with DistributedMesh use.
141  //
142  // Create a ReplicatedMesh object, with dimension to be overridden
143  // later, distributed across the default MPI communicator.
144  ReplicatedMesh mesh(init.comm());
145 
146  // Read the meshfile specified on the command line.
147  mesh.read(mesh_file);
148 
149  // Print information about the mesh to the screen.
150  mesh.print_info();
151 
152  // The node that should be monitored.
153  const unsigned int result_node = 274;
154 
155 
156  // Time stepping issues
157  //
158  // Note that the total current time is stored as a parameter
159  // in the \pEquationSystems object.
160  //
161  // the time step size
162  const Real delta_t = .0000625;
163 
164  // The number of time steps.
165  unsigned int n_time_steps = 300;
166 
167  // Create an equation systems object.
168  EquationSystems equation_systems (mesh);
169 
170  // Declare the system and its variables.
171  // Create a NewmarkSystem named "Wave"
172  equation_systems.add_system<NewmarkSystem> ("Wave");
173 
174  // Use a handy reference to this system
175  NewmarkSystem & t_system = equation_systems.get_system<NewmarkSystem> ("Wave");
176 
177  // Add the variable "p" to "Wave". "p"
178  // will be approximated using first-order approximation.
179  t_system.add_variable("p", FIRST);
180 
181  // Give the system a pointer to the matrix assembly
182  // function and the initial condition function defined
183  // below.
184  t_system.attach_assemble_function (assemble_wave);
185  t_system.attach_init_function (apply_initial);
186 
187  // Set the time step size, and optionally the
188  // Newmark parameters, so that NewmarkSystem can
189  // compute integration constants. Here we simply use
190  // pass only the time step and use default values
191  // for alpha=.25 and delta=.5.
192  t_system.set_newmark_parameters(delta_t);
193 
194  // Set the speed of sound and fluid density
195  // as EquationSystems parameter,
196  // so that assemble_wave() can access it.
197  equation_systems.parameters.set<Real>("speed") = 1000.;
198  equation_systems.parameters.set<Real>("fluid density") = 1000.;
199 
200  // Start time integration from t=0
201  t_system.time = 0.;
202 
203  // Initialize the data structures for the equation system.
204  equation_systems.init();
205 
206  // Prints information about the system to the screen.
207  equation_systems.print_info();
208 
209  // A file to store the results at certain nodes.
210  std::ofstream res_out("pressure_node.res");
211 
212  // get the dof_numbers for the nodes that
213  // should be monitored.
214  const unsigned int res_node_no = result_node;
215  const Node & res_node = mesh.node_ref(res_node_no-1);
216  unsigned int dof_no = res_node.dof_number(0, 0, 0);
217 
218  // Assemble the time independent system matrices and rhs.
219  // This function will also compute the effective system matrix
220  // K~=K+a_0*M+a_1*C and apply user specified initial
221  // conditions.
222  t_system.assemble();
223 
224  // Now solve for each time step.
225  // For convenience, use a local buffer of the
226  // current time. But once this time is updated,
227  // also update the EquationSystems parameter
228  // Start with t_time = 0 and write a short header
229  // to the nodal result file
230  res_out << "# pressure at node " << res_node_no << "\n"
231  << "# time\tpressure\n"
232  << t_system.time << "\t" << 0 << std::endl;
233 
234 
235  for (unsigned int time_step=0; time_step<n_time_steps; time_step++)
236  {
237  // Update the time. Both here and in the
238  // EquationSystems object
239  t_system.time += delta_t;
240 
241  // Update the rhs.
242  t_system.update_rhs();
243 
244  // Impose essential boundary conditions.
245  // Not that since the matrix is only assembled once,
246  // the penalty parameter should be added to the matrix
247  // only in the first time step. The applied
248  // boundary conditions may be time-dependent and hence
249  // the rhs vector is considered in each time step.
250  if (time_step == 0)
251  {
252  // The local function fill_dirichlet_bc()
253  // may also set Dirichlet boundary conditions for the
254  // matrix. When you set the flag as shown below,
255  // the flag will return true. If you want it to return
256  // false, simply do not set it.
257  equation_systems.parameters.set<bool>("Newmark set BC for Matrix") = true;
258 
259  fill_dirichlet_bc(equation_systems, "Wave");
260 
261  // unset the flag, so that it returns false
262  equation_systems.parameters.set<bool>("Newmark set BC for Matrix") = false;
263  }
264  else
265  fill_dirichlet_bc(equation_systems, "Wave");
266 
267  // Solve the system "Wave".
268  t_system.solve();
269 
270  // After solving the system, write the solution
271  // to a GMV-formatted plot file.
272  // Do only for a few time steps.
273  if (time_step == 30 || time_step == 60 ||
274  time_step == 90 || time_step == 120)
275  {
276  std::ostringstream file_name;
277 
278 #ifdef LIBMESH_HAVE_VTK
279  file_name << "out_"
280  << std::setw(3)
281  << std::setfill('0')
282  << std::right
283  << time_step
284  << ".pvtu";
285 
286  VTKIO(mesh).write_equation_systems (file_name.str(), equation_systems);
287 #else
288 
289  file_name << "out."
290  << std::setw(3)
291  << std::setfill('0')
292  << std::right
293  << time_step
294  << ".gmv";
295 
296  GMVIO(mesh).write_equation_systems (file_name.str(), equation_systems);
297 #endif
298  }
299 
300  // Update the p, v and a.
301  t_system.update_u_v_a();
302 
303  // dof_no may not be local in parallel runs, so we may need a
304  // global displacement vector
305  NumericVector<Number> & displacement = t_system.get_vector("displacement");
306  std::vector<Number> global_displacement(displacement.size());
307  displacement.localize(global_displacement);
308 
309  // Write nodal results to file. The results can then
310  // be viewed with e.g. gnuplot (run gnuplot and type
311  // 'plot "pressure_node.res" with lines' in the command line)
312  res_out << t_system.time << "\t"
313  << global_displacement[dof_no]
314  << std::endl;
315  }
316 
317  // All done.
318  return 0;
319 }
320 
321 // This function assembles the system matrix and right-hand-side
322 // for our wave equation.
324  const std::string & system_name)
325 {
326  // It is a good idea to make sure we are assembling
327  // the proper system.
328  libmesh_assert_equal_to (system_name, "Wave");
329 
330  // Get a constant reference to the mesh object.
331  const MeshBase & mesh = es.get_mesh();
332 
333  // The dimension that we are running.
334  const unsigned int dim = mesh.mesh_dimension();
335 
336  // Copy the speed of sound and fluid density
337  // to a local variable.
338  const Real speed = es.parameters.get<Real>("speed");
339  const Real rho = es.parameters.get<Real>("fluid density");
340 
341  // Get a reference to our system, as before.
342  NewmarkSystem & t_system = es.get_system<NewmarkSystem> (system_name);
343 
344  // Get a constant reference to the Finite Element type
345  // for the first (and only) variable in the system.
346  FEType fe_type = t_system.get_dof_map().variable_type(0);
347 
348  // In here, we will add the element matrices to the
349  // @e additional matrices "stiffness_mass" and "damping"
350  // and the additional vector "force", not to the members
351  // "matrix" and "rhs". Therefore, get writable
352  // references to them.
353  SparseMatrix<Number> & stiffness = t_system.get_matrix("stiffness");
354  SparseMatrix<Number> & damping = t_system.get_matrix("damping");
355  SparseMatrix<Number> & mass = t_system.get_matrix("mass");
356  NumericVector<Number> & force = t_system.get_vector("force");
357 
358  // Some solver packages (PETSc) are especially picky about
359  // allocating sparsity structure and truly assigning values
360  // to this structure. Namely, matrix additions, as performed
361  // later, exhibit acceptable performance only for identical
362  // sparsity structures. Therefore, explicitly zero the
363  // values in the collective matrix, so that matrix additions
364  // encounter identical sparsity structures.
365  SparseMatrix<Number> & matrix = *t_system.matrix;
366  DenseMatrix<Number> zero_matrix;
367 
368  // Build a Finite Element object of the specified type. Since the
369  // FEBase::build() member dynamically creates memory we will
370  // store the object as a std::unique_ptr<FEBase>. This can be thought
371  // of as a pointer that will clean up after itself.
372  std::unique_ptr<FEBase> fe (FEBase::build(dim, fe_type));
373 
374  // A 2nd order Gauss quadrature rule for numerical integration.
375  QGauss qrule (dim, SECOND);
376 
377  // Tell the finite element object to use our quadrature rule.
378  fe->attach_quadrature_rule (&qrule);
379 
380  // The element Jacobian * quadrature weight at each integration point.
381  const std::vector<Real> & JxW = fe->get_JxW();
382 
383  // The element shape functions evaluated at the quadrature points.
384  const std::vector<std::vector<Real>> & phi = fe->get_phi();
385 
386  // The element shape function gradients evaluated at the quadrature
387  // points.
388  const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
389 
390  // A reference to the DofMap object for this system. The DofMap
391  // object handles the index translation from node and element numbers
392  // to degree of freedom numbers.
393  const DofMap & dof_map = t_system.get_dof_map();
394 
395  // The element mass, damping and stiffness matrices
396  // and the element contribution to the rhs.
397  DenseMatrix<Number> Ke, Ce, Me;
399 
400  // This vector will hold the degree of freedom indices for
401  // the element. These define where in the global system
402  // the element degrees of freedom get mapped.
403  std::vector<dof_id_type> dof_indices;
404 
405  // Now we will loop over all the elements in the mesh.
406  // We will compute the element matrix and right-hand-side
407  // contribution.
408  for (const auto & elem : mesh.active_local_element_ptr_range())
409  {
410  // Get the degree of freedom indices for the
411  // current element. These define where in the global
412  // matrix and right-hand-side this element will
413  // contribute to.
414  dof_map.dof_indices (elem, dof_indices);
415 
416  // Compute the element-specific data for the current
417  // element. This involves computing the location of the
418  // quadrature points (q_point) and the shape functions
419  // (phi, dphi) for the current element.
420  fe->reinit (elem);
421 
422  // Zero the element matrices and rhs before
423  // summing them. We use the resize member here because
424  // the number of degrees of freedom might have changed from
425  // the last element. Note that this will be the case if the
426  // element type is different (i.e. the last element was HEX8
427  // and now have a PRISM6).
428  {
429  const unsigned int n_dof_indices = dof_indices.size();
430 
431  Ke.resize (n_dof_indices, n_dof_indices);
432  Ce.resize (n_dof_indices, n_dof_indices);
433  Me.resize (n_dof_indices, n_dof_indices);
434  zero_matrix.resize (n_dof_indices, n_dof_indices);
435  Fe.resize (n_dof_indices);
436  }
437 
438  // Now loop over the quadrature points. This handles
439  // the numeric integration.
440  for (unsigned int qp=0; qp<qrule.n_points(); qp++)
441  {
442  // Now we will build the element matrix. This involves
443  // a double loop to integrate the test functions (i) against
444  // the trial functions (j).
445  for (std::size_t i=0; i<phi.size(); i++)
446  for (std::size_t j=0; j<phi.size(); j++)
447  {
448  Ke(i,j) += JxW[qp]*(dphi[i][qp]*dphi[j][qp]);
449  Me(i,j) += JxW[qp]*phi[i][qp]*phi[j][qp]
450  *1./(speed*speed);
451  } // end of the matrix summation loop
452  } // end of quadrature point loop
453 
454  // Now compute the contribution to the element matrix and the
455  // right-hand-side vector if the current element lies on the
456  // boundary.
457  {
458  // In this example no natural boundary conditions will
459  // be considered. The code is left here so it can easily
460  // be extended.
461  //
462  // don't do this for any side
463  for (auto side : elem->side_index_range())
464  if (!true)
465  // if (elem->neighbor_ptr(side) == nullptr)
466  {
467  // Declare a special finite element object for
468  // boundary integration.
469  std::unique_ptr<FEBase> fe_face (FEBase::build(dim, fe_type));
470 
471  // Boundary integration requires one quadrature rule,
472  // with dimensionality one less than the dimensionality
473  // of the element.
474  QGauss qface(dim-1, SECOND);
475 
476  // Tell the finite element object to use our
477  // quadrature rule.
478  fe_face->attach_quadrature_rule (&qface);
479 
480  // The value of the shape functions at the quadrature
481  // points.
482  const std::vector<std::vector<Real>> & phi_face = fe_face->get_phi();
483 
484  // The Jacobian * Quadrature Weight at the quadrature
485  // points on the face.
486  const std::vector<Real> & JxW_face = fe_face->get_JxW();
487 
488  // Compute the shape function values on the element
489  // face.
490  fe_face->reinit(elem, side);
491 
492  // Here we consider a normal acceleration acc_n=1 applied to
493  // the whole boundary of our mesh.
494  const Real acc_n_value = 1.0;
495 
496  // Loop over the face quadrature points for integration.
497  for (unsigned int qp=0; qp<qface.n_points(); qp++)
498  {
499  // Right-hand-side contribution due to prescribed
500  // normal acceleration.
501  for (std::size_t i=0; i<phi_face.size(); i++)
502  {
503  Fe(i) += acc_n_value*rho
504  *phi_face[i][qp]*JxW_face[qp];
505  }
506  } // end face quadrature point loop
507  } // end if (elem->neighbor_ptr(side) == nullptr)
508 
509  // In this example the Dirichlet boundary conditions will be
510  // imposed via penalty method after the
511  // system is assembled.
512 
513  } // end boundary condition section
514 
515  // If this assembly program were to be used on an adaptive mesh,
516  // we would have to apply any hanging node constraint equations
517  // by uncommenting the following lines:
518  // std::vector<unsigned int> dof_indicesC = dof_indices;
519  // std::vector<unsigned int> dof_indicesM = dof_indices;
520  // dof_map.constrain_element_matrix_and_vector (Ke, Fe, dof_indices);
521  // dof_map.constrain_element_matrix (Ce, dof_indicesC);
522  // dof_map.constrain_element_matrix (Me, dof_indicesM);
523 
524  // Finally, simply add the contributions to the additional
525  // matrices and vector.
526  stiffness.add_matrix (Ke, dof_indices);
527  damping.add_matrix (Ce, dof_indices);
528  mass.add_matrix (Me, dof_indices);
529 
530  force.add_vector (Fe, dof_indices);
531 
532  // For the overall matrix, explicitly zero the entries where
533  // we added values in the other ones, so that we have
534  // identical sparsity footprints.
535  matrix.add_matrix(zero_matrix, dof_indices);
536 
537  } // end of element loop
538 }
539 
540 // This function applies the initial conditions
542  const std::string & system_name)
543 {
544  // Get a reference to our system, as before
545  NewmarkSystem & t_system = es.get_system<NewmarkSystem> (system_name);
546 
547  // Numeric vectors for the pressure, velocity and acceleration
548  // values.
549  NumericVector<Number> & pres_vec = t_system.get_vector("displacement");
550  NumericVector<Number> & vel_vec = t_system.get_vector("velocity");
551  NumericVector<Number> & acc_vec = t_system.get_vector("acceleration");
552 
553  // Assume our fluid to be at rest, which would
554  // also be the default conditions in class NewmarkSystem,
555  // but let us do it explicitly here.
556  pres_vec.zero();
557  vel_vec.zero();
558  acc_vec.zero();
559 }
560 
561 // This function applies the Dirichlet boundary conditions
563  const std::string & system_name)
564 {
565  // It is a good idea to make sure we are assembling
566  // the proper system.
567  libmesh_assert_equal_to (system_name, "Wave");
568 
569  // Get a reference to our system, as before.
570  NewmarkSystem & t_system = es.get_system<NewmarkSystem> (system_name);
571 
572  // Get writable references to the overall matrix and vector.
573  SparseMatrix<Number> & matrix = *t_system.matrix;
574  NumericVector<Number> & rhs = *t_system.rhs;
575 
576  // Get a constant reference to the mesh object.
577  const MeshBase & mesh = es.get_mesh();
578 
579  // Get libMesh's pi
580  const Real pi = libMesh::pi;
581 
582  // Ask the EquationSystems flag whether
583  // we should do this also for the matrix
584  const bool do_for_matrix =
585  es.parameters.get<bool>("Newmark set BC for Matrix");
586 
587  // Number of nodes in the mesh.
588  unsigned int n_nodes = mesh.n_nodes();
589 
590  for (unsigned int n_cnt=0; n_cnt<n_nodes; n_cnt++)
591  {
592  // Get a reference to the current node.
593  const Node & curr_node = mesh.node_ref(n_cnt);
594 
595  // Check if Dirichlet BCs should be applied to this node.
596  // Use the TOLERANCE from mesh_common.h as tolerance.
597  // Here a pressure value is applied if the z-coord.
598  // is equal to 4, which corresponds to one end of the
599  // pipe-mesh in this directory.
600  const Real z_coo = 4.;
601 
602  if (std::abs(curr_node(2)-z_coo) < TOLERANCE)
603  {
604  // The global number of the respective degree of freedom.
605  unsigned int dn = curr_node.dof_number(0, 0, 0);
606 
607  // The penalty parameter.
608  const Real penalty = 1.e10;
609 
610  // Here we apply sinusoidal pressure values for 0<t<0.002
611  // at one end of the pipe-mesh.
612  Real p_value;
613  if (t_system.time < .002)
614  p_value = sin(2*pi*t_system.time/.002);
615  else
616  p_value = .0;
617 
618  // Now add the contributions to the matrix and the rhs.
619  rhs.add(dn, p_value*penalty);
620 
621  // Add the penalty parameter to the global matrix
622  // if desired.
623  if (do_for_matrix)
624  matrix.add(dn, dn, penalty);
625  }
626  } // loop n_cnt
627 }
libMesh::NumericVector::zero
virtual void zero()=0
Set all entries to zero.
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::EquationSystems::add_system
virtual System & add_system(const std::string &system_type, const std::string &name)
Add the system of type system_type named name to the systems array.
Definition: equation_systems.C:345
libMesh::EquationSystems::get_mesh
const MeshBase & get_mesh() const
Definition: equation_systems.h:637
assemble_wave
void assemble_wave(EquationSystems &es, const std::string &system_name)
Definition: transient_ex2.C:323
libMesh::MeshBase::read
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.
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::PETSC_SOLVERS
Definition: enum_solver_package.h:36
libMesh::DofMap::dof_indices
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:1967
libMesh::MeshBase::active_local_element_ptr_range
virtual SimpleRange< element_iterator > active_local_element_ptr_range()=0
libMesh::QBase::n_points
unsigned int n_points() const
Definition: quadrature.h:126
libMesh::NewmarkSystem
This class contains a specific system class.
Definition: newmark_system.h:51
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::EquationSystems::get_system
const T_sys & get_system(const std::string &name) const
Definition: equation_systems.h:757
libMesh::TOLERANCE
static const Real TOLERANCE
Definition: libmesh_common.h:128
libMesh::DenseMatrix< Number >
libMesh::default_solver_package
SolverPackage default_solver_package()
Definition: libmesh.C:993
mesh
MeshBase & mesh
Definition: mesh_communication.C:1257
libMesh::MeshBase::mesh_dimension
unsigned int mesh_dimension() const
Definition: mesh_base.C:135
libMesh::GMVIO
This class implements writing meshes in the GMV format.
Definition: gmv_io.h:54
libMesh::VTKIO
This class implements reading and writing meshes in the VTK format.
Definition: vtk_io.h:60
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::size
virtual numeric_index_type size() const =0
libMesh::ReplicatedMesh
The ReplicatedMesh class is derived from the MeshBase class, and is used to store identical copies of...
Definition: replicated_mesh.h:47
libMesh::NumericVector< Number >
libMesh::TriangleWrapper::init
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
libMesh::DofObject::dof_number
dof_id_type dof_number(const unsigned int s, const unsigned int var, const unsigned int comp) const
Definition: dof_object.h:956
libMesh::NumericVector::add_vector
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].
Definition: numeric_vector.C:363
libMesh::MeshBase
This is the MeshBase class.
Definition: mesh_base.h:78
std::abs
MetaPhysicL::DualNumber< T, D > abs(const MetaPhysicL::DualNumber< T, D > &in)
libMesh::NumericVector::localize
virtual void localize(std::vector< T > &v_local) const =0
Creates a copy of the global vector in the local vector v_local.
libMesh::INVALID_SOLVER_PACKAGE
Definition: enum_solver_package.h:43
libMesh::FEGenericBase::build
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
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::EquationSystems::init
virtual void init()
Initialize all the systems.
Definition: equation_systems.C:96
libMesh::ImplicitSystem::matrix
SparseMatrix< Number > * matrix
The system matrix.
Definition: implicit_system.h:393
libMesh::Node
A Node is like a Point, but with more information.
Definition: node.h:52
apply_initial
void apply_initial(EquationSystems &es, const std::string &system_name)
Definition: transient_ex2.C:541
libMesh::LibMeshInit
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
Definition: libmesh.h:83
n_nodes
const dof_id_type n_nodes
Definition: tecplot_io.C:68
libMesh::EquationSystems::print_info
void print_info(std::ostream &os=libMesh::out) const
Prints information about the equation systems, by default to libMesh::out.
Definition: equation_systems.C:1314
libMesh::EquationSystems
This is the EquationSystems class.
Definition: equation_systems.h:74
libMesh::System::time
Real time
For time-dependent problems, this is the time t at the beginning of the current timestep.
Definition: system.h:1561
libMesh::Parameters::set
T & set(const std::string &)
Definition: parameters.h:460
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::MeshBase::node_ref
virtual const Node & node_ref(const dof_id_type i) const
Definition: mesh_base.h:451
main
int main(int argc, char **argv)
Definition: transient_ex2.C:100
libMesh::MeshBase::n_nodes
virtual dof_id_type n_nodes() const =0
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::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
fill_dirichlet_bc
void fill_dirichlet_bc(EquationSystems &es, const std::string &system_name)
Definition: transient_ex2.C:562
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::Parameters::get
const T & get(const std::string &) const
Definition: parameters.h:421
libMesh::out
OStreamProxy out
libMesh::System::get_vector
const NumericVector< Number > & get_vector(const std::string &vec_name) const
Definition: system.C:774
libMesh::FIRST
Definition: enum_order.h:42
libMesh::EIGEN_SOLVERS
Definition: enum_solver_package.h:40
libMesh::DenseVector< Number >
libMesh::EquationSystems::parameters
Parameters parameters
Data structure holding arbitrary parameters.
Definition: equation_systems.h:557