OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
structuralelementevaluator.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 structuralelementevaluator_h
36 #define structuralelementevaluator_h
37 
38 #include "iga/iga.h"
39 #include "matresponsemode.h"
40 
41 namespace oofem {
57 {
58 protected:
62 
65  virtual void giveCharacteristicMatrix(FloatMatrix &answer, CharType mtrx, TimeStep *tStep);
66  virtual void giveCharacteristicVector(FloatArray &answer, CharType type, ValueModeType mode, TimeStep *tStep);
67 
72  virtual IntegrationRule *giveMassMtrxIntegrationRule() { return NULL; }
81  virtual void giveMassMtrxIntegrationMask(IntArray &answer) { answer.clear(); }
92  virtual void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep);
105  virtual void computeConsistentMassMatrix(FloatMatrix &answer, TimeStep *tStep, double &mass);
106 
107 protected:
108  virtual Element *giveElement() = 0;
109 
113  virtual void computeNMatrixAt(FloatMatrix &answer, GaussPoint *gp) = 0;
114  virtual void computeBMatrixAt(FloatMatrix &answer, GaussPoint *gp) = 0;
115  virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep);
116  virtual double computeVolumeAround(GaussPoint *gp) { return 0.; }
117  virtual void giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, bool useUpdatedGpRecord = false);
119  this->giveElement()->computeVectorOf(u, tStep, answer);
120  }
121  void computeVectorOf(PrimaryField &field, const IntArray &dofIdMask, ValueModeType u, TimeStep *tStep, FloatArray &answer) {
122  this->giveElement()->computeVectorOf(field, dofIdMask, u, tStep, answer);
123  }
124  bool isActivated(TimeStep *tStep) { return true; }
125  void updateInternalState(TimeStep *tStep);
126 
134  virtual void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep) = 0;
142  virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) = 0;
150  void computeStrainVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep, FloatArray &u);
151 
152  /*
153  * Assembles the code numbers of given integration element (sub-patch)
154  * This is done by obtaining list of nonzero shape functions and
155  * by collecting the code numbers of nodes corresponding to these
156  * shape functions
157  * @returns returns nonzero if integration rule code numbers differ from element code numbers
158  */
159  //virtual int giveIntegrationElementCodeNumbers(IntArray &answer, Element *elem,
160  // IntegrationRule *ie,);
161 
169  virtual int giveIntegrationElementLocalCodeNumbers(IntArray &answer, Element *elem,
170  IntegrationRule *ie);
171 #ifdef __OOFEG
173 #endif
174 };
175 } // end namespace oofem
176 #endif //structuralelementevaluator_h
This class represent a new concept on how to define elements.
virtual IntegrationRule * giveMassMtrxIntegrationRule()
Returns the integration rule for mass matrices, if relevant.
void computeVectorOf(ValueModeType u, TimeStep *tStep, FloatArray &answer)
Returns local vector of unknowns.
Definition: element.C:86
virtual void computeNMatrixAt(FloatMatrix &answer, GaussPoint *gp)=0
Computes the matrix for which the unknown field is obtained, typically [N1, 0, N2, 0, ...; 0, N1, 0, N2, ...].
virtual void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes lumped mass matrix of receiver.
Abstract class representing field of primary variables (those, which are unknown and are typically as...
Definition: primaryfield.h:104
virtual void giveCharacteristicVector(FloatArray &answer, CharType type, ValueModeType mode, TimeStep *tStep)
ValueModeType
Type representing the mode of UnknownType or CharType, or similar types.
Definition: valuemodetype.h:78
int rotationMatrixDefined
Flag indicating if transformation matrix has been already computed.
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
Abstract base class for all finite elements.
Definition: element.h:145
virtual int giveIntegrationElementLocalCodeNumbers(IntArray &answer, Element *elem, IntegrationRule *ie)
Assembles the local element code numbers of given integration element (sub-patch) This is done by obt...
void computeStrainVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep, FloatArray &u)
Optimized version, allowing to pass element displacements as parameter.
virtual void computeConsistentMassMatrix(FloatMatrix &answer, TimeStep *tStep, double &mass)
Computes consistent mass matrix of receiver using numerical integration over element volume...
virtual double computeVolumeAround(GaussPoint *gp)
Class implementing an array of integers.
Definition: intarray.h:61
MatResponseMode
Describes the character of characteristic material matrix.
virtual Element * giveElement()=0
virtual void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep)=0
Computes the stress vector.
Abstract base class representing integration rule.
virtual void giveCharacteristicMatrix(FloatMatrix &answer, CharType mtrx, TimeStep *tStep)
virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)=0
Computes constitutive matrix of receiver.
virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
void clear()
Clears the array (zero size).
Definition: intarray.h:177
UnknownType
Type representing particular unknown (its physical meaning).
Definition: unknowntype.h:55
friend void drawIGAPatchDeformedGeometry(Element *elem, StructuralElementEvaluator *se, oofegGraphicContext &gc, TimeStep *tStep, UnknownType)
Definition: iga.C:1006
virtual void giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, bool useUpdatedGpRecord=false)
Class representing vector of real numbers.
Definition: floatarray.h:82
Implementation of matrix containing floating point numbers.
Definition: floatmatrix.h:94
void computeVectorOf(PrimaryField &field, const IntArray &dofIdMask, ValueModeType u, TimeStep *tStep, FloatArray &answer)
void computeVectorOf(ValueModeType u, TimeStep *tStep, FloatArray &answer)
CharType
Definition: chartype.h:87
virtual void giveMassMtrxIntegrationMask(IntArray &answer)
Returns mask indicating, which unknowns (their type and ordering is the same as element unknown vecto...
the oofem namespace is to define a context or scope in which all oofem names are defined.
Class representing integration point in finite element program.
Definition: gausspoint.h:93
Class representing solution step.
Definition: timestep.h:80
virtual void computeBMatrixAt(FloatMatrix &answer, GaussPoint *gp)=0

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