MAST
complex_mesh_field_function.cpp
Go to the documentation of this file.
1 /*
2  * MAST: Multidisciplinary-design Adaptation and Sensitivity Toolkit
3  * Copyright (C) 2013-2019 Manav Bhatia
4  *
5  * This library is free software; you can redistribute it and/or
6  * modify it under the terms of the GNU Lesser General Public
7  * License as published by the Free Software Foundation; either
8  * version 2.1 of the License, or (at your option) any later version.
9  *
10  * This library is distributed in the hope that it will be useful,
11  * but WITHOUT ANY WARRANTY; without even the implied warranty of
12  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
13  * Lesser General Public License for more details.
14  *
15  * You should have received a copy of the GNU Lesser General Public
16  * License along with this library; if not, write to the Free Software
17  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
18  */
19 
20 // MAST includes
23 #include "base/nonlinear_system.h"
24 
25 // libMesh includes
26 #include "libmesh/dof_map.h"
27 #include "libmesh/numeric_vector.h"
28 
29 
32  const std::string& nm):
34 _system(&sys),
35 _sol_re(nullptr),
36 _sol_im(nullptr),
37 _perturbed_sol_re(nullptr),
38 _perturbed_sol_im(nullptr),
39 _function_re(nullptr),
40 _function_im(nullptr),
41 _perturbed_function_re(nullptr),
42 _perturbed_function_im(nullptr)
43 { }
44 
45 
46 
47 
49 
50  this->clear();
51 }
52 
53 
54 
55 
56 
57 void
59 init(const libMesh::NumericVector<Real>& sol_re,
60  const libMesh::NumericVector<Real>& sol_im) {
61 
62  // first make sure that the object is not already initialized
63  libmesh_assert(!_function_re);
64 
66 
67  // next, clone this solution and localize to the sendlist
68  _sol_re = libMesh::NumericVector<Real>::build(system.comm()).release();
69  _sol_im = libMesh::NumericVector<Real>::build(system.comm()).release();
70 
71  const std::vector<libMesh::dof_id_type>& send_list =
72  system.get_dof_map().get_send_list();
73 
74  // initialize and then localize the vector with the provided solution
75  /*_sol_re->init(system.n_dofs(),
76  system.n_local_dofs(),
77  send_list,
78  false,
79  libMesh::GHOSTED);
80  sol_re.localize(*_sol_re, send_list);
81 
82  _sol_im->init(system.n_dofs(),
83  system.n_local_dofs(),
84  send_list,
85  false,
86  libMesh::GHOSTED);
87  sol_im.localize(*_sol_im, send_list);*/
88 
89  _sol_re->init(sol_re.size(), true, libMesh::SERIAL);
90  sol_re.localize(*_sol_re);
91  _sol_im->init(sol_im.size(), true, libMesh::SERIAL);
92  sol_im.localize(*_sol_im);
93 
94 
95  // finally, create the mesh interpolation function
96  _function_re = new libMesh::MeshFunction(system.get_equation_systems(),
97  *_sol_re,
98  system.get_dof_map(),
99  _system->vars());
100  _function_re->init();
101 
102  _function_im = new libMesh::MeshFunction(system.get_equation_systems(),
103  *_sol_im,
104  system.get_dof_map(),
105  _system->vars());
106  _function_im->init();
107 }
108 
109 
110 
111 
112 void
114 init_perturbation(const libMesh::NumericVector<Real>& sol_re,
115  const libMesh::NumericVector<Real>& sol_im) {
116 
117  // first make sure that the object is not already initialized
118  libmesh_assert(!_perturbed_function_re);
119 
120  MAST::NonlinearSystem& system = _system->system();
121 
122  // next, clone this solution and localize to the sendlist
123  _perturbed_sol_re = libMesh::NumericVector<Real>::build(system.comm()).release();
124  _perturbed_sol_im = libMesh::NumericVector<Real>::build(system.comm()).release();
125 
126  const std::vector<libMesh::dof_id_type>& send_list =
127  system.get_dof_map().get_send_list();
128 
129  // initialize and then localize the vector with the provided solution
130  /*_perturbed_sol_re->init(system.n_dofs(),
131  system.n_local_dofs(),
132  send_list,
133  false,
134  libMesh::GHOSTED);
135  sol_re.localize(*_sol_re, send_list);
136 
137  _perturbed_sol_im->init(system.n_dofs(),
138  system.n_local_dofs(),
139  send_list,
140  false,
141  libMesh::GHOSTED);
142  sol_im.localize(*_perturbed_sol_im, send_list);*/
143  _perturbed_sol_re->init(sol_re.size(), true, libMesh::SERIAL);
144  sol_re.localize(*_perturbed_sol_re);
145  _perturbed_sol_im->init(sol_im.size(), true, libMesh::SERIAL);
146  sol_im.localize(*_perturbed_sol_im);
147 
148 
149  // finally, create the mesh interpolation function
150  _perturbed_function_re = new libMesh::MeshFunction(system.get_equation_systems(),
152  system.get_dof_map(),
153  _system->vars());
154  _perturbed_function_re->init();
155 
156  _perturbed_function_im = new libMesh::MeshFunction(system.get_equation_systems(),
158  system.get_dof_map(),
159  _system->vars());
160  _perturbed_function_im->init();
161 }
162 
163 
164 
165 
166 void
168  const Real t,
169  ComplexVectorX& v) const {
170 
171  // make sure that the object was initialized
172  libmesh_assert(_function_re);
173 
174  DenseRealVector v_re, v_im;
175  (*_function_re)(p, t, v_re);
176  (*_function_im)(p, t, v_im);
177 
178  // make sure that the mesh function was able to find the element
179  // and a solution
180  libmesh_assert(v_re.size());
181 
182  // now copy this to the output vector
183  v = ComplexVectorX::Zero(v_re.size());
184  for (unsigned int i=0; i<v_re.size(); i++)
185  v(i) = std::complex<Real>(v_re(i), v_im(i));
186 }
187 
188 
189 
190 
191 
192 void
194  const Real t,
195  ComplexVectorX& v) const {
196 
197  // make sure that the object was initialized
198  libmesh_assert(_perturbed_function_re);
199 
200  DenseRealVector v_re, v_im;
201  (*_perturbed_function_re)(p, t, v_re);
202  (*_perturbed_function_im)(p, t, v_im);
203 
204  // make sure that the mesh function was able to find the element
205  // and a solution
206  libmesh_assert(v_re.size());
207 
208  // now copy this to the output vector
209  v = ComplexVectorX::Zero(v_re.size());
210  for (unsigned int i=0; i<v_re.size(); i++)
211  v(i) = std::complex<Real>(v_re(i), v_im(i));
212 }
213 
214 
215 
216 
217 
218 
219 
220 
221 void
223 
224  // if a pointer has been attached, then delete it and the
225  // associated vector, and clear the associated system
226  if (_function_re) {
227  delete _function_re;
228  delete _function_im;
229  _function_re = nullptr;
230  _function_im = nullptr;
231 
232  delete _sol_re;
233  delete _sol_im;
234  _sol_re = nullptr;
235  _sol_im = nullptr;
236 
237  }
238 
240  delete _perturbed_function_re;
241  delete _perturbed_function_im;
242  _perturbed_function_re = nullptr;
243  _perturbed_function_im = nullptr;
244 
245  delete _perturbed_sol_re;
246  delete _perturbed_sol_im;
247  _perturbed_sol_re = nullptr;
248  _perturbed_sol_im = nullptr;
249  }
250 }
251 
252 
253 
254 
libMesh::NumericVector< Real > * _sol_im
const std::vector< unsigned int > vars() const
MAST::NonlinearSystem & system()
libMesh::NumericVector< Real > * _perturbed_sol_im
This class implements a system for solution of nonlinear systems.
void init_perturbation(const libMesh::NumericVector< Real > &dsol_re, const libMesh::NumericVector< Real > &dsol_im)
libMesh::MeshFunction * _function_re
the MeshFunction object that performs the interpolation
Matrix< Complex, Dynamic, 1 > ComplexVectorX
libMesh::Real Real
void clear()
clear the solution and mesh function data structures
libMesh::NumericVector< Real > * _perturbed_sol_re
libMesh::DenseVector< Real > DenseRealVector
virtual void operator()(const libMesh::Point &p, const Real t, ComplexVectorX &v) const
calculates the value of the function at the specified point, p, and time, t, and returns it in v...
libMesh::MeshFunction * _perturbed_function_im
This creates the base class for functions that have a saptial and temporal dependence, and provide sensitivity operations with respect to the functions and parameters.
virtual void perturbation(const libMesh::Point &p, const Real t, ComplexVectorX &v) const
calculates the value of a perturbation in function at the specified point, p, and time...
void init(const libMesh::NumericVector< Real > &sol_re, const libMesh::NumericVector< Real > &sol_im)
libMesh::NumericVector< Real > * _sol_re
current solution that is going to be interpolated
libMesh::MeshFunction * _perturbed_function_re
MAST::SystemInitialization * _system
current system for which solution is to be interpolated
ComplexMeshFieldFunction(MAST::SystemInitialization &sys, const std::string &nm)
constructor