FEDEM Solver  R8.0
Source code of the dynamics solver
Functions
solverInterface.h File Reference

This file defines the API for the FEDEM dynamics solver core. More...

#include <stddef.h>
#include <stdbool.h>
Include dependency graph for solverInterface.h:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Functions

int solverInit (int argc, char **argv, const char *cfsi=NULL, const double *stateData=NULL, const int ndat=0, const double *gagesData=NULL, const int ngda=0, const double *extInput=NULL, int *nxinp=NULL)
 Initializes a new solver process. More...
 
int restartFromState (const double *stateData, const int ndat, const int writeToRDB=2)
 Restarts a simulation process from the provided state. More...
 
bool solveWindow (int nStep, int nInc, int nIn, int nOut, const int *fOut, const double *times, const double *inputs, double *outputs, int nDat, double *stateData, int *ierr)
 Executes the simulation for a specified time window. More...
 
int getStateSize ()
 Returns the required length of the state vector to be is used when restarting a simulation from an in-core array. More...
 
int getTransformationStateSize ()
 Returns the required length of the state vector which contains translation and rotation information for all triads and superelements. More...
 
int getPartDeformationStateSize (int bid)
 Returns the required length of the state vector which contains deformation information for the specified FE part. More...
 
int getPartStressStateSize (int bid)
 Returns the required length of the state vector which contains von Mises stress information for the specified FE part. More...
 
int getGagesSize ()
 Returns the required length of the initial strain gage vector which is used when restarting a simulation from an in-core array. More...
 
bool saveState (double *stateData, const int ndat)
 Saves current state to the provided array. More...
 
bool saveTransformationState (double *stateData, const int ndat)
 Saves current transformation state to the provided core array. More...
 
bool saveGages (double *stateData, const int ndat)
 Saves the initial gage strains to the provided core array. More...
 
bool getStrainsFromDisp (const double *disp, const int *gageIDs, double *eps, const int ndof, const int ng)
 Gets gage strain values from given displacement field. More...
 
bool getBeamForcesFromDisp (const double *disp, const int *beamIDs, double *forces, const int ndof, const int nBeam)
 Gets internal sectional forces from given displacements. More...
 
bool getRelDisp (const double *disp, const int *Ids, double *relDis, const int ndof, const int nId)
 Gets relative distances between triads from given displacements. More...
 
bool getRespVars (const double *disp, const int *Ids, double *resVar, const int ndof, const int nId)
 Gets response variable values from given displacements. More...
 
bool savePartDeformationState (int bid, double *stateData, const int ndat)
 Saves the deformation of an FE part to the provided core array. More...
 
bool savePartStressState (int bid, double *stateData, const int ndat)
 Saves the stress state of an FE part to the provided core array. More...
 
bool solveNext (int *ierr)
 Advances the solution one time/load step forward. More...
 
bool startStep (int *ierr)
 Starts a new time/load step. More...
 
bool solveIteration (int *ierr, bool all=false)
 Solves current linearized equation system and updates the state. More...
 
bool solveEigenModes (const int nmod, double *eval, double *evec, bool dofOrder, int useLaPack, int *ierr)
 Solves the eigenvalue problem at current state. More...
 
bool solveInverse (const double *x, const int *xeqs, const int *feqs, int ndis, int nfrs, int *ierr)
 Solves the inverse problem at current time/load step. More...
 
int solverDone (bool removeSingletons=true)
 Closes down current model and cleans up core memory and on disk. More...
 
void solverClose ()
 Cleans up singleton objects in core memory. More...
 
bool setExtFunc (int funcId, double value)
 Assigns a sensor value from a physical twin to the simulation model. More...
 
double getTime (bool nextStep, int *ierr)
 Returns current or next physical time of the simulation. More...
 
bool setTime (double nextTime)
 Sets the time increment size of the next time step. More...
 
double getCurrentTime (int *ierr)
 Returns the current physical time of the simulation. More...
 
double getNextTime (int *ierr)
 Returns the physical time of the next step of the simulation. More...
 
double evalFunc (int uid, const char *tag=NULL, double x=-1.0, int *ierr=NULL)
 Evaluates a specified general function in the model. More...
 
double getFunc (int uid, int *ie=NULL)
 Returns the current value of the specified general function. More...
 
int getFuncId (const char *tag)
 Returns the user ID of a tagged function. More...
 
int getEquations (int bid, int *meqn)
 Returns the equation numbers for the DOFs of a specified object. More...
 
int getStateVar (int bid, double *var)
 Returns the current state variables for a specified object. More...
 
