Teko  Version of the Day
Teko_Utilities.cpp
1 /*
2 // @HEADER
3 //
4 // ***********************************************************************
5 //
6 // Teko: A package for block and physics based preconditioning
7 // Copyright 2010 Sandia Corporation
8 //
9 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
10 // the U.S. Government retains certain rights in this software.
11 //
12 // Redistribution and use in source and binary forms, with or without
13 // modification, are permitted provided that the following conditions are
14 // met:
15 //
16 // 1. Redistributions of source code must retain the above copyright
17 // notice, this list of conditions and the following disclaimer.
18 //
19 // 2. Redistributions in binary form must reproduce the above copyright
20 // notice, this list of conditions and the following disclaimer in the
21 // documentation and/or other materials provided with the distribution.
22 //
23 // 3. Neither the name of the Corporation nor the names of the
24 // contributors may be used to endorse or promote products derived from
25 // this software without specific prior written permission.
26 //
27 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
28 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
29 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
30 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
31 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
32 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
33 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
34 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
35 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
36 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
37 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
38 //
39 // Questions? Contact Eric C. Cyr (eccyr@sandia.gov)
40 //
41 // ***********************************************************************
42 //
43 // @HEADER
44 
45 */
46 
47 #include "Teko_Config.h"
48 #include "Teko_Utilities.hpp"
49 
50 // Thyra includes
51 #include "Thyra_MultiVectorStdOps.hpp"
52 #include "Thyra_ZeroLinearOpBase.hpp"
53 #include "Thyra_DefaultDiagonalLinearOp.hpp"
54 #include "Thyra_DefaultAddedLinearOp.hpp"
55 #include "Thyra_EpetraExtDiagScaledMatProdTransformer.hpp"
56 #include "Thyra_EpetraExtDiagScalingTransformer.hpp"
57 #include "Thyra_EpetraExtAddTransformer.hpp"
58 #include "Thyra_DefaultScaledAdjointLinearOp.hpp"
59 #include "Thyra_DefaultMultipliedLinearOp.hpp"
60 #include "Thyra_DefaultZeroLinearOp.hpp"
61 #include "Thyra_DefaultProductMultiVector.hpp"
62 #include "Thyra_DefaultProductVectorSpace.hpp"
63 #include "Thyra_MultiVectorStdOps.hpp"
64 #include "Thyra_VectorStdOps.hpp"
65 #include "Thyra_SpmdVectorBase.hpp"
66 #include "Thyra_get_Epetra_Operator.hpp"
67 #include "Thyra_EpetraThyraWrappers.hpp"
68 #include "Thyra_EpetraOperatorWrapper.hpp"
69 #include "Thyra_EpetraLinearOp.hpp"
70 
71 // Teuchos includes
72 #include "Teuchos_Array.hpp"
73 
74 // Epetra includes
75 #include "Epetra_Operator.h"
76 #include "Epetra_CrsGraph.h"
77 #include "Epetra_CrsMatrix.h"
78 #include "Epetra_Vector.h"
79 #include "Epetra_Map.h"
80 
81 #include "EpetraExt_Transpose_RowMatrix.h"
82 #include "EpetraExt_MatrixMatrix.h"
83 
84 // Anasazi includes
85 #include "AnasaziBasicEigenproblem.hpp"
86 #include "AnasaziThyraAdapter.hpp"
87 #include "AnasaziBlockKrylovSchurSolMgr.hpp"
88 #include "AnasaziBlockKrylovSchur.hpp"
89 #include "AnasaziStatusTestMaxIters.hpp"
90 
91 // Isorropia includes
92 #ifdef Teko_ENABLE_Isorropia
93 #include "Isorropia_EpetraProber.hpp"
94 #endif
95 
96 // Teko includes
97 #include "Teko_EpetraHelpers.hpp"
98 #include "Teko_EpetraOperatorWrapper.hpp"
99 #include "Teko_TpetraHelpers.hpp"
100 #include "Teko_TpetraOperatorWrapper.hpp"
101 
102 // Tpetra
103 #include "Thyra_TpetraLinearOp.hpp"
104 #include "Tpetra_CrsMatrix.hpp"
105 #include "Tpetra_Vector.hpp"
106 #include "Thyra_TpetraThyraWrappers.hpp"
107 #include "TpetraExt_MatrixMatrix.hpp"
108 #include "Tpetra_RowMatrixTransposer.hpp"
109 
110 #include <cmath>
111 
112 namespace Teko {
113 
114 using Teuchos::rcp;
115 using Teuchos::rcp_dynamic_cast;
116 using Teuchos::RCP;
117 #ifdef Teko_ENABLE_Isorropia
118 using Isorropia::Epetra::Prober;
119 #endif
120 
121 const Teuchos::RCP<Teuchos::FancyOStream> getOutputStream()
122 {
123  Teuchos::RCP<Teuchos::FancyOStream> os =
124  Teuchos::VerboseObjectBase::getDefaultOStream();
125 
126  //os->setShowProcRank(true);
127  //os->setOutputToRootOnly(-1);
128  return os;
129 }
130 
131 // distance function...not parallel...entirely internal to this cpp file
132 inline double dist(int dim,double * coords,int row,int col)
133 {
134  double value = 0.0;
135  for(int i=0;i<dim;i++)
136  value += std::pow(coords[dim*row+i]-coords[dim*col+i],2.0);
137 
138  // the distance between the two
139  return std::sqrt(value);
140 }
141 
142 // distance function...not parallel...entirely internal to this cpp file
143 inline double dist(double * x,double * y,double * z,int stride,int row,int col)
144 {
145  double value = 0.0;
146  if(x!=0) value += std::pow(x[stride*row]-x[stride*col],2.0);
147  if(y!=0) value += std::pow(y[stride*row]-y[stride*col],2.0);
148  if(z!=0) value += std::pow(z[stride*row]-z[stride*col],2.0);
149 
150  // the distance between the two
151  return std::sqrt(value);
152 }
153 
172 RCP<Epetra_CrsMatrix> buildGraphLaplacian(int dim,double * coords,const Epetra_CrsMatrix & stencil)
173 {
174  // allocate a new matrix with storage for the laplacian...in case of diagonals add one extra storage
175  RCP<Epetra_CrsMatrix> gl = rcp(new Epetra_CrsMatrix(Copy,stencil.RowMap(),stencil.ColMap(),
176  stencil.MaxNumEntries()+1),true);
177 
178  // allocate an additional value for the diagonal, if neccessary
179  std::vector<double> rowData(stencil.GlobalMaxNumEntries()+1);
180  std::vector<int> rowInd(stencil.GlobalMaxNumEntries()+1);
181 
182  // loop over all the rows
183  for(int j=0;j<gl->NumMyRows();j++) {
184  int row = gl->GRID(j);
185  double diagValue = 0.0;
186  int diagInd = -1;
187  int rowSz = 0;
188 
189  // extract a copy of this row...put it in rowData, rowIndicies
190  stencil.ExtractGlobalRowCopy(row,stencil.MaxNumEntries(),rowSz,&rowData[0],&rowInd[0]);
191 
192  // loop over elements of row
193  for(int i=0;i<rowSz;i++) {
194  int col = rowInd[i];
195 
196  // is this a 0 entry masquerading as some thing else?
197  double value = rowData[i];
198  if(value==0) continue;
199 
200  // for nondiagonal entries
201  if(row!=col) {
202  double d = dist(dim,coords,row,col);
203  rowData[i] = -1.0/d;
204  diagValue += rowData[i];
205  }
206  else
207  diagInd = i;
208  }
209 
210  // handle diagonal entry
211  if(diagInd<0) { // diagonal not in row
212  rowData[rowSz] = -diagValue;
213  rowInd[rowSz] = row;
214  rowSz++;
215  }
216  else { // diagonal in row
217  rowData[diagInd] = -diagValue;
218  rowInd[diagInd] = row;
219  }
220 
221  // insert row data into graph Laplacian matrix
222  TEUCHOS_TEST_FOR_EXCEPT(gl->InsertGlobalValues(row,rowSz,&rowData[0],&rowInd[0]));
223  }
224 
225  gl->FillComplete();
226 
227  return gl;
228 }
229 
230 RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > buildGraphLaplacian(int dim,ST * coords,const Tpetra::CrsMatrix<ST,LO,GO,NT> & stencil)
231 {
232  // allocate a new matrix with storage for the laplacian...in case of diagonals add one extra storage
233  RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > gl = rcp(new Tpetra::CrsMatrix<ST,LO,GO,NT>(stencil.getRowMap(),stencil.getColMap(),
234  stencil.getGlobalMaxNumRowEntries()+1));
235 
236  // allocate an additional value for the diagonal, if neccessary
237  std::vector<ST> rowData(stencil.getGlobalMaxNumRowEntries()+1);
238  std::vector<GO> rowInd(stencil.getGlobalMaxNumRowEntries()+1);
239 
240  // loop over all the rows
241  for(LO j=0;j< (LO) gl->getNodeNumRows();j++) {
242  GO row = gl->getRowMap()->getGlobalElement(j);
243  ST diagValue = 0.0;
244  GO diagInd = -1;
245  size_t rowSz = 0;
246 
247  // extract a copy of this row...put it in rowData, rowIndicies
248  stencil.getGlobalRowCopy(row,Teuchos::ArrayView<GO>(rowInd),Teuchos::ArrayView<ST>(rowData),rowSz);
249 
250  // loop over elements of row
251  for(size_t i=0;i<rowSz;i++) {
252  GO col = rowInd[i];
253 
254  // is this a 0 entry masquerading as some thing else?
255  ST value = rowData[i];
256  if(value==0) continue;
257 
258  // for nondiagonal entries
259  if(row!=col) {
260  ST d = dist(dim,coords,row,col);
261  rowData[i] = -1.0/d;
262  diagValue += rowData[i];
263  }
264  else
265  diagInd = i;
266  }
267 
268  // handle diagonal entry
269  if(diagInd<0) { // diagonal not in row
270  rowData[rowSz] = -diagValue;
271  rowInd[rowSz] = row;
272  rowSz++;
273  }
274  else { // diagonal in row
275  rowData[diagInd] = -diagValue;
276  rowInd[diagInd] = row;
277  }
278 
279  // insert row data into graph Laplacian matrix
280  gl->replaceGlobalValues(row,Teuchos::ArrayView<GO>(rowInd),Teuchos::ArrayView<ST>(rowData));
281  }
282 
283  gl->fillComplete();
284 
285  return gl;
286 }
287 
310 RCP<Epetra_CrsMatrix> buildGraphLaplacian(double * x,double * y,double * z,int stride,const Epetra_CrsMatrix & stencil)
311 {
312  // allocate a new matrix with storage for the laplacian...in case of diagonals add one extra storage
313  RCP<Epetra_CrsMatrix> gl = rcp(new Epetra_CrsMatrix(Copy,stencil.RowMap(),stencil.ColMap(),
314  stencil.MaxNumEntries()+1),true);
315 
316  // allocate an additional value for the diagonal, if neccessary
317  std::vector<double> rowData(stencil.GlobalMaxNumEntries()+1);
318  std::vector<int> rowInd(stencil.GlobalMaxNumEntries()+1);
319 
320  // loop over all the rows
321  for(int j=0;j<gl->NumMyRows();j++) {
322  int row = gl->GRID(j);
323  double diagValue = 0.0;
324  int diagInd = -1;
325  int rowSz = 0;
326 
327  // extract a copy of this row...put it in rowData, rowIndicies
328  stencil.ExtractGlobalRowCopy(row,stencil.MaxNumEntries(),rowSz,&rowData[0],&rowInd[0]);
329 
330  // loop over elements of row
331  for(int i=0;i<rowSz;i++) {
332  int col = rowInd[i];
333 
334  // is this a 0 entry masquerading as some thing else?
335  double value = rowData[i];
336  if(value==0) continue;
337 
338  // for nondiagonal entries
339  if(row!=col) {
340  double d = dist(x,y,z,stride,row,col);
341  rowData[i] = -1.0/d;
342  diagValue += rowData[i];
343  }
344  else
345  diagInd = i;
346  }
347 
348  // handle diagonal entry
349  if(diagInd<0) { // diagonal not in row
350  rowData[rowSz] = -diagValue;
351  rowInd[rowSz] = row;
352  rowSz++;
353  }
354  else { // diagonal in row
355  rowData[diagInd] = -diagValue;
356  rowInd[diagInd] = row;
357  }
358 
359  // insert row data into graph Laplacian matrix
360  TEUCHOS_TEST_FOR_EXCEPT(gl->InsertGlobalValues(row,rowSz,&rowData[0],&rowInd[0]));
361  }
362 
363  gl->FillComplete();
364 
365  return gl;
366 }
367 
368 RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > buildGraphLaplacian(ST * x,ST * y,ST * z,GO stride,const Tpetra::CrsMatrix<ST,LO,GO,NT> & stencil)
369 {
370  // allocate a new matrix with storage for the laplacian...in case of diagonals add one extra storage
371  RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > gl = rcp(new Tpetra::CrsMatrix<ST,LO,GO,NT>(stencil.getRowMap(),stencil.getColMap(),
372  stencil.getGlobalMaxNumRowEntries()+1),true);
373 
374  // allocate an additional value for the diagonal, if neccessary
375  std::vector<ST> rowData(stencil.getGlobalMaxNumRowEntries()+1);
376  std::vector<GO> rowInd(stencil.getGlobalMaxNumRowEntries()+1);
377 
378  // loop over all the rows
379  for(LO j=0;j<(LO) gl->getNodeNumRows();j++) {
380  GO row = gl->getRowMap()->getGlobalElement(j);
381  ST diagValue = 0.0;
382  GO diagInd = -1;
383  size_t rowSz = 0;
384 
385  // extract a copy of this row...put it in rowData, rowIndicies
386  stencil.getGlobalRowCopy(row,Teuchos::ArrayView<GO>(rowInd),Teuchos::ArrayView<ST>(rowData),rowSz);
387 
388  // loop over elements of row
389  for(size_t i=0;i<rowSz;i++) {
390  GO col = rowInd[i];
391 
392  // is this a 0 entry masquerading as some thing else?
393  ST value = rowData[i];
394  if(value==0) continue;
395 
396  // for nondiagonal entries
397  if(row!=col) {
398  ST d = dist(x,y,z,stride,row,col);
399  rowData[i] = -1.0/d;
400  diagValue += rowData[i];
401  }
402  else
403  diagInd = i;
404  }
405 
406  // handle diagonal entry
407  if(diagInd<0) { // diagonal not in row
408  rowData[rowSz] = -diagValue;
409  rowInd[rowSz] = row;
410  rowSz++;
411  }
412  else { // diagonal in row
413  rowData[diagInd] = -diagValue;
414  rowInd[diagInd] = row;
415  }
416 
417  // insert row data into graph Laplacian matrix
418  gl->replaceGlobalValues(row,Teuchos::ArrayView<GO>(rowInd),Teuchos::ArrayView<ST>(rowData));
419  }
420 
421  gl->fillComplete();
422 
423  return gl;
424 }
425 
441 void applyOp(const LinearOp & A,const MultiVector & x,MultiVector & y,double alpha,double beta)
442 {
443  Thyra::apply(*A,Thyra::NOTRANS,*x,y.ptr(),alpha,beta);
444 }
445 
461 void applyTransposeOp(const LinearOp & A,const MultiVector & x,MultiVector & y,double alpha,double beta)
462 {
463  Thyra::apply(*A,Thyra::TRANS,*x,y.ptr(),alpha,beta);
464 }
465 
468 void update(double alpha,const MultiVector & x,double beta,MultiVector & y)
469 {
470  Teuchos::Array<double> scale;
471  Teuchos::Array<Teuchos::Ptr<const Thyra::MultiVectorBase<double> > > vec;
472 
473  // build arrays needed for linear combo
474  scale.push_back(alpha);
475  vec.push_back(x.ptr());
476 
477  // compute linear combination
478  Thyra::linear_combination<double>(scale,vec,beta,y.ptr());
479 }
480 
482 BlockedLinearOp getUpperTriBlocks(const BlockedLinearOp & blo,bool callEndBlockFill)
483 {
484  int rows = blockRowCount(blo);
485 
486  TEUCHOS_ASSERT(rows==blockColCount(blo));
487 
488  RCP<const Thyra::ProductVectorSpaceBase<double> > range = blo->productRange();
489  RCP<const Thyra::ProductVectorSpaceBase<double> > domain = blo->productDomain();
490 
491  // allocate new operator
492  BlockedLinearOp upper = createBlockedOp();
493 
494  // build new operator
495  upper->beginBlockFill(rows,rows);
496 
497  for(int i=0;i<rows;i++) {
498  // put zero operators on the diagonal
499  // this gurantees the vector space of
500  // the new operator are fully defined
501  RCP<const Thyra::LinearOpBase<double> > zed
502  = Thyra::zero<double>(range->getBlock(i),domain->getBlock(i));
503  upper->setBlock(i,i,zed);
504 
505  for(int j=i+1;j<rows;j++) {
506  // get block i,j
507  LinearOp uij = blo->getBlock(i,j);
508 
509  // stuff it in U
510  if(uij!=Teuchos::null)
511  upper->setBlock(i,j,uij);
512  }
513  }
514  if(callEndBlockFill)
515  upper->endBlockFill();
516 
517  return upper;
518 }
519 
521 BlockedLinearOp getLowerTriBlocks(const BlockedLinearOp & blo,bool callEndBlockFill)
522 {
523  int rows = blockRowCount(blo);
524 
525  TEUCHOS_ASSERT(rows==blockColCount(blo));
526 
527  RCP<const Thyra::ProductVectorSpaceBase<double> > range = blo->productRange();
528  RCP<const Thyra::ProductVectorSpaceBase<double> > domain = blo->productDomain();
529 
530  // allocate new operator
531  BlockedLinearOp lower = createBlockedOp();
532 
533  // build new operator
534  lower->beginBlockFill(rows,rows);
535 
536  for(int i=0;i<rows;i++) {
537  // put zero operators on the diagonal
538  // this gurantees the vector space of
539  // the new operator are fully defined
540  RCP<const Thyra::LinearOpBase<double> > zed
541  = Thyra::zero<double>(range->getBlock(i),domain->getBlock(i));
542  lower->setBlock(i,i,zed);
543 
544  for(int j=0;j<i;j++) {
545  // get block i,j
546  LinearOp uij = blo->getBlock(i,j);
547 
548  // stuff it in U
549  if(uij!=Teuchos::null)
550  lower->setBlock(i,j,uij);
551  }
552  }
553  if(callEndBlockFill)
554  lower->endBlockFill();
555 
556  return lower;
557 }
558 
578 BlockedLinearOp zeroBlockedOp(const BlockedLinearOp & blo)
579 {
580  int rows = blockRowCount(blo);
581 
582  TEUCHOS_ASSERT(rows==blockColCount(blo)); // assert that matrix is square
583 
584  RCP<const Thyra::ProductVectorSpaceBase<double> > range = blo->productRange();
585  RCP<const Thyra::ProductVectorSpaceBase<double> > domain = blo->productDomain();
586 
587  // allocate new operator
588  BlockedLinearOp zeroOp = createBlockedOp();
589 
590  // build new operator
591  zeroOp->beginBlockFill(rows,rows);
592 
593  for(int i=0;i<rows;i++) {
594  // put zero operators on the diagonal
595  // this gurantees the vector space of
596  // the new operator are fully defined
597  RCP<const Thyra::LinearOpBase<double> > zed
598  = Thyra::zero<double>(range->getBlock(i),domain->getBlock(i));
599  zeroOp->setBlock(i,i,zed);
600  }
601 
602  return zeroOp;
603 }
604 
606 bool isZeroOp(const LinearOp op)
607 {
608  // if operator is null...then its zero!
609  if(op==Teuchos::null) return true;
610 
611  // try to cast it to a zero linear operator
612  LinearOp test = rcp_dynamic_cast<const Thyra::ZeroLinearOpBase<double> >(op);
613 
614  // if it works...then its zero...otherwise its null
615  if(test!=Teuchos::null) return true;
616 
617  // See if the operator is a wrapped zero op
618  ST scalar = 0.0;
619  Thyra::EOpTransp transp = Thyra::NOTRANS;
620  RCP<const Thyra::LinearOpBase<ST> > wrapped_op;
621  Thyra::unwrap(op, &scalar, &transp, &wrapped_op);
622  test = rcp_dynamic_cast<const Thyra::ZeroLinearOpBase<double> >(wrapped_op);
623  return test!=Teuchos::null;
624 }
625 
634 ModifiableLinearOp getAbsRowSumMatrix(const LinearOp & op)
635 {
636  bool isTpetra = false;
637  RCP<const Epetra_CrsMatrix> eCrsOp;
638  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOp;
639 
640  try {
641  // get Epetra or Tpetra Operator
642  RCP<const Thyra::EpetraLinearOp > eOp = rcp_dynamic_cast<const Thyra::EpetraLinearOp >(op);
643  RCP<const Thyra::TpetraLinearOp<ST,LO,GO,NT> > tOp = rcp_dynamic_cast<const Thyra::TpetraLinearOp<ST,LO,GO,NT> >(op);
644 
645  // cast it to a CrsMatrix
646  RCP<Teuchos::FancyOStream> out = Teuchos::VerboseObjectBase::getDefaultOStream();
647  if (!eOp.is_null()){
648  eCrsOp = rcp_dynamic_cast<const Epetra_CrsMatrix>(eOp->epetra_op(),true);
649  }
650  else if (!tOp.is_null()){
651  tCrsOp = rcp_dynamic_cast<const Tpetra::CrsMatrix<ST,LO,GO,NT> >(tOp->getConstTpetraOperator(),true);
652  isTpetra = true;
653  }
654  else
655  throw std::logic_error("Neither Epetra nor Tpetra");
656  }
657  catch (std::exception & e) {
658  RCP<Teuchos::FancyOStream> out = Teuchos::VerboseObjectBase::getDefaultOStream();
659 
660  *out << "Teko: getAbsRowSumMatrix requires an Epetra_CrsMatrix or a Tpetra::CrsMatrix\n";
661  *out << " Could not extract an Epetra_Operator or a Tpetra_Operator from a \"" << op->description() << std::endl;
662  *out << " OR\n";
663  *out << " Could not cast an Epetra_Operator to a Epetra_CrsMatrix or a Tpetra_Operator to a Tpetra::CrsMatrix\n";
664  *out << std::endl;
665  *out << "*** THROWN EXCEPTION ***\n";
666  *out << e.what() << std::endl;
667  *out << "************************\n";
668 
669  throw e;
670  }
671 
672  if(!isTpetra){
673  // extract diagonal
674  const RCP<Epetra_Vector> ptrDiag = rcp(new Epetra_Vector(eCrsOp->RowMap()));
675  Epetra_Vector & diag = *ptrDiag;
676 
677  // compute absolute value row sum
678  diag.PutScalar(0.0);
679  for(int i=0;i<eCrsOp->NumMyRows();i++) {
680  double * values = 0;
681  int numEntries;
682  eCrsOp->ExtractMyRowView(i,numEntries,values);
683 
684  // build abs value row sum
685  for(int j=0;j<numEntries;j++)
686  diag[i] += std::abs(values[j]);
687  }
688 
689  // build Thyra diagonal operator
690  return Teko::Epetra::thyraDiagOp(ptrDiag,eCrsOp->RowMap(),"absRowSum( " + op->getObjectLabel() + " )");
691  }
692  else {
693  // extract diagonal
694  const RCP<Tpetra::Vector<ST,LO,GO,NT> > ptrDiag = Tpetra::createVector<ST,LO,GO,NT>(tCrsOp->getRowMap());
695  Tpetra::Vector<ST,LO,GO,NT> & diag = *ptrDiag;
696 
697  // compute absolute value row sum
698  diag.putScalar(0.0);
699  for(LO i=0;i<(LO) tCrsOp->getNodeNumRows();i++) {
700  LO numEntries = tCrsOp->getNumEntriesInLocalRow (i);
701  std::vector<LO> indices(numEntries);
702  std::vector<ST> values(numEntries);
703  Teuchos::ArrayView<const LO> indices_av(indices);
704  Teuchos::ArrayView<const ST> values_av(values);
705  tCrsOp->getLocalRowView(i,indices_av,values_av);
706 
707  // build abs value row sum
708  for(LO j=0;j<numEntries;j++)
709  diag.sumIntoLocalValue(i,std::abs(values_av[j]));
710  }
711 
712  // build Thyra diagonal operator
713  return Teko::TpetraHelpers::thyraDiagOp(ptrDiag,*tCrsOp->getRowMap(),"absRowSum( " + op->getObjectLabel() + " ))");
714 
715  }
716 
717 }
718 
727 ModifiableLinearOp getAbsRowSumInvMatrix(const LinearOp & op)
728 {
729  // if this is a blocked operator, extract diagonals block by block
730  // FIXME: this does not add in values from off-diagonal blocks
731  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_op = rcp_dynamic_cast<const Thyra::PhysicallyBlockedLinearOpBase<double> >(op);
732  if(blocked_op != Teuchos::null){
733  int numRows = blocked_op->productRange()->numBlocks();
734  TEUCHOS_ASSERT(blocked_op->productDomain()->numBlocks() == numRows);
735  RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_diag = Thyra::defaultBlockedLinearOp<double>();
736  blocked_diag->beginBlockFill(numRows,numRows);
737  for(int r = 0; r < numRows; ++r){
738  for(int c = 0; c < numRows; ++c){
739  if(r==c)
740  blocked_diag->setNonconstBlock(r,c,getAbsRowSumInvMatrix(blocked_op->getBlock(r,c)));
741  else
742  blocked_diag->setBlock(r,c,Thyra::zero<double>(blocked_op->getBlock(r,c)->range(),blocked_op->getBlock(r,c)->domain()));
743  }
744  }
745  blocked_diag->endBlockFill();
746  return blocked_diag;
747  }
748 
749  if(Teko::TpetraHelpers::isTpetraLinearOp(op)) {
750  ST scalar = 0.0;
751  bool transp = false;
752  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOp = Teko::TpetraHelpers::getTpetraCrsMatrix(op, &scalar, &transp);
753 
754  // extract diagonal
755  const RCP<Tpetra::Vector<ST,LO,GO,NT> > ptrDiag = Tpetra::createVector<ST,LO,GO,NT>(tCrsOp->getRowMap());
756  Tpetra::Vector<ST,LO,GO,NT> & diag = *ptrDiag;
757 
758  // compute absolute value row sum
759  diag.putScalar(0.0);
760  for(LO i=0;i<(LO) tCrsOp->getNodeNumRows();i++) {
761  LO numEntries = tCrsOp->getNumEntriesInLocalRow (i);
762  std::vector<LO> indices(numEntries);
763  std::vector<ST> values(numEntries);
764  Teuchos::ArrayView<const LO> indices_av(indices);
765  Teuchos::ArrayView<const ST> values_av(values);
766  tCrsOp->getLocalRowView(i,indices_av,values_av);
767 
768  // build abs value row sum
769  for(LO j=0;j<numEntries;j++)
770  diag.sumIntoLocalValue(i,std::abs(values_av[j]));
771  }
772  diag.scale(scalar);
773  diag.reciprocal(diag); // invert entries
774 
775  // build Thyra diagonal operator
776  return Teko::TpetraHelpers::thyraDiagOp(ptrDiag,*tCrsOp->getRowMap(),"absRowSum( " + op->getObjectLabel() + " ))");
777 
778  }
779  else{
780  RCP<const Thyra::EpetraLinearOp > eOp = rcp_dynamic_cast<const Thyra::EpetraLinearOp >(op,true);
781  RCP<const Epetra_CrsMatrix> eCrsOp = rcp_dynamic_cast<const Epetra_CrsMatrix>(eOp->epetra_op(),true);
782 
783  // extract diagonal
784  const RCP<Epetra_Vector> ptrDiag = rcp(new Epetra_Vector(eCrsOp->RowMap()));
785  Epetra_Vector & diag = *ptrDiag;
786 
787  // compute absolute value row sum
788  diag.PutScalar(0.0);
789  for(int i=0;i<eCrsOp->NumMyRows();i++) {
790  double * values = 0;
791  int numEntries;
792  eCrsOp->ExtractMyRowView(i,numEntries,values);
793 
794  // build abs value row sum
795  for(int j=0;j<numEntries;j++)
796  diag[i] += std::abs(values[j]);
797  }
798  diag.Reciprocal(diag); // invert entries
799 
800  // build Thyra diagonal operator
801  return Teko::Epetra::thyraDiagOp(ptrDiag,eCrsOp->RowMap(),"absRowSum( " + op->getObjectLabel() + " )");
802  }
803 
804 }
805 
813 ModifiableLinearOp getLumpedMatrix(const LinearOp & op)
814 {
815  RCP<Thyra::VectorBase<ST> > ones = Thyra::createMember(op->domain());
816  RCP<Thyra::VectorBase<ST> > diag = Thyra::createMember(op->range());
817 
818  // set to all ones
819  Thyra::assign(ones.ptr(),1.0);
820 
821  // compute lumped diagonal
822  // Thyra::apply(*op,Thyra::NONCONJ_ELE,*ones,&*diag);
823  Thyra::apply(*op,Thyra::NOTRANS,*ones,diag.ptr());
824 
825  return rcp(new Thyra::DefaultDiagonalLinearOp<ST>(diag));
826 }
827 
836 ModifiableLinearOp getInvLumpedMatrix(const LinearOp & op)
837 {
838  RCP<Thyra::VectorBase<ST> > ones = Thyra::createMember(op->domain());
839  RCP<Thyra::VectorBase<ST> > diag = Thyra::createMember(op->range());
840 
841  // set to all ones
842  Thyra::assign(ones.ptr(),1.0);
843 
844  // compute lumped diagonal
845  Thyra::apply(*op,Thyra::NOTRANS,*ones,diag.ptr());
846  Thyra::reciprocal(*diag,diag.ptr());
847 
848  return rcp(new Thyra::DefaultDiagonalLinearOp<ST>(diag));
849 }
850 
862 const ModifiableLinearOp getDiagonalOp(const LinearOp & op)
863 {
864  bool isTpetra = false;
865  RCP<const Epetra_CrsMatrix> eCrsOp;
866  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOp;
867 
868  try {
869  // get Epetra or Tpetra Operator
870  RCP<const Thyra::EpetraLinearOp > eOp = rcp_dynamic_cast<const Thyra::EpetraLinearOp >(op);
871  RCP<const Thyra::TpetraLinearOp<ST,LO,GO,NT> > tOp = rcp_dynamic_cast<const Thyra::TpetraLinearOp<ST,LO,GO,NT> >(op);
872 
873  // cast it to a CrsMatrix
874  RCP<Teuchos::FancyOStream> out = Teuchos::VerboseObjectBase::getDefaultOStream();
875  if (!eOp.is_null()){
876  eCrsOp = rcp_dynamic_cast<const Epetra_CrsMatrix>(eOp->epetra_op(),true);
877  }
878  else if (!tOp.is_null()){
879  tCrsOp = rcp_dynamic_cast<const Tpetra::CrsMatrix<ST,LO,GO,NT> >(tOp->getConstTpetraOperator(),true);
880  isTpetra = true;
881  }
882  else
883  throw std::logic_error("Neither Epetra nor Tpetra");
884  }
885  catch (std::exception & e) {
886  RCP<Teuchos::FancyOStream> out = Teuchos::VerboseObjectBase::getDefaultOStream();
887 
888  *out << "Teko: getDiagonalOp requires an Epetra_CrsMatrix or a Tpetra::CrsMatrix\n";
889  *out << " Could not extract an Epetra_Operator or a Tpetra_Operator from a \"" << op->description() << std::endl;
890  *out << " OR\n";
891  *out << " Could not cast an Epetra_Operator to a Epetra_CrsMatrix or a Tpetra_Operator to a Tpetra::CrsMatrix\n";
892  *out << std::endl;
893  *out << "*** THROWN EXCEPTION ***\n";
894  *out << e.what() << std::endl;
895  *out << "************************\n";
896 
897  throw e;
898  }
899 
900  if(!isTpetra){
901  // extract diagonal
902  const RCP<Epetra_Vector> diag = rcp(new Epetra_Vector(eCrsOp->RowMap()));
903  TEUCHOS_TEST_FOR_EXCEPT(eCrsOp->ExtractDiagonalCopy(*diag));
904 
905  // build Thyra diagonal operator
906  return Teko::Epetra::thyraDiagOp(diag,eCrsOp->RowMap(),"inv(diag( " + op->getObjectLabel() + " ))");
907  }
908  else {
909  // extract diagonal
910  const RCP<Tpetra::Vector<ST,LO,GO,NT> > diag = Tpetra::createVector<ST,LO,GO,NT>(tCrsOp->getRowMap());
911  tCrsOp->getLocalDiagCopy(*diag);
912 
913  // build Thyra diagonal operator
914  return Teko::TpetraHelpers::thyraDiagOp(diag,*tCrsOp->getRowMap(),"inv(diag( " + op->getObjectLabel() + " ))");
915 
916  }
917 }
918 
919 const MultiVector getDiagonal(const LinearOp & op)
920 {
921  bool isTpetra = false;
922  RCP<const Epetra_CrsMatrix> eCrsOp;
923  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOp;
924 
925  try {
926  // get Epetra or Tpetra Operator
927  RCP<const Thyra::EpetraLinearOp > eOp = rcp_dynamic_cast<const Thyra::EpetraLinearOp >(op);
928  RCP<const Thyra::TpetraLinearOp<ST,LO,GO,NT> > tOp = rcp_dynamic_cast<const Thyra::TpetraLinearOp<ST,LO,GO,NT> >(op);
929 
930  // cast it to a CrsMatrix
931  RCP<Teuchos::FancyOStream> out = Teuchos::VerboseObjectBase::getDefaultOStream();
932  if (!eOp.is_null()){
933  eCrsOp = rcp_dynamic_cast<const Epetra_CrsMatrix>(eOp->epetra_op(),true);
934  }
935  else if (!tOp.is_null()){
936  tCrsOp = rcp_dynamic_cast<const Tpetra::CrsMatrix<ST,LO,GO,NT> >(tOp->getConstTpetraOperator(),true);
937  isTpetra = true;
938  }
939  else
940  throw std::logic_error("Neither Epetra nor Tpetra");
941  }
942  catch (std::exception & e) {
943  RCP<Teuchos::FancyOStream> out = Teuchos::VerboseObjectBase::getDefaultOStream();
944 
945  RCP<const Thyra::EpetraLinearOp > eOp = rcp_dynamic_cast<const Thyra::EpetraLinearOp >(op);
946  RCP<const Thyra::TpetraLinearOp<ST,LO,GO,NT> > tOp = rcp_dynamic_cast<const Thyra::TpetraLinearOp<ST,LO,GO,NT> >(op);
947  *out << eOp;
948  *out << tOp;
949 
950  *out << "Teko: getDiagonal requires an Epetra_CrsMatrix or a Tpetra::CrsMatrix\n";
951  *out << " Could not extract an Epetra_Operator or a Tpetra_Operator from a \"" << op->description() << std::endl;
952  *out << " OR\n";
953  *out << " Could not cast an Epetra_Operator to a Epetra_CrsMatrix or a Tpetra_Operator to a Tpetra::CrsMatrix\n";
954  *out << std::endl;
955  *out << "*** THROWN EXCEPTION ***\n";
956  *out << e.what() << std::endl;
957  *out << "************************\n";
958 
959  throw e;
960  }
961 
962  if(!isTpetra){
963  // extract diagonal
964  const RCP<Epetra_Vector> diag = rcp(new Epetra_Vector(eCrsOp->RowMap()));
965  TEUCHOS_TEST_FOR_EXCEPT(eCrsOp->ExtractDiagonalCopy(*diag));
966 
967  return Thyra::create_Vector(diag,Thyra::create_VectorSpace(Teuchos::rcpFromRef(eCrsOp->RowMap())));
968  }
969  else {
970  // extract diagonal
971  const RCP<Tpetra::Vector<ST,LO,GO,NT> > diag = Tpetra::createVector<ST,LO,GO,NT>(tCrsOp->getRowMap());
972  tCrsOp->getLocalDiagCopy(*diag);
973 
974  // build Thyra diagonal operator
975  return Thyra::createVector<ST,LO,GO,NT>(diag,Thyra::createVectorSpace<ST,LO,GO,NT>(tCrsOp->getRowMap()));
976 
977  }
978 }
979 
980 const MultiVector getDiagonal(const Teko::LinearOp & A,const DiagonalType & dt)
981 {
982  LinearOp diagOp = Teko::getDiagonalOp(A,dt);
983 
984  Teuchos::RCP<const Thyra::MultiVectorBase<double> > v =
985  Teuchos::rcp_dynamic_cast<const Thyra::DiagonalLinearOpBase<double> >(diagOp)->getDiag();
986  return Teuchos::rcp_const_cast<Thyra::MultiVectorBase<double> >(v);
987 }
988 
1000 const ModifiableLinearOp getInvDiagonalOp(const LinearOp & op)
1001 {
1002  // if this is a diagonal linear op already, just take the reciprocal
1003  auto diagonal_op = rcp_dynamic_cast<const Thyra::DiagonalLinearOpBase<double>>(op);
1004  if(diagonal_op != Teuchos::null){
1005  auto diag = diagonal_op->getDiag();
1006  auto inv_diag = diag->clone_v();
1007  Thyra::reciprocal(*diag,inv_diag.ptr());
1008  return rcp(new Thyra::DefaultDiagonalLinearOp<double>(inv_diag));
1009  }
1010 
1011  // if this is a blocked operator, extract diagonals block by block
1012  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_op = rcp_dynamic_cast<const Thyra::PhysicallyBlockedLinearOpBase<double> >(op);
1013  if(blocked_op != Teuchos::null){
1014  int numRows = blocked_op->productRange()->numBlocks();
1015  TEUCHOS_ASSERT(blocked_op->productDomain()->numBlocks() == numRows);
1016  RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_diag = Thyra::defaultBlockedLinearOp<double>();
1017  blocked_diag->beginBlockFill(numRows,numRows);
1018  for(int r = 0; r < numRows; ++r){
1019  for(int c = 0; c < numRows; ++c){
1020  if(r==c)
1021  blocked_diag->setNonconstBlock(r,c,getInvDiagonalOp(blocked_op->getBlock(r,c)));
1022  else
1023  blocked_diag->setBlock(r,c,Thyra::zero<double>(blocked_op->getBlock(r,c)->range(),blocked_op->getBlock(r,c)->domain()));
1024  }
1025  }
1026  blocked_diag->endBlockFill();
1027  return blocked_diag;
1028  }
1029 
1030  if (Teko::TpetraHelpers::isTpetraLinearOp(op)){
1031  ST scalar = 0.0;
1032  bool transp = false;
1033  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOp = Teko::TpetraHelpers::getTpetraCrsMatrix(op, &scalar, &transp);
1034 
1035  // extract diagonal
1036  const RCP<Tpetra::Vector<ST,LO,GO,NT> > diag = Tpetra::createVector<ST,LO,GO,NT>(tCrsOp->getRowMap());
1037  diag->scale(scalar);
1038  tCrsOp->getLocalDiagCopy(*diag);
1039  diag->reciprocal(*diag);
1040 
1041  // build Thyra diagonal operator
1042  return Teko::TpetraHelpers::thyraDiagOp(diag,*tCrsOp->getRowMap(),"inv(diag( " + op->getObjectLabel() + " ))");
1043 
1044  }
1045  else {
1046  RCP<const Thyra::EpetraLinearOp > eOp = rcp_dynamic_cast<const Thyra::EpetraLinearOp >(op,true);
1047  RCP<const Epetra_CrsMatrix> eCrsOp = rcp_dynamic_cast<const Epetra_CrsMatrix>(eOp->epetra_op(),true);
1048 
1049  // extract diagonal
1050  const RCP<Epetra_Vector> diag = rcp(new Epetra_Vector(eCrsOp->RowMap()));
1051  TEUCHOS_TEST_FOR_EXCEPT(eCrsOp->ExtractDiagonalCopy(*diag));
1052  diag->Reciprocal(*diag);
1053 
1054  // build Thyra diagonal operator
1055  return Teko::Epetra::thyraDiagOp(diag,eCrsOp->RowMap(),"inv(diag( " + op->getObjectLabel() + " ))");
1056  }
1057 }
1058 
1071 const LinearOp explicitMultiply(const LinearOp & opl,const LinearOp & opm,const LinearOp & opr)
1072 {
1073  // if this is a blocked operator, multiply block by block
1074  // it is possible that not every factor in the product is blocked and these situations are handled separately
1075 
1076  bool isBlockedL = isPhysicallyBlockedLinearOp(opl);
1077  bool isBlockedM = isPhysicallyBlockedLinearOp(opm);
1078  bool isBlockedR = isPhysicallyBlockedLinearOp(opr);
1079 
1080  // all factors blocked
1081  if((isBlockedL && isBlockedM && isBlockedR)){
1082 
1083  double scalarl = 0.0;
1084  bool transpl = false;
1085  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opl = getPhysicallyBlockedLinearOp(opl,&scalarl,&transpl);
1086  double scalarm = 0.0;
1087  bool transpm = false;
1088  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opm = getPhysicallyBlockedLinearOp(opm,&scalarm,&transpm);
1089  double scalarr = 0.0;
1090  bool transpr = false;
1091  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opr = getPhysicallyBlockedLinearOp(opr,&scalarr,&transpr);
1092  double scalar = scalarl*scalarm*scalarr;
1093 
1094  int numRows = blocked_opl->productRange()->numBlocks();
1095  int numCols = blocked_opr->productDomain()->numBlocks();
1096  int numMiddle = blocked_opm->productRange()->numBlocks();
1097 
1098  // Assume that the middle block is block nxn and that it's diagonal. Otherwise use the two argument explicitMultiply twice
1099  TEUCHOS_ASSERT(blocked_opm->productDomain()->numBlocks() == numMiddle);
1100  TEUCHOS_ASSERT(blocked_opl->productDomain()->numBlocks() == numMiddle);
1101  TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == numMiddle);
1102 
1103  RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_product = Thyra::defaultBlockedLinearOp<double>();
1104  blocked_product->beginBlockFill(numRows,numCols);
1105  for(int r = 0; r < numRows; ++r){
1106  for(int c = 0; c < numCols; ++c){
1107  LinearOp product_rc = explicitMultiply(blocked_opl->getBlock(r,0),blocked_opm->getBlock(0,0),blocked_opr->getBlock(0,c));
1108  for(int m = 1; m < numMiddle; ++m){
1109  LinearOp product_m = explicitMultiply(blocked_opl->getBlock(r,m),blocked_opm->getBlock(m,m),blocked_opr->getBlock(m,c));
1110  product_rc = explicitAdd(product_rc,product_m);
1111  }
1112  blocked_product->setBlock(r,c,product_rc);
1113  }
1114  }
1115  blocked_product->endBlockFill();
1116  return Thyra::scale<double>(scalar,blocked_product.getConst());
1117  }
1118 
1119  // left and right factors blocked
1120  if(isBlockedL && !isBlockedM && isBlockedR){
1121  double scalarl = 0.0;
1122  bool transpl = false;
1123  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opl = getPhysicallyBlockedLinearOp(opl,&scalarl,&transpl);
1124  double scalarr = 0.0;
1125  bool transpr = false;
1126  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opr = getPhysicallyBlockedLinearOp(opr,&scalarr,&transpr);
1127  double scalar = scalarl*scalarr;
1128 
1129  int numRows = blocked_opl->productRange()->numBlocks();
1130  int numCols = blocked_opr->productDomain()->numBlocks();
1131  int numMiddle = 1;
1132 
1133  // Assume that the middle block is 1x1 diagonal. Left must be rx1, right 1xc
1134  TEUCHOS_ASSERT(blocked_opl->productDomain()->numBlocks() == numMiddle);
1135  TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == numMiddle);
1136 
1137  RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_product = Thyra::defaultBlockedLinearOp<double>();
1138  blocked_product->beginBlockFill(numRows,numCols);
1139  for(int r = 0; r < numRows; ++r){
1140  for(int c = 0; c < numCols; ++c){
1141  LinearOp product_rc = explicitMultiply(blocked_opl->getBlock(r,0),opm,blocked_opr->getBlock(0,c));
1142  blocked_product->setBlock(r,c,product_rc);
1143  }
1144  }
1145  blocked_product->endBlockFill();
1146  return Thyra::scale<double>(scalar,blocked_product.getConst());
1147  }
1148 
1149  // only right factor blocked
1150  if(!isBlockedL && !isBlockedM && isBlockedR){
1151  double scalarr = 0.0;
1152  bool transpr = false;
1153  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opr = getPhysicallyBlockedLinearOp(opr,&scalarr,&transpr);
1154  double scalar = scalarr;
1155 
1156  int numRows = 1;
1157  int numCols = blocked_opr->productDomain()->numBlocks();
1158  int numMiddle = 1;
1159 
1160  // Assume that the middle block is 1x1 diagonal, left is 1x1. Right must be 1xc
1161  TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == numMiddle);
1162 
1163  RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_product = Thyra::defaultBlockedLinearOp<double>();
1164  blocked_product->beginBlockFill(numRows,numCols);
1165  for(int c = 0; c < numCols; ++c){
1166  LinearOp product_c = explicitMultiply(opl,opm,blocked_opr->getBlock(0,c));
1167  blocked_product->setBlock(0,c,product_c);
1168  }
1169  blocked_product->endBlockFill();
1170  return Thyra::scale<double>(scalar,blocked_product.getConst());
1171  }
1172 
1173  //TODO: three more cases (only non-blocked - blocked - non-blocked not possible)
1174 
1175  bool isTpetral = Teko::TpetraHelpers::isTpetraLinearOp(opl);
1176  bool isTpetram = Teko::TpetraHelpers::isTpetraLinearOp(opm);
1177  bool isTpetrar = Teko::TpetraHelpers::isTpetraLinearOp(opr);
1178 
1179  if(isTpetral && isTpetram && isTpetrar){ // Both operators are Tpetra matrices so explicitly multiply them
1180 
1181  // Get left and right Tpetra crs operators
1182  ST scalarl = 0.0;
1183  bool transpl = false;
1184  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpl = Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1185  ST scalarm = 0.0;
1186  bool transpm = false;
1187  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpm = Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarm, &transpm);
1188  ST scalarr = 0.0;
1189  bool transpr = false;
1190  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpr = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1191 
1192  // Build output operator
1193  RCP<Thyra::LinearOpBase<ST> > explicitOp = rcp(new Thyra::TpetraLinearOp<ST,LO,GO,NT>());
1194  RCP<Thyra::TpetraLinearOp<ST,LO,GO,NT> > tExplicitOp = rcp_dynamic_cast<Thyra::TpetraLinearOp<ST,LO,GO,NT> >(explicitOp);
1195 
1196  // Do explicit matrix-matrix multiply
1197  RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOplm = Tpetra::createCrsMatrix<ST,LO,GO,NT>(tCrsOpl->getRowMap());
1198  RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > explicitCrsOp = Tpetra::createCrsMatrix<ST,LO,GO,NT>(tCrsOpl->getRowMap());
1199  Tpetra::MatrixMatrix::Multiply<ST,LO,GO,NT>(*tCrsOpl,transpl,*tCrsOpm,transpm,*tCrsOplm);
1200  Tpetra::MatrixMatrix::Multiply<ST,LO,GO,NT>(*tCrsOplm,false,*tCrsOpr,transpr,*explicitCrsOp);
1201  explicitCrsOp->resumeFill();
1202  explicitCrsOp->scale(scalarl*scalarm*scalarr);
1203  explicitCrsOp->fillComplete(tCrsOpr->getDomainMap(),tCrsOpl->getRangeMap());
1204  tExplicitOp->initialize(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1205  return tExplicitOp;
1206 
1207  } else if (isTpetral && !isTpetram && isTpetrar){ // Assume that the middle operator is diagonal
1208 
1209  // Get left and right Tpetra crs operators
1210  ST scalarl = 0.0;
1211  bool transpl = false;
1212  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpl = Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1213  ST scalarr = 0.0;
1214  bool transpr = false;
1215  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpr = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1216 
1217  RCP<const Tpetra::Vector<ST,LO,GO,NT> > diagPtr;
1218 
1219  // Cast middle operator as DiagonalLinearOp and extract diagonal as Vector
1220  RCP<const Thyra::DiagonalLinearOpBase<ST> > dOpm = rcp_dynamic_cast<const Thyra::DiagonalLinearOpBase<ST> >(opm);
1221  if(dOpm != Teuchos::null){
1222  RCP<const Thyra::TpetraVector<ST,LO,GO,NT> > tPtr = rcp_dynamic_cast<const Thyra::TpetraVector<ST,LO,GO,NT> >(dOpm->getDiag(),true);
1223  diagPtr = rcp_dynamic_cast<const Tpetra::Vector<ST,LO,GO,NT> >(tPtr->getConstTpetraVector(),true);
1224  }
1225  // If it's not diagonal, maybe it's zero
1226  else if(rcp_dynamic_cast<const Thyra::ZeroLinearOpBase<ST> >(opm) != Teuchos::null){
1227  diagPtr = rcp(new Tpetra::Vector<ST,LO,GO,NT>(tCrsOpl->getDomainMap()));
1228  }
1229  else
1230  TEUCHOS_ASSERT(false);
1231 
1232  RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOplm = Tpetra::importAndFillCompleteCrsMatrix<Tpetra::CrsMatrix<ST,LO,GO,NT> >(tCrsOpl, Tpetra::Import<LO,GO,NT>(tCrsOpl->getRowMap(),tCrsOpl->getRowMap()));
1233 
1234  // Do the diagonal scaling
1235  tCrsOplm->rightScale(*diagPtr);
1236 
1237  // Build output operator
1238  RCP<Thyra::LinearOpBase<ST> > explicitOp = rcp(new Thyra::TpetraLinearOp<ST,LO,GO,NT>());
1239  RCP<Thyra::TpetraLinearOp<ST,LO,GO,NT> > tExplicitOp = rcp_dynamic_cast<Thyra::TpetraLinearOp<ST,LO,GO,NT> >(explicitOp);
1240 
1241  // Do explicit matrix-matrix multiply
1242  RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > explicitCrsOp = Tpetra::createCrsMatrix<ST,LO,GO,NT>(tCrsOpl->getRowMap());
1243  Tpetra::MatrixMatrix::Multiply<ST,LO,GO,NT>(*tCrsOplm,false,*tCrsOpr,transpr,*explicitCrsOp);
1244  explicitCrsOp->resumeFill();
1245  explicitCrsOp->scale(scalarl*scalarr);
1246  explicitCrsOp->fillComplete(tCrsOpr->getDomainMap(),tCrsOpl->getRangeMap());
1247  tExplicitOp->initialize(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1248  return tExplicitOp;
1249 
1250  } else { // Assume Epetra and we can use transformers
1251 
1252  // build implicit multiply
1253  const LinearOp implicitOp = Thyra::multiply(opl,opm,opr);
1254 
1255  // build transformer
1256  const RCP<Thyra::LinearOpTransformerBase<double> > prodTrans =
1257  Thyra::epetraExtDiagScaledMatProdTransformer();
1258 
1259  // build operator and multiply
1260  const RCP<Thyra::LinearOpBase<double> > explicitOp = prodTrans->createOutputOp();
1261  prodTrans->transform(*implicitOp,explicitOp.ptr());
1262  explicitOp->setObjectLabel("explicit( " + opl->getObjectLabel() +
1263  " * " + opm->getObjectLabel() +
1264  " * " + opr->getObjectLabel() + " )");
1265 
1266  return explicitOp;
1267 
1268  }
1269 }
1270 
1285 const ModifiableLinearOp explicitMultiply(const LinearOp & opl,const LinearOp & opm,const LinearOp & opr,
1286  const ModifiableLinearOp & destOp)
1287 {
1288  bool isTpetral = Teko::TpetraHelpers::isTpetraLinearOp(opl);
1289  bool isTpetram = Teko::TpetraHelpers::isTpetraLinearOp(opm);
1290  bool isTpetrar = Teko::TpetraHelpers::isTpetraLinearOp(opr);
1291 
1292  if(isTpetral && isTpetram && isTpetrar){ // Both operators are Tpetra matrices so explicitly multiply them
1293 
1294  // Get left and right Tpetra crs operators
1295  ST scalarl = 0.0;
1296  bool transpl = false;
1297  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpl = Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1298  ST scalarm = 0.0;
1299  bool transpm = false;
1300  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpm = Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarm, &transpm);
1301  ST scalarr = 0.0;
1302  bool transpr = false;
1303  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpr = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1304 
1305  // Build output operator
1306  RCP<Thyra::LinearOpBase<ST> > explicitOp = rcp(new Thyra::TpetraLinearOp<ST,LO,GO,NT>());
1307  RCP<Thyra::TpetraLinearOp<ST,LO,GO,NT> > tExplicitOp = rcp_dynamic_cast<Thyra::TpetraLinearOp<ST,LO,GO,NT> >(explicitOp);
1308 
1309  // Do explicit matrix-matrix multiply
1310  RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOplm = Tpetra::createCrsMatrix<ST,LO,GO,NT>(tCrsOpl->getRowMap());
1311  RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > explicitCrsOp = Tpetra::createCrsMatrix<ST,LO,GO,NT>(tCrsOpl->getRowMap());
1312  Tpetra::MatrixMatrix::Multiply<ST,LO,GO,NT>(*tCrsOpl,transpl,*tCrsOpm,transpm,*tCrsOplm);
1313  Tpetra::MatrixMatrix::Multiply<ST,LO,GO,NT>(*tCrsOplm,false,*tCrsOpr,transpr,*explicitCrsOp);
1314  explicitCrsOp->resumeFill();
1315  explicitCrsOp->scale(scalarl*scalarm*scalarr);
1316  explicitCrsOp->fillComplete(tCrsOpr->getDomainMap(),tCrsOpl->getRangeMap());
1317  tExplicitOp->initialize(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1318  return tExplicitOp;
1319 
1320  } else if (isTpetral && !isTpetram && isTpetrar){ // Assume that the middle operator is diagonal
1321 
1322  // Get left and right Tpetra crs operators
1323  ST scalarl = 0.0;
1324  bool transpl = false;
1325  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpl = Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1326  ST scalarr = 0.0;
1327  bool transpr = false;
1328  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpr = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1329 
1330  // Cast middle operator as DiagonalLinearOp and extract diagonal as Vector
1331  RCP<const Thyra::DiagonalLinearOpBase<ST> > dOpm = rcp_dynamic_cast<const Thyra::DiagonalLinearOpBase<ST> >(opm,true);
1332  RCP<const Thyra::TpetraVector<ST,LO,GO,NT> > tPtr = rcp_dynamic_cast<const Thyra::TpetraVector<ST,LO,GO,NT> >(dOpm->getDiag(),true);
1333  RCP<const Tpetra::Vector<ST,LO,GO,NT> > diagPtr = rcp_dynamic_cast<const Tpetra::Vector<ST,LO,GO,NT> >(tPtr->getConstTpetraVector(),true);
1334  RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOplm = Tpetra::importAndFillCompleteCrsMatrix<Tpetra::CrsMatrix<ST,LO,GO,NT> >(tCrsOpl, Tpetra::Import<LO,GO,NT>(tCrsOpl->getRowMap(),tCrsOpl->getRowMap()));
1335 
1336  // Do the diagonal scaling
1337  tCrsOplm->rightScale(*diagPtr);
1338 
1339  // Build output operator
1340  RCP<Thyra::LinearOpBase<ST> > explicitOp = rcp(new Thyra::TpetraLinearOp<ST,LO,GO,NT>());
1341  RCP<Thyra::TpetraLinearOp<ST,LO,GO,NT> > tExplicitOp = rcp_dynamic_cast<Thyra::TpetraLinearOp<ST,LO,GO,NT> >(explicitOp);
1342 
1343  // Do explicit matrix-matrix multiply
1344  RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > explicitCrsOp = Tpetra::createCrsMatrix<ST,LO,GO,NT>(tCrsOpl->getRowMap());
1345  Tpetra::MatrixMatrix::Multiply<ST,LO,GO,NT>(*tCrsOplm,false,*tCrsOpr,transpr,*explicitCrsOp);
1346  explicitCrsOp->resumeFill();
1347  explicitCrsOp->scale(scalarl*scalarr);
1348  explicitCrsOp->fillComplete(tCrsOpr->getDomainMap(),tCrsOpl->getRangeMap());
1349  tExplicitOp->initialize(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1350  return tExplicitOp;
1351 
1352  } else { // Assume Epetra and we can use transformers
1353 
1354  // build implicit multiply
1355  const LinearOp implicitOp = Thyra::multiply(opl,opm,opr);
1356 
1357  // build transformer
1358  const RCP<Thyra::LinearOpTransformerBase<double> > prodTrans =
1359  Thyra::epetraExtDiagScaledMatProdTransformer();
1360 
1361  // build operator destination operator
1362  ModifiableLinearOp explicitOp;
1363 
1364  // if neccessary build a operator to put the explicit multiply into
1365  if(destOp==Teuchos::null)
1366  explicitOp = prodTrans->createOutputOp();
1367  else
1368  explicitOp = destOp;
1369 
1370  // perform multiplication
1371  prodTrans->transform(*implicitOp,explicitOp.ptr());
1372 
1373  // label it
1374  explicitOp->setObjectLabel("explicit( " + opl->getObjectLabel() +
1375  " * " + opm->getObjectLabel() +
1376  " * " + opr->getObjectLabel() + " )");
1377 
1378  return explicitOp;
1379 
1380  }
1381 }
1382 
1393 const LinearOp explicitMultiply(const LinearOp & opl,const LinearOp & opr)
1394 {
1395  // if this is a blocked operator, multiply block by block
1396  // it is possible that not every factor in the product is blocked and these situations are handled separately
1397 
1398  bool isBlockedL = isPhysicallyBlockedLinearOp(opl);
1399  bool isBlockedR = isPhysicallyBlockedLinearOp(opr);
1400 
1401  // both factors blocked
1402  if((isBlockedL && isBlockedR)){
1403  double scalarl = 0.0;
1404  bool transpl = false;
1405  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opl = getPhysicallyBlockedLinearOp(opl,&scalarl,&transpl);
1406  double scalarr = 0.0;
1407  bool transpr = false;
1408  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opr = getPhysicallyBlockedLinearOp(opr,&scalarr,&transpr);
1409  double scalar = scalarl*scalarr;
1410 
1411  int numRows = blocked_opl->productRange()->numBlocks();
1412  int numCols = blocked_opr->productDomain()->numBlocks();
1413  int numMiddle = blocked_opl->productDomain()->numBlocks();
1414 
1415  TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == numMiddle);
1416 
1417  RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_product = Thyra::defaultBlockedLinearOp<double>();
1418  blocked_product->beginBlockFill(numRows,numCols);
1419  for(int r = 0; r < numRows; ++r){
1420  for(int c = 0; c < numCols; ++c){
1421  LinearOp product_rc = explicitMultiply(blocked_opl->getBlock(r,0),blocked_opr->getBlock(0,c));
1422  for(int m = 1; m < numMiddle; ++m){
1423  LinearOp product_m = explicitMultiply(blocked_opl->getBlock(r,m),blocked_opr->getBlock(m,c));
1424  product_rc = explicitAdd(product_rc,product_m);
1425  }
1426  blocked_product->setBlock(r,c,Thyra::scale(scalar,product_rc));
1427  }
1428  }
1429  blocked_product->endBlockFill();
1430  return blocked_product;
1431  }
1432 
1433  // only left factor blocked
1434  if((isBlockedL && !isBlockedR)){
1435  double scalarl = 0.0;
1436  bool transpl = false;
1437  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opl = getPhysicallyBlockedLinearOp(opl,&scalarl,&transpl);
1438  double scalar = scalarl;
1439 
1440  int numRows = blocked_opl->productRange()->numBlocks();
1441  int numCols = 1;
1442  int numMiddle = 1;
1443 
1444  TEUCHOS_ASSERT(blocked_opl->productDomain()->numBlocks() == numMiddle);
1445 
1446  RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_product = Thyra::defaultBlockedLinearOp<double>();
1447  blocked_product->beginBlockFill(numRows,numCols);
1448  for(int r = 0; r < numRows; ++r){
1449  LinearOp product_r = explicitMultiply(blocked_opl->getBlock(r,0),opr);
1450  blocked_product->setBlock(r,0,Thyra::scale(scalar,product_r));
1451  }
1452  blocked_product->endBlockFill();
1453  return blocked_product;
1454  }
1455 
1456  // only right factor blocked
1457  if((!isBlockedL && isBlockedR)){
1458  double scalarr = 0.0;
1459  bool transpr = false;
1460  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opr = getPhysicallyBlockedLinearOp(opr,&scalarr,&transpr);
1461  double scalar = scalarr;
1462 
1463  int numRows = 1;
1464  int numCols = blocked_opr->productDomain()->numBlocks();
1465  int numMiddle = 1;
1466 
1467  TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == numMiddle);
1468 
1469  RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_product = Thyra::defaultBlockedLinearOp<double>();
1470  blocked_product->beginBlockFill(numRows,numCols);
1471  for(int c = 0; c < numCols; ++c){
1472  LinearOp product_c = explicitMultiply(opl,blocked_opr->getBlock(0,c));
1473  blocked_product->setBlock(0,c,Thyra::scale(scalar,product_c));
1474  }
1475  blocked_product->endBlockFill();
1476  return blocked_product;
1477  }
1478 
1479  bool isTpetral = Teko::TpetraHelpers::isTpetraLinearOp(opl);
1480  bool isTpetrar = Teko::TpetraHelpers::isTpetraLinearOp(opr);
1481 
1482  if(isTpetral && isTpetrar){ // Both operators are Tpetra matrices so explicitly multiply them
1483  // Get left and right Tpetra crs operators
1484  ST scalarl = 0.0;
1485  bool transpl = false;
1486  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpl = Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1487  ST scalarr = 0.0;
1488  bool transpr = false;
1489  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpr = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1490 
1491  // Build output operator
1492  RCP<Thyra::LinearOpBase<ST> > explicitOp = rcp(new Thyra::TpetraLinearOp<ST,LO,GO,NT>());
1493  RCP<Thyra::TpetraLinearOp<ST,LO,GO,NT> > tExplicitOp = rcp_dynamic_cast<Thyra::TpetraLinearOp<ST,LO,GO,NT> >(explicitOp);
1494 
1495  // Do explicit matrix-matrix multiply
1496  RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > explicitCrsOp = Tpetra::createCrsMatrix<ST,LO,GO,NT>(tCrsOpl->getRowMap());
1497  Tpetra::MatrixMatrix::Multiply<ST,LO,GO,NT>(*tCrsOpl,transpl,*tCrsOpr,transpr,*explicitCrsOp);
1498  explicitCrsOp->resumeFill();
1499  explicitCrsOp->scale(scalarl*scalarr);
1500  explicitCrsOp->fillComplete(tCrsOpr->getDomainMap(),tCrsOpl->getRangeMap());
1501  tExplicitOp->initialize(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1502  return tExplicitOp;
1503 
1504  } else if (isTpetral && !isTpetrar){ // Assume that the right operator is diagonal
1505 
1506  // Get left Tpetra crs operator
1507  ST scalarl = 0.0;
1508  bool transpl = false;
1509  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpl = Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1510 
1511  // Cast right operator as DiagonalLinearOp and extract diagonal as Vector
1512  RCP<const Thyra::DiagonalLinearOpBase<ST> > dOpr = rcp_dynamic_cast<const Thyra::DiagonalLinearOpBase<ST> >(opr,true);
1513  RCP<const Thyra::TpetraVector<ST,LO,GO,NT> > tPtr = rcp_dynamic_cast<const Thyra::TpetraVector<ST,LO,GO,NT> >(dOpr->getDiag(),true);
1514  RCP<const Tpetra::Vector<ST,LO,GO,NT> > diagPtr = rcp_dynamic_cast<const Tpetra::Vector<ST,LO,GO,NT> >(tPtr->getConstTpetraVector(),true);
1515  RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > explicitCrsOp = Tpetra::importAndFillCompleteCrsMatrix<Tpetra::CrsMatrix<ST,LO,GO,NT> >(tCrsOpl, Tpetra::Import<LO,GO,NT>(tCrsOpl->getRowMap(),tCrsOpl->getRowMap()));
1516 
1517  explicitCrsOp->rightScale(*diagPtr);
1518  explicitCrsOp->resumeFill();
1519  explicitCrsOp->scale(scalarl);
1520  explicitCrsOp->fillComplete(tCrsOpl->getDomainMap(),tCrsOpl->getRangeMap());
1521 
1522  return Thyra::constTpetraLinearOp<ST,LO,GO,NT>(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1523 
1524  } else if (!isTpetral && isTpetrar){ // Assume that the left operator is diagonal
1525 
1526  // Get right Tpetra crs operator
1527  ST scalarr = 0.0;
1528  bool transpr = false;
1529  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpr = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1530 
1531  RCP<const Tpetra::Vector<ST,LO,GO,NT> > diagPtr;
1532 
1533  // Cast left operator as DiagonalLinearOp and extract diagonal as Vector
1534  RCP<const Thyra::DiagonalLinearOpBase<ST> > dOpl = rcp_dynamic_cast<const Thyra::DiagonalLinearOpBase<ST> >(opl);
1535  if(dOpl != Teuchos::null){
1536  RCP<const Thyra::TpetraVector<ST,LO,GO,NT> > tPtr = rcp_dynamic_cast<const Thyra::TpetraVector<ST,LO,GO,NT> >(dOpl->getDiag(),true);
1537  diagPtr = rcp_dynamic_cast<const Tpetra::Vector<ST,LO,GO,NT> >(tPtr->getConstTpetraVector(),true);
1538  }
1539  // If it's not diagonal, maybe it's zero
1540  else if(rcp_dynamic_cast<const Thyra::ZeroLinearOpBase<ST> >(opl) != Teuchos::null){
1541  diagPtr = rcp(new Tpetra::Vector<ST,LO,GO,NT>(tCrsOpr->getRangeMap()));
1542  }
1543  else
1544  TEUCHOS_ASSERT(false);
1545 
1546  RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > explicitCrsOp = Tpetra::importAndFillCompleteCrsMatrix<Tpetra::CrsMatrix<ST,LO,GO,NT> >(tCrsOpr, Tpetra::Import<LO,GO,NT>(tCrsOpr->getRowMap(),tCrsOpr->getRowMap()));
1547 
1548  explicitCrsOp->leftScale(*diagPtr);
1549  explicitCrsOp->resumeFill();
1550  explicitCrsOp->scale(scalarr);
1551  explicitCrsOp->fillComplete(tCrsOpr->getDomainMap(),tCrsOpr->getRangeMap());
1552 
1553  return Thyra::constTpetraLinearOp<ST,LO,GO,NT>(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1554 
1555  } else { // Assume Epetra and we can use transformers
1556 
1557  // build implicit multiply
1558  const LinearOp implicitOp = Thyra::multiply(opl,opr);
1559 
1560  // build a scaling transformer
1561  RCP<Thyra::LinearOpTransformerBase<double> > prodTrans
1562  = Thyra::epetraExtDiagScalingTransformer();
1563 
1564  // check to see if a scaling transformer works: if not use the
1565  // DiagScaledMatrixProduct transformer
1566  if(not prodTrans->isCompatible(*implicitOp))
1567  prodTrans = Thyra::epetraExtDiagScaledMatProdTransformer();
1568 
1569  // build operator and multiply
1570  const RCP<Thyra::LinearOpBase<double> > explicitOp = prodTrans->createOutputOp();
1571  prodTrans->transform(*implicitOp,explicitOp.ptr());
1572  explicitOp->setObjectLabel("explicit( " + opl->getObjectLabel() +
1573  " * " + opr->getObjectLabel() + " )");
1574 
1575  return explicitOp;
1576  }
1577 }
1578 
1592 const ModifiableLinearOp explicitMultiply(const LinearOp & opl,const LinearOp & opr,
1593  const ModifiableLinearOp & destOp)
1594 {
1595  // if this is a blocked operator, multiply block by block
1596  // it is possible that not every factor in the product is blocked and these situations are handled separately
1597 
1598  bool isBlockedL = isPhysicallyBlockedLinearOp(opl);
1599  bool isBlockedR = isPhysicallyBlockedLinearOp(opr);
1600 
1601  // both factors blocked
1602  if((isBlockedL && isBlockedR)){
1603  double scalarl = 0.0;
1604  bool transpl = false;
1605  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opl = getPhysicallyBlockedLinearOp(opl,&scalarl,&transpl);
1606  double scalarr = 0.0;
1607  bool transpr = false;
1608  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opr = getPhysicallyBlockedLinearOp(opr,&scalarr,&transpr);
1609  double scalar = scalarl*scalarr;
1610 
1611  int numRows = blocked_opl->productRange()->numBlocks();
1612  int numCols = blocked_opr->productDomain()->numBlocks();
1613  int numMiddle = blocked_opl->productDomain()->numBlocks();
1614 
1615  TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == numMiddle);
1616 
1617  RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_product = Thyra::defaultBlockedLinearOp<double>();
1618  blocked_product->beginBlockFill(numRows,numCols);
1619  for(int r = 0; r < numRows; ++r){
1620  for(int c = 0; c < numCols; ++c){
1621 
1622  LinearOp product_rc = explicitMultiply(blocked_opl->getBlock(r,0),blocked_opr->getBlock(0,c));
1623  for(int m = 1; m < numMiddle; ++m){
1624  LinearOp product_m = explicitMultiply(blocked_opl->getBlock(r,m),blocked_opr->getBlock(m,c));
1625  product_rc = explicitAdd(product_rc,product_m);
1626  }
1627  blocked_product->setBlock(r,c,Thyra::scale(scalar,product_rc));
1628  }
1629  }
1630  blocked_product->endBlockFill();
1631  return blocked_product;
1632  }
1633 
1634  // only left factor blocked
1635  if((isBlockedL && !isBlockedR)){
1636  double scalarl = 0.0;
1637  bool transpl = false;
1638  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opl = getPhysicallyBlockedLinearOp(opl,&scalarl,&transpl);
1639  double scalar = scalarl;
1640 
1641  int numRows = blocked_opl->productRange()->numBlocks();
1642  int numCols = 1;
1643  int numMiddle = 1;
1644 
1645  TEUCHOS_ASSERT(blocked_opl->productDomain()->numBlocks() == numMiddle);
1646 
1647  RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_product = Thyra::defaultBlockedLinearOp<double>();
1648  blocked_product->beginBlockFill(numRows,numCols);
1649  for(int r = 0; r < numRows; ++r){
1650  LinearOp product_r = explicitMultiply(blocked_opl->getBlock(r,0),opr);
1651  blocked_product->setBlock(r,0,Thyra::scale(scalar,product_r));
1652  }
1653  blocked_product->endBlockFill();
1654  return blocked_product;
1655  }
1656 
1657  // only right factor blocked
1658  if((!isBlockedL && isBlockedR)){
1659  double scalarr = 0.0;
1660  bool transpr = false;
1661  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opr = getPhysicallyBlockedLinearOp(opr,&scalarr,&transpr);
1662  double scalar = scalarr;
1663 
1664  int numRows = 1;
1665  int numCols = blocked_opr->productDomain()->numBlocks();
1666  int numMiddle = 1;
1667 
1668  TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == numMiddle);
1669 
1670  RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_product = Thyra::defaultBlockedLinearOp<double>();
1671  blocked_product->beginBlockFill(numRows,numCols);
1672  for(int c = 0; c < numCols; ++c){
1673  LinearOp product_c = explicitMultiply(opl,blocked_opr->getBlock(0,c));
1674  blocked_product->setBlock(0,c,Thyra::scale(scalar,product_c));
1675  }
1676  blocked_product->endBlockFill();
1677  return blocked_product;
1678  }
1679 
1680  bool isTpetral = Teko::TpetraHelpers::isTpetraLinearOp(opl);
1681  bool isTpetrar = Teko::TpetraHelpers::isTpetraLinearOp(opr);
1682 
1683  if(isTpetral && isTpetrar){ // Both operators are Tpetra matrices so use the explicit Tpetra matrix-matrix multiply
1684 
1685  // Get left and right Tpetra crs operators
1686  ST scalarl = 0.0;
1687  bool transpl = false;
1688  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpl = Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1689  ST scalarr = 0.0;
1690  bool transpr = false;
1691  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpr = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1692 
1693  // Build output operator
1694  RCP<Thyra::LinearOpBase<ST> > explicitOp;
1695  if(destOp!=Teuchos::null)
1696  explicitOp = destOp;
1697  else
1698  explicitOp = rcp(new Thyra::TpetraLinearOp<ST,LO,GO,NT>());
1699  RCP<Thyra::TpetraLinearOp<ST,LO,GO,NT> > tExplicitOp = rcp_dynamic_cast<Thyra::TpetraLinearOp<ST,LO,GO,NT> >(explicitOp);
1700 
1701  // Do explicit matrix-matrix multiply
1702  RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > explicitCrsOp = Tpetra::createCrsMatrix<ST,LO,GO,NT>(tCrsOpl->getRowMap());
1703  Tpetra::MatrixMatrix::Multiply<ST,LO,GO,NT>(*tCrsOpl,transpl,*tCrsOpr,transpr,*explicitCrsOp);
1704  explicitCrsOp->resumeFill();
1705  explicitCrsOp->scale(scalarl*scalarr);
1706  explicitCrsOp->fillComplete(tCrsOpr->getDomainMap(),tCrsOpl->getRangeMap());
1707  tExplicitOp->initialize(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1708  return tExplicitOp;
1709 
1710  } else if (isTpetral && !isTpetrar){ // Assume that the right operator is diagonal
1711 
1712  // Get left Tpetra crs operator
1713  ST scalarl = 0.0;
1714  bool transpl = false;
1715  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpl = Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1716 
1717  // Cast right operator as DiagonalLinearOp and extract diagonal as Vector
1718  RCP<const Thyra::DiagonalLinearOpBase<ST> > dOpr = rcp_dynamic_cast<const Thyra::DiagonalLinearOpBase<ST> >(opr);
1719  RCP<const Thyra::TpetraVector<ST,LO,GO,NT> > tPtr = rcp_dynamic_cast<const Thyra::TpetraVector<ST,LO,GO,NT> >(dOpr->getDiag(),true);
1720  RCP<const Tpetra::Vector<ST,LO,GO,NT> > diagPtr = rcp_dynamic_cast<const Tpetra::Vector<ST,LO,GO,NT> >(tPtr->getConstTpetraVector(),true);
1721 
1722  // Scale by the diagonal operator
1723  RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > explicitCrsOp = Tpetra::importAndFillCompleteCrsMatrix<Tpetra::CrsMatrix<ST,LO,GO,NT> >(tCrsOpl, Tpetra::Import<LO,GO,NT>(tCrsOpl->getRowMap(),tCrsOpl->getRowMap()));
1724  explicitCrsOp->rightScale(*diagPtr);
1725  explicitCrsOp->resumeFill();
1726  explicitCrsOp->scale(scalarl);
1727  explicitCrsOp->fillComplete(tCrsOpl->getDomainMap(),tCrsOpl->getRangeMap());
1728  return Thyra::tpetraLinearOp<ST,LO,GO,NT>(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1729 
1730  } else if (!isTpetral && isTpetrar){ // Assume that the left operator is diagonal
1731 
1732  // Get right Tpetra crs operator
1733  ST scalarr = 0.0;
1734  bool transpr = false;
1735  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpr = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1736 
1737  // Cast leftt operator as DiagonalLinearOp and extract diagonal as Vector
1738  RCP<const Thyra::DiagonalLinearOpBase<ST> > dOpl = rcp_dynamic_cast<const Thyra::DiagonalLinearOpBase<ST> >(opl,true);
1739  RCP<const Thyra::TpetraVector<ST,LO,GO,NT> > tPtr = rcp_dynamic_cast<const Thyra::TpetraVector<ST,LO,GO,NT> >(dOpl->getDiag(),true);
1740  RCP<const Tpetra::Vector<ST,LO,GO,NT> > diagPtr = rcp_dynamic_cast<const Tpetra::Vector<ST,LO,GO,NT> >(tPtr->getConstTpetraVector(),true);
1741 
1742  // Scale by the diagonal operator
1743  RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > explicitCrsOp = Tpetra::importAndFillCompleteCrsMatrix<Tpetra::CrsMatrix<ST,LO,GO,NT> >(tCrsOpr, Tpetra::Import<LO,GO,NT>(tCrsOpr->getRowMap(),tCrsOpr->getRowMap()));
1744  explicitCrsOp->leftScale(*diagPtr);
1745  explicitCrsOp->resumeFill();
1746  explicitCrsOp->scale(scalarr);
1747  explicitCrsOp->fillComplete(tCrsOpr->getDomainMap(),tCrsOpr->getRangeMap());
1748  return Thyra::tpetraLinearOp<ST,LO,GO,NT>(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1749 
1750  } else { // Assume Epetra and we can use transformers
1751 
1752  // build implicit multiply
1753  const LinearOp implicitOp = Thyra::multiply(opl,opr);
1754 
1755  // build a scaling transformer
1756 
1757  RCP<Thyra::LinearOpTransformerBase<double> > prodTrans
1758  = Thyra::epetraExtDiagScalingTransformer();
1759 
1760  // check to see if a scaling transformer works: if not use the
1761  // DiagScaledMatrixProduct transformer
1762  if(not prodTrans->isCompatible(*implicitOp))
1763  prodTrans = Thyra::epetraExtDiagScaledMatProdTransformer();
1764 
1765  // build operator destination operator
1766  ModifiableLinearOp explicitOp;
1767 
1768  // if neccessary build a operator to put the explicit multiply into
1769  if(destOp==Teuchos::null)
1770  explicitOp = prodTrans->createOutputOp();
1771  else
1772  explicitOp = destOp;
1773 
1774  // perform multiplication
1775  prodTrans->transform(*implicitOp,explicitOp.ptr());
1776 
1777  // label it
1778  explicitOp->setObjectLabel("explicit( " + opl->getObjectLabel() +
1779  " * " + opr->getObjectLabel() + " )");
1780 
1781  return explicitOp;
1782  }
1783 }
1784 
1795 const LinearOp explicitAdd(const LinearOp & opl_in,const LinearOp & opr_in)
1796 {
1797  // if both blocked, add block by block
1798  if(isPhysicallyBlockedLinearOp(opl_in) && isPhysicallyBlockedLinearOp(opr_in)){
1799  double scalarl = 0.0;
1800  bool transpl = false;
1801  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opl = getPhysicallyBlockedLinearOp(opl_in, &scalarl, &transpl);
1802 
1803  double scalarr = 0.0;
1804  bool transpr = false;
1805  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opr = getPhysicallyBlockedLinearOp(opr_in, &scalarr, &transpr);
1806 
1807  int numRows = blocked_opl->productRange()->numBlocks();
1808  int numCols = blocked_opl->productDomain()->numBlocks();
1809  TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == numRows);
1810  TEUCHOS_ASSERT(blocked_opr->productDomain()->numBlocks() == numCols);
1811 
1812  RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_sum = Thyra::defaultBlockedLinearOp<double>();
1813  blocked_sum->beginBlockFill(numRows,numCols);
1814  for(int r = 0; r < numRows; ++r)
1815  for(int c = 0; c < numCols; ++c)
1816  blocked_sum->setBlock(r,c,explicitAdd(Thyra::scale(scalarl,blocked_opl->getBlock(r,c)),Thyra::scale(scalarr,blocked_opr->getBlock(r,c))));
1817  blocked_sum->endBlockFill();
1818  return blocked_sum;
1819  }
1820 
1821  // if only one is blocked, it must be 1x1
1822  LinearOp opl = opl_in;
1823  LinearOp opr = opr_in;
1824  if(isPhysicallyBlockedLinearOp(opl_in)){
1825  double scalarl = 0.0;
1826  bool transpl = false;
1827  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opl = getPhysicallyBlockedLinearOp(opl_in, &scalarl, &transpl);
1828  TEUCHOS_ASSERT(blocked_opl->productRange()->numBlocks() == 1);
1829  TEUCHOS_ASSERT(blocked_opl->productDomain()->numBlocks() == 1);
1830  opl = Thyra::scale(scalarl,blocked_opl->getBlock(0,0));
1831  }
1832  if(isPhysicallyBlockedLinearOp(opr_in)){
1833  double scalarr = 0.0;
1834  bool transpr = false;
1835  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opr = getPhysicallyBlockedLinearOp(opr_in, &scalarr, &transpr);
1836  TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == 1);
1837  TEUCHOS_ASSERT(blocked_opr->productDomain()->numBlocks() == 1);
1838  opr = Thyra::scale(scalarr,blocked_opr->getBlock(0,0));
1839  }
1840 
1841  bool isTpetral = Teko::TpetraHelpers::isTpetraLinearOp(opl);
1842  bool isTpetrar = Teko::TpetraHelpers::isTpetraLinearOp(opr);
1843 
1844  // if one of the operators in the sum is a thyra zero op
1845  if(isZeroOp(opl)){
1846  if(isZeroOp(opr))
1847  return opr; // return a zero op if both are zero
1848  if(isTpetrar){ // if other op is tpetra, replace this with a zero crs matrix
1849  ST scalar = 0.0;
1850  bool transp = false;
1851  auto crs_op = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalar, &transp);
1852  auto zero_crs = Tpetra::createCrsMatrix<ST,LO,GO,NT>(crs_op->getRowMap());
1853  zero_crs->fillComplete();
1854  opl = Thyra::constTpetraLinearOp<ST,LO,GO,NT>(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(crs_op->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(crs_op->getDomainMap()),zero_crs);
1855  isTpetral = true;
1856  } else
1857  return opr->clone();
1858  }
1859  if(isZeroOp(opr)){
1860  if(isTpetral){ // if other op is tpetra, replace this with a zero crs matrix
1861  ST scalar = 0.0;
1862  bool transp = false;
1863  auto crs_op = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalar, &transp);
1864  auto zero_crs = Tpetra::createCrsMatrix<ST,LO,GO,NT>(crs_op->getRowMap());
1865  zero_crs->fillComplete();
1866  opr = Thyra::constTpetraLinearOp<ST,LO,GO,NT>(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(crs_op->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(crs_op->getDomainMap()),zero_crs);
1867  isTpetrar = true;
1868  } else
1869  return opl->clone();
1870  }
1871 
1872  if(isTpetral && isTpetrar){ // Both operators are Tpetra matrices so use the explicit Tpetra matrix-matrix add
1873 
1874  // Get left and right Tpetra crs operators
1875  ST scalarl = 0.0;
1876  bool transpl = false;
1877  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpl = Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1878  ST scalarr = 0.0;
1879  bool transpr = false;
1880  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpr = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1881 
1882  // Build output operator
1883  RCP<Thyra::LinearOpBase<ST> > explicitOp = rcp(new Thyra::TpetraLinearOp<ST,LO,GO,NT>());
1884  RCP<Thyra::TpetraLinearOp<ST,LO,GO,NT> > tExplicitOp = rcp_dynamic_cast<Thyra::TpetraLinearOp<ST,LO,GO,NT> >(explicitOp);
1885 
1886  // Do explicit matrix-matrix add
1887  RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > explicitCrsOp = Tpetra::MatrixMatrix::add<ST,LO,GO,NT>(scalarl,transpl,*tCrsOpl,scalarr,transpr,*tCrsOpr);
1888  tExplicitOp->initialize(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1889  return tExplicitOp;
1890 
1891  }else{//Assume Epetra
1892  // build implicit add
1893  const LinearOp implicitOp = Thyra::add(opl,opr);
1894 
1895  // build transformer
1896  const RCP<Thyra::LinearOpTransformerBase<double> > prodTrans =
1897  Thyra::epetraExtAddTransformer();
1898 
1899  // build operator and add
1900  const RCP<Thyra::LinearOpBase<double> > explicitOp = prodTrans->createOutputOp();
1901  prodTrans->transform(*implicitOp,explicitOp.ptr());
1902  explicitOp->setObjectLabel("explicit( " + opl->getObjectLabel() +
1903  " + " + opr->getObjectLabel() + " )");
1904 
1905  return explicitOp;
1906  }
1907 }
1908 
1921 const ModifiableLinearOp explicitAdd(const LinearOp & opl,const LinearOp & opr,
1922  const ModifiableLinearOp & destOp)
1923 {
1924  // if blocked, add block by block
1925  if(isPhysicallyBlockedLinearOp(opl)){
1926  TEUCHOS_ASSERT(isPhysicallyBlockedLinearOp(opr));
1927 
1928  double scalarl = 0.0;
1929  bool transpl = false;
1930  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opl = getPhysicallyBlockedLinearOp(opl, &scalarl, &transpl);
1931 
1932  double scalarr = 0.0;
1933  bool transpr = false;
1934  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opr = getPhysicallyBlockedLinearOp(opr, &scalarr, &transpr);
1935 
1936  int numRows = blocked_opl->productRange()->numBlocks();
1937  int numCols = blocked_opl->productDomain()->numBlocks();
1938  TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == numRows);
1939  TEUCHOS_ASSERT(blocked_opr->productDomain()->numBlocks() == numCols);
1940 
1941  RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_sum = Thyra::defaultBlockedLinearOp<double>();
1942  blocked_sum->beginBlockFill(numRows,numCols);
1943  for(int r = 0; r < numRows; ++r)
1944  for(int c = 0; c < numCols; ++c)
1945  blocked_sum->setBlock(r,c,explicitAdd(Thyra::scale(scalarl,blocked_opl->getBlock(r,c)),Thyra::scale(scalarr,blocked_opr->getBlock(r,c))));
1946  blocked_sum->endBlockFill();
1947  return blocked_sum;
1948  }
1949 
1950  bool isTpetral = Teko::TpetraHelpers::isTpetraLinearOp(opl);
1951  bool isTpetrar = Teko::TpetraHelpers::isTpetraLinearOp(opr);
1952 
1953  if(isTpetral && isTpetrar){ // Both operators are Tpetra matrices so use the explicit Tpetra matrix-matrix add
1954 
1955  // Get left and right Tpetra crs operators
1956  ST scalarl = 0.0;
1957  bool transpl = false;
1958  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpl = Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1959  ST scalarr = 0.0;
1960  bool transpr = false;
1961  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpr = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1962 
1963  // Build output operator
1964  RCP<Thyra::LinearOpBase<ST> > explicitOp;
1965  if(destOp!=Teuchos::null)
1966  explicitOp = destOp;
1967  else
1968  explicitOp = rcp(new Thyra::TpetraLinearOp<ST,LO,GO,NT>());
1969  RCP<Thyra::TpetraLinearOp<ST,LO,GO,NT> > tExplicitOp = rcp_dynamic_cast<Thyra::TpetraLinearOp<ST,LO,GO,NT> >(explicitOp);
1970 
1971  // Do explicit matrix-matrix add
1972  RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > explicitCrsOp = Tpetra::MatrixMatrix::add<ST,LO,GO,NT>(scalarl,transpl,*tCrsOpl,scalarr,transpr,*tCrsOpr);
1973  tExplicitOp->initialize(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1974  return tExplicitOp;
1975 
1976  }else{ // Assume Epetra
1977 
1978  // build implicit add
1979  const LinearOp implicitOp = Thyra::add(opl,opr);
1980 
1981  // build transformer
1982  const RCP<Thyra::LinearOpTransformerBase<double> > prodTrans =
1983  Thyra::epetraExtAddTransformer();
1984 
1985  // build or reuse destination operator
1986  RCP<Thyra::LinearOpBase<double> > explicitOp;
1987  if(destOp!=Teuchos::null)
1988  explicitOp = destOp;
1989  else
1990  explicitOp = prodTrans->createOutputOp();
1991 
1992  // perform add
1993  prodTrans->transform(*implicitOp,explicitOp.ptr());
1994  explicitOp->setObjectLabel("explicit( " + opl->getObjectLabel() +
1995  " + " + opr->getObjectLabel() + " )");
1996 
1997  return explicitOp;
1998  }
1999 }
2000 
2005 const ModifiableLinearOp explicitSum(const LinearOp & op,
2006  const ModifiableLinearOp & destOp)
2007 {
2008  // convert operators to Epetra_CrsMatrix
2009  const RCP<const Epetra_CrsMatrix> epetraOp =
2010  rcp_dynamic_cast<const Epetra_CrsMatrix>(get_Epetra_Operator(*op), true);
2011 
2012  if(destOp==Teuchos::null) {
2013  Teuchos::RCP<Epetra_Operator> epetraDest = Teuchos::rcp(new Epetra_CrsMatrix(*epetraOp));
2014 
2015  return Thyra::nonconstEpetraLinearOp(epetraDest);
2016  }
2017 
2018  const RCP<Epetra_CrsMatrix> epetraDest =
2019  rcp_dynamic_cast<Epetra_CrsMatrix>(get_Epetra_Operator(*destOp), true);
2020 
2021  EpetraExt::MatrixMatrix::Add(*epetraOp,false,1.0,*epetraDest,1.0);
2022 
2023  return destOp;
2024 }
2025 
2026 const LinearOp explicitTranspose(const LinearOp & op)
2027 {
2028  if(Teko::TpetraHelpers::isTpetraLinearOp(op)){
2029 
2030  RCP<const Thyra::TpetraLinearOp<ST,LO,GO,NT> > tOp = rcp_dynamic_cast<const Thyra::TpetraLinearOp<ST,LO,GO,NT> >(op,true);
2031  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOp = rcp_dynamic_cast<const Tpetra::CrsMatrix<ST,LO,GO,NT> >(tOp->getConstTpetraOperator(),true);
2032 
2033  Tpetra::RowMatrixTransposer<ST,LO,GO,NT> transposer(tCrsOp);
2034  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > transOp = transposer.createTranspose();
2035 
2036  return Thyra::constTpetraLinearOp<ST,LO,GO,NT>(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(transOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(transOp->getDomainMap()),transOp);
2037 
2038  } else {
2039 
2040  RCP<const Epetra_Operator> eOp = Thyra::get_Epetra_Operator(*op);
2041  TEUCHOS_TEST_FOR_EXCEPTION(eOp==Teuchos::null,std::logic_error,
2042  "Teko::explicitTranspose Not an Epetra_Operator");
2043  RCP<const Epetra_RowMatrix> eRowMatrixOp
2044  = Teuchos::rcp_dynamic_cast<const Epetra_RowMatrix>(eOp);
2045  TEUCHOS_TEST_FOR_EXCEPTION(eRowMatrixOp==Teuchos::null,std::logic_error,
2046  "Teko::explicitTranspose Not an Epetra_RowMatrix");
2047 
2048  // we now have a delete transpose operator
2049  EpetraExt::RowMatrix_Transpose tranposeOp;
2050  Epetra_RowMatrix & eMat = tranposeOp(const_cast<Epetra_RowMatrix &>(*eRowMatrixOp));
2051 
2052  // this copy is because of a poor implementation of the EpetraExt::Transform
2053  // implementation
2054  Teuchos::RCP<Epetra_CrsMatrix> crsMat
2055  = Teuchos::rcp(new Epetra_CrsMatrix(dynamic_cast<Epetra_CrsMatrix &>(eMat)));
2056 
2057  return Thyra::epetraLinearOp(crsMat);
2058  }
2059 }
2060 
2061 double frobeniusNorm(const LinearOp & op_in)
2062 {
2063  LinearOp op;
2064  double scalar = 1.0;
2065 
2066  // if blocked, must be 1x1
2067  if(isPhysicallyBlockedLinearOp(op_in)){
2068  bool transp = false;
2069  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_op = getPhysicallyBlockedLinearOp(op_in,&scalar,&transp);
2070  TEUCHOS_ASSERT(blocked_op->productRange()->numBlocks() == 1);
2071  TEUCHOS_ASSERT(blocked_op->productDomain()->numBlocks() == 1);
2072  op = blocked_op->getBlock(0,0);
2073  } else
2074  op = op_in;
2075 
2076  if(Teko::TpetraHelpers::isTpetraLinearOp(op)){
2077  const RCP<const Thyra::TpetraLinearOp<ST,LO,GO,NT> > tOp = rcp_dynamic_cast<const Thyra::TpetraLinearOp<ST,LO,GO,NT> >(op);
2078  const RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > crsOp = rcp_dynamic_cast<const Tpetra::CrsMatrix<ST,LO,GO,NT> >(tOp->getConstTpetraOperator(),true);
2079  return crsOp->getFrobeniusNorm();
2080  } else {
2081  const RCP<const Epetra_Operator> epOp = Thyra::get_Epetra_Operator(*op);
2082  const RCP<const Epetra_CrsMatrix> crsOp = rcp_dynamic_cast<const Epetra_CrsMatrix>(epOp,true);
2083  return crsOp->NormFrobenius();
2084  }
2085 }
2086 
2087 double oneNorm(const LinearOp & op)
2088 {
2089  if(Teko::TpetraHelpers::isTpetraLinearOp(op)){
2090  TEUCHOS_TEST_FOR_EXCEPTION(true,std::logic_error,"One norm not currently implemented for Tpetra matrices");
2091 
2092  } else {
2093  const RCP<const Epetra_Operator> epOp = Thyra::get_Epetra_Operator(*op);
2094  const RCP<const Epetra_CrsMatrix> crsOp = rcp_dynamic_cast<const Epetra_CrsMatrix>(epOp,true);
2095  return crsOp->NormOne();
2096  }
2097 }
2098 
2099 double infNorm(const LinearOp & op)
2100 {
2101  if(Teko::TpetraHelpers::isTpetraLinearOp(op)){
2102  ST scalar = 0.0;
2103  bool transp = false;
2104  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOp = Teko::TpetraHelpers::getTpetraCrsMatrix(op, &scalar, &transp);
2105 
2106  // extract diagonal
2107  const RCP<Tpetra::Vector<ST,LO,GO,NT> > ptrDiag = Tpetra::createVector<ST,LO,GO,NT>(tCrsOp->getRowMap());
2108  Tpetra::Vector<ST,LO,GO,NT> & diag = *ptrDiag;
2109 
2110  // compute absolute value row sum
2111  diag.putScalar(0.0);
2112  for(LO i=0;i<(LO) tCrsOp->getNodeNumRows();i++) {
2113  LO numEntries = tCrsOp->getNumEntriesInLocalRow (i);
2114  std::vector<LO> indices(numEntries);
2115  std::vector<ST> values(numEntries);
2116  Teuchos::ArrayView<const LO> indices_av(indices);
2117  Teuchos::ArrayView<const ST> values_av(values);
2118  tCrsOp->getLocalRowView(i,indices_av,values_av);
2119 
2120  // build abs value row sum
2121  for(LO j=0;j<numEntries;j++)
2122  diag.sumIntoLocalValue(i,std::abs(values_av[j]));
2123  }
2124  return diag.normInf()*scalar;
2125 
2126  } else {
2127  const RCP<const Epetra_Operator> epOp = Thyra::get_Epetra_Operator(*op);
2128  const RCP<const Epetra_CrsMatrix> crsOp = rcp_dynamic_cast<const Epetra_CrsMatrix>(epOp,true);
2129  return crsOp->NormInf();
2130  }
2131 }
2132 
2133 const LinearOp buildDiagonal(const MultiVector & src,const std::string & lbl)
2134 {
2135  RCP<Thyra::VectorBase<double> > dst = Thyra::createMember(src->range());
2136  Thyra::copy(*src->col(0),dst.ptr());
2137 
2138  return Thyra::diagonal<double>(dst,lbl);
2139 }
2140 
2141 const LinearOp buildInvDiagonal(const MultiVector & src,const std::string & lbl)
2142 {
2143  const RCP<const Thyra::VectorBase<double> > srcV = src->col(0);
2144  RCP<Thyra::VectorBase<double> > dst = Thyra::createMember(srcV->range());
2145  Thyra::reciprocal<double>(*srcV,dst.ptr());
2146 
2147  return Thyra::diagonal<double>(dst,lbl);
2148 }
2149 
2151 BlockedMultiVector buildBlockedMultiVector(const std::vector<MultiVector> & mvv)
2152 {
2153  Teuchos::Array<MultiVector> mvA;
2154  Teuchos::Array<VectorSpace> vsA;
2155 
2156  // build arrays of multi vectors and vector spaces
2157  std::vector<MultiVector>::const_iterator itr;
2158  for(itr=mvv.begin();itr!=mvv.end();++itr) {
2159  mvA.push_back(*itr);
2160  vsA.push_back((*itr)->range());
2161  }
2162 
2163  // construct the product vector space
2164  const RCP<const Thyra::DefaultProductVectorSpace<double> > vs
2165  = Thyra::productVectorSpace<double>(vsA);
2166 
2167  return Thyra::defaultProductMultiVector<double>(vs,mvA);
2168 }
2169 
2180 Teuchos::RCP<Thyra::VectorBase<double> > indicatorVector(
2181  const std::vector<int> & indices,
2182  const VectorSpace & vs,
2183  double onValue, double offValue)
2184 
2185 {
2186  using Teuchos::RCP;
2187 
2188  // create a new vector
2189  RCP<Thyra::VectorBase<double> > v = Thyra::createMember<double>(vs);
2190  Thyra::put_scalar<double>(offValue,v.ptr()); // fill it with "off" values
2191 
2192  // set on values
2193  for(std::size_t i=0;i<indices.size();i++)
2194  Thyra::set_ele<double>(indices[i],onValue,v.ptr());
2195 
2196  return v;
2197 }
2198 
2223 double computeSpectralRad(const RCP<const Thyra::LinearOpBase<double> > & A, double tol,
2224  bool isHermitian,int numBlocks,int restart,int verbosity)
2225 {
2226  typedef Thyra::LinearOpBase<double> OP;
2227  typedef Thyra::MultiVectorBase<double> MV;
2228 
2229  int startVectors = 1;
2230 
2231  // construct an initial guess
2232  const RCP<MV> ivec = Thyra::createMember(A->domain());
2233  Thyra::randomize(-1.0,1.0,ivec.ptr());
2234 
2235  RCP<Anasazi::BasicEigenproblem<double,MV,OP> > eigProb
2236  = rcp(new Anasazi::BasicEigenproblem<double,MV,OP>(A,ivec));
2237  eigProb->setNEV(1);
2238  eigProb->setHermitian(isHermitian);
2239 
2240  // set the problem up
2241  if(not eigProb->setProblem()) {
2242  // big time failure!
2243  return Teuchos::ScalarTraits<double>::nan();
2244  }
2245 
2246  // we want largert magnitude eigenvalue
2247  std::string which("LM"); // largest magnitude
2248 
2249  // Create the parameter list for the eigensolver
2250  // verbosity+=Anasazi::TimingDetails;
2251  Teuchos::ParameterList MyPL;
2252  MyPL.set( "Verbosity", verbosity );
2253  MyPL.set( "Which", which );
2254  MyPL.set( "Block Size", startVectors );
2255  MyPL.set( "Num Blocks", numBlocks );
2256  MyPL.set( "Maximum Restarts", restart );
2257  MyPL.set( "Convergence Tolerance", tol );
2258 
2259  // build status test manager
2260  // RCP<Anasazi::StatusTestMaxIters<double,MV,OP> > statTest
2261  // = rcp(new Anasazi::StatusTestMaxIters<double,MV,OP>(numBlocks*(restart+1)));
2262 
2263  // Create the Block Krylov Schur solver
2264  // This takes as inputs the eigenvalue problem and the solver parameters
2265  Anasazi::BlockKrylovSchurSolMgr<double,MV,OP> MyBlockKrylovSchur(eigProb, MyPL );
2266 
2267  // Solve the eigenvalue problem, and save the return code
2268  Anasazi::ReturnType solverreturn = MyBlockKrylovSchur.solve();
2269 
2270  if(solverreturn==Anasazi::Unconverged) {
2271  double real = MyBlockKrylovSchur.getRitzValues().begin()->realpart;
2272  double comp = MyBlockKrylovSchur.getRitzValues().begin()->imagpart;
2273 
2274  return -std::abs(std::sqrt(real*real+comp*comp));
2275 
2276  // cout << "Anasazi::BlockKrylovSchur::solve() did not converge!" << std::endl;
2277  // return -std::abs(MyBlockKrylovSchur.getRitzValues().begin()->realpart);
2278  }
2279  else { // solverreturn==Anasazi::Converged
2280  double real = eigProb->getSolution().Evals.begin()->realpart;
2281  double comp = eigProb->getSolution().Evals.begin()->imagpart;
2282 
2283  return std::abs(std::sqrt(real*real+comp*comp));
2284 
2285  // cout << "Anasazi::BlockKrylovSchur::solve() converged!" << endl;
2286  // return std::abs(eigProb->getSolution().Evals.begin()->realpart);
2287  }
2288 }
2289 
2313 double computeSmallestMagEig(const RCP<const Thyra::LinearOpBase<double> > & A, double tol,
2314  bool isHermitian,int numBlocks,int restart,int verbosity)
2315 {
2316  typedef Thyra::LinearOpBase<double> OP;
2317  typedef Thyra::MultiVectorBase<double> MV;
2318 
2319  int startVectors = 1;
2320 
2321  // construct an initial guess
2322  const RCP<MV> ivec = Thyra::createMember(A->domain());
2323  Thyra::randomize(-1.0,1.0,ivec.ptr());
2324 
2325  RCP<Anasazi::BasicEigenproblem<double,MV,OP> > eigProb
2326  = rcp(new Anasazi::BasicEigenproblem<double,MV,OP>(A,ivec));
2327  eigProb->setNEV(1);
2328  eigProb->setHermitian(isHermitian);
2329 
2330  // set the problem up
2331  if(not eigProb->setProblem()) {
2332  // big time failure!
2333  return Teuchos::ScalarTraits<double>::nan();
2334  }
2335 
2336  // we want largert magnitude eigenvalue
2337  std::string which("SM"); // smallest magnitude
2338 
2339  // Create the parameter list for the eigensolver
2340  Teuchos::ParameterList MyPL;
2341  MyPL.set( "Verbosity", verbosity );
2342  MyPL.set( "Which", which );
2343  MyPL.set( "Block Size", startVectors );
2344  MyPL.set( "Num Blocks", numBlocks );
2345  MyPL.set( "Maximum Restarts", restart );
2346  MyPL.set( "Convergence Tolerance", tol );
2347 
2348  // build status test manager
2349  // RCP<Anasazi::StatusTestMaxIters<double,MV,OP> > statTest
2350  // = rcp(new Anasazi::StatusTestMaxIters<double,MV,OP>(10));
2351 
2352  // Create the Block Krylov Schur solver
2353  // This takes as inputs the eigenvalue problem and the solver parameters
2354  Anasazi::BlockKrylovSchurSolMgr<double,MV,OP> MyBlockKrylovSchur(eigProb, MyPL );
2355 
2356  // Solve the eigenvalue problem, and save the return code
2357  Anasazi::ReturnType solverreturn = MyBlockKrylovSchur.solve();
2358 
2359  if(solverreturn==Anasazi::Unconverged) {
2360  // cout << "Anasazi::BlockKrylovSchur::solve() did not converge! " << std::endl;
2361  return -std::abs(MyBlockKrylovSchur.getRitzValues().begin()->realpart);
2362  }
2363  else { // solverreturn==Anasazi::Converged
2364  // cout << "Anasazi::BlockKrylovSchur::solve() converged!" << endl;
2365  return std::abs(eigProb->getSolution().Evals.begin()->realpart);
2366  }
2367 }
2368 
2377 ModifiableLinearOp getDiagonalOp(const Teko::LinearOp & A,const DiagonalType & dt)
2378 {
2379  switch(dt) {
2380  case Diagonal:
2381  return getDiagonalOp(A);
2382  case Lumped:
2383  return getLumpedMatrix(A);
2384  case AbsRowSum:
2385  return getAbsRowSumMatrix(A);
2386  case NotDiag:
2387  default:
2388  TEUCHOS_TEST_FOR_EXCEPT(true);
2389  };
2390 
2391  return Teuchos::null;
2392 }
2393 
2402 ModifiableLinearOp getInvDiagonalOp(const Teko::LinearOp & A,const Teko::DiagonalType & dt)
2403 {
2404  switch(dt) {
2405  case Diagonal:
2406  return getInvDiagonalOp(A);
2407  case Lumped:
2408  return getInvLumpedMatrix(A);
2409  case AbsRowSum:
2410  return getAbsRowSumInvMatrix(A);
2411  case NotDiag:
2412  default:
2413  TEUCHOS_TEST_FOR_EXCEPT(true);
2414  };
2415 
2416  return Teuchos::null;
2417 }
2418 
2425 std::string getDiagonalName(const DiagonalType & dt)
2426 {
2427  switch(dt) {
2428  case Diagonal:
2429  return "Diagonal";
2430  case Lumped:
2431  return "Lumped";
2432  case AbsRowSum:
2433  return "AbsRowSum";
2434  case NotDiag:
2435  return "NotDiag";
2436  case BlkDiag:
2437  return "BlkDiag";
2438  };
2439 
2440  return "<error>";
2441 }
2442 
2451 DiagonalType getDiagonalType(std::string name)
2452 {
2453  if(name=="Diagonal")
2454  return Diagonal;
2455  if(name=="Lumped")
2456  return Lumped;
2457  if(name=="AbsRowSum")
2458  return AbsRowSum;
2459  if(name=="BlkDiag")
2460  return BlkDiag;
2461 
2462  return NotDiag;
2463 }
2464 
2465 LinearOp probe(Teuchos::RCP<const Epetra_CrsGraph> &G,const LinearOp & Op){
2466 #ifdef Teko_ENABLE_Isorropia
2467  Teuchos::ParameterList probeList;
2468  Prober prober(G,probeList,true);
2469  Teuchos::RCP<Epetra_CrsMatrix> Mat=rcp(new Epetra_CrsMatrix(Copy,*G));
2471  prober.probe(Mwrap,*Mat);
2472  return Thyra::epetraLinearOp(Mat);
2473 #else
2474  TEUCHOS_TEST_FOR_EXCEPTION(true,std::runtime_error,"Probe requires Isorropia");
2475 #endif
2476 }
2477 
2478 double norm_1(const MultiVector & v,std::size_t col)
2479 {
2480  Teuchos::Array<double> n(v->domain()->dim());
2481  Thyra::norms_1<double>(*v,n);
2482 
2483  return n[col];
2484 }
2485 
2486 double norm_2(const MultiVector & v,std::size_t col)
2487 {
2488  Teuchos::Array<double> n(v->domain()->dim());
2489  Thyra::norms_2<double>(*v,n);
2490 
2491  return n[col];
2492 }
2493 
2494 void putScalar(const ModifiableLinearOp & op,double scalar)
2495 {
2496  try {
2497  // get Epetra_Operator
2498  RCP<Epetra_Operator> eOp = Thyra::get_Epetra_Operator(*op);
2499 
2500  // cast it to a CrsMatrix
2501  RCP<Epetra_CrsMatrix> eCrsOp = rcp_dynamic_cast<Epetra_CrsMatrix>(eOp,true);
2502 
2503  eCrsOp->PutScalar(scalar);
2504  }
2505  catch (std::exception & e) {
2506  RCP<Teuchos::FancyOStream> out = Teuchos::VerboseObjectBase::getDefaultOStream();
2507 
2508  *out << "Teko: putScalar requires an Epetra_CrsMatrix\n";
2509  *out << " Could not extract an Epetra_Operator from a \"" << op->description() << std::endl;
2510  *out << " OR\n";
2511  *out << " Could not cast an Epetra_Operator to a Epetra_CrsMatrix\n";
2512  *out << std::endl;
2513  *out << "*** THROWN EXCEPTION ***\n";
2514  *out << e.what() << std::endl;
2515  *out << "************************\n";
2516 
2517  throw e;
2518  }
2519 }
2520 
2521 void clipLower(MultiVector & v,double lowerBound)
2522 {
2523  using Teuchos::RCP;
2524  using Teuchos::rcp_dynamic_cast;
2525 
2526  // cast so entries are accessible
2527  // RCP<Thyra::SpmdMultiVectorBase<double> > spmdMVec
2528  // = rcp_dynamic_cast<Thyra::DefaultSpmdMultiVector<double> >(v);
2529 
2530  for(Thyra::Ordinal i=0;i<v->domain()->dim();i++) {
2531  RCP<Thyra::SpmdVectorBase<double> > spmdVec
2532  = rcp_dynamic_cast<Thyra::SpmdVectorBase<double> >(v->col(i),true);
2533 
2534  Teuchos::ArrayRCP<double> values;
2535  // spmdMVec->getNonconstLocalData(Teuchos::ptrFromRef(values),Teuchos::ptrFromRef(i));
2536  spmdVec->getNonconstLocalData(Teuchos::ptrFromRef(values));
2537  for(Teuchos::ArrayRCP<double>::size_type j=0;j<values.size();j++) {
2538  if(values[j]<lowerBound)
2539  values[j] = lowerBound;
2540  }
2541  }
2542 }
2543 
2544 void clipUpper(MultiVector & v,double upperBound)
2545 {
2546  using Teuchos::RCP;
2547  using Teuchos::rcp_dynamic_cast;
2548 
2549  // cast so entries are accessible
2550  // RCP<Thyra::SpmdMultiVectorBase<double> > spmdMVec
2551  // = rcp_dynamic_cast<Thyra::DefaultSpmdMultiVector<double> >(v);
2552  for(Thyra::Ordinal i=0;i<v->domain()->dim();i++) {
2553  RCP<Thyra::SpmdVectorBase<double> > spmdVec
2554  = rcp_dynamic_cast<Thyra::SpmdVectorBase<double> >(v->col(i),true);
2555 
2556  Teuchos::ArrayRCP<double> values;
2557  // spmdMVec->getNonconstLocalData(Teuchos::ptrFromRef(values),Teuchos::ptrFromRef(i));
2558  spmdVec->getNonconstLocalData(Teuchos::ptrFromRef(values));
2559  for(Teuchos::ArrayRCP<double>::size_type j=0;j<values.size();j++) {
2560  if(values[j]>upperBound)
2561  values[j] = upperBound;
2562  }
2563  }
2564 }
2565 
2566 void replaceValue(MultiVector & v,double currentValue,double newValue)
2567 {
2568  using Teuchos::RCP;
2569  using Teuchos::rcp_dynamic_cast;
2570 
2571  // cast so entries are accessible
2572  // RCP<Thyra::SpmdMultiVectorBase<double> > spmdMVec
2573  // = rcp_dynamic_cast<Thyra::SpmdMultiVectorBase<double> >(v,true);
2574  for(Thyra::Ordinal i=0;i<v->domain()->dim();i++) {
2575  RCP<Thyra::SpmdVectorBase<double> > spmdVec
2576  = rcp_dynamic_cast<Thyra::SpmdVectorBase<double> >(v->col(i),true);
2577 
2578  Teuchos::ArrayRCP<double> values;
2579  // spmdMVec->getNonconstLocalData(Teuchos::ptrFromRef(values),Teuchos::ptrFromRef(i));
2580  spmdVec->getNonconstLocalData(Teuchos::ptrFromRef(values));
2581  for(Teuchos::ArrayRCP<double>::size_type j=0;j<values.size();j++) {
2582  if(values[j]==currentValue)
2583  values[j] = newValue;
2584  }
2585  }
2586 }
2587 
2588 void columnAverages(const MultiVector & v,std::vector<double> & averages)
2589 {
2590  averages.resize(v->domain()->dim());
2591 
2592  // sum over each column
2593  Thyra::sums<double>(*v,averages);
2594 
2595  // build averages
2596  Thyra::Ordinal rows = v->range()->dim();
2597  for(std::size_t i=0;i<averages.size();i++)
2598  averages[i] = averages[i]/rows;
2599 }
2600 
2601 double average(const MultiVector & v)
2602 {
2603  Thyra::Ordinal rows = v->range()->dim();
2604  Thyra::Ordinal cols = v->domain()->dim();
2605 
2606  std::vector<double> averages;
2607  columnAverages(v,averages);
2608 
2609  double sum = 0.0;
2610  for(std::size_t i=0;i<averages.size();i++)
2611  sum += averages[i] * rows;
2612 
2613  return sum/(rows*cols);
2614 }
2615 
2616 bool isPhysicallyBlockedLinearOp(const LinearOp & op)
2617 {
2618  // See if the operator is a PBLO
2619  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > pblo = rcp_dynamic_cast<const Thyra::PhysicallyBlockedLinearOpBase<double> >(op);
2620  if (!pblo.is_null())
2621  return true;
2622 
2623  // See if the operator is a wrapped PBLO
2624  ST scalar = 0.0;
2625  Thyra::EOpTransp transp = Thyra::NOTRANS;
2626  RCP<const Thyra::LinearOpBase<ST> > wrapped_op;
2627  Thyra::unwrap(op, &scalar, &transp, &wrapped_op);
2628  pblo = rcp_dynamic_cast<const Thyra::PhysicallyBlockedLinearOpBase<double> >(wrapped_op);
2629  if (!pblo.is_null())
2630  return true;
2631 
2632  return false;
2633 }
2634 
2635 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > getPhysicallyBlockedLinearOp(const LinearOp & op, ST *scalar, bool *transp)
2636 {
2637  // If the operator is a TpetraLinearOp
2638  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > pblo = rcp_dynamic_cast<const Thyra::PhysicallyBlockedLinearOpBase<double> >(op);
2639  if(!pblo.is_null()){
2640  *scalar = 1.0;
2641  *transp = false;
2642  return pblo;
2643  }
2644 
2645  // If the operator is a wrapped TpetraLinearOp
2646  RCP<const Thyra::LinearOpBase<ST> > wrapped_op;
2647  Thyra::EOpTransp eTransp = Thyra::NOTRANS;
2648  Thyra::unwrap(op, scalar, &eTransp, &wrapped_op);
2649  pblo = rcp_dynamic_cast<const Thyra::PhysicallyBlockedLinearOpBase<double> >(wrapped_op,true);
2650  if(!pblo.is_null()){
2651  *transp = true;
2652  if(eTransp == Thyra::NOTRANS)
2653  *transp = false;
2654  return pblo;
2655  }
2656 
2657  return Teuchos::null;
2658 }
2659 
2660 
2661 }
Teko::AbsRowSum
Specifies that the diagonal entry is .
Definition: Teko_Utilities.hpp:814
Teko::BlkDiag
Specifies that a block diagonal approximation is to be used.
Definition: Teko_Utilities.hpp:815
Teko::Diagonal
Specifies that just the diagonal is used.
Definition: Teko_Utilities.hpp:812
Teko::blockRowCount
int blockRowCount(const BlockedLinearOp &blo)
Get the row count in a block linear operator.
Definition: Teko_Utilities.hpp:368
Teko::scale
void scale(const double alpha, MultiVector &x)
Scale a multivector by a constant.
Definition: Teko_Utilities.hpp:593
Teko::createBlockedOp
BlockedLinearOp createBlockedOp()
Build a new blocked linear operator.
Definition: Teko_Utilities.hpp:384
Teko::blockColCount
int blockColCount(const BlockedLinearOp &blo)
Get the column count in a block linear operator.
Definition: Teko_Utilities.hpp:372
Teko::Lumped
Specifies that row sum is used to form a diagonal.
Definition: Teko_Utilities.hpp:813
Teko_Utilities.hpp
Teko::Epetra::EpetraOperatorWrapper
Implements the Epetra_Operator interface with a Thyra LinearOperator. This enables the use of absrtac...
Definition: Teko_EpetraOperatorWrapper.hpp:202
Teko::DiagonalType
DiagonalType
Type describing the type of diagonal to construct.
Definition: Teko_Utilities.hpp:812
Teko::NotDiag
For user convenience, if Teko recieves this value, exceptions will be thrown
Definition: Teko_Utilities.hpp:816