4. Clouds, meshes, scalar fields: introspection, manipulation

4.1. Clouds and scalar fields

4.1.1. cloud introspection

3D Clouds contains a set of indexed points and have several optional associated features: normals, several scalarFields, color, an octree structure, sensors…

cloud = cc.loadPointCloud("CloudComPy/Data/dataSample_5.0.xyz")

name = cloud.getName()        # 'dataSample_5 - Cloud'
cloud.setName("sampleCloud")

nbPoints = cloud.size()       # 1000000
pt = cloud.getPoint(12345)    # (-4.880000114440918, -1.5499999523162842, 0.08818499743938446)

bb = cloud.getOwnBB()         # bounding box
minCorner = bb.minCorner()    # (-5.0, -5.0, -1.0861680507659912)
maxCorner = bb.maxCorner()    # (4.989999771118164, 4.989999771118164, 5.0)

g = cloud.computeGravityCenter()

Several methods based on neighbourood computations require a radius parameter, often estimated with an average neighbour points.

radius = cc.GetPointCloudRadius(clouds=[cloud], nodes=10)

The octree is used in a lot of methods, to speed up access to the neighbourhood of the points, and it is often computed automatically.

octree = cloud.getOctree()
if octree is None:
    octree = cloud.computeOctree()

Access to scalarFields: ScalarFields are accessible by the getScalarField() method from their index (integer). It is sometimes more convenient to find them by their name from the scalarFields dictionary, which gives the index from the name.

nbSF = cloud.getNumberOfScalarFields()
if cloud.hasScalarFields()
    dic = cloud.getScalarFieldDic()

cloud.exportCoordToSF(True, True, True)
dic = cloud.getScalarFieldDic()  # {'Coord. X': 0, 'Coord. Y': 1, 'Coord. Z': 2}

SFX = cloud.getScalarField(dic['Coord. X'])
sfmin = SFX.getMin()             # -5.0
sfmax = SFX.getMax()             # 4.989999771118164
mean, variance = SFX.computeMeanAndVariance()
SFX.getValue(12345)              # -4.880000114440918

The cloud can have a color field that can be seen as a special type of scalarField. There are methods to convert scalarField to colors and vice-versa.

if not cloud.hasColors():
    cloud.setCurrentScalarField(2)
    cloud.convertCurrentScalarFieldToColors()

Normals can be computed with computeNormals():

if not cloud.hasNormals():
    cc.computeNormals([cloud])

If you need to save clouds for reopening with CloudCompare GUI, with a predefined state of what is shown (colors, normals, scalar fields), use the .bin format and define the state with the following functions:

4.1.2. cloud transformations

Basic transformations translate() and scale() allow to translate a cloud or rescale it, with separate factors along the 3 directions and an optional center (see test001.py).


# --- get a reference before transformation
g = cloud.computeGravityCenter()

print("gravityCenter: (%14.7e, %14.7e, %14.7e)" % (g[0], g[1], g[2]))
if not isCoordEqual(g, (-4.9999999e-03, -4.9999999e-03, 9.6193114e-03)):
    raise RuntimeError

# --- rescale the cloud along Z direction, center at origin

cloud.scale(1.0, 1.0, 2.0, (0., 0., 0.))

g = cloud.computeGravityCenter()
print("gravityCenter: (%14.7e, %14.7e, %14.7e)" % (g[0], g[1], g[2]))
if not isCoordEqual(g, (-4.9999999e-03, -4.9999999e-03, 1.9238623e-02)):
    raise RuntimeError

# --- inverse scaling

cloud.scale(1, 1, 0.5)

# --- a simple translation

cloud.translate((-1, -2, -3))

g = cloud.computeGravityCenter()
print("gravityCenter: (%14.7e, %14.7e, %14.7e)" % (g[0], g[1], g[2]))
if not isCoordEqual(g, (-1.0050000e+00, -2.0050001e+00, -2.9903808e+00)):
    raise RuntimeError

# --- inverse translation

cloud.translate((1, 2, 3))

More complex transformations, based on rotation, require a cloudComPy.ccGLMatrix object to define the transformation, which is applied to the cloud with the method applyRigidTransformation().

The method initFromParameters() allows to define the transformation with rotations (see Euler angles). The following code is extract from test026.py.

Define a transformation from a rotation angle and axis, plus a translation, with initFromParameters().


tr1 = cc.ccGLMatrix()
tr1.initFromParameters(0.47, (1.5, 2.9, 6.3), (5,3,-2))

Sometimes, you need to store or retreive the transformation matrix data (4x4 matrix stored in major column order), use toString(), fromString() and data().


m_data_str = tr1.toString()
m_data = tr1.data()
tr1b = cc.ccGLMatrix.fromString(m_data_str)
m_datab = tr1b.data()
for i in range(16):
    if m_data[i] != m_datab[i]:
        raise RuntimeError

