main.f90 Source File

This File Depends On

sourcefile~~main.f90~8~~EfferentGraph sourcefile~main.f90~8 main.f90 sourcefile~magneticdisorder.f90 magneticdisorder.f90 sourcefile~magneticdisorder.f90->sourcefile~main.f90~8 sourcefile~autocell.f90 autocell.f90 sourcefile~autocell.f90->sourcefile~main.f90~8
Help

Source Code


Source Code

#include "precompilerdefinitions"
program generate_structure
!!{!src/generate_structure/manual.md!}
use konstanter, only: flyt,lo_pi,lo_sqtol
use helpers, only: open_file,tochar,lo_determ,lo_mpi_helper,lo_chop,lo_get_axis_angles,lo_seed_random_numbers
use geometryfunctions, only: lo_inscribed_sphere_in_box
use type_symmetrylist, only: lo_symlist
use type_forcemap, only: lo_forcemap
use type_crystalstructure, only: lo_crystalstructure
use type_sqs, only: lo_sqs

use options, only: lo_opts
use magneticdisorder, only: lo_magdisorder
use autocell, only: return_supercellmatrix
implicit none

type(lo_opts) :: opts
type(lo_crystalstructure) :: uc,ss,p
type(lo_symlist) :: sl
type(lo_forcemap) :: map
type(lo_magdisorder) :: mag
type(lo_sqs) :: sqs
type(lo_mpi_helper) :: mw

!real(flyt), dimension(3,3) :: tm
!integer :: i,j,u

! Set some options and parse input file
init: block
    call mw%init()
    call opts%parse()
    call lo_seed_random_numbers()
    ! read positions
    call uc%readfromfile('infile.ucposcar',verbosity=opts%verbosity)
    write(*,*) '... read unitcell'
end block init

