FEDEM Solver  R8.0
Source code of the dynamics solver
solverInterface.h
Go to the documentation of this file.
1 /* SPDX-FileCopyrightText: 2023 SAP SE
2  *
3  * SPDX-License-Identifier: Apache-2.0
4  *
5  * This file is part of FEDEM - https://openfedem.org
6  */
21 #ifndef _SOLVER_INTERFACE_H
22 #define _SOLVER_INTERFACE_H
23 
24 #include <stddef.h>
25 
26 #ifdef __cplusplus
27 extern "C" {
28 #else
29 #include <stdbool.h>
30 #endif
31 
65  int solverInit(int argc, char** argv, const char* cfsi = NULL,
66  const double* stateData = NULL, const int ndat = 0,
67  const double* gagesData = NULL, const int ngda = 0,
68  const double* extInput = NULL, int* nxinp = NULL);
69 
87  int restartFromState(const double* stateData, const int ndat,
88  const int writeToRDB = 2);
89 
110  bool solveWindow(int nStep, int nInc, int nIn, int nOut, const int* fOut,
111  const double* times, const double* inputs, double* outputs,
112  int nDat, double* stateData, int* ierr);
113 
119 
125 
132 
139 
145 
156  bool saveState(double* stateData, const int ndat);
157 
163  bool saveTransformationState(double* stateData, const int ndat);
164 
170  bool saveGages(double* stateData, const int ndat);
171 
180  bool getStrainsFromDisp(const double* disp, const int* gageIDs, double* eps,
181  const int ndof, const int ng);
182 
191  bool getBeamForcesFromDisp(const double* disp, const int* beamIDs,
192  double* forces, const int ndof, const int nBeam);
193 
202  bool getRelDisp(const double* disp, const int* Ids, double* relDis,
203  const int ndof, const int nId);
204 
213  bool getRespVars(const double* disp, const int* Ids, double* resVar,
214  const int ndof, const int nId);
215 
226  bool savePartDeformationState(int bid, double* stateData, const int ndat);
227 
238  bool savePartStressState(int bid, double* stateData, const int ndat);
239 
247  bool solveNext(int* ierr);
248 
262  bool startStep(int* ierr);
263 
278  bool solveIteration(int* ierr, bool all = false);
279 
293  bool solveEigenModes(const int nmod, double* eval, double* evec,
294  bool dofOrder, int useLaPack, int* ierr);
295 
309  bool solveInverse(const double* x,
310  const int* xeqs, const int* feqs,
311  int ndis, int nfrs, int* ierr);
312 
323  int solverDone(bool removeSingletons = true);
324 
330  void solverClose();
331 
342  bool setExtFunc(int funcId, double value);
343 
349  double getTime(bool nextStep, int* ierr);
350 
355  bool setTime(double nextTime);
356 
361  inline double getCurrentTime(int* ierr) { return getTime(false,ierr); }
362 
369  inline double getNextTime(int* ierr) { return getTime(true,ierr); }
370 
383  double evalFunc(int uid, const char* tag = NULL,
384  double x = -1.0, int* ierr = NULL);
385 
393  double getFunc(int uid, int* ie = NULL) { return evalFunc(uid,NULL,-1.0,ie); }
394 
398  int getFuncId(const char* tag);
399 
407  int getEquations(int bid, int* meqn);
408 
420  int getStateVar(int bid, double* var);
421 
426  int getSystemSize(bool dofs = false);
427 
433  bool getSystemMatrix(double* Smat, int iMat = 0);
434 
439  bool getNewtonMatrix(double* Nmat);
440 
445  bool getStiffnessMatrix(double* Kmat);
446 
451  bool getMassMatrix(double* Mmat);
452 
457  bool getDampingMatrix(double* Cmat);
458 
465  bool getElementStiffnessMatrix(double* Kmat, int bid);
466 
472  bool getRhsVector(double* Rvec, int iVec = 0);
473 
478  bool setRhsVector(const double* Rvec);
479 
484  bool addRhsVector(const double* Rvec);
485 
491  bool getJointSprCoeff(double* sprCoeff, int bid);
492 
493 #ifdef __cplusplus
494 }
495 #endif
496 
497 #endif
integer, save ndof
Definition: diffractionModule.f90:16
integer(ptr), save, private x
Definition: extCtrlSysRoutinesModule.f90:16
double getNextTime(int *ierr)
Returns the physical time of the next step of the simulation.
Definition: solverInterface.h:369
bool setExtFunc(int funcId, double value)
Assigns a sensor value from a physical twin to the simulation model.
bool getNewtonMatrix(double *Nmat)
Returns current content of the system Newton matrix.
bool savePartDeformationState(int bid, double *stateData, const int ndat)
Saves the deformation of an FE part to the provided core array.
int getStateVar(int bid, double *var)
Returns the current state variables for a specified object.
int getSystemSize(bool dofs=false)
Returns the dimension (number of equations) of the system.
int getTransformationStateSize()
Returns the required length of the state vector which contains translation and rotation information f...
double getCurrentTime(int *ierr)
Returns the current physical time of the simulation.
Definition: solverInterface.h:361
bool solveIteration(int *ierr, bool all=false)
Solves current linearized equation system and updates the state.
bool getElementStiffnessMatrix(double *Kmat, int bid)
Returns the content of an element stiffness matrix.
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.
bool getStiffnessMatrix(double *Kmat)
Returns current content of the system stiffness matrix.
bool getRhsVector(double *Rvec, int iVec=0)
Returns current content of the system right-hand-side vector.
bool getStrainsFromDisp(const double *disp, const int *gageIDs, double *eps, const int ndof, const int ng)
Gets gage strain values from given displacement field.
bool addRhsVector(const double *Rvec)
Updates the content of the system right-hand-side vector.
bool getBeamForcesFromDisp(const double *disp, const int *beamIDs, double *forces, const int ndof, const int nBeam)
Gets internal sectional forces from given displacements.
int restartFromState(const double *stateData, const int ndat, const int writeToRDB=2)
Restarts a simulation process from the provided state.
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.
bool getSystemMatrix(double *Smat, int iMat=0)
Returns current content of the system Newton matrix.
double evalFunc(int uid, const char *tag=NULL, double x=-1.0, int *ierr=NULL)
Evaluates a specified general function in the model.
void solverClose()
Cleans up singleton objects in core memory.
bool savePartStressState(int bid, double *stateData, const int ndat)
Saves the stress state of an FE part to the provided core array.
double getFunc(int uid, int *ie=NULL)
Returns the current value of the specified general function.
Definition: solverInterface.h:393
int getFuncId(const char *tag)
Returns the user ID of a tagged function.
int solverDone(bool removeSingletons=true)
Closes down current model and cleans up core memory and on disk.
bool getDampingMatrix(double *Cmat)
Returns current content of the system damping matrix.
bool setTime(double nextTime)
Sets the time increment size of the next time step.
int getPartStressStateSize(int bid)
Returns the required length of the state vector which contains von Mises stress information for the s...
int getEquations(int bid, int *meqn)
Returns the equation numbers for the DOFs of a specified object.
bool startStep(int *ierr)
Starts a new time/load step.
bool getMassMatrix(double *Mmat)
Returns current content of the system mass matrix.
int getStateSize()
Returns the required length of the state vector to be is used when restarting a simulation from an in...
bool getRelDisp(const double *disp, const int *Ids, double *relDis, const int ndof, const int nId)
Gets relative distances between triads from given displacements.
int getGagesSize()
Returns the required length of the initial strain gage vector which is used when restarting a simulat...
bool saveTransformationState(double *stateData, const int ndat)
Saves current transformation state to the provided core array.
bool solveNext(int *ierr)
Advances the solution one time/load step forward.
bool saveGages(double *stateData, const int ndat)
Saves the initial gage strains to the provided core array.
bool setRhsVector(const double *Rvec)
Replaces current content of the system right-hand-side vector.
int getPartDeformationStateSize(int bid)
Returns the required length of the state vector which contains deformation information for the specif...
double getTime(bool nextStep, int *ierr)
Returns current or next physical time of the simulation.
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.
bool saveState(double *stateData, const int ndat)
Saves current state to the provided array.
bool getRespVars(const double *disp, const int *Ids, double *resVar, const int ndof, const int nId)
Gets response variable values from given displacements.
bool solveEigenModes(const int nmod, double *eval, double *evec, bool dofOrder, int useLaPack, int *ierr)
Solves the eigenvalue problem at current state.
bool getJointSprCoeff(double *sprCoeff, int bid)
Get joint spring stiffness coefficients.