example_nurbs_curve Program

Uses

  • program~~example_nurbs_curve~~UsesGraph program~example_nurbs_curve example_nurbs_curve module~forcad forcad program~example_nurbs_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 curve object to create, and finalize a NURBS curve. It sets up control points and weights, generates the curve, and exports the control points and the curve to VTK files at various stages.

Define control points for the NURBS curve

Define weights for the control points

Set control points and weights for the NURBS curve object

Deallocate local arrays

Export initial control points to a VTK file

Generate the NURBS curve with a resolution of 500

Export the generated curve to a VTK file

Show the control geometry and geometry using PyVista

Finalize the NURBS curve object


Calls

program~~example_nurbs_curve~~CallsGraph program~example_nurbs_curve example_nurbs_curve none~set~2 nurbs_curve%set program~example_nurbs_curve->none~set~2 proc~create~2 nurbs_curve%create program~example_nurbs_curve->proc~create~2 proc~export_xc~2 nurbs_curve%export_Xc program~example_nurbs_curve->proc~export_xc~2 proc~export_xg~2 nurbs_curve%export_Xg program~example_nurbs_curve->proc~export_xg~2 proc~finalize~2 nurbs_curve%finalize program~example_nurbs_curve->proc~finalize~2 proc~generate_xc~4 generate_Xc program~example_nurbs_curve->proc~generate_xc~4 proc~show~2 nurbs_curve%show program~example_nurbs_curve->proc~show~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~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 interface~elemconn_c0 elemConn_C0 proc~cmp_elem_xc_vis~2->interface~elemconn_c0 proc~cmp_elem_xg_vis~2->interface~elemconn_c0 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~get_multiplicity~2 nurbs_curve%get_multiplicity proc~cmp_degree~2->proc~get_multiplicity~2 proc~repelem repelem proc~compute_knot_vector->proc~repelem interface~compute_multiplicity compute_multiplicity proc~get_multiplicity~2->interface~compute_multiplicity proc~compute_multiplicity1 compute_multiplicity1 interface~compute_multiplicity->proc~compute_multiplicity1 proc~compute_multiplicity2 compute_multiplicity2 interface~compute_multiplicity->proc~compute_multiplicity2

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

type(nurbs_curve) :: nurbs

Declare a NURBS curve object


Functions

function generate_Xc(num_coils, radius, height, num_points_per_coil) result(control_points)

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: num_coils
real(kind=rk), intent(in) :: radius
real(kind=rk), intent(in) :: height
integer, intent(in) :: num_points_per_coil

Return Value real(kind=rk), allocatable, (:,:)


Source Code

program example_nurbs_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

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

    !> Define control points for the NURBS curve
    Xc = generate_Xc(5, 1.0_rk, 2.0_rk, 20)

    !> Define weights for the control points
    allocate(Wc(size(Xc, 1)), source = 1.0_rk)

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

    !> Deallocate local arrays
    deallocate(Xc, Wc)

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

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

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

    !> Export the generated curve to a VTK file
    call nurbs%export_Xg('vtk/demo_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/demo_curve_Xc.vtk','vtk/demo_curve_Xg.vtk')

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

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

contains

    !-----------------------------------------------------------------------------
    function generate_Xc(num_coils, radius, height, num_points_per_coil) result(control_points)
        integer, intent(in) :: num_coils, num_points_per_coil
        real(rk), intent(in) :: radius, height
        real(rk), allocatable :: control_points(:,:)
        integer :: coil, i
        real(rk) :: theta, coil_height
        allocate(control_points(num_coils * num_points_per_coil, 3))
        do coil = 1, num_coils
            coil_height = height * real(coil-1, rk) / real(num_coils-1, rk)
            theta = 0.0_rk
            do i = 1, num_points_per_coil
                theta = theta + 2.0_rk * acos(-1.0_rk) / real(num_points_per_coil, rk)
                control_points((coil - 1) * num_points_per_coil + i, 1) = radius * cos(theta)
                control_points((coil - 1) * num_points_per_coil + i, 2) = radius * sin(theta)
                control_points((coil - 1) * num_points_per_coil + i, 3) = coil_height
            end do
        end do
    end function
    !-----------------------------------------------------------------------------

end program example_nurbs_curve