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(opm, &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  auto tExplicitOp = rcp_dynamic_cast<Thyra::TpetraLinearOp<ST,LO,GO,NT> >(destOp);
1307  if(tExplicitOp.is_null())
1308  tExplicitOp = rcp(new Thyra::TpetraLinearOp<ST,LO,GO,NT>());
1309 
1310  // Do explicit matrix-matrix multiply
1311  RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOplm = Tpetra::createCrsMatrix<ST,LO,GO,NT>(tCrsOpl->getRowMap());
1312  RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > explicitCrsOp = Tpetra::createCrsMatrix<ST,LO,GO,NT>(tCrsOpl->getRowMap());
1313  Tpetra::MatrixMatrix::Multiply<ST,LO,GO,NT>(*tCrsOpl,transpl,*tCrsOpm,transpm,*tCrsOplm);
1314  Tpetra::MatrixMatrix::Multiply<ST,LO,GO,NT>(*tCrsOplm,false,*tCrsOpr,transpr,*explicitCrsOp);
1315  explicitCrsOp->resumeFill();
1316  explicitCrsOp->scale(scalarl*scalarm*scalarr);
1317  explicitCrsOp->fillComplete(tCrsOpr->getDomainMap(),tCrsOpl->getRangeMap());
1318  tExplicitOp->initialize(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1319  return tExplicitOp;
1320 
1321  } else if (isTpetral && !isTpetram && isTpetrar){ // Assume that the middle operator is diagonal
1322 
1323  // Get left and right Tpetra crs operators
1324  ST scalarl = 0.0;
1325  bool transpl = false;
1326  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpl = Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1327  ST scalarr = 0.0;
1328  bool transpr = false;
1329  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpr = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1330 
1331  // Cast middle operator as DiagonalLinearOp and extract diagonal as Vector
1332  RCP<const Thyra::DiagonalLinearOpBase<ST> > dOpm = rcp_dynamic_cast<const Thyra::DiagonalLinearOpBase<ST> >(opm,true);
1333  RCP<const Thyra::TpetraVector<ST,LO,GO,NT> > tPtr = rcp_dynamic_cast<const Thyra::TpetraVector<ST,LO,GO,NT> >(dOpm->getDiag(),true);
1334  RCP<const Tpetra::Vector<ST,LO,GO,NT> > diagPtr = rcp_dynamic_cast<const Tpetra::Vector<ST,LO,GO,NT> >(tPtr->getConstTpetraVector(),true);
1335  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()));
1336 
1337  // Do the diagonal scaling
1338  tCrsOplm->rightScale(*diagPtr);
1339 
1340  // Build output operator
1341  RCP<Thyra::LinearOpBase<ST> > explicitOp = rcp(new Thyra::TpetraLinearOp<ST,LO,GO,NT>());
1342  RCP<Thyra::TpetraLinearOp<ST,LO,GO,NT> > tExplicitOp = rcp_dynamic_cast<Thyra::TpetraLinearOp<ST,LO,GO,NT> >(explicitOp);
1343 
1344  // Do explicit matrix-matrix multiply
1345  RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > explicitCrsOp = Tpetra::createCrsMatrix<ST,LO,GO,NT>(tCrsOpl->getRowMap());
1346  Tpetra::MatrixMatrix::Multiply<ST,LO,GO,NT>(*tCrsOplm,false,*tCrsOpr,transpr,*explicitCrsOp);
1347  explicitCrsOp->resumeFill();
1348  explicitCrsOp->scale(scalarl*scalarr);
1349  explicitCrsOp->fillComplete(tCrsOpr->getDomainMap(),tCrsOpl->getRangeMap());
1350  tExplicitOp->initialize(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1351  return tExplicitOp;
1352 
1353  } else { // Assume Epetra and we can use transformers
1354 
1355  // build implicit multiply
1356  const LinearOp implicitOp = Thyra::multiply(opl,opm,opr);
1357 
1358  // build transformer
1359  const RCP<Thyra::LinearOpTransformerBase<double> > prodTrans =
1360  Thyra::epetraExtDiagScaledMatProdTransformer();
1361 
1362  // build operator destination operator
1363  ModifiableLinearOp explicitOp;
1364 
1365  // if neccessary build a operator to put the explicit multiply into
1366  if(destOp==Teuchos::null)
1367  explicitOp = prodTrans->createOutputOp();
1368  else
1369  explicitOp = destOp;
1370 
1371  // perform multiplication
1372  prodTrans->transform(*implicitOp,explicitOp.ptr());
1373 
1374  // label it
1375  explicitOp->setObjectLabel("explicit( " + opl->getObjectLabel() +
1376  " * " + opm->getObjectLabel() +
1377  " * " + opr->getObjectLabel() + " )");
1378 
1379  return explicitOp;
1380 
1381  }
1382 }
1383 
1394 const LinearOp explicitMultiply(const LinearOp & opl,const LinearOp & opr)
1395 {
1396  // if this is a blocked operator, multiply block by block
1397  // it is possible that not every factor in the product is blocked and these situations are handled separately
1398 
1399  bool isBlockedL = isPhysicallyBlockedLinearOp(opl);
1400  bool isBlockedR = isPhysicallyBlockedLinearOp(opr);
1401 
1402  // both factors blocked
1403  if((isBlockedL && isBlockedR)){
1404  double scalarl = 0.0;
1405  bool transpl = false;
1406  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opl = getPhysicallyBlockedLinearOp(opl,&scalarl,&transpl);
1407  double scalarr = 0.0;
1408  bool transpr = false;
1409  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opr = getPhysicallyBlockedLinearOp(opr,&scalarr,&transpr);
1410  double scalar = scalarl*scalarr;
1411 
1412  int numRows = blocked_opl->productRange()->numBlocks();
1413  int numCols = blocked_opr->productDomain()->numBlocks();
1414  int numMiddle = blocked_opl->productDomain()->numBlocks();
1415 
1416  TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == numMiddle);
1417 
1418  RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_product = Thyra::defaultBlockedLinearOp<double>();
1419  blocked_product->beginBlockFill(numRows,numCols);
1420  for(int r = 0; r < numRows; ++r){
1421  for(int c = 0; c < numCols; ++c){
1422  LinearOp product_rc = explicitMultiply(blocked_opl->getBlock(r,0),blocked_opr->getBlock(0,c));
1423  for(int m = 1; m < numMiddle; ++m){
1424  LinearOp product_m = explicitMultiply(blocked_opl->getBlock(r,m),blocked_opr->getBlock(m,c));
1425  product_rc = explicitAdd(product_rc,product_m);
1426  }
1427  blocked_product->setBlock(r,c,Thyra::scale(scalar,product_rc));
1428  }
1429  }
1430  blocked_product->endBlockFill();
1431  return blocked_product;
1432  }
1433 
1434  // only left factor blocked
1435  if((isBlockedL && !isBlockedR)){
1436  double scalarl = 0.0;
1437  bool transpl = false;
1438  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opl = getPhysicallyBlockedLinearOp(opl,&scalarl,&transpl);
1439  double scalar = scalarl;
1440 
1441  int numRows = blocked_opl->productRange()->numBlocks();
1442  int numCols = 1;
1443  int numMiddle = 1;
1444 
1445  TEUCHOS_ASSERT(blocked_opl->productDomain()->numBlocks() == numMiddle);
1446 
1447  RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_product = Thyra::defaultBlockedLinearOp<double>();
1448  blocked_product->beginBlockFill(numRows,numCols);
1449  for(int r = 0; r < numRows; ++r){
1450  LinearOp product_r = explicitMultiply(blocked_opl->getBlock(r,0),opr);
1451  blocked_product->setBlock(r,0,Thyra::scale(scalar,product_r));
1452  }
1453  blocked_product->endBlockFill();
1454  return blocked_product;
1455  }
1456 
1457  // only right factor blocked
1458  if((!isBlockedL && isBlockedR)){
1459  double scalarr = 0.0;
1460  bool transpr = false;
1461  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opr = getPhysicallyBlockedLinearOp(opr,&scalarr,&transpr);
1462  double scalar = scalarr;
1463 
1464  int numRows = 1;
1465  int numCols = blocked_opr->productDomain()->numBlocks();
1466  int numMiddle = 1;
1467 
1468  TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == numMiddle);
1469 
1470  RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_product = Thyra::defaultBlockedLinearOp<double>();
1471  blocked_product->beginBlockFill(numRows,numCols);
1472  for(int c = 0; c < numCols; ++c){
1473  LinearOp product_c = explicitMultiply(opl,blocked_opr->getBlock(0,c));
1474  blocked_product->setBlock(0,c,Thyra::scale(scalar,product_c));
1475  }
1476  blocked_product->endBlockFill();
1477  return blocked_product;
1478  }
1479 
1480  bool isTpetral = Teko::TpetraHelpers::isTpetraLinearOp(opl);
1481  bool isTpetrar = Teko::TpetraHelpers::isTpetraLinearOp(opr);
1482 
1483  if(isTpetral && isTpetrar){ // Both operators are Tpetra matrices so explicitly multiply them
1484  // Get left and right Tpetra crs operators
1485  ST scalarl = 0.0;
1486  bool transpl = false;
1487  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpl = Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1488  ST scalarr = 0.0;
1489  bool transpr = false;
1490  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpr = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1491 
1492  // Build output operator
1493  RCP<Thyra::LinearOpBase<ST> > explicitOp = rcp(new Thyra::TpetraLinearOp<ST,LO,GO,NT>());
1494  RCP<Thyra::TpetraLinearOp<ST,LO,GO,NT> > tExplicitOp = rcp_dynamic_cast<Thyra::TpetraLinearOp<ST,LO,GO,NT> >(explicitOp);
1495 
1496  // Do explicit matrix-matrix multiply
1497  RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > explicitCrsOp = Tpetra::createCrsMatrix<ST,LO,GO,NT>(tCrsOpl->getRowMap());
1498  Tpetra::MatrixMatrix::Multiply<ST,LO,GO,NT>(*tCrsOpl,transpl,*tCrsOpr,transpr,*explicitCrsOp);
1499  explicitCrsOp->resumeFill();
1500  explicitCrsOp->scale(scalarl*scalarr);
1501  explicitCrsOp->fillComplete(tCrsOpr->getDomainMap(),tCrsOpl->getRangeMap());
1502  tExplicitOp->initialize(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1503  return tExplicitOp;
1504 
1505  } else if (isTpetral && !isTpetrar){ // Assume that the right operator is diagonal
1506 
1507  // Get left Tpetra crs operator
1508  ST scalarl = 0.0;
1509  bool transpl = false;
1510  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpl = Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1511 
1512  // Cast right operator as DiagonalLinearOp and extract diagonal as Vector
1513  RCP<const Thyra::DiagonalLinearOpBase<ST> > dOpr = rcp_dynamic_cast<const Thyra::DiagonalLinearOpBase<ST> >(opr,true);
1514  RCP<const Thyra::TpetraVector<ST,LO,GO,NT> > tPtr = rcp_dynamic_cast<const Thyra::TpetraVector<ST,LO,GO,NT> >(dOpr->getDiag(),true);
1515  RCP<const Tpetra::Vector<ST,LO,GO,NT> > diagPtr = rcp_dynamic_cast<const Tpetra::Vector<ST,LO,GO,NT> >(tPtr->getConstTpetraVector(),true);
1516  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()));
1517 
1518  explicitCrsOp->rightScale(*diagPtr);
1519  explicitCrsOp->resumeFill();
1520  explicitCrsOp->scale(scalarl);
1521  explicitCrsOp->fillComplete(tCrsOpl->getDomainMap(),tCrsOpl->getRangeMap());
1522 
1523  return Thyra::constTpetraLinearOp<ST,LO,GO,NT>(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1524 
1525  } else if (!isTpetral && isTpetrar){ // Assume that the left operator is diagonal
1526 
1527  // Get right Tpetra crs operator
1528  ST scalarr = 0.0;
1529  bool transpr = false;
1530  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpr = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1531 
1532  RCP<const Tpetra::Vector<ST,LO,GO,NT> > diagPtr;
1533 
1534  // Cast left operator as DiagonalLinearOp and extract diagonal as Vector
1535  RCP<const Thyra::DiagonalLinearOpBase<ST> > dOpl = rcp_dynamic_cast<const Thyra::DiagonalLinearOpBase<ST> >(opl);
1536  if(dOpl != Teuchos::null){
1537  RCP<const Thyra::TpetraVector<ST,LO,GO,NT> > tPtr = rcp_dynamic_cast<const Thyra::TpetraVector<ST,LO,GO,NT> >(dOpl->getDiag(),true);
1538  diagPtr = rcp_dynamic_cast<const Tpetra::Vector<ST,LO,GO,NT> >(tPtr->getConstTpetraVector(),true);
1539  }
1540  // If it's not diagonal, maybe it's zero
1541  else if(rcp_dynamic_cast<const Thyra::ZeroLinearOpBase<ST> >(opl) != Teuchos::null){
1542  diagPtr = rcp(new Tpetra::Vector<ST,LO,GO,NT>(tCrsOpr->getRangeMap()));
1543  }
1544  else
1545  TEUCHOS_ASSERT(false);
1546 
1547  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()));
1548 
1549  explicitCrsOp->leftScale(*diagPtr);
1550  explicitCrsOp->resumeFill();
1551  explicitCrsOp->scale(scalarr);
1552  explicitCrsOp->fillComplete(tCrsOpr->getDomainMap(),tCrsOpr->getRangeMap());
1553 
1554  return Thyra::constTpetraLinearOp<ST,LO,GO,NT>(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1555 
1556  } else { // Assume Epetra and we can use transformers
1557 
1558  // build implicit multiply
1559  const LinearOp implicitOp = Thyra::multiply(opl,opr);
1560 
1561  // build a scaling transformer
1562  RCP<Thyra::LinearOpTransformerBase<double> > prodTrans
1563  = Thyra::epetraExtDiagScalingTransformer();
1564 
1565  // check to see if a scaling transformer works: if not use the
1566  // DiagScaledMatrixProduct transformer
1567  if(not prodTrans->isCompatible(*implicitOp))
1568  prodTrans = Thyra::epetraExtDiagScaledMatProdTransformer();
1569 
1570  // build operator and multiply
1571  const RCP<Thyra::LinearOpBase<double> > explicitOp = prodTrans->createOutputOp();
1572  prodTrans->transform(*implicitOp,explicitOp.ptr());
1573  explicitOp->setObjectLabel("explicit( " + opl->getObjectLabel() +
1574  " * " + opr->getObjectLabel() + " )");
1575 
1576  return explicitOp;
1577  }
1578 }
1579 
1593 const ModifiableLinearOp explicitMultiply(const LinearOp & opl,const LinearOp & opr,
1594  const ModifiableLinearOp & destOp)
1595 {
1596  // if this is a blocked operator, multiply block by block
1597  // it is possible that not every factor in the product is blocked and these situations are handled separately
1598 
1599  bool isBlockedL = isPhysicallyBlockedLinearOp(opl);
1600  bool isBlockedR = isPhysicallyBlockedLinearOp(opr);
1601 
1602  // both factors blocked
1603  if((isBlockedL && isBlockedR)){
1604  double scalarl = 0.0;
1605  bool transpl = false;
1606  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opl = getPhysicallyBlockedLinearOp(opl,&scalarl,&transpl);
1607  double scalarr = 0.0;
1608  bool transpr = false;
1609  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opr = getPhysicallyBlockedLinearOp(opr,&scalarr,&transpr);
1610  double scalar = scalarl*scalarr;
1611 
1612  int numRows = blocked_opl->productRange()->numBlocks();
1613  int numCols = blocked_opr->productDomain()->numBlocks();
1614  int numMiddle = blocked_opl->productDomain()->numBlocks();
1615 
1616  TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == numMiddle);
1617 
1618  RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_product = Thyra::defaultBlockedLinearOp<double>();
1619  blocked_product->beginBlockFill(numRows,numCols);
1620  for(int r = 0; r < numRows; ++r){
1621  for(int c = 0; c < numCols; ++c){
1622 
1623  LinearOp product_rc = explicitMultiply(blocked_opl->getBlock(r,0),blocked_opr->getBlock(0,c));
1624  for(int m = 1; m < numMiddle; ++m){
1625  LinearOp product_m = explicitMultiply(blocked_opl->getBlock(r,m),blocked_opr->getBlock(m,c));
1626  product_rc = explicitAdd(product_rc,product_m);
1627  }
1628  blocked_product->setBlock(r,c,Thyra::scale(scalar,product_rc));
1629  }
1630  }
1631  blocked_product->endBlockFill();
1632  return blocked_product;
1633  }
1634 
1635  // only left factor blocked
1636  if((isBlockedL && !isBlockedR)){
1637  double scalarl = 0.0;
1638  bool transpl = false;
1639  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opl = getPhysicallyBlockedLinearOp(opl,&scalarl,&transpl);
1640  double scalar = scalarl;
1641 
1642  int numRows = blocked_opl->productRange()->numBlocks();
1643  int numCols = 1;
1644  int numMiddle = 1;
1645 
1646  TEUCHOS_ASSERT(blocked_opl->productDomain()->numBlocks() == numMiddle);
1647 
1648  RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_product = Thyra::defaultBlockedLinearOp<double>();
1649  blocked_product->beginBlockFill(numRows,numCols);
1650  for(int r = 0; r < numRows; ++r){
1651  LinearOp product_r = explicitMultiply(blocked_opl->getBlock(r,0),opr);
1652  blocked_product->setBlock(r,0,Thyra::scale(scalar,product_r));
1653  }
1654  blocked_product->endBlockFill();
1655  return blocked_product;
1656  }
1657 
1658  // only right factor blocked
1659  if((!isBlockedL && isBlockedR)){
1660  double scalarr = 0.0;
1661  bool transpr = false;
1662  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opr = getPhysicallyBlockedLinearOp(opr,&scalarr,&transpr);
1663  double scalar = scalarr;
1664 
1665  int numRows = 1;
1666  int numCols = blocked_opr->productDomain()->numBlocks();
1667  int numMiddle = 1;
1668 
1669  TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == numMiddle);
1670 
1671  RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_product = Thyra::defaultBlockedLinearOp<double>();
1672  blocked_product->beginBlockFill(numRows,numCols);
1673  for(int c = 0; c < numCols; ++c){
1674  LinearOp product_c = explicitMultiply(opl,blocked_opr->getBlock(0,c));
1675  blocked_product->setBlock(0,c,Thyra::scale(scalar,product_c));
1676  }
1677  blocked_product->endBlockFill();
1678  return blocked_product;
1679  }
1680 
1681  bool isTpetral = Teko::TpetraHelpers::isTpetraLinearOp(opl);
1682  bool isTpetrar = Teko::TpetraHelpers::isTpetraLinearOp(opr);
1683 
1684  if(isTpetral && isTpetrar){ // Both operators are Tpetra matrices so use the explicit Tpetra matrix-matrix multiply
1685 
1686  // Get left and right Tpetra crs operators
1687  ST scalarl = 0.0;
1688  bool transpl = false;
1689  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpl = Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1690  ST scalarr = 0.0;
1691  bool transpr = false;
1692  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpr = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1693 
1694  // Build output operator
1695  RCP<Thyra::LinearOpBase<ST> > explicitOp;
1696  if(destOp!=Teuchos::null)
1697  explicitOp = destOp;
1698  else
1699  explicitOp = rcp(new Thyra::TpetraLinearOp<ST,LO,GO,NT>());
1700  RCP<Thyra::TpetraLinearOp<ST,LO,GO,NT> > tExplicitOp = rcp_dynamic_cast<Thyra::TpetraLinearOp<ST,LO,GO,NT> >(explicitOp);
1701 
1702  // Do explicit matrix-matrix multiply
1703  RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > explicitCrsOp = Tpetra::createCrsMatrix<ST,LO,GO,NT>(tCrsOpl->getRowMap());
1704  Tpetra::MatrixMatrix::Multiply<ST,LO,GO,NT>(*tCrsOpl,transpl,*tCrsOpr,transpr,*explicitCrsOp);
1705  explicitCrsOp->resumeFill();
1706  explicitCrsOp->scale(scalarl*scalarr);
1707  explicitCrsOp->fillComplete(tCrsOpr->getDomainMap(),tCrsOpl->getRangeMap());
1708  tExplicitOp->initialize(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1709  return tExplicitOp;
1710 
1711  } else if (isTpetral && !isTpetrar){ // Assume that the right operator is diagonal
1712 
1713  // Get left Tpetra crs operator
1714  ST scalarl = 0.0;
1715  bool transpl = false;
1716  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpl = Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1717 
1718  // Cast right operator as DiagonalLinearOp and extract diagonal as Vector
1719  RCP<const Thyra::DiagonalLinearOpBase<ST> > dOpr = rcp_dynamic_cast<const Thyra::DiagonalLinearOpBase<ST> >(opr);
1720  RCP<const Thyra::TpetraVector<ST,LO,GO,NT> > tPtr = rcp_dynamic_cast<const Thyra::TpetraVector<ST,LO,GO,NT> >(dOpr->getDiag(),true);
1721  RCP<const Tpetra::Vector<ST,LO,GO,NT> > diagPtr = rcp_dynamic_cast<const Tpetra::Vector<ST,LO,GO,NT> >(tPtr->getConstTpetraVector(),true);
1722 
1723  // Scale by the diagonal operator
1724  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()));
1725  explicitCrsOp->rightScale(*diagPtr);
1726  explicitCrsOp->resumeFill();
1727  explicitCrsOp->scale(scalarl);
1728  explicitCrsOp->fillComplete(tCrsOpl->getDomainMap(),tCrsOpl->getRangeMap());
1729  return Thyra::tpetraLinearOp<ST,LO,GO,NT>(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1730 
1731  } else if (!isTpetral && isTpetrar){ // Assume that the left operator is diagonal
1732 
1733  // Get right Tpetra crs operator
1734  ST scalarr = 0.0;
1735  bool transpr = false;
1736  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpr = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1737 
1738  // Cast leftt operator as DiagonalLinearOp and extract diagonal as Vector
1739  RCP<const Thyra::DiagonalLinearOpBase<ST> > dOpl = rcp_dynamic_cast<const Thyra::DiagonalLinearOpBase<ST> >(opl,true);
1740  RCP<const Thyra::TpetraVector<ST,LO,GO,NT> > tPtr = rcp_dynamic_cast<const Thyra::TpetraVector<ST,LO,GO,NT> >(dOpl->getDiag(),true);
1741  RCP<const Tpetra::Vector<ST,LO,GO,NT> > diagPtr = rcp_dynamic_cast<const Tpetra::Vector<ST,LO,GO,NT> >(tPtr->getConstTpetraVector(),true);
1742 
1743  // Scale by the diagonal operator
1744  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()));
1745  explicitCrsOp->leftScale(*diagPtr);
1746  explicitCrsOp->resumeFill();
1747  explicitCrsOp->scale(scalarr);
1748  explicitCrsOp->fillComplete(tCrsOpr->getDomainMap(),tCrsOpr->getRangeMap());
1749  return Thyra::tpetraLinearOp<ST,LO,GO,NT>(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1750 
1751  } else { // Assume Epetra and we can use transformers
1752 
1753  // build implicit multiply
1754  const LinearOp implicitOp = Thyra::multiply(opl,opr);
1755 
1756  // build a scaling transformer
1757 
1758  RCP<Thyra::LinearOpTransformerBase<double> > prodTrans
1759  = Thyra::epetraExtDiagScalingTransformer();
1760 
1761  // check to see if a scaling transformer works: if not use the
1762  // DiagScaledMatrixProduct transformer
1763  if(not prodTrans->isCompatible(*implicitOp))
1764  prodTrans = Thyra::epetraExtDiagScaledMatProdTransformer();
1765 
1766  // build operator destination operator
1767  ModifiableLinearOp explicitOp;
1768 
1769  // if neccessary build a operator to put the explicit multiply into
1770  if(destOp==Teuchos::null)
1771  explicitOp = prodTrans->createOutputOp();
1772  else
1773  explicitOp = destOp;
1774 
1775  // perform multiplication
1776  prodTrans->transform(*implicitOp,explicitOp.ptr());
1777 
1778  // label it
1779  explicitOp->setObjectLabel("explicit( " + opl->getObjectLabel() +
1780  " * " + opr->getObjectLabel() + " )");
1781 
1782  return explicitOp;
1783  }
1784 }
1785 
1796 const LinearOp explicitAdd(const LinearOp & opl_in,const LinearOp & opr_in)
1797 {
1798  // if both blocked, add block by block
1799  if(isPhysicallyBlockedLinearOp(opl_in) && isPhysicallyBlockedLinearOp(opr_in)){
1800  double scalarl = 0.0;
1801  bool transpl = false;
1802  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opl = getPhysicallyBlockedLinearOp(opl_in, &scalarl, &transpl);
1803 
1804  double scalarr = 0.0;
1805  bool transpr = false;
1806  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opr = getPhysicallyBlockedLinearOp(opr_in, &scalarr, &transpr);
1807 
1808  int numRows = blocked_opl->productRange()->numBlocks();
1809  int numCols = blocked_opl->productDomain()->numBlocks();
1810  TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == numRows);
1811  TEUCHOS_ASSERT(blocked_opr->productDomain()->numBlocks() == numCols);
1812 
1813  RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_sum = Thyra::defaultBlockedLinearOp<double>();
1814  blocked_sum->beginBlockFill(numRows,numCols);
1815  for(int r = 0; r < numRows; ++r)
1816  for(int c = 0; c < numCols; ++c)
1817  blocked_sum->setBlock(r,c,explicitAdd(Thyra::scale(scalarl,blocked_opl->getBlock(r,c)),Thyra::scale(scalarr,blocked_opr->getBlock(r,c))));
1818  blocked_sum->endBlockFill();
1819  return blocked_sum;
1820  }
1821 
1822  // if only one is blocked, it must be 1x1
1823  LinearOp opl = opl_in;
1824  LinearOp opr = opr_in;
1825  if(isPhysicallyBlockedLinearOp(opl_in)){
1826  double scalarl = 0.0;
1827  bool transpl = false;
1828  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opl = getPhysicallyBlockedLinearOp(opl_in, &scalarl, &transpl);
1829  TEUCHOS_ASSERT(blocked_opl->productRange()->numBlocks() == 1);
1830  TEUCHOS_ASSERT(blocked_opl->productDomain()->numBlocks() == 1);
1831  opl = Thyra::scale(scalarl,blocked_opl->getBlock(0,0));
1832  }
1833  if(isPhysicallyBlockedLinearOp(opr_in)){
1834  double scalarr = 0.0;
1835  bool transpr = false;
1836  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opr = getPhysicallyBlockedLinearOp(opr_in, &scalarr, &transpr);
1837  TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == 1);
1838  TEUCHOS_ASSERT(blocked_opr->productDomain()->numBlocks() == 1);
1839  opr = Thyra::scale(scalarr,blocked_opr->getBlock(0,0));
1840  }
1841 
1842  bool isTpetral = Teko::TpetraHelpers::isTpetraLinearOp(opl);
1843  bool isTpetrar = Teko::TpetraHelpers::isTpetraLinearOp(opr);
1844 
1845  // if one of the operators in the sum is a thyra zero op
1846  if(isZeroOp(opl)){
1847  if(isZeroOp(opr))
1848  return opr; // return a zero op if both are zero
1849  if(isTpetrar){ // if other op is tpetra, replace this with a zero crs matrix
1850  ST scalar = 0.0;
1851  bool transp = false;
1852  auto crs_op = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalar, &transp);
1853  auto zero_crs = Tpetra::createCrsMatrix<ST,LO,GO,NT>(crs_op->getRowMap());
1854  zero_crs->fillComplete();
1855  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);
1856  isTpetral = true;
1857  } else
1858  return opr->clone();
1859  }
1860  if(isZeroOp(opr)){
1861  if(isTpetral){ // if other op is tpetra, replace this with a zero crs matrix
1862  ST scalar = 0.0;
1863  bool transp = false;
1864  auto crs_op = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalar, &transp);
1865  auto zero_crs = Tpetra::createCrsMatrix<ST,LO,GO,NT>(crs_op->getRowMap());
1866  zero_crs->fillComplete();
1867  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);
1868  isTpetrar = true;
1869  } else
1870  return opl->clone();
1871  }
1872 
1873  if(isTpetral && isTpetrar){ // Both operators are Tpetra matrices so use the explicit Tpetra matrix-matrix add
1874 
1875  // Get left and right Tpetra crs operators
1876  ST scalarl = 0.0;
1877  bool transpl = false;
1878  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpl = Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1879  ST scalarr = 0.0;
1880  bool transpr = false;
1881  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpr = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1882 
1883  // Build output operator
1884  RCP<Thyra::LinearOpBase<ST> > explicitOp = rcp(new Thyra::TpetraLinearOp<ST,LO,GO,NT>());
1885  RCP<Thyra::TpetraLinearOp<ST,LO,GO,NT> > tExplicitOp = rcp_dynamic_cast<Thyra::TpetraLinearOp<ST,LO,GO,NT> >(explicitOp);
1886 
1887  // Do explicit matrix-matrix add
1888  RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > explicitCrsOp = Tpetra::MatrixMatrix::add<ST,LO,GO,NT>(scalarl,transpl,*tCrsOpl,scalarr,transpr,*tCrsOpr);
1889  tExplicitOp->initialize(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1890  return tExplicitOp;
1891 
1892  }else{//Assume Epetra
1893  // build implicit add
1894  const LinearOp implicitOp = Thyra::add(opl,opr);
1895 
1896  // build transformer
1897  const RCP<Thyra::LinearOpTransformerBase<double> > prodTrans =
1898  Thyra::epetraExtAddTransformer();
1899 
1900  // build operator and add
1901  const RCP<Thyra::LinearOpBase<double> > explicitOp = prodTrans->createOutputOp();
1902  prodTrans->transform(*implicitOp,explicitOp.ptr());
1903  explicitOp->setObjectLabel("explicit( " + opl->getObjectLabel() +
1904  " + " + opr->getObjectLabel() + " )");
1905 
1906  return explicitOp;
1907  }
1908 }
1909 
1922 const ModifiableLinearOp explicitAdd(const LinearOp & opl,const LinearOp & opr,
1923  const ModifiableLinearOp & destOp)
1924 {
1925  // if blocked, add block by block
1926  if(isPhysicallyBlockedLinearOp(opl)){
1927  TEUCHOS_ASSERT(isPhysicallyBlockedLinearOp(opr));
1928 
1929  double scalarl = 0.0;
1930  bool transpl = false;
1931  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opl = getPhysicallyBlockedLinearOp(opl, &scalarl, &transpl);
1932 
1933  double scalarr = 0.0;
1934  bool transpr = false;
1935  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opr = getPhysicallyBlockedLinearOp(opr, &scalarr, &transpr);
1936 
1937  int numRows = blocked_opl->productRange()->numBlocks();
1938  int numCols = blocked_opl->productDomain()->numBlocks();
1939  TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == numRows);
1940  TEUCHOS_ASSERT(blocked_opr->productDomain()->numBlocks() == numCols);
1941 
1942  RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_sum = Teuchos::rcp_dynamic_cast<Thyra::DefaultBlockedLinearOp<double>>(destOp);
1943  if(blocked_sum.is_null()) {
1944  // take care of the null case, this means we must alllocate memory
1945  blocked_sum = Thyra::defaultBlockedLinearOp<double>();
1946 
1947  blocked_sum->beginBlockFill(numRows,numCols);
1948  for(int r = 0; r < numRows; ++r) {
1949  for(int c = 0; c < numCols; ++c) {
1950  auto block = explicitAdd(Thyra::scale(scalarl,blocked_opl->getBlock(r,c)),Thyra::scale(scalarr,blocked_opr->getBlock(r,c)),Teuchos::null);
1951  blocked_sum->setNonconstBlock(r,c,block);
1952  }
1953  }
1954  blocked_sum->endBlockFill();
1955 
1956  }
1957  else {
1958  // in this case memory can be reused
1959  for(int r = 0; r < numRows; ++r)
1960  for(int c = 0; c < numCols; ++c)
1961  explicitAdd(Thyra::scale(scalarl,blocked_opl->getBlock(r,c)),Thyra::scale(scalarr,blocked_opr->getBlock(r,c)),blocked_sum->getNonconstBlock(r,c));
1962  }
1963 
1964  return blocked_sum;
1965  }
1966 
1967  bool isTpetral = Teko::TpetraHelpers::isTpetraLinearOp(opl);
1968  bool isTpetrar = Teko::TpetraHelpers::isTpetraLinearOp(opr);
1969 
1970  if(isTpetral && isTpetrar){ // Both operators are Tpetra matrices so use the explicit Tpetra matrix-matrix add
1971 
1972  // Get left and right Tpetra crs operators
1973  ST scalarl = 0.0;
1974  bool transpl = false;
1975  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpl = Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1976  ST scalarr = 0.0;
1977  bool transpr = false;
1978  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpr = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1979 
1980  // Build output operator
1981  RCP<Thyra::LinearOpBase<ST> > explicitOp;
1982  if(destOp!=Teuchos::null)
1983  explicitOp = destOp;
1984  else
1985  explicitOp = rcp(new Thyra::TpetraLinearOp<ST,LO,GO,NT>());
1986  RCP<Thyra::TpetraLinearOp<ST,LO,GO,NT> > tExplicitOp = rcp_dynamic_cast<Thyra::TpetraLinearOp<ST,LO,GO,NT> >(explicitOp);
1987 
1988  // Do explicit matrix-matrix add
1989  RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > explicitCrsOp = Tpetra::MatrixMatrix::add<ST,LO,GO,NT>(scalarl,transpl,*tCrsOpl,scalarr,transpr,*tCrsOpr);
1990  tExplicitOp->initialize(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1991  return tExplicitOp;
1992 
1993  }else{ // Assume Epetra
1994 
1995  // build implicit add
1996  const LinearOp implicitOp = Thyra::add(opl,opr);
1997 
1998  // build transformer
1999  const RCP<Thyra::LinearOpTransformerBase<double> > prodTrans =
2000  Thyra::epetraExtAddTransformer();
2001 
2002  // build or reuse destination operator
2003  RCP<Thyra::LinearOpBase<double> > explicitOp;
2004  if(destOp!=Teuchos::null)
2005  explicitOp = destOp;
2006  else
2007  explicitOp = prodTrans->createOutputOp();
2008 
2009  // perform add
2010  prodTrans->transform(*implicitOp,explicitOp.ptr());
2011  explicitOp->setObjectLabel("explicit( " + opl->getObjectLabel() +
2012  " + " + opr->getObjectLabel() + " )");
2013 
2014  return explicitOp;
2015  }
2016 }
2017 
2022 const ModifiableLinearOp explicitSum(const LinearOp & op,
2023  const ModifiableLinearOp & destOp)
2024 {
2025  // convert operators to Epetra_CrsMatrix
2026  const RCP<const Epetra_CrsMatrix> epetraOp =
2027  rcp_dynamic_cast<const Epetra_CrsMatrix>(get_Epetra_Operator(*op), true);
2028 
2029  if(destOp==Teuchos::null) {
2030  Teuchos::RCP<Epetra_Operator> epetraDest = Teuchos::rcp(new Epetra_CrsMatrix(*epetraOp));
2031 
2032  return Thyra::nonconstEpetraLinearOp(epetraDest);
2033  }
2034 
2035  const RCP<Epetra_CrsMatrix> epetraDest =
2036  rcp_dynamic_cast<Epetra_CrsMatrix>(get_Epetra_Operator(*destOp), true);
2037 
2038  EpetraExt::MatrixMatrix::Add(*epetraOp,false,1.0,*epetraDest,1.0);
2039 
2040  return destOp;
2041 }
2042 
2043 const LinearOp explicitTranspose(const LinearOp & op)
2044 {
2045  if(Teko::TpetraHelpers::isTpetraLinearOp(op)){
2046 
2047  RCP<const Thyra::TpetraLinearOp<ST,LO,GO,NT> > tOp = rcp_dynamic_cast<const Thyra::TpetraLinearOp<ST,LO,GO,NT> >(op,true);
2048  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOp = rcp_dynamic_cast<const Tpetra::CrsMatrix<ST,LO,GO,NT> >(tOp->getConstTpetraOperator(),true);
2049 
2050  Tpetra::RowMatrixTransposer<ST,LO,GO,NT> transposer(tCrsOp);
2051  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > transOp = transposer.createTranspose();
2052 
2053  return Thyra::constTpetraLinearOp<ST,LO,GO,NT>(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(transOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(transOp->getDomainMap()),transOp);
2054 
2055  } else {
2056 
2057  RCP<const Epetra_Operator> eOp = Thyra::get_Epetra_Operator(*op);
2058  TEUCHOS_TEST_FOR_EXCEPTION(eOp==Teuchos::null,std::logic_error,
2059  "Teko::explicitTranspose Not an Epetra_Operator");
2060  RCP<const Epetra_RowMatrix> eRowMatrixOp
2061  = Teuchos::rcp_dynamic_cast<const Epetra_RowMatrix>(eOp);
2062  TEUCHOS_TEST_FOR_EXCEPTION(eRowMatrixOp==Teuchos::null,std::logic_error,
2063  "Teko::explicitTranspose Not an Epetra_RowMatrix");
2064 
2065  // we now have a delete transpose operator
2066  EpetraExt::RowMatrix_Transpose tranposeOp;
2067  Epetra_RowMatrix & eMat = tranposeOp(const_cast<Epetra_RowMatrix &>(*eRowMatrixOp));
2068 
2069  // this copy is because of a poor implementation of the EpetraExt::Transform
2070  // implementation
2071  Teuchos::RCP<Epetra_CrsMatrix> crsMat
2072  = Teuchos::rcp(new Epetra_CrsMatrix(dynamic_cast<Epetra_CrsMatrix &>(eMat)));
2073 
2074  return Thyra::epetraLinearOp(crsMat);
2075  }
2076 }
2077 
2078 double frobeniusNorm(const LinearOp & op_in)
2079 {
2080  LinearOp op;
2081  double scalar = 1.0;
2082 
2083  // if blocked, must be 1x1
2084  if(isPhysicallyBlockedLinearOp(op_in)){
2085  bool transp = false;
2086  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_op = getPhysicallyBlockedLinearOp(op_in,&scalar,&transp);
2087  TEUCHOS_ASSERT(blocked_op->productRange()->numBlocks() == 1);
2088  TEUCHOS_ASSERT(blocked_op->productDomain()->numBlocks() == 1);
2089  op = blocked_op->getBlock(0,0);
2090  } else
2091  op = op_in;
2092 
2093  if(Teko::TpetraHelpers::isTpetraLinearOp(op)){
2094  const RCP<const Thyra::TpetraLinearOp<ST,LO,GO,NT> > tOp = rcp_dynamic_cast<const Thyra::TpetraLinearOp<ST,LO,GO,NT> >(op);
2095  const RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > crsOp = rcp_dynamic_cast<const Tpetra::CrsMatrix<ST,LO,GO,NT> >(tOp->getConstTpetraOperator(),true);
2096  return crsOp->getFrobeniusNorm();
2097  } else {
2098  const RCP<const Epetra_Operator> epOp = Thyra::get_Epetra_Operator(*op);
2099  const RCP<const Epetra_CrsMatrix> crsOp = rcp_dynamic_cast<const Epetra_CrsMatrix>(epOp,true);
2100  return crsOp->NormFrobenius();
2101  }
2102 }
2103 
2104 double oneNorm(const LinearOp & op)
2105 {
2106  if(Teko::TpetraHelpers::isTpetraLinearOp(op)){
2107  TEUCHOS_TEST_FOR_EXCEPTION(true,std::logic_error,"One norm not currently implemented for Tpetra matrices");
2108 
2109  } else {
2110  const RCP<const Epetra_Operator> epOp = Thyra::get_Epetra_Operator(*op);
2111  const RCP<const Epetra_CrsMatrix> crsOp = rcp_dynamic_cast<const Epetra_CrsMatrix>(epOp,true);
2112  return crsOp->NormOne();
2113  }
2114 }
2115 
2116 double infNorm(const LinearOp & op)
2117 {
2118  if(Teko::TpetraHelpers::isTpetraLinearOp(op)){
2119  ST scalar = 0.0;
2120  bool transp = false;
2121  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOp = Teko::TpetraHelpers::getTpetraCrsMatrix(op, &scalar, &transp);
2122 
2123  // extract diagonal
2124  const RCP<Tpetra::Vector<ST,LO,GO,NT> > ptrDiag = Tpetra::createVector<ST,LO,GO,NT>(tCrsOp->getRowMap());
2125  Tpetra::Vector<ST,LO,GO,NT> & diag = *ptrDiag;
2126 
2127  // compute absolute value row sum
2128  diag.putScalar(0.0);
2129  for(LO i=0;i<(LO) tCrsOp->getNodeNumRows();i++) {
2130  LO numEntries = tCrsOp->getNumEntriesInLocalRow (i);
2131  std::vector<LO> indices(numEntries);
2132  std::vector<ST> values(numEntries);
2133  Teuchos::ArrayView<const LO> indices_av(indices);
2134  Teuchos::ArrayView<const ST> values_av(values);
2135  tCrsOp->getLocalRowView(i,indices_av,values_av);
2136 
2137  // build abs value row sum
2138  for(LO j=0;j<numEntries;j++)
2139  diag.sumIntoLocalValue(i,std::abs(values_av[j]));
2140  }
2141  return diag.normInf()*scalar;
2142 
2143  } else {
2144  const RCP<const Epetra_Operator> epOp = Thyra::get_Epetra_Operator(*op);
2145  const RCP<const Epetra_CrsMatrix> crsOp = rcp_dynamic_cast<const Epetra_CrsMatrix>(epOp,true);
2146  return crsOp->NormInf();
2147  }
2148 }
2149 
2150 const LinearOp buildDiagonal(const MultiVector & src,const std::string & lbl)
2151 {
2152  RCP<Thyra::VectorBase<double> > dst = Thyra::createMember(src->range());
2153  Thyra::copy(*src->col(0),dst.ptr());
2154 
2155  return Thyra::diagonal<double>(dst,lbl);
2156 }
2157 
2158 const LinearOp buildInvDiagonal(const MultiVector & src,const std::string & lbl)
2159 {
2160  const RCP<const Thyra::VectorBase<double> > srcV = src->col(0);
2161  RCP<Thyra::VectorBase<double> > dst = Thyra::createMember(srcV->range());
2162  Thyra::reciprocal<double>(*srcV,dst.ptr());
2163 
2164  return Thyra::diagonal<double>(dst,lbl);
2165 }
2166 
2168 BlockedMultiVector buildBlockedMultiVector(const std::vector<MultiVector> & mvv)
2169 {
2170  Teuchos::Array<MultiVector> mvA;
2171  Teuchos::Array<VectorSpace> vsA;
2172 
2173  // build arrays of multi vectors and vector spaces
2174  std::vector<MultiVector>::const_iterator itr;
2175  for(itr=mvv.begin();itr!=mvv.end();++itr) {
2176  mvA.push_back(*itr);
2177  vsA.push_back((*itr)->range());
2178  }
2179 
2180  // construct the product vector space
2181  const RCP<const Thyra::DefaultProductVectorSpace<double> > vs
2182  = Thyra::productVectorSpace<double>(vsA);
2183 
2184  return Thyra::defaultProductMultiVector<double>(vs,mvA);
2185 }
2186 
2197 Teuchos::RCP<Thyra::VectorBase<double> > indicatorVector(
2198  const std::vector<int> & indices,
2199  const VectorSpace & vs,
2200  double onValue, double offValue)
2201 
2202 {
2203  using Teuchos::RCP;
2204 
2205  // create a new vector
2206  RCP<Thyra::VectorBase<double> > v = Thyra::createMember<double>(vs);
2207  Thyra::put_scalar<double>(offValue,v.ptr()); // fill it with "off" values
2208 
2209  // set on values
2210  for(std::size_t i=0;i<indices.size();i++)
2211  Thyra::set_ele<double>(indices[i],onValue,v.ptr());
2212 
2213  return v;
2214 }
2215 
2240 double computeSpectralRad(const RCP<const Thyra::LinearOpBase<double> > & A, double tol,
2241  bool isHermitian,int numBlocks,int restart,int verbosity)
2242 {
2243  typedef Thyra::LinearOpBase<double> OP;
2244  typedef Thyra::MultiVectorBase<double> MV;
2245 
2246  int startVectors = 1;
2247 
2248  // construct an initial guess
2249  const RCP<MV> ivec = Thyra::createMember(A->domain());
2250  Thyra::randomize(-1.0,1.0,ivec.ptr());
2251 
2252  RCP<Anasazi::BasicEigenproblem<double,MV,OP> > eigProb
2253  = rcp(new Anasazi::BasicEigenproblem<double,MV,OP>(A,ivec));
2254  eigProb->setNEV(1);
2255  eigProb->setHermitian(isHermitian);
2256 
2257  // set the problem up
2258  if(not eigProb->setProblem()) {
2259  // big time failure!
2260  return Teuchos::ScalarTraits<double>::nan();
2261  }
2262 
2263  // we want largert magnitude eigenvalue
2264  std::string which("LM"); // largest magnitude
2265 
2266  // Create the parameter list for the eigensolver
2267  // verbosity+=Anasazi::TimingDetails;
2268  Teuchos::ParameterList MyPL;
2269  MyPL.set( "Verbosity", verbosity );
2270  MyPL.set( "Which", which );
2271  MyPL.set( "Block Size", startVectors );
2272  MyPL.set( "Num Blocks", numBlocks );
2273  MyPL.set( "Maximum Restarts", restart );
2274  MyPL.set( "Convergence Tolerance", tol );
2275 
2276  // build status test manager
2277  // RCP<Anasazi::StatusTestMaxIters<double,MV,OP> > statTest
2278  // = rcp(new Anasazi::StatusTestMaxIters<double,MV,OP>(numBlocks*(restart+1)));
2279 
2280  // Create the Block Krylov Schur solver
2281  // This takes as inputs the eigenvalue problem and the solver parameters
2282  Anasazi::BlockKrylovSchurSolMgr<double,MV,OP> MyBlockKrylovSchur(eigProb, MyPL );
2283 
2284  // Solve the eigenvalue problem, and save the return code
2285  Anasazi::ReturnType solverreturn = MyBlockKrylovSchur.solve();
2286 
2287  if(solverreturn==Anasazi::Unconverged) {
2288  double real = MyBlockKrylovSchur.getRitzValues().begin()->realpart;
2289  double comp = MyBlockKrylovSchur.getRitzValues().begin()->imagpart;
2290 
2291  return -std::abs(std::sqrt(real*real+comp*comp));
2292 
2293  // cout << "Anasazi::BlockKrylovSchur::solve() did not converge!" << std::endl;
2294  // return -std::abs(MyBlockKrylovSchur.getRitzValues().begin()->realpart);
2295  }
2296  else { // solverreturn==Anasazi::Converged
2297  double real = eigProb->getSolution().Evals.begin()->realpart;
2298  double comp = eigProb->getSolution().Evals.begin()->imagpart;
2299 
2300  return std::abs(std::sqrt(real*real+comp*comp));
2301 
2302  // cout << "Anasazi::BlockKrylovSchur::solve() converged!" << endl;
2303  // return std::abs(eigProb->getSolution().Evals.begin()->realpart);
2304  }
2305 }
2306 
2330 double computeSmallestMagEig(const RCP<const Thyra::LinearOpBase<double> > & A, double tol,
2331  bool isHermitian,int numBlocks,int restart,int verbosity)
2332 {
2333  typedef Thyra::LinearOpBase<double> OP;
2334  typedef Thyra::MultiVectorBase<double> MV;
2335 
2336  int startVectors = 1;
2337 
2338  // construct an initial guess
2339  const RCP<MV> ivec = Thyra::createMember(A->domain());
2340  Thyra::randomize(-1.0,1.0,ivec.ptr());
2341 
2342  RCP<Anasazi::BasicEigenproblem<double,MV,OP> > eigProb
2343  = rcp(new Anasazi::BasicEigenproblem<double,MV,OP>(A,ivec));
2344  eigProb->setNEV(1);
2345  eigProb->setHermitian(isHermitian);
2346 
2347  // set the problem up
2348  if(not eigProb->setProblem()) {
2349  // big time failure!
2350  return Teuchos::ScalarTraits<double>::nan();
2351  }
2352 
2353  // we want largert magnitude eigenvalue
2354  std::string which("SM"); // smallest magnitude
2355 
2356  // Create the parameter list for the eigensolver
2357  Teuchos::ParameterList MyPL;
2358  MyPL.set( "Verbosity", verbosity );
2359  MyPL.set( "Which", which );
2360  MyPL.set( "Block Size", startVectors );
2361  MyPL.set( "Num Blocks", numBlocks );
2362  MyPL.set( "Maximum Restarts", restart );
2363  MyPL.set( "Convergence Tolerance", tol );
2364 
2365  // build status test manager
2366  // RCP<Anasazi::StatusTestMaxIters<double,MV,OP> > statTest
2367  // = rcp(new Anasazi::StatusTestMaxIters<double,MV,OP>(10));
2368 
2369  // Create the Block Krylov Schur solver
2370  // This takes as inputs the eigenvalue problem and the solver parameters
2371  Anasazi::BlockKrylovSchurSolMgr<double,MV,OP> MyBlockKrylovSchur(eigProb, MyPL );
2372 
2373  // Solve the eigenvalue problem, and save the return code
2374  Anasazi::ReturnType solverreturn = MyBlockKrylovSchur.solve();
2375 
2376  if(solverreturn==Anasazi::Unconverged) {
2377  // cout << "Anasazi::BlockKrylovSchur::solve() did not converge! " << std::endl;
2378  return -std::abs(MyBlockKrylovSchur.getRitzValues().begin()->realpart);
2379  }
2380  else { // solverreturn==Anasazi::Converged
2381  // cout << "Anasazi::BlockKrylovSchur::solve() converged!" << endl;
2382  return std::abs(eigProb->getSolution().Evals.begin()->realpart);
2383  }
2384 }
2385 
2394 ModifiableLinearOp getDiagonalOp(const Teko::LinearOp & A,const DiagonalType & dt)
2395 {
2396  switch(dt) {
2397  case Diagonal:
2398  return getDiagonalOp(A);
2399  case Lumped:
2400  return getLumpedMatrix(A);
2401  case AbsRowSum:
2402  return getAbsRowSumMatrix(A);
2403  case NotDiag:
2404  default:
2405  TEUCHOS_TEST_FOR_EXCEPT(true);
2406  };
2407 
2408  return Teuchos::null;
2409 }
2410 
2419 ModifiableLinearOp getInvDiagonalOp(const Teko::LinearOp & A,const Teko::DiagonalType & dt)
2420 {
2421  switch(dt) {
2422  case Diagonal:
2423  return getInvDiagonalOp(A);
2424  case Lumped:
2425  return getInvLumpedMatrix(A);
2426  case AbsRowSum:
2427  return getAbsRowSumInvMatrix(A);
2428  case NotDiag:
2429  default:
2430  TEUCHOS_TEST_FOR_EXCEPT(true);
2431  };
2432 
2433  return Teuchos::null;
2434 }
2435 
2442 std::string getDiagonalName(const DiagonalType & dt)
2443 {
2444  switch(dt) {
2445  case Diagonal:
2446  return "Diagonal";
2447  case Lumped:
2448  return "Lumped";
2449  case AbsRowSum:
2450  return "AbsRowSum";
2451  case NotDiag:
2452  return "NotDiag";
2453  case BlkDiag:
2454  return "BlkDiag";
2455  };
2456 
2457  return "<error>";
2458 }
2459 
2468 DiagonalType getDiagonalType(std::string name)
2469 {
2470  if(name=="Diagonal")
2471  return Diagonal;
2472  if(name=="Lumped")
2473  return Lumped;
2474  if(name=="AbsRowSum")
2475  return AbsRowSum;
2476  if(name=="BlkDiag")
2477  return BlkDiag;
2478 
2479  return NotDiag;
2480 }
2481 
2482 LinearOp probe(Teuchos::RCP<const Epetra_CrsGraph> &G,const LinearOp & Op){
2483 #ifdef Teko_ENABLE_Isorropia
2484  Teuchos::ParameterList probeList;
2485  Prober prober(G,probeList,true);
2486  Teuchos::RCP<Epetra_CrsMatrix> Mat=rcp(new Epetra_CrsMatrix(Copy,*G));
2488  prober.probe(Mwrap,*Mat);
2489  return Thyra::epetraLinearOp(Mat);
2490 #else
2491  (void)G;
2492  (void)Op;
2493  TEUCHOS_TEST_FOR_EXCEPTION(true,std::runtime_error,"Probe requires Isorropia");
2494 #endif
2495 }
2496 
2497 double norm_1(const MultiVector & v,std::size_t col)
2498 {
2499  Teuchos::Array<double> n(v->domain()->dim());
2500  Thyra::norms_1<double>(*v,n);
2501 
2502  return n[col];
2503 }
2504 
2505 double norm_2(const MultiVector & v,std::size_t col)
2506 {
2507  Teuchos::Array<double> n(v->domain()->dim());
2508  Thyra::norms_2<double>(*v,n);
2509 
2510  return n[col];
2511 }
2512 
2513 void putScalar(const ModifiableLinearOp & op,double scalar)
2514 {
2515  try {
2516  // get Epetra_Operator
2517  RCP<Epetra_Operator> eOp = Thyra::get_Epetra_Operator(*op);
2518 
2519  // cast it to a CrsMatrix
2520  RCP<Epetra_CrsMatrix> eCrsOp = rcp_dynamic_cast<Epetra_CrsMatrix>(eOp,true);
2521 
2522  eCrsOp->PutScalar(scalar);
2523  }
2524  catch (std::exception & e) {
2525  RCP<Teuchos::FancyOStream> out = Teuchos::VerboseObjectBase::getDefaultOStream();
2526 
2527  *out << "Teko: putScalar requires an Epetra_CrsMatrix\n";
2528  *out << " Could not extract an Epetra_Operator from a \"" << op->description() << std::endl;
2529  *out << " OR\n";
2530  *out << " Could not cast an Epetra_Operator to a Epetra_CrsMatrix\n";
2531  *out << std::endl;
2532  *out << "*** THROWN EXCEPTION ***\n";
2533  *out << e.what() << std::endl;
2534  *out << "************************\n";
2535 
2536  throw e;
2537  }
2538 }
2539 
2540 void clipLower(MultiVector & v,double lowerBound)
2541 {
2542  using Teuchos::RCP;
2543  using Teuchos::rcp_dynamic_cast;
2544 
2545  // cast so entries are accessible
2546  // RCP<Thyra::SpmdMultiVectorBase<double> > spmdMVec
2547  // = rcp_dynamic_cast<Thyra::DefaultSpmdMultiVector<double> >(v);
2548 
2549  for(Thyra::Ordinal i=0;i<v->domain()->dim();i++) {
2550  RCP<Thyra::SpmdVectorBase<double> > spmdVec
2551  = rcp_dynamic_cast<Thyra::SpmdVectorBase<double> >(v->col(i),true);
2552 
2553  Teuchos::ArrayRCP<double> values;
2554  // spmdMVec->getNonconstLocalData(Teuchos::ptrFromRef(values),Teuchos::ptrFromRef(i));
2555  spmdVec->getNonconstLocalData(Teuchos::ptrFromRef(values));
2556  for(Teuchos::ArrayRCP<double>::size_type j=0;j<values.size();j++) {
2557  if(values[j]<lowerBound)
2558  values[j] = lowerBound;
2559  }
2560  }
2561 }
2562 
2563 void clipUpper(MultiVector & v,double upperBound)
2564 {
2565  using Teuchos::RCP;
2566  using Teuchos::rcp_dynamic_cast;
2567 
2568  // cast so entries are accessible
2569  // RCP<Thyra::SpmdMultiVectorBase<double> > spmdMVec
2570  // = rcp_dynamic_cast<Thyra::DefaultSpmdMultiVector<double> >(v);
2571  for(Thyra::Ordinal i=0;i<v->domain()->dim();i++) {
2572  RCP<Thyra::SpmdVectorBase<double> > spmdVec
2573  = rcp_dynamic_cast<Thyra::SpmdVectorBase<double> >(v->col(i),true);
2574 
2575  Teuchos::ArrayRCP<double> values;
2576  // spmdMVec->getNonconstLocalData(Teuchos::ptrFromRef(values),Teuchos::ptrFromRef(i));
2577  spmdVec->getNonconstLocalData(Teuchos::ptrFromRef(values));
2578  for(Teuchos::ArrayRCP<double>::size_type j=0;j<values.size();j++) {
2579  if(values[j]>upperBound)
2580  values[j] = upperBound;
2581  }
2582  }
2583 }
2584 
2585 void replaceValue(MultiVector & v,double currentValue,double newValue)
2586 {
2587  using Teuchos::RCP;
2588  using Teuchos::rcp_dynamic_cast;
2589 
2590  // cast so entries are accessible
2591  // RCP<Thyra::SpmdMultiVectorBase<double> > spmdMVec
2592  // = rcp_dynamic_cast<Thyra::SpmdMultiVectorBase<double> >(v,true);
2593  for(Thyra::Ordinal i=0;i<v->domain()->dim();i++) {
2594  RCP<Thyra::SpmdVectorBase<double> > spmdVec
2595  = rcp_dynamic_cast<Thyra::SpmdVectorBase<double> >(v->col(i),true);
2596 
2597  Teuchos::ArrayRCP<double> values;
2598  // spmdMVec->getNonconstLocalData(Teuchos::ptrFromRef(values),Teuchos::ptrFromRef(i));
2599  spmdVec->getNonconstLocalData(Teuchos::ptrFromRef(values));
2600  for(Teuchos::ArrayRCP<double>::size_type j=0;j<values.size();j++) {
2601  if(values[j]==currentValue)
2602  values[j] = newValue;
2603  }
2604  }
2605 }
2606 
2607 void columnAverages(const MultiVector & v,std::vector<double> & averages)
2608 {
2609  averages.resize(v->domain()->dim());
2610 
2611  // sum over each column
2612  Thyra::sums<double>(*v,averages);
2613 
2614  // build averages
2615  Thyra::Ordinal rows = v->range()->dim();
2616  for(std::size_t i=0;i<averages.size();i++)
2617  averages[i] = averages[i]/rows;
2618 }
2619 
2620 double average(const MultiVector & v)
2621 {
2622  Thyra::Ordinal rows = v->range()->dim();
2623  Thyra::Ordinal cols = v->domain()->dim();
2624 
2625  std::vector<double> averages;
2626  columnAverages(v,averages);
2627 
2628  double sum = 0.0;
2629  for(std::size_t i=0;i<averages.size();i++)
2630  sum += averages[i] * rows;
2631 
2632  return sum/(rows*cols);
2633 }
2634 
2635 bool isPhysicallyBlockedLinearOp(const LinearOp & op)
2636 {
2637  // See if the operator is a PBLO
2638  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > pblo = rcp_dynamic_cast<const Thyra::PhysicallyBlockedLinearOpBase<double> >(op);
2639  if (!pblo.is_null())
2640  return true;
2641 
2642  // See if the operator is a wrapped PBLO
2643  ST scalar = 0.0;
2644  Thyra::EOpTransp transp = Thyra::NOTRANS;
2645  RCP<const Thyra::LinearOpBase<ST> > wrapped_op;
2646  Thyra::unwrap(op, &scalar, &transp, &wrapped_op);
2647  pblo = rcp_dynamic_cast<const Thyra::PhysicallyBlockedLinearOpBase<double> >(wrapped_op);
2648  if (!pblo.is_null())
2649  return true;
2650 
2651  return false;
2652 }
2653 
2654 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > getPhysicallyBlockedLinearOp(const LinearOp & op, ST *scalar, bool *transp)
2655 {
2656  // If the operator is a TpetraLinearOp
2657  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > pblo = rcp_dynamic_cast<const Thyra::PhysicallyBlockedLinearOpBase<double> >(op);
2658  if(!pblo.is_null()){
2659  *scalar = 1.0;
2660  *transp = false;
2661  return pblo;
2662  }
2663 
2664  // If the operator is a wrapped TpetraLinearOp
2665  RCP<const Thyra::LinearOpBase<ST> > wrapped_op;
2666  Thyra::EOpTransp eTransp = Thyra::NOTRANS;
2667  Thyra::unwrap(op, scalar, &eTransp, &wrapped_op);
2668  pblo = rcp_dynamic_cast<const Thyra::PhysicallyBlockedLinearOpBase<double> >(wrapped_op,true);
2669  if(!pblo.is_null()){
2670  *transp = true;
2671  if(eTransp == Thyra::NOTRANS)
2672  *transp = false;
2673  return pblo;
2674  }
2675 
2676  return Teuchos::null;
2677 }
2678 
2679 
2680 }
Specifies that a block diagonal approximation is to be used.
Specifies that just the diagonal is used.
int blockRowCount(const BlockedLinearOp &blo)
Get the row count in a block linear operator.
Specifies that row sum is used to form a diagonal.
void scale(const double alpha, MultiVector &x)
Scale a multivector by a constant.
Implements the Epetra_Operator interface with a Thyra LinearOperator. This enables the use of absrtac...
For user convenience, if Teko recieves this value, exceptions will be thrown.
DiagonalType
Type describing the type of diagonal to construct.
int blockColCount(const BlockedLinearOp &blo)
Get the column count in a block linear operator.
BlockedLinearOp createBlockedOp()
Build a new blocked linear operator.
Specifies that the diagonal entry is .