OOFEM 3.0
Loading...
Searching...
No Matches
quad2platesubsoil.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 "fei2dquadquad.h"
37#include "crosssection.h"
38#include "gausspoint.h"
40#include "classfactory.h"
41
42namespace oofem {
44
45FEI2dQuadQuad Quad2PlateSubSoil :: interp_quad(1, 2);
46
47Quad2PlateSubSoil :: Quad2PlateSubSoil(int n, Domain *aDomain) :
48 Quad1PlateSubSoil(n, aDomain)
49{
52}
53
54
56Quad2PlateSubSoil :: giveInterpolation() const { return & interp_quad; }
57
58
60Quad2PlateSubSoil :: giveInterpolation(DofIDItem id) const
61{
62 return & interp_quad;
63}
64
65
66void
67Quad2PlateSubSoil :: 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
80Quad2PlateSubSoil :: computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int li, int ui)
81// Returns the [3x8] 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, 8);
91 answer.zero();
92
94 for ( int i = 0; i < 8; ++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
101void
102Quad2PlateSubSoil :: SPRNodalRecoveryMI_giveSPRAssemblyPoints(IntArray &pap)
103{
104 pap.resize(8);
105 for ( int i = 1; i < 8; i++ ) {
106 pap.at(i) = this->giveNode(i)->giveNumber();
107 }
108}
109
110void
111Quad2PlateSubSoil :: SPRNodalRecoveryMI_giveDofMansDeterminedByPatch(IntArray &answer, int pap)
112{
113 int found = 0;
114 answer.resize(1);
115
116 for ( int i = 1; i <= 8; i++ ) {
117 if ( pap == this->giveNode(i)->giveNumber() ) {
118 found = 1;
119 }
120 }
121
122 if ( found ) {
123 answer.at(1) = pap;
124 } else {
125 OOFEM_ERROR("node %d not found on element %d", pap, this->giveNumber());
126 }
127}
128
129void
130Quad2PlateSubSoil ::computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer)
131// Returns the [1x8] displacement interpolation matrix {N}
132{
133 FloatArray N(8);
134 giveInterpolation()->evalN(N, iLocCoord, FEIElementGeometryWrapper(this) );
135 answer.beNMatrixOf(N, 1);
136}
137
138void
139Quad2PlateSubSoil :: computeSurfaceNMatrix(FloatMatrix &answer, int boundaryID, const FloatArray &lcoords)
140{
141 if (boundaryID == 1) {
142 this->computeNmatrixAt(lcoords, answer);
143 } else {
144 OOFEM_ERROR("computeSurfaceNMatrix: Only one surface is supported with id=1");
145 }
146}
147
148
149
150} // 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
Quad1PlateSubSoil(int n, Domain *d)
static FEI2dQuadQuad interp_quad
FEInterpolation * giveInterpolation() const override
void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer) override
#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