OOFEM 3.0
Loading...
Searching...
No Matches
engngm.h
Go to the documentation of this file.
1/*
2 *
3 * ##### ##### ###### ###### ### ###
4 * ## ## ## ## ## ## ## ### ##
5 * ## ## ## ## #### #### ## # ##
6 * ## ## ## ## ## ## ## ##
7 * ## ## ## ## ## ## ## ##
8 * ##### ##### ## ###### ## ##
9 *
10 *
11 * OOFEM : Object Oriented Finite Element Code
12 *
13 * Copyright (C) 1993 - 2025 Borek Patzak
14 *
15 *
16 *
17 * Czech Technical University, Faculty of Civil Engineering,
18 * Department of Structural Mechanics, 166 29 Prague, Czech Republic
19 *
20 * This library is free software; you can redistribute it and/or
21 * modify it under the terms of the GNU Lesser General Public
22 * License as published by the Free Software Foundation; either
23 * version 2.1 of the License, or (at your option) any later version.
24 *
25 * This program is distributed in the hope that it will be useful,
26 * but WITHOUT ANY WARRANTY; without even the implied warranty of
27 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
28 * Lesser General Public License for more details.
29 *
30 * You should have received a copy of the GNU Lesser General Public
31 * License along with this library; if not, write to the Free Software
32 * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
33 */
34
35#ifndef engngm_h
36#define engngm_h
37
38#include "oofemenv.h"
39#include "inputrecord.h"
40#include "intarray.h"
41#include "fieldmanager.h"
42#include "metastep.h"
43#include "timer.h"
44#include "assemblercallback.h"
45#include "chartype.h"
46#include "unknowntype.h"
47#include "varscaletype.h"
48#include "numericalcmpn.h"
49#include "valuemodetype.h"
50#include "internalstatetype.h"
51#include "problemmode.h"
52#include "fmode.h"
53#include "dofiditem.h"
54#include "contextoutputmode.h"
55#include "contextfilemode.h"
56#include "contextioresulttype.h"
57#include "metastep.h"
58#include "parallelcontext.h"
59#include "exportmodulemanager.h"
60#include "initmodulemanager.h"
61#include "monitormanager.h"
62#include "timestepcontroller.h"
63
64#ifdef __MPM_MODULE
65#include "../mpm/mpm.h"
66#include "../mpm/integral.h"
67#endif
68
69#ifdef __MPI_PARALLEL_MODE
70 #include "parallel.h"
71#endif
72
73#include <string>
74#include <memory>
75
77
78#define _IFT_EngngModel_nsteps "nsteps"
79#define _IFT_EngngModel_contextoutputstep "contextoutputstep"
80#define _IFT_EngngModel_renumberFlag "renumber"
81#define _IFT_EngngModel_profileOpt "profileopt"
82#define _IFT_EngngModel_nmsteps "nmsteps"
83#define _IFT_EngngModel_nonLinFormulation "nonlinform"
84#define _IFT_EngngModel_eetype "eetype"
85#define _IFT_EngngModel_parallelflag "parallelflag"
86#define _IFT_EngngModel_loadBalancingFlag "lbflag"
87#define _IFT_EngngModel_forceloadBalancingFlag "forcelb1"
88#define _IFT_EngngModel_initialGuess "initialguess"
89#define _IFT_EngngModel_referenceFile "referencefile"
90
91#define _IFT_EngngModel_lstype "lstype"
92#define _IFT_EngngModel_smtype "smtype"
93
94#define _IFT_EngngModel_suppressOutput "suppress_output" // Suppress writing to .out file
95
97
98namespace oofem {
99class Domain;
100class TimeStep;
101class Dof;
102class DofManager;
103class DataReader;
104class DataStream;
105class ErrorEstimator;
106class MetaStep;
107class MaterialInterface;
108class SparseMtrx;
109class NumericalMethod;
110class InitModuleManager;
111class ExportModuleManager;
112class FloatMatrix;
113class FloatArray;
114class LoadBalancer;
115class LoadBalancerMonitor;
116class oofegGraphicContext;
117class ProblemCommunicator;
118class ProcessCommunicatorBuff;
119class CommunicatorBuff;
120class ProcessCommunicator;
121class UnknownNumberingScheme;
122
123
134{
135protected:
138
139public:
141 FieldManager *giveFieldManager() { return & ( this->fieldManager ); }
142};
143
144
192{
193public:
201
210 //IG_Extrapolated = 2, ///< Assumes constant increment extrapolating @f$ {}^{n+1}x = {}^{n}x + \Delta t\delta{x}'@f$, where @f$ \delta x' = ({}^{n}x - {}^{n-1}x)/{}^{n}Delta t@f$.
211 };
212
213protected:
217 std :: vector< std :: unique_ptr< Domain > > domainList;
237 std :: vector< MetaStep > metaStepList;
239 std :: unique_ptr< TimeStep > stepWhenIcApply;
241 std :: unique_ptr< TimeStep > currentStep;
243 std :: unique_ptr< TimeStep > previousStep;
246
248 std :: string dataOutputFileName;
250 std :: string coreOutputFileName;
254 std :: string referenceFileName;
258
265
271 time_t startTime;
272
275
285 std::unique_ptr<ErrorEstimator> defaultErrEstimator;
287 std::unique_ptr<TimeStepController> timeStepController;
289 int rank;
294#ifdef __MPI_PARALLEL_MODE
297 #ifdef __USE_MPI
299 MPI_Comm comm;
300 #endif
303
304 std::unique_ptr<LoadBalancer> lb;
305 std::unique_ptr<LoadBalancerMonitor> lbm;
311
316
319#endif
323 std :: vector< ParallelContext > parallelContextList;
324
325
326#ifdef __MPM_MODULE
328 std :: map< std :: string, std::unique_ptr< Variable > > variableMap;
329 std :: vector < std :: unique_ptr< Term > > termList;
330 std :: vector < std :: unique_ptr< Integral > > integralList;
331#endif
332
333
334
335
338
340
341public:
345 EngngModel(int i, EngngModel * _master = NULL);
347 virtual ~EngngModel();
348 EngngModel(const EngngModel &) = delete;
349 EngngModel &operator=(const EngngModel &) = delete;
356 Domain *giveDomain(int n);
363 void setDomain(int i, Domain *ptr, bool iDeallocateOld = true);
365 int giveNumberOfDomains() { return (int)domainList.size(); }
366
367 const std :: string &giveDescription() const { return simulationDescription; }
368 const time_t &giveStartTime() { return startTime; }
369 bool giveSuppressOutput() const { return suppressOutput; }
370
374 virtual MaterialInterface *giveMaterialInterface(int n) { return nullptr; }
375 void setNumberOfEquations(int id, int neq) {
376 numberOfEquations = neq;
377 domainNeqs.at(id) = neq;
378 }
379 // input / output
381 FILE *giveOutputStream();
386 std :: string giveOutputBaseFileName() { return dataOutputFileName; }
387
391 std :: string giveReferenceFileName() { return referenceFileName;}
392
398 void letOutputBaseFileNameBe(const std :: string &src);
412 { contextOutputMode = contextMode; }
413
419 {
421 contextOutputStep = cStep;
422 }
423
424 double giveDeltaT(){return this->timeStepController->giveDeltaT();}
426 void setDeltaT(double dT){return this->timeStepController->setDeltaT(dT);}
427
428
433 void setProblemMode(problemMode pmode) { pMode = pmode; }
438 void setParallelMode(bool newParallelFlag);
445 void setProblemScale(problemScale pscale) { pScale = pscale; }
449 virtual void setRenumberFlag() { this->renumberFlag = true; }
451 virtual void resetRenumberFlag() { this->renumberFlag = false; }
452
456 double giveSolutionStepTime();
460 void giveAnalysisTime(int &rhrs, int &rmin, int &rsec, int &uhrs, int &umin, int &usec);
464 void terminateAnalysis();
465
466 // solving
473 virtual void solveYourself();
477 virtual void restartYourself(TimeStep *tS){;}
484 virtual void solveYourselfAt(TimeStep *tStep) { }
489 virtual void terminate(TimeStep *tStep);
495 virtual void doStepOutput(TimeStep *tStep);
499 void saveStepContext(TimeStep *tStep, ContextMode mode);
505 virtual void updateYourself(TimeStep *tStep);
515 virtual void initializeYourself(TimeStep *tStep);
521 virtual int initializeAdaptive(int tStepNumber) { return 0; }
522
527 virtual int giveNumberOfDomainEquations(int di, const UnknownNumberingScheme &num);
528
529 // management components
536 virtual double giveUnknownComponent(ValueModeType, TimeStep *, Domain *, Dof *) { return 0.0; }
537
545 virtual FieldPtr giveField (FieldType key, TimeStep *) { return FieldPtr();}
546 virtual FieldPtr giveField (InternalStateType key, TimeStep *) ;
547
548
550 EngngModel *giveMasterEngngModel() { return this->master; }
551
553 virtual double giveLoadLevel() { return 1.0; }
554
556 virtual double giveEigenValue(int eigNum) { return 0.0; }
558 virtual void setActiveVector(int i) { }
566 int updateSharedDofManagers(FloatArray &answer, const UnknownNumberingScheme &s, int ExchangeTag);
574 int exchangeRemoteElementData(int ExchangeTag);
579 virtual int giveCurrentNumberOfIterations() {return 1;}
580
581#ifdef __MPI_PARALLEL_MODE
583 MPI_Comm giveParallelComm() { return this->comm; }
595 int packRemoteElementData(ProcessCommunicator &processComm);
607 int unpackRemoteElementData(ProcessCommunicator &processComm);
615 int packDofManagers(ArrayWithNumbering *src, ProcessCommunicator &processComm);
623 int unpackDofManagers(ArrayWithNumbering *dest, ProcessCommunicator &processComm);
624
626 if ( t == PC_default ) {
627 return communicator;
628 } else if ( t == PC_nonlocal ) {
629 return nonlocCommunicator;
630 } else {
631 return NULL;
632 }
633 }
634#endif
635 void initializeCommMaps(bool forceInit = false);
641 virtual int instanciateYourself(DataReader &dr, InputRecord &ir, const char *outFileName, const char *desc);
648 void Instanciate_init();
654 virtual void initializeFrom(InputRecord &ir);
656 int instanciateDomains(DataReader &dr);
658 int instanciateMetaSteps(DataReader &dr);
660 virtual int instanciateDefaultMetaStep(InputRecord &ir);
661
671 virtual void updateAttributes(MetaStep *mStep);
678 this->timeStepController->initMetaStepAttributes(mStep);
679 }
680
689 virtual void saveContext(DataStream &stream, ContextMode mode);
704 virtual void restoreContext(DataStream &stream, ContextMode mode);
711 virtual void updateDomainLinks();
713 MetaStep *giveCurrentMetaStep();
717 virtual TimeStep *giveCurrentStep(bool force = false) {
718 if ( master && (!force)) {
719 return master->giveCurrentStep();
720 } else {
721 return currentStep.get();
722 //timeStepController->giveCurrentStep();
723 }
724 }
725
729 virtual void adaptTimeStep(double nIter){
730 //timeStepController->adaptTimeStep(nIter);
731 }
732
733
734
738 virtual TimeStep *givePreviousStep(bool force = false) {
739 if ( master && (!force)) {
740 return master->givePreviousStep();
741 } else {
742 return previousStep.get();
743 //timeStepController->givePreviousStep();
744 }
745 }
746
747 virtual TimeStep *giveNextStep() { return this->timeStepController->giveNextStep(); }
749 virtual void preInitializeNextStep() {}
753 virtual TimeStep *giveSolutionStepWhenIcApply(bool force = false) {
754 if ( master && (!force)) {
755 return master->giveSolutionStepWhenIcApply();
756 } else {
757 return stepWhenIcApply.get();
758 }
759 }
760
763 virtual int giveNumberOfFirstStep(bool force = false) {
764 if ( master && (!force)) {
765 return master->giveNumberOfFirstStep();
766 } else {
767 return 1;
768 }
769 }
770
771 int giveNumberOfMetaSteps() { return this->timeStepController->giveNumberOfMetaSteps(); }
773 MetaStep *giveMetaStep(int i) {return this->timeStepController->giveMetaStep(i);}
777 int giveNumberOfSteps(bool force = false) {
778 if ( master && (!force)) {
779 return master->giveNumberOfSteps();
780 } else {
781 return timeStepController->giveNumberOfSteps();
782 }
783 }
784
785 virtual double giveEndOfTimeOfInterest() { return 0.; }
789 virtual NumericalMethod *giveNumericalMethod(MetaStep *mStep) { return nullptr; }
794
795
797 virtual double giveInitialTime(){return 0.;}
801 virtual double giveFinalTime()
802 {
803 return timeStepController->giveFinalTime();
804 }
805
813 virtual int giveNewEquationNumber(int domain, DofIDItem) { return ++domainNeqs.at(domain); }
821 virtual int giveNewPrescribedEquationNumber(int domain, DofIDItem) { return ++domainPrescribedNeqs.at(domain); }
827 std :: string giveContextFileName(int tStepNumber, int stepVersion) const;
833 std :: string giveDomainFileName(int domainNum, int domainSerNum) const;
834 virtual void updateComponent(TimeStep *tStep, NumericalCmpn cmpn, Domain *d);
835 virtual void initForNewIteration(Domain *d, TimeStep *tStep, int iterationNumber, const FloatArray &solution);
843 virtual void updateSolution(FloatArray &solutionVector, TimeStep *tStep, Domain *d);
852 virtual void updateInternalRHS(FloatArray &answer, TimeStep *tStep, Domain *d, FloatArray *eNorm);
862 virtual void updateMatrix(SparseMtrx &mat, TimeStep *tStep, Domain *d);
869 virtual void initStepIncrements();
878 virtual int forceEquationNumbering(int i);
888 virtual int forceEquationNumbering();
908 virtual int requiresUnknownsDictionaryUpdate() { return false; }
915 virtual bool requiresEquationRenumbering(TimeStep *tStep) { return renumberFlag; }
916 //virtual int supportsBoundaryConditionChange () {return 0;}
930 virtual int giveUnknownDictHashIndx(ValueModeType mode, TimeStep *tStep) { return 0; }
937 virtual bool newDofHandling() { return false; }
942 virtual ParallelContext *giveParallelContext(int n);
947 virtual void initParallelContexts();
948
957 virtual void assemble(SparseMtrx &answer, TimeStep *tStep,
958 const MatrixAssembler &ma, const UnknownNumberingScheme &s, Domain *domain);
968 virtual void assemble(SparseMtrx &answer, TimeStep *tStep,
969 const MatrixAssembler &ma, const UnknownNumberingScheme &r_s, const UnknownNumberingScheme &c_s, Domain *domain);
984 void assembleVector(FloatArray &answer, TimeStep *tStep, const VectorAssembler &va, ValueModeType mode,
985 const UnknownNumberingScheme &s, Domain *domain, FloatArray *eNorms = NULL);
997 void assembleVectorFromDofManagers(FloatArray &answer, TimeStep *tStep, const VectorAssembler &va, ValueModeType mode,
998 const UnknownNumberingScheme &s, Domain *domain, FloatArray *eNorms = NULL);
1010 void assembleVectorFromElements(FloatArray &answer, TimeStep *tStep, const VectorAssembler &va, ValueModeType mode,
1011 const UnknownNumberingScheme &s, Domain *domain, FloatArray *eNorms = NULL);
1012
1023 void assembleVectorFromBC(FloatArray &answer, TimeStep *tStep, const VectorAssembler &va, ValueModeType mode,
1024 const UnknownNumberingScheme &s, Domain *domain, FloatArray *eNorms = NULL);
1025
1035 void assembleExtrapolatedForces(FloatArray &answer, TimeStep *tStep, CharType type, Domain *domain);
1036
1037 void assemblePrescribedExtrapolatedForces(FloatArray &answer, TimeStep *tStep, CharType type, Domain *domain);
1038
1039
1040#ifdef __MPM_MODULE
1042 std :: vector < std :: unique_ptr< Integral > > & giveIntegralList() {
1043 return this->integralList;
1044 }
1045 void addIntegral (std :: unique_ptr< Integral > obj) {
1046 integralList.push_back(std::move(obj));
1047 }
1048 // Needed for some of the boost-python bindings. NOTE: This takes ownership of the pointers, so it's actually completely unsafe.
1049 void py_addIntegral(Integral *obj) {
1050 std::size_t size = integralList.size();
1051 integralList.resize(size+1);
1052 integralList[size].reset(obj);
1053 }
1054 const Variable* giveVariableByName (std::string name) {
1055 if(variableMap.find(name)==variableMap.end()) OOFEM_ERROR("Unknown MPM variable '%s'",name.c_str());
1056 return variableMap[name].get();
1057
1058 }
1059 const Term* giveTerm (int indx) {
1060 // @BP: add better error handling than provided by at()
1061 if(indx<1 || (int)termList.size()<indx) OOFEM_ERROR("MPM term number %d outside of valid range 1..%d",indx,(int)termList.size());
1062 return termList[indx-1].get();
1063 }
1066 int instanciateMPM (DataReader &dr, InputRecord &ir);
1067 // end mpm experimental
1068#endif
1069
1070protected:
1077 virtual void packMigratingData(TimeStep *tStep) { }
1084 virtual void unpackMigratingData(TimeStep *tStep) { }
1085
1086public:
1091 virtual int checkConsistency() { return 1; }
1096 virtual int checkProblemConsistency();
1101 virtual void init();
1106 virtual void postInitialize();
1112 virtual void printOutputAt(FILE *file, TimeStep *tStep);
1113 virtual void printOutputAt(FILE *file, TimeStep *tStep, const IntArray &nodeSets, const IntArray &elementSets);
1121 void outputNodes(FILE *file, Domain &domain, TimeStep *tStep, int setNum);
1129 void outputElements(FILE *file, Domain &domain, TimeStep *tStep, int setNum);
1130
1131 // input / output
1133 void printYourself();
1134
1143 virtual void printDofOutputAt(FILE *stream, Dof *iDof, TimeStep *tStep);
1144
1145
1146 // identification
1148 virtual const char *giveClassName() const = 0;
1150 virtual int useNonlocalStiffnessOption() { return 0; }
1152 bool isParallel() const { return ( parallelFlag != 0 ); }
1154 int giveRank() const { return rank; }
1156 int giveNumberOfProcesses() const { return numProcs; }
1157
1158
1166 /*
1167 * Returns Load Response Mode of receiver.
1168 * This value indicates, whether nodes and elements should assemble
1169 * total or incremental load vectors.
1170 *
1171 * virtual LoadResponseMode giveLoadResponseMode () {return TotalLoad;}
1172 */
1174 EngngModelContext *giveContext() { return this->context; }
1176 virtual int giveNumberOfSlaveProblems() { return 0; }
1178 virtual EngngModel *giveSlaveProblem(int i) { return NULL; }
1179
1181 virtual bool giveEquationScalingFlag() { return false; }
1183 virtual double giveVariableScale(VarScaleType varId) { return 1.0; }
1184
1185
1195 virtual int estimateMaxPackSize(IntArray &commMap, DataStream &buff, int packUnpackType) { return 0; }
1196#ifdef __MPI_PARALLEL_MODE
1204 virtual void balanceLoad(TimeStep *tStep);
1206 virtual LoadBalancer *giveLoadBalancer() { return NULL; }
1209#endif
1211 void initParallel();
1213 EngngModel *giveEngngModel() { return this; }
1214 virtual bool isElementActivated( int elemNum ) { return true; }
1215 virtual bool isElementActivated( Element *e ) { return true; }
1217 TimeStepController* giveTimeStepController() { return this->timeStepController.get(); }
1218
1219#ifdef __OOFEG
1220 virtual void drawYourself(oofegGraphicContext &gc);
1221 virtual void drawElements(oofegGraphicContext &gc);
1222 virtual void drawNodes(oofegGraphicContext &gc);
1226 virtual void showSparseMtrxStructure(int type, oofegGraphicContext &gc, TimeStep *tStep) { }
1227#endif
1228
1229#ifdef _PYBIND_BINDINGS
1230 void setNumberOfDomains(const int& i) {this->ndomains=i;}
1231 const int& getNumberOfDomains() const {return this->ndomains;}
1232#endif
1233
1235 std :: string errorInfo(const char *func) const;
1236};
1237} // end namespace oofem
1238#endif // engngm_h
FieldManager * giveFieldManager()
Definition engngm.h:141
FieldManager fieldManager
Common fieldManager providing shared field register for the problem.
Definition engngm.h:137
int giveNumberOfTimeStepWhenIcApply()
Returns the time step number, when initial conditions should apply.
Definition engngm.h:787
EngngModelContext * context
Context.
Definition engngm.h:277
int giveContextOutputStep() const
Definition engngm.h:406
virtual void showSparseMtrxStructure(int type, oofegGraphicContext &gc, TimeStep *tStep)
Definition engngm.h:1226
int parallelFlag
Flag indicating that the receiver runs in parallel.
Definition engngm.h:281
problemScale giveProblemScale() const
Returns scale in multiscale simulation.
Definition engngm.h:447
std::unique_ptr< LoadBalancer > lb
Load Balancer.
Definition engngm.h:304
virtual int giveNewPrescribedEquationNumber(int domain, DofIDItem)
Definition engngm.h:821
std ::vector< std ::unique_ptr< Domain > > domainList
List of problem domains.
Definition engngm.h:217
enum fMode nonLinFormulation
Type of non linear formulation (total or updated formulation).
Definition engngm.h:283
void setContextOutputMode(ContextOutputMode contextMode)
Definition engngm.h:411
std::unique_ptr< TimeStepController > timeStepController
Time Step controller is responsible for collecting data from analysis, elements, and materials,...
Definition engngm.h:287
void setProblemMode(problemMode pmode)
Definition engngm.h:433
virtual EngngModel * giveSlaveProblem(int i)
Returns i-th slave problem.
Definition engngm.h:1178
virtual void preInitializeNextStep()
Does a pre-initialization of the next time step (implement if necessarry).
Definition engngm.h:749
virtual LoadBalancer * giveLoadBalancer()
Definition engngm.h:1206
bool force_load_rebalance_in_first_step
Debug flag forcing load balancing after first step.
Definition engngm.h:309
EngngModel(const EngngModel &)=delete
double giveDeltaT()
Returns time step size from the time step controlelr.
Definition engngm.h:424
int numberOfEquations
Total number of equation in current time step.
Definition engngm.h:221
virtual void setRenumberFlag()
Sets the renumber flag to true.
Definition engngm.h:449
int numProcs
Total number of collaborating processes.
Definition engngm.h:291
virtual int giveCurrentNumberOfIterations()
Definition engngm.h:579
int rank
Domain rank in a group of collaborating processes (0..groupSize-1).
Definition engngm.h:289
MPI_Comm comm
Communication object for this engineering model.
Definition engngm.h:299
virtual TimeStep * giveCurrentStep(bool force=false)
Definition engngm.h:717
std::string giveOutputBaseFileName()
Definition engngm.h:386
int giveNumberOfProcesses() const
Returns the number of collaborating processes.
Definition engngm.h:1156
std::string dataOutputFileName
Path to output stream.
Definition engngm.h:248
virtual FieldPtr giveField(FieldType key, TimeStep *)
Definition engngm.h:545
virtual LoadBalancerMonitor * giveLoadBalancerMonitor()
Definition engngm.h:1208
int giveRank() const
Returns domain rank in a group of collaborating processes (0..groupSize-1).
Definition engngm.h:1154
void setNumberOfEquations(int id, int neq)
Definition engngm.h:375
TimeStepController * giveTimeStepController()
Returns the time step controller.
Definition engngm.h:1217
ExportModuleManager * giveExportModuleManager()
Returns receiver's export module manager.
Definition engngm.h:791
ProblemCommunicator * communicator
Communicator.
Definition engngm.h:315
std::string simulationDescription
Definition engngm.h:339
virtual TimeStep * giveNextStep()
Returns next time step (next to current step) of receiver.
Definition engngm.h:747
EngngModel * giveEngngModel()
Returns reference to itself -> required by communicator.h.
Definition engngm.h:1213
virtual double giveUnknownComponent(ValueModeType, TimeStep *, Domain *, Dof *)
Definition engngm.h:536
EngngModelContext * giveContext()
Context requesting service.
Definition engngm.h:1174
bool renumberFlag
Renumbering flag (renumbers equations after each step, necessary if Dirichlet BCs change).
Definition engngm.h:229
@ IG_None
No special treatment for new iterations. Probably means ending up using for all free dofs.
Definition engngm.h:208
@ IG_Tangent
Solves an approximated tangent problem from the last iteration. Useful for changing Dirichlet boundar...
Definition engngm.h:209
EngngModel * giveMasterEngngModel()
Returns the master engnmodel.
Definition engngm.h:550
EngngModel(int i, EngngModel *_master=NULL)
Definition engngm.C:99
std ::unique_ptr< TimeStep > previousStep
Previous time step.
Definition engngm.h:243
int giveNumberOfDomains()
Returns number of domains in problem.
Definition engngm.h:365
virtual void adaptTimeStep(double nIter)
Definition engngm.h:729
ContextOutputMode giveContextOutputMode() const
Definition engngm.h:402
std::string coreOutputFileName
String with core output file name.
Definition engngm.h:250
virtual NumericalMethod * giveNumericalMethod(MetaStep *mStep)
Returns reference to receiver's numerical method.
Definition engngm.h:789
int equationNumberingCompleted
Equation numbering completed flag.
Definition engngm.h:233
bool profileOpt
Profile optimized numbering flag (using Sloan's algorithm).
Definition engngm.h:231
char processor_name[PROCESSOR_NAME_LENGTH]
Processor name.
Definition engngm.h:296
virtual void updateDofUnknownsDictionary(DofManager *, TimeStep *)
Definition engngm.h:922
std::string referenceFileName
String with reference file name.
Definition engngm.h:254
int ndomains
Number of receiver domains.
Definition engngm.h:215
virtual void solveYourselfAt(TimeStep *tStep)
Definition engngm.h:484
ProblemCommunicator * nonlocCommunicator
NonLocal Communicator. Necessary when nonlocal constitutive models are used.
Definition engngm.h:318
void setUDContextOutputMode(int cStep)
Definition engngm.h:418
int numberOfSteps
Total number of time steps.
Definition engngm.h:219
virtual int initializeAdaptive(int tStepNumber)
Definition engngm.h:521
virtual int checkConsistency()
Definition engngm.h:1091
Domain * giveDomain(int n)
Definition engngm.C:1936
std ::unique_ptr< TimeStep > currentStep
Current time step.
Definition engngm.h:241
const time_t & giveStartTime()
Definition engngm.h:368
int nonlocalExt
Flag indicating if nonlocal extension active, which will cause data to be sent between shared element...
Definition engngm.h:293
virtual void resetRenumberFlag()
Sets the renumber flag to false.
Definition engngm.h:451
@ EngngModel_Unknown_Mode
Definition engngm.h:194
bool loadBalancingFlag
If set to true, load balancing is active.
Definition engngm.h:307
void setDomain(int i, Domain *ptr, bool iDeallocateOld=true)
Definition engngm.C:1948
virtual ErrorEstimator * giveDomainErrorEstimator(int n)
Definition engngm.h:372
int giveNumberOfSteps(bool force=false)
Definition engngm.h:777
time_t startTime
Solution start time.
Definition engngm.h:271
virtual double giveLoadLevel()
Returns the current load level.
Definition engngm.h:553
MetaStep * giveMetaStep(int i)
Returns the i-th meta step.
Definition engngm.h:773
EngngModelTimer timer
E-model timer.
Definition engngm.h:279
@ InternalForcesExchangeTag
Definition engngm.h:321
@ RemoteElementExchangeTag
Definition engngm.h:321
InitModuleManager initModuleManager
Initialization module manager.
Definition engngm.h:262
std ::vector< ParallelContext > parallelContextList
List where parallel contexts are stored.
Definition engngm.h:323
virtual double giveEigenValue(int eigNum)
Only relevant for eigen value analysis. Otherwise returns zero.
Definition engngm.h:556
virtual double giveVariableScale(VarScaleType varId)
Returns the scale factor for given variable type.
Definition engngm.h:1183
virtual MaterialInterface * giveMaterialInterface(int n)
Definition engngm.h:374
std ::vector< MetaStep > metaStepList
List of problem metasteps.
Definition engngm.h:237
EngngModel & operator=(const EngngModel &)=delete
problemMode pMode
Domain mode.
Definition engngm.h:267
void setProblemScale(problemScale pscale)
Definition engngm.h:445
virtual TimeStep * givePreviousStep(bool force=false)
Definition engngm.h:738
int nMetaSteps
Number of meta steps.
Definition engngm.h:235
virtual int giveNumberOfSlaveProblems()
Returns number of slave problems.
Definition engngm.h:1176
virtual fMode giveFormulation()
Definition engngm.h:1165
bool giveSuppressOutput() const
Definition engngm.h:369
void setDeltaT(double dT)
Returns time step size through the time step controlelr.
Definition engngm.h:426
CommunicatorBuff * commBuff
Common Communicator buffer.
Definition engngm.h:313
std::string giveReferenceFileName()
Definition engngm.h:391
EngngModelTimer * giveTimer()
Returns reference to receiver timer (EngngModelTimer).
Definition engngm.h:793
const std::string & giveDescription() const
Definition engngm.h:367
ExportModuleManager exportModuleManager
Export module manager.
Definition engngm.h:260
virtual int giveNewEquationNumber(int domain, DofIDItem)
Definition engngm.h:813
ContextOutputMode contextOutputMode
Domain context output mode.
Definition engngm.h:256
EngngModel * master
Master e-model; if defined receiver is in maintained (slave) mode.
Definition engngm.h:274
virtual int useNonlocalStiffnessOption()
Returns nonzero if nonlocal stiffness option activated.
Definition engngm.h:1150
ProblemCommunicator * giveProblemCommunicator(EngngModelCommType t)
Definition engngm.h:625
MPI_Comm giveParallelComm()
Returns the communication object of reciever.
Definition engngm.h:583
virtual void unpackMigratingData(TimeStep *tStep)
Definition engngm.h:1084
std::unique_ptr< LoadBalancerMonitor > lbm
Definition engngm.h:305
void initMetaStepAttributes(MetaStep *mStep)
Definition engngm.h:677
virtual int requiresUnknownsDictionaryUpdate()
Definition engngm.h:908
virtual void restartYourself(TimeStep *tS)
Definition engngm.h:477
virtual TimeStep * giveSolutionStepWhenIcApply(bool force=false)
Definition engngm.h:753
virtual int giveUnknownDictHashIndx(ValueModeType mode, TimeStep *tStep)
Definition engngm.h:930
virtual bool giveEquationScalingFlag()
Returns the Equation scaling flag, which is used to indicate that governing equation(s) are scaled,...
Definition engngm.h:1181
int giveNumberOfMetaSteps()
Return number of meta steps.
Definition engngm.h:771
virtual double giveEndOfTimeOfInterest()
Returns end of time interest (time corresponding to end of time integration).
Definition engngm.h:785
bool suppressOutput
Flag for suppressing output to file.
Definition engngm.h:337
MonitorManager monitorManager
Monitor manager.
Definition engngm.h:264
FILE * outputStream
Output stream.
Definition engngm.h:252
virtual int giveNumberOfFirstStep(bool force=false)
Definition engngm.h:763
virtual bool isElementActivated(Element *e)
Definition engngm.h:1215
virtual void setActiveVector(int i)
Only relevant for eigen value analysis. Otherwise does noting.
Definition engngm.h:558
virtual double giveInitialTime()
return time at the begining of analysis
Definition engngm.h:797
virtual const char * giveClassName() const =0
Returns class name of the receiver.
bool isParallel() const
Returns true if receiver in parallel mode.
Definition engngm.h:1152
virtual double giveFinalTime()
Definition engngm.h:801
std ::unique_ptr< TimeStep > stepWhenIcApply
Solution step when IC (initial conditions) apply.
Definition engngm.h:239
int number
Receivers id.
Definition engngm.h:245
problemMode giveProblemMode() const
Returns domain mode.
Definition engngm.h:440
int contextOutputStep
Definition engngm.h:257
virtual bool requiresEquationRenumbering(TimeStep *tStep)
Definition engngm.h:915
IntArray domainPrescribedNeqs
Number of prescribed equations per domain.
Definition engngm.h:227
virtual bool isElementActivated(int elemNum)
Definition engngm.h:1214
std::unique_ptr< ErrorEstimator > defaultErrEstimator
Error estimator. Useful for adaptivity, or simply printing errors output.
Definition engngm.h:285
int numberOfPrescribedEquations
Total number or prescribed equations in current time step.
Definition engngm.h:223
virtual void packMigratingData(TimeStep *tStep)
Definition engngm.h:1077
problemScale pScale
Multiscale mode.
Definition engngm.h:269
virtual bool newDofHandling()
Definition engngm.h:937
virtual int estimateMaxPackSize(IntArray &commMap, DataStream &buff, int packUnpackType)
Definition engngm.h:1195
IntArray domainNeqs
Number of equations per domain.
Definition engngm.h:225
#define OOFEM_ERROR(...)
Definition error.h:79
long ContextMode
Definition contextmode.h:43
std::string errorInfo(const char *func)
Definition error.C:47
VarScaleType
Type determining the scale corresponding to particular variable.
fMode
Definition fmode.h:42
FieldType
Physical type of field.
Definition field.h:64
FloatArrayF< N > assemble(const FloatArrayF< M > &x, int const (&c)[M])
Assemble components into zero matrix.
@ COM_UserDefined
Input attribute of domain (each n-th step).
std::shared_ptr< Field > FieldPtr
Definition field.h:78
problemScale
Corresponds to macro- and micro-problem in multiscale simulations.
Definition problemmode.h:45
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
#define OOFEM_EXPORT
Definition oofemcfg.h:7
#define PROCESSOR_NAME_LENGTH
Definition parallel.h:56
Helper struct to pass array and numbering scheme as a single argument.
Definition engngm.h:197
const UnknownNumberingScheme * numbering
Definition engngm.h:199

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