FrontISTR  5.7.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) :: ncont, nstate
32 
33  do istep=1,fstrsolid%nstep_tot
34  if( fstrsolid%step_ctrl(istep)%inc_type == stepautoinc ) is_cutback_active = .true.
35  end do
36 
37  if( .not. is_cutback_active ) return
38 
39  !allocate variables to store analysis state
40  !(1)nodal values
41  allocate(fstrsolid%unode_bkup(size(fstrsolid%unode)))
42  allocate(fstrsolid%QFORCE_bkup(size(fstrsolid%QFORCE)))
43  if( fstrsolid%TEMP_ngrp_tot > 0 .or. fstrsolid%TEMP_irres > 0 ) then
44  allocate(fstrsolid%last_temp_bkup(size(fstrsolid%last_temp)))
45  endif
46 
47  !(2)elemental values
48  allocate(fstrsolid%elements_bkup(size(fstrsolid%elements)))
49  do i=1,size(fstrsolid%elements)
50  if (hecmw_is_etype_link( fstrsolid%elements(i)%etype )) cycle
51  if (hecmw_is_etype_patch( fstrsolid%elements(i)%etype )) cycle
52  ng = numofquadpoints( fstrsolid%elements(i)%etype )
53  allocate( fstrsolid%elements_bkup(i)%gausses(ng) )
54  do j=1,ng
55  fstrsolid%elements_bkup(i)%gausses(j)%pMaterial => fstrsolid%elements(i)%gausses(j)%pMaterial
56  call fstr_init_gauss( fstrsolid%elements_bkup(i)%gausses(j) )
57  end do
58  if( associated( fstrsolid%elements(i)%aux ) ) then
59  allocate( fstrsolid%elements_bkup(i)%aux(3,3) )
60  endif
61  end do
62 
63  !(3)contact values
64  ncont = fstrsolid%n_contacts
65  if( associated(fstrsolid%contacts) ) then
66  allocate(fstrsolid%contacts_bkup(ncont))
67  do i=1,ncont
68  nstate = size(fstrsolid%contacts(i)%states)
69  allocate(fstrsolid%contacts_bkup(i)%states(nstate))
70  end do
71  end if
72 
73  !(4)embed values
74  ncont = fstrsolid%n_embeds
75  if( associated(fstrsolid%embeds) ) then
76  allocate(fstrsolid%embeds_bkup(ncont))
77  do i=1,ncont
78  nstate = size(fstrsolid%embeds(i)%states)
79  allocate(fstrsolid%embeds_bkup(i)%states(nstate))
80  end do
81  end if
82 
83  end subroutine
84 
86  subroutine fstr_cutback_finalize( fstrSOLID )
87  type(fstr_solid) :: fstrSOLID
88 
89  integer(kind=kint) :: i, j
90  integer(kind=kint) :: ng
91  integer(kind=kint) :: ncont
92 
93  if( .not. is_cutback_active ) return
94 
95  !deallocate variables to store analysis state
96  !(1)nodal values
97  if( associated(fstrsolid%unode_bkup) ) deallocate(fstrsolid%unode_bkup)
98  if( associated(fstrsolid%QFORCE_bkup) ) deallocate(fstrsolid%QFORCE_bkup)
99  if( fstrsolid%TEMP_ngrp_tot > 0 .or. fstrsolid%TEMP_irres > 0 ) then
100  if( associated(fstrsolid%last_temp_bkup) ) deallocate(fstrsolid%last_temp_bkup)
101  endif
102 
103  !(2)elemental values
104  do i=1,size(fstrsolid%elements)
105  if (hecmw_is_etype_link( fstrsolid%elements(i)%etype )) cycle
106  if (hecmw_is_etype_patch( fstrsolid%elements(i)%etype )) cycle
107  ng = numofquadpoints( fstrsolid%elements(i)%etype )
108  do j=1,ng
109  call fstr_finalize_gauss( fstrsolid%elements_bkup(i)%gausses(j) )
110  end do
111  deallocate( fstrsolid%elements_bkup(i)%gausses )
112  if( associated( fstrsolid%elements_bkup(i)%aux ) ) then
113  deallocate( fstrsolid%elements_bkup(i)%aux )
114  endif
115  end do
116  deallocate(fstrsolid%elements_bkup)
117 
118  !(3)contact values
119  ncont = fstrsolid%n_contacts
120  if( associated(fstrsolid%contacts) ) then
121  do i=1,ncont
122  deallocate(fstrsolid%contacts_bkup(i)%states)
123  end do
124  deallocate(fstrsolid%contacts_bkup)
125  end if
126 
127  !(4)embed values
128  ncont = fstrsolid%n_embeds
129  if( associated(fstrsolid%embeds) ) then
130  do i=1,ncont
131  deallocate(fstrsolid%embeds_bkup(i)%states)
132  end do
133  deallocate(fstrsolid%embeds_bkup)
134  end if
135 
136 
137  end subroutine
138 
140  subroutine fstr_cutback_save( fstrSOLID, infoCTChange, infoCTChange_bak )
141  type(fstr_solid), intent(inout) :: fstrSOLID
142  type(fstr_info_contactchange), intent(inout) :: infoCTChange
143  type(fstr_info_contactchange), intent(inout) :: infoCTChange_bak
144 
145  integer(kind=kint) :: i, j
146  integer(kind=kint) :: ng
147  integer(kind=kint) :: ncont, nstate
148 
149  if( .not. is_cutback_active ) return
150 
151  !(1)nodal values
152  do i=1,size(fstrsolid%unode)
153  fstrsolid%unode_bkup(i) = fstrsolid%unode(i)
154  end do
155  do i=1,size(fstrsolid%QFORCE)
156  fstrsolid%QFORCE_bkup(i) = fstrsolid%QFORCE(i)
157  end do
158  if( fstrsolid%TEMP_ngrp_tot > 0 .or. fstrsolid%TEMP_irres > 0 ) then
159  do i=1,size(fstrsolid%last_temp)
160  fstrsolid%last_temp_bkup(i) = fstrsolid%last_temp(i)
161  end do
162  endif
163 
164  !(2)elemental values
165  do i=1,size(fstrsolid%elements)
166  if (hecmw_is_etype_link( fstrsolid%elements(i)%etype )) cycle
167  if (hecmw_is_etype_patch( fstrsolid%elements(i)%etype )) cycle
168  ng = numofquadpoints( fstrsolid%elements(i)%etype )
169  do j=1,ng
170  call fstr_copy_gauss( fstrsolid%elements(i)%gausses(j), fstrsolid%elements_bkup(i)%gausses(j) )
171  end do
172  if( associated( fstrsolid%elements(i)%aux ) ) then
173  fstrsolid%elements_bkup(i)%aux(:,:) = fstrsolid%elements(i)%aux(:,:)
174  endif
175  end do
176 
177  !(3)contact values
178  ncont = fstrsolid%n_contacts
179  if( associated(fstrsolid%contacts) ) then
180  do i=1,ncont
181  nstate = size(fstrsolid%contacts(i)%states)
182  do j=1,nstate
183  call contact_state_copy(fstrsolid%contacts(i)%states(j), fstrsolid%contacts_bkup(i)%states(j))
184  enddo
185  end do
186  infoctchange_bak = infoctchange
187  end if
188 
189  !(4)embed values
190  ncont = fstrsolid%n_embeds
191  if( associated(fstrsolid%embeds) ) then
192  do i=1,ncont
193  nstate = size(fstrsolid%embeds(i)%states)
194  do j=1,nstate
195  call contact_state_copy(fstrsolid%embeds(i)%states(j), fstrsolid%embeds_bkup(i)%states(j))
196  enddo
197  end do
198  infoctchange_bak = infoctchange
199  end if
200  end subroutine
201 
203  subroutine fstr_cutback_load( fstrSOLID, infoCTChange, infoCTChange_bak )
204  type(fstr_solid), intent(inout) :: fstrSOLID
205  type(fstr_info_contactchange), intent(inout) :: infoCTChange
206  type(fstr_info_contactchange), intent(inout) :: infoCTChange_bak
207 
208  integer(kind=kint) :: i, j
209  integer(kind=kint) :: ng
210  integer(kind=kint) :: ncont, nstate
211 
212  if( .not. is_cutback_active ) return
213 
214  !(1)nodal values
215  do i=1,size(fstrsolid%unode)
216  fstrsolid%unode(i) = fstrsolid%unode_bkup(i)
217  end do
218  do i=1,size(fstrsolid%QFORCE)
219  fstrsolid%QFORCE(i) = fstrsolid%QFORCE_bkup(i)
220  end do
221  if( fstrsolid%TEMP_ngrp_tot > 0 .or. fstrsolid%TEMP_irres > 0 ) then
222  do i=1,size(fstrsolid%last_temp)
223  fstrsolid%last_temp(i) = fstrsolid%last_temp_bkup(i)
224  end do
225  endif
226 
227  !(2)elemental values
228  do i=1,size(fstrsolid%elements)
229  if (hecmw_is_etype_link( fstrsolid%elements(i)%etype )) cycle
230  if (hecmw_is_etype_patch( fstrsolid%elements(i)%etype )) cycle
231  ng = numofquadpoints( fstrsolid%elements(i)%etype )
232  do j=1,ng
233  call fstr_copy_gauss( fstrsolid%elements_bkup(i)%gausses(j), fstrsolid%elements(i)%gausses(j) )
234  end do
235  if( associated( fstrsolid%elements(i)%aux ) ) then
236  fstrsolid%elements(i)%aux(:,:) = fstrsolid%elements_bkup(i)%aux(:,:)
237  endif
238  end do
239 
240  !(3)contact values
241  ncont = fstrsolid%n_contacts
242  if( associated(fstrsolid%contacts) ) then
243  do i=1,ncont
244  nstate = size(fstrsolid%contacts(i)%states)
245  do j=1,nstate
246  call contact_state_copy(fstrsolid%contacts_bkup(i)%states(j), fstrsolid%contacts(i)%states(j))
247  enddo
248  end do
249  infoctchange = infoctchange_bak
250  end if
251 
252  !(4)embed values
253  ncont = fstrsolid%n_embeds
254  if( associated(fstrsolid%embeds) ) then
255  do i=1,ncont
256  nstate = size(fstrsolid%embeds(i)%states)
257  do j=1,nstate
258  call contact_state_copy(fstrsolid%embeds_bkup(i)%states(j), fstrsolid%embeds(i)%states(j))
259  enddo
260  end do
261  infoctchange = infoctchange_bak
262  end if
263 
264  end subroutine
265 
266 end module m_fstr_cutback
m_fstr_cutback::fstr_cutback_init
subroutine fstr_cutback_init(hecMESH, fstrSOLID, fstrPARAM)
Initializer of cutback variables.
Definition: fstr_Cutback.f90:25
elementinfo::numofquadpoints
integer function numofquadpoints(fetype)
Obtains the number of quadrature points of the element.
Definition: element.f90:450
m_fstr_cutback::fstr_cutback_active
logical function fstr_cutback_active()
Definition: fstr_Cutback.f90:20
m_fstr_cutback::fstr_cutback_load
subroutine fstr_cutback_load(fstrSOLID, infoCTChange, infoCTChange_bak)
Load analysis status.
Definition: fstr_Cutback.f90:204
mcontactdef
This module manages the data structure for contact calculation.
Definition: fstr_contact_def.F90:12
mmechgauss::fstr_finalize_gauss
subroutine fstr_finalize_gauss(gauss)
Finializer.
Definition: mechgauss.f90:88
mmechgauss::fstr_copy_gauss
subroutine fstr_copy_gauss(gauss1, gauss2)
Copy.
Definition: mechgauss.f90:95
m_fstr::fstr_solid
Definition: m_fstr.f90:238
m_fstr::fstr_param
FSTR INNER CONTROL PARAMETERS (fstrPARAM)
Definition: m_fstr.f90:154
mmechgauss
This modules defines a structure to record history dependent parameter in static analysis.
Definition: mechgauss.f90:6
mcontactdef::fstr_info_contactchange
Definition: fstr_contact_def.F90:59
m_fstr
This module defines common data and basic structures for analysis.
Definition: m_fstr.f90:15
m_fstr_cutback::fstr_cutback_save
subroutine fstr_cutback_save(fstrSOLID, infoCTChange, infoCTChange_bak)
Save analysis status.
Definition: fstr_Cutback.f90:141
m_fstr_cutback
This module provides functions to deal with cutback.
Definition: fstr_Cutback.f90:7
elementinfo
This module encapsulate the basic functions of all elements provide by this software.
Definition: element.f90:43
mmechgauss::fstr_init_gauss
subroutine fstr_init_gauss(gauss)
Initializer.
Definition: mechgauss.f90:44
m_fstr_cutback::fstr_cutback_finalize
subroutine fstr_cutback_finalize(fstrSOLID)
Finalizer of cutback variables.
Definition: fstr_Cutback.f90:87