MAST
pressure_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 
21 // MAST includes
24 #include "numerics/utility.h"
27 #include "fluid/flight_condition.h"
28 #include "base/nonlinear_system.h"
29 
30 
31 // libMesh includes
32 #include "libmesh/numeric_vector.h"
33 #include "libmesh/system.h"
34 #include "libmesh/dof_map.h"
35 
36 
40 MAST::FieldFunction<Real>("pressure"),
41 _if_cp (false),
42 _ref_pressure (0.),
43 _system (sys),
44 _flt_cond (flt) {
45 
46 }
47 
48 
49 
50 
52 
53 }
54 
55 
56 
57 
58 void
60 init(const libMesh::NumericVector<Real>& steady_sol,
61  const libMesh::NumericVector<Real>* small_dist_sol) {
62 
64 
65  // first initialize the solution to the given vector
66  // steady state solution
67  _sol.reset(libMesh::NumericVector<Real>::build(sys.comm()).release());
68  _sol->init(steady_sol.size(), true, libMesh::SERIAL);
69 
70  // now localize the give solution to this objects's vector
71  steady_sol.localize(*_sol);
72 
73 
74  // if the mesh function has not been created so far, initialize it
75  _sol_function.reset(new libMesh::MeshFunction(sys.get_equation_systems(),
76  *_sol,
77  sys.get_dof_map(),
78  _system.vars()));
79  _sol_function->init();
80 
81 
82  if (small_dist_sol) {
83 
84  // solution real part
85  _dsol.reset(libMesh::NumericVector<Real>::build(sys.comm()).release());
86  _dsol->init(steady_sol.size(), true, libMesh::SERIAL);
87 
88  small_dist_sol->localize(*_dsol);
89 
90  _dsol_function.reset(new libMesh::MeshFunction(sys.get_equation_systems(),
91  *_dsol,
92  sys.get_dof_map(),
93  _system.vars()));
94  _dsol_function->init();
95  }
96  else {
97 
98  _dsol.reset();
99  _dsol_function.reset();
100  }
101 
102 }
103 
104 
105 
106 
107 
108 void
110 operator() (const libMesh::Point& p,
111  const Real time,
112  Real &press) const {
113 
114 
115  libmesh_assert(_sol_function.get()); // should be initialized before this call
116 
117  press = 0.;
118 
119 
120  // get the nonlinear and linearized solution
122  v;
123 
125  sol = RealVectorX::Zero(_system.system().n_vars());
126 
127 
128  // now the steady state function itself
129  (*_sol_function)(p, 0., v);
130  MAST::copy(sol, v);
131 
132 
134 
135  // now initialize the primitive variable contexts
136  p_sol.init(dynamic_cast<MAST::ConservativeFluidSystemInitialization&>(_system).dim(),
137  sol,
141 
142  if (_if_cp)
143  press = p_sol.c_pressure(_flt_cond.p0(), _flt_cond.q0());
144  else {
145  press = p_sol.p - _ref_pressure;
146  }
147 }
148 
149 
150 
151 
152 void
154 perturbation(const libMesh::Point& p,
155  const Real t,
156  Real &dpress) const {
157 
158 
159  libmesh_assert(_sol_function.get()); // should be initialized before this call
160  libmesh_assert(_dsol_function.get()); // should be initialized before this call
161 
162  dpress = 0.;
163 
164 
165  // get the nonlinear and linearized solution
167  v;
168 
170  sol = RealVectorX::Zero(_system.system().n_vars()),
171  dsol = RealVectorX::Zero(_system.system().n_vars());
172 
173 
174  // first copy the real and imaginary solutions
175  (*_dsol_function)(p, 0., v);
176  MAST::copy(sol, v);
177  dsol = sol;
178 
179  // now the steady state function itself
180  (*_sol_function)(p, 0., v);
181  MAST::copy(sol, v);
182 
183 
186 
187  // now initialize the primitive variable contexts
188  p_sol.init(dynamic_cast<MAST::ConservativeFluidSystemInitialization&>(_system).dim(),
189  sol,
193  delta_p_sol.init(p_sol,
194  dsol,
196 
197  if (_if_cp)
198  dpress = delta_p_sol.c_pressure(_flt_cond.q0());
199  else
200  dpress = delta_p_sol.dp;
201 }
202 
std::unique_ptr< libMesh::MeshFunction > _sol_function
mesh function that interpolates the solution
PressureFunction(MAST::SystemInitialization &sys, MAST::FlightCondition &flt)
const std::vector< unsigned int > vars() const
MAST::NonlinearSystem & system()
Class defines basic operations and calculation of the small disturbance primitive variables...
void init(const libMesh::NumericVector< Real > &steady_sol, const libMesh::NumericVector< Real > *small_dist_sol=nullptr)
initiate the mesh function for this solution
This class implements a system for solution of nonlinear systems.
Real _ref_pressure
the function will return pressure differential with respect to reference pressue defined in the fligh...
Class defines the conversion and some basic operations on primitive fluid variables used in calculati...
void init(const MAST::PrimitiveSolution &sol, const typename VectorType< ValType >::return_type &delta_sol, bool if_viscous)
MAST::FlightCondition & _flt_cond
flight condition
virtual void operator()(const libMesh::Point &p, const Real t, Real &press) const
provides the value of the pressure at the specified point and time
libMesh::Real Real
Real c_pressure(const Real p0, const Real q0) const
std::unique_ptr< libMesh::NumericVector< Real > > _sol
steady part of solution
libMesh::DenseVector< Real > DenseRealVector
virtual void perturbation(const libMesh::Point &p, const Real t, Real &dpress) const
provides the pressure perturbation.
Real q0() const
returns the flight dynamic pressure
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 copy(DenseRealMatrix &m1, const RealMatrixX &m2)
Definition: utility.h:167
Matrix< Real, Dynamic, 1 > RealVectorX
void init(const unsigned int dim, const RealVectorX &conservative_sol, const Real cp_val, const Real cv_val, bool if_viscous)
std::unique_ptr< libMesh::NumericVector< Real > > _dsol
small-disturbance solution
GasProperty gas_property
Ambient air properties.
std::unique_ptr< libMesh::MeshFunction > _dsol_function
MAST::SystemInitialization & _system
system associated with the mesh and solution vector
bool _if_cp
the function will return cp instead of pressure if this option is true.