OOFEM 3.0
Loading...
Searching...
No Matches
phasefieldelement.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 *
14 * Copyright (C) 1993 - 2025 Borek Patzak
15 *
16 *
17 *
18 * Czech Technical University, Faculty of Civil Engineering,
19 * Department of Structural Mechanics, 166 29 Prague, Czech Republic
20 *
21 * This library is free software; you can redistribute it and/or
22 * modify it under the terms of the GNU Lesser General Public
23 * License as published by the Free Software Foundation; either
24 * version 2.1 of the License, or (at your option) any later version.
25 *
26 * This program is distributed in the hope that it will be useful,
27 * but WITHOUT ANY WARRANTY; without even the implied warranty of
28 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
29 * Lesser General Public License for more details.
30 *
31 * You should have received a copy of the GNU Lesser General Public
32 * License along with this library; if not, write to the Free Software
33 * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
34 */
35
39#include "node.h"
40#include "material.h"
41#include "gausspoint.h"
43#include "domain.h"
44#include "mathfem.h"
45#include "timestep.h"
46#include <cstdio>
47
48
49namespace oofem {
50
51PhaseFieldElement :: PhaseFieldElement( int i, Domain *aDomain )
52{
54 internalLength = 6.0;
55 criticalEnergy = 1.0e0;
56 relaxationTime = 1.0;
57}
58
59void
60PhaseFieldElement :: computeLocationArrayOfDofIDs( const IntArray &dofIdArray, IntArray &answer )
61{
62 // Routine to extract compute the location array an element given an dofid array.
63 answer.clear();
64 NLStructuralElement *el = this->giveElement();
65 int k = 0;
66 for ( int i = 1; i <= el->giveNumberOfDofManagers(); i++ ) {
67 DofManager *dMan = el->giveDofManager( i );
68 for ( int j = 1; j <= dofIdArray.giveSize( ); j++ ) {
69
70 if ( dMan->hasDofID( (DofIDItem) dofIdArray.at( j ) ) ) {
71 Dof *d = dMan->giveDofWithID( dofIdArray.at( j ) );
72 answer.followedBy( k + d->giveNumber( ) );
73 }
74 }
75 k += dMan->giveNumberOfDofs( );
76 }
77}
78
79void
81PhaseFieldElement :: computeDisplacementUnknowns(FloatArray &answer, ValueModeType valueMode, TimeStep *stepN)
82{
83 IntArray dofIdArray;
84 this->giveDofManDofIDMask_u(dofIdArray);
85 Element :: computeVectorOf(dofIdArray, valueMode, stepN, answer);
86}
87
88void
89PhaseFieldElement :: computeDamageUnknowns(FloatArray &answer, ValueModeType valueMode, TimeStep *stepN)
90{
91 IntArray dofIdArray;
92 this->giveDofManDofIDMask_d(dofIdArray);
93 Element :: computeVectorOf(dofIdArray, valueMode, stepN, answer);
94}
95
96void
97PhaseFieldElement :: giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord)
98{
99 // Computes the internal forces corresponding to the two fields u & d
100 IntArray IdMask_u, IdMask_d;
101 this->giveDofManDofIDMask_u( IdMask_u );
102 this->giveDofManDofIDMask_d( IdMask_d );
103 this->computeLocationArrayOfDofIDs( IdMask_u, loc_u );
104 this->computeLocationArrayOfDofIDs( IdMask_d, loc_d );
105
106 int ndofs = this->giveElement()->computeNumberOfDofs();
107 answer.resize(ndofs);
108 answer.zero();
109
110 FloatArray answer_u(0);
111 FloatArray answer_d(0);
112
113 this->giveInternalForcesVector_u(answer_u, tStep, useUpdatedGpRecord);
114 this->giveInternalForcesVector_d(answer_d, tStep, useUpdatedGpRecord);
115 answer.assemble(answer_u, loc_u);
116 answer.assemble(answer_d, loc_d);
117}
118
119void
120PhaseFieldElement :: giveInternalForcesVector_u(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord)
121{
122 // computes int_V ( B_u^t * BSgima_u ) * dV
123 FloatArray NStress, BStress, vGenStress, NS;
124 FloatMatrix N, B;
125 NLStructuralElement *el = this->giveElement();
126
127 answer.clear();
128 for ( auto &gp: el->giveIntegrationRule(0) ) {
129 double dV = el->computeVolumeAround(gp);
130
131 // compute generalized stress measure
132 el->computeBmatrixAt(gp, B);
133 this->computeBStress_u(BStress, gp, tStep, useUpdatedGpRecord);
134 answer.plusProduct(B, BStress, dV);
135 }
136}
137
138void
139PhaseFieldElement :: giveInternalForcesVector_d(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord)
140{
141 // Computes int_V ( N^t *Nstress_d + B^t * g_c*l*grad(d) ) dV
142 FloatArray NStress, BStress, NS, a_d, grad_d;
143 FloatMatrix N, B;
144 NLStructuralElement *el = this->giveElement();
145 this->computeDamageUnknowns( a_d, VM_Total, tStep );
146
147 for ( auto &gp: el->giveIntegrationRule(0) ) {
148 double dV = el->computeVolumeAround(gp);
149
150 // compute generalized stress measures
151 this->computeNd_matrixAt( *gp->giveNaturalCoordinates(), N);
152 computeNStress_d(NStress, gp, tStep, useUpdatedGpRecord);
153 NS.beTProductOf(N, NStress);
154 answer.add(dV, NS);
155
156 this->computeBd_matrixAt(gp, B);
157 grad_d.beProductOf(B, a_d);
158 double l = this->giveInternalLength();
159 double g_c = this->giveCriticalEnergy();
160 BStress = grad_d * l * g_c;
161 answer.plusProduct(B, BStress, dV);
162 }
163}
164
165void
166PhaseFieldElement :: computeBStress_u(FloatArray &answer, GaussPoint *gp, TimeStep *tStep, int useUpdatedGpRecord)
167{
168 // computes G(d)*sig(u)
169 NLStructuralElement *el = this->giveElement();
170 FloatArray strain, a_u;
171 FloatMatrix B_u;
172
173 el->computeBmatrixAt(gp, B_u);
174 this->computeDisplacementUnknowns(a_u, VM_Total, tStep);
175 strain.beProductOf(B_u, a_u);
176 el->computeStressVector(answer, strain, gp, tStep);
177 answer.times( this->computeG(gp, VM_Total, tStep) );
178
179}
180
181double
182PhaseFieldElement :: computeFreeEnergy(GaussPoint *gp, TimeStep *tStep)
183{
184 StructuralMaterialStatus *matStat = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() );
185 FloatArray strain, stress;
186 stress = matStat->giveTempStressVector();
187 strain = matStat->giveTempStrainVector();
188 return 0.5 * stress.dotProduct( strain );
189}
190
191void
192PhaseFieldElement :: computeNStress_d(FloatArray &answer, GaussPoint *gp, TimeStep *tStep, int useUpdatedGpRecord)
193{
194 // computes (t*/Delta t) *(d - d_old) + g_c/l*d + G'*Psibar
195
196 //PhaseFieldCrossSection *cs = static_cast...
197 double Delta_t = tStep->giveTimeIncrement();
198 double t_star = this->giveRelaxationTime();
199 double d = computeDamageAt( gp, VM_Total, tStep );
200 double Delta_d = computeDamageAt( gp, VM_Incremental, tStep );
201
202 double l = this->giveInternalLength();
203 double g_c = this->giveCriticalEnergy();
204 double Gprim = this->computeGPrim(gp, VM_Total, tStep);
205 double Psibar = this->computeFreeEnergy( gp, tStep );
206 answer.resize( 1 );
207 answer.at( 1 ) = t_star / Delta_t * Delta_d + g_c / l *d + Gprim * Psibar;
208}
209
210double
211PhaseFieldElement :: computeDamageAt(GaussPoint *gp, ValueModeType valueMode, TimeStep *stepN)
212{
213 // d = N_d * a_d
214 NLStructuralElement *el = this->giveElement();
215 FloatArray dVec;
216 computeDamageUnknowns(dVec, valueMode, stepN);
217 FloatArray Nvec;
219 return Nvec.dotProduct(dVec);
220}
221
222double
223PhaseFieldElement :: computeG(GaussPoint *gp, ValueModeType valueMode, TimeStep *stepN)
224{
225 // computes Dg/Dd = (1-d)^2 + r0
226 double d = this->computeDamageAt(gp, valueMode, stepN);
227 double r0 = 1.0e-10;
228 return (1.0 - d) * (1.0 - d) + r0;
229}
230
231double
232PhaseFieldElement :: computeGPrim(GaussPoint *gp, ValueModeType valueMode, TimeStep *stepN)
233{
234 // compute -2*(1-d)
235 double d = this->computeDamageAt(gp, valueMode, stepN);
236 return -2.0 * (1.0 - d);
237}
238
239void
240PhaseFieldElement :: computeNd_matrixAt(const FloatArray &lCoords, FloatMatrix &N)
241{
242 NLStructuralElement *el = this->giveElement();
243 FloatArray Nvec;
244 el->giveInterpolation( )->evalN( Nvec, lCoords, FEIElementGeometryWrapper( el ) );
245 N.resize(1, Nvec.giveSize());
246 N.beNMatrixOf(Nvec,1);
247
248}
249
250void
251PhaseFieldElement :: computeBd_matrixAt(GaussPoint *gp, FloatMatrix &answer, int li, int ui)
252{
253 // Returns the [numSpaceDim x nDofs] gradient matrix {B_d} of the receiver,
254 // evaluated at gp.
255
256 NLStructuralElement *el = this->giveElement();
257 FloatMatrix dNdx;
259 answer.beTranspositionOf( dNdx );
260}
261
262void
263PhaseFieldElement :: computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
264{
265 //set displacement and nonlocal location array
267 IntArray IdMask_u, IdMask_d;
268 this->giveDofManDofIDMask_u( IdMask_u );
269 this->giveDofManDofIDMask_d( IdMask_d );
270 this->computeLocationArrayOfDofIDs( IdMask_u, loc_u );
271 this->computeLocationArrayOfDofIDs( IdMask_d, loc_d );
272
273 int nDofs = this->giveElement()->computeNumberOfDofs();
274 answer.resize( nDofs, nDofs );
275 answer.zero();
276
277 FloatMatrix answer1, answer2, answer3, answer4;
278 this->computeStiffnessMatrix_uu(answer1, rMode, tStep);
279 this->computeStiffnessMatrix_ud(answer2, rMode, tStep);
280 //this->computeStiffnessMatrix_du(answer3, rMode, tStep); //symmetric
281 answer3.beTranspositionOf( answer2 );
282 this->computeStiffnessMatrix_dd(answer4, rMode, tStep);
283
284 answer.assemble( answer1, loc_u, loc_u );
285 answer.assemble( answer2, loc_u, loc_d );
286 answer.assemble( answer3, loc_d, loc_u );
287 answer.assemble( answer4, loc_d, loc_d );
288}
289
290void
291PhaseFieldElement :: computeStiffnessMatrix_uu(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
292{
293 // This is the regular stiffness matrix times G
294 FloatMatrix B, DB, N, D_B;
295 NLStructuralElement *el = this->giveElement();
297
298 bool matStiffSymmFlag = cs->isCharacteristicMtrxSymmetric(rMode);
299 answer.clear();
300
301 for ( auto gp: el->giveIntegrationRule(0) ) {
302 double dV = el->computeVolumeAround(gp);
303 // compute int_V ( B^t * D_B * B )dV
304 el->computeBmatrixAt(gp, B );
305 cs->giveCharMaterialStiffnessMatrix(D_B, rMode, gp, tStep);
306 D_B.times( computeG(gp, VM_Total, tStep) );
307 DB.beProductOf(D_B, B);
308
309 if ( matStiffSymmFlag ) {
310 answer.plusProductSymmUpper(B, DB, dV);
311 } else {
312 answer.plusProductUnsym(B, DB, dV);
313 }
314 }
315
316 if ( matStiffSymmFlag ) {
317 answer.symmetrized();
318 }
319}
320
321void
322PhaseFieldElement :: computeStiffnessMatrix_ud(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
323{
324 FloatMatrix B, DB, N_d, DN, S(3,1);
325 FloatArray stress;
326 NLStructuralElement *el = this->giveElement();
328
329 answer.clear();
330
331 for ( auto &gp: el->giveIntegrationRule(0) ) {
332 StructuralMaterialStatus *matStat = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() );
333
334 double dV = el->computeVolumeAround(gp);
335 // compute int_V ( B^t * D_B * B )dV
336
337 this->computeNd_matrixAt(*gp->giveNaturalCoordinates(), N_d);
338
339 // stress
340 FloatArray reducedStrain, a_u;
341 FloatMatrix B_u;
342 el->computeBmatrixAt( gp, B_u );
343 stress = matStat->giveTempStressVector();
344 stress.times( this->computeGPrim( gp, VM_Total, tStep ) );
345
346 S.setColumn(stress,1);
347 DN.beProductOf(S, N_d);
348 answer.plusProductUnsym(B_u, DN, dV);
349 }
350
351}
352
353void
354PhaseFieldElement :: computeStiffnessMatrix_dd(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
355{
356 double Delta_t = tStep->giveTimeIncrement();
357 double t_star = this->giveRelaxationTime();
358 double l = this->giveInternalLength();
359 double g_c = this->giveCriticalEnergy();
360
361 FloatMatrix B_d, N_d;
362 //StructuralCrossSection *cs = dynamic_cast<StructuralCrossSection* > (this->giveCrossSection() );
363 NLStructuralElement *el = this->giveElement();
364 answer.clear();
365
366 for ( auto &gp: el->giveIntegrationRule(0) ) {
367 double dV = el->computeVolumeAround(gp);
368
369 this->computeNd_matrixAt(*gp->giveNaturalCoordinates(), N_d);
370 this->computeBd_matrixAt(gp, B_d, 1, 3);
371
372 double Gprim = this->computeGPrim(gp, VM_Total, tStep);
373 double psiBar = this->computeFreeEnergy( gp, tStep );
374 //double factorN = t_star/Delta_t + g_c/l + psiBar*Gbis;
375 double factorN = t_star / Delta_t + g_c / l + psiBar*(-2.0);
376 double factorB = g_c*l;
377
378 answer.plusProductSymmUpper(N_d, N_d, factorN*dV);
379 answer.plusProductSymmUpper(B_d, B_d, factorB*dV);
380 }
381
382 answer.symmetrized();
383}
384
385
386#if 0
387void
388PhaseFieldElement :: computeStressVectorAndLocalCumulatedStrain(FloatArray &answer, double localCumulatedStrain, GaussPoint *gp, TimeStep *stepN)
389//PhaseFieldElement :: computeStressVector(FloatArray &answer, double localCumulatedStrain, GaussPoint *gp, TimeStep *stepN)
390{
391 NLStructuralElement *elem = this->giveNLStructuralElement();
392
393 nlGeo = elem->giveGeometryMode();
394 StructuralCrossSection *cs = this->giveNLStructuralElement()->giveStructuralCrossSection();
395
396 //this->computeDamage(alpha, gp, stepN);
397 if ( nlGeo == 0 ) {
398 FloatArray Epsilon;
399 this->computeLocalStrainVector(Epsilon, gp, stepN);
400 dpmat->giveRealStressVector(answer, gp, Epsilon, nlCumulatedStrain, stepN);
401 } else if ( nlGeo == 1 ) {
402 if ( elem->giveDomain()->giveEngngModel()->giveFormulation() == TL ) {
403 FloatArray vF;
404 this->computeDeformationGradientVector(vF, gp, stepN);
405 dpmat->giveFirstPKStressVector(answer, gp, vF, nlCumulatedStrain, stepN);
406 }
407
408 }
409}
410
411
412void
413PhaseFieldElement :: giveInternalForcesVector_u(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord)
414{
415 // computes int_V ( B_u^t * BSgima_u ) * dV
416 this->giveInternalForcesVectorGen(answer, tStep, useUpdatedGpRecord,
417 NULL,
418 &this->computeBmatrixAt,
419 NULL,
420 &this->computeBStress_u,
421 &this->computeVolumeAround
422 );
423}
424
425void
426PhaseFieldElement :: giveInternalForcesVector_d(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord)
427{
428 // computes int_V ( N_d^t * NSgima_d + B_d^t * BSgima_d ) * dV
429 this->giveInternalForcesVectorGen(answer, tStep, useUpdatedGpRecord,
430 &this->computeNmatrixAt,
431 &this->computeBmatrixAt,
432 &this->computeNStress_d,
433 &this->computeBStress_d,
434 &this->computeVolumeAround
435 );
436}
437
438void
439PhaseFieldElement :: computeBStress_d(FloatArray &answer, GaussPoint *gp, TimeStep *tStep, int useUpdatedGpRecord)
440{
441 // computes g_c*l
442 // PhaseFieldCrossSection *cs = static_cast...
443 double l = this->giveInternalLength();
444 double g_c = this->giveCriticalEnergy();
445 answer.resize(1);
446 answer.at(1) = g_c*l;
447}
448
449void
450PhaseFieldElement :: computeStiffnessMatrix_du(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
451{
452
453}
454void
455PhaseFieldElement :: computeStiffnessMatrix_dd(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
456{
457 StructuralCrossSection *cs = this->giveStructuralCrossSection();
458
459 this->computeStiffnessMatrixGen(answer, rMode, tStep,
460 NULL,
461 &this->computeBmatrixAt,
462 NULL,
463 &this->Duu_B,
464 this->computeVolumeAround
465 );
466}
467
468void
469PhaseFieldElement :: computeStiffnessMatrix_uu(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
470{
471 // This is the regular stiffness matrix times G
472
473 this->computeStiffnessMatrixGen(answer, rMode, tStep,
474 NULL,
475 &this->computeBmatrixAt,
476 NULL,
477 &this->Duu_B,
478 this->computeVolumeAround
479 );
480}
481
482void
483PhaseFieldElement :: Duu_B(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
484{
485 // G*dSig/dEps
486 //StructuralCrossSection *cs = this->giveStructuralCrossSection();
487 //cs->giveCharMaterialStiffnessMatrix(answer, rMode, gp, tStep);
488
489 //answer.times( this->computeG(gp) );
490}
491#endif
492
493
494} // end namespace oofem
#define N(a, b)
bool hasDofID(DofIDItem id) const
Definition dofmanager.C:174
int giveNumberOfDofs() const
Definition dofmanager.C:287
Dof * giveDofWithID(int dofID) const
Definition dofmanager.C:127
EngngModel * giveEngngModel()
Definition domain.C:419
virtual FEInterpolation * giveInterpolation() const
Definition element.h:648
virtual IntegrationRule * giveIntegrationRule(int i)
Definition element.h:899
virtual int giveNumberOfDofManagers() const
Definition element.h:695
DofManager * giveDofManager(int i) const
Definition element.C:553
CrossSection * giveCrossSection()
Definition element.C:534
virtual double computeVolumeAround(GaussPoint *gp)
Definition element.h:530
virtual fMode giveFormulation()
Definition engngm.h:1165
virtual void evalN(FloatArray &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const =0
virtual double evaldNdx(FloatMatrix &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const =0
Domain * giveDomain() const
Definition femcmpnn.h:97
void assemble(const FloatArray &fe, const IntArray &loc)
Definition floatarray.C:616
void resize(Index s)
Definition floatarray.C:94
double & at(Index i)
Definition floatarray.h:202
void plusProduct(const FloatMatrix &b, const FloatArray &s, double dV)
Definition floatarray.C:288
double dotProduct(const FloatArray &x) const
Definition floatarray.C:524
Index giveSize() const
Returns the size of receiver.
Definition floatarray.h:261
void zero()
Zeroes all coefficients of receiver.
Definition floatarray.C:683
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Definition floatarray.C:689
void beTProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Definition floatarray.C:721
void add(const FloatArray &src)
Definition floatarray.C:218
void times(double s)
Definition floatarray.C:834
void times(double f)
void plusProductSymmUpper(const FloatMatrix &a, const FloatMatrix &b, double dV)
void resize(Index rows, Index cols)
Definition floatmatrix.C:79
*Sets size of receiver to be an empty matrix It will have zero rows and zero columns size void clear()
void plusProductUnsym(const FloatMatrix &a, const FloatMatrix &b, double dV)
void beProductOf(const FloatMatrix &a, const FloatMatrix &b)
void beTranspositionOf(const FloatMatrix &src)
void zero()
Zeroes all coefficient of receiver.
void assemble(const FloatMatrix &src, const IntArray &loc)
const FloatArray & giveNaturalCoordinates() const
Returns coordinate array of receiver.
Definition gausspoint.h:138
IntegrationPointStatus * giveMaterialStatus(IntegrationPointStatusIDType key=IPSID_Default)
Definition gausspoint.h:204
void followedBy(const IntArray &b, int allocChunk=0)
Definition intarray.C:94
int & at(std::size_t i)
Definition intarray.h:104
int giveSize() const
Definition intarray.h:211
void giveInternalForcesVector_u(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord)
void giveInternalForcesVector_d(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord)
double computeDamageAt(GaussPoint *gp, ValueModeType valueMode, TimeStep *stepN)
void computeDisplacementUnknowns(FloatArray &answer, ValueModeType valueMode, TimeStep *stepN)
void computeStiffnessMatrix_dd(FloatMatrix &, MatResponseMode, TimeStep *)
virtual NLStructuralElement * giveElement()=0
void computeNStress_d(FloatArray &answer, GaussPoint *gp, TimeStep *tStep, int useUpdatedGpRecord)
void computeDamageUnknowns(FloatArray &answer, ValueModeType valueMode, TimeStep *stepN)
virtual void giveDofManDofIDMask_u(IntArray &answer)=0
void computeLocationArrayOfDofIDs(const IntArray &dofIdArray, IntArray &answer)
virtual void computeNd_matrixAt(const FloatArray &lCoords, FloatMatrix &N)
double computeFreeEnergy(GaussPoint *gp, TimeStep *tStep)
double computeGPrim(GaussPoint *gp, ValueModeType valueMode, TimeStep *stepN)
virtual void computeBd_matrixAt(GaussPoint *, FloatMatrix &, int=1, int=ALL_STRAINS)
void computeStiffnessMatrix_ud(FloatMatrix &, MatResponseMode, TimeStep *)
void computeBStress_u(FloatArray &answer, GaussPoint *gp, TimeStep *tStep, int useUpdatedGpRecord)
virtual void giveDofManDofIDMask_d(IntArray &answer)=0
double computeG(GaussPoint *gp, ValueModeType valueMode, TimeStep *stepN)
void computeStiffnessMatrix_uu(FloatMatrix &, MatResponseMode, TimeStep *)
bool isCharacteristicMtrxSymmetric(MatResponseMode mode) const override=0
virtual void giveCharMaterialStiffnessMatrix(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)=0
virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int lowerIndx=1, int upperIndx=ALL_STRAINS)=0
virtual void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep)=0
const FloatArray & giveTempStressVector() const
Returns the const pointer to receiver's temporary stress vector.
const FloatArray & giveTempStrainVector() const
Returns the const pointer to receiver's temporary strain vector.
double giveTimeIncrement()
Returns solution step associated time increment.
Definition timestep.h:168
#define S(p)
Definition mdm.C:469
@ TL
Total Lagrange.
Definition fmode.h:44

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