IFEM  90A354
SIMbase.h
Go to the documentation of this file.
1 // $Id$
2 //==============================================================================
12 //==============================================================================
13 
14 #ifndef _SIM_BASE_H
15 #define _SIM_BASE_H
16 
17 #include "SIMadmin.h"
18 #include "SIMdependency.h"
19 #include "TimeDomain.h"
20 #include "Property.h"
21 #include "MatVec.h"
22 #include <set>
23 
24 class IntegrandBase;
25 class NormBase;
26 class ForceBase;
27 class AnaSol;
28 class SAM;
29 class AlgEqSystem;
30 class LinSolParams;
31 class SystemMatrix;
32 class SystemVector;
33 class FunctionBase;
34 class RealFunc;
35 class VecFunc;
36 class TractionFunc;
37 class ScalarFunc;
38 class Vec4;
39 class Vec3;
40 
41 
46 struct Mode
47 {
48  int eigNo;
49  double eigVal;
50  double damping;
53 
55  Mode() : eigNo(0), eigVal(0.0), damping(0.0) {}
57  bool orthonormalize(const SystemMatrix& mat);
59  bool computeDamping(const SystemMatrix& mat);
60 };
61 
62 
71 class SIMbase : public SIMadmin, public SIMdependency
72 {
73 protected:
75  explicit SIMbase(IntegrandBase* itg);
76 
77 public:
79  virtual ~SIMbase();
80 
81 
82  // Model input and pre-processing methods
83  // ======================================
84 
86  bool readModel(const char* fileName);
87 
89  virtual bool createFEMmodel(char resetNumb) = 0;
90 
94  virtual void clearProperties();
95 
97  virtual void initForSingleStep() {}
99  virtual void initForMultiStep() {}
100 
105  virtual bool preprocessC(const IntVec& ignored, bool fixDup, double time0);
106 
111  virtual bool merge(SIMbase* that,
112  const std::map<int,int>* old2new = nullptr, int poff = 0);
113 
120  bool initSystem(LinAlg::MatrixType mType,
121  size_t nMats = 1, size_t nVec = 1, size_t nScl = 0,
122  bool withRF = false);
123 
125  bool initSystem(const SIMbase* that);
126 
128  void initLHSbuffers();
129 
134  bool setAssociatedRHS(size_t iMat, size_t iVec);
135 
140  bool setMode(int mode, bool needIntegr = true, bool resetSol = false);
141 
145  void setIntegrationPrm(unsigned short int i, double prm);
146 
151  void setQuadratureRule(size_t ng, bool redimBuffers = false,
152  bool printQP = false);
153 
155  virtual bool printProblem() const;
156 
158  const IntegrandBase* getProblem() const { return myProblem; }
159 
162  { return nullptr; }
163 
167  void clearProblem() { myProblem = nullptr; }
168 
173  virtual std::string getName() const { return "SIMbase"; }
175  virtual bool mixedProblem() const { return false; }
176 
178  const LinSolParams* getSolParams() const { return mySolParams; }
179 
181  virtual unsigned short int getNoParamDim() const = 0;
183  virtual size_t getNoSpaceDim() const { return nsd; }
186  size_t getNoFields(int basis = 0) const;
188  size_t getNoDOFs(bool subSim = false) const;
191  size_t getNoNodes(int basis = 0) const;
196  size_t getNoElms(bool includeXelms = false, bool includeZelms = false) const;
198  size_t getNoSolutions(bool allocated = false) const;
200  int getNoPatches() const { return nGlPatches; }
202  size_t getNoEquations() const;
204  size_t getNoConstraints() const;
206  virtual size_t getNoRHS() const;
208  unsigned char getNoBasis() const;
209 
211  char getNodeType(int inod) const;
213  Vec4 getNodeCoord(int inod) const;
215  bool isFixed(int inod, int dof = 123) const;
217  int getGlobalNode(int node) const;
219  int getLocalNode(int node) const;
221  bool getElmNodes(IntVec& mnpc, int iel) const;
223  virtual std::vector<IntVec> getElmConnectivities() const = 0;
224 
229  void getBoundaryNodes(int pcode, IntVec& glbNodes,
230  std::vector<Vec3>* XYZ = nullptr) const;
231 
233  int findClosestNode(const Vec3&) const;
234 
236  IntVec getNodeSet(const std::string& setName) const;
237 
240  bool initDirichlet(double time = 0.0);
242  bool hasTimeDependentDirichlet() const;
243 
247  virtual bool updateDirichlet(double time = 0.0,
248  const Vector* prevSol = nullptr);
249 
251  virtual bool updateConfiguration(const Vector&) { return true; }
253  virtual bool updateRotations(const RealArray&, double = 0.0) { return true; }
254 
257  bool updateGrid(const RealArray& displ);
260  bool updateGrid(const std::string& field);
261 
264  void setRefined(int nref) { isRefined = nref; }
266  int getRefined() const { return isRefined; }
267 
269  bool hasElementActivator() const;
270 
274  void updateForNewElements(Vector& solution, const TimeDomain& time) const;
275 
276 
277  // Computational methods
278  // =====================
279 
285  virtual bool assembleSystem(const TimeDomain& time, const Vectors& prevSol,
286  bool newLHSmatrix = true, bool poorConvg = false);
287 
293  bool assembleSystem(double t0 = 0.0, const Vectors& pSol = Vectors())
294  { return this->assembleSystem(TimeDomain(t0),pSol); }
295 
300  bool extractLoadVec(Vector& loadVec, size_t idx = 0,
301  const char* hd = nullptr) const;
303  bool extractScalars(RealArray& values) const;
305  double extractScalar(size_t idx = 0) const;
306 
309  bool applyDirichlet(Vector& glbVec) const;
310 
318  bool solveEqSystem(Vector& solution, size_t idxRHS, double* rCond,
319  int printSol = 0, bool dumpEqSys = false,
320  const char* compName = "displacement");
321 
328  virtual bool solveSystem(Vector& solution, int printSol, double* rCond,
329  const char* compName = "displacement",
330  size_t idxRHS = 0)
331  {
332  return this->solveEqSystem(solution,idxRHS,rCond,printSol,true,compName);
333  }
334 
339  bool solveSystem(Vector& solution, int printSol = 0,
340  const char* compName = "displacement")
341  {
342  return this->solveSystem(solution,printSol,nullptr,compName);
343  }
344 
349  bool solveSystem(Vectors& solution, int printSol = 0,
350  const char* cmpName = "displacement");
351 
359  void getWorstDofs(const Vector& x, const Vector& r,
360  size_t nWorst, double eps, int iteNorm,
361  std::map<std::pair<int,int>,RealArray>& worst) const;
362 
369  virtual void iterationNorms(const Vector& x, const Vector& r, double& eNorm,
370  double& rNorm, double& dNorm) const;
371 
379  double solutionNorms(const Vector& x, double* inf = nullptr,
380  size_t* ind = nullptr, size_t nf = 0,
381  char type = 'D') const;
382 
390  bool solutionNorms(const TimeDomain& time,
391  const Vectors& psol, const Vectors& ssol,
392  Vectors& gNorm, Matrix* eNorm = nullptr,
393  const char* name = nullptr);
401  bool solutionNorms(const TimeDomain& time, const Vectors& psol,
402  Vectors& gNorm, Matrix* eNorm = nullptr)
403  { return this->solutionNorms(time,psol,Vectors(),gNorm,eNorm); }
412  bool solutionNorms(const Vector& psol, const Vectors& ssol,
413  Matrix& eNorm, Vectors& gNorm,
414  const char* name = nullptr)
415  {
416  return this->solutionNorms(TimeDomain(),Vectors(1,psol),ssol,
417  gNorm,&eNorm,name);
418  }
426  bool solutionNorms(const Vector& psol, Matrix& eNorm, Vectors& gNorm)
427  {
428  return this->solutionNorms(TimeDomain(),Vectors(1,psol),Vectors(),
429  gNorm,&eNorm);
430  }
439  bool solutionNorms(const Vector& psol, const Vectors& ssol, Vectors& gNorm,
440  const char* name = nullptr)
441  {
442  return this->solutionNorms(TimeDomain(),Vectors(1,psol),ssol,gNorm,
443  nullptr,name);
444  }
445 
453  virtual void printStep(int istep, const TimeDomain& time) const;
454 
460  virtual void printSolutionSummary(const Vector& solution, int printSol = 0,
461  const char* compName = nullptr,
462  std::streamsize outPrec = 0);
463 
469  bool getCurrentReactions(RealArray& RF, const Vector& psol,
470  int pcode = 0) const;
473  bool haveReactions(int pcode = 0) const;
475  virtual const RealArray* getReactionForces() const;
476 
485  bool systemModes(std::vector<Mode>& solution,
486  int nev, int ncv, int iop, double shift,
487  size_t iA = 0, size_t iB = 1);
492  bool systemModes(std::vector<Mode>& solution, size_t iA = 0, size_t iB = 1)
493  {
494  return this->systemModes(solution,opt.nev,opt.ncv,opt.eig,opt.shift,iA,iB);
495  }
496 
498  virtual bool haveBoundaryReactions(bool = false) const { return false; }
499 
505  bool assembleForces(const Vector& solution, double t0,
506  RealArray* R, Vector* S = nullptr);
507 
508 
509  // Post-processing methods
510  // =======================
511 
521  virtual bool project(Matrix& ssol, const Vector& psol,
522  SIMoptions::ProjectionMethod method = SIMoptions::GLOBAL,
523  const TimeDomain& time = TimeDomain()) const;
531  bool project(Vector& ssol, const Vector& psol,
532  SIMoptions::ProjectionMethod method = SIMoptions::GLOBAL,
533  size_t iComp = 0) const;
534 
538  bool projectAnaSol(Vector& ssol, SIMoptions::ProjectionMethod method) const;
539 
548  bool project(RealArray& values, const FunctionBase* f,
549  int basis = 1, int iField = 0, int nFields = 1,
550  SIMoptions::ProjectionMethod method = SIMoptions::GLOBAL,
551  double time = 0.0) const;
552 
556  bool evalSecondarySolution(Matrix& field, int pindx) const;
557 
559  virtual bool fieldProjections() const;
560 
562  virtual bool haveAnaSol() const { return mySol ? true : false; }
564  virtual bool haveDualSol() const { return dualField ? true : false; }
565 
570  NormBase* getNormIntegrand() const;
575  ForceBase* getBoundaryForceIntegrand(const Vec3* X0 = nullptr) const;
581 
583  const SAM* getSAM() const { return mySam; }
584 
585 protected:
589  VecFunc* getVecFunc(size_t patch, Property::Type ptype) const;
590 
600  virtual bool addConstraint(int patch, int lndx, int ldim, int dirs, int code,
601  int& ngnod, char basis = 1, bool ovrD = false);
602 
604  virtual void preprocessA() {}
606  virtual bool preprocessBeforeAsmInit(int&) { return true; }
608  virtual bool preprocessB() { return true; }
610  virtual void preprocessResultPoints() = 0;
611 
615  int renumberNodes(bool renumMNPC = false);
617  virtual bool renumberNodes(const std::map<int,int>&) { return true; }
618 
629  virtual bool extractPatchSolution(IntegrandBase* problem,
630  const Vectors& sol, size_t pindx) const;
631 
632 public:
639  void registerDependency(const std::string& name, SIMdependency* sim,
640  short int nvc = 1, unsigned char basis = 1);
641 
645  bool extractPatchSolution(const Vectors& sol, size_t pindx) const
646  { return this->extractPatchSolution(myProblem,sol,pindx); }
647 
655  size_t extractPatchSolution(const RealArray& sol, RealArray& vec,
656  const ASMbase* pch,
657  unsigned char nndof = 0,
658  unsigned char basis = 0) const;
659 
666  bool injectPatchSolution(RealArray& sol, const RealArray& vec,
667  const ASMbase* pch,
668  unsigned char nndof = 0,
669  unsigned char basis = 0) const;
670 
675  bool extractPatchElmRes(const Matrix& glbRes, Matrix& elRes, int pindx) const;
676 
682  int getLocalPatchIndex(int patchNo) const;
683 
685  const PatchVec& getFEModel() const { return myModel; }
691  ASMbase* getPatch(int idx, bool glbIndex = false) const;
692 
694  bool setPatchMaterial(size_t patch) const;
695 
697  const std::map<int,int>& getGlob2LocMap() const { return myGlb2Loc; }
698 
700  PropertyVec::const_iterator begin_prop() const { return myProps.begin(); }
702  PropertyVec::const_iterator end_prop() const { return myProps.end(); }
703 
707  SystemMatrix* getRayleighDampingMatrix(size_t iM = 1, size_t iK = 0) const;
709  SystemMatrix* getLHSmatrix(size_t idx = 0, bool copy = false) const;
711  SystemVector* getRHSvector(size_t idx = 0, bool copy = false) const;
713  void addToRHSvector(size_t idx, const SystemVector& vec, double scale = 1.0);
714 
716  RealFunc* getSclFunc(int code) const;
717 
719  void setMDflag(char flag) { mdFlag = flag; }
721  bool isFirst() const { return mdFlag <= 1; }
722 
724  void dumpEqSys(bool initialBlankLine = false);
726  void dumpSolVec(const Vector& x,bool isExpanded = true, bool expOnly = false);
727 
729  virtual int printNRforces(const IntVec& = {}) const { return 0; }
730 
731 protected:
733  char getMDflag() const { return mdFlag; }
735  virtual void shiftGlobalNums(int, int) {}
736 
738  virtual bool initMaterial(size_t) { return true; }
740  virtual bool initBodyLoad(size_t) { return true; }
742  virtual bool initNeumann(size_t) { return true; }
743 
746  const TimeDomain&) { return true; }
751  virtual bool assembleDiscreteItems(const IntegrandBase* itg,
752  const TimeDomain& time, const Vectors&)
753  { return this->assembleDiscreteTerms(itg,time); }
754 
756  virtual double externalEnergy(const Vectors& psol, const TimeDomain&) const;
757 
759  virtual bool postProcessNorms(Vectors&, Matrix*) { return true; }
760 
764  void generateThreadGroups(const Property& p, bool silence = false);
765 
767  void changeNumThreads();
768 
773  bool addMADOF(unsigned char basis, unsigned char nndof, bool other = true);
774 
776  double* theExtEnerg() { return &extEnergy; }
778  const double* getExtEnerg() const { return &extEnergy; }
779 
780 private:
784  const IntVec& getMADOF(unsigned char basis, unsigned char nndof) const;
785 
786 public:
787  static bool ignoreDirichlet;
788  static bool preserveNOrder;
789 
790 protected:
792  typedef std::map<int,RealFunc*> SclFuncMap;
794  typedef std::map<int,VecFunc*> VecFuncMap;
796  typedef std::map<int,TractionFunc*> TracFuncMap;
797 
799  typedef std::multimap<int,IntegrandBase*> IntegrandMap;
800 
801  // Model attributes
802  unsigned char nsd;
812 
813  std::vector<FunctionBase*> extrFunc;
814 
819  int nGlbNodes;
820  int isRefined;
821  bool lagMTOK;
822  bool fixZeros;
823 
825  struct DumpData
826  {
827  std::string fname;
829  std::set<int> step;
830  bool expand;
831  int count;
832  double eps;
833 
835  DumpData() : format(LinAlg::FLAT), expand(false), count(0), eps(1.0e-6) {}
836 
838  bool doDump() { return !fname.empty() && step.find(++count) != step.end(); }
839  };
840 
841  // Post-processing attributes
842  std::vector<DumpData> lhsDump;
843  std::vector<DumpData> rhsDump;
844  std::vector<DumpData> solDump;
845 
846  // Parallel computing attributes
850  std::map<int,int> myGlb2Loc;
851  const std::map<int,int>* g2l;
852  std::map<int,int> myDegenElm;
853  std::set<int> myDupNodes;
854 
855  // Equation solver attributes
860 
861 private:
862  size_t nIntGP;
863  size_t nBouGP;
864  size_t nDofS;
865  char mdFlag;
866 
868  mutable std::map<int,IntVec> extraMADOFs;
869 
870  mutable double extEnergy;
871  mutable Vector prevForces;
872 };
873 
874 #endif
std::vector< int > IntVec
General integer vector.
Definition: ASMbase.h:25
std::vector< Real > RealArray
A real-valued array without algebraic operations.
Definition: ImmersedBoundaries.h:29
Global algebraic operations on index 1-based matrices and vectors.
std::vector< Vector > Vectors
An array of real-valued vectors with algebraic operations.
Definition: MatVec.h:37
Representation of a distributed physical property.
std::vector< Property > PropertyVec
Vector of properties.
Definition: Property.h:64
Administration base class for FEM simulators.
Administration of simulators with dependencies to other simulators.
Time domain representation for time-dependent and nonlinear solvers.
Base class to check for internal boundary integrand contributions.
Definition: Interface.h:40
Base class for spline-based finite element (FE) assembly drivers.
Definition: ASMbase.h:70
Class for storage of general algebraic system of equations.
Definition: AlgEqSystem.h:33
Class for analytical solution fields (primary and secondary solution).
Definition: AnaSol.h:32
Base class representing a system level boundary force quantity.
Definition: IntegrandBase.h:523
Base class for unary spatial functions of arbitrary result type.
Definition: Function.h:147
Base class representing a system level integrated quantity.
Definition: IntegrandBase.h:42
Class for linear solver parameters.
Definition: LinSolParams.h:67
Base class representing a system level norm quantity.
Definition: IntegrandBase.h:372
Scalar-valued unary function of a spatial point.
Definition: Function.h:193
This class contains data and functions for the assembly of FE matrices.
Definition: SAM.h:39
Administration base class for FEM simulators.
Definition: SIMadmin.h:32
std::vector< int > IntVec
Convenience alias.
Definition: SIMadmin.h:34
SIMoptions & opt
Simulation control parameters.
Definition: SIMadmin.h:82
Base class for NURBS-based FEM simulators.
Definition: SIMbase.h:72
size_t getNoElms(bool includeXelms=false, bool includeZelms=false) const
Returns the model size in terms of number of elements.
Definition: SIMbase.C:835
SIMbase(IntegrandBase *itg)
The constructor initializes the pointers to dynamic data members.
Definition: SIMbase.C:48
virtual bool mixedProblem() const
Returns whether a mixed formulation is used (used by HDF5 output).
Definition: SIMbase.h:175
int findClosestNode(const Vec3 &) const
Finds the node that is closest to the given point X.
Definition: SIMbase.C:1131
virtual bool assembleDiscreteTerms(const IntegrandBase *, const TimeDomain &)
Assembles problem-dependent discrete terms, if any.
Definition: SIMbase.h:745
bool haveReactions(int pcode=0) const
Checks for total reaction forces associated with a boundary.
Definition: SIMbase.C:2105
void getWorstDofs(const Vector &x, const Vector &r, size_t nWorst, double eps, int iteNorm, std::map< std::pair< int, int >, RealArray > &worst) const
Finds the DOFs showing the worst convergence behavior.
Definition: SIMbase.C:1565
RealFunc * getSclFunc(int code) const
Returns a scalar function associated with code.
Definition: SIMbase.C:577
bool extractPatchSolution(const Vectors &sol, size_t pindx) const
Extracts all local solution vector(s) for a specified patch.
Definition: SIMbase.h:645
bool systemModes(std::vector< Mode > &solution, size_t iA=0, size_t iB=1)
Performs a generalized eigenvalue analysis of the assembled system.
Definition: SIMbase.h:492
virtual bool preprocessBeforeAsmInit(int &)
Specialized preprocessing performed before assembly initialization.
Definition: SIMbase.h:606
bool addMADOF(unsigned char basis, unsigned char nndof, bool other=true)
Adds a MADOF with an extraordinary number of DOFs on a given basis.
Definition: SIMbase.C:2628
bool solutionNorms(const TimeDomain &time, const Vectors &psol, Vectors &gNorm, Matrix *eNorm=nullptr)
Integrates some solution norm quantities.
Definition: SIMbase.h:401
const PatchVec & getFEModel() const
Returns a const reference to our FEM model.
Definition: SIMbase.h:685
virtual std::string getName() const
Returns the name of this simulator.
Definition: SIMbase.h:173
const double * getExtEnerg() const
Returns a const pointer to the external energy path integral value.
Definition: SIMbase.h:778
bool isFixed(int inod, int dof=123) const
Returns true if all DOFs in the specified global node are fixed.
Definition: SIMbase.C:1617
double extractScalar(size_t idx=0) const
Extracts an assembled global scalar quantity.
Definition: SIMbase.C:1414
virtual bool updateConfiguration(const Vector &)
Updates problem-dependent state based on the current solution.
Definition: SIMbase.h:251
bool evalSecondarySolution(Matrix &field, int pindx) const
Evaluates the secondary solution field for specified patch.
Definition: SIMbase.C:2582
std::map< int, RealFunc * > SclFuncMap
Scalar field container.
Definition: SIMbase.h:792
size_t getNoConstraints() const
Returns the number of constraint equations in the model.
Definition: SIMbase.C:876
bool solutionNorms(const Vector &psol, const Vectors &ssol, Vectors &gNorm, const char *name=nullptr)
Integrates some solution norm quantities.
Definition: SIMbase.h:439
ForceBase * getNodalForceIntegrand() const
Returns a pointer to a force integrand object for this simulator.
Definition: SIMbase.C:1760
virtual bool renumberNodes(const std::map< int, int > &)
Interface for renumbering of app-specific node number tables.
Definition: SIMbase.h:617
LinSolParams * myGl2Params
Input parameters for PETSc, for L2 projection.
Definition: SIMbase.h:859
virtual bool project(Matrix &ssol, const Vector &psol, SIMoptions::ProjectionMethod method=SIMoptions::GLOBAL, const TimeDomain &time=TimeDomain()) const
Projects the secondary solution associated with a primary solution.
Definition: SIMbase.C:2213
virtual bool haveBoundaryReactions(bool=false) const
Returns whether reaction forces are to be computed or not.
Definition: SIMbase.h:498
virtual bool preprocessB()
Preprocessing performed after the system assembly initialization.
Definition: SIMbase.h:608
virtual ASM::InterfaceChecker * getInterfaceChecker(size_t) const
Returns interface checker type for model.
Definition: SIMbase.h:161
double solutionNorms(const Vector &x, double *inf=nullptr, size_t *ind=nullptr, size_t nf=0, char type='D') const
Evaluates some norms of the primary solution vector.
Definition: SIMbase.C:1733
int getGlobalNode(int node) const
Returns the global node number from a process-local node number.
Definition: SIMbase.C:1627
std::set< int > myDupNodes
Set of duplicated nodes.
Definition: SIMbase.h:853
std::vector< DumpData > rhsDump
Right-hand-side vector dump specifications.
Definition: SIMbase.h:843
std::map< int, int > myDegenElm
Degenerated elements mapping.
Definition: SIMbase.h:852
void setMDflag(char flag)
Sets the multi-dimension simulator sequence flag.
Definition: SIMbase.h:719
bool solveSystem(Vector &solution, int printSol=0, const char *compName="displacement")
Solves the assembled linear system of equations for a given load.
Definition: SIMbase.h:339
int renumberNodes(bool renumMNPC=false)
Renumbers the global node numbers after resolving patch topology.
Definition: SIMbase.C:516
virtual bool initBodyLoad(size_t)
Initializes the body load properties for current patch.
Definition: SIMbase.h:740
bool readModel(const char *fileName)
Reads model data from the specified input file *fileName.
Definition: SIMbase.C:105
bool projectAnaSol(Vector &ssol, SIMoptions::ProjectionMethod method) const
Projects the analytical secondary solution, if any.
Definition: SIMbase.C:2478
size_t getNoSolutions(bool allocated=false) const
Returns the number of solution vectors.
Definition: SIMbase.C:864
bool extractLoadVec(Vector &loadVec, size_t idx=0, const char *hd=nullptr) const
Extracts the assembled load vector for inspection/visualization.
Definition: SIMbase.C:1371
void updateForNewElements(Vector &solution, const TimeDomain &time) const
Modifies the current solution vector when activating elements.
Definition: SIMbase.C:985
VecFuncMap myVectors
Vector property fields.
Definition: SIMbase.h:806
bool isFirst() const
Returns true, if this is the equation system owner.
Definition: SIMbase.h:721
const LinSolParams * getSolParams() const
Returns the linear equation solver parameters (for PETSc).
Definition: SIMbase.h:178
virtual bool updateRotations(const RealArray &, double=0.0)
Updates the nodal rotations for problems with rotational DOFs.
Definition: SIMbase.h:253
int getRefined() const
Returns current refinement status.
Definition: SIMbase.h:266
virtual bool assembleDiscreteItems(const IntegrandBase *itg, const TimeDomain &time, const Vectors &)
Assembles problem-dependent discrete terms, if any.
Definition: SIMbase.h:751
void addToRHSvector(size_t idx, const SystemVector &vec, double scale=1.0)
Adds a system vector to the given right-hand-side vector.
Definition: SIMbase.C:1713
bool initSystem(LinAlg::MatrixType mType, size_t nMats=1, size_t nVec=1, size_t nScl=0, bool withRF=false)
Allocates the system matrices of the FE problem to be solved.
Definition: SIMbase.C:584
size_t nIntGP
Number of interior integration points in the whole model.
Definition: SIMbase.h:862
bool systemModes(std::vector< Mode > &solution, int nev, int ncv, int iop, double shift, size_t iA=0, size_t iB=1)
Performs a generalized eigenvalue analysis of the assembled system.
Definition: SIMbase.C:2125
Vec4 getNodeCoord(int inod) const
Returns the spatial coordinates of the specified global node.
Definition: SIMbase.C:1601
std::vector< DumpData > lhsDump
Coefficient matrix dump specifications.
Definition: SIMbase.h:842
int getLocalNode(int node) const
Returns the process-local node number from a global node number.
Definition: SIMbase.C:1634
virtual void clearProperties()
Initializes the property containers of the model.
Definition: SIMbase.C:112
virtual void preprocessResultPoints()=0
Preprocesses the result sampling points.
int nGlbNodes
Total number of unique nodes in the model.
Definition: SIMbase.h:819
PatchVec myModel
The actual NURBS/spline model.
Definition: SIMbase.h:803
virtual ~SIMbase()
The destructor frees the dynamically allocated objects.
Definition: SIMbase.C:69
bool solveEqSystem(Vector &solution, size_t idxRHS, double *rCond, int printSol=0, bool dumpEqSys=false, const char *compName="displacement")
Solves the assembled linear system of equations for a given load.
Definition: SIMbase.C:1426
bool setPatchMaterial(size_t patch) const
Initializes material properties for the given patch.
Definition: SIMbase.C:2618
virtual bool addConstraint(int patch, int lndx, int ldim, int dirs, int code, int &ngnod, char basis=1, bool ovrD=false)
Preprocesses a user-defined Dirichlet boundary property.
Definition: SIMbase.C:434
bool extractPatchElmRes(const Matrix &glbRes, Matrix &elRes, int pindx) const
Extracts element results for a specified patch.
Definition: SIMbase.C:2607
SclFuncMap myScalars
Scalar property fields.
Definition: SIMbase.h:805
virtual bool haveDualSol() const
Returns whether a dual solution is available or not.
Definition: SIMbase.h:564
void generateThreadGroups(const Property &p, bool silence=false)
Generates element groups for multi-threading of boundary integrals.
Definition: SIMbase.C:548
std::multimap< int, IntegrandBase * > IntegrandMap
Property code to integrand map.
Definition: SIMbase.h:799
virtual bool initMaterial(size_t)
Initializes material properties for integration of interior terms.
Definition: SIMbase.h:738
size_t getNoDOFs(bool subSim=false) const
Returns the model size in terms of number of DOFs.
Definition: SIMbase.C:805
const std::map< int, int > * g2l
Pointer to global-to-local node mapping.
Definition: SIMbase.h:851
double * theExtEnerg()
Returns a pointer to the external energy path integral value.
Definition: SIMbase.h:776
std::map< int, VecFunc * > VecFuncMap
Vector field container.
Definition: SIMbase.h:794
PropertyVec::const_iterator end_prop() const
Returns the end of the property array.
Definition: SIMbase.h:702
static bool preserveNOrder
Set to true to preserve node ordering.
Definition: SIMbase.h:788
virtual double externalEnergy(const Vectors &psol, const TimeDomain &) const
Computes (possibly problem-dependent) external energy contribution.
Definition: SIMbase.C:2038
virtual bool assembleSystem(const TimeDomain &time, const Vectors &prevSol, bool newLHSmatrix=true, bool poorConvg=false)
Administers assembly of the linear equation system.
Definition: SIMbase.C:1168
bool setMode(int mode, bool needIntegr=true, bool resetSol=false)
Defines the solution mode before the element assembly is started.
Definition: SIMbase.C:662
PropertyVec::const_iterator begin_prop() const
Returns the beginning of the property array.
Definition: SIMbase.h:700
bool fixZeros
If true, constrain zero pivots before solving.
Definition: SIMbase.h:822
bool lagMTOK
Indicates if global multipliers is OK with multithreading.
Definition: SIMbase.h:821
void initLHSbuffers()
Initializes left-hand-side element matrix buffers for integrand.
Definition: SIMbase.C:642
virtual bool haveAnaSol() const
Returns whether an analytical solution is available or not.
Definition: SIMbase.h:562
bool assembleSystem(double t0=0.0, const Vectors &pSol=Vectors())
Administers assembly of the linear equation system.
Definition: SIMbase.h:293
virtual bool postProcessNorms(Vectors &, Matrix *)
Applies app-specific post-processing on norms.
Definition: SIMbase.h:759
virtual bool updateDirichlet(double time=0.0, const Vector *prevSol=nullptr)
Updates the time-dependent in-homogeneous Dirichlet coefficients.
Definition: SIMbase.C:934
void changeNumThreads()
Called if number of threads changes.
Definition: SIMbase.C:556
bool setAssociatedRHS(size_t iMat, size_t iVec)
Associates a system vector to a system matrix.
Definition: SIMbase.C:656
size_t getNoFields(int basis=0) const
Returns the number of primary solution fields.
Definition: SIMbase.C:799
const std::map< int, int > & getGlob2LocMap() const
Returns a const reference to our global-to-local node mapping.
Definition: SIMbase.h:697
int getNoPatches() const
Returns the total number of patches in the model.
Definition: SIMbase.h:200
bool extractScalars(RealArray &values) const
Extracts the assembled global scalar quantities.
Definition: SIMbase.C:1401
void getBoundaryNodes(int pcode, IntVec &glbNodes, std::vector< Vec3 > *XYZ=nullptr) const
Finds the list of global nodes associated with a boundary.
Definition: SIMbase.C:1052
void setIntegrationPrm(unsigned short int i, double prm)
Initializes an integration parameter for the integrand.
Definition: SIMbase.C:685
IntegrandBase * myProblem
The main integrand of this simulator.
Definition: SIMbase.h:808
char mdFlag
Sequence flag for multi-dimensional simulators.
Definition: SIMbase.h:865
std::vector< FunctionBase * > extrFunc
Extraction functions for VCP.
Definition: SIMbase.h:813
ForceBase * getBoundaryForceIntegrand(const Vec3 *X0=nullptr) const
Returns a pointer to a force integrand object for this simulator.
Definition: SIMbase.C:1754
unsigned char getNoBasis() const
Returns the number of bases in the model.
Definition: SIMbase.C:888
void clearProblem()
Clears the reference to the problem-specific data object.
Definition: SIMbase.h:167
SAM * mySam
Auxiliary data for FE assembly management.
Definition: SIMbase.h:857
static bool ignoreDirichlet
Set to true for free vibration analysis.
Definition: SIMbase.h:787
virtual bool solveSystem(Vector &solution, int printSol, double *rCond, const char *compName="displacement", size_t idxRHS=0)
Solves the assembled linear system of equations for a given load.
Definition: SIMbase.h:328
bool hasElementActivator() const
Returns true if an element activation function is specified.
Definition: SIMbase.C:978
bool solutionNorms(const Vector &psol, const Vectors &ssol, Matrix &eNorm, Vectors &gNorm, const char *name=nullptr)
Integrates some solution norm quantities.
Definition: SIMbase.h:412
bool injectPatchSolution(RealArray &sol, const RealArray &vec, const ASMbase *pch, unsigned char nndof=0, unsigned char basis=0) const
Injects a patch-wise solution vector into the global vector.
Definition: SIMbase.C:2561
std::map< int, IntVec > extraMADOFs
Additional MADOF arrays with varying DOF counts per node.
Definition: SIMbase.h:868
bool initDirichlet(double time=0.0)
Initializes time-dependent in-homogeneous Dirichlet coefficients.
Definition: SIMbase.C:910
const SAM * getSAM() const
Returns a const pointer to the SAM object of this simulator.
Definition: SIMbase.h:583
bool getCurrentReactions(RealArray &RF, const Vector &psol, int pcode=0) const
Computes the total reaction forces in the model.
Definition: SIMbase.C:2073
IntVec myLoc2Glb
Local-to-global node number mapping.
Definition: SIMbase.h:849
IntVec getNodeSet(const std::string &setName) const
Returns a predefined node set.
Definition: SIMbase.C:1114
bool applyDirichlet(Vector &glbVec) const
Applies the Dirichlet conditions to given vector.
Definition: SIMbase.C:1420
virtual bool extractPatchSolution(IntegrandBase *problem, const Vectors &sol, size_t pindx) const
Extracts all local solution vector(s) for a specified patch.
Definition: SIMbase.C:2443
virtual void iterationNorms(const Vector &x, const Vector &r, double &eNorm, double &rNorm, double &dNorm) const
Evaluates some iteration norms for convergence assessment.
Definition: SIMbase.C:1724
char getMDflag() const
Returns the multi-dimension simulator sequence flag.
Definition: SIMbase.h:733
FunctionBase * dualField
Dual solution field (extraction function)
Definition: SIMbase.h:811
std::vector< DumpData > solDump
Solution vector dump specifications.
Definition: SIMbase.h:844
AnaSol * mySol
Analytical/Exact solution.
Definition: SIMbase.h:810
virtual int printNRforces(const IntVec &={}) const
Prints out nodal reaction forces to the log stream.
Definition: SIMbase.h:729
SystemVector * getRHSvector(size_t idx=0, bool copy=false) const
Returns current system right-hand-side vector.
Definition: SIMbase.C:1704
virtual const RealArray * getReactionForces() const
Returns current reaction force container.
Definition: SIMbase.C:2119
virtual bool initNeumann(size_t)
Initializes for integration of Neumann terms for a given property.
Definition: SIMbase.h:742
Vector prevForces
Reaction forces of previous time step.
Definition: SIMbase.h:871
virtual void initForSingleStep()
Interface for app-specific single-step simulation initialisation.
Definition: SIMbase.h:97
virtual void printStep(int istep, const TimeDomain &time) const
Prints out load/time step identification.
Definition: SIMbase.C:1506
virtual bool merge(SIMbase *that, const std::map< int, int > *old2new=nullptr, int poff=0)
Merges the global equation system of that simulator with this.
Definition: SIMbase.C:473
bool updateGrid(const RealArray &displ)
Updates the grid coordinates.
Definition: SIMbase.C:946
int isRefined
Indicates if the model is adaptively refined.
Definition: SIMbase.h:820
virtual void preprocessA()
Preprocessing performed before the FEM model generation.
Definition: SIMbase.h:604
ASMbase * getPatch(int idx, bool glbIndex=false) const
Returns a pointer to a specified patch of our FEM model.
Definition: SIMbase.C:162
double extEnergy
Path integral of external forces.
Definition: SIMbase.h:870
virtual void registerDependency(SIMdependency *sim, const std::string &name, short int nvc, const PatchVec &patches, char diffBasis=0, int component=1)
Registers a dependency on a field from another SIM object.
Definition: SIMdependency.C:21
SystemMatrix * getRayleighDampingMatrix(size_t iM=1, size_t iK=0) const
Returns current Rayleigh system damping matrix.
Definition: SIMbase.C:1641
IntVec myPatches
Global patch numbers for current processor.
Definition: SIMbase.h:848
TracFuncMap myTracs
Traction property fields.
Definition: SIMbase.h:807
void setRefined(int nref)
Sets the refinement status (for restart of adaptive simulations).
Definition: SIMbase.h:264
IntegrandMap myInts
Set of all integrands involved.
Definition: SIMbase.h:809
void dumpSolVec(const Vector &x, bool isExpanded=true, bool expOnly=false)
Dumps a solution vector to file.
Definition: SIMbase.C:2754
virtual size_t getNoRHS() const
Returns the number of right-hand-side vectors.
Definition: SIMbase.C:882
virtual bool createFEMmodel(char resetNumb)=0
Creates the computational FEM model from the spline patches.
bool solutionNorms(const Vector &psol, Matrix &eNorm, Vectors &gNorm)
Integrates some solution norm quantities.
Definition: SIMbase.h:426
virtual void initForMultiStep()
Interface for app-specific multi-step simulation initialisation.
Definition: SIMbase.h:99
size_t getNoNodes(int basis=0) const
Returns the model size in terms of number of unique nodes.
Definition: SIMbase.C:811
size_t nDofS
Number of degrees of freedom in this sub-simulator.
Definition: SIMbase.h:864
virtual bool printProblem() const
Prints out problem-specific data to the log stream.
Definition: SIMbase.C:740
LinSolParams * mySolParams
Input parameters for PETSc.
Definition: SIMbase.h:858
virtual bool preprocessC(const IntVec &ignored, bool fixDup, double time0)
Performs some pre-processing tasks on the FE model.
Definition: SIMbase.C:195
void setQuadratureRule(size_t ng, bool redimBuffers=false, bool printQP=false)
Defines the spatial numerical integration scheme to use.
Definition: SIMbase.C:698
AlgEqSystem * myEqSys
The actual linear equation system.
Definition: SIMbase.h:856
int nGlPatches
Number of global patches.
Definition: SIMbase.h:847
virtual bool fieldProjections() const
Returns whether projections must be handled through fields or not.
Definition: SIMbase.C:2592
void dumpEqSys(bool initialBlankLine=false)
Dumps left-hand-side matrix and right-hand-side vector to file.
Definition: SIMbase.C:2689
PropertyVec myProps
Physical property mapping.
Definition: SIMbase.h:804
virtual void shiftGlobalNums(int, int)
Shifts global node and element numbers by constant offsets.
Definition: SIMbase.h:735
std::map< int, int > myGlb2Loc
Global-to-local node number mapping.
Definition: SIMbase.h:850
bool hasTimeDependentDirichlet() const
Checks for time-dependent in-homogeneous Dirichlet conditions.
Definition: SIMbase.C:900
unsigned char nsd
Number of spatial dimensions.
Definition: SIMbase.h:802
const IntegrandBase * getProblem() const
Returns a pointer to the problem-specific data object.
Definition: SIMbase.h:158
virtual unsigned short int getNoParamDim() const =0
Returns the number of parameter dimensions in the model.
bool getElmNodes(IntVec &mnpc, int iel) const
Finds the Matrix of Nodal Point Correspondance for element iel.
Definition: SIMbase.C:848
virtual size_t getNoSpaceDim() const
Returns the number of spatial dimensions in the model.
Definition: SIMbase.h:183
size_t nBouGP
Number of boundary integration points in the whole model.
Definition: SIMbase.h:863
char getNodeType(int inod) const
Returns the type (DOF classification) of the specified global node.
Definition: SIMbase.C:1595
std::map< int, TractionFunc * > TracFuncMap
Traction field container.
Definition: SIMbase.h:796
NormBase * getNormIntegrand() const
Returns a pointer to a norm integrand object for this simulator.
Definition: SIMbase.C:1748
size_t getNoEquations() const
Returns the number of unknowns in the linear equation system.
Definition: SIMbase.C:870
virtual void printSolutionSummary(const Vector &solution, int printSol=0, const char *compName=nullptr, std::streamsize outPrec=0)
Prints a summary of the calculated solution to std::cout.
Definition: SIMbase.C:1512
virtual std::vector< IntVec > getElmConnectivities() const =0
Obtain element-element connectivities.
const IntVec & getMADOF(unsigned char basis, unsigned char nndof) const
Returns an extraordinary MADOF array.
Definition: SIMbase.C:2668
int getLocalPatchIndex(int patchNo) const
Returns the local patch index for the given global patch number.
Definition: SIMbase.C:143
bool assembleForces(const Vector &solution, double t0, RealArray *R, Vector *S=nullptr)
Assembles reaction and interface forces for specified boundaries.
Definition: SIMbase.C:2786
VecFunc * getVecFunc(size_t patch, Property::Type ptype) const
Returns a vector function associated with given patch and property.
Definition: SIMbase.C:564
SystemMatrix * getLHSmatrix(size_t idx=0, bool copy=false) const
Returns current system left-hand-side matrix.
Definition: SIMbase.C:1695
Class administering inter-SIM field dependencies.
Definition: SIMdependency.h:30
std::vector< ASMbase * > PatchVec
Spline patch container.
Definition: SIMdependency.h:33
virtual void registerDependency(SIMdependency *sim, const std::string &name, short int nvc, const PatchVec &patches, char diffBasis=0, int component=1)
Registers a dependency on a field from another SIM object.
Definition: SIMdependency.C:21
int ncv
Number of Arnoldi vectors.
Definition: SIMoptions.h:82
int eig
Eigensolver method (1,...,5)
Definition: SIMoptions.h:80
double shift
Eigenvalue shift.
Definition: SIMoptions.h:83
int nev
Number of eigenvalues/vectors.
Definition: SIMoptions.h:81
ProjectionMethod
Enum defining the available projection methods.
Definition: SIMoptions.h:107
Scalar-valued unary function of a scalar value.
Definition: Function.h:127
Base class for representing a system matrix on different formats.
Definition: SystemMatrix.h:220
Base class for representing a system vector on different formats.
Definition: SystemMatrix.h:32
Vector-valued binary function of a spatial point and normal vector.
Definition: Function.h:292
Simple class for representing a point in 3D space.
Definition: Vec3.h:27
Simple class for representing a point in 3D space and time.
Definition: Vec3.h:209
Vector-valued unary function of a spatial point.
Definition: Function.h:242
A vector class with some added algebraic operations.
Definition: matrix.h:64
Linear algebra scope.
Definition: LinAlgenums.h:19
MatrixType
The available system matrix formats and associated solvers.
Definition: LinAlgenums.h:22
StorageFormat
Enumeration of matrix storage formats.
Definition: LinAlgenums.h:43
@ FLAT
Flat format.
Definition: LinAlgenums.h:45
Struct for storage of data associated with one mode shape.
Definition: SIMbase.h:47
Mode()
Default constructor.
Definition: SIMbase.h:55
double eigVal
Eigenvalue associated with this mode.
Definition: SIMbase.h:49
Vector eqnVec
Eigenvector associated with this mode in equation order.
Definition: SIMbase.h:52
Vector eigVec
Eigenvector associated with this mode.
Definition: SIMbase.h:51
double damping
Modal damping coefficient.
Definition: SIMbase.h:50
bool computeDamping(const SystemMatrix &mat)
Compute modal damping based on the given damping matrix.
Definition: SIMbase.C:2200
bool orthonormalize(const SystemMatrix &mat)
Orthonormalize the eigenvector w.r.t. the given matrix.
Definition: SIMbase.C:2172
int eigNo
Eigenvalue identifier.
Definition: SIMbase.h:48
Struct for representing a distributed physical property.
Definition: Property.h:26
Type
The available property types.
Definition: Property.h:32
A struct with data for system matrix/vector dumps.
Definition: SIMbase.h:826
int count
Internal step counter, dump only when step==count.
Definition: SIMbase.h:831
DumpData()
Default constructor.
Definition: SIMbase.h:835
std::set< int > step
Dump step identifiers.
Definition: SIMbase.h:829
double eps
Zero tolerance for printing small values.
Definition: SIMbase.h:832
std::string fname
File name.
Definition: SIMbase.h:827
LinAlg::StorageFormat format
File format flag.
Definition: SIMbase.h:828
bool expand
If true, dump expanded solution vectors.
Definition: SIMbase.h:830
bool doDump()
Checks if the matrix or vector should be dumped now.
Definition: SIMbase.h:838
Struct representing the time domain.
Definition: TimeDomain.h:23