FEDEM Solver  R8.0
Source code of the dynamics solver
Data Types | Functions/Subroutines
supeltypemodule Module Reference

Module with data types representing superelement objects. More...

Data Types

type  generalizeddofs
 Data type for the generalized DOFs associated with component modes. More...
 
type  nonlinforcestifftype
 Data type for the nonlinear force-displacement representation. More...
 
type  hydrodyntype
 Data type for the hydrodynamic force calculation. More...
 
type  recoverytype
 Data type for the integrated stress recovery. More...
 
type  supeltype
 Data type representing a superelement object. More...
 
type  supelptrtype
 Data type representing a superelement pointer. More...
 
interface  getptrtoid
 Returns pointer to object with specified ID. More...
 
interface  writeobject
 Standard routine for writing an object to file. More...
 
interface  updateatconvergence
 Updates the state variables pertaining to previous time step. More...
 
interface  restorefromlaststep
 Restores the state variables from the last converged time step. More...
 

Functions/Subroutines

type(supeltype) function, pointer, private getptrtoidsupel (array, id, userId)
 Returns pointer to (first) superelement with specified ID. More...
 
subroutine, private writesupeltype (sup, io, complexity)
 Standard routine for writing an object to file. More...
 
subroutine nullifysupel (sup)
 Initializes a superelement object. More...
 
subroutine initiatehydyn (hydyn, MorPar)
 Initiates the hydrodynamics quantities for a superelement. More...
 
subroutine transformoffset (sup, offset, index)
 Transforms a reference triad offset vector. More...
 
subroutine deallocatesupel (sup)
 Deallocates a superelement object. More...
 
subroutine deallocatesupels (sups)
 Deallocates an array of superelement objects. More...
 
subroutine updatesupelcorot (sup, dbgUnit, stat)
 Updates the co-rotated system for a superelement. More...
 
subroutine, private updatepreviousvalues (sup)
 Updates the state variables pertaining to the previous time step. More...
 
subroutine, private restorepreviousvalues (sup)
 Restores the state variables from the last converged time step. More...
 
subroutine buildfinit (sup)
 Establishes the deformation vector (finit) for a superelement. More...
 
real(dp) function, dimension(3, 3) getrsuptosys (sup, n)
 Returns the transformation matrix to system directions for a triad. More...
 
subroutine getsupelsvelacc (sups, velGlobal, accGlobal)
 Gets current velocity and acceleration from all superelements. More...
 
logical function isbeam (sup)
 Returns whether a superelement is a beam element or not. More...
 
character(len=5+lid_p) function getsupelid (sup)
 Returns the full ID of a superelement. More...
 

Detailed Description

Module with data types representing superelement objects.

The module also contains subroutines for accessing or manipulating the superelement data.

Function/Subroutine Documentation

◆ buildfinit()

subroutine supeltypemodule::buildfinit ( type(supeltype), intent(inout)  sup)

Establishes the deformation vector (finit) for a superelement.

Parameters
supThe supeltypemodule::supeltype object to obtain deformation for
Author
Karl Erik Thoresen
Date
3 Sep 1998
Here is the caller graph for this function:

◆ deallocatesupel()

subroutine supeltypemodule::deallocatesupel ( type(supeltype), intent(inout)  sup)

Deallocates a superelement object.

Parameters
supThe supeltypemodule::supeltype object to deallocate
Author
Knut Morten Okstad
Date
23 Jan 2017
Here is the call graph for this function:
Here is the caller graph for this function:

◆ deallocatesupels()

subroutine supeltypemodule::deallocatesupels ( type(supeltype), dimension(:), pointer  sups)

Deallocates an array of superelement objects.

Parameters
supsThe supeltypemodule::supeltype objects to deallocate
Author
Knut Morten Okstad

date 23 Jan 2017

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

◆ getptrtoidsupel()

type(supeltype) function, pointer, private supeltypemodule::getptrtoidsupel ( type(supeltype), dimension(:), intent(in), target  array,
integer, intent(in)  id,
logical, intent(in), optional  userId 
)
private

Returns pointer to (first) superelement with specified ID.

Parameters
[in]arrayArray of supeltypemodule::supeltype objects to search within
[in]idBase ID of the object to search for
[in]userIdIf .true., search for a user ID instead

If the superelement is not found, NULL is returned.

Author
Bjorn Haugen / Knut Morten Okstad
Date
13 Jan 2010

◆ getrsuptosys()

real(dp) function, dimension(3,3) supeltypemodule::getrsuptosys ( type(supeltype), intent(in)  sup,
integer, intent(in)  n 
)

Returns the transformation matrix to system directions for a triad.

Parameters
[in]supThe supeltypemodule::supeltype object to transform for
[in]n1-based index of the superelement triad to transform for

The returned 3x3 matrix is to be used to transform vectors from the superelement directions to system directions for the n'th connected triad of the superelement.

Author
Karl Erik Thoresen
Date
Dec 1998
Here is the caller graph for this function:

◆ getsupelid()

character(len=5+lid_p) function supeltypemodule::getsupelid ( type(supeltype), intent(in)  sup)

Returns the full ID of a superelement.

Parameters
[in]supThe supeltypemodule::supeltype object to get the ID for

The full ID of an object consists for of the type name of the object, the user ID and the user description.

Author
Knut Morten Okstad
Date
15 Sep 2017
Here is the caller graph for this function:

◆ getsupelsvelacc()

subroutine supeltypemodule::getsupelsvelacc ( type(supeltype), dimension(:), intent(in)  sups,
real(dp), dimension(:), intent(out)  velGlobal,
real(dp), dimension(:), intent(out)  accGlobal 
)

Gets current velocity and acceleration from all superelements.

Parameters
[in]supsAll supeltypemodule::supeltype objects in the model
[out]velGlobalSystem velocity vector
[out]accGlobalSystem acceleration vector
Author
Knut Morten Okstad
Date
1 Nov 2005
Here is the caller graph for this function:

◆ initiatehydyn()

subroutine supeltypemodule::initiatehydyn ( type(hydrodyntype), intent(out)  hydyn,
real(dp), dimension(:), intent(in), optional  MorPar 
)

Initiates the hydrodynamics quantities for a superelement.

Parameters
hydynHydrodynamics data container
[in]MorParMorison input parameters
Author
Knut Morten Okstad
Date
25 Mar 2009
Here is the caller graph for this function:

◆ isbeam()

logical function supeltypemodule::isbeam ( type(supeltype), intent(in)  sup)

Returns whether a superelement is a beam element or not.

Parameters
[in]supThe supeltypemodule::supeltype object to check
Author
Knut Morten Okstad
Date
2 Jul 2013
Here is the caller graph for this function:

◆ nullifysupel()

subroutine supeltypemodule::nullifysupel ( type(supeltype), intent(out)  sup)

Initializes a superelement object.

Parameters
supThe supeltypemodule::supeltype object to initialize
Author
Knut Morten Okstad
Date
13 Oct 2000
Here is the caller graph for this function:

◆ restorepreviousvalues()

subroutine, private supeltypemodule::restorepreviousvalues ( type(supeltype), intent(inout)  sup)
private

Restores the state variables from the last converged time step.

Parameters
supThe supeltypemodule::supeltype object to restore state for

This subroutine is invoked when doing iteration cut-back.

Author
Knut Morten Okstad
Date
28 Oct 2008

◆ transformoffset()

subroutine supeltypemodule::transformoffset ( type(supeltype), intent(inout)  sup,
real(dp), dimension(3), intent(in)  offset,
integer, intent(in)  index 
)

Transforms a reference triad offset vector.

Parameters
supThe supeltypemodule::supeltype object to transform for
[in]offsetThe reference point offset vector to transform
[in]indexWhich reference point (1, 2 or 3) to transform for

The offset vector is transformed from the superelement coordinate system to the coordinate system of the reference triad itself.

Author
Knut Morten Okstad
Date
13 Nov 2015
Here is the caller graph for this function:

◆ updatepreviousvalues()

subroutine, private supeltypemodule::updatepreviousvalues ( type(supeltype), intent(inout)  sup)
private

Updates the state variables pertaining to the previous time step.

Parameters
supThe supeltypemodule::supeltype object to update state for

This subroutine is invoked once after convergence as been reached.

Author
Knut Morten Okstad
Date
28 Oct 2008

◆ updatesupelcorot()

subroutine supeltypemodule::updatesupelcorot ( type(supeltype), intent(inout)  sup,
integer, intent(in)  dbgUnit,
integer, intent(inout)  stat 
)

Updates the co-rotated system for a superelement.

Parameters
supThe supeltypemodule::supeltype object to update for
dbgUnitFile unit number for debug print out
statStatus flag, negative value on output indicates an error

If stat equals 1 on input, the relative position matrix for the triangle is initiated (used only for shadowPosAlg == 1).

Author
Karl Erik Thoresen
Date
Sep 1999
Author
Bjorn Haugen
Date
Nov 2003
Here is the call graph for this function:
Here is the caller graph for this function:

◆ writesupeltype()

subroutine, private supeltypemodule::writesupeltype ( type(supeltype), intent(in)  sup,
integer, intent(in)  io,
integer, intent(in), optional  complexity 
)
private

Standard routine for writing an object to file.

Parameters
[in]supThe supeltypemodule::supeltype object to write
[in]ioFile unit number to write to
[in]complexityIf present, the value indicates the amount of print
Author
Karl Erik Thoresen
Date
27 Sep 1998