Module DBH
See also
opals::IDBH

Aim of module

Estimate the Diameter at Breast Height (DBH) based on 3D point cloud data (ODM).

General description

In forestry the diameter at breast height is an important measure for various estimation models and applications. Although, the model was primarily designed for stem estimation it can also be applied for modelling artificial objects like pipes or other cylindrical and conical objects. Based on given approximations (3D point, axis and radius) the module applies robust least-square fitting of cylinders or cones to neighboring points of the approximate 3d point. Hence, for each stem a separate entry within the approximation file is required. Beside the standard single fit mode, the module can also trace the stem in and opposite to the axis direction respectively (parameter trace). In trace mode the parameters patchLength, overlap and the results of the previous fit are used to compute new approximations for the next fitting step. The actual output is a shape file with 3d circles at representative position with additional statistical information of each fit. Furthermore, from python and C++ all results can be directly accessed as a vector of stem objects by requiring the output parameter outResults (see here for further details).

The provided approximations are not only used as initial values for the linearized least-square models but also for sub-selecting relevant points for the actual fit. As shown in Figure 1, Module DBH uses an extended cylinder for the data selection, where \( P_1 \) (first point of the stem approximation) is a the center of the selection cylinder, the cylinder axis \( \vec{a} \) is aligned with \( P_1 \) and \( P_2 \) (second point of the stem approximation) and the cylinder height is equal to patchLength. Whereas vertical axis approximations are appropriate for trees in general, in rare cases of titled or horizontal (e.g. dead wood) stems with many outlier points in the surrounding (e.g. undergrowth, shrubs or nearby trees) it might be necessary to have more accurate non-vertical axis approximations since this will select more relevant stem points.

Fig. 1: Data selection strategy

The radius of the data selection cylinder can be controlled by the searchRadius parameter. The definition of an appropriate formula is essential to get correct least-square fits, since too small values may excluded correct stem points from the fit. On the other hand, oversized searchRadius values may introduce too many false points precluding the implemented outlier detection to work (correctly). So depending on the input data and the accuracy of the approximations, it can be necessary to adapt the searchRadius formula. Module DBH provides the approximate point (x, y and z) radius (_r) and axis (_axisX, _axisY and _axisZ), the stem and trace id (_StemId and _TraceId) to the formula allow defining fine adoptions to the searchRadius values. In case the formula returns invalid (see Generic filter for further details), no data selection and consequently, no fit will be performed (in trace mode the iteration in the current direction is stopped).

Trace mode

As mentioned before, ModuleDBH is capable of tracing stems in and opposite to the axis direction. This feature can be used to either extract a more detailed representation of the corresponding tree or to get the stem diameter at different height(s) than the initial approximation point height. Especially in case of dense undergrowth extracting initial stem position higher than breast height can help to increase robustness. Then, using down tracing the stems allows extracting the diameter at the sought-after height. To enable trace mode in axis direction set the first value trace to true. The second value controls the opposite axis direction tracing. Hence, it is possible to activate the two tracing directions independently. In case the axis direction is always oriented upwards (P2 is always above P1 as in the approx file below) the first trace flag matches the growing direction of the trees.

As shown on the left in Figure 2 the adjusted results of the previous patch are used as approximations for the next patch. The spacing between consecutive patches are controlled by the overlap parameter but also depends on the patchLength. An overlap of 0.5 results in 50% overlapping patches, whereas an overlap of 0 leads to touching patches. The length of the shift vector between consecutive patches is computed by

\( len = (1-overlap) \cdot patchLength \)

Fig. 2: Trace mode with different overlap configurations

The module stops tracing in one direction if two consecutive fits are classified as incorrect. Either the least-square adjustment itself did fail (e.g. no convergence was achieved, not enough input data were selected, etc.) or the continuity check discard the results. This continuity check use several heuristics to decide whether the current results fits to the previous result. The result is discard if,

  • the axis changes more than 10 degrees,
  • the radius changes more than 50%,
  • the relative overlapping area of the two circles (projected on one plane) is less than 60% or,
  • the distance between two consecutive adjusted points is below a certain threshold (50% of \( len \))

