OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
idmnl1.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 idmnl1_h
36 #define idmnl1_h
37 
41 
43 
44 #define _IFT_IDNLMaterial_Name "idmnl1"
45 #define _IFT_IDNLMaterial_r "r"
46 
47 
48 namespace oofem {
49 class GaussPoint;
50 
56 {
57 protected:
60 
61  /* // Variables used to track loading/reloading
62  * public:
63  * enum LastStateType {LST_elastic, LST_loading, LST_unloading};
64  * LastStateType lst;
65  */
66 
67 public:
69  IDNLMaterialStatus(int n, Domain *d, GaussPoint *g);
71  virtual ~IDNLMaterialStatus();
72 
73  virtual void printOutputAt(FILE *file, TimeStep *tStep);
74 
78  void setLocalEquivalentStrainForAverage(double ls) { localEquivalentStrainForAverage = ls; }
79 
80  // definition
81  virtual const char *giveClassName() const { return "IDNLMaterialStatus"; }
82 
83  virtual void initTempStatus();
84  virtual void updateYourself(TimeStep *tStep);
85 
86  virtual contextIOResultType saveContext(DataStream &stream, ContextMode mode, void *obj = NULL);
87  virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj = NULL);
88 
101 };
102 
103 
110 {
111 public:
113  IDNLMaterial(int n, Domain *d);
115  virtual ~IDNLMaterial();
116 
117  // identification and auxiliary functions
118  virtual const char *giveClassName() const { return "IDNLMaterial"; }
119  virtual const char *giveInputRecordName() const { return _IFT_IDNLMaterial_Name; }
120 
122  virtual void giveInputRecord(DynamicInputRecord &input);
123 
125 
126  virtual void computeEquivalentStrain(double &kappa, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep);
136  void computeAngleAndSigmaRatio(double &nx, double &ny, double &ratio, GaussPoint *gp, bool &flag);
147  double computeStressBasedWeight(double &nx, double &ny, double &ratio, GaussPoint *gp, GaussPoint *jGp, double weight);
148  double computeStressBasedWeightForPeriodicCell(double &nx, double &ny, double &ratio, GaussPoint *gp, GaussPoint *jGp);
149 
150  void computeLocalEquivalentStrain(double &kappa, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep)
151  { IsotropicDamageMaterial1 :: computeEquivalentStrain(kappa, strain, gp, tStep); }
152 
153  virtual void updateBeforeNonlocAverage(const FloatArray &strainVector, GaussPoint *gp, TimeStep *tStep);
154 
156  virtual double giveNonlocalMetricModifierAt(GaussPoint *gp);
157 
158  virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep);
159 
160 #ifdef __OOFEG
161  virtual void NonlocalMaterialStiffnessInterface_showSparseMtrxStructure(GaussPoint *gp, oofegGraphicContext &gc, TimeStep *tStep);
163 #endif
164 
165  virtual void computeDamageParam(double &omega, double kappa, const FloatArray &strain, GaussPoint *gp);
166 
169  virtual void NonlocalMaterialStiffnessInterface_addIPContribution(SparseMtrx &dest, const UnknownNumberingScheme &s,
170  GaussPoint *gp, TimeStep *tStep);
176  virtual std :: vector< localIntegrationRecord > *NonlocalMaterialStiffnessInterface_giveIntegrationDomainList(GaussPoint *gp);
186  int giveLocalNonlocalStiffnessContribution(GaussPoint *gp, IntArray &loc, const UnknownNumberingScheme &s,
187  FloatArray &lcontrib, TimeStep *tStep);
196  void giveRemoteNonlocalStiffnessContribution(GaussPoint *gp, IntArray &rloc, const UnknownNumberingScheme &s,
197  FloatArray &rcontrib, TimeStep *tStep);
205  void giveNormalElasticStiffnessMatrix(FloatMatrix &answer,
206  MatResponseMode mode,
207  GaussPoint *gp, TimeStep *tStep);
209 
210  virtual int packUnknowns(DataStream &buff, TimeStep *tStep, GaussPoint *ip);
211  virtual int unpackAndUpdateUnknowns(DataStream &buff, TimeStep *tStep, GaussPoint *ip);
212  virtual int estimatePackSize(DataStream &buff, GaussPoint *ip);
213  virtual double predictRelativeComputationalCost(GaussPoint *gp);
214  virtual double predictRelativeRedistributionCost(GaussPoint *gp) { return 1.0; }
215 
217 
218 protected:
219  virtual void initDamaged(double kappa, FloatArray &totalStrainVector, GaussPoint *gp) { }
220 };
221 } // end namespace oofem
222 #endif // idmnl1_h
virtual const char * giveClassName() const
Definition: idmnl1.h:118
InternalStateType
Type representing the physical meaning of element or constitutive model internal variable.
Abstract base class for all nonlocal structural materials.
This class implements associated Material Status to IDNLMaterial (Nonlocal isotropic damage)...
Definition: idmnl1.h:55
GaussPoint * gp
Associated integration point.
Class and object Domain.
Definition: domain.h:115
void setLocalEquivalentStrainForAverage(double ls)
Sets the localEquivalentStrainForAverage to given value.
Definition: idmnl1.h:78
Domain * domain
Link to domain object, useful for communicating with other FEM components.
Definition: femcmpnn.h:82
Base class for all matrices stored in sparse format.
Definition: sparsemtrx.h:60
double kappa
Scalar measure of the largest strain level ever reached in material.
The purpose of DataStream abstract class is to allow to store/restore context to different streams...
Definition: datastream.h:54
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
virtual Interface * giveInterface(InterfaceType)
Interface requesting service.
Definition: idmnl1.C:930
This class implements associated Material Status to IsotropicDamageMaterial1.
Definition: idm1.h:117
virtual void giveInputRecord(DynamicInputRecord &input)
Setups the input record string of receiver.
Definition: femcmpnn.C:77
Class implementing an array of integers.
Definition: intarray.h:61
MatResponseMode
Describes the character of characteristic material matrix.
double giveLocalEquivalentStrainForAverage()
Returns the local equivalent strain to be averaged.
Definition: idmnl1.h:76
virtual const char * giveInputRecordName() const
Definition: idmnl1.h:119
double localEquivalentStrainForAverage
Equivalent strain for averaging.
Definition: idmnl1.h:59
virtual void initDamaged(double kappa, FloatArray &totalStrainVector, GaussPoint *gp)
Performs initialization, when damage first appear.
Definition: idmnl1.h:219
#define _IFT_IDNLMaterial_Name
Definition: idmnl1.h:44
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Definition: matstatus.h:140
virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Restores the receiver state previously written in stream.
Definition: idmnl1.C:912
Class Nonlocal Material Stiffness Interface.
Abstract base class allowing to control the way, how equations are assigned to individual DOFs...
virtual void initTempStatus()
Initializes the temporary internal variables, describing the current state according to previously re...
Definition: idmnl1.C:870
void computeLocalEquivalentStrain(double &kappa, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep)
Definition: idmnl1.h:150
Abstract base class representing a material status information.
Definition: matstatus.h:84
Class representing vector of real numbers.
Definition: floatarray.h:82
FloatArray strainVector
Equilibrated strain vector in reduced form.
Definition: structuralms.h:69
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
Class representing the general Input Record.
Definition: inputrecord.h:101
Base class for all nonlocal structural material statuses.
Class Interface.
Definition: interface.h:82
virtual const char * giveClassName() const
Definition: idmnl1.h:81
Class representing the a dynamic Input Record.
long ContextMode
Context mode (mask), defining the type of information written/read to/from context.
Definition: contextmode.h:43
This class implements a Nonlocal Isotropic Damage Model for Concrete in Tension Model based on nonloc...
Definition: idmnl1.h:108
This class implements a simple local isotropic damage model for concrete in tension.
Definition: idm1.h:137
virtual void updateYourself(TimeStep *tStep)
Update equilibrium history variables according to temp-variables.
Definition: idmnl1.C:882
virtual void printOutputAt(FILE *file, TimeStep *tStep)
Print receiver&#39;s output to given stream.
Definition: idmnl1.C:852
InterfaceType
Enumerative type, used to identify interface type.
Definition: interfacetype.h:43
virtual MaterialStatus * CreateStatus(GaussPoint *gp) const
Creates new copy of associated status and inserts it into given integration point.
Definition: idmnl1.h:216
the oofem namespace is to define a context or scope in which all oofem names are defined.
virtual contextIOResultType saveContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Stores receiver state to output stream.
Definition: idmnl1.C:895
virtual void computeEquivalentStrain(double &kappa, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep)
Computes the equivalent strain measure from given strain vector (full form).
Definition: idm1.C:472
Class representing integration point in finite element program.
Definition: gausspoint.h:93
virtual double predictRelativeRedistributionCost(GaussPoint *gp)
Returns the relative redistribution cost of the receiver.
Definition: idmnl1.h:214
Class representing solution step.
Definition: timestep.h:80
virtual ~IDNLMaterialStatus()
Destructor.
Definition: idmnl1.C:847
IDNLMaterialStatus(int n, Domain *d, GaussPoint *g)
Constructor.
Definition: idmnl1.C:840

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