MAST
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
heaviside_elem_homogenization_function.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"
26 #include "base/nonlinear_system.h"
27 #include "mesh/geom_elem.h"
28 
29 // libMesh includes
30 #include "libmesh/elem.h"
31 
32 
36 _width (0.1) {
37 
38 }
39 
40 
41 
44 
45 }
46 
47 
48 void
51 
52  libmesh_assert(_analysis_mesh);
53  libmesh_assert(_elem_volume_fraction.empty());
54 
55  // make sure that the system uses Lagrange shape functions for
56  // the level set.
57  libmesh_assert(_level_set_sys->system().variable_type(0).family == libMesh::LAGRANGE);
58 
59  // Next, compute the value of level-set at the nodes of the analysis mesh
60  // and obtain the value of level-set from the mesh function.
61  // These values will be provided to the elmeent for computation of
62  // the homogenized volume fraction.
63 
65  phi;
66 
67  unsigned int
68  n_nodes = 0;
69 
70  Real
71  v = 0.;
72 
73  libMesh::MeshBase::const_element_iterator
74  e_it = _analysis_mesh->active_local_elements_begin(),
75  e_end = _analysis_mesh->active_local_elements_end();
76 
77  for ( ; e_it != e_end; e_it++) {
78 
79  const libMesh::Elem* e = *e_it;
80 
81  n_nodes = e->n_nodes();
82  phi.setZero(n_nodes);
83 
84  // get the nodal values for this element
85  for (unsigned int i=0; i<n_nodes; i++) {
86  (*_level_set)(*e->node_ptr(i), 0., v);
87  phi(i) = v;
88  }
89 
90  // compute the integral of Heaviside function to approximate the
91  // volume fraction
92  MAST::GeomElem geom_elem;
93  geom_elem.init(*e, *_level_set_sys);
94  MAST::LevelSetElementBase level_set_elem(*_level_set_sys, geom_elem);
95  level_set_elem.set_solution(phi);
96  v = level_set_elem.homogenized_volume_fraction(_width);
97 
98  // store the value for this elmeent
99  _elem_volume_fraction[e] = v;
100  }
101 }
102 
103 
104 
105 void
108 
109  libmesh_assert(_analysis_mesh);
110  libmesh_assert(_elem_volume_fraction_sensitivity[&f].empty());
111  libmesh_assert(f.is_topology_parameter());
112 
113  // make sure that the system uses Lagrange shape functions for
114  // the level set.
115  libmesh_assert(_level_set_sys->system().variable_type(0).family == libMesh::LAGRANGE);
116 
117  // Next, compute the value of level-set at the nodes of the analysis mesh
118  // and obtain the value of level-set from the mesh function.
119  // These values will be provided to the elmeent for computation of
120  // the homogenized volume fraction.
121 
123  phi,
124  phi_sens;
125 
126  unsigned int
127  n_nodes = 0;
128 
129  Real
130  v = 0.;
131 
133  &p_ls = dynamic_cast<const MAST::LevelSetParameter&>(f);
134 
135  libMesh::MeshBase::const_element_iterator
136  e_it = _analysis_mesh->active_local_elements_begin(),
137  e_end = _analysis_mesh->active_local_elements_end();
138 
139  for ( ; e_it != e_end; e_it++) {
140 
141  const libMesh::Elem* e = *e_it;
142 
144 
145  n_nodes = e->n_nodes();
146  phi.setZero(n_nodes);
147  phi_sens.setZero(n_nodes);
148 
149  // get the nodal values for this element
150  for (unsigned int i=0; i<n_nodes; i++) {
151 
152  (*_level_set)(*e->node_ptr(i), 0., v);
153  phi(i) = v;
154  _level_set->derivative(f, *e->node_ptr(i), 0., v);
155  phi_sens(i) = v;
156  }
157 
158  // compute the integral of Heaviside function to approximate the
159  // volume fraction
160  MAST::GeomElem geom_elem;
161  geom_elem.init(*e, *_level_set_sys);
162  MAST::LevelSetElementBase level_set_elem(*_level_set_sys, geom_elem);
163  level_set_elem.set_solution(phi);
164  level_set_elem.set_solution(phi_sens, true);
166 
167  // store the value for this elmeent
169  }
170  else
172  }
173 }
174 
175 
MAST::NonlinearSystem & system()
This defines a parameter that is a level set function and stores a pointer to the node in the level-s...
virtual void derivative(const MAST::FunctionBase &f, ValType &v) const
calculates the value of the function derivative and returns it in v.
libMesh::Real Real
virtual bool is_topology_parameter() const
Definition: function_base.h:97
virtual void initialize_element_volume_fraction_sensitivity(const MAST::FunctionBase &f)
Matrix< Real, Dynamic, 1 > RealVectorX
std::map< const MAST::FunctionBase *, std::map< const libMesh::Elem *, Real > > _elem_volume_fraction_sensitivity
std::map< const libMesh::Elem *, Real > _elem_volume_fraction
This class acts as a wrapper around libMesh::Elem for the purpose of providing a uniform interface fo...
Definition: geom_elem.h:59
Real homogenized_volume_fraction(Real delta=0.1)
Approximates the volume fraction of the element based on integration of approximated Heaviside functi...
virtual void set_solution(const RealVectorX &vec, bool if_sens=false)
stores vec as solution for element level calculations, or its sensitivity if if_sens is true...
Definition: elem_base.cpp:60
virtual void init(const libMesh::Elem &elem, const MAST::SystemInitialization &sys_init)
initialize the object for the specified reference elem.
Definition: geom_elem.cpp:134
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 if_elem_in_domain_of_influence(const libMesh::Elem &elem, const libMesh::Node &level_set_node) const
function identifies if the given element is within the domain of influence of this specified level se...
const libMesh::Node * level_set_node() const