Mesh Object

class pylagrit.MO(name, parent)

Mesh object class

add_element_attribute(attname, keyword=None, vtype='VDOUBLE', rank='scalar', interpolate='linear', persistence='permanent', ioflag='', value=0.0)

Add a list of attributes to elements

Parameters:
  • attnames (str) – Attribute name to add
  • keyword – Keyword used by lagrit for specific attributes
  • vtype – Type of variable {‘VDOUBLE’,’VINT’,…}
add_node_attribute(attname, keyword=None, vtype='VDOUBLE', rank='scalar', interpolate='linear', persistence='permanent', ioflag='', value=0.0)

Add a list of attributes to nodes

Parameters:
  • attnames (str) – Attribute name to add
  • keyword – Keyword used by lagrit for specific attributes
  • vtype – Type of variable {‘VDOUBLE’,’VINT’,…}
addatt(attname, keyword=None, vtype='VDOUBLE', rank='scalar', length='nnodes', interpolate='linear', persistence='permanent', ioflag='', value=0.0)

Add a list of attributes

Parameters:
  • attnames (str) – Attribute name to add
  • keyword – Keyword used by lagrit for specific attributes
  • vtype – Type of variable {‘VDOUBLE’,’VINT’,…}
addatt_voronoi_volume(name='voronoi_volume')

Add voronoi volume attribute to mesh object

Parameters:name (str) – name of attribute in LaGriT
compute_distance(mo, option='distance_field', attname='dfield')

Compute distance from one mesh object to another

Parameters:
  • mo (LaGriT mesh object) – Mesh object to compute distance to base mesh from
  • option (str) – The type of distance field calculation. Available choices are ‘distance_field’ and ‘signed_distance_field’.
  • attname (str) – The name of the attribute to be created in the base mesh.

Returns: New attribute in base mesh object

Example: from pylagrit import PyLaGriT #create source mesh npts = (1,91,1) mins = (3.,0.,0.) maxs = (3.,270.,0.) src_mo = lg.create() src_mo.createpts_rtz(npts,mins,maxs,connect=False)

#create sink mesh snk_mo = lg.create() snk_mo.createpts_xyz([30,30,1],[-5.,-5.,-5.],[5.,5.,5.],connect=False)

#compute distance and store in sink mesh attribute ‘dfield’ snk_mo.compute_distance(src_mo) snk_mo.dump(‘comptest.gmv’)

compute_extrapolate(surf_mo, dir='zpos', attname='zic')
Given a 3D mesh and a 2D surface, this command will extrapolate a scalar
value from that surface onto every point of the mesh.
Parameters:
  • surf_mo (LaGriT mesh object) – Surface mesh object to extrapolate from
  • dir – The direction values are extrapolated from. Choices are one

of: ‘zpos’, ‘zneg’, ‘ypos’, ‘yneg’, ‘xpos’, ‘xneg’ :type dir: str

Parameters:attname – The name of the attribute in the surface mesh to be

extrapolated :type attname: str

Returns: New attribute in base mesh object

Example: from pylagrit import PyLaGriT #create surface mesh p1 = (-1.,-1.,-1.) p2 = (301.,-1.,-1.) p3 = (301.,301.,-1.) p4 = (-1.,301.,-1.) pts = [p1,p2,p3,p4] nnodes = (30,30,1) surf = lg.create_qua() surf.quadxy(nnodes,pts)

#make surface mesh interesting surf.math(‘sin’,’zic’,cmosrc=surf,attsrc=’xic’) surf.math(‘multiply’,’zic’,value=5.0,cmosrc=surf,attsrc=’zic’) surf.perturb(0.,0.,1.) surf.math(‘add’,’zic’,value=60.0,cmosrc=surf,attsrc=’zic’)

#create base mesh hex = lg.create_hex() hex.createpts_brick_xyz([30,30,20],[0.,0.,0.],[300.,300.,50.]) hex.resetpts_itp()

#extrapolate z values from surface mesh to base mesh hex.compute_extrapolate(surf) hex.dump(‘extrapolated.gmv’)

connect(option1='delaunay', option2=None, stride=None, big_tet_coords=[])

Connect the nodes into a Delaunay tetrahedral or triangle grid.

Parameters:
  • option1 – type of connect: delaunay, noadd, or check_interface
  • option2 – type of connect: noadd, or check_interface
  • stride (tuple(int)) – tuple of (first, last, stride) of points
connect_check_interface()

Connect the nodes into a Delaunay tetrahedral or triangle grid exhaustively checking that no edges of the mesh cross a material boundary.

connect_delaunay(option2=None, stride=None, big_tet_coords=[])

Connect the nodes into a Delaunay tetrahedral or triangle grid without adding nodes.

connect_noadd()

Connect the nodes into a Delaunay tetrahedral or triangle grid without adding nodes.

copy(name=None)

Copy mesh object

copyatt(attname_src, attname_sink=None, mo_src=None)

Add a list of attributes

Parameters:
  • attname_src (str) – Name of attribute to copy
  • attname_sink (str) – Name of sink attribute
  • mo_src (PyLaGriT Mesh Object) – Name of source mesh object
copypts(elem_type='tet', name=None)

Copy points from mesh object to new mesh object

Parameters:
  • name (str) – Name to use within lagrit for the created mesh object
  • mesh_type (str) – Mesh type for new mesh
Returns:

mesh object

create_boundary_facesets(stacked_layers=False, base_name=None, reorder=False, external=True)

Creates facesets for each boundary and writes associated avs faceset file :arg base_name: base name of faceset files :type base_name: str :arg stacked_layers: if mesh is created by stack_layers, user layertyp attr to determine top and bottom :type stacked_layers: bool :arg reorder_on_meds: reorder nodes on cell medians, usually needed for exodus file :type reorder_on_meds: bool :returns: Dictionary of facesets

createpts(crd, npts, mins, maxs, vc_switch=(1, 1, 1), rz_switch=(1, 1, 1), rz_value=(1, 1, 1), connect=False)

Create and Connect Points

Parameters:
  • crd (str) – Coordinate type of either ‘xyz’ (cartesian coordinates), ‘rtz’ (cylindrical coordinates), or ‘rtp’ (spherical coordinates).
  • npts (tuple(int)) – The number of points to create in line
  • mins (tuple(int, int, int)) – The starting value for each dimension.
  • maxs (tuple(int, int, int)) – The ending value for each dimension.
  • vc_switch (tuple(int, int, int)) – Determines if nodes represent vertices (1) or cell centers (0).
  • rz_switch (tuple(int, int, int)) – Determines true or false (1 or 0) for using ratio zoning values.
createpts_brick(crd, npts, mins, maxs, vc_switch=(1, 1, 1), rz_switch=(1, 1, 1), rz_vls=(1, 1, 1))

Create and Connect Points

Creates a grid of points in the mesh object and connects them.

