FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
fstr_solve_dynamic.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 !-------------------------------------------------------------------------------
6 
8 
9  use m_fstr
12  use fstr_frequency_analysis !Frequency analysis module
13 
14 contains
15 
16  !C================================================================C
18  !C================================================================C
19  subroutine fstr_solve_dynamic(hecMESH,hecMAT,fstrSOLID,fstrEIG &
20  ,fstrDYNAMIC,fstrRESULT,fstrPARAM &
21  ,fstrCPL,fstrFREQ, hecLagMAT &
22  ,conMAT )
24  implicit none
25  type(hecmwst_local_mesh) :: hecMESH
26  type(hecmwst_matrix) :: hecMAT
27  type(fstr_eigen) :: fstrEIG
28  type(fstr_solid) :: fstrSOLID
29  type(hecmwst_result_data) :: fstrRESULT
30  type(fstr_param) :: fstrPARAM
31  type(fstr_dynamic) :: fstrDYNAMIC
32  type(fstr_couple) :: fstrCPL !for COUPLE
33  type(fstr_freqanalysis) :: fstrFREQ
34  type(hecmwst_matrix_lagrange) :: hecLagMAT
35  type(fstr_info_contactchange) :: infoCTChange
36  type(hecmwst_matrix) :: conMAT
37  integer(kind=kint) :: i, j, num_monit, ig, is, iE, ik, in, ing, iunitS, iunit, ierror, flag, limit
38  character(len=HECMW_FILENAME_LEN) :: fname, header
39  integer(kind=kint) :: restrt_step_num, ndof
40  integer(kind=kint) :: restrt_step(1)
41 
42  num_monit = 0
43 
44  if(dabs(fstrdynamic%t_delta) < 1.0e-20) then
45  if( hecmesh%my_rank == 0 ) then
46  write(imsg,*) 'stop due to fstrDYNAMIC%t_delta = 0'
47  end if
48  call hecmw_abort( hecmw_comm_get_comm())
49  end if
50  call fstr_dynamic_alloc( hecmesh, fstrdynamic )
51 
52  !C-- file open for local use
53  !C
54  if(fstrdynamic%idx_resp == 1) then ! time history analysis
55  call hecmw_ctrl_is_subdir( flag, limit )
56  if( flag == 0)then
57  header = ""
58  else
59  header = adjustl("MONITOR/")
60  call hecmw_ctrl_make_subdir( adjustl("MONITOR/test.txt"), limit )
61  endif
62  ig = fstrdynamic%ngrp_monit
63  is = hecmesh%node_group%grp_index(ig-1)+1
64  ie = hecmesh%node_group%grp_index(ig)
65  do ik=is,ie
66  in = hecmesh%node_group%grp_item(ik)
67  if (in > hecmesh%nn_internal) cycle
68  num_monit = num_monit+1
69  ing = hecmesh%global_node_id(in)
70  iunits = 10*(num_monit-1)
71 
72  iunit = iunits + fstrdynamic%dynamic_IW4
73  write(fname,'(a,i0,a)') trim(header)//'dyna_disp_',ing,'.txt'
74  if(fstrdynamic%restart_nout < 0 ) then
75  open(iunit,file=fname, position = 'append', iostat=ierror)
76  else
77  open(iunit,file=fname, status = 'replace', iostat=ierror)
78  endif
79  if( ierror /= 0 ) then
80  write(*,*) 'stop due to file opening error:',trim(fname)
81  call hecmw_abort( hecmw_comm_get_comm())
82  end if
83 
84  iunit = iunits + fstrdynamic%dynamic_IW5
85  write(fname,'(a,i0,a)') trim(header)//'dyna_velo_',ing,'.txt'
86  if(fstrdynamic%restart_nout < 0 ) then
87  open(iunit,file=fname, position = 'append', iostat=ierror)
88  else
89  open(iunit,file=fname, status = 'replace', iostat=ierror)
90  endif
91  if( ierror /= 0 ) then
92  write(*,*) 'stop due to file opening error',trim(fname)
93  call hecmw_abort( hecmw_comm_get_comm())
94  end if
95 
96  iunit = iunits + fstrdynamic%dynamic_IW6
97  write(fname,'(a,i0,a)') trim(header)//'dyna_acce_',ing,'.txt'
98  if(fstrdynamic%restart_nout < 0 ) then
99  open(iunit,file=fname, position = 'append', iostat=ierror)
100  else
101  open(iunit,file=fname, status = 'replace', iostat=ierror)
102  endif
103  if( ierror /= 0 ) then
104  write(*,*) 'stop due to file opening error',trim(fname)
105  call hecmw_abort( hecmw_comm_get_comm())
106  end if
107 
108  iunit = iunits + fstrdynamic%dynamic_IW7
109  write(fname,'(a,i0,a)') trim(header)//'dyna_force_',ing,'.txt'
110  if(fstrdynamic%restart_nout < 0 ) then
111  open(iunit,file=fname, position = 'append', iostat=ierror)
112  else
113  open(iunit,file=fname, status = 'replace', iostat=ierror)
114  endif
115  if( ierror /= 0 ) then
116  write(*,*) 'stop due to file opening error',trim(fname)
117  call hecmw_abort( hecmw_comm_get_comm())
118  end if
119  iunit = iunits + fstrdynamic%dynamic_IW8
120  write(fname,'(a,i0,a)') trim(header)//'dyna_strain_',ing,'.txt'
121  if(fstrdynamic%restart_nout < 0 ) then
122  open(iunit,file=fname, position = 'append', iostat=ierror)
123  else
124  open(iunit,file=fname, status = 'replace', iostat=ierror)
125  endif
126  if( ierror /= 0 ) then
127  write(*,*) 'stop due to file opening error',trim(fname)
128  call hecmw_abort( hecmw_comm_get_comm())
129  end if
130 
131  iunit = iunits + fstrdynamic%dynamic_IW9
132  write(fname,'(a,i0,a)') trim(header)//'dyna_stress_',ing,'.txt'
133  if(fstrdynamic%restart_nout < 0 ) then
134  open(iunit,file=fname, position = 'append', iostat=ierror)
135  else
136  open(iunit,file=fname, status = 'replace', iostat=ierror)
137  endif
138  if( ierror /= 0 ) then
139  write(*,*) 'stop due to file opening error',trim(fname)
140  call hecmw_abort( hecmw_comm_get_comm())
141  end if
142  enddo
143  endif
144 
145  !C
146  !C-- initial condition
147  !C
148  fstrdynamic%DISP = 0.d0
149  fstrdynamic%VEL = 0.d0
150  fstrdynamic%ACC = 0.d0
151  fstrsolid%unode(:) =0.d0
152  fstrsolid%QFORCE(:) =0.d0
153  fstrdynamic%kineticEnergy=0.d0
154  fstrdynamic%strainEnergy=0.d0
155  fstrdynamic%totalEnergy=0.d0
156  call fstr_updatestate( hecmesh, fstrsolid, 0.d0 )
157 
158  ! ---- restart
159 
160  restrt_step_num = 1
161  fstrdynamic%i_step = 0
162  infoctchange%contactNode_previous = 0
163 
164  if(associated(g_initialcnd))then
165  ndof = hecmat%NDOF
166  do j = 1, size(g_initialcnd)
167  if(g_initialcnd(j)%cond_name == "velocity")then
168  do i= 1, hecmesh%n_node
169  ing = g_initialcnd(j)%intval(i)
170  if(ing <= 0) cycle
171  fstrdynamic%VEL(ndof*i-(ndof-ing),1) = g_initialcnd(j)%realval(i)
172  end do
173  elseif(g_initialcnd(j)%cond_name == "acceleration")then
174  do i = 1, hecmesh%n_node
175  ing = g_initialcnd(j)%intval(i)
176  if(ing <= 0) cycle
177  fstrdynamic%ACC(ndof*i-(ndof-ing),1) = g_initialcnd(j)%realval(i)
178  enddo
179  endif
180  enddo
181  endif
182 
183  if(fstrdynamic%restart_nout >= 0 ) then
184  call dynamic_bc_init (hecmesh, hecmat, fstrsolid, fstrdynamic)
185  call dynamic_bc_init_vl(hecmesh, hecmat, fstrsolid, fstrdynamic)
186  call dynamic_bc_init_ac(hecmesh, hecmat, fstrsolid, fstrdynamic)
187  endif
188 
189  !restart
190  if(fstrdynamic%restart_nout < 0 ) then
191  if( fstrdynamic%idx_eqa == 1 ) then
192  call fstr_read_restart_dyna_nl(restrt_step_num,hecmesh,fstrsolid,fstrdynamic,fstrparam,&
193  infoctchange%contactNode_previous)
194  elseif(fstrdynamic%idx_eqa == 11) then
195  if( .not. associated( fstrsolid%contacts ) ) then
196  call fstr_read_restart_dyna_nl(restrt_step_num,hecmesh,fstrsolid,fstrdynamic,fstrparam)
197  else
198  call fstr_read_restart_dyna_nl(restrt_step_num,hecmesh,fstrsolid,fstrdynamic,fstrparam,&
199  infoctchange%contactNode_previous)
200  endif
201  endif
202  restrt_step_num = restrt_step_num + 1
203  fstrdynamic%restart_nout = - fstrdynamic%restart_nout
204  hecmat%Iarray(98) = 1
205  end if
206 
207  if(fstrdynamic%idx_resp == 1) then ! time history analysis
208 
209  if(fstrdynamic%idx_eqa == 1) then ! implicit dynamic analysis
210  call fstr_solve_dynamic_nlimplicit_contactslag(1, hecmesh,hecmat,fstrsolid,fstreig &
211  ,fstrdynamic,fstrresult,fstrparam &
212  ,fstrcpl,heclagmat,restrt_step_num,infoctchange &
213  ,conmat )
214 
215  else if(fstrdynamic%idx_eqa == 11) then ! explicit dynamic analysis
216  call fstr_solve_dynamic_nlexplicit(hecmesh,hecmat,fstrsolid,fstreig &
217  ,fstrdynamic,fstrresult,fstrparam,infoctchange &
218  ,fstrcpl, restrt_step_num )
219  endif
220 
221  else if(fstrdynamic%idx_resp == 2) then ! frequency response analysis
222 
223  if( fstrparam%nlgeom ) then
224  if( hecmesh%my_rank .eq. 0 ) then
225  write(imsg,*) 'stop: steady-state harmonic response analysis is not yet available !'
226  end if
227  call hecmw_abort( hecmw_comm_get_comm())
228  end if
229 
230  if( hecmesh%my_rank .eq. 0 ) then
231  call fstr_solve_frequency_analysis(hecmesh, hecmat, fstrsolid, fstreig, fstrdynamic, &
232  fstrresult, fstrparam, fstrcpl, fstrfreq, heclagmat, &
233  restrt_step_num)
234  end if
235  end if
236 
237  !C-- file close for local use
238  if(fstrdynamic%idx_resp == 1) then ! time history analysis
239  do i=1,num_monit
240  iunits = 10*(i-1)
241  close(iunits + fstrdynamic%dynamic_IW4)
242  close(iunits + fstrdynamic%dynamic_IW5)
243  close(iunits + fstrdynamic%dynamic_IW6)
244  close(iunits + fstrdynamic%dynamic_IW7)
245  close(iunits + fstrdynamic%dynamic_IW8)
246  close(iunits + fstrdynamic%dynamic_IW9)
247  enddo
248  endif
249  !C-- end of finalization
250 
251  call fstr_dynamic_finalize( fstrdynamic )
252  call hecmat_finalize( hecmat )
253 
254  deallocate( fstreig%mass )
255 
256  end subroutine fstr_solve_dynamic
257 
258 end module fstr_solver_dynamic
fstr_solver_dynamic
This module contains subroutines controlling dynamic calculation.
Definition: fstr_solve_dynamic.f90:7
hecmw_ctrl_is_subdir
void hecmw_ctrl_is_subdir(int *flag, int *limit)
Definition: hecmw_control.c:2619
m_fstr_setup
This module provides functions to read in data from control file and do necessary preparation for fol...
Definition: fstr_setup.f90:7
m_fstr::fstr_eigen
Package of data used by Lanczos eigenvalue solver.
Definition: m_fstr.f90:593
fstr_frequency_analysis::fstr_solve_frequency_analysis
subroutine fstr_solve_frequency_analysis(hecMESH, hecMAT, fstrSOLID, fstrEIG, fstrDYNAMIC, fstrRESULT, fstrPARAM, fstrCPL, fstrFREQ, hecLagMAT, restart_step_num)
Definition: fstr_frequency_analysis.f90:71
fstr_dynamic_nlexplicit
This module contains subroutines for nonlinear explicit dynamic analysis.
Definition: fstr_dynamic_nlexplicit.f90:7
fstr_frequency_analysis
Definition: fstr_frequency_analysis.f90:55
m_fstr::fstr_solid
Definition: m_fstr.f90:238
fstr_dynamic_nlimplicit
This module contains subroutines for nonlinear implicit dynamic analysis.
Definition: fstr_dynamic_nlimplicit.f90:7
m_fstr::fstr_dynamic
Data for DYNAMIC ANSLYSIS (fstrDYNAMIC)
Definition: m_fstr.f90:504
m_fstr_setup::fstr_dynamic_alloc
subroutine fstr_dynamic_alloc(hecMESH, fstrDYNAMIC)
Initial setting of dynamic calculation.
Definition: fstr_setup.f90:1585
m_fstr::fstr_param
FSTR INNER CONTROL PARAMETERS (fstrPARAM)
Definition: m_fstr.f90:154
m_fstr::fstr_freqanalysis
Definition: m_fstr.f90:571
m_fstr
This module defines common data and basic structures for analysis.
Definition: m_fstr.f90:15
hecmw_ctrl_make_subdir
void hecmw_ctrl_make_subdir(char *filename, int *err, int len)
Definition: hecmw_control.c:2594
m_fstr::g_initialcnd
type(tinitialcondition), dimension(:), pointer, save g_initialcnd
Definition: m_fstr.f90:151
m_fstr::fstr_couple
Data for coupling analysis.
Definition: m_fstr.f90:611
m_fstr_setup::fstr_dynamic_finalize
subroutine fstr_dynamic_finalize(fstrDYNAMIC)
Finalizer of fstr_solid.
Definition: fstr_setup.f90:1666
fstr_dynamic_nlexplicit::fstr_solve_dynamic_nlexplicit
subroutine fstr_solve_dynamic_nlexplicit(hecMESH, hecMAT, fstrSOLID, fstrEIG, fstrDYN, fstrRESULT, fstrPARAM, infoCTChange, fstrCPL, restrt_step_num)
Definition: fstr_dynamic_nlexplicit.f90:30
fstr_dynamic_nlimplicit::fstr_solve_dynamic_nlimplicit_contactslag
subroutine fstr_solve_dynamic_nlimplicit_contactslag(cstep, hecMESH, hecMAT, fstrSOLID, fstrEIG, fstrDYNAMIC, fstrRESULT, fstrPARAM, fstrCPL, hecLagMAT, restrt_step_num, infoCTChange, conMAT)
This subroutine provides function of nonlinear implicit dynamic analysis using the Newmark method....
Definition: fstr_dynamic_nlimplicit.f90:339
fstr_solver_dynamic::fstr_solve_dynamic
subroutine fstr_solve_dynamic(hecMESH, hecMAT, fstrSOLID, fstrEIG, fstrDYNAMIC, fstrRESULT, fstrPARAM, fstrCPL, fstrFREQ, hecLagMAT, conMAT)
Master subroutine for dynamic analysis.
Definition: fstr_solve_dynamic.f90:23
m_fstr::hecmat_finalize
subroutine hecmat_finalize(hecMAT)
Definition: m_fstr.f90:922
m_fstr::imsg
integer(kind=kint), parameter imsg
Definition: m_fstr.f90:110