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

Go to the source code of this file.

Functions/Subroutines

subroutine plastics (inpc, textpart, nelcon, nmat, ntmat_, npmat_, plicon, nplicon, plkcon, nplkcon, iplas, iperturb, nstate_, ncmat_, elcon, matname, irstrt, istep, istat, n, iline, ipol, inl, ipoinp, inp, ipoinpc, ianisoplas)
 

Function/Subroutine Documentation

◆ plastics()

subroutine plastics ( character*1, dimension(*)  inpc,
character*132, dimension(16)  textpart,
integer, dimension(2,*)  nelcon,
integer  nmat,
integer  ntmat_,
integer  npmat_,
real*8, dimension(0:2*npmat_,ntmat_,*)  plicon,
integer, dimension(0:ntmat_,*)  nplicon,
real*8, dimension(0:2*npmat_,ntmat_,*)  plkcon,
integer, dimension(0:ntmat_,*)  nplkcon,
integer  iplas,
integer, dimension(*)  iperturb,
integer  nstate_,
integer  ncmat_,
real*8, dimension(0:ncmat_,ntmat_,*)  elcon,
character*80, dimension(*)  matname,
integer  irstrt,
integer  istep,
integer  istat,
integer  n,
integer  iline,
integer  ipol,
integer  inl,
integer, dimension(2,*)  ipoinp,
integer, dimension(3,*)  inp,
integer, dimension(0:*)  ipoinpc,
integer  ianisoplas 
)
23 !
24 ! reading the input deck: *PLASTIC
25 !
26  implicit none
27 !
28  logical iso
29 !
30  character*1 inpc(*)
31  character*80 matname(*)
32  character*132 textpart(16)
33 !
34  integer nelcon(2,*),nmat,ntmat_,ntmat,npmat_,npmat,istep,
35  & n,key,i,nplicon(0:ntmat_,*),nplkcon(0:ntmat_,*),ncmat_,
36  & iplas,iperturb(*),istat,nstate_,kin,itemp,ndata,ndatamax,id,
37  & irstrt,iline,ipol,inl,ipoinp(2,*),inp(3,*),ipoinpc(0:*),
38  & ianisoplas
39 !
40  real*8 plicon(0:2*npmat_,ntmat_,*),plkcon(0:2*npmat_,ntmat_,*),
41  & temperature,plconloc(802),t1l,elcon(0:ncmat_,ntmat_,*)
42 !
43  iso=.true.
44 !
45  ntmat=0
46  npmat=0
47 !
48  if((istep.gt.0).and.(irstrt.ge.0)) then
49  write(*,*) '*ERROR reading *PLASTIC: *PLASTIC should be placed'
50  write(*,*) ' before all step definitions'
51  call exit(201)
52  endif
53 !
54  if(nmat.eq.0) then
55  write(*,*)
56  & '*ERROR reading *PLASTIC: *PLASTIC should be preceded'
57  write(*,*) ' by a *MATERIAL card'
58  call exit(201)
59  endif
60 !
61  if((nelcon(1,nmat).ne.2).and.(nelcon(1,nmat).ne.9)) then
62  write(*,*)
63  & '*ERROR reading *PLASTIC: *PLASTIC should be preceded'
64  write(*,*) ' by an *ELASTIC,TYPE=ISO card or'
65  write(*,*) ' by an *ELASTIC,TYPE=ORTHO card'
66  call exit(201)
67  endif
68 !
69  iperturb(1)=3
70 c iperturb(2)=1
71 c write(*,*) '*INFO reading *PLASTIC: nonlinear geometric'
72 c write(*,*) ' effects are turned on'
73 c write(*,*)
74 !
75  if(nelcon(1,nmat).eq.2) then
76  iplas=1
77  nelcon(1,nmat)=-51
78  nstate_=max(nstate_,13)
79  else
80  ianisoplas=1
81  nelcon(1,nmat)=-114
82  nstate_=max(nstate_,14)
83  endif
84 !
85  do i=2,n
86  if(textpart(i)(1:10).eq.'HARDENING=') then
87  if(textpart(i)(11:19).eq.'KINEMATIC') then
88  iso=.false.
89  elseif(textpart(i)(11:18).eq.'COMBINED') then
90  iso=.false.
91  elseif(textpart(i)(11:14).eq.'USER') then
92  if(nelcon(1,nmat).eq.-114) then
93  write(*,*) '*ERROR reading *PLASTIC: user defined '
94  write(*,*) ' hardening is not allowed for '
95  write(*,*) ' elastically anisotropic materials'
96  call exit(201)
97  endif
98  call getnewline(inpc,textpart,istat,n,key,iline,ipol,inl,
99  & ipoinp,inp,ipoinpc)
100  return
101  endif
102  exit
103  else
104  write(*,*)
105  & '*WARNING reading *PLASTIC: parameter not recognized:'
106  write(*,*) ' ',
107  & textpart(i)(1:index(textpart(i),' ')-1)
108  call inputwarning(inpc,ipoinpc,iline,
109  &"*PLASTIC%")
110  endif
111  enddo
112 !
113  if(iso) then
114 !
115 ! isotropic hardening coefficients
116 !
117  do
118  call getnewline(inpc,textpart,istat,n,key,iline,ipol,inl,
119  & ipoinp,inp,ipoinpc)
120  if((istat.lt.0).or.(key.eq.1)) exit
121  read(textpart(3)(1:20),'(f20.0)',iostat=istat) temperature
122  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
123  &"*PLASTIC%")
124 !
125 ! first temperature
126 !
127  if(ntmat.eq.0) then
128  npmat=0
129  ntmat=ntmat+1
130  if(ntmat.gt.ntmat_) then
131  write(*,*) '*ERROR reading *PLASTIC: increase ntmat_'
132  call exit(201)
133  endif
134  nplicon(0,nmat)=ntmat
135  plicon(0,ntmat,nmat)=temperature
136 !
137 ! new temperature
138 !
139  elseif(plicon(0,ntmat,nmat).ne.temperature) then
140  npmat=0
141  ntmat=ntmat+1
142  if(ntmat.gt.ntmat_) then
143  write(*,*) '*ERROR reading *PLASTIC: increase ntmat_'
144  call exit(201)
145  endif
146  nplicon(0,nmat)=ntmat
147  plicon(0,ntmat,nmat)=temperature
148  endif
149  do i=1,2
150  read(textpart(i)(1:20),'(f20.0)',iostat=istat)
151  & plicon(2*npmat+i,ntmat,nmat)
152  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
153  &"*PLASTIC%")
154  enddo
155  npmat=npmat+1
156  if(npmat.gt.npmat_) then
157  write(*,*) '*ERROR reading *PLASTIC: increase npmat_'
158  call exit(201)
159  endif
160  nplicon(ntmat,nmat)=npmat
161  enddo
162  else
163 !
164 ! kinematic hardening coefficients
165 !
166  do
167  call getnewline(inpc,textpart,istat,n,key,iline,ipol,inl,
168  & ipoinp,inp,ipoinpc)
169  if((istat.lt.0).or.(key.eq.1)) exit
170  read(textpart(3)(1:20),'(f20.0)',iostat=istat) temperature
171  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
172  &"*PLASTIC%")
173 !
174 ! first temperature
175 !
176  if(ntmat.eq.0) then
177  npmat=0
178  ntmat=ntmat+1
179  if(ntmat.gt.ntmat_) then
180  write(*,*) '*ERROR reading *PLASTIC: increase ntmat_'
181  call exit(201)
182  endif
183  nplkcon(0,nmat)=ntmat
184  plkcon(0,ntmat,nmat)=temperature
185 !
186 ! new temperature
187 !
188  elseif(plkcon(0,ntmat,nmat).ne.temperature) then
189  npmat=0
190  ntmat=ntmat+1
191  if(ntmat.gt.ntmat_) then
192  write(*,*) '*ERROR reading *PLASTIC: increase ntmat_'
193  call exit(201)
194  endif
195  nplkcon(0,nmat)=ntmat
196  plkcon(0,ntmat,nmat)=temperature
197  endif
198  do i=1,2
199  read(textpart(i)(1:20),'(f20.0)',iostat=istat)
200  & plkcon(2*npmat+i,ntmat,nmat)
201  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
202  &"*PLASTIC%")
203  enddo
204  npmat=npmat+1
205  if(npmat.gt.npmat_) then
206  write(*,*) '*ERROR reading *PLASTIC: increase npmat_'
207  call exit(201)
208  endif
209  nplkcon(ntmat,nmat)=npmat
210  enddo
211  endif
212 !
213  if(ntmat.eq.0) then
214  write(*,*)
215  & '*ERROR reading *PLASTIC: *PLASTIC card without data'
216  call exit(201)
217  endif
218 !
219 ! elastically anisotropic materials: recasting the input data
220 ! in a format conform to the user routine umat_aniso_plas.f
221 !
222  if(nelcon(1,nmat).eq.-114) then
223  if(matname(nmat)(71:80).ne.' ') then
224  write(*,*)
225  & '*ERROR reading *PLASTIC: the material name for an'
226  write(*,*) ' elastically anisotropic material with'
227  write(*,*) ' isotropic plasticity must not exceed 70'
228  write(*,*) ' characters'
229  call exit(201)
230  else
231  do i=80,11,-1
232  matname(nmat)(i:i)=matname(nmat)(i-10:i-10)
233  enddo
234 c matname(nmat)(11:80)=matname(nmat)(1:70)
235  matname(nmat)(1:10)='ANISO_PLAS'
236  endif
237 !
238  if(iso) then
239 !
240 ! isotropic hardening
241 !
242 ! interpolating the plastic data at the elastic temperature
243 ! data points
244 !
245  ndatamax=0
246  do i=1,nelcon(2,nmat)
247  t1l=elcon(0,i,nmat)
248 c plconloc(1)=0.d0
249 c plconloc(2)=0.d0
250 c plconloc(3)=0.d0
251 c plconloc(81)=nplicon(1,nmat)+0.5d0
252 !
253  if(nplicon(0,nmat).eq.1) then
254  id=-1
255  else
256  call ident2(plicon(0,1,nmat),t1l,nplicon(0,nmat),
257  & 2*npmat_+1,id)
258  endif
259 !
260  if(nplicon(0,nmat).eq.0) then
261  continue
262  elseif((nplicon(0,nmat).eq.1).or.(id.eq.0).or.
263  & (id.eq.nplicon(0,nmat))) then
264  if(id.le.0) then
265  itemp=1
266  else
267  itemp=id
268  endif
269  kin=0
270  call plcopy(plicon,nplicon,plconloc,npmat_,ntmat_,
271  & nmat,itemp,i,kin)
272  if((id.eq.0).or.(id.eq.nplicon(0,nmat))) then
273  endif
274  else
275  kin=0
276  call plmix(plicon,nplicon,plconloc,npmat_,ntmat_,
277  & nmat,id+1,t1l,i,kin)
278  endif
279 !
280  ndata=int(plconloc(801))
281  if(ndata.eq.1) then
282  elcon(10,i,nmat)=plconloc(2)
283  elcon(11,i,nmat)=0.d0
284  elcon(12,i,nmat)=0.d0
285  elcon(13,i,nmat)=-1.d0
286  elcon(14,i,nmat)=1.d0
287  else
288  elcon(10,i,nmat)=plconloc(2)
289  elcon(11,i,nmat)=(plconloc(4)-plconloc(2))/
290  & (plconloc(3)-plconloc(1))
291  elcon(12,i,nmat)=0.d0
292  elcon(13,i,nmat)=-1.d0
293  elcon(14,i,nmat)=1.d0
294  endif
295  ndatamax=max(ndata,ndatamax)
296  enddo
297  if(ndatamax.gt.2) then
298  write(*,*)
299  & '*WARNING reading *PLASTIC: isotropic hardening'
300  write(*,*) ' curve is possibly nonlinear for'
301  write(*,*) ' the elastically anisotropic'
302  write(*,*) ' material ',matname(nmat)(11:80)
303  endif
304  else
305 !
306 ! kinematic hardening
307 !
308 ! interpolating the plastic data at the elastic temperature
309 ! data points
310 !
311  ndatamax=0
312  do i=1,nelcon(2,nmat)
313  t1l=elcon(0,i,nmat)
314 c plconloc(1)=0.d0
315 c plconloc(2)=0.d0
316 c plconloc(3)=0.d0
317 c plconloc(82)=nplkcon(1,nmat)+0.5d0
318 !
319  if(nplkcon(0,nmat).eq.1) then
320  id=-1
321  else
322  call ident2(plkcon(0,1,nmat),t1l,nplkcon(0,nmat),
323  & 2*npmat_+1,id)
324  endif
325 !
326  if(nplkcon(0,nmat).eq.0) then
327  continue
328  elseif((nplkcon(0,nmat).eq.1).or.(id.eq.0).or.
329  & (id.eq.nplkcon(0,nmat))) then
330  if(id.le.0) then
331  itemp=1
332  else
333  itemp=id
334  endif
335  kin=1
336  call plcopy(plkcon,nplkcon,plconloc,npmat_,ntmat_,
337  & nmat,itemp,i,kin)
338  if((id.eq.0).or.(id.eq.nplkcon(0,nmat))) then
339  endif
340  else
341  kin=1
342  call plmix(plkcon,nplkcon,plconloc,npmat_,ntmat_,
343  & nmat,id+1,t1l,i,kin)
344  endif
345 !
346  ndata=int(plconloc(802))
347  if(ndata.eq.1) then
348  elcon(10,i,nmat)=plconloc(402)
349  elcon(11,i,nmat)=0.d0
350  elcon(12,i,nmat)=0.d0
351  elcon(13,i,nmat)=-1.d0
352  elcon(14,i,nmat)=1.d0
353  else
354  elcon(10,i,nmat)=plconloc(402)
355  elcon(11,i,nmat)=0.d0
356  elcon(12,i,nmat)=(plconloc(404)-plconloc(402))/
357  & (plconloc(403)-plconloc(401))
358  elcon(13,i,nmat)=-1.d0
359  elcon(14,i,nmat)=1.d0
360  endif
361  ndatamax=max(ndata,ndatamax)
362  enddo
363  if(ndatamax.gt.2) then
364  write(*,*)
365  & '*WARNING reading *PLASTIC: kinematic hardening'
366  write(*,*) ' curve is possibly nonlinear for'
367  write(*,*) ' the elastically anisotropic'
368  write(*,*) ' material ',matname(nmat)(11:80)
369  endif
370  endif
371  endif
372 !
373 c if(nelcon(1,nmat).eq.-114) then
374 c write(*,*) 'anisotropic elasticity+viscoplasticity'
375 c do i=1,nelcon(2,nmat)
376 c write(*,*) (elcon(j,i,nmat),j=0,14)
377 c enddo
378 c endif
379 !
380  return
subroutine plcopy(plcon, nplcon, plconloc, npmat_, ntmat_, imat, itemp, nelem, kin)
Definition: plcopy.f:21
#define max(a, b)
Definition: cascade.c:32
subroutine inputwarning(inpc, ipoinpc, iline, text)
Definition: inputwarning.f:20
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
subroutine getnewline(inpc, textpart, istat, n, key, iline, ipol, inl, ipoinp, inp, ipoinpc)
Definition: getnewline.f:21
subroutine inputerror(inpc, ipoinpc, iline, text)
Definition: inputerror.f:20
Hosted by OpenAircraft.com, (Michigan UAV, LLC)