Details on which criteria stop the iterations will be logged, if the file or screen log level is set to debug (see Logging / Error Handling for further details).

Definition of approx file

User can either prepare the approx file in the default ASCII format or a corresponding OPALS format definition file can be specified. In the default format two points and the approximate radius must be defined for each stem

#axis pt1.x pt1.y pt1.z pt2.x pt2.y pt2.z radius
-7.581 7.944 101.467 -7.581 7.944 102.467 0.442
-9.732 8.928 101.510 -9.732 8.928 102.510 0.513
-9.797 4.386 101.446 -9.797 4.386 102.446 0.485
[...]

where the first point ( \( pt1 \)) is used as center of the selection region and the approximate axis points from the first to the second point ( \( a = pt2 - pt1 \)). Since lines with a leading '#' character are ignored, it is possible to place comment or header lines within the file but also, deactivating single stem approximations by inserting an '#' character at the beginning of the corresponding line.

The usage OPALS format definition adds some extra feature that maybe valuable in certain situations:

  • pass through extra data columns to final results: all columns that are not used for defining stem approximations (for further details see table below) will be internally stored and transfered unchanged to the output results
  • use individual filter string for each stem: although the module supports a global filter, individual stem-based filter may help to overcome tricky situation with closely located stems or dense undergrowth by e.g. using segment ids from a previously performed segmentation. In case a global and a stem-based filter are defined, both filters are combined by a logical AND operator.
  • provide point and axis information rather than two points: In situation where it is easier or more convenient to provide an approximate point and axis information, a corresponding format definition can be used (see table below).

The columns of an approximation files that are interpreted in the following way:

Interpretation of columns by ModuleDBH
Attribute Interpretation Remark
x,y,z (=default coordinates) first point
_x1, _y1, _z1 first point alternative coordinate definition
_x2, _y2, _z2 second point
_axisX, _axisY, _axisZ axis vector alternative to second point
_r radius
_radius radius alternative to attribute _r
_filter individual filter optional data filter for each stem
[...] pass through columns optional

For ModuleDHB a valid OPALS format definition files must contain an approximate radius value and

  • a first and second point or
  • a (first) point and an axis vector

for each stem entry. Individual filters and pass through columns are optional. Although the approximation will be stored in an ASCII file in general (due to the low number of stem approximation the poor IO performance of ASCII files is irrelevant), the OPALS format definition concept is flexible enough to support binary or even LAS/LAZ as long as the necessary data columns are provided. How to use pass through columns and Individual filters is shown in exampe 2.

Results as shape file and output parameter

For the cylinder and cone fit the module uses the functional models described by Lukács et. al. (1997). Using robust least-squares strategies the method can cope with gross errors and therefore classifies points into in- and outliers. To avoid extrapolation, the final resulting point ( \( P_{adj} \)) is computed by orthogonal projection of the center of gravity (CoG) of all inliers onto the adjusted axis. In other words, intersecting the adjusted cylinder/cone axis with plane e (is orthogonal to the axis and contains the inlier CoG) gives the representative point \( P_{adj} \).

Fig. 3: cylinder (left) and cone (right) fit models

The results of all successful fits are written to a shape and ASCII file (see parameter outFile). The content of both files is basically the same with one exception. As geometry object the shape file uses a visual representation of the fitted model, namely the 3d circles as shown in Figure 3, rather than just a point representation of \( P_{adj} \). Nevertheless, the coordinates of \( P_{adj} \) are attached as attributes to the circle object. When fitting a cone also the convergence angle (=half of opening angle of cone) is output. In case the cone radius is decreasing in axis direction (i.e. the cone apex is above plane e based on the axis direction, as in Figure 3), the value is positive. If the radius is increasing when advance in axis direction, the convergence angle will be set negative.

