FEDEM Solver  R8.0
Source code of the dynamics solver
Functions/Subroutines
supelroutinesmodule Module Reference

Module with subroutines for superelement calculations. More...

Functions/Subroutines

subroutine setsupelsvelacc (sups, sam, velGlobal, accGlobal)
 Extracts local velocities and accelerations for the superelements. More...
 
subroutine buildfinitvelacc (sup, beta, gamma, h)
 Calculates deformational velocities and accelerations. More...
 
subroutine incsupelsgendofs (sups, solinc, useTotalInc)
 Increments the generalized DOFs for all superelements. More...
 
subroutine updatesupels (sups, supLoads, env, beta, gamma, time, timeStep, istep, iter, newPositions, ierr)
 Updates all superelements in the model based on the computed state. More...
 
subroutine updatesupelsstatic (sups, supLoads, env, time, iter, linInc, ierr)
 Updates all superelements in the model based on the computed state. More...
 
subroutine, private updatesupelload (Q, S, supLoad, ierr)
 Adds the superelement loads to the system external load vector. More...
 
subroutine updatesupeldamping (sups, engs, alpha, ierr)
 Updates the superelement damping matrices. More...
 
subroutine updateseaenvironment (env, triads, sups, time, iter, ierr)
 Updates the current sea state. More...
 
subroutine, private calcmorisonforces (sup, env, time, istep, iter, ierr)
 Calculates Morison force contributions for a two-noded beam element. More...
 
subroutine, private calcbuoyancyforces (sup, env, time, iter, ierr)
 Calculates buoyancy forces and associated load correction stiffness. More...
 
subroutine addinsupforces (sam, sups, FSk, FDk, FIk, Qk, RFk, ierr)
 Adds superelement forces into corresponding system force vectors. More...
 
subroutine addinstaticsupforces (sam, sups, FSk, Qk, RFk, ierr)
 Adds superelement forces into corresponding system force vectors. More...
 
subroutine addinsupmat (supMat, sysMat, sup, sam, err, sysRhs, scale)
 Adds a superelement matrix into the equivalent system matrix. More...
 
subroutine comptanstiff (sup, ierr)
 Computes the tangential superelement stiffness matrix. More...
 
subroutine buildsupnewtonmat (newTangent, scaleM, scaleC, scaleK, sup, ierr)
 Computes the superelement Newton matrix. More...
 
subroutine, private scaledmatmul (m, n, alpha, A, X, Y, ldA)
 Calculates the scaled matrix-vector product Y = α*A*X. More...
 

Detailed Description

Module with subroutines for superelement calculations.

This module contains a set of subroutines for performing various computation tasks on the supeltypemodule::supeltype objects in the model during the dynamic or quasi-static simulation.

Function/Subroutine Documentation

◆ addinstaticsupforces()

subroutine supelroutinesmodule::addinstaticsupforces ( type(samtype), intent(in)  sam,
type(supeltype), dimension(:), intent(in)  sups,
real(dp), dimension(:), intent(inout)  FSk,
real(dp), dimension(:), intent(inout)  Qk,
real(dp), dimension(:), intent(inout)  RFk,
integer, intent(out)  ierr 
)

Adds superelement forces into corresponding system force vectors.

Parameters
[in]samData for managing system matrix assembly
[in]supsAll superelements in the model
FSkInternal stiffness force vector
QkExternal force vector (gravitation loads)
RFkReaction forces
[out]ierrError flag

This subroutine only considers the static forces, the stiffness and external forces.

Author
Knut Morten Okstad
Date
5 June 2002
Here is the call graph for this function:
Here is the caller graph for this function:

◆ addinsupforces()

subroutine supelroutinesmodule::addinsupforces ( type(samtype), intent(in)  sam,
type(supeltype), dimension(:), intent(in)  sups,
real(dp), dimension(:), intent(inout)  FSk,
real(dp), dimension(:), intent(inout)  FDk,
real(dp), dimension(:), intent(inout)  FIk,
real(dp), dimension(:), intent(inout)  Qk,
real(dp), dimension(:), intent(inout)  RFk,
integer, intent(out)  ierr 
)

