HighMap library (C++)
Loading...
Searching...
No Matches
hmap::gpu Namespace Reference

Functions

Array blend_gradients (const Array &array1, const Array &array2, int ir=4)
 See hmap::blend_gradients.
 
Array blend_poisson_bf (const Array &array1, const Array &array2, const int iterations=500, const Array *p_mask=nullptr)
 Blends two arrays using Poisson blending with a brute-force solver.
 
Array accumulation_curvature (const Array &z, int ir)
 See hmap::accumulation_curvature.
 
Array curvature_horizontal_cross_sectional (const Array &z, int ir)
 See hmap::curvature_horizontal_cross_sectional.
 
Array curvature_horizontal_plan (const Array &z, int ir)
 See hmap::curvature_horizontal_plan.
 
Array curvature_horizontal_tangential (const Array &z, int ir)
 See hmap::curvature_horizontal_tangential.
 
Array curvature_ring (const Array &z, int ir)
 See hmap::curvature_ring.
 
Array curvature_rotor (const Array &z, int ir)
 See hmap::curvature_rotor.
 
Array curvature_vertical_longitudinal (const Array &z, int ir)
 See hmap::curvature_vertical_longitudinal.
 
Array curvature_vertical_profile (const Array &z, int ir)
 See hmap::curvature_vertical_profile.
 
Array shape_index (const Array &z, int ir)
 See hmap::shape_index.
 
Array unsphericity (const Array &z, int ir)
 See hmap::unsphericity.
 
void hydraulic_particle (Array &z, int nparticles, int seed, Array *p_bedrock=nullptr, Array *p_moisture_map=nullptr, Array *p_erosion_map=nullptr, Array *p_deposition_map=nullptr, float c_capacity=10.f, float c_erosion=0.05f, float c_deposition=0.05f, float c_inertia=0.3f, float drag_rate=0.001f, float evap_rate=0.001f, bool post_filtering=false)
 See hmap::hydraulic_particle.
 
void hydraulic_particle (Array &z, Array *p_mask, int nparticles, int seed, Array *p_bedrock=nullptr, Array *p_moisture_map=nullptr, Array *p_erosion_map=nullptr, Array *p_deposition_map=nullptr, float c_capacity=10.f, float c_erosion=0.05f, float c_deposition=0.05f, float c_inertia=0.3f, float drag_rate=0.001f, float evap_rate=0.001f, bool post_filtering=false)
 See hmap::hydraulic_particle.
 
void hydraulic_schott (Array &z, int iterations, const Array &talus, float c_erosion=1.f, float c_thermal=0.1f, float c_deposition=0.2f, float flow_acc_exponent=0.8f, float flow_acc_exponent_depo=0.8f, float flow_routing_exponent=1.3f, float thermal_weight=1.5f, float deposition_weight=2.5f, Array *p_flow=nullptr)
 Simulates hydraulic erosion and deposition on a heightmap using the Schott method.
 
void hydraulic_schott (Array &z, int iterations, const Array &talus, Array *p_mask, float c_erosion=1.f, float c_thermal=0.1f, float c_deposition=0.2f, float flow_acc_exponent=0.8f, float flow_acc_exponent_depo=0.8f, float flow_routing_exponent=1.3f, float thermal_weight=1.5f, float deposition_weight=2.5f, Array *p_flow=nullptr)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
void hydraulic_stream_log (Array &z, float c_erosion, float talus_ref, int deposition_ir=32, float deposition_scale_ratio=1.f, float gradient_power=0.8f, float gradient_scaling_ratio=1.f, int gradient_prefilter_ir=16, float saturation_ratio=1.f, Array *p_bedrock=nullptr, Array *p_moisture_map=nullptr, Array *p_erosion_map=nullptr, Array *p_deposition_map=nullptr, Array *p_flow_map=nullptr)
 See hmap::hydraulic_stream_log.
 
void hydraulic_stream_log (Array &z, float c_erosion, float talus_ref, Array *p_mask, int deposition_ir=32, float deposition_scale_ratio=1.f, float gradient_power=0.8f, float gradient_scaling_ratio=1.f, int gradient_prefilter_ir=16, float saturation_ratio=1.f, Array *p_bedrock=nullptr, Array *p_moisture_map=nullptr, Array *p_erosion_map=nullptr, Array *p_deposition_map=nullptr, Array *p_flow_map=nullptr)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
void thermal (Array &z, const Array &talus, int iterations=10, Array *p_bedrock=nullptr, Array *p_deposition_map=nullptr)
 See hmap::thermal.
 
void thermal (Array &z, Array *p_mask, const Array &talus, int iterations=10, Array *p_bedrock=nullptr, Array *p_deposition_map=nullptr)
 See hmap::thermal.
 
void thermal (Array &z, float talus, int iterations=10, Array *p_bedrock=nullptr, Array *p_deposition_map=nullptr)
 See hmap::thermal.
 
void thermal_auto_bedrock (Array &z, const Array &talus, int iterations=10, Array *p_deposition_map=nullptr)
 See hmap::thermal_auto_bedrock.
 
void thermal_auto_bedrock (Array &z, Array *p_mask, const Array &talus, int iterations=10, Array *p_deposition_map=nullptr)
 See hmap::thermal_auto_bedrock.
 
void thermal_auto_bedrock (Array &z, float, int iterations=10, Array *p_deposition_map=nullptr)
 See hmap::thermal_auto_bedrock.
 
void thermal_inflate (Array &z, const Array &talus, int iterations=10)
 Apply thermal weathering erosion to give a scree like effect.
 
void thermal_inflate (Array &z, const Array *p_mask, const Array &talus, int iterations=10)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
void thermal_rib (Array &z, int iterations, Array *p_bedrock=nullptr)
 See hmap::thermal_rib.
 
void thermal_ridge (Array &z, const Array &talus, int iterations=10, Array *p_deposition_map=nullptr)
 Apply thermal weathering erosion to give a ridge like effect.
 
void thermal_ridge (Array &z, const Array *p_mask, const Array &talus, int iterations=10, Array *p_deposition_map=nullptr)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
void thermal_scree (Array &z, const Array &talus, const Array &zmax, int iterations=10, bool talus_constraint=true, Array *p_deposition_map=nullptr)
 Performs thermal scree erosion on a heightmap.
 
void thermal_scree (Array &z, const Array *p_mask, const Array &talus, const Array &zmax, int iterations=10, bool talus_constraint=true, Array *p_deposition_map=nullptr)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
Array local_median_deviation (const Array &array, int ir)
 See hmap::local_median_deviation.
 
Array mean_local (const Array &array, int ir)
 See hmap::mean_local.
 
Array relative_elevation (const Array &array, int ir)
 See hmap::relative_elevation.
 
Array ruggedness (const Array &array, int ir)
 See hmap::ruggedness.
 
Array rugosity (const Array &z, int ir, bool convex=true)
 See hmap::rugosity.
 
Array std_local (const Array &array, int ir)
 See hmap::std_local.
 
Array z_score (const Array &array, int ir)
 See hmap::z_score.
 
void expand (Array &array, int ir, int iterations=1)
 See hmap::expand.
 
void expand (Array &array, int ir, const Array *p_mask, int iterations=1)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
void expand (Array &array, const Array &kernel, int iterations=1)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
void expand (Array &array, const Array &kernel, const Array *p_mask, int iterations=1)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
void gamma_correction_local (Array &array, float gamma, int ir, float k=0.1f)
 See hmap::gamma_correction_local.
 
void gamma_correction_local (Array &array, float gamma, int ir, const Array *p_mask, float k=0.1f)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
void laplace (Array &array, float sigma=0.2f, int iterations=3)
 See hmap::laplace.
 
void laplace (Array &array, const Array *p_mask, float sigma=0.2f, int iterations=3)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
Array maximum_local (const Array &array, int ir)
 See hmap::maximum_local.
 
Array maximum_local_disk (const Array &array, int ir)
 See hmap::maximum_local_disk.
 
Array mean_shift (const Array &array, int ir, float talus, int iterations=1, bool talus_weighted=true)
 See hmap::mean_shift.
 
Array mean_shift (const Array &array, int ir, float talus, const Array *p_mask, int iterations=1, bool talus_weighted=true)
 
void median_3x3 (Array &array)
 See hmap::median_3x3.
 
void median_3x3 (Array &array, const Array *p_mask)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
Array median_pseudo (const Array &array, int ir)
 See hmap::median_pseudo.
 
Array minimum_local (const Array &array, int ir)
 See hmap::minimum_local.
 
Array minimum_local_disk (const Array &array, int ir)
 See hmap::minimum_local_disk.
 
void normal_displacement (Array &array, float amount=0.1f, int ir=0, bool reverse=false)
 See hmap::normal_displacement.
 
void normal_displacement (Array &array, const Array *p_mask, float amount=0.1f, int ir=0, bool reverse=false)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
void plateau (Array &array, const Array *p_mask, int ir, float factor)
 See hmap::plateau.
 
void plateau (Array &array, int ir, float factor)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
void shrink (Array &array, int ir, int iterations=1)
 See hmap::shrink.
 
void shrink (Array &array, int ir, const Array *p_mask, int iterations=1)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
void shrink (Array &array, const Array &kernel, int iterations=1)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
void shrink (Array &array, const Array &kernel, const Array *p_mask, int iterations=1)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
void smooth_cpulse (Array &array, int ir)
 See hmap::smooth_cpulse.
 
void smooth_cpulse (Array &array, int ir, const Array *p_mask)
 See hmap::smooth_cpulse.
 
void smooth_fill (Array &array, int ir, float k=0.1f, Array *p_deposition_map=nullptr)
 See hmap::smooth_fill.
 
void smooth_fill (Array &array, int ir, const Array *p_mask, float k=0.1f, Array *p_deposition_map=nullptr)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
void smooth_fill_holes (Array &array, int ir)
 See hmap::smooth_fill_holes.
 
void smooth_fill_holes (Array &array, int ir, const Array *p_mask)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
void smooth_fill_smear_peaks (Array &array, int ir)
 See hmap::smooth_fill_smear_peaks.
 
void smooth_fill_smear_peaks (Array &array, int ir, const Array *p_mask)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
Array gradient_norm (const Array &array)
 See hmap::gradient_norm.
 
Array flow_direction_d8 (const Array &z)
 See hmap::flow_direction_d8.
 
Array generate_riverbed (const Path &path, Vec2< int > shape, Vec4< float > bbox={0.f, 1.f, 0.f, 1.f}, bool bezier_smoothing=false, float depth_start=0.01f, float depth_end=1.f, float slope_start=64.f, float slope_end=32.f, float shape_exponent_start=1.f, float shape_exponent_end=10.f, float k_smoothing=0.5f, int post_filter_ir=0, Array *p_noise_x=nullptr, Array *p_noise_y=nullptr, Array *p_noise_r=nullptr)
 See hmap::generate_riverbed.
 
void interpolate_array_bicubic (const Array &source, Array &target)
 
void interpolate_array_bicubic (const Array &source, Array &target, const Vec4< float > &bbox_source, const Vec4< float > &bbox_target)
 
void interpolate_array_bilinear (const Array &source, Array &target)
 
void interpolate_array_bilinear (const Array &source, Array &target, const Vec4< float > &bbox_source, const Vec4< float > &bbox_target)
 
void interpolate_array_lagrange (const Array &source, Array &target, int order)
 
void interpolate_array_nearest (const Array &source, Array &target)
 
void interpolate_array_nearest (const Array &source, Array &target, const Vec4< float > &bbox_source, const Vec4< float > &bbox_target)
 
Array border (const Array &array, int ir)
 See hmap::border.
 
Array closing (const Array &array, int ir)
 See hmap::closing.
 
Array dilation (const Array &array, int ir)
 See hmap::dilation.
 
Array erosion (const Array &array, int ir)
 See hmap::erosion.
 
Array morphological_black_hat (const Array &array, int ir)
 See hmap::morphological_black_hat.
 
Array morphological_gradient (const Array &array, int ir)
 See hmap::morphological_gradient.
 
Array morphological_top_hat (const Array &array, int ir)
 See hmap::morphological_top_hat.
 
Array opening (const Array &array, int ir)
 See hmap::opening.
 
Array relative_distance_from_skeleton (const Array &array, int ir_search, bool zero_at_borders=true, int ir_erosion=1)
 See hmap::relative_distance_from_skeleton.
 
Array skeleton (const Array &array, bool zero_at_borders=true)
 See hmap::skeleton.
 
void helper_bind_optional_buffer (clwrapper::Run &run, const std::string &id, const Array *p_array)
 
bool init_opencl ()
 
