88       energy,qsourcesplit,active)
 
   95    integer, 
intent(in) :: ixI^L, ixO^L
 
   96    double precision, 
intent(in) :: qdt, x(ixI^S,1:ndim), wCT(ixI^S,1:nw)
 
   97    double precision, 
intent(inout) :: w(ixI^S,1:nw)
 
   98    logical, 
intent(in) :: energy,qsourcesplit
 
   99    logical, 
intent(inout) :: active
 
  101    double precision:: lambda(ixI^S,ndir,ndir),tmp(ixI^S),tmp2(ixI^S),v(ixI^S,ndir),vlambda(ixI^S,ndir)
 
  102    integer:: ix^L,idim,idir,jdir,iw
 
  106    if(qsourcesplit .eqv. 
vc_split) 
then 
  113        if({ iximin^
d>ixmin^
d .or. iximax^
d<ixmax^
d|.or.})&
 
  114          call mpistop(
"error for viscous source addition, 2 layers needed")
 
  119        if({ iximin^
d>ixmin^
d .or. iximax^
d<ixmax^
d|.or.})&
 
  120          call mpistop(
"error for viscous source addition"//&
 
  121          "requested fourth order gradients: 4 layers needed")
 
  127        v(ixi^s,idir)=wct(ixi^s,mom(idir))/wct(ixi^s,rho_)
 
  135      do idim=1,ndim; 
do idir=1,ndir
 
  138        tmp(ixi^s)=v(ixi^s,idir)
 
  142          if     (idim==
r_  .and. (idir==2 .or. idir==
phi_)) 
then 
  143            tmp(ixi^s) = tmp(ixi^s)/x(ixi^s,1)
 
  145          elseif (idim==2  .and. idir==
phi_) 
then 
  146            tmp(ixi^s)=tmp(ixi^s)/dsin(x(ixi^s,2))
 
  150        call gradient(tmp,ixi^l,ix^l,idim,tmp2)
 
  155          if (idim==
r_  .and. (idir==2 .or. idir==
phi_)) 
then 
  156            tmp2(ix^s) = tmp2(ix^s)*x(ix^s,1)
 
  158          elseif (idim==2  .and. idir==
phi_ ) 
then 
  159            tmp2(ix^s)=tmp2(ix^s)*dsin(x(ix^s,2))
 
  160          elseif (idim==2   .and. idir==2   ) 
then 
  161            tmp2(ix^s)=tmp2(ix^s)+v(ix^s,
r_)/x(ix^s,1)
 
  162          elseif (idim==
phi_.and. idir==
phi_) 
then 
  163            tmp2(ix^s)=tmp2(ix^s)+v(ix^s,
r_)/x(ix^s,1)+v(ix^s,2)/(x(ix^s,1)*dtan(x(ix^s,2)))
 
  167        lambda(ix^s,idim,idir)= lambda(ix^s,idim,idir)+ tmp2(ix^s)
 
  168        lambda(ix^s,idir,idim)= lambda(ix^s,idir,idim)+ tmp2(ix^s)
 
  172      lambda(ix^s,1:ndir,1:ndir)=lambda(ix^s,1:ndir,1:ndir)*
vc_mu*qdt
 
  180         tmp(ix^s)=tmp(ix^s)+lambda(ix^s,idir,idir)
 
  182      tmp(ix^s)=tmp(ix^s)/3.d0
 
  186         lambda(ix^s,idir,idir)=lambda(ix^s,idir,idir)-tmp(ix^s)
 
  193              tmp(ix^s)=lambda(ix^s,idir,idim)
 
  198                if (idim==
r_ .and. idir==
r_                 ) tmp(ix^s) = tmp(ix^s)*x(ix^s,1)**two
 
  199                if (idim==
r_ .and. (idir==2 .or. idir==
phi_)) tmp(ix^s) = tmp(ix^s)*x(ix^s,1)**3.d0
 
  201                if (idim==2  .and. (idir==
r_ .or. idir==2))   tmp(ix^s) = tmp(ix^s)*dsin(x(ix^s,2))
 
  202                if (idim==2  .and. idir==
phi_               ) tmp(ix^s) = tmp(ix^s)*dsin(x(ix^s,2))**two
 
  205              call gradient(tmp,ixi^l,ixo^l,idim,tmp2)
 
  210                if (idim==
r_ .and. idir==
r_                 ) tmp2(ixo^s) = tmp2(ixo^s)/(x(ixo^s,1)**two)
 
  211                if (idim==
r_ .and. (idir==2 .or. idir==
phi_)) tmp2(ixo^s) = tmp2(ixo^s)/(x(ixo^s,1)**3.d0)
 
  213                if (idim==2  .and. (idir==
r_ .or. idir==2))   tmp2(ixo^s) = tmp2(ixo^s)/(dsin(x(ixo^s,2)))
 
  214                if (idim==2  .and. idir==
phi_               ) tmp2(ixo^s) = tmp2(ixo^s)/(dsin(x(ixo^s,2))**two)
 
  217              w(ixo^s,mom(idir))=w(ixo^s,mom(idir))+tmp2(ixo^s)
 
  223        if (
coordinate==
spherical   .and. idir==2   ) w(ixo^s,mom(idir))=w(ixo^s,mom(idir))-lambda(ixo^s,
phi_,
phi_)/(x(ixo^s,1)/dtan(x(ixo^s,2)))
 
  233             vlambda(ixi^s,idim)=vlambda(ixi^s,idim)+v(ixi^s,idir)*lambda(ixi^s,idir,idim)
 
  237        w(ixo^s,e_)=w(ixo^s,e_)+tmp2(ixo^s)
 
 
  323    integer, 
intent(in)             :: ixI^L, ixO^L, idim
 
  324    double precision, 
intent(in)    :: w(ixI^S, 1:nw), x(ixI^S, 1:ndim)
 
  325    double precision, 
intent(out)   :: cross(ixI^S,ndir)
 
  327    double precision :: tmp(ixI^S), v(ixI^S)
 
  330    if (ndir/=ndim) 
call mpistop(
"This formula are probably wrong for ndim/=ndir")
 
  344        v(ixi^s)=w(ixi^s,mom(1)) 
 
  346        cross(ixi^s,2)=tmp(ixi^s)
 
  347        v(ixi^s)=w(ixi^s,mom(2))/x(ixi^s,1)  
 
  349        cross(ixi^s,2)=cross(ixi^s,2)+tmp(ixi^s)*x(ixi^s,1)
 
  351      elseif (idim==2) 
then 
  353        v(ixi^s)=w(ixi^s,mom(1)) 
 
  356        cross(ixi^s,1)=tmp(ixi^s)
 
  357        v(ixi^s)=w(ixi^s,mom(2))/x(ixi^s,1)  
 
  359        cross(ixi^s,1)=cross(ixi^s,1)+tmp(ixi^s)*x(ixi^s,1)
 
  361        v(ixi^s)=w(ixi^s,mom(2)) 
 
  363        cross(ixi^s,2)=two*tmp(ixi^s)
 
  364        v(ixi^s)=w(ixi^s,mom(1)) 
 
  365        cross(ixi^s,2)=cross(ixi^s,2)+two*v(ixi^s)/x(ixi^s,1) 
 
  367        v(ixi^s)=w(ixi^s,mom(3)) 
 
  370        cross(ixi^s,3)=tmp(ixi^s)
 
  371        v(ixi^s)=w(ixi^s,mom(2))  
 
  375        cross(ixi^s,3)=cross(ixi^s,3)+tmp(ixi^s)
 
  377      elseif (idim==3) 
then 
  382        v(ixi^s)=w(ixi^s,mom(3)) 
 
  384        cross(ixi^s,2)=tmp(ixi^s)
 
  385        v(ixi^s)=w(ixi^s,mom(2))  
 
  387        cross(ixi^s,2)=cross(ixi^s,2)+tmp(ixi^s)
 
  393        v(ixi^s)=w(ixi^s,mom(1)) 
 
  395        cross(ixi^s,1)=two*tmp(ixi^s)
 
  397        v(ixi^s)=w(ixi^s,mom(1)) 
 
  400        cross(ixi^s,2)=tmp(ixi^s)
 
  401        v(ixi^s)=w(ixi^s,mom(2))/x(ixi^s,1)  
 
  403        cross(ixi^s,2)=cross(ixi^s,2)+tmp(ixi^s)*x(ixi^s,1)
 
  407        v(ixi^s)=w(ixi^s,mom(1)) 
 
  409        cross(ixi^s,3)=tmp(ixi^s)
 
  410        v(ixi^s)=w(ixi^s,mom(3))/x(ixi^s,1) 
 
  412        cross(ixi^s,3)=cross(ixi^s,3)+tmp(ixi^s)*x(ixi^s,1)
 
  414      elseif (idim==2) 
then 
  416        v(ixi^s)=w(ixi^s,mom(1)) 
 
  419        cross(ixi^s,1)=tmp(ixi^s)
 
  420        v(ixi^s)=w(ixi^s,mom(2))/x(ixi^s,1)  
 
  422        cross(ixi^s,1)=cross(ixi^s,1)+tmp(ixi^s)*x(ixi^s,1)
 
  424        v(ixi^s)=w(ixi^s,mom(2)) 
 
  426        cross(ixi^s,2)=two*tmp(ixi^s)
 
  427        v(ixi^s)=w(ixi^s,mom(1)) 
 
  428        cross(ixi^s,2)=cross(ixi^s,2)+two*v(ixi^s)/x(ixi^s,1) 
 
  432        v(ixi^s)=w(ixi^s,mom(2)) 
 
  434        cross(ixi^s,3)=tmp(ixi^s)
 
  435        v(ixi^s)=w(ixi^s,mom(3))/dsin(x(ixi^s,2)) 
 
  437        cross(ixi^s,3)=cross(ixi^s,3)+tmp(ixi^s)*dsin(x(ixi^s,2))
 
  440      elseif (idim==3) 
then 
  442        v(ixi^s)=w(ixi^s,mom(1)) 
 
  444        cross(ixi^s,1)=tmp(ixi^s)
 
  445        v(ixi^s)=w(ixi^s,mom(3))/x(ixi^s,1) 
 
  447        cross(ixi^s,1)=cross(ixi^s,1)+tmp(ixi^s)*x(ixi^s,1)
 
  449        v(ixi^s)=w(ixi^s,mom(2)) 
 
  451        cross(ixi^s,2)=tmp(ixi^s)
 
  452        v(ixi^s)=w(ixi^s,mom(3))/dsin(x(ixi^s,2)) 
 
  454        cross(ixi^s,2)=cross(ixi^s,2)+tmp(ixi^s)*dsin(x(ixi^s,2))
 
  456        v(ixi^s)=w(ixi^s,mom(3)) 
 
  458        cross(ixi^s,3)=two*tmp(ixi^s)
 
  459        v(ixi^s)=w(ixi^s,mom(1)) 
 
  460        cross(ixi^s,3)=cross(ixi^s,3)+two*v(ixi^s)/x(ixi^s,1) 
 
  461        v(ixi^s)=w(ixi^s,mom(2)) 
 
  462        cross(ixi^s,3)=cross(ixi^s,3)+two*v(ixi^s)/(x(ixi^s,1)*dtan(x(ixi^s,2))) 
 
  466      call mpistop(
"Unknown geometry specified")
 
 
  500    integer, 
intent(in)             :: ixI^L, ixO^L
 
  501    double precision, 
intent(in)    :: qdt, x(ixI^S, 1:ndim)
 
  502    double precision, 
intent(inout) :: wCT(ixI^S, 1:nw), w(ixI^S, 1:nw)
 
  506    double precision :: vv(ixI^S), divergence(ixI^S)
 
  507    double precision :: tmp(ixI^S),tmp1(ixI^S)
 
  517      call gradient(wct(ixi^s,mom(2)),ixi^l,ixo^l,2,tmp1) 
 
  518      tmp(ixo^s)=two*(tmp1(ixo^s)+wct(ixo^s,mom(1))/x(ixo^s,1)) 
 
  520      call divvector(wct(ixi^s,mom(1:
ndir)),ixi^l,ixo^l,divergence)
 
  521      tmp(ixo^s) = tmp(ixo^s) - (2.d0/3.d0) * divergence(ixo^s)
 
  523      w(ixo^s,mom(1))=w(ixo^s,mom(1))-qdt*
