libMesh
subdomains_ex3.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>Subdomains Example 3 - Integrating discontinuous data that cuts the mesh</h1>
21 // \author Benjamin S. Kirk
22 // \date 2013
23 
24 
25 // C++ include files that we need
26 #include <iostream>
27 #include <algorithm>
28 #include <math.h>
29 
30 // Basic include file needed for the mesh functionality.
31 #include "libmesh/libmesh.h"
32 #include "libmesh/mesh.h"
33 #include "libmesh/mesh_generation.h"
34 #include "libmesh/mesh_refinement.h"
35 #include "libmesh/quadrature_gauss.h"
36 #include "libmesh/quadrature_composite.h"
37 #include "libmesh/fe.h"
38 #include "libmesh/elem.h"
39 
40 // Bring in everything from the libMesh namespace
41 using namespace libMesh;
42 
43 // declare the functions we will use
44 void integrate_function (const MeshBase & mesh);
45 
46 // signed distance function
47 const Real radius = 0.5;
48 
49 Real distance (const Point & p)
50 {
51  Point cent(0.8, 0.9);
52  return ((p-cent).norm() - radius);
53 }
54 
55 Real integrand (const Point & p)
56 {
57  return (distance(p) < 0) ? 10. : 1.;
58 }
59 
60 
61 
62 // Begin the main program.
63 int main (int argc, char ** argv)
64 {
65  // Initialize libMesh and any dependent libraries, like in example 2.
66  LibMeshInit init (argc, argv);
67 
68  // This example requires Adaptive Mesh Refinement support - although
69  // it only refines uniformly, the refinement code used is the same
70  // underneath. It also requires libmesh support for Triangle and
71  // Tetgen, which means that libmesh must be configured with
72  // --disable-strict-lgpl for this example to run.
73 #if !defined(LIBMESH_HAVE_TRIANGLE) || !defined(LIBMESH_HAVE_TETGEN) || !defined(LIBMESH_ENABLE_AMR)
74  libmesh_example_requires(false, "--disable-strict-lgpl --enable-amr");
75 #else
76 
77  // Skip this 2D example if libMesh was compiled as 1D-only.
78  libmesh_example_requires(2 <= LIBMESH_DIM, "2D support");
79 
80  // Read the mesh from file. This is the coarse mesh that will be used
81  // in example 10 to demonstrate adaptive mesh refinement. Here we will
82  // simply read it in and uniformly refine it 5 times before we compute
83  // with it.
84  Mesh mesh(init.comm());
85 
86  {
87  unsigned int dim=2;
88 
89  if (argc == 3 && std::atoi(argv[2]) == 3)
90  {
91  libmesh_here();
92  dim=3;
93  }
94 
95  mesh.read ((dim==2) ? "mesh.xda" : "hybrid_3d.xda");
96  }
97 
98  // Create a MeshRefinement object to handle refinement of our mesh.
99  // This class handles all the details of mesh refinement and coarsening.
100  MeshRefinement mesh_refinement (mesh);
101 
102  // Uniformly refine the mesh 4 times. This is the
103  // first time we use the mesh refinement capabilities
104  // of the library.
105  mesh_refinement.uniformly_refine (4);
106 
107  // Print information about the mesh to the screen.
108  mesh.print_info();
109 
110  // integrate the desired function
112 
113  // All done.
114  return 0;
115 #endif
116 }
117 
118 
119 
121 {
122 #if defined(LIBMESH_HAVE_TRIANGLE) && defined(LIBMESH_HAVE_TETGEN)
123 
124  std::vector<Real> vertex_distance;
125 
126  QComposite<QGauss> qrule (mesh.mesh_dimension(), FIRST);
127  //QGauss qrule (mesh.mesh_dimension(), FIRST);
128 
130 
131  Real int_val=0.;
132 
133  const std::vector<Point> & q_points = fe->get_xyz();
134  const std::vector<Real> & JxW = fe->get_JxW();
135 
136  for (const auto & elem : mesh.active_local_element_ptr_range())
137  {
138  vertex_distance.clear();
139 
140  for (unsigned int v=0; v<elem->n_vertices(); v++)
141  vertex_distance.push_back (distance(elem->point(v)));
142 
143  qrule.init (*elem, vertex_distance);
144 
145  fe->reinit (elem,
146  &(qrule.get_points()),
147  &(qrule.get_weights()));
148 
149 
150  // TODO: would it be valuable to have the composite quadrature rule sort
151  // from smallest to largest JxW value to help prevent
152  // ... large + small + large + large + small ...
153  // type truncation errors?
154  for (std::size_t qp=0; qp<q_points.size(); qp++)
155  int_val += JxW[qp] * integrand(q_points[qp]);
156  }
157 
158  mesh.comm().sum (int_val);
159 
160  libMesh::out << "\n***********************************\n"
161  << " int_val = " << int_val << std::endl
162  << " exact_val = " << 1*(2*2 - radius*radius*pi) + 10.*(radius*radius*pi)
163  << "\n***********************************\n"
164  << std::endl;
165 #else
166  libmesh_ignore(mesh);
167 #endif
168 }
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:178
const Real radius
unsigned int dim
MeshBase & mesh
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.
Real distance(const Point &p)
Real integrand(const Point &p)
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
void integrate_function(const MeshBase &mesh)
This is the MeshRefinement class.
void init(triangulateio &t)
Initializes the fields of t to NULL/0 as necessary.
const Parallel::Communicator & comm() const
Definition: libmesh.h:81
void libmesh_ignore(const T &)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const Parallel::Communicator & comm() const
OStreamProxy out
unsigned int mesh_dimension() const
Definition: mesh_base.C:148
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
Definition: mesh.h:50
virtual void read(const std::string &name, void *mesh_data=libmesh_nullptr, bool skip_renumber_nodes_and_elements=false, bool skip_find_neighbors=false)=0
Interfaces for reading/writing a mesh to/from a file.
void sum(T &r) const
Take a local variable and replace it with the sum of it&#39;s values on all processors.
int main(int argc, char **argv)
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:38
void print_info(std::ostream &os=libMesh::out) const
Prints relevant information about the mesh.
Definition: mesh_base.C:448
const Real pi
.
Definition: libmesh.h:172
void uniformly_refine(unsigned int n=1)
Uniformly refines the mesh n times.
static UniquePtr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.