OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
planstrssxfem.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 planstrssxfem_h
36 #define planstrssxfem_h
37 
38 #include "../sm/Elements/PlaneStress/planstrss.h"
39 #include "../sm/xfem/xfemstructuralelementinterface.h"
40 #include "vtkxmlexportmodule.h"
41 
42 #define _IFT_PlaneStress2dXfem_Name "planestress2dxfem"
43 
44 namespace oofem {
51 {
52 protected:
53  virtual void updateYourself(TimeStep *tStep);
54  virtual void postInitialize();
55 
56 public:
60  virtual ~PlaneStress2dXfem() { }
61 
63 
64  virtual const char *giveInputRecordName() const { return _IFT_PlaneStress2dXfem_Name; }
65  virtual const char *giveClassName() const { return "PlaneStress2dXfem"; }
66  virtual int computeNumberOfDofs();
67  virtual void computeGaussPoints();
68 
69  void increaseNumGP();
70 
71  virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer);
72  virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer,
73  int lowerIndx = 1, int upperIndx = ALL_STRAINS);
74  virtual void computeBHmatrixAt(GaussPoint *gp, FloatMatrix &answer);
75  virtual void giveDofManDofIDMask(int inode, IntArray &answer) const;
76  virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *, TimeStep *tStep);
77  virtual void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep);
78  virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep);
79 
80  virtual void computeDeformationGradientVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep);
81 
82  virtual void giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord);
83 
84  virtual void computeConsistentMassMatrix(FloatMatrix &answer, TimeStep *tStep, double &mass, const double *ipDensity = NULL) { XfemStructuralElementInterface :: XfemElementInterface_computeConsistentMassMatrix(answer, tStep, mass, ipDensity); }
85 
87 
88 #ifdef __OOFEG
89  virtual void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep);
90  //void drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType);
91  virtual void drawScalar(oofegGraphicContext &gc, TimeStep *tStep);
92  //virtual void drawSpecial(oofegGraphicContext &gc, TimeStep *tStep);
93 #endif
94 
97  virtual void giveInputRecord(DynamicInputRecord &input);
98 
100  virtual void giveCompositeExportData(std::vector< VTKPiece > &vtkPieces, IntArray &primaryVarsToExport, IntArray &internalVarsToExport, IntArray cellVarsToExport, TimeStep *tStep);
101 
102 };
103 } // end namespace oofem
104 #endif
virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *, TimeStep *tStep)
Computes constitutive matrix of receiver.
#define _IFT_PlaneStress2dXfem_Name
Definition: planstrssxfem.h:42
virtual void drawScalar(oofegGraphicContext &gc, TimeStep *tStep)
This class implements an isoparametric four-node quadrilateral plane- stress elasticity finite elemen...
Definition: planstrss.h:56
virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int lowerIndx=1, int upperIndx=ALL_STRAINS)
Computes the geometrical matrix of receiver in given integration point.
Class and object Domain.
Definition: domain.h:115
virtual MaterialMode giveMaterialMode()
Returns material mode for receiver integration points.
Element_Geometry_Type
Enumerative type used to classify element geometry Possible values are: EGT_point - point in space EG...
virtual void XfemElementInterface_computeConsistentMassMatrix(FloatMatrix &answer, TimeStep *tStep, double &mass, const double *ipDensity=NULL)
virtual const char * giveInputRecordName() const
Definition: planstrssxfem.h:64
virtual void giveDofManDofIDMask(int inode, IntArray &answer) const
Returns dofmanager dof mask for node.
virtual void postInitialize()
Performs post initialization steps.
Definition: planstrssxfem.C:60
Elements with geometry defined as EGT_Composite are exported using individual pieces.
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
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...
MaterialMode
Type representing material mode of integration point.
Definition: materialmode.h:89
Class implementing an array of integers.
Definition: intarray.h:61
MatResponseMode
Describes the character of characteristic material matrix.
virtual void updateYourself(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.
virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
Computes the stiffness matrix of receiver.
virtual const char * giveClassName() const
Definition: planstrssxfem.h:65
virtual int computeNumberOfDofs()
Computes or simply returns total number of element&#39;s local DOFs.
virtual void computeBHmatrixAt(GaussPoint *gp, FloatMatrix &answer)
Computes a matrix which, multiplied by the column matrix of nodal displacements, gives the displaceme...
PlaneStress2dXfem(int n, Domain *d)
Constructor.
Definition: planstrssxfem.h:58
#define ALL_STRAINS
Provides Xfem interface for a structural element.
Class representing vector of real numbers.
Definition: floatarray.h:82
Implementation of matrix containing floating point numbers.
Definition: floatmatrix.h:94
virtual Element_Geometry_Type giveGeometryType() const
Returns the element geometry type.
IRResultType
Type defining the return values of InputRecord reading operations.
Definition: irresulttype.h:47
virtual void giveInputRecord(DynamicInputRecord &input)
Setups the input record string of receiver.
Class representing the general Input Record.
Definition: inputrecord.h:101
virtual Interface * giveInterface(InterfaceType it)
Interface requesting service.
Definition: planstrssxfem.C:67
Class Interface.
Definition: interface.h:82
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
virtual void giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord)
Evaluates nodal representation of real internal forces.
Class representing the a dynamic Input Record.
virtual void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
InterfaceType
Enumerative type, used to identify interface type.
Definition: interfacetype.h:43
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.
virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer)
Computes interpolation matrix for element unknowns.
the oofem namespace is to define a context or scope in which all oofem names are defined.
virtual void computeGaussPoints()
Initializes the array of integration rules member variable.
Definition: planstrssxfem.C:80
virtual ~PlaneStress2dXfem()
Destructor.
Definition: planstrssxfem.h:60
virtual void computeConsistentMassMatrix(FloatMatrix &answer, TimeStep *tStep, double &mass, const double *ipDensity=NULL)
Computes consistent mass matrix of receiver using numerical integration over element volume...
Definition: planstrssxfem.h:84
Class representing integration point in finite element program.
Definition: gausspoint.h:93
Class representing solution step.
Definition: timestep.h:80
virtual void giveCompositeExportData(std::vector< VTKPiece > &vtkPieces, IntArray &primaryVarsToExport, IntArray &internalVarsToExport, IntArray cellVarsToExport, TimeStep *tStep)
VTK Interface.
int numberOfDofMans
Number of dofmanagers.
Definition: element.h:149
Temporary class for testing in the usual case instead of PlaneStress2dXfem there will be the standard...
Definition: planstrssxfem.h:50

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