MAST
output_assembly_elem_operations.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__output_function_base__
21 #define __mast__output_function_base__
22 
23 // C++ includes
24 #include <set>
25 #include <memory>
26 
27 // MAST includes
28 #include "base/mast_data_types.h"
30 
31 // libMesh includes
32 #include "libmesh/elem.h"
33 
34 
35 
36 namespace MAST {
37 
38  // Forward declerations
39  class FunctionBase;
40  class ElementBase;
41  class GeomElem;
42  class LevelSetIntersection;
43  template <typename ValType> class FieldFunction;
44 
56 
57  public:
58 
60 
65 
66 
75  virtual void zero_for_analysis() = 0;
76 
77 
83  virtual void zero_for_sensitivity() = 0;
84 
93  virtual Real output_for_elem() = 0;
94 
98  virtual Real output_total() = 0;
99 
108 
116  virtual Real output_sensitivity_total(const MAST::FunctionBase& p) = 0;
117 
118 
127  virtual void output_derivative_for_elem(RealVectorX& dq_dX) = 0;
128 
129 
139  virtual void evaluate() = 0;
140 
147  virtual void evaluate_sensitivity(const MAST::FunctionBase& f) = 0;
148 
155  virtual void evaluate_shape_sensitivity(const MAST::FunctionBase& f) = 0;
156 
163  virtual void
165  const MAST::FieldFunction<RealVectorX>& vel) = 0;
166 
173  void set_participating_subdomains(const std::set<libMesh::subdomain_id_type>& sids);
174 
183  void set_participating_elements(const std::set<const libMesh::Elem*>& elems);
184 
185 
190 
195  void set_participating_boundaries(const std::set<libMesh::boundary_id_type>& bids);
196 
201  const std::set<const libMesh::Elem*>&
203 
204 
209  const std::set<libMesh::subdomain_id_type>&
211 
212 
217  const std::set<libMesh::boundary_id_type>&
219 
220 
225  virtual bool if_evaluate_for_element(const MAST::GeomElem& elem) const;
226 
227 
232  virtual bool if_evaluate_for_boundary(const MAST::GeomElem& elem,
233  const unsigned int s) const;
234 
235 
236  protected:
237 
242 
247  std::set<const libMesh::Elem*> _elem_subset;
248 
249 
253  std::set<libMesh::subdomain_id_type> _sub_domain_ids;
254 
255 
259  std::set<libMesh::boundary_id_type> _bids;
260 
261  };
262 }
263 
264 #endif // __mast__output_function_base__
const std::set< libMesh::subdomain_id_type > & get_participating_subdomains()
void set_participating_subdomains(const std::set< libMesh::subdomain_id_type > &sids)
The output function can be a boundary integrated quantity, volume integrated quantity or a combinatio...
std::set< libMesh::subdomain_id_type > _sub_domain_ids
set of subdomain ids for which data will be processed.
virtual bool if_evaluate_for_element(const MAST::GeomElem &elem) const
checks to see if the object has been told about the subset of elements and if the specified element i...
void set_participating_elements(const std::set< const libMesh::Elem * > &elems)
sets the elements for which this object will evaluate and store the output data.
std::set< libMesh::boundary_id_type > _bids
set of bids for which data will be processed
virtual void evaluate_shape_sensitivity(const MAST::FunctionBase &f)=0
this evaluates all relevant shape sensitivity components on the element.
This provides the base class for definitin of element level contribution of output quantity in an ana...
virtual void zero_for_sensitivity()=0
zeroes the output quantity values stored inside this object so that assembly process can begin...
const std::set< libMesh::boundary_id_type > & get_participating_boundaries()
libMesh::Real Real
virtual void evaluate()=0
this is the abstract interface to be implemented by derived classes.
virtual ~OutputAssemblyElemOperations()
virtual destructor
virtual bool if_evaluate_for_boundary(const MAST::GeomElem &elem, const unsigned int s) const
checks to see if the specified side of the element needs evaluation of the output contribution...
Matrix< Real, Dynamic, 1 > RealVectorX
virtual Real output_sensitivity_total(const MAST::FunctionBase &p)=0
virtual void output_derivative_for_elem(RealVectorX &dq_dX)=0
returns the output quantity derivative with respect to state vector in dq_dX.
This class acts as a wrapper around libMesh::Elem for the purpose of providing a uniform interface fo...
Definition: geom_elem.h:59
virtual Real output_sensitivity_for_elem(const MAST::FunctionBase &p)=0
virtual void zero_for_analysis()=0
zeroes the output quantity values stored inside this object so that assembly process can begin...
virtual void evaluate_topology_sensitivity(const MAST::FunctionBase &f, const MAST::FieldFunction< RealVectorX > &vel)=0
this evaluates all relevant topological sensitivity components on the element.
bool _if_evaluate_on_all_elems
if true, evaluates on all elements.
const std::set< const libMesh::Elem * > & get_participating_elements() const
void set_participating_elements_to_all()
This will allow volume contribution from all elements.
void set_participating_boundaries(const std::set< libMesh::boundary_id_type > &bids)
The assembly will integration over boudnaries with ids specified in bids.
std::set< const libMesh::Elem * > _elem_subset
set of elements for which the data will be stored.
virtual void evaluate_sensitivity(const MAST::FunctionBase &f)=0
this evaluates all relevant sensitivity components on the element.