Adds superelement forces into corresponding system force vectors.

Parameters
[in]samData for managing system matrix assembly
[in]supsAll superelements in the model
FSkInternal stiffness force vector
FDkInternal damping force vector
FIkInternal internal force vector
QkExternal force vector (gravitation loads)
RFkReaction forces
[out]ierrError flag
Author
Bjorn Haugen
Date
December 2001
Here is the call graph for this function:
Here is the caller graph for this function:

◆ addinsupmat()

subroutine supelroutinesmodule::addinsupmat ( real(dp), dimension(:,:), intent(in)  supMat,
type(sysmatrixtype), intent(inout)  sysMat,
type(supeltype), intent(in)  sup,
type(samtype), intent(in)  sam,
integer, intent(out)  err,
real(dp), dimension(:), intent(inout), optional  sysRhs,
real(dp), intent(in), optional  scale 
)

Adds a superelement matrix into the equivalent system matrix.

Parameters
[in]supMatThe superelement matrix to add into the system
sysMatSystem coefficient matrix (stiffness, mass, ...)
[in]supThe superelement the provided matrix is associated with
[in]samData for managing system matrix assembly
[out]errError flag
sysRhsSystem right-hand-side vector
[in]scaleOptional scaling factor for the matrix to add

This subroutine transforms the superelement matrix supMat to the system directions for each DOF associated with it, and then adds it into the system matrix sysMat. If any of the superelement DOFs are prescribed, the associated force contributions are added into the system force vector sysRhs, if present.

Author
Karl Erik Thoresen
Date
Januar 1999
Here is the call graph for this function:
Here is the caller graph for this function:

◆ buildfinitvelacc()

subroutine supelroutinesmodule::buildfinitvelacc ( type(supeltype), intent(inout)  sup,
real(dp), intent(in)  beta,
real(dp), intent(in)  gamma,
real(dp), intent(in)  h 
)

Calculates deformational velocities and accelerations.

Parameters
supThe superelement to calculate velocity/acceleration for
[in]betaNewmark time integration parameter
[in]gammaNewmark time integration parameter
[in]hCurrent time increment size

This subroutine integrates the deformational velocities and accelerations based on the Newmark integration parameters and the state at the previous time step.

Author
Knut Morten Okstad
Date
6 Jul 2006
Here is the caller graph for this function:

◆ buildsupnewtonmat()

subroutine supelroutinesmodule::buildsupnewtonmat ( logical, intent(in)  newTangent,
real(dp), intent(in)  scaleM,
real(dp), intent(in)  scaleC,
real(dp), intent(in)  scaleK,
type(supeltype), intent(inout)  sup,
integer, intent(out)  ierr 
)

Computes the superelement Newton matrix.

Parameters
[in]newTangentIf .true., the tangential stiffness is updated
[in]scaleMMass matrix scaling factor
[in]scaleCDamping matrix scaling factor
[in]scaleKStiffness matrix scaling factor
supThe superelement to calculate tangent stiffness for
[out]ierrError flag

The Newton matrix is a linear combination of the stiffness-, damping- (if any) and stiffness matrices of the superelement.

Author
Knut Morten Okstad
Date
10 Apr 2019
Here is the call graph for this function:
Here is the caller graph for this function:

◆ calcbuoyancyforces()

subroutine, private supelroutinesmodule::calcbuoyancyforces ( type(supeltype), intent(inout)  sup,
type(environmenttype), intent(inout)  env,
real(dp), intent(in)  time,
integer, intent(in)  iter,
integer, intent(inout)  ierr 
)
private

Calculates buoyancy forces and associated load correction stiffness.

Parameters
supSuperelement to calculate buoyancy forces for
envEnvironmental data
[in]timeCurrent simulation time
[in]iterIteration counter
ierrError flag

For a (partly) submerged body.

Author
Knut Morten Okstad
Date
2 Jul 2008
Here is the call graph for this function:
Here is the caller graph for this function:

◆ calcmorisonforces()