Array basalt_field (Vec2< int > shape, Vec2< float > kw, uint seed, float warp_kw=4.f, float large_scale_warp_amp=0.2f, float large_scale_gain=6.f, float large_scale_amp=0.2f, float medium_scale_kw_ratio=3.f, float medium_scale_warp_amp=1.f, float medium_scale_gain=7.f, float medium_scale_amp=0.08f, float small_scale_kw_ratio=10.f, float small_scale_amp=0.1f, float small_scale_overlay_amp=0.002f, float rugosity_kw_ratio=1.f, float rugosity_amp=1.f, bool flatten_activate=true, float flatten_kw_ratio=1.f, float flatten_amp=0.f, const Array *p_noise_x=nullptr, const Array *p_noise_y=nullptr, Vec4< float > bbox={0.f, 1.f, 0.f, 1.f})
 Generates a synthetic procedural terrain resembling basaltic landforms.
 
Array gabor_wave (Vec2< int > shape, Vec2< float > kw, uint seed, const Array &angle, float angle_spread_ratio=1.f, Vec4< float > bbox={0.f, 1.f, 0.f, 1.f})
 Return an array filled with coherence Gabor noise.
 
Array gabor_wave (Vec2< int > shape, Vec2< float > kw, uint seed, float angle=0.f, float angle_spread_ratio=1.f, Vec4< float > bbox={0.f, 1.f, 0.f, 1.f})
 
Array gabor_wave_fbm (Vec2< int > shape, Vec2< float > kw, uint seed, const Array &angle, float angle_spread_ratio=1.f, int octaves=8, float weight=0.7f, float persistence=0.5f, float lacunarity=2.f, const Array *p_ctrl_param=nullptr, const Array *p_noise_x=nullptr, const Array *p_noise_y=nullptr, Vec4< float > bbox={0.f, 1.f, 0.f, 1.f})
 Return an array filled with coherence Gabor noise.
 
Array gabor_wave_fbm (Vec2< int > shape, Vec2< float > kw, uint seed, float angle=0.f, float angle_spread_ratio=1.f, int octaves=8, float weight=0.7f, float persistence=0.5f, float lacunarity=2.f, const Array *p_ctrl_param=nullptr, const Array *p_noise_x=nullptr, const Array *p_noise_y=nullptr, Vec4< float > bbox={0.f, 1.f, 0.f, 1.f})
 
Array gavoronoise (Vec2< int > shape, Vec2< float > kw, uint seed, const Array &angle, float amplitude=0.05f, float angle_spread_ratio=1.f, Vec2< float > kw_multiplier={4.f, 4.f}, float slope_strength=1.f, float branch_strength=2.f, float z_cut_min=0.2f, float z_cut_max=1.f, int octaves=8, float persistence=0.4f, float lacunarity=2.f, const Array *p_ctrl_param=nullptr, const Array *p_noise_x=nullptr, const Array *p_noise_y=nullptr, Vec4< float > bbox={0.f, 1.f, 0.f, 1.f})
 Generates a 2D array using the GavoroNoise algorithm, which is a procedural noise technique for terrain generation and other applications.
 
Array gavoronoise (Vec2< int > shape, Vec2< float > kw, uint seed, float angle=0.f, float amplitude=0.05f, float angle_spread_ratio=1.f, Vec2< float > kw_multiplier={4.f, 4.f}, float slope_strength=1.f, float branch_strength=2.f, float z_cut_min=0.2f, float z_cut_max=1.f, int octaves=8, float persistence=0.4f, float lacunarity=2.f, const Array *p_ctrl_param=nullptr, const Array *p_noise_x=nullptr, const Array *p_noise_y=nullptr, Vec4< float > bbox={0.f, 1.f, 0.f, 1.f})
 
Array gavoronoise (const Array &base, Vec2< float > kw, uint seed, float amplitude=0.05f, Vec2< float > kw_multiplier={4.f, 4.f}, float slope_strength=1.f, float branch_strength=2.f, float z_cut_min=0.2f, float z_cut_max=1.f, int octaves=8, float persistence=0.4f, float lacunarity=2.f, const Array *p_ctrl_param=nullptr, const Array *p_noise_x=nullptr, const Array *p_noise_y=nullptr, Vec4< float > bbox={0.f, 1.f, 0.f, 1.f})
 
Array mountain_range_radial (Vec2< int > shape, Vec2< float > kw, uint seed, float half_width=0.2f, float angle_spread_ratio=0.5f, float core_size_ratio=1.f, Vec2< float > center={0.5f, 0.5f}, int octaves=8, float weight=0.7f, float persistence=0.5f, float lacunarity=2.f, const Array *p_ctrl_param=nullptr, const Array *p_noise_x=nullptr, const Array *p_noise_y=nullptr, const Array *p_angle=nullptr, Vec4< float > bbox={0.f, 1.f, 0.f, 1.f})
 Generates a heightmap representing a radial mountain range.
 
Array noise (NoiseType noise_type, Vec2< int > shape, Vec2< float > kw, uint seed, const Array *p_noise_x=nullptr, const Array *p_noise_y=nullptr, const Array *p_stretching=nullptr, Vec4< float > bbox={0.f, 1.f, 0.f, 1.f})
 See hmap::noise.
 
Array noise_fbm (NoiseType noise_type, Vec2< int > shape, Vec2< float > kw, uint seed, int octaves=8, float weight=0.7f, float persistence=0.5f, float lacunarity=2.f, const Array *p_ctrl_param=nullptr, const Array *p_noise_x=nullptr, const Array *p_noise_y=nullptr, const Array *p_stretching=nullptr, Vec4< float > bbox={0.f, 1.f, 0.f, 1.f})
 See hmap::noise_fbm.
 
Array polygon_field (Vec2< int > shape, Vec2< float > kw, uint seed, float rmin=0.05f, float rmax=0.8f, float clamping_dist=0.1f, float clamping_k=0.1f, int n_vertices_min=3, int n_vertices_max=16, float density=0.5f, hmap::Vec2< float > jitter={0.5f, 0.5f}, float shift=0.1f, const Array *p_noise_x=nullptr, const Array *p_noise_y=nullptr, const Array *p_noise_distance=nullptr, const Array *p_density_multiplier=nullptr, const Array *p_size_multiplier=nullptr, Vec4< float > bbox={0.f, 1.f, 0.f, 1.f})
 Generates a scalar field representing the signed distance to randomly generated polygons.
 
Array polygon_field_fbm (Vec2< int > shape, Vec2< float > kw, uint seed, float rmin=0.05f, float rmax=0.8f, float clamping_dist=0.1f, float clamping_k=0.1f, int n_vertices_min=3, int n_vertices_max=16, float density=0.1f, hmap::Vec2< float > jitter={0.5f, 0.5f}, float shift=0.1f, int octaves=8, float persistence=0.5f, float lacunarity=2.f, const Array *p_noise_x=nullptr, const Array *p_noise_y=nullptr, const Array *p_noise_distance=nullptr, const Array *p_density_multiplier=nullptr, const Array *p_size_multiplier=nullptr, Vec4< float > bbox={0.f, 1.f, 0.f, 1.f})
 Generates a scalar field representing the signed distance to randomly generated polygons combined with fractal Brownian motion (fBm) noise modulation.
 
Array vorolines (Vec2< int > shape, float density, uint seed, float k_smoothing=0.f, float exp_sigma=0.f, float alpha=0.f, float alpha_span=M_PI, VoronoiReturnType return_type=VoronoiReturnType::F1_SQUARED, const Array *p_noise_x=nullptr, const Array *p_noise_y=nullptr, Vec4< float > bbox={0.f, 1.f, 0.f, 1.f}, Vec4< float > bbox_points={0.f, 1.f, 0.f, 1.f})
 Generates a Voronoi-based pattern where cells are defined by proximity to random lines.
 
Array vorolines_fbm (Vec2< int > shape, float density, uint seed, float k_smoothing=0.f, float exp_sigma=0.f, float alpha=0.f, float alpha_span=M_PI, VoronoiReturnType return_type=VoronoiReturnType::F1_SQUARED, int octaves=8, float weight=0.7f, float persistence=0.5f, float lacunarity=2.f, const Array *p_noise_x=nullptr, const Array *p_noise_y=nullptr, Vec4< float > bbox={0.f, 1.f, 0.f, 1.f}, Vec4< float > bbox_points={0.f, 1.f, 0.f, 1.f})
 Generates a Voronoi-based pattern using distances to lines defined by random points and angles, with additional fractal Brownian motion (fBm) noise modulation.
 
Array voronoi (Vec2< int > shape, Vec2< float > kw, uint seed, Vec2< float > jitter={0.5f, 0.5f}, float k_smoothing=0.f, float exp_sigma=0.f, VoronoiReturnType return_type=VoronoiReturnType::F1_SQUARED, const Array *p_ctrl_param=nullptr, const Array *p_noise_x=nullptr, const Array *p_noise_y=nullptr, Vec4< float > bbox={0.f, 1.f, 0.f, 1.f})
 Generates a Voronoi diagram in a 2D array with configurable properties.
 
Array voronoi_fbm (Vec2< int > shape, Vec2< float > kw, uint seed, Vec2< float > jitter={0.5f, 0.5f}, float k_smoothing=0.f, float exp_sigma=0.f, VoronoiReturnType return_type=VoronoiReturnType::F1_SQUARED, int octaves=8, float weight=0.7f, float persistence=0.5f, float lacunarity=2.f, const Array *p_ctrl_param=nullptr, const Array *p_noise_x=nullptr, const Array *p_noise_y=nullptr, Vec4< float > bbox={0.f, 1.f, 0.f, 1.f})
 Generates a Voronoi diagram in a 2D array with configurable properties.
 
Array voronoi_edge_distance (Vec2< int > shape, Vec2< float > kw, uint seed, Vec2< float > jitter={0.5f, 0.5f}, const Array *p_ctrl_param=nullptr, const Array *p_noise_x=nullptr, const Array *p_noise_y=nullptr, Vec4< float > bbox={0.f, 1.f, 0.f, 1.f})
 Computes the Voronoi edge distance.
 
Array voronoise (Vec2< int > shape, Vec2< float > kw, float u_param, float v_param, uint seed, const Array *p_noise_x=nullptr, const Array *p_noise_y=nullptr, Vec4< float > bbox={0.f, 1.f, 0.f, 1.f})
 Generates a 2D Voronoi noise array.
 
Array voronoise_fbm (Vec2< int > shape, Vec2< float > kw, float u_param, float v_param, uint seed, int octaves=8, float weight=0.7f, float persistence=0.5f, float lacunarity=2.f, const Array *p_ctrl_param=nullptr, const Array *p_noise_x=nullptr, const Array *p_noise_y=nullptr, Vec4< float > bbox={0.f, 1.f, 0.f, 1.f})
 Return an array filled with coherence Voronoise.
 
Array vororand (Vec2< int > shape, float density, float variability, uint seed, float k_smoothing=0.f, float exp_sigma=0.f, VoronoiReturnType return_type=VoronoiReturnType::F1_SQUARED, const Array *p_noise_x=nullptr, const Array *p_noise_y=nullptr, Vec4< float > bbox={0.f, 1.f, 0.f, 1.f}, Vec4< float > bbox_points={0.f, 1.f, 0.f, 1.f})
 Generates a 2D Voronoi-based scalar field using OpenCL.
 
Array vororand (Vec2< int > shape, const std::vector< float > &xp, const std::vector< float > &yp, float k_smoothing=0.f, float exp_sigma=0.f, VoronoiReturnType return_type=VoronoiReturnType::F1_SQUARED, const Array *p_noise_x=nullptr, const Array *p_noise_y=nullptr, Vec4< float > bbox={0.f, 1.f, 0.f, 1.f})
 
Array maximum_smooth (const Array &array1, const Array &array2, float k=0.2f)
 See hmap::maximum_smooth.
 
Array minimum_smooth (const Array &array1, const Array &array2, float k=0.2f)
 See hmap::minimum_smooth.
 
Array sdf_2d_polyline (const Path &path, Vec2< int > shape, Vec4< float > bbox={0.f, 1.f, 0.f, 1.f}, const Array *p_noise_x=nullptr, const Array *p_noise_y=nullptr)
 See hmap::sdf_2d_polyline.
 
Array sdf_2d_polyline_bezier (const Path &path, Vec2< int > shape, Vec4< float > bbox={0.f, 1.f, 0.f, 1.f}, const Array *p_noise_x=nullptr, const Array *p_noise_y=nullptr)
 See hmap::sdf_2d_polyline_bezier.
 
Array select_valley (const Array &z, int ir, bool zero_at_borders=true, bool ridge_select=false)
 See hmap::select_valley.
 
void rotate (Array &array, float angle, bool zoom_in=true)
 See hmap::rotate.
 
