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

Go to the source code of this file.

Functions/Subroutines

subroutine materialdata_sp (elcon, nelcon, imat, ntmat_, i, t1l, elconloc, kode, plicon, nplicon, npmat_, plconloc, ncmat_)
 

Function/Subroutine Documentation

◆ materialdata_sp()

subroutine materialdata_sp ( real*8, dimension(0:ncmat_,ntmat_,*), intent(in)  elcon,
integer, dimension(2,*), intent(in)  nelcon,
integer, intent(in)  imat,
integer, intent(in)  ntmat_,
integer, intent(in)  i,
real*8, intent(in)  t1l,
real*8, dimension(21), intent(inout)  elconloc,
integer, intent(in)  kode,
real*8, dimension(0:2*npmat_,ntmat_,*), intent(in)  plicon,
integer, dimension(0:ntmat_,*), intent(in)  nplicon,
integer, intent(in)  npmat_,
real*8, dimension(802), intent(inout)  plconloc,
integer, intent(in)  ncmat_ 
)
20 !
21  implicit none
22 !
23 ! determines the material data for element i
24 !
25  integer nelcon(2,*),imat,i,k,kin,ntmat_,nelconst,kode,
26  & itemp,ncmat_,id,nplicon(0:ntmat_,*),npmat_
27 !
28  real*8 elcon(0:ncmat_,ntmat_,*),t1l,elconloc(21),
29  & plicon(0:2*npmat_,ntmat_,*),plconloc(802)
30 !
31  intent(in) elcon,nelcon,imat,ntmat_,i,t1l,
32  & kode,plicon,nplicon,npmat_,ncmat_
33 !
34  intent(inout) elconloc,plconloc
35 !
36 ! nelconst: # constants read from file
37 !
38 ! calculating the hardening coefficients
39 !
40 ! for the calculation of the spring stiffness, the whole curve
41 ! has to be stored:
42 ! plconloc(2*k-1), k=1...200: displacement
43 ! plconloc(2*k),k=1...200: force
44 !
45  if(kode.lt.-50) then
46  if(npmat_.eq.0) then
47  plconloc(801)=0.5d0
48  plconloc(802)=0.5d0
49  else
50  plconloc(1)=0.d0
51  plconloc(2)=0.d0
52  plconloc(3)=0.d0
53  plconloc(801)=nplicon(1,imat)+0.5d0
54  plconloc(802)=0.5d0
55 !
56 ! nonlinear spring characteristic or gap conductance characteristic
57 !
58  if(nplicon(1,imat).ne.0) then
59 !
60  if(nplicon(0,imat).eq.1) then
61  id=-1
62  else
63  call ident2(plicon(0,1,imat),t1l,nplicon(0,imat),
64  & 2*npmat_+1,id)
65  endif
66 !
67  if(nplicon(0,imat).eq.0) then
68  continue
69  elseif((nplicon(0,imat).eq.1).or.(id.eq.0).or.
70  & (id.eq.nplicon(0,imat))) then
71  if(id.le.0) then
72  itemp=1
73  else
74  itemp=id
75  endif
76  kin=0
77  call plcopy(plicon,nplicon,plconloc,npmat_,ntmat_,
78  & imat,itemp,i,kin)
79  if((id.eq.0).or.(id.eq.nplicon(0,imat))) then
80  endif
81  else
82  kin=0
83  call plmix(plicon,nplicon,plconloc,npmat_,ntmat_,
84  & imat,id+1,t1l,i,kin)
85  endif
86  endif
87  endif
88  else
89 !
90 ! linear spring characteristic
91 !
92  nelconst=nelcon(1,imat)
93  call ident2(elcon(0,1,imat),t1l,nelcon(2,imat),ncmat_+1,id)
94  if(nelcon(2,imat).eq.0) then
95  continue
96  elseif(nelcon(2,imat).eq.1) then
97  do k=1,nelconst
98  elconloc(k)=elcon(k,1,imat)
99  enddo
100  elseif(id.eq.0) then
101  do k=1,nelconst
102  elconloc(k)=elcon(k,1,imat)
103  enddo
104  elseif(id.eq.nelcon(2,imat)) then
105  do k=1,nelconst
106  elconloc(k)=elcon(k,id,imat)
107  enddo
108  else
109  do k=1,nelconst
110  elconloc(k)=elcon(k,id,imat)+
111  & (elcon(k,id+1,imat)-elcon(k,id,imat))*
112  & (t1l-elcon(0,id,imat))/
113  & (elcon(0,id+1,imat)-elcon(0,id,imat))
114  enddo
115  endif
116  endif
117 !
118  return
subroutine plcopy(plcon, nplcon, plconloc, npmat_, ntmat_, imat, itemp, nelem, kin)
Definition: plcopy.f:21
subroutine ident2(x, px, n, ninc, id)
Definition: ident2.f:27
subroutine plmix(plcon, nplcon, plconloc, npmat_, ntmat_, imat, j, temp, nelem, kin)
Definition: plmix.f:21
Hosted by OpenAircraft.com, (Michigan UAV, LLC)