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/>.
!
!##############################################################
2009-08-31 19:43:10 +05:30
!* $Id$
2009-07-22 21:37:19 +05:30
!*****************************************************
!* Module: CONSTITUTIVE_PHENOPOWERLAW *
!*****************************************************
!* contains: *
!* - constitutive equations *
!* - parameters definition *
!*****************************************************
2009-10-21 18:40:12 +05:30
![Alu]
2012-03-12 20:13:19 +05:30
!plasticity phenopowerlaw
2009-10-21 18:40:12 +05:30
!(output) resistance_slip
!(output) shearrate_slip
!(output) resolvedstress_slip
!(output) totalshear
!(output) resistance_twin
!(output) shearrate_twin
!(output) resolvedstress_twin
!(output) totalvolfrac
!lattice_structure hex
!covera_ratio 1.587
!Nslip 3 3 6 12 # per family
!Ntwin 6 6 6 6 # per family
!
!c11 162.2e9
!c12 91.8e9
!c13 68.8e9
!c33 180.5e9
!c44 46.7e9
!
!gdot0_slip 0.001
2010-11-03 20:28:11 +05:30
!n_slip 50
2009-10-21 18:40:12 +05:30
!tau0_slip 65e6 22e6 52e6 50e6 # per family
!tausat_slip 80e6 180e6 140e6 140e6 # per family
2011-11-23 20:18:39 +05:30
!a_slip 1
2009-10-21 18:40:12 +05:30
!gdot0_twin 0.001
2010-11-03 20:28:11 +05:30
!n_twin 50
2009-10-21 18:40:12 +05:30
!tau0_twin 52e6 52e6 52e6 52e6 # per family
!s_pr 50e6 # push-up stress for slip saturation due to twinning
!twin_b 2
!twin_C 25
!twin_d 0.1
!twin_e 0.1
!h0_slipslip 10e6
!h0_sliptwin 0
!h0_twinslip 625e6
!h0_twintwin 400e6
!interaction_slipslip 5.5 5.5 1.0 52.0 5.5 5.5 1.0 52.0 27.5 0.2 72.8 1.0 72.8 72.8 27.5 1.1 1.4 5.5 7.7 7.7
!interaction_sliptwin 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0
!interaction_twinslip 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1
!interaction_twintwin 1 1 1 1 1 1 1 1 10 10 10 10 10 10 10 10 10 10 10 10
!relevantResistance 1
2009-07-22 21:37:19 +05:30
2012-03-09 01:55:28 +05:30
module constitutive_phenopowerlaw
2009-07-22 21:37:19 +05:30
use prec , only : pReal , pInt
2012-03-09 01:55:28 +05:30
implicit none
2012-04-11 19:31:02 +05:30
private
character ( len = * ) , parameter , public :: &
2012-03-09 01:55:28 +05:30
constitutive_phenopowerlaw_label = 'phenopowerlaw'
2012-04-11 19:31:02 +05:30
integer ( pInt ) , dimension ( : ) , allocatable , public :: &
2012-03-09 01:55:28 +05:30
constitutive_phenopowerlaw_sizeDotState , &
constitutive_phenopowerlaw_sizeState , &
constitutive_phenopowerlaw_sizePostResults , & ! cumulative size of post results
constitutive_phenopowerlaw_structure
2012-04-11 19:31:02 +05:30
integer ( pInt ) , dimension ( : ) , allocatable , private :: &
constitutive_phenopowerlaw_Noutput , & ! number of outputs per instance of this constitution
constitutive_phenopowerlaw_totalNslip , & ! no. of slip system used in simulation
constitutive_phenopowerlaw_totalNtwin ! no. of twin system used in simulation
2009-07-22 21:37:19 +05:30
2012-04-11 19:31:02 +05:30
integer ( pInt ) , dimension ( : , : ) , allocatable , target , public :: &
2012-03-09 01:55:28 +05:30
constitutive_phenopowerlaw_sizePostResult ! size of each post result output
2012-04-11 19:31:02 +05:30
integer ( pInt ) , dimension ( : , : ) , allocatable , private :: &
2012-03-09 01:55:28 +05:30
constitutive_phenopowerlaw_Nslip , & ! active number of slip systems per family
constitutive_phenopowerlaw_Ntwin ! active number of twin systems per family
2012-04-11 19:31:02 +05:30
character ( len = 64 ) , dimension ( : , : ) , allocatable , target , public :: &
2012-03-09 01:55:28 +05:30
constitutive_phenopowerlaw_output ! name of each post result output
2012-04-11 19:31:02 +05:30
character ( len = 32 ) , dimension ( : ) , allocatable , private :: &
2012-03-09 01:55:28 +05:30
constitutive_phenopowerlaw_structureName
2012-04-11 19:31:02 +05:30
real ( pReal ) , dimension ( : ) , allocatable , private :: &
2012-03-09 01:55:28 +05:30
constitutive_phenopowerlaw_CoverA , &
constitutive_phenopowerlaw_C11 , &
constitutive_phenopowerlaw_C12 , &
constitutive_phenopowerlaw_C13 , &
constitutive_phenopowerlaw_C33 , &
constitutive_phenopowerlaw_C44 , &
constitutive_phenopowerlaw_gdot0_slip , &
constitutive_phenopowerlaw_n_slip , &
constitutive_phenopowerlaw_n_twin , &
constitutive_phenopowerlaw_gdot0_twin
2012-04-11 19:31:02 +05:30
real ( pReal ) , dimension ( : , : ) , allocatable , private :: &
2012-03-09 01:55:28 +05:30
constitutive_phenopowerlaw_tau0_slip , &
constitutive_phenopowerlaw_tausat_slip , &
constitutive_phenopowerlaw_tau0_twin
2012-04-11 19:31:02 +05:30
real ( pReal ) , dimension ( : ) , allocatable , private :: &
2012-03-09 01:55:28 +05:30
constitutive_phenopowerlaw_spr , &
constitutive_phenopowerlaw_twinB , &
constitutive_phenopowerlaw_twinC , &
constitutive_phenopowerlaw_twinD , &
constitutive_phenopowerlaw_twinE , &
constitutive_phenopowerlaw_h0_slipslip , &
constitutive_phenopowerlaw_h0_sliptwin , &
constitutive_phenopowerlaw_h0_twinslip , &
constitutive_phenopowerlaw_h0_twintwin , &
constitutive_phenopowerlaw_a_slip , &
constitutive_phenopowerlaw_aTolResistance
2012-04-11 19:31:02 +05:30
real ( pReal ) , dimension ( : , : ) , allocatable , private :: &
2012-03-09 01:55:28 +05:30
constitutive_phenopowerlaw_interaction_slipslip , &
constitutive_phenopowerlaw_interaction_sliptwin , &
constitutive_phenopowerlaw_interaction_twinslip , &
constitutive_phenopowerlaw_interaction_twintwin
2012-04-11 19:31:02 +05:30
real ( pReal ) , dimension ( : , : , : ) , allocatable , private :: &
2012-03-09 01:55:28 +05:30
constitutive_phenopowerlaw_hardeningMatrix_slipslip , &
constitutive_phenopowerlaw_hardeningMatrix_sliptwin , &
constitutive_phenopowerlaw_hardeningMatrix_twinslip , &
constitutive_phenopowerlaw_hardeningMatrix_twintwin , &
constitutive_phenopowerlaw_Cslip_66
2012-04-11 19:31:02 +05:30
public :: &
constitutive_phenopowerlaw_init , &
constitutive_phenopowerlaw_homogenizedC , &
constitutive_phenopowerlaw_aTolState , &
constitutive_phenopowerlaw_dotState , &
2012-05-16 20:13:26 +05:30
constitutive_phenopowerlaw_deltaState , &
2012-04-11 19:31:02 +05:30
constitutive_phenopowerlaw_dotTemperature , &
constitutive_phenopowerlaw_microstructure , &
constitutive_phenopowerlaw_LpAndItsTangent , &
constitutive_phenopowerlaw_postResults , &
constitutive_phenopowerlaw_stateInit
2012-03-09 01:55:28 +05:30
contains
2009-07-22 21:37:19 +05:30
2012-02-21 21:30:00 +05:30
subroutine constitutive_phenopowerlaw_init ( myFile )
2009-07-22 21:37:19 +05:30
!**************************************
!* 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)
2012-03-09 01:55:28 +05:30
use math , only : math_Mandel3333to66 , &
math_Voigt66to3333
2009-07-22 21:37:19 +05:30
use IO
use material
2012-07-05 15:24:50 +05:30
use debug , only : debug_level , &
2012-03-09 01:55:28 +05:30
debug_constitutive , &
debug_levelBasic
2010-11-03 20:28:11 +05:30
use lattice , only : lattice_initializeStructure , lattice_symmetryType , &
lattice_maxNslipFamily , lattice_maxNtwinFamily , &
2009-07-22 21:37:19 +05:30
lattice_maxNinteraction , lattice_NslipSystem , lattice_NtwinSystem , &
2009-10-21 18:40:12 +05:30
lattice_interactionSlipSlip , &
lattice_interactionSlipTwin , &
lattice_interactionTwinSlip , &
lattice_interactionTwinTwin
2012-03-09 01:55:28 +05:30
implicit none
2012-02-21 21:30:00 +05:30
integer ( pInt ) , intent ( in ) :: myFile
2011-10-18 20:16:07 +05:30
integer ( pInt ) , parameter :: maxNchunks = lattice_maxNinteraction + 1_pInt
2009-07-22 21:37:19 +05:30
integer ( pInt ) , dimension ( 1 + 2 * maxNchunks ) :: positions
2012-02-14 20:49:59 +05:30
integer ( pInt ) section , maxNinstance , i , j , k , f , o , &
2009-10-21 18:40:12 +05:30
mySize , myStructure , index_myFamily , index_otherFamily
2012-03-09 01:55:28 +05:30
character ( len = 64 ) :: tag
character ( len = 1024 ) :: line
2009-07-22 21:37:19 +05:30
openmp parallelization working again (at least for j2 and nonlocal constitutive model).
In order to keep it like that, please follow these simple rules:
DON'T use implicit array subscripts:
example: real, dimension(3,3) :: A,B
A(:,2) = B(:,1) <--- DON'T USE
A(1:3,2) = B(1:3,1) <--- BETTER USE
In many cases the use of explicit array subscripts is inevitable for parallelization. Additionally, it is an easy means to prevent memory leaks.
Enclose all write statements with the following:
!$OMP CRITICAL (write2out)
<your write statement>
!$OMP END CRITICAL (write2out)
Whenever you change something in the code and are not sure if it affects parallelization and leads to nonconforming behavior, please ask me and/or Franz to check this.
2011-03-17 16:16:17 +05:30
!$OMP CRITICAL (write2out)
2011-03-21 16:01:17 +05:30
write ( 6 , * )
2011-08-26 19:27:29 +05:30
write ( 6 , * ) '<<<+- constitutive_' , trim ( constitutive_phenopowerlaw_label ) , ' init -+>>>'
2011-03-21 16:01:17 +05:30
write ( 6 , * ) '$Id$'
2012-02-01 00:48:55 +05:30
#include "compilation_info.f90"
openmp parallelization working again (at least for j2 and nonlocal constitutive model).
In order to keep it like that, please follow these simple rules:
DON'T use implicit array subscripts:
example: real, dimension(3,3) :: A,B
A(:,2) = B(:,1) <--- DON'T USE
A(1:3,2) = B(1:3,1) <--- BETTER USE
In many cases the use of explicit array subscripts is inevitable for parallelization. Additionally, it is an easy means to prevent memory leaks.
Enclose all write statements with the following:
!$OMP CRITICAL (write2out)
<your write statement>
!$OMP END CRITICAL (write2out)
Whenever you change something in the code and are not sure if it affects parallelization and leads to nonconforming behavior, please ask me and/or Franz to check this.
2011-03-17 16:16:17 +05:30
!$OMP END CRITICAL (write2out)
2009-07-22 21:37:19 +05:30
2012-03-12 19:39:37 +05:30
maxNinstance = int ( count ( phase_plasticity == constitutive_phenopowerlaw_label ) , pInt )
2009-07-22 21:37:19 +05:30
if ( maxNinstance == 0 ) return
2009-10-16 01:32:52 +05:30
2012-07-05 15:24:50 +05:30
if ( iand ( debug_level ( debug_constitutive ) , debug_levelBasic ) / = 0_pInt ) then
2011-03-21 16:01:17 +05:30
!$OMP CRITICAL (write2out)
2012-02-01 00:48:55 +05:30
write ( 6 , '(a16,1x,i5)' ) '# instances:' , maxNinstance
2011-03-21 16:01:17 +05:30
write ( 6 , * )
!$OMP END CRITICAL (write2out)
endif
2012-03-09 01:55:28 +05:30
allocate ( constitutive_phenopowerlaw_sizeDotState ( maxNinstance ) )
constitutive_phenopowerlaw_sizeDotState = 0_pInt
allocate ( constitutive_phenopowerlaw_sizeState ( maxNinstance ) )
constitutive_phenopowerlaw_sizeState = 0_pInt
allocate ( constitutive_phenopowerlaw_sizePostResults ( maxNinstance ) )
constitutive_phenopowerlaw_sizePostResults = 0_pInt
allocate ( constitutive_phenopowerlaw_sizePostResult ( maxval ( phase_Noutput ) , maxNinstance ) )
constitutive_phenopowerlaw_sizePostResult = 0_pInt
allocate ( constitutive_phenopowerlaw_output ( maxval ( phase_Noutput ) , maxNinstance ) )
constitutive_phenopowerlaw_output = ''
allocate ( constitutive_phenopowerlaw_Noutput ( maxNinstance ) )
constitutive_phenopowerlaw_Noutput = 0_pInt
allocate ( constitutive_phenopowerlaw_structureName ( maxNinstance ) )
constitutive_phenopowerlaw_structureName = ''
allocate ( constitutive_phenopowerlaw_structure ( maxNinstance ) )
constitutive_phenopowerlaw_structure = 0_pInt
allocate ( constitutive_phenopowerlaw_Nslip ( lattice_maxNslipFamily , maxNinstance ) )
constitutive_phenopowerlaw_Nslip = 0_pInt
allocate ( constitutive_phenopowerlaw_Ntwin ( lattice_maxNtwinFamily , maxNinstance ) )
constitutive_phenopowerlaw_Ntwin = 0_pInt
allocate ( constitutive_phenopowerlaw_totalNslip ( maxNinstance ) )
constitutive_phenopowerlaw_totalNslip = 0_pInt
allocate ( constitutive_phenopowerlaw_totalNtwin ( maxNinstance ) )
constitutive_phenopowerlaw_totalNtwin = 0_pInt
allocate ( constitutive_phenopowerlaw_CoverA ( maxNinstance ) )
constitutive_phenopowerlaw_CoverA = 0.0_pReal
allocate ( constitutive_phenopowerlaw_C11 ( maxNinstance ) )
constitutive_phenopowerlaw_C11 = 0.0_pReal
allocate ( constitutive_phenopowerlaw_C12 ( maxNinstance ) )
constitutive_phenopowerlaw_C12 = 0.0_pReal
allocate ( constitutive_phenopowerlaw_C13 ( maxNinstance ) )
constitutive_phenopowerlaw_C13 = 0.0_pReal
allocate ( constitutive_phenopowerlaw_C33 ( maxNinstance ) )
constitutive_phenopowerlaw_C33 = 0.0_pReal
allocate ( constitutive_phenopowerlaw_C44 ( maxNinstance ) )
constitutive_phenopowerlaw_C44 = 0.0_pReal
allocate ( constitutive_phenopowerlaw_Cslip_66 ( 6 , 6 , maxNinstance ) )
constitutive_phenopowerlaw_Cslip_66 = 0.0_pReal
allocate ( constitutive_phenopowerlaw_gdot0_slip ( maxNinstance ) )
constitutive_phenopowerlaw_gdot0_slip = 0.0_pReal
allocate ( constitutive_phenopowerlaw_n_slip ( maxNinstance ) )
constitutive_phenopowerlaw_n_slip = 0.0_pReal
allocate ( constitutive_phenopowerlaw_tau0_slip ( lattice_maxNslipFamily , maxNinstance ) )
constitutive_phenopowerlaw_tau0_slip = 0.0_pReal
allocate ( constitutive_phenopowerlaw_tausat_slip ( lattice_maxNslipFamily , maxNinstance ) )
constitutive_phenopowerlaw_tausat_slip = 0.0_pReal
allocate ( constitutive_phenopowerlaw_gdot0_twin ( maxNinstance ) )
constitutive_phenopowerlaw_gdot0_twin = 0.0_pReal
allocate ( constitutive_phenopowerlaw_n_twin ( maxNinstance ) )
constitutive_phenopowerlaw_n_twin = 0.0_pReal
allocate ( constitutive_phenopowerlaw_tau0_twin ( lattice_maxNtwinFamily , maxNinstance ) )
constitutive_phenopowerlaw_tau0_twin = 0.0_pReal
allocate ( constitutive_phenopowerlaw_spr ( maxNinstance ) )
constitutive_phenopowerlaw_spr = 0.0_pReal
allocate ( constitutive_phenopowerlaw_twinB ( maxNinstance ) )
constitutive_phenopowerlaw_twinB = 0.0_pReal
allocate ( constitutive_phenopowerlaw_twinC ( maxNinstance ) )
constitutive_phenopowerlaw_twinC = 0.0_pReal
allocate ( constitutive_phenopowerlaw_twinD ( maxNinstance ) )
constitutive_phenopowerlaw_twinD = 0.0_pReal
allocate ( constitutive_phenopowerlaw_twinE ( maxNinstance ) )
constitutive_phenopowerlaw_twinE = 0.0_pReal
allocate ( constitutive_phenopowerlaw_h0_slipslip ( maxNinstance ) )
constitutive_phenopowerlaw_h0_slipslip = 0.0_pReal
allocate ( constitutive_phenopowerlaw_h0_sliptwin ( maxNinstance ) )
constitutive_phenopowerlaw_h0_sliptwin = 0.0_pReal
allocate ( constitutive_phenopowerlaw_h0_twinslip ( maxNinstance ) )
constitutive_phenopowerlaw_h0_twinslip = 0.0_pReal
allocate ( constitutive_phenopowerlaw_h0_twintwin ( maxNinstance ) )
constitutive_phenopowerlaw_h0_twintwin = 0.0_pReal
2009-08-03 12:07:37 +05:30
allocate ( constitutive_phenopowerlaw_interaction_slipslip ( lattice_maxNinteraction , maxNinstance ) )
2012-03-09 01:55:28 +05:30
constitutive_phenopowerlaw_interaction_slipslip = 0.0_pReal
2009-08-03 12:07:37 +05:30
allocate ( constitutive_phenopowerlaw_interaction_sliptwin ( lattice_maxNinteraction , maxNinstance ) )
2012-03-09 01:55:28 +05:30
constitutive_phenopowerlaw_interaction_sliptwin = 0.0_pReal
2009-08-03 12:07:37 +05:30
allocate ( constitutive_phenopowerlaw_interaction_twinslip ( lattice_maxNinteraction , maxNinstance ) )
2012-03-09 01:55:28 +05:30
constitutive_phenopowerlaw_interaction_twinslip = 0.0_pReal
2009-08-03 12:07:37 +05:30
allocate ( constitutive_phenopowerlaw_interaction_twintwin ( lattice_maxNinteraction , maxNinstance ) )
2012-03-09 01:55:28 +05:30
constitutive_phenopowerlaw_interaction_twintwin = 0.0_pReal
2011-11-23 20:18:39 +05:30
allocate ( constitutive_phenopowerlaw_a_slip ( maxNinstance ) )
2012-03-09 01:55:28 +05:30
constitutive_phenopowerlaw_a_slip = 0.0_pReal
2010-10-26 18:46:37 +05:30
allocate ( constitutive_phenopowerlaw_aTolResistance ( maxNinstance ) )
2012-03-09 01:55:28 +05:30
constitutive_phenopowerlaw_aTolResistance = 0.0_pReal
2009-08-27 17:40:06 +05:30
2012-02-21 21:30:00 +05:30
rewind ( myFile )
section = 0_pInt
2009-10-16 01:32:52 +05:30
2009-07-22 21:37:19 +05:30
do while ( IO_lc ( IO_getTag ( line , '<' , '>' ) ) / = 'phase' ) ! wind forward to <phase>
2012-02-21 21:30:00 +05:30
read ( myFile , '(a1024)' , END = 100 ) line
2009-07-22 21:37:19 +05:30
enddo
do ! read thru sections of phase part
2012-02-21 21:30:00 +05:30
read ( myFile , '(a1024)' , END = 100 ) line
2009-07-22 21:37:19 +05:30
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
2009-07-22 21:37:19 +05:30
endif
2012-03-12 19:39:37 +05:30
if ( section > 0_pInt . and . phase_plasticity ( section ) == constitutive_phenopowerlaw_label ) then ! one of my sections
2012-03-12 20:13:19 +05:30
i = phase_plasticityInstance ( section ) ! which instance of my plasticity is present phase
2009-07-22 21:37:19 +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
2009-07-22 21:37:19 +05:30
select case ( tag )
2012-03-15 14:52:24 +05:30
case ( 'plasticity' , 'elasticity' )
2012-02-14 14:52:37 +05:30
cycle
2009-07-22 21:37:19 +05:30
case ( '(output)' )
2012-02-14 20:49:59 +05:30
constitutive_phenopowerlaw_Noutput ( i ) = constitutive_phenopowerlaw_Noutput ( i ) + 1_pInt
constitutive_phenopowerlaw_output ( constitutive_phenopowerlaw_Noutput ( i ) , i ) = IO_lc ( IO_stringValue ( line , positions , 2_pInt ) )
2009-07-22 21:37:19 +05:30
case ( 'lattice_structure' )
2012-02-13 19:48:07 +05:30
constitutive_phenopowerlaw_structureName ( i ) = IO_lc ( IO_stringValue ( line , positions , 2_pInt ) )
2009-07-22 21:37:19 +05:30
case ( 'covera_ratio' )
2012-02-13 19:48:07 +05:30
constitutive_phenopowerlaw_CoverA ( i ) = IO_floatValue ( line , positions , 2_pInt )
2009-07-22 21:37:19 +05:30
case ( 'c11' )
2012-02-13 19:48:07 +05:30
constitutive_phenopowerlaw_C11 ( i ) = IO_floatValue ( line , positions , 2_pInt )
2009-07-22 21:37:19 +05:30
case ( 'c12' )
2012-02-13 19:48:07 +05:30
constitutive_phenopowerlaw_C12 ( i ) = IO_floatValue ( line , positions , 2_pInt )
2009-07-22 21:37:19 +05:30
case ( 'c13' )
2012-02-13 19:48:07 +05:30
constitutive_phenopowerlaw_C13 ( i ) = IO_floatValue ( line , positions , 2_pInt )
2009-07-22 21:37:19 +05:30
case ( 'c33' )
2012-02-13 19:48:07 +05:30
constitutive_phenopowerlaw_C33 ( i ) = IO_floatValue ( line , positions , 2_pInt )
2009-07-22 21:37:19 +05:30
case ( 'c44' )
2012-02-13 19:48:07 +05:30
constitutive_phenopowerlaw_C44 ( i ) = IO_floatValue ( line , positions , 2_pInt )
2009-07-22 21:37:19 +05:30
case ( 'nslip' )
2012-02-13 19:48:07 +05:30
forall ( j = 1_pInt : lattice_maxNslipFamily ) &
constitutive_phenopowerlaw_Nslip ( j , i ) = IO_intValue ( line , positions , 1_pInt + j )
2009-07-22 21:37:19 +05:30
case ( 'gdot0_slip' )
2012-02-13 19:48:07 +05:30
constitutive_phenopowerlaw_gdot0_slip ( i ) = IO_floatValue ( line , positions , 2_pInt )
2009-07-22 21:37:19 +05:30
case ( 'n_slip' )
2012-02-13 19:48:07 +05:30
constitutive_phenopowerlaw_n_slip ( i ) = IO_floatValue ( line , positions , 2_pInt )
2009-07-22 21:37:19 +05:30
case ( 'tau0_slip' )
2012-02-13 19:48:07 +05:30
forall ( j = 1_pInt : lattice_maxNslipFamily ) &
constitutive_phenopowerlaw_tau0_slip ( j , i ) = IO_floatValue ( line , positions , 1_pInt + j )
2009-07-22 21:37:19 +05:30
case ( 'tausat_slip' )
2012-02-13 19:48:07 +05:30
forall ( j = 1_pInt : lattice_maxNslipFamily ) &
constitutive_phenopowerlaw_tausat_slip ( j , i ) = IO_floatValue ( line , positions , 1_pInt + j )
2011-11-23 20:18:39 +05:30
case ( 'a_slip' , 'w0_slip' )
2012-02-13 19:48:07 +05:30
constitutive_phenopowerlaw_a_slip ( i ) = IO_floatValue ( line , positions , 2_pInt )
2009-07-22 21:37:19 +05:30
case ( 'ntwin' )
2012-02-13 19:48:07 +05:30
forall ( j = 1_pInt : lattice_maxNtwinFamily ) &
constitutive_phenopowerlaw_Ntwin ( j , i ) = IO_intValue ( line , positions , 1_pInt + j )
2009-07-22 21:37:19 +05:30
case ( 'gdot0_twin' )
2012-02-13 19:48:07 +05:30
constitutive_phenopowerlaw_gdot0_twin ( i ) = IO_floatValue ( line , positions , 2_pInt )
2009-07-22 21:37:19 +05:30
case ( 'n_twin' )
2012-02-13 19:48:07 +05:30
constitutive_phenopowerlaw_n_twin ( i ) = IO_floatValue ( line , positions , 2_pInt )
2009-07-22 21:37:19 +05:30
case ( 'tau0_twin' )
2012-02-13 19:48:07 +05:30
forall ( j = 1_pInt : lattice_maxNtwinFamily ) &
constitutive_phenopowerlaw_tau0_twin ( j , i ) = IO_floatValue ( line , positions , 1_pInt + j )
2009-07-22 21:37:19 +05:30
case ( 's_pr' )
2012-02-13 19:48:07 +05:30
constitutive_phenopowerlaw_spr ( i ) = IO_floatValue ( line , positions , 2_pInt )
2009-07-22 21:37:19 +05:30
case ( 'twin_b' )
2012-02-13 19:48:07 +05:30
constitutive_phenopowerlaw_twinB ( i ) = IO_floatValue ( line , positions , 2_pInt )
2009-07-22 21:37:19 +05:30
case ( 'twin_c' )
2012-02-13 19:48:07 +05:30
constitutive_phenopowerlaw_twinC ( i ) = IO_floatValue ( line , positions , 2_pInt )
2009-07-22 21:37:19 +05:30
case ( 'twin_d' )
2012-02-13 19:48:07 +05:30
constitutive_phenopowerlaw_twinD ( i ) = IO_floatValue ( line , positions , 2_pInt )
2009-07-22 21:37:19 +05:30
case ( 'twin_e' )
2012-02-13 19:48:07 +05:30
constitutive_phenopowerlaw_twinE ( i ) = IO_floatValue ( line , positions , 2_pInt )
2009-07-22 21:37:19 +05:30
case ( 'h0_slipslip' )
2012-02-13 19:48:07 +05:30
constitutive_phenopowerlaw_h0_slipslip ( i ) = IO_floatValue ( line , positions , 2_pInt )
2009-07-22 21:37:19 +05:30
case ( 'h0_sliptwin' )
2012-02-13 19:48:07 +05:30
constitutive_phenopowerlaw_h0_sliptwin ( i ) = IO_floatValue ( line , positions , 2_pInt )
2009-07-22 21:37:19 +05:30
case ( 'h0_twinslip' )
2012-02-13 19:48:07 +05:30
constitutive_phenopowerlaw_h0_twinslip ( i ) = IO_floatValue ( line , positions , 2_pInt )
2009-07-22 21:37:19 +05:30
case ( 'h0_twintwin' )
2012-02-13 19:48:07 +05:30
constitutive_phenopowerlaw_h0_twintwin ( i ) = IO_floatValue ( line , positions , 2_pInt )
2010-10-26 18:46:37 +05:30
case ( 'atol_resistance' )
2012-02-13 19:48:07 +05:30
constitutive_phenopowerlaw_aTolResistance ( i ) = IO_floatValue ( line , positions , 2_pInt )
2009-07-22 21:37:19 +05:30
case ( 'interaction_slipslip' )
2012-02-13 19:48:07 +05:30
forall ( j = 1_pInt : lattice_maxNinteraction ) &
constitutive_phenopowerlaw_interaction_slipslip ( j , i ) = IO_floatValue ( line , positions , 1_pInt + j )
2009-07-22 21:37:19 +05:30
case ( 'interaction_sliptwin' )
2012-02-13 19:48:07 +05:30
forall ( j = 1_pInt : lattice_maxNinteraction ) &
constitutive_phenopowerlaw_interaction_sliptwin ( j , i ) = IO_floatValue ( line , positions , 1_pInt + j )
2009-07-22 21:37:19 +05:30
case ( 'interaction_twinslip' )
2012-02-13 19:48:07 +05:30
forall ( j = 1_pInt : lattice_maxNinteraction ) &
constitutive_phenopowerlaw_interaction_twinslip ( j , i ) = IO_floatValue ( line , positions , 1_pInt + j )
2009-07-22 21:37:19 +05:30
case ( 'interaction_twintwin' )
2012-02-13 19:48:07 +05:30
forall ( j = 1_pInt : lattice_maxNinteraction ) &
constitutive_phenopowerlaw_interaction_twintwin ( j , i ) = IO_floatValue ( line , positions , 1_pInt + j )
2012-02-14 14:52:37 +05:30
case default
call IO_error ( 220_pInt , ext_msg = tag )
2009-07-22 21:37:19 +05:30
end select
endif
enddo
2012-02-13 19:48:07 +05:30
100 do i = 1_pInt , maxNinstance
2009-10-16 01:32:52 +05:30
2009-07-22 21:37:19 +05:30
constitutive_phenopowerlaw_structure ( i ) = lattice_initializeStructure ( constitutive_phenopowerlaw_structureName ( i ) , & ! get structure
constitutive_phenopowerlaw_CoverA ( i ) )
2011-03-24 22:50:35 +05:30
constitutive_phenopowerlaw_Nslip ( 1 : lattice_maxNslipFamily , i ) = &
min ( lattice_NslipSystem ( 1 : lattice_maxNslipFamily , constitutive_phenopowerlaw_structure ( i ) ) , & ! limit active slip systems per family to min of available and requested
constitutive_phenopowerlaw_Nslip ( 1 : lattice_maxNslipFamily , i ) )
constitutive_phenopowerlaw_Ntwin ( 1 : lattice_maxNtwinFamily , i ) = &
min ( lattice_NtwinSystem ( 1 : lattice_maxNtwinFamily , constitutive_phenopowerlaw_structure ( i ) ) , & ! limit active twin systems per family to min of available and requested
constitutive_phenopowerlaw_Ntwin ( : , i ) )
2009-07-22 21:37:19 +05:30
constitutive_phenopowerlaw_totalNslip ( i ) = sum ( constitutive_phenopowerlaw_Nslip ( : , i ) ) ! how many slip systems altogether
constitutive_phenopowerlaw_totalNtwin ( i ) = sum ( constitutive_phenopowerlaw_Ntwin ( : , i ) ) ! how many twin systems altogether
2012-02-13 23:11:27 +05:30
if ( constitutive_phenopowerlaw_structure ( i ) < 1 ) call IO_error ( 205_pInt , e = i )
2009-07-22 21:37:19 +05:30
if ( any ( constitutive_phenopowerlaw_tau0_slip ( : , i ) < 0.0_pReal . and . &
2012-02-13 23:11:27 +05:30
constitutive_phenopowerlaw_Nslip ( : , i ) > 0 ) ) call IO_error ( 221_pInt , e = i , ext_msg = 'tau0_slip' )
if ( constitutive_phenopowerlaw_gdot0_slip ( i ) < = 0.0_pReal ) call IO_error ( 221_pInt , e = i , ext_msg = 'gdot0_slip' )
if ( constitutive_phenopowerlaw_n_slip ( i ) < = 0.0_pReal ) call IO_error ( 221_pInt , e = i , ext_msg = 'n_slip' )
2009-07-22 21:37:19 +05:30
if ( any ( constitutive_phenopowerlaw_tausat_slip ( : , i ) < = 0.0_pReal . and . &
2012-02-13 23:11:27 +05:30
constitutive_phenopowerlaw_Nslip ( : , i ) > 0 ) ) call IO_error ( 221_pInt , e = i , ext_msg = 'tausat_slip' )
2011-11-23 20:18:39 +05:30
if ( any ( constitutive_phenopowerlaw_a_slip ( i ) == 0.0_pReal . and . &
2012-02-13 23:11:27 +05:30
constitutive_phenopowerlaw_Nslip ( : , i ) > 0 ) ) call IO_error ( 221_pInt , e = i , ext_msg = 'a_slip' )
2009-07-22 21:37:19 +05:30
if ( any ( constitutive_phenopowerlaw_tau0_twin ( : , i ) < 0.0_pReal . and . &
2012-02-13 23:11:27 +05:30
constitutive_phenopowerlaw_Ntwin ( : , i ) > 0 ) ) call IO_error ( 221_pInt , e = i , ext_msg = 'tau0_twin' )
2009-09-18 21:07:14 +05:30
if ( constitutive_phenopowerlaw_gdot0_twin ( i ) < = 0.0_pReal . and . &
2012-02-13 23:11:27 +05:30
any ( constitutive_phenopowerlaw_Ntwin ( : , i ) > 0 ) ) call IO_error ( 221_pInt , e = i , ext_msg = 'gdot0_twin' )
2009-08-13 19:02:17 +05:30
if ( constitutive_phenopowerlaw_n_twin ( i ) < = 0.0_pReal . and . &
2012-02-13 23:11:27 +05:30
any ( constitutive_phenopowerlaw_Ntwin ( : , i ) > 0 ) ) call IO_error ( 221_pInt , e = i , ext_msg = 'n_twin' )
2010-10-26 18:46:37 +05:30
if ( constitutive_phenopowerlaw_aTolResistance ( i ) < = 0.0_pReal ) &
constitutive_phenopowerlaw_aTolResistance ( i ) = 1.0_pReal ! default absolute tolerance 1 Pa
2009-07-22 21:37:19 +05:30
enddo
2009-10-22 14:28:14 +05:30
2011-03-24 22:50:35 +05:30
allocate ( constitutive_phenopowerlaw_hardeningMatrix_slipslip ( maxval ( constitutive_phenopowerlaw_totalNslip ) , & ! slip resistance from slip activity
2009-07-22 21:37:19 +05:30
maxval ( constitutive_phenopowerlaw_totalNslip ) , &
maxNinstance ) )
2011-03-24 22:50:35 +05:30
allocate ( constitutive_phenopowerlaw_hardeningMatrix_sliptwin ( maxval ( constitutive_phenopowerlaw_totalNtwin ) , & ! slip resistance from twin activity
2009-07-22 21:37:19 +05:30
maxval ( constitutive_phenopowerlaw_totalNslip ) , &
maxNinstance ) )
2011-03-24 22:50:35 +05:30
allocate ( constitutive_phenopowerlaw_hardeningMatrix_twinslip ( maxval ( constitutive_phenopowerlaw_totalNslip ) , & ! twin resistance from slip activity
maxval ( constitutive_phenopowerlaw_totalNtwin ) , &
maxNinstance ) )
allocate ( constitutive_phenopowerlaw_hardeningMatrix_twintwin ( maxval ( constitutive_phenopowerlaw_totalNtwin ) , & ! twin resistance from twin activity
2009-07-22 21:37:19 +05:30
maxval ( constitutive_phenopowerlaw_totalNtwin ) , &
maxNinstance ) )
2009-07-23 19:03:53 +05:30
constitutive_phenopowerlaw_hardeningMatrix_slipslip = 0.0_pReal
constitutive_phenopowerlaw_hardeningMatrix_sliptwin = 0.0_pReal
constitutive_phenopowerlaw_hardeningMatrix_twinslip = 0.0_pReal
constitutive_phenopowerlaw_hardeningMatrix_twintwin = 0.0_pReal
2009-10-16 01:32:52 +05:30
2009-07-23 19:03:53 +05:30
2012-02-13 19:48:07 +05:30
do i = 1_pInt , maxNinstance
2012-02-14 20:49:59 +05:30
do j = 1_pInt , constitutive_phenopowerlaw_Noutput ( i )
2009-10-16 01:32:52 +05:30
select case ( constitutive_phenopowerlaw_output ( j , i ) )
case ( 'resistance_slip' , &
'shearrate_slip' , &
'resolvedstress_slip' &
)
mySize = constitutive_phenopowerlaw_totalNslip ( i )
case ( 'resistance_twin' , &
'shearrate_twin' , &
'resolvedstress_twin' &
)
mySize = constitutive_phenopowerlaw_totalNtwin ( i )
case ( 'totalshear' , &
'totalvolfrac' &
)
mySize = 1_pInt
case default
2012-02-13 23:11:27 +05:30
call IO_error ( 222_pInt , ext_msg = constitutive_phenopowerlaw_output ( j , i ) )
2009-10-16 01:32:52 +05:30
end select
if ( mySize > 0_pInt ) then ! any meaningful output found
constitutive_phenopowerlaw_sizePostResult ( j , i ) = mySize
constitutive_phenopowerlaw_sizePostResults ( i ) = &
constitutive_phenopowerlaw_sizePostResults ( i ) + mySize
endif
2009-07-22 21:37:19 +05:30
enddo
constitutive_phenopowerlaw_sizeDotState ( i ) = constitutive_phenopowerlaw_totalNslip ( i ) + &
2012-02-13 19:48:07 +05:30
constitutive_phenopowerlaw_totalNtwin ( i ) + 2_pInt ! s_slip, s_twin, sum(gamma), sum(f)
2009-07-22 21:37:19 +05:30
constitutive_phenopowerlaw_sizeState ( i ) = constitutive_phenopowerlaw_totalNslip ( i ) + &
2012-02-13 19:48:07 +05:30
constitutive_phenopowerlaw_totalNtwin ( i ) + 2_pInt ! s_slip, s_twin, sum(gamma), sum(f)
2009-07-22 21:37:19 +05:30
2010-11-03 20:28:11 +05:30
myStructure = constitutive_phenopowerlaw_structure ( i )
select case ( lattice_symmetryType ( myStructure ) ) ! assign elasticity tensor
2012-02-13 19:48:07 +05:30
case ( 1_pInt ) ! cubic(s)
forall ( k = 1_pInt : 3_pInt )
forall ( j = 1_pInt : 3_pInt ) &
2009-10-16 01:32:52 +05:30
constitutive_phenopowerlaw_Cslip_66 ( k , j , i ) = constitutive_phenopowerlaw_C12 ( i )
constitutive_phenopowerlaw_Cslip_66 ( k , k , i ) = constitutive_phenopowerlaw_C11 ( i )
2012-02-13 19:48:07 +05:30
constitutive_phenopowerlaw_Cslip_66 ( k + 3_pInt , k + 3_pInt , i ) = constitutive_phenopowerlaw_C44 ( i )
2009-10-16 01:32:52 +05:30
end forall
2012-02-13 19:48:07 +05:30
case ( 2_pInt ) ! hex
2009-10-16 01:32:52 +05:30
constitutive_phenopowerlaw_Cslip_66 ( 1 , 1 , i ) = constitutive_phenopowerlaw_C11 ( i )
constitutive_phenopowerlaw_Cslip_66 ( 2 , 2 , i ) = constitutive_phenopowerlaw_C11 ( i )
constitutive_phenopowerlaw_Cslip_66 ( 3 , 3 , i ) = constitutive_phenopowerlaw_C33 ( i )
constitutive_phenopowerlaw_Cslip_66 ( 1 , 2 , i ) = constitutive_phenopowerlaw_C12 ( i )
constitutive_phenopowerlaw_Cslip_66 ( 2 , 1 , i ) = constitutive_phenopowerlaw_C12 ( i )
constitutive_phenopowerlaw_Cslip_66 ( 1 , 3 , i ) = constitutive_phenopowerlaw_C13 ( i )
constitutive_phenopowerlaw_Cslip_66 ( 3 , 1 , i ) = constitutive_phenopowerlaw_C13 ( i )
constitutive_phenopowerlaw_Cslip_66 ( 2 , 3 , i ) = constitutive_phenopowerlaw_C13 ( i )
constitutive_phenopowerlaw_Cslip_66 ( 3 , 2 , i ) = constitutive_phenopowerlaw_C13 ( i )
constitutive_phenopowerlaw_Cslip_66 ( 4 , 4 , i ) = constitutive_phenopowerlaw_C44 ( i )
constitutive_phenopowerlaw_Cslip_66 ( 5 , 5 , i ) = constitutive_phenopowerlaw_C44 ( i )
constitutive_phenopowerlaw_Cslip_66 ( 6 , 6 , i ) = 0.5_pReal * ( constitutive_phenopowerlaw_C11 ( i ) - &
constitutive_phenopowerlaw_C12 ( i ) )
2009-07-22 21:37:19 +05:30
end select
constitutive_phenopowerlaw_Cslip_66 ( : , : , i ) = &
math_Mandel3333to66 ( math_Voigt66to3333 ( constitutive_phenopowerlaw_Cslip_66 ( : , : , i ) ) )
2012-02-13 19:48:07 +05:30
do f = 1_pInt , lattice_maxNslipFamily ! >>> interaction slip -- X
index_myFamily = sum ( constitutive_phenopowerlaw_Nslip ( 1 : f - 1_pInt , i ) )
do j = 1_pInt , constitutive_phenopowerlaw_Nslip ( f , i ) ! loop over (active) systems in my family (slip)
do o = 1_pInt , lattice_maxNslipFamily
index_otherFamily = sum ( constitutive_phenopowerlaw_Nslip ( 1 : o - 1_pInt , i ) )
do k = 1_pInt , constitutive_phenopowerlaw_Nslip ( o , i ) ! loop over (active) systems in other family (slip)
2009-10-21 18:40:12 +05:30
constitutive_phenopowerlaw_hardeningMatrix_slipslip ( index_otherFamily + k , index_myFamily + j , i ) = &
constitutive_phenopowerlaw_interaction_slipslip ( lattice_interactionSlipSlip ( &
sum ( lattice_NslipSystem ( 1 : o - 1 , myStructure ) ) + k , &
sum ( lattice_NslipSystem ( 1 : f - 1 , myStructure ) ) + j , &
myStructure ) , i )
enddo ; enddo
2012-02-13 19:48:07 +05:30
do o = 1_pInt , lattice_maxNtwinFamily
index_otherFamily = sum ( constitutive_phenopowerlaw_Ntwin ( 1 : o - 1_pInt , i ) )
do k = 1_pInt , constitutive_phenopowerlaw_Ntwin ( o , i ) ! loop over (active) systems in other family (twin)
2009-10-21 18:40:12 +05:30
constitutive_phenopowerlaw_hardeningMatrix_sliptwin ( index_otherFamily + k , index_myFamily + j , i ) = &
constitutive_phenopowerlaw_interaction_sliptwin ( lattice_interactionSlipTwin ( &
2012-02-13 19:48:07 +05:30
sum ( lattice_NtwinSystem ( 1 : o - 1_pInt , myStructure ) ) + k , &
sum ( lattice_NslipSystem ( 1 : f - 1_pInt , myStructure ) ) + j , &
2009-10-21 18:40:12 +05:30
myStructure ) , i )
enddo ; enddo
enddo ; enddo
2012-02-13 19:48:07 +05:30
do f = 1_pInt , lattice_maxNtwinFamily ! >>> interaction twin -- X
index_myFamily = sum ( constitutive_phenopowerlaw_Ntwin ( 1 : f - 1_pInt , i ) )
do j = 1_pInt , constitutive_phenopowerlaw_Ntwin ( f , i ) ! loop over (active) systems in my family (twin)
2009-10-21 18:40:12 +05:30
2012-02-13 19:48:07 +05:30
do o = 1_pInt , lattice_maxNslipFamily
index_otherFamily = sum ( constitutive_phenopowerlaw_Nslip ( 1 : o - 1_pInt , i ) )
do k = 1_pInt , constitutive_phenopowerlaw_Nslip ( o , i ) ! loop over (active) systems in other family (slip)
2009-10-21 18:40:12 +05:30
constitutive_phenopowerlaw_hardeningMatrix_twinslip ( index_otherFamily + k , index_myFamily + j , i ) = &
constitutive_phenopowerlaw_interaction_twinslip ( lattice_interactionTwinSlip ( &
2012-02-13 19:48:07 +05:30
sum ( lattice_NslipSystem ( 1 : o - 1_pInt , myStructure ) ) + k , &
sum ( lattice_NtwinSystem ( 1 : f - 1_pInt , myStructure ) ) + j , &
2009-10-21 18:40:12 +05:30
myStructure ) , i )
enddo ; enddo
2012-02-13 19:48:07 +05:30
do o = 1_pInt , lattice_maxNtwinFamily
index_otherFamily = sum ( constitutive_phenopowerlaw_Ntwin ( 1 : o - 1_pInt , i ) )
do k = 1_pInt , constitutive_phenopowerlaw_Ntwin ( o , i ) ! loop over (active) systems in other family (twin)
2009-10-21 18:40:12 +05:30
constitutive_phenopowerlaw_hardeningMatrix_twintwin ( index_otherFamily + k , index_myFamily + j , i ) = &
constitutive_phenopowerlaw_interaction_twintwin ( lattice_interactionTwinTwin ( &
2012-02-13 19:48:07 +05:30
sum ( lattice_NtwinSystem ( 1 : o - 1_pInt , myStructure ) ) + k , &
sum ( lattice_NtwinSystem ( 1 : f - 1_pInt , myStructure ) ) + j , &
2009-10-21 18:40:12 +05:30
myStructure ) , i )
enddo ; enddo
enddo ; enddo
2009-07-22 21:37:19 +05:30
2010-11-03 20:28:11 +05:30
! report to out file...
2009-07-22 21:37:19 +05:30
enddo
return
2012-03-09 01:55:28 +05:30
end subroutine constitutive_phenopowerlaw_init
2009-07-22 21:37:19 +05:30
function constitutive_phenopowerlaw_stateInit ( myInstance )
!*********************************************************************
!* initial microstructural state *
!*********************************************************************
use lattice , only : lattice_maxNslipFamily , lattice_maxNtwinFamily
2012-03-09 01:55:28 +05:30
2009-07-22 21:37:19 +05:30
implicit none
integer ( pInt ) , intent ( in ) :: myInstance
2012-03-09 01:55:28 +05:30
integer ( pInt ) :: i
2009-07-23 19:03:53 +05:30
real ( pReal ) , dimension ( constitutive_phenopowerlaw_sizeDotState ( myInstance ) ) :: constitutive_phenopowerlaw_stateInit
2009-07-22 21:37:19 +05:30
constitutive_phenopowerlaw_stateInit = 0.0_pReal
2012-02-21 21:30:00 +05:30
do i = 1_pInt , lattice_maxNslipFamily
2009-07-22 21:37:19 +05:30
constitutive_phenopowerlaw_stateInit ( 1 + &
sum ( constitutive_phenopowerlaw_Nslip ( 1 : i - 1 , myInstance ) ) : &
sum ( constitutive_phenopowerlaw_Nslip ( 1 : i , myInstance ) ) ) = &
constitutive_phenopowerlaw_tau0_slip ( i , myInstance )
enddo
2012-02-21 21:30:00 +05:30
do i = 1_pInt , lattice_maxNtwinFamily
2009-07-22 21:37:19 +05:30
constitutive_phenopowerlaw_stateInit ( 1 + sum ( constitutive_phenopowerlaw_Nslip ( : , myInstance ) ) + &
sum ( constitutive_phenopowerlaw_Ntwin ( 1 : i - 1 , myInstance ) ) : &
sum ( constitutive_phenopowerlaw_Nslip ( : , myInstance ) ) + &
sum ( constitutive_phenopowerlaw_Ntwin ( 1 : i , myInstance ) ) ) = &
constitutive_phenopowerlaw_tau0_twin ( i , myInstance )
enddo
return
2012-03-09 01:55:28 +05:30
end function constitutive_phenopowerlaw_stateInit
2009-07-22 21:37:19 +05:30
2009-09-18 21:07:14 +05:30
!*********************************************************************
2010-10-26 18:46:37 +05:30
!* absolute state tolerance *
2009-09-18 21:07:14 +05:30
!*********************************************************************
2010-10-26 18:46:37 +05:30
pure function constitutive_phenopowerlaw_aTolState ( myInstance )
2009-09-18 21:07:14 +05:30
implicit none
!*** input variables
2012-03-12 20:13:19 +05:30
integer ( pInt ) , intent ( in ) :: myInstance ! number specifying the current instance of the plasticity
2009-09-18 21:07:14 +05:30
!*** output variables
real ( pReal ) , dimension ( constitutive_phenopowerlaw_sizeState ( myInstance ) ) :: &
2012-03-12 20:13:19 +05:30
constitutive_phenopowerlaw_aTolState ! relevant state values for the current instance of this plasticity
2009-09-18 21:07:14 +05:30
!*** local variables
2010-10-26 18:46:37 +05:30
constitutive_phenopowerlaw_aTolState = constitutive_phenopowerlaw_aTolResistance ( myInstance )
2009-09-18 21:07:14 +05:30
2012-03-09 01:55:28 +05:30
end function constitutive_phenopowerlaw_aTolState
2009-09-18 21:07:14 +05:30
2009-07-22 21:37:19 +05:30
function constitutive_phenopowerlaw_homogenizedC ( state , ipc , ip , el )
!*********************************************************************
!* homogenized elacticity matrix *
!* INPUT: *
!* - state : state variables *
!* - ipc : component-ID of current integration point *
!* - ip : current integration point *
!* - el : current element *
!*********************************************************************
2012-03-09 01:55:28 +05:30
use prec , only : p_vec
2009-07-22 21:37:19 +05:30
use mesh , only : mesh_NcpElems , mesh_maxNips
2012-03-12 19:39:37 +05:30
use material , only : homogenization_maxNgrains , material_phase , phase_plasticityInstance
2012-03-09 01:55:28 +05:30
2009-07-22 21:37:19 +05:30
implicit none
integer ( pInt ) , intent ( in ) :: ipc , ip , el
integer ( pInt ) matID
real ( pReal ) , dimension ( 6 , 6 ) :: constitutive_phenopowerlaw_homogenizedC
type ( p_vec ) , dimension ( homogenization_maxNgrains , mesh_maxNips , mesh_NcpElems ) :: state
2012-03-12 19:39:37 +05:30
matID = phase_plasticityInstance ( material_phase ( ipc , ip , el ) )
2009-07-22 21:37:19 +05:30
constitutive_phenopowerlaw_homogenizedC = constitutive_phenopowerlaw_Cslip_66 ( : , : , matID )
return
2012-03-09 01:55:28 +05:30
end function constitutive_phenopowerlaw_homogenizedC
2009-07-22 21:37:19 +05:30
subroutine constitutive_phenopowerlaw_microstructure ( Temperature , state , ipc , ip , el )
!*********************************************************************
!* calculate derived quantities from state (not used here) *
!* INPUT: *
!* - Tp : temperature *
!* - 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
2012-03-12 19:39:37 +05:30
use material , only : homogenization_maxNgrains , material_phase , phase_plasticityInstance
2012-03-09 01:55:28 +05:30
2009-07-22 21:37:19 +05:30
implicit none
2009-08-31 19:43:10 +05:30
integer ( pInt ) ipc , ip , el , matID
2009-07-22 21:37:19 +05:30
real ( pReal ) Temperature
type ( p_vec ) , dimension ( homogenization_maxNgrains , mesh_maxNips , mesh_NcpElems ) :: state
2012-03-12 19:39:37 +05:30
matID = phase_plasticityInstance ( material_phase ( ipc , ip , el ) )
2010-11-03 20:28:11 +05:30
2012-03-09 01:55:28 +05:30
end subroutine constitutive_phenopowerlaw_microstructure
2009-07-22 21:37:19 +05:30
subroutine constitutive_phenopowerlaw_LpAndItsTangent ( Lp , dLp_dTstar , Tstar_v , Temperature , state , ipc , ip , el )
!*********************************************************************
!* plastic velocity gradient and its tangent *
!* INPUT: *
!* - 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) *
!*********************************************************************
2012-03-09 01:55:28 +05:30
use prec , only : p_vec
2009-07-22 21:37:19 +05:30
use math , only : math_Plain3333to99
use lattice , only : lattice_Sslip , lattice_Sslip_v , lattice_Stwin , lattice_Stwin_v , lattice_maxNslipFamily , lattice_maxNtwinFamily , &
2010-10-26 18:46:37 +05:30
lattice_NslipSystem , lattice_NtwinSystem
2009-07-22 21:37:19 +05:30
use mesh , only : mesh_NcpElems , mesh_maxNips
2012-03-12 19:39:37 +05:30
use material , only : homogenization_maxNgrains , material_phase , phase_plasticityInstance
2009-07-22 21:37:19 +05:30
implicit none
integer ( pInt ) ipc , ip , el
integer ( pInt ) matID , nSlip , nTwin , f , i , j , k , l , m , n , structID , index_Gamma , index_F , index_myFamily
real ( pReal ) Temperature
type ( p_vec ) , dimension ( homogenization_maxNgrains , mesh_maxNips , mesh_NcpElems ) :: state
real ( pReal ) , dimension ( 6 ) :: Tstar_v
real ( pReal ) , dimension ( 3 , 3 ) :: Lp
real ( pReal ) , dimension ( 3 , 3 , 3 , 3 ) :: dLp_dTstar3333
real ( pReal ) , dimension ( 9 , 9 ) :: dLp_dTstar
2012-03-12 19:39:37 +05:30
real ( pReal ) , dimension ( constitutive_phenopowerlaw_totalNslip ( phase_plasticityInstance ( material_phase ( ipc , ip , el ) ) ) ) :: &
2009-07-22 21:37:19 +05:30
gdot_slip , dgdot_dtauslip , tau_slip
2012-03-12 19:39:37 +05:30
real ( pReal ) , dimension ( constitutive_phenopowerlaw_totalNtwin ( phase_plasticityInstance ( material_phase ( ipc , ip , el ) ) ) ) :: &
2009-07-22 21:37:19 +05:30
gdot_twin , dgdot_dtautwin , tau_twin
2012-03-12 19:39:37 +05:30
matID = phase_plasticityInstance ( material_phase ( ipc , ip , el ) )
2009-07-22 21:37:19 +05:30
structID = constitutive_phenopowerlaw_structure ( matID )
nSlip = constitutive_phenopowerlaw_totalNslip ( matID )
nTwin = constitutive_phenopowerlaw_totalNtwin ( matID )
2012-02-21 21:30:00 +05:30
index_Gamma = nSlip + nTwin + 1_pInt
index_F = nSlip + nTwin + 2_pInt
2009-07-22 21:37:19 +05:30
Lp = 0.0_pReal
dLp_dTstar3333 = 0.0_pReal
dLp_dTstar = 0.0_pReal
2009-10-21 18:40:12 +05:30
2009-07-22 21:37:19 +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 , structID ) ) ! at which index starts my family
do i = 1_pInt , constitutive_phenopowerlaw_Nslip ( f , matID ) ! process each (active) slip system in family
2009-07-22 21:37:19 +05:30
j = j + 1_pInt
!* Calculation of Lp
openmp parallelization working again (at least for j2 and nonlocal constitutive model).
In order to keep it like that, please follow these simple rules:
DON'T use implicit array subscripts:
example: real, dimension(3,3) :: A,B
A(:,2) = B(:,1) <--- DON'T USE
A(1:3,2) = B(1:3,1) <--- BETTER USE
In many cases the use of explicit array subscripts is inevitable for parallelization. Additionally, it is an easy means to prevent memory leaks.
Enclose all write statements with the following:
!$OMP CRITICAL (write2out)
<your write statement>
!$OMP END CRITICAL (write2out)
Whenever you change something in the code and are not sure if it affects parallelization and leads to nonconforming behavior, please ask me and/or Franz to check this.
2011-03-17 16:16:17 +05:30
tau_slip ( j ) = dot_product ( Tstar_v , lattice_Sslip_v ( 1 : 6 , index_myFamily + i , structID ) )
2009-07-22 21:37:19 +05:30
gdot_slip ( j ) = constitutive_phenopowerlaw_gdot0_slip ( matID ) * ( abs ( tau_slip ( j ) ) / state ( ipc , ip , el ) % p ( j ) ) ** &
constitutive_phenopowerlaw_n_slip ( matID ) * sign ( 1.0_pReal , tau_slip ( j ) )
2009-10-22 14:28:14 +05:30
Lp = Lp + ( 1.0_pReal - state ( ipc , ip , el ) % p ( index_F ) ) * & ! 1-F
openmp parallelization working again (at least for j2 and nonlocal constitutive model).
In order to keep it like that, please follow these simple rules:
DON'T use implicit array subscripts:
example: real, dimension(3,3) :: A,B
A(:,2) = B(:,1) <--- DON'T USE
A(1:3,2) = B(1:3,1) <--- BETTER USE
In many cases the use of explicit array subscripts is inevitable for parallelization. Additionally, it is an easy means to prevent memory leaks.
Enclose all write statements with the following:
!$OMP CRITICAL (write2out)
<your write statement>
!$OMP END CRITICAL (write2out)
Whenever you change something in the code and are not sure if it affects parallelization and leads to nonconforming behavior, please ask me and/or Franz to check this.
2011-03-17 16:16:17 +05:30
gdot_slip ( j ) * lattice_Sslip ( 1 : 3 , 1 : 3 , index_myFamily + i , structID )
2009-07-22 21:37:19 +05:30
!* Calculation of the tangent of Lp
2009-10-21 18:40:12 +05:30
if ( gdot_slip ( j ) / = 0.0_pReal ) then
dgdot_dtauslip ( j ) = gdot_slip ( j ) * constitutive_phenopowerlaw_n_slip ( matID ) / tau_slip ( j )
2012-02-21 21:30:00 +05:30
forall ( k = 1_pInt : 3_pInt , l = 1_pInt : 3_pInt , m = 1_pInt : 3_pInt , n = 1_pInt : 3_pInt ) &
2009-10-21 18:40:12 +05:30
dLp_dTstar3333 ( k , l , m , n ) = dLp_dTstar3333 ( k , l , m , n ) + &
dgdot_dtauslip ( j ) * lattice_Sslip ( k , l , index_myFamily + i , structID ) * &
lattice_Sslip ( m , n , index_myFamily + i , structID )
endif
2009-07-22 21:37:19 +05:30
enddo
enddo
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 , structID ) ) ! at which index starts my family
do i = 1_pInt , constitutive_phenopowerlaw_Ntwin ( f , matID ) ! process each (active) twin system in family
2009-07-22 21:37:19 +05:30
j = j + 1_pInt
!* Calculation of Lp
openmp parallelization working again (at least for j2 and nonlocal constitutive model).
In order to keep it like that, please follow these simple rules:
DON'T use implicit array subscripts:
example: real, dimension(3,3) :: A,B
A(:,2) = B(:,1) <--- DON'T USE
A(1:3,2) = B(1:3,1) <--- BETTER USE
In many cases the use of explicit array subscripts is inevitable for parallelization. Additionally, it is an easy means to prevent memory leaks.
Enclose all write statements with the following:
!$OMP CRITICAL (write2out)
<your write statement>
!$OMP END CRITICAL (write2out)
Whenever you change something in the code and are not sure if it affects parallelization and leads to nonconforming behavior, please ask me and/or Franz to check this.
2011-03-17 16:16:17 +05:30
tau_twin ( j ) = dot_product ( Tstar_v , lattice_Stwin_v ( 1 : 6 , index_myFamily + i , structID ) )
2009-10-22 14:28:14 +05:30
gdot_twin ( j ) = ( 1.0_pReal - state ( ipc , ip , el ) % p ( index_F ) ) * & ! 1-F
constitutive_phenopowerlaw_gdot0_twin ( matID ) * &
2009-07-22 21:37:19 +05:30
( abs ( tau_twin ( j ) ) / state ( ipc , ip , el ) % p ( nSlip + j ) ) ** &
constitutive_phenopowerlaw_n_twin ( matID ) * max ( 0.0_pReal , sign ( 1.0_pReal , tau_twin ( j ) ) )
openmp parallelization working again (at least for j2 and nonlocal constitutive model).
In order to keep it like that, please follow these simple rules:
DON'T use implicit array subscripts:
example: real, dimension(3,3) :: A,B
A(:,2) = B(:,1) <--- DON'T USE
A(1:3,2) = B(1:3,1) <--- BETTER USE
In many cases the use of explicit array subscripts is inevitable for parallelization. Additionally, it is an easy means to prevent memory leaks.
Enclose all write statements with the following:
!$OMP CRITICAL (write2out)
<your write statement>
!$OMP END CRITICAL (write2out)
Whenever you change something in the code and are not sure if it affects parallelization and leads to nonconforming behavior, please ask me and/or Franz to check this.
2011-03-17 16:16:17 +05:30
Lp = Lp + gdot_twin ( j ) * lattice_Stwin ( 1 : 3 , 1 : 3 , index_myFamily + i , structID )
2009-07-22 21:37:19 +05:30
!* Calculation of the tangent of Lp
2009-10-21 18:40:12 +05:30
if ( gdot_twin ( j ) / = 0.0_pReal ) then
dgdot_dtautwin ( j ) = gdot_twin ( j ) * constitutive_phenopowerlaw_n_twin ( matID ) / tau_twin ( j )
2012-02-21 21:30:00 +05:30
forall ( k = 1_pInt : 3_pInt , l = 1_pInt : 3_pInt , m = 1_pInt : 3_pInt , n = 1_pInt : 3_pInt ) &
2009-10-21 18:40:12 +05:30
dLp_dTstar3333 ( k , l , m , n ) = dLp_dTstar3333 ( k , l , m , n ) + &
dgdot_dtautwin ( j ) * lattice_Stwin ( k , l , index_myFamily + i , structID ) * &
lattice_Stwin ( m , n , index_myFamily + i , structID )
endif
2009-07-22 21:37:19 +05:30
enddo
enddo
dLp_dTstar = math_Plain3333to99 ( dLp_dTstar3333 )
return
2012-03-09 01:55:28 +05:30
end subroutine constitutive_phenopowerlaw_LpAndItsTangent
2009-07-22 21:37:19 +05:30
function constitutive_phenopowerlaw_dotState ( Tstar_v , Temperature , state , ipc , ip , el )
!*********************************************************************
!* rate of change of microstructure *
!* INPUT: *
!* - 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 *
!*********************************************************************
2012-03-09 01:55:28 +05:30
use prec , only : p_vec
2012-02-21 21:30:00 +05:30
use lattice , only : lattice_Sslip_v , lattice_Stwin_v , lattice_maxNslipFamily , lattice_maxNtwinFamily , &
2009-07-22 21:37:19 +05:30
lattice_NslipSystem , lattice_NtwinSystem , lattice_shearTwin
use mesh , only : mesh_NcpElems , mesh_maxNips
2012-03-12 19:39:37 +05:30
use material , only : homogenization_maxNgrains , material_phase , phase_plasticityInstance
2012-03-09 01:55:28 +05:30
2009-07-22 21:37:19 +05:30
implicit none
integer ( pInt ) ipc , ip , el
2011-04-13 19:46:22 +05:30
integer ( pInt ) matID , nSlip , nTwin , f , i , j , structID , index_Gamma , index_F , index_myFamily
2009-10-22 14:28:14 +05:30
real ( pReal ) Temperature , c_slipslip , c_sliptwin , c_twinslip , c_twintwin , ssat_offset
2009-07-22 21:37:19 +05:30
type ( p_vec ) , dimension ( homogenization_maxNgrains , mesh_maxNips , mesh_NcpElems ) :: state
real ( pReal ) , dimension ( 6 ) :: Tstar_v
2012-03-12 19:39:37 +05:30
real ( pReal ) , dimension ( constitutive_phenopowerlaw_totalNslip ( phase_plasticityInstance ( material_phase ( ipc , ip , el ) ) ) ) :: &
2009-10-21 18:40:12 +05:30
gdot_slip , tau_slip , h_slipslip , h_sliptwin
2012-03-12 19:39:37 +05:30
real ( pReal ) , dimension ( constitutive_phenopowerlaw_totalNtwin ( phase_plasticityInstance ( material_phase ( ipc , ip , el ) ) ) ) :: &
2009-10-21 18:40:12 +05:30
gdot_twin , tau_twin , h_twinslip , h_twintwin
2012-03-12 19:39:37 +05:30
real ( pReal ) , dimension ( constitutive_phenopowerlaw_sizeDotState ( phase_plasticityInstance ( material_phase ( ipc , ip , el ) ) ) ) :: &
2009-07-22 21:37:19 +05:30
constitutive_phenopowerlaw_dotState
2012-03-12 19:39:37 +05:30
matID = phase_plasticityInstance ( material_phase ( ipc , ip , el ) )
2009-07-22 21:37:19 +05:30
structID = constitutive_phenopowerlaw_structure ( matID )
nSlip = constitutive_phenopowerlaw_totalNslip ( matID )
nTwin = constitutive_phenopowerlaw_totalNtwin ( matID )
2012-02-21 21:30:00 +05:30
index_Gamma = nSlip + nTwin + 1_pInt
index_F = nSlip + nTwin + 2_pInt
2009-07-22 21:37:19 +05:30
constitutive_phenopowerlaw_dotState = 0.0_pReal
!-- system-independent (nonlinear) prefactors to M_xx matrices
c_slipslip = constitutive_phenopowerlaw_h0_slipslip ( matID ) * &
( 1.0_pReal + &
constitutive_phenopowerlaw_twinC ( matID ) * state ( ipc , ip , el ) % p ( index_F ) ** constitutive_phenopowerlaw_twinB ( matID ) )
c_sliptwin = 0.0_pReal
c_twinslip = constitutive_phenopowerlaw_h0_twinslip ( matID ) * &
state ( ipc , ip , el ) % p ( index_Gamma ) ** constitutive_phenopowerlaw_twinE ( matID )
c_twintwin = constitutive_phenopowerlaw_h0_twintwin ( matID ) * &
state ( ipc , ip , el ) % p ( index_F ) ** constitutive_phenopowerlaw_twinD ( matID )
!-- add system-dependent part and calculate dot gammas
2011-02-25 14:55:53 +05:30
ssat_offset = constitutive_phenopowerlaw_spr ( matID ) * sqrt ( state ( ipc , ip , el ) % p ( index_F ) )
2009-07-22 21:37:19 +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 , structID ) ) ! at which index starts my family
do i = 1_pInt , constitutive_phenopowerlaw_Nslip ( f , matID ) ! process each (active) slip system in family
2009-07-22 21:37:19 +05:30
j = j + 1_pInt
2009-10-22 14:28:14 +05:30
h_slipslip ( j ) = c_slipslip * ( 1.0_pReal - state ( ipc , ip , el ) % p ( j ) / & ! system-dependent prefactor for slip--slip interaction
( constitutive_phenopowerlaw_tausat_slip ( f , matID ) + ssat_offset ) ) ** &
2011-11-23 20:18:39 +05:30
constitutive_phenopowerlaw_a_slip ( matID )
2009-07-22 21:37:19 +05:30
h_sliptwin ( j ) = c_sliptwin ! no system-dependent part
!* Calculation of dot gamma
openmp parallelization working again (at least for j2 and nonlocal constitutive model).
In order to keep it like that, please follow these simple rules:
DON'T use implicit array subscripts:
example: real, dimension(3,3) :: A,B
A(:,2) = B(:,1) <--- DON'T USE
A(1:3,2) = B(1:3,1) <--- BETTER USE
In many cases the use of explicit array subscripts is inevitable for parallelization. Additionally, it is an easy means to prevent memory leaks.
Enclose all write statements with the following:
!$OMP CRITICAL (write2out)
<your write statement>
!$OMP END CRITICAL (write2out)
Whenever you change something in the code and are not sure if it affects parallelization and leads to nonconforming behavior, please ask me and/or Franz to check this.
2011-03-17 16:16:17 +05:30
tau_slip ( j ) = dot_product ( Tstar_v , lattice_Sslip_v ( 1 : 6 , index_myFamily + i , structID ) )
2009-07-22 21:37:19 +05:30
gdot_slip ( j ) = constitutive_phenopowerlaw_gdot0_slip ( matID ) * ( abs ( tau_slip ( j ) ) / state ( ipc , ip , el ) % p ( j ) ) ** &
constitutive_phenopowerlaw_n_slip ( matID ) * sign ( 1.0_pReal , tau_slip ( j ) )
enddo
enddo
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 , structID ) ) ! at which index starts my family
do i = 1_pInt , constitutive_phenopowerlaw_Ntwin ( f , matID ) ! process each (active) twin system in family
2009-07-22 21:37:19 +05:30
j = j + 1_pInt
h_twinslip ( j ) = c_twinslip ! no system-dependent parts
h_twintwin ( j ) = c_twintwin
!* Calculation of dot vol frac
openmp parallelization working again (at least for j2 and nonlocal constitutive model).
In order to keep it like that, please follow these simple rules:
DON'T use implicit array subscripts:
example: real, dimension(3,3) :: A,B
A(:,2) = B(:,1) <--- DON'T USE
A(1:3,2) = B(1:3,1) <--- BETTER USE
In many cases the use of explicit array subscripts is inevitable for parallelization. Additionally, it is an easy means to prevent memory leaks.
Enclose all write statements with the following:
!$OMP CRITICAL (write2out)
<your write statement>
!$OMP END CRITICAL (write2out)
Whenever you change something in the code and are not sure if it affects parallelization and leads to nonconforming behavior, please ask me and/or Franz to check this.
2011-03-17 16:16:17 +05:30
tau_twin ( j ) = dot_product ( Tstar_v , lattice_Stwin_v ( 1 : 6 , index_myFamily + i , structID ) )
2009-10-22 14:28:14 +05:30
gdot_twin ( j ) = ( 1.0_pReal - state ( ipc , ip , el ) % p ( index_F ) ) * & ! 1-F
constitutive_phenopowerlaw_gdot0_twin ( matID ) * &
2009-07-22 21:37:19 +05:30
( abs ( tau_twin ( j ) ) / state ( ipc , ip , el ) % p ( nSlip + j ) ) ** &
constitutive_phenopowerlaw_n_twin ( matID ) * max ( 0.0_pReal , sign ( 1.0_pReal , tau_twin ( j ) ) )
enddo
enddo
!-- calculate the overall hardening based on above
j = 0_pInt
2012-02-21 21:30:00 +05:30
do f = 1_pInt , lattice_maxNslipFamily ! loop over all slip families
do i = 1_pInt , constitutive_phenopowerlaw_Nslip ( f , matID ) ! process each (active) slip system in family
2009-07-22 21:37:19 +05:30
j = j + 1_pInt
openmp parallelization working again (at least for j2 and nonlocal constitutive model).
In order to keep it like that, please follow these simple rules:
DON'T use implicit array subscripts:
example: real, dimension(3,3) :: A,B
A(:,2) = B(:,1) <--- DON'T USE
A(1:3,2) = B(1:3,1) <--- BETTER USE
In many cases the use of explicit array subscripts is inevitable for parallelization. Additionally, it is an easy means to prevent memory leaks.
Enclose all write statements with the following:
!$OMP CRITICAL (write2out)
<your write statement>
!$OMP END CRITICAL (write2out)
Whenever you change something in the code and are not sure if it affects parallelization and leads to nonconforming behavior, please ask me and/or Franz to check this.
2011-03-17 16:16:17 +05:30
constitutive_phenopowerlaw_dotState ( j ) = & ! evolution of slip resistance j
h_slipslip ( j ) * dot_product ( constitutive_phenopowerlaw_hardeningMatrix_slipslip ( 1 : nSlip , j , matID ) , abs ( gdot_slip ) ) + & ! dot gamma_slip
h_sliptwin ( j ) * dot_product ( constitutive_phenopowerlaw_hardeningMatrix_sliptwin ( 1 : nTwin , j , matID ) , gdot_twin ) ! dot gamma_twin
2009-07-22 21:37:19 +05:30
constitutive_phenopowerlaw_dotState ( index_Gamma ) = constitutive_phenopowerlaw_dotState ( index_Gamma ) + &
abs ( gdot_slip ( j ) )
enddo
enddo
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 , structID ) ) ! at which index starts my family
do i = 1_pInt , constitutive_phenopowerlaw_Ntwin ( f , matID ) ! process each (active) twin system in family
2009-07-22 21:37:19 +05:30
j = j + 1_pInt
openmp parallelization working again (at least for j2 and nonlocal constitutive model).
In order to keep it like that, please follow these simple rules:
DON'T use implicit array subscripts:
example: real, dimension(3,3) :: A,B
A(:,2) = B(:,1) <--- DON'T USE
A(1:3,2) = B(1:3,1) <--- BETTER USE
In many cases the use of explicit array subscripts is inevitable for parallelization. Additionally, it is an easy means to prevent memory leaks.
Enclose all write statements with the following:
!$OMP CRITICAL (write2out)
<your write statement>
!$OMP END CRITICAL (write2out)
Whenever you change something in the code and are not sure if it affects parallelization and leads to nonconforming behavior, please ask me and/or Franz to check this.
2011-03-17 16:16:17 +05:30
constitutive_phenopowerlaw_dotState ( j + nSlip ) = & ! evolution of twin resistance j
h_twinslip ( j ) * dot_product ( constitutive_phenopowerlaw_hardeningMatrix_twinslip ( 1 : nSlip , j , matID ) , abs ( gdot_slip ) ) + & ! dot gamma_slip
h_twintwin ( j ) * dot_product ( constitutive_phenopowerlaw_hardeningMatrix_twintwin ( 1 : nTwin , j , matID ) , gdot_twin ) ! dot gamma_twin
2009-07-22 21:37:19 +05:30
constitutive_phenopowerlaw_dotState ( index_F ) = constitutive_phenopowerlaw_dotState ( index_F ) + &
gdot_twin ( j ) / lattice_shearTwin ( index_myFamily + i , structID )
enddo
enddo
2012-03-09 01:55:28 +05:30
end function constitutive_phenopowerlaw_dotState
2009-07-22 21:37:19 +05:30
2012-05-16 20:13:26 +05:30
!*********************************************************************
!* (instantaneous) incremental change of microstructure *
!*********************************************************************
function constitutive_phenopowerlaw_deltaState ( Tstar_v , Temperature , state , g , ip , el )
use prec , only : pReal , &
pInt , &
p_vec
use mesh , only : mesh_NcpElems , &
mesh_maxNips
use material , only : homogenization_maxNgrains , &
material_phase , &
phase_plasticityInstance
implicit none
!*** input variables
integer ( pInt ) , intent ( in ) :: g , & ! current grain number
ip , & ! current integration point
el ! current element number
real ( pReal ) , intent ( in ) :: Temperature ! temperature
real ( pReal ) , dimension ( 6 ) , intent ( in ) :: Tstar_v ! current 2nd Piola-Kirchhoff stress in Mandel notation
type ( p_vec ) , dimension ( homogenization_maxNgrains , mesh_maxNips , mesh_NcpElems ) , intent ( in ) :: &
state ! current microstructural state
!*** output variables
real ( pReal ) , dimension ( constitutive_phenopowerlaw_sizeDotState ( phase_plasticityInstance ( material_phase ( g , ip , el ) ) ) ) :: &
constitutive_phenopowerlaw_deltaState ! change of state variables / microstructure
!*** local variables
constitutive_phenopowerlaw_deltaState = 0.0_pReal
endfunction
2009-07-22 21:37:19 +05:30
!****************************************************************
!* calculates the rate of change of temperature *
!****************************************************************
pure function constitutive_phenopowerlaw_dotTemperature ( Tstar_v , Temperature , state , ipc , ip , el )
!*** variables and functions from other modules ***!
use prec , only : pReal , pInt , p_vec
2012-02-21 21:30:00 +05:30
use mesh , only : mesh_NcpElems , mesh_maxNips
use material , only : homogenization_maxNgrains
2012-03-09 01:55:28 +05:30
2009-07-22 21:37:19 +05:30
implicit none
!*** input variables ***!
real ( pReal ) , dimension ( 6 ) , intent ( in ) :: Tstar_v ! 2nd Piola Kirchhoff stress tensor in Mandel notation
real ( pReal ) , intent ( in ) :: Temperature
integer ( pInt ) , intent ( in ) :: ipc , & ! grain number
ip , & ! integration point number
el ! element number
type ( p_vec ) , dimension ( homogenization_maxNgrains , mesh_maxNips , mesh_NcpElems ) , intent ( in ) :: state ! state of the current microstructure
!*** output variables ***!
real ( pReal ) constitutive_phenopowerlaw_dotTemperature ! rate of change of temparature
! calculate dotTemperature
constitutive_phenopowerlaw_dotTemperature = 0.0_pReal
2012-03-09 01:55:28 +05:30
end function constitutive_phenopowerlaw_dotTemperature
2009-07-22 21:37:19 +05:30
pure function constitutive_phenopowerlaw_postResults ( Tstar_v , Temperature , dt , state , ipc , ip , el )
!*********************************************************************
!* return array of constitutive results *
!* INPUT: *
!* - 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 lattice , only : lattice_Sslip_v , lattice_Stwin_v , lattice_maxNslipFamily , lattice_maxNtwinFamily , &
2009-10-22 14:28:14 +05:30
lattice_NslipSystem , lattice_NtwinSystem
2009-07-22 21:37:19 +05:30
use mesh , only : mesh_NcpElems , mesh_maxNips
2012-03-12 19:39:37 +05:30
use material , only : homogenization_maxNgrains , material_phase , phase_plasticityInstance , phase_Noutput
2012-03-09 01:55:28 +05:30
2009-07-22 21:37:19 +05:30
implicit none
integer ( pInt ) , intent ( in ) :: ipc , 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
integer ( pInt ) matID , o , f , i , c , nSlip , nTwin , j , structID , index_Gamma , index_F , index_myFamily
real ( pReal ) tau
2012-03-12 19:39:37 +05:30
real ( pReal ) , dimension ( constitutive_phenopowerlaw_sizePostResults ( phase_plasticityInstance ( material_phase ( ipc , ip , el ) ) ) ) :: &
2009-07-22 21:37:19 +05:30
constitutive_phenopowerlaw_postResults
2012-03-12 19:39:37 +05:30
matID = phase_plasticityInstance ( material_phase ( ipc , ip , el ) )
2009-07-22 21:37:19 +05:30
structID = constitutive_phenopowerlaw_structure ( matID )
nSlip = constitutive_phenopowerlaw_totalNslip ( matID )
nTwin = constitutive_phenopowerlaw_totalNtwin ( matID )
2012-02-21 21:30:00 +05:30
index_Gamma = nSlip + nTwin + 1_pInt
index_F = nSlip + nTwin + 2_pInt
2009-07-22 21:37:19 +05:30
constitutive_phenopowerlaw_postResults = 0.0_pReal
c = 0_pInt
2012-02-21 21:30:00 +05:30
do o = 1_pInt , phase_Noutput ( material_phase ( ipc , ip , el ) )
2009-07-22 21:37:19 +05:30
select case ( constitutive_phenopowerlaw_output ( o , matID ) )
case ( 'resistance_slip' )
2012-02-21 21:30:00 +05:30
constitutive_phenopowerlaw_postResults ( c + 1_pInt : c + nSlip ) = state ( ipc , ip , el ) % p ( 1 : nSlip )
2009-07-22 21:37:19 +05:30
c = c + nSlip
case ( 'shearrate_slip' )
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 , structID ) ) ! at which index starts my family
do i = 1_pInt , constitutive_phenopowerlaw_Nslip ( f , matID ) ! process each (active) slip system in family
2009-07-22 21:37:19 +05:30
j = j + 1_pInt
tau = dot_product ( Tstar_v , lattice_Sslip_v ( : , index_myFamily + i , structID ) )
constitutive_phenopowerlaw_postResults ( c + j ) = constitutive_phenopowerlaw_gdot0_slip ( matID ) * &
( abs ( tau ) / state ( ipc , ip , el ) % p ( j ) ) ** &
constitutive_phenopowerlaw_n_slip ( matID ) * sign ( 1.0_pReal , tau )
enddo ; enddo
c = c + nSlip
case ( 'resolvedstress_slip' )
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 , structID ) ) ! at which index starts my family
do i = 1_pInt , constitutive_phenopowerlaw_Nslip ( f , matID ) ! process each (active) slip system in family
2009-07-22 21:37:19 +05:30
j = j + 1_pInt
openmp parallelization working again (at least for j2 and nonlocal constitutive model).
In order to keep it like that, please follow these simple rules:
DON'T use implicit array subscripts:
example: real, dimension(3,3) :: A,B
A(:,2) = B(:,1) <--- DON'T USE
A(1:3,2) = B(1:3,1) <--- BETTER USE
In many cases the use of explicit array subscripts is inevitable for parallelization. Additionally, it is an easy means to prevent memory leaks.
Enclose all write statements with the following:
!$OMP CRITICAL (write2out)
<your write statement>
!$OMP END CRITICAL (write2out)
Whenever you change something in the code and are not sure if it affects parallelization and leads to nonconforming behavior, please ask me and/or Franz to check this.
2011-03-17 16:16:17 +05:30
constitutive_phenopowerlaw_postResults ( c + j ) = dot_product ( Tstar_v , lattice_Sslip_v ( 1 : 6 , index_myFamily + i , structID ) )
2009-10-21 18:40:12 +05:30
enddo ; enddo
2009-07-22 21:37:19 +05:30
c = c + nSlip
case ( 'totalshear' )
2012-02-21 21:30:00 +05:30
constitutive_phenopowerlaw_postResults ( c + 1_pInt ) = state ( ipc , ip , el ) % p ( index_Gamma )
c = c + 1_pInt
2009-07-22 21:37:19 +05:30
case ( 'resistance_twin' )
2012-02-21 21:30:00 +05:30
constitutive_phenopowerlaw_postResults ( c + 1_pInt : c + nTwin ) = state ( ipc , ip , el ) % p ( 1_pInt + nSlip : nTwin + nSlip )
2009-07-22 21:37:19 +05:30
c = c + nTwin
2009-10-21 18:40:12 +05:30
case ( 'shearrate_twin' )
2009-07-22 21:37:19 +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 , structID ) ) ! at which index starts my family
do i = 1_pInt , constitutive_phenopowerlaw_Ntwin ( f , matID ) ! process each (active) twin system in family
2009-07-22 21:37:19 +05:30
j = j + 1_pInt
openmp parallelization working again (at least for j2 and nonlocal constitutive model).
In order to keep it like that, please follow these simple rules:
DON'T use implicit array subscripts:
example: real, dimension(3,3) :: A,B
A(:,2) = B(:,1) <--- DON'T USE
A(1:3,2) = B(1:3,1) <--- BETTER USE
In many cases the use of explicit array subscripts is inevitable for parallelization. Additionally, it is an easy means to prevent memory leaks.
Enclose all write statements with the following:
!$OMP CRITICAL (write2out)
<your write statement>
!$OMP END CRITICAL (write2out)
Whenever you change something in the code and are not sure if it affects parallelization and leads to nonconforming behavior, please ask me and/or Franz to check this.
2011-03-17 16:16:17 +05:30
tau = dot_product ( Tstar_v , lattice_Stwin_v ( 1 : 6 , index_myFamily + i , structID ) )
2009-10-22 14:28:14 +05:30
constitutive_phenopowerlaw_postResults ( c + j ) = ( 1.0_pReal - state ( ipc , ip , el ) % p ( index_F ) ) * & ! 1-F
constitutive_phenopowerlaw_gdot0_twin ( matID ) * &
2009-07-22 21:37:19 +05:30
( abs ( tau ) / state ( ipc , ip , el ) % p ( j + nSlip ) ) ** &
constitutive_phenopowerlaw_n_twin ( matID ) * max ( 0.0_pReal , sign ( 1.0_pReal , tau ) )
enddo ; enddo
c = c + nTwin
case ( 'resolvedstress_twin' )
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 , structID ) ) ! at which index starts my family
do i = 1_pInt , constitutive_phenopowerlaw_Ntwin ( f , matID ) ! process each (active) twin system in family
2009-07-22 21:37:19 +05:30
j = j + 1_pInt
openmp parallelization working again (at least for j2 and nonlocal constitutive model).
In order to keep it like that, please follow these simple rules:
DON'T use implicit array subscripts:
example: real, dimension(3,3) :: A,B
A(:,2) = B(:,1) <--- DON'T USE
A(1:3,2) = B(1:3,1) <--- BETTER USE
In many cases the use of explicit array subscripts is inevitable for parallelization. Additionally, it is an easy means to prevent memory leaks.
Enclose all write statements with the following:
!$OMP CRITICAL (write2out)
<your write statement>
!$OMP END CRITICAL (write2out)
Whenever you change something in the code and are not sure if it affects parallelization and leads to nonconforming behavior, please ask me and/or Franz to check this.
2011-03-17 16:16:17 +05:30
constitutive_phenopowerlaw_postResults ( c + j ) = dot_product ( Tstar_v , lattice_Stwin_v ( 1 : 6 , index_myFamily + i , structID ) )
2009-10-21 18:40:12 +05:30
enddo ; enddo
2009-07-22 21:37:19 +05:30
c = c + nTwin
case ( 'totalvolfrac' )
2012-02-21 21:30:00 +05:30
constitutive_phenopowerlaw_postResults ( c + 1_pInt ) = state ( ipc , ip , el ) % p ( index_F )
c = c + 1_pInt
2009-07-22 21:37:19 +05:30
end select
enddo
2012-03-09 01:55:28 +05:30
end function constitutive_phenopowerlaw_postResults
2009-07-22 21:37:19 +05:30
2012-03-09 01:55:28 +05:30
end module constitutive_phenopowerlaw