Parameters:
  • crd (str) – Coordinate type of either ‘xyz’ (cartesian coordinates), ‘rtz’ (cylindrical coordinates), or ‘rtp’ (spherical coordinates).
  • npts (tuple(int, int, int)) – The number of points to create in each dimension.
  • mins (tuple(int, int, int)) – The starting value for each dimension.
  • maxs (tuple(int, int, int)) – The ending value for each dimension.
  • vc_switch (tuple(int, int, int)) – Determines if nodes represent vertices (1) or cell centers (0).
  • rz_switch (tuple(int, int, int)) – Determines true or false (1 or 0) for using ratio zmoning values.
  • rz_vls (tuple(int, int, int)) – Ratio zoning values. Each point will be multiplied by a scale of the value for that dimension.
createpts_brick_rtp(npts, mins, maxs, vc_switch=(1, 1, 1), rz_switch=(1, 1, 1), rz_vls=(1, 1, 1))

Create and connect spherical coordinates.

createpts_brick_rtz(npts, mins, maxs, vc_switch=(1, 1, 1), rz_switch=(1, 1, 1), rz_vls=(1, 1, 1))

Create and connect cylindrical coordinate points.

createpts_brick_xyz(npts, mins, maxs, vc_switch=(1, 1, 1), rz_switch=(1, 1, 1), rz_vls=(1, 1, 1))

Create and connect Cartesian coordinate points.

createpts_dxyz(dxyz, mins, maxs, clip='under', hard_bound='min', rz_switch=(1, 1, 1), rz_value=(1, 1, 1), connect=True)

Create and Connect Points to create an orthogonal hexahedral mesh. The vertex spacing is based on dxyz and the mins and maxs specified. mins (default, see hard_bound option) or maxs will be adhered to, while maxs (default) or mins will be modified based on the clip option to be truncated at the nearest value ‘under’ (default) or ‘over’ the range maxs-mins. clip and hard_bound options can be mixed by specifying tuples (see description below).

Parameters:
  • dxyz (tuple(float,float,float)) – The spacing between points in x, y, and z directions
  • mins (tuple(float,float,float)) – The starting value for each dimension.
  • maxs (tuple(float,float,float)) – The ending value for each dimension.
  • clip (string or tuple(string,string,string)) – How to handle bounds if range does not divide by dxyz, either clip ‘under’ or ‘over’ range
  • hard_bound (string or tuple(string,string,string)) – Whether to use the “min” or “max” as the hard constraint on dimension
  • rz_switch (tuple(int, int, int)) – Determines true or false (1 or 0) for using ratio zoning values.
  • connect (boolean) – Whether or not to connect points
createpts_line(npts, mins, maxs, vc_switch=(1, 1, 1), rz_switch=(1, 1, 1))

Create and Connect Points in a line

Parameters:
  • npts (int) – The number of points to create in line
  • mins (tuple(int, int, int)) – The starting value for each dimension.
  • maxs (tuple(int, int, int)) – The ending value for each dimension.
  • vc_switch (tuple(int, int, int)) – Determines if nodes represent vertices (1) or cell centers (0).
  • rz_switch (tuple(int, int, int)) – Determines true or false (1 or 0) for using ratio zoning values.
createpts_median()
createpts_rtp(npts, mins, maxs, vc_switch=(1, 1, 1), rz_switch=(1, 1, 1), rz_value=(1, 1, 1), connect=True)
createpts_rtz(npts, mins, maxs, vc_switch=(1, 1, 1), rz_switch=(1, 1, 1), rz_value=(1, 1, 1), connect=True)
createpts_xyz(npts, mins, maxs, vc_switch=(1, 1, 1), rz_switch=(1, 1, 1), rz_value=(1, 1, 1), connect=True)
delatt(attnames, force=True)

Delete a list of attributes

Parameters:
  • attnames (str or lst(str)) – Attribute names to delete
  • force (bool) – If true, delete even if the attribute permanent persistance
delete()
dump(filename=None, format=None, *args)
dump_ats_xml(filename, meshfilename, matnames={}, facenames={})

Write ats style xml file with regions :param filename: Name of xml to write :type filename: string :param meshfilename: Name of exodus file to use in xml :type meshfilename: string :param matnames: Dictionary of region names keyed by exodus material number :type matnames: dict :param facenames: Dictionary of faceset names keyed by exodus faceset number :type facenames: dict

dump_avs2(filename, points=True, elements=True, node_attr=True, element_attr=True)

Dump avs file

Parameters:
  • filename (str) – Name of avs file
  • points (bool) – Output point coordinates
  • elements (bool) – Output connectivity
  • node_attr (bool) – Output node attributes
  • element_attr (bool) – Output element attributes
dump_exo(filename, psets=False, eltsets=False, facesets=[])

Dump exo file

Parameters:
  • filename (str) – Name of exo file
  • psets (bool) – Boolean indicating that exodus will only include psets
  • eltsets (bool) – Boolean indicating that exodus will only include element sets
  • facesets (lst(FaceSet)) – Array of FaceSet objects
Example:
>>> from pylagrit import PyLaGriT
>>> l = PyLaGriT()
>>> m = l.create()
>>> m.createpts_xyz((3,3,3),(0.,0.,0.),(1.,1.,1.),rz_switch=[1,1,1],connect=True)
>>> m.status ()
>>> m.status (brief=True)
>>> fs = m.create_boundary_facesets(base_name='faceset_bounds')
>>> m.dump_exo('cube.exo',facesets=fs.values())
dump_fehm(filename, *args)
dump_gmv(filename, format='binary')
dump_lg(filename, format='binary')
dump_pflotran(filename_root, nofilter_zero=False)

Dump PFLOTRAN UGE file

Parameters:
  • filename_root (str) – root name of UGE file
  • nofilter_zero (boolean) – Set to true to write zero coefficients to file
Example:
>>> from pylagrit import PyLaGriT
>>> l = PyLaGriT()
>>> m = l.create()
>>> m.createpts_xyz((3,3,3),(0.,0.,0.),(1.,1.,1.),rz_switch=[1,1,1],connect=True)
>>> m.status ()
>>> m.status (brief=True)
>>> m.dump_pflotran('test_pflotran_dump')
dump_pset(filerootname, zonetype='zone', pset=[])

Dump zone file of psets :arg filerootname: rootname of files to create, pset name will be added to name :type filerootname: string :arg zonetype: Type of zone file to dump, ‘zone’ or ‘zonn’ :type zonetype: string :arg pset: list of psets to dump, all psets dumped if empty list :type pset: list[strings]

dump_zone_imt(filename, imt_value)
dump_zone_outside(filename, keepatt=False, keepatt_median=False, keepatt_voronoi=False)
elem_type
eltset_attribute(attribute_name, attribute_value, boolstr='eq', name=None)
eltset_bool(eset_list, boolstr='union', name=None)

Create element set from boolean operation of set of element sets

Parameters:
  • eset_list (lst(PyLaGriT element set)) – List of elements to perform boolean operation on
  • boolstr (str) – type of boolean operation to perform on element sets, one of [union,inter,not]
  • name (str) – The name to be assigned to the EltSet within LaGriT
