2013-03-22 23:05:05 +05:30
! Copyright 2011-13 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/>.
!
2012-03-09 01:55:28 +05:30
!--------------------------------------------------------------------------------------------------
2009-08-31 20:39:15 +05:30
!* $Id$
2012-03-09 01:55:28 +05:30
!--------------------------------------------------------------------------------------------------
!> @author Franz Roters, Max-Planck-Institut für Eisenforschung GmbH
!! 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 , &
p_intvec
2012-03-09 01:55:28 +05:30
implicit none
private
character ( len = 64 ) , parameter , public :: &
2013-03-28 19:20:20 +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
character ( len = 32 ) , parameter , public :: &
2013-03-28 19:20:20 +05:30
material_PARTHOMOGENIZATION = 'homogenization' , & !< keyword for homogenization part
material_PARTCRYSTALLITE = 'crystallite' , & !< keyword for crystallite part
material_PARTPHASE = 'phase' !< keyword for phase part
2012-03-09 01:55:28 +05:30
2013-01-18 17:00:52 +05:30
character ( len = 64 ) , dimension ( : ) , allocatable , public , protected :: &
2013-02-20 03:42:05 +05:30
phase_elasticity , & !< elasticity of each phase
phase_plasticity , & !< plasticity of each phase
phase_name , & !< name of each phase
homogenization_name , & !< name of each homogenization
homogenization_type , & !< type 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-01-18 17:00:52 +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
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
2013-01-18 17:00:52 +05:30
integer ( pInt ) , dimension ( : , : , : ) , allocatable , public :: &
2013-02-20 03:42:05 +05:30
material_phase !< phase (index) of each grain,IP,element
2013-01-18 17:00:52 +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-01-18 17:00:52 +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-01-18 17:00:52 +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
character ( len = 32 ) , parameter , private :: &
2013-03-28 19:20:20 +05:30
material_PARTMICROSTRUCTURE = 'microstructure' , & !< keyword for microstructure part
material_PARTTEXTURE = 'texture' !< keyword for texture part
2012-03-09 01:55:28 +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
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
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
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
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
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-01-18 17:00:52 +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
texture_Fiber !< data of each Fiber component
2012-03-09 01:55:28 +05:30
logical , dimension ( : ) , allocatable , private :: &
homogenization_active
2013-01-18 17:00:52 +05:30
public :: material_init
2012-03-09 01:55:28 +05:30
2013-01-18 17:00:52 +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
2012-03-09 01:55:28 +05:30
2009-03-04 17:18:54 +05:30
implicit none
2012-02-16 00:28:38 +05:30
integer ( pInt ) , parameter :: fileunit = 200_pInt
2013-03-28 19:20:20 +05:30
integer ( pInt ) :: m , c , h , myDebug
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$'
2013-02-25 22:04:59 +05:30
write ( 6 , '(a16,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-03-28 19:20:20 +05:30
if ( . not . IO_open_jobFile_stat ( fileunit , material_localFileExt ) ) then ! no local material configuration present...
call IO_open_file ( fileunit , material_configFile ) ! ...open material.config file
2011-08-02 15:44:16 +05:30
endif
2009-03-20 20:04:24 +05:30
call material_parseHomogenization ( fileunit , material_partHomogenization )
2012-03-09 01:55:28 +05:30
if ( iand ( myDebug , debug_levelBasic ) / = 0_pInt ) then
2013-03-28 19:20:20 +05:30
write ( 6 , '(a)' ) ' Homogenization parsed'
2011-10-18 14:48:05 +05:30
endif
2009-03-20 20:04:24 +05:30
call material_parseMicrostructure ( fileunit , material_partMicrostructure )
2012-03-09 01:55:28 +05:30
if ( iand ( myDebug , debug_levelBasic ) / = 0_pInt ) then
2013-03-28 19:20:20 +05:30
write ( 6 , '(a)' ) ' Microstructure parsed'
2011-10-18 14:48:05 +05:30
endif
2010-02-25 23:09:11 +05:30
call material_parseCrystallite ( fileunit , material_partCrystallite )
2012-03-09 01:55:28 +05:30
if ( iand ( myDebug , debug_levelBasic ) / = 0_pInt ) then
2013-03-28 19:20:20 +05:30
write ( 6 , '(a)' ) ' Crystallite parsed'
2011-10-18 14:48:05 +05:30
endif
2009-03-20 20:04:24 +05:30
call material_parseTexture ( fileunit , material_partTexture )
2012-03-09 01:55:28 +05:30
if ( iand ( myDebug , debug_levelBasic ) / = 0_pInt ) then
2013-03-28 19:20:20 +05:30
write ( 6 , '(a)' ) ' Texture parsed'
2011-10-18 14:48:05 +05:30
endif
2009-03-20 20:04:24 +05:30
call material_parsePhase ( fileunit , material_partPhase )
2012-03-09 01:55:28 +05:30
if ( iand ( myDebug , debug_levelBasic ) / = 0_pInt ) then
2013-03-28 19:20:20 +05:30
write ( 6 , '(a)' ) ' Phase parsed'
2011-10-18 14:48:05 +05:30
endif
2009-03-04 17:18:54 +05:30
close ( fileunit )
2013-03-28 19:20:20 +05:30
do m = 1_pInt , material_Nmicrostructure
if ( microstructure_crystallite ( m ) < 1_pInt . or . &
microstructure_crystallite ( m ) > material_Ncrystallite ) &
call IO_error ( 150_pInt , m )
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 ( 151_pInt , m )
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 ( 152_pInt , m )
if ( abs ( sum ( microstructure_fraction ( : , m ) ) - 1.0_pReal ) > = 1.0e-10_pReal ) then
2012-03-09 01:55:28 +05:30
if ( iand ( myDebug , debug_levelExtensive ) / = 0_pInt ) then
2013-03-28 19:20:20 +05:30
write ( 6 , '(a)' ) ' sum of microstructure fraction = ' , sum ( microstructure_fraction ( : , m ) )
2011-03-21 16:01:17 +05:30
endif
2013-03-28 19:20:20 +05:30
call IO_error ( 153_pInt , m )
2009-12-02 18:38:14 +05:30
endif
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
write ( 6 , '(1x,a32,1x,a16,1x,i4)' ) homogenization_name ( h ) , homogenization_type ( h ) , homogenization_Ngrains ( h )
2013-01-18 17:00:52 +05:30
enddo
2013-03-28 19:20:20 +05:30
write ( 6 , '(/,a32,19x,a11,1x,a12,1x,a13)' ) 'microstructure' , 'crystallite' , 'constituents' , 'homogeneous'
do m = 1_pInt , material_Nmicrostructure
write ( 6 , '(a32,4x,i4,8x,i4,8x,l1)' ) microstructure_name ( m ) , &
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
!--------------------------------------------------------------------------------------------------
2012-02-16 00:28:38 +05:30
subroutine material_parseHomogenization ( myFile , myPart )
2013-03-28 19:20:20 +05:30
use IO , only : &
IO_globalTagInPart , &
IO_countSections , &
IO_error , &
IO_countTagInPart , &
IO_lc , &
IO_getTag , &
IO_isBlank , &
IO_stringValue , &
IO_intValue , &
IO_stringPos
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
2012-03-09 01:55:28 +05:30
integer ( pInt ) , intent ( in ) :: myFile
integer ( pInt ) , parameter :: maxNchunks = 2_pInt
2009-03-04 17:18:54 +05:30
integer ( pInt ) , dimension ( 1 + 2 * maxNchunks ) :: positions
2009-03-31 14:54:12 +05:30
integer ( pInt ) Nsections , section , s
2012-03-09 01:55:28 +05:30
character ( len = 64 ) :: tag
2012-06-26 15:54:54 +05:30
character ( len = 1024 ) :: line
logical :: echo
echo = IO_globalTagInPart ( myFile , myPart , '/echo/' )
2009-03-04 17:18:54 +05:30
2012-02-16 00:28:38 +05:30
Nsections = IO_countSections ( myFile , 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-03-28 19:20:20 +05:30
allocate ( homogenization_name ( Nsections ) ) ; homogenization_name = ''
allocate ( homogenization_type ( Nsections ) ) ; homogenization_type = ''
2009-03-31 14:54:12 +05:30
allocate ( homogenization_typeInstance ( Nsections ) ) ; homogenization_typeInstance = 0_pInt
2013-03-28 19:20:20 +05:30
allocate ( homogenization_Ngrains ( Nsections ) ) ; homogenization_Ngrains = 0_pInt
allocate ( homogenization_Noutput ( Nsections ) ) ; homogenization_Noutput = 0_pInt
allocate ( homogenization_active ( Nsections ) ) ; homogenization_active = . 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-03-28 19:20:20 +05:30
homogenization_Noutput = IO_countTagInPart ( myFile , myPart , '(output)' , Nsections )
2009-03-31 14:54:12 +05:30
2012-02-16 00:28:38 +05:30
rewind ( myFile )
2009-03-04 17:18:54 +05:30
line = ''
2012-02-16 00:28:38 +05:30
section = 0_pInt
2009-03-04 17:18:54 +05:30
2013-01-18 17:00:52 +05:30
do while ( IO_lc ( IO_getTag ( line , '<' , '>' ) ) / = myPart ) ! wind forward to myPart
2012-02-16 00:28:38 +05:30
read ( myFile , '(a1024)' , END = 100 ) line
2009-03-04 17:18:54 +05:30
enddo
2013-03-28 19:20:20 +05:30
if ( echo ) write ( 6 , '(a)' ) trim ( line ) ! echo part header
2009-03-04 17:18:54 +05:30
do
2012-02-16 00:28:38 +05:30
read ( myFile , '(a1024)' , END = 100 ) line
2013-01-18 17:00:52 +05:30
if ( IO_isBlank ( line ) ) cycle ! skip empty lines
if ( IO_getTag ( line , '<' , '>' ) / = '' ) exit ! stop at next part
if ( echo ) write ( 6 , * ) trim ( line ) ! echo back read lines
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
homogenization_name ( section ) = IO_getTag ( line , '[' , ']' )
endif
2012-02-10 17:26:05 +05:30
if ( section > 0_pInt ) then
2009-03-04 17:18:54 +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-01-18 17:00:52 +05:30
homogenization_type ( section ) = IO_lc ( IO_stringValue ( line , positions , 2_pInt ) ) ! adding: IO_lc function
2012-02-10 17:26:05 +05:30
do s = 1_pInt , section
2009-03-31 14:54:12 +05:30
if ( homogenization_type ( s ) == homogenization_type ( section ) ) &
2013-01-18 17:00:52 +05:30
homogenization_typeInstance ( section ) = homogenization_typeInstance ( section ) + 1_pInt ! count instances
2009-03-31 14:54:12 +05:30
enddo
2009-03-04 17:18:54 +05:30
case ( '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
2009-07-22 21:37:19 +05:30
100 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
!--------------------------------------------------------------------------------------------------
2012-02-16 00:28:38 +05:30
subroutine material_parseMicrostructure ( myFile , 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
2012-03-09 01:55:28 +05:30
integer ( pInt ) , intent ( in ) :: myFile
2012-02-10 17:26:05 +05:30
integer ( pInt ) , parameter :: maxNchunks = 7_pInt
2012-03-09 01:55:28 +05:30
2012-02-10 17:26:05 +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
character ( len = 64 ) :: tag
character ( len = 1024 ) :: line
2012-06-26 15:54:54 +05:30
logical :: echo
echo = IO_globalTagInPart ( myFile , myPart , '/echo/' )
2011-09-13 21:24:06 +05:30
2012-02-16 00:28:38 +05:30
Nsections = IO_countSections ( myFile , 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
allocate ( microstructure_name ( Nsections ) ) ; microstructure_name = ''
allocate ( microstructure_crystallite ( Nsections ) ) ; microstructure_crystallite = 0_pInt
2009-03-04 17:18:54 +05:30
allocate ( microstructure_Nconstituents ( Nsections ) )
2009-07-22 21:37:19 +05:30
allocate ( microstructure_active ( Nsections ) )
2009-11-24 20:30:25 +05:30
allocate ( microstructure_elemhomo ( Nsections ) )
2009-07-22 21:37:19 +05:30
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
2012-02-16 00:28:38 +05:30
microstructure_Nconstituents = IO_countTagInPart ( myFile , myPart , '(constituent)' , Nsections )
2009-03-04 17:18:54 +05:30
microstructure_maxNconstituents = maxval ( microstructure_Nconstituents )
2012-02-16 00:28:38 +05:30
microstructure_elemhomo = IO_spotTagInPart ( myFile , myPart , '/elementhomogeneous/' , Nsections )
2009-11-24 20:30:25 +05:30
2013-01-18 17:00:52 +05:30
allocate ( microstructure_phase ( microstructure_maxNconstituents , Nsections ) )
microstructure_phase = 0_pInt
allocate ( microstructure_texture ( microstructure_maxNconstituents , Nsections ) )
microstructure_texture = 0_pInt
allocate ( microstructure_fraction ( microstructure_maxNconstituents , Nsections ) )
microstructure_fraction = 0.0_pReal
2009-03-04 17:18:54 +05:30
2012-02-16 00:28:38 +05:30
rewind ( myFile )
2013-03-28 19:20:20 +05:30
line = '' ! to have it initialized
section = 0_pInt ! - " -
constituent = 0_pInt ! - " -
2009-03-04 17:18:54 +05:30
2013-01-18 17:00:52 +05:30
do while ( IO_lc ( IO_getTag ( line , '<' , '>' ) ) / = myPart ) ! wind forward to myPart
2012-02-16 00:28:38 +05:30
read ( myFile , '(a1024)' , END = 100 ) line
2009-03-04 17:18:54 +05:30
enddo
2013-01-18 17:00:52 +05:30
if ( echo ) write ( 6 , * ) trim ( line ) ! echo part header
2009-03-04 17:18:54 +05:30
do
2012-02-16 00:28:38 +05:30
read ( myFile , '(a1024)' , END = 100 ) line
2013-01-18 17:00:52 +05:30
if ( IO_isBlank ( line ) ) cycle ! skip empty lines
if ( IO_getTag ( line , '<' , '>' ) / = '' ) exit ! stop at next part
if ( echo ) write ( 6 , * ) trim ( line ) ! echo back read lines
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
2009-03-04 17:18:54 +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
2012-03-09 01:55:28 +05:30
100 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
!--------------------------------------------------------------------------------------------------
2012-02-16 00:28:38 +05:30
subroutine material_parseCrystallite ( myFile , myPart )
2013-01-18 17:00:52 +05:30
use IO , only : &
IO_countSections , &
IO_error , &
IO_countTagInPart , &
IO_globalTagInPart , &
IO_getTag , &
IO_lc , &
IO_isBlank
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
2012-03-09 01:55:28 +05:30
integer ( pInt ) , intent ( in ) :: myFile
integer ( pInt ) :: Nsections , &
section
character ( len = 1024 ) :: line
2012-06-26 15:54:54 +05:30
logical :: echo
echo = IO_globalTagInPart ( myFile , myPart , '/echo/' )
2010-02-25 23:09:11 +05:30
2012-02-16 00:28:38 +05:30
Nsections = IO_countSections ( myFile , 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
allocate ( crystallite_name ( Nsections ) ) ; crystallite_name = ''
allocate ( crystallite_Noutput ( Nsections ) ) ; crystallite_Noutput = 0_pInt
2012-02-16 00:28:38 +05:30
crystallite_Noutput = IO_countTagInPart ( myFile , myPart , '(output)' , Nsections )
2010-02-25 23:09:11 +05:30
2012-02-16 00:28:38 +05:30
rewind ( myFile )
2010-02-25 23:09:11 +05:30
line = ''
2012-02-16 00:28:38 +05:30
section = 0_pInt
2010-02-25 23:09:11 +05:30
2013-01-18 17:00:52 +05:30
do while ( IO_lc ( IO_getTag ( line , '<' , '>' ) ) / = myPart ) ! wind forward to myPart
2012-02-16 00:28:38 +05:30
read ( myFile , '(a1024)' , END = 100 ) line
2010-02-25 23:09:11 +05:30
enddo
2013-01-18 17:00:52 +05:30
if ( echo ) write ( 6 , * ) trim ( line ) ! echo part header
2010-02-25 23:09:11 +05:30
do
2012-02-16 00:28:38 +05:30
read ( myFile , '(a1024)' , END = 100 ) line
2013-01-18 17:00:52 +05:30
if ( IO_isBlank ( line ) ) cycle ! skip empty lines
if ( IO_getTag ( line , '<' , '>' ) / = '' ) exit ! stop at next part
if ( echo ) write ( 6 , * ) trim ( line ) ! echo back read lines
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
2012-03-09 01:55:28 +05:30
100 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
!--------------------------------------------------------------------------------------------------
2012-02-16 00:28:38 +05:30
subroutine material_parsePhase ( myFile , myPart )
2013-03-28 19:20:20 +05:30
use IO , only : &
IO_globalTagInPart , &
IO_countSections , &
IO_error , &
IO_countTagInPart , &
IO_getTag , &
IO_spotTagInPart , &
IO_lc , &
IO_isBlank , &
IO_stringValue , &
IO_stringPos
2009-03-04 17:18:54 +05:30
implicit none
2009-03-20 20:04:24 +05:30
character ( len = * ) , intent ( in ) :: myPart
2012-03-09 01:55:28 +05:30
integer ( pInt ) , intent ( in ) :: myFile
2012-02-16 00:28:38 +05:30
integer ( pInt ) , parameter :: maxNchunks = 2_pInt
2012-03-09 01:55:28 +05:30
2009-03-04 17:18:54 +05:30
integer ( pInt ) , dimension ( 1 + 2 * maxNchunks ) :: positions
integer ( pInt ) Nsections , section , s
2012-03-09 01:55:28 +05:30
character ( len = 64 ) :: tag
character ( len = 1024 ) :: line
2012-06-26 15:54:54 +05:30
logical :: echo
echo = IO_globalTagInPart ( myFile , myPart , '/echo/' )
2009-03-04 17:18:54 +05:30
2012-02-16 00:28:38 +05:30
Nsections = IO_countSections ( myFile , 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 = ''
allocate ( phase_elasticity ( Nsections ) ) ; phase_elasticity = ''
2012-03-14 21:46:11 +05:30
allocate ( phase_elasticityInstance ( Nsections ) ) ; phase_elasticityInstance = 0_pInt
2013-03-28 19:20:20 +05:30
allocate ( phase_plasticity ( Nsections ) ) ; phase_plasticity = ''
2012-03-12 19:39:37 +05:30
allocate ( phase_plasticityInstance ( Nsections ) ) ; phase_plasticityInstance = 0_pInt
2013-03-28 19:20:20 +05:30
allocate ( phase_Noutput ( Nsections ) ) ; phase_Noutput = 0_pInt
allocate ( phase_localPlasticity ( Nsections ) ) ; phase_localPlasticity = . false .
2009-03-04 17:18:54 +05:30
2012-02-16 00:28:38 +05:30
phase_Noutput = IO_countTagInPart ( myFile , myPart , '(output)' , Nsections )
2012-03-12 19:39:37 +05:30
phase_localPlasticity = . not . IO_spotTagInPart ( myFile , myPart , '/nonlocal/' , Nsections )
2009-03-04 17:18:54 +05:30
2012-02-16 00:28:38 +05:30
rewind ( myFile )
2009-03-04 17:18:54 +05:30
line = ''
2012-02-16 00:28:38 +05:30
section = 0_pInt
2009-03-04 17:18:54 +05:30
2013-01-18 17:00:52 +05:30
do while ( IO_lc ( IO_getTag ( line , '<' , '>' ) ) / = myPart ) ! wind forward to myPart
2012-02-16 00:28:38 +05:30
read ( myFile , '(a1024)' , END = 100 ) line
2009-03-04 17:18:54 +05:30
enddo
2013-01-18 17:00:52 +05:30
if ( echo ) write ( 6 , * ) trim ( line ) ! echo part header
2009-03-04 17:18:54 +05:30
do
2012-02-16 00:28:38 +05:30
read ( myFile , '(a1024)' , END = 100 ) line
2013-01-18 17:00:52 +05:30
if ( IO_isBlank ( line ) ) cycle ! skip empty lines
if ( IO_getTag ( line , '<' , '>' ) / = '' ) exit ! stop at next part
if ( echo ) write ( 6 , * ) trim ( line ) ! echo back read lines
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
2009-03-04 17:18:54 +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' )
phase_elasticity ( section ) = IO_lc ( IO_stringValue ( line , positions , 2_pInt ) )
do s = 1_pInt , section
if ( phase_elasticity ( s ) == phase_elasticity ( section ) ) &
2013-01-18 17:00:52 +05:30
phase_elasticityInstance ( section ) = phase_elasticityInstance ( section ) + 1_pInt ! count instances
2012-03-14 21:46:11 +05:30
enddo
2012-03-12 19:39:37 +05:30
case ( 'plasticity' )
phase_plasticity ( section ) = IO_lc ( IO_stringValue ( line , positions , 2_pInt ) )
2012-02-16 00:28:38 +05:30
do s = 1_pInt , section
2012-03-12 19:39:37 +05:30
if ( phase_plasticity ( s ) == phase_plasticity ( section ) ) &
2013-01-18 17:00:52 +05:30
phase_plasticityInstance ( section ) = phase_plasticityInstance ( section ) + 1_pInt ! count instances
2009-03-04 17:18:54 +05:30
enddo
end select
endif
enddo
2012-03-09 01:55:28 +05:30
100 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
!--------------------------------------------------------------------------------------------------
2012-02-16 00:28:38 +05:30
subroutine material_parseTexture ( myFile , myPart )
2013-03-28 19:20:20 +05:30
use IO , only : &
IO_globalTagInPart , &
IO_countSections , &
IO_error , &
IO_countTagInPart , &
IO_getTag , &
IO_spotTagInPart , &
IO_lc , &
IO_isBlank , &
IO_floatValue , &
IO_stringValue , &
IO_stringPos
2013-01-18 17:00:52 +05:30
use math , only : &
inRad , &
math_sampleRandomOri
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
2012-03-09 01:55:28 +05:30
integer ( pInt ) , intent ( in ) :: myFile
integer ( pInt ) , parameter :: maxNchunks = 13_pInt
2009-03-04 17:18:54 +05:30
integer ( pInt ) , dimension ( 1 + 2 * maxNchunks ) :: positions
2013-03-28 19:20:20 +05:30
integer ( pInt ) :: Nsections , section , gauss , fiber , j
2012-03-09 01:55:28 +05:30
character ( len = 64 ) :: tag
character ( len = 1024 ) :: line
2012-06-26 15:54:54 +05:30
logical :: echo
echo = IO_globalTagInPart ( myFile , myPart , '/echo/' )
2009-03-04 17:18:54 +05:30
2012-02-16 00:28:38 +05:30
Nsections = IO_countSections ( myFile , 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
2009-03-04 17:18:54 +05:30
allocate ( texture_name ( Nsections ) ) ; texture_name = ''
allocate ( texture_ODFfile ( Nsections ) ) ; texture_ODFfile = ''
allocate ( texture_symmetry ( Nsections ) ) ; texture_symmetry = 1_pInt
allocate ( texture_Ngauss ( Nsections ) ) ; texture_Ngauss = 0_pInt
allocate ( texture_Nfiber ( Nsections ) ) ; texture_Nfiber = 0_pInt
2012-02-16 00:28:38 +05:30
texture_Ngauss = IO_countTagInPart ( myFile , myPart , '(gauss)' , Nsections ) + &
IO_countTagInPart ( myFile , myPart , '(random)' , Nsections )
texture_Nfiber = IO_countTagInPart ( myFile , myPart , '(fiber)' , Nsections )
2009-03-04 17:18:54 +05:30
texture_maxNgauss = maxval ( texture_Ngauss )
texture_maxNfiber = maxval ( texture_Nfiber )
allocate ( texture_Gauss ( 5 , texture_maxNgauss , Nsections ) ) ; texture_Gauss = 0.0_pReal
allocate ( texture_Fiber ( 6 , texture_maxNfiber , Nsections ) ) ; texture_Fiber = 0.0_pReal
2012-02-16 00:28:38 +05:30
rewind ( myFile )
2013-03-28 19:20:20 +05:30
line = '' ! to have in initialized
section = 0_pInt ! - " -
gauss = 0_pInt ! - " -
fiber = 0_pInt ! - " -
2009-03-04 17:18:54 +05:30
2013-01-18 17:00:52 +05:30
do while ( IO_lc ( IO_getTag ( line , '<' , '>' ) ) / = myPart ) ! wind forward to myPart
2012-02-16 00:28:38 +05:30
read ( myFile , '(a1024)' , END = 100 ) line
2009-03-04 17:18:54 +05:30
enddo
2013-03-28 19:20:20 +05:30
if ( echo ) write ( 6 , '(a)' ) trim ( line ) ! echo part header
2009-03-04 17:18:54 +05:30
do
2012-02-16 00:28:38 +05:30
read ( myFile , '(a1024)' , END = 100 ) line
2013-01-18 17:00:52 +05:30
if ( IO_isBlank ( line ) ) cycle ! skip empty lines
if ( IO_getTag ( line , '<' , '>' ) / = '' ) exit ! stop at next part
2013-03-28 19:20:20 +05:30
if ( echo ) write ( 6 , '(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
2009-03-04 17:18:54 +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-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
2012-03-09 01:55:28 +05:30
100 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 : &
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
2012-03-09 01:55:28 +05:30
integer ( pInt ) , dimension ( microstructure_maxNconstituents ) &
:: NgrainsOfConstituent
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
2012-03-09 01:55:28 +05:30
integer ( pInt ) :: t , e , i , g , j , m , homog , micro , sgn , hme , myDebug
integer ( pInt ) :: phaseID , textureID , dGrains , myNgrains , myNorientations , &
2012-01-12 22:01:23 +05:30
grain , constituentGrain , symExtension
2012-03-09 01:55:28 +05:30
real ( pReal ) :: extreme , rnd
2013-01-18 17:00:52 +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
2011-03-22 19:10:27 +05:30
allocate ( material_volume ( homogenization_maxNgrains , mesh_maxNips , mesh_NcpElems ) ) ; material_volume = 0.0_pReal
allocate ( material_phase ( homogenization_maxNgrains , mesh_maxNips , mesh_NcpElems ) ) ; material_phase = 0_pInt
allocate ( material_texture ( homogenization_maxNgrains , mesh_maxNips , mesh_NcpElems ) ) ; material_texture = 0_pInt
2009-03-04 17:18:54 +05:30
allocate ( material_EulerAngles ( 3 , homogenization_maxNgrains , mesh_maxNips , mesh_NcpElems ) ) ; material_EulerAngles = 0.0_pReal
allocate ( Ngrains ( material_Nhomogenization , material_Nmicrostructure ) ) ; Ngrains = 0_pInt
2012-01-12 22:01:23 +05:30
allocate ( Nelems ( material_Nhomogenization , material_Nmicrostructure ) ) ; Nelems = 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
2012-11-16 04:15:20 +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 )
2009-11-24 20:30:25 +05:30
if ( microstructure_elemhomo ( micro ) ) then
dGrains = homogenization_Ngrains ( homog )
else
2012-11-16 04:15:20 +05:30
dGrains = homogenization_Ngrains ( homog ) * FE_Nips ( t )
2009-11-24 20:30:25 +05:30
endif
Ngrains ( homog , micro ) = Ngrains ( homog , micro ) + dGrains
2012-01-12 22:01:23 +05:30
Nelems ( homog , micro ) = Nelems ( homog , micro ) + 1_pInt
2013-02-25 18:43:52 +05:30
elemsOfHomogMicro ( homog , micro ) % p ( Nelems ( homog , micro ) ) = e ! remember elements active in this homog/micro pair
2012-01-12 16:06:17 +05:30
2013-03-28 19:20:20 +05:30
enddo elementLooping
2013-01-18 17:00:52 +05:30
allocate ( volumeOfGrain ( maxval ( Ngrains ) ) ) ! reserve memory for maximum case
allocate ( phaseOfGrain ( maxval ( Ngrains ) ) ) ! reserve memory for maximum case
allocate ( textureOfGrain ( maxval ( Ngrains ) ) ) ! reserve memory for maximum case
allocate ( orientationOfGrain ( 3 , maxval ( Ngrains ) ) ) ! 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
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
2009-03-04 17:18:54 +05:30
2013-01-18 17:00:52 +05:30
!--------------------------------------------------------------------------------------------------
! calculate volume of each grain
2009-05-26 22:54:42 +05:30
volumeOfGrain = 0.0_pReal
2012-01-12 16:06:17 +05:30
grain = 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 ) ! 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 ) ) / &
2012-02-10 17:26:05 +05:30
real ( dGrains , pReal )
2013-01-18 17:00:52 +05:30
grain = grain + dGrains ! wind forward by NgrainsPerIP
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-01-18 17:00:52 +05:30
mesh_ipVolume ( i , e ) / dGrains ! assign IPvolume/Ngrains to all grains of IP
grain = grain + FE_Nips ( t ) * dGrains ! wind forward by Nips*NgrainsPerIP
2009-03-04 17:18:54 +05:30
endif
enddo
2012-01-12 16:06:17 +05:30
2013-01-18 17:00:52 +05:30
!--------------------------------------------------------------------------------------------------
! divide myNgrains as best over constituents
2009-03-04 17:18:54 +05:30
NgrainsOfConstituent = 0_pInt
2012-02-16 00:28:38 +05:30
forall ( i = 1_pInt : microstructure_Nconstituents ( micro ) ) &
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-03-28 19:20:20 +05:30
do i = 1_pInt , microstructure_Nconstituents ( micro ) ! 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
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-03-28 19:20:20 +05:30
grain = 0_pInt ! reset microstructure grain index
2009-03-04 17:18:54 +05:30
2013-03-28 19:20:20 +05:30
constituents : do i = 1_pInt , microstructure_Nconstituents ( micro ) ! loop over constituents
2009-03-04 17:18:54 +05:30
phaseID = microstructure_phase ( i , micro )
textureID = microstructure_texture ( i , micro )
2013-03-28 19:20:20 +05:30
phaseOfGrain ( grain + 1_pInt : grain + NgrainsOfConstituent ( i ) ) = phaseID ! assign resp. phase
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)
constituentGrain = 0_pInt ! constituent grain index
2013-01-18 17:00:52 +05:30
!--------------------------------------------------------------------------------------------------
! dealing with texture components
if ( texture_ODFfile ( textureID ) == '' ) then
2013-03-28 19:20:20 +05:30
do t = 1_pInt , texture_Ngauss ( textureID ) ! loop over Gauss components
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
constituentGrain = constituentGrain + int ( myNorientations * texture_Gauss ( 5 , t , textureID ) )
enddo
2013-03-28 19:20:20 +05:30
do t = 1_pInt , texture_Nfiber ( textureID ) ! loop over fiber components
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
2012-02-16 00:28:38 +05:30
constituentGrain = constituentGrain + int ( myNorientations * texture_fiber ( 6 , t , textureID ) , pInt )
2009-03-04 17:18:54 +05:30
enddo
2013-01-18 17:00:52 +05:30
do j = constituentGrain + 1_pInt , myNorientations ! fill remainder with random
2009-03-04 17:18:54 +05:30
orientationOfGrain ( : , grain + j ) = math_sampleRandomOri ( )
enddo
2013-01-18 17:00:52 +05:30
else
2013-03-28 19:20:20 +05:30
2013-01-18 17:00:52 +05:30
!--------------------------------------------------------------------------------------------------
! hybrid IA
2009-03-04 17:18:54 +05:30
orientationOfGrain ( : , grain + 1 : grain + myNorientations ) = IO_hybridIA ( myNorientations , texture_ODFfile ( textureID ) )
2012-02-13 23:11:27 +05:30
if ( all ( orientationOfGrain ( : , grain + 1 ) == - 1.0_pReal ) ) call IO_error ( 156_pInt )
2009-03-04 17:18:54 +05:30
constituentGrain = constituentGrain + myNorientations
endif
2013-01-18 17:00:52 +05:30
2009-03-04 17:18:54 +05:30
symExtension = texture_symmetry ( textureID ) - 1_pInt
2013-01-18 17:00:52 +05:30
if ( symExtension > 0_pInt ) then ! sample symmetry
constituentGrain = NgrainsOfConstituent ( i ) - myNorientations ! calc remainder of array
do j = 1_pInt , myNorientations ! loop over each "real" orientation
2009-03-04 17:18:54 +05:30
symOrientation = math_symmetricEulers ( texture_symmetry ( textureID ) , orientationOfGrain ( : , j ) ) ! get symmetric equivalents
2013-01-18 17:00:52 +05:30
e = min ( symExtension , constituentGrain ) ! are we at end of constituent grain array?
2009-03-04 17:18:54 +05:30
if ( e > 0_pInt ) then
2012-02-16 00:28:38 +05:30
orientationOfGrain ( : , grain + myNorientations + 1 + ( j - 1_pInt ) * symExtension : &
grain + myNorientations + e + ( j - 1_pInt ) * symExtension ) = &
2009-03-04 17:18:54 +05:30
symOrientation ( : , 1 : e )
2013-01-18 17:00:52 +05:30
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
grain = grain + NgrainsOfConstituent ( i ) ! advance microstructure grain index
enddo constituents
2013-03-28 19:20:20 +05:30
!--------------------------------------------------------------------------------------------------
2013-01-18 17:00:52 +05:30
! unless element homogeneous, reshuffle grains
if ( . not . microstructure_elemhomo ( micro ) ) then
do i = 1_pInt , myNgrains - 1_pInt ! walk thru grains
2011-03-22 19:10:27 +05:30
call random_number ( rnd )
2013-01-18 17:00:52 +05:30
t = nint ( rnd * ( myNgrains - i ) + i + 0.5_pReal , pInt ) ! select a grain in remaining list
m = phaseOfGrain ( t ) ! exchange current with random
2011-03-22 19:10:27 +05:30
phaseOfGrain ( t ) = phaseOfGrain ( i )
phaseOfGrain ( i ) = m
2013-01-18 17:00:52 +05:30
m = textureOfGrain ( t ) ! exchange current with random
2011-03-22 19:10:27 +05:30
textureOfGrain ( t ) = textureOfGrain ( i )
textureOfGrain ( i ) = m
orientation = orientationOfGrain ( : , t )
orientationOfGrain ( : , t ) = orientationOfGrain ( : , i )
orientationOfGrain ( : , i ) = orientation
enddo
endif
2013-01-18 17:00:52 +05:30
!--------------------------------------------------------------------------------------------------
2013-03-28 19:20:20 +05:30
! calc fraction after weighing with volumePerGrain, exchange in MC steps to improve result
2012-01-12 16:06:17 +05:30
grain = 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
forall ( i = 1_pInt : FE_Nips ( t ) , g = 1_pInt : dGrains ) ! loop over IPs and grains
2012-01-12 16:06:17 +05:30
material_volume ( g , i , e ) = volumeOfGrain ( grain + g )
material_phase ( g , i , e ) = phaseOfGrain ( grain + g )
material_texture ( g , i , e ) = textureOfGrain ( grain + g )
material_EulerAngles ( : , g , i , e ) = orientationOfGrain ( : , grain + g )
end forall
2013-01-18 17:00:52 +05:30
FEsolving_execIP ( 2 , e ) = 1_pInt ! restrict calculation to first IP only, since all other results are to be copied from this
grain = grain + dGrains ! wind forward by NgrainsPerIP
2012-01-12 16:06:17 +05:30
else
2013-01-18 17:00:52 +05:30
forall ( i = 1_pInt : FE_Nips ( t ) , g = 1_pInt : dGrains ) ! loop over IPs and grains
2012-02-16 00:28:38 +05:30
material_volume ( g , i , e ) = volumeOfGrain ( grain + ( i - 1_pInt ) * dGrains + g )
material_phase ( g , i , e ) = phaseOfGrain ( grain + ( i - 1_pInt ) * dGrains + g )
material_texture ( g , i , e ) = textureOfGrain ( grain + ( i - 1_pInt ) * dGrains + g )
material_EulerAngles ( : , g , i , e ) = orientationOfGrain ( : , grain + ( i - 1_pInt ) * dGrains + g )
2012-01-12 16:06:17 +05:30
end forall
2013-01-18 17:00:52 +05:30
grain = grain + FE_Nips ( t ) * dGrains ! wind forward by Nips*NgrainsPerIP
2009-03-04 17:18:54 +05:30
endif
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-03-28 19:20:20 +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-03-28 19:20:20 +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
2013-02-20 03:42:05 +05:30
end module material