Tag: Hatching

Export L-PBF Scan Vectors to VTK/Paraview via PySLM

The in-built visualisation for scan paths in PySLM leverages matplotlib – refer to a previous post. This is sufficient for most user’s needs when attempting to interpret and visualise the scan paths generated in PySLM, or those imported from a slice taken from an existing machine build files. Extending this beyond multiple layers or large parts becomes more tricky when factoring in visualisation of some parameters (e.g. Laser Power, effective scan speed). Admittedly, the performance of Matplotlib becomes limited to explore the intricacies and complexities embedded within the scan vectors. 

For scientific research, the fusion of scan vector geometry with volumetric datasets such as X-Ray CT during post-inspection of parts/samples, or those generated within the build process including pyrometry data, thermal-imaging offer the ability to increase our understanding and insight to observations of the effect of process on the material produced using L-PBF.  GPU based visualisation libraries such (vispy) would offer the possibility to accelerate the performance, but are not user-friendly nor offer interactivity when manipulating views and the data and are often cumbersome when processing volumetric datasets often encountered in Additive Manufacturing. Paraview is a cross-platform open-source scientific visualisation tool that is especially powerful for processing, interaction and visualisation of large-scale scientific datasets.

Paraview and the underlying VTK library offers an alternative ready-made solution to visualise this information, and are most importantly hardware accelerated with the option for raytracing provided by OSPRay and OptiX for latest RTX NVIDIA cards that include Raytracing (RT) cores. Additionally, the data can be augmented and processed using parallelised filters and tools in Paraview.

VTK File Format

Ignoring the HDF5 variations that are most useful for structured data, the underlying format within vtk that used for storing vector based data and point cloud data is the .vtp file format. The modern VTK file formats use an XML schema – unlike the legacy format, to store a structured series of geometry (volumetric data, lines, polygons, 3D elements and point clouds). The internal data format can be stored using ascii encoding or binary. Binary data can be incorporated directly within a parsable .xml format using a Base64 encoding and may additional incorporate internal compression. Alternatively data can be stored in an appended data section located at the footer of the file, which treats data section as a contiguous block of raw data. Different sub-formats exist, that are appropriate for different types of data e.g. volumetric, element based (Finite Volume / Finite Element derived) or polygon based. An approach relevant to export scan vector geometry the .vtp – format is most suitable.

The data stored in the VTK Point file consists of:

  • 3D points coordinates
  • Data attributes stored at each point location
  • Geometric elements (lines, polygons) defining connectivity with reference to the list of point coordinates

Paraview exporter implementation:

The Paraview exporter is simplistic, because the data compression is currently ignored. The process is similar to the technique used in the function pyslm.visualise.plotSequential, whereby hatch and contour vectors are merged and reprocessed in order that they represent always a series of lines (an n x 2 x 2 array). This is not the most efficient option for ContourGeometry (border scans) where scan vectors are continuously joined up, but simplifies the processing working with the data.  

Once the scan vector coordinates and the relevant data are packaged up into a single array, the data is wrote within the sub-sections of the XML file. Data is stored using floating points or integers accordingly in a binary representation. The data used to represent coordinates and indices for each vector, are stored with the ‘appended’ option within the <DataArray> element of each section. The raw data is stored and collected that are then written in the <AppendedData> element at the end of file with raw encoding option chosen. The byte offsets for the position of each ‘chunk’ of data that are referenced by the <DataArray> element are collected and stored incrementally.

For reference, the following information is provided for writing raw data, because this was difficult to obtain from the VTK documentation directly.

<AppendedData encoding=”raw”> Start of Raw Data Section
_ Underscore character is starting location for reading raw data
Section Size (Int32/Int64)Integer representing size of following section (include the size in bytes
with the offsets provided
). The integer type should match the size used in the header.
Raw data (e.g Int32, float32, float64)
….Repeated the above (two rows) for each referenced data section
</AppendedData>

Example Scan Vector Data exported to VTK

An example Aconity .ILT file was imported into PySLM and then exported to a .vtp VTK file that was processed in Paraview. The scan order is visualised by the colour map with each vertex assigned a global-id. The ‘Tube‘ filter was applied to each scan vector in order to improve their visibility.

Visualisation scan-vectors for L-PBF/SLM processed in PySLM and exported to VTK (.vtp) file format

The script excerpt can currently be found on a Gist. This will be later included in future versions of PySLM along with other import/exporters.

PySLM: Geometric Hatch Overlap Check/Visualisation

Designing scan strategies for PBF techniques, we are not entirely aware of situations that arise where the powder-bed is not fully exposed due to a mismatch when scan vectors are not sufficiently overlapped. Typically, unoptimised placement of hatch vectors lead to the creation of irregular porosity or voids in L-PBF parts.

This can arise along the intersections between the contour and interior hatches, especially a long concave regions such as sharp corner features with acute angles.

The approach is not an efficient way to examine the presence of , but provides a representative view for checking this geometrically.

Method

The approach takes advantage of the relatively new Iterator classes available within the analysis module, which vastly simplifies the generation procedure for manipulating and examining existing scan vector geometry. Firstly, generate or alternatively import the Layer and its LayerGeometry groups to examine. The group of layers are passed to the ScanVectorIterator class, which will iterate across every scan vector from both ContourGeometry and HatchGeometry objects within a Layer. Single point exposures are not considered.

import pyslm.analysis
from shapely.geometry import LineString, Polygon, MultiPolygon
from shapely.ops import cascaded_union

scanIterator =  pyslm.analysis.ScanVectorIterator([layer])

After the creation of the ScanVectorIterator, this can be readily expedited to process across all scan vectors across the Layer. The basic process relies on converting each scan vector to a Shapely polygon objects and then processing them using the geometry tools available.

For this case we use Pythonic notation to compactly operate across each scan vector and collect them. We convert each scan vector to a Shapely LineString, which has a method to then offset or buffer.

# Laser Spot Radius
laserSpotRadius = 0.04

# Iterate across each scan vector and buffer than geometry
lines = [LineString(line).buffer(laserSpotRadius) for line in scanIterator]

After offseting all the lines, this can be easily visualised by conversion to a Shapely MultiPolygon.

# Merged the offset lines into a Shapely Multi-Polygon Collection
multiPoly = MultiPolygon(lines)
pyslm.visualise.plotPolygon(multiPoly)

The geometrical result is shown below. It is relatively quick to generate and plot individual scan vectors as shown below:

Geometrical overlap of scan vectors in SLM processed by PySLM
Overlap of hatch vectors represented by geometrically offsetting the individual scan vectors within a Layer.

Each scan vector that is offset is represented by a Shapely.Polygon. It is trivial to perform a boolean operation with the Shapely library. Although, it is recommended to use the more efficient, albeit still relatively slow, shapely.ops.cascaded_union function to merge multiple geometries together:

# Cascaded union is a more efficient boolean merge for multiple polygon entities
multiPolyMerged = cascaded_union(multiPoly)

pyslm.visualise.plotPolygon(multiPolyMerged)

The combined result is shown here with a slightly smaller hatch distance to exaggerate the effect and highlight regions where the laser beam may not sufficiently provide exposure to the powder bed:

Illustration of regions that may
Regions that have insufficient coverage observed after performing a boolean merge after the scan vectors have been offset

This post shares a relatively simple example exploration using geometrical operations and the iterator class to understand potential issues related to hatch overlaps.

Final Conclusions

Potentially, this check could be extended into 3D using morphological operations. This would provide a more qualitative examination of porosity generation as a result of Furthermore, the combined use of Neural Networks or reduce order models could provide a representative exposure area to provide the geometrical offset in 2D and provide a prediction for coverage spatially in 3D.

Custom Scan Strategies in SLM / L-PBF with PySLM: Sinusoidal Scanning

Building upon the previous post that provided a detailed breakdown for creating custom island scan strategies, this further post documents a method for deploying custom ‘hatch’ infills. This is particularly desirable capability sought by researchers and has been touched upon very little in the current research. The use of unit-cell infills or in particular fractal filling curves such as the Hilbert curve have been sought for better controlling the thermal history and melt pool stability of hatch infills.

This has been previously explored in SLS [1]][2] and in SLM on a previous collaborators at the University of Nottingham investigating Fractal scanning strategy [3][4].

Typically, hatch infills are sequences of linear lines that form the the ‘hatch’ pattern. Practically, these are very efficient mechanism for infilling a 2D area by using 1D line elements when rastering a laser. Clipping of lines within polygons is intuitive. As discussed there are various scan strategies that can be employed to generate variations on this infill – i.e. stripe, checkerboard/island scan strategy and also modifying the order or sorting of the hatch vectors.

Geometrical scan strategies that adapt the infill based on the underlying geometry, i.e. lattices are acknowledges as ways for drastically improving the performance and the quality of these characteristic structures. This would be based on some medial-axis approach. This post will not specifically delve into this, rather, demonstrate an approach for custom infills on bulk regions.

Ultimately, drastically changing the behavior of the underlying hatch infill has not really been explored. This post will demonstrate an example that could be employed and explored as part of future research.

