sll_s_bsplines_new Subroutine

public subroutine sll_s_bsplines_new(bsplines, degree, periodic, xmin, xmax, ncells, breaks)

Arguments

Type IntentOptional Attributes Name
class(sll_c_bsplines), intent(inout), allocatable :: bsplines

allocatable polymorphic object

integer, intent(in) :: degree

spline degree

logical, intent(in) :: periodic

.true. if domain is periodic

real(kind=wp), intent(in) :: xmin

x coordinate of left boundary of domain

real(kind=wp), intent(in) :: xmax

x coordinate of right boundary of domain

integer, intent(in) :: ncells

number of cells in domain (one polynomial per cell)

real(kind=wp), intent(in), optional :: breaks(:)

Calls

proc~~sll_s_bsplines_new~~CallsGraph proc~sll_s_bsplines_new sll_s_bsplines_new init init proc~sll_s_bsplines_new->init

Called by

proc~~sll_s_bsplines_new~~CalledByGraph proc~sll_s_bsplines_new sll_s_bsplines_new proc~bsplprofile_new bsplProfile_new proc~bsplprofile_new->proc~sll_s_bsplines_new proc~cubspl_new cubspl_new proc~cubspl_new->proc~sll_s_bsplines_new proc~interpolate_cubic_spline interpolate_cubic_spline proc~interpolate_cubic_spline->proc~sll_s_bsplines_new proc~sbase_init t_sBase%sBase_init proc~sbase_init->proc~sll_s_bsplines_new proc~sbase_test sBase_test proc~sbase_init->proc~sbase_test interface~t_cubspl t_cubspl interface~t_cubspl->proc~cubspl_new interface~t_rprofile_bspl t_rProfile_bspl interface~t_rprofile_bspl->proc~bsplprofile_new proc~sbase_copy t_sBase%sBase_copy proc~sbase_copy->proc~sbase_init proc~sbase_new sBase_new proc~sbase_new->proc~sbase_init proc~base_copy t_base%base_copy proc~base_copy->proc~sbase_copy proc~base_new Base_new proc~base_new->proc~sbase_new proc~readstatefilefromascii ReadStateFileFromASCII proc~readstatefilefromascii->proc~sbase_new proc~readstatefilefromascii->proc~base_new proc~sbase_test->proc~sbase_new interface~readstate ReadState interface~readstate->proc~readstatefilefromascii proc~init_base Init_Base proc~init_base->proc~base_new proc~initmhd3d t_functional_mhd3d%InitMHD3D proc~initmhd3d->proc~base_new proc~transform_sfl_init t_transform_sfl%transform_SFL_init proc~transform_sfl_init->proc~base_new proc~init Init proc~init->proc~initmhd3d proc~init_gvec_to_jorek init_gvec_to_jorek proc~init_gvec_to_jorek->interface~readstate proc~init_gvec_to_jorek->proc~init_base proc~restartfromstate RestartFromState proc~restartfromstate->interface~readstate proc~rungvec rungvec proc~rungvec->proc~initmhd3d proc~transform_sfl_new transform_sfl_new proc~transform_sfl_new->proc~transform_sfl_init program~gvec_post GVEC_POST program~gvec_post->proc~initmhd3d proc~start_rungvec start_rungvec proc~start_rungvec->proc~rungvec program~gvec GVEC program~gvec->proc~rungvec

Source Code

  subroutine sll_s_bsplines_new( &
      bsplines, &
      degree  , &
      periodic, &
      xmin    , &
      xmax    , &
      ncells  , &
      breaks  )

    class(sll_c_bsplines), allocatable, intent(inout) :: bsplines
    integer                           , intent(in   ) :: degree
    logical                           , intent(in   ) :: periodic
    real(wp)                          , intent(in   ) :: xmin
    real(wp)                          , intent(in   ) :: xmax
    integer                           , intent(in   ) :: ncells
    real(wp),                 optional, intent(in   ) :: breaks(:)

    logical :: uniform

    ! Sanity checks
    SLL_ASSERT( .not. allocated( bsplines ) )
    SLL_ASSERT( degree > 0  )
    SLL_ASSERT( ncells > 0  )
    SLL_ASSERT( xmin < xmax )

    ! Determine if B-splines are uniform based on 'breaks' optional argument
    uniform = .not. present( breaks )

    ! Non-uniform case: perform sanity checks on breakpoints
    if (.not. uniform) then
      SLL_ASSERT( size( breaks ) == ncells+1     )
      SLL_ASSERT( xmin == breaks(1)              )
      SLL_ASSERT( xmax == breaks(size( breaks )) )
    end if

    ! Allocate B-splines object to correct type
    if (uniform) then
      allocate( sll_t_bsplines_uniform     :: bsplines )
    else
      allocate( sll_t_bsplines_non_uniform :: bsplines )
    end if

    ! Initialize B-splines object with type-specific constructor
    select type( bsplines )
    type is( sll_t_bsplines_uniform )
      call bsplines % init( degree, periodic, xmin, xmax, ncells )
    type is( sll_t_bsplines_non_uniform )
      call bsplines % init( degree, periodic, breaks )
    end select

  end subroutine sll_s_bsplines_new