xref: /phasta/phSolver/incompressible/hessian.f (revision 7acde132a6def0fe2daaec0d1a712dff0e5c6636)
159599516SKenneth E. Jansen        subroutine hessian ( y,         x,
259599516SKenneth E. Jansen     &                       shp,       shgl,      iBC,
359599516SKenneth E. Jansen     &                       shpb,      shglb,     iper,
459599516SKenneth E. Jansen     &                       ilwork,    uhess,     gradu  )
559599516SKenneth E. Jansen        use pointer_data  ! brings in the pointers for the blocked arrays
659599516SKenneth E. Jansen
759599516SKenneth E. Jansen        include "common.h"
859599516SKenneth E. Jansenc
959599516SKenneth E. Jansen        dimension y(nshg,ndof),
1059599516SKenneth E. Jansen     &            x(numnp,nsd),         iBC(nshg),
1159599516SKenneth E. Jansen     &            iper(nshg)
1259599516SKenneth E. Jansenc
1359599516SKenneth E. Jansen        dimension shp(MAXTOP,maxsh,MAXQPT),
1459599516SKenneth E. Jansen     &            shgl(MAXTOP,nsd,maxsh,MAXQPT),
1559599516SKenneth E. Jansen     &            shpb(MAXTOP,maxsh,MAXQPT),
1659599516SKenneth E. Jansen     &            shglb(MAXTOP,nsd,maxsh,MAXQPT)
1759599516SKenneth E. Jansenc
1859599516SKenneth E. Jansen        dimension gradu(nshg,9),     rmass(nshg),
1959599516SKenneth E. Jansen     &            uhess(nshg,27)
2059599516SKenneth E. Jansenc
2159599516SKenneth E. Jansen        dimension ilwork(nlwork)
2259599516SKenneth E. Jansen
2359599516SKenneth E. Jansenc
2459599516SKenneth E. Jansen           gradu = zero
2559599516SKenneth E. Jansen           rmass = zero
2659599516SKenneth E. Jansen
2759599516SKenneth E. Jansen           do iblk = 1, nelblk
2859599516SKenneth E. Jansen              iel    = lcblk(1,iblk)
2959599516SKenneth E. Jansen              lelCat = lcblk(2,iblk)
3059599516SKenneth E. Jansen              lcsyst = lcblk(3,iblk)
3159599516SKenneth E. Jansen              iorder = lcblk(4,iblk)
3259599516SKenneth E. Jansen              nenl   = lcblk(5,iblk) ! no. of vertices per element
3359599516SKenneth E. Jansen              nshl   = lcblk(10,iblk)
3459599516SKenneth E. Jansen              mattyp = lcblk(7,iblk)
3559599516SKenneth E. Jansen              ndofl  = lcblk(8,iblk)
3659599516SKenneth E. Jansen              nsymdl = lcblk(9,iblk)
3759599516SKenneth E. Jansen              npro   = lcblk(1,iblk+1) - iel
3859599516SKenneth E. Jansen              ngauss = nint(lcsyst)
3959599516SKenneth E. Jansenc
4059599516SKenneth E. Jansen              call velocity_gradient ( y,
4159599516SKenneth E. Jansen     &                                 x,
4259599516SKenneth E. Jansen     &                                 shp(lcsyst,1:nshl,:),
4359599516SKenneth E. Jansen     &                                 shgl(lcsyst,:,1:nshl,:),
4459599516SKenneth E. Jansen     &                                 mien(iblk)%p,
4559599516SKenneth E. Jansen     &                                 gradu,
4659599516SKenneth E. Jansen     &                                 rmass )
4759599516SKenneth E. Jansen
4859599516SKenneth E. Jansen           end do
4959599516SKenneth E. Jansen
5059599516SKenneth E. Jansenc
5159599516SKenneth E. Jansen           call reconstruct( rmass, gradu, iBC, iper, ilwork, 9 )
5259599516SKenneth E. Jansen
5359599516SKenneth E. Jansen           uhess = zero
5459599516SKenneth E. Jansen           rmass = zero
5559599516SKenneth E. Jansen
5659599516SKenneth E. Jansen           do iblk = 1, nelblk
5759599516SKenneth E. Jansen              iel    = lcblk(1,iblk)
5859599516SKenneth E. Jansen              lelCat = lcblk(2,iblk)
5959599516SKenneth E. Jansen              lcsyst = lcblk(3,iblk)
6059599516SKenneth E. Jansen              iorder = lcblk(4,iblk)
6159599516SKenneth E. Jansen              nenl   = lcblk(5,iblk) ! no. of vertices per element
6259599516SKenneth E. Jansen              nshl   = lcblk(10,iblk)
6359599516SKenneth E. Jansen              mattyp = lcblk(7,iblk)
6459599516SKenneth E. Jansen              ndofl  = lcblk(8,iblk)
6559599516SKenneth E. Jansen              nsymdl = lcblk(9,iblk)
6659599516SKenneth E. Jansen              npro   = lcblk(1,iblk+1) - iel
6759599516SKenneth E. Jansen              ngauss = nint(lcsyst)
6859599516SKenneth E. Jansenc
6959599516SKenneth E. Jansen
7059599516SKenneth E. Jansen              call velocity_hessian (  gradu,
7159599516SKenneth E. Jansen     &                                 x,
7259599516SKenneth E. Jansen     &                                 shp(lcsyst,1:nshl,:),
7359599516SKenneth E. Jansen     &                                 shgl(lcsyst,:,1:nshl,:),
7459599516SKenneth E. Jansen     &                                 mien(iblk)%p,
7559599516SKenneth E. Jansen     &                                 uhess,
7659599516SKenneth E. Jansen     &				       rmass  )
7759599516SKenneth E. Jansen           end do
7859599516SKenneth E. Jansen
7959599516SKenneth E. Jansen
8059599516SKenneth E. Jansen           call reconstruct( rmass, uhess, iBC, iper, ilwork, 27 )
8159599516SKenneth E. Jansenc
8259599516SKenneth E. Jansen      return
8359599516SKenneth E. Jansen      end
8459599516SKenneth E. Jansen
8559599516SKenneth E. Jansenc-----------------------------------------------------------------------------
8659599516SKenneth E. Jansen
8759599516SKenneth E. Jansen        subroutine velocity_gradient ( y,       x,       shp,     shgl,
8859599516SKenneth E. Jansen     &                                 ien,     gradu,   rmass    )
8959599516SKenneth E. Jansen
9059599516SKenneth E. Jansen        include "common.h"
9159599516SKenneth E. Jansenc
9259599516SKenneth E. Jansen        dimension y(nshg,ndof),               x(numnp,nsd),
9359599516SKenneth E. Jansen     &            shp(nshl,ngauss),           shgl(nsd,nshl,ngauss),
9459599516SKenneth E. Jansen     &            ien(npro,nshl),             gradu(nshg,9),
9559599516SKenneth E. Jansen     &            shdrv(npro,nsd,nshl),       shape( npro, nshl ),
9659599516SKenneth E. Jansen     &            gradul(npro,9) ,            rmass( nshg )
9759599516SKenneth E. Jansenc
9859599516SKenneth E. Jansen        dimension yl(npro,nshl,ndof),          xl(npro,nenl,nsd),
9959599516SKenneth E. Jansen     &            ql(npro,nshl,9),             dxidx(npro,nsd,nsd),
10059599516SKenneth E. Jansen     &            WdetJ(npro),		       rmassl(npro,nshl)
10159599516SKenneth E. Jansenc
10259599516SKenneth E. Jansenc
10359599516SKenneth E. Jansen        dimension sgn(npro,nshl)
10459599516SKenneth E. Jansenc
10559599516SKenneth E. Jansenc.... create the matrix of mode signs for the hierarchic basis
10659599516SKenneth E. Jansenc     functions.
10759599516SKenneth E. Jansenc
10859599516SKenneth E. Jansen        do i=1,nshl
10959599516SKenneth E. Jansen           where ( ien(:,i) < 0 )
11059599516SKenneth E. Jansen              sgn(:,i) = -one
11159599516SKenneth E. Jansen           elsewhere
11259599516SKenneth E. Jansen              sgn(:,i) = one
11359599516SKenneth E. Jansen           endwhere
11459599516SKenneth E. Jansen        enddo
11559599516SKenneth E. Jansen
11659599516SKenneth E. Jansenc
11759599516SKenneth E. Jansenc.... gather the variables
11859599516SKenneth E. Jansenc
11959599516SKenneth E. Jansen
12059599516SKenneth E. Jansen        call localy (y,    yl,     ien,    ndof,   'gather  ')
12159599516SKenneth E. Jansen        call localx (x,    xl,     ien,    nsd,    'gather  ')
12259599516SKenneth E. Jansenc
12359599516SKenneth E. Jansenc.... get the element residuals
12459599516SKenneth E. Jansenc
12559599516SKenneth E. Jansen        ql     = zero
12659599516SKenneth E. Jansen        rmassl = zero
12759599516SKenneth E. Jansen
12859599516SKenneth E. Jansen        do intp = 1, ngauss
12959599516SKenneth E. Jansen
13059599516SKenneth E. Jansen            if ( Qwt( lcsyst, intp ) .eq. zero ) cycle
13159599516SKenneth E. Jansen
13259599516SKenneth E. Jansen            gradul = zero
13359599516SKenneth E. Jansen            call getshp( shp, shgl, sgn, shape, shdrv )
13459599516SKenneth E. Jansen            call local_gradient( yl(:,:,2:4), 3,  shdrv, xl,
13559599516SKenneth E. Jansen     &                           gradul , dxidx, WdetJ )
13659599516SKenneth E. Jansen
13759599516SKenneth E. Jansenc.... assemble contribution of gradu to each element node
13859599516SKenneth E. Jansenc
13959599516SKenneth E. Jansen            do i=1,nshl
14059599516SKenneth E. Jansen                do j = 1, 9
14159599516SKenneth E. Jansen                    ql(:,i,j) = ql(:,i,j)+shape(:,i)*WdetJ*gradul(:,j)
14259599516SKenneth E. Jansen                end do
14359599516SKenneth E. Jansen
14459599516SKenneth E. Jansen                rmassl(:,i) = rmassl(:,i) + shape(:,i)*WdetJ
14559599516SKenneth E. Jansen
14659599516SKenneth E. Jansen             end do
14759599516SKenneth E. Jansen
14859599516SKenneth E. Jansen        end do
14959599516SKenneth E. Jansenc
15059599516SKenneth E. Jansenc
15159599516SKenneth E. Jansen        call local (gradu,  ql,     ien,  9,  'scatter ')
15259599516SKenneth E. Jansen        call local (rmass,  rmassl, ien,  1,  'scatter ')
15359599516SKenneth E. Jansenc
15459599516SKenneth E. Jansenc.... end
15559599516SKenneth E. Jansenc
15659599516SKenneth E. Jansen        return
15759599516SKenneth E. Jansen        end
15859599516SKenneth E. Jansen
15959599516SKenneth E. Jansen
16059599516SKenneth E. Jansenc-----------------------------------------------------------------------------
16159599516SKenneth E. Jansen
16259599516SKenneth E. Jansen        subroutine velocity_hessian ( gradu,   x,     shp,   shgl,
16359599516SKenneth E. Jansen     &                                ien,     uhess, rmass  )
16459599516SKenneth E. Jansen
16559599516SKenneth E. Jansen        include "common.h"
16659599516SKenneth E. Jansenc
16759599516SKenneth E. Jansen        dimension gradu(nshg,9),              x(numnp,nsd),
16859599516SKenneth E. Jansen     &            shp(nshl,ngauss),           shgl(nsd,nshl,ngauss),
16959599516SKenneth E. Jansen     &            ien(npro,nshl),             uhess(nshg,27),
17059599516SKenneth E. Jansen     &            shdrv(npro,nsd,nshl),       shape( npro, nshl ),
17159599516SKenneth E. Jansen     &            uhessl(npro,27),            rmass( nshg )
17259599516SKenneth E. Jansenc
17359599516SKenneth E. Jansen        dimension gradul(npro,nshl,9),          xl(npro,nenl,nsd),
17459599516SKenneth E. Jansen     &            ql(npro,nshl,27),             dxidx(npro,nsd,nsd),
17559599516SKenneth E. Jansen     &            WdetJ(npro),                  rmassl(npro, nshl)
17659599516SKenneth E. Jansenc
17759599516SKenneth E. Jansenc
17859599516SKenneth E. Jansen        dimension sgn(npro,nshl)
17959599516SKenneth E. Jansenc
18059599516SKenneth E. Jansenc.... create the matrix of mode signs for the hierarchic basis
18159599516SKenneth E. Jansenc     functions.
18259599516SKenneth E. Jansenc
18359599516SKenneth E. Jansen        do i=1,nshl
18459599516SKenneth E. Jansen           where ( ien(:,i) < 0 )
18559599516SKenneth E. Jansen              sgn(:,i) = -one
18659599516SKenneth E. Jansen           elsewhere
18759599516SKenneth E. Jansen              sgn(:,i) = one
18859599516SKenneth E. Jansen           endwhere
18959599516SKenneth E. Jansen        enddo
19059599516SKenneth E. Jansen
19159599516SKenneth E. Jansenc
19259599516SKenneth E. Jansenc.... gather the variables
19359599516SKenneth E. Jansenc
19459599516SKenneth E. Jansen
19559599516SKenneth E. Jansen        call local  (gradu,  gradul, ien,    9 ,   'gather  ')
19659599516SKenneth E. Jansen        call localx (x,      xl,     ien,    nsd,  'gather  ')
19759599516SKenneth E. Jansenc
19859599516SKenneth E. Jansenc.... get the element residuals
19959599516SKenneth E. Jansenc
20059599516SKenneth E. Jansen        ql     = zero
20159599516SKenneth E. Jansen	rmassl = zero
20259599516SKenneth E. Jansen
20359599516SKenneth E. Jansen        do intp = 1, ngauss
20459599516SKenneth E. Jansen
20559599516SKenneth E. Jansen            if ( Qwt( lcsyst, intp ) .eq. zero ) cycle
20659599516SKenneth E. Jansen
20759599516SKenneth E. Jansen            uhessl = zero
20859599516SKenneth E. Jansen            call getshp( shp, shgl, sgn, shape, shdrv )
20959599516SKenneth E. Jansen            call local_gradient( gradul, 9,  shdrv, xl,
21059599516SKenneth E. Jansen     &                           uhessl , dxidx, WdetJ )
21159599516SKenneth E. Jansen
21259599516SKenneth E. Jansenc.... assemble contribution of gradu .,
21359599516SKenneth E. Jansenc
21459599516SKenneth E. Jansen            do i=1,nshl
21559599516SKenneth E. Jansen                do j = 1,27
21659599516SKenneth E. Jansen                    ql(:,i,j)=ql(:,i,j)+shape(:,i)*WdetJ*uhessl(:,j )
21759599516SKenneth E. Jansen                end do
21859599516SKenneth E. Jansen
21959599516SKenneth E. Jansen                rmassl(:,i) = rmassl(:,i) + shape(:,i)*WdetJ
22059599516SKenneth E. Jansen             end do
22159599516SKenneth E. Jansen
22259599516SKenneth E. Jansen        end do
22359599516SKenneth E. Jansenc
22459599516SKenneth E. Jansenc
22559599516SKenneth E. Jansen        call local (uhess,  ql,     ien,  27,     'scatter ')
22659599516SKenneth E. Jansen        call local (rmass,  rmassl, ien,   1,     'scatter ')
22759599516SKenneth E. Jansenc
22859599516SKenneth E. Jansenc.... end
22959599516SKenneth E. Jansenc
23059599516SKenneth E. Jansen        return
23159599516SKenneth E. Jansen        end
23259599516SKenneth E. Jansen
23359599516SKenneth E. Jansen
23459599516SKenneth E. Jansenc--------------------------------------------------------------------
23559599516SKenneth E. Jansen      subroutine reconstruct( rmass, qres, iBC, iper, ilwork, vsize )
23659599516SKenneth E. Jansen
23759599516SKenneth E. Jansen      include "common.h"
23859599516SKenneth E. Jansen
23959599516SKenneth E. Jansen      integer vsize
24059599516SKenneth E. Jansen      dimension rmass(nshg), qres( nshg, vsize),
241*3dbc338dSCameron Smith     &          iBC(nshg), iper(nshg), ilwork(nlwork)
24259599516SKenneth E. Jansenc
24359599516SKenneth E. Jansenc
24459599516SKenneth E. Jansenc.... compute qi for node A, i.e., qres <-- qres/rmass
24559599516SKenneth E. Jansenc
24659599516SKenneth E. Jansen       if (numpe > 1) then
24759599516SKenneth E. Jansen          call commu ( qres  , ilwork,  vsize  , 'in ')
24859599516SKenneth E. Jansen          call commu ( rmass , ilwork,  1  , 'in ')
24959599516SKenneth E. Jansen       endif
25059599516SKenneth E. Jansenc
25159599516SKenneth E. Jansenc  take care of periodic boundary conditions
25259599516SKenneth E. Jansenc
25359599516SKenneth E. Jansen        do j= 1,nshg
25459599516SKenneth E. Jansen          if ((btest(iBC(j),10))) then
25559599516SKenneth E. Jansen            i = iper(j)
25659599516SKenneth E. Jansen            rmass(i) = rmass(i) + rmass(j)
25759599516SKenneth E. Jansen            qres(i,:) = qres(i,:) + qres(j,:)
25859599516SKenneth E. Jansen          endif
25959599516SKenneth E. Jansen        enddo
26059599516SKenneth E. Jansen
26159599516SKenneth E. Jansen        do j= 1,nshg
26259599516SKenneth E. Jansen          if ((btest(iBC(j),10))) then
26359599516SKenneth E. Jansen            i = iper(j)
26459599516SKenneth E. Jansen            rmass(j) = rmass(i)
26559599516SKenneth E. Jansen            qres(j,:) = qres(i,:)
26659599516SKenneth E. Jansen          endif
26759599516SKenneth E. Jansen        enddo
26859599516SKenneth E. Jansenc
26959599516SKenneth E. Jansenc.... invert the diagonal mass matrix and find q
27059599516SKenneth E. Jansenc
27159599516SKenneth E. Jansen        rmass = one/rmass
27259599516SKenneth E. Jansen
27359599516SKenneth E. Jansen       do i=1,vsize
27459599516SKenneth E. Jansen          qres(:,i) = rmass*qres(:,i)
27559599516SKenneth E. Jansen       enddo
27659599516SKenneth E. Jansen
27759599516SKenneth E. Jansen       if(numpe > 1) then
27859599516SKenneth E. Jansen          call commu (qres, ilwork, vsize, 'out')
27959599516SKenneth E. Jansen       endif
28059599516SKenneth E. Jansen
28159599516SKenneth E. Jansenc.... return
28259599516SKenneth E. Jansenc
28359599516SKenneth E. Jansen        return
28459599516SKenneth E. Jansen        end
28559599516SKenneth E. Jansen
28659599516SKenneth E. Jansenc-------------------------------------------------------------------------
28759599516SKenneth E. Jansen
28859599516SKenneth E. Jansen        subroutine local_gradient ( vector,   vsize, shgl,   xl,
28959599516SKenneth E. Jansen     &                              gradient, dxidx,   WdetJ )
29059599516SKenneth E. Jansenc
29159599516SKenneth E. Jansenc
29259599516SKenneth E. Jansen        include "common.h"
29359599516SKenneth E. Jansenc
29459599516SKenneth E. Jansenc  passed arrays
29559599516SKenneth E. Jansen
29659599516SKenneth E. Jansen        integer vsize
29759599516SKenneth E. Jansenc
29859599516SKenneth E. Jansen        dimension vector(npro,nshl,vsize),
29959599516SKenneth E. Jansen     &            shgl(npro,nsd,nshl),        xl(npro,nenl,nsd),
30059599516SKenneth E. Jansen     &            gradient(npro,vsize*3),     shg(npro,nshl,nsd),
30159599516SKenneth E. Jansen     &            dxidx(npro,nsd,nsd),        WdetJ(npro)
30259599516SKenneth E. Jansenc
30359599516SKenneth E. Jansenc  local arrays
30459599516SKenneth E. Jansenc
30559599516SKenneth E. Jansen        dimension tmp(npro),           dxdxi(npro,nsd,nsd)
30659599516SKenneth E. Jansen
30759599516SKenneth E. Jansenc
30859599516SKenneth E. Jansenc.... compute the deformation gradient
30959599516SKenneth E. Jansenc
31059599516SKenneth E. Jansen        dxdxi = zero
31159599516SKenneth E. Jansenc
31259599516SKenneth E. Jansen          do n = 1, nenl
31359599516SKenneth E. Jansen            dxdxi(:,1,1) = dxdxi(:,1,1) + xl(:,n,1) * shgl(:,1,n)
31459599516SKenneth E. Jansen            dxdxi(:,1,2) = dxdxi(:,1,2) + xl(:,n,1) * shgl(:,2,n)
31559599516SKenneth E. Jansen            dxdxi(:,1,3) = dxdxi(:,1,3) + xl(:,n,1) * shgl(:,3,n)
31659599516SKenneth E. Jansen            dxdxi(:,2,1) = dxdxi(:,2,1) + xl(:,n,2) * shgl(:,1,n)
31759599516SKenneth E. Jansen            dxdxi(:,2,2) = dxdxi(:,2,2) + xl(:,n,2) * shgl(:,2,n)
31859599516SKenneth E. Jansen            dxdxi(:,2,3) = dxdxi(:,2,3) + xl(:,n,2) * shgl(:,3,n)
31959599516SKenneth E. Jansen            dxdxi(:,3,1) = dxdxi(:,3,1) + xl(:,n,3) * shgl(:,1,n)
32059599516SKenneth E. Jansen            dxdxi(:,3,2) = dxdxi(:,3,2) + xl(:,n,3) * shgl(:,2,n)
32159599516SKenneth E. Jansen            dxdxi(:,3,3) = dxdxi(:,3,3) + xl(:,n,3) * shgl(:,3,n)
32259599516SKenneth E. Jansen          enddo
32359599516SKenneth E. Jansenc
32459599516SKenneth E. Jansenc.... compute the inverse of deformation gradient
32559599516SKenneth E. Jansenc
32659599516SKenneth E. Jansen        dxidx(:,1,1) =   dxdxi(:,2,2) * dxdxi(:,3,3)
32759599516SKenneth E. Jansen     &                 - dxdxi(:,3,2) * dxdxi(:,2,3)
32859599516SKenneth E. Jansen        dxidx(:,1,2) =   dxdxi(:,3,2) * dxdxi(:,1,3)
32959599516SKenneth E. Jansen     &                 - dxdxi(:,1,2) * dxdxi(:,3,3)
33059599516SKenneth E. Jansen        dxidx(:,1,3) =   dxdxi(:,1,2) * dxdxi(:,2,3)
33159599516SKenneth E. Jansen     &                 - dxdxi(:,1,3) * dxdxi(:,2,2)
33259599516SKenneth E. Jansen        tmp          = one / ( dxidx(:,1,1) * dxdxi(:,1,1)
33359599516SKenneth E. Jansen     &                       + dxidx(:,1,2) * dxdxi(:,2,1)
33459599516SKenneth E. Jansen     &                       + dxidx(:,1,3) * dxdxi(:,3,1) )
33559599516SKenneth E. Jansen        dxidx(:,1,1) = dxidx(:,1,1) * tmp
33659599516SKenneth E. Jansen        dxidx(:,1,2) = dxidx(:,1,2) * tmp
33759599516SKenneth E. Jansen        dxidx(:,1,3) = dxidx(:,1,3) * tmp
33859599516SKenneth E. Jansen        dxidx(:,2,1) = (dxdxi(:,2,3) * dxdxi(:,3,1)
33959599516SKenneth E. Jansen     &                - dxdxi(:,2,1) * dxdxi(:,3,3)) * tmp
34059599516SKenneth E. Jansen        dxidx(:,2,2) = (dxdxi(:,1,1) * dxdxi(:,3,3)
34159599516SKenneth E. Jansen     &                - dxdxi(:,3,1) * dxdxi(:,1,3)) * tmp
34259599516SKenneth E. Jansen        dxidx(:,2,3) = (dxdxi(:,2,1) * dxdxi(:,1,3)
34359599516SKenneth E. Jansen     &                - dxdxi(:,1,1) * dxdxi(:,2,3)) * tmp
34459599516SKenneth E. Jansen        dxidx(:,3,1) = (dxdxi(:,2,1) * dxdxi(:,3,2)
34559599516SKenneth E. Jansen     &                - dxdxi(:,2,2) * dxdxi(:,3,1)) * tmp
34659599516SKenneth E. Jansen        dxidx(:,3,2) = (dxdxi(:,3,1) * dxdxi(:,1,2)
34759599516SKenneth E. Jansen     &                - dxdxi(:,1,1) * dxdxi(:,3,2)) * tmp
34859599516SKenneth E. Jansen        dxidx(:,3,3) = (dxdxi(:,1,1) * dxdxi(:,2,2)
34959599516SKenneth E. Jansen     &                - dxdxi(:,1,2) * dxdxi(:,2,1)) * tmp
35059599516SKenneth E. Jansenc
35159599516SKenneth E. Jansen        WdetJ = Qwt(lcsyst,intp)/ tmp
35259599516SKenneth E. Jansen
35359599516SKenneth E. Jansenc
35459599516SKenneth E. Jansenc.... --------------------->  Global Gradients  <-----------------------
35559599516SKenneth E. Jansenc
35659599516SKenneth E. Jansen        gradient = zero
35759599516SKenneth E. Jansenc
35859599516SKenneth E. Jansenc
35959599516SKenneth E. Jansen        do n = 1, nshl
36059599516SKenneth E. Jansenc
36159599516SKenneth E. Jansenc.... compute the global gradient of shape-function
36259599516SKenneth E. Jansenc
36359599516SKenneth E. Jansenc            ! N_{a,x_i}= N_{a,xi_i} xi_{i,x_j}
36459599516SKenneth E. Jansenc
36559599516SKenneth E. Jansen          shg(:,n,1) = shgl(:,1,n) * dxidx(:,1,1) +
36659599516SKenneth E. Jansen     &                 shgl(:,2,n) * dxidx(:,2,1) +
36759599516SKenneth E. Jansen     &                 shgl(:,3,n) * dxidx(:,3,1)
36859599516SKenneth E. Jansen          shg(:,n,2) = shgl(:,1,n) * dxidx(:,1,2) +
36959599516SKenneth E. Jansen     &                 shgl(:,2,n) * dxidx(:,2,2) +
37059599516SKenneth E. Jansen     &                 shgl(:,3,n) * dxidx(:,3,2)
37159599516SKenneth E. Jansen          shg(:,n,3) = shgl(:,1,n) * dxidx(:,1,3) +
37259599516SKenneth E. Jansen     &                 shgl(:,2,n) * dxidx(:,2,3) +
37359599516SKenneth E. Jansen     &                 shgl(:,3,n) * dxidx(:,3,3)
37459599516SKenneth E. Jansenc
37559599516SKenneth E. Jansenc
37659599516SKenneth E. Jansenc  Y_{,x_i}=SUM_{a=1}^nenl (N_{a,x_i}(int) Ya)
37759599516SKenneth E. Jansenc
37859599516SKenneth E. Jansen          do i = 1, 3
37959599516SKenneth E. Jansen            do j = 1, vsize
38059599516SKenneth E. Jansen               k = (i-1)*vsize+j
38159599516SKenneth E. Jansen               gradient(:,k) = gradient(:,k) + shg(:,n,i)*vector(:,n,j)
38259599516SKenneth E. Jansen            end do
38359599516SKenneth E. Jansen          end do
38459599516SKenneth E. Jansen
38559599516SKenneth E. Jansen       end do
38659599516SKenneth E. Jansen
38759599516SKenneth E. Jansenc
38859599516SKenneth E. Jansenc.... return
38959599516SKenneth E. Jansenc
39059599516SKenneth E. Jansen       return
39159599516SKenneth E. Jansen       end
392