Custom Sinusoidal Approach

Sinusoidal scanning has been employed in welding research [5] and also in direct energy deposition (DED) [6][7][8] in order to improve the stability and quality of the joining or manufacturing process.

The process of generating this particular scan strategy requires some careful thought to improve the efficiency of the generation, especially given the overall increase in number of points require to essentially ‘sample’ across the sin curve.

The implementation requires subclassing the Hatcher class, by re-implementing the BaseHatcher.generateHatching and the BaseHatcher.hatch methods.

Unlike, the normal hatch vectors, the sinusoidal pattern has to be treated as a series of connected line segments, without any jumping. This requires using the ContourGeometry representation to efficiently store the discretised curve. As a result, the Hatcher.hatch method has to be re-implemented to take account of this.

The procedure builds upon previous methods to define customer behavior (see previous post). The first steps are to define a local coordinate system x' and y' for generating the individual sin curve. A sine curve y' = A \sin(k x') is generated to fill the region bounding box accordingly, given a frequency and amplitude parameter along x'.

The number of points used to discretise the sine curve is determined by \delta x. This needs to be chosen to suit the parameters for the periodicity and amplitude of the sine curve. A reasonable compromise is require as this will severely impact both the performance of clipping these curves, but also the overall file size of the build file generated.

dx = self._discretisation # num points per mm
numPoints = 2*bboxRadius * dx

x = np.arange(-bboxRadius, bboxRadius, hatchSpacing, dtype=np.float32).reshape(-1, 1)
hatches = x.copy()

"""
Generate the sinusoidal curve along the local coordinate system x' and y'. These will be later tiled and then
transformed across the entire coordinate space.
"""
xDash = np.linspace(-bboxRadius, bboxRadius, int(numPoints))
yDash = self._amplitude * np.sin(2.0*np.pi * self._frequency * xDash)

"""
We replicate and transform the sine curve along adjacent paths and transform along the y-direction
"""
y = np.tile(yDash, [x.shape[0], 1])
y += x

x = np.tile(xDash, [x.shape[0],1]).flatten()
y = y.ravel()

After generating single sine curve, numpy.tile is used to efficiently replicate the curve to fill the entire bounding box region. Each curve is then translated by an increment defined by x, to represent the effective hatch spacing or hatch distance.

The next important step is to define the sort order for scanning these. This is slightly different, in that the sort order is done per line segment used to discretise the curve. This is subtle, but very important because this ensures that the curves when clipped by the slice boundary are scanned in the same prescribed sequential order.

The increment of 1\times10^5 is used in order to potentially differentiate each curve later, if required.

# Seperate the z-order index per group
inc = np.arange(0, 10000*(xDash.shape[0]), 10000).astype(np.int64).reshape(-1,1)
zInc = np.tile(inc, [1,hatches.shape[0]]).flatten()
z += zInc

coords = np.hstack([x.reshape(-1, 1),
                    y.reshape(-1, 1),
                    z.reshape(-1, 1)])

Following the generation of these sinusoidal curves, a transformation matrix is applied accordingly, before these are clipped in the Hatcher.hatch method.

The next crucial difference, that has been implemented from PySLM version 0.3, is a new clipping method, BaseHatcher.clipContourLines. The following method is different from BaseHatcher.clipLines, in that clips ContourGeometry separately. This is important for keeping the scan vectors separate and in the correct order, which would be otherwise difficult to achieve. The clipped results are implicitly separated into contour geometry groups.

hatches = self.generateHatching(paths, self._hatchDistance, layerHatchAngle)

clippedPaths = self.clipContourLines(paths, hatches)

# Merge the lines together
if len(clippedPaths) > 0:
    for path in clippedPaths:
        clippedLines = np.vstack(path) 
        
        clippedLines = clippedLines[:,:2]
        contourGeom = ContourGeometry()

        contourGeom.coords = clippedLines.reshape(-1, 2)

        layer.geometry.append(contourGeom)

The next step is to sort the clipped paths into the right order. This is done by using the 1st value of 3rd index column accordingly sorting using sorted with a lambda function.

"""
Sort the sinusoidal vectors based on the 1st coordinate's sort id (column 3). This only sorts individual paths
rather than the contours internally.            
"""
clippedPaths = sorted(clippedPaths, key=lambda x: x[0][2])

Now, the result of the sinusoidal scan strategy can be visualised below.

Sinusoidal Hatch Scan Strategy for Selective Laser Melting - PySLM
Sinusoidal Hatch Scan Strategy for Selective Laser Melting – PySLM

This approach currently is very intensive to generate during the clipping operation, due to the number of edges along each clipping operation. Using the previous techniques with the island scan strategy in a previous post, could be use to amorotise a lot of the cost of clipping.

Example Script

The script is available on github at examples/example_custom_sinusoidal_scanning.py

References

References
1 Yang, J., Bin, H., Zhang, X., & Liu, Z. (2003). Fractal scanning path generation and control system for selective laser sintering (SLS). International Journal of Machine Tools and Manufacture, 43(3), 293–300. https://doi.org/10.1016/S0890-6955(02)00212-2
2 Ma, L., & Bin, H. (2006). Temperature and stress analysis and simulation in fractal scanning-based laser sintering. The International Journal of Advanced Manufacturing Technology, 34(9–10), 898–903. https://doi.org/10.1007/s00170-006-0665-5
3 Catchpole-Smith, S., Aboulkhair, N., Parry, L., Tuck, C., Ashcroft, I. A., & Clare, A. (2017). Fractal scan strategies for selective laser melting of ‘unweldable’ nickel superalloys. Additive Manufacturing, 15, 113–122. https://doi.org/10.1016/j.addma.2017.02.002
4 Sebastian, R., Catchpole-Smith, S., Simonelli, M., Rushworth, A., Chen, H., & Clare, A. (2020). ‘Unit cell’ type scan strategies for powder bed fusion: The Hilbert fractal. Additive Manufacturing, 36(July), 101588. https://doi.org/10.1016/j.addma.2020.101588
5 Tongtong Liu, Zhongyan Mu, Renzhi Hu, Shengyong Pang,
Sinusoidal oscillating laser welding of 7075 aluminum alloy: Hydrodynamics, porosity formation and optimization, International Journal of Heat and Mass Transfer, Volume 140, 2019, Pages 346-358, ISSN 0017-9310, https://doi.org/10.1016/j.ijheatmasstransfer.2019.05.111
6 Cao, Y., Zhu, S., Liang, X., & Wang, W. (2011). Overlapping model of beads and curve fitting of bead section for rapid manufacturing by robotic MAG welding process. Robotics and Computer-Integrated Manufacturing, 27(3), 641–645. https://doi.org/10.1016/j.rcim.2010.11.002
7 Zhang, W., Tong, M., & Harrison, N. M. (2020). Scanning strategies effect on temperature, residual stress and deformation by multi-laser beam powder bed fusion manufacturing. Additive Manufacturing, 36(June), 101507. https://doi.org/10.1016/j.addma.2020.101507
8 Ding, D., Pan, Z., Cuiuri, D., & Li, H. (2015). A multi-bead overlapping model for robotic wire and arc additive manufacturing (WAAM). Robotics and Computer-Integrated Manufacturing, 31, 101–110. https://doi.org/10.1016/j.rcim.2014.08.008

Custom Island Scan Strategies for L-PBF/SLM using PySLM

The fact that most island scan strategies employed in SLM are nearly always square raised the question whether we could do more. I recently came across this ability to define ‘hexagon’ island regions advertised in the 2020 release of Autodesk Netfabb. Unfortunately this is a commercial tool and not always available. The practical reasons for implementing a hexagon island scanning strategy are largely unclear, but this prompted to create an example to illustrate how one would create custom island regions using PySLM. This in future could open some interesting ideas of tuning the scan strategy spatially across a layer.

Structural materials in cells - OpenLearn - Open University - T356_3
Honeycombs or heaxgonal lattices observed in nature are a popular structure used in composites engineering. Could the same be applied in Additive Manufacturing?

The user needs to customise the behaviour they desire by deriving subclasses from:

These classes serve the purpose for defining a ‘regular’ tessellated sub-region containing hatches. Regular regions that share the same shape characteristics for using the infill optimises the overall clipping performance outlined in the previous post.

PySLM: Checkerboard Island Scan Strategy Implementation used for L-PBF (Selective Laser Melting)
Illustration of Checkerboard Island Scan Strategy Implementation

Theoretically, we could build 2D unstructured cells e.g. Voronoi patterns, however, internally hatches for each region will require individual clipping and penalised with a significant performance hit during the hatching process.

Voronoi Diagram --
Example of a Voronoi diagram: regions are dibi based on the boundaries between.

The Island subclass region is the most important part to re-define the behavior. If we want to change the island regions to become regular tessellated polygons, the localBoundary method should be re-defined. In this example, it will generate a hexagon region, but the implementation below should be generic to cover other N-gon primitives:

   def localBoundary(self) -> np.ndarray:
    # Redefine the local boundary to be the hexagon shape

    if HexIsland._boundary is None:
        # Simple approach is to use a radius to define the overall island size
        #radius = np.sqrt(2*(self._islandWidth*0.5 + self._islandOverlap)**2)

        numPoints = 6

        radius = self._islandWidth / np.cos(np.pi/numPoints)  / 2 + self._islandOverlap

        print('island', radius, self._islandWidth)

        # Generate polygon island
        coords = np.zeros((numPoints+1, 2))

        for i in np.arange(0,numPoints):
            # Subtracting -0.5 orientates the polygon along its face
            angle = (i-0.5)/numPoints*2*np.pi
            coords[i] = [np.cos(angle), np.sin(angle)]

        # Close the polygon
        coords[-1] = coords[0]

        # Scale the polygon
        coords *= radius

        # Assign to the static class attribute
        HexIsland._boundary = coords

    return HexIsland._boundary

The polygon shape is defined by numPoints, so this can be changed to another polygon if desired. The polygon boundary is defined using a radius for the island region and from this a regular polygon is constructed on the outside. The polygon points are rotated by adjusting the start angle so there is a vertical edge on the RHS.

PySLM SLM Additive Manufacturing Scan Stragies: Hexagonal Island Tessellation
The Polygon is constructed around the island size (radius) and is orientated with the RHS edge vertically

This is generated once as a static class attribute, stored in _boundary to remove the overhead when generating the boundary.

The next step is to generate the internal hatch, which in this occasion needs to be clipped with the local boundary. First, the hatch vectors are generated covering the exterior region using the same radius as the polygon. This ensures that for any rotation transformation of the hatch vectors within the island are fully covered. This is relatively familiar to other code which generates these.

def generateInternalHatch(self, isOdd = True) -> np.ndarray:
    """
    Generates a set of hatches orthogonal to the island's coordinate system :math:`(x\\prime, y\\prime)`.

    :param isOdd: The chosen orientation of the hatching
    :return: (nx3) Set of sorted hatch coordinates
    """

    numPoints = 6

    radius = self._islandWidth / np.cos(np.pi / numPoints) / 2 + self._islandOverlap

    startX = -radius
    startY = -radius

    endX = radius
    endY = radius

    # Generate the basic hatch lines to fill the island region
    x = np.tile(np.arange(startX, endX, self._hatchDistance).reshape(-1, 1), 2).flatten()
    y = np.array([startY, endY])
    y = np.resize(y, x.shape)

    z = np.arange(0, y.shape[0] / 2, 0.5).astype(np.int64)

    coords =  np.hstack([x.reshape(-1, 1),
                            y.reshape(-1, 1),
                            z.reshape(-1,1)])

    # Toggle the hatch angle
    theta_h = np.deg2rad(90.0) if isOdd else np.deg2rad(0.0)

    # Create the 2D rotation matrix with an additional row, column to preserve the hatch order
    c, s = np.cos(theta_h), np.sin(theta_h)
    R = np.array([(c, -s, 0),
                  (s, c, 0),
                  (0, 0, 1.0)])

    # Apply the rotation matrix and translate to bounding box centre
    coords = np.matmul(R, coords.T).T

The next stage is to clip the hatch vectors with the local boundary. This is achieved using the static class method hatching.BaseHatcher.clipLines. The clipped hatches need to be sorted using the ‘z’ index or 2nd column of the clippedLines.

# Clip the hatch fill to the boundary
boundary = [[self.localBoundary()]]
clippedLines = np.array(hatching.BaseHatcher.clipLines(boundary, coords))

# Sort the hatches
clippedLines = clippedLines[:, :, :3]
id = np.argsort(clippedLines[:, 0, 2])
clippedLines = clippedLines[id, :, :]

# Convert to a flat 2D array of hatches and resort the indices
coordsUp = clippedLines.reshape(-1,3)
coordsUp[:,2] = np.arange(0, coordsUp.shape[0] / 2, 0.5).astype(np.int64)
return coordsUp

After sorting, the ‘z’ indexes need to the be condensed or flattened by re-building the ‘z’ index into sequential order. This is done to ensure when the hatches for islands are merged, we simply increment the index of the island using the length of the hatch array rather than performing np.max each time. This is later seen in the method hatching.IslandHatcher.hatch

# Generate the hatches for all the islands
idx = 0
for island in sortedIslands:

    # Generate the hatches for each island subregion
    coords = island.hatch()

    # Note for sorting later the order of the hatch vector is updated based on the sortedIsland
    coords[:, 2] += idx
    ...
    
    ...
    # 
    idx += coords.shape[0] / 2

clippedCoords = np.vstack(clippedCoords)
unclippedCoords = np.vstack(unclippedCoords).reshape(-1,2,3)

HexIslandHatcher

The final stage, is to re-implement hatching.IslandHatcher as a subclass. In this class, at a minimum, the generateIsland method needs to be redefined to correctly positioned the islands so that they tessellate correctly.

def generateIslands(self, paths, hatchAngle: float = 90.0):
    """
    Generate a series of tessellating Hex Islands to fill the region. For now this requires re-implementing because
    the boundaries of the island may be different shapes and require a specific placement in order to correctly
    tessellate within a region.
    """

    # Hatch angle
    theta_h = np.radians(hatchAngle)  # 'rad'

    # Get the bounding box of the boundary
    bbox = self.boundaryBoundingBox(paths)

    print('bounding box bbox', bbox)
    # Expand the bounding box
    bboxCentre = np.mean(bbox.reshape(2, 2), axis=0)

    # Calculates the diagonal length for which is the longest
    diagonal = bbox[2:] - bboxCentre
    bboxRadius = np.sqrt(diagonal.dot(diagonal))

    # Number of sides of the polygon island
    numPoints = 6

    # Construct a square which wraps the radius
    numIslandsX = int(2 * bboxRadius / self._islandWidth) + 1
    numIslandsY = int(2 * bboxRadius / ((self._islandWidth + self._islandOverlap) * np.sin(2*np.pi/numPoints)) )+ 1

The key difference here is defining the number of islands in the y-direction to account for the tessellation of the polygons. This is a simple geometry problem. The y-offset for the islands is simply the vertical component of the 2 x island radius at the angular increment to form the polygon.

Example of tesselation of hexagon islands

The HexIsland are generated with the offsets and appended to the list. These are then treat internally by the parent class IslandHatcher.

...

...

for i in np.arange(0, numIslandsX):
    for j in np.arange(0, numIslandsY):

        # gGenerate the island position
        startX = -bboxRadius + i * self._islandWidth + np.mod(j, 2) * self._islandWidth / 2
        startY = -bboxRadius + j * (self._islandWidth) * np.sin(2*np.pi/numPoints)

        pos = np.array([(startX, startY)])

        # Apply the rotation matrix and translate to bounding box centre
        pos = np.matmul(R, pos.T)
        pos = pos.T + bboxCentre

        # Generate a HexIsland and append to the island
        island = HexIsland(origin=pos, orientation=theta_h,
                            islandWidth=self._islandWidth, islandOverlap=self._islandOverlap,
                            hatchDistance=self._hatchDistance)

        island.posId = (i, j)
        island.id = id
        islands.append(island)

        id += 1

return islands

The island tessellation generated is shown below, with the an offset between islands applied by modifying the radius.

PySLM - Additive Manufacturing Library for Selective Laser Melting. The figure shows the generation of hexagonal hatch island regions.
Hexagon Island Boundaries generated across the entire region. The boundaries of the layer are shown, which are used for the intersection test.

The fully clipped scan strategy is shown below with the scanning ordered in the Y-direction.

PySLM - Additive Manufacturing Library for Selective Laser Melting. Figure shows the fully clipped hexagon islands in a custom island scan strategy
Hexagonal Island Scan Strategy: Consists of 5 mm Island (radius) with an offset at the boundaries of 0.1 mm.

Conclusions

This post illustrates how one can effectively decompose a layer region into a series of repeatable ‘island’ units which can be processed in an efficient manner, by only clipping hatches at boundary regions. This potentially has the ability to define spatially aware island regions; for example this could be redefining island sizes or parameters towards the boundary of a part. It could be used to alter the scan strategies within the region too, with the effect of changing the thermal behavior.

The full excerpt of the example can be found on github at examples/example_custom_island_hatcher.py.

PySLM: Exposure Map Generation & Visualisation

A curious feature was seeing the possibility to both generate and visualise the exposure points generated in along scan vectors. This also includes generating a pseudo ‘heatmap’ or exposure map of the energy deposited by the laser scanning across the surface. Variations of this idea, I have seen across a variety of commercial software – in particular Renishaw’s QuantAM, which is perhaps the source of inspiration for coming up with an approach.

Usually, there is not much indication as the relative exposure of energy from the laser into a layer based both the chosen laser parameters applied to the scan vectors of a specific LayerGeometry group and also the overall spatial distribution of scan vectors in a region. For instance, overlaps between scan vectors will in theory deposit more energy across time into a localised region . Ideally, the exposure of energy from the laser is very uniform across the regions.

