SPEED
MAKE_STRESS_TENSOR_DAMPED.f90 File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine make_stress_tensor_damped (nn, lambda, mu, duxdx, duydx, duzdx, duxdy, duydy, duzdy, duxdz, duydz, duzdz, strain_visc_xx, strain_visc_xy, strain_vi
 Computes the stress tensor.
 

Function/Subroutine Documentation

◆ make_stress_tensor_damped()

subroutine make_stress_tensor_damped ( integer*4  nn,
real*8, dimension(nn,nn,nn)  lambda,
real*8, dimension(nn,nn,nn)  mu,
real*8, dimension(nn,nn,nn)  duxdx,
real*8, dimension(nn,nn,nn)  duydx,
  duzdx,
real*8, dimension(nn,nn,nn)  duxdy,
real*8, dimension(nn,nn,nn)  duydy,
  duzdy,
real*8, dimension(nn,nn,nn)  duxdz,
real*8, dimension(nn,nn,nn)  duydz,
  duzdz,
real*8, dimension(nn,nn,nn,n_sls)  strain_visc_xx,
  strain_visc_xy,
  strain_vi 
)

Computes the stress tensor.

Author
Ilario Mazzieri
Date
November, 2014
Version
1.0
Parameters
[in]nnnumber of 1-D Legendre nodes
[in]lambdanodal values of Lame coefficient lambda
[in]munodal values of Lame coefficient mu
[in]duxdxnodal values for spatial derivatives of the displacement
[in]duydxnodal values for spatial derivatives of the displacement
[in]duzdxnodal values for spatial derivatives of the displacement
[in]duxdynodal values for spatial derivatives of the displacement
[in]duydynodal values for spatial derivatives of the displacement
[in]duzdynodal values for spatial derivatives of the displacement
[in]duxdznodal values for spatial derivatives of the displacement
[in]duydznodal values for spatial derivatives of the displacement
[in]duzdznodal values for spatial derivatives of the displacement
[in]strain_visc_xxnodal values for viscoelastic stress tensor
[in]strain_visc_xynodal values for viscoelastic stress tensor
[in]strain_visc_xznodal values for viscoelastic stress tensor
[in]strain_visc_yynodal values for viscoelastic stress tensor
[in]strain_visc_yznodal values for viscoelastic stress tensor
[in]strain_visc_zznodal values for viscoelastic stress tensor
[in]Ylambdaanelastic coefficients (1st Lame modulus)
[in]Ymuanelastic coefficients (2nd Lame modulus)
[in]ielocal element id
[out]sxxnodal values for the stress tensor
[out]syynodal values for the stress tensor
[out]szznodal values for the stress tensor
[out]syznodal values for the stress tensor
[out]szxnodal values for the stress tensor
[out]sxynodal values for the stress tensor

Definition at line 52 of file MAKE_STRESS_TENSOR_DAMPED.f90.

57 strain_visc_yy, strain_visc_yz, strain_visc_zz, &
58 ylambda,ymu, n_sls, &
59 sxx,syy,szz,&
60 syz,szx,sxy,ie)
61
62
63 implicit none
64
65 integer*4 :: nn,N_SLS
66 integer*4 :: p,q,r,ie
67
68 real*8, dimension(N_SLS) :: ylambda, ymu, temp
69 real*8, dimension(nn,nn,nn) :: lambda,mu
70 real*8, dimension(nn,nn,nn) :: sxx,syy,szz,syz,szx,sxy
71 real*8, dimension(nn,nn,nn) :: duxdx,duxdy,duxdz,duydx,duydy,duydz,duzdx,duzdy,duzdz
72 real*8, dimension(nn,nn,nn,N_SLS) :: strain_visc_xx, strain_visc_xy, strain_visc_xz, &
73 strain_visc_yy, strain_visc_yz, strain_visc_zz
74
75
76 do r = 1,nn
77 do q = 1,nn
78 do p = 1,nn
79
80
81 temp = strain_visc_xx(p,q,r,:) + strain_visc_yy(p,q,r,:) + strain_visc_zz(p,q,r,:)
82
83 sxx(p,q,r) = lambda(p,q,r) * (duxdx(p,q,r) +duydy(p,q,r) +duzdz(p,q,r) &
84 - dot_product(ylambda,temp)) &
85 +2.0d0*mu(p,q,r) * (duxdx(p,q,r) - dot_product(ymu,strain_visc_xx(p,q,r,:)))
86
87 syy(p,q,r) = lambda(p,q,r) * (duxdx(p,q,r) +duydy(p,q,r) +duzdz(p,q,r) &
88 - dot_product(ylambda,temp)) &
89 +2.0d0*mu(p,q,r) * (duydy(p,q,r) - dot_product(ymu,strain_visc_yy(p,q,r,:)))
90
91 szz(p,q,r) = lambda(p,q,r) * (duxdx(p,q,r) +duydy(p,q,r) +duzdz(p,q,r) &
92 - dot_product(ylambda,temp)) &
93 +2.0d0*mu(p,q,r) * (duzdz(p,q,r) - dot_product(ymu,strain_visc_zz(p,q,r,:)))
94
95 syz(p,q,r) = mu(p,q,r) * (duydz(p,q,r) + duzdy(p,q,r) - dot_product(ymu,strain_visc_yz(p,q,r,:)))
96 szx(p,q,r) = mu(p,q,r) * (duzdx(p,q,r) + duxdz(p,q,r) - dot_product(ymu,strain_visc_xz(p,q,r,:)))
97 sxy(p,q,r) = mu(p,q,r) * (duxdy(p,q,r) + duydx(p,q,r) - dot_product(ymu,strain_visc_xy(p,q,r,:)))
98
99 ! write(*,*) p,q,r
100 ! write(*,*) temp
101 ! write(*,*) DOT_PRODUCT(Ylambda,temp)
102 ! write(*,*) DOT_PRODUCT(Ymu,strain_visc_xx(p,q,r,:))
103 ! write(*,*) DOT_PRODUCT(Ymu,strain_visc_yy(p,q,r,:))
104 ! write(*,*) DOT_PRODUCT(Ymu,strain_visc_zz(p,q,r,:))
105 ! write(*,*) DOT_PRODUCT(Ymu,strain_visc_yz(p,q,r,:))
106 ! write(*,*) DOT_PRODUCT(Ymu,strain_visc_xz(p,q,r,:))
107 ! write(*,*) DOT_PRODUCT(Ymu,strain_visc_xy(p,q,r,:))
108 ! write(*,*) '--------------------------------------'
109 enddo
110 enddo
111 enddo
112
113 return
114