You can not select more than 25 topics Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.
CACmb/src/elements.f90

187 lines
7.1 KiB

module elements
!This module contains the elements data structures, structures needed for building regions
!and operations that are done on elements
use parameters
use subroutines
implicit none
!Data structures used to represent the CAC elements. Each index represents an element
integer,allocatable :: tag_ele(:) !Element tag (used to keep track of id's
character(len=100), allocatable :: type_ele(:) !Element type
integer, allocatable :: size_ele(:), lat_ele(:) !Element siz
real(kind=dp), allocatable :: r_node(:,:,:,:) !Nodal position array
integer :: ele_num=0 !Number of elements
!Data structure used to represent atoms
integer, allocatable :: tag_atom(:) !atom id
character(len=100), allocatable:: type_atom(:) !atom type
real(kind =dp),allocatable :: r_atom(:,:) !atom position
integer :: atom_num=0 !Number of atoms
!Variables for creating elements based on primitive cells
real(kind=dp) :: cubic_cell(3,8), fcc_cell(3,8), fcc_mat(3,3), fcc_inv(3,3)
integer :: cubic_faces(4,6)
!Below are lattice type arrays which provide information on the general form of the elements.
!We currently have a limit of 10 lattice types for simplicities sake but this can be easily increased.
integer :: max_ng_node, ng_node(10) !Max number of nodes per element and number of nodes per element for each lattice type
!These variables contain information on the basis, for simplicities sake we limit
!the user to the definition of 10 lattice types with 10 basis atoms at each lattice point.
!This can be easily increased with no change to efficiency
integer :: max_basisnum, basisnum(10)!Max basis atom number, number of basis atoms in each lattice type
character(len=2) :: basis_type(10,10) !Atom type of each basis
real(kind=dp) :: basis_pos(3,10,10) !Basis atom positions
!Simulation cell parameters
real(kind=dp) :: box_bd(6)
public
contains
subroutine lattice_init
!This subroutine just intializes variables needed for building the different finite
!element types
!First initialize the cubic cell
cubic_cell = reshape((/ 0.0_dp, 0.0_dp, 0.0_dp, &
1.0_dp, 0.0_dp, 0.0_dp, &
1.0_dp, 1.0_dp, 0.0_dp, &
0.0_dp, 1.0_dp, 0.0_dp, &
0.0_dp, 0.0_dp, 1.0_dp, &
1.0_dp, 0.0_dp, 1.0_dp, &
1.0_dp, 1.0_dp, 1.0_dp, &
0.0_dp, 1.0_dp, 1.0_dp /), &
shape(fcc_cell))
!Now we create a list containing the list of vertices needed to describe the 6 cube faces
cubic_faces(:,1) = (/ 1, 4, 8, 5 /)
cubic_faces(:,2) = (/ 2, 3, 7, 6 /)
cubic_faces(:,3) = (/ 1, 2, 6, 5 /)
cubic_faces(:,4) = (/ 3, 4, 8, 7 /)
cubic_faces(:,5) = (/ 1, 2, 3, 4 /)
cubic_faces(:,6) = (/ 5, 6, 7, 8 /)
!!Now initialize the fcc primitive cell
fcc_cell = reshape((/ 0.0_dp, 0.0_dp, 0.0_dp, &
0.5_dp, 0.5_dp, 0.0_dp, &
0.5_dp, 1.0_dp, 0.5_dp, &
0.0_dp, 0.5_dp, 0.5_dp, &
0.5_dp, 0.0_dp, 0.5_dp, &
1.0_dp, 0.5_dp, 0.5_dp, &
1.0_dp, 1.0_dp, 1.0_dp, &
0.5_dp, 0.5_dp, 1.0_dp /), &
shape(fcc_cell))
fcc_mat = reshape((/ 0.5_dp, 0.5_dp, 0.0_dp, &
0.0_dp, 0.5_dp, 0.5_dp, &
0.5_dp, 0.0_dp, 0.5_dp /), &
shape(fcc_mat))
call matrix_inverse(fcc_mat,3,fcc_inv)
max_basisnum = 0
basisnum(:) = 0
basis_pos(:,:,:) = 0.0_dp
ng_node(:) = 0
end subroutine lattice_init
subroutine cell_init(lapa,esize,ele_type, orient_mat, cell_mat)
!This subroutine uses the user provided information to transform the finite element cell to the correct
!size, orientation, and dimensions using the esize, lattice parameter, element_type, and orientation provided
!by the user
integer, intent(in) :: esize
real(kind=dp), intent(in) :: lapa, orient_mat(3,3)
character(len=100), intent(in) :: ele_type
real(kind=dp), dimension(3,max_ng_node), intent(out) :: cell_mat
select case(trim(ele_type))
case('fcc')
cell_mat(:,1:8) = lapa * ((esize-1)*matmul(orient_mat, fcc_cell))
case default
print *, "Element type ", trim(ele_type), " currently not accepted"
stop
end select
end subroutine cell_init
subroutine alloc_ele_arrays(n,m)
!This subroutine used to provide initial allocation for the atom and element arrays
integer, intent(in) :: n, m !n-size of element arrays, m-size of atom arrays
integer :: allostat
!Allocate element arrays
if(n > 0) then
allocate(tag_ele(n), type_ele(n), size_ele(n), lat_ele(n), r_node(3,max_basisnum, max_ng_node,n), &
stat=allostat)
if(allostat > 0) then
print *, "Error allocating element arrays in elements.f90 because of: ", allostat
stop
end if
end if
if(m > 0) then
!Allocate atom arrays
allocate(tag_atom(m), type_atom(m), r_atom(3,m), stat=allostat)
if(allostat > 0) then
print *, "Error allocating atom arrays in elements.f90 because of: ", allostat
stop
end if
end if
end subroutine
subroutine add_element(type, size, lat, r)
!Subroutine which adds an element to the element arrays
integer, intent(in) :: size, lat
character(len=100), intent(in) :: type
real(kind=dp), intent(in) :: r(3, max_basisnum, max_ng_node)
ele_num = ele_num + 1
tag_ele(ele_num) = ele_num
type_ele(ele_num) = type
size_ele(ele_num) = size
lat_ele(ele_num) = lat
r_node(:,:,:,ele_num) = r(:,:,:)
end subroutine add_element
subroutine add_atom(type, r)
!Subroutine which adds an atom to the atom arrays
character(len=2), intent(in) :: type
real(kind=dp), intent(in), dimension(3) :: r
atom_num = atom_num+1
tag_atom(atom_num) = atom_num
type_atom(atom_num) = type
r_atom(:,atom_num) = r(:)
end subroutine add_atom
subroutine def_ng_node(n, element_types)
!This subroutine defines the maximum node number among n element types
integer, intent(in) :: n !Number of element types
character(len=100), dimension(n) :: element_types !Array of element type strings
integer :: i
max_ng_node = 0
do i=1, n
select case(trim(adjustl(element_types(i))))
case('fcc')
ng_node(i) = 8
end select
if(ng_node(i) > max_ng_node) max_ng_node = ng_node(i)
end do
end subroutine def_ng_node
end module elements