OOFEM 3.0
Loading...
Searching...
No Matches
oofem::NlDEIDynamic Class Reference

#include <nldeidynamic.h>

Inheritance diagram for oofem::NlDEIDynamic:
Collaboration diagram for oofem::NlDEIDynamic:

Public Member Functions

 NlDEIDynamic (int i, EngngModel *master=nullptr)
virtual ~NlDEIDynamic ()
void solveYourself () override
void solveYourselfAt (TimeStep *tStep) override
void updateYourself (TimeStep *tStep) override
double giveUnknownComponent (ValueModeType type, TimeStep *tStep, Domain *d, Dof *dof) override
void initializeFrom (InputRecord &ir) override
TimeStepgiveNextStep () override
 Returns next time step (next to current step) of receiver.
NumericalMethodgiveNumericalMethod (MetaStep *mStep) override
 Returns reference to receiver's numerical method.
void saveContext (DataStream &stream, ContextMode mode) override
void restoreContext (DataStream &stream, ContextMode mode) override
void printOutputAt (FILE *file, TimeStep *tStep) override
void printDofOutputAt (FILE *stream, Dof *iDof, TimeStep *tStep) override
const char * giveInputRecordName () const
const char * giveClassName () const override
 Returns class name of the receiver.
fMode giveFormulation () override
int giveNumberOfFirstStep (bool force=false) override
int estimateMaxPackSize (IntArray &commMap, DataStream &buff, int packUnpackType) override
Public Member Functions inherited from oofem::StructuralEngngModel
 StructuralEngngModel (int i, EngngModel *master=nullptr)
 Creates new StructuralEngngModel with number i, associated to domain d.
virtual ~StructuralEngngModel ()
 Destructor.
void updateYourself (TimeStep *tStep) override
int checkConsistency () override
void computeReaction (FloatArray &answer, TimeStep *tStep, int di)
void terminate (TimeStep *tStep) override
void buildReactionTable (IntArray &restrDofMans, IntArray &restrDofs, IntArray &eqn, TimeStep *tStep, int di)
void updateInternalRHS (FloatArray &answer, TimeStep *tStep, Domain *d, FloatArray *eNorm) override
void showSparseMtrxStructure (int type, oofegGraphicContext &gc, TimeStep *tStep) override
 EngngModel (int i, EngngModel *_master=NULL)
virtual ~EngngModel ()
 Destructor.
 EngngModel (const EngngModel &)=delete
EngngModeloperator= (const EngngModel &)=delete
DomaingiveDomain (int n)
void setDomain (int i, Domain *ptr, bool iDeallocateOld=true)
int giveNumberOfDomains ()
 Returns number of domains in problem.
const std::string & giveDescription () const
const time_t & giveStartTime ()
bool giveSuppressOutput () const
virtual ErrorEstimatorgiveDomainErrorEstimator (int n)
virtual MaterialInterfacegiveMaterialInterface (int n)
void setNumberOfEquations (int id, int neq)
FILE * giveOutputStream ()
 Returns file descriptor of output file.
std::string giveOutputBaseFileName ()
std::string giveReferenceFileName ()
void letOutputBaseFileNameBe (const std ::string &src)
ContextOutputMode giveContextOutputMode () const
int giveContextOutputStep () const
void setContextOutputMode (ContextOutputMode contextMode)
void setUDContextOutputMode (int cStep)
double giveDeltaT ()
 Returns time step size from the time step controlelr.
void setDeltaT (double dT)
 Returns time step size through the time step controlelr.
void setProblemMode (problemMode pmode)
void setParallelMode (bool newParallelFlag)
problemMode giveProblemMode () const
 Returns domain mode.
void setProblemScale (problemScale pscale)
problemScale giveProblemScale () const
 Returns scale in multiscale simulation.
virtual void setRenumberFlag ()
 Sets the renumber flag to true.
virtual void resetRenumberFlag ()
 Sets the renumber flag to false.
double giveSolutionStepTime ()
void giveAnalysisTime (int &rhrs, int &rmin, int &rsec, int &uhrs, int &umin, int &usec)
void terminateAnalysis ()
virtual void restartYourself (TimeStep *tS)
virtual void doStepOutput (TimeStep *tStep)
void saveStepContext (TimeStep *tStep, ContextMode mode)
virtual void initializeYourself (TimeStep *tStep)
virtual int initializeAdaptive (int tStepNumber)
virtual int giveNumberOfDomainEquations (int di, const UnknownNumberingScheme &num)
virtual FieldPtr giveField (FieldType key, TimeStep *)
virtual FieldPtr giveField (InternalStateType key, TimeStep *)
EngngModelgiveMasterEngngModel ()
 Returns the master engnmodel.