If you need to keep precision while iterating on transformations, for instance, it is better to work with double precision transformations cloudComPy.ccGLMatrixd


tr2 = cc.ccGLMatrixd()
tr2.initFromParameters(0.47, (1.5, 2.9, 6.3), (5,3,-2))

For these transformation examples, we create an ellipsoid cloud:


sphere = cc.ccSphere(1.0)
cloud = sphere.samplePoints(False, 100000)
cloud.scale(1.0, 3.0, 9.0)

And we rotate it, with a transformation and the method applyRigidTransformation()


transform1 = cc.ccGLMatrix()
transform1.initFromParameters(0.25, (1.5, 2.9, 6.3), (0,0,0))
cloud.applyRigidTransformation(transform1)

You can get the rotation parameters of a transformation in the form of a rotation angle and a vector with getParameters1()


matrixParams1 = transform1.getParameters1()

if not math.isclose(matrixParams1.alpha_rad, 0.25, rel_tol=1.e-5):
    raise RuntimeError
if not isCoordEqual(matrixParams1.axis3D, (0.211393, 0.408694, 0.887852), 1.e-5):
    raise RuntimeError

You can also get the rotation parameters of a transformation in the form of the 3 angles phi, psi, theta with getParameters2()


matrixParams2 = transform1.getParameters2()

if not math.isclose(matrixParams2.phi_rad, 0.225260, rel_tol=1.e-5):
    raise RuntimeError
if not math.isclose(matrixParams2.psi_rad, 0.0639141, rel_tol=1.e-5):
    raise RuntimeError
if not math.isclose(matrixParams2.theta_rad, 0.0954226, rel_tol=1.e-5):
    raise RuntimeError
if not isCoordEqual(matrixParams2.t3D, (0., 0., 0.)):
    raise RuntimeError

The product of two transformations gives a transformation


x =transform1*transform1

mpx = x.getParameters1()
if not math.isclose(mpx.alpha_rad, 0.5, rel_tol=1.e-5):
    raise RuntimeError
if not isCoordEqual(mpx.axis3D, (0.211393, 0.408694, 0.887852), 1.e-5):
    raise RuntimeError

To get the inverse of a transformation, use inverse()


rotinv = transform1.inverse()

mpi = rotinv.getParameters1()
if not math.isclose(mpi.alpha_rad, 0.25, rel_tol=1.e-5):
    raise RuntimeError
if not isCoordEqual(mpi.axis3D, (-0.211393, -0.408694, -0.887852), 1.e-5):
    raise RuntimeError

Let’s play with the 3 rotations phi, psi, theta:


rotphi = cc.ccGLMatrix()
rotphi.initFromParameters(0.1, (0., 0., 1.), (0,0,0))

a1 = rotphi.getParameters2()
if not math.isclose(a1.phi_rad, 0.1, rel_tol=1.e-7):
    raise RuntimeError

rotpsi = cc.ccGLMatrix()
rotpsi.initFromParameters(0.2, (1., 0., 0.), (0,0,0))

a2 = rotpsi.getParameters2()
if not math.isclose(a2.psi_rad, 0.2, rel_tol=1.e-7):
    raise RuntimeError

rottheta = cc.ccGLMatrix()
rottheta.initFromParameters(0.3, (0., 1., 0.), (0,0,0))

a3 = rottheta.getParameters2()
if not math.isclose(a3.theta_rad, 0.3, rel_tol=1.e-7):
    raise RuntimeError

rotTaitBryan = rotphi*rottheta*rotpsi

a = rotTaitBryan.getParameters2()
if not isCoordEqual((a.phi_rad, a.psi_rad, a.theta_rad), (0.1, 0.2, 0.3), 1.e-7):
    raise RuntimeError

The transposed version of a transformation gives an inverse rotation:


mat1 = rotTaitBryan.transposed()
p1 = rotTaitBryan.getParameters1()
pt1 = mat1.getParameters1()
if not math.isclose(p1.alpha_rad, pt1.alpha_rad, rel_tol=1.e-7):
    raise RuntimeError
angle1 = p1.axis3D
angle2 = pt1.axis3D
for i in range(3):
    if not math.isclose(angle1[i], - angle2[i], rel_tol=1.e-7):
        raise RuntimeError

Another way to check:


mat = (rotTaitBryan.transposed()).inverse()

# --- access to transformation data from numpy

v= np.array(rotTaitBryan.data()) - np.array(mat.data()) # should be null

d2 = np.inner(v, v)
if d2>1e-12:
    raise RuntimeError

An interpolation between transformations can be computed with Interpolate():


rot1 = cc.ccGLMatrix()
rot1.initFromParameters(math.pi/3., (0., 0., 1.), (0,0,0))
rot2 = cc.ccGLMatrix()
rot2.initFromParameters(2*math.pi/3., (0., 0., 1.), (0,0,0))
rotation = cc.ccGLMatrix.Interpolate(0.5, rot1, rot2)

