OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
AbaqusUserElement.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 abaqususerelement_h
36 #define abaqususerelement_h
37 
38 #include "structuralelement.h"
39 #include "nlstructuralelement.h"
40 #include "floatarray.h"
41 #include "floatmatrix.h"
42 #include "timestep.h"
43 
45 
46 #define _IFT_AbaqusUserElement_Name "abaqususerelement"
47 #define _IFT_AbaqusUserElement_userElement "uel"
48 #define _IFT_AbaqusUserElement_numcoords "coords"
49 #define _IFT_AbaqusUserElement_dofs "dofs"
50 #define _IFT_AbaqusUserElement_numsvars "numsvars"
51 #define _IFT_AbaqusUserElement_properties "properties"
52 #define _IFT_AbaqusUserElement_type "type"
53 #define _IFT_AbaqusUserElement_name "name"
54 
55 
56 namespace oofem {
57 
78 {
79 private:
81  void *uelobj;
82 
84  int nCoords;
87 
89  int numSvars;
90 
93 
96 
98  int jtype;
99 
101  int nrhs = 2;
103 
106 
108  int ndofel = 0;
109 
111  int mcrd = 2;
112 
114  int mlvarx = 1;
115 
119 
122 
124  int kinc = 1, kstep = 1;
125 
127  int npredef = 1;
129 
132 
134  int ndLoad = 0;
135  int mdLoad = 0;
136 
138  int *jdltype = NULL; // Temporary init.
139 
141  double params [ 3 ];
142 
145 
148 
150  void (*uel)(double *rhs, double *amatrx, double *svars, double energy [ 8 ], int *ndofel, // 5
151  int *nrhs, int *nsvars, double *props, int *nprops, double *coords, int *mcrd, // 6
152  int *nnode, double *u, double *du, double *v, double *a, int *jtype, // 6
153  double time [ 2 ], double *dtime, int *kstep, int *kinc, int *jelem, // 5
154  double params [ 3 ], int *ndload, int *jdltyp, double *adlmag, double *predef, // 5
155  int *npredef, int *lflags, int *mvarx, double *ddlmag, int *mdload, // 5
156  double *pnewdt, int *jprops, int *njprop, double *period); // 4 - tot 36
157 
159  std :: string filename;
160 
161 public:
163  AbaqusUserElement(int n, Domain *d);
164 
166  virtual ~AbaqusUserElement();
167 
168  virtual void computeConsistentMassMatrix(FloatMatrix &answer, TimeStep *tStep, double &mass, const double *ipDensity = NULL);
169  //virtual void computeInitialStressMatrix(FloatMatrix &answer, TimeStep *tStep);
170  virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep);
171  virtual void giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord = 0);
172  virtual void giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, FloatArray &U, FloatMatrix &DU, int useUpdatedGpRecord);
173  virtual int computeNumberOfDofs() { return this->ndofel; }
174  virtual void giveDofManDofIDMask(int inode, IntArray &answer) const;
175  virtual void computeField(ValueModeType mode, TimeStep *tStep, const FloatArray &lcoords, FloatArray &answer)
176  { OOFEM_ERROR("Abaqus user element cannot support computing local unknown vector\n"); }
177  virtual void updateYourself(TimeStep *tStep);
178  virtual void updateInternalState(TimeStep *tStep);
179 
180  bool hasTangent() const {
181  return hasTangentFlag;
182  }
183  virtual void letTempTangentBe(FloatMatrix &src) {
184  tempAmatrx = src;
185  hasTangentFlag = true;
186  }
188  return tempRHS = src;
189  }
191  return tempSvars = src;
192  }
193  virtual const FloatArray &giveStateVector() const {
194  return svars;
195  }
196  virtual const FloatArray &giveTempStateVector() const {
197  return tempSvars;
198  }
199  virtual const FloatMatrix &giveTempTangent() {
200  return tempAmatrx;
201  }
202 
204 
206  virtual void giveInputRecord(DynamicInputRecord &input);
207  virtual void postInitialize();
208 
209  // definition & identification
210  virtual const char *giveClassName() const { return "AbaqusUserElement"; }
211  virtual const char *giveInputRecordName() const { return _IFT_AbaqusUserElement_Name; }
214  }
216  return EGT_line_1;
217  }
218 
219 protected:
220  virtual void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep) {
221  OOFEM_ERROR("Function not defined for AbaqusUserElement and should never be called. This is a bug.");
222  }
223  virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) {
224  OOFEM_ERROR("Function not defined for AbaqusUserElement and should never be called. This is a bug.");
225  }
226 
227  virtual void computeBmatrixAt(GaussPoint *, FloatMatrix &, int = 1, int = ALL_STRAINS) {
228  OOFEM_ERROR("Function not defined for AbaqusUserElement and should never be called. This is a bug.");
229  }
230 
231  virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) {
232  OOFEM_ERROR("Function not defined for AbaqusUserElement and should never be called. This is a bug.");
233  return 0;
234  }
235 };
236 
237 }// namespace oofem
238 
239 #endif // abaqususerelement_h
void(* uel)(double *rhs, double *amatrx, double *svars, double energy[8], int *ndofel, int *nrhs, int *nsvars, double *props, int *nprops, double *coords, int *mcrd, int *nnode, double *u, double *du, double *v, double *a, int *jtype, double time[2], double *dtime, int *kstep, int *kinc, int *jelem, double params[3], int *ndload, int *jdltyp, double *adlmag, double *predef, int *npredef, int *lflags, int *mvarx, double *ddlmag, int *mdload, double *pnewdt, int *jprops, int *njprop, double *period)
Pointer to the dynamically loaded uel-function (translated to C)
InternalStateType
Type representing the physical meaning of element or constitutive model internal variable.
integrationDomain
Used by integrator class to supply integration points for proper domain to be integrated (Area...
virtual const FloatMatrix & giveTempTangent()
virtual void giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord=0)
Evaluates nodal representation of real internal forces.
Class and object Domain.
Definition: domain.h:115
UEL interface from Abaqus user elements.
Element_Geometry_Type
Enumerative type used to classify element geometry Possible values are: EGT_point - point in space EG...
virtual Element_Geometry_Type giveGeometryType() const
Returns the element geometry type.
Abstract base class for "structural" finite elements with geometrical nonlinearities.
virtual const FloatArray & giveTempStateVector() const
std::string filename
File containing the uel function.
virtual void postInitialize()
Performs post initialization steps.
ValueModeType
Type representing the mode of UnknownType or CharType, or similar types.
Definition: valuemodetype.h:78
FloatArray props
Element properties.
virtual void updateInternalState(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.
virtual FloatArray & letTempSvarsBe(FloatArray &src)
virtual void computeBmatrixAt(GaussPoint *, FloatMatrix &, int=1, int=ALL_STRAINS)
Computes the geometrical matrix of receiver in given integration point.
virtual const char * giveInputRecordName() const
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.
int numSvars
Number of status variables.
Class implementing an array of integers.
Definition: intarray.h:61
MatResponseMode
Describes the character of characteristic material matrix.
AbaqusUserElement(int n, Domain *d)
Constructor.
virtual Interface * giveInterface(InterfaceType it)
Interface requesting service.
#define _IFT_AbaqusUserElement_Name
virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
Computes constitutive matrix of receiver.
virtual void giveDofManDofIDMask(int inode, IntArray &answer) const
Returns dofmanager dof mask for node.
void * uelobj
Dynamically loaded uel.
virtual const char * giveClassName() const
FloatMatrix amatrx
Element amatrx.
#define OOFEM_ERROR(...)
Definition: error.h:61
FloatArray U
Inputs to element routines. Velocity and Acceleration currently ignored.
virtual void letTempTangentBe(FloatMatrix &src)
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...
int nCoords
Coord transf.
virtual ~AbaqusUserElement()
Destructor.
#define ALL_STRAINS
virtual void updateYourself(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.
bool hasTangentFlag
Keeps track of whether the tangent has been obtained already.
Class representing vector of real numbers.
Definition: floatarray.h:82
virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
Computes the stiffness matrix of receiver.
virtual void computeField(ValueModeType mode, TimeStep *tStep, const FloatArray &lcoords, FloatArray &answer)
Computes the unknown vector interpolated at the specified local coordinates.
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 int computeNumberOfDofs()
Computes or simply returns total number of element's local DOFs.
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
Class representing the general Input Record.
Definition: inputrecord.h:101
Class Interface.
Definition: interface.h:82
Class representing the a dynamic Input Record.
FloatArray svars
Status variables.
InterfaceType
Enumerative type, used to identify interface type.
Definition: interfacetype.h:43
virtual FloatMatrix & letTempRhsBe(FloatMatrix &src)
the oofem namespace is to define a context or scope in which all oofem names are defined.
virtual const FloatArray & giveStateVector() const
virtual integrationDomain giveIntegrationDomain() const
Returns integration domain for receiver, used to initialize integration point over receiver volume...
virtual void giveInputRecord(DynamicInputRecord &input)
Setups the input record string of receiver.
Class representing integration point in finite element program.
Definition: gausspoint.h:93
Class representing solution step.
Definition: timestep.h:80

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