smbx_t Derived Type

type, public :: smbx_t


Contents

Source Code


Components

TypeVisibility AttributesNameInitial
integer, public :: imcon
real(kind=rp), public, dimension(3,3):: basis
real(kind=rp), public, dimension(3,3):: dl_basis
real(kind=rp), public :: volume
logical, public :: is_deforming
logical, public :: is_aligned

Whether the basis vectors are aligned with the laboratory frame


Type-Bound Procedures

procedure, public :: init => smbx_init

  • public subroutine smbx_init(this, imcon)

    Initializes an instance of smbx_t. Can also be called to reset.

    Arguments

    Type IntentOptional AttributesName
    class(smbx_t), intent(out) :: this

    An instance of smbx_t.

    integer, intent(in) :: imcon

    Flag specifying boundary conditions on the simulation box.

    Read more…

procedure, public :: set_basis => smbx_set_basis

  • public subroutine smbx_set_basis(this, bv)

    Sets all three basis vectors.

    Arguments

    Type IntentOptional AttributesName
    class(smbx_t), intent(inout) :: this
    real(kind=rp), intent(in), dimension(3,3):: bv

procedure, public :: freeze => smbx_freeze

  • public subroutine smbx_freeze(this)

    Specifies this as non-deforming.

    Arguments

    Type IntentOptional AttributesName
    class(smbx_t), intent(inout) :: this

procedure, public :: unfreeze => smbx_unfreeze

  • public subroutine smbx_unfreeze(this)

    Specifies this as deforming.

    Arguments

    Type IntentOptional AttributesName
    class(smbx_t), intent(inout) :: this

procedure, public :: get_image => smbx_get_image

  • public subroutine smbx_get_image(this, r)

    Returns the image of r under PBC.

    Arguments

    Type IntentOptional AttributesName
    class(smbx_t), intent(in) :: this
    real(kind=rp), intent(inout), dimension(3):: r

procedure, public :: wrap_all => smbx_wrap_all

  • public subroutine smbx_wrap_all(this, coords)

    Wraps atom positions w.r.t. periodic boundary conditions.

    Arguments

    Type IntentOptional AttributesName
    class(smbx_t), intent(in) :: this
    real(kind=rp), intent(inout), dimension(:,:):: coords

procedure, public :: to_center => smbx_to_center

  • public subroutine smbx_to_center(this, coords, com)

    Adjusts atom positions such that the c.o.m. of the atoms is at the center of the box. Assumes all atoms to have the same mass and aligned axis. Optionally returns the original c.o.m.

    Arguments

    Type IntentOptional AttributesName
    class(smbx_t), intent(in) :: this
    real(kind=rp), intent(inout), dimension(:,:):: coords
    real(kind=rp), intent(out), optional dimension(3):: com

procedure, public :: to_origin => smbx_to_origin

  • public subroutine smbx_to_origin(this, coords, com)

    Adjusts atom positions such that the c.o.m. of the atoms is at the origin. Assumes all atoms to have the same mass and aligned axis. Optionally returns the original c.o.m.

    Arguments

    Type IntentOptional AttributesName
    class(smbx_t), intent(in) :: this
    real(kind=rp), intent(inout), dimension(:,:):: coords
    real(kind=rp), intent(out), optional dimension(3):: com

procedure, public :: get_rnd_points => smbx_get_rnd_points

  • public subroutine smbx_get_rnd_points(this, coords)

    Returns uniformly distributed points within the box.

    Arguments

    Type IntentOptional AttributesName
    class(smbx_t), intent(in) :: this
    real(kind=rp), intent(out), dimension(:,:):: coords

Source Code

type smbx_t
    integer :: imcon
    real(rp), dimension(3,3) :: basis
    real(rp), dimension(3,3) :: dl_basis
    real(rp) :: volume
    logical :: is_deforming
    logical :: is_aligned
        !! Whether the basis vectors are aligned with the laboratory frame

    contains
        procedure :: init => smbx_init
        procedure :: set_basis => smbx_set_basis
        procedure :: freeze => smbx_freeze
        procedure :: unfreeze => smbx_unfreeze
        procedure :: get_image => smbx_get_image
        procedure :: wrap_all => smbx_wrap_all
        procedure :: to_center => smbx_to_center
        procedure :: to_origin => smbx_to_origin
        procedure :: get_rnd_points => smbx_get_rnd_points
end type smbx_t