OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
htselement.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 htselement_h
36 #define htselement_h
37 
38 #include "../sm/Elements/structuralelement.h"
39 #include "gaussintegrationrule.h"
40 
41 #define _IFT_HTSelement_Name "htselement"
42 
43 namespace oofem {
50 {
51 protected:
53  //debug
54  double lambda, mu;
55  double cgX, cgY;
58 
59 public:
60  HTSelement(int n, Domain * d);
61  virtual ~HTSelement() { }
62 
64 
65  virtual const char *giveInputRecordName() const { return _IFT_HTSelement_Name; }
66  virtual const char *giveClassName() const { return "HTSelement"; }
67 
68 protected:
69  virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int, int) {; }
70  virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer) {; }
71  virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep);
72  double computeVolumeAroundSide(GaussPoint *gp, int elemSideNumber);
73  Node *giveSideNode(int elementSideNumber, int nodeNumber);
74  double giveSideLength(int sideNumber);
75  virtual int computeNumberOfDofs() { return 4 * numberOfEdges; }
76  virtual void computeGaussPoints();
77  virtual void giveDofManDofIDMask(int inode, IntArray &) const;
78  virtual StructuralElement *giveStructuralElement() { return this; }
79  //jak se pocita deformace???
80  virtual void computeStrainVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep) { answer.resize(numberOfStressDofs); }
81  //dodelat vypocet napeti!!!
82  virtual void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep) { answer.resize(numberOfStressDofs); }
83  virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) {
84  OOFEM_ERROR("Function not defined for this element and should never be called. This is a bug.");
85  }
86  //dodelat internal forces, budou potreba pro nelinearni vypocet
87  virtual void giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord) { answer.resize(numberOfDofs); }
88 
89  void computePuVectorAt(FloatArray &answer, FloatMatrix N, FloatArray u, GaussPoint *gp, int sideNumber);
90  void computePsVectorAt(FloatArray &answer, FloatArray t, GaussPoint *gp);
92  virtual int testElementExtension(ElementExtension ext) { return ( ext == Element_EdgeLoadSupport ); }
93 
94 
95  void computeFMatrixAt(FloatMatrix &answer, FloatMatrix N, GaussPoint *gp, int sideNumber);
96  void computeAMatrixAt(FloatMatrix &answer, FloatMatrix N, GaussPoint *gp, int sideNumber);
97  void computeUvMatrixAt(FloatMatrix &answer, GaussPoint *gp, int sideNubmer);
98  void computeSvMatrixAt(FloatMatrix &answer, GaussPoint *gp, int sideNumber);
99  void computeUgammaMatrixAt(FloatMatrix &answer, GaussPoint *gp);
100  void computeOutwardNormalMatrix(FloatMatrix &answer, int sideNumber);
101 
102 
103  void computeCenterOfGravity();
104  virtual int giveNumberOfNodes() const { return numberOfEdges; }
105  //uv functions
106  void uv1(FloatArray &answer, double x, double y);
107  void uv2(FloatArray &answer, double x, double y);
108  void uv3(FloatArray &answer, double x, double y);
109  void uv4(FloatArray &answer, double x, double y);
110  void uv5(FloatArray &answer, double x, double y);
111  void uv6(FloatArray &answer, double x, double y);
112  void uv7(FloatArray &answer, double x, double y);
113  void uv8(FloatArray &answer, double x, double y);
114  void uv9(FloatArray &answer, double x, double y);
115  void uv10(FloatArray &answer, double x, double y);
116  void uv11(FloatArray &answer, double x, double y);
117  void uv12(FloatArray &answer, double x, double y);
118  void uv25_4(FloatArray &answer, double x, double y);
119 
120  //sv functions
121  void sv1(FloatArray &answer, double x, double y);
122  void sv2(FloatArray &answer, double x, double y);
123  void sv3(FloatArray &answer, double x, double y);
124  void sv4(FloatArray &answer, double x, double y);
125  void sv5(FloatArray &answer, double x, double y);
126  void sv6(FloatArray &answer, double x, double y);
127  void sv7(FloatArray &answer, double x, double y);
128  void sv8(FloatArray &answer, double x, double y);
129  void sv9(FloatArray &answer, double x, double y);
130  void sv10(FloatArray &answer, double x, double y);
131  void sv11(FloatArray &answer, double x, double y);
132  void sv12(FloatArray &answer, double x, double y);
133  void sv25_4(FloatArray &answer, double x, double y);
134 
135  //u_gamma functions
136  double u_gammaConst(GaussPoint *gp);
137  double u_gammaLin(GaussPoint *gp);
138 };
139 } // end namespace oofem
140 #endif // htselement_h
virtual int giveNumberOfNodes() const
Returns number of nodes of receiver.
Definition: htselement.h:104
double computeVolumeAroundSide(GaussPoint *gp, int elemSideNumber)
Definition: htselement.C:148
void uv5(FloatArray &answer, double x, double y)
Definition: htselement.C:576
void uv12(FloatArray &answer, double x, double y)
Definition: htselement.C:657
void uv3(FloatArray &answer, double x, double y)
Definition: htselement.C:530
#define _IFT_HTSelement_Name
Definition: htselement.h:41
Class and object Domain.
Definition: domain.h:115
void uv4(FloatArray &answer, double x, double y)
Definition: htselement.C:537
void sv1(FloatArray &answer, double x, double y)
Definition: htselement.C:544
void computePrescribedDisplacementLoadVectorAt(FloatArray &answer, TimeStep *tStep, ValueModeType mode)
Definition: htselement.C:235
void uv8(FloatArray &answer, double x, double y)
Definition: htselement.C:597
void computeFMatrixAt(FloatMatrix &answer, FloatMatrix N, GaussPoint *gp, int sideNumber)
Definition: htselement.C:301
void uv6(FloatArray &answer, double x, double y)
Definition: htselement.C:583
ValueModeType
Type representing the mode of UnknownType or CharType, or similar types.
Definition: valuemodetype.h:78
void uv7(FloatArray &answer, double x, double y)
Definition: htselement.C:590
void sv11(FloatArray &answer, double x, double y)
Definition: htselement.C:680
void sv7(FloatArray &answer, double x, double y)
Definition: htselement.C:620
void uv10(FloatArray &answer, double x, double y)
Definition: htselement.C:643
Implements a Hybrid-Trefftz element See http://en.wikipedia.org/wiki/Trefftz_method for description...
Definition: htselement.h:49
void computeCenterOfGravity()
Definition: htselement.C:104
void sv12(FloatArray &answer, double x, double y)
Definition: htselement.C:688
Class implementing an array of integers.
Definition: intarray.h:61
MatResponseMode
Describes the character of characteristic material matrix.
void computeSvMatrixAt(FloatMatrix &answer, GaussPoint *gp, int sideNumber)
Definition: htselement.C:402
void uv11(FloatArray &answer, double x, double y)
Definition: htselement.C:650
virtual StructuralElement * giveStructuralElement()
Definition: htselement.h:78
void uv9(FloatArray &answer, double x, double y)
Definition: htselement.C:636
virtual const char * giveInputRecordName() const
Definition: htselement.h:65
virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int, int)
Computes the geometrical matrix of receiver in given integration point.
Definition: htselement.h:69
Abstract base class for all "structural" finite elements.
void computeAMatrixAt(FloatMatrix &answer, FloatMatrix N, GaussPoint *gp, int sideNumber)
Definition: htselement.C:317
virtual void giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord)
Returns equivalent nodal forces vectors.
Definition: htselement.h:87
void sv5(FloatArray &answer, double x, double y)
Definition: htselement.C:604
void sv2(FloatArray &answer, double x, double y)
Definition: htselement.C:552
#define OOFEM_ERROR(...)
Definition: error.h:61
ElementExtension
Type representing element extension.
void sv8(FloatArray &answer, double x, double y)
Definition: htselement.C:628
void sv3(FloatArray &answer, double x, double y)
Definition: htselement.C:560
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Definition: htselement.C:75
Node * giveSideNode(int elementSideNumber, int nodeNumber)
Definition: htselement.C:129
void sv10(FloatArray &answer, double x, double y)
Definition: htselement.C:672
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.
Definition: htselement.h:82
void computeOutwardNormalMatrix(FloatMatrix &answer, int sideNumber)
Definition: htselement.C:282
#define N(p, q)
Definition: mdm.C:367
void uv25_4(FloatArray &answer, double x, double y)
Definition: htselement.C:880
void sv9(FloatArray &answer, double x, double y)
Definition: htselement.C:664
void sv25_4(FloatArray &answer, double x, double y)
Definition: htselement.C:887
void computePsVectorAt(FloatArray &answer, FloatArray t, GaussPoint *gp)
Definition: htselement.C:224
double u_gammaLin(GaussPoint *gp)
Definition: htselement.C:903
void sv4(FloatArray &answer, double x, double y)
Definition: htselement.C:568
virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
Computes numerically stiffness matrix of receiver.
Definition: htselement.C:170
virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer)
Computes interpolation matrix for element unknowns.
Definition: htselement.h:70
Class representing vector of real numbers.
Definition: floatarray.h:82
void computePuVectorAt(FloatArray &answer, FloatMatrix N, FloatArray u, GaussPoint *gp, int sideNumber)
Definition: htselement.C:211
Implementation of matrix containing floating point numbers.
Definition: floatmatrix.h:94
void computeUgammaMatrixAt(FloatMatrix &answer, GaussPoint *gp)
Definition: htselement.C:491
IRResultType
Type defining the return values of InputRecord reading operations.
Definition: irresulttype.h:47
virtual int computeNumberOfDofs()
Computes or simply returns total number of element's local DOFs.
Definition: htselement.h:75
virtual void giveDofManDofIDMask(int inode, IntArray &) const
Returns dofmanager dof mask for node.
Definition: htselement.C:63
HTSelement(int n, Domain *d)
Definition: htselement.C:49
Class representing the general Input Record.
Definition: inputrecord.h:101
virtual void computeGaussPoints()
Initializes the array of integration rules member variable.
Definition: htselement.C:91
void uv2(FloatArray &answer, double x, double y)
Definition: htselement.C:523
virtual const char * giveClassName() const
Definition: htselement.h:66
void computeUvMatrixAt(FloatMatrix &answer, GaussPoint *gp, int sideNubmer)
Definition: htselement.C:331
double u_gammaConst(GaussPoint *gp)
Definition: htselement.C:897
virtual int testElementExtension(ElementExtension ext)
Tests if the element implements required extension.
Definition: htselement.h:92
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...
Definition: htselement.h:80
void uv1(FloatArray &answer, double x, double y)
Definition: htselement.C:516
virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
Computes constitutive matrix of receiver.
Definition: htselement.h:83
the oofem namespace is to define a context or scope in which all oofem names are defined.
void sv6(FloatArray &answer, double x, double y)
Definition: htselement.C:612
Class implementing node in finite element mesh.
Definition: node.h:87
double giveSideLength(int sideNumber)
Definition: htselement.C:156
Class representing integration point in finite element program.
Definition: gausspoint.h:93
Class representing solution step.
Definition: timestep.h:80
virtual ~HTSelement()
Definition: htselement.h:61
Element extension for edge loads.
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