virtual double giveLoadLevel ()
 Returns the current load level.
virtual double giveEigenValue (int eigNum)
 Only relevant for eigen value analysis. Otherwise returns zero.
virtual void setActiveVector (int i)
 Only relevant for eigen value analysis. Otherwise does noting.
int updateSharedDofManagers (FloatArray &answer, const UnknownNumberingScheme &s, int ExchangeTag)
int exchangeRemoteElementData (int ExchangeTag)
virtual int giveCurrentNumberOfIterations ()
MPI_Comm giveParallelComm ()
 Returns the communication object of reciever.
int packRemoteElementData (ProcessCommunicator &processComm)
int unpackRemoteElementData (ProcessCommunicator &processComm)
int packDofManagers (ArrayWithNumbering *src, ProcessCommunicator &processComm)
int unpackDofManagers (ArrayWithNumbering *dest, ProcessCommunicator &processComm)
ProblemCommunicatorgiveProblemCommunicator (EngngModelCommType t)
void initializeCommMaps (bool forceInit=false)
virtual int instanciateYourself (DataReader &dr, InputRecord &ir, const char *outFileName, const char *desc)
void Instanciate_init ()
int instanciateDomains (DataReader &dr)
 Instanciate problem domains by calling their instanciateYourself() service.
int instanciateMetaSteps (DataReader &dr)
 Instanciate problem meta steps by calling their instanciateYourself() service.
virtual int instanciateDefaultMetaStep (InputRecord &ir)
 Instanciate default metastep, if nmsteps is zero.
virtual void updateAttributes (MetaStep *mStep)
void initMetaStepAttributes (MetaStep *mStep)
virtual void updateDomainLinks ()
MetaStepgiveCurrentMetaStep ()
 Returns current meta step.
virtual TimeStepgiveCurrentStep (bool force=false)
virtual void adaptTimeStep (double nIter)
virtual TimeStepgivePreviousStep (bool force=false)
virtual void preInitializeNextStep ()
 Does a pre-initialization of the next time step (implement if necessarry).
virtual TimeStepgiveSolutionStepWhenIcApply (bool force=false)
int giveNumberOfMetaSteps ()
 Return number of meta steps.
MetaStepgiveMetaStep (int i)
 Returns the i-th meta step.
int giveNumberOfSteps (bool force=false)
virtual double giveEndOfTimeOfInterest ()
 Returns end of time interest (time corresponding to end of time integration).
int giveNumberOfTimeStepWhenIcApply ()
 Returns the time step number, when initial conditions should apply.
ExportModuleManagergiveExportModuleManager ()
 Returns receiver's export module manager.
EngngModelTimergiveTimer ()
 Returns reference to receiver timer (EngngModelTimer).
virtual double giveFinalTime ()
virtual int giveNewEquationNumber (int domain, DofIDItem)
virtual int giveNewPrescribedEquationNumber (int domain, DofIDItem)
std::string giveContextFileName (int tStepNumber, int stepVersion) const
std::string giveDomainFileName (int domainNum, int domainSerNum) const
virtual void updateComponent (TimeStep *tStep, NumericalCmpn cmpn, Domain *d)
virtual void initForNewIteration (Domain *d, TimeStep *tStep, int iterationNumber, const FloatArray &solution)
virtual void updateSolution (FloatArray &solutionVector, TimeStep *tStep, Domain *d)
virtual void updateMatrix (SparseMtrx &mat, TimeStep *tStep, Domain *d)
virtual void initStepIncrements ()
virtual int forceEquationNumbering (int i)
virtual int forceEquationNumbering ()
virtual int requiresUnknownsDictionaryUpdate ()
virtual bool requiresEquationRenumbering (TimeStep *tStep)
virtual void updateDofUnknownsDictionary (DofManager *, TimeStep *)
virtual int giveUnknownDictHashIndx (ValueModeType mode, TimeStep *tStep)
virtual bool newDofHandling ()
virtual ParallelContextgiveParallelContext (int n)
virtual void initParallelContexts ()
virtual void assemble (SparseMtrx &answer, TimeStep *tStep, const MatrixAssembler &ma, const UnknownNumberingScheme &s, Domain *domain)
virtual void assemble (SparseMtrx &answer, TimeStep *tStep, const MatrixAssembler &ma, const UnknownNumberingScheme &r_s, const UnknownNumberingScheme &c_s, Domain *domain)
void assembleVector (FloatArray &answer, TimeStep *tStep, const VectorAssembler &va, ValueModeType mode, const UnknownNumberingScheme &s, Domain *domain, FloatArray *eNorms=NULL)
void assembleVectorFromDofManagers (FloatArray &answer, TimeStep *tStep, const VectorAssembler &va, ValueModeType mode, const UnknownNumberingScheme &s, Domain *domain, FloatArray *eNorms=NULL)
void assembleVectorFromElements (FloatArray &answer, TimeStep *tStep, const VectorAssembler &va, ValueModeType mode, const UnknownNumberingScheme &s, Domain *domain, FloatArray *eNorms=NULL)
void assembleVectorFromBC (FloatArray &answer, TimeStep *tStep, const VectorAssembler &va, ValueModeType mode, const UnknownNumberingScheme &s, Domain *domain, FloatArray *eNorms=NULL)
void assembleExtrapolatedForces (FloatArray &answer, TimeStep *tStep, CharType type, Domain *domain)
void assemblePrescribedExtrapolatedForces (FloatArray &answer, TimeStep *tStep, CharType type, Domain *domain)
virtual int checkProblemConsistency ()
virtual void init ()
virtual void postInitialize ()
virtual void printOutputAt (FILE *file, TimeStep *tStep, const IntArray &nodeSets, const IntArray &elementSets)
void outputNodes (FILE *file, Domain &domain, TimeStep *tStep, int setNum)
void outputElements (FILE *file, Domain &domain, TimeStep *tStep, int setNum)
void printYourself ()
 Prints state of receiver. Useful for debugging.
