example_nurbs_volume Program

Uses

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

Define control points for the NURBS volume

Define weights for the control points

Set control points and weights for the NURBS volume object

Deallocate local arrays

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

Generate the NURBS volume with a resolution of 15X15X15

Export the generated volume to a VTK file

Export the parameter space to a VTK file

Show the control geometry, geometry and parameters using PyVista

Finalize the NURBS volume object


Calls

program~~example_nurbs_volume~~CallsGraph program~example_nurbs_volume example_nurbs_volume none~set~2 nurbs_volume%set program~example_nurbs_volume->none~set~2 proc~create~2 nurbs_volume%create program~example_nurbs_volume->proc~create~2 proc~elevate_degree~2 nurbs_volume%elevate_degree program~example_nurbs_volume->proc~elevate_degree~2 proc~export_xc~2 nurbs_volume%export_Xc program~example_nurbs_volume->proc~export_xc~2 proc~export_xg~2 nurbs_volume%export_Xg program~example_nurbs_volume->proc~export_xg~2 proc~export_xth_in_xg~2 nurbs_volume%export_Xth_in_Xg program~example_nurbs_volume->proc~export_xth_in_xg~2 proc~finalize~2 nurbs_volume%finalize program~example_nurbs_volume->proc~finalize~2 proc~generate_xc~7 generate_Xc program~example_nurbs_volume->proc~generate_xc~7 proc~insert_knots~2 nurbs_volume%insert_knots program~example_nurbs_volume->proc~insert_knots~2 proc~show~2 nurbs_volume%show program~example_nurbs_volume->proc~show~2 proc~set1~2 nurbs_volume%set1 none~set~2->proc~set1~2 proc~set2~2 nurbs_volume%set2 none~set~2->proc~set2~2 proc~set3~2 nurbs_volume%set3 none~set~2->proc~set3~2 proc~set4~2 nurbs_volume%set4 none~set~2->proc~set4~2 interface~compute_xg~2 compute_Xg proc~create~2->interface~compute_xg~2 interface~ndgrid ndgrid proc~create~2->interface~ndgrid proc~is_rational~2 nurbs_volume%is_rational proc~create~2->proc~is_rational~2 set set proc~create~2->set proc~elevate_degree~2->none~set~2 none~get_knot~2 nurbs_volume%get_knot proc~elevate_degree~2->none~get_knot~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_volume%cmp_elem_Xc_vis proc~export_xc~2->proc~cmp_elem_xc_vis~2 proc~export_vtk_legacy export_vtk_legacy proc~export_xc~2->proc~export_vtk_legacy proc~export_xc~2->set proc~cmp_elem_xg_vis~2 nurbs_volume%cmp_elem_Xg_vis proc~export_xg~2->proc~cmp_elem_xg_vis~2 proc~export_xg~2->proc~export_vtk_legacy proc~export_xg~2->set proc~export_xth_in_xg~2->interface~compute_xg~2 interface~unique unique proc~export_xth_in_xg~2->interface~unique proc~export_xth_in_xg~2->proc~export_vtk_legacy proc~export_xth_in_xg~2->proc~is_rational~2 proc~export_xth_in_xg~2->set interface~compute_multiplicity compute_multiplicity proc~insert_knots~2->interface~compute_multiplicity proc~insert_knots~2->none~get_knot~2 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 s_loc s_loc proc~insert_knots~2->s_loc proc~insert_knots~2->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_3d compute_Xg_bspline_3d interface~compute_xg~2->proc~compute_xg_bspline_3d proc~compute_xg_bspline_3d_1point compute_Xg_bspline_3d_1point interface~compute_xg~2->proc~compute_xg_bspline_3d_1point proc~compute_xg_nurbs_3d compute_Xg_nurbs_3d interface~compute_xg~2->proc~compute_xg_nurbs_3d proc~compute_xg_nurbs_3d_1point compute_Xg_nurbs_3d_1point interface~compute_xg~2->proc~compute_xg_nurbs_3d_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~2 nurbs_volume%get_knot_all none~get_knot~2->proc~get_knot_all~2 proc~get_knoti~2 nurbs_volume%get_knoti none~get_knot~2->proc~get_knoti~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~elevate_degree_a_5_9->interface~compute_multiplicity proc~bincoeff bincoeff proc~elevate_degree_a_5_9->proc~bincoeff proc~set1~2->set proc~cmp_degree~2 nurbs_volume%cmp_degree proc~set1~2->proc~cmp_degree~2 proc~cmp_nc~2 nurbs_volume%cmp_nc proc~set1~2->proc~cmp_nc~2 proc~set2~2->set proc~set2~2->proc~cmp_nc~2 proc~compute_knot_vector compute_knot_vector proc~set2~2->proc~compute_knot_vector proc~set3~2->set proc~set3~2->proc~cmp_degree~2 proc~set4~2->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~2->set proc~get_multiplicity~2 nurbs_volume%get_multiplicity proc~cmp_degree~2->proc~get_multiplicity~2 proc~cmp_nc~2->interface~compute_multiplicity proc~cmp_nc~2->set proc~repelem repelem proc~compute_knot_vector->proc~repelem interface~kron kron proc~compute_xg_bspline_3d->interface~kron proc~basis_bspline basis_bspline proc~compute_xg_bspline_3d->proc~basis_bspline proc~compute_xg_bspline_3d_1point->interface~kron proc~compute_xg_bspline_3d_1point->proc~basis_bspline proc~cmp_tgc_3d cmp_Tgc_3d proc~compute_xg_nurbs_3d->proc~cmp_tgc_3d proc~compute_xg_nurbs_3d_1point->interface~kron proc~compute_xg_nurbs_3d_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_3d->interface~kron proc~cmp_tgc_3d->proc~basis_bspline 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

