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

Go to the source code of this file.

Functions/Subroutines

subroutine restarts (istep, nset, nload, nforc, nboun, nk, ne, nmpc, nalset, nmat, ntmat_, npmat_, norien, nam, nprint, mi, ntrans, ncs_, namtot_, ncmat_, mpcfree, maxlenmpc, ne1d, ne2d, nflow, nlabel, iplas, nkon, ithermal, nmethod, iperturb, nstate_, nener, set, istartset, iendset, ialset, co, kon, ipkon, lakon, nodeboun, ndirboun, iamboun, xboun, ikboun, ilboun, ipompc, nodempc, coefmpc, labmpc, ikmpc, ilmpc, nodeforc, ndirforc, iamforc, xforc, ikforc, ilforc, nelemload, iamload, sideload, xload, elcon, nelcon, rhcon, nrhcon, alcon, nalcon, alzero, plicon, nplicon, plkcon, nplkcon, orname, orab, ielorien, trab, inotr, amname, amta, namta, t0, t1, iamt1, veold, ielmat, matname, prlab, prset, filab, vold, nodebounold, ndirbounold, xbounold, xforcold, xloadold, t1old, eme, iponor, xnor, knor, thicke, offset, iponoel, inoel, rig, shcon, nshcon, cocon, ncocon, ics, sti, ener, xstate, jobnamec, infree, irstrt, inpc, textpart, istat, n, key, prestr, iprestr, cbody, ibody, xbody, nbody, xbodyold, ttime, qaold, cs, mcs, output, physcon, ctrl, typeboun, iline, ipol, inl, ipoinp, inp, fmpc, tieset, ntie, tietol, ipoinpc, nslavs, t0g, t1g, nprop, ielprop, prop, mortar, nintpoint, ifacecount, islavsurf, pslavsurf, clearini)
 

Function/Subroutine Documentation

◆ restarts()

