!------------------------------------------------------------------------------|
!------------------------------------------------------------------------------|
!                                                                              |
!              ||===\\                                                         | 
!              ||    \\                                                        |
!              ||     ||   //==\\   ||  ||   //==||  ||/==\\                   |
!              ||     ||  ||    ||  ||  ||  ||   ||  ||    ||                  |
!              ||    //   ||    ||  ||  ||  ||   ||  ||                        |
!              ||===//     \\==//    \\==\\  \\==\\  ||                        |
!                                                                              |
!------------------------------------------------------------------------------|
!------------------------------------------------------------------------------|
!                                                                              |
!              MAKE_MATRIX    Jun. 2007                                        |
!                                                                              |
!------------------------------------------------------------------------------|
!------------------------------------------------------------------------------|

subroutine make_matrix (params,ndof,ael,bel,icon,xg,yg,zg,kfix,viscosity0,     &
                        density,penalty,expon,activationenergy,expansion,      &
                        diffusivity,heat,plasticity_type,                      &
                        plasticity_ss_type_coh,plasticity_ss_type_phi,         &
                        plasticity_parameters,fviscosity,unode,vnode,wnode,    &
                        temp,pressure,strain,e2dp,is_plastic,nnode,ffff,r0,s0, &
                        t0,rst,ileaves,eviscosity,vbounded,yield_ratio,        &
                        frict_angle,threadinfo,weight)

!------------------------------------------------------------------------------|
!(((((((((((((((( Purpose of the routine  ))))))))))))))))))))))))))))))))))))))
!------------------------------------------------------------------------------|
! note that this routine is called to create the FE matrix and rhs vector
! for both the Stokes (ndof=3) and Energy equations (ndof=1)

! mpe       : number of nodes per element (8)
! ndof      : number of degrees of freedom per node (3 or 1)
! ael       : computed finite element matrix
! bel       : computed rhs vector
! icon      : connectivity array for the current element
! xg,yg,zg  : global coordinate arrays of length nnode
! kfix      : bc array of length ndof*nnode 
!             (kfix=1 means the dof is fixed to the value stored in velo)
! viscosity : viscosity
! density   : density at zero temperature (in C)
! penalty   : incompressibility penalty factor
! penaltyg  : penalty factor to impose bc
! expon     : nonlinear viscosity exponent
! unode,vnode,wnode : velocity array (obtained from previous time step or at 
!             least containing the proper velocity at the fixed dofs)
! strain    : total accumulated strain as interpolated from the 3D cloud
! e2dp      : effective strain rate from previous step as interpolated from the
!             3D cloud
! nnode     : number of nodes
! f         : global rhs vector

!------------------------------------------------------------------------------|
!((((((((((((((((  declaration of the subroutine arguments  ))))))))))))))))))))
!------------------------------------------------------------------------------|

use constants
use definitions
use gauss
use invariants
!use mpi
use threads

implicit none

include 'mpif.h'

type(parameters) params
integer ndof
double precision ael(params%mpe*ndof,params%mpe*ndof)
double precision bel(params%mpe*ndof)
integer icon(params%mpe) 
double precision xg(nnode),yg(nnode),zg(nnode)
integer kfix(nnode*ndof)
double precision viscosity0
double precision density
double precision penalty
double precision expon
double precision activationenergy
double precision expansion
double precision diffusivity
double precision heat
character (len=8) plasticity_type
character (len=16) plasticity_ss_type_coh,plasticity_ss_type_phi
double precision plasticity_parameters(9)
double precision :: fviscosity
double precision unode(nnode),vnode(nnode),wnode(nnode)
double precision temp(nnode)
double precision pressure
double precision strain(nnode)
double precision e2dp(nnode)
integer nnode
double precision ffff(nnode*ndof)
double precision r0,s0,t0,rst
integer ileaves
double precision eviscosity
logical vbounded
double precision :: yield_ratio
double precision :: frict_angle
logical is_plastic
type (thread) threadinfo
double precision weight

