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

Go to the source code of this file.

Functions/Subroutines

subroutine stiff2mat (elas, ckl, vj, cauchy)
 

Function/Subroutine Documentation

◆ stiff2mat()

subroutine stiff2mat ( real*8, dimension(21)  elas,
real*8, dimension(3,3)  ckl,
real*8  vj,
integer  cauchy 
)
20 !
21 ! elas(21): stiffness constants in the spatial description, i.e.
22 ! the derivative of the Cauchy stress or the Kirchhoff
23 ! stress with respect to the Eulerian strain
24 ! ckl(3,3): inverse deformation gradient
25 ! vj: Jacobian determinant
26 ! cauchy: if 1: elas is written in terms of Cauchy stress
27 ! if 0: elas is written in terms of Kirchhoff stress
28 !
29 ! OUTPUT:
30 !
31 ! elas(21): stiffness constants in the material description,i.e.
32 ! the derivative of the second Piola-Kirchhoff stress (PK2)
33 ! with respect to the Lagrangian strain
34 !
35  implicit none
36 !
37  integer cauchy
38 !
39  integer kk(84),i,nt,k,l,m,n
40 !
41  real*8 elas(21),e(21),ckl(3,3),vj
42 !
43  data 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,
44  & 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,
45  & 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,
46  & 2,3,2,3/
47 !
48  nt=0
49  do i=1,21
50  k=kk(nt+1)
51  l=kk(nt+2)
52  m=kk(nt+3)
53  n=kk(nt+4)
54  nt=nt+4
55  e(i)=elas(1)*ckl(k,1)*ckl(l,1)*ckl(m,1)*ckl(n,1)
56  & +elas(2)*(ckl(k,2)*ckl(l,2)*ckl(m,1)*ckl(n,1)+
57  & ckl(k,1)*ckl(l,1)*ckl(m,2)*ckl(n,2))
58  & +elas(3)*ckl(k,2)*ckl(l,2)*ckl(m,2)*ckl(n,2)
59  & +elas(4)*(ckl(k,3)*ckl(l,3)*ckl(m,1)*ckl(n,1)+
60  & ckl(k,1)*ckl(l,1)*ckl(m,3)*ckl(n,3))
61  & +elas(5)*(ckl(k,3)*ckl(l,3)*ckl(m,2)*ckl(n,2)+
62  & ckl(k,2)*ckl(l,2)*ckl(m,3)*ckl(n,3))
63  & +elas(6)*ckl(k,3)*ckl(l,3)*ckl(m,3)*ckl(n,3)
64  & +elas(7)*(ckl(k,2)*ckl(l,1)*ckl(m,1)*ckl(n,1)+
65  & ckl(k,1)*ckl(l,2)*ckl(m,1)*ckl(n,1)+
66  & ckl(k,1)*ckl(l,1)*ckl(m,2)*ckl(n,1)+
67  & ckl(k,1)*ckl(l,1)*ckl(m,1)*ckl(n,2))
68  & +elas(8)*(ckl(k,2)*ckl(l,2)*ckl(m,2)*ckl(n,1)+
69  & ckl(k,2)*ckl(l,2)*ckl(m,1)*ckl(n,2)+
70  & ckl(k,2)*ckl(l,1)*ckl(m,2)*ckl(n,2)+
71  & ckl(k,1)*ckl(l,2)*ckl(m,2)*ckl(n,2))
72  & +elas(9)*(ckl(k,3)*ckl(l,3)*ckl(m,2)*ckl(n,1)+
73  & ckl(k,3)*ckl(l,3)*ckl(m,1)*ckl(n,2)+
74  & ckl(k,2)*ckl(l,1)*ckl(m,3)*ckl(n,3)+
75  & ckl(k,1)*ckl(l,2)*ckl(m,3)*ckl(n,3))
76  & +elas(10)*(ckl(k,2)*ckl(l,1)*ckl(m,2)*ckl(n,1)+
77  & ckl(k,1)*ckl(l,2)*ckl(m,2)*ckl(n,1)+
78  & ckl(k,2)*ckl(l,1)*ckl(m,1)*ckl(n,2)+
79  & ckl(k,1)*ckl(l,2)*ckl(m,1)*ckl(n,2))
80  & +elas(11)*(ckl(k,3)*ckl(l,1)*ckl(m,1)*ckl(n,1)+
81  & ckl(k,1)*ckl(l,3)*ckl(m,1)*ckl(n,1)+
82  & ckl(k,1)*ckl(l,1)*ckl(m,3)*ckl(n,1)+
83  & ckl(k,1)*ckl(l,1)*ckl(m,1)*ckl(n,3))
84  & +elas(12)*(ckl(k,2)*ckl(l,2)*ckl(m,3)*ckl(n,1)+
85  & ckl(k,3)*ckl(l,1)*ckl(m,2)*ckl(n,2)+
86  & ckl(k,1)*ckl(l,3)*ckl(m,2)*ckl(n,2)+
87  & ckl(k,2)*ckl(l,2)*ckl(m,3)*ckl(n,1))
88  & +elas(13)*(ckl(k,3)*ckl(l,3)*ckl(m,3)*ckl(n,1)+
89  & ckl(k,3)*ckl(l,3)*ckl(m,1)*ckl(n,3)+
90  & ckl(k,3)*ckl(l,1)*ckl(m,3)*ckl(n,3)+
91  & ckl(k,1)*ckl(l,3)*ckl(m,3)*ckl(n,3))
92  & +elas(14)*(ckl(k,3)*ckl(l,1)*ckl(m,2)*ckl(n,1)+
93  & ckl(k,1)*ckl(l,3)*ckl(m,2)*ckl(n,1)+
94  & ckl(k,2)*ckl(l,1)*ckl(m,3)*ckl(n,1)+
95  & ckl(k,1)*ckl(l,2)*ckl(m,3)*ckl(n,1)+
96  & ckl(k,3)*ckl(l,1)*ckl(m,1)*ckl(n,2)+
97  & ckl(k,1)*ckl(l,3)*ckl(m,1)*ckl(n,2)+
98  & ckl(k,2)*ckl(l,1)*ckl(m,1)*ckl(n,3)+
99  & ckl(k,1)*ckl(l,2)*ckl(m,1)*ckl(n,3))
100  & +elas(15)*(ckl(k,3)*ckl(l,1)*ckl(m,3)*ckl(n,1)+
101  & ckl(k,1)*ckl(l,3)*ckl(m,3)*ckl(n,1)+
102  & ckl(k,3)*ckl(l,1)*ckl(m,1)*ckl(n,3)+
103  & ckl(k,1)*ckl(l,3)*ckl(m,1)*ckl(n,3))
104  & +elas(16)*(ckl(k,3)*ckl(l,2)*ckl(m,1)*ckl(n,1)+
105  & ckl(k,2)*ckl(l,3)*ckl(m,1)*ckl(n,1)+
106  & ckl(k,1)*ckl(l,1)*ckl(m,3)*ckl(n,2)+
107  & ckl(k,1)*ckl(l,1)*ckl(m,2)*ckl(n,3))
108  & +elas(17)*(ckl(k,3)*ckl(l,2)*ckl(m,2)*ckl(n,2)+
109  & ckl(k,2)*ckl(l,3)*ckl(m,2)*ckl(n,2)+
110  & ckl(k,2)*ckl(l,2)*ckl(m,3)*ckl(n,2)+
111  & ckl(k,2)*ckl(l,2)*ckl(m,2)*ckl(n,3))
112  & +elas(18)*(ckl(k,3)*ckl(l,3)*ckl(m,3)*ckl(n,2)+
113  & ckl(k,3)*ckl(l,3)*ckl(m,2)*ckl(n,3)+
114  & ckl(k,3)*ckl(l,2)*ckl(m,3)*ckl(n,3)+
115  & ckl(k,2)*ckl(l,3)*ckl(m,3)*ckl(n,3))
116  & +elas(19)*(ckl(k,3)*ckl(l,2)*ckl(m,2)*ckl(n,1)+
117  & ckl(k,2)*ckl(l,3)*ckl(m,2)*ckl(n,1)+
118  & ckl(k,3)*ckl(l,2)*ckl(m,1)*ckl(n,2)+
119  & ckl(k,2)*ckl(l,3)*ckl(m,1)*ckl(n,2)+
120  & ckl(k,2)*ckl(l,1)*ckl(m,3)*ckl(n,2)+
121  & ckl(k,1)*ckl(l,2)*ckl(m,3)*ckl(n,2)+
122  & ckl(k,2)*ckl(l,1)*ckl(m,2)*ckl(n,3)+
123  & ckl(k,1)*ckl(l,2)*ckl(m,2)*ckl(n,3))
124  & +elas(20)*(ckl(k,3)*ckl(l,2)*ckl(m,3)*ckl(n,1)+
125  & ckl(k,2)*ckl(l,3)*ckl(m,3)*ckl(n,1)+
126  & ckl(k,3)*ckl(l,1)*ckl(m,3)*ckl(n,2)+
127  & ckl(k,1)*ckl(l,3)*ckl(m,3)*ckl(n,2)+
128  & ckl(k,3)*ckl(l,2)*ckl(m,1)*ckl(n,3)+
129  & ckl(k,2)*ckl(l,3)*ckl(m,1)*ckl(n,3)+
130  & ckl(k,3)*ckl(l,1)*ckl(m,2)*ckl(n,3)+
131  & ckl(k,1)*ckl(l,3)*ckl(m,2)*ckl(n,3))
132  & +elas(21)*(ckl(k,3)*ckl(l,2)*ckl(m,3)*ckl(n,2)+
133  & ckl(k,2)*ckl(l,3)*ckl(m,3)*ckl(n,2)+
134  & ckl(k,3)*ckl(l,2)*ckl(m,2)*ckl(n,3)+
135  & ckl(k,2)*ckl(l,3)*ckl(m,2)*ckl(n,3))
136  enddo
137 !
138  if(cauchy.eq.1) then
139  do i=1,21
140  elas(i)=e(i)*vj
141  enddo
142  else
143  do i=1,21
144  elas(i)=e(i)
145  enddo
146  endif
147 !
148  return
subroutine cauchy(n, x, l, u, nbd, g, iorder, iwhere, t, d, xcp, m, wy, ws, sy, wt, theta, col, head, p, c, wbp, v, nseg, iprint, sbgnrm, info, epsmch)
Definition: lbfgsb.f:1225
Hosted by OpenAircraft.com, (Michigan UAV, LLC)