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

Go to the source code of this file.

Functions/Subroutines

subroutine cyclichardenings (inpc, textpart, nelcon, nmat, ntmat_, npmat_, plicon, nplicon, ncmat_, elcon, matname, irstrt, istep, istat, n, iline, ipol, inl, ipoinp, inp, ipoinpc)
 

Function/Subroutine Documentation

◆ cyclichardenings()

subroutine cyclichardenings ( 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,
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 
)
22 !
23 ! reading the input deck: *CYCLIC HARDENING
24 !
25  implicit none
26 !
27  character*1 inpc(*)
28  character*80 matname(*)
29  character*132 textpart(16)
30 !
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,*)
34 !
35  real*8 plicon(0:2*npmat_,ntmat_,*),temperature,
36  & elcon(0:ncmat_,ntmat_,*),plconloc(802),t1l
37 !
38  ntmat=0
39  npmat=0
40 !
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'
45  call exit(201)
46  endif
47 !
48  if(nmat.eq.0) then
49  write(*,*) '*ERROR in cychards: *CYCLIC HARDENING'
50  write(*,*) ' should be preceded'
51  write(*,*) ' by a *MATERIAL card'
52  call exit(201)
53  endif
54 !
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'
60  call exit(201)
61  endif
62 !
63  do i=2,n
64  write(*,*)
65  & '*WARNING in cychards: parameter not recognized:'
66  write(*,*) ' ',
67  & textpart(i)(1:index(textpart(i),' ')-1)
68  call inputwarning(inpc,ipoinpc,iline,
69  &"*CYCLIC HARDENING%")
70  enddo
71 !
72 ! isotropic hardening coefficients
73 !
74  do
75  call getnewline(inpc,textpart,istat,n,key,iline,ipol,inl,
76  & ipoinp,inp,ipoinpc)
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%")
81 !
82 ! first temperature
83 !
84  if(ntmat.eq.0) then
85  npmat=0
86  ntmat=ntmat+1
87  if(ntmat.gt.ntmat_) then
88  write(*,*) '*ERROR in cychards: increase ntmat_'
89  call exit(201)
90  endif
91  nplicon(0,nmat)=ntmat
92  plicon(0,ntmat,nmat)=temperature
93 !
94 ! new temperature
95 !
96  elseif(plicon(0,ntmat,nmat).ne.temperature) then
97  npmat=0
98  ntmat=ntmat+1
99  if(ntmat.gt.ntmat_) then
100  write(*,*) '*ERROR in cychards: increase ntmat_'
101  call exit(201)
102  endif
103  nplicon(0,nmat)=ntmat
104  plicon(0,ntmat,nmat)=temperature
105  endif
106  do i=1,2
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%")
111  enddo
112  npmat=npmat+1
113  if(npmat.gt.npmat_) then
114  write(*,*) '*ERROR in cychards: increase npmat_'
115  call exit(201)
116  endif
117  nplicon(ntmat,nmat)=npmat
118  enddo
119 !
120  if(ntmat.eq.0) then
121  write(*,*) '*ERROR in cychards: *CYCLIC HARDENING card'
122  write(*,*) ' without data encountered'
123  call exit(201)
124  endif
125 !
126 ! elastically anisotropic materials: recasting the input data
127 ! in a format conform to the user routine umat_aniso_plas.f
128 !
129  if(nelcon(1,nmat).eq.-114) then
130 !
131 ! isotropic hardening
132 !
133 ! interpolating the plastic data at the elastic temperature
134 ! data points
135 !
136  ndatamax=0
137  do i=1,nelcon(2,nmat)
138  t1l=elcon(0,i,nmat)
139 !
140  if(nplicon(0,nmat).eq.1) then
141  id=-1
142  else
143  call ident2(plicon(0,1,nmat),t1l,nplicon(0,nmat),
144  & 2*npmat_+1,id)
145  endif
146 !
147  if(nplicon(0,nmat).eq.0) then
148  continue
149  elseif((nplicon(0,nmat).eq.1).or.(id.eq.0).or.
150  & (id.eq.nplicon(0,nmat))) then
151  if(id.le.0) then
152  itemp=1
153  else
154  itemp=id
155  endif
156  kin=0
157  call plcopy(plicon,nplicon,plconloc,npmat_,ntmat_,
158  & nmat,itemp,i,kin)
159  if((id.eq.0).or.(id.eq.nplicon(0,nmat))) then
160  endif
161  else
162  kin=0
163  call plmix(plicon,nplicon,plconloc,npmat_,ntmat_,
164  & nmat,id+1,t1l,i,kin)
165  endif
166 !
167  ndata=int(plconloc(801))
168  if(ndata.eq.1) then
169  elcon(10,i,nmat)=plconloc(2)
170  elcon(11,i,nmat)=0.d0
171  else
172  elcon(10,i,nmat)=plconloc(2)
173  elcon(11,i,nmat)=(plconloc(4)-plconloc(2))/
174  & (plconloc(3)-plconloc(1))
175  endif
176  ndatamax=max(ndata,ndatamax)
177  enddo
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)
183  endif
184  endif
185 !
186  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)