OOFEM 3.0
Loading...
Searching...
No Matches
oofem::PFEM Class Reference

#include <pfem.h>

Inheritance diagram for oofem::PFEM:
Collaboration diagram for oofem::PFEM:

Public Member Functions

 PFEM (int i, EngngModel *_master=nullptr)
 ~PFEM ()
void solveYourselfAt (TimeStep *) override
void updateYourself (TimeStep *tStep) override
double giveUnknownComponent (ValueModeType mode, TimeStep *tStep, Domain *d, Dof *dof) override
void saveContext (DataStream &stream, ContextMode mode) override
void restoreContext (DataStream &stream, ContextMode mode) override
TimeStepgiveNextStep () override
 Returns next time step (next to current step) of receiver.
TimeStepgiveSolutionStepWhenIcApply (bool force=false) override
NumericalMethodgiveNumericalMethod (MetaStep *) override
 Returns reference to receiver's numerical method.
void preInitializeNextStep () override
int forceEquationNumbering (int id) override
int requiresUnknownsDictionaryUpdate () override
void initializeFrom (InputRecord &ir) override
 Initialization from given input record.
int checkConsistency () override
const char * giveClassName () const override
 Returns class name of the receiver.
fMode giveFormulation () override
void printDofOutputAt (FILE *stream, Dof *iDof, TimeStep *atTime) override
int giveNumberOfDomainEquations (int, const UnknownNumberingScheme &num) override
int giveNewEquationNumber (int domain, DofIDItem) override
int giveNewPrescribedEquationNumber (int domain, DofIDItem) override
int giveUnknownDictHashIndx (ValueModeType mode, TimeStep *stepN) override
void updateDofUnknownsDictionary (DofManager *inode, TimeStep *tStep) override
void updateDofUnknownsDictionaryPressure (DofManager *inode, TimeStep *tStep)
 Writes pressures into the dof unknown dictionaries.
void updateDofUnknownsDictionaryVelocities (DofManager *inode, TimeStep *tStep)
 Writes velocities into the dof unknown dictionaries.
void resetEquationNumberings ()
 Resets the equation numberings as the mesh is recreated in each time step.
int giveAssociatedMaterialNumber ()
 Returns number of material to be associated with created elements.
int giveAssociatedCrossSectionNumber ()
 Returns number of cross section to be associated with created elements.
int giveAssociatedPressureBC ()
 Returns number of zero pressure boundary condition to be prescribed on free surface.
 EngngModel (int i, EngngModel *_master=NULL)
virtual ~EngngModel ()
 Destructor.
 EngngModel (const EngngModel &)=delete
EngngModeloperator= (const EngngModel &)=delete
DomaingiveDomain (int n)
void setDomain (int i, Domain *ptr, bool iDeallocateOld=true)
int giveNumberOfDomains ()
 Returns number of domains in problem.
const std::string & giveDescription () const
const time_t & giveStartTime ()
bool giveSuppressOutput () const
virtual ErrorEstimatorgiveDomainErrorEstimator (int n)
virtual MaterialInterfacegiveMaterialInterface (int n)
void setNumberOfEquations (int id, int neq)
FILE * giveOutputStream ()
 Returns file descriptor of output file.
std::string giveOutputBaseFileName ()
std::string giveReferenceFileName ()
void letOutputBaseFileNameBe (const std ::string &src)
ContextOutputMode giveContextOutputMode () const
int giveContextOutputStep () const
void setContextOutputMode (ContextOutputMode contextMode)
void setUDContextOutputMode (int cStep)
double giveDeltaT ()
 Returns time step size from the time step controlelr.
void setDeltaT (double dT)
 Returns time step size through the time step controlelr.
void setProblemMode (problemMode pmode)
void setParallelMode (bool newParallelFlag)
problemMode giveProblemMode () const
 Returns domain mode.
void setProblemScale (problemScale pscale)
problemScale giveProblemScale () const
 Returns scale in multiscale simulation.
virtual void setRenumberFlag ()
 Sets the renumber flag to true.
virtual void resetRenumberFlag ()
 Sets the renumber flag to false.
double giveSolutionStepTime ()
void giveAnalysisTime (int &rhrs, int &rmin, int &rsec, int &uhrs, int &umin, int &usec)
void terminateAnalysis ()
virtual void solveYourself ()
virtual void restartYourself (TimeStep *tS)
virtual void terminate (TimeStep *tStep)
virtual void doStepOutput (TimeStep *tStep)
void saveStepContext (TimeStep *tStep, ContextMode mode)
virtual void initializeYourself (TimeStep *tStep)
virtual int initializeAdaptive (int tStepNumber)
virtual FieldPtr giveField (FieldType key, TimeStep *)
virtual FieldPtr giveField (InternalStateType key, TimeStep *)
EngngModelgiveMasterEngngModel ()
 Returns the master engnmodel.