type(nurbs_volume) :: nurbs

Declare a NURBS volume object


Functions

pure function generate_Xc(L) result(control_points)

Arguments

Type IntentOptional Attributes Name
real(kind=rk), intent(in) :: L

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


Source Code

program example_nurbs_volume

    use forcad, only: rk, nurbs_volume

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

    !-----------------------------------------------------------------------------
    ! Setting up the NURBS volume
    !-----------------------------------------------------------------------------

    !> Define control points for the NURBS volume
    Xc = generate_Xc(1.0_rk)

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

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

    !> Deallocate local arrays
    deallocate(Xc, Wc)

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

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

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


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

    !-----------------------------------------------------------------------------
    ! Creating the NURBS volume
    !-----------------------------------------------------------------------------

    !> Generate the NURBS volume with a resolution of 15X15X15
    call nurbs%create(15, 15, 15)


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

    !> Export the parameter space to a VTK file
    call nurbs%export_Xth_in_Xg('vtk/demo_volume_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_volume_Xc.vtk', 'vtk/demo_volume_Xg.vtk', 'vtk/demo_volume_Xth_in_Xg.vtk')

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

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

contains

    !-----------------------------------------------------------------------------
    pure function generate_Xc(L) result(control_points)
        implicit none
        real(rk), intent(in) :: L
        real(rk), allocatable :: control_points(:,:)
        real(rk) :: L2
        L2 = L / 2.0_rk
        allocate(control_points(8, 3))
        control_points(1,:) = [ L2, -L2,  L2]
        control_points(2,:) = [ L2, -L2, -L2]
        control_points(3,:) = [-L2, -L2,  L2]
        control_points(4,:) = [-L2, -L2, -L2]
        control_points(5,:) = [ L2,  L2,  L2]
        control_points(6,:) = [ L2,  L2, -L2]
        control_points(7,:) = [-L2,  L2,  L2]
        control_points(8,:) = [-L2,  L2, -L2]
    end function
    !-----------------------------------------------------------------------------

end program example_nurbs_volume