example_nurbs_surface Program

Uses

  • program~~example_nurbs_surface~~UsesGraph program~example_nurbs_surface example_nurbs_surface module~forcad forcad program~example_nurbs_surface->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_nurbs_curve->module~forcad_kinds fordebug fordebug module~forcad_nurbs_curve->fordebug module~forcad_utils forcad_utils module~forcad_nurbs_curve->module~forcad_utils module~forcad_nurbs_surface->module~forcad_kinds module~forcad_nurbs_surface->fordebug module~forcad_nurbs_surface->module~forcad_utils module~forcad_nurbs_volume->module~forcad_kinds module~forcad_nurbs_volume->fordebug module~forcad_nurbs_volume->module~forcad_utils module~forcad_utils->module~forcad_kinds

This program demonstrates the usage of a NURBS surface object to create, and finalize a NURBS surface. It sets up control points and weights, generates the surface, and exports the control points and the surface to VTK files at various stages.

Define control points for the NURBS surface

Define weights for the control points

Set control points and weights for the NURBS surface object

Deallocate local arrays

Elevate the degree of the NURBS surface in the first and second directions Insert knots into the NURBS surface in the first and second directions Export control points to a VTK file

Generate the NURBS surface with a resolution of 30x30

Export the generated surface to a VTK file

Export the parameter space to a VTK file

Show the control geometry, geometry and parameters using PyVista

Finalize the NURBS surface object


Calls

program~~example_nurbs_surface~~CallsGraph program~example_nurbs_surface example_nurbs_surface none~set nurbs_surface%set program~example_nurbs_surface->none~set proc~create nurbs_surface%create program~example_nurbs_surface->proc~create proc~elevate_degree nurbs_surface%elevate_degree program~example_nurbs_surface->proc~elevate_degree proc~export_xc nurbs_surface%export_Xc program~example_nurbs_surface->proc~export_xc proc~export_xg nurbs_surface%export_Xg program~example_nurbs_surface->proc~export_xg proc~export_xth_in_xg nurbs_surface%export_Xth_in_Xg program~example_nurbs_surface->proc~export_xth_in_xg proc~finalize nurbs_surface%finalize program~example_nurbs_surface->proc~finalize proc~generate_xc~8 generate_Xc program~example_nurbs_surface->proc~generate_xc~8 proc~insert_knots nurbs_surface%insert_knots program~example_nurbs_surface->proc~insert_knots proc~show nurbs_surface%show program~example_nurbs_surface->proc~show proc~set1 nurbs_surface%set1 none~set->proc~set1 proc~set2 nurbs_surface%set2 none~set->proc~set2 proc~set3 nurbs_surface%set3 none~set->proc~set3 proc~set4 nurbs_surface%set4 none~set->proc~set4 interface~compute_xg compute_Xg proc~create->interface~compute_xg interface~ndgrid ndgrid proc~create->interface~ndgrid proc~is_rational nurbs_surface%is_rational proc~create->proc~is_rational set set proc~create->set proc~elevate_degree->none~set none~get_knot nurbs_surface%get_knot proc~elevate_degree->none~get_knot 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 proc~cmp_elem_xc_vis nurbs_surface%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~export_xc->set proc~cmp_elem_xg_vis nurbs_surface%cmp_elem_Xg_vis proc~export_xg->proc~cmp_elem_xg_vis proc~export_xg->proc~export_vtk_legacy proc~export_xg->set proc~export_xth_in_xg->interface~compute_xg interface~unique unique proc~export_xth_in_xg->interface~unique proc~export_xth_in_xg->proc~export_vtk_legacy proc~export_xth_in_xg->proc~is_rational proc~export_xth_in_xg->set interface~compute_multiplicity compute_multiplicity proc~insert_knots->interface~compute_multiplicity proc~insert_knots->none~get_knot 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 s_loc s_loc proc~insert_knots->s_loc proc~insert_knots->set 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_2d compute_Xg_bspline_2d interface~compute_xg->proc~compute_xg_bspline_2d proc~compute_xg_bspline_2d_1point compute_Xg_bspline_2d_1point interface~compute_xg->proc~compute_xg_bspline_2d_1point proc~compute_xg_nurbs_2d compute_Xg_nurbs_2d interface~compute_xg->proc~compute_xg_nurbs_2d proc~compute_xg_nurbs_2d_1point compute_Xg_nurbs_2d_1point interface~compute_xg->proc~compute_xg_nurbs_2d_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~get_knot_all nurbs_surface%get_knot_all none~get_knot->proc~get_knot_all proc~get_knoti nurbs_surface%get_knoti none~get_knot->proc~get_knoti 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 proc~set1->set proc~cmp_degree nurbs_surface%cmp_degree proc~set1->proc~cmp_degree proc~cmp_nc nurbs_surface%cmp_nc proc~set1->proc~cmp_nc proc~set2->set proc~set2->proc~cmp_nc proc~compute_knot_vector compute_knot_vector proc~set2->proc~compute_knot_vector proc~set3->set proc~set3->proc~cmp_degree proc~set4->set 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~cmp_degree->set proc~get_multiplicity nurbs_surface%get_multiplicity proc~cmp_degree->proc~get_multiplicity proc~cmp_nc->interface~compute_multiplicity proc~cmp_nc->set proc~repelem repelem proc~compute_knot_vector->proc~repelem interface~kron kron proc~compute_xg_bspline_2d->interface~kron proc~basis_bspline basis_bspline proc~compute_xg_bspline_2d->proc~basis_bspline proc~compute_xg_bspline_2d_1point->interface~kron proc~compute_xg_bspline_2d_1point->proc~basis_bspline proc~cmp_tgc_2d cmp_Tgc_2d proc~compute_xg_nurbs_2d->proc~cmp_tgc_2d proc~compute_xg_nurbs_2d_1point->interface~kron proc~compute_xg_nurbs_2d_1point->proc~basis_bspline proc~kron3 kron3 interface~kron->proc~kron3 proc~kron_t1_t1 kron_t1_t1 interface~kron->proc~kron_t1_t1 proc~kron_t1_t2 kron_t1_t2 interface~kron->proc~kron_t1_t2 proc~cmp_tgc_2d->interface~kron proc~cmp_tgc_2d->proc~basis_bspline proc~get_multiplicity->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

