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

Go to the source code of this file.

Functions/Subroutines

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)
 

Function/Subroutine Documentation

◆ restartread()

subroutine restartread ( 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  irestartstep,
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(*)  cs,
integer  mcs,
character*3  output,
real*8, dimension(*)  physcon,
real*8, dimension(*)  ctrl,
character*1, dimension(*)  typeboun,
real*8, dimension(*)  fmpc,
character*81, dimension(*)  tieset,
integer  ntie,
real*8, dimension(*)  tietol,
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 
)
39 !
40  implicit none
41 !
42  character*1 typeboun(*)
43  character*3 output
44  character*6 prlab(*)
45  character*8 lakon(*)
46  character*20 labmpc(*),sideload(*)
47  character*80 orname(*),amname(*),matname(*),version
48  character*81 set(*),prset(*),tieset(*),cbody(*)
49  character*87 filab(*)
50  character*132 fnrstrt,jobnamec(*)
51 !
52  integer istep,nset,nload,nforc,nboun,nk,ne,nmpc,nalset,nmat,
53  & ntmat_,npmat_,norien,nam,nprint,mi(*),ntrans,ncs_,
54  & namtot_,ncmat_,mpcfree,ne1d,ne2d,nflow,nlabel,iplas,nkon,
55  & ithermal,nmethod,iperturb(*),nstate_,istartset(*),iendset(*),
56  & ialset(*),kon(*),ipkon(*),nodeboun(*),ndirboun(*),iamboun(*),
57  & ikboun(*),ilboun(*),ipompc(*),nodempc(*),ikmpc(*),ilmpc(*),
58  & nodeforc(*),ndirforc(*),iamforc(*),ikforc(*),ilforc(*),
59  & nelemload(*),iamload(*),nelcon(*),mt,nprop,ielprop(*),
60  & nrhcon(*),nalcon(*),nplicon(*),nplkcon(*),ielorien(*),
61  & inotr(*),mortar,nintpoint,ifacecount,islavsurf(*),
62  & namta(*),iamt1(*),ielmat(*),nodebounold(*),ndirbounold(*),
63  & iponor(*),knor(*),iponoel(*),inoel(*),rig(*),
64  & nshcon(*),ncocon(*),ics(*),infree(*),i,ipos,
65  & nener,irestartstep,istat,iprestr,
66  & maxlenmpc,mcs,mpcend,ntie,ibody(*),nbody,nslavs
67 !
68  real*8 co(*),xboun(*),coefmpc(*),xforc(*),xload(*),elcon(*),
69  & rhcon(*),alcon(*),alzero(*),plicon(*),plkcon(*),orab(*),
70  & trab(*),amta(*),t0(*),t1(*),veold(*),tietol(*),pslavsurf(*),
71  & vold(*),xbounold(*),xforcold(*),xloadold(*),t1old(*),eme(*),
72  & xnor(*),thicke(*),offset(*),t0g(*),t1g(*),clearini(*),
73  & shcon(*),cocon(*),sti(*),ener(*),xstate(*),prestr(*),ttime,
74  & qaold(2),physcon(*),ctrl(*),cs(*),fmpc(*),xbody(*),
75  & xbodyold(*),prop(*)
76 !
77  ipos=index(jobnamec(1),char(0))
78  fnrstrt(1:ipos-1)=jobnamec(1)(1:ipos-1)
79  fnrstrt(ipos:ipos+3)=".rin"
80  do i=ipos+4,132
81  fnrstrt(i:i)=' '
82  enddo
83 !
84  open(15,file=fnrstrt,access='SEQUENTIAL',form='UNFORMATTED',
85  & err=15)
86 !
87  do
88 !
89  read(15,iostat=istat) version
90 !
91  if(istat.lt.0) then
92  if(irestartstep.eq.0) then
93 !
94 ! reading the last step
95 !
96  irestartstep=istep
97  close(15)
98  open(15,file=fnrstrt,access='SEQUENTIAL',
99  & form='UNFORMATTED',err=15)
100  read(15) version
101  else
102  write(*,*) '*ERROR reading *RESTART,READ: requested step'
103  write(*,*) ' is not in the restart file'
104  call exit(201)
105  endif
106  endif
107  read(15)istep
108 !
109 ! set size
110 !
111  read(15)nset
112  read(15)nalset
113 !
114 ! load size
115 !
116  read(15)nload
117  read(15)nbody
118  read(15)nforc
119  read(15)nboun
120  read(15)nflow
121 !
122 ! mesh size
123 !
124  read(15)nk
125  read(15)ne
126  read(15)nkon
127  read(15)(mi(i),i=1,3)
128  mt=mi(2)+1
129 !
130 ! constraint size
131 !
132  read(15)nmpc
133  read(15)mpcend
134  read(15)maxlenmpc
135 !
136 ! material size
137 !
138  read(15)nmat
139  read(15)ntmat_
140  read(15)npmat_
141  read(15)ncmat_
142 !
143 ! property info
144 !
145  read(15)nprop
146 !
147 ! transformation size
148 !
149  read(15)norien
150  read(15)ntrans
151 !
152 ! amplitude size
153 !
154  read(15)nam
155  read(15)namtot_
156 !
157 ! print size
158 !
159  read(15)nprint
160  read(15)nlabel
161 !
162 ! tie size
163 !
164  read(15)ntie
165 !
166 ! cyclic symmetry size
167 !
168  read(15)ncs_
169  read(15)mcs
170 !
171 ! 1d and 2d element size
172 !
173  read(15)ne1d
174  read(15)ne2d
175  read(15)(infree(i),i=1,4)
176 !
177 ! procedure info
178 !
179  read(15)nmethod
180  read(15)(iperturb(i),i=1,2)
181  read(15)nener
182  read(15)iplas
183  read(15)ithermal
184  read(15)nstate_
185  read(15)nslavs
186  read(15)iprestr
187  read(15)mortar
188  if(mortar.eq.1) then
189  read(15)ifacecount
190  read(15)nintpoint
191  endif
192 !
193  if(istep.eq.irestartstep) exit
194 !
195 ! skipping the next entries until the requested step is found
196 !
197  call skip(nset,nalset,nload,nbody,nforc,nboun,nk,ne,nkon,
198  & mi(1),nmpc,mpcend,nmat,ntmat_,npmat_,ncmat_,norien,ntrans,
199  & nam,nprint,nlabel,ncs_,ne1d,ne2d,infree,nmethod,
200  & iperturb,nener,ithermal,nstate_,iprestr,mcs,ntie,
201  & nslavs,nprop,mortar,ifacecount,nintpoint)
202 !
203  enddo
204 !
205 ! sets
206 !
207  read(15)(set(i),i=1,nset)
208  read(15)(istartset(i),i=1,nset)
209  read(15)(iendset(i),i=1,nset)
210  do i=1,nalset
211  read(15)ialset(i)
212  enddo
213 !
214 ! mesh
215 !
216  read(15)(co(i),i=1,3*nk)
217  read(15)(kon(i),i=1,nkon)
218  read(15)(ipkon(i),i=1,ne)
219  read(15)(lakon(i),i=1,ne)
220 !
221 ! single point constraints
222 !
223  read(15)(nodeboun(i),i=1,nboun)
224  read(15)(ndirboun(i),i=1,nboun)
225  read(15)(typeboun(i),i=1,nboun)
226  read(15)(xboun(i),i=1,nboun)
227  read(15)(ikboun(i),i=1,nboun)
228  read(15)(ilboun(i),i=1,nboun)
229  if(nam.gt.0) read(15)(iamboun(i),i=1,nboun)
230  read(15)(nodebounold(i),i=1,nboun)
231  read(15)(ndirbounold(i),i=1,nboun)
232  read(15)(xbounold(i),i=1,nboun)
233 !
234 ! multiple point constraints
235 !
236  read(15)(ipompc(i),i=1,nmpc)
237  read(15)(labmpc(i),i=1,nmpc)
238  read(15)(ikmpc(i),i=1,nmpc)
239  read(15)(ilmpc(i),i=1,nmpc)
240  read(15)(fmpc(i),i=1,nmpc)
241  read(15)(nodempc(i),i=1,3*mpcend)
242  read(15)(coefmpc(i),i=1,mpcend)
243  mpcfree=mpcend+1
244 !
245 ! point forces
246 !
247  read(15)(nodeforc(i),i=1,2*nforc)
248  read(15)(ndirforc(i),i=1,nforc)
249  read(15)(xforc(i),i=1,nforc)
250  read(15)(ikforc(i),i=1,nforc)
251  read(15)(ilforc(i),i=1,nforc)
252  if(nam.gt.0) read(15)(iamforc(i),i=1,nforc)
253  read(15)(xforcold(i),i=1,nforc)
254 !
255 ! distributed loads
256 !
257  read(15)(nelemload(i),i=1,2*nload)
258  read(15)(sideload(i),i=1,nload)
259  read(15)(xload(i),i=1,2*nload)
260  if(nam.gt.0) read(15)(iamload(i),i=1,2*nload)
261  read(15)(xloadold(i),i=1,2*nload)
262  read(15)(cbody(i),i=1,nbody)
263  read(15)(ibody(i),i=1,3*nbody)
264  read(15)(xbody(i),i=1,7*nbody)
265  read(15)(xbodyold(i),i=1,7*nbody)
266 !
267 ! prestress
268 !
269  if(iprestr.gt.0) read(15) (prestr(i),i=1,6*mi(1)*ne)
270 !
271 ! labels
272 !
273  read(15)(prlab(i),i=1,nprint)
274  read(15)(prset(i),i=1,nprint)
275  read(15)(filab(i),i=1,nlabel)
276 !
277 ! elastic constants
278 !
279  read(15)(elcon(i),i=1,(ncmat_+1)*ntmat_*nmat)
280  read(15)(nelcon(i),i=1,2*nmat)
281 !
282 ! density
283 !
284  read(15)(rhcon(i),i=1,2*ntmat_*nmat)
285  read(15)(nrhcon(i),i=1,nmat)
286 !
287 ! specific heat
288 !
289  read(15)(shcon(i),i=1,4*ntmat_*nmat)
290  read(15)(nshcon(i),i=1,nmat)
291 !
292 ! conductivity
293 !
294  read(15)(cocon(i),i=1,7*ntmat_*nmat)
295  read(15)(ncocon(i),i=1,2*nmat)
296 !
297 ! expansion coefficients
298 !
299  read(15)(alcon(i),i=1,7*ntmat_*nmat)
300  read(15)(nalcon(i),i=1,2*nmat)
301  read(15)(alzero(i),i=1,nmat)
302 !
303 ! physical constants
304 !
305  read(15)(physcon(i),i=1,10)
306 !
307 ! plastic data
308 !
309  if(npmat_.ne.0)then
310  read(15)(plicon(i),i=1,(2*npmat_+1)*ntmat_*nmat)
311  read(15)(nplicon(i),i=1,(ntmat_+1)*nmat)
312  read(15)(plkcon(i),i=1,(2*npmat_+1)*ntmat_*nmat)
313  read(15)(nplkcon(i),i=1,(ntmat_+1)*nmat)
314  endif
315 !
316 ! material orientation
317 !
318  if(norien.ne.0)then
319  read(15)(orname(i),i=1,norien)
320  read(15)(orab(i),i=1,7*norien)
321  read(15)(ielorien(i),i=1,mi(3)*ne)
322  endif
323 !
324 ! fluid section properties
325 !
326  if(nprop.ne.0) then
327  read(15)(ielprop(i),i=1,ne)
328  read(15)(prop(i),i=1,nprop)
329  endif
330 !
331 ! transformations
332 !
333  if(ntrans.ne.0)then
334  read(15)(trab(i),i=1,7*ntrans)
335  read(15)(inotr(i),i=1,2*nk)
336  endif
337 !
338 ! amplitudes
339 !
340  if(nam.gt.0)then
341  read(15)(amname(i),i=1,nam)
342  read(15)(namta(i),i=1,3*nam-1)
343  read(15) namta(3*nam)
344  read(15)(amta(i),i=1,2*namta(3*nam-1))
345  endif
346 !
347 ! temperatures
348 !
349  if(ithermal.gt.0)then
350  read(15)(t0(i),i=1,nk)
351  read(15)(t1(i),i=1,nk)
352  if((ne1d.gt.0).or.(ne2d.gt.0))then
353  read(15)(t0g(i),i=1,2*nk)
354  read(15)(t1g(i),i=1,2*nk)
355  endif
356  if(nam.gt.0) read(15)(iamt1(i),i=1,nk)
357  read(15)(t1old(i),i=1,nk)
358  endif
359 !
360 ! materials
361 !
362  read(15)(matname(i),i=1,nmat)
363  read(15)(ielmat(i),i=1,mi(3)*ne)
364 !
365 ! temperature, displacement, static pressure, velocity and acceleration
366 !
367  read(15)(vold(i),i=1,mt*nk)
368  if((nmethod.eq.4).or.((nmethod.eq.1).and.(iperturb(1).ge.2))) then
369  read(15)(veold(i),i=1,mt*nk)
370  endif
371 !
372 ! 1d and 2d elements
373 !
374  if((ne1d.gt.0).or.(ne2d.gt.0))then
375  read(15)(iponor(i),i=1,2*nkon)
376  read(15)(xnor(i),i=1,infree(1))
377  read(15)(knor(i),i=1,infree(2))
378  read(15)(thicke(i),i=1,mi(3)*nkon)
379  read(15)(offset(i),i=1,2*ne)
380  read(15)(iponoel(i),i=1,infree(4))
381  read(15)(inoel(i),i=1,3*(infree(3)-1))
382  read(15)(rig(i),i=1,infree(4))
383  endif
384 !
385 ! tie constraints
386 !
387  if(ntie.gt.0) then
388  read(15)(tieset(i),i=1,3*ntie)
389  read(15)(tietol(i),i=1,3*ntie)
390  endif
391 !
392 ! cyclic symmetry
393 !
394  if(ncs_.gt.0)then
395  read(15)(ics(i),i=1,ncs_)
396  endif
397  if(mcs.gt.0) then
398  read(15)(cs(i),i=1,17*mcs)
399  endif
400 !
401 ! integration point variables
402 !
403  read(15)(sti(i),i=1,6*mi(1)*ne)
404  read(15)(eme(i),i=1,6*mi(1)*ne)
405  if(nener.eq.1) then
406  read(15)(ener(i),i=1,mi(1)*ne)
407  endif
408  if(nstate_.gt.0)then
409  if(mortar.eq.0) then
410  read(15)(xstate(i),i=1,nstate_*mi(1)*(ne+nslavs))
411  elseif(mortar.eq.1) then
412  read(15)(xstate(i),i=1,nstate_*mi(1)*(ne+nintpoint))
413  endif
414  endif
415 !
416 ! face-to-face penalty contact variables
417 !
418  if(mortar.eq.1) then
419  read(15) (islavsurf(i),i=1,2*ifacecount+2)
420  read(15) (pslavsurf(i),i=1,3*nintpoint)
421  read(15) (clearini(i),i=1,3*9*ifacecount)
422  endif
423 !
424 ! control parameters
425 !
426  read(15) (ctrl(i),i=1,39)
427  read(15) (qaold(i),i=1,2)
428  read(15) output
429  read(15) ttime
430 !
431  close(15)
432 !
433  return
434 !
435  15 write(*,*) '*ERROR reading *RESTART,READ: could not open file ',
436  & fnrstrt
437  call exit(201)
subroutine skip(nset, nalset, nload, nbody, nforc, nboun, nk, ne, nkon, mi, nmpc, memmpc_, nmat, ntmat_, npmat_, ncmat_, norien, ntrans, nam, nprint, nlabel, ncs_, ne1d, ne2d, infree, nmethod, iperturb, nener, ithermal, nstate_, iprestr, mcs, ntie, nslavs, nprop, mortar, ifacecount, nintpoint)
Definition: skip.f:25
Hosted by OpenAircraft.com, (Michigan UAV, LLC)