Returns:

PyLaGriT element set object

eltset_inter(eset_list, name=None)
eltset_not(eset_list, name=None)
eltset_object(mo, name=None)

Create element set from the intersecting elements with another mesh object

eltset_region(region, name=None)
eltset_union(eset_list, name=None)
eltset_write(filename_root, eset_name=None, ascii=True)

Write element set(s) to a file in ascii or binary format

Parameters:
  • filename_root (str) – root name of file
  • eset_name (EltSet object) – name of eltset to write; if blank, all eltsets in mesh object are written
  • ascii – Switch to indicate ascii [True] or binary [False]
Example:
>>> from pylagrit import PyLaGriT
>>> import numpy as np
>>> import sys
>>> 
>>> lg = PyLaGriT()
>>> 
>>> dxyz = np.array([0.1,0.25,0.25])
>>> mins = np.array([0.,0.,0.])
>>> maxs = np.array([1.,1.,1.])
>>> mqua = lg.createpts_dxyz(dxyz,mins,maxs,'quad',hard_bound=('min','max','min'),connect=True)
>>> 
>>> example_pset1 = mqua.pset_geom_xyz(mins,maxs-(maxs-mins)/2)
>>> example_eset1 = example_pset1.eltset()
>>> example_pset2 = mqua.pset_geom_xyz(mins+maxs/2,maxs)
>>> example_eset2 = example_pset2.eltset()
>>> # to write one specific eltset
>>> mqua.eltset_write('test_specific',eset_name=example_eset1)
>>> # to write all eltsets
>>> mqua.eltset_write('test_all')
extract_surfmesh(name=None, stride=[1, 0, 0], reorder=False, resetpts_itp=True, external=False)
extrude(offset, offset_type='const', return_type='volume', direction=[], name=None)

Extrude mesh object to new mesh object This command takes the current mesh object (topologically 1d or 2d mesh (a line, a set of line segments, or a planar or non-planar surface)) and extrudes it into three dimensions along either the normal to the curve or surface (default), along a user defined vector, or to a set of points that the user has specified. If the extrusion was along the normal of the surface or along a user defined vector, the command can optionally find the external surface of the volume created and return that to the user. Refer to http://lagrit.lanl.gov/docs/commands/extrude.html for more details on arguments.

Parameters:
  • name (str) – Name to use within lagrit for the created mesh object
  • offset (float) – Distance to extrude
  • offset_type (str) – either const or min (interp will be handled in the PSET class in the future)
  • return_type (str) – either volume for entire mesh or bubble for just the external surface
  • direction (lst[float,float,float]) – Direction to extrude in, defaults to normal of the object
Returns:

mesh object

filter(stride=[1, 0, 0], tolerance=None, boolean=None, attribute=None)
gmv(exe=None, filename=None)
grid2grid(ioption, name=None)

Convert a mesh with one element type to a mesh with another

Parameters:
  • ioption – type of conversion: quadtotri2 quad to 2 triangles, no new points. prismtotet3 prism to 3 tets, no new points. quadtotri4 quad to 4 triangles, with one new point. pyrtotet4 pyramid to 4 tets, with one new point. hextotet5 hex to 5 tets, no new points. hextotet6 hex to 6 tets, no new points. prismtotet14 prism to 14 tets, four new points (1 + 3 faces). prismtotet18 prism to 18 tets, six new points (1 + 5 faces). hextotet24 hex to 24 tets, seven new points (1 + 6 faces). tree_to_fe quadtree or octree grid to grid with no parent-type elements.
  • name (str) – Internal Lagrit name of new mesh object, automatically created if None

Returns MO object

grid2grid_hextotet24(name=None)

Hex to 24 tets, seven new points (1 + 6 faces) :arg name: Internal Lagrit name of new mesh object, automatically created if None :type name: str

Returns MO object

grid2grid_hextotet5(name=None)

Hex to 5 tets, no new points :arg name: Internal Lagrit name of new mesh object, automatically created if None :type name: str

Returns MO object

grid2grid_hextotet6(name=None)

Hex to 6 tets, no new points :arg name: Internal Lagrit name of new mesh object, automatically created if None :type name: str

Returns MO object

grid2grid_prismtotet14(name=None)

Prism to 14 tets, four new points (1 + 3 faces) :arg name: Internal Lagrit name of new mesh object, automatically created if None :type name: str

Returns MO object

grid2grid_prismtotet18(name=None)

Prism to 18 tets, four new points (1 + 3 faces) :arg name: Internal Lagrit name of new mesh object, automatically created if None :type name: str

Returns MO object

grid2grid_prismtotet3(name=None)

Quad to 2 triangles, no new points. Prism to 3 tets, no new points. :arg name: Internal Lagrit name of new mesh object, automatically created if None :type name: str

Returns MO object

grid2grid_pyrtotet4(name=None)

Pyramid to 4 tets, with one new point :arg name: Internal Lagrit name of new mesh object, automatically created if None :type name: str

Returns MO object

grid2grid_quadtotri2(name=None)

Quad to 2 triangles, no new points. :arg name: Internal Lagrit name of new mesh object, automatically created if None :type name: str

Returns MO object

grid2grid_quadtotri4(name=None)

Quad to 4 triangles, with one new point :arg name: Internal Lagrit name of new mesh object, automatically created if None :type name: str

Returns MO object

grid2grid_tree_to_fe(name=None)

Quadtree or octree grid to grid with no parent-type elements. :arg name: Internal Lagrit name of new mesh object, automatically created if None :type name: str

Returns MO object

information()

Returns a formatted dictionary with mesh information.

Information is that found in cmo/status/MO

interpolate(method, attsink, cmosrc, attsrc, stride=[1, 0, 0], tie_option=None, flag_option=None, keep_option=None, interp_function=None)

Interpolate values from attribute attsrc from mesh object cmosrc to current mesh object

interpolate_continuous(attsink, cmosrc, attsrc, stride=[1, 0, 0], interp_function=None, nearest=None)
interpolate_default(attsink, cmosrc, attsrc, stride=[1, 0, 0], tie_option='tiemax', flag_option='plus1', keep_option='delatt', interp_function=None)
interpolate_map(attsink, cmosrc, attsrc, stride=[1, 0, 0], tie_option=None, flag_option=None, keep_option=None, interp_function=None)
interpolate_voronoi(attsink, cmosrc, attsrc, stride=[1, 0, 0], interp_function=None)
intersect_elements(mo, attr_name='attr00')

This command takes two meshes and creates an element-based attribute in mesh1 that contains the number of elements in mesh2 that intersected the respective element in mesh1. We define intersection as two elements sharing any common point.

Parameters:
  • mo (PyLaGriT mesh object) – Mesh object to intersect with current mesh object to determine where to refine
  • attr_name (str) – Name to give created attribute
Returns:

attr_name

