FrontISTR  5.7.1
Large-scale structural analysis program with finit element method
fstr_ctrl_material.f90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
2 ! Copyright (c) 2019 FrontISTR Commons
3 ! This software is released under the MIT License, see LICENSE.txt
4 !-------------------------------------------------------------------------------
7  use hecmw
8  use mmaterial
9  use m_table
11  implicit none
12 
13  private :: read_user_matl
14 
15 contains
16 
17 
18  !----------------------------------------------------------------------
20  integer function fstr_ctrl_get_material( ctrl, matname )
21  integer(kind=kint), intent(in) :: ctrl
22  character(len=*), intent(out) :: matname
23 
24  matname=""
25  fstr_ctrl_get_material = fstr_ctrl_get_param_ex( ctrl, 'NAME ', '# ', 0, 'S', matname )
26  end function fstr_ctrl_get_material
27 
28  !----------------------------------------------------------------------
30  integer function fstr_ctrl_get_usermaterial( ctrl, mattype, nlgeom, nstatus, matval )
31  integer(kind=kint), intent(in) :: ctrl
32  integer(kind=kint), intent(inout) :: mattype
33  integer(kind=kint), intent(out) :: nlgeom
34  integer(kind=kint), intent(out) :: nstatus
35  real(kind=kreal),intent(out) :: matval(:)
36 
37  integer(kind=kint) :: ipt
38  character(len=HECMW_NAME_LEN) :: data_fmt
39  character(len=256) :: s, fname
40 
42  mattype = usermaterial
43  nlgeom = updatelag !default value
44  nstatus = 1
45  if( fstr_ctrl_get_param_ex( ctrl, 'NSTATUS ', '# ', 0, 'I', nstatus )/= 0) return
46  if( fstr_ctrl_get_param_ex( ctrl, 'KIRCHHOFF ', '# ', 0, 'E', ipt )/= 0) return
47  if( ipt/=0 ) nlgeom = totallag
48 
49  fstr_ctrl_get_usermaterial = read_user_matl( ctrl, matval )
50  end function fstr_ctrl_get_usermaterial
51 
52  !----------------------------------------------------------------------
54  integer function fstr_ctrl_get_elasticity( ctrl, mattype, nlgeom, matval, dict )
55  integer(kind=kint), intent(in) :: ctrl
56  integer(kind=kint), intent(inout) :: mattype
57  integer(kind=kint), intent(out) :: nlgeom
58  real(kind=kreal),intent(out) :: matval(:)
59  type(dict_struct), pointer :: dict
60 
61  integer(kind=kint) :: i,j, rcode, depends, ipt, n
62  real(kind=kreal), allocatable :: fval(:,:)
63  character(len=HECMW_NAME_LEN) :: data_fmt
64  type( ttable ) :: mattable
65  logical :: isok
66  character(len=256) :: s
67 
69  depends = 0
70  rcode = fstr_ctrl_get_param_ex( ctrl, 'DEPENDENCIES ', '# ', 0, 'I', depends )
71  if( depends>1 ) depends=1 ! temperature depends only currently
72  if( depends > 3 ) stop "We cannot read dependencies>3 right now"
73  nlgeom = totallag !default value
74 
75  if( fstr_ctrl_get_param_ex( ctrl, 'CAUCHY ', '# ', 0, 'E', ipt )/= 0) return
76  if( ipt/=0 ) nlgeom = updatelag
77  if( fstr_ctrl_get_param_ex( ctrl, 'INFINITESIMAL ', '# ', 0, 'E', ipt )/= 0) return
78  if( ipt/=0 ) nlgeom = infinitesimal
79 
80  ! for backward compatibility
81  if( fstr_ctrl_get_param_ex( ctrl, 'INFINITE ', '# ', 0, 'E', ipt )/= 0) return
82  if( ipt/=0 ) then
83  write(*,*) "Warning : !ELASTIC : parameter 'INFINITE' is deprecated." &
84  & // " Please use the replacement parameter 'INFINITESIMAL'"
85  nlgeom = infinitesimal
86  endif
87 
88  ipt=1
89  s = 'ISOTROPIC,ORTHOTROPIC,USER '
90  if( fstr_ctrl_get_param_ex( ctrl, 'TYPE ', s, 0, 'P', ipt ) /= 0 ) return
91 
92  n = fstr_ctrl_get_data_line_n( ctrl )
93  ! ISOTROPIC
94  if( ipt==1 ) then
95  allocate( fval(2+depends,n) )
96  if( depends==0 ) then
97  data_fmt = "RR "
98  block
99  real(kind=kreal) :: fval1(n), fval2(n)
100  fval1 = fval(1,:)
101  fval2 = fval(2,:)
103  fstr_ctrl_get_data_array_ex( ctrl, data_fmt, fval1, fval2 )
104  fval(1,:) = fval1
105  fval(2,:) = fval2
106  end block
107  endif
108  if( depends==1 ) then
109  data_fmt = "RRR "
110  block
111  real(kind=kreal) :: fval1(n), fval2(n), fval3(n)
112  fval1 = fval(1,:)
113  fval2 = fval(2,:)
114  fval3 = fval(3,:)
116  fstr_ctrl_get_data_array_ex( ctrl, data_fmt, fval1, fval2, fval3 )
117  fval(1,:) = fval1
118  fval(2,:) = fval2
119  fval(3,:) = fval3
120  end block
121  endif
122  if( fstr_ctrl_get_elasticity ==0 ) then
123  matval(m_youngs) = fval(1,1)
124  matval(m_poisson) = fval(2,1)
125  call init_table( mattable, depends, 2+depends,n, fval )
126  call dict_add_key( dict, mc_isoelastic, mattable )
127  ! call print_table( mattable, 6 ); pause
128  endif
129  mattype = elastic
130 
131  ! ORTHOTROPIC
132  else if( ipt==2 ) then
133  allocate( fval(9+depends,n) )
134  if( depends==0 ) then
135  data_fmt = "RRRRRRRRR "
136  block
137  real(kind=kreal) :: fval1(n), fval2(n), fval3(n), fval4(n), fval5(n), &
138  fval6(n), fval7(n), fval8(n), fval9(n)
139  fval1 = fval(1,:)
140  fval2 = fval(2,:)
141  fval3 = fval(3,:)
142  fval4 = fval(4,:)
143  fval5 = fval(5,:)
144  fval6 = fval(6,:)
145  fval7 = fval(7,:)
146  fval8 = fval(8,:)
147  fval9 = fval(9,:)
149  fstr_ctrl_get_data_array_ex( ctrl, data_fmt, &
150  fval1, fval2, fval3, fval4, fval5, fval6, fval7, fval8, fval9 )
151  fval(1,:) = fval1
152  fval(2,:) = fval2
153  fval(3,:) = fval3
154  fval(4,:) = fval4
155  fval(5,:) = fval5
156  fval(6,:) = fval6
157  fval(7,:) = fval7
158  fval(8,:) = fval8
159  fval(9,:) = fval9
160  end block
161  else if( depends==1 ) then
162  data_fmt = "RRRRRRRRRR "
163  block
164  real(kind=kreal) :: fval1(n), fval2(n), fval3(n), fval4(n), fval5(n), &
165  fval6(n), fval7(n), fval8(n), fval9(n), fval10(n)
166  fval1 = fval(1,:)
167  fval2 = fval(2,:)
168  fval3 = fval(3,:)
169  fval4 = fval(4,:)
170  fval5 = fval(5,:)
171  fval6 = fval(6,:)
172  fval7 = fval(7,:)
173  fval8 = fval(8,:)
174  fval9 = fval(9,:)
175  fval10 = fval(10,:)
177  fstr_ctrl_get_data_array_ex( ctrl, data_fmt, &
178  fval1, fval2, fval3, fval4, fval5, fval6, fval7, fval8, fval9, fval10 )
179  fval(1,:) = fval1
180  fval(2,:) = fval2
181  fval(3,:) = fval3
182  fval(4,:) = fval4
183  fval(5,:) = fval5
184  fval(6,:) = fval6
185  fval(7,:) = fval7
186  fval(8,:) = fval8
187  fval(9,:) = fval9
188  fval(10,:) = fval10
189  end block
190  endif
191  if( fstr_ctrl_get_elasticity ==0 ) then
192  isok = .true.
193  do i=1,n
194  if( fval(1,i)<=0.d0 .or. fval(2,i)<=0.d0 .or. fval(3,i)<=0.d0 .or. &
195  fval(7,i)<=0.d0 .or. fval(8,i)<=0.d0 .or. fval(9,i)<=0.d0 ) then
196  isok = .false.; fstr_ctrl_get_elasticity=1; exit
197  endif
198  enddo
199  if( isok ) then
200  call init_table( mattable, depends, 9+depends,n, fval )
201  call dict_add_key( dict, mc_orthoelastic, mattable )
202  mattype = mn_orthoelastic
203  endif
204  endif
205 
206  else if( ipt==3 ) then
207  fstr_ctrl_get_elasticity = read_user_matl( ctrl, matval(101:200))
208  mattype = userelastic
209  nlgeom = infinitesimal
210 
211  else
212  stop "ERROR: Material type not supported"
213 
214  endif
215 
216  call finalize_table( mattable )
217  if( allocated(fval) ) deallocate(fval)
218  end function fstr_ctrl_get_elasticity
219 
220 
221  !----------------------------------------------------------------------
223  integer function fstr_ctrl_get_hyperelastic( ctrl, mattype, nlgeom, matval )
224  integer(kind=kint), intent(in) :: ctrl
225  integer(kind=kint), intent(inout) :: mattype
226  integer(kind=kint), intent(out) :: nlgeom
227  real(kind=kreal),intent(out) :: matval(:)
228 
229  integer(kind=kint) :: i,j, rcode, depends, ipt
230  real(kind=kreal), allocatable :: fval(:,:)
231  character(len=HECMW_NAME_LEN) :: data_fmt
232  character(len=256) :: s
233 
235  depends = 0
236  rcode = fstr_ctrl_get_param_ex( ctrl, 'DEPENDENCIES ', '# ', 0, 'I', depends )
237  if( depends > 3 ) stop "We cannot read dependencies>3 right now"
238  nlgeom = totallag !default value
239  if( fstr_ctrl_get_param_ex( ctrl, 'CAUCHY ', '# ', 0, 'E', ipt )/= 0) return
240  if( ipt/=0 ) nlgeom = updatelag
241 
242  ipt=1
243  s = 'NEOHOOKE,MOONEY-RIVLIN,ARRUDA-BOYCE,USER,MOONEY-RIVLIN-ANISO '
244  if( fstr_ctrl_get_param_ex( ctrl, 'TYPE ', s, 0, 'P', ipt ) /= 0 ) return
245 
246  ! NEOHOOKE
247  if( ipt==1 ) then
248  allocate( fval(2,depends+1) )
249  fval =0.0d0
250  if( depends==0 ) then
251  data_fmt = "RR "
252  block
253  real(kind=kreal) :: fval1(depends+1), fval2(depends+1)
254  fval1 = fval(1,:)
255  fval2 = fval(2,:)
257  fstr_ctrl_get_data_array_ex( ctrl, data_fmt, fval1, fval2 )
258  fval(1,:) = fval1
259  fval(2,:) = fval2
260  end block
261  endif
262  if( fstr_ctrl_get_hyperelastic ==0 ) then
263  if( fval(2,1)==0.d0 ) stop "We cannot deal with incompressible material currently"
264  matval(m_plconst1) = fval(1,1)
265  matval(m_plconst2) = 0.d0
266  matval(m_plconst3) = fval(2,1)
267  ! matval(M_YOUNGS) = fval(1,1)
268  ! matval(M_POISSON) = fval(2,1)
269  endif
270  mattype = neohooke
271 
272  ! MOONEY
273  else if( ipt==2 ) then
274  allocate( fval(3,depends+1) )
275  fval =0.0d0
276  if( depends==0 ) then
277  data_fmt = "RRR "
278  block
279  real(kind=kreal) :: fval1(depends+1), fval2(depends+1), fval3(depends+1)
280  fval1 = fval(1,:)
281  fval2 = fval(2,:)
282  fval3 = fval(3,:)
284  fstr_ctrl_get_data_array_ex( ctrl, data_fmt, fval1, fval2, fval3 )
285  fval(1,:) = fval1
286  fval(2,:) = fval2
287  fval(3,:) = fval3
288  end block
289  endif
290  if( fstr_ctrl_get_hyperelastic ==0 ) then
291  matval(m_plconst1) = fval(1,1)
292  matval(m_plconst2) = fval(2,1)
293  matval(m_plconst3) = fval(3,1)
294  endif
295  mattype = mooneyrivlin
296 
297  ! ARRUDA
298  else if( ipt==3 ) then
299  allocate( fval(3,depends+1) )
300  fval =0.0d0
301  if( depends==0 ) then
302  data_fmt = "RRR "
303  block
304  real(kind=kreal) :: fval1(depends+1), fval2(depends+1), fval3(depends+1)
305  fval1 = fval(1,:)
306  fval2 = fval(2,:)
307  fval3 = fval(3,:)
309  fstr_ctrl_get_data_array_ex( ctrl, data_fmt, fval1, fval2, fval3 )
310  fval(1,:) = fval1
311  fval(2,:) = fval2
312  fval(3,:) = fval3
313  end block
314  endif
315  if( fstr_ctrl_get_hyperelastic ==0 ) then
316  matval(m_plconst1) = fval(1,1)
317  matval(m_plconst2) = fval(2,1)
318  matval(m_plconst3) = fval(3,1)
319  endif
320  mattype = arrudaboyce
321 
322  else if( ipt==4 ) then !User
323  fstr_ctrl_get_hyperelastic = read_user_matl( ctrl, matval(101:200))
324  mattype = userhyperelastic
325 
326  ! MOONEY-ORTHO
327  else if( ipt==5 ) then
328  allocate( fval(10,depends+1) )
329  fval =0.0d0
330  if( depends==0 ) then
331  data_fmt = "RRRRRrrrrr "
332  block
333  real(kind=kreal) :: fval1(depends+1), fval2(depends+1), fval3(depends+1), fval4(depends+1), fval5(depends+1), &
334  fval6(depends+1), fval7(depends+1), fval8(depends+1), fval9(depends+1), fval10(depends+1)
335  fval1 = fval(1,:)
336  fval2 = fval(2,:)
337  fval3 = fval(3,:)
338  fval4 = fval(4,:)
339  fval5 = fval(5,:)
340  fval6 = fval(6,:)
341  fval7 = fval(7,:)
342  fval8 = fval(8,:)
343  fval9 = fval(9,:)
344  fval10 = fval(10,:)
346  fstr_ctrl_get_data_array_ex( ctrl, data_fmt, &
347  fval1, fval2, fval3, fval4, fval5, fval6, fval7, fval8, fval9, fval10 )
348  fval(1,:) = fval1
349  fval(2,:) = fval2
350  fval(3,:) = fval3
351  fval(4,:) = fval4
352  fval(5,:) = fval5
353  fval(6,:) = fval6
354  fval(7,:) = fval7
355  fval(8,:) = fval8
356  fval(9,:) = fval9
357  fval(10,:) = fval10
358  end block
359  endif
360  if( fstr_ctrl_get_hyperelastic ==0 ) then
361  matval(m_plconst1) = fval(1,1)
362  matval(m_plconst2) = fval(2,1)
363  matval(m_plconst3) = fval(3,1)
364  matval(m_plconst4) = fval(4,1)
365  matval(m_plconst5) = fval(5,1)
366  matval(m_plconst6) = fval(6,1)
367  matval(m_plconst7) = fval(7,1)
368  matval(m_plconst8) = fval(8,1)
369  matval(m_plconst9) = fval(9,1)
370  matval(m_plconst10) = fval(10,1)
371  endif
372  mattype = mooneyrivlin_aniso
373 
374  endif
375 
376  if( allocated(fval) ) deallocate(fval)
377  end function fstr_ctrl_get_hyperelastic
378 
379 
380  !----------------------------------------------------------------------
382  integer function fstr_ctrl_get_viscoelasticity( ctrl, mattype, nlgeom, dict )
383  integer(kind=kint), intent(in) :: ctrl
384  integer(kind=kint), intent(inout) :: mattype
385  integer(kind=kint), intent(out) :: nlgeom
386  type(dict_struct), pointer :: dict
387 
388  integer(kind=kint) :: i,j, rcode, depends, ipt, n
389  real(kind=kreal), allocatable :: fval(:,:)
390  character(len=HECMW_NAME_LEN) :: data_fmt
391  type( ttable ) :: mattable
392  character(len=256) :: s
393 
395  depends = 0
396  rcode = fstr_ctrl_get_param_ex( ctrl, 'DEPENDENCIES ', '# ', 0, 'I', depends )
397  if( depends>1 ) depends=1 ! temperature depends only currently
398  !depends = 0
399  nlgeom = totallag !default value
400  if( fstr_ctrl_get_param_ex( ctrl, 'INFINITESIMAL ', '# ', 0, 'E', ipt )/= 0) return
401  if( ipt/=0 ) nlgeom = infinitesimal
402 
403  ! for backward compatibility
404  if( fstr_ctrl_get_param_ex( ctrl, 'INFINITE ', '# ', 0, 'E', ipt )/= 0) return
405  if( ipt/=0 ) then
406  write(*,*) "Warning : !VISCOELASTIC : parameter 'INFINITE' is deprecated." &
407  & // " Please use the replacement parameter 'INFINITESIMAL'"
408  nlgeom = infinitesimal
409  endif
410 
411  ipt=1
412  s = 'ISOTROPIC,USER '
413  if( fstr_ctrl_get_param_ex( ctrl, 'TYPE ', s, 0, 'P', ipt ) /= 0 ) return
414  ipt = 1
415 
416  ! ISOTROPIC
417  if( ipt==1 ) then
418  n = fstr_ctrl_get_data_line_n( ctrl )
419  allocate( fval(2+depends,n) )
420  if( depends==0 ) then
421  data_fmt = "RR "
422  block
423  real(kind=kreal) :: fval1(n), fval2(n)
424  fval1 = fval(1,:)
425  fval2 = fval(2,:)
427  fstr_ctrl_get_data_array_ex( ctrl, data_fmt, fval1, fval2 )
428  fval(1,:) = fval1
429  fval(2,:) = fval2
430  end block
431  if( fval(2,1)==0.d0 ) stop "Error in defining viscoelasticity: Relaxation time cannot be zero!"
432  endif
433  if( depends==1 ) then
434  data_fmt = "RRR "
435  block
436  real(kind=kreal) :: fval1(n), fval2(n), fval3(n)
437  fval1 = fval(1,:)
438  fval2 = fval(2,:)
439  fval3 = fval(3,:)
441  fstr_ctrl_get_data_array_ex( ctrl, data_fmt, fval1, fval2, fval3 )
442  fval(1,:) = fval1
443  fval(2,:) = fval2
444  fval(3,:) = fval3
445  end block
446  endif
447  if( fstr_ctrl_get_viscoelasticity ==0 ) then
448  call init_table( mattable, 1, 2+depends,n, fval )
449  call dict_add_key( dict, mc_viscoelastic, mattable )
450  ! call print_table( mattable, 6 ); pause
451  endif
452  mattype = viscoelastic
453 
454  else
455  stop "ERROR: Material type not supported"
456 
457  endif
458 
459  call finalize_table( mattable )
460  if( allocated(fval) ) deallocate(fval)
461  end function fstr_ctrl_get_viscoelasticity
462 
463 
465  integer function fstr_ctrl_get_trs( ctrl, mattype, matval )
466  integer(kind=kint), intent(in) :: ctrl
467  integer(kind=kint), intent(inout) :: mattype
468  real(kind=kreal),intent(out) :: matval(:)
469 
470  integer :: ipt
471  character(len=256) :: s
472 
473  ipt=1
474  s = 'WLF,ARRHENIUS '
475  if( fstr_ctrl_get_param_ex( ctrl, 'DEFINITION ', s, 0, 'P', ipt ) /= 0 ) return
476 
478  = fstr_ctrl_get_data_ex( ctrl, 1, "RRR ", matval(1), matval(2), matval(3) )
479  if( fstr_ctrl_get_trs/=0 ) return
480  mattype = mattype+ipt
481 
482  end function fstr_ctrl_get_trs
483 
484 
485  !----------------------------------------------------------------------
487  integer function fstr_ctrl_get_plasticity( ctrl, mattype, nlgeom, matval, mattable, dict )
488  integer(kind=kint), intent(in) :: ctrl
489  integer(kind=kint), intent(inout) :: mattype
490  integer(kind=kint), intent(out) :: nlgeom
491  real(kind=kreal),intent(out) :: matval(:)
492  real(kind=kreal), pointer :: mattable(:)
493  type(dict_struct), pointer :: dict
494 
495  integer(kind=kint) :: i, n, rcode, depends, ipt, hipt
496  real(kind=kreal), allocatable :: fval(:,:)
497  real(kind=kreal) :: phi, psi
498  character(len=HECMW_NAME_LEN) :: data_fmt
499  character(len=256) :: s
500  type( ttable ) :: mttable
501  real(kind=kreal), parameter :: pi=3.14159265358979d0
502 
504  ipt = 0; hipt = 0
505 
506  depends = 0
507  rcode = fstr_ctrl_get_param_ex( ctrl, 'DEPENDENCIES ', '# ', 0, 'I', depends )
508  if( depends>1 ) depends = 1 ! we consider temperature dependence only currently
509  if( depends > 3 ) stop "We cannot read dependencies>3 right now"
510  nlgeom = updatelag !default value
511  if( fstr_ctrl_get_param_ex( ctrl, 'KIRCHHOFF ', '# ', 0, 'E', ipt )/= 0) return
512  ! rcode = fstr_ctrl_get_param_ex( ctrl, 'FILE ', '# ', 0, 'S', fname )
513  if( ipt/=0 ) nlgeom = totallag
514  if( fstr_ctrl_get_param_ex( ctrl, 'INFINITESIMAL ', '# ', 0, 'E', ipt )/= 0) return
515  if( ipt/=0 ) nlgeom = infinitesimal
516 
517  ! for backward compatibility
518  if( fstr_ctrl_get_param_ex( ctrl, 'INFINITE ', '# ', 0, 'E', ipt )/= 0) return
519  if( ipt/=0 ) then
520  write(*,*) "Warning : !PLASTIC : parameter 'INFINITE' is deprecated." &
521  & // " Please use the replacement parameter 'INFINITESIMAL'"
522  nlgeom = infinitesimal
523  endif
524 
525  call setdigit( 1, 1, mattype )
526  call setdigit( 2, 2, mattype )
527 
528  ! hardening
529  s = 'BILINEAR,MULTILINEAR,SWIFT,RAMBERG-OSGOOD,KINEMATIC,COMBINED '
530  if( fstr_ctrl_get_param_ex( ctrl, 'HARDEN ', s , 0, 'P', hipt ) /= 0 ) return
531  if( hipt==0 ) hipt=1 ! default: linear hardening
532  call setdigit( 5, hipt-1, mattype )
533 
534  ! yield function
535  s = 'MISES,MOHR-COULOMB,DRUCKER-PRAGER,USER '
536  call setdigit( 2, 2, mattype )
537  if( fstr_ctrl_get_param_ex( ctrl, 'YIELD ', s , 0, 'P', ipt ) /= 0 ) return
538  if( ipt==0 ) ipt=1 ! default: mises yield function
539  call setdigit( 4, ipt-1, mattype )
540 
541  n = fstr_ctrl_get_data_line_n( ctrl )
542  if( n == 0 ) return ! fail in reading plastic
543  if( hipt==2 .and. n<2 ) return ! not enough data
544  if( ( ipt==2 .or. ipt==3 ) .and. hipt>2 ) hipt = 1
545 
546  select case (ipt)
547  case (1) !Mises
548  select case (hipt)
549  case (1,5) ! linear hardening, kinematic hardening
550  allocate( fval(2,n) )
551  data_fmt = "RR "
552  block
553  real(kind=kreal) :: fval1(n), fval2(n)
554  fval1 = fval(1,:)
555  fval2 = fval(2,:)
557  fstr_ctrl_get_data_array_ex( ctrl, data_fmt, fval1, fval2 )
558  fval(1,:) = fval1
559  fval(2,:) = fval2
560  end block
561  if( fstr_ctrl_get_plasticity ==0 ) then
562  matval(m_plconst1) = fval(1,1)
563  if(hipt==1) then
564  matval(m_plconst2) = fval(2,1)
565  else
566  matval(m_plconst2) = 0.d0
567  matval(m_plconst3) = fval(2,1)
568  endif
569  endif
570  case (2) ! multilinear approximation
571  allocate( fval(depends+2,n) )
572  if( depends==0 ) then
573  data_fmt = "RR "
574  block
575  real(kind=kreal) :: fval1(n), fval2(n)
576  fval1 = fval(1,:)
577  fval2 = fval(2,:)
579  fstr_ctrl_get_data_array_ex( ctrl, data_fmt, fval1, fval2 )
580  fval(1,:) = fval1
581  fval(2,:) = fval2
582  end block
583  if( fstr_ctrl_get_plasticity ==0 ) then
584  if( fval(2,1)/=0.d0 ) then
585  print *, "Multilinear hardening: First plastic strain must be zero"
586  stop
587  endif
588  do i=1,n
589  if( fval(2,i)<0.0 ) &
590  stop "Multilinear hardening: Error in plastic strain definition"
591  enddo
592  call init_table( mttable,1, 2+depends, n, fval )
593  call dict_add_key( dict, mc_yield, mttable )
594 
595  endif
596  else ! depends==1
597  data_fmt = "RRR "
598  block
599  real(kind=kreal) :: fval1(n), fval2(n), fval3(n)
600  fval1 = fval(1,:)
601  fval2 = fval(2,:)
602  fval3 = fval(3,:)
604  fstr_ctrl_get_data_array_ex( ctrl, data_fmt, fval1, fval2, fval3 )
605  fval(1,:) = fval1
606  fval(2,:) = fval2
607  fval(3,:) = fval3
608  end block
609  if( fstr_ctrl_get_plasticity ==0 ) then
610  call init_table( mttable,2, 2+depends,n, fval )
611  call dict_add_key( dict, mc_yield, mttable )
612  endif
613  endif
614  case (3, 4, 6) ! swift, Ramberg-Osgood, Combined
615  allocate( fval(3,1) )
616  data_fmt = "RRR "
617  block
618  real(kind=kreal) :: fval1(1), fval2(1), fval3(1)
619  fval1 = fval(1,:)
620  fval2 = fval(2,:)
621  fval3 = fval(3,:)
623  fstr_ctrl_get_data_array_ex( ctrl, data_fmt, fval1, fval2, fval3 )
624  fval(1,:) = fval1
625  fval(2,:) = fval2
626  fval(3,:) = fval3
627  end block
628  if( fstr_ctrl_get_plasticity ==0 ) then
629  matval(m_plconst1) = fval(1,1)
630  matval(m_plconst2) = fval(2,1)
631  matval(m_plconst3) = fval(3,1)
632  endif
633  case default
634  print *, "Error in hardening definition!"
635  stop
636  end select
637  case (2, 3) ! Mohr-Coulomb, Drucker-Prager
638  select case (hipt)
639  case (1) ! linear hardening
640  allocate( fval(4,n) )
641  data_fmt = "RRrr "
642  fval(4,:) = -1.d0
643  block
644  real(kind=kreal) :: fval1(n), fval2(n), fval3(n), fval4(n)
645  fval1 = fval(1,:)
646  fval2 = fval(2,:)
647  fval3 = fval(3,:)
648  fval4 = fval(4,:)
650  = fstr_ctrl_get_data_array_ex( ctrl, data_fmt, fval1, fval2, fval3, fval4 )
651  fval(1,:) = fval1
652  fval(2,:) = fval2
653  fval(3,:) = fval3
654  fval(4,:) = fval4
655  end block
656  if( fstr_ctrl_get_plasticity ==0 ) then
657  matval(m_plconst1) = fval(1,1) ! c
658  matval(m_plconst2) = fval(3,1) ! H
659  phi = fval(2,1)*pi/180.d0
660  if( fval(4,1) >= 0.d0 ) then
661  psi = fval(4,1)*pi/180.d0
662  else
663  psi = phi
664  endif
665  endif
666  case (2) ! multilinear hardening
667  if( depends>0 ) then
668  stop "Mohr-Coulomb and Drucker-Prager do not support temperature dependency"
669  endif
670  allocate( fval(2,n) )
671  data_fmt = "Rr "
672  fval(2,:) = -1.d0
673  block
674  real(kind=kreal) :: fval1(n), fval2(n)
675  fval1 = fval(1,:)
676  fval2 = fval(2,:)
678  fstr_ctrl_get_data_array_ex( ctrl, data_fmt, fval1, fval2 )
679  fval(1,:) = fval1
680  fval(2,:) = fval2
681  end block
682  if( fstr_ctrl_get_plasticity ==0 ) then
683  phi =fval(1,1)*pi/180.d0
684  if( fval(2,1) >= 0.d0 ) then
685  psi = fval(2,1)*pi/180.d0
686  else
687  psi = phi
688  endif
689  if( fval(2,2)/=0.d0 ) then
690  print *, "Multilinear hardening: First plastic strain must be zero"
691  stop
692  endif
693  do i=2,n
694  if( fval(2,i)<0.0 ) &
695  stop "Multilinear hardening: Error in plastic strain definition"
696  enddo
697  call init_table( mttable,1, 2, n-1, fval(1:2,2:n) )
698  call dict_add_key( dict, mc_yield, mttable )
699  endif
700  end select
701  if( ipt==3 ) then ! Drucker-Prager
702  matval(m_plconst3) = 2.d0*sin(phi)/ ( sqrt(3.d0)*(3.d0+sin(phi)) ) ! eta
703  matval(m_plconst4) = 6.d0*cos(phi)/ ( sqrt(3.d0)*(3.d0+sin(phi)) ) ! xi
704  matval(m_plconst5) = 2.d0*sin(psi)/ ( sqrt(3.d0)*(3.d0+sin(psi)) ) ! etabar
705  else ! Mohr-Coulomb
706  matval(m_plconst3) = phi
707  matval(m_plconst4) = psi
708  endif
709 
710  case(4)
711  fstr_ctrl_get_plasticity = read_user_matl( ctrl, matval(101:200) )
712 
713  case default
714  stop "Yield function not supported"
715  end select
716 
717  if( allocated(fval) ) deallocate(fval)
718  call finalize_table( mttable )
719  end function fstr_ctrl_get_plasticity
720 
721 
722  !----------------------------------------------------------------------
724  integer function fstr_ctrl_get_viscoplasticity( ctrl, mattype, nlgeom, dict )
725  integer(kind=kint), intent(in) :: ctrl
726  integer(kind=kint), intent(inout) :: mattype
727  integer(kind=kint), intent(out) :: nlgeom
728  type(dict_struct), pointer :: dict
729 
730  integer(kind=kint) :: i,j, rcode, depends, ipt, n
731  real(kind=kreal), allocatable :: fval(:,:)
732  character(len=HECMW_NAME_LEN) :: data_fmt
733  type( ttable ) :: mattable
734  character(len=256) :: s
735 
737  depends = 0
738  rcode = fstr_ctrl_get_param_ex( ctrl, 'DEPENDENCIES ', '# ', 0, 'I', depends )
739  if( depends>1 ) depends=1 ! temperature depends only currently
740  nlgeom = updatelag !default value
741  if( fstr_ctrl_get_param_ex( ctrl, 'KIRCHHOFF ', '# ', 0, 'E', ipt )/= 0) return
742  if( ipt/=0 ) nlgeom = totallag
743 
744  ipt=1
745  s = 'NORTON,USER '
746  if( fstr_ctrl_get_param_ex( ctrl, 'TYPE ', s, 0, 'P', ipt ) /= 0 ) return
747  ipt = 1
748 
749  ! NORTON
750  if( ipt==1 ) then
751  n = fstr_ctrl_get_data_line_n( ctrl )
752  allocate( fval(3+depends,n) )
753  if( depends==0 ) then
754  data_fmt = "RRR "
755  block
756  real(kind=kreal) :: fval1(n), fval2(n), fval3(n)
757  fval1 = fval(1,:)
758  fval2 = fval(2,:)
759  fval3 = fval(3,:)
761  fstr_ctrl_get_data_array_ex( ctrl, data_fmt, fval1, fval2, fval3 )
762  fval(1,:) = fval1
763  fval(2,:) = fval2
764  fval(3,:) = fval3
765  end block
766  endif
767  if( depends==1 ) then
768  data_fmt = "RRRR "
769  block
770  real(kind=kreal) :: fval1(n), fval2(n), fval3(n), fval4(n)
771  fval1 = fval(1,:)
772  fval2 = fval(2,:)
773  fval3 = fval(3,:)
774  fval4 = fval(4,:)
776  fstr_ctrl_get_data_array_ex( ctrl, data_fmt, fval1, fval2, fval3, fval4 )
777  fval(1,:) = fval1
778  fval(2,:) = fval2
779  fval(3,:) = fval3
780  fval(4,:) = fval4
781  end block
782  endif
783  if( fstr_ctrl_get_viscoplasticity ==0 ) then
784  call init_table( mattable, depends, 3+depends,n, fval )
785  call dict_add_key( dict, mc_norton, mattable )
786  ! call print_table( mattable, 6 ); pause
787  endif
788  mattype = norton
789 
790  else
791  stop "ERROR: Material type not supported"
792 
793  endif
794 
795  call finalize_table( mattable )
796  if( allocated(fval) ) deallocate(fval)
797  end function fstr_ctrl_get_viscoplasticity
798 
799  !----------------------------------------------------------------------
801  integer function fstr_ctrl_get_density( ctrl, matval )
802  integer(kind=kint), intent(in) :: ctrl
803  real(kind=kreal),intent(out) :: matval(:)
804 
805  integer(kind=kint) :: i, rcode, depends
806  real(kind=kreal), allocatable :: fval(:,:)
807  character(len=HECMW_NAME_LEN) :: data_fmt
808 
809  data_fmt = "R "
810 
812 
813  depends = 0
814  rcode = fstr_ctrl_get_param_ex( ctrl, 'DEPENDENCIES ', '# ', 0, 'I', depends )
815  if( depends>1 ) depends = 1 ! we consider temperature dependence only currently
816 
817  allocate( fval(1,depends+1) )
818  do i=2,1+depends
819  data_fmt = data_fmt //"R "
820  enddo
821  block
822  real(kind=kreal) :: fval1(depends+1)
823  fval1 = fval(1,:)
825  = fstr_ctrl_get_data_array_ex( ctrl, data_fmt, fval1 )
826  fval(1,:) = fval1
827  end block
828  if( fstr_ctrl_get_density==0 ) matval(m_density) = fval(1,1)
829 
830  if( allocated(fval) ) deallocate(fval)
831 
832  end function fstr_ctrl_get_density
833 
834 
835  !----------------------------------------------------------------------
837  integer function fstr_ctrl_get_expansion_coeff( ctrl, matval, dict )
838  integer(kind=kint), intent(in) :: ctrl
839  real(kind=kreal),intent(out) :: matval(:)
840  type(dict_struct), pointer :: dict
841 
842  integer(kind=kint) :: i, n, rcode, depends, ipt
843  real(kind=kreal), allocatable :: fval(:,:)
844  type( ttable ) :: mttable
845  character(len=HECMW_NAME_LEN) :: data_fmt, ss
846 
847  data_fmt = "R "
848 
850  n = fstr_ctrl_get_data_line_n( ctrl )
851  if( n == 0 ) return ! fail in reading plastic
852 
853  ss = 'ISOTROPIC,ORTHOTROPIC '
854  ipt = 1 !default
855  if( fstr_ctrl_get_param_ex( ctrl, 'TYPE ', ss, 0, 'P', ipt ) /= 0 ) return
856 
857  depends = 0
858  rcode = fstr_ctrl_get_param_ex( ctrl, 'DEPENDENCIES ', '# ', 0, 'I', depends )
859  if( depends>1 ) depends = 1 ! we consider temperature dependence only currently
860 
861  if( ipt==1 ) then
862  allocate( fval(depends+1, n) )
863  do i=2,1+depends
864  data_fmt = data_fmt //"R "
865  enddo
866  if( depends==0 ) then
867  block
868  real(kind=kreal) :: fval1(n)
869  fval1 = fval(1,:)
871  fstr_ctrl_get_data_array_ex( ctrl, "R ", fval1 )
872  fval(1,:) = fval1
873  end block
874  else
875  block
876  real(kind=kreal) :: fval1(n), fval2(n)
877  fval1 = fval(1,:)
878  fval2 = fval(2,:)
880  fstr_ctrl_get_data_array_ex( ctrl, "RR ", fval1, fval2 )
881  fval(1,:) = fval1
882  fval(2,:) = fval2
883  end block
884  endif
885  if( fstr_ctrl_get_expansion_coeff==0 ) then
886  matval(m_exapnsion) = fval(1,1)
887  call init_table( mttable,depends, 1+depends, n, fval )
888  call dict_add_key( dict, mc_themoexp, mttable )
889  endif
890  else
891  allocate( fval(3+depends,n) )
892  do i=2,3+depends
893  data_fmt = trim(data_fmt) //"R "
894  enddo
895  if( depends==0 ) then
896  block
897  real(kind=kreal) :: fval1(n), fval2(n), fval3(n)
898  fval1 = fval(1,:)
899  fval2 = fval(2,:)
900  fval3 = fval(3,:)
902  fstr_ctrl_get_data_array_ex( ctrl, data_fmt, fval1, fval2, fval3 )
903  fval(1,:) = fval1
904  fval(2,:) = fval2
905  fval(3,:) = fval3
906  end block
907  elseif( depends==1 ) then
908  block
909  real(kind=kreal) :: fval1(n), fval2(n), fval3(n), fval4(n)
910  fval1 = fval(1,:)
911  fval2 = fval(2,:)
912  fval3 = fval(3,:)
913  fval4 = fval(4,:)
915  fstr_ctrl_get_data_array_ex( ctrl, data_fmt, fval1, fval2, fval3, fval4 )
916  fval(1,:) = fval1
917  fval(2,:) = fval2
918  fval(3,:) = fval3
919  fval(4,:) = fval4
920  end block
921  endif
922  if( fstr_ctrl_get_expansion_coeff==0 ) then
923  call init_table( mttable, depends, 3+depends,n, fval )
924  if( fstr_ctrl_get_expansion_coeff==0 ) call dict_add_key( dict, mc_orthoexp, mttable )
925  endif
926  endif
927 
928  call finalize_table( mttable )
929  if( allocated(fval) ) deallocate(fval)
930  end function fstr_ctrl_get_expansion_coeff
931 
932 
933  integer function read_user_matl( ctrl, matval )
934  integer(kind=kint), intent(in) :: ctrl
935  real(kind=kreal),intent(out) :: matval(:)
936 
937  integer(kind=kint) :: n, i, j
938  real(kind=kreal) :: fval(10,10)
939 
940  read_user_matl = -1
941 
942  n = fstr_ctrl_get_data_line_n( ctrl )
943  if( n > 10 ) stop "Num of data lines for user-defined material exceeds 10"
944  fval =0.d0
945  block
946  real(kind=kreal) :: fval1(10), fval2(10), fval3(10), fval4(10), fval5(10), &
947  fval6(10), fval7(10), fval8(10), fval9(10), fval10(10)
948  fval1 = fval(1,:)
949  fval2 = fval(2,:)
950  fval3 = fval(3,:)
951  fval4 = fval(4,:)
952  fval5 = fval(5,:)
953  fval6 = fval(6,:)
954  fval7 = fval(7,:)
955  fval8 = fval(8,:)
956  fval9 = fval(9,:)
957  fval10 = fval(10,:)
958  if( fstr_ctrl_get_data_array_ex( ctrl, 'rrrrrrrrrr ', fval1, fval2, fval3, &
959  fval4, fval5, fval6, fval7, fval8, fval9, fval10 ) /= 0 ) return
960  fval(1,:) = fval1
961  fval(2,:) = fval2
962  fval(3,:) = fval3
963  fval(4,:) = fval4
964  fval(5,:) = fval5
965  fval(6,:) = fval6
966  fval(7,:) = fval7
967  fval(8,:) = fval8
968  fval(9,:) = fval9
969  fval(10,:) = fval10
970  end block
971  do i=1,10
972  do j=1,10
973  matval((i-1)*10+j)=fval(j,i)
974  enddo
975  enddo
976 
977  read_user_matl = 0
978  end function read_user_matl
979 
980  !----------------------------------------------------------------------
982  integer function fstr_ctrl_get_fluid( ctrl, mattype, nlgeom, matval, dict )
983  integer(kind=kint), intent(in) :: ctrl
984  integer(kind=kint), intent(inout) :: mattype
985  integer(kind=kint), intent(out) :: nlgeom
986  real(kind=kreal),intent(out) :: matval(:)
987  type(dict_struct), pointer :: dict
988 
989  integer(kind=kint) :: i,j, rcode, depends, ipt, n
990  real(kind=kreal), allocatable :: fval(:,:)
991  character(len=HECMW_NAME_LEN) :: data_fmt
992  type( ttable ) :: mattable
993  logical :: isok
994  character(len=256) :: s
995 
997  depends = 0
998  rcode = fstr_ctrl_get_param_ex( ctrl, 'DEPENDENCIES ', '# ', 0, 'I', depends )
999  if( depends>1 ) depends=1 ! temperature depends only currently
1000  if( depends > 3 ) stop "We cannot read dependencies>3 right now"
1001  nlgeom = totallag !default value
1002 
1003  ipt=1
1004  s = 'INCOMP_NEWTONIAN '
1005  if( fstr_ctrl_get_param_ex( ctrl, 'TYPE ', s, 0, 'P', ipt ) /= 0 ) return
1006 
1007  n = fstr_ctrl_get_data_line_n( ctrl )
1008  ! ISOTROPIC
1009  if( ipt==1 ) then
1010  allocate( fval(1+depends,n) )
1011  if( depends==0 ) then
1012  data_fmt = "R "
1013  block
1014  real(kind=kreal) :: fval1(n)
1015  fval1 = fval(1,:)
1017  fstr_ctrl_get_data_array_ex( ctrl, data_fmt, fval1 )
1018  fval(1,:) = fval1
1019  end block
1020  endif
1021  if( depends==1 ) then
1022  data_fmt = "RR "
1023  block
1024  real(kind=kreal) :: fval1(n), fval2(n)
1025  fval1 = fval(1,:)
1026  fval2 = fval(2,:)
1028  fstr_ctrl_get_data_array_ex( ctrl, data_fmt, fval1, fval2 )
1029  fval(1,:) = fval1
1030  fval(2,:) = fval2
1031  end block
1032  endif
1033  if( fstr_ctrl_get_fluid ==0 ) then
1034  matval(m_viscocity) = fval(1,1)
1035  call init_table( mattable, depends, 1+depends,n, fval )
1036  call dict_add_key( dict, mc_incomp_newtonian, mattable )
1037  ! call print_table( mattable, 6 ); pause
1038  endif
1039  mattype = incomp_newtonian
1040 
1041  else
1042  stop "ERROR: Material type not supported"
1043 
1044  endif
1045 
1046  call finalize_table( mattable )
1047  if( allocated(fval) ) deallocate(fval)
1048  end function fstr_ctrl_get_fluid
1049 
1050  !----------------------------------------------------------------------
1052  integer function fstr_ctrl_get_spring_d( ctrl, mattype, nlgeom, matval_i, dict )
1053  integer(kind=kint), intent(in) :: ctrl
1054  integer(kind=kint), intent(inout) :: mattype
1055  integer(kind=kint), intent(out) :: nlgeom
1056  integer(kind=kint), intent(out) :: matval_i(:)
1057  type(dict_struct), pointer :: dict
1058 
1059  integer(kind=kint) :: i,j, rcode, depends, ipt, n, dof1, dof2
1060  real(kind=kreal), allocatable :: fval(:,:)
1061  integer(kind=kint), allocatable :: ival(:,:)
1062  character(len=HECMW_NAME_LEN) :: data_fmt
1063  type( ttable ) :: mattable
1064  logical :: isok
1065  character(len=DICT_KEY_LENGTH) :: spkey
1066 
1068  depends = 0
1069  !rcode = fstr_ctrl_get_param_ex( ctrl, 'DEPENDENCIES ', '# ', 0, 'I', depends )
1070  !if( depends>1 ) depends=1 ! temperature depends only currently
1071  nlgeom = infinitesimal !default value
1072 
1073  n = fstr_ctrl_get_data_line_n( ctrl )
1074  allocate( fval(1+depends,n), ival(2,n) )
1075  if( depends==0 ) then
1076  data_fmt = "IIR "
1077  block
1078  integer(kind=kint) :: ival1(n), ival2(n)
1079  real(kind=kreal) :: fval1(n)
1080  ival1 = ival(1,:)
1081  ival2 = ival(2,:)
1082  fval1 = fval(1,:)
1084  fstr_ctrl_get_data_array_ex( ctrl, data_fmt, ival1, ival2, fval1 )
1085  ival(1,:) = ival1
1086  ival(2,:) = ival2
1087  fval(1,:) = fval1
1088  end block
1089  !else if( depends==1 ) then
1090  ! data_fmt = "IIRR "
1091  ! fstr_ctrl_get_SPRING_D = &
1092  ! fstr_ctrl_get_data_array_ex( ctrl, data_fmt, ival(1,:), ival(2,:), fval(1,:), fval(2,:) )
1093  endif
1094  if( fstr_ctrl_get_spring_d ==0 ) then
1095  matval_i(m_spring_d_ndoffset+1) = n ! number of dof pairs
1096  do i=1,n
1097  matval_i(m_spring_d_ndoffset+2*i ) = ival(1,i) ! first dof
1098  matval_i(m_spring_d_ndoffset+2*i+1) = ival(2,i) ! second dof
1099  call init_table( mattable, depends, 1+depends, 1, fval(1,i:i) )
1100  write(spkey,'(A,I0)') trim(mc_spring),i
1101  call dict_add_key( dict, trim(spkey), mattable )
1102  call finalize_table( mattable )
1103  enddo
1104  endif
1105 
1106  mattype = connector
1107  if( allocated(fval) ) deallocate(fval)
1108  if( allocated(ival) ) deallocate(ival)
1109  end function fstr_ctrl_get_spring_d
1110 
1111  !----------------------------------------------------------------------
1113  integer function fstr_ctrl_get_spring_a( ctrl, mattype, nlgeom, matval_i, dict )
1114  integer(kind=kint), intent(in) :: ctrl
1115  integer(kind=kint), intent(inout) :: mattype
1116  integer(kind=kint), intent(out) :: nlgeom
1117  integer(kind=kint), intent(out) :: matval_i(:)
1118  type(dict_struct), pointer :: dict
1119 
1120  integer(kind=kint) :: i,j, rcode, depends, ipt, n, dof1, dof2
1121  real(kind=kreal), allocatable :: fval(:,:)
1122  character(len=HECMW_NAME_LEN) :: data_fmt
1123  type( ttable ) :: mattable
1124  logical :: isok
1125  character(len=DICT_KEY_LENGTH) :: spkey
1126 
1128  depends = 0
1129  nlgeom = infinitesimal !default value
1130 
1131  n = fstr_ctrl_get_data_line_n( ctrl )
1132  allocate( fval(1+depends,n) )
1133  if( depends==0 ) then
1134  data_fmt = "R "
1135  block
1136  real(kind=kreal) :: fval1(n)
1137  fval1 = fval(1,:)
1139  fstr_ctrl_get_data_array_ex( ctrl, data_fmt, fval1 )
1140  fval(1,:) = fval1
1141  end block
1142  !else if( depends==1 ) then
1143  ! data_fmt = "IIRR "
1144  ! fstr_ctrl_get_SPRING_A = &
1145  ! fstr_ctrl_get_data_array_ex( ctrl, data_fmt, fval(1,:), fval(2,:) )
1146  endif
1147  if( fstr_ctrl_get_spring_a ==0 ) then
1148  matval_i(m_spring_a_ndoffset+1) = 1 ! number of dof pairs
1149  do i=1,n
1150  call init_table( mattable, depends, 1+depends, 1, fval(1,i:i) )
1151  write(spkey,'(A,A)') trim(mc_spring),'_A'
1152  call dict_add_key( dict, trim(spkey), mattable )
1153  call finalize_table( mattable )
1154  enddo
1155  endif
1156 
1157  mattype = connector
1158  if( allocated(fval) ) deallocate(fval)
1159  end function fstr_ctrl_get_spring_a
1160 
1161  !----------------------------------------------------------------------
1163  integer function fstr_ctrl_get_dashpot_d( ctrl, mattype, nlgeom, matval_i, dict )
1164  integer(kind=kint), intent(in) :: ctrl
1165  integer(kind=kint), intent(inout) :: mattype
1166  integer(kind=kint), intent(out) :: nlgeom
1167  integer(kind=kint), intent(out) :: matval_i(:)
1168  type(dict_struct), pointer :: dict
1169 
1170  integer(kind=kint) :: i,j, rcode, depends, ipt, n, dof1, dof2
1171  real(kind=kreal), allocatable :: fval(:,:)
1172  integer(kind=kint), allocatable :: ival(:,:)
1173  character(len=HECMW_NAME_LEN) :: data_fmt
1174  type( ttable ) :: mattable
1175  logical :: isok
1176  character(len=DICT_KEY_LENGTH) :: spkey
1177 
1179  depends = 0
1180  !rcode = fstr_ctrl_get_param_ex( ctrl, 'DEPENDENCIES ', '# ', 0, 'I', depends )
1181  !if( depends>1 ) depends=1 ! temperature depends only currently
1182  nlgeom = infinitesimal !default value
1183 
1184  n = fstr_ctrl_get_data_line_n( ctrl )
1185  allocate( fval(1+depends,n), ival(2,n) )
1186  if( depends==0 ) then
1187  data_fmt = "IIR "
1188  block
1189  integer(kind=kint) :: ival1(n), ival2(n)
1190  real(kind=kreal) :: fval1(n)
1191  ival1 = ival(1,:)
1192  ival2 = ival(2,:)
1193  fval1 = fval(1,:)
1195  fstr_ctrl_get_data_array_ex( ctrl, data_fmt, ival1, ival2, fval1 )
1196  ival(1,:) = ival1
1197  ival(2,:) = ival2
1198  fval(1,:) = fval1
1199  end block
1200  !else if( depends==1 ) then
1201  ! data_fmt = "IIRR "
1202  ! fstr_ctrl_get_DASHPOT_D = &
1203  ! fstr_ctrl_get_data_array_ex( ctrl, data_fmt, ival(1,:), ival(2,:), fval(1,:), fval(2,:) )
1204  endif
1205  if( fstr_ctrl_get_dashpot_d ==0 ) then
1206  matval_i(m_dashpot_d_ndoffset+1) = n
1207  do i=1,n
1208  matval_i(m_dashpot_d_ndoffset+2*i ) = ival(1,i) ! first dof
1209  matval_i(m_dashpot_d_ndoffset+2*i+1) = ival(2,i) ! second dof
1210  call init_table( mattable, depends, 1+depends, 1, fval(1,i:i) )
1211  write(spkey,'(A,I0)') trim(mc_dashpot),i
1212  call dict_add_key( dict, trim(spkey), mattable )
1213  call finalize_table( mattable )
1214  enddo
1215  endif
1216 
1217  mattype = connector
1218  if( allocated(fval) ) deallocate(fval)
1219  if( allocated(ival) ) deallocate(ival)
1220  end function fstr_ctrl_get_dashpot_d
1221 
1222  !----------------------------------------------------------------------
1224  integer function fstr_ctrl_get_dashpot_a( ctrl, mattype, nlgeom, matval_i, dict )
1225  integer(kind=kint), intent(in) :: ctrl
1226  integer(kind=kint), intent(inout) :: mattype
1227  integer(kind=kint), intent(out) :: nlgeom
1228  integer(kind=kint), intent(out) :: matval_i(:)
1229  type(dict_struct), pointer :: dict
1230 
1231  integer(kind=kint) :: i,j, rcode, depends, ipt, n, dof1, dof2
1232  real(kind=kreal), allocatable :: fval(:,:)
1233  character(len=HECMW_NAME_LEN) :: data_fmt
1234  type( ttable ) :: mattable
1235  logical :: isok
1236  character(len=DICT_KEY_LENGTH) :: spkey
1237 
1239  depends = 0
1240  nlgeom = infinitesimal !default value
1241 
1242  n = fstr_ctrl_get_data_line_n( ctrl )
1243  allocate( fval(1+depends,n) )
1244  if( depends==0 ) then
1245  data_fmt = "R "
1246  block
1247  real(kind=kreal) :: fval1(n)
1248  fval1 = fval(1,:)
1250  fstr_ctrl_get_data_array_ex( ctrl, data_fmt, fval1 )
1251  fval(1,:) = fval1
1252  end block
1253  !else if( depends==1 ) then
1254  ! data_fmt = "IIRR "
1255  ! fstr_ctrl_get_DASHPOT_A = &
1256  ! fstr_ctrl_get_data_array_ex( ctrl, data_fmt, fval(1,:), fval(2,:) )
1257  endif
1258  if( fstr_ctrl_get_dashpot_a ==0 ) then
1259  matval_i(m_dashpot_a_ndoffset+1) = 1 ! number of dof pairs
1260  do i=1,n
1261  call init_table( mattable, depends, 1+depends, 1, fval(1,i:i) )
1262  write(spkey,'(A,A)') trim(mc_dashpot),'_A'
1263  call dict_add_key( dict, trim(spkey), mattable )
1264  call finalize_table( mattable )
1265  enddo
1266  endif
1267 
1268  mattype = connector
1269  if( allocated(fval) ) deallocate(fval)
1270  end function fstr_ctrl_get_dashpot_a
1271 
1272 end module fstr_ctrl_material
int fstr_ctrl_get_param_ex(int *ctrl, const char *param_name, const char *value_list, int *necessity, char *type, void *val)
int fstr_ctrl_get_data_line_n(int *ctrl)
int fstr_ctrl_get_data_array_ex(int *ctrl, const char *format,...)
int fstr_ctrl_get_data_ex(int *ctrl, int *line_no, const char *format,...)
This module manages read in of various material properties.
integer function fstr_ctrl_get_dashpot_d(ctrl, mattype, nlgeom, matval_i, dict)
Read in !DASHPOT_D.
integer function fstr_ctrl_get_hyperelastic(ctrl, mattype, nlgeom, matval)
Read in !HYPERELASTIC.
integer function fstr_ctrl_get_viscoelasticity(ctrl, mattype, nlgeom, dict)
Read in !VISCOELASTIC.
integer function fstr_ctrl_get_viscoplasticity(ctrl, mattype, nlgeom, dict)
Read in !CREEP.
integer function fstr_ctrl_get_usermaterial(ctrl, mattype, nlgeom, nstatus, matval)
Read in !USER_MATERIAL.
integer function fstr_ctrl_get_expansion_coeff(ctrl, matval, dict)
Read in !EXPANSION_COEFF.
integer function fstr_ctrl_get_trs(ctrl, mattype, matval)
Read in !TRS.
integer function fstr_ctrl_get_elasticity(ctrl, mattype, nlgeom, matval, dict)
Read in !ELASTIC.
integer function fstr_ctrl_get_plasticity(ctrl, mattype, nlgeom, matval, mattable, dict)
Read in !PLASTIC.
integer function fstr_ctrl_get_dashpot_a(ctrl, mattype, nlgeom, matval_i, dict)
Read in !DASHPOT_A.
integer function fstr_ctrl_get_material(ctrl, matname)
Read in !MATERIAL.
integer function fstr_ctrl_get_density(ctrl, matval)
Read in !DENSITY.
integer function fstr_ctrl_get_spring_a(ctrl, mattype, nlgeom, matval_i, dict)
Read in !SPRING_A.
integer function fstr_ctrl_get_fluid(ctrl, mattype, nlgeom, matval, dict)
Read in !FLUID.
integer function fstr_ctrl_get_spring_d(ctrl, mattype, nlgeom, matval_i, dict)
Read in !SPRING_D.
Definition: hecmw.f90:6
This module provides data structure table which would be dictionaried afterwards.
Definition: ttable.f90:7
This module summarizes all information of material properties.
Definition: material.f90:6
integer(kind=kint), parameter totallag
Definition: material.f90:14
integer(kind=kint), parameter infinitesimal
Definition: material.f90:13
integer(kind=kint), parameter usermaterial
Definition: material.f90:63
integer(kind=kint), parameter updatelag
Definition: material.f90:15