TABLE OF CONTENTS


ABINIT/cont35 [ Functions ]

[ Top ] [ Functions ]

NAME

 cont35

FUNCTION

 Contract symmetric rank3 tensor gxa with rank5 symmetric tensor to
 produce symmetric rank2 tensor.

COPYRIGHT

 Copyright (C) 1998-2018 ABINIT group (DCA, XG, GMR)
 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

  gxa(2,10)=rank 3 symmetric complex tensor in order
  rank5(2,21)=rank 5 complex tensor (symmetric storage)

OUTPUT

  rank2(6)=rank 2 real tensor (symmetric storage) 11 22 33 32 31 21.

NOTES

 Tensors are in "symmetric" storage mode.
 For rank 3 tensor gxa this is
     111 221 331 321 311 211 222 332 322 333;
 For rank 5 tensor rank5 this is
     11111 22111 33111 32111 31111 21111 22211 33211 32211 33311
     22221 33221 32221 33321 33331 22222 33222 32222 33322 33332 33333;
 For rank 2 tensor rank2 this is 11, 22, 33, 32, 31, 21;
 gxa and rank5 are complex; rank2 is real.
 Want $2 Re[contraction]$.
 $rank2(a,b)=2 Re[gxa(i,j,k)^"*" rank5(a,b,i,j,k)]$.

 Note that the input gxa is typically the result of
 gxa(i,j,k)=[{5 \over 2} gmet(i,l) gmet(j,m) gmet(k,n) - {3 \over 2} gmet(i,j) gmet(l,m) gmet(k,n)] gxa_old(l,m)
 where the subroutine "metcon" already includes weights in the definition
 of gxa for off-diagonal elements.

PARENTS

      nonlop_pl

CHILDREN

