DAMASK with grid solvers  Revision: v2.0.3-2204-gdb1f2151
The Düsseldorf Advanced Material Simulation Kit with Grid Solvers
lattice.f90
Go to the documentation of this file.
1 # 1 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/lattice.f90"
2 # 1 "<built-in>"
3 # 1 "<command-line>"
4 # 1 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/lattice.f90"
5 !--------------------------------------------------------------------------------------------------
11 ! and cleavage as well as interaction among the various systems
12 !--------------------------------------------------------------------------------------------------
13 module lattice
14  use prec
15  use io
16  use config
17  use math
18  use rotations
19 
20  implicit none
21  private
22 
23 !--------------------------------------------------------------------------------------------------
24 ! face centered cubic
25  integer, dimension(2), parameter :: &
26  fcc_nslipsystem = [12, 6]
27 
28  integer, dimension(1), parameter :: &
29  fcc_ntwinsystem = [12]
30 
31  integer, dimension(1), parameter :: &
32  fcc_ntranssystem = [12]
33 
34  integer, dimension(1), parameter :: &
36 
37  integer, parameter :: &
38 
39  fcc_nslip = sum(fcc_nslipsystem), &
40  fcc_ntwin = sum(fcc_ntwinsystem), &
41  fcc_ntrans = sum(fcc_ntranssystem), &
43 
44 
45 
46 
47 
48 
49 
50  real(preal), dimension(3+3,FCC_NSLIP), parameter :: &
51  fcc_systemslip = reshape(real([&
52  ! Slip direction Plane normal ! SCHMID-BOAS notation
53  0, 1,-1, 1, 1, 1, & ! B2
54  -1, 0, 1, 1, 1, 1, & ! B4
55  1,-1, 0, 1, 1, 1, & ! B5
56  0,-1,-1, -1,-1, 1, & ! C1
57  1, 0, 1, -1,-1, 1, & ! C3
58  -1, 1, 0, -1,-1, 1, & ! C5
59  0,-1, 1, 1,-1,-1, & ! A2
60  -1, 0,-1, 1,-1,-1, & ! A3
61  1, 1, 0, 1,-1,-1, & ! A6
62  0, 1, 1, -1, 1,-1, & ! D1
63  1, 0,-1, -1, 1,-1, & ! D4
64  -1,-1, 0, -1, 1,-1, & ! D6
65  ! Slip system <110>{110}
66  1, 1, 0, 1,-1, 0, &
67  1,-1, 0, 1, 1, 0, &
68  1, 0, 1, 1, 0,-1, &
69  1, 0,-1, 1, 0, 1, &
70  0, 1, 1, 0, 1,-1, &
71  0, 1,-1, 0, 1, 1 &
72  ],preal),shape(fcc_systemslip))
73 
74  real(preal), dimension(3+3,FCC_NTWIN), parameter :: &
75  fcc_systemtwin = reshape(real( [&
76  -2, 1, 1, 1, 1, 1, &
77  1,-2, 1, 1, 1, 1, &
78  1, 1,-2, 1, 1, 1, &
79  2,-1, 1, -1,-1, 1, &
80  -1, 2, 1, -1,-1, 1, &
81  -1,-1,-2, -1,-1, 1, &
82  -2,-1,-1, 1,-1,-1, &
83  1, 2,-1, 1,-1,-1, &
84  1,-1, 2, 1,-1,-1, &
85  2, 1,-1, -1, 1,-1, &
86  -1,-2,-1, -1, 1,-1, &
87  -1, 1, 2, -1, 1,-1 &
88  ],preal),shape(fcc_systemtwin))
89 
90  integer, dimension(2,FCC_NTWIN), parameter, public :: &
92  2,3, &
93  1,3, &
94  1,2, &
95  5,6, &
96  4,6, &
97  4,5, &
98  8,9, &
99  7,9, &
100  7,8, &
101  11,12, &
102  10,12, &
103  10,11 &
105 
106  real(preal), dimension(3+3,FCC_NCLEAVAGE), parameter :: &
107  fcc_systemcleavage = reshape(real([&
108  ! Cleavage direction Plane normal
109  0, 1, 0, 1, 0, 0, &
110  0, 0, 1, 0, 1, 0, &
111  1, 0, 0, 0, 0, 1 &
112  ],preal),shape(fcc_systemcleavage))
113 
114 !--------------------------------------------------------------------------------------------------
115 ! body centered cubic
116  integer, dimension(2), parameter :: &
117  bcc_nslipsystem = [12, 12]
118 
119  integer, dimension(1), parameter :: &
120  bcc_ntwinsystem = [12]
121 
122  integer, dimension(1), parameter :: &
123  bcc_ncleavagesystem = [3]
124 
125  integer, parameter :: &
126 
127  bcc_nslip = sum(bcc_nslipsystem), &
128  bcc_ntwin = sum(bcc_ntwinsystem), &
130 
131 
132 
133 
134 
135 
136  real(preal), dimension(3+3,BCC_NSLIP), parameter :: &
137  bcc_systemslip = reshape(real([&
138  ! Slip direction Plane normal
139  ! Slip system <111>{110}
140  1,-1, 1, 0, 1, 1, &
141  -1,-1, 1, 0, 1, 1, &
142  1, 1, 1, 0,-1, 1, &
143  -1, 1, 1, 0,-1, 1, &
144  -1, 1, 1, 1, 0, 1, &
145  -1,-1, 1, 1, 0, 1, &
146  1, 1, 1, -1, 0, 1, &
147  1,-1, 1, -1, 0, 1, &
148  -1, 1, 1, 1, 1, 0, &
149  -1, 1,-1, 1, 1, 0, &
150  1, 1, 1, -1, 1, 0, &
151  1, 1,-1, -1, 1, 0, &
152  ! Slip system <111>{112}
153  -1, 1, 1, 2, 1, 1, &
154  1, 1, 1, -2, 1, 1, &
155  1, 1,-1, 2,-1, 1, &
156  1,-1, 1, 2, 1,-1, &
157  1,-1, 1, 1, 2, 1, &
158  1, 1,-1, -1, 2, 1, &
159  1, 1, 1, 1,-2, 1, &
160  -1, 1, 1, 1, 2,-1, &
161  1, 1,-1, 1, 1, 2, &
162  1,-1, 1, -1, 1, 2, &
163  -1, 1, 1, 1,-1, 2, &
164  1, 1, 1, 1, 1,-2 &
165  ],preal),shape(bcc_systemslip))
166 
167  real(preal), dimension(3+3,BCC_NTWIN), parameter :: &
168  bcc_systemtwin = reshape(real([&
169  ! Twin system <111>{112}
170  -1, 1, 1, 2, 1, 1, &
171  1, 1, 1, -2, 1, 1, &
172  1, 1,-1, 2,-1, 1, &
173  1,-1, 1, 2, 1,-1, &
174  1,-1, 1, 1, 2, 1, &
175  1, 1,-1, -1, 2, 1, &
176  1, 1, 1, 1,-2, 1, &
177  -1, 1, 1, 1, 2,-1, &
178  1, 1,-1, 1, 1, 2, &
179  1,-1, 1, -1, 1, 2, &
180  -1, 1, 1, 1,-1, 2, &
181  1, 1, 1, 1, 1,-2 &
182  ],preal),shape(bcc_systemtwin))
183 
184  real(preal), dimension(3+3,BCC_NCLEAVAGE), parameter :: &
185  bcc_systemcleavage = reshape(real([&
186  ! Cleavage direction Plane normal
187  0, 1, 0, 1, 0, 0, &
188  0, 0, 1, 0, 1, 0, &
189  1, 0, 0, 0, 0, 1 &
190  ],preal),shape(bcc_systemcleavage))
191 
192 !--------------------------------------------------------------------------------------------------
193 ! hexagonal
194  integer, dimension(6), parameter :: &
195  hex_nslipsystem = [3, 3, 3, 6, 12, 6]
196 
197  integer, dimension(4), parameter :: &
198  hex_ntwinsystem = [6, 6, 6, 6]
199 
200  integer, parameter :: &
201 
202  hex_nslip = sum(hex_nslipsystem), &
203  hex_ntwin = sum(hex_ntwinsystem)
204 
205 
206 
207 
208 
209  real(preal), dimension(4+4,HEX_NSLIP), parameter :: &
210  hex_systemslip = reshape(real([&
211  ! Slip direction Plane normal
212  ! Basal systems <-1-1.0>{00.1} (independent of c/a-ratio, Bravais notation (4 coordinate base))
213  2, -1, -1, 0, 0, 0, 0, 1, &
214  -1, 2, -1, 0, 0, 0, 0, 1, &
215  -1, -1, 2, 0, 0, 0, 0, 1, &
216  ! 1st type prismatic systems <-1-1.0>{1-1.0} (independent of c/a-ratio)
217  2, -1, -1, 0, 0, 1, -1, 0, &
218  -1, 2, -1, 0, -1, 0, 1, 0, &
219  -1, -1, 2, 0, 1, -1, 0, 0, &
220  ! 2nd type prismatic systems <-11.0>{11.0} -- a slip; plane normals independent of c/a-ratio
221  -1, 1, 0, 0, 1, 1, -2, 0, &
222  0, -1, 1, 0, -2, 1, 1, 0, &
223  1, 0, -1, 0, 1, -2, 1, 0, &
224  ! 1st type 1st order pyramidal systems <-1-1.0>{-11.1} -- plane normals depend on the c/a-ratio
225  -1, 2, -1, 0, 1, 0, -1, 1, &
226  -2, 1, 1, 0, 0, 1, -1, 1, &
227  -1, -1, 2, 0, -1, 1, 0, 1, &
228  1, -2, 1, 0, -1, 0, 1, 1, &
229  2, -1, -1, 0, 0, -1, 1, 1, &
230  1, 1, -2, 0, 1, -1, 0, 1, &
231  ! pyramidal system: c+a slip <11.3>{-10.1} -- plane normals depend on the c/a-ratio
232  -2, 1, 1, 3, 1, 0, -1, 1, &
233  -1, -1, 2, 3, 1, 0, -1, 1, &
234  -1, -1, 2, 3, 0, 1, -1, 1, &
235  1, -2, 1, 3, 0, 1, -1, 1, &
236  1, -2, 1, 3, -1, 1, 0, 1, &
237  2, -1, -1, 3, -1, 1, 0, 1, &
238  2, -1, -1, 3, -1, 0, 1, 1, &
239  1, 1, -2, 3, -1, 0, 1, 1, &
240  1, 1, -2, 3, 0, -1, 1, 1, &
241  -1, 2, -1, 3, 0, -1, 1, 1, &
242  -1, 2, -1, 3, 1, -1, 0, 1, &
243  -2, 1, 1, 3, 1, -1, 0, 1, &
244  ! pyramidal system: c+a slip <11.3>{-1-1.2} -- as for hexagonal ice (Castelnau et al. 1996, similar to twin system found below)
245  -1, -1, 2, 3, 1, 1, -2, 2, & ! <11.3>{-1-1.2} shear = 2((c/a)^2-2)/(3 c/a)
246  1, -2, 1, 3, -1, 2, -1, 2, &
247  2, -1, -1, 3, -2, 1, 1, 2, &
248  1, 1, -2, 3, -1, -1, 2, 2, &
249  -1, 2, -1, 3, 1, -2, 1, 2, &
250  -2, 1, 1, 3, 2, -1, -1, 2 &
251  ],preal),shape(hex_systemslip))
252 
253  real(preal), dimension(4+4,HEX_NTWIN), parameter :: &
254  hex_systemtwin = reshape(real([&
255  ! Compression or Tension = f(twinning shear=f(c/a)) for each metal ! (according to Yoo 1981)
256  -1, 0, 1, 1, 1, 0, -1, 2, & ! <-10.1>{10.2} shear = (3-(c/a)^2)/(sqrt(3) c/a)
257  0, -1, 1, 1, 0, 1, -1, 2, &
258  1, -1, 0, 1, -1, 1, 0, 2, &
259  1, 0, -1, 1, -1, 0, 1, 2, &
260  0, 1, -1, 1, 0, -1, 1, 2, &
261  -1, 1, 0, 1, 1, -1, 0, 2, &
262 !
263  -1, -1, 2, 6, 1, 1, -2, 1, & ! <11.6>{-1-1.1} shear = 1/(c/a)
264  1, -2, 1, 6, -1, 2, -1, 1, &
265  2, -1, -1, 6, -2, 1, 1, 1, &
266  1, 1, -2, 6, -1, -1, 2, 1, &
267  -1, 2, -1, 6, 1, -2, 1, 1, &
268  -2, 1, 1, 6, 2, -1, -1, 1, &
269 !
270  1, 0, -1, -2, 1, 0, -1, 1, & ! <10.-2>{10.1} shear = (4(c/a)^2-9)/(4 sqrt(3) c/a)
271  0, 1, -1, -2, 0, 1, -1, 1, &
272  -1, 1, 0, -2, -1, 1, 0, 1, &
273  -1, 0, 1, -2, -1, 0, 1, 1, &
274  0, -1, 1, -2, 0, -1, 1, 1, &
275  1, -1, 0, -2, 1, -1, 0, 1, &
276 !
277  1, 1, -2, -3, 1, 1, -2, 2, & ! <11.-3>{11.2} shear = 2((c/a)^2-2)/(3 c/a)
278  -1, 2, -1, -3, -1, 2, -1, 2, &
279  -2, 1, 1, -3, -2, 1, 1, 2, &
280  -1, -1, 2, -3, -1, -1, 2, 2, &
281  1, -2, 1, -3, 1, -2, 1, 2, &
282  2, -1, -1, -3, 2, -1, -1, 2 &
283  ],preal),shape(hex_systemtwin))
284 
285 !--------------------------------------------------------------------------------------------------
286 ! body centered tetragonal
287  integer, dimension(13), parameter :: &
288  bct_nslipsystem = [2, 2, 2, 4, 2, 4, 2, 2, 4, 8, 4, 8, 8 ]
289 
290  integer, parameter :: &
291 
292  bct_nslip = sum(bct_nslipsystem)
293 
294 
295 
296 
297  real(preal), dimension(3+3,BCT_NSLIP), parameter :: &
298  bct_systemslip = reshape(real([&
299  ! Slip direction Plane normal
300  ! Slip family 1 {100)<001] (Bravais notation {hkl)<uvw] for bct c/a = 0.5456)
301  0, 0, 1, 1, 0, 0, &
302  0, 0, 1, 0, 1, 0, &
303  ! Slip family 2 {110)<001]
304  0, 0, 1, 1, 1, 0, &
305  0, 0, 1, -1, 1, 0, &
306  ! slip family 3 {100)<010]
307  0, 1, 0, 1, 0, 0, &
308  1, 0, 0, 0, 1, 0, &
309  ! Slip family 4 {110)<1-11]/2
310  1,-1, 1, 1, 1, 0, &
311  1,-1,-1, 1, 1, 0, &
312  -1,-1,-1, -1, 1, 0, &
313  -1,-1, 1, -1, 1, 0, &
314  ! Slip family 5 {110)<1-10]
315  1, -1, 0, 1, 1, 0, &
316  1, 1, 0, 1,-1, 0, &
317  ! Slip family 6 {100)<011]
318  0, 1, 1, 1, 0, 0, &
319  0,-1, 1, 1, 0, 0, &
320  -1, 0, 1, 0, 1, 0, &
321  1, 0, 1, 0, 1, 0, &
322  ! Slip family 7 {001)<010]
323  0, 1, 0, 0, 0, 1, &
324  1, 0, 0, 0, 0, 1, &
325  ! Slip family 8 {001)<110]
326  1, 1, 0, 0, 0, 1, &
327  -1, 1, 0, 0, 0, 1, &
328  ! Slip family 9 {011)<01-1]
329  0, 1,-1, 0, 1, 1, &
330  0,-1,-1, 0,-1, 1, &
331  -1, 0,-1, -1, 0, 1, &
332  1, 0,-1, 1, 0, 1, &
333  ! Slip family 10 {011)<1-11]/2
334  1,-1, 1, 0, 1, 1, &
335  1, 1,-1, 0, 1, 1, &
336  1, 1, 1, 0, 1,-1, &
337  -1, 1, 1, 0, 1,-1, &
338  1,-1,-1, 1, 0, 1, &
339  -1,-1, 1, 1, 0, 1, &
340  1, 1, 1, 1, 0,-1, &
341  1,-1, 1, 1, 0,-1, &
342  ! Slip family 11 {011)<100]
343  1, 0, 0, 0, 1, 1, &
344  1, 0, 0, 0, 1,-1, &
345  0, 1, 0, 1, 0, 1, &
346  0, 1, 0, 1, 0,-1, &
347  ! Slip family 12 {211)<01-1]
348  0, 1,-1, 2, 1, 1, &
349  0,-1,-1, 2,-1, 1, &
350  1, 0,-1, 1, 2, 1, &
351  -1, 0,-1, -1, 2, 1, &
352  0, 1,-1, -2, 1, 1, &
353  0,-1,-1, -2,-1, 1, &
354  -1, 0,-1, -1,-2, 1, &
355  1, 0,-1, 1,-2, 1, &
356  ! Slip family 13 {211)<-111]/2
357  -1, 1, 1, 2, 1, 1, &
358  -1,-1, 1, 2,-1, 1, &
359  1,-1, 1, 1, 2, 1, &
360  -1,-1, 1, -1, 2, 1, &
361  1, 1, 1, -2, 1, 1, &
362  1,-1, 1, -2,-1, 1, &
363  -1, 1, 1, -1,-2, 1, &
364  1, 1, 1, 1,-2, 1 &
365  ],preal),shape(bct_systemslip))
366 
367 !--------------------------------------------------------------------------------------------------
368 ! orthorhombic
369  integer, dimension(3), parameter :: &
370  ort_ncleavagesystem = [1, 1, 1]
371 
372  integer, parameter :: &
373 
375 
376 
377 
378 
379  real(preal), dimension(3+3,ORT_NCLEAVAGE), parameter :: &
380  ort_systemcleavage = reshape(real([&
381  ! Cleavage direction Plane normal
382  0, 1, 0, 1, 0, 0, &
383  0, 0, 1, 0, 1, 0, &
384  1, 0, 0, 0, 0, 1 &
385  ],preal),shape(ort_systemcleavage))
386 
387 
388  enum, bind(c); enumerator :: &
390  lattice_iso_id, &
391  lattice_fcc_id, &
392  lattice_bcc_id, &
393  lattice_bct_id, &
394  lattice_hex_id, &
396  end enum
397 
398 ! SHOULD NOT BE PART OF LATTICE BEGIN
399  real(preal), dimension(:), allocatable, public, protected :: &
404  real(preal), dimension(:,:,:), allocatable, public, protected :: &
405  lattice_c66, &
408  integer(kind(lattice_UNDEFINED_ID)), dimension(:), allocatable, public, protected :: &
410 ! SHOULD NOT BE PART OF LATTICE END
411 
413  module procedure slipprojection_transverse
414  end interface lattice_forestprojection_edge
415 
417  module procedure slipprojection_direction
418  end interface lattice_forestprojection_screw
419 
420  public :: &
421  lattice_init, &
422  lattice_iso_id, &
423  lattice_fcc_id, &
424  lattice_bcc_id, &
425  lattice_bct_id, &
426  lattice_hex_id, &
427  lattice_ort_id, &
450 
451 contains
452 
453 !--------------------------------------------------------------------------------------------------
455 !--------------------------------------------------------------------------------------------------
456 subroutine lattice_init
457 
458  integer :: nphases, p,i
459  character(len=pStringLen) :: structure = ''
460 
461  write(6,'(/,a)') ' <<<+- lattice init -+>>>'; flush(6)
462 
463  nphases = size(config_phase)
464 
465  allocate(lattice_structure(nphases),source = lattice_undefined_id)
466  allocate(lattice_c66(6,6,nphases), source=0.0_preal)
467 
468  allocate(lattice_thermalconductivity(3,3,nphases), source=0.0_preal)
469  allocate(lattice_damagediffusion(3,3,nphases), source=0.0_preal)
470 
471  allocate(lattice_damagemobility,&
474  source=[(0.0_preal,i=1,nphases)])
475 
476  do p = 1, size(config_phase)
477 
478  lattice_c66(1,1,p) = config_phase(p)%getFloat('c11')
479  lattice_c66(1,2,p) = config_phase(p)%getFloat('c12')
480 
481  lattice_c66(1,3,p) = config_phase(p)%getFloat('c13',defaultval=0.0_preal)
482  lattice_c66(2,2,p) = config_phase(p)%getFloat('c22',defaultval=0.0_preal)
483  lattice_c66(2,3,p) = config_phase(p)%getFloat('c23',defaultval=0.0_preal)
484  lattice_c66(3,3,p) = config_phase(p)%getFloat('c33',defaultval=0.0_preal)
485  lattice_c66(4,4,p) = config_phase(p)%getFloat('c44',defaultval=0.0_preal)
486  lattice_c66(5,5,p) = config_phase(p)%getFloat('c55',defaultval=0.0_preal)
487  lattice_c66(6,6,p) = config_phase(p)%getFloat('c66',defaultval=0.0_preal)
488 
489  structure = config_phase(p)%getString('lattice_structure')
490  select case(trim(structure))
491  case('iso')
493  case('fcc')
495  case('bcc')
497  case('hex')
499  case('bct')
501  case('ort')
503  case default
504  call io_error(130,ext_msg='lattice_init: '//trim(structure))
505  end select
506 
507  lattice_c66(1:6,1:6,p) = applylatticesymmetryc66(lattice_c66(1:6,1:6,p),structure)
508 
509  lattice_mu(p) = equivalent_mu(lattice_c66(1:6,1:6,p),'voigt')
510  lattice_nu(p) = equivalent_nu(lattice_c66(1:6,1:6,p),'voigt')
511 
512  lattice_c66(1:6,1:6,p) = math_sym3333to66(math_voigt66to3333(lattice_c66(1:6,1:6,p))) ! Literature data is in Voigt notation
513  do i = 1, 6
514  if (abs(lattice_c66(i,i,p))<tol_math_check) &
515  call io_error(135,el=i,ip=p,ext_msg='matrix diagonal "el"ement of phase "ip"')
516  enddo
517 
518 
519  ! SHOULD NOT BE PART OF LATTICE BEGIN
520  lattice_thermalconductivity(1,1,p) = config_phase(p)%getFloat('thermal_conductivity11',defaultval=0.0_preal)
521  lattice_thermalconductivity(2,2,p) = config_phase(p)%getFloat('thermal_conductivity22',defaultval=0.0_preal)
522  lattice_thermalconductivity(3,3,p) = config_phase(p)%getFloat('thermal_conductivity33',defaultval=0.0_preal)
524 
525  lattice_specificheat(p) = config_phase(p)%getFloat('specific_heat',defaultval=0.0_preal)
526  lattice_massdensity(p) = config_phase(p)%getFloat('mass_density', defaultval=0.0_preal)
527 
528  lattice_damagediffusion(1,1,p) = config_phase(p)%getFloat('damage_diffusion11',defaultval=0.0_preal)
529  lattice_damagediffusion(2,2,p) = config_phase(p)%getFloat('damage_diffusion22',defaultval=0.0_preal)
530  lattice_damagediffusion(3,3,p) = config_phase(p)%getFloat('damage_diffusion33',defaultval=0.0_preal)
532 
533  lattice_damagemobility(p) = config_phase(p)%getFloat('damage_mobility',defaultval=0.0_preal)
534  ! SHOULD NOT BE PART OF LATTICE END
535 
536  call unittest
537 
538  enddo
539 
540 end subroutine lattice_init
541 
542 
543 !--------------------------------------------------------------------------------------------------
545 !--------------------------------------------------------------------------------------------------
546 function lattice_characteristicshear_twin(Ntwin,structure,CoverA) result(characteristicShear)
547 
548  integer, dimension(:), intent(in) :: ntwin
549  character(len=*), intent(in) :: structure
550  real(preal), intent(in) :: covera
551  real(preal), dimension(sum(Ntwin)) :: characteristicshear
552 
553  integer :: &
554  a, & !< index of active system
555  p, & !< index in potential system list
556  f, & !< index of my family
557  s
558 
559  integer, dimension(HEX_NTWIN), parameter :: &
560  hex_sheartwin = reshape( [&
561  1, & ! <-10.1>{10.2}
562  1, &
563  1, &
564  1, &
565  1, &
566  1, &
567  2, & ! <11.6>{-1-1.1}
568  2, &
569  2, &
570  2, &
571  2, &
572  2, &
573  3, & ! <10.-2>{10.1}
574  3, &
575  3, &
576  3, &
577  3, &
578  3, &
579  4, & ! <11.-3>{11.2}
580  4, &
581  4, &
582  4, &
583  4, &
584  4 &
585  ],[hex_ntwin]) ! indicator to formulas below
586 
587  if (len_trim(structure) /= 3) &
588  call io_error(137,ext_msg='lattice_characteristicShear_Twin: '//trim(structure))
589 
590  a = 0
591  myfamilies: do f = 1,size(ntwin,1)
592  mysystems: do s = 1,ntwin(f)
593  a = a + 1
594  select case(structure)
595  case('fcc','bcc')
596  characteristicshear(a) = 0.5_preal*sqrt(2.0_preal)
597  case('hex')
598  if (covera < 1.0_preal .or. covera > 2.0_preal) &
599  call io_error(131,ext_msg='lattice_characteristicShear_Twin')
600  p = sum(hex_ntwinsystem(1:f-1))+s
601  select case(hex_sheartwin(p)) ! from Christian & Mahajan 1995 p.29
602  case (1) ! <-10.1>{10.2}
603  characteristicshear(a) = (3.0_preal-covera**2.0_preal)/sqrt(3.0_preal)/covera
604  case (2) ! <11.6>{-1-1.1}
605  characteristicshear(a) = 1.0_preal/covera
606  case (3) ! <10.-2>{10.1}
607  characteristicshear(a) = (4.0_preal*covera**2.0_preal-9.0_preal)/sqrt(48.0_preal)/covera
608  case (4) ! <11.-3>{11.2}
609  characteristicshear(a) = 2.0_preal*(covera**2.0_preal-2.0_preal)/3.0_preal/covera
610  end select
611  case default
612  call io_error(137,ext_msg='lattice_characteristicShear_Twin: '//trim(structure))
613  end select
614  enddo mysystems
615  enddo myfamilies
616 
618 
619 
620 !--------------------------------------------------------------------------------------------------
622 !--------------------------------------------------------------------------------------------------
623 function lattice_c66_twin(Ntwin,C66,structure,CoverA)
624 
625  integer, dimension(:), intent(in) :: ntwin
626  character(len=*), intent(in) :: structure
627  real(preal), dimension(6,6), intent(in) :: c66
628  real(preal), intent(in) :: covera
629  real(preal), dimension(6,6,sum(Ntwin)) :: lattice_c66_twin
630 
631  real(preal), dimension(3,3,sum(Ntwin)):: coordinatesystem
632  type(rotation) :: r
633  integer :: i
634 
635  if (len_trim(structure) /= 3) &
636  call io_error(137,ext_msg='lattice_C66_twin: '//trim(structure))
637 
638  select case(structure)
639  case('fcc')
640  coordinatesystem = buildcoordinatesystem(ntwin,fcc_nslipsystem,fcc_systemtwin,&
641  trim(structure),0.0_preal)
642  case('bcc')
643  coordinatesystem = buildcoordinatesystem(ntwin,bcc_nslipsystem,bcc_systemtwin,&
644  trim(structure),0.0_preal)
645  case('hex')
646  coordinatesystem = buildcoordinatesystem(ntwin,hex_nslipsystem,hex_systemtwin,&
647  'hex',covera)
648  case default
649  call io_error(137,ext_msg='lattice_C66_twin: '//trim(structure))
650  end select
651 
652  do i = 1, sum(ntwin)
653  call r%fromAxisAngle([coordinatesystem(1:3,2,i),pi],p=1) ! ToDo: Why always 180 deg?
654  lattice_c66_twin(1:6,1:6,i) = r%rotTensor4sym(c66)
655  enddo
656 
657 end function lattice_c66_twin
658 
659 
660 !--------------------------------------------------------------------------------------------------
662 !--------------------------------------------------------------------------------------------------
663 function lattice_c66_trans(Ntrans,C_parent66,structure_target, &
664  cOverA_trans,a_bcc,a_fcc)
665 
666  integer, dimension(:), intent(in) :: ntrans
667  character(len=*), intent(in) :: structure_target
668  real(preal), dimension(6,6), intent(in) :: c_parent66
669  real(preal), dimension(6,6,sum(Ntrans)) :: lattice_c66_trans
670 
671  real(preal), dimension(6,6) :: c_bar66, c_target_unrotated66
672  real(preal), dimension(3,3,sum(Ntrans)) :: q,s
673  type(rotation) :: r
674  real(preal) :: a_bcc, a_fcc, covera_trans
675  integer :: i
676 
677  if (len_trim(structure_target) /= 3) &
678  call io_error(137,ext_msg='lattice_C66_trans (target): '//trim(structure_target))
679 
680  !--------------------------------------------------------------------------------------------------
681  ! elasticity matrix of the target phase in cube orientation
682  if (structure_target(1:3) == 'hex') then
683  if (covera_trans < 1.0_preal .or. covera_trans > 2.0_preal) &
684  call io_error(131,ext_msg='lattice_C66_trans: '//trim(structure_target))
685  c_bar66(1,1) = (c_parent66(1,1) + c_parent66(1,2) + 2.0_preal*c_parent66(4,4))/2.0_preal
686  c_bar66(1,2) = (c_parent66(1,1) + 5.0_preal*c_parent66(1,2) - 2.0_preal*c_parent66(4,4))/6.0_preal
687  c_bar66(3,3) = (c_parent66(1,1) + 2.0_preal*c_parent66(1,2) + 4.0_preal*c_parent66(4,4))/3.0_preal
688  c_bar66(1,3) = (c_parent66(1,1) + 2.0_preal*c_parent66(1,2) - 2.0_preal*c_parent66(4,4))/3.0_preal
689  c_bar66(4,4) = (c_parent66(1,1) - c_parent66(1,2) + c_parent66(4,4))/3.0_preal
690  c_bar66(1,4) = (c_parent66(1,1) - c_parent66(1,2) - 2.0_preal*c_parent66(4,4)) /(3.0_preal*sqrt(2.0_preal))
691 
692  c_target_unrotated66 = 0.0_preal
693  c_target_unrotated66(1,1) = c_bar66(1,1) - c_bar66(1,4)**2.0_preal/c_bar66(4,4)
694  c_target_unrotated66(1,2) = c_bar66(1,2) + c_bar66(1,4)**2.0_preal/c_bar66(4,4)
695  c_target_unrotated66(1,3) = c_bar66(1,3)
696  c_target_unrotated66(3,3) = c_bar66(3,3)
697  c_target_unrotated66(4,4) = c_bar66(4,4) - c_bar66(1,4)**2.0_preal/(0.5_preal*(c_bar66(1,1) - c_bar66(1,2)))
698  c_target_unrotated66 = applylatticesymmetryc66(c_target_unrotated66,'hex')
699  elseif (structure_target(1:3) == 'bcc') then
700  if (a_bcc <= 0.0_preal .or. a_fcc <= 0.0_preal) &
701  call io_error(134,ext_msg='lattice_C66_trans: '//trim(structure_target))
702  c_target_unrotated66 = c_parent66
703  else
704  call io_error(137,ext_msg='lattice_C66_trans : '//trim(structure_target))
705  endif
706 
707  do i = 1, 6
708  if (abs(c_target_unrotated66(i,i))<tol_math_check) &
709  call io_error(135,el=i,ext_msg='matrix diagonal "el"ement in transformation')
710  enddo
711 
712  call buildtransformationsystem(q,s,ntrans,covera_trans,a_fcc,a_bcc)
713 
714  do i = 1, sum(ntrans)
715  call r%fromMatrix(q(1:3,1:3,i))
716  lattice_c66_trans(1:6,1:6,i) = r%rotTensor4sym(c_target_unrotated66)
717  enddo
718 
719  end function lattice_c66_trans
720 
721 
722 !--------------------------------------------------------------------------------------------------
724 ! Koester et al. 2012, Acta Materialia 60 (2012) 3894–3901, eq. (17)
725 ! Gröger et al. 2008, Acta Materialia 56 (2008) 5412–5425, table 1
726 !--------------------------------------------------------------------------------------------------
727 function lattice_nonschmidmatrix(Nslip,nonSchmidCoefficients,sense) result(nonSchmidMatrix)
728 
729  integer, dimension(:), intent(in) :: nslip
730  real(preal), dimension(:), intent(in) :: nonschmidcoefficients
731  integer, intent(in) :: sense
732  real(preal), dimension(1:3,1:3,sum(Nslip)) :: nonschmidmatrix
733 
734  real(preal), dimension(1:3,1:3,sum(Nslip)) :: coordinatesystem
735  real(preal), dimension(3) :: direction, normal, np
736  type(rotation) :: r
737  integer :: i
738 
739  if (abs(sense) /= 1) call io_error(0,ext_msg='lattice_nonSchmidMatrix')
740 
741  coordinatesystem = buildcoordinatesystem(nslip,bcc_nslipsystem,bcc_systemslip,&
742  'bcc',0.0_preal)
743  coordinatesystem(1:3,1,1:sum(nslip)) = coordinatesystem(1:3,1,1:sum(nslip))*real(sense,preal) ! convert unidirectional coordinate system
744  nonschmidmatrix = lattice_schmidmatrix_slip(nslip,'bcc',0.0_preal) ! Schmid contribution
745 
746  do i = 1,sum(nslip)
747  direction = coordinatesystem(1:3,1,i)
748  normal = coordinatesystem(1:3,2,i)
749  call r%fromAxisAngle([direction,60.0_preal],degrees=.true.,p=1)
750  np = r%rotate(normal)
751 
752  if (size(nonschmidcoefficients)>0) nonschmidmatrix(1:3,1:3,i) = nonschmidmatrix(1:3,1:3,i) &
753  + nonschmidcoefficients(1) * math_outer(direction, np)
754  if (size(nonschmidcoefficients)>1) nonschmidmatrix(1:3,1:3,i) = nonschmidmatrix(1:3,1:3,i) &
755  + nonschmidcoefficients(2) * math_outer(math_cross(normal, direction), normal)
756  if (size(nonschmidcoefficients)>2) nonschmidmatrix(1:3,1:3,i) = nonschmidmatrix(1:3,1:3,i) &
757  + nonschmidcoefficients(3) * math_outer(math_cross(np, direction), np)
758  if (size(nonschmidcoefficients)>3) nonschmidmatrix(1:3,1:3,i) = nonschmidmatrix(1:3,1:3,i) &
759  + nonschmidcoefficients(4) * math_outer(normal, normal)
760  if (size(nonschmidcoefficients)>4) nonschmidmatrix(1:3,1:3,i) = nonschmidmatrix(1:3,1:3,i) &
761  + nonschmidcoefficients(5) * math_outer(math_cross(normal, direction), &
762  math_cross(normal, direction))
763  if (size(nonschmidcoefficients)>5) nonschmidmatrix(1:3,1:3,i) = nonschmidmatrix(1:3,1:3,i) &
764  + nonschmidcoefficients(6) * math_outer(direction, direction)
765  enddo
766 
767 end function lattice_nonschmidmatrix
768 
769 
770 !--------------------------------------------------------------------------------------------------
773 !--------------------------------------------------------------------------------------------------
774 function lattice_interaction_slipbyslip(Nslip,interactionValues,structure) result(interactionMatrix)
775 
776  integer, dimension(:), intent(in) :: nslip
777  real(preal), dimension(:), intent(in) :: interactionvalues
778  character(len=*), intent(in) :: structure
779  real(preal), dimension(sum(Nslip),sum(Nslip)) :: interactionmatrix
780 
781  integer, dimension(:), allocatable :: nslipmax
782  integer, dimension(:,:), allocatable :: interactiontypes
783 
784  integer, dimension(FCC_NSLIP,FCC_NSLIP), parameter :: &
785  fcc_interactionslipslip = reshape( [&
786  1, 2, 2, 4, 6, 5, 3, 5, 5, 4, 5, 6, 9,10, 9,10,11,12, & ! -----> acting
787  2, 1, 2, 6, 4, 5, 5, 4, 6, 5, 3, 5, 9,10,11,12, 9,10, & ! |
788  2, 2, 1, 5, 5, 3, 5, 6, 4, 6, 5, 4, 11,12, 9,10, 9,10, & ! |
789  4, 6, 5, 1, 2, 2, 4, 5, 6, 3, 5, 5, 9,10,10, 9,12,11, & ! v
790  6, 4, 5, 2, 1, 2, 5, 3, 5, 5, 4, 6, 9,10,12,11,10, 9, & ! reacting
791  5, 5, 3, 2, 2, 1, 6, 5, 4, 5, 6, 4, 11,12,10, 9,10, 9, &
792  3, 5, 5, 4, 5, 6, 1, 2, 2, 4, 6, 5, 10, 9,10, 9,11,12, &
793  5, 4, 6, 5, 3, 5, 2, 1, 2, 6, 4, 5, 10, 9,12,11, 9,10, &
794  5, 6, 4, 6, 5, 4, 2, 2, 1, 5, 5, 3, 12,11,10, 9, 9,10, &
795  4, 5, 6, 3, 5, 5, 4, 6, 5, 1, 2, 2, 10, 9, 9,10,12,11, &
796  5, 3, 5, 5, 4, 6, 6, 4, 5, 2, 1, 2, 10, 9,11,12,10, 9, &
797  6, 5, 4, 5, 6, 4, 5, 5, 3, 2, 2, 1, 12,11, 9,10,10, 9, &
798 
799  9, 9,11, 9, 9,11,10,10,12,10,10,12, 1, 7, 8, 8, 8, 8, &
800  10,10,12,10,10,12, 9, 9,11, 9, 9,11, 7, 1, 8, 8, 8, 8, &
801  9,11, 9,10,12,10,10,12,10, 9,11, 9, 8, 8, 1, 7, 8, 8, &
802  10,12,10, 9,11, 9, 9,11, 9,10,12,10, 8, 8, 7, 1, 8, 8, &
803  11, 9, 9,12,10,10,11, 9, 9,12,10,10, 8, 8, 8, 8, 1, 7, &
804  12,10,10,11, 9, 9,12,10,10,11, 9, 9, 8, 8, 8, 8, 7, 1 &
805  ],shape(fcc_interactionslipslip))
818 
819  integer, dimension(BCC_NSLIP,BCC_NSLIP), parameter :: &
820  bcc_interactionslipslip = reshape( [&
821  1,2,6,6,5,4,4,3,4,3,5,4, 6,6,4,3,3,4,6,6,4,3,6,6, & ! -----> acting
822  2,1,6,6,4,3,5,4,5,4,4,3, 6,6,3,4,4,3,6,6,3,4,6,6, & ! |
823  6,6,1,2,4,5,3,4,4,5,3,4, 4,3,6,6,6,6,3,4,6,6,4,3, & ! |
824  6,6,2,1,3,4,4,5,3,4,4,5, 3,4,6,6,6,6,4,3,6,6,3,4, & ! v
825  5,4,4,3,1,2,6,6,3,4,5,4, 3,6,4,6,6,4,6,3,4,6,3,6, & ! reacting
826  4,3,5,4,2,1,6,6,4,5,4,3, 4,6,3,6,6,3,6,4,3,6,4,6, &
827  4,5,3,4,6,6,1,2,5,4,3,4, 6,3,6,4,4,6,3,6,6,4,6,3, &
828  3,4,4,5,6,6,2,1,4,3,4,5, 6,4,6,3,3,6,4,6,6,3,6,4, &
829  4,5,4,3,3,4,5,4,1,2,6,6, 3,6,6,4,4,6,6,3,6,4,3,6, &
830  3,4,5,4,4,5,4,3,2,1,6,6, 4,6,6,3,3,6,6,4,6,3,4,6, &
831  5,4,3,4,5,4,3,4,6,6,1,2, 6,3,4,6,6,4,3,6,4,6,6,3, &
832  4,3,4,5,4,3,4,5,6,6,2,1, 6,4,3,6,6,3,4,6,3,6,6,4, &
833  !
834  6,6,4,3,3,4,6,6,3,4,6,6, 1,5,6,6,5,6,6,3,5,6,3,6, &
835  6,6,3,4,6,6,3,4,6,6,3,4, 5,1,6,6,6,5,3,6,6,5,6,3, &
836  4,3,6,6,4,3,6,6,6,6,4,3, 6,6,1,5,6,3,5,6,3,6,5,6, &
837  3,4,6,6,6,6,4,3,4,3,6,6, 6,6,5,1,3,6,6,5,6,3,6,5, &
838  3,4,6,6,6,6,4,3,4,3,6,6, 5,6,6,3,1,6,5,6,5,3,6,6, &
839  4,3,6,6,4,3,6,6,6,6,4,3, 6,5,3,6,6,1,6,5,3,5,6,6, &
840  6,6,3,4,6,6,3,4,6,6,3,4, 6,3,5,6,5,6,1,6,6,6,5,3, &
841  6,6,4,3,3,4,6,6,3,4,6,6, 3,6,6,5,6,5,6,1,6,6,3,5, &
842  4,3,6,6,4,3,6,6,6,6,4,3, 5,6,3,6,5,3,6,6,1,6,6,5, &
843  3,4,6,6,6,6,4,3,4,3,6,6, 6,5,6,3,3,5,6,6,6,1,5,6, &
844  6,6,4,3,3,4,6,6,3,4,6,6, 3,6,5,6,6,6,5,3,6,5,1,6, &
845  6,6,3,4,6,6,3,4,6,6,3,4, 6,3,6,5,6,6,3,5,5,6,6,1 &
846  ],shape(bcc_interactionslipslip))
853 
854  integer, dimension(HEX_NSLIP,HEX_NSLIP), parameter :: &
855  hex_interactionslipslip = reshape( [&
856  1, 2, 2, 3, 3, 3, 7, 7, 7, 13,13,13,13,13,13, 21,21,21,21,21,21,21,21,21,21,21,21, 31,31,31,31,31,31, & ! -----> acting
857  2, 1, 2, 3, 3, 3, 7, 7, 7, 13,13,13,13,13,13, 21,21,21,21,21,21,21,21,21,21,21,21, 31,31,31,31,31,31, & ! |
858  2, 2, 1, 3, 3, 3, 7, 7, 7, 13,13,13,13,13,13, 21,21,21,21,21,21,21,21,21,21,21,21, 31,31,31,31,31,31, & ! |
859  ! ! v
860  6, 6, 6, 4, 5, 5, 8, 8, 8, 14,14,14,14,14,14, 22,22,22,22,22,22,22,22,22,22,22,22, 32,32,32,32,32,32, & ! reacting
861  6, 6, 6, 5, 4, 5, 8, 8, 8, 14,14,14,14,14,14, 22,22,22,22,22,22,22,22,22,22,22,22, 32,32,32,32,32,32, &
862  6, 6, 6, 5, 5, 4, 8, 8, 8, 14,14,14,14,14,14, 22,22,22,22,22,22,22,22,22,22,22,22, 32,32,32,32,32,32, &
863  !
864  12,12,12, 11,11,11, 9,10,10, 15,15,15,15,15,15, 23,23,23,23,23,23,23,23,23,23,23,23, 33,33,33,33,33,33, &
865  12,12,12, 11,11,11, 10, 9,10, 15,15,15,15,15,15, 23,23,23,23,23,23,23,23,23,23,23,23, 33,33,33,33,33,33, &
866  12,12,12, 11,11,11, 10,10, 9, 15,15,15,15,15,15, 23,23,23,23,23,23,23,23,23,23,23,23, 33,33,33,33,33,33, &
867  !
868  20,20,20, 19,19,19, 18,18,18, 16,17,17,17,17,17, 24,24,24,24,24,24,24,24,24,24,24,24, 34,34,34,34,34,34, &
869  20,20,20, 19,19,19, 18,18,18, 17,16,17,17,17,17, 24,24,24,24,24,24,24,24,24,24,24,24, 34,34,34,34,34,34, &
870  20,20,20, 19,19,19, 18,18,18, 17,17,16,17,17,17, 24,24,24,24,24,24,24,24,24,24,24,24, 34,34,34,34,34,34, &
871  20,20,20, 19,19,19, 18,18,18, 17,17,17,16,17,17, 24,24,24,24,24,24,24,24,24,24,24,24, 34,34,34,34,34,34, &
872  20,20,20, 19,19,19, 18,18,18, 17,17,17,17,16,17, 24,24,24,24,24,24,24,24,24,24,24,24, 34,34,34,34,34,34, &
873  20,20,20, 19,19,19, 18,18,18, 17,17,17,17,17,16, 24,24,24,24,24,24,24,24,24,24,24,24, 34,34,34,34,34,34, &
874  !
875  30,30,30, 29,29,29, 28,28,28, 27,27,27,27,27,27, 25,26,26,26,26,26,26,26,26,26,26,26, 35,35,35,35,35,35, &
876  30,30,30, 29,29,29, 28,28,28, 27,27,27,27,27,27, 26,25,26,26,26,26,26,26,26,26,26,26, 35,35,35,35,35,35, &
877  30,30,30, 29,29,29, 28,28,28, 27,27,27,27,27,27, 26,26,25,26,26,26,26,26,26,26,26,26, 35,35,35,35,35,35, &
878  30,30,30, 29,29,29, 28,28,28, 27,27,27,27,27,27, 26,26,26,25,26,26,26,26,26,26,26,26, 35,35,35,35,35,35, &
879  30,30,30, 29,29,29, 28,28,28, 27,27,27,27,27,27, 26,26,26,26,25,26,26,26,26,26,26,26, 35,35,35,35,35,35, &
880  30,30,30, 29,29,29, 28,28,28, 27,27,27,27,27,27, 26,26,26,26,26,25,26,26,26,26,26,26, 35,35,35,35,35,35, &
881  30,30,30, 29,29,29, 28,28,28, 27,27,27,27,27,27, 26,26,26,26,26,26,25,26,26,26,26,26, 35,35,35,35,35,35, &
882  30,30,30, 29,29,29, 28,28,28, 27,27,27,27,27,27, 26,26,26,26,26,26,26,25,26,26,26,26, 35,35,35,35,35,35, &
883  30,30,30, 29,29,29, 28,28,28, 27,27,27,27,27,27, 26,26,26,26,26,26,26,26,25,26,26,26, 35,35,35,35,35,35, &
884  30,30,30, 29,29,29, 28,28,28, 27,27,27,27,27,27, 26,26,26,26,26,26,26,26,26,25,26,26, 35,35,35,35,35,35, &
885  30,30,30, 29,29,29, 28,28,28, 27,27,27,27,27,27, 26,26,26,26,26,26,26,26,26,26,25,26, 35,35,35,35,35,35, &
886  30,30,30, 29,29,29, 28,28,28, 27,27,27,27,27,27, 26,26,26,26,26,26,26,26,26,26,26,25, 35,35,35,35,35,35, &
887  !
888  42,42,42, 41,41,41, 40,40,40, 39,39,39,39,39,39, 38,38,38,38,38,38,38,38,38,38,38,38, 36,37,37,37,37,37, &
889  42,42,42, 41,41,41, 40,40,40, 39,39,39,39,39,39, 38,38,38,38,38,38,38,38,38,38,38,38, 37,36,37,37,37,37, &
890  42,42,42, 41,41,41, 40,40,40, 39,39,39,39,39,39, 38,38,38,38,38,38,38,38,38,38,38,38, 37,37,36,37,37,37, &
891  42,42,42, 41,41,41, 40,40,40, 39,39,39,39,39,39, 38,38,38,38,38,38,38,38,38,38,38,38, 37,37,37,36,37,37, &
892  42,42,42, 41,41,41, 40,40,40, 39,39,39,39,39,39, 38,38,38,38,38,38,38,38,38,38,38,38, 37,37,37,37,36,37, &
893  42,42,42, 41,41,41, 40,40,40, 39,39,39,39,39,39, 38,38,38,38,38,38,38,38,38,38,38,38, 37,37,37,37,37,36 &
894  ],shape(hex_interactionslipslip))
895 
896  integer, dimension(BCT_NSLIP,BCT_NSLIP), parameter :: &
897  bct_interactionslipslip = reshape( [&
898  1, 2, 3, 3, 7, 7, 13, 13, 13, 13, 21, 21, 31, 31, 31, 31, 43, 43, 57, 57, 73, 73, 73, 73, 91, 91, 91, 91, 91, 91, 91, 91, 111, 111, 111, 111, 133,133,133,133,133,133,133,133, 157,157,157,157,157,157,157,157, & ! -----> acting
899  2, 1, 3, 3, 7, 7, 13, 13, 13, 13, 21, 21, 31, 31, 31, 31, 43, 43, 57, 57, 73, 73, 73, 73, 91, 91, 91, 91, 91, 91, 91, 91, 111, 111, 111, 111, 133,133,133,133,133,133,133,133, 157,157,157,157,157,157,157,157, & ! |
900  ! |
901  6, 6, 4, 5, 8, 8, 14, 14, 14, 14, 22, 22, 32, 32, 32, 32, 44, 44, 58, 58, 74, 74, 74, 74, 92, 92, 92, 92, 92, 92, 92, 92, 112, 112, 112, 112, 134,134,134,134,134,134,134,134, 158,158,158,158,158,158,158,158, & ! v
902  6, 6, 5, 4, 8, 8, 14, 14, 14, 14, 22, 22, 32, 32, 32, 32, 44, 44, 58, 58, 74, 74, 74, 74, 92, 92, 92, 92, 92, 92, 92, 92, 112, 112, 112, 112, 134,134,134,134,134,134,134,134, 158,158,158,158,158,158,158,158, & ! reacting
903  !
904  12, 12, 11, 11, 9, 10, 15, 15, 15, 15, 23, 23, 33, 33, 33, 33, 45, 45, 59, 59, 75, 75, 75, 75, 93, 93, 93, 93, 93, 93, 93, 93, 113, 113, 113, 113, 135,135,135,135,135,135,135,135, 159,159,159,159,159,159,159,159, &
905  12, 12, 11, 11, 10, 9, 15, 15, 15, 15, 23, 23, 33, 33, 33, 33, 45, 45, 59, 59, 75, 75, 75, 75, 93, 93, 93, 93, 93, 93, 93, 93, 113, 113, 113, 113, 135,135,135,135,135,135,135,135, 159,159,159,159,159,159,159,159, &
906  !
907  20, 20, 19, 19, 18, 18, 16, 17, 17, 17, 24, 24, 34, 34, 34, 34, 46, 46, 60, 60, 76, 76, 76, 76, 94, 94, 94, 94, 94, 94, 94, 94, 114, 114, 114, 114, 136,136,136,136,136,136,136,136, 160,160,160,160,160,160,160,160, &
908  20, 20, 19, 19, 18, 18, 17, 16, 17, 17, 24, 24, 34, 34, 34, 34, 46, 46, 60, 60, 76, 76, 76, 76, 94, 94, 94, 94, 94, 94, 94, 94, 114, 114, 114, 114, 136,136,136,136,136,136,136,136, 160,160,160,160,160,160,160,160, &
909  20, 20, 19, 19, 18, 18, 17, 17, 16, 17, 24, 24, 34, 34, 34, 34, 46, 46, 60, 60, 76, 76, 76, 76, 94, 94, 94, 94, 94, 94, 94, 94, 114, 114, 114, 114, 136,136,136,136,136,136,136,136, 160,160,160,160,160,160,160,160, &
910  20, 20, 19, 19, 18, 18, 17, 17, 17, 16, 24, 24, 34, 34, 34, 34, 46, 46, 60, 60, 76, 76, 76, 76, 94, 94, 94, 94, 94, 94, 94, 94, 114, 114, 114, 114, 136,136,136,136,136,136,136,136, 160,160,160,160,160,160,160,160, &
911  !
912  30, 30, 29, 29, 28, 28, 27, 27, 27, 27, 25, 26, 35, 35, 35, 35, 47, 47, 61, 61, 77, 77, 77, 77, 95, 95, 95, 95, 95, 95, 95, 95, 115, 115, 115, 115, 137,137,137,137,137,137,137,137, 161,161,161,161,161,161,161,161, &
913  30, 30, 29, 29, 28, 28, 27, 27, 27, 27, 26, 25, 35, 35, 35, 35, 47, 47, 61, 61, 77, 77, 77, 77, 95, 95, 95, 95, 95, 95, 95, 95, 115, 115, 115, 115, 137,137,137,137,137,137,137,137, 161,161,161,161,161,161,161,161, &
914  !
915  42, 42, 41, 41, 40, 40, 39, 39, 39, 39, 38, 38, 36, 37, 37, 37, 48, 48, 62, 62, 78, 78, 78, 78, 96, 96, 96, 96, 96, 96, 96, 96, 116, 116, 116, 116, 138,138,138,138,138,138,138,138, 162,162,162,162,162,162,162,162, &
916  42, 42, 41, 41, 40, 40, 39, 39, 39, 39, 38, 38, 37, 36, 37, 37, 48, 48, 62, 62, 78, 78, 78, 78, 96, 96, 96, 96, 96, 96, 96, 96, 116, 116, 116, 116, 138,138,138,138,138,138,138,138, 162,162,162,162,162,162,162,162, &
917  42, 42, 41, 41, 40, 40, 39, 39, 39, 39, 38, 38, 37, 37, 36, 37, 48, 48, 62, 62, 78, 78, 78, 78, 96, 96, 96, 96, 96, 96, 96, 96, 116, 116, 116, 116, 138,138,138,138,138,138,138,138, 162,162,162,162,162,162,162,162, &
918  42, 42, 41, 41, 40, 40, 39, 39, 39, 39, 38, 38, 37, 37, 37, 36, 48, 48, 62, 62, 78, 78, 78, 78, 96, 96, 96, 96, 96, 96, 96, 96, 116, 116, 116, 116, 138,138,138,138,138,138,138,138, 162,162,162,162,162,162,162,162, &
919  !
920  56, 56, 55, 55, 54, 54, 53, 53, 53, 53, 52, 52, 51, 51, 51, 51, 49, 50, 63, 63, 79, 79, 79, 79, 97, 97, 97, 97, 97, 97, 97, 97, 117, 117, 117, 117, 139,139,139,139,139,139,139,139, 163,163,163,163,163,163,163,163, &
921  56, 56, 55, 55, 54, 54, 53, 53, 53, 53, 52, 52, 51, 51, 51, 51, 50, 49, 63, 63, 79, 79, 79, 79, 97, 97, 97, 97, 97, 97, 97, 97, 117, 117, 117, 117, 139,139,139,139,139,139,139,139, 163,163,163,163,163,163,163,163, &
922  !
923  72, 72, 71, 71, 70, 70, 69, 69, 69, 69, 68, 68, 67, 67, 67, 67, 66, 66, 64, 65, 80, 80, 80, 80, 98, 98, 98, 98, 98, 98, 98, 98, 118, 118, 118, 118, 140,140,140,140,140,140,140,140, 164,164,164,164,164,164,164,164, &
924  72, 72, 71, 71, 70, 70, 69, 69, 69, 69, 68, 68, 67, 67, 67, 67, 66, 66, 65, 64, 80, 80, 80, 80, 98, 98, 98, 98, 98, 98, 98, 98, 118, 118, 118, 118, 140,140,140,140,140,140,140,140, 164,164,164,164,164,164,164,164, &
925  !
926  90, 90, 89, 89, 88, 88, 87, 87, 87, 87, 86, 86, 85, 85, 85, 85, 84, 84, 83, 83, 81, 82, 82, 82, 99, 99, 99, 99, 99, 99, 99, 99, 119, 119, 119, 119, 141,141,141,141,141,141,141,141, 165,165,165,165,165,165,165,165, &
927  90, 90, 89, 89, 88, 88, 87, 87, 87, 87, 86, 86, 85, 85, 85, 85, 84, 84, 83, 83, 82, 81, 82, 82, 99, 99, 99, 99, 99, 99, 99, 99, 119, 119, 119, 119, 141,141,141,141,141,141,141,141, 165,165,165,165,165,165,165,165, &
928  90, 90, 89, 89, 88, 88, 87, 87, 87, 87, 86, 86, 85, 85, 85, 85, 84, 84, 83, 83, 82, 82, 81, 82, 99, 99, 99, 99, 99, 99, 99, 99, 119, 119, 119, 119, 141,141,141,141,141,141,141,141, 165,165,165,165,165,165,165,165, &
929  90, 90, 89, 89, 88, 88, 87, 87, 87, 87, 86, 86, 85, 85, 85, 85, 84, 84, 83, 83, 82, 82, 82, 81, 99, 99, 99, 99, 99, 99, 99, 99, 119, 119, 119, 119, 141,141,141,141,141,141,141,141, 165,165,165,165,165,165,165,165, &
930  !
931  110,110, 109,109, 108,108, 107,107,107,107, 106,106, 105,105,105,105, 104,104, 103,103, 102,102,102,102, 100,101,101,101,101,101,101,101, 120, 120, 120, 120, 142,142,142,142,142,142,142,142, 166,166,166,166,166,166,166,166, &
932  110,110, 109,109, 108,108, 107,107,107,107, 106,106, 105,105,105,105, 104,104, 103,103, 102,102,102,102, 101,100,101,101,101,101,101,101, 120, 120, 120, 120, 142,142,142,142,142,142,142,142, 166,166,166,166,166,166,166,166, &
933  110,110, 109,109, 108,108, 107,107,107,107, 106,106, 105,105,105,105, 104,104, 103,103, 102,102,102,102, 101,101,100,101,101,101,101,101, 120, 120, 120, 120, 142,142,142,142,142,142,142,142, 166,166,166,166,166,166,166,166, &
934  110,110, 109,109, 108,108, 107,107,107,107, 106,106, 105,105,105,105, 104,104, 103,103, 102,102,102,102, 101,101,101,100,101,101,101,101, 120, 120, 120, 120, 142,142,142,142,142,142,142,142, 166,166,166,166,166,166,166,166, &
935  110,110, 109,109, 108,108, 107,107,107,107, 106,106, 105,105,105,105, 104,104, 103,103, 102,102,102,102, 101,101,101,101,100,101,101,101, 120, 120, 120, 120, 142,142,142,142,142,142,142,142, 166,166,166,166,166,166,166,166, &
936  110,110, 109,109, 108,108, 107,107,107,107, 106,106, 105,105,105,105, 104,104, 103,103, 102,102,102,102, 101,101,101,101,101,100,101,101, 120, 120, 120, 120, 142,142,142,142,142,142,142,142, 166,166,166,166,166,166,166,166, &
937  110,110, 109,109, 108,108, 107,107,107,107, 106,106, 105,105,105,105, 104,104, 103,103, 102,102,102,102, 101,101,101,101,101,101,100,101, 120, 120, 120, 120, 142,142,142,142,142,142,142,142, 166,166,166,166,166,166,166,166, &
938  110,110, 109,109, 108,108, 107,107,107,107, 106,106, 105,105,105,105, 104,104, 103,103, 102,102,102,102, 101,101,101,101,101,101,101,100, 120, 120, 120, 120, 142,142,142,142,142,142,142,142, 166,166,166,166,166,166,166,166, &
939  !
940  132,132, 131,131, 130,130, 129,129,129,129, 128,128, 127,127,127,127, 126,126, 125,125, 124,124,124,124, 123,123,123,123,123,123,123,123, 121, 122, 122, 122, 143,143,143,143,143,143,143,143, 167,167,167,167,167,167,167,167, &
941  132,132, 131,131, 130,130, 129,129,129,129, 128,128, 127,127,127,127, 126,126, 125,125, 124,124,124,124, 123,123,123,123,123,123,123,123, 121, 121, 122, 122, 143,143,143,143,143,143,143,143, 167,167,167,167,167,167,167,167, &
942  132,132, 131,131, 130,130, 129,129,129,129, 128,128, 127,127,127,127, 126,126, 125,125, 124,124,124,124, 123,123,123,123,123,123,123,123, 121, 122, 121, 122, 143,143,143,143,143,143,143,143, 167,167,167,167,167,167,167,167, &
943  132,132, 131,131, 130,130, 129,129,129,129, 128,128, 127,127,127,127, 126,126, 125,125, 124,124,124,124, 123,123,123,123,123,123,123,123, 121, 122, 122, 121, 143,143,143,143,143,143,143,143, 167,167,167,167,167,167,167,167, &
944  !
945  156,156, 155,155, 154,154, 153,153,153,153, 152,152, 151,151,151,151, 150,150, 149,149, 148,148,148,148, 147,147,147,147,147,147,147,147, 146, 146, 146, 146, 144,145,145,145,145,145,145,145, 168,168,168,168,168,168,168,168, &
946  156,156, 155,155, 154,154, 153,153,153,153, 152,152, 151,151,151,151, 150,150, 149,149, 148,148,148,148, 147,147,147,147,147,147,147,147, 146, 146, 146, 146, 145,144,145,145,145,145,145,145, 168,168,168,168,168,168,168,168, &
947  156,156, 155,155, 154,154, 153,153,153,153, 152,152, 151,151,151,151, 150,150, 149,149, 148,148,148,148, 147,147,147,147,147,147,147,147, 146, 146, 146, 146, 145,145,144,145,145,145,145,145, 168,168,168,168,168,168,168,168, &
948  156,156, 155,155, 154,154, 153,153,153,153, 152,152, 151,151,151,151, 150,150, 149,149, 148,148,148,148, 147,147,147,147,147,147,147,147, 146, 146, 146, 146, 145,145,145,144,145,145,145,145, 168,168,168,168,168,168,168,168, &
949  156,156, 155,155, 154,154, 153,153,153,153, 152,152, 151,151,151,151, 150,150, 149,149, 148,148,148,148, 147,147,147,147,147,147,147,147, 146, 146, 146, 146, 145,145,145,145,144,145,145,145, 168,168,168,168,168,168,168,168, &
950  156,156, 155,155, 154,154, 153,153,153,153, 152,152, 151,151,151,151, 150,150, 149,149, 148,148,148,148, 147,147,147,147,147,147,147,147, 146, 146, 146, 146, 145,145,145,145,145,144,145,145, 168,168,168,168,168,168,168,168, &
951  156,156, 155,155, 154,154, 153,153,153,153, 152,152, 151,151,151,151, 150,150, 149,149, 148,148,148,148, 147,147,147,147,147,147,147,147, 146, 146, 146, 146, 145,145,145,145,145,145,144,145, 168,168,168,168,168,168,168,168, &
952  156,156, 155,155, 154,154, 153,153,153,153, 152,152, 151,151,151,151, 150,150, 149,149, 148,148,148,148, 147,147,147,147,147,147,147,147, 146, 146, 146, 146, 145,145,145,145,145,145,145,144, 168,168,168,168,168,168,168,168, &
953  !
954  182,182, 181,181, 180,180, 179,179,179,179, 178,178, 177,177,177,177, 176,176, 175,175, 174,174,174,174, 173,173,173,173,173,173,173,173, 172, 172, 172, 172, 171,171,171,171,171,171,171,171, 169,170,170,170,170,170,170,170, &
955  182,182, 181,181, 180,180, 179,179,179,179, 178,178, 177,177,177,177, 176,176, 175,175, 174,174,174,174, 173,173,173,173,173,173,173,173, 172, 172, 172, 172, 171,171,171,171,171,171,171,171, 170,169,170,170,170,170,170,170, &
956  182,182, 181,181, 180,180, 179,179,179,179, 178,178, 177,177,177,177, 176,176, 175,175, 174,174,174,174, 173,173,173,173,173,173,173,173, 172, 172, 172, 172, 171,171,171,171,171,171,171,171, 170,170,169,170,170,170,170,170, &
957  182,182, 181,181, 180,180, 179,179,179,179, 178,178, 177,177,177,177, 176,176, 175,175, 174,174,174,174, 173,173,173,173,173,173,173,173, 172, 172, 172, 172, 171,171,171,171,171,171,171,171, 170,170,170,169,170,170,170,170, &
958  182,182, 181,181, 180,180, 179,179,179,179, 178,178, 177,177,177,177, 176,176, 175,175, 174,174,174,174, 173,173,173,173,173,173,173,173, 172, 172, 172, 172, 171,171,171,171,171,171,171,171, 170,170,170,170,169,170,170,170, &
959  182,182, 181,181, 180,180, 179,179,179,179, 178,178, 177,177,177,177, 176,176, 175,175, 174,174,174,174, 173,173,173,173,173,173,173,173, 172, 172, 172, 172, 171,171,171,171,171,171,171,171, 169,170,170,170,170,169,170,170, &
960  182,182, 181,181, 180,180, 179,179,179,179, 178,178, 177,177,177,177, 176,176, 175,175, 174,174,174,174, 173,173,173,173,173,173,173,173, 172, 172, 172, 172, 171,171,171,171,171,171,171,171, 169,170,170,170,170,170,169,170, &
961  182,182, 181,181, 180,180, 179,179,179,179, 178,178, 177,177,177,177, 176,176, 175,175, 174,174,174,174, 173,173,173,173,173,173,173,173, 172, 172, 172, 172, 171,171,171,171,171,171,171,171, 169,170,170,170,170,170,170,169 &
962  ],shape(bct_interactionslipslip))
963 
964 
965  if (len_trim(structure) /= 3) &
966  call io_error(137,ext_msg='lattice_interaction_SlipBySlip: '//trim(structure))
967 
968  select case(structure)
969  case('fcc')
970  interactiontypes = fcc_interactionslipslip
971  nslipmax = fcc_nslipsystem
972  case('bcc')
973  interactiontypes = bcc_interactionslipslip
974  nslipmax = bcc_nslipsystem
975  case('hex')
976  interactiontypes = hex_interactionslipslip
977  nslipmax = hex_nslipsystem
978  case('bct')
979  interactiontypes = bct_interactionslipslip
980  nslipmax = bct_nslipsystem
981  case default
982  call io_error(137,ext_msg='lattice_interaction_SlipBySlip: '//trim(structure))
983  end select
984 
985  interactionmatrix = buildinteraction(nslip,nslip,nslipmax,nslipmax,interactionvalues,interactiontypes)
986 
988 
989 
990 !--------------------------------------------------------------------------------------------------
993 !--------------------------------------------------------------------------------------------------
994 function lattice_interaction_twinbytwin(Ntwin,interactionValues,structure) result(interactionMatrix)
995 
996  integer, dimension(:), intent(in) :: ntwin
997  real(preal), dimension(:), intent(in) :: interactionvalues
998  character(len=*), intent(in) :: structure
999  real(preal), dimension(sum(Ntwin),sum(Ntwin)) :: interactionmatrix
1000 
1001  integer, dimension(:), allocatable :: ntwinmax
1002  integer, dimension(:,:), allocatable :: interactiontypes
1003 
1004  integer, dimension(FCC_NTWIN,FCC_NTWIN), parameter :: &
1005  fcc_interactiontwintwin = reshape( [&
1006  1,1,1,2,2,2,2,2,2,2,2,2, & ! -----> acting
1007  1,1,1,2,2,2,2,2,2,2,2,2, & ! |
1008  1,1,1,2,2,2,2,2,2,2,2,2, & ! |
1009  2,2,2,1,1,1,2,2,2,2,2,2, & ! v
1010  2,2,2,1,1,1,2,2,2,2,2,2, & ! reacting
1011  2,2,2,1,1,1,2,2,2,2,2,2, &
1012  2,2,2,2,2,2,1,1,1,2,2,2, &
1013  2,2,2,2,2,2,1,1,1,2,2,2, &
1014  2,2,2,2,2,2,1,1,1,2,2,2, &
1015  2,2,2,2,2,2,2,2,2,1,1,1, &
1016  2,2,2,2,2,2,2,2,2,1,1,1, &
1017  2,2,2,2,2,2,2,2,2,1,1,1 &
1018  ],shape(fcc_interactiontwintwin))
1019 
1020  integer, dimension(BCC_NTWIN,BCC_NTWIN), parameter :: &
1021  bcc_interactiontwintwin = reshape( [&
1022  1,3,3,3,3,3,3,2,3,3,2,3, & ! -----> acting
1023  3,1,3,3,3,3,2,3,3,3,3,2, & ! |
1024  3,3,1,3,3,2,3,3,2,3,3,3, & ! |
1025  3,3,3,1,2,3,3,3,3,2,3,3, & ! v
1026  3,3,3,2,1,3,3,3,3,2,3,3, & ! reacting
1027  3,3,2,3,3,1,3,3,2,3,3,3, &
1028  3,2,3,3,3,3,1,3,3,3,3,2, &
1029  2,3,3,3,3,3,3,1,3,3,2,3, &
1030  3,3,2,3,3,2,3,3,1,3,3,3, &
1031  3,3,3,2,2,3,3,3,3,1,3,3, &
1032  2,3,3,3,3,3,3,2,3,3,1,3, &
1033  3,2,3,3,3,3,2,3,3,3,3,1 &
1034  ],shape(bcc_interactiontwintwin))
1038  integer, dimension(HEX_NTWIN,HEX_NTWIN), parameter :: &
1039  hex_interactiontwintwin = reshape( [&
1040  1, 2, 2, 2, 2, 2, 3, 3, 3, 3, 3, 3, 7, 7, 7, 7, 7, 7, 13,13,13,13,13,13, & ! -----> acting
1041  2, 1, 2, 2, 2, 2, 3, 3, 3, 3, 3, 3, 7, 7, 7, 7, 7, 7, 13,13,13,13,13,13, & ! |
1042  2, 2, 1, 2, 2, 2, 3, 3, 3, 3, 3, 3, 7, 7, 7, 7, 7, 7, 13,13,13,13,13,13, & ! |
1043  2, 2, 2, 1, 2, 2, 3, 3, 3, 3, 3, 3, 7, 7, 7, 7, 7, 7, 13,13,13,13,13,13, & ! v
1044  2, 2, 2, 2, 1, 2, 3, 3, 3, 3, 3, 3, 7, 7, 7, 7, 7, 7, 13,13,13,13,13,13, & ! reacting
1045  2, 2, 2, 2, 2, 1, 3, 3, 3, 3, 3, 3, 7, 7, 7, 7, 7, 7, 13,13,13,13,13,13, &
1046  !
1047  6, 6, 6, 6, 6, 6, 4, 5, 5, 5, 5, 5, 8, 8, 8, 8, 8, 8, 14,14,14,14,14,14, &
1048  6, 6, 6, 6, 6, 6, 5, 4, 5, 5, 5, 5, 8, 8, 8, 8, 8, 8, 14,14,14,14,14,14, &
1049  6, 6, 6, 6, 6, 6, 5, 5, 4, 5, 5, 5, 8, 8, 8, 8, 8, 8, 14,14,14,14,14,14, &
1050  6, 6, 6, 6, 6, 6, 5, 5, 5, 4, 5, 5, 8, 8, 8, 8, 8, 8, 14,14,14,14,14,14, &
1051  6, 6, 6, 6, 6, 6, 5, 5, 5, 5, 4, 5, 8, 8, 8, 8, 8, 8, 14,14,14,14,14,14, &
1052  6, 6, 6, 6, 6, 6, 5, 5, 5, 5, 5, 4, 8, 8, 8, 8, 8, 8, 14,14,14,14,14,14, &
1053  !
1054  12,12,12,12,12,12, 11,11,11,11,11,11, 9,10,10,10,10,10, 15,15,15,15,15,15, &
1055  12,12,12,12,12,12, 11,11,11,11,11,11, 10, 9,10,10,10,10, 15,15,15,15,15,15, &
1056  12,12,12,12,12,12, 11,11,11,11,11,11, 10,10, 9,10,10,10, 15,15,15,15,15,15, &
1057  12,12,12,12,12,12, 11,11,11,11,11,11, 10,10,10, 9,10,10, 15,15,15,15,15,15, &
1058  12,12,12,12,12,12, 11,11,11,11,11,11, 10,10,10,10, 9,10, 15,15,15,15,15,15, &
1059  12,12,12,12,12,12, 11,11,11,11,11,11, 10,10,10,10,10, 9, 15,15,15,15,15,15, &
1060  !
1061  20,20,20,20,20,20, 19,19,19,19,19,19, 18,18,18,18,18,18, 16,17,17,17,17,17, &
1062  20,20,20,20,20,20, 19,19,19,19,19,19, 18,18,18,18,18,18, 17,16,17,17,17,17, &
1063  20,20,20,20,20,20, 19,19,19,19,19,19, 18,18,18,18,18,18, 17,17,16,17,17,17, &
1064  20,20,20,20,20,20, 19,19,19,19,19,19, 18,18,18,18,18,18, 17,17,17,16,17,17, &
1065  20,20,20,20,20,20, 19,19,19,19,19,19, 18,18,18,18,18,18, 17,17,17,17,16,17, &
1066  20,20,20,20,20,20, 19,19,19,19,19,19, 18,18,18,18,18,18, 17,17,17,17,17,16 &
1067  ],shape(hex_interactiontwintwin))
1068 
1069  if (len_trim(structure) /= 3) &
1070  call io_error(137,ext_msg='lattice_interaction_TwinByTwin: '//trim(structure))
1071 
1072  select case(structure)
1073  case('fcc')
1074  interactiontypes = fcc_interactiontwintwin
1075  ntwinmax = fcc_ntwinsystem
1076  case('bcc')
1077  interactiontypes = bcc_interactiontwintwin
1078  ntwinmax = bcc_ntwinsystem
1079  case('hex')
1080  interactiontypes = hex_interactiontwintwin
1081  ntwinmax = hex_ntwinsystem
1082  case default
1083  call io_error(137,ext_msg='lattice_interaction_TwinByTwin: '//trim(structure))
1084  end select
1085 
1086  interactionmatrix = buildinteraction(ntwin,ntwin,ntwinmax,ntwinmax,interactionvalues,interactiontypes)
1087 
1088 end function lattice_interaction_twinbytwin
1089 
1090 
1091 !--------------------------------------------------------------------------------------------------
1094 !--------------------------------------------------------------------------------------------------
1095 function lattice_interaction_transbytrans(Ntrans,interactionValues,structure) result(interactionMatrix)
1097  integer, dimension(:), intent(in) :: ntrans
1098  real(preal), dimension(:), intent(in) :: interactionvalues
1099  character(len=*), intent(in) :: structure
1100  real(preal), dimension(sum(Ntrans),sum(Ntrans)) :: interactionmatrix
1101 
1102  integer, dimension(:), allocatable :: ntransmax
1103  integer, dimension(:,:), allocatable :: interactiontypes
1104 
1105  integer, dimension(FCC_NTRANS,FCC_NTRANS), parameter :: &
1106  fcc_interactiontranstrans = reshape( [&
1107  1,1,1,2,2,2,2,2,2,2,2,2, & ! -----> acting
1108  1,1,1,2,2,2,2,2,2,2,2,2, & ! |
1109  1,1,1,2,2,2,2,2,2,2,2,2, & ! |
1110  2,2,2,1,1,1,2,2,2,2,2,2, & ! v
1111  2,2,2,1,1,1,2,2,2,2,2,2, & ! reacting
1112  2,2,2,1,1,1,2,2,2,2,2,2, &
1113  2,2,2,2,2,2,1,1,1,2,2,2, &
1114  2,2,2,2,2,2,1,1,1,2,2,2, &
1115  2,2,2,2,2,2,1,1,1,2,2,2, &
1116  2,2,2,2,2,2,2,2,2,1,1,1, &
1117  2,2,2,2,2,2,2,2,2,1,1,1, &
1118  2,2,2,2,2,2,2,2,2,1,1,1 &
1119  ],shape(fcc_interactiontranstrans))
1120 
1121  if (len_trim(structure) /= 3) &
1122  call io_error(137,ext_msg='lattice_interaction_TransByTrans: '//trim(structure))
1123 
1124  if(structure == 'fcc') then
1125  interactiontypes = fcc_interactiontranstrans
1126  ntransmax = fcc_ntranssystem
1127  else
1128  call io_error(137,ext_msg='lattice_interaction_TransByTrans: '//trim(structure))
1129  end if
1130 
1131  interactionmatrix = buildinteraction(ntrans,ntrans,ntransmax,ntransmax,interactionvalues,interactiontypes)
1132 
1134 
1135 
1136 !--------------------------------------------------------------------------------------------------
1139 !--------------------------------------------------------------------------------------------------
1140 function lattice_interaction_slipbytwin(Nslip,Ntwin,interactionValues,structure) result(interactionMatrix)
1142  integer, dimension(:), intent(in) :: nslip, & !< number of active slip systems per family
1143  ntwin
1144  real(preal), dimension(:), intent(in) :: interactionvalues
1145  character(len=*), intent(in) :: structure
1146  real(preal), dimension(sum(Nslip),sum(Ntwin)) :: interactionmatrix
1147 
1148  integer, dimension(:), allocatable :: nslipmax, &
1149  ntwinmax
1150  integer, dimension(:,:), allocatable :: interactiontypes
1151 
1152  integer, dimension(FCC_NTWIN,FCC_NSLIP), parameter :: &
1153  fcc_interactionsliptwin = reshape( [&
1154  1,1,1,3,3,3,2,2,2,3,3,3, & ! -----> twin (acting)
1155  1,1,1,3,3,3,3,3,3,2,2,2, & ! |
1156  1,1,1,2,2,2,3,3,3,3,3,3, & ! |
1157  3,3,3,1,1,1,3,3,3,2,2,2, & ! v
1158  3,3,3,1,1,1,2,2,2,3,3,3, & ! slip (reacting)
1159  2,2,2,1,1,1,3,3,3,3,3,3, &
1160  2,2,2,3,3,3,1,1,1,3,3,3, &
1161  3,3,3,2,2,2,1,1,1,3,3,3, &
1162  3,3,3,3,3,3,1,1,1,2,2,2, &
1163  3,3,3,2,2,2,3,3,3,1,1,1, &
1164  2,2,2,3,3,3,3,3,3,1,1,1, &
1165  3,3,3,3,3,3,2,2,2,1,1,1, &
1166 
1167  4,4,4,4,4,4,4,4,4,4,4,4, &
1168  4,4,4,4,4,4,4,4,4,4,4,4, &
1169  4,4,4,4,4,4,4,4,4,4,4,4, &
1170  4,4,4,4,4,4,4,4,4,4,4,4, &
1171  4,4,4,4,4,4,4,4,4,4,4,4, &
1172  4,4,4,4,4,4,4,4,4,4,4,4 &
1173  ],shape(fcc_interactionsliptwin))
1177  integer, dimension(BCC_NTWIN,BCC_NSLIP), parameter :: &
1178  bcc_interactionsliptwin = reshape( [&
1179  3,3,3,2,2,3,3,3,3,2,3,3, & ! -----> twin (acting)
1180  3,3,2,3,3,2,3,3,2,3,3,3, & ! |
1181  3,2,3,3,3,3,2,3,3,3,3,2, & ! |
1182  2,3,3,3,3,3,3,2,3,3,2,3, & ! v
1183  2,3,3,3,3,3,3,2,3,3,2,3, & ! slip (reacting)
1184  3,3,2,3,3,2,3,3,2,3,3,3, &
1185  3,2,3,3,3,3,2,3,3,3,3,2, &
1186  3,3,3,2,2,3,3,3,3,2,3,3, &
1187  2,3,3,3,3,3,3,2,3,3,2,3, &
1188  3,3,3,2,2,3,3,3,3,2,3,3, &
1189  3,2,3,3,3,3,2,3,3,3,3,2, &
1190  3,3,2,3,3,2,3,3,2,3,3,3, &
1191  !
1192  1,3,3,3,3,3,3,2,3,3,2,3, &
1193  3,1,3,3,3,3,2,3,3,3,3,2, &
1194  3,3,1,3,3,2,3,3,2,3,3,3, &
1195  3,3,3,1,2,3,3,3,3,2,3,3, &
1196  3,3,3,2,1,3,3,3,3,2,3,3, &
1197  3,3,2,3,3,1,3,3,2,3,3,3, &
1198  3,2,3,3,3,3,1,3,3,3,3,2, &
1199  2,3,3,3,3,3,3,1,3,3,2,3, &
1200  3,3,2,3,3,2,3,3,1,3,3,3, &
1201  3,3,3,2,2,3,3,3,3,1,3,3, &
1202  2,3,3,3,3,3,3,2,3,3,1,3, &
1203  3,2,3,3,3,3,2,3,3,3,3,1 &
1204  ],shape(bcc_interactionsliptwin))
1208  integer, dimension(HEX_NTWIN,HEX_NSLIP), parameter :: &
1209  hex_interactionsliptwin = reshape( [&
1210  1, 1, 1, 1, 1, 1, 2, 2, 2, 2, 2, 2, 3, 3, 3, 3, 3, 3, 4, 4, 4, 4, 4, 4, & ! ----> twin (acting)
1211  1, 1, 1, 1, 1, 1, 2, 2, 2, 2, 2, 2, 3, 3, 3, 3, 3, 3, 4, 4, 4, 4, 4, 4, & ! |
1212  1, 1, 1, 1, 1, 1, 2, 2, 2, 2, 2, 2, 3, 3, 3, 3, 3, 3, 4, 4, 4, 4, 4, 4, & ! |
1213  ! v
1214  5, 5, 5, 5, 5, 5, 6, 6, 6, 6, 6, 6, 7, 7, 7, 7, 7, 7, 8, 8, 8, 8, 8, 8, & ! slip (reacting)
1215  5, 5, 5, 5, 5, 5, 6, 6, 6, 6, 6, 6, 7, 7, 7, 7, 7, 7, 8, 8, 8, 8, 8, 8, &
1216  5, 5, 5, 5, 5, 5, 6, 6, 6, 6, 6, 6, 7, 7, 7, 7, 7, 7, 8, 8, 8, 8, 8, 8, &
1217  !
1218  9, 9, 9, 9, 9, 9, 10,10,10,10,10,10, 11,11,11,11,11,11, 12,12,12,12,12,12, &
1219  9, 9, 9, 9, 9, 9, 10,10,10,10,10,10, 11,11,11,11,11,11, 12,12,12,12,12,12, &
1220  9, 9, 9, 9, 9, 9, 10,10,10,10,10,10, 11,11,11,11,11,11, 12,12,12,12,12,12, &
1221  !
1222  13,13,13,13,13,13, 14,14,14,14,14,14, 15,15,15,15,15,15, 16,16,16,16,16,16, &
1223  13,13,13,13,13,13, 14,14,14,14,14,14, 15,15,15,15,15,15, 16,16,16,16,16,16, &
1224  13,13,13,13,13,13, 14,14,14,14,14,14, 15,15,15,15,15,15, 16,16,16,16,16,16, &
1225  13,13,13,13,13,13, 14,14,14,14,14,14, 15,15,15,15,15,15, 16,16,16,16,16,16, &
1226  13,13,13,13,13,13, 14,14,14,14,14,14, 15,15,15,15,15,15, 16,16,16,16,16,16, &
1227  13,13,13,13,13,13, 14,14,14,14,14,14, 15,15,15,15,15,15, 16,16,16,16,16,16, &
1228  !
1229  17,17,17,17,17,17, 18,18,18,18,18,18, 19,19,19,19,19,19, 20,20,20,20,20,20, &
1230  17,17,17,17,17,17, 18,18,18,18,18,18, 19,19,19,19,19,19, 20,20,20,20,20,20, &
1231  17,17,17,17,17,17, 18,18,18,18,18,18, 19,19,19,19,19,19, 20,20,20,20,20,20, &
1232  17,17,17,17,17,17, 18,18,18,18,18,18, 19,19,19,19,19,19, 20,20,20,20,20,20, &
1233  17,17,17,17,17,17, 18,18,18,18,18,18, 19,19,19,19,19,19, 20,20,20,20,20,20, &
1234  17,17,17,17,17,17, 18,18,18,18,18,18, 19,19,19,19,19,19, 20,20,20,20,20,20, &
1235  17,17,17,17,17,17, 18,18,18,18,18,18, 19,19,19,19,19,19, 20,20,20,20,20,20, &
1236  17,17,17,17,17,17, 18,18,18,18,18,18, 19,19,19,19,19,19, 20,20,20,20,20,20, &
1237  17,17,17,17,17,17, 18,18,18,18,18,18, 19,19,19,19,19,19, 20,20,20,20,20,20, &
1238  17,17,17,17,17,17, 18,18,18,18,18,18, 19,19,19,19,19,19, 20,20,20,20,20,20, &
1239  17,17,17,17,17,17, 18,18,18,18,18,18, 19,19,19,19,19,19, 20,20,20,20,20,20, &
1240  17,17,17,17,17,17, 18,18,18,18,18,18, 19,19,19,19,19,19, 20,20,20,20,20,20, &
1241  !
1242  21,21,21,21,21,21, 22,22,22,22,22,22, 23,23,23,23,23,23, 24,24,24,24,24,24, &
1243  21,21,21,21,21,21, 22,22,22,22,22,22, 23,23,23,23,23,23, 24,24,24,24,24,24, &
1244  21,21,21,21,21,21, 22,22,22,22,22,22, 23,23,23,23,23,23, 24,24,24,24,24,24, &
1245  21,21,21,21,21,21, 22,22,22,22,22,22, 23,23,23,23,23,23, 24,24,24,24,24,24, &
1246  21,21,21,21,21,21, 22,22,22,22,22,22, 23,23,23,23,23,23, 24,24,24,24,24,24, &
1247  21,21,21,21,21,21, 22,22,22,22,22,22, 23,23,23,23,23,23, 24,24,24,24,24,24 &
1248  !
1249  ],shape(hex_interactionsliptwin))
1250 
1251  if (len_trim(structure) /= 3) &
1252  call io_error(137,ext_msg='lattice_interaction_SlipByTwin: '//trim(structure))
1253 
1254  select case(structure)
1255  case('fcc')
1256  interactiontypes = fcc_interactionsliptwin
1257  nslipmax = fcc_nslipsystem
1258  ntwinmax = fcc_ntwinsystem
1259  case('bcc')
1260  interactiontypes = bcc_interactionsliptwin
1261  nslipmax = bcc_nslipsystem
1262  ntwinmax = bcc_ntwinsystem
1263  case('hex')
1264  interactiontypes = hex_interactionsliptwin
1265  nslipmax = hex_nslipsystem
1266  ntwinmax = hex_ntwinsystem
1267  case default
1268  call io_error(137,ext_msg='lattice_interaction_SlipByTwin: '//trim(structure))
1269  end select
1270 
1271  interactionmatrix = buildinteraction(nslip,ntwin,nslipmax,ntwinmax,interactionvalues,interactiontypes)
1272 
1273 end function lattice_interaction_slipbytwin
1274 
1275 
1276 !--------------------------------------------------------------------------------------------------
1279 !--------------------------------------------------------------------------------------------------
1280 function lattice_interaction_slipbytrans(Nslip,Ntrans,interactionValues,structure) result(interactionMatrix)
1282  integer, dimension(:), intent(in) :: nslip, & !< number of active slip systems per family
1283  ntrans
1284  real(preal), dimension(:), intent(in) :: interactionvalues
1285  character(len=*), intent(in) :: structure
1286  real(preal), dimension(sum(Nslip),sum(Ntrans)) :: interactionmatrix
1287 
1288  integer, dimension(:), allocatable :: nslipmax, &
1289  ntransmax
1290  integer, dimension(:,:), allocatable :: interactiontypes
1291 
1292  integer, dimension(FCC_NTRANS,FCC_NSLIP), parameter :: &
1293  fcc_interactionsliptrans = reshape( [&
1294  1,1,1,3,3,3,2,2,2,3,3,3, & ! -----> trans (acting)
1295  1,1,1,3,3,3,3,3,3,2,2,2, & ! |
1296  1,1,1,2,2,2,3,3,3,3,3,3, & ! |
1297  3,3,3,1,1,1,3,3,3,2,2,2, & ! v
1298  3,3,3,1,1,1,2,2,2,3,3,3, & ! slip (reacting)
1299  2,2,2,1,1,1,3,3,3,3,3,3, &
1300  2,2,2,3,3,3,1,1,1,3,3,3, &
1301  3,3,3,2,2,2,1,1,1,3,3,3, &
1302  3,3,3,3,3,3,1,1,1,2,2,2, &
1303  3,3,3,2,2,2,3,3,3,1,1,1, &
1304  2,2,2,3,3,3,3,3,3,1,1,1, &
1305  3,3,3,3,3,3,2,2,2,1,1,1, &
1306 
1307  4,4,4,4,4,4,4,4,4,4,4,4, &
1308  4,4,4,4,4,4,4,4,4,4,4,4, &
1309  4,4,4,4,4,4,4,4,4,4,4,4, &
1310  4,4,4,4,4,4,4,4,4,4,4,4, &
1311  4,4,4,4,4,4,4,4,4,4,4,4, &
1312  4,4,4,4,4,4,4,4,4,4,4,4 &
1313  ],shape(fcc_interactionsliptrans))
1314 
1315  if (len_trim(structure) /= 3) &
1316  call io_error(137,ext_msg='lattice_interaction_SlipByTrans: '//trim(structure))
1317 
1318  select case(structure)
1319  case('fcc')
1320  interactiontypes = fcc_interactionsliptrans
1321  nslipmax = fcc_nslipsystem
1322  ntransmax = fcc_ntranssystem
1323  case default
1324  call io_error(137,ext_msg='lattice_interaction_SlipByTrans: '//trim(structure))
1325  end select
1326 
1327  interactionmatrix = buildinteraction(nslip,ntrans,nslipmax,ntransmax,interactionvalues,interactiontypes)
1328 
1329  end function lattice_interaction_slipbytrans
1330 
1331 
1332 !--------------------------------------------------------------------------------------------------
1335 !--------------------------------------------------------------------------------------------------
1336 function lattice_interaction_twinbyslip(Ntwin,Nslip,interactionValues,structure) result(interactionMatrix)
1338  integer, dimension(:), intent(in) :: ntwin, & !< number of active twin systems per family
1339  nslip
1340  real(preal), dimension(:), intent(in) :: interactionvalues
1341  character(len=*), intent(in) :: structure
1342  real(preal), dimension(sum(Ntwin),sum(Nslip)) :: interactionmatrix
1343 
1344  integer, dimension(:), allocatable :: ntwinmax, &
1345  nslipmax
1346  integer, dimension(:,:), allocatable :: interactiontypes
1347 
1348  integer, dimension(FCC_NSLIP,FCC_NTWIN), parameter :: &
1349  fcc_interactiontwinslip = 1
1350 
1351  integer, dimension(BCC_NSLIP,BCC_NTWIN), parameter :: &
1352  bcc_interactiontwinslip = 1
1353 
1354  integer, dimension(HEX_NSLIP,HEX_NTWIN), parameter :: &
1355  hex_interactiontwinslip = reshape( [&
1356  1, 1, 1, 5, 5, 5, 9, 9, 9, 13,13,13,13,13,13, 17,17,17,17,17,17,17,17,17,17,17,17, 21,21,21,21,21,21, & ! ----> slip (acting)
1357  1, 1, 1, 5, 5, 5, 9, 9, 9, 13,13,13,13,13,13, 17,17,17,17,17,17,17,17,17,17,17,17, 21,21,21,21,21,21, & ! |
1358  1, 1, 1, 5, 5, 5, 9, 9, 9, 13,13,13,13,13,13, 17,17,17,17,17,17,17,17,17,17,17,17, 21,21,21,21,21,21, & ! |
1359  1, 1, 1, 5, 5, 5, 9, 9, 9, 13,13,13,13,13,13, 17,17,17,17,17,17,17,17,17,17,17,17, 21,21,21,21,21,21, & ! v
1360  1, 1, 1, 5, 5, 5, 9, 9, 9, 13,13,13,13,13,13, 17,17,17,17,17,17,17,17,17,17,17,17, 21,21,21,21,21,21, & ! twin (reacting)
1361  1, 1, 1, 5, 5, 5, 9, 9, 9, 13,13,13,13,13,13, 17,17,17,17,17,17,17,17,17,17,17,17, 21,21,21,21,21,21, &
1362  !
1363  2, 2, 2, 6, 6, 6, 10,10,10, 14,14,14,14,14,14, 18,18,18,18,18,18,18,18,18,18,18,18, 22,22,22,22,22,22, &
1364  2, 2, 2, 6, 6, 6, 10,10,10, 14,14,14,14,14,14, 18,18,18,18,18,18,18,18,18,18,18,18, 22,22,22,22,22,22, &
1365  2, 2, 2, 6, 6, 6, 10,10,10, 14,14,14,14,14,14, 18,18,18,18,18,18,18,18,18,18,18,18, 22,22,22,22,22,22, &
1366  2, 2, 2, 6, 6, 6, 10,10,10, 14,14,14,14,14,14, 18,18,18,18,18,18,18,18,18,18,18,18, 22,22,22,22,22,22, &
1367  2, 2, 2, 6, 6, 6, 10,10,10, 14,14,14,14,14,14, 18,18,18,18,18,18,18,18,18,18,18,18, 22,22,22,22,22,22, &
1368  2, 2, 2, 6, 6, 6, 10,10,10, 14,14,14,14,14,14, 18,18,18,18,18,18,18,18,18,18,18,18, 22,22,22,22,22,22, &
1369  !
1370  3, 3, 3, 7, 7, 7, 11,11,11, 15,15,15,15,15,15, 19,19,19,19,19,19,19,19,19,19,19,19, 23,23,23,23,23,23, &
1371  3, 3, 3, 7, 7, 7, 11,11,11, 15,15,15,15,15,15, 19,19,19,19,19,19,19,19,19,19,19,19, 23,23,23,23,23,23, &
1372  3, 3, 3, 7, 7, 7, 11,11,11, 15,15,15,15,15,15, 19,19,19,19,19,19,19,19,19,19,19,19, 23,23,23,23,23,23, &
1373  3, 3, 3, 7, 7, 7, 11,11,11, 15,15,15,15,15,15, 19,19,19,19,19,19,19,19,19,19,19,19, 23,23,23,23,23,23, &
1374  3, 3, 3, 7, 7, 7, 11,11,11, 15,15,15,15,15,15, 19,19,19,19,19,19,19,19,19,19,19,19, 23,23,23,23,23,23, &
1375  3, 3, 3, 7, 7, 7, 11,11,11, 15,15,15,15,15,15, 19,19,19,19,19,19,19,19,19,19,19,19, 23,23,23,23,23,23, &
1376  !
1377  4, 4, 4, 8, 8, 8, 12,12,12, 16,16,16,16,16,16, 20,20,20,20,20,20,20,20,20,20,20,20, 24,24,24,24,24,24, &
1378  4, 4, 4, 8, 8, 8, 12,12,12, 16,16,16,16,16,16, 20,20,20,20,20,20,20,20,20,20,20,20, 24,24,24,24,24,24, &
1379  4, 4, 4, 8, 8, 8, 12,12,12, 16,16,16,16,16,16, 20,20,20,20,20,20,20,20,20,20,20,20, 24,24,24,24,24,24, &
1380  4, 4, 4, 8, 8, 8, 12,12,12, 16,16,16,16,16,16, 20,20,20,20,20,20,20,20,20,20,20,20, 24,24,24,24,24,24, &
1381  4, 4, 4, 8, 8, 8, 12,12,12, 16,16,16,16,16,16, 20,20,20,20,20,20,20,20,20,20,20,20, 24,24,24,24,24,24, &
1382  4, 4, 4, 8, 8, 8, 12,12,12, 16,16,16,16,16,16, 20,20,20,20,20,20,20,20,20,20,20,20, 24,24,24,24,24,24 &
1383  ],shape(hex_interactiontwinslip))
1384 
1385  if (len_trim(structure) /= 3) &
1386  call io_error(137,ext_msg='lattice_interaction_TwinBySlip: '//trim(structure))
1387 
1388  select case(structure)
1389  case('fcc')
1390  interactiontypes = fcc_interactiontwinslip
1391  ntwinmax = fcc_ntwinsystem
1392  nslipmax = fcc_nslipsystem
1393  case('bcc')
1394  interactiontypes = bcc_interactiontwinslip
1395  ntwinmax = bcc_ntwinsystem
1396  nslipmax = bcc_nslipsystem
1397  case('hex')
1398  interactiontypes = hex_interactiontwinslip
1399  ntwinmax = hex_ntwinsystem
1400  nslipmax = hex_nslipsystem
1401  case default
1402  call io_error(137,ext_msg='lattice_interaction_TwinBySlip: '//trim(structure))
1403  end select
1404 
1405  interactionmatrix = buildinteraction(ntwin,nslip,ntwinmax,nslipmax,interactionvalues,interactiontypes)
1406 
1407 end function lattice_interaction_twinbyslip
1408 
1409 
1410 !--------------------------------------------------------------------------------------------------
1413 !--------------------------------------------------------------------------------------------------
1414 function lattice_schmidmatrix_slip(Nslip,structure,cOverA) result(SchmidMatrix)
1416  integer, dimension(:), intent(in) :: nslip
1417  character(len=*), intent(in) :: structure
1418  real(preal), intent(in) :: covera
1419  real(preal), dimension(3,3,sum(Nslip)) :: schmidmatrix
1420 
1421  real(preal), dimension(3,3,sum(Nslip)) :: coordinatesystem
1422  real(preal), dimension(:,:), allocatable :: slipsystems
1423  integer, dimension(:), allocatable :: nslipmax
1424  integer :: i
1425 
1426  if (len_trim(structure) /= 3) &
1427  call io_error(137,ext_msg='lattice_SchmidMatrix_slip: '//trim(structure))
1428 
1429  select case(structure)
1430  case('fcc')
1431  nslipmax = fcc_nslipsystem
1432  slipsystems = fcc_systemslip
1433  case('bcc')
1434  nslipmax = bcc_nslipsystem
1435  slipsystems = bcc_systemslip
1436  case('hex')
1437  nslipmax = hex_nslipsystem
1438  slipsystems = hex_systemslip
1439  case('bct')
1440  nslipmax = bct_nslipsystem
1441  slipsystems = bct_systemslip
1442  case default
1443  call io_error(137,ext_msg='lattice_SchmidMatrix_slip: '//trim(structure))
1444  end select
1445 
1446  if (any(nslipmax(1:size(nslip)) - nslip < 0)) &
1447  call io_error(145,ext_msg='Nslip '//trim(structure))
1448  if (any(nslip < 0)) &
1449  call io_error(144,ext_msg='Nslip '//trim(structure))
1450 
1451  coordinatesystem = buildcoordinatesystem(nslip,nslipmax,slipsystems,structure,covera)
1452 
1453  do i = 1, sum(nslip)
1454  schmidmatrix(1:3,1:3,i) = math_outer(coordinatesystem(1:3,1,i),coordinatesystem(1:3,2,i))
1455  if (abs(math_trace33(schmidmatrix(1:3,1:3,i))) > tol_math_check) &
1456  call io_error(0,i,ext_msg = 'dilatational Schmid matrix for slip')
1457  enddo
1458 
1459 end function lattice_schmidmatrix_slip
1460 
1461 
1462 !--------------------------------------------------------------------------------------------------
1465 !--------------------------------------------------------------------------------------------------
1466 function lattice_schmidmatrix_twin(Ntwin,structure,cOverA) result(SchmidMatrix)
1468  integer, dimension(:), intent(in) :: ntwin
1469  character(len=*), intent(in) :: structure
1470  real(preal), intent(in) :: covera
1471  real(preal), dimension(3,3,sum(Ntwin)) :: schmidmatrix
1472 
1473  real(preal), dimension(3,3,sum(Ntwin)) :: coordinatesystem
1474  real(preal), dimension(:,:), allocatable :: twinsystems
1475  integer, dimension(:), allocatable :: ntwinmax
1476  integer :: i
1477 
1478  if (len_trim(structure) /= 3) &
1479  call io_error(137,ext_msg='lattice_SchmidMatrix_twin: '//trim(structure))
1480 
1481  select case(structure)
1482  case('fcc')
1483  ntwinmax = fcc_ntwinsystem
1484  twinsystems = fcc_systemtwin
1485  case('bcc')
1486  ntwinmax = bcc_ntwinsystem
1487  twinsystems = bcc_systemtwin
1488  case('hex')
1489  ntwinmax = hex_ntwinsystem
1490  twinsystems = hex_systemtwin
1491  case default
1492  call io_error(137,ext_msg='lattice_SchmidMatrix_twin: '//trim(structure))
1493  end select
1494 
1495  if (any(ntwinmax(1:size(ntwin)) - ntwin < 0)) &
1496  call io_error(145,ext_msg='Ntwin '//trim(structure))
1497  if (any(ntwin < 0)) &
1498  call io_error(144,ext_msg='Ntwin '//trim(structure))
1499 
1500  coordinatesystem = buildcoordinatesystem(ntwin,ntwinmax,twinsystems,structure,covera)
1501 
1502  do i = 1, sum(ntwin)
1503  schmidmatrix(1:3,1:3,i) = math_outer(coordinatesystem(1:3,1,i),coordinatesystem(1:3,2,i))
1504  if (abs(math_trace33(schmidmatrix(1:3,1:3,i))) > tol_math_check) &
1505  call io_error(0,i,ext_msg = 'dilatational Schmid matrix for twin')
1506  enddo
1507 
1508 end function lattice_schmidmatrix_twin
1509 
1510 
1511 !--------------------------------------------------------------------------------------------------
1514 !--------------------------------------------------------------------------------------------------
1515 function lattice_schmidmatrix_trans(Ntrans,structure_target,cOverA,a_bcc,a_fcc) result(SchmidMatrix)
1517  integer, dimension(:), intent(in) :: ntrans
1518  character(len=*), intent(in) :: structure_target
1519  real(preal), intent(in) :: covera
1520  real(preal), dimension(3,3,sum(Ntrans)) :: schmidmatrix
1521 
1522  real(preal), dimension(3,3,sum(Ntrans)) :: devnull
1523  real(preal) :: a_bcc, a_fcc
1524 
1525  if (len_trim(structure_target) /= 3) &
1526  call io_error(137,ext_msg='lattice_SchmidMatrix_trans: '//trim(structure_target))
1527  if (structure_target(1:3) /= 'bcc' .and. structure_target(1:3) /= 'hex') &
1528  call io_error(137,ext_msg='lattice_SchmidMatrix_trans: '//trim(structure_target))
1529 
1530  if (structure_target(1:3) == 'hex' .and. (covera < 1.0_preal .or. covera > 2.0_preal)) &
1531  call io_error(131,ext_msg='lattice_SchmidMatrix_trans: '//trim(structure_target))
1532 
1533  if (structure_target(1:3) == 'bcc' .and. (a_bcc <= 0.0_preal .or. a_fcc <= 0.0_preal)) &
1534  call io_error(134,ext_msg='lattice_SchmidMatrix_trans: '//trim(structure_target))
1535 
1536  call buildtransformationsystem(devnull,schmidmatrix,ntrans,covera,a_fcc,a_bcc)
1537 
1538 end function lattice_schmidmatrix_trans
1539 
1540 
1541 !--------------------------------------------------------------------------------------------------
1544 !--------------------------------------------------------------------------------------------------
1545 function lattice_schmidmatrix_cleavage(Ncleavage,structure,cOverA) result(SchmidMatrix)
1547  integer, dimension(:), intent(in) :: ncleavage
1548  character(len=*), intent(in) :: structure
1549  real(preal), intent(in) :: covera
1550  real(preal), dimension(3,3,3,sum(Ncleavage)) :: schmidmatrix
1551 
1552  real(preal), dimension(3,3,sum(Ncleavage)) :: coordinatesystem
1553  real(preal), dimension(:,:), allocatable :: cleavagesystems
1554  integer, dimension(:), allocatable :: ncleavagemax
1555  integer :: i
1556 
1557  if (len_trim(structure) /= 3) &
1558  call io_error(137,ext_msg='lattice_SchmidMatrix_cleavage: '//trim(structure))
1559 
1560  select case(structure)
1561  case('ort')
1562  ncleavagemax = ort_ncleavagesystem
1563  cleavagesystems = ort_systemcleavage
1564  case('fcc')
1565  ncleavagemax = fcc_ncleavagesystem
1566  cleavagesystems = fcc_systemcleavage
1567  case('bcc')
1568  ncleavagemax = bcc_ncleavagesystem
1569  cleavagesystems = bcc_systemcleavage
1570  case default
1571  call io_error(137,ext_msg='lattice_SchmidMatrix_cleavage: '//trim(structure))
1572  end select
1573 
1574  if (any(ncleavagemax(1:size(ncleavage)) - ncleavage < 0)) &
1575  call io_error(145,ext_msg='Ncleavage '//trim(structure))
1576  if (any(ncleavage < 0)) &
1577  call io_error(144,ext_msg='Ncleavage '//trim(structure))
1578 
1579  coordinatesystem = buildcoordinatesystem(ncleavage,ncleavagemax,cleavagesystems,structure,covera)
1580 
1581  do i = 1, sum(ncleavage)
1582  schmidmatrix(1:3,1:3,1,i) = math_outer(coordinatesystem(1:3,1,i),coordinatesystem(1:3,2,i))
1583  schmidmatrix(1:3,1:3,2,i) = math_outer(coordinatesystem(1:3,3,i),coordinatesystem(1:3,2,i))
1584  schmidmatrix(1:3,1:3,3,i) = math_outer(coordinatesystem(1:3,2,i),coordinatesystem(1:3,2,i))
1585  enddo
1586 
1587 end function lattice_schmidmatrix_cleavage
1588 
1589 
1590 !--------------------------------------------------------------------------------------------------
1592 !--------------------------------------------------------------------------------------------------
1593 function lattice_slip_direction(Nslip,structure,cOverA) result(d)
1595  integer, dimension(:), intent(in) :: nslip
1596  character(len=*), intent(in) :: structure
1597  real(preal), intent(in) :: covera
1598  real(preal), dimension(3,sum(Nslip)) :: d
1599 
1600  real(preal), dimension(3,3,sum(Nslip)) :: coordinatesystem
1601 
1602  coordinatesystem = coordinatesystem_slip(nslip,structure,covera)
1603  d = coordinatesystem(1:3,1,1:sum(nslip))
1604 
1605 end function lattice_slip_direction
1606 
1607 
1608 !--------------------------------------------------------------------------------------------------
1610 !--------------------------------------------------------------------------------------------------
1611 function lattice_slip_normal(Nslip,structure,cOverA) result(n)
1613  integer, dimension(:), intent(in) :: nslip
1614  character(len=*), intent(in) :: structure
1615  real(preal), intent(in) :: covera
1616  real(preal), dimension(3,sum(Nslip)) :: n
1617 
1618  real(preal), dimension(3,3,sum(Nslip)) :: coordinatesystem
1619 
1620  coordinatesystem = coordinatesystem_slip(nslip,structure,covera)
1621  n = coordinatesystem(1:3,2,1:sum(nslip))
1622 
1623 end function lattice_slip_normal
1624 
1625 
1626 !--------------------------------------------------------------------------------------------------
1628 !--------------------------------------------------------------------------------------------------
1629 function lattice_slip_transverse(Nslip,structure,cOverA) result(t)
1631  integer, dimension(:), intent(in) :: nslip
1632  character(len=*), intent(in) :: structure
1633  real(preal), intent(in) :: covera
1634  real(preal), dimension(3,sum(Nslip)) :: t
1635 
1636  real(preal), dimension(3,3,sum(Nslip)) :: coordinatesystem
1637 
1638  coordinatesystem = coordinatesystem_slip(nslip,structure,covera)
1639  t = coordinatesystem(1:3,3,1:sum(nslip))
1640 
1641 end function lattice_slip_transverse
1642 
1643 
1644 !--------------------------------------------------------------------------------------------------
1647 !--------------------------------------------------------------------------------------------------
1648 function lattice_labels_slip(Nslip,structure) result(labels)
1650  integer, dimension(:), intent(in) :: nslip
1651  character(len=*), intent(in) :: structure
1652 
1653  character(len=:), dimension(:), allocatable :: labels
1654 
1655  real(preal), dimension(:,:), allocatable :: slipsystems
1656  integer, dimension(:), allocatable :: nslipmax
1657 
1658  if (len_trim(structure) /= 3) &
1659  call io_error(137,ext_msg='lattice_labels_slip: '//trim(structure))
1660 
1661  select case(structure)
1662  case('fcc')
1663  nslipmax = fcc_nslipsystem
1664  slipsystems = fcc_systemslip
1665  case('bcc')
1666  nslipmax = bcc_nslipsystem
1667  slipsystems = bcc_systemslip
1668  case('hex')
1669  nslipmax = hex_nslipsystem
1670  slipsystems = hex_systemslip
1671  case('bct')
1672  nslipmax = bct_nslipsystem
1673  slipsystems = bct_systemslip
1674  case default
1675  call io_error(137,ext_msg='lattice_labels_slip: '//trim(structure))
1676  end select
1677 
1678  if (any(nslipmax(1:size(nslip)) - nslip < 0)) &
1679  call io_error(145,ext_msg='Nslip '//trim(structure))
1680  if (any(nslip < 0)) &
1681  call io_error(144,ext_msg='Nslip '//trim(structure))
1682 
1683  labels = getlabels(nslip,nslipmax,slipsystems)
1684 
1685 end function lattice_labels_slip
1686 
1687 
1688 !--------------------------------------------------------------------------------------------------
1690 !--------------------------------------------------------------------------------------------------
1691 function lattice_applylatticesymmetry33(T,structure) result(T_sym)
1693  real(preal), dimension(3,3) :: t_sym
1694 
1695  real(preal), dimension(3,3), intent(in) :: t
1696  character(len=*), intent(in) :: structure
1697 
1698  integer :: k
1699 
1700  t_sym = 0.0_preal
1701 
1702  if (len_trim(structure) /= 3) &
1703  call io_error(137,ext_msg='lattice_applyLatticeSymmetry33: '//trim(structure))
1704 
1705  select case(structure)
1706  case('iso','fcc','bcc')
1707  do k=1,3
1708  t_sym(k,k) = t(1,1)
1709  enddo
1710  case('hex')
1711  t_sym(1,1) = t(1,1)
1712  t_sym(2,2) = t(1,1)
1713  t_sym(3,3) = t(3,3)
1714  case('ort','bct')
1715  t_sym(1,1) = t(1,1)
1716  t_sym(2,2) = t(2,2)
1717  t_sym(3,3) = t(3,3)
1718  case default
1719  call io_error(137,ext_msg='lattice_applyLatticeSymmetry33: '//trim(structure))
1720  end select
1721 
1722 end function lattice_applylatticesymmetry33
1723 
1724 
1725 !--------------------------------------------------------------------------------------------------
1728 !--------------------------------------------------------------------------------------------------
1729 function applylatticesymmetryc66(C66,structure) result(C66_sym)
1731  real(preal), dimension(6,6) :: c66_sym
1732 
1733  real(preal), dimension(6,6), intent(in) :: c66
1734  character(len=*), intent(in) :: structure
1735 
1736  integer :: j,k
1737 
1738  c66_sym = 0.0_preal
1739 
1740  if (len_trim(structure) /= 3) &
1741  call io_error(137,ext_msg='applyLatticeSymmetryC66: '//trim(structure))
1742 
1743  select case(structure)
1744  case ('iso')
1745  do k=1,3
1746  do j=1,3
1747  c66_sym(k,j) = c66(1,2)
1748  enddo
1749  c66_sym(k,k) = c66(1,1)
1750  c66_sym(k+3,k+3) = 0.5_preal*(c66(1,1)-c66(1,2))
1751  enddo
1752  case ('fcc','bcc')
1753  do k=1,3
1754  do j=1,3
1755  c66_sym(k,j) = c66(1,2)
1756  enddo
1757  c66_sym(k,k) = c66(1,1)
1758  c66_sym(k+3,k+3) = c66(4,4)
1759  enddo
1760  case ('hex')
1761  c66_sym(1,1) = c66(1,1)
1762  c66_sym(2,2) = c66(1,1)
1763  c66_sym(3,3) = c66(3,3)
1764  c66_sym(1,2) = c66(1,2)
1765  c66_sym(2,1) = c66(1,2)
1766  c66_sym(1,3) = c66(1,3)
1767  c66_sym(3,1) = c66(1,3)
1768  c66_sym(2,3) = c66(1,3)
1769  c66_sym(3,2) = c66(1,3)
1770  c66_sym(4,4) = c66(4,4)
1771  c66_sym(5,5) = c66(4,4)
1772  c66_sym(6,6) = 0.5_preal*(c66(1,1)-c66(1,2))
1773  case ('ort')
1774  c66_sym(1,1) = c66(1,1)
1775  c66_sym(2,2) = c66(2,2)
1776  c66_sym(3,3) = c66(3,3)
1777  c66_sym(1,2) = c66(1,2)
1778  c66_sym(2,1) = c66(1,2)
1779  c66_sym(1,3) = c66(1,3)
1780  c66_sym(3,1) = c66(1,3)
1781  c66_sym(2,3) = c66(2,3)
1782  c66_sym(3,2) = c66(2,3)
1783  c66_sym(4,4) = c66(4,4)
1784  c66_sym(5,5) = c66(5,5)
1785  c66_sym(6,6) = c66(6,6)
1786  case ('bct')
1787  c66_sym(1,1) = c66(1,1)
1788  c66_sym(2,2) = c66(1,1)
1789  c66_sym(3,3) = c66(3,3)
1790  c66_sym(1,2) = c66(1,2)
1791  c66_sym(2,1) = c66(1,2)
1792  c66_sym(1,3) = c66(1,3)
1793  c66_sym(3,1) = c66(1,3)
1794  c66_sym(2,3) = c66(1,3)
1795  c66_sym(3,2) = c66(1,3)
1796  c66_sym(4,4) = c66(4,4)
1797  c66_sym(5,5) = c66(4,4)
1798  c66_sym(6,6) = c66(6,6)
1799  case default
1800  call io_error(137,ext_msg='applyLatticeSymmetryC66: '//trim(structure))
1801  end select
1802 
1803 end function applylatticesymmetryc66
1804 
1805 
1806 !--------------------------------------------------------------------------------------------------
1809 !--------------------------------------------------------------------------------------------------
1810 function lattice_labels_twin(Ntwin,structure) result(labels)
1812  integer, dimension(:), intent(in) :: ntwin
1813  character(len=*), intent(in) :: structure
1814 
1815  character(len=:), dimension(:), allocatable :: labels
1816 
1817  real(preal), dimension(:,:), allocatable :: twinsystems
1818  integer, dimension(:), allocatable :: ntwinmax
1819 
1820  if (len_trim(structure) /= 3) &
1821  call io_error(137,ext_msg='lattice_labels_twin: '//trim(structure))
1822 
1823  select case(structure)
1824  case('fcc')
1825  ntwinmax = fcc_ntwinsystem
1826  twinsystems = fcc_systemtwin
1827  case('bcc')
1828  ntwinmax = bcc_ntwinsystem
1829  twinsystems = bcc_systemtwin
1830  case('hex')
1831  ntwinmax = hex_ntwinsystem
1832  twinsystems = hex_systemtwin
1833  case default
1834  call io_error(137,ext_msg='lattice_labels_twin: '//trim(structure))
1835  end select
1836 
1837  if (any(ntwinmax(1:size(ntwin)) - ntwin < 0)) &
1838  call io_error(145,ext_msg='Ntwin '//trim(structure))
1839  if (any(ntwin < 0)) &
1840  call io_error(144,ext_msg='Ntwin '//trim(structure))
1841 
1842  labels = getlabels(ntwin,ntwinmax,twinsystems)
1843 
1844 end function lattice_labels_twin
1845 
1846 
1847 !--------------------------------------------------------------------------------------------------
1850 !--------------------------------------------------------------------------------------------------
1851 function slipprojection_transverse(Nslip,structure,cOverA) result(projection)
1853  integer, dimension(:), intent(in) :: nslip
1854  character(len=*), intent(in) :: structure
1855  real(preal), intent(in) :: covera
1856  real(preal), dimension(sum(Nslip),sum(Nslip)) :: projection
1857 
1858  real(preal), dimension(3,sum(Nslip)) :: n, t
1859  integer :: i, j
1860 
1861  n = lattice_slip_normal(nslip,structure,covera)
1862  t = lattice_slip_transverse(nslip,structure,covera)
1863 
1864  do i=1, sum(nslip); do j=1, sum(nslip)
1865  projection(i,j) = abs(math_inner(n(:,i),t(:,j)))
1866  enddo; enddo
1867 
1868 end function slipprojection_transverse
1869 
1870 
1871 !--------------------------------------------------------------------------------------------------
1874 !--------------------------------------------------------------------------------------------------
1875 function slipprojection_direction(Nslip,structure,cOverA) result(projection)
1877  integer, dimension(:), intent(in) :: nslip
1878  character(len=*), intent(in) :: structure
1879  real(preal), intent(in) :: covera
1880  real(preal), dimension(sum(Nslip),sum(Nslip)) :: projection
1881 
1882  real(preal), dimension(3,sum(Nslip)) :: n, d
1883  integer :: i, j
1884 
1885  n = lattice_slip_normal(nslip,structure,covera)
1886  d = lattice_slip_direction(nslip,structure,covera)
1887 
1888  do i=1, sum(nslip); do j=1, sum(nslip)
1889  projection(i,j) = abs(math_inner(n(:,i),d(:,j)))
1890  enddo; enddo
1891 
1892 end function slipprojection_direction
1893 
1894 
1895 !--------------------------------------------------------------------------------------------------
1898 !--------------------------------------------------------------------------------------------------
1899 function coordinatesystem_slip(Nslip,structure,cOverA) result(coordinateSystem)
1901  integer, dimension(:), intent(in) :: nslip
1902  character(len=*), intent(in) :: structure
1903  real(preal), intent(in) :: covera
1904  real(preal), dimension(3,3,sum(Nslip)) :: coordinatesystem
1905 
1906  real(preal), dimension(:,:), allocatable :: slipsystems
1907  integer, dimension(:), allocatable :: nslipmax
1908 
1909  if (len_trim(structure) /= 3) &
1910  call io_error(137,ext_msg='coordinateSystem_slip: '//trim(structure))
1911 
1912  select case(structure)
1913  case('fcc')
1914  nslipmax = fcc_nslipsystem
1915  slipsystems = fcc_systemslip
1916  case('bcc')
1917  nslipmax = bcc_nslipsystem
1918  slipsystems = bcc_systemslip
1919  case('hex')
1920  nslipmax = hex_nslipsystem
1921  slipsystems = hex_systemslip
1922  case('bct')
1923  nslipmax = bct_nslipsystem
1924  slipsystems = bct_systemslip
1925  case default
1926  call io_error(137,ext_msg='coordinateSystem_slip: '//trim(structure))
1927  end select
1928 
1929  if (any(nslipmax(1:size(nslip)) - nslip < 0)) &
1930  call io_error(145,ext_msg='Nslip '//trim(structure))
1931  if (any(nslip < 0)) &
1932  call io_error(144,ext_msg='Nslip '//trim(structure))
1933 
1934  coordinatesystem = buildcoordinatesystem(nslip,nslipmax,slipsystems,structure,covera)
1935 
1936 end function coordinatesystem_slip
1937 
1938 
1939 !--------------------------------------------------------------------------------------------------
1941 !--------------------------------------------------------------------------------------------------
1942 function buildinteraction(reacting_used,acting_used,reacting_max,acting_max,values,matrix)
1944  integer, dimension(:), intent(in) :: &
1945  reacting_used, & !< # of reacting systems per family as specified in material.config
1946  acting_used, & !< # of acting systems per family as specified in material.config
1947  reacting_max, & !< max # of reacting systems per family for given lattice
1948  acting_max
1949  real(preal), dimension(:), intent(in) :: values
1950  integer, dimension(:,:), intent(in) :: matrix
1951  real(preal), dimension(sum(reacting_used),sum(acting_used)) :: buildinteraction
1952 
1953  integer :: &
1954  acting_family_index, acting_family, acting_system, &
1955  reacting_family_index, reacting_family, reacting_system, &
1956  i,j,k,l
1957 
1958  do acting_family = 1,size(acting_used,1)
1959  acting_family_index = sum(acting_used(1:acting_family-1))
1960  do acting_system = 1,acting_used(acting_family)
1961 
1962  do reacting_family = 1,size(reacting_used,1)
1963  reacting_family_index = sum(reacting_used(1:reacting_family-1))
1964  do reacting_system = 1,reacting_used(reacting_family)
1965 
1966  i = sum( acting_max(1: acting_family-1)) + acting_system
1967  j = sum(reacting_max(1:reacting_family-1)) + reacting_system
1968 
1969  k = acting_family_index + acting_system
1970  l = reacting_family_index + reacting_system
1971 
1972  if (matrix(i,j) > size(values)) call io_error(138,ext_msg='buildInteraction')
1973 
1974  buildinteraction(l,k) = values(matrix(i,j))
1975 
1976  enddo; enddo
1977  enddo; enddo
1978 
1979 end function buildinteraction
1980 
1981 
1982 !--------------------------------------------------------------------------------------------------
1985 !--------------------------------------------------------------------------------------------------
1986 function buildcoordinatesystem(active,potential,system,structure,cOverA)
1988  integer, dimension(:), intent(in) :: &
1989  active, & !< # of active systems per family
1990  potential
1991  real(preal), dimension(:,:), intent(in) :: &
1992  system
1993  character(len=*), intent(in) :: &
1994  structure
1995  real(preal), intent(in) :: &
1996  covera
1997  real(preal), dimension(3,3,sum(active)) :: &
1999 
2000  real(preal), dimension(3) :: &
2001  direction, normal
2002  integer :: &
2003  a, & !< index of active system
2004  p, & !< index in potential system matrix
2005  f, & !< index of my family
2006  s
2007 
2008  if (len_trim(structure) /= 3) &
2009  call io_error(137,ext_msg='buildCoordinateSystem: '//trim(structure))
2010  if (trim(structure) == 'bct' .and. covera > 2.0_preal) &
2011  call io_error(131,ext_msg='buildCoordinateSystem:'//trim(structure))
2012  if (trim(structure) == 'hex' .and. (covera < 1.0_preal .or. covera > 2.0_preal)) &
2013  call io_error(131,ext_msg='buildCoordinateSystem:'//trim(structure))
2014 
2015  a = 0
2016  activefamilies: do f = 1,size(active,1)
2017  activesystems: do s = 1,active(f)
2018  a = a + 1
2019  p = sum(potential(1:f-1))+s
2020 
2021  select case(trim(structure))
2022 
2023  case ('fcc','bcc','iso','ort','bct')
2024  direction = system(1:3,p)
2025  normal = system(4:6,p)
2026 
2027  case ('hex')
2028  direction = [ system(1,p)*1.5_preal, &
2029  (system(1,p)+2.0_preal*system(2,p))*sqrt(0.75_preal), &
2030  system(4,p)*covera ] ! direction [uvtw]->[3u/2 (u+2v)*sqrt(3)/2 w*(p/a)])
2031  normal = [ system(5,p), &
2032  (system(5,p)+2.0_preal*system(6,p))/sqrt(3.0_preal), &
2033  system(8,p)/covera ] ! plane (hkil)->(h (h+2k)/sqrt(3) l/(p/a))
2034 
2035  case default
2036  call io_error(137,ext_msg='buildCoordinateSystem: '//trim(structure))
2037 
2038  end select
2039 
2040  buildcoordinatesystem(1:3,1,a) = direction/norm2(direction)
2041  buildcoordinatesystem(1:3,2,a) = normal /norm2(normal)
2042  buildcoordinatesystem(1:3,3,a) = math_cross(direction/norm2(direction),&
2043  normal /norm2(normal))
2044 
2045  enddo activesystems
2046  enddo activefamilies
2047 
2048 end function buildcoordinatesystem
2049 
2050 
2051 !--------------------------------------------------------------------------------------------------
2053 ! Needed to calculate Schmid matrix and rotated stiffness matrices.
2054 ! @details: set c/a = 0.0 for fcc -> bcc transformation
2055 ! set a_Xcc = 0.0 for fcc -> hex transformation
2056 !--------------------------------------------------------------------------------------------------
2057 subroutine buildtransformationsystem(Q,S,Ntrans,cOverA,a_fcc,a_bcc)
2059  integer, dimension(:), intent(in) :: &
2060  Ntrans
2061  real(pReal), dimension(3,3,sum(Ntrans)), intent(out) :: &
2062  Q, & !< Total rotation: Q = R*B
2063  S
2064  real(pReal), intent(in) :: &
2065  cOverA, & !< c/a for target hex structure
2066  a_bcc, & !< lattice parameter a for target bcc structure
2067  a_fcc
2068 
2069  type(rotation) :: &
2070  R, & !< Pitsch rotation
2071  B
2072  real(pReal), dimension(3,3) :: &
2073  U, & !< Bain deformation
2074  ss, sd
2075  real(pReal), dimension(3) :: &
2076  x, y, z
2077  integer :: &
2078  i
2079  real(pReal), dimension(3+3,FCC_NTRANS), parameter :: &
2080  FCCTOHEX_SYSTEMTRANS = reshape(real( [&
2081  -2, 1, 1, 1, 1, 1, &
2082  1,-2, 1, 1, 1, 1, &
2083  1, 1,-2, 1, 1, 1, &
2084  2,-1, 1, -1,-1, 1, &
2085  -1, 2, 1, -1,-1, 1, &
2086  -1,-1,-2, -1,-1, 1, &
2087  -2,-1,-1, 1,-1,-1, &
2088  1, 2,-1, 1,-1,-1, &
2089  1,-1, 2, 1,-1,-1, &
2090  2, 1,-1, -1, 1,-1, &
2091  -1,-2,-1, -1, 1,-1, &
2092  -1, 1, 2, -1, 1,-1 &
2093  ],preal),shape(fcctohex_systemtrans))
2094  real(pReal), dimension(4,fcc_Ntrans), parameter :: &
2095  FCCTOBCC_SYSTEMTRANS = reshape([&
2096  0.0, 1.0, 0.0, 10.26, & ! Pitsch OR (Ma & Hartmaier 2014, Table 3)
2097  0.0,-1.0, 0.0, 10.26, &
2098  0.0, 0.0, 1.0, 10.26, &
2099  0.0, 0.0,-1.0, 10.26, &
2100  1.0, 0.0, 0.0, 10.26, &
2101  -1.0, 0.0, 0.0, 10.26, &
2102  0.0, 0.0, 1.0, 10.26, &
2103  0.0, 0.0,-1.0, 10.26, &
2104  1.0, 0.0, 0.0, 10.26, &
2105  -1.0, 0.0, 0.0, 10.26, &
2106  0.0, 1.0, 0.0, 10.26, &
2107  0.0,-1.0, 0.0, 10.26 &
2108  ],shape(fcctobcc_systemtrans))
2109 
2110  integer, dimension(9,fcc_Ntrans), parameter :: &
2111  FCCTOBCC_BAINVARIANT = reshape( [&
2112  1, 0, 0, 0, 1, 0, 0, 0, 1, & ! Pitsch OR (Ma & Hartmaier 2014, Table 3)
2113  1, 0, 0, 0, 1, 0, 0, 0, 1, &
2114  1, 0, 0, 0, 1, 0, 0, 0, 1, &
2115  1, 0, 0, 0, 1, 0, 0, 0, 1, &
2116  0, 1, 0, 1, 0, 0, 0, 0, 1, &
2117  0, 1, 0, 1, 0, 0, 0, 0, 1, &
2118  0, 1, 0, 1, 0, 0, 0, 0, 1, &
2119  0, 1, 0, 1, 0, 0, 0, 0, 1, &
2120  0, 0, 1, 1, 0, 0, 0, 1, 0, &
2121  0, 0, 1, 1, 0, 0, 0, 1, 0, &
2122  0, 0, 1, 1, 0, 0, 0, 1, 0, &
2123  0, 0, 1, 1, 0, 0, 0, 1, 0 &
2124  ],shape(fcctobcc_bainvariant))
2125 
2126  real(pReal), dimension(4,fcc_Ntrans), parameter :: &
2127  FCCTOBCC_BAINROT = reshape([&
2128  1.0, 0.0, 0.0, 45.0, & ! Rotate fcc austensite to bain variant
2129  1.0, 0.0, 0.0, 45.0, &
2130  1.0, 0.0, 0.0, 45.0, &
2131  1.0, 0.0, 0.0, 45.0, &
2132  0.0, 1.0, 0.0, 45.0, &
2133  0.0, 1.0, 0.0, 45.0, &
2134  0.0, 1.0, 0.0, 45.0, &
2135  0.0, 1.0, 0.0, 45.0, &
2136  0.0, 0.0, 1.0, 45.0, &
2137  0.0, 0.0, 1.0, 45.0, &
2138  0.0, 0.0, 1.0, 45.0, &
2139  0.0, 0.0, 1.0, 45.0 &
2140  ],shape(fcctobcc_bainrot))
2141 
2142  if (a_bcc > 0.0_preal .and. a_fcc > 0.0_preal .and. deq0(covera)) then ! fcc -> bcc transformation
2143  do i = 1,sum(ntrans)
2144  call r%fromAxisAngle(fcctobcc_systemtrans(:,i),degrees=.true.,p=1)
2145  call b%fromAxisAngle(fcctobcc_bainrot(:,i), degrees=.true.,p=1)
2146  x = real(fcctobcc_bainvariant(1:3,i),preal)
2147  y = real(fcctobcc_bainvariant(4:6,i),preal)
2148  z = real(fcctobcc_bainvariant(7:9,i),preal)
2149 
2150  u = (a_bcc/a_fcc)*math_outer(x,x) &
2151  + (a_bcc/a_fcc)*math_outer(y,y) * sqrt(2.0_preal) &
2152  + (a_bcc/a_fcc)*math_outer(z,z) * sqrt(2.0_preal)
2153  q(1:3,1:3,i) = matmul(r%asMatrix(),b%asMatrix())
2154  s(1:3,1:3,i) = matmul(r%asMatrix(),u) - math_i3
2155  enddo
2156  elseif (covera > 0.0_preal .and. deq0(a_bcc)) then ! fcc -> hex transformation
2157  ss = math_i3
2158  sd = math_i3
2159  ss(1,3) = sqrt(2.0_preal)/4.0_preal
2160  sd(3,3) = covera/sqrt(8.0_preal/3.0_preal)
2161 
2162  do i = 1,sum(ntrans)
2163  x = fcctohex_systemtrans(1:3,i)/norm2(fcctohex_systemtrans(1:3,i))
2164  z = fcctohex_systemtrans(4:6,i)/norm2(fcctohex_systemtrans(4:6,i))
2165  y = -math_cross(x,z)
2166  q(1:3,1,i) = x
2167  q(1:3,2,i) = y
2168  q(1:3,3,i) = z
2169  s(1:3,1:3,i) = matmul(q(1:3,1:3,i), matmul(matmul(sd,ss), transpose(q(1:3,1:3,i)))) - math_i3 ! ToDo: This is of interest for the Schmid matrix only
2170  enddo
2171  else
2172  call io_error(132,ext_msg='buildTransformationSystem')
2173  endif
2174 
2175 end subroutine buildtransformationsystem
2176 
2177 
2178 !--------------------------------------------------------------------------------------------------
2180 !--------------------------------------------------------------------------------------------------
2181 function getlabels(active,potential,system) result(labels)
2183  integer, dimension(:), intent(in) :: &
2184  active, & !< # of active systems per family
2185  potential
2186  real(preal), dimension(:,:), intent(in) :: &
2187  system
2188 
2189  character(len=:), dimension(:), allocatable :: labels
2190  character(len=:), allocatable :: label
2191 
2192  integer :: i,j
2193  integer :: &
2194  a, & !< index of active system
2195  p, & !< index in potential system matrix
2196  f, & !< index of my family
2197  s
2198 
2199  i = 2*size(system,1) + (size(system,1) - 2) + 4 ! 2 letters per index + spaces + brackets
2200  allocate(character(len=i) :: labels(sum(active)), label)
2201 
2202  a = 0
2203  activefamilies: do f = 1,size(active,1)
2204  activesystems: do s = 1,active(f)
2205  a = a + 1
2206  p = sum(potential(1:f-1))+s
2207 
2208  i = 1
2209  label(i:i) = '['
2210  direction: do j = 1, size(system,1)/2
2211  write(label(i+1:i+2),'(I2.1)') int(system(j,p))
2212  label(i+3:i+3) = ' '
2213  i = i + 3
2214  enddo direction
2215  label(i:i) = ']'
2216 
2217  i = i +1
2218  label(i:i) = '('
2219  normal: do j = size(system,1)/2+1, size(system,1)
2220  write(label(i+1:i+2),'(I2.1)') int(system(j,p))
2221  label(i+3:i+3) = ' '
2222  i = i + 3
2223  enddo normal
2224  label(i:i) = ')'
2225 
2226  labels(s) = label
2227 
2228  enddo activesystems
2229  enddo activefamilies
2230 
2231 end function getlabels
2232 
2233 
2234 !--------------------------------------------------------------------------------------------------
2237 !--------------------------------------------------------------------------------------------------
2238 function equivalent_nu(C,assumption) result(nu)
2240  real(preal), dimension(6,6), intent(in) :: c
2241  character(len=*), intent(in) :: assumption
2242 
2243  real(preal) :: k, mu, nu
2244  logical :: error
2245  real(preal), dimension(6,6) :: s
2246 
2247  if (io_lc(assumption) == 'voigt') then
2248  k = (c(1,1)+c(2,2)+c(3,3) +2.0_preal*(c(1,2)+c(2,3)+c(1,3))) &
2249  / 9.0_preal
2250  elseif(io_lc(assumption) == 'reuss') then
2251  call math_invert(s,error,c)
2252  if(error) call io_error(0)
2253  k = 1.0_preal &
2254  / (s(1,1)+s(2,2)+s(3,3) +2.0_preal*(s(1,2)+s(2,3)+s(1,3)))
2255  else
2256  call io_error(0)
2257  k = 0.0_preal
2258  endif
2259 
2260  mu = equivalent_mu(c,assumption)
2261  nu = (1.5_preal*k -mu)/(3.0_preal*k+mu)
2262 
2263 end function equivalent_nu
2264 
2265 
2266 !--------------------------------------------------------------------------------------------------
2269 !--------------------------------------------------------------------------------------------------
2270 function equivalent_mu(C,assumption) result(mu)
2272  real(preal), dimension(6,6), intent(in) :: c
2273  character(len=*), intent(in) :: assumption
2274 
2275  real(preal) :: mu
2276  logical :: error
2277  real(preal), dimension(6,6) :: s
2278 
2279  if (io_lc(assumption) == 'voigt') then
2280  mu = (1.0_preal*(c(1,1)+c(2,2)+c(3,3)) -1.0_preal*(c(1,2)+c(2,3)+c(1,3)) +3.0_preal*(c(4,4)+c(5,5)+c(6,6))) &
2281  / 15.0_preal
2282  elseif(io_lc(assumption) == 'reuss') then
2283  call math_invert(s,error,c)
2284  if(error) call io_error(0)
2285  mu = 15.0_preal &
2286  / (4.0_preal*(s(1,1)+s(2,2)+s(3,3)) -4.0_preal*(s(1,2)+s(2,3)+s(1,3)) +3.0_preal*(s(4,4)+s(5,5)+s(6,6)))
2287  else
2288  call io_error(0)
2289  mu = 0.0_preal
2290  endif
2291 
2292 end function equivalent_mu
2293 
2294 
2295 !--------------------------------------------------------------------------------------------------
2297 !--------------------------------------------------------------------------------------------------
2298 subroutine unittest
2300  real(pReal), dimension(:,:,:), allocatable :: CoSy
2301  real(pReal), dimension(:,:), allocatable :: system
2302 
2303  real(pReal), dimension(6,6) :: C
2304  real(pReal), dimension(2) :: r
2305  real(pReal) :: lambda
2306 
2307  call random_number(r)
2308 
2309  system = reshape([1.0_preal+r(1),0.0_preal,0.0_preal, 0.0_preal,1.0_preal+r(2),0.0_preal],[6,1])
2310  cosy = buildcoordinatesystem([1],[1],system,'fcc',0.0_preal)
2311 
2312  if(any(dneq(cosy(1:3,1:3,1),math_i3))) &
2313  call io_error(0)
2314 
2315  call random_number(c)
2316  c(1,1) = c(1,1) + 1.0_preal
2317  c = applylatticesymmetryc66(c,'iso')
2318  if(dneq(c(6,6),equivalent_mu(c,'voigt'),1.0e-12_preal)) &
2319  call io_error(0,ext_msg='equivalent_mu/voigt')
2320  if(dneq(c(6,6),equivalent_mu(c,'voigt'),1.0e-12_preal)) &
2321  call io_error(0,ext_msg='equivalent_mu/reuss')
2322  lambda = c(1,2)
2323  if(dneq(lambda*0.5_preal/(lambda+equivalent_mu(c,'voigt')),equivalent_nu(c,'voigt'),1.0e-12_preal)) &
2324  call io_error(0,ext_msg='equivalent_nu/voigt')
2325  if(dneq(lambda*0.5_preal/(lambda+equivalent_mu(c,'reuss')),equivalent_nu(c,'reuss'),1.0e-12_preal)) &
2326  call io_error(0,ext_msg='equivalent_nu/reuss')
2327 
2328 end subroutine unittest
2329 
2330 end module lattice
lattice::fcc_ncleavagesystem
integer, dimension(1), parameter fcc_ncleavagesystem
Definition: lattice.f90:34
lattice::buildtransformationsystem
subroutine buildtransformationsystem(Q, S, Ntrans, cOverA, a_fcc, a_bcc)
Helper function to define transformation systems.
Definition: lattice.f90:2058
lattice::lattice_interaction_twinbyslip
real(preal) function, dimension(sum(ntwin), sum(nslip)), public lattice_interaction_twinbyslip(Ntwin, Nslip, interactionValues, structure)
Twin-slip interaction matrix details only active twin and slip systems are considered.
Definition: lattice.f90:1337
lattice::lattice_interaction_slipbytwin
real(preal) function, dimension(sum(nslip), sum(ntwin)), public lattice_interaction_slipbytwin(Nslip, Ntwin, interactionValues, structure)
Slip-twin interaction matrix details only active slip and twin systems are considered.
Definition: lattice.f90:1141
lattice::buildcoordinatesystem
real(preal) function, dimension(3, 3, sum(active)) buildcoordinatesystem(active, potential, system, structure, cOverA)
Build a local coordinate system on slip, twin, trans, cleavage systems.
Definition: lattice.f90:1987
rotations
rotation storage and conversion
Definition: rotations.f90:53
lattice::bct_nslipsystem
integer, dimension(13), parameter bct_nslipsystem
Definition: lattice.f90:287
lattice::lattice_bcc_id
@, public lattice_bcc_id
Definition: lattice.f90:395
prec::tol_math_check
real(preal), parameter tol_math_check
tolerance for internal math self-checks (rotation)
Definition: prec.f90:30
lattice::lattice_damagediffusion
real(preal), dimension(:,:,:), allocatable, public, protected lattice_damagediffusion
Definition: lattice.f90:404
lattice::lattice_c66
real(preal), dimension(:,:,:), allocatable, public, protected lattice_c66
Definition: lattice.f90:404
lattice::lattice_schmidmatrix_twin
real(preal) function, dimension(3, 3, sum(ntwin)), public lattice_schmidmatrix_twin(Ntwin, structure, cOverA)
Schmid matrix for twinning details only active twin systems are considered.
Definition: lattice.f90:1467
lattice::fcc_ntrans
integer, parameter fcc_ntrans
total # of transformation systems for fcc
Definition: lattice.f90:37
math::math_voigt66to3333
pure real(preal) function, dimension(3, 3, 3, 3) math_voigt66to3333(m66)
convert 66 Voigt matrix into symmetric 3x3x3x3 matrix
Definition: math.f90:834
io::io_error
subroutine, public io_error(error_ID, el, ip, g, instance, ext_msg)
write error statements to standard out and terminate the Marc/spectral run with exit #9xxx
Definition: IO.f90:305
rotations::rotation
Definition: rotations.f90:63
lattice::bct_nslip
integer, parameter bct_nslip
total # of slip systems for bct
Definition: lattice.f90:290
lattice::bcc_ncleavagesystem
integer, dimension(1), parameter bcc_ncleavagesystem
Definition: lattice.f90:122
config::config_phase
type(tpartitionedstringlist), dimension(:), allocatable, public, protected config_phase
Definition: config.f90:23
lattice::bct_systemslip
real(preal), dimension(3+3, bct_nslip), parameter bct_systemslip
slip systems for bct sorted by Bieler
Definition: lattice.f90:297
math::math_outer
pure real(preal) function, dimension(size(a, 1), size(b, 1)) math_outer(A, B)
outer product of arbitrary sized vectors (A ⊗ B / i,j)
Definition: math.f90:327
lattice::lattice_specificheat
real(preal), dimension(:), allocatable, public, protected lattice_specificheat
Definition: lattice.f90:399
lattice::hex_ntwin
integer, parameter hex_ntwin
total # of twin systems for hex
Definition: lattice.f90:200
lattice::equivalent_nu
real(preal) function equivalent_nu(C, assumption)
Equivalent Poisson's ratio (ν)
Definition: lattice.f90:2239
lattice::fcc_nslip
integer, parameter fcc_nslip
total # of slip systems for fcc
Definition: lattice.f90:37
lattice::ort_systemcleavage
real(preal), dimension(3+3, ort_ncleavage), parameter ort_systemcleavage
Definition: lattice.f90:379
lattice::bcc_ncleavage
integer, parameter bcc_ncleavage
total # of cleavage systems for bcc
Definition: lattice.f90:125
lattice::lattice_forestprojection_edge
Definition: lattice.f90:412
math::math_sym3333to66
pure real(preal) function, dimension(6, 6) math_sym3333to66(m3333, weighted)
convert symmetric 3x3x3x3 matrix into 6x6 matrix
Definition: math.f90:778
lattice::getlabels
character(len=:) function, dimension(:), allocatable getlabels(active, potential, system)
select active systems as strings
Definition: lattice.f90:2182
lattice::bcc_systemslip
real(preal), dimension(3+3, bcc_nslip), parameter bcc_systemslip
Definition: lattice.f90:136
lattice::lattice_iso_id
@, public lattice_iso_id
Definition: lattice.f90:395
lattice::lattice_massdensity
real(preal), dimension(:), allocatable, public, protected lattice_massdensity
Definition: lattice.f90:399
config
Reads in the material configuration from file.
Definition: config.f90:13
lattice::slipprojection_transverse
real(preal) function, dimension(sum(nslip), sum(nslip)) slipprojection_transverse(Nslip, structure, cOverA)
Projection of the transverse direction onto the slip plane.
Definition: lattice.f90:1852
lattice::lattice_interaction_twinbytwin
real(preal) function, dimension(sum(ntwin), sum(ntwin)), public lattice_interaction_twinbytwin(Ntwin, interactionValues, structure)
Twin-twin interaction matrix details only active twin systems are considered.
Definition: lattice.f90:995
lattice::lattice_nu
real(preal), dimension(:), allocatable, public, protected lattice_nu
Definition: lattice.f90:399
lattice::hex_systemslip
real(preal), dimension(4+4, hex_nslip), parameter hex_systemslip
slip systems for hex, sorted by P. Eisenlohr CCW around starting next to a_1 axis
Definition: lattice.f90:209
lattice::lattice_applylatticesymmetry33
real(preal) function, dimension(3, 3), public lattice_applylatticesymmetry33(T, structure)
Return 3x3 tensor with symmetry according to given crystal structure.
Definition: lattice.f90:1692
lattice::lattice_characteristicshear_twin
real(preal) function, dimension(sum(ntwin)), public lattice_characteristicshear_twin(Ntwin, structure, CoverA)
Characteristic shear for twinning.
Definition: lattice.f90:547
lattice::slipprojection_direction
real(preal) function, dimension(sum(nslip), sum(nslip)) slipprojection_direction(Nslip, structure, cOverA)
Projection of the slip direction onto the slip plane.
Definition: lattice.f90:1876
lattice::bcc_nslipsystem
integer, dimension(2), parameter bcc_nslipsystem
Definition: lattice.f90:116
lattice::lattice_schmidmatrix_trans
real(preal) function, dimension(3, 3, sum(ntrans)), public lattice_schmidmatrix_trans(Ntrans, structure_target, cOverA, a_bcc, a_fcc)
Schmid matrix for twinning details only active twin systems are considered.
Definition: lattice.f90:1516
prec
setting precision for real and int type
Definition: prec.f90:13
lattice::applylatticesymmetryc66
real(preal) function, dimension(6, 6) applylatticesymmetryc66(C66, structure)
Return stiffness matrix in 6x6 notation with symmetry according to given crystal structure.
Definition: lattice.f90:1730
lattice::fcc_systemcleavage
real(preal), dimension(3+3, fcc_ncleavage), parameter fcc_systemcleavage
Definition: lattice.f90:106
lattice::fcc_ncleavage
integer, parameter fcc_ncleavage
total # of cleavage systems for fcc
Definition: lattice.f90:37
lattice::lattice_init
subroutine, public lattice_init
Module initialization.
Definition: lattice.f90:457
lattice::lattice_c66_twin
real(preal) function, dimension(6, 6, sum(ntwin)), public lattice_c66_twin(Ntwin, C66, structure, CoverA)
Rotated elasticity matrices for twinning in 66-vector notation.
Definition: lattice.f90:624
lattice::lattice_slip_direction
real(preal) function, dimension(3, sum(nslip)), public lattice_slip_direction(Nslip, structure, cOverA)
Slip direction of slip systems (|| b)
Definition: lattice.f90:1594
lattice::lattice_fcc_twinnucleationslippair
integer, dimension(2, fcc_ntwin), parameter, public lattice_fcc_twinnucleationslippair
Definition: lattice.f90:90
lattice::lattice_hex_id
@, public lattice_hex_id
Definition: lattice.f90:395
io
input/output functions, partly depending on chosen solver
Definition: IO.f90:12
lattice::lattice_mu
real(preal), dimension(:), allocatable, public, protected lattice_mu
Definition: lattice.f90:399
prec::preal
integer, parameter preal
number with 15 significant digits, up to 1e+-307 (typically 64 bit)
Definition: prec.f90:20
lattice::lattice_interaction_slipbyslip
real(preal) function, dimension(sum(nslip), sum(nslip)), public lattice_interaction_slipbyslip(Nslip, interactionValues, structure)
Slip-slip interaction matrix details only active slip systems are considered.
Definition: lattice.f90:775
lattice::lattice_slip_normal
real(preal) function, dimension(3, sum(nslip)), public lattice_slip_normal(Nslip, structure, cOverA)
Normal direction of slip systems (|| n)
Definition: lattice.f90:1612
lattice::lattice_schmidmatrix_slip
real(preal) function, dimension(3, 3, sum(nslip)), public lattice_schmidmatrix_slip(Nslip, structure, cOverA)
Schmid matrix for slip details only active slip systems are considered.
Definition: lattice.f90:1415
math::math_inner
real(preal) pure function math_inner(A, B)
inner product of arbitrary sized vectors (A · B / i,i)
Definition: math.f90:343
lattice::ort_ncleavagesystem
integer, dimension(3), parameter ort_ncleavagesystem
Definition: lattice.f90:369
prec::dneq
logical elemental pure function dneq(a, b, tol)
inequality comparison for float with double precision
Definition: prec.f90:146
lattice::fcc_systemslip
real(preal), dimension(3+3, fcc_nslip), parameter fcc_systemslip
Slip system <110>{111} directions. Sorted according to Eisenlohr & Hantcherli.
Definition: lattice.f90:50
lattice::coordinatesystem_slip
real(preal) function, dimension(3, 3, sum(nslip)) coordinatesystem_slip(Nslip, structure, cOverA)
build a local coordinate system on slip systems
Definition: lattice.f90:1900
lattice::lattice_fcc_id
@, public lattice_fcc_id
Definition: lattice.f90:395
lattice::lattice_forestprojection_screw
Definition: lattice.f90:416
lattice::lattice_c66_trans
real(preal) function, dimension(6, 6, sum(ntrans)), public lattice_c66_trans(Ntrans, C_parent66, structure_target, cOverA_trans, a_bcc, a_fcc)
Rotated elasticity matrices for transformation in 66-vector notation.
Definition: lattice.f90:665
lattice::lattice_bct_id
@, public lattice_bct_id
Definition: lattice.f90:395
lattice::lattice_labels_twin
character(len=:) function, dimension(:), allocatable, public lattice_labels_twin(Ntwin, structure)
Labels for twin systems details only active twin systems are considered.
Definition: lattice.f90:1811
lattice::lattice_slip_transverse
real(preal) function, dimension(3, sum(nslip)), public lattice_slip_transverse(Nslip, structure, cOverA)
Transverse direction of slip systems ( || t = b x n)
Definition: lattice.f90:1630
lattice::bcc_ntwinsystem
integer, dimension(1), parameter bcc_ntwinsystem
Definition: lattice.f90:119
lattice::lattice_labels_slip
character(len=:) function, dimension(:), allocatable, public lattice_labels_slip(Nslip, structure)
Labels for slip systems details only active slip systems are considered.
Definition: lattice.f90:1649
lattice
contains lattice structure definitions including Schmid matrices for slip, twin, trans,
Definition: lattice.f90:13
lattice::equivalent_mu
real(preal) function equivalent_mu(C, assumption)
Equivalent shear modulus (μ)
Definition: lattice.f90:2271
math::math_invert
subroutine math_invert(InvA, error, A)
invert quadratic matrix of arbitrary dimension
Definition: math.f90:521
lattice::hex_ntwinsystem
integer, dimension(4), parameter hex_ntwinsystem
Definition: lattice.f90:197
lattice::lattice_schmidmatrix_cleavage
real(preal) function, dimension(3, 3, 3, sum(ncleavage)), public lattice_schmidmatrix_cleavage(Ncleavage, structure, cOverA)
Schmid matrix for cleavage details only active cleavage systems are considered.
Definition: lattice.f90:1546
math::math_trace33
real(preal) pure function math_trace33(m)
trace of a 3x3 matrix
Definition: math.f90:611
lattice::ort_ncleavage
integer, parameter ort_ncleavage
total # of cleavage systems for ortho
Definition: lattice.f90:372
lattice::hex_systemtwin
real(preal), dimension(4+4, hex_ntwin), parameter hex_systemtwin
twin systems for hex, sorted by P. Eisenlohr CCW around starting next to a_1 axis
Definition: lattice.f90:253
lattice::fcc_systemtwin
real(preal), dimension(3+3, fcc_ntwin), parameter fcc_systemtwin
Twin system <112>{111} directions. Sorted according to Eisenlohr & Hantcherli.
Definition: lattice.f90:74
math
Mathematical library, including random number generation and tensor representations.
Definition: math.f90:12
lattice::lattice_interaction_transbytrans
real(preal) function, dimension(sum(ntrans), sum(ntrans)), public lattice_interaction_transbytrans(Ntrans, interactionValues, structure)
Trans-trans interaction matrix details only active trans systems are considered.
Definition: lattice.f90:1096
lattice::fcc_ntwinsystem
integer, dimension(1), parameter fcc_ntwinsystem
Definition: lattice.f90:28
lattice::bcc_systemcleavage
real(preal), dimension(3+3, bcc_ncleavage), parameter bcc_systemcleavage
Definition: lattice.f90:184
lattice::lattice_interaction_slipbytrans
real(preal) function, dimension(sum(nslip), sum(ntrans)), public lattice_interaction_slipbytrans(Nslip, Ntrans, interactionValues, structure)
Slip-trans interaction matrix details only active slip and trans systems are considered.
Definition: lattice.f90:1281
lattice::lattice_undefined_id
@ lattice_undefined_id
Definition: lattice.f90:395
lattice::bcc_ntwin
integer, parameter bcc_ntwin
total # of twin systems for bcc
Definition: lattice.f90:125
io::unittest
subroutine unittest
check correctness of some IO functions
Definition: IO.f90:662
lattice::lattice_thermalconductivity
real(preal), dimension(:,:,:), allocatable, public, protected lattice_thermalconductivity
Definition: lattice.f90:404
lattice::buildinteraction
real(preal) function, dimension(sum(reacting_used), sum(acting_used)) buildinteraction(reacting_used, acting_used, reacting_max, acting_max, values, matrix)
Populate reduced interaction matrix.
Definition: lattice.f90:1943
io::io_lc
pure character(len=len(string)) function, public io_lc(string)
changes characters in string to lower case
Definition: IO.f90:280
lattice::bcc_nslip
integer, parameter bcc_nslip
total # of slip systems for bcc
Definition: lattice.f90:125
math::pi
real(preal), parameter pi
ratio of a circle's circumference to its diameter
Definition: math.f90:27
lattice::lattice_damagemobility
real(preal), dimension(:), allocatable, public, protected lattice_damagemobility
Definition: lattice.f90:399
math::math_cross
pure real(preal) function, dimension(3) math_cross(A, B)
cross product a x b
Definition: math.f90:312
lattice::hex_nslip
integer, parameter hex_nslip
total # of slip systems for hex
Definition: lattice.f90:200
lattice::hex_nslipsystem
integer, dimension(6), parameter hex_nslipsystem
Definition: lattice.f90:194
lattice::lattice_structure
integer(kind(lattice_undefined_id)), dimension(:), allocatable, public, protected lattice_structure
Definition: lattice.f90:408
lattice::fcc_nslipsystem
integer, dimension(2), parameter fcc_nslipsystem
Definition: lattice.f90:25
lattice::fcc_ntranssystem
integer, dimension(1), parameter fcc_ntranssystem
Definition: lattice.f90:31
lattice::lattice_ort_id
@, public lattice_ort_id
Definition: lattice.f90:395
prec::deq0
logical elemental pure function deq0(a, tol)
equality to 0 comparison for float with double precision
Definition: prec.f90:166
lattice::lattice_nonschmidmatrix
real(preal) function, dimension(1:3, 1:3, sum(nslip)), public lattice_nonschmidmatrix(Nslip, nonSchmidCoefficients, sense)
Non-schmid projections for bcc with up to 6 coefficients.
Definition: lattice.f90:728
lattice::fcc_ntwin
integer, parameter fcc_ntwin
total # of twin systems for fcc
Definition: lattice.f90:37
lattice::bcc_systemtwin
real(preal), dimension(3+3, bcc_ntwin), parameter bcc_systemtwin
Definition: lattice.f90:167
math::math_i3
real(preal), dimension(3, 3), parameter math_i3
3x3 Identity
Definition: math.f90:32