FrontISTR  5.9.0
Large-scale structural analysis program with finit element method
m_fstr_creatematrix_and_dampingforce Module Reference

This module assembles the tangent stiffness matrix and, in the implicit dynamic case, also the effective dynamic system matrix A = c1*K + c2*M (+ b3*C for connectors) together with the corresponding RHS contribution stored in fstrSOLIDDFORCE. More...

Functions/Subroutines

subroutine, public fstr_creatematrix_and_dampingforce (hecMESH, hecMAT, fstrSOLID, time, tincr, fstrDYNAMIC, coef)
 Assemble the system matrix and, optionally, the dynamic damping force. More...
 
subroutine calc_stiff_and_mass_elem (ic_type, nn, ndof, ecoord, u, u_prev, tt, time, tincr, is_dynamic, fstrSOLID, hecMESH, icel, nodLOCAL, stiff_mat, mass_mat, damp_mat, lumped)
 Compute the element tangent stiffness (and, in the dynamic case, the element mass and connector damping matrices) for the given element type. More...
 
subroutine calc_damping_mat_and_force_elem (ic_type, nn, ndof, material, fstrDYNAMIC, a3, b3, vecA, vecB, stiff_mat, mass_mat, damp_mat, mat, df)
 Combine the element K, M (and connector C) into the effective dynamic system matrix mat = c1*K + c2*M (+ b3*C) and compute the corresponding RHS contribution df. More...
 

Detailed Description

This module assembles the tangent stiffness matrix and, in the implicit dynamic case, also the effective dynamic system matrix A = c1*K + c2*M (+ b3*C for connectors) together with the corresponding RHS contribution stored in fstrSOLIDDFORCE.

Function/Subroutine Documentation

◆ calc_damping_mat_and_force_elem()

subroutine m_fstr_creatematrix_and_dampingforce::calc_damping_mat_and_force_elem ( integer(kind=kint), intent(in)  ic_type,
integer(kind=kint), intent(in)  nn,
integer(kind=kint), intent(in)  ndof,
type(tmaterial), intent(in), pointer  material,
type(fstr_dynamic), intent(in)  fstrDYNAMIC,
real(kind=kreal), intent(in)  a3,
real(kind=kreal), intent(in)  b3,
real(kind=kreal), dimension(:), intent(in)  vecA,
real(kind=kreal), dimension(:), intent(in)  vecB,
real(kind=kreal), dimension(:,:), intent(in)  stiff_mat,
real(kind=kreal), dimension(:,:), intent(in)  mass_mat,
real(kind=kreal), dimension(:,:), intent(in)  damp_mat,
real(kind=kreal), dimension(:,:), intent(out)  mat,
real(kind=kreal), dimension(:), intent(out)  df 
)

Combine the element K, M (and connector C) into the effective dynamic system matrix mat = c1*K + c2*M (+ b3*C) and compute the corresponding RHS contribution df.

Per-material Rayleigh damping (is_elem_Rayleigh_damping) overrides the global fstrDYNAMICray_{m,k} for the current element.

Definition at line 363 of file fstr_CreateMatrix.f90.

Here is the caller graph for this function:

◆ calc_stiff_and_mass_elem()

subroutine m_fstr_creatematrix_and_dampingforce::calc_stiff_and_mass_elem ( integer(kind=kint), intent(in)  ic_type,
integer(kind=kint), intent(inout)  nn,
integer(kind=kint), intent(in)  ndof,
real(kind=kreal), dimension(:,:), intent(in)  ecoord,
real(kind=kreal), dimension(:,:), intent(in)  u,
real(kind=kreal), dimension(:,:), intent(in)  u_prev,
real(kind=kreal), dimension(:), intent(in)  tt,
real(kind=kreal), intent(in)  time,
real(kind=kreal), intent(in)  tincr,
logical, intent(in)  is_dynamic,
type(fstr_solid), intent(inout)  fstrSOLID,
type(hecmwst_local_mesh), intent(in)  hecMESH,
integer(kind=kint), intent(in)  icel,
integer(kind=kint), dimension(:), intent(inout)  nodLOCAL,
real(kind=kreal), dimension(:,:), intent(inout)  stiff_mat,
real(kind=kreal), dimension(:,:), intent(inout)  mass_mat,
real(kind=kreal), dimension(:,:), intent(inout)  damp_mat,
real(kind=kreal), dimension(:), intent(inout)  lumped 
)

Compute the element tangent stiffness (and, in the dynamic case, the element mass and connector damping matrices) for the given element type.

Special case: 341 elements with the selective ES/NS smoothed FEM option (kel341SESNS) leave stiff_mat as zero here because their stiffness is assembled separately through the 881/891 element-type path. The mass matrix is still computed normally so that the dynamic system matrix contains the proper inertia contribution.

Definition at line 177 of file fstr_CreateMatrix.f90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ fstr_creatematrix_and_dampingforce()

subroutine, public m_fstr_creatematrix_and_dampingforce::fstr_creatematrix_and_dampingforce ( type(hecmwst_local_mesh hecMESH,
type(hecmwst_matrix hecMAT,
type(fstr_solid fstrSOLID,
real(kind=kreal), intent(in)  time,
real(kind=kreal), intent(in)  tincr,
type(fstr_dynamic), intent(in), optional  fstrDYNAMIC,
real(kind=kreal), dimension(6), intent(in), optional  coef 
)

Assemble the system matrix and, optionally, the dynamic damping force.

Static / eigen / frequency (fstrDYNAMIC and coef absent): Assemble only the tangent stiffness matrix K into hecMAT.

Implicit dynamic (fstrDYNAMIC and coef present): Assemble the effective dynamic system matrix A = c1*K + c2*M (+ b3*C for connector elements) into hecMAT, and accumulate the mass / Rayleigh damping / connector damping contribution to the RHS into fstrSOLIDDFORCE.

Parameters
hecmeshmesh information
hecmatsystem matrix
fstrsolidsolid state
[in]timecurrent time
[in]tincrtime increment
[in]fstrdynamicdynamic info (omit for static)
[in]coefNewmark time-integration coefficients

Definition at line 30 of file fstr_CreateMatrix.f90.

Here is the call graph for this function:
Here is the caller graph for this function: