OOFEM 3.0
Loading...
Searching...
No Matches
trabbonenl3d.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#include "trabbonenl3d.h"
36#include "gausspoint.h"
37#include "floatmatrix.h"
38#include "floatarray.h"
39#include "mathfem.h"
40#include "sparsemtrx.h"
41#include "nonlocalmaterialext.h"
42#include "classfactory.h"
43#include "dynamicinputrecord.h"
45
46#ifdef __OOFEG
47 #include "oofeggraphiccontext.h"
48 #include "connectivitytable.h"
49#endif
50
51#include <cstdlib>
52
53namespace oofem {
55
57{
58}
59
60
61void
62TrabBoneNL3D :: updateBeforeNonlocAverage(const FloatArray &strainVector, GaussPoint *gp, TimeStep *tStep) const
63{
64 auto nlStatus = static_cast< TrabBoneNL3DStatus * >( this->giveStatus(gp) );
65 FloatArray SDstrainVector;
66
67 this->initTempStatus(gp);
68 this->giveStressDependentPartOfStrainVector(SDstrainVector, gp, strainVector, tStep, VM_Total);
69
70 nlStatus->letTempStrainVectorBe(strainVector);
71
72 this->performPlasticityReturn(gp, strainVector, tStep);
73 double cumPlastStrain = this->computeLocalCumPlastStrain(strainVector, gp, tStep);
74 nlStatus->setLocalCumPlastStrainForAverage(cumPlastStrain);
75}
76
77
79TrabBoneNL3D :: give3dMaterialStiffnessMatrix(MatResponseMode mode, GaussPoint *gp,
80 TimeStep *tStep) const
81{
82 auto nlStatus = static_cast< TrabBoneNL3DStatus * >( this->giveStatus(gp) );
83
84 if ( mode == ElasticStiffness ) {
85 auto compliance = this->constructAnisoComplTensor();
86 auto elasticity = inv(compliance);
87 return elasticity;
88 } else if ( mode == SecantStiffness ) {
89 auto compliance = this->constructAnisoComplTensor();
90 auto elasticity = inv(compliance);
91 auto tempDam = nlStatus->giveTempDam();
92 return elasticity * (1.0 - tempDam);
93 } else /*if ( mode == TangentStiffness )*/ {
94 double kappa = nlStatus->giveKappa();
95 double tempKappa = nlStatus->giveTempKappa();
96 double dKappa = tempKappa - kappa;
97 if ( dKappa < 10.e-9 ) {
98 dKappa = 0;
99 }
100
101 if ( dKappa > 0.0 ) {
102 // Imports
103 auto &tempEffectiveStress = nlStatus->giveTempEffectiveStress();
104 double nlKappa = this->computeCumPlastStrain(gp, tStep);
105 double tempDam = nlStatus->giveTempDam();
106 double dam = nlStatus->giveDam();
107 auto &plasFlowDirec = nlStatus->givePlasFlowDirec();
108 auto &SSaTensor = nlStatus->giveSSaTensor();
109 double beta = nlStatus->giveBeta();
110
111 // Construction of the dyadic product tensor
112 auto prodTensor = Tdot(SSaTensor, plasFlowDirec);
113 // Construction of the tangent stiffness second term
114 auto secondTerm = dyad(dot(SSaTensor, plasFlowDirec), prodTensor) * (-( 1.0 - tempDam ) / beta);
115
116 auto tangentMatrix = SSaTensor * (1.0 - tempDam) + secondTerm;
117
118 // Construction of the tangent stiffness third term
119 if ( tempDam - dam > 0 ) {
120 tangentMatrix += dyad(tempEffectiveStress, prodTensor) * (-expDam * critDam * exp(-expDam * nlKappa) * ( 1.0 - mParam ) / beta);
121
122 }
123 return tangentMatrix;
124 } else {
125 // Import of state variables
126 double tempDam = nlStatus->giveTempDam();
127 // Construction of the tangent stiffness
128 auto compliance = this->constructAnisoComplTensor();
129 auto elasticity = inv(compliance);
130 return elasticity * (1.0 - tempDam);
131 }
132 }
133}
134
135
136void
137TrabBoneNL3D :: NonlocalMaterialStiffnessInterface_addIPContribution(SparseMtrx &dest, const UnknownNumberingScheme &s, GaussPoint *gp, TimeStep *tStep)
138{
139 auto nlStatus = static_cast< TrabBoneNL3DStatus * >( this->giveStatus(gp) );
140 auto list = nlStatus->giveIntegrationDomainList();
141
142 FloatArray rcontrib, lcontrib;
143 IntArray loc, rloc;
144
145 FloatMatrix contrib;
146
147 if ( this->giveLocalNonlocalStiffnessContribution(gp, loc, s, lcontrib, tStep) == 0 ) {
148 return;
149 }
150
151 for ( auto &lir: *list ) {
152 auto rmat = dynamic_cast< TrabBoneNL3D * >( lir.nearGp->giveMaterial() );
153 if ( rmat ) {
154 rmat->giveRemoteNonlocalStiffnessContribution(lir.nearGp, rloc, s, rcontrib, tStep);
155 double coeff = gp->giveElement()->computeVolumeAround(gp) * lir.weight / nlStatus->giveIntegrationScale();
156
157 contrib.clear();
158 contrib.plusDyadUnsym(lcontrib, rcontrib, - 1.0 * coeff);
159 dest.assemble(loc, rloc, contrib);
160 }
161 }
162}
163
164std :: vector< localIntegrationRecord > *
165TrabBoneNL3D :: NonlocalMaterialStiffnessInterface_giveIntegrationDomainList(GaussPoint *gp)
166{
167 auto nlStatus = static_cast< TrabBoneNL3DStatus * >( this->giveStatus(gp) );
168 this->buildNonlocalPointTable(gp);
169 return nlStatus->giveIntegrationDomainList();
170}
171
172int
173TrabBoneNL3D :: giveLocalNonlocalStiffnessContribution(GaussPoint *gp, IntArray &loc, const UnknownNumberingScheme &s,
174 FloatArray &lcontrib, TimeStep *tStep)
175{
176 auto nlStatus = static_cast< TrabBoneNL3DStatus * >( this->giveStatus(gp) );
177 auto elem = static_cast< StructuralElement * >( gp->giveElement() );
178
179 double nlKappa = this->computeCumPlastStrain(gp, tStep);
180 double dam = nlStatus->giveDam();
181 double tempDam = nlStatus->giveTempDam();
182
183 if ( ( tempDam - dam ) > 0.0 ) {
184 FloatMatrix b;
185
186 elem->giveLocationArray(loc, s);
187 auto &localNu = nlStatus->giveTempEffectiveStress();
188
189 elem->giveLocationArray(loc, EModelDefaultEquationNumbering() );
190 elem->computeBmatrixAt(gp, b);
191 double dDamFunc = expDam * critDam * exp(-expDam * nlKappa);
192
193 int nrows = b.giveNumberOfColumns();
194 int nsize = localNu.giveSize();
195 lcontrib.resize(nrows);
196 for ( int i = 1; i <= nrows; i++ ) {
197 double sum = 0.0;
198 for ( int j = 1; j <= nsize; j++ ) {
199 sum += b.at(j, i) * localNu.at(j);
200 }
201
202 lcontrib.at(i) = mParam * dDamFunc * sum;
203 }
204
205 return 1;
206 } else {
207 loc.clear();
208 return 0;
209 }
210}
211
212void
213TrabBoneNL3D :: giveRemoteNonlocalStiffnessContribution(GaussPoint *gp, IntArray &rloc, const UnknownNumberingScheme &s,
214 FloatArray &rcontrib, TimeStep *tStep)
215{
216 auto nlStatus = static_cast< TrabBoneNL3DStatus * >( this->giveStatus(gp) );
217 auto elem = static_cast< StructuralElement * >( gp->giveElement() );
218
219 FloatMatrix b;
220
221 elem->giveLocationArray(rloc, s);
222 elem->computeBmatrixAt(gp, b);
223
224 double kappa = nlStatus->giveKappa();
225 double tempKappa = nlStatus->giveTempKappa();
226 double dKappa = tempKappa - kappa;
227 if ( dKappa < 10.e-9 ) {
228 dKappa = 0;
229 }
230
231 if ( dKappa > 0.0 ) {
232 FloatArray remoteNu, prodTensor;
233 const FloatArray &plasFlowDirec = nlStatus->givePlasFlowDirec();
234 const FloatMatrix &SSaTensor = nlStatus->giveSSaTensor();
235 double beta = nlStatus->giveBeta();
236
237 prodTensor.beTProductOf(SSaTensor, plasFlowDirec);
238 remoteNu = 1 / beta * prodTensor;
239 rcontrib.beTProductOf(b, remoteNu);
240 } else {
241 rcontrib.resize(b.giveNumberOfColumns());
242 rcontrib.zero();
243 }
244}
245
246
248TrabBoneNL3D :: giveRealStressVector_3d(const FloatArrayF<6> &strain, GaussPoint *gp,
249 TimeStep *tStep) const
250{
251 auto nlStatus = static_cast< TrabBoneNL3DStatus * >( this->giveStatus(gp) );
252
253 this->initTempStatus(gp);
254
255 performPlasticityReturn(gp, strain, tStep);
256 auto tempDam = computeDamage(gp, tStep);
257 auto &effStress = nlStatus->giveTempEffectiveStress();
258
259 auto stress = ( 1 - tempDam ) * effStress;
260
261 for ( int i = 1; i <= 6; i++ ) {
262 if ( sqrt( stress.at(i) * stress.at(i) ) < 1e-8 ) {
263 stress.at(i) = 0.;
264 }
265 }
266
267 computePlasStrainEnerDensity(gp, strain, stress);
268
269 if ( densCrit != 0. ) {
270 stress += computeDensificationStress(gp, strain, tStep);
271 }
272
273 nlStatus->setTempDam(tempDam);
274 nlStatus->letTempStrainVectorBe(strain);
275 nlStatus->letTempStressVectorBe(stress);
276 return stress;
277}
278
279
280double
281TrabBoneNL3D :: computeCumPlastStrain(GaussPoint *gp, TimeStep *tStep) const
282{
283 auto nlStatus = static_cast< TrabBoneNL3DStatus * >( this->giveStatus(gp) );
284
285 this->buildNonlocalPointTable(gp);
287
288 auto list = nlStatus->giveIntegrationDomainList();
289
290 double nonlocalCumPlastStrain = 0.0;
291 for ( auto &lir: *list ) {
292 auto nonlocStatus = static_cast< TrabBoneNL3DStatus * >( this->giveStatus(lir.nearGp) );
293 double nonlocalContribution = nonlocStatus->giveLocalCumPlastStrainForAverage();
294 nonlocalContribution *= lir.weight;
295 nonlocalCumPlastStrain += nonlocalContribution;
296 }
297
298 nonlocalCumPlastStrain *= 1. / nlStatus->giveIntegrationScale();
299
300 double localCumPlastStrain = nlStatus->giveLocalCumPlastStrainForAverage();
301 return mParam * nonlocalCumPlastStrain + ( 1 - mParam ) * localCumPlastStrain;
302}
303
304
305Interface *
306TrabBoneNL3D :: giveInterface(InterfaceType type)
307{
309 return static_cast< StructuralNonlocalMaterialExtensionInterface * >(this);
310 } else if ( type == NonlocalMaterialStiffnessInterfaceType ) {
311 return static_cast< NonlocalMaterialStiffnessInterface * >(this);
312 } else {
313 return nullptr;
314 }
315}
316
317
318void
319TrabBoneNL3D :: initializeFrom(InputRecord &ir)
320{
321 TrabBone3D :: initializeFrom(ir);
322 StructuralNonlocalMaterialExtensionInterface :: initializeFrom(ir);
323
325 if ( R < 0.0 ) {
326 R = 0.0;
327 }
328
329 mParam = 2.;
331}
332
333
334void
335TrabBoneNL3D :: giveInputRecord(DynamicInputRecord &input)
336{
337 TrabBone3D :: giveInputRecord(input);
338 input.setField(this->R, _IFT_TrabBoneNL3D_r);
339 input.setField(this->mParam, _IFT_TrabBoneNL3D_m);
340}
341
342
343
344double
345TrabBoneNL3D :: computeWeightFunction(const double R, const FloatArray &src, const FloatArray &coord) const
346{
347 double dist = distance(src, coord);
348
349 if ( ( dist >= 0. ) && ( dist <= this->R ) ) {
350 double help = ( 1. - dist * dist / ( R * R ) );
351 return help * help;
352 }
353
354 return 0.0;
355}
356
357
361
362TrabBoneNL3DStatus :: TrabBoneNL3DStatus(GaussPoint *g) :
364{
365}
366
367
368void
369TrabBoneNL3DStatus :: printOutputAt(FILE *file, TimeStep *tStep) const
370{
371 StructuralMaterialStatus :: printOutputAt(file, tStep);
372 fprintf(file, "status {");
373 fprintf(file, " plastrains ");
374 for ( auto &val : plasDef ) {
375 fprintf( file, " %.4e", val );
376 }
377
378 fprintf(file, " ,");
379 fprintf(file, " kappa %.4e ,", kappa);
380 fprintf(file, " dam %.4e ,", tempDam);
381 fprintf(file, " esed %.4e ,", this->tempTSED - this->tempPSED);
382 fprintf(file, " psed %.4e ,", this->tempPSED);
383 fprintf(file, " tsed %.4e", this->tempTSED);
384 fprintf(file, "}\n");
385}
386
387
388void
389TrabBoneNL3DStatus :: initTempStatus()
390{
391 TrabBone3DStatus :: initTempStatus();
392}
393
394
395void
396TrabBoneNL3DStatus :: updateYourself(TimeStep *tStep)
397{
398 TrabBone3DStatus :: updateYourself(tStep);
399}
400
401
402Interface *
403TrabBoneNL3DStatus :: giveInterface(InterfaceType type)
404{
406 return this;
407 } else {
408 return nullptr;
409 }
410}
411
412
413int
414TrabBoneNL3D :: packUnknowns(DataStream &buff, TimeStep *tStep, GaussPoint *ip)
415{
416 abort();
417 #if 0
418 IDNLMaterialStatus *nlStatus = static_cast< IDNLMaterialStatus * >( this->giveStatus(ip) );
419
420 this->buildNonlocalPointTable(ip);
422
423 return buff.write( nlStatus->giveLocalEquivalentStrainForAverage() );
424
425 #endif
426}
427
428int
429TrabBoneNL3D :: unpackAndUpdateUnknowns(DataStream &buff, TimeStep *tStep, GaussPoint *ip)
430{
431 abort();
432 #if 0
433 int result;
434 IDNLMaterialStatus *nlStatus = static_cast< IDNLMaterialStatus * >( this->giveStatus(ip) );
435 double localEquivalentStrainForAverage;
436
437 result = buff.read(localEquivalentStrainForAverage);
438 nlStatus->setLocalEquivalentStrainForAverage(localEquivalentStrainForAverage);
439 return result;
440
441 #endif
442}
443
444int
445TrabBoneNL3D :: estimatePackSize(DataStream &buff, GaussPoint *ip)
446{
447 abort();
448 #if 0
449 // Note: nlStatus localStrainVectorForAverage memeber must be properly sized!
450 // IDNLMaterialStatus *nlStatus = (IDNLMaterialStatus*) this -> giveStatus (ip);
451 return buff.givePackSizeOfDouble(1);
452
453 #endif
454}
455
456//
457// END: PARALLEL MODE OPTION
459} // end namespace oofem
#define REGISTER_Material(class)
virtual int read(int *data, std::size_t count)=0
Reads count integer values into array pointed by data.
virtual int givePackSizeOfDouble(std::size_t count)=0
virtual int write(const int *data, std::size_t count)=0
Writes count integer values from array pointed by data.
void setField(int item, InputFieldType id)
virtual double computeVolumeAround(GaussPoint *gp)
Definition element.h:530
void resize(Index s)
Definition floatarray.C:94
double & at(Index i)
Definition floatarray.h:202
void zero()
Zeroes all coefficients of receiver.
Definition floatarray.C:683
void beTProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Definition floatarray.C:721
*Sets size of receiver to be an empty matrix It will have zero rows and zero columns size void clear()
int giveNumberOfColumns() const
Returns number of columns of receiver.
void plusDyadUnsym(const FloatArray &a, const FloatArray &b, double dV)
double at(std::size_t i, std::size_t j) const
Element * giveElement()
Returns corresponding element to receiver.
Definition gausspoint.h:187
double giveLocalEquivalentStrainForAverage()
Returns the local equivalent strain to be averaged.
Definition idmnl1.h:74
void setLocalEquivalentStrainForAverage(double ls)
Sets the localEquivalentStrainForAverage to given value.
Definition idmnl1.h:76
virtual MaterialStatus * giveStatus(GaussPoint *gp) const
Definition material.C:206
virtual void initTempStatus(GaussPoint *gp) const
Definition material.C:221
void buildNonlocalPointTable(GaussPoint *gp) const
void updateDomainBeforeNonlocAverage(TimeStep *tStep) const
std ::vector< localIntegrationRecord > * giveIntegrationDomainList()
virtual int assemble(const IntArray &loc, const FloatMatrix &mat)=0
void giveStressDependentPartOfStrainVector(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrainVector, TimeStep *tStep, ValueModeType mode) const
FloatArrayF< 6 > plasDef
Definition trabbone3d.h:111
TrabBone3DStatus(GaussPoint *g)
Definition trabbone3d.C:979
FloatMatrixF< 6, 6 > constructAnisoComplTensor() const
Construct anisotropic compliance tensor.
Definition trabbone3d.C:580
void computePlasStrainEnerDensity(GaussPoint *gp, const FloatArrayF< 6 > &strain, const FloatArrayF< 6 > &stress) const
Definition trabbone3d.C:54
TrabBone3D(int n, Domain *d)
Definition trabbone3d.C:50
void performPlasticityReturn(GaussPoint *gp, const FloatArrayF< 6 > &strain, TimeStep *tStep) const
Definition trabbone3d.C:175
double computeDamage(GaussPoint *gp, TimeStep *tStep) const
Definition trabbone3d.C:510
FloatArrayF< 6 > computeDensificationStress(GaussPoint *gp, const FloatArrayF< 6 > &totalStrain, TimeStep *tStep) const
Definition trabbone3d.C:530
double giveLocalCumPlastStrainForAverage() const
double computeCumPlastStrain(GaussPoint *gp, TimeStep *tStep) const override
TrabBoneNL3D(int n, Domain *d)
int giveLocalNonlocalStiffnessContribution(GaussPoint *gp, IntArray &loc, const UnknownNumberingScheme &s, FloatArray &lcontrib, TimeStep *tStep)
double computeLocalCumPlastStrain(const FloatArrayF< 6 > &strain, GaussPoint *gp, TimeStep *tStep) const
#define IR_GIVE_OPTIONAL_FIELD(__ir, __value, __id)
Definition inputrecord.h:75
#define IR_GIVE_FIELD(__ir, __value, __id)
Definition inputrecord.h:67
FloatMatrixF< N, P > Tdot(const FloatMatrixF< M, N > &a, const FloatMatrixF< M, P > &b)
Computes .
FloatMatrixF< N, M > dyad(const FloatArrayF< N > &a, const FloatArrayF< M > &b)
Computes the dyadic product .
double sum(const FloatArray &x)
double dot(const FloatArray &x, const FloatArray &y)
FloatMatrixF< N, N > inv(const FloatMatrixF< N, N > &mat, double zeropiv=1e-24)
Computes the inverse.
@ NonlocalMaterialStiffnessInterfaceType
@ NonlocalMaterialStatusExtensionInterfaceType
@ NonlocalMaterialExtensionInterfaceType
double distance(const FloatArray &x, const FloatArray &y)
#define _IFT_TrabBoneNL3D_m
#define _IFT_TrabBoneNL3D_r

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