2012-03-09 01:55:28 +05:30
!--------------------------------------------------------------------------------------------------
2014-03-29 13:50:36 +05:30
! $Id$
2012-03-09 01:55:28 +05:30
!--------------------------------------------------------------------------------------------------
!> @author Franz Roters, Max-Planck-Institut für Eisenforschung GmbH
2013-06-11 22:05:04 +05:30
!> @author Philip Eisenlohr, Max-Planck-Institut für Eisenforschung GmbH
2013-01-18 17:00:52 +05:30
!> @brief Parses material config file, either solverJobName.materialConfig or material.config
!> @details reads the material configuration file, where solverJobName.materialConfig takes
!! precedence over material.config and parses the sections 'homogenization', 'crystallite',
!! 'phase', 'texture', and 'microstucture'
2012-03-09 01:55:28 +05:30
!--------------------------------------------------------------------------------------------------
module material
2013-01-18 17:00:52 +05:30
use prec , only : &
pReal , &
2013-02-25 18:43:52 +05:30
pInt , &
2014-04-15 14:50:38 +05:30
#ifdef NEWSTATE
tState , &
#endif
2013-02-25 18:43:52 +05:30
p_intvec
2012-03-09 01:55:28 +05:30
implicit none
private
2013-12-16 16:26:56 +05:30
character ( len = * ) , parameter , public :: &
2013-11-27 13:35:23 +05:30
ELASTICITY_HOOKE_label = 'hooke' , &
PLASTICITY_NONE_label = 'none' , &
PLASTICITY_J2_label = 'j2' , &
PLASTICITY_PHENOPOWERLAW_label = 'phenopowerlaw' , &
PLASTICITY_DISLOTWIN_label = 'dislotwin' , &
PLASTICITY_TITANMOD_label = 'titanmod' , &
PLASTICITY_NONLOCAL_label = 'nonlocal' , &
2014-06-23 00:28:29 +05:30
DAMAGE_NONE_label = 'none' , &
DAMAGE_LOCAL_label = 'local' , &
DAMAGE_GRADIENT_label = 'gradient' , &
THERMAL_NONE_label = 'none' , &
THERMAL_ISO_label = 'isothermal' , &
THERMAL_CONDUCTION_label = 'conduction' , &
THERMAL_ADIABATIC_label = 'adiabatic' , &
2014-03-14 04:50:50 +05:30
HOMOGENIZATION_NONE_label = 'none' , &
2013-11-27 13:35:23 +05:30
HOMOGENIZATION_ISOSTRAIN_label = 'isostrain' , &
HOMOGENIZATION_RGC_label = 'rgc'
enum , bind ( c )
2013-12-12 22:39:59 +05:30
enumerator :: ELASTICITY_undefined_ID , &
ELASTICITY_hooke_ID
end enum
enum , bind ( c )
enumerator :: PLASTICITY_undefined_ID , &
PLASTICITY_none_ID , &
2013-11-27 13:35:23 +05:30
PLASTICITY_J2_ID , &
PLASTICITY_phenopowerlaw_ID , &
PLASTICITY_dislotwin_ID , &
PLASTICITY_titanmod_ID , &
PLASTICITY_nonlocal_ID
2013-12-12 22:39:59 +05:30
end enum
2014-06-23 00:28:29 +05:30
enum , bind ( c )
enumerator :: DAMAGE_undefined_ID , &
DAMAGE_none_ID , &
DAMAGE_local_ID , &
DAMAGE_gradient_ID
end enum
enum , bind ( c )
enumerator :: THERMAL_undefined_ID , &
THERMAL_none_ID , &
THERMAL_iso_ID , &
THERMAL_conduction_ID , &
THERMAL_adiabatic_ID
end enum
2013-12-12 22:39:59 +05:30
enum , bind ( c )
enumerator :: HOMOGENIZATION_undefined_ID , &
2014-03-14 04:50:50 +05:30
HOMOGENIZATION_none_ID , &
2013-12-12 22:39:59 +05:30
HOMOGENIZATION_isostrain_ID , &
2013-11-27 13:35:23 +05:30
HOMOGENIZATION_RGC_ID
end enum
2013-12-16 16:26:56 +05:30
character ( len = * ) , parameter , public :: &
2013-06-27 00:49:00 +05:30
MATERIAL_configFile = 'material.config' , & !< generic name for material configuration file
MATERIAL_localFileExt = 'materialConfig' !< extension of solver job name depending material configuration file
2012-03-09 01:55:28 +05:30
2013-12-16 16:26:56 +05:30
character ( len = * ) , parameter , public :: &
2013-06-27 00:49:00 +05:30
MATERIAL_partHomogenization = 'homogenization' , & !< keyword for homogenization part
MATERIAL_partCrystallite = 'crystallite' , & !< keyword for crystallite part
MATERIAL_partPhase = 'phase' !< keyword for phase part
2013-11-27 13:35:23 +05:30
2013-12-13 04:33:37 +05:30
integer ( kind ( ELASTICITY_undefined_ID ) ) , dimension ( : ) , allocatable , public , protected :: &
2013-12-13 19:44:17 +05:30
phase_elasticity !< elasticity of each phase
integer ( kind ( PLASTICITY_undefined_ID ) ) , dimension ( : ) , allocatable , public , protected :: &
phase_plasticity !< plasticity of each phase
2014-06-23 00:28:29 +05:30
integer ( kind ( DAMAGE_undefined_ID ) ) , dimension ( : ) , allocatable , public , protected :: &
phase_damage !< damage of each phase
integer ( kind ( THERMAL_undefined_ID ) ) , dimension ( : ) , allocatable , public , protected :: &
phase_thermal !< thermal of each phase
2013-12-13 19:44:17 +05:30
integer ( kind ( HOMOGENIZATION_undefined_ID ) ) , dimension ( : ) , allocatable , public , protected :: &
2013-11-27 13:35:23 +05:30
homogenization_type !< type of each homogenization
2013-12-16 16:26:56 +05:30
character ( len = 64 ) , dimension ( : ) , allocatable , public , protected :: &
2013-02-20 03:42:05 +05:30
phase_name , & !< name of each phase
homogenization_name , & !< name of each homogenization
crystallite_name !< name of each crystallite setting
2012-03-09 01:55:28 +05:30
2013-01-18 17:00:52 +05:30
integer ( pInt ) , public , protected :: &
2013-02-20 03:42:05 +05:30
homogenization_maxNgrains , & !< max number of grains in any USED homogenization
material_Nphase , & !< number of phases
material_Nhomogenization , & !< number of homogenizations
material_Nmicrostructure , & !< number of microstructures
material_Ncrystallite !< number of crystallite settings
2012-03-09 01:55:28 +05:30
2013-12-16 16:26:56 +05:30
integer ( pInt ) , dimension ( : ) , allocatable , public , protected :: &
2013-02-20 03:42:05 +05:30
homogenization_Ngrains , & !< number of grains in each homogenization
homogenization_Noutput , & !< number of '(output)' items per homogenization
phase_Noutput , & !< number of '(output)' items per phase
phase_elasticityInstance , & !< instance of particular elasticity of each phase
phase_plasticityInstance , & !< instance of particular plasticity of each phase
2014-06-23 00:28:29 +05:30
phase_damageInstance , & !< instance of particular plasticity of each phase
phase_thermalInstance , & !< instance of particular plasticity of each phase
2013-02-20 03:42:05 +05:30
crystallite_Noutput , & !< number of '(output)' items per crystallite setting
homogenization_typeInstance , & !< instance of particular type of each homogenization
microstructure_crystallite !< crystallite setting ID of each microstructure
2012-03-09 01:55:28 +05:30
2014-03-12 13:03:51 +05:30
integer ( pInt ) , dimension ( : , : , : ) , allocatable , public :: &
2014-05-08 20:25:19 +05:30
material_phase !< phase (index) of each grain,IP,element
2014-04-15 14:50:38 +05:30
#ifdef NEWSTATE
2014-05-08 20:25:19 +05:30
type ( tState ) , allocatable , dimension ( : ) , public :: &
2014-04-15 14:50:38 +05:30
plasticState , &
2014-06-23 00:28:29 +05:30
elasticState , &
damageState , &
thermalState
2014-04-15 14:50:38 +05:30
#endif
2014-05-08 20:25:19 +05:30
2013-12-16 16:26:56 +05:30
integer ( pInt ) , dimension ( : , : , : ) , allocatable , public , protected :: &
2013-02-20 03:42:05 +05:30
material_texture !< texture (index) of each grain,IP,element
2012-03-09 01:55:28 +05:30
2013-12-16 16:26:56 +05:30
real ( pReal ) , dimension ( : , : , : , : ) , allocatable , public , protected :: &
2013-02-20 03:42:05 +05:30
material_EulerAngles !< initial orientation of each grain,IP,element
2012-03-09 01:55:28 +05:30
2013-12-16 16:26:56 +05:30
logical , dimension ( : ) , allocatable , public , protected :: &
2012-03-09 01:55:28 +05:30
microstructure_active , &
2013-02-20 03:42:05 +05:30
microstructure_elemhomo , & !< flag to indicate homogeneous microstructure distribution over element's IPs
phase_localPlasticity !< flags phases with local constitutive law
2012-03-09 01:55:28 +05:30
2013-12-16 16:26:56 +05:30
character ( len = * ) , parameter , private :: &
2013-06-27 00:49:00 +05:30
MATERIAL_partMicrostructure = 'microstructure' , & !< keyword for microstructure part
MATERIAL_partTexture = 'texture' !< keyword for texture part
2012-03-09 01:55:28 +05:30
2013-12-16 16:26:56 +05:30
character ( len = 64 ) , dimension ( : ) , allocatable , private :: &
2013-02-20 03:42:05 +05:30
microstructure_name , & !< name of each microstructure
texture_name !< name of each texture
2012-03-09 01:55:28 +05:30
2013-12-16 16:26:56 +05:30
character ( len = 256 ) , dimension ( : ) , allocatable , private :: &
2013-02-20 03:42:05 +05:30
texture_ODFfile !< name of each ODF file
2012-03-09 01:55:28 +05:30
2013-12-16 16:26:56 +05:30
integer ( pInt ) , private :: &
2013-02-20 03:42:05 +05:30
material_Ntexture , & !< number of textures
microstructure_maxNconstituents , & !< max number of constituents in any phase
texture_maxNgauss , & !< max number of Gauss components in any texture
texture_maxNfiber !< max number of Fiber components in any texture
2012-03-09 01:55:28 +05:30
2013-12-16 16:26:56 +05:30
integer ( pInt ) , dimension ( : ) , allocatable , private :: &
2013-02-20 03:42:05 +05:30
microstructure_Nconstituents , & !< number of constituents in each microstructure
texture_symmetry , & !< number of symmetric orientations per texture
texture_Ngauss , & !< number of Gauss components per texture
texture_Nfiber !< number of Fiber components per texture
2012-03-09 01:55:28 +05:30
2013-12-16 16:26:56 +05:30
integer ( pInt ) , dimension ( : , : ) , allocatable , private :: &
2013-02-20 03:42:05 +05:30
microstructure_phase , & !< phase IDs of each microstructure
microstructure_texture !< texture IDs of each microstructure
2012-03-09 01:55:28 +05:30
2013-12-16 16:26:56 +05:30
real ( pReal ) , dimension ( : , : ) , allocatable , private :: &
2013-02-20 03:42:05 +05:30
microstructure_fraction !< vol fraction of each constituent in microstructure
2012-03-09 01:55:28 +05:30
2013-12-16 16:26:56 +05:30
real ( pReal ) , dimension ( : , : , : ) , allocatable , private :: &
2013-02-20 03:42:05 +05:30
material_volume , & !< volume of each grain,IP,element
texture_Gauss , & !< data of each Gauss component
2013-05-02 14:05:37 +05:30
texture_Fiber , & !< data of each Fiber component
2013-07-24 16:39:39 +05:30
texture_transformation !< transformation for each texture
2012-03-09 01:55:28 +05:30
2013-12-16 16:26:56 +05:30
logical , dimension ( : ) , allocatable , private :: &
2012-03-09 01:55:28 +05:30
homogenization_active
2014-06-23 00:28:29 +05:30
#if defined(HDF) || defined(NEWSTATE)
integer ( pInt ) , dimension ( : , : , : , : ) , allocatable , public , protected :: mappingConstitutive
integer ( pInt ) , dimension ( : , : , : ) , allocatable , public , protected :: mappingCrystallite
integer ( pInt ) , dimension ( : ) , allocatable :: ConstitutivePosition
integer ( pInt ) , dimension ( : ) , allocatable :: CrystallitePosition
#endif
2012-03-09 01:55:28 +05:30
2013-10-16 18:34:59 +05:30
public :: &
2013-11-27 13:35:23 +05:30
material_init , &
2013-12-16 16:26:56 +05:30
ELASTICITY_hooke_ID , &
2013-11-27 13:35:23 +05:30
PLASTICITY_none_ID , &
PLASTICITY_J2_ID , &
PLASTICITY_phenopowerlaw_ID , &
PLASTICITY_dislotwin_ID , &
PLASTICITY_titanmod_ID , &
PLASTICITY_nonlocal_ID , &
2014-06-23 00:28:29 +05:30
DAMAGE_undefined_ID , &
DAMAGE_none_ID , &
DAMAGE_local_ID , &
DAMAGE_gradient_ID , &
THERMAL_undefined_ID , &
THERMAL_none_ID , &
THERMAL_iso_ID , &
THERMAL_conduction_ID , &
THERMAL_adiabatic_ID , &
2014-03-14 04:50:50 +05:30
HOMOGENIZATION_none_ID , &
2013-11-27 13:35:23 +05:30
HOMOGENIZATION_isostrain_ID , &
2014-03-25 22:51:47 +05:30
#ifdef HDF
material_NconstituentsPhase , &
#endif
2013-11-27 13:35:23 +05:30
HOMOGENIZATION_RGC_ID
2013-10-16 18:34:59 +05:30
private :: &
material_parseHomogenization , &
material_parseMicrostructure , &
material_parseCrystallite , &
material_parsePhase , &
material_parseTexture , &
material_populateGrains
2012-03-09 01:55:28 +05:30
contains
2009-03-04 17:18:54 +05:30
2013-01-18 17:00:52 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief parses material configuration file
!> @details figures out if solverJobName.materialConfig is present, if not looks for
!> material.config
!--------------------------------------------------------------------------------------------------
2012-03-09 01:55:28 +05:30
subroutine material_init
2013-01-18 17:00:52 +05:30
use , intrinsic :: iso_fortran_env ! to get compiler_version and compiler_options (at least for gfortran 4.6 at the moment)
use IO , only : &
IO_error , &
IO_open_file , &
2013-02-25 22:04:59 +05:30
IO_open_jobFile_stat , &
IO_timeStamp
2013-01-18 17:00:52 +05:30
use debug , only : &
debug_level , &
debug_material , &
debug_levelBasic , &
debug_levelExtensive
2014-06-23 00:28:29 +05:30
use mesh , only : &
mesh_maxNips , &
mesh_NcpElems , &
mesh_element , &
FE_Nips , &
FE_geomtype
2012-03-09 01:55:28 +05:30
2009-03-04 17:18:54 +05:30
implicit none
2013-12-12 22:39:59 +05:30
integer ( pInt ) , parameter :: FILEUNIT = 200_pInt
2013-03-28 19:20:20 +05:30
integer ( pInt ) :: m , c , h , myDebug
2014-06-23 00:28:29 +05:30
integer ( pInt ) :: &
g , & !< grain number
i , & !< integration point number
e , & !< element number
phase
2012-07-05 15:24:50 +05:30
myDebug = debug_level ( debug_material )
2009-06-18 19:58:02 +05:30
2013-01-18 17:00:52 +05:30
write ( 6 , '(/,a)' ) ' <<<+- material init -+>>>'
write ( 6 , '(a)' ) ' $Id$'
2014-05-15 15:10:43 +05:30
write ( 6 , '(a15,a)' ) ' Current time: ' , IO_timeStamp ( )
2012-02-01 00:48:55 +05:30
#include "compilation_info.f90"
2013-01-18 17:00:52 +05:30
2013-12-12 22:39:59 +05:30
if ( . not . IO_open_jobFile_stat ( FILEUNIT , material_localFileExt ) ) & ! no local material configuration present...
call IO_open_file ( FILEUNIT , material_configFile ) ! ...open material.config file
call material_parseHomogenization ( FILEUNIT , material_partHomogenization )
2013-11-27 13:35:23 +05:30
if ( iand ( myDebug , debug_levelBasic ) / = 0_pInt ) write ( 6 , '(a)' ) ' Homogenization parsed'
2013-12-12 22:39:59 +05:30
call material_parseMicrostructure ( FILEUNIT , material_partMicrostructure )
2013-11-27 13:35:23 +05:30
if ( iand ( myDebug , debug_levelBasic ) / = 0_pInt ) write ( 6 , '(a)' ) ' Microstructure parsed'
2013-12-12 22:39:59 +05:30
call material_parseCrystallite ( FILEUNIT , material_partCrystallite )
2013-11-27 13:35:23 +05:30
if ( iand ( myDebug , debug_levelBasic ) / = 0_pInt ) write ( 6 , '(a)' ) ' Crystallite parsed'
2013-12-12 22:39:59 +05:30
call material_parseTexture ( FILEUNIT , material_partTexture )
2013-11-27 13:35:23 +05:30
if ( iand ( myDebug , debug_levelBasic ) / = 0_pInt ) write ( 6 , '(a)' ) ' Texture parsed'
2013-12-12 22:39:59 +05:30
call material_parsePhase ( FILEUNIT , material_partPhase )
2013-11-27 13:35:23 +05:30
if ( iand ( myDebug , debug_levelBasic ) / = 0_pInt ) write ( 6 , '(a)' ) ' Phase parsed'
2013-12-12 22:39:59 +05:30
close ( FILEUNIT )
2014-04-15 14:50:38 +05:30
#ifdef NEWSTATE
allocate ( plasticState ( material_Nphase ) )
2014-05-08 20:25:19 +05:30
allocate ( elasticState ( material_Nphase ) )
2014-06-23 00:28:29 +05:30
allocate ( damageState ( material_Nphase ) )
allocate ( thermalState ( material_Nphase ) )
#endif
#if defined(HDF) || defined(NEWSTATE)
allocate ( mappingConstitutive ( 2 , homogenization_maxNgrains , mesh_maxNips , mesh_NcpElems ) , source = 0_pInt )
allocate ( mappingCrystallite ( 2 , homogenization_maxNgrains , mesh_NcpElems ) , source = 0_pInt )
allocate ( ConstitutivePosition ( material_Nphase ) , source = 0_pInt )
allocate ( CrystallitePosition ( material_Nphase ) , source = 0_pInt )
ElemLoop : do e = 1_pInt , mesh_NcpElems ! loop over elements
IPloop : do i = 1_pInt , FE_Nips ( FE_geomtype ( mesh_element ( 2 , e ) ) ) ! loop over IPs
GrainLoop : do g = 1_pInt , homogenization_Ngrains ( mesh_element ( 3 , e ) ) ! loop over grains
phase = material_phase ( g , i , e )
ConstitutivePosition ( phase ) = ConstitutivePosition ( phase ) + 1_pInt ! not distinguishing between instances of same phase
mappingConstitutive ( 1 : 2 , g , i , e ) = [ ConstitutivePosition ( phase ) , phase ]
enddo GrainLoop
enddo IPloop
enddo ElemLoop
2014-04-15 14:50:38 +05:30
#endif
2014-05-08 20:25:19 +05:30
2013-03-28 19:20:20 +05:30
do m = 1_pInt , material_Nmicrostructure
2013-11-27 13:35:23 +05:30
if ( microstructure_crystallite ( m ) < 1_pInt . or . &
microstructure_crystallite ( m ) > material_Ncrystallite ) &
call IO_error ( 150_pInt , m , ext_msg = 'crystallite' )
if ( minval ( microstructure_phase ( 1 : microstructure_Nconstituents ( m ) , m ) ) < 1_pInt . or . &
maxval ( microstructure_phase ( 1 : microstructure_Nconstituents ( m ) , m ) ) > material_Nphase ) &
call IO_error ( 150_pInt , m , ext_msg = 'phase' )
if ( minval ( microstructure_texture ( 1 : microstructure_Nconstituents ( m ) , m ) ) < 1_pInt . or . &
maxval ( microstructure_texture ( 1 : microstructure_Nconstituents ( m ) , m ) ) > material_Ntexture ) &
call IO_error ( 150_pInt , m , ext_msg = 'texture' )
if ( microstructure_Nconstituents ( m ) < 1_pInt ) &
call IO_error ( 151_pInt , m )
2009-03-04 17:18:54 +05:30
enddo
2013-03-28 19:20:20 +05:30
debugOut : if ( iand ( myDebug , debug_levelExtensive ) / = 0_pInt ) then
2013-01-18 17:00:52 +05:30
write ( 6 , '(/,a,/)' ) ' MATERIAL configuration'
write ( 6 , '(a32,1x,a16,1x,a6)' ) 'homogenization ' , 'type ' , 'grains'
2013-03-28 19:20:20 +05:30
do h = 1_pInt , material_Nhomogenization
2013-05-29 22:53:49 +05:30
write ( 6 , '(1x,a32,1x,a16,1x,i6)' ) homogenization_name ( h ) , homogenization_type ( h ) , homogenization_Ngrains ( h )
2013-01-18 17:00:52 +05:30
enddo
2013-05-29 22:53:49 +05:30
write ( 6 , '(/,a14,18x,1x,a11,1x,a12,1x,a13)' ) 'microstructure' , 'crystallite' , 'constituents' , 'homogeneous'
2013-03-28 19:20:20 +05:30
do m = 1_pInt , material_Nmicrostructure
2013-05-29 22:53:49 +05:30
write ( 6 , '(1x,a32,1x,i11,1x,i12,1x,l13)' ) microstructure_name ( m ) , &
2013-03-28 19:20:20 +05:30
microstructure_crystallite ( m ) , &
microstructure_Nconstituents ( m ) , &
microstructure_elemhomo ( m )
if ( microstructure_Nconstituents ( m ) > 0_pInt ) then
do c = 1_pInt , microstructure_Nconstituents ( m )
write ( 6 , '(a1,1x,a32,1x,a32,1x,f7.4)' ) '>' , phase_name ( microstructure_phase ( c , m ) ) , &
texture_name ( microstructure_texture ( c , m ) ) , &
microstructure_fraction ( c , m )
2013-01-18 17:00:52 +05:30
enddo
write ( 6 , * )
endif
enddo
2013-03-28 19:20:20 +05:30
endif debugOut
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
2012-03-09 01:55:28 +05:30
call material_populateGrains
2009-03-04 17:18:54 +05:30
2012-03-09 01:55:28 +05:30
end subroutine material_init
2009-03-04 17:18:54 +05:30
2013-01-18 17:00:52 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief parses the homogenization part in the material configuration file
!--------------------------------------------------------------------------------------------------
2013-12-12 22:39:59 +05:30
subroutine material_parseHomogenization ( fileUnit , myPart )
2013-03-28 19:20:20 +05:30
use IO , only : &
2013-06-27 00:49:00 +05:30
IO_read , &
2013-03-28 19:20:20 +05:30
IO_globalTagInPart , &
IO_countSections , &
IO_error , &
IO_countTagInPart , &
IO_lc , &
IO_getTag , &
IO_isBlank , &
IO_stringValue , &
IO_intValue , &
2013-12-12 22:39:59 +05:30
IO_stringPos , &
IO_EOF
2013-01-18 17:00:52 +05:30
use mesh , only : &
mesh_element
2012-03-09 01:55:28 +05:30
2009-03-04 17:18:54 +05:30
implicit none
2009-03-20 20:04:24 +05:30
character ( len = * ) , intent ( in ) :: myPart
2013-12-12 22:39:59 +05:30
integer ( pInt ) , intent ( in ) :: fileUnit
2012-03-09 01:55:28 +05:30
2013-12-12 22:39:59 +05:30
integer ( pInt ) , parameter :: MAXNCHUNKS = 2_pInt
2012-03-09 01:55:28 +05:30
2013-12-12 22:39:59 +05:30
integer ( pInt ) , dimension ( 1 + 2 * MAXNCHUNKS ) :: positions
2013-11-27 13:35:23 +05:30
integer ( pInt ) :: Nsections , section , s
2013-12-12 22:39:59 +05:30
character ( len = 65536 ) :: &
tag , line
2013-06-27 00:49:00 +05:30
logical :: echo
2013-12-16 16:26:56 +05:30
2013-12-12 22:39:59 +05:30
echo = IO_globalTagInPart ( fileUnit , myPart , '/echo/' )
Nsections = IO_countSections ( fileUnit , myPart )
2009-03-04 17:18:54 +05:30
material_Nhomogenization = Nsections
2012-02-13 23:11:27 +05:30
if ( Nsections < 1_pInt ) call IO_error ( 160_pInt , ext_msg = myPart )
2009-05-07 21:57:36 +05:30
2013-12-16 16:26:56 +05:30
allocate ( homogenization_name ( Nsections ) ) ; homogenization_name = ''
2014-03-14 04:50:50 +05:30
allocate ( homogenization_type ( Nsections ) , source = HOMOGENIZATION_undefined_ID )
allocate ( homogenization_typeInstance ( Nsections ) , source = 0_pInt )
allocate ( homogenization_Ngrains ( Nsections ) , source = 0_pInt )
allocate ( homogenization_Noutput ( Nsections ) , source = 0_pInt )
allocate ( homogenization_active ( Nsections ) , source = . false . )
2009-07-22 21:37:19 +05:30
2013-01-18 17:00:52 +05:30
forall ( s = 1_pInt : Nsections ) homogenization_active ( s ) = any ( mesh_element ( 3 , : ) == s ) ! current homogenization used in model? Homogenization view, maximum operations depend on maximum number of homog schemes
2013-12-12 22:39:59 +05:30
homogenization_Noutput = IO_countTagInPart ( fileUnit , myPart , '(output)' , Nsections )
2009-03-04 17:18:54 +05:30
2013-12-12 22:39:59 +05:30
rewind ( fileUnit )
line = '' ! to have it initialized
section = 0_pInt ! - " -
2013-12-19 14:19:47 +05:30
do while ( trim ( line ) / = IO_EOF . and . IO_lc ( IO_getTag ( line , '<' , '>' ) ) / = myPart ) ! wind forward to <homogenization>
2013-12-12 22:39:59 +05:30
line = IO_read ( fileUnit )
2009-03-04 17:18:54 +05:30
enddo
2013-07-08 21:18:13 +05:30
if ( echo ) write ( 6 , '(/,1x,a)' ) trim ( line ) ! echo part header
2009-03-04 17:18:54 +05:30
2013-12-12 22:39:59 +05:30
do while ( trim ( line ) / = IO_EOF ) ! read through sections of material part
line = IO_read ( fileUnit )
2013-01-18 17:00:52 +05:30
if ( IO_isBlank ( line ) ) cycle ! skip empty lines
2013-12-12 22:39:59 +05:30
if ( IO_getTag ( line , '<' , '>' ) / = '' ) then ! stop at next part
line = IO_read ( fileUnit , . true . ) ! reset IO_read
exit
endif
2013-07-08 21:18:13 +05:30
if ( echo ) write ( 6 , '(2x,a)' ) trim ( line ) ! echo back read lines
2013-01-18 17:00:52 +05:30
if ( IO_getTag ( line , '[' , ']' ) / = '' ) then ! next section
2012-02-16 00:28:38 +05:30
section = section + 1_pInt
2013-12-16 16:26:56 +05:30
homogenization_name ( section ) = IO_getTag ( line , '[' , ']' )
2009-03-04 17:18:54 +05:30
endif
2012-02-10 17:26:05 +05:30
if ( section > 0_pInt ) then
2013-12-12 22:39:59 +05:30
positions = IO_stringPos ( line , MAXNCHUNKS )
2013-01-18 17:00:52 +05:30
tag = IO_lc ( IO_stringValue ( line , positions , 1_pInt ) ) ! extract key
2009-03-04 17:18:54 +05:30
select case ( tag )
case ( 'type' )
2013-11-27 13:35:23 +05:30
select case ( IO_lc ( IO_stringValue ( line , positions , 2_pInt ) ) )
2014-03-14 04:50:50 +05:30
case ( HOMOGENIZATION_NONE_label )
homogenization_type ( section ) = HOMOGENIZATION_NONE_ID
homogenization_Ngrains ( section ) = 1_pInt
2013-11-27 13:35:23 +05:30
case ( HOMOGENIZATION_ISOSTRAIN_label )
homogenization_type ( section ) = HOMOGENIZATION_ISOSTRAIN_ID
case ( HOMOGENIZATION_RGC_label )
homogenization_type ( section ) = HOMOGENIZATION_RGC_ID
case default
call IO_error ( 500_pInt , ext_msg = trim ( IO_stringValue ( line , positions , 2_pInt ) ) )
end select
homogenization_typeInstance ( section ) = &
count ( homogenization_type == homogenization_type ( section ) ) ! count instances
2013-12-12 05:12:33 +05:30
case ( 'nconstituents' , 'ngrains' )
2012-02-10 17:26:05 +05:30
homogenization_Ngrains ( section ) = IO_intValue ( line , positions , 2_pInt )
2009-03-04 17:18:54 +05:30
end select
endif
enddo
2013-06-27 00:49:00 +05:30
homogenization_maxNgrains = maxval ( homogenization_Ngrains , homogenization_active )
2009-03-04 17:18:54 +05:30
2013-03-28 19:20:20 +05:30
end subroutine material_parseHomogenization
2009-03-04 17:18:54 +05:30
2013-01-18 17:00:52 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief parses the microstructure part in the material configuration file
!--------------------------------------------------------------------------------------------------
2013-12-12 22:39:59 +05:30
subroutine material_parseMicrostructure ( fileUnit , myPart )
2009-03-04 17:18:54 +05:30
use IO
2013-01-18 17:00:52 +05:30
use mesh , only : &
mesh_element , &
mesh_NcpElems
2012-03-09 01:55:28 +05:30
2009-03-04 17:18:54 +05:30
implicit none
2009-03-20 20:04:24 +05:30
character ( len = * ) , intent ( in ) :: myPart
2013-12-12 22:39:59 +05:30
integer ( pInt ) , intent ( in ) :: fileUnit
2012-03-09 01:55:28 +05:30
2013-12-12 22:39:59 +05:30
integer ( pInt ) , parameter :: MAXNCHUNKS = 7_pInt
2012-03-09 01:55:28 +05:30
2013-12-12 22:39:59 +05:30
integer ( pInt ) , dimension ( 1_pInt + 2_pInt * MAXNCHUNKS ) :: positions
2012-03-09 01:55:28 +05:30
integer ( pInt ) :: Nsections , section , constituent , e , i
2013-12-12 22:39:59 +05:30
character ( len = 65536 ) :: &
tag , line
2013-06-27 00:49:00 +05:30
logical :: echo
2012-06-26 15:54:54 +05:30
2013-12-12 22:39:59 +05:30
echo = IO_globalTagInPart ( fileUnit , myPart , '/echo/' )
2011-09-13 21:24:06 +05:30
2013-12-12 22:39:59 +05:30
Nsections = IO_countSections ( fileUnit , myPart )
2009-03-04 17:18:54 +05:30
material_Nmicrostructure = Nsections
2012-02-13 23:11:27 +05:30
if ( Nsections < 1_pInt ) call IO_error ( 160_pInt , ext_msg = myPart )
2010-02-25 23:09:11 +05:30
2013-12-16 16:26:56 +05:30
allocate ( microstructure_name ( Nsections ) ) ; microstructure_name = ''
2013-12-13 04:33:37 +05:30
allocate ( microstructure_crystallite ( Nsections ) , source = 0_pInt )
allocate ( microstructure_Nconstituents ( Nsections ) , source = 0_pInt )
allocate ( microstructure_active ( Nsections ) , source = . false . )
allocate ( microstructure_elemhomo ( Nsections ) , source = . false . )
2009-07-22 21:37:19 +05:30
2014-05-15 15:10:43 +05:30
if ( any ( mesh_element ( 4 , 1 : mesh_NcpElems ) > Nsections ) ) &
call IO_error ( 155_pInt , ext_msg = 'Microstructure in geometry > Sections in material.config' )
2013-01-18 17:00:52 +05:30
forall ( e = 1_pInt : mesh_NcpElems ) microstructure_active ( mesh_element ( 4 , e ) ) = . true . ! current microstructure used in model? Elementwise view, maximum N operations for N elements
2011-10-18 14:48:05 +05:30
2013-12-12 22:39:59 +05:30
microstructure_Nconstituents = IO_countTagInPart ( fileUnit , myPart , '(constituent)' , Nsections )
2009-03-04 17:18:54 +05:30
microstructure_maxNconstituents = maxval ( microstructure_Nconstituents )
2013-12-12 22:39:59 +05:30
microstructure_elemhomo = IO_spotTagInPart ( fileUnit , myPart , '/elementhomogeneous/' , Nsections )
2009-11-24 20:30:25 +05:30
2013-12-13 04:33:37 +05:30
allocate ( microstructure_phase ( microstructure_maxNconstituents , Nsections ) , source = 0_pInt )
allocate ( microstructure_texture ( microstructure_maxNconstituents , Nsections ) , source = 0_pInt )
allocate ( microstructure_fraction ( microstructure_maxNconstituents , Nsections ) , source = 0.0_pReal )
2013-12-16 16:26:56 +05:30
2013-12-12 22:39:59 +05:30
rewind ( fileUnit )
2013-03-28 19:20:20 +05:30
line = '' ! to have it initialized
section = 0_pInt ! - " -
constituent = 0_pInt ! - " -
2013-12-16 16:26:56 +05:30
2013-12-19 14:19:47 +05:30
do while ( trim ( line ) / = IO_EOF . and . IO_lc ( IO_getTag ( line , '<' , '>' ) ) / = myPart ) ! wind forward to <microstructure>
2013-12-12 22:39:59 +05:30
line = IO_read ( fileUnit )
2009-03-04 17:18:54 +05:30
enddo
2013-07-08 21:18:13 +05:30
if ( echo ) write ( 6 , '(/,1x,a)' ) trim ( line ) ! echo part header
2009-03-04 17:18:54 +05:30
2013-12-12 22:39:59 +05:30
do while ( trim ( line ) / = IO_EOF ) ! read through sections of material part
line = IO_read ( fileUnit )
2013-01-18 17:00:52 +05:30
if ( IO_isBlank ( line ) ) cycle ! skip empty lines
2013-12-12 22:39:59 +05:30
if ( IO_getTag ( line , '<' , '>' ) / = '' ) then ! stop at next part
line = IO_read ( fileUnit , . true . ) ! reset IO_read
exit
endif
2013-07-08 21:18:13 +05:30
if ( echo ) write ( 6 , '(2x,a)' ) trim ( line ) ! echo back read lines
2013-01-18 17:00:52 +05:30
if ( IO_getTag ( line , '[' , ']' ) / = '' ) then ! next section
2012-02-16 00:28:38 +05:30
section = section + 1_pInt
constituent = 0_pInt
2009-03-04 17:18:54 +05:30
microstructure_name ( section ) = IO_getTag ( line , '[' , ']' )
endif
2012-02-16 00:28:38 +05:30
if ( section > 0_pInt ) then
2013-12-12 22:39:59 +05:30
positions = IO_stringPos ( line , MAXNCHUNKS )
2013-01-18 17:00:52 +05:30
tag = IO_lc ( IO_stringValue ( line , positions , 1_pInt ) ) ! extract key
2009-03-04 17:18:54 +05:30
select case ( tag )
2010-02-25 23:09:11 +05:30
case ( 'crystallite' )
2012-02-10 17:26:05 +05:30
microstructure_crystallite ( section ) = IO_intValue ( line , positions , 2_pInt )
2009-03-04 17:18:54 +05:30
case ( '(constituent)' )
2012-02-16 00:28:38 +05:30
constituent = constituent + 1_pInt
do i = 2_pInt , 6_pInt , 2_pInt
2009-03-04 17:18:54 +05:30
tag = IO_lc ( IO_stringValue ( line , positions , i ) )
select case ( tag )
case ( 'phase' )
2012-02-10 17:26:05 +05:30
microstructure_phase ( constituent , section ) = IO_intValue ( line , positions , i + 1_pInt )
2009-03-04 17:18:54 +05:30
case ( 'texture' )
2012-02-10 17:26:05 +05:30
microstructure_texture ( constituent , section ) = IO_intValue ( line , positions , i + 1_pInt )
2009-03-04 17:18:54 +05:30
case ( 'fraction' )
2012-02-10 17:26:05 +05:30
microstructure_fraction ( constituent , section ) = IO_floatValue ( line , positions , i + 1_pInt )
2009-03-04 17:18:54 +05:30
end select
enddo
end select
endif
enddo
2013-06-27 00:49:00 +05:30
end subroutine material_parseMicrostructure
2009-03-04 17:18:54 +05:30
2013-01-18 17:00:52 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief parses the crystallite part in the material configuration file
!--------------------------------------------------------------------------------------------------
2013-12-12 22:39:59 +05:30
subroutine material_parseCrystallite ( fileUnit , myPart )
2013-01-18 17:00:52 +05:30
use IO , only : &
2013-06-27 00:49:00 +05:30
IO_read , &
2013-01-18 17:00:52 +05:30
IO_countSections , &
IO_error , &
IO_countTagInPart , &
IO_globalTagInPart , &
IO_getTag , &
IO_lc , &
2013-12-12 22:39:59 +05:30
IO_isBlank , &
IO_EOF
2010-02-25 23:09:11 +05:30
2012-03-09 01:55:28 +05:30
implicit none
2010-02-25 23:09:11 +05:30
character ( len = * ) , intent ( in ) :: myPart
2013-12-12 22:39:59 +05:30
integer ( pInt ) , intent ( in ) :: fileUnit
2012-03-09 01:55:28 +05:30
2013-06-27 00:49:00 +05:30
integer ( pInt ) :: Nsections , &
section
character ( len = 65536 ) :: line
logical :: echo
2012-06-26 15:54:54 +05:30
2013-12-12 22:39:59 +05:30
echo = IO_globalTagInPart ( fileUnit , myPart , '/echo/' )
2010-02-25 23:09:11 +05:30
2013-12-12 22:39:59 +05:30
Nsections = IO_countSections ( fileUnit , myPart )
2010-02-25 23:09:11 +05:30
material_Ncrystallite = Nsections
2012-02-13 23:11:27 +05:30
if ( Nsections < 1_pInt ) call IO_error ( 160_pInt , ext_msg = myPart )
2010-02-25 23:09:11 +05:30
2013-12-16 16:26:56 +05:30
allocate ( crystallite_name ( Nsections ) ) ; crystallite_name = ''
2013-12-13 04:33:37 +05:30
allocate ( crystallite_Noutput ( Nsections ) , source = 0_pInt )
2010-02-25 23:09:11 +05:30
2013-12-12 22:39:59 +05:30
crystallite_Noutput = IO_countTagInPart ( fileUnit , myPart , '(output)' , Nsections )
2010-02-25 23:09:11 +05:30
2013-12-12 22:39:59 +05:30
rewind ( fileUnit )
line = '' ! to have it initialized
section = 0_pInt ! - " -
2013-12-19 14:19:47 +05:30
do while ( trim ( line ) / = IO_EOF . and . IO_lc ( IO_getTag ( line , '<' , '>' ) ) / = myPart ) ! wind forward to <Crystallite>
2013-12-12 22:39:59 +05:30
line = IO_read ( fileUnit )
2010-02-25 23:09:11 +05:30
enddo
2013-07-08 21:18:13 +05:30
if ( echo ) write ( 6 , '(/,1x,a)' ) trim ( line ) ! echo part header
2010-02-25 23:09:11 +05:30
2013-12-12 22:39:59 +05:30
do while ( trim ( line ) / = IO_EOF ) ! read through sections of material part
line = IO_read ( fileUnit )
2013-01-18 17:00:52 +05:30
if ( IO_isBlank ( line ) ) cycle ! skip empty lines
2013-12-12 22:39:59 +05:30
if ( IO_getTag ( line , '<' , '>' ) / = '' ) then ! stop at next part
line = IO_read ( fileUnit , . true . ) ! reset IO_read
exit
endif
2013-07-08 21:18:13 +05:30
if ( echo ) write ( 6 , '(2x,a)' ) trim ( line ) ! echo back read lines
2013-01-18 17:00:52 +05:30
if ( IO_getTag ( line , '[' , ']' ) / = '' ) then ! next section
2012-02-16 00:28:38 +05:30
section = section + 1_pInt
2010-02-25 23:09:11 +05:30
crystallite_name ( section ) = IO_getTag ( line , '[' , ']' )
endif
enddo
2013-06-27 00:49:00 +05:30
end subroutine material_parseCrystallite
2010-02-25 23:09:11 +05:30
2013-01-18 17:00:52 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief parses the phase part in the material configuration file
!--------------------------------------------------------------------------------------------------
2013-12-12 22:39:59 +05:30
subroutine material_parsePhase ( fileUnit , myPart )
2013-03-28 19:20:20 +05:30
use IO , only : &
2013-06-27 00:49:00 +05:30
IO_read , &
2013-03-28 19:20:20 +05:30
IO_globalTagInPart , &
IO_countSections , &
IO_error , &
IO_countTagInPart , &
IO_getTag , &
IO_spotTagInPart , &
IO_lc , &
IO_isBlank , &
IO_stringValue , &
2013-12-12 22:39:59 +05:30
IO_stringPos , &
IO_EOF
2013-03-28 19:20:20 +05:30
2009-03-04 17:18:54 +05:30
implicit none
2009-03-20 20:04:24 +05:30
character ( len = * ) , intent ( in ) :: myPart
2013-12-12 22:39:59 +05:30
integer ( pInt ) , intent ( in ) :: fileUnit
2012-03-09 01:55:28 +05:30
2013-12-12 22:39:59 +05:30
integer ( pInt ) , parameter :: MAXNCHUNKS = 2_pInt
2012-03-09 01:55:28 +05:30
2013-12-12 22:39:59 +05:30
integer ( pInt ) , dimension ( 1 + 2 * MAXNCHUNKS ) :: positions
integer ( pInt ) :: Nsections , section
character ( len = 65536 ) :: &
tag , line
2013-06-27 00:49:00 +05:30
logical :: echo
2012-06-26 15:54:54 +05:30
2013-12-12 22:39:59 +05:30
echo = IO_globalTagInPart ( fileUnit , myPart , '/echo/' )
2009-03-04 17:18:54 +05:30
2013-12-12 22:39:59 +05:30
Nsections = IO_countSections ( fileUnit , myPart )
2009-03-04 17:18:54 +05:30
material_Nphase = Nsections
2012-02-13 23:11:27 +05:30
if ( Nsections < 1_pInt ) call IO_error ( 160_pInt , ext_msg = myPart )
2010-02-25 23:09:11 +05:30
2013-03-28 19:20:20 +05:30
allocate ( phase_name ( Nsections ) ) ; phase_name = ''
2013-12-13 04:33:37 +05:30
allocate ( phase_elasticity ( Nsections ) , source = ELASTICITY_undefined_ID )
allocate ( phase_elasticityInstance ( Nsections ) , source = 0_pInt )
allocate ( phase_plasticity ( Nsections ) , source = PLASTICITY_undefined_ID )
allocate ( phase_plasticityInstance ( Nsections ) , source = 0_pInt )
2014-06-23 00:28:29 +05:30
allocate ( phase_damage ( Nsections ) , source = DAMAGE_undefined_ID )
allocate ( phase_damageInstance ( Nsections ) , source = 0_pInt )
allocate ( phase_thermal ( Nsections ) , source = THERMAL_undefined_ID )
allocate ( phase_thermalInstance ( Nsections ) , source = 0_pInt )
2013-12-13 04:33:37 +05:30
allocate ( phase_Noutput ( Nsections ) , source = 0_pInt )
allocate ( phase_localPlasticity ( Nsections ) , source = . false . )
2009-03-04 17:18:54 +05:30
2013-12-12 22:39:59 +05:30
phase_Noutput = IO_countTagInPart ( fileUnit , myPart , '(output)' , Nsections )
phase_localPlasticity = . not . IO_spotTagInPart ( fileUnit , myPart , '/nonlocal/' , Nsections )
2009-03-04 17:18:54 +05:30
2013-12-12 22:39:59 +05:30
rewind ( fileUnit )
line = '' ! to have it initialized
section = 0_pInt ! - " -
2013-12-19 14:19:47 +05:30
do while ( trim ( line ) / = IO_EOF . and . IO_lc ( IO_getTag ( line , '<' , '>' ) ) / = myPart ) ! wind forward to <Phase>
2013-12-12 22:39:59 +05:30
line = IO_read ( fileUnit )
2009-03-04 17:18:54 +05:30
enddo
2013-07-08 21:18:13 +05:30
if ( echo ) write ( 6 , '(/,1x,a)' ) trim ( line ) ! echo part header
2009-03-04 17:18:54 +05:30
2013-12-12 22:39:59 +05:30
do while ( trim ( line ) / = IO_EOF ) ! read through sections of material part
line = IO_read ( fileUnit )
2013-01-18 17:00:52 +05:30
if ( IO_isBlank ( line ) ) cycle ! skip empty lines
2013-12-12 22:39:59 +05:30
if ( IO_getTag ( line , '<' , '>' ) / = '' ) then ! stop at next part
line = IO_read ( fileUnit , . true . ) ! reset IO_read
exit
endif
2013-07-08 21:18:13 +05:30
if ( echo ) write ( 6 , '(2x,a)' ) trim ( line ) ! echo back read lines
2013-01-18 17:00:52 +05:30
if ( IO_getTag ( line , '[' , ']' ) / = '' ) then ! next section
2012-02-16 00:28:38 +05:30
section = section + 1_pInt
2009-03-04 17:18:54 +05:30
phase_name ( section ) = IO_getTag ( line , '[' , ']' )
endif
2012-02-10 17:26:05 +05:30
if ( section > 0_pInt ) then
2013-12-12 22:39:59 +05:30
positions = IO_stringPos ( line , MAXNCHUNKS )
2013-01-18 17:00:52 +05:30
tag = IO_lc ( IO_stringValue ( line , positions , 1_pInt ) ) ! extract key
2009-03-04 17:18:54 +05:30
select case ( tag )
2012-03-14 21:46:11 +05:30
case ( 'elasticity' )
2013-11-27 13:35:23 +05:30
select case ( IO_lc ( IO_stringValue ( line , positions , 2_pInt ) ) )
case ( ELASTICITY_HOOKE_label )
phase_elasticity ( section ) = ELASTICITY_HOOKE_ID
case default
call IO_error ( 200_pInt , ext_msg = trim ( IO_stringValue ( line , positions , 2_pInt ) ) )
end select
phase_elasticityInstance ( section ) = count ( phase_elasticity == phase_elasticity ( section ) ) ! count instances
2012-03-12 19:39:37 +05:30
case ( 'plasticity' )
2013-11-27 13:35:23 +05:30
select case ( IO_lc ( IO_stringValue ( line , positions , 2_pInt ) ) )
case ( PLASTICITY_NONE_label )
phase_plasticity ( section ) = PLASTICITY_NONE_ID
case ( PLASTICITY_J2_label )
phase_plasticity ( section ) = PLASTICITY_J2_ID
case ( PLASTICITY_PHENOPOWERLAW_label )
phase_plasticity ( section ) = PLASTICITY_PHENOPOWERLAW_ID
case ( PLASTICITY_DISLOTWIN_label )
phase_plasticity ( section ) = PLASTICITY_DISLOTWIN_ID
case ( PLASTICITY_TITANMOD_label )
phase_plasticity ( section ) = PLASTICITY_TITANMOD_ID
case ( PLASTICITY_NONLOCAL_label )
phase_plasticity ( section ) = PLASTICITY_NONLOCAL_ID
case default
call IO_error ( 201_pInt , ext_msg = trim ( IO_stringValue ( line , positions , 2_pInt ) ) )
end select
phase_plasticityInstance ( section ) = count ( phase_plasticity == phase_plasticity ( section ) ) ! count instances
2014-06-23 00:28:29 +05:30
case ( 'damage' )
select case ( IO_lc ( IO_stringValue ( line , positions , 2_pInt ) ) )
case ( DAMAGE_NONE_label )
phase_damage ( section ) = DAMAGE_none_ID
case ( DAMAGE_LOCAL_label )
phase_damage ( section ) = DAMAGE_local_ID
case ( DAMAGE_GRADIENT_label )
phase_damage ( section ) = DAMAGE_gradient_ID
case default
call IO_error ( 200_pInt , ext_msg = trim ( IO_stringValue ( line , positions , 2_pInt ) ) )
end select
phase_damageInstance ( section ) = count ( phase_damage == phase_damage ( section ) ) ! count instances
case ( 'thermal' )
select case ( IO_lc ( IO_stringValue ( line , positions , 2_pInt ) ) )
case ( THERMAL_NONE_label )
phase_thermal ( section ) = THERMAL_none_ID
case ( THERMAL_ISO_label )
phase_thermal ( section ) = THERMAL_iso_ID
case ( THERMAL_CONDUCTION_label )
phase_thermal ( section ) = THERMAL_conduction_ID
case ( THERMAL_ADIABATIC_label )
phase_thermal ( section ) = THERMAL_adiabatic_ID
case default
call IO_error ( 200_pInt , ext_msg = trim ( IO_stringValue ( line , positions , 2_pInt ) ) )
end select
phase_thermalInstance ( section ) = count ( phase_thermal == phase_thermal ( section ) ) ! count instances
2009-03-04 17:18:54 +05:30
end select
endif
enddo
2013-06-27 00:49:00 +05:30
end subroutine material_parsePhase
2009-03-04 17:18:54 +05:30
2013-01-18 17:00:52 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief parses the texture part in the material configuration file
!--------------------------------------------------------------------------------------------------
2013-12-12 22:39:59 +05:30
subroutine material_parseTexture ( fileUnit , myPart )
2013-03-28 19:20:20 +05:30
use IO , only : &
2013-06-27 00:49:00 +05:30
IO_read , &
2013-03-28 19:20:20 +05:30
IO_globalTagInPart , &
IO_countSections , &
IO_error , &
IO_countTagInPart , &
IO_getTag , &
IO_spotTagInPart , &
IO_lc , &
IO_isBlank , &
IO_floatValue , &
IO_stringValue , &
2013-12-12 22:39:59 +05:30
IO_stringPos , &
IO_EOF
2013-01-18 17:00:52 +05:30
use math , only : &
inRad , &
2013-05-02 14:05:37 +05:30
math_sampleRandomOri , &
math_I3 , &
math_inv33
2012-03-09 01:55:28 +05:30
2009-03-04 17:18:54 +05:30
implicit none
2009-03-20 20:04:24 +05:30
character ( len = * ) , intent ( in ) :: myPart
2013-12-12 22:39:59 +05:30
integer ( pInt ) , intent ( in ) :: fileUnit
2012-03-09 01:55:28 +05:30
2013-12-12 22:39:59 +05:30
integer ( pInt ) , parameter :: MAXNCHUNKS = 13_pInt
2012-03-09 01:55:28 +05:30
2013-12-12 22:39:59 +05:30
integer ( pInt ) , dimension ( 1 + 2 * MAXNCHUNKS ) :: positions
2013-03-28 19:20:20 +05:30
integer ( pInt ) :: Nsections , section , gauss , fiber , j
2013-06-27 00:49:00 +05:30
character ( len = 65536 ) :: tag
character ( len = 65536 ) :: line
logical :: echo
2012-06-26 15:54:54 +05:30
2013-12-12 22:39:59 +05:30
echo = IO_globalTagInPart ( fileUnit , myPart , '/echo/' )
2009-03-04 17:18:54 +05:30
2013-12-12 22:39:59 +05:30
Nsections = IO_countSections ( fileUnit , myPart )
2009-03-04 17:18:54 +05:30
material_Ntexture = Nsections
2012-02-13 23:11:27 +05:30
if ( Nsections < 1_pInt ) call IO_error ( 160_pInt , ext_msg = myPart )
2010-02-25 23:09:11 +05:30
2013-12-13 04:33:37 +05:30
allocate ( texture_name ( Nsections ) ) ; texture_name = ''
allocate ( texture_ODFfile ( Nsections ) ) ; texture_ODFfile = ''
2013-12-16 16:26:56 +05:30
allocate ( texture_symmetry ( Nsections ) , source = 1_pInt )
2013-12-13 04:33:37 +05:30
allocate ( texture_Ngauss ( Nsections ) , source = 0_pInt )
allocate ( texture_Nfiber ( Nsections ) , source = 0_pInt )
2009-03-04 17:18:54 +05:30
2013-12-12 22:39:59 +05:30
texture_Ngauss = IO_countTagInPart ( fileUnit , myPart , '(gauss)' , Nsections ) + &
IO_countTagInPart ( fileUnit , myPart , '(random)' , Nsections )
texture_Nfiber = IO_countTagInPart ( fileUnit , myPart , '(fiber)' , Nsections )
2009-03-04 17:18:54 +05:30
texture_maxNgauss = maxval ( texture_Ngauss )
texture_maxNfiber = maxval ( texture_Nfiber )
2013-12-13 04:33:37 +05:30
allocate ( texture_Gauss ( 5 , texture_maxNgauss , Nsections ) , source = 0.0_pReal )
allocate ( texture_Fiber ( 6 , texture_maxNfiber , Nsections ) , source = 0.0_pReal )
allocate ( texture_transformation ( 3 , 3 , Nsections ) , source = 0.0_pReal )
texture_transformation = spread ( math_I3 , 3 , Nsections )
2013-12-16 16:26:56 +05:30
2013-12-12 22:39:59 +05:30
rewind ( fileUnit )
2013-03-28 19:20:20 +05:30
line = '' ! to have in initialized
section = 0_pInt ! - " -
gauss = 0_pInt ! - " -
fiber = 0_pInt ! - " -
2013-12-19 14:19:47 +05:30
do while ( trim ( line ) / = IO_EOF . and . IO_lc ( IO_getTag ( line , '<' , '>' ) ) / = myPart ) ! wind forward to <texture>
2013-12-12 22:39:59 +05:30
line = IO_read ( fileUnit )
2009-03-04 17:18:54 +05:30
enddo
2013-07-08 21:18:13 +05:30
if ( echo ) write ( 6 , '(/,1x,a)' ) trim ( line ) ! echo part header
2009-03-04 17:18:54 +05:30
2013-12-12 22:39:59 +05:30
do while ( trim ( line ) / = IO_EOF )
line = IO_read ( fileUnit )
2013-01-18 17:00:52 +05:30
if ( IO_isBlank ( line ) ) cycle ! skip empty lines
2013-12-12 22:39:59 +05:30
if ( IO_getTag ( line , '<' , '>' ) / = '' ) then ! stop at next part
line = IO_read ( fileUnit , . true . ) ! reset IO_read
exit
endif
2013-07-08 21:18:13 +05:30
if ( echo ) write ( 6 , '(2x,a)' ) trim ( line ) ! echo back read lines
2013-01-18 17:00:52 +05:30
if ( IO_getTag ( line , '[' , ']' ) / = '' ) then ! next section
2012-02-16 00:28:38 +05:30
section = section + 1_pInt
gauss = 0_pInt
fiber = 0_pInt
2009-03-04 17:18:54 +05:30
texture_name ( section ) = IO_getTag ( line , '[' , ']' )
endif
2012-02-16 00:28:38 +05:30
if ( section > 0_pInt ) then
2013-12-12 22:39:59 +05:30
positions = IO_stringPos ( line , MAXNCHUNKS )
2013-01-18 17:00:52 +05:30
tag = IO_lc ( IO_stringValue ( line , positions , 1_pInt ) ) ! extract key
2013-03-28 19:20:20 +05:30
textureType : select case ( tag )
2009-03-04 17:18:54 +05:30
2013-07-24 16:39:39 +05:30
case ( 'axes' , 'rotation' ) textureType
2013-05-29 22:53:49 +05:30
do j = 1_pInt , 3_pInt ! look for "x", "y", and "z" entries
2013-05-02 14:05:37 +05:30
tag = IO_lc ( IO_stringValue ( line , positions , j + 1_pInt ) )
select case ( tag )
case ( 'x' , '+x' )
2013-12-16 16:26:56 +05:30
texture_transformation ( j , 1 : 3 , section ) = [ 1.0_pReal , 0.0_pReal , 0.0_pReal ] ! original axis is now +x-axis
2013-05-02 14:05:37 +05:30
case ( '-x' )
2013-12-16 16:26:56 +05:30
texture_transformation ( j , 1 : 3 , section ) = [ - 1.0_pReal , 0.0_pReal , 0.0_pReal ] ! original axis is now -x-axis
2013-05-02 14:05:37 +05:30
case ( 'y' , '+y' )
2013-12-16 16:26:56 +05:30
texture_transformation ( j , 1 : 3 , section ) = [ 0.0_pReal , 1.0_pReal , 0.0_pReal ] ! original axis is now +y-axis
2013-05-02 14:05:37 +05:30
case ( '-y' )
2013-12-16 16:26:56 +05:30
texture_transformation ( j , 1 : 3 , section ) = [ 0.0_pReal , - 1.0_pReal , 0.0_pReal ] ! original axis is now -y-axis
2013-05-02 14:05:37 +05:30
case ( 'z' , '+z' )
2013-12-16 16:26:56 +05:30
texture_transformation ( j , 1 : 3 , section ) = [ 0.0_pReal , 0.0_pReal , 1.0_pReal ] ! original axis is now +z-axis
2013-05-02 14:05:37 +05:30
case ( '-z' )
2013-12-16 16:26:56 +05:30
texture_transformation ( j , 1 : 3 , section ) = [ 0.0_pReal , 0.0_pReal , - 1.0_pReal ] ! original axis is now -z-axis
2013-05-02 14:05:37 +05:30
case default
2013-05-29 22:53:49 +05:30
call IO_error ( 157_pInt , section )
2013-05-02 14:05:37 +05:30
end select
enddo
2013-03-28 19:20:20 +05:30
case ( 'hybridia' ) textureType
2012-02-10 17:26:05 +05:30
texture_ODFfile ( section ) = IO_stringValue ( line , positions , 2_pInt )
2009-03-04 17:18:54 +05:30
2013-03-28 19:20:20 +05:30
case ( 'symmetry' ) textureType
2012-02-10 17:26:05 +05:30
tag = IO_lc ( IO_stringValue ( line , positions , 2_pInt ) )
2009-03-04 17:18:54 +05:30
select case ( tag )
case ( 'orthotropic' )
2012-02-16 00:28:38 +05:30
texture_symmetry ( section ) = 4_pInt
2009-03-04 17:18:54 +05:30
case ( 'monoclinic' )
2012-02-16 00:28:38 +05:30
texture_symmetry ( section ) = 2_pInt
2009-03-04 17:18:54 +05:30
case default
2012-02-16 00:28:38 +05:30
texture_symmetry ( section ) = 1_pInt
2009-03-04 17:18:54 +05:30
end select
2013-03-28 19:20:20 +05:30
case ( '(random)' ) textureType
2012-02-16 00:28:38 +05:30
gauss = gauss + 1_pInt
2011-05-28 15:14:43 +05:30
texture_Gauss ( 1 : 3 , gauss , section ) = math_sampleRandomOri ( )
2013-03-28 19:20:20 +05:30
do j = 2_pInt , 4_pInt , 2_pInt
tag = IO_lc ( IO_stringValue ( line , positions , j ) )
2011-05-28 15:14:43 +05:30
select case ( tag )
case ( 'scatter' )
2013-03-28 19:20:20 +05:30
texture_Gauss ( 4 , gauss , section ) = IO_floatValue ( line , positions , j + 1_pInt ) * inRad
2011-05-28 15:14:43 +05:30
case ( 'fraction' )
2013-03-28 19:20:20 +05:30
texture_Gauss ( 5 , gauss , section ) = IO_floatValue ( line , positions , j + 1_pInt )
2011-05-28 15:14:43 +05:30
end select
enddo
2013-03-28 19:20:20 +05:30
case ( '(gauss)' ) textureType
2012-02-16 00:28:38 +05:30
gauss = gauss + 1_pInt
2013-03-28 19:20:20 +05:30
do j = 2_pInt , 10_pInt , 2_pInt
tag = IO_lc ( IO_stringValue ( line , positions , j ) )
2009-03-04 17:18:54 +05:30
select case ( tag )
case ( 'phi1' )
2013-03-28 19:20:20 +05:30
texture_Gauss ( 1 , gauss , section ) = IO_floatValue ( line , positions , j + 1_pInt ) * inRad
2009-03-04 17:18:54 +05:30
case ( 'phi' )
2013-03-28 19:20:20 +05:30
texture_Gauss ( 2 , gauss , section ) = IO_floatValue ( line , positions , j + 1_pInt ) * inRad
2009-03-04 17:18:54 +05:30
case ( 'phi2' )
2013-03-28 19:20:20 +05:30
texture_Gauss ( 3 , gauss , section ) = IO_floatValue ( line , positions , j + 1_pInt ) * inRad
2009-03-04 17:18:54 +05:30
case ( 'scatter' )
2013-03-28 19:20:20 +05:30
texture_Gauss ( 4 , gauss , section ) = IO_floatValue ( line , positions , j + 1_pInt ) * inRad
2009-03-04 17:18:54 +05:30
case ( 'fraction' )
2013-03-28 19:20:20 +05:30
texture_Gauss ( 5 , gauss , section ) = IO_floatValue ( line , positions , j + 1_pInt )
2009-03-04 17:18:54 +05:30
end select
enddo
2011-05-28 15:14:43 +05:30
2013-03-28 19:20:20 +05:30
case ( '(fiber)' ) textureType
2012-02-16 00:28:38 +05:30
fiber = fiber + 1_pInt
2013-03-28 19:20:20 +05:30
do j = 2_pInt , 12_pInt , 2_pInt
tag = IO_lc ( IO_stringValue ( line , positions , j ) )
2009-03-04 17:18:54 +05:30
select case ( tag )
2009-04-28 15:15:52 +05:30
case ( 'alpha1' )
2013-03-28 19:20:20 +05:30
texture_Fiber ( 1 , fiber , section ) = IO_floatValue ( line , positions , j + 1_pInt ) * inRad
2009-03-04 17:18:54 +05:30
case ( 'alpha2' )
2013-03-28 19:20:20 +05:30
texture_Fiber ( 2 , fiber , section ) = IO_floatValue ( line , positions , j + 1_pInt ) * inRad
2009-03-04 17:18:54 +05:30
case ( 'beta1' )
2013-03-28 19:20:20 +05:30
texture_Fiber ( 3 , fiber , section ) = IO_floatValue ( line , positions , j + 1_pInt ) * inRad
2009-03-04 17:18:54 +05:30
case ( 'beta2' )
2013-03-28 19:20:20 +05:30
texture_Fiber ( 4 , fiber , section ) = IO_floatValue ( line , positions , j + 1_pInt ) * inRad
2009-03-04 17:18:54 +05:30
case ( 'scatter' )
2013-03-28 19:20:20 +05:30
texture_Fiber ( 5 , fiber , section ) = IO_floatValue ( line , positions , j + 1_pInt ) * inRad
2009-03-04 17:18:54 +05:30
case ( 'fraction' )
2013-03-28 19:20:20 +05:30
texture_Fiber ( 6 , fiber , section ) = IO_floatValue ( line , positions , j + 1_pInt )
2009-03-04 17:18:54 +05:30
end select
enddo
2013-03-28 19:20:20 +05:30
end select textureType
2009-03-04 17:18:54 +05:30
endif
enddo
2013-06-27 00:49:00 +05:30
end subroutine material_parseTexture
2009-03-04 17:18:54 +05:30
2013-01-18 17:00:52 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief populates the grains
!> @details populates the grains by identifying active microstructure/homogenization pairs,
!! calculates the volume of the grains and deals with texture components and hybridIA
!--------------------------------------------------------------------------------------------------
2012-03-09 01:55:28 +05:30
subroutine material_populateGrains
2013-01-18 17:00:52 +05:30
use math , only : &
2013-05-29 22:53:49 +05:30
math_RtoEuler , &
math_EulerToR , &
math_mul33x33 , &
math_range , &
2013-01-18 17:00:52 +05:30
math_sampleRandomOri , &
math_sampleGaussOri , &
math_sampleFiberOri , &
math_symmetricEulers
use mesh , only : &
mesh_element , &
mesh_maxNips , &
mesh_NcpElems , &
mesh_ipVolume , &
FE_Nips , &
FE_geomtype
use IO , only : &
IO_error , &
IO_hybridIA
use FEsolving , only : &
FEsolving_execIP
use debug , only : &
debug_level , &
debug_material , &
debug_levelBasic
2012-03-09 01:55:28 +05:30
2009-03-04 17:18:54 +05:30
implicit none
integer ( pInt ) , dimension ( : , : ) , allocatable :: Ngrains
2013-05-29 22:53:49 +05:30
integer ( pInt ) , dimension ( microstructure_maxNconstituents ) :: &
NgrainsOfConstituent , &
currentGrainOfConstituent , &
randomOrder
real ( pReal ) , dimension ( microstructure_maxNconstituents ) :: &
rndArray
2010-02-18 21:24:10 +05:30
real ( pReal ) , dimension ( : ) , allocatable :: volumeOfGrain
2009-03-04 17:18:54 +05:30
real ( pReal ) , dimension ( : , : ) , allocatable :: orientationOfGrain
2012-03-09 01:55:28 +05:30
real ( pReal ) , dimension ( 3 ) :: orientation
real ( pReal ) , dimension ( 3 , 3 ) :: symOrientation
2011-03-22 19:10:27 +05:30
integer ( pInt ) , dimension ( : ) , allocatable :: phaseOfGrain , textureOfGrain
2013-05-29 22:53:49 +05:30
integer ( pInt ) :: t , e , i , g , j , m , c , r , homog , micro , sgn , hme , myDebug , &
phaseID , textureID , dGrains , myNgrains , myNorientations , myNconstituents , &
grain , constituentGrain , ipGrain , symExtension , ip
2012-03-09 01:55:28 +05:30
real ( pReal ) :: extreme , rnd
2014-05-08 20:25:19 +05:30
integer ( pInt ) , dimension ( : , : ) , allocatable :: Nelems ! counts number of elements in homog, micro array
2013-02-25 18:43:52 +05:30
type ( p_intvec ) , dimension ( : , : ) , allocatable :: elemsOfHomogMicro ! lists element number in homog, micro array
2009-03-04 17:18:54 +05:30
2012-07-05 15:24:50 +05:30
myDebug = debug_level ( debug_material )
2012-03-09 01:55:28 +05:30
2013-12-13 04:33:37 +05:30
allocate ( material_volume ( homogenization_maxNgrains , mesh_maxNips , mesh_NcpElems ) , source = 0.0_pReal )
allocate ( material_phase ( homogenization_maxNgrains , mesh_maxNips , mesh_NcpElems ) , source = 0_pInt )
allocate ( material_texture ( homogenization_maxNgrains , mesh_maxNips , mesh_NcpElems ) , source = 0_pInt )
allocate ( material_EulerAngles ( 3 , homogenization_maxNgrains , mesh_maxNips , mesh_NcpElems ) , source = 0.0_pReal )
2009-03-04 17:18:54 +05:30
2013-12-13 04:33:37 +05:30
allocate ( Ngrains ( material_Nhomogenization , material_Nmicrostructure ) , source = 0_pInt )
allocate ( Nelems ( material_Nhomogenization , material_Nmicrostructure ) , source = 0_pInt )
2012-01-12 16:06:17 +05:30
2013-01-18 17:00:52 +05:30
!--------------------------------------------------------------------------------------------------
2012-01-12 22:01:23 +05:30
! precounting of elements for each homog/micro pair
2012-02-16 00:28:38 +05:30
do e = 1_pInt , mesh_NcpElems
2012-01-12 22:01:23 +05:30
homog = mesh_element ( 3 , e )
micro = mesh_element ( 4 , e )
Nelems ( homog , micro ) = Nelems ( homog , micro ) + 1_pInt
2012-01-12 16:06:17 +05:30
enddo
2013-02-25 18:43:52 +05:30
allocate ( elemsOfHomogMicro ( material_Nhomogenization , material_Nmicrostructure ) )
do homog = 1 , material_Nhomogenization
do micro = 1 , material_Nmicrostructure
if ( Nelems ( homog , micro ) > 0_pInt ) then
allocate ( elemsOfHomogMicro ( homog , micro ) % p ( Nelems ( homog , micro ) ) )
elemsOfHomogMicro ( homog , micro ) % p = 0_pInt
endif
enddo
enddo
2012-01-12 22:01:23 +05:30
2013-01-18 17:00:52 +05:30
!--------------------------------------------------------------------------------------------------
2012-01-12 22:01:23 +05:30
! identify maximum grain count per IP (from element) and find grains per homog/micro pair
2013-01-18 17:00:52 +05:30
Nelems = 0_pInt ! reuse as counter
2013-03-28 19:20:20 +05:30
elementLooping : do e = 1_pInt , mesh_NcpElems
2013-05-29 22:53:49 +05:30
t = FE_geomtype ( mesh_element ( 2 , e ) )
2009-03-04 17:18:54 +05:30
homog = mesh_element ( 3 , e )
micro = mesh_element ( 4 , e )
2013-01-18 17:00:52 +05:30
if ( homog < 1_pInt . or . homog > material_Nhomogenization ) & ! out of bounds
2012-02-13 23:11:27 +05:30
call IO_error ( 154_pInt , e , 0_pInt , 0_pInt )
2013-01-18 17:00:52 +05:30
if ( micro < 1_pInt . or . micro > material_Nmicrostructure ) & ! out of bounds
2012-02-13 23:11:27 +05:30
call IO_error ( 155_pInt , e , 0_pInt , 0_pInt )
2013-05-29 22:53:49 +05:30
if ( microstructure_elemhomo ( micro ) ) then ! how many grains are needed at this element?
dGrains = homogenization_Ngrains ( homog ) ! only one set of Ngrains (other IPs are plain copies)
2009-11-24 20:30:25 +05:30
else
2013-05-29 22:53:49 +05:30
dGrains = homogenization_Ngrains ( homog ) * FE_Nips ( t ) ! each IP has Ngrains
2009-11-24 20:30:25 +05:30
endif
2013-05-29 22:53:49 +05:30
Ngrains ( homog , micro ) = Ngrains ( homog , micro ) + dGrains ! total grain count
Nelems ( homog , micro ) = Nelems ( homog , micro ) + 1_pInt ! total element count
2013-02-25 18:43:52 +05:30
elemsOfHomogMicro ( homog , micro ) % p ( Nelems ( homog , micro ) ) = e ! remember elements active in this homog/micro pair
2013-03-28 19:20:20 +05:30
enddo elementLooping
2013-05-29 22:53:49 +05:30
2013-12-13 04:33:37 +05:30
allocate ( volumeOfGrain ( maxval ( Ngrains ) ) , source = 0.0_pReal ) ! reserve memory for maximum case
allocate ( phaseOfGrain ( maxval ( Ngrains ) ) , source = 0_pInt ) ! reserve memory for maximum case
allocate ( textureOfGrain ( maxval ( Ngrains ) ) , source = 0_pInt ) ! reserve memory for maximum case
allocate ( orientationOfGrain ( 3 , maxval ( Ngrains ) ) , source = 0.0_pReal ) ! reserve memory for maximum case
2009-03-04 17:18:54 +05:30
2012-03-09 01:55:28 +05:30
if ( iand ( myDebug , debug_levelBasic ) / = 0_pInt ) then
2011-03-21 16:01:17 +05:30
!$OMP CRITICAL (write2out)
2013-03-28 19:20:20 +05:30
write ( 6 , '(/,a/)' ) ' MATERIAL grain population'
2012-10-02 18:23:25 +05:30
write ( 6 , '(a32,1x,a32,1x,a6)' ) 'homogenization_name' , 'microstructure_name' , 'grain#'
2011-03-21 16:01:17 +05:30
!$OMP END CRITICAL (write2out)
endif
2013-01-18 17:00:52 +05:30
do homog = 1_pInt , material_Nhomogenization ! loop over homogenizations
dGrains = homogenization_Ngrains ( homog ) ! grain number per material point
do micro = 1_pInt , material_Nmicrostructure ! all pairs of homog and micro
if ( Ngrains ( homog , micro ) > 0_pInt ) then ! an active pair of homog and micro
myNgrains = Ngrains ( homog , micro ) ! assign short name for total number of grains to populate
2013-05-29 22:53:49 +05:30
myNconstituents = microstructure_Nconstituents ( micro ) ! assign short name for number of constituents
2012-03-09 01:55:28 +05:30
if ( iand ( myDebug , debug_levelBasic ) / = 0_pInt ) then
2011-03-21 16:01:17 +05:30
!$OMP CRITICAL (write2out)
2013-03-28 19:20:20 +05:30
write ( 6 , '(/,a32,1x,a32,1x,i6)' ) homogenization_name ( homog ) , microstructure_name ( micro ) , myNgrains
2011-03-21 16:01:17 +05:30
!$OMP END CRITICAL (write2out)
endif
2013-05-29 22:53:49 +05:30
2013-01-18 17:00:52 +05:30
!--------------------------------------------------------------------------------------------------
! calculate volume of each grain
2013-05-29 22:53:49 +05:30
2009-05-26 22:54:42 +05:30
volumeOfGrain = 0.0_pReal
2012-01-12 16:06:17 +05:30
grain = 0_pInt
2013-05-29 22:53:49 +05:30
2012-01-12 22:01:23 +05:30
do hme = 1_pInt , Nelems ( homog , micro )
2013-05-29 22:53:49 +05:30
e = elemsOfHomogMicro ( homog , micro ) % p ( hme ) ! my combination of homog and micro, only perform calculations for elements with homog, micro combinations which is indexed in cpElemsindex
2012-11-16 04:15:20 +05:30
t = FE_geomtype ( mesh_element ( 2 , e ) )
2013-01-18 17:00:52 +05:30
if ( microstructure_elemhomo ( micro ) ) then ! homogeneous distribution of grains over each element's IPs
2012-11-16 04:15:20 +05:30
volumeOfGrain ( grain + 1_pInt : grain + dGrains ) = sum ( mesh_ipVolume ( 1 : FE_Nips ( t ) , e ) ) / &
2013-05-29 22:53:49 +05:30
real ( dGrains , pReal ) ! each grain combines size of all IPs in that element
grain = grain + dGrains ! wind forward by Ngrains@IP
2012-01-12 16:06:17 +05:30
else
2013-01-18 17:00:52 +05:30
forall ( i = 1_pInt : FE_Nips ( t ) ) & ! loop over IPs
2012-02-16 00:28:38 +05:30
volumeOfGrain ( grain + ( i - 1 ) * dGrains + 1_pInt : grain + i * dGrains ) = &
2013-05-29 22:53:49 +05:30
mesh_ipVolume ( i , e ) / dGrains ! assign IPvolume/Ngrains@IP to all grains of IP
grain = grain + FE_Nips ( t ) * dGrains ! wind forward by Nips*Ngrains@IP
2009-03-04 17:18:54 +05:30
endif
enddo
2013-05-29 22:53:49 +05:30
if ( grain / = myNgrains ) &
2013-09-18 19:37:55 +05:30
call IO_error ( 0 , el = homog , ip = micro , ext_msg = 'inconsistent grain count after volume calc' )
2013-05-29 22:53:49 +05:30
2013-01-18 17:00:52 +05:30
!--------------------------------------------------------------------------------------------------
! divide myNgrains as best over constituents
2013-05-29 22:53:49 +05:30
!
! example: three constituents with fractions of 0.25, 0.25, and 0.5 distributed over 20 (microstructure) grains
!
! ***** ***** **********
! NgrainsOfConstituent: 5, 5, 10
! counters:
! |-----> grain (if constituent == 2)
! |--> constituentGrain (of constituent 2)
!
NgrainsOfConstituent = 0_pInt ! reset counter of grains per constituent
forall ( i = 1_pInt : myNconstituents ) &
2013-03-28 19:20:20 +05:30
NgrainsOfConstituent ( i ) = nint ( microstructure_fraction ( i , micro ) * myNgrains , pInt ) ! do rounding integer conversion
do while ( sum ( NgrainsOfConstituent ) / = myNgrains ) ! total grain count over constituents wrong?
sgn = sign ( 1_pInt , myNgrains - sum ( NgrainsOfConstituent ) ) ! direction of required change
2009-03-04 17:18:54 +05:30
extreme = 0.0_pReal
t = 0_pInt
2013-05-29 22:53:49 +05:30
do i = 1_pInt , myNconstituents ! find largest deviator
2012-02-10 17:26:05 +05:30
if ( real ( sgn , pReal ) * log ( NgrainsOfConstituent ( i ) / myNgrains / microstructure_fraction ( i , micro ) ) > extreme ) then
extreme = real ( sgn , pReal ) * log ( NgrainsOfConstituent ( i ) / myNgrains / microstructure_fraction ( i , micro ) )
2009-03-04 17:18:54 +05:30
t = i
endif
enddo
2013-03-28 19:20:20 +05:30
NgrainsOfConstituent ( t ) = NgrainsOfConstituent ( t ) + sgn ! change that by one
2009-06-15 18:41:21 +05:30
enddo
2013-01-18 17:00:52 +05:30
2013-05-29 22:53:49 +05:30
!--------------------------------------------------------------------------------------------------
! assign phase and texture info
2009-03-04 17:18:54 +05:30
phaseOfGrain = 0_pInt
2011-03-22 19:10:27 +05:30
textureOfGrain = 0_pInt
2009-03-04 17:18:54 +05:30
orientationOfGrain = 0.0_pReal
2013-05-29 22:53:49 +05:30
texture : do i = 1_pInt , myNconstituents ! loop over constituents
grain = sum ( NgrainsOfConstituent ( 1_pInt : i - 1_pInt ) ) ! set microstructure grain index of current constituent
! "grain" points to start of this constituent's grain population
constituentGrain = 0_pInt ! constituent grain index
2009-03-04 17:18:54 +05:30
phaseID = microstructure_phase ( i , micro )
textureID = microstructure_texture ( i , micro )
2013-05-29 22:53:49 +05:30
phaseOfGrain ( grain + 1_pInt : grain + NgrainsOfConstituent ( i ) ) = phaseID ! assign resp. phase
2013-03-28 19:20:20 +05:30
textureOfGrain ( grain + 1_pInt : grain + NgrainsOfConstituent ( i ) ) = textureID ! assign resp. texture
2009-03-04 17:18:54 +05:30
2012-02-10 17:26:05 +05:30
myNorientations = ceiling ( real ( NgrainsOfConstituent ( i ) , pReal ) / &
2013-03-28 19:20:20 +05:30
real ( texture_symmetry ( textureID ) , pReal ) , pInt ) ! max number of unique orientations (excl. symmetry)
2013-01-18 17:00:52 +05:30
!--------------------------------------------------------------------------------------------------
2013-05-29 22:53:49 +05:30
! ...has texture components
2013-01-18 17:00:52 +05:30
if ( texture_ODFfile ( textureID ) == '' ) then
2013-05-29 22:53:49 +05:30
gauss : do t = 1_pInt , texture_Ngauss ( textureID ) ! loop over Gauss components
2013-03-28 19:20:20 +05:30
do g = 1_pInt , int ( myNorientations * texture_Gauss ( 5 , t , textureID ) , pInt ) ! loop over required grain count
2009-03-04 17:18:54 +05:30
orientationOfGrain ( : , grain + constituentGrain + g ) = &
math_sampleGaussOri ( texture_Gauss ( 1 : 3 , t , textureID ) , &
texture_Gauss ( 4 , t , textureID ) )
enddo
2013-05-29 22:53:49 +05:30
constituentGrain = &
constituentGrain + int ( myNorientations * texture_Gauss ( 5 , t , textureID ) ) ! advance counter for grains of current constituent
enddo gauss
2009-03-04 17:18:54 +05:30
2013-05-29 22:53:49 +05:30
fiber : do t = 1_pInt , texture_Nfiber ( textureID ) ! loop over fiber components
2013-03-28 19:20:20 +05:30
do g = 1_pInt , int ( myNorientations * texture_Fiber ( 6 , t , textureID ) , pInt ) ! loop over required grain count
2009-03-04 17:18:54 +05:30
orientationOfGrain ( : , grain + constituentGrain + g ) = &
math_sampleFiberOri ( texture_Fiber ( 1 : 2 , t , textureID ) , &
texture_Fiber ( 3 : 4 , t , textureID ) , &
texture_Fiber ( 5 , t , textureID ) )
enddo
2013-05-29 22:53:49 +05:30
constituentGrain = &
constituentGrain + int ( myNorientations * texture_fiber ( 6 , t , textureID ) , pInt ) ! advance counter for grains of current constituent
enddo fiber
2009-03-04 17:18:54 +05:30
2013-05-29 22:53:49 +05:30
random : do constituentGrain = constituentGrain + 1_pInt , myNorientations ! fill remainder with random
orientationOfGrain ( : , grain + constituentGrain ) = math_sampleRandomOri ( )
enddo random
!--------------------------------------------------------------------------------------------------
! ...has hybrid IA
2013-01-18 17:00:52 +05:30
else
2013-05-29 22:53:49 +05:30
orientationOfGrain ( 1 : 3 , grain + 1_pInt : grain + myNorientations ) = &
IO_hybridIA ( myNorientations , texture_ODFfile ( textureID ) )
if ( all ( orientationOfGrain ( 1 : 3 , grain + 1_pInt ) == - 1.0_pReal ) ) call IO_error ( 156_pInt )
endif
2013-03-28 19:20:20 +05:30
2013-01-18 17:00:52 +05:30
!--------------------------------------------------------------------------------------------------
2013-07-24 16:39:39 +05:30
! ...texture transformation
2013-05-29 22:53:49 +05:30
do j = 1_pInt , myNorientations ! loop over each "real" orientation
orientationOfGrain ( 1 : 3 , grain + j ) = math_RtoEuler ( & ! translate back to Euler angles
math_mul33x33 ( & ! pre-multiply
math_EulertoR ( orientationOfGrain ( 1 : 3 , grain + j ) ) , & ! face-value orientation
2013-07-24 16:39:39 +05:30
texture_transformation ( 1 : 3 , 1 : 3 , textureID ) & ! and transformation matrix
2013-05-29 22:53:49 +05:30
) &
)
enddo
!--------------------------------------------------------------------------------------------------
! ...sample symmetry
2013-01-18 17:00:52 +05:30
2009-03-04 17:18:54 +05:30
symExtension = texture_symmetry ( textureID ) - 1_pInt
2013-05-29 22:53:49 +05:30
if ( symExtension > 0_pInt ) then ! sample symmetry (number of additional equivalent orientations)
constituentGrain = myNorientations ! start right after "real" orientations
2013-01-18 17:00:52 +05:30
do j = 1_pInt , myNorientations ! loop over each "real" orientation
2013-05-29 22:53:49 +05:30
symOrientation = math_symmetricEulers ( texture_symmetry ( textureID ) , &
orientationOfGrain ( 1 : 3 , grain + j ) ) ! get symmetric equivalents
e = min ( symExtension , NgrainsOfConstituent ( i ) - constituentGrain ) ! do not overshoot end of constituent grain array
2009-03-04 17:18:54 +05:30
if ( e > 0_pInt ) then
2013-05-29 22:53:49 +05:30
orientationOfGrain ( 1 : 3 , grain + constituentGrain + 1 : &
grain + constituentGrain + e ) = &
symOrientation ( 1 : 3 , 1 : e )
constituentGrain = constituentGrain + e ! remainder shrinks by e
2009-03-04 17:18:54 +05:30
endif
enddo
endif
2013-01-18 17:00:52 +05:30
2013-03-28 19:20:20 +05:30
!--------------------------------------------------------------------------------------------------
2013-05-29 22:53:49 +05:30
! shuffle grains within current constituent
do j = 1_pInt , NgrainsOfConstituent ( i ) - 1_pInt ! walk thru grains of current constituent
2011-03-22 19:10:27 +05:30
call random_number ( rnd )
2013-05-29 22:53:49 +05:30
t = nint ( rnd * ( NgrainsOfConstituent ( i ) - j ) + j + 0.5_pReal , pInt ) ! select a grain in remaining list
m = phaseOfGrain ( grain + t ) ! exchange current with random
phaseOfGrain ( grain + t ) = phaseOfGrain ( grain + j )
phaseOfGrain ( grain + j ) = m
m = textureOfGrain ( grain + t ) ! exchange current with random
textureOfGrain ( grain + t ) = textureOfGrain ( grain + j )
textureOfGrain ( grain + j ) = m
orientation = orientationOfGrain ( 1 : 3 , grain + t ) ! exchange current with random
orientationOfGrain ( 1 : 3 , grain + t ) = orientationOfGrain ( 1 : 3 , grain + j )
orientationOfGrain ( 1 : 3 , grain + j ) = orientation
2011-03-22 19:10:27 +05:30
enddo
2013-05-29 22:53:49 +05:30
enddo texture
2013-06-11 12:58:08 +05:30
!< @todo calc fraction after weighing with volumePerGrain, exchange in MC steps to improve result (humbug at the moment)
2013-05-29 22:53:49 +05:30
!--------------------------------------------------------------------------------------------------
! distribute grains of all constituents as accurately as possible to given constituent fractions
ip = 0_pInt
currentGrainOfConstituent = 0_pInt
2012-01-12 22:01:23 +05:30
do hme = 1_pInt , Nelems ( homog , micro )
2013-02-25 18:43:52 +05:30
e = elemsOfHomogMicro ( homog , micro ) % p ( hme ) ! only perform calculations for elements with homog, micro combinations which is indexed in cpElemsindex
2012-11-21 22:27:57 +05:30
t = FE_geomtype ( mesh_element ( 2 , e ) )
2013-01-18 17:00:52 +05:30
if ( microstructure_elemhomo ( micro ) ) then ! homogeneous distribution of grains over each element's IPs
2013-05-29 22:53:49 +05:30
m = 1_pInt ! process only first IP
2012-01-12 16:06:17 +05:30
else
2013-05-29 22:53:49 +05:30
m = FE_Nips ( t ) ! process all IPs
2009-03-04 17:18:54 +05:30
endif
2013-05-29 22:53:49 +05:30
do i = 1_pInt , m ! loop over necessary IPs
ip = ip + 1_pInt ! keep track of total ip count
ipGrain = 0_pInt ! count number of grains assigned at this IP
randomOrder = math_range ( microstructure_maxNconstituents ) ! start out with ordered sequence of constituents
call random_number ( rndArray ) ! as many rnd numbers as (max) constituents
do j = 1_pInt , myNconstituents - 1_pInt ! loop over constituents ...
r = nint ( rndArray ( j ) * ( myNconstituents - j ) + j + 0.5_pReal , pInt ) ! ... select one in remaining list
c = randomOrder ( r ) ! ... call it "c"
randomOrder ( r ) = randomOrder ( j ) ! ... and exchange with present position in constituent list
grain = sum ( NgrainsOfConstituent ( 1 : c - 1_pInt ) ) ! figure out actual starting index in overall/consecutive grain population
do g = 1_pInt , min ( dGrains - ipGrain , & ! leftover number of grains at this IP
max ( 0_pInt , & ! no negative values
nint ( real ( ip * dGrains * NgrainsOfConstituent ( c ) ) / & ! fraction of grains scaled to this constituent...
real ( myNgrains ) , pInt ) - & ! ...minus those already distributed
currentGrainOfConstituent ( c ) ) )
ipGrain = ipGrain + 1_pInt ! advance IP grain counter
currentGrainOfConstituent ( c ) = currentGrainOfConstituent ( c ) + 1_pInt ! advance index of grain population for constituent c
material_volume ( ipGrain , i , e ) = volumeOfGrain ( grain + currentGrainOfConstituent ( c ) ) ! assign properties
material_phase ( ipGrain , i , e ) = phaseOfGrain ( grain + currentGrainOfConstituent ( c ) )
material_texture ( ipGrain , i , e ) = textureOfGrain ( grain + currentGrainOfConstituent ( c ) )
material_EulerAngles ( 1 : 3 , ipGrain , i , e ) = orientationOfGrain ( 1 : 3 , grain + currentGrainOfConstituent ( c ) )
enddo ; enddo
c = randomOrder ( microstructure_Nconstituents ( micro ) ) ! look up constituent remaining after random shuffling
grain = sum ( NgrainsOfConstituent ( 1 : c - 1_pInt ) ) ! figure out actual starting index in overall/consecutive grain population
do ipGrain = ipGrain + 1_pInt , dGrains ! ensure last constituent fills up to dGrains
currentGrainOfConstituent ( c ) = currentGrainOfConstituent ( c ) + 1_pInt
material_volume ( ipGrain , i , e ) = volumeOfGrain ( grain + currentGrainOfConstituent ( c ) )
material_phase ( ipGrain , i , e ) = phaseOfGrain ( grain + currentGrainOfConstituent ( c ) )
material_texture ( ipGrain , i , e ) = textureOfGrain ( grain + currentGrainOfConstituent ( c ) )
material_EulerAngles ( 1 : 3 , ipGrain , i , e ) = orientationOfGrain ( 1 : 3 , grain + currentGrainOfConstituent ( c ) )
enddo
enddo
do i = i , FE_Nips ( t ) ! loop over IPs to (possibly) distribute copies from first IP
material_volume ( 1_pInt : dGrains , i , e ) = material_volume ( 1_pInt : dGrains , 1 , e )
material_phase ( 1_pInt : dGrains , i , e ) = material_phase ( 1_pInt : dGrains , 1 , e )
material_texture ( 1_pInt : dGrains , i , e ) = material_texture ( 1_pInt : dGrains , 1 , e )
material_EulerAngles ( 1 : 3 , 1_pInt : dGrains , i , e ) = material_EulerAngles ( 1 : 3 , 1_pInt : dGrains , 1 , e )
enddo
2009-03-04 17:18:54 +05:30
enddo
2013-01-18 17:00:52 +05:30
endif ! active homog,micro pair
2009-03-04 17:18:54 +05:30
enddo
enddo
2009-05-26 22:54:42 +05:30
deallocate ( volumeOfGrain )
2009-03-04 17:18:54 +05:30
deallocate ( phaseOfGrain )
2011-03-22 19:10:27 +05:30
deallocate ( textureOfGrain )
2009-03-04 17:18:54 +05:30
deallocate ( orientationOfGrain )
2012-01-12 22:01:23 +05:30
deallocate ( Nelems )
2013-06-11 22:05:04 +05:30
!> @todo - causing segmentation fault: needs looking into
2013-02-25 18:43:52 +05:30
!do homog = 1,material_Nhomogenization
! do micro = 1,material_Nmicrostructure
2013-05-29 22:53:49 +05:30
! if (Nelems(homog,micro) > 0_pInt) deallocate(elemsOfHomogMicro(homog,micro)%p)
2013-02-25 18:43:52 +05:30
! enddo
!enddo
2012-01-12 22:01:23 +05:30
deallocate ( elemsOfHomogMicro )
2009-03-04 17:18:54 +05:30
2012-03-09 01:55:28 +05:30
end subroutine material_populateGrains
2009-03-04 17:18:54 +05:30
2014-05-08 20:25:19 +05:30
#ifdef HDF
2014-03-25 22:51:47 +05:30
integer ( pInt ) pure function material_NconstituentsPhase ( matID )
implicit none
integer ( pInt ) , intent ( in ) :: matID
material_NconstituentsPhase = count ( microstructure_phase == matID )
end function
2014-05-08 20:25:19 +05:30
#endif
2014-03-25 22:51:47 +05:30
2013-02-20 03:42:05 +05:30
end module material