virtual double giveLoadLevel ()
 Returns the current load level.
virtual double giveEigenValue (int eigNum)
 Only relevant for eigen value analysis. Otherwise returns zero.
virtual void setActiveVector (int i)
 Only relevant for eigen value analysis. Otherwise does noting.
int updateSharedDofManagers (FloatArray &answer, const UnknownNumberingScheme &s, int ExchangeTag)
int exchangeRemoteElementData (int ExchangeTag)
virtual int giveCurrentNumberOfIterations ()
MPI_Comm giveParallelComm ()
 Returns the communication object of reciever.
int packRemoteElementData (ProcessCommunicator &processComm)
int unpackRemoteElementData (ProcessCommunicator &processComm)
int packDofManagers (ArrayWithNumbering *src, ProcessCommunicator &processComm)
int unpackDofManagers (ArrayWithNumbering *dest, ProcessCommunicator &processComm)
ProblemCommunicatorgiveProblemCommunicator (EngngModelCommType t)
void initializeCommMaps (bool forceInit=false)
virtual int instanciateYourself (DataReader &dr, InputRecord &ir, const char *outFileName, const char *desc)
void Instanciate_init ()
int instanciateDomains (DataReader &dr)
 Instanciate problem domains by calling their instanciateYourself() service.
int instanciateMetaSteps (DataReader &dr)
 Instanciate problem meta steps by calling their instanciateYourself() service.
virtual int instanciateDefaultMetaStep (InputRecord &ir)
 Instanciate default metastep, if nmsteps is zero.
virtual void updateAttributes (MetaStep *mStep)
void initMetaStepAttributes (MetaStep *mStep)
virtual void updateDomainLinks ()
MetaStepgiveCurrentMetaStep ()
 Returns current meta step.
virtual TimeStepgiveCurrentStep (bool force=false)
virtual void adaptTimeStep (double nIter)
virtual TimeStepgivePreviousStep (bool force=false)
virtual int giveNumberOfFirstStep (bool force=false)
int giveNumberOfMetaSteps ()
 Return number of meta steps.
MetaStepgiveMetaStep (int i)
 Returns the i-th meta step.
int giveNumberOfSteps (bool force=false)
virtual double giveEndOfTimeOfInterest ()
 Returns end of time interest (time corresponding to end of time integration).
int giveNumberOfTimeStepWhenIcApply ()
 Returns the time step number, when initial conditions should apply.
ExportModuleManagergiveExportModuleManager ()
 Returns receiver's export module manager.
EngngModelTimergiveTimer ()
 Returns reference to receiver timer (EngngModelTimer).
virtual double giveInitialTime ()
 return time at the begining of analysis
