example3_volume Program

Uses

  • program~~example3_volume~~UsesGraph program~example3_volume example3_volume module~forcad forcad program~example3_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_utils forcad_utils module~forcad->module~forcad_utils module~forcad_nurbs_curve->module~forcad_kinds module~forcad_nurbs_curve->module~forcad_utils module~forcad_nurbs_surface->module~forcad_kinds module~forcad_nurbs_surface->module~forcad_utils module~forcad_nurbs_volume->module~forcad_kinds module~forcad_nurbs_volume->module~forcad_utils module~forcad_utils->module~forcad_kinds

This program demonstrates the usage of a NURBS (Non-Uniform Rational B-Spline) volume object to create and finalize a NURBS volume. It sets up control points, weights, and knot vectors for all three dimensions, generates the volume, and exports the control points and the volume to VTK files.

Define the control points for the NURBS volume

Define weights for the control points (optional) Define knot vectors for all three dimensions Set knot vectors, control points, and weights for the NURBS volume object Wc is optional.

Deallocate local arrays

Export parameter space to a VTK file

Export the control points to a VTK file

Generate the NURBS volume with resolutions of 20, 20, and 20 in the three dimensions

Export the generated volume to a VTK file

Export the NURBS volume to an IGES file Not supported for volumes. Show the control geometry and geometry using PyVista

Print size of knot vectors Insert knots 0.25 and 0.75 in all three directions Print size of knot vectors after inserting knots Print degrees

Elevate degree by 2 in all three directions Print degrees after elevating

Print size of knot vectors Print size of knot vectors after removing knots Generate the refined NURBS volume with resolutions of 40, 40, and 40 in the three dimensions

Export refined parameter space to a VTK file

Export updated control points to a VTK file

Export the refined generated volume to a VTK file

Export the NURBS volume to an IGES file Not supported for volumes. 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 parameter space to a VTK file

Export the transformed control points to a VTK file

Export the transformed generated volume to a VTK file

Export the NURBS volume to an IGES file Not supported for volumes. Show the control geometry and geometry using PyVista

first compute and set the connectivities of volume elements

get the connectivity of the face1 of the first element get the degree of the faces Finalize the NURBS volume object


Calls

