NeoPZ
TPZFrontSym.h
Go to the documentation of this file.
1 
5 template<class TVar>
6 class TPZEqnArray;
7 
8 #ifndef TPZFRONTSYM_H
9 #define TPZFRONTSYM_H
10 
11 #ifdef USING_ATLAS
12 extern "C"{
13 #include <cblas.h>
14 };
15 #endif
16 
17 #include <iostream>
18 #include <math.h>
19 #include <stdlib.h>
20 #include <stdio.h>
21 #include <fstream>
22 
23 #include "pzmatrix.h"
24 #include "pzstack.h"
25 #include "pzvec.h"
26 #include "TPZFront.h"
27 #include "TPZFileEqnStorage.h"
28 #include "TPZStackEqnStorage.h"
29 
39 template <class TVar>
40 class TPZFrontSym : public TPZFront<TVar> {
41 public:
43  std::string GetMatrixType();
44 
46  static void main();
48  ~TPZFrontSym();
50  TPZFrontSym();
51 
53  TPZFront<TVar>(cp)
54  {
55  }
57  TPZFrontSym(int64_t GlobalSize);
58 
60  virtual void SetDecomposeType(DecomposeType dectype) override
61  {
62  if (dectype == ECholesky || dectype == ELDLt) {
63  this->fDecomposeType = dectype;
64  }
65  else
66  {
67  DebugStop();
68  }
69  }
70 
71 
79  void DecomposeEquations(int64_t mineq, int64_t maxeq, TPZEqnArray<TVar> & result);
80 
86  void SymbolicDecomposeEquations(int64_t mineq, int64_t maxeq);
87 
89  void SymbolicAddKel(TPZVec < int64_t > & destinationindex);
90 
92  void Compress();
93 
95  void Expand(int largefrontsize);
96 
98  TVar & Element(int64_t i, int64_t j){
99  if(i>j){
100  int64_t i_temp=i;
101  i=j;
102  j=i_temp;
103  }
104  return this->fData[(j*(j+1))/2+i];
105  }
107  const TVar & Element(int64_t i, int64_t j) const {
108  if(i>j){
109  int64_t i_temp=i;
110  i=j;
111  j=i_temp;
112  }
113  return this->fData[(j*(j+1))/2+i];
114  }
116  void AddKel(TPZFMatrix<TVar> &elmat, TPZVec<int64_t> &destinationindex);
117 
119  virtual void AddKel(TPZFMatrix<TVar> &elmat, TPZVec<int64_t> &sourceindex, TPZVec<int64_t> &destinationindex);
120 
122  virtual void ExtractFrontMatrix(TPZFMatrix<TVar> &front) override;
123 
124  public:
125 int ClassId() const override;
126 
127 private:
128 
129  TVar & Element4JGreatEqualI(int64_t i, int64_t j){
130 #ifdef PZDEBUG
131  if(i>j){
132  DebugStop();
133  }
134 #endif
135  return this->fData[(j*(j+1))/2+i];
136  }
137 
138 
144  void DecomposeOneEquation(int64_t ieq, TPZEqnArray<TVar> &eqnarray);
145 
152  void FreeGlobal(int64_t global);
154  int Local(int64_t global);
155 public:
157  virtual int64_t NFree() override;
159  void Reset(int64_t GlobalSize=0);
161  void AllocData();
162 
164  void Print(const char *name, std::ostream& out=std::cout) const;
165  void PrintGlobal(const char *name, std::ostream& out = std::cout);
166 
169 
171  virtual void TensorProductIJ(int ithread, typename TPZFront<TVar>::STensorProductMTData *data) override;
172 
174  /*# TPZFileEqnStorage lnkTPZFileEqnStorage; */
175 
177  /*# TPZStackEqnStorage lnkTPZStackEqnStorage; */
178 };
179 
180 #endif //TPZFRONTSYM_H
void Expand(int largefrontsize)
Expand the front matrix.
void Reset(int64_t GlobalSize=0)
Resets data structure.
Definition: pzmatrix.h:52
int Local(int64_t global)
return a local index corresponding to a global equation number
TPZVec< TVar > fData
Frontal matrix data.
Definition: TPZFront.h:164
void SymbolicDecomposeEquations(int64_t mineq, int64_t maxeq)
Decompose these equations in a symbolic way and store freed indexes in fFree.
It is an equation array, generally in its decomposed form. Frontal.
Definition: tpzeqnarray.h:36
Contains the TPZFront class which implements decomposition process of the frontal matrix...
Templated vector implementation.
void AddKel(TPZFMatrix< TVar > &elmat, TPZVec< int64_t > &destinationindex)
Add a contribution of a stiffness matrix.
struct para paralelizar a decomposicao da matriz
Definition: TPZFront.h:178
void SymbolicAddKel(TPZVec< int64_t > &destinationindex)
Add a contribution of a stiffness matrix using the indexes to compute the frontwidth.
~TPZFrontSym()
Simple destructor.
Abstract class implements storage and decomposition process of the frontal matrix involving simmetry ...
const TVar & Element(int64_t i, int64_t j) const
Returns ith, jth element of matrix. .
Definition: TPZFrontSym.h:107
Abstract class implements storage and decomposition process of the frontal matrix. Frontal.
virtual void TensorProductIJ(int ithread, typename TPZFront< TVar >::STensorProductMTData *data) override
Does the tensor product betweem two vectors in the positions dependent of ithread.
void DecomposeEquations(int64_t mineq, int64_t maxeq, TPZEqnArray< TVar > &result)
Decompose these equations and put the result in eqnarray Default decompose method is Cholesky...
DecomposeType fDecomposeType
Used Decomposition method.
Definition: TPZFront.h:172
#define DebugStop()
Returns a message to user put a breakpoint in.
Definition: pzerror.h:20
TVar & Element4JGreatEqualI(int64_t i, int64_t j)
Definition: TPZFrontSym.h:129
TVar & Element(int64_t i, int64_t j)
Returns ith, jth element of matrix. .
Definition: TPZFrontSym.h:98
virtual void SetDecomposeType(DecomposeType dectype) override
Set the decomposition type.
Definition: TPZFrontSym.h:60
TPZFrontSym(const TPZFrontSym< TVar > &cp)
Definition: TPZFrontSym.h:52
void DecomposeOneEquation(int64_t ieq, TPZEqnArray< TVar > &eqnarray)
Decomposes ieq equation and add the result to EqnArray.
DecomposeType GetDecomposeType() const
Returns decomposition type.
Definition: TPZFrontSym.cpp:28
Full matrix class. Matrix.
Definition: pzfmatrix.h:32
Contains the TPZFileEqnStorage class which implements an equation array and stores the EqnArrays...
A simple stack.
Contains TPZMatrix<TVar>class, root matrix class.
virtual void ExtractFrontMatrix(TPZFMatrix< TVar > &front) override
Reorders the elements of the frontmatrix into the full matrix.
virtual int64_t NFree() override
Returns the number of free equations.
int ClassId() const override
Define the class id associated with the class.
Contains the TPZStackEqnStorage class responsible for storing arrays of equations.
void Print(const char *name, std::ostream &out=std::cout) const
Prints TPZFront data.
Definition: TPZFrontSym.cpp:51
void FreeGlobal(int64_t global)
Sets the global equation as freed, allowing the space.
void AllocData()
Allocates data for Front.
Definition: TPZFrontSym.cpp:91
std::string GetMatrixType()
Returns its type.
void PrintGlobal(const char *name, std::ostream &out=std::cout)
Definition: TPZFrontSym.cpp:33
static void main()
TPZFrontSym()
Simple constructor.
DecomposeType
Defines decomposition type for any matrix classes.
Definition: pzmatrix.h:52
void Compress()
Compress data structure.