It is debatable whether there is particular usefulness of visualising the exposure map of the laser – but in practice can be used to differentiate laser scan parameters used on models in parametric studies. With further thought, a cheap numerical simulation could potentially capture the transient heat distribution generated by the sequential deposition of these exposure points.

A Light Background on Laser Exposure:

Hopefully someone can further comment on of this and perhaps correct me…

SLM systems are built using fiber laser systems offering a very fine, high quality controlled laser beam. Two types of laser operation exist: continuous wave (CW) and pulsed.

In practice, with pulsed lasers used in AM platforms, the scan vector is decomposed across a number of exposure points at a set ‘pointDistance’ from each other along the scan vector. Each exposure irradiates energy into the powder-bed surface at a fixed energy (inferred from the average laserPower) for a set time period referred as a ‘pointExposureTime‘. After exposure, the galvo-mirrors instantaneously move to direct the beam to the next exposure point. The relative speed this occurs optically gives the illusion that the laser is continuously moving. Continuous wave, as the name describes, emits continuous irradiation whilst the beam traverses the surface. The laser speed parameter is used here.

Different systems use either pulsed, continuous, or even both modes of operation, offering different benefits. This in principle is the reason for BuildStyle class containing the following parameters to cover both situations:

  • Point Exposure Time
  • Point Distance
  • Laser Speed
  • Laser Power – used across both modes

Ultimately, this does not matter, but the proposed method obviously requires the user to set at minimum the pointDistance parameter in order for the following methodology to work.

Proposed Methodology

Taking an existing Layer containing LayerGeometry, the exposure points are generated individually for each scan vector. Given the potential number of scan vectors within a layer, it would be beneficial to perform this efficiently.

PySLM: Island Hatching for SLM
Example part that is hatched with a series of overlapping 5mm island. The scan order is indicated by the blue line.

The basic process to do this is based on the following. The scan vector \textbf{v} has a particular length. Based on the point exposure distance, we can equally divide along this a number of points. Using both v_o and its direction ,we can then project a given number of exposure points along this line.

PySLM: definition of the laser hatch vectors and their decomposition into exposure points.
Definition of a hath vector and decomposing each scan vector into a series of exposure points.

An attempt to do this more efficiently and exploit vectorisation built into numpy was achieved by the following procedure. The following properties of each scan vector (\textbf{v}) is obtained:

  • start point (v_o)
  • normalised scan vector direction (\hat{\textbf{v}})
  • scan vector length \lVert \textbf{v} \rVert

For each LayerGeometry object, these are obtained across each scan vector as follows assuming that it is a HatchGeometry object:

# Calculate the length of the hatch vector and the direction
coords = layerGeom.coords.reshape(-1, 2, 2)
delta = np.diff(coords, axis=1).reshape(-1, 2)
lineDist = np.hypot(delta[:, 0], delta[:, 1]).reshape(-1, 1)

# Take the first coordinate
v0 = coords[:, 1, :].reshape(-1, 2)

# Normalise each scan vector direction
vhat = -1.0 * delta / lineDist

The number of point exposure across each vector is calculated by simple division. Unfortunately there is rounding, so the exposure may be missing at the end of the scan vector.

# Calculate the number of exposure points across the hatch vector based on its length
numPoints = np.ceil(lineDist / pointDistance).astype(np.int)

The total number of exposure points is calculated and is used to pre-populate some arrays, where the terms are used to extrapolate the exposure point distance (Line 16) are generated. Line 10 is key here and idxArray stores a range to account for each exposure point along the scan vector

Unfortunately, there isn’t a convenient way vectorise this so it is performed across each scan vector. Finally, the exposure points are generated by extrapolating them from v_o along \hat{\textbf{v}}.

# Pre-populate some arrays to extrapolate the exposure points from
totalPoints = int(np.sum(numPoints))
idxArray = np.zeros([totalPoints, 1])
pntsArray = np.zeros([totalPoints, 2])
dirArray = np.zeros([totalPoints, 2])

idx = 0
for i in range(len(numPoints)):
    j = int(numPoints[i])
    idxArray[idx:idx + j, 0] = np.arange(0, j)
    pntsArray[idx:idx + j] = p0[i]
    dirArray[idx:idx + j] = vhat[i]
    idx += j

# Calculate the hatch exposure points
hatchExposurePoints = pntsArray + pointDistance * idxArray * dirArray

This generates a set of exposure points, which can be seen in the figure below:

PySLM: Exposure points generated across the scan vectors between adjacent islands.
Exposure points plotted along each hatch vector are shown. Note the overlapping area between adjacent islands.

Once the exposure points are generated, the deposited energy (laserPower x pointExposureTime) is added to each point and stored as a third column for later use.

Exposure Map Generation

Generating the exposure map is trivial. For a selected resolution chosen by the user, the bitmap slice is used currently to get the image to process with. The exposure coordinates are simply mapped, given an offset and resolution onto the image. The deposited energy is then summed accordingly across each pixel. The energy per area (latex]Q[/latex]), is calculated by diving by the aerial coverage of a single pixel.

PySLM: A pseudo 'heatmap' or exposure Map showing the deposited energy for an island hatch region processed using SLM.
Exposure / Heat Map showing the deposited energy applied to the position and laser parameters used. Notice the overlap region between islands show a higher energy deposition.

The resolution is arbitrarily selected. The above figure is generated with a resolution of 0.2 mm. It can be seen in this example that the there is a greater energy deposition across the the overlap regions between the islands. There are some ‘aliasing’ effects due to discretion onto a finer grid which is dependent on the resolution chosen.

Conclusions

The above exposure map generation might have little scientific utility by itself. Rather it can offer a method to potentially visualise the energy deposition across a layer.

However, many thermo-mechanical AM build simulation that incorporate thermal effects in addition to the inherent strain approach could potentially incorporate a reference exposure map across the layer rather than assuming a single aerial heat-flux applied to the layer. This could actually be carried out across a ‘packet’ or number of layers to better correspond with the geometry and scan strategies employed.

Documented Example

The code for this can be found in examples/example_heatmap.py on the github repository.

Multi-threading Slicing & Hatching in PySLM

In PySLM, the slicing and hatching problem is inherently parallelisable simply because we can discretise areas the geometry into disrcete layers that for most situations behave independent. However, the actual underlying algorithms for slicing, offsetting the boundaries, clipping the hatch vectors is serial (single threaded). In order to significantly reduce the process time, multi-threaded processing is very desirable

Multi-threading in Python

Unfortunately, Python like most scripting or interpreter languages of the past are not inherently designed or destined to be multi-threaded. Perhaps, this may change in the future, but other scripting languages may fill this computational void (Rust, Julia, Go). Python, by intentions limits any multi-threaded use in scripts by using a construct known as the GIL – Global Interpreter Lock. This is a shared situation in other common scripting languages Matlab (ParPool), Javascript (Worker) where the parallel computing capability of multi-core CPUs cannot be exploited in a straightforward manner.

To some extent special distributions such as via the Anaconda distribution, core processing libraries such as numpy, scipy, scikit libraries internally are generally multi-threaded and support vectorisation via native CPU extensions. More computational mathematical operations and algorithms can to some extent be optimised to run in parallel automatically using numba, numexpr, and however, this cannot cover more broad multi-functional algorithms, such as those used in PySLM.

Python has the thread module for multi-threaded processing, however, for CPU bound processing it has very limited use This is because Python uses the global interpreter lock – GIL and this only allows one programming thread (i.e. one line) to be executed at any instance. It is mainly used for asynchronous IO (network or filesystem) operations which can be processed in the background.

Use of Multiprocessing Library in PySLM

The other option is to use the multiprocessing library built into the core Python distribution. Without going into too much of the formalities, multi-processing spawns multiple python processes and assign batches of work. The following programming article I found as a useful reference to the pitfalls of using the library.

In this implementation, the Pool and Manager modules are used to more optimally process the geometry. The most important section is to initialise the multiprocessing library with the ‘spawn‘ method, which stops random interruptions during the operation as discussed in the previous article.

from multiprocessing import Manager
from multiprocessing.pool import Pool
from multiprocessing import set_start_method

set_start_method("spawn")

The Manager.dict acts as a ‘proxy‘ object used to more efficiently store data which is shared between each process that is launched. Without using manager, for each process launch, a copy of the objects passed are made. This is important for the geometry or Part object, which if it were to contain a lattice of a set ofs complex surface would become expensive to copy.

d = Manager().dict()
d['part'] = solidPart
d['layerThickness'] = layerThickness # [mm]

A Pool object is used to create a set number of processes via setting the parameter processes=8 (typically one per CPU core). This is a fixed number re-used across a batch through the entire computation which removes the cost of additional copying and initialising many process instances. A series of z slice levels are created representing the layer z-id. These are then merged into a list of tuple pairs with the Manager dictionary and is stored in processList.

Pool.map is used to perform the slice function (calculateLayer) and collect all computed layers following the computation.

p = Pool(processes=8)

