2011-04-07 12:50:28 +05:30
! Copyright 2011 Max-Planck-Institut für Eisenforschung GmbH
2011-04-04 19:39:54 +05:30
!
! This file is part of DAMASK,
2011-04-07 12:50:28 +05:30
! the Düsseldorf Advanced MAterial Simulation Kit.
2011-04-04 19:39:54 +05:30
!
! DAMASK is free software: you can redistribute it and/or modify
! it under the terms of the GNU General Public License as published by
! the Free Software Foundation, either version 3 of the License, or
! (at your option) any later version.
!
! DAMASK is distributed in the hope that it will be useful,
! but WITHOUT ANY WARRANTY; without even the implied warranty of
! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
! GNU General Public License for more details.
!
! You should have received a copy of the GNU General Public License
! along with DAMASK. If not, see <http://www.gnu.org/licenses/>.
!
!##############################################################
2010-09-23 13:35:50 +05:30
!* $Id$
2010-07-13 13:49:25 +05:30
2011-03-31 14:51:43 +05:30
! states for titanmod
! Basic states
! rhoedge
! rhoscrew
! shear_system
! Dependent states
! segment_edge
! segment_screw
! resistance_edge
! resistance_screw
! tau_slip
! gdot_slip
! velocity_edge
! velocity_screw
! gdot_slip_edge
! gdot_slip_screw
! stressratio_edge_p
! stressratio_screw_p
! shear_basal
! shear_prism
! shear_pyra
! shear_pyrca
2010-07-13 13:49:25 +05:30
MODULE constitutive_titanmod
!* Include other modules
use prec , only : pReal , pInt
implicit none
!* Lists of states and physical parameters
character ( len = * ) , parameter :: constitutive_titanmod_label = 'titanmod'
2011-08-01 15:41:32 +05:30
character ( len = 18 ) , dimension ( 3 ) , parameter :: constitutive_titanmod_listBasicSlipStates = ( / 'rho_edge ' , &
'rho_screw ' , &
2011-03-31 14:51:43 +05:30
'shear_system' / )
2010-09-22 14:30:40 +05:30
character ( len = 18 ) , dimension ( 1 ) , parameter :: constitutive_titanmod_listBasicTwinStates = ( / 'gdot_twin' / )
2010-09-01 15:37:52 +05:30
2012-02-10 17:26:05 +05:30
character ( len = 19 ) , dimension ( 11 ) , parameter :: constitutive_titanmod_listDependentSlipStates = ( / 'segment_edge ' , &
2011-08-01 15:41:32 +05:30
'segment_screw ' , &
'resistance_edge ' , &
'resistance_screw ' , &
'tau_slip ' , &
'velocity_edge ' , &
'velocity_screw ' , &
'gdot_slip_edge ' , &
'gdot_slip_screw ' , &
'stressratio_edge_p ' , &
'stressratio_screw_p' &
2010-09-07 20:14:37 +05:30
/ )
2010-09-22 14:30:40 +05:30
character ( len = 18 ) , dimension ( 2 ) , parameter :: constitutive_titanmod_listDependentTwinStates = ( / 'twin_fraction' , &
2011-08-01 15:41:32 +05:30
'tau_twin ' &
2010-09-22 14:30:40 +05:30
/ )
2010-07-13 13:49:25 +05:30
real ( pReal ) , parameter :: kB = 1.38e-23_pReal ! Boltzmann constant in J/Kelvin
!* Definition of global variables
integer ( pInt ) , dimension ( : ) , allocatable :: constitutive_titanmod_sizeDotState , & ! number of dotStates
constitutive_titanmod_sizeState , & ! total number of microstructural state variables
constitutive_titanmod_sizePostResults ! cumulative size of post results
integer ( pInt ) , dimension ( : , : ) , allocatable , target :: constitutive_titanmod_sizePostResult ! size of each post result output
2012-02-14 20:49:59 +05:30
character ( len = 64 ) , dimension ( : , : ) , allocatable , target :: constitutive_titanmod_output ! name of each post result output
integer ( pInt ) , dimension ( : ) , allocatable :: constitutive_titanmod_Noutput ! number of outputs per instance of this constitution
2010-07-13 13:49:25 +05:30
character ( len = 32 ) , dimension ( : ) , allocatable :: constitutive_titanmod_structureName ! name of the lattice structure
integer ( pInt ) , dimension ( : ) , allocatable :: constitutive_titanmod_structure , & ! number representing the kind of lattice structure
constitutive_titanmod_totalNslip , & ! total number of active slip systems for each instance
constitutive_titanmod_totalNtwin ! total number of active twin systems for each instance
integer ( pInt ) , dimension ( : , : ) , allocatable :: constitutive_titanmod_Nslip , & ! number of active slip systems for each family and instance
constitutive_titanmod_Ntwin , & ! number of active twin systems for each family and instance
constitutive_titanmod_slipFamily , & ! lookup table relating active slip system to slip family for each instance
constitutive_titanmod_twinFamily , & ! lookup table relating active twin system to twin family for each instance
constitutive_titanmod_slipSystemLattice , & ! lookup table relating active slip system index to lattice slip system index for each instance
constitutive_titanmod_twinSystemLattice ! lookup table relating active twin system index to lattice twin system index for each instance
real ( pReal ) , dimension ( : ) , allocatable :: constitutive_titanmod_CoverA , & ! c/a ratio for hex type lattice
constitutive_titanmod_C11 , & ! C11 element in elasticity matrix
constitutive_titanmod_C12 , & ! C12 element in elasticity matrix
constitutive_titanmod_C13 , & ! C13 element in elasticity matrix
constitutive_titanmod_C33 , & ! C33 element in elasticity matrix
constitutive_titanmod_C44 , & ! C44 element in elasticity matrix
2010-09-29 12:05:08 +05:30
constitutive_titanmod_debyefrequency , & !Debye frequency
constitutive_titanmod_kinkf0 , & !Debye frequency
2010-07-13 13:49:25 +05:30
constitutive_titanmod_Gmod , & ! shear modulus
constitutive_titanmod_CAtomicVolume , & ! atomic volume in Bugers vector unit
2010-09-07 20:14:37 +05:30
constitutive_titanmod_dc , & ! prefactor for self-diffusion coefficient
2010-09-22 14:30:40 +05:30
constitutive_titanmod_twinhpconstant , & ! activation energy for dislocation climb
2010-07-15 12:46:15 +05:30
constitutive_titanmod_GrainSize , & ! grain size - Not being used
2010-07-13 13:49:25 +05:30
constitutive_titanmod_MaxTwinFraction , & ! maximum allowed total twin volume fraction
constitutive_titanmod_r , & ! r-exponent in twin nucleation rate
2010-07-15 12:46:15 +05:30
constitutive_titanmod_CEdgeDipMinDistance , & ! Not being used
constitutive_titanmod_Cmfptwin , & ! Not being used
constitutive_titanmod_Cthresholdtwin , & ! Not being used
2010-10-26 18:46:37 +05:30
constitutive_titanmod_aTolRho ! absolute tolerance for integration of dislocation density
2010-07-13 13:49:25 +05:30
real ( pReal ) , dimension ( : , : , : ) , allocatable :: constitutive_titanmod_Cslip_66 ! elasticity matrix in Mandel notation for each instance
real ( pReal ) , dimension ( : , : , : , : ) , allocatable :: constitutive_titanmod_Ctwin_66 ! twin elasticity matrix in Mandel notation for each instance
real ( pReal ) , dimension ( : , : , : , : , : ) , allocatable :: constitutive_titanmod_Cslip_3333 ! elasticity matrix for each instance
real ( pReal ) , dimension ( : , : , : , : , : , : ) , allocatable :: constitutive_titanmod_Ctwin_3333 ! twin elasticity matrix for each instance
2010-07-15 12:46:15 +05:30
real ( pReal ) , dimension ( : , : ) , allocatable :: constitutive_titanmod_rho_edge0 , & ! initial edge dislocation density per slip system for each family and instance
constitutive_titanmod_rho_screw0 , & ! initial screw dislocation density per slip system for each family and instance
2011-03-31 14:51:43 +05:30
constitutive_titanmod_shear_system0 , & ! accumulated shear on each system
2010-07-13 13:49:25 +05:30
constitutive_titanmod_burgersPerSlipFamily , & ! absolute length of burgers vector [m] for each slip family and instance
constitutive_titanmod_burgersPerSlipSystem , & ! absolute length of burgers vector [m] for each slip system and instance
constitutive_titanmod_burgersPerTwinFamily , & ! absolute length of burgers vector [m] for each twin family and instance
constitutive_titanmod_burgersPerTwinSystem , & ! absolute length of burgers vector [m] for each twin system and instance
2010-07-15 12:46:15 +05:30
constitutive_titanmod_f0_PerSlipFamily , & ! activation energy for glide [J] for each slip family and instance
constitutive_titanmod_f0_PerSlipSystem , & ! activation energy for glide [J] for each slip system and instance
2010-09-01 15:37:52 +05:30
constitutive_titanmod_twinf0_PerTwinFamily , & ! activation energy for glide [J] for each twin family and instance
constitutive_titanmod_twinf0_PerTwinSystem , & ! activation energy for glide [J] for each twin system and instance
2010-09-07 20:14:37 +05:30
constitutive_titanmod_twinshearconstant_PerTwinFamily , & ! activation energy for glide [J] for each twin family and instance
constitutive_titanmod_twinshearconstant_PerTwinSystem , & ! activation energy for glide [J] for each twin system and instance
2010-08-06 21:23:45 +05:30
constitutive_titanmod_tau0e_PerSlipFamily , & ! Initial yield stress for edge dislocations per slip family
constitutive_titanmod_tau0e_PerSlipSystem , & ! Initial yield stress for edge dislocations per slip system
constitutive_titanmod_tau0s_PerSlipFamily , & ! Initial yield stress for screw dislocations per slip family
constitutive_titanmod_tau0s_PerSlipSystem , & ! Initial yield stress for screw dislocations per slip system
2010-09-22 14:30:40 +05:30
constitutive_titanmod_twintau0_PerTwinFamily , & ! Initial yield stress for edge dislocations per twin family
constitutive_titanmod_twintau0_PerTwinSystem , & ! Initial yield stress for edge dislocations per twin system
2010-07-15 12:46:15 +05:30
constitutive_titanmod_capre_PerSlipFamily , & ! Capture radii for edge dislocations per slip family
constitutive_titanmod_capre_PerSlipSystem , & ! Capture radii for edge dislocations per slip system
constitutive_titanmod_caprs_PerSlipFamily , & ! Capture radii for screw dislocations per slip family
constitutive_titanmod_caprs_PerSlipSystem , & ! Capture radii for screw dislocations per slip system
constitutive_titanmod_pe_PerSlipFamily , & ! p-exponent in glide velocity
constitutive_titanmod_ps_PerSlipFamily , & ! p-exponent in glide velocity
constitutive_titanmod_qe_PerSlipFamily , & ! q-exponent in glide velocity
constitutive_titanmod_qs_PerSlipFamily , & ! q-exponent in glide velocity
constitutive_titanmod_pe_PerSlipSystem , & ! p-exponent in glide velocity
constitutive_titanmod_ps_PerSlipSystem , & ! p-exponent in glide velocity
constitutive_titanmod_qe_PerSlipSystem , & ! q-exponent in glide velocity
constitutive_titanmod_qs_PerSlipSystem , & ! q-exponent in glide velocity
2010-09-22 14:30:40 +05:30
constitutive_titanmod_twinp_PerTwinFamily , & ! p-exponent in glide velocity
constitutive_titanmod_twinq_PerTwinFamily , & ! q-exponent in glide velocity
constitutive_titanmod_twinp_PerTwinSystem , & ! p-exponent in glide velocity
constitutive_titanmod_twinq_PerTwinSystem , & ! p-exponent in glide velocity
2010-08-06 21:23:45 +05:30
constitutive_titanmod_v0e_PerSlipFamily , & ! edge dislocation velocity prefactor [m/s] for each family and instance
constitutive_titanmod_v0e_PerSlipSystem , & ! screw dislocation velocity prefactor [m/s] for each slip system and instance
constitutive_titanmod_v0s_PerSlipFamily , & ! edge dislocation velocity prefactor [m/s] for each family and instance
constitutive_titanmod_v0s_PerSlipSystem , & ! screw dislocation velocity prefactor [m/s] for each slip system and instance
2010-09-22 14:30:40 +05:30
constitutive_titanmod_twingamma0_PerTwinFamily , & ! edge dislocation velocity prefactor [m/s] for each family and instance
constitutive_titanmod_twingamma0_PerTwinSystem , & ! screw dislocation velocity prefactor [m/s] for each slip system and instance
constitutive_titanmod_kinkcriticallength_PerSlipFamily , & ! screw dislocation mobility prefactor for kink-pairs per slip family
constitutive_titanmod_kinkcriticallength_PerSlipSystem , & ! screw dislocation mobility prefactor for kink-pairs per slip system
2010-07-13 13:49:25 +05:30
constitutive_titanmod_twinsizePerTwinFamily , & ! twin thickness [m] for each twin family and instance
constitutive_titanmod_twinsizePerTwinSystem , & ! twin thickness [m] for each twin system and instance
2010-07-15 12:46:15 +05:30
constitutive_titanmod_CeLambdaSlipPerSlipFamily , & ! Adj. parameter for distance between 2 forest dislocations for each slip family and instance
constitutive_titanmod_CeLambdaSlipPerSlipSystem , & ! Adj. parameter for distance between 2 forest dislocations for each slip system and instance
constitutive_titanmod_CsLambdaSlipPerSlipFamily , & ! Adj. parameter for distance between 2 forest dislocations for each slip family and instance
constitutive_titanmod_CsLambdaSlipPerSlipSystem , & ! Adj. parameter for distance between 2 forest dislocations for each slip system and instance
2010-09-22 14:30:40 +05:30
constitutive_titanmod_twinLambdaSlipPerTwinFamily , & ! Adj. parameter for distance between 2 forest dislocations for each slip family and instance
constitutive_titanmod_twinLambdaSlipPerTwinSystem , & ! Adj. parameter for distance between 2 forest dislocations for each slip system and instance
2010-07-13 13:49:25 +05:30
constitutive_titanmod_interactionSlipSlip , & ! coefficients for slip-slip interaction for each interaction type and instance
2010-08-23 17:06:51 +05:30
constitutive_titanmod_interaction_ee , & ! coefficients for e-e interaction for each interaction type and instance
constitutive_titanmod_interaction_ss , & ! coefficients for s-s interaction for each interaction type and instance
constitutive_titanmod_interaction_es , & ! coefficients for e-s-twin interaction for each interaction type and instance
constitutive_titanmod_interactionSlipTwin , & ! coefficients for twin-slip interaction for each interaction type and instance
2010-07-13 13:49:25 +05:30
constitutive_titanmod_interactionTwinSlip , & ! coefficients for twin-slip interaction for each interaction type and instance
constitutive_titanmod_interactionTwinTwin ! coefficients for twin-twin interaction for each interaction type and instance
real ( pReal ) , dimension ( : , : , : ) , allocatable :: constitutive_titanmod_interactionMatrixSlipSlip , & ! interaction matrix of the different slip systems for each instance
2010-08-23 17:06:51 +05:30
constitutive_titanmod_interactionMatrix_ee , & ! interaction matrix of e-e for each instance
constitutive_titanmod_interactionMatrix_ss , & ! interaction matrix of s-s for each instance
constitutive_titanmod_interactionMatrix_es , & ! interaction matrix of e-s for each instance
2010-07-13 13:49:25 +05:30
constitutive_titanmod_interactionMatrixSlipTwin , & ! interaction matrix of slip systems with twin systems for each instance
constitutive_titanmod_interactionMatrixTwinSlip , & ! interaction matrix of twin systems with slip systems for each instance
constitutive_titanmod_interactionMatrixTwinTwin , & ! interaction matrix of the different twin systems for each instance
2010-08-06 21:23:45 +05:30
constitutive_titanmod_forestProjectionEdge , & ! matrix of forest projections of edge dislocations for each instance
2010-09-01 15:37:52 +05:30
constitutive_titanmod_forestProjectionScrew , & ! matrix of forest projections of screw dislocations for each instance
constitutive_titanmod_TwinforestProjectionEdge , & ! matrix of forest projections of edge dislocations in twin system for each instance
constitutive_titanmod_TwinforestProjectionScrew ! matrix of forest projections of screw dislocations in twin system for each instance
2010-07-13 13:49:25 +05:30
CONTAINS
!****************************************
!* - constitutive_titanmod_init
!* - constitutive_titanmod_stateInit
!* - constitutive_titanmod_relevantState
!* - constitutive_titanmod_homogenizedC
!* - constitutive_titanmod_microstructure
!* - constitutive_titanmod_LpAndItsTangent
!* - consistutive_titanmod_dotState
!* - constitutive_titanmod_dotTemperature
!* - consistutive_titanmod_postResults
!****************************************
2010-08-06 21:23:45 +05:30
2010-07-13 13:49:25 +05:30
subroutine constitutive_titanmod_init ( file )
!**************************************
!* Module initialization *
!**************************************
2012-02-21 21:30:00 +05:30
use , intrinsic :: iso_fortran_env ! to get compiler_version and compiler_options (at least for gfortran 4.6 at the moment)
2010-07-13 13:49:25 +05:30
use prec , only : pInt , pReal
use math , only : math_Mandel3333to66 , math_Voigt66to3333 , math_mul3x3
use IO
use material
use lattice
!* Input variables
integer ( pInt ) , intent ( in ) :: file
!* Local variables
2012-02-21 21:30:00 +05:30
integer ( pInt ) , parameter :: maxNchunks = 21_pInt
integer ( pInt ) , dimension ( 1_pInt + 2_pInt * maxNchunks ) :: positions
integer ( pInt ) section , f , i , j , k , l , m , n , o , p , q , r , s , s1 , s2 , t , t1 , t2 , ns , nt , mySize , myStructure , maxTotalNslip , &
2010-08-23 17:06:51 +05:30
maxTotalNtwin
2012-02-21 21:30:00 +05:30
integer :: maxNinstance !no pInt
2010-07-13 13:49:25 +05:30
character ( len = 64 ) tag
character ( len = 1024 ) line
2010-09-22 14:30:40 +05:30
write ( 6 , * )
2011-08-26 19:27:29 +05:30
write ( 6 , * ) '<<<+- constitutive_' , trim ( constitutive_titanmod_label ) , ' init -+>>>'
2010-09-23 13:35:50 +05:30
write ( 6 , * ) '$Id$'
2012-02-01 00:48:55 +05:30
#include "compilation_info.f90"
2010-07-13 13:49:25 +05:30
maxNinstance = count ( phase_constitution == constitutive_titanmod_label )
if ( maxNinstance == 0 ) return
!* Space allocation for global variables
2012-02-21 21:30:00 +05:30
allocate ( constitutive_titanmod_sizeDotState ( maxNinstance ) )
constitutive_titanmod_sizeDotState = 0_pInt
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_sizeState ( maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_sizeState = 0_pInt
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_sizePostResults ( maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_sizePostResults = 0_pInt
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_sizePostResult ( maxval ( phase_Noutput ) , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_sizePostResult = 0_pInt
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_output ( maxval ( phase_Noutput ) , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_output = ''
2012-02-14 20:49:59 +05:30
allocate ( constitutive_titanmod_Noutput ( maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_Noutput = 0_pInt
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_structureName ( maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_structureName = ''
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_structure ( maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_structure = 0_pInt
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_Nslip ( lattice_maxNslipFamily , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_Nslip = 0_pInt
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_Ntwin ( lattice_maxNtwinFamily , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_Ntwin = 0_pInt
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_slipFamily ( lattice_maxNslip , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_slipFamily = 0_pInt
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_twinFamily ( lattice_maxNtwin , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_twinFamily = 0_pInt
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_slipSystemLattice ( lattice_maxNslip , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_slipSystemLattice = 0_pInt
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_twinSystemLattice ( lattice_maxNtwin , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_twinSystemLattice = 0_pInt
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_totalNslip ( maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_totalNslip = 0_pInt
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_totalNtwin ( maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_totalNtwin = 0_pInt
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_CoverA ( maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_CoverA = 0.0_pReal
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_C11 ( maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_C11 = 0.0_pReal
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_C12 ( maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_C12 = 0.0_pReal
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_C13 ( maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_C13 = 0.0_pReal
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_C33 ( maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_C33 = 0.0_pReal
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_C44 ( maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_C44 = 0.0_pReal
2010-09-29 12:05:08 +05:30
allocate ( constitutive_titanmod_debyefrequency ( maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_debyefrequency = 0.0_pReal
2010-09-29 12:05:08 +05:30
allocate ( constitutive_titanmod_kinkf0 ( maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_kinkf0 = 0.0_pReal
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_Gmod ( maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_Gmod = 0.0_pReal
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_CAtomicVolume ( maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_CAtomicVolume = 0.0_pReal
2010-09-07 20:14:37 +05:30
allocate ( constitutive_titanmod_dc ( maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_dc = 0.0_pReal
2010-09-07 20:14:37 +05:30
allocate ( constitutive_titanmod_twinhpconstant ( maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_twinhpconstant = 0.0_pReal
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_GrainSize ( maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_GrainSize = 0.0_pReal
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_MaxTwinFraction ( maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_MaxTwinFraction = 0.0_pReal
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_r ( maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_r = 0.0_pReal
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_CEdgeDipMinDistance ( maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_CEdgeDipMinDistance = 0.0_pReal
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_Cmfptwin ( maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_Cmfptwin = 0.0_pReal
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_Cthresholdtwin ( maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_Cthresholdtwin = 0.0_pReal
2010-10-26 18:46:37 +05:30
allocate ( constitutive_titanmod_aTolRho ( maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_aTolRho = 0.0_pReal
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_Cslip_66 ( 6 , 6 , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_Cslip_66 = 0.0_pReal
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_Cslip_3333 ( 3 , 3 , 3 , 3 , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_Cslip_3333 = 0.0_pReal
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_rho_edge0 ( lattice_maxNslipFamily , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_rho_edge0 = 0.0_pReal
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_rho_screw0 ( lattice_maxNslipFamily , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_rho_screw0 = 0.0_pReal
2011-03-31 14:51:43 +05:30
allocate ( constitutive_titanmod_shear_system0 ( lattice_maxNslipFamily , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_shear_system0 = 0.0_pReal
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_burgersPerSlipFamily ( lattice_maxNslipFamily , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_burgersPerSlipFamily = 0.0_pReal
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_burgersPerTwinFamily ( lattice_maxNtwinFamily , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_burgersPerTwinFamily = 0.0_pReal
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_f0_PerSlipFamily ( lattice_maxNslipFamily , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_f0_PerSlipFamily = 0.0_pReal
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_tau0e_PerSlipFamily ( lattice_maxNslipFamily , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_tau0e_PerSlipFamily = 0.0_pReal
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_tau0s_PerSlipFamily ( lattice_maxNslipFamily , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_tau0s_PerSlipFamily = 0.0_pReal
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_capre_PerSlipFamily ( lattice_maxNslipFamily , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_capre_PerSlipFamily = 0.0_pReal
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_caprs_PerSlipFamily ( lattice_maxNslipFamily , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_caprs_PerSlipFamily = 0.0_pReal
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_pe_PerSlipFamily ( lattice_maxNslipFamily , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_pe_PerSlipFamily = 0.0_pReal
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_ps_PerSlipFamily ( lattice_maxNslipFamily , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_ps_PerSlipFamily = 0.0_pReal
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_qe_PerSlipFamily ( lattice_maxNslipFamily , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_qe_PerSlipFamily = 0.0_pReal
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_qs_PerSlipFamily ( lattice_maxNslipFamily , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_qs_PerSlipFamily = 0.0_pReal
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_v0e_PerSlipFamily ( lattice_maxNslipFamily , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_v0e_PerSlipFamily = 0.0_pReal
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_v0s_PerSlipFamily ( lattice_maxNslipFamily , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_v0s_PerSlipFamily = 0.0_pReal
2010-09-22 14:30:40 +05:30
allocate ( constitutive_titanmod_kinkcriticallength_PerSlipFamily ( lattice_maxNslipFamily , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_kinkcriticallength_PerSlipFamily = 0.0_pReal
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_twinsizePerTwinFamily ( lattice_maxNtwinFamily , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_twinsizePerTwinFamily = 0.0_pReal
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_CeLambdaSlipPerSlipFamily ( lattice_maxNslipFamily , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_CeLambdaSlipPerSlipFamily = 0.0_pReal
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_CsLambdaSlipPerSlipFamily ( lattice_maxNslipFamily , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_CsLambdaSlipPerSlipFamily = 0.0_pReal
2010-09-01 15:37:52 +05:30
allocate ( constitutive_titanmod_twinf0_PerTwinFamily ( lattice_maxNTwinFamily , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_twinf0_PerTwinFamily = 0.0_pReal
2010-09-07 20:14:37 +05:30
allocate ( constitutive_titanmod_twinshearconstant_PerTwinFamily ( lattice_maxNTwinFamily , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_twinshearconstant_PerTwinFamily = 0.0_pReal
2010-09-22 14:30:40 +05:30
allocate ( constitutive_titanmod_twintau0_PerTwinFamily ( lattice_maxNTwinFamily , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_twintau0_PerTwinFamily = 0.0_pReal
2010-09-22 14:30:40 +05:30
allocate ( constitutive_titanmod_twinp_PerTwinFamily ( lattice_maxNTwinFamily , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_twingamma0_PerTwinFamily = 0.0_pReal
2010-09-22 14:30:40 +05:30
allocate ( constitutive_titanmod_twinq_PerTwinFamily ( lattice_maxNTwinFamily , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_twinLambdaSlipPerTwinFamily = 0.0_pReal
2010-09-22 14:30:40 +05:30
allocate ( constitutive_titanmod_twingamma0_PerTwinFamily ( lattice_maxNTwinFamily , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_twinp_PerTwinFamily = 0.0_pReal
2010-09-22 14:30:40 +05:30
allocate ( constitutive_titanmod_twinLambdaSlipPerTwinFamily ( lattice_maxNTwinFamily , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_twinq_PerTwinFamily = 0.0_pReal
2010-09-01 15:37:52 +05:30
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_interactionSlipSlip ( lattice_maxNinteraction , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_interactionSlipSlip = 0.0_pReal
2010-08-23 17:06:51 +05:30
allocate ( constitutive_titanmod_interaction_ee ( lattice_maxNinteraction , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_interaction_ee = 0.0_pReal
2010-08-23 17:06:51 +05:30
allocate ( constitutive_titanmod_interaction_ss ( lattice_maxNinteraction , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_interaction_ss = 0.0_pReal
2010-08-23 17:06:51 +05:30
allocate ( constitutive_titanmod_interaction_es ( lattice_maxNinteraction , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_interaction_ss = 0.0_pReal
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_interactionSlipTwin ( lattice_maxNinteraction , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_interactionSlipTwin = 0.0_pReal
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_interactionTwinSlip ( lattice_maxNinteraction , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_interactionTwinSlip = 0.0_pReal
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_interactionTwinTwin ( lattice_maxNinteraction , maxNinstance ) )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_interactionTwinTwin = 0.0_pReal
2010-07-13 13:49:25 +05:30
!* Readout data from material.config file
rewind ( file )
line = ''
2012-02-21 21:30:00 +05:30
section = 0_pInt
2010-07-13 13:49:25 +05:30
2010-09-22 14:30:40 +05:30
write ( 6 , * ) 'titanmod: Reading material parameters from material config file'
2010-07-13 13:49:25 +05:30
do while ( IO_lc ( IO_getTag ( line , '<' , '>' ) ) / = 'phase' ) ! wind forward to <phase>
read ( file , '(a1024)' , END = 100 ) line
enddo
2012-02-14 05:00:59 +05:30
do ! read thru sections of phase part
2010-07-13 13:49:25 +05:30
read ( file , '(a1024)' , END = 100 ) line
if ( IO_isBlank ( line ) ) cycle ! skip empty lines
if ( IO_getTag ( line , '<' , '>' ) / = '' ) exit ! stop at next part
if ( IO_getTag ( line , '[' , ']' ) / = '' ) then ! next section
2012-02-14 05:00:59 +05:30
section = section + 1_pInt ! advance section counter
cycle ! skip to next line
2010-07-13 13:49:25 +05:30
endif
2012-02-13 19:48:07 +05:30
if ( section > 0_pInt . and . phase_constitution ( section ) == constitutive_titanmod_label ) then ! one of my sections
2012-02-14 05:00:59 +05:30
i = phase_constitutionInstance ( section ) ! which instance of my constitution is present phase
2010-07-13 13:49:25 +05:30
positions = IO_stringPos ( line , maxNchunks )
2012-02-14 05:00:59 +05:30
tag = IO_lc ( IO_stringValue ( line , positions , 1_pInt ) ) ! extract key
2010-07-13 13:49:25 +05:30
select case ( tag )
2012-02-14 14:52:37 +05:30
case ( 'constitution' )
cycle
2010-07-13 13:49:25 +05:30
case ( '(output)' )
2012-02-14 20:49:59 +05:30
constitutive_titanmod_Noutput ( i ) = constitutive_titanmod_Noutput ( i ) + 1_pInt
constitutive_titanmod_output ( constitutive_titanmod_Noutput ( i ) , i ) = IO_lc ( IO_stringValue ( line , positions , 2_pInt ) )
write ( 6 , * ) tag , constitutive_titanmod_output ( constitutive_titanmod_Noutput ( i ) , i )
2010-09-01 15:37:52 +05:30
case ( 'lattice_structure' )
2012-02-13 19:48:07 +05:30
constitutive_titanmod_structureName ( i ) = IO_lc ( IO_stringValue ( line , positions , 2_pInt ) )
2010-09-22 14:30:40 +05:30
write ( 6 , * ) tag , constitutive_titanmod_structureName ( i )
2010-07-13 13:49:25 +05:30
case ( 'covera_ratio' )
2012-02-13 19:48:07 +05:30
constitutive_titanmod_CoverA ( i ) = IO_floatValue ( line , positions , 2_pInt )
2010-09-22 14:30:40 +05:30
write ( 6 , * ) tag , constitutive_titanmod_CoverA ( i )
2010-07-13 13:49:25 +05:30
case ( 'c11' )
2012-02-13 19:48:07 +05:30
constitutive_titanmod_C11 ( i ) = IO_floatValue ( line , positions , 2_pInt )
2010-08-23 17:06:51 +05:30
write ( 6 , * ) tag , constitutive_titanmod_C11 ( i )
2010-07-13 13:49:25 +05:30
case ( 'c12' )
2012-02-13 19:48:07 +05:30
constitutive_titanmod_C12 ( i ) = IO_floatValue ( line , positions , 2_pInt )
2010-08-23 17:06:51 +05:30
write ( 6 , * ) tag , constitutive_titanmod_C12 ( i )
2010-07-13 13:49:25 +05:30
case ( 'c13' )
2012-02-13 19:48:07 +05:30
constitutive_titanmod_C13 ( i ) = IO_floatValue ( line , positions , 2_pInt )
2010-08-23 17:06:51 +05:30
write ( 6 , * ) tag , constitutive_titanmod_C13 ( i )
2010-07-13 13:49:25 +05:30
case ( 'c33' )
2012-02-13 19:48:07 +05:30
constitutive_titanmod_C33 ( i ) = IO_floatValue ( line , positions , 2_pInt )
2010-08-23 17:06:51 +05:30
write ( 6 , * ) tag , constitutive_titanmod_C33 ( i )
2010-07-13 13:49:25 +05:30
case ( 'c44' )
2012-02-13 19:48:07 +05:30
constitutive_titanmod_C44 ( i ) = IO_floatValue ( line , positions , 2_pInt )
2010-08-23 17:06:51 +05:30
write ( 6 , * ) tag , constitutive_titanmod_C44 ( i )
2010-09-29 12:05:08 +05:30
case ( 'debyefrequency' )
2012-02-13 19:48:07 +05:30
constitutive_titanmod_debyefrequency ( i ) = IO_floatValue ( line , positions , 2_pInt )
2010-09-29 12:05:08 +05:30
write ( 6 , * ) tag , constitutive_titanmod_debyefrequency ( i )
case ( 'kinkf0' )
2012-02-13 19:48:07 +05:30
constitutive_titanmod_kinkf0 ( i ) = IO_floatValue ( line , positions , 2_pInt )
2010-09-29 12:05:08 +05:30
write ( 6 , * ) tag , constitutive_titanmod_kinkf0 ( i )
2010-07-13 13:49:25 +05:30
case ( 'nslip' )
2012-02-13 19:48:07 +05:30
forall ( j = 1_pInt : lattice_maxNslipFamily ) &
constitutive_titanmod_Nslip ( j , i ) = IO_intValue ( line , positions , 1_pInt + j )
2010-09-22 14:30:40 +05:30
write ( 6 , * ) tag , constitutive_titanmod_Nslip ( 1 : 4 , i )
2010-07-13 13:49:25 +05:30
case ( 'ntwin' )
2012-02-13 19:48:07 +05:30
forall ( j = 1_pInt : lattice_maxNtwinFamily ) &
constitutive_titanmod_Ntwin ( j , i ) = IO_intValue ( line , positions , 1_pInt + j )
2010-09-22 14:30:40 +05:30
write ( 6 , * ) tag , constitutive_titanmod_Ntwin ( 1 : 4 , i )
2010-07-13 13:49:25 +05:30
case ( 'rho_edge0' )
2012-02-13 19:48:07 +05:30
forall ( j = 1_pInt : lattice_maxNslipFamily ) &
constitutive_titanmod_rho_edge0 ( j , i ) = IO_floatValue ( line , positions , 1_pInt + j )
2010-09-22 14:30:40 +05:30
write ( 6 , * ) tag , constitutive_titanmod_rho_edge0 ( 1 : 4 , i )
2010-07-13 13:49:25 +05:30
case ( 'rho_screw0' )
2012-02-21 21:30:00 +05:30
forall ( j = 1_pInt : lattice_maxNslipFamily ) &
2012-02-13 19:48:07 +05:30
constitutive_titanmod_rho_screw0 ( j , i ) = IO_floatValue ( line , positions , 1_pInt + j )
2010-09-22 14:30:40 +05:30
write ( 6 , * ) tag , constitutive_titanmod_rho_screw0 ( 1 : 4 , i )
2010-07-13 13:49:25 +05:30
case ( 'slipburgers' )
2012-02-13 19:48:07 +05:30
forall ( j = 1_pInt : lattice_maxNslipFamily ) &
constitutive_titanmod_burgersPerSlipFamily ( j , i ) = IO_floatValue ( line , positions , 1_pInt + j )
2010-09-22 14:30:40 +05:30
write ( 6 , * ) tag , constitutive_titanmod_burgersPerSlipFamily ( 1 : 4 , i )
2010-07-13 13:49:25 +05:30
case ( 'twinburgers' )
2012-02-13 19:48:07 +05:30
forall ( j = 1_pInt : lattice_maxNtwinFamily ) &
constitutive_titanmod_burgersPerTwinFamily ( j , i ) = IO_floatValue ( line , positions , 1_pInt + j )
2010-09-22 14:30:40 +05:30
write ( 6 , * ) tag , constitutive_titanmod_burgersPerTwinFamily ( 1 : 4 , i )
2010-07-13 13:49:25 +05:30
case ( 'f0' )
2012-02-13 19:48:07 +05:30
forall ( j = 1_pInt : lattice_maxNslipFamily ) &
constitutive_titanmod_f0_PerSlipFamily ( j , i ) = IO_floatValue ( line , positions , 1_pInt + j )
2010-09-22 14:30:40 +05:30
write ( 6 , * ) tag , constitutive_titanmod_f0_PerSlipFamily ( 1 : 4 , i )
2010-09-01 15:37:52 +05:30
case ( 'twinf0' )
2012-02-13 19:48:07 +05:30
forall ( j = 1_pInt : lattice_maxNtwinFamily ) &
constitutive_titanmod_twinf0_PerTwinFamily ( j , i ) = IO_floatValue ( line , positions , 1_pInt + j )
2010-09-22 14:30:40 +05:30
write ( 6 , * ) tag , constitutive_titanmod_twinf0_PerTwinFamily ( 1 : 4 , i )
2010-07-13 13:49:25 +05:30
case ( 'tau0e' )
2012-02-13 19:48:07 +05:30
forall ( j = 1_pInt : lattice_maxNslipFamily ) &
constitutive_titanmod_tau0e_PerSlipFamily ( j , i ) = IO_floatValue ( line , positions , 1_pInt + j )
2010-09-22 14:30:40 +05:30
write ( 6 , * ) tag , constitutive_titanmod_tau0e_PerSlipFamily ( 1 : 4 , i )
case ( 'twintau0' )
2012-02-13 19:48:07 +05:30
forall ( j = 1_pInt : lattice_maxNtwinFamily ) &
constitutive_titanmod_twintau0_PerTwinFamily ( j , i ) = IO_floatValue ( line , positions , 1_pInt + j )
2010-09-22 14:30:40 +05:30
write ( 6 , * ) tag , constitutive_titanmod_twintau0_PerTwinFamily ( 1 : 4 , i )
2010-07-13 13:49:25 +05:30
case ( 'tau0s' )
2012-02-13 19:48:07 +05:30
forall ( j = 1_pInt : lattice_maxNslipFamily ) &
constitutive_titanmod_tau0s_PerSlipFamily ( j , i ) = IO_floatValue ( line , positions , 1_pInt + j )
2010-09-22 14:30:40 +05:30
write ( 6 , * ) tag , constitutive_titanmod_tau0s_PerSlipFamily ( 1 : 4 , i )
2010-07-13 13:49:25 +05:30
case ( 'capre' )
2012-02-13 19:48:07 +05:30
forall ( j = 1_pInt : lattice_maxNslipFamily ) &
constitutive_titanmod_capre_PerSlipFamily ( j , i ) = IO_floatValue ( line , positions , 1_pInt + j )
2010-09-22 14:30:40 +05:30
write ( 6 , * ) tag , constitutive_titanmod_capre_PerSlipFamily ( 1 : 4 , i )
2010-07-13 13:49:25 +05:30
case ( 'caprs' )
2012-02-13 19:48:07 +05:30
forall ( j = 1_pInt : lattice_maxNslipFamily ) &
constitutive_titanmod_caprs_PerSlipFamily ( j , i ) = IO_floatValue ( line , positions , 1_pInt + j )
2010-09-22 14:30:40 +05:30
write ( 6 , * ) tag , constitutive_titanmod_caprs_PerSlipFamily ( 1 : 4 , i )
2010-07-13 13:49:25 +05:30
case ( 'v0e' )
2012-02-13 19:48:07 +05:30
forall ( j = 1_pInt : lattice_maxNslipFamily ) &
constitutive_titanmod_v0e_PerSlipFamily ( j , i ) = IO_floatValue ( line , positions , 1_pInt + j )
2010-09-22 14:30:40 +05:30
write ( 6 , * ) tag , constitutive_titanmod_v0e_PerSlipFamily ( 1 : 4 , i )
case ( 'twingamma0' )
2012-02-13 19:48:07 +05:30
forall ( j = 1_pInt : lattice_maxNtwinFamily ) &
constitutive_titanmod_twingamma0_PerTwinFamily ( j , i ) = IO_floatValue ( line , positions , 1_pInt + j )
2010-09-22 14:30:40 +05:30
write ( 6 , * ) tag , constitutive_titanmod_twingamma0_PerTwinFamily ( 1 : 4 , i )
2010-07-13 13:49:25 +05:30
case ( 'v0s' )
2012-02-13 19:48:07 +05:30
forall ( j = 1_pInt : lattice_maxNslipFamily ) &
constitutive_titanmod_v0s_PerSlipFamily ( j , i ) = IO_floatValue ( line , positions , 1_pInt + j )
2010-09-22 14:30:40 +05:30
write ( 6 , * ) tag , constitutive_titanmod_v0s_PerSlipFamily ( 1 : 4 , i )
case ( 'kinkcriticallength' )
2012-02-13 19:48:07 +05:30
forall ( j = 1_pInt : lattice_maxNslipFamily ) &
constitutive_titanmod_kinkcriticallength_PerSlipFamily ( j , i ) = IO_floatValue ( line , positions , 1_pInt + j )
2010-09-22 14:30:40 +05:30
write ( 6 , * ) tag , constitutive_titanmod_kinkcriticallength_PerSlipFamily ( 1 : 4 , i )
2010-07-13 13:49:25 +05:30
case ( 'twinsize' )
2012-02-13 19:48:07 +05:30
forall ( j = 1_pInt : lattice_maxNtwinFamily ) &
constitutive_titanmod_twinsizePerTwinFamily ( j , i ) = IO_floatValue ( line , positions , 1_pInt + j )
2010-08-23 17:06:51 +05:30
write ( 6 , * ) tag
2010-07-13 13:49:25 +05:30
case ( 'celambdaslip' )
2012-02-13 19:48:07 +05:30
forall ( j = 1_pInt : lattice_maxNslipFamily ) &
constitutive_titanmod_CeLambdaSlipPerSlipFamily ( j , i ) = IO_floatValue ( line , positions , 1_pInt + j )
2010-08-23 17:06:51 +05:30
write ( 6 , * ) tag
2010-09-22 14:30:40 +05:30
case ( 'twinlambdaslip' )
2012-02-13 19:48:07 +05:30
forall ( j = 1_pInt : lattice_maxNtwinFamily ) &
constitutive_titanmod_twinlambdaslipPerTwinFamily ( j , i ) = IO_floatValue ( line , positions , 1_pInt + j )
2010-09-22 14:30:40 +05:30
write ( 6 , * ) tag , constitutive_titanmod_twinlambdaslipPerTwinFamily ( 1 : 4 , i )
2010-07-13 13:49:25 +05:30
case ( 'cslambdaslip' )
2012-02-13 19:48:07 +05:30
forall ( j = 1_pInt : lattice_maxNslipFamily ) &
constitutive_titanmod_CsLambdaSlipPerSlipFamily ( j , i ) = IO_floatValue ( line , positions , 1_pInt + j )
2010-08-23 17:06:51 +05:30
write ( 6 , * ) tag
2010-07-13 13:49:25 +05:30
case ( 'grainsize' )
2012-02-13 19:48:07 +05:30
constitutive_titanmod_GrainSize ( i ) = IO_floatValue ( line , positions , 2_pInt )
2010-08-23 17:06:51 +05:30
write ( 6 , * ) tag
2010-07-13 13:49:25 +05:30
case ( 'maxtwinfraction' )
2012-02-13 19:48:07 +05:30
constitutive_titanmod_MaxTwinFraction ( i ) = IO_floatValue ( line , positions , 2_pInt )
2010-08-23 17:06:51 +05:30
write ( 6 , * ) tag
2010-07-13 13:49:25 +05:30
case ( 'pe' )
2012-02-13 19:48:07 +05:30
forall ( j = 1_pInt : lattice_maxNslipFamily ) &
constitutive_titanmod_pe_PerSlipFamily ( j , i ) = IO_floatValue ( line , positions , 1_pInt + j )
2010-09-22 14:30:40 +05:30
write ( 6 , * ) tag , constitutive_titanmod_pe_PerSlipFamily ( 1 : 4 , i )
case ( 'twinp' )
2012-02-13 19:48:07 +05:30
forall ( j = 1_pInt : lattice_maxNtwinFamily ) &
constitutive_titanmod_twinp_PerTwinFamily ( j , i ) = IO_floatValue ( line , positions , 1_pInt + j )
2010-09-22 14:30:40 +05:30
write ( 6 , * ) tag , constitutive_titanmod_twinp_PerTwinFamily ( 1 : 4 , i )
2010-07-13 13:49:25 +05:30
case ( 'ps' )
2012-02-13 19:48:07 +05:30
forall ( j = 1_pInt : lattice_maxNslipFamily ) &
constitutive_titanmod_ps_PerSlipFamily ( j , i ) = IO_floatValue ( line , positions , 1_pInt + j )
2010-09-22 14:30:40 +05:30
write ( 6 , * ) tag , constitutive_titanmod_ps_PerSlipFamily ( 1 : 4 , i )
2010-07-13 13:49:25 +05:30
case ( 'qe' )
2012-02-13 19:48:07 +05:30
forall ( j = 1_pInt : lattice_maxNslipFamily ) &
constitutive_titanmod_qe_PerSlipFamily ( j , i ) = IO_floatValue ( line , positions , 1_pInt + j )
2010-09-22 14:30:40 +05:30
write ( 6 , * ) tag , constitutive_titanmod_qe_PerSlipFamily ( 1 : 4 , i )
case ( 'twinq' )
2012-02-13 19:48:07 +05:30
forall ( j = 1_pInt : lattice_maxNtwinFamily ) &
constitutive_titanmod_twinq_PerTwinFamily ( j , i ) = IO_floatValue ( line , positions , 1_pInt + j )
2010-09-22 14:30:40 +05:30
write ( 6 , * ) tag , constitutive_titanmod_twinq_PerTwinFamily ( 1 : 4 , i )
2010-07-13 13:49:25 +05:30
case ( 'qs' )
2012-02-13 19:48:07 +05:30
forall ( j = 1_pInt : lattice_maxNslipFamily ) &
constitutive_titanmod_qs_PerSlipFamily ( j , i ) = IO_floatValue ( line , positions , 1_pInt + j )
2010-09-22 14:30:40 +05:30
write ( 6 , * ) tag , constitutive_titanmod_qs_PerSlipFamily ( 1 : 4 , i )
2010-09-07 20:14:37 +05:30
case ( 'twinshearconstant' )
2012-02-13 19:48:07 +05:30
forall ( j = 1_pInt : lattice_maxNtwinFamily ) &
constitutive_titanmod_twinshearconstant_PerTwinFamily ( j , i ) = IO_floatValue ( line , positions , 1_pInt + j )
2010-09-22 14:30:40 +05:30
write ( 6 , * ) tag , constitutive_titanmod_twinshearconstant_PerTwinFamily ( 1 : 4 , i )
2010-09-07 20:14:37 +05:30
case ( 'dc' )
2012-02-13 19:48:07 +05:30
constitutive_titanmod_dc ( i ) = IO_floatValue ( line , positions , 2_pInt )
2010-08-23 17:06:51 +05:30
write ( 6 , * ) tag
2010-09-07 20:14:37 +05:30
case ( 'twinhpconstant' )
2012-02-13 19:48:07 +05:30
constitutive_titanmod_twinhpconstant ( i ) = IO_floatValue ( line , positions , 2_pInt )
2010-08-23 17:06:51 +05:30
write ( 6 , * ) tag
2010-10-26 18:46:37 +05:30
case ( 'atol_rho' )
2012-02-13 19:48:07 +05:30
constitutive_titanmod_aTolRho ( i ) = IO_floatValue ( line , positions , 2_pInt )
2010-08-23 17:06:51 +05:30
write ( 6 , * ) tag
2010-07-13 13:49:25 +05:30
case ( 'interactionslipslip' )
2012-02-13 19:48:07 +05:30
forall ( j = 1_pInt : lattice_maxNinteraction ) &
constitutive_titanmod_interactionSlipSlip ( j , i ) = IO_floatValue ( line , positions , 1_pInt + j )
2010-08-23 17:06:51 +05:30
write ( 6 , * ) tag
case ( 'interactionee' )
2012-02-13 19:48:07 +05:30
forall ( j = 1_pInt : lattice_maxNinteraction ) &
constitutive_titanmod_interaction_ee ( j , i ) = IO_floatValue ( line , positions , 1_pInt + j )
2010-08-23 17:06:51 +05:30
write ( 6 , * ) tag
case ( 'interactionss' )
2012-02-13 19:48:07 +05:30
forall ( j = 1_pInt : lattice_maxNinteraction ) &
constitutive_titanmod_interaction_ss ( j , i ) = IO_floatValue ( line , positions , 1_pInt + j )
2010-08-23 17:06:51 +05:30
write ( 6 , * ) tag
case ( 'interactiones' )
2012-02-13 19:48:07 +05:30
forall ( j = 1_pInt : lattice_maxNinteraction ) &
constitutive_titanmod_interaction_es ( j , i ) = IO_floatValue ( line , positions , 1_pInt + j )
2010-08-23 17:06:51 +05:30
write ( 6 , * ) tag
2010-07-13 13:49:25 +05:30
case ( 'interactionsliptwin' )
2012-02-13 19:48:07 +05:30
forall ( j = 1_pInt : lattice_maxNinteraction ) &
constitutive_titanmod_interactionSlipTwin ( j , i ) = IO_floatValue ( line , positions , 1_pInt + j )
2010-08-23 17:06:51 +05:30
write ( 6 , * ) tag
2010-07-13 13:49:25 +05:30
case ( 'interactiontwinslip' )
2012-02-13 19:48:07 +05:30
forall ( j = 1_pInt : lattice_maxNinteraction ) &
constitutive_titanmod_interactionTwinSlip ( j , i ) = IO_floatValue ( line , positions , 1_pInt + j )
2010-08-23 17:06:51 +05:30
write ( 6 , * ) tag
2010-07-13 13:49:25 +05:30
case ( 'interactiontwintwin' )
2012-02-13 19:48:07 +05:30
forall ( j = 1_pInt : lattice_maxNinteraction ) &
constitutive_titanmod_interactionTwinTwin ( j , i ) = IO_floatValue ( line , positions , 1_pInt + j )
2010-08-23 17:06:51 +05:30
write ( 6 , * ) tag
2012-02-14 14:52:37 +05:30
case default
call IO_error ( 230_pInt , ext_msg = tag )
2010-07-13 13:49:25 +05:30
end select
endif
enddo
write ( 6 , * ) 'Material Property reading done'
2012-02-13 19:48:07 +05:30
100 do i = 1_pInt , maxNinstance
2010-07-13 13:49:25 +05:30
constitutive_titanmod_structure ( i ) = &
lattice_initializeStructure ( constitutive_titanmod_structureName ( i ) , constitutive_titanmod_CoverA ( i ) )
myStructure = constitutive_titanmod_structure ( i )
!* Sanity checks
2012-02-13 23:11:27 +05:30
if ( myStructure < 1_pInt . or . myStructure > 3_pInt ) call IO_error ( 205_pInt , e = i )
if ( sum ( constitutive_titanmod_Nslip ( : , i ) ) < = 0_pInt ) call IO_error ( 231_pInt , e = i , ext_msg = 'nslip' )
if ( sum ( constitutive_titanmod_Ntwin ( : , i ) ) < 0_pInt ) call IO_error ( 231_pInt , e = i , ext_msg = 'ntwin' )
2012-02-21 21:30:00 +05:30
do f = 1_pInt , lattice_maxNslipFamily
2010-07-13 13:49:25 +05:30
if ( constitutive_titanmod_Nslip ( f , i ) > 0_pInt ) then
2012-02-13 23:11:27 +05:30
if ( constitutive_titanmod_rho_edge0 ( f , i ) < 0.0_pReal ) call IO_error ( 231_pInt , e = i , ext_msg = 'rho_edge0' )
if ( constitutive_titanmod_rho_screw0 ( f , i ) < 0.0_pReal ) call IO_error ( 231_pInt , e = i , ext_msg = 'rho_screw0' )
if ( constitutive_titanmod_burgersPerSlipFamily ( f , i ) < = 0.0_pReal ) call IO_error ( 231_pInt , e = i , ext_msg = 'slipburgers' )
if ( constitutive_titanmod_f0_PerSlipFamily ( f , i ) < = 0.0_pReal ) call IO_error ( 231_pInt , e = i , ext_msg = 'f0' )
if ( constitutive_titanmod_tau0e_PerSlipFamily ( f , i ) < = 0.0_pReal ) call IO_error ( 231_pInt , e = i , ext_msg = 'tau0e' )
if ( constitutive_titanmod_tau0s_PerSlipFamily ( f , i ) < = 0.0_pReal ) call IO_error ( 231_pInt , e = i , ext_msg = 'tau0s' )
if ( constitutive_titanmod_capre_PerSlipFamily ( f , i ) < = 0.0_pReal ) call IO_error ( 231_pInt , e = i , ext_msg = 'capre' )
if ( constitutive_titanmod_caprs_PerSlipFamily ( f , i ) < = 0.0_pReal ) call IO_error ( 231_pInt , e = i , ext_msg = 'carrs' )
if ( constitutive_titanmod_v0e_PerSlipFamily ( f , i ) < = 0.0_pReal ) call IO_error ( 231_pInt , e = i , ext_msg = 'v0e' )
if ( constitutive_titanmod_v0s_PerSlipFamily ( f , i ) < = 0.0_pReal ) call IO_error ( 231_pInt , e = i , ext_msg = 'v0s' )
if ( constitutive_titanmod_kinkcriticallength_PerSlipFamily ( f , i ) < = 0.0_pReal ) &
call IO_error ( 231_pInt , e = i , ext_msg = 'kinkCriticalLength' )
2010-07-13 13:49:25 +05:30
endif
enddo
2012-02-21 21:30:00 +05:30
do f = 1_pInt , lattice_maxNtwinFamily
2010-09-01 15:37:52 +05:30
if ( constitutive_titanmod_Ntwin ( f , i ) > 0_pInt ) then
2012-02-13 23:11:27 +05:30
if ( constitutive_titanmod_burgersPerTwinFamily ( f , i ) < = 0.0_pReal ) call IO_error ( 231_pInt , e = i , ext_msg = 'twinburgers' )
if ( constitutive_titanmod_twinf0_PerTwinFamily ( f , i ) < = 0.0_pReal ) call IO_error ( 231_pInt , e = i , ext_msg = 'twinf0' )
if ( constitutive_titanmod_twinshearconstant_PerTwinFamily ( f , i ) < = 0.0_pReal ) &
call IO_error ( 231_pInt , e = i , ext_msg = 'twinshearconstant' )
if ( constitutive_titanmod_twintau0_PerTwinFamily ( f , i ) < = 0.0_pReal ) call IO_error ( 231_pInt , e = i , ext_msg = 'twintau0' )
if ( constitutive_titanmod_twingamma0_PerTwinFamily ( f , i ) < = 0.0_pReal ) call IO_error ( 231_pInt , e = i , ext_msg = 'twingamma0' )
2010-07-13 13:49:25 +05:30
endif
enddo
2012-02-13 23:11:27 +05:30
if ( constitutive_titanmod_dc ( i ) < = 0.0_pReal ) call IO_error ( 231_pInt , e = i , ext_msg = 'dc' )
if ( constitutive_titanmod_twinhpconstant ( i ) < = 0.0_pReal ) call IO_error ( 231_pInt , e = i , ext_msg = 'twinhpconstant' )
if ( constitutive_titanmod_aTolRho ( i ) < = 0.0_pReal ) call IO_error ( 231_pInt , e = i , ext_msg = 'aTolRho' )
2010-07-13 13:49:25 +05:30
!* Determine total number of active slip or twin systems
constitutive_titanmod_Nslip ( : , i ) = min ( lattice_NslipSystem ( : , myStructure ) , constitutive_titanmod_Nslip ( : , i ) )
constitutive_titanmod_Ntwin ( : , i ) = min ( lattice_NtwinSystem ( : , myStructure ) , constitutive_titanmod_Ntwin ( : , i ) )
constitutive_titanmod_totalNslip ( i ) = sum ( constitutive_titanmod_Nslip ( : , i ) )
constitutive_titanmod_totalNtwin ( i ) = sum ( constitutive_titanmod_Ntwin ( : , i ) )
write ( 6 , * ) 'Sanity Checks done !'
enddo
!* Allocation of variables whose size depends on the total number of active slip systems
maxTotalNslip = maxval ( constitutive_titanmod_totalNslip )
2012-02-13 23:11:27 +05:30
maxTotalNtwin = maxval ( constitutive_titanmod_totalNtwin )
2010-07-13 13:49:25 +05:30
write ( 6 , * ) 'maxTotalNslip' , maxTotalNslip
write ( 6 , * ) 'maxTotalNtwin' , maxTotalNtwin
allocate ( constitutive_titanmod_burgersPerSlipSystem ( maxTotalNslip , maxNinstance ) )
allocate ( constitutive_titanmod_burgersPerTwinSystem ( maxTotalNtwin , maxNinstance ) )
allocate ( constitutive_titanmod_f0_PerSlipSystem ( maxTotalNslip , maxNinstance ) )
allocate ( constitutive_titanmod_tau0e_PerSlipSystem ( maxTotalNslip , maxNinstance ) )
allocate ( constitutive_titanmod_tau0s_PerSlipSystem ( maxTotalNslip , maxNinstance ) )
allocate ( constitutive_titanmod_capre_PerSlipSystem ( maxTotalNslip , maxNinstance ) )
allocate ( constitutive_titanmod_caprs_PerSlipSystem ( maxTotalNslip , maxNinstance ) )
allocate ( constitutive_titanmod_pe_PerSlipSystem ( maxTotalNslip , maxNinstance ) )
allocate ( constitutive_titanmod_ps_PerSlipSystem ( maxTotalNslip , maxNinstance ) )
allocate ( constitutive_titanmod_qe_PerSlipSystem ( maxTotalNslip , maxNinstance ) )
allocate ( constitutive_titanmod_qs_PerSlipSystem ( maxTotalNslip , maxNinstance ) )
allocate ( constitutive_titanmod_v0e_PerSlipSystem ( maxTotalNslip , maxNinstance ) )
allocate ( constitutive_titanmod_v0s_PerSlipSystem ( maxTotalNslip , maxNinstance ) )
2010-09-22 14:30:40 +05:30
allocate ( constitutive_titanmod_kinkcriticallength_PerSlipSystem ( maxTotalNslip , maxNinstance ) )
allocate ( constitutive_titanmod_CeLambdaSlipPerSlipSystem ( maxTotalNslip , maxNinstance ) )
allocate ( constitutive_titanmod_CsLambdaSlipPerSlipSystem ( maxTotalNslip , maxNinstance ) )
2010-07-13 13:49:25 +05:30
2010-09-01 15:37:52 +05:30
allocate ( constitutive_titanmod_twinf0_PerTwinSystem ( maxTotalNTwin , maxNinstance ) )
2010-09-07 20:14:37 +05:30
allocate ( constitutive_titanmod_twinshearconstant_PerTwinSystem ( maxTotalNTwin , maxNinstance ) )
2010-09-22 14:30:40 +05:30
allocate ( constitutive_titanmod_twintau0_PerTwinSystem ( maxTotalNTwin , maxNinstance ) )
allocate ( constitutive_titanmod_twinp_PerTwinSystem ( maxTotalNTwin , maxNinstance ) )
allocate ( constitutive_titanmod_twinq_PerTwinSystem ( maxTotalNTwin , maxNinstance ) )
allocate ( constitutive_titanmod_twingamma0_PerTwinSystem ( maxTotalNTwin , maxNinstance ) )
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_twinsizePerTwinSystem ( maxTotalNtwin , maxNinstance ) )
2010-09-22 14:30:40 +05:30
allocate ( constitutive_titanmod_twinLambdaSlipPerTwinSystem ( maxTotalNtwin , maxNinstance ) )
2010-09-01 15:37:52 +05:30
2010-07-13 13:49:25 +05:30
constitutive_titanmod_burgersPerSlipSystem = 0.0_pReal
constitutive_titanmod_burgersPerTwinSystem = 0.0_pReal
constitutive_titanmod_f0_PerSlipSystem = 0.0_pReal
constitutive_titanmod_tau0e_PerSlipSystem = 0.0_pReal
constitutive_titanmod_tau0s_PerSlipSystem = 0.0_pReal
constitutive_titanmod_capre_PerSlipSystem = 0.0_pReal
constitutive_titanmod_caprs_PerSlipSystem = 0.0_pReal
constitutive_titanmod_v0e_PerSlipSystem = 0.0_pReal
constitutive_titanmod_v0s_PerSlipSystem = 0.0_pReal
2010-09-22 14:30:40 +05:30
constitutive_titanmod_kinkcriticallength_PerSlipSystem = 0.0_pReal
2010-08-23 17:06:51 +05:30
constitutive_titanmod_pe_PerSlipSystem = 0.0_pReal
constitutive_titanmod_ps_PerSlipSystem = 0.0_pReal
constitutive_titanmod_qe_PerSlipSystem = 0.0_pReal
constitutive_titanmod_qs_PerSlipSystem = 0.0_pReal
2010-07-13 13:49:25 +05:30
2010-09-01 15:37:52 +05:30
constitutive_titanmod_twinf0_PerTwinSystem = 0.0_pReal
2010-09-07 20:14:37 +05:30
constitutive_titanmod_twinshearconstant_PerTwinSystem = 0.0_pReal
2010-09-22 14:30:40 +05:30
constitutive_titanmod_twintau0_PerTwinSystem = 0.0_pReal
constitutive_titanmod_twingamma0_PerTwinSystem = 0.0_pReal
constitutive_titanmod_twinp_PerTwinSystem = 0.0_pReal
constitutive_titanmod_twinq_PerTwinSystem = 0.0_pReal
2010-07-13 13:49:25 +05:30
constitutive_titanmod_twinsizePerTwinSystem = 0.0_pReal
constitutive_titanmod_CeLambdaSlipPerSlipSystem = 0.0_pReal
constitutive_titanmod_CsLambdaSlipPerSlipSystem = 0.0_pReal
2010-09-22 14:30:40 +05:30
constitutive_titanmod_twinLambdaSlipPerTwinSystem = 0.0_pReal
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_interactionMatrixSlipSlip ( maxTotalNslip , maxTotalNslip , maxNinstance ) )
2010-08-23 17:06:51 +05:30
allocate ( constitutive_titanmod_interactionMatrix_ee ( maxTotalNslip , maxTotalNslip , maxNinstance ) )
allocate ( constitutive_titanmod_interactionMatrix_ss ( maxTotalNslip , maxTotalNslip , maxNinstance ) )
allocate ( constitutive_titanmod_interactionMatrix_es ( maxTotalNslip , maxTotalNslip , maxNinstance ) )
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_interactionMatrixSlipTwin ( maxTotalNslip , maxTotalNtwin , maxNinstance ) )
allocate ( constitutive_titanmod_interactionMatrixTwinSlip ( maxTotalNtwin , maxTotalNslip , maxNinstance ) )
allocate ( constitutive_titanmod_interactionMatrixTwinTwin ( maxTotalNtwin , maxTotalNtwin , maxNinstance ) )
allocate ( constitutive_titanmod_forestProjectionEdge ( maxTotalNslip , maxTotalNslip , maxNinstance ) )
2010-08-06 21:23:45 +05:30
allocate ( constitutive_titanmod_forestProjectionScrew ( maxTotalNslip , maxTotalNslip , maxNinstance ) )
2010-09-07 20:14:37 +05:30
allocate ( constitutive_titanmod_TwinforestProjectionEdge ( maxTotalNtwin , maxTotalNtwin , maxNinstance ) )
allocate ( constitutive_titanmod_TwinforestProjectionScrew ( maxTotalNtwin , maxTotalNtwin , maxNinstance ) )
2010-07-13 13:49:25 +05:30
constitutive_titanmod_interactionMatrixSlipSlip = 0.0_pReal
2010-08-23 17:06:51 +05:30
constitutive_titanmod_interactionMatrix_ee = 0.0_pReal
constitutive_titanmod_interactionMatrix_ss = 0.0_pReal
constitutive_titanmod_interactionMatrix_es = 0.0_pReal
2010-07-13 13:49:25 +05:30
constitutive_titanmod_interactionMatrixSlipTwin = 0.0_pReal
constitutive_titanmod_interactionMatrixTwinSlip = 0.0_pReal
constitutive_titanmod_interactionMatrixTwinTwin = 0.0_pReal
constitutive_titanmod_forestProjectionEdge = 0.0_pReal
2010-08-06 21:23:45 +05:30
constitutive_titanmod_forestProjectionScrew = 0.0_pReal
2010-09-01 15:37:52 +05:30
constitutive_titanmod_TwinforestProjectionEdge = 0.0_pReal
constitutive_titanmod_TwinforestProjectionScrew = 0.0_pReal
2010-07-13 13:49:25 +05:30
allocate ( constitutive_titanmod_Ctwin_66 ( 6 , 6 , maxTotalNtwin , maxNinstance ) )
allocate ( constitutive_titanmod_Ctwin_3333 ( 3 , 3 , 3 , 3 , maxTotalNtwin , maxNinstance ) )
constitutive_titanmod_Ctwin_66 = 0.0_pReal
constitutive_titanmod_Ctwin_3333 = 0.0_pReal
2010-09-01 15:37:52 +05:30
2010-07-13 13:49:25 +05:30
write ( 6 , * ) 'Allocated slip system variables'
2012-02-13 19:48:07 +05:30
do i = 1_pInt , maxNinstance
2010-07-13 13:49:25 +05:30
myStructure = constitutive_titanmod_structure ( i )
2010-09-07 20:14:37 +05:30
!* Inverse lookup of slip system family
2010-07-13 13:49:25 +05:30
l = 0_pInt
2012-02-13 19:48:07 +05:30
do f = 1_pInt , lattice_maxNslipFamily
do k = 1_pInt , constitutive_titanmod_Nslip ( f , i )
l = l + 1_pInt
2010-07-13 13:49:25 +05:30
constitutive_titanmod_slipFamily ( l , i ) = f
2012-02-13 19:48:07 +05:30
constitutive_titanmod_slipSystemLattice ( l , i ) = sum ( lattice_NslipSystem ( 1 : f - 1_pInt , myStructure ) ) + k
2010-07-13 13:49:25 +05:30
enddo ; enddo
2010-09-22 14:30:40 +05:30
2010-09-07 20:14:37 +05:30
!* Inverse lookup of twin system family
2010-07-13 13:49:25 +05:30
l = 0_pInt
2012-02-13 19:48:07 +05:30
do f = 1_pInt , lattice_maxNtwinFamily
do k = 1_pInt , constitutive_titanmod_Ntwin ( f , i )
l = l + 1_pInt
2010-07-13 13:49:25 +05:30
constitutive_titanmod_twinFamily ( l , i ) = f
2012-02-13 19:48:07 +05:30
constitutive_titanmod_twinSystemLattice ( l , i ) = sum ( lattice_NtwinSystem ( 1 : f - 1_pInt , myStructure ) ) + k
2010-07-13 13:49:25 +05:30
enddo ; enddo
!* Determine size of state array
ns = constitutive_titanmod_totalNslip ( i )
nt = constitutive_titanmod_totalNtwin ( i )
constitutive_titanmod_sizeDotState ( i ) = &
size ( constitutive_titanmod_listBasicSlipStates ) * ns + size ( constitutive_titanmod_listBasicTwinStates ) * nt
constitutive_titanmod_sizeState ( i ) = &
constitutive_titanmod_sizeDotState ( i ) + &
size ( constitutive_titanmod_listDependentSlipStates ) * ns + size ( constitutive_titanmod_listDependentTwinStates ) * nt
write ( 6 , * ) 'Determined size of state and dot state'
2010-09-29 12:05:08 +05:30
!* Determine size of postResults array
2012-02-21 21:30:00 +05:30
do o = 1_pInt , constitutive_titanmod_Noutput ( i )
2010-07-13 13:49:25 +05:30
select case ( constitutive_titanmod_output ( o , i ) )
case ( 'rhoedge' , &
'rhoscrew' , &
2010-08-06 21:23:45 +05:30
'segment_edge' , &
'segment_screw' , &
2010-09-01 15:37:52 +05:30
'resistance_edge' , &
'resistance_screw' , &
2010-09-07 20:14:37 +05:30
'velocity_edge' , &
'velocity_screw' , &
2011-03-31 14:51:43 +05:30
'tau_slip' , &
2010-09-29 12:05:08 +05:30
'gdot_slip_edge' , &
'gdot_slip_screw' , &
2011-03-31 14:51:43 +05:30
'gdot_slip' , &
2010-09-29 12:05:08 +05:30
'stressratio_edge_p' , &
2011-03-31 14:51:43 +05:30
'stressratio_screw_p' , &
'shear_system' &
2010-07-13 13:49:25 +05:30
)
mySize = constitutive_titanmod_totalNslip ( i )
2010-09-22 14:30:40 +05:30
case ( 'twin_fraction' , &
2010-09-07 20:14:37 +05:30
'gdot_twin' , &
2010-09-22 14:30:40 +05:30
'tau_twin' &
2010-07-13 13:49:25 +05:30
)
mySize = constitutive_titanmod_totalNtwin ( i )
2011-03-31 14:51:43 +05:30
case ( 'shear_basal' , & ! use only if all 4 slip familiies in hex are considered
'shear_prism' , & ! use only if all 4 slip familiies in hex are considered
'shear_pyra' , & ! use only if all 4 slip familiies in hex are considered
'shear_pyrca' , & ! use only if all 4 slip familiies in hex are considered
'rhoedge_basal' , &
'rhoedge_prism' , &
'rhoedge_pyra' , &
'rhoedge_pyrca' , &
'rhoscrew_basal' , &
'rhoscrew_prism' , &
'rhoscrew_pyra' , &
'rhoscrew_pyrca' , &
'shear_total' &
)
mySize = 1_pInt
2010-07-13 13:49:25 +05:30
case default
2012-02-13 23:11:27 +05:30
call IO_error ( 232_pInt , ext_msg = constitutive_titanmod_output ( o , i ) )
2010-07-13 13:49:25 +05:30
end select
if ( mySize > 0_pInt ) then ! any meaningful output found
constitutive_titanmod_sizePostResult ( o , i ) = mySize
constitutive_titanmod_sizePostResults ( i ) = constitutive_titanmod_sizePostResults ( i ) + mySize
endif
enddo
write ( 6 , * ) 'Determining elasticity matrix'
!* Elasticity matrix and shear modulus according to material.config
select case ( myStructure )
2012-02-13 19:48:07 +05:30
case ( 1_pInt : 2_pInt ) ! cubic(s)
forall ( k = 1_pInt : 3_pInt )
forall ( j = 1_pInt : 3_pInt ) &
2010-07-13 13:49:25 +05:30
constitutive_titanmod_Cslip_66 ( k , j , i ) = constitutive_titanmod_C12 ( i )
constitutive_titanmod_Cslip_66 ( k , k , i ) = constitutive_titanmod_C11 ( i )
2012-02-13 19:48:07 +05:30
constitutive_titanmod_Cslip_66 ( k + 3_pInt , k + 3_pInt , i ) = constitutive_titanmod_C44 ( i )
2010-07-13 13:49:25 +05:30
end forall
2012-02-13 19:48:07 +05:30
case ( 3_pInt : ) ! all hex
2010-07-13 13:49:25 +05:30
constitutive_titanmod_Cslip_66 ( 1 , 1 , i ) = constitutive_titanmod_C11 ( i )
constitutive_titanmod_Cslip_66 ( 2 , 2 , i ) = constitutive_titanmod_C11 ( i )
constitutive_titanmod_Cslip_66 ( 3 , 3 , i ) = constitutive_titanmod_C33 ( i )
constitutive_titanmod_Cslip_66 ( 1 , 2 , i ) = constitutive_titanmod_C12 ( i )
constitutive_titanmod_Cslip_66 ( 2 , 1 , i ) = constitutive_titanmod_C12 ( i )
constitutive_titanmod_Cslip_66 ( 1 , 3 , i ) = constitutive_titanmod_C13 ( i )
constitutive_titanmod_Cslip_66 ( 3 , 1 , i ) = constitutive_titanmod_C13 ( i )
constitutive_titanmod_Cslip_66 ( 2 , 3 , i ) = constitutive_titanmod_C13 ( i )
constitutive_titanmod_Cslip_66 ( 3 , 2 , i ) = constitutive_titanmod_C13 ( i )
constitutive_titanmod_Cslip_66 ( 4 , 4 , i ) = constitutive_titanmod_C44 ( i )
constitutive_titanmod_Cslip_66 ( 5 , 5 , i ) = constitutive_titanmod_C44 ( i )
constitutive_titanmod_Cslip_66 ( 6 , 6 , i ) = 0.5_pReal * ( constitutive_titanmod_C11 ( i ) - constitutive_titanmod_C12 ( i ) )
end select
2010-09-22 14:30:40 +05:30
constitutive_titanmod_Cslip_66 ( : , : , i ) = math_Mandel3333to66 ( math_Voigt66to3333 ( constitutive_titanmod_Cslip_66 ( : , : , i ) ) )
2010-07-13 13:49:25 +05:30
constitutive_titanmod_Cslip_3333 ( : , : , : , : , i ) = math_Voigt66to3333 ( constitutive_titanmod_Cslip_66 ( : , : , i ) )
constitutive_titanmod_Gmod ( i ) = &
0.2_pReal * ( constitutive_titanmod_C11 ( i ) - constitutive_titanmod_C12 ( i ) ) + 0.3_pReal * constitutive_titanmod_C44 ( i )
!* Construction of the twin elasticity matrices
2012-02-13 19:48:07 +05:30
do j = 1_pInt , lattice_maxNtwinFamily
do k = 1_pInt , constitutive_titanmod_Ntwin ( j , i )
do l = 1_pInt , 3_pInt ; do m = 1_pInt , 3_pInt ; do n = 1_pInt , 3_pInt ; do o = 1_pInt , 3_pInt
do p = 1_pInt , 3_pInt ; do q = 1_pInt , 3_pInt ; do r = 1_pInt , 3_pInt ; do s = 1_pInt , 3_pInt
2012-02-21 21:30:00 +05:30
constitutive_titanmod_Ctwin_3333 ( l , m , n , o , sum ( constitutive_titanmod_Nslip ( 1 : j - 1_pInt , i ) ) + k , i ) = &
constitutive_titanmod_Ctwin_3333 ( l , m , n , o , sum ( constitutive_titanmod_Nslip ( 1 : j - 1_pInt , i ) ) + k , i ) + &
2010-07-13 13:49:25 +05:30
constitutive_titanmod_Cslip_3333 ( p , q , r , s , i ) * &
2012-02-21 21:30:00 +05:30
lattice_Qtwin ( l , p , sum ( lattice_NslipSystem ( 1 : j - 1_pInt , myStructure ) ) + k , myStructure ) * &
lattice_Qtwin ( m , q , sum ( lattice_NslipSystem ( 1 : j - 1_pInt , myStructure ) ) + k , myStructure ) * &
lattice_Qtwin ( n , r , sum ( lattice_NslipSystem ( 1 : j - 1_pInt , myStructure ) ) + k , myStructure ) * &
lattice_Qtwin ( o , s , sum ( lattice_NslipSystem ( 1 : j - 1_pInt , myStructure ) ) + k , myStructure )
2010-07-13 13:49:25 +05:30
enddo ; enddo ; enddo ; enddo ; enddo ; enddo ; enddo ; enddo
constitutive_titanmod_Ctwin_66 ( : , : , k , i ) = math_Mandel3333to66 ( constitutive_titanmod_Ctwin_3333 ( : , : , : , : , k , i ) )
enddo
enddo
2010-09-07 20:14:37 +05:30
!* Burgers vector, dislocation velocity prefactor for each slip system
2012-02-21 21:30:00 +05:30
do s = 1_pInt , constitutive_titanmod_totalNslip ( i )
2010-07-13 13:49:25 +05:30
f = constitutive_titanmod_slipFamily ( s , i )
constitutive_titanmod_burgersPerSlipSystem ( s , i ) = constitutive_titanmod_burgersPerSlipFamily ( f , i )
constitutive_titanmod_f0_PerSlipSystem ( s , i ) = constitutive_titanmod_f0_PerSlipFamily ( f , i )
constitutive_titanmod_tau0e_PerSlipSystem ( s , i ) = constitutive_titanmod_tau0e_PerSlipFamily ( f , i )
constitutive_titanmod_tau0s_PerSlipSystem ( s , i ) = constitutive_titanmod_tau0s_PerSlipFamily ( f , i )
constitutive_titanmod_capre_PerSlipSystem ( s , i ) = constitutive_titanmod_capre_PerSlipFamily ( f , i )
constitutive_titanmod_caprs_PerSlipSystem ( s , i ) = constitutive_titanmod_caprs_PerSlipFamily ( f , i )
constitutive_titanmod_v0e_PerSlipSystem ( s , i ) = constitutive_titanmod_v0e_PerSlipFamily ( f , i )
constitutive_titanmod_v0s_PerSlipSystem ( s , i ) = constitutive_titanmod_v0s_PerSlipFamily ( f , i )
2010-09-22 14:30:40 +05:30
constitutive_titanmod_kinkcriticallength_PerSlipSystem ( s , i ) = constitutive_titanmod_kinkcriticallength_PerSlipFamily ( f , i )
2010-07-13 13:49:25 +05:30
constitutive_titanmod_pe_PerSlipSystem ( s , i ) = constitutive_titanmod_pe_PerSlipFamily ( f , i )
constitutive_titanmod_ps_PerSlipSystem ( s , i ) = constitutive_titanmod_ps_PerSlipFamily ( f , i )
constitutive_titanmod_qe_PerSlipSystem ( s , i ) = constitutive_titanmod_qe_PerSlipFamily ( f , i )
constitutive_titanmod_qs_PerSlipSystem ( s , i ) = constitutive_titanmod_qs_PerSlipFamily ( f , i )
constitutive_titanmod_CeLambdaSlipPerSlipSystem ( s , i ) = constitutive_titanmod_CeLambdaSlipPerSlipFamily ( f , i )
constitutive_titanmod_CsLambdaSlipPerSlipSystem ( s , i ) = constitutive_titanmod_CsLambdaSlipPerSlipFamily ( f , i )
enddo
!* Burgers vector, nucleation rate prefactor and twin size for each twin system
2012-02-21 21:30:00 +05:30
do t = 1_pInt , constitutive_titanmod_totalNtwin ( i )
2010-09-22 14:30:40 +05:30
f = constitutive_titanmod_twinFamily ( t , i )
constitutive_titanmod_burgersPerTwinSystem ( t , i ) = constitutive_titanmod_burgersPerTwinFamily ( f , i )
constitutive_titanmod_twinsizePerTwinSystem ( t , i ) = constitutive_titanmod_twinsizePerTwinFamily ( f , i )
constitutive_titanmod_twinf0_PerTwinSystem ( t , i ) = constitutive_titanmod_twinf0_PerTwinFamily ( f , i )
constitutive_titanmod_twinshearconstant_PerTwinSystem ( t , i ) = constitutive_titanmod_twinshearconstant_PerTwinFamily ( f , i )
constitutive_titanmod_twintau0_PerTwinSystem ( t , i ) = constitutive_titanmod_twintau0_PerTwinFamily ( f , i )
constitutive_titanmod_twingamma0_PerTwinSystem ( t , i ) = constitutive_titanmod_twingamma0_PerTwinFamily ( f , i )
constitutive_titanmod_twinp_PerTwinSystem ( t , i ) = constitutive_titanmod_twinp_PerTwinFamily ( f , i )
constitutive_titanmod_twinq_PerTwinSystem ( t , i ) = constitutive_titanmod_twinq_PerTwinFamily ( f , i )
constitutive_titanmod_twinLambdaSlipPerTwinSystem ( t , i ) = constitutive_titanmod_twinLambdaSlipPerTwinFamily ( f , i )
2010-09-01 15:37:52 +05:30
enddo
2010-07-13 13:49:25 +05:30
!* Construction of interaction matrices
2012-02-21 21:30:00 +05:30
do s1 = 1_pInt , constitutive_titanmod_totalNslip ( i )
do s2 = 1_pInt , constitutive_titanmod_totalNslip ( i )
2010-07-13 13:49:25 +05:30
constitutive_titanmod_interactionMatrixSlipSlip ( s1 , s2 , i ) = &
constitutive_titanmod_interactionSlipSlip ( lattice_interactionSlipSlip ( constitutive_titanmod_slipSystemLattice ( s1 , i ) , &
constitutive_titanmod_slipSystemLattice ( s2 , i ) , &
myStructure ) , i )
2010-08-23 17:06:51 +05:30
constitutive_titanmod_interactionMatrix_ee ( s1 , s2 , i ) = &
constitutive_titanmod_interaction_ee ( lattice_interactionSlipSlip ( constitutive_titanmod_slipSystemLattice ( s1 , i ) , &
constitutive_titanmod_slipSystemLattice ( s2 , i ) , &
myStructure ) , i )
constitutive_titanmod_interactionMatrix_ss ( s1 , s2 , i ) = &
constitutive_titanmod_interaction_ss ( lattice_interactionSlipSlip ( constitutive_titanmod_slipSystemLattice ( s1 , i ) , &
constitutive_titanmod_slipSystemLattice ( s2 , i ) , &
myStructure ) , i )
constitutive_titanmod_interactionMatrix_es ( s1 , s2 , i ) = &
constitutive_titanmod_interaction_es ( lattice_interactionSlipSlip ( constitutive_titanmod_slipSystemLattice ( s1 , i ) , &
constitutive_titanmod_slipSystemLattice ( s2 , i ) , &
myStructure ) , i )
enddo ; enddo
2010-07-13 13:49:25 +05:30
2012-02-21 21:30:00 +05:30
do s1 = 1_pInt , constitutive_titanmod_totalNslip ( i )
do t2 = 1_pInt , constitutive_titanmod_totalNtwin ( i )
2010-07-13 13:49:25 +05:30
constitutive_titanmod_interactionMatrixSlipTwin ( s1 , t2 , i ) = &
constitutive_titanmod_interactionSlipTwin ( lattice_interactionSlipTwin ( constitutive_titanmod_slipSystemLattice ( s1 , i ) , &
constitutive_titanmod_twinSystemLattice ( t2 , i ) , &
myStructure ) , i )
enddo ; enddo
2012-02-21 21:30:00 +05:30
do t1 = 1_pInt , constitutive_titanmod_totalNtwin ( i )
do s2 = 1_pInt , constitutive_titanmod_totalNslip ( i )
2010-07-13 13:49:25 +05:30
constitutive_titanmod_interactionMatrixTwinSlip ( t1 , s2 , i ) = &
constitutive_titanmod_interactionTwinSlip ( lattice_interactionTwinSlip ( constitutive_titanmod_twinSystemLattice ( t1 , i ) , &
constitutive_titanmod_slipSystemLattice ( s2 , i ) , &
myStructure ) , i )
enddo ; enddo
2012-02-21 21:30:00 +05:30
do t1 = 1_pInt , constitutive_titanmod_totalNtwin ( i )
do t2 = 1_pInt , constitutive_titanmod_totalNtwin ( i )
2010-07-13 13:49:25 +05:30
constitutive_titanmod_interactionMatrixTwinTwin ( t1 , t2 , i ) = &
constitutive_titanmod_interactionTwinTwin ( lattice_interactionTwinTwin ( constitutive_titanmod_twinSystemLattice ( t1 , i ) , &
constitutive_titanmod_twinSystemLattice ( t2 , i ) , &
myStructure ) , i )
enddo ; enddo
!* Calculation of forest projections for edge dislocations
2012-02-21 21:30:00 +05:30
do s1 = 1_pInt , constitutive_titanmod_totalNslip ( i )
do s2 = 1_pInt , constitutive_titanmod_totalNslip ( i )
2010-07-13 13:49:25 +05:30
constitutive_titanmod_forestProjectionEdge ( s1 , s2 , i ) = &
abs ( math_mul3x3 ( lattice_sn ( : , constitutive_titanmod_slipSystemLattice ( s1 , i ) , myStructure ) , &
lattice_st ( : , constitutive_titanmod_slipSystemLattice ( s2 , i ) , myStructure ) ) )
2010-08-06 21:23:45 +05:30
!* Calculation of forest projections for screw dislocations
constitutive_titanmod_forestProjectionScrew ( s1 , s2 , i ) = &
abs ( math_mul3x3 ( lattice_sn ( : , constitutive_titanmod_slipSystemLattice ( s1 , i ) , myStructure ) , &
lattice_sd ( : , constitutive_titanmod_slipSystemLattice ( s2 , i ) , myStructure ) ) )
enddo ; enddo
2010-07-13 13:49:25 +05:30
2010-09-01 15:37:52 +05:30
!* Calculation of forest projections for edge dislocations in twin system
2012-02-21 21:30:00 +05:30
do t1 = 1_pInt , constitutive_titanmod_totalNtwin ( i )
do t2 = 1_pInt , constitutive_titanmod_totalNtwin ( i )
2010-09-01 15:37:52 +05:30
constitutive_titanmod_TwinforestProjectionEdge ( t1 , t2 , i ) = &
abs ( math_mul3x3 ( lattice_tn ( : , constitutive_titanmod_twinSystemLattice ( t1 , i ) , myStructure ) , &
lattice_tt ( : , constitutive_titanmod_twinSystemLattice ( t2 , i ) , myStructure ) ) )
!* Calculation of forest projections for screw dislocations in twin system
constitutive_titanmod_TwinforestProjectionScrew ( t1 , t2 , i ) = &
abs ( math_mul3x3 ( lattice_tn ( : , constitutive_titanmod_twinSystemLattice ( t1 , i ) , myStructure ) , &
lattice_td ( : , constitutive_titanmod_twinSystemLattice ( t2 , i ) , myStructure ) ) )
enddo ; enddo
enddo
2010-07-13 13:49:25 +05:30
write ( 6 , * ) 'Init All done'
return
end subroutine
function constitutive_titanmod_stateInit ( myInstance )
!*********************************************************************
!* initial microstructural state *
!*********************************************************************
use prec , only : pReal , pInt
use lattice , only : lattice_maxNslipFamily , lattice_maxNtwinFamily
implicit none
!* Input-Output variables
integer ( pInt ) :: myInstance
real ( pReal ) , dimension ( constitutive_titanmod_sizeState ( myInstance ) ) :: constitutive_titanmod_stateInit
!* Local variables
2010-09-01 15:37:52 +05:30
integer ( pInt ) s0 , s1 , s , t , f , ns , nt , ts0 , ts1 , tf , ts
2010-08-06 21:23:45 +05:30
real ( pReal ) , dimension ( constitutive_titanmod_totalNslip ( myInstance ) ) :: rho_edge0 , &
2010-07-13 13:49:25 +05:30
rho_screw0 , &
2011-03-31 14:51:43 +05:30
shear_system0 , &
2010-08-06 21:23:45 +05:30
segment_edge0 , &
segment_screw0 , &
resistance_edge0 , &
2010-08-23 17:06:51 +05:30
resistance_screw0
2010-09-22 14:30:40 +05:30
real ( pReal ) , dimension ( constitutive_titanmod_totalNtwin ( myInstance ) ) :: twingamma_dot0 , &
resistance_twin0
2010-07-13 13:49:25 +05:30
ns = constitutive_titanmod_totalNslip ( myInstance )
nt = constitutive_titanmod_totalNtwin ( myInstance )
constitutive_titanmod_stateInit = 0.0_pReal
!* Initialize basic slip state variables
2010-09-01 15:37:52 +05:30
! For slip
2010-07-13 13:49:25 +05:30
s1 = 0_pInt
2012-02-21 21:30:00 +05:30
do f = 1_pInt , lattice_maxNslipFamily
2010-07-13 13:49:25 +05:30
s0 = s1 + 1_pInt
s1 = s0 + constitutive_titanmod_Nslip ( f , myInstance ) - 1_pInt
do s = s0 , s1
rho_edge0 ( s ) = constitutive_titanmod_rho_edge0 ( f , myInstance )
rho_screw0 ( s ) = constitutive_titanmod_rho_screw0 ( f , myInstance )
2011-03-31 14:51:43 +05:30
shear_system0 ( s ) = 0.0_pReal
2010-07-13 13:49:25 +05:30
enddo
enddo
2010-09-01 15:37:52 +05:30
!* Initialize basic slip state variables
! For twin
ts1 = 0_pInt
2012-02-21 21:30:00 +05:30
do tf = 1_pInt , lattice_maxNtwinFamily
2010-09-01 15:37:52 +05:30
ts0 = ts1 + 1_pInt
ts1 = ts0 + constitutive_titanmod_Ntwin ( tf , myInstance ) - 1_pInt
do ts = ts0 , ts1
2010-09-07 20:14:37 +05:30
twingamma_dot0 ( ts ) = 0.0_pReal
2010-09-01 15:37:52 +05:30
enddo
enddo
2010-07-13 13:49:25 +05:30
constitutive_titanmod_stateInit ( 1 : ns ) = rho_edge0
2012-02-21 21:30:00 +05:30
constitutive_titanmod_stateInit ( ns + 1_pInt : 2_pInt * ns ) = rho_screw0
constitutive_titanmod_stateInit ( 2_pInt * ns + 1_pInt : 3_pInt * ns ) = shear_system0
constitutive_titanmod_stateInit ( 3_pInt * ns + 1_pInt : 3_pInt * ns + nt ) = twingamma_dot0
2010-07-13 13:49:25 +05:30
!* Initialize dependent slip microstructural variables
2012-02-21 21:30:00 +05:30
forall ( s = 1_pInt : ns ) &
2010-08-06 21:23:45 +05:30
segment_edge0 ( s ) = constitutive_titanmod_CeLambdaSlipPerSlipSystem ( s , myInstance ) / &
2010-09-22 14:30:40 +05:30
sqrt ( dot_product ( ( rho_edge0 ) , constitutive_titanmod_forestProjectionEdge ( 1 : ns , s , myInstance ) ) + &
dot_product ( ( rho_screw0 ) , constitutive_titanmod_forestProjectionScrew ( 1 : ns , s , myInstance ) ) )
2010-08-06 21:23:45 +05:30
2012-02-21 21:30:00 +05:30
constitutive_titanmod_stateInit ( 3_pInt * ns + nt + 1_pInt : 4_pInt * ns + nt ) = segment_edge0
2010-07-13 13:49:25 +05:30
2012-02-21 21:30:00 +05:30
forall ( s = 1_pInt : ns ) &
2010-08-06 21:23:45 +05:30
segment_screw0 ( s ) = constitutive_titanmod_CsLambdaSlipPerSlipSystem ( s , myInstance ) / &
2010-09-22 14:30:40 +05:30
sqrt ( dot_product ( ( rho_edge0 ) , constitutive_titanmod_forestProjectionEdge ( 1 : ns , s , myInstance ) ) + &
dot_product ( ( rho_screw0 ) , constitutive_titanmod_forestProjectionScrew ( 1 : ns , s , myInstance ) ) )
2010-08-06 21:23:45 +05:30
2012-02-21 21:30:00 +05:30
constitutive_titanmod_stateInit ( 4_pInt * ns + nt + 1_pInt : 5_pInt * ns + nt ) = segment_screw0
2010-07-13 13:49:25 +05:30
2012-02-21 21:30:00 +05:30
forall ( s = 1_pInt : ns ) &
2010-08-06 21:23:45 +05:30
resistance_edge0 ( s ) = &
2010-07-13 13:49:25 +05:30
constitutive_titanmod_Gmod ( myInstance ) * constitutive_titanmod_burgersPerSlipSystem ( s , myInstance ) * &
2010-08-23 17:06:51 +05:30
sqrt ( dot_product ( ( rho_edge0 ) , constitutive_titanmod_interactionMatrix_ee ( 1 : ns , s , myInstance ) ) + dot_product ( ( rho_screw0 ) , &
constitutive_titanmod_interactionMatrix_es ( 1 : ns , s , myInstance ) ) )
2011-03-31 14:51:43 +05:30
2012-02-21 21:30:00 +05:30
constitutive_titanmod_stateInit ( 5_pInt * ns + nt + 1_pInt : 6_pInt * ns + nt ) = resistance_edge0
2010-07-13 13:49:25 +05:30
2012-02-21 21:30:00 +05:30
forall ( s = 1_pInt : ns ) &
2010-08-06 21:23:45 +05:30
resistance_screw0 ( s ) = &
2010-07-13 13:49:25 +05:30
constitutive_titanmod_Gmod ( myInstance ) * constitutive_titanmod_burgersPerSlipSystem ( s , myInstance ) * &
2010-08-23 17:06:51 +05:30
sqrt ( dot_product ( ( rho_edge0 ) , constitutive_titanmod_interactionMatrix_es ( 1 : ns , s , myInstance ) ) + dot_product ( ( rho_screw0 ) , &
constitutive_titanmod_interactionMatrix_ss ( 1 : ns , s , myInstance ) ) )
2011-03-31 14:51:43 +05:30
2012-02-21 21:30:00 +05:30
constitutive_titanmod_stateInit ( 6_pInt * ns + nt + 1_pInt : 7_pInt * ns + nt ) = resistance_screw0
2010-07-13 13:49:25 +05:30
2012-02-21 21:30:00 +05:30
forall ( t = 1_pInt : nt ) &
2010-09-22 14:30:40 +05:30
resistance_twin0 ( t ) = 0.0_pReal
2012-02-21 21:30:00 +05:30
constitutive_titanmod_stateInit ( 7_pInt * ns + nt + 1_pInt : 7_pInt * ns + 2_pInt * nt ) = resistance_twin0
2010-07-13 13:49:25 +05:30
return
end function
2010-10-26 18:46:37 +05:30
pure function constitutive_titanmod_aTolState ( myInstance )
2010-07-13 13:49:25 +05:30
!*********************************************************************
2010-10-26 18:46:37 +05:30
!* absolute state tolerance *
2010-07-13 13:49:25 +05:30
!*********************************************************************
use prec , only : pReal , pInt
implicit none
!* Input-Output variables
integer ( pInt ) , intent ( in ) :: myInstance
2010-10-26 18:46:37 +05:30
real ( pReal ) , dimension ( constitutive_titanmod_sizeState ( myInstance ) ) :: constitutive_titanmod_aTolState
2010-07-13 13:49:25 +05:30
2010-10-26 18:46:37 +05:30
constitutive_titanmod_aTolState = constitutive_titanmod_aTolRho ( myInstance )
2010-07-13 13:49:25 +05:30
return
endfunction
pure function constitutive_titanmod_homogenizedC ( state , g , ip , el )
!*********************************************************************
!* calculates homogenized elacticity matrix *
!* - state : microstructure quantities *
!* - g : component-ID of current integration point *
!* - ip : current integration point *
!* - el : current element *
!*********************************************************************
use prec , only : pReal , pInt , p_vec
use mesh , only : mesh_NcpElems , mesh_maxNips
use material , only : homogenization_maxNgrains , material_phase , phase_constitutionInstance
implicit none
!* Input-Output variables
integer ( pInt ) , intent ( in ) :: g , ip , el
type ( p_vec ) , dimension ( homogenization_maxNgrains , mesh_maxNips , mesh_NcpElems ) , intent ( in ) :: state
real ( pReal ) , dimension ( 6 , 6 ) :: constitutive_titanmod_homogenizedC
2010-09-07 20:14:37 +05:30
real ( pReal ) , dimension ( constitutive_titanmod_totalNtwin ( phase_constitutionInstance ( material_phase ( g , ip , el ) ) ) ) :: &
volumefraction_pertwinsystem
2010-07-13 13:49:25 +05:30
!* Local variables
integer ( pInt ) myInstance , ns , nt , i
real ( pReal ) sumf
!* Shortened notation
myInstance = phase_constitutionInstance ( material_phase ( g , ip , el ) )
ns = constitutive_titanmod_totalNslip ( myInstance )
nt = constitutive_titanmod_totalNtwin ( myInstance )
!* Total twin volume fraction
2012-02-21 21:30:00 +05:30
do i = 1_pInt , nt
volumefraction_pertwinsystem ( i ) = state ( g , ip , el ) % p ( 3_pInt * ns + i ) / &
2010-09-07 20:14:37 +05:30
constitutive_titanmod_twinshearconstant_PerTwinSystem ( i , myInstance )
enddo
!sumf = sum(state(g,ip,el)%p((6*ns+7*nt+1):(6*ns+8*nt))) ! safe for nt == 0
sumf = sum ( abs ( volumefraction_pertwinsystem ( 1 : nt ) ) ) ! safe for nt == 0
2010-07-13 13:49:25 +05:30
!* Homogenized elasticity matrix
constitutive_titanmod_homogenizedC = ( 1.0_pReal - sumf ) * constitutive_titanmod_Cslip_66 ( : , : , myInstance )
2012-02-21 21:30:00 +05:30
do i = 1_pInt , nt
2010-07-13 13:49:25 +05:30
constitutive_titanmod_homogenizedC = &
2010-09-07 20:14:37 +05:30
! constitutive_titanmod_homogenizedC + state(g,ip,el)%p(6*ns+7*nt+i)*constitutive_titanmod_Ctwin_66(:,:,i,myInstance)
constitutive_titanmod_homogenizedC + volumefraction_pertwinsystem ( i ) * constitutive_titanmod_Ctwin_66 ( : , : , i , myInstance )
2010-07-13 13:49:25 +05:30
enddo
return
end function
subroutine constitutive_titanmod_microstructure ( Temperature , state , g , ip , el )
!*********************************************************************
!* calculates quantities characterizing the microstructure *
!* - Temperature : temperature *
!* - state : microstructure quantities *
!* - ipc : component-ID of current integration point *
!* - ip : current integration point *
!* - el : current element *
!*********************************************************************
use prec , only : pReal , pInt , p_vec
use mesh , only : mesh_NcpElems , mesh_maxNips
use material , only : homogenization_maxNgrains , material_phase , phase_constitutionInstance
!use debug, only: debugger
implicit none
!* Input-Output variables
integer ( pInt ) , intent ( in ) :: g , ip , el
real ( pReal ) , intent ( in ) :: Temperature
type ( p_vec ) , dimension ( homogenization_maxNgrains , mesh_maxNips , mesh_NcpElems ) , intent ( inout ) :: state
!* Local variables
integer ( pInt ) myInstance , myStructure , ns , nt , s , t , i
real ( pReal ) sumf , sfe
2010-09-07 20:14:37 +05:30
real ( pReal ) , dimension ( constitutive_titanmod_totalNtwin ( phase_constitutionInstance ( material_phase ( g , ip , el ) ) ) ) :: &
2011-04-13 19:46:22 +05:30
volumefraction_pertwinsystem
2010-07-13 13:49:25 +05:30
!* Shortened notation
myInstance = phase_constitutionInstance ( material_phase ( g , ip , el ) )
myStructure = constitutive_titanmod_structure ( myInstance )
ns = constitutive_titanmod_totalNslip ( myInstance )
nt = constitutive_titanmod_totalNtwin ( myInstance )
2010-09-01 15:37:52 +05:30
! Need to update this list
2010-07-13 13:49:25 +05:30
!* State: 1 : ns rho_edge
!* State: ns+1 : 2*ns rho_screw
2011-03-31 14:51:43 +05:30
!* State: 2*ns+1 : 3*ns shear_system
!* State: 3*ns+1 : 3*ns+nt gamma_twin
!* State: 3*ns+nt+1 : 4*ns+nt segment_edge
!* State: 4*ns+nt+1 : 5*ns+nt segment_screw
!* State: 5*ns+nt+1 : 6*ns+nt resistance_edge
!* State: 6*ns+nt+1 : 7*ns+nt resistance_screw
!* State: 7*ns+nt+1 : 7*ns+2*nt resistance_twin
!* State: 7*ns+2*nt+1 : 8*ns+2*nt velocity_edge
!* State: 8*ns+2*nt+1 : 9*ns+2*nt velocity_screw
!* State: 9*ns+2*nt+1 : 10*ns+2*nt tau_slip
!* State: 10*ns+2*nt+1 : 11*ns+2*nt gdot_slip_edge
!* State: 11*ns+2*nt+1 : 12*ns+2*nt gdot_slip_screw
!* State: 12*ns+2*nt+1 : 13*ns+2*nt StressRatio_edge_p
!* State: 13*ns+2*nt+1 : 14*ns+2*nt StressRatio_screw_p
2010-07-13 13:49:25 +05:30
!* Total twin volume fraction
2012-02-21 21:30:00 +05:30
do i = 1_pInt , nt
volumefraction_pertwinsystem ( i ) = state ( g , ip , el ) % p ( 3_pInt * ns + i ) / &
2010-09-07 20:14:37 +05:30
constitutive_titanmod_twinshearconstant_PerTwinSystem ( i , myInstance )
enddo
!sumf = sum(state(g,ip,el)%p((6*ns+7*nt+1):(6*ns+8*nt))) ! safe for nt == 0
sumf = sum ( abs ( volumefraction_pertwinsystem ( 1 : nt ) ) ) ! safe for nt == 0
2010-07-13 13:49:25 +05:30
!* Stacking fault energy
sfe = 0.0002_pReal * Temperature - 0.0396_pReal
!* rescaled twin volume fraction for topology
2010-09-01 15:37:52 +05:30
!forall (t = 1:nt) &
! fOverStacksize(t) = &
! state(g,ip,el)%p(2*ns+t)/constitutive_titanmod_twinsizePerTwinSystem(t,myInstance)
2010-08-23 17:06:51 +05:30
2010-09-01 15:37:52 +05:30
! average segment length for edge dislocations in matrix
2012-02-21 21:30:00 +05:30
forall ( s = 1_pInt : ns ) &
state ( g , ip , el ) % p ( 3_pInt * ns + nt + s ) = constitutive_titanmod_CeLambdaSlipPerSlipSystem ( s , myInstance ) / &
2010-09-22 14:30:40 +05:30
sqrt ( dot_product ( state ( g , ip , el ) % p ( 1 : ns ) , &
constitutive_titanmod_forestProjectionEdge ( 1 : ns , s , myInstance ) ) + &
2012-02-21 21:30:00 +05:30
dot_product ( state ( g , ip , el ) % p ( ns + 1_pInt : 2_pInt * ns ) , &
2010-09-22 14:30:40 +05:30
constitutive_titanmod_forestProjectionScrew ( 1 : ns , s , myInstance ) ) )
! average segment length for screw dislocations in matrix
2012-02-21 21:30:00 +05:30
forall ( s = 1_pInt : ns ) &
state ( g , ip , el ) % p ( 4_pInt * ns + nt + s ) = constitutive_titanmod_CsLambdaSlipPerSlipSystem ( s , myInstance ) / &
2010-09-22 14:30:40 +05:30
sqrt ( dot_product ( state ( g , ip , el ) % p ( 1 : ns ) , &
constitutive_titanmod_forestProjectionEdge ( 1 : ns , s , myInstance ) ) + &
2012-02-21 21:30:00 +05:30
dot_product ( state ( g , ip , el ) % p ( ns + 1_pInt : 2_pInt * ns ) , &
2010-09-22 14:30:40 +05:30
constitutive_titanmod_forestProjectionScrew ( 1 : ns , s , myInstance ) ) )
2010-07-13 13:49:25 +05:30
2010-08-06 21:23:45 +05:30
!* threshold stress or slip resistance for edge dislocation motion
2012-02-21 21:30:00 +05:30
forall ( s = 1_pInt : ns ) &
state ( g , ip , el ) % p ( 5_pInt * ns + nt + s ) = &
2010-07-13 13:49:25 +05:30
constitutive_titanmod_Gmod ( myInstance ) * constitutive_titanmod_burgersPerSlipSystem ( s , myInstance ) * &
2010-08-23 17:06:51 +05:30
sqrt ( dot_product ( ( state ( g , ip , el ) % p ( 1 : ns ) ) , &
constitutive_titanmod_interactionMatrix_ee ( 1 : ns , s , myInstance ) ) + &
2012-02-21 21:30:00 +05:30
dot_product ( ( state ( g , ip , el ) % p ( ns + 1_pInt : 2_pInt * ns ) ) , &
2010-08-23 17:06:51 +05:30
constitutive_titanmod_interactionMatrix_es ( 1 : ns , s , myInstance ) ) )
2010-07-13 13:49:25 +05:30
2010-08-06 21:23:45 +05:30
!* threshold stress or slip resistance for screw dislocation motion
2012-02-21 21:30:00 +05:30
forall ( s = 1_pInt : ns ) &
state ( g , ip , el ) % p ( 6_pInt * ns + nt + s ) = &
2010-07-13 13:49:25 +05:30
constitutive_titanmod_Gmod ( myInstance ) * constitutive_titanmod_burgersPerSlipSystem ( s , myInstance ) * &
2010-08-23 17:06:51 +05:30
sqrt ( dot_product ( ( state ( g , ip , el ) % p ( 1 : ns ) ) , &
constitutive_titanmod_interactionMatrix_es ( 1 : ns , s , myInstance ) ) + &
2012-02-21 21:30:00 +05:30
dot_product ( ( state ( g , ip , el ) % p ( ns + 1_pInt : 2_pInt * ns ) ) , &
2010-08-23 17:06:51 +05:30
constitutive_titanmod_interactionMatrix_ss ( 1 : ns , s , myInstance ) ) )
2010-09-01 15:37:52 +05:30
2010-09-22 14:30:40 +05:30
!* threshold stress or slip resistance for dislocation motion in twin
2012-02-21 21:30:00 +05:30
forall ( t = 1_pInt : nt ) &
state ( g , ip , el ) % p ( 7_pInt * ns + nt + t ) = &
2010-09-01 15:37:52 +05:30
constitutive_titanmod_Gmod ( myInstance ) * constitutive_titanmod_burgersPerTwinSystem ( t , myInstance ) * &
2012-02-21 21:30:00 +05:30
( dot_product ( ( abs ( state ( g , ip , el ) % p ( 2_pInt * ns + 1_pInt : 2_pInt * ns + nt ) ) ) , &
2010-09-01 15:37:52 +05:30
constitutive_titanmod_interactionMatrixTwinTwin ( 1 : nt , t , myInstance ) ) )
2010-07-13 13:49:25 +05:30
return
end subroutine
subroutine constitutive_titanmod_LpAndItsTangent ( Lp , dLp_dTstar , Tstar_v , Temperature , state , g , ip , el )
!*********************************************************************
!* calculates plastic velocity gradient and its tangent *
!* INPUT: *
!* - Temperature : temperature *
!* - state : microstructure quantities *
!* - Tstar_v : 2nd Piola Kirchhoff stress tensor (Mandel) *
!* - ipc : component-ID at current integration point *
!* - ip : current integration point *
!* - el : current element *
!* OUTPUT: *
!* - Lp : plastic velocity gradient *
!* - dLp_dTstar : derivative of Lp (4th-rank tensor) *
!*********************************************************************
use prec , only : pReal , pInt , p_vec
use math , only : math_Plain3333to99
use mesh , only : mesh_NcpElems , mesh_maxNips
use material , only : homogenization_maxNgrains , material_phase , phase_constitutionInstance
2012-02-21 21:30:00 +05:30
use lattice , only : lattice_Sslip , lattice_Sslip_v , lattice_Stwin_v , lattice_maxNslipFamily , lattice_maxNtwinFamily , &
lattice_NslipSystem , lattice_NtwinSystem , lattice_Stwin
2010-07-13 13:49:25 +05:30
implicit none
!* Input-Output variables
integer ( pInt ) , intent ( in ) :: g , ip , el
real ( pReal ) , intent ( in ) :: Temperature
real ( pReal ) , dimension ( 6 ) , intent ( in ) :: Tstar_v
type ( p_vec ) , dimension ( homogenization_maxNgrains , mesh_maxNips , mesh_NcpElems ) , intent ( inout ) :: state
real ( pReal ) , dimension ( 3 , 3 ) , intent ( out ) :: Lp
real ( pReal ) , dimension ( 9 , 9 ) , intent ( out ) :: dLp_dTstar
!* Local variables
integer ( pInt ) myInstance , myStructure , ns , nt , f , i , j , k , l , m , n , index_myFamily
2010-08-06 21:23:45 +05:30
real ( pReal ) sumf , StressRatio_edge_p , minusStressRatio_edge_p , StressRatio_edge_pminus1 , StressRatio_screw_p , &
2011-04-13 19:46:22 +05:30
StressRatio_screw_pminus1 , BoltzmannRatioedge , minusStressRatio_screw_p , &
2010-09-22 14:30:40 +05:30
screwvelocity_prefactor , twinStressRatio_p , twinminusStressRatio_p , twinStressRatio_pminus1 , &
2011-04-13 19:46:22 +05:30
twinDotGamma0 , BoltzmannRatioscrew , BoltzmannRatiotwin , bottomstress_edge , bottomstress_screw
2010-07-13 13:49:25 +05:30
real ( pReal ) , dimension ( 3 , 3 , 3 , 3 ) :: dLp_dTstar3333
real ( pReal ) , dimension ( constitutive_titanmod_totalNslip ( phase_constitutionInstance ( material_phase ( g , ip , el ) ) ) ) :: &
2010-09-29 12:05:08 +05:30
gdot_slip , dgdot_dtauslip , tau_slip , edge_velocity , screw_velocity , gdot_slip_edge , gdot_slip_screw
2010-07-13 13:49:25 +05:30
real ( pReal ) , dimension ( constitutive_titanmod_totalNtwin ( phase_constitutionInstance ( material_phase ( g , ip , el ) ) ) ) :: &
2011-04-13 19:46:22 +05:30
gdot_twin , dgdot_dtautwin , tau_twin , volumefraction_pertwinsystem
2010-07-13 13:49:25 +05:30
!* Shortened notation
myInstance = phase_constitutionInstance ( material_phase ( g , ip , el ) )
myStructure = constitutive_titanmod_structure ( myInstance )
ns = constitutive_titanmod_totalNslip ( myInstance )
nt = constitutive_titanmod_totalNtwin ( myInstance )
2012-02-21 21:30:00 +05:30
do i = 1_pInt , nt
volumefraction_pertwinsystem ( i ) = state ( g , ip , el ) % p ( 3_pInt * ns + i ) / &
2010-09-07 20:14:37 +05:30
constitutive_titanmod_twinshearconstant_PerTwinSystem ( i , myInstance )
enddo
sumf = sum ( abs ( volumefraction_pertwinsystem ( 1 : nt ) ) ) ! safe for nt == 0
2010-07-13 13:49:25 +05:30
Lp = 0.0_pReal
dLp_dTstar3333 = 0.0_pReal
dLp_dTstar = 0.0_pReal
!* Dislocation glide part
gdot_slip = 0.0_pReal
2010-09-29 12:05:08 +05:30
gdot_slip_edge = 0.0_pReal
gdot_slip_screw = 0.0_pReal
2010-07-13 13:49:25 +05:30
dgdot_dtauslip = 0.0_pReal
j = 0_pInt
2012-02-21 21:30:00 +05:30
do f = 1_pInt , lattice_maxNslipFamily ! loop over all slip families
index_myFamily = sum ( lattice_NslipSystem ( 1 : f - 1_pInt , myStructure ) ) ! at which index starts my family
do i = 1_pInt , constitutive_titanmod_Nslip ( f , myInstance ) ! process each (active) slip system in family
2010-07-13 13:49:25 +05:30
j = j + 1_pInt
!* Calculation of Lp
!* Resolved shear stress on slip system
tau_slip ( j ) = dot_product ( Tstar_v , lattice_Sslip_v ( : , index_myFamily + i , myStructure ) )
2010-09-22 14:30:40 +05:30
!*************************************************
2011-03-31 14:51:43 +05:30
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! if(myStructure>=3.and.j>3) then ! for all non-basal slip systems
2012-02-21 21:30:00 +05:30
if ( myStructure == 3_pInt ) then ! only for prismatic and pyr <a> systems in hex
2010-09-29 12:05:08 +05:30
screwvelocity_prefactor = constitutive_titanmod_debyefrequency ( myInstance ) * &
2012-02-21 21:30:00 +05:30
state ( g , ip , el ) % p ( 4_pInt * ns + nt + j ) * ( constitutive_titanmod_burgersPerSlipSystem ( j , myInstance ) / &
2010-09-22 14:30:40 +05:30
constitutive_titanmod_kinkcriticallength_PerSlipSystem ( j , myInstance ) ) ** 2
2010-09-29 12:05:08 +05:30
!* Stress ratio for screw ! No slip resistance for screw dislocations, only Peierls stress
2011-03-31 14:51:43 +05:30
bottomstress_screw = constitutive_titanmod_tau0s_PerSlipSystem ( j , myInstance )
2010-08-23 17:06:51 +05:30
StressRatio_screw_p = ( ( abs ( tau_slip ( j ) ) ) / &
2011-03-31 14:51:43 +05:30
( bottomstress_screw ) &
2010-08-06 21:23:45 +05:30
) ** constitutive_titanmod_ps_PerSlipSystem ( j , myInstance )
2011-03-31 14:51:43 +05:30
2010-09-01 15:37:52 +05:30
if ( ( 1.0_pReal - StressRatio_screw_p ) > 0.001_pReal ) then
2010-08-23 17:06:51 +05:30
minusStressRatio_screw_p = 1.0_pReal - StressRatio_screw_p
else
2010-09-01 15:37:52 +05:30
minusStressRatio_screw_p = 0.001_pReal
2010-08-23 17:06:51 +05:30
endif
2011-03-31 14:51:43 +05:30
bottomstress_screw = constitutive_titanmod_tau0s_PerSlipSystem ( j , myInstance )
2010-08-06 21:23:45 +05:30
StressRatio_screw_pminus1 = ( ( abs ( tau_slip ( j ) ) ) / &
2011-03-31 14:51:43 +05:30
( bottomstress_screw ) &
2010-08-06 21:23:45 +05:30
) ** ( constitutive_titanmod_ps_PerSlipSystem ( j , myInstance ) - 1.0_pReal )
2010-07-13 13:49:25 +05:30
2010-09-29 12:05:08 +05:30
!* Boltzmann ratio for screw
BoltzmannRatioscrew = constitutive_titanmod_kinkf0 ( myInstance ) / ( kB * Temperature )
2010-07-13 13:49:25 +05:30
2011-03-31 14:51:43 +05:30
else ! if the structure is not hex or the slip family is basal
screwvelocity_prefactor = constitutive_titanmod_v0s_PerSlipSystem ( j , myInstance )
bottomstress_screw = constitutive_titanmod_tau0s_PerSlipSystem ( j , myInstance ) + state ( g , ip , el ) % p ( 6 * ns + nt + j )
StressRatio_screw_p = ( ( abs ( tau_slip ( j ) ) ) / ( bottomstress_screw ) ) ** constitutive_titanmod_ps_PerSlipSystem ( j , myInstance )
if ( ( 1.0_pReal - StressRatio_screw_p ) > 0.001_pReal ) then
minusStressRatio_screw_p = 1.0_pReal - StressRatio_screw_p
else
minusStressRatio_screw_p = 0.001_pReal
endif
StressRatio_screw_pminus1 = ( ( abs ( tau_slip ( j ) ) ) / ( bottomstress_screw ) ) ** &
( constitutive_titanmod_ps_PerSlipSystem ( j , myInstance ) - 1.0_pReal )
!* Boltzmann ratio for screw
BoltzmannRatioscrew = constitutive_titanmod_f0_PerSlipSystem ( j , myInstance ) / ( kB * Temperature )
endif
!* Stress ratio for edge
bottomstress_edge = constitutive_titanmod_tau0e_PerSlipSystem ( j , myInstance ) + state ( g , ip , el ) % p ( 5 * ns + nt + j )
StressRatio_edge_p = ( ( abs ( tau_slip ( j ) ) ) / &
( bottomstress_edge ) &
) ** constitutive_titanmod_pe_PerSlipSystem ( j , myInstance )
if ( ( 1.0_pReal - StressRatio_edge_p ) > 0.001_pReal ) then
minusStressRatio_edge_p = 1.0_pReal - StressRatio_edge_p
else
minusStressRatio_edge_p = 0.001_pReal
endif
StressRatio_edge_pminus1 = ( ( abs ( tau_slip ( j ) ) ) / ( bottomstress_edge ) ) ** &
( constitutive_titanmod_pe_PerSlipSystem ( j , myInstance ) - 1.0_pReal )
!* Boltzmann ratio for edge. For screws it is defined above
BoltzmannRatioedge = constitutive_titanmod_f0_PerSlipSystem ( j , myInstance ) / ( kB * Temperature )
2010-09-22 14:30:40 +05:30
screw_velocity ( j ) = screwvelocity_prefactor * & ! there is no v0 for screw now because it is included in the prefactor
2010-09-29 12:05:08 +05:30
exp ( - BoltzmannRatioscrew * ( minusStressRatio_screw_p ) ** &
2010-08-06 21:23:45 +05:30
constitutive_titanmod_qs_PerSlipSystem ( j , myInstance ) )
2011-03-31 14:51:43 +05:30
edge_velocity ( j ) = constitutive_titanmod_v0e_PerSlipSystem ( j , myInstance ) * exp ( - BoltzmannRatioedge * &
( minusStressRatio_edge_p ) ** &
constitutive_titanmod_qe_PerSlipSystem ( j , myInstance ) )
2010-09-29 12:05:08 +05:30
!* Shear rates due to edge slip
gdot_slip_edge ( j ) = constitutive_titanmod_burgersPerSlipSystem ( j , myInstance ) * ( state ( g , ip , el ) % p ( j ) * &
edge_velocity ( j ) ) * sign ( 1.0_pReal , tau_slip ( j ) )
!* Shear rates due to screw slip
gdot_slip_screw ( j ) = constitutive_titanmod_burgersPerSlipSystem ( j , myInstance ) * ( state ( g , ip , el ) % p ( ns + j ) * &
screw_velocity ( j ) ) * sign ( 1.0_pReal , tau_slip ( j ) )
!Total shear rate
gdot_slip ( j ) = gdot_slip_edge ( j ) + gdot_slip_screw ( j )
2010-08-23 17:06:51 +05:30
2011-03-31 14:51:43 +05:30
state ( g , ip , el ) % p ( 7 * ns + 2 * nt + j ) = edge_velocity ( j )
state ( g , ip , el ) % p ( 8 * ns + 2 * nt + j ) = screw_velocity ( j )
state ( g , ip , el ) % p ( 9 * ns + 2 * nt + j ) = tau_slip ( j )
state ( g , ip , el ) % p ( 10 * ns + 2 * nt + j ) = gdot_slip_edge ( j )
state ( g , ip , el ) % p ( 11 * ns + 2 * nt + j ) = gdot_slip_screw ( j )
state ( g , ip , el ) % p ( 12 * ns + 2 * nt + j ) = StressRatio_edge_p
state ( g , ip , el ) % p ( 13 * ns + 2 * nt + j ) = StressRatio_screw_p
2010-09-22 14:30:40 +05:30
2010-07-13 13:49:25 +05:30
!* Derivatives of shear rates
2010-09-22 14:30:40 +05:30
dgdot_dtauslip ( j ) = constitutive_titanmod_burgersPerSlipSystem ( j , myInstance ) * ( ( &
2010-08-23 17:06:51 +05:30
( &
( &
( &
2010-09-22 14:30:40 +05:30
( edge_velocity ( j ) * state ( g , ip , el ) % p ( j ) ) ) * &
2011-03-31 14:51:43 +05:30
BoltzmannRatioedge * &
2010-07-15 12:46:15 +05:30
constitutive_titanmod_pe_PerSlipSystem ( j , myInstance ) * &
2010-08-23 17:06:51 +05:30
constitutive_titanmod_qe_PerSlipSystem ( j , myInstance ) &
) / &
2011-03-31 14:51:43 +05:30
bottomstress_edge &
2010-08-23 17:06:51 +05:30
) * &
2010-08-06 21:23:45 +05:30
StressRatio_edge_pminus1 * ( minusStressRatio_edge_p ) ** &
2010-08-23 17:06:51 +05:30
( constitutive_titanmod_qe_PerSlipSystem ( j , myInstance ) - 1.0_pReal ) &
) + &
( &
( &
( &
2010-09-22 14:30:40 +05:30
( state ( g , ip , el ) % p ( ns + j ) * screw_velocity ( j ) ) * &
2011-03-31 14:51:43 +05:30
BoltzmannRatioscrew * &
2010-07-15 12:46:15 +05:30
constitutive_titanmod_ps_PerSlipSystem ( j , myInstance ) * &
2010-08-23 17:06:51 +05:30
constitutive_titanmod_qs_PerSlipSystem ( j , myInstance ) &
) / &
2011-03-31 14:51:43 +05:30
bottomstress_screw &
2010-08-23 17:06:51 +05:30
) * &
2010-08-06 21:23:45 +05:30
StressRatio_screw_pminus1 * ( minusStressRatio_screw_p ) ** ( constitutive_titanmod_qs_PerSlipSystem ( j , myInstance ) - 1.0_pReal ) &
2010-08-23 17:06:51 +05:30
) &
) !* sign(1.0_pReal,tau_slip(j))
2011-03-31 14:51:43 +05:30
2010-08-23 17:06:51 +05:30
!*************************************************
2010-09-07 20:14:37 +05:30
!sumf=0.0_pReal
2010-07-13 13:49:25 +05:30
!* Plastic velocity gradient for dislocation glide
Lp = Lp + ( 1.0_pReal - sumf ) * gdot_slip ( j ) * lattice_Sslip ( : , : , index_myFamily + i , myStructure )
!* Calculation of the tangent of Lp
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 ) &
2010-07-13 13:49:25 +05:30
dLp_dTstar3333 ( k , l , m , n ) = &
dLp_dTstar3333 ( k , l , m , n ) + dgdot_dtauslip ( j ) * &
lattice_Sslip ( k , l , index_myFamily + i , myStructure ) * &
lattice_Sslip ( m , n , index_myFamily + i , myStructure )
enddo
enddo
!* Mechanical twinning part
gdot_twin = 0.0_pReal
dgdot_dtautwin = 0.0_pReal
j = 0_pInt
2012-02-21 21:30:00 +05:30
do f = 1_pInt , lattice_maxNtwinFamily ! loop over all slip families
index_myFamily = sum ( lattice_NtwinSystem ( 1 : f - 1_pInt , myStructure ) ) ! at which index starts my family
do i = 1_pInt , constitutive_titanmod_Ntwin ( f , myInstance ) ! process each (active) slip system in family
2010-07-13 13:49:25 +05:30
j = j + 1_pInt
!* Calculation of Lp
!* Resolved shear stress on twin system
tau_twin ( j ) = dot_product ( Tstar_v , lattice_Stwin_v ( : , index_myFamily + i , myStructure ) )
2010-09-01 15:37:52 +05:30
!**************************************************************************************
!* Stress ratios
! StressRatio_r = (state(g,ip,el)%p(6*ns+3*nt+j)/tau_twin(j))**constitutive_titanmod_r(myInstance)
2010-07-13 13:49:25 +05:30
2010-08-23 17:06:51 +05:30
!* Shear rates and their derivatives due to twin
2010-08-06 21:23:45 +05:30
! if ( tau_twin(j) > 0.0_pReal ) !then
2010-09-01 15:37:52 +05:30
! gdot_twin(j) = 0.0_pReal!&
2010-08-06 21:23:45 +05:30
! (constitutive_titanmod_MaxTwinFraction(myInstance)-sumf)*lattice_shearTwin(index_myFamily+i,myStructure)*&
! state(g,ip,el)%p(6*ns+4*nt+j)*constitutive_titanmod_Ndot0PerTwinSystem(f,myInstance)*exp(-StressRatio_r)
! dgdot_dtautwin(j) = ((gdot_twin(j)*constitutive_titanmod_r(myInstance))/tau_twin(j))*StressRatio_r
! endif
2010-09-01 15:37:52 +05:30
!**************************************************************************************
!* Stress ratio for edge
2010-09-22 14:30:40 +05:30
twinStressRatio_p = ( ( abs ( tau_twin ( j ) ) ) / &
2011-03-31 14:51:43 +05:30
( constitutive_titanmod_twintau0_PerTwinSystem ( j , myInstance ) + state ( g , ip , el ) % p ( 7 * ns + nt + j ) ) &
2010-09-22 14:30:40 +05:30
) ** constitutive_titanmod_twinp_PerTwinSystem ( j , myInstance )
2011-03-31 14:51:43 +05:30
2010-09-22 14:30:40 +05:30
if ( ( 1.0_pReal - twinStressRatio_p ) > 0.001_pReal ) then
twinminusStressRatio_p = 1.0_pReal - twinStressRatio_p
2010-09-01 15:37:52 +05:30
else
2010-09-22 14:30:40 +05:30
twinminusStressRatio_p = 0.001_pReal
2010-09-01 15:37:52 +05:30
endif
2011-03-31 14:51:43 +05:30
2010-09-22 14:30:40 +05:30
twinStressRatio_pminus1 = ( ( abs ( tau_twin ( j ) ) ) / &
2011-03-31 14:51:43 +05:30
( constitutive_titanmod_twintau0_PerTwinSystem ( j , myInstance ) + state ( g , ip , el ) % p ( 7 * ns + nt + j ) ) &
2010-09-22 14:30:40 +05:30
) ** ( constitutive_titanmod_twinp_PerTwinSystem ( j , myInstance ) - 1.0_pReal )
2010-09-01 15:37:52 +05:30
!* Boltzmann ratio
2011-03-31 14:51:43 +05:30
BoltzmannRatiotwin = constitutive_titanmod_twinf0_PerTwinSystem ( j , myInstance ) / ( kB * Temperature )
2010-09-01 15:37:52 +05:30
2010-09-22 14:30:40 +05:30
!* Initial twin shear rates
2010-09-01 15:37:52 +05:30
TwinDotGamma0 = &
2010-09-22 14:30:40 +05:30
constitutive_titanmod_twingamma0_PerTwinSystem ( j , myInstance )
!* Shear rates due to twin
gdot_twin ( j ) = sign ( 1.0_pReal , tau_twin ( j ) ) * constitutive_titanmod_twingamma0_PerTwinSystem ( j , myInstance ) * &
2011-03-31 14:51:43 +05:30
exp ( - BoltzmannRatiotwin * ( twinminusStressRatio_p ) ** constitutive_titanmod_twinq_PerTwinSystem ( j , myInstance ) )
2010-09-22 14:30:40 +05:30
2011-03-31 14:51:43 +05:30
2010-09-01 15:37:52 +05:30
!* Derivatives of shear rates in twin
dgdot_dtautwin ( j ) = ( &
( &
( &
( abs ( gdot_twin ( j ) ) ) * &
2011-03-31 14:51:43 +05:30
BoltzmannRatiotwin * &
2010-09-22 14:30:40 +05:30
constitutive_titanmod_twinp_PerTwinSystem ( j , myInstance ) * &
constitutive_titanmod_twinq_PerTwinSystem ( j , myInstance ) &
2010-09-01 15:37:52 +05:30
) / &
2010-09-22 14:30:40 +05:30
constitutive_titanmod_twintau0_PerTwinSystem ( j , myInstance ) &
2010-09-01 15:37:52 +05:30
) * &
2010-09-22 14:30:40 +05:30
twinStressRatio_pminus1 * ( twinminusStressRatio_p ) ** &
( constitutive_titanmod_twinq_PerTwinSystem ( j , myInstance ) - 1.0_pReal ) &
2010-09-01 15:37:52 +05:30
) !* sign(1.0_pReal,tau_slip(j))
2010-08-23 17:06:51 +05:30
!* Plastic velocity gradient for mechanical twinning
2010-09-13 14:59:03 +05:30
! Lp = Lp + sumf*gdot_twin(j)*lattice_Stwin(:,:,index_myFamily+i,myStructure)
Lp = Lp + gdot_twin ( j ) * lattice_Stwin ( : , : , index_myFamily + i , myStructure )
2010-07-13 13:49:25 +05:30
!* Calculation of the tangent of Lp
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 ) &
2010-07-13 13:49:25 +05:30
dLp_dTstar3333 ( k , l , m , n ) = &
dLp_dTstar3333 ( k , l , m , n ) + dgdot_dtautwin ( j ) * &
lattice_Stwin ( k , l , index_myFamily + i , myStructure ) * &
lattice_Stwin ( m , n , index_myFamily + i , myStructure )
enddo
enddo
dLp_dTstar = math_Plain3333to99 ( dLp_dTstar3333 )
!if ((ip==1).and.(el==1)) then
! write(6,*) '#LP/TANGENT#'
! write(6,*)
! write(6,*) 'Tstar_v', Tstar_v
! write(6,*) 'tau_slip', tau_slip
2012-02-02 01:50:05 +05:30
! write(6,'(a10,/,4(3(e20.8,1x),/))') 'state',state(1,1,1)%p
! write(6,'(a,/,3(3(f10.4,1x)/))') 'Lp',Lp
! write(6,'(a,/,9(9(f10.4,1x)/))') 'dLp_dTstar',dLp_dTstar
2010-07-13 13:49:25 +05:30
!endif
return
end subroutine
function constitutive_titanmod_dotState ( Tstar_v , Temperature , state , g , ip , el )
!*********************************************************************
!* rate of change of microstructure *
!* INPUT: *
!* - Temperature : temperature *
!* - state : microstructure quantities *
!* - Tstar_v : 2nd Piola Kirchhoff stress tensor (Mandel) *
!* - ipc : component-ID at current integration point *
!* - ip : current integration point *
!* - el : current element *
!* OUTPUT: *
!* - constitutive_dotState : evolution of state variable *
!*********************************************************************
use prec , only : pReal , pInt , p_vec
use mesh , only : mesh_NcpElems , mesh_maxNips
use material , only : homogenization_maxNgrains , material_phase , phase_constitutionInstance
2012-02-21 21:30:00 +05:30
use lattice , only : lattice_maxNslipFamily , lattice_maxNtwinFamily , &
lattice_NslipSystem , lattice_NtwinSystem , lattice_Stwin_v
2010-07-13 13:49:25 +05:30
implicit none
!* Input-Output variables
integer ( pInt ) , intent ( in ) :: g , ip , el
real ( pReal ) , intent ( in ) :: Temperature
real ( pReal ) , dimension ( 6 ) , intent ( in ) :: Tstar_v
type ( p_vec ) , dimension ( homogenization_maxNgrains , mesh_maxNips , mesh_NcpElems ) , intent ( in ) :: state
real ( pReal ) , dimension ( constitutive_titanmod_sizeDotState ( phase_constitutionInstance ( material_phase ( g , ip , el ) ) ) ) :: &
constitutive_titanmod_dotState
!* Local variables
2011-04-13 19:46:22 +05:30
integer ( pInt ) MyInstance , MyStructure , ns , nt , f , i , j , index_myFamily
real ( pReal ) sumf , BoltzmannRatio , &
twinStressRatio_p , twinminusStressRatio_p
2010-07-13 13:49:25 +05:30
real ( pReal ) , dimension ( constitutive_titanmod_totalNslip ( phase_constitutionInstance ( material_phase ( g , ip , el ) ) ) ) :: &
2011-04-13 19:46:22 +05:30
DotRhoEdgeGeneration , DotRhoEdgeAnnihilation , DotRhoScrewAnnihilation , &
DotRhoScrewGeneration
2010-08-23 17:06:51 +05:30
real ( pReal ) , dimension ( constitutive_titanmod_totalNtwin ( phase_constitutionInstance ( material_phase ( g , ip , el ) ) ) ) :: gdot_twin , &
2011-04-13 19:46:22 +05:30
tau_twin , &
volumefraction_pertwinsystem
2010-07-13 13:49:25 +05:30
!* Shortened notation
myInstance = phase_constitutionInstance ( material_phase ( g , ip , el ) )
MyStructure = constitutive_titanmod_structure ( myInstance )
ns = constitutive_titanmod_totalNslip ( myInstance )
nt = constitutive_titanmod_totalNtwin ( myInstance )
2012-02-21 21:30:00 +05:30
do i = 1_pInt , nt
volumefraction_pertwinsystem ( i ) = state ( g , ip , el ) % p ( 3_pInt * ns + i ) / &
2010-09-07 20:14:37 +05:30
constitutive_titanmod_twinshearconstant_PerTwinSystem ( i , myInstance )
enddo
2012-02-21 21:30:00 +05:30
sumf = sum ( abs ( volumefraction_pertwinsystem ( 1_pInt : nt ) ) ) ! safe for nt == 0
2010-07-13 13:49:25 +05:30
constitutive_titanmod_dotState = 0.0_pReal
2010-09-01 15:37:52 +05:30
j = 0_pInt
2012-02-21 21:30:00 +05:30
do f = 1_pInt , lattice_maxNslipFamily ! loop over all slip families
index_myFamily = sum ( lattice_NslipSystem ( 1 : f - 1_pInt , myStructure ) ) ! at which index starts my family
do i = 1_pInt , constitutive_titanmod_Nslip ( f , myInstance ) ! process each (active) slip system in family
2010-07-13 13:49:25 +05:30
j = j + 1_pInt
!* Multiplication of edge dislocations
2011-03-31 14:51:43 +05:30
DotRhoEdgeGeneration ( j ) = ( state ( g , ip , el ) % p ( ns + j ) * state ( g , ip , el ) % p ( 8 * ns + 2 * nt + j ) / state ( g , ip , el ) % p ( 4 * ns + nt + j ) )
2010-07-13 13:49:25 +05:30
!* Multiplication of screw dislocations
2011-03-31 14:51:43 +05:30
DotRhoScrewGeneration ( j ) = ( state ( g , ip , el ) % p ( j ) * state ( g , ip , el ) % p ( 7 * ns + 2 * nt + j ) / state ( g , ip , el ) % p ( 3 * ns + nt + j ) )
2010-07-13 13:49:25 +05:30
!* Annihilation of edge dislocations
2010-08-06 21:23:45 +05:30
DotRhoEdgeAnnihilation ( j ) = - ( ( state ( g , ip , el ) % p ( j ) ) ** 2 ) * &
2011-03-31 14:51:43 +05:30
constitutive_titanmod_capre_PerSlipSystem ( j , myInstance ) * state ( g , ip , el ) % p ( 7 * ns + 2 * nt + j ) / 2.0_pReal
2010-07-13 13:49:25 +05:30
!* Annihilation of screw dislocations
2010-08-06 21:23:45 +05:30
DotRhoScrewAnnihilation ( j ) = - ( ( state ( g , ip , el ) % p ( ns + j ) ) ** 2 ) * &
2011-03-31 14:51:43 +05:30
constitutive_titanmod_caprs_PerSlipSystem ( j , myInstance ) * state ( g , ip , el ) % p ( 8 * ns + 2 * nt + j ) / 2.0_pReal
2010-07-13 13:49:25 +05:30
!* Edge dislocation density rate of change
constitutive_titanmod_dotState ( j ) = &
DotRhoEdgeGeneration ( j ) + DotRhoEdgeAnnihilation ( j )
!* Screw dislocation density rate of change
constitutive_titanmod_dotState ( ns + j ) = &
DotRhoScrewGeneration ( j ) + DotRhoScrewAnnihilation ( j )
2010-07-15 12:46:15 +05:30
2011-03-31 14:51:43 +05:30
constitutive_titanmod_dotState ( 2 * ns + j ) = &
state ( g , ip , el ) % p ( 10 * ns + 2 * nt + j ) + state ( g , ip , el ) % p ( 11 * ns + 2 * nt + j ) ! sum of shear due to edge and screw
2010-07-13 13:49:25 +05:30
enddo
enddo
2010-08-06 21:23:45 +05:30
2010-09-01 15:37:52 +05:30
!* Twin fraction evolution
2010-07-13 13:49:25 +05:30
j = 0_pInt
2012-02-21 21:30:00 +05:30
do f = 1_pInt , lattice_maxNtwinFamily ! loop over all twin families
index_myFamily = sum ( lattice_NtwinSystem ( 1 : f - 1_pInt , MyStructure ) ) ! at which index starts my family
do i = 1_pInt , constitutive_titanmod_Ntwin ( f , myInstance ) ! process each (active) twin system in family
2010-07-13 13:49:25 +05:30
j = j + 1_pInt
2010-09-01 15:37:52 +05:30
!*************************************************************************
!This was in dislotwin - keeping it for safety
!*************************************************************************
! !* Resolved shear stress on twin system
! tau_twin(j) = dot_product(Tstar_v,lattice_Stwin_v(:,index_myFamily+i,myStructure))
! !* Stress ratios
! StressRatio_r = (state(g,ip,el)%p(6*ns+3*nt+j)/tau_twin(j))**constitutive_titanmod_r(myInstance)
!
! !* Shear rates and their derivatives due to twin
! if ( tau_twin(j) > 0.0_pReal ) then
! constitutive_titanmod_dotState(2*ns+j) = &
! (constitutive_titanmod_MaxTwinFraction(myInstance)-sumf)*&
! state(g,ip,el)%p(6*ns+4*nt+j)*constitutive_titanmod_Ndot0PerTwinSystem(f,myInstance)*exp(-StressRatio_r)
! endif
!*************************************************************************
2010-07-13 13:49:25 +05:30
!* Resolved shear stress on twin system
tau_twin ( j ) = dot_product ( Tstar_v , lattice_Stwin_v ( : , index_myFamily + i , myStructure ) )
2010-09-01 15:37:52 +05:30
!* Stress ratio for edge
2010-09-22 14:30:40 +05:30
twinStressRatio_p = ( ( abs ( tau_twin ( j ) ) ) / &
2011-03-31 14:51:43 +05:30
( constitutive_titanmod_twintau0_PerTwinSystem ( j , myInstance ) + state ( g , ip , el ) % p ( 7 * ns + nt + j ) ) &
2010-09-22 14:30:40 +05:30
) ** ( constitutive_titanmod_twinp_PerTwinSystem ( j , myInstance ) )
2010-09-01 15:37:52 +05:30
2010-09-22 14:30:40 +05:30
if ( ( 1.0_pReal - twinStressRatio_p ) > 0.001_pReal ) then
twinminusStressRatio_p = 1.0_pReal - twinStressRatio_p
2010-09-01 15:37:52 +05:30
else
2010-09-22 14:30:40 +05:30
twinminusStressRatio_p = 0.001_pReal
2010-09-01 15:37:52 +05:30
endif
!* Boltzmann ratio
BoltzmannRatio = constitutive_titanmod_twinf0_PerTwinSystem ( j , myInstance ) / ( kB * Temperature )
2010-09-29 12:05:08 +05:30
2010-09-22 14:30:40 +05:30
gdot_twin ( j ) = constitutive_titanmod_twingamma0_PerTwinSystem ( j , myInstance ) * exp ( - BoltzmannRatio * &
( twinminusStressRatio_p ) ** &
constitutive_titanmod_twinq_PerTwinSystem ( j , myInstance ) ) * sign ( 1.0_pReal , tau_twin ( j ) )
2011-03-31 14:51:43 +05:30
constitutive_titanmod_dotState ( 3 * ns + j ) = gdot_twin ( j )
2010-09-07 20:14:37 +05:30
enddo
2010-07-13 13:49:25 +05:30
enddo
!write(6,*) '#DOTSTATE#'
!write(6,*)
2012-02-02 01:50:05 +05:30
!write(6,'(a,/,4(3(f30.20,1x)/))') 'EdgeGeneration',DotRhoEdgeGeneration
!write(6,'(a,/,4(3(f30.20,1x)/))') 'ScrewGeneration',DotRhoScrewGeneration
!write(6,'(a,/,4(3(f30.20,1x)/))') 'EdgeAnnihilation',DotRhoEdgeAnnihilation
!write(6,'(a,/,4(3(f30.20,1x)/))') 'ScrewAnnihilation',DotRhoScrewAnnihilation
2010-09-29 12:05:08 +05:30
2010-07-13 13:49:25 +05:30
return
end function
pure function constitutive_titanmod_dotTemperature ( Tstar_v , Temperature , state , g , ip , el )
!*********************************************************************
!* rate of change of microstructure *
!* INPUT: *
!* - Temperature : temperature *
!* - Tstar_v : 2nd Piola Kirchhoff stress tensor (Mandel) *
!* - ipc : component-ID at current integration point *
!* - ip : current integration point *
!* - el : current element *
!* OUTPUT: *
!* - constitutive_dotTemperature : evolution of Temperature *
!*********************************************************************
use prec , only : pReal , pInt , p_vec
use mesh , only : mesh_NcpElems , mesh_maxNips
use material , only : homogenization_maxNgrains
implicit none
!* Input-Output variables
integer ( pInt ) , intent ( in ) :: g , ip , el
real ( pReal ) , intent ( in ) :: Temperature
real ( pReal ) , dimension ( 6 ) , intent ( in ) :: Tstar_v
type ( p_vec ) , dimension ( homogenization_maxNgrains , mesh_maxNips , mesh_NcpElems ) , intent ( in ) :: state
real ( pReal ) constitutive_titanmod_dotTemperature
constitutive_titanmod_dotTemperature = 0.0_pReal
return
end function
pure function constitutive_titanmod_postResults ( Tstar_v , Temperature , dt , state , g , ip , el )
!*********************************************************************
!* return array of constitutive results *
!* INPUT: *
!* - Temperature : temperature *
!* - Tstar_v : 2nd Piola Kirchhoff stress tensor (Mandel) *
!* - dt : current time increment *
!* - ipc : component-ID at current integration point *
!* - ip : current integration point *
!* - el : current element *
!*********************************************************************
use prec , only : pReal , pInt , p_vec
use mesh , only : mesh_NcpElems , mesh_maxNips
use material , only : homogenization_maxNgrains , material_phase , phase_constitutionInstance , phase_Noutput
implicit none
!* Definition of variables
integer ( pInt ) , intent ( in ) :: g , ip , el
real ( pReal ) , intent ( in ) :: dt , Temperature
real ( pReal ) , dimension ( 6 ) , intent ( in ) :: Tstar_v
type ( p_vec ) , dimension ( homogenization_maxNgrains , mesh_maxNips , mesh_NcpElems ) , intent ( in ) :: state
2011-04-13 19:46:22 +05:30
integer ( pInt ) myInstance , myStructure , ns , nt , o , i , c
2010-09-29 12:05:08 +05:30
real ( pReal ) sumf
2010-07-13 13:49:25 +05:30
real ( pReal ) , dimension ( constitutive_titanmod_sizePostResults ( phase_constitutionInstance ( material_phase ( g , ip , el ) ) ) ) :: &
constitutive_titanmod_postResults
2010-09-07 20:14:37 +05:30
real ( pReal ) , dimension ( constitutive_titanmod_totalNtwin ( phase_constitutionInstance ( material_phase ( g , ip , el ) ) ) ) :: &
volumefraction_pertwinsystem
2010-07-13 13:49:25 +05:30
!* Shortened notation
myInstance = phase_constitutionInstance ( material_phase ( g , ip , el ) )
myStructure = constitutive_titanmod_structure ( myInstance )
ns = constitutive_titanmod_totalNslip ( myInstance )
nt = constitutive_titanmod_totalNtwin ( myInstance )
2012-02-21 21:30:00 +05:30
do i = 1_pInt , nt
volumefraction_pertwinsystem ( i ) = state ( g , ip , el ) % p ( 3_pInt * ns + i ) / &
2010-09-07 20:14:37 +05:30
constitutive_titanmod_twinshearconstant_PerTwinSystem ( i , myInstance )
enddo
!sumf = sum(state(g,ip,el)%p((6*ns+7*nt+1):(6*ns+8*nt))) ! safe for nt == 0
sumf = sum ( abs ( volumefraction_pertwinsystem ( 1 : nt ) ) ) ! safe for nt == 0
2010-07-13 13:49:25 +05:30
!* Required output
c = 0_pInt
constitutive_titanmod_postResults = 0.0_pReal
2012-02-21 21:30:00 +05:30
do o = 1_pInt , phase_Noutput ( material_phase ( g , ip , el ) )
2010-07-13 13:49:25 +05:30
select case ( constitutive_titanmod_output ( o , myInstance ) )
case ( 'rhoedge' )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_postResults ( c + 1_pInt : c + ns ) = state ( g , ip , el ) % p ( 1_pInt : ns )
2010-07-13 13:49:25 +05:30
c = c + ns
case ( 'rhoscrew' )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_postResults ( c + 1_pInt : c + ns ) = state ( g , ip , el ) % p ( ns + 1_pInt : 2_pInt * ns )
2010-07-13 13:49:25 +05:30
c = c + ns
2010-08-06 21:23:45 +05:30
case ( 'segment_edge' )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_postResults ( c + 1_pInt : c + ns ) = state ( g , ip , el ) % p ( ( 3_pInt * ns + nt + 1_pInt ) : ( 4_pInt * ns + nt ) )
2010-08-06 21:23:45 +05:30
c = c + ns
2011-03-31 14:51:43 +05:30
case ( 'segment_screw' )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_postResults ( c + 1_pInt : c + ns ) = state ( g , ip , el ) % p ( ( 4_pInt * ns + nt + 1_pInt ) : ( 5_pInt * ns + nt ) )
2010-08-06 21:23:45 +05:30
c = c + ns
2011-03-31 14:51:43 +05:30
case ( 'resistance_edge' )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_postResults ( c + 1_pInt : c + ns ) = state ( g , ip , el ) % p ( ( 5_pInt * ns + nt + 1_pInt ) : ( 6_pInt * ns + nt ) )
2010-08-06 21:23:45 +05:30
c = c + ns
2011-03-31 14:51:43 +05:30
case ( 'resistance_screw' )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_postResults ( c + 1_pInt : c + ns ) = state ( g , ip , el ) % p ( ( 6_pInt * ns + nt + 1_pInt ) : ( 7_pInt * ns + nt ) )
2010-09-29 12:05:08 +05:30
c = c + ns
case ( 'velocity_edge' )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_postResults ( c + 1_pInt : c + ns ) = state ( g , ip , el ) % p ( ( 7 * ns + 2 * nt + 1 ) : ( 8 * ns + 2 * nt ) )
2010-09-29 12:05:08 +05:30
c = c + ns
case ( 'velocity_screw' )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_postResults ( c + 1_pInt : c + ns ) = state ( g , ip , el ) % p ( ( 8 * ns + 2 * nt + 1 ) : ( 9 * ns + 2 * nt ) )
2011-03-31 14:51:43 +05:30
c = c + ns
case ( 'tau_slip' )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_postResults ( c + 1_pInt : c + ns ) = abs ( state ( g , ip , el ) % p ( ( 9 * ns + 2 * nt + 1 ) : ( 10 * ns + 2 * nt ) ) )
2010-09-29 12:05:08 +05:30
c = c + ns
case ( 'gdot_slip_edge' )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_postResults ( c + 1_pInt : c + ns ) = abs ( state ( g , ip , el ) % p ( ( 10 * ns + 2 * nt + 1 ) : ( 11 * ns + 2 * nt ) ) )
2010-09-29 12:05:08 +05:30
c = c + ns
case ( 'gdot_slip_screw' )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_postResults ( c + 1_pInt : c + ns ) = abs ( state ( g , ip , el ) % p ( ( 11 * ns + 2 * nt + 1 ) : ( 12 * ns + 2 * nt ) ) )
2011-03-31 14:51:43 +05:30
c = c + ns
case ( 'gdot_slip' )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_postResults ( c + 1_pInt : c + ns ) = abs ( state ( g , ip , el ) % p ( ( 10 * ns + 2 * nt + 1 ) : ( 11 * ns + 2 * nt ) ) ) + &
2011-03-31 14:51:43 +05:30
abs ( state ( g , ip , el ) % p ( ( 11 * ns + 2 * nt + 1 ) : ( 12 * ns + 2 * nt ) ) )
2010-09-29 12:05:08 +05:30
c = c + ns
case ( 'stressratio_edge_p' )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_postResults ( c + 1_pInt : c + ns ) = abs ( state ( g , ip , el ) % p ( ( 12 * ns + 2 * nt + 1 ) : ( 13 * ns + 2 * nt ) ) )
2010-09-29 12:05:08 +05:30
c = c + ns
case ( 'stressratio_screw_p' )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_postResults ( c + 1_pInt : c + ns ) = abs ( state ( g , ip , el ) % p ( ( 13 * ns + 2 * nt + 1 ) : ( 14 * ns + 2 * nt ) ) )
2010-09-29 12:05:08 +05:30
c = c + ns
2011-03-31 14:51:43 +05:30
case ( 'shear_system' )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_postResults ( c + 1_pInt : c + ns ) = abs ( state ( g , ip , el ) % p ( ( 2 * ns + 1 ) : ( 3 * ns ) ) )
2011-03-31 14:51:43 +05:30
c = c + ns
case ( 'shear_basal' )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_postResults ( c + 1_pInt : c + 1_pInt ) = sum ( abs ( state ( g , ip , el ) % p ( ( 2 * ns + 1 ) : ( 2 * ns + 3 ) ) ) )
c = c + 1_pInt
2011-03-31 14:51:43 +05:30
case ( 'shear_prism' )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_postResults ( c + 1_pInt : c + 1_pInt ) = sum ( abs ( state ( g , ip , el ) % p ( ( 2 * ns + 4 ) : ( 2 * ns + 6 ) ) ) )
c = c + 1_pInt
2011-03-31 14:51:43 +05:30
case ( 'shear_pyra' )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_postResults ( c + 1_pInt : c + 1_pInt ) = sum ( abs ( state ( g , ip , el ) % p ( ( 2 * ns + 7 ) : ( 2 * ns + 12 ) ) ) )
c = c + 1_pInt
2011-03-31 14:51:43 +05:30
case ( 'shear_pyrca' )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_postResults ( c + 1_pInt : c + 1_pInt ) = sum ( abs ( state ( g , ip , el ) % p ( ( 2 * ns + 13 ) : ( 2 * ns + 24 ) ) ) )
c = c + 1_pInt
2011-03-31 14:51:43 +05:30
case ( 'rhoedge_basal' )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_postResults ( c + 1_pInt : c + 1_pInt ) = sum ( state ( g , ip , el ) % p ( ( 1 ) : ( 3 ) ) )
c = c + 1_pInt
2011-03-31 14:51:43 +05:30
case ( 'rhoedge_prism' )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_postResults ( c + 1_pInt : c + 1_pInt ) = sum ( state ( g , ip , el ) % p ( ( 4 ) : ( 6 ) ) )
c = c + 1_pInt
2011-03-31 14:51:43 +05:30
case ( 'rhoedge_pyra' )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_postResults ( c + 1_pInt : c + 1_pInt ) = sum ( state ( g , ip , el ) % p ( ( 7 ) : ( 12 ) ) )
c = c + 1_pInt
2011-03-31 14:51:43 +05:30
case ( 'rhoedge_pyrca' )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_postResults ( c + 1_pInt : c + 1_pInt ) = sum ( state ( g , ip , el ) % p ( ( 13 ) : ( 24 ) ) )
c = c + 1_pInt
2011-03-31 14:51:43 +05:30
case ( 'rhoscrew_basal' )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_postResults ( c + 1_pInt : c + 1_pInt ) = sum ( state ( g , ip , el ) % p ( ( ns + 1 ) : ( ns + 3 ) ) )
c = c + 1_pInt
2011-03-31 14:51:43 +05:30
case ( 'rhoscrew_prism' )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_postResults ( c + 1_pInt : c + 1_pInt ) = sum ( state ( g , ip , el ) % p ( ( ns + 4 ) : ( ns + 6 ) ) )
c = c + 1_pInt
2011-03-31 14:51:43 +05:30
case ( 'rhoscrew_pyra' )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_postResults ( c + 1_pInt : c + 1_pInt ) = sum ( state ( g , ip , el ) % p ( ( ns + 7 ) : ( ns + 12 ) ) )
c = c + 1_pInt
2011-03-31 14:51:43 +05:30
case ( 'rhoscrew_pyrca' )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_postResults ( c + 1_pInt : c + 1_pInt ) = sum ( state ( g , ip , el ) % p ( ( ns + 13 ) : ( ns + 24 ) ) )
c = c + 1_pInt
2011-03-31 14:51:43 +05:30
case ( 'shear_total' )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_postResults ( c + 1_pInt : c + 1_pInt ) = sum ( abs ( state ( g , ip , el ) % p ( ( 2 * ns + 1 ) : ( 3 * ns ) ) ) )
c = c + 1_pInt
2010-07-13 13:49:25 +05:30
case ( 'twin_fraction' )
2012-02-21 21:30:00 +05:30
constitutive_titanmod_postResults ( c + 1_pInt : c + nt ) = abs ( volumefraction_pertwinsystem ( 1 : nt ) )
2010-07-13 13:49:25 +05:30
c = c + nt
2010-09-29 12:05:08 +05:30
! 'rhoedge', &
! 'rhoscrew', &
! 'segment_edge', &
! 'segment_screw', &
! 'resistance_edge', &
! 'resistance_screw', &
! 'velocity_edge', &
! 'velocity_screw', &
2011-03-31 14:51:43 +05:30
! 'tau_slip', &
2010-09-29 12:05:08 +05:30
! 'gdot_slip_edge' &
! 'gdot_slip_screw' &
2011-03-31 14:51:43 +05:30
! 'gdot_slip', &
2010-09-29 12:05:08 +05:30
! 'StressRatio_edge_p' &
! 'StressRatio_screw_p'
2011-03-31 14:51:43 +05:30
! 'shear_total', &
! 'shear_basal' &
! 'shear_prism', &
! 'shear_pyra', &
! 'shear_pyrca', &
! 'twin_fraction', &
2010-07-13 13:49:25 +05:30
end select
enddo
return
end function
END MODULE