|
OOFEM 3.0
|
#include <stokesflow.h>
Public Member Functions | |
| StokesFlow (int i, EngngModel *_master=nullptr) | |
| ~StokesFlow () | |
| void | solveYourselfAt (TimeStep *tStep) override |
| void | updateYourself (TimeStep *tStep) override |
| double | giveUnknownComponent (ValueModeType mode, TimeStep *tStep, Domain *domain, Dof *dof) override |
| bool | newDofHandling () override |
| double | giveReynoldsNumber () override |
| int | forceEquationNumbering (int id) override |
| void | initializeFrom (InputRecord &ir) override |
| int | checkConsistency () override |
| void | doStepOutput (TimeStep *tStep) override |
| void | updateInternalState (TimeStep *tStep) |
| void | updateComponent (TimeStep *tStep, NumericalCmpn cmpn, Domain *d) override |
| void | updateSolution (FloatArray &solutionVector, TimeStep *tStep, Domain *d) override |
| void | updateInternalRHS (FloatArray &answer, TimeStep *tStep, Domain *d, FloatArray *eNorm) override |
| void | updateMatrix (SparseMtrx &mat, TimeStep *tStep, Domain *d) override |
| NumericalMethod * | giveNumericalMethod (MetaStep *mStep) override |
| Returns reference to receiver's numerical method. | |
| TimeStep * | giveNextStep () override |
| Returns next time step (next to current step) of receiver. | |
| const char * | giveClassName () const override |
| Returns class name of the receiver. | |
| const char * | giveInputRecordName () const |
| Public Member Functions inherited from oofem::FluidModel | |
| FluidModel (int i, EngngModel *master) | |
| int | forceEquationNumbering () override |
| EngngModel (int i, EngngModel *_master=NULL) | |
| virtual | ~EngngModel () |
| Destructor. | |
| EngngModel (const EngngModel &)=delete | |
| EngngModel & | operator= (const EngngModel &)=delete |
| Domain * | giveDomain (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 ErrorEstimator * | giveDomainErrorEstimator (int n) |
| virtual MaterialInterface * | giveMaterialInterface (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) |
| void | saveStepContext (TimeStep *tStep, ContextMode mode) |
| virtual void | initializeYourself (TimeStep *tStep) |
| virtual int | initializeAdaptive (int tStepNumber) |
| virtual int | giveNumberOfDomainEquations (int di, const UnknownNumberingScheme &num) |
| virtual FieldPtr | giveField (FieldType key, TimeStep *) |
| virtual FieldPtr | giveField (InternalStateType key, TimeStep *) |
| EngngModel * | giveMasterEngngModel () |
| 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) |
| ProblemCommunicator * | giveProblemCommunicator (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 | saveContext (DataStream &stream, ContextMode mode) |
| virtual void | restoreContext (DataStream &stream, ContextMode mode) |
| virtual void | updateDomainLinks () |
| MetaStep * | giveCurrentMetaStep () |
| Returns current meta step. | |
| virtual TimeStep * | giveCurrentStep (bool force=false) |
| virtual void | adaptTimeStep (double nIter) |
| virtual TimeStep * | givePreviousStep (bool force=false) |
| virtual void | preInitializeNextStep () |
| Does a pre-initialization of the next time step (implement if necessarry). | |
| virtual TimeStep * | giveSolutionStepWhenIcApply (bool force=false) |
| virtual int | giveNumberOfFirstStep (bool force=false) |
| int | giveNumberOfMetaSteps () |
| Return number of meta steps. | |
| MetaStep * | giveMetaStep (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. | |
| ExportModuleManager * | giveExportModuleManager () |
| Returns receiver's export module manager. | |
| EngngModelTimer * | giveTimer () |
| Returns reference to receiver timer (EngngModelTimer). | |
| virtual double | giveInitialTime () |
| return time at the begining of analysis | |
| virtual double | giveFinalTime () |
| virtual int | giveNewEquationNumber (int domain, DofIDItem) |
| virtual int | giveNewPrescribedEquationNumber (int domain, DofIDItem) |
| std::string | giveContextFileName (int tStepNumber, int stepVersion) const |
| std::string | giveDomainFileName (int domainNum, int domainSerNum) const |
| virtual void | initForNewIteration (Domain *d, TimeStep *tStep, int iterationNumber, const FloatArray &solution) |
| virtual void | initStepIncrements () |
| virtual int | requiresUnknownsDictionaryUpdate () |
| virtual bool | requiresEquationRenumbering (TimeStep *tStep) |
| virtual void | updateDofUnknownsDictionary (DofManager *, TimeStep *) |
| virtual int | giveUnknownDictHashIndx (ValueModeType mode, TimeStep *tStep) |
| virtual ParallelContext * | giveParallelContext (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 void | printDofOutputAt (FILE *stream, Dof *iDof, TimeStep *tStep) |
| 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. | |
| virtual fMode | giveFormulation () |
| EngngModelContext * | giveContext () |
| Context requesting service. | |
| virtual int | giveNumberOfSlaveProblems () |
| Returns number of slave problems. | |
| virtual EngngModel * | giveSlaveProblem (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 LoadBalancer * | giveLoadBalancer () |
| virtual LoadBalancerMonitor * | giveLoadBalancerMonitor () |
| void | initParallel () |
| Request domain rank and problem size. | |
| EngngModel * | giveEngngModel () |
| Returns reference to itself -> required by communicator.h. | |
| virtual bool | isElementActivated (int elemNum) |
| virtual bool | isElementActivated (Element *e) |
| TimeStepController * | giveTimeStepController () |
| 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 Attributes | |
| double | deltaT = 1. |
| Time increment read from input record. | |
| std ::unique_ptr< PrimaryField > | velocityPressureField |
| Primary unknowns. | |
| SparseMtrxType | sparseMtrxType |
| Sparse matrix type. | |
| std ::unique_ptr< SparseNonLinearSystemNM > | nMethod |
| Numerical method. | |
| LinSystSolverType | solverType = ST_Direct |
| Linear solver type. | |
| FloatArray | eNorm |
| Element norm for nonlinear analysis (squared). | |
| std ::unique_ptr< MeshQualityErrorEstimator > | meshqualityee |
| Used for determining if a new mesh must be created. | |
| double | maxdef = 1. |
| Maximum deformation allowed. | |
| std ::unique_ptr< SparseMtrx > | stiffnessMatrix |
| FloatArray | solutionVector |
| FloatArray | internalForces |
| TopologyState | ts |
| Topology state, most notably used for determining if there is a need to remesh. | |
| 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< MetaStep > | metaStepList |
| List of problem metasteps. | |
| std ::unique_ptr< TimeStep > | stepWhenIcApply |
| Solution step when IC (initial conditions) apply. | |
| std ::unique_ptr< TimeStep > | currentStep |
| Current time step. | |
| std ::unique_ptr< TimeStep > | previousStep |
| 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. | |
| EngngModel * | master |
| Master e-model; if defined receiver is in maintained (slave) mode. | |
| EngngModelContext * | context |
| 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< ErrorEstimator > | defaultErrEstimator |
| Error estimator. Useful for adaptivity, or simply printing errors output. | |
| std::unique_ptr< TimeStepController > | timeStepController |
| 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< LoadBalancer > | lb |
| Load Balancer. | |
| std::unique_ptr< LoadBalancerMonitor > | lbm |
| 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. | |
| CommunicatorBuff * | commBuff |
| Common Communicator buffer. | |
| ProblemCommunicator * | communicator |
| Communicator. | |
| ProblemCommunicator * | nonlocCommunicator |
| NonLocal Communicator. Necessary when nonlocal constitutive models are used. | |
| std ::vector< ParallelContext > | parallelContextList |
| 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... | |
| virtual void | packMigratingData (TimeStep *tStep) |
| virtual void | unpackMigratingData (TimeStep *tStep) |
Implements the engineering model to solve incompressible Stokes flow. Stokes flow means acceleration is ignored.
Definition at line 63 of file stokesflow.h.
| oofem::StokesFlow::StokesFlow | ( | int | i, |
| EngngModel * | _master = nullptr ) |
Definition at line 53 of file stokesflow.C.
References oofem::FluidModel::FluidModel(), and oofem::EngngModel::ndomains.
Referenced by oofem::StokesFlowVelocityHomogenization::StokesFlowVelocityHomogenization().
| oofem::StokesFlow::~StokesFlow | ( | ) |
Definition at line 59 of file stokesflow.C.
|
overridevirtual |
Allows programmer to test some receiver's internal data, before computation begins.
Reimplemented from oofem::EngngModel.
Definition at line 257 of file stokesflow.C.
References oofem::EngngModel::giveDomain(), oofem::Domain::giveElements(), and OOFEM_WARNING.
|
overridevirtual |
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.
Reimplemented from oofem::EngngModel.
Definition at line 286 of file stokesflow.C.
References oofem::TopologyDescription::doOutput(), and oofem::EngngModel::giveDomain().
|
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::FluidModel.
Definition at line 236 of file stokesflow.C.
References oofem::EngngModel::equationNumberingCompleted, and stiffnessMatrix.
|
inlineoverridevirtual |
Returns class name of the receiver.
Implements oofem::EngngModel.
Reimplemented in oofem::StokesFlowVelocityHomogenization.
Definition at line 123 of file stokesflow.h.
|
inline |
Definition at line 124 of file stokesflow.h.
References _IFT_StokesFlow_Name.
|
overridevirtual |
Returns next time step (next to current step) of receiver.
Reimplemented from oofem::EngngModel.
Definition at line 304 of file stokesflow.C.
References oofem::EngngModel::currentStep, deltaT, oofem::EngngModel::giveNumberOfTimeStepWhenIcApply(), and oofem::EngngModel::previousStep.
|
overridevirtual |
Returns reference to receiver's numerical method.
Reimplemented from oofem::EngngModel.
Definition at line 296 of file stokesflow.C.
References oofem::EngngModel::giveDomain(), and nMethod.
Referenced by solveYourselfAt().
|
overridevirtual |
Implements oofem::FluidModel.
Definition at line 251 of file stokesflow.C.
|
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.
Reimplemented from oofem::EngngModel.
Definition at line 245 of file stokesflow.C.
References velocityPressureField.
|
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 62 of file stokesflow.C.
References _IFT_EngngModel_lstype, _IFT_EngngModel_smtype, _IFT_StokesFlow_deltat, deltaT, IR_GIVE_OPTIONAL_FIELD, maxdef, meshqualityee, oofem::SMT_PetscMtrx, solverType, sparseMtrxType, oofem::ST_Petsc, stiffnessMatrix, ts, oofem::TS_OK, and velocityPressureField.
|
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
Reimplemented from oofem::EngngModel.
Definition at line 105 of file stokesflow.h.
|
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.
Reimplemented from oofem::EngngModel.
Definition at line 88 of file stokesflow.C.
References oofem::EngngModel::assembleVector(), oofem::classFactory, oofem::TimeStep::convergedReason, oofem::CR_CONVERGED, eNorm, oofem::EngngModel::giveCurrentMetaStep(), oofem::EngngModel::giveDomain(), oofem::EngngModel::giveExportModuleManager(), oofem::TimeStep::giveMetaStepNumber(), oofem::TimeStep::giveNumber(), oofem::EngngModel::giveNumberOfDomainEquations(), oofem::Domain::giveNumberOfElements(), giveNumericalMethod(), oofem::EngngModel::giveProblemScale(), oofem::Domain::giveTopology(), oofem::globalErrorEEV, oofem::EngngModel::initMetaStepAttributes(), internalForces, oofem::EngngModel::LoadExchangeTag, oofem::macroScale, maxdef, meshqualityee, nMethod, oofem::TimeStep::numberOfIterations, OOFEM_ERROR, OOFEM_LOG_INFO, oofem::TopologyDescription::replaceFEMesh(), oofem::FloatArray::resize(), solutionVector, sparseMtrxType, stiffnessMatrix, ts, oofem::TS_NeedsRemeshing, updateInternalRHS(), updateMatrix(), oofem::EngngModel::updateSharedDofManagers(), updateSolution(), velocityPressureField, and oofem::FloatArray::zero().
|
overridevirtual |
Reimplemented from oofem::EngngModel.
Definition at line 177 of file stokesflow.C.
References oofem::assemble(), oofem::EngngModel::assembleVector(), eNorm, oofem::Domain::giveElements(), internalForces, oofem::EngngModel::InternalForcesExchangeTag, oofem::InternalRhs, oofem::NonLinearLhs, OOFEM_ERROR, solutionVector, stiffnessMatrix, oofem::EngngModel::updateSharedDofManagers(), oofem::FMElement::updateStabilizationCoeffs(), and velocityPressureField.
|
overridevirtual |
Updates the solution (guess) according to the new values. Callback for nonlinear solvers (e.g. Newton-Raphson).
| solutionVector | New solution. |
| tStep | Time when component is updated. |
| d | Domain. |
| eNorm | Optional per-element norm (for normalization). |
Reimplemented from oofem::EngngModel.
Definition at line 212 of file stokesflow.C.
References oofem::EngngModel::assembleVector(), oofem::EngngModel::InternalForcesExchangeTag, oofem::EngngModel::updateSharedDofManagers(), and oofem::FloatArray::zero().
Referenced by solveYourselfAt().
| void oofem::StokesFlow::updateInternalState | ( | TimeStep * | tStep | ) |
Definition at line 272 of file stokesflow.C.
References oofem::EngngModel::domainList, and ts.
Referenced by updateYourself().
|
overridevirtual |
Updates the solution (guess) according to the new values. Callback for nonlinear solvers (e.g. Newton-Raphson).
| solutionVector | New solution. |
| tStep | Time when component is updated. |
| d | Domain. |
Reimplemented from oofem::EngngModel.
Definition at line 221 of file stokesflow.C.
References oofem::assemble(), and oofem::SparseMtrx::zero().
Referenced by solveYourselfAt().
|
overridevirtual |
Updates the solution (guess) according to the new values. Callback for nonlinear solvers (e.g. Newton-Raphson), and are called before new internal forces are computed.
| solutionVector | New solution. |
| tStep | Time when component is updated. |
| d | Domain. |
Reimplemented from oofem::EngngModel.
Definition at line 202 of file stokesflow.C.
References oofem::Domain::giveElements(), solutionVector, oofem::FMElement::updateStabilizationCoeffs(), and velocityPressureField.
Referenced by solveYourselfAt().
|
overridevirtual |
Updates everything for the problem. Updates the internal state for the elements. Also calls updateNodalPositions.
Reimplemented from oofem::EngngModel.
Definition at line 229 of file stokesflow.C.
References updateInternalState().
|
protected |
Time increment read from input record.
Definition at line 67 of file stokesflow.h.
Referenced by giveNextStep(), and initializeFrom().
|
protected |
Element norm for nonlinear analysis (squared).
Definition at line 77 of file stokesflow.h.
Referenced by solveYourselfAt(), and updateComponent().
|
protected |
Definition at line 86 of file stokesflow.h.
Referenced by solveYourselfAt(), and updateComponent().
|
protected |
Maximum deformation allowed.
Definition at line 82 of file stokesflow.h.
Referenced by initializeFrom(), and solveYourselfAt().
|
protected |
Used for determining if a new mesh must be created.
Definition at line 80 of file stokesflow.h.
Referenced by initializeFrom(), and solveYourselfAt().
|
protected |
Numerical method.
Definition at line 73 of file stokesflow.h.
Referenced by giveNumericalMethod(), and solveYourselfAt().
|
protected |
Definition at line 85 of file stokesflow.h.
Referenced by solveYourselfAt(), updateComponent(), and updateSolution().
|
protected |
Linear solver type.
Definition at line 75 of file stokesflow.h.
Referenced by oofem::StokesFlowVelocityHomogenization::computeTangent(), and initializeFrom().
|
protected |
Sparse matrix type.
Definition at line 71 of file stokesflow.h.
Referenced by initializeFrom(), and solveYourselfAt().
|
protected |
Definition at line 84 of file stokesflow.h.
Referenced by forceEquationNumbering(), initializeFrom(), solveYourselfAt(), and updateComponent().
|
protected |
Topology state, most notably used for determining if there is a need to remesh.
Definition at line 89 of file stokesflow.h.
Referenced by initializeFrom(), solveYourselfAt(), and updateInternalState().
|
protected |
Primary unknowns.
Definition at line 69 of file stokesflow.h.
Referenced by giveUnknownComponent(), initializeFrom(), solveYourselfAt(), updateComponent(), and updateSolution().