SOURCE

 50 #if defined HAVE_CONFIG_H
 51 #include "config.h"
 52 #endif
 53 
 54 #include "abi_common.h"
 55 
 56 
 57 subroutine cont35(gxa,rank5,rank2)
 58 
 59  use defs_basis
 60 
 61 !This section has been created automatically by the script Abilint (TD).
 62 !Do not modify the following lines by hand.
 63 #undef ABI_FUNC
 64 #define ABI_FUNC 'cont35'
 65 !End of the abilint section
 66 
 67  implicit none
 68 
 69 !Arguments ------------------------------------
 70 !arrays
 71  real(dp),intent(in) :: gxa(2,10),rank5(2,21)
 72  real(dp),intent(out) :: rank2(6)
 73 
 74 !Local variables-------------------------------
 75 !scalars
 76  integer,parameter :: im=2,re=1
 77 
 78 ! *************************************************************************
 79 
 80 !Simply write out index summations
 81 
 82 !a=1, b=1 in rank2(a,b) --> maps to index 1
 83  rank2(1)=2.0d0*(&
 84 & gxa(re, 1)*rank5(re, 1)+gxa(im, 1)*rank5(im, 1)+&
 85 & gxa(re, 7)*rank5(re, 7)+gxa(im, 7)*rank5(im, 7)+&
 86 & gxa(re,10)*rank5(re,10)+gxa(im,10)*rank5(im,10)+&
 87 & gxa(re, 2)*rank5(re, 2)+gxa(im, 2)*rank5(im, 2)+&
 88 & gxa(re, 3)*rank5(re, 3)+gxa(im, 3)*rank5(im, 3)+&
 89 & gxa(re, 5)*rank5(re, 5)+gxa(im, 5)*rank5(im, 5)+&
 90 & gxa(re, 6)*rank5(re, 6)+gxa(im, 6)*rank5(im, 6)+&
 91 & gxa(re, 8)*rank5(re, 8)+gxa(im, 8)*rank5(im, 8)+&
 92 & gxa(re, 9)*rank5(re, 9)+gxa(im, 9)*rank5(im, 9)+&
 93 & gxa(re, 4)*rank5(re, 4)+gxa(im, 4)*rank5(im, 4))
 94 
 95 
 96 !a=2, b=2 in rank2(a,b) --> maps to index 2
 97  rank2(2)=2.0d0*(&
 98 & gxa(re, 1)*rank5(re, 2)+gxa(im, 1)*rank5(im, 2)+&
 99 & gxa(re, 7)*rank5(re,16)+gxa(im, 7)*rank5(im,16)+&
100 & gxa(re,10)*rank5(re,19)+gxa(im,10)*rank5(im,19)+&
101 & gxa(re, 2)*rank5(re,11)+gxa(im, 2)*rank5(im,11)+&
102 & gxa(re, 3)*rank5(re,12)+gxa(im, 3)*rank5(im,12)+&
103 & gxa(re, 5)*rank5(re, 9)+gxa(im, 5)*rank5(im, 9)+&
104 & gxa(re, 6)*rank5(re, 7)+gxa(im, 6)*rank5(im, 7)+&
105 & gxa(re, 8)*rank5(re,17)+gxa(im, 8)*rank5(im,17)+&
106 & gxa(re, 9)*rank5(re,18)+gxa(im, 9)*rank5(im,18)+&
107 & gxa(re, 4)*rank5(re,13)+gxa(im, 4)*rank5(im,13))
108 
109 !a=3, b=3 in rank2(a,b) --> maps to index 3
110  rank2(3)=2.0d0*(&
111 & gxa(re, 1)*rank5(re, 3)+gxa(im, 1)*rank5(im, 3)+&
112 & gxa(re, 7)*rank5(re,17)+gxa(im, 7)*rank5(im,17)+&
113 & gxa(re,10)*rank5(re,21)+gxa(im,10)*rank5(im,21)+&
114 & gxa(re, 2)*rank5(re,12)+gxa(im, 2)*rank5(im,12)+&
115 & gxa(re, 3)*rank5(re,15)+gxa(im, 3)*rank5(im,15)+&
116 & gxa(re, 5)*rank5(re,10)+gxa(im, 5)*rank5(im,10)+&
117 & gxa(re, 6)*rank5(re, 8)+gxa(im, 6)*rank5(im, 8)+&
118 & gxa(re, 8)*rank5(re,20)+gxa(im, 8)*rank5(im,20)+&
119 & gxa(re, 9)*rank5(re,19)+gxa(im, 9)*rank5(im,19)+&
120 & gxa(re, 4)*rank5(re,14)+gxa(im, 4)*rank5(im,14))
121 
122 !a=3, b=2 in rank2(a,b) --> maps to index 4
123  rank2(4)=2.0d0*(&
124 & gxa(re, 1)*rank5(re, 4)+gxa(im, 1)*rank5(im, 4)+&
125 & gxa(re, 7)*rank5(re,18)+gxa(im, 7)*rank5(im,18)+&
126 & gxa(re,10)*rank5(re,20)+gxa(im,10)*rank5(im,20)+&
127 & gxa(re, 2)*rank5(re,13)+gxa(im, 2)*rank5(im,13)+&
128 & gxa(re, 3)*rank5(re,14)+gxa(im, 3)*rank5(im,14)+&
129 & gxa(re, 5)*rank5(re, 8)+gxa(im, 5)*rank5(im, 8)+&
130 & gxa(re, 6)*rank5(re, 9)+gxa(im, 6)*rank5(im, 9)+&
131 & gxa(re, 8)*rank5(re,19)+gxa(im, 8)*rank5(im,19)+&
132 & gxa(re, 9)*rank5(re,17)+gxa(im, 9)*rank5(im,17)+&
133 & gxa(re, 4)*rank5(re,12)+gxa(im, 4)*rank5(im,12))
134 
135 !a=3, b=1 in rank2(a,b) --> maps to index 5
136  rank2(5)=2.0d0*(&
137 & gxa(re, 1)*rank5(re, 5)+gxa(im, 1)*rank5(im, 5)+&
138 & gxa(re, 7)*rank5(re,13)+gxa(im, 7)*rank5(im,13)+&
139 & gxa(re,10)*rank5(re,15)+gxa(im,10)*rank5(im,15)+&
140 & gxa(re, 2)*rank5(re, 9)+gxa(im, 2)*rank5(im, 9)+&
141 & gxa(re, 3)*rank5(re,10)+gxa(im, 3)*rank5(im,10)+&
142 & gxa(re, 5)*rank5(re, 3)+gxa(im, 5)*rank5(im, 3)+&
143 & gxa(re, 6)*rank5(re, 4)+gxa(im, 6)*rank5(im, 4)+&
144 & gxa(re, 8)*rank5(re,14)+gxa(im, 8)*rank5(im,14)+&
145 & gxa(re, 9)*rank5(re,12)+gxa(im, 9)*rank5(im,12)+&
146 & gxa(re, 4)*rank5(re, 8)+gxa(im, 4)*rank5(im, 8))
147 
148 !a=2, b=1 in rank2(a,b) --> maps to index 6
149  rank2(6)=2.0d0*(&
150 & gxa(re, 1)*rank5(re, 6)+gxa(im, 1)*rank5(im, 6)+&
151 & gxa(re, 7)*rank5(re,11)+gxa(im, 7)*rank5(im,11)+&
152 & gxa(re,10)*rank5(re,14)+gxa(im,10)*rank5(im,14)+&
153 & gxa(re, 2)*rank5(re, 7)+gxa(im, 2)*rank5(im, 7)+&
154 & gxa(re, 3)*rank5(re, 8)+gxa(im, 3)*rank5(im, 8)+&
155 & gxa(re, 5)*rank5(re, 4)+gxa(im, 5)*rank5(im, 4)+&
156 & gxa(re, 6)*rank5(re, 2)+gxa(im, 6)*rank5(im, 2)+&
157 & gxa(re, 8)*rank5(re,12)+gxa(im, 8)*rank5(im,12)+&
158 & gxa(re, 9)*rank5(re,13)+gxa(im, 9)*rank5(im,13)+&
159 & gxa(re, 4)*rank5(re, 9)+gxa(im, 4)*rank5(im, 9))
160 
161 end subroutine cont35