29 subroutine matlmatrix( gauss, sectType, matrix, time, dtime, cdsys, temperature, isEp, hdflag )
31 integer,
intent(in) :: sectType
32 real(kind=kreal),
intent(out) :: matrix(:,:)
33 real(kind=kreal),
intent(in) :: time
34 real(kind=kreal),
intent(in) :: dtime
35 real(kind=kreal),
intent(in) :: cdsys(3,3)
36 real(kind=kreal),
intent(in) :: temperature
37 integer(kind=kint),
intent(in),
optional :: isEp
38 integer(kind=kint),
intent(in),
optional :: hdflag
41 integer :: flag, hdflag_in
42 real(kind=kreal) :: cijkl(3,3,3,3)
47 if(
present(isep) )
then
48 if( isep == 1 )flag = 1
52 if(
present(hdflag) ) hdflag_in = hdflag
58 elseif(
iselastic(matl%mtype) .or. flag==1 )
then
66 call calelasticmatrix( matl, secttype, matrix, temperature, hdflag=hdflag_in )
70 print *,
"Elasticity type", matl%mtype,
"not supported"
76 call mat_c2d( cijkl, matrix, secttype )
79 call mat_c2d( cijkl, matrix, secttype )
82 call mat_c2d( cijkl, matrix, secttype )
87 gauss%istatus(1), gauss%fstatus, gauss%plstrain, matrix, temperature, hdflag=hdflag_in )
89 call umatlmatrix( matl%name, matl%variables(101:), gauss%strain, &
90 gauss%stress, gauss%fstatus, matrix, dtime, time )
91 elseif( matl%mtype==
norton )
then
92 call iso_creep( matl, secttype, gauss%stress, gauss%strain, gauss%fstatus, &
93 gauss%plstrain, dtime, time, matrix, temperature, hdflag=hdflag_in )
95 stop
"Material type not supported!"
102 subroutine stressupdate( gauss, sectType, strain, stress, cdsys, time, dtime, temp, tempn, hdflag )
104 integer,
intent(in) :: sectType
105 real(kind=kreal),
intent(in) :: strain(6)
106 real(kind=kreal),
intent(out) :: stress(6)
107 real(kind=kreal),
intent(in) :: cdsys(3,3)
108 real(kind=kreal),
intent(in),
optional :: time
109 real(kind=kreal),
intent(in),
optional :: dtime
110 real(kind=kreal),
intent(in) :: temp
111 real(kind=kreal),
intent(in) :: tempn
112 integer(kind=kint),
intent(in) :: hdflag
121 call uelasticupdate( gauss%pMaterial%variables(101:), strain, stress )
123 if( .not.
present(dtime) ) stop
"error in viscoelastic update!"
124 call updateviscoelastic( gauss%pMaterial, secttype, strain, stress, gauss%fstatus, dtime, temp, tempn )
125 elseif ( gauss%pMaterial%mtype==
norton )
then
126 if( .not.
present(dtime) ) stop
"error in viscoelastic update!"
127 call update_iso_creep( gauss%pMaterial, secttype, strain, stress, gauss%fstatus, &
128 & gauss%plstrain, dtime, time, temp, hdflag=hdflag )
130 call uupdate( gauss%pMaterial%name, gauss%pMaterial%variables(101:), &
131 strain, stress, gauss%fstatus, dtime, time )
138 real(kind=kreal),
intent(in) :: cijkl(3,3,3,3)
139 real(kind=kreal),
intent(out) :: dij(6,6)
140 integer,
intent(in) :: itype
145 dij(1,1) = cijkl(1,1,1,1)
146 dij(1,2) = cijkl(1,1,2,2)
147 dij(1,3) = cijkl(1,1,3,3)
148 dij(1,4) = cijkl(1,1,1,2)
149 dij(1,5) = cijkl(1,1,2,3)
150 dij(1,6) = cijkl(1,1,3,1)
151 dij(2,1) = cijkl(2,2,1,1)
152 dij(2,2) = cijkl(2,2,2,2)
153 dij(2,3) = cijkl(2,2,3,3)
154 dij(2,4) = cijkl(2,2,1,2)
155 dij(2,5) = cijkl(2,2,2,3)
156 dij(2,6) = cijkl(2,2,3,1)
157 dij(3,1) = cijkl(3,3,1,1)
158 dij(3,2) = cijkl(3,3,2,2)
159 dij(3,3) = cijkl(3,3,3,3)
160 dij(3,4) = cijkl(3,3,1,2)
161 dij(3,5) = cijkl(3,3,2,3)
162 dij(3,6) = cijkl(3,3,3,1)
163 dij(4,1) = cijkl(1,2,1,1)
164 dij(4,2) = cijkl(1,2,2,2)
165 dij(4,3) = cijkl(1,2,3,3)
166 dij(4,4) = cijkl(1,2,1,2)
167 dij(4,5) = cijkl(1,2,2,3)
168 dij(4,6) = cijkl(1,2,3,1)
169 dij(5,1) = cijkl(2,3,1,1)
170 dij(5,2) = cijkl(2,3,2,2)
171 dij(5,3) = cijkl(2,3,3,3)
172 dij(5,4) = cijkl(2,3,1,2)
173 dij(5,5) = cijkl(2,3,2,3)
174 dij(5,6) = cijkl(2,3,3,1)
175 dij(6,1) = cijkl(3,1,1,1)
176 dij(6,2) = cijkl(3,1,2,2)
177 dij(6,3) = cijkl(3,1,3,3)
178 dij(6,4) = cijkl(3,1,1,2)
179 dij(6,5) = cijkl(3,1,2,3)
180 dij(6,6) = cijkl(3,1,3,1)
183 dij(1,1) = cijkl(1,1,1,1)
184 dij(1,2) = cijkl(1,1,2,2)
185 dij(1,3) = cijkl(1,1,1,2)
186 dij(2,1) = cijkl(2,2,1,1)
187 dij(2,2) = cijkl(2,2,2,2)
188 dij(2,3) = cijkl(2,2,1,2)
189 dij(3,1) = cijkl(1,2,1,1)
190 dij(3,2) = cijkl(1,2,2,2)
191 dij(3,3) = cijkl(1,2,1,2)
193 dij(1,1) = cijkl(1,1,1,1)
194 dij(1,2) = cijkl(1,1,2,2)
195 dij(1,3) = cijkl(1,1,1,2)
196 dij(1,4) = cijkl(1,1,3,3)
197 dij(2,1) = cijkl(2,2,1,1)
198 dij(2,2) = cijkl(2,2,2,2)
199 dij(2,3) = cijkl(2,2,1,2)
200 dij(2,4) = cijkl(2,2,3,3)
201 dij(3,1) = cijkl(1,2,1,1)
202 dij(3,2) = cijkl(1,2,2,2)
203 dij(3,3) = cijkl(1,2,1,2)
204 dij(3,4) = cijkl(1,2,3,3)
205 dij(4,1) = cijkl(3,3,1,1)
206 dij(4,2) = cijkl(3,3,2,2)
207 dij(4,3) = cijkl(3,3,1,2)
208 dij(4,4) = cijkl(3,3,3,3)
217 (gauss, secttype, d, &
218 e1_hat, e2_hat, e3_hat, cg1, cg2, cg3, &
223 integer,
intent(in) :: sectType, n_layer
224 real(kind = kreal),
intent(out) :: d(:, :)
225 real(kind = kreal),
intent(in) :: e1_hat(3), e2_hat(3), e3_hat(3)
226 real(kind = kreal),
intent(in) :: cg1(3), cg2(3), cg3(3)
227 real(kind = kreal),
intent(out) :: alpha
231 real(kind = kreal) :: c(3, 3, 3, 3)
236 matl => gauss%pMaterial
242 (matl, secttype, c, &
243 e1_hat, e2_hat, e3_hat, cg1, cg2, cg3, &
248 stop
"Material type not supported!"
266 real(kind = kreal),
intent(in) :: c(:, :, :, :)
267 real(kind = kreal),
intent(out) :: d(:, :)
268 integer,
intent(in) :: itype
272 integer :: index_i(5), index_j(5), &
273 index_k(5), index_l(5)
274 integer :: i, j, k, l
321 d(is, js) = c(i, j, k, l)
340 integer(kind=kint),
intent(in) :: ctype
341 integer(kind=kint),
intent(in) :: dofid
342 real(kind=kreal),
intent(out) :: rparams(:)
343 integer(kind=kint),
intent(out) :: iparams(:)
344 real(kind=kreal),
intent(in),
optional :: stretch
346 real(kind=kreal) ina(1), outa(4)
348 character(len=DICT_KEY_LENGTH) :: cnkey
354 write(cnkey,
'(A,I0)') trim(
mc_spring),dofid
358 write(cnkey,
'(A,A)') trim(
mc_spring),
'_A'
366 stop
"CONNECTOR ctype is not defined"
369 if(
present(stretch) )
then
371 call fetch_tabledata( cnkey, gauss%pMaterial%dict, outa(1:1), ierr, ina )
373 call fetch_tabledata( cnkey, gauss%pMaterial%dict, outa(1:1), ierr )
This module provides functions for elastic material.
subroutine calelasticmatrix(matl, sectType, D, temp, hdflag)
Calculate isotropic elastic matrix.
subroutine linearelastic_shell(matl, sectType, c, e1_hat, e2_hat, e3_hat, cg1, cg2, cg3, alpha, n_layer)
subroutine calelasticmatrix_ortho(matl, sectType, bij, DMAT, temp)
Calculate orthotropic elastic matrix.
This module provide functions for elastoplastic calculation.
subroutine, public calelastoplasticmatrix(matl, sectType, stress, istat, extval, plstrain, D, temperature, hdflag)
This subroutine calculates elastoplastic constitutive relation.
This module manages calculation relates with materials.
subroutine getconnectorproperty(gauss, ctype, dofid, rparams, iparams, stretch)
subroutine matlmatrix(gauss, sectType, matrix, time, dtime, cdsys, temperature, isEp, hdflag)
Calculate constituive matrix.
subroutine mat_c2d_shell(c, D, itype)
integer function getnlgeomflag(gauss)
Fetch the nlgeom flag of the material.
subroutine matlmatrix_shell(gauss, sectType, D, e1_hat, e2_hat, e3_hat, cg1, cg2, cg3, alpha, n_layer)
subroutine mat_c2d(cijkl, dij, itype)
Transfer rank 4 constituive matrix to rank 2 form.
subroutine stressupdate(gauss, sectType, strain, stress, cdsys, time, dtime, temp, tempn, hdflag)
Update strain and stress for elastic and hyperelastic materials.
This module provides functions for creep calculation.
subroutine update_iso_creep(matl, sectType, strain, stress, extval, plstrain, dtime, ttime, temp, hdflag)
This subroutine calculates stresses and creep status for an elastically isotropic material with isotr...
subroutine iso_creep(matl, sectType, stress, strain, extval, plstrain, dtime, ttime, stiffness, temp, hdflag)
This subroutine calculates stiffness for elastically isotropic materials with isotropic creep.
This module provides functions for hyperelastic calculation.
subroutine calupdateelasticarrudaboyce(matl, sectType, dstrain, dstress)
This subroutine provides to update stress and strain for Arrude-Royce material.
subroutine calelasticarrudaboyce(matl, sectType, cijkl, strain)
This subroutine provides elastic tangent coefficient for Arruda-Boyce hyperelastic material.
subroutine calelasticmooneyrivlinaniso(matl, sectType, cijkl, strain, cdsys, hdflag)
This subroutine provides elastic tangent coefficient for anisotropic Mooney-Rivlin hyperelastic mater...
subroutine calupdateelasticmooneyrivlinaniso(matl, sectType, strain, stress, cdsys, hdflag)
This subroutine provides to update stress and strain for anisotropic Mooney-Rivlin material.
subroutine calelasticmooneyrivlin(matl, sectType, cijkl, strain, hdflag)
This subroutine provides elastic tangent coefficient for Mooney-Rivlin hyperelastic material.
subroutine calupdateelasticmooneyrivlin(matl, sectType, strain, stress, strain_energy, hdflag)
This subroutine provides to update stress and strain for Mooney-Rivlin material.
This module summarizes all information of material properties.
integer(kind=kint), parameter m_dashpot_axial
integer(kind=kint), parameter mooneyrivlin
integer(kind=kint), parameter mooneyrivlin_aniso
integer(kind=kint), parameter planestress
integer(kind=kint), parameter m_spring_axial
integer function getelastictype(mtype)
Get elastic type.
character(len=dict_key_length) mc_spring
integer(kind=kint), parameter d3
integer(kind=kint), parameter planestrain
integer(kind=kint), parameter elastic
integer(kind=kint), parameter arrudaboyce
integer(kind=kint), parameter shell
integer(kind=kint), parameter norton
character(len=dict_key_length) mc_dashpot
integer(kind=kint), parameter axissymetric
integer(kind=kint), parameter m_dashpot_dof
integer(kind=kint), parameter userelastic
integer(kind=kint), parameter neohooke
integer(kind=kint), parameter usermaterial
integer(kind=kint), parameter m_spring_dof
logical function iselastic(mtype)
If it is an elastic material?
integer(kind=kint), parameter userhyperelastic
integer(kind=kint), parameter m_spring_d_ndoffset
integer(kind=kint), parameter m_dashpot_d_ndoffset
logical function isviscoelastic(mtype)
If it is an viscoelastic material?
logical function iselastoplastic(mtype)
If it is an elastoplastic material?
This modules defines a structure to record history dependent parameter in static analysis.
This module provides interface for elastic or hyperelastic calculation.
subroutine uelasticmatrix(matl, strain, D)
This subroutine calculates constitutive relation.
subroutine uelasticupdate(matl, strain, stress)
This subroutine calculate updated strain and stress.
This subroutine read in used-defined material properties tangent.
subroutine uupdate(mname, matl, strain, stress, fstat, dtime, ttime, temperature)
This subroutine calculate strain and stress increment.
subroutine umatlmatrix(mname, matl, strain, stress, fstat, D, dtime, ttime, temperature)
This subroutine calculates constitutive matrix.
This module provides functions for viscoelastic calculation.
subroutine updateviscoelastic(matl, sectType, eps, sig, vsig, dt, temp, tempn)
This subroutine provides to update stress for viscoelastic material.
subroutine calviscoelasticmatrix(matl, sectType, dt, D, temp)
This subroutine calculates tangent moduli for isotropic viscoelastic material.
Structure to manage all material related data.
All data should be recorded in every quadrature points.