File: genstrain.f90

package info (click to toggle)
elkcode 5.4.24-2
  • links: PTS, VCS
  • area: main
  • in suites: buster, sid
  • size: 12,840 kB
  • sloc: f90: 48,415; fortran: 22,457; perl: 965; makefile: 384; sh: 369; python: 105; ansic: 67
file content (70 lines) | stat: -rw-r--r-- 1,714 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70

! Copyright (C) 2013 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
! This file is distributed under the terms of the GNU General Public License.
! See the file COPYING for license details.

subroutine genstrain
use modmain
implicit none
! local variables
integer i,j,k
real(8) a(3,3),b(3,3),t1
! external functions
real(8) ddot,dnrm2
external ddot,dnrm2
! set first strain equal to isotropic scaling
t1=dnrm2(9,avec,1)
strain(:,:,1)=avec(:,:)/t1
nstrain=1
do i=1,3
  do j=1,3
! set strain tensor in lattice coordinates to delta_ij
    a(:,:)=0.d0
    a(i,j)=1.d0
! symmetrise strain tensor
    call symmat(a)
! convert to mixed Cartesian-lattice coordinates
    call r3mtm(ainv,a,b)
! orthogonalise strain tensor to previous tensors
    do k=1,nstrain
      t1=-ddot(9,strain(:,:,k),1,b,1)
      call daxpy(9,t1,strain(:,:,k),1,b,1)
    end do
! compute the norm
    t1=dnrm2(9,b,1)
    if (t1.lt.epslat) cycle
! normalise tensor and store in global array
    nstrain=nstrain+1
    strain(:,:,nstrain)=b(:,:)/t1
  end do
end do
! zero small components
do k=1,nstrain
  do i=1,3
    do j=1,3
      if (abs(strain(i,j,k)).lt.epslat) strain(i,j,k)=0.d0
    end do
  end do
end do
! lattice optimisation case
if ((task.eq.2).or.(task.eq.3)) then
  if (latvopt.eq.2) then
! remove isotropic scaling when latvopt=2
    strain(:,:,1)=strain(:,:,nstrain)
    nstrain=nstrain-1
  else if (latvopt.lt.0) then
! optimise over particular strain when latvopt < 0
    i=abs(latvopt)
    if (i.gt.nstrain) then
      write(*,*)
      write(*,'("Error(genstrain): |latvopt| > nstrain : ",2I8)') i,nstrain
      write(*,*)
      stop
    end if
    strain(:,:,1)=strain(:,:,i)
    nstrain=1
  end if
end if
return
end subroutine