XPalm.VPalm submodule API

XPalm.VPalm.AbstractGeometryModelType

geometry process abstract model.

All models implemented to simulate the geometry process must be a subtype of this type, e.g. struct MyGeometryModel <: AbstractGeometryModel end.

You can list all models implementing this process using subtypes:

Examples

subtypes(AbstractGeometryModel)
source
XPalm.VPalm.LeafGeometryModelType
LeafGeometryModel(;vpalm_parameters, rng=Random.MersenneTwister())

A PlantSimEngine model that builds the 3D geometry for a leaf, including the petiole, rachis, and leaflets. This model operates at the phytomer scale and modifies the MTG directly.

Arguments

  • vpalm_parameters::Dict{String,Any}: VPalm model parameters.
  • rng::Random.AbstractRNG: Random number generator for stochastic processes.

Inputs

  • height_internodes: Internode height (from plantsimengine_status)
  • radius_internodes: Internode radius (from plantsimengine_status)
  • biomass_leaves: Leaf biomass (from plantsimengine_status)
  • rank_leaves: Leaf rank (from plantsimengine_status)

Outputs

This model has no outputs as it modifies the MTG directly by adding geometric properties and child nodes.

Notes

The model requires access to the VPalm parameters via the parameters dictionary under the "vpalm" key.

source
PlantSimEngine.run!Method
run!(model, models, status, meteo, constants, node)

Builds the 3D geometry for a leaf by adding internode properties and creating child nodes for petiole, rachis, and leaflets.

Arguments

  • model::LeafGeometryModel: The leaf geometry model
  • models: A ModelList struct holding the parameters for the model
  • status: The status of the model with inputs (height, radius, biomass, rank)
  • meteo: Meteorology structure (not used by this model)
  • constants: Physical constants (not used by this model)
  • node: MTG node of the phytomer

Notes

The model expects node to be the phytomer MTG node and accesses VPalm parameters from model.vpalm_parameters.

source
XPalm.VPalm.add_geometry!Method
add_geometry!(
    mtg, refmesh_cylinder, refmesh_snag, refmesh_plane;
    snag_width=0.20u"m", # see defaultOrthotropyAttribute in the trunk in the java implementation
    snag_height=0.15u"m",
    snag_length=3.0u"m",
)

Adds geometry to the MTG (Multiscale Tree Graph) for the oil palm plant architecture, i.e. compute the meshes.

source
XPalm.VPalm.add_leaflet_geometry!Method
add_leaflet_geometry!(
    leaflet_node,
    internode_width,
    internode_height,
    rachis_position,
    rachis_orientation,
    rachis_rotation,
    stem_bending,
    refmesh_plane
)

Create the leaflet geometry based on its segments.

Arguments

  • leaflet_node: The MTG node of the leaflet
  • internode_width: Width of the internode (used for positioning)
  • internode_height: Height of the internode (used for positioning)
  • rachis_position: Position of the rachis section where the leaflet is attached
  • rachis_orientation: Orientation angles [zenithal, azimuthal, torsion] of the rachis section
  • rachis_rotation: Rotation of the rachis due to phyllotaxy (degrees)
  • stem_bending: Bending of the stem (degrees)
  • refmesh_plane: Reference mesh used for the planar leaflet segments

Returns

  • Nothing (the geometry is added directly to the leaflet node and its segments)
source
XPalm.VPalm.add_section_geometry!Function
add_section_geometry!(
    node, internode_width, internode_height, internode_phyllotaxy, stem_bending, 
    refmesh_cylinder, position_section=Ref(Meshes.Point(0.0, 0.0, 0.0)), angles=[0.0, 0.0, 0.0],
    type::String,
)

Create the petiole/rachis sections geometry based on their dimensions.

Arguments

  • node: the MTG node of the petiole/rachis
  • refmesh_cylinder: the reference mesh used for a cylinder (PlantGeom.RefMesh)
  • internode_width: the width of the internode on the stipe (m)
  • internode_height: the height of the internode on the stipe (m)
  • internode_phyllotaxy: the phyllotaxy of the internode on the stipe (°)
  • stem_bending: the bending of the stipe (°)
  • type::String: the type of the section ("PetioleSegment" or "RachisSegment")
  • position_section=Ref(Meshes.Point(0.0, 0.0, 0.0)): the position of the section relative to the first one.
source
XPalm.VPalm.bendMethod
bend(
    type, width_bend, height_bend, init_torsion, x, y, z, mass_rachis, mass_leaflets_right, mass_leaflets_left,
    distance_application, elastic_modulus, shear_modulus, step, npoints, nsegments;
    all_points=false,
    angle_max=deg2rad(21),
    force=true,
    verbose=true
)

Compute the deformation of the rachis by applying both bending and torsion.

Arguments

  • type: Vector of section types (1: triangle bottom, 2: rectangle, 3: triangle top, 4: ellipse, 5: circle).
  • width_bend: Vector of segment section widths (m).
  • height_bend: Vector of segment section heights (m).
  • init_torsion: Vector of initial torsion angles (degrees).
  • x: Vector of x coordinates of the segments.
  • y: Vector of y coordinates of the segments.
  • z: Vector of z coordinates of the segments.
  • mass_rachis: Vector of rachis segment masses (kg).
  • mass_leaflets_right: Vector of leaflet masses carried by the segment, on the right side (kg).
  • mass_leaflets_left: Vector of leaflet masses carried by the segment, on the left side (kg).
  • distance_application: Vector of application distances for the left and right weights (m).
  • elastic_modulus: Vector of elasticity moduli (bending, MPa).
  • shear_modulus: Vector of shear moduli (torsion, MPa).
  • step: Length of the segments that discretize the object (m).
  • npoints: Number of points used in the grid discretizing the section.
  • nsegments: Number of segments dividing the rachis to compute the torsion and bending.
  • all_points=false: return all points used in the computation (true), or only the input points corresponding to x, y and z coordinates (false, default).
  • angle_max=deg2rad(21): Maximum angle for testing the small displacement hypothesis (radians).
  • force=true: Check if verify the small displacements hypothesis and bounds the values to be at maximum angle_max
  • verbose=true: Provide information during computation.

Returns

Named tuple with geometrical fields describing the rachis bended and with torsion applied

  • x: x coordinates of the points.
  • y: y coordinates of the points.
  • z: z coordinates of the points.
  • length: length of the segments.
  • angle_xy: angle between the xy-plan and the segment.
  • angle_xz: angle between the xz-plan and the segment.
  • torsion: torsion angle of the segment.

All these fields are vectors of the same length as the input vectors (i.e. number of segments).

Details

The bending and torsion are applied to the sections of the rachis defined by 5 segments.

source
XPalm.VPalm.beta_distribution_normMethod
beta_distribution_norm(x, xm, ym)

Calculate the normalized beta distribution value at point x. This is the exact implementation from the Java version.

Arguments

  • x: Position [0 to 1].
  • xm: Mode of the beta distribution.
  • ym: Maximum value of the beta distribution.

Returns

  • Normalized beta distribution value.
source
XPalm.VPalm.beta_distribution_norm_integralMethod
beta_distribution_norm_integral(xm, ym)

Calculate the integral (area) of the normalized beta distribution. Equivalent to betaDistributionNormIntegral in the Java version.

Arguments

  • xm: Mode of the beta distribution.
  • ym: Value of the function at the mode.

Returns

  • Approximate area under the beta distribution curve.
source
XPalm.VPalm.biomechanical_properties_rachisMethod
biomechanical_properties_rachis(
    rachis_twist_initial_angle, rachis_twist_initial_angle_sdp,
    elastic_modulus, shear_modulus, rachis_length,
    leaflet_length_at_b_intercept, leaflet_length_at_b_slope, relative_position_bpoint,
    relative_position_bpoint_sd, relative_length_first_leaflet, relative_length_last_leaflet, relative_position_leaflet_max_length,
    rachis_fresh_weight, rank, height_cpoint, zenithal_cpoint_angle, nb_sections,
    height_rachis_tappering,
    points, iterations, angle_max;
    verbose, rng
)

Use of the biomechanical model to compute the properties of the rachis.

Arguments

  • rachis_twist_initial_angle: initial twist angle of the rachis (°)
  • rachis_twist_initial_angle_sdp: standard deviation of the initial twist angle of the rachis (°)
  • elastic_modulus: elastic modulus of the rachis (Pa)
  • shear_modulus: shear modulus of the rachis (Pa)
  • rachis_length: length of the rachis (m)
  • leaflet_length_at_b_intercept: intercept of the linear function for the leaflet length at the B point (m)
  • leaflet_length_at_b_slope: slope of the linear function for the leaflet length at the B point (m)
  • relative_position_bpoint: relative position of the B point on the rachis (0: base to 1: tip)
  • relative_position_bpoint_sd: standard deviation of the relative position of the B point on the rachis
  • relative_length_first_leaflet: relative length of the first leaflet on the rachis (0 to 1)
  • relative_length_last_leaflet: relative length of the last leaflet on the rachis (0 to 1)
  • relative_position_leaflet_max_length: relative position of the longest leaflet on the rachis (0.111 to 0.999)
  • rachis_fresh_weight: fresh weight of the rachis (kg)
  • rank: rank of the rachis
  • height_cpoint: height of the C point (m)
  • zenithal_cpoint_angle: zenithal angle of the C point (°)
  • nb_sections: number of sections to compute the bending
  • height_rachis_tappering: tappering factor for the rachis height
  • npoints_computed: number of points to compute the bending
  • iterations: number of iterations to compute the bending
  • angle_max: maximum angle to compute the bending (°)
  • verbose: display information about the computation (e.g. checks on the units)
  • rng: the random number generator

Returns

