OOFEM 3.0
Loading...
Searching...
No Matches
quad1platesubsoil.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 "fei2dquadlin.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
54FEI2dQuadLin Quad1PlateSubSoil :: interp_lin(1, 2);
55
56Quad1PlateSubSoil :: Quad1PlateSubSoil(int n, Domain *aDomain) :
59{
62}
63
64
66Quad1PlateSubSoil :: giveInterpolation() const { return & interp_lin; }
67
68
70Quad1PlateSubSoil :: giveInterpolation(DofIDItem id) const
71{
72 return & interp_lin;
73}
74
75
76void
77Quad1PlateSubSoil :: 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
89Quad1PlateSubSoil :: 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
96Quad1PlateSubSoil :: computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int li, int ui)
97// Returns the [3x4] 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, 4);
107 answer.zero();
108
110 for ( int i = 0; i < 4; ++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
119Quad1PlateSubSoil :: computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep)
120{
121 answer = this->giveStructuralCrossSection()->giveGeneralizedStress_PlateSubSoil(strain, gp, tStep);
122}
123
124
125void
126Quad1PlateSubSoil :: computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
127{
128 answer = this->giveStructuralCrossSection()->give2dPlateSubSoilStiffMtrx(rMode, gp, tStep);
129}
130
131void
132Quad1PlateSubSoil :: giveDofManDofIDMask(int inode, IntArray &answer) const
133{
134 answer = {D_w};
135}
136
137
138void
139Quad1PlateSubSoil :: computeMidPlaneNormal(FloatArray &answer, const GaussPoint *gp)
140{
141 FloatArray u, v;
142 u.beDifferenceOf( this->giveNode(2)->giveCoordinates(), this->giveNode(1)->giveCoordinates() );
143 v.beDifferenceOf( this->giveNode(3)->giveCoordinates(), this->giveNode(1)->giveCoordinates() );
144
145 answer.beVectorProductOf(u, v);
146 answer.normalize();
147}
148
149
150double
151Quad1PlateSubSoil :: giveCharacteristicLength(const FloatArray &normalToCrackPlane)
152//
153// returns receiver's characteristic length for crack band models
154// for a crack formed in the plane with normal normalToCrackPlane.
155//
156{
157 return this->giveCharacteristicLengthForPlaneElements(normalToCrackPlane);
158}
159
160
161double
162Quad1PlateSubSoil :: computeVolumeAround(GaussPoint *gp)
163{
164 double weight = gp->giveWeight();
165 double detJ = fabs( this->interp_lin.giveTransformationJacobian( gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(this) ) );
166 return detJ * weight;
167}
168
169
170void
171Quad1PlateSubSoil :: computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep)
172// Returns the lumped mass matrix of the receiver.
173{
174 OOFEM_ERROR("Mass matrix not provided");
175}
176
177
178int
179Quad1PlateSubSoil :: giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
180{
181 if ( type == IST_ShellForceTensor ){
182 FloatArray help;
183 answer.resize(6);
184 help = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() )->giveStressVector();
185 answer.at(1) = 0.0; // nx
186 answer.at(2) = 0.0; // ny
187 answer.at(3) = help.at(1); // nz
188 answer.at(4) = help.at(3); // vyz
189 answer.at(5) = help.at(2); // vxz
190 answer.at(6) = 0.0; // vxy
191 return 1;
192 }
193 return StructuralElement :: giveIPValue(answer, gp, type, tStep);
194}
195
196Interface *
197Quad1PlateSubSoil :: giveInterface(InterfaceType interface)
198{
199 if ( interface == ZZNodalRecoveryModelInterfaceType ) {
200 return static_cast< ZZNodalRecoveryModelInterface * >(this);
201 } else if ( interface == SPRNodalRecoveryModelInterfaceType ) {
202 return static_cast< SPRNodalRecoveryModelInterface * >(this);
203 }
204
205 return NULL;
206}
207
208void
209Quad1PlateSubSoil :: SPRNodalRecoveryMI_giveSPRAssemblyPoints(IntArray &pap)
210{
211 pap.resize(4);
212 for ( int i = 1; i < 5; i++ ) {
213 pap.at(i) = this->giveNode(i)->giveNumber();
214 }
215}
216
217void
218Quad1PlateSubSoil :: SPRNodalRecoveryMI_giveDofMansDeterminedByPatch(IntArray &answer, int pap)
219{
220 int found = 0;
221 answer.resize(1);
222
223 for ( int i = 1; i < 5; i++ ) {
224 if ( pap == this->giveNode(i)->giveNumber() ) {
225 found = 1;
226 }
227 }
228
229 if ( found ) {
230 answer.at(1) = pap;
231 } else {
232 OOFEM_ERROR("node unknown");
233 }
234}
235
236void
237Quad1PlateSubSoil ::computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer)
238// Returns the [1x4] displacement interpolation matrix {N}
239{
240 FloatArray N(4);
241 giveInterpolation()->evalN(N, iLocCoord, FEIElementGeometryWrapper(this) );
242 answer.beNMatrixOf(N, 1);
243}
244
245
246void
247Quad1PlateSubSoil ::computeSurfaceNMatrix(FloatMatrix &answer, int boundaryID, const FloatArray &lcoords)
248{
249 if (boundaryID == 1) {
250 this->computeNmatrixAt(lcoords, answer);
251 } else {
252 OOFEM_ERROR("computeSurfaceNMatrix: Only one surface is supported with id=1");
253 }
254}
255
256
257
258// void
259// Quad1PlateSubSoil :: giveSurfaceDofMapping(IntArray &answer, int iSurf) const
260// {
261// answer.resize(4);
262// answer.zero();
263// if ( iSurf == 1 ) {
264// for (int i = 1; i <= 4; i++) {
265// answer.at(i) = i;
266// }
267// } else {
268// OOFEM_ERROR("wrong surface number");
269// }
270// }
271
272
273// double
274// Quad1PlateSubSoil :: computeSurfaceVolumeAround(GaussPoint *gp, int iSurf)
275// {
276// return this->computeVolumeAround(gp);
277// }
278
279
280// int
281// Quad1PlateSubSoil :: computeLoadLSToLRotationMatrix(FloatMatrix &answer, int isurf, GaussPoint *gp)
282// {
283// return 0;
284// }
285
286
287} // 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 resize(Index s)
Definition floatarray.C:94
double & at(Index i)
Definition floatarray.h:202
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
IntegrationPointStatus * giveMaterialStatus(IntegrationPointStatusIDType key=IPSID_Default)
Definition gausspoint.h:204
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
static FEI2dQuadLin interp_lin
FEInterpolation * giveInterpolation() const override
void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer) override
StructuralCrossSection * giveStructuralCrossSection()
Helper function which returns the structural cross-section for the element.
StructuralElement(int n, Domain *d)
const FloatArray & giveStressVector() const
Returns the const pointer to receiver's stress vector.
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