OOFEM 3.0
Loading...
Searching...
No Matches
tria1platesubsoil.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
38#include "fei2dtrlin.h"
39#include "node.h"
40#include "material.h"
41#include "crosssection.h"
42#include "gausspoint.h"
44#include "floatmatrix.h"
45#include "floatarray.h"
46#include "intarray.h"
47#include "load.h"
48#include "mathfem.h"
49#include "classfactory.h"
50
51namespace oofem {
53
54FEI2dTrLin Tria1PlateSubSoil :: interp_lin(1, 2);
55
56Tria1PlateSubSoil :: Tria1PlateSubSoil(int n, Domain *aDomain) :
59{
62}
63
64
66Tria1PlateSubSoil :: giveInterpolation(DofIDItem id) const
67{
68 return & interp_lin;
69}
70
71
73Tria1PlateSubSoil :: giveInterpolation() const { return & interp_lin; }
74
75
76void
77Tria1PlateSubSoil :: computeGaussPoints()
78// Sets up the array containing the four Gauss points of the receiver.
79{
80 if ( integrationRulesArray.size() == 0 ) {
81 integrationRulesArray.resize( 1 );
82 integrationRulesArray [ 0 ] = std::make_unique<GaussIntegrationRule>(1, this, 1, 5);
83 this->giveCrossSection()->setupIntegrationPoints(* integrationRulesArray [ 0 ], numberOfGaussPoints, this);
84 }
85}
86
87
88void
89Tria1PlateSubSoil :: computeBodyLoadVectorAt(FloatArray &answer, Load *forLoad, TimeStep *tStep, ValueModeType mode)
90{
91 OOFEM_ERROR("Body load not supported, use surface load instead");
92}
93
94
95void
96Tria1PlateSubSoil :: computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int li, int ui)
97// Returns the [3x3] strain-displacement matrix {B} of the receiver,
98// evaluated at gp.
99{
100 FloatArray n;
101 FloatMatrix dn;
102
103 this->interp_lin.evaldNdx( dn, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(this) );
104 this->interp_lin.evalN( n, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(this) );
105
106 answer.resize(3, 3);
107 answer.zero();
108
110 for ( int i = 0; i < 3; ++i ) {
111 answer(0, i) = n(i); // eps_z
112 answer(1, i) = dn(i, 0); // gamma_xz
113 answer(2, i) = dn(i, 1); // gamma_yz
114 }
115}
116
117
118void
119Tria1PlateSubSoil :: computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep)
120{
121 answer = this->giveStructuralCrossSection()->giveGeneralizedStress_PlateSubSoil(strain, gp, tStep);
122}
123
124
125void
126Tria1PlateSubSoil :: computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
127{
128 answer = this->giveStructuralCrossSection()->give2dPlateSubSoilStiffMtrx(rMode, gp, tStep);
129}
130
131
132void
133Tria1PlateSubSoil :: giveDofManDofIDMask(int inode, IntArray &answer) const
134{
135 answer = {D_w};
136}
137
138
139void
140Tria1PlateSubSoil :: computeMidPlaneNormal(FloatArray &answer, const GaussPoint *gp)
141{
142 FloatArray u, v;
143 u.beDifferenceOf( this->giveNode(2)->giveCoordinates(), this->giveNode(1)->giveCoordinates() );
144 v.beDifferenceOf( this->giveNode(3)->giveCoordinates(), this->giveNode(1)->giveCoordinates() );
145
146 answer.beVectorProductOf(u, v);
147 answer.normalize();
148}
149
150
151double
152Tria1PlateSubSoil :: giveCharacteristicLength(const FloatArray &normalToCrackPlane)
153//
154// returns receiver's characteristic length for crack band models
155// for a crack formed in the plane with normal normalToCrackPlane.
156//
157{
158 return this->giveCharacteristicLengthForPlaneElements(normalToCrackPlane);
159}
160
161
162double
163Tria1PlateSubSoil :: computeVolumeAround(GaussPoint *gp)
164{
165 double detJ, weight;
166
167 weight = gp->giveWeight();
168 detJ = fabs( this->interp_lin.giveTransformationJacobian( gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(this) ) );
169 return detJ * weight;
170}
171
172
173void
174Tria1PlateSubSoil :: computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep)
175// Returns the lumped mass matrix of the receiver.
176{
177 OOFEM_ERROR("Mass matrix not provided");
178}
179
180
181int
182Tria1PlateSubSoil :: giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
183{
184 return StructuralElement :: giveIPValue(answer, gp, type, tStep);
185}
186
187Interface *
188Tria1PlateSubSoil :: giveInterface(InterfaceType interface)
189{
190 if ( interface == ZZNodalRecoveryModelInterfaceType ) {
191 return static_cast< ZZNodalRecoveryModelInterface * >(this);
192 } else if ( interface == SPRNodalRecoveryModelInterfaceType ) {
193 return static_cast< SPRNodalRecoveryModelInterface * >(this);
194 }
195
196 return NULL;
197}
198
199void
200Tria1PlateSubSoil :: SPRNodalRecoveryMI_giveSPRAssemblyPoints(IntArray &pap)
201{
202 pap.resize(3);
203 for ( int i = 1; i < 4; i++ ) {
204 pap.at(i) = this->giveNode(i)->giveNumber();
205 }
206}
207
208void
209Tria1PlateSubSoil :: SPRNodalRecoveryMI_giveDofMansDeterminedByPatch(IntArray &answer, int pap)
210{
211 int found = 0;
212 answer.resize(1);
213
214 for ( int i = 1; i < 4; i++ ) {
215 if ( pap == this->giveNode(i)->giveNumber() ) {
216 found = 1;
217 }
218 }
219
220 if ( found ) {
221 answer.at(1) = pap;
222 } else {
223 OOFEM_ERROR("node unknown");
224 }
225}
226
227void
228Tria1PlateSubSoil ::computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer)
229// Returns the [1x3] displacement interpolation matrix {N}
230{
231 FloatArray N(3);
232 giveInterpolation()->evalN(N, iLocCoord, FEIElementGeometryWrapper(this) );
233 answer.beNMatrixOf(N, 1);
234}
235
236
237void
238Tria1PlateSubSoil :: computeSurfaceNMatrix(FloatMatrix &answer, int boundaryID, const FloatArray &lcoords)
239{
240 if (boundaryID == 1) {
241 this->computeNmatrixAt(lcoords, answer);
242 } else {
243 OOFEM_ERROR("computeSurfaceNMatrix: Only one surface is supported with id=1");
244 }
245}
246
247// void
248// Tria1PlateSubSoil :: giveSurfaceDofMapping(IntArray &answer, int iSurf) const
249// {
250// answer.resize(3);
251// answer.zero();
252// if ( iSurf == 1 ) {
253// for (int i = 1; i<=3; i++) {
254// answer.at(i) = i;
255// }
256// } else {
257// OOFEM_ERROR("wrong surface number");
258// }
259// }
260//
261//
262// double
263// Tria1PlateSubSoil :: computeSurfaceVolumeAround(GaussPoint *gp, int iSurf)
264// {
265// return this->computeVolumeAround(gp);
266// }
267//
268//
269// int
270// Tria1PlateSubSoil :: computeLoadLSToLRotationMatrix(FloatMatrix &answer, int isurf, GaussPoint *gp)
271// {
272// return 0;
273// }
274
275
276} // 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
double giveCharacteristicLengthForPlaneElements(const FloatArray &normalToCrackPlane)
Definition element.C:1188
CrossSection * giveCrossSection()
Definition element.C:534
int giveNumber() const
Definition femcmpnn.h:104
void beDifferenceOf(const FloatArray &a, const FloatArray &b)
Definition floatarray.C:403
void beVectorProductOf(const FloatArray &v1, const FloatArray &v2)
Definition floatarray.C:476
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
double giveWeight()
Returns integration weight of receiver.
Definition gausspoint.h:180
void resize(int n)
Definition intarray.C:73
int & at(std::size_t i)
Definition intarray.h:104
StructuralCrossSection * giveStructuralCrossSection()
Helper function which returns the structural cross-section for the element.
StructuralElement(int n, Domain *d)
FEInterpolation * giveInterpolation() const override
void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer) override
static FEI2dTrLin interp_lin
ZZNodalRecoveryModelInterface(Element *element)
Constructor.
#define OOFEM_ERROR(...)
Definition error.h:79
@ SPRNodalRecoveryModelInterfaceType
@ ZZNodalRecoveryModelInterfaceType

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