example3_surface Program

Uses

  • program~~example3_surface~~UsesGraph program~example3_surface example3_surface module~forcad forcad program~example3_surface->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) surface object to create and finalize a NURBS surface. It sets up control points, weights, and knot vectors for all three dimensions, generates the surface, and exports the control points and the surface to VTK files.

Define control points for the NURBS surface

Define weights for the control points Define knot vectors for both dimensions Set knot vectors, control points, and weights for the NURBS surface object

Deallocate local arrays

Export the control points to a VTK file

Generate the NURBS surface with resolutions of 30 in both dimensions

Export the generated surface to a VTK file

Show the control geometry and geometry using PyVista

Print size of the knot vectors Insert knots 0.25, twice and 0.75, once in both directions Print size of the knot vectors after inserting knots Print the degrees

Elevate degree by 2 in both directions Print the degrees after elevating

Print size of the knot vectors Remove knots 0.25, twice and 0.75, once in both directions Print size of the knot vectors after removing knots Generate the refined NURBS surface with resolutions of 30 in both dimensions

Export updated control points to a VTK file

Export the refined generated surface 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 surface object


Calls

program~~example3_surface~~CallsGraph program~example3_surface example3_surface none~get_degree~3 nurbs_surface%get_degree program~example3_surface->none~get_degree~3 none~get_knot~3 nurbs_surface%get_knot program~example3_surface->none~get_knot~3 none~set~3 nurbs_surface%set program~example3_surface->none~set~3 proc~create~3 nurbs_surface%create program~example3_surface->proc~create~3 proc~elevate_degree~3 nurbs_surface%elevate_degree program~example3_surface->proc~elevate_degree~3 proc~export_xc~3 nurbs_surface%export_Xc program~example3_surface->proc~export_xc~3 proc~export_xg~3 nurbs_surface%export_Xg program~example3_surface->proc~export_xg~3 proc~finalize~3 nurbs_surface%finalize program~example3_surface->proc~finalize~3 proc~generate_xc~5 generate_Xc program~example3_surface->proc~generate_xc~5 proc~insert_knots~3 nurbs_surface%insert_knots program~example3_surface->proc~insert_knots~3 proc~remove_knots~3 nurbs_surface%remove_knots program~example3_surface->proc~remove_knots~3 proc~rotate_xc~3 nurbs_surface%rotate_Xc program~example3_surface->proc~rotate_xc~3 proc~rotate_xg~3 nurbs_surface%rotate_Xg program~example3_surface->proc~rotate_xg~3 proc~show~3 nurbs_surface%show program~example3_surface->proc~show~3 proc~translate_xc~3 nurbs_surface%translate_Xc program~example3_surface->proc~translate_xc~3 proc~translate_xg~3 nurbs_surface%translate_Xg program~example3_surface->proc~translate_xg~3 proc~get_degree_all~2 nurbs_surface%get_degree_all none~get_degree~3->proc~get_degree_all~2 proc~get_degree_dir~2 nurbs_surface%get_degree_dir none~get_degree~3->proc~get_degree_dir~2 proc~get_knot_all~3 nurbs_surface%get_knot_all none~get_knot~3->proc~get_knot_all~3 proc~get_knoti~3 nurbs_surface%get_knoti none~get_knot~3->proc~get_knoti~3 proc~set1~3 nurbs_surface%set1 none~set~3->proc~set1~3 proc~set2~3 nurbs_surface%set2 none~set~3->proc~set2~3 proc~set3~3 nurbs_surface%set3 none~set~3->proc~set3~3 proc~set4~3 nurbs_surface%set4 none~set~3->proc~set4~3 interface~compute_xg~3 compute_Xg proc~create~3->interface~compute_xg~3 interface~ndgrid ndgrid proc~create~3->interface~ndgrid proc~is_rational~3 nurbs_surface%is_rational proc~create~3->proc~is_rational~3 proc~elevate_degree~3->none~get_knot~3 proc~elevate_degree~3->none~set~3 proc~elevate_degree_a_5_9 elevate_degree_A_5_9 proc~elevate_degree~3->proc~elevate_degree_a_5_9 proc~elevate_degree~3->proc~is_rational~3 proc~cmp_elem_xc_vis~3 nurbs_surface%cmp_elem_Xc_vis proc~export_xc~3->proc~cmp_elem_xc_vis~3 proc~cmp_elem_xg_vis~3 nurbs_surface%cmp_elem_Xg_vis proc~export_xg~3->proc~cmp_elem_xg_vis~3 proc~insert_knots~3->none~get_knot~3 proc~insert_knots~3->none~set~3 interface~compute_multiplicity compute_multiplicity proc~insert_knots~3->interface~compute_multiplicity proc~findspan findspan proc~insert_knots~3->proc~findspan proc~insert_knot_a_5_1 insert_knot_A_5_1 proc~insert_knots~3->proc~insert_knot_a_5_1 proc~insert_knots~3->proc~is_rational~3 proc~remove_knots~3->none~get_knot~3 proc~remove_knots~3->none~set~3 proc~remove_knots~3->interface~compute_multiplicity proc~remove_knots~3->proc~findspan proc~remove_knots~3->proc~is_rational~3 proc~remove_knots_a_5_8 remove_knots_A_5_8 proc~remove_knots~3->proc~remove_knots_a_5_8 proc~rotation rotation proc~rotate_xc~3->proc~rotation proc~rotate_xg~3->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~ndgrid2 ndgrid2 interface~ndgrid->proc~ndgrid2 proc~ndgrid3 ndgrid3 interface~ndgrid->proc~ndgrid3 interface~elemconn_c0 elemConn_C0 proc~cmp_elem_xc_vis~3->interface~elemconn_c0 proc~cmp_elem_xg_vis~3->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~3 nurbs_surface%cmp_degree proc~set1~3->proc~cmp_degree~3 proc~cmp_nc~3 nurbs_surface%cmp_nc proc~set1~3->proc~cmp_nc~3 proc~set2~3->proc~cmp_nc~3 proc~compute_knot_vector compute_knot_vector proc~set2~3->proc~compute_knot_vector proc~set3~3->proc~cmp_degree~3 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~3 nurbs_surface%get_multiplicity proc~cmp_degree~3->proc~get_multiplicity~3 proc~cmp_nc~3->interface~compute_multiplicity proc~repelem repelem proc~compute_knot_vector->proc~repelem proc~get_multiplicity~3->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) :: knot1(6)

