def create_mesh(self, extra_nodes=True): """ Create a mesh from the field region, optionally including the field extra nodes. """ mesh = self.domain.mesh if self.approx_order != 0: conns, mat_ids, descs = [], [], [] for ig, ap in self.aps.iteritems(): group = self.domain.groups[ig] if extra_nodes: conn = ap.econn else: offset = group.shape.n_ep conn = ap.econn[:,:offset] conns.append(conn) mat_ids.append(mesh.mat_ids[ig]) descs.append(mesh.descs[ig]) if extra_nodes: coors = self.coors else: coors = self.coors[:self.n_vertex_dof] mesh = Mesh.from_data(self.name, coors, None, conns, mat_ids, descs) return mesh
def gen_block_mesh(dims, shape, centre, mat_id=0, name='block', coors=None, verbose=True): """ Generate a 2D or 3D block mesh. The dimension is determined by the lenght of the shape argument. Parameters ---------- dims : array of 2 or 3 floats Dimensions of the block. shape : array of 2 or 3 ints Shape (counts of nodes in x, y, z) of the block mesh. centre : array of 2 or 3 floats Centre of the block. mat_id : int, optional The material id of all elements. name : string Mesh name. verbose : bool If True, show progress of the mesh generation. Returns ------- mesh : Mesh instance """ dims = nm.asarray(dims, dtype=nm.float64) shape = nm.asarray(shape, dtype=nm.int32) centre = nm.asarray(centre, dtype=nm.float64) dim = shape.shape[0] centre = centre[:dim] dims = dims[:dim] n_nod = nm.prod(shape) output('generating %d vertices...' % n_nod, verbose=verbose) x0 = centre - 0.5 * dims dd = dims / (shape - 1) ngrid = nm.mgrid[[slice(ii) for ii in shape]] ngrid.shape = (dim, n_nod) coors = x0 + ngrid.T * dd output('...done', verbose=verbose) n_el = nm.prod(shape - 1) output('generating %d cells...' % n_el, verbose=verbose) mat_ids = nm.empty((n_el,), dtype=nm.int32) mat_ids.fill(mat_id) conn, desc = get_tensor_product_conn(shape) output('...done', verbose=verbose) mesh = Mesh.from_data(name, coors, None, [conn], [mat_ids], [desc]) return mesh
def create_mesh_from_control_points( self ): offset = 0 dim = self.spbs[0].cxyz.shape[1] coors = nm.empty((0, dim), dtype=nm.float64) conns = [] mat_ids = [] descs = [] for ib, spb in enumerate( self.spbs ): n_nod = spb.cxyz.shape[0] coors = nm.concatenate( (coors, spb.cxyz), 0 ) descs.append( '3_2' ) conn = [] for ij in xrange( spb.cpi.shape[1] ): for ik in xrange( spb.cpi.shape[2] ): inx = spb.cpi[:,ij,ik] row = [[p1, p2] for p1, p2 in zip( inx[:-1], inx[1:] )] conn.extend( row ) for ij in xrange( spb.cpi.shape[0] ): for ik in xrange( spb.cpi.shape[2] ): inx = spb.cpi[ij,:,ik] row = [[p1, p2] for p1, p2 in zip( inx[:-1], inx[1:] )] conn.extend( row ) for ij in xrange( spb.cpi.shape[0] ): for ik in xrange( spb.cpi.shape[1] ): inx = spb.cpi[ij,ik,:] row = [[p1, p2] for p1, p2 in zip( inx[:-1], inx[1:] )] conn.extend( row ) aux = nm.empty(len(conn), dtype=nm.int32) aux.fill(ib) mat_ids.append(aux) conns.append( offset + nm.array( conn, dtype = nm.int32 ) ) offset += n_nod mesh = Mesh.from_data('control_points', coors, None, conns, mat_ids, descs) return mesh
def gen_cylinder_mesh(dims, shape, centre, axis='x', force_hollow=False, is_open=False, open_angle=0.0, non_uniform=False, name='cylinder', verbose=True): """ Generate a cylindrical mesh along an axis. Its cross-section can be ellipsoidal. Parameters ---------- dims : array of 5 floats Dimensions of the cylinder: inner surface semi-axes a1, b1, outer surface semi-axes a2, b2, length. shape : array of 3 ints Shape (counts of nodes in radial, circumferential and longitudinal directions) of the cylinder mesh. centre : array of 3 floats Centre of the cylinder. axis: one of 'x', 'y', 'z' The axis of the cylinder. force_hollow : boolean Force hollow mesh even if inner radii a1 = b1 = 0. is_open : boolean Generate an open cylinder segment. open_angle : float Opening angle in radians. non_uniform : boolean If True, space the mesh nodes in radial direction so that the element volumes are (approximately) the same, making thus the elements towards the outer surface thinner. name : string Mesh name. verbose : bool If True, show progress of the mesh generation. Returns ------- mesh : Mesh instance """ dims = nm.asarray(dims, dtype=nm.float64) shape = nm.asarray(shape, dtype=nm.int32) centre = nm.asarray(centre, dtype=nm.float64) a1, b1, a2, b2, length = dims nr, nfi, nl = shape origin = centre - nm.array([0.5 * length, 0.0, 0.0]) dfi = 2.0 * (nm.pi - open_angle) / nfi if is_open: nnfi = nfi + 1 else: nnfi = nfi is_hollow = force_hollow or not (max(abs(a1), abs(b1)) < 1e-15) if is_hollow: mr = 0 else: mr = (nnfi - 1) * nl grid = nm.zeros((nr, nnfi, nl), dtype=nm.int32) n_nod = nr * nnfi * nl - mr coors = nm.zeros((n_nod, 3), dtype=nm.float64) angles = nm.linspace(open_angle, open_angle + (nfi) * dfi, nfi + 1) xs = nm.linspace(0.0, length, nl) if non_uniform: ras = nm.zeros((nr, ), dtype=nm.float64) rbs = nm.zeros_like(ras) advol = (a2**2 - a1**2) / (nr - 1) bdvol = (b2**2 - b1**2) / (nr - 1) ras[0], rbs[0] = a1, b1 for ii in range(1, nr): ras[ii] = nm.sqrt(advol + ras[ii - 1]**2) rbs[ii] = nm.sqrt(bdvol + rbs[ii - 1]**2) else: ras = nm.linspace(a1, a2, nr) rbs = nm.linspace(b1, b2, nr) # This is 3D only... bar = MyBar(" nodes:", verbose=verbose) bar.init(n_nod) ii = 0 for ix in range(nr): a, b = ras[ix], rbs[ix] for iy, fi in enumerate(angles[:nnfi]): for iz, x in enumerate(xs): grid[ix, iy, iz] = ii coors[ii] = origin + [x, a * nm.cos(fi), b * nm.sin(fi)] if not (ii % 100): bar.update(ii) ii += 1 if not is_hollow and (ix == 0): if iy > 0: grid[ix, iy, iz] = grid[ix, 0, iz] ii -= 1 assert_(ii == n_nod) n_el = (nr - 1) * nnfi * (nl - 1) conn = nm.zeros((n_el, 8), dtype=nm.int32) bar = MyBar(" elements:", verbose=verbose) bar.init(n_el) ii = 0 for (ix, iy, iz) in cycle([nr - 1, nnfi, nl - 1]): if iy < (nnfi - 1): conn[ii, :] = [ grid[ix, iy, iz], grid[ix + 1, iy, iz], grid[ix + 1, iy + 1, iz], grid[ix, iy + 1, iz], grid[ix, iy, iz + 1], grid[ix + 1, iy, iz + 1], grid[ix + 1, iy + 1, iz + 1], grid[ix, iy + 1, iz + 1] ] ii += 1 elif not is_open: conn[ii, :] = [ grid[ix, iy, iz], grid[ix + 1, iy, iz], grid[ix + 1, 0, iz], grid[ix, 0, iz], grid[ix, iy, iz + 1], grid[ix + 1, iy, iz + 1], grid[ix + 1, 0, iz + 1], grid[ix, 0, iz + 1] ] ii += 1 if not (ii % 100): bar.update(ii) mat_id = nm.zeros((n_el, ), dtype=nm.int32) desc = '3_8' assert_(n_nod == (conn.max() + 1)) if axis == 'z': coors = coors[:, [1, 2, 0]] elif axis == 'y': coors = coors[:, [2, 0, 1]] mesh = Mesh.from_data(name, coors, None, [conn], [mat_id], [desc]) return mesh
def gen_mesh_from_voxels(voxels, dims, etype='q'): """ Generate FE mesh from voxels (volumetric data). Parameters ---------- voxels : array Voxel matrix, 1=material. dims : array Size of one voxel. etype : integer, optional 'q' - quadrilateral or hexahedral elements 't' - triangular or tetrahedral elements Returns ------- mesh : Mesh instance Finite element mesh. """ dims = dims.squeeze() dim = len(dims) nddims = nm.array(voxels.shape) + 2 nodemtx = nm.zeros(nddims, dtype=nm.int32) if dim == 2: #iy, ix = nm.where(voxels.transpose()) iy, ix = nm.where(voxels) nel = ix.shape[0] if etype == 'q': nodemtx[ix, iy] += 1 nodemtx[ix + 1, iy] += 1 nodemtx[ix + 1, iy + 1] += 1 nodemtx[ix, iy + 1] += 1 elif etype == 't': nodemtx[ix, iy] += 2 nodemtx[ix + 1, iy] += 1 nodemtx[ix + 1, iy + 1] += 2 nodemtx[ix, iy + 1] += 1 nel *= 2 elif dim == 3: #iy, ix, iz = nm.where(voxels.transpose(1, 0, 2)) iy, ix, iz = nm.where(voxels) nel = ix.shape[0] if etype == 'q': nodemtx[ix, iy, iz] += 1 nodemtx[ix + 1, iy, iz] += 1 nodemtx[ix + 1, iy + 1, iz] += 1 nodemtx[ix, iy + 1, iz] += 1 nodemtx[ix, iy, iz + 1] += 1 nodemtx[ix + 1, iy, iz + 1] += 1 nodemtx[ix + 1, iy + 1, iz + 1] += 1 nodemtx[ix, iy + 1, iz + 1] += 1 elif etype == 't': nodemtx[ix, iy, iz] += 6 nodemtx[ix + 1, iy, iz] += 2 nodemtx[ix + 1, iy + 1, iz] += 2 nodemtx[ix, iy + 1, iz] += 2 nodemtx[ix, iy, iz + 1] += 2 nodemtx[ix + 1, iy, iz + 1] += 2 nodemtx[ix + 1, iy + 1, iz + 1] += 6 nodemtx[ix, iy + 1, iz + 1] += 2 nel *= 6 else: msg = 'incorrect voxel dimension! (%d)' % dim raise ValueError(msg) ndidx = nm.where(nodemtx) coors = nm.array(ndidx).transpose() * dims nnod = coors.shape[0] nodeid = -nm.ones(nddims, dtype=nm.int32) nodeid[ndidx] = nm.arange(nnod) # generate elements if dim == 2: elems = nm.array([ nodeid[ix, iy], nodeid[ix + 1, iy], nodeid[ix + 1, iy + 1], nodeid[ix, iy + 1] ]).transpose() elif dim == 3: elems = nm.array([ nodeid[ix, iy, iz], nodeid[ix + 1, iy, iz], nodeid[ix + 1, iy + 1, iz], nodeid[ix, iy + 1, iz], nodeid[ix, iy, iz + 1], nodeid[ix + 1, iy, iz + 1], nodeid[ix + 1, iy + 1, iz + 1], nodeid[ix, iy + 1, iz + 1] ]).transpose() if etype == 't': elems = elems_q2t(elems) eid = etype + str(dim) eltab = {'q2': 4, 'q3': 8, 't2': 3, 't3': 4} mesh = Mesh.from_data('voxel_data', coors, nm.ones( (nnod, ), dtype=nm.int32), {0: nm.ascontiguousarray(elems)}, {0: nm.ones((nel, ), dtype=nm.int32)}, {0: '%d_%d' % (dim, eltab[eid])}) return mesh
def gen_tiled_mesh(mesh, grid=None, scale=1.0, eps=1e-6, ret_ndmap=False): """ Generate a new mesh by repeating a given periodic element along each axis. Parameters ---------- mesh : Mesh instance The input periodic FE mesh. grid : array Number of repetition along each axis. scale : float, optional Scaling factor. eps : float, optional Tolerance for boundary detection. ret_ndmap : bool, optional If True, return global node map. Returns ------- mesh_out : Mesh instance FE mesh. ndmap : array Maps: actual node id --> node id in the reference cell. """ bbox = mesh.get_bounding_box() if grid is None: iscale = max(int(1.0 / scale), 1) grid = [iscale] * mesh.dim conns = mesh.conns[0] for ii in mesh.conns[1:]: conns = nm.vstack((conns, ii)) mat_ids = mesh.mat_ids[0] for ii in mesh.mat_ids[1:]: mat_ids = nm.hstack((mat_ids, ii)) coors = mesh.coors ngrps = mesh.ngroups nrep = nm.prod(grid) ndmap = None bar = MyBar(" repeating:") bar.init(nrep) nblk = 1 for ii, gr in enumerate(grid): if ret_ndmap: (conns, coors, ngrps, ndmap0) = tiled_mesh1d(conns, coors, ngrps, ii, gr, bbox.transpose()[ii], eps=eps, mybar=(bar, nblk), ndmap=ndmap) ndmap = ndmap0 else: conns, coors, ngrps = tiled_mesh1d(conns, coors, ngrps, ii, gr, bbox.transpose()[ii], eps=eps, mybar=(bar, nblk)) nblk *= gr bar.update(nblk) mat_ids = nm.tile(mat_ids, (nrep, )) mesh_out = Mesh.from_data('tiled mesh', coors * scale, ngrps, [conns], [mat_ids], [mesh.descs[0]]) if ret_ndmap: return mesh_out, ndmap else: return mesh_out
def gen_block_mesh(dims, shape, centre, mat_id=0, name='block', verbose=True): """ Generate a 2D or 3D block mesh. The dimension is determined by the lenght of the shape argument. Parameters ---------- dims : array of 2 or 3 floats Dimensions of the block. shape : array of 2 or 3 ints Shape (counts of nodes in x, y, z) of the block mesh. centre : array of 2 or 3 floats Centre of the block. mat_id : int, optional The material id of all elements. name : string Mesh name. verbose : bool If True, show progress of the mesh generation. Returns ------- mesh : Mesh instance """ dims = nm.asarray(dims, dtype=nm.float64) shape = nm.asarray(shape, dtype=nm.int32) centre = nm.asarray(centre, dtype=nm.float64) dim = shape.shape[0] centre = centre[:dim] dims = dims[:dim] x0 = centre - 0.5 * dims dd = dims / (shape - 1) grid = nm.zeros(shape, dtype=nm.int32) n_nod = nm.prod(shape) coors = nm.zeros((n_nod, dim), dtype=nm.float64) bar = MyBar(" nodes:", verbose=verbose) bar.init(n_nod) for ii, ic in enumerate(cycle(shape)): grid[tuple(ic)] = ii coors[ii] = x0 + ic * dd if not (ii % 100): bar.update(ii) bar.update(ii + 1) n_el = nm.prod(shape - 1) mat_ids = nm.empty((n_el, ), dtype=nm.int32) mat_ids.fill(mat_id) if (dim == 2): conn = nm.zeros((n_el, 4), dtype=nm.int32) bar = MyBar(" elements:", verbose=verbose) bar.init(n_el) for ii, (ix, iy) in enumerate(cycle(shape - 1)): conn[ii, :] = [ grid[ix, iy], grid[ix + 1, iy], grid[ix + 1, iy + 1], grid[ix, iy + 1] ] if not (ii % 100): bar.update(ii) bar.update(ii + 1) desc = '2_4' else: conn = nm.zeros((n_el, 8), dtype=nm.int32) bar = MyBar(" elements:", verbose=verbose) bar.init(n_el) for ii, (ix, iy, iz) in enumerate(cycle(shape - 1)): conn[ii, :] = [ grid[ix, iy, iz], grid[ix + 1, iy, iz], grid[ix + 1, iy + 1, iz], grid[ix, iy + 1, iz], grid[ix, iy, iz + 1], grid[ix + 1, iy, iz + 1], grid[ix + 1, iy + 1, iz + 1], grid[ix, iy + 1, iz + 1] ] if not (ii % 100): bar.update(ii) bar.update(ii + 1) desc = '3_8' mesh = Mesh.from_data(name, coors, None, [conn], [mat_ids], [desc]) return mesh
def gen_block_mesh(dims, shape, centre, name='block'): """ Generate a 2D or 3D block mesh. The dimension is determined by the lenght of the shape argument. Parameters ---------- dims : array of 2 or 3 floats Dimensions of the block. shape : array of 2 or 3 ints Shape (counts of nodes in x, y, z) of the block mesh. centre : array of 2 or 3 floats Centre of the block. name : string Mesh name. Returns ------- mesh : Mesh instance """ dims = nm.asarray(dims, dtype=nm.float64) shape = nm.asarray(shape, dtype=nm.int32) centre = nm.asarray(centre, dtype=nm.float64) dim = shape.shape[0] centre = centre[:dim] dims = dims[:dim] x0 = centre - 0.5 * dims dd = dims / (shape - 1) grid = nm.zeros(shape, dtype = nm.int32) n_nod = nm.prod(shape) coors = nm.zeros((n_nod, dim), dtype = nm.float64) bar = MyBar(" nodes:") bar.init(n_nod) for ii, ic in enumerate(cycle(shape)): grid[tuple(ic)] = ii coors[ii] = x0 + ic * dd if not (ii % 100): bar.update(ii) bar.update(ii + 1) n_el = nm.prod(shape - 1) mat_id = nm.zeros((n_el,), dtype = nm.int32) if (dim == 2): conn = nm.zeros((n_el, 4), dtype = nm.int32) bar = MyBar(" elements:") bar.init(n_el) for ii, (ix, iy) in enumerate(cycle(shape - 1)): conn[ii,:] = [grid[ix ,iy], grid[ix+1,iy ], grid[ix+1,iy+1], grid[ix ,iy+1]] if not (ii % 100): bar.update(ii) bar.update(ii + 1) desc = '2_4' else: conn = nm.zeros((n_el, 8), dtype = nm.int32) bar = MyBar(" elements:") bar.init(n_el) for ii, (ix, iy, iz) in enumerate(cycle(shape - 1)): conn[ii,:] = [grid[ix ,iy ,iz ], grid[ix+1,iy ,iz ], grid[ix+1,iy+1,iz ], grid[ix ,iy+1,iz ], grid[ix ,iy ,iz+1], grid[ix+1,iy ,iz+1], grid[ix+1,iy+1,iz+1], grid[ix ,iy+1,iz+1]] if not (ii % 100): bar.update(ii) bar.update(ii + 1) desc = '3_8' mesh = Mesh.from_data(name, coors, None, [conn], [mat_id], [desc]) return mesh
def gen_cylinder_mesh(dims, shape, centre, axis='x', force_hollow=False, is_open=False, open_angle=0.0, non_uniform=False, name='cylinder'): """ Generate a cylindrical mesh along an axis. Its cross-section can be ellipsoidal. Parameters ---------- axis: one of 'x', 'y', 'z' The axis of the cylinder. dims : array of 5 floats Dimensions of the cylinder: inner surface semi-axes a1, b1, outer surface semi-axes a2, b2, length. shape : array of 3 ints Shape (counts of nodes in radial, circumferential and longitudinal directions) of the cylinder mesh. centre : array of 3 floats Centre of the cylinder. force_hollow : boolean Force hollow mesh even if inner radii a1 = b1 = 0. is_open : boolean Generate an open cylinder segment. open_angle : float Opening angle in radians. non_uniform : boolean If True, space the mesh nodes in radial direction so that the element volumes are (approximately) the same, making thus the elements towards the outer surface thinner. name : string Mesh name. Returns ------- mesh : Mesh instance """ dims = nm.asarray(dims, dtype=nm.float64) shape = nm.asarray(shape, dtype=nm.int32) centre = nm.asarray(centre, dtype=nm.float64) a1, b1, a2, b2, length = dims nr, nfi, nl = shape origin = centre - nm.array([0.5 * length, 0.0, 0.0]) dfi = 2.0 * (nm.pi - open_angle) / nfi if is_open: nnfi = nfi + 1 else: nnfi = nfi is_hollow = force_hollow or not (max(abs(a1), abs(b1)) < 1e-15) if is_hollow: mr = 0 else: mr = (nnfi - 1) * nl grid = nm.zeros((nr, nnfi, nl), dtype=nm.int32) n_nod = nr * nnfi * nl - mr coors = nm.zeros((n_nod, 3), dtype=nm.float64) angles = nm.linspace(open_angle, open_angle+(nfi)*dfi, nfi+1) xs = nm.linspace(0.0, length, nl) if non_uniform: ras = nm.zeros((nr,), dtype=nm.float64) rbs = nm.zeros_like(ras) advol = (a2**2 - a1**2) / (nr - 1) bdvol = (b2**2 - b1**2) / (nr - 1) ras[0], rbs[0] = a1, b1 for ii in range(1, nr): ras[ii] = nm.sqrt(advol + ras[ii-1]**2) rbs[ii] = nm.sqrt(bdvol + rbs[ii-1]**2) else: ras = nm.linspace(a1, a2, nr) rbs = nm.linspace(b1, b2, nr) # This is 3D only... bar = MyBar(" nodes:") bar.init(n_nod) ii = 0 for ix in range(nr): a, b = ras[ix], rbs[ix] for iy, fi in enumerate(angles[:nnfi]): for iz, x in enumerate(xs): grid[ix,iy,iz] = ii coors[ii] = origin + [x, a * nm.cos(fi), b * nm.sin(fi)] if not (ii % 100): bar.update(ii) ii += 1 if not is_hollow and (ix == 0): if iy > 0: grid[ix,iy,iz] = grid[ix,0,iz] ii -= 1 print assert_(ii == n_nod) n_el = (nr - 1) * nnfi * (nl - 1) conn = nm.zeros((n_el, 8), dtype=nm.int32) bar = MyBar(" elements:") bar.init(n_el) ii = 0 for (ix, iy, iz) in cycle([nr-1, nnfi, nl-1]): if iy < (nnfi - 1): conn[ii,:] = [grid[ix ,iy ,iz ], grid[ix+1,iy ,iz ], grid[ix+1,iy+1,iz ], grid[ix ,iy+1,iz ], grid[ix ,iy ,iz+1], grid[ix+1,iy ,iz+1], grid[ix+1,iy+1,iz+1], grid[ix ,iy+1,iz+1]] ii += 1 elif not is_open: conn[ii,:] = [grid[ix ,iy ,iz ], grid[ix+1,iy ,iz ], grid[ix+1,0,iz ], grid[ix ,0,iz ], grid[ix ,iy ,iz+1], grid[ix+1,iy ,iz+1], grid[ix+1,0,iz+1], grid[ix ,0,iz+1]] ii += 1 if not (ii % 100): bar.update(ii) print mat_id = nm.zeros((n_el,), dtype = nm.int32) desc = '3_8' assert_(n_nod == (conn.max() + 1)) if axis == 'z': coors = coors[:,[1,2,0]] elif axis == 'y': coors = coors[:,[2,0,1]] mesh = Mesh.from_data(name, coors, None, [conn], [mat_id], [desc]) return mesh
def gen_mesh_from_voxels(voxels, dims, etype='q'): """ Generate FE mesh from voxels (volumetric data). Parameters ---------- voxels : array Voxel matrix, 1=material. dims : array Size of one voxel. etype : integer, optional 'q' - quadrilateral or hexahedral elements 't' - triangular or tetrahedral elements Returns ------- mesh : Mesh instance Finite element mesh. """ dims = dims.squeeze() dim = len(dims) nddims = nm.array(voxels.shape) + 1 nodemtx = nm.zeros(nddims, dtype=nm.int32) if dim == 2: #iy, ix = nm.where(voxels.transpose()) iy, ix = nm.where(voxels) nel = ix.shape[0] if etype == 'q': nodemtx[ix,iy] += 1 nodemtx[ix + 1,iy] += 1 nodemtx[ix + 1,iy + 1] += 1 nodemtx[ix,iy + 1] += 1 elif etype == 't': nodemtx[ix,iy] += 2 nodemtx[ix + 1,iy] += 1 nodemtx[ix + 1,iy + 1] += 2 nodemtx[ix,iy + 1] += 1 nel *= 2 elif dim == 3: #iy, ix, iz = nm.where(voxels.transpose(1, 0, 2)) iy, ix, iz = nm.where(voxels) nel = ix.shape[0] if etype == 'q': nodemtx[ix,iy,iz] += 1 nodemtx[ix + 1,iy,iz] += 1 nodemtx[ix + 1,iy + 1,iz] += 1 nodemtx[ix,iy + 1,iz] += 1 nodemtx[ix,iy,iz + 1] += 1 nodemtx[ix + 1,iy,iz + 1] += 1 nodemtx[ix + 1,iy + 1,iz + 1] += 1 nodemtx[ix,iy + 1,iz + 1] += 1 elif etype == 't': nodemtx[ix,iy,iz] += 6 nodemtx[ix + 1,iy,iz] += 2 nodemtx[ix + 1,iy + 1,iz] += 2 nodemtx[ix,iy + 1,iz] += 2 nodemtx[ix,iy,iz + 1] += 2 nodemtx[ix + 1,iy,iz + 1] += 2 nodemtx[ix + 1,iy + 1,iz + 1] += 6 nodemtx[ix,iy + 1,iz + 1] += 2 nel *= 6 else: msg = 'incorrect voxel dimension! (%d)' % dim raise ValueError(msg) ndidx = nm.where(nodemtx) coors = nm.array(ndidx).transpose() * dims nnod = coors.shape[0] nodeid = -nm.ones(nddims, dtype=nm.int32) nodeid[ndidx] = nm.arange(nnod) # generate elements if dim == 2: elems = nm.array([nodeid[ix,iy], nodeid[ix + 1,iy], nodeid[ix + 1,iy + 1], nodeid[ix,iy + 1]]).transpose() elif dim == 3: elems = nm.array([nodeid[ix,iy,iz], nodeid[ix + 1,iy,iz], nodeid[ix + 1,iy + 1,iz], nodeid[ix,iy + 1,iz], nodeid[ix,iy,iz + 1], nodeid[ix + 1,iy,iz + 1], nodeid[ix + 1,iy + 1,iz + 1], nodeid[ix,iy + 1,iz + 1]]).transpose() if etype == 't': elems = elems_q2t(elems) eid = etype + str(dim) eltab = {'q2': 4, 'q3': 8, 't2': 3, 't3': 4} mesh = Mesh.from_data('voxel_data', coors, nm.ones((nnod,), dtype=nm.int32), {0: nm.ascontiguousarray(elems)}, {0: nm.ones((nel,), dtype=nm.int32)}, {0: '%d_%d' % (dim, eltab[eid])}) return mesh
def gen_tiled_mesh(mesh, grid=None, scale=1.0, eps=1e-6, ret_ndmap=False): """ Generate a new mesh by repeating a given periodic element along each axis. Parameters ---------- mesh : Mesh instance The input periodic FE mesh. grid : array Number of repetition along each axis. scale : float, optional Scaling factor. eps : float, optional Tolerance for boundary detection. ret_ndmap : bool, optional If True, return global node map. Returns ------- mesh_out : Mesh instance FE mesh. ndmap : array Maps: actual node id --> node id in the reference cell. """ bbox = mesh.get_bounding_box() if grid is None: iscale = max(int(1.0 / scale), 1) grid = [iscale] * mesh.dim conns = mesh.conns[0] for ii in mesh.conns[1:]: conns = nm.vstack((conns, ii)) mat_ids = mesh.mat_ids[0] for ii in mesh.mat_ids[1:]: mat_ids = nm.hstack((mat_ids, ii)) coors = mesh.coors ngrps = mesh.ngroups nrep = nm.prod(grid) ndmap = None bar = MyBar(" repeating:") bar.init(nrep) nblk = 1 for ii, gr in enumerate(grid): if ret_ndmap: (conns, coors, ngrps, ndmap0) = tiled_mesh1d(conns, coors, ngrps, ii, gr, bbox.transpose()[ii], eps=eps, mybar=(bar, nblk), ndmap=ndmap) ndmap = ndmap0 else: conns, coors, ngrps = tiled_mesh1d(conns, coors, ngrps, ii, gr, bbox.transpose()[ii], eps=eps, mybar=(bar, nblk)) nblk *= gr bar.update(nblk) mat_ids = nm.tile(mat_ids, (nrep,)) mesh_out = Mesh.from_data('tiled mesh', coors * scale, ngrps, [conns], [mat_ids], [mesh.descs[0]]) if ret_ndmap: return mesh_out, ndmap else: return mesh_out
def linearize(self, dofs, min_level=0, max_level=1, eps=1e-4): """ Linearize the solution for post-processing. Parameters ---------- dofs : array, shape (n_nod, n_component) The array of DOFs reshaped so that each column corresponds to one component. min_level : int The minimum required level of mesh refinement. max_level : int The maximum level of mesh refinement. eps : float The relative tolerance parameter of mesh adaptivity. Returns ------- mesh : Mesh instance The adapted, nonconforming, mesh. vdofs : array The DOFs defined in vertices of `mesh`. levels : array of ints The refinement level used for each element group. """ assert_(dofs.ndim == 2) n_nod, dpn = dofs.shape assert_(n_nod == self.n_nod) assert_(dpn == self.shape[0]) vertex_coors = self.coors[:self.n_vertex_dof, :] coors = [] vdofs = [] conns = [] mat_ids = [] levels = [] offset = 0 for ig, ap in self.aps.iteritems(): ps = ap.interp.poly_spaces['v'] gps = ap.interp.gel.interp.poly_spaces['v'] group = self.domain.groups[ig] vertex_conn = ap.econn[:, :group.shape.n_ep] eval_dofs = get_eval_dofs(dofs, ap.econn, ps, ori=ap.ori) eval_coors = get_eval_coors(vertex_coors, vertex_conn, gps) (level, _coors, conn, _vdofs, _mat_ids) = create_output(eval_dofs, eval_coors, group.shape.n_el, ps, min_level=min_level, max_level=max_level, eps=eps) _mat_ids[:] = self.domain.mesh.mat_ids[ig][0] coors.append(_coors) vdofs.append(_vdofs) conns.append(conn + offset) mat_ids.append(_mat_ids) levels.append(level) offset += _coors.shape[0] coors = nm.concatenate(coors, axis=0) vdofs = nm.concatenate(vdofs, axis=0) mesh = Mesh.from_data('linearized_mesh', coors, None, conns, mat_ids, self.domain.mesh.descs) return mesh, vdofs, levels
def create_expression_output(expression, name, primary_field_name, fields, materials, variables, functions=None, mode='eval', term_mode=None, extra_args=None, verbose=True, kwargs=None, min_level=0, max_level=1, eps=1e-4): """ Create output mesh and data for the expression using the adaptive linearizer. Parameters ---------- expression : str The expression to evaluate. name : str The name of the data. primary_field_name : str The name of field that defines the element groups and polynomial spaces. fields : dict The dictionary of fields used in `variables`. materials : Materials instance The materials used in the expression. variables : Variables instance The variables used in the expression. functions : Functions instance, optional The user functions for materials etc. mode : one of 'eval', 'el_avg', 'qp' The evaluation mode - 'qp' requests the values in quadrature points, 'el_avg' element averages and 'eval' means integration over each term region. term_mode : str The term call mode - some terms support different call modes and depending on the call mode different values are returned. extra_args : dict, optional Extra arguments to be passed to terms in the expression. verbose : bool If False, reduce verbosity. kwargs : dict, optional The variables (dictionary of (variable name) : (Variable instance)) to be used in the expression. min_level : int The minimum required level of mesh refinement. max_level : int The maximum level of mesh refinement. eps : float The relative tolerance parameter of mesh adaptivity. Returns ------- out : dict The output dictionary. """ field = fields[primary_field_name] vertex_coors = field.coors[:field.n_vertex_dof, :] coors = [] vdofs = [] conns = [] mat_ids = [] levels = [] offset = 0 for ig, ap in field.aps.iteritems(): ps = ap.interp.poly_spaces['v'] gps = ap.interp.gel.interp.poly_spaces['v'] group = field.domain.groups[ig] vertex_conn = ap.econn[:, :group.shape.n_ep] eval_dofs = get_eval_expression(expression, ig, fields, materials, variables, functions=functions, mode=mode, extra_args=extra_args, verbose=verbose, kwargs=kwargs) eval_coors = get_eval_coors(vertex_coors, vertex_conn, gps) (level, _coors, conn, _vdofs, _mat_ids) = create_output(eval_dofs, eval_coors, group.shape.n_el, ps, min_level=min_level, max_level=max_level, eps=eps) _mat_ids[:] = field.domain.mesh.mat_ids[ig][0] coors.append(_coors) vdofs.append(_vdofs) conns.append(conn + offset) mat_ids.append(_mat_ids) levels.append(level) offset += _coors.shape[0] coors = nm.concatenate(coors, axis=0) vdofs = nm.concatenate(vdofs, axis=0) mesh = Mesh.from_data('linearized_mesh', coors, None, conns, mat_ids, field.domain.mesh.descs) out = {} out[name] = Struct(name='output_data', mode='vertex', data=vdofs, var_name=name, dofs=None, mesh=mesh, levels=levels) out = convert_complex_output(out) return out