list(attname=None, stride=[1, 0, 0], pset=None)
massage(bisection_len, merge_len, toldamage, tolroughness=None, stride=None, nosmooth=False, norecon=False, strictmergelength=False, checkaxy=False, semiexclusive=False, ignoremats=False, lite=False)

MASSAGE creates, annihilates, and moves nodes and swaps connections in a 2D or 3D mesh in order to improve element aspect ratios and establish user-desired edge lengths.

The actions of MASSAGE are controlled by values of these four parameters:

bisection_length - edge length that will trigger bisection. merge_length - edge length that will trigger merging. toldamage - maximum grid deformation of interfaces and external boundaries

allowed in a single merge, smooth or reconnection event.
tolroughness - (for 2D surface grids only) measure of grid roughness
(deviation from average surface normal) that triggers refinement.

The final, optional keywork argument(s) can be one or more of nosmooth, norecon, lite, ignoremats, strictmergelength, checkaxy, semiexclusive, and exclusive.

Specifying nosmooth will turn off the ‘smooth’ step by skipping the call to SGD. Specifying norecon will turn off all ‘recon’ steps. If lite is specified, only one iteration of the merging/reconnection/smoothing loop is executed, and a reconnection after edge refinement is omitted. This is suitable for applications, such as Gradient Weighted Moving Finite Elements, where MASSAGE is called repeatedly.

The optional argument ignoremats causes MASSAGE to process the multimaterial mesh in a single material mode; it ignores the material interfaces.

The optional argument strictmergelength forces strict interpretation of merge_length so that there is no merging along the edges of flat elements. This is important if ignoremats is specified to avoid losing the interfaces.

If checkaxy is given, then we insure that for 2D meshes, the output mesh will have positive xy-projected triangle areas, provided that the input mesh had them in the first place.

If exclusive is given, then edge refinement operations will only be performed on edges whose endpoints are both in the PSET that MASSAGE is working on. (As usual, new nodes created by refinement are added to the PSET so that MASSAGE can refine edges recursively.) The default behavior is ‘inclusive’, where only ONE edge endpoint has to belong to the PSET for the edge to be eligible for refinement.

If semiexclusive is given, refinement will only be triggered by edges with both endpoints in the PSET, but some edges with less than two endpoints in the PSET might be refined as part of a ‘Rivara chain’ triggered by the refinement of an edge with both endpoints in the PSET. This represents an intermediate case between ‘inclusive’ and exclusive

massage2(filename, min_scale, bisection_len, merge_len, toldamage, tolroughness=None, stride=None, nosmooth=False, norecon=False, strictmergelength=False, checkaxy=False, semiexclusive=False, ignoremats=False, lite=False)

MASSAGE2 iteratively calls MASSAGE to refine adaptively according to a gradient field. Thus, the bisection_length option must be a field.

file_name is a file which contains a set of LaGriT commands that calculates the gradient field based on the distance field. In other words, the gradient field is a function of the distance field. It is necessary to have this file when using this routine, as the field must be updated after each refinement iteration.

Use this function in conjunction with PyLaGriT.define(**kwargs) for best results.

See MASSAGE for other arguments.

math(operation, attsink, value=None, stride=[1, 0, 0], cmosrc=None, attsrc=None)
maxs
minmax(attname=None, stride=[1, 0, 0])
minmax_xyz(stride=[1, 0, 0], verbose=True)
mins
mregion(boolstr, name=None)

Create mregion using boolean string

Parameters:
  • boolstr (str) – String of boolean operations
  • name (string) – Internal lagrit name for mesh object
Returns:

MRegion

ndim_geo
ndim_topo
nelems
nnodes
paraview(exe=None, filename=None)
perturb(xfactor, yfactor, zfactor, stride=(1, 0, 0))

This command moves node coordinates in the following manner.

Three pairs of random numbers between 0 and 1 are generated. These pairs refer to the x, y and z coordinates of the nodes respectively. The first random number of each pair is multiplied by the factor given in the command. The second random number is used to determine if the calculated offset is to be added or subtracted from the coordinate.

printatt(attname=None, stride=[1, 0, 0], pset=None, eltset=None, ptype='value')
pset_attribute(attribute, value, comparison='eq', stride=(1, 0, 0), name=None)

Define PSet by attribute

Parameters:
  • attribute (str) – Nodes defined by attribute ID.
  • value (integer) – attribute ID value.
  • comparison (can use default without specifiy anything, or list[lt|le|gt|ge|eq|ne]) – attribute comparison, default is eq.
  • stride (list[int, int, int]) – Nodes defined by ifirst, ilast, and istride.
  • name (str) – The name to be assigned to the PSet created.

Returns: PSet object

pset_bool(pset_list, boolean='union', name=None)

Return PSet from boolean operation on list of psets

Defines and returns a PSet from points that are not inside the PSet, ps.

pset_geom(mins, maxs, ctr=(0, 0, 0), geom='xyz', stride=(1, 0, 0), name=None)

Define PSet by Geometry

Selects points from geometry specified by string geom and returns a PSet.

Parameters:
  • mins (tuple(int, int, int)) – Coordinate of one of the shape’s defining points. xyz (Cartesian): (x1, y1, z1); rtz (Cylindrical): (radius1, theta1, z1); rtp (Spherical): (radius1, theta1, phi1);
  • maxs (tuple(int, int, int)) – Coordinate of one of the shape’s defining points. xyz (Cartesian): (x2, y2, z2); rtz (Cylindrical): (radius2, theta2, z2); rtp (Spherical): (radius2, theta2, phi2);
  • ctr (tuple(int, int, int)) – Coordinate of the relative center.
  • geom (str) – Type of geometric shape: ‘xyz’ (spherical), ‘rtz’ (cylindrical), ‘rtp’ (spherical)
  • stride (list[int, int, int]) – Nodes defined by ifirst, ilast, and istride.
  • name (str) – The name to be assigned to the PSet created.

Returns: PSet object

pset_geom_rtp(mins, maxs, ctr=(0, 0, 0), stride=(1, 0, 0), name=None)

Forms a pset of nodes within the sphere, sperical shell or sperical section given by radius1 to radius2, and angles theta1 to theta2 (0 - 180) and angles phi1 to phi2 (0 - 360). Refer to http://lagrit.lanl.gov/docs/conventions.html for an explanation of angles

Parameters:
  • mins (tuple(int, int, int)) – Defines radius1, theta1, and phi1.
  • maxs (tuple(int, int, int)) – Defines radius2, theta2, and phi2.
  • stride (list[int, int, int]) – Nodes defined by ifirst, ilast, and istride.
  • name (str) – The name to be assigned to the PSet created.
  • ctr (tuple(int, int, int)) – Coordinate of the relative center.
  • stride – Nodes defined by ifirst, ilast, and istride.
  • name – The name to be assigned to the PSet created.

Returns: PSet object

pset_geom_rtz(mins, maxs, ctr=(0, 0, 0), stride=(1, 0, 0), name=None)

