OOFEM 3.0
Loading...
Searching...
No Matches
structuralelement.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
44#include "feinterpol.h"
45#include "domain.h"
46#include "material.h"
47#include "nonlocalmaterialext.h"
48#include "load.h"
49#include "boundaryload.h"
50#include "pointload.h"
51#include "gausspoint.h"
53#include "intarray.h"
54#include "floatarray.h"
55#include "floatmatrix.h"
57#include "mathfem.h"
60#include "set.h"
61
62#ifdef __OOFEG
63 #include "oofeggraphiccontext.h"
64 #include "connectivitytable.h"
65#endif
66
67
68namespace oofem {
69 StructuralElement :: StructuralElement(int n, Domain *aDomain) :
70 Element(n, aDomain)
71{}
72
73
74StructuralElement :: ~StructuralElement()
75{}
76
77#if 0
78void
79StructuralElement :: computeConstitutiveMatrixAt(FloatMatrix &answer,
80 MatResponseMode rMode, GaussPoint *gp,
81 TimeStep *tStep)
82// Returns the material matrix {E} of the receiver.
83// type of matrix is determined by this->giveMaterialMode()
84// rMode parameter determines type of stiffness matrix to be requested
85// (tangent, secant, ...)
86{
87 this->giveStructuralCrossSection()->giveCharMaterialStiffnessMatrix(answer, rMode, gp, tStep);
88}
89#endif
90
91void StructuralElement :: computeLoadVector(FloatArray &answer, BodyLoad *load, CharType type, ValueModeType mode, TimeStep *tStep)
92{
93 if ( type != ExternalForcesVector ) {
94 answer.clear();
95 return;
96 }
97 // Just a wrapper for the deadweight body load computations:
98 PointLoad *p = dynamic_cast< PointLoad * >( load );
99 if ( p ) {
100 FloatArray lcoords;
101 if ( this->computeLocalCoordinates(lcoords, p->giveCoordinates() ) ) {
102 this->computePointLoadVectorAt(answer, load, tStep, mode);
103 }
104 } else {
106 this->computeBodyLoadVectorAt(answer, load, tStep, mode);
107 }
108}
109
110void StructuralElement :: computeBoundarySurfaceLoadVector(FloatArray &answer, BoundaryLoad *load, int boundary, CharType type, ValueModeType mode, TimeStep *tStep, bool global)
111{
112 answer.clear();
113 if ( type != ExternalForcesVector ) {
114 return;
115 }
116
117 FEInterpolation *fei = this->giveInterpolation();
118 if ( !fei ) {
119 OOFEM_ERROR("No interpolator available");
120 }
121
122 FloatArray n_vec;
123 FloatMatrix n, T;
124 FloatArray force, globalIPcoords;
125 //int nsd = fei->giveNsd();
126
127 std :: unique_ptr< IntegrationRule >iRule(this->giveBoundarySurfaceIntegrationRule(load->giveApproxOrder(), boundary) );
128
129 for ( GaussPoint *gp : * iRule ) {
130 const FloatArray &lcoords = gp->giveNaturalCoordinates();
131
132 if ( load->giveFormulationType() == Load :: FT_Entity ) {
133 load->computeValueAt(force, tStep, lcoords, mode);
134 } else {
135 fei->boundaryLocal2Global(globalIPcoords, boundary, lcoords, FEIElementGeometryWrapper(this) );
136 load->computeValueAt(force, tStep, globalIPcoords, mode);
137 }
138
140 // We always want the global values in the end, so we might as well compute them here directly:
141 // transform force
142 if ( load->giveCoordSystMode() == Load :: CST_Global ) {
143 // then just keep it in global c.s
144 } else {
146 // transform from local boundary to element local c.s
147 // uncommented since the other (now commented) approach did not work correctly
148 if ( this->computeLoadLSToLRotationMatrix(T, boundary, gp) ) {
149 force.rotatedWith(T, 'n');
150 }
151 // then to global c.s
152 //if ( this->computeLoadGToLRotationMtrx(T) ) {
153 // force.rotatedWith(T, 't');
154 //}
155 }
156
157 // Construct n-matrix
158 this->computeSurfaceNMatrix(n, boundary, lcoords); // to allow adaptation on element level
159
161 double thickness = 1.0; // Should be the circumference for axisymm-elements.
162 double dV = thickness * this->computeSurfaceVolumeAround(gp, boundary);
163 answer.plusProduct(n, force, dV);
164 }
165}
166
167double
168StructuralElement :: computeSurfaceVolumeAround(GaussPoint *gp, int iSurf)
169{
170 FEInterpolation *fei = this->giveInterpolation();
171 const FloatArray &lcoords = gp->giveNaturalCoordinates();
172 double J = fei->boundarySurfaceGiveTransformationJacobian(iSurf, lcoords, FEIElementGeometryWrapper(this) );
173
174 return ( gp->giveWeight() * J );
175}
176
177void
178StructuralElement :: computeSurfaceNMatrix(FloatMatrix &answer, int boundaryID, const FloatArray &lcoords)
179{
180 FloatArray n_vec;
181 this->giveInterpolation()->boundarySurfaceEvalN(n_vec, boundaryID, lcoords, FEIElementGeometryWrapper(this) );
182 answer.beNMatrixOf(n_vec, this->giveInterpolation()->giveNsd(this->giveGeometryType()) );
183}
184
185
186void StructuralElement :: computeBoundaryEdgeLoadVector(FloatArray &answer, BoundaryLoad *load, int boundary, CharType type, ValueModeType mode, TimeStep *tStep, bool global)
187{
188 answer.clear();
189 if ( type != ExternalForcesVector ) {
190 return;
191 }
192
193 FEInterpolation *fei = this->giveInterpolation();
194 if ( !fei ) {
195 OOFEM_ERROR("No interpolator available");
196 }
197
198 FloatMatrix n, T;
199 FloatArray force, globalIPcoords;
200
201 std :: unique_ptr< IntegrationRule >iRule(this->giveBoundaryEdgeIntegrationRule(load->giveApproxOrder(), boundary) );
202
203 for ( GaussPoint *gp : * iRule ) {
204 const FloatArray &lcoords = gp->giveNaturalCoordinates();
205
206 if ( load->giveFormulationType() == Load :: FT_Entity ) {
207 load->computeValueAt(force, tStep, lcoords, mode);
208 } else {
209 fei->boundaryEdgeLocal2Global(globalIPcoords, boundary, lcoords, FEIElementGeometryWrapper(this) );
210 load->computeValueAt(force, tStep, globalIPcoords, mode);
211 }
212
214 // We always want the global values in the end, so we might as well compute them here directly:
215 // transform force
216 if ( load->giveCoordSystMode() == Load :: CST_Global ) {
217 // then just keep it in global c.s
218 } else {
220 // transform from local boundary to element local c.s
221 if ( this->computeLoadLEToLRotationMatrix(T, boundary, gp) ) {
222 force.rotatedWith(T, 'n');
223 }
224 // then to global c.s
225 if ( this->computeLoadGToLRotationMtrx(T) ) {
226 force.rotatedWith(T, 't');
227 }
228 }
229
230 // Construct n-matrix
231 //fei->boundaryEdgeEvalN( n_vec, boundary, lcoords, FEIElementGeometryWrapper(this) );
232 //n.beNMatrixOf(n_vec, nsd);
233 this->computeEdgeNMatrix(n, boundary, lcoords); // to allow adapttation on element level
234
235 double dV = this->computeEdgeVolumeAround(gp, boundary);
236 answer.plusProduct(n, force, dV);
237 }
238}
239
240double
241StructuralElement :: computeEdgeVolumeAround(GaussPoint *gp, int iEdge)
242{
243 FEInterpolation *fei = this->giveInterpolation();
244 const FloatArray &lcoords = gp->giveNaturalCoordinates();
245 double J = fei->boundaryEdgeGiveTransformationJacobian(iEdge, lcoords, FEIElementGeometryWrapper(this) );
246
247 return ( gp->giveWeight() * J );
248}
249
250void
251StructuralElement :: computeEdgeNMatrix(FloatMatrix &answer, int boundaryID, const FloatArray &lcoords)
252{
253 FloatArray n_vec;
254 this->giveInterpolation()->boundaryEdgeEvalN(n_vec, boundaryID, lcoords, FEIElementGeometryWrapper(this) );
255 answer.beNMatrixOf(n_vec, this->giveInterpolation()->giveNsd(this->giveGeometryType()) );
256}
257
258
259
260
261void
262StructuralElement :: computeBodyLoadVectorAt(FloatArray &answer, Load *forLoad, TimeStep *tStep, ValueModeType mode)
263// Computes numerically the load vector of the receiver due to the body
264// loads, at tStep.
265// load is assumed to be in global cs.
266// load vector is then transformed to coordinate system in each node.
267// (should be global coordinate system, but there may be defined
268// different coordinate system in each node)
269{
270 double dens, dV;
271 FloatArray force, ntf;
272 FloatMatrix n, T;
273
274 if ( ( forLoad->giveBCGeoType() != BodyLoadBGT ) || ( forLoad->giveBCValType() != ForceLoadBVT ) ) {
275 OOFEM_ERROR("unknown load type");
276 }
277
278 // note: force is assumed to be in global coordinate system.
279 forLoad->computeComponentArrayAt(force, tStep, mode);
280 // transform from global to element local c.s
281 if ( this->computeLoadGToLRotationMtrx(T) ) {
282 force.rotatedWith(T, 'n');
283 }
284
285 answer.clear();
286
287 if ( force.giveSize() ) {
288 for ( GaussPoint *gp : * this->giveDefaultIntegrationRulePtr() ) {
290 dV = this->computeVolumeAround(gp);
291 dens = this->giveCrossSection()->give('d', gp);
292 ntf.beTProductOf(n, force);
293 answer.add(dV * dens, ntf);
294 }
295 } else {
296 return;
297 }
298}
299
300void
301StructuralElement :: computePointLoadVectorAt(FloatArray &answer, Load *load, TimeStep *tStep, ValueModeType mode, bool global)
302{
303 FloatArray force, lcoords;
304 FloatMatrix T, n;
305
306 PointLoad *pointLoad = dynamic_cast< PointLoad * >( load );
307 FloatArray coords = pointLoad->giveCoordinates();
308 pointLoad->computeValueAt(force, tStep, coords, mode);
309 if ( this->computeLocalCoordinates(lcoords, pointLoad->giveCoordinates() ) ) {
310 this->computeNmatrixAt(lcoords, n);
311 answer.beTProductOf(n, force);
312 } else {
313 OOFEM_WARNING("point load outside element");
314 }
315
316 // transform force
317 if ( pointLoad->giveCoordSystMode() == Load :: CST_Global ) {
318 // transform from global to element local c.s
319 if ( this->computeLoadGToLRotationMtrx(T) ) {
320 answer.rotatedWith(T, 'n');
321 }
322 }
323}
324
325int
326StructuralElement :: giveNumberOfIPForMassMtrxIntegration()
327{
328 IntegrationRule *iRule = this->giveIntegrationRule(0);
329 // returns necessary number of ip to integrate the mass matrix
330 // \int_V N^T*N dV => (order of the approximation)*2 (constant density assumed)
332 int order = this->giveInterpolation()->giveInterpolationOrder();
334}
335
336void
337StructuralElement :: computeConsistentMassMatrix(FloatMatrix &answer, TimeStep *tStep, double &mass, const double *ipDensity)
338// Computes numerically the consistent (full) mass matrix of the receiver.
339{
340 int ndofs = computeNumberOfDofs();
341 double density, dV;
342 FloatMatrix n;
343 GaussIntegrationRule iRule(1, this, 1, 1);
344 IntArray mask;
345
346 answer.resize(ndofs, ndofs);
347 answer.zero();
348 if ( !this->isActivated(tStep) ) {
349 return;
350 }
351
353
355
356 mass = 0.;
357
358 for ( GaussPoint *gp : iRule ) {
360 density = this->giveCrossSection()->give('d', gp);
361
362 if ( ipDensity != NULL ) {
363 // Override density if desired
364 density = * ipDensity;
365 }
366
367 dV = this->computeVolumeAround(gp);
368 mass += density * dV;
369
370 if ( mask.isEmpty() ) {
371 answer.plusProductSymmUpper(n, n, density * dV);
372 } else {
373 for ( int i = 1; i <= ndofs; i++ ) {
374 for ( int j = i; j <= ndofs; j++ ) {
375 double summ = 0.;
376 for ( int k = 1; k <= n.giveNumberOfRows(); k++ ) {
377 if ( mask.at(k) == 0 ) {
378 continue;
379 }
380
381 summ += n.at(k, i) * n.at(k, j);
382 }
383
384 answer.at(i, j) += summ * density * dV;
385 }
386 }
387 }
388 }
389
390 answer.symmetrized();
391}
392
393void
394StructuralElement :: setupIRForMassMtrxIntegration(IntegrationRule &iRule)
395{
396 int nip;
397 if ( ( nip = this->giveNumberOfIPForMassMtrxIntegration() ) == 0 ) {
398 OOFEM_ERROR("no integration points available");
399 }
401}
402
403void
404StructuralElement :: computeMassMatrix(FloatMatrix &answer, TimeStep *tStep)
405// Returns the lumped mass matrix of the receiver.
406{
407 double mass;
408
409 this->computeConsistentMassMatrix(answer, tStep, mass);
410}
411
412void
413StructuralElement :: computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep)
414// Returns the lumped mass matrix of the receiver.
415{
416 double mass = 0.;
417
418 IntArray nodeDofIDMask, dimFlag(3);
419 int indx = 0, ldofs, dim;
420 double summ;
421
422 if ( !this->isActivated(tStep) ) {
423 int ndofs = computeNumberOfDofs();
424 answer.resize(ndofs, ndofs);
425 answer.zero();
426 return;
427 }
428
429 this->computeConsistentMassMatrix(answer, tStep, mass);
430 ldofs = answer.giveNumberOfRows();
431
432 for ( int i = 1; i <= this->giveNumberOfDofManagers(); i++ ) {
433 this->giveDofManDofIDMask(i, nodeDofIDMask);
434 //this->giveDofManager(i)->giveLocationArray(nodeDofIDMask, nodalArray);
435 for ( int j = 1; j <= nodeDofIDMask.giveSize(); j++ ) {
436 indx++;
437 // zero all off-diagonal terms
438 for ( int k = 1; k <= ldofs; k++ ) {
439 if ( k != indx ) {
440 answer.at(indx, k) = 0.;
441 answer.at(k, indx) = 0.;
442 }
443 }
444
445 if ( ( nodeDofIDMask.at(j) != D_u ) && ( nodeDofIDMask.at(j) != D_v ) && ( nodeDofIDMask.at(j) != D_w ) ) {
446 // zero corresponding diagonal member too <= no displacement dof
447 answer.at(indx, indx) = 0.;
448 } else if ( nodeDofIDMask.at(j) == D_u ) {
449 dimFlag.at(1) = 1;
450 } else if ( nodeDofIDMask.at(j) == D_v ) {
451 dimFlag.at(2) = 1;
452 } else if ( nodeDofIDMask.at(j) == D_w ) {
453 dimFlag.at(3) = 1;
454 }
455 }
456 }
457
458 if ( indx != ldofs ) {
459 OOFEM_ERROR("internal consistency check failed");
460 }
461
462 dim = dimFlag.at(1) + dimFlag.at(2) + dimFlag.at(3);
463 summ = 0.;
464 for ( int k = 1; k <= ldofs; k++ ) {
465 summ += answer.at(k, k);
466 }
467
468 answer.times(dim * mass / summ);
469}
470
471
472void
473StructuralElement :: computeResultingIPTemperatureAt(FloatArray &answer, TimeStep *tStep, GaussPoint *gp, ValueModeType mode)
474// Computes at tStep the resulting force due to all temperature loads that act
475// on the receiver. This force is used by the element for computing its
476// body load vector.
477{
478 int n, nLoads;
479 Load *load;
480 FloatArray gCoords, temperature;
481 //int nbc = domain->giveNumberOfBoundaryConditions();
482
483 if ( this->computeGlobalCoordinates(gCoords, gp->giveNaturalCoordinates() ) == 0 ) {
484 OOFEM_ERROR("computeGlobalCoordinates failed");
485 }
486
487 answer.clear();
488 nLoads = this->giveBodyLoadArray()->giveSize();
489 for ( int i = 1; i <= nLoads; i++ ) {
490 n = bodyLoadArray.at(i);
491 load = domain->giveLoad(n);
492 if ( load->giveBCValType() == TemperatureBVT ) {
493 load->computeValueAt(temperature, tStep, gCoords, mode);
494 answer.add(temperature);
495 }
496 }
497
498
499 // new approach using sets
500 // add exact end forces due to nonnodal loading applied indirectly (via sets)
501 BCTracker *bct = this->domain->giveBCTracker();
502 const BCTracker::entryListType& bcList = bct->getElementRecords(this->number);
503
504 for (BCTracker::entryListType::const_iterator it = bcList.begin(); it != bcList.end(); ++it) {
505 GeneralBoundaryCondition *bc = this->domain->giveBc((*it).bcNumber);
506 if ( ( load = dynamic_cast< StructuralTemperatureLoad * >( bc ) ) ) {
507 if ( bc->giveSetNumber() && bc->isImposed(tStep) ) {
508 if ( load->giveBCValType() == TemperatureBVT ) {
509 load->computeValueAt(temperature, tStep, gCoords, mode);
510 answer.add(temperature);
511 }
512 }
513 }
514 }
515}
516
517void
518StructuralElement :: computeResultingIPEigenstrainAt(FloatArray &answer, TimeStep *tStep, GaussPoint *gp, ValueModeType mode)
519// Computes at tStep all eigenstrains that act on the receiver. Eigenstrains are used by the element for computing its body load vector.
520{
521 int n, nLoads;
522 Load *load;
523 FloatArray gCoords, eigenstrain;
524 //int nbc = domain->giveNumberOfBoundaryConditions();
525
526 if ( this->computeGlobalCoordinates(gCoords, gp->giveNaturalCoordinates() ) == 0 ) {
527 OOFEM_ERROR("computeGlobalCoordinates failed");
528 }
529
530 answer.clear();
531
532 // old approach preserved for compatibility
533
534 nLoads = this->giveBodyLoadArray()->giveSize();
535 for ( int i = 1; i <= nLoads; i++ ) {
536 n = bodyLoadArray.at(i);
537 load = domain->giveLoad(n);
538 if ( load->giveBCValType() == EigenstrainBVT ) {
539 load->computeValueAt(eigenstrain, tStep, gCoords, mode);
540 answer.add(eigenstrain);
541 }
542 }
543
544 // new approach using sets
545 // add exact end forces due to nonnodal loading applied indirectly (via sets)
546 BCTracker *bct = this->domain->giveBCTracker();
548
549 for (BCTracker::entryListType::iterator it = bcList.begin(); it != bcList.end(); ++it) {
550 GeneralBoundaryCondition *bc = this->domain->giveBc((*it).bcNumber);
551
552 if ( ( load = dynamic_cast< StructuralEigenstrainLoad * >( bc ) ) ) {
553 if ( bc->giveSetNumber() && bc->isImposed(tStep) ) {
554 if ( load->giveBCValType() == EigenstrainBVT ) {
555 load->computeValueAt(eigenstrain, tStep, gCoords, mode);
556 answer.add(eigenstrain);
557 }
558 }
559 }
560 }
561}
562
563
564void
565StructuralElement :: computeField(ValueModeType mode, TimeStep *tStep, const FloatArray &lcoords, FloatArray &answer)
566{
567 FloatArray u;
568 FloatMatrix n;
569 this->computeNmatrixAt(lcoords, n);
570 this->computeVectorOf(mode, tStep, u);
571 answer.beProductOf(n, u);
572}
573
574void
575StructuralElement :: computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode,
576 TimeStep *tStep)
577// Computes numerically the stiffness matrix of the receiver.
578{
579 int iStartIndx, iEndIndx, jStartIndx, jEndIndx;
580 double dV;
581 FloatMatrix d, bi, bj, dbj, dij;
582 bool matStiffSymmFlag = this->giveCrossSection()->isCharacteristicMtrxSymmetric(rMode);
583
584 answer.clear();
585
586 if ( !this->isActivated(tStep) ) {
587 return;
588 }
589
590 if ( integrationRulesArray.size() > 1 ) {
591 for ( int i = 0; i < ( int ) integrationRulesArray.size(); i++ ) {
592 iStartIndx = integrationRulesArray [ i ]->getStartIndexOfLocalStrainWhereApply();
593 iEndIndx = integrationRulesArray [ i ]->getEndIndexOfLocalStrainWhereApply();
594 for ( int j = 0; j < ( int ) integrationRulesArray.size(); j++ ) {
595 IntegrationRule *iRule;
596 jStartIndx = integrationRulesArray [ j ]->getStartIndexOfLocalStrainWhereApply();
597 jEndIndx = integrationRulesArray [ j ]->getEndIndexOfLocalStrainWhereApply();
598 if ( i == j ) {
599 iRule = integrationRulesArray [ i ].get();
600 } else if ( integrationRulesArray [ i ]->giveNumberOfIntegrationPoints() < integrationRulesArray [ j ]->giveNumberOfIntegrationPoints() ) {
601 iRule = integrationRulesArray [ i ].get();
602 } else {
603 iRule = integrationRulesArray [ j ].get();
604 }
605
606 for ( GaussPoint *gp : * iRule ) {
607 this->computeBmatrixAt(gp, bi, iStartIndx, iEndIndx);
608 this->computeConstitutiveMatrixAt(d, rMode, gp, tStep);
609 dij.beSubMatrixOf(d, iStartIndx, iEndIndx, jStartIndx, jEndIndx);
610 if ( i != j ) {
611 this->computeBmatrixAt(gp, bj, jStartIndx, jEndIndx);
612 } else {
613 bj = bi;
614 }
615
616 dV = this->computeVolumeAround(gp);
617 dbj.beProductOf(dij, bj);
618 if ( matStiffSymmFlag ) {
619 answer.plusProductSymmUpper(bi, dbj, dV);
620 } else {
621 answer.plusProductUnsym(bi, dbj, dV);
622 }
623 }
624 }
625 }
626 } else {
627 for ( GaussPoint *gp : * this->giveDefaultIntegrationRulePtr() ) {
628 this->computeBmatrixAt(gp, bj);
629 this->computeConstitutiveMatrixAt(d, rMode, gp, tStep);
630 dV = this->computeVolumeAround(gp);
631 dbj.beProductOf(d, bj);
632 if ( matStiffSymmFlag ) {
633 answer.plusProductSymmUpper(bj, dbj, dV);
634 } else {
635 answer.plusProductUnsym(bj, dbj, dV);
636 }
637 }
638 }
639
640 if ( matStiffSymmFlag ) {
641 answer.symmetrized();
642 }
643}
644
645void StructuralElement :: computeStiffnessMatrix_withIRulesAsSubcells(FloatMatrix &answer,
646 MatResponseMode rMode, TimeStep *tStep)
647{
648 FloatMatrix temp, bj, d, dbj;
649 int ndofs = this->computeNumberOfDofs();
650 bool matStiffSymmFlag = this->giveCrossSection()->isCharacteristicMtrxSymmetric(rMode);
651 IntArray irlocnum;
652
653 answer.resize(ndofs, ndofs);
654 answer.zero();
655
656 FloatMatrix *m = & answer;
657 if ( this->giveInterpolation() && this->giveInterpolation()->hasSubPatchFormulation() ) {
658 m = & temp;
659 }
660
661 // loop over individual integration rules
662 for ( auto &iRule : integrationRulesArray ) {
663 // loop over individual integration points
664 for ( GaussPoint *gp : * iRule ) {
665 double dV = this->computeVolumeAround(gp);
666 this->computeBmatrixAt(gp, bj);
667 this->computeConstitutiveMatrixAt(d, rMode, gp, tStep);
668
669 dbj.beProductOf(d, bj);
670 if ( matStiffSymmFlag ) {
671 m->plusProductSymmUpper(bj, dbj, dV);
672 } else {
673 m->plusProductUnsym(bj, dbj, dV);
674 }
675 }
676
677 // localize irule contribution into element matrix
678 if ( this->giveIntegrationRuleLocalCodeNumbers(irlocnum, * iRule) ) {
679 answer.assemble(* m, irlocnum);
680 m->clear();
681 }
682 } // end loop over irules
683
684 if ( matStiffSymmFlag ) {
685 answer.symmetrized();
686 }
687}
688
689
690void
691StructuralElement :: computeStrainVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep)
692// Computes the vector containing the strains at the Gauss point gp of
693// the receiver, at time step tStep. The nature of these strains depends
694// on the element's type.
695{
696 FloatMatrix b;
697 FloatArray u;
698
699 if ( !this->isActivated(tStep) ) {
700 answer.resize(StructuralMaterial :: giveSizeOfVoigtSymVector(gp->giveMaterialMode() ) );
701 answer.zero();
702 return;
703 }
704
705 this->computeBmatrixAt(gp, b);
706 this->computeVectorOf(VM_Total, tStep, u);
707
708 // subtract initial displacements, if defined
709 if ( initialDisplacements ) {
711 }
712 answer.beProductOf(b, u);
713}
714
715#if 0
716void
717StructuralElement :: computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep)
718{
719 this->giveStructuralCrossSection()->giveRealStresses(answer, gp, strain, tStep);
720}
721#endif
722
723void
724StructuralElement :: giveInternalForcesVector(FloatArray &answer,
725 TimeStep *tStep, int useUpdatedGpRecord)
726//
727// returns nodal representation of real internal forces - necessary only for
728// non-linear analysis.
729// if useGpRecord == 1 then data stored in gp->giveStressVector() are used
730// instead computing stressVector through this->ComputeStressVector();
731// this must be done after you want internal forces after element->updateYourself()
732// has been called for the same time step.
733//
734{
735 FloatMatrix b;
736 FloatArray u, stress, strain;
737
738 // This function can be quite costly to do inside the loops when one has many slave dofs.
739 this->computeVectorOf(VM_Total, tStep, u);
740 // subtract initial displacements, if defined
741 if ( initialDisplacements ) {
743 }
744
745 // zero answer will resize accordingly when adding first contribution
746 answer.clear();
747
748 for ( GaussPoint *gp : * this->giveDefaultIntegrationRulePtr() ) {
749 this->computeBmatrixAt(gp, b);
750
751 if ( useUpdatedGpRecord == 1 ) {
752 auto status = gp->giveMaterialStatus();
753 StructuralMaterialStatus *matStat = dynamic_cast< StructuralMaterialStatus * >( status );
754 if ( matStat ) {
755 stress = matStat->giveStressVector();
756 } else {
757 LatticeMaterialStatus *lmatStat = dynamic_cast< LatticeMaterialStatus * >( status );
758 if ( lmatStat ) {
759 stress = lmatStat->giveLatticeStress();
760 } else {
762 stress = ms->giveTraction();
763 }
764 }
765 } else {
766 if ( !this->isActivated(tStep) ) {
767 strain.resize(StructuralMaterial :: giveSizeOfVoigtSymVector(gp->giveMaterialMode() ) );
768 strain.zero();
769 }
770 strain.beProductOf(b, u);
771 this->computeStressVector(stress, strain, gp, tStep);
772 }
773
774 // updates gp stress and strain record acording to current
775 // increment of displacement
776 if ( stress.giveSize() == 0 ) {
777 break;
778 }
779
780 // now every gauss point has real stress vector
781 // compute nodal representation of internal forces using f = B^T*Sigma dV
782 double dV = this->computeVolumeAround(gp);
783 if ( stress.giveSize() == 6 ) {
784 // It may happen that e.g. plane strain is computed
785 // using the default 3D implementation. If so,
786 // the stress needs to be reduced.
787 // (Note that no reduction will take place if
788 // the simulation is actually 3D.)
789 FloatArray stressTemp;
790 StructuralMaterial :: giveReducedSymVectorForm(stressTemp, stress, gp->giveMaterialMode() );
791 answer.plusProduct(b, stressTemp, dV);
792 } else {
793 answer.plusProduct(b, stress, dV);
794 }
795 }
796
797 // if inactive update state, but no contribution to global system
798 if ( !this->isActivated(tStep) ) {
799 answer.zero();
800 return;
801 }
802}
803
804
805
806void
807StructuralElement :: giveInternalForcesVector_withIRulesAsSubcells(FloatArray &answer,
808 TimeStep *tStep, int useUpdatedGpRecord)
809//
810// returns nodal representation of real internal forces - necessary only for
811// non-linear analysis.
812// if useGpRecord == 1 then data stored in gp->giveStressVector() are used
813// instead computing stressVector through this->ComputeStressVector();
814// this must be done after you want internal forces after element->updateYourself()
815// has been called for the same time step.
816//
817{
818 FloatMatrix b;
819 FloatArray temp, u, stress, strain;
820 IntArray irlocnum;
821
822 // This function can be quite costly to do inside the loops when one has many slave dofs.
823 this->computeVectorOf(VM_Total, tStep, u);
824 // subtract initial displacements, if defined
825 if ( initialDisplacements ) {
827 }
828
829 // zero answer will resize accordingly when adding first contribution
830 answer.clear();
831
832 FloatArray *m = & answer;
833 if ( this->giveInterpolation() && this->giveInterpolation()->hasSubPatchFormulation() ) {
834 m = & temp;
835 }
836
837 // loop over individual integration rules
838 for ( auto &iRule : integrationRulesArray ) {
839 for ( GaussPoint *gp : * iRule ) {
840 StructuralMaterialStatus *matStat = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() );
841 this->computeBmatrixAt(gp, b);
842
843 if ( useUpdatedGpRecord == 1 ) {
844 stress = matStat->giveStressVector();
845 } else {
846 if ( !this->isActivated(tStep) ) {
847 strain.resize(StructuralMaterial :: giveSizeOfVoigtSymVector(gp->giveMaterialMode() ) );
848 strain.zero();
849 }
850 strain.beProductOf(b, u);
851 this->computeStressVector(stress, strain, gp, tStep);
852 }
853
854 //
855 // updates gp stress and strain record acording to current
856 // increment of displacement
857 //
858 if ( stress.giveSize() == 0 ) {
859 break;
860 }
861
862 //
863 // now every gauss point has real stress vector
864 //
865 // compute nodal representation of internal forces using f = B^T*Sigma dV
866 //
867 double dV = this->computeVolumeAround(gp);
868 m->plusProduct(b, stress, dV);
869
870 // localize irule contribution into element matrix
871 if ( this->giveIntegrationRuleLocalCodeNumbers(irlocnum, * iRule) ) {
872 answer.assemble(* m, irlocnum);
873 m->clear();
874 }
875 }
876 }
877
878 // if inactive update state, but no contribution to global system
879 if ( !this->isActivated(tStep) ) {
880 answer.zero();
881 return;
882 }
883}
884
885
886void
887StructuralElement :: giveCharacteristicMatrix(FloatMatrix &answer,
888 CharType mtrx, TimeStep *tStep)
889//
890// returns characteristics matrix of receiver according to mtrx
891//
892{
893 if ( mtrx == TangentStiffnessMatrix ) {
894 this->computeStiffnessMatrix(answer, TangentStiffness, tStep);
895 } else if ( mtrx == SecantStiffnessMatrix ) {
896 this->computeStiffnessMatrix(answer, SecantStiffness, tStep);
897 } else if ( mtrx == ElasticStiffnessMatrix ) {
898 this->computeStiffnessMatrix(answer, ElasticStiffness, tStep);
899 } else if ( mtrx == MassMatrix ) {
900 this->computeMassMatrix(answer, tStep);
901 } else if ( mtrx == LumpedMassMatrix ) {
902 this->computeLumpedMassMatrix(answer, tStep);
903 } else if ( mtrx == InitialStressMatrix ) {
904 this->computeInitialStressMatrix(answer, tStep);
905 } else if ( mtrx == LumpedInitialStressMatrix ) {
906 this->computeLumpedInitialStressMatrix(answer, tStep);
907 } else {
908 OOFEM_ERROR("Unknown Type of characteristic mtrx (%s)", __CharTypeToString(mtrx) );
909 }
910}
911
912
913
914void
915StructuralElement :: giveCharacteristicVector(FloatArray &answer, CharType mtrx, ValueModeType mode,
916 TimeStep *tStep)
917//
918// returns characteristics vector of receiver according to mtrx
919//
920{
921 if ( mtrx == ExternalForcesVector ) {
922 //this->computeForceLoadVector(answer, tStep, mode); // bp: assembled by emodel
923 answer.resize(0);
924 } else if ( ( mtrx == InternalForcesVector ) && ( mode == VM_Total ) ) {
925 this->giveInternalForcesVector(answer, tStep);
926 } else if ( ( mtrx == LastEquilibratedInternalForcesVector ) && ( mode == VM_Total ) ) {
927 /* here tstep is not relevant, we set useUpdatedGpRecord = 1
928 * and this will cause to integrate internal forces using existing (nontemp, equlibrated) stresses in
929 * statuses. Mainly used to compute reaction forces */
930 this->giveInternalForcesVector(answer, tStep, 1);
931 } else if ( mtrx == LumpedMassMatrix ) {
932 FloatMatrix M;
933 this->computeLumpedMassMatrix(M, tStep);
934 answer.resize(M.giveNumberOfColumns() );
935 for ( int i = 0; i < M.giveNumberOfColumns(); ++i ) {
936 answer [ i ] = M(i, i);
937 }
938 } else {
939 OOFEM_ERROR("Unknown Type of characteristic mtrx (%s)", __CharTypeToString(mtrx) );
940 }
941}
942
943void
944StructuralElement :: updateYourself(TimeStep *tStep)
945{
946 Element :: updateYourself(tStep);
947
948 // record initial displacement if element not active
949 if ( activityTimeFunction && !isActivated(tStep) ) {
950 if ( !initialDisplacements ) {
951 initialDisplacements = std :: make_unique< FloatArray >();
952 }
953
954 this->computeVectorOf(VM_Total, tStep, * initialDisplacements);
955 }
956}
957
958
959void
960StructuralElement :: updateInternalState(TimeStep *tStep)
961// Updates the receiver at end of step.
962{
963 FloatArray stress, strain;
964
965 // force updating strains & stresses
966 for ( auto &iRule : integrationRulesArray ) {
967 for ( GaussPoint *gp : * iRule ) {
968 this->computeStrainVector(strain, gp, tStep);
969 this->computeStressVector(stress, strain, gp, tStep);
970 }
971 }
972}
973
974void
975StructuralElement :: updateBeforeNonlocalAverage(TimeStep *tStep)
976// Updates the local material quantities before nonlocal averaging
977{
978 /*
979 * Nonlocal material, related to receiver is also passed, because it is not possible to
980 * ask receiver for its material model pointer using giveMaterial() service (returns Material type)
981 * and then to cast this Material pointer into NonlocalMaterial pointer type,
982 * because it is possible to cast only the pointer of derived class to pointer to base class.
983 */
984 FloatArray epsilon;
985
986 if ( this->giveParallelMode() == Element_remote ) {
987 return;
988 }
989
990 // force updating local quantities
991 for ( auto &iRule : integrationRulesArray ) {
992 for ( GaussPoint *gp : * iRule ) {
993 this->computeStrainVector(epsilon, gp, tStep);
994 // provide material local strain increment - as is provided to computeRealStresVector
995 // allows to update internal vars to be averaged to new state
996
997 // not possible - produces wrong result
999 materialExt = static_cast< StructuralNonlocalMaterialExtensionInterface * >( this->giveStructuralCrossSection()->
1000 giveMaterialInterface(NonlocalMaterialExtensionInterfaceType, gp) );
1001
1002 if ( !materialExt ) {
1003 return; //_error("updateBeforeNonlocalAverage: material with no StructuralNonlocalMaterial support");
1004 }
1005 materialExt->updateBeforeNonlocAverage(epsilon, gp, tStep);
1006 }
1007 }
1008}
1009
1010
1011int
1012StructuralElement :: checkConsistency()
1013//
1014// check internal consistency
1015// mainly tests, whether crossSection data
1016// are safe for conversion to "Structural" version
1017//
1018{
1019 int result = 1;
1020 if ( !this->giveCrossSection()->testCrossSectionExtension(CS_StructuralCapability) ) {
1021 OOFEM_WARNING("cross-section %s without structural support", this->giveCrossSection()->giveClassName() );
1022 result = 0;
1023 }
1024
1025 return result;
1026}
1027
1028void
1029StructuralElement :: computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer)
1030{
1031 int numNodes = this->giveNumberOfDofManagers();
1032 FloatArray N(numNodes);
1033
1034 int dim = this->giveSpatialDimension();
1035
1036 answer.resize(dim, dim * numNodes);
1037 answer.zero();
1038 giveInterpolation()->evalN(N, iLocCoord, FEIElementGeometryWrapper(this) );
1039
1040 answer.beNMatrixOf(N, dim);
1041}
1042
1043void
1044StructuralElement :: condense(FloatMatrix *stiff, FloatMatrix *mass, FloatArray *load, IntArray *what)
1045{
1046 /*
1047 * function for condensation of stiffness matrix and if requested of load vector and
1048 * initial stress or mass matrices (if mass and load arguments are nonzero (not NULL) then
1049 * their condensation is done).
1050 * Based on Rayleigh-Ritz method
1051 */
1052 int i, ii, j, k;
1053 int nkon = what->giveSize();
1054 int size = stiff->giveNumberOfRows();
1055 double coeff, dii, lii = 0;
1056 FloatArray gaussCoeff;
1057 // stored gauss coefficient
1058 // for mass condensation
1059
1060 // check
1061 if ( !stiff->isSquare() ) {
1062 OOFEM_ERROR("stiffness size mismatch");
1063 }
1064
1065 if ( mass ) {
1066 if ( !( mass->isSquare() && mass->giveNumberOfRows() == size ) ) {
1067 OOFEM_ERROR("mass size mismatch");
1068 }
1069 }
1070
1071 if ( load ) {
1072 if ( !( load->giveSize() == size ) ) {
1073 OOFEM_ERROR("load size mismatch");
1074 }
1075 }
1076
1077 // create gauss coeff array if mass condensation requested
1078 if ( mass ) {
1079 gaussCoeff.resize(size);
1080 }
1081
1082 for ( i = 1; i <= nkon; i++ ) {
1083 ii = what->at(i);
1084 if ( ( ii > size ) || ( ii <= 0 ) ) {
1085 OOFEM_ERROR("wrong dof number");
1086 }
1087
1088 dii = stiff->at(ii, ii);
1089 if ( load ) {
1090 lii = load->at(ii);
1091 }
1092
1093 // stiffness matrix condensation
1094 for ( j = 1; j <= size; j++ ) {
1095 coeff = -stiff->at(j, ii) / dii;
1096 if ( ii != j ) {
1097 for ( k = 1; k <= size; k++ ) {
1098 stiff->at(j, k) += stiff->at(ii, k) * coeff;
1099 }
1100 }
1101
1102 if ( load ) {
1103 load->at(j) += coeff * lii;
1104 }
1105
1106 if ( mass ) {
1107 gaussCoeff.at(j) = coeff;
1108 }
1109 }
1110
1111 for ( k = 1; k <= size; k++ ) {
1112 stiff->at(ii, k) = 0.;
1113 stiff->at(k, ii) = 0.;
1114 }
1115
1116 // mass or initial stress matrix condensation
1117 // uses gauss coefficients stored in gaussCoeff.
1118 if ( mass ) {
1119 for ( j = 1; j <= size; j++ ) {
1120 for ( k = 1; k <= size; k++ ) {
1121 if ( ( ii != j ) && ( ii != k ) ) {
1122 mass->at(j, k) += mass->at(j, ii) * gaussCoeff.at(k) +
1123 mass->at(ii, k) * gaussCoeff.at(j) +
1124 mass->at(ii, ii) * gaussCoeff.at(j) * gaussCoeff.at(k);
1125 }
1126 }
1127 }
1128
1129 for ( k = 1; k <= size; k++ ) {
1130 mass->at(ii, k) = 0.;
1131 mass->at(k, ii) = 0.;
1132 }
1133 }
1134 }
1135}
1136
1137
1138int
1139StructuralElement :: giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
1140{
1141 if ( type == IST_DisplacementVector ) {
1142 FloatArray u;
1143 FloatMatrix N;
1144 this->computeVectorOf(VM_Total, tStep, u);
1146 answer.beProductOf(N, u);
1147 return 1;
1148 }
1149 return Element :: giveIPValue(answer, gp, type, tStep);
1150}
1151
1152
1153void
1154StructuralElement :: giveNonlocalLocationArray(IntArray &locationArray, const UnknownNumberingScheme &s)
1155{
1156 IntArray elemLocArry;
1157
1158 locationArray.clear();
1159 // loop over element IP
1160 for ( IntegrationPoint *ip : * this->giveDefaultIntegrationRulePtr() ) {
1161 auto interface = static_cast< NonlocalMaterialStiffnessInterface * >( this->giveStructuralCrossSection()->
1162 giveMaterialInterface(NonlocalMaterialStiffnessInterfaceType, ip) );
1163
1164
1165 if ( interface == NULL ) {
1166 locationArray.clear();
1167 return;
1168 }
1169
1170 auto integrationDomainList = interface->
1171 NonlocalMaterialStiffnessInterface_giveIntegrationDomainList(ip);
1172 // loop over IP influencing IPs, extract corresponding element numbers and their code numbers
1173 for ( auto &lir : * integrationDomainList ) {
1174 lir.nearGp->giveElement()->giveLocationArray(elemLocArry, s);
1175 /*
1176 * Currently no care given to multiple occurences of code number in locationArray.
1177 */
1178 locationArray.followedBy(elemLocArry, 20);
1179 }
1180 } // end loop over IPs
1181}
1182
1183
1184void
1185StructuralElement :: addNonlocalStiffnessContributions(SparseMtrx &dest, const UnknownNumberingScheme &s, TimeStep *tStep)
1186{
1188
1189 if ( !this->isActivated(tStep) ) {
1190 return;
1191 }
1192
1193 // loop over element IP
1194 for ( IntegrationPoint *ip : * this->giveDefaultIntegrationRulePtr() ) {
1195 auto interface = static_cast< NonlocalMaterialStiffnessInterface * >( this->giveStructuralCrossSection()->
1196 giveMaterialInterface(NonlocalMaterialStiffnessInterfaceType, ip) );
1197 if ( interface == NULL ) {
1198 return;
1199 }
1200
1201 interface->NonlocalMaterialStiffnessInterface_addIPContribution(dest, s, ip, tStep);
1202 }
1203}
1204
1205
1206int
1207StructuralElement :: adaptiveUpdate(TimeStep *tStep)
1208{
1209 int result = 1;
1210 FloatArray strain;
1211
1212 for ( auto &iRule : integrationRulesArray ) {
1213 for ( IntegrationPoint *ip : * iRule ) {
1214 auto interface = static_cast< MaterialModelMapperInterface * >( this->giveStructuralCrossSection()->
1215 giveMaterialInterface(MaterialModelMapperInterfaceType, ip) );
1216
1217 if ( interface == NULL ) {
1218 return 0;
1219 }
1220 this->computeStrainVector(strain, ip, tStep);
1221 result &= interface->MMI_update(ip, tStep, & strain);
1222 }
1223 }
1224
1225 return result;
1226}
1227
1228void StructuralElement :: giveInputRecord(DynamicInputRecord &input)
1229{
1230 Element :: giveInputRecord(input);
1231
1233}
1234
1235
1236StructuralCrossSection *StructuralElement :: giveStructuralCrossSection()
1237{
1238 return static_cast< StructuralCrossSection * >( this->giveCrossSection() );
1239}
1240
1241void StructuralElement :: createMaterialStatus()
1242{
1244 for ( auto &iRule : integrationRulesArray ) {
1245 for ( GaussPoint *gp : * iRule ) {
1246 cs->createMaterialStatus(* gp);
1247 }
1248 }
1249}
1250
1251
1252#ifdef __OOFEG
1253
1254//
1255int
1256StructuralElement :: giveInternalStateAtNode(FloatArray &answer, InternalStateType type, InternalStateMode mode, int node, TimeStep *tStep)
1257{
1258 if ( type == IST_DisplacementVector ) {
1259 Node *n = this->giveNode(node);
1260 answer.resize(3);
1261 answer.at(1) = n->giveUpdatedCoordinate(1, tStep) - n->giveCoordinate(1);
1262 answer.at(2) = n->giveUpdatedCoordinate(2, tStep) - n->giveCoordinate(2);
1263 answer.at(3) = n->giveUpdatedCoordinate(3, tStep) - n->giveCoordinate(3);
1264 return 1;
1265 } else {
1266 return Element :: giveInternalStateAtNode(answer, type, mode, node, tStep);
1267 }
1268}
1269
1270
1271
1272void
1273StructuralElement :: showSparseMtrxStructure(CharType mtrx, oofegGraphicContext &gc, TimeStep *tStep)
1274{
1275 if ( mtrx == TangentStiffnessMatrix ||
1276 mtrx == SecantStiffnessMatrix || mtrx == ElasticStiffnessMatrix ) {
1277 int i, j, n;
1278 IntArray loc;
1280
1281 WCRec p [ 4 ];
1282 GraphicObj *go;
1283
1284 EASValsSetLineWidth(OOFEG_SPARSE_PROFILE_WIDTH);
1285 EASValsSetColor(gc.getStandardSparseProfileColor() );
1286 EASValsSetLayer(OOFEG_SPARSE_PROFILE_LAYER);
1287 EASValsSetFillStyle(FILL_SOLID);
1288
1289 n = loc.giveSize();
1290 for ( i = 1; i <= n; i++ ) {
1291 if ( loc.at(i) == 0 ) {
1292 continue;
1293 }
1294
1295 for ( j = i; j <= n; j++ ) {
1296 if ( loc.at(j) == 0 ) {
1297 continue;
1298 }
1299
1300 if ( gc.getSparseProfileMode() == 0 ) {
1301 p [ 0 ].x = ( FPNum ) loc.at(i) - 0.5;
1302 p [ 0 ].y = ( FPNum ) loc.at(j) - 0.5;
1303 p [ 0 ].z = 0.;
1304 p [ 1 ].x = ( FPNum ) loc.at(i) + 0.5;
1305 p [ 1 ].y = ( FPNum ) loc.at(j) - 0.5;
1306 p [ 1 ].z = 0.;
1307 p [ 2 ].x = ( FPNum ) loc.at(i) + 0.5;
1308 p [ 2 ].y = ( FPNum ) loc.at(j) + 0.5;
1309 p [ 2 ].z = 0.;
1310 p [ 3 ].x = ( FPNum ) loc.at(i) - 0.5;
1311 p [ 3 ].y = ( FPNum ) loc.at(j) + 0.5;
1312 p [ 3 ].z = 0.;
1313 go = CreateQuad3D(p);
1314 EGWithMaskChangeAttributes(WIDTH_MASK | FILL_MASK | COLOR_MASK | LAYER_MASK, go);
1315 EMAddGraphicsToModel(ESIModel(), go);
1316
1317 p [ 0 ].x = ( FPNum ) loc.at(j) - 0.5;
1318 p [ 0 ].y = ( FPNum ) loc.at(i) - 0.5;
1319 p [ 0 ].z = 0.;
1320 p [ 1 ].x = ( FPNum ) loc.at(j) + 0.5;
1321 p [ 1 ].y = ( FPNum ) loc.at(i) - 0.5;
1322 p [ 1 ].z = 0.;
1323 p [ 2 ].x = ( FPNum ) loc.at(j) + 0.5;
1324 p [ 2 ].y = ( FPNum ) loc.at(i) + 0.5;
1325 p [ 2 ].z = 0.;
1326 p [ 3 ].x = ( FPNum ) loc.at(j) - 0.5;
1327 p [ 3 ].y = ( FPNum ) loc.at(i) + 0.5;
1328 p [ 3 ].z = 0.;
1329 go = CreateQuad3D(p);
1330 EGWithMaskChangeAttributes(WIDTH_MASK | FILL_MASK | COLOR_MASK | LAYER_MASK, go);
1331 EMAddGraphicsToModel(ESIModel(), go);
1332 } else {
1333 p [ 0 ].x = ( FPNum ) loc.at(i);
1334 p [ 0 ].y = ( FPNum ) loc.at(j);
1335 p [ 0 ].z = 0.;
1336
1337 EASValsSetMType(SQUARE_MARKER);
1338 go = CreateMarker3D(p);
1339 EGWithMaskChangeAttributes(COLOR_MASK | LAYER_MASK | VECMTYPE_MASK, go);
1340 EMAddGraphicsToModel(ESIModel(), go);
1341
1342 p [ 0 ].x = ( FPNum ) loc.at(j);
1343 p [ 0 ].y = ( FPNum ) loc.at(i);
1344 p [ 0 ].z = 0.;
1345
1346 EASValsSetMType(SQUARE_MARKER);
1347 go = CreateMarker3D(p);
1348 EGWithMaskChangeAttributes(COLOR_MASK | LAYER_MASK | VECMTYPE_MASK, go);
1349 EMAddGraphicsToModel(ESIModel(), go);
1350 }
1351 }
1352 }
1353 }
1354}
1355
1356void
1357StructuralElement :: showExtendedSparseMtrxStructure(CharType mtrx, oofegGraphicContext &gc, TimeStep *tStep)
1358{
1359 if ( mtrx == TangentStiffnessMatrix ) {
1360 // loop over element IP
1361 for ( IntegrationPoint *ip : * this->giveDefaultIntegrationRulePtr() ) {
1362 auto interface = static_cast< NonlocalMaterialStiffnessInterface * >( this->giveStructuralCrossSection()->
1363 giveMaterialInterface(NonlocalMaterialStiffnessInterfaceType, ip) );
1364
1365 if ( interface == NULL ) {
1366 return;
1367 }
1369 }
1370 }
1371}
1372
1373#endif
1374} // end namespace oofem
#define N(a, b)
const entryListType & getElementRecords(int elem)
Definition bctracker.C:133
std::list< Entry > entryListType
Definition bctracker.h:67
void computeValueAt(FloatArray &answer, TimeStep *tStep, const FloatArray &coords, ValueModeType mode) override
int giveApproxOrder() override=0
CoordSystType giveCoordSystMode() override
double giveCoordinate(int i) const
Definition dofmanager.h:383
Node * giveNode(int i) const
Definition element.h:629
virtual int computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords)
Definition element.C:1225
virtual bool isActivated(TimeStep *tStep)
Definition element.C:838
virtual FEInterpolation * giveInterpolation() const
Definition element.h:648
virtual int giveSpatialDimension()
Definition element.C:1391
virtual std::unique_ptr< IntegrationRule > giveBoundaryEdgeIntegrationRule(int order, int boundary)
Definition element.C:911
virtual int computeNumberOfDofs()
Definition element.h:436
virtual void giveDofManDofIDMask(int inode, IntArray &answer) const
Definition element.h:488
IntArray * giveBodyLoadArray()
Returns array containing load numbers of loads acting on element.
Definition element.C:411
virtual MaterialMode giveMaterialMode()
Definition element.h:738
virtual std::unique_ptr< IntegrationRule > giveBoundarySurfaceIntegrationRule(int order, int boundary)
Definition element.C:917
void giveLocationArray(IntArray &locationArray, const UnknownNumberingScheme &s, IntArray *dofIds=NULL) const
Definition element.C:429
virtual int giveIntegrationRuleLocalCodeNumbers(IntArray &answer, IntegrationRule &ie)
Definition element.h:747
virtual integrationDomain giveIntegrationDomain() const
Definition element.C:1693
int activityTimeFunction
Element activity time function. If defined, nonzero value indicates active receiver,...
Definition element.h:163
void computeVectorOf(ValueModeType u, TimeStep *tStep, FloatArray &answer)
Definition element.C:103
virtual IntegrationRule * giveIntegrationRule(int i)
Definition element.h:899
std::vector< std ::unique_ptr< IntegrationRule > > integrationRulesArray
Definition element.h:157
elementParallelMode giveParallelMode() const
Definition element.h:1139
virtual int giveNumberOfDofManagers() const
Definition element.h:695
CrossSection * giveCrossSection()
Definition element.C:534
virtual IntegrationRule * giveDefaultIntegrationRulePtr()
Definition element.h:886
Element(int n, Domain *aDomain)
Definition element.C:86
virtual bool computeLocalCoordinates(FloatArray &answer, const FloatArray &gcoords)
Definition element.C:1240
virtual double computeVolumeAround(GaussPoint *gp)
Definition element.h:530
IntArray bodyLoadArray
Definition element.h:147
virtual Element_Geometry_Type giveGeometryType() const =0
virtual void boundaryEdgeLocal2Global(FloatArray &answer, int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const =0
virtual double boundarySurfaceGiveTransformationJacobian(int isurf, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const =0
virtual double boundaryEdgeGiveTransformationJacobian(int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const =0
virtual void boundaryLocal2Global(FloatArray &answer, int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const =0
Domain * domain
Link to domain object, useful for communicating with other FEM components.
Definition femcmpnn.h:79
int number
Component number.
Definition femcmpnn.h:77
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
Index giveSize() const
Returns the size of receiver.
Definition floatarray.h:261
void zero()
Zeroes all coefficients of receiver.
Definition floatarray.C:683
void rotatedWith(FloatMatrix &r, char mode)
Definition floatarray.C:814
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 subtract(const FloatArray &src)
Definition floatarray.C:320
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 beSubMatrixOf(const FloatMatrix &src, Index topRow, Index bottomRow, Index topCol, Index bottomCol)
void plusProductUnsym(const FloatMatrix &a, const FloatMatrix &b, double dV)
int giveNumberOfColumns() const
Returns number of columns of receiver.
void beProductOf(const FloatMatrix &a, const FloatMatrix &b)
void beNMatrixOf(const FloatArray &n, int nsd)
void zero()
Zeroes all coefficient of receiver.
int giveNumberOfRows() const
Returns number of rows of receiver.
void assemble(const FloatMatrix &src, const IntArray &loc)
double at(std::size_t i, std::size_t j) const
bool isSquare() const
Returns nonzero if receiver is square matrix.
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
const FloatArray & giveSubPatchCoordinates() const
Returns local sub-patch coordinates of the receiver.
Definition gausspoint.h:142
IntegrationPointStatus * giveMaterialStatus(IntegrationPointStatusIDType key=IPSID_Default)
Definition gausspoint.h:204
double giveWeight()
Returns integration weight of receiver.
Definition gausspoint.h:180
virtual bcGeomType giveBCGeoType() const
virtual bool isImposed(TimeStep *tStep)
void followedBy(const IntArray &b, int allocChunk=0)
Definition intarray.C:94
bool isEmpty() const
Definition intarray.h:217
int & at(std::size_t i)
Definition intarray.h:104
int giveSize() const
Definition intarray.h:211
virtual int getRequiredNumberOfIntegrationPoints(integrationDomain dType, int approxOrder)
int setUpIntegrationPoints(integrationDomain intdomain, int nPoints, MaterialMode matMode)
const FloatArrayF< 6 > & giveLatticeStress() const
Returns lattice stress.
virtual void computeValueAt(FloatArray &answer, TimeStep *tStep, const FloatArray &coords, ValueModeType mode)=0
virtual void computeComponentArrayAt(FloatArray &answer, TimeStep *tStep, ValueModeType mode)
Definition load.C:84
virtual FormulationType giveFormulationType()
Definition load.h:170
virtual double giveUpdatedCoordinate(int ic, TimeStep *tStep, double scale=1.)
Definition node.C:234
virtual void NonlocalMaterialStiffnessInterface_addIPContribution(SparseMtrx &dest, const UnknownNumberingScheme &s, GaussPoint *gp, TimeStep *tStep)=0
Computes and adds IP contributions to destination matrix.
virtual void NonlocalMaterialStiffnessInterface_showSparseMtrxStructure(GaussPoint *gp, oofegGraphicContext &gc, TimeStep *)
const FloatArray & giveCoordinates() const
Definition pointload.h:85
void computeValueAt(FloatArray &answer, TimeStep *tStep, const FloatArray &coords, ValueModeType mode) override
Definition pointload.C:45
CoordSystType giveCoordSystMode() override
Definition pointload.h:87
virtual void createMaterialStatus(GaussPoint &iGP)=0
virtual int computeLoadLSToLRotationMatrix(FloatMatrix &answer, int iSurf, GaussPoint *gp)
virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)=0
virtual void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep)
virtual void giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord=0)
virtual void computeInitialStressMatrix(FloatMatrix &answer, TimeStep *tStep)
virtual void computeSurfaceNMatrix(FloatMatrix &answer, int boundaryID, const FloatArray &lcoords)
StructuralCrossSection * giveStructuralCrossSection()
Helper function which returns the structural cross-section for the element.
virtual void computeConsistentMassMatrix(FloatMatrix &answer, TimeStep *tStep, double &mass, const double *ipDensity=NULL)
virtual void computePointLoadVectorAt(FloatArray &answer, Load *load, TimeStep *tStep, ValueModeType mode, bool global=true)
virtual void setupIRForMassMtrxIntegration(IntegrationRule &iRule)
virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
virtual void giveMassMtrxIntegrationgMask(IntArray &answer)
virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer)
virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int lowerIndx=1, int upperIndx=ALL_STRAINS)=0
virtual void computeLumpedInitialStressMatrix(FloatMatrix &answer, TimeStep *tStep)
virtual double computeSurfaceVolumeAround(GaussPoint *gp, int iSurf)
virtual void computeBodyLoadVectorAt(FloatArray &answer, Load *load, TimeStep *tStep, ValueModeType mode)
std::unique_ptr< FloatArray > initialDisplacements
Initial displacement vector, describes the initial nodal displacements when element has been casted.
virtual void computeStrainVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep)
virtual double computeEdgeVolumeAround(GaussPoint *gp, int iEdge)
virtual int giveNumberOfIPForMassMtrxIntegration()
virtual int computeLoadLEToLRotationMatrix(FloatMatrix &answer, int iEdge, GaussPoint *gp)
virtual void computeMassMatrix(FloatMatrix &answer, TimeStep *tStep)
virtual int computeLoadGToLRotationMtrx(FloatMatrix &answer)
virtual void computeEdgeNMatrix(FloatMatrix &answer, int boundaryID, const FloatArray &lcoords)
computes edge interpolation matrix
const char * giveClassName() const override
virtual void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep)=0
const FloatArrayF< 3 > & giveTraction() const
Returns the const pointer to receiver's traction vector.
const FloatArray & giveStressVector() const
Returns the const pointer to receiver's stress vector.
virtual void updateBeforeNonlocAverage(const FloatArray &strainVector, GaussPoint *gp, TimeStep *tStep) const =0
#define OOFEM_WARNING(...)
Definition error.h:80
#define OOFEM_ERROR(...)
Definition error.h:79
@ CS_StructuralCapability
Structural capability.
@ Element_remote
Element in active domain is only mirror of some remote element.
Definition element.h:89
GaussPoint IntegrationPoint
Definition gausspoint.h:272
@ BodyLoadBGT
Distributed body load.
Definition bcgeomtype.h:43
const char * __CharTypeToString(CharType _value)
Definition cltypes.C:338
InternalStateMode
Determines the mode of internal variable.
@ ForceLoadBVT
Definition bcvaltype.h:43
@ EigenstrainBVT
Definition bcvaltype.h:48
@ TemperatureBVT
Definition bcvaltype.h:42
@ NonlocalMaterialStiffnessInterfaceType
@ MaterialModelMapperInterfaceType
@ NonlocalMaterialExtensionInterfaceType
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
#define OOFEG_SPARSE_PROFILE_WIDTH
#define OOFEG_SPARSE_PROFILE_LAYER

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