List of columns in the output files
Attribute Interpretation Remark
Id output id unique identifier needed to fulfill the shape file specification
StemId input stem id matches the approx input stem id
TraceId trace id 0 for the initial fit, increasing while forward tracing and decreasing (i.e. negative numbers) while backward tracing
x, y, z coordinates of \( P_{adj} \)
r adjusted radius
ax, ay, az adjusted axis vector
convAngle convergence angle set for cone fits only: it is the angle between the cone generatrix and the cone axis. The value is set positive if the cone radius is decreasing when advancing in axis direction and negative if it is increasing. This columns allow differentiating between cylinder and cone fits.
offsetX, offsetY, offsetZ offset vector from the approximate point to the adjusted axis
dr correction of approximate radius value
RadialDev radial deviation based on inlier points
Redundancy redundancy of adjustment
nObs number of adjustment observations (i.e. points)
nUsed number of used inlier observations (i.e active points)
[...] pass through columns optional

Beside the output files it is possible to directly access the result from Python and C++ by querying the output parameter outResults. The returned object is a vector of opals::StemModel objects. Each stem object in turn stores at least one successful model fit (or multiple fits in trace mode). Each successful fit is represented as opals::StemEntry. Therefore, it is possible to further process or analyse the full results without parsing any result file.

Parameter description

-inFileinput ODM file
Type: opals::Path
Remarks: mandatory
ODM file containing 3D point cloud of forest area.
-approxFilefile containing tree approximations
Type: opals::Path
Remarks: mandatory
White space separated ASCII file containing two 3d axis points and cylinder radius (x1,y1,z1,x2,y2,z2,r).
-aFormatfile format for tree approximations
Type: opals::String
Remarks: default=auto
More details to come
-filtermodify the input using a (tree of) filter(s)
Type: opals::String
Remarks: optional
A filter string in EBNF syntax as described in section 'Filter syntax' can be passed to restrict the set of input points used within the ICP algorithm (e.g. only last echoes or only smooth areas). Note that this parameter does not affect the final ODM files in which (regardless of the filter settings) always all points are included.
-outFileoutput (shape) file
Type: opals::Path
Remarks: estimable
Name of Shape file to which the successfully estimated circles are written to. The file contains the outline of the 3d circles + attributes.
Estimation rule: If not specified, the file is composed of: current working directory + name (body) of the input file + '.shp'.
-weightFuncpoint weight formula
Type: opals::String
Remarks: optional
An algebraic formula in generic filter syntax can be specified to calculate an individual weight for each point. Both, the points' attributes and the distance w.r.t the grid location can be used to determine the point weight. Inverse Distance Weighting (IDW) can be realized by using the built-in "Dist2D(n[0],n[1])" or "SqrDist2D(n[0],n[1])" function returning the linear or squared distance between the current grid location n[0] and the considered neighbour point n[1]. The following predefined formulae can be used:
IDW1.... inverse linear distance (i.e., 1/Dist2D(n[0],n[1])+epsilon)
IDW2.... inverse squared distance (i.e., 1/SqrDist2Dn[0],n[1])+epsilon)
IDW:i... inverse distance to the power of i ( i.e., 1/pow(Dist2D(n[0],n[1])+epsilon),i)
where "epsilon" in the above formula is replaced with an appropriate infinitesimal value during run-time in order to avoid divisions-by-zero. Please note, that a valid grid height can only be obtained if the weight function can successfully be evaluated for all k nearest neighbours of the respective grid location.
-geometryModelstem fitting geometry
Type: opals::DBHModel
Remarks: default=cylinder
Possible values:  
  auto_ ...... automatic model selection
  cylinder ... cylinder model
  cone ....... cone model
