DAMASK with grid solvers  Revision: v2.0.3-2204-gdb1f2151
The Düsseldorf Advanced Material Simulation Kit with Grid Solvers
homogenization_mech_isostrain.f90
Go to the documentation of this file.
1 # 1 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/homogenization_mech_isostrain.f90"
2 # 1 "<built-in>"
3 # 1 "<command-line>"
4 # 1 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/homogenization_mech_isostrain.f90"
5 !--------------------------------------------------------------------------------------------------
10 !--------------------------------------------------------------------------------------------------
11 submodule(homogenization) homogenization_mech_isostrain
12 
13  enum, bind(c); enumerator :: &
14  parallel_id, &
15  average_id
16  end enum
17 
18  type :: tparameters
19  integer :: &
20  Nconstituents
21  integer(kind(average_ID)) :: &
22  mapping
23  end type
24 
25  type(tParameters), dimension(:), allocatable :: param
26 
27 
28 contains
29 
30 !--------------------------------------------------------------------------------------------------
32 !--------------------------------------------------------------------------------------------------
33 module subroutine mech_isostrain_init
34 
35  integer :: &
36  Ninstance, &
37  h, &
38  NofMyHomog
39  character(len=pStringLen) :: &
40  tag = ''
41 
42  write(6,'(/,a)') ' <<<+- homogenization_'//homogenization_isostrain_label//' init -+>>>'
43 
44  ninstance = count(homogenization_type == homogenization_isostrain_id)
45  if (iand(debug_level(debug_homogenization),debug_levelbasic) /= 0) &
46  write(6,'(a16,1x,i5,/)') '# instances:',ninstance
47 
48  allocate(param(ninstance)) ! one container of parameters per instance
49 
50  do h = 1, size(homogenization_type)
51  if (homogenization_type(h) /= homogenization_isostrain_id) cycle
52 
53  associate(prm => param(homogenization_typeinstance(h)),&
54  config => config_homogenization(h))
55 
56  prm%Nconstituents = config_homogenization(h)%getInt('nconstituents')
57  tag = 'sum'
58  select case(trim(config%getString('mapping',defaultval = tag)))
59  case ('sum')
60  prm%mapping = parallel_id
61  case ('avg')
62  prm%mapping = average_id
63  case default
64  call io_error(211,ext_msg=trim(tag)//' ('//homogenization_isostrain_label//')')
65  end select
66 
67  nofmyhomog = count(material_homogenizationat == h)
68  homogstate(h)%sizeState = 0
69  allocate(homogstate(h)%state0 (0,nofmyhomog))
70  allocate(homogstate(h)%subState0(0,nofmyhomog))
71  allocate(homogstate(h)%state (0,nofmyhomog))
72 
73  end associate
74 
75  enddo
76 
77 end subroutine mech_isostrain_init
78 
79 
80 !--------------------------------------------------------------------------------------------------
82 !--------------------------------------------------------------------------------------------------
83 module subroutine mech_isostrain_partitiondeformation(f,avgf)
84 
85  real(pReal), dimension (:,:,:), intent(out) :: F
86 
87  real(pReal), dimension (3,3), intent(in) :: avgF
88 
89  f = spread(avgf,3,size(f,3))
90 
91 end subroutine mech_isostrain_partitiondeformation
92 
93 
94 !--------------------------------------------------------------------------------------------------
96 !--------------------------------------------------------------------------------------------------
97 module subroutine mech_isostrain_averagestressanditstangent(avgp,davgpdavgf,p,dpdf,instance)
98 
99  real(pReal), dimension (3,3), intent(out) :: avgP
100  real(pReal), dimension (3,3,3,3), intent(out) :: dAvgPdAvgF
101 
102  real(pReal), dimension (:,:,:), intent(in) :: P
103  real(pReal), dimension (:,:,:,:,:), intent(in) :: dPdF
104  integer, intent(in) :: instance
105 
106  associate(prm => param(instance))
107 
108  select case (prm%mapping)
109  case (parallel_id)
110  avgp = sum(p,3)
111  davgpdavgf = sum(dpdf,5)
112  case (average_id)
113  avgp = sum(p,3) /real(prm%Nconstituents,preal)
114  davgpdavgf = sum(dpdf,5)/real(prm%Nconstituents,preal)
115  end select
116 
117  end associate
118 
119 end subroutine mech_isostrain_averagestressanditstangent
120 
121 end submodule homogenization_mech_isostrain
config
Reads in the material configuration from file.
Definition: config.f90:13
homogenization
homogenization manager, organizing deformation partitioning and stress homogenization
Definition: homogenization.f90:11