From 03f69c6df7942360f466b1dc7c068043ac5f314b Mon Sep 17 00:00:00 2001 From: Alex Selimov Date: Thu, 5 Dec 2019 11:03:18 -0500 Subject: [PATCH] Removed extra variables from mode_create.f90, added a new module to contain simulation box information and changed code accordingly, new grow subroutine in elements. --- src/box.f90 | 48 +++++++++++++++++++++++++++++++++++ src/elements.f90 | 61 ++++++++++++++++++++++++++++++++++++++++----- src/main.f90 | 4 ++- src/mode_create.f90 | 10 +++++--- 4 files changed, 113 insertions(+), 10 deletions(-) create mode 100644 src/box.f90 diff --git a/src/box.f90 b/src/box.f90 new file mode 100644 index 0000000..fd1586a --- /dev/null +++ b/src/box.f90 @@ -0,0 +1,48 @@ +module box + !This module contains information on the properties of the current box. + use parameters + + implicit none + + real(kind=dp) :: box_bd(6) !Global box boundaries + + !The subbox variables contain values for each subbox, being the boxes read in through some + !command. Currently only mode_merge will require sub_boxes, for mode_create it will always + !allocate to only 1 sub_box + integer :: sub_box_num = 0 + real(kind=dp), allocatable :: sub_box_ori(:,:,:) + real(kind=dp), allocatable :: sub_box_bd(:,:) + + public + contains + + subroutine box_init + !Initialize some box functions + box_bd(:) = 0.0_dp + end subroutine box_init + + subroutine alloc_sub_box(n) + !Allocate the sub_box variables + + integer, intent(in) :: n + + sub_box_num = n + allocate(sub_box_ori(3,3,n), sub_box_bd(6,n)) + end subroutine alloc_sub_box + + subroutine grow_box(temp_box_bd) + !This function takes in a temporary box boundary and adjusts the overall box boundaries + !to include it + + real(kind=dp), dimension(6), intent(in) :: temp_box_bd + + integer :: i + + do i = 1, 3 + if(temp_box_bd(2*i-1) < box_bd(2*i-1)) box_bd(2*i-1) = temp_box_bd(2*i-1) + if(temp_box_bd(2*i) > box_bd(2*i)) box_bd(2*i) = temp_box_bd(2*i) + end do + return + end subroutine grow_box + +end module box \ No newline at end of file diff --git a/src/elements.f90 b/src/elements.f90 index 3b66afb..633a50a 100644 --- a/src/elements.f90 +++ b/src/elements.f90 @@ -7,7 +7,6 @@ module elements 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 @@ -16,7 +15,7 @@ module elements integer :: node_num=0 !Total number of nodes !Data structure used to represent atoms - integer, allocatable :: tag_atom(:), type_atom(:)!atom id + integer, allocatable :: type_atom(:)!atom type real(kind =dp),allocatable :: r_atom(:,:) !atom position integer :: atom_num=0 !Number of atoms @@ -119,7 +118,7 @@ module elements !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), & + allocate(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 @@ -129,7 +128,7 @@ module elements if(m > 0) then !Allocate atom arrays - allocate(tag_atom(m), type_atom(m), r_atom(3,m), stat=allostat) + allocate(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 @@ -137,6 +136,58 @@ module elements end if end subroutine + subroutine grow_ele_arrays(n, m) + + integer, intent(in) :: n, m + + integer :: ele_size, atom_size, buffer_size + integer, allocatable :: temp_int(:) + real(kind=dp), allocatable :: temp_ele_real(:,:,:,:), temp_real(:,:) + character(len=100), allocatable :: char_temp(:) + + !The default size we grow the + buffer_size = 1024 + !Figure out the size of the atom and element arrays + ele_size = size(size_ele) + atom_size = size(type_atom) + + !Check if we need to grow the ele_size, if so grow all the variables + if ( n > size(size_ele)) then + + allocate(temp_int(n+buffer_size)) + temp_int(1:ele_size) = lat_ele + temp_int(ele_size+1:) = 0 + call move_alloc(temp_int(1:ele_size), lat_ele) + + allocate(temp_int(n+buffer_size)) + temp_int(1:ele_size) = size_ele + temp_int(ele_size+1:) = 0 + call move_alloc(temp_int(1:ele_size), size_ele) + + allocate(char_temp(n+buffer_size)) + char_temp(1:ele_size) = type_ele + call move_alloc(char_temp, type_ele) + + allocate(temp_ele_real(3, max_basisnum, max_ng_node, n+buffer_size)) + temp_ele_real(:,:,:,1:ele_size) = r_node + temp_ele_real(:,:,:,ele_size+1:) = 0.0_dp + call move_alloc(temp_ele_real, r_node) + end if + + !Now grow atom arrays if needed + if (m > atom_size) then + allocate(temp_int(m+buffer_size)) + temp_int(1:atom_size) = type_atom + temp_int(atom_size+1:) = 0 + call move_alloc(temp_int, type_atom) + + allocate(temp_real(3,m+buffer_size)) + temp_real(:,1:atom_size) = r_atom + temp_real(:, atom_size+1:) = 0.0_dp + call move_alloc(temp_real, r_atom) + 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 @@ -144,7 +195,6 @@ module elements 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 @@ -160,7 +210,6 @@ module elements 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(:) diff --git a/src/main.f90 b/src/main.f90 index 7f76b6b..c7ac008 100644 --- a/src/main.f90 +++ b/src/main.f90 @@ -20,8 +20,10 @@ program main integer :: arg_num character(len=100) :: mode + !Call initialization functions call lattice_init - + call box_init + ! Command line parsing arg_num = command_argument_count() diff --git a/src/mode_create.f90 b/src/mode_create.f90 index 2d41ae8..50b5f4c 100644 --- a/src/mode_create.f90 +++ b/src/mode_create.f90 @@ -135,6 +135,10 @@ module mode_create end if end if + !The last thing we do is setup the sub_box_boundaries + call alloc_sub_box(1) + sub_box_ori(:,:,1) = orient + sub_box_bd(:,1) = box_bd end subroutine create !This subroutine parses the command and pulls out information needed for mode_create subroutine parse_command() @@ -280,9 +284,9 @@ module mode_create integer, dimension(3,8), intent(in) :: box_in_lat !The box vertices transformed to lattice space real(kind=dp), dimension(3,3), intent(in) :: transform_matrix !The transformation matrix from lattice_space to real space !Internal variables - integer :: i, inod, bd_in_lat(6), bd_in_array(6), ix, iy, iz, numlatpoints, templatpoints, ele(3,8), rzero(3), ilat, & - type_interp(basisnum(1)*esize**3), vlat(3), temp_lat(3,8), m, n, o - real(kind=dp) :: v(3), temp_nodes(3,1,8), ele_atoms(3,esize**3), r_interp(3,basisnum(1)*esize**3) + integer :: i, inod, bd_in_lat(6), bd_in_array(6), ix, iy, iz, numlatpoints, ele(3,8), rzero(3), & + vlat(3), temp_lat(3,8), m, n, o + real(kind=dp) :: v(3), temp_nodes(3,1,8) real(kind=dp), allocatable :: resize_lat_array(:,:) logical, allocatable :: lat_points(:,:,:) logical :: node_in_bd(8)