numLayers = solidPart.boundingBox[5] / layerThickness
z = np.arange(0, numLayers).tolist()

processList = list(zip([d] * len(z), z))

# Run the pro
layers = p.map(calculateLayer, processList)

The slicing function is fairly straightforward and just unpacks the arguments and performs the slicing and hatching operation. Note: each layer needs to currently initialise its own instance of a Hatcher class because this is not shared across all the processes. This carries a small cost, but means each layer can process entirely independently; in this example the change is the hatchAngle across layers. The layer position is calculated using the layer position (zid) and layerThickness.

def calculateLayer(input):
    # Typically the hatch angle is globally rotated per layer by usually 66.7 degrees per layer
    d = input[0]
    zid= input[1]

    layerThickness = d['layerThickness']
    solidPart = d['part']

    # Create a StripeHatcher object for performing any hatching operations
    myHatcher = hatching.Hatcher()

    # Set the base hatching parameters which are generated within Hatcher
    layerAngleOffset = 66.7
    myHatcher.hatchAngle = 10 + zid * 66.7
    myHatcher.volumeOffsetHatch = 0.08
    myHatcher.spotCompensation = 0.06
    myHatcher.numInnerContours = 2
    myHatcher.numOuterContours = 1
    myHatcher.hatchSortMethod = hatching.AlternateSort()

    #myHatcher.hatchAngle += 10

    # Slice the boundary
    geomSlice = solidPart.getVectorSlice(zid*layerThickness)

    # Hatch the boundary using myHatcher
    layer = myHatcher.hatch(geomSlice)

    # The layer height is set in integer increment of microns to ensure no rounding error during manufacturing
    layer.z = int(zid*layerThickness * 1000)
    layer.layerId = int(zid)

    return zid

The final step to use multiprocessing in Python is the inclusion of the python __name__ guard i.e:

if __name__ == '__main__':
   main()

The above is unfortunate because it makes debugging slightly more tedious in editors, but is the price for extra performance.

Performance Improvement

The performance improvement using the multiprocesssing library is shown in the table below for a modest 4 core laptop (my budget doesn’t stretch that far).

PySLM: A matplotlib showing hatching and slicing across multiple layers for the a cubic geometry using the python multi-processing library.
Matplotlib Figure showing every subsequent 10 layers of hatching for the geometry but is shown reduced scale.

This was performed on the examples/inversePyramid.stl geometry with an overall bounding box size [90 \times 90 \times 60] mm, hatch distance h_d=0.08mm and the layer thickness set at 40 μm.

Number of ProcessesRun Time [s]
Base-line (Simple For loop)121
1108
265.4
442
637.1
831.8
Approximate timings for 4 core CPU i7 Processors Using the Multi-Processing Library.

Given these are approximate timings, it is nearly linear performance improvement for the simple example. However, it can be seen choosing more processes beyond cores, does squeeze some extra performance out – perhaps due to Intel’s hyperthreading. Below shows that the CPU is fully utilised.

PySLM: Multi-threading options

Conclusions:

This post shows how one can adapt existing routines to generate multi-processing slicing and hatching with PySLM. In the future, it is desirable to explore a more integrated class structure for hooking functions onto. Other areas that are of interest to explore are potentially the use of GPU computing to parallelise some of the fundamental algorithms.

Example

The example showing this that can be run is example_3d_multithread.py in the github repo.

Improving Performance of Island Checkerboard Scan Strategy Hatching in PySLM

The hatching performance of PySLM using ClipperLib via PyClipper is reasonably good considering the age of the library using the Vatti polygon clipping algorithm. Without attempting to optimise the underlying library and clipping algorithm for most scenarios, the hatch clipping process should be sufficient for most use case. Future investigation will explore alternative clipping algorithms to further improve the performance of this intensive computational process

For the unfamiliar with the basic hatching process of a single layer, the laser or electron beam (a 1D single point source) must scan across an aerial (2D) region. This is done by creating a series of lines/vectors which infill or raster across the surface.

The most basic form of hatch infill for bulk regions is an alternating, meander, or in some locales referred to a serpentine scan strategy. This tends to be undesirable in SLM due to the creation of localised heat build-up [1] resulting in porosity, poor surface finish [2], residual stress and resultant distortion and anisotropy due to preferential grain growth [3]. Stripe or Island scan strategies are employed in attempt to mitigate these by limiting the length of scan vectors used across a region [4][5][6]. Within the layer hatch vectors for each island are oriented orthogonal to each other and the scan vector length can be precisely controlled in order to reduce the magnitude of residual stresses generated [7].

However, when the user desires a stripe or an island scan strategy, the number of clipping operations for the individual hatch vectors increases drastically. The increase in number of clipping operations increases due to division of the area into fixed size regions corresponding to the desired scan vector length (typically 5 mm)]:

  • Standard Meander Scan Strategy: n_{clip} \propto \frac{A}{hatchDistance(h_d)}
  • Stripe Scan Strategy: n_{clip} \propto \frac{A}{StripeWidth}
  • Island Scan Strategy: n_{clip} \propto \frac{A}{IslandWidth^2}

As can be observed, the performance of hatching with an island scan strategy degrades rapidly when using the island scan due to reciprocal square. As a result, using a naive approach, hatching a very large planar region using an island scan strategy could quickly result in 100,000+ clipping operations for a single layer for a large flat. In addition, this is irrespective of the sparsity of the layer geometry. The way the hatch filling approach works in PySLM, the maximum extent of a contour/polygon region is found. A circle is projected based on this maximum extent, and an outer bounding box is covered. This is explained in a previous post.

The scan vectors are tiled across the region. The reason behind this is to guarantee complete coverage irrespective of the chosen hatch angle, \theta_h, across the layer and largely simplifies the computation. The issue is that many regions will be outside the boundary of the part. Sparse regions both void and solid will not require additional clipping.

The Proposed Technique:

In summary, the proposed technique takes advantage that each island is regular, and therefore each island can be used to discretise the region. This can be used to perform intersection tests for region that may be clipped, whilst recycling existing hatch vectors for those within the interior boundary.

Given that use an island scan strategy provides essentially structured grid, this can be easily transformed into a a method for selecting regions. Using the shapely library, each island boundary consisting of 4 edges can be quickly tested to check if it overlaps internally with the solid part and also intersected with the boundary. This is an efficient operation to perform, despite shapely (libGEOS) being not as efficient as PyClipper.

from shapely.geometry.polygon import LinearRing, Polygon

intersectIslands = []
overlapIslands = []

intersectIslandsSet = set()
overlapIslandsSet= set()

for i in range(len(islands)):
    
    island = islands[i]
    s = Polygon(LinearRing(island[:-1]))

    if poly.overlaps(s):
        overlapIslandsSet.add(i) # id
        overlapIslands.append(island)

    if poly.intersects(s):
        intersectIslandsSet.add(i)  # id
        intersectIslands.append(island)


# Perform difference between the python sets
unTouchedIslandSet = intersectIslandsSet-overlapIslandsSet
unTouchedIslands = [islands[i] for i in unTouchedIslandSet]

This library is used because the user may re-test the same polygon consecutively, unlike re-building the polygon state in ClipperLib. Ultimately, this presents three unique cases:

  1. Non-Intersecting (shapely.polygon.intersects(island) == False) – The Island resides outside of the boundary and is discarded,
  2. Intersecting (shapely.polygon.intersects(island) == True) – The Island is in an internal region, but may be also clipped by the boundary,
  3. Clipped (shapely.polygon.intersects(island) == True) – The island intersects with the boundary and requires clipping.

PySLM - Clipping of island regions when generating Island Scan Strategies for Selective Laser Melting
The result is shown here for a simple 200 mm square filled with 5 mm islands:

Taking the difference between cases 2) and 3), the islands with hatch scan vectors can be generated without requiring unnecessary clipping of the interior scan vectors. As a result this significantly reduces the computational effort required.

Although extreme, the previous example generated a total number of 2209 5 mm islands to cover the entire region. The breakdown of the island intersections are:

  1. Non-intersecting islands: 1591 (72%),
  2. Non-clipped islands: 419 (19%),
  3. Clipped islands: 199 (9%).

With respect to solid regions, the number of clipped islands account for 32% of the total area. The overall result is shown below. The total area of the hatch region that was hatched is 1.97 \times 10^3 \ mm^2, which is equivalent to a square length of 445 mm, significantly larger than what is capable on most commercial SLM systems. Using an island size of 5 mm with an 80 μm hatch spacing, the approximate hatching time is 6.5 s on a modest laptop system. For this example, 780 000 hatch vectors were generated.

PySLM - A close-up view showing the clipped scan vectors using the more efficient island scan strategy.
A close up view showing the 5mm Island Hatching with 0.8 mm Hatch Distance. Blue Lines show the overall path traversed by the laser beam. The total time taken for hatching was approximately 8 seconds.

