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

#include <cbs.h>

Inheritance diagram for oofem::CBS:
Collaboration diagram for oofem::CBS:

Public Member Functions

 CBS (int i, EngngModel *_master=NULL)
virtual ~CBS ()
void solveYourselfAt (TimeStep *tStep) override
void updateYourself (TimeStep *tStep) override
int giveUnknownDictHashIndx (ValueModeType mode, TimeStep *tStep) override
double giveUnknownComponent (ValueModeType type, TimeStep *tStep, Domain *d, Dof *dof) override
bool newDofHandling () override
double giveReynoldsNumber () override
double giveTheta1 ()
double giveTheta2 ()
double giveTractionPressure (Dof *dof)
void saveContext (DataStream &stream, ContextMode mode) override
void restoreContext (DataStream &stream, ContextMode mode) override
void updateDomainLinks () override
TimeStepgiveNextStep () override
 Returns next time step (next to current step) of receiver.
TimeStepgiveSolutionStepWhenIcApply (bool force=false) override
NumericalMethodgiveNumericalMethod (MetaStep *mStep) override
 Returns reference to receiver's numerical method.
void initializeFrom (InputRecord &ir) override
int checkConsistency () override
const char * giveClassName () const override
 Returns class name of the receiver.
const char * giveInputRecordName () const
fMode giveFormulation () override
void printDofOutputAt (FILE *stream, Dof *iDof, TimeStep *tStep) override
int giveNumberOfDomainEquations (int, const UnknownNumberingScheme &num) override
int giveNewEquationNumber (int domain, DofIDItem) override
int giveNewPrescribedEquationNumber (int domain, DofIDItem) override
bool giveEquationScalingFlag () override
 Returns the Equation scaling flag, which is used to indicate that governing equation(s) are scaled, or non-dimensionalized.
double giveVariableScale (VarScaleType varId) override
 Returns the scale factor for given variable type.
Public Member Functions inherited from oofem::FluidModel
 FluidModel (int i, EngngModel *master)
int forceEquationNumbering (int id) override
int forceEquationNumbering () override
 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)
MetaStepgiveCurrentMetaStep ()
 Returns current meta step.
virtual TimeStepgiveCurrentStep (bool force=false)
virtual void adaptTimeStep (double nIter)
virtual TimeStepgivePreviousStep (bool force=false)
virtual void preInitializeNextStep ()
 Does a pre-initialization of the next time step (implement if necessarry).
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 requiresUnknownsDictionaryUpdate ()
virtual bool requiresEquationRenumbering (TimeStep *tStep)
virtual void updateDofUnknownsDictionary (DofManager *, TimeStep *)
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 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 *tStep)
void applyIC (TimeStep *tStep)
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
SparseMtrxType sparseMtrxType
std ::unique_ptr< SparseMtrxlhs
DofDistributedPrimaryField pressureField
 Pressure field.
DofDistributedPrimaryField velocityField
 Velocity field.
FloatArray deltaAuxVelocity
FloatArray prescribedTractionPressure
FloatArray nodalPrescribedTractionPressureConnectivity
FloatArray mm
 Lumped mass matrix.
std ::unique_ptr< SparseMtrxmss
 Sparse consistent mass.
double deltaT
 Time step and its minimal value.
double minDeltaT
double theta1
 Integration constants.
double theta2
int initFlag
int consistentMassFlag
 Consistent mass flag.
VelocityEquationNumbering vnum
VelocityEquationNumbering vnumPrescribed
PressureEquationNumbering pnum
PressureEquationNumbering pnumPrescribed
bool equationScalingFlag
double lscale
 Length scale.
double uscale
 Velocity scale.
double dscale
 Density scale.
double Re
 Reynolds number.
std ::unique_ptr< MaterialInterfacematerialInterface
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 CBS algorithm for solving incompressible Navier-Stokes equations

Definition at line 171 of file cbs.h.

Constructor & Destructor Documentation

◆ CBS()

◆ ~CBS()

oofem::CBS::~CBS ( )
virtual

Definition at line 119 of file cbs.C.

Member Function Documentation

◆ applyIC()

void oofem::CBS::applyIC ( TimeStep * tStep)
protected

◆ checkConsistency()

int oofem::CBS::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 499 of file cbs.C.

References equationScalingFlag, oofem::ForceLoadBVT, oofem::Domain::giveBcs(), oofem::EngngModel::giveDomain(), oofem::Domain::giveElements(), oofem::Domain::giveIcs(), giveVariableScale(), OOFEM_WARNING, oofem::PressureBVT, uscale, oofem::VelocityBVT, oofem::VST_Force, and oofem::VST_Pressure.

◆ giveClassName()

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

Returns class name of the receiver.

Implements oofem::EngngModel.

Definition at line 251 of file cbs.h.

◆ giveEquationScalingFlag()

bool oofem::CBS::giveEquationScalingFlag ( )
inlineoverridevirtual

Returns the Equation scaling flag, which is used to indicate that governing equation(s) are scaled, or non-dimensionalized.

Reimplemented from oofem::EngngModel.

Definition at line 262 of file cbs.h.

References equationScalingFlag.

◆ giveFormulation()

fMode oofem::CBS::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 253 of file cbs.h.

References oofem::TL.

◆ giveInputRecordName()

const char * oofem::CBS::giveInputRecordName ( ) const
inline

Definition at line 252 of file cbs.h.

References _IFT_CBS_Name.

◆ giveNewEquationNumber()

int oofem::CBS::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 590 of file cbs.C.

References OOFEM_ERROR, pnum, and vnum.

◆ giveNewPrescribedEquationNumber()

int oofem::CBS::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 605 of file cbs.C.

References OOFEM_ERROR, pnumPrescribed, and vnumPrescribed.

◆ giveNextStep()

◆ giveNumberOfDomainEquations()

int oofem::CBS::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 619 of file cbs.C.

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

Referenced by solveYourselfAt().

◆ giveNumericalMethod()

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

Returns reference to receiver's numerical method.

Reimplemented from oofem::EngngModel.

Definition at line 123 of file cbs.C.

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

Referenced by solveYourselfAt(), and updateDomainLinks().

◆ giveReynoldsNumber()

double oofem::CBS::giveReynoldsNumber ( )
overridevirtual

Implements oofem::FluidModel.

Definition at line 198 of file cbs.C.

References equationScalingFlag, and Re.

◆ giveSolutionStepWhenIcApply()

TimeStep * oofem::CBS::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 224 of file cbs.C.

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

Referenced by giveNextStep().

◆ giveTheta1()

double oofem::CBS::giveTheta1 ( )

◆ giveTheta2()

double oofem::CBS::giveTheta2 ( )

Definition at line 206 of file cbs.C.

References theta2.

◆ giveTractionPressure()

double oofem::CBS::giveTractionPressure ( Dof * dof)
Todo
This should just disappear completely.

Definition at line 209 of file cbs.C.

References oofem::Dof::__givePrescribedEquationNumber(), OOFEM_ERROR, and prescribedTractionPressure.

Referenced by oofem::TractionPressureBC::give().

◆ giveUnknownComponent()

double oofem::CBS::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 187 of file cbs.C.

References oofem::Dof::giveDofID(), pressureField, and velocityField.

◆ giveUnknownDictHashIndx()

int oofem::CBS::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 443 of file cbs.C.

References oofem::TimeStep::giveNumber().

◆ giveVariableScale()

double oofem::CBS::giveVariableScale ( VarScaleType varId)
overridevirtual

Returns the scale factor for given variable type.

Reimplemented from oofem::EngngModel.

Definition at line 629 of file cbs.C.

References dscale, lscale, OOFEM_ERROR, uscale, oofem::VST_Density, oofem::VST_Force, oofem::VST_Length, oofem::VST_Pressure, oofem::VST_Time, oofem::VST_Velocity, and oofem::VST_Viscosity.

Referenced by checkConsistency(), and giveNextStep().

◆ initializeFrom()

void oofem::CBS::initializeFrom ( InputRecord & ir)
overridevirtual

Initializes receiver according to object description in input reader. InitString 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.

Reimplemented from oofem::EngngModel.

