example1_curve Program

Uses

  • program~~example1_curve~~UsesGraph program~example1_curve example1_curve module~forcad forcad program~example1_curve->module~forcad 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_utils module~forcad_nurbs_surface->module~forcad_utils module~forcad_nurbs_volume->module~forcad_utils stdlib_quadrature stdlib_quadrature module~forcad_utils->stdlib_quadrature

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 control points to a VTK file

Generate the NURBS curve with a resolution of 20

Export the generated curve to a VTK 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 updated control points to a VTK file

Export the refined generated curve to a VTK 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 the transformed control points to a VTK file

Export the transformed generated volume to a VTK 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~2 nurbs_curve%get_knot program~example1_curve->none~get_knot~2 none~set~2 nurbs_curve%set program~example1_curve->none~set~2 proc~create~2 nurbs_curve%create program~example1_curve->proc~create~2 proc~elevate_degree~2 nurbs_curve%elevate_degree program~example1_curve->proc~elevate_degree~2 proc~export_xc~2 nurbs_curve%export_Xc program~example1_curve->proc~export_xc~2 proc~export_xg~2 nurbs_curve%export_Xg program~example1_curve->proc~export_xg~2 proc~finalize~2 nurbs_curve%finalize program~example1_curve->proc~finalize~2 proc~get_degree nurbs_curve%get_degree program~example1_curve->proc~get_degree proc~insert_knots~2 nurbs_curve%insert_knots program~example1_curve->proc~insert_knots~2 proc~remove_knots~2 nurbs_curve%remove_knots program~example1_curve->proc~remove_knots~2 proc~rotate_xc~2 nurbs_curve%rotate_Xc program~example1_curve->proc~rotate_xc~2 proc~rotate_xg~2 nurbs_curve%rotate_Xg program~example1_curve->proc~rotate_xg~2 proc~show~2 nurbs_curve%show program~example1_curve->proc~show~2 proc~translate_xc~2 nurbs_curve%translate_Xc program~example1_curve->proc~translate_xc~2 proc~translate_xg~2 nurbs_curve%translate_Xg program~example1_curve->proc~translate_xg~2 proc~get_knot_all~2 nurbs_curve%get_knot_all none~get_knot~2->proc~get_knot_all~2 proc~get_knoti~2 nurbs_curve%get_knoti none~get_knot~2->proc~get_knoti~2 proc~set1a nurbs_curve%set1a none~set~2->proc~set1a proc~set1~2 nurbs_curve%set1 none~set~2->proc~set1~2 proc~set2~2 nurbs_curve%set2 none~set~2->proc~set2~2 proc~set3~2 nurbs_curve%set3 none~set~2->proc~set3~2 proc~set4~2 nurbs_curve%set4 none~set~2->proc~set4~2 interface~compute_xg~2 compute_Xg proc~create~2->interface~compute_xg~2 proc~is_rational~2 nurbs_curve%is_rational proc~create~2->proc~is_rational~2 proc~elevate_degree~2->none~set~2 proc~elevate_degree_a_5_9 elevate_degree_A_5_9 proc~elevate_degree~2->proc~elevate_degree_a_5_9 proc~elevate_degree~2->proc~is_rational~2 proc~cmp_elem_xc_vis~2 nurbs_curve%cmp_elem_Xc_vis proc~export_xc~2->proc~cmp_elem_xc_vis~2 proc~cmp_elem_xg_vis~2 nurbs_curve%cmp_elem_Xg_vis proc~export_xg~2->proc~cmp_elem_xg_vis~2 proc~insert_knots~2->none~set~2 interface~compute_multiplicity compute_multiplicity proc~insert_knots~2->interface~compute_multiplicity proc~findspan findspan proc~insert_knots~2->proc~findspan proc~insert_knot_a_5_1 insert_knot_A_5_1 proc~insert_knots~2->proc~insert_knot_a_5_1 proc~insert_knots~2->proc~is_rational~2 proc~remove_knots~2->none~set~2 proc~remove_knots~2->interface~compute_multiplicity proc~remove_knots~2->proc~findspan proc~remove_knots~2->proc~is_rational~2 proc~remove_knots_a_5_8 remove_knots_A_5_8 proc~remove_knots~2->proc~remove_knots_a_5_8 proc~rotation rotation proc~rotate_xc~2->proc~rotation proc~rotate_xg~2->proc~rotation proc~compute_multiplicity1 compute_multiplicity1 interface~compute_multiplicity->proc~compute_multiplicity1 proc~compute_multiplicity2 compute_multiplicity2 interface~compute_multiplicity->proc~compute_multiplicity2 interface~elemconn_c0 elemConn_C0 proc~cmp_elem_xc_vis~2->interface~elemconn_c0 proc~cmp_elem_xg_vis~2->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~2 nurbs_curve%cmp_degree proc~set1a->proc~cmp_degree~2 proc~set1~2->proc~cmp_degree~2 proc~compute_knot_vector compute_knot_vector proc~set2~2->proc~compute_knot_vector proc~set3~2->proc~cmp_degree~2 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~factln factln proc~bincoeff->proc~factln proc~get_multiplicity~2 nurbs_curve%get_multiplicity proc~cmp_degree~2->proc~get_multiplicity~2 proc~repelem repelem proc~compute_knot_vector->proc~repelem proc~get_multiplicity~2->interface~compute_multiplicity

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 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')


    !-----------------------------------------------------------------------------
    ! 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 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')

    !-----------------------------------------------------------------------------
    ! 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 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')

    !-----------------------------------------------------------------------------
    ! 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