OOFEM
2.4
OOFEM.org - Object Oriented Finite Element Solver
|
This class represents PFEM method for solving incompressible Navier-Stokes equations. More...
#include <pfem.h>
Public Member Functions | |
PFEM (int i, EngngModel *_master=NULL) | |
~PFEM () | |
void | solveYourselfAt (TimeStep *) |
Solves problem for given time step. More... | |
virtual void | updateYourself (TimeStep *tStep) |
Updates nodal values (calls also this->updateDofUnknownsDictionary for updating dofs unknowns dictionaries if model supports changes of static system). More... | |
double | giveUnknownComponent (ValueModeType mode, TimeStep *tStep, Domain *d, Dof *dof) |
Returns requested unknown. More... | |
virtual contextIOResultType | saveContext (DataStream &stream, ContextMode mode) |
Stores the state of model to output stream. More... | |
virtual contextIOResultType | restoreContext (DataStream &stream, ContextMode mode) |
Restores the state of model from output stream. More... | |
TimeStep * | giveNextStep () |
Returns next time step (next to current step) of receiver. More... | |
TimeStep * | giveSolutionStepWhenIcApply () |
NumericalMethod * | giveNumericalMethod (MetaStep *) |
Returns reference to receiver's numerical method. More... | |
virtual void | preInitializeNextStep () |
Removes all elements and call DelaunayTriangulator to build up new mesh with new recognized boundary. More... | |
virtual int | forceEquationNumbering (int id) |
Forces equation renumbering on given domain. More... | |
virtual int | requiresUnknownsDictionaryUpdate () |
Indicates if EngngModel requires Dofs dictionaries to be updated. More... | |
IRResultType | initializeFrom (InputRecord *ir) |
Initialization from given input record. More... | |
virtual int | checkConsistency () |
Allows programmer to test some receiver's internal data, before computation begins. More... | |
const char * | giveClassName () const |
Returns class name of the receiver. More... | |
fMode | giveFormulation () |
Indicates type of non linear computation (total or updated formulation). More... | |
virtual void | printDofOutputAt (FILE *stream, Dof *iDof, TimeStep *atTime) |
DOF printing routine. More... | |
virtual int | giveNumberOfDomainEquations (int, const UnknownNumberingScheme &num) |
Returns number of equations for given domain in active (current time step) time step. More... | |
virtual int | giveNewEquationNumber (int domain, DofIDItem) |
Increases number of equations of receiver's domain and returns newly created equation number. More... | |
virtual int | giveNewPrescribedEquationNumber (int domain, DofIDItem) |
Increases number of prescribed equations of receiver's domain and returns newly created equation number. More... | |
virtual int | giveUnknownDictHashIndx (ValueModeType mode, TimeStep *stepN) |
This method is responsible for computing unique dictionary id (ie hash value) from given valueModeType and time step. More... | |
virtual void | updateDofUnknownsDictionary (DofManager *inode, TimeStep *tStep) |
Updates necessary values in Dofs unknown dictionaries. More... | |
void | updateDofUnknownsDictionaryPressure (DofManager *inode, TimeStep *tStep) |
Writes pressures into the dof unknown dictionaries. More... | |
void | updateDofUnknownsDictionaryVelocities (DofManager *inode, TimeStep *tStep) |
Writes velocities into the dof unknown dictionaries. More... | |
void | resetEquationNumberings () |
Resets the equation numberings as the mesh is recreated in each time step. More... | |
int | giveAssociatedMaterialNumber () |
Returns number of material to be associated with created elements. More... | |
int | giveAssociatedCrossSectionNumber () |
Returns number of cross section to be associated with created elements. More... | |
int | giveAssociatedPressureBC () |
Returns number of zero pressure boundary condition to be prescribed on free surface. More... | |
Public Member Functions inherited from oofem::EngngModel | |
EngngModel (int i, EngngModel *_master=NULL) | |
Constructor. More... | |
virtual | ~EngngModel () |
Destructor. More... | |
EngngModel (const EngngModel &)=delete | |
EngngModel & | operator= (const EngngModel &)=delete |
Domain * | giveDomain (int n) |
Service for accessing particular problem domain. More... | |
void | setDomain (int i, Domain *ptr, bool iDeallocateOld=true) |
Sets i-th domain of receiver. More... | |
int | giveNumberOfDomains () |
Returns number of domains in problem. More... | |
const std::string & | giveDescription () const |
const time_t & | giveStartTime () |
bool | giveSuppressOutput () const |
virtual ErrorEstimator * | giveDomainErrorEstimator (int n) |
Service for accessing ErrorEstimator corresponding to particular domain. More... | |
virtual MaterialInterface * | giveMaterialInterface (int n) |
Returns material interface representation for given domain. More... | |
void | setNumberOfEquations (int id, int neq) |
FILE * | giveOutputStream () |
Returns file descriptor of output file. More... | |
std::string | giveOutputBaseFileName () |
Returns base output file name to which extensions, like .out .vtu .osf should be added. More... | |
std::string | giveReferenceFileName () |
Returns reference file name. More... | |
void | letOutputBaseFileNameBe (const std::string &src) |
Sets the base output file name. More... | |
ContextOutputMode | giveContextOutputMode () |
Returns domain context output mode. More... | |
int | giveContextOutputStep () |
Returns domain context output step. More... | |
void | setContextOutputMode (ContextOutputMode contextMode) |
Sets context output mode of receiver. More... | |
void | setUDContextOutputMode (int cStep) |
Sets user defined context output mode (it sets contextOutputMode to contextOutputMode), setting contextOutputStep to given value. More... | |
void | setProblemMode (problemMode pmode) |
Sets domain mode to given mode. More... | |
void | setParallelMode (bool newParallelFlag) |
Sets the problem to run in parallel (or not). More... | |
problemMode | giveProblemMode () |
Returns domain mode. More... | |
void | setProblemScale (problemScale pscale) |
Sets scale in multiscale simulation. More... | |
problemScale | giveProblemScale () |
Returns scale in multiscale simulation. More... | |
virtual void | setRenumberFlag () |
Sets the renumber flag to true. More... | |
virtual void | resetRenumberFlag () |
Sets the renumber flag to false. More... | |
double | giveSolutionStepTime () |
Returns the user time of the current simulation step in seconds. More... | |
void | giveAnalysisTime (int &rhrs, int &rmin, int &rsec, int &uhrs, int &umin, int &usec) |
Returns the real and user time for the analysis. More... | |
void | terminateAnalysis () |
Performs analysis termination after finishing analysis. More... | |
virtual void | solveYourself () |
Starts solution process. More... | |
virtual void | terminate (TimeStep *tStep) |
Terminates the solution of time step. More... | |
virtual void | doStepOutput (TimeStep *tStep) |
Prints the ouput of the solution step (using virtual this->printOutputAtservice) to the stream detemined using this->giveOutputStream() method and calls exportModuleManager to do output. More... | |
void | saveStepContext (TimeStep *tStep, ContextMode mode) |
Saves context of given solution step, if required (determined using this->giveContextOutputMode() method). More... | |
virtual void | initializeYourself (TimeStep *tStep) |
Provides the opportunity to initialize state variables stored in element integration points according to initial conditions using function initializeYourself() on element level. More... | |
virtual int | initializeAdaptive (int tStepNumber) |
Initializes the newly generated discretization state according to previous solution. More... | |
virtual FieldPtr | giveField (FieldType key, TimeStep *) |
Returns the smart pointer to requested field, Null otherwise. More... | |
EngngModel * | giveMasterEngngModel () |
Returns the master engnmodel. More... | |
virtual double | giveLoadLevel () |
Returns the current load level. More... | |
virtual double | giveEigenValue (int eigNum) |
Only relevant for eigen value analysis. Otherwise returns zero. More... | |
virtual void | setActiveVector (int i) |
Only relevant for eigen value analysis. Otherwise does noting. More... | |
int | updateSharedDofManagers (FloatArray &answer, const UnknownNumberingScheme &s, int ExchangeTag) |
Exchanges necessary remote DofManagers data. More... | |
int | exchangeRemoteElementData (int ExchangeTag) |
Exchanges necessary remote element data with remote partitions. More... | |
virtual int | giveCurrentNumberOfIterations () |
Returns number of iterations that was required to reach equilibrium - used for adaptive step length in staggered problem. More... | |
MPI_Comm | giveParallelComm () |
Returns the communication object of reciever. More... | |
int | packRemoteElementData (ProcessCommunicator &processComm) |
Packs data of local element to be received by their remote counterpart on remote partitions. More... | |
int | unpackRemoteElementData (ProcessCommunicator &processComm) |
Unpacks data for remote elements (which are mirrors of remote partition's local elements). More... | |
int | packDofManagers (ArrayWithNumbering *src, ProcessCommunicator &processComm) |
Packing function for vector values of DofManagers. More... | |
int | unpackDofManagers (ArrayWithNumbering *dest, ProcessCommunicator &processComm) |
Unpacking function for vector values of DofManagers . More... | |
ProblemCommunicator * | giveProblemCommunicator (EngngModelCommType t) |
void | initializeCommMaps (bool forceInit=false) |
virtual int | instanciateYourself (DataReader &dr, InputRecord *ir, const char *outFileName, const char *desc) |
Initializes whole problem according to its description stored in inputStream. More... | |
void | Instanciate_init () |
Initialization of the receiver state (opening the default output stream, empty domain creation, initialization of parallel context, etc) before Initialization form DataReader. More... | |
int | instanciateDomains (DataReader &dr) |
Instanciate problem domains by calling their instanciateYourself() service. More... | |
int | instanciateMetaSteps (DataReader &dr) |
Instanciate problem meta steps by calling their instanciateYourself() service. More... | |
virtual int | instanciateDefaultMetaStep (InputRecord *ir) |
Instanciate default metastep, if nmsteps is zero. More... | |
virtual void | updateAttributes (MetaStep *mStep) |
Update receiver attributes according to step metaStep attributes. More... | |
void | initMetaStepAttributes (MetaStep *mStep) |
Update e-model attributes attributes according to step metaStep attributes. More... | |
virtual void | updateDomainLinks () |
Updates domain links after the domains of receiver have changed. More... | |
MetaStep * | giveCurrentMetaStep () |
Returns current meta step. More... | |
virtual TimeStep * | giveCurrentStep (bool force=false) |
Returns current time step. More... | |
virtual TimeStep * | givePreviousStep (bool force=false) |
Returns previous time step. More... | |
TimeStep * | generateNextStep () |
Generate new time step (and associate metastep). More... | |
virtual TimeStep * | giveSolutionStepWhenIcApply (bool force=false) |
Returns the solution step when Initial Conditions (IC) apply. More... | |
virtual int | giveNumberOfFirstStep (bool force=false) |
Returns number of first time step used by receiver. More... | |
int | giveNumberOfMetaSteps () |
Return number of meta steps. More... | |
MetaStep * | giveMetaStep (int i) |
Returns the i-th meta step. More... | |
int | giveNumberOfSteps (bool force=false) |
Returns total number of steps. More... | |
virtual double | giveEndOfTimeOfInterest () |
Returns end of time interest (time corresponding to end of time integration). More... | |
int | giveNumberOfTimeStepWhenIcApply () |
Returns the time step number, when initial conditions should apply. More... | |
ExportModuleManager * | giveExportModuleManager () |
Returns receiver's export module manager. More... | |
EngngModelTimer * | giveTimer () |
Returns reference to receiver timer (EngngModelTimer). More... | |
std::string | giveContextFileName (int tStepNumber, int stepVersion) const |
Returns the filename for the context file for the given step and version. More... | |
std::string | giveDomainFileName (int domainNum, int domainSerNum) const |
Returns the filename for the given domain (used by adaptivity and restore) More... | |
virtual void | updateComponent (TimeStep *tStep, NumericalCmpn cmpn, Domain *d) |
Updates components mapped to numerical method if necessary during solution process. More... | |
virtual void | initStepIncrements () |
Initializes solution of new time step. More... | |
virtual int | forceEquationNumbering () |
Forces equation renumbering on all domains associated to engng model. More... | |
virtual bool | requiresEquationRenumbering (TimeStep *tStep) |
Returns true if equation renumbering is required for given solution step. More... | |
virtual ParallelContext * | giveParallelContext (int n) |
Returns the parallel context corresponding to given domain (n) and unknown type Default implementation returns i-th context from parallelContextList. More... | |
virtual void | initParallelContexts () |
Creates parallel contexts. More... | |
virtual void | assemble (SparseMtrx &answer, TimeStep *tStep, const MatrixAssembler &ma, const UnknownNumberingScheme &s, Domain *domain) |
Assembles characteristic matrix of required type into given sparse matrix. More... | |
virtual void | assemble (SparseMtrx &answer, TimeStep *tStep, const MatrixAssembler &ma, const UnknownNumberingScheme &r_s, const UnknownNumberingScheme &c_s, Domain *domain) |
Assembles characteristic matrix of required type into given sparse matrix. More... | |
void | assembleVector (FloatArray &answer, TimeStep *tStep, const VectorAssembler &va, ValueModeType mode, const UnknownNumberingScheme &s, Domain *domain, FloatArray *eNorms=NULL) |
Assembles characteristic vector of required type from dofManagers, element, and active boundary conditions, into given vector. More... | |
void | assembleVectorFromDofManagers (FloatArray &answer, TimeStep *tStep, const VectorAssembler &va, ValueModeType mode, const UnknownNumberingScheme &s, Domain *domain, FloatArray *eNorms=NULL) |
Assembles characteristic vector of required type from dofManagers into given vector. More... | |
void | assembleVectorFromElements (FloatArray &answer, TimeStep *tStep, const VectorAssembler &va, ValueModeType mode, const UnknownNumberingScheme &s, Domain *domain, FloatArray *eNorms=NULL) |
Assembles characteristic vector of required type from elements into given vector. More... | |
void | assembleVectorFromBC (FloatArray &answer, TimeStep *tStep, const VectorAssembler &va, ValueModeType mode, const UnknownNumberingScheme &s, Domain *domain, FloatArray *eNorms=NULL) |
Assembles characteristic vector of required type from boundary conditions. More... | |
void | assembleExtrapolatedForces (FloatArray &answer, TimeStep *tStep, CharType type, Domain *domain) |
Assembles the extrapolated internal forces vector, useful for obtaining a good initial guess in nonlinear analysis with Dirichlet boundary conditions. More... | |
void | assemblePrescribedExtrapolatedForces (FloatArray &answer, TimeStep *tStep, CharType type, Domain *domain) |
void | assembleVectorFromContacts (FloatArray &answer, TimeStep *tStep, CharType type, ValueModeType mode, const UnknownNumberingScheme &s, Domain *domain, FloatArray *eNorms=NULL) |
virtual int | checkProblemConsistency () |
Allows programmer to test problem its internal data, before computation begins. More... | |
virtual void | init () |
Initializes the receiver state. More... | |
virtual void | postInitialize () |
Performs post-initialization for all the problem contents (which is called after initializeFrom). More... | |
virtual void | printOutputAt (FILE *file, TimeStep *tStep) |
Prints output of receiver to output domain stream, for given time step. More... | |
virtual void | printOutputAt (FILE *file, TimeStep *tStep, const IntArray &nodeSets, const IntArray &elementSets) |
void | outputNodes (FILE *file, Domain &domain, TimeStep *tStep, int setNum) |
Outputs all nodes in the given set. More... | |
void | outputElements (FILE *file, Domain &domain, TimeStep *tStep, int setNum) |
Outputs all elements in the given set. More... | |
void | printYourself () |
Prints state of receiver. Useful for debugging. More... | |
virtual int | useNonlocalStiffnessOption () |
Returns nonzero if nonlocal stiffness option activated. More... | |
bool | isParallel () const |
Returns true if receiver in parallel mode. More... | |
int | giveRank () const |
Returns domain rank in a group of collaborating processes (0..groupSize-1) More... | |
int | giveNumberOfProcesses () const |
Returns the number of collaborating processes. More... | |
EngngModelContext * | giveContext () |
Context requesting service. More... | |
virtual int | giveNumberOfSlaveProblems () |
Returns number of slave problems. More... | |
virtual EngngModel * | giveSlaveProblem (int i) |
Returns i-th slave problem. More... | |
virtual bool | giveEquationScalingFlag () |
Returns the Equation scaling flag, which is used to indicate that governing equation(s) are scaled, or non-dimensionalized. More... | |
virtual double | giveVariableScale (VarScaleType varId) |
Returns the scale factor for given variable type. More... | |
virtual int | estimateMaxPackSize (IntArray &commMap, DataStream &buff, int packUnpackType) |
Determines the space necessary for send/receive buffer. More... | |
virtual void | balanceLoad (TimeStep *tStep) |
Recovers the load balance between processors, if needed. More... | |
virtual LoadBalancer * | giveLoadBalancer () |
Returns reference to receiver's load balancer. More... | |
virtual LoadBalancerMonitor * | giveLoadBalancerMonitor () |
Returns reference to receiver's load balancer monitor. More... | |
void | initParallel () |
Request domain rank and problem size. More... | |
EngngModel * | giveEngngModel () |
Returns reference to itself -> required by communicator.h. More... | |
virtual bool | isElementActivated (int elemNum) |
virtual bool | isElementActivated (Element *e) |
virtual void | drawYourself (oofegGraphicContext &gc) |
virtual void | drawElements (oofegGraphicContext &gc) |
virtual void | drawNodes (oofegGraphicContext &gc) |
virtual void | showSparseMtrxStructure (int type, oofegGraphicContext &gc, TimeStep *tStep) |
Shows the sparse structure of required matrix, type == 1 stiffness. More... | |
std::string | errorInfo (const char *func) const |
Returns string for prepending output (used by error reporting macros). More... | |
Protected Member Functions | |
void | updateInternalState (TimeStep *) |
Updates nodal values (calls also this->updateDofUnknownsDictionary for updating dofs unknowns dictionaries if model supports changes of static system). More... | |
void | applyIC (TimeStep *) |
Initializes velocity and pressure fields in the first step with prescribed values. More... | |
void | deactivateTooCloseParticles () |
Deactivates particles upon the particalRemovalRatio. More... | |
Protected Member Functions inherited from oofem::EngngModel | |
virtual void | packMigratingData (TimeStep *tStep) |
Packs receiver data when rebalancing load. More... | |
virtual void | unpackMigratingData (TimeStep *tStep) |
Unpacks receiver data when rebalancing load. More... | |
Protected Attributes | |
SparseLinearSystemNM * | nMethod |
Numerical method used to solve the problem. More... | |
LinSystSolverType | solverType |
Used solver type for linear system of equations. More... | |
SparseMtrxType | sparseMtrxType |
Used type of sparse matrix. More... | |
FloatArray | avLhs |
Left-hand side matrix for the auxiliary velocity equations. More... | |
std::unique_ptr< SparseMtrx > | pLhs |
Left-hand side matrix for the pressure equations. More... | |
FloatArray | vLhs |
Left-hand side matrix for the velocity equations. More... | |
PrimaryField | PressureField |
Pressure field. More... | |
PrimaryField | VelocityField |
Velocity field. More... | |
FloatArray | AuxVelocity |
Array of auxiliary velocities used during computation. More... | |
double | deltaT |
Time step length. More... | |
double | minDeltaT |
Minimal value of time step. More... | |
double | alphaShapeCoef |
Value of alpha coefficient for the boundary recognition. More... | |
double | particleRemovalRatio |
Element side ratio for the removal of the close partices. More... | |
double | rtolv |
Convergence tolerance. More... | |
double | rtolp |
int | maxiter |
Max number of iterations. More... | |
double | domainVolume |
Area or volume of the fluid domain, which can be controlled. More... | |
bool | printVolumeReport |
Flag for volume report. More... | |
int | discretizationScheme |
Explicit or implicit time discretization. More... | |
int | associatedCrossSection |
Number of cross section to associate with created elements. More... | |
int | associatedMaterial |
Number of material to associate with created elements. More... | |
int | associatedPressureBC |
Number of pressure boundary condition to be prescribed on free surface. More... | |
PressureNumberingScheme | pns |
Pressure numbering. More... | |
AuxVelocityNumberingScheme | avns |
Auxiliary Velocity numbering. More... | |
VelocityNumberingScheme | vns |
Velocity numbering. More... | |
VelocityNumberingScheme | prescribedVns |
Prescribed velocity numbering. More... | |
Protected Attributes inherited from oofem::EngngModel | |
int | ndomains |
Number of receiver domains. More... | |
std::vector< std::unique_ptr< Domain > > | domainList |
List of problem domains. More... | |
int | numberOfSteps |
Total number of time steps. More... | |
int | numberOfEquations |
Total number of equation in current time step. More... | |
int | numberOfPrescribedEquations |
Total number or prescribed equations in current time step. More... | |
IntArray | domainNeqs |
Number of equations per domain. More... | |
IntArray | domainPrescribedNeqs |
Number of prescribed equations per domain. More... | |
bool | renumberFlag |
Renumbering flag (renumbers equations after each step, necessary if Dirichlet BCs change). More... | |
bool | profileOpt |
Profile optimized numbering flag (using Sloan's algorithm). More... | |
int | equationNumberingCompleted |
Equation numbering completed flag. More... | |
int | nMetaSteps |
Number of meta steps. More... | |
std::vector< MetaStep > | metaStepList |
List of problem metasteps. More... | |
std::unique_ptr< TimeStep > | stepWhenIcApply |
Solution step when IC (initial conditions) apply. More... | |
std::unique_ptr< TimeStep > | currentStep |
Current time step. More... | |
std::unique_ptr< TimeStep > | previousStep |
Previous time step. More... | |
int | number |
Receivers id. More... | |
std::string | dataOutputFileName |
Path to output stream. More... | |
std::string | coreOutputFileName |
String with core output file name. More... | |
FILE * | outputStream |
Output stream. More... | |
std::string | referenceFileName |
String with reference file name. More... | |
ContextOutputMode | contextOutputMode |
Domain context output mode. More... | |
int | contextOutputStep |
ExportModuleManager * | exportModuleManager |
Export module manager. More... | |
InitModuleManager * | initModuleManager |
Initialization module manager. More... | |
problemMode | pMode |
Domain mode. More... | |
problemScale | pScale |
Multiscale mode. More... | |
time_t | startTime |
Solution start time. More... | |
EngngModel * | master |
Master e-model; if defined receiver is in maintained (slave) mode. More... | |
EngngModelContext * | context |
Context. More... | |
EngngModelTimer | timer |
E-model timer. More... | |
int | parallelFlag |
Flag indicating that the receiver runs in parallel. More... | |
enum fMode | nonLinFormulation |
Type of non linear formulation (total or updated formulation). More... | |
ErrorEstimator * | defaultErrEstimator |
Error estimator. Useful for adaptivity, or simply printing errors output. More... | |
int | rank |
Domain rank in a group of collaborating processes (0..groupSize-1). More... | |
int | numProcs |
Total number of collaborating processes. More... | |
int | nonlocalExt |
Flag indicating if nonlocal extension active, which will cause data to be sent between shared elements before computing the internal forces. More... | |
char | processor_name [PROCESSOR_NAME_LENGTH] |
Processor name. More... | |
MPI_Comm | comm |
Communication object for this engineering model. More... | |
CommunicatorBuff * | commBuff |
Common Communicator buffer. More... | |
ProblemCommunicator * | communicator |
Communicator. More... | |
ProblemCommunicator * | nonlocCommunicator |
NonLocal Communicator. Necessary when nonlocal constitutive models are used. More... | |
std::vector< ParallelContext > | parallelContextList |
List where parallel contexts are stored. More... | |
bool | suppressOutput |
Flag for suppressing output to file. More... | |
std::string | simulationDescription |
LoadBalancer * | lb |
Load Balancer. More... | |
LoadBalancerMonitor * | lbm |
bool | loadBalancingFlag |
If set to true, load balancing is active. More... | |
bool | force_load_rebalance_in_first_step |
Debug flag forcing load balancing after first step. More... | |
Additional Inherited Members | |
Public Types inherited from oofem::EngngModel | |
enum | EngngModel_UpdateMode { EngngModel_Unknown_Mode, EngngModel_SUMM_Mode, EngngModel_SET_Mode } |
enum | EngngModelCommType { PC_default, PC_nonlocal } |
enum | InitialGuess { IG_None = 0, IG_Tangent = 1 } |
Means to choose methods for finding a good initial guess. More... | |
Protected Types inherited from oofem::EngngModel | |
enum | { InternalForcesExchangeTag, MassExchangeTag, LoadExchangeTag, ReactionExchangeTag, RemoteElementExchangeTag } |
Message tags. More... | |
This class represents PFEM method for solving incompressible Navier-Stokes equations.
|
inline |
|
protected |
Initializes velocity and pressure fields in the first step with prescribed values.
Definition at line 800 of file pfem.C.
References oofem::FloatArray::at(), oofem::Element::computeArea(), oofem::Domain::giveDofManagers(), oofem::Domain::giveElements(), OOFEM_LOG_INFO, oofem::FloatArray::resize(), oofem::PFEMElement::updateInternalState(), oofem::Element::updateYourself(), and oofem::FloatArray::zero().
|
virtual |
Allows programmer to test some receiver's internal data, before computation begins.
Reimplemented from oofem::EngngModel.
Definition at line 753 of file pfem.C.
References oofem::EngngModel::checkConsistency(), oofem::Domain::giveElements(), and OOFEM_WARNING.
|
protected |
Deactivates particles upon the particalRemovalRatio.
Definition at line 909 of file pfem.C.
References oofem::PFEMParticle::deactivate(), oofem::FloatArray::distance(), oofem::Node::giveCoordinates(), oofem::Domain::giveElement(), oofem::Element::giveNode(), oofem::Domain::giveNumberOfElements(), oofem::PFEMParticle::isActive(), oofem::max(), oofem::min(), and OOFEM_ERROR.
|
virtual |
Forces equation renumbering on given domain.
All equation numbers in all dofManagers are invalidated, and new equation numbers are generated starting from domainNeqs entry corresponding to given domain. It will update numberOfEquations variable accordingly. Should be used at startup to force equation numbering and therefore sets numberOfEquations. Must be used if model supports changes of static system to assign new valid equation numbers to dofManagers.
Reimplemented from oofem::EngngModel.
Definition at line 178 of file pfem.C.
References oofem::Domain::giveDofManagers().
|
inline |
Returns number of cross section to be associated with created elements.
Definition at line 275 of file pfem.h.
Referenced by oofem::DelaunayTriangulator::writeMesh().
|
inline |
Returns number of material to be associated with created elements.
Definition at line 273 of file pfem.h.
Referenced by oofem::DelaunayTriangulator::writeMesh().
|
inline |
Returns number of zero pressure boundary condition to be prescribed on free surface.
Definition at line 277 of file pfem.h.
Referenced by oofem::DelaunayTriangulator::writeMesh().
|
inlinevirtual |
Returns class name of the receiver.
Implements oofem::EngngModel.
|
inlinevirtual |
Indicates type of non linear computation (total or updated formulation).
This is used for example on Nodal level to update coordinates if updated formulation is done, or on element level, when non linear contributions are computed.
Reimplemented from oofem::EngngModel.
Definition at line 245 of file pfem.h.
References oofem::AL.
|
virtual |
Increases number of equations of receiver's domain and returns newly created equation number.
Used mainly by DofManagers to allocate their corresponding equation number if it is not currently allocated. The DofIDItem parameter allows to distinguish between several possible governing equations, that can be numbered separately.
Reimplemented from oofem::EngngModel.
Definition at line 861 of file pfem.C.
References OOFEM_ERROR.
|
virtual |
Increases number of prescribed equations of receiver's domain and returns newly created equation number.
Used mainly by DofManagers to allocate their corresponding equation number if it is not currently allocated. The DofIDItem parameter allows to distinguish between several possible governing equations, that can be numbered separately.
Reimplemented from oofem::EngngModel.
Definition at line 873 of file pfem.C.
References OOFEM_ERROR.
|
virtual |
Returns next time step (next to current step) of receiver.
Reimplemented from oofem::EngngModel.
Definition at line 343 of file pfem.C.
References oofem::Element::computeArea(), oofem::TR1_2D_PFEM::computeCriticalTimeStep(), oofem::EngngModel::forceEquationNumbering(), oofem::Domain::giveElement(), oofem::Domain::giveNumberOfElements(), oofem::min(), OOFEM_LOG_INFO, and oofem::VST_Time.
|
virtual |
Returns number of equations for given domain in active (current time step) time step.
The numbering scheme determines which system the result is requested for.
Reimplemented from oofem::EngngModel.
Definition at line 884 of file pfem.C.
References oofem::EngngModel::forceEquationNumbering(), and oofem::UnknownNumberingScheme::giveRequiredNumberOfDomainEquation().
|
virtual |
Returns reference to receiver's numerical method.
Reimplemented from oofem::EngngModel.
Definition at line 162 of file pfem.C.
References oofem::classFactory, oofem::ClassFactory::createSparseLinSolver(), and OOFEM_ERROR.
|
virtual |
Returns requested unknown.
Unknown at give time step is characterized by its type and mode and by its equation number. This function is used by Dofs, when they are requested for their associated unknowns.
Reimplemented from oofem::EngngModel.
Definition at line 263 of file pfem.C.
References oofem::__ValueModeTypeToString(), oofem::Dof::giveUnknowns(), and OOFEM_ERROR.
Referenced by oofem::InteractionLoad::computeValueAt(), and oofem::FluidStructureProblem::solveYourselfAt().
|
virtual |
This method is responsible for computing unique dictionary id (ie hash value) from given valueModeType and time step.
This function is used by particular dofs to access unknown identified by given parameters from its dictionary using computed index. Usually the hash algorithm should produce index that depend on time step relatively to actual one to avoid storage of complete history.
Reimplemented from oofem::EngngModel.
Definition at line 629 of file pfem.C.
References oofem::TimeStep::giveNumber(), and OOFEM_ERROR.
|
virtual |
Initialization from given input record.
Reimplemented from oofem::EngngModel.
Definition at line 215 of file pfem.C.
References _IFT_EngngModel_lstype, _IFT_EngngModel_smtype, _IFT_PFEM_alphashapecoef, _IFT_PFEM_associatedCrossSection, _IFT_PFEM_associatedMaterial, _IFT_PFEM_deltat, _IFT_PFEM_discretizationScheme, _IFT_PFEM_maxiter, _IFT_PFEM_mindeltat, _IFT_PFEM_particalRemovalRatio, _IFT_PFEM_pressureBC, _IFT_PFEM_printVolumeReport, _IFT_PFEM_rtolp, _IFT_PFEM_rtolv, oofem::EngngModel::initializeFrom(), IR_GIVE_FIELD, IR_GIVE_OPTIONAL_FIELD, and oofem::IRRT_OK.
|
virtual |
Removes all elements and call DelaunayTriangulator to build up new mesh with new recognized boundary.
Non-meshed particles are set free and move according Newton laws of motion
Reimplemented from oofem::EngngModel.
Definition at line 304 of file pfem.C.
References oofem::Domain::clearElements(), oofem::DelaunayTriangulator::generateMesh(), oofem::Domain::giveDofManagers(), oofem::Domain::giveElements(), oofem::Domain::giveNumberOfElements(), oofem::PFEMParticle::setFree(), and VERBOSE_PRINTS.
DOF printing routine.
Called by DofManagers to print Dof specific part. Dof class provides component printing routines, but emodel is responsible for what will be printed at DOF level.
stream | output stream |
iDof | dof to be processed |
atTime | solution step |
Reimplemented from oofem::EngngModel.
Definition at line 770 of file pfem.C.
References oofem::DofManager::giveCoordinate(), oofem::Dof::giveDofID(), oofem::Dof::giveDofManager(), OOFEM_ERROR, and oofem::Dof::printSingleOutputAt().
|
inlinevirtual |
Indicates if EngngModel requires Dofs dictionaries to be updated.
If EngngModel does not support changes of static system, the dof forwards the requests for its unknowns to EngngModel, where unknowns are naturally kept. This is possible, because dof equation number is same during whole solution. But when changes of static system are allowed, several problem arise. For example by solving simple incremental static with allowed static changes, the incremental displacement vector of structure can not be added to total displacement vector of structure, because equation numbers may have changed, and one can not simply add these vector to obtain new total displacement vector, because incompatible displacement will be added. To solve this problem, unknown dictionary at dof level has been assumed. Dof then keeps its unknowns in its own private dictionary. After computing increment of solution, engngModel updates for each dof its unknowns in its dictionary (using updateUnknownsDictionary function). For aforementioned example engngModel updates incremental values but also total value by asking dof for previous total value (dof will use its dictionary, does not asks back EngngModel) adds corresponding increment and updates total value in dictionary.
Reimplemented from oofem::EngngModel.
void oofem::PFEM::resetEquationNumberings | ( | ) |
|
virtual |
Restores the state of model from output stream.
Restores not only the receiver state, but also same function is invoked for all DofManagers and Elements in associated domain. Note that by restoring element context also contexts of all associated integration points (and material statuses) are restored. Each context is associated with unique time step. Only one context per time step is allowed. Restore context function will restore such context, which is related (through its step number) to time step number and version given in obj parameter. Restoring context will change current time step in order to correspond to newly restored context.
stream | Context file. |
mode | Determines amount of info in stream. |
ContextIOERR | exception if error encountered. |
Reimplemented from oofem::EngngModel.
Definition at line 732 of file pfem.C.
References oofem::CIO_OK, oofem::EngngModel::restoreContext(), and THROW_CIOERR.
|
virtual |
Stores the state of model to output stream.
Stores not only the receiver state, but also same function is invoked for all DofManagers and Elements in associated domain. Note that by storing element context also contexts of all associated integration points (and material statuses) are stored.
stream | Context stream. |
mode | Determines amount of info in stream. |
ContextIOERR | If error encountered. |
Reimplemented from oofem::EngngModel.
Definition at line 710 of file pfem.C.
References oofem::CIO_OK, oofem::EngngModel::saveContext(), and THROW_CIOERR.
|
virtual |
Solves problem for given time step.
Should assemble characteristic matrices and vectors if necessary and solve problem using appropriate numerical method. After finishing solution, this->updateYourself function for updating solution state and then this->terminate function (for updating nodal and element values) should be called.
assumming the mass matrix is always lumped, PFEMMassVelocityAssembler could be replaced by:
Reimplemented from oofem::EngngModel.
Definition at line 392 of file pfem.C.
References oofem::FloatArray::add(), oofem::FloatArray::at(), oofem::classFactory, oofem::Load::computeComponentArrayAt(), oofem::FloatArray::computeNorm(), oofem::ClassFactory::createSparseMtrx(), oofem::Domain::giveDofManagers(), oofem::Domain::giveLoad(), oofem::TimeStep::giveMetaStepNumber(), oofem::TimeStep::giveNumber(), oofem::TimeStep::givePreviousStep(), oofem::FloatArray::giveSize(), oofem::TimeStep::giveTimeIncrement(), oofem::TimeStep::incrementStateCounter(), oofem::PFEMParticle::isActive(), oofem::PFEMParticle::isFree(), oofem::TimeStep::isTheFirstStep(), oofem::max(), oofem::FloatArray::negated(), OOFEM_ERROR, OOFEM_LOG_INFO, oofem::FloatArray::resize(), oofem::FloatArray::subtract(), oofem::FloatArray::times(), and oofem::FloatArray::zero().
|
virtual |
Updates necessary values in Dofs unknown dictionaries.
Reimplemented from oofem::EngngModel.
Definition at line 642 of file pfem.C.
References oofem::FloatArray::at(), and oofem::FloatArray::giveSize().
void oofem::PFEM::updateDofUnknownsDictionaryPressure | ( | DofManager * | inode, |
TimeStep * | tStep | ||
) |
Writes pressures into the dof unknown dictionaries.
Definition at line 670 of file pfem.C.
References oofem::FloatArray::at(), oofem::Dof::giveBcValue(), oofem::DofManager::giveDofWithID(), oofem::Dof::hasBc(), and oofem::Dof::updateUnknownsDictionary().
void oofem::PFEM::updateDofUnknownsDictionaryVelocities | ( | DofManager * | inode, |
TimeStep * | tStep | ||
) |
Writes velocities into the dof unknown dictionaries.
Definition at line 686 of file pfem.C.
References oofem::FloatArray::at(), and oofem::FloatArray::giveSize().
|
protected |
|
virtual |
Updates nodal values (calls also this->updateDofUnknownsDictionary for updating dofs unknowns dictionaries if model supports changes of static system).
The element internal state update is also forced using updateInternalState service.giv
Reimplemented from oofem::EngngModel.
Definition at line 602 of file pfem.C.
References oofem::EngngModel::updateYourself().
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |