queso-0.53.0
LibMeshNegativeLaplacianOperator.C
Go to the documentation of this file.
1 //-----------------------------------------------------------------------bl-
2 //--------------------------------------------------------------------------
3 //
4 // QUESO - a library to support the Quantification of Uncertainty
5 // for Estimation, Simulation and Optimization
6 //
7 // Copyright (C) 2008-2015 The PECOS Development Team
8 //
9 // This library is free software; you can redistribute it and/or
10 // modify it under the terms of the Version 2.1 GNU Lesser General
11 // Public License as published by the Free Software Foundation.
12 //
13 // This library is distributed in the hope that it will be useful,
14 // but WITHOUT ANY WARRANTY; without even the implied warranty of
15 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
16 // Lesser General Public License for more details.
17 //
18 // You should have received a copy of the GNU Lesser General Public
19 // License along with this library; if not, write to the Free Software
20 // Foundation, Inc. 51 Franklin Street, Fifth Floor,
21 // Boston, MA 02110-1301 USA
22 //
23 //-----------------------------------------------------------------------el-
24 
25 #include <queso/Defines.h>
26 
27 #ifdef QUESO_HAVE_LIBMESH
28 
29 #include <set>
30 #include <vector>
31 
32 #include <boost/shared_ptr.hpp>
33 
34 #include <queso/FunctionOperatorBuilder.h>
35 #include <queso/LibMeshNegativeLaplacianOperator.h>
36 #include <libmesh/libmesh_common.h>
37 #include <libmesh/mesh.h>
38 #include <libmesh/equation_systems.h>
39 #include <libmesh/condensed_eigen_system.h>
40 #include <libmesh/auto_ptr.h>
41 #include <libmesh/fe_base.h>
42 #include <libmesh/quadrature_gauss.h>
43 #include <libmesh/dense_matrix.h>
44 #include <libmesh/sparse_matrix.h>
45 #include <libmesh/dof_map.h>
46 #include <libmesh/id_types.h>
47 #include <libmesh/elem.h>
48 #include <libmesh/zero_function.h>
49 #include <libmesh/dirichlet_boundaries.h>
50 #include <libmesh/utility.h>
51 #include <libmesh/string_to_enum.h>
52 #include <libmesh/enum_order.h>
53 #include <libmesh/enum_fe_family.h>
54 
55 namespace QUESO {
56 
57 LibMeshNegativeLaplacianOperator::LibMeshNegativeLaplacianOperator(
58  const FunctionOperatorBuilder & builder, libMesh::MeshBase & m)
59  : LibMeshOperatorBase(builder, m)
60 {
61  boost::shared_ptr<libMesh::EquationSystems> es(this->equation_systems);
62 
63  // Give the system a pointer to the matrix assembly
64  // function defined below.
65  libMesh::CondensedEigenSystem & eigen_system =
66  es->get_system<libMesh::CondensedEigenSystem>("Eigensystem");
67 
68  // Declare the system variables.
69  // Adds the variable "u" to "Eigensystem". "u"
70  // will be approximated using second-order approximation.
71  unsigned int u_var = eigen_system.add_variable("u",
72  libMesh::Utility::string_to_enum<libMeshEnums::Order>(this->builder.order),
73  libMesh::Utility::string_to_enum<libMeshEnums::FEFamily>(this->builder.family));
74 
75  // This works because *this is a subclass of System::Assembly
76  // and requires the class to implement an 'assemble'
77  eigen_system.attach_assemble_object(*this);
78 
79  // Set the number of requested eigenpairs \p n_evals and the number
80  // of basis vectors used in the solution algorithm.
81  es->parameters.set<unsigned int>("eigenpairs") =
82  this->builder.num_req_eigenpairs;
83  es->parameters.set<unsigned int>("basis vectors") =
84  this->builder.num_req_eigenpairs * 3;
85 
86  // Set the solver tolerance and the maximum number of iterations.
87  es->parameters.set<libMesh::Real>("linear solver tolerance") =
88  pow(libMesh::TOLERANCE, 5./3.);
89  es->parameters.set<unsigned int>("linear solver maximum iterations") = 1000;
90 
91  // Set the type of the problem, here we deal with
92  // a generalized Hermitian problem.
93  eigen_system.set_eigenproblem_type(libMeshEnums::GHEP);
94 
95  // Order the eigenvalues "smallest first"
96  // This hoses performance?
97  eigen_system.eigen_solver->set_position_of_spectrum
98  (libMeshEnums::SMALLEST_MAGNITUDE);
99 
100  // Set up the boundary (only works if this->m is a square)
101  // We'll just the whole boundary to be Dirichlet, because why not
102  std::set<libMesh::boundary_id_type> boundary_ids;
103  boundary_ids.insert(0);
104  boundary_ids.insert(1);
105  boundary_ids.insert(2);
106  boundary_ids.insert(3);
107 
108  // Assign which variables must satisfy the boundary constraints
109  std::vector<unsigned int> vars;
110  vars.push_back(u_var);
111 
112  // Add the boundary object to the eigensystem
113  libMesh::ZeroFunction<> zf;
114  libMesh::DirichletBoundary dirichlet_bc(boundary_ids, vars, &zf);
115  eigen_system.get_dof_map().add_dirichlet_boundary(dirichlet_bc);
116 
117  // Initialize the data structures for the equation system.
118  es->init();
119 
120  // Here we obtain the ids of the Dirichlet dofs
121  std::set<unsigned int> global_dirichlet_dof_ids;
122  libMesh::DofConstraints::const_iterator it =
123  eigen_system.get_dof_map().constraint_rows_begin();
124 
125  for ( ; it != eigen_system.get_dof_map().constraint_rows_end(); ++it) {
126  global_dirichlet_dof_ids.insert(it->first);
127  }
128 
129  eigen_system.initialize_condensed_dofs(global_dirichlet_dof_ids);
130 
131  // Solve the system "Eigensystem".
132  eigen_system.solve();
133 
134  // Get the number of converged eigen pairs.
135  this->nconv = eigen_system.get_n_converged();
136 }
137 
138 LibMeshNegativeLaplacianOperator::~LibMeshNegativeLaplacianOperator()
139 {
140 }
141 
142 void LibMeshNegativeLaplacianOperator::assemble()
143 {
144 #ifdef LIBMESH_HAVE_SLEPC
145 
146  boost::shared_ptr<libMesh::EquationSystems> es(this->equation_systems);
147 
148  // Get a constant reference to the mesh object.
149  const libMesh::MeshBase& mesh = es->get_mesh();
150 
151  // The dimension that we are running.
152  const unsigned int dim = mesh.mesh_dimension();
153 
154  // Get a reference to our system.
155  libMesh::EigenSystem & eigen_system = es->get_system<libMesh::EigenSystem>("Eigensystem");
156 
157  // Get a constant reference to the Finite Element type
158  // for the first (and only) variable in the system.
159  libMesh::FEType fe_type = eigen_system.get_dof_map().variable_type(0);
160 
161  // A reference to the two system matrices
162  libMesh::SparseMatrix<libMesh::Number>& matrix_A = *eigen_system.matrix_A;
163  libMesh::SparseMatrix<libMesh::Number>& matrix_B = *eigen_system.matrix_B;
164 
165  // Build a Finite Element object of the specified type. Since the
166  // \p FEBase::build() member dynamically creates memory we will
167  // store the object as an \p AutoPtr<FEBase>. This can be thought
168  // of as a pointer that will clean up after itself.
169  libMesh::AutoPtr<libMesh::FEBase> fe (libMesh::FEBase::build(dim, fe_type));
170 
171  // A Gauss quadrature rule for numerical integration.
172  // Use the default quadrature order.
173  libMesh::QGauss qrule(dim, fe_type.default_quadrature_order());
174 
175  // Tell the finite element object to use our quadrature rule.
176  fe->attach_quadrature_rule (&qrule);
177 
178  // The element Jacobian * quadrature weight at each integration point.
179  const std::vector<libMesh::Real>& JxW = fe->get_JxW();
180 
181  // The element shape functions evaluated at the quadrature points.
182  const std::vector<std::vector<libMesh::Real> >& phi = fe->get_phi();
183  const std::vector<std::vector<libMesh::RealGradient> >& dphi = fe->get_dphi();
184 
185  // A reference to the \p DofMap object for this system. The \p DofMap
186  // object handles the index translation from node and element numbers
187  // to degree of freedom numbers.
188  const libMesh::DofMap& dof_map = eigen_system.get_dof_map();
189 
190  // The element stiffness matrix.
191  libMesh::DenseMatrix<libMesh::Number> Me;
192  libMesh::DenseMatrix<libMesh::Number> Ke;
193 
194  // This vector will hold the degree of freedom indices for
195  // the element. These define where in the global system
196  // the element degrees of freedom get mapped.
197  std::vector<libMesh::dof_id_type> dof_indices;
198 
199  // Now we will loop over all the elements in the mesh that
200  // live on the local processor. We will compute the element
201  // matrix and right-hand-side contribution. In case users
202  // later modify this program to include refinement, we will
203  // be safe and will only consider the active elements;
204  // hence we use a variant of the \p active_elem_iterator.
205  libMesh::MeshBase::const_element_iterator el = mesh.active_local_elements_begin();
206  const libMesh::MeshBase::const_element_iterator end_el =
207  mesh.active_local_elements_end();
208 
209  for ( ; el != end_el; ++el) {
210  // Store a pointer to the element we are currently
211  // working on. This allows for nicer syntax later.
212  const libMesh::Elem* elem = *el;
213 
214  // Get the degree of freedom indices for the
215  // current element. These define where in the global
216  // matrix and right-hand-side this element will
217  // contribute to.
218  dof_map.dof_indices(elem, dof_indices);
219 
220  // Compute the element-specific data for the current
221  // element. This involves computing the location of the
222  // quadrature points (q_point) and the shape functions
223  // (phi, dphi) for the current element.
224  fe->reinit(elem);
225 
226  // Zero the element matrices before
227  // summing them. We use the resize member here because
228  // the number of degrees of freedom might have changed from
229  // the last element. Note that this will be the case if the
230  // element type is different (i.e. the last element was a
231  // triangle, now we are on a quadrilateral).
232  Ke.resize(dof_indices.size(), dof_indices.size());
233  Me.resize(dof_indices.size(), dof_indices.size());
234 
235  // Now loop over the quadrature points. This handles
236  // the numeric integration.
237  //
238  // We will build the element matrix. This involves
239  // a double loop to integrate the test funcions (i) against
240  // the trial functions (j).
241  for (unsigned int qp=0; qp < qrule.n_points(); qp++)
242  for (unsigned int i=0; i < phi.size(); i++)
243  for (unsigned int j=0; j < phi.size(); j++)
244  {
245  Me(i, j) += JxW[qp] * phi[i][qp] * phi[j][qp];
246  Ke(i, j) += JxW[qp] * (dphi[i][qp] * dphi[j][qp]);
247  }
248 
249  // On an unrefined mesh, constrain_element_matrix does
250  // nothing. If this assembly function is ever repurposed to
251  // run on a refined mesh, getting the hanging node constraints
252  // right will be important. Note that, even with
253  // asymmetric_constraint_rows = false, the constrained dof
254  // diagonals still exist in the matrix, with diagonal entries
255  // that are there to ensure non-singular matrices for linear
256  // solves but which would generate positive non-physical
257  // eigenvalues for eigensolves.
258  // dof_map.constrain_element_matrix(Ke, dof_indices, false);
259  // dof_map.constrain_element_matrix(Me, dof_indices, false);
260 
261  // Finally, simply add the element contribution to the
262  // overall matrices A and B.
263  matrix_A.add_matrix (Ke, dof_indices);
264  matrix_B.add_matrix (Me, dof_indices);
265 
266  } // end of element loop
267 
268 #endif // LIBMESH_HAVE_SLEPC
269 }
270 
271 void LibMeshNegativeLaplacianOperator::print_info() const
272 {
273  // Prints information about the system to the screen.
274  this->equation_systems->get_mesh().print_info();
275  this->equation_systems->print_info();
276 }
277 
278 } // End namespace QUESO
279 
280 #endif // QUESO_HAVE_LIBMESH
int dim
Definition: ann2fig.cpp:81
that you receive source code or can get it if you want it
Definition: License.txt:41

Generated on Thu Jun 11 2015 13:52:31 for queso-0.53.0 by  doxygen 1.8.5