virtual double giveFinalTime ()
std::string giveContextFileName (int tStepNumber, int stepVersion) const
std::string giveDomainFileName (int domainNum, int domainSerNum) const
virtual void updateComponent (TimeStep *tStep, NumericalCmpn cmpn, Domain *d)
virtual void initForNewIteration (Domain *d, TimeStep *tStep, int iterationNumber, const FloatArray &solution)
virtual void updateSolution (FloatArray &solutionVector, TimeStep *tStep, Domain *d)
virtual void updateInternalRHS (FloatArray &answer, TimeStep *tStep, Domain *d, FloatArray *eNorm)
virtual void updateMatrix (SparseMtrx &mat, TimeStep *tStep, Domain *d)
virtual void initStepIncrements ()
virtual int forceEquationNumbering ()
virtual bool requiresEquationRenumbering (TimeStep *tStep)
virtual bool newDofHandling ()
virtual ParallelContextgiveParallelContext (int n)
virtual void initParallelContexts ()
virtual void assemble (SparseMtrx &answer, TimeStep *tStep, const MatrixAssembler &ma, const UnknownNumberingScheme &s, Domain *domain)
virtual void assemble (SparseMtrx &answer, TimeStep *tStep, const MatrixAssembler &ma, const UnknownNumberingScheme &r_s, const UnknownNumberingScheme &c_s, Domain *domain)
void assembleVector (FloatArray &answer, TimeStep *tStep, const VectorAssembler &va, ValueModeType mode, const UnknownNumberingScheme &s, Domain *domain, FloatArray *eNorms=NULL)
void assembleVectorFromDofManagers (FloatArray &answer, TimeStep *tStep, const VectorAssembler &va, ValueModeType mode, const UnknownNumberingScheme &s, Domain *domain, FloatArray *eNorms=NULL)
void assembleVectorFromElements (FloatArray &answer, TimeStep *tStep, const VectorAssembler &va, ValueModeType mode, const UnknownNumberingScheme &s, Domain *domain, FloatArray *eNorms=NULL)
void assembleVectorFromBC (FloatArray &answer, TimeStep *tStep, const VectorAssembler &va, ValueModeType mode, const UnknownNumberingScheme &s, Domain *domain, FloatArray *eNorms=NULL)
void assembleExtrapolatedForces (FloatArray &answer, TimeStep *tStep, CharType type, Domain *domain)
void assemblePrescribedExtrapolatedForces (FloatArray &answer, TimeStep *tStep, CharType type, Domain *domain)
virtual int checkProblemConsistency ()
virtual void init ()
virtual void postInitialize ()
virtual void printOutputAt (FILE *file, TimeStep *tStep)
virtual void printOutputAt (FILE *file, TimeStep *tStep, const IntArray &nodeSets, const IntArray &elementSets)
void outputNodes (FILE *file, Domain &domain, TimeStep *tStep, int setNum)
void outputElements (FILE *file, Domain &domain, TimeStep *tStep, int setNum)
void printYourself ()
 Prints state of receiver. Useful for debugging.
virtual int useNonlocalStiffnessOption ()
 Returns nonzero if nonlocal stiffness option activated.
bool isParallel () const
 Returns true if receiver in parallel mode.
int giveRank () const
 Returns domain rank in a group of collaborating processes (0..groupSize-1).
int giveNumberOfProcesses () const
 Returns the number of collaborating processes.
EngngModelContextgiveContext ()
 Context requesting service.
virtual int giveNumberOfSlaveProblems ()
 Returns number of slave problems.
virtual EngngModelgiveSlaveProblem (int i)
 Returns i-th slave problem.
virtual bool giveEquationScalingFlag ()
 Returns the Equation scaling flag, which is used to indicate that governing equation(s) are scaled, or non-dimensionalized.
virtual double giveVariableScale (VarScaleType varId)
 Returns the scale factor for given variable type.
virtual int estimateMaxPackSize (IntArray &commMap, DataStream &buff, int packUnpackType)
virtual void balanceLoad (TimeStep *tStep)
virtual LoadBalancergiveLoadBalancer ()
virtual LoadBalancerMonitorgiveLoadBalancerMonitor ()
void initParallel ()
 Request domain rank and problem size.
EngngModelgiveEngngModel ()
 Returns reference to itself -> required by communicator.h.
virtual bool isElementActivated (int elemNum)
virtual bool isElementActivated (Element *e)
TimeStepControllergiveTimeStepController ()
 Returns the time step controller.
virtual void drawYourself (oofegGraphicContext &gc)
virtual void drawElements (oofegGraphicContext &gc)
virtual void drawNodes (oofegGraphicContext &gc)
virtual void showSparseMtrxStructure (int type, oofegGraphicContext &gc, TimeStep *tStep)
std::string errorInfo (const char *func) const
 Returns string for prepending output (used by error reporting macros).

Protected Member Functions

void updateInternalState (TimeStep *)
void applyIC (TimeStep *)
 Initializes velocity and pressure fields in the first step with prescribed values.
void deactivateTooCloseParticles ()
 Deactivates particles upon the particalRemovalRatio.
virtual void packMigratingData (TimeStep *tStep)
virtual void unpackMigratingData (TimeStep *tStep)

Protected Attributes

std ::unique_ptr< SparseLinearSystemNMnMethod
 Numerical method used to solve the problem.
LinSystSolverType solverType
 Used solver type for linear system of equations.
SparseMtrxType sparseMtrxType
 Used type of sparse matrix.
FloatArray avLhs
 Left-hand side matrix for the auxiliary velocity equations.
std ::unique_ptr< SparseMtrxpLhs
 Left-hand side matrix for the pressure equations.
FloatArray vLhs
 Left-hand side matrix for the velocity equations.
PrimaryField PressureField
 Pressure field.
PrimaryField VelocityField
 Velocity field.
FloatArray AuxVelocity
 Array of auxiliary velocities used during computation.
double deltaT
 Time step length.
double minDeltaT
 Minimal value of time step.