The order of hatching scanned is shown by the blue lines, which trace the midpoints of the vectors. Hatches inside the island are scanned sequentially. The order of scanning in this case is chosen to go vertically upwards and then horizontally across using the in-built Python 3 sorting function with a lambda expression Remarkably, all performed using one line:

sortedIslands = sorted(islands, key=lambda island: (island.posId[0], island.posId[1]) )

A future post will elaborate further methods for sorting hatch vectors and island groups.

Comparison to Original Implementation:

The following is a non-scientific benchmark performed to illustrate the performance profile of the proposed method in PySLM.

Island Size [mm]Original Method Time [s]Proposed Method Time [s]
34665.3
52586.5
101217.9
20758.23
Approximate benchmark comparing Island Hatching Techniques in PySLM

It is clearly evident that the proposed method reduces the overall time by 1-2 orders for hatching a region. What is strange is that with the new proposed method, the overall time increases with the island size.

Generally it is expected that the number of clipping operations n_{clip} to be the following:

n_{clip} \propto \frac{Perimiter}{IslandWidth}

Potentially, this allows bespoke complex ‘sub-island’ scan strategies to be employed without a significant additional cost because scan vectors within un-clipped island regions can be very quickly replicated across the layer.

Other Benefits

The other benefits of taking approach is making a more modular object orientated approach for generating island based strategies, which don’t arbitrarily follow regular structured patterns. A future article will illustrate further explain the procedures for generating these.

The example can be seen and run in examples/example_island_hatcher.py in the Github repository.

References

References
1 Parry, L. A., Ashcroft, I. A., & Wildman, R. D. (2019). Geometrical effects on residual stress in selective laser melting. Additive Manufacturing, 25. https://doi.org/10.1016/j.addma.2018.09.026
2 Valente, E. H., Gundlach, C., Christiansen, T. L., & Somers, M. A. J. (2019). Effect of scanning strategy during selective laser melting on surface topography, porosity, and microstructure of additively manufactured Ti-6Al-4V. Applied Sciences (Switzerland), 9(24). https://doi.org/10.3390/app9245554
3, 4 Zhang, W., Tong, M., & Harrison, N. M. (2020). Scanning strategies effect on temperature, residual stress and deformation by multi-laser beam powder bed fusion manufacturing. Additive Manufacturing, 36(June), 101507. https://doi.org/10.1016/j.addma.2020.101507
5 Ali, H., Ghadbeigi, H., & Mumtaz, K. (2018). Effect of scanning strategies on residual stress and mechanical properties of Selective Laser Melted Ti6Al4V. Materials Science and Engineering A, 712(October 2017), 175–187. https://doi.org/10.1016/j.msea.2017.11.103
6 Robinson, J., Ashton, I., Fox, P., Jones, E., & Sutcliffe, C. (2018). Determination of the effect of scan strategy on residual stress in laser powder bed fusion additive manufacturing. Additive Manufacturing, 23(February), 13–24. https://doi.org/10.1016/j.addma.2018.07.001
7 Mercelis, P., & Kruth, J.-P. (2006). Residual stresses in selective laser sintering and selective laser melting. Rapid Prototyping Journal, 12(5), 254–265. https://doi.org/10.1108/13552540610707013

Optimising and Improving PyClipper (ClipperLib) Hatch Clipping Performance for AM

In previous incarnations throughout my time, the renowned ClipperLib has been used to provide the polygon and line clipping and offsetting operations for generating the hatching used in Selective Laser Melting . This library is also used as part of the slicing engine in the popular Cura software used in the Ultimakr filament 3D printers. his opensource c++ polygon clipping library at present is very efficient, robust and arguably does its specific job very well.

The First Attempt

The first implementation was done through a custom mex extension for Matlab – infact there is actually now a mex extension available on mathworks. The clipping of hatch vectors with a polygon boundary was done on an individual basis; one hatch vector at a time. This is necessary because returned clipped lines in ClipperLib are not guaranteed to be returned in the same sequential order as they were passed due to the underlying Vatti algorithm which works using a ‘scan-line’ approach.

Unfortunately, this was very inefficient operation requiring the re-initialisation of the boundary polygon for each clipped line. Through bruteforce and persistence, this method did a satisfactory job hatching simple geometries for research purposes.

The Second Attempt

The next progression was transitioning to Python, using PyClipper. Appreciating the challenges of trying to run hatching on embedded hardware, clipping each hatch vector didn’t fair well. The second approach was to clip all the hatch vectors together to remove the overhead of re-initializing the ClipperLib state. After clipping the hatch vectors were sorted using the following method.

The inner product or dot project between the midpoint of each hatch vector and the effective vector of the hatch angle is made. The projected distance is the used to sort the relative line position. This can be seen in LinearSort class or in the following code excerpt

theta_h = np.deg2rad(hatchAngle)

# Find the unit vector normal based on the hatch angle
norm = np.array([np.cos(theta_h), np.sin(theta_h)])

midPoints = np.mean(scanVectors, axis=1)

# Project the midpoint onto the hatch angle vector
idx2 = norm.dot(midPoints.T)

# Sort the scan vectors
idx3 = np.argsort(idx2)
sortIdx = np.arange(len(midPoints))[idx3]

return scanVectors[sortIdx]

This worked fairly well. However, it has severe drawback in that it is limited to just hatching and sorting across a single region and relies on a constant hatch direction. The sorting mechanism severely constrains what approaches can be used for scan strategy design, i.e. island scan strategies.

The current approach in PySLM

The current implementation builds on a clever trick built into ClipperLib that unfortunately is not implemented in PyClipper. Within ClipperLib there is an option to include the z coordinate as part of the Point Struct by defining the use_xyz preprocessor directive. This allows the user to pass an additional 4 byte integer for each point. Infact, this can represent any information, and this can infact be the unique hatch id, which defines the order of scanning. Note: We could further modify the library to use any data-type but this is unnecessary.

Each hatch vector coordinate is given an individual id.During the clipping operation with the boundary, the intersection point is assigned the id from the original ‘z’ coordinate.

During clipping, the z-component is not actually used for clipping. ClipperLib does not know which z-component should be used at the intersection between edges; the boundary or the hatch line. A function maxZFillFunc is defined in ClipperLib.cpp which resolves this ambiguity. To resolve this, a specific function finds the maximum z component from all intersecting edge points, and stores it in the new point pt generated. As a result the original index used for sorting is stored in the clipped result.


void maxZFillFunc(const IntPoint& e1bot, IntPoint& e1top, IntPoint& e2bot, IntPoint& e2top, IntPoint& pt) {
        // Find the maximum z value from all points. This provides a non ambiguious cases, as for PySLM the background
        // contour has a value of zero.
        long maxZ = -1;
        if(e1bot.Z > maxZ)
            maxZ = e1bot.Z;

        if(e1top.Z > maxZ)
            maxZ = e1top.Z;

        if(e2top.Z > maxZ)
            maxZ = e2top.Z;

        if(e2bot.Z > maxZ)
            maxZ = e2bot.Z;

        // Assign the z value to pt
        pt.Z = maxZ;
  };

void Clipper::ZFillFunction(ZFillCallback zFillFunc)
{
  m_ZFill = zFillFunc;
}

Further changes in the cython pyrex definition – pyclipper.pyx was needed to enable this functionality, hence, this needs to be compiled internally for this project.

Back in Python, the order index is assigned to each pairs of coordinates representing a hatch vector, resulting in an 2n\times3 array where n is the number of hatches:

idx = np.arange(len(x) / 2)
hatchArray = np.hstack([x,y,id])

After clipping the hatch vectors generated requiring a quick sort based on their id. Typically this isn’t expensive because most of the clipped hatch vectors are returned in the same order.

clippedLines = self.clipperToHatchArray(clippedPaths)

# Extract only x-y coordinates and sort based on the pseudo-order stored in the z component.
clippedLines = clippedLines[:, :, :3]
id = np.argsort(clippedLines[:, 0, 2])
clippedLines = clippedLines[id, :, :]

Conclusions

Throughout PySLM, the user simply has to generate a sequence of hatch vectors as the same order they wish them to be scanned. An increasing unique index is applied to each hatch vector.

This proposed technique is remarkbly simple, efficient and effective for clipping . It is very trivial but given the limit number of polygon clipping libraries available this was particularly useful.

In the future, I plan to explore storing additional data associated with the hatch vectors. This itself may not be necessary because this data can be looked up from the global hatch data array, but it may result in simplification of the code.

Basics of Hatching in PySLM

The generation of hatch infills in PySLM is relatively straightforward and can be adapted to suit the user’s needs. Anecdotally, this was a remarkable improvement on performance and simplicity over preexisting code that was used previously during University, which used a combination of Matlab and a ClipperLib library.

The generation of hatch vectors are nearly always a series of adjacent parallel scan vectors that an infill a bounded area. There has been some unusual approaches for infills that have been tried such as spiral [1][2] and fractal [3][4][5] scan strategies. Generally, exploring different raster patterns across an SLM part remains largely unexplored due to the unavailability of opensource tools, code and open access to machine platforms now in research environments. However, potentially this could unlock the ability to micro-structural development in some metal alloys. Attard et al. demonstrates this by modifying the island size throughout the part to promote different microstructures [6].

