MAST
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
homogenized_density_function_base.cpp
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 // MAST includes
23 #include "level_set/filter_base.h"
24 
25 
27 HomogenizedDensityFunctionBase(const std::string& nm):
28 MAST::FieldFunction<Real> (nm),
29 _level_set_sys (nullptr),
30 _analysis_mesh (nullptr),
31 _level_set (nullptr),
32 _filter (nullptr) {
33 
34 }
35 
36 
39 
40 }
41 
42 
43 void
45  libMesh::MeshBase& analysis_mesh,
46  MAST::FieldFunction<Real>& level_set,
47  MAST::FilterBase& filter) {
48 
49  libmesh_assert(!_level_set_sys);
50 
51  _level_set_sys = &level_set_sys;
52  _analysis_mesh = &analysis_mesh;
53  _level_set = &level_set;
54  _filter = &filter;
55  _sub_point_locator.reset(_analysis_mesh->sub_point_locator().release());
56 }
57 
58 
59 void
61 operator() (const libMesh::Point& p, const Real t, Real& v) const {
62 
63  libmesh_assert(_analysis_mesh);
64  libmesh_assert(!_elem_volume_fraction.empty());
65 
66  // identify the element in the analysis mesh that this point belongs to
67  // and return the value of the homogenized function
68  const libMesh::Elem* e = (*_sub_point_locator)(p);
69 
70  libmesh_assert(e);
71 
72  std::map<const libMesh::Elem*, Real>::const_iterator
73  it = _elem_volume_fraction.find(e);
74 
75  libmesh_assert(it != _elem_volume_fraction.end());
76 
77  v = it->second;
78 }
79 
80 
81 
82 void
85  const libMesh::Point& p, const Real t, Real& v) const {
86 
87  libmesh_assert(_analysis_mesh);
88  libmesh_assert(!_elem_volume_fraction_sensitivity.empty());
89 
90  // identify the element in the analysis mesh that this point belongs to
91  // and return the value of the homogenized function
92  const libMesh::Elem* e = (*_sub_point_locator)(p);
93 
94  libmesh_assert(e);
95 
96  v = this->get_elem_volume_fraction_sensitivity(f, *e);
97 }
98 
99 
100 Real
102 get_elem_volume_fraction(const libMesh::Elem& e) const {
103 
104  std::map<const libMesh::Elem*, Real>::const_iterator
105  it = _elem_volume_fraction.find(&e);
106 
107  libmesh_assert(it != _elem_volume_fraction.end());
108 
109  return it->second;
110 }
111 
112 
113 
114 Real
117  const libMesh::Elem& e) const {
118 
119  std::map<const MAST::FunctionBase*, std::map<const libMesh::Elem*, Real>>::const_iterator
120  it_f = _elem_volume_fraction_sensitivity.find(&f);
121 
122  libmesh_assert (it_f != _elem_volume_fraction_sensitivity.end());
123 
124  std::map<const libMesh::Elem*, Real>::const_iterator
125  it_e = it_f->second.find(&e);
126 
127  libmesh_assert (it_e != it_f->second.end());
128 
129  return it_e->second;
130 }
131 
132 
133 const std::map<const libMesh::Elem*, Real>*
136 {
137  std::map<const MAST::FunctionBase*, std::map<const libMesh::Elem*, Real>>::const_iterator
138  it = _elem_volume_fraction_sensitivity.find(&f);
139 
140  if (it != _elem_volume_fraction_sensitivity.end())
141  return &(it->second);
142  else
143  return nullptr;
144 }
const std::map< const libMesh::Elem *, Real > * get_elem_volume_fraction_sensitivity_map(const MAST::FunctionBase &f) const
libMesh::Real Real
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.
Real get_elem_volume_fraction(const libMesh::Elem &e) const
Creates a geometric filter for the level-set design variables.
Definition: filter_base.h:39
virtual void init(MAST::SystemInitialization &level_set_sys, libMesh::MeshBase &analysis_mesh, MAST::FieldFunction< Real > &level_set, MAST::FilterBase &filter)
std::map< const MAST::FunctionBase *, std::map< const libMesh::Elem *, Real > > _elem_volume_fraction_sensitivity
virtual void derivative(const MAST::FunctionBase &f, const libMesh::Point &p, const Real t, Real &v) const
calculates the value of the derivative of function with respect to the function f at the specified po...
std::map< const libMesh::Elem *, Real > _elem_volume_fraction
std::unique_ptr< libMesh::PointLocatorBase > _sub_point_locator
Real get_elem_volume_fraction_sensitivity(const MAST::FunctionBase &f, const libMesh::Elem &e) const
virtual void operator()(const libMesh::Point &p, const Real t, Real &v) const
calculates the value of the function at the specified point, p, and time, t, and returns it in v...