double alphaShapeCoef
 Value of alpha coefficient for the boundary recognition.
double particleRemovalRatio
 Element side ratio for the removal of the close partices.
double rtolv
 Convergence tolerance.
double rtolp
int maxiter
 Max number of iterations.
double domainVolume
 Area or volume of the fluid domain, which can be controlled.
bool printVolumeReport
 Flag for volume report.
int discretizationScheme
 Explicit or implicit time discretization.
int associatedCrossSection
 Number of cross section to associate with created elements.
int associatedMaterial
 Number of material to associate with created elements.
int associatedPressureBC
 Number of pressure boundary condition to be prescribed on free surface.
PressureNumberingScheme pns
 Pressure numbering.
AuxVelocityNumberingScheme avns
 Auxiliary Velocity numbering.
VelocityNumberingScheme vns
 Velocity numbering.
VelocityNumberingScheme prescribedVns
 Prescribed velocity numbering.
Protected Attributes inherited from oofem::EngngModel
int ndomains
 Number of receiver domains.
std ::vector< std ::unique_ptr< Domain > > domainList
 List of problem domains.
int numberOfSteps
 Total number of time steps.
int numberOfEquations
 Total number of equation in current time step.
int numberOfPrescribedEquations
 Total number or prescribed equations in current time step.
IntArray domainNeqs
 Number of equations per domain.
IntArray domainPrescribedNeqs
 Number of prescribed equations per domain.
bool renumberFlag
 Renumbering flag (renumbers equations after each step, necessary if Dirichlet BCs change).
bool profileOpt
 Profile optimized numbering flag (using Sloan's algorithm).
int equationNumberingCompleted
 Equation numbering completed flag.
int nMetaSteps
 Number of meta steps.
std ::vector< MetaStepmetaStepList
 List of problem metasteps.
std ::unique_ptr< TimeStepstepWhenIcApply
 Solution step when IC (initial conditions) apply.
std ::unique_ptr< TimeStepcurrentStep
 Current time step.
std ::unique_ptr< TimeSteppreviousStep
 Previous time step.
int number
 Receivers id.
std::string dataOutputFileName
 Path to output stream.
std::string coreOutputFileName
 String with core output file name.
FILE * outputStream
 Output stream.
std::string referenceFileName
 String with reference file name.
ContextOutputMode contextOutputMode
 Domain context output mode.
int contextOutputStep
ExportModuleManager exportModuleManager
 Export module manager.
InitModuleManager initModuleManager
 Initialization module manager.
MonitorManager monitorManager
 Monitor manager.
problemMode pMode
 Domain mode.
problemScale pScale
 Multiscale mode.
time_t startTime
 Solution start time.
EngngModelmaster
 Master e-model; if defined receiver is in maintained (slave) mode.
EngngModelContextcontext
 Context.
EngngModelTimer timer
 E-model timer.
int parallelFlag
 Flag indicating that the receiver runs in parallel.
enum fMode nonLinFormulation
 Type of non linear formulation (total or updated formulation).
std::unique_ptr< ErrorEstimatordefaultErrEstimator
 Error estimator. Useful for adaptivity, or simply printing errors output.
std::unique_ptr< TimeStepControllertimeStepController
 Time Step controller is responsible for collecting data from analysis, elements, and materials, and select the appropriate timestep size for the next step, or reduce the step in case of convergence problems.
int rank
 Domain rank in a group of collaborating processes (0..groupSize-1).
int numProcs
 Total number of collaborating processes.
int nonlocalExt
 Flag indicating if nonlocal extension active, which will cause data to be sent between shared elements before computing the internal forces.
char processor_name [PROCESSOR_NAME_LENGTH]
 Processor name.
MPI_Comm comm
 Communication object for this engineering model.
std::unique_ptr< LoadBalancerlb
 Load Balancer.
std::unique_ptr< LoadBalancerMonitorlbm
bool loadBalancingFlag
 If set to true, load balancing is active.
bool force_load_rebalance_in_first_step
 Debug flag forcing load balancing after first step.
CommunicatorBuffcommBuff
 Common Communicator buffer.
ProblemCommunicatorcommunicator
 Communicator.
ProblemCommunicatornonlocCommunicator
 NonLocal Communicator. Necessary when nonlocal constitutive models are used.
std ::vector< ParallelContextparallelContextList
 List where parallel contexts are stored.
bool suppressOutput
 Flag for suppressing output to file.
std::string simulationDescription

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 }
Protected Types inherited from oofem::EngngModel
enum  {
  InternalForcesExchangeTag , MassExchangeTag , LoadExchangeTag , ReactionExchangeTag ,
  RemoteElementExchangeTag
}
 Message tags. More...

