NeoPZ
TPZFrontNonSym.h
Go to the documentation of this file.
1 
5 /* Generated by Together */
6 template<class TVar>
7 class TPZEqnArray;
8 
9 #ifndef TPZFRONTNONSYM_H
10 #define TPZFRONTNONSYM_H
11 
12 
13 #include "pzmatrix.h"
14 #include "pzstack.h"
15 #include "pzvec.h"
16 #include "TPZFront.h"
17 
18 #include <math.h>
19 #include <iostream>
20 #include <fstream>
21 #include "TPZStackEqnStorage.h"
22 #include "TPZFileEqnStorage.h"
23 
24 #ifdef USING_BLAS
25 #ifdef USING_MKL
26 #include <mkl.h>
27 #elif MACOSX
28 #include <Accelerate/Accelerate.h>
29 #else
30 #ifdef MACOSX
31 #include <Accelerate/Accelerate.h>
32 #else
33 extern "C"{
34 #include "cblas.h"
35  //#include "g2c.h"
36  //#include "fblaswr.h"
37 };
38 #endif
39 #endif
40 #endif
41 
42 #ifdef USING_ATLAS
43 extern "C"{
44 #include <cblas.h>
45  //#include "g2c.h"
46  //#include "fblaswr.h"
47 };
48 #endif
49 
50 
60 template<class TVar>
61 class TPZFrontNonSym : public TPZFront<TVar> {
62 public:
64  std::string GetMatrixType();
65 
67  static void main();
73  TPZFrontNonSym(int64_t GlobalSize);
74 
75  TPZFrontNonSym(const TPZFrontNonSym &cp);
76 
78  virtual void SetDecomposeType(DecomposeType dectype) override
79  {
80  if (dectype == ELU) {
81  this->fDecomposeType = dectype;
82  }
83  else
84  {
85  DebugStop();
86  }
87  }
88 
89 
97  void DecomposeEquations(int64_t mineq, int64_t maxeq, TPZEqnArray<TVar> & result);
98 
100  void SymbolicDecomposeEquations(int64_t mineq, int64_t maxeq);
101 
103  void SymbolicAddKel(TPZVec < int64_t > & destinationindex);
104 
106  void Compress();
107 
112  void Expand(int largefrontsize);
113 
118  TVar & Element(int64_t i, int64_t j) {
119  return this->fData[this->fMaxFront*j + i];
120  }
121 
126  const TVar & Element(int64_t i, int64_t j) const{
127  return this->fData[this->fMaxFront*j + i];
128  }
129 
135  void AddKel(TPZFMatrix<TVar> &elmat, TPZVec<int64_t> &destinationindex);
136 
143  virtual void AddKel(TPZFMatrix<TVar> &elmat, TPZVec<int64_t> &sourceindex, TPZVec<int64_t> &destinationindex);
144 
146  virtual void ExtractFrontMatrix(TPZFMatrix<TVar> &front) override;
147 
148  public:
149 int ClassId() const override;
150 
151 private:
152 
158  void DecomposeOneEquation(int64_t ieq, TPZEqnArray<TVar> &eqnarray);
159 
165  void FreeGlobal(int64_t global);
170  int Local(int64_t global);
171 
172 public:
174  virtual int64_t NFree() override;
179  void Reset(int64_t GlobalSize=0);
181  void AllocData();
182 
184  void Print(const char *name, std::ostream& out) const;
185  void PrintGlobal(const char *name, std::ostream& out);
186 
187 
188 private:
189 
191  /*# TPZStackEqnStorage lnkTPZStackEqnStorage; */
192 
194  /*# TPZFileEqnStorage lnkTPZFileEqnStorage; */
195 
196 public:
197 
198  virtual void TensorProductIJ(int ithread, typename TPZFront<TVar>::STensorProductMTData *data) override;
199 
200 };
201 
202 template<class TVar>
204  return Hash("TPZFrontNonSym") ^ TPZFront<TVar>::ClassId() << 1;
205 }
206 
207 
208 #endif //TPZFRONTNONSYM_H
void AddKel(TPZFMatrix< TVar > &elmat, TPZVec< int64_t > &destinationindex)
Add a contribution of a stiffness matrix.
TPZVec< TVar > fData
Frontal matrix data.
Definition: TPZFront.h:164
void Expand(int largefrontsize)
Expand the front matrix.
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...
TPZFrontNonSym()
Simple constructor.
Templated vector implementation.
static void main()
virtual void TensorProductIJ(int ithread, typename TPZFront< TVar >::STensorProductMTData *data) override
struct para paralelizar a decomposicao da matriz
Definition: TPZFront.h:178
void Compress()
Compress data structure.
void Reset(int64_t GlobalSize=0)
Resets data structure.
Abstract class implements storage and decomposition process of the frontal matrix. Frontal.
void AllocData()
Allocates data for Front.
DecomposeType fDecomposeType
Used Decomposition method.
Definition: TPZFront.h:172
void DecomposeEquations(int64_t mineq, int64_t maxeq, TPZEqnArray< TVar > &result)
Decompose these equations and put the result in eqnarray. Default decompose method is LU...
const TVar & Element(int64_t i, int64_t j) const
Returns the ith,jth element of the matrix. .
int ClassId() const override
Define the class id associated with the class.
void FreeGlobal(int64_t global)
Sets the global equation as freed, allowing the space used by this equation to be used by future ass...
#define DebugStop()
Returns a message to user put a breakpoint in.
Definition: pzerror.h:20
virtual void ExtractFrontMatrix(TPZFMatrix< TVar > &front) override
Extract the front matrix.
Definition: pzmatrix.h:52
~TPZFrontNonSym()
Simple destructor.
Full matrix class. Matrix.
Definition: pzfmatrix.h:32
Contains the TPZFileEqnStorage class which implements an equation array and stores the EqnArrays...
A simple stack.
int32_t Hash(std::string str)
Definition: TPZHash.cpp:10
TVar & Element(int64_t i, int64_t j)
Returns the ith,jth element of the matrix. .
Contains TPZMatrix<TVar>class, root matrix class.
virtual void SetDecomposeType(DecomposeType dectype) override
Set the decomposition type.
void SymbolicAddKel(TPZVec< int64_t > &destinationindex)
Add a contribution of a stiffness matrix using the indexes to compute the frontwidth.
Contains the TPZStackEqnStorage class responsible for storing arrays of equations.
int fMaxFront
Maximum size of the front.
Definition: TPZFront.h:136
void SymbolicDecomposeEquations(int64_t mineq, int64_t maxeq)
Decompose these equations in a symbolic way and store freed indexes in fFree.
virtual int64_t NFree() override
Returns the number of free equations.
void DecomposeOneEquation(int64_t ieq, TPZEqnArray< TVar > &eqnarray)
Decomposes ieq equation and add the result to EqnArray.
void Print(const char *name, std::ostream &out) const
It prints TPZFront data.
int Local(int64_t global)
Returns a local index corresponding to a global equation number.
std::string GetMatrixType()
Type of matrix.
void PrintGlobal(const char *name, std::ostream &out)
Abstract class implements storage and decomposition process of the frontal matrix involving non-simet...
int ClassId() const override
Define the class id associated with the class.
Definition: TPZFront.h:317
DecomposeType
Defines decomposition type for any matrix classes.
Definition: pzmatrix.h:52