FrontISTR  5.9.0
Large-scale structural analysis program with finit element method
fstr_Cutback.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  use m_fstr
9  use elementinfo
10  use mmechgauss
11  use mcontactdef
12 
13  implicit none
14 
15  logical, private :: is_cutback_active = .false.
16 
17 contains
18 
19  logical function fstr_cutback_active()
20  fstr_cutback_active = is_cutback_active
21  end function
22 
24  subroutine fstr_cutback_init( hecMESH, fstrSOLID, fstrPARAM )
25  type(hecmwst_local_mesh) :: hecMESH
26  type(fstr_param) :: fstrPARAM
27  type(fstr_solid) :: fstrSOLID
28 
29  integer(kind=kint) :: istep, i, j
30  integer(kind=kint) :: ng
31  integer(kind=kint) :: nthick
32  integer(kind=kint) :: ncont, nstate
33 
34  do istep=1,fstrsolid%nstep_tot
35  if( fstrsolid%step_ctrl(istep)%inc_type == stepautoinc ) is_cutback_active = .true.
36  end do
37 
38  if( .not. is_cutback_active ) return
39 
40  !allocate variables to store analysis state
41  !(1)nodal values
42  allocate(fstrsolid%unode_bkup(size(fstrsolid%unode)))
43  allocate(fstrsolid%QFORCE_bkup(size(fstrsolid%QFORCE)))
44  if( associated(fstrsolid%shell_triad) ) allocate(fstrsolid%shell_triad_bkup(size(fstrsolid%shell_triad)))
45  if( associated(fstrsolid%shell_drill) ) allocate(fstrsolid%shell_drill_bkup(size(fstrsolid%shell_drill)))
46  if( fstrsolid%TEMP_ngrp_tot > 0 .or. fstrsolid%TEMP_irres > 0 ) then
47  allocate(fstrsolid%last_temp_bkup(size(fstrsolid%last_temp)))
48  endif
49 
50  !(2)elemental values
51  allocate(fstrsolid%elements_bkup(size(fstrsolid%elements)))
52  do i=1,size(fstrsolid%elements)
53  if (hecmw_is_etype_link( fstrsolid%elements(i)%etype )) cycle
54  if (hecmw_is_etype_patch( fstrsolid%elements(i)%etype )) cycle
55  ng = numofquadpoints( fstrsolid%elements(i)%etype )
56  allocate( fstrsolid%elements_bkup(i)%gausses(ng) )
57  do j=1,ng
58  fstrsolid%elements_bkup(i)%gausses(j)%pMaterial => fstrsolid%elements(i)%gausses(j)%pMaterial
59  call fstr_init_gauss( fstrsolid%elements_bkup(i)%gausses(j) )
60  end do
61  nthick = fstr_shell_num_thickness_points( fstrsolid%elements(i)%etype )
62  if( nthick > 0 ) call fstr_init_shell_layer_gausses( fstrsolid%elements_bkup(i), ng, &
63  fstrsolid%elements(i)%shell_nlayer, nthick )
64  if( associated( fstrsolid%elements(i)%aux ) ) then
65  allocate( fstrsolid%elements_bkup(i)%aux(3,3) )
66  endif
67  end do
68 
69  !(3)contact values
70  ncont = fstrsolid%n_contacts
71  if( associated(fstrsolid%contacts) ) then
72  allocate(fstrsolid%contacts_bkup(ncont))
73  do i=1,ncont
74  nstate = size(fstrsolid%contacts(i)%states)
75  allocate(fstrsolid%contacts_bkup(i)%states(nstate))
76  end do
77  end if
78 
79  !(4)embed values
80  ncont = fstrsolid%n_embeds
81  if( associated(fstrsolid%embeds) ) then
82  allocate(fstrsolid%embeds_bkup(ncont))
83  do i=1,ncont
84  nstate = size(fstrsolid%embeds(i)%states)
85  allocate(fstrsolid%embeds_bkup(i)%states(nstate))
86  end do
87  end if
88 
89  end subroutine
90 
92  subroutine fstr_cutback_finalize( fstrSOLID )
93  type(fstr_solid) :: fstrSOLID
94 
95  integer(kind=kint) :: i, j
96  integer(kind=kint) :: ng
97  integer(kind=kint) :: ncont
98 
99  if( .not. is_cutback_active ) return
100 
101  !deallocate variables to store analysis state
102  !(1)nodal values
103  if( associated(fstrsolid%unode_bkup) ) deallocate(fstrsolid%unode_bkup)
104  if( associated(fstrsolid%QFORCE_bkup) ) deallocate(fstrsolid%QFORCE_bkup)
105  if( associated(fstrsolid%shell_triad_bkup) ) deallocate(fstrsolid%shell_triad_bkup)
106  if( associated(fstrsolid%shell_drill_bkup) ) deallocate(fstrsolid%shell_drill_bkup)
107  if( fstrsolid%TEMP_ngrp_tot > 0 .or. fstrsolid%TEMP_irres > 0 ) then
108  if( associated(fstrsolid%last_temp_bkup) ) deallocate(fstrsolid%last_temp_bkup)
109  endif
110 
111  !(2)elemental values
112  do i=1,size(fstrsolid%elements)
113  if (hecmw_is_etype_link( fstrsolid%elements(i)%etype )) cycle
114  if (hecmw_is_etype_patch( fstrsolid%elements(i)%etype )) cycle
115  ng = numofquadpoints( fstrsolid%elements(i)%etype )
116  do j=1,ng
117  call fstr_finalize_gauss( fstrsolid%elements_bkup(i)%gausses(j) )
118  end do
119  deallocate( fstrsolid%elements_bkup(i)%gausses )
120  call fstr_finalize_shell_layer_gausses( fstrsolid%elements_bkup(i) )
121  if( associated( fstrsolid%elements_bkup(i)%aux ) ) then
122  deallocate( fstrsolid%elements_bkup(i)%aux )
123  endif
124  end do
125  deallocate(fstrsolid%elements_bkup)
126 
127  !(3)contact values
128  ncont = fstrsolid%n_contacts
129  if( associated(fstrsolid%contacts) ) then
130  do i=1,ncont
131  deallocate(fstrsolid%contacts_bkup(i)%states)
132  end do
133  deallocate(fstrsolid%contacts_bkup)
134  end if
135 
136  !(4)embed values
137  ncont = fstrsolid%n_embeds
138  if( associated(fstrsolid%embeds) ) then
139  do i=1,ncont
140  deallocate(fstrsolid%embeds_bkup(i)%states)
141  end do
142  deallocate(fstrsolid%embeds_bkup)
143  end if
144 
145 
146  end subroutine
147 
149  subroutine fstr_cutback_save( fstrSOLID, infoCTChange, infoCTChange_bak )
150  type(fstr_solid), intent(inout) :: fstrSOLID
151  type(fstr_info_contactchange), intent(inout) :: infoCTChange
152  type(fstr_info_contactchange), intent(inout) :: infoCTChange_bak
153 
154  integer(kind=kint) :: i, j
155  integer(kind=kint) :: ng
156  integer(kind=kint) :: ncont, nstate
157 
158  if( .not. is_cutback_active ) return
159 
160  !(1)nodal values
161  do i=1,size(fstrsolid%unode)
162  fstrsolid%unode_bkup(i) = fstrsolid%unode(i)
163  end do
164  do i=1,size(fstrsolid%QFORCE)
165  fstrsolid%QFORCE_bkup(i) = fstrsolid%QFORCE(i)
166  end do
167  if( associated(fstrsolid%shell_triad_bkup) ) then
168  do i=1,size(fstrsolid%shell_triad)
169  fstrsolid%shell_triad_bkup(i) = fstrsolid%shell_triad(i)
170  end do
171  endif
172  if( associated(fstrsolid%shell_drill_bkup) ) then
173  do i=1,size(fstrsolid%shell_drill)
174  fstrsolid%shell_drill_bkup(i) = fstrsolid%shell_drill(i)
175  end do
176  endif
177  if( fstrsolid%TEMP_ngrp_tot > 0 .or. fstrsolid%TEMP_irres > 0 ) then
178  do i=1,size(fstrsolid%last_temp)
179  fstrsolid%last_temp_bkup(i) = fstrsolid%last_temp(i)
180  end do
181  endif
182 
183  !(2)elemental values
184  do i=1,size(fstrsolid%elements)
185  if (hecmw_is_etype_link( fstrsolid%elements(i)%etype )) cycle
186  if (hecmw_is_etype_patch( fstrsolid%elements(i)%etype )) cycle
187  ng = numofquadpoints( fstrsolid%elements(i)%etype )
188  do j=1,ng
189  call fstr_copy_gauss( fstrsolid%elements(i)%gausses(j), fstrsolid%elements_bkup(i)%gausses(j) )
190  end do
191  call fstr_copy_shell_layer_gausses( fstrsolid%elements(i), fstrsolid%elements_bkup(i) )
192  if( associated( fstrsolid%elements(i)%aux ) ) then
193  fstrsolid%elements_bkup(i)%aux(:,:) = fstrsolid%elements(i)%aux(:,:)
194  endif
195  end do
196 
197  !(3)contact values
198  ncont = fstrsolid%n_contacts
199  if( associated(fstrsolid%contacts) ) then
200  do i=1,ncont
201  nstate = size(fstrsolid%contacts(i)%states)
202  do j=1,nstate
203  call contact_state_copy(fstrsolid%contacts(i)%states(j), fstrsolid%contacts_bkup(i)%states(j))
204  enddo
205  end do
206  infoctchange_bak = infoctchange
207  end if
208 
209  !(4)embed values
210  ncont = fstrsolid%n_embeds
211  if( associated(fstrsolid%embeds) ) then
212  do i=1,ncont
213  nstate = size(fstrsolid%embeds(i)%states)
214  do j=1,nstate
215  call contact_state_copy(fstrsolid%embeds(i)%states(j), fstrsolid%embeds_bkup(i)%states(j))
216  enddo
217  end do
218  infoctchange_bak = infoctchange
219  end if
220  end subroutine
221 
223  subroutine fstr_cutback_load( fstrSOLID, infoCTChange, infoCTChange_bak )
224  type(fstr_solid), intent(inout) :: fstrSOLID
225  type(fstr_info_contactchange), intent(inout) :: infoCTChange
226  type(fstr_info_contactchange), intent(inout) :: infoCTChange_bak
227 
228  integer(kind=kint) :: i, j
229  integer(kind=kint) :: ng
230  integer(kind=kint) :: ncont, nstate
231 
232  if( .not. is_cutback_active ) return
233 
234  !(1)nodal values
235  do i=1,size(fstrsolid%unode)
236  fstrsolid%unode(i) = fstrsolid%unode_bkup(i)
237  end do
238  do i=1,size(fstrsolid%QFORCE)
239  fstrsolid%QFORCE(i) = fstrsolid%QFORCE_bkup(i)
240  end do
241  if( associated(fstrsolid%shell_triad_bkup) ) then
242  do i=1,size(fstrsolid%shell_triad)
243  fstrsolid%shell_triad(i) = fstrsolid%shell_triad_bkup(i)
244  fstrsolid%shell_dtriad(i) = fstrsolid%shell_triad_bkup(i)
245  end do
246  endif
247  if( associated(fstrsolid%shell_drill_bkup) ) then
248  do i=1,size(fstrsolid%shell_drill)
249  fstrsolid%shell_drill(i) = fstrsolid%shell_drill_bkup(i)
250  fstrsolid%shell_ddrill(i) = fstrsolid%shell_drill_bkup(i)
251  end do
252  endif
253  if( fstrsolid%TEMP_ngrp_tot > 0 .or. fstrsolid%TEMP_irres > 0 ) then
254  do i=1,size(fstrsolid%last_temp)
255  fstrsolid%last_temp(i) = fstrsolid%last_temp_bkup(i)
256  end do
257  endif
258 
259  !(2)elemental values
260  do i=1,size(fstrsolid%elements)
261  if (hecmw_is_etype_link( fstrsolid%elements(i)%etype )) cycle
262  if (hecmw_is_etype_patch( fstrsolid%elements(i)%etype )) cycle
263  ng = numofquadpoints( fstrsolid%elements(i)%etype )
264  do j=1,ng
265  call fstr_copy_gauss( fstrsolid%elements_bkup(i)%gausses(j), fstrsolid%elements(i)%gausses(j) )
266  end do
267  call fstr_copy_shell_layer_gausses( fstrsolid%elements_bkup(i), fstrsolid%elements(i) )
268  if( associated( fstrsolid%elements(i)%aux ) ) then
269  fstrsolid%elements(i)%aux(:,:) = fstrsolid%elements_bkup(i)%aux(:,:)
270  endif
271  end do
272 
273  !(3)contact values
274  ncont = fstrsolid%n_contacts
275  if( associated(fstrsolid%contacts) ) then
276  do i=1,ncont
277  nstate = size(fstrsolid%contacts(i)%states)
278  do j=1,nstate
279  call contact_state_copy(fstrsolid%contacts_bkup(i)%states(j), fstrsolid%contacts(i)%states(j))
280  enddo
281  end do
282  infoctchange = infoctchange_bak
283  end if
284 
285  !(4)embed values
286  ncont = fstrsolid%n_embeds
287  if( associated(fstrsolid%embeds) ) then
288  do i=1,ncont
289  nstate = size(fstrsolid%embeds(i)%states)
290  do j=1,nstate
291  call contact_state_copy(fstrsolid%embeds_bkup(i)%states(j), fstrsolid%embeds(i)%states(j))
292  enddo
293  end do
294  infoctchange = infoctchange_bak
295  end if
296 
297  end subroutine
298 
299 end module m_fstr_cutback
This module encapsulate the basic functions of all elements provide by this software.
Definition: element.f90:43
integer function numofquadpoints(fetype)
Obtains the number of quadrature points of the element.
Definition: element.f90:450
This module provides functions to deal with cutback.
Definition: fstr_Cutback.f90:7
subroutine fstr_cutback_save(fstrSOLID, infoCTChange, infoCTChange_bak)
Save analysis status.
subroutine fstr_cutback_load(fstrSOLID, infoCTChange, infoCTChange_bak)
Load analysis status.
subroutine fstr_cutback_init(hecMESH, fstrSOLID, fstrPARAM)
Initializer of cutback variables.
subroutine fstr_cutback_finalize(fstrSOLID)
Finalizer of cutback variables.
logical function fstr_cutback_active()
This module defines common data and basic structures for analysis.
Definition: m_fstr.F90:15
This module manages the data structure for contact calculation.
subroutine contact_state_copy(cstate1, cstate2)
Copy.
This modules defines a structure to record history dependent parameter in static analysis.
Definition: mechgauss.f90:6
integer(kind=kint) function fstr_shell_num_thickness_points(etype)
Number of through-thickness quadrature points used by shell stiffness.
Definition: mechgauss.f90:103
subroutine fstr_init_gauss(gauss)
Initializer.
Definition: mechgauss.f90:52
subroutine fstr_init_shell_layer_gausses(element, ng, nlayer, nthick)
Allocate shell history for every surface Gauss point, layer, and thickness point.
Definition: mechgauss.f90:117
subroutine fstr_finalize_gauss(gauss)
Finializer.
Definition: mechgauss.f90:96
subroutine fstr_finalize_shell_layer_gausses(element)
Release shell layer/thickness history.
Definition: mechgauss.f90:255
subroutine fstr_copy_gauss(gauss1, gauss2)
Copy.
Definition: mechgauss.f90:270
subroutine fstr_copy_shell_layer_gausses(element1, element2)
Copy shell layer/thickness history.
Definition: mechgauss.f90:295
FSTR INNER CONTROL PARAMETERS (fstrPARAM)
Definition: m_fstr.F90:156