20 #ifndef __mast_mesh_generation_truss_2d_h__ 21 #define __mast_mesh_generation_truss_2d_h__ 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> 45 namespace Generation {
50 template <
typename ScalarType>
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.;
66 template <
typename ContextType,
typename ScalarFieldType>
68 const ScalarFieldType& f)
const {
77 template <
typename ScalarType>
80 template <
typename Context>
85 length = c.input(
"length",
"length of domain along x-axis", 0.24),
86 height = c.input(
"height",
"length of domain along y-axis", 0.04);
88 return length * height;
106 return i + (2*nx+1)*j;
110 Error(
false,
"Invalid element type");
113 return libMesh::invalid_uint;
123 const real_t dirichlet_length_fraction,
124 const libMesh::ElemType type) {
126 Assert0(type == libMesh::QUAD4 || type == libMesh::QUAD9,
127 "Method only implemented for Quad4/Quad9");
132 libMesh::BoundaryInfo & boundary_info = mesh.get_boundary_info();
134 mesh.set_mesh_dimension(3);
135 mesh.set_spatial_dimension(3);
136 mesh.reserve_elem(nx*ny);
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));
148 std::map<uint_t, libMesh::Node*> nodes;
156 case libMesh::QUAD4: {
158 for (
uint_t j=0; j<=ny; j++)
159 for (
uint_t i=0; i<=nx; i++) {
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,
172 case libMesh::QUAD9: {
174 for (
uint_t j=0; j<=(2*ny); j++)
175 for (
uint_t i=0; i<=(2*nx); i++) {
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,
189 Assert0(
false,
"ERROR: Unrecognized 2D element type.");
197 case libMesh::QUAD4: {
199 for (
uint_t j=0; j<ny; j++)
200 for (
uint_t i=0; i<nx; i++) {
204 elem->set_id(elem_id++);
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) ];
213 boundary_info.add_side(elem, 0, 0);
216 boundary_info.add_side(elem, 2, 2);
219 boundary_info.add_side(elem, 3, 3);
222 boundary_info.add_side(elem, 1, 1);
224 if (j == 0 && i <= dirichlet_length_fraction * nx)
225 boundary_info.add_side(elem, 0, 6);
227 if (j == 0 && i >= (1.-dirichlet_length_fraction)* nx)
228 boundary_info.add_side(elem, 0, 7);
230 if (j == 0 && i == 0)
231 boundary_info.add_side(elem, 0, 8);
237 case libMesh::QUAD9: {
239 for (
uint_t j=0; j<(2*ny); j += 2)
240 for (
uint_t i=0; i<(2*nx); i += 2) {
244 elem->set_id(elem_id++);
247 elem->set_node(0) = nodes[
idx(type,nx,i, j) ];
248 elem->set_node(1) = nodes[
idx(type,nx,i+2,j) ];
249 elem->set_node(2) = nodes[
idx(type,nx,i+2,j+2) ];
250 elem->set_node(3) = nodes[
idx(type,nx,i, j+2) ];
251 elem->set_node(4) = nodes[
idx(type,nx,i+1,j) ];
252 elem->set_node(5) = nodes[
idx(type,nx,i+2,j+1) ];
253 elem->set_node(6) = nodes[
idx(type,nx,i+1,j+2) ];
254 elem->set_node(7) = nodes[
idx(type,nx,i, j+1) ];
255 elem->set_node(8) = nodes[
idx(type,nx,i+1,j+1) ];
258 boundary_info.add_side(elem, 0, 0);
261 boundary_info.add_side(elem, 2, 2);
264 boundary_info.add_side(elem, 3, 3);
267 boundary_info.add_side(elem, 1, 1);
269 if (j == 0 && i <= dirichlet_length_fraction * 2*nx)
270 boundary_info.add_side(elem, 0, 6);
272 if (j == 0 && i >= (1.-dirichlet_length_fraction)* 2*nx)
273 boundary_info.add_side(elem, 0, 7);
275 if (j == 0 && i == 0)
276 boundary_info.add_side(elem, 0, 8);
282 Assert0(
false,
"ERROR: Unrecognized 2D element type.");
286 boundary_info.sideset_name(0) =
"bottom";
287 boundary_info.sideset_name(1) =
"right";
288 boundary_info.sideset_name(2) =
"top";
289 boundary_info.sideset_name(3) =
"left";
290 boundary_info.sideset_name(6) =
"left_dirichlet";
291 boundary_info.sideset_name(7) =
"right_dirichlet";
294 boundary_info.nodeset_name(0) =
"bottom";
295 boundary_info.nodeset_name(1) =
"right";
296 boundary_info.nodeset_name(2) =
"top";
297 boundary_info.nodeset_name(3) =
"left";
300 mesh.prepare_for_use ();
304 template <
typename Context>
307 libMesh::UnstructuredMesh& mesh) {
310 length = c.input(
"length",
"length of domain along x-axis", 0.24),
311 height = c.input(
"height",
"length of domain along y-axis", 0.04),
312 dirichlet_length_fraction = c.input
313 (
"dirichlet_length_fraction",
314 "length fraction of the truss boundary where dirichlet condition is applied",
318 nx_divs = c.input(
"nx_divs",
"number of elements along x-axis", 30),
319 ny_divs = c.input(
"ny_divs",
"number of elements along y-axis", 5);
322 t = c.input(
"elem_type",
"type of geometric element in the mesh",
"quad4");
325 e_type = libMesh::Utility::string_to_enum<libMesh::ElemType>(t);
331 if (c.sol_fe_order > 1 && e_type == libMesh::QUAD4)
332 e_type = libMesh::QUAD9;
341 dirichlet_length_fraction,
347 template <
typename Context>
351 c.sys->get_dof_map().add_dirichlet_boundary
352 (libMesh::DirichletBoundary({6, 7}, {0, 1}, libMesh::ZeroFunction<real_t>()));
357 template <
typename ScalarType,
typename InitType>
358 std::unique_ptr<pressure_t<ScalarType>>
362 length = c.input(
"length",
"length of domain along x-axis", 0.24),
363 frac = c.input(
"loadlength_fraction",
"fraction of boundary length on which pressure will act", 0.2),
364 p_val = c.input(
"pressure",
"pressure on side of domain", 2.e4);
367 std::unique_ptr<pressure_t<ScalarType>>
375 template <
typename ScalarType,
typename Context>
384 Assert2(c.rho_fe_family == libMesh::LAGRANGE,
385 c.rho_fe_family, libMesh::LAGRANGE,
386 "Method assumes Lagrange interpolation function for density");
387 Assert2(c.rho_fe_order == libMesh::FIRST,
388 c.rho_fe_order, libMesh::FIRST,
389 "Method assumes Lagrange interpolation function for density");
393 length = c.input(
"length",
"length of domain along x-axis", 0.24),
394 height = c.input(
"height",
"length of domain along y-axis", 0.04),
395 frac = c.input(
"loadlength_fraction",
"fraction of boundary length on which pressure will act", 0.2),
396 filter_radius = c.input(
"filter_radius",
"radius of geometric filter for level set field", 0.008),
397 vf = c.input(
"volume_fraction",
"upper limit for the volume fraction", 0.2);
400 sys_num = c.rho_sys->number(),
401 first_dof = c.rho_sys->get_dof_map().first_dof(c.rho_sys->comm().rank()),
402 end_dof = c.rho_sys->get_dof_map().end_dof(c.rho_sys->comm().rank()),
408 libMesh::MeshBase::const_element_iterator
409 e_it = c.mesh->active_local_elements_begin(),
410 e_end = c.mesh->active_local_elements_end();
412 std::set<const libMesh::Node*> nodes;
414 for ( ; e_it != e_end; e_it++) {
416 const libMesh::Elem* e = *e_it;
418 Assert0(e->type() == libMesh::QUAD4 ||
419 e->type() == libMesh::QUAD9,
420 "Method requires Quad4/Quad9 element");
423 for (
uint_t i=0; i<4; i++) {
425 const libMesh::Node& n = *e->node_ptr(i);
436 dof_id = n.dof_number(sys_num, 0, 0);
455 if (dof_id >= first_dof &&
465 c.rho_sys->solution->close();
474 #endif // __mast_mesh_generation_truss_2d_h__
void init_analysis_mesh(Context &c, libMesh::UnstructuredMesh &mesh)
uint_t idx(const libMesh::ElemType type, const uint_t nx, const uint_t i, const uint_t j)
MAST::Base::ParameterData & add_ghosted_topology_parameter(MAST::Optimization::DesignParameter< ScalarType > &p, const uint_t id)
void init_analysis_dirichlet_conditions(Context &c)
std::unique_ptr< pressure_t< ScalarType > > build_pressure_load(InitType &c)
MAST::Base::ParameterData & add_topology_parameter(MAST::Optimization::DesignParameter< ScalarType > &p, const uint_t id)
ScalarType value(ContextType &c) const
ScalarType derivative(ContextType &c, const ScalarFieldType &f) const
real_t reference_volume(Context &c)
#define Assert0(cond, msg)
#define Assert2(cond, v1, v2, msg)
void init_simp_dvs(Context &c, MAST::Optimization::DesignParameterVector< ScalarType > &dvs)
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)
void set_point(real_t x, real_t y=0., real_t z=0.)
void synchronize(const libMesh::DofMap &dof_map)
Pressure(real_t p, real_t l1, real_t frac)
std::unique_ptr< ValType > build(const libMesh::System &sys)