MAST3
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
heat_sink2d.hpp
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 #ifndef __mast_mesh_generation_heat_sink_2d_h__
21 #define __mast_mesh_generation_heat_sink_2d_h__
22 
23 // MAST includes
28 
29 // libMesh includes
30 #include <libmesh/system.h>
31 #include <libmesh/unstructured_mesh.h>
32 #include <libmesh/fe_type.h>
33 #include <libmesh/string_to_enum.h>
34 #include <libmesh/mesh_generation.h>
35 #include <libmesh/elem.h>
36 #include <libmesh/node.h>
37 #include <libmesh/boundary_info.h>
38 #include <libmesh/dirichlet_boundaries.h>
39 #include <libmesh/zero_function.h>
40 
41 
42 
43 namespace MAST {
44 namespace Mesh {
45 namespace Generation {
46 
47 
48 struct HeatSink2D {
49 
50  static const uint_t dim = 2;
51 
52  template <typename Context>
53  inline real_t
54  reference_volume(Context& c) {
55 
56  real_t
57  length = c.input("length", "length of domain along x-axis", 1.0),
58  height = c.input("height", "length of domain along y-axis", 1.0);
59 
60  return length * height;
61  }
62 
63 
64  inline uint_t idx(const libMesh::ElemType type,
65  const uint_t nx,
66  const uint_t i,
67  const uint_t j)
68  {
69  switch(type)
70  {
71  case libMesh::QUAD4:
72  {
73  return i + (nx+1)*j;
74  }
75 
76  case libMesh::QUAD9:
77  {
78  return i + (2*nx+1)*j;
79  }
80 
81  default:
82  Error(false, "Invalid element type");
83  }
84 
85  return libMesh::invalid_uint;
86  }
87 
88 
89 
90  inline void build_mesh(libMesh::UnstructuredMesh & mesh,
91  const uint_t nx,
92  const uint_t ny,
93  const real_t length,
94  const real_t height,
95  const real_t dirichlet_length_fraction,
96  const libMesh::ElemType type) {
97 
98  Assert0(type == libMesh::QUAD4 || type == libMesh::QUAD9,
99  "Method only implemented for Quad4/Quad9");
100 
101  // Clear the mesh and start from scratch
102  mesh.clear();
103 
104  libMesh::BoundaryInfo & boundary_info = mesh.get_boundary_info();
105 
106  mesh.set_mesh_dimension(3);
107  mesh.set_spatial_dimension(3);
108  mesh.reserve_elem(nx*ny);
109 
110  if (type == libMesh::QUAD4)
111  mesh.reserve_nodes( (nx+1)*(ny+1));
112  else if (type == libMesh::QUAD9)
113  mesh.reserve_nodes( (2*nx+1)*(2*ny+1));
114 
115  real_t
116  xmax = length,
117  ymax = height;
118 
119 
120  std::map<uint_t, libMesh::Node*> nodes;
121 
122  // Build the nodes.
123  uint_t
124  node_id = 0,
125  n = 0;
126  switch (type)
127  {
128  case libMesh::QUAD4: {
129 
130  for (uint_t j=0; j<=ny; j++)
131  for (uint_t i=0; i<=nx; i++) {
132  nodes[node_id] =
133  mesh.add_point(libMesh::Point(static_cast<real_t>(i)/static_cast<real_t>(nx)*length,
134  static_cast<real_t>(j)/static_cast<real_t>(ny)*height,
135  0.),
136  n++);
137  node_id++;
138  }
139 
140 
141  break;
142  }
143 
144  case libMesh::QUAD9: {
145 
146  for (uint_t j=0; j<=(2*ny); j++)
147  for (uint_t i=0; i<=(2*nx); i++) {
148  nodes[node_id] =
149  mesh.add_point(libMesh::Point(static_cast<real_t>(i)/static_cast<real_t>(2*nx)*length,
150  static_cast<real_t>(j)/static_cast<real_t>(2*ny)*height,
151  0.),
152  n++);
153  node_id++;
154  }
155 
156  break;
157  }
158 
159 
160  default:
161  Assert0(false, "ERROR: Unrecognized 2D element type.");
162  }
163 
164  // Build the elements.
165  uint_t
166  elem_id = 0;
167  switch (type) {
168 
169  case libMesh::QUAD4: {
170 
171  for (uint_t j=0; j<ny; j++)
172  for (uint_t i=0; i<nx; i++) {
173 
174  libMesh::Elem
175  *elem = libMesh::Elem::build(libMesh::QUAD4).release();
176  elem->set_id(elem_id++);
177  mesh.add_elem(elem);
178 
179  elem->set_node(0) = nodes[idx(type,nx,i,j) ];
180  elem->set_node(1) = nodes[idx(type,nx,i+1,j) ];
181  elem->set_node(2) = nodes[idx(type,nx,i+1,j+1) ];
182  elem->set_node(3) = nodes[idx(type,nx,i,j+1) ];
183 
184  if (j == 0)
185  boundary_info.add_side(elem, 0, 0);
186 
187  if (j == (ny-1))
188  boundary_info.add_side(elem, 2, 2);
189 
190  if (i == 0)
191  boundary_info.add_side(elem, 3, 3);
192 
193  if (i == (nx-1))
194  boundary_info.add_side(elem, 1, 1);
195 
196  if ((j == 0) &&
197  (i <= .5*(1.+dirichlet_length_fraction) * nx) &&
198  (i >= .5*(1.-dirichlet_length_fraction) * nx))
199  boundary_info.add_side(elem, 0, 6);
200  }
201  break;
202  }
203 
204 
205  case libMesh::QUAD9: {
206 
207  for (uint_t j=0; j<(2*ny); j += 2)
208  for (uint_t i=0; i<(2*nx); i += 2) {
209 
210  libMesh::Elem
211  *elem = libMesh::Elem::build(libMesh::QUAD9).release();
212  elem->set_id(elem_id++);
213  mesh.add_elem(elem);
214 
215  elem->set_node(0) = nodes[idx(type,nx,i, j) ];
216  elem->set_node(1) = nodes[idx(type,nx,i+2,j) ];
217  elem->set_node(2) = nodes[idx(type,nx,i+2,j+2) ];
218  elem->set_node(3) = nodes[idx(type,nx,i, j+2) ];
219  elem->set_node(4) = nodes[idx(type,nx,i+1,j) ];
220  elem->set_node(5) = nodes[idx(type,nx,i+2,j+1) ];
221  elem->set_node(6) = nodes[idx(type,nx,i+1,j+2) ];
222  elem->set_node(7) = nodes[idx(type,nx,i, j+1) ];
223  elem->set_node(8) = nodes[idx(type,nx,i+1,j+1) ];
224 
225  if (j == 0)
226  boundary_info.add_side(elem, 0, 0);
227 
228  if (j == 2*(ny-1))
229  boundary_info.add_side(elem, 2, 2);
230 
231  if (i == 0)
232  boundary_info.add_side(elem, 3, 3);
233 
234  if (i == 2*(nx-1))
235  boundary_info.add_side(elem, 1, 1);
236 
237  if ((j == 0) &&
238  (i <= .5*(1.+dirichlet_length_fraction) * 2*nx) &&
239  (i >= .5*(1.-dirichlet_length_fraction) * 2*nx))
240  boundary_info.add_side(elem, 0, 6);
241  }
242  break;
243  }
244 
245  default:
246  Assert0(false, "ERROR: Unrecognized 2D element type.");
247  }
248 
249  // Add sideset names to boundary info (Z axis out of the screen)
250  boundary_info.sideset_name(0) = "bottom";
251  boundary_info.sideset_name(1) = "right";
252  boundary_info.sideset_name(2) = "top";
253  boundary_info.sideset_name(3) = "left";
254  boundary_info.sideset_name(6) = "dirichlet";
255 
256  // Add nodeset names to boundary info
257  boundary_info.nodeset_name(0) = "bottom";
258  boundary_info.nodeset_name(1) = "right";
259  boundary_info.nodeset_name(2) = "top";
260  boundary_info.nodeset_name(3) = "left";
261 
262  // Done building the mesh. Now prepare it for use.
263  mesh.prepare_for_use ();
264  }
265 
266 
267  template <typename Context>
268  inline void
269  init_analysis_mesh(Context& c,
270  libMesh::UnstructuredMesh& mesh) {
271 
272  real_t
273  length = c.input("length", "length of domain along x-axis", 1.),
274  height = c.input("height", "length of domain along y-axis", 1.),
275  dirichlet_length_fraction = c.input
276  ("dirichlet_length_fraction",
277  "length fraction of the boundary where dirichlet condition is applied",
278  0.1);
279 
280  uint_t
281  nx_divs = c.input("nx_divs", "number of elements along x-axis", 20),
282  ny_divs = c.input("ny_divs", "number of elements along y-axis", 20);
283 
284  std::string
285  t = c.input("elem_type", "type of geometric element in the mesh", "quad4");
286 
287  libMesh::ElemType
288  e_type = libMesh::Utility::string_to_enum<libMesh::ElemType>(t);
289 
290  //
291  // if high order FE is used, libMesh requires atleast a second order
292  // geometric element.
293  //
294  if (c.fe_order > 1 && e_type == libMesh::QUAD4)
295  e_type = libMesh::QUAD9;
296 
297  //
298  // initialize the mesh with one element
299  //
300  build_mesh(mesh,
301  nx_divs, ny_divs,
302  length,
303  height,
304  dirichlet_length_fraction,
305  e_type);
306  }
307 
308 
309 
310  template <typename Context>
311  inline void
313 
314  c.sys->get_dof_map().add_dirichlet_boundary
315  (libMesh::DirichletBoundary({6}, {0}, libMesh::ZeroFunction<real_t>()));
316  }
317 
318 
319 
320 
321  template <typename ScalarType, typename Context>
322  inline void
324  (Context &c,
326 
327  //
328  // this assumes that density variable has a constant value per element
329  //
330  Assert2(c.fe_family == libMesh::LAGRANGE,
331  c.fe_family, libMesh::LAGRANGE,
332  "Method assumes Lagrange interpolation function for density");
333 
334  real_t
335  tol = 1.e-12,
336  length = c.input("length", "length of domain along x-axis", 1.),
337  height = c.input("height", "length of domain along y-axis", 1.),
338  filter_radius = c.input("filter_radius", "radius of geometric filter for level set field", 0.1),
339  vf = c.input("volume_fraction", "upper limit for the volume fraction", 0.3);
340 
341  uint_t
342  sys_num = c.rho_sys->number(),
343  first_dof = c.rho_sys->get_dof_map().first_dof(c.rho_sys->comm().rank()),
344  end_dof = c.rho_sys->get_dof_map().end_dof(c.rho_sys->comm().rank()),
345  dof_id = 0;
346 
347  real_t
348  val = 0.;
349 
350  libMesh::MeshBase::const_element_iterator
351  e_it = c.mesh->active_local_elements_begin(),
352  e_end = c.mesh->active_local_elements_end();
353 
354  std::set<const libMesh::Node*> nodes;
355 
356  for ( ; e_it != e_end; e_it++) {
357 
358  const libMesh::Elem* e = *e_it;
359 
360  for (uint_t i=0; i<e->n_nodes(); i++) {
361 
362  const libMesh::Node& n = *e->node_ptr(i);
363 
364  // if we have alredy operated on this node, then
365  // we skip it
366  if (nodes.count(&n))
367  continue;
368 
369  // otherwise, we add it to the set of operated nodes and
370  // check if a design parameter should be computed for this
371  nodes.insert(&n);
372 
373  dof_id = n.dof_number(sys_num, 0, 0);
374 
377  dv->set_point(n(0), n(1), n(2));
378 
379  if (dof_id >= first_dof &&
380  dof_id < end_dof)
381  dvs.add_topology_parameter(*dv, dof_id);
382  else
383  dvs.add_ghosted_topology_parameter(*dv, dof_id);
384  }
385  }
386 
387  dvs.synchronize(c.rho_sys->get_dof_map());
388  c.rho_sys->solution->close();
389  }
390 };
391 
392 } // namespace Generation
393 } // namespace Mesh
394 } // namespace MAST
395 
396 
397 #endif // __mast_mesh_generation_heat_sink_2d_h__
#define Error(cond, msg)
Definition: exceptions.hpp:166
MAST::Base::ParameterData & add_ghosted_topology_parameter(MAST::Optimization::DesignParameter< ScalarType > &p, const uint_t id)
void init_analysis_mesh(Context &c, libMesh::UnstructuredMesh &mesh)
void init_analysis_dirichlet_conditions(Context &c)
real_t reference_volume(Context &c)
Definition: heat_sink2d.hpp:54
MAST::Base::ParameterData & add_topology_parameter(MAST::Optimization::DesignParameter< ScalarType > &p, const uint_t id)
#define Assert0(cond, msg)
Definition: exceptions.hpp:134
unsigned int uint_t
#define Assert2(cond, v1, v2, msg)
Definition: exceptions.hpp:152
void build_mesh(libMesh::UnstructuredMesh &mesh, const uint_t nx, const uint_t ny, const real_t length, const real_t height, const real_t dirichlet_length_fraction, const libMesh::ElemType type)
Definition: heat_sink2d.hpp:90
void set_point(real_t x, real_t y=0., real_t z=0.)
void init_simp_dvs(Context &c, MAST::Optimization::DesignParameterVector< ScalarType > &dvs)
double real_t
void synchronize(const libMesh::DofMap &dof_map)
std::unique_ptr< ValType > build(const libMesh::System &sys)
Definition: utility.hpp:284
uint_t idx(const libMesh::ElemType type, const uint_t nx, const uint_t i, const uint_t j)
Definition: heat_sink2d.hpp:64