The first step of generating the hatch is finding the region to ensure scan vector generated guarantees full coverage across the boundary for the chosen hatch angle, h_d. This is needed irrespective of scan strategy. My previous approach was rather cumbersome and required trigonometry to be used to calculate the start and end points according the chosen (\theta_h) to cover the boundary. A far simpler approach was conceived.

Basic method:

The approach simply relies on covering the entire boundary irrespective of hatch angle and letting the polygon clipping library do the actual heavy lifting.

In this method, the bounding box for a closed region is found. A rather simple calculation made by taking the min and max of the boundary coordinates – easily obtained via shapely.

The maximum extent or diagonal length, r, from the bounding box corner to the bounding centre is then calculated. In the diagram, it is observed that this forms a circular locus, which upon choosing any hatch angle, all hatch lines after rotation will be guaranteed to fully cover boundary. Remarkably simple.

The basic code can be found in hatching.py:

# Expand the bounding box
bboxCentre = np.mean(bbox.reshape(2,2), axis=0)

# Calculates the diagonal length for which is the longest
diagonal = bbox[2:] - bboxCentre
bboxRadius = np.sqrt(diagonal.dot(diagonal))

From the locus, we assume the hatch coordinate local coordinate system with an origin (-r,-r). A series of parallel hatch vectors with hatch distance, h_d, can be generated conveniently covering the entire square region . This can be done using numpy, using np.tile:

# Construct a square which wraps the radius
x = np.tile(np.arange(-bboxRadius, bboxRadius, hatchSpacing, dtype=np.float32).reshape(-1, 1), (2)).flatten()
y = np.array([-bboxRadius, bboxRadius]);
y = np.resize(y, x.shape)

Unfortunately, ClipperLib by default does not guarantee the order lines are clipped in. The naive approach is to clip these individually in a sequential order, which introduces a significant overhead. Alternatively the clipped hatch lines can be sorted. For achieving this, the PyClipper library was customized to guarantee the sequential order of the scan , which is explained in this post.

The approach to guarantee the z-order is to provide a z-coordinate which defines the order of scanning. This should be in pairs of ascending order which group the coordinates in each scan vector: i.e. 0,0, 1,1,..,2,2 etc. Finally before being clipped, an n\times3 coordinate matrix is formed.

# Generate a linear sorting according to the order of scanning.
# i.e. 0,0,1,1,..,n-1,n-1,n,n 
z = np.arange(0, x.shape[0] / 2, 0.5).astype(np.int64)

# Group the XY coordinates and z order 
coords = np.hstack([x.reshape(-1, 1),
                    y.reshape(-1, 1),
                    z.reshape(-1, 1)])

From this, a 2D rotation matrix, R(\theta_h), and translation to the bounding box centre can be applied to the coordinates of the unclipped hatch lines.

# Create the 2D rotation matrix
c, s = np.cos(theta_h), np.sin(theta_h)
R = np.array([(c, -s, 0),
              (s, c, 0),
              (0, 0, 1.0)])

# Apply the rotation matrix and translate to bounding box centre
coords = np.matmul(R, coords.T)
coords = coords.T + np.hstack([bboxCentre, 0.0])

There is one issue with the method is that using region / island based scan strategies. Many of the regions are not clipped, therefore, for large regions it can become inefficient. However, assuming a normal hatch in-fill is used, the clipping operation mostly amortises the additional cost.

References

References
1 Zhang, W., Tong, M., & Harrison, N. M. (2020). Scanning strategies effect on temperature, residual stress and deformation by multi-laser beam powder bed fusion manufacturing. Additive Manufacturing, 36(June), 101507. https://doi.org/10.1016/j.addma.2020.101507
2 Qian, B., Shi, Y. S., Wei, Q. S., & Wang, H. B. (2012). The helix scan strategy applied to the selective laser melting. International Journal of Advanced Manufacturing Technology, 63(5–8), 631–640. https://doi.org/10.1007/s00170-012-3922-9
3 Ma, L., & Bin, H. (2006). Temperature and stress analysis and simulation in fractal scanning-based laser sintering. The International Journal of Advanced Manufacturing Technology, 34(9–10), 898–903. https://doi.org/10.1007/s00170-006-0665-5
4 Catchpole-Smith, S., Aboulkhair, N., Parry, L., Tuck, C., Ashcroft, I., & Clare, A. (2017). Fractal scan strategies for selective laser melting of ‘unweldable’ nickel superalloys. Additive Manufacturing, 15, 113–122. https://doi.org/10.1016/j.addma.2017.02.002
5 Sebastian, R., Catchpole-Smith, S., Simonelli, M., Rushworth, A., Chen, H., & Clare, A. (2020). ‘Unit cell’ type scan strategies for powder bed fusion: The Hilbert fractal. Additive Manufacturing, 36(July), 101588. https://doi.org/10.1016/j.addma.2020.101588
6 Attard, B., Cruchley, S., Beetz, C., Megahed, M., Chiu, Y. L., & Attallah, M. M. (2020). Microstructural control during laser powder fusion to create graded microstructure Ni-superalloy components. Additive Manufacturing, 36, 101432. https://doi.org/10.1016/j.addma.2020.101432

Slicing and Hatching for Selective Laser Melting (L-PBF)

Much of slicing and hatching process is already taken for granted in commercial software mostly offered by the OEMs of these systems rarely discussed amongst academic research. Already we observe practically the implications direct control over laser parameters and scan strategy on the quality of the bulk material – reduction in defects, minimising distortion due to residual stress, and the surface quality of parts manufactured using these process. Additionally, it can have a profound impact the the metallic phase generation, micro-structural texture driven via physics-informed models [1], grading of the bulk properties and offer precise control over manufacturing intricate features such as thin-wall or lattice structures [2].

This post hopefully highlights to those unfamiliar some of the basis process encountered in the generation of machine build files used in AM systems and get a better understanding to the operation behind PySLM. I have tried my best to generalise this as much as possible, but I imagine there are subtleties I have not come across.

This post is to provide some reference into the generation of hatches or scan vectors are created for use in AM processes such as selective laser melting (SLM), which uses a point energy source to raster across a medium. Some people prefer to more generally to classify the family of processes using the technical ASTM F42 committee standards 52900 and 52911 – Powder Bed Fusion (PBF). I won’t go into the basic process of the manufacturing processes such as EBM, SLM, SLA, BJF, as there are many excellent articles already that explain these in far greater detail.

Machine Build Files

AM processes require a digital representation to manufacture an object. These tend to be computed offline – separate from the 3D Printer, using specialist or dedicated pre-processing software. I expect this will become a closed-loop system in the future, such that the manufacturing integrated directly into the machine.

For some AM process families, the control operations may be exceedingly granular – i.e. G-code. G-code formats state specific instructions or functional commands for the 3D printer to sequentially or linearly execute. These tend to fit with deposition methods such as Filament Extrusion, Direct-Ink-Writing (robo-casting) and direct energy deposition (DED) methods. Typically, these tends to be for deposition with a machine systems, which requires coordination of physical motion in-conjunction with some mechanised actuation to deposit/fuse material.

Machine Build File Formats for L-PBF

For exposure (laser, electron-beam) based AM processes, commercial systems use a compact notation solely for representing the scan path the exposure source will traverse . The formats are often binary to aid their compactness.

To summarise, within these build files, an intermediate representation consists of index-based referenceable parameters for the build. The remainder consists of a series of layers, that contain geometric entities (points, vectors) that are used to to control the exposure for the border or contour or raster or infill the interior region. For L-PBF processes, the digital files, commonly referred as “machine build file” comes in various flavours dependent on the machine manufacture:

  • Renshaw .mtt,
  • SLM Solution .slm,
  • DMG Mori Realizer .rea
  • EOS .sli
  • Aconity .cli+ or .ilt wrapper

Some file formats, such as Open Beam Path format can specify bezier curves [3]. Another recently proposed open source format created by RWTH Aachen in 2022 called OpenVector Format based on Google’s Protobuf schema. The format aims to offer a specification universally compatible across a swathe of PBF processes and supplement existing commercial formats with additional build-process meta-data (e.g. build, platform temperature, dosing) and detailed definition with further advancements in the process, such as multi-beam builds.

Build-File Formats

Higher level representations that describe the distribution of material(s) defining geometry – this could be bitmap slices or even a 3D model. Processes such as Jetting, BJF, High Speed Sintering, DLP Vat-polymerisation currently available offer this a reality. With time, polymer and metal processes will evolve to become 2D:, diode aerial melting [4] or more aerial based scanning based on holographic additive manufacturing methods, such as those proposed by Seurat AM [5] based off research at LLNL, and recently at University of Cambridge [6] . In the future, we can already observe the exciting prospect of new processes such as computed axial lithograph [7] that will provide us near instantaneous volumetric additive manufacturing.

For now, single and multi point exposure systems for the imminent future will remain with us as the currently available process. PySLM uses an intermediate representation – specifying a set of points and lines to control the exposure of energy into a layer.

The Slicing and Hatching Process in L-PBF

With nearly most conventional 3D printing process, it begins with a 3D representation of a solid volume or geometry. 2D planar slices or layers are extracted from a 3D mesh or B-Rep surface in CAD by taking cross-sections from a geometry. Each slice layer consist of a set of boundaries and holes describing the cross-section of an object. Note: non-planar deposition does exist for DED/Filament processes, such as this Curved Layer Fused Deposition Modeling [ref] and a spherical slicing technique [8].

For consolidating material, an exposure beam must raster across the surface medium (metal or polymer powder, or a photo-polymer resin) depending on the process. Currently this is a single or multiple point which moves at a velocity vwith a power P across the surface. The designated exposure or energy deposited into the medium is principally a function of these two parameters, depending on the type of laser:

  • (Quasi)-Continious Wave: The laser remains active switched on (typically modulated using a form of PWM) across the entire length of the scan vector
  • Pulsed Mode (Q-Switched): Laser is pulsed at set distances and exposure times across the scan vector

Numerous experiments often tend to result in parametric power/speed maps to the achieved part bulk density, that result in usually optimal processing windows that produce stable and consistent melt-tracks [9][10]. Recently, process maps are based on a non-dimensional parameter such as the normalised enthalpy approach, that more reliably assist selecting a suitable process windows [11].

Illustration of a scan vector commonly used in Laser Powder-Bed Fusion (SLM)

However, the complexity of the process extends further and is related to a many additional variables dependent on the process such as layer thickness, absorption coefficient (powder and material), exposure beam profile etc.. Additionally, the cumulative energy deposited spatially over a period of time must consider overlap of scan vectors within an area.

Scan Vector Generation

Each boundary polygon is offset initially to account for the the radius of the beam exposure, which is termed a ‘spot compensation factor‘. Some processes such as SLS or BJF account for global part shrinkage volumetrically throughout the part by having a global scale factor or deformed mesh to compensate to non-uniform shrinkage across the part.

The composition of laser scan vectors used in a slice or layer for L-PBF or Selective Laser Melting. The boundary is offset multiple times, with the interior or core filled with hatch vectors.
The typical composition of a layer used for scanning in exposure based processes. This consists of outer and inner contours, with the core interior filled with hatches.

This first initial offset is the outer-contour which would be visible on the exterior of the part. This contour will have a different set of laser parameters in order to optimise and improve the surface roughness of the part obtained. A further offset is applied to generate a set of inner-contours before hatching begins.

Depending on the orientation of the surface (e.g. up-skin or down-skin), the boundary and interior region may be intersected to fine-tune the laser parameters to provide better surface texture, or surface roughness – typically varying between Ra = 3-13 μm [12] primarily determined by the surface angle and a combination of the process variables including,

  • the powder feedstock (bulk material, powder size distribution)
  • laser parameters
  • layer thickness (pre-dominantly fixed or constant for most AM processes)

Overhang regions and surfaces with a low overhang angles tend to be susceptible to high surface-roughness. Roller re-coater L-PBF systems – available only on 3DSystems or AddUp system,, tend to offer far superior surface quality on low inclined or overhang regions. Additionally, progressive advancement and maturity of laser parameter optimisation, and those computationally driven using part geometry [13] are able to further enhance the quality and potentially eliminate the need for support structures. Depending on the machine platform, these regions are identified by sampling across two-three layers. Overhang regions obviously require support geometry, which is an entirely different topic discussed in this post.

Laser parameters in SLM (L-PBF) can be optimised based on the adjacent surface regions. Special regions, include the upskin, downskin and overhang regions
Laser parameters can be optimised based on the adjacent surface regions. Special regions, include the upskin, downskin and overhang regions needed to improve the surface roughness and reduce density in regions.

Following the generation of the contours, the inner core region requires filling with hatches. Hatches are a series of parallel scan vectors placed adjacent at a set hatch distance, h_d. This parameter is optimized according to the material processed, but is essentially related to the spot radius of the exposure point r_s in order to reduce inter-track and inter layer porosity. Across each layer these tend to be placed at a particular orientation \theta_h, which is is then incrementally rotated globally for subsequent layers, typically 66.6°. This rotation aims to smooth out the build process in order to minimise inter-track porosity, and generate homogeneous material, and in the case of SLM mitigate the effects of anisotropic residual stress generation.

The composition and terminology (hatch distance, hatch spacing, hatch angle) used in L-PBF. The Layer Geometry objects used to scan across a Layer in Selective Laser Melting (L-PBF). The various parameters such as the hatch distance and hatch angle are shown.
A general composition of the various LayerGeometry objects used to scan across a Layer. The various parameters such as the hatch distance, spacing and hatch angle are shown.

The distribution (position, length, rotation) of these hatch vectors are arranged using a laser scan strategy. The most common include a simple alternating hatch, stripe and island or checkerboard scan strategy.

Each set or group of scan vectors is stored together in a LayerGeometry, depending on the type (either a set of point exposures, contour or hatch vectors). These LayerGeometry groups usually share a set of exposure parameters – power, laser scan speed (point exposure time, point distance for a pulsed laser), focus position).

Some systems offer a greater degree of control and can control individual power across the scan vectors. Other can fine tune the acceleration and modulate the power along the scan vectors to support techniques known as ‘skywriting‘. For instance in SLM, it has been proposed that careful tuning of the laser parameters towards the end of the scan vector, i.e. turning can reduce porosity by preventing premature collapse of key holing phenomena [14]. In theory, PySLM could be extended to provide greater control of the electro-optic systems used in the process if so desired.

Hopefully, this provides enough background for those who are interested and engaged in working with developing scan strategies and material development using PySLM.

References

References
1 Plotkowski, A., Ferguson, J., Stump, B., Halsey, W., Paquit, V., Joslin, C., Babu, S. S., Marquez Rossy, A., Kirka, M. M., & Dehoff, R. R. (2021). A stochastic scan strategy for grain structure control in complex geometries using electron beam powder bed fusion. Additive Manufacturing46. https://doi.org/10.1016/j.addma.2021.102092
2 Ghouse, S., Babu, S., van Arkel, R. J., Nai, K., Hooper, P. A., & Jeffers, J. R. T. (2017). The influence of laser parameters and scanning strategies on the mechanical properties of a stochastic porous material. Materials and Design131, 498–508. https://doi.org/10.1016/j.matdes.2017.06.041
3 Open Beam Path – Freemelt, https://gitlab.com/freemelt/openmelt/obplib-python
4 Zavala Arredondo, Miguel Angel (2017) Diode Area Melting Use of High Power Diode Lasers in Additive Manufacturing of Metallic Components. PhD thesis, University of Sheffield.
5 Seurat AM. https://www.seuratech.com/
6 https://www.theengineer.co.uk/holographic-additive-manufacturing-lasers/
7 Kelly, B., Bhattacharya, I., Shusteff, M., Panas, R. M., Taylor, H. K., & Spadaccini, C. M. (2017). Computed Axial Lithography (CAL): Toward Single Step 3D Printing of Arbitrary Geometries. Retrieved from http://arxiv.org/abs/1705.05893
8 Yigit, I. E., & Lazoglu, I. (2020). Spherical slicing method and its application on robotic additive manufacturing. Progress in Additive Manufacturing, 5(4), 387–394. https://doi.org/10.1007/s40964-020-00135-5
9 Yadroitsev, I., & Smurov, I. (2010). Selective laser melting technology: From the single laser melted track stability to 3D parts of complex shape. Physics Procedia, 5(Part 2), 551–560. https://doi.org/10.1016/j.phpro.2010.08.083
10 Maamoun, A. H., Xue, Y. F., Elbestawi, M. A., & Veldhuis, S. C. (2018). Effect of selective laser melting process parameters on the quality of al alloy parts: Powder characterization, density, surface roughness, and dimensional accuracy. Materials, 11(12). https://doi.org/10.3390/ma11122343
11 Ferro, P., Meneghello, R., Savio, G., & Berto, F. (2020). A modified volumetric energy density–based approach for porosity assessment in additive manufacturing process design. International Journal of Advanced Manufacturing Technology, 110(7–8), 1911–1921. https://doi.org/10.1007/s00170-020-05949-9
12 Ni, C., Shi, Y., & Liu, J. (2019). Effects of inclination angle on surface roughness and corrosion properties of selective laser melted 316L stainless steel. Materials Research Express, 6(3). https://doi.org/10.1088/2053-1591/aaf2d3
13 Velo3D Sapphire Printer – SupportFree Technology. https://blog.velo3d.com/blog/supportfree-what-does-it-mean-why-is-it-important
14 Martin, A. A., Calta, N. P., Khairallah, S. A., Wang, J., Depond, P. J., Fong, A. Y., … Matthews, M. J. (2019). Dynamics of pore formation during laser powder bed fusion additive manufacturing. Nature Communications, 10(1), 1–10. https://doi.org/10.1038/s41467-019-10009-2