a = rotation.getParameters1()
if not isCoordEqual(a.axis3D, (0., 0., 1.)):
    raise RuntimeError
if not math.isclose(a.alpha_rad, math.pi/2., rel_tol=1.e-7):
    raise RuntimeError

There is also the method FromToRotation() to get a transformation matrix that rotates a vector to another:


rot1 = cc.ccGLMatrix.FromToRotation((1., 0., 0.),(0., 0., 1.))

a = rot1.getParameters1()
if not isCoordEqual(a.axis3D, (0., -1., 0.)):
    raise RuntimeError
if not math.isclose(a.alpha_rad, math.pi/2., rel_tol=1.e-7):
    raise RuntimeError

You can generate a ‘viewing’ matrix from a looking vector and a ‘up’ direction with FromViewDirAndUpDir()


r2 = math.sqrt(2.)
rot1 = cc.ccGLMatrix.FromViewDirAndUpDir((0., r2, r2), (0., -r2, r2))

a = rot1.getParameters1()
if not isCoordEqual(a.axis3D, (-1., 0., 0.)):
    raise RuntimeError
if not math.isclose(a.alpha_rad, 3*math.pi/4., rel_tol=1.e-7):
    raise RuntimeError

4.1.3. cloud copy, destruction

A cloud can be cloned with all its features except from the octree with cloneThis():


cloned = cloud.cloneThis()

To free memory, a cloud can be deleted with cloudComPy.deleteEntity() (WARNING be sure to have no more Python objects referencing the deleted object):


cc.deleteEntity(cloned) # delete the cloud copy
cloned = None

The above code snippets are from test026.py.

Some methods give a selection from a cloud as a cloudComPy.ReferenceCloud which is a light structure referencing the selected nodes in the original cloud. To convert this selection in a new cloud, use the method partialClone():


cloud = cc.loadPointCloud(getSampleCloud(5.0))

refCloud = cc.CloudSamplingTools.noiseFilter(cloud, 0.04, 1.0) # selection on cloud
origCloud = refCloud.getAssociatedCloud()            # the original cloud ~ cloud
(noiseCloud, res) = origCloud.partialClone(refCloud) # ccPointCloud from selection, status
noiseCloud.setName("noiseCloud")

if refCloud.__class__ != cc.ReferenceCloud:
    raise RuntimeError
if refCloud.size() < 7470 or refCloud.size() > 7570:
    raise RuntimeError
if res != 0:
    raise RuntimeError

The above code snippet is from test019.py.

4.1.4. cloud normals

Once computed (see computeNormals()), normals can be exported to scalar fields with exportNormalToSF():


cloud.exportNormalToSF(True, True, True)

sf=cloud.getScalarField(2)
if sf.getName() != 'Nz':
    raise RuntimeError

Normals can also be converted to color with convertNormalToRGB():

The colors obtained are an HSV color field, with H = dip direction, S = dip and V = 1.


cloud.convertNormalToRGB()

if not cloud.hasColors():
    raise RuntimeError

To compute strike and dip fields from normals, use convertNormalToDipDirSFs() (see definition for Strike and dip):


cloud.convertNormalToDipDirSFs()

dicsf = cloud.getScalarFieldDic()
sfdip = cloud.getScalarField(dicsf['Dip (degrees)'])
sfdipd = cloud.getScalarField(dicsf['Dip direction (degrees)'])

Normals can be inverted with cloudComPy.invertNormals(), reoriented with a Fast Marching method or a Minimum Spanning Tree method, with orientNormalsWithFM() and orientNormalsWithMST().


if not cloud.orientNormalsWithFM():
    raise RuntimeError

if not cloud.orientNormalsWithMST():
    raise RuntimeError

if not cc.invertNormals([cloud]):
    raise RuntimeError

The above code snippets are from test014.py.

A cloud can be shifted along the normals with shiftPointsAlongNormals(). For instance:


cloud2.shiftPointsAlongNormals(0.5)

The above code snippet is from test055.py.

4.1.5. cloud colors

The cloudComPy.QColor class wraps the Qt5 QColor class and provides provides colors based on RGB, HSV or CMYK values. See test028.py for an example of use of cloudComPy.QColor methods.

Several methods of cloudComPy.ccPointCloud allow to colorize the cloud:


cloud.colorize(0.2, 0.3, 0.4, 1.0) 
if not cloud.hasColors():
    raise RuntimeError

col = cc.QColor.fromRgb(32,48,64) 
if not cloud.setColor(col):
    raise RuntimeError

c1=cc.QColor.fromRgbF(1.0, 0., 0.)
c2=cc.QColor.fromRgbF(0., 0., 1.0)
if not cloud.setColorGradient(0, c1, c2):
    raise RuntimeError

