libMesh
miscellaneous_ex1.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2017 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>Miscellaneous Example 1 - Infinite Elements for the Wave Equation</h1>
21 // \author Daniel Dreyer
22 // \date 2003
23 //
24 // This is the sixth example program. It builds on
25 // the previous examples, and introduces the Infinite
26 // Element class. Note that the library must be compiled
27 // with Infinite Elements enabled. Otherwise, this
28 // example will abort.
29 // This example intends to demonstrate the similarities
30 // between the FE and the InfFE classes in libMesh.
31 // The matrices are assembled according to the wave equation.
32 // However, for practical applications a time integration
33 // scheme (as introduced in subsequent examples) should be
34 // used.
35 
36 // C++ include files that we need
37 #include <iostream>
38 #include <algorithm>
39 #include <math.h>
40 
41 // Basic include file needed for the mesh functionality.
42 #include "libmesh/exodusII_io.h"
43 #include "libmesh/libmesh.h"
44 #include "libmesh/replicated_mesh.h"
45 #include "libmesh/mesh_generation.h"
46 #include "libmesh/linear_implicit_system.h"
47 #include "libmesh/equation_systems.h"
48 
49 // Define the Finite and Infinite Element object.
50 #include "libmesh/fe.h"
51 #include "libmesh/inf_fe.h"
52 #include "libmesh/inf_elem_builder.h"
53 
54 // Define Gauss quadrature rules.
55 #include "libmesh/quadrature_gauss.h"
56 
57 // Define useful datatypes for finite element
58 // matrix and vector components.
59 #include "libmesh/sparse_matrix.h"
60 #include "libmesh/numeric_vector.h"
61 #include "libmesh/dense_matrix.h"
62 #include "libmesh/dense_vector.h"
63 
64 // Define the DofMap, which handles degree of freedom
65 // indexing.
66 #include "libmesh/dof_map.h"
67 
68 // The definition of a vertex associated with a Mesh.
69 #include "libmesh/node.h"
70 
71 // The definition of a geometric element
72 #include "libmesh/elem.h"
73 
74 // Bring in everything from the libMesh namespace
75 using namespace libMesh;
76 
77 // Function prototype. This is similar to the Poisson
78 // assemble function of example 4.
80  const std::string & system_name);
81 
82 // Begin the main program.
83 int main (int argc, char ** argv)
84 {
85  // Initialize libMesh, like in example 2.
86  LibMeshInit init (argc, argv);
87 
88  // This example requires Infinite Elements
89 #ifndef LIBMESH_ENABLE_INFINITE_ELEMENTS
90  libmesh_example_requires(false, "--enable-ifem");
91 #else
92 
93  // Skip this 3D example if libMesh was compiled as 1D/2D-only.
94  libmesh_example_requires(3 <= LIBMESH_DIM, "3D support");
95 
96  // Tell the user what we are doing.
97  libMesh::out << "Running ex6 with dim = 3" << std::endl << std::endl;
98 
99  // Create a serialized mesh, distributed across the default MPI
100  // communicator.
101  // InfElemBuilder still requires some updates to be DistributedMesh
102  // compatible
103 
104  ReplicatedMesh mesh(init.comm());
105 
106  // Use the internal mesh generator to create elements
107  // on the square [-1,1]^3, of type Hex8.
109  4, 4, 4,
110  -1., 1.,
111  -1., 1.,
112  -1., 1.,
113  HEX8);
114 
115  // Print information about the mesh to the screen.
116  mesh.print_info();
117 
118  // Write the mesh before the infinite elements are added
119 #ifdef LIBMESH_HAVE_EXODUS_API
120  ExodusII_IO(mesh).write ("orig_mesh.e");
121 #endif
122 
123  // Normally, when a mesh is imported or created in
124  // libMesh, only conventional elements exist. The infinite
125  // elements used here, however, require prescribed
126  // nodal locations (with specified distances from an imaginary
127  // origin) and configurations that a conventional mesh creator
128  // in general does not offer. Therefore, an efficient method
129  // for building infinite elements is offered. It can account
130  // for symmetry planes and creates infinite elements in a fully
131  // automatic way.
132  //
133  // Right now, the simplified interface is used, automatically
134  // determining the origin. Check MeshBase for a generalized
135  // method that can even return the element faces of interior
136  // vibrating surfaces. The bool determines whether to be
137  // verbose.
138  InfElemBuilder builder(mesh);
139  builder.build_inf_elem(true);
140 
141  // Reassign subdomain_id() of all infinite elements.
142  // Otherwise, the exodus-api will fail on the mesh.
143  for (auto & elem : mesh.element_ptr_range())
144  if (elem->infinite())
145  elem->subdomain_id() = 1;
146 
147  // Print information about the mesh to the screen.
148  mesh.print_info();
149 
150  // Write the mesh with the infinite elements added.
151  // Compare this to the original mesh.
152 #ifdef LIBMESH_HAVE_EXODUS_API
153  ExodusII_IO(mesh).write ("ifems_added.e");
154 #endif
155 
156  // After building infinite elements, we have to let
157  // the elements find their neighbors again.
159 
160  // Create an equation systems object, where ThinSystem
161  // offers only the crucial functionality for solving a
162  // system. Use ThinSystem when you want the sleekest
163  // system possible.
164  EquationSystems equation_systems (mesh);
165 
166  // Declare the system and its variables.
167  // Create a system named "Wave". This can
168  // be a simple, steady system
169  equation_systems.add_system<LinearImplicitSystem> ("Wave");
170 
171  // Create an FEType describing the approximation
172  // characteristics of the InfFE object. Note that
173  // the constructor automatically defaults to some
174  // sensible values. But use FIRST order
175  // approximation.
176  FEType fe_type(FIRST);
177 
178  // Add the variable "p" to "Wave". Note that there exist
179  // various approaches in adding variables. In example 3,
180  // add_variable took the order of approximation and used
181  // default values for the FEFamily, while here the FEType
182  // is used.
183  equation_systems.get_system("Wave").add_variable("p", fe_type);
184 
185  // Give the system a pointer to the matrix assembly
186  // function.
187  equation_systems.get_system("Wave").attach_assemble_function (assemble_wave);
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") = 1.;
193  equation_systems.parameters.set<Real>("fluid density") = 1.;
194 
195  // Initialize the data structures for the equation system.
196  equation_systems.init();
197 
198  // Prints information about the system to the screen.
199  equation_systems.print_info();
200 
201  // Solve the system "Wave".
202  equation_systems.get_system("Wave").solve();
203 
204  // Write the whole EquationSystems object to file.
205  // For infinite elements, the concept of nodal_soln()
206  // is not applicable. Therefore, writing the mesh in
207  // some format @e always gives all-zero results at
208  // the nodes of the infinite elements. Instead,
209  // use the FEInterface::compute_data() methods to
210  // determine physically correct results within an
211  // infinite element.
212  equation_systems.write ("eqn_sys.dat", WRITE);
213 
214  // All done.
215  return 0;
216 
217 #endif // else part of ifndef LIBMESH_ENABLE_INFINITE_ELEMENTS
218 }
219 
220 // This function assembles the system matrix and right-hand-side
221 // for the discrete form of our wave equation.
223  const std::string & libmesh_dbg_var(system_name))
224 {
225  // It is a good idea to make sure we are assembling
226  // the proper system.
227  libmesh_assert_equal_to (system_name, "Wave");
228 
229  // Avoid unused variable warnings when compiling without infinite
230  // elements enabled.
231  libmesh_ignore(es);
232 
233 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
234 
235  // Get a constant reference to the mesh object.
236  const MeshBase & mesh = es.get_mesh();
237 
238  // Get a reference to the system we are solving.
240 
241  // A reference to the DofMap object for this system. The DofMap
242  // object handles the index translation from node and element numbers
243  // to degree of freedom numbers.
244  const DofMap & dof_map = system.get_dof_map();
245 
246  // The dimension that we are running.
247  const unsigned int dim = mesh.mesh_dimension();
248 
249  // Copy the speed of sound to a local variable.
250  const Real speed = es.parameters.get<Real>("speed");
251 
252  // Get a constant reference to the Finite Element type
253  // for the first (and only) variable in the system.
254  const FEType & fe_type = dof_map.variable_type(0);
255 
256  // Build a Finite Element object of the specified type. Since the
257  // FEBase::build() member dynamically creates memory we will
258  // store the object as a UniquePtr<FEBase>.
259  UniquePtr<FEBase> fe (FEBase::build(dim, fe_type));
260 
261  // Do the same for an infinite element.
262  UniquePtr<FEBase> inf_fe (FEBase::build_InfFE(dim, fe_type));
263 
264  // A 2nd order Gauss quadrature rule for numerical integration.
265  QGauss qrule (dim, SECOND);
266 
267  // Tell the finite element object to use our quadrature rule.
268  fe->attach_quadrature_rule (&qrule);
269 
270  // Due to its internal structure, the infinite element handles
271  // quadrature rules differently. It takes the quadrature
272  // rule which has been initialized for the FE object, but
273  // creates suitable quadrature rules by @e itself. The user
274  // need not worry about this.
275  inf_fe->attach_quadrature_rule (&qrule);
276 
277  // Define data structures to contain the element matrix
278  // and right-hand-side vector contribution. Following
279  // basic finite element terminology we will denote these
280  // "Ke", "Ce", "Me", and "Fe" for the stiffness, damping
281  // and mass matrices, and the load vector. Note that in
282  // Acoustics, these descriptors though do @e not match the
283  // true physical meaning of the projectors. The final
284  // overall system, however, resembles the conventional
285  // notation again.
290 
291  // This vector will hold the degree of freedom indices for
292  // the element. These define where in the global system
293  // the element degrees of freedom get mapped.
294  std::vector<dof_id_type> dof_indices;
295 
296  // Now we will loop over all the elements in the mesh.
297  // We will compute the element matrix and right-hand-side
298  // contribution.
299  for (const auto & elem : mesh.active_local_element_ptr_range())
300  {
301  // Get the degree of freedom indices for the
302  // current element. These define where in the global
303  // matrix and right-hand-side this element will
304  // contribute to.
305  dof_map.dof_indices (elem, dof_indices);
306 
307  // The mesh contains both finite and infinite elements. These
308  // elements are handled through different classes, namely
309  // FE and InfFE, respectively. However, since both
310  // are derived from FEBase, they share the same interface,
311  // and overall burden of coding is @e greatly reduced through
312  // using a pointer, which is adjusted appropriately to the
313  // current element type.
314  FEBase * cfe = libmesh_nullptr;
315 
316  // This here is almost the only place where we need to
317  // distinguish between finite and infinite elements.
318  // For faster computation, however, different approaches
319  // may be feasible.
320  //
321  // Up to now, we do not know what kind of element we
322  // have. Aske the element of what type it is:
323  if (elem->infinite())
324  {
325  // We have an infinite element. Let cfe point
326  // to our InfFE object. This is handled through
327  // a UniquePtr. Through the UniquePtr::get() we "borrow"
328  // the pointer, while the UniquePtr inf_fe is
329  // still in charge of memory management.
330  cfe = inf_fe.get();
331  }
332  else
333  {
334  // This is a conventional finite element. Let fe handle it.
335  cfe = fe.get();
336 
337  // Boundary conditions.
338  // Here we just zero the rhs-vector. For natural boundary
339  // conditions check e.g. previous examples.
340  {
341  // Zero the RHS for this element.
342  Fe.resize (dof_indices.size());
343 
344  system.rhs->add_vector (Fe, dof_indices);
345  } // end boundary condition section
346  } // else if (elem->infinite())
347 
348  // This is slightly different from the Poisson solver:
349  // Since the finite element object may change, we have to
350  // initialize the constant references to the data fields
351  // each time again, when a new element is processed.
352  //
353  // The element Jacobian * quadrature weight at each integration point.
354  const std::vector<Real> & JxW = cfe->get_JxW();
355 
356  // The element shape functions evaluated at the quadrature points.
357  const std::vector<std::vector<Real>> & phi = cfe->get_phi();
358 
359  // The element shape function gradients evaluated at the quadrature
360  // points.
361  const std::vector<std::vector<RealGradient>> & dphi = cfe->get_dphi();
362 
363  // The infinite elements need more data fields than conventional FE.
364  // These are the gradients of the phase term dphase, an additional
365  // radial weight for the test functions Sobolev_weight, and its
366  // gradient.
367  //
368  // Note that these data fields are also initialized appropriately by
369  // the FE method, so that the weak form (below) is valid for @e both
370  // finite and infinite elements.
371  const std::vector<RealGradient> & dphase = cfe->get_dphase();
372  const std::vector<Real> & weight = cfe->get_Sobolev_weight();
373  const std::vector<RealGradient> & dweight = cfe->get_Sobolev_dweight();
374 
375  // Now this is all independent of whether we use an FE
376  // or an InfFE. Nice, hm? ;-)
377  //
378  // Compute the element-specific data, as described
379  // in previous examples.
380  cfe->reinit (elem);
381 
382  // Zero the element matrices. Boundary conditions were already
383  // processed in the FE-only section, see above.
384  Ke.resize (dof_indices.size(), dof_indices.size());
385  Ce.resize (dof_indices.size(), dof_indices.size());
386  Me.resize (dof_indices.size(), dof_indices.size());
387 
388  // The total number of quadrature points for infinite elements
389  // @e has to be determined in a different way, compared to
390  // conventional finite elements. This type of access is also
391  // valid for finite elements, so this can safely be used
392  // anytime, instead of asking the quadrature rule, as
393  // seen in previous examples.
394  unsigned int max_qp = cfe->n_quadrature_points();
395 
396  // Loop over the quadrature points.
397  for (unsigned int qp=0; qp<max_qp; qp++)
398  {
399  // Similar to the modified access to the number of quadrature
400  // points, the number of shape functions may also be obtained
401  // in a different manner. This offers the great advantage
402  // of being valid for both finite and infinite elements.
403  const unsigned int n_sf = cfe->n_shape_functions();
404 
405  // Now we will build the element matrices. Since the infinite
406  // elements are based on a Petrov-Galerkin scheme, the
407  // resulting system matrices are non-symmetric. The additional
408  // weight, described before, is part of the trial space.
409  //
410  // For the finite elements, though, these matrices are symmetric
411  // just as we know them, since the additional fields dphase,
412  // weight, and dweight are initialized appropriately.
413  //
414  // test functions: weight[qp]*phi[i][qp]
415  // trial functions: phi[j][qp]
416  // phase term: phase[qp]
417  //
418  // derivatives are similar, but note that these are of type
419  // Point, not of type Real.
420  for (unsigned int i=0; i<n_sf; i++)
421  for (unsigned int j=0; j<n_sf; j++)
422  {
423  // (ndt*Ht + nHt*d) * nH
424  Ke(i,j) +=
425  (
426  (dweight[qp] * phi[i][qp] // Point * Real = Point
427  + // +
428  dphi[i][qp] * weight[qp] // Point * Real = Point
429  ) * dphi[j][qp]
430  ) * JxW[qp];
431 
432  // (d*Ht*nmut*nH - ndt*nmu*Ht*H - d*nHt*nmu*H)
433  Ce(i,j) +=
434  (
435  (dphase[qp] * dphi[j][qp]) // (Point * Point) = Real
436  * weight[qp] * phi[i][qp] // * Real * Real = Real
437  - // -
438  (dweight[qp] * dphase[qp]) // (Point * Point) = Real
439  * phi[i][qp] * phi[j][qp] // * Real * Real = Real
440  - // -
441  (dphi[i][qp] * dphase[qp]) // (Point * Point) = Real
442  * weight[qp] * phi[j][qp] // * Real * Real = Real
443  ) * JxW[qp];
444 
445  // (d*Ht*H * (1 - nmut*nmu))
446  Me(i,j) +=
447  (
448  (1. - (dphase[qp] * dphase[qp])) // (Real - (Point * Point)) = Real
449  * phi[i][qp] * phi[j][qp] * weight[qp] // * Real * Real * Real = Real
450  ) * JxW[qp];
451 
452  } // end of the matrix summation loop
453  } // end of quadrature point loop
454 
455  // The element matrices are now built for this element.
456  // Collect them in Ke, and then add them to the global matrix.
457  // The SparseMatrix::add_matrix() member does this for us.
458  Ke.add(1./speed , Ce);
459  Ke.add(1./(speed*speed), Me);
460 
461  // If this assembly program were to be used on an adaptive mesh,
462  // we would have to apply any hanging node constraint equations
463  dof_map.constrain_element_matrix(Ke, dof_indices);
464 
465  system.matrix->add_matrix (Ke, dof_indices);
466  } // end of element loop
467 
468  // Note that we have not applied any boundary conditions so far.
469  // Here we apply a unit load at the node located at (0,0,0).
470  for (const auto & node : mesh.local_node_ptr_range())
471  if (std::abs((*node)(0)) < TOLERANCE &&
472  std::abs((*node)(1)) < TOLERANCE &&
473  std::abs((*node)(2)) < TOLERANCE)
474  {
475  // The global number of the respective degree of freedom.
476  unsigned int dn = node->dof_number(0,0,0);
477 
478  system.rhs->add (dn, 1.);
479  }
480 
481 #else
482 
483  // dummy assert
484  libmesh_assert_not_equal_to (es.get_mesh().mesh_dimension(), 1);
485 
486 #endif //ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
487 }
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:178
double abs(double a)
This is the EquationSystems class.
int main(int argc, char **argv)
The ReplicatedMesh class is derived from the MeshBase class, and is used to store identical copies of...
void write(const std::string &name, const XdrMODE, const unsigned int write_flags=(WRITE_DATA), bool partition_agnostic=true) const
Write the systems to disk using the XDR data format.
virtual SimpleRange< node_iterator > local_node_ptr_range()=0
unsigned int dim
virtual unsigned int n_quadrature_points() const =0
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
Definition: exodusII_io.h:52
This class provides a specific system class.
void resize(const unsigned int n)
Resize the vector.
Definition: dense_vector.h:350
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 class libmesh_nullptr_t libmesh_nullptr
NumericVector< Number > * rhs
The system matrix.
const T & get(const std::string &) const
Definition: parameters.h:430
static const Real TOLERANCE
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
Definition: libmesh.h:62
The libMesh namespace provides an interface to certain functionality in the library.
This is the MeshBase class.
Definition: mesh_base.h:68
std::unique_ptr< T > UniquePtr
Definition: auto_ptr.h:46
virtual SimpleRange< element_iterator > active_local_element_ptr_range()=0
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:167
virtual unsigned int n_shape_functions() const =0
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 find_neighbors(const bool reset_remote_elements=false, const bool reset_current_list=true)=0
Locate element face (edge in 2D) neighbors.
const DofMap & get_dof_map() const
Definition: system.h:2030
const std::vector< Real > & get_JxW() const
Definition: fe_abstract.h:237
const std::vector< Real > & get_Sobolev_weight() const
Definition: fe_base.h:420
virtual SimpleRange< element_iterator > element_ptr_range()=0
dof_id_type weight(const MeshBase &mesh, const processor_id_type pid)
Definition: mesh_tools.C:245
void print_info(std::ostream &os=libMesh::out) const
Prints information about the equation systems, by default to libMesh::out.
boostcopy::enable_if_c< ScalarTraits< T2 >::value, void >::type add(const T2 factor, const DenseMatrix< T3 > &mat)
Adds factor times mat to this matrix.
Definition: dense_matrix.h:883
const std::vector< std::vector< OutputShape > > & get_phi() const
Definition: fe_base.h:208
void init(triangulateio &t)
Initializes the fields of t to NULL/0 as necessary.
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.
const Parallel::Communicator & comm() const
Definition: libmesh.h:81
PetscErrorCode Vec Mat libmesh_dbg_var(j)
static UniquePtr< FEGenericBase > build_InfFE(const unsigned int dim, const FEType &type)
Builds a specific infinite element type.
void libmesh_ignore(const T &)
const std::vector< RealGradient > & get_Sobolev_dweight() const
Definition: fe_base.h:428
This class is used to build infinite elements on top of an existing mesh.
void assemble_wave(EquationSystems &es, const std::string &system_name)
const std::vector< std::vector< OutputGradient > > & get_dphi() const
Definition: fe_base.h:216
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const Point build_inf_elem(const bool be_verbose=false)
Build infinite elements atop a volume-based mesh, determine origin automatically. ...
virtual void write(const std::string &fname) libmesh_override
This method implements writing a mesh to a specified file.
Definition: exodusII_io.C:789
T & set(const std::string &)
Definition: parameters.h:469
OStreamProxy out
SparseMatrix< Number > * matrix
The system matrix.
unsigned int mesh_dimension() const
Definition: mesh_base.C:148
void resize(const unsigned int new_m, const unsigned int new_n)
Resize the matrix.
Definition: dense_matrix.h:776
This class implements specific orders of Gauss quadrature.
Parameters parameters
Data structure holding arbitrary parameters.
const MeshBase & get_mesh() const
const T_sys & get_system(const std::string &name) const
virtual void init()
Initialize all the systems.
virtual void add(const numeric_index_type i, const T value)=0
Adds value to each entry of the vector.
void print_info(std::ostream &os=libMesh::out) const
Prints relevant information about the mesh.
Definition: mesh_base.C:448
void build_cube(UnstructuredMesh &mesh, const unsigned int nx=0, const unsigned int ny=0, const unsigned int nz=0, const Real xmin=0., const Real xmax=1., const Real ymin=0., const Real ymax=1., const Real zmin=0., const Real zmax=1., const ElemType type=INVALID_ELEM, const bool gauss_lobatto_grid=false)
Builds a (elements) cube.
This class forms the foundation from which generic finite elements may be derived.
virtual void reinit(const Elem *elem, const std::vector< Point > *const pts=libmesh_nullptr, const std::vector< Real > *const weights=libmesh_nullptr)=0
This is at the core of this class.
const std::vector< OutputGradient > & get_dphase() const
Definition: fe_base.h:404
static UniquePtr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.