virtual int useNonlocalStiffnessOption ()
 Returns nonzero if nonlocal stiffness option activated.
bool isParallel () const
 Returns true if receiver in parallel mode.
int giveRank () const
 Returns domain rank in a group of collaborating processes (0..groupSize-1).
int giveNumberOfProcesses () const
 Returns the number of collaborating processes.
EngngModelContextgiveContext ()
 Context requesting service.
virtual int giveNumberOfSlaveProblems ()
 Returns number of slave problems.
virtual EngngModelgiveSlaveProblem (int i)
 Returns i-th slave problem.
virtual bool giveEquationScalingFlag ()
 Returns the Equation scaling flag, which is used to indicate that governing equation(s) are scaled, or non-dimensionalized.
virtual double giveVariableScale (VarScaleType varId)
 Returns the scale factor for given variable type.
virtual void balanceLoad (TimeStep *tStep)
virtual LoadBalancergiveLoadBalancer ()
virtual LoadBalancerMonitorgiveLoadBalancerMonitor ()
void initParallel ()
 Request domain rank and problem size.
EngngModelgiveEngngModel ()
 Returns reference to itself -> required by communicator.h.
virtual bool isElementActivated (int elemNum)
virtual bool isElementActivated (Element *e)
TimeStepControllergiveTimeStepController ()
 Returns the time step controller.
virtual void drawYourself (oofegGraphicContext &gc)
virtual void drawElements (oofegGraphicContext &gc)
virtual void drawNodes (oofegGraphicContext &gc)
std::string errorInfo (const char *func) const
 Returns string for prepending output (used by error reporting macros).

Protected Member Functions

void computeLoadVector (FloatArray &answer, ValueModeType mode, TimeStep *tStep)
void computeMassMtrx (FloatArray &mass, double &maxOm, TimeStep *tStep)
void computeMassMtrx2 (FloatMatrix &mass, double &maxOm, TimeStep *tStep)
double giveInitialTime () override
 return time at the begining of analysis
Protected Member Functions inherited from oofem::StructuralEngngModel
void printReactionForces (TimeStep *tStep, int id, FILE *out)
virtual void computeExternalLoadReactionContribution (FloatArray &reactions, TimeStep *tStep, int di)
void updateInternalState (TimeStep *tStep)
void printOutputAt (FILE *file, TimeStep *tStep) override
virtual void packMigratingData (TimeStep *tStep)
virtual void unpackMigratingData (TimeStep *tStep)

Protected Attributes

FloatArray massMatrix
 Mass matrix.
FloatArray loadVector
 Load vector.
FloatArray previousIncrementOfDisplacementVector
 Vector storing displacement increments.
FloatArray displacementVector
 Displacement, velocity and acceleration vectors.
FloatArray velocityVector
FloatArray accelerationVector
FloatArray internalForces
 Vector of real nodal forces.
double dumpingCoef
 Dumping coefficient (C = dumpingCoef * MassMtrx).
double deltaT
 Time step.
int initFlag
 Flag indicating the need for initialization.
double reductionFactor
 Optional reduction factor for time step deltaT.
int drFlag
 Flag indicating whether dynamic relaxation takes place.
FloatArray loadRefVector
 Reference load vector.
double c
 Parameter determining rate of the loading process.
double pt
 Load level.
double Tau
 End of time interval.
double pyEstimate
 Estimate of loadRefVector^T*displacementVector(Tau).
double pMp
 Product of p^tM^(-1)p; where p is reference load vector.
LinSystSolverType solverType
SparseMtrxType sparseMtrxType
std::unique_ptr< SparseLinearSystemNMnMethod
Protected Attributes inherited from oofem::StructuralEngngModel
StateCounterType internalVarUpdateStamp
FloatArray internalForcesEBENorm
 Norm of nodal internal forces evaluated on element by element basis (squared).
Protected Attributes inherited from oofem::EngngModel
int ndomains
 Number of receiver domains.
std ::vector< std ::unique_ptr< Domain > > domainList
 List of problem domains.
int numberOfSteps
 Total number of time steps.
int numberOfEquations
 Total number of equation in current time step.
int numberOfPrescribedEquations
 Total number or prescribed equations in current time step.
IntArray domainNeqs
 Number of equations per domain.
IntArray domainPrescribedNeqs
 Number of prescribed equations per domain.
bool renumberFlag
 Renumbering flag (renumbers equations after each step, necessary if Dirichlet BCs change).
bool profileOpt
 Profile optimized numbering flag (using Sloan's algorithm).
int equationNumberingCompleted
 Equation numbering completed flag.
int nMetaSteps
 Number of meta steps.
std ::vector< MetaStepmetaStepList
 List of problem metasteps.
std ::unique_ptr< TimeStepstepWhenIcApply
 Solution step when IC (initial conditions) apply.
std ::unique_ptr< TimeStepcurrentStep
 Current time step.
std ::unique_ptr< TimeSteppreviousStep
 Previous time step.
int number
 Receivers id.
std::string dataOutputFileName
 Path to output stream.
std::string coreOutputFileName
 String with core output file name.
FILE * outputStream
 Output stream.
std::string referenceFileName
 String with reference file name.
ContextOutputMode contextOutputMode
 Domain context output mode.
int contextOutputStep
ExportModuleManager exportModuleManager
 Export module manager.
InitModuleManager initModuleManager
 Initialization module manager.
MonitorManager monitorManager
 Monitor manager.
problemMode pMode
 Domain mode.
problemScale pScale
 Multiscale mode.
time_t startTime
 Solution start time.
EngngModelmaster
 Master e-model; if defined receiver is in maintained (slave) mode.
EngngModelContextcontext
 Context.
EngngModelTimer timer
 E-model timer.
int parallelFlag
 Flag indicating that the receiver runs in parallel.
enum fMode nonLinFormulation
 Type of non linear formulation (total or updated formulation).
std::unique_ptr< ErrorEstimatordefaultErrEstimator
 Error estimator. Useful for adaptivity, or simply printing errors output.
std::unique_ptr< TimeStepControllertimeStepController
 Time Step controller is responsible for collecting data from analysis, elements, and materials, and select the appropriate timestep size for the next step, or reduce the step in case of convergence problems.
int rank
 Domain rank in a group of collaborating processes (0..groupSize-1).
int numProcs
 Total number of collaborating processes.
int nonlocalExt
 Flag indicating if nonlocal extension active, which will cause data to be sent between shared elements before computing the internal forces.
char processor_name [PROCESSOR_NAME_LENGTH]
 Processor name.
MPI_Comm comm
 Communication object for this engineering model.
std::unique_ptr< LoadBalancerlb
 Load Balancer.
std::unique_ptr< LoadBalancerMonitorlbm
bool loadBalancingFlag
 If set to true, load balancing is active.
bool force_load_rebalance_in_first_step
 Debug flag forcing load balancing after first step.
CommunicatorBuffcommBuff
 Common Communicator buffer.
ProblemCommunicatorcommunicator
 Communicator.
ProblemCommunicatornonlocCommunicator
 NonLocal Communicator. Necessary when nonlocal constitutive models are used.
std ::vector< ParallelContextparallelContextList
 List where parallel contexts are stored.
bool suppressOutput
 Flag for suppressing output to file.
std::string simulationDescription

Additional Inherited Members

Public Types inherited from oofem::EngngModel
enum  EngngModel_UpdateMode { EngngModel_Unknown_Mode , EngngModel_SUMM_Mode , EngngModel_SET_Mode }
enum  EngngModelCommType { PC_default , PC_nonlocal }
enum  InitialGuess { IG_None = 0 , IG_Tangent = 1 }
Protected Types inherited from oofem::EngngModel
enum  {
  InternalForcesExchangeTag , MassExchangeTag , LoadExchangeTag , ReactionExchangeTag ,
  RemoteElementExchangeTag
}
 Message tags. More...

Detailed Description

This class implements NonLinear (- may be changed) solution of dynamic problems using Direct Explicit Integration scheme - Central Difference Method. For efficiency reasons it uses diagonal mass matrix. It is formulated in increments of displacements rather than in total variables.

Solution of this problem is series of loading cases, maintained as sequence of time-steps.

Analysis starts assembling governing equations at time step 0 ( 0 given by boundary and initial cond.) they result in response at time step 1. For time step 0 we need special start code. Because this method is explicit, when solving equations for step t, we obtain solution in step t+dt. But printing is performed for step t. So, when analyst specifies initial conditions, then he/she specifies them in time step 0.

Current implementation supports parallel processing. Both node- and element cut strategies can be used.

  • In node cut strategy, partitions are divided using cut, which goes through nodes. These cut nodes are called "shared" nodes. Generally, unknown values in shared nodes are composed from local partition contributions as well as from contributions from remote partitions sharing this node. Particularly, masses and real nodal forces have to be exchanged for shared nodes.
  • In element cut strategy, partitions are divided using cut running through elements. The cut elements are replicated on neighbouring partitions. The nodes belonging to replicated elements belonging to remote partitions are called remote nodes. The are mirrors or remote copies of corresponding nodes on neighbouring partition.
  • Additional mode has been introduced remote element mode. It introduces the "remote" elements, the exact local mirrors of remote counterparts. Introduced to support general nonlocal constitutive models, in order to provide efficient way, how to average local data without need of fine grain communication.

Definition at line 91 of file nldeidynamic.h.

Constructor & Destructor Documentation

◆ NlDEIDynamic()

◆ ~NlDEIDynamic()

oofem::NlDEIDynamic::~NlDEIDynamic ( )
virtual

Definition at line 68 of file nldeidynamic.C.

Member Function Documentation

◆ computeLoadVector()

void oofem::NlDEIDynamic::computeLoadVector ( FloatArray & answer,
ValueModeType mode,
TimeStep * tStep )
protected

Assembles the load vector. If in parallel mode, the loads of shared/remote nodes are exchanged and remote contributions are taken into account.

Parameters
answerLoad vector.
modeValue type mode of load vector.
tStepSolution step.

Definition at line 473 of file nldeidynamic.C.

References oofem::EngngModel::assembleVector(), oofem::EngngModel::giveDomain(), oofem::EngngModel::giveNumberOfDomainEquations(), oofem::EngngModel::LoadExchangeTag, oofem::FloatArray::resize(), oofem::EngngModel::updateSharedDofManagers(), and oofem::FloatArray::zero().

Referenced by solveYourselfAt().

◆ computeMassMtrx()

void oofem::NlDEIDynamic::computeMassMtrx ( FloatArray & mass,
double & maxOm,
TimeStep * tStep )
protected

Assembles the diagonal mass matrix of receiver. Local or Global variant of zero mass elements replacement is performed. If runs in parallel, the masses of shared nodes are exchanged and remote contributions are added accordingly.

Parameters
massAssembled mass matrix.
maxOmEstimate of eigenfrequency.
tStepTime step.
Todo
This rotation matrix is not flexible enough.. it can only work with full size matrices and doesn't allow for flexibility in the matrixassembler.
Todo
This rotation matrix is not flexible enough.. it can only work with full size matrices and doesn't allow for flexibility in the matrixassembler.

Definition at line 492 of file nldeidynamic.C.

References oofem::FloatArray::at(), oofem::FloatMatrix::at(), oofem::IntArray::at(), oofem::Element_remote, oofem::Element::giveCharacteristicMatrix(), oofem::EngngModel::giveDomain(), oofem::Domain::giveElement(), oofem::Domain::giveElements(), oofem::Element::giveLocationArray(), oofem::EngngModel::giveNumberOfDomainEquations(), oofem::Domain::giveNumberOfElements(), oofem::FloatMatrix::giveNumberOfRows(), oofem::Element::giveParallelMode(), oofem::EngngModel::giveRank(), oofem::Element::giveRotationMatrix(), oofem::IntArray::giveSize(), oofem::FloatMatrix::isNotEmpty(), oofem::EngngModel::MassExchangeTag, massMatrix, oofem::max(), OOFEM_ERROR, OOFEM_WARNING, oofem::FloatMatrix::rotatedWith(), oofem::EngngModel::updateSharedDofManagers(), VERBOSEPARALLEL_PRINT, and ZERO_REL_MASS.

Referenced by solveYourselfAt().

◆ computeMassMtrx2()

void oofem::NlDEIDynamic::computeMassMtrx2 ( FloatMatrix & mass,
double & maxOm,
TimeStep * tStep )
protected

◆ estimateMaxPackSize()

int oofem::NlDEIDynamic::estimateMaxPackSize ( IntArray & commMap,
DataStream & buff,
int packUnpackType )
overridevirtual

Determines the space necessary for send/receive buffer. It uses related communication map pattern to determine the maximum size needed.

Parameters
commMapCommunication map used to send/receive messages.
buffCommunication buffer.
packUnpackTypeDetermines the type of packed quantity, used by receiver to estimate the size of pack/unpack buffer accordingly.
Returns
Upper bound of space needed.
Todo
Fix this old ProblemCommMode__NODE_CUT value

Reimplemented from oofem::EngngModel.

Definition at line 665 of file nldeidynamic.C.

References oofem::Element::estimatePackSize(), oofem::Domain::giveDofManager(), oofem::EngngModel::giveDomain(), oofem::Domain::giveElement(), oofem::DataStream::givePackSizeOfDouble(), and oofem::max().

◆ giveClassName()

const char * oofem::NlDEIDynamic::giveClassName ( ) const
inlineoverridevirtual

Returns class name of the receiver.

Implements oofem::EngngModel.

Definition at line 156 of file nldeidynamic.h.

◆ giveFormulation()

fMode oofem::NlDEIDynamic::giveFormulation ( )
inlineoverridevirtual

Indicates type of non linear computation (total or updated formulation). This is used for example on Nodal level to update coordinates if updated formulation is done, or on element level, when non linear contributions are computed.

Reimplemented from oofem::EngngModel.

Definition at line 157 of file nldeidynamic.h.

References oofem::TL.

◆ giveInitialTime()

double oofem::NlDEIDynamic::giveInitialTime ( )
inlineoverrideprotectedvirtual

return time at the begining of analysis

Reimplemented from oofem::EngngModel.

Definition at line 182 of file nldeidynamic.h.

◆ giveInputRecordName()

const char * oofem::NlDEIDynamic::giveInputRecordName ( ) const
inline

Definition at line 155 of file nldeidynamic.h.

References _IFT_NlDEIDynamic_Name.

◆ giveNextStep()

TimeStep * oofem::NlDEIDynamic::giveNextStep ( )
overridevirtual

Returns next time step (next to current step) of receiver.

Reimplemented from oofem::EngngModel.

Definition at line 150 of file nldeidynamic.C.

References oofem::EngngModel::currentStep, deltaT, and oofem::EngngModel::previousStep.

◆ giveNumberOfFirstStep()

int oofem::NlDEIDynamic::giveNumberOfFirstStep ( bool force = false)
inlineoverridevirtual

Returns number of first time step used by receiver.

Parameters
forcewhen set to true then receiver reply is returned instead of master (default)

Reimplemented from oofem::EngngModel.

Definition at line 159 of file nldeidynamic.h.

◆ giveNumericalMethod()

NumericalMethod * oofem::NlDEIDynamic::giveNumericalMethod ( MetaStep * mStep)
overridevirtual

Returns reference to receiver's numerical method.

Reimplemented from oofem::EngngModel.

Definition at line 72 of file nldeidynamic.C.

References oofem::classFactory, oofem::EngngModel::giveDomain(), nMethod, and solverType.

◆ giveUnknownComponent()

double oofem::NlDEIDynamic::giveUnknownComponent ( ValueModeType ,
TimeStep * ,
Domain * ,
Dof *  )
overridevirtual

Returns requested unknown. Unknown at give time step is characterized by its type and mode and by its equation number. This function is used by Dofs, when they are requested for their associated unknowns.

See also
Dof::giveUnknown

Reimplemented from oofem::EngngModel.

Definition at line 114 of file nldeidynamic.C.

References oofem::Dof::__giveEquationNumber(), accelerationVector, displacementVector, oofem::EngngModel::giveCurrentStep(), OOFEM_ERROR, previousIncrementOfDisplacementVector, and velocityVector.

◆ initializeFrom()

void oofem::NlDEIDynamic::initializeFrom ( InputRecord & ir)
overridevirtual

Initializes receiver according to object description in input reader. InitString can be imagined as data record in component database belonging to receiver. Receiver may use value-name extracting functions to extract particular field from record.

Reimplemented from oofem::EngngModel.

Definition at line 82 of file nldeidynamic.C.

References _IFT_NlDEIDynamic_deltat, _IFT_NlDEIDynamic_drflag, _IFT_NlDEIDynamic_dumpcoef, _IFT_NlDEIDynamic_nonlocalext, _IFT_NlDEIDynamic_py, _IFT_NlDEIDynamic_reduct, _IFT_NlDEIDynamic_tau, oofem::EngngModel::commBuff, oofem::EngngModel::communicator, deltaT, drFlag, dumpingCoef, oofem::EngngModel::giveNumberOfProcesses(), oofem::EngngModel::giveRank(), oofem::InputRecord::hasField(), IR_GIVE_FIELD, IR_GIVE_OPTIONAL_FIELD, oofem::EngngModel::nonlocalExt, oofem::EngngModel::nonlocCommunicator, pyEstimate, reductionFactor, and Tau.

◆ printDofOutputAt()

void oofem::NlDEIDynamic::printDofOutputAt ( FILE * stream,
Dof * iDof,
TimeStep * tStep )
overridevirtual

DOF printing routine. Called by DofManagers to print Dof specific part. Dof class provides component printing routines, but emodel is responsible for what will be printed at DOF level.

Parameters
streamoutput stream
iDofdof to be processed
tStepsolution step

Reimplemented from oofem::EngngModel.

Definition at line 751 of file nldeidynamic.C.

References oofem::Dof::printMultipleOutputAt().

◆ printOutputAt()

void oofem::NlDEIDynamic::printOutputAt ( FILE * file,
TimeStep * tStep )
overridevirtual

Prints output of receiver to output domain stream, for given time step. Corresponding function for element gauss points is invoked (gaussPoint::printOutputAt).

Reimplemented from oofem::EngngModel.

Definition at line 763 of file nldeidynamic.C.

References drFlag, oofem::EngngModel::giveDomain(), oofem::TimeStep::giveNumber(), oofem::TimeStep::giveTargetTime(), oofem::StructuralEngngModel::printReactionForces(), and pt.

◆ restoreContext()

void oofem::NlDEIDynamic::restoreContext ( DataStream & stream,
ContextMode mode )
overridevirtual

Restores the state of model from output stream. Restores not only the receiver state, but also same function is invoked for all DofManagers and Elements in associated domain. Note that by restoring element context also contexts of all associated integration points (and material statuses) are restored. Each context is associated with unique time step. Only one context per time step is allowed. Restore context function will restore such context, which is related (through its step number) to time step number and version given in obj parameter. Restoring context will change current time step in order to correspond to newly restored context.

Parameters
streamContext file.
modeDetermines amount of info in stream.
Exceptions
ContextIOERRexception if error encountered.

Reimplemented from oofem::EngngModel.

Definition at line 722 of file nldeidynamic.C.

References accelerationVector, oofem::CIO_IOERR, oofem::CIO_OK, deltaT, displacementVector, previousIncrementOfDisplacementVector, oofem::DataStream::read(), THROW_CIOERR, and velocityVector.

◆ saveContext()

void oofem::NlDEIDynamic::saveContext ( DataStream & stream,
ContextMode mode )
overridevirtual

Stores the state of model to output stream. Stores not only the receiver state, but also same function is invoked for all DofManagers and Elements in associated domain. Note that by storing element context also contexts of all associated integration points (and material statuses) are stored.

Parameters
streamContext stream.
modeDetermines amount of info in stream.
Exceptions
ContextIOERRIf error encountered.

Reimplemented from oofem::EngngModel.

Definition at line 694 of file nldeidynamic.C.

References accelerationVector, oofem::CIO_IOERR, oofem::CIO_OK, deltaT, displacementVector, previousIncrementOfDisplacementVector, THROW_CIOERR, velocityVector, and oofem::DataStream::write().

◆ solveYourself()

void oofem::NlDEIDynamic::solveYourself ( )
overridevirtual

Starts solution process. Implementation should invoke for each time step solveYourselfAt function with time step as parameter. Time steps are created using giveNextStep function (this will set current time step to newly created, and updates previous step).

Reimplemented from oofem::EngngModel.

Definition at line 169 of file nldeidynamic.C.

References oofem::EngngModel::giveNumberOfDomainEquations(), oofem::EngngModel::giveRank(), oofem::EngngModel::initializeCommMaps(), oofem::EngngModel::isParallel(), and OOFEM_LOG_INFO.

◆ solveYourselfAt()

void oofem::NlDEIDynamic::solveYourselfAt ( TimeStep * tStep)
overridevirtual

◆ updateYourself()

void oofem::NlDEIDynamic::updateYourself ( TimeStep * tStep)
overridevirtual

Updates internal state after finishing time step. (for example total values may be updated according to previously solved increments). Then element values are also updated (together with related integration points and material statuses).

Reimplemented from oofem::EngngModel.

Definition at line 461 of file nldeidynamic.C.

Member Data Documentation

◆ accelerationVector

FloatArray oofem::NlDEIDynamic::accelerationVector
protected

◆ c

double oofem::NlDEIDynamic::c
protected

Parameter determining rate of the loading process.

Definition at line 118 of file nldeidynamic.h.

Referenced by solveYourselfAt().

◆ deltaT

double oofem::NlDEIDynamic::deltaT
protected

Time step.

Definition at line 107 of file nldeidynamic.h.

Referenced by giveNextStep(), initializeFrom(), restoreContext(), saveContext(), and solveYourselfAt().

◆ displacementVector

FloatArray oofem::NlDEIDynamic::displacementVector
protected

Displacement, velocity and acceleration vectors.

Definition at line 101 of file nldeidynamic.h.

Referenced by giveUnknownComponent(), NlDEIDynamic(), restoreContext(), saveContext(), and solveYourselfAt().

◆ drFlag

int oofem::NlDEIDynamic::drFlag
protected

Flag indicating whether dynamic relaxation takes place.

Definition at line 114 of file nldeidynamic.h.

Referenced by initializeFrom(), printOutputAt(), and solveYourselfAt().

◆ dumpingCoef

double oofem::NlDEIDynamic::dumpingCoef
protected

Dumping coefficient (C = dumpingCoef * MassMtrx).

Definition at line 105 of file nldeidynamic.h.

Referenced by initializeFrom(), and solveYourselfAt().

◆ initFlag

int oofem::NlDEIDynamic::initFlag
protected

Flag indicating the need for initialization.

Definition at line 109 of file nldeidynamic.h.

Referenced by NlDEIDynamic(), and solveYourselfAt().

◆ internalForces

FloatArray oofem::NlDEIDynamic::internalForces
protected

Vector of real nodal forces.

Definition at line 103 of file nldeidynamic.h.

Referenced by NlDEIDynamic(), and solveYourselfAt().

◆ loadRefVector

FloatArray oofem::NlDEIDynamic::loadRefVector
protected

Reference load vector.

Definition at line 116 of file nldeidynamic.h.

Referenced by solveYourselfAt().

◆ loadVector

FloatArray oofem::NlDEIDynamic::loadVector
protected

Load vector.

Definition at line 97 of file nldeidynamic.h.

Referenced by NlDEIDynamic(), and solveYourselfAt().

◆ massMatrix

FloatArray oofem::NlDEIDynamic::massMatrix
protected

Mass matrix.

Definition at line 95 of file nldeidynamic.h.

Referenced by computeMassMtrx(), NlDEIDynamic(), and solveYourselfAt().

◆ nMethod

std::unique_ptr<SparseLinearSystemNM> oofem::NlDEIDynamic::nMethod
protected

Definition at line 130 of file nldeidynamic.h.

Referenced by giveNumericalMethod().

◆ pMp

double oofem::NlDEIDynamic::pMp
protected

Product of p^tM^(-1)p; where p is reference load vector.

Definition at line 126 of file nldeidynamic.h.

Referenced by solveYourselfAt().

◆ previousIncrementOfDisplacementVector

FloatArray oofem::NlDEIDynamic::previousIncrementOfDisplacementVector
protected

Vector storing displacement increments.

Definition at line 99 of file nldeidynamic.h.

Referenced by giveUnknownComponent(), NlDEIDynamic(), restoreContext(), saveContext(), and solveYourselfAt().

◆ pt

double oofem::NlDEIDynamic::pt
protected

Load level.

Definition at line 120 of file nldeidynamic.h.

Referenced by printOutputAt(), and solveYourselfAt().

◆ pyEstimate

double oofem::NlDEIDynamic::pyEstimate
protected

Estimate of loadRefVector^T*displacementVector(Tau).

Definition at line 124 of file nldeidynamic.h.

Referenced by initializeFrom(), and solveYourselfAt().

◆ reductionFactor

double oofem::NlDEIDynamic::reductionFactor
protected

Optional reduction factor for time step deltaT.

Definition at line 111 of file nldeidynamic.h.

Referenced by initializeFrom(), and solveYourselfAt().

◆ solverType

LinSystSolverType oofem::NlDEIDynamic::solverType
protected

Definition at line 128 of file nldeidynamic.h.

Referenced by giveNumericalMethod().

◆ sparseMtrxType

SparseMtrxType oofem::NlDEIDynamic::sparseMtrxType
protected

Definition at line 129 of file nldeidynamic.h.

◆ Tau

double oofem::NlDEIDynamic::Tau
protected

End of time interval.

Definition at line 122 of file nldeidynamic.h.

Referenced by initializeFrom(), and solveYourselfAt().

◆ velocityVector

FloatArray oofem::NlDEIDynamic::velocityVector
protected

The documentation for this class was generated from the following files:

This page is part of the OOFEM-3.0 documentation. Copyright Copyright (C) 1994-2025 Borek Patzak Bořek Patzák
Project e-mail: oofem@fsv.cvut.cz
Generated at for OOFEM by doxygen 1.15.0 written by Dimitri van Heesch, © 1997-2011