13 private :: read_user_matl
21 integer(kind=kint),
intent(in) :: ctrl
22 character(len=*),
intent(out) :: matname
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(:)
37 integer(kind=kint) :: ipt
38 character(len=HECMW_NAME_LEN) :: data_fmt
39 character(len=256) :: s, fname
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
61 integer(kind=kint) :: i,j, rcode, depends, ipt, n
62 real(kind=kreal),
allocatable :: fval(:,:)
63 character(len=HECMW_NAME_LEN) :: data_fmt
66 character(len=256) :: s
71 if( depends>1 ) depends=1
72 if( depends > 3 ) stop
"We cannot read dependencies>3 right now"
83 write(*,*)
"Warning : !ELASTIC : parameter 'INFINITE' is deprecated." &
84 & //
" Please use the replacement parameter 'INFINITESIMAL'"
89 s =
'ISOTROPIC,ORTHOTROPIC,USER '
95 allocate( fval(2+depends,n) )
99 real(kind=kreal) :: fval1(n), fval2(n)
108 if( depends==1 )
then
111 real(kind=kreal) :: fval1(n), fval2(n), fval3(n)
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 )
132 else if( ipt==2 )
then
133 allocate( fval(9+depends,n) )
134 if( depends==0 )
then
135 data_fmt =
"RRRRRRRRR "
137 real(kind=kreal) :: fval1(n), fval2(n), fval3(n), fval4(n), fval5(n), &
138 fval6(n), fval7(n), fval8(n), fval9(n)
150 fval1, fval2, fval3, fval4, fval5, fval6, fval7, fval8, fval9 )
161 else if( depends==1 )
then
162 data_fmt =
"RRRRRRRRRR "
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)
178 fval1, fval2, fval3, fval4, fval5, fval6, fval7, fval8, fval9, fval10 )
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
200 call init_table( mattable, depends, 9+depends,n, fval )
201 call dict_add_key( dict, mc_orthoelastic, mattable )
202 mattype = mn_orthoelastic
206 else if( ipt==3 )
then
208 mattype = userelastic
209 nlgeom = infinitesimal
212 stop
"ERROR: Material type not supported"
216 call finalize_table( mattable )
217 if(
allocated(fval) )
deallocate(fval)
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(:)
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
237 if( depends > 3 ) stop
"We cannot read dependencies>3 right now"
240 if( ipt/=0 ) nlgeom = updatelag
243 s =
'NEOHOOKE,MOONEY-RIVLIN,ARRUDA-BOYCE,USER,MOONEY-RIVLIN-ANISO '
248 allocate( fval(2,depends+1) )
250 if( depends==0 )
then
253 real(kind=kreal) :: fval1(depends+1), fval2(depends+1)
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)
273 else if( ipt==2 )
then
274 allocate( fval(3,depends+1) )
276 if( depends==0 )
then
279 real(kind=kreal) :: fval1(depends+1), fval2(depends+1), fval3(depends+1)
291 matval(m_plconst1) = fval(1,1)
292 matval(m_plconst2) = fval(2,1)
293 matval(m_plconst3) = fval(3,1)
295 mattype = mooneyrivlin
298 else if( ipt==3 )
then
299 allocate( fval(3,depends+1) )
301 if( depends==0 )
then
304 real(kind=kreal) :: fval1(depends+1), fval2(depends+1), fval3(depends+1)
316 matval(m_plconst1) = fval(1,1)
317 matval(m_plconst2) = fval(2,1)
318 matval(m_plconst3) = fval(3,1)
320 mattype = arrudaboyce
322 else if( ipt==4 )
then
324 mattype = userhyperelastic
327 else if( ipt==5 )
then
328 allocate( fval(10,depends+1) )
330 if( depends==0 )
then
331 data_fmt =
"RRRRRrrrrr "
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)
347 fval1, fval2, fval3, fval4, fval5, fval6, fval7, fval8, fval9, fval10 )
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)
372 mattype = mooneyrivlin_aniso
376 if(
allocated(fval) )
deallocate(fval)
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
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
397 if( depends>1 ) depends=1
401 if( ipt/=0 ) nlgeom = infinitesimal
406 write(*,*)
"Warning : !VISCOELASTIC : parameter 'INFINITE' is deprecated." &
407 & //
" Please use the replacement parameter 'INFINITESIMAL'"
408 nlgeom = infinitesimal
412 s =
'ISOTROPIC,USER '
419 allocate( fval(2+depends,n) )
420 if( depends==0 )
then
423 real(kind=kreal) :: fval1(n), fval2(n)
431 if( fval(2,1)==0.d0 ) stop
"Error in defining viscoelasticity: Relaxation time cannot be zero!"
433 if( depends==1 )
then
436 real(kind=kreal) :: fval1(n), fval2(n), fval3(n)
448 call init_table( mattable, 1, 2+depends,n, fval )
449 call dict_add_key( dict, mc_viscoelastic, mattable )
452 mattype = viscoelastic
455 stop
"ERROR: Material type not supported"
459 call finalize_table( mattable )
460 if(
allocated(fval) )
deallocate(fval)
466 integer(kind=kint),
intent(in) :: ctrl
467 integer(kind=kint),
intent(inout) :: mattype
468 real(kind=kreal),
intent(out) :: matval(:)
471 character(len=256) :: s
480 mattype = mattype+ipt
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
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
508 if( depends>1 ) depends = 1
509 if( depends > 3 ) stop
"We cannot read dependencies>3 right now"
513 if( ipt/=0 ) nlgeom = totallag
515 if( ipt/=0 ) nlgeom = infinitesimal
520 write(*,*)
"Warning : !PLASTIC : parameter 'INFINITE' is deprecated." &
521 & //
" Please use the replacement parameter 'INFINITESIMAL'"
522 nlgeom = infinitesimal
525 call setdigit( 1, 1, mattype )
526 call setdigit( 2, 2, mattype )
529 s =
'BILINEAR,MULTILINEAR,SWIFT,RAMBERG-OSGOOD,KINEMATIC,COMBINED '
532 call setdigit( 5, hipt-1, mattype )
535 s =
'MISES,MOHR-COULOMB,DRUCKER-PRAGER,USER '
536 call setdigit( 2, 2, mattype )
539 call setdigit( 4, ipt-1, mattype )
543 if( hipt==2 .and. n<2 )
return
544 if( ( ipt==2 .or. ipt==3 ) .and. hipt>2 ) hipt = 1
550 allocate( fval(2,n) )
553 real(kind=kreal) :: fval1(n), fval2(n)
562 matval(m_plconst1) = fval(1,1)
564 matval(m_plconst2) = fval(2,1)
566 matval(m_plconst2) = 0.d0
567 matval(m_plconst3) = fval(2,1)
571 allocate( fval(depends+2,n) )
572 if( depends==0 )
then
575 real(kind=kreal) :: fval1(n), fval2(n)
584 if( fval(2,1)/=0.d0 )
then
585 print *,
"Multilinear hardening: First plastic strain must be zero"
589 if( fval(2,i)<0.0 ) &
590 stop
"Multilinear hardening: Error in plastic strain definition"
592 call init_table( mttable,1, 2+depends, n, fval )
593 call dict_add_key( dict, mc_yield, mttable )
599 real(kind=kreal) :: fval1(n), fval2(n), fval3(n)
610 call init_table( mttable,2, 2+depends,n, fval )
611 call dict_add_key( dict, mc_yield, mttable )
615 allocate( fval(3,1) )
618 real(kind=kreal) :: fval1(1), fval2(1), fval3(1)
629 matval(m_plconst1) = fval(1,1)
630 matval(m_plconst2) = fval(2,1)
631 matval(m_plconst3) = fval(3,1)
634 print *,
"Error in hardening definition!"
640 allocate( fval(4,n) )
644 real(kind=kreal) :: fval1(n), fval2(n), fval3(n), fval4(n)
657 matval(m_plconst1) = fval(1,1)
658 matval(m_plconst2) = fval(3,1)
659 phi = fval(2,1)*pi/180.d0
660 if( fval(4,1) >= 0.d0 )
then
661 psi = fval(4,1)*pi/180.d0
668 stop
"Mohr-Coulomb and Drucker-Prager do not support temperature dependency"
670 allocate( fval(2,n) )
674 real(kind=kreal) :: fval1(n), fval2(n)
683 phi =fval(1,1)*pi/180.d0
684 if( fval(2,1) >= 0.d0 )
then
685 psi = fval(2,1)*pi/180.d0
689 if( fval(2,2)/=0.d0 )
then
690 print *,
"Multilinear hardening: First plastic strain must be zero"
694 if( fval(2,i)<0.0 ) &
695 stop
"Multilinear hardening: Error in plastic strain definition"
697 call init_table( mttable,1, 2, n-1, fval(1:2,2:n) )
698 call dict_add_key( dict, mc_yield, mttable )
702 matval(m_plconst3) = 2.d0*sin(phi)/ ( sqrt(3.d0)*(3.d0+sin(phi)) )
703 matval(m_plconst4) = 6.d0*cos(phi)/ ( sqrt(3.d0)*(3.d0+sin(phi)) )
704 matval(m_plconst5) = 2.d0*sin(psi)/ ( sqrt(3.d0)*(3.d0+sin(psi)) )
706 matval(m_plconst3) = phi
707 matval(m_plconst4) = psi
714 stop
"Yield function not supported"
717 if(
allocated(fval) )
deallocate(fval)
718 call finalize_table( mttable )
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
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
739 if( depends>1 ) depends=1
742 if( ipt/=0 ) nlgeom = totallag
752 allocate( fval(3+depends,n) )
753 if( depends==0 )
then
756 real(kind=kreal) :: fval1(n), fval2(n), fval3(n)
767 if( depends==1 )
then
770 real(kind=kreal) :: fval1(n), fval2(n), fval3(n), fval4(n)
784 call init_table( mattable, depends, 3+depends,n, fval )
785 call dict_add_key( dict, mc_norton, mattable )
791 stop
"ERROR: Material type not supported"
795 call finalize_table( mattable )
796 if(
allocated(fval) )
deallocate(fval)
802 integer(kind=kint),
intent(in) :: ctrl
803 real(kind=kreal),
intent(out) :: matval(:)
805 integer(kind=kint) :: i, rcode, depends
806 real(kind=kreal),
allocatable :: fval(:,:)
807 character(len=HECMW_NAME_LEN) :: data_fmt
815 if( depends>1 ) depends = 1
817 allocate( fval(1,depends+1) )
819 data_fmt = data_fmt //
"R "
822 real(kind=kreal) :: fval1(depends+1)
830 if(
allocated(fval) )
deallocate(fval)
838 integer(kind=kint),
intent(in) :: ctrl
839 real(kind=kreal),
intent(out) :: matval(:)
840 type(dict_struct),
pointer :: dict
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
853 ss =
'ISOTROPIC,ORTHOTROPIC '
859 if( depends>1 ) depends = 1
862 allocate( fval(depends+1, n) )
864 data_fmt = data_fmt //
"R "
866 if( depends==0 )
then
868 real(kind=kreal) :: fval1(n)
876 real(kind=kreal) :: fval1(n), fval2(n)
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 )
891 allocate( fval(3+depends,n) )
893 data_fmt = trim(data_fmt) //
"R "
895 if( depends==0 )
then
897 real(kind=kreal) :: fval1(n), fval2(n), fval3(n)
907 elseif( depends==1 )
then
909 real(kind=kreal) :: fval1(n), fval2(n), fval3(n), fval4(n)
923 call init_table( mttable, depends, 3+depends,n, fval )
928 call finalize_table( mttable )
929 if(
allocated(fval) )
deallocate(fval)
933 integer function read_user_matl( ctrl, matval )
934 integer(kind=kint),
intent(in) :: ctrl
935 real(kind=kreal),
intent(out) :: matval(:)
937 integer(kind=kint) :: n, i, j
938 real(kind=kreal) :: fval(10,10)
943 if( n > 10 ) stop
"Num of data lines for user-defined material exceeds 10"
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)
959 fval4, fval5, fval6, fval7, fval8, fval9, fval10 ) /= 0 )
return
973 matval((i-1)*10+j)=fval(j,i)
978 end function read_user_matl
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
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
994 character(len=256) :: s
999 if( depends>1 ) depends=1
1000 if( depends > 3 ) stop
"We cannot read dependencies>3 right now"
1004 s =
'INCOMP_NEWTONIAN '
1010 allocate( fval(1+depends,n) )
1011 if( depends==0 )
then
1014 real(kind=kreal) :: fval1(n)
1021 if( depends==1 )
then
1024 real(kind=kreal) :: fval1(n), fval2(n)
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 )
1039 mattype = incomp_newtonian
1042 stop
"ERROR: Material type not supported"
1046 call finalize_table( mattable )
1047 if(
allocated(fval) )
deallocate(fval)
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
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
1065 character(len=DICT_KEY_LENGTH) :: spkey
1071 nlgeom = infinitesimal
1074 allocate( fval(1+depends,n), ival(2,n) )
1075 if( depends==0 )
then
1078 integer(kind=kint) :: ival1(n), ival2(n)
1079 real(kind=kreal) :: fval1(n)
1095 matval_i(m_spring_d_ndoffset+1) = n
1097 matval_i(m_spring_d_ndoffset+2*i ) = ival(1,i)
1098 matval_i(m_spring_d_ndoffset+2*i+1) = ival(2,i)
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 )
1107 if(
allocated(fval) )
deallocate(fval)
1108 if(
allocated(ival) )
deallocate(ival)
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
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
1125 character(len=DICT_KEY_LENGTH) :: spkey
1129 nlgeom = infinitesimal
1132 allocate( fval(1+depends,n) )
1133 if( depends==0 )
then
1136 real(kind=kreal) :: fval1(n)
1148 matval_i(m_spring_a_ndoffset+1) = 1
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 )
1158 if(
allocated(fval) )
deallocate(fval)
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
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
1176 character(len=DICT_KEY_LENGTH) :: spkey
1182 nlgeom = infinitesimal
1185 allocate( fval(1+depends,n), ival(2,n) )
1186 if( depends==0 )
then
1189 integer(kind=kint) :: ival1(n), ival2(n)
1190 real(kind=kreal) :: fval1(n)
1206 matval_i(m_dashpot_d_ndoffset+1) = n
1208 matval_i(m_dashpot_d_ndoffset+2*i ) = ival(1,i)
1209 matval_i(m_dashpot_d_ndoffset+2*i+1) = ival(2,i)
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 )
1218 if(
allocated(fval) )
deallocate(fval)
1219 if(
allocated(ival) )
deallocate(ival)
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
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
1236 character(len=DICT_KEY_LENGTH) :: spkey
1240 nlgeom = infinitesimal
1243 allocate( fval(1+depends,n) )
1244 if( depends==0 )
then
1247 real(kind=kreal) :: fval1(n)
1259 matval_i(m_dashpot_a_ndoffset+1) = 1
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 )
1269 if(
allocated(fval) )
deallocate(fval)
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.
This module provides data structure table which would be dictionaried afterwards.
This module summarizes all information of material properties.
integer(kind=kint), parameter totallag
integer(kind=kint), parameter infinitesimal
integer(kind=kint), parameter usermaterial
integer(kind=kint), parameter updatelag