14 fNTVarVar(0),fRowBlock(),fColBlock(),fColPosition(0),fNumberofColumnBlocks(0),fColumnBlockPosition(0),
15 fColumnBlockNumber(0),fColumnBlockLastUsed(0),fDoubleValues(0),fDoubleValLastUsed(0) {
27 SetBlocks(row,col,nvar,nrowblocks,ncolblocks);
36 out <<
"TPZTransfer<TVar>::Print : ";
38 out <<
"rows : " << this->
Rows() <<
" cols : " << this->
Cols() << endl;
40 out << this->
Rows() <<
' ' << this->
Cols() << endl;
43 for(irb=0 ; irb <
fRowBlock.MaxBlockSize(); irb++) {
45 out <<
"block row number : " << irb << endl;
47 out <<
"position of first column block : " <<
fColPosition[irb] << endl;
52 for(icbcounter=0; icbcounter < numcolbl; icbcounter++) {
56 out <<
"column block counter : " << icbcounter <<
57 " column block number " << icb << endl;
64 out <<
"row position number : " << rownumber <<
65 " column position number " << colnumber << endl;
66 out <<
"block sizes : row : " << rowsize <<
" col " << colsize << endl;
70 loc.
Print(0,out,form);
73 for(ir=0; ir<rowsize; ir++) {
74 for(ic=0; ic<colsize; ic++) {
75 out << (rownumber+ir) <<
' ' << (colnumber+ic) <<
' ' << loc(ir,ic) << endl;
123 cout <<
"TPZTransfer:SetBlocks called for an already defined row = " <<
132 for(;ic < lastic; ic++) {
160 if(colpos == -1 || numcolblocks == -1) {
161 cout <<
"TPZTransfer<TVar>::SetBlockMatrix called for ilegal parameters : " 162 <<
" row = " << row <<
" col = " << col <<
" colpos = " << colpos <<
163 " numcolblocks = " << numcolblocks << endl;
166 int ic = colpos, lastic = colpos+numcolblocks;
167 for(;ic<lastic;ic++) {
171 cout <<
"TPZTransfer<TVar>::SetBlockMatrix column not found for row = " << row <<
172 " col = " << col << endl;
177 if(nblrows != mat.
Rows() || nblcols != mat.
Cols()) {
178 cout <<
"TPZTransfer<TVar>::SetBlockMatrix matrix has incompatible dimensions : " 179 " nblrows = " << nblrows <<
" nblcols = " << nblcols <<
" mat.rows = " <<
180 mat.
Rows() <<
" mat.cols " << mat.
Cols() << endl;
187 fDoubleValLastUsed += nblrows*nblcols;
203 TVar alpha, TVar beta,
int opt)
const{
205 if ((!opt && this->
Cols() != x.
Rows()) || (opt && this->
Rows() != x.
Rows()))
206 this->
Error(
"TPZTransfer<TVar>::MultAdd <matrices with incompatible dimensions>" );
208 this->
Error (
"TPZTransfer<TVar>::MultiplyAdd incompatible dimensions\n");
211 int xcols = x.
Cols();
225 for (
int i=0; i<thisr; i++) {
226 for (
int c=0; c<ncc; c++) {
227 tempcoarse(i,c) = x.
GetVal(iv+i*fNTVarVar,c);
230 MultAddScalar(tempcoarse, tempfine, tempfine, alpha, 0., opt);
231 for (
int i=0; i<thisr; i++) {
232 for (
int c=0; c<ncc; c++) {
233 z(iv+i*fNTVarVar,c) = tempfine.
GetVal(i,c);
242 TVar alpha, TVar beta,
int opt)
const{
244 if ((!opt && this->
Cols() != x.
Rows()) || (opt && this->
Rows() != x.
Rows()))
245 this->
Error(
"TPZTransfer<TVar>::MultAdd <matrices with incompatible dimensions>" );
247 this->
Error (
"TPZTransfer<TVar>::MultiplyAdd incompatible dimensions\n");
250 int xcols = x.
Cols();
253 for (ic = 0; ic < xcols; ic++) {
255 for ( r=0; r < rows; r++) {
258 if(!rowblocksize)
continue;
260 if(colpos == -1)
continue;
262 if(numcolbl == 0)
continue;
264 for( c=0; c < numcolbl; c++) {
266 int colblockpos =
fColBlock.Position(col);
271 aloc.
MultAdd(xloc,zloc,zloc,alpha,1.,opt);
275 for ( r=0; r < rows; r++) {
278 if(!rowblocksize)
continue;
281 if(colpos == -1)
continue;
283 if(numcolbl == 0)
continue;
284 for( c=0; c < numcolbl; c++) {
286 int colblockpos =
fColBlock.Position(col);
291 aloc.
MultAdd(xloc,zloc,zloc,alpha,1.,opt);
307 int nrf = finesol.
Rows();
308 int ncf = finesol.
Cols();
309 int nrc = coarsesol.
Rows();
310 int ncc = coarsesol.
Cols();
311 int thisr = this->
Rows();
312 int thisc = this->
Cols();
316 finesol.
Redim(nrf,ncf);
325 for (
int i=0; i<thisr; i++) {
326 for (
int c=0; c<ncf; c++) {
327 tempcoarse(i,c) = coarsesol.
GetVal(iv+i*fNTVarVar,c);
332 for (
int i=0; i<thisr; i++) {
333 for (
int c=0; c<ncf; c++) {
334 finesol(iv+i*fNTVarVar,c) = tempfine.
GetVal(i,c);
349 int nrf = fine.
Rows();
350 int ncf = fine.
Cols();
351 int nrc = coarse.
Rows();
352 int ncc = coarse.
Cols();
353 int thisr = this->
Rows();
354 int thisc = this->
Cols();
358 coarse.
Redim(nrc,ncf);
367 for (
int i=0; i<thisc; i++) {
368 for (
int c=0; c<ncf; c++) {
369 tempfine(i,c) = fine.
GetVal(iv+i*fNTVarVar,c);
374 for (
int i=0; i<thisr; i++) {
375 for (
int c=0; c<ncf; c++) {
376 coarse(iv+i*fNTVarVar,c) = tempcoarse(i,c);
387 if ((opt==0 && this->
Cols() != A.
Rows()) || (opt ==1 && this->
Rows() != A.
Rows()))
389 this->
Error(
"TPZTransfer<TVar>::Multiply incompatible dimensions" );
TPZManVector< int > fColumnBlockPosition
Vector indicating the starting point of each column block.
int ClassId() const override
Define the class id associated with the class.
void SetBlockMatrix(int row, int col, TPZFMatrix< TVar > &mat)
Sets the row,col block equal to matrix mat if row col was not specified by AddBlockNumbers, an error will be issued and exit.
TPZManVector< TVar > fDoubleValues
Storage space for the matrix blocks.
static int Error(const char *msg, const char *msg2=0)
Returns error messages.
TPZBlock< TVar > fColBlock
Block sizes of the columns.
void AddBlockNumbers(int row, TPZVec< int > &colnumbers)
Will specify the sparsity pattern of row.
MatrixOutputFormat
Defines output format.
void SetBlocks(TPZBlock< TVar > &row, TPZBlock< TVar > &col, int nvar, int nrowblocks, int ncolblocks)
This operation will reset the matrix to zero with no rows defined.
int64_t fRow
Number of rows in matrix.
int fNTVarVar
Number of variables associated with each shape function.
void PrepareZ(const TPZFMatrix< TVar > &y, TPZFMatrix< TVar > &z, const TVar beta, const int opt) const
Is an auxiliar method used by MultiplyAdd.
int fDoubleValLastUsed
Indicates the next free position of fDoubleValues.
virtual void Print(const char *name=NULL, std::ostream &out=std::cout, const MatrixOutputFormat form=EFormatted) const override
It prints the matrix data in a MatrixFormat Rows X Cols.
int NAlloc() const
Returns number of elements allocated for this object.
void Expand(const int64_t newsize)
Expands the allocated storage to fit the newsize parameter.
virtual void Resize(const int64_t newsize, const T &object)
Resizes the vector object.
Contains the TPZTransfer class which implements a rectangular sparse block matrix.
void ExpandColumnVectorEntries(int numcol)
Increases the storage allocated int fColPosition to include numcol more values.
virtual void Resize(const int64_t newsize, const T &object)
Resizes the vector object reallocating the necessary storage, copying the existing objects to the new...
TPZVec< int > fNumberofColumnBlocks
Vector indicating the number of column blocks associated with each row.
void MultiplyScalar(const TPZFMatrix< TVar > &A, TPZFMatrix< TVar > &B, int opt) const
int HasRowDefinition(int row)
Returns 1 if the row is defined (i.e. has column entries)
Contains TPZMatrixclass which implements full matrix (using column major representation).
#define DebugStop()
Returns a message to user put a breakpoint in.
TPZTransfer()
Default constructor.
void TransferResidual(const TPZFMatrix< TVar > &fine, TPZFMatrix< TVar > &coarse)
Will transfer the residual, taking into acount there may be more than one TVar variable.
int64_t Rows() const
Returns number of rows.
TPZVec< int > fColPosition
Vector indicating the starting column block for each row.
Full matrix class. Matrix.
int Redim(const int64_t newRows, const int64_t newCols) override
Redimension a matrix and ZERO your elements.
void TransferSolution(const TPZFMatrix< TVar > &coarsesol, TPZFMatrix< TVar > &finesol)
Will transfer the solution, taking into acount there may be more than one TVar variable.
int64_t fCol
Number of cols in matrix.
int fColumnBlockLastUsed
Indicates the next free position.
Implements block matrices. Matrix utility.
TPZBlock< TVar > fRowBlock
Block sizes of the rows.
int Size(const int block_diagonal) const
Returns block dimension.
void Fill(const T ©, const int64_t from=0, const int64_t numelem=-1)
Will fill the elements of the vector with a copy object.
virtual void MultAdd(const TPZFMatrix< TVar > &x, const TPZFMatrix< TVar > &y, TPZFMatrix< TVar > &z, const TVar alpha=1., const TVar beta=0., const int opt=0) const override
It computes z = beta * y + alpha * opt(this)*x but z and x can not overlap in memory.
int64_t Cols() const
Returns number of cols.
void MultAdd(const TPZFMatrix< TVar > &x, const TPZFMatrix< TVar > &y, TPZFMatrix< TVar > &z, const TVar alpha, const TVar beta, const int opt=0) const override
Multiplies the transfer matrix and puts the result in z.
virtual void Print(std::ostream &out) const
int64_t NElements() const
Returns the number of elements of the vector.
int MaxBlockSize() const
Returns the max number of blocks on diagonal.
void ExpandDoubleValueEntries(int numval)
Increases the storage space available in the fDoubleValues vector to include numval entries...
Implements rectangular matrix which extends a solution vector of the coarse mesh to a solution vector...
TPZManVector< int > fColumnBlockNumber
Vector indicating the number of the column corresponding to the block.
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.
TVar & g(const int64_t row, const int64_t col) const
void MultAddScalar(const TPZFMatrix< TVar > &x, const TPZFMatrix< TVar > &y, TPZFMatrix< TVar > &z, const TVar alpha, const TVar beta, const int opt=0) const
Multiplies the transfer matrix and puts the result in z.
Root matrix class (abstract). Matrix.