Forms a pset of nodes within the cylinder or cylindrical shell section given by radius1 to radius2, and angles theta1 to theta2 and height z1 to z2. Refer to http://lagrit.lanl.gov/docs/conventions.html for an explanation of angles

Parameters:
  • mins (tuple(int, int, int)) – Defines radius1, theta1, and z1.
  • maxs (tuple(int, int, int)) – Defines radius2, theta2, and z2.
  • stride (list[int, int, int]) – Nodes defined by ifirst, ilast, and istride.
  • name (str) – The name to be assigned to the PSet created.
  • ctr (tuple(int, int, int)) – Coordinate of the relative center.
  • stride – Nodes defined by ifirst, ilast, and istride.
  • name – The name to be assigned to the PSet created.

Returns: PSet object

pset_geom_xyz(mins, maxs, ctr=(0, 0, 0), stride=(1, 0, 0), name=None)

Define PSet by Tetrahedral Geometry

Selects points from a Tetrahedral region.

Parameters:
  • mins (tuple(int, int, int)) – Coordinate point of 1 of the tetrahedral’s corners.
  • maxs (tuple(int, int, int)) – Coordinate point of 1 of the tetrahedral’s corners.
  • ctr (tuple(int, int, int)) – Coordinate of the relative center.
  • stride (list[int, int, int]) – Nodes defined by ifirst, ilast, and istride.
  • name (str) – The name to be assigned to the PSet created.

Returns: PSet object

pset_inter(pset_list, name=None)
pset_not(pset_list, name=None)
pset_region(region, stride=(1, 0, 0), name=None)

Define PSet by region

Parameters:
  • region – region to create pset
  • stride (list[int, int, int]) – Nodes defined by ifirst, ilast, and istride.
  • name (str) – The name to be assigned to the PSet created.

Returns: PSet object

pset_surface(surface, stride=(1, 0, 0), name=None)

Define PSet by surface

Parameters:
  • surface – surface to create pset
  • stride (list[int, int, int]) – Nodes defined by ifirst, ilast, and istride.
  • name (str) – The name to be assigned to the PSet created.

Returns: PSet object

pset_union(pset_list, name=None)
quadxy(nnodes, pts, connect=True)

Define an arbitrary, logical quad of points in 3D space with nnodes(x,y,z) nodes. By default, the nodes will be connected.

Parameters:
  • nnodes (tuple(int, int, int)) – The number of nodes to create in each dimension. One value must == 1 and the other two must be > 1.
  • pts (list of four 3-tuples (float)) – The four corners of the quad surface, defined in counter clockwise order (the normal to the quad points is defined using the right hand rule and the order of the points).
  • connect (bool) – connect points
Example:
>>> #To use pylagrit, import the module.
>>> import pylagrit
>>> #Start the lagrit session.
>>> lg = pylagrit.PyLaGriT()
>>> #Create a mesh object.
>>> qua = lg.create_qua()
>>> #Define 4 points in correct order
>>> p1 = (0.0,200.0,-400.0)
>>> p2 = (0.0,-200.0,-400.0)
>>> p3 = (140.0,-200.0,0.0)
>>> p4 = (118.0,200.0,0.0)
>>> pts = [p1,p2,p3,p4]
>>> #Define nnodes
>>> nnodes = (29,1,82)
>>> #Create and connect skewed plane
>>> qua.quadxy(nnodes,pts)
quadxyz(nnodes, pts, connect=True)

Define an arbitrary and logical set of points in 3D (xyz) space. The set of points will be connected into hexahedrons by default. Set ‘connect=False’ to prevent connection.

arg nnodes:The number of nodes including the 1st and last point along each X, Y, Z axis. The number of points will be 1 more than the number of elements in each dimension.
type nnodes:tuple(int, int, int)
arg pts:The eight corners of the hexahedron. The four bottom corners are listed first,

then the four top corners. Each set of corners (bottom and top) are defined in counter-clockwise order (the normal to the quad points is defined using the right hand rule and the order of the points). :arg pts: list of eight 3-tuples (float)

arg connect:connect points
type connect:bool
Example:
>>> #To use pylagrit, import the module.
>>> import pylagrit
>>> #Start the lagrit session.
>>> lg = pylagrit.PyLaGriT()
>>> #Create a mesh object.
>>> hex = lg.create()
>>> #Define 4 bottom points in correct order
>>> p1 = (0.0,0.0,0.0)
>>> p2 = (1.0,0.0,0.02)
>>> p3 = (1.0,1.0,0.0)
>>> p4 = (0.0,1.0,0.1)
>>> #Define 4 top points in correct order
>>> p5 = (0.0,0.0,1.0)
>>> p6 = (1.0,0.0,1.0)
>>> p7 = (1.0,1.0,1.0)
>>> p8 = (0.0,1.0,1.1)
>>> pts = [p1,p2,p3,p4,p5,p6,p7,p8]
>>> #Define nnodes
>>> nnodes = (3,3,3)
>>> #Create and connect skewed hex mesh
>>> hex.quadxyz(nnodes,pts)
>>> #Dump mesh
>>> hex.dump('quadxyz_test.gmv')
quality(*args, quality_type=None, save_att=False)
quality_angle(value, boolean='gt', save_att=False)
quality_aspect(save_att=False)
quality_edge_max(save_att=False)
quality_edge_min(save_att=False)
quality_edge_ratio(save_att=False)
quality_pcc()
read(filename, filetype=None)
recon(option='', damage='', checkaxy=False)
refine(refine_option='constant', field=' ', interpolation=' ', refine_type='element', stride=[1, 0, 0], values=[1.0], inclusive_flag='exclusive', prd_choice=None)
refine_to_object(mo, level=None, imt=None, prd_choice=None)

Refine mesh at locations that intersect another mesh object

Parameters:
  • mo (PyLaGriT mesh object) – Mesh object to intersect with current mesh object to determine where to refine
  • level (int) – max level of refinement
  • imt (int) – Value to assign to imt (LaGriT material type attribute)
  • prd_choice (int) – directions of refinement
region(boolstr, name=None)

Create region using boolean string

Parameters:
  • boolstr (str) – String of boolean operations
  • name (string) – Internal lagrit name for mesh object
Returns:

Region

Example:
>>> from pylagrit import PyLaGriT
>>> import numpy
>>> lg = PyLaGriT()
>>> mesh = lg.create()
>>> mins = (0,0,0)
>>> maxs = (5,5,5)
>>> eighth = mesh.surface_box(mins,maxs)
>>> boolstr1 = 'le '+eighth.name
>>> boolstr2 = 'gt '+eighth.name
>>> reg1 = mesh.region(boolstr1)
>>> reg2 = mesh.region(boolstr2)
>>> mreg1 = mesh.mregion(boolstr1)
>>> mreg2 = mesh.mregion(boolstr2)
>>> mesh.createpts_brick_xyz((10,10,10), (0,0,0), (10,10,10))
>>> mesh.rmregion(reg1)
>>> mesh.dump('reg_test.gmv')
region_bool(bool, name=None)

This method is deprecated and will be replaced by the MO.region() method in future releases.

regnpts(geom, ray_points, region, ptdist, stride=[1, 0, 0], irratio=0, rrz=0, maxpenetr=None)
regnpts_rtp(ray_points, region, ptdist, stride=[1, 0, 0], irratio=0, rrz=0, maxpenetr=None)

Generates points in a region previously defined by the region command. The points are generated by shooting rays through a user specified set of points from an origin point and finding the intersection of each ray with the surfaces that define the region.

Parameters:
  • ray_points (float 3-tuple) – single (x,y,z) point that defines center of spher which rays emante from
  • region (Region) – region to generate points within
  • ptdist (int float or str) – parameter that determines point distribution pattern
  • stride (int or PSet) – points to shoot rays through
  • irratio (int) – parameter that determines point distribution pattern
  • rrz (int or float) – ratio zoning value
  • maxpenetr (int or float) – maximum distance along ray that points will be distributed
regnpts_rtz(ray_points, region, ptdist, stride=[1, 0, 0], irratio=0, rrz=0, maxpenetr=None)

Generates points in a region previously defined by the region command. The points are generated by shooting rays through a user specified set of points from a line and finding the intersection of each ray with the surfaces that define the region.

Parameters:
  • ray_points (2-tuple of float 3-tuples) – two points that define cylinder which rays emante from
  • region (Region) – region to generate points within
  • ptdist (int float or str) – parameter that determines point distribution pattern
  • stride (int or PSet) – points to shoot rays through
  • irratio (int) – parameter that determines point distribution pattern
  • rrz (int or float) – ratio zoning value
  • maxpenetr (int or float) – maximum distance along ray that points will be distributed
regnpts_xyz(ray_points, region, ptdist, stride=[1, 0, 0], irratio=0, rrz=0, maxpenetr=None)

Generates points in a region previously defined by the region command. The points are generated by shooting rays through a user specified set of points from a plane and finding the intersection of each ray with the surfaces that define the region.

Parameters:
  • ray_points (3-tuple of float 3-tuples) – three points that define plane which rays emante from
  • region (Region) – region to generate points within
  • ptdist (int float or str) – parameter that determines point distribution pattern
  • stride (int or PSet) – points to shoot rays through
  • irratio (int) – parameter that determines point distribution pattern
  • rrz (int or float) – ratio zoning value
  • maxpenetr – maximum distance along ray that points will be distributed
reorder_nodes(order='ascending', cycle='zic yic xic')
resetpts_itp()

set node type from connectivity of mesh

rmmat(material_number, option='', exclusive=False)

This routine is used to remove points that are of a specified material value (itetclr for elements or imt for nodes). Elements with the specified material value are flagged by setting the element material type negative. They are not removed from the mesh object. :param material_number: Number of material :type material_number: int :param option: {‘’,’node’,’element’,’all’}, ‘node’ removes nodes with imt=material_number, ‘element’ removes elements with itetclr=material_number, ‘all’ or ‘’ removes nodes and elements with material_number equal to imt and itetclr, respectively :type option: str :param exclusive: if True, removes everything except nodes with imt=material and removes everything except elements with itetclr= material number. :type exclusive: bool

rmmat_element(material_number, exclusive=False)

This routine is used to remove elements that are of a specified material value (itetclr for elements). Elements with the specified material value are flagged by setting the element material type negative. They are not removed from the mesh object. :param material_number: Number of material :type material_number: int :param exclusive: if True, removes everything except elements with itetclr=material number. :type exclusive: bool

rmmat_node(material_number, exclusive=False)

This routine is used to remove points that are of a specified material value (imt). :param material_number: Number of material (imt) :type material_number: int :param exclusive: if True, removes everything except nodes with imt=material_number :type exclusive: bool

rmpoint_compress(filter_bool=False, resetpts_itp=True)

remove all marked nodes and correct the itet array

Parameters:resetpts_itp (bool) – set node type from connectivity of mesh
rmpoint_eltset(eltset, compress=True, resetpts_itp=True)
rmpoint_pset(pset, itype='exclusive', compress=True, resetpts_itp=True)
rmregion(region, rmpoints=True, filter_bool=False, resetpts_itp=True)

Remove points that lie inside region

Parameters:region (Region) – name of region points will be removed from
rotateln(coord1, coord2, theta, center=[0, 0, 0], copy=False, stride=(1, 0, 0))

Rotates a point distribution (specified by ifirst,ilast,istride) about a line. The copy option allows the user to make a copy of the original points as well as the rotated points, while copy=False just keeps the rotated points themselves. The line of rotation defined by coord1 and coord2 needs to be defined such that the endpoints extend beyond the point distribution being rotated. theta (in degrees) is the angle of rotation whose positive direction is determined by the right-hand-rule, that is, if the thumb of your right hand points in the direction of the line (1 to 2), then your fingers will curl in the direction of rotation. center is the point where the line can be shifted to before rotation takes place. If the copy option is chosen, the new points will have only coordinate values (xic, yic, zic); no values will be set for any other mesh object attribute for these points. Note: The end points of the line segment must extend well beyond the point set being rotated.

Example 1:
>>> from pylagrit import PyLaGriT
>>> import numpy
>>> x = numpy.arange(0,10.1,1)
>>> y = x
>>> z = [0,1]
>>> lg = PyLaGriT()
>>> mqua = lg.gridder(x,y,z,elem_type='hex',connect=True)
>>> mqua.rotateln([mqua.xmin-0.1,0,0],[mqua.xmax+0.1,0,0],25)
>>> mqua.dump_exo('rotated.exo')
>>> mqua.dump_ats_xml('rotated.xml','rotated.exo')
>>> mqua.paraview()
Example 2:
>>> from pylagrit import PyLaGriT
>>> import numpy
>>> x = numpy.arange(0,10.1,1)
>>> y = [0,1]
>>> #z = [0,1]
>>> lg = PyLaGriT()
>>> layer = lg.gridder(x=x,y=y,elem_type='quad',connect=True)
>>> layer.rotateln([0,layer.ymin-0.10,0],[0,layer.ymax+0.1,0],25)
>>> layer.dump('tmp_lay_top.inp')
>>> # Layer depths?
>>> #           1   2   3    4    5    6    7   8    9   10
>>> layers = [ .1, 1.]
>>> addnum = [  4, 2]
>>> #matnum = [2]*len(layers)
>>> matnum = [2, 1]
>>> layer_interfaces = numpy.cumsum(layers)
>>> mtop = layer.copy()
>>> stack_files = ['tmp_lay_top.inp 1,9']
>>> #stack_files.append('tmp_lay_peat_bot.inp 1,33')
>>> i = 1
>>> for li,m,a in zip(layer_interfaces,matnum,addnum):
>>>     layer.math('sub',li,'zic',cmosrc=mtop)
>>>     stack_files.append('tmp_lay'+str(i)+'.inp '+str(int(m))+', '+str(a))
>>>     layer.dump('tmp_lay'+str(i)+'.inp')
>>>     i += 1
>>> layer.math('sub',2,'zic',cmosrc=mtop)
>>> #layer.setatt('zic',-2.)
>>> layer.dump('tmp_lay_bot.inp')
>>> stack_files.append('tmp_lay_bot.inp 2')
>>> stack_files.reverse()
>>> # Create stacked layer mesh and fill
>>> stack = lg.create()
>>> stack.stack_layers('avs',stack_files,flip_opt=True)
>>> stack_hex = stack.stack_fill()
>>> stack_hex.dump_exo('rotated.exo')
>>> stack_hex.dump_ats_xml('rotated.xml','rotated.exo')
>>> stack_hex.paraview()
rzbrick(n_ijk, connect=True, stride=(1, 0, 0), coordinate_space='xyz')

