TABLE OF CONTENTS


ABINIT/ddkten [ Functions ]

[ Top ] [ Functions ]

NAME

 ddkten

FUNCTION

 Compact or decompact the tensors related to the ffnl(:,1,...)
 part of the ddk operator, taking into account the direction
 of the ddk perturbation.

COPYRIGHT

 Copyright (C) 1999-2018 ABINIT group (XG)
 This file is distributed under the terms of the
 GNU General Public License, see ~abinit/COPYING
 or http://www.gnu.org/copyleft/gpl.txt .
 For the initials of contributors, see ~abinit/doc/developers/contributors.txt.

INPUTS

  compact= if 1, compact from tmpfac
  idir=direction of the ddk perturbation
  rank=0,1,2, or 3 = rank of tmpfac tensor, also angular momentum (=l)

OUTPUT

  (see side effects)

SIDE EFFECTS

 Input/Output:
  temp(2,(rank*(rank+1))/2)=compacted tensor
    for l=1, just a scalar
    for l=2, a vector
  tmpfac(2,(rank+1)*(rank+2)/2)=decompacted tensor
    for l=1, a vector
    for l=2, a symmetric matrix, stored as
     (1 . .)
     (6 2 .)
     (5 4 3)

NOTES

 For l=0, there is no contribution.

PARENTS

      nonlop_pl

CHILDREN

SOURCE

 48 #if defined HAVE_CONFIG_H
 49 #include "config.h"
 50 #endif
 51 
 52 #include "abi_common.h"
 53 
 54 
 55 subroutine ddkten(compact,idir,rank,temp,tmpfac)
 56 
 57  use defs_basis
 58  use m_errors
 59 
 60 !This section has been created automatically by the script Abilint (TD).
 61 !Do not modify the following lines by hand.
 62 #undef ABI_FUNC
 63 #define ABI_FUNC 'ddkten'
 64 !End of the abilint section
 65 
 66  implicit none
 67 
 68 !Arguments ------------------------------------
 69 !scalars
 70  integer,intent(in) :: compact,idir,rank
 71 !arrays
 72  real(dp),intent(inout) :: temp(2,(rank*(rank+1))/2)
 73  real(dp),intent(inout) :: tmpfac(2,((rank+1)*(rank+2))/2)
 74 
 75 !Local variables-------------------------------
 76 !scalars
 77  character(len=500) :: message
 78 
 79 ! *************************************************************************
 80 
 81  if(rank/=1 .and. rank/=2 .and. rank/=3)then
 82    write(message, '(a,i10,a,a,a)' )&
 83 &   'Input rank=',rank,' not allowed.',ch10,&
 84 &   'Possible values are 1,2,3 only.'
 85    MSG_BUG(message)
 86  end if
 87 
 88 !Take care of p angular momentum
 89  if(rank==1)then
 90 
 91 !  Compaction tmpfac -> temp
 92    if(compact==1)then
 93      temp(:,1)=tmpfac(:,idir)
 94 
 95 !    Decompaction temp -> tmpfac
 96    else
 97      tmpfac(:,1:3)=0.0d0
 98      tmpfac(:,idir)=temp(:,1)
 99    end if
100 
101 !  Take care of d angular momentum
102 !  rank=2 11->1 22->2 33->3 32->4 31->5 21->6
103 
104  else if(rank==2)then
105 
106 !  Compaction tmpfac -> temp
107    if(compact==1)then
108      if(idir==1)then
109 !      Count the number of non-zero derivatives with respect to k(idir)
110 !      The factor of 2 on the diagonal comes from the derivative with
111 !      respect to the first K then to the second K
112        temp(:,1)=2.0d0*tmpfac(:,1); temp(:,2)=tmpfac(:,6); temp(:,3)=tmpfac(:,5)
113      else if(idir==2)then
114        temp(:,2)=2.0d0*tmpfac(:,2); temp(:,1)=tmpfac(:,6); temp(:,3)=tmpfac(:,4)
115      else if(idir==3)then
116        temp(:,3)=2.0d0*tmpfac(:,3); temp(:,1)=tmpfac(:,5); temp(:,2)=tmpfac(:,4)
117      end if
118 !    Decompaction temp -> tmpfac
119    else
120      tmpfac(:,1:6)=0.0d0
121      tmpfac(:,idir)=2.0d0*temp(:,idir)
122      if(idir==1)then
123        tmpfac(:,5)=temp(:,3); tmpfac(:,6)=temp(:,2)
124      else if(idir==2)then
125        tmpfac(:,4)=temp(:,3); tmpfac(:,6)=temp(:,1)
126      else if(idir==3)then
127        tmpfac(:,4)=temp(:,2); tmpfac(:,5)=temp(:,1)
128      end if
129    end if
130 
131 !  Take care of f angular momentum
132  else if(rank==3)then
133 !  rank=3 111->1 221->2 331->3 321->4 311->5 211->6 222->7 332->8 322->9 333->10
134 !  rank=2 11->1 22->2 33->3 32->4 31->5 21->6
135 
136 !  Compaction tmpfac -> temp
137    if(compact==1)then
138      if(idir==1)then
139 !      Count the number of non-zero derivatives with respect to k(idir)
140        temp(:,1)=3.0d0*tmpfac(:,1)
141        temp(:,2:4)=tmpfac(:,2:4)
142        temp(:,5:6)=2.0d0*tmpfac(:,5:6)
143      else if(idir==2)then
144        temp(:,6)=2.0d0*tmpfac(:,2)
145        temp(:,4)=2.0d0*tmpfac(:,9)
146        temp(:,5)=tmpfac(:,4)
147        temp(:,1)=tmpfac(:,6)
148        temp(:,3)=tmpfac(:,8)
149        temp(:,2)=3.0d0*tmpfac(:,7)
150      else if(idir==3)then
151        temp(:,3)=3.0d0*tmpfac(:,10)
152        temp(:,5)=2.0d0*tmpfac(:,3)
153        temp(:,4)=2.0d0*tmpfac(:,8)
154        temp(:,6)=tmpfac(:,4)
155        temp(:,1)=tmpfac(:,5)
156        temp(:,2)=tmpfac(:,9)
157      end if
158 !    Decompaction temp -> tmpfac
159    else
160      tmpfac(:,1:10)=0.0d0
161      if(idir==1)then
162        tmpfac(:,1)=3.0d0*temp(:,1)
163        tmpfac(:,2:4)=temp(:,2:4)
164        tmpfac(:,5:6)=2.0d0*temp(:,5:6)
165      else if(idir==2)then
166        tmpfac(:,2)=2.0d0*temp(:,6)
167        tmpfac(:,9)=2.0d0*temp(:,4)
168        tmpfac(:,4)=temp(:,5)
169        tmpfac(:,6)=temp(:,1)
170        tmpfac(:,8)=temp(:,3)
171        tmpfac(:,7)=3.0d0*temp(:,2)
172      else if(idir==3)then
173        tmpfac(:,10)=3.0d0*temp(:,3)
174        tmpfac(:,3)=2.0d0*temp(:,5)
175        tmpfac(:,8)=2.0d0*temp(:,4)
176        tmpfac(:,4)=temp(:,6)
177        tmpfac(:,5)=temp(:,1)
178        tmpfac(:,9)=temp(:,2)
179      end if
180    end if
181 
182  end if
183 
184 end subroutine ddkten