2012-10-11 20:19:12 +05:30
!--------------------------------------------------------------------------------------------------
! $Id$
!--------------------------------------------------------------------------------------------------
!> @author Franz Roters, Max-Planck-Institut für Eisenforschung GmbH
!> @author Philip Eisenlohr, Max-Planck-Institut für Eisenforschung GmbH
2013-07-01 11:40:42 +05:30
!> @brief material subroutine for phenomenological crystal plasticity formulation using a powerlaw
!! fitting
2012-10-11 20:19:12 +05:30
!--------------------------------------------------------------------------------------------------
2012-03-09 01:55:28 +05:30
module constitutive_phenopowerlaw
2013-07-01 11:40:42 +05:30
use prec , only : &
pReal , &
2013-09-18 19:37:55 +05:30
pInt
2009-07-22 21:37:19 +05:30
2012-03-09 01:55:28 +05:30
implicit none
2012-04-11 19:31:02 +05:30
private
2013-12-12 03:33:09 +05:30
integer ( pInt ) , dimension ( : ) , allocatable , public , protected :: &
2014-05-22 20:46:05 +05:30
#ifndef NEWSTATE
2012-03-09 01:55:28 +05:30
constitutive_phenopowerlaw_sizeDotState , &
constitutive_phenopowerlaw_sizeState , &
2014-05-22 20:46:05 +05:30
#endif
2014-03-09 02:20:31 +05:30
constitutive_phenopowerlaw_sizePostResults !< cumulative size of post results
2012-04-11 19:31:02 +05:30
2013-12-12 03:33:09 +05:30
integer ( pInt ) , dimension ( : , : ) , allocatable , target , public :: &
2013-07-01 11:40:42 +05:30
constitutive_phenopowerlaw_sizePostResult !< size of each post result output
2013-12-12 03:33:09 +05:30
character ( len = 64 ) , dimension ( : , : ) , allocatable , target , public :: &
2013-07-01 11:40:42 +05:30
constitutive_phenopowerlaw_output !< name of each post result output
2013-12-12 03:33:09 +05:30
integer ( pInt ) , dimension ( : ) , allocatable , private :: &
2012-10-11 20:19:12 +05:30
constitutive_phenopowerlaw_Noutput , & !< number of outputs per instance of this constitution
constitutive_phenopowerlaw_totalNslip , & !< no. of slip system used in simulation
constitutive_phenopowerlaw_totalNtwin !< no. of twin system used in simulation
2012-03-09 01:55:28 +05:30
2013-12-12 03:33:09 +05:30
integer ( pInt ) , dimension ( : , : ) , allocatable , private :: &
2012-10-11 20:19:12 +05:30
constitutive_phenopowerlaw_Nslip , & !< active number of slip systems per family (input parameter, per family)
constitutive_phenopowerlaw_Ntwin !< active number of twin systems per family (input parameter, per family)
2012-03-09 01:55:28 +05:30
2013-12-12 03:33:09 +05:30
real ( pReal ) , dimension ( : ) , allocatable , private :: &
2012-10-11 20:19:12 +05:30
constitutive_phenopowerlaw_gdot0_slip , & !< reference shear strain rate for slip (input parameter)
constitutive_phenopowerlaw_gdot0_twin , & !< reference shear strain rate for twin (input parameter)
2012-10-22 20:25:07 +05:30
constitutive_phenopowerlaw_n_slip , & !< stress exponent for slip (input parameter)
2013-07-01 11:40:42 +05:30
constitutive_phenopowerlaw_n_twin , & !< stress exponent for twin (input parameter)
2012-03-09 01:55:28 +05:30
2012-10-11 20:19:12 +05:30
constitutive_phenopowerlaw_spr , & !< push-up factor for slip saturation due to twinning
2012-03-09 01:55:28 +05:30
constitutive_phenopowerlaw_twinB , &
constitutive_phenopowerlaw_twinC , &
constitutive_phenopowerlaw_twinD , &
constitutive_phenopowerlaw_twinE , &
2012-11-14 15:52:34 +05:30
constitutive_phenopowerlaw_h0_SlipSlip , & !< reference hardening slip - slip (input parameter)
constitutive_phenopowerlaw_h0_SlipTwin , & !< reference hardening slip - twin (input parameter, no effect at the moment)
constitutive_phenopowerlaw_h0_TwinSlip , & !< reference hardening twin - slip (input parameter)
constitutive_phenopowerlaw_h0_TwinTwin , & !< reference hardening twin - twin (input parameter)
2012-03-09 01:55:28 +05:30
constitutive_phenopowerlaw_a_slip , &
2012-10-22 20:25:07 +05:30
constitutive_phenopowerlaw_aTolResistance , &
constitutive_phenopowerlaw_aTolShear , &
constitutive_phenopowerlaw_aTolTwinfrac
2012-03-09 01:55:28 +05:30
2013-12-12 03:33:09 +05:30
real ( pReal ) , dimension ( : , : ) , allocatable , private :: &
2013-07-01 11:40:42 +05:30
constitutive_phenopowerlaw_tau0_slip , & !< initial critical shear stress for slip (input parameter, per family)
constitutive_phenopowerlaw_tau0_twin , & !< initial critical shear stress for twin (input parameter, per family)
constitutive_phenopowerlaw_tausat_slip , & !< maximum critical shear stress for slip (input parameter, per family)
constitutive_phenopowerlaw_nonSchmidCoeff , &
2012-11-14 15:52:34 +05:30
constitutive_phenopowerlaw_interaction_SlipSlip , & !< interaction factors slip - slip (input parameter)
constitutive_phenopowerlaw_interaction_SlipTwin , & !< interaction factors slip - twin (input parameter)
constitutive_phenopowerlaw_interaction_TwinSlip , & !< interaction factors twin - slip (input parameter)
constitutive_phenopowerlaw_interaction_TwinTwin !< interaction factors twin - twin (input parameter)
2012-03-09 01:55:28 +05:30
2013-12-12 03:33:09 +05:30
real ( pReal ) , dimension ( : , : , : ) , allocatable , private :: &
2012-11-14 15:52:34 +05:30
constitutive_phenopowerlaw_hardeningMatrix_SlipSlip , &
constitutive_phenopowerlaw_hardeningMatrix_SlipTwin , &
constitutive_phenopowerlaw_hardeningMatrix_TwinSlip , &
2014-03-09 02:20:31 +05:30
constitutive_phenopowerlaw_hardeningMatrix_TwinTwin
2013-11-27 17:09:28 +05:30
enum , bind ( c )
2013-12-12 03:33:09 +05:30
enumerator :: undefined_ID , &
resistance_slip_ID , &
2013-11-27 17:09:28 +05:30
accumulatedshear_slip_ID , &
shearrate_slip_ID , &
resolvedstress_slip_ID , &
totalshear_ID , &
resistance_twin_ID , &
accumulatedshear_twin_ID , &
shearrate_twin_ID , &
resolvedstress_twin_ID , &
totalvolfrac_ID
end enum
2013-12-12 03:33:09 +05:30
integer ( kind ( undefined_ID ) ) , dimension ( : , : ) , allocatable , private :: &
2013-11-27 17:09:28 +05:30
constitutive_phenopowerlaw_outputID !< ID of each post result output
2012-04-11 19:31:02 +05:30
public :: &
constitutive_phenopowerlaw_init , &
2014-05-22 20:46:05 +05:30
#ifndef NEWSTATE
2013-07-01 11:40:42 +05:30
constitutive_phenopowerlaw_stateInit , &
2012-04-11 19:31:02 +05:30
constitutive_phenopowerlaw_aTolState , &
2014-05-22 20:46:05 +05:30
#endif
2013-07-01 11:40:42 +05:30
constitutive_phenopowerlaw_LpAndItsTangent , &
2012-04-11 19:31:02 +05:30
constitutive_phenopowerlaw_dotState , &
2013-07-01 11:40:42 +05:30
constitutive_phenopowerlaw_postResults
2012-04-11 19:31:02 +05:30
2014-05-22 20:46:05 +05:30
#ifdef NEWSTATE
private :: &
constitutive_phenopowerlaw_aTolState , &
constitutive_phenopowerlaw_stateInit
#endif
2012-03-09 01:55:28 +05:30
contains
2009-07-22 21:37:19 +05:30
2013-07-01 11:40:42 +05:30
2012-10-11 20:19:12 +05:30
!--------------------------------------------------------------------------------------------------
2013-07-01 11:40:42 +05:30
!> @brief module initialization
!> @details reads in material parameters, allocates arrays, and does sanity checks
2012-10-11 20:19:12 +05:30
!--------------------------------------------------------------------------------------------------
2013-12-12 03:33:09 +05:30
subroutine constitutive_phenopowerlaw_init ( fileUnit )
2012-10-11 20:19:12 +05:30
use , intrinsic :: iso_fortran_env ! to get compiler_version and compiler_options (at least for gfortran 4.6 at the moment)
2013-09-18 19:37:55 +05:30
use prec , only : &
tol_math_check
2013-07-01 11:40:42 +05:30
use debug , only : &
debug_level , &
debug_constitutive , &
debug_levelBasic
2013-12-19 14:19:47 +05:30
use math , only : &
math_Mandel3333to66 , &
math_Voigt66to3333
use IO , only : &
IO_read , &
IO_lc , &
IO_getTag , &
IO_isBlank , &
IO_stringPos , &
IO_stringValue , &
IO_floatValue , &
IO_intValue , &
IO_warning , &
IO_error , &
IO_timeStamp , &
IO_EOF
use material , only : &
homogenization_maxNgrains , &
phase_plasticity , &
phase_plasticityInstance , &
phase_Noutput , &
PLASTICITY_PHENOPOWERLAW_label , &
PLASTICITY_PHENOPOWERLAW_ID , &
2014-05-22 20:46:05 +05:30
material_phase , &
#ifdef NEWSTATE
plasticState , &
#endif
2013-12-19 14:19:47 +05:30
MATERIAL_partPhase
2013-02-08 21:25:53 +05:30
use lattice
2014-05-22 20:46:05 +05:30
use numerics , only : &
numerics_integrator
2009-10-21 18:40:12 +05:30
2012-03-09 01:55:28 +05:30
implicit none
2013-12-12 03:33:09 +05:30
integer ( pInt ) , intent ( in ) :: fileUnit
2013-07-01 11:40:42 +05:30
2013-10-09 11:42:16 +05:30
integer ( pInt ) , parameter :: MAXNCHUNKS = LATTICE_maxNinteraction + 1_pInt
2013-07-01 11:40:42 +05:30
integer ( pInt ) , dimension ( 1_pInt + 2_pInt * MAXNCHUNKS ) :: positions
integer ( pInt ) :: &
maxNinstance , &
2014-03-09 02:20:31 +05:30
instance , phase , j , k , f , o , &
2013-07-01 11:40:42 +05:30
Nchunks_SlipSlip , Nchunks_SlipTwin , Nchunks_TwinSlip , Nchunks_TwinTwin , &
2013-09-12 20:17:09 +05:30
Nchunks_SlipFamilies , Nchunks_TwinFamilies , Nchunks_nonSchmid , &
2014-03-09 02:20:31 +05:30
index_myFamily , index_otherFamily , &
2014-05-27 17:40:16 +05:30
mySize = 0_pInt , sizeState , sizeDotState
2013-07-01 11:40:42 +05:30
character ( len = 65536 ) :: &
tag = '' , &
2014-05-22 20:46:05 +05:30
line = ''
integer ( pInt ) :: NofMyPhase
2014-04-29 23:20:59 +05:30
real ( pReal ) , dimension ( : ) , allocatable :: tempPerSlip
2013-01-09 03:41:59 +05:30
2013-11-27 13:34:05 +05:30
write ( 6 , '(/,a)' ) ' <<<+- constitutive_' / / PLASTICITY_PHENOPOWERLAW_label / / ' init -+>>>'
2013-05-28 23:01:55 +05:30
write ( 6 , '(a)' ) ' $Id$'
write ( 6 , '(a15,a)' ) ' Current time: ' , IO_timeStamp ( )
2012-02-01 00:48:55 +05:30
#include "compilation_info.f90"
2009-07-22 21:37:19 +05:30
2013-11-27 13:34:05 +05:30
maxNinstance = int ( count ( phase_plasticity == PLASTICITY_PHENOPOWERLAW_ID ) , pInt )
2013-07-01 11:40:42 +05:30
if ( maxNinstance == 0_pInt ) return
2009-10-16 01:32:52 +05:30
2013-07-01 11:40:42 +05:30
if ( iand ( debug_level ( debug_constitutive ) , debug_levelBasic ) / = 0_pInt ) &
write ( 6 , '(a16,1x,i5,/)' ) '# instances:' , maxNinstance
2014-05-22 20:46:05 +05:30
#ifndef NEWSTATE
2013-12-12 03:33:09 +05:30
allocate ( constitutive_phenopowerlaw_sizeDotState ( maxNinstance ) , source = 0_pInt )
allocate ( constitutive_phenopowerlaw_sizeState ( maxNinstance ) , source = 0_pInt )
2014-05-22 20:46:05 +05:30
#endif
allocate ( constitutive_phenopowerlaw_sizePostResults ( maxNinstance ) , source = 0_pInt )
2013-12-12 03:33:09 +05:30
allocate ( constitutive_phenopowerlaw_sizePostResult ( maxval ( phase_Noutput ) , maxNinstance ) , &
source = 0_pInt )
2012-03-09 01:55:28 +05:30
allocate ( constitutive_phenopowerlaw_output ( maxval ( phase_Noutput ) , maxNinstance ) )
constitutive_phenopowerlaw_output = ''
2013-12-12 03:33:09 +05:30
allocate ( constitutive_phenopowerlaw_outputID ( maxval ( phase_Noutput ) , maxNinstance ) , source = undefined_ID )
allocate ( constitutive_phenopowerlaw_Noutput ( maxNinstance ) , source = 0_pInt )
allocate ( constitutive_phenopowerlaw_Nslip ( lattice_maxNslipFamily , maxNinstance ) , source = 0_pInt )
allocate ( constitutive_phenopowerlaw_Ntwin ( lattice_maxNtwinFamily , maxNinstance ) , source = 0_pInt )
allocate ( constitutive_phenopowerlaw_totalNslip ( maxNinstance ) , source = 0_pInt )
allocate ( constitutive_phenopowerlaw_totalNtwin ( maxNinstance ) , source = 0_pInt )
allocate ( constitutive_phenopowerlaw_gdot0_slip ( maxNinstance ) , source = 0.0_pReal )
allocate ( constitutive_phenopowerlaw_n_slip ( maxNinstance ) , source = 0.0_pReal )
allocate ( constitutive_phenopowerlaw_tau0_slip ( lattice_maxNslipFamily , maxNinstance ) , &
source = 0.0_pReal )
allocate ( constitutive_phenopowerlaw_tausat_slip ( lattice_maxNslipFamily , maxNinstance ) , &
source = 0.0_pReal )
allocate ( constitutive_phenopowerlaw_gdot0_twin ( maxNinstance ) , source = 0.0_pReal )
allocate ( constitutive_phenopowerlaw_n_twin ( maxNinstance ) , source = 0.0_pReal )
allocate ( constitutive_phenopowerlaw_tau0_twin ( lattice_maxNtwinFamily , maxNinstance ) , &
source = 0.0_pReal )
allocate ( constitutive_phenopowerlaw_spr ( maxNinstance ) , source = 0.0_pReal )
allocate ( constitutive_phenopowerlaw_twinB ( maxNinstance ) , source = 0.0_pReal )
allocate ( constitutive_phenopowerlaw_twinC ( maxNinstance ) , source = 0.0_pReal )
allocate ( constitutive_phenopowerlaw_twinD ( maxNinstance ) , source = 0.0_pReal )
allocate ( constitutive_phenopowerlaw_twinE ( maxNinstance ) , source = 0.0_pReal )
allocate ( constitutive_phenopowerlaw_h0_SlipSlip ( maxNinstance ) , source = 0.0_pReal )
allocate ( constitutive_phenopowerlaw_h0_SlipTwin ( maxNinstance ) , source = 0.0_pReal )
allocate ( constitutive_phenopowerlaw_h0_TwinSlip ( maxNinstance ) , source = 0.0_pReal )
allocate ( constitutive_phenopowerlaw_h0_TwinTwin ( maxNinstance ) , source = 0.0_pReal )
allocate ( constitutive_phenopowerlaw_interaction_SlipSlip ( lattice_maxNinteraction , maxNinstance ) , &
source = 0.0_pReal )
allocate ( constitutive_phenopowerlaw_interaction_SlipTwin ( lattice_maxNinteraction , maxNinstance ) , &
source = 0.0_pReal )
allocate ( constitutive_phenopowerlaw_interaction_TwinSlip ( lattice_maxNinteraction , maxNinstance ) , &
source = 0.0_pReal )
allocate ( constitutive_phenopowerlaw_interaction_TwinTwin ( lattice_maxNinteraction , maxNinstance ) , &
source = 0.0_pReal )
allocate ( constitutive_phenopowerlaw_a_slip ( maxNinstance ) , source = 0.0_pReal )
allocate ( constitutive_phenopowerlaw_aTolResistance ( maxNinstance ) , source = 0.0_pReal )
allocate ( constitutive_phenopowerlaw_aTolShear ( maxNinstance ) , source = 0.0_pReal )
allocate ( constitutive_phenopowerlaw_aTolTwinfrac ( maxNinstance ) , source = 0.0_pReal )
allocate ( constitutive_phenopowerlaw_nonSchmidCoeff ( lattice_maxNnonSchmid , maxNinstance ) , &
source = 0.0_pReal )
rewind ( fileUnit )
2014-03-09 02:20:31 +05:30
phase = 0_pInt
2013-12-19 14:19:47 +05:30
do while ( trim ( line ) / = IO_EOF . and . IO_lc ( IO_getTag ( line , '<' , '>' ) ) / = material_partPhase ) ! wind forward to <phase>
2013-12-12 03:33:09 +05:30
line = IO_read ( fileUnit )
2009-07-22 21:37:19 +05:30
enddo
2013-12-12 22:39:59 +05:30
2014-03-09 02:20:31 +05:30
parsingFile : do while ( trim ( line ) / = IO_EOF ) ! read through sections of phase part
2013-12-12 03:33:09 +05:30
line = IO_read ( fileUnit )
2012-10-11 20:19:12 +05:30
if ( IO_isBlank ( line ) ) cycle ! skip empty lines
2013-12-12 22:39:59 +05:30
if ( IO_getTag ( line , '<' , '>' ) / = '' ) then ! stop at next part
line = IO_read ( fileUnit , . true . ) ! reset IO_read
exit
endif
2014-03-09 02:20:31 +05:30
if ( IO_getTag ( line , '[' , ']' ) / = '' ) then ! next phase
phase = phase + 1_pInt ! advance phase section counter
if ( phase_plasticity ( phase ) == PLASTICITY_PHENOPOWERLAW_ID ) then
Nchunks_SlipFamilies = count ( lattice_NslipSystem ( : , phase ) > 0_pInt )
Nchunks_TwinFamilies = count ( lattice_NtwinSystem ( : , phase ) > 0_pInt )
Nchunks_SlipSlip = maxval ( lattice_interactionSlipSlip ( : , : , phase ) )
Nchunks_SlipTwin = maxval ( lattice_interactionSlipTwin ( : , : , phase ) )
Nchunks_TwinSlip = maxval ( lattice_interactionTwinSlip ( : , : , phase ) )
Nchunks_TwinTwin = maxval ( lattice_interactionTwinTwin ( : , : , phase ) )
Nchunks_nonSchmid = lattice_NnonSchmid ( phase )
2014-04-29 23:20:59 +05:30
if ( allocated ( tempPerSlip ) ) deallocate ( tempPerSlip )
allocate ( tempPerSlip ( Nchunks_SlipFamilies ) )
2013-06-12 01:46:40 +05:30
endif
2014-02-10 20:01:19 +05:30
cycle ! skip to next line
2009-07-22 21:37:19 +05:30
endif
2014-03-09 02:20:31 +05:30
if ( phase > 0_pInt ) then ; if ( phase_plasticity ( phase ) == PLASTICITY_PHENOPOWERLAW_ID ) then ! one of my phases. Do not short-circuit here (.and. between if-statements), it's not safe in Fortran
instance = phase_plasticityInstance ( phase ) ! which instance of my plasticity is present phase
2014-02-10 20:01:19 +05:30
positions = IO_stringPos ( line , MAXNCHUNKS )
tag = IO_lc ( IO_stringValue ( line , positions , 1_pInt ) ) ! extract key
select case ( tag )
case ( '(output)' )
select case ( IO_lc ( IO_stringValue ( line , positions , 2_pInt ) ) )
case ( 'resistance_slip' )
2014-02-28 15:48:40 +05:30
constitutive_phenopowerlaw_outputID ( constitutive_phenopowerlaw_Noutput ( instance ) , instance ) = resistance_slip_ID
2014-06-25 04:23:25 +05:30
constitutive_phenopowerlaw_Noutput ( instance ) = constitutive_phenopowerlaw_Noutput ( instance ) + 1_pInt
constitutive_phenopowerlaw_output ( constitutive_phenopowerlaw_Noutput ( instance ) , instance ) = &
IO_lc ( IO_stringValue ( line , positions , 2_pInt ) )
2014-02-10 20:01:19 +05:30
case ( 'accumulatedshear_slip' )
2014-02-28 15:48:40 +05:30
constitutive_phenopowerlaw_outputID ( constitutive_phenopowerlaw_Noutput ( instance ) , instance ) = accumulatedshear_slip_ID
2014-06-25 04:23:25 +05:30
constitutive_phenopowerlaw_Noutput ( instance ) = constitutive_phenopowerlaw_Noutput ( instance ) + 1_pInt
constitutive_phenopowerlaw_output ( constitutive_phenopowerlaw_Noutput ( instance ) , instance ) = &
IO_lc ( IO_stringValue ( line , positions , 2_pInt ) )
2014-02-10 20:01:19 +05:30
case ( 'shearrate_slip' )
2014-02-28 15:48:40 +05:30
constitutive_phenopowerlaw_outputID ( constitutive_phenopowerlaw_Noutput ( instance ) , instance ) = shearrate_slip_ID
2014-06-25 04:23:25 +05:30
constitutive_phenopowerlaw_Noutput ( instance ) = constitutive_phenopowerlaw_Noutput ( instance ) + 1_pInt
constitutive_phenopowerlaw_output ( constitutive_phenopowerlaw_Noutput ( instance ) , instance ) = &
IO_lc ( IO_stringValue ( line , positions , 2_pInt ) )
2014-02-10 20:01:19 +05:30
case ( 'resolvedstress_slip' )
2014-02-28 15:48:40 +05:30
constitutive_phenopowerlaw_outputID ( constitutive_phenopowerlaw_Noutput ( instance ) , instance ) = resolvedstress_slip_ID
2014-06-25 04:23:25 +05:30
constitutive_phenopowerlaw_Noutput ( instance ) = constitutive_phenopowerlaw_Noutput ( instance ) + 1_pInt
constitutive_phenopowerlaw_output ( constitutive_phenopowerlaw_Noutput ( instance ) , instance ) = &
IO_lc ( IO_stringValue ( line , positions , 2_pInt ) )
2014-02-10 20:01:19 +05:30
case ( 'totalshear' )
2014-02-28 15:48:40 +05:30
constitutive_phenopowerlaw_outputID ( constitutive_phenopowerlaw_Noutput ( instance ) , instance ) = totalshear_ID
2014-06-25 04:23:25 +05:30
constitutive_phenopowerlaw_Noutput ( instance ) = constitutive_phenopowerlaw_Noutput ( instance ) + 1_pInt
constitutive_phenopowerlaw_output ( constitutive_phenopowerlaw_Noutput ( instance ) , instance ) = &
IO_lc ( IO_stringValue ( line , positions , 2_pInt ) )
2014-02-10 20:01:19 +05:30
case ( 'resistance_twin' )
2014-02-28 15:48:40 +05:30
constitutive_phenopowerlaw_outputID ( constitutive_phenopowerlaw_Noutput ( instance ) , instance ) = resistance_twin_ID
2014-06-25 04:23:25 +05:30
constitutive_phenopowerlaw_Noutput ( instance ) = constitutive_phenopowerlaw_Noutput ( instance ) + 1_pInt
constitutive_phenopowerlaw_output ( constitutive_phenopowerlaw_Noutput ( instance ) , instance ) = &
IO_lc ( IO_stringValue ( line , positions , 2_pInt ) )
2014-02-10 20:01:19 +05:30
case ( 'accumulatedshear_twin' )
2014-02-28 15:48:40 +05:30
constitutive_phenopowerlaw_outputID ( constitutive_phenopowerlaw_Noutput ( instance ) , instance ) = accumulatedshear_twin_ID
2014-06-25 04:23:25 +05:30
constitutive_phenopowerlaw_Noutput ( instance ) = constitutive_phenopowerlaw_Noutput ( instance ) + 1_pInt
constitutive_phenopowerlaw_output ( constitutive_phenopowerlaw_Noutput ( instance ) , instance ) = &
IO_lc ( IO_stringValue ( line , positions , 2_pInt ) )
2014-02-10 20:01:19 +05:30
case ( 'shearrate_twin' )
2014-02-28 15:48:40 +05:30
constitutive_phenopowerlaw_outputID ( constitutive_phenopowerlaw_Noutput ( instance ) , instance ) = shearrate_twin_ID
2014-06-25 04:23:25 +05:30
constitutive_phenopowerlaw_Noutput ( instance ) = constitutive_phenopowerlaw_Noutput ( instance ) + 1_pInt
constitutive_phenopowerlaw_output ( constitutive_phenopowerlaw_Noutput ( instance ) , instance ) = &
IO_lc ( IO_stringValue ( line , positions , 2_pInt ) )
2014-02-10 20:01:19 +05:30
case ( 'resolvedstress_twin' )
2014-02-28 15:48:40 +05:30
constitutive_phenopowerlaw_outputID ( constitutive_phenopowerlaw_Noutput ( instance ) , instance ) = resolvedstress_twin_ID
2014-06-25 04:23:25 +05:30
constitutive_phenopowerlaw_Noutput ( instance ) = constitutive_phenopowerlaw_Noutput ( instance ) + 1_pInt
constitutive_phenopowerlaw_output ( constitutive_phenopowerlaw_Noutput ( instance ) , instance ) = &
IO_lc ( IO_stringValue ( line , positions , 2_pInt ) )
2014-02-10 20:01:19 +05:30
case ( 'totalvolfrac' )
2014-02-28 15:48:40 +05:30
constitutive_phenopowerlaw_outputID ( constitutive_phenopowerlaw_Noutput ( instance ) , instance ) = totalvolfrac_ID
2014-06-25 04:23:25 +05:30
constitutive_phenopowerlaw_Noutput ( instance ) = constitutive_phenopowerlaw_Noutput ( instance ) + 1_pInt
constitutive_phenopowerlaw_output ( constitutive_phenopowerlaw_Noutput ( instance ) , instance ) = &
IO_lc ( IO_stringValue ( line , positions , 2_pInt ) )
2014-02-10 20:01:19 +05:30
case default
2014-06-25 04:23:25 +05:30
2014-02-10 20:01:19 +05:30
end select
2014-04-29 23:20:59 +05:30
!--------------------------------------------------------------------------------------------------
2014-05-22 20:46:05 +05:30
! parameters depending on number of slip families
2014-02-10 20:01:19 +05:30
case ( 'nslip' )
2014-04-29 23:20:59 +05:30
if ( positions ( 1 ) < Nchunks_SlipFamilies + 1_pInt ) &
2014-02-10 20:01:19 +05:30
call IO_warning ( 50_pInt , ext_msg = trim ( tag ) / / ' (' / / PLASTICITY_PHENOPOWERLAW_label / / ')' )
2014-04-29 23:20:59 +05:30
if ( positions ( 1 ) > Nchunks_SlipFamilies + 1_pInt ) &
call IO_error ( 150_pInt , ext_msg = trim ( tag ) / / ' (' / / PLASTICITY_PHENOPOWERLAW_label / / ')' )
2014-02-10 20:01:19 +05:30
Nchunks_SlipFamilies = positions ( 1 ) - 1_pInt
do j = 1_pInt , Nchunks_SlipFamilies
2014-04-29 23:20:59 +05:30
constitutive_phenopowerlaw_Nslip ( j , instance ) = IO_intValue ( line , positions , 1_pInt + j )
2014-02-10 20:01:19 +05:30
enddo
2014-04-29 23:20:59 +05:30
case ( 'tausat_slip' , 'tau0_slip' )
2014-02-10 20:01:19 +05:30
do j = 1_pInt , Nchunks_SlipFamilies
2014-04-29 23:20:59 +05:30
tempPerSlip ( j ) = IO_floatValue ( line , positions , 1_pInt + j )
2014-02-10 20:01:19 +05:30
enddo
2014-04-29 23:20:59 +05:30
select case ( tag )
case ( 'tausat_slip' )
constitutive_phenopowerlaw_tausat_slip ( 1 : Nchunks_SlipFamilies , instance ) = tempPerSlip ( 1 : Nchunks_SlipFamilies )
case ( 'tau0_slip' )
constitutive_phenopowerlaw_tau0_slip ( 1 : Nchunks_SlipFamilies , instance ) = tempPerSlip ( 1 : Nchunks_SlipFamilies )
end select
!--------------------------------------------------------------------------------------------------
2014-05-22 20:46:05 +05:30
! parameters depending on number of twin families
2014-02-10 20:01:19 +05:30
case ( 'ntwin' )
2014-04-29 23:20:59 +05:30
if ( positions ( 1 ) < Nchunks_TwinFamilies + 1_pInt ) &
2014-06-17 20:54:44 +05:30
call IO_warning ( 51_pInt , ext_msg = trim ( tag ) / / ' (' / / PLASTICITY_PHENOPOWERLAW_label / / ')' )
2014-04-29 23:20:59 +05:30
if ( positions ( 1 ) > Nchunks_TwinFamilies + 1_pInt ) &
call IO_error ( 150_pInt , ext_msg = trim ( tag ) / / ' (' / / PLASTICITY_PHENOPOWERLAW_label / / ')' )
2014-02-10 20:01:19 +05:30
Nchunks_TwinFamilies = positions ( 1 ) - 1_pInt
do j = 1_pInt , Nchunks_TwinFamilies
2014-04-29 23:20:59 +05:30
constitutive_phenopowerlaw_Ntwin ( j , instance ) = IO_intValue ( line , positions , 1_pInt + j )
2014-02-10 20:01:19 +05:30
enddo
case ( 'tau0_twin' )
do j = 1_pInt , Nchunks_TwinFamilies
2014-02-28 15:48:40 +05:30
constitutive_phenopowerlaw_tau0_twin ( j , instance ) = IO_floatValue ( line , positions , 1_pInt + j )
2014-02-10 20:01:19 +05:30
enddo
2014-04-29 23:20:59 +05:30
!--------------------------------------------------------------------------------------------------
! parameters depending on number of interactions
case ( 'interaction_sliptwin' )
if ( positions ( 1 ) < 1_pInt + Nchunks_SlipTwin ) &
call IO_warning ( 52_pInt , ext_msg = trim ( tag ) / / ' (' / / PLASTICITY_PHENOPOWERLAW_label / / ')' )
do j = 1_pInt , Nchunks_SlipTwin
constitutive_phenopowerlaw_interaction_SlipTwin ( j , instance ) = IO_floatValue ( line , positions , 1_pInt + j )
enddo
case ( 'interaction_twinslip' )
if ( positions ( 1 ) < 1_pInt + Nchunks_TwinSlip ) &
call IO_warning ( 52_pInt , ext_msg = trim ( tag ) / / ' (' / / PLASTICITY_PHENOPOWERLAW_label / / ')' )
do j = 1_pInt , Nchunks_TwinSlip
constitutive_phenopowerlaw_interaction_TwinSlip ( j , instance ) = IO_floatValue ( line , positions , 1_pInt + j )
enddo
case ( 'interaction_twintwin' )
if ( positions ( 1 ) < 1_pInt + Nchunks_TwinTwin ) &
call IO_warning ( 52_pInt , ext_msg = trim ( tag ) / / ' (' / / PLASTICITY_PHENOPOWERLAW_label / / ')' )
do j = 1_pInt , Nchunks_TwinTwin
constitutive_phenopowerlaw_interaction_TwinTwin ( j , instance ) = IO_floatValue ( line , positions , 1_pInt + j )
enddo
case ( 'nonschmid_coefficients' )
if ( positions ( 1 ) < 1_pInt + Nchunks_nonSchmid ) &
call IO_warning ( 52_pInt , ext_msg = trim ( tag ) / / ' (' / / PLASTICITY_PHENOPOWERLAW_label / / ')' )
do j = 1_pInt , Nchunks_nonSchmid
constitutive_phenopowerlaw_nonSchmidCoeff ( j , instance ) = IO_floatValue ( line , positions , 1_pInt + j )
enddo
!--------------------------------------------------------------------------------------------------
! parameters independent of number of slip/twin systems
case ( 'gdot0_slip' )
constitutive_phenopowerlaw_gdot0_slip ( instance ) = IO_floatValue ( line , positions , 2_pInt )
case ( 'n_slip' )
constitutive_phenopowerlaw_n_slip ( instance ) = IO_floatValue ( line , positions , 2_pInt )
case ( 'a_slip' , 'w0_slip' )
constitutive_phenopowerlaw_a_slip ( instance ) = IO_floatValue ( line , positions , 2_pInt )
case ( 'gdot0_twin' )
constitutive_phenopowerlaw_gdot0_twin ( instance ) = IO_floatValue ( line , positions , 2_pInt )
case ( 'n_twin' )
constitutive_phenopowerlaw_n_twin ( instance ) = IO_floatValue ( line , positions , 2_pInt )
2014-02-10 20:01:19 +05:30
case ( 's_pr' )
2014-02-28 15:48:40 +05:30
constitutive_phenopowerlaw_spr ( instance ) = IO_floatValue ( line , positions , 2_pInt )
2014-02-10 20:01:19 +05:30
case ( 'twin_b' )
2014-02-28 15:48:40 +05:30
constitutive_phenopowerlaw_twinB ( instance ) = IO_floatValue ( line , positions , 2_pInt )
2014-02-10 20:01:19 +05:30
case ( 'twin_c' )
2014-02-28 15:48:40 +05:30
constitutive_phenopowerlaw_twinC ( instance ) = IO_floatValue ( line , positions , 2_pInt )
2014-02-10 20:01:19 +05:30
case ( 'twin_d' )
2014-02-28 15:48:40 +05:30
constitutive_phenopowerlaw_twinD ( instance ) = IO_floatValue ( line , positions , 2_pInt )
2014-02-10 20:01:19 +05:30
case ( 'twin_e' )
2014-02-28 15:48:40 +05:30
constitutive_phenopowerlaw_twinE ( instance ) = IO_floatValue ( line , positions , 2_pInt )
2014-02-10 20:01:19 +05:30
case ( 'h0_slipslip' )
2014-02-28 15:48:40 +05:30
constitutive_phenopowerlaw_h0_SlipSlip ( instance ) = IO_floatValue ( line , positions , 2_pInt )
2014-02-10 20:01:19 +05:30
case ( 'h0_sliptwin' )
2014-02-28 15:48:40 +05:30
constitutive_phenopowerlaw_h0_SlipTwin ( instance ) = IO_floatValue ( line , positions , 2_pInt )
2014-02-10 20:01:19 +05:30
call IO_warning ( 42_pInt , ext_msg = trim ( tag ) / / ' (' / / PLASTICITY_PHENOPOWERLAW_label / / ')' )
case ( 'h0_twinslip' )
2014-02-28 15:48:40 +05:30
constitutive_phenopowerlaw_h0_TwinSlip ( instance ) = IO_floatValue ( line , positions , 2_pInt )
2014-02-10 20:01:19 +05:30
case ( 'h0_twintwin' )
2014-02-28 15:48:40 +05:30
constitutive_phenopowerlaw_h0_TwinTwin ( instance ) = IO_floatValue ( line , positions , 2_pInt )
2014-02-10 20:01:19 +05:30
case ( 'atol_resistance' )
2014-02-28 15:48:40 +05:30
constitutive_phenopowerlaw_aTolResistance ( instance ) = IO_floatValue ( line , positions , 2_pInt )
2014-02-10 20:01:19 +05:30
case ( 'atol_shear' )
2014-02-28 15:48:40 +05:30
constitutive_phenopowerlaw_aTolShear ( instance ) = IO_floatValue ( line , positions , 2_pInt )
2014-02-10 20:01:19 +05:30
case ( 'atol_twinfrac' )
2014-02-28 15:48:40 +05:30
constitutive_phenopowerlaw_aTolTwinfrac ( instance ) = IO_floatValue ( line , positions , 2_pInt )
2014-02-10 20:01:19 +05:30
case ( 'interaction_slipslip' )
if ( positions ( 1 ) < 1_pInt + Nchunks_SlipSlip ) &
call IO_warning ( 52_pInt , ext_msg = trim ( tag ) / / ' (' / / PLASTICITY_PHENOPOWERLAW_label / / ')' )
do j = 1_pInt , Nchunks_SlipSlip
2014-02-28 15:48:40 +05:30
constitutive_phenopowerlaw_interaction_SlipSlip ( j , instance ) = IO_floatValue ( line , positions , 1_pInt + j )
2014-02-10 20:01:19 +05:30
enddo
case default
2014-06-25 04:23:25 +05:30
2014-02-10 20:01:19 +05:30
end select
endif ; endif
2014-03-09 02:20:31 +05:30
enddo parsingFile
2009-07-22 21:37:19 +05:30
2014-03-09 02:20:31 +05:30
sanityChecks : do phase = 1_pInt , size ( phase_plasticity )
myPhase : if ( phase_plasticity ( phase ) == PLASTICITY_phenopowerlaw_ID ) then
instance = phase_plasticityInstance ( phase )
constitutive_phenopowerlaw_Nslip ( 1 : lattice_maxNslipFamily , instance ) = &
2014-06-18 14:40:16 +05:30
min ( lattice_NslipSystem ( 1 : lattice_maxNslipFamily , phase ) , & ! limit active slip systems per family to min of available and requested
constitutive_phenopowerlaw_Nslip ( 1 : lattice_maxNslipFamily , instance ) )
2014-03-09 02:20:31 +05:30
constitutive_phenopowerlaw_Ntwin ( 1 : lattice_maxNtwinFamily , instance ) = &
2014-06-18 14:40:16 +05:30
min ( lattice_NtwinSystem ( 1 : lattice_maxNtwinFamily , phase ) , & ! limit active twin systems per family to min of available and requested
constitutive_phenopowerlaw_Ntwin ( : , instance ) )
2014-03-09 02:20:31 +05:30
constitutive_phenopowerlaw_totalNslip ( instance ) = sum ( constitutive_phenopowerlaw_Nslip ( : , instance ) ) ! how many slip systems altogether
constitutive_phenopowerlaw_totalNtwin ( instance ) = sum ( constitutive_phenopowerlaw_Ntwin ( : , instance ) ) ! how many twin systems altogether
if ( any ( constitutive_phenopowerlaw_tau0_slip ( : , instance ) < 0.0_pReal . and . &
constitutive_phenopowerlaw_Nslip ( : , instance ) > 0 ) ) &
call IO_error ( 211_pInt , el = instance , ext_msg = 'tau0_slip (' / / PLASTICITY_PHENOPOWERLAW_label / / ')' )
if ( constitutive_phenopowerlaw_gdot0_slip ( instance ) < = 0.0_pReal ) &
call IO_error ( 211_pInt , el = instance , ext_msg = 'gdot0_slip (' / / PLASTICITY_PHENOPOWERLAW_label / / ')' )
if ( constitutive_phenopowerlaw_n_slip ( instance ) < = 0.0_pReal ) &
call IO_error ( 211_pInt , el = instance , ext_msg = 'n_slip (' / / PLASTICITY_PHENOPOWERLAW_label / / ')' )
if ( any ( constitutive_phenopowerlaw_tausat_slip ( : , instance ) < = 0.0_pReal . and . &
constitutive_phenopowerlaw_Nslip ( : , instance ) > 0 ) ) &
call IO_error ( 211_pInt , el = instance , ext_msg = 'tausat_slip (' / / PLASTICITY_PHENOPOWERLAW_label / / ')' )
if ( any ( constitutive_phenopowerlaw_a_slip ( instance ) == 0.0_pReal . and . &
constitutive_phenopowerlaw_Nslip ( : , instance ) > 0 ) ) &
call IO_error ( 211_pInt , el = instance , ext_msg = 'a_slip (' / / PLASTICITY_PHENOPOWERLAW_label / / ')' )
if ( any ( constitutive_phenopowerlaw_tau0_twin ( : , instance ) < 0.0_pReal . and . &
constitutive_phenopowerlaw_Ntwin ( : , instance ) > 0 ) ) &
call IO_error ( 211_pInt , el = instance , ext_msg = 'tau0_twin (' / / PLASTICITY_PHENOPOWERLAW_label / / ')' )
if ( constitutive_phenopowerlaw_gdot0_twin ( instance ) < = 0.0_pReal . and . &
any ( constitutive_phenopowerlaw_Ntwin ( : , instance ) > 0 ) ) &
call IO_error ( 211_pInt , el = instance , ext_msg = 'gdot0_twin (' / / PLASTICITY_PHENOPOWERLAW_label / / ')' )
if ( constitutive_phenopowerlaw_n_twin ( instance ) < = 0.0_pReal . and . &
any ( constitutive_phenopowerlaw_Ntwin ( : , instance ) > 0 ) ) &
call IO_error ( 211_pInt , el = instance , ext_msg = 'n_twin (' / / PLASTICITY_PHENOPOWERLAW_label / / ')' )
if ( constitutive_phenopowerlaw_aTolResistance ( instance ) < = 0.0_pReal ) &
constitutive_phenopowerlaw_aTolResistance ( instance ) = 1.0_pReal ! default absolute tolerance 1 Pa
if ( constitutive_phenopowerlaw_aTolShear ( instance ) < = 0.0_pReal ) &
constitutive_phenopowerlaw_aTolShear ( instance ) = 1.0e-6_pReal ! default absolute tolerance 1e-6
if ( constitutive_phenopowerlaw_aTolTwinfrac ( instance ) < = 0.0_pReal ) &
constitutive_phenopowerlaw_aTolTwinfrac ( instance ) = 1.0e-6_pReal ! default absolute tolerance 1e-6
endif myPhase
2013-07-01 11:40:42 +05:30
enddo sanityChecks
2009-10-22 14:28:14 +05:30
2013-07-12 12:27:15 +05:30
!--------------------------------------------------------------------------------------------------
! allocation of variables whose size depends on the total number of active slip systems
2012-11-14 15:52:34 +05:30
allocate ( constitutive_phenopowerlaw_hardeningMatrix_SlipSlip ( maxval ( constitutive_phenopowerlaw_totalNslip ) , & ! slip resistance from slip activity
2009-07-22 21:37:19 +05:30
maxval ( constitutive_phenopowerlaw_totalNslip ) , &
2013-12-12 22:39:59 +05:30
maxNinstance ) , source = 0.0_pReal )
2012-11-14 15:52:34 +05:30
allocate ( constitutive_phenopowerlaw_hardeningMatrix_SlipTwin ( maxval ( constitutive_phenopowerlaw_totalNslip ) , & ! slip resistance from twin activity
2011-03-24 22:50:35 +05:30
maxval ( constitutive_phenopowerlaw_totalNtwin ) , &
2013-12-12 22:39:59 +05:30
maxNinstance ) , source = 0.0_pReal )
2012-11-14 15:52:34 +05:30
allocate ( constitutive_phenopowerlaw_hardeningMatrix_TwinSlip ( maxval ( constitutive_phenopowerlaw_totalNtwin ) , & ! twin resistance from slip activity
maxval ( constitutive_phenopowerlaw_totalNslip ) , &
2013-12-12 22:39:59 +05:30
maxNinstance ) , source = 0.0_pReal )
2012-11-14 15:52:34 +05:30
allocate ( constitutive_phenopowerlaw_hardeningMatrix_TwinTwin ( maxval ( constitutive_phenopowerlaw_totalNtwin ) , & ! twin resistance from twin activity
2009-07-22 21:37:19 +05:30
maxval ( constitutive_phenopowerlaw_totalNtwin ) , &
2013-12-12 22:39:59 +05:30
maxNinstance ) , source = 0.0_pReal )
2009-10-16 01:32:52 +05:30
2014-03-09 02:20:31 +05:30
initializeInstances : do phase = 1_pInt , size ( phase_plasticity )
if ( phase_plasticity ( phase ) == PLASTICITY_phenopowerlaw_ID ) then
2014-05-22 20:46:05 +05:30
NofMyPhase = count ( material_phase == phase )
2014-03-09 02:20:31 +05:30
instance = phase_plasticityInstance ( phase )
outputsLoop : do o = 1_pInt , constitutive_phenopowerlaw_Noutput ( instance )
select case ( constitutive_phenopowerlaw_outputID ( o , instance ) )
case ( resistance_slip_ID , &
shearrate_slip_ID , &
accumulatedshear_slip_ID , &
resolvedstress_slip_ID &
)
mySize = constitutive_phenopowerlaw_totalNslip ( instance )
case ( resistance_twin_ID , &
shearrate_twin_ID , &
accumulatedshear_twin_ID , &
resolvedstress_twin_ID &
)
mySize = constitutive_phenopowerlaw_totalNtwin ( instance )
case ( totalshear_ID , &
totalvolfrac_ID &
)
mySize = 1_pInt
case default
end select
outputFound : if ( mySize > 0_pInt ) then
constitutive_phenopowerlaw_sizePostResult ( o , instance ) = mySize
constitutive_phenopowerlaw_sizePostResults ( instance ) = constitutive_phenopowerlaw_sizePostResults ( instance ) + mySize
endif outputFound
enddo outputsLoop
2014-05-22 20:46:05 +05:30
#ifndef NEWSTATE
2014-03-09 02:20:31 +05:30
constitutive_phenopowerlaw_sizeDotState ( instance ) = constitutive_phenopowerlaw_totalNslip ( instance ) + &
constitutive_phenopowerlaw_totalNtwin ( instance ) + &
2_pInt + &
constitutive_phenopowerlaw_totalNslip ( instance ) + &
constitutive_phenopowerlaw_totalNtwin ( instance ) ! s_slip, s_twin, sum(gamma), sum(f), accshear_slip, accshear_twin
constitutive_phenopowerlaw_sizeState ( instance ) = constitutive_phenopowerlaw_sizeDotState ( instance )
2014-05-22 20:46:05 +05:30
#else
2014-05-27 17:40:16 +05:30
sizeState = constitutive_phenopowerlaw_totalNslip ( instance ) + &
2014-05-22 20:46:05 +05:30
constitutive_phenopowerlaw_totalNtwin ( instance ) + &
2_pInt + &
constitutive_phenopowerlaw_totalNslip ( instance ) + &
constitutive_phenopowerlaw_totalNtwin ( instance ) ! s_slip, s_twin, sum(gamma), sum(f), accshear_slip, accshear_twin
2014-05-27 20:16:03 +05:30
plasticState ( phase ) % sizeState = sizeState
2014-06-11 22:22:18 +05:30
sizeDotState = sizeState
2014-05-27 20:16:03 +05:30
plasticState ( phase ) % sizeDotState = sizeState
2014-05-27 17:40:16 +05:30
allocate ( plasticState ( phase ) % aTolState ( sizeState ) , source = 0.0_pReal )
allocate ( plasticState ( phase ) % state0 ( sizeState , NofMyPhase ) , source = 0.0_pReal )
allocate ( plasticState ( phase ) % partionedState0 ( sizeState , NofMyPhase ) , source = 0.0_pReal )
allocate ( plasticState ( phase ) % subState0 ( sizeState , NofMyPhase ) , source = 0.0_pReal )
allocate ( plasticState ( phase ) % state ( sizeState , NofMyPhase ) , source = 0.0_pReal )
allocate ( plasticState ( phase ) % state_backup ( sizeState , NofMyPhase ) , source = 0.0_pReal )
2014-05-22 20:46:05 +05:30
2014-06-11 22:22:18 +05:30
allocate ( plasticState ( phase ) % dotState ( sizeDotState , NofMyPhase ) , source = 0.0_pReal )
allocate ( plasticState ( phase ) % dotState_backup ( sizeDotState , NofMyPhase ) , source = 0.0_pReal )
2014-05-22 20:46:05 +05:30
if ( any ( numerics_integrator == 1_pInt ) ) then
2014-06-11 22:22:18 +05:30
allocate ( plasticState ( phase ) % previousDotState ( sizeDotState , NofMyPhase ) , source = 0.0_pReal )
allocate ( plasticState ( phase ) % previousDotState2 ( sizeDotState , NofMyPhase ) , source = 0.0_pReal )
2014-05-22 20:46:05 +05:30
endif
if ( any ( numerics_integrator == 4_pInt ) ) &
2014-06-11 22:22:18 +05:30
allocate ( plasticState ( phase ) % RK4dotState ( sizeDotState , NofMyPhase ) , source = 0.0_pReal )
2014-05-22 20:46:05 +05:30
if ( any ( numerics_integrator == 5_pInt ) ) &
2014-06-11 22:22:18 +05:30
allocate ( plasticState ( phase ) % RKCK45dotState ( 6 , sizeDotState , NofMyPhase ) , source = 0.0_pReal )
2014-05-22 20:46:05 +05:30
#endif
2014-05-08 20:25:19 +05:30
do f = 1_pInt , lattice_maxNslipFamily ! >>> interaction slip -- X
2014-03-09 02:20:31 +05:30
index_myFamily = sum ( constitutive_phenopowerlaw_Nslip ( 1 : f - 1_pInt , instance ) )
do j = 1_pInt , constitutive_phenopowerlaw_Nslip ( f , instance ) ! loop over (active) systems in my family (slip)
do o = 1_pInt , lattice_maxNslipFamily
index_otherFamily = sum ( constitutive_phenopowerlaw_Nslip ( 1 : o - 1_pInt , instance ) )
do k = 1_pInt , constitutive_phenopowerlaw_Nslip ( o , instance ) ! loop over (active) systems in other family (slip)
constitutive_phenopowerlaw_hardeningMatrix_SlipSlip ( index_myFamily + j , index_otherFamily + k , instance ) = &
constitutive_phenopowerlaw_interaction_SlipSlip ( lattice_interactionSlipSlip ( &
sum ( lattice_NslipSystem ( 1 : f - 1 , phase ) ) + j , &
sum ( lattice_NslipSystem ( 1 : o - 1 , phase ) ) + k , &
phase ) , instance )
enddo ; enddo
do o = 1_pInt , lattice_maxNtwinFamily
index_otherFamily = sum ( constitutive_phenopowerlaw_Ntwin ( 1 : o - 1_pInt , instance ) )
do k = 1_pInt , constitutive_phenopowerlaw_Ntwin ( o , instance ) ! loop over (active) systems in other family (twin)
constitutive_phenopowerlaw_hardeningMatrix_SlipTwin ( index_myFamily + j , index_otherFamily + k , instance ) = &
constitutive_phenopowerlaw_interaction_SlipTwin ( lattice_interactionSlipTwin ( &
sum ( lattice_NslipSystem ( 1 : f - 1_pInt , phase ) ) + j , &
sum ( lattice_NtwinSystem ( 1 : o - 1_pInt , phase ) ) + k , &
phase ) , instance )
enddo ; enddo
enddo ; enddo
2014-05-08 20:25:19 +05:30
do f = 1_pInt , lattice_maxNtwinFamily ! >>> interaction twin -- X
2014-03-09 02:20:31 +05:30
index_myFamily = sum ( constitutive_phenopowerlaw_Ntwin ( 1 : f - 1_pInt , instance ) )
do j = 1_pInt , constitutive_phenopowerlaw_Ntwin ( f , instance ) ! loop over (active) systems in my family (twin)
do o = 1_pInt , lattice_maxNslipFamily
index_otherFamily = sum ( constitutive_phenopowerlaw_Nslip ( 1 : o - 1_pInt , instance ) )
do k = 1_pInt , constitutive_phenopowerlaw_Nslip ( o , instance ) ! loop over (active) systems in other family (slip)
constitutive_phenopowerlaw_hardeningMatrix_TwinSlip ( index_myFamily + j , index_otherFamily + k , instance ) = &
constitutive_phenopowerlaw_interaction_TwinSlip ( lattice_interactionTwinSlip ( &
sum ( lattice_NtwinSystem ( 1 : f - 1_pInt , phase ) ) + j , &
sum ( lattice_NslipSystem ( 1 : o - 1_pInt , phase ) ) + k , &
phase ) , instance )
enddo ; enddo
do o = 1_pInt , lattice_maxNtwinFamily
index_otherFamily = sum ( constitutive_phenopowerlaw_Ntwin ( 1 : o - 1_pInt , instance ) )
do k = 1_pInt , constitutive_phenopowerlaw_Ntwin ( o , instance ) ! loop over (active) systems in other family (twin)
constitutive_phenopowerlaw_hardeningMatrix_TwinTwin ( index_myFamily + j , index_otherFamily + k , instance ) = &
constitutive_phenopowerlaw_interaction_TwinTwin ( lattice_interactionTwinTwin ( &
sum ( lattice_NtwinSystem ( 1 : f - 1_pInt , phase ) ) + j , &
sum ( lattice_NtwinSystem ( 1 : o - 1_pInt , phase ) ) + k , &
phase ) , instance )
enddo ; enddo
enddo ; enddo
2014-05-22 20:46:05 +05:30
#ifdef NEWSTATE
call constitutive_phenopowerlaw_stateInit ( phase , instance )
call constitutive_phenopowerlaw_aTolState ( phase , instance )
#endif
2014-03-09 02:20:31 +05:30
endif
enddo initializeInstances
2009-07-22 21:37:19 +05:30
2012-03-09 01:55:28 +05:30
end subroutine constitutive_phenopowerlaw_init
2009-07-22 21:37:19 +05:30
2014-05-22 20:46:05 +05:30
#ifdef NEWSTATE
!--------------------------------------------------------------------------------------------------
!> @brief sets the initial microstructural state for a given instance of this plasticity
!--------------------------------------------------------------------------------------------------
subroutine constitutive_phenopowerlaw_stateInit ( phase , instance )
use lattice , only : &
lattice_maxNslipFamily , &
lattice_maxNtwinFamily
use material , only : &
plasticState
implicit none
integer ( pInt ) , intent ( in ) :: &
instance , & !< number specifying the instance of the plasticity
phase
integer ( pInt ) :: &
i
real ( pReal ) , dimension ( size ( plasticState ( phase ) % state ( : , 1 ) ) ) :: tempState
tempState = 0.0_pReal
do i = 1_pInt , lattice_maxNslipFamily
tempState ( 1 + &
sum ( constitutive_phenopowerlaw_Nslip ( 1 : i - 1 , instance ) ) : &
sum ( constitutive_phenopowerlaw_Nslip ( 1 : i , instance ) ) ) = &
constitutive_phenopowerlaw_tau0_slip ( i , instance )
enddo
do i = 1_pInt , lattice_maxNtwinFamily
tempState ( 1 + sum ( constitutive_phenopowerlaw_Nslip ( : , instance ) ) + &
sum ( constitutive_phenopowerlaw_Ntwin ( 1 : i - 1 , instance ) ) : &
sum ( constitutive_phenopowerlaw_Nslip ( : , instance ) ) + &
sum ( constitutive_phenopowerlaw_Ntwin ( 1 : i , instance ) ) ) = &
constitutive_phenopowerlaw_tau0_twin ( i , instance )
enddo
plasticState ( phase ) % state = spread ( tempState , 2 , size ( plasticState ( phase ) % state ( 1 , : ) ) )
plasticState ( phase ) % state0 = plasticState ( phase ) % state
plasticState ( phase ) % partionedState0 = plasticState ( phase ) % state
end subroutine constitutive_phenopowerlaw_stateInit
!--------------------------------------------------------------------------------------------------
!> @brief sets the relevant state values for a given instance of this plasticity
!--------------------------------------------------------------------------------------------------
subroutine constitutive_phenopowerlaw_aTolState ( phase , instance )
use material , only : &
plasticState
implicit none
integer ( pInt ) , intent ( in ) :: &
instance , & !< number specifying the instance of the plasticity
phase
real ( pReal ) , dimension ( size ( plasticState ( phase ) % aTolState ( : ) ) ) :: tempTol
tempTol = 0.0_pReal
tempTol ( 1 : constitutive_phenopowerlaw_totalNslip ( instance ) + &
constitutive_phenopowerlaw_totalNtwin ( instance ) ) = &
constitutive_phenopowerlaw_aTolResistance ( instance )
tempTol ( 1 + constitutive_phenopowerlaw_totalNslip ( instance ) + &
constitutive_phenopowerlaw_totalNtwin ( instance ) ) = &
constitutive_phenopowerlaw_aTolShear ( instance )
tempTol ( 2 + constitutive_phenopowerlaw_totalNslip ( instance ) + &
constitutive_phenopowerlaw_totalNtwin ( instance ) ) = &
constitutive_phenopowerlaw_aTolTwinFrac ( instance )
tempTol ( 3 + constitutive_phenopowerlaw_totalNslip ( instance ) + &
constitutive_phenopowerlaw_totalNtwin ( instance ) : &
2 + 2 * ( constitutive_phenopowerlaw_totalNslip ( instance ) + &
constitutive_phenopowerlaw_totalNtwin ( instance ) ) ) = &
constitutive_phenopowerlaw_aTolShear ( instance )
plasticState ( phase ) % aTolState = tempTol
end subroutine constitutive_phenopowerlaw_aTolState
#else
2012-10-11 20:19:12 +05:30
!--------------------------------------------------------------------------------------------------
2013-07-01 11:40:42 +05:30
!> @brief sets the initial microstructural state for a given instance of this plasticity
2012-10-11 20:19:12 +05:30
!--------------------------------------------------------------------------------------------------
2014-02-28 15:48:40 +05:30
pure function constitutive_phenopowerlaw_stateInit ( instance )
2013-07-01 11:40:42 +05:30
use lattice , only : &
lattice_maxNslipFamily , &
lattice_maxNtwinFamily
2012-03-09 01:55:28 +05:30
2009-07-22 21:37:19 +05:30
implicit none
2013-07-01 11:40:42 +05:30
integer ( pInt ) , intent ( in ) :: &
2014-03-09 02:20:31 +05:30
instance !< number specifying the instance of the plasticity
2014-02-28 15:48:40 +05:30
real ( pReal ) , dimension ( constitutive_phenopowerlaw_sizeDotState ( instance ) ) :: &
2013-07-01 11:40:42 +05:30
constitutive_phenopowerlaw_stateInit
integer ( pInt ) :: &
i
2009-07-22 21:37:19 +05:30
constitutive_phenopowerlaw_stateInit = 0.0_pReal
2012-02-21 21:30:00 +05:30
do i = 1_pInt , lattice_maxNslipFamily
2009-07-22 21:37:19 +05:30
constitutive_phenopowerlaw_stateInit ( 1 + &
2014-02-28 15:48:40 +05:30
sum ( constitutive_phenopowerlaw_Nslip ( 1 : i - 1 , instance ) ) : &
sum ( constitutive_phenopowerlaw_Nslip ( 1 : i , instance ) ) ) = &
constitutive_phenopowerlaw_tau0_slip ( i , instance )
2009-07-22 21:37:19 +05:30
enddo
2012-02-21 21:30:00 +05:30
do i = 1_pInt , lattice_maxNtwinFamily
2014-02-28 15:48:40 +05:30
constitutive_phenopowerlaw_stateInit ( 1 + sum ( constitutive_phenopowerlaw_Nslip ( : , instance ) ) + &
sum ( constitutive_phenopowerlaw_Ntwin ( 1 : i - 1 , instance ) ) : &
sum ( constitutive_phenopowerlaw_Nslip ( : , instance ) ) + &
sum ( constitutive_phenopowerlaw_Ntwin ( 1 : i , instance ) ) ) = &
constitutive_phenopowerlaw_tau0_twin ( i , instance )
2009-07-22 21:37:19 +05:30
enddo
2012-03-09 01:55:28 +05:30
end function constitutive_phenopowerlaw_stateInit
2009-07-22 21:37:19 +05:30
2012-10-11 20:19:12 +05:30
!--------------------------------------------------------------------------------------------------
2013-07-01 11:40:42 +05:30
!> @brief sets the relevant state values for a given instance of this plasticity
2012-10-11 20:19:12 +05:30
!--------------------------------------------------------------------------------------------------
2014-02-28 15:48:40 +05:30
pure function constitutive_phenopowerlaw_aTolState ( instance )
2013-07-01 11:40:42 +05:30
2012-10-11 20:19:12 +05:30
implicit none
2014-03-09 02:20:31 +05:30
integer ( pInt ) , intent ( in ) :: instance !< number specifying the instance of the plasticity
2013-07-01 11:40:42 +05:30
2014-03-09 02:20:31 +05:30
real ( pReal ) , dimension ( constitutive_phenopowerlaw_sizeState ( instance ) ) :: &
2013-07-01 11:40:42 +05:30
constitutive_phenopowerlaw_aTolState
2012-10-22 20:25:07 +05:30
2014-02-28 15:48:40 +05:30
constitutive_phenopowerlaw_aTolState ( 1 : constitutive_phenopowerlaw_totalNslip ( instance ) + &
constitutive_phenopowerlaw_totalNtwin ( instance ) ) = &
constitutive_phenopowerlaw_aTolResistance ( instance )
constitutive_phenopowerlaw_aTolState ( 1 + constitutive_phenopowerlaw_totalNslip ( instance ) + &
constitutive_phenopowerlaw_totalNtwin ( instance ) ) = &
constitutive_phenopowerlaw_aTolShear ( instance )
constitutive_phenopowerlaw_aTolState ( 2 + constitutive_phenopowerlaw_totalNslip ( instance ) + &
constitutive_phenopowerlaw_totalNtwin ( instance ) ) = &
constitutive_phenopowerlaw_aTolTwinFrac ( instance )
constitutive_phenopowerlaw_aTolState ( 3 + constitutive_phenopowerlaw_totalNslip ( instance ) + &
constitutive_phenopowerlaw_totalNtwin ( instance ) : &
2 + 2 * ( constitutive_phenopowerlaw_totalNslip ( instance ) + &
constitutive_phenopowerlaw_totalNtwin ( instance ) ) ) = &
constitutive_phenopowerlaw_aTolShear ( instance )
2009-09-18 21:07:14 +05:30
2012-03-09 01:55:28 +05:30
end function constitutive_phenopowerlaw_aTolState
2009-09-18 21:07:14 +05:30
2014-05-22 20:46:05 +05:30
#endif
2012-10-11 20:19:12 +05:30
!--------------------------------------------------------------------------------------------------
2013-07-01 11:40:42 +05:30
!> @brief calculates plastic velocity gradient and its tangent
2012-10-11 20:19:12 +05:30
!--------------------------------------------------------------------------------------------------
2013-10-16 18:34:59 +05:30
pure subroutine constitutive_phenopowerlaw_LpAndItsTangent ( Lp , dLp_dTstar99 , Tstar_v , state , ipc , ip , el )
2013-07-01 11:40:42 +05:30
use prec , only : &
p_vec
use math , only : &
math_Plain3333to99 , &
2014-01-22 21:17:49 +05:30
math_Mandel6to33
2013-07-01 11:40:42 +05:30
use lattice , only : &
lattice_Sslip , &
lattice_Sslip_v , &
lattice_Stwin , &
lattice_Stwin_v , &
lattice_maxNslipFamily , &
lattice_maxNtwinFamily , &
lattice_NslipSystem , &
lattice_NtwinSystem , &
2013-08-05 14:53:21 +05:30
lattice_NnonSchmid
2013-07-01 11:40:42 +05:30
use mesh , only : &
mesh_NcpElems , &
mesh_maxNips
use material , only : &
homogenization_maxNgrains , &
material_phase , &
phase_plasticityInstance
2009-07-22 21:37:19 +05:30
implicit none
2014-03-13 12:13:49 +05:30
real ( pReal ) , dimension ( 3 , 3 ) , intent ( out ) :: &
2013-07-01 11:40:42 +05:30
Lp !< plastic velocity gradient
2014-03-13 12:13:49 +05:30
real ( pReal ) , dimension ( 9 , 9 ) , intent ( out ) :: &
2013-10-08 21:57:26 +05:30
dLp_dTstar99 !< derivative of Lp with respect to 2nd Piola Kirchhoff stress
2013-07-01 11:40:42 +05:30
2014-03-13 12:13:49 +05:30
real ( pReal ) , dimension ( 6 ) , intent ( in ) :: &
2013-07-01 11:40:42 +05:30
Tstar_v !< 2nd Piola Kirchhoff stress tensor in Mandel notation
2014-03-13 12:13:49 +05:30
integer ( pInt ) , intent ( in ) :: &
2013-07-01 11:40:42 +05:30
ipc , & !< component-ID of integration point
ip , & !< integration point
el !< element
2014-05-22 20:46:05 +05:30
#ifdef NEWSTATE
real ( pReal ) , dimension ( : ) , intent ( in ) :: &
state
#else
type ( p_vec ) , intent ( in ) :: &
2013-07-01 11:40:42 +05:30
state !< microstructure state
2014-05-22 20:46:05 +05:30
#endif
2013-07-01 11:40:42 +05:30
integer ( pInt ) :: &
2014-02-28 15:48:40 +05:30
instance , &
2013-07-01 11:40:42 +05:30
nSlip , &
2014-03-09 02:20:31 +05:30
nTwin , phase , index_Gamma , index_F , index_myFamily , &
2013-07-01 11:40:42 +05:30
f , i , j , k , l , m , n
real ( pReal ) , dimension ( 3 , 3 , 3 , 3 ) :: &
dLp_dTstar3333 !< derivative of Lp with respect to Tstar as 4th order tensor
real ( pReal ) , dimension ( 3 , 3 , 2 ) :: &
nonSchmid_tensor
2012-03-12 19:39:37 +05:30
real ( pReal ) , dimension ( constitutive_phenopowerlaw_totalNslip ( phase_plasticityInstance ( material_phase ( ipc , ip , el ) ) ) ) :: &
2013-01-22 04:41:16 +05:30
gdot_slip_pos , gdot_slip_neg , dgdot_dtauslip_pos , dgdot_dtauslip_neg , tau_slip_pos , tau_slip_neg
2012-03-12 19:39:37 +05:30
real ( pReal ) , dimension ( constitutive_phenopowerlaw_totalNtwin ( phase_plasticityInstance ( material_phase ( ipc , ip , el ) ) ) ) :: &
2009-07-22 21:37:19 +05:30
gdot_twin , dgdot_dtautwin , tau_twin
2014-03-09 02:20:31 +05:30
phase = material_phase ( ipc , ip , el )
instance = phase_plasticityInstance ( phase )
2014-02-28 15:48:40 +05:30
nSlip = constitutive_phenopowerlaw_totalNslip ( instance )
nTwin = constitutive_phenopowerlaw_totalNtwin ( instance )
2012-02-21 21:30:00 +05:30
index_Gamma = nSlip + nTwin + 1_pInt
index_F = nSlip + nTwin + 2_pInt
2009-07-22 21:37:19 +05:30
Lp = 0.0_pReal
dLp_dTstar3333 = 0.0_pReal
2013-07-12 12:27:15 +05:30
dLp_dTstar99 = 0.0_pReal
2009-10-21 18:40:12 +05:30
2009-07-22 21:37:19 +05:30
j = 0_pInt
2013-07-01 11:40:42 +05:30
slipFamiliesLoop : do f = 1_pInt , lattice_maxNslipFamily
2014-03-09 02:20:31 +05:30
index_myFamily = sum ( lattice_NslipSystem ( 1 : f - 1_pInt , phase ) ) ! at which index starts my family
2014-02-28 15:48:40 +05:30
do i = 1_pInt , constitutive_phenopowerlaw_Nslip ( f , instance ) ! process each (active) slip system in family
2009-07-22 21:37:19 +05:30
j = j + 1_pInt
2012-10-11 20:19:12 +05:30
!--------------------------------------------------------------------------------------------------
! Calculation of Lp
2014-03-09 02:20:31 +05:30
tau_slip_pos ( j ) = dot_product ( Tstar_v , lattice_Sslip_v ( 1 : 6 , 1 , index_myFamily + i , phase ) )
2014-01-22 21:17:49 +05:30
tau_slip_neg ( j ) = tau_slip_pos ( j )
2014-03-09 02:20:31 +05:30
nonSchmid_tensor ( 1 : 3 , 1 : 3 , 1 ) = lattice_Sslip ( 1 : 3 , 1 : 3 , 1 , index_myFamily + i , phase )
2013-08-05 14:53:21 +05:30
nonSchmid_tensor ( 1 : 3 , 1 : 3 , 2 ) = nonSchmid_tensor ( 1 : 3 , 1 : 3 , 1 )
2014-03-09 02:20:31 +05:30
do k = 1 , lattice_NnonSchmid ( phase )
2014-02-28 15:48:40 +05:30
tau_slip_pos ( j ) = tau_slip_pos ( j ) + constitutive_phenopowerlaw_nonSchmidCoeff ( k , instance ) * &
2014-03-09 02:20:31 +05:30
dot_product ( Tstar_v , lattice_Sslip_v ( 1 : 6 , 2 * k , index_myFamily + i , phase ) )
2014-02-28 15:48:40 +05:30
tau_slip_neg ( j ) = tau_slip_neg ( j ) + constitutive_phenopowerlaw_nonSchmidCoeff ( k , instance ) * &
2014-03-09 02:20:31 +05:30
dot_product ( Tstar_v , lattice_Sslip_v ( 1 : 6 , 2 * k + 1 , index_myFamily + i , phase ) )
2014-02-28 15:48:40 +05:30
nonSchmid_tensor ( 1 : 3 , 1 : 3 , 1 ) = nonSchmid_tensor ( 1 : 3 , 1 : 3 , 1 ) + constitutive_phenopowerlaw_nonSchmidCoeff ( k , instance ) * &
2014-03-09 02:20:31 +05:30
lattice_Sslip ( 1 : 3 , 1 : 3 , 2 * k , index_myFamily + i , phase )
2014-02-28 15:48:40 +05:30
nonSchmid_tensor ( 1 : 3 , 1 : 3 , 2 ) = nonSchmid_tensor ( 1 : 3 , 1 : 3 , 2 ) + constitutive_phenopowerlaw_nonSchmidCoeff ( k , instance ) * &
2014-03-09 02:20:31 +05:30
lattice_Sslip ( 1 : 3 , 1 : 3 , 2 * k + 1 , index_myFamily + i , phase )
2013-01-22 04:41:16 +05:30
enddo
2014-05-22 20:46:05 +05:30
#ifdef NEWSTATE
gdot_slip_pos ( j ) = 0.5_pReal * constitutive_phenopowerlaw_gdot0_slip ( instance ) * &
( ( abs ( tau_slip_pos ( j ) ) / state ( j ) ) ** constitutive_phenopowerlaw_n_slip ( instance ) ) * &
sign ( 1.0_pReal , tau_slip_pos ( j ) )
gdot_slip_neg ( j ) = 0.5_pReal * constitutive_phenopowerlaw_gdot0_slip ( instance ) * &
( ( abs ( tau_slip_neg ( j ) ) / state ( j ) ) ** constitutive_phenopowerlaw_n_slip ( instance ) ) * &
sign ( 1.0_pReal , tau_slip_neg ( j ) )
Lp = Lp + ( 1.0_pReal - state ( index_F ) ) * & ! 1-F
( gdot_slip_pos ( j ) + gdot_slip_neg ( j ) ) * lattice_Sslip ( 1 : 3 , 1 : 3 , 1 , index_myFamily + i , phase )
#else
2014-02-28 15:48:40 +05:30
gdot_slip_pos ( j ) = 0.5_pReal * constitutive_phenopowerlaw_gdot0_slip ( instance ) * &
2014-03-13 12:13:49 +05:30
( ( abs ( tau_slip_pos ( j ) ) / state % p ( j ) ) ** constitutive_phenopowerlaw_n_slip ( instance ) ) * &
2014-01-22 21:17:49 +05:30
sign ( 1.0_pReal , tau_slip_pos ( j ) )
2014-05-22 20:46:05 +05:30
2014-02-28 15:48:40 +05:30
gdot_slip_neg ( j ) = 0.5_pReal * constitutive_phenopowerlaw_gdot0_slip ( instance ) * &
2014-03-13 12:13:49 +05:30
( ( abs ( tau_slip_neg ( j ) ) / state % p ( j ) ) ** constitutive_phenopowerlaw_n_slip ( instance ) ) * &
2014-01-22 21:17:49 +05:30
sign ( 1.0_pReal , tau_slip_neg ( j ) )
2014-05-22 20:46:05 +05:30
Lp = Lp + ( 1.0_pReal - state % p ( index_F ) ) * & ! 1-F
2014-05-12 14:58:32 +05:30
( gdot_slip_pos ( j ) + gdot_slip_neg ( j ) ) * lattice_Sslip ( 1 : 3 , 1 : 3 , 1 , index_myFamily + i , phase )
2014-05-22 20:46:05 +05:30
#endif
2009-07-22 21:37:19 +05:30
2012-10-11 20:19:12 +05:30
!--------------------------------------------------------------------------------------------------
! Calculation of the tangent of Lp
2013-01-22 04:41:16 +05:30
if ( gdot_slip_pos ( j ) / = 0.0_pReal ) then
2014-02-28 15:48:40 +05:30
dgdot_dtauslip_pos ( j ) = gdot_slip_pos ( j ) * constitutive_phenopowerlaw_n_slip ( instance ) / tau_slip_pos ( j )
2013-01-22 04:41:16 +05:30
forall ( k = 1_pInt : 3_pInt , l = 1_pInt : 3_pInt , m = 1_pInt : 3_pInt , n = 1_pInt : 3_pInt ) &
dLp_dTstar3333 ( k , l , m , n ) = dLp_dTstar3333 ( k , l , m , n ) + &
2014-03-09 02:20:31 +05:30
dgdot_dtauslip_pos ( j ) * lattice_Sslip ( k , l , 1 , index_myFamily + i , phase ) * &
2013-01-22 04:41:16 +05:30
nonSchmid_tensor ( m , n , 1 )
endif
if ( gdot_slip_neg ( j ) / = 0.0_pReal ) then
2014-02-28 15:48:40 +05:30
dgdot_dtauslip_neg ( j ) = gdot_slip_neg ( j ) * constitutive_phenopowerlaw_n_slip ( instance ) / tau_slip_neg ( j )
2012-02-21 21:30:00 +05:30
forall ( k = 1_pInt : 3_pInt , l = 1_pInt : 3_pInt , m = 1_pInt : 3_pInt , n = 1_pInt : 3_pInt ) &
2009-10-21 18:40:12 +05:30
dLp_dTstar3333 ( k , l , m , n ) = dLp_dTstar3333 ( k , l , m , n ) + &
2014-03-09 02:20:31 +05:30
dgdot_dtauslip_neg ( j ) * lattice_Sslip ( k , l , 1 , index_myFamily + i , phase ) * &
2013-01-22 04:41:16 +05:30
nonSchmid_tensor ( m , n , 2 )
2009-10-21 18:40:12 +05:30
endif
2009-07-22 21:37:19 +05:30
enddo
2013-07-01 11:40:42 +05:30
enddo slipFamiliesLoop
2009-07-22 21:37:19 +05:30
j = 0_pInt
2013-07-01 11:40:42 +05:30
twinFamiliesLoop : do f = 1_pInt , lattice_maxNtwinFamily
2014-03-13 12:13:49 +05:30
index_myFamily = sum ( lattice_NtwinSystem ( 1 : f - 1_pInt , phase ) ) ! at which index starts my family
do i = 1_pInt , constitutive_phenopowerlaw_Ntwin ( f , instance ) ! process each (active) twin system in family
2009-07-22 21:37:19 +05:30
j = j + 1_pInt
2012-10-11 20:19:12 +05:30
!--------------------------------------------------------------------------------------------------
! Calculation of Lp
2014-03-09 02:20:31 +05:30
tau_twin ( j ) = dot_product ( Tstar_v , lattice_Stwin_v ( 1 : 6 , index_myFamily + i , phase ) )
2014-05-22 20:46:05 +05:30
#ifdef NEWSTATE
gdot_twin ( j ) = ( 1.0_pReal - state ( index_F ) ) * & ! 1-F
constitutive_phenopowerlaw_gdot0_twin ( instance ) * &
( abs ( tau_twin ( j ) ) / state ( nSlip + j ) ) ** &
constitutive_phenopowerlaw_n_twin ( instance ) * max ( 0.0_pReal , sign ( 1.0_pReal , tau_twin ( j ) ) )
#else
2014-03-13 12:13:49 +05:30
gdot_twin ( j ) = ( 1.0_pReal - state % p ( index_F ) ) * & ! 1-F
2014-02-28 15:48:40 +05:30
constitutive_phenopowerlaw_gdot0_twin ( instance ) * &
2014-03-13 12:13:49 +05:30
( abs ( tau_twin ( j ) ) / state % p ( nSlip + j ) ) ** &
2014-02-28 15:48:40 +05:30
constitutive_phenopowerlaw_n_twin ( instance ) * max ( 0.0_pReal , sign ( 1.0_pReal , tau_twin ( j ) ) )
2014-05-22 20:46:05 +05:30
#endif
2014-03-09 02:20:31 +05:30
Lp = Lp + gdot_twin ( j ) * lattice_Stwin ( 1 : 3 , 1 : 3 , index_myFamily + i , phase )
2009-07-22 21:37:19 +05:30
2012-10-11 20:19:12 +05:30
!--------------------------------------------------------------------------------------------------
! Calculation of the tangent of Lp
2009-10-21 18:40:12 +05:30
if ( gdot_twin ( j ) / = 0.0_pReal ) then
2014-02-28 15:48:40 +05:30
dgdot_dtautwin ( j ) = gdot_twin ( j ) * constitutive_phenopowerlaw_n_twin ( instance ) / tau_twin ( j )
2012-02-21 21:30:00 +05:30
forall ( k = 1_pInt : 3_pInt , l = 1_pInt : 3_pInt , m = 1_pInt : 3_pInt , n = 1_pInt : 3_pInt ) &
2009-10-21 18:40:12 +05:30
dLp_dTstar3333 ( k , l , m , n ) = dLp_dTstar3333 ( k , l , m , n ) + &
2014-03-09 02:20:31 +05:30
dgdot_dtautwin ( j ) * lattice_Stwin ( k , l , index_myFamily + i , phase ) * &
lattice_Stwin ( m , n , index_myFamily + i , phase )
2009-10-21 18:40:12 +05:30
endif
2009-07-22 21:37:19 +05:30
enddo
2013-07-01 11:40:42 +05:30
enddo twinFamiliesLoop
2009-07-22 21:37:19 +05:30
2013-07-12 12:27:15 +05:30
dLp_dTstar99 = math_Plain3333to99 ( dLp_dTstar3333 )
2009-07-22 21:37:19 +05:30
2012-03-09 01:55:28 +05:30
end subroutine constitutive_phenopowerlaw_LpAndItsTangent
2009-07-22 21:37:19 +05:30
2012-10-11 20:19:12 +05:30
!--------------------------------------------------------------------------------------------------
2013-07-01 11:40:42 +05:30
!> @brief calculates the rate of change of microstructure
2012-10-11 20:19:12 +05:30
!--------------------------------------------------------------------------------------------------
2013-10-16 18:34:59 +05:30
function constitutive_phenopowerlaw_dotState ( Tstar_v , state , ipc , ip , el )
2013-07-01 11:40:42 +05:30
use prec , only : &
p_vec
use lattice , only : &
lattice_Sslip_v , &
lattice_Stwin_v , &
lattice_maxNslipFamily , &
lattice_maxNtwinFamily , &
lattice_NslipSystem , &
lattice_NtwinSystem , &
lattice_shearTwin , &
2013-08-05 14:53:21 +05:30
lattice_NnonSchmid
2013-07-01 11:40:42 +05:30
use mesh , only : &
mesh_NcpElems , &
mesh_maxNips
use material , only : &
homogenization_maxNgrains , &
material_phase , &
phase_plasticityInstance
2012-03-09 01:55:28 +05:30
2009-07-22 21:37:19 +05:30
implicit none
2014-03-13 12:13:49 +05:30
real ( pReal ) , dimension ( 6 ) , intent ( in ) :: &
2013-07-01 11:40:42 +05:30
Tstar_v !< 2nd Piola Kirchhoff stress tensor in Mandel notation
2014-03-13 12:13:49 +05:30
integer ( pInt ) , intent ( in ) :: &
2013-07-01 11:40:42 +05:30
ipc , & !< component-ID of integration point
ip , & !< integration point
2014-05-22 20:46:05 +05:30
el !< element !< microstructure state
#ifdef NEWSTATE
real ( pReal ) , dimension ( : ) , intent ( in ) :: &
state
real ( pReal ) , dimension ( size ( state ) ) :: &
constitutive_phenopowerlaw_dotState
#else
type ( p_vec ) , intent ( in ) :: &
2013-07-01 11:40:42 +05:30
state !< microstructure state
real ( pReal ) , dimension ( constitutive_phenopowerlaw_sizeDotState ( phase_plasticityInstance ( material_phase ( ipc , ip , el ) ) ) ) :: &
constitutive_phenopowerlaw_dotState
2014-05-22 20:46:05 +05:30
#endif
2013-07-01 11:40:42 +05:30
2013-10-08 21:57:26 +05:30
integer ( pInt ) :: &
2014-03-09 02:20:31 +05:30
instance , phase , &
2013-10-08 21:57:26 +05:30
nSlip , nTwin , &
f , i , j , k , &
index_Gamma , index_F , index_myFamily , &
offset_accshear_slip , offset_accshear_twin
real ( pReal ) :: &
c_SlipSlip , c_SlipTwin , c_TwinSlip , c_TwinTwin , &
ssat_offset
2013-07-01 11:40:42 +05:30
2012-03-12 19:39:37 +05:30
real ( pReal ) , dimension ( constitutive_phenopowerlaw_totalNslip ( phase_plasticityInstance ( material_phase ( ipc , ip , el ) ) ) ) :: &
2013-01-22 04:41:16 +05:30
gdot_slip , tau_slip_pos , tau_slip_neg , left_SlipSlip , left_SlipTwin , right_SlipSlip , right_TwinSlip
2012-03-12 19:39:37 +05:30
real ( pReal ) , dimension ( constitutive_phenopowerlaw_totalNtwin ( phase_plasticityInstance ( material_phase ( ipc , ip , el ) ) ) ) :: &
2012-11-14 15:52:34 +05:30
gdot_twin , tau_twin , left_TwinSlip , left_TwinTwin , right_SlipTwin , right_TwinTwin
2014-03-09 02:20:31 +05:30
phase = material_phase ( ipc , ip , el )
instance = phase_plasticityInstance ( phase )
2009-07-22 21:37:19 +05:30
2014-02-28 15:48:40 +05:30
nSlip = constitutive_phenopowerlaw_totalNslip ( instance )
nTwin = constitutive_phenopowerlaw_totalNtwin ( instance )
2009-07-22 21:37:19 +05:30
2012-02-21 21:30:00 +05:30
index_Gamma = nSlip + nTwin + 1_pInt
index_F = nSlip + nTwin + 2_pInt
2013-02-08 19:03:25 +05:30
offset_accshear_slip = nSlip + nTwin + 2_pInt
offset_accshear_twin = nSlip + nTwin + 2_pInt + nSlip
2009-07-22 21:37:19 +05:30
constitutive_phenopowerlaw_dotState = 0.0_pReal
2012-10-11 20:19:12 +05:30
!--------------------------------------------------------------------------------------------------
2012-10-22 20:25:07 +05:30
! system-independent (nonlinear) prefactors to M_Xx (X influenced by x) matrices
2014-05-22 20:46:05 +05:30
#ifdef NEWSTATE
c_SlipSlip = constitutive_phenopowerlaw_h0_SlipSlip ( instance ) * &
( 1.0_pReal + constitutive_phenopowerlaw_twinC ( instance ) * state ( index_F ) ** &
constitutive_phenopowerlaw_twinB ( instance ) )
c_SlipTwin = 0.0_pReal
c_TwinSlip = constitutive_phenopowerlaw_h0_TwinSlip ( instance ) * &
state ( index_Gamma ) ** constitutive_phenopowerlaw_twinE ( instance )
c_TwinTwin = constitutive_phenopowerlaw_h0_TwinTwin ( instance ) * &
state ( index_F ) ** constitutive_phenopowerlaw_twinD ( instance )
#else
2014-02-28 15:48:40 +05:30
c_SlipSlip = constitutive_phenopowerlaw_h0_SlipSlip ( instance ) * &
2014-03-13 12:13:49 +05:30
( 1.0_pReal + constitutive_phenopowerlaw_twinC ( instance ) * state % p ( index_F ) ** &
2014-02-28 15:48:40 +05:30
constitutive_phenopowerlaw_twinB ( instance ) )
2012-11-14 15:52:34 +05:30
c_SlipTwin = 0.0_pReal
2014-02-28 15:48:40 +05:30
c_TwinSlip = constitutive_phenopowerlaw_h0_TwinSlip ( instance ) * &
2014-03-13 12:13:49 +05:30
state % p ( index_Gamma ) ** constitutive_phenopowerlaw_twinE ( instance )
2014-02-28 15:48:40 +05:30
c_TwinTwin = constitutive_phenopowerlaw_h0_TwinTwin ( instance ) * &
2014-03-13 12:13:49 +05:30
state % p ( index_F ) ** constitutive_phenopowerlaw_twinD ( instance )
2014-05-22 20:46:05 +05:30
#endif
2013-07-01 11:40:42 +05:30
!--------------------------------------------------------------------------------------------------
! calculate left and right vectors and calculate dot gammas
2014-05-22 20:46:05 +05:30
#ifdef NEWSTATE
ssat_offset = constitutive_phenopowerlaw_spr ( instance ) * sqrt ( state ( index_F ) )
#else
ssat_offset = constitutive_phenopowerlaw_spr ( instance ) * sqrt ( state % p ( index_F ) ) !< microstructure state
#endif
2009-07-22 21:37:19 +05:30
j = 0_pInt
2013-07-01 11:40:42 +05:30
slipFamiliesLoop1 : do f = 1_pInt , lattice_maxNslipFamily
2014-03-09 02:20:31 +05:30
index_myFamily = sum ( lattice_NslipSystem ( 1 : f - 1_pInt , phase ) ) ! at which index starts my family
do i = 1_pInt , constitutive_phenopowerlaw_Nslip ( f , instance ) ! process each (active) slip system in family
2009-07-22 21:37:19 +05:30
j = j + 1_pInt
2013-07-01 11:40:42 +05:30
left_SlipSlip ( j ) = 1.0_pReal ! no system-dependent left part
left_SlipTwin ( j ) = 1.0_pReal ! no system-dependent left part
2014-05-22 20:46:05 +05:30
#ifdef NEWSTATE
right_SlipSlip ( j ) = abs ( 1.0_pReal - state ( j ) / &
( constitutive_phenopowerlaw_tausat_slip ( f , instance ) + ssat_offset ) ) &
** constitutive_phenopowerlaw_a_slip ( instance ) &
* sign ( 1.0_pReal , 1.0_pReal - state ( j ) / &
( constitutive_phenopowerlaw_tausat_slip ( f , instance ) + ssat_offset ) )
#else
2014-03-13 12:13:49 +05:30
right_SlipSlip ( j ) = abs ( 1.0_pReal - state % p ( j ) / &
2014-02-28 15:48:40 +05:30
( constitutive_phenopowerlaw_tausat_slip ( f , instance ) + ssat_offset ) ) &
** constitutive_phenopowerlaw_a_slip ( instance ) &
2014-03-13 12:13:49 +05:30
* sign ( 1.0_pReal , 1.0_pReal - state % p ( j ) / &
2014-05-22 20:46:05 +05:30
( constitutive_phenopowerlaw_tausat_slip ( f , instance ) + ssat_offset ) ) !< microstructure state
#endif
2013-07-01 11:40:42 +05:30
right_TwinSlip ( j ) = 1.0_pReal ! no system-dependent part
2009-07-22 21:37:19 +05:30
2012-10-11 20:19:12 +05:30
!--------------------------------------------------------------------------------------------------
! Calculation of dot gamma
2014-03-09 02:20:31 +05:30
tau_slip_pos ( j ) = dot_product ( Tstar_v , lattice_Sslip_v ( 1 : 6 , 1 , index_myFamily + i , phase ) )
2013-01-22 04:41:16 +05:30
tau_slip_neg ( j ) = tau_slip_pos ( j )
2014-03-09 02:20:31 +05:30
do k = 1 , lattice_NnonSchmid ( phase )
2014-02-28 15:48:40 +05:30
tau_slip_pos ( j ) = tau_slip_pos ( j ) + constitutive_phenopowerlaw_nonSchmidCoeff ( k , instance ) * &
2014-03-09 02:20:31 +05:30
dot_product ( Tstar_v , lattice_Sslip_v ( 1 : 6 , 2 * k , index_myFamily + i , phase ) )
2014-02-28 15:48:40 +05:30
tau_slip_neg ( j ) = tau_slip_neg ( j ) + constitutive_phenopowerlaw_nonSchmidCoeff ( k , instance ) * &
2014-03-09 02:20:31 +05:30
dot_product ( Tstar_v , lattice_Sslip_v ( 1 : 6 , 2 * k + 1 , index_myFamily + i , phase ) )
2013-01-22 04:41:16 +05:30
enddo
2014-05-22 20:46:05 +05:30
#ifdef NEWSTATE
gdot_slip ( j ) = constitutive_phenopowerlaw_gdot0_slip ( instance ) * 0.5_pReal * &
( ( abs ( tau_slip_pos ( j ) ) / state ( j ) ) ** constitutive_phenopowerlaw_n_slip ( instance ) &
+ ( abs ( tau_slip_neg ( j ) ) / state ( j ) ) ** constitutive_phenopowerlaw_n_slip ( instance ) ) &
* sign ( 1.0_pReal , tau_slip_pos ( j ) )
#else
2014-02-28 15:48:40 +05:30
gdot_slip ( j ) = constitutive_phenopowerlaw_gdot0_slip ( instance ) * 0.5_pReal * &
2014-03-13 12:13:49 +05:30
( ( abs ( tau_slip_pos ( j ) ) / state % p ( j ) ) ** constitutive_phenopowerlaw_n_slip ( instance ) &
+ ( abs ( tau_slip_neg ( j ) ) / state % p ( j ) ) ** constitutive_phenopowerlaw_n_slip ( instance ) ) &
2014-05-22 20:46:05 +05:30
* sign ( 1.0_pReal , tau_slip_pos ( j ) )
#endif
2013-07-01 11:40:42 +05:30
enddo
enddo slipFamiliesLoop1
2009-07-22 21:37:19 +05:30
2014-05-22 20:46:05 +05:30
2009-07-22 21:37:19 +05:30
j = 0_pInt
2013-07-01 11:40:42 +05:30
twinFamiliesLoop1 : do f = 1_pInt , lattice_maxNtwinFamily
2014-03-09 02:20:31 +05:30
index_myFamily = sum ( lattice_NtwinSystem ( 1 : f - 1_pInt , phase ) ) ! at which index starts my family
do i = 1_pInt , constitutive_phenopowerlaw_Ntwin ( f , instance ) ! process each (active) twin system in family
2009-07-22 21:37:19 +05:30
j = j + 1_pInt
2012-11-14 15:52:34 +05:30
left_TwinSlip ( j ) = 1.0_pReal ! no system-dependent right part
left_TwinTwin ( j ) = 1.0_pReal ! no system-dependent right part
right_SlipTwin ( j ) = 1.0_pReal ! no system-dependent right part
right_TwinTwin ( j ) = 1.0_pReal ! no system-dependent right part
2012-10-22 20:25:07 +05:30
2013-07-01 11:40:42 +05:30
!--------------------------------------------------------------------------------------------------
! Calculation of dot vol frac
2014-03-09 02:20:31 +05:30
tau_twin ( j ) = dot_product ( Tstar_v , lattice_Stwin_v ( 1 : 6 , index_myFamily + i , phase ) )
2014-05-22 20:46:05 +05:30
#ifdef NEWSTATE
gdot_twin ( j ) = ( 1.0_pReal - state ( index_F ) ) * & ! 1-F
constitutive_phenopowerlaw_gdot0_twin ( instance ) * &
( abs ( tau_twin ( j ) ) / state ( nSlip + j ) ) ** &
constitutive_phenopowerlaw_n_twin ( instance ) * max ( 0.0_pReal , sign ( 1.0_pReal , tau_twin ( j ) ) )
#else
2014-03-13 12:13:49 +05:30
gdot_twin ( j ) = ( 1.0_pReal - state % p ( index_F ) ) * & ! 1-F
2014-02-28 15:48:40 +05:30
constitutive_phenopowerlaw_gdot0_twin ( instance ) * &
2014-03-13 12:13:49 +05:30
( abs ( tau_twin ( j ) ) / state % p ( nSlip + j ) ) ** &
2014-02-28 15:48:40 +05:30
constitutive_phenopowerlaw_n_twin ( instance ) * max ( 0.0_pReal , sign ( 1.0_pReal , tau_twin ( j ) ) )
2014-05-22 20:46:05 +05:30
#endif
2009-07-22 21:37:19 +05:30
enddo
2013-07-01 11:40:42 +05:30
enddo twinFamiliesLoop1
2009-07-22 21:37:19 +05:30
2012-10-11 20:19:12 +05:30
!--------------------------------------------------------------------------------------------------
! calculate the overall hardening based on above
2009-07-22 21:37:19 +05:30
j = 0_pInt
2013-07-01 11:40:42 +05:30
slipFamiliesLoop2 : do f = 1_pInt , lattice_maxNslipFamily
2014-03-09 02:20:31 +05:30
do i = 1_pInt , constitutive_phenopowerlaw_Nslip ( f , instance ) ! process each (active) slip system in family
2009-07-22 21:37:19 +05:30
j = j + 1_pInt
2012-10-22 20:25:07 +05:30
constitutive_phenopowerlaw_dotState ( j ) = & ! evolution of slip resistance j
2012-11-14 15:52:34 +05:30
c_SlipSlip * left_SlipSlip ( j ) * &
2014-02-28 15:48:40 +05:30
dot_product ( constitutive_phenopowerlaw_hardeningMatrix_SlipSlip ( j , 1 : nSlip , instance ) , &
2012-11-14 15:52:34 +05:30
right_SlipSlip * abs ( gdot_slip ) ) + & ! dot gamma_slip modulated by right-side slip factor
c_SlipTwin * left_SlipTwin ( j ) * &
2014-02-28 15:48:40 +05:30
dot_product ( constitutive_phenopowerlaw_hardeningMatrix_SlipTwin ( j , 1 : nTwin , instance ) , &
2012-11-14 15:52:34 +05:30
right_SlipTwin * gdot_twin ) ! dot gamma_twin modulated by right-side twin factor
2009-07-22 21:37:19 +05:30
constitutive_phenopowerlaw_dotState ( index_Gamma ) = constitutive_phenopowerlaw_dotState ( index_Gamma ) + &
abs ( gdot_slip ( j ) )
2013-02-08 19:03:25 +05:30
constitutive_phenopowerlaw_dotState ( offset_accshear_slip + j ) = abs ( gdot_slip ( j ) )
2009-07-22 21:37:19 +05:30
enddo
2013-07-01 11:40:42 +05:30
enddo slipFamiliesLoop2
2014-05-22 20:46:05 +05:30
2009-07-22 21:37:19 +05:30
j = 0_pInt
2013-07-01 11:40:42 +05:30
twinFamiliesLoop2 : do f = 1_pInt , lattice_maxNtwinFamily
2014-03-09 02:20:31 +05:30
index_myFamily = sum ( lattice_NtwinSystem ( 1 : f - 1_pInt , phase ) ) ! at which index starts my family
do i = 1_pInt , constitutive_phenopowerlaw_Ntwin ( f , instance ) ! process each (active) twin system in family
2009-07-22 21:37:19 +05:30
j = j + 1_pInt
2012-10-22 20:25:07 +05:30
constitutive_phenopowerlaw_dotState ( j + nSlip ) = & ! evolution of twin resistance j
2012-11-14 15:52:34 +05:30
c_TwinSlip * left_TwinSlip ( j ) * &
2014-02-28 15:48:40 +05:30
dot_product ( constitutive_phenopowerlaw_hardeningMatrix_TwinSlip ( j , 1 : nSlip , instance ) , &
2012-11-14 15:52:34 +05:30
right_TwinSlip * abs ( gdot_slip ) ) + & ! dot gamma_slip modulated by right-side slip factor
c_TwinTwin * left_TwinTwin ( j ) * &
2014-02-28 15:48:40 +05:30
dot_product ( constitutive_phenopowerlaw_hardeningMatrix_TwinTwin ( j , 1 : nTwin , instance ) , &
2012-11-14 15:52:34 +05:30
right_TwinTwin * gdot_twin ) ! dot gamma_twin modulated by right-side twin factor
2014-05-22 20:46:05 +05:30
#ifndef NEWSTATE
2014-03-13 12:13:49 +05:30
if ( state % p ( index_F ) < 0.98_pReal ) & ! ensure twin volume fractions stays below 1.0
2014-01-22 21:17:49 +05:30
constitutive_phenopowerlaw_dotState ( index_F ) = constitutive_phenopowerlaw_dotState ( index_F ) + &
2014-03-09 02:20:31 +05:30
gdot_twin ( j ) / lattice_shearTwin ( index_myFamily + i , phase )
2014-05-22 20:46:05 +05:30
#else
if ( state ( index_F ) < 0.98_pReal ) & ! ensure twin volume fractions stays below 1.0
constitutive_phenopowerlaw_dotState ( index_F ) = constitutive_phenopowerlaw_dotState ( index_F ) + &
gdot_twin ( j ) / lattice_shearTwin ( index_myFamily + i , phase )
#endif
2013-02-08 19:03:25 +05:30
constitutive_phenopowerlaw_dotState ( offset_accshear_twin + j ) = abs ( gdot_twin ( j ) )
2009-07-22 21:37:19 +05:30
enddo
2013-07-01 11:40:42 +05:30
enddo twinFamiliesLoop2
2014-05-22 20:46:05 +05:30
2012-03-09 01:55:28 +05:30
end function constitutive_phenopowerlaw_dotState
2009-07-22 21:37:19 +05:30
2012-10-11 20:19:12 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief return array of constitutive results
!--------------------------------------------------------------------------------------------------
2013-10-16 18:34:59 +05:30
pure function constitutive_phenopowerlaw_postResults ( Tstar_v , state , ipc , ip , el )
2013-07-01 11:40:42 +05:30
use prec , only : &
p_vec
use mesh , only : &
mesh_NcpElems , &
mesh_maxNips
use material , only : &
homogenization_maxNgrains , &
material_phase , &
phase_plasticityInstance , &
phase_Noutput
use lattice , only : &
lattice_Sslip_v , &
lattice_Stwin_v , &
lattice_maxNslipFamily , &
lattice_maxNtwinFamily , &
lattice_NslipSystem , &
2013-08-05 14:53:21 +05:30
lattice_NtwinSystem , &
lattice_NnonSchmid
2013-07-01 11:40:42 +05:30
use mesh , only : &
mesh_NcpElems , &
mesh_maxNips
2009-07-22 21:37:19 +05:30
implicit none
2014-03-13 12:13:49 +05:30
real ( pReal ) , dimension ( 6 ) , intent ( in ) :: &
2013-07-01 11:40:42 +05:30
Tstar_v !< 2nd Piola Kirchhoff stress tensor in Mandel notation
2014-03-13 12:13:49 +05:30
integer ( pInt ) , intent ( in ) :: &
2013-07-01 11:40:42 +05:30
ipc , & !< component-ID of integration point
ip , & !< integration point
2014-05-22 20:46:05 +05:30
el !< element !< microstructure state
#ifdef NEWSTATE
real ( pReal ) , dimension ( : ) , intent ( in ) :: &
state
#else
type ( p_vec ) , intent ( in ) :: &
2013-07-01 11:40:42 +05:30
state !< microstructure state
2014-05-22 20:46:05 +05:30
#endif
2013-07-01 11:40:42 +05:30
2012-03-12 19:39:37 +05:30
real ( pReal ) , dimension ( constitutive_phenopowerlaw_sizePostResults ( phase_plasticityInstance ( material_phase ( ipc , ip , el ) ) ) ) :: &
2009-07-22 21:37:19 +05:30
constitutive_phenopowerlaw_postResults
2013-07-01 11:40:42 +05:30
integer ( pInt ) :: &
2014-03-09 02:20:31 +05:30
instance , phase , &
2013-07-01 11:40:42 +05:30
nSlip , nTwin , &
o , f , i , c , j , k , &
index_Gamma , index_F , index_accshear_slip , index_accshear_twin , index_myFamily
real ( pReal ) :: &
tau_slip_pos , tau_slip_neg , tau
2014-03-09 02:20:31 +05:30
phase = material_phase ( ipc , ip , el )
instance = phase_plasticityInstance ( phase )
2009-07-22 21:37:19 +05:30
2014-02-28 15:48:40 +05:30
nSlip = constitutive_phenopowerlaw_totalNslip ( instance )
nTwin = constitutive_phenopowerlaw_totalNtwin ( instance )
2009-07-22 21:37:19 +05:30
2012-02-21 21:30:00 +05:30
index_Gamma = nSlip + nTwin + 1_pInt
index_F = nSlip + nTwin + 2_pInt
2013-02-06 23:39:11 +05:30
index_accshear_slip = nSlip + nTwin + 3_pInt
index_accshear_twin = nSlip + nTwin + 3_pInt + nSlip
2009-07-22 21:37:19 +05:30
constitutive_phenopowerlaw_postResults = 0.0_pReal
c = 0_pInt
2014-06-25 04:23:25 +05:30
outputsLoop : do o = 1_pInt , constitutive_phenopowerlaw_Noutput ( instance )
2014-02-28 15:48:40 +05:30
select case ( constitutive_phenopowerlaw_outputID ( o , instance ) )
2013-11-27 17:09:28 +05:30
case ( resistance_slip_ID )
2014-05-22 20:46:05 +05:30
#ifdef NEWSTATE
constitutive_phenopowerlaw_postResults ( c + 1_pInt : c + nSlip ) = state ( 1 : nSlip )
#else
2014-03-13 12:13:49 +05:30
constitutive_phenopowerlaw_postResults ( c + 1_pInt : c + nSlip ) = state % p ( 1 : nSlip )
2014-05-22 20:46:05 +05:30
#endif
2009-07-22 21:37:19 +05:30
c = c + nSlip
2013-11-27 17:09:28 +05:30
case ( accumulatedshear_slip_ID )
2014-05-22 20:46:05 +05:30
#ifdef NEWSTATE
constitutive_phenopowerlaw_postResults ( c + 1_pInt : c + nSlip ) = state ( index_accshear_slip : &
index_accshear_slip + nSlip - 1_pInt )
#else
2014-03-13 12:13:49 +05:30
constitutive_phenopowerlaw_postResults ( c + 1_pInt : c + nSlip ) = state % p ( index_accshear_slip : &
2014-05-22 20:46:05 +05:30
index_accshear_slip + nSlip - 1_pInt )
#endif
2013-02-06 23:39:11 +05:30
c = c + nSlip
2013-11-27 17:09:28 +05:30
case ( shearrate_slip_ID )
2009-07-22 21:37:19 +05:30
j = 0_pInt
2013-07-01 11:40:42 +05:30
slipFamiliesLoop1 : do f = 1_pInt , lattice_maxNslipFamily
2014-03-09 02:20:31 +05:30
index_myFamily = sum ( lattice_NslipSystem ( 1 : f - 1_pInt , phase ) ) ! at which index starts my family
do i = 1_pInt , constitutive_phenopowerlaw_Nslip ( f , instance ) ! process each (active) slip system in family
2009-07-22 21:37:19 +05:30
j = j + 1_pInt
2014-03-09 02:20:31 +05:30
tau_slip_pos = dot_product ( Tstar_v , lattice_Sslip_v ( 1 : 6 , 1 , index_myFamily + i , phase ) )
2013-01-22 04:41:16 +05:30
tau_slip_neg = tau_slip_pos
2014-03-09 02:20:31 +05:30
do k = 1 , lattice_NnonSchmid ( phase )
2014-02-28 15:48:40 +05:30
tau_slip_pos = tau_slip_pos + constitutive_phenopowerlaw_nonSchmidCoeff ( k , instance ) * &
2014-03-09 02:20:31 +05:30
dot_product ( Tstar_v , lattice_Sslip_v ( 1 : 6 , 2 * k , index_myFamily + i , phase ) )
2014-02-28 15:48:40 +05:30
tau_slip_neg = tau_slip_neg + constitutive_phenopowerlaw_nonSchmidCoeff ( k , instance ) * &
2014-03-09 02:20:31 +05:30
dot_product ( Tstar_v , lattice_Sslip_v ( 1 : 6 , 2 * k + 1 , index_myFamily + i , phase ) )
2013-01-22 04:41:16 +05:30
enddo
2014-05-22 20:46:05 +05:30
#ifdef NEWSTATE
constitutive_phenopowerlaw_postResults ( c + j ) = constitutive_phenopowerlaw_gdot0_slip ( instance ) * 0.5_pReal * &
( ( abs ( tau_slip_pos ) / state ( j ) ) ** constitutive_phenopowerlaw_n_slip ( instance ) &
+ ( abs ( tau_slip_neg ) / state ( j ) ) ** constitutive_phenopowerlaw_n_slip ( instance ) ) &
* sign ( 1.0_pReal , tau_slip_pos )
#else
2014-02-28 15:48:40 +05:30
constitutive_phenopowerlaw_postResults ( c + j ) = constitutive_phenopowerlaw_gdot0_slip ( instance ) * 0.5_pReal * &
2014-03-13 12:13:49 +05:30
( ( abs ( tau_slip_pos ) / state % p ( j ) ) ** constitutive_phenopowerlaw_n_slip ( instance ) &
+ ( abs ( tau_slip_neg ) / state % p ( j ) ) ** constitutive_phenopowerlaw_n_slip ( instance ) ) &
2013-01-22 04:41:16 +05:30
* sign ( 1.0_pReal , tau_slip_pos )
2014-05-22 20:46:05 +05:30
#endif
2013-07-01 11:40:42 +05:30
enddo
enddo slipFamiliesLoop1
2009-07-22 21:37:19 +05:30
c = c + nSlip
2013-11-27 17:09:28 +05:30
case ( resolvedstress_slip_ID )
2009-07-22 21:37:19 +05:30
j = 0_pInt
2013-07-01 11:40:42 +05:30
slipFamiliesLoop2 : do f = 1_pInt , lattice_maxNslipFamily
2014-03-09 02:20:31 +05:30
index_myFamily = sum ( lattice_NslipSystem ( 1 : f - 1_pInt , phase ) ) ! at which index starts my family
do i = 1_pInt , constitutive_phenopowerlaw_Nslip ( f , instance ) ! process each (active) slip system in family
2009-07-22 21:37:19 +05:30
j = j + 1_pInt
2013-07-01 11:40:42 +05:30
constitutive_phenopowerlaw_postResults ( c + j ) = &
2014-03-09 02:20:31 +05:30
dot_product ( Tstar_v , lattice_Sslip_v ( 1 : 6 , 1 , index_myFamily + i , phase ) )
2013-07-01 11:40:42 +05:30
enddo
enddo slipFamiliesLoop2
2009-07-22 21:37:19 +05:30
c = c + nSlip
2013-11-27 17:09:28 +05:30
case ( totalshear_ID )
2014-05-22 20:46:05 +05:30
#ifdef NEWSTATE
constitutive_phenopowerlaw_postResults ( c + 1_pInt ) = &
state ( index_Gamma )
#else
2013-07-01 11:40:42 +05:30
constitutive_phenopowerlaw_postResults ( c + 1_pInt ) = &
2014-03-13 12:13:49 +05:30
state % p ( index_Gamma )
2014-05-22 20:46:05 +05:30
#endif
2012-02-21 21:30:00 +05:30
c = c + 1_pInt
2009-07-22 21:37:19 +05:30
2013-11-27 17:09:28 +05:30
case ( resistance_twin_ID )
2014-05-22 20:46:05 +05:30
#ifdef NEWSTATE
constitutive_phenopowerlaw_postResults ( c + 1_pInt : c + nTwin ) = &
state ( 1_pInt + nSlip : nTwin + nSlip - 1_pInt )
#else
2013-07-01 11:40:42 +05:30
constitutive_phenopowerlaw_postResults ( c + 1_pInt : c + nTwin ) = &
2014-05-22 20:46:05 +05:30
state % p ( 1_pInt + nSlip : nTwin + nSlip - 1_pInt )
#endif
2009-07-22 21:37:19 +05:30
c = c + nTwin
2009-10-21 18:40:12 +05:30
2013-11-27 17:09:28 +05:30
case ( accumulatedshear_twin_ID )
2014-05-22 20:46:05 +05:30
#ifdef NEWSTATE
2013-07-01 11:40:42 +05:30
constitutive_phenopowerlaw_postResults ( c + 1_pInt : c + nTwin ) = &
2014-05-22 20:46:05 +05:30
state ( index_accshear_twin : index_accshear_twin + nTwin - 1_pInt )
#else
constitutive_phenopowerlaw_postResults ( c + 1_pInt : c + nTwin ) = &
state % p ( index_accshear_twin : index_accshear_twin + nTwin - 1_pInt )
#endif
2013-02-06 23:39:11 +05:30
c = c + nTwin
2013-11-27 17:09:28 +05:30
case ( shearrate_twin_ID )
2009-07-22 21:37:19 +05:30
j = 0_pInt
2013-07-01 11:40:42 +05:30
twinFamiliesLoop1 : do f = 1_pInt , lattice_maxNtwinFamily
2014-03-09 02:20:31 +05:30
index_myFamily = sum ( lattice_NtwinSystem ( 1 : f - 1_pInt , phase ) ) ! at which index starts my family
do i = 1_pInt , constitutive_phenopowerlaw_Ntwin ( f , instance ) ! process each (active) twin system in family
2009-07-22 21:37:19 +05:30
j = j + 1_pInt
2014-03-09 02:20:31 +05:30
tau = dot_product ( Tstar_v , lattice_Stwin_v ( 1 : 6 , index_myFamily + i , phase ) )
2014-05-22 20:46:05 +05:30
#ifdef NEWSTATE
constitutive_phenopowerlaw_postResults ( c + j ) = ( 1.0_pReal - state ( index_F ) ) * & ! 1-F
constitutive_phenopowerlaw_gdot0_twin ( instance ) * &
( abs ( tau ) / state ( j + nSlip ) ) ** &
constitutive_phenopowerlaw_n_twin ( instance ) * max ( 0.0_pReal , sign ( 1.0_pReal , tau ) )
#else
2014-03-13 12:13:49 +05:30
constitutive_phenopowerlaw_postResults ( c + j ) = ( 1.0_pReal - state % p ( index_F ) ) * & ! 1-F
2014-02-28 15:48:40 +05:30
constitutive_phenopowerlaw_gdot0_twin ( instance ) * &
2014-03-13 12:13:49 +05:30
( abs ( tau ) / state % p ( j + nSlip ) ) ** &
2014-02-28 15:48:40 +05:30
constitutive_phenopowerlaw_n_twin ( instance ) * max ( 0.0_pReal , sign ( 1.0_pReal , tau ) )
2014-05-22 20:46:05 +05:30
#endif
2013-07-01 11:40:42 +05:30
enddo
enddo twinFamiliesLoop1
2009-07-22 21:37:19 +05:30
c = c + nTwin
2013-11-27 17:09:28 +05:30
case ( resolvedstress_twin_ID )
2009-07-22 21:37:19 +05:30
j = 0_pInt
2013-07-01 11:40:42 +05:30
twinFamiliesLoop2 : do f = 1_pInt , lattice_maxNtwinFamily
2014-03-09 02:20:31 +05:30
index_myFamily = sum ( lattice_NtwinSystem ( 1 : f - 1_pInt , phase ) ) ! at which index starts my family
do i = 1_pInt , constitutive_phenopowerlaw_Ntwin ( f , instance ) ! process each (active) twin system in family
2009-07-22 21:37:19 +05:30
j = j + 1_pInt
2013-07-01 11:40:42 +05:30
constitutive_phenopowerlaw_postResults ( c + j ) = &
2014-03-09 02:20:31 +05:30
dot_product ( Tstar_v , lattice_Stwin_v ( 1 : 6 , index_myFamily + i , phase ) )
2013-07-01 11:40:42 +05:30
enddo
enddo twinFamiliesLoop2
2009-07-22 21:37:19 +05:30
c = c + nTwin
2013-11-27 17:09:28 +05:30
case ( totalvolfrac_ID )
2014-05-22 20:46:05 +05:30
#ifdef NEWSTATE
constitutive_phenopowerlaw_postResults ( c + 1_pInt ) = state ( index_F )
#else
2014-03-13 12:13:49 +05:30
constitutive_phenopowerlaw_postResults ( c + 1_pInt ) = state % p ( index_F )
2014-05-22 20:46:05 +05:30
#endif
2012-02-21 21:30:00 +05:30
c = c + 1_pInt
2009-07-22 21:37:19 +05:30
end select
2013-07-01 11:40:42 +05:30
enddo outputsLoop
2009-07-22 21:37:19 +05:30
2012-03-09 01:55:28 +05:30
end function constitutive_phenopowerlaw_postResults
2009-07-22 21:37:19 +05:30
2014-05-22 20:46:05 +05:30
2012-03-09 01:55:28 +05:30
end module constitutive_phenopowerlaw