Skip to content
Snippets Groups Projects
vrm.f90 14.6 KiB
Newer Older
  • Learn to ignore specific revisions
  • !------------------------------------------------------------------------------|
    !------------------------------------------------------------------------------|
    !                                                                              |
    !              ||===\\                                                         | 
    !              ||    \\                                                        |
    !              ||     ||   //==\\   ||  ||   //==||  ||/==\\                   |
    !              ||     ||  ||    ||  ||  ||  ||   ||  ||    ||                  |
    !              ||    //   ||    ||  ||  ||  ||   ||  ||                        |
    !              ||===//     \\==//    \\==\\  \\==\\  ||                        |
    !                                                                              |
    !------------------------------------------------------------------------------|
    !------------------------------------------------------------------------------|
    !                                                                              |
    !              VRM    Feb. 2007                                                |
    !                                                                              |
    !------------------------------------------------------------------------------|
    !------------------------------------------------------------------------------|
    
    
    subroutine vrm (params,J2d,J3d,compressibility,viscosity,pressure,             &
    
                   friction_angle,plasticity_parameters,plasticity_type,           &
                   plasticity_ss_type,is_plastic,flag_vrm_pb,straintot,yield)
    
    
    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
    
    character (len=8) plasticity_type
    
    character (len=12) plasticity_ss_type
    
    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,e2d,zeta,alpha,k
    
    double precision :: sin_theta,cos_theta,sin_phi,cos_phi
    double precision :: strain_soft_in,strain_soft_out,strain_soft_phi,phi,fact
    
    double precision :: strain_soft_yield,eps,strain_soft_ref
    
    !-------------------------------------------------------------------------------
    !-------------------------------------------------------------------------------
    
    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)
        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
        strain_soft_in=plasticity_parameters(3)
        if (strain_soft_in.gt.0.d0) then
          strain_soft_out=plasticity_parameters(4)
          strain_soft_phi=plasticity_parameters(5)*pi/180.d0
    
          ss_type: select case (trim(plasticity_ss_type))
            case ('tot_strain')
              strain_soft_ref=straintot
            case ('strain_rate')
              strain_soft_ref=e2d
            case default
              write (*,*) 'plasticity_ss_type: ',plasticity_ss_type
              call stop_run('error in vrm$')
              flag_vrm_pb = .true.
          end select ss_type
          if (strain_soft_ref.gt.strain_soft_in) then
            fact=(strain_soft_ref-strain_soft_in)/(strain_soft_out-strain_soft_in)
    
            fact=min(fact,1.d0)
            phi=phi+(strain_soft_phi-phi)*fact
    
        friction_angle = phi*180.d0/pi
    
        c     = plasticity_parameters(2)
        alpha = (2.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 ('DPII')
        ! Strain softening/hardening
        phi = plasticity_parameters(1)*pi/180.d0
        strain_soft_in=plasticity_parameters(3)
        if (strain_soft_in.gt.0.d0) then
          strain_soft_out=plasticity_parameters(4)
          strain_soft_phi=plasticity_parameters(5)*pi/180.d0
    
          ss_type: select case (trim(plasticity_ss_type))
            case ('tot_strain')
              strain_soft_ref=straintot
            case ('strain_rate')
              strain_soft_ref=e2d
            case default
              write (*,*) 'plasticity_ss_type: ',plasticity_ss_type
              call stop_run('error in vrm$')
              flag_vrm_pb = .true.
          end select ss_type
          if (strain_soft_ref.gt.strain_soft_in) then
            fact=(strain_soft_ref-strain_soft_in)/(strain_soft_out-strain_soft_in)
    
            fact=min(fact,1.d0)
            phi=phi+(strain_soft_phi-phi)*fact
          endif
        endif
    
        friction_angle = phi*180.d0/pi
    
        c     = plasticity_parameters(2)
        alpha = (2.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
        strain_soft_in=plasticity_parameters(3)
        if (strain_soft_in.gt.0.d0) then
          strain_soft_out=plasticity_parameters(4)
          strain_soft_phi=plasticity_parameters(5)*pi/180.d0
    
          ss_type: select case (trim(plasticity_ss_type))
            case ('tot_strain')
              strain_soft_ref=straintot
            case ('strain_rate')
              strain_soft_ref=e2d
            case default
              write (*,*) 'plasticity_ss_type: ',plasticity_ss_type
              call stop_run('error in vrm$')
              flag_vrm_pb = .true.
          end select ss_type
          if (strain_soft_ref.gt.strain_soft_in) then
            fact=(strain_soft_ref-strain_soft_in)/(strain_soft_out-strain_soft_in)
    
            fact=min(fact,1.d0)
            phi=phi+(strain_soft_phi-phi)*fact
    
        friction_angle = phi*180.d0/pi
    
        c     = plasticity_parameters(2)
        alpha = (2.d0*sin(phi))  /(sqrt(3.d0)*3.d0)
        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
        strain_soft_in=plasticity_parameters(3)
        if (strain_soft_in.gt.0.d0) then
    
          strain_soft_out=plasticity_parameters(4)
    
          strain_soft_phi=plasticity_parameters(5)*pi/180.d0
    
          ss_type: select case (trim(plasticity_ss_type))
            case ('tot_strain')
              strain_soft_ref=straintot
            case ('strain_rate')
              strain_soft_ref=e2d
            case default
              write (*,*) 'plasticity_ss_type: ',plasticity_ss_type
              call stop_run('error in vrm$')
              flag_vrm_pb = .true.
          end select ss_type
          if (strain_soft_ref.gt.strain_soft_in) then
            fact=(strain_soft_ref-strain_soft_in)/(strain_soft_out-strain_soft_in)
    
            fact=min(fact,1.d0)
            phi=phi+(strain_soft_phi-phi)*fact
          endif
    
        friction_angle = phi*180.d0/pi
    
        c     = plasticity_parameters(2)
        alpha = 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
        strain_soft_in=plasticity_parameters(3)
        if (strain_soft_in.gt.0.d0) then
          strain_soft_out=plasticity_parameters(4)
          strain_soft_phi=plasticity_parameters(5)*pi/180.d0
    
          ss_type: select case (trim(plasticity_ss_type))
            case ('tot_strain')
              strain_soft_ref=straintot
            case ('strain_rate')
              strain_soft_ref=e2d
            case default
              write (*,*) 'plasticity_ss_type: ',plasticity_ss_type
              call stop_run('error in vrm$')
              flag_vrm_pb = .true.
          end select ss_type
          if (strain_soft_ref.gt.strain_soft_in) then
            fact=(strain_soft_ref-strain_soft_in)/(strain_soft_out-strain_soft_in)
    
            fact=min(fact,1.d0)
            phi=phi+(strain_soft_phi-phi)*fact
    
        friction_angle = phi*180.d0/pi
    
        c     = plasticity_parameters(2)
    
        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
    
        strain_soft_in=plasticity_parameters(3)
        if (strain_soft_in.gt.0.d0) then
          strain_soft_out=plasticity_parameters(4)
          strain_soft_phi=plasticity_parameters(5)*pi/180.d0
    
          ss_type: select case (trim(plasticity_ss_type))
            case ('tot_strain')
              strain_soft_ref=straintot
            case ('strain_rate')
              strain_soft_ref=e2d
            case default
              write (*,*) 'plasticity_ss_type: ',plasticity_ss_type
              call stop_run('error in vrm$')
              flag_vrm_pb = .true.
          end select ss_type
          if (strain_soft_ref.gt.strain_soft_in) then
            fact=(strain_soft_ref-strain_soft_in)/(strain_soft_out-strain_soft_in)
    
            fact=min(fact,1.d0)
            phi=phi+(strain_soft_phi-phi)*fact
          endif
        endif
    
        friction_angle = phi*180.d0/pi
    
        c = plasticity_parameters(2)
    
        alpha = sin(phi)
        k     = c*cos(phi)
    
        if (params%plastic_stress_min.lt.eps) then
    
          yield = max(plasticity_parameters(2)/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 ('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
    
        strain_soft_in=plasticity_parameters(3)
        if (strain_soft_in.gt.0.d0) then
          strain_soft_out=plasticity_parameters(4)
          strain_soft_phi=plasticity_parameters(5)*pi/180.d0
    
          ss_type: select case (trim(plasticity_ss_type))
            case ('tot_strain')
              strain_soft_ref=straintot
            case ('strain_rate')
              strain_soft_ref=e2d
            case default
              write (*,*) 'plasticity_ss_type: ',plasticity_ss_type
              call stop_run('error in vrm$')
              flag_vrm_pb = .true.
          end select ss_type
          if (strain_soft_ref.gt.strain_soft_in) then
            fact=(strain_soft_ref-strain_soft_in)/(strain_soft_out-strain_soft_in)
    
            fact=min(fact,1.d0)
            phi=phi+(strain_soft_phi-phi)*fact
          endif
        endif
    
        friction_angle = phi*180.d0/pi
    
        c = plasticity_parameters(2)
        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(plasticity_parameters(2)/compressibility,(-pressure*sin_phi+c*cos_phi)/zeta)
        else
    
          yield = max(params%plastic_stress_min,(-pressure*sin_phi+c*cos_phi)/zeta)
    
        if (fail.gt.0.d0) then
    
    Dave Whipp's avatar
    Dave Whipp committed
          !if (viscosity.le.params%viscositymin) then
          !  write(*,*) 'VISCOSITY BELOW MINIMUM'
          !  write(*,*) 'viscosity: ',viscosity
          !  write(*,*) 'yield: ',yield
          !  write(*,*) 'fail: ',fail
          !  write(*,*) '-pressure: ',-pressure
          !  write(*,*) 'phi: ',phi
          !  write(*,*) 'sin_phi: ',sin_phi
          !  write(*,*) 'c: ',c
          !  write(*,*) 'cos_phi: ',cos_phi
          !  write(*,*) 'zeta: ',zeta
          !  write(*,*) 'e2d: ',e2d
          !endif
    
          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)
        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.
    
    return
    end
    !-------------------------------------------------------------------------------
    
    !-------------------------------------------------------------------------------