OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
linedistributedspring.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 - 2014 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 
35 #include "../sm/Elements/linedistributedspring.h"
36 #include "../sm/Materials/structuralms.h"
37 #include "../sm/CrossSections/structuralcrosssection.h"
38 #include "fei3dlinelin.h"
39 #include "node.h"
40 #include "material.h"
41 #include "crosssection.h"
42 #include "gausspoint.h"
43 #include "gaussintegrationrule.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 
51 namespace oofem {
52 REGISTER_Element(LineDistributedSpring);
53 
55 
59 {
61  numberOfDofMans = 2;
62 }
63 
64 
67 {
68  return & interp_lin;
69 }
70 
71 
74 
75 
76 void
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 ].reset( new GaussIntegrationRule(1, this, 1, 5) );
84  }
85 }
86 
87 
88 void
90 {
91  OOFEM_ERROR("Body load not supported, use surface load instead");
92 }
93 
94 
95 void
97 // Returns the [3x3] strain-displacement matrix {B} of the receiver,
98 // evaluated at gp.
99 {
100  FloatArray n;
101  int ndofs = this->dofs.giveSize();
102 
104 
105  answer.resize(ndofs, ndofs*2);
106  answer.zero();
107 
108  for (int idof=1; idof<=ndofs; idof++) {
109  answer.at(idof, idof) = n.at(1);
110  answer.at(idof, ndofs+idof) = n.at(2);
111  }
112 }
113 
114 
115 void
117 {
118  int ndofs = this->dofs.giveSize();
119  answer.resize(ndofs);
120 
121  for (int idof=1; idof<=ndofs; idof++) {
122  answer.at(idof) = strain.at(idof)*this->springStiffnesses.at(idof);
123  }
124  //this->giveStructuralCrossSection()->giveGeneralizedStress_DistributedSpring(answer, gp, strain, tStep);
125 }
126 
127 
128 void
130 {
131  //this->giveStructuralCrossSection()->give2dPlateSubSoilStiffMtrx(answer, rMode, gp, tStep);
132  int ndofs = this->dofs.giveSize();
133  answer.resize(ndofs, ndofs);
134  answer.beDiagonal(this->springStiffnesses);
135 }
136 
137 
138 void
140  TimeStep *tStep, int useUpdatedGpRecord)
141 {
142  FloatArray u;
143  FloatMatrix k;
144 
145  this->computeVectorOf(VM_Total, tStep, u);
146  this->computeStiffnessMatrix(k, TangentStiffness, tStep);
147  answer.beProductOf(k, u);
148 
149 }
150 
151 
152 
155 {
156  IRResultType result; // Required by IR_GIVE_FIELD macro
157 
160 
162  OOFEM_ERROR ("dofs and k params size mismatch");
163  }
164  // from element
166 }
167 
168 int
170 {
171  // skip StructuralElement consistency as there is checjk for structurak cross section capability
172  return Element::checkConsistency();
173 }
174 
176 {
177  FloatArray stress, strain;
178 
179  fprintf(File, "element %d (%8d) :\n", this->giveLabel(), this->giveNumber() );
180 
181  for ( GaussPoint *gp: *this->giveDefaultIntegrationRulePtr() ) {
182  fprintf(File, "\tGP 1.%d :\n", gp->giveNumber());
183  this->computeStrainVector(strain, gp, tStep);
184  this->computeStressVector(stress, strain, gp, tStep);
185  fprintf(File, "\t\tStrain");
186  for (int i=1; i<=strain.giveSize(); i++) fprintf (File, " %e", strain.at(i));
187  fprintf(File, "\n\t\tStress");
188  for (int i=1; i<=stress.giveSize(); i++) fprintf (File, " %e", stress.at(i));
189  fprintf(File, "\n");
190  }
191 
192 }
193 
194 
195 
196 void
198 {
199  answer = this->dofs;
200 }
201 
202 
203 double
205 {
206  double detJ, weight;
207 
208  weight = gp->giveWeight();
210  return detJ * weight;
211 }
212 
213 
214 void
216 // Returns the lumped mass matrix of the receiver.
217 {
218  OOFEM_ERROR("Mass matrix not provided");
219 }
220 
221 
222 int
224 {
225  return StructuralElement :: giveIPValue(answer, gp, type, tStep);
226 }
227 
228 Interface *
230 {
231  if ( interface == ZZNodalRecoveryModelInterfaceType ) {
232  return static_cast< ZZNodalRecoveryModelInterface * >(this);
233  } else if ( interface == SPRNodalRecoveryModelInterfaceType ) {
234  return static_cast< SPRNodalRecoveryModelInterface * >(this);
235  }
236 
237  return NULL;
238 }
239 
240 void
242 {
243  pap.resize(2);
244  for ( int i = 1; i < 3; i++ ) {
245  pap.at(i) = this->giveNode(i)->giveNumber();
246  }
247 }
248 
249 void
251 {
252  int found = 0;
253  answer.resize(1);
254 
255  for ( int i = 1; i < 3; i++ ) {
256  if ( pap == this->giveNode(i)->giveNumber() ) {
257  found = 1;
258  }
259  }
260 
261  if ( found ) {
262  answer.at(1) = pap;
263  } else {
264  OOFEM_ERROR("node unknown");
265  }
266 }
267 
268 
269 } // end namespace oofem
CrossSection * giveCrossSection()
Definition: element.C:495
InternalStateType
Type representing the physical meaning of element or constitutive model internal variable.
The element interface required by ZZNodalRecoveryModel.
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
Class and object Domain.
Definition: domain.h:115
void computeVectorOf(ValueModeType u, TimeStep *tStep, FloatArray &answer)
Returns local vector of unknowns.
Definition: element.C:86
virtual IntegrationRule * giveDefaultIntegrationRulePtr()
Access method for default integration rule.
Definition: element.h:822
int checkConsistency()
Performs consistency check.
The element interface required by ZZNodalRecoveryModel.
virtual void giveDofManDofIDMask(int inode, IntArray &) const
Returns dofmanager dof mask for node.
double & at(int i)
Coefficient access function.
Definition: floatarray.h:131
ValueModeType
Type representing the mode of UnknownType or CharType, or similar types.
Definition: valuemodetype.h:78
virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
Computes constitutive matrix of receiver.
#define _IFT_LineDistributedSpring_Stifnesses
void beDiagonal(const FloatArray &diag)
Modifies receiver to be a diagonal matrix with the components specified in diag.
Definition: floatmatrix.C:1433
virtual double giveTransformationJacobian(const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the determinant of the transformation.
Definition: fei3dlinelin.C:112
Class implementing an array of integers.
Definition: intarray.h:61
int & at(int i)
Coefficient access function.
Definition: intarray.h:103
MatResponseMode
Describes the character of characteristic material matrix.
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
virtual void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep)
Computes the stress vector of receiver at given integration point, at time step tStep.
Class representing a general abstraction for finite element interpolation class.
Definition: feinterpol.h:132
Abstract base class for all "structural" finite elements.
virtual void SPRNodalRecoveryMI_giveSPRAssemblyPoints(IntArray &pap)
#define OOFEM_ERROR(...)
Definition: error.h:61
REGISTER_Element(LSpace)
DofIDItem
Type representing particular dof type.
Definition: dofiditem.h:86
virtual void evalN(FloatArray &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the array of interpolation functions (shape functions) at given point.
Definition: fei3dlinelin.C:49
virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int=1, int=ALL_STRAINS)
Computes the geometrical matrix of receiver in given integration point.
virtual double giveWeight()
Returns integration weight of receiver.
Definition: gausspoint.h:181
virtual Interface * giveInterface(InterfaceType it)
Interface requesting service.
Wrapper around element definition to provide FEICellGeometry interface.
Definition: feinterpol.h:95
virtual int checkConsistency()
Performs consistency check.
Definition: element.h:763
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Receiver becomes the result of the product of aMatrix and anArray.
Definition: floatarray.C:676
void giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord)
Returns equivalent nodal forces vectors.
virtual void computeStrainVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep)
Compute strain vector of receiver evaluated at given integration point at given time step from elemen...
double at(int i, int j) const
Coefficient access function.
Definition: floatmatrix.h:176
virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
Computes numerically stiffness matrix of receiver.
void resize(int n)
Checks size of receiver towards requested bounds.
Definition: intarray.C:124
virtual int setupIntegrationPoints(IntegrationRule &irule, int npoints, Element *element)
Sets up integration rule for the given element.
Definition: crosssection.C:54
int numberOfGaussPoints
Number of integration points as specified by nip.
Definition: element.h:188
Class representing vector of real numbers.
Definition: floatarray.h:82
Implementation of matrix containing floating point numbers.
Definition: floatmatrix.h:94
IRResultType
Type defining the return values of InputRecord reading operations.
Definition: irresulttype.h:47
virtual FEInterpolation * giveInterpolation() const
void resize(int rows, int cols)
Checks size of receiver towards requested bounds.
Definition: floatmatrix.C:1358
Class representing the general Input Record.
Definition: inputrecord.h:101
Class Interface.
Definition: interface.h:82
std::vector< std::unique_ptr< IntegrationRule > > integrationRulesArray
List of integration rules of receiver (each integration rule contains associated integration points a...
Definition: element.h:170
void zero()
Zeroes all coefficient of receiver.
Definition: floatmatrix.C:1326
virtual void computeBodyLoadVectorAt(FloatArray &answer, Load *load, TimeStep *tStep, ValueModeType mode)
Computes the load vector due to body load acting on receiver, at given time step. ...
void printOutputAt(FILE *File, TimeStep *tStep)
Prints output of receiver to stream, for given time step.
InterfaceType
Enumerative type, used to identify interface type.
Definition: interfacetype.h:43
#define _IFT_LineDistributedSpring_Dofs
virtual void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes lumped mass matrix of receiver.
Load is base abstract class for all loads.
Definition: load.h:61
int giveSize() const
Definition: intarray.h:203
int giveSize() const
Returns the size of receiver.
Definition: floatarray.h:218
int giveLabel() const
Definition: element.h:1055
the oofem namespace is to define a context or scope in which all oofem names are defined.
#define IR_GIVE_FIELD(__ir, __value, __id)
Macro facilitating the use of input record reading methods.
Definition: inputrecord.h:69
int giveNumber() const
Definition: femcmpnn.h:107
virtual void computeGaussPoints()
Initializes the array of integration rules member variable.
Node * giveNode(int i) const
Returns reference to the i-th node of element.
Definition: element.h:610
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Class representing integration point in finite element program.
Definition: gausspoint.h:93
virtual double computeVolumeAround(GaussPoint *gp)
Returns volume related to given integration point.
Class representing solution step.
Definition: timestep.h:80
virtual void SPRNodalRecoveryMI_giveDofMansDeterminedByPatch(IntArray &answer, int pap)
int numberOfDofMans
Number of dofmanagers.
Definition: element.h:149
const FloatArray & giveNaturalCoordinates()
Returns coordinate array of receiver.
Definition: gausspoint.h:138
Class representing Gaussian-quadrature integration rule.
void resize(int s)
Resizes receiver towards requested size.
Definition: floatarray.C:631

This page is part of the OOFEM documentation. Copyright (c) 2011 Borek Patzak
Project e-mail: info@oofem.org
Generated at Tue Jan 2 2018 20:07:29 for OOFEM by doxygen 1.8.11 written by Dimitri van Heesch, © 1997-2011