FEDEM Solver
R8.0
Source code of the dynamics solver
|
This file defines the API for the FEDEM dynamics solver core. More...
#include <stddef.h>
#include <stdbool.h>
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... | |
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.
bool addRhsVector | ( | const double * | Rvec | ) |
Updates the content of the system right-hand-side vector.
[in] | Rvec | Vector to be added to the current right-hand-side vector |
double evalFunc | ( | int | uid, |
const char * | tag = NULL , |
||
double | x = -1.0 , |
||
int * | ierr = NULL |
||
) |
Evaluates a specified general function in the model.
[in] | uid | User ID of the general function to be evaluated |
[in] | tag | Alternative function identification |
[in] | x | If 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] | ierr | Equal to zero on a successful call, a non-zero value indicates that the function could not be evaluated |
bool getBeamForcesFromDisp | ( | const double * | disp, |
const int * | beamIDs, | ||
double * | forces, | ||
const int | ndof, | ||
const int | nBeam | ||
) |
Gets internal sectional forces from given displacements.
[in] | disp | Displacement vector |
[in] | beamIDs | Array of beam IDs |
[out] | forces | Force vector at the sections |
[in] | ndof | Number of degrees of freedom |
[in] | nBeam | Number of beams |
|
inline |
Returns the current physical time of the simulation.
[out] | ierr | Always equal to zero |
bool getDampingMatrix | ( | double * | Cmat | ) |
Returns current content of the system damping matrix.
[out] | Cmat | The system damping matrix in full format |
bool getElementStiffnessMatrix | ( | double * | Kmat, |
int | bid | ||
) |
Returns the content of an element stiffness matrix.
[out] | Kmat | The element stiffness matrix |
[in] | bid | Base ID of the superelement (beam, reduced FE part or generic part) to return the stiffness matrix for |
int getEquations | ( | int | bid, |
int * | meqn | ||
) |
Returns the equation numbers for the DOFs of a specified object.
[in] | bid | Base ID of the mechanism object to get equations for |
[out] | meqn | List of equation numbers (1-based) |
double getFunc | ( | int | uid, |
int * | ie = NULL |
||
) |
Returns the current value of the specified general function.
[in] | uid | User ID of the general function to be evaluated |
[out] | ie | Equal to zero on a successful call, a non-zero value indicates that the function could not be evaluated |
int getFuncId | ( | const char * | tag | ) |
Returns the user ID of a tagged function.
int getGagesSize | ( | ) |
Returns the required length of the initial strain gage vector which is used when restarting a simulation from an in-core array.
bool getJointSprCoeff | ( | double * | sprCoeff, |
int | bid | ||
) |
Get joint spring stiffness coefficients.
[out] | sprCoeff | spring coefficient |
[in] | bid | Base ID of joint |
bool getMassMatrix | ( | double * | Mmat | ) |
Returns current content of the system mass matrix.
[out] | Mmat | The system mass matrix in full format |
bool getNewtonMatrix | ( | double * | Nmat | ) |
Returns current content of the system Newton matrix.
[out] | Nmat | The system Newton matrix in full format |
|
inline |
Returns the physical time of the next step of the simulation.
[out] | ierr | Equal 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. |
int getPartDeformationStateSize | ( | int | bid | ) |
Returns the required length of the state vector which contains deformation information for the specified FE part.
[in] | bid | Base ID of the FE part to consider |
int getPartStressStateSize | ( | int | bid | ) |
Returns the required length of the state vector which contains von Mises stress information for the specified FE part.
[in] | bid | Base ID of the FE part to consider |
bool getRelDisp | ( | const double * | disp, |
const int * | Ids, | ||
double * | relDis, | ||
const int | ndof, | ||
const int | nId | ||
) |
Gets relative distances between triads from given displacements.
[in] | disp | Displacement vector |
[in] | Ids | Array of engine IDs |
[out] | relDis | Relative distance change |
[in] | ndof | Number of degrees of freedom |
[in] | nId | Number of engine IDs |
bool getRespVars | ( | const double * | disp, |
const int * | Ids, | ||
double * | resVar, | ||
const int | ndof, | ||
const int | nId | ||
) |
Gets response variable values from given displacements.
[in] | disp | Displacement vector |
[in] | Ids | Array of engine IDs |
[out] | resVar | Response variable values |
[in] | ndof | Number of degrees of freedom |
[in] | nId | Number of engine IDs |
bool getRhsVector | ( | double * | Rvec, |
int | iVec = 0 |
||
) |
Returns current content of the system right-hand-side vector.
[out] | Rvec | The system right-hand-side vector |
[in] | iVec | System vector type flag (0=residual or 5=external) |
int getStateSize | ( | ) |
Returns the required length of the state vector to be is used when restarting a simulation from an in-core array.
int getStateVar | ( | int | bid, |
double * | var | ||
) |
Returns the current state variables for a specified object.
[in] | bid | Base ID of the mechanism object to get statevariables for |
[out] | var | List of variable values |
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).
bool getStiffnessMatrix | ( | double * | Kmat | ) |
Returns current content of the system stiffness matrix.
[out] | Kmat | The system stiffness matrix in full format |
bool getStrainsFromDisp | ( | const double * | disp, |
const int * | gageIDs, | ||
double * | eps, | ||
const int | ndof, | ||
const int | ng | ||
) |
Gets gage strain values from given displacement field.
[in] | disp | Displacement vector |
[in] | gageIDs | Array of strain gage ID numbers |
[out] | eps | Strain tensor at the strain gages |
[in] | ndof | Length of disp (number of degrees of freedom) |
[in] | ng | Length of gageIDs (number of gages) |
bool getSystemMatrix | ( | double * | Smat, |
int | iMat = 0 |
||
) |
Returns current content of the system Newton matrix.
[out] | Smat | The system matrix in full format |
[in] | iMat | System matrix type flag (0, 1, 2 or 3) |
int getSystemSize | ( | bool | dofs = false | ) |
Returns the dimension (number of equations) of the system.
[in] | dofs | If true, return the total number of DOFs in the system |
double getTime | ( | bool | nextStep, |
int * | ierr | ||
) |
Returns current or next physical time of the simulation.
[in] | nextStep | If true, return next time, otherwise current |
[out] | ierr | Error flag |
int getTransformationStateSize | ( | ) |
Returns the required length of the state vector which contains translation and rotation information for all triads and superelements.
int restartFromState | ( | const double * | stateData, |
const int | ndat, | ||
const int | writeToRDB = 2 |
||
) |
Restarts a simulation process from the provided state.
[in] | stateData | Complete state vector to restart simulation from |
[in] | ndat | Length of the restart state vector |
[in] | writeToRDB | Flag for saving response variables (see below) |
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:
bool saveGages | ( | double * | stateData, |
const int | ndat | ||
) |
Saves the initial gage strains to the provided core array.
[out] | stateData | Gage strain state vector |
[in] | ndat | Length of the state vector |
bool savePartDeformationState | ( | int | bid, |
double * | stateData, | ||
const int | ndat | ||
) |
Saves the deformation of an FE part to the provided core array.
[in] | bid | Base ID of the FE part to consider |
[out] | stateData | Deformation state vector |
[in] | ndat | Length 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.
bool savePartStressState | ( | int | bid, |
double * | stateData, | ||
const int | ndat | ||
) |
Saves the stress state of an FE part to the provided core array.
[in] | bid | Base ID of the FE part to consider |
[out] | stateData | Stress state vector |
[in] | ndat | Length 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.
bool saveState | ( | double * | stateData, |
const int | ndat | ||
) |
Saves current state to the provided array.
[out] | stateData | Complete state vector |
[in] | ndat | Length 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.
bool saveTransformationState | ( | double * | stateData, |
const int | ndat | ||
) |
Saves current transformation state to the provided core array.
[out] | stateData | Transformation state vector |
[in] | ndat | Length of the state vector |
bool setExtFunc | ( | int | funcId, |
double | value | ||
) |
Assigns a sensor value from a physical twin to the simulation model.
[in] | funcId | ID of the external function to receive the sensor value |
[in] | value | The actual sensor value to be assigned |
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.
bool setRhsVector | ( | const double * | Rvec | ) |
Replaces current content of the system right-hand-side vector.
[in] | Rvec | The new system right-hand-side vector |
bool setTime | ( | double | nextTime | ) |
Sets the time increment size of the next time step.
[in] | nextTime | Physical time to calculate next time increment size from |
bool solveEigenModes | ( | const int | nmod, |
double * | eval, | ||
double * | evec, | ||
bool | dofOrder, | ||
int | useLaPack, | ||
int * | ierr | ||
) |
Solves the eigenvalue problem at current state.
[in] | nmod | Number of eigenmodes to calculate |
[out] | eval | Computed eigenvalues |
[out] | evec | Computed eigenvectors |
[in] | dofOrder | If true, return eigenvectors in DOF-order. Otherwise, they are returned in equation order |
[in] | useLaPack | Flag usage of LAPACK eigensolvers (1=DSYGVX, 2=DGGEVX) |
[out] | ierr | Equal 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 |
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.
[in] | x | Specified displacement values at a set of degrees of freedom |
[in] | xeqs | Equation numbers (1-based) for the specified displacements |
[in] | feqs | Equation numbers (1-based) for the DOFs with unknown external forces |
[in] | ndis | Number of DOFs with specified displacements |
[in] | nfrs | Number of DOFs with unknown external forces |
[out] | ierr | Equal to zero on a successful computation, a non-zero value indicates some error that requires the simulation to terminate |
bool solveIteration | ( | int * | ierr, |
bool | all = false |
||
) |
Solves current linearized equation system and updates the state.
[out] | ierr | Equal to zero on a successful computation, a non-zero value indicates some error that requires to simulation to terminate |
[in] | all | If true, repeat until 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.
bool solveNext | ( | int * | ierr | ) |
Advances the solution one time/load step forward.
[out] | ierr | Equal to zero on a successful computation, a non-zero value indicates some error that requires to simulation to terminate |
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.
int solverDone | ( | bool | removeSingletons = true | ) |
Closes down current model and cleans up core memory and on disk.
[in] | removeSingletons | If false, singleton objects are not deleted |
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.
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.
[in] | argc | Number of command-line arguments passed to the solver |
[in] | argv | List of command-line arguments passed to the solver |
[in] | cfsi | Content of the solver input file describing the model |
[in] | stateData | Complete state vector to restart simulation from |
[in] | ndat | Length of the restart state vector |
[in] | gagesData | Initial strain gage values for restart |
[in] | ngda | Length of the gages array |
[in] | extInput | Initial external function values |
[in,out] | nxinp | Length of the extInput array |
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.
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.
[in] | nStep | Number of time/load steps to solve for from current state |
[in] | nInc | Number of explicit time increments |
[in] | nIn | Number of input sensor values |
[in] | nOut | Number of output sensor values |
[in] | fOut | List of user IDs identifying the output sensors in the model |
[in] | times | Explicit times to solver for (dimension nInc) |
[in] | inputs | List of input sensor values (dimension nIn*nStep) |
[out] | outputs | List of output sensor values (dimension nOut*nStep) |
[in] | nDat | Length of the restart state vector |
[out] | stateData | Complete state vector at the end of the time window |
[out] | ierr | A non-zero value indicates some error condition that requires the simulation to be terminated |
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.
bool startStep | ( | int * | ierr | ) |
Starts a new time/load step.
[out] | ierr | Equal to zero on a successful computation, a non-zero value indicates some error that requires to simulation to terminate |
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.