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

Go to the source code of this file.

Functions/Subroutines

subroutine str2mat (str, ckl, vj, cauchy)
 

Function/Subroutine Documentation

◆ str2mat()

subroutine str2mat ( real*8, dimension(6)  str,
real*8, dimension(3,3)  ckl,
real*8  vj,
integer  cauchy 
)
20 !
21 ! converts the stress in spatial coordinates into material coordinates
22 ! or the strain in material coordinates into spatial coordinates.
23 !
24 ! INPUT:
25 !
26 ! str(6): Cauchy stress, Kirchhoff stress or Lagrange strain
27 ! component order: 11,22,33,12,13,23
28 ! ckl(3,3): the inverse deformation gradient
29 ! vj: Jakobian determinant
30 ! cauchy: integer variable
31 ! if 1: str contains the Cauchy stress
32 ! if 0: str contains the Kirchhoff stress or
33 ! Lagrange strain
34 !
35 ! OUTPUT:
36 !
37 ! str(6): Piola-Kirchhoff stress of the second kind (PK2) or
38 ! Euler strain
39 !
40  implicit none
41 !
42  integer cauchy
43 !
44  integer i,m1,m2
45 !
46  real*8 str(6),s(6),ckl(3,3),vj
47 !
48  do i=1,6
49  if(i.eq.1) then
50  m1=1
51  m2=1
52  elseif(i.eq.2) then
53  m1=2
54  m2=2
55  elseif(i.eq.3) then
56  m1=3
57  m2=3
58  elseif(i.eq.4) then
59  m1=2
60  m2=1
61  elseif(i.eq.5) then
62  m1=3
63  m2=1
64  else
65  m1=3
66  m2=2
67  endif
68 !
69  s(i)=(str(1)*ckl(m1,1)*ckl(m2,1)+
70  & str(2)*ckl(m1,2)*ckl(m2,2)+
71  & str(3)*ckl(m1,3)*ckl(m2,3)+
72  & str(4)*(ckl(m1,1)*ckl(m2,2)+ckl(m1,2)*ckl(m2,1))+
73  & str(5)*(ckl(m1,1)*ckl(m2,3)+ckl(m1,3)*ckl(m2,1))+
74  & str(6)*(ckl(m1,2)*ckl(m2,3)+ckl(m1,3)*ckl(m2,2)))
75 !
76  enddo
77 !
78  if(cauchy.eq.1) then
79  do i=1,6
80  str(i)=s(i)*vj
81  enddo
82  else
83  do i=1,6
84  str(i)=s(i)
85  enddo
86  endif
87 !
88  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)