Builds a brick mesh and generates a nearest neighbor connectivity matrix

Currently only configured for this flavor of syntax:

rzbrick/xyz|rtz|rtp/ni,nj,nk/pset,get,name/connect/

Use this option with quadxyz to connect logically rectangular grids.

Parameters:
  • n_ijk (tuple) – number of points to be created in each direction.
  • connect (bool) – connect points
  • stride (tuple) – Stride to select
  • coordinate_space (str) – xyz,rtz,or rtp coordinate spaces
select()
sendline(cmd, verbose=True, expectstr='Enter a command')
set_id(option, node_attname='id_node', elem_attname='id_elem')

This command creates integer attributes that contain the node and/or element number. If later operations delete nodes or elements causing renumbering, these attributes will contain the original node or element number.

Parameters:
  • option (str) – create attribute for nodes, elements, or both {‘both’,’node’,’element’}
  • node_attname (str) – name for new node attribute
  • elem_attname (str) – name for new element attribute

Example: from pylagrit import PyLaGriT #instantiate PyLaGriT lg = PyLaGriT() #create source mesh npts = (11,11,11) mins = (0.,0.,0.) maxs = (1.,1.,1.) mesh = lg.create() mesh.createpts_brick_xyz(npts,mins,maxs) #write node and element attribute numbers mesh.set_id(‘both’,node_attname=’node_att1’,elem_attname=’elem_att1’) #select and remove points p_mins = (0.5,0.,0.) p_maxs = (1.,1.,1.) points = mesh.pset_geom_xyz(p_mins,p_maxs) mesh.rmpoint_pset(points) #dump mesh with original node and element numbering saved mesh.dump(‘set_id_test.gmv’)

setatt(attname, value, stride=[1, 0, 0])
setpts(no_interface=False, closed_surfaces=False)

Set point types and imt material by calling surfset and regset routines.

Parameters:
  • ray_points (float 3-tuple) – single (x,y,z) point that defines center of spher which rays emante from
  • region (Region) – region to generate points within
  • ptdist (int float or str) – parameter that determines point distribution pattern
  • stride (int or PSet) – points to shoot rays through
  • irratio (int) – parameter that determines point distribution pattern
  • rrz (int or float) – ratio zoning value
  • maxpenetr – maximum distance along ray that points will be distributed
settets(method=None)
settets_color_points()
settets_color_tets()
settets_geometry()
settets_newtets()
settets_normal()
settets_parents()
smooth(*args, **kwargs)
stack_fill(name=None)
stack_layers(filelist, file_type='avs', nlayers=None, matids=None, xy_subset=None, buffer_opt=None, truncate_opt=None, pinchout_opt=None, flip_opt=False, fill=False)
status(brief=False, verbose=True)
subset(mins, maxs, geom='xyz')

Return Mesh Object Subset

Creates a new mesh object that contains only a geometric subset defined by mins and maxs.

Parameters:
  • mins – Coordinate of one of the shape’s defining points. xyz (Cartesian): (x1, y1, z1); rtz (Cylindrical): (radius1, theta1, z1); rtp (Spherical): (radius1, theta1, phi1);
  • maxs (tuple(int, int, int)) – Coordinate of one of the shape’s defining points. xyz (Cartesian): (x2, y2, z2); rtz (Cylindrical): (radius2, theta2, z2); rtp (Spherical): (radius2, theta2, phi2);
  • geom (str) – Type of geometric shape: ‘xyz’ (spherical), ‘rtz’ (cylindrical), ‘rtp’ (spherical)
Typep mins:

tuple(int, int, int)

Returns: MO object

Example:
>>> #To use pylagrit, import the module.
>>> import pylagrit
>>> #Start the lagrit session.
>>> lg = pylagrit.PyLaGriT()
>>> #Create a mesh object.
>>> mo = lg.create()
>>> mo.createpts_brick_xyz((5,5,5), (0,0,0), (5,5,5))
>>> #Take the subset from (3,3,3)
>>> mo.subset((3,3,3),(5,5,5))
subset_rtp(mins, maxs)

Return Spherical MO Subset

Creates a new mesh object that contains only a spherical subset defined by mins and maxs.

Parameters:
  • mins (tuple(int, int, int)) – Defines radius1, theta1, and phi1.
  • maxs (tuple(int, int, int)) – Defines radius2, theta2, and phi2.

Returns: MO object

subset_rtz(mins, maxs)

Return Cylindrical MO Subset

Creates a new mesh object that contains only a cylindrical subset defined by mins and maxs.

Parameters:
  • mins (tuple(int, int, int)) – Defines radius1, theta1, and z1.
  • maxs (tuple(int, int, int)) – Defines radius2, theta2, and z2.

Returns: MO object

subset_xyz(mins, maxs)

Return Tetrehedral MO Subset

Creates a new mesh object that contains only a tetrehedral subset defined by mins and maxs.

Parameters:
  • mins (tuple(int, int, int)) – Coordinate point of 1 of the tetrahedral’s corners.
  • maxs (tuple(int, int, int)) – Coordinate point of 1 of the tetrahedral’s corners.

Returns: MO object

surface(name=None, ibtype='reflect')
surface_box(mins, maxs, name=None, ibtype='reflect')
surface_cylinder(coord1, coord2, radius, name=None, ibtype='reflect')
surface_plane(coord1, coord2, coord3, name=None, ibtype='reflect')
trans(xold, xnew, stride=(1, 0, 0))

Translate mesh according to old coordinates “xold” to new coordinates “xnew”

Parameters:
  • xold (tuple(float,float,float)) – old position
  • xnew (tuple(float,float,float)) – new position
  • stride (tuple(int,int,int)) – tuple of (first, last, stride) of points
tri_mesh_output_prep()

Prepare tri mesh for output, remove dudded points, ensure delaunay volumes, etc. Combination of lagrit commands: filter/1 0 0 rmpoint/compress recon/1 resetpts/itp

triangulate(order='clockwise')

