OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
abaqususermaterial.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 abaqususermaterial_h
36 #define abaqususermaterial_h
37 
38 #include "../sm/Materials/structuralmaterial.h"
39 #include "../sm/Materials/structuralms.h"
40 #include "floatarray.h"
41 #include "floatmatrix.h"
42 
44 
45 #define _IFT_AbaqusUserMaterial_Name "abaqususermaterial"
46 #define _IFT_AbaqusUserMaterial_numState "numstate"
47 #define _IFT_AbaqusUserMaterial_properties "properties"
48 #define _IFT_AbaqusUserMaterial_initialStress "initialstress"
49 #define _IFT_AbaqusUserMaterial_userMaterial "umat"
50 #define _IFT_AbaqusUserMaterial_name "name"
51 #define _IFT_AbaqusUserMaterial_numericalTangent "numericaltangent"
52 #define _IFT_AbaqusUserMaterial_numericalTangentPerturbation "perturbation"
53 
54 
55 namespace oofem {
81 {
82 private:
84  static int n;
86  void *umatobj;
87 
89  void ( * umat )(double *stress, double *statev, double *ddsdde, double *sse, double *spd, // 5
90  double *scd, double *rpl, double *ddsddt, double *drplde, double *drpldt, // 5
91  double *stran, double *dstran, double time [ 2 ], double *dtime, double *temp, // 4
92  double *dtemp, double predef [ 1 ], double dpred [ 1 ], char cmname [ 80 ], int *ndi, // 5
93  int *nshr, int *ntens, int *nstatv, double *props, int *nprops, double coords [ 3 ], // 6
94  double *drot, double *pnewdt, double *celent, double *dfgrd0, double *dfgrd1, // 5
95  int *noel, int *npt, int *layer, int *kspt, int *kstep, int *kinc); // 6
97  char cmname [ 80 ];
99  int numState;
109 
114 
119 
121  std :: string filename;
122 
123 public:
125  AbaqusUserMaterial(int n, Domain * d);
127  virtual ~AbaqusUserMaterial();
128 
137  virtual void giveInputRecord(DynamicInputRecord &input);
138 
139  virtual MaterialStatus *CreateStatus(GaussPoint *gp) const;
140 
141  virtual void give3dMaterialStiffnessMatrix(FloatMatrix &answer,
142  MatResponseMode mode, GaussPoint *gp, TimeStep *tStep);
143 
144  virtual void give3dMaterialStiffnessMatrix_dPdF(FloatMatrix &answer,
145  MatResponseMode mode,
146  GaussPoint *gp,
147  TimeStep *tStep);
148 
149  virtual void givePlaneStrainStiffMtrx_dPdF(FloatMatrix &answer,
150  MatResponseMode mmode, GaussPoint *gp,
151  TimeStep *tStep);
152 
153  virtual void giveRealStressVector_3d(FloatArray &answer, GaussPoint *gp,
154  const FloatArray &reducedStrain, TimeStep *tStep);
155 
156  virtual void giveFirstPKStressVector_3d(FloatArray &answer, GaussPoint *gp,
157  const FloatArray &reducedF, TimeStep *tStep);
158 
159  virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep);
160 
161  virtual int hasNonLinearBehaviour() { return true; }
162  virtual const char *giveClassName() const { return "AbaqusUserMaterial"; }
163  virtual const char *giveInputRecordName() const { return _IFT_AbaqusUserMaterial_Name; }
164 };
165 
167 {
168 protected:
170  int numState;
177 
180 
181 public:
183  AbaqusUserMaterialStatus(int n, Domain * d, GaussPoint * gp, int numState);
186 
187  virtual void initTempStatus();
188  virtual void updateYourself(TimeStep *tStep);
189 
190  bool hasTangent() const { return hasTangentFlag; }
191 
192  const FloatArray &giveStateVector() const { return stateVector; }
193  FloatArray &letStateVectorBe(FloatArray &s) { return stateVector = s; }
194  const FloatArray &giveTempStateVector() const { return tempStateVector; }
195  FloatArray &letTempStateVectorBe(FloatArray &s) { return tempStateVector = s; }
196  const FloatMatrix &giveTempTangent() { return tempTangent; }
198  tempTangent = std :: move(t);
199  hasTangentFlag = true;
200  }
201 
202 
203 
204  virtual void printOutputAt(FILE *file, TimeStep *tStep);
205 
206  virtual const char *giveClassName() const { return "AbaqusUserMaterialStatus"; }
207 };
208 } // end namespace oofem
209 #endif // abaqususermaterial_h
InternalStateType
Type representing the physical meaning of element or constitutive model internal variable.
int numState
Number of state variables.
virtual void giveRealStressVector_3d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep)
Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress.
Class and object Domain.
Definition: domain.h:115
double mPerturbation
Size of perturbation if numerical tangent is used.
FloatArray & letTempStateVectorBe(FloatArray &s)
AbaqusUserMaterial(int n, Domain *d)
Constructor.
FloatArray tempStateVector
Temporary state vector.
int numState
Size of the state vector.
This class implements a structural material status information.
Definition: structuralms.h:65
char cmname[80]
Name for material routine.
virtual void give3dMaterialStiffnessMatrix_dPdF(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in Reduced form.
virtual void giveInputRecord(DynamicInputRecord &input)
Setups the input record string of receiver.
const FloatArray & giveTempStateVector() const
MatResponseMode
Describes the character of characteristic material matrix.
virtual const char * giveClassName() const
virtual void give3dMaterialStiffnessMatrix(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Computes full 3d material stiffness matrix at given integration point, time, respecting load history ...
static int n
Gausspoint counter.
FloatArray stateVector
General state vector.
virtual ~AbaqusUserMaterialStatus()
Destructor.
bool mUseNumericalTangent
Flag to determine if numerical tangent should be used.
const FloatArray & giveStateVector() const
virtual void giveFirstPKStressVector_3d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress.
std::string filename
Name of the file that contains the umat function.
virtual ~AbaqusUserMaterial()
Destructor.
virtual IRResultType initializeFrom(InputRecord *ir)
Reads the following values;.
Abstract base class representing a material status information.
Definition: matstatus.h:84
Class representing vector of real numbers.
Definition: floatarray.h:82
int mStressInterpretation
Flag to determine how the stress and Jacobian are interpreted.
FloatArray initialStress
Initial stress.
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
FloatArray properties
Material properties.
void * umatobj
Dynamically loaded umat.
Class representing the general Input Record.
Definition: inputrecord.h:101
virtual const char * giveInputRecordName() const
virtual const char * giveClassName() const
FloatArray & letStateVectorBe(FloatArray &s)
virtual void printOutputAt(FILE *file, TimeStep *tStep)
Prints output of receiver to stream, for given time step.
Definition: femcmpnn.h:171
#define _IFT_AbaqusUserMaterial_Name
Class representing the a dynamic Input Record.
Abstract base class for all "structural" constitutive models.
This class allows for custom user materials from Abaqus (UMAT).
FloatMatrix tempTangent
Temporary elastic tangent.
the oofem namespace is to define a context or scope in which all oofem names are defined.
void(* umat)(double *stress, double *statev, double *ddsdde, double *sse, double *spd, double *scd, double *rpl, double *ddsddt, double *drplde, double *drpldt, double *stran, double *dstran, double time[2], double *dtime, double *temp, double *dtemp, double predef[1], double dpred[1], char cmname[80], int *ndi, int *nshr, int *ntens, int *nstatv, double *props, int *nprops, double coords[3], double *drot, double *pnewdt, double *celent, double *dfgrd0, double *dfgrd1, int *noel, int *npt, int *layer, int *kspt, int *kstep, int *kinc)
Pointer to the dynamically loaded umat-function (translated to C)
virtual void initTempStatus(GaussPoint *gp)
Initializes temporary variables stored in integration point status at the beginning of new time step...
Definition: material.C:267
bool hasTangentFlag
Checker to see if tangent has been computed.
virtual MaterialStatus * CreateStatus(GaussPoint *gp) const
Creates new copy of associated status and inserts it into given integration point.
Class representing integration point in finite element program.
Definition: gausspoint.h:93
Class representing solution step.
Definition: timestep.h:80
virtual void givePlaneStrainStiffMtrx_dPdF(FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
virtual int hasNonLinearBehaviour()
Returns nonzero if receiver is non linear.
const FloatMatrix & giveTempTangent()

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