20 #ifndef __mast_mesh_generation_heat_sink_2d_h__ 21 #define __mast_mesh_generation_heat_sink_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 {
52 template <
typename Context>
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);
60 return length * height;
78 return i + (2*nx+1)*j;
82 Error(
false,
"Invalid element type");
85 return libMesh::invalid_uint;
95 const real_t dirichlet_length_fraction,
96 const libMesh::ElemType type) {
98 Assert0(type == libMesh::QUAD4 || type == libMesh::QUAD9,
99 "Method only implemented for Quad4/Quad9");
104 libMesh::BoundaryInfo & boundary_info = mesh.get_boundary_info();
106 mesh.set_mesh_dimension(3);
107 mesh.set_spatial_dimension(3);
108 mesh.reserve_elem(nx*ny);
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));
120 std::map<uint_t, libMesh::Node*> nodes;
128 case libMesh::QUAD4: {
130 for (
uint_t j=0; j<=ny; j++)
131 for (
uint_t i=0; i<=nx; i++) {
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,
144 case libMesh::QUAD9: {
146 for (
uint_t j=0; j<=(2*ny); j++)
147 for (
uint_t i=0; i<=(2*nx); i++) {
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,
161 Assert0(
false,
"ERROR: Unrecognized 2D element type.");
169 case libMesh::QUAD4: {
171 for (
uint_t j=0; j<ny; j++)
172 for (
uint_t i=0; i<nx; i++) {
176 elem->set_id(elem_id++);
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) ];
185 boundary_info.add_side(elem, 0, 0);
188 boundary_info.add_side(elem, 2, 2);
191 boundary_info.add_side(elem, 3, 3);
194 boundary_info.add_side(elem, 1, 1);
197 (i <= .5*(1.+dirichlet_length_fraction) * nx) &&
198 (i >= .5*(1.-dirichlet_length_fraction) * nx))
199 boundary_info.add_side(elem, 0, 6);
205 case libMesh::QUAD9: {
207 for (
uint_t j=0; j<(2*ny); j += 2)
208 for (
uint_t i=0; i<(2*nx); i += 2) {
212 elem->set_id(elem_id++);
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) ];
226 boundary_info.add_side(elem, 0, 0);
229 boundary_info.add_side(elem, 2, 2);
232 boundary_info.add_side(elem, 3, 3);
235 boundary_info.add_side(elem, 1, 1);
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);
246 Assert0(
false,
"ERROR: Unrecognized 2D element type.");
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";
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";
263 mesh.prepare_for_use ();
267 template <
typename Context>
270 libMesh::UnstructuredMesh& mesh) {
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",
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);
285 t = c.input(
"elem_type",
"type of geometric element in the mesh",
"quad4");
288 e_type = libMesh::Utility::string_to_enum<libMesh::ElemType>(t);
294 if (c.fe_order > 1 && e_type == libMesh::QUAD4)
295 e_type = libMesh::QUAD9;
304 dirichlet_length_fraction,
310 template <
typename Context>
314 c.sys->get_dof_map().add_dirichlet_boundary
315 (libMesh::DirichletBoundary({6}, {0}, libMesh::ZeroFunction<real_t>()));
321 template <
typename ScalarType,
typename Context>
330 Assert2(c.fe_family == libMesh::LAGRANGE,
331 c.fe_family, libMesh::LAGRANGE,
332 "Method assumes Lagrange interpolation function for density");
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);
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()),
350 libMesh::MeshBase::const_element_iterator
351 e_it = c.mesh->active_local_elements_begin(),
352 e_end = c.mesh->active_local_elements_end();
354 std::set<const libMesh::Node*> nodes;
356 for ( ; e_it != e_end; e_it++) {
358 const libMesh::Elem* e = *e_it;
360 for (
uint_t i=0; i<e->n_nodes(); i++) {
362 const libMesh::Node& n = *e->node_ptr(i);
373 dof_id = n.dof_number(sys_num, 0, 0);
379 if (dof_id >= first_dof &&
388 c.rho_sys->solution->close();
397 #endif // __mast_mesh_generation_heat_sink_2d_h__
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)
MAST::Base::ParameterData & add_topology_parameter(MAST::Optimization::DesignParameter< ScalarType > &p, const uint_t id)
#define Assert0(cond, msg)
#define Assert2(cond, v1, v2, msg)
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 init_simp_dvs(Context &c, MAST::Optimization::DesignParameterVector< ScalarType > &dvs)
void synchronize(const libMesh::DofMap &dof_map)
std::unique_ptr< ValType > build(const libMesh::System &sys)
uint_t idx(const libMesh::ElemType type, const uint_t nx, const uint_t i, const uint_t j)