int getSystemSize (bool dofs=false)
 Returns the dimension (number of equations) of the system. More...
 
bool getSystemMatrix (double *Smat, int iMat=0)
 Returns current content of the system Newton matrix. More...
 
bool getNewtonMatrix (double *Nmat)
 Returns current content of the system Newton matrix. More...
 
bool getStiffnessMatrix (double *Kmat)
 Returns current content of the system stiffness matrix. More...
 
bool getMassMatrix (double *Mmat)
 Returns current content of the system mass matrix. More...
 
bool getDampingMatrix (double *Cmat)
 Returns current content of the system damping matrix. More...
 
bool getElementStiffnessMatrix (double *Kmat, int bid)
 Returns the content of an element stiffness matrix. More...
 
bool getRhsVector (double *Rvec, int iVec=0)
 Returns current content of the system right-hand-side vector. More...
 
bool setRhsVector (const double *Rvec)
 Replaces current content of the system right-hand-side vector. More...
 
bool addRhsVector (const double *Rvec)
 Updates the content of the system right-hand-side vector. More...
 
bool getJointSprCoeff (double *sprCoeff, int bid)
 Get joint spring stiffness coefficients. More...
 

Detailed Description

This file defines the API for the FEDEM dynamics solver core.

Each function of the dynamics solver shared library that is supposed to be accessible for outside applications need to have their declaration in this file.

Author
Knut Morten Okstad, Fedem Technology AS
Date
20 Dec 2016

Function Documentation

◆ addRhsVector()

bool addRhsVector ( const double *  Rvec)

Updates the content of the system right-hand-side vector.

Parameters
[in]RvecVector to be added to the current right-hand-side vector

◆ evalFunc()

double evalFunc ( int  uid,
const char *  tag = NULL,
double  x = -1.0,
int *  ierr = NULL 
)

Evaluates a specified general function in the model.

Parameters
[in]uidUser ID of the general function to be evaluated
[in]tagAlternative function identification
[in]xIf zero or positive, this is the function argument value. If negative, the value is ignored and the function is evaluated for the current state of the sensor object that is defined as the function argument. If tag is specified, the value is ignored (same behavior as if negative).
[out]ierrEqual to zero on a successful call, a non-zero value indicates that the function could not be evaluated
Returns
The evaluated function value

◆ getBeamForcesFromDisp()

bool getBeamForcesFromDisp ( const double *  disp,
const int *  beamIDs,
double *  forces,
const int  ndof,
const int  nBeam 
)

Gets internal sectional forces from given displacements.

Parameters
[in]dispDisplacement vector
[in]beamIDsArray of beam IDs
[out]forcesForce vector at the sections
[in]ndofNumber of degrees of freedom
[in]nBeamNumber of beams

◆ getCurrentTime()

double getCurrentTime ( int *  ierr)
inline

Returns the current physical time of the simulation.

Parameters
[out]ierrAlways equal to zero

◆ getDampingMatrix()

bool getDampingMatrix ( double *  Cmat)

Returns current content of the system damping matrix.

Parameters
[out]CmatThe system damping matrix in full format

◆ getElementStiffnessMatrix()

bool getElementStiffnessMatrix ( double *  Kmat,
int  bid 
)

Returns the content of an element stiffness matrix.

Parameters
[out]KmatThe element stiffness matrix
[in]bidBase ID of the superelement (beam, reduced FE part or generic part) to return the stiffness matrix for

◆ getEquations()

int getEquations ( int  bid,
int *  meqn 
)

Returns the equation numbers for the DOFs of a specified object.

Parameters
[in]bidBase ID of the mechanism object to get equations for
[out]meqnList of equation numbers (1-based)
Returns
   0 : Non-existing object or object with no DOFs
> 0 : Number of DOFs for the specified object (length of meqn)

◆ getFunc()

double getFunc ( int  uid,
int *  ie = NULL 
)

Returns the current value of the specified general function.

Parameters
[in]uidUser ID of the general function to be evaluated
[out]ieEqual to zero on a successful call, a non-zero value indicates that the function could not be evaluated
Returns
The evaluated function value

◆ getFuncId()

int getFuncId ( const char *  tag)

Returns the user ID of a tagged function.

◆ getGagesSize()

int getGagesSize ( )

Returns the required length of the initial strain gage vector which is used when restarting a simulation from an in-core array.

◆ getJointSprCoeff()

bool getJointSprCoeff ( double *  sprCoeff,
int  bid 
)

Get joint spring stiffness coefficients.

Parameters
[out]sprCoeffspring coefficient
[in]bidBase ID of joint

◆ getMassMatrix()

bool getMassMatrix ( double *  Mmat)

Returns current content of the system mass matrix.

Parameters
[out]MmatThe system mass matrix in full format

◆ getNewtonMatrix()

bool getNewtonMatrix ( double *  Nmat)

Returns current content of the system Newton matrix.

Parameters
[out]NmatThe system Newton matrix in full format

◆ getNextTime()

double getNextTime ( int *  ierr)
inline

Returns the physical time of the next step of the simulation.

Parameters
[out]ierrEqual to zero on a successful call, a non-zero value may occur only if the time step size in the model is defined via a general function that could not be evaluated.

◆ getPartDeformationStateSize()

int getPartDeformationStateSize ( int  bid)

Returns the required length of the state vector which contains deformation information for the specified FE part.

Parameters
[in]bidBase ID of the FE part to consider

◆ getPartStressStateSize()

int getPartStressStateSize ( int  bid)

Returns the required length of the state vector which contains von Mises stress information for the specified FE part.

Parameters
[in]bidBase ID of the FE part to consider

◆ getRelDisp()

bool getRelDisp ( const double *  disp,
const int *  Ids,
double *  relDis,
const int  ndof,
const int  nId 
)

Gets relative distances between triads from given displacements.

Parameters
[in]dispDisplacement vector
[in]IdsArray of engine IDs
[out]relDisRelative distance change
[in]ndofNumber of degrees of freedom
[in]nIdNumber of engine IDs

◆ getRespVars()

bool getRespVars ( const double *  disp,
const int *  Ids,
double *  resVar,
const int  ndof,
const int  nId 
)

Gets response variable values from given displacements.

Parameters
[in]dispDisplacement vector
[in]IdsArray of engine IDs
[out]resVarResponse variable values
[in]ndofNumber of degrees of freedom
[in]nIdNumber of engine IDs

◆ getRhsVector()

bool getRhsVector ( double *  Rvec,
int  iVec = 0 
)

Returns current content of the system right-hand-side vector.

Parameters
[out]RvecThe system right-hand-side vector
[in]iVecSystem vector type flag (0=residual or 5=external)

◆ getStateSize()

int getStateSize ( )

Returns the required length of the state vector to be is used when restarting a simulation from an in-core array.

◆ getStateVar()

int getStateVar ( int  bid,
double *  var 
)

Returns the current state variables for a specified object.

Parameters
[in]bidBase ID of the mechanism object to get statevariables for
[out]varList of variable values
Returns
   0 : Non-existing object or object with no DOFs
> 0 : Number of state variables for the specified object

This function is provided to facilitate unit testing where we want to assert on certain state variables. For triads, only the updated positions are returned (no rotation quantities).

◆ getStiffnessMatrix()

bool getStiffnessMatrix ( double *  Kmat)

Returns current content of the system stiffness matrix.

Parameters
[out]KmatThe system stiffness matrix in full format

◆ getStrainsFromDisp()

bool getStrainsFromDisp ( const double *  disp,
const int *  gageIDs,
double *  eps,
const int  ndof,
const int  ng 
)

Gets gage strain values from given displacement field.

Parameters
[in]dispDisplacement vector
[in]gageIDsArray of strain gage ID numbers
[out]epsStrain tensor at the strain gages
[in]ndofLength of disp (number of degrees of freedom)
[in]ngLength of gageIDs (number of gages)

◆ getSystemMatrix()

bool getSystemMatrix ( double *  Smat,
int  iMat = 0 
)

Returns current content of the system Newton matrix.

Parameters
[out]SmatThe system matrix in full format
[in]iMatSystem matrix type flag (0, 1, 2 or 3)

◆ getSystemSize()

int getSystemSize ( bool  dofs = false)

Returns the dimension (number of equations) of the system.

Parameters
[in]dofsIf true, return the total number of DOFs in the system

◆ getTime()

double getTime ( bool  nextStep,
int *  ierr 
)

Returns current or next physical time of the simulation.

Parameters
[in]nextStepIf true, return next time, otherwise current
[out]ierrError flag

◆ getTransformationStateSize()

int getTransformationStateSize ( )

Returns the required length of the state vector which contains translation and rotation information for all triads and superelements.

◆ restartFromState()

int restartFromState ( const double *  stateData,
const int  ndat,
const int  writeToRDB = 2 
)

Restarts a simulation process from the provided state.

Parameters
[in]stateDataComplete state vector to restart simulation from
[in]ndatLength of the restart state vector
[in]writeToRDBFlag for saving response variables (see below)
Returns
   0 : Success
< 0 : An error occurred

This function re-initializes the mechanism elements with data from the provided state array, such that the simulation can continue from that state.

The value of the input flag writeToRDB is interpreted as follows:

  • 0: No results saving
  • 1: Append results to the already opened results database
  • 2: Increment the results database and write to new files

◆ saveGages()

bool saveGages ( double *  stateData,
const int  ndat 
)

Saves the initial gage strains to the provided core array.

Parameters
[out]stateDataGage strain state vector
[in]ndatLength of the state vector

◆ savePartDeformationState()

bool savePartDeformationState ( int  bid,
double *  stateData,
const int  ndat 
)

Saves the deformation of an FE part to the provided core array.

Parameters
[in]bidBase ID of the FE part to consider
[out]stateDataDeformation state vector
[in]ndatLength of the state vector

This function is used to save current deformation state for the specified FE part to the provided core array after a successful solveNext call, for access by external visualisation modules.

◆ savePartStressState()

bool savePartStressState ( int  bid,
double *  stateData,
const int  ndat 
)

Saves the stress state of an FE part to the provided core array.

Parameters
[in]bidBase ID of the FE part to consider
[out]stateDataStress state vector
[in]ndatLength of the state vector

This function is used to save current von Mises stress state for the specified FE part to the provided core array after a successful solveNext call, for access by external visualisation modules.

◆ saveState()

bool saveState ( double *  stateData,
const int  ndat 
)

Saves current state to the provided array.

Parameters
[out]stateDataComplete state vector
[in]ndatLength of the state vector

This function can be used after a succesful solveWindow call, or after a solveNext call returning false indicating the end of the simulation, to save current state to an in-core array that can be used to restart the simulation later on.

◆ saveTransformationState()

bool saveTransformationState ( double *  stateData,
const int  ndat 
)

Saves current transformation state to the provided core array.

Parameters
[out]stateDataTransformation state vector
[in]ndatLength of the state vector

◆ setExtFunc()

bool setExtFunc ( int  funcId,
double  value 
)

Assigns a sensor value from a physical twin to the simulation model.

Parameters
[in]funcIdID of the external function to receive the sensor value
[in]valueThe actual sensor value to be assigned
Returns
false if funcId is out of range, otherwise true

This function may be used prior to the solveNext call, to assign a sensor value from a physical twin to the specified actuator or load in the digital twin model.

◆ setRhsVector()

bool setRhsVector ( const double *  Rvec)

Replaces current content of the system right-hand-side vector.

Parameters
[in]RvecThe new system right-hand-side vector

◆ setTime()

bool setTime ( double  nextTime)

Sets the time increment size of the next time step.

Parameters
[in]nextTimePhysical time to calculate next time increment size from

◆ solveEigenModes()

bool solveEigenModes ( const int  nmod,
double *  eval,
double *  evec,
bool  dofOrder,
int  useLaPack,
int *  ierr 
)

Solves the eigenvalue problem at current state.

Parameters
[in]nmodNumber of eigenmodes to calculate
[out]evalComputed eigenvalues
[out]evecComputed eigenvectors
[in]dofOrderIf true, return eigenvectors in DOF-order. Otherwise, they are returned in equation order
[in]useLaPackFlag usage of LAPACK eigensolvers (1=DSYGVX, 2=DGGEVX)
[out]ierrEqual to zero on a successful computation, a negative value indicates memory allocation error or similar issues, positive value indicate that the eigenvalue solver failed, but the execution may continue
Returns
true, unless the simulation has to stop due to an error condition

◆ solveInverse()

bool solveInverse ( const double *  x,
const int *  xeqs,
const int *  feqs,
int  ndis,
int  nfrs,
int *  ierr 
)

Solves the inverse problem at current time/load step.

Parameters
[in]xSpecified displacement values at a set of degrees of freedom
[in]xeqsEquation numbers (1-based) for the specified displacements
[in]feqsEquation numbers (1-based) for the DOFs with unknown external forces
[in]ndisNumber of DOFs with specified displacements
[in]nfrsNumber of DOFs with unknown external forces
[out]ierrEqual to zero on a successful computation, a non-zero value indicates some error that requires the simulation to terminate
Returns
true, unless the simulation has to stop due to some error condition, or the end of the simulation has been reached

◆ solveIteration()

bool solveIteration ( int *  ierr,
bool  all = false 
)

Solves current linearized equation system and updates the state.

Parameters
[out]ierrEqual to zero on a successful computation, a non-zero value indicates some error that requires to simulation to terminate
[in]allIf true, repeat until convergence
Returns
true, unless the simulation has to stop, either because of an error condition, or (if all is false) the current time/load step has reached convergence

This function solves the current linearized equation system and updates all state variables in the model. Then it assembles the linearized system of equations for next iteration, unless convergence has been reached. if all is true, this is repeated until convergence is achieved.

◆ solveNext()

bool solveNext ( int *  ierr)

Advances the solution one time/load step forward.

Parameters
[out]ierrEqual to zero on a successful computation, a non-zero value indicates some error that requires to simulation to terminate
Returns
true, unless current time/load step failed to converge or the end time of the simulation has been reached

◆ solverClose()

void solverClose ( )

Cleans up singleton objects in core memory.

This function needs to be invoked only if the solveDone() call was performed with removeSingletons = false, before finishing.

◆ solverDone()

int solverDone ( bool  removeSingletons = true)

Closes down current model and cleans up core memory and on disk.

Parameters
[in]removeSingletonsIf false, singleton objects are not deleted
Returns
Zero on success, a non-zero value indicates some error

This function should be invoked only when the time (or load) step loop is finished. If the process is going to be restarted later on the same model, the removeSingletons argument should be true, to ensure that objects allocated on the heap only on program startup are not deleted.

◆ solverInit()

int solverInit ( int  argc,
char **  argv,
const char *  cfsi = NULL,
const double *  stateData = NULL,
const int  ndat = 0,
const double *  gagesData = NULL,
const int  ngda = 0,
const double *  extInput = NULL,
int *  nxinp = NULL 
)

Initializes a new solver process.

Parameters
[in]argcNumber of command-line arguments passed to the solver
[in]argvList of command-line arguments passed to the solver
[in]cfsiContent of the solver input file describing the model
[in]stateDataComplete state vector to restart simulation from
[in]ndatLength of the restart state vector
[in]gagesDataInitial strain gage values for restart
[in]ngdaLength of the gages array
[in]extInputInitial external function values
[in,out]nxinpLength of the extInput array
Returns
≥ 0 : Success
< 0 : An error has occured

This function processes the input and sets up necessary data structures prior to the time integration loop. It also performs the license checks. See the Fedem R7.4 Users Guide, Appendix C.3 for a complete list of all command-line arguments that may be specified and their default values.

If the argument cfsi is NULL, the model is assumed defined in the file specified via the command-line argument -fsifile instead.

The arguments extInput and nxinp are used to transfer initial values of the external functions to the solver before the initial equilibrirum iterations is performed. The argument nxinp is reset to zero on exit, if initial equilibrium iterations is not performed, and to a negative valiue if the array extInput is too short.

The return value is zero on success, a negative value indicates some error. If cfsi is NULL and stateData is NULL, the return value is the required minimum length of the array to be used as state vector for restart.

◆ solveWindow()

bool solveWindow ( int  nStep,
int  nInc,
int  nIn,
int  nOut,
const int *  fOut,
const double *  times,
const double *  inputs,
double *  outputs,
int  nDat,
double *  stateData,
int *  ierr 
)

Executes the simulation for a specified time window.

Parameters
[in]nStepNumber of time/load steps to solve for from current state
[in]nIncNumber of explicit time increments
[in]nInNumber of input sensor values
[in]nOutNumber of output sensor values
[in]fOutList of user IDs identifying the output sensors in the model
[in]timesExplicit times to solver for (dimension nInc)
[in]inputsList of input sensor values (dimension nIn*nStep)
[out]outputsList of output sensor values (dimension nOut*nStep)
[in]nDatLength of the restart state vector
[out]stateDataComplete state vector at the end of the time window
[out]ierrA non-zero value indicates some error condition that requires the simulation to be terminated
Returns
true, unless the end of the simulation has been reached

This function solves the problem for a time/load step window, with given values for the external functions (physical sensor readings), and extraction of results from another set of general functions in the model.

◆ startStep()

bool startStep ( int *  ierr)

Starts a new time/load step.

Parameters
[out]ierrEqual to zero on a successful computation, a non-zero value indicates some error that requires to simulation to terminate
Returns
true, unless the simulation has to stop, either because of an error condition, or because the end time of the simulation has been reached

This function starts a new time (or load) step, by calculating the predicted response and the coefficient matrix and right-hand-side vector of the first nonlinear iteration. It has to be followed up by a series of solveIteration calls in order to continue the simulation, but the linear equation system can be manipulated in between.