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

Go to the source code of this file.

Functions/Subroutines

subroutine plmix (plcon, nplcon, plconloc, npmat_, ntmat_, imat, j, temp, nelem, kin)
 

Function/Subroutine Documentation

◆ plmix()

subroutine plmix ( real*8, dimension(0:2*npmat_,ntmat_,*), intent(in)  plcon,
integer, dimension(0:ntmat_,*), intent(in)  nplcon,
real*8, dimension(802), intent(inout)  plconloc,
integer, intent(in)  npmat_,
integer, intent(in)  ntmat_,
integer, intent(in)  imat,
integer, intent(in)  j,
real*8, intent(in)  temp,
integer, intent(in)  nelem,
integer, intent(in)  kin 
)
21 !
22 ! interpolates the hardening data for material imat and temperature
23 ! j and j-1 to obtain data for temperature temp. The data is taken
24 ! from plcon and stored in plconloc.
25 ! The Von Mises stress is interpolated for a given equivalent
26 ! plastic strain. If the equivalent strain data points for
27 ! temperature j and j-1 do not coincide, the union of both is
28 ! taken. If this union exceeds 200 (ierror=1), the equivalent plastic
29 ! strain range is divided into 199 intervals yielding 200 new
30 ! equivalent strain data points, for which the Von Mises stress
31 ! is interpolated.
32 ! Attention: in plcon the odd storage spaces contain the Von
33 ! Mises stress, the even ones the equivalent plastic
34 ! strain. For plconloc, this order is reversed.
35 !
36  implicit none
37 !
38  integer imat,ndata,ntmat_,npmat_,nplcon(0:ntmat_,*),nelem,
39  & kin,k,j,k1,k2,ierror,ndata1,ndata2,itemp
40 !
41  real*8 eplmin,eplmax,depl,epla,plcon(0:2*npmat_,ntmat_,*),
42  & plconloc(802),dummy,temp,ep1,ep2,t1,t2,s1,s2,ratio
43 !
44  intent(in) plcon,nplcon,npmat_,ntmat_,
45  & imat,j,temp,nelem,kin
46 !
47  intent(inout) plconloc
48 !
49  ndata=0
50  ierror=0
51 !
52  ndata1=nplcon(j-1,imat)
53  ndata2=nplcon(j,imat)
54  t1=plcon(0,j-1,imat)
55  t2=plcon(0,j,imat)
56  ratio=(temp-t1)/(t2-t1)
57 !
58 ! the interval on which the stress interpolation is performed
59 ! is the intersection of the domain of the two curves
60 !
61  k1=1
62  k2=1
63  ep1=plcon(2,j-1,imat)
64  ep2=plcon(2,j,imat)
65  if(ep1.gt.ep2) then
66  do k2=1,ndata2
67  ep2=plcon(2*k2,j,imat)
68  if(ep2.gt.ep1) exit
69  enddo
70  if(k2.gt.ndata2) then
71  write(*,*) '*ERROR in plmix: there exist two temperatures'
72  write(*,*) ' for which the hardening curves are'
73  write(*,*) ' disjunct'
74  call exit(201)
75  endif
76  elseif(ep2.gt.ep1) then
77  do k1=1,ndata1
78  ep1=plcon(2*k1,j-1,imat)
79  if(ep1.gt.ep2) exit
80  enddo
81  if(k1.gt.ndata1) then
82  write(*,*) '*ERROR in plmix: there exist two temperatures'
83  write(*,*) ' for which the hardening curves are'
84  write(*,*) ' disjunct'
85  call exit(201)
86  endif
87  endif
88 !
89  do
90  s1=plcon(2*k1-1,j-1,imat)
91  s2=plcon(2*k2-1,j,imat)
92  ep1=plcon(2*k1,j-1,imat)
93  ep2=plcon(2*k2,j,imat)
94 !
95  if(dabs(ep1-ep2).lt.1.d-10) then
96  if(k2.lt.ndata2) then
97  k2=k2+1
98  elseif(k1.lt.ndata1) then
99  k1=k1+1
100  else
101  ndata=ndata+1
102  if(ndata.gt.200) then
103  ierror=1
104  exit
105  endif
106  if(kin.eq.0) then
107  plconloc(2*ndata-1)=ep1+ratio*(ep2-ep1)
108  plconloc(2*ndata)=s1+ratio*(s2-s1)
109  else
110  plconloc(399+2*ndata)=ep1+ratio*(ep2-ep1)
111  plconloc(400+2*ndata)=s1+ratio*(s2-s1)
112  endif
113  exit
114  endif
115  cycle
116  endif
117  if(ep1.lt.ep2) then
118  ndata=ndata+1
119  if(ndata.gt.200) then
120  ierror=1
121  exit
122  endif
123  call plinterpol(plcon,nplcon,j,s2,dummy,npmat_,ntmat_,
124  & imat,nelem,ep1)
125  if(kin.eq.0) then
126  plconloc(2*ndata-1)=ep1
127  plconloc(2*ndata)=s1+ratio*(s2-s1)
128  else
129  plconloc(399+2*ndata)=ep1
130  plconloc(400+2*ndata)=s1+ratio*(s2-s1)
131  endif
132  if(k1.lt.ndata1) then
133  k1=k1+1
134  cycle
135  else
136  exit
137  endif
138  else
139  ndata=ndata+1
140  if(ndata.gt.200) then
141  ierror=1
142  exit
143  endif
144  call plinterpol(plcon,nplcon,j-1,s1,dummy,npmat_,ntmat_,
145  & imat,nelem,ep2)
146  if(kin.eq.0) then
147  plconloc(2*ndata-1)=ep2
148  plconloc(2*ndata)=s1+ratio*(s2-s1)
149  else
150  plconloc(399+2*ndata)=ep2
151  plconloc(400+2*ndata)=s1+ratio*(s2-s1)
152  endif
153  if(k2.lt.ndata2) then
154  k2=k2+1
155  cycle
156  else
157  exit
158  endif
159  endif
160  enddo
161 !
162 ! if more than 200 data points result, the interval is divided into
163 ! 199 equidistant intervals
164 !
165  if(ierror.eq.0) then
166  if(kin.eq.0) then
167  plconloc(801)=real(ndata)+0.5d0
168  else
169  plconloc(802)=real(ndata)+0.5d0
170  endif
171  else
172  if(kin.eq.0) then
173  eplmin=max(plcon(2,j-1,imat),plcon(2,j,imat))
174  eplmax=min(plcon(2*ndata1,j-1,imat),plcon(2*ndata2,j,imat))
175  & -1.d-10
176  depl=(eplmax-eplmin)/199.d0
177  do k=1,200
178  epla=eplmin+(k-1)*depl
179  itemp=j-1
180  call plinterpol(plcon,nplcon,itemp,s1,
181  & dummy,npmat_,ntmat_,imat,nelem,epla)
182  itemp=j
183  call plinterpol(plcon,nplcon,itemp,s2,
184  & dummy,npmat_,ntmat_,imat,nelem,epla)
185  plconloc(2*k-1)=epla
186  plconloc(2*k)=s1+ratio*(s2-s1)
187  enddo
188  plconloc(801)=200.5d0
189  else
190  eplmin=max(plcon(2,j-1,imat),plcon(2,j,imat))
191  eplmax=min(plcon(2*ndata1,j-1,imat),plcon(2*ndata2,j,imat))
192  & -1.d-10
193  depl=(eplmax-eplmin)/199.d0
194  do k=1,200
195  epla=eplmin+(k-1)*depl
196  itemp=j-1
197  call plinterpol(plcon,nplcon,itemp,s1,
198  & dummy,npmat_,ntmat_,imat,nelem,epla)
199  itemp=j
200  call plinterpol(plcon,nplcon,itemp,s2,
201  & dummy,npmat_,ntmat_,imat,nelem,epla)
202  plconloc(399+2*k)=epla
203  plconloc(400+2*k)=s1+ratio*(s2-s1)
204  enddo
205  plconloc(802)=200.5d0
206  endif
207  endif
208 !
209  return
#define max(a, b)
Definition: cascade.c:32
#define min(a, b)
Definition: cascade.c:31
subroutine plinterpol(plcon, nplcon, itemp, f, df, npmat_, ntmat_, imat, nelem, epl)
Definition: plinterpol.f:20
Hosted by OpenAircraft.com, (Michigan UAV, LLC)