mimic_cells.f90 Source File


Contents

Source Code


Source Code

!
!    MiMiC: A Framework for Multiscale Modeling in Computational Chemistry
!    Copyright (C) 2015-2022  The MiMiC Contributors (see CONTRIBUTORS file for details).
!
!    This file is part of MiMiC.
!
!    MiMiC is free software: you can redistribute it and/or modify
!    it under the terms of the GNU Lesser General Public License as
!    published by the Free Software Foundation, either version 3 of
!    the License, or (at your option) any later version.
!
!    MiMiC is distributed in the hope that it will be useful, but
!    WITHOUT ANY WARRANTY; without even the implied warranty of
!    MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
!    GNU Lesser General Public License for more details.
!
!    You should have received a copy of the GNU Lesser General Public License
!    along with this program.  If not, see <http://www.gnu.org/licenses/>.
!

!> cell types
module mimic_cells

    use mimic_precision, only: dp

    implicit none

    private

    public :: cell_type

    !> simulation cell type
    type :: cell_type
        private
        integer, public :: id
        !> number of points along each direction
        integer, dimension(3), public :: num_points
        !> real number of points along each direction
        integer, dimension(3), public :: num_points_r
        !> volume of the cell
        real(dp), public :: volume
        !> coordinates of the point of origin
        real(dp), dimension(3), public :: origin
        !> lattice vectors
        real(dp), dimension(3,3), public :: lattice
        !> inverse lattice vectors
        real(dp), dimension(3,3), public :: lattice_inverse
    contains
        private
        procedure, public :: init => init_cell
    end type cell_type

contains

subroutine init_cell(this, id, num_points, num_points_r, origin, lattice)

    class(cell_type), intent(inout) :: this
    !> ID of a cell
    integer, intent(in) :: id
    !> number of points in each dimension
    integer, dimension(3), intent(in) :: num_points
    !> Real dimensions of arrays (in CPMD that would be +1 along each dim)
    integer, dimension(3), intent(in) :: num_points_r
    !> Coordinates of a point of origin
    real(dp), dimension(3), intent(in) :: origin
    !> lattice vectors
    real(dp), dimension(3,3), intent(in) :: lattice

    real(dp), dimension(3,3) :: cofactor

    this%id = id
    this%num_points = num_points
    this%num_points_r = num_points_r
    this%origin = origin
    this%lattice = lattice

    cofactor(1,1) = lattice(2,2) * lattice(3,3) - lattice(3,2) * lattice(2,3)
    cofactor(1,2) = lattice(1,3) * lattice(3,2) - lattice(3,3) * lattice(1,2)
    cofactor(1,3) = lattice(1,2) * lattice(2,3) - lattice(2,2) * lattice(1,3)
    cofactor(2,1) = lattice(2,3) * lattice(3,1) - lattice(3,3) * lattice(2,1)
    cofactor(2,2) = lattice(1,1) * lattice(3,3) - lattice(3,1) * lattice(1,3)
    cofactor(2,3) = lattice(1,3) * lattice(2,1) - lattice(2,3) * lattice(1,1)
    cofactor(3,1) = lattice(2,1) * lattice(3,2) - lattice(3,1) * lattice(2,2)
    cofactor(3,2) = lattice(1,2) * lattice(3,1) - lattice(3,2) * lattice(1,1)
    cofactor(3,3) = lattice(1,1) * lattice(2,2) - lattice(2,1) * lattice(1,2)

    this%volume = lattice(1,1) * cofactor(1,1) + &
                & lattice(2,1) * cofactor(1,2) + &
                & lattice(3,1) * cofactor(1,3)

    this%lattice_inverse = cofactor / this%volume

end subroutine init_cell

end module mimic_cells