Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(cell_type), | intent(inout) | :: | this | |||
integer, | intent(in) | :: | id | ID of a cell |
||
integer, | intent(in), | dimension(3) | :: | num_points | number of points in each dimension |
|
integer, | intent(in), | dimension(3) | :: | num_points_r | Real dimensions of arrays (in CPMD that would be +1 along each dim) |
|
real(kind=dp), | intent(in), | dimension(3) | :: | origin | Coordinates of a point of origin |
|
real(kind=dp), | intent(in), | dimension(3,3) | :: | lattice | lattice vectors |
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