OOFEM 3.0
Loading...
Searching...
No Matches
deidynamic.C
Go to the documentation of this file.
1/*
2 *
3 * ##### ##### ###### ###### ### ###
4 * ## ## ## ## ## ## ## ### ##
5 * ## ## ## ## #### #### ## # ##
6 * ## ## ## ## ## ## ## ##
7 * ## ## ## ## ## ## ## ##
8 * ##### ##### ## ###### ## ##
9 *
10 *
11 * OOFEM : Object Oriented Finite Element Code
12 *
13 * Copyright (C) 1993 - 2025 Borek Patzak
14 *
15 *
16 *
17 * Czech Technical University, Faculty of Civil Engineering,
18 * Department of Structural Mechanics, 166 29 Prague, Czech Republic
19 *
20 * This library is free software; you can redistribute it and/or
21 * modify it under the terms of the GNU Lesser General Public
22 * License as published by the Free Software Foundation; either
23 * version 2.1 of the License, or (at your option) any later version.
24 *
25 * This program is distributed in the hope that it will be useful,
26 * but WITHOUT ANY WARRANTY; without even the implied warranty of
27 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
28 * Lesser General Public License for more details.
29 *
30 * You should have received a copy of the GNU Lesser General Public
31 * License along with this library; if not, write to the Free Software
32 * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
33 */
34
36#include "timestep.h"
37#include "dofmanager.h"
38#include "element.h"
39#include "dof.h"
40#include "verbose.h"
41#include "mathfem.h"
42#include "classfactory.h"
44
45namespace oofem {
46#define ZERO_MASS 1.E-10 // unit dependent !!!!
47
49
50DEIDynamic :: ~DEIDynamic() { }
51
52NumericalMethod *DEIDynamic :: giveNumericalMethod(MetaStep *mStep)
53// only one has reason for DEIDynamic
54// - SolutionOfLinearEquations
55
56{
57 return NULL; // not necessary here - diagonal matrix is used-simple inversion
58}
59
60
61void
62DEIDynamic :: initializeFrom(InputRecord &ir)
63{
64 StructuralEngngModel :: initializeFrom(ir);
65
66 IR_GIVE_FIELD(ir, dumpingCoef, _IFT_DEIDynamic_dumpcoef); // C = dumpingCoef * M
68}
69
70
71double DEIDynamic :: giveUnknownComponent(ValueModeType mode, TimeStep *tStep, Domain *d, Dof *dof)
72// returns unknown quantity like displacement, velocity of equation eq in time t
73// This function translates this request to numerical method language
74{
75 int eq = dof->__giveEquationNumber();
76#ifdef DEBUG
77 if ( eq == 0 ) {
78 OOFEM_ERROR("invalid equation number");
79 }
80#endif
81
82 if ( tStep != this->giveCurrentStep() ) {
83 OOFEM_ERROR("unknown time step encountered");
84 }
85
86 switch ( mode ) {
87 case VM_Total:
88 return displacementVector.at(eq);
89
90 case VM_Velocity:
91 return velocityVector.at(eq);
92
93 case VM_Acceleration:
94 return accelerationVector.at(eq);
95
96 default:
97 OOFEM_ERROR("Unknown is of undefined ValueModeType for this problem");
98 }
99
100 // return 0.0;
101}
102
103
104TimeStep *DEIDynamic :: giveNextStep()
105{
106 int istep = giveNumberOfFirstStep();
107 double totalTime = 0.;
108 StateCounterType counter = 1;
109
110 if ( currentStep ) {
111 totalTime = currentStep->giveTargetTime() + deltaT;
112 istep = currentStep->giveNumber() + 1;
113 counter = currentStep->giveSolutionStateCounter() + 1;
114 }
115
116 previousStep = std :: move(currentStep);
117 currentStep = std::make_unique<TimeStep>(istep, this, 1, totalTime, deltaT, counter);
118
119 return currentStep.get();
120}
121
122
123void DEIDynamic :: solveYourselfAt(TimeStep *tStep)
124{
125 //
126 // creates system of governing eq's and solves them at given time step
127 //
128 // this is an explicit problem: we assemble governing equating at time t
129 // and solution is obtained for time t+dt
130 //
131 // first assemble problem at current time step to obtain results in following
132 // time step.
133 // and then print results for this step also.
134 // for first time step we need special start code
135 Domain *domain = this->giveDomain(1);
136 int nelem = domain->giveNumberOfElements();
137 IntArray loc;
138 Element *element;
139 int neq;
140 int n, init = 0;
141 double coeff, maxDt, maxOmi, maxOm = 0., maxOmEl, c1, c2, c3;
142 FloatMatrix charMtrx, charMtrx2, R;
143 FloatArray previousDisplacementVector;
144
145
147 if ( tStep->giveNumber() == giveNumberOfFirstStep() ) {
148 init = 1;
149#ifdef VERBOSE
150 OOFEM_LOG_INFO("Assembling mass matrix\n");
151#endif
152
153 //
154 // first step assemble mass Matrix
155 //
156
157 massMatrix.resize(neq);
158 massMatrix.zero();
160 for ( int i = 1; i <= nelem; i++ ) {
161 element = domain->giveElement(i);
162 element->giveLocationArray(loc, dn);
163 element->giveCharacteristicMatrix(charMtrx, LumpedMassMatrix, tStep);
164 // charMtrx.beLumpedOf(fullCharMtrx);
165 element->giveCharacteristicMatrix(charMtrx2, TangentStiffnessMatrix, tStep);
166 if ( charMtrx2.isNotEmpty() ) {
168 if ( element->giveRotationMatrix(R) ) {
169 charMtrx2.rotatedWith(R);
170 charMtrx.rotatedWith(R);
171 }
172 }
173
174 //
175 // assemble it manually
176 //
177#ifdef DEBUG
178 if ( loc.giveSize() != charMtrx.giveNumberOfRows() ) {
179 OOFEM_ERROR("dimension mismatch");
180 }
181
182#endif
183
184 n = loc.giveSize();
185
186 maxOmEl = 0.;
187 for ( int j = 1; j <= n; j++ ) {
188 if ( charMtrx.at(j, j) > ZERO_MASS ) {
189 maxOmi = charMtrx2.at(j, j) / charMtrx.at(j, j);
190 if ( init ) {
191 maxOmEl = ( maxOmEl > maxOmi ) ? ( maxOmEl ) : ( maxOmi );
192 }
193 }
194 }
195
196 maxOm = ( maxOm > maxOmEl ) ? ( maxOm ) : ( maxOmEl );
197
198 if (maxOmEl > ZERO_MASS) {
199 // update zero mass diagonal entries ; but skip massless elements
200 for ( int j = 1; j <= n; j++ ) {
201 int jj = loc.at(j);
202 if ( ( jj ) && ( charMtrx.at(j, j) <= ZERO_MASS ) ) {
203 charMtrx.at(j, j) = charMtrx2.at(j, j) / maxOmEl;
204 }
205 }
206 }
207
208 for ( int j = 1; j <= n; j++ ) {
209 int jj = loc.at(j);
210 if ( jj ) {
211 massMatrix.at(jj) += charMtrx.at(j, j);
212 }
213 }
214 }
215
216 // if init - try to determine the best deltaT
217 if ( init ) {
218 maxDt = 2 / sqrt(maxOm);
219 if ( deltaT > maxDt ) {
220 OOFEM_LOG_RELEVANT("DEIDynamic: deltaT reduced to %e\n", maxDt);
221 deltaT = maxDt;
222 tStep->setTimeIncrement(deltaT);
223 }
224 }
225
226
227 //
228 // special init step - compute displacements at tstep 0
229 //
230 displacementVector.resize(neq);
231 displacementVector.zero();
232 nextDisplacementVector.resize(neq);
234 velocityVector.resize(neq);
235 velocityVector.zero();
236 accelerationVector.resize(neq);
237 accelerationVector.zero();
238
239
240 for ( auto &node : domain->giveDofManagers() ) {
241 for ( Dof *iDof: *node ) {
242 // ask for initial values obtained from
243 // bc (boundary conditions) and ic (initial conditions)
244 // now we are setting initial cond. for step -1.
245 if ( !iDof->isPrimaryDof() ) {
246 continue;
247 }
248
249 int jj = iDof->__giveEquationNumber();
250 if ( jj ) {
251 nextDisplacementVector.at(jj) = iDof->giveUnknown(VM_Total, tStep);
252 // become displacementVector after init
253 velocityVector.at(jj) = iDof->giveUnknown(VM_Velocity, tStep);
254 // accelerationVector = iDof->giveUnknown(AccelerartionVector,tStep) ;
255 }
256 }
257 }
258
260
261 return;
262 } // end of init step
263
264#ifdef VERBOSE
265 OOFEM_LOG_INFO("Assembling right hand side\n");
266#endif
267
268
269 c1 = ( 1. / ( deltaT * deltaT ) );
270 c2 = ( 1. / ( 2. * deltaT ) );
271 c3 = ( 2. / ( deltaT * deltaT ) );
272
273 previousDisplacementVector = displacementVector;
275
276 //
277 // assembling the element part of load vector
278 //
280 loadVector.zero();
282 VM_Total, EModelDefaultEquationNumbering(), domain);
284
285
286 //
287 // assembling additional parts of right hand side
288 //
290 for ( int i = 1; i <= nelem; i++ ) {
291 element = domain->giveElement(i);
292 element->giveLocationArray(loc, dn);
293 element->giveCharacteristicMatrix(charMtrx, TangentStiffnessMatrix, tStep);
294 if ( charMtrx.isNotEmpty() ) {
296 if ( element->giveRotationMatrix(R) ) {
297 charMtrx.rotatedWith(R);
298 }
299 }
300
301 n = loc.giveSize();
302 for ( int j = 1; j <= n; j++ ) {
303 int jj = loc.at(j);
304 if ( jj ) {
305 for ( int k = 1; k <= n; k++ ) {
306 int kk = loc.at(k);
307 if ( kk ) {
308 loadVector.at(jj) -= charMtrx.at(j, k) * displacementVector.at(kk);
309 }
310 }
311
312 //
313 // if init step - find minimum period of vibration in order to
314 // determine maximal admissible time step
315 //
316 //maxOmi = charMtrx.at(j,j)/massMatrix.at(jj) ;
317 //if (init) maxOm = (maxOm > maxOmi) ? (maxOm) : (maxOmi) ;
318 }
319 }
320 }
321
322
323
324 for ( int j = 1; j <= neq; j++ ) {
325 coeff = massMatrix.at(j);
326 loadVector.at(j) += coeff * c3 * displacementVector.at(j) -
327 coeff * ( c1 - dumpingCoef * c2 ) *
328 previousDisplacementVector.at(j);
329 }
330
331 //
332 // set-up numerical model
333 //
334 /* it is not necessary to call numerical method
335 * approach used here is not good, but effective enough
336 * inverse of diagonal mass matrix is done here
337 */
338 //
339 // call numerical model to solve arose problem - done locally here
340 //
341#ifdef VERBOSE
342 OOFEM_LOG_RELEVANT( "Solving [step number %8d, time %15e]\n", tStep->giveNumber(), tStep->giveTargetTime() );
343#endif
344 double prevD;
345
346 for ( int i = 1; i <= neq; i++ ) {
347 prevD = previousDisplacementVector.at(i);
348 nextDisplacementVector.at(i) = loadVector.at(i) /
349 ( massMatrix.at(i) * ( c1 + dumpingCoef * c2 ) );
350 velocityVector.at(i) = nextDisplacementVector.at(i) - prevD;
351 accelerationVector.at(i) =
353 2. * displacementVector.at(i) + prevD;
354 }
355
356 accelerationVector.times(c1);
357 velocityVector.times(c2);
358}
359
360
361void DEIDynamic :: printDofOutputAt(FILE *stream, Dof *iDof, TimeStep *tStep)
362{
363 static char dofchar[] = "dva";
364 static ValueModeType dofmodes[] = {
365 VM_Total, VM_Velocity, VM_Acceleration
366 };
367
368 iDof->printMultipleOutputAt(stream, tStep, dofchar, dofmodes, 3);
369}
370} // end namespace oofem
#define REGISTER_EngngModel(class)
FloatArray massMatrix
Definition deidynamic.h:76
FloatArray loadVector
Definition deidynamic.h:77
FloatArray displacementVector
Definition deidynamic.h:79
FloatArray velocityVector
Definition deidynamic.h:79
FloatArray nextDisplacementVector
Definition deidynamic.h:78
FloatArray accelerationVector
Definition deidynamic.h:79
int giveNumberOfFirstStep(bool force=false) override
Definition deidynamic.h:102
virtual void printMultipleOutputAt(FILE *File, TimeStep *tStep, char *ch, ValueModeType *mode, int nite)
Definition dof.C:92
virtual int __giveEquationNumber() const =0
int giveNumberOfElements() const
Returns number of elements in domain.
Definition domain.h:463
std ::vector< std ::unique_ptr< DofManager > > & giveDofManagers()
Definition domain.h:427
Element * giveElement(int n)
Definition domain.C:165
virtual bool giveRotationMatrix(FloatMatrix &answer)
Definition element.C:298
virtual void giveCharacteristicMatrix(FloatMatrix &answer, CharType type, TimeStep *tStep)
Definition element.C:620
void giveLocationArray(IntArray &locationArray, const UnknownNumberingScheme &s, IntArray *dofIds=NULL) const
Definition element.C:429
virtual TimeStep * giveCurrentStep(bool force=false)
Definition engngm.h:717
virtual int giveNumberOfDomainEquations(int di, const UnknownNumberingScheme &num)
Definition engngm.C:452
virtual void init()
Definition engngm.C:2082
std ::unique_ptr< TimeStep > previousStep
Previous time step.
Definition engngm.h:243
Domain * giveDomain(int n)
Definition engngm.C:1936
std ::unique_ptr< TimeStep > currentStep
Current time step.
Definition engngm.h:241
int updateSharedDofManagers(FloatArray &answer, const UnknownNumberingScheme &s, int ExchangeTag)
Definition engngm.C:2162
void assembleVector(FloatArray &answer, TimeStep *tStep, const VectorAssembler &va, ValueModeType mode, const UnknownNumberingScheme &s, Domain *domain, FloatArray *eNorms=NULL)
Definition engngm.C:1106
double & at(Index i)
Definition floatarray.h:202
void rotatedWith(const FloatMatrix &r, char mode='n')
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
int & at(std::size_t i)
Definition intarray.h:104
int giveSize() const
Definition intarray.h:211
double giveTargetTime()
Returns target time.
Definition timestep.h:164
void setTimeIncrement(double newDt)
Sets solution step time increment.
Definition timestep.h:170
int giveNumber()
Returns receiver's number.
Definition timestep.h:144
#define ZERO_MASS
Definition deidynamic.C:46
#define _IFT_DEIDynamic_deltat
Definition deidynamic.h:44
#define _IFT_DEIDynamic_dumpcoef
Definition deidynamic.h:43
#define OOFEM_ERROR(...)
Definition error.h:79
#define IR_GIVE_FIELD(__ir, __value, __id)
Definition inputrecord.h:67
#define OOFEM_LOG_INFO(...)
Definition logger.h:143
#define OOFEM_LOG_RELEVANT(...)
Definition logger.h:142
long StateCounterType
StateCounterType type used to indicate solution state.

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