CalculiX  2.13
A Free Software Three-Dimensional Structural Finite Element Program
umat_elastic_fiber.f File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine umat_elastic_fiber (amat, iel, iint, kode, elconloc, emec, emec0, beta, xokl, voj, xkl, vj, ithermal, t1l, dtime, time, ttime, icmd, ielas, mi, nstate_, xstateini, xstate, stre, stiff, iorien, pgauss, orab)
 

Function/Subroutine Documentation

◆ umat_elastic_fiber()

subroutine umat_elastic_fiber ( character*80  amat,
integer  iel,
integer  iint,
integer  kode,
real*8, dimension(21)  elconloc,
real*8, dimension(6)  emec,
real*8, dimension(6)  emec0,
real*8, dimension(6)  beta,
real*8, dimension(3,3)  xokl,
real*8  voj,
real*8, dimension(3,3)  xkl,
real*8  vj,
integer  ithermal,
real*8  t1l,
real*8  dtime,
real*8  time,
real*8  ttime,
integer  icmd,
integer  ielas,
integer, dimension(*)  mi,
integer  nstate_,
real*8, dimension(nstate_,mi(1),*)  xstateini,
real*8, dimension(nstate_,mi(1),*)  xstate,
real*8, dimension(6)  stre,
real*8, dimension(21)  stiff,
integer  iorien,
real*8, dimension(3)  pgauss,
real*8, dimension(7,*)  orab 
)
24 !
25 ! calculates stiffness and stresses for a user defined material
26 ! law
27 !
28 ! icmd=3: calcutates stress at mechanical strain
29 ! else: calculates stress at mechanical strain and the stiffness
30 ! matrix
31 !
32 ! INPUT:
33 !
34 ! amat material name
35 ! iel element number
36 ! iint integration point number
37 !
38 ! kode material type (-100-#of constants entered
39 ! under *USER MATERIAL): can be used for materials
40 ! with varying number of constants
41 !
42 ! elconloc(21) user defined constants defined by the keyword
43 ! card *USER MATERIAL (max. 21, actual # =
44 ! -kode-100), interpolated for the
45 ! actual temperature t1l
46 !
47 ! emec(6) Lagrange mechanical strain tensor (component order:
48 ! 11,22,33,12,13,23) at the end of the increment
49 ! (thermal strains are subtracted)
50 ! emec0(6) Lagrange mechanical strain tensor at the start of the
51 ! increment (thermal strains are subtracted)
52 ! beta(6) residual stress tensor (the stress entered under
53 ! the keyword *INITIAL CONDITIONS,TYPE=STRESS)
54 !
55 ! xokl(3,3) deformation gradient at the start of the increment
56 ! voj Jacobian at the start of the increment
57 ! xkl(3,3) deformation gradient at the end of the increment
58 ! vj Jacobian at the end of the increment
59 !
60 ! ithermal 0: no thermal effects are taken into account
61 ! >0: thermal effects are taken into account (triggered
62 ! by the keyword *INITIAL CONDITIONS,TYPE=TEMPERATURE)
63 ! t1l temperature at the end of the increment
64 ! dtime time length of the increment
65 ! time step time at the end of the current increment
66 ! ttime total time at the start of the current step
67 !
68 ! icmd not equal to 3: calculate stress and stiffness
69 ! 3: calculate only stress
70 ! ielas 0: no elastic iteration: irreversible effects
71 ! are allowed
72 ! 1: elastic iteration, i.e. no irreversible
73 ! deformation allowed
74 !
75 ! mi(1) max. # of integration points per element in the
76 ! model
77 ! nstate_ max. # of state variables in the model
78 !
79 ! xstateini(nstate_,mi(1),# of elements)
80 ! state variables at the start of the increment
81 ! xstate(nstate_,mi(1),# of elements)
82 ! state variables at the end of the increment
83 !
84 ! stre(6) Piola-Kirchhoff stress of the second kind
85 ! at the start of the increment
86 !
87 ! iorien number of the local coordinate axis system
88 ! in the integration point at stake (takes the value
89 ! 0 if no local system applies)
90 ! pgauss(3) global coordinates of the integration point
91 ! orab(7,*) description of all local coordinate systems.
92 ! If a local coordinate system applies the global
93 ! tensors can be obtained by premultiplying the local
94 ! tensors with skl(3,3). skl is determined by calling
95 ! the subroutine transformatrix:
96 ! call transformatrix(orab(1,iorien),pgauss,skl)
97 !
98 !
99 ! OUTPUT:
100 !
101 ! xstate(nstate_,mi(1),# of elements)
102 ! updated state variables at the end of the increment
103 ! stre(6) Piola-Kirchhoff stress of the second kind at the
104 ! end of the increment
105 ! stiff(21): consistent tangent stiffness matrix in the material
106 ! frame of reference at the end of the increment. In
107 ! other words: the derivative of the PK2 stress with
108 ! respect to the Lagrangian strain tensor. The matrix
109 ! is supposed to be symmetric, only the upper half is
110 ! to be given in the same order as for a fully
111 ! anisotropic elastic material (*ELASTIC,TYPE=ANISO).
112 ! Notice that the matrix is an integral part of the
113 ! fourth order material tensor, i.e. the Voigt notation
114 ! is not used.
115 !
116  implicit none
117 !
118  character*80 amat
119 !
120  integer ithermal,icmd,kode,ielas,iel,iint,nstate_,mi(*),nfiber,
121  & i,
122  & j,k,l,m,n,ioffset,nt,kk(84),iorien
123 !
124  real*8 elconloc(21),stiff(21),emec0(6),beta(6),stre(6),
125  & vj,t1l,dtime,xkl(3,3),xokl(3,3),voj,c(3,3),a(3),pgauss(3),
126  & orab(7,*),skl(3,3),aa(3),emec(6),time,ttime
127 !
128  real*8 xstate(nstate_,mi(1),*),xstateini(nstate_,mi(1),*),
129  & constant(21),dd,dm(3,3,4),djdc(3,3,4),d2jdc2(3,3,3,3,4),
130  & v1,v1b,v3,v3bi,v4(4),v4br(4),djbdc(3,3,4),d2jbdc2(3,3,3,3,4),
131  & didc(3,3,3),d2idc2(3,3,3,3,3),dibdc(3,3,3),d2ibdc2(3,3,3,3,3),
132  & dudc(3,3),d2udc2(3,3,3,3),v33,cinv(3,3),xk1,xk2,d(3,3),term
133 !
134  kk=(/1,1,1,1,1,1,2,2,2,2,2,2,1,1,3,3,2,2,3,3,3,3,3,3,
135  & 1,1,1,2,2,2,1,2,3,3,1,2,1,2,1,2,1,1,1,3,2,2,1,3,3,3,1,3,
136  & 1,2,1,3,1,3,1,3,1,1,2,3,2,2,2,3,3,3,2,3,1,2,2,3,1,3,2,3,
137  & 2,3,2,3/)
138 !
139 ! calculating the transformation matrix
140 !
141  if(iorien.gt.0) then
142  call transformatrix(orab(1,iorien),pgauss,skl)
143  endif
144 !
145 ! # of fibers
146 !
147  nfiber=(-kode-102)/4
148  do i=1,-kode-100
149  constant(i)=elconloc(i)
150  enddo
151  if(dabs(constant(2)).lt.1.d-10) then
152  constant(2)=1.d0/(20.d0*constant(1))
153  endif
154 !
155 ! calculation of the Green deformation tensor for the
156 ! mechanical strain
157 !
158  do i=1,3
159  c(i,i)=emec(i)*2.d0+1.d0
160  enddo
161  c(1,2)=2.d0*emec(4)
162  c(1,3)=2.d0*emec(5)
163  c(2,3)=2.d0*emec(6)
164 !
165 ! creation of the delta Dirac matrix d
166 !
167  do i=1,3
168  d(i,i)=1.d0
169  enddo
170  d(1,2)=0.d0
171  d(1,3)=0.d0
172  d(2,3)=0.d0
173 !
174 ! calculation of the structural tensors
175 !
176  do k=1,nfiber
177  ioffset=4*k-1
178  a(1)=constant(ioffset)
179  a(2)=constant(ioffset+1)
180  dd=a(1)*a(1)+a(2)*a(2)
181  if(dd.gt.1.d0) then
182  write(*,*) '*ERROR in umat_el_fiber: components of'
183  write(*,*) ' direction vector ',k,' are too big'
184  call exit(201)
185  endif
186  a(3)=dsqrt(1.d0-dd)
187 !
188 ! check for local coordinate systems
189 !
190  if(iorien.gt.0) then
191  do j=1,3
192  aa(j)=a(j)
193  enddo
194  do j=1,3
195  a(j)=skl(j,1)*aa(1)+skl(j,2)*aa(2)+skl(j,3)*aa(3)
196  enddo
197  endif
198 !
199  do j=1,3
200  do i=1,j
201  dm(i,j,k)=a(i)*a(j)
202  enddo
203  enddo
204  enddo
205 !
206 ! calculation of the invariants
207 !
208  v1=c(1,1)+c(2,2)+c(3,3)
209  v3=c(1,1)*(c(2,2)*c(3,3)-c(2,3)*c(2,3))
210  & -c(1,2)*(c(1,2)*c(3,3)-c(1,3)*c(2,3))
211  & +c(1,3)*(c(1,2)*c(2,3)-c(1,3)*c(2,2))
212  do j=1,nfiber
213  v4(j)=dm(1,1,j)*c(1,1)+dm(2,2,j)*c(2,2)+dm(3,3,j)*c(3,3)+
214  & 2.d0*(dm(1,2,j)*c(1,2)+dm(1,3,j)*c(1,3)+dm(2,3,j)*c(2,3))
215  enddo
216 !
217  v33=v3**(-1.d0/3.d0)
218 !
219 ! inversion of c
220 !
221  cinv(1,1)=(c(2,2)*c(3,3)-c(2,3)*c(2,3))/v3
222  cinv(2,2)=(c(1,1)*c(3,3)-c(1,3)*c(1,3))/v3
223  cinv(3,3)=(c(1,1)*c(2,2)-c(1,2)*c(1,2))/v3
224  cinv(1,2)=(c(1,3)*c(2,3)-c(1,2)*c(3,3))/v3
225  cinv(1,3)=(c(1,2)*c(2,3)-c(2,2)*c(1,3))/v3
226  cinv(2,3)=(c(1,2)*c(1,3)-c(1,1)*c(2,3))/v3
227  cinv(2,1)=cinv(1,2)
228  cinv(3,1)=cinv(1,3)
229  cinv(3,2)=cinv(2,3)
230 !
231 ! first derivative of the invariants with respect to c(k,l)
232 !
233  do l=1,3
234  do k=1,l
235  didc(k,l,1)=d(k,l)
236  didc(k,l,3)=v3*cinv(k,l)
237  do j=1,nfiber
238  djdc(k,l,j)=dm(k,l,j)
239  enddo
240  enddo
241  enddo
242 !
243 ! second derivative of the invariants w.r.t. c(k,l)
244 ! and c(m,n)
245 !
246  if(icmd.ne.3) then
247  nt=0
248  do i=1,21
249  k=kk(nt+1)
250  l=kk(nt+2)
251  m=kk(nt+3)
252  n=kk(nt+4)
253  nt=nt+4
254  d2idc2(k,l,m,n,1)=0.d0
255  d2idc2(k,l,m,n,3)=v3*(cinv(m,n)*cinv(k,l)-
256  & (cinv(k,m)*cinv(n,l)+cinv(k,n)*cinv(m,l))/2.d0)
257  do j=1,nfiber
258  d2jdc2(k,l,m,n,j)=0.d0
259  enddo
260  enddo
261  endif
262 !
263 ! derivatives for the reduced invariants
264 !
265  v1b=v1*v33
266  v3bi=1.d0/dsqrt(v3)
267  do j=1,nfiber
268  v4br(j)=v4(j)*v33-1.d0
269  enddo
270 !
271 ! first derivative of the reduced c-invariants w.r.t. c(k,l)
272 !
273  do l=1,3
274  do k=1,l
275  dibdc(k,l,1)=-v33**4*v1*didc(k,l,3)/3.d0
276  & +v33*didc(k,l,1)
277  do j=1,nfiber
278  djbdc(k,l,j)=-v33**4*v4(j)*didc(k,l,3)/3.d0
279  & +v33*djdc(k,l,j)
280  enddo
281  enddo
282  enddo
283 !
284 ! second derivative of the reduced c-invariants w.r.t. c(k,l)
285 ! and c(m,n)
286 !
287  if(icmd.ne.3) then
288  nt=0
289  do i=1,21
290  k=kk(nt+1)
291  l=kk(nt+2)
292  m=kk(nt+3)
293  n=kk(nt+4)
294  nt=nt+4
295  d2ibdc2(k,l,m,n,1)=4.d0/9.d0*v33**7*v1*didc(k,l,3)
296  & *didc(m,n,3)-v33**4/3.d0*(didc(m,n,1)*didc(k,l,3)
297  & +didc(k,l,1)*didc(m,n,3))-v33**4/3.d0*v1*
298  & d2idc2(k,l,m,n,3)+v33*d2idc2(k,l,m,n,1)
299  do j=1,nfiber
300  d2jbdc2(k,l,m,n,j)=4.d0/9.d0*v33**7*v4(j)*didc(k,l,3)
301  & *didc(m,n,3)-v33**4/3.d0*(djdc(m,n,j)*didc(k,l,3)
302  & +djdc(k,l,j)*didc(m,n,3))-v33**4/3.d0*v4(j)*
303  & d2idc2(k,l,m,n,3)+v33*d2jdc2(k,l,m,n,j)
304  enddo
305  enddo
306  endif
307 !
308 ! calculation of the stress
309 ! the anisotropy is only taken into account for v4br(j)>=0
310 !
311  do l=1,3
312  do k=1,l
313  dudc(k,l)=constant(1)*dibdc(k,l,1)+
314  & (1.d0-v3bi)*didc(k,l,3)/constant(2)
315  do j=1,nfiber
316  if(v4br(j).lt.0.d0) cycle
317  if(xk2*v4br(j)**2.gt.227.d0) then
318  write(*,*) '*ERROR in umat_elastic_fiber'
319  write(*,*) ' fiber extension is too large'
320  write(*,*) ' for exponential function'
321  call exit(201)
322  endif
323  ioffset=4*j
324  xk1=constant(ioffset+1)
325  xk2=constant(ioffset+2)
326  dudc(k,l)=dudc(k,l)+xk1*v4br(j)*
327  & dexp(xk2*v4br(j)**2)*djbdc(k,l,j)
328  enddo
329  enddo
330  enddo
331 !
332 ! calculation of the stiffness matrix
333 ! the anisotropy is only taken into account for v4br(j)>=0
334 !
335  if(icmd.ne.3) then
336  nt=0
337  do i=1,21
338  k=kk(nt+1)
339  l=kk(nt+2)
340  m=kk(nt+3)
341  n=kk(nt+4)
342  nt=nt+4
343  term=constant(1)*d2ibdc2(k,l,m,n,1)+
344  & v3bi**3*didc(k,l,3)*didc(m,n,3)/(2.d0*constant(2))
345  & +(1.d0-v3bi)*d2idc2(k,l,m,n,3)/constant(2)
346  do j=1,nfiber
347  if(v4br(j).lt.0.d0) cycle
348  ioffset=4*j
349  xk1=constant(ioffset+1)
350  xk2=constant(ioffset+2)
351  term=term+xk1*dexp(xk2*v4br(j)**2)*
352  & (djbdc(k,l,j)*djbdc(m,n,j)*(1.d0+2.d0*xk2*v4br(j)**2)+
353  & v4br(j)*d2jbdc2(k,l,m,n,j))
354  enddo
355  d2udc2(k,l,m,n)=term
356  enddo
357  endif
358 !
359 ! storing the stiffness matrix and/or the stress
360 !
361  if(icmd.ne.3) then
362 !
363 ! storing the stiffness matrix
364 !
365  nt=0
366  do i=1,21
367  k=kk(nt+1)
368  l=kk(nt+2)
369  m=kk(nt+3)
370  n=kk(nt+4)
371  nt=nt+4
372  stiff(i)=4.d0*d2udc2(k,l,m,n)
373  enddo
374  endif
375 !
376 ! store the stress at mechanical strain
377 !
378  stre(1)=2.d0*dudc(1,1)
379  stre(2)=2.d0*dudc(2,2)
380  stre(3)=2.d0*dudc(3,3)
381  stre(4)=2.d0*dudc(1,2)
382  stre(5)=2.d0*dudc(1,3)
383  stre(6)=2.d0*dudc(2,3)
384 !
385  return
subroutine transformatrix(xab, p, a)
Definition: transformatrix.f:20
static double * v1
Definition: mafillsmmain_se.c:40
Hosted by OpenAircraft.com, (Michigan UAV, LLC)