OOFEM 3.0
Loading...
Searching...
No Matches
trplanestressrotallman.C
Go to the documentation of this file.
1/*
2 *
3 * ##### ##### ###### ###### ### ###
4 * ## ## ## ## ## ## ## ### ##
5 * ## ## ## ## #### #### ## # ##
6 * ## ## ## ## ## ## ## ##
7 * ## ## ## ## ## ## ## ##
8 * ##### ##### ## ###### ## ##
9 *
10 *
11 * OOFEM : Object Oriented Finite Element Code
12 *
13 * Copyright (C) 1993 - 2025 Borek Patzak
14 *
15 *
16 *
17 * Czech Technical University, Faculty of Civil Engineering,
18 * Department of Structural Mechanics, 166 29 Prague, Czech Republic
19 *
20 * This library is free software; you can redistribute it and/or
21 * modify it under the terms of the GNU Lesser General Public
22 * License as published by the Free Software Foundation; either
23 * version 2.1 of the License, or (at your option) any later version.
24 *
25 * This program is distributed in the hope that it will be useful,
26 * but WITHOUT ANY WARRANTY; without even the implied warranty of
27 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
28 * Lesser General Public License for more details.
29 *
30 * You should have received a copy of the GNU Lesser General Public
31 * License along with this library; if not, write to the Free Software
32 * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
33 */
34
36#include "fei2dtrquad.h"
37#include "fei2dtrlin.h"
38#include "node.h"
39#include "gausspoint.h"
41#include "floatmatrix.h"
42#include "floatarray.h"
43#include "intarray.h"
44#include "mathfem.h"
45#include "boundaryload.h"
46#include "classfactory.h"
47
48#ifdef __OOFEG
49 #include "oofeggraphiccontext.h"
50 #include "oofegutils.h"
51#endif
52
53namespace oofem {
55
56FEI2dTrQuad TrPlanestressRotAllman :: qinterpolation(1, 2);
57
58TrPlanestressRotAllman :: TrPlanestressRotAllman(int n, Domain *aDomain) :
59 TrPlaneStress2d(n, aDomain)
60{
63}
64
66TrPlanestressRotAllman :: giveInterface(InterfaceType interface)
67{
68 if ( interface == LayeredCrossSectionInterfaceType ) {
69 return static_cast< LayeredCrossSectionInterface * >(this);
70 } else if ( interface == ZZNodalRecoveryModelInterfaceType ) {
71 return static_cast< ZZNodalRecoveryModelInterface * >(this);
72 } else if ( interface == SPRNodalRecoveryModelInterfaceType ) {
73 return static_cast< SPRNodalRecoveryModelInterface * >(this);
74 } else if ( interface == SpatialLocalizerInterfaceType ) {
75 return static_cast< SpatialLocalizerInterface * >(this);
76
77 }
78 return NULL;
79}
80
81void
82TrPlanestressRotAllman :: computeLocalNodalCoordinates(std::vector< FloatArray > &lxy)
83{
84 lxy.resize(6);
85 for ( int i = 0; i < 3; i++ ) {
86 lxy [ i ] = this->giveNode(i + 1)->giveCoordinates();
87 }
88 lxy [ 3 ].resize(2);
89 lxy [ 4 ].resize(2);
90 lxy [ 5 ].resize(2);
91 for ( int i = 1; i <= 2; i++ ) {
92 lxy [ 3 ].at(i) = 0.5 * ( lxy [ 0 ].at(i) + lxy [ 1 ].at(i) );
93 lxy [ 4 ].at(i) = 0.5 * ( lxy [ 1 ].at(i) + lxy [ 2 ].at(i) );
94 lxy [ 5 ].at(i) = 0.5 * ( lxy [ 2 ].at(i) + lxy [ 0 ].at(i) );
95 }
96}
97
98
99
100void
101TrPlanestressRotAllman :: computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer)
102// Returns the displacement interpolation matrix {N} of the receiver, eva-
103// luated at gp.
104{
105 FloatArray L, n;
106 std::vector< FloatArray > lxy;
107
108 answer.resize(3, 9);
109 answer.zero();
110
111 this->computeLocalNodalCoordinates(lxy); // get ready for tranformation into 3d
112 this->qinterpolation.evalN( n, iLocCoord, FEIVertexListGeometryWrapper(lxy, this->qinterpolation.giveGeometryType()) );
113 this->interp.evalN( L, iLocCoord, FEIVertexListGeometryWrapper(lxy, this->interp.giveGeometryType()));
114
115 answer.at(1, 1) = answer.at(2, 2) = n.at(1) + n.at(4) / 2. + n.at(6) / 2.;
116 answer.at(1, 4) = answer.at(2, 5) = n.at(2) + n.at(4) / 2. + n.at(5) / 2.;
117 answer.at(1, 7) = answer.at(2, 8) = n.at(3) + n.at(5) / 2. + n.at(6) / 2.;
118 answer.at(1, 3) = n.at(6) * ( lxy [ 0 ].at(2) - lxy [ 2 ].at(2) ) / 8.0 - n.at(4) * ( lxy [ 1 ].at(2) - lxy [ 0 ].at(2) ) / 8.0;
119 answer.at(1, 6) = n.at(4) * ( lxy [ 1 ].at(2) - lxy [ 0 ].at(2) ) / 8.0 - n.at(5) * ( lxy [ 2 ].at(2) - lxy [ 1 ].at(2) ) / 8.0;
120 answer.at(1, 9) = n.at(5) * ( lxy [ 2 ].at(2) - lxy [ 1 ].at(2) ) / 8.0 - n.at(6) * ( lxy [ 0 ].at(2) - lxy [ 2 ].at(2) ) / 8.0;
121 answer.at(2, 3) = -n.at(6) * ( lxy [ 0 ].at(1) - lxy [ 2 ].at(1) ) / 8.0 + n.at(4) * ( lxy [ 1 ].at(1) - lxy [ 0 ].at(1) ) / 8.0;
122 answer.at(2, 6) = -n.at(4) * ( lxy [ 1 ].at(1) - lxy [ 0 ].at(1) ) / 8.0 + n.at(5) * ( lxy [ 2 ].at(1) - lxy [ 1 ].at(1) ) / 8.0;
123 answer.at(2, 9) = -n.at(5) * ( lxy [ 2 ].at(1) - lxy [ 1 ].at(1) ) / 8.0 + n.at(6) * ( lxy [ 0 ].at(1) - lxy [ 2 ].at(1) ) / 8.0;
124 // linear approx for rotations
125 answer.at(3, 3) = L.at(1);
126 answer.at(3, 6) = L.at(2);
127 answer.at(3, 9) = L.at(3);
128}
129
130void
131TrPlanestressRotAllman :: computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int li, int ui)
132// Returns the [3x12] strain-displacement matrix {B} of the receiver, eva-
133// luated at gp.
134{
135 FloatMatrix dnx;
136 std::vector< FloatArray > lxy;
137
138 this->computeLocalNodalCoordinates(lxy); // get ready for tranformation into 3d
139 this->qinterpolation.evaldNdx( dnx, gp->giveNaturalCoordinates(), FEIVertexListGeometryWrapper(lxy, this->qinterpolation.giveGeometryType()) );
140
141 answer.resize(3, 9);
142 answer.zero();
143
144 // epsilon_x
145 answer.at(1, 1) = dnx.at(1, 1) + 0.5 * dnx.at(4, 1) + 0.5 * dnx.at(6, 1);
146 answer.at(1, 4) = dnx.at(2, 1) + 0.5 * dnx.at(4, 1) + 0.5 * dnx.at(5, 1);
147 answer.at(1, 7) = dnx.at(3, 1) + 0.5 * dnx.at(5, 1) + 0.5 * dnx.at(6, 1);
148 answer.at(1, 3) =+dnx.at(6, 1) * ( lxy [ 0 ].at(2) - lxy [ 2 ].at(2) ) / 8.0 - dnx.at(4, 1) * ( lxy [ 1 ].at(2) - lxy [ 0 ].at(2) ) / 8.0;
149 answer.at(1, 6) =+dnx.at(4, 1) * ( lxy [ 1 ].at(2) - lxy [ 0 ].at(2) ) / 8.0 - dnx.at(5, 1) * ( lxy [ 2 ].at(2) - lxy [ 1 ].at(2) ) / 8.0;
150 answer.at(1, 9) =+dnx.at(5, 1) * ( lxy [ 2 ].at(2) - lxy [ 1 ].at(2) ) / 8.0 - dnx.at(6, 1) * ( lxy [ 0 ].at(2) - lxy [ 2 ].at(2) ) / 8.0;
151
152 // epsilon_y
153 answer.at(2, 2) = dnx.at(1, 2) + 0.5 * dnx.at(4, 2) + 0.5 * dnx.at(6, 2);
154 answer.at(2, 5) = dnx.at(2, 2) + 0.5 * dnx.at(4, 2) + 0.5 * dnx.at(5, 2);
155 answer.at(2, 8) = dnx.at(3, 2) + 0.5 * dnx.at(5, 2) + 0.5 * dnx.at(6, 2);
156 answer.at(2, 3) =-dnx.at(6, 2) * ( lxy [ 0 ].at(1) - lxy [ 2 ].at(1) ) / 8.0 + dnx.at(4, 2) * ( lxy [ 1 ].at(1) - lxy [ 0 ].at(1) ) / 8.0;
157 answer.at(2, 6) =-dnx.at(4, 2) * ( lxy [ 1 ].at(1) - lxy [ 0 ].at(1) ) / 8.0 + dnx.at(5, 2) * ( lxy [ 2 ].at(1) - lxy [ 1 ].at(1) ) / 8.0;
158 answer.at(2, 9) =-dnx.at(5, 2) * ( lxy [ 2 ].at(1) - lxy [ 1 ].at(1) ) / 8.0 + dnx.at(6, 2) * ( lxy [ 0 ].at(1) - lxy [ 2 ].at(1) ) / 8.0;
159
160 // gamma_xy (shear)
161 answer.at(3, 1) = dnx.at(1, 2) + 0.5 * dnx.at(4, 2) + 0.5 * dnx.at(6, 2);
162 answer.at(3, 2) = dnx.at(1, 1) + 0.5 * dnx.at(4, 1) + 0.5 * dnx.at(6, 1);
163 answer.at(3, 4) = dnx.at(2, 2) + 0.5 * dnx.at(4, 2) + 0.5 * dnx.at(5, 2);
164 answer.at(3, 5) = dnx.at(2, 1) + 0.5 * dnx.at(4, 1) + 0.5 * dnx.at(5, 1);
165 answer.at(3, 7) = dnx.at(3, 2) + 0.5 * dnx.at(5, 2) + 0.5 * dnx.at(6, 2);
166 answer.at(3, 8) = dnx.at(3, 1) + 0.5 * dnx.at(5, 1) + 0.5 * dnx.at(6, 1);
167
168 answer.at(3, 3) = dnx.at(6, 2) * ( lxy [ 0 ].at(2) - lxy [ 2 ].at(2) ) / 8.0 - dnx.at(4, 2) * ( lxy [ 1 ].at(2) - lxy [ 0 ].at(2) ) / 8.0;
169 answer.at(3, 3) += -dnx.at(6, 1) * ( lxy [ 0 ].at(1) - lxy [ 2 ].at(1) ) / 8.0 + dnx.at(4, 1) * ( lxy [ 1 ].at(1) - lxy [ 0 ].at(1) ) / 8.0;
170 answer.at(3, 6) = dnx.at(4, 2) * ( lxy [ 1 ].at(2) - lxy [ 0 ].at(2) ) / 8.0 - dnx.at(5, 2) * ( lxy [ 2 ].at(2) - lxy [ 1 ].at(2) ) / 8.0;
171 answer.at(3, 6) += -dnx.at(4, 1) * ( lxy [ 1 ].at(1) - lxy [ 0 ].at(1) ) / 8.0 + dnx.at(5, 1) * ( lxy [ 2 ].at(1) - lxy [ 1 ].at(1) ) / 8.0;
172 answer.at(3, 9) = dnx.at(5, 2) * ( lxy [ 2 ].at(2) - lxy [ 1 ].at(2) ) / 8.0 - dnx.at(6, 2) * ( lxy [ 0 ].at(2) - lxy [ 2 ].at(2) ) / 8.0;
173 answer.at(3, 9) += -dnx.at(5, 1) * ( lxy [ 2 ].at(1) - lxy [ 1 ].at(1) ) / 8.0 + dnx.at(6, 1) * ( lxy [ 0 ].at(1) - lxy [ 2 ].at(1) ) / 8.0;
174}
175
176
177void
178TrPlanestressRotAllman :: computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
179{
180 // compute standard stiffness matrix
181 TrPlaneStress2d :: computeStiffnessMatrix(answer, rMode, tStep);
182 // add zero energy mode stabilization
183 FloatMatrix ks;
184 this->computeStiffnessMatrixZeroEnergyStabilization(ks, rMode, tStep);
185 answer.add(ks);
186}
187
188void
189TrPlanestressRotAllman :: computeStiffnessMatrixZeroEnergyStabilization(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
190{
191 FloatMatrix b(1, 9);
192 FloatMatrix dnx;
193 FloatArray lec = Vec3(0.333333333333, 0.333333333333, 0.333333333333); // element center in local coordinates
194 std::vector< FloatArray > lxy;
195
196 this->computeLocalNodalCoordinates(lxy); // get ready for tranformation into 3d
197 this->qinterpolation.evaldNdx( dnx, lec, FEIVertexListGeometryWrapper(lxy, this->qinterpolation.giveGeometryType()) );
198
199 // evaluate (dv/dx-du/dy)/2. at element center
200 b.at(1, 1) = -1.0 * ( dnx.at(1, 2) + 0.5 * dnx.at(4, 2) + 0.5 * dnx.at(6, 2) );
201 b.at(1, 2) = dnx.at(1, 1) + 0.5 * dnx.at(4, 1) + 0.5 * dnx.at(6, 1);
202 b.at(1, 4) = -1.0 * ( dnx.at(2, 2) + 0.5 * dnx.at(4, 2) + 0.5 * dnx.at(5, 2) );
203 b.at(1, 5) = dnx.at(2, 1) + 0.5 * dnx.at(4, 1) + 0.5 * dnx.at(5, 1);
204 b.at(1, 7) = -1.0 * ( dnx.at(3, 2) + 0.5 * dnx.at(5, 2) + 0.5 * dnx.at(6, 2) );
205 b.at(1, 8) = dnx.at(3, 1) + 0.5 * dnx.at(5, 1) + 0.5 * dnx.at(6, 1);
206
207 b.at(1, 3) = -dnx.at(6, 2) * ( lxy [ 0 ].at(2) - lxy [ 2 ].at(2) ) / 8.0 + dnx.at(4, 2) * ( lxy [ 1 ].at(2) - lxy [ 0 ].at(2) ) / 8.0;
208 b.at(1, 3) += -dnx.at(6, 1) * ( lxy [ 0 ].at(1) - lxy [ 2 ].at(1) ) / 8.0 + dnx.at(4, 1) * ( lxy [ 1 ].at(1) - lxy [ 0 ].at(1) ) / 8.0;
209 b.at(1, 6) = -dnx.at(4, 2) * ( lxy [ 1 ].at(2) - lxy [ 0 ].at(2) ) / 8.0 + dnx.at(5, 2) * ( lxy [ 2 ].at(2) - lxy [ 1 ].at(2) ) / 8.0;
210 b.at(1, 6) += -dnx.at(4, 1) * ( lxy [ 1 ].at(1) - lxy [ 0 ].at(1) ) / 8.0 + dnx.at(5, 1) * ( lxy [ 2 ].at(1) - lxy [ 1 ].at(1) ) / 8.0;
211 b.at(1, 9) = -dnx.at(5, 2) * ( lxy [ 2 ].at(2) - lxy [ 1 ].at(2) ) / 8.0 + dnx.at(6, 2) * ( lxy [ 0 ].at(2) - lxy [ 2 ].at(2) ) / 8.0;
212 b.at(1, 9) += -dnx.at(5, 1) * ( lxy [ 2 ].at(1) - lxy [ 1 ].at(1) ) / 8.0 + dnx.at(6, 1) * ( lxy [ 0 ].at(1) - lxy [ 2 ].at(1) ) / 8.0;
213 b.times(0.5);
214 // add -1.0*sum(r_w)/3.0
215 b.at(1, 3) -= 1.0 / 3.0;
216 b.at(1, 6) -= 1.0 / 3.0;
217 b.at(1, 9) -= 1.0 / 3.0;
218 // add alpha*Volume*B^T[G]B to element stiffness matrix
219 double G = this->giveStructuralCrossSection()->give( Gxy, this->giveDefaultIntegrationRulePtr()->getIntegrationPoint(0) );
220 double coeff = G * this->giveArea() * this->giveCrossSection()->give( CS_Thickness, this->giveDefaultIntegrationRulePtr()->getIntegrationPoint(0) ) * 1.e-6;
221 answer.beTProductOf(b, b);
222 answer.times(coeff);
223}
224
225
226
227void
228TrPlanestressRotAllman :: giveDofManDofIDMask(int inode, IntArray &answer) const
229{
230 answer = {D_u, D_v, R_w};
231}
232
233
234
235double
236TrPlanestressRotAllman :: giveArea()
237// returns the area occupied by the receiver
238{
239 if ( area > 0 ) {
240 return area; // check if previously computed
241 }
242
243 return ( area = fabs( this->interp.giveArea( FEIElementGeometryWrapper(this) ) ) );
244}
245
246void TrPlanestressRotAllman :: computeGaussPoints()
247// Sets up the array containing the four Gauss points of the receiver.
248{
249 if ( integrationRulesArray.size() == 0 ) {
250 integrationRulesArray.resize( 1 );
251 integrationRulesArray [ 0 ] = std::make_unique<GaussIntegrationRule>(1, this, 1, 3);
252 this->giveCrossSection()->setupIntegrationPoints(* integrationRulesArray [ 0 ], numberOfGaussPoints, this);
253 }
254}
255
256void
257TrPlanestressRotAllman :: computeEgdeNMatrixAt(FloatMatrix &answer, int iedge, GaussPoint *gp)
258{
259 std::vector< FloatArray > lxy;
260 FloatArray l, n;
261 FEI2dTrQuad qi(1, 2);
262
263 this->computeLocalNodalCoordinates(lxy); // get ready for tranformation into 3d
265 const auto &en = qi.computeLocalEdgeMapping(iedge); // get edge mapping
266 this->interp.edgeEvalN( l, iedge, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(this) );
267 answer.resize(3, 6);
268
269 answer.at(1, 1) = answer.at(2, 2) = n.at(1) + n.at(3) / 2.0;
270 answer.at(1, 4) = answer.at(2, 5) = n.at(2) + n.at(3) / 2.0;
271 answer.at(1, 3) = n.at(3) * ( lxy [ en.at(2) - 1 ].at(2) - lxy [ en.at(1) - 1 ].at(2) ) / 8.0;
272 answer.at(1, 6) = -answer.at(1, 3);
273 answer.at(2, 3) = n.at(3) * ( lxy [ en.at(2) - 1 ].at(1) - lxy [ en.at(1) - 1 ].at(1) ) / 8.0;
274 answer.at(2, 6) = -answer.at(2, 3);
275 answer.at(3, 3) = l.at(1);
276 answer.at(3, 6) = l.at(2);
277}
278
279void
280TrPlanestressRotAllman :: giveEdgeDofMapping(IntArray &answer, int iEdge) const
281{
282 /*
283 * provides dof mapping of local edge dofs (only nonzero are taken into account)
284 * to global element dofs
285 */
286
287 answer.resize(6);
288 if ( iEdge == 1 ) { // edge between nodes 1,2
289 answer.at(1) = 1;
290 answer.at(2) = 2;
291 answer.at(3) = 3;
292 answer.at(4) = 4;
293 answer.at(5) = 5;
294 answer.at(6) = 6;
295 } else if ( iEdge == 2 ) { // edge between nodes 2 3
296 answer.at(1) = 4;
297 answer.at(2) = 5;
298 answer.at(3) = 6;
299 answer.at(4) = 7;
300 answer.at(5) = 8;
301 answer.at(6) = 9;
302 } else if ( iEdge == 3 ) { // edge between nodes 3 1
303 answer.at(1) = 7;
304 answer.at(2) = 8;
305 answer.at(3) = 9;
306 answer.at(4) = 1;
307 answer.at(5) = 2;
308 answer.at(6) = 3;
309 } else {
310 OOFEM_ERROR("wrong edge number");
311 }
312}
313
314void TrPlanestressRotAllman :: computeBoundaryEdgeLoadVector(FloatArray &answer, BoundaryLoad *load, int boundary, CharType type, ValueModeType mode, TimeStep *tStep, bool global)
315{
316 answer.clear();
317 if ( type != ExternalForcesVector ) {
318 return;
319 }
320
321 double dV;
322 FloatMatrix T;
323 FloatArray globalIPcoords;
324
325 EdgeLoad *edgeLoad = dynamic_cast< EdgeLoad * >(load);
326 if ( edgeLoad ) {
327 int approxOrder = edgeLoad->giveApproxOrder() + this->giveInterpolation()->giveInterpolationOrder();
328 int numberOfGaussPoints = ( int ) ceil( ( approxOrder + 1. ) / 2. );
329 GaussIntegrationRule iRule(1, this, 1, 1);
330 iRule.SetUpPointsOnLine(numberOfGaussPoints, _Unknown);
331 FloatArray reducedAnswer, force, ntf;
332 IntArray mask;
333 FloatMatrix n;
334
335 for ( GaussPoint *gp: iRule ) {
336 this->computeEgdeNMatrixAt(n, boundary, gp);
337 dV = this->computeEdgeVolumeAround(gp, boundary);
338
339 if ( edgeLoad->giveFormulationType() == Load :: FT_Entity ) {
340 edgeLoad->computeValueAt(force, tStep, gp->giveNaturalCoordinates(), mode);
341 } else {
342 //this->computeEdgeIpGlobalCoords(globalIPcoords, gp, boundary);
343 this->giveInterpolation()->boundaryEdgeLocal2Global( globalIPcoords, boundary, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(this) );
344 edgeLoad->computeValueAt(force, tStep, globalIPcoords, mode);
345 }
346
347 // transform force
348 if ( edgeLoad->giveCoordSystMode() == Load :: CST_Global ) {
349 } else {
350 // transform from local boundary to element local c.s
351 if ( this->computeLoadLEToLRotationMatrix(T, boundary, gp) ) {
352 force.rotatedWith(T, 'n');
353 }
354 // then to global c.s
355 if ( this->computeLoadGToLRotationMtrx(T) ) {
356 force.rotatedWith(T, 't');
357 }
358 }
359
360 ntf.beTProductOf(n, force);
361 answer.add(dV, ntf);
362 }
363
364
365 return;
366 } else {
367 OOFEM_ERROR("incompatible load");
368 return;
369 }
370
371
372}
373
374
375
376/*
377 * double
378 * TrPlanestressRotAllman :: computeEdgeVolumeAround(GaussPoint *gp, int iEdge)
379 * {
380 * // edge with linear geometry -> one can use linear interpolation safely
381 * double detJ = this->interp.edgeGiveTransformationJacobian( iEdge, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(this) );
382 * return detJ *gp->giveWeight();
383 * }
384 *
385 *
386 * void
387 * TrPlanestressRotAllman :: computeEdgeIpGlobalCoords(FloatArray &answer, GaussPoint *gp, int iEdge)
388 * {
389 * // edge with linear geometry -> one can use linear interpolation safely
390 * this->interp.edgeLocal2global( answer, iEdge, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(this) );
391 * }
392 *
393 *
394 * int
395 * TrPlanestressRotAllman :: computeLoadLEToLRotationMatrix(FloatMatrix &answer, int iEdge, GaussPoint *gp)
396 * {
397 * // returns transformation matrix from
398 * // edge local coordinate system
399 * // to element local c.s
400 * // (same as global c.s in this case)
401 * //
402 * // i.e. f(element local) = T * f(edge local)
403 * //
404 * double dx, dy, length;
405 * Node *nodeA, *nodeB;
406 * int aNode = 0, bNode = 0;
407 *
408 * answer.resize(2, 2);
409 * answer.zero();
410 *
411 * if ( iEdge == 1 ) { // edge between nodes 1 2
412 * aNode = 1;
413 * bNode = 2;
414 * } else if ( iEdge == 2 ) { // edge between nodes 2 3
415 * aNode = 2;
416 * bNode = 3;
417 * } else if ( iEdge == 3 ) { // edge between nodes 2 3
418 * aNode = 3;
419 * bNode = 1;
420 * } else {
421 * OOFEM_ERROR("wrong egde number");
422 * }
423 *
424 * nodeA = this->giveNode(aNode);
425 * nodeB = this->giveNode(bNode);
426 *
427 * dx = nodeB->giveCoordinate(1) - nodeA->giveCoordinate(1);
428 * dy = nodeB->giveCoordinate(2) - nodeA->giveCoordinate(2);
429 * length = sqrt(dx * dx + dy * dy);
430 *
431 * answer.at(1, 1) = dx / length;
432 * answer.at(1, 2) = -dy / length;
433 * answer.at(2, 1) = dy / length;
434 * answer.at(2, 2) = dx / length;
435 *
436 * return 1;
437 * }
438 *
439 *
440 *
441 *
442 * double TrPlanestressRotAllman :: computeVolumeAround(GaussPoint *gp)
443 * // Returns the portion of the receiver which is attached to gp.
444 * {
445 * double detJ, weight;
446 *
447 * weight = gp->giveWeight();
448 * // safe to use linear interpolation here (geometry is linear)
449 * detJ = fabs( this->interp.giveTransformationJacobian( gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(this) ) );
450 *
451 * return detJ *weight *this->giveCrossSection()->give(CS_Thickness), gp;
452 * }
453 */
454
455
456
457int
458TrPlanestressRotAllman :: SPRNodalRecoveryMI_giveNumberOfIP()
459{
460 return this->numberOfGaussPoints;
461}
462} // end namespace oofem
#define REGISTER_Element(class)
void computeValueAt(FloatArray &answer, TimeStep *tStep, const FloatArray &coords, ValueModeType mode) override
int giveApproxOrder() override=0
CoordSystType giveCoordSystMode() override
Node * giveNode(int i) const
Definition element.h:629
int numberOfDofMans
Number of dofmanagers.
Definition element.h:136
std::vector< std ::unique_ptr< IntegrationRule > > integrationRulesArray
Definition element.h:157
int numberOfGaussPoints
Definition element.h:175
CrossSection * giveCrossSection()
Definition element.C:534
virtual IntegrationRule * giveDefaultIntegrationRulePtr()
Definition element.h:886
void edgeEvalN(FloatArray &answer, int iedge, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const override
const Element_Geometry_Type giveGeometryType() const override
Definition fei2dtrquad.h:50
IntArray computeLocalEdgeMapping(int iedge) const override
double & at(Index i)
Definition floatarray.h:202
void rotatedWith(FloatMatrix &r, char mode)
Definition floatarray.C:814
void beTProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Definition floatarray.C:721
void add(const FloatArray &src)
Definition floatarray.C:218
void times(double f)
void add(const FloatMatrix &a)
void resize(Index rows, Index cols)
Definition floatmatrix.C:79
void zero()
Zeroes all coefficient of receiver.
double at(std::size_t i, std::size_t j) const
void beTProductOf(const FloatMatrix &a, const FloatMatrix &b)
int SetUpPointsOnLine(int nPoints, MaterialMode mode) override
const FloatArray & giveNaturalCoordinates() const
Returns coordinate array of receiver.
Definition gausspoint.h:138
void resize(int n)
Definition intarray.C:73
int & at(std::size_t i)
Definition intarray.h:104
virtual FormulationType giveFormulationType()
Definition load.h:170
int computeLoadLEToLRotationMatrix(FloatMatrix &answer, int iEdge, GaussPoint *gp) override
double computeEdgeVolumeAround(GaussPoint *gp, int iEdge) override
StructuralCrossSection * giveStructuralCrossSection()
Helper function which returns the structural cross-section for the element.
virtual int computeLoadGToLRotationMtrx(FloatMatrix &answer)
TrPlaneStress2d(int n, Domain *d)
Definition trplanstrss.C:58
static FEI2dTrLin interp
Definition trplanstrss.h:70
FEInterpolation * giveInterpolation() const override
Definition trplanstrss.C:70
void computeEgdeNMatrixAt(FloatMatrix &answer, int iedge, GaussPoint *gp)
void computeStiffnessMatrixZeroEnergyStabilization(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
virtual void computeLocalNodalCoordinates(std::vector< FloatArray > &lxy)
#define OOFEM_ERROR(...)
Definition error.h:79
#define Gxy
Definition matconst.h:72
@ CS_Thickness
Thickness.
@ SPRNodalRecoveryModelInterfaceType
@ ZZNodalRecoveryModelInterfaceType
@ LayeredCrossSectionInterfaceType
@ SpatialLocalizerInterfaceType
static FloatArray Vec3(const double &a, const double &b, const double &c)
Definition floatarray.h:607

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