Nodes of different colours represent the following:
Solid arrows point from a submodule to the (sub)module which it is
descended from. Dashed arrows point from a module or program unit to
modules which it uses.
Where possible, edges connecting nodes are
given different colours to make them easier to distinguish in
large graphs.
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
Export initial control points to a VTK file
Generate the NURBS surface with a resolution of 30x30
Export the generated surface to a VTK file
Show the control geometry and geometry using PyVista
Nodes of different colours represent the following:
Solid arrows point from a procedure to one which it calls. Dashed
arrows point from an interface to procedures which implement that interface.
This could include the module procedures in a generic interface or the
implementation in a submodule of an interface in a parent module.
Where possible, edges connecting nodes are
given different colours to make them easier to distinguish in
large graphs.
function generate_Xc(num_rows, num_cols, peak_height) result(control_points)
Arguments
Type
Intent
Optional
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_surfaceuse forcad,only:rk,nurbs_surfaceimplicit none type(nurbs_surface)::nurbs!! Declare a NURBS surface objectreal(rk),allocatable::Xc(:,:),Wc(:)!! Arrays for control points and weights!-----------------------------------------------------------------------------! Setting up the NURBS surface!-----------------------------------------------------------------------------!> Define control points for the NURBS surfaceXc=generate_Xc(10,10,1.5_rk)!> Define weights for the control pointsallocate(Wc(size(Xc,1)),source=1.0_rk)!> Set control points and weights for the NURBS surface objectcall nurbs%set([10,10],Xc,Wc)!> Deallocate local arraysdeallocate(Xc,Wc)!> Export initial control points to a VTK filecall nurbs%export_Xc('vtk/demo_surface_Xc.vtk')!-----------------------------------------------------------------------------! Creating the NURBS surface!-----------------------------------------------------------------------------!> Generate the NURBS surface with a resolution of 30x30call nurbs%create(res1=30,res2=30)!> Export the generated surface to a VTK filecall nurbs%export_Xg('vtk/demo_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 PyVistacall nurbs%show('vtk/demo_surface_Xc.vtk','vtk/demo_surface_Xg.vtk')!-----------------------------------------------------------------------------! Finalizing!-----------------------------------------------------------------------------!> Finalize the NURBS surface objectcall nurbs%finalize()contains!-----------------------------------------------------------------------------function generate_Xc(num_rows,num_cols,peak_height)result(control_points)integer,intent(in)::num_rows,num_colsreal(rk),intent(in)::peak_heightreal(rk),allocatable::control_points(:,:)integer::i,jreal(rk)::x_spacing,y_spacing,x_offset,y_offsetx_spacing=1.0_rk/real(num_cols-1,rk)y_spacing=1.0_rk/real(num_rows-1,rk)x_offset=-0.5_rky_offset=-0.5_rkallocate(control_points(num_rows*num_cols,3))do i=1,num_rowsdo j=1,num_colscontrol_points((i-1)*num_cols+j,1)=x_offset+real(j-1,rk)*x_spacingcontrol_points((i-1)*num_cols+j,2)=y_offset+real(i-1,rk)*y_spacingcontrol_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_rkend do end do end function!-----------------------------------------------------------------------------end program example_nurbs_surface