FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
fstr_Restart.f90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
2 ! Copyright (c) 2019 FrontISTR Commons
3 ! This software is released under the MIT License, see LICENSE.txt
4 !-------------------------------------------------------------------------------
7 
9  use m_utilities
10  use m_fstr
11  implicit none
12 
13 contains
14 
16  !----------------------------------------------------------------------*
17  subroutine fstr_read_restart(cstep,substep,step_count,ctime,dtime,hecMESH,fstrSOLID,fstrPARAM,contactNode)
18  !----------------------------------------------------------------------*
19  integer, intent(out) :: cstep
20  integer, intent(out) :: substep
21  integer, intent(out) :: step_count
22  real(kind=kreal), intent(out) :: ctime
23  real(kind=kreal), intent(out) :: dtime
24  integer, intent(out) :: contactNode
25  type (hecmwST_local_mesh), intent(in) :: hecMESH
26  type (fstr_solid),intent(inout) :: fstrSOLID
27  type(fstr_param), intent(in) :: fstrPARAM
28 
29  integer :: i,j,restrt_step(3),nif(2),istat(1),nload_prev(1),naux(2)
30  real(kind=kreal) :: times(3)
31 
32  call hecmw_restart_open()
33 
34  call hecmw_restart_read_int(restrt_step)
35  if( fstrparam%restart_version >= 5 ) then
36  if( myrank == 0 ) write(*,*) 'Reading restart file as new format(>=ver5.0)'
37  call hecmw_restart_read_real(times)
38  call hecmw_restart_read_int(fstrsolid%NRstat_i)
39  call hecmw_restart_read_real(fstrsolid%NRstat_r)
40  call hecmw_restart_read_int(istat)
41  else
42  if( myrank == 0 ) write(*,*) 'Reading restart file as old format(<ver5.0)'
43  endif
44  call hecmw_restart_read_int(nload_prev) !load info at previous step
45  if( nload_prev(1)>0 ) then
46  allocate(fstrsolid%step_ctrl_restart%Load(nload_prev(1)))
47  call hecmw_restart_read_int(fstrsolid%step_ctrl_restart%Load)
48  endif
49 
50  call hecmw_restart_read_real(fstrsolid%unode)
51  call hecmw_restart_read_real(fstrsolid%unode_bak)
52  call hecmw_restart_read_real(fstrsolid%QFORCE)
53 
54  do i= 1, hecmesh%n_elem
55  if (hecmw_is_etype_link( fstrsolid%elements(i)%etype )) cycle
56  if (hecmw_is_etype_patch( fstrsolid%elements(i)%etype )) cycle
57  do j= 1, size(fstrsolid%elements(i)%gausses)
58  call hecmw_restart_read_int(nif)
59  call hecmw_restart_read_real(fstrsolid%elements(i)%gausses(j)%strain)
60  call hecmw_restart_read_real(fstrsolid%elements(i)%gausses(j)%strain_bak)
61  call hecmw_restart_read_real(fstrsolid%elements(i)%gausses(j)%stress)
62  call hecmw_restart_read_real(fstrsolid%elements(i)%gausses(j)%stress_bak)
63  if( nif(1)>0 ) call hecmw_restart_read_int(fstrsolid%elements(i)%gausses(j)%istatus)
64  if( nif(2)>0 ) call hecmw_restart_read_real(fstrsolid%elements(i)%gausses(j)%fstatus)
65  enddo
66  call hecmw_restart_read_int(naux)
67  do j= 1, naux(2)
68  call hecmw_restart_read_real(fstrsolid%elements(i)%aux(:,j))
69  enddo
70  enddo
71 
72  if( associated( fstrsolid%contacts ) ) then
73  call hecmw_restart_read_int(nif)
74  contactnode = nif(1)
75  do i= 1, fstrsolid%n_contacts
76  do j= 1, size(fstrsolid%contacts(i)%slave)
77  call hecmw_restart_read_int(nif)
78  fstrsolid%contacts(i)%states(j)%surface = nif(1)
79  fstrsolid%contacts(i)%states(j)%state = nif(2)
80  call hecmw_restart_read_real(fstrsolid%contacts(i)%states(j)%lpos)
81  call hecmw_restart_read_real(fstrsolid%contacts(i)%states(j)%direction)
82  call hecmw_restart_read_real(fstrsolid%contacts(i)%states(j)%multiplier)
83  call hecmw_restart_read_real(fstrsolid%contacts(i)%states(j)%tangentForce_trial)
84  call hecmw_restart_read_real(fstrsolid%contacts(i)%states(j)%tangentForce_final)
85  enddo
86  enddo
87  endif
88 
89  call hecmw_restart_close()
90 
91  cstep = restrt_step(1)
92  substep = restrt_step(2) + 1
93  step_count = restrt_step(3)
94  if( fstrparam%restart_version >= 5 ) then
95  ctime = times(1)
96  dtime = times(2)
97  fstrsolid%AutoINC_stat = istat(1)
98  if( dabs(times(1)-times(3)) < 1.d-10 ) then
99  cstep = cstep + 1
100  substep = 1
101  endif
102  do i=1,size(fstrsolid%step_ctrl) !shift start time
103  fstrsolid%step_ctrl(i)%starttime = fstrsolid%step_ctrl(i)%starttime + times(3)
104  end do
105  else
106  ctime = fstrsolid%step_ctrl(cstep)%starttime
107  ctime = ctime + dble(substep-1)*fstrsolid%step_ctrl(cstep)%initdt
108  dtime = fstrsolid%step_ctrl(cstep)%initdt
109  if( dabs(ctime-fstrsolid%step_ctrl(cstep)%starttime-fstrsolid%step_ctrl(cstep)%elapsetime) < 1.d-10 ) then
110  cstep = cstep + 1
111  substep = 1
112  endif
113  endif
114 
115  end subroutine fstr_read_restart
116 
118  !----------------------------------------------------------------------*
119  subroutine fstr_write_restart(cstep,cstep_ext,substep,step_count,ctime,dtime,hecMESH, &
120  & fstrSOLID,fstrPARAM,is_StepFinished,contactNode)
121  !----------------------------------------------------------------------*
122  integer, intent(in) :: cstep
123  integer, intent(in) :: cstep_ext
124  integer, intent(in) :: substep
125  integer, intent(in) :: step_count
126  real(kind=kreal), intent(in) :: ctime
127  real(kind=kreal), intent(in) :: dtime
128  logical, intent(in) :: is_StepFinished
129  integer, intent(in) :: contactNode
130  type (hecmwST_local_mesh), intent(in) :: hecMESH
131  type (fstr_solid), intent(in) :: fstrSOLID
132  type(fstr_param), intent(in) :: fstrPARAM
133 
134  integer :: i,j,restrt_step(3),nif(2),istat(1),nload_prev(1),naux(2)
135  real(kind=kreal) :: times(3)
136 
137  restrt_step(1) = cstep_ext
138  restrt_step(2) = substep
139  restrt_step(3) = step_count
140  times(1) = ctime
141  times(2) = dtime
142  if( is_stepfinished ) then
143  times(3) = ctime
144  else
145  times(3) = fstrsolid%step_ctrl(cstep)%starttime
146  end if
147  istat(1) = fstrsolid%AutoINC_stat
148  call hecmw_restart_add_int(restrt_step,size(restrt_step))
149  if( fstrparam%restart_version >= 5 ) then
150  call hecmw_restart_add_real(times,size(times))
151  call hecmw_restart_add_int(fstrsolid%NRstat_i,size(fstrsolid%NRstat_i))
152  call hecmw_restart_add_real(fstrsolid%NRstat_r,size(fstrsolid%NRstat_r))
153  call hecmw_restart_add_int(istat,1)
154  endif
155  nload_prev(1) = 0
156  if( is_stepfinished ) then
157  if( associated(fstrsolid%step_ctrl(cstep)%Load) ) nload_prev(1) = size(fstrsolid%step_ctrl(cstep)%Load)
158  call hecmw_restart_add_int(nload_prev,1)
159  if( nload_prev(1)>0 ) call hecmw_restart_add_int(fstrsolid%step_ctrl(cstep)%Load,nload_prev(1))
160  else
161  if( cstep>1 ) then
162  if( associated(fstrsolid%step_ctrl(cstep-1)%Load) ) nload_prev(1) = size(fstrsolid%step_ctrl(cstep-1)%Load)
163  call hecmw_restart_add_int(nload_prev,1)
164  if( nload_prev(1)>0 ) call hecmw_restart_add_int(fstrsolid%step_ctrl(cstep-1)%Load,nload_prev(1))
165  else
166  if( associated(fstrsolid%step_ctrl_restart%Load) ) nload_prev(1) = size(fstrsolid%step_ctrl_restart%Load)
167  call hecmw_restart_add_int(nload_prev,1)
168  if( nload_prev(1)>0 ) call hecmw_restart_add_int(fstrsolid%step_ctrl_restart%Load,nload_prev(1))
169  endif
170  end if
171 
172  call hecmw_restart_add_real(fstrsolid%unode,size(fstrsolid%unode))
173  call hecmw_restart_add_real(fstrsolid%unode_bak,size(fstrsolid%unode_bak))
174  call hecmw_restart_add_real(fstrsolid%QFORCE,size(fstrsolid%QFORCE))
175 
176  do i= 1, hecmesh%n_elem
177  if (hecmw_is_etype_link( fstrsolid%elements(i)%etype )) cycle
178  if (hecmw_is_etype_patch( fstrsolid%elements(i)%etype )) cycle
179  do j= 1, size(fstrsolid%elements(i)%gausses)
180  nif = 0
181  if( associated(fstrsolid%elements(i)%gausses(j)%istatus) ) nif(1)=size(fstrsolid%elements(i)%gausses(j)%istatus)
182  if( associated(fstrsolid%elements(i)%gausses(j)%fstatus) ) nif(2)=size(fstrsolid%elements(i)%gausses(j)%fstatus)
183  call hecmw_restart_add_int(nif,size(nif))
184  call hecmw_restart_add_real(fstrsolid%elements(i)%gausses(j)%strain,size(fstrsolid%elements(i)%gausses(j)%strain))
185  call hecmw_restart_add_real(fstrsolid%elements(i)%gausses(j)%strain_bak,size(fstrsolid%elements(i)%gausses(j)%strain_bak))
186  call hecmw_restart_add_real(fstrsolid%elements(i)%gausses(j)%stress,size(fstrsolid%elements(i)%gausses(j)%stress))
187  call hecmw_restart_add_real(fstrsolid%elements(i)%gausses(j)%stress_bak,size(fstrsolid%elements(i)%gausses(j)%stress_bak))
188  if( nif(1)>0 ) then
189  call hecmw_restart_add_int(fstrsolid%elements(i)%gausses(j)%istatus,size(fstrsolid%elements(i)%gausses(j)%istatus))
190  endif
191  if( nif(2)>0 ) then
192  call hecmw_restart_add_real(fstrsolid%elements(i)%gausses(j)%fstatus,size(fstrsolid%elements(i)%gausses(j)%fstatus))
193  endif
194  enddo
195  naux = 0
196  if( associated(fstrsolid%elements(i)%aux) ) naux=shape(fstrsolid%elements(i)%aux)
197  call hecmw_restart_add_int(naux,size(naux))
198  do j= 1, naux(2)
199  call hecmw_restart_add_real(fstrsolid%elements(i)%aux(:,j),naux(1))
200  enddo
201  enddo
202 
203  if( associated( fstrsolid%contacts ) ) then
204  nif(1) = contactnode
205  call hecmw_restart_add_int(nif,size(nif))
206  do i= 1, fstrsolid%n_contacts
207  do j= 1, size(fstrsolid%contacts(i)%slave)
208  nif(1) = fstrsolid%contacts(i)%states(j)%surface
209  nif(2) = fstrsolid%contacts(i)%states(j)%state
210  call hecmw_restart_add_int(nif,size(nif))
211  call hecmw_restart_add_real(fstrsolid%contacts(i)%states(j)%lpos,size(fstrsolid%contacts(i)%states(j)%lpos))
212  call hecmw_restart_add_real(fstrsolid%contacts(i)%states(j)%direction,size(fstrsolid%contacts(i)%states(j)%direction))
213  call hecmw_restart_add_real(fstrsolid%contacts(i)%states(j)%multiplier,size(fstrsolid%contacts(i)%states(j)%multiplier))
214  call hecmw_restart_add_real(fstrsolid%contacts(i)%states(j)%tangentForce_trial, &
215  size(fstrsolid%contacts(i)%states(j)%tangentForce_trial))
216  call hecmw_restart_add_real(fstrsolid%contacts(i)%states(j)%tangentForce_final, &
217  size(fstrsolid%contacts(i)%states(j)%tangentForce_final))
218  enddo
219  enddo
220  endif
221 
222  call hecmw_restart_write()
223 
224  end subroutine fstr_write_restart
225 
227  !----------------------------------------------------------------------*
228  subroutine fstr_read_restart_dyna_nl(cstep,hecMESH,fstrSOLID,fstrDYNAMIC,fstrPARAM,contactNode)
229  !----------------------------------------------------------------------*
230  integer, intent(out) :: cstep
231  integer, intent(out), optional :: contactNode
232  type (hecmwST_local_mesh), intent(in) :: hecMESH
233  type (fstr_solid),intent(inout) :: fstrSOLID
234  type ( fstr_dynamic), intent(inout) :: fstrDYNAMIC
235  type(fstr_param), intent(in) :: fstrPARAM
236 
237  integer :: i,j,restrt_step(1),nif(2),naux(2)
238  real(kind=kreal) :: data(2)
239 
240  call hecmw_restart_open()
241 
242  call hecmw_restart_read_int(restrt_step)
243  cstep = restrt_step(1)
244  call hecmw_restart_read_real(fstrsolid%unode)
245  call hecmw_restart_read_real(fstrsolid%QFORCE)
246 
247  do i= 1, hecmesh%n_elem
248  if (hecmw_is_etype_link( fstrsolid%elements(i)%etype )) cycle
249  if (hecmw_is_etype_patch( fstrsolid%elements(i)%etype )) cycle
250  do j= 1, size(fstrsolid%elements(i)%gausses)
251  call hecmw_restart_read_int(nif)
252  call hecmw_restart_read_real(fstrsolid%elements(i)%gausses(j)%strain)
253  call hecmw_restart_read_real(fstrsolid%elements(i)%gausses(j)%strain_bak)
254  call hecmw_restart_read_real(fstrsolid%elements(i)%gausses(j)%stress)
255  call hecmw_restart_read_real(fstrsolid%elements(i)%gausses(j)%stress_bak)
256  if( nif(1)>0 ) call hecmw_restart_read_int(fstrsolid%elements(i)%gausses(j)%istatus)
257  if( nif(2)>0 ) call hecmw_restart_read_real(fstrsolid%elements(i)%gausses(j)%fstatus)
258  enddo
259  call hecmw_restart_read_int(naux)
260  do j= 1, naux(2)
261  call hecmw_restart_read_real(fstrsolid%elements(i)%aux(:,j))
262  enddo
263  enddo
264 
265  if(present(contactnode)) then
266  call hecmw_restart_read_int(nif)
267  contactnode = nif(1)
268  do i= 1, fstrsolid%n_contacts
269  do j= 1, size(fstrsolid%contacts(i)%slave)
270  call hecmw_restart_read_int(nif)
271  fstrsolid%contacts(i)%states(j)%surface = nif(1)
272  fstrsolid%contacts(i)%states(j)%state = nif(2)
273  call hecmw_restart_read_real(fstrsolid%contacts(i)%states(j)%lpos)
274  call hecmw_restart_read_real(fstrsolid%contacts(i)%states(j)%direction)
275  call hecmw_restart_read_real(fstrsolid%contacts(i)%states(j)%multiplier)
276  call hecmw_restart_read_real(fstrsolid%contacts(i)%states(j)%tangentForce_trial)
277  call hecmw_restart_read_real(fstrsolid%contacts(i)%states(j)%tangentForce_final)
278  enddo
279  enddo
280  endif
281 
282  call hecmw_restart_read_int(restrt_step)
283  fstrdynamic%idx_eqa = restrt_step(1)
284  call hecmw_restart_read_real(data)
285  fstrdynamic%t_curr = data(1)
286  fstrdynamic%strainEnergy = data(2)
287  if( fstrdynamic%idx_eqa == 1 ) then
288  call hecmw_restart_read_real(fstrdynamic%DISP(:,1))
289  call hecmw_restart_read_real(fstrdynamic%VEL(:,1))
290  call hecmw_restart_read_real(fstrdynamic%ACC(:,1))
291  else
292  call hecmw_restart_read_real(fstrdynamic%DISP(:,1))
293  call hecmw_restart_read_real(fstrdynamic%DISP(:,3))
294  endif
295  do i= 1, hecmesh%n_elem
296  if (hecmw_is_etype_link( fstrsolid%elements(i)%etype )) cycle
297  if (hecmw_is_etype_patch( fstrsolid%elements(i)%etype )) cycle
298  call hecmw_restart_read_real(fstrsolid%elements(i)%equiForces)
299  enddo
300 
301  call hecmw_restart_close()
302 
303  end subroutine fstr_read_restart_dyna_nl
304 
306  !----------------------------------------------------------------------*
307  subroutine fstr_write_restart_dyna_nl(cstep,hecMESH,fstrSOLID,fstrDYNAMIC,fstrPARAM,contactNode)
308  !----------------------------------------------------------------------*
309  integer, intent(in) :: cstep
310  integer, intent(in), optional :: contactNode
311  type (hecmwST_local_mesh), intent(in) :: hecMESH
312  type (fstr_solid), intent(in) :: fstrSOLID
313  type ( fstr_dynamic), intent(in) :: fstrDYNAMIC
314  type(fstr_param), intent(in) :: fstrPARAM
315 
316  integer :: i,j,restrt_step(1),nif(2),naux(2)
317  real(kind=kreal) :: data(2)
318 
319  restrt_step(1) = cstep
320  call hecmw_restart_add_int(restrt_step,size(restrt_step))
321  call hecmw_restart_add_real(fstrsolid%unode,size(fstrsolid%unode))
322  call hecmw_restart_add_real(fstrsolid%QFORCE,size(fstrsolid%QFORCE))
323 
324  do i= 1, hecmesh%n_elem
325  if (hecmw_is_etype_link( fstrsolid%elements(i)%etype )) cycle
326  if (hecmw_is_etype_patch( fstrsolid%elements(i)%etype )) cycle
327  do j= 1, size(fstrsolid%elements(i)%gausses)
328  nif = 0
329  if( associated(fstrsolid%elements(i)%gausses(j)%istatus) ) nif(1)=size(fstrsolid%elements(i)%gausses(j)%istatus)
330  if( associated(fstrsolid%elements(i)%gausses(j)%fstatus) ) nif(2)=size(fstrsolid%elements(i)%gausses(j)%fstatus)
331  call hecmw_restart_add_int(nif,size(nif))
332  call hecmw_restart_add_real(fstrsolid%elements(i)%gausses(j)%strain,size(fstrsolid%elements(i)%gausses(j)%strain))
333  call hecmw_restart_add_real(fstrsolid%elements(i)%gausses(j)%strain_bak,size(fstrsolid%elements(i)%gausses(j)%strain_bak))
334  call hecmw_restart_add_real(fstrsolid%elements(i)%gausses(j)%stress,size(fstrsolid%elements(i)%gausses(j)%stress))
335  call hecmw_restart_add_real(fstrsolid%elements(i)%gausses(j)%stress_bak,size(fstrsolid%elements(i)%gausses(j)%stress_bak))
336  if( nif(1)>0 ) then
337  call hecmw_restart_add_int(fstrsolid%elements(i)%gausses(j)%istatus,size(fstrsolid%elements(i)%gausses(j)%istatus))
338  endif
339  if( nif(2)>0 ) then
340  call hecmw_restart_add_real(fstrsolid%elements(i)%gausses(j)%fstatus,size(fstrsolid%elements(i)%gausses(j)%fstatus))
341  endif
342  enddo
343  naux = 0
344  if( associated(fstrsolid%elements(i)%aux) ) naux=shape(fstrsolid%elements(i)%aux)
345  call hecmw_restart_add_int(naux,size(naux))
346  do j= 1, naux(2)
347  call hecmw_restart_add_real(fstrsolid%elements(i)%aux(:,j),naux(1))
348  enddo
349  enddo
350 
351  if(present(contactnode)) then
352  nif(1) = contactnode
353  call hecmw_restart_add_int(nif,size(nif))
354  do i= 1, fstrsolid%n_contacts
355  do j= 1, size(fstrsolid%contacts(i)%slave)
356  nif(1) = fstrsolid%contacts(i)%states(j)%surface
357  nif(2) = fstrsolid%contacts(i)%states(j)%state
358  call hecmw_restart_add_int(nif,size(nif))
359  call hecmw_restart_add_real(fstrsolid%contacts(i)%states(j)%lpos,size(fstrsolid%contacts(i)%states(j)%lpos))
360  call hecmw_restart_add_real(fstrsolid%contacts(i)%states(j)%direction,size(fstrsolid%contacts(i)%states(j)%direction))
361  call hecmw_restart_add_real(fstrsolid%contacts(i)%states(j)%multiplier,size(fstrsolid%contacts(i)%states(j)%multiplier))
362  call hecmw_restart_add_real(fstrsolid%contacts(i)%states(j)%tangentForce_trial, &
363  size(fstrsolid%contacts(i)%states(j)%tangentForce_trial))
364  call hecmw_restart_add_real(fstrsolid%contacts(i)%states(j)%tangentForce_final, &
365  size(fstrsolid%contacts(i)%states(j)%tangentForce_final))
366  enddo
367  enddo
368  endif
369 
370  restrt_step(1) = fstrdynamic%idx_eqa
371  call hecmw_restart_add_int(restrt_step,size(restrt_step))
372  data(1) = fstrdynamic%t_curr
373  data(2) = fstrdynamic%strainEnergy
374  call hecmw_restart_add_real(data,size(data))
375  if( fstrdynamic%idx_eqa == 1 ) then
376  call hecmw_restart_add_real(fstrdynamic%DISP(:,1),size(fstrdynamic%DISP(:,1)))
377  call hecmw_restart_add_real(fstrdynamic%VEL(:,1),size(fstrdynamic%VEL(:,1)))
378  call hecmw_restart_add_real(fstrdynamic%ACC(:,1),size(fstrdynamic%ACC(:,1)))
379  else
380  call hecmw_restart_add_real(fstrdynamic%DISP(:,1),size(fstrdynamic%DISP(:,1)))
381  call hecmw_restart_add_real(fstrdynamic%DISP(:,3),size(fstrdynamic%DISP(:,3)))
382  endif
383  do i= 1, hecmesh%n_elem
384  if (hecmw_is_etype_link( fstrsolid%elements(i)%etype )) cycle
385  if (hecmw_is_etype_patch( fstrsolid%elements(i)%etype )) cycle
386  call hecmw_restart_add_real(fstrsolid%elements(i)%equiForces,size(fstrsolid%elements(i)%equiForces))
387  enddo
388 
389  call hecmw_restart_write()
390 
391  end subroutine fstr_write_restart_dyna_nl
392 
393 end module m_fstr_restart
m_utilities
This module provides aux functions.
Definition: utilities.f90:6
m_fstr_restart::fstr_read_restart_dyna_nl
subroutine fstr_read_restart_dyna_nl(cstep, hecMESH, fstrSOLID, fstrDYNAMIC, fstrPARAM, contactNode)
Read in restart file for nonlinear dynamic analysis.
Definition: fstr_Restart.f90:229
m_fstr::myrank
integer(kind=kint) myrank
PARALLEL EXECUTION.
Definition: m_fstr.f90:96
m_fstr_restart::fstr_read_restart
subroutine fstr_read_restart(cstep, substep, step_count, ctime, dtime, hecMESH, fstrSOLID, fstrPARAM, contactNode)
Read in restart file.
Definition: fstr_Restart.f90:18
m_fstr::fstr_param
FSTR INNER CONTROL PARAMETERS (fstrPARAM)
Definition: m_fstr.f90:154
m_fstr_restart::fstr_write_restart
subroutine fstr_write_restart(cstep, cstep_ext, substep, step_count, ctime, dtime, hecMESH, fstrSOLID, fstrPARAM, is_StepFinished, contactNode)
write out restart file
Definition: fstr_Restart.f90:121
m_fstr
This module defines common data and basic structures for analysis.
Definition: m_fstr.f90:15
m_fstr_restart::fstr_write_restart_dyna_nl
subroutine fstr_write_restart_dyna_nl(cstep, hecMESH, fstrSOLID, fstrDYNAMIC, fstrPARAM, contactNode)
write out restart file for nonlinear dynamic analysis
Definition: fstr_Restart.f90:308
m_fstr_restart
This module provides functions to read in and write out restart files.
Definition: fstr_Restart.f90:8