example1_curve Program

Uses

  • program~~example1_curve~~UsesGraph program~example1_curve example1_curve module~forcad forcad program~example1_curve->module~forcad module~forcad_kinds forcad_kinds module~forcad->module~forcad_kinds module~forcad_nurbs_curve forcad_nurbs_curve module~forcad->module~forcad_nurbs_curve module~forcad_nurbs_surface forcad_nurbs_surface module~forcad->module~forcad_nurbs_surface module~forcad_nurbs_volume forcad_nurbs_volume module~forcad->module~forcad_nurbs_volume module~forcad_utils forcad_utils module~forcad->module~forcad_utils module~forcad_nurbs_curve->module~forcad_kinds module~forcad_nurbs_curve->module~forcad_utils module~forcad_nurbs_surface->module~forcad_kinds module~forcad_nurbs_surface->module~forcad_utils module~forcad_nurbs_volume->module~forcad_kinds module~forcad_nurbs_volume->module~forcad_utils module~forcad_utils->module~forcad_kinds

This program demonstrates the usage of a NURBS (Non-Uniform Rational B-Spline) curve object to create and finalize a NURBS curve. It sets up control points, weights, and knot vectors for all three dimensions, generates the curve, and exports the control points and the curve to VTK files.

Define control points for the NURBS curve Define weights for the control points (optional) Define knot vector

Set knot vector, control points, and weights for the NURBS curve object. Wc is optional

Deallocate local arrays

Export parameter space to a VTK file

Export control points to a VTK file

Generate the NURBS curve with a resolution of 20

Export the generated curve to a VTK file

Export the NURBS curve to an IGES file

Show the control geometry and geometry using PyVista

Print size of the knot vector

Insert knots 0.25, twice and 0.75, once

Print size of the updated knot vector

Print the degree of the curve

Elevate the degree of the curve (2 times)

Print the updated degree of the curve

Print size of the knot vector

Remove knots 0.25, twice and 0.75, once

Print size of the updated knot vector

Generate the refined curve with a resolution of 20

Export refined parameter space to a VTK file

Export updated control points to a VTK file

Export the refined generated curve to a VTK file

Export the refined NURBS curve to an IGES file

Show the control geometry and geometry using PyVista

Rotate the control points

Rotate the generated curve

Translate the control points

Translate the generated curve

Export parameter space to a VTK file

Export the transformed control points to a VTK file

Export the transformed generated volume to a VTK file

Export the transformed NURBS curve to an IGES file

Show the control geometry and geometry using PyVista

Finalize the NURBS curve object


Calls

