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