OOFEM
2.4
OOFEM.org - Object Oriented Finite Element Solver
|
This class is a base class for microproblem. More...
#include <micromaterial.h>
Public Member Functions | |
MicroMaterial (int n, Domain *d) | |
Constructor. More... | |
virtual | ~MicroMaterial () |
Destructor. More... | |
virtual IRResultType | initializeFrom (InputRecord *ir) |
Initializes receiver according to object description stored in input record. More... | |
virtual const char * | giveInputRecordName () const |
virtual const char * | giveClassName () const |
virtual void | giveRealStressVector_3d (FloatArray &answer, GaussPoint *, const FloatArray &, TimeStep *) |
Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress. More... | |
virtual MaterialStatus * | CreateStatus (GaussPoint *gp) const |
Creates new copy of associated status and inserts it into given integration point. More... | |
void | giveMacroStiffnessMatrix (FloatMatrix &answer, TimeStep *tStep, MatResponseMode rMode, const IntArray µMasterNodes, const IntArray µBoundaryNodes) |
void | setMacroProperties (Domain *macroDomain, MacroLSpace *macroLSpaceElement, const IntArray µMasterNodes, const IntArray µBoundaryNodes) |
void | init (void) |
Related to numbering scheme. More... | |
int | giveDofEquationNumber (Dof *dof) const |
Returns the equation number for corresponding DOF. More... | |
virtual bool | isDefault () const |
Returns true, if receiver is the default engngModel equation numbering scheme; This is useful for some components (typically elements), that cache their code numbers for default numbering to avoid repeated evaluation. More... | |
virtual int | giveRequiredNumberOfDomainEquation () const |
Returns required number of domain equation. More... | |
Public Member Functions inherited from oofem::StructuralMaterial | |
StructuralMaterial (int n, Domain *d) | |
Constructor. More... | |
virtual | ~StructuralMaterial () |
Destructor. More... | |
virtual int | hasMaterialModeCapability (MaterialMode mode) |
Tests if material supports material mode. More... | |
virtual void | giveInputRecord (DynamicInputRecord &input) |
Setups the input record string of receiver. More... | |
virtual void | giveStiffnessMatrix (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) |
Computes the stiffness matrix for giveRealStressVector of receiver in given integration point, respecting its history. More... | |
virtual void | giveRealStressVector (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep) |
Computes the real stress vector for given total strain and integration point. More... | |
virtual void | giveRealStressVector_PlaneStrain (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep) |
Default implementation relies on giveRealStressVector_3d. More... | |
virtual void | giveRealStressVector_StressControl (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, const IntArray &strainControl, TimeStep *tStep) |
Iteratively calls giveRealStressVector_3d to find the stress controlled equal to zero· More... | |
virtual void | giveRealStressVector_ShellStressControl (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, const IntArray &strainControl, TimeStep *tStep) |
virtual void | giveRealStressVector_PlaneStress (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep) |
Default implementation relies on giveRealStressVector_StressControl. More... | |
virtual void | giveRealStressVector_1d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep) |
Default implementation relies on giveRealStressVector_StressControl. More... | |
virtual void | giveRealStressVector_Warping (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep) |
Default implementation relies on giveRealStressVector_StressControl. More... | |
virtual void | giveRealStressVector_2dBeamLayer (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep) |
Default implementation relies on giveRealStressVector_StressControl. More... | |
virtual void | giveRealStressVector_PlateLayer (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep) |
Default implementation relies on giveRealStressVector_StressControl. More... | |
virtual void | giveRealStressVector_Fiber (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep) |
Default implementation relies on giveRealStressVector_StressControl. More... | |
virtual void | giveRealStressVector_Lattice2d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep) |
virtual void | giveRealStressVector_Lattice3d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep) |
virtual void | giveRealStressVector_2dPlateSubSoil (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep) |
Default implementation is not provided. More... | |
virtual void | giveRealStressVector_3dBeamSubSoil (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep) |
virtual void | giveEshelbyStressVector_PlaneStrain (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
Prototype for computation of Eshelby stress. More... | |
void | give_dPdF_from (const FloatMatrix &dSdE, FloatMatrix &answer, GaussPoint *gp) |
void | convert_dSdE_2_dPdF (FloatMatrix &answer, const FloatMatrix &dSdE, const FloatArray &S, const FloatArray &F, MaterialMode matMode) |
virtual void | giveThermalDilatationVector (FloatArray &answer, GaussPoint *gp, TimeStep *tStep) |
Returns a vector of coefficients of thermal dilatation in direction of each material principal (local) axis. More... | |
double | giveReferenceTemperature () |
Returns the reference temperature of receiver. More... | |
virtual void | computeStressIndependentStrainVector (FloatArray &answer, GaussPoint *gp, TimeStep *tStep, ValueModeType mode) |
Computes reduced strain vector in given integration point, generated by internal processes in material, which are independent on loading in particular integration point. More... | |
virtual void | computeStressIndependentStrainVector_3d (FloatArray &answer, GaussPoint *gp, TimeStep *tStep, ValueModeType mode) |
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 in integration point. More... | |
virtual void | give3dMaterialStiffnessMatrix_dPdF (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) |
virtual void | give3dMaterialStiffnessMatrix_dCde (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) |
void | giveStressDependentPartOfStrainVector (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrainVector, TimeStep *tStep, ValueModeType mode) |
Method for subtracting from reduced space strain vector its stress-independent parts (caused by temperature, shrinkage, creep and possibly by other phenomena). More... | |
void | giveStressDependentPartOfStrainVector_3d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrainVector, TimeStep *tStep, ValueModeType mode) |
virtual int | setIPValue (const FloatArray &value, GaussPoint *gp, InternalStateType type) |
Sets the value of a certain variable at a given integration point to the given value. More... | |
virtual int | giveIPValue (FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) |
Returns the integration point corresponding value in Reduced form. More... | |
virtual void | give2dBeamLayerStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
Method for computing 2d beam layer stiffness matrix of receiver. More... | |
virtual void | givePlateLayerStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
Method for computing 2d plate layer stiffness matrix of receiver. More... | |
virtual void | giveFiberStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
Method for computing 1d fiber stiffness matrix of receiver. More... | |
virtual void | give2dLatticeStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
Method for computing 2d lattice stiffness matrix of receiver. More... | |
virtual void | give3dLatticeStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
Method for computing 3d lattice stiffness matrix of receiver. More... | |
virtual void | give2dPlateSubSoilStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
Method for computing stiffness matrix of plate subsoil model. More... | |
virtual void | give3dBeamSubSoilStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
Method for computing stiffness matrix of beam3d subsoil model. More... | |
virtual void | giveFirstPKStressVector_3d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress. More... | |
virtual void | giveFirstPKStressVector_PlaneStrain (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
Default implementation relies on giveFirstPKStressVector_3d. More... | |
virtual void | giveFirstPKStressVector_PlaneStress (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
Default implementation relies on giveFirstPKStressVector_3d. More... | |
virtual void | giveFirstPKStressVector_1d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
Default implementation relies on giveFirstPKStressVector_3d. More... | |
virtual void | giveCauchyStressVector_3d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
virtual void | giveCauchyStressVector_PlaneStrain (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
virtual void | giveCauchyStressVector_PlaneStress (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
virtual void | giveCauchyStressVector_1d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
virtual void | givePlaneStressStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
Method for computing plane stress stiffness matrix of receiver. More... | |
virtual void | givePlaneStressStiffMtrx_dPdF (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
virtual void | givePlaneStressStiffMtrx_dCde (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
virtual void | givePlaneStrainStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
Method for computing plane strain stiffness matrix of receiver. More... | |
virtual void | givePlaneStrainStiffMtrx_dPdF (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
virtual void | givePlaneStrainStiffMtrx_dCde (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
virtual void | give1dStressStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
Method for computing 1d stiffness matrix of receiver. More... | |
virtual void | give1dStressStiffMtrx_dPdF (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
virtual void | give1dStressStiffMtrx_dCde (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
Public Member Functions inherited from oofem::Material | |
Material (int n, Domain *d) | |
Constructor. More... | |
virtual | ~Material () |
Destructor. More... | |
virtual bool | isCharacteristicMtrxSymmetric (MatResponseMode rMode) |
Returns true if stiffness matrix of receiver is symmetric Default implementation returns true. More... | |
virtual double | give (int aProperty, GaussPoint *gp) |
Returns the value of material property 'aProperty'. More... | |
virtual bool | hasProperty (int aProperty, GaussPoint *gp) |
Returns true if 'aProperty' exists on material. More... | |
virtual void | modifyProperty (int aProperty, double value, GaussPoint *gp) |
Modify 'aProperty', which already exists on material. More... | |
double | giveCastingTime () |
virtual bool | isActivated (TimeStep *tStep) |
virtual int | hasNonLinearBehaviour () |
Returns nonzero if receiver is non linear. More... | |
virtual int | hasCastingTimeSupport () |
Tests if material supports casting time. More... | |
virtual void | printYourself () |
Prints receiver state on stdout. Useful for debugging. More... | |
virtual contextIOResultType | saveIPContext (DataStream &stream, ContextMode mode, GaussPoint *gp) |
Stores integration point state to output stream. More... | |
virtual contextIOResultType | restoreIPContext (DataStream &stream, ContextMode mode, GaussPoint *gp) |
Reads integration point state to output stream. More... | |
virtual int | checkConsistency () |
Allows programmer to test some internal data, before computation begins. More... | |
virtual int | initMaterial (Element *element) |
Optional function to call specific procedures when initializing a material. More... | |
virtual MaterialStatus * | giveStatus (GaussPoint *gp) const |
Returns material status of receiver in given integration point. More... | |
virtual int | packUnknowns (DataStream &buff, TimeStep *tStep, GaussPoint *ip) |
Pack all necessary data of integration point (according to element parallel_mode) into given communication buffer. More... | |
virtual int | unpackAndUpdateUnknowns (DataStream &buff, TimeStep *tStep, GaussPoint *ip) |
Unpack and updates all necessary data of given integration point (according to element parallel_mode) into given communication buffer. More... | |
virtual int | estimatePackSize (DataStream &buff, GaussPoint *ip) |
Estimates the necessary pack size to hold all packed data of receiver. More... | |
virtual double | predictRelativeComputationalCost (GaussPoint *gp) |
Returns the weight representing relative computational cost of receiver The reference material model is linear isotropic material - its weight is set to 1.0 The other material models should compare to this reference model. More... | |
virtual double | predictRelativeRedistributionCost (GaussPoint *gp) |
Returns the relative redistribution cost of the receiver. More... | |
virtual void | initTempStatus (GaussPoint *gp) |
Initializes temporary variables stored in integration point status at the beginning of new time step. More... | |
Public Member Functions inherited from oofem::FEMComponent | |
FEMComponent (int n, Domain *d) | |
Regular constructor, creates component with given number and belonging to given domain. More... | |
virtual | ~FEMComponent () |
Virtual destructor. More... | |
Domain * | giveDomain () const |
virtual void | setDomain (Domain *d) |
Sets associated Domain. More... | |
int | giveNumber () const |
void | setNumber (int num) |
Sets number of receiver. More... | |
virtual void | updateLocalNumbering (EntityRenumberingFunctor &f) |
Local renumbering support. More... | |
virtual contextIOResultType | saveContext (DataStream &stream, ContextMode mode, void *obj=NULL) |
Stores receiver state to output stream. More... | |
virtual contextIOResultType | restoreContext (DataStream &stream, ContextMode mode, void *obj=NULL) |
Restores the receiver state previously written in stream. More... | |
virtual void | printOutputAt (FILE *file, TimeStep *tStep) |
Prints output of receiver to stream, for given time step. More... | |
virtual Interface * | giveInterface (InterfaceType t) |
Interface requesting service. More... | |
std::string | errorInfo (const char *func) const |
Returns string for prepending output (used by error reporting macros). More... | |
Public Member Functions inherited from oofem::UnknownNumberingScheme | |
UnknownNumberingScheme (void) | |
virtual | ~UnknownNumberingScheme () |
Public Attributes | |
std::string | inputFileNameMicro |
EngngModel * | problemMicro |
Pointer to the underlying micro problem. More... | |
Domain * | macroDomain |
Pointer to the macroscale domain. More... | |
MacroLSpace * | macroLSpaceElement |
Pointer to the macroscale element. More... | |
std::vector< FloatArray > | microMasterCoords |
Array containing coordinates of 8 master nodes of microproblem. More... | |
int ** | microBoundaryDofs |
Array containing equation numbers for boundary nodes [DofManagerNumber][DOF]. More... | |
IntArray | microBoundaryDofsArr |
Array of equation numbers associated to boundary nodes. More... | |
int ** | microInternalDofs |
Array containing equation numbers for internal nodes to be condensed out [DofManagerNumber][DOF]. More... | |
IntArray | microInternalDofsArr |
Array of equation numbers associated to internal nodes. More... | |
int ** | microDefaultDofs |
Array containing default equation numbers for all nodes [DofManagerNumber][DOF]. More... | |
bool | microMatIsUsed |
Flag signalizing whether micromaterial is used by other element. More... | |
Protected Types | |
enum | EquationNumbering { AllNodes, BoundaryNodes, InteriorNodes } |
Protected Attributes | |
bool | isDefaultNumbering |
int | maxNumberOfDomainEquation |
The maximum DOFs corresponding to released all of the boundary conditions. More... | |
int | reqNumberOfDomainEquation |
Required number of domain equations. More... | |
int | NumberOfDofManagers |
Number of DOF Managers. More... | |
EquationNumbering | DofEquationNumbering |
int | totalBoundaryDofs |
Number of equations associated with boundary nodes. More... | |
int | totalInternalDofs |
Number of equations associated with boundary nodes. More... | |
Protected Attributes inherited from oofem::StructuralMaterial | |
double | referenceTemperature |
Reference temperature (temperature, when material has been built into structure). More... | |
Protected Attributes inherited from oofem::Material | |
Dictionary | propertyDictionary |
Property dictionary. More... | |
double | castingTime |
Casting time. More... | |
Protected Attributes inherited from oofem::FEMComponent | |
int | number |
Component number. More... | |
Domain * | domain |
Link to domain object, useful for communicating with other FEM components. More... | |
Additional Inherited Members | |
Static Public Member Functions inherited from oofem::StructuralMaterial | |
static int | giveSymVI (int ind1, int ind2) |
static int | giveVI (int ind1, int ind2) |
static int | giveVoigtVectorMask (IntArray &answer, MaterialMode mmode) |
Returns a mask of the vector indicies corresponding to components in a general (non-symmetric) second order tensor of some stress/strain/deformation measure that performes work. More... | |
static int | giveVoigtSymVectorMask (IntArray &answer, MaterialMode mmode) |
The same as giveVoigtVectorMask but returns a mask corresponding to a symmetric second order tensor. More... | |
static void | giveInvertedVoigtVectorMask (IntArray &answer, MaterialMode mmode) |
Gives the inverted version of giveVoigtVectorMask. More... | |
static int | giveSizeOfVoigtVector (MaterialMode mmode) |
Returns the size of reduced stress/strain vector according to given mode. More... | |
static int | giveSizeOfVoigtSymVector (MaterialMode mmode) |
Returns the size of symmetric part of a reduced stress/strain vector according to given mode. More... | |
static void | giveFullVectorForm (FloatArray &answer, const FloatArray &strainVector, MaterialMode matMode) |
Converts the reduced symmetric Voigt vector (2nd order tensor) to full form. More... | |
static void | giveFullVectorFormF (FloatArray &answer, const FloatArray &strainVector, MaterialMode matMode) |
Converts the reduced deformation gradient Voigt vector (2nd order tensor). More... | |
static void | giveFullSymVectorForm (FloatArray &answer, const FloatArray &vec, MaterialMode matMode) |
Converts the reduced unsymmetric Voigt vector (2nd order tensor) to full form. More... | |
static void | giveReducedVectorForm (FloatArray &answer, const FloatArray &vec, MaterialMode matMode) |
Converts the full symmetric Voigt vector (2nd order tensor) to reduced form. More... | |
static void | giveReducedSymVectorForm (FloatArray &answer, const FloatArray &vec, MaterialMode matMode) |
Converts the full unsymmetric Voigt vector (2nd order tensor) to reduced form. More... | |
static void | giveFullSymMatrixForm (FloatMatrix &answer, const FloatMatrix &red, MaterialMode matMode) |
Converts the full unsymmetric Voigt matrix (4th order tensor) to reduced form. More... | |
static void | giveReducedMatrixForm (FloatMatrix &answer, const FloatMatrix &full, MaterialMode matMode) |
Converts the full symmetric Voigt matrix (4th order tensor) to reduced form. More... | |
static void | giveReducedSymMatrixForm (FloatMatrix &answer, const FloatMatrix &full, MaterialMode matMode) |
Converts the full unsymmetric Voigt matrix (4th order tensor) to reduced form. More... | |
static void | transformStrainVectorTo (FloatArray &answer, const FloatMatrix &base, const FloatArray &strainVector, bool transpose=false) |
Transforms 3d strain vector into another coordinate system. More... | |
static void | transformStressVectorTo (FloatArray &answer, const FloatMatrix &base, const FloatArray &stressVector, bool transpose=false) |
Transforms 3d stress vector into another coordinate system. More... | |
static double | computeVonMisesStress (const FloatArray *currentStress) |
Computes equivalent of von Mises stress. More... | |
static void | giveStrainVectorTranformationMtrx (FloatMatrix &answer, const FloatMatrix &base, bool transpose=false) |
Computes 3d strain vector transformation matrix from standard vector transformation matrix. More... | |
static void | give2DStrainVectorTranformationMtrx (FloatMatrix &answer, const FloatMatrix &base, bool transpose=false) |
Computes 2d strain vector transformation matrix from standard vector transformation matrix. More... | |
static void | giveStressVectorTranformationMtrx (FloatMatrix &answer, const FloatMatrix &base, bool transpose=false) |
Computes 3d stress vector transformation matrix from standard vector transformation matrix. More... | |
static void | givePlaneStressVectorTranformationMtrx (FloatMatrix &answer, const FloatMatrix &base, bool transpose=false) |
Computes 2d stress vector transformation matrix from standard vector transformation matrix. More... | |
static void | sortPrincDirAndValCloseTo (FloatArray *pVal, FloatMatrix *pDir, FloatMatrix *toPDir) |
Method for sorting newly computed principal values (pVal) and corresponding principal directions (pDir) to be closed to some (often previous) principal directions (toPDir). More... | |
static void | computePrincipalValues (FloatArray &answer, const FloatArray &s, stressStrainPrincMode mode) |
Common functions for convenience. More... | |
static void | computePrincipalValDir (FloatArray &answer, FloatMatrix &dir, const FloatArray &s, stressStrainPrincMode mode) |
Computes principal values and directions of stress or strain vector. More... | |
static double | computeDeviatoricVolumetricSplit (FloatArray &dev, const FloatArray &s) |
Computes split of receiver into deviatoric and volumetric part. More... | |
static void | computeDeviatoricVolumetricSum (FloatArray &s, const FloatArray &dev, double mean) |
static void | applyDeviatoricElasticCompliance (FloatArray &strain, const FloatArray &stress, double EModulus, double nu) |
static void | applyDeviatoricElasticCompliance (FloatArray &strain, const FloatArray &stress, double GModulus) |
static void | applyDeviatoricElasticStiffness (FloatArray &stress, const FloatArray &strain, double EModulus, double nu) |
static void | applyDeviatoricElasticStiffness (FloatArray &stress, const FloatArray &strain, double GModulus) |
static void | applyElasticStiffness (FloatArray &stress, const FloatArray &strain, double EModulus, double nu) |
static void | applyElasticCompliance (FloatArray &strain, const FloatArray &stress, double EModulus, double nu) |
static double | computeStressNorm (const FloatArray &stress) |
static double | computeFirstInvariant (const FloatArray &s) |
static double | computeSecondStressInvariant (const FloatArray &s) |
static double | computeThirdStressInvariant (const FloatArray &s) |
static double | computeFirstCoordinate (const FloatArray &s) |
static double | computeSecondCoordinate (const FloatArray &s) |
static double | computeThirdCoordinate (const FloatArray &s) |
Static Public Attributes inherited from oofem::StructuralMaterial | |
static std::vector< std::vector< int > > | vIindex |
Voigt index map. More... | |
static std::vector< std::vector< int > > | svIndex |
Symmetric Voigt index map. More... | |
This class is a base class for microproblem.
The microproblem represents itself a problem which is solved separately from the macroproblem with appropriate boundary conditions. MacroLspace needs stiffness matrix derived from this microproblem. For this purpose, natural boundary conditions on microproblem have to be excluded. Stiffness matrix of microproblem is condensed to provide stiffness matrix for macroelement.
Definition at line 89 of file micromaterial.h.
|
protected |
Enumerator | |
---|---|
AllNodes | |
BoundaryNodes | |
InteriorNodes |
Definition at line 150 of file micromaterial.h.
oofem::MicroMaterial::MicroMaterial | ( | int | n, |
Domain * | d | ||
) |
Constructor.
Definition at line 95 of file micromaterial.C.
References AllNodes, oofem::IntArray::clear(), DofEquationNumbering, isDefaultNumbering, microBoundaryDofs, microBoundaryDofsArr, microDefaultDofs, microInternalDofs, microInternalDofsArr, microMatIsUsed, and problemMicro.
|
virtual |
Destructor.
Definition at line 108 of file micromaterial.C.
References microBoundaryDofs, microDefaultDofs, microInternalDofs, NumberOfDofManagers, and problemMicro.
|
virtual |
Creates new copy of associated status and inserts it into given integration point.
gp | Integration point where newly created status will be stored. |
Reimplemented from oofem::Material.
Definition at line 227 of file micromaterial.C.
References oofem::FEMComponent::giveDomain().
|
inlinevirtual |
Reimplemented from oofem::StructuralMaterial.
Definition at line 102 of file micromaterial.h.
References oofem::IntegrationPointStatus::gp.
|
virtual |
Returns the equation number for corresponding DOF.
The numbering should return nonzero value if the equation is assigned to the given DOF, zero otherwise.
Implements oofem::UnknownNumberingScheme.
Definition at line 608 of file micromaterial.C.
References AllNodes, BoundaryNodes, DofEquationNumbering, oofem::Dof::giveDofID(), oofem::Dof::giveDofManNumber(), InteriorNodes, microBoundaryDofs, microDefaultDofs, microInternalDofs, and OOFEM_ERROR.
Referenced by giveMacroStiffnessMatrix().
|
inlinevirtual |
Implements oofem::FEMComponent.
Definition at line 101 of file micromaterial.h.
References _IFT_MicroMaterial_Name.
void oofem::MicroMaterial::giveMacroStiffnessMatrix | ( | FloatMatrix & | answer, |
TimeStep * | tStep, | ||
MatResponseMode | rMode, | ||
const IntArray & | microMasterNodes, | ||
const IntArray & | microBoundaryNodes | ||
) |
Definition at line 235 of file micromaterial.C.
References AllNodes, oofem::EngngModel::assemble(), oofem::IntArray::at(), oofem::FloatArray::at(), oofem::FloatMatrix::at(), oofem::SparseMtrx::at(), oofem::SparseMtrx::backSubstitutionWith(), oofem::FloatMatrix::beProductOf(), oofem::FloatMatrix::beTProductOf(), oofem::SparseMtrx::buildInternalStructure(), oofem::classFactory, oofem::ClassFactory::createSparseMtrx(), DofEquationNumbering, oofem::MacroLSpace::evalInterpolation(), oofem::SparseMtrx::factorized(), oofem::IntArray::findFirstIndexOf(), oofem::DofManager::giveCoordinates(), giveDofEquationNumber(), oofem::Domain::giveDofManager(), oofem::DofManager::giveDofWithID(), oofem::EngngModel::giveDomain(), oofem::Domain::giveEngngModel(), oofem::DofManager::giveGlobalNumber(), oofem::SparseMtrx::giveNumberOfColumns(), oofem::FloatMatrix::giveNumberOfColumns(), oofem::Element::giveNumberOfNodes(), oofem::FloatMatrix::giveNumberOfRows(), giveRequiredNumberOfDomainEquation(), oofem::IntArray::giveSize(), InteriorNodes, oofem::SparseMtrx::isAllocatedAt(), isDefaultNumbering, macroLSpaceElement, maxNumberOfDomainEquation, microBoundaryDofsArr, microInternalDofsArr, microMasterCoords, OOFEM_ERROR, OOFEM_LOG_INFO, problemMicro, reqNumberOfDomainEquation, oofem::FloatMatrix::resize(), totalBoundaryDofs, totalInternalDofs, oofem::SparseMtrx::zero(), and oofem::FloatMatrix::zero().
Referenced by oofem::MacroLSpace::computeStiffnessMatrix().
|
virtual |
Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress.
Reimplemented from oofem::StructuralMaterial.
Definition at line 162 of file micromaterial.C.
References OOFEM_ERROR.
|
virtual |
Returns required number of domain equation.
Number is always less or equal to the sum of all DOFs gathered from all nodes.
Reimplemented from oofem::UnknownNumberingScheme.
Definition at line 637 of file micromaterial.C.
References reqNumberOfDomainEquation.
Referenced by giveMacroStiffnessMatrix().
|
virtual |
Related to numbering scheme.
Reimplemented from oofem::UnknownNumberingScheme.
Definition at line 597 of file micromaterial.C.
References oofem::EngngModel::giveNumberOfDomains(), microMatIsUsed, OOFEM_ERROR, and problemMicro.
Referenced by oofem::MacroLSpace::computeStiffnessMatrix().
|
virtual |
Initializes receiver according to object description stored in input record.
This function is called immediately after creating object using constructor. Input record can be imagined as data record in component database belonging to receiver. Receiver may use value-name extracting functions to extract particular field from record.
ir | Input record to initialize from. |
Reimplemented from oofem::StructuralMaterial.
Definition at line 141 of file micromaterial.C.
References _IFT_MicroMaterial_fileName, oofem::_processor, inputFileNameMicro, oofem::InstanciateProblem(), IR_GIVE_FIELD, oofem::IRRT_OK, OOFEM_LOG_INFO, and problemMicro.
|
inlinevirtual |
Returns true, if receiver is the default engngModel equation numbering scheme; This is useful for some components (typically elements), that cache their code numbers for default numbering to avoid repeated evaluation.
Reimplemented from oofem::UnknownNumberingScheme.
Definition at line 124 of file micromaterial.h.
void oofem::MicroMaterial::setMacroProperties | ( | Domain * | macroDomain, |
MacroLSpace * | macroLSpaceElement, | ||
const IntArray & | microMasterNodes, | ||
const IntArray & | microBoundaryNodes | ||
) |
Definition at line 511 of file micromaterial.C.
References oofem::IntArray::contains(), oofem::IntArray::followedBy(), oofem::Node::giveCoordinates(), oofem::EngngModel::giveCurrentMetaStep(), oofem::Domain::giveDofManager(), oofem::EngngModel::giveDomain(), oofem::Domain::giveEngngModel(), oofem::DofManager::giveGlobalNumber(), oofem::EngngModel::giveMetaStep(), oofem::EngngModel::giveNextStep(), oofem::Domain::giveNode(), oofem::Domain::giveNumberOfDofManagers(), oofem::DofManager::giveNumberOfDofs(), oofem::MetaStep::giveNumberOfSteps(), oofem::IntArray::giveSize(), macroDomain, macroLSpaceElement, maxNumberOfDomainEquation, microBoundaryDofs, microBoundaryDofsArr, microDefaultDofs, microInternalDofs, microInternalDofsArr, microMasterCoords, NumberOfDofManagers, OOFEM_ERROR, problemMicro, oofem::MetaStep::setNumberOfSteps(), totalBoundaryDofs, and totalInternalDofs.
Referenced by oofem::MacroLSpace::computeStiffnessMatrix().
|
protected |
Definition at line 151 of file micromaterial.h.
Referenced by giveDofEquationNumber(), giveMacroStiffnessMatrix(), and MicroMaterial().
std :: string oofem::MicroMaterial::inputFileNameMicro |
Definition at line 97 of file micromaterial.h.
Referenced by initializeFrom().
|
protected |
Definition at line 143 of file micromaterial.h.
Referenced by giveMacroStiffnessMatrix(), and MicroMaterial().
Domain* oofem::MicroMaterial::macroDomain |
Pointer to the macroscale domain.
Definition at line 116 of file micromaterial.h.
Referenced by setMacroProperties().
MacroLSpace* oofem::MicroMaterial::macroLSpaceElement |
Pointer to the macroscale element.
Definition at line 119 of file micromaterial.h.
Referenced by giveMacroStiffnessMatrix(), and setMacroProperties().
|
protected |
The maximum DOFs corresponding to released all of the boundary conditions.
Definition at line 145 of file micromaterial.h.
Referenced by giveMacroStiffnessMatrix(), and setMacroProperties().
int** oofem::MicroMaterial::microBoundaryDofs |
Array containing equation numbers for boundary nodes [DofManagerNumber][DOF].
Definition at line 130 of file micromaterial.h.
Referenced by giveDofEquationNumber(), MicroMaterial(), setMacroProperties(), and ~MicroMaterial().
IntArray oofem::MicroMaterial::microBoundaryDofsArr |
Array of equation numbers associated to boundary nodes.
Definition at line 132 of file micromaterial.h.
Referenced by giveMacroStiffnessMatrix(), MicroMaterial(), and setMacroProperties().
int** oofem::MicroMaterial::microDefaultDofs |
Array containing default equation numbers for all nodes [DofManagerNumber][DOF].
Definition at line 138 of file micromaterial.h.
Referenced by giveDofEquationNumber(), MicroMaterial(), setMacroProperties(), and ~MicroMaterial().
int** oofem::MicroMaterial::microInternalDofs |
Array containing equation numbers for internal nodes to be condensed out [DofManagerNumber][DOF].
Definition at line 134 of file micromaterial.h.
Referenced by giveDofEquationNumber(), MicroMaterial(), setMacroProperties(), and ~MicroMaterial().
IntArray oofem::MicroMaterial::microInternalDofsArr |
Array of equation numbers associated to internal nodes.
Definition at line 136 of file micromaterial.h.
Referenced by giveMacroStiffnessMatrix(), MicroMaterial(), and setMacroProperties().
std::vector< FloatArray > oofem::MicroMaterial::microMasterCoords |
Array containing coordinates of 8 master nodes of microproblem.
Definition at line 128 of file micromaterial.h.
Referenced by oofem::MacroLSpace::changeMicroBoundaryConditions(), oofem::MacroLSpace::giveInternalForcesVector(), giveMacroStiffnessMatrix(), and setMacroProperties().
bool oofem::MicroMaterial::microMatIsUsed |
Flag signalizing whether micromaterial is used by other element.
Definition at line 140 of file micromaterial.h.
Referenced by oofem::MacroLSpace::computeStiffnessMatrix(), init(), and MicroMaterial().
|
protected |
Number of DOF Managers.
Definition at line 149 of file micromaterial.h.
Referenced by setMacroProperties(), and ~MicroMaterial().
EngngModel* oofem::MicroMaterial::problemMicro |
Pointer to the underlying micro problem.
Definition at line 113 of file micromaterial.h.
Referenced by oofem::MacroLSpace::computeStiffnessMatrix(), oofem::MacroLSpace::giveInternalForcesVector(), giveMacroStiffnessMatrix(), init(), initializeFrom(), MicroMaterial(), setMacroProperties(), and ~MicroMaterial().
|
protected |
Required number of domain equations.
Definition at line 147 of file micromaterial.h.
Referenced by giveMacroStiffnessMatrix(), and giveRequiredNumberOfDomainEquation().
|
protected |
Number of equations associated with boundary nodes.
Definition at line 153 of file micromaterial.h.
Referenced by giveMacroStiffnessMatrix(), and setMacroProperties().
|
protected |
Number of equations associated with boundary nodes.
Definition at line 155 of file micromaterial.h.
Referenced by giveMacroStiffnessMatrix(), and setMacroProperties().