Definition at line 135 of file cbs.C.

References _IFT_CBS_cmflag, _IFT_CBS_deltat, _IFT_CBS_dscale, _IFT_CBS_lscale, _IFT_CBS_miflag, _IFT_CBS_mindeltat, _IFT_CBS_scaleflag, _IFT_CBS_theta1, _IFT_CBS_theta2, _IFT_CBS_uscale, _IFT_EngngModel_lstype, _IFT_EngngModel_smtype, consistentMassFlag, deltaT, dscale, equationScalingFlag, oofem::EngngModel::giveContext(), oofem::EngngModel::giveDomain(), IR_GIVE_FIELD, IR_GIVE_OPTIONAL_FIELD, lscale, materialInterface, minDeltaT, Re, oofem::FieldManager::registerField(), solverType, sparseMtrxType, theta1, theta2, uscale, and velocityField.

◆ newDofHandling()

bool oofem::CBS::newDofHandling ( )
inlineoverridevirtual

Temporary method for allowing code to seamlessly convert from the old to new way of handling DOF values. (the new way expects the field to store all values, regardless of if they are computed, from BC, or IC.) This is used by MasterDof

Todo
When all models have converted to using a field, this should be removed.

Reimplemented from oofem::EngngModel.

Definition at line 232 of file cbs.h.

◆ printDofOutputAt()

void oofem::CBS::printDofOutputAt ( FILE * stream,
Dof * iDof,
TimeStep * tStep )
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
tStepsolution step

Reimplemented from oofem::EngngModel.

Definition at line 556 of file cbs.C.

References dscale, oofem::Dof::giveDofID(), OOFEM_ERROR, oofem::Dof::printSingleOutputAt(), and uscale.

◆ restoreContext()

void oofem::CBS::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 487 of file cbs.C.

References oofem::CIO_OK, prescribedTractionPressure, and THROW_CIOERR.

◆ saveContext()

void oofem::CBS::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 475 of file cbs.C.

References oofem::CIO_OK, prescribedTractionPressure, and THROW_CIOERR.

◆ solveYourselfAt()

◆ updateDomainLinks()

void oofem::CBS::updateDomainLinks ( )
overridevirtual

Updates domain links after the domains of receiver have changed. Used mainly after restoring context - the domains may change and this service is then used to update domain variables in all components belonging to receiver like error estimators, solvers, etc, having domains as attributes.

Reimplemented from oofem::EngngModel.

Definition at line 548 of file cbs.C.

References oofem::EngngModel::giveCurrentMetaStep(), oofem::EngngModel::giveDomain(), and giveNumericalMethod().

◆ updateInternalState()

void oofem::CBS::updateInternalState ( TimeStep * tStep)
protected

Updates the IP values for the new solution velocities

Definition at line 464 of file cbs.C.

References oofem::EngngModel::domainList.

Referenced by solveYourselfAt().

◆ updateYourself()

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

Updates internal state after finishing time step. (for example total values may be updated according to previously solved increments). Then element values are also updated (together with related integration points and material statuses).

Reimplemented from oofem::EngngModel.

Definition at line 450 of file cbs.C.

References materialInterface.

Member Data Documentation

◆ consistentMassFlag

int oofem::CBS::consistentMassFlag
protected

Consistent mass flag.

Definition at line 201 of file cbs.h.

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

◆ deltaAuxVelocity

FloatArray oofem::CBS::deltaAuxVelocity
protected

Definition at line 185 of file cbs.h.

Referenced by solveYourselfAt().

◆ deltaT

double oofem::CBS::deltaT
protected

Time step and its minimal value.

Definition at line 195 of file cbs.h.

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

◆ dscale

double oofem::CBS::dscale
protected

Density scale.

Definition at line 214 of file cbs.h.

Referenced by CBS(), giveVariableScale(), initializeFrom(), and printDofOutputAt().

◆ equationScalingFlag

bool oofem::CBS::equationScalingFlag
protected

◆ initFlag

int oofem::CBS::initFlag
protected

Definition at line 199 of file cbs.h.

Referenced by CBS(), and solveYourselfAt().

◆ lhs