vc_mu*tmp(ixo^s)/x(ixo^s,1)
 
  525      call gradient(wct(ixi^s,mom(1)),ixi^l,ixo^l,2,tmp) 
 
  526      vv(ixi^s)=wct(ixi^s,mom(2))/x(ixi^s,1)  
 
  527      call gradient(vv,ixi^l,ixo^l,1,tmp1) 
 
  528      tmp(ixo^s)=tmp(ixo^s)+tmp1(ixo^s)*x(ixo^s,1)
 
  530      w(ixo^s,mom(2))=w(ixo^s,mom(2))+qdt*
vc_mu*tmp(ixo^s)/x(ixo^s,1)
 
  536      call gradient(wct(ixi^s,mom(2)),ixi^l,ixo^l,2,tmp1) 
 
  537      tmp(ixo^s)=two*(tmp1(ixo^s)+wct(ixo^s,mom(1))/x(ixo^s,1)) 
 
  539      call divvector(wct(ixi^s,mom(1:
ndir)),ixi^l,ixo^l,divergence)
 
  540      tmp(ixo^s) = tmp(ixo^s) - (2.d0/3.d0) * divergence(ixo^s)
 
  542      w(ixo^s,mom(1))=w(ixo^s,mom(1))-qdt*
vc_mu*tmp(ixo^s)/x(ixo^s,1)
 
  547      call gradient(wct(ixi^s,mom(3)),ixi^l,ixo^l,3,tmp1) 
 
  548      tmp(ixo^s)=two*(tmp1(ixo^s)+wct(ixo^s,mom(1))/x(ixo^s,1)+wct(ixo^s,mom(2))/(x(ixo^s,1)*dtan(x(ixo^s,2)))) 
 
  551      tmp(ixo^s) = tmp(ixo^s) - (2.d0/3.d0) * divergence(ixo^s)
 
  553      w(ixo^s,mom(1))=w(ixo^s,mom(1))-qdt*
vc_mu*tmp(ixo^s)/x(ixo^s,1)
 
  556      w(ixo^s,mom(2))=w(ixo^s,mom(2))-qdt*
vc_mu*tmp(ixo^s)/(x(ixo^s,1)*dtan(x(ixo^s,2)))
 
  559      call gradient(wct(ixi^s,mom(1)),ixi^l,ixo^l,2,tmp) 
 
  560      vv(ixi^s)=wct(ixi^s,mom(2))/x(ixi^s,1)  
 
  561      call gradient(vv,ixi^l,ixo^l,1,tmp1) 
 
  562      tmp(ixo^s)=tmp(ixo^s)+tmp1(ixo^s)*x(ixo^s,1)
 
  564      w(ixo^s,mom(2))=w(ixo^s,mom(2))+qdt*
vc_mu*tmp(ixo^s)/x(ixo^s,1)
 
  567      call gradient(wct(ixi^s,mom(1)),ixi^l,ixo^l,3,tmp) 
 
  569      vv(ixi^s)=wct(ixi^s,mom(3))/x(ixi^s,1) 
 
  570      call gradient(vv,ixi^l,ixo^l,1,tmp1) 
 
  571      tmp(ixo^s)=tmp(ixo^s)+tmp1(ixo^s)*x(ixo^s,1)
 
  573      w(ixo^s,mom(3))=w(ixo^s,mom(3))+qdt*
vc_mu*tmp(ixo^s)/x(ixo^s,1)
 
  576      call gradient(wct(ixi^s,mom(2)),ixi^l,ixo^l,3,tmp) 
 
  579      vv(ixi^s)=wct(ixi^s,mom(3))/dsin(x(ixi^s,2)) 
 
  580      call gradient(vv,ixi^l,ixo^l,2,tmp1) 
 
  581      tmp(ixo^s)=tmp(ixo^s)+tmp1(ixo^s)*dsin(x(ixo^s,2))
 
  583      w(ixo^s,mom(3))=w(ixo^s,mom(3))+qdt*
vc_mu*tmp(ixo^s)/(x(ixo^s,1)*dtan(x(ixo^s,2)))