Detailed Description

This class represents PFEM method for solving incompressible Navier-Stokes equations

Author
David Krybus

Definition at line 124 of file pfem.h.

Constructor & Destructor Documentation

◆ PFEM()

◆ ~PFEM()

oofem::PFEM::~PFEM ( )
inline

Definition at line 207 of file pfem.h.

Member Function Documentation

◆ applyIC()

◆ checkConsistency()

int oofem::PFEM::checkConsistency ( )
overridevirtual

Allows programmer to test some receiver's internal data, before computation begins.

Returns
Nonzero if receiver check is o.k.

Reimplemented from oofem::EngngModel.

Definition at line 726 of file pfem.C.

References oofem::EngngModel::giveDomain(), oofem::Domain::giveElements(), and OOFEM_WARNING.

◆ deactivateTooCloseParticles()

◆ forceEquationNumbering()

int oofem::PFEM::forceEquationNumbering ( int i)
overridevirtual

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 176 of file pfem.C.

References avns, oofem::EngngModel::domainNeqs, oofem::EngngModel::domainPrescribedNeqs, oofem::EngngModel::giveCurrentStep(), oofem::Domain::giveDofManagers(), oofem::EngngModel::giveDomain(), pns, and resetEquationNumberings().

◆ giveAssociatedCrossSectionNumber()

int oofem::PFEM::giveAssociatedCrossSectionNumber ( )
inline

Returns number of cross section to be associated with created elements.

Definition at line 275 of file pfem.h.

References associatedCrossSection.

Referenced by oofem::DelaunayTriangulator::writeMesh().

◆ giveAssociatedMaterialNumber()

int oofem::PFEM::giveAssociatedMaterialNumber ( )
inline

Returns number of material to be associated with created elements.

Definition at line 273 of file pfem.h.

References associatedMaterial.

Referenced by oofem::DelaunayTriangulator::writeMesh().

◆ giveAssociatedPressureBC()

int oofem::PFEM::giveAssociatedPressureBC ( )
inline

Returns number of zero pressure boundary condition to be prescribed on free surface.

Definition at line 277 of file pfem.h.

References associatedPressureBC.

Referenced by oofem::DelaunayTriangulator::writeMesh().

◆ giveClassName()

const char * oofem::PFEM::giveClassName ( ) const
inlineoverridevirtual

Returns class name of the receiver.

Implements oofem::EngngModel.

Definition at line 244 of file pfem.h.

◆ giveFormulation()

fMode oofem::PFEM::giveFormulation ( )
inlineoverridevirtual

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.

◆ giveNewEquationNumber()

int oofem::PFEM::giveNewEquationNumber ( int domain,
DofIDItem  )
overridevirtual

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 834 of file pfem.C.

References OOFEM_ERROR, and vns.

◆ giveNewPrescribedEquationNumber()

int oofem::PFEM::giveNewPrescribedEquationNumber ( int domain,
DofIDItem  )
overridevirtual

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 846 of file pfem.C.

References OOFEM_ERROR, and prescribedVns.

◆ giveNextStep()

◆ giveNumberOfDomainEquations()

int oofem::PFEM::giveNumberOfDomainEquations ( int di,
const UnknownNumberingScheme & num )
overridevirtual

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 857 of file pfem.C.

References oofem::EngngModel::equationNumberingCompleted, and oofem::UnknownNumberingScheme::giveRequiredNumberOfDomainEquation().

Referenced by applyIC(), and solveYourselfAt().

◆ giveNumericalMethod()

NumericalMethod * oofem::PFEM::giveNumericalMethod ( MetaStep * mStep)
overridevirtual

Returns reference to receiver's numerical method.

Reimplemented from oofem::EngngModel.

Definition at line 163 of file pfem.C.

References oofem::classFactory, oofem::EngngModel::giveDomain(), nMethod, OOFEM_ERROR, and solverType.

Referenced by solveYourselfAt().

◆ giveSolutionStepWhenIcApply()

TimeStep * oofem::PFEM::giveSolutionStepWhenIcApply ( bool force = false)
overridevirtual

Returns the solution step when Initial Conditions (IC) apply.

Parameters
forcewhen set to true then receiver reply is returned instead of master (default)

Reimplemented from oofem::EngngModel.

Definition at line 288 of file pfem.C.

