MAST
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
mesh_field_function.h
Go to the documentation of this file.
1 /*
2  * MAST: Multidisciplinary-design Adaptation and Sensitivity Toolkit
3  * Copyright (C) 2013-2020 Manav Bhatia and MAST authors
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 #ifndef __mast__mesh_field_function__
21 #define __mast__mesh_field_function__
22 
23 // MAST includes
25 
26 
27 // libMesh includes
28 #include "libmesh/numeric_vector.h"
29 #include "libmesh/mesh_function.h"
30 #include "libmesh/system.h"
31 
32 
33 
34 namespace MAST {
35 
36  // Forward declerations
37  class SystemInitialization;
38 
39 
45  public MAST::FieldFunction<RealVectorX> {
46 
47  public:
52  const std::string& nm,
53  libMesh::ParallelType p_type);
54 
55 
59  MeshFieldFunction(libMesh::System& sys,
60  const std::string& nm,
61  libMesh::ParallelType p_type);
62 
63 
67  virtual ~MeshFieldFunction();
68 
69 
74  virtual void operator() (const libMesh::Point& p,
75  const Real t,
76  RealVectorX& v) const;
77 
83  virtual void gradient(const libMesh::Point& p,
84  const Real t,
85  RealMatrixX& g) const;
86 
92  virtual void perturbation (const libMesh::Point& p,
93  const Real t,
94  RealVectorX& v) const;
95 
96 
101  virtual void perturbation_gradient (const libMesh::Point& p,
102  const Real t,
103  RealMatrixX& v) const;
104 
105 
110  virtual void derivative (const MAST::FunctionBase& f,
111  const libMesh::Point& p,
112  const Real t,
113  RealVectorX& v) const;
114 
119  virtual void derivative_gradient (const MAST::FunctionBase& f,
120  const libMesh::Point& p,
121  const Real t,
122  RealMatrixX& v) const;
123 
132  void init(const libMesh::NumericVector<Real>& sol, bool reuse_vector);
133 
134 
138  void init_sens(const MAST::FunctionBase& f,
139  const libMesh::NumericVector<Real>& dsol,
140  bool reuse_vector);
141 
142 
146  libMesh::MeshFunction& get_function() {
147 
148  libmesh_assert(_function);
149  return *_function->_func;
150  }
151 
156  libMesh::MeshFunction& get_perturbed_function() {
157 
158  libmesh_assert(_perturbed_function);
159  return *_perturbed_function->_func;
160  }
161 
162 
171 
172 
178 
179 
183  void clear();
184 
185  protected:
186 
187  struct SolFunc {
189  _sol (nullptr),
190  _cloned_sol (nullptr),
191  _func (nullptr) {}
193  if (_cloned_sol) delete _cloned_sol;
194  delete _func;
195  }
196 
198  const libMesh::NumericVector<Real>* _sol;
199  libMesh::NumericVector<Real>* _cloned_sol;
200  libMesh::MeshFunction* _func;
201  };
202 
203 
204  void _init_sol_func(bool reuse_sol,
205  const libMesh::NumericVector<Real>& sol,
207 
213 
217  libMesh::ParallelType _p_type;
218 
223 
227  libMesh::System* _sys;
228 
233 
234 
239 
243  std::map<const MAST::FunctionBase*, MAST::MeshFieldFunction::SolFunc*> _function_sens;
244  };
245 }
246 
247 #endif // __mast__mesh_field_function__
248 
virtual void perturbation_gradient(const libMesh::Point &p, const Real t, RealMatrixX &v) const
calculates the value of the function at the specified point, p, and time, t, and returns it in v...
virtual void gradient(const libMesh::Point &p, const Real t, RealMatrixX &g) const
calculates the gradient of value of the function at the specified point, p, and time, t, and returns it in g.
virtual void perturbation(const libMesh::Point &p, const Real t, RealVectorX &v) const
calculates the value of perturbation in the function at the specified point, p, and time...
void init(const libMesh::NumericVector< Real > &sol, bool reuse_vector)
initializes the data structures to perform the interpolation function of sol.
void _init_sol_func(bool reuse_sol, const libMesh::NumericVector< Real > &sol, MAST::MeshFieldFunction::SolFunc &sol_func)
SolFunc * _function
current solution that is going to be interpolated
MeshFieldFunction(MAST::SystemInitialization &sys, const std::string &nm, libMesh::ParallelType p_type)
constructor
SolFunc * _perturbed_function
current perturbation solution that is going to be interpolated
This provides a wrapper FieldFunction compatible class that interpolates the solution using libMesh&#39;s...
virtual void derivative(const MAST::FunctionBase &f, const libMesh::Point &p, const Real t, RealVectorX &v) const
calculates the value of the function at the specified point, p, and time, t, and returns it in v...
libMesh::ParallelType _p_type
type of parallel vector required for this mesh function.
virtual void clear_element_quadrature_point_solution()
clears the quadrature point solution provided by the corresponding set method above.
libMesh::Real Real
virtual void derivative_gradient(const MAST::FunctionBase &f, const libMesh::Point &p, const Real t, RealMatrixX &v) const
calculates the value of the function at the specified point, p, and time, t, and returns it in v...
libMesh::NumericVector< Real > * _cloned_sol
RealVectorX _qp_sol
quadrature point solution of the element
virtual void set_element_quadrature_point_solution(RealVectorX &sol)
When a mesh field function is attached to an assembly routine during system assembly, then the current solution can be provided by the element quadrature point update.
libMesh::MeshFunction & get_perturbed_function()
Matrix< Real, Dynamic, Dynamic > RealMatrixX
virtual void operator()(const libMesh::Point &p, const Real t, RealVectorX &v) const
calculates the value of the function at the specified point, p, and time, t, and returns it in v...
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.
void clear()
clear the solution
Matrix< Real, Dynamic, 1 > RealVectorX
libMesh::MeshFunction & get_function()
const libMesh::NumericVector< Real > * _sol
bool _use_qp_sol
flag is set to true when the quadrature point solution is provided by an element
void init_sens(const MAST::FunctionBase &f, const libMesh::NumericVector< Real > &dsol, bool reuse_vector)
initializes the the data structures for computation of sensitivity for the specified function...
libMesh::System * _sys
current system for which solution is to be interpolated
virtual ~MeshFieldFunction()
destructor
std::map< const MAST::FunctionBase *, MAST::MeshFieldFunction::SolFunc * > _function_sens
solution sensitivity for specified value