IFEM  90A354
SAM.h
Go to the documentation of this file.
1 // $Id$
2 //==============================================================================
12 //==============================================================================
13 
14 #ifndef _SAM_H
15 #define _SAM_H
16 
17 #include "MatVec.h"
18 #include <set>
19 #include <map>
20 
21 class SystemMatrix;
22 class SystemVector;
23 
24 using IntVec = std::vector<int>;
25 using IntSet = std::set<int>;
26 
27 
38 class SAM
39 {
40 public:
42  SAM();
44  virtual ~SAM();
45 
47  void print(std::ostream& os) const;
48 
50  template<class T> void printVector(T& os, const RealArray& vec,
51  const char* heading = nullptr) const
52  {
53  if (heading)
54  os <<"\n"<< heading <<":";
55  for (int inod = 1; inod <= nnod; inod++)
56  for (int jdof = madof[inod-1]; jdof < madof[inod]; jdof++)
57  if (utl::trunc(vec[jdof-1]) != Real(0))
58  {
59  os <<"\nNode "<< inod <<":";
60  for (int idof = madof[inod-1]; idof < madof[inod]; idof++)
61  os <<" "<< utl::trunc(vec[idof-1]);
62  break;
63  }
64  os << std::endl;
65  }
66 
68  int getNoElms() const { return nel; }
71  int getNoNodes(char dofType = 'A') const;
73  int getNoDOFs() const { return ndof; }
75  int getNoSpecifiedDOFs() const { return nspdof; }
77  int getNoConstraints() const { return nceq; }
79  int getNoEquations() const { return neq; }
83  IntSet getEquations(char dType, int dof = 0) const;
85  const int* getMADOF() const { return madof; }
87  const int* getMEQN() const { return meqn; }
88 
92  void getDofCouplings(IntVec& irow, IntVec& jcol) const;
94  void getDofCouplings(std::vector<IntSet>& dofc) const;
95 
106  bool assembleSystem(SystemMatrix& sysK, SystemVector& sysRHS,
107  const Matrix& eK, int iel = 0,
108  RealArray* reactionForces = nullptr) const;
109 
115  bool assembleSystem(SystemMatrix& sysM, const Matrix& eM, int iel = 0) const;
116 
126  bool assembleSystem(SystemVector& sysRHS,
127  const Matrix& eK, int iel = 0,
128  RealArray* reactionForces = nullptr) const;
129 
136  bool assembleSystem(SystemVector& sysRHS,
137  const RealArray& eS, int iel = 0,
138  RealArray* reactionForces = nullptr) const;
139 
146  bool assembleSystem(SystemVector& sysRHS,
147  const Real* nS, int inod = 0,
148  RealArray* reactionForces = nullptr) const;
149 
155  bool assembleSystem(SystemVector& sysRHS, Real S,
156  const std::pair<int,int>& dof) const;
157 
161  void addToRHS(SystemVector& sysRHS, const RealArray& S) const;
162 
166  bool getElmNodes(IntVec& mnpc, int iel) const;
172  bool getElmEqns(IntVec& meen, int iel, size_t nedof = 0) const;
176  void getUniqueEqns(IntSet& meen, int iel) const;
177 
181  bool getNodeEqns(IntVec& mnen, int inod) const;
182 
185  char getNodeType(int inod) const
186  { return inod-- > 0 && inod < (int)nodeType.size() ? nodeType[inod] : ' '; }
187 
190  std::pair<int,int> getNodeDOFs(int inod) const;
191 
197  std::pair<int,int> getNodeAndLocalDof(int idof, bool eqno = false) const;
198 
203  int getEquation(int inod, int ldof) const;
204 
210  Real getDofVal(const RealArray& dofVec, int inod, int ldof) const;
211 
223  virtual bool expandSolution(const SystemVector& solVec,
224  Vector& dofVec, Real scaleSD = Real(1)) const;
225 
232  bool expandVector(const Vector& solVec, Vector& dofVec) const;
233 
240  bool getSolVec(RealArray& solVec, const RealArray& dofVec) const;
241 
246  bool applyDirichlet(Vector& dofVec) const;
247 
252  virtual Real dot(const Vector& x, const Vector& y, char dofType = 'D') const;
256  Real norm2(const Vector& x, char dofType = 'D') const
257  { return sqrt(this->dot(x,x,dofType)); }
261  virtual Real normL2(const Vector& x, char dofType = 'D') const;
266  virtual Real normInf(const Vector& x, size_t& comp, char dofType = 'D') const;
267 
271  virtual Real normReact(const RealArray& u, const RealArray& rf) const;
276  Real getReaction(int dir, const RealArray& rf,
277  const IntVec* nodes = nullptr) const;
281  bool haveReaction(int dir, const IntVec* nodes = nullptr) const;
287  bool getNodalReactions(int inod, const RealArray& rf, RealArray& nrf) const;
288 
290  virtual bool merge(const SAM*, const std::map<int,int>*) { return false; }
291 
292 protected:
294  bool initSystemEquations();
295 
304  void assembleRHS(Real* RHS, Real value, int ieq) const;
305 
310  void assembleReactions(RealArray& rf, const RealArray& eS, int iel) const;
311 
316  bool expandVector(const Real* solVec, Vector& dofVec, Real scaleSD) const;
317 
319  void printStatusCodes(std::ostream& os) const;
321  void printCEQ(std::ostream& os, int iceq) const;
323  void printMEQN(std::ostream& os) const;
324 
325 private:
326  int mpar[50];
327 
328 protected:
329  // The following parameters are pointers to specific locations in MPAR
330  int& nnod;
331  int& nel;
332  int& ndof;
333  int& nspdof;
334  int& nceq;
335  int& neq;
336  int& nmmnpc;
337  int& nmmceq;
338 
339  // The standard %SAM arrays (see K. Bell's reports for detailed explanation).
340  // We are using plane C-pointers for these items such that they more easily
341  // can be passed directly as arguments to FORTRAN subroutines.
342  int* mpmnpc;
343  int* mmnpc;
344  int* madof;
345  int* msc;
346  int* mpmceq;
347  int* mmceq;
349  int* minex;
350  int* meqn;
351 
352  std::vector<char> nodeType;
353  std::vector<char> dof_type;
354 
355  friend class DenseMatrix;
356  friend class SPRMatrix;
357  friend class SparseMatrix;
358  friend class DiagMatrix;
359  friend class PETScMatrix;
360 };
361 
362 #endif
std::vector< int > IntVec
General integer vector.
Definition: ASMbase.h:25
std::set< int > IntSet
General integer set.
Definition: ASMunstruct.h:22
#define Real
The floating point type to use.
Definition: ImmersedBoundaries.h:18
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.
Class for representing a dense system matrix.
Definition: DenseMatrix.h:27
Class for representing a diagonal system matrix.
Definition: DiagMatrix.h:25
Class for representing the system matrix in PETSc format.
Definition: PETScMatrix.h:172
This class contains data and functions for the assembly of FE matrices.
Definition: SAM.h:39
void addToRHS(SystemVector &sysRHS, const RealArray &S) const
Adds a global load vector into the system load vector.
Definition: SAM.C:593
virtual Real dot(const Vector &x, const Vector &y, char dofType='D') const
Computes the dot-product of two vectors of length NDOF.
Definition: SAM.C:842
virtual bool expandSolution(const SystemVector &solVec, Vector &dofVec, Real scaleSD=Real(1)) const
Expands a solution vector from equation-ordering to DOF-ordering.
Definition: SAM.C:768
void printVector(T &os, const RealArray &vec, const char *heading=nullptr) const
Prints out a nodal DOF vector to the given stream.
Definition: SAM.h:50
void assembleReactions(RealArray &rf, const RealArray &eS, int iel) const
Assembles reaction forces for the fixed and prescribed DOFs.
Definition: SAM.C:620
int & nnod
Number of nodes.
Definition: SAM.h:330
void printCEQ(std::ostream &os, int iceq) const
Prints out a constraint equation to the given stream.
Definition: SAM.C:127
bool assembleSystem(SystemMatrix &sysK, SystemVector &sysRHS, const Matrix &eK, int iel=0, RealArray *reactionForces=nullptr) const
Adds an element stiffness matrix into the system stiffness matrix.
Definition: SAM.C:430
SAM()
The constructor initializes an empty object.
Definition: SAM.C:66
int getNoSpecifiedDOFs() const
Returns the number of specified DOFs in the model.
Definition: SAM.h:75
Real norm2(const Vector &x, char dofType='D') const
Computes the l2-norm of a vector of length NDOF.
Definition: SAM.h:256
int mpar[50]
Matrix of parameters.
Definition: SAM.h:326
std::pair< int, int > getNodeAndLocalDof(int idof, bool eqno=false) const
Returns the internal node number and local index for a global DOF.
Definition: SAM.C:101
int * madof
Matrix of accumulated DOFs.
Definition: SAM.h:344
Real * ttcc
Table of tables of constraint equation coefficients.
Definition: SAM.h:348
int & nspdof
Number of specified DOFs.
Definition: SAM.h:333
bool getNodalReactions(int inod, const RealArray &rf, RealArray &nrf) const
Returns a vector of reaction forces for a given node.
Definition: SAM.C:976
virtual ~SAM()
The destructor frees the dynamically allocated arrays.
Definition: SAM.C:86
int * minex
Matrix of internal to external node numbers.
Definition: SAM.h:349
int & neq
Number of system equations.
Definition: SAM.h:335
virtual Real normInf(const Vector &x, size_t &comp, char dofType='D') const
Computes the L_infinity-norm of a vector of length NDOF.
Definition: SAM.C:883
bool getElmNodes(IntVec &mnpc, int iel) const
Finds the matrix of nodal point correspondance for an element.
Definition: SAM.C:635
void printStatusCodes(std::ostream &os) const
Prints out the DOF status codes to the given stream.
Definition: SAM.C:208
std::vector< char > nodeType
Nodal DOF classification.
Definition: SAM.h:352
int & nel
Number of elements.
Definition: SAM.h:331
int getEquation(int inod, int ldof) const
Finds the equation number corresponding to a local nodal DOF.
Definition: SAM.C:746
int * msc
Matrix of status codes.
Definition: SAM.h:345
std::vector< char > dof_type
Individual DOF classification.
Definition: SAM.h:353
int getNoElms() const
Returns the number of elements in the model.
Definition: SAM.h:68
IntSet getEquations(char dType, int dof=0) const
Returns the equations numbers for a given DOF type.
Definition: SAM.C:1002
int * mmnpc
Matrix of matrices of nodal point correspondances.
Definition: SAM.h:343
int & nmmceq
Number of elements in MMCEQ.
Definition: SAM.h:337
int * mpmceq
Matrix of pointers to MCEQs in MMCEQ.
Definition: SAM.h:346
char getNodeType(int inod) const
Returns the DOF classification of a given node.
Definition: SAM.h:185
Real getReaction(int dir, const RealArray &rf, const IntVec *nodes=nullptr) const
Returns the total reaction force in the given coordinate direction.
Definition: SAM.C:956
void getDofCouplings(IntVec &irow, IntVec &jcol) const
Computes the sparse structure (DOF couplings) in the system matrix.
Definition: SAM.C:354
const int * getMADOF() const
Returns the Matrix of Accumulated DOFs.
Definition: SAM.h:85
bool getElmEqns(IntVec &meen, int iel, size_t nedof=0) const
Finds the matrix of equation numbers for an element.
Definition: SAM.C:655
bool getNodeEqns(IntVec &mnen, int inod) const
Finds the matrix of equation numbers for a node.
Definition: SAM.C:721
int getNoNodes(char dofType='A') const
Returns the number of FE nodes in the model.
Definition: SAM.C:339
void getUniqueEqns(IntSet &meen, int iel) const
Finds the set of unique equation numbers for an element.
Definition: SAM.C:699
int & nmmnpc
Number of elements in MMNPC.
Definition: SAM.h:336
int getNoEquations() const
Returns the number of equations (free DOFs) in the model.
Definition: SAM.h:79
int * mpmnpc
Matrix of pointers to MNPCs in MMNPC.
Definition: SAM.h:342
bool haveReaction(int dir, const IntVec *nodes=nullptr) const
Checks for total reaction force in the given coordinate direction.
Definition: SAM.C:943
Real getDofVal(const RealArray &dofVec, int inod, int ldof) const
Returns specified solution component from a global solution vector.
Definition: SAM.C:758
std::pair< int, int > getNodeDOFs(int inod) const
Returns the first and last DOFs for a node.
Definition: SAM.C:738
bool expandVector(const Vector &solVec, Vector &dofVec) const
Expands a solution vector from equation-ordering to DOF-ordering.
Definition: SAM.C:777
int * mmceq
Matrix of matrices of constraint equation definitions.
Definition: SAM.h:347
int getNoConstraints() const
Returns the total number of constraint equations in the model.
Definition: SAM.h:77
void assembleRHS(Real *RHS, Real value, int ieq) const
Adds a scalar value into a system right hand-side vector.
Definition: SAM.C:605
bool getSolVec(RealArray &solVec, const RealArray &dofVec) const
Extracts equation-ordered solution vector from DOF-ordered vector.
Definition: SAM.C:812
bool initSystemEquations()
Initializes the DOF-to-equation connectivity array MEQN.
Definition: SAM.C:225
int getNoDOFs() const
Returns the total number of DOFs in the model.
Definition: SAM.h:73
void printMEQN(std::ostream &os) const
Prints out the equation number mapping to the given stream.
Definition: SAM.C:180
void print(std::ostream &os) const
Prints out the key data to the given stream.
Definition: SAM.C:144
bool applyDirichlet(Vector &dofVec) const
Applies the non-homogenous Dirichlet BCs to the given vector.
Definition: SAM.C:825
int & nceq
Number of constraint equations.
Definition: SAM.h:334
const int * getMEQN() const
Returns the Matrix of EQuation Numbers.
Definition: SAM.h:87
virtual Real normL2(const Vector &x, char dofType='D') const
Computes the L2-norm of a vector of length NDOF.
Definition: SAM.C:860
virtual Real normReact(const RealArray &u, const RealArray &rf) const
Computes the energy norm contributions from nodal reaction forces.
Definition: SAM.C:923
int * meqn
Matrix of equation numbers.
Definition: SAM.h:350
virtual bool merge(const SAM *, const std::map< int, int > *)
Merges the assembly data from another SIM object with this.
Definition: SAM.h:290
int & ndof
Number of DOFs.
Definition: SAM.h:332
Class for representing the system matrix on the SPR-format.
Definition: SPRMatrix.h:37
Class for representing a system matrix on an unstructured sparse form.
Definition: SparseMatrix.h:38
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
A vector class with some added algebraic operations.
Definition: matrix.h:64
T trunc(T v)
Truncate a value to zero when it is less than a given threshold.
Definition: matrix.h:1765