A named tuple with the following fields:

  • length: vector with the length of each segment
  • points_positions: the position of the points along the rachis
  • bending: the bending angle of the rachis
  • deviation: the deviation of the rachis (angle in the xz plane)
  • torsion: the torsion of the rachis
  • x: the x coordinates of the rachis
  • y: the y coordinates of the rachis
  • z: the z coordinates of the rachis

Details

Split the rachis into 5 segments defined by remarkable points (C, C-B, B, B-A, A). Each segment has a particular shape, a mass, and the leaflets on both sides of the rachis have a mass. Coefficents are used to compute the mass distribution and relative lengths of segments. The rachis is bent using the bend function.

source
XPalm.VPalm.build_mockupMethod
build_mockup(parameters; merge_scale=:leaflet, rng=Random.MersenneTwister(parameters["seed"]))

Construct a mockup of an oil palm plant architecture using the specified parameters.

Arguments

  • parameters::Dict: Dictionary containing model parameters for the oil palm plant architecture.
  • merge_scale::Symbol: (optional) The scale at which to merge geometry.
    • :node: Geometry is not merged, each node has its own mesh (finer scale is leaflet segments).
    • :leaflet (default): Geometry is merged at the leaflet level.
    • :leaf: All geometry for a leaf is merged into a single mesh.
    • :plant: All plant geometry is merged into a single mesh.
  • rng: (optional) A random number generator for stochastic processes. Defaults to a Mersenne Twister seeded with the value in parameters["seed"]. If set to nothing, randomness is disabled (useful for testing).

Description

The merge_scale argument controls how the geometry is structured within the Multiscale Tree Graph (MTG). The resulting mesh is identical in all cases, but its organization differs.

  • Using :leaflet retains the finest detail, with each leaflet having its own mesh. This is best for analyses like light interception at the organ level.
  • Using :leaf or :plant merges geometry into larger components. A single mesh for the whole plant (:plant) is the most performant for rendering, but it prevents querying information for individual organs from the mesh (e.g., which part of the mesh is a given leaflet).

Returns

  • mtg: An MTG (Multiscale Tree Graph) representing the oil palm plant architecture, including geometry at the specified merge scale.

Example

using XPalm.VPalm
file = joinpath(dirname(dirname(pathof(XPalm))), "test", "references", "vpalm-parameter_file.yml")
parameters = read_parameters(file)
mtg = build_mockup(parameters; merge_scale=:plant)
source
XPalm.VPalm.c_point_angleMethod
c_point_angle(leaf_rank, cpoint_decli_intercept, cpoint_decli_slope, cpoint_angle_SDP; rng)

Compute the angle at the C point of the leaf.

Arguments

  • leaf_rank: Rank of the leaf
  • cpoint_decli_intercept: Intercept of the linear relationship between leaf rank and C point declination
  • cpoint_decli_slope: Slope of the linear relationship
  • cpoint_angle_SDP: Standard deviation of the C point angle
  • rng: Random number generator

Returns

The zenithal angle at the C point of the leaf (°)

source
XPalm.VPalm.calculate_segmentFunction
calculate_segment(relative_position, num_segments=10)

Calculate the segment index for a given relative position along the rachis.

Arguments

  • relative_position: Relative position along the rachis [0 to 1), where 0 is the base and 1 is the tip.
  • num_segments: Number of segments the rachis is divided into (default: 10).

Details

We divide the rachis into segments to capture variations in properties along its length. This function:

  1. Converts a continuous relative position (0-1) into a discrete segment index
  2. Ensures the segment index is within valid bounds (1 to num_segments)

Biological Context

The palm rachis exhibits changing properties along its length, including:

  • Leaflet grouping patterns
  • Leaflet sizes and angles

Dividing the rachis into discrete segments allows the model to represent these gradual changes in a computationally efficient manner. Each segment can have different parameter values that together create the characteristic patterns seen in real palms.

Returns

The segment index (starts at 1 in Julia).

source
XPalm.VPalm.calculate_segment_anglesMethod
calculate_segment_angles(young_modulus, initial_angle, leaflet_length, tapering, segment_positions)

Calculate the global angles for each segment of a bent leaflet based on the Young's modulus model.

Arguments

  • young_modulus: Value of Young's modulus
  • initial_angle: Initial angle from vertical in radians
  • leaflet_length: Total length of the leaflet
  • tapering: Tapering factor
  • segment_positions: Array of segment boundary positions (normalized 0-1)

Returns

  • Array of segment angles in radians
source
XPalm.VPalm.compute_leaf_rankFunction
compute_leaf_rank(nb_internodes, index_leaf)

Compute the rank of a leaf based on the total number of internodes and the index of the leaf.

Arguments

  • nb_internodes: The total number of internodes until leaf of rank 1.
  • index_leaf: The index of the leaf.
  • leaves_in_sheath: The number of leaves in the sheath, i.e. with rank < 1 (default is 0).

Note

This is a simple leaf rank, not considering the leaves of rank <= 0.

Returns

The leaf rank, i.e. 1 for the first opened leaf, 2 for the second leaf, etc.

source
XPalm.VPalm.compute_leaflet_type_frequenciesMethod
compute_leaflet_type_frequencies(leaflet_frequency_high, leaflet_frequency_low)

Compute the frequency of leaflet type within the sub-sections of a rachis.

Arguments

  • leaflet_frequency_high: Vector of frequency values for the +1 leaflet types (high) along the rachis sub-sections.
  • leaflet_frequency_low: Vector of frequency values for the -1 leaflet types (low) along the rachis sub-sections..

Note that the length of the two vectors must be the same. It will define how many sub-sections the rachis is divided into for this computation.

Returns

A vector of NamedTuples representing the (;high, medium, low) frequencies for each sub-section.

source
XPalm.VPalm.compute_number_of_leafletsMethod
compute_number_of_leaflets(rachis_final_length, nb_max, nb_slope, nb_infl, nbLeaflets_SDP; rng)

Compute the number of leaflets based on the logistic function, a standard deviation and a minimum value allowed.

Arguments

  • rachis_final_length: Final length of the rachis (m).
  • nb_max: Maximum number of leaflets.
  • nb_min: Minimum number of leaflets.
  • nb_slope: Slope parameter for the logistic function (leaflet m⁻¹).
  • nb_infl: Inflection point parameter for the logistic function (m).
  • nbLeaflets_SDP: Standard deviation of the normal distribution for the number of leaflets.
  • rng: Random number generator.

Returns

The computed number of leaflets (integer).

source
XPalm.VPalm.compute_properties_internode!Method
compute_properties_internode!(node, index, nb_internodes, rank, stem_height, stem_diameter, parameters, rng)

Computes the mtg properties of an internode.

Arguments

  • node: the internode node
  • index: the index of the internode
  • nb_internodes: the total number of internodes
  • rank: the rank of the internode
  • stem_height: the height of the stem (m)
  • stem_diameter: the diameter of the stem (m)
  • parameters: the parameters of the model
  • rng: the random number generator

Returns

The internode node updated with properties.

Details

The internode dimensions are computed based on the dimensions of the stem and the parameters of the model:

  • width: width of the internode (m)
  • diameter: diameter of the internode (m)
  • length: length of the internode (m)
  • rank: rank of the internode
  • Orthotropy: orthotropy of the internode (set as a constant value)
  • XEuler: Euler / phyllotactic angle of the internode (rad)

Examples

using XPalm.VPalm
using Unitful

parameters = VPalm.default_parameters()
nb_internodes = parameters["nb_leaves_emitted"] + parameters["nb_internodes_before_planting"] # The number of internodes emitted since the seed
# Plant / Scale 1
plant = Node(NodeMTG("/", "Plant", 1, 1))
# Stem (& Roots) / Scale 2
stem = Node(plant, NodeMTG("+", "Stem", 1, 2))
compute_properties_stem!(stem, parameters, 3.0u"m"; rng=rng)
stem_height = stem[:stem_height]
stem_diameter = stem[:stem_diameter]
# Phytomer / Scale 3
phytomer = Node(stem, NodeMTG("/", "Phytomer", 1, 3))
# Internode & Leaf / Scale 4
internode = Node(phytomer, NodeMTG("/", "Internode", 1, 4))
compute_properties_internode!(internode, 1, nb_internodes, stem_height, stem_diameter, parameters, rng)
source
XPalm.VPalm.compute_properties_leaf!Method
compute_properties_leaf!(node, leaf_rank, final_length, parameters, rng)

Compute the properties of a leaf node:

  • zenithalinsertionangle: the zenithal insertion angle of the leaf (rad)
  • rachis_length: the length of the rachis (m)
  • zenithalcpointangle: the zenithal angle at C-point (rad)

Arguments

  • node: the leaf node
  • leaf_rank: the rank of the leaf
  • final_length: the final length of the leaf (m)
  • parameters: the parameters of the model
  • rng: the random number generator

Returns

The leaf node updated with properties.

Details

The leaf dimensions are computed based on the dimensions of the stem and the parameters of the model:

  • zenithalinsertionangle: the zenithal insertion angle of the leaf (rad). Uses the VPalm.leaf_insertion_angle function.
  • rachislength: the length of the rachis (m). Uses the `rachisexpansion` function.
  • zenithalcpointangle: the zenithal angle at C-point (rad). Uses the c_point_angle function.

Examples

using XPalm.VPalm
using Unitful