program~~example1_curve~~CallsGraph program~example1_curve example1_curve none~get_knot nurbs_curve%get_knot program~example1_curve->none~get_knot none~set nurbs_curve%set program~example1_curve->none~set proc~create nurbs_curve%create program~example1_curve->proc~create proc~elevate_degree nurbs_curve%elevate_degree program~example1_curve->proc~elevate_degree proc~export_iges nurbs_curve%export_iges program~example1_curve->proc~export_iges proc~export_xc nurbs_curve%export_Xc program~example1_curve->proc~export_xc proc~export_xg nurbs_curve%export_Xg program~example1_curve->proc~export_xg proc~export_xth nurbs_curve%export_Xth program~example1_curve->proc~export_xth proc~finalize nurbs_curve%finalize program~example1_curve->proc~finalize proc~get_degree nurbs_curve%get_degree program~example1_curve->proc~get_degree proc~insert_knots nurbs_curve%insert_knots program~example1_curve->proc~insert_knots proc~remove_knots nurbs_curve%remove_knots program~example1_curve->proc~remove_knots proc~rotate_xc nurbs_curve%rotate_Xc program~example1_curve->proc~rotate_xc proc~rotate_xg nurbs_curve%rotate_Xg program~example1_curve->proc~rotate_xg proc~show nurbs_curve%show program~example1_curve->proc~show proc~translate_xc nurbs_curve%translate_Xc program~example1_curve->proc~translate_xc proc~translate_xg nurbs_curve%translate_Xg program~example1_curve->proc~translate_xg proc~get_knot_all nurbs_curve%get_knot_all none~get_knot->proc~get_knot_all proc~get_knoti nurbs_curve%get_knoti none~get_knot->proc~get_knoti proc~set1 nurbs_curve%set1 none~set->proc~set1 proc~set1a nurbs_curve%set1a none~set->proc~set1a proc~set2 nurbs_curve%set2 none~set->proc~set2 proc~set3 nurbs_curve%set3 none~set->proc~set3 proc~set4 nurbs_curve%set4 none~set->proc~set4 interface~compute_xg compute_Xg proc~create->interface~compute_xg proc~is_rational nurbs_curve%is_rational proc~create->proc~is_rational proc~elevate_degree->none~set proc~elevate_degree_a_5_9 elevate_degree_A_5_9 proc~elevate_degree->proc~elevate_degree_a_5_9 proc~elevate_degree->proc~is_rational append append proc~export_iges->append delete delete proc~export_iges->delete init init proc~export_iges->init makedpsections makedpsections proc~export_iges->makedpsections makegsection makegsection proc~export_iges->makegsection makessection makessection proc~export_iges->makessection proc~export_iges->proc~is_rational writeigesfile writeigesfile proc~export_iges->writeigesfile proc~cmp_elem_xc_vis nurbs_curve%cmp_elem_Xc_vis proc~export_xc->proc~cmp_elem_xc_vis proc~export_vtk_legacy export_vtk_legacy proc~export_xc->proc~export_vtk_legacy proc~cmp_elem_xg_vis nurbs_curve%cmp_elem_Xg_vis proc~export_xg->proc~cmp_elem_xg_vis proc~export_xg->proc~export_vtk_legacy proc~export_xth->none~set interface~ndgrid ndgrid proc~export_xth->interface~ndgrid interface~unique unique proc~export_xth->interface~unique proc~cmp_elem nurbs_curve%cmp_elem proc~export_xth->proc~cmp_elem proc~export_xth->proc~export_vtk_legacy proc~insert_knots->none~set interface~compute_multiplicity compute_multiplicity proc~insert_knots->interface~compute_multiplicity proc~findspan findspan proc~insert_knots->proc~findspan proc~insert_knot_a_5_1 insert_knot_A_5_1 proc~insert_knots->proc~insert_knot_a_5_1 proc~insert_knots->proc~is_rational proc~remove_knots->none~set proc~remove_knots->interface~compute_multiplicity proc~remove_knots->proc~findspan proc~remove_knots->proc~is_rational proc~remove_knots_a_5_8 remove_knots_A_5_8 proc~remove_knots->proc~remove_knots_a_5_8 proc~rotation rotation proc~rotate_xc->proc~rotation proc~rotate_xg->proc~rotation proc~compute_multiplicity1 compute_multiplicity1 interface~compute_multiplicity->proc~compute_multiplicity1 proc~compute_multiplicity2 compute_multiplicity2 interface~compute_multiplicity->proc~compute_multiplicity2 proc~compute_xg_bspline_1d compute_Xg_bspline_1d interface~compute_xg->proc~compute_xg_bspline_1d proc~compute_xg_bspline_1d_1point compute_Xg_bspline_1d_1point interface~compute_xg->proc~compute_xg_bspline_1d_1point proc~compute_xg_nurbs_1d compute_Xg_nurbs_1d interface~compute_xg->proc~compute_xg_nurbs_1d proc~compute_xg_nurbs_1d_1point compute_Xg_nurbs_1d_1point interface~compute_xg->proc~compute_xg_nurbs_1d_1point proc~ndgrid2 ndgrid2 interface~ndgrid->proc~ndgrid2 proc~ndgrid3 ndgrid3 interface~ndgrid->proc~ndgrid3 proc~unique_integer unique_integer interface~unique->proc~unique_integer proc~unique_real unique_real interface~unique->proc~unique_real proc~cmp_elem->interface~unique interface~elemconn_cn elemConn_Cn proc~cmp_elem->interface~elemconn_cn proc~get_multiplicity nurbs_curve%get_multiplicity proc~cmp_elem->proc~get_multiplicity interface~elemconn_c0 elemConn_C0 proc~cmp_elem_xc_vis->interface~elemconn_c0 proc~cmp_elem_xg_vis->interface~elemconn_c0 proc~elevate_degree_a_5_9->interface~compute_multiplicity proc~bincoeff bincoeff proc~elevate_degree_a_5_9->proc~bincoeff cosd cosd proc~rotation->cosd sind sind proc~rotation->sind proc~cmp_degree nurbs_curve%cmp_degree proc~set1->proc~cmp_degree proc~set1a->proc~cmp_degree proc~compute_knot_vector compute_knot_vector proc~set2->proc~compute_knot_vector proc~set3->proc~cmp_degree proc~cmp_elemconn_c0_l cmp_elemConn_C0_L interface~elemconn_c0->proc~cmp_elemconn_c0_l proc~cmp_elemconn_c0_s cmp_elemConn_C0_S interface~elemconn_c0->proc~cmp_elemconn_c0_s proc~cmp_elemconn_c0_v cmp_elemConn_C0_V interface~elemconn_c0->proc~cmp_elemconn_c0_v proc~cmp_elemconn_cn_l cmp_elemConn_Cn_L interface~elemconn_cn->proc~cmp_elemconn_cn_l proc~cmp_elemconn_cn_s cmp_elemConn_Cn_S interface~elemconn_cn->proc~cmp_elemconn_cn_s proc~cmp_elemconn_cn_v cmp_elemConn_Cn_V interface~elemconn_cn->proc~cmp_elemconn_cn_v proc~factln factln proc~bincoeff->proc~factln proc~cmp_degree->proc~get_multiplicity proc~repelem repelem proc~compute_knot_vector->proc~repelem proc~basis_bspline basis_bspline proc~compute_xg_bspline_1d->proc~basis_bspline proc~compute_xg_bspline_1d_1point->proc~basis_bspline proc~cmp_tgc_1d cmp_Tgc_1d proc~compute_xg_nurbs_1d->proc~cmp_tgc_1d proc~compute_xg_nurbs_1d_1point->proc~basis_bspline proc~get_multiplicity->interface~compute_multiplicity proc~cmp_tgc_1d->proc~basis_bspline