std :: unique_ptr< SparseMtrx > oofem::CBS::lhs
protected

Definition at line 180 of file cbs.h.

Referenced by solveYourselfAt().

◆ lscale

double oofem::CBS::lscale
protected

Length scale.

Definition at line 210 of file cbs.h.

Referenced by CBS(), giveVariableScale(), and initializeFrom().

◆ materialInterface

std :: unique_ptr< MaterialInterface > oofem::CBS::materialInterface
protected

Definition at line 220 of file cbs.h.

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

◆ minDeltaT

double oofem::CBS::minDeltaT
protected

Definition at line 195 of file cbs.h.

Referenced by giveNextStep(), and initializeFrom().

◆ mm

FloatArray oofem::CBS::mm
protected

Lumped mass matrix.

Definition at line 190 of file cbs.h.

Referenced by solveYourselfAt().

◆ mss

std :: unique_ptr< SparseMtrx > oofem::CBS::mss
protected

Sparse consistent mass.

Definition at line 192 of file cbs.h.

Referenced by solveYourselfAt().

◆ nMethod

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

Numerical method used to solve the problem.

Definition at line 175 of file cbs.h.

Referenced by giveNumericalMethod(), and solveYourselfAt().

◆ nodalPrescribedTractionPressureConnectivity

FloatArray oofem::CBS::nodalPrescribedTractionPressureConnectivity
protected

Definition at line 187 of file cbs.h.

Referenced by solveYourselfAt().

◆ pnum

PressureEquationNumbering oofem::CBS::pnum
protected

Definition at line 205 of file cbs.h.

Referenced by CBS(), giveNewEquationNumber(), and solveYourselfAt().

◆ pnumPrescribed

PressureEquationNumbering oofem::CBS::pnumPrescribed
protected

Definition at line 206 of file cbs.h.

Referenced by CBS(), giveNewPrescribedEquationNumber(), and solveYourselfAt().

◆ prescribedTractionPressure

FloatArray oofem::CBS::prescribedTractionPressure
protected

Definition at line 186 of file cbs.h.

Referenced by giveTractionPressure(), restoreContext(), saveContext(), and solveYourselfAt().

◆ pressureField

DofDistributedPrimaryField oofem::CBS::pressureField
protected

Pressure field.

Definition at line 182 of file cbs.h.

Referenced by applyIC(), CBS(), giveUnknownComponent(), and solveYourselfAt().

◆ Re

double oofem::CBS::Re
protected

Reynolds number.

Definition at line 216 of file cbs.h.

Referenced by CBS(), giveReynoldsNumber(), and initializeFrom().

◆ solverType

LinSystSolverType oofem::CBS::solverType
protected

Definition at line 177 of file cbs.h.

Referenced by giveNumericalMethod(), and initializeFrom().

◆ sparseMtrxType

SparseMtrxType oofem::CBS::sparseMtrxType
protected

Definition at line 178 of file cbs.h.

Referenced by initializeFrom(), and solveYourselfAt().

◆ theta1

double oofem::CBS::theta1
protected

Integration constants.

Definition at line 197 of file cbs.h.

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

◆ theta2

double oofem::CBS::theta2
protected

Definition at line 197 of file cbs.h.

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

◆ uscale

double oofem::CBS::uscale
protected

Velocity scale.

Definition at line 212 of file cbs.h.

Referenced by CBS(), checkConsistency(), giveVariableScale(), initializeFrom(), and printDofOutputAt().

◆ velocityField

DofDistributedPrimaryField oofem::CBS::velocityField
protected

Velocity field.

Definition at line 184 of file cbs.h.

Referenced by applyIC(), CBS(), giveUnknownComponent(), initializeFrom(), and solveYourselfAt().

◆ vnum

VelocityEquationNumbering oofem::CBS::vnum
protected

Definition at line 203 of file cbs.h.

Referenced by CBS(), giveNewEquationNumber(), and solveYourselfAt().

◆ vnumPrescribed

VelocityEquationNumbering oofem::CBS::vnumPrescribed
protected

Definition at line 204 of file cbs.h.

Referenced by CBS(), and giveNewPrescribedEquationNumber().


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