discretize.CylMesh¶
-
class
discretize.
CylMesh
(h=None, x0=None, **kwargs)[source]¶ Bases:
discretize.base.base_tensor_mesh.BaseTensorMesh
,discretize.base.base_mesh.BaseRectangularMesh
,discretize.InnerProducts.InnerProducts
,discretize.View.CylView
,discretize.DiffOperators.DiffOperators
CylMesh is a mesh class for cylindrical problems. It supports both cylindrically symmetric and 3D cylindrical meshes that include an azimuthal discretization.
For a cylindrically symmetric mesh use
h = [hx, 1, hz]
. For example:import discretize from discretize import utils cs, nc, npad = 20., 30, 8 hx = utils.meshTensor([(cs, npad+10, -0.7), (cs, nc), (cs, npad, 1.3)]) hz = utils.meshTensor([(cs, npad ,-1.3), (cs, nc), (cs, npad, 1.3)]) mesh = discretize.CylMesh([hx, 1, hz], x0=[0, 0, -hz.sum()/2]) mesh.plotGrid()
To create a 3D cylindrical mesh, we also include an azimuthal discretization
import discretize from discretize import utils cs, nc, npad = 20., 30, 8 nc_theta = 8 hx = utils.meshTensor([(cs, npad+10, -0.7), (cs, nc), (cs, npad, 1.3)]) hy = 2 * np.pi/nc_theta * np.ones(nc_theta) hz = utils.meshTensor([(cs,npad, -1.3), (cs,nc), (cs, npad, 1.3)]) mesh = discretize.CylMesh([hx, hy, hz], x0=[0, 0, -hz.sum()/2]) mesh.plotGrid()
Required Properties:
- axis_u (
Vector3
): Vector orientation of u-direction. For more details see the docs for therotation_matrix
property., a 3D Vector of <class ‘float’> with shape (3), Default: X - axis_v (
Vector3
): Vector orientation of v-direction. For more details see the docs for therotation_matrix
property., a 3D Vector of <class ‘float’> with shape (3), Default: Y - axis_w (
Vector3
): Vector orientation of w-direction. For more details see the docs for therotation_matrix
property., a 3D Vector of <class ‘float’> with shape (3), Default: Z - cartesianOrigin (
Array
): Cartesian origin of the mesh, a list or numpy array of <class ‘float’> with shape (*) - h (a list of
Array
): h is a list containing the cell widths of the tensor mesh in each dimension., a list (each item is a list or numpy array of <class ‘float’> with shape (*)) with length between 0 and 3 - reference_system (
String
): The type of coordinate reference frame. Can take on the values cartesian, cylindrical, or spherical. Abbreviations of these are allowed., a unicode string, Default: cartesian - x0 (
Array
): origin of the mesh (dim, ), a list or numpy array of <class ‘float’> with shape (*)
- axis_u (
Attributes¶
-
CylMesh.
area
¶ Face areas
For a 3D cyl mesh: [radial, azimuthal, vertical], while a cylindrically symmetric mesh doesn’t have y-Faces, so it returns [radial, vertical]
Returns: face areas Return type: numpy.ndarray
-
CylMesh.
areaFx
¶ Area of the x-faces (radial faces). Radial faces exist on all cylindrical meshes
\[A_x = r \theta h_z\]Returns: area of x-faces Return type: numpy.ndarray
-
CylMesh.
areaFy
¶ Area of y-faces (Azimuthal faces). Azimuthal faces exist only on 3D cylindrical meshes.
\[A_y = h_x h_z\]Returns: area of y-faces Return type: numpy.ndarray
-
CylMesh.
areaFz
¶ Area of z-faces.
\[A_z = \frac{\theta}{2} (r_2^2 - r_1^2)z\]Returns: area of the z-faces Return type: numpy.ndarray
-
CylMesh.
aveCC2F
¶ Construct the averaging operator on cell centers to faces.
-
CylMesh.
aveCCV2F
¶ Construct the averaging operator on cell centers to faces as a vector.
-
CylMesh.
aveE2CC
¶ averaging operator of edges to cell centers
Returns: matrix that averages from edges to cell centers Return type: scipy.sparse.csr_matrix
-
CylMesh.
aveE2CCV
¶ averaging operator of edges to a cell centered vector
Returns: matrix that averages from edges to cell centered vectors Return type: scipy.sparse.csr_matrix
-
CylMesh.
aveEx2CC
¶ averaging operator of x-edges (radial) to cell centers
Returns: matrix that averages from x-edges to cell centers Return type: scipy.sparse.csr_matrix
-
CylMesh.
aveEy2CC
¶ averaging operator of y-edges (azimuthal) to cell centers
Returns: matrix that averages from y-edges to cell centers Return type: scipy.sparse.csr_matrix
-
CylMesh.
aveEz2CC
¶ averaging operator of z-edges to cell centers
Returns: matrix that averages from z-edges to cell centers Return type: scipy.sparse.csr_matrix
-
CylMesh.
aveF2CC
¶ averaging operator of faces to cell centers
Returns: matrix that averages from faces to cell centers Return type: scipy.sparse.csr_matrix
-
CylMesh.
aveF2CCV
¶ averaging operator of x-faces (radial) to cell centered vectors
Returns: matrix that averages from faces to cell centered vectors Return type: scipy.sparse.csr_matrix
-
CylMesh.
aveFx2CC
¶ averaging operator of x-faces (radial) to cell centers
Returns: matrix that averages from x-faces to cell centers Return type: scipy.sparse.csr_matrix
-
CylMesh.
aveFy2CC
¶ averaging operator of y-faces (azimuthal) to cell centers
Returns: matrix that averages from y-faces to cell centers Return type: scipy.sparse.csr_matrix
-
CylMesh.
aveFz2CC
¶ averaging operator of z-faces (vertical) to cell centers
Returns: matrix that averages from z-faces to cell centers Return type: scipy.sparse.csr_matrix
-
CylMesh.
aveN2CC
¶ Construct the averaging operator on cell nodes to cell centers.
-
CylMesh.
aveN2E
¶ Construct the averaging operator on cell nodes to cell edges, keeping each dimension separate.
-
CylMesh.
aveN2F
¶ Construct the averaging operator on cell nodes to cell faces, keeping each dimension separate.
-
CylMesh.
axis_u
¶ X
Type: axis_u ( Vector3
)Type: Vector orientation of u-direction. For more details see the docs for the rotation_matrix
property., a 3D Vector of <class ‘float’> with shape (3), Default
-
CylMesh.
axis_v
¶ Y
Type: axis_v ( Vector3
)Type: Vector orientation of v-direction. For more details see the docs for the rotation_matrix
property., a 3D Vector of <class ‘float’> with shape (3), Default
-
CylMesh.
axis_w
¶ Z
Type: axis_w ( Vector3
)Type: Vector orientation of w-direction. For more details see the docs for the rotation_matrix
property., a 3D Vector of <class ‘float’> with shape (3), Default
-
CylMesh.
cartesianOrigin
¶ Cartesian origin of the mesh, a list or numpy array of <class ‘float’> with shape (*)
Type: cartesianOrigin ( Array
)
-
CylMesh.
cellGrad
¶ The cell centered Gradient, takes you to cell faces.
-
CylMesh.
cellGradBC
¶ The cell centered Gradient boundary condition matrix
-
CylMesh.
cellGradx
¶ Cell centered Gradient in the x dimension. Has neumann boundary conditions.
-
CylMesh.
cellGrady
¶
-
CylMesh.
cellGradz
¶ Cell centered Gradient in the x dimension. Has neumann boundary conditions.
-
CylMesh.
edge
¶ Edge lengths
Returns: vector of edge lengths \((r, \theta, z)\) Return type: numpy.ndarray
-
CylMesh.
edgeCurl
¶ The edgeCurl (edges to faces)
Returns: edge curl operator Return type: scipy.sparse.csr_matrix
-
CylMesh.
edgeEx
¶ x-edge lengths - these are the radial edges. Radial edges only exist for a 3D cyl mesh.
Returns: vector of radial edge lengths Return type: numpy.ndarray
-
CylMesh.
edgeEy
¶ y-edge lengths - these are the azimuthal edges. Azimuthal edges exist for all cylindrical meshes. These are arc-lengths (\(\theta r\))
Returns: vector of the azimuthal edges Return type: numpy.ndarray
-
CylMesh.
edgeEz
¶ z-edge lengths - these are the vertical edges. Vertical edges only exist for a 3D cyl mesh.
Returns: vector of the vertical edges Return type: numpy.ndarray
-
CylMesh.
faceDiv
¶ Construct divergence operator (faces to cell-centres).
-
CylMesh.
faceDivx
¶ Construct divergence operator in the x component (faces to cell-centres).
-
CylMesh.
faceDivy
¶ Construct divergence operator in the y component (faces to cell-centres).
-
CylMesh.
faceDivz
¶ Construct divergence operator in the z component (faces to cell-centres).
-
CylMesh.
gridCC
¶ Cell-centered grid.
-
CylMesh.
gridEx
¶ Edge staggered grid in the x direction.
-
CylMesh.
gridEy
¶ Grid of y-edges (azimuthal-faces) in cylindrical coordinates \((r, \theta, z)\).
Returns: grid locations of azimuthal faces Return type: numpy.ndarray
-
CylMesh.
gridEz
¶ Grid of z-faces (vertical-faces) in cylindrical coordinates \((r, \theta, z)\).
Returns: grid locations of radial faces Return type: numpy.ndarray
-
CylMesh.
gridFx
¶ Grid of x-faces (radial-faces) in cylindrical coordinates \((r, \theta, z)\).
Returns: grid locations of radial faces Return type: numpy.ndarray
-
CylMesh.
gridFy
¶ Face staggered grid in the y direction.
-
CylMesh.
gridFz
¶ Face staggered grid in the z direction.
-
CylMesh.
gridN
¶ Nodal grid in cylindrical coordinates \((r, \theta, z)\). Nodes do not exist in a cylindrically symmetric mesh.
Returns: grid locations of nodes Return type: numpy.ndarray
-
CylMesh.
h
¶ h is a list containing the cell widths of the tensor mesh in each dimension., a list (each item is a list or numpy array of <class ‘float’> with shape (*)) with length between 0 and 3
Type: h (a list of Array
)
-
CylMesh.
h_gridded
¶ Returns an (nC, dim) numpy array with the widths of all cells in order
-
CylMesh.
hx
¶ Width of cells in the x direction
-
CylMesh.
hy
¶ Width of cells in the y direction
-
CylMesh.
hz
¶ Width of cells in the z direction
-
CylMesh.
isSymmetric
¶ Is the mesh cylindrically symmetric?
Returns: True if the mesh is cylindrically symmetric, False otherwise Return type: bool
-
CylMesh.
nE
¶ Total number of edges.
Returns: nE Return type: int = sum([nEx, nEy, nEz])
-
CylMesh.
nEz
¶ returns: Number of z-edges :rtype: int
-
CylMesh.
nN
¶ returns: Total number of nodes :rtype: int
-
CylMesh.
nNx
¶ returns: Number of nodes in the x-direction :rtype: int
-
CylMesh.
nNy
¶ returns: Number of nodes in the y-direction :rtype: int
-
CylMesh.
nodalGrad
¶ Construct gradient operator (nodes to edges).
-
CylMesh.
nodalLaplacian
¶ Construct laplacian operator (nodes to edges).
-
CylMesh.
normals
¶ Face Normals
Return type: numpy.ndarray Returns: normals, (sum(nF), dim)
-
CylMesh.
reference_is_rotated
¶ True if the axes are rotated from the traditional <X,Y,Z> system with vectors of \((1,0,0)\), \((0,1,0)\), and \((0,0,1)\)
-
CylMesh.
reference_system
¶ cartesian
Type: reference_system ( String
)Type: The type of coordinate reference frame. Can take on the values cartesian, cylindrical, or spherical. Abbreviations of these are allowed., a unicode string, Default
-
CylMesh.
rotation_matrix
¶ Builds a rotation matrix to transform coordinates from their coordinate system into a conventional cartesian system. This is built off of the three axis_u, axis_v, and axis_w properties; these mapping coordinates use the letters U, V, and W (the three letters preceding X, Y, and Z in the alphabet) to define the projection of the X, Y, and Z durections. These UVW vectors describe the placement and transformation of the mesh’s coordinate sytem assuming at most 3 directions.
Why would you want to use these UVW mapping vectors the this rotation_matrix property? They allow us to define the realationship between local and global coordinate systems and provide a tool for switching between the two while still maintaing the connectivity of the mesh’s cells. For a visual example of this, please see the figure in the docs for the
vtkInterface
.
-
CylMesh.
tangents
¶ Edge Tangents
Return type: numpy.ndarray Returns: normals, (sum(nE), dim)
-
CylMesh.
vectorCCx
¶ Cell-centered grid vector (1D) in the x direction.
-
CylMesh.
vectorCCy
¶ Cell-centered grid vector (1D) in the y direction.
-
CylMesh.
vectorCCz
¶ Cell-centered grid vector (1D) in the z direction.
-
CylMesh.
vectorNx
¶ Nodal grid vector (1D) in the x direction.
-
CylMesh.
vectorNy
¶ Nodal grid vector (1D) in the y direction.
-
CylMesh.
vectorNz
¶ Nodal grid vector (1D) in the z direction.
-
CylMesh.
vnC
¶ Total number of cells in each direction
Return type: numpy.ndarray Returns: [nCx, nCy, nCz]
-
CylMesh.
vnE
¶ Total number of edges in each direction
Returns: - vnE (numpy.ndarray = [nEx, nEy, nEz], (dim, ))
- .. plot:: – :include-source:
import discretize import numpy as np M = discretize.TensorMesh([np.ones(n) for n in [2,3]]) M.plotGrid(edges=True, showIt=True)
-
CylMesh.
vnEx
¶ Number of x-edges in each direction
Return type: numpy.ndarray Returns: vnEx
-
CylMesh.
vnEy
¶ Number of y-edges in each direction
Returns: vnEy or None if dim < 2, (dim, ) Return type: numpy.ndarray
-
CylMesh.
vnEz
¶ returns: Number of z-edges in each direction or None if nCy > 1, (dim, ) :rtype: numpy.ndarray
-
CylMesh.
vnF
¶ Total number of faces in each direction
Return type: numpy.ndarray Returns: [nFx, nFy, nFz], (dim, ) import discretize import numpy as np M = discretize.TensorMesh([np.ones(n) for n in [2,3]]) M.plotGrid(faces=True, showIt=True)
(Source code, png, hires.png, pdf)
-
CylMesh.
vnFx
¶ returns: Number of x-faces in each direction, (dim, ) :rtype: numpy.ndarray
-
CylMesh.
vnFy
¶ Number of y-faces in each direction
Return type: numpy.ndarray Returns: vnFy or None if dim < 2
-
CylMesh.
vnFz
¶ Number of z-faces in each direction
Return type: numpy.ndarray Returns: vnFz or None if dim < 3
-
CylMesh.
vnN
¶ Total number of nodes in each direction
Return type: numpy.ndarray Returns: [nNx, nNy, nNz]
-
CylMesh.
vol
¶ Volume of each cell
Returns: cell volumes Return type: numpy.ndarray
Methods¶
-
CylMesh.
cartesianGrid
(locType='CC', theta_shift=None)[source]¶ Takes a grid location (‘CC’, ‘N’, ‘Ex’, ‘Ey’, ‘Ez’, ‘Fx’, ‘Fy’, ‘Fz’) and returns that grid in cartesian coordinates
Parameters: locType (str) – grid location Returns: cartesian coordinates for the cylindrical grid Return type: numpy.ndarray
-
CylMesh.
copy
()¶ Make a copy of the current mesh
-
classmethod
CylMesh.
deserialize
(value, trusted=False, strict=False, assert_valid=False, **kwargs)¶ Creates HasProperties instance from serialized dictionary
This uses the Property deserializers to deserialize all JSON-compatible dictionary values into their corresponding Property values on a new instance of a HasProperties class. Extra keys in the dictionary that do not correspond to Properties will be ignored.
Parameters:
- value - Dictionary to deserialize new instance from.
- trusted - If True (and if the input dictionary has
'__class__'
keyword and this class is in the registry), the new HasProperties class will come from the dictionary. If False (the default), only the HasProperties class this method is called on will be constructed. - strict - Requires
'__class__'
, if present on the input dictionary, to match the deserialized instance’s class. Also disallows unused properties in the input dictionary. Default is False. - assert_valid - Require deserialized instance to be valid. Default is False.
- Any other keyword arguments will be passed through to the Property deserializers.
-
CylMesh.
equal
(other)¶ Determine if two HasProperties instances are equivalent
Equivalence is determined by checking if all Property values on two instances are equal, using
Property.equal
.
-
CylMesh.
getBCProjWF
(BC, discretization='CC')¶ The weak form boundary condition projection matrices.
Example
# Neumann in all directions BC = 'neumann' # 3D, Dirichlet in y Neumann else BC = ['neumann', 'dirichlet', 'neumann'] # 3D, Neumann in x on bottom of domain, Dirichlet else BC = [['neumann', 'dirichlet'], 'dirichlet', 'dirichlet']
-
CylMesh.
getBCProjWF_simple
(discretization='CC')¶ The weak form boundary condition projection matrices when mixed boundary condition is used
-
CylMesh.
getEdgeInnerProduct
(prop=None, invProp=False, invMat=False, doFast=True)¶ Generate the edge inner product matrix
Parameters: - prop (numpy.ndarray) – material property (tensor properties are possible) at each cell center (nC, (1, 3, or 6))
- invProp (bool) – inverts the material property
- invMat (bool) – inverts the matrix
- doFast (bool) – do a faster implementation if available.
Returns: M, the inner product matrix (nE, nE)
Return type:
-
CylMesh.
getEdgeInnerProductDeriv
(prop, doFast=True, invProp=False, invMat=False)¶ Parameters: - prop (numpy.ndarray) – material property (tensor properties are possible) at each cell center (nC, (1, 3, or 6))
- doFast (bool) – do a faster implementation if available.
- invProp (bool) – inverts the material property
- invMat (bool) – inverts the matrix
Returns: dMdm, the derivative of the inner product matrix (nE, nC*nA)
Return type:
-
CylMesh.
getFaceInnerProduct
(prop=None, invProp=False, invMat=False, doFast=True)¶ Generate the face inner product matrix
Parameters: - prop (numpy.ndarray) – material property (tensor properties are possible) at each cell center (nC, (1, 3, or 6))
- invProp (bool) – inverts the material property
- invMat (bool) – inverts the matrix
- doFast (bool) – do a faster implementation if available.
Returns: M, the inner product matrix (nF, nF)
Return type:
-
CylMesh.
getFaceInnerProductDeriv
(prop, doFast=True, invProp=False, invMat=False)¶ Parameters: - prop (numpy.ndarray) – material property (tensor properties are possible) at each cell center (nC, (1, 3, or 6))
- doFast – bool do a faster implementation if available.
- invProp (bool) – inverts the material property
- invMat (bool) – inverts the matrix
Returns: dMdmu(u), the derivative of the inner product matrix for a certain u
Return type:
-
CylMesh.
getInterpolationMat
(loc, locType='CC', zerosOutside=False)[source]¶ Produces interpolation matrix
Parameters: - loc (numpy.ndarray) – Location of points to interpolate to
- locType (str) –
What to interpolate locType can be:
'Ex' -> x-component of field defined on edges 'Ey' -> y-component of field defined on edges 'Ez' -> z-component of field defined on edges 'Fx' -> x-component of field defined on faces 'Fy' -> y-component of field defined on faces 'Fz' -> z-component of field defined on faces 'N' -> scalar field defined on nodes 'CC' -> scalar field defined on cell centers 'CCVx' -> x-component of vector field defined on cell centers 'CCVy' -> y-component of vector field defined on cell centers 'CCVz' -> z-component of vector field defined on cell centers
Returns: M, the interpolation matrix
Return type:
-
CylMesh.
getInterpolationMatCartMesh
(Mrect, locType='CC', locTypeTo=None)[source]¶ Takes a cartesian mesh and returns a projection to translate onto the cartesian grid.
Parameters: - Mrect (discretize.base.BaseMesh) – the mesh to interpolate on to
- locType (str) – grid location (‘CC’, ‘N’, ‘Ex’, ‘Ey’, ‘Ez’, ‘Fx’, ‘Fy’, ‘Fz’)
- locTypeTo (str) – grid location to interpolate to. If None, the same grid type as locType will be assumed
Returns: M, the interpolation matrix
Return type:
-
CylMesh.
getTensor
(key)¶ Returns a tensor list.
Parameters: key (str) – Which tensor (see below)
key can be:
'CC' -> scalar field defined on cell centers 'N' -> scalar field defined on nodes 'Fx' -> x-component of field defined on faces 'Fy' -> y-component of field defined on faces 'Fz' -> z-component of field defined on faces 'Ex' -> x-component of field defined on edges 'Ey' -> y-component of field defined on edges 'Ez' -> z-component of field defined on edges
Returns: list of the tensors that make up the mesh. Return type: list
-
CylMesh.
isInside
(pts, locType='N')¶ Determines if a set of points are inside a mesh.
Parameters: pts (numpy.ndarray) – Location of points to test Return type: numpy.ndarray Returns: inside, numpy array of booleans
-
CylMesh.
plotGrid
(*args, **kwargs)¶
-
CylMesh.
plotImage
(*args, **kwargs)¶
-
CylMesh.
projectEdgeVector
(eV)¶ Given a vector, eV, in cartesian coordinates, this will project it onto the mesh using the tangents
Parameters: eV (numpy.ndarray) – edge vector with shape (nE, dim) Return type: numpy.ndarray Returns: projected edge vector, (nE, )
-
CylMesh.
projectFaceVector
(fV)¶ Given a vector, fV, in cartesian coordinates, this will project it onto the mesh using the normals
Parameters: fV (numpy.ndarray) – face vector with shape (nF, dim) Return type: numpy.ndarray Returns: projected face vector, (nF, )
-
CylMesh.
r
(x, xType='CC', outType='CC', format='V')¶ r is a quick reshape command that will do the best it can at giving you what you want.
For example, you have a face variable, and you want the x component of it reshaped to a 3D matrix.
r can fulfil your dreams:
mesh.r(V, 'F', 'Fx', 'M') | | | | | | | { | | | How: 'M' or ['V'] for a matrix | | | (ndgrid style) or a vector (n x dim) | | | } | | { | | What you want: ['CC'], 'N', | | 'F', 'Fx', 'Fy', 'Fz', | | 'E', 'Ex', 'Ey', or 'Ez' | | } | { | What is it: ['CC'], 'N', | 'F', 'Fx', 'Fy', 'Fz', | 'E', 'Ex', 'Ey', or 'Ez' | } { The input: as a list or ndarray }
For example:
# Separates each component of the Ex grid into 3 matrices Xex, Yex, Zex = r(mesh.gridEx, 'Ex', 'Ex', 'M') # Given an edge vector, return just the x edges as a vector XedgeVector = r(edgeVector, 'E', 'Ex', 'V') # Separates each component of the edgeVector into 3 vectors eX, eY, eZ = r(edgeVector, 'E', 'E', 'V')
-
CylMesh.
save
(filename='mesh.json', verbose=False)¶ Save the mesh to json :param str file: filename for saving the casing properties :param str directory: working directory for saving the file
-
CylMesh.
serialize
(include_class=True, save_dynamic=False, **kwargs)¶ Serializes a HasProperties instance to dictionary
This uses the Property serializers to serialize all Property values to a JSON-compatible dictionary. Properties that are undefined are not included. If the HasProperties instance contains a reference to itself, a
properties.SelfReferenceError
will be raised.Parameters:
- include_class - If True (the default), the name of the class
will also be saved to the serialized dictionary under key
'__class__'
- save_dynamic - If True, dynamic properties are written to the serialized dict (default: False).
- Any other keyword arguments will be passed through to the Property serializers.
- include_class - If True (the default), the name of the class
will also be saved to the serialized dictionary under key
-
CylMesh.
setCellGradBC
(BC)¶ Function that sets the boundary conditions for cell-centred derivative operators.
Example
..code:: python
# Neumann in all directions BC = ‘neumann’
# 3D, Dirichlet in y Neumann else BC = [‘neumann’, ‘dirichlet’, ‘neumann’]
# 3D, Neumann in x on bottom of domain, Dirichlet else BC = [[‘neumann’, ‘dirichlet’], ‘dirichlet’, ‘dirichlet’]
-
CylMesh.
toVTK
(models=None)¶ Convert this mesh object to it’s proper
vtki
data object with the given model dictionary as the cell data of that dataset.Parameters: models (dict(numpy.ndarray)) – Name(‘s) and array(‘s). Match number of cells
-
CylMesh.
validate
()¶ Call all registered class validator methods
These are all methods decorated with
@properties.validator
. Validator methods are expected to raise a ValidationError if they fail.
-
CylMesh.
writeVTK
(fileName, models=None, directory='')¶ Makes and saves a VTK object from this mesh and given models
Parameters: