OOFEM 3.0
Loading...
Searching...
No Matches
tria2platesubsoil.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 "crosssection.h"
38#include "gausspoint.h"
40#include "classfactory.h"
41
42namespace oofem {
44
45FEI2dTrQuad Tria2PlateSubSoil :: interp_quad(1, 2);
46
47Tria2PlateSubSoil :: Tria2PlateSubSoil(int n, Domain *aDomain) :
48 Tria1PlateSubSoil(n, aDomain)
49{
52}
53
54
56Tria2PlateSubSoil :: giveInterpolation(DofIDItem id) const
57{
58 return & interp_quad;
59}
60
61
63Tria2PlateSubSoil :: giveInterpolation() const { return & interp_quad; }
64
65
66void
67Tria2PlateSubSoil :: computeGaussPoints()
68// Sets up the array containing the four Gauss points of the receiver.
69{
70 if ( integrationRulesArray.size() == 0 ) {
71 integrationRulesArray.resize( 1 );
72 integrationRulesArray [ 0 ] = std::make_unique<GaussIntegrationRule>(1, this, 1, 5);
73 this->giveCrossSection()->setupIntegrationPoints(* integrationRulesArray [ 0 ], numberOfGaussPoints, this);
74 }
75}
76
77
78
79void
80Tria2PlateSubSoil :: computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int li, int ui)
81// Returns the [3x6] strain-displacement matrix {B} of the receiver,
82// evaluated at gp.
83{
84 FloatArray n;
85 FloatMatrix dn;
86
87 this->interp_quad.evaldNdx( dn, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(this) );
89
90 answer.resize(3, 6);
91 answer.zero();
92
94 for ( int i = 0; i < 6; ++i ) {
95 answer(0, i) = n(i); // eps_z
96 answer(1, i) = dn(i, 0); // gamma_xz
97 answer(2, i) = dn(i, 1); // gamma_yz
98 }
99}
100
101
102//TODO ZZNodalRecoveryModel can not determine some values. This is caused by sum of row entries is zero for (N^T)N matrix for vertices, yielding zero entries in the lumped form.
103
104void
105Tria2PlateSubSoil :: SPRNodalRecoveryMI_giveSPRAssemblyPoints(IntArray &pap)
106{
107 pap.resize(3);
108 pap.at(1) = this->giveNode(1)->giveNumber();
109 pap.at(2) = this->giveNode(2)->giveNumber();
110 pap.at(3) = this->giveNode(3)->giveNumber();
111}
112
113void
114Tria2PlateSubSoil :: SPRNodalRecoveryMI_giveDofMansDeterminedByPatch(IntArray &answer, int pap)
115{
116 answer.resize(3);
117 if ( pap == this->giveNode(1)->giveNumber() ) {
118 answer.at(1) = pap;
119 answer.at(2) = this->giveNode(4)->giveNumber();
120 answer.at(3) = this->giveNode(6)->giveNumber();
121 } else if ( pap == this->giveNode(2)->giveNumber() ) {
122 answer.at(1) = pap;
123 answer.at(2) = this->giveNode(5)->giveNumber();
124 answer.at(3) = this->giveNode(4)->giveNumber();
125 } else if ( pap == this->giveNode(3)->giveNumber() ) {
126 answer.at(1) = pap;
127 answer.at(2) = this->giveNode(6)->giveNumber();
128 answer.at(3) = this->giveNode(5)->giveNumber();
129 } else {
130 OOFEM_ERROR("node unknown");
131 }
132}
133
135Tria2PlateSubSoil :: SPRNodalRecoveryMI_givePatchType()
136{
138}
139
140void
141Tria2PlateSubSoil ::computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer)
142// Returns the [1x6] displacement interpolation matrix {N}
143{
144 FloatArray N(6);
145 giveInterpolation()->evalN(N, iLocCoord, FEIElementGeometryWrapper(this) );
146 answer.beNMatrixOf(N, 1);
147}
148
149void
150Tria2PlateSubSoil :: computeSurfaceNMatrix(FloatMatrix &answer, int boundaryID, const FloatArray &lcoords)
151{
152 if (boundaryID == 1) {
153 this->computeNmatrixAt(lcoords, answer);
154 } else {
155 OOFEM_ERROR("computeSurfaceNMatrix: Only one surface is supported with id=1");
156 }
157}
158
159} // end namespace oofem
#define N(a, b)
#define REGISTER_Element(class)
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
int giveNumber() const
Definition femcmpnn.h:104
void resize(Index rows, Index cols)
Definition floatmatrix.C:79
void beNMatrixOf(const FloatArray &n, int nsd)
void zero()
Zeroes all coefficient of receiver.
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
Tria1PlateSubSoil(int n, Domain *d)
void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer) override
FEInterpolation * giveInterpolation() const override
static FEI2dTrQuad interp_quad
#define OOFEM_ERROR(...)
Definition error.h:79

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