void warp (Array &array, const Array *p_dx, const Array *p_dy)
 See hmap::warp.
 
Vec4< float > helper_transform_bbox (const Vec4< float > &bbox_source, const Vec4< float > &bbox_target)
 
void helper_bind_optional_buffers (clwrapper::Run &run, const Array *p_noise_x, const Array *p_noise_y)
 

Function Documentation

◆ blend_gradients()

Array hmap::gpu::blend_gradients ( const Array array1,
const Array array2,
int  ir = 4 
)

◆ blend_poisson_bf()

Array hmap::gpu::blend_poisson_bf ( const Array array1,
const Array array2,
const int  iterations = 500,
const Array p_mask = nullptr 
)

Blends two arrays using Poisson blending with a brute-force solver.

This function performs Poisson blending between array1 and array2 over a specified number of iterations. Optionally, a mask can be provided to control the blending regions.

Parameters
array1The first input array.
array2The second input array.
iterationsThe number of iterations for the blending process (default: 500).
p_maskOptional pointer to an array defining the blending mask. If null, blending is applied globally.
Returns
The blended array resulting from the Poisson blending operation.

Example

#include "highmap.hpp"
int main(void)
{
hmap::Vec2<int> shape = {256, 256};
hmap::Vec2<float> kw = {2.f, 2.f};
int seed = 2;
shape,
2.f * kw,
++seed);
int iterations = 5000;
hmap::Array z3 = hmap::gpu::blend_poisson_bf(z1, z2, iterations);
hmap::export_banner_png("ex_blend_poisson_bf.png",
{z1, z2, z3},
}
Array class, helper to manipulate 2D float array with "(i, j)" indexing.
Definition array.hpp:32
Array blend_poisson_bf(const Array &array1, const Array &array2, const int iterations=500, const Array *p_mask=nullptr)
Blends two arrays using Poisson blending with a brute-force solver.
Definition blending_gpu.cpp:30
bool init_opencl()
Definition gpu_opencl.cpp:24
void export_banner_png(const std::string &fname, const std::vector< Array > &arrays, int cmap, bool hillshading=false)
Exports a set of arrays as a banner PNG image file.
Definition export_banner_png.cpp:11
void remap(Array &array, float vmin, float vmax, float from_min, float from_max)
Remap array elements from a starting range to a target range.
Definition range.cpp:374
Array noise_fbm(NoiseType noise_type, Vec2< int > shape, Vec2< float > kw, uint seed, int octaves=8, float weight=0.7f, float persistence=0.5f, float lacunarity=2.f, const Array *p_ctrl_param=nullptr, const Array *p_noise_x=nullptr, const Array *p_noise_y=nullptr, const Array *p_stretching=nullptr, Vec4< float > bbox={0.f, 1.f, 0.f, 1.f})
Return an array filled with coherence fbm noise.
Definition noise.cpp:41
@ WORLEY
Worley.
Definition functions.hpp:73
@ PERLIN
Perlin.
Definition functions.hpp:64
@ JET
Definition colormaps.hpp:86
Vec2 class for basic manipulation of 2D vectors.
Definition algebra.hpp:40

Result

◆ accumulation_curvature()

Array hmap::gpu::accumulation_curvature ( const Array z,
int  ir 
)

◆ curvature_horizontal_cross_sectional()

Array hmap::gpu::curvature_horizontal_cross_sectional ( const Array z,
int  ir 
)

◆ curvature_horizontal_plan()

Array hmap::gpu::curvature_horizontal_plan ( const Array z,
int  ir 
)

◆ curvature_horizontal_tangential()

Array hmap::gpu::curvature_horizontal_tangential ( const Array z,
int  ir 
)

◆ curvature_ring()

Array hmap::gpu::curvature_ring ( const Array z,
int  ir 
)

◆ curvature_rotor()

Array hmap::gpu::curvature_rotor ( const Array z,
int  ir 
)

◆ curvature_vertical_longitudinal()

Array hmap::gpu::curvature_vertical_longitudinal ( const Array z,
int  ir 
)

◆ curvature_vertical_profile()

Array hmap::gpu::curvature_vertical_profile ( const Array z,
int  ir 
)

◆ shape_index()

Array hmap::gpu::shape_index ( const Array z,
int  ir 
)

◆ unsphericity()

Array hmap::gpu::unsphericity ( const Array z,
int  ir 
)

◆ hydraulic_particle() [1/2]

void hmap::gpu::hydraulic_particle ( Array z,
int  nparticles,
int  seed,
Array p_bedrock = nullptr,
Array p_moisture_map = nullptr,
Array p_erosion_map = nullptr,
Array p_deposition_map = nullptr,
float  c_capacity = 10.f,
float  c_erosion = 0.05f,
float  c_deposition = 0.05f,
float  c_inertia = 0.3f,
float  drag_rate = 0.001f,
float  evap_rate = 0.001f,
bool  post_filtering = false 
)

◆ hydraulic_particle() [2/2]

void hmap::gpu::hydraulic_particle ( Array z,
Array p_mask,
int  nparticles,
int  seed,
Array p_bedrock = nullptr,
Array p_moisture_map = nullptr,
Array p_erosion_map = nullptr,
Array p_deposition_map = nullptr,
float  c_capacity = 10.f,
float  c_erosion = 0.05f,
float  c_deposition = 0.05f,
float  c_inertia = 0.3f,
float  drag_rate = 0.001f,
float  evap_rate = 0.001f,
bool  post_filtering = false 
)

◆ hydraulic_schott() [1/2]

void hmap::gpu::hydraulic_schott ( Array z,
int  iterations,
const Array talus,
float  c_erosion = 1.f,
float  c_thermal = 0.1f,
float  c_deposition = 0.2f,
float  flow_acc_exponent = 0.8f,
float  flow_acc_exponent_depo = 0.8f,
float  flow_routing_exponent = 1.3f,
float  thermal_weight = 1.5f,
float  deposition_weight = 2.5f,
Array p_flow = nullptr 
)

Simulates hydraulic erosion and deposition on a heightmap using the Schott method.

This function performs hydraulic erosion on the given heightmap z over a specified number of iterations. It includes parameters for controlling erosion, deposition, and flow routing. Optional flow accumulation can also be computed and stored in the p_flow array.

Parameters
[in,out]zThe heightmap array to be modified. Heights are updated in-place.
[in]iterationsThe number of iterations for the hydraulic erosion process.
[in]talusAn array defining the slope threshold for erosion.
[in]c_erosionErosion coefficient (default: 1.0f).
[in]c_thermalThermal erosion coefficient (default: 0.1f).
[in]c_depositionDeposition coefficient (default: 0.2f).
[in]flow_acc_exponentExponent controlling the influence of flow accumulation on erosion (default: 0.8f).
[in]flow_acc_exponent_depoExponent controlling the influence of flow accumulation on deposition (default: 0.8f).
[in]flow_routing_exponentExponent controlling flow routing behavior (default: 1.3f).
[in]thermal_weightWeight of thermal erosion effects (default: 1.5f).
[in]deposition_weightWeight of deposition effects (default: 2.5f).
[out]p_flowOptional pointer to an array for storing flow accumulation data. If null, flow data is not returned (default: nullptr).
Note
Taken from https://hal.science/hal-04565030v1/document
Only available if OpenCL is enabled.

Example

#include "highmap.hpp"
int main(void)
{
hmap::Vec2<int> shape = {256, 256};
shape = {512, 512};
hmap::Vec2<float> kw = {2.f, 2.f};
int seed = 2;
auto z0 = z;
int iterations = 400;
hmap::Array talus(shape, 2.f / (float)shape.x); // for thermal
hmap::gpu::hydraulic_schott(z, iterations, talus);
hmap::export_banner_png("ex_hydraulic_schott.png",
{z0, z},
true);
}
void hydraulic_schott(Array &z, int iterations, const Array &talus, float c_erosion=1.f, float c_thermal=0.1f, float c_deposition=0.2f, float flow_acc_exponent=0.8f, float flow_acc_exponent_depo=0.8f, float flow_routing_exponent=1.3f, float thermal_weight=1.5f, float deposition_weight=2.5f, Array *p_flow=nullptr)
Simulates hydraulic erosion and deposition on a heightmap using the Schott method.
Definition hydraulic_schott_gpu.cpp:10
@ TERRAIN
Definition colormaps.hpp:90
T x
Definition algebra.hpp:41

Result

◆ hydraulic_schott() [2/2]

void hmap::gpu::hydraulic_schott ( Array z,
int  iterations,
const Array talus,
Array p_mask,
float  c_erosion = 1.f,
float  c_thermal = 0.1f,
float  c_deposition = 0.2f,
float  flow_acc_exponent = 0.8f,
float  flow_acc_exponent_depo = 0.8f,
float  flow_routing_exponent = 1.3f,
float  thermal_weight = 1.5f,
float  deposition_weight = 2.5f,
Array p_flow = nullptr 
)

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

◆ hydraulic_stream_log() [1/2]

void hmap::gpu::hydraulic_stream_log ( Array z,
float  c_erosion,
float  talus_ref,
int  deposition_ir = 32,
float  deposition_scale_ratio = 1.f,
float  gradient_power = 0.8f,
float  gradient_scaling_ratio = 1.f,
int  gradient_prefilter_ir = 16,
float  saturation_ratio = 1.f,
Array p_bedrock = nullptr,
Array p_moisture_map = nullptr,
Array p_erosion_map = nullptr,
Array p_deposition_map = nullptr,
Array p_flow_map = nullptr 
)

◆ hydraulic_stream_log() [2/2]

void hmap::gpu::hydraulic_stream_log ( Array z,
float  c_erosion,
float  talus_ref,
Array p_mask,
int  deposition_ir = 32,
float  deposition_scale_ratio = 1.f,
float  gradient_power = 0.8f,
float  gradient_scaling_ratio = 1.f,
int  gradient_prefilter_ir = 16,
float  saturation_ratio = 1.f,
Array p_bedrock = nullptr,
Array p_moisture_map = nullptr,
Array p_erosion_map = nullptr,
Array p_deposition_map = nullptr,
Array p_flow_map = nullptr 
)

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

◆ thermal() [1/3]

void hmap::gpu::thermal ( Array z,
const Array talus,
int  iterations = 10,
Array p_bedrock = nullptr,
Array p_deposition_map = nullptr 
)

◆ thermal() [2/3]

void hmap::gpu::thermal ( Array z,
Array p_mask,
const Array talus,
int  iterations = 10,
Array p_bedrock = nullptr,
Array p_deposition_map = nullptr 
)

◆ thermal() [3/3]

void hmap::gpu::thermal ( Array z,
float  talus,
int  iterations = 10,
Array p_bedrock = nullptr,
Array p_deposition_map = nullptr 
)

◆ thermal_auto_bedrock() [1/3]

void hmap::gpu::thermal_auto_bedrock ( Array z,
const Array talus,
int  iterations = 10,
Array p_deposition_map = nullptr 
)

◆ thermal_auto_bedrock() [2/3]

void hmap::gpu::thermal_auto_bedrock ( Array z,
Array p_mask,
const Array talus,
int  iterations = 10,
Array p_deposition_map = nullptr 
)

◆ thermal_auto_bedrock() [3/3]

void hmap::gpu::thermal_auto_bedrock ( Array z,
float  talus,
int  iterations = 10,
Array p_deposition_map = nullptr 
)

◆ thermal_inflate() [1/2]

void hmap::gpu::thermal_inflate ( Array z,
const Array talus,
int  iterations = 10 
)

Apply thermal weathering erosion to give a scree like effect.

Note
Only available if OpenCL is enabled.
Parameters
zInput array.
talusTalus limit.
p_deposition_map[out] Reference to the deposition map, provided as an output field.
iterationsNumber of iterations.

Example

Result

◆ thermal_inflate() [2/2]

void hmap::gpu::thermal_inflate ( Array z,
const Array p_mask,
const Array talus,
int  iterations = 10 
)

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

◆ thermal_rib()

void hmap::gpu::thermal_rib ( Array z,
int  iterations,
Array p_bedrock = nullptr 
)

◆ thermal_ridge() [1/2]

void hmap::gpu::thermal_ridge ( Array z,
const Array talus,
int  iterations = 10,
Array p_deposition_map = nullptr 
)

Apply thermal weathering erosion to give a ridge like effect.

Note
Based on https://www.fractal-landscapes.co.uk/maths.html
Only available if OpenCL is enabled.
Parameters
zInput array.
talusTalus limit.
p_deposition_map[out] Reference to the deposition map, provided as an output field.
iterationsNumber of iterations.

Example

Result

◆ thermal_ridge() [2/2]