if not cloud.setColorGradientDefault(2):
    raise RuntimeError

if not cloud.setColorGradientBanded(1, 3.0):
    raise RuntimeError

if not cloud.changeColorLevels(60, 180, 20, 220, 1, 1, 1):
    raise RuntimeError

if not cloud.convertRGBToGreyScale():
    raise RuntimeError

ScalarFields can be used to define or modify colors, colors can be used to define scalarFields:


cloud.exportCoordToSF(True, True, True)
cloud.setCurrentDisplayedScalarField(0)

if not cloud.convertCurrentScalarFieldToColors(mixWithExistingColor=False):
    raise RuntimeError

if not cloud.enhanceRGBWithIntensitySF(sfIdx=1):
    raise RuntimeError

n1 = cloud.getNumberOfScalarFields()

cloud.sfFromColor(True, True, True, False, False) # export R, G, B, not alpha, not composite

n2 = cloud.getNumberOfScalarFields()
if (n2-n1) != 3:
    raise RuntimeError

If you don’t need any more colors, it is possible to free some memory with unallocateColors() (WARNING be sure to have no more Python objects referencing the deleted object):


cloud.unallocateColors()
if cloud.hasColors():
    raise RuntimeError

With two clouds “sharing a same region”, it is possible to define the color of one cloud by interpolation from the other cloud, with interpolateColorsFrom()


cloud1 = cc.loadPointCloud(getSampleCloud(1.0))

if not cloud1.interpolateColorsFrom(cloud):
    raise RuntimeError

The above code snippets are from test029.py.

4.1.6. scalar fields

when a scalar field is modified with Numpy (see Read, modify or create a scalar field with Numpy), you must reinitialise the min and max value of the scalar field with computeMinAndMax().


dic = cloud1.getScalarFieldDic()
sf1 = cloud1.getScalarField(dic['Coord. Z'])
max1 = sf1.getMax()
asf1 = sf1.toNpArray()   # access to Numpy array, without copy
asf1[0] = 2*max1         # modification in place
sf1.computeMinAndMax()

General information and statistics are available with:


sfname = sf1.getName()
sfmin = sf1.getMin()
sfmax = sf1.getMax()
mean, var = sf1.computeMeanAndVariance()
val = sf1.getValue(23)

The above code snippets are from test002.py.

Some scalar fields may be shifted to prevent a loss of accuracy. This is the case, for instance, for GPS time in some lidar files in las format.

To get the global shift, use:


cloud=cc.loadPointCloud(os.path.join(dataExtDir,"PTS_LAMB93_IGN69_extract.las"))
dic = cloud.getScalarFieldDic()
sf = cloud.getScalarField(dic['Gps Time'])
timeShift = sf.getGlobalShift()

The above code snippets is from test020.py.

To change the scalar field name, set a value on a point, fill the scalar field with a uniform value, use:


sf1.setName("aNewName")
sf1.fill(0.5)
sf1.setValue(3, 2.25) # index, value

The above code snippet is from test002.py.

Scalar fields can be built from normals or colors, and can be used to define colors: see cloud normals and cloud colors.

To get a scalar field gradient, use cloudComPy.ccPointCloud.computeScalarFieldGradient() with an appropriate radius:


radius = cc.GetPointCloudRadius([cloud], 12) # number of nodes wanted within the radius
if not cloud.computeScalarFieldGradient(1, radius, True):
    raise RuntimeError

The previous code is extract from test003.py.

Some methods of cloudComPy.ccPointCloud, directly wrapped from CloudCompare GUI, are used to select a current scalar field. This is required by some other methods, working with the current selected scalar field:

For instance, the filterPointsByScalarValue() method relies on setCurrentDisplayedScalarField() to work. The points of the input cloud are filtered by keeping (or excluding) the points for which the value of the scalar field belongs to a (min,max) interval.


dupSFindex = cloud.addScalarField("DuplicateFlags")
cloud.setCurrentScalarField(dupSFindex)
ret = cc.GeometricalAnalysisTools.FlagDuplicatePoints(cloud) # identify duplicated points
if ret != cc.ErrorCode.NoError:
    raise RuntimeError
noDuplCloud = cloud.filterPointsByScalarValue(0., 0., outside=False) # remove duplicated pts

The previous code is extract from test019.py.

If you need to free some memory and do not need any more some scalar fields, you can use the following methods of cloudComPy.ccPointCloud (WARNING be sure to have no more Python objects referencing the deleted object):

4.1.7. sensors

Acquisition data files sometimes contain information on the sensors used to capture data.

A short introduction to sensors in CloudCompare is given in CloudCompare wiki - Sensors.

CloudComPy gives access to basic information on sensors:

  • type of sensor,

  • position and orientation,

  • graphic scale,

  • uncertainty.