Specifies which geometry model will be fitted to the stem points.
-patchLengthdata selection patch length
Type: double
Remarks: default=1
Defines the point selection length based on the approximate axis point and axis direction. The approximate axis point and axis direction can be interpreted as reference plane. Points that have a larger orthogonal distance to the reference plane as patchLength/2 are not used for the model fit.
-overlappatch overlap in trace mode
Type: double
Remarks: default=0.5
In trace mode this factor defines how much consecutive patches should overlap. E.g. a value of 0.8 results in a patch overlap of 80%, 0 means that patches are just touching, and -1 results to a 'patch length' gap between consecutive patches.
-tracetracing stems up and/or downwards
Type: opals::Vector<bool>
Remarks: default=0 0
Based on the approximate axis values, this option enables (=true) or disabled (=false) tracing the stem in axis direction (first value) and in the inverse axis direction (second value). If only one value is provided the inverse axis direction will never be traced.
-searchRadiuscalculates radius for data selection
Type: opals::Calculator<DM::ICalculator::ReadAccess::attributes, DM::ICalculator::WriteAccess::none>
Remarks: default=_r < 0.08 ? 0.1 : _r*1.25
The module provide the approximate point (x,y,z), axis (_axisX,_axisY,_axisZ) and radius (_r), as well as, the stem and trace id (_StemId, _TraceId) to the calculator allows defining a fine adapted search radius. If the formula returns invalid no fit will be performed.
Calculator<attributes, none>: has read access to: attributes; has write access to: none;
See Calculator
-outResultsestimated stem models
Type: opals::Vector<opals::StemModel>
Remarks: output
A list of successfully estimated stem models is returned.
StemModel: A StemModel object stores a vector of fitted stem objects. ModuleDBH creates this object for each input stem and attaches all corresponding model fits to it. Furthermore, some administrative information is stored:
  • stem ID of the input stem approximation;
  • optional pass through transfer attribute object that was provided by the user;
  • optional stem based filter object;
-byproductstem points datamanager file
Type: opals::Path
Remarks: optional
This optional datamanager output can be used to further visualize or analyse all points that have been used for the stem fit. Each points contains several attributes including the stem id, the residual vetcor.
-tileSizetile size length for stem points ODM
Type: double
Remarks: optional
This parameter is only relevant if a stem points datamanager is written and the automatic tile size estimation should be overruled. Inappropriate tile sizes can lead to memory problems or slow performance in later processing steps. On average, 100,000 to 200,000 points should be stored within one tile (check index statistics after export).

Examples

Basic example

The data used in the following example are located in the $OPALS_ROOT/demo/dbh directory. This simple example shows a how the module can be used based on given stem approximation.

As a prerequisite, the TLS point cloud of four trees must be imported into the ODM. To achieve that, change to the demo/dbh directory and type:

opalsImport -inf trees.laz

Now, run the following command

opalsDBH -inf trees.odm -approx stem_approx.txt -geometryModel cylinder -outf trees_results.shp -patchlen 0.5

which outputs the shape file trees_results.shp and trees_results.txt. The ASCII file contains the same results but in human readable form. Using opalsView a quick visual 3D check of the results can be performed by runing

opalsView -inf trees.odm -filter "generic[z>101 and z<102]"

and adding the shape file to the view:

Fig. 4: Height section of points including stem fits (colored circles)

Example using individual stem filters and pass through columns

The module allows passing through columns from the approximation to the result files. This feature is useful for preserving pre-computed information or properties of individual stems in the resulting products as shown below. Furthermore, the user can provide individual filter for each stem which is a high-end feature that is only rarely required. Nevertheless, its mode of operation is shown in this example. First, a segmentation of the stem points is performed by the following command,

opalsSegmentation -inf trees.odm -searchRadius 0.07 -minSegSize 25 -method condClustering -sort descendingPSize -filter "generic[z>101 and z<102]"

which segments the 4 stems and 6 branchs within the specified height range. By sorting the segments by their point size, it is secured that the segment ids (attribute segmentId) result from 0 to 3. Assuming known tree heights (pass through attribute), approximation of the stem centers and their corresponding segement id, an approximation file for the four trees can be prepared:

#pt1.x pt1.y pt1.z pt2.x pt2.y pt2.z radius treeHeight filter
-7.581 7.944 101.467 -7.581 7.944 102.467 0.442 12.40 "generic[SegmentId == 0]"
-9.732 8.928 101.510 -9.732 8.928 102.510 0.513 22.40 "generic[SegmentId == 1]"
-9.797 4.386 101.446 -9.797 4.386 102.446 0.485 32.40 "generic[SegmentId == 2]"
-6.122 4.157 101.463 -6.122 4.157 102.463 0.262 42.40 "generic[SegmentId == 3]"

The individual filters secure that only the corresponding segmented stem points will be used for stem fitting, which eliminats outliers and non-stem points for the actual adjustment step. To feed this information into the Module DBH, a corresponding OPALS Format Definition is required.