void hmap::gpu::thermal_ridge ( Array z,
const Array p_mask,
const Array talus,
int  iterations = 10,
Array p_deposition_map = nullptr 
)

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

◆ thermal_scree() [1/2]

void hmap::gpu::thermal_scree ( Array z,
const Array talus,
const Array zmax,
int  iterations = 10,
bool  talus_constraint = true,
Array p_deposition_map = nullptr 
)

Performs thermal scree erosion on a heightmap.

This function applies a thermal erosion process that redistributes material from steeper slopes to flatter areas, simulating talus formation. The process iterates a given number of times to achieve a more stable terrain profile.

Parameters
[out]zThe heightmap to be modified in-place by the erosion process.
[in]talusThe threshold slope angles that determine where material is moved.
[in]zmaxThe maximum allowed elevation for erosion effects.
[in]iterationsThe number of erosion iterations to apply (default: 10).
[in]talus_constraintWhether to enforce a constraint on the talus slope (default: true).
[out]p_deposition_mapOptional pointer to an array that stores the deposited material per cell.

◆ thermal_scree() [2/2]

void hmap::gpu::thermal_scree ( Array z,
const Array p_mask,
const Array talus,
const Array zmax,
int  iterations = 10,
bool  talus_constraint = true,
Array p_deposition_map = nullptr 
)

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

◆ local_median_deviation()

Array hmap::gpu::local_median_deviation ( const Array array,
int  ir 
)

◆ mean_local()

Array hmap::gpu::mean_local ( const Array array,
int  ir 
)

◆ relative_elevation()

Array hmap::gpu::relative_elevation ( const Array array,
int  ir 
)

◆ ruggedness()

Array hmap::gpu::ruggedness ( const Array array,
int  ir 
)

◆ rugosity()

Array hmap::gpu::rugosity ( const Array z,
int  ir,
bool  convex = true 
)

◆ std_local()

Array hmap::gpu::std_local ( const Array array,
int  ir 
)

◆ z_score()

Array hmap::gpu::z_score ( const Array array,
int  ir 
)

◆ expand() [1/4]

void hmap::gpu::expand ( Array array,
int  ir,
int  iterations = 1 
)

◆ expand() [2/4]

void hmap::gpu::expand ( Array array,
int  ir,
const Array p_mask,
int  iterations = 1 
)

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

◆ expand() [3/4]

void hmap::gpu::expand ( Array array,
const Array kernel,
int  iterations = 1 
)

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

◆ expand() [4/4]

void hmap::gpu::expand ( Array array,
const Array kernel,
const Array p_mask,
int  iterations = 1 
)

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

◆ gamma_correction_local() [1/2]

void hmap::gpu::gamma_correction_local ( Array array,
float  gamma,
int  ir,
float  k = 0.1f 
)

◆ gamma_correction_local() [2/2]

void hmap::gpu::gamma_correction_local ( Array array,
float  gamma,
int  ir,
const Array p_mask,
float  k = 0.1f 
)

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

◆ laplace() [1/2]

void hmap::gpu::laplace ( Array array,
float  sigma = 0.2f,
int  iterations = 3 
)

◆ laplace() [2/2]

void hmap::gpu::laplace ( Array array,
const Array p_mask,
float  sigma = 0.2f,
int  iterations = 3 
)

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

◆ maximum_local()

Array hmap::gpu::maximum_local ( const Array array,
int  ir 
)

◆ maximum_local_disk()

Array hmap::gpu::maximum_local_disk ( const Array array,
int  ir 
)

◆ mean_shift() [1/2]

Array hmap::gpu::mean_shift ( const Array array,
int  ir,
float  talus,
int  iterations = 1,
bool  talus_weighted = true 
)

◆ mean_shift() [2/2]

Array hmap::gpu::mean_shift ( const Array array,
int  ir,
float  talus,
const Array p_mask,
int  iterations = 1,
bool  talus_weighted = true 
)

◆ median_3x3() [1/2]

void hmap::gpu::median_3x3 ( Array array)

◆ median_3x3() [2/2]

void hmap::gpu::median_3x3 ( Array array,
const Array p_mask 
)

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

◆ median_pseudo()

Array hmap::gpu::median_pseudo ( const Array array,
int  ir 
)

◆ minimum_local()

Array hmap::gpu::minimum_local ( const Array array,
int  ir 
)

◆ minimum_local_disk()

Array hmap::gpu::minimum_local_disk ( const Array array,
int  ir 
)

◆ normal_displacement() [1/2]

void hmap::gpu::normal_displacement ( Array array,
float  amount = 0.1f,
int  ir = 0,
bool  reverse = false 
)

◆ normal_displacement() [2/2]

void hmap::gpu::normal_displacement ( Array array,
const Array p_mask,
float  amount = 0.1f,
int  ir = 0,
bool  reverse = false 
)

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

◆ plateau() [1/2]

void hmap::gpu::plateau ( Array array,
const Array p_mask,
int  ir,
float  factor 
)

◆ plateau() [2/2]

void hmap::gpu::plateau ( Array array,
int  ir,
float  factor 
)

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

◆ shrink() [1/4]

void hmap::gpu::shrink ( Array array,
int  ir,
int  iterations = 1 
)

◆ shrink() [2/4]

void hmap::gpu::shrink ( Array array,
int  ir,
const Array p_mask,
int  iterations = 1 
)

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

◆ shrink() [3/4]

void hmap::gpu::shrink ( Array array,
const Array kernel,
int  iterations = 1 
)

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

◆ shrink() [4/4]

void hmap::gpu::shrink ( Array array,
const Array kernel,
const Array p_mask,
int  iterations = 1 
)

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

◆ smooth_cpulse() [1/2]

void hmap::gpu::smooth_cpulse ( Array array,
int  ir 
)

◆ smooth_cpulse() [2/2]

void hmap::gpu::smooth_cpulse ( Array array,
int  ir,
const Array p_mask 
)

◆ smooth_fill() [1/2]

void hmap::gpu::smooth_fill ( Array array,
int  ir,
float  k = 0.1f,
Array p_deposition_map = nullptr 
)

◆ smooth_fill() [2/2]

void hmap::gpu::smooth_fill ( Array array,
int  ir,
const Array p_mask,
float  k = 0.1f,
Array p_deposition_map = nullptr 
)

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

◆ smooth_fill_holes() [1/2]

void hmap::gpu::smooth_fill_holes ( Array array,
int  ir 
)

◆ smooth_fill_holes() [2/2]

void hmap::gpu::smooth_fill_holes ( Array array,
int  ir,
const Array p_mask 
)

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

◆ smooth_fill_smear_peaks() [1/2]

void hmap::gpu::smooth_fill_smear_peaks ( Array array,
int  ir 
)

◆ smooth_fill_smear_peaks() [2/2]

void hmap::gpu::smooth_fill_smear_peaks ( Array array,
int  ir,
const Array p_mask 
)

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

◆ gradient_norm()

Array hmap::gpu::gradient_norm ( const Array array)

◆ flow_direction_d8()

Array hmap::gpu::flow_direction_d8 ( const Array z)

◆ generate_riverbed()

Array hmap::gpu::generate_riverbed ( const Path path,
Vec2< int >  shape,
Vec4< float >  bbox = {0.f, 1.f, 0.f, 1.f},
bool  bezier_smoothing = false,
float  depth_start = 0.01f,
float  depth_end = 1.f,
float  slope_start = 64.f,
float  slope_end = 32.f,
float  shape_exponent_start = 1.f,
float  shape_exponent_end = 10.f,
float  k_smoothing = 0.5f,
int  post_filter_ir = 0,
Array p_noise_x = nullptr,
Array p_noise_y = nullptr,
Array p_noise_r = nullptr 
)

◆ interpolate_array_bicubic() [1/2]

void hmap::gpu::interpolate_array_bicubic ( const Array source,
Array target 
)

◆ interpolate_array_bicubic() [2/2]

void hmap::gpu::interpolate_array_bicubic ( const Array source,
Array target,
const Vec4< float > &  bbox_source,
const Vec4< float > &  bbox_target 
)

◆ interpolate_array_bilinear() [1/2]

void hmap::gpu::interpolate_array_bilinear ( const Array source,
Array target 
)

◆ interpolate_array_bilinear() [2/2]

void hmap::gpu::interpolate_array_bilinear ( const Array source,
Array target,
const Vec4< float > &  bbox_source,
const Vec4< float > &  bbox_target 
)

◆ interpolate_array_lagrange()

void hmap::gpu::interpolate_array_lagrange ( const Array source,
Array target,
int  order 
)

◆ interpolate_array_nearest() [1/2]

void hmap::gpu::interpolate_array_nearest ( const Array source,
Array target 
)

◆ interpolate_array_nearest() [2/2]

void hmap::gpu::interpolate_array_nearest ( const Array source,
Array target,
const Vec4< float > &  bbox_source,
const Vec4< float > &  bbox_target 
)

◆ border()

Array hmap::gpu::border ( const Array array,
int  ir 
)

◆ closing()

Array hmap::gpu::closing ( const Array array,
int  ir 
)

◆ dilation()

Array hmap::gpu::dilation ( const Array array,
int  ir 
)

◆ erosion()

Array hmap::gpu::erosion ( const Array array,
int  ir 
)

◆ morphological_black_hat()

Array hmap::gpu::morphological_black_hat ( const Array array,
int  ir 
)

◆ morphological_gradient()

Array hmap::gpu::morphological_gradient ( const Array array,
int  ir 
)

◆ morphological_top_hat()

Array hmap::gpu::morphological_top_hat ( const Array array,
int  ir 
)

◆ opening()

Array hmap::gpu::opening ( const Array array,
int  ir 
)

◆ relative_distance_from_skeleton()

Array hmap::gpu::relative_distance_from_skeleton ( const Array array,
int  ir_search,
bool  zero_at_borders = true,
int  ir_erosion = 1 
)

◆ skeleton()

Array hmap::gpu::skeleton ( const Array array,
bool  zero_at_borders = true 
)

◆ helper_bind_optional_buffer()

void hmap::gpu::helper_bind_optional_buffer ( clwrapper::Run &  run,
const std::string &  id,
const Array p_array 
)

◆ init_opencl()

bool hmap::gpu::init_opencl ( )

◆ basalt_field()

Array hmap::gpu::basalt_field ( Vec2< int >  shape,
Vec2< float >  kw,
uint  seed,
float  warp_kw = 4.f,
float  large_scale_warp_amp = 0.2f,
float  large_scale_gain = 6.f,
float  large_scale_amp = 0.2f,
float  medium_scale_kw_ratio = 3.f,
float  medium_scale_warp_amp = 1.f,
float  medium_scale_gain = 7.f,
float  medium_scale_amp = 0.08f,
float  small_scale_kw_ratio = 10.f,
float  small_scale_amp = 0.1f,
float  small_scale_overlay_amp = 0.002f,
float  rugosity_kw_ratio = 1.f,
float  rugosity_amp = 1.f,
bool  flatten_activate = true,
float  flatten_kw_ratio = 1.f,
float  flatten_amp = 0.f,
const Array p_noise_x = nullptr,
const Array p_noise_y = nullptr,
Vec4< float >  bbox = {0.f, 1.f, 0.f, 1.f} 
)

Generates a synthetic procedural terrain resembling basaltic landforms.

This function creates a multi-scale procedural field combining large, medium, and small-scale Voronoi-based patterns, noise warping, and optional flattening, simulating the morphology of fractured basalt or volcanic terrains. The terrain is constructed using a combination of Voronoi diagrams (via voronoi_fbm) and fractal noise (noise_fbm), layered with frequency-domain manipulations and amplitude/gain controls at each scale.

The final output is a heightmap represented as an Array, normalized and composed of:

  • Large-scale cellular patterns with smoothed Voronoi edge distances.
  • Medium and small-scale structures introducing finer surface variation.
  • Optional rugosity (fine detail) and flattening to simulate erosion or flow effects.
Parameters
shapeOutput resolution (width x height) of the field.
kwBase wave numbers (frequency) for the terrain features.
seedInitial seed used for deterministic random generation.
warp_kwFrequency of the warping noise that displaces Voronoi positions.
large_scale_warp_ampAmplitude of displacement for large-scale Voronoi warping.
large_scale_gainGain adjustment applied to the large-scale features.
large_scale_ampFinal amplitude of the large-scale height contribution.
medium_scale_kw_ratioScaling factor for the frequency of the medium-scale patterns.
medium_scale_warp_ampAmplitude of warping for the medium-scale displacement.
medium_scale_gainGain control for medium-scale modulation.
medium_scale_ampAmplitude of the medium-scale heightmap.
small_scale_kw_ratioFrequency ratio for small-scale details.
small_scale_ampAmplitude of small-scale pattern contribution.
small_scale_overlay_ampAdditional overlay strength for repeating the small-scale pattern.
rugosity_kw_ratioFrequency ratio for high-frequency noise applied as fine roughness.
rugosity_ampStrength of the rugosity (high-frequency modulation).
flatten_activateEnables or disables the final flattening operation.
flatten_kw_ratioFrequency scaling of the flattening noise field.
flatten_ampAmplitude control of the flattening operation.
p_noise_xOptional pointer to a noise field used to displace grid coordinates in X.
p_noise_yOptional pointer to a noise field used to displace grid coordinates in Y.
bboxThe 2D bounding box ({xmin, xmax, ymin, ymax}) over which the terrain is generated.
Returns
A procedurally generated Array representing the synthetic basalt-like terrain field.
Note
  • This function relies on OpenCL-based kernels via the gpu:: namespace.
  • The returned field is normalized in amplitude but may require rescaling to match specific physical units.
  • Adjusting seed, warp_kw, and the gain/amplitude values can produce a wide variety of terrain features.

Example

#include "highmap.hpp"
int main(void)
{
hmap::Vec2<int> shape = {256, 256};
hmap::Vec2<float> kw = {4.f, 4.f};
int seed = 0;
hmap::Array z = hmap::gpu::basalt_field(shape, kw, seed);
hmap::export_banner_png("ex_basalt_field.png",
{z},
true);
}
Array basalt_field(Vec2< int > shape, Vec2< float > kw, uint seed, float warp_kw=4.f, float large_scale_warp_amp=0.2f, float large_scale_gain=6.f, float large_scale_amp=0.2f, float medium_scale_kw_ratio=3.f, float medium_scale_warp_amp=1.f, float medium_scale_gain=7.f, float medium_scale_amp=0.08f, float small_scale_kw_ratio=10.f, float small_scale_amp=0.1f, float small_scale_overlay_amp=0.002f, float rugosity_kw_ratio=1.f, float rugosity_amp=1.f, bool flatten_activate=true, float flatten_kw_ratio=1.f, float flatten_amp=0.f, const Array *p_noise_x=nullptr, const Array *p_noise_y=nullptr, Vec4< float > bbox={0.f, 1.f, 0.f, 1.f})
Generates a synthetic procedural terrain resembling basaltic landforms.
Definition basalt_field.cpp:13
@ INFERNO
Definition colormaps.hpp:85

Result

◆ gabor_wave() [1/2]

Array hmap::gpu::gabor_wave ( Vec2< int >  shape,
Vec2< float >  kw,
uint  seed,
const Array angle,
float  angle_spread_ratio = 1.f,
Vec4< float >  bbox = {0.f, 1.f, 0.f, 1.f} 
)

Return an array filled with coherence Gabor noise.

Parameters
shapeArray shape.
kwNoise wavenumbers {kx, ky} for each directions.
seedRandom seed number.
angleBase orientation angle for the Gabor wavelets (in radians). Defaults to 0.
angle_spread_ratioRatio that controls the spread of wave orientations around the base angle. Defaults to 1.
bboxDomain bounding box.
Returns
Array Fractal noise.
Note
Taken from https://www.shadertoy.com/view/clGyWm
Only available if OpenCL is enabled.

Example

#include "highmap.hpp"
int main(void)
{
hmap::Vec2<int> shape = {256, 256};
hmap::Vec2<float> kw = {2.f, 2.f};
int seed = 1;
// --- base
hmap::Array z = hmap::gpu::gabor_wave(shape, kw, seed);
hmap::Array z_fbm = hmap::gpu::gabor_wave_fbm(shape, kw, seed);
// --- angle control
float angle = 45.f;
float angle_spread_ratio = 0.f;
{16.f, 16.f},
seed,
angle_spread_ratio);
hmap::Array za1 = hmap::gpu::gabor_wave(shape, {16.f, 16.f}, seed, 0.f, 0.5f);
{16.f, 16.f},
seed,
0.f,
0.1f);
// --- local angle
hmap::Array field = hmap::noise(hmap::NoiseType::PERLIN, shape, kw, seed);
hmap::Array array_angle = hmap::gradient_angle(field) * 180.f / 3.14159f;
angle_spread_ratio = 0.f;
{16.f, 16.f},
seed,
array_angle,
angle_spread_ratio);
{16.f, 16.f},
seed,
array_angle,
angle_spread_ratio);
hmap::export_banner_png("ex_gabor_wave.png",
{z, z_fbm, za0, za1, za2, zr1, zr2},
true);
}
Array gabor_wave(Vec2< int > shape, Vec2< float > kw, uint seed, const Array &angle, float angle_spread_ratio=1.f, Vec4< float > bbox={0.f, 1.f, 0.f, 1.f})
Return an array filled with coherence Gabor noise.
Definition primitives_gpu.cpp:42
Array gabor_wave_fbm(Vec2< int > shape, Vec2< float > kw, uint seed, const Array &angle, float angle_spread_ratio=1.f, int octaves=8, float weight=0.7f, float persistence=0.5f, float lacunarity=2.f, const Array *p_ctrl_param=nullptr, const Array *p_noise_x=nullptr, const Array *p_noise_y=nullptr, Vec4< float > bbox={0.f, 1.f, 0.f, 1.f})
Return an array filled with coherence Gabor noise.
Definition primitives_gpu.cpp:86
Array noise(NoiseType noise_type, Vec2< int > shape, Vec2< float > kw, uint seed, const Array *p_noise_x=nullptr, const Array *p_noise_y=nullptr, const Array *p_stretching=nullptr, Vec4< float > bbox={0.f, 1.f, 0.f, 1.f})
Return an array filled with coherence noise.
Definition noise.cpp:16
float angle(const Point &p1, const Point &p2)
Computes the angle between two points relative to the x-axis.
Definition points.cpp:42
Array gradient_angle(const Array &array, bool downward=false)
Compute the polar angle of the gradient of a 2D array.
Definition gradient.cpp:60

Result

◆ gabor_wave() [2/2]

Array hmap::gpu::gabor_wave ( Vec2< int >  shape,
Vec2< float >  kw,
uint  seed,
float  angle = 0.f,
float  angle_spread_ratio = 1.f,
Vec4< float >  bbox = {0.f, 1.f, 0.f, 1.f} 
)

◆ gabor_wave_fbm() [1/2]

Array hmap::gpu::gabor_wave_fbm ( Vec2< int >  shape,
Vec2< float >  kw,
uint  seed,
const Array angle,
float  angle_spread_ratio = 1.f,
int  octaves = 8,
float  weight = 0.7f,
float  persistence = 0.5f,
float  lacunarity = 2.f,
const Array p_ctrl_param = nullptr,
const Array p_noise_x = nullptr,
const Array p_noise_y = nullptr,
Vec4< float >  bbox = {0.f, 1.f, 0.f, 1.f} 
)

Return an array filled with coherence Gabor noise.

Parameters
shapeArray shape.
kwNoise wavenumbers {kx, ky} for each directions.
seedRandom seed number.
angleBase orientation angle for the Gabor wavelets (in radians). Defaults to 0.
angle_spread_ratioRatio that controls the spread of wave orientations around the base angle. Defaults to 1.
octavesNumber of octaves.
weigthOctave weighting.
persistenceOctave persistence.
lacunarityDefines the wavenumber ratio between each octaves.
p_ctrl_paramReference to the control parameter array (acts as a multiplier for the weight parameter).
p_noise_x,p_noise_yReference to the input noise arrays.
bboxDomain bounding box.
Returns
Array Fractal noise.
Note
Taken from https://www.shadertoy.com/view/clGyWm
Only available if OpenCL is enabled.

Example

#include "highmap.hpp"
int main(void)
{
hmap::Vec2<int> shape = {256, 256};
hmap::Vec2<float> kw = {2.f, 2.f};
int seed = 1;
// --- base
hmap::Array z = hmap::gpu::gabor_wave(shape, kw, seed);
hmap::Array z_fbm = hmap::gpu::gabor_wave_fbm(shape, kw, seed);
// --- angle control
float angle = 45.f;
float angle_spread_ratio = 0.f;
{16.f, 16.f},
seed,
angle_spread_ratio);
hmap::Array za1 = hmap::gpu::gabor_wave(shape, {16.f, 16.f}, seed, 0.f, 0.5f);
{16.f, 16.f},
seed,
0.f,
0.1f);
// --- local angle
hmap::Array field = hmap::noise(hmap::NoiseType::PERLIN, shape, kw, seed);
hmap::Array array_angle = hmap::gradient_angle(field) * 180.f / 3.14159f;
angle_spread_ratio = 0.f;
{16.f, 16.f},
seed,
array_angle,
angle_spread_ratio);
{16.f, 16.f},
seed,
array_angle,
angle_spread_ratio);
hmap::export_banner_png("ex_gabor_wave.png",
{z, z_fbm, za0, za1, za2, zr1, zr2},
true);
}

Result

◆ gabor_wave_fbm() [2/2]

Array hmap::gpu::gabor_wave_fbm ( Vec2< int >  shape,
Vec2< float >  kw,
uint  seed,
float  angle = 0.f,
float  angle_spread_ratio = 1.f,
int  octaves = 8,
float  weight = 0.7f,
float  persistence = 0.5f,
float  lacunarity = 2.f,
const Array p_ctrl_param = nullptr,
const Array p_noise_x = nullptr,
const Array p_noise_y = nullptr,
Vec4< float >  bbox = {0.f, 1.f, 0.f, 1.f} 
)

◆ gavoronoise() [1/3]

Array hmap::gpu::gavoronoise ( Vec2< int >  shape,
Vec2< float >  kw,
uint  seed,
const Array angle,
float  amplitude = 0.05f,
float  angle_spread_ratio = 1.f,
Vec2< float >  kw_multiplier = {4.f, 4.f},
float  slope_strength = 1.f,
float  branch_strength = 2.f,
float  z_cut_min = 0.2f,
float  z_cut_max = 1.f,
int  octaves = 8,
float  persistence = 0.4f,
float  lacunarity = 2.f,
const Array p_ctrl_param = nullptr,
const Array p_noise_x = nullptr,
const Array p_noise_y = nullptr,
Vec4< float >  bbox = {0.f, 1.f, 0.f, 1.f} 
)

Generates a 2D array using the GavoroNoise algorithm, which is a procedural noise technique for terrain generation and other applications.

Parameters
shapeDimensions of the output array.
kwWave number vector controlling the noise frequency.
seedSeed value for random number generation.
amplitudeAmplitude of the noise.
kw_multiplierMultiplier for wave numbers in the noise function.
slope_strengthStrength of slope-based directional erosion in the noise.
branch_strengthStrength of branch-like structures in the generated noise.
z_cut_minMinimum cutoff for Z-value in the noise.
z_cut_maxMaximum cutoff for Z-value in the noise.
octavesNumber of octaves for fractal Brownian motion (fBm).
persistenceAmplitude scaling factor between noise octaves.
lacunarityFrequency scaling factor between noise octaves.
p_ctrl_paramOptional array for control parameters, can modify the Z cutoff dynamically.
p_noise_xOptional array for X-axis noise perturbation.
p_noise_yOptional array for Y-axis noise perturbation.
bboxBounding box for mapping grid coordinates to world space.
Returns
A 2D array containing the generated GavoroNoise values.
Note
Taken from https://www.shadertoy.com/view/MtGcWh
Only available if OpenCL is enabled.

This function leverages an OpenCL kernel to compute the GavoroNoise values on the GPU, allowing for efficient large-scale generation. The kernel applies a combination of fractal Brownian motion (fBm), directional erosion, and other procedural techniques to generate intricate noise patterns.

The optional p_ctrl_param, p_noise_x, and p_noise_y buffers provide additional flexibility for dynamically adjusting noise parameters and perturbations.

Example

#include "highmap.hpp"
int main(void)
{
hmap::Vec2<int> shape = {256, 256};
hmap::Vec2<float> kw = {2.f, 2.f};
int seed = 1;
// --- base usage
hmap::Array z1 = hmap::gpu::gavoronoise(shape, kw, seed);
float amp = 0.05f;
float angle = 45.f;
float angle_spread_ratio = 0.f;
kw,
seed,
amp,
angle_spread_ratio);
// --- with input base noise
int octaves = 2;
shape,
kw,
++seed,
octaves);
// base amplitude amplitude expected to be in [-1, 1] (approx.)
hmap::remap(base, -1.f, 1.f);
hmap::Array z3 = hmap::gpu::gavoronoise(base, kw, seed, amp);
// --- local angle
hmap::Array field = hmap::noise(hmap::NoiseType::PERLIN, shape, kw, seed);
hmap::Array array_angle = hmap::gradient_angle(field) * 180.f / 3.14159f;
angle_spread_ratio = 0.f;
{8.f, 8.f},
seed,
array_angle,
amp,
angle_spread_ratio);
hmap::export_banner_png("ex_gavoronoise.png",
{z1, z2, z3, z4},
true);
}
Array gavoronoise(Vec2< int > shape, Vec2< float > kw, uint seed, const Array &angle, float amplitude=0.05f, float angle_spread_ratio=1.f, Vec2< float > kw_multiplier={4.f, 4.f}, float slope_strength=1.f, float branch_strength=2.f, float z_cut_min=0.2f, float z_cut_max=1.f, int octaves=8, float persistence=0.4f, float lacunarity=2.f, const Array *p_ctrl_param=nullptr, const Array *p_noise_x=nullptr, const Array *p_noise_y=nullptr, Vec4< float > bbox={0.f, 1.f, 0.f, 1.f})
Generates a 2D array using the GavoroNoise algorithm, which is a procedural noise technique for terra...
Definition primitives_gpu.cpp:169

Result

◆ gavoronoise() [2/3]

Array hmap::gpu::gavoronoise ( Vec2< int >  shape,
Vec2< float >  kw,
uint  seed,
float  angle = 0.f,
float  amplitude = 0.05f,
float  angle_spread_ratio = 1.f,
Vec2< float >  kw_multiplier = {4.f, 4.f},
float  slope_strength = 1.f,
float  branch_strength = 2.f,
float  z_cut_min = 0.2f,
float  z_cut_max = 1.f,
int  octaves = 8,
float  persistence = 0.4f,
float  lacunarity = 2.f,
const Array p_ctrl_param = nullptr,
const Array p_noise_x = nullptr,
const Array p_noise_y = nullptr,
Vec4< float >  bbox = {0.f, 1.f, 0.f, 1.f} 
)

◆ gavoronoise() [3/3]

Array hmap::gpu::gavoronoise ( const Array base,
Vec2< float >  kw,
uint  seed,
float  amplitude = 0.05f,
Vec2< float >  kw_multiplier = {4.f, 4.f},
float  slope_strength = 1.f,
float  branch_strength = 2.f,
float  z_cut_min = 0.2f,
float  z_cut_max = 1.f,
int  octaves = 8,
float  persistence = 0.4f,
float  lacunarity = 2.f,
const Array p_ctrl_param = nullptr,
const Array p_noise_x = nullptr,
const Array p_noise_y = nullptr,
Vec4< float >  bbox = {0.f, 1.f, 0.f, 1.f} 
)

◆ mountain_range_radial()

Array hmap::gpu::mountain_range_radial ( Vec2< int >  shape,
Vec2< float >  kw,
uint  seed,
float  half_width = 0.2f,
float  angle_spread_ratio = 0.5f,
float  core_size_ratio = 1.f,
Vec2< float >  center = {0.5f, 0.5f},
int  octaves = 8,
float  weight = 0.7f,
float  persistence = 0.5f,
float  lacunarity = 2.f,
const Array p_ctrl_param = nullptr,
const Array p_noise_x = nullptr,
const Array p_noise_y = nullptr,
const Array p_angle = nullptr,
Vec4< float >  bbox = {0.f, 1.f, 0.f, 1.f} 
)

Generates a heightmap representing a radial mountain range.

This function creates a heightmap that simulates a mountain range emanating radially from a specified center. The mountain range is influenced by various noise parameters and control attributes.

Parameters
shapeThe dimensions of the output heightmap as a 2D vector.
kwThe wave numbers (frequency components) as a 2D vector.
seedThe seed for random noise generation.
half_widthThe half-width of the radial mountain range, controlling its spread. Default is 0.2f.
angle_spread_ratioThe ratio controlling the angular spread of the mountain range. Default is 0.5f.
centerThe center point of the radial mountain range as normalized coordinates within [0, 1]. Default is {0.5f, 0.5f}.
octavesThe number of octaves for fractal noise generation. Default is 8.
weightThe initial weight for noise contribution. Default is 0.7f.
persistenceThe amplitude scaling factor for subsequent noise octaves. Default is 0.5f.
lacunarityThe frequency scaling factor for subsequent noise octaves. Default is 2.0f.
p_ctrl_paramOptional pointer to an array of control parameters influencing the terrain generation.
p_noise_xOptional pointer to a precomputed noise array for the X-axis.
p_noise_yOptional pointer to a precomputed noise array for the Y-axis.
p_angleOptional pointer to an array to output the angle.
bboxThe bounding box of the output heightmap in normalized coordinates [xmin, xmax, ymin, ymax]. Default is {0.0f, 1.0f, 0.0f, 1.0f}.
Returns
Array The generated heightmap representing the radial mountain range.
Note
Only available if OpenCL is enabled.

Example

#include "highmap.hpp"
int main(void)
{
hmap::Vec2<int> shape = {256, 256};
shape = {1024, 1024};
hmap::Vec2<float> kw = {8.f, 8.f};
int seed = 1;
hmap::export_banner_png("ex_mountain_range_radial.png",
{z1},
true);
}
Array mountain_range_radial(Vec2< int > shape, Vec2< float > kw, uint seed, float half_width=0.2f, float angle_spread_ratio=0.5f, float core_size_ratio=1.f, Vec2< float > center={0.5f, 0.5f}, int octaves=8, float weight=0.7f, float persistence=0.5f, float lacunarity=2.f, const Array *p_ctrl_param=nullptr, const Array *p_noise_x=nullptr, const Array *p_noise_y=nullptr, const Array *p_angle=nullptr, Vec4< float > bbox={0.f, 1.f, 0.f, 1.f})
Generates a heightmap representing a radial mountain range.
Definition primitives_gpu.cpp:328

Result

◆ noise()

Array hmap::gpu::noise ( NoiseType  noise_type,
Vec2< int >  shape,
Vec2< float >  kw,
uint  seed,
const Array p_noise_x = nullptr,
const Array p_noise_y = nullptr,
const Array p_stretching = nullptr,
Vec4< float >  bbox = {0.f, 1.f, 0.f, 1.f} 
)

◆ noise_fbm()

Array hmap::gpu::noise_fbm ( NoiseType  noise_type,
Vec2< int >  shape,
Vec2< float >  kw,
uint  seed,
int  octaves = 8,
float  weight = 0.7f,
float  persistence = 0.5f,
float  lacunarity = 2.f,
const Array p_ctrl_param = nullptr,
const Array p_noise_x = nullptr,
const Array p_noise_y = nullptr,
const Array p_stretching = nullptr,
Vec4< float >  bbox = {0.f, 1.f, 0.f, 1.f} 
)

◆ polygon_field()

Array hmap::gpu::polygon_field ( Vec2< int >  shape,
Vec2< float >  kw,
uint  seed,
float  rmin = 0.05f,
float  rmax = 0.8f,
float  clamping_dist = 0.1f,
float  clamping_k = 0.1f,
int  n_vertices_min = 3,
int  n_vertices_max = 16,
float  density = 0.5f,
hmap::Vec2< float >  jitter = {0.5f, 0.5f},
float  shift = 0.1f,
const Array p_noise_x = nullptr,
const Array p_noise_y = nullptr,
const Array p_noise_distance = nullptr,
const Array p_density_multiplier = nullptr,
const Array p_size_multiplier = nullptr,
Vec4< float >  bbox = {0.f, 1.f, 0.f, 1.f} 
)

Generates a scalar field representing the signed distance to randomly generated polygons.

This function procedurally generates a field where each pixel stores the signed distance to the nearest polygon, using a smooth signed distance function (SDF) with optional clamping. The polygons are generated randomly within the specified bounding box, with configurable vertex counts, sizes, jitter, and density. Optionally, displacement noise fields can be applied to polygon positions.

Parameters
shapeResolution of the output field (width, height).
kwScaling factor for field coordinates (world units per pixel).
seedRandom seed for polygon generation.
rminMinimum polygon radius (relative to bbox size).
rmaxMaximum polygon radius (relative to bbox size).
clamping_distDistance threshold for clamping the SDF value (used to soften edges).
clamping_kSmoothness parameter for clamping.
n_vertices_minMinimum number of vertices per polygon.
n_vertices_maxMaximum number of vertices per polygon.
densityFraction of pixels covered by polygons (approximate).
jitterRandom displacement factor applied to polygon vertices.
shiftRandom position shift applied to polygon center.
p_noise_xOptional displacement noise field in the X direction (nullptr to disable).
p_noise_yOptional displacement noise field in the Y direction (nullptr to disable).
bboxBounding box in normalized coordinates {xmin, xmax, ymin, ymax}.
Returns
Array 2D array containing the signed distance field.
Note
Polygons are randomly generated per call and are not guaranteed to be convex.
The sign of the SDF is negative inside polygons and positive outside.

Example

Result

◆ polygon_field_fbm()

Array hmap::gpu::polygon_field_fbm ( Vec2< int >  shape,
Vec2< float >  kw,
uint  seed,
float  rmin = 0.05f,
float  rmax = 0.8f,
float  clamping_dist = 0.1f,
float  clamping_k = 0.1f,
int  n_vertices_min = 3,
int  n_vertices_max = 16,
float  density = 0.1f,
hmap::Vec2< float >  jitter = {0.5f, 0.5f},
float  shift = 0.1f,
int  octaves = 8,
float  persistence = 0.5f,
float  lacunarity = 2.f,
const Array p_noise_x = nullptr,
const Array p_noise_y = nullptr,
const Array p_noise_distance = nullptr,
const Array p_density_multiplier = nullptr,
const Array p_size_multiplier = nullptr,
Vec4< float >  bbox = {0.f, 1.f, 0.f, 1.f} 
)

Generates a scalar field representing the signed distance to randomly generated polygons combined with fractal Brownian motion (fBm) noise modulation.

Similar to polygon_field(), but the generated SDF is modulated by an fBm noise function to create more natural, irregular shapes. The fBm parameters allow control over the noise persistence, lacunarity, and number of octaves.

Parameters
shapeResolution of the output field (width, height).
kwScaling factor for field coordinates (world units per pixel).
seedRandom seed for polygon generation and fBm noise.
rminMinimum polygon radius (relative to bbox size).
rmaxMaximum polygon radius (relative to bbox size).
clamping_distDistance threshold for clamping the SDF value (used to soften edges).
clamping_kSmoothness parameter for clamping.
n_vertices_minMinimum number of vertices per polygon.
n_vertices_maxMaximum number of vertices per polygon.
densityFraction of pixels covered by polygons (approximate).
jitterRandom displacement factor applied to polygon vertices.
shiftRandom position shift applied to polygon center.
octavesNumber of fBm octaves.
persistenceAmplitude decay per octave in fBm.
lacunarityFrequency multiplier per octave in fBm.
p_noise_xOptional displacement noise field in the X direction (nullptr to disable).
p_noise_yOptional displacement noise field in the Y direction (nullptr to disable).
bboxBounding box in normalized coordinates {xmin, xmax, ymin, ymax}.
Returns
Array 2D array containing the signed distance field modulated by fBm.
Note
The sign of the SDF is negative inside polygons and positive outside.

Example

Result

◆ vorolines()

Array hmap::gpu::vorolines ( Vec2< int >  shape,
float  density,
uint  seed,
float  k_smoothing = 0.f,
float  exp_sigma = 0.f,
float  alpha = 0.f,
float  alpha_span = M_PI,
VoronoiReturnType  return_type = VoronoiReturnType::F1_SQUARED,
const Array p_noise_x = nullptr,
const Array p_noise_y = nullptr,
Vec4< float >  bbox = {0.f, 1.f, 0.f, 1.f},
Vec4< float >  bbox_points = {0.f, 1.f, 0.f, 1.f} 
)

Generates a Voronoi-based pattern where cells are defined by proximity to random lines.

This function generates an OpenCL-accelerated Voronoi-like pattern based on the distance from each pixel to a set of randomly oriented lines. Each line is defined by a random point and a direction sampled from a uniform distribution around a given angle.

Parameters
shapeThe resolution of the resulting 2D array (width, height).
densityNumber of base points per unit area used to define lines.
seedSeed for the random number generator used to generate base points and directions.
k_smoothingKernel smoothing factor; controls how sharp or soft the distance fields are.
exp_sigmaExponential smoothing parameter applied to the computed distance field.
alphaBase angle (in radians) used to orient the generated lines.
alpha_spanMaximum angular deviation from alpha; controls line orientation variability.
return_typeType of Voronoi output to return (e.g., F1, F2, edge distance, smoothed field, etc.).
p_noise_xOptional pointer to an input noise field applied to the X coordinates (can be nullptr).
p_noise_yOptional pointer to an input noise field applied to the Y coordinates (can be nullptr).
bboxBounding box in normalized coordinates (min_x, max_x, min_y, max_y) of the final array.
bbox_pointsBounding box within which random base points are sampled.
Returns
A 2D array (of type Array) containing the computed distance field based on line proximity.
Note
Each line is defined from a point (x, y) to a direction offset using angle theta = alpha + rand * alpha_span.
The OpenCL kernel "vorolines" must be defined and compiled beforehand.