Variables

Type Attributes Name Initial
real(kind=rk), allocatable :: Wc(:)

Arrays for control points and weights

real(kind=rk), allocatable :: Xc(:,:)

Arrays for control points and weights

real(kind=rk) :: knot(6)

Array for knot vector

type(nurbs_curve) :: nurbs

Declare a NURBS curve object


Source Code

program example1_curve

    use forcad, only: rk, nurbs_curve

    implicit none
    type(nurbs_curve) :: nurbs              !! Declare a NURBS curve object
    real(rk), allocatable :: Xc(:,:), Wc(:) !! Arrays for control points and weights
    real(rk) :: knot(6)                     !! Array for knot vector

    !-----------------------------------------------------------------------------
    ! Setting up the NURBS curve
    !-----------------------------------------------------------------------------

    !> Define control points for the NURBS curve
    allocate(Xc(3, 3))
    Xc(1,:) = [0.0_rk, 0.0_rk, 0.0_rk]
    Xc(2,:) = [0.0_rk, 5.0_rk, 0.0_rk]
    Xc(3,:) = [5.0_rk, 5.0_rk, 0.0_rk]

    !> Define weights for the control points (optional)
    allocate(Wc(3))
    Wc = [1.0_rk, 2.0_rk, 0.3_rk]

    !> Define knot vector
    knot = [0.0_rk, 0.0_rk, 0.0_rk, 1.0_rk, 1.0_rk, 1.0_rk]

    !> Set knot vector, control points, and weights for the NURBS curve object.
    !> Wc is optional
    call nurbs%set(knot, Xc, Wc)

    !> Deallocate local arrays
    deallocate(Xc, Wc)

    !> Export parameter space to a VTK file
    call nurbs%export_Xth('vtk/nurbs_curve_Xth.vtk')

    !> Export control points to a VTK file
    call nurbs%export_Xc('vtk/nurbs_curve_Xc.vtk')

    !-----------------------------------------------------------------------------
    ! Creating the NURBS curve
    !-----------------------------------------------------------------------------

    !> Generate the NURBS curve with a resolution of 20
    call nurbs%create(res = 20)

    !> Export the generated curve to a VTK file
    call nurbs%export_Xg('vtk/nurbs_curve_Xg.vtk')

    !> Export the NURBS curve to an IGES file
    call nurbs%export_iges('iges/nurbs_curve.iges')

    !-----------------------------------------------------------------------------
    ! Visualization using PyVista
    ! Note: PyVista is required for visualization. Install it using `pip install pyvista`
    !-----------------------------------------------------------------------------

    !> Show the control geometry and geometry using PyVista
    call nurbs%show('vtk/nurbs_curve_Xc.vtk','vtk/nurbs_curve_Xg.vtk')

    !-----------------------------------------------------------------------------
    ! Refinements
    !-----------------------------------------------------------------------------

    !> Print size of the knot vector
    print*, size(nurbs%get_knot())

    !> Insert knots 0.25, twice and 0.75, once
    call nurbs%insert_knots([0.25_rk, 0.75_rk], [2,1])

    !> Print size of the updated knot vector
    print*, size(nurbs%get_knot())

    !> Print the degree of the curve
    print*, nurbs%get_degree()

    !> Elevate the degree of the curve (2 times)
    call nurbs%elevate_degree(2)

    !> Print the updated degree of the curve
    print*, nurbs%get_degree()

    !> Print size of the knot vector
    print*, size(nurbs%get_knot())

    !> Remove knots 0.25, twice and 0.75, once
    call nurbs%remove_knots([0.25_rk, 0.75_rk], [2,1])

    !> Print size of the updated knot vector
    print*, size(nurbs%get_knot())

    !> Generate the refined curve with a resolution of 20
    call nurbs%create()

    !> Export refined parameter space to a VTK file
    call nurbs%export_Xth('vtk/nurbs_curve_Xth2.vtk')

    !> Export updated control points to a VTK file
    call nurbs%export_Xc('vtk/nurbs_curve_Xc2.vtk')

    !> Export the refined generated curve to a VTK file
    call nurbs%export_Xg('vtk/nurbs_curve_Xg2.vtk')

    !> Export the refined NURBS curve to an IGES file
    call nurbs%export_iges('iges/nurbs_curve2.iges')

    !-----------------------------------------------------------------------------
    ! Visualization using PyVista
    ! Note: PyVista is required for visualization. Install it using `pip install pyvista`
    !-----------------------------------------------------------------------------

    !> Show the control geometry and geometry using PyVista
    call nurbs%show('vtk/nurbs_curve_Xc2.vtk','vtk/nurbs_curve_Xg2.vtk')

    !-----------------------------------------------------------------------------
    ! Transformations
    !-----------------------------------------------------------------------------

    !> Rotate the control points
    call nurbs%rotate_Xc(alpha=45.0_rk, beta=0.0_rk, theta=90.0_rk)

    !> Rotate the generated curve
    call nurbs%rotate_Xg(alpha=-45.0_rk, beta=0.0_rk, theta=-90.0_rk)

    !> Translate the control points
    call nurbs%translate_Xc([1.0_rk, 2.0_rk, -3.0_rk])

    !> Translate the generated curve
    call nurbs%translate_Xg([-1.0_rk, -2.0_rk, 3.0_rk])

    !> Export parameter space to a VTK file
    call nurbs%export_Xth('vtk/nurbs_curve_Xth3.vtk')

    !> Export the transformed control points to a VTK file
    call nurbs%export_Xc('vtk/nurbs_curve_Xc3.vtk')

    !> Export the transformed generated volume to a VTK file
    call nurbs%export_Xg('vtk/nurbs_curve_Xg3.vtk')

    !> Export the transformed NURBS curve to an IGES file
    call nurbs%export_iges('iges/nurbs_curve3.iges')

    !-----------------------------------------------------------------------------
    ! Visualization using PyVista
    ! Note: PyVista is required for visualization. Install it using `pip install pyvista`
    !-----------------------------------------------------------------------------

    !> Show the control geometry and geometry using PyVista
    call nurbs%show('vtk/nurbs_curve_Xc3.vtk','vtk/nurbs_curve_Xg3.vtk')

    !-----------------------------------------------------------------------------
    ! Finalizing
    !-----------------------------------------------------------------------------

    !> Finalize the NURBS curve object
    call nurbs%finalize()

end program example1_curve