References deltaT, oofem::EngngModel::giveNumberOfTimeStepWhenIcApply(), oofem::EngngModel::master, and oofem::EngngModel::stepWhenIcApply.

Referenced by giveNextStep().

◆ giveUnknownComponent()

double oofem::PFEM::giveUnknownComponent ( ValueModeType ,
TimeStep * ,
Domain * ,
Dof *  )
overridevirtual

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.

See also
Dof::giveUnknown

Reimplemented from oofem::EngngModel.

Definition at line 257 of file pfem.C.

References AuxVelocity, avns, giveUnknownDictHashIndx(), oofem::Dof::giveUnknowns(), OOFEM_ERROR, and requiresUnknownsDictionaryUpdate().

Referenced by oofem::InteractionLoad::computeValueAt(), and oofem::FluidStructureProblem::solveYourselfAt().

◆ giveUnknownDictHashIndx()

int oofem::PFEM::giveUnknownDictHashIndx ( ValueModeType mode,
TimeStep * tStep )
overridevirtual

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 626 of file pfem.C.

References oofem::EngngModel::giveCurrentStep(), oofem::TimeStep::giveNumber(), oofem::EngngModel::givePreviousStep(), and OOFEM_ERROR.

Referenced by giveUnknownComponent().

◆ initializeFrom()

◆ preInitializeNextStep()

void oofem::PFEM::preInitializeNextStep ( )
overridevirtual

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 301 of file pfem.C.

References alphaShapeCoef, oofem::Domain::clearElements(), oofem::DelaunayTriangulator::generateMesh(), oofem::Domain::giveDofManagers(), oofem::EngngModel::giveDomain(), oofem::Domain::giveElements(), oofem::Domain::giveNumberOfElements(), oofem::PFEMParticle::setFree(), and VERBOSE_PRINTS.

◆ printDofOutputAt()

void oofem::PFEM::printDofOutputAt ( FILE * stream,
Dof * iDof,
TimeStep * atTime )
overridevirtual

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.

Parameters
streamoutput stream
iDofdof to be processed
atTimesolution step

Reimplemented from oofem::EngngModel.

Definition at line 743 of file pfem.C.

References oofem::DofManager::giveCoordinate(), oofem::Dof::giveDofID(), oofem::Dof::giveDofManager(), OOFEM_ERROR, and oofem::Dof::printSingleOutputAt().

◆ requiresUnknownsDictionaryUpdate()

int oofem::PFEM::requiresUnknownsDictionaryUpdate ( )
inlineoverridevirtual

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.

Definition at line 235 of file pfem.h.

Referenced by giveUnknownComponent(), and updateInternalState().

◆ resetEquationNumberings()

void oofem::PFEM::resetEquationNumberings ( )

Resets the equation numberings as the mesh is recreated in each time step.

Definition at line 873 of file pfem.C.

References avns, pns, prescribedVns, and vns.

Referenced by forceEquationNumbering().

◆ restoreContext()

void oofem::PFEM::restoreContext ( DataStream & stream,
ContextMode mode )
overridevirtual

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.

Parameters
streamContext file.
modeDetermines amount of info in stream.
Exceptions
ContextIOERRexception if error encountered.

Reimplemented from oofem::EngngModel.

Definition at line 717 of file pfem.C.

References PressureField, and VelocityField.

◆ saveContext()

void oofem::PFEM::saveContext ( DataStream & stream,
ContextMode mode )
overridevirtual

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.

Parameters
streamContext stream.
modeDetermines amount of info in stream.
Exceptions
ContextIOERRIf error encountered.

Reimplemented from oofem::EngngModel.

Definition at line 707 of file pfem.C.

References PressureField, and VelocityField.

◆ solveYourselfAt()

void oofem::PFEM::solveYourselfAt ( TimeStep * tStep)
overridevirtual

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 389 of file pfem.C.

