28 character*80 matname(*)
29 character*132 textpart(16)
31 integer nelcon(2,*),nmat,ntmat_,ntmat,npmat_,npmat,istep,
32 & n,key,i,nplicon(0:ntmat_,*),istat,ncmat_,itemp,id,ipoinpc(0:*),
33 & ndata,ndatamax,kin,irstrt,iline,ipol,inl,ipoinp(2,*),inp(3,*)
35 real*8 plicon(0:2*npmat_,ntmat_,*),temperature,
36 & elcon(0:ncmat_,ntmat_,*),plconloc(802),t1l
41 if((istep.gt.0).and.(irstrt.ge.0))
then 42 write(*,*)
'*ERROR in cychards: *CYCLIC HARDENING' 43 write(*,*)
' should be placed before all step' 44 write(*,*)
' definitions' 49 write(*,*)
'*ERROR in cychards: *CYCLIC HARDENING' 50 write(*,*)
' should be preceded' 51 write(*,*)
' by a *MATERIAL card' 55 if(((nelcon(1,nmat).ne.-51).and.(nelcon(1,nmat).ne.-114)).or.
56 & (nplicon(0,nmat).ne.0))
then 57 write(*,*)
'*ERROR in cychards: *CYCLIC HARDENING' 58 write(*,*)
' should be preceded' 59 write(*,*)
' by an *PLASTIC,HARDENING=COMBINED card' 65 &
'*WARNING in cychards: parameter not recognized:' 67 & textpart(i)(1:index(textpart(i),
' ')-1)
69 &
"*CYCLIC HARDENING%")
75 call getnewline(inpc,textpart,istat,n,key,iline,ipol,inl,
77 if((istat.lt.0).or.(key.eq.1))
exit 78 read(textpart(3)(1:20),
'(f20.0)',iostat=istat) temperature
79 if(istat.gt.0)
call inputerror(inpc,ipoinpc,iline,
80 &
"*CYCLIC HARDENING%")
87 if(ntmat.gt.ntmat_)
then 88 write(*,*)
'*ERROR in cychards: increase ntmat_' 92 plicon(0,ntmat,nmat)=temperature
96 elseif(plicon(0,ntmat,nmat).ne.temperature)
then 99 if(ntmat.gt.ntmat_)
then 100 write(*,*)
'*ERROR in cychards: increase ntmat_' 103 nplicon(0,nmat)=ntmat
104 plicon(0,ntmat,nmat)=temperature
107 read(textpart(i)(1:20),
'(f20.0)',iostat=istat)
108 & plicon(2*npmat+i,ntmat,nmat)
109 if(istat.gt.0)
call inputerror(inpc,ipoinpc,iline,
110 &
"*CYCLIC HARDENING%")
113 if(npmat.gt.npmat_)
then 114 write(*,*)
'*ERROR in cychards: increase npmat_' 117 nplicon(ntmat,nmat)=npmat
121 write(*,*)
'*ERROR in cychards: *CYCLIC HARDENING card' 122 write(*,*)
' without data encountered' 129 if(nelcon(1,nmat).eq.-114)
then 137 do i=1,nelcon(2,nmat)
140 if(nplicon(0,nmat).eq.1)
then 143 call ident2(plicon(0,1,nmat),t1l,nplicon(0,nmat),
147 if(nplicon(0,nmat).eq.0)
then 149 elseif((nplicon(0,nmat).eq.1).or.(id.eq.0).or.
150 & (id.eq.nplicon(0,nmat)))
then 157 call plcopy(plicon,nplicon,plconloc,npmat_,ntmat_,
159 if((id.eq.0).or.(id.eq.nplicon(0,nmat)))
then 163 call plmix(plicon,nplicon,plconloc,npmat_,ntmat_,
164 & nmat,id+1,t1l,i,kin)
167 ndata=int(plconloc(801))
169 elcon(10,i,nmat)=plconloc(2)
170 elcon(11,i,nmat)=0.d0
172 elcon(10,i,nmat)=plconloc(2)
173 elcon(11,i,nmat)=(plconloc(4)-plconloc(2))/
174 & (plconloc(3)-plconloc(1))
176 ndatamax=
max(ndata,ndatamax)
178 if(ndatamax.gt.2)
then 179 write(*,*)
'*WARNING in plastics: isotropic hardening' 180 write(*,*)
' curve is possibly nonlinear for' 181 write(*,*)
' the elastically anisotropic' 182 write(*,*)
' material ',matname(nmat)(71:80)
subroutine plcopy(plcon, nplcon, plconloc, npmat_, ntmat_, imat, itemp, nelem, kin)
Definition: plcopy.f:21
#define max(a, b)
Definition: cascade.c:32
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