FrontISTR  5.7.1
Large-scale structural analysis program with finit element method
fstr_Residual.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 hecmw
10  implicit none
11 
12  public :: fstr_update_ndforce
13  public :: fstr_update_ndforce_spc
14  public :: fstr_get_residual
15  public :: fstr_get_norm_contact
17  public :: fstr_get_x_norm_contact
18  public :: fstr_get_potential
19 
20  private :: fstr_update_ndforce_mpc
21 
22 contains
23 
24  !C---------------------------------------------------------------------*
25  subroutine fstr_update_ndforce(cstep,hecMESH,hecMAT,fstrSOLID,conMAT)
26  !C---------------------------------------------------------------------*
32  use m_fstr
33  use muload
34  use m_fstr_spring
35  integer(kind=kint), intent(in) :: cstep
36  type(hecmwst_local_mesh), intent(in) :: hecmesh
37  type(hecmwst_matrix), intent(inout) :: hecmat
38  type(fstr_solid), intent(inout) :: fstrsolid
39  type(hecmwst_matrix), intent(inout), optional :: conmat
40  ! Local variables
41  integer(kind=kint) :: ndof, idof
42  real(kind=kreal) :: factor
43 
44  factor = fstrsolid%factor(2)
45 
46  ! Set residual load
47  do idof=1, hecmesh%n_node* hecmesh%n_dof
48  hecmat%B(idof)=fstrsolid%GL(idof)-fstrsolid%QFORCE(idof)
49  end do
50  ndof = hecmat%NDOF
51 
52  call fstr_update_ndforce_spring( cstep, hecmesh, fstrsolid, hecmat%B )
53 
54  ! Consider Uload
55  call uresidual( cstep, factor, hecmat%B )
56 
57  ! Consider EQUATION condition
58  call fstr_update_ndforce_mpc( hecmesh, hecmat%B )
59 
60  ! Consider SPC condition
61  call fstr_update_ndforce_spc( cstep, hecmesh, fstrsolid, hecmat%B )
62  if(present(conmat)) call fstr_update_ndforce_spc( cstep, hecmesh, fstrsolid, conmat%B )
63 
64  !
65  call hecmw_update_r(hecmesh,hecmat%B,hecmesh%n_node, ndof)
66  end subroutine fstr_update_ndforce
67 
68  subroutine fstr_update_ndforce_mpc( hecMESH, B )
69  use m_fstr
70  type(hecmwst_local_mesh), intent(in) :: hecmesh
71  real(kind=kreal), intent(inout) :: b(:)
72  ! Local variables
73  integer(kind=kint) ndof, ig0, is0, ie0, ik, in, idof
74  real(kind=kreal) :: rhs, lambda
75 
76  ndof = hecmesh%n_dof
77  outer: do ig0=1,hecmesh%mpc%n_mpc
78  is0= hecmesh%mpc%mpc_index(ig0-1)+1
79  ie0= hecmesh%mpc%mpc_index(ig0)
80  do ik= is0, ie0
81  if (hecmesh%mpc%mpc_dof(ik) > ndof) cycle outer
82  enddo
83  ! Suppose the lagrange multiplier= first dof of first node
84  in = hecmesh%mpc%mpc_item(is0)
85  idof = hecmesh%mpc%mpc_dof(is0)
86  rhs = hecmesh%mpc%mpc_val(is0)
87  lambda = b(ndof*(in-1)+idof)/rhs
88  ! update nodal residual
89  do ik= is0, ie0
90  in = hecmesh%mpc%mpc_item(ik)
91  idof = hecmesh%mpc%mpc_dof(ik)
92  rhs = hecmesh%mpc%mpc_val(ik)
93  b(ndof*(in-1)+idof) = b(ndof*(in-1)+idof) - rhs*lambda
94  enddo
95  enddo outer
96  end subroutine fstr_update_ndforce_mpc
97 
98  subroutine fstr_update_ndforce_spc( cstep, hecMESH, fstrSOLID, B )
99  use m_fstr
100  integer(kind=kint), intent(in) :: cstep
101  type(hecmwst_local_mesh), intent(in) :: hecmesh
102  type(fstr_solid), intent(in) :: fstrsolid
103  real(kind=kreal), intent(inout) :: b(:)
104  ! Local variables
105  integer(kind=kint) ndof, ig0, ig, ityp, is0, ie0, ik, in, idof1, idof2, idof
106  integer(kind=kint) :: grpid
107  real(kind=kreal) :: rhs
108 
109  ndof = hecmesh%n_dof
110  fstrsolid%REACTION = 0.d0
111 
112  do ig0= 1, fstrsolid%BOUNDARY_ngrp_tot
113  grpid = fstrsolid%BOUNDARY_ngrp_GRPID(ig0)
114  if( .not. fstr_isboundaryactive( fstrsolid, grpid, cstep ) ) cycle
115  ig= fstrsolid%BOUNDARY_ngrp_ID(ig0)
116  rhs= fstrsolid%BOUNDARY_ngrp_val(ig0)
117  ityp= fstrsolid%BOUNDARY_ngrp_type(ig0)
118  is0= hecmesh%node_group%grp_index(ig-1) + 1
119  ie0= hecmesh%node_group%grp_index(ig )
120  do ik= is0, ie0
121  in = hecmesh%node_group%grp_item(ik)
122  idof1 = ityp/10
123  idof2 = ityp - idof1*10
124  if( fstrsolid%BOUNDARY_ngrp_rotID(ig0) > 0 ) then
125  idof1 = 1
126  idof2 = ndof
127  end if
128  do idof=idof1,idof2
129  b( ndof*(in-1) + idof ) = 0.d0
130  !for output reaction force
131  fstrsolid%REACTION(ndof*(in-1)+idof) = fstrsolid%QFORCE(ndof*(in-1)+idof)
132  !count embed force as reaction force
133  if( associated(fstrsolid%EMBED_NFORCE) ) fstrsolid%REACTION(ndof*(in-1)+idof) = &
134  & fstrsolid%QFORCE(ndof*(in-1)+idof) - fstrsolid%EMBED_NFORCE(ndof*(in-1)+idof)
135  enddo
136  enddo
137  enddo
138  end subroutine fstr_update_ndforce_spc
139 
140 
142  subroutine fstr_calc_residual_vector(hecMESH, hecMAT, fstrSOLID, ctime, tincr, iter, cstep, dtime, fstrPARAM)
143  use m_fstr_update
144  use m_fstr_ass_load
145  implicit none
146  integer, intent(in) :: cstep
147  type (hecmwST_local_mesh) :: hecMESH
148  type (hecmwST_matrix) :: hecMAT
149  type (fstr_solid) :: fstrSOLID
150  real(kind=kreal), intent(in) :: ctime
151  real(kind=kreal), intent(in) :: dtime
152  type (fstr_param) :: fstrPARAM
153  real(kind=kreal), intent(in) :: tincr
154  integer(kind=kint) :: iter
155 
156  ! ----- update the strain, stress, and internal force
157  call fstr_updatenewton(hecmesh, hecmat, fstrsolid, ctime, tincr, iter)
158 
159  ! ----- Set residual
160  if( fstrsolid%DLOAD_follow /= 0 .or. fstrsolid%CLOAD_ngrp_rot /= 0 ) &
161  & call fstr_ass_load(cstep, ctime+dtime, hecmesh, hecmat, fstrsolid, fstrparam )
162 
163 
164  if( fstrsolid%elemact%ELEMACT_egrp_tot > 0 ) then
165  call fstr_update_elemact_solid_by_value( hecmesh, fstrsolid, cstep, ctime )
166  endif
167 
168  call fstr_update_ndforce(cstep, hecmesh, hecmat, fstrsolid)
169 
170  end subroutine fstr_calc_residual_vector
171 
173  subroutine fstr_calc_residual_vector_with_x(hecMESH, hecMAT, fstrSOLID, ctime, tincr, iter, cstep, dtime, fstrPARAM)
174  use m_fstr_update
175  implicit none
176  integer, intent(in) :: cstep
177  type (hecmwST_local_mesh) :: hecMESH
178  type (hecmwST_matrix) :: hecMAT
179  type (fstr_solid) :: fstrSOLID
180  real(kind=kreal), intent(in) :: ctime
181  real(kind=kreal), intent(in) :: dtime
182  type (fstr_param) :: fstrPARAM
183  real(kind=kreal), intent(in) :: tincr
184  integer(kind=kint) :: iter
185 
186  integer :: i
187  real(kind=kreal) :: dunode_bak(hecmat%ndof*hecmesh%n_node)
188 
189  do i=1, hecmat%ndof*hecmesh%n_node
190  dunode_bak(i) = fstrsolid%dunode(i)
191  fstrsolid%dunode(i) = fstrsolid%dunode(i) + hecmat%X(i)
192  enddo
193  call fstr_calc_residual_vector(hecmesh, hecmat, fstrsolid, ctime, tincr, iter, cstep, dtime, fstrparam)
194  do i=1, hecmat%ndof*hecmesh%n_node
195  fstrsolid%dunode(i) = dunode_bak(i)
196  enddo
197  end subroutine fstr_calc_residual_vector_with_x
198 
200  real(kind=kreal) function fstr_get_residual( force, hecMESH )
201  use m_fstr
202  real(kind=kreal), intent(in) :: force(:)
203  type(hecmwst_local_mesh), intent(in) :: hecmesh
204  integer :: ndof
205  ndof = hecmesh%n_dof
206  call hecmw_innerproduct_r(hecmesh,ndof,force,force,fstr_get_residual)
207  end function
208 
210  real(kind=kreal) function fstr_get_energy( force, displacement, hecMESH )
211  use m_fstr
212  real(kind=kreal), intent(in) :: force(:), displacement(:)
213  type(hecmwst_local_mesh), intent(in) :: hecmesh
214  integer :: ndof
215  ndof = hecmesh%n_dof
216  call hecmw_innerproduct_r(hecmesh, ndof, force, displacement, fstr_get_energy)
217  end function
218 
220  real(kind=kreal) function fstr_get_norm_contact(flag,hecMESH,hecMAT,fstrSOLID,hecLagMAT)
221  use m_fstr
222  type(hecmwst_local_mesh), intent(in) :: hecmesh
223  type(hecmwst_matrix), intent(in) :: hecmat
224  type(fstr_solid), intent(in) :: fstrsolid
225  type(hecmwst_matrix_lagrange), intent(in) :: heclagmat
226  character(len=13) :: flag
227  real(kind=kreal) :: tmp1, tmp2, bi
228  integer :: i, i0, ndof
229  if( flag=='residualForce' )then
230  ndof = hecmesh%n_dof
231  call hecmw_innerproduct_r(hecmesh,ndof,hecmat%B,hecmat%B,tmp1)
232  tmp2 = 0.0d0
233  i0 = hecmesh%n_node*ndof
234  do i=1,heclagmat%num_lagrange
235  bi = hecmat%B(i0+i)
236  tmp2 = tmp2 + bi*bi
237  enddo
238  call hecmw_allreduce_r1(hecmesh,tmp2,hecmw_sum)
239  fstr_get_norm_contact = tmp1 + tmp2
240  elseif( flag==' force' )then
241  call hecmw_innerproduct_r(hecmesh,ndof,fstrsolid%QFORCE,fstrsolid%QFORCE,fstr_get_norm_contact)
242  endif
243  end function
244 
245  !
246  function fstr_get_norm_para_contact(hecMAT,hecLagMAT,conMAT,hecMESH) result(rhsB)
247  use m_fstr
248  implicit none
249  type(hecmwst_matrix), intent(in) :: hecmat
250  type(hecmwst_matrix_lagrange), intent(in) :: heclagmat
251  type(hecmwst_matrix), intent(in) :: conmat
252  type(hecmwst_local_mesh), intent(in) :: hecmesh
253  !
254  real(kind=kreal) :: rhsb
255  integer(kind=kint) :: i,ndof,nndof,npndof,num_lagrange
256  real(kind=kreal), allocatable :: rhs_con(:)
257  real(kind=kreal), pointer :: rhs_lag(:)
258 
259  ndof = conmat%ndof
260  nndof = conmat%N * ndof
261  npndof = conmat%NP * ndof
262  num_lagrange = heclagmat%num_lagrange
263 
264  allocate(rhs_con(npndof))
265  do i=1,npndof
266  rhs_con(i) = conmat%B(i)
267  enddo
268  call hecmw_assemble_r(hecmesh, rhs_con, conmat%NP, conmat%NDOF)
269 
270  do i=1,nndof
271  rhs_con(i) = rhs_con(i) + hecmat%B(i)
272  enddo
273 
274  rhs_lag => conmat%B(npndof+1:npndof+num_lagrange)
275 
276  rhsb = dot_product(rhs_con(1:nndof), rhs_con(1:nndof)) + dot_product(rhs_lag(:), rhs_lag(:))
277  call hecmw_allreduce_r1(hecmesh, rhsb, hecmw_sum)
278  deallocate(rhs_con)
279 
280  end function fstr_get_norm_para_contact
281 
282  function fstr_get_x_norm_contact(hecMAT,hecLagMAT,hecMESH) result(rhsX)
283  use m_fstr
284  implicit none
285  type(hecmwst_matrix), intent(in) :: hecmat
286  type(hecmwst_matrix_lagrange), intent(in) :: heclagmat
287  type(hecmwst_local_mesh), intent(in) :: hecmesh
288  real(kind=kreal) :: rhsx
289  integer(kind=kint) :: nndof, npndof, i
290 
291  nndof = hecmat%N * hecmat%NDOF
292  npndof = hecmat%NP * hecmat%NDOF
293  rhsx = 0.d0
294  do i=1,nndof
295  rhsx = rhsx + hecmat%X(i) ** 2
296  end do
297  do i=1,heclagmat%num_lagrange
298  rhsx = rhsx + hecmat%X(npndof+i) ** 2
299  end do
300  call hecmw_allreduce_r1(hecmesh, rhsx, hecmw_sum)
301 
302  end function fstr_get_x_norm_contact
303 
304  !C---------------------------------------------------------------------*
305  function fstr_get_potential(cstep,hecMESH,hecMAT,fstrSOLID,ptype) result(potential)
306  use m_fstr
307  use muload
308  use m_fstr_spring
309  integer(kind=kint), intent(in) :: cstep
310  type(hecmwst_local_mesh), intent(in) :: hecmesh
311  type(hecmwst_matrix), intent(inout) :: hecmat
312  type(fstr_solid), intent(inout) :: fstrsolid
313  integer(kind=kint), intent(in) :: ptype
314  real(kind=kreal) :: potential
315  ! Local variables
316  integer(kind=kint) :: ndof, i, icel, icel0, ig
317  real(kind=kreal), allocatable :: totdisp(:), totload(:)
318  real(kind=kreal) :: factor
319  real(kind=kreal) :: extenal_work, internal_work
320 
321  ndof = hecmat%NDOF
322  factor = fstrsolid%factor(2)
323  potential = 0.d0
324 
325  allocate(totload(ndof*hecmesh%n_node))
326 
327  ! Set external load
328  do i=1,ndof*hecmesh%n_node
329  totload(i) = 0.5d0*(fstrsolid%GL(i)+fstrsolid%GL0(i))
330  end do
331  ! Consiter Spring
332  call fstr_update_ndforce_spring( cstep, hecmesh, fstrsolid, totload )
333  ! Consider Uload
334  call uresidual( cstep, factor, totload )
335  ! Consider EQUATION condition
336  call fstr_update_ndforce_mpc( hecmesh, totload )
337 
338  if( ptype == 1 ) then
339  ! calc internal work
340  internal_work = 0.d0
341  do icel0=1,hecmesh%ne_internal
342  icel=hecmesh%elem_internal_list(icel0)
343  do ig=1,size(fstrsolid%elements(icel)%gausses)
344  internal_work = internal_work + fstrsolid%elements(icel)%gausses(ig)%strain_energy
345  enddo
346  enddo
347  call hecmw_allreduce_r1(hecmesh, internal_work, hecmw_sum)
348  ! calc external work
349  call hecmw_innerproduct_r(hecmesh,ndof,totload,fstrsolid%dunode,extenal_work)
350  potential = internal_work - extenal_work
351  else if( ptype == 2 ) then ! multiply internal force with displacement
352  ! calc internal work
353  call hecmw_innerproduct_r(hecmesh,ndof,0.5d0*(fstrsolid%QFORCE+fstrsolid%QFORCE_bak),fstrsolid%dunode,internal_work)
354  ! calc external work
355  call hecmw_innerproduct_r(hecmesh,ndof,totload,fstrsolid%dunode,extenal_work)
356  potential = internal_work - extenal_work
357  else if( ptype == 3 ) then ! norm of residual
358  ! Consider internal force
359  totload(:) = fstrsolid%QFORCE(:) - totload(:)
360  ! Consider SPC condition
361  call fstr_update_ndforce_spc( cstep, hecmesh, fstrsolid, totload )
362  ! calc norm of residual
363  call hecmw_innerproduct_r(hecmesh,ndof,totload,totload,potential)
364  potential = dsqrt(potential)
365  else
366  stop "ptype error in fstr_get_potential"
367  endif
368  !write(IMSG,*) "internal_work,extenal_work",internal_work,extenal_work
369  end function fstr_get_potential
370 
371  function fstr_get_potential_with_x(cstep,hecMESH,hecMAT,fstrSOLID,ptype) result(potential)
372  use m_fstr
373  use muload
374  use m_fstr_spring
375  implicit none
376  integer(kind=kint), intent(in) :: cstep
377  type(hecmwst_local_mesh), intent(in) :: hecmesh
378  type(hecmwst_matrix), intent(inout) :: hecmat
379  type(fstr_solid), intent(inout) :: fstrsolid
380  integer(kind=kint), intent(in) :: ptype
381  real(kind=kreal) :: potential
382 
383  integer :: i
384  real(kind=kreal) :: dunode_bak(hecmat%ndof*hecmesh%n_node)
385 
386  do i=1, hecmat%ndof*hecmesh%n_node
387  dunode_bak(i) = fstrsolid%dunode(i)
388  fstrsolid%dunode(i) = fstrsolid%dunode(i) + hecmat%X(i)
389  enddo
390  potential = fstr_get_potential(cstep,hecmesh,hecmat,fstrsolid,ptype)
391  do i=1, hecmat%ndof*hecmesh%n_node
392  fstrsolid%dunode(i) = dunode_bak(i)
393  enddo
394  end function fstr_get_potential_with_x
395 
396  subroutine plot_potential_graph( cstep, hecMESH, hecMAT, fstrSOLID, fstrPARAM, &
397  restrt_step_num, sub_step, ctime, dtime, iter, tincr )
398  use m_fstr
399  integer, intent(in) :: cstep
400  type (hecmwST_local_mesh) :: hecMESH
401  type (hecmwST_matrix) :: hecMAT
402  type (fstr_solid) :: fstrSOLID
403  integer, intent(in) :: sub_step
404  real(kind=kreal), intent(in) :: ctime
405  real(kind=kreal), intent(in) :: dtime
406  type (fstr_param) :: fstrPARAM
407  type (hecmwST_matrix_lagrange) :: hecLagMAT
408  integer(kind=kint) :: restrt_step_num
409 
410  integer(kind=kint) :: iter, i, j
411  integer(kind=kint), parameter :: ntot = 30
412  real(kind=kreal) :: tincr, alpha, dulen
413  real(kind=kreal) :: pot(3)
414  real(kind=kreal), allocatable :: dunode_bak(:)
415 
416  dulen = dsqrt(dot_product(fstrsolid%dunode(:),fstrsolid%dunode(:)))
417  write(imsg,*) "dulen",dulen
418 
419  allocate(dunode_bak(size(fstrsolid%dunode)))
420  do i=1, hecmat%ndof*hecmesh%n_node
421  dunode_bak(i) = fstrsolid%dunode(i)
422  enddo
423 
424  do i=-ntot,ntot
425  alpha = 5.d-3*dble(i)/dble(ntot)
426  hecmat%X(:) = alpha*fstrsolid%dunode(:)
427  do j=1, size(fstrsolid%dunode)
428  fstrsolid%dunode(j) = dunode_bak(j)+hecmat%X(j)
429  enddo
430  call fstr_calc_residual_vector(hecmesh, hecmat, fstrsolid, ctime, tincr, iter, cstep, dtime, fstrparam)
431  pot(1) = fstr_get_potential(cstep,hecmesh,hecmat,fstrsolid,1)
432  pot(2) = fstr_get_potential(cstep,hecmesh,hecmat,fstrsolid,2)
433  pot(3) = fstr_get_potential(cstep,hecmesh,hecmat,fstrsolid,3)
434  write(imsg,'(I4,5(",",1pE16.9))') i,alpha,alpha*dulen,pot(1:3)
435  enddo
436 
437  do i=1, size(fstrsolid%dunode)
438  fstrsolid%dunode(i) = dunode_bak(i)
439  enddo
440  hecmat%X(:) = 0.d0
441  call fstr_calc_residual_vector_with_x(hecmesh, hecmat, fstrsolid, ctime, tincr, iter, cstep, dtime, fstrparam)
442  end subroutine
443 
444 end module m_fstr_residual
Definition: hecmw.f90:6
This module provides functions to take into account external load.
subroutine fstr_ass_load(cstep, ctime, hecMESH, hecMAT, fstrSOLID, fstrPARAM)
This subroutine assmble following external force into fstrSOLIDGL and hecMATB afterwards.
This module provide a function to elemact elements.
subroutine fstr_update_elemact_solid_by_value(hecMESH, fstrSOLID, cstep, ctime)
This module provides function to calculate residual of nodal force.
subroutine, public fstr_update_ndforce_spc(cstep, hecMESH, fstrSOLID, B)
subroutine, public fstr_update_ndforce(cstep, hecMESH, hecMAT, fstrSOLID, conMAT)
subroutine plot_potential_graph(cstep, hecMESH, hecMAT, fstrSOLID, fstrPARAM, restrt_step_num, sub_step, ctime, dtime, iter, tincr)
real(kind=kreal) function, public fstr_get_norm_para_contact(hecMAT, hecLagMAT, conMAT, hecMESH)
real(kind=kreal) function, public fstr_get_residual(force, hecMESH)
Calculate magnitude of a real vector.
real(kind=kreal) function, public fstr_get_potential(cstep, hecMESH, hecMAT, fstrSOLID, ptype)
subroutine fstr_calc_residual_vector(hecMESH, hecMAT, fstrSOLID, ctime, tincr, iter, cstep, dtime, fstrPARAM)
\breaf This subroutine calculate residual vector
subroutine fstr_calc_residual_vector_with_x(hecMESH, hecMAT, fstrSOLID, ctime, tincr, iter, cstep, dtime, fstrPARAM)
\breaf This subroutine calculate residual vector
real(kind=kreal) function fstr_get_potential_with_x(cstep, hecMESH, hecMAT, fstrSOLID, ptype)
real(kind=kreal) function, public fstr_get_norm_contact(flag, hecMESH, hecMAT, fstrSOLID, hecLagMAT)
Calculate square norm.
real(kind=kreal) function fstr_get_energy(force, displacement, hecMESH)
Calculate magnitude of a real vector.
real(kind=kreal) function, public fstr_get_x_norm_contact(hecMAT, hecLagMAT, hecMESH)
This module provides functions to deal with spring force.
Definition: fstr_Spring.f90:7
subroutine fstr_update_ndforce_spring(cstep, hecMESH, fstrSOLID, B)
Definition: fstr_Spring.f90:46
This module provides function to calculate to do updates.
Definition: fstr_Update.f90:6
subroutine fstr_updatenewton(hecMESH, hecMAT, fstrSOLID, time, tincr, iter, strainEnergy)
Update displacement, stress, strain and internal forces.
Definition: fstr_Update.f90:26
This module defines common data and basic structures for analysis.
Definition: m_fstr.f90:15
integer(kind=kint), parameter imsg
Definition: m_fstr.f90:111
logical function fstr_isboundaryactive(fstrSOLID, nbc, cstep)
Definition: m_fstr.f90:1036
This subroutine read in used-defined loading tangent.
Definition: uload.f90:7
subroutine uresidual(cstep, factor, residual)
This subroutine take consider of user-defined external loading.
Definition: uload.f90:39