program~~example3_volume~~CallsGraph program~example3_volume example3_volume none~get_degree~3 nurbs_volume%get_degree program~example3_volume->none~get_degree~3 none~get_knot~3 nurbs_volume%get_knot program~example3_volume->none~get_knot~3 none~set~3 nurbs_volume%set program~example3_volume->none~set~3 proc~cmp_degreeface nurbs_volume%cmp_degreeFace program~example3_volume->proc~cmp_degreeface proc~cmp_elemface nurbs_volume%cmp_elemFace program~example3_volume->proc~cmp_elemface proc~cmp_elem~3 nurbs_volume%cmp_elem program~example3_volume->proc~cmp_elem~3 proc~create~3 nurbs_volume%create program~example3_volume->proc~create~3 proc~elevate_degree~3 nurbs_volume%elevate_degree program~example3_volume->proc~elevate_degree~3 proc~export_xc~3 nurbs_volume%export_Xc program~example3_volume->proc~export_xc~3 proc~export_xg~3 nurbs_volume%export_Xg program~example3_volume->proc~export_xg~3 proc~export_xth~3 nurbs_volume%export_Xth program~example3_volume->proc~export_xth~3 proc~finalize~3 nurbs_volume%finalize program~example3_volume->proc~finalize~3 proc~generate_xc~3 generate_Xc program~example3_volume->proc~generate_xc~3 proc~insert_knots~3 nurbs_volume%insert_knots program~example3_volume->proc~insert_knots~3 proc~remove_knots~3 nurbs_volume%remove_knots program~example3_volume->proc~remove_knots~3 proc~rotate_xc~3 nurbs_volume%rotate_Xc program~example3_volume->proc~rotate_xc~3 proc~rotate_xg~3 nurbs_volume%rotate_Xg program~example3_volume->proc~rotate_xg~3 proc~set_elem~3 nurbs_volume%set_elem program~example3_volume->proc~set_elem~3 proc~show~3 nurbs_volume%show program~example3_volume->proc~show~3 proc~translate_xc~3 nurbs_volume%translate_Xc program~example3_volume->proc~translate_xc~3 proc~translate_xg~3 nurbs_volume%translate_Xg program~example3_volume->proc~translate_xg~3 proc~get_degree_all~2 nurbs_volume%get_degree_all none~get_degree~3->proc~get_degree_all~2 proc~get_degree_dir~2 nurbs_volume%get_degree_dir none~get_degree~3->proc~get_degree_dir~2 proc~get_knot_all~3 nurbs_volume%get_knot_all none~get_knot~3->proc~get_knot_all~3 proc~get_knoti~3 nurbs_volume%get_knoti none~get_knot~3->proc~get_knoti~3 proc~set1~3 nurbs_volume%set1 none~set~3->proc~set1~3 proc~set2~3 nurbs_volume%set2 none~set~3->proc~set2~3 proc~set3~3 nurbs_volume%set3 none~set~3->proc~set3~3 proc~set4~3 nurbs_volume%set4 none~set~3->proc~set4~3 interface~elemconn_cn elemConn_Cn proc~cmp_elem~3->interface~elemconn_cn interface~unique unique proc~cmp_elem~3->interface~unique proc~get_multiplicity~3 nurbs_volume%get_multiplicity proc~cmp_elem~3->proc~get_multiplicity~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_volume%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~cmp_elem_xc_vis~3 nurbs_volume%cmp_elem_Xc_vis proc~export_xc~3->proc~cmp_elem_xc_vis~3 proc~export_vtk_legacy export_vtk_legacy proc~export_xc~3->proc~export_vtk_legacy proc~cmp_elem_xg_vis~3 nurbs_volume%cmp_elem_Xg_vis proc~export_xg~3->proc~cmp_elem_xg_vis~3 proc~export_xg~3->proc~export_vtk_legacy proc~export_xth~3->none~set~3 proc~export_xth~3->proc~cmp_elem~3 proc~export_xth~3->interface~ndgrid proc~export_xth~3->interface~unique proc~export_xth~3->proc~export_vtk_legacy 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~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_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~compute_xg_bspline_3d compute_Xg_bspline_3d interface~compute_xg~3->proc~compute_xg_bspline_3d proc~compute_xg_bspline_3d_1point compute_Xg_bspline_3d_1point interface~compute_xg~3->proc~compute_xg_bspline_3d_1point proc~compute_xg_nurbs_3d compute_Xg_nurbs_3d interface~compute_xg~3->proc~compute_xg_nurbs_3d proc~compute_xg_nurbs_3d_1point compute_Xg_nurbs_3d_1point interface~compute_xg~3->proc~compute_xg_nurbs_3d_1point proc~cmp_elemconn_cn_l cmp_elemConn_Cn_L interface~elemconn_cn->proc~cmp_elemconn_cn_l proc~cmp_elemconn_cn_s cmp_elemConn_Cn_S interface~elemconn_cn->proc~cmp_elemconn_cn_s proc~cmp_elemconn_cn_v cmp_elemConn_Cn_V interface~elemconn_cn->proc~cmp_elemconn_cn_v 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 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 proc~get_multiplicity~3->interface~compute_multiplicity cosd cosd proc~rotation->cosd sind sind proc~rotation->sind proc~cmp_degree~3 nurbs_volume%cmp_degree proc~set1~3->proc~cmp_degree~3 proc~cmp_nc~3 nurbs_volume%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~cmp_degree~3->proc~get_multiplicity~3 proc~cmp_nc~3->interface~compute_multiplicity proc~repelem repelem proc~compute_knot_vector->proc~repelem proc~basis_bspline basis_bspline proc~compute_xg_bspline_3d->proc~basis_bspline proc~kron kron proc~compute_xg_bspline_3d->proc~kron proc~compute_xg_bspline_3d_1point->proc~basis_bspline proc~compute_xg_bspline_3d_1point->proc~kron proc~cmp_tgc_3d cmp_Tgc_3d proc~compute_xg_nurbs_3d->proc~cmp_tgc_3d proc~compute_xg_nurbs_3d_1point->proc~basis_bspline proc~compute_xg_nurbs_3d_1point->proc~kron proc~cmp_tgc_3d->proc~basis_bspline proc~cmp_tgc_3d->proc~kron

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(4)