file = joinpath(dirname(dirname(pathof(VPalm))), "test", "files", "parameter_file.yml")
parameters = read_parameters(file)
nb_internodes = parameters["nb_leaves_emitted"] + parameters["nb_internodes_before_planting"] # The number of internodes emitted since the seed
nb_leaves_alive = floor(Int, mean_and_sd(parameters["nb_leaves_mean"], parameters["nb_leaves_sd"], rng))
nb_leaves_alive = min(nb_leaves_alive, nb_internodes)
# Plant / Scale 1
plant = Node(NodeMTG("/", "Plant", 1, 1))
# Stem (& Roots) / Scale 2
stem = Node(plant, NodeMTG("+", "Stem", 1, 2))
compute_properties_stem!(stem, parameters, 3.0u"m"; rng=rng)
stem_height = stem[:stem_height]
stem_diameter = stem[:stem_diameter]
# Phytomer / Scale 3
phytomer = Node(stem, NodeMTG("/", "Phytomer", 1, 3))
# Internode & Leaf / Scale 4
internode = Node(phytomer, NodeMTG("/", "Internode", 1, 4))
leaf = Node(internode, NodeMTG("+", "Leaf", 1, 4))
compute_properties_leaf!(leaf, 1, nb_internodes, nb_leaves_alive, parameters, rng)
source
XPalm.VPalm.compute_properties_petiole!Method
compute_properties_petiole!(
    petiole_node,
    insertion_angle, rachis_length, zenithal_cpoint_angle,
    width_base, height_base, cpoint_width_intercept,
    cpoint_width_slope, cpoint_height_width_ratio,
    petiole_rachis_ratio_mean,
    petiole_rachis_ratio_sd, nb_sections;
    rng=Random.MersenneTwister(1)
)

Compute the dimensional properties of a petiole.

Arguments

  • petiole_node: the MTG Node of the petiole
  • insertion_angle: the angle of insertion of the petiole on the stem (°)
  • rachis_length: the length of the rachis (m)
  • zenithal_cpoint_angle: the zenithal angle of the C point of the petiole, i.e. the tip (°)
  • width_base: the width of the petiole at its base (m)
  • height_base: the height of the petiole at its base (m)
  • cpoint_width_intercept: the intercept of the linear function for the width at the C point (m)
  • cpoint_width_slope: the slope of the linear function for the width at the C point
  • cpoint_height_width_ratio: the ratio of the height to width at the C point
  • petiole_rachis_ratio_mean: the mean ratio of the petiole to rachis length
  • petiole_rachis_ratio_sd: the standard deviation of the ratio of the petiole to rachis length
  • nb_sections: the number of sections discretizing the petiole
  • rng=Random.MersenneTwister(1): the random number generator

Returns

The petiole node updated with properties.

Details

Properties are computed based on the allometries of the petiole and the rachis:

  • length: the length of the petiole (m)
  • azimuthal_angle: the azimuthal angle of the petiole (°)
  • width_base: the width of the petiole at its base (m)
  • height_base: the height of the petiole at its base (m)
  • width_cpoint: the width of the petiole at the C point (m)
  • height_cpoint: the height of the petiole at the C point (m)
  • zenithalinsertionangle: the zenithal angle of insertion of the petiole on the stem (°)
  • zenithalcpointangle: the zenithal angle of the C point of the petiole (°)
  • section_length: the length of the petiole sections (m)
  • sectioninsertionangle: the zenithal angle of insertion between the petioles sections (°)
source
XPalm.VPalm.compute_properties_petiole_section!Method
compute_properties_petiole_section!(petiole_node, section_node, index, nb_sections)

Compute the dimension of a petiole section based on the dimensions of the petiole.

Arguments

  • petiole_node: the MTG Node of the petiole
  • section_node: the MTG Node of the section to be computed
  • index: the index of the section on the petiole, from 1 at the base to nb_sections.
  • nb_sections: the number of sections discretizing the petiole
  • section_insertion_angle: the zenithal angle of the petioles sections (global angle, °)

Returns

The section node updated with dimensional properties.

Details

The petiole_node should have the following attributes:

  • width_base: the width of the petiole at its base (m)
  • height_base: the height of the petiole at its base (m)
  • width_cpoint: the width of the petiole at the C point (m)
  • height_cpoint: the height of the petiole at the C point (m)
  • section_length: the length of the petiole sections (m)
  • insertion_angle: the angle of insertion of the petiole on the stem (°)
  • section_insertion_angle: the zenithal angle of insertion between the petioles sections (°)
  • azimuthal_angle: the azimuthal angle at the insertion (°)
source
XPalm.VPalm.compute_properties_stem!Method
compute_properties_stem!(node, parameters, length_reference_leaf; rng=MersenneTwister(1234))

Compute the properties of the stem node.

Arguments

  • node: the stem node
  • parameters: the parameters of the Vpalm model
  • length_reference_leaf: the length of the reference leaf (usually, rank 17)
  • rng=MersenneTwister(1234): the random number generator

Returns

The stem node updated with properties.

Details

The stem dimensions are computed based on the parameters of the model:

  • stembending: the bending of the stem. Uses the `VPalm.stembending` function.
  • stemheight: the height of the stem. Uses the `VPalm.stemheight` function.
  • stemdiameter: the diameter of the stem. Uses the `VPalm.stemdiameter` function.

Examples

using XPalm.VPalm
using Unitful

file = joinpath(dirname(dirname(pathof(VPalm))), "test", "references", "vpalm-parameter_file.yml")
parameters = read_parameters(file)
nb_internodes = parameters["nb_leaves_emitted"] + parameters["nb_internodes_before_planting"] # The number of internodes emitted since the seed
nb_leaves_alive = floor(Int, mean_and_sd(parameters["nb_leaves_mean"], parameters["nb_leaves_sd"], rng))
nb_leaves_alive = min(nb_leaves_alive, nb_internodes)
# Plant / Scale 1
plant = Node(NodeMTG("/", "Plant", 1, 1))
# Stem (& Roots) / Scale 2
stem = Node(plant, NodeMTG("+", "Stem", 1, 2))
compute_properties_stem!(stem, parameters, 3.0u"m"; rng=rng)
source
XPalm.VPalm.create_leaflet_segments!Method
create_leaflet_segments!(
    unique_mtg_id,
    leaflet_node,
    scale,
    leaflet_length,
    width_max,
    stiffness,
    tapering,
    leaflet_relative_pos;
    xm_intercept, xm_slope,
    ym_intercept, ym_slope  
)

Create the segments that make up a leaflet with proper shape and bending properties.