Example

#include "highmap.hpp"
int main(void)
{
hmap::Vec2<int> shape = {256, 256};
int seed = 1;
float density = 8.f;
float k_smoothing = 0.005f;
float exp_sigma = 0.01f;
float alpha = 0.f;
float alpha_span = 0.5f * M_PI;
std::vector<hmap::VoronoiReturnType> types = {
std::vector<hmap::Array> zs = {};
for (auto type : types)
{
density,
seed,
k_smoothing,
exp_sigma,
alpha,
alpha_span,
type);
z = sqrt(z);
zs.push_back(z);
}
hmap::export_banner_png("ex_vorolines.png", zs, hmap::Cmap::INFERNO);
// FBM
zs = {};
for (auto type : types)
{
density,
seed,
k_smoothing,
exp_sigma,
alpha,
alpha_span,
type);
z = sqrt(z);
zs.push_back(z);
}
hmap::export_banner_png("ex_vorolines_fbm.png", zs, hmap::Cmap::INFERNO);
}
Array vorolines(Vec2< int > shape, float density, uint seed, float k_smoothing=0.f, float exp_sigma=0.f, float alpha=0.f, float alpha_span=M_PI, VoronoiReturnType return_type=VoronoiReturnType::F1_SQUARED, const Array *p_noise_x=nullptr, const Array *p_noise_y=nullptr, Vec4< float > bbox={0.f, 1.f, 0.f, 1.f}, Vec4< float > bbox_points={0.f, 1.f, 0.f, 1.f})
Generates a Voronoi-based pattern where cells are defined by proximity to random lines.
Definition primitives_gpu.cpp:602
Array vorolines_fbm(Vec2< int > shape, float density, uint seed, float k_smoothing=0.f, float exp_sigma=0.f, float alpha=0.f, float alpha_span=M_PI, VoronoiReturnType return_type=VoronoiReturnType::F1_SQUARED, int octaves=8, float weight=0.7f, float persistence=0.5f, float lacunarity=2.f, const Array *p_noise_x=nullptr, const Array *p_noise_y=nullptr, Vec4< float > bbox={0.f, 1.f, 0.f, 1.f}, Vec4< float > bbox_points={0.f, 1.f, 0.f, 1.f})
Generates a Voronoi-based pattern using distances to lines defined by random points and angles,...
Definition primitives_gpu.cpp:681
Array sqrt(const Array &array)
Return the square root of the array elements.
Definition math.cpp:389
@ EDGE_DISTANCE_SQUARED
Definition primitives.hpp:38
@ F1_SQUARED
Definition primitives.hpp:32
@ CONSTANT_F2MF1_SQUARED
Definition primitives.hpp:40
@ F2_SQUARED
Definition primitives.hpp:33
@ CONSTANT
Definition primitives.hpp:39
@ F1TF2_SQUARED
Definition primitives.hpp:34
@ F2MF1_SQUARED
Definition primitives.hpp:36
@ F1DF2_SQUARED
Definition primitives.hpp:35
@ EDGE_DISTANCE_EXP
Definition primitives.hpp:37

Result

◆ vorolines_fbm()

Array hmap::gpu::vorolines_fbm ( Vec2< int >  shape,
float  density,
uint  seed,
float  k_smoothing = 0.f,
float  exp_sigma = 0.f,
float  alpha = 0.f,
float  alpha_span = M_PI,
VoronoiReturnType  return_type = VoronoiReturnType::F1_SQUARED,
int  octaves = 8,
float  weight = 0.7f,
float  persistence = 0.5f,
float  lacunarity = 2.f,
const Array p_noise_x = nullptr,
const Array p_noise_y = nullptr,
Vec4< float >  bbox = {0.f, 1.f, 0.f, 1.f},
Vec4< float >  bbox_points = {0.f, 1.f, 0.f, 1.f} 
)

Generates a Voronoi-based pattern using distances to lines defined by random points and angles, with additional fractal Brownian motion (fBm) noise modulation.

This function extends the standard vorolines generation by introducing fBm-based warping of the coordinate space, resulting in more organic and fractal-like structures. It creates a Voronoi distance field based on proximity to oriented line segments and distorts the result using multi-octave procedural noise.

Parameters
shapeOutput resolution of the 2D array (width, height).
densityNumber of base points per unit area used to define lines.
seedSeed value for the random number generator.
k_smoothingKernel smoothing coefficient to soften distance values (e.g., for blending).
exp_sigmaSigma value for optional exponential smoothing on the final field.
alphaBase orientation angle (in radians) of lines generated from random points.
alpha_spanMaximum angle deviation from alpha, determining directional randomness of lines.
return_typeType of output to return (e.g., F1, F2, distance to edge, smoothed version).
octavesNumber of noise octaves used in the fBm modulation.
weightWeight of each octave's contribution to the total noise.
persistenceAmplitude decay factor for each successive octave (commonly 0.5–0.8).
lacunarityFrequency multiplier for each successive octave (commonly 2.0).
p_noise_xOptional pointer to an external noise field applied to X coordinates (can be nullptr).
p_noise_yOptional pointer to an external noise field applied to Y coordinates (can be nullptr).
bboxBounding box for the final image domain (min_x, max_x, min_y, max_y).
bbox_pointsBounding box from which the initial set of points are sampled.
Returns
A 2D Array representing the Voronoi-fBm field, distorted by noise and influenced by distance to random lines.
Note
This version uses internally computed fBm noise unless external fields (p_noise_x, p_noise_y) are provided.
This function requires an OpenCL kernel named "vorolines_fbm" to be compiled and accessible.

Example

#include "highmap.hpp"
int main(void)
{
hmap::Vec2<int> shape = {256, 256};
int seed = 1;
float density = 8.f;
float k_smoothing = 0.005f;
float exp_sigma = 0.01f;
float alpha = 0.f;
float alpha_span = 0.5f * M_PI;
std::vector<hmap::VoronoiReturnType> types = {
std::vector<hmap::Array> zs = {};
for (auto type : types)
{
density,
seed,
k_smoothing,
exp_sigma,
alpha,
alpha_span,
type);
z = sqrt(z);
zs.push_back(z);
}
hmap::export_banner_png("ex_vorolines.png", zs, hmap::Cmap::INFERNO);
// FBM
zs = {};
for (auto type : types)
{
density,
seed,
k_smoothing,
exp_sigma,
alpha,
alpha_span,
type);
z = sqrt(z);
zs.push_back(z);
}
hmap::export_banner_png("ex_vorolines_fbm.png", zs, hmap::Cmap::INFERNO);
}

Result

◆ voronoi()

Array hmap::gpu::voronoi ( Vec2< int >  shape,
Vec2< float >  kw,
uint  seed,
Vec2< float >  jitter = {0.5f, 0.5f},
float  k_smoothing = 0.f,
float  exp_sigma = 0.f,
VoronoiReturnType  return_type = VoronoiReturnType::F1_SQUARED,
const Array p_ctrl_param = nullptr,
const Array p_noise_x = nullptr,
const Array p_noise_y = nullptr,
Vec4< float >  bbox = {0.f, 1.f, 0.f, 1.f} 
)

Generates a Voronoi diagram in a 2D array with configurable properties.

Parameters
shapeThe dimensions of the output array as a 2D vector of integers.
kwThe frequency scale factors for the Voronoi cells, given as a 2D vector of floats.
seedA seed value for random number generation, ensuring reproducibility.
jitter(Optional) The amount of random variation in the positions of Voronoi cell sites, given as a 2D vector of floats. Defaults to {0.5f, 0.5f}.
return_type(Optional) The type of value to compute for the Voronoi diagram. Defaults to VoronoiReturnType::F1_SQUARED.
p_ctrl_param(Optional) A pointer to an Array used to control the Voronoi computation. Used here as a multiplier for the jitter. If nullptr, no control is applied.
p_noise_x(Optional) A pointer to an Array providing additional noise in the x-direction for cell positions. If nullptr, no x-noise is applied.
p_noise_y(Optional) A pointer to an Array providing additional noise in the y-direction for cell positions. If nullptr, no y-noise is applied.
bbox(Optional) The bounding box for the Voronoi computation, given as a 4D vector of floats representing {min_x, max_x, min_y, max_y}. Defaults to {0.f, 1.f, 0.f, 1.f}.
Returns
A 2D array representing the generated Voronoi diagram.
Note
Only available if OpenCL is enabled.

Example

#include "highmap.hpp"
int main(void)
{
hmap::Vec2<int> shape = {256, 256};
hmap::Vec2<float> kw = {8.f, 8.f};
int seed = 1;
hmap::Vec2<float> jitter = {1.f, 1.f};
float k_smoothing = 0.f;
float exp_sigma = 0.05f;
std::vector<hmap::VoronoiReturnType> types = {
std::vector<hmap::Array> zs = {};
for (auto type : types)
{
kw,
seed,
jitter,
k_smoothing,
exp_sigma,
type);
zs.push_back(z);
}
for (auto type : types)
{
kw,
seed,
jitter,
k_smoothing,
exp_sigma,
type);
zs.push_back(z);
}
}
Array voronoi(Vec2< int > shape, Vec2< float > kw, uint seed, Vec2< float > jitter={0.5f, 0.5f}, float k_smoothing=0.f, float exp_sigma=0.f, VoronoiReturnType return_type=VoronoiReturnType::F1_SQUARED, const Array *p_ctrl_param=nullptr, const Array *p_noise_x=nullptr, const Array *p_noise_y=nullptr, Vec4< float > bbox={0.f, 1.f, 0.f, 1.f})
Generates a Voronoi diagram in a 2D array with configurable properties.
Definition primitives_gpu.cpp:726
Array voronoi_fbm(Vec2< int > shape, Vec2< float > kw, uint seed, Vec2< float > jitter={0.5f, 0.5f}, float k_smoothing=0.f, float exp_sigma=0.f, VoronoiReturnType return_type=VoronoiReturnType::F1_SQUARED, int octaves=8, float weight=0.7f, float persistence=0.5f, float lacunarity=2.f, const Array *p_ctrl_param=nullptr, const Array *p_noise_x=nullptr, const Array *p_noise_y=nullptr, Vec4< float > bbox={0.f, 1.f, 0.f, 1.f})
Generates a Voronoi diagram in a 2D array with configurable properties.
Definition primitives_gpu.cpp:768

Result

◆ voronoi_fbm()

Array hmap::gpu::voronoi_fbm ( Vec2< int >  shape,
Vec2< float >  kw,
uint  seed,
Vec2< float >  jitter = {0.5f, 0.5f},
float  k_smoothing = 0.f,
float  exp_sigma = 0.f,
VoronoiReturnType  return_type = VoronoiReturnType::F1_SQUARED,
int  octaves = 8,
float  weight = 0.7f,
float  persistence = 0.5f,
float  lacunarity = 2.f,
const Array p_ctrl_param = nullptr,
const Array p_noise_x = nullptr,
const Array p_noise_y = nullptr,
Vec4< float >  bbox = {0.f, 1.f, 0.f, 1.f} 
)

Generates a Voronoi diagram in a 2D array with configurable properties.

Parameters
shapeThe dimensions of the output array as a 2D vector of integers.
kwThe frequency scale factors for the base Voronoi cells, given as a 2D vector of floats.
seedA seed value for random number generation, ensuring reproducibility.
jitter(Optional) The amount of random variation in the positions of Voronoi cell sites, given as a 2D vector of floats. Defaults to {0.5f, 0.5f}.
return_type(Optional) The type of value to compute for the Voronoi diagram. Defaults to VoronoiReturnType::F1_SQUARED.
octaves(Optional) The number of layers (octaves) in the fractal Brownian motion. Defaults to 8.
weight(Optional) The initial weight of the base layer in the FBM computation. Defaults to 0.7f.
persistence(Optional) The persistence factor that controls the amplitude reduction between octaves. Defaults to 0.5f.
lacunarity(Optional) The lacunarity factor that controls the frequency increase between octaves. Defaults to 2.f.
p_ctrl_param(Optional) A pointer to an Array used to control the Voronoi computation. If nullptr, no control is applied.
p_noise_x(Optional) A pointer to an Array providing additional noise in the x-direction for cell positions. If nullptr, no x-noise is applied.
p_noise_y(Optional) A pointer to an Array providing additional noise in the y-direction for cell positions. If nullptr, no y-noise is applied.
bbox(Optional) The bounding box for the Voronoi computation, given as a 4D vector of floats representing {min_x, max_x, min_y, max_y}. Defaults to {0.f, 1.f, 0.f, 1.f}.
Returns
A 2D array representing the generated Voronoi diagram.
Note
Only available if OpenCL is enabled.