For the test, we need a data file to download here.


# example data available here: http://sourceforge.net/projects/e57-3d-imgfmt/files/E57Example-data/manitouNoInvalidPoints.e57/download
if not os.path.isfile(os.path.join(dataExtDir,"manitouNoInvalidPoints.e57")):
    if not os.path.exists(dataExtDir):
        os.makedirs(dataExtDir)
    url = "https://www.simulation.openfields.fr/phocadownload/manitouNoInvalidPoints.e57"
    r = requests.get(url)
    with open(os.path.join(dataExtDir,"manitouNoInvalidPoints.e57"), 'wb') as f:
        f.write(r.content)
        
entities = cc.importFile(os.path.join(dataExtDir,"manitouNoInvalidPoints.e57"))

There are several point clouds in the example, with one sensor per cloud. The getSensors() method gives a list of sensors associated to the cloud. Generic sensor cloudComPy.ccSensor provides the following methods:

The Ground Based Lidar (GBL) sensor cloudComPy.ccGBLSensor provides also:


for entity in entities[1]:
    entity.getName()
    sensors = entity.getSensors()
    if len(sensors) < 1:
        raise RuntimeError
    sensor= sensors[0]
    if sensor.getType() != cc.CC_SENSOR_TYPE.GROUND_BASED_LIDAR:
        raise RuntimeError
    if sensor.getClassID() != cc.CC_TYPES.GBL_SENSOR:
        raise RuntimeError
    grscale = sensor.getGraphicScale()
    print("graphic scale:", grscale)
    tr=sensor.getRigidTransformation()
    ret1=tr.getParameters1()
    alpha = ret1.alpha_rad*180/math.pi
    print("alpha:", alpha)
    axis = ret1.axis3D
    print("axis", axis)
    trans = ret1.t3D
    print("translation:", trans)
    uncertainty = sensor.getUncertainty()
    print("uncertainty", uncertainty)

Sensors are also used to compute scattering angles with ComputeScatteringAngles().


cc.computeNormals(entities[1])
for entity in entities[1]:
    sensors = entity.getSensors()
    if len(sensors) < 1:
        raise RuntimeError
    sensor= sensors[0]
    sensor.ComputeScatteringAngles()

The above code snippets are from test041.py.

4.2. meshes

4.2.1. create a mesh from a cloud

When a cloud represents a kind of 2.5D elevated surface, a mesh can be built using the nodes of a the cloud with cloudComPy.ccMesh.triangulate().


cloud1 = cc.loadPointCloud(getSampleCloud2(3.0, 0, 0.1))
cloud1.setName("cloud1")

mesh1 = cc.ccMesh.triangulate(cloud1, cc.TRIANGULATION_TYPES.DELAUNAY_2D_AXIS_ALIGNED, dim=2)
mesh1.setName("mesh1")

The above code snippet is from test011.py.

4.2.2. meshes introspection

In CloudCompare and CloudComPy 3D meshes are triangular meshes (ccMesh), built on a set of vertices, which constitute an associated cloud. To get this cloud, use getAssociatedCloud() which gives a ccPointCloud.


meshply = cc.loadMesh(os.path.join(dataDir, "mesh.ply"))
cloud1 = meshply.getAssociatedCloud()

The above code snippet is from test020.py.

A mesh have a name, distinct of the associated cloud name, a size which is the number of triangles.

A mesh can have normals, computed either on vertices or on triangles:


cloud1 = cc.loadPointCloud(getSampleCloud2(3.0,0, 0.1))
cloud1.setName("cloud1")
plane = cc.ccPlane.Fit(cloud1)
mesh1 = cc.ccMesh.triangulate(cloud1, cc.TRIANGULATION_TYPES.DELAUNAY_2D_AXIS_ALIGNED)
mesh1.setName("mesh1")
cc.computeNormals([mesh1], computePerVertexNormals=False)

The above code snippet is from test014.py.

Normals on triangles can be deleted with clearTriNormals().

If you need to iterate through the triangles and their vertices, use getTriangleVertIndexes():


# --- access to triangle nodes, per triangle indice
cloud = mesh1.getAssociatedCloud()
indexes = mesh1.getTriangleVertIndexes(453)
p0 = cloud.getPoint(indexes[0])
p1 = cloud.getPoint(indexes[1])
p2 = cloud.getPoint(indexes[2])

The above code snippet is from test011.py.

If you want to save meshes to reopen them with the CloudCompare GUI, with a predefined state of what is displayed (colors, normals, scalar fields), use the .bin format and define the state with the following functions:

It is possible to compute the surface of a mesh (sum of the triangles surfaces), or the volume inside, with computeMeshArea() or computeMeshVolume().

