parent
8733dfe2e0
commit
66645c8dfc
@ -0,0 +1,116 @@
|
||||
module opt_bubble
|
||||
!This module contains the bubble option which can be used to create voids with specific pressures of certain atoms
|
||||
|
||||
use parameters
|
||||
use elements
|
||||
use box
|
||||
use opt_group
|
||||
|
||||
implicit none
|
||||
|
||||
real(kind=dp), private :: br, bt, bp, c(3)
|
||||
character(len=2), private :: species
|
||||
|
||||
public
|
||||
contains
|
||||
subroutine bubble(arg_pos)
|
||||
integer, intent(inout) :: arg_pos
|
||||
|
||||
integer :: new_type, n, j, i
|
||||
real(kind = dp) :: p(3), rand, factor, per, vol
|
||||
|
||||
print *, '------------------------------------------------------------'
|
||||
print *, 'Option Bubble'
|
||||
print *, '------------------------------------------------------------'
|
||||
!First we parse the bubble command
|
||||
call parse_bubble(arg_pos)
|
||||
|
||||
!Now we use the existing group code to delete a sphere which represents the bubble
|
||||
centroid=c
|
||||
radius = br
|
||||
type='both'
|
||||
gshape='sphere'
|
||||
group_nodes = .true.
|
||||
|
||||
call get_group
|
||||
call refine_group
|
||||
call get_group
|
||||
call delete_group
|
||||
|
||||
!Now we create a new atom type with the desired species
|
||||
call add_atom_type(species, new_type)
|
||||
|
||||
!Now we calculate the number of atoms we need to add for the desired pressure
|
||||
print *, "Creating a bubble with pressure", bp, " at temperature ", bt, " with radius ", br
|
||||
|
||||
factor = 1.0e24/6.02214e23
|
||||
if (bp <= 10.0) then
|
||||
per=factor*(3.29674113+4.51777872*bp**(-0.80473167))
|
||||
else if (bp .le. 20.0) then
|
||||
per=factor*(4.73419689-0.072919483*bp)
|
||||
else
|
||||
per=factor*(4.73419689-0.072919483*bp)
|
||||
print *, 'warning: pressure is too high'
|
||||
print *, 'equation of state is only valid for < 20 GPa'
|
||||
endif
|
||||
vol = 4.0*pi/3.0*br**3.0
|
||||
|
||||
n = vol/per
|
||||
print *, "Adding ", n, " atoms of species ", species
|
||||
|
||||
!Now add n atoms randomly within the sphere
|
||||
do i = 1, n
|
||||
do while(.true.)
|
||||
do j = 1, 3
|
||||
call random_number(rand)
|
||||
p(j) = rand*(2*br) + c(j)-br
|
||||
end do
|
||||
if (norm2(p-c) < br) exit
|
||||
end do
|
||||
|
||||
call add_atom(0, new_type, 1, p)
|
||||
end do
|
||||
|
||||
end subroutine bubble
|
||||
|
||||
subroutine parse_bubble(arg_pos)
|
||||
integer, intent(inout) :: arg_pos
|
||||
|
||||
integer :: i, arglen
|
||||
character(len=100) :: tmptxt
|
||||
|
||||
!First read in the bubble centroid
|
||||
do i = 1, 3
|
||||
arg_pos = arg_pos + 1
|
||||
call get_command_argument(arg_pos, tmptxt, arglen)
|
||||
call parse_pos(i, tmptxt, c(i))
|
||||
end do
|
||||
|
||||
!Now the bubble radius
|
||||
arg_pos = arg_pos + 1
|
||||
call get_command_argument(arg_pos, tmptxt, arglen)
|
||||
if(arglen == 0) stop "Missing bubble radius"
|
||||
read(tmptxt, *) br
|
||||
|
||||
!Now bubble pressure
|
||||
arg_pos = arg_pos + 1
|
||||
call get_command_argument(arg_pos, tmptxt, arglen)
|
||||
if(arglen == 0) stop "Missing bubble pressure"
|
||||
read(tmptxt, *) bp
|
||||
|
||||
!Now bubble temperature
|
||||
arg_pos = arg_pos + 1
|
||||
call get_command_argument(arg_pos, tmptxt, arglen)
|
||||
if(arglen == 0) stop "Missing bubble temperature"
|
||||
read(tmptxt, *) bt
|
||||
|
||||
!Now the bubble species
|
||||
arg_pos = arg_pos + 1
|
||||
call get_command_argument(arg_pos, species, arglen)
|
||||
if(arglen == 0) stop "Missing bubble species"
|
||||
|
||||
arg_pos = arg_pos +1
|
||||
return
|
||||
end subroutine parse_bubble
|
||||
end module opt_bubble
|
||||
|
Loading…
Reference in new issue