triangulate will take an ordered set of nodes in the current 2d mesh object that define a perimeter of a polygon and create a trangulation of the polygon. The nodes are assumed to lie in the xy plane; the z coordinate is ignored. No checks are performed to verify that the nodes define a legal perimeter (i.e. that segments of the perimeter do not cross). The code will connect the last node to the first node to complete the perimeter.

This code support triangulation of self-intersecting polygons (polygon with holes), assuming that the order of the nodes are correct. Moreover the connectivity of the polyline must also be defined correctly. No checks are made.

One disadvantage of the algorithm for triangulating self-intersecting polygons is that it does not always work. For example, if the holes have complicated shapes, with many concave vertices, the code might fail. In this case, the user may try to rotate the order of the nodes:
NODE_ID:
1 -> 2 2 -> 3 … N -> 1
param order:direction of point ordering
type order:string
Example:
>>> from pylagrit import PyLaGriT
>>> # Create pylagrit object
>>> lg = PyLaGriT()
>>> # Define polygon points in clockwise direction
>>> # and create tri mesh object
>>> coords = [[0.0, 0.0, 0.0],
>>>           [0.0, 1000.0, 0.0],
>>>           [2200.0, 200.0, 0.0],
>>>           [2200.0, 0.0, 0.0]]
>>> motri = lg.tri_mo_from_polyline(coords)
>>> # Triangulate polygon
>>> motri.triangulate()
>>> motri.setatt('imt',1)
>>> motri.setatt('itetclr',1)
>>> # refine mesh with successively smaller edge length constraints
>>> edge_length = [1000,500,250,125,75,40,20,15]
>>> for i,l in enumerate(edge_length):
>>>     motri.resetpts_itp()
>>>     motri.refine(refine_option='rivara',refine_type='edge',values=[l],inclusive_flag='inclusive')
>>>     motri.smooth()
>>>     motri.recon(0)
>>> # provide additional smoothing after the last refine
>>> for i in range(5):
>>>     motri.smooth()
>>>     motri.recon(0)
>>> # create delaunay mesh and clean up
>>> motri.tri_mesh_output_prep()
>>> # dump fehm files
>>> motri.dump_fehm('nk_mesh00')
>>> # view results
>>> motri.paraview() 
upscale(method, attsink, cmosrc, attsrc=None, stride=(1, 0, 0), boundary_choice=None, keepatt=False, set_id=False)

The upscale command is used to interpolate attribute values from nodes of a fine source mesh to node attributes of a coarse sink mesh. The subroutine finds nodes of the fine source mesh within the Voronoi cell of every node in the coarser sink mesh. Nodes on cell boundaries are assigned to two or more sink nodes. Then the attributes of all the source nodes within a source node’s cell are upscaled into a single value based on the chosen method. Mesh elements and connectivity are ignored and only node values are used to upscale values on to the sink mesh nodes.

Parameters:
  • method (str) – Type of upscaling: sum, min, max, and averages ariave, harave, geoave
  • attsink (str) – attribute sink
  • cmosrc (PyLaGriT Mesh Object) – PyLaGriT mesh object source
  • attsrc (str) – attribute src, defaults to name of attsink
  • stride (tuple(int)) – tuple of (first, last, stride) of points
  • boundary_choice (str) – method of choice when source nodes are found on the boundary of multiple Voronoi volumes of sink nodes: single, divide, or multiple
upscale_ariave(attsink, cmosrc, attsrc=None, stride=(1, 0, 0), boundary_choice=None, keepatt=False, set_id=False)

Upscale using arithmetic average of cmosrc points within Voronoi volumes of current mesh

Parameters:
  • attsink (str) – attribute sink
  • cmosrc (PyLaGriT Mesh Object) – PyLaGriT mesh object source
  • attsrc (str) – attribute src
  • stride (tuple(int)) – tuple of (first, last, stride) of points
  • boundary_choice (str) – method of choice when source nodes are found on the boundary of multiple Voronoi volumes of sink nodes: single, divide, or multiple
upscale_geoave(attsink, cmosrc, attsrc=None, stride=(1, 0, 0), boundary_choice=None, keepatt=False, set_id=False)

Upscale using geometric average of cmosrc points within Voronoi volumes of current mesh

Parameters:
  • attsink (str) – attribute sink
  • cmosrc (PyLaGriT Mesh Object) – PyLaGriT mesh object source
  • attsrc (str) – attribute src
  • stride (tuple(int)) – tuple of (first, last, stride) of points
  • boundary_choice (str) – method of choice when source nodes are found on the boundary of multiple Voronoi volumes of sink nodes: single, divide, or multiple
upscale_harave(attsink, cmosrc, attsrc=None, stride=(1, 0, 0), boundary_choice=None, keepatt=False, set_id=False)

Upscale using harmonic average of cmosrc points within Voronoi volumes of current mesh

Parameters:
  • attsink (str) – attribute sink
  • cmosrc (PyLaGriT Mesh Object) – PyLaGriT mesh object source
  • attsrc (str) – attribute src
  • stride (tuple(int)) – tuple of (first, last, stride) of points
  • boundary_choice (str) – method of choice when source nodes are found on the boundary of multiple Voronoi volumes of sink nodes: single, divide, or multiple
upscale_max(attsink, cmosrc, attsrc=None, stride=(1, 0, 0), boundary_choice=None, keepatt=False, set_id=False)

Upscale using maximum of cmosrc points within Voronoi volumes of current mesh

Parameters:
  • attsink (str) – attribute sink
  • cmosrc (PyLaGriT Mesh Object) – PyLaGriT mesh object source
  • attsrc (str) – attribute src
  • stride (tuple(int)) – tuple of (first, last, stride) of points
  • boundary_choice (str) – method of choice when source nodes are found on the boundary of multiple Voronoi volumes of sink nodes: single, divide, or multiple
upscale_min(attsink, cmosrc, attsrc=None, stride=(1, 0, 0), boundary_choice=None, keepatt=False, set_id=False)

Upscale using minimum of cmosrc points within Voronoi volumes of current mesh

Parameters:
  • attsink (str) – attribute sink
  • cmosrc (PyLaGriT Mesh Object) – PyLaGriT mesh object source
  • attsrc (str) – attribute src
  • stride (tuple(int)) – tuple of (first, last, stride) of points
  • boundary_choice (str) – method of choice when source nodes are found on the boundary of multiple Voronoi volumes of sink nodes: single, divide, or multiple
upscale_sum(attsink, cmosrc, attsrc=None, stride=(1, 0, 0), boundary_choice=None, keepatt=False, set_id=False)

Upscale using sum of cmosrc points within Voronoi volumes of current mesh

Parameters:
  • attsink (str) – attribute sink
  • cmosrc (PyLaGriT Mesh Object) – PyLaGriT mesh object source
  • attsrc (str) – attribute src
  • stride (tuple(int)) – tuple of (first, last, stride) of points
  • boundary_choice (str) – method of choice when source nodes are found on the boundary of multiple Voronoi volumes of sink nodes: single, divide, or multiple
xlength
xmax
xmin
ylength
ymax
ymin
zlength
zmax
zmin