<opalsFormatDefinition>
<description>X1,Y1,Z1,R1,X2,Y2,Z2,R2 Text Format</description>
<text>
<header text='#pt1.x pt1.y pt1.z pt2.x pt2.y pt2.z radius treeHeight filter' />
<column name='x' type='double' format='7.3' />
<column name='y' type='double' format='7.3' />
<column name='z' type='double' format='9.3' />
<column name='_x2' type='double' format='7.3' />
<column name='_y2' type='double' format='7.3' />
<column name='_z2' type='double' format='9.3' />
<column name='_r' type='double' format='7.3' />
<column name='_treeHeight' type='float' format='10.2' />
<column name='_filter' type='string' typeFile='quotedString' />
</text>
</opalsFormatDefinition>

Columns that do not serve as stem approximations, as described above, will be passed through to the result output table (treeHeight in this case). The OFD-xml file needs to be set as aFormat parameter.

opalsDBH -inf trees.odm -approx approx_with_info.txt -aformat approx_with_info_def.xml -geometryModel cylinder -outf results_with_info.shp -patchlen 0.5

The command above will create results_with_info.shp which containes the following attribute table

Fig. 5: Attribute table of resuting shape file. Right column (treeHeight) was passed through from the approx file

Extended example including stem approximation extraction

The examples above are based on known/given stem approximations. But in general this information is not known in advance. In this extended example, a complete workflow starting from the raw 3D point cloud to the final DBH results including a method for extracting approximations is shown.

Assuming that the 3D point cloud has been imported into an ODM, one needs to first compute a DTM. This helps to subselect the relevant stem points but is also required to determine the breast height for each stem (typically 1.3m above ground). A simple method that results in approximated DTM is to selected the lowest point within 2D cells. The cell size must be selected that large that at least one ground point lies within each cell. Especially, in steep terrain the resulting model will be of low quality, but for the stem fitting it might be good enough.

opalsCell -inf trees.odm -outf final_dtm.tif -cells 0.25 -feat min

To compute a high quality DTM one can use Module RobFilter or Module TerrainFilter (not included in the stable release yet) to classify ground points first, and then compute a detailed and accurate DTM, as shown be the following command sequence

opalsTerrainFilter -inf trees.odm -pyramidLevels 3 -filterThresholds 0.1 0.2 1 -grids 0.25
opalsDTM -inf trees.odm -grids 0.1 -filter "class[ground]"
opalsFillGaps -inf trees_dtm.tif -outf final_dtm.tif -method adaptive

Using Module FillGaps unwanted holes can be filled to get a continuous model as shown in Figure 6.

Fig. 6: Shaded DTM. Left: after opalsDTM. Right: after opalsFillGaps

Next, we computed the attribute normalizedZ for all points which helps to subselect relevant stem points and for finding robust DBH approximations.

opalsAddInfo -inf trees.odm -gridf final_dtm.tif -attr "normalizedZ=z-r[0]"

In the following steps, we pre-filter stem point candidates and export them into a new ODM (trees_sub.odm) to increase robustness and processing speed. We consider echos that have a normalized height between 0.4 and 3m and which are last echos (in TLS solid bark echos are typically last echos). To further increase robustness also the vertical distribution of the points is considered using Module Cell. The vertical range and point count within 5cm 2D cells are computed and attached as attributes to the points. Echos are not consider (as stem points) if the vertical range and point count is less than 1m and 50, respectively.

opalsCell -inf trees.odm -cells 0.05 -feat pcount range -filter "generic[normalizedZ>0.4 and normalizedZ<3] and echo[last]"
opalsAddInfo -inf trees.odm -grid trees_z_pcount.tif trees_z_range.tif -attribute "_zrange=r[1]; _pcount=r[0]" -filter "generic[normalizedZ>0.4 and normalizedZ<3] and echo[last]"
opalsPointStats -inf trees.odm -searchR 0.05 -feat pcount range -filter "generic[normalizedZ>0.4 and normalizedZ<3] and echo[last]"
opalsExport -inf trees.odm -outf trees_sub.odm -filter "generic[normalizedz>0.4 and normalizedZ<3 and _zrange>1.0 and _pcount>50] and echo[last]"