type(nurbs_surface) :: nurbs

Declare a NURBS surface object


Functions

pure 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 example_nurbs_surface

    use forcad, only: rk, nurbs_surface

    implicit none
    type(nurbs_surface) :: nurbs             !! Declare a NURBS surface object
    real(rk), allocatable :: Xc(:,:), Wc(:)  !! Arrays for control points and weights

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

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

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

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

    !> Deallocate local arrays
    deallocate(Xc, Wc)

    !-----------------------------------------------------------------------------
    ! Refinement
    !-----------------------------------------------------------------------------

    !> Elevate the degree of the NURBS surface in the first and second directions
    call nurbs%elevate_degree(1,3)
    call nurbs%elevate_degree(2,3)

    !> Insert knots into the NURBS surface in the first and second directions
    call nurbs%insert_knots(1,[0.5_rk], [1])
    call nurbs%insert_knots(2,[0.5_rk], [1])


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

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

    !> Generate the NURBS surface with a resolution of 30x30
    call nurbs%create(res1=30, res2=30)


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

    !> Export the parameter space to a VTK file
    call nurbs%export_Xth_in_Xg('vtk/demo_surface_Xth_in_Xg.vtk')

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

    !> Show the control geometry, geometry and parameters using PyVista
    call nurbs%show('vtk/demo_surface_Xc.vtk', 'vtk/demo_surface_Xg.vtk', 'vtk/demo_surface_Xth_in_Xg.vtk')

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

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

contains

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