Example

#include "highmap.hpp"
int main(void)
{
hmap::Vec2<int> shape = {256, 256};
hmap::Vec2<float> kw = {8.f, 8.f};
int seed = 1;
hmap::Vec2<float> jitter = {1.f, 1.f};
float k_smoothing = 0.f;
float exp_sigma = 0.05f;
std::vector<hmap::VoronoiReturnType> types = {
std::vector<hmap::Array> zs = {};
for (auto type : types)
{
kw,
seed,
jitter,
k_smoothing,
exp_sigma,
type);
zs.push_back(z);
}
for (auto type : types)
{
kw,
seed,
jitter,
k_smoothing,
exp_sigma,
type);
zs.push_back(z);
}
}

Result

◆ voronoi_edge_distance()

Array hmap::gpu::voronoi_edge_distance ( Vec2< int >  shape,
Vec2< float >  kw,
uint  seed,
Vec2< float >  jitter = {0.5f, 0.5f},
const Array p_ctrl_param = nullptr,
const Array p_noise_x = nullptr,
const Array p_noise_y = nullptr,
Vec4< float >  bbox = {0.f, 1.f, 0.f, 1.f} 
)

Computes the Voronoi edge distance.

Parameters
shapeThe shape of the grid as a 2D vector (width, height).
kwThe weights for the Voronoi kernel as a 2D vector.
seedThe random seed used for generating Voronoi points.
jitterOptional parameter for controlling jitter in Voronoi point placement (default is {0.5f, 0.5f}).
p_ctrl_paramOptional pointer to an Array specifying control parameters for Voronoi grid jitter (default is nullptr).
p_noise_x,p_noise_yReference to the input noise arrays.
bboxThe bounding box for the Voronoi diagram as {x_min, x_max, y_min, y_max} (default is {0.f, 1.f, 0.f, 1.f}).
Note
Taken from https://www.shadertoy.com/view/llG3zy
Only available if OpenCL is enabled.
The resulting Array has the same dimensions as the input shape.

◆ voronoise()

Array hmap::gpu::voronoise ( Vec2< int >  shape,
Vec2< float >  kw,
float  u_param,
float  v_param,
uint  seed,
const Array p_noise_x = nullptr,
const Array p_noise_y = nullptr,
Vec4< float >  bbox = {0.f, 1.f, 0.f, 1.f} 
)

Generates a 2D Voronoi noise array.

This function computes a Voronoi noise pattern based on the input parameters and returns it as a 2D array. The noise is calculated in the OpenCL kernel noise_voronoise, which uses a combination of hashing and smoothstep functions to generate a weighted Voronoi noise field.

Note
Taken from https://www.shadertoy.com/view/Xd23Dh
Only available if OpenCL is enabled.
Parameters
shapeThe dimensions of the 2D output array as a vector (width and height).
kwWave numbers for scaling the noise pattern, represented as a 2D vector.
u_paramA control parameter for the noise, adjusting the contribution of random offsets.
v_paramA control parameter for the noise, affecting the smoothness of the pattern.
p_ctrl_paramOptional pointer to an Array specifying control parameters for Voronoi grid jitter (default is nullptr).
p_noise_x,p_noise_yReference to the input noise arrays.
seedA seed value for random number generation, ensuring reproducibility.
Returns
An Array object containing the generated 2D Voronoi noise values.

Example

#include "highmap.hpp"
int main(void)
{
hmap::Vec2<int> shape = {256, 256};
hmap::Vec2<float> kw = {4.f, 4.f};
int seed = 1;
hmap::Array z00 = hmap::gpu::voronoise(shape, kw, 0.f, 0.f, seed);
hmap::Array z10 = hmap::gpu::voronoise(shape, kw, 1.f, 0.f, seed);
hmap::Array z01 = hmap::gpu::voronoise(shape, kw, 0.f, 1.f, seed);
hmap::Array z11 = hmap::gpu::voronoise(shape, kw, 1.f, 1.f, seed);
hmap::Array zfbm = hmap::gpu::voronoise_fbm(shape, kw, 1.f, 0.3f, seed);
zfbm.infos();
hmap::export_banner_png("ex_voronoise.png",
{z00, z10, z01, z11, zfbm},
}
void infos(std::string msg="") const
Display various information about the array.
Definition io.cpp:58
Array voronoise_fbm(Vec2< int > shape, Vec2< float > kw, float u_param, float v_param, uint seed, int octaves=8, float weight=0.7f, float persistence=0.5f, float lacunarity=2.f, const Array *p_ctrl_param=nullptr, const Array *p_noise_x=nullptr, const Array *p_noise_y=nullptr, Vec4< float > bbox={0.f, 1.f, 0.f, 1.f})
Return an array filled with coherence Voronoise.
Definition primitives_gpu.cpp:853
Array voronoise(Vec2< int > shape, Vec2< float > kw, float u_param, float v_param, uint seed, const Array *p_noise_x=nullptr, const Array *p_noise_y=nullptr, Vec4< float > bbox={0.f, 1.f, 0.f, 1.f})
Generates a 2D Voronoi noise array.
Definition primitives_gpu.cpp:818

Result

◆ voronoise_fbm()

Array hmap::gpu::voronoise_fbm ( Vec2< int >  shape,
Vec2< float >  kw,
float  u_param,
float  v_param,
uint  seed,
int  octaves = 8,
float  weight = 0.7f,
float  persistence = 0.5f,
float  lacunarity = 2.f,
const Array p_ctrl_param = nullptr,
const Array p_noise_x = nullptr,
const Array p_noise_y = nullptr,
Vec4< float >  bbox = {0.f, 1.f, 0.f, 1.f} 
)

Return an array filled with coherence Voronoise.

Parameters
shapeArray shape.
kwNoise wavenumbers {kx, ky} for each directions.
seedRandom seed number.
octavesNumber of octaves.
weigthOctave weighting.
persistenceOctave persistence.
lacunarityDefines the wavenumber ratio between each octaves.
p_ctrl_paramReference to the control parameter array (acts as a multiplier for the weight parameter).
p_noise_x,p_noise_yReference to the input noise arrays.
bboxDomain bounding box.
Returns
Array Fractal noise.
Note
Taken from https://www.shadertoy.com/view/clGyWm
Only available if OpenCL is enabled.

Example

#include "highmap.hpp"
int main(void)
{
hmap::Vec2<int> shape = {256, 256};
hmap::Vec2<float> kw = {4.f, 4.f};
int seed = 1;
hmap::Array z00 = hmap::gpu::voronoise(shape, kw, 0.f, 0.f, seed);
hmap::Array z10 = hmap::gpu::voronoise(shape, kw, 1.f, 0.f, seed);
hmap::Array z01 = hmap::gpu::voronoise(shape, kw, 0.f, 1.f, seed);
hmap::Array z11 = hmap::gpu::voronoise(shape, kw, 1.f, 1.f, seed);
hmap::Array zfbm = hmap::gpu::voronoise_fbm(shape, kw, 1.f, 0.3f, seed);
zfbm.infos();
hmap::export_banner_png("ex_voronoise.png",
{z00, z10, z01, z11, zfbm},
}

Result

◆ vororand() [1/2]

Array hmap::gpu::vororand ( Vec2< int >  shape,
float  density,
float  variability,
uint  seed,
float  k_smoothing = 0.f,
float  exp_sigma = 0.f,
VoronoiReturnType  return_type = VoronoiReturnType::F1_SQUARED,
const Array p_noise_x = nullptr,
const Array p_noise_y = nullptr,
Vec4< float >  bbox = {0.f, 1.f, 0.f, 1.f},
Vec4< float >  bbox_points = {0.f, 1.f, 0.f, 1.f} 
)

Generates a 2D Voronoi-based scalar field using OpenCL.

This function computes a Voronoi diagram or derived metric (such as F1, F2, or edge distances) on a grid of given shape. A set of random points is generated within an extended bounding box, based on the desired density and variability, to reduce edge artifacts. Optionally, per-pixel displacement can be applied through noise fields.

The result is stored in an Array object representing a 2D scalar field.

Parameters
shapeDimensions of the output array (width x height).
densityNumber of random points per unit area for Voronoi diagram.
variabilityAmount of randomness added to the point generation bounding box.
seedSeed for random number generation used in point sampling.
k_smoothingSmoothing factor used in soft minimum/maximum Voronoi distance computations.
exp_sigmaStandard deviation used in exponential falloff for edge distance computation.
return_typeType of Voronoi computation to perform (e.g., F1, F2, F2-F1, edge distance).
p_noise_xOptional pointer to a noise field applied to X coordinates of grid points.
p_noise_yOptional pointer to a noise field applied to Y coordinates of grid points.
bboxBounding box of the domain in which the field is computed: {xmin, xmax, ymin, ymax}.
bbox_pointsBounding box for point generation, usually larger than bbox to avoid edge effects.
Returns
An Array object containing the computed scalar field.
Note
  • The kernel "vororand" must be compiled and available in the OpenCL context.
  • If p_noise_x or p_noise_y are provided, they must match the shape of the output array.
  • The generated point cloud will be larger than bbox to reduce border artifacts.

Example

#include "highmap.hpp"
int main(void)
{
hmap::Vec2<int> shape = {256, 256};
float density = 8.f;
float variability = 4.f;
int seed = 1;
float k_smoothing = 0.05f;
float exp_sigma = 0.01f;
std::vector<hmap::VoronoiReturnType> types = {
std::vector<hmap::Array> zs = {};
for (auto type : types)
{
density,
variability,
seed,
k_smoothing,
exp_sigma,
type);
zs.push_back(z);
}
}
Array vororand(Vec2< int > shape, float density, float variability, uint seed, float k_smoothing=0.f, float exp_sigma=0.f, VoronoiReturnType return_type=VoronoiReturnType::F1_SQUARED, const Array *p_noise_x=nullptr, const Array *p_noise_y=nullptr, Vec4< float > bbox={0.f, 1.f, 0.f, 1.f}, Vec4< float > bbox_points={0.f, 1.f, 0.f, 1.f})
Generates a 2D Voronoi-based scalar field using OpenCL.
Definition primitives_gpu.cpp:935

Result

◆ vororand() [2/2]

Array hmap::gpu::vororand ( Vec2< int >  shape,
const std::vector< float > &  xp,
const std::vector< float > &  yp,
float  k_smoothing = 0.f,
float  exp_sigma = 0.f,
VoronoiReturnType  return_type = VoronoiReturnType::F1_SQUARED,
const Array p_noise_x = nullptr,
const Array p_noise_y = nullptr,
Vec4< float >  bbox = {0.f, 1.f, 0.f, 1.f} 
)

◆ maximum_smooth()

Array hmap::gpu::maximum_smooth ( const Array array1,
const Array array2,
float  k = 0.2f 
)

◆ minimum_smooth()

Array hmap::gpu::minimum_smooth ( const Array array1,
const Array array2,
float  k = 0.2f 
)

◆ sdf_2d_polyline()

Array hmap::gpu::sdf_2d_polyline ( const Path path,
Vec2< int >  shape,
Vec4< float >  bbox = {0.f, 1.f, 0.f, 1.f},
const Array p_noise_x = nullptr,
const Array p_noise_y = nullptr 
)

◆ sdf_2d_polyline_bezier()

Array hmap::gpu::sdf_2d_polyline_bezier ( const Path path,
Vec2< int >  shape,
Vec4< float >  bbox = {0.f, 1.f, 0.f, 1.f},
const Array p_noise_x = nullptr,
const Array p_noise_y = nullptr 
)

◆ select_valley()

Array hmap::gpu::select_valley ( const Array z,
int  ir,
bool  zero_at_borders = true,
bool  ridge_select = false 
)

◆ rotate()

void hmap::gpu::rotate ( Array array,
float  angle,
bool  zoom_in = true 
)

◆ warp()

void hmap::gpu::warp ( Array array,
const Array p_dx,
const Array p_dy 
)

See hmap::warp.

◆ helper_transform_bbox()

Vec4< float > hmap::gpu::helper_transform_bbox ( const Vec4< float > &  bbox_source,
const Vec4< float > &  bbox_target 
)

◆ helper_bind_optional_buffers()

void hmap::gpu::helper_bind_optional_buffers ( clwrapper::Run &  run,
const Array p_noise_x,
const Array p_noise_y 
)