As shown in Figure 7 the previous filter step results in a good pre-selection of possible stem points.

Fig. 7: Extracted stem point candidates

Now its time to extract the actual stem approximation which can be done by using a small python script. The script utilises Module Segmentation to segment each stem using conditional clustering. The CoG of each segment and an estimated radius derived from the x- and y-extend of the segments are stored within a new approximation ODM. Although, the segment CoG can be shifted compared to stem center (especially if the stem has not been fully covered for all sides), this approximation error is usually uncritical for the subsequent process.

from __future__ import print_function # print function syntax as in python 3
from opals import pyDM, Segmentation
import opals
import sys
def treeSeg(input_odm, output_odm):
# perform segmentation using conditional clustering
mod = Segmentation.Segmentation()
mod.inFile = input_odm
mod.searchRadius = 0.15
mod.minSegSize = 25
mod.method = opals.Types.SegmentationMethod.condClustering
mod.sort = opals.Types.SegmentSort.descendingPSize
mod.run()
# get segmentation manager (gives access to the segment objects)
segManager = mod.segments
print("The segmentation process has found " + str(segManager.sizeSegments()) + " tree segment(s)")
odm = pyDM.Datamanager.create(output_odm)
layouter.addColumn(pyDM.ColumnType.uint32,"_pointID")
layouter.addColumn(pyDM.ColumnType.float_,"_r")
layout = layouter.getLayout()
# iterate over all segments and write relevant segments properties
for seg in segManager.getSegments():
pt = pyDM.Point( seg.getCoG().x, seg.getCoG().y, seg.getCoG().z, layout ) # create pt based on segment CoG
dx = (seg.getBBox().xmax - seg.getBBox().xmin) # x extend of segment
dy = (seg.getBBox().ymax - seg.getBBox().ymin) # y extend of segment
pt.info().set(0, seg.getID())
pt.info().set(1, max([dx, dy])/2. ) # store radius estimate
odm.addPoint(pt)
odm.save()
print("'%s' odm with %d stem approximation was been created" % (output_odm, segManager.sizeSegments()))
if __name__ == "__main__":
treeSeg(sys.argv[1], sys.argv[2])

The height value, however, needs to be taken from the DTM rather from the segment itself. Since the DBH measure usually refers to the highest surface point in the stem vicinity, the height determination requires some additional processing. Therefore, the DTM is imported as additional point cloud into the approximation ODM. Then, the highest DTM point is determined within in a 0.25m radius to the stem approximation and attached as attribute _ZMax using Module PointStats. Finally, the actual DBH height is computed by adding 1.3m to _ZMax and stored in attribute _z1 by a corresponding Module AddInfo command. In this data set the stems are nearly vertical, which is why the second point of each stem approximation is simply set 1m above the base point (see attribute _z2). For death wood or titled stems better axis approximation might be necessary.

python treeSeg.py trees_sub.odm stem_approx.odm
opalsImport -inf stem_approx.odm final_dtm.tif -outf stem_approx.odm
opalsPointStats -inf stem_approx.odm -feat max -searchr 0.25 -searchm d2 -refmodel zeroPlane -filter "generic[fileId==0]" "generic[fileId==1]"
opalsAddInfo -inf stem_approx.odm -gridf final_dtm.tif -attr "_z1=_ZMax+1.3; _z2=_ZMax+2.3; _dtm=r[0]" -filter "generic[fileId==0]"

After the command sequence above the ODM contains all necessary information (see Figure 8) to write the actual approximation file.

Fig. 8: Attribute statistics of stem approximations
opalsExport -inf stem_approx.odm -oformat approx_extracted_def.xml -outf approx_extracted.txt -filter "generic[fileId==0]"

Once, the approximation file has been written, Module DBH can be called. When testing and evaluating different parameters and geometry models as shown below, it can be useful to additionally output fitting ODMs (see parameter byproduct). Those ODMs contain all points the where selected for fitting from the input ODM. Pointwise values computed during the robust least-squares adjustment are attached to each point, so e.g. detected outliers, the adjusted point sigma and the residual value can be easily visualised.