Arrays for knot vectors in all three dimensions

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

Arrays for knot vectors in all three dimensions

real(kind=rk) :: knot3(4)

Arrays for knot vectors in all three dimensions

type(nurbs_volume) :: nurbs

Declare a NURBS volume object


Functions

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 example3_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
    real(rk) :: knot1(4), knot2(4), knot3(4) !! Arrays for knot vectors in all three dimensions

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

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

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

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

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

    !> Deallocate local arrays
    deallocate(Xc, Wc)

    !> Export parameter space to a VTK file
    call nurbs%export_Xth('vtk/nurbs_volume_Xth.vtk')

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

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

    !> Generate the NURBS volume with resolutions of 20, 20, and 20 in the three dimensions
    call nurbs%create(20, 20, 20)

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

    !> Export the NURBS volume to an IGES file
    !> Not supported for volumes.

    !-----------------------------------------------------------------------------
    ! 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_volume_Xc.vtk','vtk/nurbs_volume_Xg.vtk')

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

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

    !> Insert knots 0.25 and 0.75 in all three directions
    call nurbs%insert_knots(1, [0.25_rk, 0.75_rk], [1,1]) ! direction 1
    call nurbs%insert_knots(2, [0.25_rk, 0.75_rk], [1,1]) ! direction 2
    call nurbs%insert_knots(3, [0.25_rk, 0.75_rk], [1,1]) ! direction 3

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

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

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

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

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

    call nurbs%remove_knots(1, [0.25_rk, 0.75_rk], [1,1]) ! direction 1
    call nurbs%remove_knots(2, [0.25_rk, 0.75_rk], [1,1]) ! direction 2
    call nurbs%remove_knots(3, [0.25_rk, 0.75_rk], [1,1]) ! direction 3

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

    !> Generate the refined NURBS volume with resolutions of 40, 40, and 40 in the three dimensions
    call nurbs%create()

    !> Export refined parameter space to a VTK file
    call nurbs%export_Xth('vtk/nurbs_volume_Xth2.vtk')

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

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

    !> Export the NURBS volume to an IGES file
    !> Not supported for volumes.

    !-----------------------------------------------------------------------------
    ! 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_volume_Xc2.vtk','vtk/nurbs_volume_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 parameter space to a VTK file
    call nurbs%export_Xth('vtk/nurbs_volume_Xth3.vtk')

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

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

    !> Export the NURBS volume to an IGES file
    !> Not supported for volumes.

    !-----------------------------------------------------------------------------
    ! 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_volume_Xc3.vtk','vtk/nurbs_volume_Xg3.vtk')

    !-----------------------------------------------------------------------------
    ! Extract faces
    !-----------------------------------------------------------------------------

    !> first compute and set the connectivities of volume elements
    call nurbs%set_elem(nurbs%cmp_elem())

    !> get the connectivity of the face1 of the first element
    print*, 'Face 1 of element 1:', nurbs%cmp_elemFace(elem=1, face=1)
    print*, 'Face 2 of element 1:', nurbs%cmp_elemFace(elem=1, face=2)
    print*, 'Face 3 of element 1:', nurbs%cmp_elemFace(elem=1, face=3)
    print*, 'Face 4 of element 1:', nurbs%cmp_elemFace(elem=1, face=4)
    print*, 'Face 5 of element 1:', nurbs%cmp_elemFace(elem=1, face=5)
    print*, 'Face 6 of element 1:', nurbs%cmp_elemFace(elem=1, face=6)

    !> get the degree of the faces
    print*, 'Degree of face 1:', nurbs%cmp_degreeFace(face=1)
    print*, 'Degree of face 2:', nurbs%cmp_degreeFace(face=2)
    print*, 'Degree of face 3:', nurbs%cmp_degreeFace(face=3)
    print*, 'Degree of face 4:', nurbs%cmp_degreeFace(face=4)
    print*, 'Degree of face 5:', nurbs%cmp_degreeFace(face=5)
    print*, 'Degree of face 6:', nurbs%cmp_degreeFace(face=6)

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

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

contains

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