MAST
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
level_set_elem_base.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 
21 #ifndef __mast__level_set_elem_base__
22 #define __mast__level_set_elem_base__
23 
24 // MAST includes
25 #include "base/elem_base.h"
26 
27 namespace MAST {
28 
29  // Forward declerations
30  class ElementPropertyCardBase;
31  class BoundaryConditionBase;
32  class FEMOperatorMatrix;
33  template <typename ValType> class FieldFunction;
34 
35 
37  public MAST::ElementBase
38  {
39  public:
44  const MAST::GeomElem& elem);
45 
46 
47  virtual ~LevelSetElementBase();
48 
49 
50  void
52 
53  libmesh_assert(!_phi_vel);
54  _phi_vel = &vel;
55  }
56 
62  void set_propagation_mode(bool f) {
63  _if_propagation = f;
64  }
65 
66 
73 
74 
78  virtual bool
79  internal_residual (bool request_jacobian,
80  RealVectorX& f,
81  RealMatrixX& jac);
82 
83 
87  virtual bool
88  velocity_residual (bool request_jacobian,
89  RealVectorX& f,
90  RealMatrixX& jac_xdot,
91  RealMatrixX& jac);
92 
96  bool
97  side_external_residual (bool request_jacobian,
98  RealVectorX& f,
99  RealMatrixX& jac,
100  std::multimap<libMesh::boundary_id_type, MAST::BoundaryConditionBase*>& bc);
101 
105  bool
106  volume_external_residual (bool request_jacobian,
107  RealVectorX& f,
108  RealMatrixX& jac,
109  std::multimap<libMesh::subdomain_id_type, MAST::BoundaryConditionBase*>& bc);
110 
114  virtual bool
115  internal_residual_sensitivity (bool request_jacobian,
116  RealVectorX& f,
117  RealMatrixX& jac);
121  virtual bool
122  velocity_residual_sensitivity (bool request_jacobian,
123  RealVectorX& f,
124  RealMatrixX& jac);
125 
129  bool
130  side_external_residual_sensitivity (bool request_jacobian,
131  RealVectorX& f,
132  RealMatrixX& jac,
133  std::multimap<libMesh::boundary_id_type, MAST::BoundaryConditionBase*>& bc);
134 
138  bool
139  volume_external_residual_sensitivity (bool request_jacobian,
140  RealVectorX& f,
141  RealMatrixX& jac,
142  std::multimap<libMesh::subdomain_id_type, MAST::BoundaryConditionBase*>& bc);
143 
148  Real
149  volume();
150 
151 
157  Real
159 
160 
165  Real
167 
168 
182  Real
183  perimeter(Real delta=0.1);
184 
190  Real
191  perimeter_sensitivity(Real delta=0.1);
192 
193 
199  Real
200  volume_boundary_velocity_on_side(unsigned int s);
201 
202  protected:
203 
207  void _velocity_and_source(const unsigned int qp,
208  const libMesh::Point& p,
209  const Real t,
210  const MAST::FEMOperatorMatrix& Bmat,
211  const std::vector<MAST::FEMOperatorMatrix>& dBmat,
212  RealVectorX& vel,
213  Real& source);
214 
218  void _tau(const MAST::FEBase& fe,
219  unsigned int qp,
220  const MAST::FEMOperatorMatrix& Bmat,
221  const std::vector<MAST::FEMOperatorMatrix>& dBmat,
222  const RealVectorX& vel,
223  RealMatrixX& tau);
224 
225 
226 
227 
228  void
229  _dc_operator(const MAST::FEBase& fe,
230  const unsigned int qp,
231  const std::vector<MAST::FEMOperatorMatrix>& dB_mat,
232  const RealVectorX& vel,
233  Real& dc);
234 
235  void
236  _calculate_dxidX (const MAST::FEBase& fe,
237  const unsigned int qp,
238  RealMatrixX& dxi_dX);
239 
252  void _initialize_fem_operators(const unsigned int qp,
253  const MAST::FEBase& fe,
255  std::vector<MAST::FEMOperatorMatrix>& dBmat);
256 
257 
262 
268 
273  };
274 
275 }
276 
277 
278 #endif // __mast__level_set_elem_base__
279 
virtual bool internal_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac)
internal force contribution to system residual
void _calculate_dxidX(const MAST::FEBase &fe, const unsigned int qp, RealMatrixX &dxi_dX)
const MAST::FieldFunction< Real > * _phi_vel
element property
bool volume_external_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, std::multimap< libMesh::subdomain_id_type, MAST::BoundaryConditionBase *> &bc)
volume external force contribution to system residual
const RealVectorX & sol(bool if_sens=false) const
Definition: elem_base.cpp:51
LevelSetElementBase(MAST::SystemInitialization &sys, const MAST::GeomElem &elem)
Constructor.
virtual bool velocity_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xdot, RealMatrixX &jac)
inertial force contribution to system residual
void _tau(const MAST::FEBase &fe, unsigned int qp, const MAST::FEMOperatorMatrix &Bmat, const std::vector< MAST::FEMOperatorMatrix > &dBmat, const RealVectorX &vel, RealMatrixX &tau)
initializes the tau operator
libMesh::Real Real
RealVectorX _ref_sol
reference solution for reinitialization of the level set
Real volume_boundary_velocity_on_side(unsigned int s)
void set_reference_solution_for_initialization(const RealVectorX &sol)
For reinitialization to , the solution before initialization is used to calculate the source and velo...
bool volume_external_residual_sensitivity(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, std::multimap< libMesh::subdomain_id_type, MAST::BoundaryConditionBase *> &bc)
sensitivity of the volume external force contribution to system residual
void set_propagation_mode(bool f)
This can operate in one of two modes: propagation of level set given Vn, or reinitialization of level...
Matrix< Real, Dynamic, Dynamic > RealMatrixX
virtual bool velocity_residual_sensitivity(bool request_jacobian, RealVectorX &f, RealMatrixX &jac)
sensitivity of the damping force contribution to system residual
Matrix< Real, Dynamic, 1 > RealVectorX
void set_velocity_function(const MAST::FieldFunction< Real > &vel)
Real perimeter(Real delta=0.1)
Approximates the integral of the Dirac delta function to approximate the perimeter.
This class acts as a wrapper around libMesh::Elem for the purpose of providing a uniform interface fo...
Definition: geom_elem.h:59
const MAST::GeomElem & elem() const
Definition: elem_base.h:109
Real homogenized_volume_fraction(Real delta=0.1)
Approximates the volume fraction of the element based on integration of approximated Heaviside functi...
bool side_external_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, std::multimap< libMesh::boundary_id_type, MAST::BoundaryConditionBase *> &bc)
side external force contribution to system residual
bool _if_propagation
this can operate in one of two modes: propagation of level set given Vn, or reinitialization of level...
void _initialize_fem_operators(const unsigned int qp, const MAST::FEBase &fe, MAST::FEMOperatorMatrix &Bmat, std::vector< MAST::FEMOperatorMatrix > &dBmat)
When mass = false, initializes the FEM operator matrix to the shape functions as ...
Real homogenized_volume_fraction_sensitivity(Real delta=0.1)
Sensitivity of the homogenized volume fraction for the specified level set value and its sensitivity...
bool side_external_residual_sensitivity(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, std::multimap< libMesh::boundary_id_type, MAST::BoundaryConditionBase *> &bc)
sensitivity of the side external force contribution to system residual
void _velocity_and_source(const unsigned int qp, const libMesh::Point &p, const Real t, const MAST::FEMOperatorMatrix &Bmat, const std::vector< MAST::FEMOperatorMatrix > &dBmat, RealVectorX &vel, Real &source)
calculates the velocity at the quadrature point
Real perimeter_sensitivity(Real delta=0.1)
computes the partial derivative of the integral of the Dirac delta function using the solution and se...
virtual bool internal_residual_sensitivity(bool request_jacobian, RealVectorX &f, RealMatrixX &jac)
sensitivity of the internal force contribution to system residual
void _dc_operator(const MAST::FEBase &fe, const unsigned int qp, const std::vector< MAST::FEMOperatorMatrix > &dB_mat, const RealVectorX &vel, Real &dc)
This is the base class for elements that implement calculation of finite element quantities over the ...
Definition: elem_base.h:72