References oofem::FloatArray::add(), applyIC(), oofem::assemble(), oofem::EngngModel::assembleVector(), oofem::EngngModel::assembleVectorFromElements(), oofem::FloatArray::at(), AuxVelocity, avLhs, avns, oofem::classFactory, oofem::Load::computeComponentArrayAt(), oofem::FloatArray::computeNorm(), deltaT, discretizationScheme, oofem::Domain::giveDofManagers(), oofem::EngngModel::giveDomain(), oofem::Domain::giveLoad(), oofem::EngngModel::giveMetaStep(), oofem::TimeStep::giveMetaStepNumber(), oofem::TimeStep::giveNumber(), giveNumberOfDomainEquations(), oofem::EngngModel::giveNumberOfFirstStep(), giveNumericalMethod(), oofem::TimeStep::givePreviousStep(), oofem::FloatArray::giveSize(), oofem::TimeStep::giveTimeIncrement(), oofem::TimeStep::incrementStateCounter(), oofem::PFEMParticle::isActive(), oofem::PFEMParticle::isFree(), oofem::TimeStep::isTheFirstStep(), oofem::max(), maxiter, oofem::FloatArray::negated(), nMethod, OOFEM_ERROR, OOFEM_LOG_INFO, pLhs, pns, PressureField, oofem::FloatArray::resize(), rtolp, rtolv, sparseMtrxType, oofem::EngngModel::stepWhenIcApply, oofem::FloatArray::subtract(), oofem::FloatArray::times(), updateDofUnknownsDictionary(), updateDofUnknownsDictionaryPressure(), updateDofUnknownsDictionaryVelocities(), VelocityField, vLhs, vns, and oofem::FloatArray::zero().

◆ updateDofUnknownsDictionary()

void oofem::PFEM::updateDofUnknownsDictionary ( DofManager * ,
TimeStep *  )
overridevirtual

Updates necessary values in Dofs unknown dictionaries.

See also
EngngModel::requiresUnknownsDictionaryUpdate
Dof::updateUnknownsDictionary

Reimplemented from oofem::EngngModel.

Definition at line 639 of file pfem.C.

References oofem::FloatArray::at(), oofem::FloatArray::giveSize(), pns, PressureField, and VelocityField.

Referenced by solveYourselfAt(), and updateInternalState().

◆ updateDofUnknownsDictionaryPressure()

void oofem::PFEM::updateDofUnknownsDictionaryPressure ( DofManager * inode,
TimeStep * tStep )

Writes pressures into the dof unknown dictionaries.

Definition at line 667 of file pfem.C.

References oofem::FloatArray::at(), oofem::Dof::giveBcValue(), oofem::DofManager::giveDofWithID(), oofem::Dof::hasBc(), pns, PressureField, and oofem::Dof::updateUnknownsDictionary().

Referenced by solveYourselfAt().

◆ updateDofUnknownsDictionaryVelocities()

void oofem::PFEM::updateDofUnknownsDictionaryVelocities ( DofManager * inode,
TimeStep * tStep )

Writes velocities into the dof unknown dictionaries.

Definition at line 683 of file pfem.C.

References oofem::FloatArray::at(), oofem::FloatArray::giveSize(), and VelocityField.

Referenced by solveYourselfAt().

◆ updateInternalState()

void oofem::PFEM::updateInternalState ( TimeStep * stepN)
protected

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.

Definition at line 610 of file pfem.C.

References oofem::EngngModel::domainList, requiresUnknownsDictionaryUpdate(), and updateDofUnknownsDictionary().

Referenced by updateYourself().

◆ updateYourself()

void oofem::PFEM::updateYourself ( TimeStep * tStep)
overridevirtual

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 599 of file pfem.C.

References deactivateTooCloseParticles(), and updateInternalState().

Member Data Documentation

◆ alphaShapeCoef

double oofem::PFEM::alphaShapeCoef
protected

Value of alpha coefficient for the boundary recognition.

Definition at line 154 of file pfem.h.

Referenced by initializeFrom(), and preInitializeNextStep().

◆ associatedCrossSection

int oofem::PFEM::associatedCrossSection
protected

Number of cross section to associate with created elements.

Definition at line 171 of file pfem.h.

Referenced by giveAssociatedCrossSectionNumber(), initializeFrom(), and PFEM().

◆ associatedMaterial

int oofem::PFEM::associatedMaterial
protected

Number of material to associate with created elements.

Definition at line 173 of file pfem.h.

Referenced by giveAssociatedMaterialNumber(), initializeFrom(), and PFEM().

◆ associatedPressureBC

int oofem::PFEM::associatedPressureBC
protected

Number of pressure boundary condition to be prescribed on free surface.

Definition at line 175 of file pfem.h.

Referenced by giveAssociatedPressureBC(), initializeFrom(), and PFEM().

◆ AuxVelocity

FloatArray oofem::PFEM::AuxVelocity
protected

Array of auxiliary velocities used during computation.

Definition at line 147 of file pfem.h.

Referenced by giveUnknownComponent(), and solveYourselfAt().

◆ avLhs

FloatArray oofem::PFEM::avLhs
protected

Left-hand side matrix for the auxiliary velocity equations.

Definition at line 135 of file pfem.h.

