TABLE OF CONTENTS
ABINIT/contistr01 [ Functions ]
NAME
contistr01
FUNCTION
Carries out specialized metric tensor operations needed for contraction of the 2nd strain derivative of the l=0,1,2,3 nonlocal Kleinman-Bylander pseudopotential operation. Derivatives are wrt a pair of cartesian strain components. Full advantage is taken of the full permutational symmetry of these tensors.
COPYRIGHT
Copyright (C) 1998-2018 ABINIT group (DRH) 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
istr=1,...6 specifies cartesian strain component 11,22,33,32,31,21 rank=angular momentum gm(3,3)=metric tensor (array is symmetric but stored as 3x3) gprimd(3,3)=reciprocal space dimensional primitive translations aa(2,*)=unique elements of complex right-hand tensor bb(2,*)=unique elements of complex left-hand tensor
OUTPUT
eisnl(3)=contraction for nonlocal internal strain derivative energy
NOTES
All tensors are stored in a compressed storage mode defined below; input and output conform to this scheme. When tensor elements occur repeatedly due to symmetry, the WEIGHT IS INCLUDED in the output tensor element to simplify later contractions with other tensors of the same rank and form, i.e. the next contraction is then simply a dot product over the unique elements.
PARENTS
nonlop_pl
CHILDREN
SOURCE
47 #if defined HAVE_CONFIG_H 48 #include "config.h" 49 #endif 50 51 #include "abi_common.h" 52 53 subroutine contistr01(istr,rank,gm,gprimd,eisnl,aa,bb) 54 55 use defs_basis 56 use m_profiling_abi 57 58 !This section has been created automatically by the script Abilint (TD). 59 !Do not modify the following lines by hand. 60 #undef ABI_FUNC 61 #define ABI_FUNC 'contistr01' 62 !End of the abilint section 63 64 implicit none 65 66 !Arguments ------------------------------------ 67 !scalars 68 integer,intent(in) :: istr,rank 69 !arrays 70 real(dp),intent(in) :: aa(2,((rank+1)*(rank+2))/2),bb(2,((rank+2)*(rank+3))/2) 71 real(dp),intent(in) :: gm(3,3),gprimd(3,3) 72 real(dp),intent(out) :: eisnl(3) 73 74 !Local variables------------------------------- 75 !scalars 76 integer :: ii,jj,ka,kb 77 !arrays 78 integer,parameter :: idx(12)=(/1,1,2,2,3,3,3,2,3,1,2,1/) 79 real(dp) :: dgm(3,3),tmp(2,3) 80 real(dp),allocatable :: cm(:,:,:) 81 82 ! ************************************************************************* 83 84 ABI_ALLOCATE(cm,(3,((rank+1)*(rank+2))/2,((rank+2)*(rank+3))/2)) 85 86 ka=idx(2*istr-1);kb=idx(2*istr) 87 88 do ii = 1,3 89 dgm(:,ii)=-(gprimd(ka,:)*gprimd(kb,ii)+gprimd(kb,:)*gprimd(ka,ii)) 90 end do 91 92 cm(:,:,:)=0.d0 93 ! 94 !The code below was written by a Mathematica program and formatted by 95 !a combination of editing scripts. It is not intended to be read 96 !by human beings, and certainly not to be modified by one. Conceivably 97 !it could be shortened somewhat by identifying common subexpressions. 98 ! 99 if(rank==1)then 100 cm(1,1,1)=dgm(1,1) 101 cm(1,2,1)=dgm(1,2) 102 cm(1,3,1)=dgm(1,3) 103 cm(1,1,5)=dgm(1,3) 104 cm(1,2,5)=dgm(2,3) 105 cm(1,3,5)=dgm(3,3) 106 cm(1,1,6)=dgm(1,2) 107 cm(1,2,6)=dgm(2,2) 108 cm(1,3,6)=dgm(2,3) 109 cm(2,1,2)=dgm(1,2) 110 cm(2,2,2)=dgm(2,2) 111 cm(2,3,2)=dgm(2,3) 112 cm(2,1,4)=dgm(1,3) 113 cm(2,2,4)=dgm(2,3) 114 cm(2,3,4)=dgm(3,3) 115 cm(2,1,6)=dgm(1,1) 116 cm(2,2,6)=dgm(1,2) 117 cm(2,3,6)=dgm(1,3) 118 cm(3,1,3)=dgm(1,3) 119 cm(3,2,3)=dgm(2,3) 120 cm(3,3,3)=dgm(3,3) 121 cm(3,1,4)=dgm(1,2) 122 cm(3,2,4)=dgm(2,2) 123 cm(3,3,4)=dgm(2,3) 124 cm(3,1,5)=dgm(1,1) 125 cm(3,2,5)=dgm(1,2) 126 cm(3,3,5)=dgm(1,3) 127 elseif(rank==2)then 128 cm(1,1,1)=2*gm(1,1)*dgm(1,1) 129 cm(1,2,1)=-0.5d0*gm(2,2)*dgm(1,1)+3*gm(1,2)*dgm(1,2)-0.5d0*gm(1,1)& 130 & *dgm(2,2) 131 cm(1,3,1)=-0.5d0*gm(3,3)*dgm(1,1)+3*gm(1,3)*dgm(1,3)-0.5d0*gm(1,1)& 132 & *dgm(3,3) 133 cm(1,4,1)=-gm(2,3)*dgm(1,1)+3*gm(1,3)*dgm(1,2)+3*gm(1,2)*dgm(1,3)& 134 & -gm(1,1)*dgm(2,3) 135 cm(1,5,1)=2*(gm(1,3)*dgm(1,1)+gm(1,1)*dgm(1,3)) 136 cm(1,6,1)=2*(gm(1,2)*dgm(1,1)+gm(1,1)*dgm(1,2)) 137 cm(1,1,2)=-0.5d0*gm(2,2)*dgm(1,1)+3*gm(1,2)*dgm(1,2)-0.5d0*gm(1,1)& 138 & *dgm(2,2) 139 cm(1,2,2)=2*gm(2,2)*dgm(2,2) 140 cm(1,3,2)=-0.5d0*gm(3,3)*dgm(2,2)+3*gm(2,3)*dgm(2,3)-0.5d0*gm(2,2)& 141 & *dgm(3,3) 142 cm(1,4,2)=2*(gm(2,3)*dgm(2,2)+gm(2,2)*dgm(2,3)) 143 cm(1,5,2)=3*gm(2,3)*dgm(1,2)-gm(2,2)*dgm(1,3)-gm(1,3)*dgm(2,2)& 144 & +3*gm(1,2)*dgm(2,3) 145 cm(1,6,2)=2*(gm(2,2)*dgm(1,2)+gm(1,2)*dgm(2,2)) 146 cm(1,1,3)=-0.5d0*gm(3,3)*dgm(1,1)+3*gm(1,3)*dgm(1,3)-0.5d0*gm(1,1)& 147 & *dgm(3,3) 148 cm(1,2,3)=-0.5d0*gm(3,3)*dgm(2,2)+3*gm(2,3)*dgm(2,3)-0.5d0*gm(2,2)& 149 & *dgm(3,3) 150 cm(1,3,3)=2*gm(3,3)*dgm(3,3) 151 cm(1,4,3)=2*(gm(3,3)*dgm(2,3)+gm(2,3)*dgm(3,3)) 152 cm(1,5,3)=2*(gm(3,3)*dgm(1,3)+gm(1,3)*dgm(3,3)) 153 cm(1,6,3)=-gm(3,3)*dgm(1,2)+3*gm(2,3)*dgm(1,3)+3*gm(1,3)*dgm(2,3)& 154 & -gm(1,2)*dgm(3,3) 155 cm(1,1,4)=-gm(2,3)*dgm(1,1)+3*gm(1,3)*dgm(1,2)+3*gm(1,2)*dgm(1,3)& 156 & -gm(1,1)*dgm(2,3) 157 cm(1,2,4)=2*(gm(2,3)*dgm(2,2)+gm(2,2)*dgm(2,3)) 158 cm(1,3,4)=2*(gm(3,3)*dgm(2,3)+gm(2,3)*dgm(3,3)) 159 cm(1,4,4)=3*gm(3,3)*dgm(2,2)+2*gm(2,3)*dgm(2,3)+3*gm(2,2)*dgm(3,3) 160 cm(1,5,4)=3*gm(3,3)*dgm(1,2)+gm(2,3)*dgm(1,3)+gm(1,3)*dgm(2,3)& 161 & +3*gm(1,2)*dgm(3,3) 162 cm(1,6,4)=gm(2,3)*dgm(1,2)+3*gm(2,2)*dgm(1,3)+3*gm(1,3)*dgm(2,2)& 163 & +gm(1,2)*dgm(2,3) 164 cm(1,1,5)=2*(gm(1,3)*dgm(1,1)+gm(1,1)*dgm(1,3)) 165 cm(1,2,5)=3*gm(2,3)*dgm(1,2)-gm(2,2)*dgm(1,3)-gm(1,3)*dgm(2,2)& 166 & +3*gm(1,2)*dgm(2,3) 167 cm(1,3,5)=2*(gm(3,3)*dgm(1,3)+gm(1,3)*dgm(3,3)) 168 cm(1,4,5)=3*gm(3,3)*dgm(1,2)+gm(2,3)*dgm(1,3)+gm(1,3)*dgm(2,3)& 169 & +3*gm(1,2)*dgm(3,3) 170 cm(1,5,5)=3*gm(3,3)*dgm(1,1)+2*gm(1,3)*dgm(1,3)+3*gm(1,1)*dgm(3,3) 171 cm(1,6,5)=3*gm(2,3)*dgm(1,1)+gm(1,3)*dgm(1,2)+gm(1,2)*dgm(1,3)& 172 & +3*gm(1,1)*dgm(2,3) 173 cm(1,1,6)=2*(gm(1,2)*dgm(1,1)+gm(1,1)*dgm(1,2)) 174 cm(1,2,6)=2*(gm(2,2)*dgm(1,2)+gm(1,2)*dgm(2,2)) 175 cm(1,3,6)=-gm(3,3)*dgm(1,2)+3*gm(2,3)*dgm(1,3)+3*gm(1,3)*dgm(2,3)& 176 & -gm(1,2)*dgm(3,3) 177 cm(1,4,6)=gm(2,3)*dgm(1,2)+3*gm(2,2)*dgm(1,3)+3*gm(1,3)*dgm(2,2)& 178 & +gm(1,2)*dgm(2,3) 179 cm(1,5,6)=3*gm(2,3)*dgm(1,1)+gm(1,3)*dgm(1,2)+gm(1,2)*dgm(1,3)& 180 & +3*gm(1,1)*dgm(2,3) 181 cm(1,6,6)=3*gm(2,2)*dgm(1,1)+2*gm(1,2)*dgm(1,2)+3*gm(1,1)*dgm(2,2) 182 cm(2,1,2)=2*(gm(1,2)*dgm(1,1)+gm(1,1)*dgm(1,2)) 183 cm(2,2,2)=2*(gm(2,2)*dgm(1,2)+gm(1,2)*dgm(2,2)) 184 cm(2,3,2)=-gm(3,3)*dgm(1,2)+3*gm(2,3)*dgm(1,3)+3*gm(1,3)*dgm(2,3)& 185 & -gm(1,2)*dgm(3,3) 186 cm(2,4,2)=gm(2,3)*dgm(1,2)+3*gm(2,2)*dgm(1,3)+3*gm(1,3)*dgm(2,2)& 187 & +gm(1,2)*dgm(2,3) 188 cm(2,5,2)=3*gm(2,3)*dgm(1,1)+gm(1,3)*dgm(1,2)+gm(1,2)*dgm(1,3)& 189 & +3*gm(1,1)*dgm(2,3) 190 cm(2,6,2)=3*gm(2,2)*dgm(1,1)+2*gm(1,2)*dgm(1,2)+3*gm(1,1)*dgm(2,2) 191 cm(2,1,4)=2*(gm(1,3)*dgm(1,1)+gm(1,1)*dgm(1,3)) 192 cm(2,2,4)=3*gm(2,3)*dgm(1,2)-gm(2,2)*dgm(1,3)-gm(1,3)*dgm(2,2)& 193 & +3*gm(1,2)*dgm(2,3) 194 cm(2,3,4)=2*(gm(3,3)*dgm(1,3)+gm(1,3)*dgm(3,3)) 195 cm(2,4,4)=3*gm(3,3)*dgm(1,2)+gm(2,3)*dgm(1,3)+gm(1,3)*dgm(2,3)& 196 & +3*gm(1,2)*dgm(3,3) 197 cm(2,5,4)=3*gm(3,3)*dgm(1,1)+2*gm(1,3)*dgm(1,3)+3*gm(1,1)*dgm(3,3) 198 cm(2,6,4)=3*gm(2,3)*dgm(1,1)+gm(1,3)*dgm(1,2)+gm(1,2)*dgm(1,3)& 199 & +3*gm(1,1)*dgm(2,3) 200 cm(2,1,6)=2*gm(1,1)*dgm(1,1) 201 cm(2,2,6)=-0.5d0*gm(2,2)*dgm(1,1)+3*gm(1,2)*dgm(1,2)-0.5d0*gm(1,1)& 202 & *dgm(2,2) 203 cm(2,3,6)=-0.5d0*gm(3,3)*dgm(1,1)+3*gm(1,3)*dgm(1,3)-0.5d0*gm(1,1)& 204 & *dgm(3,3) 205 cm(2,4,6)=-gm(2,3)*dgm(1,1)+3*gm(1,3)*dgm(1,2)+3*gm(1,2)*dgm(1,3)& 206 & -gm(1,1)*dgm(2,3) 207 cm(2,5,6)=2*(gm(1,3)*dgm(1,1)+gm(1,1)*dgm(1,3)) 208 cm(2,6,6)=2*(gm(1,2)*dgm(1,1)+gm(1,1)*dgm(1,2)) 209 cm(2,1,7)=-0.5d0*gm(2,2)*dgm(1,1)+3*gm(1,2)*dgm(1,2)-0.5d0*gm(1,1)& 210 & *dgm(2,2) 211 cm(2,2,7)=2*gm(2,2)*dgm(2,2) 212 cm(2,3,7)=-0.5d0*gm(3,3)*dgm(2,2)+3*gm(2,3)*dgm(2,3)-0.5d0*gm(2,2)& 213 & *dgm(3,3) 214 cm(2,4,7)=2*(gm(2,3)*dgm(2,2)+gm(2,2)*dgm(2,3)) 215 cm(2,5,7)=3*gm(2,3)*dgm(1,2)-gm(2,2)*dgm(1,3)-gm(1,3)*dgm(2,2)& 216 & +3*gm(1,2)*dgm(2,3) 217 cm(2,6,7)=2*(gm(2,2)*dgm(1,2)+gm(1,2)*dgm(2,2)) 218 cm(2,1,8)=-0.5d0*gm(3,3)*dgm(1,1)+3*gm(1,3)*dgm(1,3)-0.5d0*gm(1,1)& 219 & *dgm(3,3) 220 cm(2,2,8)=-0.5d0*gm(3,3)*dgm(2,2)+3*gm(2,3)*dgm(2,3)-0.5d0*gm(2,2)& 221 & *dgm(3,3) 222 cm(2,3,8)=2*gm(3,3)*dgm(3,3) 223 cm(2,4,8)=2*(gm(3,3)*dgm(2,3)+gm(2,3)*dgm(3,3)) 224 cm(2,5,8)=2*(gm(3,3)*dgm(1,3)+gm(1,3)*dgm(3,3)) 225 cm(2,6,8)=-gm(3,3)*dgm(1,2)+3*gm(2,3)*dgm(1,3)+3*gm(1,3)*dgm(2,3)& 226 & -gm(1,2)*dgm(3,3) 227 cm(2,1,9)=-gm(2,3)*dgm(1,1)+3*gm(1,3)*dgm(1,2)+3*gm(1,2)*dgm(1,3)& 228 & -gm(1,1)*dgm(2,3) 229 cm(2,2,9)=2*(gm(2,3)*dgm(2,2)+gm(2,2)*dgm(2,3)) 230 cm(2,3,9)=2*(gm(3,3)*dgm(2,3)+gm(2,3)*dgm(3,3)) 231 cm(2,4,9)=3*gm(3,3)*dgm(2,2)+2*gm(2,3)*dgm(2,3)+3*gm(2,2)*dgm(3,3) 232 cm(2,5,9)=3*gm(3,3)*dgm(1,2)+gm(2,3)*dgm(1,3)+gm(1,3)*dgm(2,3)& 233 & +3*gm(1,2)*dgm(3,3) 234 cm(2,6,9)=gm(2,3)*dgm(1,2)+3*gm(2,2)*dgm(1,3)+3*gm(1,3)*dgm(2,2)& 235 & +gm(1,2)*dgm(2,3) 236 cm(3,1,3)=2*(gm(1,3)*dgm(1,1)+gm(1,1)*dgm(1,3)) 237 cm(3,2,3)=3*gm(2,3)*dgm(1,2)-gm(2,2)*dgm(1,3)-gm(1,3)*dgm(2,2)& 238 & +3*gm(1,2)*dgm(2,3) 239 cm(3,3,3)=2*(gm(3,3)*dgm(1,3)+gm(1,3)*dgm(3,3)) 240 cm(3,4,3)=3*gm(3,3)*dgm(1,2)+gm(2,3)*dgm(1,3)+gm(1,3)*dgm(2,3)& 241 & +3*gm(1,2)*dgm(3,3) 242 cm(3,5,3)=3*gm(3,3)*dgm(1,1)+2*gm(1,3)*dgm(1,3)+3*gm(1,1)*dgm(3,3) 243 cm(3,6,3)=3*gm(2,3)*dgm(1,1)+gm(1,3)*dgm(1,2)+gm(1,2)*dgm(1,3)& 244 & +3*gm(1,1)*dgm(2,3) 245 cm(3,1,4)=2*(gm(1,2)*dgm(1,1)+gm(1,1)*dgm(1,2)) 246 cm(3,2,4)=2*(gm(2,2)*dgm(1,2)+gm(1,2)*dgm(2,2)) 247 cm(3,3,4)=-gm(3,3)*dgm(1,2)+3*gm(2,3)*dgm(1,3)+3*gm(1,3)*dgm(2,3)& 248 & -gm(1,2)*dgm(3,3) 249 cm(3,4,4)=gm(2,3)*dgm(1,2)+3*gm(2,2)*dgm(1,3)+3*gm(1,3)*dgm(2,2)& 250 & +gm(1,2)*dgm(2,3) 251 cm(3,5,4)=3*gm(2,3)*dgm(1,1)+gm(1,3)*dgm(1,2)+gm(1,2)*dgm(1,3)& 252 & +3*gm(1,1)*dgm(2,3) 253 cm(3,6,4)=3*gm(2,2)*dgm(1,1)+2*gm(1,2)*dgm(1,2)+3*gm(1,1)*dgm(2,2) 254 cm(3,1,5)=2*gm(1,1)*dgm(1,1) 255 cm(3,2,5)=-0.5d0*gm(2,2)*dgm(1,1)+3*gm(1,2)*dgm(1,2)-0.5d0*gm(1,1)& 256 & *dgm(2,2) 257 cm(3,3,5)=-0.5d0*gm(3,3)*dgm(1,1)+3*gm(1,3)*dgm(1,3)-0.5d0*gm(1,1)& 258 & *dgm(3,3) 259 cm(3,4,5)=-gm(2,3)*dgm(1,1)+3*gm(1,3)*dgm(1,2)+3*gm(1,2)*dgm(1,3)& 260 & -gm(1,1)*dgm(2,3) 261 cm(3,5,5)=2*(gm(1,3)*dgm(1,1)+gm(1,1)*dgm(1,3)) 262 cm(3,6,5)=2*(gm(1,2)*dgm(1,1)+gm(1,1)*dgm(1,2)) 263 cm(3,1,8)=-gm(2,3)*dgm(1,1)+3*gm(1,3)*dgm(1,2)+3*gm(1,2)*dgm(1,3)& 264 & -gm(1,1)*dgm(2,3) 265 cm(3,2,8)=2*(gm(2,3)*dgm(2,2)+gm(2,2)*dgm(2,3)) 266 cm(3,3,8)=2*(gm(3,3)*dgm(2,3)+gm(2,3)*dgm(3,3)) 267 cm(3,4,8)=3*gm(3,3)*dgm(2,2)+2*gm(2,3)*dgm(2,3)+3*gm(2,2)*dgm(3,3) 268 cm(3,5,8)=3*gm(3,3)*dgm(1,2)+gm(2,3)*dgm(1,3)+gm(1,3)*dgm(2,3)& 269 & +3*gm(1,2)*dgm(3,3) 270 cm(3,6,8)=gm(2,3)*dgm(1,2)+3*gm(2,2)*dgm(1,3)+3*gm(1,3)*dgm(2,2)& 271 & +gm(1,2)*dgm(2,3) 272 cm(3,1,9)=-0.5d0*gm(2,2)*dgm(1,1)+3*gm(1,2)*dgm(1,2)-0.5d0*gm(1,1)& 273 & *dgm(2,2) 274 cm(3,2,9)=2*gm(2,2)*dgm(2,2) 275 cm(3,3,9)=-0.5d0*gm(3,3)*dgm(2,2)+3*gm(2,3)*dgm(2,3)-0.5d0*gm(2,2)& 276 & *dgm(3,3) 277 cm(3,4,9)=2*(gm(2,3)*dgm(2,2)+gm(2,2)*dgm(2,3)) 278 cm(3,5,9)=3*gm(2,3)*dgm(1,2)-gm(2,2)*dgm(1,3)-gm(1,3)*dgm(2,2)& 279 & +3*gm(1,2)*dgm(2,3) 280 cm(3,6,9)=2*(gm(2,2)*dgm(1,2)+gm(1,2)*dgm(2,2)) 281 cm(3,1,10)=-0.5d0*gm(3,3)*dgm(1,1)+3*gm(1,3)*dgm(1,3)-0.5d0*gm(1,1)& 282 & *dgm(3,3) 283 cm(3,2,10)=-0.5d0*gm(3,3)*dgm(2,2)+3*gm(2,3)*dgm(2,3)-0.5d0*gm(2,2)& 284 & *dgm(3,3) 285 cm(3,3,10)=2*gm(3,3)*dgm(3,3) 286 cm(3,4,10)=2*(gm(3,3)*dgm(2,3)+gm(2,3)*dgm(3,3)) 287 cm(3,5,10)=2*(gm(3,3)*dgm(1,3)+gm(1,3)*dgm(3,3)) 288 cm(3,6,10)=-gm(3,3)*dgm(1,2)+3*gm(2,3)*dgm(1,3)+3*gm(1,3)*dgm(2,3)& 289 & -gm(1,2)*dgm(3,3) 290 elseif(rank==3)then 291 cm(1,1,1)=3*gm(1,1)**2*dgm(1,1) 292 cm(1,2,1)=4.5d0*gm(1,2)**2*dgm(1,1)+9*gm(1,1)*gm(1,2)*dgm(1,2)& 293 & +gm(1,1)*(-3*gm(2,2)*dgm(1,1)-1.5d0*gm(1,1)*dgm(2,2)) 294 cm(1,3,1)=4.5d0*gm(1,3)**2*dgm(1,1)+9*gm(1,1)*gm(1,3)*dgm(1,3)& 295 & +gm(1,1)*(-3*gm(3,3)*dgm(1,1)-1.5d0*gm(1,1)*dgm(3,3)) 296 cm(1,4,1)=9*gm(1,2)*(gm(1,3)*dgm(1,1)+gm(1,1)*dgm(1,3))+gm(1,1)& 297 & *(-6*gm(2,3)*dgm(1,1)+9*gm(1,3)*dgm(1,2)-3*gm(1,1)*dgm(2,3)) 298 cm(1,5,1)=gm(1,1)*(6*gm(1,3)*dgm(1,1)+3*gm(1,1)*dgm(1,3)) 299 cm(1,6,1)=gm(1,1)*(6*gm(1,2)*dgm(1,1)+3*gm(1,1)*dgm(1,2)) 300 cm(1,7,1)=7.5d0*gm(1,2)**2*dgm(1,2)-1.5d0*gm(1,1)*gm(2,2)*dgm(1,2)& 301 & -1.5d0*gm(1,2)*(gm(2,2)*dgm(1,1)+gm(1,1)*dgm(2,2)) 302 cm(1,8,1)=7.5d0*gm(1,3)**2*dgm(1,2)+gm(1,1)*(-1.5d0*gm(3,3)*dgm(1,2)& 303 & -3*gm(2,3)*dgm(1,3))+gm(1,3)*(-3*gm(2,3)*dgm(1,1)+15*gm(1,2)& 304 & *dgm(1,3)-3*gm(1,1)*dgm(2,3))-1.5d0*gm(1,2)*(gm(3,3)*dgm(1,1)& 305 & +gm(1,1)*dgm(3,3)) 306 cm(1,9,1)=7.5d0*gm(1,2)**2*dgm(1,3)+gm(1,1)*(-3*gm(2,3)*dgm(1,2)& 307 & -1.5d0*gm(2,2)*dgm(1,3))+gm(1,3)*(-1.5d0*gm(2,2)*dgm(1,1)+15*gm(1,2)& 308 & *dgm(1,2)-1.5d0*gm(1,1)*dgm(2,2))-3*gm(1,2)*(gm(2,3)*dgm(1,1)& 309 & +gm(1,1)*dgm(2,3)) 310 cm(1,10,1)=7.5d0*gm(1,3)**2*dgm(1,3)-1.5d0*gm(1,1)*gm(3,3)*dgm(1,3)& 311 & -1.5d0*gm(1,3)*(gm(3,3)*dgm(1,1)+gm(1,1)*dgm(3,3)) 312 cm(1,1,2)=4.5d0*gm(1,2)**2*dgm(1,1)+9*gm(1,1)*gm(1,2)*dgm(1,2)& 313 & +gm(1,1)*(-3*gm(2,2)*dgm(1,1)-1.5d0*gm(1,1)*dgm(2,2)) 314 cm(1,2,2)=6*gm(2,2)**2*dgm(1,1)+3*gm(1,2)**2*dgm(2,2)+gm(2,2)& 315 & *(6*gm(1,2)*dgm(1,2)+12*gm(1,1)*dgm(2,2)) 316 cm(1,3,2)=7.5d0*gm(2,3)**2*dgm(1,1)-6*gm(1,2)*gm(3,3)*dgm(1,2)& 317 & -3*gm(1,3)**2*dgm(2,2)-1.5d0*gm(1,1)*gm(3,3)*dgm(2,2)+9*gm(1,2)& 318 & *gm(1,3)*dgm(2,3)+gm(2,3)*(9*gm(1,3)*dgm(1,2)+9*gm(1,2)*dgm(1,3)& 319 & +15*gm(1,1)*dgm(2,3))-3*gm(1,2)**2*dgm(3,3)+gm(2,2)*(-1.5d0*gm(3,3)& 320 & *dgm(1,1)-6*gm(1,3)*dgm(1,3)-1.5d0*gm(1,1)*dgm(3,3)) 321 cm(1,4,2)=12*gm(1,1)*gm(2,3)*dgm(2,2)+gm(1,2)*(6*gm(2,3)*dgm(1,2)& 322 & +3*gm(1,3)*dgm(2,2))+3*gm(1,2)**2*dgm(2,3)+gm(2,2)*(12*gm(2,3)& 323 & *dgm(1,1)+3*gm(1,3)*dgm(1,2)+3*gm(1,2)*dgm(1,3)+12*gm(1,1)*dgm(2,3)) 324 cm(1,5,2)=1.5d0*gm(1,2)**2*dgm(1,3)+gm(1,1)*(12*gm(2,3)*dgm(1,2)& 325 & -4.5d0*gm(2,2)*dgm(1,3))+gm(1,3)*(-4.5d0*gm(2,2)*dgm(1,1)+3*gm(1,2)& 326 & *dgm(1,2)-4.5d0*gm(1,1)*dgm(2,2))+12*gm(1,2)*(gm(2,3)*dgm(1,1)& 327 & +gm(1,1)*dgm(2,3)) 328 cm(1,6,2)=4.5d0*gm(1,2)**2*dgm(1,2)+7.5d0*gm(1,1)*gm(2,2)*dgm(1,2)& 329 & +7.5d0*gm(1,2)*(gm(2,2)*dgm(1,1)+gm(1,1)*dgm(2,2)) 330 cm(1,7,2)=gm(2,2)*(3*gm(2,2)*dgm(1,2)+6*gm(1,2)*dgm(2,2)) 331 cm(1,8,2)=1.5d0*gm(2,3)**2*dgm(1,2)-4.5d0*gm(1,2)*gm(3,3)*dgm(2,2)& 332 & +gm(2,3)*(12*gm(2,2)*dgm(1,3)+12*gm(1,3)*dgm(2,2)+3*gm(1,2)*dgm(2,3))& 333 & +gm(2,2)*(-4.5d0*gm(3,3)*dgm(1,2)+12*gm(1,3)*dgm(2,3)-4.5d0*gm(1,2)& 334 & *dgm(3,3)) 335 cm(1,9,2)=6*gm(2,2)**2*dgm(1,3)+3*gm(1,2)*gm(2,3)*dgm(2,2)+gm(2,2)& 336 & *(3*gm(2,3)*dgm(1,2)+12*gm(1,3)*dgm(2,2)+3*gm(1,2)*dgm(2,3)) 337 cm(1,10,2)=7.5d0*gm(2,3)**2*dgm(1,3)+gm(3,3)*(-1.5d0*gm(1,3)*dgm(2,2)& 338 & -3*gm(1,2)*dgm(2,3))+gm(2,3)*(-3*gm(3,3)*dgm(1,2)+15*gm(1,3)& 339 & *dgm(2,3)-3*gm(1,2)*dgm(3,3))-1.5d0*gm(2,2)*(gm(3,3)*dgm(1,3)& 340 & +gm(1,3)*dgm(3,3)) 341 cm(1,1,3)=4.5d0*gm(1,3)**2*dgm(1,1)+9*gm(1,1)*gm(1,3)*dgm(1,3)& 342 & +gm(1,1)*(-3*gm(3,3)*dgm(1,1)-1.5d0*gm(1,1)*dgm(3,3)) 343 cm(1,2,3)=7.5d0*gm(2,3)**2*dgm(1,1)-6*gm(1,2)*gm(3,3)*dgm(1,2)& 344 & -3*gm(1,3)**2*dgm(2,2)-1.5d0*gm(1,1)*gm(3,3)*dgm(2,2)+9*gm(1,2)& 345 & *gm(1,3)*dgm(2,3)+gm(2,3)*(9*gm(1,3)*dgm(1,2)+9*gm(1,2)*dgm(1,3)& 346 & +15*gm(1,1)*dgm(2,3))-3*gm(1,2)**2*dgm(3,3)+gm(2,2)*(-1.5d0*gm(3,3)& 347 & *dgm(1,1)-6*gm(1,3)*dgm(1,3)-1.5d0*gm(1,1)*dgm(3,3)) 348 cm(1,3,3)=6*gm(3,3)**2*dgm(1,1)+3*gm(1,3)**2*dgm(3,3)+gm(3,3)& 349 & *(6*gm(1,3)*dgm(1,3)+12*gm(1,1)*dgm(3,3)) 350 cm(1,4,3)=3*gm(1,3)**2*dgm(2,3)+gm(3,3)*(3*gm(1,2)*dgm(1,3)+12*gm(1,1)& 351 & *dgm(2,3))+gm(2,3)*(12*gm(3,3)*dgm(1,1)+6*gm(1,3)*dgm(1,3)+12*gm(1,1)& 352 & *dgm(3,3))+3*gm(1,3)*(gm(3,3)*dgm(1,2)+gm(1,2)*dgm(3,3)) 353 cm(1,5,3)=4.5d0*gm(1,3)**2*dgm(1,3)+7.5d0*gm(1,1)*gm(3,3)*dgm(1,3)& 354 & +7.5d0*gm(1,3)*(gm(3,3)*dgm(1,1)+gm(1,1)*dgm(3,3)) 355 cm(1,6,3)=1.5d0*gm(1,3)**2*dgm(1,2)+gm(1,1)*(-4.5d0*gm(3,3)*dgm(1,2)& 356 & +12*gm(2,3)*dgm(1,3))+gm(1,3)*(12*gm(2,3)*dgm(1,1)+3*gm(1,2)& 357 & *dgm(1,3)+12*gm(1,1)*dgm(2,3))-4.5d0*gm(1,2)*(gm(3,3)*dgm(1,1)& 358 & +gm(1,1)*dgm(3,3)) 359 cm(1,7,3)=7.5d0*gm(2,3)**2*dgm(1,2)-1.5d0*gm(1,2)*gm(3,3)*dgm(2,2)& 360 & +gm(2,3)*(-3*gm(2,2)*dgm(1,3)-3*gm(1,3)*dgm(2,2)+15*gm(1,2)*dgm(2,3))& 361 & +gm(2,2)*(-1.5d0*gm(3,3)*dgm(1,2)-3*gm(1,3)*dgm(2,3)-1.5d0*gm(1,2)& 362 & *dgm(3,3)) 363 cm(1,8,3)=6*gm(3,3)**2*dgm(1,2)+3*gm(1,3)*gm(2,3)*dgm(3,3)+gm(3,3)& 364 & *(3*gm(2,3)*dgm(1,3)+3*gm(1,3)*dgm(2,3)+12*gm(1,2)*dgm(3,3)) 365 cm(1,9,3)=1.5d0*gm(2,3)**2*dgm(1,3)+gm(3,3)*(-4.5d0*gm(1,3)*dgm(2,2)& 366 & +12*gm(1,2)*dgm(2,3))+gm(2,3)*(12*gm(3,3)*dgm(1,2)+3*gm(1,3)& 367 & *dgm(2,3)+12*gm(1,2)*dgm(3,3))-4.5d0*gm(2,2)*(gm(3,3)*dgm(1,3)& 368 & +gm(1,3)*dgm(3,3)) 369 cm(1,10,3)=gm(3,3)*(3*gm(3,3)*dgm(1,3)+6*gm(1,3)*dgm(3,3)) 370 cm(1,1,4)=9*gm(1,2)*(gm(1,3)*dgm(1,1)+gm(1,1)*dgm(1,3))+gm(1,1)& 371 & *(-6*gm(2,3)*dgm(1,1)+9*gm(1,3)*dgm(1,2)-3*gm(1,1)*dgm(2,3)) 372 cm(1,2,4)=12*gm(1,1)*gm(2,3)*dgm(2,2)+gm(1,2)*(6*gm(2,3)*dgm(1,2)& 373 & +3*gm(1,3)*dgm(2,2))+3*gm(1,2)**2*dgm(2,3)+gm(2,2)*(12*gm(2,3)& 374 & *dgm(1,1)+3*gm(1,3)*dgm(1,2)+3*gm(1,2)*dgm(1,3)+12*gm(1,1)*dgm(2,3)) 375 cm(1,3,4)=3*gm(1,3)**2*dgm(2,3)+gm(3,3)*(3*gm(1,2)*dgm(1,3)+12*gm(1,1)& 376 & *dgm(2,3))+gm(2,3)*(12*gm(3,3)*dgm(1,1)+6*gm(1,3)*dgm(1,3)+12*gm(1,1)& 377 & *dgm(3,3))+3*gm(1,3)*(gm(3,3)*dgm(1,2)+gm(1,2)*dgm(3,3)) 378 cm(1,4,4)=9*gm(2,3)**2*dgm(1,1)+18*gm(1,2)*gm(3,3)*dgm(1,2)+9*gm(1,3)& 379 & **2*dgm(2,2)+15*gm(1,1)*gm(3,3)*dgm(2,2)-6*gm(1,2)*gm(1,3)*dgm(2,3)& 380 & +gm(2,3)*(-6*gm(1,3)*dgm(1,2)-6*gm(1,2)*dgm(1,3)+18*gm(1,1)*dgm(2,3))& 381 & +9*gm(1,2)**2*dgm(3,3)+gm(2,2)*(15*gm(3,3)*dgm(1,1)+18*gm(1,3)& 382 & *dgm(1,3)+15*gm(1,1)*dgm(3,3)) 383 cm(1,5,4)=3*gm(1,3)**2*dgm(1,2)+gm(1,1)*(12*gm(3,3)*dgm(1,2)+3*gm(2,3)& 384 & *dgm(1,3))+gm(1,3)*(3*gm(2,3)*dgm(1,1)+6*gm(1,2)*dgm(1,3)+3*gm(1,1)& 385 & *dgm(2,3))+12*gm(1,2)*(gm(3,3)*dgm(1,1)+gm(1,1)*dgm(3,3)) 386 cm(1,6,4)=3*gm(1,2)**2*dgm(1,3)+gm(1,1)*(3*gm(2,3)*dgm(1,2)+12*gm(2,2)& 387 & *dgm(1,3))+gm(1,3)*(12*gm(2,2)*dgm(1,1)+6*gm(1,2)*dgm(1,2)+12*gm(1,1)& 388 & *dgm(2,2))+3*gm(1,2)*(gm(2,3)*dgm(1,1)+gm(1,1)*dgm(2,3)) 389 cm(1,7,4)=-3*gm(2,2)**2*dgm(1,3)+9*gm(1,2)*gm(2,3)*dgm(2,2)+gm(2,2)& 390 & *(9*gm(2,3)*dgm(1,2)-6*gm(1,3)*dgm(2,2)+9*gm(1,2)*dgm(2,3)) 391 cm(1,8,4)=3*gm(2,3)**2*dgm(1,3)+gm(3,3)*(12*gm(1,3)*dgm(2,2)+3*gm(1,2)& 392 & *dgm(2,3))+gm(2,3)*(3*gm(3,3)*dgm(1,2)+6*gm(1,3)*dgm(2,3)+3*gm(1,2)& 393 & *dgm(3,3))+12*gm(2,2)*(gm(3,3)*dgm(1,3)+gm(1,3)*dgm(3,3)) 394 cm(1,9,4)=3*gm(2,3)**2*dgm(1,2)+12*gm(1,2)*gm(3,3)*dgm(2,2)+gm(2,3)& 395 & *(3*gm(2,2)*dgm(1,3)+3*gm(1,3)*dgm(2,2)+6*gm(1,2)*dgm(2,3))+gm(2,2)& 396 & *(12*gm(3,3)*dgm(1,2)+3*gm(1,3)*dgm(2,3)+12*gm(1,2)*dgm(3,3)) 397 cm(1,10,4)=-3*gm(3,3)**2*dgm(1,2)+9*gm(1,3)*gm(2,3)*dgm(3,3)+gm(3,3)& 398 & *(9*gm(2,3)*dgm(1,3)+9*gm(1,3)*dgm(2,3)-6*gm(1,2)*dgm(3,3)) 399 cm(1,1,5)=gm(1,1)*(6*gm(1,3)*dgm(1,1)+3*gm(1,1)*dgm(1,3)) 400 cm(1,2,5)=1.5d0*gm(1,2)**2*dgm(1,3)+gm(1,1)*(12*gm(2,3)*dgm(1,2)& 401 & -4.5d0*gm(2,2)*dgm(1,3))+gm(1,3)*(-4.5d0*gm(2,2)*dgm(1,1)+3*gm(1,2)& 402 & *dgm(1,2)-4.5d0*gm(1,1)*dgm(2,2))+12*gm(1,2)*(gm(2,3)*dgm(1,1)& 403 & +gm(1,1)*dgm(2,3)) 404 cm(1,3,5)=4.5d0*gm(1,3)**2*dgm(1,3)+7.5d0*gm(1,1)*gm(3,3)*dgm(1,3)& 405 & +7.5d0*gm(1,3)*(gm(3,3)*dgm(1,1)+gm(1,1)*dgm(3,3)) 406 cm(1,4,5)=3*gm(1,3)**2*dgm(1,2)+gm(1,1)*(12*gm(3,3)*dgm(1,2)+3*gm(2,3)& 407 & *dgm(1,3))+gm(1,3)*(3*gm(2,3)*dgm(1,1)+6*gm(1,2)*dgm(1,3)+3*gm(1,1)& 408 & *dgm(2,3))+12*gm(1,2)*(gm(3,3)*dgm(1,1)+gm(1,1)*dgm(3,3)) 409 cm(1,5,5)=3*gm(1,3)**2*dgm(1,1)+6*gm(1,1)*gm(1,3)*dgm(1,3)+gm(1,1)& 410 & *(12*gm(3,3)*dgm(1,1)+6*gm(1,1)*dgm(3,3)) 411 cm(1,6,5)=3*gm(1,2)*(gm(1,3)*dgm(1,1)+gm(1,1)*dgm(1,3))+gm(1,1)& 412 & *(12*gm(2,3)*dgm(1,1)+3*gm(1,3)*dgm(1,2)+6*gm(1,1)*dgm(2,3)) 413 cm(1,7,5)=-1.5d0*gm(1,1)*gm(2,3)*dgm(2,2)+gm(1,2)*(15*gm(2,3)& 414 & *dgm(1,2)-3*gm(1,3)*dgm(2,2))+7.5d0*gm(1,2)**2*dgm(2,3)+gm(2,2)& 415 & *(-1.5d0*gm(2,3)*dgm(1,1)-3*gm(1,3)*dgm(1,2)-3*gm(1,2)*dgm(1,3)& 416 & -1.5d0*gm(1,1)*dgm(2,3)) 417 cm(1,8,5)=1.5d0*gm(1,3)**2*dgm(2,3)+gm(3,3)*(12*gm(1,2)*dgm(1,3)& 418 & -4.5d0*gm(1,1)*dgm(2,3))+gm(2,3)*(-4.5d0*gm(3,3)*dgm(1,1)+3*gm(1,3)& 419 & *dgm(1,3)-4.5d0*gm(1,1)*dgm(3,3))+12*gm(1,3)*(gm(3,3)*dgm(1,2)& 420 & +gm(1,2)*dgm(3,3)) 421 cm(1,9,5)=-3*gm(2,3)**2*dgm(1,1)+15*gm(1,2)*gm(3,3)*dgm(1,2)-3*gm(1,3)& 422 & **2*dgm(2,2)-1.5d0*gm(1,1)*gm(3,3)*dgm(2,2)+9*gm(1,2)*gm(1,3)& 423 & *dgm(2,3)+gm(2,3)*(9*gm(1,3)*dgm(1,2)+9*gm(1,2)*dgm(1,3)-6*gm(1,1)& 424 & *dgm(2,3))+7.5d0*gm(1,2)**2*dgm(3,3)+gm(2,2)*(-1.5d0*gm(3,3)& 425 & *dgm(1,1)-6*gm(1,3)*dgm(1,3)-1.5d0*gm(1,1)*dgm(3,3)) 426 cm(1,10,5)=-1.5d0*gm(3,3)**2*dgm(1,1)+4.5d0*gm(1,3)**2*dgm(3,3)& 427 & +gm(3,3)*(9*gm(1,3)*dgm(1,3)-3*gm(1,1)*dgm(3,3)) 428 cm(1,1,6)=gm(1,1)*(6*gm(1,2)*dgm(1,1)+3*gm(1,1)*dgm(1,2)) 429 cm(1,2,6)=4.5d0*gm(1,2)**2*dgm(1,2)+7.5d0*gm(1,1)*gm(2,2)*dgm(1,2)& 430 & +7.5d0*gm(1,2)*(gm(2,2)*dgm(1,1)+gm(1,1)*dgm(2,2)) 431 cm(1,3,6)=1.5d0*gm(1,3)**2*dgm(1,2)+gm(1,1)*(-4.5d0*gm(3,3)*dgm(1,2)& 432 & +12*gm(2,3)*dgm(1,3))+gm(1,3)*(12*gm(2,3)*dgm(1,1)+3*gm(1,2)& 433 & *dgm(1,3)+12*gm(1,1)*dgm(2,3))-4.5d0*gm(1,2)*(gm(3,3)*dgm(1,1)& 434 & +gm(1,1)*dgm(3,3)) 435 cm(1,4,6)=3*gm(1,2)**2*dgm(1,3)+gm(1,1)*(3*gm(2,3)*dgm(1,2)+12*gm(2,2)& 436 & *dgm(1,3))+gm(1,3)*(12*gm(2,2)*dgm(1,1)+6*gm(1,2)*dgm(1,2)+12*gm(1,1)& 437 & *dgm(2,2))+3*gm(1,2)*(gm(2,3)*dgm(1,1)+gm(1,1)*dgm(2,3)) 438 cm(1,5,6)=3*gm(1,2)*(gm(1,3)*dgm(1,1)+gm(1,1)*dgm(1,3))+gm(1,1)& 439 & *(12*gm(2,3)*dgm(1,1)+3*gm(1,3)*dgm(1,2)+6*gm(1,1)*dgm(2,3)) 440 cm(1,6,6)=3*gm(1,2)**2*dgm(1,1)+6*gm(1,1)*gm(1,2)*dgm(1,2)+gm(1,1)& 441 & *(12*gm(2,2)*dgm(1,1)+6*gm(1,1)*dgm(2,2)) 442 cm(1,7,6)=-1.5d0*gm(2,2)**2*dgm(1,1)+4.5d0*gm(1,2)**2*dgm(2,2)& 443 & +gm(2,2)*(9*gm(1,2)*dgm(1,2)-3*gm(1,1)*dgm(2,2)) 444 cm(1,8,6)=-3*gm(2,3)**2*dgm(1,1)-6*gm(1,2)*gm(3,3)*dgm(1,2)+7.5d0*gm(1,3)& 445 & **2*dgm(2,2)-1.5d0*gm(1,1)*gm(3,3)*dgm(2,2)+9*gm(1,2)*gm(1,3)& 446 & *dgm(2,3)+gm(2,3)*(9*gm(1,3)*dgm(1,2)+9*gm(1,2)*dgm(1,3)-6*gm(1,1)& 447 & *dgm(2,3))-3*gm(1,2)**2*dgm(3,3)+gm(2,2)*(-1.5d0*gm(3,3)*dgm(1,1)& 448 & +15*gm(1,3)*dgm(1,3)-1.5d0*gm(1,1)*dgm(3,3)) 449 cm(1,9,6)=-4.5d0*gm(1,1)*gm(2,3)*dgm(2,2)+gm(1,2)*(3*gm(2,3)*dgm(1,2)& 450 & +12*gm(1,3)*dgm(2,2))+1.5d0*gm(1,2)**2*dgm(2,3)+gm(2,2)*(-4.5d0*gm(2,3)& 451 & *dgm(1,1)+12*gm(1,3)*dgm(1,2)+12*gm(1,2)*dgm(1,3)-4.5d0*gm(1,1)& 452 & *dgm(2,3)) 453 cm(1,10,6)=7.5d0*gm(1,3)**2*dgm(2,3)+gm(3,3)*(-3*gm(1,2)*dgm(1,3)& 454 & -1.5d0*gm(1,1)*dgm(2,3))+gm(2,3)*(-1.5d0*gm(3,3)*dgm(1,1)+15*gm(1,3)& 455 & *dgm(1,3)-1.5d0*gm(1,1)*dgm(3,3))-3*gm(1,3)*(gm(3,3)*dgm(1,2)& 456 & +gm(1,2)*dgm(3,3)) 457 cm(1,1,7)=7.5d0*gm(1,2)**2*dgm(1,2)-1.5d0*gm(1,1)*gm(2,2)*dgm(1,2)& 458 & -1.5d0*gm(1,2)*(gm(2,2)*dgm(1,1)+gm(1,1)*dgm(2,2)) 459 cm(1,2,7)=gm(2,2)*(3*gm(2,2)*dgm(1,2)+6*gm(1,2)*dgm(2,2)) 460 cm(1,3,7)=7.5d0*gm(2,3)**2*dgm(1,2)-1.5d0*gm(1,2)*gm(3,3)*dgm(2,2)& 461 & +gm(2,3)*(-3*gm(2,2)*dgm(1,3)-3*gm(1,3)*dgm(2,2)+15*gm(1,2)*dgm(2,3))& 462 & +gm(2,2)*(-1.5d0*gm(3,3)*dgm(1,2)-3*gm(1,3)*dgm(2,3)-1.5d0*gm(1,2)& 463 & *dgm(3,3)) 464 cm(1,4,7)=-3*gm(2,2)**2*dgm(1,3)+9*gm(1,2)*gm(2,3)*dgm(2,2)+gm(2,2)& 465 & *(9*gm(2,3)*dgm(1,2)-6*gm(1,3)*dgm(2,2)+9*gm(1,2)*dgm(2,3)) 466 cm(1,5,7)=-1.5d0*gm(1,1)*gm(2,3)*dgm(2,2)+gm(1,2)*(15*gm(2,3)& 467 & *dgm(1,2)-3*gm(1,3)*dgm(2,2))+7.5d0*gm(1,2)**2*dgm(2,3)+gm(2,2)& 468 & *(-1.5d0*gm(2,3)*dgm(1,1)-3*gm(1,3)*dgm(1,2)-3*gm(1,2)*dgm(1,3)& 469 & -1.5d0*gm(1,1)*dgm(2,3)) 470 cm(1,6,7)=-1.5d0*gm(2,2)**2*dgm(1,1)+4.5d0*gm(1,2)**2*dgm(2,2)& 471 & +gm(2,2)*(9*gm(1,2)*dgm(1,2)-3*gm(1,1)*dgm(2,2)) 472 cm(1,7,7)=3*gm(2,2)**2*dgm(2,2) 473 cm(1,8,7)=4.5d0*gm(2,3)**2*dgm(2,2)+9*gm(2,2)*gm(2,3)*dgm(2,3)& 474 & +gm(2,2)*(-3*gm(3,3)*dgm(2,2)-1.5d0*gm(2,2)*dgm(3,3)) 475 cm(1,9,7)=gm(2,2)*(6*gm(2,3)*dgm(2,2)+3*gm(2,2)*dgm(2,3)) 476 cm(1,10,7)=7.5d0*gm(2,3)**2*dgm(2,3)-1.5d0*gm(2,2)*gm(3,3)*dgm(2,3)& 477 & -1.5d0*gm(2,3)*(gm(3,3)*dgm(2,2)+gm(2,2)*dgm(3,3)) 478 cm(1,1,8)=7.5d0*gm(1,3)**2*dgm(1,2)+gm(1,1)*(-1.5d0*gm(3,3)*dgm(1,2)& 479 & -3*gm(2,3)*dgm(1,3))+gm(1,3)*(-3*gm(2,3)*dgm(1,1)+15*gm(1,2)& 480 & *dgm(1,3)-3*gm(1,1)*dgm(2,3))-1.5d0*gm(1,2)*(gm(3,3)*dgm(1,1)& 481 & +gm(1,1)*dgm(3,3)) 482 cm(1,2,8)=1.5d0*gm(2,3)**2*dgm(1,2)-4.5d0*gm(1,2)*gm(3,3)*dgm(2,2)& 483 & +gm(2,3)*(12*gm(2,2)*dgm(1,3)+12*gm(1,3)*dgm(2,2)+3*gm(1,2)*dgm(2,3))& 484 & +gm(2,2)*(-4.5d0*gm(3,3)*dgm(1,2)+12*gm(1,3)*dgm(2,3)-4.5d0*gm(1,2)& 485 & *dgm(3,3)) 486 cm(1,3,8)=6*gm(3,3)**2*dgm(1,2)+3*gm(1,3)*gm(2,3)*dgm(3,3)+gm(3,3)& 487 & *(3*gm(2,3)*dgm(1,3)+3*gm(1,3)*dgm(2,3)+12*gm(1,2)*dgm(3,3)) 488 cm(1,4,8)=3*gm(2,3)**2*dgm(1,3)+gm(3,3)*(12*gm(1,3)*dgm(2,2)+3*gm(1,2)& 489 & *dgm(2,3))+gm(2,3)*(3*gm(3,3)*dgm(1,2)+6*gm(1,3)*dgm(2,3)+3*gm(1,2)& 490 & *dgm(3,3))+12*gm(2,2)*(gm(3,3)*dgm(1,3)+gm(1,3)*dgm(3,3)) 491 cm(1,5,8)=1.5d0*gm(1,3)**2*dgm(2,3)+gm(3,3)*(12*gm(1,2)*dgm(1,3)& 492 & -4.5d0*gm(1,1)*dgm(2,3))+gm(2,3)*(-4.5d0*gm(3,3)*dgm(1,1)+3*gm(1,3)& 493 & *dgm(1,3)-4.5d0*gm(1,1)*dgm(3,3))+12*gm(1,3)*(gm(3,3)*dgm(1,2)& 494 & +gm(1,2)*dgm(3,3)) 495 cm(1,6,8)=-3*gm(2,3)**2*dgm(1,1)-6*gm(1,2)*gm(3,3)*dgm(1,2)+7.5d0*gm(1,3)& 496 & **2*dgm(2,2)-1.5d0*gm(1,1)*gm(3,3)*dgm(2,2)+9*gm(1,2)*gm(1,3)& 497 & *dgm(2,3)+gm(2,3)*(9*gm(1,3)*dgm(1,2)+9*gm(1,2)*dgm(1,3)-6*gm(1,1)& 498 & *dgm(2,3))-3*gm(1,2)**2*dgm(3,3)+gm(2,2)*(-1.5d0*gm(3,3)*dgm(1,1)& 499 & +15*gm(1,3)*dgm(1,3)-1.5d0*gm(1,1)*dgm(3,3)) 500 cm(1,7,8)=4.5d0*gm(2,3)**2*dgm(2,2)+9*gm(2,2)*gm(2,3)*dgm(2,3)& 501 & +gm(2,2)*(-3*gm(3,3)*dgm(2,2)-1.5d0*gm(2,2)*dgm(3,3)) 502 cm(1,8,8)=6*gm(3,3)**2*dgm(2,2)+3*gm(2,3)**2*dgm(3,3)+gm(3,3)& 503 & *(6*gm(2,3)*dgm(2,3)+12*gm(2,2)*dgm(3,3)) 504 cm(1,9,8)=4.5d0*gm(2,3)**2*dgm(2,3)+7.5d0*gm(2,2)*gm(3,3)*dgm(2,3)& 505 & +7.5d0*gm(2,3)*(gm(3,3)*dgm(2,2)+gm(2,2)*dgm(3,3)) 506 cm(1,10,8)=gm(3,3)*(3*gm(3,3)*dgm(2,3)+6*gm(2,3)*dgm(3,3)) 507 cm(1,1,9)=7.5d0*gm(1,2)**2*dgm(1,3)+gm(1,1)*(-3*gm(2,3)*dgm(1,2)& 508 & -1.5d0*gm(2,2)*dgm(1,3))+gm(1,3)*(-1.5d0*gm(2,2)*dgm(1,1)+15*gm(1,2)& 509 & *dgm(1,2)-1.5d0*gm(1,1)*dgm(2,2))-3*gm(1,2)*(gm(2,3)*dgm(1,1)& 510 & +gm(1,1)*dgm(2,3)) 511 cm(1,2,9)=6*gm(2,2)**2*dgm(1,3)+3*gm(1,2)*gm(2,3)*dgm(2,2)+gm(2,2)& 512 & *(3*gm(2,3)*dgm(1,2)+12*gm(1,3)*dgm(2,2)+3*gm(1,2)*dgm(2,3)) 513 cm(1,3,9)=1.5d0*gm(2,3)**2*dgm(1,3)+gm(3,3)*(-4.5d0*gm(1,3)*dgm(2,2)& 514 & +12*gm(1,2)*dgm(2,3))+gm(2,3)*(12*gm(3,3)*dgm(1,2)+3*gm(1,3)& 515 & *dgm(2,3)+12*gm(1,2)*dgm(3,3))-4.5d0*gm(2,2)*(gm(3,3)*dgm(1,3)& 516 & +gm(1,3)*dgm(3,3)) 517 cm(1,4,9)=3*gm(2,3)**2*dgm(1,2)+12*gm(1,2)*gm(3,3)*dgm(2,2)+gm(2,3)& 518 & *(3*gm(2,2)*dgm(1,3)+3*gm(1,3)*dgm(2,2)+6*gm(1,2)*dgm(2,3))+gm(2,2)& 519 & *(12*gm(3,3)*dgm(1,2)+3*gm(1,3)*dgm(2,3)+12*gm(1,2)*dgm(3,3)) 520 cm(1,5,9)=-3*gm(2,3)**2*dgm(1,1)+15*gm(1,2)*gm(3,3)*dgm(1,2)-3*gm(1,3)& 521 & **2*dgm(2,2)-1.5d0*gm(1,1)*gm(3,3)*dgm(2,2)+9*gm(1,2)*gm(1,3)& 522 & *dgm(2,3)+gm(2,3)*(9*gm(1,3)*dgm(1,2)+9*gm(1,2)*dgm(1,3)-6*gm(1,1)& 523 & *dgm(2,3))+7.5d0*gm(1,2)**2*dgm(3,3)+gm(2,2)*(-1.5d0*gm(3,3)& 524 & *dgm(1,1)-6*gm(1,3)*dgm(1,3)-1.5d0*gm(1,1)*dgm(3,3)) 525 cm(1,6,9)=-4.5d0*gm(1,1)*gm(2,3)*dgm(2,2)+gm(1,2)*(3*gm(2,3)*dgm(1,2)& 526 & +12*gm(1,3)*dgm(2,2))+1.5d0*gm(1,2)**2*dgm(2,3)+gm(2,2)*(-4.5d0*gm(2,3)& 527 & *dgm(1,1)+12*gm(1,3)*dgm(1,2)+12*gm(1,2)*dgm(1,3)-4.5d0*gm(1,1)& 528 & *dgm(2,3)) 529 cm(1,7,9)=gm(2,2)*(6*gm(2,3)*dgm(2,2)+3*gm(2,2)*dgm(2,3)) 530 cm(1,8,9)=4.5d0*gm(2,3)**2*dgm(2,3)+7.5d0*gm(2,2)*gm(3,3)*dgm(2,3)& 531 & +7.5d0*gm(2,3)*(gm(3,3)*dgm(2,2)+gm(2,2)*dgm(3,3)) 532 cm(1,9,9)=3*gm(2,3)**2*dgm(2,2)+6*gm(2,2)*gm(2,3)*dgm(2,3)+gm(2,2)& 533 & *(12*gm(3,3)*dgm(2,2)+6*gm(2,2)*dgm(3,3)) 534 cm(1,10,9)=-1.5d0*gm(3,3)**2*dgm(2,2)+4.5d0*gm(2,3)**2*dgm(3,3)& 535 & +gm(3,3)*(9*gm(2,3)*dgm(2,3)-3*gm(2,2)*dgm(3,3)) 536 cm(1,1,10)=7.5d0*gm(1,3)**2*dgm(1,3)-1.5d0*gm(1,1)*gm(3,3)*dgm(1,3)& 537 & -1.5d0*gm(1,3)*(gm(3,3)*dgm(1,1)+gm(1,1)*dgm(3,3)) 538 cm(1,2,10)=7.5d0*gm(2,3)**2*dgm(1,3)+gm(3,3)*(-1.5d0*gm(1,3)*dgm(2,2)& 539 & -3*gm(1,2)*dgm(2,3))+gm(2,3)*(-3*gm(3,3)*dgm(1,2)+15*gm(1,3)& 540 & *dgm(2,3)-3*gm(1,2)*dgm(3,3))-1.5d0*gm(2,2)*(gm(3,3)*dgm(1,3)& 541 & +gm(1,3)*dgm(3,3)) 542 cm(1,3,10)=gm(3,3)*(3*gm(3,3)*dgm(1,3)+6*gm(1,3)*dgm(3,3)) 543 cm(1,4,10)=-3*gm(3,3)**2*dgm(1,2)+9*gm(1,3)*gm(2,3)*dgm(3,3)+gm(3,3)& 544 & *(9*gm(2,3)*dgm(1,3)+9*gm(1,3)*dgm(2,3)-6*gm(1,2)*dgm(3,3)) 545 cm(1,5,10)=-1.5d0*gm(3,3)**2*dgm(1,1)+4.5d0*gm(1,3)**2*dgm(3,3)& 546 & +gm(3,3)*(9*gm(1,3)*dgm(1,3)-3*gm(1,1)*dgm(3,3)) 547 cm(1,6,10)=7.5d0*gm(1,3)**2*dgm(2,3)+gm(3,3)*(-3*gm(1,2)*dgm(1,3)& 548 & -1.5d0*gm(1,1)*dgm(2,3))+gm(2,3)*(-1.5d0*gm(3,3)*dgm(1,1)+15*gm(1,3)& 549 & *dgm(1,3)-1.5d0*gm(1,1)*dgm(3,3))-3*gm(1,3)*(gm(3,3)*dgm(1,2)& 550 & +gm(1,2)*dgm(3,3)) 551 cm(1,7,10)=7.5d0*gm(2,3)**2*dgm(2,3)-1.5d0*gm(2,2)*gm(3,3)*dgm(2,3)& 552 & -1.5d0*gm(2,3)*(gm(3,3)*dgm(2,2)+gm(2,2)*dgm(3,3)) 553 cm(1,8,10)=gm(3,3)*(3*gm(3,3)*dgm(2,3)+6*gm(2,3)*dgm(3,3)) 554 cm(1,9,10)=-1.5d0*gm(3,3)**2*dgm(2,2)+4.5d0*gm(2,3)**2*dgm(3,3)& 555 & +gm(3,3)*(9*gm(2,3)*dgm(2,3)-3*gm(2,2)*dgm(3,3)) 556 cm(1,10,10)=3*gm(3,3)**2*dgm(3,3) 557 cm(2,1,2)=gm(1,1)*(6*gm(1,2)*dgm(1,1)+3*gm(1,1)*dgm(1,2)) 558 cm(2,2,2)=4.5d0*gm(1,2)**2*dgm(1,2)+7.5d0*gm(1,1)*gm(2,2)*dgm(1,2)& 559 & +7.5d0*gm(1,2)*(gm(2,2)*dgm(1,1)+gm(1,1)*dgm(2,2)) 560 cm(2,3,2)=1.5d0*gm(1,3)**2*dgm(1,2)+gm(1,1)*(-4.5d0*gm(3,3)*dgm(1,2)& 561 & +12*gm(2,3)*dgm(1,3))+gm(1,3)*(12*gm(2,3)*dgm(1,1)+3*gm(1,2)& 562 & *dgm(1,3)+12*gm(1,1)*dgm(2,3))-4.5d0*gm(1,2)*(gm(3,3)*dgm(1,1)& 563 & +gm(1,1)*dgm(3,3)) 564 cm(2,4,2)=3*gm(1,2)**2*dgm(1,3)+gm(1,1)*(3*gm(2,3)*dgm(1,2)+12*gm(2,2)& 565 & *dgm(1,3))+gm(1,3)*(12*gm(2,2)*dgm(1,1)+6*gm(1,2)*dgm(1,2)+12*gm(1,1)& 566 & *dgm(2,2))+3*gm(1,2)*(gm(2,3)*dgm(1,1)+gm(1,1)*dgm(2,3)) 567 cm(2,5,2)=3*gm(1,2)*(gm(1,3)*dgm(1,1)+gm(1,1)*dgm(1,3))+gm(1,1)& 568 & *(12*gm(2,3)*dgm(1,1)+3*gm(1,3)*dgm(1,2)+6*gm(1,1)*dgm(2,3)) 569 cm(2,6,2)=3*gm(1,2)**2*dgm(1,1)+6*gm(1,1)*gm(1,2)*dgm(1,2)+gm(1,1)& 570 & *(12*gm(2,2)*dgm(1,1)+6*gm(1,1)*dgm(2,2)) 571 cm(2,7,2)=-1.5d0*gm(2,2)**2*dgm(1,1)+4.5d0*gm(1,2)**2*dgm(2,2)& 572 & +gm(2,2)*(9*gm(1,2)*dgm(1,2)-3*gm(1,1)*dgm(2,2)) 573 cm(2,8,2)=-3*gm(2,3)**2*dgm(1,1)-6*gm(1,2)*gm(3,3)*dgm(1,2)+7.5d0*gm(1,3)& 574 & **2*dgm(2,2)-1.5d0*gm(1,1)*gm(3,3)*dgm(2,2)+9*gm(1,2)*gm(1,3)& 575 & *dgm(2,3)+gm(2,3)*(9*gm(1,3)*dgm(1,2)+9*gm(1,2)*dgm(1,3)-6*gm(1,1)& 576 & *dgm(2,3))-3*gm(1,2)**2*dgm(3,3)+gm(2,2)*(-1.5d0*gm(3,3)*dgm(1,1)& 577 & +15*gm(1,3)*dgm(1,3)-1.5d0*gm(1,1)*dgm(3,3)) 578 cm(2,9,2)=-4.5d0*gm(1,1)*gm(2,3)*dgm(2,2)+gm(1,2)*(3*gm(2,3)*dgm(1,2)& 579 & +12*gm(1,3)*dgm(2,2))+1.5d0*gm(1,2)**2*dgm(2,3)+gm(2,2)*(-4.5d0*gm(2,3)& 580 & *dgm(1,1)+12*gm(1,3)*dgm(1,2)+12*gm(1,2)*dgm(1,3)-4.5d0*gm(1,1)& 581 & *dgm(2,3)) 582 cm(2,10,2)=7.5d0*gm(1,3)**2*dgm(2,3)+gm(3,3)*(-3*gm(1,2)*dgm(1,3)& 583 & -1.5d0*gm(1,1)*dgm(2,3))+gm(2,3)*(-1.5d0*gm(3,3)*dgm(1,1)+15*gm(1,3)& 584 & *dgm(1,3)-1.5d0*gm(1,1)*dgm(3,3))-3*gm(1,3)*(gm(3,3)*dgm(1,2)& 585 & +gm(1,2)*dgm(3,3)) 586 cm(2,1,4)=gm(1,1)*(6*gm(1,3)*dgm(1,1)+3*gm(1,1)*dgm(1,3)) 587 cm(2,2,4)=1.5d0*gm(1,2)**2*dgm(1,3)+gm(1,1)*(12*gm(2,3)*dgm(1,2)& 588 & -4.5d0*gm(2,2)*dgm(1,3))+gm(1,3)*(-4.5d0*gm(2,2)*dgm(1,1)+3*gm(1,2)& 589 & *dgm(1,2)-4.5d0*gm(1,1)*dgm(2,2))+12*gm(1,2)*(gm(2,3)*dgm(1,1)& 590 & +gm(1,1)*dgm(2,3)) 591 cm(2,3,4)=4.5d0*gm(1,3)**2*dgm(1,3)+7.5d0*gm(1,1)*gm(3,3)*dgm(1,3)& 592 & +7.5d0*gm(1,3)*(gm(3,3)*dgm(1,1)+gm(1,1)*dgm(3,3)) 593 cm(2,4,4)=3*gm(1,3)**2*dgm(1,2)+gm(1,1)*(12*gm(3,3)*dgm(1,2)+3*gm(2,3)& 594 & *dgm(1,3))+gm(1,3)*(3*gm(2,3)*dgm(1,1)+6*gm(1,2)*dgm(1,3)+3*gm(1,1)& 595 & *dgm(2,3))+12*gm(1,2)*(gm(3,3)*dgm(1,1)+gm(1,1)*dgm(3,3)) 596 cm(2,5,4)=3*gm(1,3)**2*dgm(1,1)+6*gm(1,1)*gm(1,3)*dgm(1,3)+gm(1,1)& 597 & *(12*gm(3,3)*dgm(1,1)+6*gm(1,1)*dgm(3,3)) 598 cm(2,6,4)=3*gm(1,2)*(gm(1,3)*dgm(1,1)+gm(1,1)*dgm(1,3))+gm(1,1)& 599 & *(12*gm(2,3)*dgm(1,1)+3*gm(1,3)*dgm(1,2)+6*gm(1,1)*dgm(2,3)) 600 cm(2,7,4)=-1.5d0*gm(1,1)*gm(2,3)*dgm(2,2)+gm(1,2)*(15*gm(2,3)& 601 & *dgm(1,2)-3*gm(1,3)*dgm(2,2))+7.5d0*gm(1,2)**2*dgm(2,3)+gm(2,2)& 602 & *(-1.5d0*gm(2,3)*dgm(1,1)-3*gm(1,3)*dgm(1,2)-3*gm(1,2)*dgm(1,3)& 603 & -1.5d0*gm(1,1)*dgm(2,3)) 604 cm(2,8,4)=1.5d0*gm(1,3)**2*dgm(2,3)+gm(3,3)*(12*gm(1,2)*dgm(1,3)& 605 & -4.5d0*gm(1,1)*dgm(2,3))+gm(2,3)*(-4.5d0*gm(3,3)*dgm(1,1)+3*gm(1,3)& 606 & *dgm(1,3)-4.5d0*gm(1,1)*dgm(3,3))+12*gm(1,3)*(gm(3,3)*dgm(1,2)& 607 & +gm(1,2)*dgm(3,3)) 608 cm(2,9,4)=-3*gm(2,3)**2*dgm(1,1)+15*gm(1,2)*gm(3,3)*dgm(1,2)-3*gm(1,3)& 609 & **2*dgm(2,2)-1.5d0*gm(1,1)*gm(3,3)*dgm(2,2)+9*gm(1,2)*gm(1,3)& 610 & *dgm(2,3)+gm(2,3)*(9*gm(1,3)*dgm(1,2)+9*gm(1,2)*dgm(1,3)-6*gm(1,1)& 611 & *dgm(2,3))+7.5d0*gm(1,2)**2*dgm(3,3)+gm(2,2)*(-1.5d0*gm(3,3)& 612 & *dgm(1,1)-6*gm(1,3)*dgm(1,3)-1.5d0*gm(1,1)*dgm(3,3)) 613 cm(2,10,4)=-1.5d0*gm(3,3)**2*dgm(1,1)+4.5d0*gm(1,3)**2*dgm(3,3)& 614 & +gm(3,3)*(9*gm(1,3)*dgm(1,3)-3*gm(1,1)*dgm(3,3)) 615 cm(2,1,6)=3*gm(1,1)**2*dgm(1,1) 616 cm(2,2,6)=4.5d0*gm(1,2)**2*dgm(1,1)+9*gm(1,1)*gm(1,2)*dgm(1,2)& 617 & +gm(1,1)*(-3*gm(2,2)*dgm(1,1)-1.5d0*gm(1,1)*dgm(2,2)) 618 cm(2,3,6)=4.5d0*gm(1,3)**2*dgm(1,1)+9*gm(1,1)*gm(1,3)*dgm(1,3)& 619 & +gm(1,1)*(-3*gm(3,3)*dgm(1,1)-1.5d0*gm(1,1)*dgm(3,3)) 620 cm(2,4,6)=9*gm(1,2)*(gm(1,3)*dgm(1,1)+gm(1,1)*dgm(1,3))+gm(1,1)& 621 & *(-6*gm(2,3)*dgm(1,1)+9*gm(1,3)*dgm(1,2)-3*gm(1,1)*dgm(2,3)) 622 cm(2,5,6)=gm(1,1)*(6*gm(1,3)*dgm(1,1)+3*gm(1,1)*dgm(1,3)) 623 cm(2,6,6)=gm(1,1)*(6*gm(1,2)*dgm(1,1)+3*gm(1,1)*dgm(1,2)) 624 cm(2,7,6)=7.5d0*gm(1,2)**2*dgm(1,2)-1.5d0*gm(1,1)*gm(2,2)*dgm(1,2)& 625 & -1.5d0*gm(1,2)*(gm(2,2)*dgm(1,1)+gm(1,1)*dgm(2,2)) 626 cm(2,8,6)=7.5d0*gm(1,3)**2*dgm(1,2)+gm(1,1)*(-1.5d0*gm(3,3)*dgm(1,2)& 627 & -3*gm(2,3)*dgm(1,3))+gm(1,3)*(-3*gm(2,3)*dgm(1,1)+15*gm(1,2)& 628 & *dgm(1,3)-3*gm(1,1)*dgm(2,3))-1.5d0*gm(1,2)*(gm(3,3)*dgm(1,1)& 629 & +gm(1,1)*dgm(3,3)) 630 cm(2,9,6)=7.5d0*gm(1,2)**2*dgm(1,3)+gm(1,1)*(-3*gm(2,3)*dgm(1,2)& 631 & -1.5d0*gm(2,2)*dgm(1,3))+gm(1,3)*(-1.5d0*gm(2,2)*dgm(1,1)+15*gm(1,2)& 632 & *dgm(1,2)-1.5d0*gm(1,1)*dgm(2,2))-3*gm(1,2)*(gm(2,3)*dgm(1,1)& 633 & +gm(1,1)*dgm(2,3)) 634 cm(2,10,6)=7.5d0*gm(1,3)**2*dgm(1,3)-1.5d0*gm(1,1)*gm(3,3)*dgm(1,3)& 635 & -1.5d0*gm(1,3)*(gm(3,3)*dgm(1,1)+gm(1,1)*dgm(3,3)) 636 cm(2,1,7)=4.5d0*gm(1,2)**2*dgm(1,1)+9*gm(1,1)*gm(1,2)*dgm(1,2)& 637 & +gm(1,1)*(-3*gm(2,2)*dgm(1,1)-1.5d0*gm(1,1)*dgm(2,2)) 638 cm(2,2,7)=6*gm(2,2)**2*dgm(1,1)+3*gm(1,2)**2*dgm(2,2)+gm(2,2)& 639 & *(6*gm(1,2)*dgm(1,2)+12*gm(1,1)*dgm(2,2)) 640 cm(2,3,7)=7.5d0*gm(2,3)**2*dgm(1,1)-6*gm(1,2)*gm(3,3)*dgm(1,2)& 641 & -3*gm(1,3)**2*dgm(2,2)-1.5d0*gm(1,1)*gm(3,3)*dgm(2,2)+9*gm(1,2)& 642 & *gm(1,3)*dgm(2,3)+gm(2,3)*(9*gm(1,3)*dgm(1,2)+9*gm(1,2)*dgm(1,3)& 643 & +15*gm(1,1)*dgm(2,3))-3*gm(1,2)**2*dgm(3,3)+gm(2,2)*(-1.5d0*gm(3,3)& 644 & *dgm(1,1)-6*gm(1,3)*dgm(1,3)-1.5d0*gm(1,1)*dgm(3,3)) 645 cm(2,4,7)=12*gm(1,1)*gm(2,3)*dgm(2,2)+gm(1,2)*(6*gm(2,3)*dgm(1,2)& 646 & +3*gm(1,3)*dgm(2,2))+3*gm(1,2)**2*dgm(2,3)+gm(2,2)*(12*gm(2,3)& 647 & *dgm(1,1)+3*gm(1,3)*dgm(1,2)+3*gm(1,2)*dgm(1,3)+12*gm(1,1)*dgm(2,3)) 648 cm(2,5,7)=1.5d0*gm(1,2)**2*dgm(1,3)+gm(1,1)*(12*gm(2,3)*dgm(1,2)& 649 & -4.5d0*gm(2,2)*dgm(1,3))+gm(1,3)*(-4.5d0*gm(2,2)*dgm(1,1)+3*gm(1,2)& 650 & *dgm(1,2)-4.5d0*gm(1,1)*dgm(2,2))+12*gm(1,2)*(gm(2,3)*dgm(1,1)& 651 & +gm(1,1)*dgm(2,3)) 652 cm(2,6,7)=4.5d0*gm(1,2)**2*dgm(1,2)+7.5d0*gm(1,1)*gm(2,2)*dgm(1,2)& 653 & +7.5d0*gm(1,2)*(gm(2,2)*dgm(1,1)+gm(1,1)*dgm(2,2)) 654 cm(2,7,7)=gm(2,2)*(3*gm(2,2)*dgm(1,2)+6*gm(1,2)*dgm(2,2)) 655 cm(2,8,7)=1.5d0*gm(2,3)**2*dgm(1,2)-4.5d0*gm(1,2)*gm(3,3)*dgm(2,2)& 656 & +gm(2,3)*(12*gm(2,2)*dgm(1,3)+12*gm(1,3)*dgm(2,2)+3*gm(1,2)*dgm(2,3))& 657 & +gm(2,2)*(-4.5d0*gm(3,3)*dgm(1,2)+12*gm(1,3)*dgm(2,3)-4.5d0*gm(1,2)& 658 & *dgm(3,3)) 659 cm(2,9,7)=6*gm(2,2)**2*dgm(1,3)+3*gm(1,2)*gm(2,3)*dgm(2,2)+gm(2,2)& 660 & *(3*gm(2,3)*dgm(1,2)+12*gm(1,3)*dgm(2,2)+3*gm(1,2)*dgm(2,3)) 661 cm(2,10,7)=7.5d0*gm(2,3)**2*dgm(1,3)+gm(3,3)*(-1.5d0*gm(1,3)*dgm(2,2)& 662 & -3*gm(1,2)*dgm(2,3))+gm(2,3)*(-3*gm(3,3)*dgm(1,2)+15*gm(1,3)& 663 & *dgm(2,3)-3*gm(1,2)*dgm(3,3))-1.5d0*gm(2,2)*(gm(3,3)*dgm(1,3)& 664 & +gm(1,3)*dgm(3,3)) 665 cm(2,1,8)=4.5d0*gm(1,3)**2*dgm(1,1)+9*gm(1,1)*gm(1,3)*dgm(1,3)& 666 & +gm(1,1)*(-3*gm(3,3)*dgm(1,1)-1.5d0*gm(1,1)*dgm(3,3)) 667 cm(2,2,8)=7.5d0*gm(2,3)**2*dgm(1,1)-6*gm(1,2)*gm(3,3)*dgm(1,2)& 668 & -3*gm(1,3)**2*dgm(2,2)-1.5d0*gm(1,1)*gm(3,3)*dgm(2,2)+9*gm(1,2)& 669 & *gm(1,3)*dgm(2,3)+gm(2,3)*(9*gm(1,3)*dgm(1,2)+9*gm(1,2)*dgm(1,3)& 670 & +15*gm(1,1)*dgm(2,3))-3*gm(1,2)**2*dgm(3,3)+gm(2,2)*(-1.5d0*gm(3,3)& 671 & *dgm(1,1)-6*gm(1,3)*dgm(1,3)-1.5d0*gm(1,1)*dgm(3,3)) 672 cm(2,3,8)=6*gm(3,3)**2*dgm(1,1)+3*gm(1,3)**2*dgm(3,3)+gm(3,3)& 673 & *(6*gm(1,3)*dgm(1,3)+12*gm(1,1)*dgm(3,3)) 674 cm(2,4,8)=3*gm(1,3)**2*dgm(2,3)+gm(3,3)*(3*gm(1,2)*dgm(1,3)+12*gm(1,1)& 675 & *dgm(2,3))+gm(2,3)*(12*gm(3,3)*dgm(1,1)+6*gm(1,3)*dgm(1,3)+12*gm(1,1)& 676 & *dgm(3,3))+3*gm(1,3)*(gm(3,3)*dgm(1,2)+gm(1,2)*dgm(3,3)) 677 cm(2,5,8)=4.5d0*gm(1,3)**2*dgm(1,3)+7.5d0*gm(1,1)*gm(3,3)*dgm(1,3)& 678 & +7.5d0*gm(1,3)*(gm(3,3)*dgm(1,1)+gm(1,1)*dgm(3,3)) 679 cm(2,6,8)=1.5d0*gm(1,3)**2*dgm(1,2)+gm(1,1)*(-4.5d0*gm(3,3)*dgm(1,2)& 680 & +12*gm(2,3)*dgm(1,3))+gm(1,3)*(12*gm(2,3)*dgm(1,1)+3*gm(1,2)& 681 & *dgm(1,3)+12*gm(1,1)*dgm(2,3))-4.5d0*gm(1,2)*(gm(3,3)*dgm(1,1)& 682 & +gm(1,1)*dgm(3,3)) 683 cm(2,7,8)=7.5d0*gm(2,3)**2*dgm(1,2)-1.5d0*gm(1,2)*gm(3,3)*dgm(2,2)& 684 & +gm(2,3)*(-3*gm(2,2)*dgm(1,3)-3*gm(1,3)*dgm(2,2)+15*gm(1,2)*dgm(2,3))& 685 & +gm(2,2)*(-1.5d0*gm(3,3)*dgm(1,2)-3*gm(1,3)*dgm(2,3)-1.5d0*gm(1,2)& 686 & *dgm(3,3)) 687 cm(2,8,8)=6*gm(3,3)**2*dgm(1,2)+3*gm(1,3)*gm(2,3)*dgm(3,3)+gm(3,3)& 688 & *(3*gm(2,3)*dgm(1,3)+3*gm(1,3)*dgm(2,3)+12*gm(1,2)*dgm(3,3)) 689 cm(2,9,8)=1.5d0*gm(2,3)**2*dgm(1,3)+gm(3,3)*(-4.5d0*gm(1,3)*dgm(2,2)& 690 & +12*gm(1,2)*dgm(2,3))+gm(2,3)*(12*gm(3,3)*dgm(1,2)+3*gm(1,3)& 691 & *dgm(2,3)+12*gm(1,2)*dgm(3,3))-4.5d0*gm(2,2)*(gm(3,3)*dgm(1,3)& 692 & +gm(1,3)*dgm(3,3)) 693 cm(2,10,8)=gm(3,3)*(3*gm(3,3)*dgm(1,3)+6*gm(1,3)*dgm(3,3)) 694 cm(2,1,9)=9*gm(1,2)*(gm(1,3)*dgm(1,1)+gm(1,1)*dgm(1,3))+gm(1,1)& 695 & *(-6*gm(2,3)*dgm(1,1)+9*gm(1,3)*dgm(1,2)-3*gm(1,1)*dgm(2,3)) 696 cm(2,2,9)=12*gm(1,1)*gm(2,3)*dgm(2,2)+gm(1,2)*(6*gm(2,3)*dgm(1,2)& 697 & +3*gm(1,3)*dgm(2,2))+3*gm(1,2)**2*dgm(2,3)+gm(2,2)*(12*gm(2,3)& 698 & *dgm(1,1)+3*gm(1,3)*dgm(1,2)+3*gm(1,2)*dgm(1,3)+12*gm(1,1)*dgm(2,3)) 699 cm(2,3,9)=3*gm(1,3)**2*dgm(2,3)+gm(3,3)*(3*gm(1,2)*dgm(1,3)+12*gm(1,1)& 700 & *dgm(2,3))+gm(2,3)*(12*gm(3,3)*dgm(1,1)+6*gm(1,3)*dgm(1,3)+12*gm(1,1)& 701 & *dgm(3,3))+3*gm(1,3)*(gm(3,3)*dgm(1,2)+gm(1,2)*dgm(3,3)) 702 cm(2,4,9)=9*gm(2,3)**2*dgm(1,1)+18*gm(1,2)*gm(3,3)*dgm(1,2)+9*gm(1,3)& 703 & **2*dgm(2,2)+15*gm(1,1)*gm(3,3)*dgm(2,2)-6*gm(1,2)*gm(1,3)*dgm(2,3)& 704 & +gm(2,3)*(-6*gm(1,3)*dgm(1,2)-6*gm(1,2)*dgm(1,3)+18*gm(1,1)*dgm(2,3))& 705 & +9*gm(1,2)**2*dgm(3,3)+gm(2,2)*(15*gm(3,3)*dgm(1,1)+18*gm(1,3)& 706 & *dgm(1,3)+15*gm(1,1)*dgm(3,3)) 707 cm(2,5,9)=3*gm(1,3)**2*dgm(1,2)+gm(1,1)*(12*gm(3,3)*dgm(1,2)+3*gm(2,3)& 708 & *dgm(1,3))+gm(1,3)*(3*gm(2,3)*dgm(1,1)+6*gm(1,2)*dgm(1,3)+3*gm(1,1)& 709 & *dgm(2,3))+12*gm(1,2)*(gm(3,3)*dgm(1,1)+gm(1,1)*dgm(3,3)) 710 cm(2,6,9)=3*gm(1,2)**2*dgm(1,3)+gm(1,1)*(3*gm(2,3)*dgm(1,2)+12*gm(2,2)& 711 & *dgm(1,3))+gm(1,3)*(12*gm(2,2)*dgm(1,1)+6*gm(1,2)*dgm(1,2)+12*gm(1,1)& 712 & *dgm(2,2))+3*gm(1,2)*(gm(2,3)*dgm(1,1)+gm(1,1)*dgm(2,3)) 713 cm(2,7,9)=-3*gm(2,2)**2*dgm(1,3)+9*gm(1,2)*gm(2,3)*dgm(2,2)+gm(2,2)& 714 & *(9*gm(2,3)*dgm(1,2)-6*gm(1,3)*dgm(2,2)+9*gm(1,2)*dgm(2,3)) 715 cm(2,8,9)=3*gm(2,3)**2*dgm(1,3)+gm(3,3)*(12*gm(1,3)*dgm(2,2)+3*gm(1,2)& 716 & *dgm(2,3))+gm(2,3)*(3*gm(3,3)*dgm(1,2)+6*gm(1,3)*dgm(2,3)+3*gm(1,2)& 717 & *dgm(3,3))+12*gm(2,2)*(gm(3,3)*dgm(1,3)+gm(1,3)*dgm(3,3)) 718 cm(2,9,9)=3*gm(2,3)**2*dgm(1,2)+12*gm(1,2)*gm(3,3)*dgm(2,2)+gm(2,3)& 719 & *(3*gm(2,2)*dgm(1,3)+3*gm(1,3)*dgm(2,2)+6*gm(1,2)*dgm(2,3))+gm(2,2)& 720 & *(12*gm(3,3)*dgm(1,2)+3*gm(1,3)*dgm(2,3)+12*gm(1,2)*dgm(3,3)) 721 cm(2,10,9)=-3*gm(3,3)**2*dgm(1,2)+9*gm(1,3)*gm(2,3)*dgm(3,3)+gm(3,3)& 722 & *(9*gm(2,3)*dgm(1,3)+9*gm(1,3)*dgm(2,3)-6*gm(1,2)*dgm(3,3)) 723 cm(2,1,11)=7.5d0*gm(1,2)**2*dgm(1,2)-1.5d0*gm(1,1)*gm(2,2)*dgm(1,2)& 724 & -1.5d0*gm(1,2)*(gm(2,2)*dgm(1,1)+gm(1,1)*dgm(2,2)) 725 cm(2,2,11)=gm(2,2)*(3*gm(2,2)*dgm(1,2)+6*gm(1,2)*dgm(2,2)) 726 cm(2,3,11)=7.5d0*gm(2,3)**2*dgm(1,2)-1.5d0*gm(1,2)*gm(3,3)*dgm(2,2)& 727 & +gm(2,3)*(-3*gm(2,2)*dgm(1,3)-3*gm(1,3)*dgm(2,2)+15*gm(1,2)*dgm(2,3))& 728 & +gm(2,2)*(-1.5d0*gm(3,3)*dgm(1,2)-3*gm(1,3)*dgm(2,3)-1.5d0*gm(1,2)& 729 & *dgm(3,3)) 730 cm(2,4,11)=-3*gm(2,2)**2*dgm(1,3)+9*gm(1,2)*gm(2,3)*dgm(2,2)+gm(2,2)& 731 & *(9*gm(2,3)*dgm(1,2)-6*gm(1,3)*dgm(2,2)+9*gm(1,2)*dgm(2,3)) 732 cm(2,5,11)=-1.5d0*gm(1,1)*gm(2,3)*dgm(2,2)+gm(1,2)*(15*gm(2,3)& 733 & *dgm(1,2)-3*gm(1,3)*dgm(2,2))+7.5d0*gm(1,2)**2*dgm(2,3)+gm(2,2)& 734 & *(-1.5d0*gm(2,3)*dgm(1,1)-3*gm(1,3)*dgm(1,2)-3*gm(1,2)*dgm(1,3)& 735 & -1.5d0*gm(1,1)*dgm(2,3)) 736 cm(2,6,11)=-1.5d0*gm(2,2)**2*dgm(1,1)+4.5d0*gm(1,2)**2*dgm(2,2)& 737 & +gm(2,2)*(9*gm(1,2)*dgm(1,2)-3*gm(1,1)*dgm(2,2)) 738 cm(2,7,11)=3*gm(2,2)**2*dgm(2,2) 739 cm(2,8,11)=4.5d0*gm(2,3)**2*dgm(2,2)+9*gm(2,2)*gm(2,3)*dgm(2,3)& 740 & +gm(2,2)*(-3*gm(3,3)*dgm(2,2)-1.5d0*gm(2,2)*dgm(3,3)) 741 cm(2,9,11)=gm(2,2)*(6*gm(2,3)*dgm(2,2)+3*gm(2,2)*dgm(2,3)) 742 cm(2,10,11)=7.5d0*gm(2,3)**2*dgm(2,3)-1.5d0*gm(2,2)*gm(3,3)*dgm(2,3)& 743 & -1.5d0*gm(2,3)*(gm(3,3)*dgm(2,2)+gm(2,2)*dgm(3,3)) 744 cm(2,1,12)=7.5d0*gm(1,3)**2*dgm(1,2)+gm(1,1)*(-1.5d0*gm(3,3)*dgm(1,2)& 745 & -3*gm(2,3)*dgm(1,3))+gm(1,3)*(-3*gm(2,3)*dgm(1,1)+15*gm(1,2)& 746 & *dgm(1,3)-3*gm(1,1)*dgm(2,3))-1.5d0*gm(1,2)*(gm(3,3)*dgm(1,1)& 747 & +gm(1,1)*dgm(3,3)) 748 cm(2,2,12)=1.5d0*gm(2,3)**2*dgm(1,2)-4.5d0*gm(1,2)*gm(3,3)*dgm(2,2)& 749 & +gm(2,3)*(12*gm(2,2)*dgm(1,3)+12*gm(1,3)*dgm(2,2)+3*gm(1,2)*dgm(2,3))& 750 & +gm(2,2)*(-4.5d0*gm(3,3)*dgm(1,2)+12*gm(1,3)*dgm(2,3)-4.5d0*gm(1,2)& 751 & *dgm(3,3)) 752 cm(2,3,12)=6*gm(3,3)**2*dgm(1,2)+3*gm(1,3)*gm(2,3)*dgm(3,3)+gm(3,3)& 753 & *(3*gm(2,3)*dgm(1,3)+3*gm(1,3)*dgm(2,3)+12*gm(1,2)*dgm(3,3)) 754 cm(2,4,12)=3*gm(2,3)**2*dgm(1,3)+gm(3,3)*(12*gm(1,3)*dgm(2,2)& 755 & +3*gm(1,2)*dgm(2,3))+gm(2,3)*(3*gm(3,3)*dgm(1,2)+6*gm(1,3)*dgm(2,3)& 756 & +3*gm(1,2)*dgm(3,3))+12*gm(2,2)*(gm(3,3)*dgm(1,3)+gm(1,3)*dgm(3,3)) 757 cm(2,5,12)=1.5d0*gm(1,3)**2*dgm(2,3)+gm(3,3)*(12*gm(1,2)*dgm(1,3)& 758 & -4.5d0*gm(1,1)*dgm(2,3))+gm(2,3)*(-4.5d0*gm(3,3)*dgm(1,1)+3*gm(1,3)& 759 & *dgm(1,3)-4.5d0*gm(1,1)*dgm(3,3))+12*gm(1,3)*(gm(3,3)*dgm(1,2)& 760 & +gm(1,2)*dgm(3,3)) 761 cm(2,6,12)=-3*gm(2,3)**2*dgm(1,1)-6*gm(1,2)*gm(3,3)*dgm(1,2)+7.5d0*gm(1,3)& 762 & **2*dgm(2,2)-1.5d0*gm(1,1)*gm(3,3)*dgm(2,2)+9*gm(1,2)*gm(1,3)& 763 & *dgm(2,3)+gm(2,3)*(9*gm(1,3)*dgm(1,2)+9*gm(1,2)*dgm(1,3)-6*gm(1,1)& 764 & *dgm(2,3))-3*gm(1,2)**2*dgm(3,3)+gm(2,2)*(-1.5d0*gm(3,3)*dgm(1,1)& 765 & +15*gm(1,3)*dgm(1,3)-1.5d0*gm(1,1)*dgm(3,3)) 766 cm(2,7,12)=4.5d0*gm(2,3)**2*dgm(2,2)+9*gm(2,2)*gm(2,3)*dgm(2,3)& 767 & +gm(2,2)*(-3*gm(3,3)*dgm(2,2)-1.5d0*gm(2,2)*dgm(3,3)) 768 cm(2,8,12)=6*gm(3,3)**2*dgm(2,2)+3*gm(2,3)**2*dgm(3,3)+gm(3,3)& 769 & *(6*gm(2,3)*dgm(2,3)+12*gm(2,2)*dgm(3,3)) 770 cm(2,9,12)=4.5d0*gm(2,3)**2*dgm(2,3)+7.5d0*gm(2,2)*gm(3,3)*dgm(2,3)& 771 & +7.5d0*gm(2,3)*(gm(3,3)*dgm(2,2)+gm(2,2)*dgm(3,3)) 772 cm(2,10,12)=gm(3,3)*(3*gm(3,3)*dgm(2,3)+6*gm(2,3)*dgm(3,3)) 773 cm(2,1,13)=7.5d0*gm(1,2)**2*dgm(1,3)+gm(1,1)*(-3*gm(2,3)*dgm(1,2)& 774 & -1.5d0*gm(2,2)*dgm(1,3))+gm(1,3)*(-1.5d0*gm(2,2)*dgm(1,1)+15*gm(1,2)& 775 & *dgm(1,2)-1.5d0*gm(1,1)*dgm(2,2))-3*gm(1,2)*(gm(2,3)*dgm(1,1)& 776 & +gm(1,1)*dgm(2,3)) 777 cm(2,2,13)=6*gm(2,2)**2*dgm(1,3)+3*gm(1,2)*gm(2,3)*dgm(2,2)+gm(2,2)& 778 & *(3*gm(2,3)*dgm(1,2)+12*gm(1,3)*dgm(2,2)+3*gm(1,2)*dgm(2,3)) 779 cm(2,3,13)=1.5d0*gm(2,3)**2*dgm(1,3)+gm(3,3)*(-4.5d0*gm(1,3)*dgm(2,2)& 780 & +12*gm(1,2)*dgm(2,3))+gm(2,3)*(12*gm(3,3)*dgm(1,2)+3*gm(1,3)& 781 & *dgm(2,3)+12*gm(1,2)*dgm(3,3))-4.5d0*gm(2,2)*(gm(3,3)*dgm(1,3)& 782 & +gm(1,3)*dgm(3,3)) 783 cm(2,4,13)=3*gm(2,3)**2*dgm(1,2)+12*gm(1,2)*gm(3,3)*dgm(2,2)+gm(2,3)& 784 & *(3*gm(2,2)*dgm(1,3)+3*gm(1,3)*dgm(2,2)+6*gm(1,2)*dgm(2,3))+gm(2,2)& 785 & *(12*gm(3,3)*dgm(1,2)+3*gm(1,3)*dgm(2,3)+12*gm(1,2)*dgm(3,3)) 786 cm(2,5,13)=-3*gm(2,3)**2*dgm(1,1)+15*gm(1,2)*gm(3,3)*dgm(1,2)& 787 & -3*gm(1,3)**2*dgm(2,2)-1.5d0*gm(1,1)*gm(3,3)*dgm(2,2)+9*gm(1,2)& 788 & *gm(1,3)*dgm(2,3)+gm(2,3)*(9*gm(1,3)*dgm(1,2)+9*gm(1,2)*dgm(1,3)& 789 & -6*gm(1,1)*dgm(2,3))+7.5d0*gm(1,2)**2*dgm(3,3)+gm(2,2)*(-1.5d0*gm(3,3)& 790 & *dgm(1,1)-6*gm(1,3)*dgm(1,3)-1.5d0*gm(1,1)*dgm(3,3)) 791 cm(2,6,13)=-4.5d0*gm(1,1)*gm(2,3)*dgm(2,2)+gm(1,2)*(3*gm(2,3)& 792 & *dgm(1,2)+12*gm(1,3)*dgm(2,2))+1.5d0*gm(1,2)**2*dgm(2,3)+gm(2,2)& 793 & *(-4.5d0*gm(2,3)*dgm(1,1)+12*gm(1,3)*dgm(1,2)+12*gm(1,2)*dgm(1,3)& 794 & -4.5d0*gm(1,1)*dgm(2,3)) 795 cm(2,7,13)=gm(2,2)*(6*gm(2,3)*dgm(2,2)+3*gm(2,2)*dgm(2,3)) 796 cm(2,8,13)=4.5d0*gm(2,3)**2*dgm(2,3)+7.5d0*gm(2,2)*gm(3,3)*dgm(2,3)& 797 & +7.5d0*gm(2,3)*(gm(3,3)*dgm(2,2)+gm(2,2)*dgm(3,3)) 798 cm(2,9,13)=3*gm(2,3)**2*dgm(2,2)+6*gm(2,2)*gm(2,3)*dgm(2,3)+gm(2,2)& 799 & *(12*gm(3,3)*dgm(2,2)+6*gm(2,2)*dgm(3,3)) 800 cm(2,10,13)=-1.5d0*gm(3,3)**2*dgm(2,2)+4.5d0*gm(2,3)**2*dgm(3,3)& 801 & +gm(3,3)*(9*gm(2,3)*dgm(2,3)-3*gm(2,2)*dgm(3,3)) 802 cm(2,1,14)=7.5d0*gm(1,3)**2*dgm(1,3)-1.5d0*gm(1,1)*gm(3,3)*dgm(1,3)& 803 & -1.5d0*gm(1,3)*(gm(3,3)*dgm(1,1)+gm(1,1)*dgm(3,3)) 804 cm(2,2,14)=7.5d0*gm(2,3)**2*dgm(1,3)+gm(3,3)*(-1.5d0*gm(1,3)*dgm(2,2)& 805 & -3*gm(1,2)*dgm(2,3))+gm(2,3)*(-3*gm(3,3)*dgm(1,2)+15*gm(1,3)& 806 & *dgm(2,3)-3*gm(1,2)*dgm(3,3))-1.5d0*gm(2,2)*(gm(3,3)*dgm(1,3)& 807 & +gm(1,3)*dgm(3,3)) 808 cm(2,3,14)=gm(3,3)*(3*gm(3,3)*dgm(1,3)+6*gm(1,3)*dgm(3,3)) 809 cm(2,4,14)=-3*gm(3,3)**2*dgm(1,2)+9*gm(1,3)*gm(2,3)*dgm(3,3)+gm(3,3)& 810 & *(9*gm(2,3)*dgm(1,3)+9*gm(1,3)*dgm(2,3)-6*gm(1,2)*dgm(3,3)) 811 cm(2,5,14)=-1.5d0*gm(3,3)**2*dgm(1,1)+4.5d0*gm(1,3)**2*dgm(3,3)& 812 & +gm(3,3)*(9*gm(1,3)*dgm(1,3)-3*gm(1,1)*dgm(3,3)) 813 cm(2,6,14)=7.5d0*gm(1,3)**2*dgm(2,3)+gm(3,3)*(-3*gm(1,2)*dgm(1,3)& 814 & -1.5d0*gm(1,1)*dgm(2,3))+gm(2,3)*(-1.5d0*gm(3,3)*dgm(1,1)+15*gm(1,3)& 815 & *dgm(1,3)-1.5d0*gm(1,1)*dgm(3,3))-3*gm(1,3)*(gm(3,3)*dgm(1,2)& 816 & +gm(1,2)*dgm(3,3)) 817 cm(2,7,14)=7.5d0*gm(2,3)**2*dgm(2,3)-1.5d0*gm(2,2)*gm(3,3)*dgm(2,3)& 818 & -1.5d0*gm(2,3)*(gm(3,3)*dgm(2,2)+gm(2,2)*dgm(3,3)) 819 cm(2,8,14)=gm(3,3)*(3*gm(3,3)*dgm(2,3)+6*gm(2,3)*dgm(3,3)) 820 cm(2,9,14)=-1.5d0*gm(3,3)**2*dgm(2,2)+4.5d0*gm(2,3)**2*dgm(3,3)& 821 & +gm(3,3)*(9*gm(2,3)*dgm(2,3)-3*gm(2,2)*dgm(3,3)) 822 cm(2,10,14)=3*gm(3,3)**2*dgm(3,3) 823 cm(3,1,3)=gm(1,1)*(6*gm(1,3)*dgm(1,1)+3*gm(1,1)*dgm(1,3)) 824 cm(3,2,3)=1.5d0*gm(1,2)**2*dgm(1,3)+gm(1,1)*(12*gm(2,3)*dgm(1,2)& 825 & -4.5d0*gm(2,2)*dgm(1,3))+gm(1,3)*(-4.5d0*gm(2,2)*dgm(1,1)+3*gm(1,2)& 826 & *dgm(1,2)-4.5d0*gm(1,1)*dgm(2,2))+12*gm(1,2)*(gm(2,3)*dgm(1,1)& 827 & +gm(1,1)*dgm(2,3)) 828 cm(3,3,3)=4.5d0*gm(1,3)**2*dgm(1,3)+7.5d0*gm(1,1)*gm(3,3)*dgm(1,3)& 829 & +7.5d0*gm(1,3)*(gm(3,3)*dgm(1,1)+gm(1,1)*dgm(3,3)) 830 cm(3,4,3)=3*gm(1,3)**2*dgm(1,2)+gm(1,1)*(12*gm(3,3)*dgm(1,2)+3*gm(2,3)& 831 & *dgm(1,3))+gm(1,3)*(3*gm(2,3)*dgm(1,1)+6*gm(1,2)*dgm(1,3)+3*gm(1,1)& 832 & *dgm(2,3))+12*gm(1,2)*(gm(3,3)*dgm(1,1)+gm(1,1)*dgm(3,3)) 833 cm(3,5,3)=3*gm(1,3)**2*dgm(1,1)+6*gm(1,1)*gm(1,3)*dgm(1,3)+gm(1,1)& 834 & *(12*gm(3,3)*dgm(1,1)+6*gm(1,1)*dgm(3,3)) 835 cm(3,6,3)=3*gm(1,2)*(gm(1,3)*dgm(1,1)+gm(1,1)*dgm(1,3))+gm(1,1)& 836 & *(12*gm(2,3)*dgm(1,1)+3*gm(1,3)*dgm(1,2)+6*gm(1,1)*dgm(2,3)) 837 cm(3,7,3)=-1.5d0*gm(1,1)*gm(2,3)*dgm(2,2)+gm(1,2)*(15*gm(2,3)& 838 & *dgm(1,2)-3*gm(1,3)*dgm(2,2))+7.5d0*gm(1,2)**2*dgm(2,3)+gm(2,2)& 839 & *(-1.5d0*gm(2,3)*dgm(1,1)-3*gm(1,3)*dgm(1,2)-3*gm(1,2)*dgm(1,3)& 840 & -1.5d0*gm(1,1)*dgm(2,3)) 841 cm(3,8,3)=1.5d0*gm(1,3)**2*dgm(2,3)+gm(3,3)*(12*gm(1,2)*dgm(1,3)& 842 & -4.5d0*gm(1,1)*dgm(2,3))+gm(2,3)*(-4.5d0*gm(3,3)*dgm(1,1)+3*gm(1,3)& 843 & *dgm(1,3)-4.5d0*gm(1,1)*dgm(3,3))+12*gm(1,3)*(gm(3,3)*dgm(1,2)& 844 & +gm(1,2)*dgm(3,3)) 845 cm(3,9,3)=-3*gm(2,3)**2*dgm(1,1)+15*gm(1,2)*gm(3,3)*dgm(1,2)-3*gm(1,3)& 846 & **2*dgm(2,2)-1.5d0*gm(1,1)*gm(3,3)*dgm(2,2)+9*gm(1,2)*gm(1,3)& 847 & *dgm(2,3)+gm(2,3)*(9*gm(1,3)*dgm(1,2)+9*gm(1,2)*dgm(1,3)-6*gm(1,1)& 848 & *dgm(2,3))+7.5d0*gm(1,2)**2*dgm(3,3)+gm(2,2)*(-1.5d0*gm(3,3)& 849 & *dgm(1,1)-6*gm(1,3)*dgm(1,3)-1.5d0*gm(1,1)*dgm(3,3)) 850 cm(3,10,3)=-1.5d0*gm(3,3)**2*dgm(1,1)+4.5d0*gm(1,3)**2*dgm(3,3)& 851 & +gm(3,3)*(9*gm(1,3)*dgm(1,3)-3*gm(1,1)*dgm(3,3)) 852 cm(3,1,4)=gm(1,1)*(6*gm(1,2)*dgm(1,1)+3*gm(1,1)*dgm(1,2)) 853 cm(3,2,4)=4.5d0*gm(1,2)**2*dgm(1,2)+7.5d0*gm(1,1)*gm(2,2)*dgm(1,2)& 854 & +7.5d0*gm(1,2)*(gm(2,2)*dgm(1,1)+gm(1,1)*dgm(2,2)) 855 cm(3,3,4)=1.5d0*gm(1,3)**2*dgm(1,2)+gm(1,1)*(-4.5d0*gm(3,3)*dgm(1,2)& 856 & +12*gm(2,3)*dgm(1,3))+gm(1,3)*(12*gm(2,3)*dgm(1,1)+3*gm(1,2)& 857 & *dgm(1,3)+12*gm(1,1)*dgm(2,3))-4.5d0*gm(1,2)*(gm(3,3)*dgm(1,1)& 858 & +gm(1,1)*dgm(3,3)) 859 cm(3,4,4)=3*gm(1,2)**2*dgm(1,3)+gm(1,1)*(3*gm(2,3)*dgm(1,2)+12*gm(2,2)& 860 & *dgm(1,3))+gm(1,3)*(12*gm(2,2)*dgm(1,1)+6*gm(1,2)*dgm(1,2)+12*gm(1,1)& 861 & *dgm(2,2))+3*gm(1,2)*(gm(2,3)*dgm(1,1)+gm(1,1)*dgm(2,3)) 862 cm(3,5,4)=3*gm(1,2)*(gm(1,3)*dgm(1,1)+gm(1,1)*dgm(1,3))+gm(1,1)& 863 & *(12*gm(2,3)*dgm(1,1)+3*gm(1,3)*dgm(1,2)+6*gm(1,1)*dgm(2,3)) 864 cm(3,6,4)=3*gm(1,2)**2*dgm(1,1)+6*gm(1,1)*gm(1,2)*dgm(1,2)+gm(1,1)& 865 & *(12*gm(2,2)*dgm(1,1)+6*gm(1,1)*dgm(2,2)) 866 cm(3,7,4)=-1.5d0*gm(2,2)**2*dgm(1,1)+4.5d0*gm(1,2)**2*dgm(2,2)& 867 & +gm(2,2)*(9*gm(1,2)*dgm(1,2)-3*gm(1,1)*dgm(2,2)) 868 cm(3,8,4)=-3*gm(2,3)**2*dgm(1,1)-6*gm(1,2)*gm(3,3)*dgm(1,2)+7.5d0*gm(1,3)& 869 & **2*dgm(2,2)-1.5d0*gm(1,1)*gm(3,3)*dgm(2,2)+9*gm(1,2)*gm(1,3)& 870 & *dgm(2,3)+gm(2,3)*(9*gm(1,3)*dgm(1,2)+9*gm(1,2)*dgm(1,3)-6*gm(1,1)& 871 & *dgm(2,3))-3*gm(1,2)**2*dgm(3,3)+gm(2,2)*(-1.5d0*gm(3,3)*dgm(1,1)& 872 & +15*gm(1,3)*dgm(1,3)-1.5d0*gm(1,1)*dgm(3,3)) 873 cm(3,9,4)=-4.5d0*gm(1,1)*gm(2,3)*dgm(2,2)+gm(1,2)*(3*gm(2,3)*dgm(1,2)& 874 & +12*gm(1,3)*dgm(2,2))+1.5d0*gm(1,2)**2*dgm(2,3)+gm(2,2)*(-4.5d0*gm(2,3)& 875 & *dgm(1,1)+12*gm(1,3)*dgm(1,2)+12*gm(1,2)*dgm(1,3)-4.5d0*gm(1,1)& 876 & *dgm(2,3)) 877 cm(3,10,4)=7.5d0*gm(1,3)**2*dgm(2,3)+gm(3,3)*(-3*gm(1,2)*dgm(1,3)& 878 & -1.5d0*gm(1,1)*dgm(2,3))+gm(2,3)*(-1.5d0*gm(3,3)*dgm(1,1)+15*gm(1,3)& 879 & *dgm(1,3)-1.5d0*gm(1,1)*dgm(3,3))-3*gm(1,3)*(gm(3,3)*dgm(1,2)& 880 & +gm(1,2)*dgm(3,3)) 881 cm(3,1,5)=3*gm(1,1)**2*dgm(1,1) 882 cm(3,2,5)=4.5d0*gm(1,2)**2*dgm(1,1)+9*gm(1,1)*gm(1,2)*dgm(1,2)& 883 & +gm(1,1)*(-3*gm(2,2)*dgm(1,1)-1.5d0*gm(1,1)*dgm(2,2)) 884 cm(3,3,5)=4.5d0*gm(1,3)**2*dgm(1,1)+9*gm(1,1)*gm(1,3)*dgm(1,3)& 885 & +gm(1,1)*(-3*gm(3,3)*dgm(1,1)-1.5d0*gm(1,1)*dgm(3,3)) 886 cm(3,4,5)=9*gm(1,2)*(gm(1,3)*dgm(1,1)+gm(1,1)*dgm(1,3))+gm(1,1)& 887 & *(-6*gm(2,3)*dgm(1,1)+9*gm(1,3)*dgm(1,2)-3*gm(1,1)*dgm(2,3)) 888 cm(3,5,5)=gm(1,1)*(6*gm(1,3)*dgm(1,1)+3*gm(1,1)*dgm(1,3)) 889 cm(3,6,5)=gm(1,1)*(6*gm(1,2)*dgm(1,1)+3*gm(1,1)*dgm(1,2)) 890 cm(3,7,5)=7.5d0*gm(1,2)**2*dgm(1,2)-1.5d0*gm(1,1)*gm(2,2)*dgm(1,2)& 891 & -1.5d0*gm(1,2)*(gm(2,2)*dgm(1,1)+gm(1,1)*dgm(2,2)) 892 cm(3,8,5)=7.5d0*gm(1,3)**2*dgm(1,2)+gm(1,1)*(-1.5d0*gm(3,3)*dgm(1,2)& 893 & -3*gm(2,3)*dgm(1,3))+gm(1,3)*(-3*gm(2,3)*dgm(1,1)+15*gm(1,2)& 894 & *dgm(1,3)-3*gm(1,1)*dgm(2,3))-1.5d0*gm(1,2)*(gm(3,3)*dgm(1,1)& 895 & +gm(1,1)*dgm(3,3)) 896 cm(3,9,5)=7.5d0*gm(1,2)**2*dgm(1,3)+gm(1,1)*(-3*gm(2,3)*dgm(1,2)& 897 & -1.5d0*gm(2,2)*dgm(1,3))+gm(1,3)*(-1.5d0*gm(2,2)*dgm(1,1)+15*gm(1,2)& 898 & *dgm(1,2)-1.5d0*gm(1,1)*dgm(2,2))-3*gm(1,2)*(gm(2,3)*dgm(1,1)& 899 & +gm(1,1)*dgm(2,3)) 900 cm(3,10,5)=7.5d0*gm(1,3)**2*dgm(1,3)-1.5d0*gm(1,1)*gm(3,3)*dgm(1,3)& 901 & -1.5d0*gm(1,3)*(gm(3,3)*dgm(1,1)+gm(1,1)*dgm(3,3)) 902 cm(3,1,8)=9*gm(1,2)*(gm(1,3)*dgm(1,1)+gm(1,1)*dgm(1,3))+gm(1,1)& 903 & *(-6*gm(2,3)*dgm(1,1)+9*gm(1,3)*dgm(1,2)-3*gm(1,1)*dgm(2,3)) 904 cm(3,2,8)=12*gm(1,1)*gm(2,3)*dgm(2,2)+gm(1,2)*(6*gm(2,3)*dgm(1,2)& 905 & +3*gm(1,3)*dgm(2,2))+3*gm(1,2)**2*dgm(2,3)+gm(2,2)*(12*gm(2,3)& 906 & *dgm(1,1)+3*gm(1,3)*dgm(1,2)+3*gm(1,2)*dgm(1,3)+12*gm(1,1)*dgm(2,3)) 907 cm(3,3,8)=3*gm(1,3)**2*dgm(2,3)+gm(3,3)*(3*gm(1,2)*dgm(1,3)+12*gm(1,1)& 908 & *dgm(2,3))+gm(2,3)*(12*gm(3,3)*dgm(1,1)+6*gm(1,3)*dgm(1,3)+12*gm(1,1)& 909 & *dgm(3,3))+3*gm(1,3)*(gm(3,3)*dgm(1,2)+gm(1,2)*dgm(3,3)) 910 cm(3,4,8)=9*gm(2,3)**2*dgm(1,1)+18*gm(1,2)*gm(3,3)*dgm(1,2)+9*gm(1,3)& 911 & **2*dgm(2,2)+15*gm(1,1)*gm(3,3)*dgm(2,2)-6*gm(1,2)*gm(1,3)*dgm(2,3)& 912 & +gm(2,3)*(-6*gm(1,3)*dgm(1,2)-6*gm(1,2)*dgm(1,3)+18*gm(1,1)*dgm(2,3))& 913 & +9*gm(1,2)**2*dgm(3,3)+gm(2,2)*(15*gm(3,3)*dgm(1,1)+18*gm(1,3)& 914 & *dgm(1,3)+15*gm(1,1)*dgm(3,3)) 915 cm(3,5,8)=3*gm(1,3)**2*dgm(1,2)+gm(1,1)*(12*gm(3,3)*dgm(1,2)+3*gm(2,3)& 916 & *dgm(1,3))+gm(1,3)*(3*gm(2,3)*dgm(1,1)+6*gm(1,2)*dgm(1,3)+3*gm(1,1)& 917 & *dgm(2,3))+12*gm(1,2)*(gm(3,3)*dgm(1,1)+gm(1,1)*dgm(3,3)) 918 cm(3,6,8)=3*gm(1,2)**2*dgm(1,3)+gm(1,1)*(3*gm(2,3)*dgm(1,2)+12*gm(2,2)& 919 & *dgm(1,3))+gm(1,3)*(12*gm(2,2)*dgm(1,1)+6*gm(1,2)*dgm(1,2)+12*gm(1,1)& 920 & *dgm(2,2))+3*gm(1,2)*(gm(2,3)*dgm(1,1)+gm(1,1)*dgm(2,3)) 921 cm(3,7,8)=-3*gm(2,2)**2*dgm(1,3)+9*gm(1,2)*gm(2,3)*dgm(2,2)+gm(2,2)& 922 & *(9*gm(2,3)*dgm(1,2)-6*gm(1,3)*dgm(2,2)+9*gm(1,2)*dgm(2,3)) 923 cm(3,8,8)=3*gm(2,3)**2*dgm(1,3)+gm(3,3)*(12*gm(1,3)*dgm(2,2)+3*gm(1,2)& 924 & *dgm(2,3))+gm(2,3)*(3*gm(3,3)*dgm(1,2)+6*gm(1,3)*dgm(2,3)+3*gm(1,2)& 925 & *dgm(3,3))+12*gm(2,2)*(gm(3,3)*dgm(1,3)+gm(1,3)*dgm(3,3)) 926 cm(3,9,8)=3*gm(2,3)**2*dgm(1,2)+12*gm(1,2)*gm(3,3)*dgm(2,2)+gm(2,3)& 927 & *(3*gm(2,2)*dgm(1,3)+3*gm(1,3)*dgm(2,2)+6*gm(1,2)*dgm(2,3))+gm(2,2)& 928 & *(12*gm(3,3)*dgm(1,2)+3*gm(1,3)*dgm(2,3)+12*gm(1,2)*dgm(3,3)) 929 cm(3,10,8)=-3*gm(3,3)**2*dgm(1,2)+9*gm(1,3)*gm(2,3)*dgm(3,3)+gm(3,3)& 930 & *(9*gm(2,3)*dgm(1,3)+9*gm(1,3)*dgm(2,3)-6*gm(1,2)*dgm(3,3)) 931 cm(3,1,9)=4.5d0*gm(1,2)**2*dgm(1,1)+9*gm(1,1)*gm(1,2)*dgm(1,2)& 932 & +gm(1,1)*(-3*gm(2,2)*dgm(1,1)-1.5d0*gm(1,1)*dgm(2,2)) 933 cm(3,2,9)=6*gm(2,2)**2*dgm(1,1)+3*gm(1,2)**2*dgm(2,2)+gm(2,2)& 934 & *(6*gm(1,2)*dgm(1,2)+12*gm(1,1)*dgm(2,2)) 935 cm(3,3,9)=7.5d0*gm(2,3)**2*dgm(1,1)-6*gm(1,2)*gm(3,3)*dgm(1,2)& 936 & -3*gm(1,3)**2*dgm(2,2)-1.5d0*gm(1,1)*gm(3,3)*dgm(2,2)+9*gm(1,2)& 937 & *gm(1,3)*dgm(2,3)+gm(2,3)*(9*gm(1,3)*dgm(1,2)+9*gm(1,2)*dgm(1,3)& 938 & +15*gm(1,1)*dgm(2,3))-3*gm(1,2)**2*dgm(3,3)+gm(2,2)*(-1.5d0*gm(3,3)& 939 & *dgm(1,1)-6*gm(1,3)*dgm(1,3)-1.5d0*gm(1,1)*dgm(3,3)) 940 cm(3,4,9)=12*gm(1,1)*gm(2,3)*dgm(2,2)+gm(1,2)*(6*gm(2,3)*dgm(1,2)& 941 & +3*gm(1,3)*dgm(2,2))+3*gm(1,2)**2*dgm(2,3)+gm(2,2)*(12*gm(2,3)& 942 & *dgm(1,1)+3*gm(1,3)*dgm(1,2)+3*gm(1,2)*dgm(1,3)+12*gm(1,1)*dgm(2,3)) 943 cm(3,5,9)=1.5d0*gm(1,2)**2*dgm(1,3)+gm(1,1)*(12*gm(2,3)*dgm(1,2)& 944 & -4.5d0*gm(2,2)*dgm(1,3))+gm(1,3)*(-4.5d0*gm(2,2)*dgm(1,1)+3*gm(1,2)& 945 & *dgm(1,2)-4.5d0*gm(1,1)*dgm(2,2))+12*gm(1,2)*(gm(2,3)*dgm(1,1)& 946 & +gm(1,1)*dgm(2,3)) 947 cm(3,6,9)=4.5d0*gm(1,2)**2*dgm(1,2)+7.5d0*gm(1,1)*gm(2,2)*dgm(1,2)& 948 & +7.5d0*gm(1,2)*(gm(2,2)*dgm(1,1)+gm(1,1)*dgm(2,2)) 949 cm(3,7,9)=gm(2,2)*(3*gm(2,2)*dgm(1,2)+6*gm(1,2)*dgm(2,2)) 950 cm(3,8,9)=1.5d0*gm(2,3)**2*dgm(1,2)-4.5d0*gm(1,2)*gm(3,3)*dgm(2,2)& 951 & +gm(2,3)*(12*gm(2,2)*dgm(1,3)+12*gm(1,3)*dgm(2,2)+3*gm(1,2)*dgm(2,3))& 952 & +gm(2,2)*(-4.5d0*gm(3,3)*dgm(1,2)+12*gm(1,3)*dgm(2,3)-4.5d0*gm(1,2)& 953 & *dgm(3,3)) 954 cm(3,9,9)=6*gm(2,2)**2*dgm(1,3)+3*gm(1,2)*gm(2,3)*dgm(2,2)+gm(2,2)& 955 & *(3*gm(2,3)*dgm(1,2)+12*gm(1,3)*dgm(2,2)+3*gm(1,2)*dgm(2,3)) 956 cm(3,10,9)=7.5d0*gm(2,3)**2*dgm(1,3)+gm(3,3)*(-1.5d0*gm(1,3)*dgm(2,2)& 957 & -3*gm(1,2)*dgm(2,3))+gm(2,3)*(-3*gm(3,3)*dgm(1,2)+15*gm(1,3)& 958 & *dgm(2,3)-3*gm(1,2)*dgm(3,3))-1.5d0*gm(2,2)*(gm(3,3)*dgm(1,3)& 959 & +gm(1,3)*dgm(3,3)) 960 cm(3,1,10)=4.5d0*gm(1,3)**2*dgm(1,1)+9*gm(1,1)*gm(1,3)*dgm(1,3)& 961 & +gm(1,1)*(-3*gm(3,3)*dgm(1,1)-1.5d0*gm(1,1)*dgm(3,3)) 962 cm(3,2,10)=7.5d0*gm(2,3)**2*dgm(1,1)-6*gm(1,2)*gm(3,3)*dgm(1,2)& 963 & -3*gm(1,3)**2*dgm(2,2)-1.5d0*gm(1,1)*gm(3,3)*dgm(2,2)+9*gm(1,2)& 964 & *gm(1,3)*dgm(2,3)+gm(2,3)*(9*gm(1,3)*dgm(1,2)+9*gm(1,2)*dgm(1,3)& 965 & +15*gm(1,1)*dgm(2,3))-3*gm(1,2)**2*dgm(3,3)+gm(2,2)*(-1.5d0*gm(3,3)& 966 & *dgm(1,1)-6*gm(1,3)*dgm(1,3)-1.5d0*gm(1,1)*dgm(3,3)) 967 cm(3,3,10)=6*gm(3,3)**2*dgm(1,1)+3*gm(1,3)**2*dgm(3,3)+gm(3,3)& 968 & *(6*gm(1,3)*dgm(1,3)+12*gm(1,1)*dgm(3,3)) 969 cm(3,4,10)=3*gm(1,3)**2*dgm(2,3)+gm(3,3)*(3*gm(1,2)*dgm(1,3)+12*gm(1,1)& 970 & *dgm(2,3))+gm(2,3)*(12*gm(3,3)*dgm(1,1)+6*gm(1,3)*dgm(1,3)+12*gm(1,1)& 971 & *dgm(3,3))+3*gm(1,3)*(gm(3,3)*dgm(1,2)+gm(1,2)*dgm(3,3)) 972 cm(3,5,10)=4.5d0*gm(1,3)**2*dgm(1,3)+7.5d0*gm(1,1)*gm(3,3)*dgm(1,3)& 973 & +7.5d0*gm(1,3)*(gm(3,3)*dgm(1,1)+gm(1,1)*dgm(3,3)) 974 cm(3,6,10)=1.5d0*gm(1,3)**2*dgm(1,2)+gm(1,1)*(-4.5d0*gm(3,3)*dgm(1,2)& 975 & +12*gm(2,3)*dgm(1,3))+gm(1,3)*(12*gm(2,3)*dgm(1,1)+3*gm(1,2)& 976 & *dgm(1,3)+12*gm(1,1)*dgm(2,3))-4.5d0*gm(1,2)*(gm(3,3)*dgm(1,1)& 977 & +gm(1,1)*dgm(3,3)) 978 cm(3,7,10)=7.5d0*gm(2,3)**2*dgm(1,2)-1.5d0*gm(1,2)*gm(3,3)*dgm(2,2)& 979 & +gm(2,3)*(-3*gm(2,2)*dgm(1,3)-3*gm(1,3)*dgm(2,2)+15*gm(1,2)*dgm(2,3))& 980 & +gm(2,2)*(-1.5d0*gm(3,3)*dgm(1,2)-3*gm(1,3)*dgm(2,3)-1.5d0*gm(1,2)& 981 & *dgm(3,3)) 982 cm(3,8,10)=6*gm(3,3)**2*dgm(1,2)+3*gm(1,3)*gm(2,3)*dgm(3,3)+gm(3,3)& 983 & *(3*gm(2,3)*dgm(1,3)+3*gm(1,3)*dgm(2,3)+12*gm(1,2)*dgm(3,3)) 984 cm(3,9,10)=1.5d0*gm(2,3)**2*dgm(1,3)+gm(3,3)*(-4.5d0*gm(1,3)*dgm(2,2)& 985 & +12*gm(1,2)*dgm(2,3))+gm(2,3)*(12*gm(3,3)*dgm(1,2)+3*gm(1,3)& 986 & *dgm(2,3)+12*gm(1,2)*dgm(3,3))-4.5d0*gm(2,2)*(gm(3,3)*dgm(1,3)& 987 & +gm(1,3)*dgm(3,3)) 988 cm(3,10,10)=gm(3,3)*(3*gm(3,3)*dgm(1,3)+6*gm(1,3)*dgm(3,3)) 989 cm(3,1,12)=7.5d0*gm(1,2)**2*dgm(1,3)+gm(1,1)*(-3*gm(2,3)*dgm(1,2)& 990 & -1.5d0*gm(2,2)*dgm(1,3))+gm(1,3)*(-1.5d0*gm(2,2)*dgm(1,1)+15*gm(1,2)& 991 & *dgm(1,2)-1.5d0*gm(1,1)*dgm(2,2))-3*gm(1,2)*(gm(2,3)*dgm(1,1)& 992 & +gm(1,1)*dgm(2,3)) 993 cm(3,2,12)=6*gm(2,2)**2*dgm(1,3)+3*gm(1,2)*gm(2,3)*dgm(2,2)+gm(2,2)& 994 & *(3*gm(2,3)*dgm(1,2)+12*gm(1,3)*dgm(2,2)+3*gm(1,2)*dgm(2,3)) 995 cm(3,3,12)=1.5d0*gm(2,3)**2*dgm(1,3)+gm(3,3)*(-4.5d0*gm(1,3)*dgm(2,2)& 996 & +12*gm(1,2)*dgm(2,3))+gm(2,3)*(12*gm(3,3)*dgm(1,2)+3*gm(1,3)& 997 & *dgm(2,3)+12*gm(1,2)*dgm(3,3))-4.5d0*gm(2,2)*(gm(3,3)*dgm(1,3)& 998 & +gm(1,3)*dgm(3,3)) 999 cm(3,4,12)=3*gm(2,3)**2*dgm(1,2)+12*gm(1,2)*gm(3,3)*dgm(2,2)+gm(2,3)& 1000 & *(3*gm(2,2)*dgm(1,3)+3*gm(1,3)*dgm(2,2)+6*gm(1,2)*dgm(2,3))+gm(2,2)& 1001 & *(12*gm(3,3)*dgm(1,2)+3*gm(1,3)*dgm(2,3)+12*gm(1,2)*dgm(3,3)) 1002 cm(3,5,12)=-3*gm(2,3)**2*dgm(1,1)+15*gm(1,2)*gm(3,3)*dgm(1,2)& 1003 & -3*gm(1,3)**2*dgm(2,2)-1.5d0*gm(1,1)*gm(3,3)*dgm(2,2)+9*gm(1,2)& 1004 & *gm(1,3)*dgm(2,3)+gm(2,3)*(9*gm(1,3)*dgm(1,2)+9*gm(1,2)*dgm(1,3)& 1005 & -6*gm(1,1)*dgm(2,3))+7.5d0*gm(1,2)**2*dgm(3,3)+gm(2,2)*(-1.5d0*gm(3,3)& 1006 & *dgm(1,1)-6*gm(1,3)*dgm(1,3)-1.5d0*gm(1,1)*dgm(3,3)) 1007 cm(3,6,12)=-4.5d0*gm(1,1)*gm(2,3)*dgm(2,2)+gm(1,2)*(3*gm(2,3)& 1008 & *dgm(1,2)+12*gm(1,3)*dgm(2,2))+1.5d0*gm(1,2)**2*dgm(2,3)+gm(2,2)& 1009 & *(-4.5d0*gm(2,3)*dgm(1,1)+12*gm(1,3)*dgm(1,2)+12*gm(1,2)*dgm(1,3)& 1010 & -4.5d0*gm(1,1)*dgm(2,3)) 1011 cm(3,7,12)=gm(2,2)*(6*gm(2,3)*dgm(2,2)+3*gm(2,2)*dgm(2,3)) 1012 cm(3,8,12)=4.5d0*gm(2,3)**2*dgm(2,3)+7.5d0*gm(2,2)*gm(3,3)*dgm(2,3)& 1013 & +7.5d0*gm(2,3)*(gm(3,3)*dgm(2,2)+gm(2,2)*dgm(3,3)) 1014 cm(3,9,12)=3*gm(2,3)**2*dgm(2,2)+6*gm(2,2)*gm(2,3)*dgm(2,3)+gm(2,2)& 1015 & *(12*gm(3,3)*dgm(2,2)+6*gm(2,2)*dgm(3,3)) 1016 cm(3,10,12)=-1.5d0*gm(3,3)**2*dgm(2,2)+4.5d0*gm(2,3)**2*dgm(3,3)& 1017 & +gm(3,3)*(9*gm(2,3)*dgm(2,3)-3*gm(2,2)*dgm(3,3)) 1018 cm(3,1,13)=7.5d0*gm(1,2)**2*dgm(1,2)-1.5d0*gm(1,1)*gm(2,2)*dgm(1,2)& 1019 & -1.5d0*gm(1,2)*(gm(2,2)*dgm(1,1)+gm(1,1)*dgm(2,2)) 1020 cm(3,2,13)=gm(2,2)*(3*gm(2,2)*dgm(1,2)+6*gm(1,2)*dgm(2,2)) 1021 cm(3,3,13)=7.5d0*gm(2,3)**2*dgm(1,2)-1.5d0*gm(1,2)*gm(3,3)*dgm(2,2)& 1022 & +gm(2,3)*(-3*gm(2,2)*dgm(1,3)-3*gm(1,3)*dgm(2,2)+15*gm(1,2)*dgm(2,3))& 1023 & +gm(2,2)*(-1.5d0*gm(3,3)*dgm(1,2)-3*gm(1,3)*dgm(2,3)-1.5d0*gm(1,2)& 1024 & *dgm(3,3)) 1025 cm(3,4,13)=-3*gm(2,2)**2*dgm(1,3)+9*gm(1,2)*gm(2,3)*dgm(2,2)+gm(2,2)& 1026 & *(9*gm(2,3)*dgm(1,2)-6*gm(1,3)*dgm(2,2)+9*gm(1,2)*dgm(2,3)) 1027 cm(3,5,13)=-1.5d0*gm(1,1)*gm(2,3)*dgm(2,2)+gm(1,2)*(15*gm(2,3)& 1028 & *dgm(1,2)-3*gm(1,3)*dgm(2,2))+7.5d0*gm(1,2)**2*dgm(2,3)+gm(2,2)& 1029 & *(-1.5d0*gm(2,3)*dgm(1,1)-3*gm(1,3)*dgm(1,2)-3*gm(1,2)*dgm(1,3)& 1030 & -1.5d0*gm(1,1)*dgm(2,3)) 1031 cm(3,6,13)=-1.5d0*gm(2,2)**2*dgm(1,1)+4.5d0*gm(1,2)**2*dgm(2,2)& 1032 & +gm(2,2)*(9*gm(1,2)*dgm(1,2)-3*gm(1,1)*dgm(2,2)) 1033 cm(3,7,13)=3*gm(2,2)**2*dgm(2,2) 1034 cm(3,8,13)=4.5d0*gm(2,3)**2*dgm(2,2)+9*gm(2,2)*gm(2,3)*dgm(2,3)& 1035 & +gm(2,2)*(-3*gm(3,3)*dgm(2,2)-1.5d0*gm(2,2)*dgm(3,3)) 1036 cm(3,9,13)=gm(2,2)*(6*gm(2,3)*dgm(2,2)+3*gm(2,2)*dgm(2,3)) 1037 cm(3,10,13)=7.5d0*gm(2,3)**2*dgm(2,3)-1.5d0*gm(2,2)*gm(3,3)*dgm(2,3)& 1038 & -1.5d0*gm(2,3)*(gm(3,3)*dgm(2,2)+gm(2,2)*dgm(3,3)) 1039 cm(3,1,14)=7.5d0*gm(1,3)**2*dgm(1,2)+gm(1,1)*(-1.5d0*gm(3,3)*dgm(1,2)& 1040 & -3*gm(2,3)*dgm(1,3))+gm(1,3)*(-3*gm(2,3)*dgm(1,1)+15*gm(1,2)& 1041 & *dgm(1,3)-3*gm(1,1)*dgm(2,3))-1.5d0*gm(1,2)*(gm(3,3)*dgm(1,1)& 1042 & +gm(1,1)*dgm(3,3)) 1043 cm(3,2,14)=1.5d0*gm(2,3)**2*dgm(1,2)-4.5d0*gm(1,2)*gm(3,3)*dgm(2,2)& 1044 & +gm(2,3)*(12*gm(2,2)*dgm(1,3)+12*gm(1,3)*dgm(2,2)+3*gm(1,2)*dgm(2,3))& 1045 & +gm(2,2)*(-4.5d0*gm(3,3)*dgm(1,2)+12*gm(1,3)*dgm(2,3)-4.5d0*gm(1,2)& 1046 & *dgm(3,3)) 1047 cm(3,3,14)=6*gm(3,3)**2*dgm(1,2)+3*gm(1,3)*gm(2,3)*dgm(3,3)+gm(3,3)& 1048 & *(3*gm(2,3)*dgm(1,3)+3*gm(1,3)*dgm(2,3)+12*gm(1,2)*dgm(3,3)) 1049 cm(3,4,14)=3*gm(2,3)**2*dgm(1,3)+gm(3,3)*(12*gm(1,3)*dgm(2,2)& 1050 & +3*gm(1,2)*dgm(2,3))+gm(2,3)*(3*gm(3,3)*dgm(1,2)+6*gm(1,3)*dgm(2,3)& 1051 & +3*gm(1,2)*dgm(3,3))+12*gm(2,2)*(gm(3,3)*dgm(1,3)+gm(1,3)*dgm(3,3)) 1052 cm(3,5,14)=1.5d0*gm(1,3)**2*dgm(2,3)+gm(3,3)*(12*gm(1,2)*dgm(1,3)& 1053 & -4.5d0*gm(1,1)*dgm(2,3))+gm(2,3)*(-4.5d0*gm(3,3)*dgm(1,1)+3*gm(1,3)& 1054 & *dgm(1,3)-4.5d0*gm(1,1)*dgm(3,3))+12*gm(1,3)*(gm(3,3)*dgm(1,2)& 1055 & +gm(1,2)*dgm(3,3)) 1056 cm(3,6,14)=-3*gm(2,3)**2*dgm(1,1)-6*gm(1,2)*gm(3,3)*dgm(1,2)+7.5d0*gm(1,3)& 1057 & **2*dgm(2,2)-1.5d0*gm(1,1)*gm(3,3)*dgm(2,2)+9*gm(1,2)*gm(1,3)& 1058 & *dgm(2,3)+gm(2,3)*(9*gm(1,3)*dgm(1,2)+9*gm(1,2)*dgm(1,3)-6*gm(1,1)& 1059 & *dgm(2,3))-3*gm(1,2)**2*dgm(3,3)+gm(2,2)*(-1.5d0*gm(3,3)*dgm(1,1)& 1060 & +15*gm(1,3)*dgm(1,3)-1.5d0*gm(1,1)*dgm(3,3)) 1061 cm(3,7,14)=4.5d0*gm(2,3)**2*dgm(2,2)+9*gm(2,2)*gm(2,3)*dgm(2,3)& 1062 & +gm(2,2)*(-3*gm(3,3)*dgm(2,2)-1.5d0*gm(2,2)*dgm(3,3)) 1063 cm(3,8,14)=6*gm(3,3)**2*dgm(2,2)+3*gm(2,3)**2*dgm(3,3)+gm(3,3)& 1064 & *(6*gm(2,3)*dgm(2,3)+12*gm(2,2)*dgm(3,3)) 1065 cm(3,9,14)=4.5d0*gm(2,3)**2*dgm(2,3)+7.5d0*gm(2,2)*gm(3,3)*dgm(2,3)& 1066 & +7.5d0*gm(2,3)*(gm(3,3)*dgm(2,2)+gm(2,2)*dgm(3,3)) 1067 cm(3,10,14)=gm(3,3)*(3*gm(3,3)*dgm(2,3)+6*gm(2,3)*dgm(3,3)) 1068 cm(3,1,15)=7.5d0*gm(1,3)**2*dgm(1,3)-1.5d0*gm(1,1)*gm(3,3)*dgm(1,3)& 1069 & -1.5d0*gm(1,3)*(gm(3,3)*dgm(1,1)+gm(1,1)*dgm(3,3)) 1070 cm(3,2,15)=7.5d0*gm(2,3)**2*dgm(1,3)+gm(3,3)*(-1.5d0*gm(1,3)*dgm(2,2)& 1071 & -3*gm(1,2)*dgm(2,3))+gm(2,3)*(-3*gm(3,3)*dgm(1,2)+15*gm(1,3)& 1072 & *dgm(2,3)-3*gm(1,2)*dgm(3,3))-1.5d0*gm(2,2)*(gm(3,3)*dgm(1,3)& 1073 & +gm(1,3)*dgm(3,3)) 1074 cm(3,3,15)=gm(3,3)*(3*gm(3,3)*dgm(1,3)+6*gm(1,3)*dgm(3,3)) 1075 cm(3,4,15)=-3*gm(3,3)**2*dgm(1,2)+9*gm(1,3)*gm(2,3)*dgm(3,3)+gm(3,3)& 1076 & *(9*gm(2,3)*dgm(1,3)+9*gm(1,3)*dgm(2,3)-6*gm(1,2)*dgm(3,3)) 1077 cm(3,5,15)=-1.5d0*gm(3,3)**2*dgm(1,1)+4.5d0*gm(1,3)**2*dgm(3,3)& 1078 & +gm(3,3)*(9*gm(1,3)*dgm(1,3)-3*gm(1,1)*dgm(3,3)) 1079 cm(3,6,15)=7.5d0*gm(1,3)**2*dgm(2,3)+gm(3,3)*(-3*gm(1,2)*dgm(1,3)& 1080 & -1.5d0*gm(1,1)*dgm(2,3))+gm(2,3)*(-1.5d0*gm(3,3)*dgm(1,1)+15*gm(1,3)& 1081 & *dgm(1,3)-1.5d0*gm(1,1)*dgm(3,3))-3*gm(1,3)*(gm(3,3)*dgm(1,2)& 1082 & +gm(1,2)*dgm(3,3)) 1083 cm(3,7,15)=7.5d0*gm(2,3)**2*dgm(2,3)-1.5d0*gm(2,2)*gm(3,3)*dgm(2,3)& 1084 & -1.5d0*gm(2,3)*(gm(3,3)*dgm(2,2)+gm(2,2)*dgm(3,3)) 1085 cm(3,8,15)=gm(3,3)*(3*gm(3,3)*dgm(2,3)+6*gm(2,3)*dgm(3,3)) 1086 cm(3,9,15)=-1.5d0*gm(3,3)**2*dgm(2,2)+4.5d0*gm(2,3)**2*dgm(3,3)& 1087 & +gm(3,3)*(9*gm(2,3)*dgm(2,3)-3*gm(2,2)*dgm(3,3)) 1088 cm(3,10,15)=3*gm(3,3)**2*dgm(3,3) 1089 end if 1090 ! 1091 !Section for volume derivative term for diagonal strains. This is 1092 !equivalent to calculating the reduced gradients but it is easier 1093 !to embed it here as a separate calculation than cope with the logic 1094 !and ffkg index rearrangements necessary to do a separate force 1095 !calculation. 1096 ! 1097 if(istr<=3) then 1098 if(rank==0)then 1099 cm(1,1,1)=cm(1,1,1)-(1) 1100 cm(2,1,2)=cm(2,1,2)-(1) 1101 cm(3,1,3)=cm(3,1,3)-(1) 1102 elseif(rank==1)then 1103 cm(1,1,1)=cm(1,1,1)-(gm(1,1)) 1104 cm(1,2,1)=cm(1,2,1)-(gm(1,2)) 1105 cm(1,3,1)=cm(1,3,1)-(gm(1,3)) 1106 cm(1,1,5)=cm(1,1,5)-(gm(1,3)) 1107 cm(1,2,5)=cm(1,2,5)-(gm(2,3)) 1108 cm(1,3,5)=cm(1,3,5)-(gm(3,3)) 1109 cm(1,1,6)=cm(1,1,6)-(gm(1,2)) 1110 cm(1,2,6)=cm(1,2,6)-(gm(2,2)) 1111 cm(1,3,6)=cm(1,3,6)-(gm(2,3)) 1112 cm(2,1,2)=cm(2,1,2)-(gm(1,2)) 1113 cm(2,2,2)=cm(2,2,2)-(gm(2,2)) 1114 cm(2,3,2)=cm(2,3,2)-(gm(2,3)) 1115 cm(2,1,4)=cm(2,1,4)-(gm(1,3)) 1116 cm(2,2,4)=cm(2,2,4)-(gm(2,3)) 1117 cm(2,3,4)=cm(2,3,4)-(gm(3,3)) 1118 cm(2,1,6)=cm(2,1,6)-(gm(1,1)) 1119 cm(2,2,6)=cm(2,2,6)-(gm(1,2)) 1120 cm(2,3,6)=cm(2,3,6)-(gm(1,3)) 1121 cm(3,1,3)=cm(3,1,3)-(gm(1,3)) 1122 cm(3,2,3)=cm(3,2,3)-(gm(2,3)) 1123 cm(3,3,3)=cm(3,3,3)-(gm(3,3)) 1124 cm(3,1,4)=cm(3,1,4)-(gm(1,2)) 1125 cm(3,2,4)=cm(3,2,4)-(gm(2,2)) 1126 cm(3,3,4)=cm(3,3,4)-(gm(2,3)) 1127 cm(3,1,5)=cm(3,1,5)-(gm(1,1)) 1128 cm(3,2,5)=cm(3,2,5)-(gm(1,2)) 1129 cm(3,3,5)=cm(3,3,5)-(gm(1,3)) 1130 elseif(rank==2)then 1131 cm(1,1,1)=cm(1,1,1)-(gm(1,1)**2) 1132 cm(1,2,1)=cm(1,2,1)-(1.5d0*gm(1,2)**2-0.5d0*gm(1,1)*gm(2,2)) 1133 cm(1,3,1)=cm(1,3,1)-(1.5d0*gm(1,3)**2-0.5d0*gm(1,1)*gm(3,3)) 1134 cm(1,4,1)=cm(1,4,1)-(3*gm(1,2)*gm(1,3)-gm(1,1)*gm(2,3)) 1135 cm(1,5,1)=cm(1,5,1)-(2*gm(1,1)*gm(1,3)) 1136 cm(1,6,1)=cm(1,6,1)-(2*gm(1,1)*gm(1,2)) 1137 cm(1,1,2)=cm(1,1,2)-(1.5d0*gm(1,2)**2-0.5d0*gm(1,1)*gm(2,2)) 1138 cm(1,2,2)=cm(1,2,2)-(gm(2,2)**2) 1139 cm(1,3,2)=cm(1,3,2)-(1.5d0*gm(2,3)**2-0.5d0*gm(2,2)*gm(3,3)) 1140 cm(1,4,2)=cm(1,4,2)-(2*gm(2,2)*gm(2,3)) 1141 cm(1,5,2)=cm(1,5,2)-(-gm(1,3)*gm(2,2)+3*gm(1,2)*gm(2,3)) 1142 cm(1,6,2)=cm(1,6,2)-(2*gm(1,2)*gm(2,2)) 1143 cm(1,1,3)=cm(1,1,3)-(1.5d0*gm(1,3)**2-0.5d0*gm(1,1)*gm(3,3)) 1144 cm(1,2,3)=cm(1,2,3)-(1.5d0*gm(2,3)**2-0.5d0*gm(2,2)*gm(3,3)) 1145 cm(1,3,3)=cm(1,3,3)-(gm(3,3)**2) 1146 cm(1,4,3)=cm(1,4,3)-(2*gm(2,3)*gm(3,3)) 1147 cm(1,5,3)=cm(1,5,3)-(2*gm(1,3)*gm(3,3)) 1148 cm(1,6,3)=cm(1,6,3)-(3*gm(1,3)*gm(2,3)-gm(1,2)*gm(3,3)) 1149 cm(1,1,4)=cm(1,1,4)-(3*gm(1,2)*gm(1,3)-gm(1,1)*gm(2,3)) 1150 cm(1,2,4)=cm(1,2,4)-(2*gm(2,2)*gm(2,3)) 1151 cm(1,3,4)=cm(1,3,4)-(2*gm(2,3)*gm(3,3)) 1152 cm(1,4,4)=cm(1,4,4)-(gm(2,3)**2+3*gm(2,2)*gm(3,3)) 1153 cm(1,5,4)=cm(1,5,4)-(gm(1,3)*gm(2,3)+3*gm(1,2)*gm(3,3)) 1154 cm(1,6,4)=cm(1,6,4)-(3*gm(1,3)*gm(2,2)+gm(1,2)*gm(2,3)) 1155 cm(1,1,5)=cm(1,1,5)-(2*gm(1,1)*gm(1,3)) 1156 cm(1,2,5)=cm(1,2,5)-(-gm(1,3)*gm(2,2)+3*gm(1,2)*gm(2,3)) 1157 cm(1,3,5)=cm(1,3,5)-(2*gm(1,3)*gm(3,3)) 1158 cm(1,4,5)=cm(1,4,5)-(gm(1,3)*gm(2,3)+3*gm(1,2)*gm(3,3)) 1159 cm(1,5,5)=cm(1,5,5)-(gm(1,3)**2+3*gm(1,1)*gm(3,3)) 1160 cm(1,6,5)=cm(1,6,5)-(gm(1,2)*gm(1,3)+3*gm(1,1)*gm(2,3)) 1161 cm(1,1,6)=cm(1,1,6)-(2*gm(1,1)*gm(1,2)) 1162 cm(1,2,6)=cm(1,2,6)-(2*gm(1,2)*gm(2,2)) 1163 cm(1,3,6)=cm(1,3,6)-(3*gm(1,3)*gm(2,3)-gm(1,2)*gm(3,3)) 1164 cm(1,4,6)=cm(1,4,6)-(3*gm(1,3)*gm(2,2)+gm(1,2)*gm(2,3)) 1165 cm(1,5,6)=cm(1,5,6)-(gm(1,2)*gm(1,3)+3*gm(1,1)*gm(2,3)) 1166 cm(1,6,6)=cm(1,6,6)-(gm(1,2)**2+3*gm(1,1)*gm(2,2)) 1167 cm(2,1,2)=cm(2,1,2)-(2*gm(1,1)*gm(1,2)) 1168 cm(2,2,2)=cm(2,2,2)-(2*gm(1,2)*gm(2,2)) 1169 cm(2,3,2)=cm(2,3,2)-(3*gm(1,3)*gm(2,3)-gm(1,2)*gm(3,3)) 1170 cm(2,4,2)=cm(2,4,2)-(3*gm(1,3)*gm(2,2)+gm(1,2)*gm(2,3)) 1171 cm(2,5,2)=cm(2,5,2)-(gm(1,2)*gm(1,3)+3*gm(1,1)*gm(2,3)) 1172 cm(2,6,2)=cm(2,6,2)-(gm(1,2)**2+3*gm(1,1)*gm(2,2)) 1173 cm(2,1,4)=cm(2,1,4)-(2*gm(1,1)*gm(1,3)) 1174 cm(2,2,4)=cm(2,2,4)-(-gm(1,3)*gm(2,2)+3*gm(1,2)*gm(2,3)) 1175 cm(2,3,4)=cm(2,3,4)-(2*gm(1,3)*gm(3,3)) 1176 cm(2,4,4)=cm(2,4,4)-(gm(1,3)*gm(2,3)+3*gm(1,2)*gm(3,3)) 1177 cm(2,5,4)=cm(2,5,4)-(gm(1,3)**2+3*gm(1,1)*gm(3,3)) 1178 cm(2,6,4)=cm(2,6,4)-(gm(1,2)*gm(1,3)+3*gm(1,1)*gm(2,3)) 1179 cm(2,1,6)=cm(2,1,6)-(gm(1,1)**2) 1180 cm(2,2,6)=cm(2,2,6)-(1.5d0*gm(1,2)**2-0.5d0*gm(1,1)*gm(2,2)) 1181 cm(2,3,6)=cm(2,3,6)-(1.5d0*gm(1,3)**2-0.5d0*gm(1,1)*gm(3,3)) 1182 cm(2,4,6)=cm(2,4,6)-(3*gm(1,2)*gm(1,3)-gm(1,1)*gm(2,3)) 1183 cm(2,5,6)=cm(2,5,6)-(2*gm(1,1)*gm(1,3)) 1184 cm(2,6,6)=cm(2,6,6)-(2*gm(1,1)*gm(1,2)) 1185 cm(2,1,7)=cm(2,1,7)-(1.5d0*gm(1,2)**2-0.5d0*gm(1,1)*gm(2,2)) 1186 cm(2,2,7)=cm(2,2,7)-(gm(2,2)**2) 1187 cm(2,3,7)=cm(2,3,7)-(1.5d0*gm(2,3)**2-0.5d0*gm(2,2)*gm(3,3)) 1188 cm(2,4,7)=cm(2,4,7)-(2*gm(2,2)*gm(2,3)) 1189 cm(2,5,7)=cm(2,5,7)-(-gm(1,3)*gm(2,2)+3*gm(1,2)*gm(2,3)) 1190 cm(2,6,7)=cm(2,6,7)-(2*gm(1,2)*gm(2,2)) 1191 cm(2,1,8)=cm(2,1,8)-(1.5d0*gm(1,3)**2-0.5d0*gm(1,1)*gm(3,3)) 1192 cm(2,2,8)=cm(2,2,8)-(1.5d0*gm(2,3)**2-0.5d0*gm(2,2)*gm(3,3)) 1193 cm(2,3,8)=cm(2,3,8)-(gm(3,3)**2) 1194 cm(2,4,8)=cm(2,4,8)-(2*gm(2,3)*gm(3,3)) 1195 cm(2,5,8)=cm(2,5,8)-(2*gm(1,3)*gm(3,3)) 1196 cm(2,6,8)=cm(2,6,8)-(3*gm(1,3)*gm(2,3)-gm(1,2)*gm(3,3)) 1197 cm(2,1,9)=cm(2,1,9)-(3*gm(1,2)*gm(1,3)-gm(1,1)*gm(2,3)) 1198 cm(2,2,9)=cm(2,2,9)-(2*gm(2,2)*gm(2,3)) 1199 cm(2,3,9)=cm(2,3,9)-(2*gm(2,3)*gm(3,3)) 1200 cm(2,4,9)=cm(2,4,9)-(gm(2,3)**2+3*gm(2,2)*gm(3,3)) 1201 cm(2,5,9)=cm(2,5,9)-(gm(1,3)*gm(2,3)+3*gm(1,2)*gm(3,3)) 1202 cm(2,6,9)=cm(2,6,9)-(3*gm(1,3)*gm(2,2)+gm(1,2)*gm(2,3)) 1203 cm(3,1,3)=cm(3,1,3)-(2*gm(1,1)*gm(1,3)) 1204 cm(3,2,3)=cm(3,2,3)-(-gm(1,3)*gm(2,2)+3*gm(1,2)*gm(2,3)) 1205 cm(3,3,3)=cm(3,3,3)-(2*gm(1,3)*gm(3,3)) 1206 cm(3,4,3)=cm(3,4,3)-(gm(1,3)*gm(2,3)+3*gm(1,2)*gm(3,3)) 1207 cm(3,5,3)=cm(3,5,3)-(gm(1,3)**2+3*gm(1,1)*gm(3,3)) 1208 cm(3,6,3)=cm(3,6,3)-(gm(1,2)*gm(1,3)+3*gm(1,1)*gm(2,3)) 1209 cm(3,1,4)=cm(3,1,4)-(2*gm(1,1)*gm(1,2)) 1210 cm(3,2,4)=cm(3,2,4)-(2*gm(1,2)*gm(2,2)) 1211 cm(3,3,4)=cm(3,3,4)-(3*gm(1,3)*gm(2,3)-gm(1,2)*gm(3,3)) 1212 cm(3,4,4)=cm(3,4,4)-(3*gm(1,3)*gm(2,2)+gm(1,2)*gm(2,3)) 1213 cm(3,5,4)=cm(3,5,4)-(gm(1,2)*gm(1,3)+3*gm(1,1)*gm(2,3)) 1214 cm(3,6,4)=cm(3,6,4)-(gm(1,2)**2+3*gm(1,1)*gm(2,2)) 1215 cm(3,1,5)=cm(3,1,5)-(gm(1,1)**2) 1216 cm(3,2,5)=cm(3,2,5)-(1.5d0*gm(1,2)**2-0.5d0*gm(1,1)*gm(2,2)) 1217 cm(3,3,5)=cm(3,3,5)-(1.5d0*gm(1,3)**2-0.5d0*gm(1,1)*gm(3,3)) 1218 cm(3,4,5)=cm(3,4,5)-(3*gm(1,2)*gm(1,3)-gm(1,1)*gm(2,3)) 1219 cm(3,5,5)=cm(3,5,5)-(2*gm(1,1)*gm(1,3)) 1220 cm(3,6,5)=cm(3,6,5)-(2*gm(1,1)*gm(1,2)) 1221 cm(3,1,8)=cm(3,1,8)-(3*gm(1,2)*gm(1,3)-gm(1,1)*gm(2,3)) 1222 cm(3,2,8)=cm(3,2,8)-(2*gm(2,2)*gm(2,3)) 1223 cm(3,3,8)=cm(3,3,8)-(2*gm(2,3)*gm(3,3)) 1224 cm(3,4,8)=cm(3,4,8)-(gm(2,3)**2+3*gm(2,2)*gm(3,3)) 1225 cm(3,5,8)=cm(3,5,8)-(gm(1,3)*gm(2,3)+3*gm(1,2)*gm(3,3)) 1226 cm(3,6,8)=cm(3,6,8)-(3*gm(1,3)*gm(2,2)+gm(1,2)*gm(2,3)) 1227 cm(3,1,9)=cm(3,1,9)-(1.5d0*gm(1,2)**2-0.5d0*gm(1,1)*gm(2,2)) 1228 cm(3,2,9)=cm(3,2,9)-(gm(2,2)**2) 1229 cm(3,3,9)=cm(3,3,9)-(1.5d0*gm(2,3)**2-0.5d0*gm(2,2)*gm(3,3)) 1230 cm(3,4,9)=cm(3,4,9)-(2*gm(2,2)*gm(2,3)) 1231 cm(3,5,9)=cm(3,5,9)-(-gm(1,3)*gm(2,2)+3*gm(1,2)*gm(2,3)) 1232 cm(3,6,9)=cm(3,6,9)-(2*gm(1,2)*gm(2,2)) 1233 cm(3,1,10)=cm(3,1,10)-(1.5d0*gm(1,3)**2-0.5d0*gm(1,1)*gm(3,3)) 1234 cm(3,2,10)=cm(3,2,10)-(1.5d0*gm(2,3)**2-0.5d0*gm(2,2)*gm(3,3)) 1235 cm(3,3,10)=cm(3,3,10)-(gm(3,3)**2) 1236 cm(3,4,10)=cm(3,4,10)-(2*gm(2,3)*gm(3,3)) 1237 cm(3,5,10)=cm(3,5,10)-(2*gm(1,3)*gm(3,3)) 1238 cm(3,6,10)=cm(3,6,10)-(3*gm(1,3)*gm(2,3)-gm(1,2)*gm(3,3)) 1239 elseif(rank==3)then 1240 cm(1,1,1)=cm(1,1,1)-(gm(1,1)**3) 1241 cm(1,2,1)=cm(1,2,1)-(gm(1,1)*(4.5d0*gm(1,2)**2-1.5d0*gm(1,1)*gm(2,2))) 1242 cm(1,3,1)=cm(1,3,1)-(gm(1,1)*(4.5d0*gm(1,3)**2-1.5d0*gm(1,1)*gm(3,3))) 1243 cm(1,4,1)=cm(1,4,1)-(gm(1,1)*(9*gm(1,2)*gm(1,3)-3*gm(1,1)*gm(2,3))) 1244 cm(1,5,1)=cm(1,5,1)-(3*gm(1,1)**2*gm(1,3)) 1245 cm(1,6,1)=cm(1,6,1)-(3*gm(1,1)**2*gm(1,2)) 1246 cm(1,7,1)=cm(1,7,1)-(2.5d0*gm(1,2)**3-1.5d0*gm(1,1)*gm(1,2)*gm(2,2)) 1247 cm(1,8,1)=cm(1,8,1)-(-3*gm(1,1)*gm(1,3)*gm(2,3)+gm(1,2)*(7.5d0*gm(1,3)& 1248 & **2-1.5d0*gm(1,1)*gm(3,3))) 1249 cm(1,9,1)=cm(1,9,1)-(7.5d0*gm(1,2)**2*gm(1,3)-1.5d0*gm(1,1)*gm(1,3)& 1250 & *gm(2,2)-3*gm(1,1)*gm(1,2)*gm(2,3)) 1251 cm(1,10,1)=cm(1,10,1)-(2.5d0*gm(1,3)**3-1.5d0*gm(1,1)*gm(1,3)& 1252 & *gm(3,3)) 1253 cm(1,1,2)=cm(1,1,2)-(gm(1,1)*(4.5d0*gm(1,2)**2-1.5d0*gm(1,1)*gm(2,2))) 1254 cm(1,2,2)=cm(1,2,2)-(gm(2,2)*(3*gm(1,2)**2+6*gm(1,1)*gm(2,2))) 1255 cm(1,3,2)=cm(1,3,2)-(-3*gm(1,3)**2*gm(2,2)+9*gm(1,2)*gm(1,3)*gm(2,3)& 1256 & -3*gm(1,2)**2*gm(3,3)+gm(1,1)*(7.5d0*gm(2,3)**2-1.5d0*gm(2,2)& 1257 & *gm(3,3))) 1258 cm(1,4,2)=cm(1,4,2)-(3*gm(1,2)*gm(1,3)*gm(2,2)+3*gm(1,2)**2*gm(2,3)& 1259 & +12*gm(1,1)*gm(2,2)*gm(2,3)) 1260 cm(1,5,2)=cm(1,5,2)-(1.5d0*gm(1,2)**2*gm(1,3)-4.5d0*gm(1,1)*gm(1,3)& 1261 & *gm(2,2)+12*gm(1,1)*gm(1,2)*gm(2,3)) 1262 cm(1,6,2)=cm(1,6,2)-(1.5d0*gm(1,2)**3+7.5d0*gm(1,1)*gm(1,2)*gm(2,2)) 1263 cm(1,7,2)=cm(1,7,2)-(3*gm(1,2)*gm(2,2)**2) 1264 cm(1,8,2)=cm(1,8,2)-(12*gm(1,3)*gm(2,2)*gm(2,3)+gm(1,2)*(1.5d0*gm(2,3)& 1265 & **2-4.5d0*gm(2,2)*gm(3,3))) 1266 cm(1,9,2)=cm(1,9,2)-(gm(2,2)*(6*gm(1,3)*gm(2,2)+3*gm(1,2)*gm(2,3))) 1267 cm(1,10,2)=cm(1,10,2)-(-3*gm(1,2)*gm(2,3)*gm(3,3)+gm(1,3)*(7.5d0*gm(2,3)& 1268 & **2-1.5d0*gm(2,2)*gm(3,3))) 1269 cm(1,1,3)=cm(1,1,3)-(gm(1,1)*(4.5d0*gm(1,3)**2-1.5d0*gm(1,1)*gm(3,3))) 1270 cm(1,2,3)=cm(1,2,3)-(-3*gm(1,3)**2*gm(2,2)+9*gm(1,2)*gm(1,3)*gm(2,3)& 1271 & -3*gm(1,2)**2*gm(3,3)+gm(1,1)*(7.5d0*gm(2,3)**2-1.5d0*gm(2,2)& 1272 & *gm(3,3))) 1273 cm(1,3,3)=cm(1,3,3)-(gm(3,3)*(3*gm(1,3)**2+6*gm(1,1)*gm(3,3))) 1274 cm(1,4,3)=cm(1,4,3)-(3*gm(1,3)**2*gm(2,3)+3*gm(1,2)*gm(1,3)*gm(3,3)& 1275 & +12*gm(1,1)*gm(2,3)*gm(3,3)) 1276 cm(1,5,3)=cm(1,5,3)-(1.5d0*gm(1,3)**3+7.5d0*gm(1,1)*gm(1,3)*gm(3,3)) 1277 cm(1,6,3)=cm(1,6,3)-(12*gm(1,1)*gm(1,3)*gm(2,3)+gm(1,2)*(1.5d0*gm(1,3)& 1278 & **2-4.5d0*gm(1,1)*gm(3,3))) 1279 cm(1,7,3)=cm(1,7,3)-(-3*gm(1,3)*gm(2,2)*gm(2,3)+gm(1,2)*(7.5d0*gm(2,3)& 1280 & **2-1.5d0*gm(2,2)*gm(3,3))) 1281 cm(1,8,3)=cm(1,8,3)-(gm(3,3)*(3*gm(1,3)*gm(2,3)+6*gm(1,2)*gm(3,3))) 1282 cm(1,9,3)=cm(1,9,3)-(12*gm(1,2)*gm(2,3)*gm(3,3)+gm(1,3)*(1.5d0*gm(2,3)& 1283 & **2-4.5d0*gm(2,2)*gm(3,3))) 1284 cm(1,10,3)=cm(1,10,3)-(3*gm(1,3)*gm(3,3)**2) 1285 cm(1,1,4)=cm(1,1,4)-(gm(1,1)*(9*gm(1,2)*gm(1,3)-3*gm(1,1)*gm(2,3))) 1286 cm(1,2,4)=cm(1,2,4)-(3*gm(1,2)*gm(1,3)*gm(2,2)+3*gm(1,2)**2*gm(2,3)& 1287 & +12*gm(1,1)*gm(2,2)*gm(2,3)) 1288 cm(1,3,4)=cm(1,3,4)-(3*gm(1,3)**2*gm(2,3)+3*gm(1,2)*gm(1,3)*gm(3,3)& 1289 & +12*gm(1,1)*gm(2,3)*gm(3,3)) 1290 cm(1,4,4)=cm(1,4,4)-(9*gm(1,3)**2*gm(2,2)-6*gm(1,2)*gm(1,3)*gm(2,3)& 1291 & +9*gm(1,2)**2*gm(3,3)+gm(1,1)*(9*gm(2,3)**2+15*gm(2,2)*gm(3,3))) 1292 cm(1,5,4)=cm(1,5,4)-(3*gm(1,1)*gm(1,3)*gm(2,3)+gm(1,2)*(3*gm(1,3)& 1293 & **2+12*gm(1,1)*gm(3,3))) 1294 cm(1,6,4)=cm(1,6,4)-(3*gm(1,2)**2*gm(1,3)+12*gm(1,1)*gm(1,3)*gm(2,2)& 1295 & +3*gm(1,1)*gm(1,2)*gm(2,3)) 1296 cm(1,7,4)=cm(1,7,4)-(gm(2,2)*(-3*gm(1,3)*gm(2,2)+9*gm(1,2)*gm(2,3))) 1297 cm(1,8,4)=cm(1,8,4)-(3*gm(1,2)*gm(2,3)*gm(3,3)+gm(1,3)*(3*gm(2,3)& 1298 & **2+12*gm(2,2)*gm(3,3))) 1299 cm(1,9,4)=cm(1,9,4)-(3*gm(1,3)*gm(2,2)*gm(2,3)+gm(1,2)*(3*gm(2,3)& 1300 & **2+12*gm(2,2)*gm(3,3))) 1301 cm(1,10,4)=cm(1,10,4)-(gm(3,3)*(9*gm(1,3)*gm(2,3)-3*gm(1,2)*gm(3,3))) 1302 cm(1,1,5)=cm(1,1,5)-(3*gm(1,1)**2*gm(1,3)) 1303 cm(1,2,5)=cm(1,2,5)-(1.5d0*gm(1,2)**2*gm(1,3)-4.5d0*gm(1,1)*gm(1,3)& 1304 & *gm(2,2)+12*gm(1,1)*gm(1,2)*gm(2,3)) 1305 cm(1,3,5)=cm(1,3,5)-(1.5d0*gm(1,3)**3+7.5d0*gm(1,1)*gm(1,3)*gm(3,3)) 1306 cm(1,4,5)=cm(1,4,5)-(3*gm(1,1)*gm(1,3)*gm(2,3)+gm(1,2)*(3*gm(1,3)& 1307 & **2+12*gm(1,1)*gm(3,3))) 1308 cm(1,5,5)=cm(1,5,5)-(gm(1,1)*(3*gm(1,3)**2+6*gm(1,1)*gm(3,3))) 1309 cm(1,6,5)=cm(1,6,5)-(gm(1,1)*(3*gm(1,2)*gm(1,3)+6*gm(1,1)*gm(2,3))) 1310 cm(1,7,5)=cm(1,7,5)-(-3*gm(1,2)*gm(1,3)*gm(2,2)+7.5d0*gm(1,2)& 1311 & **2*gm(2,3)-1.5d0*gm(1,1)*gm(2,2)*gm(2,3)) 1312 cm(1,8,5)=cm(1,8,5)-(1.5d0*gm(1,3)**2*gm(2,3)+12*gm(1,2)*gm(1,3)& 1313 & *gm(3,3)-4.5d0*gm(1,1)*gm(2,3)*gm(3,3)) 1314 cm(1,9,5)=cm(1,9,5)-(-3*gm(1,3)**2*gm(2,2)+9*gm(1,2)*gm(1,3)*gm(2,3)& 1315 & +7.5d0*gm(1,2)**2*gm(3,3)+gm(1,1)*(-3*gm(2,3)**2-1.5d0*gm(2,2)& 1316 & *gm(3,3))) 1317 cm(1,10,5)=cm(1,10,5)-(gm(3,3)*(4.5d0*gm(1,3)**2-1.5d0*gm(1,1)& 1318 & *gm(3,3))) 1319 cm(1,1,6)=cm(1,1,6)-(3*gm(1,1)**2*gm(1,2)) 1320 cm(1,2,6)=cm(1,2,6)-(1.5d0*gm(1,2)**3+7.5d0*gm(1,1)*gm(1,2)*gm(2,2)) 1321 cm(1,3,6)=cm(1,3,6)-(12*gm(1,1)*gm(1,3)*gm(2,3)+gm(1,2)*(1.5d0*gm(1,3)& 1322 & **2-4.5d0*gm(1,1)*gm(3,3))) 1323 cm(1,4,6)=cm(1,4,6)-(3*gm(1,2)**2*gm(1,3)+12*gm(1,1)*gm(1,3)*gm(2,2)& 1324 & +3*gm(1,1)*gm(1,2)*gm(2,3)) 1325 cm(1,5,6)=cm(1,5,6)-(gm(1,1)*(3*gm(1,2)*gm(1,3)+6*gm(1,1)*gm(2,3))) 1326 cm(1,6,6)=cm(1,6,6)-(gm(1,1)*(3*gm(1,2)**2+6*gm(1,1)*gm(2,2))) 1327 cm(1,7,6)=cm(1,7,6)-(gm(2,2)*(4.5d0*gm(1,2)**2-1.5d0*gm(1,1)*gm(2,2))) 1328 cm(1,8,6)=cm(1,8,6)-(7.5d0*gm(1,3)**2*gm(2,2)+9*gm(1,2)*gm(1,3)& 1329 & *gm(2,3)-3*gm(1,2)**2*gm(3,3)+gm(1,1)*(-3*gm(2,3)**2-1.5d0*gm(2,2)& 1330 & *gm(3,3))) 1331 cm(1,9,6)=cm(1,9,6)-(12*gm(1,2)*gm(1,3)*gm(2,2)+1.5d0*gm(1,2)& 1332 & **2*gm(2,3)-4.5d0*gm(1,1)*gm(2,2)*gm(2,3)) 1333 cm(1,10,6)=cm(1,10,6)-(7.5d0*gm(1,3)**2*gm(2,3)-3*gm(1,2)*gm(1,3)& 1334 & *gm(3,3)-1.5d0*gm(1,1)*gm(2,3)*gm(3,3)) 1335 cm(1,1,7)=cm(1,1,7)-(2.5d0*gm(1,2)**3-1.5d0*gm(1,1)*gm(1,2)*gm(2,2)) 1336 cm(1,2,7)=cm(1,2,7)-(3*gm(1,2)*gm(2,2)**2) 1337 cm(1,3,7)=cm(1,3,7)-(-3*gm(1,3)*gm(2,2)*gm(2,3)+gm(1,2)*(7.5d0*gm(2,3)& 1338 & **2-1.5d0*gm(2,2)*gm(3,3))) 1339 cm(1,4,7)=cm(1,4,7)-(gm(2,2)*(-3*gm(1,3)*gm(2,2)+9*gm(1,2)*gm(2,3))) 1340 cm(1,5,7)=cm(1,5,7)-(-3*gm(1,2)*gm(1,3)*gm(2,2)+7.5d0*gm(1,2)& 1341 & **2*gm(2,3)-1.5d0*gm(1,1)*gm(2,2)*gm(2,3)) 1342 cm(1,6,7)=cm(1,6,7)-(gm(2,2)*(4.5d0*gm(1,2)**2-1.5d0*gm(1,1)*gm(2,2))) 1343 cm(1,7,7)=cm(1,7,7)-(gm(2,2)**3) 1344 cm(1,8,7)=cm(1,8,7)-(gm(2,2)*(4.5d0*gm(2,3)**2-1.5d0*gm(2,2)*gm(3,3))) 1345 cm(1,9,7)=cm(1,9,7)-(3*gm(2,2)**2*gm(2,3)) 1346 cm(1,10,7)=cm(1,10,7)-(2.5d0*gm(2,3)**3-1.5d0*gm(2,2)*gm(2,3)& 1347 & *gm(3,3)) 1348 cm(1,1,8)=cm(1,1,8)-(-3*gm(1,1)*gm(1,3)*gm(2,3)+gm(1,2)*(7.5d0*gm(1,3)& 1349 & **2-1.5d0*gm(1,1)*gm(3,3))) 1350 cm(1,2,8)=cm(1,2,8)-(12*gm(1,3)*gm(2,2)*gm(2,3)+gm(1,2)*(1.5d0*gm(2,3)& 1351 & **2-4.5d0*gm(2,2)*gm(3,3))) 1352 cm(1,3,8)=cm(1,3,8)-(gm(3,3)*(3*gm(1,3)*gm(2,3)+6*gm(1,2)*gm(3,3))) 1353 cm(1,4,8)=cm(1,4,8)-(3*gm(1,2)*gm(2,3)*gm(3,3)+gm(1,3)*(3*gm(2,3)& 1354 & **2+12*gm(2,2)*gm(3,3))) 1355 cm(1,5,8)=cm(1,5,8)-(1.5d0*gm(1,3)**2*gm(2,3)+12*gm(1,2)*gm(1,3)& 1356 & *gm(3,3)-4.5d0*gm(1,1)*gm(2,3)*gm(3,3)) 1357 cm(1,6,8)=cm(1,6,8)-(7.5d0*gm(1,3)**2*gm(2,2)+9*gm(1,2)*gm(1,3)& 1358 & *gm(2,3)-3*gm(1,2)**2*gm(3,3)+gm(1,1)*(-3*gm(2,3)**2-1.5d0*gm(2,2)& 1359 & *gm(3,3))) 1360 cm(1,7,8)=cm(1,7,8)-(gm(2,2)*(4.5d0*gm(2,3)**2-1.5d0*gm(2,2)*gm(3,3))) 1361 cm(1,8,8)=cm(1,8,8)-(gm(3,3)*(3*gm(2,3)**2+6*gm(2,2)*gm(3,3))) 1362 cm(1,9,8)=cm(1,9,8)-(1.5d0*gm(2,3)**3+7.5d0*gm(2,2)*gm(2,3)*gm(3,3)) 1363 cm(1,10,8)=cm(1,10,8)-(3*gm(2,3)*gm(3,3)**2) 1364 cm(1,1,9)=cm(1,1,9)-(7.5d0*gm(1,2)**2*gm(1,3)-1.5d0*gm(1,1)*gm(1,3)& 1365 & *gm(2,2)-3*gm(1,1)*gm(1,2)*gm(2,3)) 1366 cm(1,2,9)=cm(1,2,9)-(gm(2,2)*(6*gm(1,3)*gm(2,2)+3*gm(1,2)*gm(2,3))) 1367 cm(1,3,9)=cm(1,3,9)-(12*gm(1,2)*gm(2,3)*gm(3,3)+gm(1,3)*(1.5d0*gm(2,3)& 1368 & **2-4.5d0*gm(2,2)*gm(3,3))) 1369 cm(1,4,9)=cm(1,4,9)-(3*gm(1,3)*gm(2,2)*gm(2,3)+gm(1,2)*(3*gm(2,3)& 1370 & **2+12*gm(2,2)*gm(3,3))) 1371 cm(1,5,9)=cm(1,5,9)-(-3*gm(1,3)**2*gm(2,2)+9*gm(1,2)*gm(1,3)*gm(2,3)& 1372 & +7.5d0*gm(1,2)**2*gm(3,3)+gm(1,1)*(-3*gm(2,3)**2-1.5d0*gm(2,2)& 1373 & *gm(3,3))) 1374 cm(1,6,9)=cm(1,6,9)-(12*gm(1,2)*gm(1,3)*gm(2,2)+1.5d0*gm(1,2)& 1375 & **2*gm(2,3)-4.5d0*gm(1,1)*gm(2,2)*gm(2,3)) 1376 cm(1,7,9)=cm(1,7,9)-(3*gm(2,2)**2*gm(2,3)) 1377 cm(1,8,9)=cm(1,8,9)-(1.5d0*gm(2,3)**3+7.5d0*gm(2,2)*gm(2,3)*gm(3,3)) 1378 cm(1,9,9)=cm(1,9,9)-(gm(2,2)*(3*gm(2,3)**2+6*gm(2,2)*gm(3,3))) 1379 cm(1,10,9)=cm(1,10,9)-(gm(3,3)*(4.5d0*gm(2,3)**2-1.5d0*gm(2,2)& 1380 & *gm(3,3))) 1381 cm(1,1,10)=cm(1,1,10)-(2.5d0*gm(1,3)**3-1.5d0*gm(1,1)*gm(1,3)& 1382 & *gm(3,3)) 1383 cm(1,2,10)=cm(1,2,10)-(-3*gm(1,2)*gm(2,3)*gm(3,3)+gm(1,3)*(7.5d0*gm(2,3)& 1384 & **2-1.5d0*gm(2,2)*gm(3,3))) 1385 cm(1,3,10)=cm(1,3,10)-(3*gm(1,3)*gm(3,3)**2) 1386 cm(1,4,10)=cm(1,4,10)-(gm(3,3)*(9*gm(1,3)*gm(2,3)-3*gm(1,2)*gm(3,3))) 1387 cm(1,5,10)=cm(1,5,10)-(gm(3,3)*(4.5d0*gm(1,3)**2-1.5d0*gm(1,1)& 1388 & *gm(3,3))) 1389 cm(1,6,10)=cm(1,6,10)-(7.5d0*gm(1,3)**2*gm(2,3)-3*gm(1,2)*gm(1,3)& 1390 & *gm(3,3)-1.5d0*gm(1,1)*gm(2,3)*gm(3,3)) 1391 cm(1,7,10)=cm(1,7,10)-(2.5d0*gm(2,3)**3-1.5d0*gm(2,2)*gm(2,3)& 1392 & *gm(3,3)) 1393 cm(1,8,10)=cm(1,8,10)-(3*gm(2,3)*gm(3,3)**2) 1394 cm(1,9,10)=cm(1,9,10)-(gm(3,3)*(4.5d0*gm(2,3)**2-1.5d0*gm(2,2)& 1395 & *gm(3,3))) 1396 cm(1,10,10)=cm(1,10,10)-(gm(3,3)**3) 1397 cm(2,1,2)=cm(2,1,2)-(3*gm(1,1)**2*gm(1,2)) 1398 cm(2,2,2)=cm(2,2,2)-(1.5d0*gm(1,2)**3+7.5d0*gm(1,1)*gm(1,2)*gm(2,2)) 1399 cm(2,3,2)=cm(2,3,2)-(12*gm(1,1)*gm(1,3)*gm(2,3)+gm(1,2)*(1.5d0*gm(1,3)& 1400 & **2-4.5d0*gm(1,1)*gm(3,3))) 1401 cm(2,4,2)=cm(2,4,2)-(3*gm(1,2)**2*gm(1,3)+12*gm(1,1)*gm(1,3)*gm(2,2)& 1402 & +3*gm(1,1)*gm(1,2)*gm(2,3)) 1403 cm(2,5,2)=cm(2,5,2)-(gm(1,1)*(3*gm(1,2)*gm(1,3)+6*gm(1,1)*gm(2,3))) 1404 cm(2,6,2)=cm(2,6,2)-(gm(1,1)*(3*gm(1,2)**2+6*gm(1,1)*gm(2,2))) 1405 cm(2,7,2)=cm(2,7,2)-(gm(2,2)*(4.5d0*gm(1,2)**2-1.5d0*gm(1,1)*gm(2,2))) 1406 cm(2,8,2)=cm(2,8,2)-(7.5d0*gm(1,3)**2*gm(2,2)+9*gm(1,2)*gm(1,3)& 1407 & *gm(2,3)-3*gm(1,2)**2*gm(3,3)+gm(1,1)*(-3*gm(2,3)**2-1.5d0*gm(2,2)& 1408 & *gm(3,3))) 1409 cm(2,9,2)=cm(2,9,2)-(12*gm(1,2)*gm(1,3)*gm(2,2)+1.5d0*gm(1,2)& 1410 & **2*gm(2,3)-4.5d0*gm(1,1)*gm(2,2)*gm(2,3)) 1411 cm(2,10,2)=cm(2,10,2)-(7.5d0*gm(1,3)**2*gm(2,3)-3*gm(1,2)*gm(1,3)& 1412 & *gm(3,3)-1.5d0*gm(1,1)*gm(2,3)*gm(3,3)) 1413 cm(2,1,4)=cm(2,1,4)-(3*gm(1,1)**2*gm(1,3)) 1414 cm(2,2,4)=cm(2,2,4)-(1.5d0*gm(1,2)**2*gm(1,3)-4.5d0*gm(1,1)*gm(1,3)& 1415 & *gm(2,2)+12*gm(1,1)*gm(1,2)*gm(2,3)) 1416 cm(2,3,4)=cm(2,3,4)-(1.5d0*gm(1,3)**3+7.5d0*gm(1,1)*gm(1,3)*gm(3,3)) 1417 cm(2,4,4)=cm(2,4,4)-(3*gm(1,1)*gm(1,3)*gm(2,3)+gm(1,2)*(3*gm(1,3)& 1418 & **2+12*gm(1,1)*gm(3,3))) 1419 cm(2,5,4)=cm(2,5,4)-(gm(1,1)*(3*gm(1,3)**2+6*gm(1,1)*gm(3,3))) 1420 cm(2,6,4)=cm(2,6,4)-(gm(1,1)*(3*gm(1,2)*gm(1,3)+6*gm(1,1)*gm(2,3))) 1421 cm(2,7,4)=cm(2,7,4)-(-3*gm(1,2)*gm(1,3)*gm(2,2)+7.5d0*gm(1,2)& 1422 & **2*gm(2,3)-1.5d0*gm(1,1)*gm(2,2)*gm(2,3)) 1423 cm(2,8,4)=cm(2,8,4)-(1.5d0*gm(1,3)**2*gm(2,3)+12*gm(1,2)*gm(1,3)& 1424 & *gm(3,3)-4.5d0*gm(1,1)*gm(2,3)*gm(3,3)) 1425 cm(2,9,4)=cm(2,9,4)-(-3*gm(1,3)**2*gm(2,2)+9*gm(1,2)*gm(1,3)*gm(2,3)& 1426 & +7.5d0*gm(1,2)**2*gm(3,3)+gm(1,1)*(-3*gm(2,3)**2-1.5d0*gm(2,2)& 1427 & *gm(3,3))) 1428 cm(2,10,4)=cm(2,10,4)-(gm(3,3)*(4.5d0*gm(1,3)**2-1.5d0*gm(1,1)& 1429 & *gm(3,3))) 1430 cm(2,1,6)=cm(2,1,6)-(gm(1,1)**3) 1431 cm(2,2,6)=cm(2,2,6)-(gm(1,1)*(4.5d0*gm(1,2)**2-1.5d0*gm(1,1)*gm(2,2))) 1432 cm(2,3,6)=cm(2,3,6)-(gm(1,1)*(4.5d0*gm(1,3)**2-1.5d0*gm(1,1)*gm(3,3))) 1433 cm(2,4,6)=cm(2,4,6)-(gm(1,1)*(9*gm(1,2)*gm(1,3)-3*gm(1,1)*gm(2,3))) 1434 cm(2,5,6)=cm(2,5,6)-(3*gm(1,1)**2*gm(1,3)) 1435 cm(2,6,6)=cm(2,6,6)-(3*gm(1,1)**2*gm(1,2)) 1436 cm(2,7,6)=cm(2,7,6)-(2.5d0*gm(1,2)**3-1.5d0*gm(1,1)*gm(1,2)*gm(2,2)) 1437 cm(2,8,6)=cm(2,8,6)-(-3*gm(1,1)*gm(1,3)*gm(2,3)+gm(1,2)*(7.5d0*gm(1,3)& 1438 & **2-1.5d0*gm(1,1)*gm(3,3))) 1439 cm(2,9,6)=cm(2,9,6)-(7.5d0*gm(1,2)**2*gm(1,3)-1.5d0*gm(1,1)*gm(1,3)& 1440 & *gm(2,2)-3*gm(1,1)*gm(1,2)*gm(2,3)) 1441 cm(2,10,6)=cm(2,10,6)-(2.5d0*gm(1,3)**3-1.5d0*gm(1,1)*gm(1,3)& 1442 & *gm(3,3)) 1443 cm(2,1,7)=cm(2,1,7)-(gm(1,1)*(4.5d0*gm(1,2)**2-1.5d0*gm(1,1)*gm(2,2))) 1444 cm(2,2,7)=cm(2,2,7)-(gm(2,2)*(3*gm(1,2)**2+6*gm(1,1)*gm(2,2))) 1445 cm(2,3,7)=cm(2,3,7)-(-3*gm(1,3)**2*gm(2,2)+9*gm(1,2)*gm(1,3)*gm(2,3)& 1446 & -3*gm(1,2)**2*gm(3,3)+gm(1,1)*(7.5d0*gm(2,3)**2-1.5d0*gm(2,2)& 1447 & *gm(3,3))) 1448 cm(2,4,7)=cm(2,4,7)-(3*gm(1,2)*gm(1,3)*gm(2,2)+3*gm(1,2)**2*gm(2,3)& 1449 & +12*gm(1,1)*gm(2,2)*gm(2,3)) 1450 cm(2,5,7)=cm(2,5,7)-(1.5d0*gm(1,2)**2*gm(1,3)-4.5d0*gm(1,1)*gm(1,3)& 1451 & *gm(2,2)+12*gm(1,1)*gm(1,2)*gm(2,3)) 1452 cm(2,6,7)=cm(2,6,7)-(1.5d0*gm(1,2)**3+7.5d0*gm(1,1)*gm(1,2)*gm(2,2)) 1453 cm(2,7,7)=cm(2,7,7)-(3*gm(1,2)*gm(2,2)**2) 1454 cm(2,8,7)=cm(2,8,7)-(12*gm(1,3)*gm(2,2)*gm(2,3)+gm(1,2)*(1.5d0*gm(2,3)& 1455 & **2-4.5d0*gm(2,2)*gm(3,3))) 1456 cm(2,9,7)=cm(2,9,7)-(gm(2,2)*(6*gm(1,3)*gm(2,2)+3*gm(1,2)*gm(2,3))) 1457 cm(2,10,7)=cm(2,10,7)-(-3*gm(1,2)*gm(2,3)*gm(3,3)+gm(1,3)*(7.5d0*gm(2,3)& 1458 & **2-1.5d0*gm(2,2)*gm(3,3))) 1459 cm(2,1,8)=cm(2,1,8)-(gm(1,1)*(4.5d0*gm(1,3)**2-1.5d0*gm(1,1)*gm(3,3))) 1460 cm(2,2,8)=cm(2,2,8)-(-3*gm(1,3)**2*gm(2,2)+9*gm(1,2)*gm(1,3)*gm(2,3)& 1461 & -3*gm(1,2)**2*gm(3,3)+gm(1,1)*(7.5d0*gm(2,3)**2-1.5d0*gm(2,2)& 1462 & *gm(3,3))) 1463 cm(2,3,8)=cm(2,3,8)-(gm(3,3)*(3*gm(1,3)**2+6*gm(1,1)*gm(3,3))) 1464 cm(2,4,8)=cm(2,4,8)-(3*gm(1,3)**2*gm(2,3)+3*gm(1,2)*gm(1,3)*gm(3,3)& 1465 & +12*gm(1,1)*gm(2,3)*gm(3,3)) 1466 cm(2,5,8)=cm(2,5,8)-(1.5d0*gm(1,3)**3+7.5d0*gm(1,1)*gm(1,3)*gm(3,3)) 1467 cm(2,6,8)=cm(2,6,8)-(12*gm(1,1)*gm(1,3)*gm(2,3)+gm(1,2)*(1.5d0*gm(1,3)& 1468 & **2-4.5d0*gm(1,1)*gm(3,3))) 1469 cm(2,7,8)=cm(2,7,8)-(-3*gm(1,3)*gm(2,2)*gm(2,3)+gm(1,2)*(7.5d0*gm(2,3)& 1470 & **2-1.5d0*gm(2,2)*gm(3,3))) 1471 cm(2,8,8)=cm(2,8,8)-(gm(3,3)*(3*gm(1,3)*gm(2,3)+6*gm(1,2)*gm(3,3))) 1472 cm(2,9,8)=cm(2,9,8)-(12*gm(1,2)*gm(2,3)*gm(3,3)+gm(1,3)*(1.5d0*gm(2,3)& 1473 & **2-4.5d0*gm(2,2)*gm(3,3))) 1474 cm(2,10,8)=cm(2,10,8)-(3*gm(1,3)*gm(3,3)**2) 1475 cm(2,1,9)=cm(2,1,9)-(gm(1,1)*(9*gm(1,2)*gm(1,3)-3*gm(1,1)*gm(2,3))) 1476 cm(2,2,9)=cm(2,2,9)-(3*gm(1,2)*gm(1,3)*gm(2,2)+3*gm(1,2)**2*gm(2,3)& 1477 & +12*gm(1,1)*gm(2,2)*gm(2,3)) 1478 cm(2,3,9)=cm(2,3,9)-(3*gm(1,3)**2*gm(2,3)+3*gm(1,2)*gm(1,3)*gm(3,3)& 1479 & +12*gm(1,1)*gm(2,3)*gm(3,3)) 1480 cm(2,4,9)=cm(2,4,9)-(9*gm(1,3)**2*gm(2,2)-6*gm(1,2)*gm(1,3)*gm(2,3)& 1481 & +9*gm(1,2)**2*gm(3,3)+gm(1,1)*(9*gm(2,3)**2+15*gm(2,2)*gm(3,3))) 1482 cm(2,5,9)=cm(2,5,9)-(3*gm(1,1)*gm(1,3)*gm(2,3)+gm(1,2)*(3*gm(1,3)& 1483 & **2+12*gm(1,1)*gm(3,3))) 1484 cm(2,6,9)=cm(2,6,9)-(3*gm(1,2)**2*gm(1,3)+12*gm(1,1)*gm(1,3)*gm(2,2)& 1485 & +3*gm(1,1)*gm(1,2)*gm(2,3)) 1486 cm(2,7,9)=cm(2,7,9)-(gm(2,2)*(-3*gm(1,3)*gm(2,2)+9*gm(1,2)*gm(2,3))) 1487 cm(2,8,9)=cm(2,8,9)-(3*gm(1,2)*gm(2,3)*gm(3,3)+gm(1,3)*(3*gm(2,3)& 1488 & **2+12*gm(2,2)*gm(3,3))) 1489 cm(2,9,9)=cm(2,9,9)-(3*gm(1,3)*gm(2,2)*gm(2,3)+gm(1,2)*(3*gm(2,3)& 1490 & **2+12*gm(2,2)*gm(3,3))) 1491 cm(2,10,9)=cm(2,10,9)-(gm(3,3)*(9*gm(1,3)*gm(2,3)-3*gm(1,2)*gm(3,3))) 1492 cm(2,1,11)=cm(2,1,11)-(2.5d0*gm(1,2)**3-1.5d0*gm(1,1)*gm(1,2)& 1493 & *gm(2,2)) 1494 cm(2,2,11)=cm(2,2,11)-(3*gm(1,2)*gm(2,2)**2) 1495 cm(2,3,11)=cm(2,3,11)-(-3*gm(1,3)*gm(2,2)*gm(2,3)+gm(1,2)*(7.5d0*gm(2,3)& 1496 & **2-1.5d0*gm(2,2)*gm(3,3))) 1497 cm(2,4,11)=cm(2,4,11)-(gm(2,2)*(-3*gm(1,3)*gm(2,2)+9*gm(1,2)*gm(2,3))) 1498 cm(2,5,11)=cm(2,5,11)-(-3*gm(1,2)*gm(1,3)*gm(2,2)+7.5d0*gm(1,2)& 1499 & **2*gm(2,3)-1.5d0*gm(1,1)*gm(2,2)*gm(2,3)) 1500 cm(2,6,11)=cm(2,6,11)-(gm(2,2)*(4.5d0*gm(1,2)**2-1.5d0*gm(1,1)& 1501 & *gm(2,2))) 1502 cm(2,7,11)=cm(2,7,11)-(gm(2,2)**3) 1503 cm(2,8,11)=cm(2,8,11)-(gm(2,2)*(4.5d0*gm(2,3)**2-1.5d0*gm(2,2)& 1504 & *gm(3,3))) 1505 cm(2,9,11)=cm(2,9,11)-(3*gm(2,2)**2*gm(2,3)) 1506 cm(2,10,11)=cm(2,10,11)-(2.5d0*gm(2,3)**3-1.5d0*gm(2,2)*gm(2,3)& 1507 & *gm(3,3)) 1508 cm(2,1,12)=cm(2,1,12)-(-3*gm(1,1)*gm(1,3)*gm(2,3)+gm(1,2)*(7.5d0*gm(1,3)& 1509 & **2-1.5d0*gm(1,1)*gm(3,3))) 1510 cm(2,2,12)=cm(2,2,12)-(12*gm(1,3)*gm(2,2)*gm(2,3)+gm(1,2)*(1.5d0*gm(2,3)& 1511 & **2-4.5d0*gm(2,2)*gm(3,3))) 1512 cm(2,3,12)=cm(2,3,12)-(gm(3,3)*(3*gm(1,3)*gm(2,3)+6*gm(1,2)*gm(3,3))) 1513 cm(2,4,12)=cm(2,4,12)-(3*gm(1,2)*gm(2,3)*gm(3,3)+gm(1,3)*(3*gm(2,3)& 1514 & **2+12*gm(2,2)*gm(3,3))) 1515 cm(2,5,12)=cm(2,5,12)-(1.5d0*gm(1,3)**2*gm(2,3)+12*gm(1,2)*gm(1,3)& 1516 & *gm(3,3)-4.5d0*gm(1,1)*gm(2,3)*gm(3,3)) 1517 cm(2,6,12)=cm(2,6,12)-(7.5d0*gm(1,3)**2*gm(2,2)+9*gm(1,2)*gm(1,3)& 1518 & *gm(2,3)-3*gm(1,2)**2*gm(3,3)+gm(1,1)*(-3*gm(2,3)**2-1.5d0*gm(2,2)& 1519 & *gm(3,3))) 1520 cm(2,7,12)=cm(2,7,12)-(gm(2,2)*(4.5d0*gm(2,3)**2-1.5d0*gm(2,2)& 1521 & *gm(3,3))) 1522 cm(2,8,12)=cm(2,8,12)-(gm(3,3)*(3*gm(2,3)**2+6*gm(2,2)*gm(3,3))) 1523 cm(2,9,12)=cm(2,9,12)-(1.5d0*gm(2,3)**3+7.5d0*gm(2,2)*gm(2,3)& 1524 & *gm(3,3)) 1525 cm(2,10,12)=cm(2,10,12)-(3*gm(2,3)*gm(3,3)**2) 1526 cm(2,1,13)=cm(2,1,13)-(7.5d0*gm(1,2)**2*gm(1,3)-1.5d0*gm(1,1)& 1527 & *gm(1,3)*gm(2,2)-3*gm(1,1)*gm(1,2)*gm(2,3)) 1528 cm(2,2,13)=cm(2,2,13)-(gm(2,2)*(6*gm(1,3)*gm(2,2)+3*gm(1,2)*gm(2,3))) 1529 cm(2,3,13)=cm(2,3,13)-(12*gm(1,2)*gm(2,3)*gm(3,3)+gm(1,3)*(1.5d0*gm(2,3)& 1530 & **2-4.5d0*gm(2,2)*gm(3,3))) 1531 cm(2,4,13)=cm(2,4,13)-(3*gm(1,3)*gm(2,2)*gm(2,3)+gm(1,2)*(3*gm(2,3)& 1532 & **2+12*gm(2,2)*gm(3,3))) 1533 cm(2,5,13)=cm(2,5,13)-(-3*gm(1,3)**2*gm(2,2)+9*gm(1,2)*gm(1,3)& 1534 & *gm(2,3)+7.5d0*gm(1,2)**2*gm(3,3)+gm(1,1)*(-3*gm(2,3)**2-1.5d0*gm(2,2)& 1535 & *gm(3,3))) 1536 cm(2,6,13)=cm(2,6,13)-(12*gm(1,2)*gm(1,3)*gm(2,2)+1.5d0*gm(1,2)& 1537 & **2*gm(2,3)-4.5d0*gm(1,1)*gm(2,2)*gm(2,3)) 1538 cm(2,7,13)=cm(2,7,13)-(3*gm(2,2)**2*gm(2,3)) 1539 cm(2,8,13)=cm(2,8,13)-(1.5d0*gm(2,3)**3+7.5d0*gm(2,2)*gm(2,3)& 1540 & *gm(3,3)) 1541 cm(2,9,13)=cm(2,9,13)-(gm(2,2)*(3*gm(2,3)**2+6*gm(2,2)*gm(3,3))) 1542 cm(2,10,13)=cm(2,10,13)-(gm(3,3)*(4.5d0*gm(2,3)**2-1.5d0*gm(2,2)& 1543 & *gm(3,3))) 1544 cm(2,1,14)=cm(2,1,14)-(2.5d0*gm(1,3)**3-1.5d0*gm(1,1)*gm(1,3)& 1545 & *gm(3,3)) 1546 cm(2,2,14)=cm(2,2,14)-(-3*gm(1,2)*gm(2,3)*gm(3,3)+gm(1,3)*(7.5d0*gm(2,3)& 1547 & **2-1.5d0*gm(2,2)*gm(3,3))) 1548 cm(2,3,14)=cm(2,3,14)-(3*gm(1,3)*gm(3,3)**2) 1549 cm(2,4,14)=cm(2,4,14)-(gm(3,3)*(9*gm(1,3)*gm(2,3)-3*gm(1,2)*gm(3,3))) 1550 cm(2,5,14)=cm(2,5,14)-(gm(3,3)*(4.5d0*gm(1,3)**2-1.5d0*gm(1,1)& 1551 & *gm(3,3))) 1552 cm(2,6,14)=cm(2,6,14)-(7.5d0*gm(1,3)**2*gm(2,3)-3*gm(1,2)*gm(1,3)& 1553 & *gm(3,3)-1.5d0*gm(1,1)*gm(2,3)*gm(3,3)) 1554 cm(2,7,14)=cm(2,7,14)-(2.5d0*gm(2,3)**3-1.5d0*gm(2,2)*gm(2,3)& 1555 & *gm(3,3)) 1556 cm(2,8,14)=cm(2,8,14)-(3*gm(2,3)*gm(3,3)**2) 1557 cm(2,9,14)=cm(2,9,14)-(gm(3,3)*(4.5d0*gm(2,3)**2-1.5d0*gm(2,2)& 1558 & *gm(3,3))) 1559 cm(2,10,14)=cm(2,10,14)-(gm(3,3)**3) 1560 cm(3,1,3)=cm(3,1,3)-(3*gm(1,1)**2*gm(1,3)) 1561 cm(3,2,3)=cm(3,2,3)-(1.5d0*gm(1,2)**2*gm(1,3)-4.5d0*gm(1,1)*gm(1,3)& 1562 & *gm(2,2)+12*gm(1,1)*gm(1,2)*gm(2,3)) 1563 cm(3,3,3)=cm(3,3,3)-(1.5d0*gm(1,3)**3+7.5d0*gm(1,1)*gm(1,3)*gm(3,3)) 1564 cm(3,4,3)=cm(3,4,3)-(3*gm(1,1)*gm(1,3)*gm(2,3)+gm(1,2)*(3*gm(1,3)& 1565 & **2+12*gm(1,1)*gm(3,3))) 1566 cm(3,5,3)=cm(3,5,3)-(gm(1,1)*(3*gm(1,3)**2+6*gm(1,1)*gm(3,3))) 1567 cm(3,6,3)=cm(3,6,3)-(gm(1,1)*(3*gm(1,2)*gm(1,3)+6*gm(1,1)*gm(2,3))) 1568 cm(3,7,3)=cm(3,7,3)-(-3*gm(1,2)*gm(1,3)*gm(2,2)+7.5d0*gm(1,2)& 1569 & **2*gm(2,3)-1.5d0*gm(1,1)*gm(2,2)*gm(2,3)) 1570 cm(3,8,3)=cm(3,8,3)-(1.5d0*gm(1,3)**2*gm(2,3)+12*gm(1,2)*gm(1,3)& 1571 & *gm(3,3)-4.5d0*gm(1,1)*gm(2,3)*gm(3,3)) 1572 cm(3,9,3)=cm(3,9,3)-(-3*gm(1,3)**2*gm(2,2)+9*gm(1,2)*gm(1,3)*gm(2,3)& 1573 & +7.5d0*gm(1,2)**2*gm(3,3)+gm(1,1)*(-3*gm(2,3)**2-1.5d0*gm(2,2)& 1574 & *gm(3,3))) 1575 cm(3,10,3)=cm(3,10,3)-(gm(3,3)*(4.5d0*gm(1,3)**2-1.5d0*gm(1,1)& 1576 & *gm(3,3))) 1577 cm(3,1,4)=cm(3,1,4)-(3*gm(1,1)**2*gm(1,2)) 1578 cm(3,2,4)=cm(3,2,4)-(1.5d0*gm(1,2)**3+7.5d0*gm(1,1)*gm(1,2)*gm(2,2)) 1579 cm(3,3,4)=cm(3,3,4)-(12*gm(1,1)*gm(1,3)*gm(2,3)+gm(1,2)*(1.5d0*gm(1,3)& 1580 & **2-4.5d0*gm(1,1)*gm(3,3))) 1581 cm(3,4,4)=cm(3,4,4)-(3*gm(1,2)**2*gm(1,3)+12*gm(1,1)*gm(1,3)*gm(2,2)& 1582 & +3*gm(1,1)*gm(1,2)*gm(2,3)) 1583 cm(3,5,4)=cm(3,5,4)-(gm(1,1)*(3*gm(1,2)*gm(1,3)+6*gm(1,1)*gm(2,3))) 1584 cm(3,6,4)=cm(3,6,4)-(gm(1,1)*(3*gm(1,2)**2+6*gm(1,1)*gm(2,2))) 1585 cm(3,7,4)=cm(3,7,4)-(gm(2,2)*(4.5d0*gm(1,2)**2-1.5d0*gm(1,1)*gm(2,2))) 1586 cm(3,8,4)=cm(3,8,4)-(7.5d0*gm(1,3)**2*gm(2,2)+9*gm(1,2)*gm(1,3)& 1587 & *gm(2,3)-3*gm(1,2)**2*gm(3,3)+gm(1,1)*(-3*gm(2,3)**2-1.5d0*gm(2,2)& 1588 & *gm(3,3))) 1589 cm(3,9,4)=cm(3,9,4)-(12*gm(1,2)*gm(1,3)*gm(2,2)+1.5d0*gm(1,2)& 1590 & **2*gm(2,3)-4.5d0*gm(1,1)*gm(2,2)*gm(2,3)) 1591 cm(3,10,4)=cm(3,10,4)-(7.5d0*gm(1,3)**2*gm(2,3)-3*gm(1,2)*gm(1,3)& 1592 & *gm(3,3)-1.5d0*gm(1,1)*gm(2,3)*gm(3,3)) 1593 cm(3,1,5)=cm(3,1,5)-(gm(1,1)**3) 1594 cm(3,2,5)=cm(3,2,5)-(gm(1,1)*(4.5d0*gm(1,2)**2-1.5d0*gm(1,1)*gm(2,2))) 1595 cm(3,3,5)=cm(3,3,5)-(gm(1,1)*(4.5d0*gm(1,3)**2-1.5d0*gm(1,1)*gm(3,3))) 1596 cm(3,4,5)=cm(3,4,5)-(gm(1,1)*(9*gm(1,2)*gm(1,3)-3*gm(1,1)*gm(2,3))) 1597 cm(3,5,5)=cm(3,5,5)-(3*gm(1,1)**2*gm(1,3)) 1598 cm(3,6,5)=cm(3,6,5)-(3*gm(1,1)**2*gm(1,2)) 1599 cm(3,7,5)=cm(3,7,5)-(2.5d0*gm(1,2)**3-1.5d0*gm(1,1)*gm(1,2)*gm(2,2)) 1600 cm(3,8,5)=cm(3,8,5)-(-3*gm(1,1)*gm(1,3)*gm(2,3)+gm(1,2)*(7.5d0*gm(1,3)& 1601 & **2-1.5d0*gm(1,1)*gm(3,3))) 1602 cm(3,9,5)=cm(3,9,5)-(7.5d0*gm(1,2)**2*gm(1,3)-1.5d0*gm(1,1)*gm(1,3)& 1603 & *gm(2,2)-3*gm(1,1)*gm(1,2)*gm(2,3)) 1604 cm(3,10,5)=cm(3,10,5)-(2.5d0*gm(1,3)**3-1.5d0*gm(1,1)*gm(1,3)& 1605 & *gm(3,3)) 1606 cm(3,1,8)=cm(3,1,8)-(gm(1,1)*(9*gm(1,2)*gm(1,3)-3*gm(1,1)*gm(2,3))) 1607 cm(3,2,8)=cm(3,2,8)-(3*gm(1,2)*gm(1,3)*gm(2,2)+3*gm(1,2)**2*gm(2,3)& 1608 & +12*gm(1,1)*gm(2,2)*gm(2,3)) 1609 cm(3,3,8)=cm(3,3,8)-(3*gm(1,3)**2*gm(2,3)+3*gm(1,2)*gm(1,3)*gm(3,3)& 1610 & +12*gm(1,1)*gm(2,3)*gm(3,3)) 1611 cm(3,4,8)=cm(3,4,8)-(9*gm(1,3)**2*gm(2,2)-6*gm(1,2)*gm(1,3)*gm(2,3)& 1612 & +9*gm(1,2)**2*gm(3,3)+gm(1,1)*(9*gm(2,3)**2+15*gm(2,2)*gm(3,3))) 1613 cm(3,5,8)=cm(3,5,8)-(3*gm(1,1)*gm(1,3)*gm(2,3)+gm(1,2)*(3*gm(1,3)& 1614 & **2+12*gm(1,1)*gm(3,3))) 1615 cm(3,6,8)=cm(3,6,8)-(3*gm(1,2)**2*gm(1,3)+12*gm(1,1)*gm(1,3)*gm(2,2)& 1616 & +3*gm(1,1)*gm(1,2)*gm(2,3)) 1617 cm(3,7,8)=cm(3,7,8)-(gm(2,2)*(-3*gm(1,3)*gm(2,2)+9*gm(1,2)*gm(2,3))) 1618 cm(3,8,8)=cm(3,8,8)-(3*gm(1,2)*gm(2,3)*gm(3,3)+gm(1,3)*(3*gm(2,3)& 1619 & **2+12*gm(2,2)*gm(3,3))) 1620 cm(3,9,8)=cm(3,9,8)-(3*gm(1,3)*gm(2,2)*gm(2,3)+gm(1,2)*(3*gm(2,3)& 1621 & **2+12*gm(2,2)*gm(3,3))) 1622 cm(3,10,8)=cm(3,10,8)-(gm(3,3)*(9*gm(1,3)*gm(2,3)-3*gm(1,2)*gm(3,3))) 1623 cm(3,1,9)=cm(3,1,9)-(gm(1,1)*(4.5d0*gm(1,2)**2-1.5d0*gm(1,1)*gm(2,2))) 1624 cm(3,2,9)=cm(3,2,9)-(gm(2,2)*(3*gm(1,2)**2+6*gm(1,1)*gm(2,2))) 1625 cm(3,3,9)=cm(3,3,9)-(-3*gm(1,3)**2*gm(2,2)+9*gm(1,2)*gm(1,3)*gm(2,3)& 1626 & -3*gm(1,2)**2*gm(3,3)+gm(1,1)*(7.5d0*gm(2,3)**2-1.5d0*gm(2,2)& 1627 & *gm(3,3))) 1628 cm(3,4,9)=cm(3,4,9)-(3*gm(1,2)*gm(1,3)*gm(2,2)+3*gm(1,2)**2*gm(2,3)& 1629 & +12*gm(1,1)*gm(2,2)*gm(2,3)) 1630 cm(3,5,9)=cm(3,5,9)-(1.5d0*gm(1,2)**2*gm(1,3)-4.5d0*gm(1,1)*gm(1,3)& 1631 & *gm(2,2)+12*gm(1,1)*gm(1,2)*gm(2,3)) 1632 cm(3,6,9)=cm(3,6,9)-(1.5d0*gm(1,2)**3+7.5d0*gm(1,1)*gm(1,2)*gm(2,2)) 1633 cm(3,7,9)=cm(3,7,9)-(3*gm(1,2)*gm(2,2)**2) 1634 cm(3,8,9)=cm(3,8,9)-(12*gm(1,3)*gm(2,2)*gm(2,3)+gm(1,2)*(1.5d0*gm(2,3)& 1635 & **2-4.5d0*gm(2,2)*gm(3,3))) 1636 cm(3,9,9)=cm(3,9,9)-(gm(2,2)*(6*gm(1,3)*gm(2,2)+3*gm(1,2)*gm(2,3))) 1637 cm(3,10,9)=cm(3,10,9)-(-3*gm(1,2)*gm(2,3)*gm(3,3)+gm(1,3)*(7.5d0*gm(2,3)& 1638 & **2-1.5d0*gm(2,2)*gm(3,3))) 1639 cm(3,1,10)=cm(3,1,10)-(gm(1,1)*(4.5d0*gm(1,3)**2-1.5d0*gm(1,1)& 1640 & *gm(3,3))) 1641 cm(3,2,10)=cm(3,2,10)-(-3*gm(1,3)**2*gm(2,2)+9*gm(1,2)*gm(1,3)& 1642 & *gm(2,3)-3*gm(1,2)**2*gm(3,3)+gm(1,1)*(7.5d0*gm(2,3)**2-1.5d0*gm(2,2)& 1643 & *gm(3,3))) 1644 cm(3,3,10)=cm(3,3,10)-(gm(3,3)*(3*gm(1,3)**2+6*gm(1,1)*gm(3,3))) 1645 cm(3,4,10)=cm(3,4,10)-(3*gm(1,3)**2*gm(2,3)+3*gm(1,2)*gm(1,3)& 1646 & *gm(3,3)+12*gm(1,1)*gm(2,3)*gm(3,3)) 1647 cm(3,5,10)=cm(3,5,10)-(1.5d0*gm(1,3)**3+7.5d0*gm(1,1)*gm(1,3)& 1648 & *gm(3,3)) 1649 cm(3,6,10)=cm(3,6,10)-(12*gm(1,1)*gm(1,3)*gm(2,3)+gm(1,2)*(1.5d0*gm(1,3)& 1650 & **2-4.5d0*gm(1,1)*gm(3,3))) 1651 cm(3,7,10)=cm(3,7,10)-(-3*gm(1,3)*gm(2,2)*gm(2,3)+gm(1,2)*(7.5d0*gm(2,3)& 1652 & **2-1.5d0*gm(2,2)*gm(3,3))) 1653 cm(3,8,10)=cm(3,8,10)-(gm(3,3)*(3*gm(1,3)*gm(2,3)+6*gm(1,2)*gm(3,3))) 1654 cm(3,9,10)=cm(3,9,10)-(12*gm(1,2)*gm(2,3)*gm(3,3)+gm(1,3)*(1.5d0*gm(2,3)& 1655 & **2-4.5d0*gm(2,2)*gm(3,3))) 1656 cm(3,10,10)=cm(3,10,10)-(3*gm(1,3)*gm(3,3)**2) 1657 cm(3,1,12)=cm(3,1,12)-(7.5d0*gm(1,2)**2*gm(1,3)-1.5d0*gm(1,1)& 1658 & *gm(1,3)*gm(2,2)-3*gm(1,1)*gm(1,2)*gm(2,3)) 1659 cm(3,2,12)=cm(3,2,12)-(gm(2,2)*(6*gm(1,3)*gm(2,2)+3*gm(1,2)*gm(2,3))) 1660 cm(3,3,12)=cm(3,3,12)-(12*gm(1,2)*gm(2,3)*gm(3,3)+gm(1,3)*(1.5d0*gm(2,3)& 1661 & **2-4.5d0*gm(2,2)*gm(3,3))) 1662 cm(3,4,12)=cm(3,4,12)-(3*gm(1,3)*gm(2,2)*gm(2,3)+gm(1,2)*(3*gm(2,3)& 1663 & **2+12*gm(2,2)*gm(3,3))) 1664 cm(3,5,12)=cm(3,5,12)-(-3*gm(1,3)**2*gm(2,2)+9*gm(1,2)*gm(1,3)& 1665 & *gm(2,3)+7.5d0*gm(1,2)**2*gm(3,3)+gm(1,1)*(-3*gm(2,3)**2-1.5d0*gm(2,2)& 1666 & *gm(3,3))) 1667 cm(3,6,12)=cm(3,6,12)-(12*gm(1,2)*gm(1,3)*gm(2,2)+1.5d0*gm(1,2)& 1668 & **2*gm(2,3)-4.5d0*gm(1,1)*gm(2,2)*gm(2,3)) 1669 cm(3,7,12)=cm(3,7,12)-(3*gm(2,2)**2*gm(2,3)) 1670 cm(3,8,12)=cm(3,8,12)-(1.5d0*gm(2,3)**3+7.5d0*gm(2,2)*gm(2,3)& 1671 & *gm(3,3)) 1672 cm(3,9,12)=cm(3,9,12)-(gm(2,2)*(3*gm(2,3)**2+6*gm(2,2)*gm(3,3))) 1673 cm(3,10,12)=cm(3,10,12)-(gm(3,3)*(4.5d0*gm(2,3)**2-1.5d0*gm(2,2)& 1674 & *gm(3,3))) 1675 cm(3,1,13)=cm(3,1,13)-(2.5d0*gm(1,2)**3-1.5d0*gm(1,1)*gm(1,2)& 1676 & *gm(2,2)) 1677 cm(3,2,13)=cm(3,2,13)-(3*gm(1,2)*gm(2,2)**2) 1678 cm(3,3,13)=cm(3,3,13)-(-3*gm(1,3)*gm(2,2)*gm(2,3)+gm(1,2)*(7.5d0*gm(2,3)& 1679 & **2-1.5d0*gm(2,2)*gm(3,3))) 1680 cm(3,4,13)=cm(3,4,13)-(gm(2,2)*(-3*gm(1,3)*gm(2,2)+9*gm(1,2)*gm(2,3))) 1681 cm(3,5,13)=cm(3,5,13)-(-3*gm(1,2)*gm(1,3)*gm(2,2)+7.5d0*gm(1,2)& 1682 & **2*gm(2,3)-1.5d0*gm(1,1)*gm(2,2)*gm(2,3)) 1683 cm(3,6,13)=cm(3,6,13)-(gm(2,2)*(4.5d0*gm(1,2)**2-1.5d0*gm(1,1)& 1684 & *gm(2,2))) 1685 cm(3,7,13)=cm(3,7,13)-(gm(2,2)**3) 1686 cm(3,8,13)=cm(3,8,13)-(gm(2,2)*(4.5d0*gm(2,3)**2-1.5d0*gm(2,2)& 1687 & *gm(3,3))) 1688 cm(3,9,13)=cm(3,9,13)-(3*gm(2,2)**2*gm(2,3)) 1689 cm(3,10,13)=cm(3,10,13)-(2.5d0*gm(2,3)**3-1.5d0*gm(2,2)*gm(2,3)& 1690 & *gm(3,3)) 1691 cm(3,1,14)=cm(3,1,14)-(-3*gm(1,1)*gm(1,3)*gm(2,3)+gm(1,2)*(7.5d0*gm(1,3)& 1692 & **2-1.5d0*gm(1,1)*gm(3,3))) 1693 cm(3,2,14)=cm(3,2,14)-(12*gm(1,3)*gm(2,2)*gm(2,3)+gm(1,2)*(1.5d0*gm(2,3)& 1694 & **2-4.5d0*gm(2,2)*gm(3,3))) 1695 cm(3,3,14)=cm(3,3,14)-(gm(3,3)*(3*gm(1,3)*gm(2,3)+6*gm(1,2)*gm(3,3))) 1696 cm(3,4,14)=cm(3,4,14)-(3*gm(1,2)*gm(2,3)*gm(3,3)+gm(1,3)*(3*gm(2,3)& 1697 & **2+12*gm(2,2)*gm(3,3))) 1698 cm(3,5,14)=cm(3,5,14)-(1.5d0*gm(1,3)**2*gm(2,3)+12*gm(1,2)*gm(1,3)& 1699 & *gm(3,3)-4.5d0*gm(1,1)*gm(2,3)*gm(3,3)) 1700 cm(3,6,14)=cm(3,6,14)-(7.5d0*gm(1,3)**2*gm(2,2)+9*gm(1,2)*gm(1,3)& 1701 & *gm(2,3)-3*gm(1,2)**2*gm(3,3)+gm(1,1)*(-3*gm(2,3)**2-1.5d0*gm(2,2)& 1702 & *gm(3,3))) 1703 cm(3,7,14)=cm(3,7,14)-(gm(2,2)*(4.5d0*gm(2,3)**2-1.5d0*gm(2,2)& 1704 & *gm(3,3))) 1705 cm(3,8,14)=cm(3,8,14)-(gm(3,3)*(3*gm(2,3)**2+6*gm(2,2)*gm(3,3))) 1706 cm(3,9,14)=cm(3,9,14)-(1.5d0*gm(2,3)**3+7.5d0*gm(2,2)*gm(2,3)& 1707 & *gm(3,3)) 1708 cm(3,10,14)=cm(3,10,14)-(3*gm(2,3)*gm(3,3)**2) 1709 cm(3,1,15)=cm(3,1,15)-(2.5d0*gm(1,3)**3-1.5d0*gm(1,1)*gm(1,3)& 1710 & *gm(3,3)) 1711 cm(3,2,15)=cm(3,2,15)-(-3*gm(1,2)*gm(2,3)*gm(3,3)+gm(1,3)*(7.5d0*gm(2,3)& 1712 & **2-1.5d0*gm(2,2)*gm(3,3))) 1713 cm(3,3,15)=cm(3,3,15)-(3*gm(1,3)*gm(3,3)**2) 1714 cm(3,4,15)=cm(3,4,15)-(gm(3,3)*(9*gm(1,3)*gm(2,3)-3*gm(1,2)*gm(3,3))) 1715 cm(3,5,15)=cm(3,5,15)-(gm(3,3)*(4.5d0*gm(1,3)**2-1.5d0*gm(1,1)& 1716 & *gm(3,3))) 1717 cm(3,6,15)=cm(3,6,15)-(7.5d0*gm(1,3)**2*gm(2,3)-3*gm(1,2)*gm(1,3)& 1718 & *gm(3,3)-1.5d0*gm(1,1)*gm(2,3)*gm(3,3)) 1719 cm(3,7,15)=cm(3,7,15)-(2.5d0*gm(2,3)**3-1.5d0*gm(2,2)*gm(2,3)& 1720 & *gm(3,3)) 1721 cm(3,8,15)=cm(3,8,15)-(3*gm(2,3)*gm(3,3)**2) 1722 cm(3,9,15)=cm(3,9,15)-(gm(3,3)*(4.5d0*gm(2,3)**2-1.5d0*gm(2,2)& 1723 & *gm(3,3))) 1724 cm(3,10,15)=cm(3,10,15)-(gm(3,3)**3) 1725 end if 1726 end if 1727 ! 1728 !contraction to output 3-vector 1729 ! 1730 eisnl(:)=0.d0 1731 do jj=1,((rank+2)*(rank+3))/2 1732 tmp(:,:)=0.d0 1733 do ii=1,((rank+1)*(rank+2))/2 1734 tmp(:,1)=tmp(:,1)+aa(:,ii)*cm(1,ii,jj) 1735 tmp(:,2)=tmp(:,2)+aa(:,ii)*cm(2,ii,jj) 1736 tmp(:,3)=tmp(:,3)+aa(:,ii)*cm(3,ii,jj) 1737 end do 1738 eisnl(:)=eisnl(:)+tmp(1,:)*bb(1,jj)+tmp(2,:)*bb(2,jj) 1739 end do 1740 eisnl(:)=2.0d0*eisnl(:) 1741 !factor of 2 multiplied in to drop call to conjugate contraction 1742 !eisnl(:)=0.5d0*eisnl(:) 1743 1744 ABI_DEALLOCATE(cm) 1745 1746 end subroutine contistr01