OOFEM 3.0
Loading...
Searching...
No Matches
simplecrosssection.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
39#include "gausspoint.h"
40#include "floatarray.h"
41#include "floatarrayf.h"
42#include "floatmatrixf.h"
43#include "classfactory.h"
44#include "dynamicinputrecord.h"
45#include "datastream.h"
46#include "contextioerr.h"
47#include "engngm.h"
48
49namespace oofem {
51
52
55{
56 auto mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
57 return mat->giveRealStressVector_3d(strain, gp, tStep);
58}
59
62{
63 auto mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
64 IntArray strainControl = {
65 1, 2, 4, 5, 6
66 };
67 return mat->giveRealStressVector_ShellStressControl(strain, strainControl, gp, tStep);
68}
69
70
73{
74 auto mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
75 return mat->giveRealStressVector_PlaneStrain(strain, gp, tStep);
76}
77
78
81{
82 auto mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
83 return mat->giveRealStressVector_PlaneStress(strain, gp, tStep);
84}
85
86
89{
90 auto mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
91 return mat->giveRealStressVector_1d(strain, gp, tStep);
92}
93
94
97{
98 auto mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
99 return mat->giveRealStressVector_Warping(strain, gp, tStep);
100}
101
102
104SimpleCrossSection::giveStiffnessMatrix_3d(MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) const
105{
106 auto mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
107 return mat->give3dMaterialStiffnessMatrix(rMode, gp, tStep);
108}
109
110
113{
114 auto mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
115 return mat->givePlaneStressStiffMtrx(rMode, gp, tStep);
116}
117
118
121{
122 auto mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
123 return mat->givePlaneStrainStiffMtrx(rMode, gp, tStep);
124}
125
126
128SimpleCrossSection::giveStiffnessMatrix_1d(MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) const
129{
130 auto mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
131 return mat->give1dStressStiffMtrx(rMode, gp, tStep);
132}
133
134
137{
145 auto mat = static_cast< StructuralMaterial * >( this->giveMaterial(gp) );
146 auto elasticStrain = strain;
147 FloatArray et;
148 this->giveTemperatureVector(et, gp, tStep);
149 if ( et.giveSize() > 0 ) {
150 auto e0 = mat->giveThermalDilatationVector(gp, tStep);
151 double thick = this->give(CS_Thickness, gp);
152 elasticStrain.at(1) -= e0.at(1) * ( et.at(1) - mat->giveReferenceTemperature() );
153 if ( et.giveSize() > 1 ) {
154 elasticStrain.at(2) -= e0.at(1) * et.at(2) / thick; // kappa_x
155 }
156 }
157 auto tangent = this->give2dBeamStiffMtrx(ElasticStiffness, gp, tStep);
158 auto answer = dot(tangent, elasticStrain);
159
160 auto status = static_cast< StructuralMaterialStatus * >( mat->giveStatus(gp) );
161 status->letTempStrainVectorBe(strain);
162 status->letTempStressVectorBe(answer);
163
164 return answer;
165}
166
167
170{
176 auto mat = static_cast< StructuralMaterial * >( this->giveMaterial(gp) );
177 auto elasticStrain = strain;
178 FloatArray et;
179 this->giveTemperatureVector(et, gp, tStep);
180 if ( et.giveSize() > 0 ) {
181 double thick = this->give(CS_Thickness, gp);
182 double width = this->give(CS_Width, gp);
183 auto e0 = mat->giveThermalDilatationVector(gp, tStep);
184 elasticStrain.at(1) -= e0.at(1) * ( et.at(1) - mat->giveReferenceTemperature() );
185 if ( et.giveSize() > 1 ) {
186 elasticStrain.at(5) -= e0.at(1) * et.at(2) / thick; // kappa_y
187 if ( et.giveSize() > 2 ) {
188 elasticStrain.at(6) -= e0.at(1) * et.at(3) / width; // kappa_z
189 }
190 }
191 }
192 auto tangent = this->give3dBeamStiffMtrx(ElasticStiffness, gp, tStep);
193 auto answer = dot(tangent, elasticStrain);
194
195 auto status = static_cast< StructuralMaterialStatus * >( mat->giveStatus(gp) );
196 status->letTempStrainVectorBe(strain);
197 status->letTempStressVectorBe(answer);
198
199 return answer;
200}
201
202
205{
212 auto mat = static_cast< StructuralMaterial * >( this->giveMaterial(gp) );
213 auto elasticStrain = strain;
214 FloatArray et;
215 this->giveTemperatureVector(et, gp, tStep); // FIXME use return channel, maybe fixed size/std::pair?
216 if ( et.giveSize() > 0 ) {
217 double thick = this->give(CS_Thickness, gp);
218 auto e0 = mat->giveThermalDilatationVector(gp, tStep);
219 if ( et.giveSize() > 1 ) {
220 elasticStrain.at(1) -= e0.at(1) * et.at(2) / thick; // kappa_x
221 elasticStrain.at(2) -= e0.at(2) * et.at(2) / thick; // kappa_y
222 }
223 }
224 auto tangent = this->give2dPlateStiffMtrx(ElasticStiffness, gp, tStep);
225 auto answer = dot(tangent, elasticStrain);
226
227 auto status = static_cast< StructuralMaterialStatus * >( mat->giveStatus(gp) );
228 status->letTempStrainVectorBe(strain);
229 status->letTempStressVectorBe(answer);
230
231 return answer;
232}
233
234
237{
244 auto mat = static_cast< StructuralMaterial * >( this->giveMaterial(gp) );
245 auto elasticStrain = strain;
246 FloatArray et;
247 this->giveTemperatureVector(et, gp, tStep);
248 if ( et.giveSize() ) {
249 double thick = this->give(CS_Thickness, gp);
250 auto e0 = mat->giveThermalDilatationVector(gp, tStep);
251 elasticStrain.at(1) -= e0.at(1) * ( et.at(1) - mat->giveReferenceTemperature() );
252 elasticStrain.at(2) -= e0.at(2) * ( et.at(1) - mat->giveReferenceTemperature() );
253 if ( et.giveSize() > 1 ) {
254 elasticStrain.at(4) -= e0.at(1) * et.at(2) / thick; // kappa_x
255 elasticStrain.at(5) -= e0.at(2) * et.at(2) / thick; // kappa_y
256 }
257 }
258 auto tangent = this->give3dShellStiffMtrx(ElasticStiffness, gp, tStep);
259 auto answer = dot(tangent, elasticStrain);
260
261 auto status = static_cast< StructuralMaterialStatus * >( mat->giveStatus(gp) );
262 status->letTempStrainVectorBe(strain);
263 status->letTempStressVectorBe(answer);
264
265 return answer;
266}
267
268
271{
278
279 FloatArrayF< 9 >answer;
280 FloatArrayF< 8 >rstrain;
281 for ( int i = 1; i <= 8; i++ ) {
282 rstrain.at(i) = strain.at(i);
283 }
284 FloatArray ra = this->giveGeneralizedStress_Shell(rstrain, gp, tStep);
285 for ( int i = 1; i <= 8; i++ ) {
286 answer.at(i) = ra.at(i);
287 }
288 answer.at(9) = this->give(CS_DrillingStiffness, gp) * strain.at(9);
289
290 auto mat = static_cast< StructuralMaterial * >( this->giveMaterial(gp) );
291 auto status = static_cast< StructuralMaterialStatus * >( mat->giveStatus(gp) );
292 status->letTempStrainVectorBe(strain);
293 status->letTempStressVectorBe(answer);
294
295 return answer;
296}
297
298
301{
302 auto tangent = this->giveMembraneRotStiffMtrx(ElasticStiffness, gp, tStep);
303 auto answer = dot(tangent, strain);
304
305 auto status = static_cast< StructuralMaterialStatus * >( this->giveMaterial(gp)->giveStatus(gp) );
306 status->letTempStrainVectorBe(strain);
307 status->letTempStressVectorBe(answer);
308
311 return answer;
312}
313
316{
317 auto mat = static_cast< StructuralMaterial * >( this->giveMaterial(gp) );
318 return mat->giveRealStressVector_2dPlateSubSoil(generalizedStrain, gp, tStep);
319}
320
321
322void
324{
325 MaterialMode mode = gp->giveMaterialMode();
326 if ( mode == _2dBeam ) {
327 answer = this->give2dBeamStiffMtrx(rMode, gp, tStep);
328 } else if ( mode == _3dBeam ) {
329 answer = this->give3dBeamStiffMtrx(rMode, gp, tStep);
330 } else if ( mode == _2dPlate ) {
331 answer = this->give2dPlateStiffMtrx(rMode, gp, tStep);
332 } else if ( mode == _3dShell ) {
333 answer = this->give3dShellStiffMtrx(rMode, gp, tStep);
334 } else if ( mode == _3dDegeneratedShell ) {
335 answer = this->give3dDegeneratedShellStiffMtrx(rMode, gp, tStep);
336 } else {
337 auto mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
338
339 if ( mode == _3dMat ) {
340 answer = mat->give3dMaterialStiffnessMatrix(rMode, gp, tStep);
341 } else if ( mode == _PlaneStress ) {
342 answer = mat->givePlaneStressStiffMtrx(rMode, gp, tStep);
343 } else if ( mode == _PlaneStrain ) {
344 answer = mat->givePlaneStrainStiffMtrx(rMode, gp, tStep);
345 } else if ( mode == _1dMat ) {
346 answer = mat->give1dStressStiffMtrx(rMode, gp, tStep);
347 } else {
348 mat->giveStiffnessMatrix(answer, rMode, gp, tStep);
349 }
350 }
351}
352
353
355SimpleCrossSection::give2dBeamStiffMtrx(MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) const
356{
357 auto mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
358
359 auto mat1d = mat->give1dStressStiffMtrx(rMode, gp, tStep);
360 auto area = this->give(CS_Area, gp);
361 auto Iy = this->give(CS_InertiaMomentY, gp);
362 auto shearAreaz = this->give(CS_ShearAreaZ, gp);
363
365 answer.at(1, 1) = mat1d.at(1, 1) * area;
366 answer.at(2, 2) = mat1d.at(1, 1) * Iy;
367 answer.at(3, 3) = shearAreaz * mat->give('G', gp);
368 return answer;
369}
370
371
373SimpleCrossSection::give3dBeamStiffMtrx(MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) const
374{
375 auto mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
376
377 auto mat1d = mat->give1dStressStiffMtrx(rMode, gp, tStep);
378 double E = mat1d.at(1, 1);
379 double G = mat->give('G', gp);
380 double area = this->give(CS_Area, gp);
381 double Iy = this->give(CS_InertiaMomentY, gp);
382 double Iz = this->give(CS_InertiaMomentZ, gp);
383 double Ik = this->give(CS_TorsionConstantX, gp); // St. Venant torsional constant
384
385 //shearCoeff = this->give(CS_BeamShearCoeff);
386 double shearAreay = this->give(CS_ShearAreaY, gp);
387 double shearAreaz = this->give(CS_ShearAreaZ, gp);
388
390
391 answer.at(1, 1) = E * area;
393 answer.at(2, 2) = shearAreay * G;
394 answer.at(3, 3) = shearAreaz * G;
395 //answer.at(2, 2) = shearCoeff * G * area;
396 //answer.at(3, 3) = shearCoeff * G * area;
397 answer.at(4, 4) = G * Ik;
398 answer.at(5, 5) = E * Iy;
399 answer.at(6, 6) = E * Iz;
400 return answer;
401}
402
403
405SimpleCrossSection::give2dPlateStiffMtrx(MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) const
406{
407 auto mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
408
409 auto mat2d = mat->givePlaneStressStiffMtrx(rMode, gp, tStep);
410 double thickness = this->give(CS_Thickness, gp);
411 double thickness3 = thickness * thickness * thickness;
412
414
415 for ( int i = 1; i <= 2; i++ ) {
416 for ( int j = 1; j <= 2; j++ ) {
417 answer.at(i, j) = mat2d.at(i, j) * thickness3 / 12.;
418 }
419 }
420
421 answer.at(3, 3) = mat2d.at(3, 3) * thickness3 / 12.;
422 answer.at(4, 4) = mat2d.at(3, 3) * thickness * ( 5. / 6. );
423 answer.at(5, 5) = answer.at(4, 4);
424 return answer;
425}
426
427
429SimpleCrossSection::give3dShellStiffMtrx(MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) const
430{
431 auto mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
432
433 double thickness = this->give(CS_Thickness, gp);
434 double thickness3 = thickness * thickness * thickness;
435
436 auto mat2d = mat->givePlaneStressStiffMtrx(rMode, gp, tStep);
437
439
440 for ( int i = 1; i <= 3; i++ ) {
441 for ( int j = 1; j <= 3; j++ ) {
442 answer.at(i, j) = mat2d.at(i, j) * thickness;
443 }
444 }
445 for ( int i = 1; i <= 3; i++ ) {
446 for ( int j = 1; j <= 3; j++ ) {
447 answer.at(i + 3, j + 3) = mat2d.at(i, j) * thickness3 / 12.0;
448 }
449 }
450
451 answer.at(8, 8) = answer.at(7, 7) = mat2d.at(3, 3) * thickness * ( 5. / 6. );
452 return answer;
453}
454
455
457SimpleCrossSection::give3dShellRotStiffMtrx(MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) const
458{
459 FloatMatrix d, answer;
460 d = this->give3dShellStiffMtrx(rMode, gp, tStep);
461 answer.resize(9, 9);
462 answer.zero();
463 answer.assemble(d, { 1, 2, 3, 4, 5, 6, 7, 8 });
464 answer.at(9, 9) = this->give(CS_DrillingStiffness, gp);
465 return answer;
466}
467
468
471{
472 auto mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
473 auto answer = mat->give3dMaterialStiffnessMatrix(rMode, gp, tStep);
474
475 answer.at(1, 1) -= answer.at(1, 3) * answer.at(3, 1) / answer.at(3, 3);
476 answer.at(2, 1) -= answer.at(2, 3) * answer.at(3, 1) / answer.at(3, 3);
477 answer.at(1, 2) -= answer.at(1, 3) * answer.at(3, 2) / answer.at(3, 3);
478 answer.at(2, 2) -= answer.at(2, 3) * answer.at(3, 2) / answer.at(3, 3);
479
480 answer.at(3, 1) = 0.0;
481 answer.at(3, 2) = 0.0;
482 answer.at(3, 3) = 0.0;
483 answer.at(2, 3) = 0.0;
484 answer.at(1, 3) = 0.0;
485 return answer;
486}
487
489SimpleCrossSection::giveMembraneRotStiffMtrx(MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) const
490{
491 auto mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
492 auto d = mat->givePlaneStressStiffMtrx(ElasticStiffness, gp, tStep);
493 auto ds = assemble< 4, 4 >(d, { 0, 1, 2 }, { 0, 1, 2 });
494 ds.at(4, 4) = this->give(CS_DrillingStiffness, gp);
495 return ds;
496}
497
500{
501 auto mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
502 return mat->give2dPlateSubSoilStiffMtrx(ElasticStiffness, gp, tStep);
503}
504
505
506void
508{
510
511 double value;
512
513 double thick = 0.0;
517 }
518
519 double width = 0.0;
522 propertyDictionary.add(CS_Width, width);
523 }
524
525 double area = 0.0;
528 } else {
529 area = thick * width;
530 }
531 propertyDictionary.add(CS_Area, area);
532
533 value = 0.0;
536
537 value = 0.0;
540
541 value = 0.0;
544
545 double beamshearcoeff = 0.0;
547 propertyDictionary.add(CS_BeamShearCoeff, beamshearcoeff);
548
549 value = 0.0;
551 if ( value == 0.0 ) {
552 value = beamshearcoeff * area;
553 }
555
556 value = 0.0;
558 if ( value == 0.0 ) {
559 value = beamshearcoeff * area;
560 }
562
563 value = 0.0;
566
567 value = 0.0;
570
571 value = 0.0;
574
575 this->materialNumber = 0;
577
579 value = 0.0;
582
583 value = 0.0;
586
587 value = 1.0;
590 }
591}
592
593
594void
596{
598
599 if ( this->propertyDictionary.includes(CS_Thickness) ) {
601 }
602
603 if ( this->propertyDictionary.includes(CS_Width) ) {
605 }
606
607 if ( this->propertyDictionary.includes(CS_Area) ) {
609 }
610
611 if ( this->propertyDictionary.includes(CS_TorsionConstantX) ) {
613 }
614
615 if ( this->propertyDictionary.includes(CS_InertiaMomentY) ) {
617 }
618
619 if ( this->propertyDictionary.includes(CS_InertiaMomentZ) ) {
621 }
622
623 if ( this->propertyDictionary.includes(CS_ShearAreaY) ) {
625 }
626
627 if ( this->propertyDictionary.includes(CS_ShearAreaZ) ) {
629 }
630
631 if ( this->propertyDictionary.includes(CS_BeamShearCoeff) ) {
633 }
634
636
637 if ( this->propertyDictionary.includes(CS_DirectorVectorX) ) {
639 }
640
641 if ( this->propertyDictionary.includes(CS_DirectorVectorY) ) {
643 }
644
645 if ( this->propertyDictionary.includes(CS_DirectorVectorZ) ) {
647 }
648}
649
650void
656
657
658bool
660{
661 if ( this->giveMaterialNumber() ) {
662 return this->domain->giveMaterial(this->giveMaterialNumber() )->isCharacteristicMtrxSymmetric(rMode);
663 } else {
664 return false; // Bet false...
665 }
666}
667
668
669Material *
671{
672 if ( this->giveMaterialNumber() ) {
673 return this->giveDomain()->giveMaterial(this->giveMaterialNumber() );
674 } else {
675 return ip->giveElement()->giveMaterial();
676 }
677}
678
679
680double
681SimpleCrossSection::give(int aProperty, GaussPoint *gp) const
682{
683 return this->giveMaterial(gp)->give(aProperty, gp);
684}
685
686
687int
689{
690 if ( type == IST_CrossSectionNumber ) {
691 answer.resize(1);
692 answer.at(1) = this->giveNumber();
693 return 1;
694 }
695 return this->giveMaterial(ip)->giveIPValue(answer, ip, type, tStep);
696}
697
698
699
700int
702//
703// check internal consistency
704// mainly tests, whether material and crossSection data
705// are safe for conversion to "Structural" versions
706//
707{
708 int result = 1;
709 Material *mat = this->giveDomain()->giveMaterial(this->materialNumber);
710 if ( !dynamic_cast< StructuralMaterial * >( mat ) ) {
711 OOFEM_WARNING("material %s without structural support", mat->giveClassName() );
712 result = 0;
713 }
714
715 return result;
716}
717
718
724
725
726
727
728
729void
731{
732 // This function returns the Cauchy stress in vector format
733 // corresponding to a given deformation gradient according to the stress-deformation
734 // mode stored in the each gp.
735
736 MaterialMode mode = gp->giveMaterialMode();
737 StructuralMaterial *mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
738
739 if ( mode == _3dMat ) {
740 mat->giveCauchyStressVector_3d(answer, gp, reducedvF, tStep);
741 } else if ( mode == _PlaneStrain ) {
742 mat->giveCauchyStressVector_PlaneStrain(answer, gp, reducedvF, tStep);
743 } else if ( mode == _PlaneStress ) {
744 mat->giveCauchyStressVector_PlaneStress(answer, gp, reducedvF, tStep);
745 } else if ( mode == _1dMat ) {
746 mat->giveCauchyStressVector_1d(answer, gp, reducedvF, tStep);
747 }
748}
749
750void
752{
753 MaterialMode mode = gp->giveMaterialMode();
754 StructuralMaterial *mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
755
756 if ( mode == _3dMat ) {
757 // mat->giveCauchyStressVector_3d(answer, gp, reducedvF, tStep);
758 } else if ( mode == _PlaneStrain ) {
759 mat->giveEshelbyStressVector_PlaneStrain(answer, gp, reducedvF, tStep);
760 } else if ( mode == _PlaneStress ) {
761 // mat->giveCauchyStressVector_PlaneStress(answer, gp, reducedvF, tStep);
762 } else if ( mode == _1dMat ) {
763 // mat->giveCauchyStressVector_1d(answer, gp, reducedvF, tStep);
764 }
765}
766
767
768
769
772{
773 auto mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
774 return mat->giveFirstPKStressVector_3d(vF, gp, tStep);
775}
776
777
780{
781 auto mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
782 return mat->giveFirstPKStressVector_PlaneStrain(vF, gp, tStep);
783}
784
785
788{
789 auto mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
790 return mat->giveFirstPKStressVector_PlaneStress(strain, gp, tStep);
791}
792
793
796{
797 auto mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
798 return mat->giveFirstPKStressVector_1d(vF, gp, tStep);
799}
800
801
802
803void
805{
806 MaterialMode mode = gp->giveMaterialMode();
807 auto mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
808 if ( mode == _3dMat ) {
809 answer = mat->give3dMaterialStiffnessMatrix_dPdF(rMode, gp, tStep);
810 } else if ( mode == _PlaneStress ) {
811 answer = mat->givePlaneStressStiffnessMatrix_dPdF(rMode, gp, tStep);
812 } else if ( mode == _PlaneStrain ) {
813 answer = mat->givePlaneStrainStiffnessMatrix_dPdF(rMode, gp, tStep);
814 } else if ( mode == _1dMat ) {
815 answer = mat->give1dStressStiffnessMatrix_dPdF(rMode, gp, tStep);
816 } else {
817 OOFEM_ERROR("Unsupported material mode for large strain analysis");
818 }
819}
820
821
824{
825 auto mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
826 return mat->give3dMaterialStiffnessMatrix_dPdF(rMode, gp, tStep);
827}
828
829
832{
833 auto mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
834 return mat->givePlaneStrainStiffnessMatrix_dPdF(rMode, gp, tStep);
835}
836
837
838
841{
842 auto mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
843 return mat->givePlaneStressStiffnessMatrix_dPdF(rMode, gp, tStep);
844}
845
846
849{
850 auto mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
851 return mat->give1dStressStiffnessMatrix_dPdF(rMode, gp, tStep);
852}
853
854
855
856
857void
859 MatResponseMode rMode, GaussPoint *gp,
860 TimeStep *tStep)
861{
862 StructuralMaterial *mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
863
864 MaterialMode mode = gp->giveMaterialMode();
865 if ( mode == _3dMat ) {
866 mat->give3dMaterialStiffnessMatrix_dCde(answer, rMode, gp, tStep);
867 } else if ( mode == _PlaneStress ) {
868 mat->givePlaneStressStiffMtrx_dCde(answer, rMode, gp, tStep);
869 } else if ( mode == _PlaneStrain ) {
870 mat->givePlaneStrainStiffMtrx_dCde(answer, rMode, gp, tStep);
871 } else if ( mode == _1dMat ) {
872 mat->give1dStressStiffMtrx_dCde(answer, rMode, gp, tStep);
873 } else {
874 OOFEM_ERROR("unknown mode (%s)", __MaterialModeToString(mode) );
875 }
876}
877
878
879void
881{
882 Element *elem = gp->giveElement();
883 answer.clear();
884 //sum up all prescribed temperatures over an element
885 StructuralElement *selem = dynamic_cast< StructuralElement * >( elem );
886 selem->computeResultingIPTemperatureAt(answer, tStep, gp, VM_Total);
887
888 /* add external source, if provided */
889 FieldManager *fm = this->domain->giveEngngModel()->giveContext()->giveFieldManager();
890 FieldPtr tf;
891
892 if ( ( tf = fm->giveField(FT_Temperature) ) ) {
893 // temperature field registered
894 FloatArray gcoords, et2;
895 int err;
897 if ( ( err = tf->evaluateAt(et2, gcoords, VM_Total, tStep) ) ) {
898 OOFEM_ERROR("tf->evaluateAt failed, element %d, error code %d", elem->giveNumber(), err);
899 }
900 if ( et2.isNotEmpty() ) {
901 if ( answer.isEmpty() ) {
902 answer = et2;
903 } else {
904 answer.at(1) += et2.at(1);
905 }
906 }
907 }
908}
909
910int
912{
913 return this->giveMaterial(gp)->packUnknowns(buff, tStep, gp);
914}
915
916int
918{
919 return this->giveMaterial(gp)->unpackAndUpdateUnknowns(buff, tStep, gp);
920}
921
922int
924{
925 return this->giveMaterial(gp)->estimatePackSize(buff, gp);
926}
927
928void
930{
932
933 if ( mode & CM_Definition ) {
934 if ( !stream.write(materialNumber) ) {
936 }
937 if ( !stream.write(czMaterialNumber) ) {
939 }
940 }
941}
942
943void
945{
947
948 if ( mode & CM_Definition ) {
949 if ( !stream.read(materialNumber) ) {
951 }
952
953 if ( !stream.read(czMaterialNumber) ) {
955 }
956 }
957}
958} // end namespace oofem
#define E(a, b)
#define REGISTER_CrossSection(class)
void giveInputRecord(DynamicInputRecord &input) override
Dictionary propertyDictionary
void restoreContext(DataStream &stream, ContextMode mode) override
void saveContext(DataStream &stream, ContextMode mode) override
void initializeFrom(InputRecord &ir) override
virtual int read(int *data, std::size_t count)=0
Reads count integer values into array pointed by data.
virtual int write(const int *data, std::size_t count)=0
Writes count integer values from array pointed by data.
Material * giveMaterial(int n)
Definition domain.C:284
void setField(int item, InputFieldType id)
virtual int computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords)
Definition element.C:1225
virtual Material * giveMaterial()
Definition element.C:523
Domain * giveDomain() const
Definition femcmpnn.h:97
virtual Interface * giveInterface(InterfaceType t)
Definition femcmpnn.h:181
Domain * domain
Link to domain object, useful for communicating with other FEM components.
Definition femcmpnn.h:79
virtual const char * giveClassName() const =0
int giveNumber() const
Definition femcmpnn.h:104
FieldPtr giveField(FieldType key)
double & at(std::size_t i)
void resize(Index s)
Definition floatarray.C:94
double & at(Index i)
Definition floatarray.h:202
Index giveSize() const
Returns the size of receiver.
Definition floatarray.h:261
bool isEmpty() const
Returns true if receiver is empty.
Definition floatarray.h:265
bool isNotEmpty() const
Returns true if receiver is not empty.
Definition floatarray.h:263
double at(std::size_t i, std::size_t j) const
void resize(Index rows, Index cols)
Definition floatmatrix.C:79
void zero()
Zeroes all coefficient of receiver.
void assemble(const FloatMatrix &src, const IntArray &loc)
double at(std::size_t i, std::size_t j) const
IntegrationPointStatus * setMaterialStatus(std::unique_ptr< IntegrationPointStatus > ptr, IntegrationPointStatusIDType key=IPSID_Default)
Definition gausspoint.h:213
const FloatArray & giveNaturalCoordinates() const
Returns coordinate array of receiver.
Definition gausspoint.h:138
MaterialMode giveMaterialMode()
Returns corresponding material mode of receiver.
Definition gausspoint.h:190
Element * giveElement()
Returns corresponding element to receiver.
Definition gausspoint.h:187
virtual bool hasField(InputFieldType id)=0
Returns true if record contains field identified by idString keyword.
virtual std::unique_ptr< MaterialStatus > CreateStatus(GaussPoint *gp) const
Definition material.h:332
virtual int estimatePackSize(DataStream &buff, GaussPoint *ip)
Definition material.h:315
virtual MaterialStatus * giveStatus(GaussPoint *gp) const
Definition material.C:206
virtual int unpackAndUpdateUnknowns(DataStream &buff, TimeStep *tStep, GaussPoint *ip)
Definition material.h:311
virtual double give(int aProperty, GaussPoint *gp) const
Definition material.C:51
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Definition material.C:138
virtual int packUnknowns(DataStream &buff, TimeStep *tStep, GaussPoint *ip)
Definition material.h:302
Interface * giveMaterialInterface(InterfaceType t, IntegrationPoint *ip) override
void giveCharMaterialStiffnessMatrix_dPdF(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) override
FloatMatrixF< 4, 4 > giveStiffnessMatrix_dPdF_PlaneStress(MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) const override
int estimatePackSize(DataStream &buff, GaussPoint *gp) override
void saveContext(DataStream &stream, ContextMode mode) override
int materialNumber
Material number.
virtual FloatArrayF< 9 > giveFirstPKStress_3d(const FloatArrayF< 9 > &reducedF, GaussPoint *gp, TimeStep *tStep) const override
FloatArrayF< 6 > giveRealStress_3d(const FloatArrayF< 6 > &reducedStrain, GaussPoint *gp, TimeStep *tStep) const override
FloatMatrixF< 6, 6 > giveStiffnessMatrix_3d(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override
int czMaterialNumber
Cohesive zone material number.
FloatMatrixF< 4, 4 > giveStiffnessMatrix_PlaneStrain(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override
Material * giveMaterial(IntegrationPoint *ip) const override
hidden by virtual oofem::Material* TransportCrossSection::giveMaterial() const
void giveInputRecord(DynamicInputRecord &input) override
int packUnknowns(DataStream &buff, TimeStep *tStep, GaussPoint *gp) override
void initializeFrom(InputRecord &ir) override
FloatArrayF< 6 > giveGeneralizedStress_Beam3d(const FloatArrayF< 6 > &generalizedStrain, GaussPoint *gp, TimeStep *tStep) const override
virtual FloatArrayF< 1 > giveFirstPKStress_1d(const FloatArrayF< 1 > &reducedF, GaussPoint *gp, TimeStep *tStep) const override
int giveIPValue(FloatArray &answer, GaussPoint *ip, InternalStateType type, TimeStep *tStep) override
FloatArrayF< 8 > giveGeneralizedStress_Shell(const FloatArrayF< 8 > &generalizedStrain, GaussPoint *gp, TimeStep *tStep) const override
FloatArrayF< 5 > giveGeneralizedStress_Plate(const FloatArrayF< 5 > &generalizedStrain, GaussPoint *gp, TimeStep *tStep) const override
FloatMatrixF< 9, 9 > give3dShellRotStiffMtrx(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override
FloatArrayF< 6 > giveRealStress_3dDegeneratedShell(const FloatArrayF< 6 > &reducedStrain, GaussPoint *gp, TimeStep *tStep) const override
FloatMatrixF< 5, 5 > give2dPlateStiffMtrx(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override
FloatMatrixF< 5, 5 > giveStiffnessMatrix_dPdF_PlaneStrain(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override
FloatMatrixF< 9, 9 > giveStiffnessMatrix_dPdF_3d(MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) const override
virtual FloatArrayF< 4 > giveFirstPKStress_PlaneStress(const FloatArrayF< 4 > &reducedF, GaussPoint *gp, TimeStep *tStep) const override
double give(int aProperty, GaussPoint *gp) const override
FloatMatrixF< 1, 1 > giveStiffnessMatrix_dPdF_1d(MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) const override
FloatArrayF< 2 > giveRealStress_Warping(const FloatArrayF< 2 > &reducedStrain, GaussPoint *gp, TimeStep *tStep) const override
void giveTemperatureVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep) const
void giveCharMaterialStiffnessMatrix(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) override
FloatMatrixF< 3, 3 > give2dBeamStiffMtrx(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override
FloatMatrixF< 1, 1 > giveStiffnessMatrix_1d(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override
int unpackAndUpdateUnknowns(DataStream &buff, TimeStep *tStep, GaussPoint *gp) override
FloatMatrixF< 6, 6 > give3dBeamStiffMtrx(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override
void createMaterialStatus(GaussPoint &iGP) override
FloatArrayF< 4 > giveGeneralizedStress_MembraneRot(const FloatArrayF< 4 > &generalizedStrain, GaussPoint *gp, TimeStep *tStep) const override
FloatArrayF< 3 > giveGeneralizedStress_PlateSubSoil(const FloatArrayF< 3 > &generalizedStrain, GaussPoint *gp, TimeStep *tStep) const override
void giveStiffnessMatrix_dCde(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) override
FloatArrayF< 9 > giveGeneralizedStress_ShellRot(const FloatArrayF< 9 > &generalizedStrain, GaussPoint *gp, TimeStep *tStep) const override
FloatArrayF< 3 > giveGeneralizedStress_Beam2d(const FloatArrayF< 3 > &generalizedStrain, GaussPoint *gp, TimeStep *tStep) const override
FloatMatrixF< 3, 3 > giveStiffnessMatrix_PlaneStress(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override
void giveEshelbyStresses(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedvF, TimeStep *tStep) override
FloatArrayF< 1 > giveRealStress_1d(const FloatArrayF< 1 > &reducedStrain, GaussPoint *gp, TimeStep *tStep) const override
FloatMatrixF< 4, 4 > giveMembraneRotStiffMtrx(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override
FloatMatrixF< 8, 8 > give3dShellStiffMtrx(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override
FloatMatrixF< 6, 6 > give3dDegeneratedShellStiffMtrx(MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) const override
FloatArrayF< 4 > giveRealStress_PlaneStrain(const FloatArrayF< 4 > &reducedStrain, GaussPoint *gp, TimeStep *tStep) const override
void restoreContext(DataStream &stream, ContextMode mode) override
bool isCharacteristicMtrxSymmetric(MatResponseMode mode) const override
FloatMatrixF< 3, 3 > give2dPlateSubSoilStiffMtrx(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override
void giveCauchyStresses(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedFIncrement, TimeStep *tStep) override
virtual FloatArrayF< 5 > giveFirstPKStress_PlaneStrain(const FloatArrayF< 5 > &reducedF, GaussPoint *gp, TimeStep *tStep) const override
FloatArrayF< 3 > giveRealStress_PlaneStress(const FloatArrayF< 3 > &reducedStrain, GaussPoint *gp, TimeStep *tStep) const override
virtual void computeResultingIPTemperatureAt(FloatArray &answer, TimeStep *tStep, GaussPoint *gp, ValueModeType mode)
void letTempStrainVectorBe(const FloatArray &v)
Assigns tempStrainVector to given vector v.
virtual void give1dStressStiffMtrx_dCde(FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
virtual FloatMatrixF< 3, 3 > givePlaneStressStiffMtrx(MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const
virtual FloatArrayF< 3 > giveRealStressVector_2dPlateSubSoil(const FloatArrayF< 3 > &reducedE, GaussPoint *gp, TimeStep *tStep) const
Default implementation is not provided.
virtual void giveCauchyStressVector_PlaneStress(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
virtual void giveEshelbyStressVector_PlaneStrain(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
virtual FloatArrayF< 1 > giveRealStressVector_1d(const FloatArrayF< 1 > &reducedE, GaussPoint *gp, TimeStep *tStep) const
Default implementation relies on giveRealStressVector_StressControl.
virtual FloatArrayF< 1 > giveFirstPKStressVector_1d(const FloatArrayF< 1 > &vF, GaussPoint *gp, TimeStep *tStep) const
Default implementation relies on giveFirstPKStressVector_StressControl.
virtual FloatArrayF< 2 > giveRealStressVector_Warping(const FloatArrayF< 2 > &reducedE, GaussPoint *gp, TimeStep *tStep) const
Default implementation relies on giveRealStressVector_StressControl.
virtual FloatMatrixF< 6, 6 > give3dMaterialStiffnessMatrix(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const
virtual FloatMatrixF< 3, 3 > give2dPlateSubSoilStiffMtrx(MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const
virtual FloatMatrixF< 1, 1 > give1dStressStiffnessMatrix_dPdF(MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const
virtual void give3dMaterialStiffnessMatrix_dCde(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
virtual void givePlaneStrainStiffMtrx_dCde(FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
virtual FloatMatrixF< 5, 5 > givePlaneStrainStiffnessMatrix_dPdF(MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const
virtual FloatMatrixF< 4, 4 > givePlaneStrainStiffMtrx(MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const
virtual void giveCauchyStressVector_1d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
virtual FloatArrayF< 4 > giveRealStressVector_PlaneStrain(const FloatArrayF< 4 > &strain, GaussPoint *gp, TimeStep *tStep) const
Default implementation relies on giveRealStressVector_3d.
virtual FloatArrayF< 5 > giveFirstPKStressVector_PlaneStrain(const FloatArrayF< 5 > &vF, GaussPoint *gp, TimeStep *tStep) const
Default implementation relies on giveFirstPKStressVector_3d.
virtual void giveCauchyStressVector_PlaneStrain(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
virtual FloatArrayF< 6 > giveRealStressVector_3d(const FloatArrayF< 6 > &strain, GaussPoint *gp, TimeStep *tStep) const
Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress.
virtual FloatMatrixF< 4, 4 > givePlaneStressStiffnessMatrix_dPdF(MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const
virtual void giveCauchyStressVector_3d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
virtual FloatMatrixF< 1, 1 > give1dStressStiffMtrx(MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const
virtual FloatArrayF< 4 > giveFirstPKStressVector_PlaneStress(const FloatArrayF< 4 > &vF, GaussPoint *gp, TimeStep *tStep) const
Default implementation relies on giveFirstPKStressVector_StressControl.
virtual FloatArrayF< 3 > giveRealStressVector_PlaneStress(const FloatArrayF< 3 > &reducedE, GaussPoint *gp, TimeStep *tStep) const
Default implementation relies on giveRealStressVector_StressControl.
virtual FloatMatrixF< 9, 9 > give3dMaterialStiffnessMatrix_dPdF(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const
virtual FloatArrayF< 9 > giveFirstPKStressVector_3d(const FloatArrayF< 9 > &vF, GaussPoint *gp, TimeStep *tStep) const
Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress.
virtual void givePlaneStressStiffMtrx_dCde(FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
#define THROW_CIOERR(e)
#define CM_Definition
Definition contextmode.h:47
#define OOFEM_WARNING(...)
Definition error.h:80
#define OOFEM_ERROR(...)
Definition error.h:79
#define IR_GIVE_OPTIONAL_FIELD(__ir, __value, __id)
Definition inputrecord.h:75
#define IR_GIVE_FIELD(__ir, __value, __id)
Definition inputrecord.h:67
long ContextMode
Definition contextmode.h:43
@ CS_ShearAreaY
Shear area in y direction.
@ CS_DrillingType
Type of artificially added drilling stiffness for drilling DOFs.
@ CS_InertiaMomentZ
Moment of inertia around z-axis.
@ CS_DirectorVectorY
Director vector component in y-axis.
@ CS_DrillingStiffness
Penalty stiffness for drilling DOFs.
@ CS_TorsionConstantX
Saint-Venant torsional constant (J).
@ CS_RelDrillingStiffness
Relative penalty stiffness for drilling DOFs.
@ CS_BeamShearCoeff
Shear coefficient of beam.
@ CS_DirectorVectorX
Director vector component in x-axis.
@ CS_InertiaMomentY
Moment of inertia around y-axis.
@ CS_Area
Area.
@ CS_Thickness
Thickness.
@ CS_Width
Width.
@ CS_ShearAreaZ
Shear area in z direction.
@ CS_DirectorVectorZ
Director vector component in z-axis.
GaussPoint IntegrationPoint
Definition gausspoint.h:272
FloatArrayF< N > assemble(const FloatArrayF< M > &x, int const (&c)[M])
Assemble components into zero matrix.
double dot(const FloatArray &x, const FloatArray &y)
std::shared_ptr< Field > FieldPtr
Definition field.h:78
@ CIO_IOERR
General IO error.
#define _IFT_SimpleCrossSection_iz
Inertia moment z.
#define _IFT_SimpleCrossSection_thick
#define _IFT_SimpleCrossSection_drillType
Type of artificially added stiffnes for drilling DOFs.
#define _IFT_SimpleCrossSection_directorx
#define _IFT_SimpleCrossSection_shearareaz
Shear area z direction.
#define _IFT_SimpleCrossSection_shearcoeff
#define _IFT_SimpleCrossSection_area
#define _IFT_SimpleCrossSection_relDrillStiffness
Relative penalty term for drilling stiffness.
#define _IFT_SimpleCrossSection_width
#define _IFT_SimpleCrossSection_iy
Inertia moment y.
#define _IFT_SimpleCrossSection_directorz
#define _IFT_SimpleCrossSection_MaterialNumber
Material number for the bulk material.
#define _IFT_SimpleCrossSection_ik
Saint-Venant torsional constant.
#define _IFT_SimpleCrossSection_directory
#define _IFT_SimpleCrossSection_drillStiffness
Penalty term for drilling stiffness.
#define _IFT_SimpleCrossSection_shearareay
Shear area y direction.

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