TABLE OF CONTENTS
ABINIT/contstr26 [ Functions ]
NAME
contstr26
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
istr1=1,...6 specifies cartesian strain component 11,22,33,32,31,21 istr2=seccond strain component 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
e2nl=contraction for nonlocal 2nd-order 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
48 #if defined HAVE_CONFIG_H 49 #include "config.h" 50 #endif 51 52 #include "abi_common.h" 53 54 subroutine contstr26(istr1,istr2,rank,gm,gprimd,e2nl,aa,bb) 55 56 use m_profiling_abi 57 58 use defs_basis 59 60 !This section has been created automatically by the script Abilint (TD). 61 !Do not modify the following lines by hand. 62 #undef ABI_FUNC 63 #define ABI_FUNC 'contstr26' 64 !End of the abilint section 65 66 implicit none 67 68 !Arguments ------------------------------------ 69 !scalars 70 integer,intent(in) :: istr1,istr2,rank 71 real(dp),intent(out) :: e2nl 72 !arrays 73 real(dp),intent(in) :: aa(2,((rank+1)*(rank+2))/2),bb(2,((rank+1)*(rank+2))/2) 74 real(dp),intent(in) :: gm(3,3),gprimd(3,3) 75 76 !Local variables------------------------------- 77 !scalars 78 integer,parameter :: mrank=3 79 integer :: ii,jj,ka,kb,kd,kg 80 !arrays 81 integer,save :: idx(12)=(/1,1,2,2,3,3,3,2,3,1,2,1/) 82 real(dp) :: d2gm(3,3),dgm01(3,3),dgm10(3,3),tmp(2) 83 real(dp),allocatable :: cm(:,:) 84 85 ! ************************************************************************* 86 ABI_ALLOCATE(cm,(((mrank+1)*(mrank+2))/2,((mrank+1)*(mrank+2))/2)) 87 88 ka=idx(2*istr1-1);kb=idx(2*istr1);kg=idx(2*istr2-1);kd=idx(2*istr2) 89 90 do ii = 1,3 91 dgm01(:,ii)=-(gprimd(ka,:)*gprimd(kb,ii)+gprimd(kb,:)*gprimd(ka,ii)) 92 dgm10(:,ii)=-(gprimd(kg,:)*gprimd(kd,ii)+gprimd(kd,:)*gprimd(kg,ii)) 93 end do 94 95 d2gm(:,:)=0.d0 96 do ii = 1,3 97 if(ka==kg) d2gm(:,ii)=d2gm(:,ii)& 98 & +gprimd(kb,:)*gprimd(kd,ii)+gprimd(kd,:)*gprimd(kb,ii) 99 if(ka==kd) d2gm(:,ii)=d2gm(:,ii)& 100 & +gprimd(kb,:)*gprimd(kg,ii)+gprimd(kg,:)*gprimd(kb,ii) 101 if(kb==kg) d2gm(:,ii)=d2gm(:,ii)& 102 & +gprimd(ka,:)*gprimd(kd,ii)+gprimd(kd,:)*gprimd(ka,ii) 103 if(kb==kd) d2gm(:,ii)=d2gm(:,ii)& 104 & +gprimd(ka,:)*gprimd(kg,ii)+gprimd(kg,:)*gprimd(ka,ii) 105 end do 106 d2gm(:,:)=0.5d0*d2gm(:,:) 107 108 ! 109 !The code below was written by a Mathematica program and formatted by 110 !a combination of editing scripts. It is not intended to be read 111 !by human beings, and certainly not to be modified by one. Conceivably 112 !it could be shortened somewhat by identifying common subexpressions. 113 ! 114 if(rank==0)then 115 cm(1,1)=0.0d0 116 elseif(rank==1)then 117 cm(1,1)=d2gm(1,1) 118 cm(2,1)=d2gm(1,2) 119 cm(3,1)=d2gm(1,3) 120 cm(1,2)=d2gm(1,2) 121 cm(2,2)=d2gm(2,2) 122 cm(3,2)=d2gm(2,3) 123 cm(1,3)=d2gm(1,3) 124 cm(2,3)=d2gm(2,3) 125 cm(3,3)=d2gm(3,3) 126 elseif(rank==2)then 127 cm(1,1)=2*(dgm01(1,1)*dgm10(1,1)+gm(1,1)*d2gm(1,1)) 128 cm(2,1)=-0.5d0*dgm01(2,2)*dgm10(1,1)+3*dgm01(1,2)*dgm10(1,2)& 129 & -0.5d0*dgm01(1,1)*dgm10(2,2)-0.5d0*gm(2,2)*d2gm(1,1)+3*gm(1,2)& 130 & *d2gm(1,2)-0.5d0*gm(1,1)*d2gm(2,2) 131 cm(3,1)=-0.5d0*dgm01(3,3)*dgm10(1,1)+3*dgm01(1,3)*dgm10(1,3)& 132 & -0.5d0*dgm01(1,1)*dgm10(3,3)-0.5d0*gm(3,3)*d2gm(1,1)+3*gm(1,3)& 133 & *d2gm(1,3)-0.5d0*gm(1,1)*d2gm(3,3) 134 cm(4,1)=-dgm01(2,3)*dgm10(1,1)+3*dgm01(1,3)*dgm10(1,2)+3*dgm01(1,2)& 135 & *dgm10(1,3)-dgm01(1,1)*dgm10(2,3)-gm(2,3)*d2gm(1,1)+3*gm(1,3)& 136 & *d2gm(1,2)+3*gm(1,2)*d2gm(1,3)-gm(1,1)*d2gm(2,3) 137 cm(5,1)=2*(dgm01(1,3)*dgm10(1,1)+dgm01(1,1)*dgm10(1,3)+gm(1,3)& 138 & *d2gm(1,1)+gm(1,1)*d2gm(1,3)) 139 cm(6,1)=2*(dgm01(1,2)*dgm10(1,1)+dgm01(1,1)*dgm10(1,2)+gm(1,2)& 140 & *d2gm(1,1)+gm(1,1)*d2gm(1,2)) 141 cm(1,2)=-0.5d0*dgm01(2,2)*dgm10(1,1)+3*dgm01(1,2)*dgm10(1,2)& 142 & -0.5d0*dgm01(1,1)*dgm10(2,2)-0.5d0*gm(2,2)*d2gm(1,1)+3*gm(1,2)& 143 & *d2gm(1,2)-0.5d0*gm(1,1)*d2gm(2,2) 144 cm(2,2)=2*(dgm01(2,2)*dgm10(2,2)+gm(2,2)*d2gm(2,2)) 145 cm(3,2)=-0.5d0*dgm01(3,3)*dgm10(2,2)+3*dgm01(2,3)*dgm10(2,3)& 146 & -0.5d0*dgm01(2,2)*dgm10(3,3)-0.5d0*gm(3,3)*d2gm(2,2)+3*gm(2,3)& 147 & *d2gm(2,3)-0.5d0*gm(2,2)*d2gm(3,3) 148 cm(4,2)=2*(dgm01(2,3)*dgm10(2,2)+dgm01(2,2)*dgm10(2,3)+gm(2,3)& 149 & *d2gm(2,2)+gm(2,2)*d2gm(2,3)) 150 cm(5,2)=3*dgm01(2,3)*dgm10(1,2)-dgm01(2,2)*dgm10(1,3)-dgm01(1,3)& 151 & *dgm10(2,2)+3*dgm01(1,2)*dgm10(2,3)+3*gm(2,3)*d2gm(1,2)-gm(2,2)& 152 & *d2gm(1,3)-gm(1,3)*d2gm(2,2)+3*gm(1,2)*d2gm(2,3) 153 cm(6,2)=2*(dgm01(2,2)*dgm10(1,2)+dgm01(1,2)*dgm10(2,2)+gm(2,2)& 154 & *d2gm(1,2)+gm(1,2)*d2gm(2,2)) 155 cm(1,3)=-0.5d0*dgm01(3,3)*dgm10(1,1)+3*dgm01(1,3)*dgm10(1,3)& 156 & -0.5d0*dgm01(1,1)*dgm10(3,3)-0.5d0*gm(3,3)*d2gm(1,1)+3*gm(1,3)& 157 & *d2gm(1,3)-0.5d0*gm(1,1)*d2gm(3,3) 158 cm(2,3)=-0.5d0*dgm01(3,3)*dgm10(2,2)+3*dgm01(2,3)*dgm10(2,3)& 159 & -0.5d0*dgm01(2,2)*dgm10(3,3)-0.5d0*gm(3,3)*d2gm(2,2)+3*gm(2,3)& 160 & *d2gm(2,3)-0.5d0*gm(2,2)*d2gm(3,3) 161 cm(3,3)=2*(dgm01(3,3)*dgm10(3,3)+gm(3,3)*d2gm(3,3)) 162 cm(4,3)=2*(dgm01(3,3)*dgm10(2,3)+dgm01(2,3)*dgm10(3,3)+gm(3,3)& 163 & *d2gm(2,3)+gm(2,3)*d2gm(3,3)) 164 cm(5,3)=2*(dgm01(3,3)*dgm10(1,3)+dgm01(1,3)*dgm10(3,3)+gm(3,3)& 165 & *d2gm(1,3)+gm(1,3)*d2gm(3,3)) 166 cm(6,3)=-dgm01(3,3)*dgm10(1,2)+3*dgm01(2,3)*dgm10(1,3)+3*dgm01(1,3)& 167 & *dgm10(2,3)-dgm01(1,2)*dgm10(3,3)-gm(3,3)*d2gm(1,2)+3*gm(2,3)& 168 & *d2gm(1,3)+3*gm(1,3)*d2gm(2,3)-gm(1,2)*d2gm(3,3) 169 cm(1,4)=-dgm01(2,3)*dgm10(1,1)+3*dgm01(1,3)*dgm10(1,2)+3*dgm01(1,2)& 170 & *dgm10(1,3)-dgm01(1,1)*dgm10(2,3)-gm(2,3)*d2gm(1,1)+3*gm(1,3)& 171 & *d2gm(1,2)+3*gm(1,2)*d2gm(1,3)-gm(1,1)*d2gm(2,3) 172 cm(2,4)=2*(dgm01(2,3)*dgm10(2,2)+dgm01(2,2)*dgm10(2,3)+gm(2,3)& 173 & *d2gm(2,2)+gm(2,2)*d2gm(2,3)) 174 cm(3,4)=2*(dgm01(3,3)*dgm10(2,3)+dgm01(2,3)*dgm10(3,3)+gm(3,3)& 175 & *d2gm(2,3)+gm(2,3)*d2gm(3,3)) 176 cm(4,4)=3*dgm01(3,3)*dgm10(2,2)+2*dgm01(2,3)*dgm10(2,3)+3*dgm01(2,2)& 177 & *dgm10(3,3)+3*gm(3,3)*d2gm(2,2)+2*gm(2,3)*d2gm(2,3)+3*gm(2,2)& 178 & *d2gm(3,3) 179 cm(5,4)=3*dgm01(3,3)*dgm10(1,2)+dgm01(2,3)*dgm10(1,3)+dgm01(1,3)& 180 & *dgm10(2,3)+3*dgm01(1,2)*dgm10(3,3)+3*gm(3,3)*d2gm(1,2)+gm(2,3)& 181 & *d2gm(1,3)+gm(1,3)*d2gm(2,3)+3*gm(1,2)*d2gm(3,3) 182 cm(6,4)=dgm01(2,3)*dgm10(1,2)+3*dgm01(2,2)*dgm10(1,3)+3*dgm01(1,3)& 183 & *dgm10(2,2)+dgm01(1,2)*dgm10(2,3)+gm(2,3)*d2gm(1,2)+3*gm(2,2)& 184 & *d2gm(1,3)+3*gm(1,3)*d2gm(2,2)+gm(1,2)*d2gm(2,3) 185 cm(1,5)=2*(dgm01(1,3)*dgm10(1,1)+dgm01(1,1)*dgm10(1,3)+gm(1,3)& 186 & *d2gm(1,1)+gm(1,1)*d2gm(1,3)) 187 cm(2,5)=3*dgm01(2,3)*dgm10(1,2)-dgm01(2,2)*dgm10(1,3)-dgm01(1,3)& 188 & *dgm10(2,2)+3*dgm01(1,2)*dgm10(2,3)+3*gm(2,3)*d2gm(1,2)-gm(2,2)& 189 & *d2gm(1,3)-gm(1,3)*d2gm(2,2)+3*gm(1,2)*d2gm(2,3) 190 cm(3,5)=2*(dgm01(3,3)*dgm10(1,3)+dgm01(1,3)*dgm10(3,3)+gm(3,3)& 191 & *d2gm(1,3)+gm(1,3)*d2gm(3,3)) 192 cm(4,5)=3*dgm01(3,3)*dgm10(1,2)+dgm01(2,3)*dgm10(1,3)+dgm01(1,3)& 193 & *dgm10(2,3)+3*dgm01(1,2)*dgm10(3,3)+3*gm(3,3)*d2gm(1,2)+gm(2,3)& 194 & *d2gm(1,3)+gm(1,3)*d2gm(2,3)+3*gm(1,2)*d2gm(3,3) 195 cm(5,5)=3*dgm01(3,3)*dgm10(1,1)+2*dgm01(1,3)*dgm10(1,3)+3*dgm01(1,1)& 196 & *dgm10(3,3)+3*gm(3,3)*d2gm(1,1)+2*gm(1,3)*d2gm(1,3)+3*gm(1,1)& 197 & *d2gm(3,3) 198 cm(6,5)=3*dgm01(2,3)*dgm10(1,1)+dgm01(1,3)*dgm10(1,2)+dgm01(1,2)& 199 & *dgm10(1,3)+3*dgm01(1,1)*dgm10(2,3)+3*gm(2,3)*d2gm(1,1)+gm(1,3)& 200 & *d2gm(1,2)+gm(1,2)*d2gm(1,3)+3*gm(1,1)*d2gm(2,3) 201 cm(1,6)=2*(dgm01(1,2)*dgm10(1,1)+dgm01(1,1)*dgm10(1,2)+gm(1,2)& 202 & *d2gm(1,1)+gm(1,1)*d2gm(1,2)) 203 cm(2,6)=2*(dgm01(2,2)*dgm10(1,2)+dgm01(1,2)*dgm10(2,2)+gm(2,2)& 204 & *d2gm(1,2)+gm(1,2)*d2gm(2,2)) 205 cm(3,6)=-dgm01(3,3)*dgm10(1,2)+3*dgm01(2,3)*dgm10(1,3)+3*dgm01(1,3)& 206 & *dgm10(2,3)-dgm01(1,2)*dgm10(3,3)-gm(3,3)*d2gm(1,2)+3*gm(2,3)& 207 & *d2gm(1,3)+3*gm(1,3)*d2gm(2,3)-gm(1,2)*d2gm(3,3) 208 cm(4,6)=dgm01(2,3)*dgm10(1,2)+3*dgm01(2,2)*dgm10(1,3)+3*dgm01(1,3)& 209 & *dgm10(2,2)+dgm01(1,2)*dgm10(2,3)+gm(2,3)*d2gm(1,2)+3*gm(2,2)& 210 & *d2gm(1,3)+3*gm(1,3)*d2gm(2,2)+gm(1,2)*d2gm(2,3) 211 cm(5,6)=3*dgm01(2,3)*dgm10(1,1)+dgm01(1,3)*dgm10(1,2)+dgm01(1,2)& 212 & *dgm10(1,3)+3*dgm01(1,1)*dgm10(2,3)+3*gm(2,3)*d2gm(1,1)+gm(1,3)& 213 & *d2gm(1,2)+gm(1,2)*d2gm(1,3)+3*gm(1,1)*d2gm(2,3) 214 cm(6,6)=3*dgm01(2,2)*dgm10(1,1)+2*dgm01(1,2)*dgm10(1,2)+3*dgm01(1,1)& 215 & *dgm10(2,2)+3*gm(2,2)*d2gm(1,1)+2*gm(1,2)*d2gm(1,2)+3*gm(1,1)& 216 & *d2gm(2,2) 217 elseif(rank==3)then 218 cm(1,1)=gm(1,1)*(6*dgm01(1,1)*dgm10(1,1)+3*gm(1,1)*d2gm(1,1)) 219 cm(2,1)=4.5d0*gm(1,2)**2*d2gm(1,1)-3*gm(2,2)*(dgm01(1,1)*dgm10(1,1)& 220 & +gm(1,1)*d2gm(1,1))+9*gm(1,2)*(dgm01(1,2)*dgm10(1,1)+dgm01(1,1)& 221 & *dgm10(1,2)+gm(1,1)*d2gm(1,2))+gm(1,1)*(-3*dgm01(2,2)*dgm10(1,1)& 222 & +9*dgm01(1,2)*dgm10(1,2)-3*dgm01(1,1)*dgm10(2,2)-1.5d0*gm(1,1)& 223 & *d2gm(2,2)) 224 cm(3,1)=4.5d0*gm(1,3)**2*d2gm(1,1)-3*gm(3,3)*(dgm01(1,1)*dgm10(1,1)& 225 & +gm(1,1)*d2gm(1,1))+9*gm(1,3)*(dgm01(1,3)*dgm10(1,1)+dgm01(1,1)& 226 & *dgm10(1,3)+gm(1,1)*d2gm(1,3))+gm(1,1)*(-3*dgm01(3,3)*dgm10(1,1)& 227 & +9*dgm01(1,3)*dgm10(1,3)-3*dgm01(1,1)*dgm10(3,3)-1.5d0*gm(1,1)& 228 & *d2gm(3,3)) 229 cm(4,1)=9*gm(1,2)*dgm01(1,3)*dgm10(1,1)-6*gm(1,1)*dgm01(2,3)& 230 & *dgm10(1,1)+9*gm(1,1)*dgm01(1,3)*dgm10(1,2)+9*gm(1,2)*dgm01(1,1)& 231 & *dgm10(1,3)+9*gm(1,1)*dgm01(1,2)*dgm10(1,3)-6*gm(1,1)*dgm01(1,1)& 232 & *dgm10(2,3)-6*gm(2,3)*(dgm01(1,1)*dgm10(1,1)+gm(1,1)*d2gm(1,1))& 233 & +9*gm(1,3)*(dgm01(1,2)*dgm10(1,1)+dgm01(1,1)*dgm10(1,2)+gm(1,2)& 234 & *d2gm(1,1)+gm(1,1)*d2gm(1,2))+9*gm(1,1)*gm(1,2)*d2gm(1,3)-3*gm(1,1)& 235 & **2*d2gm(2,3) 236 cm(5,1)=6*gm(1,3)*(dgm01(1,1)*dgm10(1,1)+gm(1,1)*d2gm(1,1))& 237 & +gm(1,1)*(6*dgm01(1,3)*dgm10(1,1)+6*dgm01(1,1)*dgm10(1,3)+3*gm(1,1)& 238 & *d2gm(1,3)) 239 cm(6,1)=6*gm(1,2)*(dgm01(1,1)*dgm10(1,1)+gm(1,1)*d2gm(1,1))& 240 & +gm(1,1)*(6*dgm01(1,2)*dgm10(1,1)+6*dgm01(1,1)*dgm10(1,2)+3*gm(1,1)& 241 & *d2gm(1,2)) 242 cm(7,1)=-1.5d0*gm(1,1)*(dgm01(2,2)*dgm10(1,2)+dgm01(1,2)*dgm10(2,2))& 243 & +7.5d0*gm(1,2)**2*d2gm(1,2)-1.5d0*gm(2,2)*(dgm01(1,2)*dgm10(1,1)& 244 & +dgm01(1,1)*dgm10(1,2)+gm(1,2)*d2gm(1,1)+gm(1,1)*d2gm(1,2))+gm(1,2)& 245 & *(-1.5d0*dgm01(2,2)*dgm10(1,1)+15*dgm01(1,2)*dgm10(1,2)-1.5d0*dgm01(1,1)& 246 & *dgm10(2,2)-1.5d0*gm(1,1)*d2gm(2,2)) 247 cm(8,1)=-3*gm(1,3)*dgm01(2,3)*dgm10(1,1)-1.5d0*gm(1,2)*dgm01(3,3)& 248 & *dgm10(1,1)+15*gm(1,3)*dgm01(1,3)*dgm10(1,2)-1.5d0*gm(1,1)*dgm01(3,3)& 249 & *dgm10(1,2)+15*gm(1,3)*dgm01(1,2)*dgm10(1,3)+15*gm(1,2)*dgm01(1,3)& 250 & *dgm10(1,3)-3*gm(1,1)*dgm01(2,3)*dgm10(1,3)-3*gm(1,3)*dgm01(1,1)& 251 & *dgm10(2,3)-3*gm(1,1)*dgm01(1,3)*dgm10(2,3)-1.5d0*gm(1,2)*dgm01(1,1)& 252 & *dgm10(3,3)-1.5d0*gm(1,1)*dgm01(1,2)*dgm10(3,3)+7.5d0*gm(1,3)& 253 & **2*d2gm(1,2)-1.5d0*gm(3,3)*(dgm01(1,2)*dgm10(1,1)+dgm01(1,1)& 254 & *dgm10(1,2)+gm(1,2)*d2gm(1,1)+gm(1,1)*d2gm(1,2))+15*gm(1,2)*gm(1,3)& 255 & *d2gm(1,3)-3*gm(2,3)*(dgm01(1,3)*dgm10(1,1)+dgm01(1,1)*dgm10(1,3)& 256 & +gm(1,3)*d2gm(1,1)+gm(1,1)*d2gm(1,3))-3*gm(1,1)*gm(1,3)*d2gm(2,3)& 257 & -1.5d0*gm(1,1)*gm(1,2)*d2gm(3,3) 258 cm(9,1)=-1.5d0*gm(1,3)*dgm01(2,2)*dgm10(1,1)-3*gm(1,2)*dgm01(2,3)& 259 & *dgm10(1,1)+15*gm(1,3)*dgm01(1,2)*dgm10(1,2)+15*gm(1,2)*dgm01(1,3)& 260 & *dgm10(1,2)-3*gm(1,1)*dgm01(2,3)*dgm10(1,2)+15*gm(1,2)*dgm01(1,2)& 261 & *dgm10(1,3)-1.5d0*gm(1,1)*dgm01(2,2)*dgm10(1,3)-1.5d0*gm(1,3)& 262 & *dgm01(1,1)*dgm10(2,2)-1.5d0*gm(1,1)*dgm01(1,3)*dgm10(2,2)-3*gm(1,2)& 263 & *dgm01(1,1)*dgm10(2,3)-3*gm(1,1)*dgm01(1,2)*dgm10(2,3)+15*gm(1,2)& 264 & *gm(1,3)*d2gm(1,2)-3*gm(2,3)*(dgm01(1,2)*dgm10(1,1)+dgm01(1,1)& 265 & *dgm10(1,2)+gm(1,2)*d2gm(1,1)+gm(1,1)*d2gm(1,2))+7.5d0*gm(1,2)& 266 & **2*d2gm(1,3)-1.5d0*gm(2,2)*(dgm01(1,3)*dgm10(1,1)+dgm01(1,1)& 267 & *dgm10(1,3)+gm(1,3)*d2gm(1,1)+gm(1,1)*d2gm(1,3))-1.5d0*gm(1,1)& 268 & *gm(1,3)*d2gm(2,2)-3*gm(1,1)*gm(1,2)*d2gm(2,3) 269 cm(10,1)=-1.5d0*gm(1,1)*(dgm01(3,3)*dgm10(1,3)+dgm01(1,3)*dgm10(3,3))& 270 & +7.5d0*gm(1,3)**2*d2gm(1,3)-1.5d0*gm(3,3)*(dgm01(1,3)*dgm10(1,1)& 271 & +dgm01(1,1)*dgm10(1,3)+gm(1,3)*d2gm(1,1)+gm(1,1)*d2gm(1,3))+gm(1,3)& 272 & *(-1.5d0*dgm01(3,3)*dgm10(1,1)+15*dgm01(1,3)*dgm10(1,3)-1.5d0*dgm01(1,1)& 273 & *dgm10(3,3)-1.5d0*gm(1,1)*d2gm(3,3)) 274 cm(1,2)=4.5d0*gm(1,2)**2*d2gm(1,1)-3*gm(2,2)*(dgm01(1,1)*dgm10(1,1)& 275 & +gm(1,1)*d2gm(1,1))+9*gm(1,2)*(dgm01(1,2)*dgm10(1,1)+dgm01(1,1)& 276 & *dgm10(1,2)+gm(1,1)*d2gm(1,2))+gm(1,1)*(-3*dgm01(2,2)*dgm10(1,1)& 277 & +9*dgm01(1,2)*dgm10(1,2)-3*dgm01(1,1)*dgm10(2,2)-1.5d0*gm(1,1)& 278 & *d2gm(2,2)) 279 cm(2,2)=12*gm(1,1)*dgm01(2,2)*dgm10(2,2)+6*gm(1,2)*(dgm01(2,2)& 280 & *dgm10(1,2)+dgm01(1,2)*dgm10(2,2))+6*gm(2,2)**2*d2gm(1,1)+3*gm(1,2)& 281 & **2*d2gm(2,2)+gm(2,2)*(12*dgm01(2,2)*dgm10(1,1)+6*dgm01(1,2)& 282 & *dgm10(1,2)+12*dgm01(1,1)*dgm10(2,2)+6*gm(1,2)*d2gm(1,2)+12*gm(1,1)& 283 & *d2gm(2,2)) 284 cm(3,2)=-1.5d0*gm(2,2)*dgm01(3,3)*dgm10(1,1)+9*gm(1,3)*dgm01(2,3)& 285 & *dgm10(1,2)-6*gm(1,2)*dgm01(3,3)*dgm10(1,2)-6*gm(2,2)*dgm01(1,3)& 286 & *dgm10(1,3)-6*gm(1,3)*dgm01(2,2)*dgm10(1,3)+9*gm(1,2)*dgm01(2,3)& 287 & *dgm10(1,3)-6*gm(1,3)*dgm01(1,3)*dgm10(2,2)-1.5d0*gm(1,1)*dgm01(3,3)& 288 & *dgm10(2,2)+9*gm(1,3)*dgm01(1,2)*dgm10(2,3)+9*gm(1,2)*dgm01(1,3)& 289 & *dgm10(2,3)+15*gm(1,1)*dgm01(2,3)*dgm10(2,3)-1.5d0*gm(2,2)*dgm01(1,1)& 290 & *dgm10(3,3)-6*gm(1,2)*dgm01(1,2)*dgm10(3,3)-1.5d0*gm(1,1)*dgm01(2,2)& 291 & *dgm10(3,3)+7.5d0*gm(2,3)**2*d2gm(1,1)-6*gm(1,3)*gm(2,2)*d2gm(1,3)& 292 & -3*gm(1,3)**2*d2gm(2,2)+gm(3,3)*(-1.5d0*dgm01(2,2)*dgm10(1,1)& 293 & -6*dgm01(1,2)*dgm10(1,2)-1.5d0*dgm01(1,1)*dgm10(2,2)-1.5d0*gm(2,2)& 294 & *d2gm(1,1)-6*gm(1,2)*d2gm(1,2)-1.5d0*gm(1,1)*d2gm(2,2))+9*gm(1,2)& 295 & *gm(1,3)*d2gm(2,3)+gm(2,3)*(15*dgm01(2,3)*dgm10(1,1)+9*dgm01(1,3)& 296 & *dgm10(1,2)+9*dgm01(1,2)*dgm10(1,3)+15*dgm01(1,1)*dgm10(2,3)& 297 & +9*gm(1,3)*d2gm(1,2)+9*gm(1,2)*d2gm(1,3)+15*gm(1,1)*d2gm(2,3))& 298 & -3*gm(1,2)**2*d2gm(3,3)-1.5d0*gm(1,1)*gm(2,2)*d2gm(3,3) 299 cm(4,2)=3*gm(1,3)*dgm01(2,2)*dgm10(1,2)+6*gm(1,2)*dgm01(2,3)& 300 & *dgm10(1,2)+3*gm(1,2)*dgm01(2,2)*dgm10(1,3)+3*gm(1,3)*dgm01(1,2)& 301 & *dgm10(2,2)+3*gm(1,2)*dgm01(1,3)*dgm10(2,2)+12*gm(1,1)*dgm01(2,3)& 302 & *dgm10(2,2)+6*gm(1,2)*dgm01(1,2)*dgm10(2,3)+12*gm(1,1)*dgm01(2,2)& 303 & *dgm10(2,3)+3*gm(1,2)*gm(1,3)*d2gm(2,2)+gm(2,3)*(12*dgm01(2,2)& 304 & *dgm10(1,1)+6*dgm01(1,2)*dgm10(1,2)+12*dgm01(1,1)*dgm10(2,2)& 305 & +12*gm(2,2)*d2gm(1,1)+6*gm(1,2)*d2gm(1,2)+12*gm(1,1)*d2gm(2,2))& 306 & +3*gm(1,2)**2*d2gm(2,3)+gm(2,2)*(12*dgm01(2,3)*dgm10(1,1)+3*dgm01(1,3)& 307 & *dgm10(1,2)+3*dgm01(1,2)*dgm10(1,3)+12*dgm01(1,1)*dgm10(2,3)& 308 & +3*gm(1,3)*d2gm(1,2)+3*gm(1,2)*d2gm(1,3)+12*gm(1,1)*d2gm(2,3)) 309 cm(5,2)=-4.5d0*gm(1,3)*dgm01(2,2)*dgm10(1,1)+12*gm(1,2)*dgm01(2,3)& 310 & *dgm10(1,1)+3*gm(1,3)*dgm01(1,2)*dgm10(1,2)+3*gm(1,2)*dgm01(1,3)& 311 & *dgm10(1,2)+12*gm(1,1)*dgm01(2,3)*dgm10(1,2)+3*gm(1,2)*dgm01(1,2)& 312 & *dgm10(1,3)-4.5d0*gm(1,1)*dgm01(2,2)*dgm10(1,3)-4.5d0*gm(1,3)& 313 & *dgm01(1,1)*dgm10(2,2)-4.5d0*gm(1,1)*dgm01(1,3)*dgm10(2,2)+12*gm(1,2)& 314 & *dgm01(1,1)*dgm10(2,3)+12*gm(1,1)*dgm01(1,2)*dgm10(2,3)+3*gm(1,2)& 315 & *gm(1,3)*d2gm(1,2)+12*gm(2,3)*(dgm01(1,2)*dgm10(1,1)+dgm01(1,1)& 316 & *dgm10(1,2)+gm(1,2)*d2gm(1,1)+gm(1,1)*d2gm(1,2))+1.5d0*gm(1,2)& 317 & **2*d2gm(1,3)-4.5d0*gm(2,2)*(dgm01(1,3)*dgm10(1,1)+dgm01(1,1)& 318 & *dgm10(1,3)+gm(1,3)*d2gm(1,1)+gm(1,1)*d2gm(1,3))-4.5d0*gm(1,1)& 319 & *gm(1,3)*d2gm(2,2)+12*gm(1,1)*gm(1,2)*d2gm(2,3) 320 cm(6,2)=7.5d0*gm(1,1)*(dgm01(2,2)*dgm10(1,2)+dgm01(1,2)*dgm10(2,2))& 321 & +4.5d0*gm(1,2)**2*d2gm(1,2)+7.5d0*gm(2,2)*(dgm01(1,2)*dgm10(1,1)& 322 & +dgm01(1,1)*dgm10(1,2)+gm(1,2)*d2gm(1,1)+gm(1,1)*d2gm(1,2))+gm(1,2)& 323 & *(7.5d0*dgm01(2,2)*dgm10(1,1)+9*dgm01(1,2)*dgm10(1,2)+7.5d0*dgm01(1,1)& 324 & *dgm10(2,2)+7.5d0*gm(1,1)*d2gm(2,2)) 325 cm(7,2)=6*gm(1,2)*dgm01(2,2)*dgm10(2,2)+3*gm(2,2)**2*d2gm(1,2)& 326 & +6*gm(2,2)*(dgm01(2,2)*dgm10(1,2)+dgm01(1,2)*dgm10(2,2)+gm(1,2)& 327 & *d2gm(2,2)) 328 cm(8,2)=-4.5d0*gm(2,2)*dgm01(3,3)*dgm10(1,2)+12*gm(2,2)*dgm01(2,3)& 329 & *dgm10(1,3)+12*gm(1,3)*dgm01(2,3)*dgm10(2,2)-4.5d0*gm(1,2)*dgm01(3,3)& 330 & *dgm10(2,2)+12*gm(2,2)*dgm01(1,3)*dgm10(2,3)+12*gm(1,3)*dgm01(2,2)& 331 & *dgm10(2,3)+3*gm(1,2)*dgm01(2,3)*dgm10(2,3)-4.5d0*gm(2,2)*dgm01(1,2)& 332 & *dgm10(3,3)-4.5d0*gm(1,2)*dgm01(2,2)*dgm10(3,3)+1.5d0*gm(2,3)& 333 & **2*d2gm(1,2)-4.5d0*gm(3,3)*(dgm01(2,2)*dgm10(1,2)+dgm01(1,2)& 334 & *dgm10(2,2)+gm(2,2)*d2gm(1,2)+gm(1,2)*d2gm(2,2))+12*gm(1,3)*gm(2,2)& 335 & *d2gm(2,3)+gm(2,3)*(3*dgm01(2,3)*dgm10(1,2)+12*dgm01(2,2)*dgm10(1,3)& 336 & +12*dgm01(1,3)*dgm10(2,2)+3*dgm01(1,2)*dgm10(2,3)+12*gm(2,2)& 337 & *d2gm(1,3)+12*gm(1,3)*d2gm(2,2)+3*gm(1,2)*d2gm(2,3))-4.5d0*gm(1,2)& 338 & *gm(2,2)*d2gm(3,3) 339 cm(9,2)=12*gm(1,3)*dgm01(2,2)*dgm10(2,2)+3*gm(1,2)*dgm01(2,3)& 340 & *dgm10(2,2)+3*gm(1,2)*dgm01(2,2)*dgm10(2,3)+6*gm(2,2)**2*d2gm(1,3)& 341 & +3*gm(2,3)*(dgm01(2,2)*dgm10(1,2)+dgm01(1,2)*dgm10(2,2)+gm(2,2)& 342 & *d2gm(1,2)+gm(1,2)*d2gm(2,2))+gm(2,2)*(3*dgm01(2,3)*dgm10(1,2)& 343 & +12*dgm01(2,2)*dgm10(1,3)+12*dgm01(1,3)*dgm10(2,2)+3*dgm01(1,2)& 344 & *dgm10(2,3)+12*gm(1,3)*d2gm(2,2)+3*gm(1,2)*d2gm(2,3)) 345 cm(10,2)=-1.5d0*gm(2,2)*dgm01(3,3)*dgm10(1,3)-1.5d0*gm(1,3)& 346 & *dgm01(3,3)*dgm10(2,2)+15*gm(1,3)*dgm01(2,3)*dgm10(2,3)-3*gm(1,2)& 347 & *dgm01(3,3)*dgm10(2,3)-1.5d0*gm(2,2)*dgm01(1,3)*dgm10(3,3)-1.5d0*gm(1,3)& 348 & *dgm01(2,2)*dgm10(3,3)-3*gm(1,2)*dgm01(2,3)*dgm10(3,3)+7.5d0*gm(2,3)& 349 & **2*d2gm(1,3)+gm(3,3)*(-3*dgm01(2,3)*dgm10(1,2)-1.5d0*dgm01(2,2)& 350 & *dgm10(1,3)-1.5d0*dgm01(1,3)*dgm10(2,2)-3*dgm01(1,2)*dgm10(2,3)& 351 & -3*gm(2,3)*d2gm(1,2)-1.5d0*gm(2,2)*d2gm(1,3)-1.5d0*gm(1,3)*d2gm(2,2)& 352 & -3*gm(1,2)*d2gm(2,3))-1.5d0*gm(1,3)*gm(2,2)*d2gm(3,3)+gm(2,3)& 353 & *(-3*dgm01(3,3)*dgm10(1,2)+15*dgm01(2,3)*dgm10(1,3)+15*dgm01(1,3)& 354 & *dgm10(2,3)-3*dgm01(1,2)*dgm10(3,3)+15*gm(1,3)*d2gm(2,3)-3*gm(1,2)& 355 & *d2gm(3,3)) 356 cm(1,3)=4.5d0*gm(1,3)**2*d2gm(1,1)-3*gm(3,3)*(dgm01(1,1)*dgm10(1,1)& 357 & +gm(1,1)*d2gm(1,1))+9*gm(1,3)*(dgm01(1,3)*dgm10(1,1)+dgm01(1,1)& 358 & *dgm10(1,3)+gm(1,1)*d2gm(1,3))+gm(1,1)*(-3*dgm01(3,3)*dgm10(1,1)& 359 & +9*dgm01(1,3)*dgm10(1,3)-3*dgm01(1,1)*dgm10(3,3)-1.5d0*gm(1,1)& 360 & *d2gm(3,3)) 361 cm(2,3)=-1.5d0*gm(2,2)*dgm01(3,3)*dgm10(1,1)+9*gm(1,3)*dgm01(2,3)& 362 & *dgm10(1,2)-6*gm(1,2)*dgm01(3,3)*dgm10(1,2)-6*gm(2,2)*dgm01(1,3)& 363 & *dgm10(1,3)-6*gm(1,3)*dgm01(2,2)*dgm10(1,3)+9*gm(1,2)*dgm01(2,3)& 364 & *dgm10(1,3)-6*gm(1,3)*dgm01(1,3)*dgm10(2,2)-1.5d0*gm(1,1)*dgm01(3,3)& 365 & *dgm10(2,2)+9*gm(1,3)*dgm01(1,2)*dgm10(2,3)+9*gm(1,2)*dgm01(1,3)& 366 & *dgm10(2,3)+15*gm(1,1)*dgm01(2,3)*dgm10(2,3)-1.5d0*gm(2,2)*dgm01(1,1)& 367 & *dgm10(3,3)-6*gm(1,2)*dgm01(1,2)*dgm10(3,3)-1.5d0*gm(1,1)*dgm01(2,2)& 368 & *dgm10(3,3)+7.5d0*gm(2,3)**2*d2gm(1,1)-6*gm(1,3)*gm(2,2)*d2gm(1,3)& 369 & -3*gm(1,3)**2*d2gm(2,2)+gm(3,3)*(-1.5d0*dgm01(2,2)*dgm10(1,1)& 370 & -6*dgm01(1,2)*dgm10(1,2)-1.5d0*dgm01(1,1)*dgm10(2,2)-1.5d0*gm(2,2)& 371 & *d2gm(1,1)-6*gm(1,2)*d2gm(1,2)-1.5d0*gm(1,1)*d2gm(2,2))+9*gm(1,2)& 372 & *gm(1,3)*d2gm(2,3)+gm(2,3)*(15*dgm01(2,3)*dgm10(1,1)+9*dgm01(1,3)& 373 & *dgm10(1,2)+9*dgm01(1,2)*dgm10(1,3)+15*dgm01(1,1)*dgm10(2,3)& 374 & +9*gm(1,3)*d2gm(1,2)+9*gm(1,2)*d2gm(1,3)+15*gm(1,1)*d2gm(2,3))& 375 & -3*gm(1,2)**2*d2gm(3,3)-1.5d0*gm(1,1)*gm(2,2)*d2gm(3,3) 376 cm(3,3)=12*gm(1,1)*dgm01(3,3)*dgm10(3,3)+6*gm(1,3)*(dgm01(3,3)& 377 & *dgm10(1,3)+dgm01(1,3)*dgm10(3,3))+6*gm(3,3)**2*d2gm(1,1)+3*gm(1,3)& 378 & **2*d2gm(3,3)+gm(3,3)*(12*dgm01(3,3)*dgm10(1,1)+6*dgm01(1,3)& 379 & *dgm10(1,3)+12*dgm01(1,1)*dgm10(3,3)+6*gm(1,3)*d2gm(1,3)+12*gm(1,1)& 380 & *d2gm(3,3)) 381 cm(4,3)=3*gm(1,3)*dgm01(3,3)*dgm10(1,2)+6*gm(1,3)*dgm01(2,3)& 382 & *dgm10(1,3)+3*gm(1,2)*dgm01(3,3)*dgm10(1,3)+6*gm(1,3)*dgm01(1,3)& 383 & *dgm10(2,3)+12*gm(1,1)*dgm01(3,3)*dgm10(2,3)+3*gm(1,3)*dgm01(1,2)& 384 & *dgm10(3,3)+3*gm(1,2)*dgm01(1,3)*dgm10(3,3)+12*gm(1,1)*dgm01(2,3)& 385 & *dgm10(3,3)+3*gm(1,3)**2*d2gm(2,3)+gm(3,3)*(12*dgm01(2,3)*dgm10(1,1)& 386 & +3*dgm01(1,3)*dgm10(1,2)+3*dgm01(1,2)*dgm10(1,3)+12*dgm01(1,1)& 387 & *dgm10(2,3)+12*gm(2,3)*d2gm(1,1)+3*gm(1,3)*d2gm(1,2)+3*gm(1,2)& 388 & *d2gm(1,3)+12*gm(1,1)*d2gm(2,3))+3*gm(1,2)*gm(1,3)*d2gm(3,3)& 389 & +gm(2,3)*(12*dgm01(3,3)*dgm10(1,1)+6*dgm01(1,3)*dgm10(1,3)+12*dgm01(1,1)& 390 & *dgm10(3,3)+6*gm(1,3)*d2gm(1,3)+12*gm(1,1)*d2gm(3,3)) 391 cm(5,3)=7.5d0*gm(1,1)*(dgm01(3,3)*dgm10(1,3)+dgm01(1,3)*dgm10(3,3))& 392 & +4.5d0*gm(1,3)**2*d2gm(1,3)+7.5d0*gm(3,3)*(dgm01(1,3)*dgm10(1,1)& 393 & +dgm01(1,1)*dgm10(1,3)+gm(1,3)*d2gm(1,1)+gm(1,1)*d2gm(1,3))+gm(1,3)& 394 & *(7.5d0*dgm01(3,3)*dgm10(1,1)+9*dgm01(1,3)*dgm10(1,3)+7.5d0*dgm01(1,1)& 395 & *dgm10(3,3)+7.5d0*gm(1,1)*d2gm(3,3)) 396 cm(6,3)=12*gm(1,3)*dgm01(2,3)*dgm10(1,1)-4.5d0*gm(1,2)*dgm01(3,3)& 397 & *dgm10(1,1)+3*gm(1,3)*dgm01(1,3)*dgm10(1,2)-4.5d0*gm(1,1)*dgm01(3,3)& 398 & *dgm10(1,2)+3*gm(1,3)*dgm01(1,2)*dgm10(1,3)+3*gm(1,2)*dgm01(1,3)& 399 & *dgm10(1,3)+12*gm(1,1)*dgm01(2,3)*dgm10(1,3)+12*gm(1,3)*dgm01(1,1)& 400 & *dgm10(2,3)+12*gm(1,1)*dgm01(1,3)*dgm10(2,3)-4.5d0*gm(1,2)*dgm01(1,1)& 401 & *dgm10(3,3)-4.5d0*gm(1,1)*dgm01(1,2)*dgm10(3,3)+1.5d0*gm(1,3)& 402 & **2*d2gm(1,2)-4.5d0*gm(3,3)*(dgm01(1,2)*dgm10(1,1)+dgm01(1,1)& 403 & *dgm10(1,2)+gm(1,2)*d2gm(1,1)+gm(1,1)*d2gm(1,2))+3*gm(1,2)*gm(1,3)& 404 & *d2gm(1,3)+12*gm(2,3)*(dgm01(1,3)*dgm10(1,1)+dgm01(1,1)*dgm10(1,3)& 405 & +gm(1,3)*d2gm(1,1)+gm(1,1)*d2gm(1,3))+12*gm(1,1)*gm(1,3)*d2gm(2,3)& 406 & -4.5d0*gm(1,1)*gm(1,2)*d2gm(3,3) 407 cm(7,3)=-1.5d0*gm(2,2)*dgm01(3,3)*dgm10(1,2)-3*gm(2,2)*dgm01(2,3)& 408 & *dgm10(1,3)-3*gm(1,3)*dgm01(2,3)*dgm10(2,2)-1.5d0*gm(1,2)*dgm01(3,3)& 409 & *dgm10(2,2)-3*gm(2,2)*dgm01(1,3)*dgm10(2,3)-3*gm(1,3)*dgm01(2,2)& 410 & *dgm10(2,3)+15*gm(1,2)*dgm01(2,3)*dgm10(2,3)-1.5d0*gm(2,2)*dgm01(1,2)& 411 & *dgm10(3,3)-1.5d0*gm(1,2)*dgm01(2,2)*dgm10(3,3)+7.5d0*gm(2,3)& 412 & **2*d2gm(1,2)-1.5d0*gm(3,3)*(dgm01(2,2)*dgm10(1,2)+dgm01(1,2)& 413 & *dgm10(2,2)+gm(2,2)*d2gm(1,2)+gm(1,2)*d2gm(2,2))-3*gm(1,3)*gm(2,2)& 414 & *d2gm(2,3)+gm(2,3)*(15*dgm01(2,3)*dgm10(1,2)-3*dgm01(2,2)*dgm10(1,3)& 415 & -3*dgm01(1,3)*dgm10(2,2)+15*dgm01(1,2)*dgm10(2,3)-3*gm(2,2)*d2gm(1,3)& 416 & -3*gm(1,3)*d2gm(2,2)+15*gm(1,2)*d2gm(2,3))-1.5d0*gm(1,2)*gm(2,2)& 417 & *d2gm(3,3) 418 cm(8,3)=3*gm(1,3)*dgm01(3,3)*dgm10(2,3)+3*gm(1,3)*dgm01(2,3)& 419 & *dgm10(3,3)+12*gm(1,2)*dgm01(3,3)*dgm10(3,3)+6*gm(3,3)**2*d2gm(1,2)& 420 & +gm(3,3)*(12*dgm01(3,3)*dgm10(1,2)+3*dgm01(2,3)*dgm10(1,3)+3*dgm01(1,3)& 421 & *dgm10(2,3)+12*dgm01(1,2)*dgm10(3,3)+3*gm(2,3)*d2gm(1,3)+3*gm(1,3)& 422 & *d2gm(2,3)+12*gm(1,2)*d2gm(3,3))+3*gm(2,3)*(dgm01(3,3)*dgm10(1,3)& 423 & +dgm01(1,3)*dgm10(3,3)+gm(1,3)*d2gm(3,3)) 424 cm(9,3)=-4.5d0*gm(2,2)*dgm01(3,3)*dgm10(1,3)-4.5d0*gm(1,3)*dgm01(3,3)& 425 & *dgm10(2,2)+3*gm(1,3)*dgm01(2,3)*dgm10(2,3)+12*gm(1,2)*dgm01(3,3)& 426 & *dgm10(2,3)-4.5d0*gm(2,2)*dgm01(1,3)*dgm10(3,3)-4.5d0*gm(1,3)& 427 & *dgm01(2,2)*dgm10(3,3)+12*gm(1,2)*dgm01(2,3)*dgm10(3,3)+1.5d0*gm(2,3)& 428 & **2*d2gm(1,3)+gm(3,3)*(12*dgm01(2,3)*dgm10(1,2)-4.5d0*dgm01(2,2)& 429 & *dgm10(1,3)-4.5d0*dgm01(1,3)*dgm10(2,2)+12*dgm01(1,2)*dgm10(2,3)& 430 & +12*gm(2,3)*d2gm(1,2)-4.5d0*gm(2,2)*d2gm(1,3)-4.5d0*gm(1,3)*d2gm(2,2)& 431 & +12*gm(1,2)*d2gm(2,3))-4.5d0*gm(1,3)*gm(2,2)*d2gm(3,3)+gm(2,3)& 432 & *(12*dgm01(3,3)*dgm10(1,2)+3*dgm01(2,3)*dgm10(1,3)+3*dgm01(1,3)& 433 & *dgm10(2,3)+12*dgm01(1,2)*dgm10(3,3)+3*gm(1,3)*d2gm(2,3)+12*gm(1,2)& 434 & *d2gm(3,3)) 435 cm(10,3)=6*gm(1,3)*dgm01(3,3)*dgm10(3,3)+3*gm(3,3)**2*d2gm(1,3)& 436 & +6*gm(3,3)*(dgm01(3,3)*dgm10(1,3)+dgm01(1,3)*dgm10(3,3)+gm(1,3)& 437 & *d2gm(3,3)) 438 cm(1,4)=9*gm(1,2)*dgm01(1,3)*dgm10(1,1)-6*gm(1,1)*dgm01(2,3)& 439 & *dgm10(1,1)+9*gm(1,1)*dgm01(1,3)*dgm10(1,2)+9*gm(1,2)*dgm01(1,1)& 440 & *dgm10(1,3)+9*gm(1,1)*dgm01(1,2)*dgm10(1,3)-6*gm(1,1)*dgm01(1,1)& 441 & *dgm10(2,3)-6*gm(2,3)*(dgm01(1,1)*dgm10(1,1)+gm(1,1)*d2gm(1,1))& 442 & +9*gm(1,3)*(dgm01(1,2)*dgm10(1,1)+dgm01(1,1)*dgm10(1,2)+gm(1,2)& 443 & *d2gm(1,1)+gm(1,1)*d2gm(1,2))+9*gm(1,1)*gm(1,2)*d2gm(1,3)-3*gm(1,1)& 444 & **2*d2gm(2,3) 445 cm(2,4)=3*gm(1,3)*dgm01(2,2)*dgm10(1,2)+6*gm(1,2)*dgm01(2,3)& 446 & *dgm10(1,2)+3*gm(1,2)*dgm01(2,2)*dgm10(1,3)+3*gm(1,3)*dgm01(1,2)& 447 & *dgm10(2,2)+3*gm(1,2)*dgm01(1,3)*dgm10(2,2)+12*gm(1,1)*dgm01(2,3)& 448 & *dgm10(2,2)+6*gm(1,2)*dgm01(1,2)*dgm10(2,3)+12*gm(1,1)*dgm01(2,2)& 449 & *dgm10(2,3)+3*gm(1,2)*gm(1,3)*d2gm(2,2)+gm(2,3)*(12*dgm01(2,2)& 450 & *dgm10(1,1)+6*dgm01(1,2)*dgm10(1,2)+12*dgm01(1,1)*dgm10(2,2)& 451 & +12*gm(2,2)*d2gm(1,1)+6*gm(1,2)*d2gm(1,2)+12*gm(1,1)*d2gm(2,2))& 452 & +3*gm(1,2)**2*d2gm(2,3)+gm(2,2)*(12*dgm01(2,3)*dgm10(1,1)+3*dgm01(1,3)& 453 & *dgm10(1,2)+3*dgm01(1,2)*dgm10(1,3)+12*dgm01(1,1)*dgm10(2,3)& 454 & +3*gm(1,3)*d2gm(1,2)+3*gm(1,2)*d2gm(1,3)+12*gm(1,1)*d2gm(2,3)) 455 cm(3,4)=3*gm(1,3)*dgm01(3,3)*dgm10(1,2)+6*gm(1,3)*dgm01(2,3)& 456 & *dgm10(1,3)+3*gm(1,2)*dgm01(3,3)*dgm10(1,3)+6*gm(1,3)*dgm01(1,3)& 457 & *dgm10(2,3)+12*gm(1,1)*dgm01(3,3)*dgm10(2,3)+3*gm(1,3)*dgm01(1,2)& 458 & *dgm10(3,3)+3*gm(1,2)*dgm01(1,3)*dgm10(3,3)+12*gm(1,1)*dgm01(2,3)& 459 & *dgm10(3,3)+3*gm(1,3)**2*d2gm(2,3)+gm(3,3)*(12*dgm01(2,3)*dgm10(1,1)& 460 & +3*dgm01(1,3)*dgm10(1,2)+3*dgm01(1,2)*dgm10(1,3)+12*dgm01(1,1)& 461 & *dgm10(2,3)+12*gm(2,3)*d2gm(1,1)+3*gm(1,3)*d2gm(1,2)+3*gm(1,2)& 462 & *d2gm(1,3)+12*gm(1,1)*d2gm(2,3))+3*gm(1,2)*gm(1,3)*d2gm(3,3)& 463 & +gm(2,3)*(12*dgm01(3,3)*dgm10(1,1)+6*dgm01(1,3)*dgm10(1,3)+12*dgm01(1,1)& 464 & *dgm10(3,3)+6*gm(1,3)*d2gm(1,3)+12*gm(1,1)*d2gm(3,3)) 465 cm(4,4)=15*gm(2,2)*dgm01(3,3)*dgm10(1,1)-6*gm(1,3)*dgm01(2,3)& 466 & *dgm10(1,2)+18*gm(1,2)*dgm01(3,3)*dgm10(1,2)+18*gm(2,2)*dgm01(1,3)& 467 & *dgm10(1,3)+18*gm(1,3)*dgm01(2,2)*dgm10(1,3)-6*gm(1,2)*dgm01(2,3)& 468 & *dgm10(1,3)+18*gm(1,3)*dgm01(1,3)*dgm10(2,2)+15*gm(1,1)*dgm01(3,3)& 469 & *dgm10(2,2)-6*gm(1,3)*dgm01(1,2)*dgm10(2,3)-6*gm(1,2)*dgm01(1,3)& 470 & *dgm10(2,3)+18*gm(1,1)*dgm01(2,3)*dgm10(2,3)+15*gm(2,2)*dgm01(1,1)& 471 & *dgm10(3,3)+18*gm(1,2)*dgm01(1,2)*dgm10(3,3)+15*gm(1,1)*dgm01(2,2)& 472 & *dgm10(3,3)+9*gm(2,3)**2*d2gm(1,1)+18*gm(1,3)*gm(2,2)*d2gm(1,3)& 473 & +9*gm(1,3)**2*d2gm(2,2)+gm(3,3)*(15*dgm01(2,2)*dgm10(1,1)+18*dgm01(1,2)& 474 & *dgm10(1,2)+15*dgm01(1,1)*dgm10(2,2)+15*gm(2,2)*d2gm(1,1)+18*gm(1,2)& 475 & *d2gm(1,2)+15*gm(1,1)*d2gm(2,2))-6*gm(1,2)*gm(1,3)*d2gm(2,3)& 476 & +gm(2,3)*(18*dgm01(2,3)*dgm10(1,1)-6*dgm01(1,3)*dgm10(1,2)-6*dgm01(1,2)& 477 & *dgm10(1,3)+18*dgm01(1,1)*dgm10(2,3)-6*gm(1,3)*d2gm(1,2)-6*gm(1,2)& 478 & *d2gm(1,3)+18*gm(1,1)*d2gm(2,3))+9*gm(1,2)**2*d2gm(3,3)+15*gm(1,1)& 479 & *gm(2,2)*d2gm(3,3) 480 cm(5,4)=3*gm(1,3)*dgm01(2,3)*dgm10(1,1)+12*gm(1,2)*dgm01(3,3)& 481 & *dgm10(1,1)+6*gm(1,3)*dgm01(1,3)*dgm10(1,2)+12*gm(1,1)*dgm01(3,3)& 482 & *dgm10(1,2)+6*gm(1,3)*dgm01(1,2)*dgm10(1,3)+6*gm(1,2)*dgm01(1,3)& 483 & *dgm10(1,3)+3*gm(1,1)*dgm01(2,3)*dgm10(1,3)+3*gm(1,3)*dgm01(1,1)& 484 & *dgm10(2,3)+3*gm(1,1)*dgm01(1,3)*dgm10(2,3)+12*gm(1,2)*dgm01(1,1)& 485 & *dgm10(3,3)+12*gm(1,1)*dgm01(1,2)*dgm10(3,3)+3*gm(1,3)**2*d2gm(1,2)& 486 & +12*gm(3,3)*(dgm01(1,2)*dgm10(1,1)+dgm01(1,1)*dgm10(1,2)+gm(1,2)& 487 & *d2gm(1,1)+gm(1,1)*d2gm(1,2))+6*gm(1,2)*gm(1,3)*d2gm(1,3)+3*gm(2,3)& 488 & *(dgm01(1,3)*dgm10(1,1)+dgm01(1,1)*dgm10(1,3)+gm(1,3)*d2gm(1,1)& 489 & +gm(1,1)*d2gm(1,3))+3*gm(1,1)*gm(1,3)*d2gm(2,3)+12*gm(1,1)*gm(1,2)& 490 & *d2gm(3,3) 491 cm(6,4)=12*gm(1,3)*dgm01(2,2)*dgm10(1,1)+3*gm(1,2)*dgm01(2,3)& 492 & *dgm10(1,1)+6*gm(1,3)*dgm01(1,2)*dgm10(1,2)+6*gm(1,2)*dgm01(1,3)& 493 & *dgm10(1,2)+3*gm(1,1)*dgm01(2,3)*dgm10(1,2)+6*gm(1,2)*dgm01(1,2)& 494 & *dgm10(1,3)+12*gm(1,1)*dgm01(2,2)*dgm10(1,3)+12*gm(1,3)*dgm01(1,1)& 495 & *dgm10(2,2)+12*gm(1,1)*dgm01(1,3)*dgm10(2,2)+3*gm(1,2)*dgm01(1,1)& 496 & *dgm10(2,3)+3*gm(1,1)*dgm01(1,2)*dgm10(2,3)+6*gm(1,2)*gm(1,3)& 497 & *d2gm(1,2)+3*gm(2,3)*(dgm01(1,2)*dgm10(1,1)+dgm01(1,1)*dgm10(1,2)& 498 & +gm(1,2)*d2gm(1,1)+gm(1,1)*d2gm(1,2))+3*gm(1,2)**2*d2gm(1,3)& 499 & +12*gm(2,2)*(dgm01(1,3)*dgm10(1,1)+dgm01(1,1)*dgm10(1,3)+gm(1,3)& 500 & *d2gm(1,1)+gm(1,1)*d2gm(1,3))+12*gm(1,1)*gm(1,3)*d2gm(2,2)+3*gm(1,1)& 501 & *gm(1,2)*d2gm(2,3) 502 cm(7,4)=-6*gm(1,3)*dgm01(2,2)*dgm10(2,2)+9*gm(1,2)*dgm01(2,3)& 503 & *dgm10(2,2)+9*gm(1,2)*dgm01(2,2)*dgm10(2,3)-3*gm(2,2)**2*d2gm(1,3)& 504 & +9*gm(2,3)*(dgm01(2,2)*dgm10(1,2)+dgm01(1,2)*dgm10(2,2)+gm(2,2)& 505 & *d2gm(1,2)+gm(1,2)*d2gm(2,2))+gm(2,2)*(9*dgm01(2,3)*dgm10(1,2)& 506 & -6*dgm01(2,2)*dgm10(1,3)-6*dgm01(1,3)*dgm10(2,2)+9*dgm01(1,2)& 507 & *dgm10(2,3)-6*gm(1,3)*d2gm(2,2)+9*gm(1,2)*d2gm(2,3)) 508 cm(8,4)=12*gm(2,2)*dgm01(3,3)*dgm10(1,3)+12*gm(1,3)*dgm01(3,3)& 509 & *dgm10(2,2)+6*gm(1,3)*dgm01(2,3)*dgm10(2,3)+3*gm(1,2)*dgm01(3,3)& 510 & *dgm10(2,3)+12*gm(2,2)*dgm01(1,3)*dgm10(3,3)+12*gm(1,3)*dgm01(2,2)& 511 & *dgm10(3,3)+3*gm(1,2)*dgm01(2,3)*dgm10(3,3)+3*gm(2,3)**2*d2gm(1,3)& 512 & +gm(3,3)*(3*dgm01(2,3)*dgm10(1,2)+12*dgm01(2,2)*dgm10(1,3)+12*dgm01(1,3)& 513 & *dgm10(2,2)+3*dgm01(1,2)*dgm10(2,3)+3*gm(2,3)*d2gm(1,2)+12*gm(2,2)& 514 & *d2gm(1,3)+12*gm(1,3)*d2gm(2,2)+3*gm(1,2)*d2gm(2,3))+12*gm(1,3)& 515 & *gm(2,2)*d2gm(3,3)+gm(2,3)*(3*dgm01(3,3)*dgm10(1,2)+6*dgm01(2,3)& 516 & *dgm10(1,3)+6*dgm01(1,3)*dgm10(2,3)+3*dgm01(1,2)*dgm10(3,3)+6*gm(1,3)& 517 & *d2gm(2,3)+3*gm(1,2)*d2gm(3,3)) 518 cm(9,4)=12*gm(2,2)*dgm01(3,3)*dgm10(1,2)+3*gm(2,2)*dgm01(2,3)& 519 & *dgm10(1,3)+3*gm(1,3)*dgm01(2,3)*dgm10(2,2)+12*gm(1,2)*dgm01(3,3)& 520 & *dgm10(2,2)+3*gm(2,2)*dgm01(1,3)*dgm10(2,3)+3*gm(1,3)*dgm01(2,2)& 521 & *dgm10(2,3)+6*gm(1,2)*dgm01(2,3)*dgm10(2,3)+12*gm(2,2)*dgm01(1,2)& 522 & *dgm10(3,3)+12*gm(1,2)*dgm01(2,2)*dgm10(3,3)+3*gm(2,3)**2*d2gm(1,2)& 523 & +12*gm(3,3)*(dgm01(2,2)*dgm10(1,2)+dgm01(1,2)*dgm10(2,2)+gm(2,2)& 524 & *d2gm(1,2)+gm(1,2)*d2gm(2,2))+3*gm(1,3)*gm(2,2)*d2gm(2,3)+gm(2,3)& 525 & *(6*dgm01(2,3)*dgm10(1,2)+3*dgm01(2,2)*dgm10(1,3)+3*dgm01(1,3)& 526 & *dgm10(2,2)+6*dgm01(1,2)*dgm10(2,3)+3*gm(2,2)*d2gm(1,3)+3*gm(1,3)& 527 & *d2gm(2,2)+6*gm(1,2)*d2gm(2,3))+12*gm(1,2)*gm(2,2)*d2gm(3,3) 528 cm(10,4)=9*gm(1,3)*dgm01(3,3)*dgm10(2,3)+9*gm(1,3)*dgm01(2,3)& 529 & *dgm10(3,3)-6*gm(1,2)*dgm01(3,3)*dgm10(3,3)-3*gm(3,3)**2*d2gm(1,2)& 530 & +gm(3,3)*(-6*dgm01(3,3)*dgm10(1,2)+9*dgm01(2,3)*dgm10(1,3)+9*dgm01(1,3)& 531 & *dgm10(2,3)-6*dgm01(1,2)*dgm10(3,3)+9*gm(2,3)*d2gm(1,3)+9*gm(1,3)& 532 & *d2gm(2,3)-6*gm(1,2)*d2gm(3,3))+9*gm(2,3)*(dgm01(3,3)*dgm10(1,3)& 533 & +dgm01(1,3)*dgm10(3,3)+gm(1,3)*d2gm(3,3)) 534 cm(1,5)=6*gm(1,3)*(dgm01(1,1)*dgm10(1,1)+gm(1,1)*d2gm(1,1))& 535 & +gm(1,1)*(6*dgm01(1,3)*dgm10(1,1)+6*dgm01(1,1)*dgm10(1,3)+3*gm(1,1)& 536 & *d2gm(1,3)) 537 cm(2,5)=-4.5d0*gm(1,3)*dgm01(2,2)*dgm10(1,1)+12*gm(1,2)*dgm01(2,3)& 538 & *dgm10(1,1)+3*gm(1,3)*dgm01(1,2)*dgm10(1,2)+3*gm(1,2)*dgm01(1,3)& 539 & *dgm10(1,2)+12*gm(1,1)*dgm01(2,3)*dgm10(1,2)+3*gm(1,2)*dgm01(1,2)& 540 & *dgm10(1,3)-4.5d0*gm(1,1)*dgm01(2,2)*dgm10(1,3)-4.5d0*gm(1,3)& 541 & *dgm01(1,1)*dgm10(2,2)-4.5d0*gm(1,1)*dgm01(1,3)*dgm10(2,2)+12*gm(1,2)& 542 & *dgm01(1,1)*dgm10(2,3)+12*gm(1,1)*dgm01(1,2)*dgm10(2,3)+3*gm(1,2)& 543 & *gm(1,3)*d2gm(1,2)+12*gm(2,3)*(dgm01(1,2)*dgm10(1,1)+dgm01(1,1)& 544 & *dgm10(1,2)+gm(1,2)*d2gm(1,1)+gm(1,1)*d2gm(1,2))+1.5d0*gm(1,2)& 545 & **2*d2gm(1,3)-4.5d0*gm(2,2)*(dgm01(1,3)*dgm10(1,1)+dgm01(1,1)& 546 & *dgm10(1,3)+gm(1,3)*d2gm(1,1)+gm(1,1)*d2gm(1,3))-4.5d0*gm(1,1)& 547 & *gm(1,3)*d2gm(2,2)+12*gm(1,1)*gm(1,2)*d2gm(2,3) 548 cm(3,5)=7.5d0*gm(1,1)*(dgm01(3,3)*dgm10(1,3)+dgm01(1,3)*dgm10(3,3))& 549 & +4.5d0*gm(1,3)**2*d2gm(1,3)+7.5d0*gm(3,3)*(dgm01(1,3)*dgm10(1,1)& 550 & +dgm01(1,1)*dgm10(1,3)+gm(1,3)*d2gm(1,1)+gm(1,1)*d2gm(1,3))+gm(1,3)& 551 & *(7.5d0*dgm01(3,3)*dgm10(1,1)+9*dgm01(1,3)*dgm10(1,3)+7.5d0*dgm01(1,1)& 552 & *dgm10(3,3)+7.5d0*gm(1,1)*d2gm(3,3)) 553 cm(4,5)=3*gm(1,3)*dgm01(2,3)*dgm10(1,1)+12*gm(1,2)*dgm01(3,3)& 554 & *dgm10(1,1)+6*gm(1,3)*dgm01(1,3)*dgm10(1,2)+12*gm(1,1)*dgm01(3,3)& 555 & *dgm10(1,2)+6*gm(1,3)*dgm01(1,2)*dgm10(1,3)+6*gm(1,2)*dgm01(1,3)& 556 & *dgm10(1,3)+3*gm(1,1)*dgm01(2,3)*dgm10(1,3)+3*gm(1,3)*dgm01(1,1)& 557 & *dgm10(2,3)+3*gm(1,1)*dgm01(1,3)*dgm10(2,3)+12*gm(1,2)*dgm01(1,1)& 558 & *dgm10(3,3)+12*gm(1,1)*dgm01(1,2)*dgm10(3,3)+3*gm(1,3)**2*d2gm(1,2)& 559 & +12*gm(3,3)*(dgm01(1,2)*dgm10(1,1)+dgm01(1,1)*dgm10(1,2)+gm(1,2)& 560 & *d2gm(1,1)+gm(1,1)*d2gm(1,2))+6*gm(1,2)*gm(1,3)*d2gm(1,3)+3*gm(2,3)& 561 & *(dgm01(1,3)*dgm10(1,1)+dgm01(1,1)*dgm10(1,3)+gm(1,3)*d2gm(1,1)& 562 & +gm(1,1)*d2gm(1,3))+3*gm(1,1)*gm(1,3)*d2gm(2,3)+12*gm(1,1)*gm(1,2)& 563 & *d2gm(3,3) 564 cm(5,5)=3*gm(1,3)**2*d2gm(1,1)+12*gm(3,3)*(dgm01(1,1)*dgm10(1,1)& 565 & +gm(1,1)*d2gm(1,1))+6*gm(1,3)*(dgm01(1,3)*dgm10(1,1)+dgm01(1,1)& 566 & *dgm10(1,3)+gm(1,1)*d2gm(1,3))+gm(1,1)*(12*dgm01(3,3)*dgm10(1,1)& 567 & +6*dgm01(1,3)*dgm10(1,3)+12*dgm01(1,1)*dgm10(3,3)+6*gm(1,1)*d2gm(3,3)) 568 cm(6,5)=3*gm(1,2)*dgm01(1,3)*dgm10(1,1)+12*gm(1,1)*dgm01(2,3)& 569 & *dgm10(1,1)+3*gm(1,1)*dgm01(1,3)*dgm10(1,2)+3*gm(1,2)*dgm01(1,1)& 570 & *dgm10(1,3)+3*gm(1,1)*dgm01(1,2)*dgm10(1,3)+12*gm(1,1)*dgm01(1,1)& 571 & *dgm10(2,3)+12*gm(2,3)*(dgm01(1,1)*dgm10(1,1)+gm(1,1)*d2gm(1,1))& 572 & +3*gm(1,3)*(dgm01(1,2)*dgm10(1,1)+dgm01(1,1)*dgm10(1,2)+gm(1,2)& 573 & *d2gm(1,1)+gm(1,1)*d2gm(1,2))+3*gm(1,1)*gm(1,2)*d2gm(1,3)+6*gm(1,1)& 574 & **2*d2gm(2,3) 575 cm(7,5)=(-36*gm(1,3)*dgm01(2,2)*dgm10(1,2)+180*gm(1,2)*dgm01(2,3)& 576 & *dgm10(1,2)-36*gm(1,2)*dgm01(2,2)*dgm10(1,3)-36*gm(1,3)*dgm01(1,2)& 577 & *dgm10(2,2)-36*gm(1,2)*dgm01(1,3)*dgm10(2,2)-18*gm(1,1)*dgm01(2,3)& 578 & *dgm10(2,2)+180*gm(1,2)*dgm01(1,2)*dgm10(2,3)-18*gm(1,1)*dgm01(2,2)& 579 & *dgm10(2,3)-36*gm(1,2)*gm(1,3)*d2gm(2,2)+gm(2,3)*(-18*dgm01(2,2)& 580 & *dgm10(1,1)+180*dgm01(1,2)*dgm10(1,2)-18*dgm01(1,1)*dgm10(2,2)& 581 & -18*gm(2,2)*d2gm(1,1)+180*gm(1,2)*d2gm(1,2)-18*gm(1,1)*d2gm(2,2))& 582 & +90*gm(1,2)**2*d2gm(2,3)+gm(2,2)*(-18*dgm01(2,3)*dgm10(1,1)-36*dgm01(1,3)& 583 & *dgm10(1,2)-36*dgm01(1,2)*dgm10(1,3)-18*dgm01(1,1)*dgm10(2,3)& 584 & -36*gm(1,3)*d2gm(1,2)-36*gm(1,2)*d2gm(1,3)-18*gm(1,1)*d2gm(2,3)))& 585 & /12.d0 586 cm(8,5)=12*gm(1,3)*dgm01(3,3)*dgm10(1,2)+3*gm(1,3)*dgm01(2,3)& 587 & *dgm10(1,3)+12*gm(1,2)*dgm01(3,3)*dgm10(1,3)+3*gm(1,3)*dgm01(1,3)& 588 & *dgm10(2,3)-4.5d0*gm(1,1)*dgm01(3,3)*dgm10(2,3)+12*gm(1,3)*dgm01(1,2)& 589 & *dgm10(3,3)+12*gm(1,2)*dgm01(1,3)*dgm10(3,3)-4.5d0*gm(1,1)*dgm01(2,3)& 590 & *dgm10(3,3)+1.5d0*gm(1,3)**2*d2gm(2,3)+gm(3,3)*(-4.5d0*dgm01(2,3)& 591 & *dgm10(1,1)+12*dgm01(1,3)*dgm10(1,2)+12*dgm01(1,2)*dgm10(1,3)& 592 & -4.5d0*dgm01(1,1)*dgm10(2,3)-4.5d0*gm(2,3)*d2gm(1,1)+12*gm(1,3)& 593 & *d2gm(1,2)+12*gm(1,2)*d2gm(1,3)-4.5d0*gm(1,1)*d2gm(2,3))+12*gm(1,2)& 594 & *gm(1,3)*d2gm(3,3)+gm(2,3)*(-4.5d0*dgm01(3,3)*dgm10(1,1)+3*dgm01(1,3)& 595 & *dgm10(1,3)-4.5d0*dgm01(1,1)*dgm10(3,3)+3*gm(1,3)*d2gm(1,3)-4.5d0*gm(1,1)& 596 & *d2gm(3,3)) 597 cm(9,5)=-1.5d0*gm(2,2)*dgm01(3,3)*dgm10(1,1)+9*gm(1,3)*dgm01(2,3)& 598 & *dgm10(1,2)+15*gm(1,2)*dgm01(3,3)*dgm10(1,2)-6*gm(2,2)*dgm01(1,3)& 599 & *dgm10(1,3)-6*gm(1,3)*dgm01(2,2)*dgm10(1,3)+9*gm(1,2)*dgm01(2,3)& 600 & *dgm10(1,3)-6*gm(1,3)*dgm01(1,3)*dgm10(2,2)-1.5d0*gm(1,1)*dgm01(3,3)& 601 & *dgm10(2,2)+9*gm(1,3)*dgm01(1,2)*dgm10(2,3)+9*gm(1,2)*dgm01(1,3)& 602 & *dgm10(2,3)-6*gm(1,1)*dgm01(2,3)*dgm10(2,3)-1.5d0*gm(2,2)*dgm01(1,1)& 603 & *dgm10(3,3)+15*gm(1,2)*dgm01(1,2)*dgm10(3,3)-1.5d0*gm(1,1)*dgm01(2,2)& 604 & *dgm10(3,3)-3*gm(2,3)**2*d2gm(1,1)-6*gm(1,3)*gm(2,2)*d2gm(1,3)& 605 & -3*gm(1,3)**2*d2gm(2,2)+gm(3,3)*(-1.5d0*dgm01(2,2)*dgm10(1,1)& 606 & +15*dgm01(1,2)*dgm10(1,2)-1.5d0*dgm01(1,1)*dgm10(2,2)-1.5d0*gm(2,2)& 607 & *d2gm(1,1)+15*gm(1,2)*d2gm(1,2)-1.5d0*gm(1,1)*d2gm(2,2))+9*gm(1,2)& 608 & *gm(1,3)*d2gm(2,3)+gm(2,3)*(-6*dgm01(2,3)*dgm10(1,1)+9*dgm01(1,3)& 609 & *dgm10(1,2)+9*dgm01(1,2)*dgm10(1,3)-6*dgm01(1,1)*dgm10(2,3)+9*gm(1,3)& 610 & *d2gm(1,2)+9*gm(1,2)*d2gm(1,3)-6*gm(1,1)*d2gm(2,3))+7.5d0*gm(1,2)& 611 & **2*d2gm(3,3)-1.5d0*gm(1,1)*gm(2,2)*d2gm(3,3) 612 cm(10,5)=-3*gm(1,1)*dgm01(3,3)*dgm10(3,3)+9*gm(1,3)*(dgm01(3,3)& 613 & *dgm10(1,3)+dgm01(1,3)*dgm10(3,3))-1.5d0*gm(3,3)**2*d2gm(1,1)& 614 & +4.5d0*gm(1,3)**2*d2gm(3,3)+gm(3,3)*(-3*dgm01(3,3)*dgm10(1,1)& 615 & +9*dgm01(1,3)*dgm10(1,3)-3*dgm01(1,1)*dgm10(3,3)+9*gm(1,3)*d2gm(1,3)& 616 & -3*gm(1,1)*d2gm(3,3)) 617 cm(1,6)=6*gm(1,2)*(dgm01(1,1)*dgm10(1,1)+gm(1,1)*d2gm(1,1))& 618 & +gm(1,1)*(6*dgm01(1,2)*dgm10(1,1)+6*dgm01(1,1)*dgm10(1,2)+3*gm(1,1)& 619 & *d2gm(1,2)) 620 cm(2,6)=7.5d0*gm(1,1)*(dgm01(2,2)*dgm10(1,2)+dgm01(1,2)*dgm10(2,2))& 621 & +4.5d0*gm(1,2)**2*d2gm(1,2)+7.5d0*gm(2,2)*(dgm01(1,2)*dgm10(1,1)& 622 & +dgm01(1,1)*dgm10(1,2)+gm(1,2)*d2gm(1,1)+gm(1,1)*d2gm(1,2))+gm(1,2)& 623 & *(7.5d0*dgm01(2,2)*dgm10(1,1)+9*dgm01(1,2)*dgm10(1,2)+7.5d0*dgm01(1,1)& 624 & *dgm10(2,2)+7.5d0*gm(1,1)*d2gm(2,2)) 625 cm(3,6)=12*gm(1,3)*dgm01(2,3)*dgm10(1,1)-4.5d0*gm(1,2)*dgm01(3,3)& 626 & *dgm10(1,1)+3*gm(1,3)*dgm01(1,3)*dgm10(1,2)-4.5d0*gm(1,1)*dgm01(3,3)& 627 & *dgm10(1,2)+3*gm(1,3)*dgm01(1,2)*dgm10(1,3)+3*gm(1,2)*dgm01(1,3)& 628 & *dgm10(1,3)+12*gm(1,1)*dgm01(2,3)*dgm10(1,3)+12*gm(1,3)*dgm01(1,1)& 629 & *dgm10(2,3)+12*gm(1,1)*dgm01(1,3)*dgm10(2,3)-4.5d0*gm(1,2)*dgm01(1,1)& 630 & *dgm10(3,3)-4.5d0*gm(1,1)*dgm01(1,2)*dgm10(3,3)+1.5d0*gm(1,3)& 631 & **2*d2gm(1,2)-4.5d0*gm(3,3)*(dgm01(1,2)*dgm10(1,1)+dgm01(1,1)& 632 & *dgm10(1,2)+gm(1,2)*d2gm(1,1)+gm(1,1)*d2gm(1,2))+3*gm(1,2)*gm(1,3)& 633 & *d2gm(1,3)+12*gm(2,3)*(dgm01(1,3)*dgm10(1,1)+dgm01(1,1)*dgm10(1,3)& 634 & +gm(1,3)*d2gm(1,1)+gm(1,1)*d2gm(1,3))+12*gm(1,1)*gm(1,3)*d2gm(2,3)& 635 & -4.5d0*gm(1,1)*gm(1,2)*d2gm(3,3) 636 cm(4,6)=12*gm(1,3)*dgm01(2,2)*dgm10(1,1)+3*gm(1,2)*dgm01(2,3)& 637 & *dgm10(1,1)+6*gm(1,3)*dgm01(1,2)*dgm10(1,2)+6*gm(1,2)*dgm01(1,3)& 638 & *dgm10(1,2)+3*gm(1,1)*dgm01(2,3)*dgm10(1,2)+6*gm(1,2)*dgm01(1,2)& 639 & *dgm10(1,3)+12*gm(1,1)*dgm01(2,2)*dgm10(1,3)+12*gm(1,3)*dgm01(1,1)& 640 & *dgm10(2,2)+12*gm(1,1)*dgm01(1,3)*dgm10(2,2)+3*gm(1,2)*dgm01(1,1)& 641 & *dgm10(2,3)+3*gm(1,1)*dgm01(1,2)*dgm10(2,3)+6*gm(1,2)*gm(1,3)& 642 & *d2gm(1,2)+3*gm(2,3)*(dgm01(1,2)*dgm10(1,1)+dgm01(1,1)*dgm10(1,2)& 643 & +gm(1,2)*d2gm(1,1)+gm(1,1)*d2gm(1,2))+3*gm(1,2)**2*d2gm(1,3)& 644 & +12*gm(2,2)*(dgm01(1,3)*dgm10(1,1)+dgm01(1,1)*dgm10(1,3)+gm(1,3)& 645 & *d2gm(1,1)+gm(1,1)*d2gm(1,3))+12*gm(1,1)*gm(1,3)*d2gm(2,2)+3*gm(1,1)& 646 & *gm(1,2)*d2gm(2,3) 647 cm(5,6)=3*gm(1,2)*dgm01(1,3)*dgm10(1,1)+12*gm(1,1)*dgm01(2,3)& 648 & *dgm10(1,1)+3*gm(1,1)*dgm01(1,3)*dgm10(1,2)+3*gm(1,2)*dgm01(1,1)& 649 & *dgm10(1,3)+3*gm(1,1)*dgm01(1,2)*dgm10(1,3)+12*gm(1,1)*dgm01(1,1)& 650 & *dgm10(2,3)+12*gm(2,3)*(dgm01(1,1)*dgm10(1,1)+gm(1,1)*d2gm(1,1))& 651 & +3*gm(1,3)*(dgm01(1,2)*dgm10(1,1)+dgm01(1,1)*dgm10(1,2)+gm(1,2)& 652 & *d2gm(1,1)+gm(1,1)*d2gm(1,2))+3*gm(1,1)*gm(1,2)*d2gm(1,3)+6*gm(1,1)& 653 & **2*d2gm(2,3) 654 cm(6,6)=3*gm(1,2)**2*d2gm(1,1)+12*gm(2,2)*(dgm01(1,1)*dgm10(1,1)& 655 & +gm(1,1)*d2gm(1,1))+6*gm(1,2)*(dgm01(1,2)*dgm10(1,1)+dgm01(1,1)& 656 & *dgm10(1,2)+gm(1,1)*d2gm(1,2))+gm(1,1)*(12*dgm01(2,2)*dgm10(1,1)& 657 & +6*dgm01(1,2)*dgm10(1,2)+12*dgm01(1,1)*dgm10(2,2)+6*gm(1,1)*d2gm(2,2)) 658 cm(7,6)=-3*gm(1,1)*dgm01(2,2)*dgm10(2,2)+9*gm(1,2)*(dgm01(2,2)& 659 & *dgm10(1,2)+dgm01(1,2)*dgm10(2,2))-1.5d0*gm(2,2)**2*d2gm(1,1)& 660 & +4.5d0*gm(1,2)**2*d2gm(2,2)+gm(2,2)*(-3*dgm01(2,2)*dgm10(1,1)& 661 & +9*dgm01(1,2)*dgm10(1,2)-3*dgm01(1,1)*dgm10(2,2)+9*gm(1,2)*d2gm(1,2)& 662 & -3*gm(1,1)*d2gm(2,2)) 663 cm(8,6)=-1.5d0*gm(2,2)*dgm01(3,3)*dgm10(1,1)+9*gm(1,3)*dgm01(2,3)& 664 & *dgm10(1,2)-6*gm(1,2)*dgm01(3,3)*dgm10(1,2)+15*gm(2,2)*dgm01(1,3)& 665 & *dgm10(1,3)+15*gm(1,3)*dgm01(2,2)*dgm10(1,3)+9*gm(1,2)*dgm01(2,3)& 666 & *dgm10(1,3)+15*gm(1,3)*dgm01(1,3)*dgm10(2,2)-1.5d0*gm(1,1)*dgm01(3,3)& 667 & *dgm10(2,2)+9*gm(1,3)*dgm01(1,2)*dgm10(2,3)+9*gm(1,2)*dgm01(1,3)& 668 & *dgm10(2,3)-6*gm(1,1)*dgm01(2,3)*dgm10(2,3)-1.5d0*gm(2,2)*dgm01(1,1)& 669 & *dgm10(3,3)-6*gm(1,2)*dgm01(1,2)*dgm10(3,3)-1.5d0*gm(1,1)*dgm01(2,2)& 670 & *dgm10(3,3)-3*gm(2,3)**2*d2gm(1,1)+15*gm(1,3)*gm(2,2)*d2gm(1,3)& 671 & +7.5d0*gm(1,3)**2*d2gm(2,2)+gm(3,3)*(-1.5d0*dgm01(2,2)*dgm10(1,1)& 672 & -6*dgm01(1,2)*dgm10(1,2)-1.5d0*dgm01(1,1)*dgm10(2,2)-1.5d0*gm(2,2)& 673 & *d2gm(1,1)-6*gm(1,2)*d2gm(1,2)-1.5d0*gm(1,1)*d2gm(2,2))+9*gm(1,2)& 674 & *gm(1,3)*d2gm(2,3)+gm(2,3)*(-6*dgm01(2,3)*dgm10(1,1)+9*dgm01(1,3)& 675 & *dgm10(1,2)+9*dgm01(1,2)*dgm10(1,3)-6*dgm01(1,1)*dgm10(2,3)+9*gm(1,3)& 676 & *d2gm(1,2)+9*gm(1,2)*d2gm(1,3)-6*gm(1,1)*d2gm(2,3))-3*gm(1,2)& 677 & **2*d2gm(3,3)-1.5d0*gm(1,1)*gm(2,2)*d2gm(3,3) 678 cm(9,6)=12*gm(1,3)*dgm01(2,2)*dgm10(1,2)+3*gm(1,2)*dgm01(2,3)& 679 & *dgm10(1,2)+12*gm(1,2)*dgm01(2,2)*dgm10(1,3)+12*gm(1,3)*dgm01(1,2)& 680 & *dgm10(2,2)+12*gm(1,2)*dgm01(1,3)*dgm10(2,2)-4.5d0*gm(1,1)*dgm01(2,3)& 681 & *dgm10(2,2)+3*gm(1,2)*dgm01(1,2)*dgm10(2,3)-4.5d0*gm(1,1)*dgm01(2,2)& 682 & *dgm10(2,3)+12*gm(1,2)*gm(1,3)*d2gm(2,2)+gm(2,3)*(-4.5d0*dgm01(2,2)& 683 & *dgm10(1,1)+3*dgm01(1,2)*dgm10(1,2)-4.5d0*dgm01(1,1)*dgm10(2,2)& 684 & -4.5d0*gm(2,2)*d2gm(1,1)+3*gm(1,2)*d2gm(1,2)-4.5d0*gm(1,1)*d2gm(2,2))& 685 & +1.5d0*gm(1,2)**2*d2gm(2,3)+gm(2,2)*(-4.5d0*dgm01(2,3)*dgm10(1,1)& 686 & +12*dgm01(1,3)*dgm10(1,2)+12*dgm01(1,2)*dgm10(1,3)-4.5d0*dgm01(1,1)& 687 & *dgm10(2,3)+12*gm(1,3)*d2gm(1,2)+12*gm(1,2)*d2gm(1,3)-4.5d0*gm(1,1)& 688 & *d2gm(2,3)) 689 cm(10,6)=(-36*gm(1,3)*dgm01(3,3)*dgm10(1,2)+180*gm(1,3)*dgm01(2,3)& 690 & *dgm10(1,3)-36*gm(1,2)*dgm01(3,3)*dgm10(1,3)+180*gm(1,3)*dgm01(1,3)& 691 & *dgm10(2,3)-18*gm(1,1)*dgm01(3,3)*dgm10(2,3)-36*gm(1,3)*dgm01(1,2)& 692 & *dgm10(3,3)-36*gm(1,2)*dgm01(1,3)*dgm10(3,3)-18*gm(1,1)*dgm01(2,3)& 693 & *dgm10(3,3)+90*gm(1,3)**2*d2gm(2,3)+gm(3,3)*(-18*dgm01(2,3)*dgm10(1,1)& 694 & -36*dgm01(1,3)*dgm10(1,2)-36*dgm01(1,2)*dgm10(1,3)-18*dgm01(1,1)& 695 & *dgm10(2,3)-18*gm(2,3)*d2gm(1,1)-36*gm(1,3)*d2gm(1,2)-36*gm(1,2)& 696 & *d2gm(1,3)-18*gm(1,1)*d2gm(2,3))-36*gm(1,2)*gm(1,3)*d2gm(3,3)& 697 & +gm(2,3)*(-18*dgm01(3,3)*dgm10(1,1)+180*dgm01(1,3)*dgm10(1,3)& 698 & -18*dgm01(1,1)*dgm10(3,3)+180*gm(1,3)*d2gm(1,3)-18*gm(1,1)*d2gm(3,3)))& 699 & /12.d0 700 cm(1,7)=-1.5d0*gm(1,1)*(dgm01(2,2)*dgm10(1,2)+dgm01(1,2)*dgm10(2,2))& 701 & +7.5d0*gm(1,2)**2*d2gm(1,2)-1.5d0*gm(2,2)*(dgm01(1,2)*dgm10(1,1)& 702 & +dgm01(1,1)*dgm10(1,2)+gm(1,2)*d2gm(1,1)+gm(1,1)*d2gm(1,2))+gm(1,2)& 703 & *(-1.5d0*dgm01(2,2)*dgm10(1,1)+15*dgm01(1,2)*dgm10(1,2)-1.5d0*dgm01(1,1)& 704 & *dgm10(2,2)-1.5d0*gm(1,1)*d2gm(2,2)) 705 cm(2,7)=6*gm(1,2)*dgm01(2,2)*dgm10(2,2)+3*gm(2,2)**2*d2gm(1,2)& 706 & +6*gm(2,2)*(dgm01(2,2)*dgm10(1,2)+dgm01(1,2)*dgm10(2,2)+gm(1,2)& 707 & *d2gm(2,2)) 708 cm(3,7)=-1.5d0*gm(2,2)*dgm01(3,3)*dgm10(1,2)-3*gm(2,2)*dgm01(2,3)& 709 & *dgm10(1,3)-3*gm(1,3)*dgm01(2,3)*dgm10(2,2)-1.5d0*gm(1,2)*dgm01(3,3)& 710 & *dgm10(2,2)-3*gm(2,2)*dgm01(1,3)*dgm10(2,3)-3*gm(1,3)*dgm01(2,2)& 711 & *dgm10(2,3)+15*gm(1,2)*dgm01(2,3)*dgm10(2,3)-1.5d0*gm(2,2)*dgm01(1,2)& 712 & *dgm10(3,3)-1.5d0*gm(1,2)*dgm01(2,2)*dgm10(3,3)+7.5d0*gm(2,3)& 713 & **2*d2gm(1,2)-1.5d0*gm(3,3)*(dgm01(2,2)*dgm10(1,2)+dgm01(1,2)& 714 & *dgm10(2,2)+gm(2,2)*d2gm(1,2)+gm(1,2)*d2gm(2,2))-3*gm(1,3)*gm(2,2)& 715 & *d2gm(2,3)+gm(2,3)*(15*dgm01(2,3)*dgm10(1,2)-3*dgm01(2,2)*dgm10(1,3)& 716 & -3*dgm01(1,3)*dgm10(2,2)+15*dgm01(1,2)*dgm10(2,3)-3*gm(2,2)*d2gm(1,3)& 717 & -3*gm(1,3)*d2gm(2,2)+15*gm(1,2)*d2gm(2,3))-1.5d0*gm(1,2)*gm(2,2)& 718 & *d2gm(3,3) 719 cm(4,7)=-6*gm(1,3)*dgm01(2,2)*dgm10(2,2)+9*gm(1,2)*dgm01(2,3)& 720 & *dgm10(2,2)+9*gm(1,2)*dgm01(2,2)*dgm10(2,3)-3*gm(2,2)**2*d2gm(1,3)& 721 & +9*gm(2,3)*(dgm01(2,2)*dgm10(1,2)+dgm01(1,2)*dgm10(2,2)+gm(2,2)& 722 & *d2gm(1,2)+gm(1,2)*d2gm(2,2))+gm(2,2)*(9*dgm01(2,3)*dgm10(1,2)& 723 & -6*dgm01(2,2)*dgm10(1,3)-6*dgm01(1,3)*dgm10(2,2)+9*dgm01(1,2)& 724 & *dgm10(2,3)-6*gm(1,3)*d2gm(2,2)+9*gm(1,2)*d2gm(2,3)) 725 cm(5,7)=-3*gm(1,3)*dgm01(2,2)*dgm10(1,2)+15*gm(1,2)*dgm01(2,3)& 726 & *dgm10(1,2)-3*gm(1,2)*dgm01(2,2)*dgm10(1,3)-3*gm(1,3)*dgm01(1,2)& 727 & *dgm10(2,2)-3*gm(1,2)*dgm01(1,3)*dgm10(2,2)-1.5d0*gm(1,1)*dgm01(2,3)& 728 & *dgm10(2,2)+15*gm(1,2)*dgm01(1,2)*dgm10(2,3)-1.5d0*gm(1,1)*dgm01(2,2)& 729 & *dgm10(2,3)-3*gm(1,2)*gm(1,3)*d2gm(2,2)+gm(2,3)*(-1.5d0*dgm01(2,2)& 730 & *dgm10(1,1)+15*dgm01(1,2)*dgm10(1,2)-1.5d0*dgm01(1,1)*dgm10(2,2)& 731 & -1.5d0*gm(2,2)*d2gm(1,1)+15*gm(1,2)*d2gm(1,2)-1.5d0*gm(1,1)*d2gm(2,2))& 732 & +7.5d0*gm(1,2)**2*d2gm(2,3)+gm(2,2)*(-1.5d0*dgm01(2,3)*dgm10(1,1)& 733 & -3*dgm01(1,3)*dgm10(1,2)-3*dgm01(1,2)*dgm10(1,3)-1.5d0*dgm01(1,1)& 734 & *dgm10(2,3)-3*gm(1,3)*d2gm(1,2)-3*gm(1,2)*d2gm(1,3)-1.5d0*gm(1,1)& 735 & *d2gm(2,3)) 736 cm(6,7)=-3*gm(1,1)*dgm01(2,2)*dgm10(2,2)+9*gm(1,2)*(dgm01(2,2)& 737 & *dgm10(1,2)+dgm01(1,2)*dgm10(2,2))-1.5d0*gm(2,2)**2*d2gm(1,1)& 738 & +4.5d0*gm(1,2)**2*d2gm(2,2)+gm(2,2)*(-3*dgm01(2,2)*dgm10(1,1)& 739 & +9*dgm01(1,2)*dgm10(1,2)-3*dgm01(1,1)*dgm10(2,2)+9*gm(1,2)*d2gm(1,2)& 740 & -3*gm(1,1)*d2gm(2,2)) 741 cm(7,7)=gm(2,2)*(6*dgm01(2,2)*dgm10(2,2)+3*gm(2,2)*d2gm(2,2)) 742 cm(8,7)=4.5d0*gm(2,3)**2*d2gm(2,2)-3*gm(3,3)*(dgm01(2,2)*dgm10(2,2)& 743 & +gm(2,2)*d2gm(2,2))+9*gm(2,3)*(dgm01(2,3)*dgm10(2,2)+dgm01(2,2)& 744 & *dgm10(2,3)+gm(2,2)*d2gm(2,3))+gm(2,2)*(-3*dgm01(3,3)*dgm10(2,2)& 745 & +9*dgm01(2,3)*dgm10(2,3)-3*dgm01(2,2)*dgm10(3,3)-1.5d0*gm(2,2)& 746 & *d2gm(3,3)) 747 cm(9,7)=6*gm(2,3)*(dgm01(2,2)*dgm10(2,2)+gm(2,2)*d2gm(2,2))& 748 & +gm(2,2)*(6*dgm01(2,3)*dgm10(2,2)+6*dgm01(2,2)*dgm10(2,3)+3*gm(2,2)& 749 & *d2gm(2,3)) 750 cm(10,7)=-1.5d0*gm(2,2)*(dgm01(3,3)*dgm10(2,3)+dgm01(2,3)*dgm10(3,3))& 751 & +7.5d0*gm(2,3)**2*d2gm(2,3)-1.5d0*gm(3,3)*(dgm01(2,3)*dgm10(2,2)& 752 & +dgm01(2,2)*dgm10(2,3)+gm(2,3)*d2gm(2,2)+gm(2,2)*d2gm(2,3))+gm(2,3)& 753 & *(-1.5d0*dgm01(3,3)*dgm10(2,2)+15*dgm01(2,3)*dgm10(2,3)-1.5d0*dgm01(2,2)& 754 & *dgm10(3,3)-1.5d0*gm(2,2)*d2gm(3,3)) 755 cm(1,8)=-3*gm(1,3)*dgm01(2,3)*dgm10(1,1)-1.5d0*gm(1,2)*dgm01(3,3)& 756 & *dgm10(1,1)+15*gm(1,3)*dgm01(1,3)*dgm10(1,2)-1.5d0*gm(1,1)*dgm01(3,3)& 757 & *dgm10(1,2)+15*gm(1,3)*dgm01(1,2)*dgm10(1,3)+15*gm(1,2)*dgm01(1,3)& 758 & *dgm10(1,3)-3*gm(1,1)*dgm01(2,3)*dgm10(1,3)-3*gm(1,3)*dgm01(1,1)& 759 & *dgm10(2,3)-3*gm(1,1)*dgm01(1,3)*dgm10(2,3)-1.5d0*gm(1,2)*dgm01(1,1)& 760 & *dgm10(3,3)-1.5d0*gm(1,1)*dgm01(1,2)*dgm10(3,3)+7.5d0*gm(1,3)& 761 & **2*d2gm(1,2)-1.5d0*gm(3,3)*(dgm01(1,2)*dgm10(1,1)+dgm01(1,1)& 762 & *dgm10(1,2)+gm(1,2)*d2gm(1,1)+gm(1,1)*d2gm(1,2))+15*gm(1,2)*gm(1,3)& 763 & *d2gm(1,3)-3*gm(2,3)*(dgm01(1,3)*dgm10(1,1)+dgm01(1,1)*dgm10(1,3)& 764 & +gm(1,3)*d2gm(1,1)+gm(1,1)*d2gm(1,3))-3*gm(1,1)*gm(1,3)*d2gm(2,3)& 765 & -1.5d0*gm(1,1)*gm(1,2)*d2gm(3,3) 766 cm(2,8)=-4.5d0*gm(2,2)*dgm01(3,3)*dgm10(1,2)+12*gm(2,2)*dgm01(2,3)& 767 & *dgm10(1,3)+12*gm(1,3)*dgm01(2,3)*dgm10(2,2)-4.5d0*gm(1,2)*dgm01(3,3)& 768 & *dgm10(2,2)+12*gm(2,2)*dgm01(1,3)*dgm10(2,3)+12*gm(1,3)*dgm01(2,2)& 769 & *dgm10(2,3)+3*gm(1,2)*dgm01(2,3)*dgm10(2,3)-4.5d0*gm(2,2)*dgm01(1,2)& 770 & *dgm10(3,3)-4.5d0*gm(1,2)*dgm01(2,2)*dgm10(3,3)+1.5d0*gm(2,3)& 771 & **2*d2gm(1,2)-4.5d0*gm(3,3)*(dgm01(2,2)*dgm10(1,2)+dgm01(1,2)& 772 & *dgm10(2,2)+gm(2,2)*d2gm(1,2)+gm(1,2)*d2gm(2,2))+12*gm(1,3)*gm(2,2)& 773 & *d2gm(2,3)+gm(2,3)*(3*dgm01(2,3)*dgm10(1,2)+12*dgm01(2,2)*dgm10(1,3)& 774 & +12*dgm01(1,3)*dgm10(2,2)+3*dgm01(1,2)*dgm10(2,3)+12*gm(2,2)& 775 & *d2gm(1,3)+12*gm(1,3)*d2gm(2,2)+3*gm(1,2)*d2gm(2,3))-4.5d0*gm(1,2)& 776 & *gm(2,2)*d2gm(3,3) 777 cm(3,8)=3*gm(1,3)*dgm01(3,3)*dgm10(2,3)+3*gm(1,3)*dgm01(2,3)& 778 & *dgm10(3,3)+12*gm(1,2)*dgm01(3,3)*dgm10(3,3)+6*gm(3,3)**2*d2gm(1,2)& 779 & +gm(3,3)*(12*dgm01(3,3)*dgm10(1,2)+3*dgm01(2,3)*dgm10(1,3)+3*dgm01(1,3)& 780 & *dgm10(2,3)+12*dgm01(1,2)*dgm10(3,3)+3*gm(2,3)*d2gm(1,3)+3*gm(1,3)& 781 & *d2gm(2,3)+12*gm(1,2)*d2gm(3,3))+3*gm(2,3)*(dgm01(3,3)*dgm10(1,3)& 782 & +dgm01(1,3)*dgm10(3,3)+gm(1,3)*d2gm(3,3)) 783 cm(4,8)=12*gm(2,2)*dgm01(3,3)*dgm10(1,3)+12*gm(1,3)*dgm01(3,3)& 784 & *dgm10(2,2)+6*gm(1,3)*dgm01(2,3)*dgm10(2,3)+3*gm(1,2)*dgm01(3,3)& 785 & *dgm10(2,3)+12*gm(2,2)*dgm01(1,3)*dgm10(3,3)+12*gm(1,3)*dgm01(2,2)& 786 & *dgm10(3,3)+3*gm(1,2)*dgm01(2,3)*dgm10(3,3)+3*gm(2,3)**2*d2gm(1,3)& 787 & +gm(3,3)*(3*dgm01(2,3)*dgm10(1,2)+12*dgm01(2,2)*dgm10(1,3)+12*dgm01(1,3)& 788 & *dgm10(2,2)+3*dgm01(1,2)*dgm10(2,3)+3*gm(2,3)*d2gm(1,2)+12*gm(2,2)& 789 & *d2gm(1,3)+12*gm(1,3)*d2gm(2,2)+3*gm(1,2)*d2gm(2,3))+12*gm(1,3)& 790 & *gm(2,2)*d2gm(3,3)+gm(2,3)*(3*dgm01(3,3)*dgm10(1,2)+6*dgm01(2,3)& 791 & *dgm10(1,3)+6*dgm01(1,3)*dgm10(2,3)+3*dgm01(1,2)*dgm10(3,3)+6*gm(1,3)& 792 & *d2gm(2,3)+3*gm(1,2)*d2gm(3,3)) 793 cm(5,8)=12*gm(1,3)*dgm01(3,3)*dgm10(1,2)+3*gm(1,3)*dgm01(2,3)& 794 & *dgm10(1,3)+12*gm(1,2)*dgm01(3,3)*dgm10(1,3)+3*gm(1,3)*dgm01(1,3)& 795 & *dgm10(2,3)-4.5d0*gm(1,1)*dgm01(3,3)*dgm10(2,3)+12*gm(1,3)*dgm01(1,2)& 796 & *dgm10(3,3)+12*gm(1,2)*dgm01(1,3)*dgm10(3,3)-4.5d0*gm(1,1)*dgm01(2,3)& 797 & *dgm10(3,3)+1.5d0*gm(1,3)**2*d2gm(2,3)+gm(3,3)*(-4.5d0*dgm01(2,3)& 798 & *dgm10(1,1)+12*dgm01(1,3)*dgm10(1,2)+12*dgm01(1,2)*dgm10(1,3)& 799 & -4.5d0*dgm01(1,1)*dgm10(2,3)-4.5d0*gm(2,3)*d2gm(1,1)+12*gm(1,3)& 800 & *d2gm(1,2)+12*gm(1,2)*d2gm(1,3)-4.5d0*gm(1,1)*d2gm(2,3))+12*gm(1,2)& 801 & *gm(1,3)*d2gm(3,3)+gm(2,3)*(-4.5d0*dgm01(3,3)*dgm10(1,1)+3*dgm01(1,3)& 802 & *dgm10(1,3)-4.5d0*dgm01(1,1)*dgm10(3,3)+3*gm(1,3)*d2gm(1,3)-4.5d0*gm(1,1)& 803 & *d2gm(3,3)) 804 cm(6,8)=-1.5d0*gm(2,2)*dgm01(3,3)*dgm10(1,1)+9*gm(1,3)*dgm01(2,3)& 805 & *dgm10(1,2)-6*gm(1,2)*dgm01(3,3)*dgm10(1,2)+15*gm(2,2)*dgm01(1,3)& 806 & *dgm10(1,3)+15*gm(1,3)*dgm01(2,2)*dgm10(1,3)+9*gm(1,2)*dgm01(2,3)& 807 & *dgm10(1,3)+15*gm(1,3)*dgm01(1,3)*dgm10(2,2)-1.5d0*gm(1,1)*dgm01(3,3)& 808 & *dgm10(2,2)+9*gm(1,3)*dgm01(1,2)*dgm10(2,3)+9*gm(1,2)*dgm01(1,3)& 809 & *dgm10(2,3)-6*gm(1,1)*dgm01(2,3)*dgm10(2,3)-1.5d0*gm(2,2)*dgm01(1,1)& 810 & *dgm10(3,3)-6*gm(1,2)*dgm01(1,2)*dgm10(3,3)-1.5d0*gm(1,1)*dgm01(2,2)& 811 & *dgm10(3,3)-3*gm(2,3)**2*d2gm(1,1)+15*gm(1,3)*gm(2,2)*d2gm(1,3)& 812 & +7.5d0*gm(1,3)**2*d2gm(2,2)+gm(3,3)*(-1.5d0*dgm01(2,2)*dgm10(1,1)& 813 & -6*dgm01(1,2)*dgm10(1,2)-1.5d0*dgm01(1,1)*dgm10(2,2)-1.5d0*gm(2,2)& 814 & *d2gm(1,1)-6*gm(1,2)*d2gm(1,2)-1.5d0*gm(1,1)*d2gm(2,2))+9*gm(1,2)& 815 & *gm(1,3)*d2gm(2,3)+gm(2,3)*(-6*dgm01(2,3)*dgm10(1,1)+9*dgm01(1,3)& 816 & *dgm10(1,2)+9*dgm01(1,2)*dgm10(1,3)-6*dgm01(1,1)*dgm10(2,3)+9*gm(1,3)& 817 & *d2gm(1,2)+9*gm(1,2)*d2gm(1,3)-6*gm(1,1)*d2gm(2,3))-3*gm(1,2)& 818 & **2*d2gm(3,3)-1.5d0*gm(1,1)*gm(2,2)*d2gm(3,3) 819 cm(7,8)=4.5d0*gm(2,3)**2*d2gm(2,2)-3*gm(3,3)*(dgm01(2,2)*dgm10(2,2)& 820 & +gm(2,2)*d2gm(2,2))+9*gm(2,3)*(dgm01(2,3)*dgm10(2,2)+dgm01(2,2)& 821 & *dgm10(2,3)+gm(2,2)*d2gm(2,3))+gm(2,2)*(-3*dgm01(3,3)*dgm10(2,2)& 822 & +9*dgm01(2,3)*dgm10(2,3)-3*dgm01(2,2)*dgm10(3,3)-1.5d0*gm(2,2)& 823 & *d2gm(3,3)) 824 cm(8,8)=12*gm(2,2)*dgm01(3,3)*dgm10(3,3)+6*gm(2,3)*(dgm01(3,3)& 825 & *dgm10(2,3)+dgm01(2,3)*dgm10(3,3))+6*gm(3,3)**2*d2gm(2,2)+3*gm(2,3)& 826 & **2*d2gm(3,3)+gm(3,3)*(12*dgm01(3,3)*dgm10(2,2)+6*dgm01(2,3)& 827 & *dgm10(2,3)+12*dgm01(2,2)*dgm10(3,3)+6*gm(2,3)*d2gm(2,3)+12*gm(2,2)& 828 & *d2gm(3,3)) 829 cm(9,8)=7.5d0*gm(2,2)*(dgm01(3,3)*dgm10(2,3)+dgm01(2,3)*dgm10(3,3))& 830 & +4.5d0*gm(2,3)**2*d2gm(2,3)+7.5d0*gm(3,3)*(dgm01(2,3)*dgm10(2,2)& 831 & +dgm01(2,2)*dgm10(2,3)+gm(2,3)*d2gm(2,2)+gm(2,2)*d2gm(2,3))+gm(2,3)& 832 & *(7.5d0*dgm01(3,3)*dgm10(2,2)+9*dgm01(2,3)*dgm10(2,3)+7.5d0*dgm01(2,2)& 833 & *dgm10(3,3)+7.5d0*gm(2,2)*d2gm(3,3)) 834 cm(10,8)=6*gm(2,3)*dgm01(3,3)*dgm10(3,3)+3*gm(3,3)**2*d2gm(2,3)& 835 & +6*gm(3,3)*(dgm01(3,3)*dgm10(2,3)+dgm01(2,3)*dgm10(3,3)+gm(2,3)& 836 & *d2gm(3,3)) 837 cm(1,9)=(-18*gm(1,3)*dgm01(2,2)*dgm10(1,1)-36*gm(1,2)*dgm01(2,3)& 838 & *dgm10(1,1)+180*gm(1,3)*dgm01(1,2)*dgm10(1,2)+180*gm(1,2)*dgm01(1,3)& 839 & *dgm10(1,2)-36*gm(1,1)*dgm01(2,3)*dgm10(1,2)+180*gm(1,2)*dgm01(1,2)& 840 & *dgm10(1,3)-18*gm(1,1)*dgm01(2,2)*dgm10(1,3)-18*gm(1,3)*dgm01(1,1)& 841 & *dgm10(2,2)-18*gm(1,1)*dgm01(1,3)*dgm10(2,2)-36*gm(1,2)*dgm01(1,1)& 842 & *dgm10(2,3)-36*gm(1,1)*dgm01(1,2)*dgm10(2,3)+180*gm(1,2)*gm(1,3)& 843 & *d2gm(1,2)-36*gm(2,3)*(dgm01(1,2)*dgm10(1,1)+dgm01(1,1)*dgm10(1,2)& 844 & +gm(1,2)*d2gm(1,1)+gm(1,1)*d2gm(1,2))+90*gm(1,2)**2*d2gm(1,3)& 845 & -18*gm(2,2)*(dgm01(1,3)*dgm10(1,1)+dgm01(1,1)*dgm10(1,3)+gm(1,3)& 846 & *d2gm(1,1)+gm(1,1)*d2gm(1,3))-18*gm(1,1)*gm(1,3)*d2gm(2,2)-36*gm(1,1)& 847 & *gm(1,2)*d2gm(2,3))/12.d0 848 cm(2,9)=12*gm(1,3)*dgm01(2,2)*dgm10(2,2)+3*gm(1,2)*dgm01(2,3)& 849 & *dgm10(2,2)+3*gm(1,2)*dgm01(2,2)*dgm10(2,3)+6*gm(2,2)**2*d2gm(1,3)& 850 & +3*gm(2,3)*(dgm01(2,2)*dgm10(1,2)+dgm01(1,2)*dgm10(2,2)+gm(2,2)& 851 & *d2gm(1,2)+gm(1,2)*d2gm(2,2))+gm(2,2)*(3*dgm01(2,3)*dgm10(1,2)& 852 & +12*dgm01(2,2)*dgm10(1,3)+12*dgm01(1,3)*dgm10(2,2)+3*dgm01(1,2)& 853 & *dgm10(2,3)+12*gm(1,3)*d2gm(2,2)+3*gm(1,2)*d2gm(2,3)) 854 cm(3,9)=-4.5d0*gm(2,2)*dgm01(3,3)*dgm10(1,3)-4.5d0*gm(1,3)*dgm01(3,3)& 855 & *dgm10(2,2)+3*gm(1,3)*dgm01(2,3)*dgm10(2,3)+12*gm(1,2)*dgm01(3,3)& 856 & *dgm10(2,3)-4.5d0*gm(2,2)*dgm01(1,3)*dgm10(3,3)-4.5d0*gm(1,3)& 857 & *dgm01(2,2)*dgm10(3,3)+12*gm(1,2)*dgm01(2,3)*dgm10(3,3)+1.5d0*gm(2,3)& 858 & **2*d2gm(1,3)+gm(3,3)*(12*dgm01(2,3)*dgm10(1,2)-4.5d0*dgm01(2,2)& 859 & *dgm10(1,3)-4.5d0*dgm01(1,3)*dgm10(2,2)+12*dgm01(1,2)*dgm10(2,3)& 860 & +12*gm(2,3)*d2gm(1,2)-4.5d0*gm(2,2)*d2gm(1,3)-4.5d0*gm(1,3)*d2gm(2,2)& 861 & +12*gm(1,2)*d2gm(2,3))-4.5d0*gm(1,3)*gm(2,2)*d2gm(3,3)+gm(2,3)& 862 & *(12*dgm01(3,3)*dgm10(1,2)+3*dgm01(2,3)*dgm10(1,3)+3*dgm01(1,3)& 863 & *dgm10(2,3)+12*dgm01(1,2)*dgm10(3,3)+3*gm(1,3)*d2gm(2,3)+12*gm(1,2)& 864 & *d2gm(3,3)) 865 cm(4,9)=12*gm(2,2)*dgm01(3,3)*dgm10(1,2)+3*gm(2,2)*dgm01(2,3)& 866 & *dgm10(1,3)+3*gm(1,3)*dgm01(2,3)*dgm10(2,2)+12*gm(1,2)*dgm01(3,3)& 867 & *dgm10(2,2)+3*gm(2,2)*dgm01(1,3)*dgm10(2,3)+3*gm(1,3)*dgm01(2,2)& 868 & *dgm10(2,3)+6*gm(1,2)*dgm01(2,3)*dgm10(2,3)+12*gm(2,2)*dgm01(1,2)& 869 & *dgm10(3,3)+12*gm(1,2)*dgm01(2,2)*dgm10(3,3)+3*gm(2,3)**2*d2gm(1,2)& 870 & +12*gm(3,3)*(dgm01(2,2)*dgm10(1,2)+dgm01(1,2)*dgm10(2,2)+gm(2,2)& 871 & *d2gm(1,2)+gm(1,2)*d2gm(2,2))+3*gm(1,3)*gm(2,2)*d2gm(2,3)+gm(2,3)& 872 & *(6*dgm01(2,3)*dgm10(1,2)+3*dgm01(2,2)*dgm10(1,3)+3*dgm01(1,3)& 873 & *dgm10(2,2)+6*dgm01(1,2)*dgm10(2,3)+3*gm(2,2)*d2gm(1,3)+3*gm(1,3)& 874 & *d2gm(2,2)+6*gm(1,2)*d2gm(2,3))+12*gm(1,2)*gm(2,2)*d2gm(3,3) 875 cm(5,9)=-1.5d0*gm(2,2)*dgm01(3,3)*dgm10(1,1)+9*gm(1,3)*dgm01(2,3)& 876 & *dgm10(1,2)+15*gm(1,2)*dgm01(3,3)*dgm10(1,2)-6*gm(2,2)*dgm01(1,3)& 877 & *dgm10(1,3)-6*gm(1,3)*dgm01(2,2)*dgm10(1,3)+9*gm(1,2)*dgm01(2,3)& 878 & *dgm10(1,3)-6*gm(1,3)*dgm01(1,3)*dgm10(2,2)-1.5d0*gm(1,1)*dgm01(3,3)& 879 & *dgm10(2,2)+9*gm(1,3)*dgm01(1,2)*dgm10(2,3)+9*gm(1,2)*dgm01(1,3)& 880 & *dgm10(2,3)-6*gm(1,1)*dgm01(2,3)*dgm10(2,3)-1.5d0*gm(2,2)*dgm01(1,1)& 881 & *dgm10(3,3)+15*gm(1,2)*dgm01(1,2)*dgm10(3,3)-1.5d0*gm(1,1)*dgm01(2,2)& 882 & *dgm10(3,3)-3*gm(2,3)**2*d2gm(1,1)-6*gm(1,3)*gm(2,2)*d2gm(1,3)& 883 & -3*gm(1,3)**2*d2gm(2,2)+gm(3,3)*(-1.5d0*dgm01(2,2)*dgm10(1,1)& 884 & +15*dgm01(1,2)*dgm10(1,2)-1.5d0*dgm01(1,1)*dgm10(2,2)-1.5d0*gm(2,2)& 885 & *d2gm(1,1)+15*gm(1,2)*d2gm(1,2)-1.5d0*gm(1,1)*d2gm(2,2))+9*gm(1,2)& 886 & *gm(1,3)*d2gm(2,3)+gm(2,3)*(-6*dgm01(2,3)*dgm10(1,1)+9*dgm01(1,3)& 887 & *dgm10(1,2)+9*dgm01(1,2)*dgm10(1,3)-6*dgm01(1,1)*dgm10(2,3)+9*gm(1,3)& 888 & *d2gm(1,2)+9*gm(1,2)*d2gm(1,3)-6*gm(1,1)*d2gm(2,3))+7.5d0*gm(1,2)& 889 & **2*d2gm(3,3)-1.5d0*gm(1,1)*gm(2,2)*d2gm(3,3) 890 cm(6,9)=12*gm(1,3)*dgm01(2,2)*dgm10(1,2)+3*gm(1,2)*dgm01(2,3)& 891 & *dgm10(1,2)+12*gm(1,2)*dgm01(2,2)*dgm10(1,3)+12*gm(1,3)*dgm01(1,2)& 892 & *dgm10(2,2)+12*gm(1,2)*dgm01(1,3)*dgm10(2,2)-4.5d0*gm(1,1)*dgm01(2,3)& 893 & *dgm10(2,2)+3*gm(1,2)*dgm01(1,2)*dgm10(2,3)-4.5d0*gm(1,1)*dgm01(2,2)& 894 & *dgm10(2,3)+12*gm(1,2)*gm(1,3)*d2gm(2,2)+gm(2,3)*(-4.5d0*dgm01(2,2)& 895 & *dgm10(1,1)+3*dgm01(1,2)*dgm10(1,2)-4.5d0*dgm01(1,1)*dgm10(2,2)& 896 & -4.5d0*gm(2,2)*d2gm(1,1)+3*gm(1,2)*d2gm(1,2)-4.5d0*gm(1,1)*d2gm(2,2))& 897 & +1.5d0*gm(1,2)**2*d2gm(2,3)+gm(2,2)*(-4.5d0*dgm01(2,3)*dgm10(1,1)& 898 & +12*dgm01(1,3)*dgm10(1,2)+12*dgm01(1,2)*dgm10(1,3)-4.5d0*dgm01(1,1)& 899 & *dgm10(2,3)+12*gm(1,3)*d2gm(1,2)+12*gm(1,2)*d2gm(1,3)-4.5d0*gm(1,1)& 900 & *d2gm(2,3)) 901 cm(7,9)=6*gm(2,3)*(dgm01(2,2)*dgm10(2,2)+gm(2,2)*d2gm(2,2))& 902 & +gm(2,2)*(6*dgm01(2,3)*dgm10(2,2)+6*dgm01(2,2)*dgm10(2,3)+3*gm(2,2)& 903 & *d2gm(2,3)) 904 cm(8,9)=7.5d0*gm(2,2)*(dgm01(3,3)*dgm10(2,3)+dgm01(2,3)*dgm10(3,3))& 905 & +4.5d0*gm(2,3)**2*d2gm(2,3)+7.5d0*gm(3,3)*(dgm01(2,3)*dgm10(2,2)& 906 & +dgm01(2,2)*dgm10(2,3)+gm(2,3)*d2gm(2,2)+gm(2,2)*d2gm(2,3))+gm(2,3)& 907 & *(7.5d0*dgm01(3,3)*dgm10(2,2)+9*dgm01(2,3)*dgm10(2,3)+7.5d0*dgm01(2,2)& 908 & *dgm10(3,3)+7.5d0*gm(2,2)*d2gm(3,3)) 909 cm(9,9)=3*gm(2,3)**2*d2gm(2,2)+12*gm(3,3)*(dgm01(2,2)*dgm10(2,2)& 910 & +gm(2,2)*d2gm(2,2))+6*gm(2,3)*(dgm01(2,3)*dgm10(2,2)+dgm01(2,2)& 911 & *dgm10(2,3)+gm(2,2)*d2gm(2,3))+gm(2,2)*(12*dgm01(3,3)*dgm10(2,2)& 912 & +6*dgm01(2,3)*dgm10(2,3)+12*dgm01(2,2)*dgm10(3,3)+6*gm(2,2)*d2gm(3,3)) 913 cm(10,9)=-3*gm(2,2)*dgm01(3,3)*dgm10(3,3)+9*gm(2,3)*(dgm01(3,3)& 914 & *dgm10(2,3)+dgm01(2,3)*dgm10(3,3))-1.5d0*gm(3,3)**2*d2gm(2,2)& 915 & +4.5d0*gm(2,3)**2*d2gm(3,3)+gm(3,3)*(-3*dgm01(3,3)*dgm10(2,2)& 916 & +9*dgm01(2,3)*dgm10(2,3)-3*dgm01(2,2)*dgm10(3,3)+9*gm(2,3)*d2gm(2,3)& 917 & -3*gm(2,2)*d2gm(3,3)) 918 cm(1,10)=-1.5d0*gm(1,1)*(dgm01(3,3)*dgm10(1,3)+dgm01(1,3)*dgm10(3,3))& 919 & +7.5d0*gm(1,3)**2*d2gm(1,3)-1.5d0*gm(3,3)*(dgm01(1,3)*dgm10(1,1)& 920 & +dgm01(1,1)*dgm10(1,3)+gm(1,3)*d2gm(1,1)+gm(1,1)*d2gm(1,3))+gm(1,3)& 921 & *(-1.5d0*dgm01(3,3)*dgm10(1,1)+15*dgm01(1,3)*dgm10(1,3)-1.5d0*dgm01(1,1)& 922 & *dgm10(3,3)-1.5d0*gm(1,1)*d2gm(3,3)) 923 cm(2,10)=-1.5d0*gm(2,2)*dgm01(3,3)*dgm10(1,3)-1.5d0*gm(1,3)& 924 & *dgm01(3,3)*dgm10(2,2)+15*gm(1,3)*dgm01(2,3)*dgm10(2,3)-3*gm(1,2)& 925 & *dgm01(3,3)*dgm10(2,3)-1.5d0*gm(2,2)*dgm01(1,3)*dgm10(3,3)-1.5d0*gm(1,3)& 926 & *dgm01(2,2)*dgm10(3,3)-3*gm(1,2)*dgm01(2,3)*dgm10(3,3)+7.5d0*gm(2,3)& 927 & **2*d2gm(1,3)+gm(3,3)*(-3*dgm01(2,3)*dgm10(1,2)-1.5d0*dgm01(2,2)& 928 & *dgm10(1,3)-1.5d0*dgm01(1,3)*dgm10(2,2)-3*dgm01(1,2)*dgm10(2,3)& 929 & -3*gm(2,3)*d2gm(1,2)-1.5d0*gm(2,2)*d2gm(1,3)-1.5d0*gm(1,3)*d2gm(2,2)& 930 & -3*gm(1,2)*d2gm(2,3))-1.5d0*gm(1,3)*gm(2,2)*d2gm(3,3)+gm(2,3)& 931 & *(-3*dgm01(3,3)*dgm10(1,2)+15*dgm01(2,3)*dgm10(1,3)+15*dgm01(1,3)& 932 & *dgm10(2,3)-3*dgm01(1,2)*dgm10(3,3)+15*gm(1,3)*d2gm(2,3)-3*gm(1,2)& 933 & *d2gm(3,3)) 934 cm(3,10)=6*gm(1,3)*dgm01(3,3)*dgm10(3,3)+3*gm(3,3)**2*d2gm(1,3)& 935 & +6*gm(3,3)*(dgm01(3,3)*dgm10(1,3)+dgm01(1,3)*dgm10(3,3)+gm(1,3)& 936 & *d2gm(3,3)) 937 cm(4,10)=9*gm(1,3)*dgm01(3,3)*dgm10(2,3)+9*gm(1,3)*dgm01(2,3)& 938 & *dgm10(3,3)-6*gm(1,2)*dgm01(3,3)*dgm10(3,3)-3*gm(3,3)**2*d2gm(1,2)& 939 & +gm(3,3)*(-6*dgm01(3,3)*dgm10(1,2)+9*dgm01(2,3)*dgm10(1,3)+9*dgm01(1,3)& 940 & *dgm10(2,3)-6*dgm01(1,2)*dgm10(3,3)+9*gm(2,3)*d2gm(1,3)+9*gm(1,3)& 941 & *d2gm(2,3)-6*gm(1,2)*d2gm(3,3))+9*gm(2,3)*(dgm01(3,3)*dgm10(1,3)& 942 & +dgm01(1,3)*dgm10(3,3)+gm(1,3)*d2gm(3,3)) 943 cm(5,10)=-3*gm(1,1)*dgm01(3,3)*dgm10(3,3)+9*gm(1,3)*(dgm01(3,3)& 944 & *dgm10(1,3)+dgm01(1,3)*dgm10(3,3))-1.5d0*gm(3,3)**2*d2gm(1,1)& 945 & +4.5d0*gm(1,3)**2*d2gm(3,3)+gm(3,3)*(-3*dgm01(3,3)*dgm10(1,1)& 946 & +9*dgm01(1,3)*dgm10(1,3)-3*dgm01(1,1)*dgm10(3,3)+9*gm(1,3)*d2gm(1,3)& 947 & -3*gm(1,1)*d2gm(3,3)) 948 cm(6,10)=-3*gm(1,3)*dgm01(3,3)*dgm10(1,2)+15*gm(1,3)*dgm01(2,3)& 949 & *dgm10(1,3)-3*gm(1,2)*dgm01(3,3)*dgm10(1,3)+15*gm(1,3)*dgm01(1,3)& 950 & *dgm10(2,3)-1.5d0*gm(1,1)*dgm01(3,3)*dgm10(2,3)-3*gm(1,3)*dgm01(1,2)& 951 & *dgm10(3,3)-3*gm(1,2)*dgm01(1,3)*dgm10(3,3)-1.5d0*gm(1,1)*dgm01(2,3)& 952 & *dgm10(3,3)+7.5d0*gm(1,3)**2*d2gm(2,3)+gm(3,3)*(-1.5d0*dgm01(2,3)& 953 & *dgm10(1,1)-3*dgm01(1,3)*dgm10(1,2)-3*dgm01(1,2)*dgm10(1,3)-1.5d0*dgm01(1,1)& 954 & *dgm10(2,3)-1.5d0*gm(2,3)*d2gm(1,1)-3*gm(1,3)*d2gm(1,2)-3*gm(1,2)& 955 & *d2gm(1,3)-1.5d0*gm(1,1)*d2gm(2,3))-3*gm(1,2)*gm(1,3)*d2gm(3,3)& 956 & +gm(2,3)*(-1.5d0*dgm01(3,3)*dgm10(1,1)+15*dgm01(1,3)*dgm10(1,3)& 957 & -1.5d0*dgm01(1,1)*dgm10(3,3)+15*gm(1,3)*d2gm(1,3)-1.5d0*gm(1,1)& 958 & *d2gm(3,3)) 959 cm(7,10)=-1.5d0*gm(2,2)*(dgm01(3,3)*dgm10(2,3)+dgm01(2,3)*dgm10(3,3))& 960 & +7.5d0*gm(2,3)**2*d2gm(2,3)-1.5d0*gm(3,3)*(dgm01(2,3)*dgm10(2,2)& 961 & +dgm01(2,2)*dgm10(2,3)+gm(2,3)*d2gm(2,2)+gm(2,2)*d2gm(2,3))+gm(2,3)& 962 & *(-1.5d0*dgm01(3,3)*dgm10(2,2)+15*dgm01(2,3)*dgm10(2,3)-1.5d0*dgm01(2,2)& 963 & *dgm10(3,3)-1.5d0*gm(2,2)*d2gm(3,3)) 964 cm(8,10)=6*gm(2,3)*dgm01(3,3)*dgm10(3,3)+3*gm(3,3)**2*d2gm(2,3)& 965 & +6*gm(3,3)*(dgm01(3,3)*dgm10(2,3)+dgm01(2,3)*dgm10(3,3)+gm(2,3)& 966 & *d2gm(3,3)) 967 cm(9,10)=-3*gm(2,2)*dgm01(3,3)*dgm10(3,3)+9*gm(2,3)*(dgm01(3,3)& 968 & *dgm10(2,3)+dgm01(2,3)*dgm10(3,3))-1.5d0*gm(3,3)**2*d2gm(2,2)& 969 & +4.5d0*gm(2,3)**2*d2gm(3,3)+gm(3,3)*(-3*dgm01(3,3)*dgm10(2,2)& 970 & +9*dgm01(2,3)*dgm10(2,3)-3*dgm01(2,2)*dgm10(3,3)+9*gm(2,3)*d2gm(2,3)& 971 & -3*gm(2,2)*d2gm(3,3)) 972 cm(10,10)=gm(3,3)*(6*dgm01(3,3)*dgm10(3,3)+3*gm(3,3)*d2gm(3,3)) 973 974 975 end if 976 ! 977 !contraction to output scalar 978 ! 979 e2nl=0.d0 980 do jj=1,((rank+1)*(rank+2))/2 981 tmp(:)=0.d0 982 do ii=1,((rank+1)*(rank+2))/2 983 tmp(:)=tmp(:)+aa(:,ii)*cm(ii,jj) 984 end do 985 e2nl=e2nl+tmp(1)*bb(1,jj)+tmp(2)*bb(2,jj) 986 end do 987 988 ABI_DEALLOCATE(cm) 989 990 end subroutine contstr26