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-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 #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 
54 
58  MeshFieldFunction(libMesh::System& sys, const std::string& nm);
59 
60 
64  virtual ~MeshFieldFunction();
65 
66 
71  virtual void operator() (const libMesh::Point& p,
72  const Real t,
73  RealVectorX& v) const;
74 
80  virtual void gradient(const libMesh::Point& p,
81  const Real t,
82  RealMatrixX& g) const;
83 
89  virtual void perturbation (const libMesh::Point& p,
90  const Real t,
91  RealVectorX& v) const;
92 
93 
98  virtual void derivative (const MAST::FunctionBase& f,
99  const libMesh::Point& p,
100  const Real t,
101  RealVectorX& v) const;
102 
103 
109  void init(const libMesh::NumericVector<Real>& sol,
110  const libMesh::NumericVector<Real>* dsol = nullptr);
111 
112 
116  libMesh::MeshFunction& get_function() {
117 
118  libmesh_assert(_function);
119  return *_function;
120  }
121 
126  libMesh::MeshFunction& get_perturbed_function() {
127 
128  libmesh_assert(_perturbed_function);
129  return *_perturbed_function;
130  }
131 
132 
141 
142 
148 
149 
153  void clear();
154 
155  protected:
156 
162 
163 
168 
172  libMesh::System* _sys;
173 
177  libMesh::NumericVector<Real> *_sol, *_dsol;
178 
182  libMesh::MeshFunction *_function, *_perturbed_function;
183  };
184 }
185 
186 #endif // __mast__mesh_field_function__
187 
MeshFieldFunction(MAST::SystemInitialization &sys, const std::string &nm)
constructor
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::NumericVector< Real > * _sol
current solution that is going to be interpolated
libMesh::NumericVector< Real > * _dsol
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.
libMesh::MeshFunction * _function
the MeshFunction object that performs the interpolation
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...
This provides a wrapper FieldFunction compatible class that interpolates the solution using libMesh&#39;s...
libMesh::MeshFunction * _perturbed_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 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...
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
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()
bool _use_qp_sol
flag is set to true when the quadrature point solution is provided by an element
void init(const libMesh::NumericVector< Real > &sol, const libMesh::NumericVector< Real > *dsol=nullptr)
initializes the data structures to perform the interpolation function of sol.
libMesh::System * _sys
current system for which solution is to be interpolated
virtual ~MeshFieldFunction()
destructor