Arrays for knot vectors in both dimensions

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

Arrays for knot vectors in both dimensions

type(nurbs_surface) :: nurbs

Declare a NURBS surface object


Functions

function generate_Xc(num_rows, num_cols, peak_height) result(control_points)

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: num_rows
integer, intent(in) :: num_cols
real(kind=rk), intent(in) :: peak_height

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


Source Code

program example3_surface

    use forcad

    implicit none

    type(nurbs_surface) :: nurbs            !! Declare a NURBS surface object
    real(rk), allocatable :: Xc(:,:), Wc(:) !! Arrays for control points and weights
    real(rk) :: knot1(6), knot2(6)          !! Arrays for knot vectors in both dimensions

    !-----------------------------------------------------------------------------
    ! Setting up the NURBS surface
    !-----------------------------------------------------------------------------

    !> Define control points for the NURBS surface
    Xc = generate_Xc(3, 3, 1.0_rk)

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

    !> Define knot vectors for both dimensions
    knot1 = [0.0_rk, 0.0_rk, 0.0_rk, 1.0_rk, 1.0_rk, 1.0_rk]
    knot2 = [0.0_rk, 0.0_rk, 0.0_rk, 1.0_rk, 1.0_rk, 1.0_rk]

    !> Set knot vectors, control points, and weights for the NURBS surface object
    call nurbs%set(knot1, knot2, Xc, Wc)

    !> Deallocate local arrays
    deallocate(Xc, Wc)

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

    !-----------------------------------------------------------------------------
    ! Creating the NURBS surface
    !-----------------------------------------------------------------------------

    !> Generate the NURBS surface with resolutions of 30 in both dimensions
    call nurbs%create(30, 30)

    !> Export the generated surface to a VTK file
    call nurbs%export_Xg('vtk/nurbs_surface_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_surface_Xc.vtk','vtk/nurbs_surface_Xg.vtk')

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

    !> Print size of the knot vectors
    print*, size(nurbs%get_knot(1))
    print*, size(nurbs%get_knot(2))

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

    !> Print size of the knot vectors after inserting knots
    print*, size(nurbs%get_knot(1))
    print*, size(nurbs%get_knot(2))

    !> Print the degrees
    print*, nurbs%get_degree()

    !> Elevate degree by 2 in both directions
    call nurbs%elevate_degree(1, 2) ! direction 1
    call nurbs%elevate_degree(2, 2) ! direction 2

    !> Print the degrees after elevating
    print*, nurbs%get_degree()

    !> Print size of the knot vectors
    print*, size(nurbs%get_knot(1))
    print*, size(nurbs%get_knot(2))

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

    !> Print size of the knot vectors after removing knots
    print*, size(nurbs%get_knot(1))
    print*, size(nurbs%get_knot(2))

    !> Generate the refined NURBS surface with resolutions of 30 in both dimensions
    call nurbs%create()

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

    !> Export the refined generated surface to a VTK file
    call nurbs%export_Xg('vtk/nurbs_surface_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_surface_Xc2.vtk','vtk/nurbs_surface_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_surface_Xc3.vtk')

    !> Export the transformed generated volume to a VTK file
    call nurbs%export_Xg('vtk/nurbs_surface_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_surface_Xc3.vtk','vtk/nurbs_surface_Xg3.vtk')

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

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

contains

    !-----------------------------------------------------------------------------
    function generate_Xc(num_rows, num_cols, peak_height) result(control_points)
        integer, intent(in) :: num_rows, num_cols
        real(rk), intent(in) :: peak_height
        real(rk), allocatable :: control_points(:,:)
        integer :: i, j
        real(rk) :: x_spacing, y_spacing, x_offset, y_offset
        x_spacing = 1.0_rk / real(num_cols - 1, rk)
        y_spacing = 1.0_rk / real(num_rows - 1, rk)
        x_offset = -0.5_rk
        y_offset = -0.5_rk
        allocate(control_points(num_rows * num_cols, 3))
        do i = 1, num_rows
            do j = 1, num_cols
                control_points((i - 1) * num_cols + j, 1) = x_offset + real(j - 1, rk) * x_spacing
                control_points((i - 1) * num_cols + j, 2) = y_offset + real(i - 1, rk) * y_spacing
                control_points((i - 1) * num_cols + j, 3) = &
                    peak_height * exp(-((control_points((i - 1) * num_cols + j, 1) ** 2) &
                    + (control_points((i - 1) * num_cols + j, 2) ** 2))) + 0.5_rk * peak_height * 0.2_rk
            end do
        end do
    end function
    !-----------------------------------------------------------------------------

end program example3_surface