!------------------------------------------------------------------------------|
!(((((((((((((((( declaration of the subroutine internal variables )))))))))))))
!------------------------------------------------------------------------------|
logical is_plastic_temp,flag_vrm_pb
integer i,j,k,k1,k2,k3,ii,jj,ij,iint,nint,nb,err,ic
integer iproc,nproc,ierr,mpe
double precision eps,yield,fail,fixt,yieldnow
double precision r,s,t,w,volume,vol
double precision tr,straintot,e2dprev
double precision xcond,xmass,aref
double precision exx,eyy,ezz,exy,eyz,ezx,e2d,e3d
double precision f1(3,3),defgrad(3,3)
double precision detf1,detf,e2dref,viscomean,frictmean,friction_angle
double precision viscosity,compressibility,temperature
double precision velox,veloy,veloz,tau,uvwnorm,xmin,xmax,ymin,ymax
double precision zmin,zmax,dx,dy,dz,alpha,J2d,J3d
double precision jcb(3,3),jcbi(3,3),jcbp(3,3),jcbip(3,3),d(6,6),dl(6,6)
double precision,dimension(:),allocatable :: x,y,z,tnode,be
double precision,dimension(:),allocatable :: h,dhdr,dhds,dhdt,ht
double precision,dimension(:),allocatable :: dhdx,dhdy,dhdz
double precision,dimension(:),allocatable :: vel
double precision,dimension(:,:),allocatable :: b,bd,bl,aelp

!-------------------------------------------------------------------------------
!-------------------------------------------------------------------------------

call mpi_comm_size (mpi_comm_world,nproc,ierr)
call mpi_comm_rank (mpi_comm_world,iproc,ierr)

compressibility=penalty

mpe=params%mpe

eps=tiny(eps)

nint=8

nb=6
if (ndof.eq.1) nb=3

allocate (x(mpe),stat=err)                  ; if (err.ne.0) call stop_run ('Error alloc x in make_matrix$')
allocate (y(mpe),stat=err)                  ; if (err.ne.0) call stop_run ('Error alloc y in make_matrix$')
allocate (z(mpe),stat=err)                  ; if (err.ne.0) call stop_run ('Error alloc z in make_matrix$')
allocate (tnode(mpe),stat=err)              ; if (err.ne.0) call stop_run ('Error alloc tnode in make_matrix$')
allocate (be(mpe*ndof),stat=err)            ; if (err.ne.0) call stop_run ('Error alloc be in make_matrix$')
allocate (h(mpe),stat=err)                  ; if (err.ne.0) call stop_run ('Error alloc h in make_matrix$')
allocate (dhdr(mpe),stat=err)               ; if (err.ne.0) call stop_run ('Error alloc dhdr in make_matrix$')
allocate (dhds(mpe),stat=err)               ; if (err.ne.0) call stop_run ('Error alloc dhds in make_matrix$')
allocate (dhdt(mpe),stat=err)               ; if (err.ne.0) call stop_run ('Error alloc dhdt in make_matrix$')
allocate (ht(mpe),stat=err)                 ; if (err.ne.0) call stop_run ('Error alloc ht in make_matrix$')
allocate (dhdx(mpe),stat=err)               ; if (err.ne.0) call stop_run ('Error alloc dhdx in make_matrix$')
allocate (dhdy(mpe),stat=err)               ; if (err.ne.0) call stop_run ('Error alloc dhdy in make_matrix$')
allocate (dhdz(mpe),stat=err)               ; if (err.ne.0) call stop_run ('Error alloc dhdz in make_matrix$')
allocate (b(ndof*mpe,nb),stat=err)          ; if (err.ne.0) call stop_run ('Error alloc b in make_matrix$')
allocate (bd(ndof*mpe,nb),stat=err)         ; if (err.ne.0) call stop_run ('Error alloc bd in make_matrix$')
allocate (bl(ndof*mpe,nb),stat=err)         ; if (err.ne.0) call stop_run ('Error alloc bl in make_matrix$')
allocate (aelp(mpe*ndof,mpe*ndof),stat=err) ; if (err.ne.0) call stop_run ('Error alloc aelp in make_matrix$')
allocate (vel(mpe*ndof),stat=err)           ; if (err.ne.0) call stop_run ('Error alloc vel in make_matrix$')

do k=1,mpe
   x(k)=xg(icon(k))
   y(k)=yg(icon(k))
   z(k)=zg(icon(k))*params%vex
   tnode(k)=temp(icon(k))
enddo

ael=0.d0
aelp=0.d0
be=0.d0
bel=0.d0
bl=0.d0
vol=0.d0

!-------------------------------------------------------------------------------
! ndof=3
!-------------------------------------------------------------------------------

if (ndof.eq.1) goto 1111

weight=0.d0

! first add viscosity term (8 point integration)

is_plastic=.false.
is_plastic_temp=.false.
viscomean=0.d0
yield_ratio=0.d0
frict_angle=0.d0
frictmean=0.d0
yield=-1.0

do iint=1,nint
   r=r0+rst*(rr(iint)+1.d0)/2.d0
   s=s0+rst*(ss(iint)+1.d0)/2.d0
   t=t0+rst*(tt(iint)+1.d0)/2.d0
   w=ww(iint)
   call compute_dhdx_dhdy_dhdz (mpe,r,s,t,x,y,z,dhdx,dhdy,dhdz,volume)

   h(1)=(1.d0-r)*(1.d0-s)*(1.d0-t)/8.d0
   h(2)=(1.d0+r)*(1.d0-s)*(1.d0-t)/8.d0
   h(3)=(1.d0-r)*(1.d0+s)*(1.d0-t)/8.d0
   h(4)=(1.d0+r)*(1.d0+s)*(1.d0-t)/8.d0
   h(5)=(1.d0-r)*(1.d0-s)*(1.d0+t)/8.d0
   h(6)=(1.d0+r)*(1.d0-s)*(1.d0+t)/8.d0
   h(7)=(1.d0-r)*(1.d0+s)*(1.d0+t)/8.d0
   h(8)=(1.d0+r)*(1.d0+s)*(1.d0+t)/8.d0
   
   temperature=0.d0
   do k=1,mpe
      temperature=temperature+h(k)*tnode(k)
   enddo

   exx=0.d0
   eyy=0.d0
   ezz=0.d0
   exy=0.d0
   eyz=0.d0
   ezx=0.d0
   do k=1,mpe
      ic=icon(k)
      exx=exx+dhdx(k)*unode(ic)
      eyy=eyy+dhdy(k)*vnode(ic)
      ezz=ezz+dhdz(k)*wnode(ic)
      exy=exy+(dhdx(k)*vnode(ic)+dhdy(k)*unode(ic))/2.d0
      eyz=eyz+(dhdy(k)*wnode(ic)+dhdz(k)*vnode(ic))/2.d0
      ezx=ezx+(dhdz(k)*unode(ic)+dhdx(k)*wnode(ic))/2.d0
   enddo

   call deviatoric (params%invariants_2d,exx,eyy,ezz,exy,eyz,ezx)

   J2d=second_invariant(exx,eyy,ezz,exy,eyz,ezx)
   J3d=third_invariant (exx,eyy,ezz,exy,eyz,ezx)

   e2d=0.d0
   if (J2d.gt.eps) e2d=sqrt(J2d)

   straintot=0.d0
   do k=1,mpe
      straintot=straintot+h(k)*strain(icon(k))
   enddo

   e2dprev=0.d0
   do k=1,mpe
      e2dprev=e2dprev+h(k)*e2dp(icon(k))
   enddo

!   if (straintot.gt.plasticity_parameters(3)) &
!   write(*,'(4es15.6)') sum(x)/8.d0,sum(y)/8.d0,sum(z)/8.d0,straintot

   e2dref=1.d0
   viscosity=viscosity0      
   if (e2d.gt.eps) then
     if (expon.gt.1.d0+eps .or. expon.lt.1.d0-eps) viscosity=viscosity*(e2d/e2dref)**(1.d0/expon-1.d0)
   endif
   viscosity=viscosity*fviscosity*exp(activationenergy/(expon*Rgas*(temperature*params%tempscale+273.15d0)))

   friction_angle=0.d0
   if (trim(plasticity_type)/='No') &
   call vrm (params,J2d,J3d,compressibility,viscosity,pressure,friction_angle, &
            plasticity_parameters,plasticity_type,plasticity_ss_type_coh,      &
            plasticity_ss_type_phi,is_plastic_temp,flag_vrm_pb,straintot,      &
            e2dprev,yield)

   is_plastic=is_plastic.or.is_plastic_temp

   if (params%viscositymax.ge.0.d0) then
      viscosity=min(params%viscositymax,viscosity)
   end if
   if (params%viscositymin.ge.0.d0) then
      viscosity=max(params%viscositymin,viscosity)
   end if

   viscomean=viscomean+viscosity
   yieldnow=2.d0*viscosity*e2d
   yield_ratio=yield_ratio+yieldnow/yield
   frictmean=frictmean+friction_angle

   d=0.d0
   d(1,1)=viscosity*2.d0
   d(2,2)=viscosity*2.d0
   d(3,3)=viscosity*2.d0
   d(4,4)=viscosity
   d(5,5)=viscosity
   d(6,6)=viscosity

   do k=1,mpe
      k1=(k-1)*3+1
      k2=(k-1)*3+2
      k3=(k-1)*3+3
      b(k1,1)=dhdx(k)
      b(k2,1)=0.d0
      b(k3,1)=0.d0
      b(k1,2)=0.d0
      b(k2,2)=dhdy(k)
      b(k3,2)=0.d0
      b(k1,3)=0.d0
      b(k2,3)=0.d0
      b(k3,3)=dhdz(k)
      b(k1,4)=0.d0
      b(k2,4)=dhdz(k)
      b(k3,4)=dhdy(k)
      b(k1,5)=dhdz(k)
      b(k2,5)=0.d0
      b(k3,5)=dhdx(k)
      b(k1,6)=dhdy(k)
      b(k2,6)=dhdx(k)
      b(k3,6)=0.d0
   enddo

   do k=1,mpe
      k1=(k-1)*3+1
      k2=(k-1)*3+2
      k3=(k-1)*3+3
      bl(k1,1)=bl(k1,1)+dhdx(k)*volume*w
      bl(k2,2)=bl(k2,2)+dhdy(k)*volume*w
      bl(k3,3)=bl(k3,3)+dhdz(k)*volume*w
      bl(k2,4)=bl(k2,4)+dhdz(k)*volume*w
      bl(k3,4)=bl(k3,4)+dhdy(k)*volume*w
      bl(k1,5)=bl(k1,5)+dhdz(k)*volume*w
      bl(k3,5)=bl(k3,5)+dhdx(k)*volume*w
      bl(k1,6)=bl(k1,6)+dhdy(k)*volume*w
      bl(k2,6)=bl(k2,6)+dhdx(k)*volume*w
   enddo

   vol=vol+volume*w
   weight=weight+volume*w*density

   do j=1,mpe*ndof
      do i=1,6
         bd(j,i)=0.d0
         do k=1,6
            bd(j,i)=bd(j,i)+b(j,k)*d(k,i)
         enddo
      enddo
   enddo

   do j=1,mpe*ndof
      do i=1,mpe*ndof
         do k=1,6
            ael(i,j)=ael(i,j)+b(i,k)*bd(j,k)*volume*w
         enddo
      enddo
   enddo

   do j=1,mpe
      jj=(j-1)*ndof+ndof
      !bel(jj)=bel(jj)+h(j)*volume*w*density*(1.d0-expansion*temperature*params%tempscale)
      bel(jj)=bel(jj)+h(j)*volume*w*density
   enddo

enddo

! second add compressibility term (1 point integration)

viscomean=viscomean/nint
yield_ratio=yield_ratio/nint
frictmean=frictmean/nint

eviscosity=viscomean
frict_angle=frictmean

!viscomean=1.d0

dl=0.d0
if (params%bulkvisc) then
  dl(1,1)=compressibility
  dl(1,2)=compressibility
  dl(1,3)=compressibility
  dl(2,1)=compressibility
  dl(2,2)=compressibility
  dl(2,3)=compressibility
  dl(3,1)=compressibility
  dl(3,2)=compressibility
  dl(3,3)=compressibility
else
  dl(1,1)=compressibility*viscomean
  dl(1,2)=compressibility*viscomean
  dl(1,3)=compressibility*viscomean
  dl(2,1)=compressibility*viscomean
  dl(2,2)=compressibility*viscomean
  dl(2,3)=compressibility*viscomean
  dl(3,1)=compressibility*viscomean
  dl(3,2)=compressibility*viscomean
  dl(3,3)=compressibility*viscomean
endif

r=0.d0
s=0.d0
t=0.d0
w=8.d0
!added by Jean on 13/2/08
!   r=r0+rst*(r+1.d0)/2.d0
!   s=s0+rst*(s+1.d0)/2.d0
!   t=t0+rst*(t+1.d0)/2.d0
!   w=w
! up to here
! was wrong...
call compute_dhdx_dhdy_dhdz (mpe,r,s,t,x,y,z,dhdx,dhdy,dhdz,volume)

!----------effective viscosity--------------------------

!exx=0.d0
!eyy=0.d0
!ezz=0.d0
!exy=0.d0
!eyz=0.d0
!ezx=0.d0
!do k=1,mpe
!   ic=icon(k)
!   exx=exx+dhdx(k)*unode(ic)
!   eyy=eyy+dhdy(k)*vnode(ic)
!   ezz=ezz+dhdz(k)*wnode(ic)
!   exy=exy+(dhdx(k)*vnode(ic)+dhdy(k)*unode(ic))/2.d0
!   eyz=eyz+(dhdy(k)*wnode(ic)+dhdz(k)*vnode(ic))/2.d0
!   ezx=ezx+(dhdz(k)*unode(ic)+dhdx(k)*wnode(ic))/2.d0
!enddo

!call deviatoric (exx,eyy,ezz,exy,eyz,ezx)
!
!J2d=second_invariant (exx,eyy,ezz,exy,eyz,ezx)
!J3d=third_invariant (exx,eyy,ezz,exy,eyz,ezx)
!
!e2d=0.d0
!if (J2d.ne.0.d0) e2d=sqrt(J2d)
!
!e2dref=1.d0
!eviscosity=viscosity0      
!if (e2d.ne.0.d0 .and. expon.ne.1.d0) eviscosity=eviscosity*(e2d/e2dref)**(1.d0/expon-1.d0)
!viscosity=viscosity*exp(activationenergy/(expon*Rgas*(temperature*params%tempscale+273.15d0)))
!
!if (trim(plasticity_type)/='No') &
!call vrm(J2d,J3d,eviscosity,pressure,plasticity_parameters,plasticity_type,is_plastic_dummy,flag_vrm_pb)

!----------effective viscosity--------------------------


do k=1,mpe
   k1=(k-1)*3+1
   k2=(k-1)*3+2
   k3=(k-1)*3+3
   b(k1,1)=dhdx(k)
   b(k2,1)=0.d0
   b(k3,1)=0.d0
   b(k1,2)=0.d0
   b(k2,2)=dhdy(k)
   b(k3,2)=0.d0
   b(k1,3)=0.d0
   b(k2,3)=0.d0
   b(k3,3)=dhdz(k)
   b(k1,4)=0.d0
   b(k2,4)=dhdz(k)
   b(k3,4)=dhdy(k)
   b(k1,5)=dhdz(k)
   b(k2,5)=0.d0
   b(k3,5)=dhdx(k)
   b(k1,6)=dhdy(k)
   b(k2,6)=dhdx(k)
   b(k3,6)=0.d0
enddo

do j=1,mpe*ndof
   do i=1,6
      bd(j,i)=0.d0
      do k=1,6
         bd(j,i)=bd(j,i)+b(j,k)*dl(k,i)
      enddo
   enddo
enddo

do j=1,mpe*ndof
   do i=1,mpe*ndof
      do k=1,6
         ael(i,j)=ael(i,j)+bd(i,k)*b(j,k)*volume*w
      enddo
   enddo
enddo

!do i=1,mpe*ndof
!   do j=1,mpe*ndof
!      if (abs(ael(i,j)-ael(j,i))>eps) stop 'Ael po symmetric'
!   end do
!end do

   ! third add fixed velocity boundary conditions
   do ii=1,mpe
   do k=1,ndof
      ij=(icon(ii)-1)*ndof+k
      if (kfix(ij).eq.1) then
         if (k.eq.1) fixt=unode(icon(ii))
         if (k.eq.2) fixt=vnode(icon(ii))
         if (k.eq.3) fixt=wnode(icon(ii))
         i=(ii-1)*ndof+k
         aref=ael(i,i)
         do j=1,mpe*ndof
            bel(j)=bel(j)-ael(j,i)*fixt
            ael(i,j)=0.d0
            ael(j,i)=0.d0
         enddo
         ael(i,i)=aref
         bel(i)=aref*fixt

!         ael(i,i)=1.d0*penaltyg
!         bel(i)=fixt*penaltyg
      endif
   enddo
   enddo

goto 1112

1111  continue

!-------------------------------------------------------------------------------
! temperature calculations
!-------------------------------------------------------------------------------

alpha=1.d0

do iint=1,nint

   r=r0+rst*(rr(iint)+1.d0)/2.d0
   s=s0+rst*(ss(iint)+1.d0)/2.d0
   t=t0+rst*(tt(iint)+1.d0)/2.d0
   w=ww(iint)

   h(1)=(1.d0-r)*(1.d0-s)*(1.d0-t)/8.d0
   h(2)=(1.d0+r)*(1.d0-s)*(1.d0-t)/8.d0
   h(3)=(1.d0-r)*(1.d0+s)*(1.d0-t)/8.d0
   h(4)=(1.d0+r)*(1.d0+s)*(1.d0-t)/8.d0
   h(5)=(1.d0-r)*(1.d0-s)*(1.d0+t)/8.d0
   h(6)=(1.d0+r)*(1.d0-s)*(1.d0+t)/8.d0
   h(7)=(1.d0-r)*(1.d0+s)*(1.d0+t)/8.d0
   h(8)=(1.d0+r)*(1.d0+s)*(1.d0+t)/8.d0

   velox=0.d0
   veloy=0.d0
   veloz=0.d0
   do k=1,mpe
      ic=icon(k)
      velox=velox+h(k)*unode(ic)
      veloy=veloy+h(k)*vnode(ic)
      veloz=veloz+h(k)*wnode(ic)
   enddo
 
   dhdr(1)=-(1.d0-s)*(1.d0-t)/8.d0
   dhdr(2)=(1.d0-s)*(1.d0-t)/8.d0
   dhdr(3)=-(1.d0+s)*(1.d0-t)/8.d0
   dhdr(4)=(1.d0+s)*(1.d0-t)/8.d0
   dhdr(5)=-(1.d0-s)*(1.d0+t)/8.d0
   dhdr(6)=(1.d0-s)*(1.d0+t)/8.d0
   dhdr(7)=-(1.d0+s)*(1.d0+t)/8.d0
   dhdr(8)=(1.d0+s)*(1.d0+t)/8.d0
   dhds(1)=-(1.d0-r)*(1.d0-t)/8.d0
   dhds(2)=-(1.d0+r)*(1.d0-t)/8.d0
   dhds(3)=(1.d0-r)*(1.d0-t)/8.d0
   dhds(4)=(1.d0+r)*(1.d0-t)/8.d0
   dhds(5)=-(1.d0-r)*(1.d0+t)/8.d0
   dhds(6)=-(1.d0+r)*(1.d0+t)/8.d0
   dhds(7)=(1.d0-r)*(1.d0+t)/8.d0
   dhds(8)=(1.d0+r)*(1.d0+t)/8.d0
   dhdt(1)=-(1.d0-r)*(1.d0-s)/8.d0
   dhdt(2)=-(1.d0+r)*(1.d0-s)/8.d0
   dhdt(3)=-(1.d0-r)*(1.d0+s)/8.d0
   dhdt(4)=-(1.d0+r)*(1.d0+s)/8.d0
   dhdt(5)=(1.d0-r)*(1.d0-s)/8.d0
   dhdt(6)=(1.d0+r)*(1.d0-s)/8.d0
   dhdt(7)=(1.d0-r)*(1.d0+s)/8.d0
   dhdt(8)=(1.d0+r)*(1.d0+s)/8.d0

   jcb=0.
   do k=1,mpe
      jcb(1,1)=jcb(1,1)+dhdr(k)*x(k)
      jcb(1,2)=jcb(1,2)+dhdr(k)*y(k)
      jcb(1,3)=jcb(1,3)+dhdr(k)*z(k)
      jcb(2,1)=jcb(2,1)+dhds(k)*x(k)
      jcb(2,2)=jcb(2,2)+dhds(k)*y(k)
      jcb(2,3)=jcb(2,3)+dhds(k)*z(k)
      jcb(3,1)=jcb(3,1)+dhdt(k)*x(k)
      jcb(3,2)=jcb(3,2)+dhdt(k)*y(k)
      jcb(3,3)=jcb(3,3)+dhdt(k)*z(k)
   enddo

   volume=jcb(1,1)*jcb(2,2)*jcb(3,3) &
         +jcb(1,2)*jcb(2,3)*jcb(3,1) &
         +jcb(2,1)*jcb(3,2)*jcb(1,3) &
         -jcb(1,3)*jcb(2,2)*jcb(3,1) &
         -jcb(1,2)*jcb(2,1)*jcb(3,3) &
         -jcb(2,3)*jcb(3,2)*jcb(1,1)

   if (volume.le.eps) then
      if (iproc.eq.0) call write_error_vtk (ileaves,volume,icon,x,y,z)
      call stop_run('Element bow-tied or collapsed in make_matrix temperature$')
   endif

   jcbi(1,1)=(jcb(2,2)*jcb(3,3)-jcb(2,3)*jcb(3,2))/volume
   jcbi(2,1)=(jcb(2,3)*jcb(3,1)-jcb(2,1)*jcb(3,3))/volume
   jcbi(3,1)=(jcb(2,1)*jcb(3,2)-jcb(2,2)*jcb(3,1))/volume
   jcbi(1,2)=(jcb(1,3)*jcb(3,2)-jcb(1,2)*jcb(3,3))/volume
   jcbi(2,2)=(jcb(1,1)*jcb(3,3)-jcb(1,3)*jcb(3,1))/volume
   jcbi(3,2)=(jcb(1,2)*jcb(3,1)-jcb(1,1)*jcb(3,2))/volume
   jcbi(1,3)=(jcb(1,2)*jcb(2,3)-jcb(1,3)*jcb(2,2))/volume
   jcbi(2,3)=(jcb(1,3)*jcb(2,1)-jcb(1,1)*jcb(2,3))/volume
   jcbi(3,3)=(jcb(1,1)*jcb(2,2)-jcb(1,2)*jcb(2,1))/volume

   do k=1,mpe
      do i=1,3
         b(k,i)=jcbi(i,1)*dhdr(k)+jcbi(i,2)*dhds(k)+jcbi(i,3)*dhdt(k)
      enddo
   enddo

   vol=vol+volume*w

   uvwnorm=velox**2+veloy**2+veloz**2
   xmin=minval(x)
   xmax=maxval(x)
   dx=xmax-xmin
   ymin=minval(y)
   ymax=maxval(y)
   dy=ymax-ymin
   zmin=minval(z)
   zmax=maxval(z)
   dz=zmax-zmin
   if (uvwnorm.lt.tiny(uvwnorm)) then
      tau=0.d0
   else
      tau=(abs(velox)*dx+abs(veloy)*dy+abs(veloz)*dz)/sqrt(15.d0)/uvwnorm
   endif

   do k=1,mpe
      ht(k)=h(k)+tau*(velox*b(k,1)+veloy*b(k,2)+veloz*b(k,3))
   enddo

   do j=1,mpe
      bel(j)=bel(j)+h(j)*params%dt*heat*volume*w
      do i=1,mpe
         ! xcond=diffusivity*(b(i,1)*b(j,1)+b(i,2)*b(j,2)+b(i,3)*b(j,3))*volume
         xcond=diffusivity*(b(i,1)*b(j,1)+b(i,2)*b(j,2)+b(i,3)*b(j,3))*volume &
              +ht(i)*(velox*b(j,1)+veloy*b(j,2)+veloz*b(j,3))*volume
         xmass=h(i)*h(j)*volume
         ael(i,j)=ael(i,j)+(xmass+params%dt*xcond*alpha)*w
         aelp(i,j)=aelp(i,j)+(xmass-params%dt*xcond*(1.d0-alpha))*w
      enddo
   enddo
   
enddo

do i=1,mpe
   do j=1,mpe
      bel(i)=bel(i)+aelp(i,j)*temp(icon(j))
   enddo
enddo

! third add fixed temperature boundary conditions

do i=1,mpe
   ic=icon(i)
   if (kfix(ic).eq.1) then
      fixt=temp(ic)
      do j=1,mpe
         bel(j)=bel(j)-ael(j,i)*fixt
         ael(i,j)=0.d0
         ael(j,i)=0.d0
      enddo
      ael(i,i)=1.d0*vol
      bel(i)=fixt*vol
   endif
enddo

1112  continue

deallocate (x)
deallocate (y)
deallocate (z)
deallocate (be)
deallocate (tnode)
deallocate (h)
deallocate (dhdr)
deallocate (dhds)
deallocate (dhdt)
deallocate (ht)
deallocate (dhdx)
deallocate (dhdy)
deallocate (dhdz)
deallocate (b)
deallocate (bd)
deallocate (bl)
deallocate (aelp)
deallocate (vel)

return
end
!-------------------------------------------------------------------------------
!-------------------------------------------------------------------------------