migflow.scontact
Model for Immersed Granular Flow: Fluid user interface
Contact: jonathan.lambrechts@uclouvain.be Webpage: www.migflow.be
MigFlow computes immersed granular flows using an unresolved FEM-DEM model. The continuous phase (mixture) is computed by solving averaged Navier-Stokes equations on unstructured meshes with the continuous finite element method.
- class migflow.scontact.ParticleProblem(dim)
Bases:
object
Creates the numerical structure containing all the physical particles that appear in the problem
Builds the particle problem structure.
- Parameters
dim (
int
) – Dimension of the problem- Raises
ValueError – If the dimension differs from 2 or 3
NameError – If C builder cannot be found
- __del__()
Deletes the particle solver structure.
- add_body(x, invertM, invertI)
Adds a body in the particle problem. — Only in 2D
- Parameters
x (
ndarray
) – array of the position of the bodyInvertM – inverse mass of the body
InvertI – inverse inertia of the body
- Returns
body id
- Return type
int
- add_boundary_periodic_entity(dim, transform)
Adds a periodic entity.
- Parameters
dim (
int
) – dimension of the entitytransform (
ndarray
) – array of the transformation to applied to the periodic entity
- add_boundary_periodic_point(x, tag)
Adds a boundary periodic point.
- Parameters
x (
ndarray
) – array of the coordinates of the pointtag (
int
) – entity tag
- add_boundary_periodic_segment(x0, x1, tag)
Adds a boundary periodic segment.
- Parameters
x0 (
ndarray
) – array of the coordinates of the first endpointx1 (
ndarray
) – array of the coordinates of the second endpointtag (
int
) – entity tag
- add_boundary_periodic_triangle(x0, x1, x2, tag)
Adds a boundary periodic triangle.
- Parameters
x0 (
ndarray
) – array of the coordinates of the first summitx1 (
ndarray
) – array of the coordinates of the second summitx2 (
ndarray
) – array of the coordinates of the third summittag (
int
) – tag of the periodic entity
- add_boundary_segment(x0, x1, tag, material='default')
Adds a boundary segment.
- Parameters
x0 (
ndarray
) – array of the coordinates of the first endpointx1 (
ndarray
) – array of the coordinates of the second endpointtag (
str
) – segment tagmaterial (
str
) – segment material
- add_boundary_triangle(x0, x1, x2, tag, material='default')
Adds a boundary triangle.
- Parameters
x0 (
ndarray
) – array of the coordinates of the first summitx1 (
ndarray
) – array of the coordinates of the second summitx2 (
ndarray
) – array of the coordinates of the third summittag (
str
) – triangle tagmaterial (
str
) – triangle material
- add_entities(xbody, invertM=0.0, invertI=0.0, particles_coordinates=None, particles_radius=None, particles_materials='default', segments_coordinates=None, segments_tags=None, segments_materials='default', triangles_coordinates=None, triangles_tags=None, triangles_materials='default', is_relative_position=True)
Adds entities (either particles, segments or triangles) to a body described by its position, inverse of mass and inverse of inertia.
- Parameters
xbody (np.ndarray) – position of the body
invertM (float, optional) – inverse mass of the body. Defaults to 0.
invertI (float, optional) – inverse inertia of the body. Defaults to 0.
particles_coordinates (np.ndarray, optional) – position of the particles. Defaults to None.
particles_radius (Union[np.ndarray, float], optional) – array of the particles radii. Defaults to None.
particles_materials (Union[str, List[str]], optional) – array of the particles material. Defaults to “default”.
segments_coordinates (np.ndarray, optional) – array of the segments coordinates as n_segments x 2 x dimension. Defaults to None.
segments_tags (Union[str, List[str]], optional) – array of the segments tags. Defaults to None.
segments_materials (Union[str, List[str]], optional) – array of the segments materials. Defaults to “default”.
triangles_coordinates (np.ndarray, optional) – array of the segments coordinates as n_segments x 3 x dimension. Defaults to None.
triangles_tags (str, optional) – array of the triangles tags. Defaults to None.
triangles_materials (Union[str, List[str]], optional) – array of the triangles materials. Defaults to “default”.
is_relative_position (bool, optional) – flag if the given position is the relative position towards the body or the absolute position. Defaults to True.
- Returns
body id
- Return type
int
- add_particle(x, r, m, tag='default', forced=False, inverse_inertia=None)
Adds a particle in the particle container.
- Parameters
x (
ndarray
) – array of the particle positionr (
float
) – particle radiusm (
float
) – particle masstag (
str
) – particle material
- Returns
body id
- Return type
int
- add_particle_body(positions, radii, masses, material='default', forced=False)
Adds a body in the particle problem with a list of constituting particles — Only in 2D
- Parameters
positions (
ndarray
) – positions of the particles constituting the bodyradii (
ndarray
) – radii of the particlesmasses (
ndarray
) – masses of the particlesmaterial (
str
) – material of the particlesforced (
bool
) – Flag indicating whether the body is forced or not
- Returns
body id
- Return type
int
- add_particle_to_body(x, r, body, material='default')
Adds a particle to a body.
- Parameters
x (
ndarray
) – array of the particle positon with respect to the body’s centre of massr (
float
) – Particle radiusbody (
int
) – number of the body to which the particle will be addedmaterial (
str
) – Particle material
- add_particles(x, r, m, tag='default', forced=False, inverse_inertia=None)
Adds particles in the particle container.
- Parameters
x (
ndarray
) – array of tuple to set the centre particle positionr (
ndarray
) – array of particle radiusm (
ndarray
) – array of particle masstag (
str
) – particle material
- Returns
id of the added bodies
- Return type
np.ndarray
- add_segment_to_body(x0, x1, body, tag, material='default')
Adds a segment to a body.
- Parameters
x0 (np.ndarray) – array of the coordinates of the first endpoint
x1 (np.ndarray) – array of the coordinates of the second endpoint
body (int) – body id
tag (str) – boundary id associated to the segment
material (str) – Segment material
- add_triangle_to_body(x0, x1, x2, body, tag, material='default')
Adds a segment to a body.
- Parameters
x0 (np.ndarray) – array of the coordinates of the first endpoint
x1 (np.ndarray) – array of the coordinates of the second endpoint
body (int) – body id
tag (str) – boundary id associated to the segment
material (str) – Segment material
- body_constraint()
Returns the invert mases of the bodies. :returns: Array of the body invert masses :rtype: np.ndarray
- body_invert_inertia()
Returns the invert inertias of the bodies. :returns: Array of the body invert inertias :rtype: np.ndarray
- body_invert_mass()
Returns the invert mases of the bodies. :returns: Array of the body invert masses :rtype: np.ndarray
- body_omega()
Returns the angular velocities of the bodies. :returns: Array of the body angular velocities :rtype: np.ndarray
- body_position()
Returns the positions of the bodies. :returns: Array of the body positions :rtype: np.ndarray
- body_theta()
Returns the array of the body angles. :returns: Array of the body angles :rtype: np.ndarray
- body_velocity()
Returns the velocities of the bodies. :returns: Array of the body velocities :rtype: np.ndarray
- bounding_box()
Returns the bounding box of the particle problem. :returns: position of the lower left corner, position of the upper right corner :rtype: Tuple[np.ndarray, np.ndarray]
- compute_relative_position(mass_center, particles_position, theta)
- compute_stress_tensor(nodes, radius)
Computes the stress tensor of the granular material :type nodes:
ndarray
:param nodes: Array of nodes at which to compute the stress tensor :param r: Radius within which the contact forces will be averaged- Returns
Array of stresses on the nodes of the grid
- Return type
np.ndarray
- contact_forces()
Returns the contact forces on each particle. :returns: Array of the particle contact forces :rtype: np.ndarray
- contacts()
Gives the contacts between the bodies. :returns: Array of contact forces between the particles :rtype: np.ndarray
- delassus()
Returns the delassus operators of the particles. :returns: Array of the particle delassus operators :rtype: np.ndarray
- dim()
Returns the dimension of the particle problem. :returns: Dimension of the problem :rtype: int
- get_boundary_forces(tag='default')
Returns the forces acting on a boundary because of the contacts in the global basis. :type tag:
str
:param tag: Name of the boundary- Returns
Array of forces on boundary named tag
- Return type
np.ndarray
- get_material_id(material)
Returns the number associated to a string material. :type material:
str
:param material: string material.- Returns
number associated to the string material
- Return type
int
- get_tag_id(tag='default')
Returns the number associated to a string tag. :type tag:
str
:param tag: string tag.- Returns
number associated to the string tag
- Return type
int
- id()
Returns the ids of the particles. :returns: Array of the particle ids :rtype: np.ndarray
- init_body(reference_position, relative_position, particle_masses, particle_radii, material='default')
- Return type
int
- iterate(dt, forces, torques=None, body_forces=None, tol=1e-08, force_motion=0)
- Computes iteratively the set of velocities that respect the non-interpenetration constrain.
Returns 1 if the computation converged, 0 otherwise.
- Parameters
dt (
float
) – Numerical time stepforces (
ndarray
) – List of vectors containing the forces applied on the particlesbody_forces (
Optional
[ndarray
]) – List of forces applied to the bodiestol (
float
) – Optional argument defining the interpenetration tolerance to stop the NLGS iterations of the NSCDforce_motion (
int
) – Optional argument allowing to move the grains if convergence has not been reached when set to 1
- load_msh_boundaries(filename, tags=None, shift=None, material='default')
Loads the numerical domain and set the physical boundaries the particles cannot go through.
- Parameters
filename (
str
) – name of the msh file. If None the current gmsh model is loaded without the need to a msh file.tags (
Optional
[str
]) – tags of physical boundary defined in the msh fileshift (
Optional
[ndarray
]) – optional argument to shift the numerical domainmaterial (
str
) – material of the boundary
- material()
Returns the materials of the particles. :returns: Array of the particle materials :rtype: np.ndarray
- mu()
Returns the array of the friction coefficients between materials. :returns: Array of the friction coefficients between materials. :rtype: np.ndarray
- n_bodies()
Returns the number of bodies. :returns: number of bodies :rtype: int
- n_particles()
Returns the number of particles. :returns: number of particles :rtype: int
- omega()
Returns the omega of the particles. :returns: Array of the particle omega :rtype: np.ndarray
- parse_type(tname)
- particle_body()
Returns the bodies associated to each particle
- Returns
body id
- Return type
np.ndarray
- particle_relative_position()
Returns the relative positions of the particles
- Returns
Array of the relative positions
- Return type
np.ndarray
- periodic_entity()
Returns the array of periodic entities. :returns: Array of periodic entities :rtype: np.ndarray
- periodic_points()
Returns the array of periodic points. :returns: Array of periodic points :rtype: np.ndarray
- periodic_segments()
Returns the array of periodic segments.
- Return type
ndarray
- periodic_triangles()
Returns the array of periodic triangles. :returns: Array of periodic segments :rtype: np.ndarray
- position()
Returns the positions of the particles. :returns: Array of the particle positions :rtype: np.ndarray
- r()
Returns the raddi of the particles. :returns: Array of the particle radii :rtype: np.ndarray
- read_mig(odir, t=None, iteration=None)
Reads output files. :type odir:
str
:param odir: Directory in which to read the file :type t:Optional
[float
] :param t: Time at which to read the file :param iteration; iteration to read: :param -1 for last iteration:
- remove_bodies_flag(flag)
Removes particles based on given flag. :type flag:
ndarray
:param flag: flag with the same size as the number of bodies, with 0 for bodies to be removed and 1 for bodies to be kept
- restore_state()
Restores the saved state of the particle problem.
- save_state()
Saves the current state of the particle problem.
- segments()
Returns the array of boundary segments. :returns: Array of boundary segments :rtype: np.ndarray
- set_boundary_velocity(tag, v)
Sets the velocity of a boundary to a given value. :type tag:
str
:param tag: Name of the boundary :type v:ndarray
:param v: Velocity vector to be given to the boundary
- set_compute_contact_forces(compute_contact_forces=1)
Enables the computation of the contact forces on each particle if 1 and disables it if 0.
- set_fixed_contact_geometry(flag)
Sets the management of the contact geometry, either fixed (default) or changing each time the contact is solved.
- Parameters
flag (
bool
) – True for fixed False for changing.
- set_friction_coefficient(mu, mat0='default', mat1='default')
Sets the friction coefficient between two materials.
- Parameters
mu (
float
) – Value of the friction coefficientmat0 (
str
) – First materialmat1 (
str
) – Second material
- set_predict_basis(flag)
Set the management of the contact geometry, if true, a contact basis estimation is used in the contact resolution
- Parameters
flag (bool) – True for the estimation
- set_use_queue(use_queue=1)
Enables the use of the queue if 1 and disables it if 0.
- triangles()
Returns the array of boundary triangles. :returns: Array of boundary triangles :rtype: np.ndarray
- velocity()
Returns the velocities of the particles. :returns: Array of the particle velocities :rtype: np.ndarray
- volume()
Returns the volumes of the particles. :returns: Array of the particle volumes :rtype: np.ndarray
- write_mig(output_dir, t, extra_fields={}, write_contacts=True)
Writes output files for post-visualization. :type output_dir:
str
:param output_dir: Directory in which to write the file (expected to end in “.mig”) :type t:float
:param t: Time at which the simulation is :type extra_fields:Dict
[str
,ndarray
] :param extra_fields: extra field as a dictionnary {name : values}
- migflow.scontact.fastdeposit_2d(boundary, r, svgfilename=None)
Fast 2d deposit of particles in a simple domain. :type boundary:
ndarray
:param boundary: points describing the boundary (single loop counter-clockwise oriented, no holes). :type boundary: np.ndarray [n,2] :type r:ndarray
:param r: radius of the particles. :type r: np.ndarray [m] :type svgfilename:Optional
[str
] :param svgfilename: if not None, a svg file is generated illustrating the algorithm- Returns
positions of the particles.
- Return type
np.ndarray [m,2]