1*59599516SKenneth E. Jansen module lhsGkeep 2*59599516SKenneth E. Jansen 3*59599516SKenneth E. Jansen real*8, allocatable :: lhsG(:) 4*59599516SKenneth E. Jansen 5*59599516SKenneth E. Jansen end module 6*59599516SKenneth E. Jansen 7*59599516SKenneth E. Jansenc---------------------------------------------------------------------------- 8*59599516SKenneth E. Jansen 9*59599516SKenneth E. Jansen subroutine keeplhsG 10*59599516SKenneth E. Jansen 11*59599516SKenneth E. Jansen use lhsGkeep 12*59599516SKenneth E. Jansen 13*59599516SKenneth E. Jansen include "common.h" 14*59599516SKenneth E. Jansen 15*59599516SKenneth E. Jansen allocate ( lhsG(nnz*nshg) ) 16*59599516SKenneth E. Jansen 17*59599516SKenneth E. Jansen return 18*59599516SKenneth E. Jansen 19*59599516SKenneth E. Jansen end 20*59599516SKenneth E. Jansen subroutine cmass (shp, shgl, xl, em) 21*59599516SKenneth E. Jansenc 22*59599516SKenneth E. Jansenc---------------------------------------------------------------------- 23*59599516SKenneth E. Jansenc 24*59599516SKenneth E. Jansenc This subroutine computes the consistent mass matrices 25*59599516SKenneth E. Jansenc 26*59599516SKenneth E. Jansenc Ken Jansen, Spring 2000 27*59599516SKenneth E. Jansenc---------------------------------------------------------------------- 28*59599516SKenneth E. Jansenc 29*59599516SKenneth E. Jansenc 30*59599516SKenneth E. Jansen include "common.h" 31*59599516SKenneth E. Jansenc 32*59599516SKenneth E. Jansen integer ne, na, nb, nodlcla, nodlclb, iel 33*59599516SKenneth E. Jansen dimension 34*59599516SKenneth E. Jansen & shp(nshl,MAXQPT), shgl(nsd,nshl,MAXQPT), 35*59599516SKenneth E. Jansen & em(npro,nshl,nshl), 36*59599516SKenneth E. Jansen & xl(npro,nenl,nsd) 37*59599516SKenneth E. Jansenc 38*59599516SKenneth E. Jansen dimension shape(npro,nshl), shdrv(npro,nsd,nshl), 39*59599516SKenneth E. Jansen & sgn(npro,nshl), dxidx(npro,nsd,nsd), 40*59599516SKenneth E. Jansen & shg(npro,nshl,nsd), 41*59599516SKenneth E. Jansen & WdetJ(npro) 42*59599516SKenneth E. Jansenc 43*59599516SKenneth E. Jansen em = zero 44*59599516SKenneth E. Jansenc 45*59599516SKenneth E. Jansenc.... loop through the integration points 46*59599516SKenneth E. Jansenc 47*59599516SKenneth E. Jansen do intp = 1, ngauss ! (these are in common.h) 48*59599516SKenneth E. Jansenc 49*59599516SKenneth E. Jansenc.... get the hierarchic shape functions at this int point 50*59599516SKenneth E. Jansenc 51*59599516SKenneth E. Jansen call getshp(shp, shgl, sgn, 52*59599516SKenneth E. Jansen & shape, shdrv, intp) 53*59599516SKenneth E. Jansenc 54*59599516SKenneth E. Jansenc.... calculate the determinant of the jacobian and weight it 55*59599516SKenneth E. Jansenc 56*59599516SKenneth E. Jansen call e3metric( xl, shdrv,dxidx,shg,WdetJ) 57*59599516SKenneth E. Jansenc 58*59599516SKenneth E. Jansen do iel = 1, npro 59*59599516SKenneth E. Jansen do na = 1, nshl 60*59599516SKenneth E. Jansen do nb = 1, nshl 61*59599516SKenneth E. Jansen shp2 = shape(iel,na) * shape(iel,nb) 62*59599516SKenneth E. Jansen em(iel,na,nb) = em(iel,na,nb) + shp2*WdetJ(iel) 63*59599516SKenneth E. Jansen enddo 64*59599516SKenneth E. Jansen enddo 65*59599516SKenneth E. Jansen enddo 66*59599516SKenneth E. Jansen enddo 67*59599516SKenneth E. Jansenc 68*59599516SKenneth E. Jansenc.... return 69*59599516SKenneth E. Jansenc 70*59599516SKenneth E. Jansen return 71*59599516SKenneth E. Jansen end 72*59599516SKenneth E. Jansen 73*59599516SKenneth E. Jansen subroutine cmassl (shp, shgl, shpf, shglf, xl, em, Qwtf) 74*59599516SKenneth E. Jansenc 75*59599516SKenneth E. Jansenc---------------------------------------------------------------------- 76*59599516SKenneth E. Jansenc 77*59599516SKenneth E. Jansenc This subroutine computes the consistent mass matrices 78*59599516SKenneth E. Jansenc 79*59599516SKenneth E. Jansenc Ken Jansen, Spring 2000 80*59599516SKenneth E. Jansenc---------------------------------------------------------------------- 81*59599516SKenneth E. Jansenc 82*59599516SKenneth E. Jansenc 83*59599516SKenneth E. Jansen include "common.h" 84*59599516SKenneth E. Jansenc 85*59599516SKenneth E. Jansen integer ne, na, nb, nodlcla, nodlclb, iel 86*59599516SKenneth E. Jansen dimension 87*59599516SKenneth E. Jansen & shp(nshl,MAXQPT), shgl(nsd,nshl,MAXQPT), 88*59599516SKenneth E. Jansen & shpf(nshl,MAXQPT), shglf(nsd,nshl,MAXQPT), 89*59599516SKenneth E. Jansen & em(npro,nshl,nshl), eml(npro,nshl), 90*59599516SKenneth E. Jansen & xl(npro,nenl,nsd) 91*59599516SKenneth E. Jansenc 92*59599516SKenneth E. Jansen dimension shape(npro,nshl), shdrv(npro,nsd,nshl), 93*59599516SKenneth E. Jansen & sgn(npro,nshl), dxidx(npro,nsd,nsd), 94*59599516SKenneth E. Jansen & shg(npro,nshl,nsd), Qwtf(ngaussf), 95*59599516SKenneth E. Jansen & WdetJ(npro) 96*59599516SKenneth E. Jansenc 97*59599516SKenneth E. Jansen em = zero 98*59599516SKenneth E. Jansen eml= zero 99*59599516SKenneth E. Jansen 100*59599516SKenneth E. Jansen if (ifproj.eq.1)then 101*59599516SKenneth E. Jansen nods = nshl 102*59599516SKenneth E. Jansen else 103*59599516SKenneth E. Jansen nods = nenl 104*59599516SKenneth E. Jansen endif 105*59599516SKenneth E. Jansen 106*59599516SKenneth E. Jansenc----------------> Get the lumped mass matrix <----------------------- 107*59599516SKenneth E. Jansen 108*59599516SKenneth E. Jansenc 109*59599516SKenneth E. Jansenc.... loop through the integration points 110*59599516SKenneth E. Jansenc 111*59599516SKenneth E. Jansen do intp = 1, ngaussf ! (these are in common.h) 112*59599516SKenneth E. Jansenc 113*59599516SKenneth E. Jansenc.... get the hierarchic shape functions at this int point 114*59599516SKenneth E. Jansenc 115*59599516SKenneth E. Jansen call getshp(shpf, shglf, sgn, 116*59599516SKenneth E. Jansen & shape, shdrv, intp) 117*59599516SKenneth E. Jansenc 118*59599516SKenneth E. Jansenc.... calculate the determinant of the jacobian and weight it 119*59599516SKenneth E. Jansenc 120*59599516SKenneth E. Jansen call e3metricf( xl, shdrv,dxidx,shg,WdetJ,Qwtf) 121*59599516SKenneth E. Jansenc 122*59599516SKenneth E. Jansen do i=1,nods !nenl !nshl 123*59599516SKenneth E. Jansen eml(:,i) = eml(:,i) + shape(:,i)*WdetJ(:) 124*59599516SKenneth E. Jansen enddo 125*59599516SKenneth E. Jansen 126*59599516SKenneth E. Jansen enddo ! End loop over quad points. 127*59599516SKenneth E. Jansen 128*59599516SKenneth E. Jansen 129*59599516SKenneth E. Jansenc--------------> Get the consistent mass matrix <------------------------ 130*59599516SKenneth E. Jansen 131*59599516SKenneth E. Jansen 132*59599516SKenneth E. Jansen shape = zero 133*59599516SKenneth E. Jansen shdrv = zero 134*59599516SKenneth E. Jansen dxidx = zero 135*59599516SKenneth E. Jansen WdetJ = zero 136*59599516SKenneth E. Jansen shg = zero 137*59599516SKenneth E. Jansen 138*59599516SKenneth E. Jansenc 139*59599516SKenneth E. Jansenc.... loop through the integration points 140*59599516SKenneth E. Jansenc 141*59599516SKenneth E. Jansen do intp = 1, ngauss ! (these are in common.h) 142*59599516SKenneth E. Jansen 143*59599516SKenneth E. Jansenc.... get the hierarchic shape functions at this int point 144*59599516SKenneth E. Jansenc 145*59599516SKenneth E. Jansenc 146*59599516SKenneth E. Jansenc.... for the mass matrix to be consistent shp and shgl must be 147*59599516SKenneth E. Jansenc.... evaluated with at least higher quadrature than one-pt. quad. 148*59599516SKenneth E. Jansen 149*59599516SKenneth E. Jansen call getshp(shp, shgl, sgn, 150*59599516SKenneth E. Jansen & shape, shdrv, intp) 151*59599516SKenneth E. Jansen 152*59599516SKenneth E. Jansenc 153*59599516SKenneth E. Jansenc.... calculate the determinant of the jacobian and weight it 154*59599516SKenneth E. Jansenc 155*59599516SKenneth E. Jansen call e3metric( xl, shdrv,dxidx,shg,WdetJ) 156*59599516SKenneth E. Jansenc 157*59599516SKenneth E. Jansen 158*59599516SKenneth E. Jansen do iel = 1, npro 159*59599516SKenneth E. Jansen do na = 1, nods !nenl !nshl 160*59599516SKenneth E. Jansen do nb = 1, nods !nenl !nshl 161*59599516SKenneth E. Jansen shp2 = shape(iel,na) * shape(iel,nb) 162*59599516SKenneth E. Jansen em(iel,na,nb) = em(iel,na,nb) + shp2*WdetJ(iel) 163*59599516SKenneth E. Jansen enddo 164*59599516SKenneth E. Jansen enddo 165*59599516SKenneth E. Jansen enddo 166*59599516SKenneth E. Jansen 167*59599516SKenneth E. Jansen enddo ! End loop over quadrature points 168*59599516SKenneth E. Jansen 169*59599516SKenneth E. Jansen 170*59599516SKenneth E. Jansen 171*59599516SKenneth E. Jansenc----------> Obtain a mixed (lumped/consistent) mass matrix <------------ 172*59599516SKenneth E. Jansen 173*59599516SKenneth E. Jansenc... Different combinations of the lump and mass matrices yield 174*59599516SKenneth E. Jansenc... filters of varying widths. In the limiting case were 175*59599516SKenneth E. Jansenc... the entire matrix is lumped, we obtain the same filter as 176*59599516SKenneth E. Jansenc... in getdmc.f. Note that in these equivalent ways of 177*59599516SKenneth E. Jansenc... filtering one-point quadrature is used for shpf and shgl.. 178*59599516SKenneth E. Jansen 179*59599516SKenneth E. Jansen 180*59599516SKenneth E. Jansen em = (one-flump)*em 181*59599516SKenneth E. Jansen 182*59599516SKenneth E. Jansen do iel = 1, npro 183*59599516SKenneth E. Jansen do na = 1, nods !nenl !nshl 184*59599516SKenneth E. Jansen em(iel,na,na) = em(iel,na,na)+flump*eml(iel,na) 185*59599516SKenneth E. Jansen enddo 186*59599516SKenneth E. Jansen enddo 187*59599516SKenneth E. Jansen 188*59599516SKenneth E. Jansenc 189*59599516SKenneth E. Jansenc.... return 190*59599516SKenneth E. Jansenc 191*59599516SKenneth E. Jansen return 192*59599516SKenneth E. Jansen end 193*59599516SKenneth E. Jansen 194*59599516SKenneth E. Jansen 195*59599516SKenneth E. Jansen 196*59599516SKenneth E. Jansen subroutine cmasstl (shp, shgl, shpf, shglf, xl, em, Qwtf) 197*59599516SKenneth E. Jansenc 198*59599516SKenneth E. Jansenc---------------------------------------------------------------------- 199*59599516SKenneth E. Jansenc 200*59599516SKenneth E. Jansenc This subroutine computes the consistent mass matrices 201*59599516SKenneth E. Jansenc 202*59599516SKenneth E. Jansenc Ken Jansen, Spring 2000 203*59599516SKenneth E. Jansenc---------------------------------------------------------------------- 204*59599516SKenneth E. Jansenc 205*59599516SKenneth E. Jansenc 206*59599516SKenneth E. Jansen include "common.h" 207*59599516SKenneth E. Jansenc 208*59599516SKenneth E. Jansen integer ne, na, nb, nodlcla, nodlclb, iel 209*59599516SKenneth E. Jansen dimension 210*59599516SKenneth E. Jansen & shp(nshl,MAXQPT), shgl(nsd,nshl,MAXQPT), 211*59599516SKenneth E. Jansen & shpf(nshl,MAXQPT), shglf(nsd,nshl,MAXQPT), 212*59599516SKenneth E. Jansen & em(npro,nshl,nshl), eml(npro,nshl), 213*59599516SKenneth E. Jansen & xl(npro,nenl,nsd) 214*59599516SKenneth E. Jansenc 215*59599516SKenneth E. Jansen dimension shape(npro,nshl), shdrv(npro,nsd,nshl), 216*59599516SKenneth E. Jansen & sgn(npro,nshl), dxidx(npro,nsd,nsd), 217*59599516SKenneth E. Jansen & shg(npro,nshl,nsd), Qwtf(ngaussf), 218*59599516SKenneth E. Jansen & WdetJ(npro) 219*59599516SKenneth E. Jansenc 220*59599516SKenneth E. Jansen em = zero 221*59599516SKenneth E. Jansen eml= zero 222*59599516SKenneth E. Jansen 223*59599516SKenneth E. Jansenc----------------> Get the lumped mass matrix <----------------------- 224*59599516SKenneth E. Jansen 225*59599516SKenneth E. Jansenc 226*59599516SKenneth E. Jansenc.... loop through the integration points 227*59599516SKenneth E. Jansenc 228*59599516SKenneth E. Jansen do intp = 1, ngaussf ! (these are in common.h) 229*59599516SKenneth E. Jansenc 230*59599516SKenneth E. Jansenc.... get the hierarchic shape functions at this int point 231*59599516SKenneth E. Jansenc 232*59599516SKenneth E. Jansen call getshp(shpf, shglf, sgn, 233*59599516SKenneth E. Jansen & shape, shdrv, intp) 234*59599516SKenneth E. Jansenc 235*59599516SKenneth E. Jansenc.... calculate the determinant of the jacobian and weight it 236*59599516SKenneth E. Jansenc 237*59599516SKenneth E. Jansen call e3metricf( xl, shdrv,dxidx,shg,WdetJ,Qwtf) 238*59599516SKenneth E. Jansenc 239*59599516SKenneth E. Jansen do i=1,nshl 240*59599516SKenneth E. Jansen eml(:,i) = eml(:,i) + shape(:,i)*WdetJ(:) 241*59599516SKenneth E. Jansen enddo 242*59599516SKenneth E. Jansen 243*59599516SKenneth E. Jansen enddo ! End loop over quad points. 244*59599516SKenneth E. Jansen 245*59599516SKenneth E. Jansen 246*59599516SKenneth E. Jansenc--------------> Get the consistent mass matrix <------------------------ 247*59599516SKenneth E. Jansen 248*59599516SKenneth E. Jansen 249*59599516SKenneth E. Jansen shape = zero 250*59599516SKenneth E. Jansen shdrv = zero 251*59599516SKenneth E. Jansen dxidx = zero 252*59599516SKenneth E. Jansen WdetJ = zero 253*59599516SKenneth E. Jansen shg = zero 254*59599516SKenneth E. Jansen 255*59599516SKenneth E. Jansenc 256*59599516SKenneth E. Jansenc.... loop through the integration points 257*59599516SKenneth E. Jansenc 258*59599516SKenneth E. Jansen do intp = 1, ngauss ! (these are in common.h) 259*59599516SKenneth E. Jansen 260*59599516SKenneth E. Jansenc.... get the hierarchic shape functions at this int point 261*59599516SKenneth E. Jansenc 262*59599516SKenneth E. Jansenc 263*59599516SKenneth E. Jansenc.... for the mass matrix to be consistent shp and shgl must be 264*59599516SKenneth E. Jansenc.... evaluated with at least higher quadrature than one-pt. quad. 265*59599516SKenneth E. Jansen 266*59599516SKenneth E. Jansen call getshp(shp, shgl, sgn, 267*59599516SKenneth E. Jansen & shape, shdrv, intp) 268*59599516SKenneth E. Jansen 269*59599516SKenneth E. Jansenc 270*59599516SKenneth E. Jansenc.... calculate the determinant of the jacobian and weight it 271*59599516SKenneth E. Jansenc 272*59599516SKenneth E. Jansen call e3metric( xl, shdrv,dxidx,shg,WdetJ) 273*59599516SKenneth E. Jansenc 274*59599516SKenneth E. Jansen 275*59599516SKenneth E. Jansen do iel = 1, npro 276*59599516SKenneth E. Jansen do na = 1, nshl 277*59599516SKenneth E. Jansen do nb = 1, nshl 278*59599516SKenneth E. Jansen shp2 = shape(iel,na) * shape(iel,nb) 279*59599516SKenneth E. Jansen em(iel,na,nb) = em(iel,na,nb) + shp2*WdetJ(iel) 280*59599516SKenneth E. Jansen enddo 281*59599516SKenneth E. Jansen enddo 282*59599516SKenneth E. Jansen enddo 283*59599516SKenneth E. Jansen 284*59599516SKenneth E. Jansen enddo ! End loop over quadrature points 285*59599516SKenneth E. Jansen 286*59599516SKenneth E. Jansen 287*59599516SKenneth E. Jansen 288*59599516SKenneth E. Jansenc----------> Obtain a mixed (lumped/consistent) mass matrix <------------ 289*59599516SKenneth E. Jansen 290*59599516SKenneth E. Jansenc... Different combinations of the lump and mass matrices yield 291*59599516SKenneth E. Jansenc... filters of varying widths. In the limiting case were 292*59599516SKenneth E. Jansenc... the entire matrix is lumped, we obtain the same filter as 293*59599516SKenneth E. Jansenc... in getdmc.f. Note that in these equivalent ways of 294*59599516SKenneth E. Jansenc... filtering one-point quadrature is used for shpf and shgl.. 295*59599516SKenneth E. Jansen 296*59599516SKenneth E. Jansen 297*59599516SKenneth E. Jansen do iel = 1, npro 298*59599516SKenneth E. Jansen do na = 1, nshl 299*59599516SKenneth E. Jansen em(iel,na,na) = eml(iel,na) 300*59599516SKenneth E. Jansen enddo 301*59599516SKenneth E. Jansen enddo 302*59599516SKenneth E. Jansen 303*59599516SKenneth E. Jansenc 304*59599516SKenneth E. Jansenc.... return 305*59599516SKenneth E. Jansenc 306*59599516SKenneth E. Jansen return 307*59599516SKenneth E. Jansen end 308*59599516SKenneth E. Jansen 309*59599516SKenneth E. Jansenc----------------------------------------------------------------------- 310*59599516SKenneth E. Jansenc 311*59599516SKenneth E. Jansenc compute the metrics of the mapping from global to local 312*59599516SKenneth E. Jansenc coordinates and the jacobian of the mapping (weighted by 313*59599516SKenneth E. Jansenc the quadrature weight 314*59599516SKenneth E. Jansenc 315*59599516SKenneth E. Jansenc----------------------------------------------------------------------- 316*59599516SKenneth E. Jansen subroutine e3metricf( xl, shgl, dxidx, 317*59599516SKenneth E. Jansen & shg, WdetJ, Qwtf) 318*59599516SKenneth E. Jansen 319*59599516SKenneth E. Jansen include "common.h" 320*59599516SKenneth E. Jansen 321*59599516SKenneth E. Jansen real*8 xl(npro,nenl,nsd), shgl(npro,nsd,nshl), 322*59599516SKenneth E. Jansen & dxidx(npro,nsd,nsd), shg(npro,nshl,nsd), 323*59599516SKenneth E. Jansen & WdetJ(npro), Qwtf(ngaussf) 324*59599516SKenneth E. Jansen 325*59599516SKenneth E. Jansen real*8 dxdxi(npro,nsd,nsd), tmp(npro) 326*59599516SKenneth E. Jansen 327*59599516SKenneth E. Jansenc 328*59599516SKenneth E. Jansenc.... compute the deformation gradient 329*59599516SKenneth E. Jansenc 330*59599516SKenneth E. Jansen dxdxi = zero 331*59599516SKenneth E. Jansenc 332*59599516SKenneth E. Jansen do n = 1, nenl 333*59599516SKenneth E. Jansen dxdxi(:,1,1) = dxdxi(:,1,1) + xl(:,n,1) * shgl(:,1,n) 334*59599516SKenneth E. Jansen dxdxi(:,1,2) = dxdxi(:,1,2) + xl(:,n,1) * shgl(:,2,n) 335*59599516SKenneth E. Jansen dxdxi(:,1,3) = dxdxi(:,1,3) + xl(:,n,1) * shgl(:,3,n) 336*59599516SKenneth E. Jansen dxdxi(:,2,1) = dxdxi(:,2,1) + xl(:,n,2) * shgl(:,1,n) 337*59599516SKenneth E. Jansen dxdxi(:,2,2) = dxdxi(:,2,2) + xl(:,n,2) * shgl(:,2,n) 338*59599516SKenneth E. Jansen dxdxi(:,2,3) = dxdxi(:,2,3) + xl(:,n,2) * shgl(:,3,n) 339*59599516SKenneth E. Jansen dxdxi(:,3,1) = dxdxi(:,3,1) + xl(:,n,3) * shgl(:,1,n) 340*59599516SKenneth E. Jansen dxdxi(:,3,2) = dxdxi(:,3,2) + xl(:,n,3) * shgl(:,2,n) 341*59599516SKenneth E. Jansen dxdxi(:,3,3) = dxdxi(:,3,3) + xl(:,n,3) * shgl(:,3,n) 342*59599516SKenneth E. Jansen enddo 343*59599516SKenneth E. Jansenc 344*59599516SKenneth E. Jansenc.... compute the inverse of deformation gradient 345*59599516SKenneth E. Jansenc 346*59599516SKenneth E. Jansen dxidx(:,1,1) = dxdxi(:,2,2) * dxdxi(:,3,3) 347*59599516SKenneth E. Jansen & - dxdxi(:,3,2) * dxdxi(:,2,3) 348*59599516SKenneth E. Jansen dxidx(:,1,2) = dxdxi(:,3,2) * dxdxi(:,1,3) 349*59599516SKenneth E. Jansen & - dxdxi(:,1,2) * dxdxi(:,3,3) 350*59599516SKenneth E. Jansen dxidx(:,1,3) = dxdxi(:,1,2) * dxdxi(:,2,3) 351*59599516SKenneth E. Jansen & - dxdxi(:,1,3) * dxdxi(:,2,2) 352*59599516SKenneth E. Jansen tmp = one / ( dxidx(:,1,1) * dxdxi(:,1,1) 353*59599516SKenneth E. Jansen & + dxidx(:,1,2) * dxdxi(:,2,1) 354*59599516SKenneth E. Jansen & + dxidx(:,1,3) * dxdxi(:,3,1) ) 355*59599516SKenneth E. Jansen dxidx(:,1,1) = dxidx(:,1,1) * tmp 356*59599516SKenneth E. Jansen dxidx(:,1,2) = dxidx(:,1,2) * tmp 357*59599516SKenneth E. Jansen dxidx(:,1,3) = dxidx(:,1,3) * tmp 358*59599516SKenneth E. Jansen dxidx(:,2,1) = (dxdxi(:,2,3) * dxdxi(:,3,1) 359*59599516SKenneth E. Jansen & - dxdxi(:,2,1) * dxdxi(:,3,3)) * tmp 360*59599516SKenneth E. Jansen dxidx(:,2,2) = (dxdxi(:,1,1) * dxdxi(:,3,3) 361*59599516SKenneth E. Jansen & - dxdxi(:,3,1) * dxdxi(:,1,3)) * tmp 362*59599516SKenneth E. Jansen dxidx(:,2,3) = (dxdxi(:,2,1) * dxdxi(:,1,3) 363*59599516SKenneth E. Jansen & - dxdxi(:,1,1) * dxdxi(:,2,3)) * tmp 364*59599516SKenneth E. Jansen dxidx(:,3,1) = (dxdxi(:,2,1) * dxdxi(:,3,2) 365*59599516SKenneth E. Jansen & - dxdxi(:,2,2) * dxdxi(:,3,1)) * tmp 366*59599516SKenneth E. Jansen dxidx(:,3,2) = (dxdxi(:,3,1) * dxdxi(:,1,2) 367*59599516SKenneth E. Jansen & - dxdxi(:,1,1) * dxdxi(:,3,2)) * tmp 368*59599516SKenneth E. Jansen dxidx(:,3,3) = (dxdxi(:,1,1) * dxdxi(:,2,2) 369*59599516SKenneth E. Jansen & - dxdxi(:,1,2) * dxdxi(:,2,1)) * tmp 370*59599516SKenneth E. Jansenc 371*59599516SKenneth E. Jansenc WdetJ = Qwt(lcsyst,intp) / tmp 372*59599516SKenneth E. Jansen 373*59599516SKenneth E. Jansen WdetJ = Qwtf(intp) / tmp 374*59599516SKenneth E. Jansenc 375*59599516SKenneth E. Jansenc.... compute the global gradient of shape-functions 376*59599516SKenneth E. Jansenc 377*59599516SKenneth E. Jansen do n = 1, nshl 378*59599516SKenneth E. Jansen shg(:,n,1) = shgl(:,1,n) * dxidx(:,1,1) + 379*59599516SKenneth E. Jansen & shgl(:,2,n) * dxidx(:,2,1) + 380*59599516SKenneth E. Jansen & shgl(:,3,n) * dxidx(:,3,1) 381*59599516SKenneth E. Jansen shg(:,n,2) = shgl(:,1,n) * dxidx(:,1,2) + 382*59599516SKenneth E. Jansen & shgl(:,2,n) * dxidx(:,2,2) + 383*59599516SKenneth E. Jansen & shgl(:,3,n) * dxidx(:,3,2) 384*59599516SKenneth E. Jansen shg(:,n,3) = shgl(:,1,n) * dxidx(:,1,3) + 385*59599516SKenneth E. Jansen & shgl(:,2,n) * dxidx(:,2,3) + 386*59599516SKenneth E. Jansen & shgl(:,3,n) * dxidx(:,3,3) 387*59599516SKenneth E. Jansen enddo 388*59599516SKenneth E. Jansen 389*59599516SKenneth E. Jansen return 390*59599516SKenneth E. Jansen end 391*59599516SKenneth E. Jansen 392*59599516SKenneth E. Jansen 393*59599516SKenneth E. Jansen 394