MAST3
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
panel2d.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_panel_2d_h__
21 #define __mast_mesh_generation_panel_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 Panel2D {
49 
50  template <typename ScalarType>
51  class Pressure {
52  public:
54  real_t l1,
55  real_t frac):
56  _p(p), _l1(l1), _frac(frac)
57  {}
58  virtual ~Pressure() {}
59 
60  template <typename ContextType>
61  inline ScalarType value(ContextType& c) const {
62  ScalarType v=(fabs(c.qp_location(0)-_l1*0.5) <= 0.5*_frac*_l1)?_p:0.;
63  return v;
64  }
65 
66  template <typename ContextType, typename ScalarFieldType>
67  inline ScalarType derivative(ContextType& c,
68  const ScalarFieldType& f) const {
69  return 0.;
70  }
71 
72  private:
74  };
75 
76  static const uint_t dim = 2;
77  template <typename ScalarType>
79 
80  template <typename Context>
81  inline real_t
82  reference_volume(Context& c) {
83 
84  real_t
85  length = c.input("length", "length of domain along x-axis", 0.3),
86  height = c.input("height", "length of domain along y-axis", 0.03);
87 
88  return length * height;
89  }
90 
91 
92  inline uint_t idx(const libMesh::ElemType type,
93  const uint_t nx,
94  const uint_t i,
95  const uint_t j)
96  {
97  switch(type)
98  {
99  case libMesh::QUAD4:
100  {
101  return i + (nx+1)*j;
102  }
103 
104  case libMesh::QUAD9:
105  {
106  return i + (2*nx+1)*j;
107  }
108 
109  default:
110  Error(false, "Invalid element type");
111  }
112 
113  return libMesh::invalid_uint;
114  }
115 
116 
117 
118  inline void build_mesh(libMesh::UnstructuredMesh & mesh,
119  const uint_t nx,
120  const uint_t ny,
121  const real_t length,
122  const real_t height,
123  const real_t dirichlet_length_fraction,
124  const libMesh::ElemType type) {
125 
126  Assert0(type == libMesh::QUAD4 || type == libMesh::QUAD9,
127  "Method only implemented for Quad4/Quad9");
128 
129  // Clear the mesh and start from scratch
130  mesh.clear();
131 
132  libMesh::BoundaryInfo & boundary_info = mesh.get_boundary_info();
133 
134  mesh.set_mesh_dimension(3);
135  mesh.set_spatial_dimension(3);
136  mesh.reserve_elem(nx*ny);
137 
138  if (type == libMesh::QUAD4)
139  mesh.reserve_nodes( (nx+1)*(ny+1));
140  else if (type == libMesh::QUAD9)
141  mesh.reserve_nodes( (2*nx+1)*(2*ny+1));
142 
143  real_t
144  xmax = length,
145  ymax = height;
146 
147 
148  std::map<uint_t, libMesh::Node*> nodes;
149 
150  // Build the nodes.
151  uint_t
152  node_id = 0,
153  n = 0;
154  switch (type)
155  {
156  case libMesh::QUAD4: {
157 
158  for (uint_t j=0; j<=ny; j++)
159  for (uint_t i=0; i<=nx; i++) {
160  nodes[node_id] =
161  mesh.add_point(libMesh::Point(static_cast<real_t>(i)/static_cast<real_t>(nx)*length,
162  static_cast<real_t>(j)/static_cast<real_t>(ny)*height,
163  0.),
164  n++);
165  node_id++;
166  }
167 
168 
169  break;
170  }
171 
172  case libMesh::QUAD9: {
173 
174  for (uint_t j=0; j<=(2*ny); j++)
175  for (uint_t i=0; i<=(2*nx); i++) {
176  nodes[node_id] =
177  mesh.add_point(libMesh::Point(static_cast<real_t>(i)/static_cast<real_t>(2*nx)*length,
178  static_cast<real_t>(j)/static_cast<real_t>(2*ny)*height,
179  0.),
180  n++);
181  node_id++;
182  }
183 
184  break;
185  }
186 
187 
188  default:
189  Assert0(false, "ERROR: Unrecognized 2D element type.");
190  }
191 
192  // Build the elements.
193  uint_t
194  elem_id = 0;
195  switch (type) {
196 
197  case libMesh::QUAD4: {
198 
199  for (uint_t j=0; j<ny; j++)
200  for (uint_t i=0; i<nx; i++) {
201 
202  libMesh::Elem
203  *elem = libMesh::Elem::build(libMesh::QUAD4).release();
204  elem->set_id(elem_id++);
205  mesh.add_elem(elem);
206 
207  elem->set_node(0) = nodes[idx(type,nx,i,j) ];
208  elem->set_node(1) = nodes[idx(type,nx,i+1,j) ];
209  elem->set_node(2) = nodes[idx(type,nx,i+1,j+1) ];
210  elem->set_node(3) = nodes[idx(type,nx,i,j+1) ];
211 
212  if (j == 0)
213  boundary_info.add_side(elem, 0, 0);
214 
215  if (j == (ny-1))
216  boundary_info.add_side(elem, 2, 2);
217 
218  if (i == 0)
219  boundary_info.add_side(elem, 3, 3);
220 
221  if (i == (nx-1))
222  boundary_info.add_side(elem, 1, 1);
223 
224  if (i == 0 && j >= (1.-dirichlet_length_fraction) * ny)
225  boundary_info.add_side(elem, 3, 6);
226 
227  if (i == (nx-1) && j >= (1.-dirichlet_length_fraction)* ny)
228  boundary_info.add_side(elem, 1, 7);
229  }
230  break;
231  }
232 
233 
234  case libMesh::QUAD9: {
235 
236  for (uint_t j=0; j<(2*ny); j += 2)
237  for (uint_t i=0; i<(2*nx); i += 2) {
238 
239  libMesh::Elem
240  *elem = libMesh::Elem::build(libMesh::QUAD9).release();
241  elem->set_id(elem_id++);
242  mesh.add_elem(elem);
243 
244  elem->set_node(0) = nodes[idx(type,nx,i, j) ];
245  elem->set_node(1) = nodes[idx(type,nx,i+2,j) ];
246  elem->set_node(2) = nodes[idx(type,nx,i+2,j+2) ];
247  elem->set_node(3) = nodes[idx(type,nx,i, j+2) ];
248  elem->set_node(4) = nodes[idx(type,nx,i+1,j) ];
249  elem->set_node(5) = nodes[idx(type,nx,i+2,j+1) ];
250  elem->set_node(6) = nodes[idx(type,nx,i+1,j+2) ];
251  elem->set_node(7) = nodes[idx(type,nx,i, j+1) ];
252  elem->set_node(8) = nodes[idx(type,nx,i+1,j+1) ];
253 
254  if (j == 0)
255  boundary_info.add_side(elem, 0, 0);
256 
257  if (j == 2*(ny-1))
258  boundary_info.add_side(elem, 2, 2);
259 
260  if (i == 0)
261  boundary_info.add_side(elem, 3, 3);
262 
263  if (i == 2*(nx-1))
264  boundary_info.add_side(elem, 1, 1);
265 
266  if (i == 0 && j >= (1.-dirichlet_length_fraction)* 2*ny)
267  boundary_info.add_side(elem, 3, 6);
268 
269  if (i == 2*(nx-1) && j >= (1.-dirichlet_length_fraction)* 2*ny)
270  boundary_info.add_side(elem, 1, 7);
271  }
272  break;
273  }
274 
275  default:
276  Assert0(false, "ERROR: Unrecognized 2D element type.");
277  }
278 
279  // Add sideset names to boundary info (Z axis out of the screen)
280  boundary_info.sideset_name(0) = "bottom";
281  boundary_info.sideset_name(1) = "right";
282  boundary_info.sideset_name(2) = "top";
283  boundary_info.sideset_name(3) = "left";
284  boundary_info.sideset_name(6) = "left_dirichlet";
285  boundary_info.sideset_name(7) = "right_dirichlet";
286 
287  // Add nodeset names to boundary info
288  boundary_info.nodeset_name(0) = "bottom";
289  boundary_info.nodeset_name(1) = "right";
290  boundary_info.nodeset_name(2) = "top";
291  boundary_info.nodeset_name(3) = "left";
292 
293  // Done building the mesh. Now prepare it for use.
294  mesh.prepare_for_use ();
295  }
296 
297 
298  template <typename Context>
299  inline void
300  init_analysis_mesh(Context& c,
301  libMesh::UnstructuredMesh& mesh) {
302 
303  real_t
304  length = c.input("length", "length of domain along x-axis", 0.3),
305  height = c.input("height", "length of domain along y-axis", 0.03),
306  dirichlet_length_fraction = c.input
307  ("dirichlet_length_fraction",
308  "length fraction of the truss boundary where dirichlet condition is applied",
309  0.1);
310 
311  uint_t
312  nx_divs = c.input("nx_divs", "number of elements along x-axis", 100),
313  ny_divs = c.input("ny_divs", "number of elements along y-axis", 10);
314 
315  std::string
316  t = c.input("elem_type", "type of geometric element in the mesh", "quad4");
317 
318  libMesh::ElemType
319  e_type = libMesh::Utility::string_to_enum<libMesh::ElemType>(t);
320 
321  //
322  // if high order FE is used, libMesh requires atleast a second order
323  // geometric element.
324  //
325  if (c.sol_fe_order > 1 && e_type == libMesh::QUAD4)
326  e_type = libMesh::QUAD9;
327 
328  //
329  // initialize the mesh with one element
330  //
331  build_mesh(mesh,
332  nx_divs, ny_divs,
333  length,
334  height,
335  dirichlet_length_fraction,
336  e_type);
337  }
338 
339 
340 
341  template <typename Context>
342  inline void
344 
345  c.sys->get_dof_map().add_dirichlet_boundary
346  (libMesh::DirichletBoundary({6, 7}, {0, 1}, libMesh::ZeroFunction<real_t>()));
347  }
348 
349 
350 
351  template <typename ScalarType, typename InitType>
352  std::unique_ptr<pressure_t<ScalarType>>
353  build_pressure_load(InitType& c) {
354 
355  real_t
356  length = c.input("length", "length of domain along x-axis", 0.3),
357  frac = c.input("loadlength_fraction", "fraction of boundary length on which pressure will act", 1.0),
358  p_val = c.input("pressure", "pressure on side of domain", 2.e6);
359  c.p_side_id = 2;
360 
361  std::unique_ptr<pressure_t<ScalarType>>
362  press(new pressure_t<ScalarType>(p_val, length, frac));
363 
364  return press;
365  }
366 
367 
368 
369  template <typename ScalarType, typename Context>
370  inline void
372  (Context &c,
374 
375  //
376  // this assumes that density variable has a constant value per element
377  //
378  Assert2(c.rho_fe_family == libMesh::LAGRANGE,
379  c.rho_fe_family, libMesh::LAGRANGE,
380  "Method assumes Lagrange interpolation function for density");
381  Assert2(c.rho_fe_order == libMesh::FIRST,
382  c.rho_fe_order, libMesh::FIRST,
383  "Method assumes Lagrange interpolation function for density");
384 
385  real_t
386  vf = c.input("volume_fraction", "upper limit for the volume fraction", 0.2);
387 
388  uint_t
389  sys_num = c.rho_sys->number(),
390  first_dof = c.rho_sys->get_dof_map().first_dof(c.rho_sys->comm().rank()),
391  end_dof = c.rho_sys->get_dof_map().end_dof(c.rho_sys->comm().rank()),
392  dof_id = 0;
393 
394  real_t
395  val = 0.;
396 
397  libMesh::MeshBase::const_element_iterator
398  e_it = c.mesh->active_local_elements_begin(),
399  e_end = c.mesh->active_local_elements_end();
400 
401  std::set<const libMesh::Node*> nodes;
402 
403  for ( ; e_it != e_end; e_it++) {
404 
405  const libMesh::Elem* e = *e_it;
406 
407  Assert0(e->type() == libMesh::QUAD4 ||
408  e->type() == libMesh::QUAD9,
409  "Method requires Quad4/Quad9 element");
410 
411  // only the first four nodes of the quad element are used
412  for (uint_t i=0; i<4; i++) {
413 
414  const libMesh::Node& n = *e->node_ptr(i);
415 
416  // if we have alredy operated on this node, then
417  // we skip it
418  if (nodes.count(&n))
419  continue;
420 
421  // otherwise, we add it to the set of operated nodes and
422  // check if a design parameter should be computed for this
423  nodes.insert(&n);
424 
425  dof_id = n.dof_number(sys_num, 0, 0);
426 
429  dv->set_point(n(0), n(1), n(2));
430 
431  if (dof_id >= first_dof &&
432  dof_id < end_dof)
433  dvs.add_topology_parameter(*dv, dof_id);
434  else
435  dvs.add_ghosted_topology_parameter(*dv, dof_id);
436  }
437  }
438 
439  dvs.synchronize(c.rho_sys->get_dof_map());
440  c.rho_sys->solution->close();
441  }
442 };
443 
444 } // namespace Generation
445 } // namespace Mesh
446 } // namespace MAST
447 
448 
449 #endif // __mast_mesh_generation_truss_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_dirichlet_conditions(Context &c)
Definition: panel2d.hpp:343
real_t reference_volume(Context &c)
Definition: panel2d.hpp:82
ScalarType value(ContextType &c) const
Definition: panel2d.hpp:61
void init_simp_dvs(Context &c, MAST::Optimization::DesignParameterVector< ScalarType > &dvs)
Definition: panel2d.hpp:372
uint_t idx(const libMesh::ElemType type, const uint_t nx, const uint_t i, const uint_t j)
Definition: panel2d.hpp:92
MAST::Base::ParameterData & add_topology_parameter(MAST::Optimization::DesignParameter< ScalarType > &p, const uint_t id)
Pressure(real_t p, real_t l1, real_t frac)
Definition: panel2d.hpp:53
#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: panel2d.hpp:118
void init_analysis_mesh(Context &c, libMesh::UnstructuredMesh &mesh)
Definition: panel2d.hpp:300
void set_point(real_t x, real_t y=0., real_t z=0.)
double real_t
void synchronize(const libMesh::DofMap &dof_map)
static const uint_t dim
Definition: panel2d.hpp:76
std::unique_ptr< ValType > build(const libMesh::System &sys)
Definition: utility.hpp:284
ScalarType derivative(ContextType &c, const ScalarFieldType &f) const
Definition: panel2d.hpp:67
std::unique_ptr< pressure_t< ScalarType > > build_pressure_load(InitType &c)
Definition: panel2d.hpp:353