Arguments

  • unique_mtg_id: Reference to the unique ID counter
  • leaflet_node: Parent leaflet node
  • scale: MTG scale for the segments
  • leaflet_length: Total length of the leaflet in meters
  • width_max: Maximum width of the leaflet in meters
  • stiffness: Stiffness value (Young's modulus) for biomechanical bending
  • tapering: Tapering factor (how width decreases along length)
  • leaflet_relative_pos: Relative position of the leaflet on the rachis (0-1)
  • xm_intercept, xm_slope: Parameters for defining maximum leaflet width position
  • ym_intercept, ym_slope: Parameters for defining maximum leaflet width value

Returns

Nothing (segments are added directly to the leaflet node as children)

source
XPalm.VPalm.create_leaflets_for_side!Method
create_leaflets_for_side!(
    unique_mtg_id,
    rachis_node,
    scale,
    leaf_rank,
    rachis_length,
    nb_rachis_sections,
    leaflets_position,
    leaflets,
    leaflet_max_length,
    leaflet_max_width,
    side,
    parameters;
    last_rank_unfolding=2,
    rng=Random.MersenneTwister(1234)
)

Create leaflets for one side of the palm frond rachis.

Arguments

  • unique_mtg_id: Reference to the unique ID counter for MTG nodes
  • rachis_node: Root node of the rachis
  • scale: MTG scale for leaflets
  • leaf_rank: Rank of the leaf (affects unfolding for young fronds)
  • rachis_length: Total length of the rachis in meters
  • nb_rachis_sections: Number of segments dividing the rachis
  • leaflets_position: Array of positions along the rachis for each leaflet
  • leaflets: NamedTuple with leaflet grouping information (group, group_size, plane)
  • leaflet_max_length: Maximum length of leaflets (length of the longest leaflet)
  • leaflet_max_width: Maximum width of leaflets (width of the widest leaflet)
  • side: Side of rachis (1=right, -1=left)
  • parameters: Model parameters
  • last_rank_unfolding=2: Rank at which leaflets are fully unfolded (default is 2)
  • rng=Random.MersenneTwister(1234): Random number generator

Returns

Nothing (leaflets are attached directly to the rachis node in the MTG structure)

source
XPalm.VPalm.create_sectionMethod
create_section(section, section_type)

Fill in the matrix according to the section shape.

Arguments

  • section: Section matrix.
  • section_type: Section type (1: triangle bottom, 2: rectangle, 3: triangle top, 4: ellipse, 5: circle).

Returns

  • The filled section matrix with 1s for cells inside the shape and 0s outside.
source
XPalm.VPalm.create_single_leafletMethod
create_single_leaflet(
    unique_mtg_id,
    index,
    scale,
    leaf_rank,
    leaflet_relative_pos,
    norm_leaflet_rank,
    plane,
    side,
    leaflet_max_length,
    leaflet_max_width,
    parameters;
    offset=0.0,
    last_rank_unfolding=2,
    rng=Random.MersenneTwister(1234)
)

Create a single leaflet with properly computed angles, dimensions and segments.

Arguments

  • unique_mtg_id: Reference to the unique ID counter
  • index: Index for the leaflet node (for identification in MTG)
  • scale: MTG scale level for the leaflet
  • leaf_rank: Rank of the leaf (affects unfolding for young leaves)
  • leaflet_relative_pos: Relative position of leaflet on rachis (0 to 1)
  • norm_leaflet_rank: Normalized rank of the leaflet (0 to 1)
  • plane: Plane type of leaflet (1=high/upward, 0=medium/horizontal, -1=low/downward)
  • side: Side of the leaf (1=right, -1=left)
  • leaflet_max_length: Maximum leaflet length in meters (length of the longest leaflet)
  • leaflet_max_width: Maximum leaflet width in meters (width of the widest leaflet)
  • parameters: Model parameters dictionary
  • offset: Offset from the start of parent node (when applicable)
  • last_rank_unfolding: Rank at which leaflets are fully unfolded (default is 2)
  • rng: Random number generator

Returns

The created leaflet node with all its segment children

source
XPalm.VPalm.cylinderMethod
cylinder()
cylinder(r, l)

Returns a normalized cylinder mesh, or a cylinder with radius r and length l.

Arguments

  • r: The radius of the cylinder.
  • l: The length of the cylinder.
source
XPalm.VPalm.default_parametersMethod
default_parameters(; type="static")

Returns a dictionary of default parameters for the VPalm model.

Arguments

  • type: The type of parameters to return, either "static" or "dynamic". Default is "static".

Details

VPalm can be used in two modes:

  • "static": For static plant architecture, where the plant structure does not change over time. The parameters are measured from one or several real oil palm plants and are used to build mockups of the plant architecture,

which can then be used for simulations or visualizations around this age.

  • "dynamic": For dynamic plant architecture, where the plant structure can change over time (e.g., growth, environmental effects). This is typically used for simulations that involve plant growth over time (like XPalm), or for digital twins of oil palm plants.

Example

default_params = default_parameters()
source
XPalm.VPalm.dist_and_angles_to_xyzMethod
dist_and_angles_to_xyz(dist_p2p1, vangle_xy, vangle_xz)

Transform distances and angles into point coordinates.

Arguments

  • dist_p2p1: Vector of segment lengths (m).
  • vangle_xy: Vector of angles between the segment and the XY plane (radians).
  • vangle_xz: Vector of angles between the segment and the XZ plane (radians).

Returns

The points as a vector of Meshes.Point.

source
XPalm.VPalm.draw_group_sizeMethod
draw_group_size(group, leaflet_type_frequencies, rng)

Determine the size of a leaflet group based on the relative position along the rachis and frequency patterns.

Arguments

  • group: Index of the leaflet group based on its relative position on the rachis (1 to length(leaflet_type_frequencies)).
  • leaflet_type_frequencies: Vector of NamedTuples representing frequency distributions for each rachis segment, with fields:
    • high: Frequency of plane=+1 leaflets (first leaflet in each group), i.e. leaflets on "high" position
    • medium: Frequency of plane=0 leaflets (intermediate leaflets in groups), i.e. leaflets on "medium" position, horizontally inserted on the rachis
    • low: Frequency of plane=-1 leaflets (terminal leaflets in groups), i.e. leaflets on "low" position
  • rng: Random number generator for stochastic determination, or nothing for deterministic behavior.

Details

This function implements an inverse relationship between the frequency of high (plane=1) leaflets and group size, modeling a fundamental biological pattern in palm frond architecture:

  • Segments with high frequency of high leaflets produce many small groups of leaflets
  • Segments with low frequency of high leaflets produce fewer, larger groups of leaflets

When rng is provided, the calculation uses a probabilistic rounding mechanism to ensure proper statistical distribution of group sizes. When rng is nothing, it uses deterministic rounding to the nearest integer. This creates the natural variation in leaflet grouping patterns seen along real palm fronds, where clustering patterns change systematically from base to tip.

Returns

An integer representing the number of leaflets in the group.

source
XPalm.VPalm.draw_planeMethod
draw_plane(frequencies, rng::Nothing)

Draw the plane position for a leaflet based on its frequencies.

Arguments

  • frequencies: A NamedTuple with fields high, medium, and low representing the frequency of each leaflet type in the current rachis segment.
  • rng: A random number generator (optional, can be Nothing for deterministic behavior

Returns

  • The plane position for the leaflet, which can be:
    • 1 for high position (upward)
    • 0 for medium position (horizontal)
    • -1 for low position (downward)
source
XPalm.VPalm.elliptical_cylinderMethod
elliptical_cylinder(r1, r2, l)

Create an elliptical cylinder mesh.

Arguments

  • r1: The radius of the cylinder in the x direction.
  • r2: The radius of the cylinder in the y direction.
  • l: The length of the cylinder.
source
XPalm.VPalm.exponetialMethod
exponetial(x, a, b)

Compute an exponential function at given x value.

Arguments

  • x: The input value.
  • a: The coefficient a of the exponential function.
  • b: The coefficient b of the exponential function.

Note

The exponential function is defined as a * exp(b * x).

source
XPalm.VPalm.final_angleMethod
final_angle(young_modulus, z_angle, length, tapering)

Calculate the maximal deformation angle of a beam.

Arguments

  • young_modulus: Value of Young's modulus
  • z_angle: Angle from vertical (upright) in radians
  • length: Length of the beam where the load is applied
  • tapering: Tapering factor of the beam

Returns

  • The final angle from vertical at the cantilever extremity (in radians)
source
XPalm.VPalm.group_leafletsMethod
group_leaflets(leaflets_relative_position, leaflets_type_frequency, rng)

Compute the group, group size and plane positions of each leaflet along the rachis.

Arguments

  • leaflets_relative_position: Array of relative positions for the leaflets along the rachis (see relative_leaflet_position()).
  • leaflets_type_frequency: Vector of NamedTuples representing frequency distributions along the rachis (if e.g. 10 values are provided, it means the rachis is divided into 10 sub-sections), with fields:
    • high: Frequency of plane=+1 leaflets (first leaflet in each group), i.e. leaflets on "high" position
    • medium: Frequency of plane=0 leaflets (intermediate leaflets in groups), i.e. leaflets on "medium" position, horizontally inserted on the rachis
    • low: Frequency of plane=-1 leaflets (terminal leaflets in groups), i.e. leaflets on "low" position
  • rng: Random number generator.

Details

This function:

1. Organizes leaflets into groups based on position-dependent size distributions
2. Assigns a spatial plane to each leaflet within a group:
    - The first leaflet in each group is always placed on the high position (plane=1)
    - Subsequent leaflets are positioned on medium (plane=0) or low (plane=-1) positions based on their frequency distribution at that rachis segment

Biological Context

Grouping of leaflets is a key morphological feature in palm species, particularly in oil palm (Elaeis guineensis). Unlike some palms with regularly spaced leaflets, oil palms exhibit distinctive clustering patterns where:

  1. Leaflets occur in groups of variable sizes, but typically around 3 leaflets per group
  2. Within each group, leaflets emerge at different angles:
    • The first leaflet points upward (high position)
    • Others point horizontally or downward (medium and low positions)
  3. The pattern of grouping changes along the rachis:
    • Closer to the base: typically larger groups with more leaflets
    • Toward the tip: smaller groups or single leaflets

The model uses an inverse relationship between high-position leaflet frequency and group size to recreate the natural variation in leaflet insertion angle - sections with many high-position leaflets have smaller groups (but more of them), while sections with few high-position leaflets form larger groups.

The grouping pattern changes along the rachis, creating the characteristic appearance of palm fronds with varying leaflet arrangement patterns from base to tip.

Returns

A NamedTuple containing arrays for:

  • group: Group identifier for each leaflet
  • group_size: Size of the group that each leaflet belongs to
  • plane: Spatial position/orientation of each leaflet (1=high, 0=medium, -1=low)
source
XPalm.VPalm.height_to_width_ratioMethod
height_to_width_ratio(x, ratio_point_c, ratio_point_a, pos_ratio_max, ratio_max)

Computes the relative width along the rachis.

Arguments

  • x: relative position on the rachis
  • ratio_point_c: ratio at point C
  • ratio_point_a: ratio at point A
  • pos_ratio_max: relative position of the maximum value of the ratio
  • ratio_max: maximum ratio value
source
XPalm.VPalm.inertia_flex_rotaFunction
inertia_flex_rota(base_width, height, orientation_angle, section_type, grid_size = 100)

Compute the inertia of bending and torsion, and the cross-section area.

Arguments

  • base_width: Dimension of the base.
  • height: Dimension of the height.
  • orientation_angle: Section orientation angle (torsion, in radians).
  • section_type: Section type (see details).
  • grid_size: Number of discretizations (default to 100).

Details

For the section type, possible values are:

  • section_type = 1: triangle (bottom-oriented)
  • section_type = 2: rectangle
  • section_type = 3: triangle (top-oriented)
  • section_type = 4: ellipse
  • section_type = 5: circle

Returns

  • A NamedTuple with fields:
    • ig_flex: Bending inertia.
    • ig_tor: Torsion inertia.
    • sr: Cross-section surface.
source
XPalm.VPalm.init_attributes_seed!Method
init_attributes_seed!(plant, parameters; rng=Random.MersenneTwister(parameters["seed"]))

Initialize the attributes of a palm plant seed (one internode with one leaf), based on the provided parameters.

source
XPalm.VPalm.internode_diameterMethod
internode_diameter(internode_index, rank, stem_diameter, stem_base_shrinkage, stem_top_shrinkage)

Computes the diameter of an internode at a given rank.

Arguments

  • internode_index: The index of the internode.
  • rank: The rank of the internode.
  • stem_diameter: The diameter of the stem at the base.
  • stem_base_shrinkage: The shrinkage coefficient at the stem base.
  • stem_top_shrinkage: The shrinkage coefficient at the stem top.

Returns

The diameter of the internode (m).

Details

A shrinking function is applied to the stem base and top to compute the diameter of the internode.

source
XPalm.VPalm.internode_lengthMethod
Internode length model

Computes the length of an internode at a given rank.

Arguments

  • i/ internode_index: The index of the internode.
  • Nbl / nb_internodes: The total number of internodes == number of leaves emitted since planting.
  • sh / stem_height: The height of the stem.
  • R / internode_rank_no_expansion: The rank of the internode that will not expand.
  • N / nb_internodes_before_planting: The number of internodes before planting.
  • l_0 / internode_min_height: The minimal length of the internode.

Returns

The length of the internode (m).

Details

The internode length is computed using a quadratic function. The objective is to have a internodes that are short and growing for the first emitted leaves (before nb_internodes_before_planting), and then getting to a stable "constant" height, and at the end for the youngest leaves, having nodes currently growing (smaller).

The internode length is computed as follows : Internode length ^ l | _____________________ | /| || / | | l0 | / | | |-|–-|––––––––––-|–-|––> Internode number 1 N N + N + Nbl Nbl - R where : - l0 is internode_min_height (m), the minimum height of the internode. - l is internode_heigth_final (m), the maximum height of the internode.- N isnbinternodesbeforeplanting, the number of internodes before planting. - R isinternoderanknoexpansion, the number of internodes not in expansion. - Nbl is the number of leaves emitted since planting. with the conditions that : - the sum of the areas of the first triangle, the rectangle and the last triangle is equal tostemheight. - if the equation of the first line isa * x + b: -a = (l - l0) / (N - 1)-b = l0 - a- the area of the first triangle isa * N * (N + 1) / 2 + b * Nand after development :l * N/2 + l0 * N/2- the area of the rectangle (between N + 1 and N + Nbl - R - 1) is(Nbl - R - 1) * l- if the equation of the last line isc * x + d, then: -c = (l0 - l) / R-d = l0 - c * (Nbl + N)- the area of the last triangle is(R + 1) * (c * (2N + 2Nbl - R) / 2 + d)and after development :l * ((R + 1)/ 2) + l_0 * (-(R + 1) / 2 + R + 1)reminder: - the sum of integers from m to n isn * (n + 1) / 2 - m * (m - 1) / 2- the sum of cx + d from m to n isc * (n * (n + 1) / 2 - m * (m - 1) / 2) + d * (n - m + 1)or(n - m + 1) * (c * (n + m) / 2 + d)`

source
XPalm.VPalm.interp_pointsMethod
interp_points(points, step)

Interpolate points along a path to have equidistant points.

Arguments

  • points: Vector of Meshes.Point objects defining the original path.
  • step: Distance between interpolated points.

Returns

  • vec_points: Vector of interpolated Meshes.Point objects.
  • i_discret_pts_exp: Indices of the original points in the interpolated path.
  • vec_dist_p2p1: Vector of distances between consecutive points.
  • vec_angle_xy: Vector of angles between segments and the XY plane.
  • vec_angle_xz: Vector of angles between segments and the XZ plane.
source
XPalm.VPalm.leaf_insertion_angleFunction
leaf_insertion_angle(rank, leaf_max_angle=90, leaf_slope_angle=0.05, leaf_inflection_angle=40)

Compute the insertion angle of the leaf on the internode.

Note: The insertion angle is computed using a logistic function.

Arguments

  • rank: The rank of the leaf.
  • leaf_max_angle: The maximum angle of the leaf.
  • leaf_slope_angle: The slope of the logistic function.
  • leaf_inflection_angle: The inflection point of the logistic function.
source
XPalm.VPalm.leaflet_azimuthal_angleMethod
leaflet_azimuthal_angle(relative_pos, side, angle_c, angle_slope, angle_a, angle_sdp, rng)

Calculate the leaflet insertion angle in the horizontal plane (in degrees).

Arguments

  • relative_pos: Relative position of the leaflet on the rachis [0 to 1].
  • side: Side of the leaf (1 for right, -1 for left).
  • angle_c: Constant parameter for the axial angle calculation (°).
  • angle_slope: Slope parameter for the axial angle calculation (°).
  • angle_a: Amplitude parameter for the axial angle calculation (°).
  • angle_sdp: Standard deviation percentage for random variation (°).
  • rng: Random number generator.

Returns

  • Horizontal insertion angle in degrees.
source
XPalm.VPalm.leaflet_length_at_bpointMethod
leaflet_length_at_bpoint(rachis_length, intercept, slope)

Compute the length of leaflets at the B point of the rachis using a linear relationship.

Arguments

  • rachis_length: The total length of the rachis (m).
  • intercept: The intercept parameter of the linear function (m).
  • slope: The slope parameter of the linear function (dimensionless).

Details

This function uses a linear model to determine leaflet length at the B point:

leaflet_length = intercept + slope * rachis_length

The B point is a key reference point on the rachis that marks the transition from an oval to a round shape of the rachis. The leaflet length at this point serves as a reference for calculating the distribution of leaflet lengths along the entire rachis.

Returns

The length of leaflets at the B point position (m).

source
XPalm.VPalm.leaflet_length_maxMethod
leaflet_length_max(
    leaflet_length_at_b, 
    relative_position_bpoint, 
    relative_length_first_leaflet, 
    relative_length_last_leaflet, 
    relative_position_leaflet_max_length, 
    relative_position_bpoint_sd, 
    rng
)

Calculate the maximum leaflet length for the rachis, used to scale the relative length profile.

Arguments

  • leaflet_length_at_b: Length of leaflets at the B point on the rachis (m).
  • relative_position_bpoint: Relative position of the B point along the rachis (0 to 1).
  • relative_length_first_leaflet: Relative length of the first leaflet at rachis base [0 to 1].
  • relative_length_last_leaflet: Relative length of the last leaflet at rachis tip [0 to 1].
  • relative_position_leaflet_max_length: Relative position where leaflets reach maximum length [0 to 1].
  • relative_position_bpoint_sd: Standard deviation for stochastic variation in B point position.
  • rng: Random number generator.

Details

This function calculates the maximum leaflet length that would result in the specified leaflet length at the B point, considering the shape of the length profile along the rachis.

The calculation uses the inverse of the relative length function at the B point position to determine what maximum value would yield the desired length at that specific position.

Biological Context

In palm fronds, leaflet length typically follows a bell-shaped distribution along the rachis:

  1. Leaflets are short at the base (petiole end)
  2. They increase in length to reach a maximum somewhere close to the middle of the rachis
  3. They decrease in length toward the tip

The B point is a key morphological reference point where the rachis cross-section transitions from oval to round. By knowing the leaflet length at this specific point, we can calculate the maximum leaflet length for the entire frond, which serves as a scaling factor for all other leaflets.

The stochastic variation in B point position reflects natural biological variability between individual palms or fronds.

Returns

The maximum leaflet length for the rachis (m).

source
XPalm.VPalm.leaflet_width_at_bpointMethod
leaflet_width_at_bpoint(rachis_length, intercept, slope)

Calculate leaflet width at B point (reference point).

Arguments

  • rachis_length: The total length of the rachis (m).
  • intercept: The intercept parameter of the linear function (m).
  • slope: The slope parameter of the linear function (dimensionless).

Details

This function uses a linear model to determine leaflet width at the B point:

leaflet_width = intercept + slope * rachis_length

The B point is a key reference point on the rachis that marks the transition between different architectural zones. The leaflet width at this point serves as a reference for calculating the distribution of leaflet widths along the entire rachis.

Returns

  • The width of leaflets at the B point position (m).
source
XPalm.VPalm.leaflet_width_maxMethod
leaflet_width_max(
    leaflet_width_at_b,
    relative_position_bpoint,
    width_first,
    width_last,
    pos_width_max,
    relative_position_bpoint_sd,
    rng
)

leaflet_width_max(
    leaflet_width_at_b,
    relative_position_bpoint,
    width_first,
    width_last,
    pos_width_max,
)

Calculate the maximum leaflet width for the rachis, used to scale the width profile.

Arguments

  • leaflet_width_at_b: Width of leaflets at the B point on the rachis (m).
  • relative_position_bpoint: Mean relative position of the B point along the rachis [0 to 1].
  • width_first: Relative width of the first leaflet at rachis base [0 to 1].
  • width_last: Relative width of the last leaflet at rachis tip [0 to 1].
  • pos_width_max: Relative position where leaflets reach maximum width [0 to 1].
  • relative_position_bpoint_sd: Standard deviation for stochastic variation in B point position (optional).
  • rng: Random number generator (optional).

Details

This function calculates the maximum leaflet width that would result in the specified width at the B point, considering the shape of the width profile along the rachis.

The calculation uses the inverse of the relative width function at the B point position to determine what maximum value would yield the desired width at that specific position.

Biological Context

In palm fronds, leaflet width typically varies along the rachis:

  1. Narrow leaflets at the base (petiole end)
  2. Wider leaflets in the middle region
  3. Narrowing again toward the tip

By knowing the leaflet width at the B point, we can calculate the maximum leaflet width for the entire frond, which serves as a scaling factor for all other leaflets.

Returns

The maximum leaflet width for the rachis (m).

source
XPalm.VPalm.leaflet_zenithal_angleMethod
leaflet_zenithal_angle(relative_pos, leaflet_type, side, high_a0_sup, high_amax_sup, high_a0_inf, high_amax_inf, 
             low_a0_sup, low_amax_sup, low_a0_inf, low_amax_inf, rng)

Calculate the leaflet insertion angle in the vertical plane (in degrees).

Arguments

  • relative_pos: Relative position of the leaflet on the rachis [0 to 1].
  • leaflet_type: Type of leaflet (-1=down, 0=medium, 1=up).
  • side: Side of the leaf (1 for right, -1 for left).
  • high_a0_sup: Upper bound of angle at position 0 for high position leaflets.
  • high_amax_sup: Upper bound of maximum angle for high position leaflets.
  • high_a0_inf: Lower bound of angle at position 0 for high position leaflets.
  • high_amax_inf: Lower bound of maximum angle for high position leaflets.
  • low_a0_sup: Upper bound of angle at position 0 for low position leaflets.
  • low_amax_sup: Upper bound of maximum angle for low position leaflets.
  • low_a0_inf: Lower bound of angle at position 0 for low position leaflets.
  • low_amax_inf: Lower bound of maximum angle for low position leaflets.
  • rng: Random number generator.

Returns

  • Vertical insertion angle in degrees.
source
XPalm.VPalm.leaflet_zenithal_angle_boundariesFunction
leaflet_zenithal_angle_boundaries(rel_pos, a0, a_max, xm=0.5)

Calculate the boundaries of the radial angle based on position along the rachis.

Arguments

  • rel_pos: Relative position on the rachis [0 to 1].
  • a0: Radial angle around C point.
  • a_max: Maximum value of radial angle (in degrees).
  • xm: Relative position on rachis of the maximum radial angle (default: 0.5).

Returns

  • Radial angle in degrees.
source
XPalm.VPalm.leaflets!Method
leaflets!(unique_mtg_id, rachis_node, scale, leaf_rank, rachis_length, parameters; rng=Random.MersenneTwister())

Create leaflets for a given rachis node, computing their positions, types, and dimensions.

Arguments

  • unique_mtg_id: Reference to a unique identifier for the MTG nodes.
  • rachis_node: The parent node of the rachis where leaflets will be attached.
  • scale: The scale of the leaflets in the MTG.
  • leaf_rank: The rank of the leaf associated with the rachis.
  • rachis_length: The total length of the rachis in meters.
  • height_cpoint: The height of the central point of the rachis in meters.
  • width_cpoint: The width of the central point of the rachis in meters.
  • zenithal_cpoint_angle: The zenithal angle of the central point of the rachis in degrees.
  • parameters: A dictionary containing biomechanical parameters for the leaflets.
  • rng: A random number generator for stochastic processes (default is a new MersenneTwister).

Note

The parameters is a Dict{String} containing the following keys:

  • "leaflets_nb_max": Maximum number of leaflets per rachis.
  • "leaflets_nb_min": Minimum number of leaflets per rachis.
  • "leaflets_nb_slope": Slope for the number of leaflets distribution.
  • "leaflets_nb_inflexion": Inflexion point for the number of leaflets distribution.
  • "nbLeaflets_SDP": Standard deviation for the number of leaflets.
  • "leaflet_position_shape_coefficient": Shape coefficient for the relative positions of leaflets.
  • "leaflet_frequency_high": Frequency of high-position leaflets.
  • "leaflet_frequency_low": Frequency of low-position leaflets.
  • "leaflet_frequency_shape_coefficient": Shape coefficient for the frequency distribution of leaflets.
  • "leaflet_between_to_within_group_ratio": Ratio of spacing between groups to within groups.
  • "leaflet_length_at_b_intercept": Intercept for the length of leaflets at point B.
  • "leaflet_length_at_b_slope": Slope for the length of leaflets at point B.
  • "relative_position_bpoint": Relative position of point B along the rachis.
  • "relative_position_bpoint_sd": Standard deviation of the relative position of point B.
  • "relative_length_first_leaflet": Relative length of the first leaflet.
  • "relative_length_last_leaflet": Relative length of the last leaflet.
  • "relative_position_leaflet_max_length": Relative position of the leaflet with maximum length.
  • "leaflet_width_at_b_intercept": Intercept for the width of leaflets at point B.
  • "leaflet_width_at_b_slope": Slope for the width of leaflets at point B.
  • "relative_width_first_leaflet": Relative width of the first leaflet.
  • "relative_width_last_leaflet": Relative width of the last leaflet.
  • "relative_position_leaflet_max_width": Relative position of the leaflet with maximum width.
  • "rachis_nb_segments": Number of segments in the rachis.
source
XPalm.VPalm.linearMethod
linear(x, intercept, slope)

Compute a linear function at given x value.

Arguments

  • x: The input value.
  • intercept: The intercept of the linear function.
  • slope: The slope of the linear function.
source
XPalm.VPalm.local_flexionMethod
local_flexion(current_angle, final_angle, young_modulus, tapering, relative_position)

Calculate the local bending angle at a specific position along the beam.

Arguments

  • current_angle: Current angle in radians
  • final_angle: Final angle of the beam in radians
  • young_modulus: Value of Young's modulus
  • tapering: Tapering factor of the beam
  • relative_position: Relative position along the beam (0 to 1)

Returns

  • Flexion angle at the current position (in radians)
source
XPalm.VPalm.logisticMethod
logistic(x, max, slope, inflection)

Compute a logistic function.

Arguments

  • x: The input value.
  • max: The maximum value of the logistic function.
  • slope: The slope of the logistic function.
  • inflection: The inflection point of the logistic function.
source
XPalm.VPalm.mean_and_sdMethod
mean_and_sd(mean, sd, rng)

Compute a random value from a normal distribution with a given mean and standard deviation.

Arguments

  • mean: The mean of the normal distribution.
  • sd: The standard deviation of the normal distribution.
  • rng: The random number generator.
source
XPalm.VPalm.mtg_skeletonMethod
mtg_skeleton(parameters; rng=Random.MersenneTwister(parameters["seed"]))

Makes an MTG skeleton with nb_leaves_emitted leaves, including all intermediate organs:

  • Plant: the whole palm
  • Stem: the stem of the plant, i.e. the remaining part of the plant after the leaves have been removed
  • Phytomer: the part that includes the leaf and the internode
  • Internodes: the part of the phytomer that is between two leaves
  • Leaf: the leaf of the plant, also called frond

Note: this skeleton does not include reproductive organs (inflorescences, fruits) or the scales that decompose the leaf (petiole, rachis, leaflets).

Arguments

  • parameters: The parameters for the MTG skeleton. See VPalm.default_parameters() for the default parameters.
  • rng: (optional) The random number generator to use for stochastic processes. Defaults to Random.MersenneTwister(parameters["seed"]), but can be set to nothing to disable randomness (for testing).

Examples

file = joinpath(dirname(dirname(pathof(VPalm))), "test", "files", "parameter_file.yml")
parameters = read_parameters(file)
mtg_skeleton(parameters)
source
XPalm.VPalm.normal_deviation_drawMethod
normal_deviation_draw(sd, rng)

Draw a random value from a normal distribution with a given standard deviation.

Arguments

  • sd: The standard deviation of the normal distribution.

Optional arguments

  • rng: The random number generator.
source
XPalm.VPalm.normal_deviation_percent_drawMethod
normal_deviation_percent_draw(value, sdp, rng)

Calculate a normally distributed random deviation based on a percentage of the value.

Arguments

  • value: Base value.
  • sd: Standard deviation in %.
  • rng: Random number generator.

Returns

  • The random deviation.
source
XPalm.VPalm.normalize_positions!Method
normalize_positions!(positions, rachis_length)

Scale and offset positions to span the full rachis length.

Arguments

  • positions: Vector of positions to be modified in place.
  • rachis_length: Total length of rachis in meters.

Details

This function:

  1. Offsets positions so the first leaflet is at position 0
  2. Scales all positions to ensure the last leaflet is exactly at rachis_length
  3. Maintains the relative spacing pattern established by previous processing

This ensures leaflets are properly distributed along the entire rachis while preserving the characteristic grouping patterns.

source
XPalm.VPalm.petioleMethod
petiole(parent_node, index, scale, rachis_length, zenithal_insertion_angle, zenithal_cpoint_angle, parameters)

Make a leaf petiole.

Arguments

  • unique_mtg_id: a next free unique id for the MTG nodes
  • parent_node: the parent node on which the petiole will be attached
  • index: the MTG index of the petiole
  • scale: the MTG scale of the petiole
  • rachis_length: the rachis length, used to feed allometries to compute the petiole dimensions
  • zenithal_insertion_angle: petiole insertion angle
  • zenithal_cpoint_angle: angle at the C point (tip of the petiole, starting point of the rachis)
  • parameters: a list of parameters as a Dict{String}:
    • "leafbasewidth": the base width of the petiole (m)
    • "leafbaseheight": the base heigth of the petiole (m)
    • "cpointwidthintercept": petiole width at the c-point intercept for linear interpolation (m)
    • "cpointwidthslope": petiole width at the c-point slope for linear interpolation
    • "cpointheightwidth_ratio": height to width ratio at the C point
    • "petiolerachisratio_mean": the average value of the ratio between rachis length and petiole length
    • "petiolerachisratio_sd": its standard deviation
    • "petiolenbsegments": the number of segments used to discretize the petiole
source
XPalm.VPalm.petiole_azimuthal_angleMethod
petiole_azimuthal_angle(; rng=Random.MersenneTwister(1))

Compute the azimuthal angle of the petiole based on the petiole/rachis length ratio.

Arguments

  • petiole_rachis_ratio_mean: Average value of the petiole/rachis length ratio
  • petiole_rachis_ratio_sd: Standard deviation of the petiole/rachis length ratio
  • rng: Random number generator

Returns

The azimuthal angle of the petiole (°)

source
XPalm.VPalm.petiole_dimensions_at_cpointFunction
petiole_dimensions_at_cpoint(rachis_length, cpoint_width_intercept, cpoint_width_slope, cpoint_height_width_ratio)

Compute the width and height of the petiole at the C point (end-point).

Arguments

  • rachis_length: Length of the rachis (m)
  • cpoint_width_intercept=0.0098u"m": Intercept of the linear relationship between rachis width at C point and rachis length (m)
  • cpoint_width_slope=0.012: Slope of the linear relationship
  • cpoint_height_width_ratio=0.568: Ratio between the height and width of the leaf at C point

Returns

A named tuple with the following keys:

  • width_cpoint: Width at the C point of the petiole (m)
  • height_cpoint: Height at the C point of the petiole (m)
source
XPalm.VPalm.petiole_heightMethod
petiole_height(relative_position, height_cpoint, height_base)

Compute the height profile along the petiole (m).

Arguments

  • relative_position: Position along the petiole (0-1)
  • height_base: Height at the base of the leaf
  • height_cpoint: Height of the leaf section at C point
source
XPalm.VPalm.petiole_lengthFunction
petiole_length(petiole_rachis_ratio_mean, petiole_rachis_ratio_sd, rachis_length; rng=Random.MersenneTwister(1))

Compute the length of the petiole based on the rachis length and the petiole/rachis length ratio.

Arguments

  • rachis_length: Length of the rachis (m)
  • petiole_rachis_ratio_mean=0.25: Average value of the petiole/rachis length ratio
  • petiole_rachis_ratio_sd=0.034: Standard deviation of the petiole/rachis length ratio
  • rng: Random number generator

Returns

The length of the petiole (m)

source
XPalm.VPalm.petiole_sections!Method
petiole_sections!(petiole_node, petiole_nb_segments, unique_mtg_id)

Create the sections of a petiole.

Arguments

  • petiole_node: the MTG Node of the petiole
  • petiole_nb_segments: the number of segments used to discretize the petiole
  • unique_mtg_id: a next free unique id for the MTG nodes, given as a Ref{Int}

Returns

Nothing, the petiole node is updated in-place with its sections.

source
XPalm.VPalm.petiole_widthMethod
petiole_width(relative_position, width_cpoint, width_base)

Compute the width profile along the petiole (m).

Arguments

  • relative_position: Position along the petiole (0-1)
  • width_base: Width at base of leaf
  • width_cpoint: Width of the leaf at C point
source
XPalm.VPalm.phyllotactic_angleMethod
phyllotactic_angle(phyllotactic_angle_mean, phyllotactic_angle_sd; rng)

Computes the phyllotactic angle (°) using an average angle and a standard deviation (random draw from a normal distribution).

Arguments

  • phyllotactic_angle_mean: The average phyllotactic angle (°).
  • phyllotactic_angle_sd: The standard deviation of the phyllotactic angle (°).

Optional arguments

  • rng: The random number generator.
source
XPalm.VPalm.piecewise_linear_areaMethod
piecewise_linear_area(x, y)

Calculate the area under a piecewise linear function. Equivalent to PiecewiseFunctionArea in the Java version.

Arguments

  • x: Array of x-coordinates of the control points.
  • y: Array of y-coordinates of the control points.

Returns

  • Area under the piecewise linear function.
source
XPalm.VPalm.planeMethod
create_plane_mesh()

Create a simple rectangular plane mesh that will be used as a reference for leaflet segments. The plane is created in the XZ plane with width along X and length along Z.

Returns

A Meshes.SimpleMesh object representing a simple rectangular plane mesh

source
XPalm.VPalm.properties_petiole_sectionMethod
properties_petiole_section(
    index, nb_sections, width_base, height_base,
    width_cpoint, height_cpoint, petiole_section_length,
    petiole_insertion_angle, petiole_section_insertion_angle,
    azimuthal_angle
)

Compute the properties of each section of the petiole.

Arguments

  • index: The index of the section within all sections (1-nb_sections)
  • nb_sections: The number of sections discretizing the petiole
  • width_base: Width of the petiole at its base (m)
  • heigth_base: Height of the petiole at its base (m)
  • width_cpoint: Width of the petiole at the C point (tip of the petiole, i.e. transition point to rachis, m)
  • height_cpoint: Height at the C point (m)
  • petiole_section_length: The length of the petiole sections (m)
  • petiole_insertion_angle: Zenithal angle of insertion between the petiole and the stipe (local angle, relative to the stipe, °)
  • petiole_section_insertion_angle: The zenithal angle of insertion between the petioles sections (°)
  • azimuthal_angle: Azimuthal angle at the insertion (°)

Returns

A vector of dimensions for each section, given as a named tuple:

  • width: width of the section (m)
  • height: height of the section (m)
  • length: length of the section (m)
  • zenithal_angle: zenithal angle of the section (global angle, °)
  • azimuthal_angle: azimuthal angle of the section (global angle, °)
  • torsion_angle: torsion angle of the section (°)
source
XPalm.VPalm.rachisMethod
rachis(unique_mtg_id, index, scale, leaf_rank, rachis_length, height_cpoint, width_cpoint, zenithal_cpoint_angle, fresh_biomass, parameters; rng)

Builds a rachis node in the MTG structure.

Arguments

  • unique_mtg_id: A reference to a unique identifier for the MTG nodes.
  • index: The index of the rachis segment.
  • scale: The scale of the rachis segment.
  • leaf_rank: The rank of the leaf associated with the rachis.
  • rachis_length: The length of the rachis (m).
  • height_cpoint: The height of the central point of the rachis (m).
  • width_cpoint: The width of the central point of the rachis (m).
  • zenithal_cpoint_angle: The zenithal angle of the central point of the rachis (°).
  • fresh_biomass: The fresh biomass of the rachis (kg).
  • parameters: A dictionary containing biomechanical parameters for the rachis (see below).
  • rng: A random number generator for stochastic processes.

Note

The parameters is a Dict{String} containing the following keys:

  • "rachis_nb_segments": The number of segments in the rachis.
  • "rachis_twist_initial_angle": The initial twist angle of the rachis (°).
  • "rachis_twist_initial_angle_sdp": The standard deviation of the initial twist angle (°).
  • "elastic_modulus": The elastic modulus of the rachis material (Pa).
  • "shear_modulus": The shear modulus of the rachis material (Pa).
  • "leaflet_length_at_b_intercept": The intercept for the leaflet length at point B (m).
  • "leaflet_length_at_b_slope": The slope for the leaflet length at point B (m).
  • "relative_position_bpoint": The relative position of point B along the rachis (m).
  • "relative_position_bpoint_sd": The standard deviation of the relative position of point B (m).
  • "relative_length_first_leaflet": The relative length of the first leaflet (m).
  • "relative_length_last_leaflet": The relative length of the last leaflet (m).
  • "relative_position_leaflet_max_length": The relative position of the leaflet with maximum length (m).
  • "height_rachis_tappering": The height at which the rachis tapers (m).
  • "rachis_width_tip": The width of the rachis tip (m).
source
XPalm.VPalm.rachis_expansionMethod
rachis_expansion(leaf_rank, rachis_final_length)

Simple function to compute the rachis expansion (using an expansion factor)
    based on the leaf rank.

# Arguments

- `leaf_rank`: The rank of the leaf.
- `rachis_final_length`: The final length of the rachis.
source
XPalm.VPalm.rachis_heightMethod
rachis_height(relative_position, cpoint_height, rachis_height_tappering)

Computes the rachis height (m) at a given relative position using a the height at C Point and rachis tappering.

Arguments

  • relative_position: The relative position along the rachis (0: base to 1: tip).
  • cpoint_height: The height of the rachis at the C point, i.e. rachis base (m).
  • rachis_height_tappering: The tappering factor for the rachis height.
source
XPalm.VPalm.rachis_length_from_biomassMethod
rachis_length_from_biomass(rachis_biomass, leaf_length_intercept, leaf_length_slope)

Compute the length of the rachis based on its biomass using a linear relationship.

Arguments

  • rachis_biomass: The biomass of the rachis (g).
  • leaf_length_intercept: The intercept of the linear relationship for leaf length.
  • leaf_length_slope: The slope of the linear relationship for leaf length.

Returns

The length of the rachis (m).

source
XPalm.VPalm.rachis_widthMethod
rachis_width(relative_position, cpoint_width, rachis_width_tip)

Computes the rachis width (m) at a given relative position using the width at C Point and rachis width at the tip.

Arguments

  • relative_position: The relative position along the rachis (0: base to 1: tip).
  • cpoint_width: The width of the rachis at the C point, i.e. rachis base (m).
  • rachis_width_tip: The width of the rachis at the tip (m).
source
XPalm.VPalm.read_parametersMethod
read_parameters(file; verbose=true)

Reads a parameter file and returns the contents as an ordered dictionary.

Arguments

  • file: The path to the parameter file.
  • verbose: Whether to show warnings for units (default: true)

Returns

An ordered dictionary containing the contents of the parameter file with appropriate units.

Example

file = joinpath(dirname(dirname(pathof(VPalm))),"test","files","parameter_file.yml")
read_parameters(file)
source
XPalm.VPalm.read_plyMethod
read_ply(fname)

Reads a PLY file and returns a Meshes.SimpleMesh object.

Arguments

  • fname: The path to the PLY file.

Returns

A Meshes.SimpleMesh object.

source
XPalm.VPalm.relative_leaflet_lengthMethod
relative_leaflet_length(x, relative_length_first_leaflet, relative_length_last_leaflet, relative_position_leaflet_max_length)

Relative leaflet length given by their relative position along the rachis.

Arguments

  • x: relative leaflet position on the rachis (0: base to 1: tip)
  • relative_length_first_leaflet: relative length of the first leaflet on the rachis (0 to 1)
  • relative_length_last_leaflet: relative length of the last leaflet on the rachis (0 to 1)
  • relative_position_leaflet_max_length: relative position of the longest leaflet on the rachis (0.111 to 0.999)
source
XPalm.VPalm.relative_leaflet_positionMethod
relative_leaflet_position(relative_rank, shape_coefficient)

Compute the relative leaflet position on the rachis.

Arguments

  • relative_rank: Relative leaflet rank, usually in the form of (0 to 1].
  • shape_coefficient: Shape coefficient (around 0).

Returns

The relative leaflet position, in the same form as relative_rank, usually (0 to 1].

source
XPalm.VPalm.relative_leaflet_widthMethod
relative_leaflet_width(x, width_first, width_last, pos_width_max)

Calculate the relative leaflet width at a given position along the rachis.

Arguments

  • x: Relative position of the leaflet on the rachis [0 to 1].
  • width_first: Relative width of the first leaflet (at rachis base).
  • width_last: Relative width of the last leaflet (at rachis tip).
  • pos_width_max: Relative position where leaflets reach maximum width [0 to 1].

Details

This function uses a piecewise linear model to calculate relative leaflet width:

  1. From base to maximum width position: Linear increase from width_first to 1.0
  2. From maximum width position to tip: Linear decrease from 1.0 to width_last

Biological Context

The width of leaflets along a palm frond typically follows a pattern where:

  • Leaflets start relatively narrow at the base
  • Widen to reach maximum width at some point along the rachis
  • Narrow again toward the tip

Returns

The relative width at position x [0 to 1].

source
XPalm.VPalm.shrink_leaflets_in_groups!Function
shrink_leaflets_in_groups!(positions, leaflets, ratio=2.0)

Adjust the spacing between leaflets to create appropriate within-group and between-group distances.

Arguments

  • positions: Vector of current leaflet positions along the rachis.
  • leaflets: A NamedTuple containing arrays for leaflet properties (group, group_size, plane).
  • ratio=2.0: Ratio of inter-group to intra-group spacing.

Details

This function implements a biological principle where leaflets within the same group are positioned closer together than leaflets in different groups. It:

  1. Uses a fixed ratio (2:1) between inter-group and intra-group spacing
  2. Preserves the overall distribution pattern while creating distinct groups
  3. Processes each group sequentially, adjusting positions based on group size

Biological Context

In many palm species, particularly oil palm, leaflets appear in distinct groups along the rachis. This grouping pattern is characterized by:

  1. Consistent, smaller spacing between leaflets within the same group
  2. Larger spacing between adjacent groups
  3. The ratio between these spacings is typically species-specific

This spacing pattern is essential for the palm's characteristic appearance and affects light interception patterns along the frond.

source
XPalm.VPalm.snagMethod
snag()
snag(l, w, h)

Returns a normalized snag mesh, or a snag mesh with given dimensions in m.

source
XPalm.VPalm.stem_bendingMethod
stem_bending(stem_bending_mean, stem_bending_sd; rng)

Computes the stem bending (°) using an average bending and a standard deviation (random draw from a normal distribution).

Arguments

  • stem_bending_mean: The average stem bending (°).
  • stem_bending_sd: The standard deviation of the stem bending (°).

Optional arguments

  • rng: The random number generator.
source
XPalm.VPalm.stem_diameterMethod
stem_diameter(rachis_length_reference, stem_diameter_max, stem_diameter_slope, stem_diameter_inflection, stem_diameter_residual; rng)

Computes the stem diameter (m) at a given rachis length reference (m).

Arguments

  • rachis_length_reference: The rachis length reference (m). Taken as the rachis length of the first leaf.
  • stem_diameter_max: The maximum stem diameter (m).
  • stem_diameter_slope: The slope of the logistic function.
  • stem_diameter_inflection: The inflection point of the logistic function.
  • stem_diameter_residual: The residual of the stem diameter (m).
  • stem_diameter_snag: The diameter estimation due to snags (m).
  • rng: The random number generator.

Details

The stem diameter is computed using a logistic function, and then some variability is added to simulate natural variations that might occur in real-world scenarios.

source
XPalm.VPalm.stem_heightMethod
stem_height(nb_leaves_emitted, initial_stem_height, stem_height_coefficient, internode_length_at_maturity, stem_growth_start; rng)

Computes the stem height (m) at a given number of leaves emitted.

Arguments

  • nb_leaves_emitted: The number of leaves emitted from planting.
  • initial_stem_height: The initial stem height at planting (m).
  • stem_height_coefficient: The coefficient of the exponential function.
  • internode_length_at_maturity: The internode length when the plant is mature (m).
  • stem_growth_start: The number of leaves emitted at which the stem starts to grow (m). This is because the stem does not grow at the same rate at the beginning of the plant's life,

because it first grows more in diameter than in height.

  • stem_height_variation: The variation of the stem height (m) due to the random draw from a normal distribution.
  • rng: The random number generator.

Details

The stem height is computed using an exponential function for the first stem_growth_start leaves emitted, and then a linear function for the remaining leaves emitted.

Note that the stem height can also be subject to some variability using stem_height_variation, simulating natural variations that might occur in real-world scenarios, but this variability will never make the stem height go below 30% of the intial computed height.

source
XPalm.VPalm.unbendMethod
unbend(distance, inclination)

Removes torsion and bending of a bent beam by transforming it into a straight line, while preserving its insertion angle (inclination angle of the first segment).

Arguments

  • distance: Vector of distances between consecutive points (in meters)
  • inclination: Vector of inclination angles (in degrees), only the first value is used

Returns

  • A vector of Meshes.Point objects representing the unbent positions

Details

This function creates a straight line with the same cumulative length as the input distances, while maintaining the insertion angle specified by the first inclination value. The output points represent the unbent state of a curved structure.

Note

Mainly used to compute the input coordinates for bend() from experimental points.

Example

using VPalm, Unitful, Meshes
distances = [0.001, 1.340, 1.340, 0.770, 0.770]u"m";
inclinations = [48.8, 48.8, 48.8, 48.8, 48.8];  # degrees
points = VPalm.unbend(distances, inclinations)
source
XPalm.VPalm.update_leaflet_angles!Method
update_leaflet_angles!(
    leaflet, leaf_rank; 
    last_rank_unfolding=2, unique_mtg_id=new_id(leaflet), 
    xm_intercept=0.176, xm_slope=0.08, 
    ym_intercept=0.51, ym_slope=-0.025
)

Update the angles and stiffness of a leaflet based on its position, side, and leaf rank.

Arguments

  • leaflet: The leaflet node to update
  • leaf_rank: The rank of the leaf (affects unfolding for young leaves)
  • last_rank_unfolding: Rank at which leaflets are fully unfolded (default is 2)
  • unique_mtg_id: Reference to the unique ID counter for MTG nodes (default is the maximum ID in the MTG)
  • xm_intercept, xm_slope: Parameters for defining maximum leaflet width position
  • ym_intercept, ym_slope: Parameters for defining maximum leaflet width value
source
XPalm.VPalm.update_petiole_angles!Method
update_petiole_angles!(petiole_node)

Update the angles of the petiole segments based on the petiole node properties.

Arguments

  • petiole_node: the MTG Node of the petiole

Returns

Nothing, the petiole node is updated in-place with its segments angles.

source
XPalm.VPalm.update_rachis_angles!Method
update_rachis_angles!(rachis_node, leaf_rank, rachis_length, height_cpoint, width_cpoint, zenithal_cpoint_angle, fresh_biomass, parameters; rng)

Update the angles and dimensions of the rachis segments based on biomechanical properties. This function is called when the rachis exists already, but needs updating its angles and dimensions following a change in the leaf rank, length, or fresh biomass.

source
XPalm.VPalm.update_segment_angles!Method
update_segment_angles!(leaflet, young_modulus, initial_angle, leaflet_length, tapering)

Update the zenithal angles of each segment in a leaflet based on the Young's modulus model.

Arguments

  • leaflet: The leaflet MTG node containing segments
  • young_modulus: Value of Young's modulus
  • initial_angle: Initial angle from vertical in radians
  • leaflet_length: Total length of the leaflet
  • tapering: Tapering factor

Returns

  • Nothing (the angles are updated in place)
source
XPalm.VPalm.width_at_cpointMethod
width_at_cpoint(rachis_length, cpoint_width_intercept, cpoint_width_slope)

Compute width at C point based on rachis length.

Arguments

  • rachis_length: Length of rachis (m)
  • cpoint_width_intercept: Intercept of linear function (m)
  • cpoint_width_slope: Slope of linear function
source
XPalm.VPalm.write_parametersMethod
write_parameters(file, params)

Write the given parameters to a file using YAML format.

Arguments

  • file: The file path to write the parameters to.
  • params: The parameters to be written.

Example

file = joinpath(dirname(dirname(pathof(VPalm))),"test","files","parameter_file.yml")
params = read_parameters(file)
write_parameters(tempname(), params)
source
XPalm.VPalm.xyz_to_dist_and_anglesMethod
xyz_to_dist_and_angles(points)

Compute segment lengths and angles from point coordinates.

Arguments

  • points: A vector of Meshes.Points.

Returns

  • A NamedTuple with fields:
    • dist_p2p1: Vector of segment lengths (m).
    • vangle_xy: Vector of angles between the segment and the XY plane (radians).
    • vangle_xz: Vector of angles between the segment and the XZ plane (radians).
source
XPalm.VPalm.@check_unitMacro
@check_unit(variable, expected_unit, [verbose=true])

Check if a variable has the expected unit type. If no unit is found, assign the expected unit. If a different unit is found, try to convert to the expected unit.

Arguments

  • variable: The variable to check
  • expected_unit: The expected unit (e.g., u"m")
  • verbose: Whether to show warnings (default: true)

Examples

```julia rachislength = 10 @checkunit rachis_length u"m" # Will add u"m" and warn

petiolelength = 10 @checkunit petiole_length u"m" false # Will add u"m" without a warning

mass = 5.0u"g" @check_unit mass u"kg" # Will convert g to kg

source