The mesh is not necessarily closed: the object returned by computeMeshVolume() is a tuple, containing the volume, a boolean indicating if the mesh is not closed, and structure giving the number of the different types of edges: (total, shared, not shared, and shared by more of 2 triangles).

An example with a closed mesh:


sphere = cc.ccSphere(radius=2, precision=128)
area = sphere.computeMeshArea()
res = sphere.computeMeshVolume()
if res[1]:                # there are edges not shared or shared by more than 2 triangles
    raise RuntimeError
vol = res[0]
sphereArea = 16*math.pi            # geometric values 
sphereVolume = (4/3.)*math.pi*8
if not math.isclose(area, sphereArea, rel_tol=1.e-3):
    raise RuntimeError
if not math.isclose(vol, sphereVolume, rel_tol=1.e-3):
    raise RuntimeError

An example with an open mesh:


sphere = cc.ccSphere(radius=2, precision=128)
area = sphere.computeMeshArea()
res = sphere.computeMeshVolume()
if res[1]:                # there are edges not shared or shared by more than 2 triangles
    raise RuntimeError
vol = res[0]
sphereArea = 16*math.pi            # geometric values 
sphereVolume = (4/3.)*math.pi*8
if not math.isclose(area, sphereArea, rel_tol=1.e-3):
    raise RuntimeError
if not math.isclose(vol, sphereVolume, rel_tol=1.e-3):
    raise RuntimeError

The above code snippets are from test055.py.

4.2.3. meshes modifications

A mesh can be refined using subdivide() to force all triangles to have an area smaller than a given maximum. The result is a new mesh.


mesh3 = mesh2.subdivide(0.001)
mesh3.setName("mesh3")

The mesh can be “smoothed” with laplacianSmooth() by moving the vertices slightly over several iterations:


mesh3.laplacianSmooth(nbIteration=20, factor=0.2)

The above code snippets are from test011.py.

The mesh triangles can be reversed (invert the normals) with flipTriangles().


dish.flipTriangles()

The above code snippet is from test055.py.

4.2.4. Generate a cloud from a mesh

A mesh can be used to create a cloud with samplePoints(), with a target either of number of points or of density of points:


cloud2=mesh1.samplePoints(densityBased=True, samplingParameter=50, withNormals=True)
cloud2.setName("cloud2")

The above code snippet is from test015.py.

4.2.5. mesh copy, destruction

A mesh can be cloned with all its features, except from the octree, with :py:meth`~.cloudComPy.ccMesh.cloneMesh`:


mesh2 = mesh1.cloneMesh()
if mesh2.getName() != "mesh1.clone":
    raise RuntimeError

To free memory, a mesh can be deleted with cloudComPy.deleteEntity() (WARNING be sure to have no more Python objects referencing the deleted object):


cc.deleteEntity(mesh3)
mesh3=None

The above code snippets are from test011.py.

4.3. primitives

Primitives regroups all the geometric generators provided by CloudCompare:

All the primitives derives from ccGenericPrimitive which is derived from ccMesh. Thus, all the mesh methods apply here.

All primitive constructors use as argument: geometrical features, an optional drawing precision, an optional transformation (always possible afterwards).

The primitives have some more methods than the meshes:


cylinder = cc.ccCylinder(0.5, 6.0)
if cylinder.getTypeName() != 'Cylinder':
    raise RuntimeError
if cylinder.size() != 96:
    raise RuntimeError

The above code snippet is from test011.py.

4.3.1. Box


tr1 = cc.ccGLMatrix()
tr1.initFromParameters(0.1, 0.2, 0.3, (8.0, 0.0, 0.0))
box = cc.ccBox((1., 2., 3.), tr1, "aBox")
if box.getName() != 'aBox':
    raise RuntimeError
if box.size() != 12:
    raise RuntimeError

4.3.2. Cone

Cone axis corresponds to the ‘Z’ dimension by default.


tr2 = cc.ccGLMatrix()
tr2.initFromParameters(0.5, (0., 1., 0.), (5.0, 6.0, 3.0))
cone = cc.ccCone(3., 1., 2., 0., 0., tr2, "aCone", 12)
if cone.getName() != 'aCone':
    raise RuntimeError
if cone.size() != 48:
    raise RuntimeError

A ccCone offers several methods to retreive its geometric features:

4.3.3. Cylinder


tr3 = cc.ccGLMatrix()
tr3.initFromParameters(0., (0., 0., 0.), (3.0, 0.0, 4.0))
cylinder = cc.ccCylinder(0.5, 3.0, tr3, 'aCylinder', 48)
if cylinder.getName() != 'aCylinder':
    raise RuntimeError
if cylinder.size() != 192:
    raise RuntimeError

A ccCylinder derives from a ccCone. It is internally represented by a cone with the same top and bottom radius.

4.3.4. Plane


