queso-0.53.0
LibMeshOperatorBase.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 <iostream>
30 #include <fstream>
31 #include <vector>
32 #include <queso/FunctionOperatorBuilder.h>
33 #include <queso/LibMeshFunction.h>
34 #include <queso/LibMeshOperatorBase.h>
35 
36 #include <libmesh/libmesh.h>
37 #include <libmesh/mesh.h>
38 #include <libmesh/mesh_generation.h>
39 #include <libmesh/equation_systems.h>
40 #include <libmesh/explicit_system.h>
41 #include <libmesh/condensed_eigen_system.h>
42 #include <libmesh/numeric_vector.h>
43 #include <libmesh/exodusII_io.h>
44 #include <libmesh/parallel.h>
45 #include <libmesh/parallel_implementation.h>
46 
47 namespace QUESO {
48 
49 LibMeshOperatorBase::LibMeshOperatorBase(
50  const FunctionOperatorBuilder & builder, libMesh::MeshBase & m)
51  : OperatorBase(),
52  equation_systems(new libMesh::EquationSystems(m)),
53  builder(builder)
54 {
55 #ifndef LIBMESH_HAVE_SLEPC
56  if (m.processor_id() == 0)
57  std::cerr << "ERROR: This example requires libMesh to be\n"
58  << "compiled with SLEPc eigen solvers support!"
59  << std::endl;
60 #else
61 
62 #ifdef LIBMESH_DEFAULT_SINGLE_PRECISION
63  // SLEPc currently gives us a nasty crash with Real==float
64  libmesh_example_assert(false, "--disable-singleprecision");
65 #endif
66 
67 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
68  // SLEPc currently gives us an "inner product not well defined" with
69  // Number==complex
70  libmesh_example_assert(false, "--disable-complex");
71 #endif
72 
73  // Create a CondensedEigenSystem named "Eigensystem" and (for convenience)
74  // use a reference to the system we create.
75  this->equation_systems->add_system<libMesh::CondensedEigenSystem>("Eigensystem");
76 
77 #endif // LIBMESH_HAVE_SLEPC
78 }
79 
80 LibMeshOperatorBase::~LibMeshOperatorBase()
81 {
82 }
83 
84 void LibMeshOperatorBase::save_converged_evals(const std::string & filename) const
85 {
86  unsigned int i;
87  std::ofstream evals_file(filename.c_str());
88 
89  std::pair<libMesh::Real, libMesh::Real> eval;
90  for (i = 0; i < this->nconv; i++) {
91  eval = this->equation_systems
92  ->get_system<libMesh::EigenSystem>("Eigensystem").get_eigenpair(i);
93  if (equation_systems->processor_id() == 0) {
94  evals_file << eval.first << " " << eval.second << std::endl;
95  }
96  }
97  if (this->equation_systems->processor_id() == 0) {
98  evals_file.close();
99  }
100 }
101 
102 void LibMeshOperatorBase::save_converged_evec(const std::string & filename,
103  unsigned int i) const
104 {
105  if (i < this->nconv) {
106  boost::shared_ptr<libMesh::EquationSystems> es(this->equation_systems);
107  es->get_system<libMesh::EigenSystem>("Eigensystem").get_eigenpair(i);
108  libMesh::ExodusII_IO(es->get_mesh()).write_equation_systems(filename, *es);
109  }
110  else {
111  std::cerr << "Warning: eigenpair " << i
112  << " did not converge. Not saving."
113  << std::endl;
114  }
115 }
116 
117 unsigned int LibMeshOperatorBase::get_num_converged() const {
118  return this->nconv;
119 }
120 
121 double LibMeshOperatorBase::get_eigenvalue(unsigned int i) const
122 {
123  if (i < this->nconv) {
124  std::pair<libMesh::Real, libMesh::Real> eval;
125  boost::shared_ptr<libMesh::EquationSystems> es(this->equation_systems);
126  eval = es->get_system<libMesh::EigenSystem>("Eigensystem").get_eigenpair(i);
127  return eval.first;
128  }
129  else {
130  return -1;
131  }
132 }
133 
134 double LibMeshOperatorBase::get_inverted_eigenvalue(unsigned int i) const
135 {
136  return 1.0 / this->get_eigenvalue(i);
137 }
138 
139 libMesh::EquationSystems & LibMeshOperatorBase::get_equation_systems() const
140 {
141  return *this->equation_systems;
142 }
143 
144 boost::shared_ptr<FunctionBase>
145 LibMeshOperatorBase::inverse_kl_transform(std::vector<double> & xi,
146  double alpha) const
147 {
148  unsigned int i;
149  boost::shared_ptr<libMesh::EquationSystems> es(this->equation_systems);
150  LibMeshFunction *kl = new LibMeshFunction(this->builder, es->get_mesh());
151 
152  // Make sure all procs in libmesh mpi communicator all have the same xi. No,
153  // I can't set the seed in QUESO. That would mess with the QUESO
154  // communicator.
155  this->equation_systems->comm().broadcast(xi);
156 
157  boost::shared_ptr<libMesh::EquationSystems> kl_eq_sys(kl->get_equation_systems());
158 
159  std::pair<libMesh::Real, libMesh::Real> eval;
160  for (i = 0; i < this->get_num_converged(); i++) {
161  eval = es->get_system<libMesh::EigenSystem>("Eigensystem").get_eigenpair(i);
162  kl_eq_sys->get_system<libMesh::ExplicitSystem>("Function").solution->add(
163  xi[i] / pow(eval.first, alpha / 2.0),
164  *es->get_system<libMesh::EigenSystem>("Eigensystem").solution);
165  }
166 
167  boost::shared_ptr<FunctionBase> ap(kl);
168  return ap;
169 }
170 
171 } // End namespace QUESO
172 
173 #endif // QUESO_HAVE_LIBMESH

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