opalsDBH -inf trees_sub.odm -approxfile approx_extracted.txt -geometryModel cone -outf dbh_con_20cm.shp -byproduct dbh_con_20cm.odm -patchlength 0.2
opalsDBH -inf trees_sub.odm -approxfile approx_extracted.txt -geometryModel cylinder -outf dbh_cyl_60cm.shp -byproduct dbh_cyl_60cm.odm -patchlength 0.6

The signed fitting residuals of the first stem of the two module runs from above are visualised in Figure 9. Whereas the different patch length is easly visible, the effect of the different fitting model is hardly recognisable since the stem is narrowing very little over the seletect section.

Fig. 9: Byproduct ODM colored by residual. Left: cone model and 20cm patch length. Right: cylinder model and 60cm patch length ODM

GM: Idee für Abschlußsatz?

References

Bruggisser M., Hollaus M., Kükenbrink D., Pfeifer N., 2019, Comparison of forest structure metrics derived from uav lidar and als data, ISPRS Ann. Photogramm. Remote Sens. Spat. Inf. Sci., IV-2/W5, pp. 325-332, 10.5194/isprs-annals-IV-2-W5-325-2019

Bruggisser M, Hollaus M., Otepka J., Pfeifer N., 2020, Influence of ULS acquisition characteristics on tree stem parameter estimation, ISPRS Journal of Photogrammetry and Remote Sensing, Volume 168, pp. 28-40, ISSN 0924-2716, 10.1016/j.isprsjprs.2020.08.002

Lukács, G., Marshall, A. D., Martin, R. R., 1997, Geometric least-squares fitting of spheres,cylinders,cones and tori.

Wieser M., Hollaus M., Mandlburger G., Glira P., Pfeifer N., 2016, ULS LiDAR supported analyses of laser beam penetration from different ALS systems into vegetation, ISPRS Ann. Photogramm. Remote Sens. Spat. Inf. Sci., III-3, pp. 233-239, 10.5194/isprsannals-III-3-233-2016

Author
jo, gm
Date
10.03.2020
@ pcount
number of points used for interpolation of grid point
Definition: GridFeature.hpp:14
@ cone
cone model
@ d2
Search based on 2D coordinates (x and y) only.
Interface to a factory object for creating AddInfo layouts.
Definition: pyDM.py:1332
opalsAddInfo is the executable file of Module AddInfo
Definition: ModuleExecutables.hpp:8
@ condClustering
Parameter group 'condClustering' containing options for conditional clustering (opalsSegmentation)
opalsSegmentation is the executable file of Module Segmentation
Definition: ModuleExecutables.hpp:193
opalsImport is the executable file of Module Import
Definition: ModuleExecutables.hpp:113
3d point object
Definition: pyDM.py:1541
@ method
method of computation (opalsFillGaps, opalsLineTopology)
opalsView is the executable file of Module View
Definition: ModuleExecutables.hpp:248
@ pyramidLevels
Number of data/image pyramid levels (opalsTerrainFilter)
opalsExport is the executable file of Module Export
Definition: ModuleExecutables.hpp:73
@ byproduct
optional output that is not the central result of a module run (opalsSegmentation: segment odm repres...
This is the fake namespace of all opals Python scripts.
Definition: doc/temp_doxy/python/__init__.py:1
opalsFillGaps is the executable file of Module FillGaps
Definition: ModuleExecutables.hpp:78
opalsPointStats is the executable file of Module PointStats
Definition: ModuleExecutables.hpp:163
@ descendingPSize
sort descending by point per segment
opalsCell is the executable file of Module Cell
Definition: ModuleExecutables.hpp:23
@ cylinder
cylinder model
@ adaptive
Uses adaptive plane fit with inverse distance weighting.
@ filter
string to be parsed in construction of DM::IFilter (various modules)
opalsDBH is the executable file of Module DBH
Definition: ModuleExecutables.hpp:38
@ range
scanner range group (opalsStripAdjust)
opalsTerrainFilter is the executable file of Module TerrainFilter
Definition: ModuleExecutables.hpp:228
@ sort
sorting method (opalsSegmentation: sorting of segments)
@ filterThresholds
filter thresholds for hierarchical levels (opalsTerrainFilter)
@ grid
Grid data file.