subroutine, private supelroutinesmodule::calcmorisonforces ( type(supeltype), intent(inout)  sup,
type(environmenttype), intent(inout)  env,
real(dp), intent(in)  time,
integer, intent(in)  istep,
integer, intent(in)  iter,
integer, intent(inout)  ierr 
)
private

Calculates Morison force contributions for a two-noded beam element.

Parameters
supThe beam element to calculate Morison forces for
envEnvironmental data
[in]timeCurrent simulation time
[in]istepTime increment counter
[in]iterIteration counter
ierrError flag
Author
Knut Morten Okstad
Date
24 Sep 2019
Here is the call graph for this function:
Here is the caller graph for this function:

◆ comptanstiff()

subroutine supelroutinesmodule::comptanstiff ( type(supeltype), intent(inout)  sup,
integer, intent(out)  ierr 
)

Computes the tangential superelement stiffness matrix.

Parameters
supThe superelement to calculate tangent stiffness for
[out]ierrError flag
Author
Dag Rune Christensen
Date
23 Oct 1997
Author
Bjorn Haugen
Date
3 Oct 2003
Author
Knut Morten Okstad
Date
10 Apr 2019
Here is the call graph for this function:
Here is the caller graph for this function:

◆ incsupelsgendofs()

subroutine supelroutinesmodule::incsupelsgendofs ( type(supeltype), dimension(:), intent(inout)  sups,
real(dp), dimension(:), intent(in)  solinc,
logical, intent(in), optional  useTotalInc 
)

Increments the generalized DOFs for all superelements.

Parameters
supsAll superelements in the model
[in]solincCurrent global solution increment
[in]useTotalIncIf .true., update from previous solution state
Author
Bjorn Haugen
Date
16 Nov 2001
Here is the caller graph for this function:

◆ scaledmatmul()

subroutine, private supelroutinesmodule::scaledmatmul ( integer, intent(in)  m,
integer, intent(in)  n,
real(dp), intent(in)  alpha,
real(dp), dimension(:,:), intent(in)  A,
real(dp), dimension(:), intent(in)  X,
real(dp), dimension(:), intent(out)  Y,
integer, intent(in), optional  ldA 
)
private

Calculates the scaled matrix-vector product Y = α*A*X.

Parameters
[in]mNumber of rows in matrix A
[in]nNumber of columns in matrix A and length of vector X.
[in]alphaScaling factor, α
[in]AThe matrix to multiply with
[in]XThe vector to multiply with
[out]YThe output vector
[in]ldALeading dimension of matrix A

This is just a convenience wrapper over the BLAS subroutine DGEMV.

Author
Knut Morten Okstad
Date
21 Oct 2019
Here is the caller graph for this function:

◆ setsupelsvelacc()

subroutine supelroutinesmodule::setsupelsvelacc ( type(supeltype), dimension(:), intent(inout)  sups,
type(samtype), intent(in)  sam,
real(dp), dimension(:), intent(in)  velGlobal,
real(dp), dimension(:), intent(in)  accGlobal 
)

Extracts local velocities and accelerations for the superelements.

Parameters
supsAll superelements in the model
[in]samData for managing system matrix assembly
[in]velGlobalGlobal velocity vector
[in]accGlobalGlobal acceleration vector

This subroutine extracts the local velocities and accelerations (uld and uldd, respectively) from the given global velocity and acceleration vectors for all superelements.

Author
Bjorn Haugen
Date
16 Nov 2001
Here is the caller graph for this function:

◆ updateseaenvironment()

subroutine supelroutinesmodule::updateseaenvironment ( type(environmenttype), intent(inout)  env,
type(triadtype), dimension(:), intent(in)  triads,
type(supeltype), dimension(:), intent(in)  sups,
real(dp), intent(in)  time,
integer, intent(in)  iter,
integer, intent(out)  ierr 
)

Updates the current sea state.

Parameters
envEnvironmental data
[in]triadsAll triads in the model
[in]supsAll superelements in the model
[in]timeCurrent physical time
[in]iterIteration counter
[out]ierrError flag
Author
Knut Morten Okstad
Date
23 Oct 2008
Here is the call graph for this function:
Here is the caller graph for this function:

◆ updatesupeldamping()

subroutine supelroutinesmodule::updatesupeldamping ( type(supeltype), dimension(:), intent(inout)  sups,
type(enginetype), dimension(:), intent(inout)  engs,
real(dp), dimension(2), intent(in)  alpha,
integer, intent(out)  ierr 
)

Updates the superelement damping matrices.

Parameters
supsAll superelements in the model
engsAll general functions in the model
[in]alphaGlobal structural damping parameters
[out]ierrError flag

Effective with time-dependent structural damping coefficients. This subroutine also updates the time-dependent stiffness scaling factor.

Author
Knut Morten Okstad
Date
1 Dec 2010
Here is the call graph for this function:
Here is the caller graph for this function:

◆ updatesupelload()

subroutine, private supelroutinesmodule::updatesupelload ( real(dp), dimension(:), intent(inout)  Q,
real(dp), dimension(:,:), intent(in)  S,
type(supelloadtype), intent(inout)  supLoad,
integer, intent(inout)  ierr 
)
private

Adds the superelement loads to the system external load vector.

Parameters
QSystem external load vector
SNormalized superelement load vectors
supLoadThe superelement load object to assemble
ierrError flag

This subroutine updates the amplitude of the superelement load and adds its contribution to the external force vector of the global system.

Author
Knut Morten okstad
Date
8 Apr 2008
Here is the call graph for this function:
Here is the caller graph for this function:

◆ updatesupels()

subroutine supelroutinesmodule::updatesupels ( type(supeltype), dimension(:), intent(inout)  sups,
type(supelloadtype), dimension(:), intent(inout)  supLoads,
type(environmenttype), intent(inout)  env,
real(dp), intent(in)  beta,
real(dp), intent(in)  gamma,
real(dp), intent(in)  time,
real(dp), intent(in)  timeStep,
integer, intent(in)  istep,
integer, intent(in)  iter,
logical, intent(in)  newPositions,
integer, intent(inout)  ierr 
)

Updates all superelements in the model based on the computed state.

Parameters
supsAll superelements in the model
supLoadsAll superelement loads in the model
envEnvironmental data
[in]betaNewmark time integration parameter
[in]gammaNewmark time integration parameter
[in]timeCurrent physical time
[in]timeStepCurrent time increment size
[in]istepTime increment counter
[in]iterIteration counter
[in]newPositionsIf .true., the position variables have been updated
ierrError flag

This subroutine updates the superelement objects:

  • the deflection vector, Finit, is established
  • superelement forces are calculated and transformed to system directions
  • the total nodal forces are calculated in system directions
Author
Karl Erik Thoresen
Date
Sep 1999
Author
Knut Morten Okstad
Date
Nov 2001
Author
Knut Morten Okstad
Date
Feb 2002
Author
Knut Morten Okstad
Date
Apr 2008
Author
Knut Morten Okstad
Date
Aug 2011
Here is the call graph for this function:
Here is the caller graph for this function:

◆ updatesupelsstatic()

subroutine supelroutinesmodule::updatesupelsstatic ( type(supeltype), dimension(:), intent(inout)  sups,
type(supelloadtype), dimension(:), intent(inout)  supLoads,
type(environmenttype), intent(inout)  env,
real(dp), intent(in)  time,
integer, intent(in)  iter,
integer, intent(in), optional  linInc,
integer, intent(inout)  ierr 
)

Updates all superelements in the model based on the computed state.

Parameters
supsAll superelements in the model
supLoadsAll superelement loads in the model
envEnvironmental data
[in]timeCurrent physical time
[in]iterIteration counter
[in]linIncOption for linear (1) or incremental (2) deformations
ierrError flag

This subroutine updates the superelement objects statically. It is used in place of supelroutinesmodule::updatesupels in a quasi-static simulation.

Author
Knut Morten Okstad
Date
5 June 2002
Here is the call graph for this function:
Here is the caller graph for this function: