OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
nlstructuralelement.h
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 - 2013 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 #ifndef nlstructuralelement_h
36 #define nlstructuralelement_h
37 
38 #include "../sm/Elements/structuralelement.h"
39 
41 
42 #define _IFT_NLStructuralElement_nlgeoflag "nlgeo"
43 
44 
45 namespace oofem {
46 class TimeStep;
47 class Node;
48 class Material;
49 class GaussPoint;
50 class FloatMatrix;
51 class FloatArray;
52 class IntArray;
53 
87 {
88 protected:
91 
92 public:
98  NLStructuralElement(int n, Domain * d);
100  virtual ~NLStructuralElement() { }
101 
108  int giveGeometryMode() { return nlGeometry; }
109 
120  void computeFirstPKStressVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep);
121 
131  void computeCauchyStressVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep);
132 
144  virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep);
145 
146 
155  virtual void computeInitialStressMatrix(FloatMatrix &answer, TimeStep *tStep);
156 
172 
182  virtual void giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord = 0);
183 
196  void giveInternalForcesVector_withIRulesAsSubcells(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord = 0);
197 
206  virtual void computeDeformationGradientVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep);
207 
211  double computeCurrentVolume(TimeStep *tStep);
212 
213  // data management
215  virtual void giveInputRecord(DynamicInputRecord &input);
216 
217  // definition
218  virtual const char *giveClassName() const { return "NLStructuralElement"; }
219 
220 protected:
221  virtual int checkConsistency();
231  virtual void computeBHmatrixAt(GaussPoint *gp, FloatMatrix &answer) {
232  OOFEM_ERROR("method not implemented for this element");
233  return;
234  }
235  friend class GradDpElement;
236  friend class PhaseFieldElement;
238 };
239 } // end namespace oofem
240 #endif // nlstructuralelement_h
virtual const char * giveClassName() const
int nlGeometry
Flag indicating if geometrical nonlinearities apply.
Class and object Domain.
Definition: domain.h:115
Abstract base class for "structural" finite elements with geometrical nonlinearities.
virtual ~NLStructuralElement()
Destructor.
Abstract class for phase field formulation.
void giveInternalForcesVector_withIRulesAsSubcells(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord=0)
Evaluates nodal representation of real internal forces.
virtual int checkConsistency()
Performs consistency check.
MatResponseMode
Describes the character of characteristic material matrix.
Abstract class for gradient formulation of coupled damage-plasticity model(GradDp).
Definition: graddpelement.h:48
void computeFirstPKStressVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep)
Computes the first Piola-Kirchhoff stress tensor on Voigt format.
Abstract base class for all "structural" finite elements.
double computeCurrentVolume(TimeStep *tStep)
Computes the current volume of element.
#define OOFEM_ERROR(...)
Definition: error.h:61
void computeCauchyStressVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep)
Computes the Cauchy stress tensor on Voigt format.
int giveGeometryMode()
Returns the geometry mode describing the formulation used in the internal work 0 - Engineering (small...
Provides Xfem interface for a structural element.
virtual void giveInputRecord(DynamicInputRecord &input)
Setups the input record string of receiver.
Class representing vector of real numbers.
Definition: floatarray.h:82
virtual void computeDeformationGradientVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep)
Computes the deformation gradient in Voigt form at integration point ip and at time step tStep...
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 void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
Computes the stiffness matrix of receiver.
Class representing the general Input Record.
Definition: inputrecord.h:101
NLStructuralElement(int n, Domain *d)
Constructor.
Class representing the a dynamic Input Record.
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
virtual void computeBHmatrixAt(GaussPoint *gp, FloatMatrix &answer)
Computes a matrix which, multiplied by the column matrix of nodal displacements, gives the displaceme...
virtual void giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord=0)
Evaluates nodal representation of real internal forces.
the oofem namespace is to define a context or scope in which all oofem names are defined.
virtual void computeInitialStressMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes the initial stiffness matrix of receiver.
Class representing integration point in finite element program.
Definition: gausspoint.h:93
Class representing solution step.
Definition: timestep.h:80
void computeStiffnessMatrix_withIRulesAsSubcells(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
Computes the stiffness matrix of receiver.

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:30 for OOFEM by doxygen 1.8.11 written by Dimitri van Heesch, © 1997-2011