tr4 = cc.ccGLMatrix()
tr4.initFromParameters(0.5, (0., 1., 0.), (-3.0, 0.0, 4.0))
plane = cc.ccPlane(3.0, 4.0, tr4, "aPlane")
if plane.getName() != 'aPlane':
    raise RuntimeError
if plane.size() != 2:
    raise RuntimeError

By default, plane normal corresponds to ‘Z’ dimension.

A ccPlane offers several methods to retreive its geometric features:

The method Fit() allows to adjust a plane primitive on a cloud:


cloud1 = cc.loadPointCloud(getSampleCloud2(3.0,0, 0.1))
cloud1.setName("cloud1")
tr1 = cc.ccGLMatrix()
tr1.initFromParameters(0.2*math.pi, (1., 1., 1.), (0.0, 0.0, 10.0))
cloud1.applyRigidTransformation(tr1)

plane = cc.ccPlane.Fit(cloud1)
equation = plane.getEquation()
tr2 = plane.getTransformation()

The above code snippet is from test012.py.

4.3.5. Quadric


tr7 = cc.ccGLMatrix()
tr7.initFromParameters(0.5*math.pi, (1., 0., 0.), (-5.0, -2.0, -2.0))
quadric = cc.ccQuadric((-1., -1.), (1., 1.), (1., 2., 1., 1., 2., 2.), dims=(0,1,2), transMat=tr7, name="aQuadric", precision=60)
if quadric.getName() != 'aQuadric':
    raise RuntimeError
if quadric.size() != 6962:
    raise RuntimeError

Quadric orthogonal dimension is ‘Z’ by default.

4.3.6. Sphere


tr5 = cc.ccGLMatrix()
tr5.initFromParameters(0.0, (0., 0., 0.), (-7.0, 5.0, 1.0))
sphere = cc.ccSphere(1.5, tr5, "aSphere", 72)
if sphere.getName() != 'aSphere':
    raise RuntimeError
if sphere.size() != 10224:
    raise RuntimeError

A ccSphere offers the method getRadius().

4.3.7. Torus


tr6 = cc.ccGLMatrix()
tr6.initFromParameters(0.2, (1., 2., 0.), (-0.0, -5.0, -9.0))
torus = cc.ccTorus(5.0, 7.0, math.pi, True, 3.0, tr6, "aTorus", 60)
if torus.getName() != 'aTorus':
    raise RuntimeError
if torus.size() != 964:
    raise RuntimeError

4.3.8. Dish


tr8 = cc.ccGLMatrix()
tr8.initFromParameters(0.5*math.pi, (0., 1., 0.), (-0.0, -0.0, -4.0))
dish = cc.ccDish(2.0, 1.0, 3.0, tr8, "aDish", 72)
if dish.getName() != 'aDish':
    raise RuntimeError
if dish.size() != 2520:
    raise RuntimeError

Unless Otherwise noted, the above code snippets are from test009.py.

4.4. polylines

A ccPolyline is considered as a cloud of points (in a specific order) with a open/closed state information. Polylines are often used as tools to cut clouds (see Cut a cloud or a mesh with a polyline).

A polyline can be created from a cloud (the points have to be in the good order):


cloud1 = cc.ccPointCloud("boundingBox1")
cloud1.coordsFromNPArray_copy(coords)
cloud1.applyRigidTransformation(transform1)

poly1 = cc.ccPolyline(cloud1)
poly1.setName("poly1")
poly1.addChild(cloud1)
poly1.addPointIndex(0, cloud1.size())
poly1.setClosed(True)

Another way (here, we simulate the cloud creation):


cloud2 = cc.ccPointCloud("boundingBox2")
poly2 = cc.ccPolyline(cloud2)
poly2.setName("poly2")
poly2.addChild(child=cloud2, dependencyFlags = cc.DEPENDENCY_FLAGS.DP_NONE, insertIndex=-1)
cloud2.reserve(cloud1.size())   # fill the cloud with ordered points
for i in range(cloud1.size()):
    cloud2.addPoint(cloud1.getPoint(i))
cloud2.shrinkToFit()
poly2.addPointIndex(0, cloud2.size())
poly2.setClosed(True)

It is also possible tor create a polyline directly from a list of 2d or 3d coordinates. With 2d coordinates, the Z coordinate is set to 0.


poly2d = cc.ccPolyline([(0., 0.), (0., 1.), (1., 1.), (1., 0.)], True)
poly2d.setName("poly2d")
poly3d = cc.ccPolyline([(0., 0., 0.), (0., 1., 0.1), (1., 1., 0.2), (1., 0., 0.3)], False)
poly3d.setName("poly3d")
# get the associated cloud and add it as child, in order to save the poly...
cloud2d = poly2d.getAssociatedCloud()
poly2d.addChild(cloud2d)
cloud3d = poly3d.getAssociatedCloud()
poly3d.addChild(cloud3d)

The above code snippets are from test026.py.