subroutine restarts ( integer  istep,
integer  nset,
integer  nload,
integer  nforc,
integer  nboun,
integer  nk,
integer  ne,
integer  nmpc,
integer  nalset,
integer  nmat,
integer  ntmat_,
integer  npmat_,
integer  norien,
integer  nam,
integer  nprint,
integer, dimension(*)  mi,
integer  ntrans,
integer  ncs_,
integer  namtot_,
integer  ncmat_,
integer  mpcfree,
integer  maxlenmpc,
integer  ne1d,
integer  ne2d,
integer  nflow,
integer  nlabel,
integer  iplas,
integer  nkon,
integer  ithermal,
integer  nmethod,
integer, dimension(*)  iperturb,
integer  nstate_,
integer  nener,
character*81, dimension(*)  set,
integer, dimension(*)  istartset,
integer, dimension(*)  iendset,
integer, dimension(*)  ialset,
real*8, dimension(*)  co,
integer, dimension(*)  kon,
integer, dimension(*)  ipkon,
character*8, dimension(*)  lakon,
integer, dimension(*)  nodeboun,
integer, dimension(*)  ndirboun,
integer, dimension(*)  iamboun,
real*8, dimension(*)  xboun,
integer, dimension(*)  ikboun,
integer, dimension(*)  ilboun,
integer, dimension(*)  ipompc,
integer, dimension(*)  nodempc,
real*8, dimension(*)  coefmpc,
character*20, dimension(*)  labmpc,
integer, dimension(*)  ikmpc,
integer, dimension(*)  ilmpc,
integer, dimension(*)  nodeforc,
integer, dimension(*)  ndirforc,
integer, dimension(*)  iamforc,
real*8, dimension(*)  xforc,
integer, dimension(*)  ikforc,
integer, dimension(*)  ilforc,
integer, dimension(*)  nelemload,
integer, dimension(*)  iamload,
character*20, dimension(*)  sideload,
real*8, dimension(*)  xload,
real*8, dimension(*)  elcon,
integer, dimension(*)  nelcon,
real*8, dimension(*)  rhcon,
integer, dimension(*)  nrhcon,
real*8, dimension(*)  alcon,
integer, dimension(*)  nalcon,
real*8, dimension(*)  alzero,
real*8, dimension(*)  plicon,
integer, dimension(*)  nplicon,
real*8, dimension(*)  plkcon,
integer, dimension(*)  nplkcon,
character*80, dimension(*)  orname,
real*8, dimension(*)  orab,
integer, dimension(*)  ielorien,
real*8, dimension(*)  trab,
integer, dimension(*)  inotr,
character*80, dimension(*)  amname,
real*8, dimension(*)  amta,
integer, dimension(*)  namta,
real*8, dimension(*)  t0,
real*8, dimension(*)  t1,
integer, dimension(*)  iamt1,
real*8, dimension(*)  veold,
integer, dimension(*)  ielmat,
character*80, dimension(*)  matname,
character*6, dimension(*)  prlab,
character*81, dimension(*)  prset,
character*87, dimension(*)  filab,
real*8, dimension(*)  vold,
integer, dimension(*)  nodebounold,
integer, dimension(*)  ndirbounold,
real*8, dimension(*)  xbounold,
real*8, dimension(*)  xforcold,
real*8, dimension(*)  xloadold,
real*8, dimension(*)  t1old,
real*8, dimension(*)  eme,
integer, dimension(*)  iponor,
real*8, dimension(*)  xnor,
integer, dimension(*)  knor,
real*8, dimension(*)  thicke,
real*8, dimension(*)  offset,
integer, dimension(*)  iponoel,
integer, dimension(*)  inoel,
integer, dimension(*)  rig,
real*8, dimension(*)  shcon,
integer, dimension(*)  nshcon,
real*8, dimension(*)  cocon,
integer, dimension(*)  ncocon,
integer, dimension(*)  ics,
real*8, dimension(*)  sti,
real*8, dimension(*)  ener,
real*8, dimension(*)  xstate,
character*132, dimension(*)  jobnamec,
integer, dimension(*)  infree,
integer  irstrt,
character*1, dimension(*)  inpc,
character*132, dimension(16)  textpart,
integer  istat,
integer  n,
integer  key,
real*8, dimension(*)  prestr,
integer  iprestr,
character*81, dimension(*)  cbody,
integer, dimension(*)  ibody,
real*8, dimension(*)  xbody,
integer  nbody,
real*8, dimension(*)  xbodyold,
real*8  ttime,
real*8, dimension(2)  qaold,
real*8, dimension(17,*)  cs,
integer  mcs,
character*3  output,
real*8, dimension(*)  physcon,
real*8, dimension(*)  ctrl,
character*1, dimension(*)  typeboun,
integer  iline,
integer  ipol,
integer  inl,
integer, dimension(2,*)  ipoinp,
integer, dimension(3,*)  inp,
real*8, dimension(*)  fmpc,
character*81, dimension(3,*)  tieset,
integer  ntie,
real*8, dimension(*)  tietol,
integer, dimension(0:*)  ipoinpc,
integer  nslavs,
real*8, dimension(*)  t0g,
real*8, dimension(*)  t1g,
integer  nprop,
integer, dimension(*)  ielprop,
real*8, dimension(*)  prop,
integer  mortar,
integer  nintpoint,
integer  ifacecount,
integer, dimension(*)  islavsurf,
real*8, dimension(*)  pslavsurf,
real*8, dimension(*)  clearini 
)
40 !
41  implicit none
42 !
43  character*1 typeboun(*),inpc(*)
44  character*3 output
45  character*6 prlab(*)
46  character*8 lakon(*)
47  character*20 labmpc(*),sideload(*)
48  character*80 orname(*),amname(*),matname(*)
49  character*81 set(*),prset(*),tieset(3,*),cbody(*)
50  character*87 filab(*)
51  character*132 jobnamec(*),textpart(16)
52 !
53  integer istep,nset,nload,nforc,nboun,nk,ne,nmpc,nalset,nmat,
54  & ntmat_,npmat_,norien,nam,nprint,mi(*),ntrans,ncs_,
55  & namtot_,ncmat_,mpcfree,ne1d,ne2d,nflow,nlabel,iplas,nkon,
56  & ithermal,nmethod,iperturb(*),nstate_,istartset(*),iendset(*),
57  & ialset(*),kon(*),ipkon(*),nodeboun(*),ndirboun(*),iamboun(*),
58  & ikboun(*),ilboun(*),ipompc(*),nodempc(*),ikmpc(*),ilmpc(*),
59  & nodeforc(*),ndirforc(*),iamforc(*),ikforc(*),ilforc(*),
60  & nelemload(*),iamload(*),nelcon(*),ipoinpc(0:*),
61  & nrhcon(*),nalcon(*),nplicon(*),nplkcon(*),ielorien(*),
62  & inotr(*),nprop,ielprop(*),mortar,nintpoint,ifacecount,
63  & namta(*),iamt1(*),ielmat(*),nodebounold(*),ndirbounold(*),
64  & iponor(*),knor(*),iponoel(*),inoel(*),rig(*),islavsurf(*),
65  & nshcon(*),ncocon(*),ics(*),infree(*),
66  & nener,irestartstep,irestartread,irstrt,istat,n,i,key,
67  & iprestr,mcs,maxlenmpc,iline,ipol,inl,
68  & ipoinp(2,*),inp(3,*),ntie,ibody(*),nbody,nslavs
69 !
70  real*8 co(*),xboun(*),coefmpc(*),xforc(*),xload(*),elcon(*),
71  & rhcon(*),alcon(*),alzero(*),plicon(*),plkcon(*),orab(*),
72  & trab(*),amta(*),t0(*),t1(*),prestr(*),veold(*),tietol(*),
73  & vold(*),xbounold(*),xforcold(*),xloadold(*),t1old(*),eme(*),
74  & xnor(*),thicke(*),offset(*),t0g(*),t1g(*),clearini(*),
75  & shcon(*),cocon(*),sti(*),ener(*),xstate(*),prop(*),
76  & ttime,qaold(2),cs(17,*),physcon(*),pslavsurf(*),
77  & ctrl(*),fmpc(*),xbody(*),xbodyold(*)
78 !
79  irestartread=0
80  irestartstep=0
81 !
82  do i=2,n
83  if(textpart(i)(1:4).eq.'READ') then
84  irestartread=1
85 c if(irestartstep.eq.0) irestartstep=1
86  elseif(textpart(i)(1:5).eq.'STEP=') then
87  read(textpart(i)(6:15),'(i10)',iostat=istat) irestartstep
88  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
89  &"*RESTART%")
90  elseif(textpart(i)(1:5).eq.'WRITE') then
91  irstrt=1
92  elseif(textpart(i)(1:10).eq.'FREQUENCY=') then
93  read(textpart(i)(11:20),'(i10)',iostat=istat) irstrt
94  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
95  &"*RESTART%")
96  else
97  write(*,*)
98  & '*WARNING in restarts: parameter not recognized:'
99  write(*,*) ' ',
100  & textpart(i)(1:index(textpart(i),' ')-1)
101  call inputwarning(inpc,ipoinpc,iline,
102  &"*RESTART%")
103  endif
104  enddo
105 !
106  if(irestartread.eq.1) then
107  call restartread(istep,nset,nload,nforc, nboun,nk,ne,
108  & nmpc,nalset,nmat,ntmat_,npmat_,norien,nam,nprint,mi,
109  & ntrans,ncs_,namtot_,ncmat_,mpcfree,maxlenmpc,
110  & ne1d,ne2d,nflow,nlabel,iplas,
111  & nkon,ithermal,nmethod,iperturb,nstate_,nener,set,istartset,
112  & iendset,ialset,co,kon,ipkon,lakon,nodeboun,ndirboun,iamboun,
113  & xboun,ikboun,ilboun,ipompc,nodempc,coefmpc,labmpc,ikmpc,ilmpc,
114  & nodeforc,ndirforc,iamforc,xforc,ikforc,ilforc,nelemload,iamload,
115  & sideload,xload,elcon,nelcon,rhcon,nrhcon,
116  & alcon,nalcon,alzero,plicon,nplicon,plkcon,nplkcon,orname,orab,
117  & ielorien,trab,inotr,amname,amta,namta,t0,t1,iamt1,veold,
118  & ielmat,matname,prlab,prset,filab,vold,nodebounold,
119  & ndirbounold,xbounold,xforcold,xloadold,t1old,eme,
120  & iponor,xnor,knor,thicke,offset,iponoel,inoel,rig,
121  & shcon,nshcon,cocon,ncocon,ics,sti,
122  & ener,xstate,jobnamec,infree,irestartstep,prestr,iprestr,
123  & cbody,ibody,xbody,nbody,xbodyold,ttime,qaold,cs,mcs,
124  & output,physcon,ctrl,typeboun,fmpc,tieset,ntie,tietol,nslavs,
125  & t0g,t1g,nprop,ielprop,prop,mortar,nintpoint,ifacecount,
126  & islavsurf,pslavsurf,clearini)
127  endif
128 !
129  call getnewline(inpc,textpart,istat,n,key,iline,ipol,inl,
130  & ipoinp,inp,ipoinpc)
131 !
132  return
subroutine inputwarning(inpc, ipoinpc, iline, text)
Definition: inputwarning.f:20
subroutine restartread(istep, nset, nload, nforc, nboun, nk, ne, nmpc, nalset, nmat, ntmat_, npmat_, norien, nam, nprint, mi, ntrans, ncs_, namtot_, ncmat_, mpcfree, maxlenmpc, ne1d, ne2d, nflow, nlabel, iplas, nkon, ithermal, nmethod, iperturb, nstate_, nener, set, istartset, iendset, ialset, co, kon, ipkon, lakon, nodeboun, ndirboun, iamboun, xboun, ikboun, ilboun, ipompc, nodempc, coefmpc, labmpc, ikmpc, ilmpc, nodeforc, ndirforc, iamforc, xforc, ikforc, ilforc, nelemload, iamload, sideload, xload, elcon, nelcon, rhcon, nrhcon, alcon, nalcon, alzero, plicon, nplicon, plkcon, nplkcon, orname, orab, ielorien, trab, inotr, amname, amta, namta, t0, t1, iamt1, veold, ielmat, matname, prlab, prset, filab, vold, nodebounold, ndirbounold, xbounold, xforcold, xloadold, t1old, eme, iponor, xnor, knor, thicke, offset, iponoel, inoel, rig, shcon, nshcon, cocon, ncocon, ics, sti, ener, xstate, jobnamec, infree, irestartstep, prestr, iprestr, cbody, ibody, xbody, nbody, xbodyold, ttime, qaold, cs, mcs, output, physcon, ctrl, typeboun, fmpc, tieset, ntie, tietol, nslavs, t0g, t1g, nprop, ielprop, prop, mortar, nintpoint, ifacecount, islavsurf, pslavsurf, clearini)
Definition: restartread.f:39
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)