Skip to content
Snippets Groups Projects
vrm.f90 13.3 KiB
Newer Older
!------------------------------------------------------------------------------|
!------------------------------------------------------------------------------|
!                                                                              |
!              ||===\\                                                         | 
!              ||    \\                                                        |
!              ||     ||   //==\\   ||  ||   //==||  ||/==\\                   |
!              ||     ||  ||    ||  ||  ||  ||   ||  ||    ||                  |
!              ||    //   ||    ||  ||  ||  ||   ||  ||                        |
!              ||===//     \\==//    \\==\\  \\==\\  ||                        |
!                                                                              |
!------------------------------------------------------------------------------|
!------------------------------------------------------------------------------|
!                                                                              |
!              VRM    Feb. 2007                                                |
!                                                                              |
!------------------------------------------------------------------------------|
!------------------------------------------------------------------------------|

subroutine vrm (params,J2d,J3d,compressibility,viscosity,pressure,             &
               friction_angle,plasticity_parameters,plasticity_type,           &
Dave Whipp's avatar
Dave Whipp committed
               plasticity_ss_type_coh,plasticity_ss_type_phi,is_plastic,       &

use constants
use invariants

!------------------------------------------------------------------------------|
!(((((((((((((((( Purpose of the routine  ))))))))))))))))))))))))))))))))))))))
!------------------------------------------------------------------------------|
!This routine computes the rescaled viscosity in the element in the case it is
!above the failure criterion

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

implicit none

double precision :: J2d,J3d,viscosity,pressure,compressibility,straintot
double precision :: e2dprev,friction_angle
character (len=8) plasticity_type
character (len=16) plasticity_ss_type_coh
character (len=16) plasticity_ss_type_phi
double precision plasticity_parameters(9)
logical ::  is_plastic,flag_vrm_pb

!------------------------------------------------------------------------------|
!(((((((((((((((( declaration of the subroutine internal variables )))))))))))))
!------------------------------------------------------------------------------|

Dave Whipp's avatar
Dave Whipp committed
double precision :: yield,fail,theta,c,phi,e2d,zeta,alpha,k,eps
double precision :: sin_theta,cos_theta,sin_phi,cos_phi
Dave Whipp's avatar
Dave Whipp committed

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

flag_vrm_pb = .false.
if (params%init_e2d) then
  e2d=params%e2d0
else
  e2d=sqrt(J2d)
endif

select case (trim(plasticity_type))

  case ('vM')
    yield = plasticity_parameters(1)
    ! Call strain softening routine if initial cohesion (yield) is positive
Dave Whipp's avatar
Dave Whipp committed
    if (plasticity_parameters(6).gt.0.d0) then
      yield = ss_new(yield,plasticity_ss_type_coh,plasticity_parameters(6),    &
                     plasticity_parameters(7),plasticity_parameters(8),        &
                     straintot,e2d,e2dprev,flag_vrm_pb)
    fail = 2.d0*viscosity*e2d-yield
    if (fail.gt.0.d0) then
      viscosity=0.5d0*yield/e2d
      is_plastic=.true.
    else
      is_plastic=.false.
    end if
  case ('DPI')
    ! Strain softening/hardening
    phi = plasticity_parameters(1)*pi/180.d0
Dave Whipp's avatar
Dave Whipp committed
    if (plasticity_parameters(6).gt.0.d0) then
      c = ss_new(c,plasticity_ss_type_coh,plasticity_parameters(6),            &
                 plasticity_parameters(7),plasticity_parameters(8),            &
                 straintot,e2d,e2dprev,flag_vrm_pb)
    endif
Dave Whipp's avatar
Dave Whipp committed
    if (plasticity_parameters(3).gt.0.d0) then
      phi = ss_new(phi,plasticity_ss_type_phi,plasticity_parameters(3),        &
                   plasticity_parameters(4),plasticity_parameters(5)*pi/180.d0,&
                   straintot,e2d,e2dprev,flag_vrm_pb)
    friction_angle = phi*180.d0/pi
    alpha = (6.d0*sin(phi))  /(sqrt(3.d0)*(3.d0-sin(phi)))
    k     = (6.d0*c*cos(phi))/(sqrt(3.d0)*(3.d0-sin(phi)))
    if (params%plastic_stress_min.lt.eps) then
      yield = max(c/compressibility,k-alpha*pressure)
    else
      yield = max(params%plastic_stress_min,k-alpha*pressure)
    if (fail.gt.0.d0) then
      viscosity=0.5d0*yield/e2d
      is_plastic=.true.
    else
      is_plastic=.false.
    end if
  case ('DPII')
    ! Strain softening/hardening
    phi = plasticity_parameters(1)*pi/180.d0
Dave Whipp's avatar
Dave Whipp committed
    if (plasticity_parameters(6).gt.0.d0) then
      c = ss_new(c,plasticity_ss_type_coh,plasticity_parameters(6),            &
                 plasticity_parameters(7),plasticity_parameters(8),            &
                 straintot,e2d,e2dprev,flag_vrm_pb)
    endif
Dave Whipp's avatar
Dave Whipp committed
    if (plasticity_parameters(3).gt.0.d0) then
      phi = ss_new(phi,plasticity_ss_type_phi,plasticity_parameters(3),        &
                   plasticity_parameters(4),plasticity_parameters(5)*pi/180.d0,&
                   straintot,e2d,e2dprev,flag_vrm_pb)
    friction_angle = phi*180.d0/pi
    alpha = (6.d0*sin(phi))  /(sqrt(3.d0)*(3.d0+sin(phi)))
    k     = (6.d0*c*cos(phi))/(sqrt(3.d0)*(3.d0+sin(phi)))
    if (params%plastic_stress_min.lt.eps) then
      yield = max(c/compressibility,k-alpha*pressure)
    else
      yield = max(params%plastic_stress_min,k-alpha*pressure)
    fail  = 2.d0*viscosity*e2d-yield 
    if (fail.gt.0.d0) then
      viscosity=0.5d0*yield/e2d
      is_plastic=.true.
    else
      is_plastic=.false.
    end if
  case ('DPIII')
    ! Strain softening/hardening
    phi = plasticity_parameters(1)*pi/180.d0
Dave Whipp's avatar
Dave Whipp committed
    if (plasticity_parameters(6).gt.0.d0) then
      c = ss_new(c,plasticity_ss_type_coh,plasticity_parameters(6),            &
                 plasticity_parameters(7),plasticity_parameters(8),            &
                 straintot,e2d,e2dprev,flag_vrm_pb)
    endif
Dave Whipp's avatar
Dave Whipp committed
    if (plasticity_parameters(3).gt.0.d0) then
      phi = ss_new(phi,plasticity_ss_type_phi,plasticity_parameters(3),        &
                   plasticity_parameters(4),plasticity_parameters(5)*pi/180.d0,&
                   straintot,e2d,e2dprev,flag_vrm_pb)
    friction_angle = phi*180.d0/pi
    k     = (6.d0*c*cos(phi))/(sqrt(3.d0)*3.d0)
    if (params%plastic_stress_min.lt.eps) then
	  yield = max(c/compressibility,k-alpha*pressure)
	else
	  yield = max(params%plastic_stress_min,k-alpha*pressure)
    fail  = 2.d0*viscosity*e2d-yield 
    if (fail.gt.0.d0) then
      viscosity=0.5d0*yield/e2d
      is_plastic=.true.
    else
      is_plastic=.false.
    end if

  case ('DPIV')
    ! Strain softening/hardening
    phi = plasticity_parameters(1)*pi/180.d0
Dave Whipp's avatar
Dave Whipp committed
    if (plasticity_parameters(6).gt.0.d0) then
      c = ss_new(c,plasticity_ss_type_coh,plasticity_parameters(6),            &
                 plasticity_parameters(7),plasticity_parameters(8),            &
                 straintot,e2d,e2dprev,flag_vrm_pb)
    endif
Dave Whipp's avatar
Dave Whipp committed
    if (plasticity_parameters(3).gt.0.d0) then
      phi = ss_new(phi,plasticity_ss_type_phi,plasticity_parameters(3),        &
                   plasticity_parameters(4),plasticity_parameters(5)*pi/180.d0,&
                   straintot,e2d,e2dprev,flag_vrm_pb)
    friction_angle = phi*180.d0/pi
    alpha = (3.d0*tan(phi))/(sqrt(9.d0+12.d0*(tan(phi))**2.d0))
    k     = (3.d0*c)/(sqrt(9.d0+12.d0*(tan(phi))**2.d0))
    if (params%plastic_stress_min.lt.eps) then
      yield = max(c/compressibility,k-alpha*pressure)
    else
      yield = max(params%plastic_stress_min,k-alpha*pressure)
    fail  = 2.d0*viscosity*e2d-yield 
    if (fail.gt.0.d0) then
      viscosity=0.5d0*yield/e2d
      is_plastic=.true.
    else
      is_plastic=.false.
    end if
  case ('DPV') 
    alpha = plasticity_parameters(1)
    k     = plasticity_parameters(2)
    if (params%plastic_stress_min.lt.eps) then
      yield = max(k/compressibility,k-alpha*pressure)
    else
      yield = max(params%plastic_stress_min,k-alpha*pressure)
    fail  = 2.d0*viscosity*e2d-yield 
    if (fail.gt.0.d0) then
      viscosity=0.5d0*yield/e2d
      is_plastic=.true.
    else
      is_plastic=.false.
    end if
    
  case ('DPVI')
    ! Strain softening/hardening
    phi = plasticity_parameters(1)*pi/180.d0
Dave Whipp's avatar
Dave Whipp committed
    if (plasticity_parameters(6).gt.0.d0) then
      c = ss_new(c,plasticity_ss_type_coh,plasticity_parameters(6),            &
                 plasticity_parameters(7),plasticity_parameters(8),            &
                 straintot,e2d,e2dprev,flag_vrm_pb)
    endif
Dave Whipp's avatar
Dave Whipp committed
    if (plasticity_parameters(3).gt.0.d0) then
      phi = ss_new(phi,plasticity_ss_type_phi,plasticity_parameters(3),        &
                   plasticity_parameters(4),plasticity_parameters(5)*pi/180.d0,&
                   straintot,e2d,e2dprev,flag_vrm_pb)
    friction_angle = phi*180.d0/pi
    if (params%plastic_stress_min.lt.eps) then
      yield = max(c/compressibility,k-alpha*pressure)
    else
      yield = max(params%plastic_stress_min,k-alpha*pressure)
    fail  = 2.d0*viscosity*e2d-yield
    if (fail.gt.0.d0) then
      viscosity=0.5d0*yield/e2d
      is_plastic=.true.
    else
      is_plastic=.false.
    end if

  case ('DPVII')
    phi = plasticity_parameters(1)*pi/180.d0
    c = plasticity_parameters(2)
Dave Whipp's avatar
Dave Whipp committed
    if (plasticity_parameters(6).gt.0.d0) then
      c = ss_new(c,plasticity_ss_type_coh,plasticity_parameters(6),            &
                 plasticity_parameters(7),plasticity_parameters(8),            &
                 straintot,e2d,e2dprev,flag_vrm_pb)
    endif
Dave Whipp's avatar
Dave Whipp committed
    if (plasticity_parameters(3).gt.0.d0) then
      phi = ss_new(phi,plasticity_ss_type_phi,plasticity_parameters(3),        &
                   plasticity_parameters(4),plasticity_parameters(5)*pi/180.d0,&
                   straintot,e2d,e2dprev,flag_vrm_pb)
    friction_angle = phi*180.d0/pi
    alpha = sin(phi)
    k     = c*cos(phi)
    if (params%plastic_stress_min.lt.eps) then
      yield = max(c/compressibility,k-alpha*pressure)
      yield = max(params%plastic_stress_min,k-alpha*pressure)
    fail  = 2.d0*viscosity*e2d-yield 
    if (fail.gt.0.d0) then
      viscosity=0.5d0*yield/e2d
      is_plastic=.true.
    else
      is_plastic=.false.
    end if

  case ('MC')
    if (params%init_e2d) then
      theta = 0.d0
    else
      theta = lode_angle(J2d,J3d)
    endif
    sin_theta = sin(theta)
    cos_theta = cos(theta)
    phi = plasticity_parameters(1)*pi/180.d0
Dave Whipp's avatar
Dave Whipp committed
    if (plasticity_parameters(6).gt.0.d0) then
      c = ss_new(c,plasticity_ss_type_coh,plasticity_parameters(6),            &
                 plasticity_parameters(7),plasticity_parameters(8),            &
                 straintot,e2d,e2dprev,flag_vrm_pb)
    endif
Dave Whipp's avatar
Dave Whipp committed
    if (plasticity_parameters(3).gt.0.d0) then
      phi = ss_new(phi,plasticity_ss_type_phi,plasticity_parameters(3),        &
                   plasticity_parameters(4),plasticity_parameters(5)*pi/180.d0,&
                   straintot,e2d,e2dprev,flag_vrm_pb)
    friction_angle = phi*180.d0/pi
    sin_phi = sin(phi)
    cos_phi = cos(phi)
    zeta    = cos_theta-sin_phi*sin_theta/sqrt3
    if (params%plastic_stress_min.lt.eps) then
      yield = max(c/compressibility,(-pressure*sin_phi+c*cos_phi)/zeta)
      yield = max(params%plastic_stress_min,(-pressure*sin_phi+c*cos_phi)/zeta)
    if (fail.gt.0.d0) then
      is_plastic=.true.
    else
      is_plastic=.false.
    end if

  case ('Tresca')
    if (params%init_e2d) then
      theta = 0.d0
    else
      theta = lode_angle(J2d,J3d)
    endif
    yield = plasticity_parameters(1)
    ! Call strain softening routine if initial cohesion (yield) is positive
Dave Whipp's avatar
Dave Whipp committed
    if (plasticity_parameters(6).gt.0.d0) then
      yield = ss_new(yield,plasticity_ss_type_coh,plasticity_parameters(6),    &
                 plasticity_parameters(7),plasticity_parameters(8),            &
                 straintot,e2d,e2dprev,flag_vrm_pb)
    cos_theta = cos(theta)
    fail  = 2.d0*viscosity*e2d*cos_theta-yield 
    if (fail.gt.0.d0) then
      viscosity=0.5d0*yield/e2d/cos_theta
      is_plastic=.true.
    else
      is_plastic=.false.
    end if

  case default
    print *,plasticity_type
    call stop_run('error in vrm$')
    flag_vrm_pb = .true.
!-------------------------------------------------------------------------------
Dave Whipp's avatar
Dave Whipp committed
!-------------------------------------------------------------------------------