TABLE OF CONTENTS


ABINIT/contstr26 [ Functions ]

[ Top ] [ 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-2024 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.

SOURCE

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