OOFEM
2.4
OOFEM.org - Object Oriented Finite Element Solver
|
Solves general nonlinear transient transport problems. More...
#include <transienttransportproblem.h>
Public Member Functions | |
TransientTransportProblem (int i, EngngModel *_master) | |
Constructor. More... | |
virtual | ~TransientTransportProblem () |
Destructor. More... | |
virtual void | solveYourselfAt (TimeStep *tStep) |
Solves problem for given time step. More... | |
virtual void | updateComponent (TimeStep *tStep, NumericalCmpn cmpn, Domain *d) |
Updates components mapped to numerical method if necessary during solution process. More... | |
virtual double | giveUnknownComponent (ValueModeType mode, TimeStep *tStep, Domain *d, Dof *dof) |
Returns requested unknown. More... | |
virtual contextIOResultType | saveContext (DataStream &stream, ContextMode mode) |
Stores the state of model to output stream. More... | |
virtual contextIOResultType | restoreContext (DataStream &stream, ContextMode mode) |
Restores the state of model from output stream. More... | |
virtual void | applyIC () |
virtual int | requiresUnknownsDictionaryUpdate () |
Indicates if EngngModel requires Dofs dictionaries to be updated. More... | |
virtual int | giveUnknownDictHashIndx (ValueModeType mode, TimeStep *tStep) |
This method is responsible for computing unique dictionary id (ie hash value) from given valueModeType and time step. More... | |
virtual void | updateDomainLinks () |
Updates domain links after the domains of receiver have changed. More... | |
Function * | giveDtFunction () |
double | giveDeltaT (int n) |
double | giveDiscreteTime (int iStep) |
virtual TimeStep * | giveNextStep () |
Returns next time step (next to current step) of receiver. More... | |
virtual TimeStep * | giveSolutionStepWhenIcApply (bool force=false) |
Returns the solution step when Initial Conditions (IC) apply. More... | |
virtual NumericalMethod * | giveNumericalMethod (MetaStep *mStep) |
Returns reference to receiver's numerical method. More... | |
virtual IRResultType | initializeFrom (InputRecord *ir) |
Initializes receiver according to object description in input reader. More... | |
virtual bool | requiresEquationRenumbering (TimeStep *tStep) |
Returns true if equation renumbering is required for given solution step. More... | |
virtual int | forceEquationNumbering () |
Forces equation renumbering on all domains associated to engng model. More... | |
virtual void | updateYourself (TimeStep *tStep) |
Updates internal state after finishing time step. More... | |
virtual int | checkConsistency () |
Allows programmer to test some receiver's internal data, before computation begins. More... | |
virtual FieldPtr | giveField (FieldType key, TimeStep *) |
Returns the smart pointer to requested field, Null otherwise. More... | |
virtual const char * | giveInputRecordName () const |
virtual const char * | giveClassName () const |
Returns class name of the receiver. More... | |
virtual fMode | giveFormulation () |
Indicates type of non linear computation (total or updated formulation). More... | |
Public Member Functions inherited from oofem::EngngModel | |
EngngModel (int i, EngngModel *_master=NULL) | |
Constructor. More... | |
virtual | ~EngngModel () |
Destructor. More... | |
EngngModel (const EngngModel &)=delete | |
EngngModel & | operator= (const EngngModel &)=delete |
Domain * | giveDomain (int n) |
Service for accessing particular problem domain. More... | |
void | setDomain (int i, Domain *ptr, bool iDeallocateOld=true) |
Sets i-th domain of receiver. More... | |
int | giveNumberOfDomains () |
Returns number of domains in problem. More... | |
const std::string & | giveDescription () const |
const time_t & | giveStartTime () |
bool | giveSuppressOutput () const |
virtual ErrorEstimator * | giveDomainErrorEstimator (int n) |
Service for accessing ErrorEstimator corresponding to particular domain. More... | |
virtual MaterialInterface * | giveMaterialInterface (int n) |
Returns material interface representation for given domain. More... | |
void | setNumberOfEquations (int id, int neq) |
FILE * | giveOutputStream () |
Returns file descriptor of output file. More... | |
std::string | giveOutputBaseFileName () |
Returns base output file name to which extensions, like .out .vtu .osf should be added. More... | |
std::string | giveReferenceFileName () |
Returns reference file name. More... | |
void | letOutputBaseFileNameBe (const std::string &src) |
Sets the base output file name. More... | |
ContextOutputMode | giveContextOutputMode () |
Returns domain context output mode. More... | |
int | giveContextOutputStep () |
Returns domain context output step. More... | |
void | setContextOutputMode (ContextOutputMode contextMode) |
Sets context output mode of receiver. More... | |
void | setUDContextOutputMode (int cStep) |
Sets user defined context output mode (it sets contextOutputMode to contextOutputMode), setting contextOutputStep to given value. More... | |
void | setProblemMode (problemMode pmode) |
Sets domain mode to given mode. More... | |
void | setParallelMode (bool newParallelFlag) |
Sets the problem to run in parallel (or not). More... | |
problemMode | giveProblemMode () |
Returns domain mode. More... | |
void | setProblemScale (problemScale pscale) |
Sets scale in multiscale simulation. More... | |
problemScale | giveProblemScale () |
Returns scale in multiscale simulation. More... | |
virtual void | setRenumberFlag () |
Sets the renumber flag to true. More... | |
virtual void | resetRenumberFlag () |
Sets the renumber flag to false. More... | |
double | giveSolutionStepTime () |
Returns the user time of the current simulation step in seconds. More... | |
void | giveAnalysisTime (int &rhrs, int &rmin, int &rsec, int &uhrs, int &umin, int &usec) |
Returns the real and user time for the analysis. More... | |
void | terminateAnalysis () |
Performs analysis termination after finishing analysis. More... | |
virtual void | solveYourself () |
Starts solution process. More... | |
virtual void | terminate (TimeStep *tStep) |
Terminates the solution of time step. More... | |
virtual void | doStepOutput (TimeStep *tStep) |
Prints the ouput of the solution step (using virtual this->printOutputAtservice) to the stream detemined using this->giveOutputStream() method and calls exportModuleManager to do output. More... | |
void | saveStepContext (TimeStep *tStep, ContextMode mode) |
Saves context of given solution step, if required (determined using this->giveContextOutputMode() method). More... | |
virtual void | initializeYourself (TimeStep *tStep) |
Provides the opportunity to initialize state variables stored in element integration points according to initial conditions using function initializeYourself() on element level. More... | |
virtual int | initializeAdaptive (int tStepNumber) |
Initializes the newly generated discretization state according to previous solution. More... | |
virtual int | giveNumberOfDomainEquations (int di, const UnknownNumberingScheme &num) |
Returns number of equations for given domain in active (current time step) time step. More... | |
EngngModel * | giveMasterEngngModel () |
Returns the master engnmodel. More... | |
virtual double | giveLoadLevel () |
Returns the current load level. More... | |
virtual double | giveEigenValue (int eigNum) |
Only relevant for eigen value analysis. Otherwise returns zero. More... | |
virtual void | setActiveVector (int i) |
Only relevant for eigen value analysis. Otherwise does noting. More... | |
int | updateSharedDofManagers (FloatArray &answer, const UnknownNumberingScheme &s, int ExchangeTag) |
Exchanges necessary remote DofManagers data. More... | |
int | exchangeRemoteElementData (int ExchangeTag) |
Exchanges necessary remote element data with remote partitions. More... | |
virtual int | giveCurrentNumberOfIterations () |
Returns number of iterations that was required to reach equilibrium - used for adaptive step length in staggered problem. More... | |
MPI_Comm | giveParallelComm () |
Returns the communication object of reciever. More... | |
int | packRemoteElementData (ProcessCommunicator &processComm) |
Packs data of local element to be received by their remote counterpart on remote partitions. More... | |
int | unpackRemoteElementData (ProcessCommunicator &processComm) |
Unpacks data for remote elements (which are mirrors of remote partition's local elements). More... | |
int | packDofManagers (ArrayWithNumbering *src, ProcessCommunicator &processComm) |
Packing function for vector values of DofManagers. More... | |
int | unpackDofManagers (ArrayWithNumbering *dest, ProcessCommunicator &processComm) |
Unpacking function for vector values of DofManagers . More... | |
ProblemCommunicator * | giveProblemCommunicator (EngngModelCommType t) |
void | initializeCommMaps (bool forceInit=false) |
virtual int | instanciateYourself (DataReader &dr, InputRecord *ir, const char *outFileName, const char *desc) |
Initializes whole problem according to its description stored in inputStream. More... | |
void | Instanciate_init () |
Initialization of the receiver state (opening the default output stream, empty domain creation, initialization of parallel context, etc) before Initialization form DataReader. More... | |
int | instanciateDomains (DataReader &dr) |
Instanciate problem domains by calling their instanciateYourself() service. More... | |
int | instanciateMetaSteps (DataReader &dr) |
Instanciate problem meta steps by calling their instanciateYourself() service. More... | |
virtual int | instanciateDefaultMetaStep (InputRecord *ir) |
Instanciate default metastep, if nmsteps is zero. More... | |
virtual void | updateAttributes (MetaStep *mStep) |
Update receiver attributes according to step metaStep attributes. More... | |
void | initMetaStepAttributes (MetaStep *mStep) |
Update e-model attributes attributes according to step metaStep attributes. More... | |
MetaStep * | giveCurrentMetaStep () |
Returns current meta step. More... | |
virtual TimeStep * | giveCurrentStep (bool force=false) |
Returns current time step. More... | |
virtual TimeStep * | givePreviousStep (bool force=false) |
Returns previous time step. More... | |
TimeStep * | generateNextStep () |
Generate new time step (and associate metastep). More... | |
virtual void | preInitializeNextStep () |
Does a pre-initialization of the next time step (implement if necessarry) More... | |
virtual int | giveNumberOfFirstStep (bool force=false) |
Returns number of first time step used by receiver. More... | |
int | giveNumberOfMetaSteps () |
Return number of meta steps. More... | |
MetaStep * | giveMetaStep (int i) |
Returns the i-th meta step. More... | |
int | giveNumberOfSteps (bool force=false) |
Returns total number of steps. More... | |
virtual double | giveEndOfTimeOfInterest () |
Returns end of time interest (time corresponding to end of time integration). More... | |
int | giveNumberOfTimeStepWhenIcApply () |
Returns the time step number, when initial conditions should apply. More... | |
ExportModuleManager * | giveExportModuleManager () |
Returns receiver's export module manager. More... | |
EngngModelTimer * | giveTimer () |
Returns reference to receiver timer (EngngModelTimer). More... | |
virtual int | giveNewEquationNumber (int domain, DofIDItem) |
Increases number of equations of receiver's domain and returns newly created equation number. More... | |
virtual int | giveNewPrescribedEquationNumber (int domain, DofIDItem) |
Increases number of prescribed equations of receiver's domain and returns newly created equation number. More... | |
std::string | giveContextFileName (int tStepNumber, int stepVersion) const |
Returns the filename for the context file for the given step and version. More... | |
std::string | giveDomainFileName (int domainNum, int domainSerNum) const |
Returns the filename for the given domain (used by adaptivity and restore) More... | |
virtual void | initStepIncrements () |
Initializes solution of new time step. More... | |
virtual int | forceEquationNumbering (int i) |
Forces equation renumbering on given domain. More... | |
virtual void | updateDofUnknownsDictionary (DofManager *, TimeStep *) |
Updates necessary values in Dofs unknown dictionaries. More... | |
virtual ParallelContext * | giveParallelContext (int n) |
Returns the parallel context corresponding to given domain (n) and unknown type Default implementation returns i-th context from parallelContextList. More... | |
virtual void | initParallelContexts () |
Creates parallel contexts. More... | |
virtual void | assemble (SparseMtrx &answer, TimeStep *tStep, const MatrixAssembler &ma, const UnknownNumberingScheme &s, Domain *domain) |
Assembles characteristic matrix of required type into given sparse matrix. More... | |
virtual void | assemble (SparseMtrx &answer, TimeStep *tStep, const MatrixAssembler &ma, const UnknownNumberingScheme &r_s, const UnknownNumberingScheme &c_s, Domain *domain) |
Assembles characteristic matrix of required type into given sparse matrix. More... | |
void | assembleVector (FloatArray &answer, TimeStep *tStep, const VectorAssembler &va, ValueModeType mode, const UnknownNumberingScheme &s, Domain *domain, FloatArray *eNorms=NULL) |
Assembles characteristic vector of required type from dofManagers, element, and active boundary conditions, into given vector. More... | |
void | assembleVectorFromDofManagers (FloatArray &answer, TimeStep *tStep, const VectorAssembler &va, ValueModeType mode, const UnknownNumberingScheme &s, Domain *domain, FloatArray *eNorms=NULL) |
Assembles characteristic vector of required type from dofManagers into given vector. More... | |
void | assembleVectorFromElements (FloatArray &answer, TimeStep *tStep, const VectorAssembler &va, ValueModeType mode, const UnknownNumberingScheme &s, Domain *domain, FloatArray *eNorms=NULL) |
Assembles characteristic vector of required type from elements into given vector. More... | |
void | assembleVectorFromBC (FloatArray &answer, TimeStep *tStep, const VectorAssembler &va, ValueModeType mode, const UnknownNumberingScheme &s, Domain *domain, FloatArray *eNorms=NULL) |
Assembles characteristic vector of required type from boundary conditions. More... | |
void | assembleExtrapolatedForces (FloatArray &answer, TimeStep *tStep, CharType type, Domain *domain) |
Assembles the extrapolated internal forces vector, useful for obtaining a good initial guess in nonlinear analysis with Dirichlet boundary conditions. More... | |
void | assemblePrescribedExtrapolatedForces (FloatArray &answer, TimeStep *tStep, CharType type, Domain *domain) |
void | assembleVectorFromContacts (FloatArray &answer, TimeStep *tStep, CharType type, ValueModeType mode, const UnknownNumberingScheme &s, Domain *domain, FloatArray *eNorms=NULL) |
virtual int | checkProblemConsistency () |
Allows programmer to test problem its internal data, before computation begins. More... | |
virtual void | init () |
Initializes the receiver state. More... | |
virtual void | postInitialize () |
Performs post-initialization for all the problem contents (which is called after initializeFrom). More... | |
virtual void | printOutputAt (FILE *file, TimeStep *tStep) |
Prints output of receiver to output domain stream, for given time step. More... | |
virtual void | printOutputAt (FILE *file, TimeStep *tStep, const IntArray &nodeSets, const IntArray &elementSets) |
void | outputNodes (FILE *file, Domain &domain, TimeStep *tStep, int setNum) |
Outputs all nodes in the given set. More... | |
void | outputElements (FILE *file, Domain &domain, TimeStep *tStep, int setNum) |
Outputs all elements in the given set. More... | |
void | printYourself () |
Prints state of receiver. Useful for debugging. More... | |
virtual void | printDofOutputAt (FILE *stream, Dof *iDof, TimeStep *tStep) |
DOF printing routine. More... | |
virtual int | useNonlocalStiffnessOption () |
Returns nonzero if nonlocal stiffness option activated. More... | |
bool | isParallel () const |
Returns true if receiver in parallel mode. More... | |
int | giveRank () const |
Returns domain rank in a group of collaborating processes (0..groupSize-1) More... | |
int | giveNumberOfProcesses () const |
Returns the number of collaborating processes. More... | |
EngngModelContext * | giveContext () |
Context requesting service. More... | |
virtual int | giveNumberOfSlaveProblems () |
Returns number of slave problems. More... | |
virtual EngngModel * | giveSlaveProblem (int i) |
Returns i-th slave problem. More... | |
virtual bool | giveEquationScalingFlag () |
Returns the Equation scaling flag, which is used to indicate that governing equation(s) are scaled, or non-dimensionalized. More... | |
virtual double | giveVariableScale (VarScaleType varId) |
Returns the scale factor for given variable type. More... | |
virtual int | estimateMaxPackSize (IntArray &commMap, DataStream &buff, int packUnpackType) |
Determines the space necessary for send/receive buffer. More... | |
virtual void | balanceLoad (TimeStep *tStep) |
Recovers the load balance between processors, if needed. More... | |
virtual LoadBalancer * | giveLoadBalancer () |
Returns reference to receiver's load balancer. More... | |
virtual LoadBalancerMonitor * | giveLoadBalancerMonitor () |
Returns reference to receiver's load balancer monitor. More... | |
void | initParallel () |
Request domain rank and problem size. More... | |
EngngModel * | giveEngngModel () |
Returns reference to itself -> required by communicator.h. More... | |
virtual bool | isElementActivated (int elemNum) |
virtual bool | isElementActivated (Element *e) |
virtual void | drawYourself (oofegGraphicContext &gc) |
virtual void | drawElements (oofegGraphicContext &gc) |
virtual void | drawNodes (oofegGraphicContext &gc) |
virtual void | showSparseMtrxStructure (int type, oofegGraphicContext &gc, TimeStep *tStep) |
Shows the sparse structure of required matrix, type == 1 stiffness. More... | |
std::string | errorInfo (const char *func) const |
Returns string for prepending output (used by error reporting macros). More... | |
Protected Attributes | |
SparseMtrxType | sparseMtrxType |
std::unique_ptr< PrimaryField > | field |
std::unique_ptr< SparseMtrx > | effectiveMatrix |
FloatArray | solution |
FloatArray | internalForces |
FloatArray | eNorm |
std::unique_ptr< SparseNonLinearSystemNM > | nMethod |
Numerical method used to solve the problem. More... | |
double | alpha |
int | dtFunction |
FloatArray | prescribedTimes |
double | deltaT |
bool | keepTangent |
bool | hasTangent |
bool | lumped |
IntArray | exportFields |
Protected Attributes inherited from oofem::EngngModel | |
int | ndomains |
Number of receiver domains. More... | |
std::vector< std::unique_ptr< Domain > > | domainList |
List of problem domains. More... | |
int | numberOfSteps |
Total number of time steps. More... | |
int | numberOfEquations |
Total number of equation in current time step. More... | |
int | numberOfPrescribedEquations |
Total number or prescribed equations in current time step. More... | |
IntArray | domainNeqs |
Number of equations per domain. More... | |
IntArray | domainPrescribedNeqs |
Number of prescribed equations per domain. More... | |
bool | renumberFlag |
Renumbering flag (renumbers equations after each step, necessary if Dirichlet BCs change). More... | |
bool | profileOpt |
Profile optimized numbering flag (using Sloan's algorithm). More... | |
int | equationNumberingCompleted |
Equation numbering completed flag. More... | |
int | nMetaSteps |
Number of meta steps. More... | |
std::vector< MetaStep > | metaStepList |
List of problem metasteps. More... | |
std::unique_ptr< TimeStep > | stepWhenIcApply |
Solution step when IC (initial conditions) apply. More... | |
std::unique_ptr< TimeStep > | currentStep |
Current time step. More... | |
std::unique_ptr< TimeStep > | previousStep |
Previous time step. More... | |
int | number |
Receivers id. More... | |
std::string | dataOutputFileName |
Path to output stream. More... | |
std::string | coreOutputFileName |
String with core output file name. More... | |
FILE * | outputStream |
Output stream. More... | |
std::string | referenceFileName |
String with reference file name. More... | |
ContextOutputMode | contextOutputMode |
Domain context output mode. More... | |
int | contextOutputStep |
ExportModuleManager * | exportModuleManager |
Export module manager. More... | |
InitModuleManager * | initModuleManager |
Initialization module manager. More... | |
problemMode | pMode |
Domain mode. More... | |
problemScale | pScale |
Multiscale mode. More... | |
time_t | startTime |
Solution start time. More... | |
EngngModel * | master |
Master e-model; if defined receiver is in maintained (slave) mode. More... | |
EngngModelContext * | context |
Context. More... | |
EngngModelTimer | timer |
E-model timer. More... | |
int | parallelFlag |
Flag indicating that the receiver runs in parallel. More... | |
enum fMode | nonLinFormulation |
Type of non linear formulation (total or updated formulation). More... | |
ErrorEstimator * | defaultErrEstimator |
Error estimator. Useful for adaptivity, or simply printing errors output. More... | |
int | rank |
Domain rank in a group of collaborating processes (0..groupSize-1). More... | |
int | numProcs |
Total number of collaborating processes. More... | |
int | nonlocalExt |
Flag indicating if nonlocal extension active, which will cause data to be sent between shared elements before computing the internal forces. More... | |
char | processor_name [PROCESSOR_NAME_LENGTH] |
Processor name. More... | |
MPI_Comm | comm |
Communication object for this engineering model. More... | |
CommunicatorBuff * | commBuff |
Common Communicator buffer. More... | |
ProblemCommunicator * | communicator |
Communicator. More... | |
ProblemCommunicator * | nonlocCommunicator |
NonLocal Communicator. Necessary when nonlocal constitutive models are used. More... | |
std::vector< ParallelContext > | parallelContextList |
List where parallel contexts are stored. More... | |
bool | suppressOutput |
Flag for suppressing output to file. More... | |
std::string | simulationDescription |
LoadBalancer * | lb |
Load Balancer. More... | |
LoadBalancerMonitor * | lbm |
bool | loadBalancingFlag |
If set to true, load balancing is active. More... | |
bool | force_load_rebalance_in_first_step |
Debug flag forcing load balancing after first step. More... | |
Additional Inherited Members | |
Public Types inherited from oofem::EngngModel | |
enum | EngngModel_UpdateMode { EngngModel_Unknown_Mode, EngngModel_SUMM_Mode, EngngModel_SET_Mode } |
enum | EngngModelCommType { PC_default, PC_nonlocal } |
enum | InitialGuess { IG_None = 0, IG_Tangent = 1 } |
Means to choose methods for finding a good initial guess. More... | |
Protected Types inherited from oofem::EngngModel | |
enum | { InternalForcesExchangeTag, MassExchangeTag, LoadExchangeTag, ReactionExchangeTag, RemoteElementExchangeTag } |
Message tags. More... | |
Protected Member Functions inherited from oofem::EngngModel | |
virtual void | packMigratingData (TimeStep *tStep) |
Packs receiver data when rebalancing load. More... | |
virtual void | unpackMigratingData (TimeStep *tStep) |
Unpacks receiver data when rebalancing load. More... | |
Solves general nonlinear transient transport problems.
Definition at line 65 of file transienttransportproblem.h.
oofem::TransientTransportProblem::TransientTransportProblem | ( | int | i, |
EngngModel * | _master = NULL |
||
) |
Constructor.
Definition at line 56 of file transienttransportproblem.C.
References oofem::EngngModel::ndomains.
|
virtual |
Destructor.
Definition at line 67 of file transienttransportproblem.C.
|
virtual |
Definition at line 339 of file transienttransportproblem.C.
References field, oofem::EngngModel::giveDomain(), oofem::Domain::giveElements(), giveSolutionStepWhenIcApply(), OOFEM_LOG_INFO, oofem::TransportElement::updateInternalState(), and oofem::Element::updateYourself().
Referenced by solveYourselfAt().
|
virtual |
Allows programmer to test some receiver's internal data, before computation begins.
Reimplemented from oofem::EngngModel.
Definition at line 441 of file transienttransportproblem.C.
References oofem::EngngModel::checkConsistency(), oofem::EngngModel::giveDomain(), and OOFEM_WARNING.
|
virtual |
Forces equation renumbering on all domains associated to engng model.
All equation numbers in all domains for all dofManagers are invalidated, and new equation numbers are generated starting from 1 on each 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 381 of file transienttransportproblem.C.
References effectiveMatrix, and oofem::EngngModel::forceEquationNumbering().
|
inlinevirtual |
Returns class name of the receiver.
Implements oofem::EngngModel.
Definition at line 126 of file transienttransportproblem.h.
double oofem::TransientTransportProblem::giveDeltaT | ( | int | n | ) |
Definition at line 149 of file transienttransportproblem.C.
References deltaT, dtFunction, oofem::Function::evaluateAtTime(), giveDiscreteTime(), oofem::EngngModel::giveDomain(), oofem::Domain::giveFunction(), oofem::FloatArray::giveSize(), and prescribedTimes.
Referenced by giveNextStep(), and giveSolutionStepWhenIcApply().
double oofem::TransientTransportProblem::giveDiscreteTime | ( | int | iStep | ) |
Definition at line 161 of file transienttransportproblem.C.
References oofem::FloatArray::at(), oofem::FloatArray::giveSize(), OOFEM_ERROR, and prescribedTimes.
Referenced by giveDeltaT().
Function* oofem::TransientTransportProblem::giveDtFunction | ( | ) |
Returns the smart pointer to requested field, Null otherwise.
The return value uses shared_ptr, as some registered fields may be owned (and maintained) by emodel, while some may be created on demand and thus reliable reference counting mechanism is essential.
Reimplemented from oofem::EngngModel.
Definition at line 462 of file transienttransportproblem.C.
References field, oofem::EngngModel::giveCurrentStep(), and OOFEM_ERROR.
|
inlinevirtual |
Indicates type of non linear computation (total or updated formulation).
This is used for example on Nodal level to update coordinates if updated formulation is done, or on element level, when non linear contributions are computed.
Reimplemented from oofem::EngngModel.
Definition at line 127 of file transienttransportproblem.h.
References oofem::TL.
|
inlinevirtual |
Definition at line 125 of file transienttransportproblem.h.
References _IFT_TransientTransportProblem_Name.
|
virtual |
Returns next time step (next to current step) of receiver.
Reimplemented from oofem::EngngModel.
Definition at line 174 of file transienttransportproblem.C.
References alpha, oofem::EngngModel::currentStep, giveDeltaT(), giveSolutionStepWhenIcApply(), and oofem::EngngModel::previousStep.
|
virtual |
Returns reference to receiver's numerical method.
Reimplemented from oofem::EngngModel.
Definition at line 70 of file transienttransportproblem.C.
References oofem::EngngModel::giveDomain(), and nMethod.
Referenced by solveYourselfAt(), and updateDomainLinks().
|
virtual |
Returns the solution step when Initial Conditions (IC) apply.
force | when set to true then receiver reply is returned instead of master (default) |
Reimplemented from oofem::EngngModel.
Definition at line 189 of file transienttransportproblem.C.
References alpha, giveDeltaT(), oofem::EngngModel::giveNumberOfTimeStepWhenIcApply(), oofem::EngngModel::giveSolutionStepWhenIcApply(), oofem::EngngModel::master, and oofem::EngngModel::stepWhenIcApply.
Referenced by applyIC(), and giveNextStep().
|
virtual |
Returns requested unknown.
Unknown at give time step is characterized by its type and mode and by its equation number. This function is used by Dofs, when they are requested for their associated unknowns.
Reimplemented from oofem::EngngModel.
Definition at line 126 of file transienttransportproblem.C.
References alpha, field, oofem::TimeStep::givePreviousStep(), oofem::TimeStep::giveTimeIncrement(), and OOFEM_ERROR.
|
virtual |
This method is responsible for computing unique dictionary id (ie hash value) from given valueModeType and time step.
This function is used by particular dofs to access unknown identified by given parameters from its dictionary using computed index. Usually the hash algorithm should produce index that depend on time step relatively to actual one to avoid storage of complete history.
Reimplemented from oofem::EngngModel.
Definition at line 428 of file transienttransportproblem.C.
References oofem::TimeStep::giveNumber().
|
virtual |
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 80 of file transienttransportproblem.C.
References _IFT_EngngModel_smtype, _IFT_TransientTransportProblem_alpha, _IFT_TransientTransportProblem_deltaT, _IFT_TransientTransportProblem_dtFunction, _IFT_TransientTransportProblem_exportFields, _IFT_TransientTransportProblem_keepTangent, _IFT_TransientTransportProblem_lumped, _IFT_TransientTransportProblem_prescribedTimes, alpha, oofem::IntArray::at(), oofem::IntArray::clear(), oofem::FloatArray::clear(), deltaT, dtFunction, exportFields, field, oofem::EngngModel::giveContext(), oofem::EngngModelContext::giveFieldManager(), oofem::IntArray::giveSize(), oofem::InputRecord::hasField(), oofem::EngngModel::initializeFrom(), IR_GIVE_FIELD, IR_GIVE_OPTIONAL_FIELD, keepTangent, lumped, prescribedTimes, oofem::FieldManager::registerField(), oofem::SMT_Skyline, and sparseMtrxType.
|
virtual |
Returns true if equation renumbering is required for given solution step.
This may of course change the number of equation and in general there is no guarantee that for a certain dof the same equation will be assigned. So the use of DOF unknowns dictionaries is generally recommended.
Reimplemented from oofem::EngngModel.
Definition at line 358 of file transienttransportproblem.C.
References oofem::Domain::giveBcs(), oofem::EngngModel::giveDomain(), oofem::GeneralBoundaryCondition::giveNumberOfInternalDofManagers(), oofem::TimeStep::givePreviousStep(), oofem::TimeStep::isTheFirstStep(), and oofem::ActiveBoundaryCondition::requiresActiveDofs().
|
virtual |
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 435 of file transienttransportproblem.C.
|
virtual |
Restores the state of model from output stream.
Restores not only the receiver state, but also same function is invoked for all DofManagers and Elements in associated domain. Note that by restoring element context also contexts of all associated integration points (and material statuses) are restored. Each context is associated with unique time step. Only one context per time step is allowed. Restore context function will restore such context, which is related (through its step number) to time step number and version given in obj parameter. Restoring context will change current time step in order to correspond to newly restored context.
stream | Context file. |
mode | Determines amount of info in stream. |
ContextIOERR | exception if error encountered. |
Reimplemented from oofem::EngngModel.
Definition at line 411 of file transienttransportproblem.C.
References oofem::CIO_OK, field, oofem::EngngModel::restoreContext(), and THROW_CIOERR.
|
virtual |
Stores the state of model to output stream.
Stores not only the receiver state, but also same function is invoked for all DofManagers and Elements in associated domain. Note that by storing element context also contexts of all associated integration points (and material statuses) are stored.
stream | Context stream. |
mode | Determines amount of info in stream. |
ContextIOERR | If error encountered. |
Reimplemented from oofem::EngngModel.
Definition at line 394 of file transienttransportproblem.C.
References oofem::CIO_OK, field, oofem::EngngModel::saveContext(), and THROW_CIOERR.
|
virtual |
Solves problem for given time step.
Should assemble characteristic matrices and vectors if necessary and solve problem using appropriate numerical method. After finishing solution, this->updateYourself function for updating solution state and then this->terminate function (for updating nodal and element values) should be called.
Reimplemented from oofem::EngngModel.
Definition at line 206 of file transienttransportproblem.C.
References applyIC(), oofem::EngngModel::assembleVector(), oofem::classFactory, oofem::ClassFactory::createSparseMtrx(), effectiveMatrix, eNorm, field, oofem::EngngModel::giveCurrentMetaStep(), oofem::EngngModel::giveDomain(), oofem::TimeStep::giveNumber(), oofem::EngngModel::giveNumberOfDomainEquations(), giveNumericalMethod(), oofem::TimeStep::givePreviousStep(), oofem::TimeStep::giveTargetTime(), internalForces, oofem::InternalRhs, oofem::TimeStep::isTheFirstStep(), oofem::EngngModel::LoadExchangeTag, nMethod, OOFEM_LOG_INFO, oofem::FloatArray::resize(), oofem::SparseNonLinearSystemNM::rlm_total, solution, sparseMtrxType, updateComponent(), oofem::EngngModel::updateSharedDofManagers(), and oofem::FloatArray::zero().
|
virtual |
Updates components mapped to numerical method if necessary during solution process.
Some numerical methods may require updating mapped components during solution process (e.g., updating of tangent stiffness when using updated Newton-Raphson method).
tStep | Time when component is updated. |
cmpn | Numerical component to update. |
d | Domain. |
Reimplemented from oofem::EngngModel.
Definition at line 283 of file transienttransportproblem.C.
References oofem::FloatArray::add(), alpha, oofem::EngngModel::assemble(), oofem::EngngModel::assembleVector(), oofem::FloatArray::beDifferenceOf(), effectiveMatrix, eNorm, field, oofem::TimeStep::givePreviousStep(), oofem::FloatArray::giveSize(), oofem::TimeStep::giveTimeIncrement(), hasTangent, internalForces, oofem::EngngModel::InternalForcesExchangeTag, oofem::InternalRhs, keepTangent, lumped, oofem::NonLinearLhs, OOFEM_ERROR, solution, oofem::FloatArray::times(), oofem::EngngModel::updateSharedDofManagers(), and oofem::FloatArray::zero().
Referenced by solveYourselfAt().
|
virtual |
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 456 of file transienttransportproblem.C.
References oofem::EngngModel::giveCurrentMetaStep(), oofem::EngngModel::giveDomain(), giveNumericalMethod(), oofem::NumericalMethod::setDomain(), and oofem::EngngModel::updateDomainLinks().
|
virtual |
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 388 of file transienttransportproblem.C.
References oofem::EngngModel::updateYourself().
|
protected |
Definition at line 80 of file transienttransportproblem.h.
Referenced by giveNextStep(), giveSolutionStepWhenIcApply(), giveUnknownComponent(), initializeFrom(), and updateComponent().
|
protected |
Definition at line 83 of file transienttransportproblem.h.
Referenced by giveDeltaT(), and initializeFrom().
|
protected |
Definition at line 81 of file transienttransportproblem.h.
Referenced by giveDeltaT(), and initializeFrom().
|
protected |
Definition at line 71 of file transienttransportproblem.h.
Referenced by forceEquationNumbering(), solveYourselfAt(), and updateComponent().
|
protected |
Definition at line 75 of file transienttransportproblem.h.
Referenced by solveYourselfAt(), and updateComponent().
|
protected |
Definition at line 87 of file transienttransportproblem.h.
Referenced by initializeFrom().
|
protected |
Definition at line 69 of file transienttransportproblem.h.
Referenced by applyIC(), giveField(), giveUnknownComponent(), initializeFrom(), restoreContext(), saveContext(), solveYourselfAt(), and updateComponent().
|
protected |
Definition at line 84 of file transienttransportproblem.h.
Referenced by updateComponent().
|
protected |
Definition at line 74 of file transienttransportproblem.h.
Referenced by solveYourselfAt(), and updateComponent().
|
protected |
Definition at line 84 of file transienttransportproblem.h.
Referenced by initializeFrom(), and updateComponent().
|
protected |
Definition at line 85 of file transienttransportproblem.h.
Referenced by initializeFrom(), and updateComponent().
|
protected |
Numerical method used to solve the problem.
Definition at line 78 of file transienttransportproblem.h.
Referenced by giveNumericalMethod(), and solveYourselfAt().
|
protected |
Definition at line 82 of file transienttransportproblem.h.
Referenced by giveDeltaT(), giveDiscreteTime(), and initializeFrom().
|
protected |
Definition at line 73 of file transienttransportproblem.h.
Referenced by solveYourselfAt(), and updateComponent().
|
protected |
Definition at line 68 of file transienttransportproblem.h.
Referenced by initializeFrom(), and solveYourselfAt().