A polyline can be imported from several formats (see Known formats).

The .poly format is a very simple ASCII format, a line per point, x y z. for example:

-4.0  4.0  0.0
 0.0  3.5  0.0
 4.0  4.0  0.0
 3.5  0.0  0.0
 4.0 -4.0  0.0
 0.0 -3.5  0.0
-4.0 -4.0  0.0
-3.5  0.0  0.0

The class ccPolyline offers several information methods:

The status of the polyline can be modified:

When creating the polyline from a cloud (see examples above) you will need:

The method smoothChaikin() is used to smooth the polyline by adding nodes over several iterations:


p2 = poly.smoothChaikin(0.2, 4)

The above code snippet is from test007.py.

4.5. 2D polygons (facets)

The class ccFacet is a tool used to fit a 2D polygon to a point cloud. It is very similar to the Fit() method but the extents of the fitted plane follows the contour of the cloud (on its 2D projection on the plane).

The 2D polygon is created using the static method Create(). If you want a convex contour, leave the argument maxEdgeLength to its default, 0. A positive value will gives a concave contour respecting maxEdgeLength.

It is possible to provide a plane equation as optional argument (see cloudComPy.ccPlane.getEquation()).

The class ccFacet provides several useful methods:

getCenter()

the facet plane center

getNormal()

the facet plane normal

getPlaneEquation()

the facet plane equation: [a, b, c, d] as ax+by+cz=d

getContour()

the contour polyline

getContourVertices()

the point cloud with the contour vertices

getSurface()

the area of the polygon


polygon3 = cc.ccFacet.Create(cloudCropZ, maxEdgeLength=0.5)
center3 = polygon3.getCenter()
normal3 = polygon3.getNormal()
eq3 = polygon3.getPlaneEquation()
contour3 = polygon3.getContour()
vert3 = polygon3.getContourVertices()
surface3 = polygon3.getSurface()
if not math.isclose(surface3, 56, rel_tol=2e-02):
    raise RuntimeError

With the plane equation:


polygon4 = cc.ccFacet.Create(cloud=cloudCropZ, maxEdgeLength=0.05, planeEquation=(0., 0., 1., 0.))
center4 = polygon4.getCenter()
normal4 = polygon4.getNormal()
eq4 = polygon4.getPlaneEquation()
contour4 = polygon4.getContour()
vert4 = polygon4.getContourVertices()
surface4 = polygon4.getSurface()
if not math.isclose(surface4, 56, rel_tol=1e-02):
    raise RuntimeError

The above code snippets are from test021.py.

4.6. octrees

The Octree structure of CloudCompare is very efficient for nearest neighbour extraction and used throughout CloudCompare. The CloudCompare wiki provides a good introduction to the CloudCompare octree.

Most of the time, you don’t need to manipulate explicitely the octree, but some methods have been wrapped in CloudComPy, to allow for some very specific processing.

If you need to search points in particular neighbourhood, ccOctree provides useful methods illustrated in the following example.

Firstly, compute the octree. Define your neighbourhood either with a radius or the number of nodes you want in this neighbourhood.


# --- get a sample cloud, build the octree

cloud = cc.loadPointCloud(getSampleCloud(5.0))
octree = cloud.computeOctree(progressCb=None, autoAddChild=True)
if (octree.getNumberOfProjectedPoints() != 1000000):
    raise RuntimeError

# --- search points in neighbourhood, within a given radius r

r=0.05
level = octree.findBestLevelForAGivenNeighbourhoodSizeExtraction(r)
level12 = octree.findBestLevelForAGivenPopulationPerCell(12)

Search neighbours in a sphere:


# --- search in a sphere

neighbours = octree.getPointsInSphericalNeighbourhood((-1.5, 2.0, -0.026529), r, level)
if len(neighbours) != 37:
    raise RuntimeError

Get information on the points found:


# --- get information on some points

pt = neighbours[0].point
id = neighbours[0].pointIndex
sd = neighbours[0].squareDistd
if (math.sqrt(sd) > r):
    raise RuntimeError

It is also possible to search in a cylinder:


# --- search in a cylinder

params = cc.CylindricalNeighbourhood()
params.center = (-1.5, 2.0, 0.0)
params.level = level
params.radius = r
params.maxHalfLength = 1

neighboursCyl = octree.getPointsInCylindricalNeighbourhood(params)
if len(neighboursCyl) != 81:
    raise RuntimeError

or in a box:


# --- search in a box

params = cc.BoxNeighbourhood()
params.center=(-1.500000, 2.000000, -0.026529)
params.level = level
params.dimensions = (0.05, 0.04, 0.03)

neighboursBox = octree.getPointsInBoxNeighbourhood(params)
if len(neighboursBox) != 9:
    raise RuntimeError

The above code snippets are from test016.py where you will find examples of how to use the other methods in the octree.