NeoPZ
tpzdohrassembly.cpp
Go to the documentation of this file.
1 
8 #include "tpzdohrassembly.h"
9 #include "pzfmatrix.h"
10 #include "pzvec.h"
11 #include "pzlog.h"
12 
13 #ifdef LOG4CXX
14 static LoggerPtr logger(Logger::getLogger("substruct.dohrassembly"));
15 #endif
16 
17 // sum the values in the local matrix into the global matrix
18 template<class TVar>
19 void TPZDohrAssembly<TVar>::Assemble(int isub, const TPZFMatrix<TVar> &local, TPZFMatrix<TVar> &global) const
20 {
21  TPZVec<int> &avec = fFineEqs[isub];
22  int neq = avec.NElements();
23  int ncols = local.Cols();
24  int ieq;
25  for (int ic=0; ic<ncols; ic++)
26  {
27  for(ieq=0; ieq<neq; ieq++)
28  {
29  global(avec[ieq],ic) += local.GetVal(ieq,ic);
30  }
31  }
32 #ifdef LOG4CXX
33  {
34  std::stringstream sout;
35  sout << "Assembling destination indices " << avec << std::endl;
36  local.Print("Input vector",sout);
37  global.Print("Resulting vector",sout);
38  if (logger->isDebugEnabled())
39  {
40  LOGPZ_DEBUG(logger, sout.str());
41  }
42  }
43 #endif
44 }
45 
46 // extract the values from the global matrix into the local matrix
47 template<class TVar>
48 void TPZDohrAssembly<TVar>::Extract(int isub, const TPZFMatrix<TVar> &global, TPZFMatrix<TVar> &local) const
49 {
50  TPZVec<int> &avec = fFineEqs[isub];
51  int neq = avec.NElements();
52  int ncols = global.Cols();
53  int ieq;
54  local.Resize(neq,ncols);
55  for (int ic=0; ic<ncols; ic++)
56  {
57  for(ieq=0; ieq<neq; ieq++)
58  {
59  local(ieq,ic) = global.GetVal(avec[ieq],ic);
60  }
61  }
62 #ifdef LOG4CXX
63  {
64  std::stringstream sout;
65  sout << "sub structure " << isub << " Extracting destination indices " << avec << std::endl;
66  local.Print("extracted vector",sout);
67  global.Print("Global vector",sout);
68  if (logger->isDebugEnabled())
69  {
70  LOGPZ_DEBUG(logger, sout.str());
71  }
72  }
73 #endif
74 }
75 
76 // sum the values in the local matrix into the global matrix
77 template<class TVar>
78 void TPZDohrAssembly<TVar>::AssembleCoarse(int isub, const TPZFMatrix<TVar> &local, TPZFMatrix<TVar> &global) const
79 {
80  TPZVec<int> &avec = fCoarseEqs[isub];
81  int neq = avec.NElements();
82  int ieq;
83  int ncols = local.Cols();
84  for (int ic=0; ic<ncols; ic++)
85  {
86  for(ieq=0; ieq<neq; ieq++)
87  {
88  global(avec[ieq],ic) += local.GetVal(ieq,ic);
89  }
90  }
91 #ifdef LOG4CXX
92  {
93  std::stringstream sout;
94  sout << "Assembling destination indices " << avec << std::endl;
95  local.Print("Input vector",sout);
96  global.Print("Resulting vector",sout);
97  if (logger->isDebugEnabled())
98  {
99  LOGPZ_DEBUG(logger, sout.str());
100  }
101  }
102 #endif
103 }
104 
105 // extract the values from the global matrix into the local matrix
106 template<class TVar>
107 void TPZDohrAssembly<TVar>::ExtractCoarse(int isub, const TPZFMatrix<TVar> &global, TPZFMatrix<TVar> &local) const
108 {
109  TPZVec<int> &avec = fCoarseEqs[isub];
110  int neq = avec.NElements();
111  int ieq;
112  int ncols = global.Cols();
113  local.Resize(neq,ncols);
114  for (int ic=0; ic<ncols; ic++)
115  {
116  for(ieq=0; ieq<neq; ieq++)
117  {
118  local(ieq,ic) = global.GetVal(avec[ieq],ic);
119  }
120  }
121 }
122 
123 template<class TVar>
125  return Hash("TPZDohrAssembly") ^ ClassIdOrHash<TVar>() << 1;
126 }
127 
129 template<class TVar>
130 void TPZDohrAssembly<TVar>::Write(TPZStream &buf, int withclassid) const
131 {
132  int nfine = fFineEqs.size();
133  buf.Write(&nfine,1);
134  for (int f=0; f<nfine; f++) {
135  buf.Write( fFineEqs[f]);
136  }
137  int ncoarse = fCoarseEqs.size();
138  buf.Write(&ncoarse,1);
139  for (int nc=0; nc<ncoarse; nc++) {
140  buf.Write( fCoarseEqs[nc]);
141  }
142 }
143 
145 template<class TVar>
146 void TPZDohrAssembly<TVar>::Read(TPZStream &buf, void *context)
147 {
148  int nfine;
149  buf.Read(&nfine);
150  fFineEqs.resize(nfine);
151  for (int f=0; f<nfine; f++) {
152  buf.Read( fFineEqs[f]);
153  }
154  int ncoarse;
155  buf.Read(&ncoarse);
156  fCoarseEqs.resize(ncoarse);
157  for (int nc=0; nc<ncoarse; nc++) {
158  buf.Read( fCoarseEqs[nc]);
159  }
160 }
161 
162 template class TPZDohrAssembly<float>;
163 template class TPZDohrAssembly<double>;
164 template class TPZDohrAssembly<long double>;
165 
166 template class TPZDohrAssembly<std::complex<float> >;
167 template class TPZDohrAssembly<std::complex<double> >;
void Extract(int isub, const TPZFMatrix< TVar > &global, TPZFMatrix< TVar > &local) const
Extract the values from the global matrix into the local matrix.
void Assemble(int isub, const TPZFMatrix< TVar > &local, TPZFMatrix< TVar > &global) const
Sum the values in the local matrix into the global matrix.
Assembling using Dohrmann algorithm. Sub structure.
int ClassId() const override
Define the class id associated with the class.
void Read(TPZStream &buf, void *context) override
method for reading the object for a stream
Contains definitions to LOGPZ_DEBUG, LOGPZ_INFO, LOGPZ_WARN, LOGPZ_ERROR and LOGPZ_FATAL, and the implementation of the inline InitializePZLOG(string) function using log4cxx library or not. It must to be called out of "#ifdef LOG4CXX" scope.
Templated vector implementation.
void AssembleCoarse(int isub, const TPZFMatrix< TVar > &local, TPZFMatrix< TVar > &global) const
Sum the values in the local matrix into the global matrix.
Contains the TPZDohrAssembly class which implements assembling using Dohrmann algorithm.
void ExtractCoarse(int isub, const TPZFMatrix< TVar > &global, TPZFMatrix< TVar > &local) const
Extract the values from the global matrix into the local matrix.
f
Definition: test.py:287
virtual void Write(const bool val)
Definition: TPZStream.cpp:8
Contains TPZMatrixclass which implements full matrix (using column major representation).
#define LOGPZ_DEBUG(A, B)
Define log for debug info.
Definition: pzlog.h:87
void Write(TPZStream &buf, int withclassid) const override
method for streaming the object to a stream
Full matrix class. Matrix.
Definition: pzfmatrix.h:32
int32_t Hash(std::string str)
Definition: TPZHash.cpp:10
int64_t Cols() const
Returns number of cols.
Definition: pzmatrix.h:809
virtual void Print(std::ostream &out) const
Definition: pzmatrix.h:253
int Resize(const int64_t newRows, const int64_t wCols) override
Redimension a matrix, but maintain your elements.
Definition: pzfmatrix.cpp:1016
Defines the interface for saving and reading data. Persistency.
Definition: TPZStream.h:50
int64_t NElements() const
Returns the number of elements of the vector.
Definition: pzvec.h:190
const TVar & GetVal(const int64_t row, const int64_t col) const override
Get values without bounds checking This method is faster than "Get" if DEBUG is defined.
Definition: pzfmatrix.h:566
virtual void Read(bool &val)
Definition: TPZStream.cpp:91