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

Go to the source code of this file.

Functions/Subroutines

subroutine plinterpol (plcon, nplcon, itemp, f, df, npmat_, ntmat_, imat, nelem, epl)
 

Function/Subroutine Documentation

◆ plinterpol()

subroutine plinterpol ( real*8, dimension(0:2*npmat_,ntmat_,*), intent(in)  plcon,
integer, dimension(0:ntmat_,*), intent(in)  nplcon,
integer, intent(in)  itemp,
real*8, intent(inout)  f,
real*8, intent(inout)  df,
integer, intent(in)  npmat_,
integer, intent(in)  ntmat_,
integer, intent(in)  imat,
integer, intent(in)  nelem,
real*8, intent(in)  epl 
)
20 !
21  implicit none
22 !
23 ! interpolation of isotropic or kinematic hardening data
24 ! input: hardening data plcon and nplcon, temperature itemp,
25 ! size parameters npmat_ and ntmat_, material number imat
26 ! and equivalent plastic strain at which the coefficients
27 ! are to be determined
28 ! output: hardening coefficient and its local derivative f and df
29 !
30  integer npmat_,ntmat_,nplcon(0:ntmat_,*),itemp,ndata,imat,j,
31  & nelem
32 !
33  real*8 plcon(0:2*npmat_,ntmat_,*),f,df,epl
34 !
35  intent(in) plcon,nplcon,itemp,npmat_,ntmat_,
36  & imat,nelem,epl
37 !
38  intent(inout) f,df
39 !
40  ndata=nplcon(itemp,imat)
41 !
42  do j=1,ndata
43  if(epl.lt.plcon(2*j,itemp,imat)) exit
44  enddo
45 !
46  if((j.eq.1).or.(j.gt.ndata)) then
47  if(j.eq.1) then
48  f=plcon(1,itemp,imat)
49  df=0.d0
50  else
51  f=plcon(2*ndata-1,itemp,imat)
52  df=0.d0
53  endif
54  write(*,*) '*WARNING in plinterpol: plastic strain ',epl
55  write(*,*) ' outside material plastic strain range'
56  write(*,*) ' in element ',nelem,' and material ',imat
57  write(*,*) ' for temperature ',plcon(0,itemp,imat)
58  else
59  df=(plcon(2*j-1,itemp,imat)-plcon(2*j-3,itemp,imat))/
60  & (plcon(2*j,itemp,imat)-plcon(2*j-2,itemp,imat))
61  f=plcon(2*j-3,itemp,imat)+
62  & df*(epl-plcon(2*j-2,itemp,imat))
63  endif
64 !
65  return
subroutine df(x, u, uprime, rpar, nev)
Definition: subspace.f:133
Hosted by OpenAircraft.com, (Michigan UAV, LLC)