! Create the supercell
getsupercell: block
    real(flyt), dimension(3,3) :: tm
    real(flyt), parameter :: perfect_fill=0.523598775598299_flyt
    real(flyt) :: a,b,c,al,be,gm,fillratio,r0
    integer, dimension(3,3) :: supercellmatrix
    integer :: i,j,u
    character(len=2000) :: dumstr

    ! decide how to build a supercell, a few different routines are provided.
    if ( sum(abs(opts%ndssdim)) .gt. 0.0_flyt ) then
        ! non-diagonal supercell
        call uc%build_supercell(ss,nondiagdimensions=opts%ndssdim)
        supercellmatrix=opts%ndssdim
    elseif ( opts%desired_na .gt. 0 ) then
        ! number of atoms instead of dimensions are defined
        call return_supercellmatrix(uc,opts%desired_na,supercellmatrix)
        call uc%build_supercell(ss,nondiagdimensions=supercellmatrix)
    else
        ! normal, diagonal supercell
        call uc%build_supercell(ss,opts%ssdim)
        supercellmatrix=0
        do i=1,3
            supercellmatrix(i,i)=opts%ssdim(i)
        enddo
    endif

    ! build a supercell
    write(*,*) '... built supercell'

    r0=lo_inscribed_sphere_in_box(ss%latticevectors)
    fillratio=100*(4*lo_pi*(r0**3)/(ss%volume*3.0_flyt))/perfect_fill
    call lo_get_axis_angles(ss%latticevectors,a,b,c,al,be,gm)

    write(*,*) ' Supercellmatrix:'
    do i=1,3
        write(*,*) tochar(supercellmatrix(:,i))
    enddo
    write(*,*) '    Filling ratio: ',tochar(fillratio),'% of the ideal cube'
    write(*,"(1X,'            a,b,c:',3(2X,F14.7))") a,b,c
    write(*,"(1X,' alpha,beta,gamma:',3(2X,F14.7))") al*180/lo_pi,be*180/lo_pi,gm*180/lo_pi
    write(*,*) '  number of atoms: ',lo_determ(supercellmatrix)*uc%na
    write(*,*) ''

    select case(opts%outputformat)
    case(1) ! VASP
        call ss%writetofile('outfile.ssposcar',opts%outputformat)
        write(*,*) '... wrote supercell in VASP format'
    case(2) ! Abinit
        call ss%writetofile('outfile.supercell_abinit',opts%outputformat)
        write(*,*) '... wrote supercell in Abinit format'
    case(3) ! LAMMPS
        call ss%writetofile('outfile.supercell_lammps',opts%outputformat,transformationmatrix=tm)
        write(*,*) '... wrote supercell in LAMMPS format'
        uc%latticevectors=matmul(tm,uc%latticevectors)
        ss%latticevectors=matmul(tm,ss%latticevectors)
        call uc%writetofile('outfile.uc_lammps',1)
        call ss%writetofile('outfile.ss_lammps',1)
    case(4) ! FHI-Aims
        call ss%writetofile('outfile.supercell_aims',opts%outputformat,transformationmatrix=tm)
        write(*,*) '... wrote supercell in FHI-Aims format'
    case(5) ! xyz format, for i-pi
        ! create that weird comment line that I-PI likes:
        write(*,*) 'FIXME ATOMIC UNITS IPI'
        stop
        dumstr="# CELL{H}: "
        do i=1,3
        do j=1,3
            dumstr=trim(dumstr)//" "//tochar(ss%latticevectors(j,i),ndecimals=10)
        enddo
        enddo
        dumstr=trim(dumstr)//" cell{angstrom} positions{angstrom}"
        u=open_file('out','outfile.cell_ipi.xyz')
            write(u,*) ss%na
            write(u,*) trim(dumstr)
            do i=1,ss%na
                write(u,"(2X,A4,4X,3(1X,E19.12))") trim(ss%atomic_symbol(ss%species(i))),lo_chop(ss%rcart(:,i),lo_sqtol)
            enddo
        close(u)
        ! and print the normal output stuff
        call uc%writetofile('outfile.uc_ipi',1)
        call ss%writetofile('outfile.ss_ipi',1)
        ! and the lammps file
        ss%latticevectors=ss%latticevectors
        ss%inv_latticevectors=ss%inv_latticevectors
        call ss%writetofile('outfile.supercell_lammps_ipi',3,transformationmatrix=tm)
    end select
end block getsupercell

! If I want to do something disordered, this is where I generate that.
disorderthings: block
    ! Now stop if there is no alloy stuff going on.
    if ( uc%info%alloy .or. uc%info%collmag .or. uc%info%noncollmag ) then
        ! Best to generate a symmetry table
        if ( opts%cutoff2 .lt. 0 ) opts%cutoff2=ss%maxcutoff()
        call sl%generate(uc,ss,opts%cutoff2,-1.0_flyt,-1.0_flyt,verbosity=opts%verbosity,polar=.false.,wraparound=.true.,magcutoff2=opts%cutoff2)
        call map%generate(sl,uc,ss,polarcorrectiontype=0,verbosity=opts%verbosity)
    else
        ! If neither an alloy or magnetically disordered, stop here.
        write(*,*) '... done'
        stop
    endif

    ! Generate magnetically disordered things
    if ( uc%info%collmag .or. uc%info%noncollmag ) then
        ! Get a bunch of magnetic configurations
        call mag%generate(sl,uc,ss)
        call mag%optimize(ss,opts%magnconf,opts%magnbin)
        ! Dump them
        call mag%dump_configurations(map,uc,ss,30)
    endif

    if ( uc%info%alloy ) then
        ! Generate an sqs, optimize on randomness and so on
        call sqs%generate(uc,ss,sl,opts%verbosity)
        ! Dump a configuration
        call sqs%returnstructure(ss,p,opts%verbosity)
        call p%writetofile('outfile.sqs',1)
    endif
end block disorderthings

call mw%destroy()

end program