Referenced by PFEM(), and solveYourselfAt().

◆ avns

AuxVelocityNumberingScheme oofem::PFEM::avns
protected

Auxiliary Velocity numbering.

Definition at line 180 of file pfem.h.

Referenced by forceEquationNumbering(), giveUnknownComponent(), resetEquationNumberings(), and solveYourselfAt().

◆ deltaT

double oofem::PFEM::deltaT
protected

Time step length.

Definition at line 150 of file pfem.h.

Referenced by giveNextStep(), giveSolutionStepWhenIcApply(), initializeFrom(), and solveYourselfAt().

◆ discretizationScheme

int oofem::PFEM::discretizationScheme
protected

Explicit or implicit time discretization.

Definition at line 168 of file pfem.h.

Referenced by initializeFrom(), PFEM(), and solveYourselfAt().

◆ domainVolume

double oofem::PFEM::domainVolume
protected

Area or volume of the fluid domain, which can be controlled.

Definition at line 163 of file pfem.h.

Referenced by applyIC(), giveNextStep(), and PFEM().

◆ maxiter

int oofem::PFEM::maxiter
protected

Max number of iterations.

Definition at line 160 of file pfem.h.

Referenced by initializeFrom(), and solveYourselfAt().

◆ minDeltaT

double oofem::PFEM::minDeltaT
protected

Minimal value of time step.

Definition at line 152 of file pfem.h.

Referenced by initializeFrom().

◆ nMethod

std :: unique_ptr< SparseLinearSystemNM > oofem::PFEM::nMethod
protected

Numerical method used to solve the problem.

Definition at line 128 of file pfem.h.

Referenced by giveNumericalMethod(), PFEM(), and solveYourselfAt().

◆ particleRemovalRatio

double oofem::PFEM::particleRemovalRatio
protected

Element side ratio for the removal of the close partices.

Definition at line 156 of file pfem.h.

Referenced by deactivateTooCloseParticles(), and initializeFrom().

◆ pLhs

std :: unique_ptr< SparseMtrx > oofem::PFEM::pLhs
protected

Left-hand side matrix for the pressure equations.

Definition at line 137 of file pfem.h.

Referenced by PFEM(), and solveYourselfAt().

◆ pns

◆ prescribedVns

VelocityNumberingScheme oofem::PFEM::prescribedVns
protected

Prescribed velocity numbering.

Definition at line 184 of file pfem.h.

Referenced by giveNewPrescribedEquationNumber(), PFEM(), and resetEquationNumberings().

◆ PressureField

PrimaryField oofem::PFEM::PressureField
protected

◆ printVolumeReport

bool oofem::PFEM::printVolumeReport
protected

Flag for volume report.

Definition at line 165 of file pfem.h.

Referenced by giveNextStep(), initializeFrom(), and PFEM().

◆ rtolp

double oofem::PFEM::rtolp
protected

Definition at line 158 of file pfem.h.

Referenced by initializeFrom(), and solveYourselfAt().

◆ rtolv

double oofem::PFEM::rtolv
protected

Convergence tolerance.

Definition at line 158 of file pfem.h.

Referenced by initializeFrom(), and solveYourselfAt().

◆ solverType

LinSystSolverType oofem::PFEM::solverType
protected

Used solver type for linear system of equations.

Definition at line 130 of file pfem.h.

Referenced by giveNumericalMethod(), and initializeFrom().

◆ sparseMtrxType

SparseMtrxType oofem::PFEM::sparseMtrxType
protected

Used type of sparse matrix.

Definition at line 132 of file pfem.h.

Referenced by initializeFrom(), and solveYourselfAt().

◆ VelocityField

PrimaryField oofem::PFEM::VelocityField
protected

◆ vLhs

FloatArray oofem::PFEM::vLhs
protected

Left-hand side matrix for the velocity equations.

Definition at line 139 of file pfem.h.

Referenced by PFEM(), and solveYourselfAt().

◆ vns

VelocityNumberingScheme oofem::PFEM::vns
protected

Velocity numbering.

Definition at line 182 of file pfem.h.

Referenced by applyIC(), giveNewEquationNumber(), PFEM(), resetEquationNumberings(), and solveYourselfAt().


The documentation for this class was generated from the following files:

This page is part of the OOFEM-3.0 documentation. Copyright Copyright (C) 1994-2025 Borek Patzak Bořek Patzák
Project e-mail: oofem@fsv.cvut.cz
Generated at for OOFEM by doxygen 1.15.0 written by Dimitri van Heesch, © 1997-2011