61DIIDynamic :: ~DIIDynamic()
85 StructuralEngngModel :: initializeFrom(ir);
120 if (
theta < 1.37 ) {
121 OOFEM_WARNING(
"Found theta < 1.37. Performing correction, theta = 1.37");
132double DIIDynamic :: giveUnknownComponent(ValueModeType mode,
TimeStep *tStep,
Domain *d,
Dof *dof)
154 case VM_Acceleration:
158 OOFEM_ERROR(
"Unknown is of undefined ValueModeType for this problem");
168 double totalTime =
deltaT;
175 counter =
currentStep->giveSolutionStateCounter() + 1;
185 currentStep = std::make_unique<TimeStep>(istep,
this, 1, totalTime,
deltaT, counter, td);
190void DIIDynamic :: initializeYourself(
TimeStep *tStep)
213 for (
Dof *iDof: *node ) {
217 if ( !iDof->isPrimaryDof() ) {
221 int jj = iDof->__giveEquationNumber();
234void DIIDynamic :: solveYourself()
236 StructuralEngngModel :: solveYourself();
302 for (
int i = 1; i <= neq; i++ ) {
319 for (
int i = 1; i <= neq; i++ ) {
330 for (
int i = 1; i <= neq; i++ ) {
336 for (
int i = 1; i <= neq; i++ ) {
341 for (
int i = 1; i <= neq; i++ ) {
356 OOFEM_LOG_INFO(
"TD_Wilson: Updating acceleration, velocity and displacement.\n");
357 for (
int i = 1; i <= neq; i++ ) {
371 for (
int i = 1; i <= neq; i++ ) {
379 for (
int i = 1; i <= neq; i++ ) {
400 StructuralEngngModel :: updateYourself(tStep);
405DIIDynamic :: printDofOutputAt(FILE *stream,
Dof *iDof,
TimeStep *tStep)
407 static char dofchar[] =
"dva";
408 static ValueModeType dofmodes[] = {
409 VM_Total, VM_Velocity, VM_Acceleration
421 int i, j, k, jj, kk, n;
430 for ( i = 1; i <= nelem; i++ ) {
454 for ( j = 1; j <= n; j++ ) {
457 for ( k = 1; k <= n; k++ ) {
460 answer.
at(jj) += charMtrx.
at(j, k) * vec.
at(kk);
494 for (
int ielem = 1; ielem <= nelem; ielem++ ) {
560 a3 = 1 / ( 2 *
beta ) - 1;
565 a8 = dt2 * ( ( 1 / 2 ) -
beta );
620 EngngModel :: saveContext(stream, mode);
647 EngngModel :: restoreContext(stream, mode);
#define REGISTER_EngngModel(class)
void timesMtrx(FloatArray &answer, FloatArray &vec, CharType type, Domain *domain, TimeStep *tStep)
FloatArray previousDisplacementVector
void determineConstants(TimeStep *tStep)
void assembleDirichletBcRhsVector(FloatArray &answer, Domain *d, TimeStep *tStep)
std ::unique_ptr< SparseLinearSystemNM > nMethod
Numerical method used to solve the problem.
FloatArray previousVelocityVector
LinSystSolverType solverType
void assembleLoadVector(FloatArray &_loadVector, Domain *domain, ValueModeType mode, TimeStep *tStep)
TimeDiscretizationType initialTimeDiscretization
FloatArray velocityVector
FloatArray accelerationVector
FloatArray previousAccelerationVector
FloatArray previousIncrementOfDisplacement
FloatArray previousLoadVector
FloatArray displacementVector
std ::unique_ptr< SparseMtrx > stiffnessMatrix
SparseMtrxType sparseMtrxType
virtual void printMultipleOutputAt(FILE *File, TimeStep *tStep, char *ch, ValueModeType *mode, int nite)
virtual int __giveEquationNumber() const =0
int giveNumber()
Returns domain number.
int giveNumberOfElements() const
Returns number of elements in domain.
std ::vector< std ::unique_ptr< DofManager > > & giveDofManagers()
Element * giveElement(int n)
virtual void giveElementDofIDMask(IntArray &answer) const
virtual bool giveRotationMatrix(FloatMatrix &answer)
virtual void giveCharacteristicMatrix(FloatMatrix &answer, CharType type, TimeStep *tStep)
void computeVectorOfPrescribed(ValueModeType u, TimeStep *tStep, FloatArray &answer)
void giveLocationArray(IntArray &locationArray, const UnknownNumberingScheme &s, IntArray *dofIds=NULL) const
elementParallelMode giveParallelMode() const
int giveNumberOfTimeStepWhenIcApply()
Returns the time step number, when initial conditions should apply.
virtual TimeStep * giveCurrentStep(bool force=false)
virtual int giveNumberOfDomainEquations(int di, const UnknownNumberingScheme &num)
std ::unique_ptr< TimeStep > previousStep
Previous time step.
int ndomains
Number of receiver domains.
Domain * giveDomain(int n)
std ::unique_ptr< TimeStep > currentStep
Current time step.
virtual int giveNumberOfFirstStep(bool force=false)
int updateSharedDofManagers(FloatArray &answer, const UnknownNumberingScheme &s, int ExchangeTag)
std ::unique_ptr< TimeStep > stepWhenIcApply
Solution step when IC (initial conditions) apply.
void assembleVector(FloatArray &answer, TimeStep *tStep, const VectorAssembler &va, ValueModeType mode, const UnknownNumberingScheme &s, Domain *domain, FloatArray *eNorms=NULL)
bool containsOnlyZeroes() const
void assemble(const FloatArray &fe, const IntArray &loc)
void zero()
Zeroes all coefficients of receiver.
void rotatedWith(FloatMatrix &r, char mode)
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
void rotatedWith(const FloatMatrix &r, char mode='n')
void add(const FloatMatrix &a)
bool isNotEmpty() const
Tests for empty matrix.
int giveNumberOfRows() const
Returns number of rows of receiver.
double at(std::size_t i, std::size_t j) const
StructuralEngngModel(int i, EngngModel *master=nullptr)
Creates new StructuralEngngModel with number i, associated to domain d.
double giveTimeIncrement()
Returns solution step associated time increment.
double giveTargetTime()
Returns target time.
TimeDiscretizationType giveTimeDiscretization()
Returns time discretization.
int giveNumber()
Returns receiver's number.
TimeStep * givePreviousStep()
Returns pointer to previous solution step.
#define _IFT_DIIDynamic_eta
#define _IFT_DIIDynamic_theta
#define _IFT_DIIDynamic_gamma
#define _IFT_DIIDynamic_deltat
#define _IFT_DIIDynamic_delta
#define _IFT_DIIDynamic_ddtScheme
#define _IFT_DIIDynamic_beta
#define _IFT_EngngModel_lstype
#define _IFT_EngngModel_smtype
#define OOFEM_WARNING(...)
#define OOFEM_LOG_INFO(...)
#define OOFEM_LOG_RELEVANT(...)
#define OOFEM_LOG_DEBUG(...)
@ Element_remote
Element in active domain is only mirror of some remote element.
FloatArrayF< N > assemble(const FloatArrayF< M > &x, int const (&c)[M])
Assemble components into zero matrix.
long StateCounterType
StateCounterType type used to indicate solution state.
ClassFactory & classFactory
TimeDiscretizationType
Time discretization used by transient solvers.
@ TD_TwoPointBackward
Two-point Backward Euler method.
@ TD_Newmark
Newmark-beta method.
@ TD_ThreePointBackward
Three-point Backward Euler method.
@ TD_Wilson
Wilson-theta method.