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  const LinearOp test = rcp_dynamic_cast<const Thyra::ZeroLinearOpBase<double> >(op);
613 
614  // if it works...then its zero...otherwise its null
615  return test!=Teuchos::null;
616 }
617 
626 ModifiableLinearOp getAbsRowSumMatrix(const LinearOp & op)
627 {
628  bool isTpetra = false;
629  RCP<const Epetra_CrsMatrix> eCrsOp;
630  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOp;
631 
632  try {
633  // get Epetra or Tpetra Operator
634  RCP<const Thyra::EpetraLinearOp > eOp = rcp_dynamic_cast<const Thyra::EpetraLinearOp >(op);
635  RCP<const Thyra::TpetraLinearOp<ST,LO,GO,NT> > tOp = rcp_dynamic_cast<const Thyra::TpetraLinearOp<ST,LO,GO,NT> >(op);
636 
637  // cast it to a CrsMatrix
638  RCP<Teuchos::FancyOStream> out = Teuchos::VerboseObjectBase::getDefaultOStream();
639  if (!eOp.is_null()){
640  eCrsOp = rcp_dynamic_cast<const Epetra_CrsMatrix>(eOp->epetra_op(),true);
641  }
642  else if (!tOp.is_null()){
643  tCrsOp = rcp_dynamic_cast<const Tpetra::CrsMatrix<ST,LO,GO,NT> >(tOp->getConstTpetraOperator(),true);
644  isTpetra = true;
645  }
646  else
647  throw std::logic_error("Neither Epetra nor Tpetra");
648  }
649  catch (std::exception & e) {
650  RCP<Teuchos::FancyOStream> out = Teuchos::VerboseObjectBase::getDefaultOStream();
651 
652  *out << "Teko: getAbsRowSumMatrix requires an Epetra_CrsMatrix or a Tpetra::CrsMatrix\n";
653  *out << " Could not extract an Epetra_Operator or a Tpetra_Operator from a \"" << op->description() << std::endl;
654  *out << " OR\n";
655  *out << " Could not cast an Epetra_Operator to a Epetra_CrsMatrix or a Tpetra_Operator to a Tpetra::CrsMatrix\n";
656  *out << std::endl;
657  *out << "*** THROWN EXCEPTION ***\n";
658  *out << e.what() << std::endl;
659  *out << "************************\n";
660 
661  throw e;
662  }
663 
664  if(!isTpetra){
665  // extract diagonal
666  const RCP<Epetra_Vector> ptrDiag = rcp(new Epetra_Vector(eCrsOp->RowMap()));
667  Epetra_Vector & diag = *ptrDiag;
668 
669  // compute absolute value row sum
670  diag.PutScalar(0.0);
671  for(int i=0;i<eCrsOp->NumMyRows();i++) {
672  double * values = 0;
673  int numEntries;
674  eCrsOp->ExtractMyRowView(i,numEntries,values);
675 
676  // build abs value row sum
677  for(int j=0;j<numEntries;j++)
678  diag[i] += std::abs(values[j]);
679  }
680 
681  // build Thyra diagonal operator
682  return Teko::Epetra::thyraDiagOp(ptrDiag,eCrsOp->RowMap(),"absRowSum( " + op->getObjectLabel() + " )");
683  }
684  else {
685  // extract diagonal
686  const RCP<Tpetra::Vector<ST,LO,GO,NT> > ptrDiag = Tpetra::createVector<ST,LO,GO,NT>(tCrsOp->getRowMap());
687  Tpetra::Vector<ST,LO,GO,NT> & diag = *ptrDiag;
688 
689  // compute absolute value row sum
690  diag.putScalar(0.0);
691  for(LO i=0;i<(LO) tCrsOp->getNodeNumRows();i++) {
692  LO numEntries = tCrsOp->getNumEntriesInLocalRow (i);
693  std::vector<LO> indices(numEntries);
694  std::vector<ST> values(numEntries);
695  Teuchos::ArrayView<const LO> indices_av(indices);
696  Teuchos::ArrayView<const ST> values_av(values);
697  tCrsOp->getLocalRowView(i,indices_av,values_av);
698 
699  // build abs value row sum
700  for(LO j=0;j<numEntries;j++)
701  diag.sumIntoLocalValue(i,std::abs(values_av[j]));
702  }
703 
704  // build Thyra diagonal operator
705  return Teko::TpetraHelpers::thyraDiagOp(ptrDiag,*tCrsOp->getRowMap(),"absRowSum( " + op->getObjectLabel() + " ))");
706 
707  }
708 
709 }
710 
719 ModifiableLinearOp getAbsRowSumInvMatrix(const LinearOp & op)
720 {
721  // if this is a blocked operator, extract diagonals block by block
722  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_op = rcp_dynamic_cast<const Thyra::PhysicallyBlockedLinearOpBase<double> >(op);
723  if(blocked_op != Teuchos::null){
724  int numRows = blocked_op->productRange()->numBlocks();
725  TEUCHOS_ASSERT(blocked_op->productDomain()->numBlocks() == numRows);
726  RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_diag = Thyra::defaultBlockedLinearOp<double>();
727  blocked_diag->beginBlockFill(numRows,numRows);
728  for(int r = 0; r < numRows; ++r){
729  for(int c = 0; c < numRows; ++c){
730  if(r==c)
731  blocked_diag->setNonconstBlock(r,c,getAbsRowSumInvMatrix(blocked_op->getBlock(r,c)));
732  else
733  blocked_diag->setBlock(r,c,Thyra::zero<double>(blocked_op->getBlock(r,c)->range(),blocked_op->getBlock(r,c)->domain()));
734  }
735  }
736  blocked_diag->endBlockFill();
737  return blocked_diag;
738  }
739 
740  if(Teko::TpetraHelpers::isTpetraLinearOp(op)) {
741  ST scalar = 0.0;
742  bool transp = false;
743  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOp = Teko::TpetraHelpers::getTpetraCrsMatrix(op, &scalar, &transp);
744 
745  // extract diagonal
746  const RCP<Tpetra::Vector<ST,LO,GO,NT> > ptrDiag = Tpetra::createVector<ST,LO,GO,NT>(tCrsOp->getRowMap());
747  Tpetra::Vector<ST,LO,GO,NT> & diag = *ptrDiag;
748 
749  // compute absolute value row sum
750  diag.putScalar(0.0);
751  for(LO i=0;i<(LO) tCrsOp->getNodeNumRows();i++) {
752  LO numEntries = tCrsOp->getNumEntriesInLocalRow (i);
753  std::vector<LO> indices(numEntries);
754  std::vector<ST> values(numEntries);
755  Teuchos::ArrayView<const LO> indices_av(indices);
756  Teuchos::ArrayView<const ST> values_av(values);
757  tCrsOp->getLocalRowView(i,indices_av,values_av);
758 
759  // build abs value row sum
760  for(LO j=0;j<numEntries;j++)
761  diag.sumIntoLocalValue(i,std::abs(values_av[j]));
762  }
763  diag.scale(scalar);
764  diag.reciprocal(diag); // invert entries
765 
766  // build Thyra diagonal operator
767  return Teko::TpetraHelpers::thyraDiagOp(ptrDiag,*tCrsOp->getRowMap(),"absRowSum( " + op->getObjectLabel() + " ))");
768 
769  }
770  else{
771  RCP<const Thyra::EpetraLinearOp > eOp = rcp_dynamic_cast<const Thyra::EpetraLinearOp >(op,true);
772  RCP<const Epetra_CrsMatrix> eCrsOp = rcp_dynamic_cast<const Epetra_CrsMatrix>(eOp->epetra_op(),true);
773 
774  // extract diagonal
775  const RCP<Epetra_Vector> ptrDiag = rcp(new Epetra_Vector(eCrsOp->RowMap()));
776  Epetra_Vector & diag = *ptrDiag;
777 
778  // compute absolute value row sum
779  diag.PutScalar(0.0);
780  for(int i=0;i<eCrsOp->NumMyRows();i++) {
781  double * values = 0;
782  int numEntries;
783  eCrsOp->ExtractMyRowView(i,numEntries,values);
784 
785  // build abs value row sum
786  for(int j=0;j<numEntries;j++)
787  diag[i] += std::abs(values[j]);
788  }
789  diag.Reciprocal(diag); // invert entries
790 
791  // build Thyra diagonal operator
792  return Teko::Epetra::thyraDiagOp(ptrDiag,eCrsOp->RowMap(),"absRowSum( " + op->getObjectLabel() + " )");
793  }
794 
795 }
796 
804 ModifiableLinearOp getLumpedMatrix(const LinearOp & op)
805 {
806  RCP<Thyra::VectorBase<ST> > ones = Thyra::createMember(op->domain());
807  RCP<Thyra::VectorBase<ST> > diag = Thyra::createMember(op->range());
808 
809  // set to all ones
810  Thyra::assign(ones.ptr(),1.0);
811 
812  // compute lumped diagonal
813  // Thyra::apply(*op,Thyra::NONCONJ_ELE,*ones,&*diag);
814  Thyra::apply(*op,Thyra::NOTRANS,*ones,diag.ptr());
815 
816  return rcp(new Thyra::DefaultDiagonalLinearOp<ST>(diag));
817 }
818 
827 ModifiableLinearOp getInvLumpedMatrix(const LinearOp & op)
828 {
829  RCP<Thyra::VectorBase<ST> > ones = Thyra::createMember(op->domain());
830  RCP<Thyra::VectorBase<ST> > diag = Thyra::createMember(op->range());
831 
832  // set to all ones
833  Thyra::assign(ones.ptr(),1.0);
834 
835  // compute lumped diagonal
836  Thyra::apply(*op,Thyra::NOTRANS,*ones,diag.ptr());
837  Thyra::reciprocal(*diag,diag.ptr());
838 
839  return rcp(new Thyra::DefaultDiagonalLinearOp<ST>(diag));
840 }
841 
853 const ModifiableLinearOp getDiagonalOp(const LinearOp & op)
854 {
855  bool isTpetra = false;
856  RCP<const Epetra_CrsMatrix> eCrsOp;
857  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOp;
858 
859  try {
860  // get Epetra or Tpetra Operator
861  RCP<const Thyra::EpetraLinearOp > eOp = rcp_dynamic_cast<const Thyra::EpetraLinearOp >(op);
862  RCP<const Thyra::TpetraLinearOp<ST,LO,GO,NT> > tOp = rcp_dynamic_cast<const Thyra::TpetraLinearOp<ST,LO,GO,NT> >(op);
863 
864  // cast it to a CrsMatrix
865  RCP<Teuchos::FancyOStream> out = Teuchos::VerboseObjectBase::getDefaultOStream();
866  if (!eOp.is_null()){
867  eCrsOp = rcp_dynamic_cast<const Epetra_CrsMatrix>(eOp->epetra_op(),true);
868  }
869  else if (!tOp.is_null()){
870  tCrsOp = rcp_dynamic_cast<const Tpetra::CrsMatrix<ST,LO,GO,NT> >(tOp->getConstTpetraOperator(),true);
871  isTpetra = true;
872  }
873  else
874  throw std::logic_error("Neither Epetra nor Tpetra");
875  }
876  catch (std::exception & e) {
877  RCP<Teuchos::FancyOStream> out = Teuchos::VerboseObjectBase::getDefaultOStream();
878 
879  *out << "Teko: getDiagonalOp requires an Epetra_CrsMatrix or a Tpetra::CrsMatrix\n";
880  *out << " Could not extract an Epetra_Operator or a Tpetra_Operator from a \"" << op->description() << std::endl;
881  *out << " OR\n";
882  *out << " Could not cast an Epetra_Operator to a Epetra_CrsMatrix or a Tpetra_Operator to a Tpetra::CrsMatrix\n";
883  *out << std::endl;
884  *out << "*** THROWN EXCEPTION ***\n";
885  *out << e.what() << std::endl;
886  *out << "************************\n";
887 
888  throw e;
889  }
890 
891  if(!isTpetra){
892  // extract diagonal
893  const RCP<Epetra_Vector> diag = rcp(new Epetra_Vector(eCrsOp->RowMap()));
894  TEUCHOS_TEST_FOR_EXCEPT(eCrsOp->ExtractDiagonalCopy(*diag));
895 
896  // build Thyra diagonal operator
897  return Teko::Epetra::thyraDiagOp(diag,eCrsOp->RowMap(),"inv(diag( " + op->getObjectLabel() + " ))");
898  }
899  else {
900  // extract diagonal
901  const RCP<Tpetra::Vector<ST,LO,GO,NT> > diag = Tpetra::createVector<ST,LO,GO,NT>(tCrsOp->getRowMap());
902  tCrsOp->getLocalDiagCopy(*diag);
903 
904  // build Thyra diagonal operator
905  return Teko::TpetraHelpers::thyraDiagOp(diag,*tCrsOp->getRowMap(),"inv(diag( " + op->getObjectLabel() + " ))");
906 
907  }
908 }
909 
910 const MultiVector getDiagonal(const LinearOp & op)
911 {
912  bool isTpetra = false;
913  RCP<const Epetra_CrsMatrix> eCrsOp;
914  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOp;
915 
916  try {
917  // get Epetra or Tpetra Operator
918  RCP<const Thyra::EpetraLinearOp > eOp = rcp_dynamic_cast<const Thyra::EpetraLinearOp >(op);
919  RCP<const Thyra::TpetraLinearOp<ST,LO,GO,NT> > tOp = rcp_dynamic_cast<const Thyra::TpetraLinearOp<ST,LO,GO,NT> >(op);
920 
921  // cast it to a CrsMatrix
922  RCP<Teuchos::FancyOStream> out = Teuchos::VerboseObjectBase::getDefaultOStream();
923  if (!eOp.is_null()){
924  eCrsOp = rcp_dynamic_cast<const Epetra_CrsMatrix>(eOp->epetra_op(),true);
925  }
926  else if (!tOp.is_null()){
927  tCrsOp = rcp_dynamic_cast<const Tpetra::CrsMatrix<ST,LO,GO,NT> >(tOp->getConstTpetraOperator(),true);
928  isTpetra = true;
929  }
930  else
931  throw std::logic_error("Neither Epetra nor Tpetra");
932  }
933  catch (std::exception & e) {
934  RCP<Teuchos::FancyOStream> out = Teuchos::VerboseObjectBase::getDefaultOStream();
935 
936  RCP<const Thyra::EpetraLinearOp > eOp = rcp_dynamic_cast<const Thyra::EpetraLinearOp >(op);
937  RCP<const Thyra::TpetraLinearOp<ST,LO,GO,NT> > tOp = rcp_dynamic_cast<const Thyra::TpetraLinearOp<ST,LO,GO,NT> >(op);
938  *out << eOp;
939  *out << tOp;
940 
941  *out << "Teko: getDiagonal requires an Epetra_CrsMatrix or a Tpetra::CrsMatrix\n";
942  *out << " Could not extract an Epetra_Operator or a Tpetra_Operator from a \"" << op->description() << std::endl;
943  *out << " OR\n";
944  *out << " Could not cast an Epetra_Operator to a Epetra_CrsMatrix or a Tpetra_Operator to a Tpetra::CrsMatrix\n";
945  *out << std::endl;
946  *out << "*** THROWN EXCEPTION ***\n";
947  *out << e.what() << std::endl;
948  *out << "************************\n";
949 
950  throw e;
951  }
952 
953  if(!isTpetra){
954  // extract diagonal
955  const RCP<Epetra_Vector> diag = rcp(new Epetra_Vector(eCrsOp->RowMap()));
956  TEUCHOS_TEST_FOR_EXCEPT(eCrsOp->ExtractDiagonalCopy(*diag));
957 
958  return Thyra::create_Vector(diag,Thyra::create_VectorSpace(Teuchos::rcpFromRef(eCrsOp->RowMap())));
959  }
960  else {
961  // extract diagonal
962  const RCP<Tpetra::Vector<ST,LO,GO,NT> > diag = Tpetra::createVector<ST,LO,GO,NT>(tCrsOp->getRowMap());
963  tCrsOp->getLocalDiagCopy(*diag);
964 
965  // build Thyra diagonal operator
966  return Thyra::createVector<ST,LO,GO,NT>(diag,Thyra::createVectorSpace<ST,LO,GO,NT>(tCrsOp->getRowMap()));
967 
968  }
969 }
970 
971 const MultiVector getDiagonal(const Teko::LinearOp & A,const DiagonalType & dt)
972 {
973  LinearOp diagOp = Teko::getDiagonalOp(A,dt);
974 
975  Teuchos::RCP<const Thyra::MultiVectorBase<double> > v =
976  Teuchos::rcp_dynamic_cast<const Thyra::DiagonalLinearOpBase<double> >(diagOp)->getDiag();
977  return Teuchos::rcp_const_cast<Thyra::MultiVectorBase<double> >(v);
978 }
979 
991 const ModifiableLinearOp getInvDiagonalOp(const LinearOp & op)
992 {
993  if (Teko::TpetraHelpers::isTpetraLinearOp(op)){
994  ST scalar = 0.0;
995  bool transp = false;
996  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOp = Teko::TpetraHelpers::getTpetraCrsMatrix(op, &scalar, &transp);
997 
998  // extract diagonal
999  const RCP<Tpetra::Vector<ST,LO,GO,NT> > diag = Tpetra::createVector<ST,LO,GO,NT>(tCrsOp->getRowMap());
1000  diag->scale(scalar);
1001  tCrsOp->getLocalDiagCopy(*diag);
1002  diag->reciprocal(*diag);
1003 
1004  // build Thyra diagonal operator
1005  return Teko::TpetraHelpers::thyraDiagOp(diag,*tCrsOp->getRowMap(),"inv(diag( " + op->getObjectLabel() + " ))");
1006 
1007  }
1008  else {
1009  RCP<const Thyra::EpetraLinearOp > eOp = rcp_dynamic_cast<const Thyra::EpetraLinearOp >(op,true);
1010  RCP<const Epetra_CrsMatrix> eCrsOp = rcp_dynamic_cast<const Epetra_CrsMatrix>(eOp->epetra_op(),true);
1011 
1012  // extract diagonal
1013  const RCP<Epetra_Vector> diag = rcp(new Epetra_Vector(eCrsOp->RowMap()));
1014  TEUCHOS_TEST_FOR_EXCEPT(eCrsOp->ExtractDiagonalCopy(*diag));
1015  diag->Reciprocal(*diag);
1016 
1017  // build Thyra diagonal operator
1018  return Teko::Epetra::thyraDiagOp(diag,eCrsOp->RowMap(),"inv(diag( " + op->getObjectLabel() + " ))");
1019  }
1020 }
1021 
1034 const LinearOp explicitMultiply(const LinearOp & opl,const LinearOp & opm,const LinearOp & opr)
1035 {
1036  // if this is a blocked operator, multiply block by block
1037  // it is possible that not every factor in the product is blocked and these situations are handled separately
1038 
1039  bool isBlockedL = isPhysicallyBlockedLinearOp(opl);
1040  bool isBlockedM = isPhysicallyBlockedLinearOp(opm);
1041  bool isBlockedR = isPhysicallyBlockedLinearOp(opr);
1042 
1043  // all factors blocked
1044  if((isBlockedL && isBlockedM && isBlockedR)){
1045 
1046  double scalarl = 0.0;
1047  bool transpl = false;
1048  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opl = getPhysicallyBlockedLinearOp(opl,&scalarl,&transpl);
1049  double scalarm = 0.0;
1050  bool transpm = false;
1051  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opm = getPhysicallyBlockedLinearOp(opm,&scalarm,&transpm);
1052  double scalarr = 0.0;
1053  bool transpr = false;
1054  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opr = getPhysicallyBlockedLinearOp(opr,&scalarr,&transpr);
1055  double scalar = scalarl*scalarm*scalarr;
1056 
1057  int numRows = blocked_opl->productRange()->numBlocks();
1058  int numCols = blocked_opr->productDomain()->numBlocks();
1059  int numMiddle = blocked_opm->productRange()->numBlocks();
1060 
1061  // Assume that the middle block is block nxn and that it's diagonal. Otherwise use the two argument explicitMultiply twice
1062  TEUCHOS_ASSERT(blocked_opm->productDomain()->numBlocks() == numMiddle);
1063  TEUCHOS_ASSERT(blocked_opl->productDomain()->numBlocks() == numMiddle);
1064  TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == numMiddle);
1065 
1066  RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_product = Thyra::defaultBlockedLinearOp<double>();
1067  blocked_product->beginBlockFill(numRows,numCols);
1068  for(int r = 0; r < numRows; ++r){
1069  for(int c = 0; c < numCols; ++c){
1070  LinearOp product_rc = explicitMultiply(blocked_opl->getBlock(r,0),blocked_opm->getBlock(0,0),blocked_opr->getBlock(0,c));
1071  for(int m = 1; m < numMiddle; ++m){
1072  LinearOp product_m = explicitMultiply(blocked_opl->getBlock(r,m),blocked_opm->getBlock(m,m),blocked_opr->getBlock(m,c));
1073  product_rc = explicitAdd(product_rc,product_m);
1074  }
1075  blocked_product->setBlock(r,c,product_rc);
1076  }
1077  }
1078  blocked_product->endBlockFill();
1079  return Thyra::scale<double>(scalar,blocked_product.getConst());
1080  }
1081 
1082  // left and right factors blocked
1083  if(isBlockedL && !isBlockedM && isBlockedR){
1084  double scalarl = 0.0;
1085  bool transpl = false;
1086  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opl = getPhysicallyBlockedLinearOp(opl,&scalarl,&transpl);
1087  double scalarr = 0.0;
1088  bool transpr = false;
1089  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opr = getPhysicallyBlockedLinearOp(opr,&scalarr,&transpr);
1090  double scalar = scalarl*scalarr;
1091 
1092  int numRows = blocked_opl->productRange()->numBlocks();
1093  int numCols = blocked_opr->productDomain()->numBlocks();
1094  int numMiddle = 1;
1095 
1096  // Assume that the middle block is 1x1 diagonal. Left must be rx1, right 1xc
1097  TEUCHOS_ASSERT(blocked_opl->productDomain()->numBlocks() == numMiddle);
1098  TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == numMiddle);
1099 
1100  RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_product = Thyra::defaultBlockedLinearOp<double>();
1101  blocked_product->beginBlockFill(numRows,numCols);
1102  for(int r = 0; r < numRows; ++r){
1103  for(int c = 0; c < numCols; ++c){
1104  LinearOp product_rc = explicitMultiply(blocked_opl->getBlock(r,0),opm,blocked_opr->getBlock(0,c));
1105  blocked_product->setBlock(r,c,product_rc);
1106  }
1107  }
1108  blocked_product->endBlockFill();
1109  return Thyra::scale<double>(scalar,blocked_product.getConst());
1110  }
1111 
1112  // only right factor blocked
1113  if(!isBlockedL && !isBlockedM && isBlockedR){
1114  double scalarr = 0.0;
1115  bool transpr = false;
1116  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opr = getPhysicallyBlockedLinearOp(opr,&scalarr,&transpr);
1117  double scalar = scalarr;
1118 
1119  int numRows = 1;
1120  int numCols = blocked_opr->productDomain()->numBlocks();
1121  int numMiddle = 1;
1122 
1123  // Assume that the middle block is 1x1 diagonal, left is 1x1. Right must be 1xc
1124  TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == numMiddle);
1125 
1126  RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_product = Thyra::defaultBlockedLinearOp<double>();
1127  blocked_product->beginBlockFill(numRows,numCols);
1128  for(int c = 0; c < numCols; ++c){
1129  LinearOp product_c = explicitMultiply(opl,opm,blocked_opr->getBlock(0,c));
1130  blocked_product->setBlock(0,c,product_c);
1131  }
1132  blocked_product->endBlockFill();
1133  return Thyra::scale<double>(scalar,blocked_product.getConst());
1134  }
1135 
1136  //TODO: three more cases (only non-blocked - blocked - non-blocked not possible)
1137 
1138  bool isTpetral = Teko::TpetraHelpers::isTpetraLinearOp(opl);
1139  bool isTpetram = Teko::TpetraHelpers::isTpetraLinearOp(opm);
1140  bool isTpetrar = Teko::TpetraHelpers::isTpetraLinearOp(opr);
1141 
1142  if(isTpetral && isTpetram && isTpetrar){ // Both operators are Tpetra matrices so explicitly multiply them
1143 
1144  // Get left and right Tpetra crs operators
1145  ST scalarl = 0.0;
1146  bool transpl = false;
1147  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpl = Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1148  ST scalarm = 0.0;
1149  bool transpm = false;
1150  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpm = Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarm, &transpm);
1151  ST scalarr = 0.0;
1152  bool transpr = false;
1153  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpr = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1154 
1155  // Build output operator
1156  RCP<Thyra::LinearOpBase<ST> > explicitOp = rcp(new Thyra::TpetraLinearOp<ST,LO,GO,NT>());
1157  RCP<Thyra::TpetraLinearOp<ST,LO,GO,NT> > tExplicitOp = rcp_dynamic_cast<Thyra::TpetraLinearOp<ST,LO,GO,NT> >(explicitOp);
1158 
1159  // Do explicit matrix-matrix multiply
1160  RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOplm = Tpetra::createCrsMatrix<ST,LO,GO,NT>(tCrsOpl->getRowMap());
1161  RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > explicitCrsOp = Tpetra::createCrsMatrix<ST,LO,GO,NT>(tCrsOpl->getRowMap());
1162  Tpetra::MatrixMatrix::Multiply<ST,LO,GO,NT>(*tCrsOpl,transpl,*tCrsOpm,transpm,*tCrsOplm);
1163  Tpetra::MatrixMatrix::Multiply<ST,LO,GO,NT>(*tCrsOplm,false,*tCrsOpr,transpr,*explicitCrsOp);
1164  explicitCrsOp->resumeFill();
1165  explicitCrsOp->scale(scalarl*scalarm*scalarr);
1166  explicitCrsOp->fillComplete(tCrsOpr->getDomainMap(),tCrsOpl->getRangeMap());
1167  tExplicitOp->initialize(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1168  return tExplicitOp;
1169 
1170  } else if (isTpetral && !isTpetram && isTpetrar){ // Assume that the middle operator is diagonal
1171 
1172  // Get left and right Tpetra crs operators
1173  ST scalarl = 0.0;
1174  bool transpl = false;
1175  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpl = Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1176  ST scalarr = 0.0;
1177  bool transpr = false;
1178  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpr = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1179 
1180  RCP<const Tpetra::Vector<ST,LO,GO,NT> > diagPtr;
1181 
1182  // Cast middle operator as DiagonalLinearOp and extract diagonal as Vector
1183  RCP<const Thyra::DiagonalLinearOpBase<ST> > dOpm = rcp_dynamic_cast<const Thyra::DiagonalLinearOpBase<ST> >(opm);
1184  if(dOpm != Teuchos::null){
1185  RCP<const Thyra::TpetraVector<ST,LO,GO,NT> > tPtr = rcp_dynamic_cast<const Thyra::TpetraVector<ST,LO,GO,NT> >(dOpm->getDiag(),true);
1186  diagPtr = rcp_dynamic_cast<const Tpetra::Vector<ST,LO,GO,NT> >(tPtr->getConstTpetraVector(),true);
1187  }
1188  // If it's not diagonal, maybe it's zero
1189  else if(rcp_dynamic_cast<const Thyra::ZeroLinearOpBase<ST> >(opm) != Teuchos::null){
1190  diagPtr = rcp(new Tpetra::Vector<ST,LO,GO,NT>(tCrsOpl->getDomainMap()));
1191  }
1192  else
1193  TEUCHOS_ASSERT(false);
1194 
1195  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()));
1196 
1197  // Do the diagonal scaling
1198  tCrsOplm->rightScale(*diagPtr);
1199 
1200  // Build output operator
1201  RCP<Thyra::LinearOpBase<ST> > explicitOp = rcp(new Thyra::TpetraLinearOp<ST,LO,GO,NT>());
1202  RCP<Thyra::TpetraLinearOp<ST,LO,GO,NT> > tExplicitOp = rcp_dynamic_cast<Thyra::TpetraLinearOp<ST,LO,GO,NT> >(explicitOp);
1203 
1204  // Do explicit matrix-matrix multiply
1205  RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > explicitCrsOp = Tpetra::createCrsMatrix<ST,LO,GO,NT>(tCrsOpl->getRowMap());
1206  Tpetra::MatrixMatrix::Multiply<ST,LO,GO,NT>(*tCrsOplm,false,*tCrsOpr,transpr,*explicitCrsOp);
1207  explicitCrsOp->resumeFill();
1208  explicitCrsOp->scale(scalarl*scalarr);
1209  explicitCrsOp->fillComplete(tCrsOpr->getDomainMap(),tCrsOpl->getRangeMap());
1210  tExplicitOp->initialize(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1211  return tExplicitOp;
1212 
1213  } else { // Assume Epetra and we can use transformers
1214 
1215  // build implicit multiply
1216  const LinearOp implicitOp = Thyra::multiply(opl,opm,opr);
1217 
1218  // build transformer
1219  const RCP<Thyra::LinearOpTransformerBase<double> > prodTrans =
1220  Thyra::epetraExtDiagScaledMatProdTransformer();
1221 
1222  // build operator and multiply
1223  const RCP<Thyra::LinearOpBase<double> > explicitOp = prodTrans->createOutputOp();
1224  prodTrans->transform(*implicitOp,explicitOp.ptr());
1225  explicitOp->setObjectLabel("explicit( " + opl->getObjectLabel() +
1226  " * " + opm->getObjectLabel() +
1227  " * " + opr->getObjectLabel() + " )");
1228 
1229  return explicitOp;
1230 
1231  }
1232 }
1233 
1248 const ModifiableLinearOp explicitMultiply(const LinearOp & opl,const LinearOp & opm,const LinearOp & opr,
1249  const ModifiableLinearOp & destOp)
1250 {
1251  bool isTpetral = Teko::TpetraHelpers::isTpetraLinearOp(opl);
1252  bool isTpetram = Teko::TpetraHelpers::isTpetraLinearOp(opm);
1253  bool isTpetrar = Teko::TpetraHelpers::isTpetraLinearOp(opr);
1254 
1255  if(isTpetral && isTpetram && isTpetrar){ // Both operators are Tpetra matrices so explicitly multiply them
1256 
1257  // Get left and right Tpetra crs operators
1258  ST scalarl = 0.0;
1259  bool transpl = false;
1260  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpl = Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1261  ST scalarm = 0.0;
1262  bool transpm = false;
1263  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpm = Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarm, &transpm);
1264  ST scalarr = 0.0;
1265  bool transpr = false;
1266  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpr = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1267 
1268  // Build output operator
1269  RCP<Thyra::LinearOpBase<ST> > explicitOp = rcp(new Thyra::TpetraLinearOp<ST,LO,GO,NT>());
1270  RCP<Thyra::TpetraLinearOp<ST,LO,GO,NT> > tExplicitOp = rcp_dynamic_cast<Thyra::TpetraLinearOp<ST,LO,GO,NT> >(explicitOp);
1271 
1272  // Do explicit matrix-matrix multiply
1273  RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOplm = Tpetra::createCrsMatrix<ST,LO,GO,NT>(tCrsOpl->getRowMap());
1274  RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > explicitCrsOp = Tpetra::createCrsMatrix<ST,LO,GO,NT>(tCrsOpl->getRowMap());
1275  Tpetra::MatrixMatrix::Multiply<ST,LO,GO,NT>(*tCrsOpl,transpl,*tCrsOpm,transpm,*tCrsOplm);
1276  Tpetra::MatrixMatrix::Multiply<ST,LO,GO,NT>(*tCrsOplm,false,*tCrsOpr,transpr,*explicitCrsOp);
1277  explicitCrsOp->resumeFill();
1278  explicitCrsOp->scale(scalarl*scalarm*scalarr);
1279  explicitCrsOp->fillComplete(tCrsOpr->getDomainMap(),tCrsOpl->getRangeMap());
1280  tExplicitOp->initialize(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1281  return tExplicitOp;
1282 
1283  } else if (isTpetral && !isTpetram && isTpetrar){ // Assume that the middle operator is diagonal
1284 
1285  // Get left and right Tpetra crs operators
1286  ST scalarl = 0.0;
1287  bool transpl = false;
1288  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpl = Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1289  ST scalarr = 0.0;
1290  bool transpr = false;
1291  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpr = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1292 
1293  // Cast middle operator as DiagonalLinearOp and extract diagonal as Vector
1294  RCP<const Thyra::DiagonalLinearOpBase<ST> > dOpm = rcp_dynamic_cast<const Thyra::DiagonalLinearOpBase<ST> >(opm,true);
1295  RCP<const Thyra::TpetraVector<ST,LO,GO,NT> > tPtr = rcp_dynamic_cast<const Thyra::TpetraVector<ST,LO,GO,NT> >(dOpm->getDiag(),true);
1296  RCP<const Tpetra::Vector<ST,LO,GO,NT> > diagPtr = rcp_dynamic_cast<const Tpetra::Vector<ST,LO,GO,NT> >(tPtr->getConstTpetraVector(),true);
1297  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()));
1298 
1299  // Do the diagonal scaling
1300  tCrsOplm->rightScale(*diagPtr);
1301 
1302  // Build output operator
1303  RCP<Thyra::LinearOpBase<ST> > explicitOp = rcp(new Thyra::TpetraLinearOp<ST,LO,GO,NT>());
1304  RCP<Thyra::TpetraLinearOp<ST,LO,GO,NT> > tExplicitOp = rcp_dynamic_cast<Thyra::TpetraLinearOp<ST,LO,GO,NT> >(explicitOp);
1305 
1306  // Do explicit matrix-matrix multiply
1307  RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > explicitCrsOp = Tpetra::createCrsMatrix<ST,LO,GO,NT>(tCrsOpl->getRowMap());
1308  Tpetra::MatrixMatrix::Multiply<ST,LO,GO,NT>(*tCrsOplm,false,*tCrsOpr,transpr,*explicitCrsOp);
1309  explicitCrsOp->resumeFill();
1310  explicitCrsOp->scale(scalarl*scalarr);
1311  explicitCrsOp->fillComplete(tCrsOpr->getDomainMap(),tCrsOpl->getRangeMap());
1312  tExplicitOp->initialize(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1313  return tExplicitOp;
1314 
1315  } else { // Assume Epetra and we can use transformers
1316 
1317  // build implicit multiply
1318  const LinearOp implicitOp = Thyra::multiply(opl,opm,opr);
1319 
1320  // build transformer
1321  const RCP<Thyra::LinearOpTransformerBase<double> > prodTrans =
1322  Thyra::epetraExtDiagScaledMatProdTransformer();
1323 
1324  // build operator destination operator
1325  ModifiableLinearOp explicitOp;
1326 
1327  // if neccessary build a operator to put the explicit multiply into
1328  if(destOp==Teuchos::null)
1329  explicitOp = prodTrans->createOutputOp();
1330  else
1331  explicitOp = destOp;
1332 
1333  // perform multiplication
1334  prodTrans->transform(*implicitOp,explicitOp.ptr());
1335 
1336  // label it
1337  explicitOp->setObjectLabel("explicit( " + opl->getObjectLabel() +
1338  " * " + opm->getObjectLabel() +
1339  " * " + opr->getObjectLabel() + " )");
1340 
1341  return explicitOp;
1342 
1343  }
1344 }
1345 
1356 const LinearOp explicitMultiply(const LinearOp & opl,const LinearOp & opr)
1357 {
1358  // if this is a blocked operator, multiply block by block
1359  // it is possible that not every factor in the product is blocked and these situations are handled separately
1360 
1361  bool isBlockedL = isPhysicallyBlockedLinearOp(opl);
1362  bool isBlockedR = isPhysicallyBlockedLinearOp(opr);
1363 
1364  // both factors blocked
1365  if((isBlockedL && isBlockedR)){
1366  double scalarl = 0.0;
1367  bool transpl = false;
1368  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opl = getPhysicallyBlockedLinearOp(opl,&scalarl,&transpl);
1369  double scalarr = 0.0;
1370  bool transpr = false;
1371  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opr = getPhysicallyBlockedLinearOp(opr,&scalarr,&transpr);
1372  double scalar = scalarl*scalarr;
1373 
1374  int numRows = blocked_opl->productRange()->numBlocks();
1375  int numCols = blocked_opr->productDomain()->numBlocks();
1376  int numMiddle = blocked_opl->productDomain()->numBlocks();
1377 
1378  TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == numMiddle);
1379 
1380  RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_product = Thyra::defaultBlockedLinearOp<double>();
1381  blocked_product->beginBlockFill(numRows,numCols);
1382  for(int r = 0; r < numRows; ++r){
1383  for(int c = 0; c < numCols; ++c){
1384  LinearOp product_rc = explicitMultiply(blocked_opl->getBlock(r,0),blocked_opr->getBlock(0,c));
1385  for(int m = 1; m < numMiddle; ++m){
1386  LinearOp product_m = explicitMultiply(blocked_opl->getBlock(r,m),blocked_opr->getBlock(m,c));
1387  product_rc = explicitAdd(product_rc,product_m);
1388  }
1389  blocked_product->setBlock(r,c,Thyra::scale(scalar,product_rc));
1390  }
1391  }
1392  blocked_product->endBlockFill();
1393  return blocked_product;
1394  }
1395 
1396  // only left factor blocked
1397  if((isBlockedL && !isBlockedR)){
1398  double scalarl = 0.0;
1399  bool transpl = false;
1400  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opl = getPhysicallyBlockedLinearOp(opl,&scalarl,&transpl);
1401  double scalar = scalarl;
1402 
1403  int numRows = blocked_opl->productRange()->numBlocks();
1404  int numCols = 1;
1405  int numMiddle = 1;
1406 
1407  TEUCHOS_ASSERT(blocked_opl->productDomain()->numBlocks() == numMiddle);
1408 
1409  RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_product = Thyra::defaultBlockedLinearOp<double>();
1410  blocked_product->beginBlockFill(numRows,numCols);
1411  for(int r = 0; r < numRows; ++r){
1412  LinearOp product_r = explicitMultiply(blocked_opl->getBlock(r,0),opr);
1413  blocked_product->setBlock(r,0,Thyra::scale(scalar,product_r));
1414  }
1415  blocked_product->endBlockFill();
1416  return blocked_product;
1417  }
1418 
1419  // only right factor blocked
1420  if((!isBlockedL && isBlockedR)){
1421  double scalarr = 0.0;
1422  bool transpr = false;
1423  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opr = getPhysicallyBlockedLinearOp(opr,&scalarr,&transpr);
1424  double scalar = scalarr;
1425 
1426  int numRows = 1;
1427  int numCols = blocked_opr->productDomain()->numBlocks();
1428  int numMiddle = 1;
1429 
1430  TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == numMiddle);
1431 
1432  RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_product = Thyra::defaultBlockedLinearOp<double>();
1433  blocked_product->beginBlockFill(numRows,numCols);
1434  for(int c = 0; c < numCols; ++c){
1435  LinearOp product_c = explicitMultiply(opl,blocked_opr->getBlock(0,c));
1436  blocked_product->setBlock(0,c,Thyra::scale(scalar,product_c));
1437  }
1438  blocked_product->endBlockFill();
1439  return blocked_product;
1440  }
1441 
1442  bool isTpetral = Teko::TpetraHelpers::isTpetraLinearOp(opl);
1443  bool isTpetrar = Teko::TpetraHelpers::isTpetraLinearOp(opr);
1444 
1445  if(isTpetral && isTpetrar){ // Both operators are Tpetra matrices so explicitly multiply them
1446  // Get left and right Tpetra crs operators
1447  ST scalarl = 0.0;
1448  bool transpl = false;
1449  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpl = Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1450  ST scalarr = 0.0;
1451  bool transpr = false;
1452  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpr = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1453 
1454  // Build output operator
1455  RCP<Thyra::LinearOpBase<ST> > explicitOp = rcp(new Thyra::TpetraLinearOp<ST,LO,GO,NT>());
1456  RCP<Thyra::TpetraLinearOp<ST,LO,GO,NT> > tExplicitOp = rcp_dynamic_cast<Thyra::TpetraLinearOp<ST,LO,GO,NT> >(explicitOp);
1457 
1458  // Do explicit matrix-matrix multiply
1459  RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > explicitCrsOp = Tpetra::createCrsMatrix<ST,LO,GO,NT>(tCrsOpl->getRowMap());
1460  Tpetra::MatrixMatrix::Multiply<ST,LO,GO,NT>(*tCrsOpl,transpl,*tCrsOpr,transpr,*explicitCrsOp);
1461  explicitCrsOp->resumeFill();
1462  explicitCrsOp->scale(scalarl*scalarr);
1463  explicitCrsOp->fillComplete(tCrsOpr->getDomainMap(),tCrsOpl->getRangeMap());
1464  tExplicitOp->initialize(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1465  return tExplicitOp;
1466 
1467  } else if (isTpetral && !isTpetrar){ // Assume that the right operator is diagonal
1468 
1469  // Get left Tpetra crs operator
1470  ST scalarl = 0.0;
1471  bool transpl = false;
1472  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpl = Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1473 
1474  // Cast right operator as DiagonalLinearOp and extract diagonal as Vector
1475  RCP<const Thyra::DiagonalLinearOpBase<ST> > dOpr = rcp_dynamic_cast<const Thyra::DiagonalLinearOpBase<ST> >(opr,true);
1476  RCP<const Thyra::TpetraVector<ST,LO,GO,NT> > tPtr = rcp_dynamic_cast<const Thyra::TpetraVector<ST,LO,GO,NT> >(dOpr->getDiag(),true);
1477  RCP<const Tpetra::Vector<ST,LO,GO,NT> > diagPtr = rcp_dynamic_cast<const Tpetra::Vector<ST,LO,GO,NT> >(tPtr->getConstTpetraVector(),true);
1478  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()));
1479 
1480  explicitCrsOp->rightScale(*diagPtr);
1481  explicitCrsOp->resumeFill();
1482  explicitCrsOp->scale(scalarl);
1483  explicitCrsOp->fillComplete(tCrsOpl->getDomainMap(),tCrsOpl->getRangeMap());
1484 
1485  return Thyra::constTpetraLinearOp<ST,LO,GO,NT>(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1486 
1487  } else if (!isTpetral && isTpetrar){ // Assume that the left operator is diagonal
1488 
1489  // Get right Tpetra crs operator
1490  ST scalarr = 0.0;
1491  bool transpr = false;
1492  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpr = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1493 
1494  RCP<const Tpetra::Vector<ST,LO,GO,NT> > diagPtr;
1495 
1496  // Cast left operator as DiagonalLinearOp and extract diagonal as Vector
1497  RCP<const Thyra::DiagonalLinearOpBase<ST> > dOpl = rcp_dynamic_cast<const Thyra::DiagonalLinearOpBase<ST> >(opl);
1498  if(dOpl != Teuchos::null){
1499  RCP<const Thyra::TpetraVector<ST,LO,GO,NT> > tPtr = rcp_dynamic_cast<const Thyra::TpetraVector<ST,LO,GO,NT> >(dOpl->getDiag(),true);
1500  diagPtr = rcp_dynamic_cast<const Tpetra::Vector<ST,LO,GO,NT> >(tPtr->getConstTpetraVector(),true);
1501  }
1502  // If it's not diagonal, maybe it's zero
1503  else if(rcp_dynamic_cast<const Thyra::ZeroLinearOpBase<ST> >(opl) != Teuchos::null){
1504  diagPtr = rcp(new Tpetra::Vector<ST,LO,GO,NT>(tCrsOpr->getRangeMap()));
1505  }
1506  else
1507  TEUCHOS_ASSERT(false);
1508 
1509  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()));
1510 
1511  explicitCrsOp->leftScale(*diagPtr);
1512  explicitCrsOp->resumeFill();
1513  explicitCrsOp->scale(scalarr);
1514  explicitCrsOp->fillComplete(tCrsOpr->getDomainMap(),tCrsOpr->getRangeMap());
1515 
1516  return Thyra::constTpetraLinearOp<ST,LO,GO,NT>(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1517 
1518  } else { // Assume Epetra and we can use transformers
1519 
1520  // build implicit multiply
1521  const LinearOp implicitOp = Thyra::multiply(opl,opr);
1522 
1523  // build a scaling transformer
1524  RCP<Thyra::LinearOpTransformerBase<double> > prodTrans
1525  = Thyra::epetraExtDiagScalingTransformer();
1526 
1527  // check to see if a scaling transformer works: if not use the
1528  // DiagScaledMatrixProduct transformer
1529  if(not prodTrans->isCompatible(*implicitOp))
1530  prodTrans = Thyra::epetraExtDiagScaledMatProdTransformer();
1531 
1532  // build operator and multiply
1533  const RCP<Thyra::LinearOpBase<double> > explicitOp = prodTrans->createOutputOp();
1534  prodTrans->transform(*implicitOp,explicitOp.ptr());
1535  explicitOp->setObjectLabel("explicit( " + opl->getObjectLabel() +
1536  " * " + opr->getObjectLabel() + " )");
1537 
1538  return explicitOp;
1539  }
1540 }
1541 
1555 const ModifiableLinearOp explicitMultiply(const LinearOp & opl,const LinearOp & opr,
1556  const ModifiableLinearOp & destOp)
1557 {
1558  // if this is a blocked operator, multiply block by block
1559  // it is possible that not every factor in the product is blocked and these situations are handled separately
1560 
1561  bool isBlockedL = isPhysicallyBlockedLinearOp(opl);
1562  bool isBlockedR = isPhysicallyBlockedLinearOp(opr);
1563 
1564  // both factors blocked
1565  if((isBlockedL && isBlockedR)){
1566  double scalarl = 0.0;
1567  bool transpl = false;
1568  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opl = getPhysicallyBlockedLinearOp(opl,&scalarl,&transpl);
1569  double scalarr = 0.0;
1570  bool transpr = false;
1571  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opr = getPhysicallyBlockedLinearOp(opr,&scalarr,&transpr);
1572  double scalar = scalarl*scalarr;
1573 
1574  int numRows = blocked_opl->productRange()->numBlocks();
1575  int numCols = blocked_opr->productDomain()->numBlocks();
1576  int numMiddle = blocked_opl->productDomain()->numBlocks();
1577 
1578  TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == numMiddle);
1579 
1580 std::cout << "rows, cols, middle = " << numRows << ", " << numCols << ", " << numMiddle << std::endl;
1581 
1582  RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_product = Thyra::defaultBlockedLinearOp<double>();
1583  blocked_product->beginBlockFill(numRows,numCols);
1584  for(int r = 0; r < numRows; ++r){
1585  for(int c = 0; c < numCols; ++c){
1586 
1587 std::cout << "Block multiply " << r << ", " << c << std::endl;
1588 std::cout << blocked_opl->getBlock(r,0)->range()->dim() << " x " << blocked_opl->getBlock(r,0)->domain()->dim() << " times ";
1589 std::cout << blocked_opr->getBlock(0,c)->range()->dim() << " x " << blocked_opr->getBlock(0,c)->domain()->dim() << std::endl;
1590 
1591  LinearOp product_rc = explicitMultiply(blocked_opl->getBlock(r,0),blocked_opr->getBlock(0,c));
1592  for(int m = 1; m < numMiddle; ++m){
1593 std::cout << blocked_opl->getBlock(r,m)->range()->dim() << " x " << blocked_opl->getBlock(r,m)->domain()->dim() << " times ";
1594 std::cout << blocked_opr->getBlock(m,c)->range()->dim() << " x " << blocked_opr->getBlock(m,c)->domain()->dim() << std::endl;
1595  LinearOp product_m = explicitMultiply(blocked_opl->getBlock(r,m),blocked_opr->getBlock(m,c));
1596  product_rc = explicitAdd(product_rc,product_m);
1597  }
1598  blocked_product->setBlock(r,c,Thyra::scale(scalar,product_rc));
1599  }
1600  }
1601  blocked_product->endBlockFill();
1602  return blocked_product;
1603  }
1604 
1605  // only left factor blocked
1606  if((isBlockedL && !isBlockedR)){
1607  double scalarl = 0.0;
1608  bool transpl = false;
1609  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opl = getPhysicallyBlockedLinearOp(opl,&scalarl,&transpl);
1610  double scalar = scalarl;
1611 
1612  int numRows = blocked_opl->productRange()->numBlocks();
1613  int numCols = 1;
1614  int numMiddle = 1;
1615 
1616  TEUCHOS_ASSERT(blocked_opl->productDomain()->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  LinearOp product_r = explicitMultiply(blocked_opl->getBlock(r,0),opr);
1622  blocked_product->setBlock(r,0,Thyra::scale(scalar,product_r));
1623  }
1624  blocked_product->endBlockFill();
1625  return blocked_product;
1626  }
1627 
1628  // only right factor blocked
1629  if((!isBlockedL && isBlockedR)){
1630  double scalarr = 0.0;
1631  bool transpr = false;
1632  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opr = getPhysicallyBlockedLinearOp(opr,&scalarr,&transpr);
1633  double scalar = scalarr;
1634 
1635  int numRows = 1;
1636  int numCols = blocked_opr->productDomain()->numBlocks();
1637  int numMiddle = 1;
1638 
1639  TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == numMiddle);
1640 
1641  RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_product = Thyra::defaultBlockedLinearOp<double>();
1642  blocked_product->beginBlockFill(numRows,numCols);
1643  for(int c = 0; c < numCols; ++c){
1644  LinearOp product_c = explicitMultiply(opl,blocked_opr->getBlock(0,c));
1645  blocked_product->setBlock(0,c,Thyra::scale(scalar,product_c));
1646  }
1647  blocked_product->endBlockFill();
1648  return blocked_product;
1649  }
1650 
1651  bool isTpetral = Teko::TpetraHelpers::isTpetraLinearOp(opl);
1652  bool isTpetrar = Teko::TpetraHelpers::isTpetraLinearOp(opr);
1653 
1654  if(isTpetral && isTpetrar){ // Both operators are Tpetra matrices so use the explicit Tpetra matrix-matrix multiply
1655 
1656  // Get left and right Tpetra crs operators
1657  ST scalarl = 0.0;
1658  bool transpl = false;
1659  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpl = Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1660  ST scalarr = 0.0;
1661  bool transpr = false;
1662  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpr = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1663 
1664  // Build output operator
1665  RCP<Thyra::LinearOpBase<ST> > explicitOp;
1666  if(destOp!=Teuchos::null)
1667  explicitOp = destOp;
1668  else
1669  explicitOp = rcp(new Thyra::TpetraLinearOp<ST,LO,GO,NT>());
1670  RCP<Thyra::TpetraLinearOp<ST,LO,GO,NT> > tExplicitOp = rcp_dynamic_cast<Thyra::TpetraLinearOp<ST,LO,GO,NT> >(explicitOp);
1671 
1672  // Do explicit matrix-matrix multiply
1673  RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > explicitCrsOp = Tpetra::createCrsMatrix<ST,LO,GO,NT>(tCrsOpl->getRowMap());
1674  Tpetra::MatrixMatrix::Multiply<ST,LO,GO,NT>(*tCrsOpl,transpl,*tCrsOpr,transpr,*explicitCrsOp);
1675  explicitCrsOp->resumeFill();
1676  explicitCrsOp->scale(scalarl*scalarr);
1677  explicitCrsOp->fillComplete(tCrsOpr->getDomainMap(),tCrsOpl->getRangeMap());
1678  tExplicitOp->initialize(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1679  return tExplicitOp;
1680 
1681  } else if (isTpetral && !isTpetrar){ // Assume that the right operator is diagonal
1682 
1683  // Get left Tpetra crs operator
1684  ST scalarl = 0.0;
1685  bool transpl = false;
1686  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpl = Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1687 
1688  // Cast right operator as DiagonalLinearOp and extract diagonal as Vector
1689  RCP<const Thyra::DiagonalLinearOpBase<ST> > dOpr = rcp_dynamic_cast<const Thyra::DiagonalLinearOpBase<ST> >(opr);
1690  RCP<const Thyra::TpetraVector<ST,LO,GO,NT> > tPtr = rcp_dynamic_cast<const Thyra::TpetraVector<ST,LO,GO,NT> >(dOpr->getDiag(),true);
1691  RCP<const Tpetra::Vector<ST,LO,GO,NT> > diagPtr = rcp_dynamic_cast<const Tpetra::Vector<ST,LO,GO,NT> >(tPtr->getConstTpetraVector(),true);
1692 
1693  // Scale by the diagonal operator
1694  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()));
1695  explicitCrsOp->rightScale(*diagPtr);
1696  explicitCrsOp->resumeFill();
1697  explicitCrsOp->scale(scalarl);
1698  explicitCrsOp->fillComplete(tCrsOpl->getDomainMap(),tCrsOpl->getRangeMap());
1699  return Thyra::tpetraLinearOp<ST,LO,GO,NT>(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1700 
1701  } else if (!isTpetral && isTpetrar){ // Assume that the left operator is diagonal
1702 
1703  // Get right Tpetra crs operator
1704  ST scalarr = 0.0;
1705  bool transpr = false;
1706  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpr = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1707 
1708  // Cast leftt operator as DiagonalLinearOp and extract diagonal as Vector
1709  RCP<const Thyra::DiagonalLinearOpBase<ST> > dOpl = rcp_dynamic_cast<const Thyra::DiagonalLinearOpBase<ST> >(opl,true);
1710  RCP<const Thyra::TpetraVector<ST,LO,GO,NT> > tPtr = rcp_dynamic_cast<const Thyra::TpetraVector<ST,LO,GO,NT> >(dOpl->getDiag(),true);
1711  RCP<const Tpetra::Vector<ST,LO,GO,NT> > diagPtr = rcp_dynamic_cast<const Tpetra::Vector<ST,LO,GO,NT> >(tPtr->getConstTpetraVector(),true);
1712 
1713  // Scale by the diagonal operator
1714  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()));
1715  explicitCrsOp->leftScale(*diagPtr);
1716  explicitCrsOp->resumeFill();
1717  explicitCrsOp->scale(scalarr);
1718  explicitCrsOp->fillComplete(tCrsOpr->getDomainMap(),tCrsOpr->getRangeMap());
1719  return Thyra::tpetraLinearOp<ST,LO,GO,NT>(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1720 
1721  } else { // Assume Epetra and we can use transformers
1722 
1723  // build implicit multiply
1724  const LinearOp implicitOp = Thyra::multiply(opl,opr);
1725 
1726  // build a scaling transformer
1727 
1728  RCP<Thyra::LinearOpTransformerBase<double> > prodTrans
1729  = Thyra::epetraExtDiagScalingTransformer();
1730 
1731  // check to see if a scaling transformer works: if not use the
1732  // DiagScaledMatrixProduct transformer
1733  if(not prodTrans->isCompatible(*implicitOp))
1734  prodTrans = Thyra::epetraExtDiagScaledMatProdTransformer();
1735 
1736  // build operator destination operator
1737  ModifiableLinearOp explicitOp;
1738 
1739  // if neccessary build a operator to put the explicit multiply into
1740  if(destOp==Teuchos::null)
1741  explicitOp = prodTrans->createOutputOp();
1742  else
1743  explicitOp = destOp;
1744 
1745  // perform multiplication
1746  prodTrans->transform(*implicitOp,explicitOp.ptr());
1747 
1748  // label it
1749  explicitOp->setObjectLabel("explicit( " + opl->getObjectLabel() +
1750  " * " + opr->getObjectLabel() + " )");
1751 
1752  return explicitOp;
1753  }
1754 }
1755 
1766 const LinearOp explicitAdd(const LinearOp & opl,const LinearOp & opr)
1767 {
1768  // if blocked, add block by block
1769  if(isPhysicallyBlockedLinearOp(opl)){
1770  TEUCHOS_ASSERT(isPhysicallyBlockedLinearOp(opr));
1771 
1772  double scalarl = 0.0;
1773  bool transpl = false;
1774  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opl = getPhysicallyBlockedLinearOp(opl, &scalarl, &transpl);
1775 
1776  double scalarr = 0.0;
1777  bool transpr = false;
1778  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opr = getPhysicallyBlockedLinearOp(opr, &scalarr, &transpr);
1779 
1780  int numRows = blocked_opl->productRange()->numBlocks();
1781  int numCols = blocked_opl->productDomain()->numBlocks();
1782  TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == numRows);
1783  TEUCHOS_ASSERT(blocked_opr->productDomain()->numBlocks() == numCols);
1784 
1785  RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_sum = Thyra::defaultBlockedLinearOp<double>();
1786  blocked_sum->beginBlockFill(numRows,numCols);
1787  for(int r = 0; r < numRows; ++r)
1788  for(int c = 0; c < numCols; ++c)
1789  blocked_sum->setBlock(r,c,explicitAdd(Thyra::scale(scalarl,blocked_opl->getBlock(r,c)),Thyra::scale(scalarr,blocked_opr->getBlock(r,c))));
1790  blocked_sum->endBlockFill();
1791  return blocked_sum;
1792  }
1793 
1794  bool isTpetral = Teko::TpetraHelpers::isTpetraLinearOp(opl);
1795  bool isTpetrar = Teko::TpetraHelpers::isTpetraLinearOp(opr);
1796 
1797  if(isTpetral && isTpetrar){ // Both operators are Tpetra matrices so use the explicit Tpetra matrix-matrix add
1798 
1799  // Get left and right Tpetra crs operators
1800  ST scalarl = 0.0;
1801  bool transpl = false;
1802  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpl = Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1803  ST scalarr = 0.0;
1804  bool transpr = false;
1805  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpr = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1806 
1807  // Build output operator
1808  RCP<Thyra::LinearOpBase<ST> > explicitOp = rcp(new Thyra::TpetraLinearOp<ST,LO,GO,NT>());
1809  RCP<Thyra::TpetraLinearOp<ST,LO,GO,NT> > tExplicitOp = rcp_dynamic_cast<Thyra::TpetraLinearOp<ST,LO,GO,NT> >(explicitOp);
1810 
1811  // Do explicit matrix-matrix add
1812  RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > explicitCrsOp = Tpetra::MatrixMatrix::add<ST,LO,GO,NT>(scalarl,transpl,*tCrsOpl,scalarr,transpr,*tCrsOpr);
1813  tExplicitOp->initialize(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1814  return tExplicitOp;
1815 
1816  }else{//Assume Epetra
1817  // build implicit add
1818  const LinearOp implicitOp = Thyra::add(opl,opr);
1819 
1820  // build transformer
1821  const RCP<Thyra::LinearOpTransformerBase<double> > prodTrans =
1822  Thyra::epetraExtAddTransformer();
1823 
1824  // build operator and add
1825  const RCP<Thyra::LinearOpBase<double> > explicitOp = prodTrans->createOutputOp();
1826  prodTrans->transform(*implicitOp,explicitOp.ptr());
1827  explicitOp->setObjectLabel("explicit( " + opl->getObjectLabel() +
1828  " + " + opr->getObjectLabel() + " )");
1829 
1830  return explicitOp;
1831  }
1832 }
1833 
1846 const ModifiableLinearOp explicitAdd(const LinearOp & opl,const LinearOp & opr,
1847  const ModifiableLinearOp & destOp)
1848 {
1849  // if blocked, add block by block
1850  if(isPhysicallyBlockedLinearOp(opl)){
1851  TEUCHOS_ASSERT(isPhysicallyBlockedLinearOp(opr));
1852 
1853  double scalarl = 0.0;
1854  bool transpl = false;
1855  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opl = getPhysicallyBlockedLinearOp(opl, &scalarl, &transpl);
1856 
1857  double scalarr = 0.0;
1858  bool transpr = false;
1859  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opr = getPhysicallyBlockedLinearOp(opr, &scalarr, &transpr);
1860 
1861  int numRows = blocked_opl->productRange()->numBlocks();
1862  int numCols = blocked_opl->productDomain()->numBlocks();
1863  TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == numRows);
1864  TEUCHOS_ASSERT(blocked_opr->productDomain()->numBlocks() == numCols);
1865 
1866  RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_sum = Thyra::defaultBlockedLinearOp<double>();
1867  blocked_sum->beginBlockFill(numRows,numCols);
1868  for(int r = 0; r < numRows; ++r)
1869  for(int c = 0; c < numCols; ++c)
1870  blocked_sum->setBlock(r,c,explicitAdd(Thyra::scale(scalarl,blocked_opl->getBlock(r,c)),Thyra::scale(scalarr,blocked_opr->getBlock(r,c))));
1871  blocked_sum->endBlockFill();
1872  return blocked_sum;
1873  }
1874 
1875  bool isTpetral = Teko::TpetraHelpers::isTpetraLinearOp(opl);
1876  bool isTpetrar = Teko::TpetraHelpers::isTpetraLinearOp(opr);
1877 
1878  if(isTpetral && isTpetrar){ // Both operators are Tpetra matrices so use the explicit Tpetra matrix-matrix add
1879 
1880  // Get left and right Tpetra crs operators
1881  ST scalarl = 0.0;
1882  bool transpl = false;
1883  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpl = Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1884  ST scalarr = 0.0;
1885  bool transpr = false;
1886  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpr = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1887 
1888  // Build output operator
1889  RCP<Thyra::LinearOpBase<ST> > explicitOp;
1890  if(destOp!=Teuchos::null)
1891  explicitOp = destOp;
1892  else
1893  explicitOp = rcp(new Thyra::TpetraLinearOp<ST,LO,GO,NT>());
1894  RCP<Thyra::TpetraLinearOp<ST,LO,GO,NT> > tExplicitOp = rcp_dynamic_cast<Thyra::TpetraLinearOp<ST,LO,GO,NT> >(explicitOp);
1895 
1896  // Do explicit matrix-matrix add
1897  RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > explicitCrsOp = Tpetra::MatrixMatrix::add<ST,LO,GO,NT>(scalarl,transpl,*tCrsOpl,scalarr,transpr,*tCrsOpr);
1898  tExplicitOp->initialize(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1899  return tExplicitOp;
1900 
1901  }else{ // Assume Epetra
1902 
1903  // build implicit add
1904  const LinearOp implicitOp = Thyra::add(opl,opr);
1905 
1906  // build transformer
1907  const RCP<Thyra::LinearOpTransformerBase<double> > prodTrans =
1908  Thyra::epetraExtAddTransformer();
1909 
1910  // build or reuse destination operator
1911  RCP<Thyra::LinearOpBase<double> > explicitOp;
1912  if(destOp!=Teuchos::null)
1913  explicitOp = destOp;
1914  else
1915  explicitOp = prodTrans->createOutputOp();
1916 
1917  // perform add
1918  prodTrans->transform(*implicitOp,explicitOp.ptr());
1919  explicitOp->setObjectLabel("explicit( " + opl->getObjectLabel() +
1920  " + " + opr->getObjectLabel() + " )");
1921 
1922  return explicitOp;
1923  }
1924 }
1925 
1930 const ModifiableLinearOp explicitSum(const LinearOp & op,
1931  const ModifiableLinearOp & destOp)
1932 {
1933  // convert operators to Epetra_CrsMatrix
1934  const RCP<const Epetra_CrsMatrix> epetraOp =
1935  rcp_dynamic_cast<const Epetra_CrsMatrix>(get_Epetra_Operator(*op), true);
1936 
1937  if(destOp==Teuchos::null) {
1938  Teuchos::RCP<Epetra_Operator> epetraDest = Teuchos::rcp(new Epetra_CrsMatrix(*epetraOp));
1939 
1940  return Thyra::nonconstEpetraLinearOp(epetraDest);
1941  }
1942 
1943  const RCP<Epetra_CrsMatrix> epetraDest =
1944  rcp_dynamic_cast<Epetra_CrsMatrix>(get_Epetra_Operator(*destOp), true);
1945 
1946  EpetraExt::MatrixMatrix::Add(*epetraOp,false,1.0,*epetraDest,1.0);
1947 
1948  return destOp;
1949 }
1950 
1951 const LinearOp explicitTranspose(const LinearOp & op)
1952 {
1953  if(Teko::TpetraHelpers::isTpetraLinearOp(op)){
1954 
1955  RCP<const Thyra::TpetraLinearOp<ST,LO,GO,NT> > tOp = rcp_dynamic_cast<const Thyra::TpetraLinearOp<ST,LO,GO,NT> >(op,true);
1956  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOp = rcp_dynamic_cast<const Tpetra::CrsMatrix<ST,LO,GO,NT> >(tOp->getConstTpetraOperator(),true);
1957 
1958  Tpetra::RowMatrixTransposer<ST,LO,GO,NT> transposer(tCrsOp);
1959  RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > transOp = transposer.createTranspose();
1960 
1961  return Thyra::constTpetraLinearOp<ST,LO,GO,NT>(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(transOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(transOp->getDomainMap()),transOp);
1962 
1963  } else {
1964 
1965  RCP<const Epetra_Operator> eOp = Thyra::get_Epetra_Operator(*op);
1966  TEUCHOS_TEST_FOR_EXCEPTION(eOp==Teuchos::null,std::logic_error,
1967  "Teko::explicitTranspose Not an Epetra_Operator");
1968  RCP<const Epetra_RowMatrix> eRowMatrixOp
1969  = Teuchos::rcp_dynamic_cast<const Epetra_RowMatrix>(eOp);
1970  TEUCHOS_TEST_FOR_EXCEPTION(eRowMatrixOp==Teuchos::null,std::logic_error,
1971  "Teko::explicitTranspose Not an Epetra_RowMatrix");
1972 
1973  // we now have a delete transpose operator
1974  EpetraExt::RowMatrix_Transpose tranposeOp;
1975  Epetra_RowMatrix & eMat = tranposeOp(const_cast<Epetra_RowMatrix &>(*eRowMatrixOp));
1976 
1977  // this copy is because of a poor implementation of the EpetraExt::Transform
1978  // implementation
1979  Teuchos::RCP<Epetra_CrsMatrix> crsMat
1980  = Teuchos::rcp(new Epetra_CrsMatrix(dynamic_cast<Epetra_CrsMatrix &>(eMat)));
1981 
1982  return Thyra::epetraLinearOp(crsMat);
1983  }
1984 }
1985 
1986 double frobeniusNorm(const LinearOp & op)
1987 {
1988  if(Teko::TpetraHelpers::isTpetraLinearOp(op)){
1989  const RCP<const Thyra::TpetraLinearOp<ST,LO,GO,NT> > tOp = rcp_dynamic_cast<const Thyra::TpetraLinearOp<ST,LO,GO,NT> >(op);
1990  const RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > crsOp = rcp_dynamic_cast<const Tpetra::CrsMatrix<ST,LO,GO,NT> >(tOp->getConstTpetraOperator(),true);
1991  return crsOp->getFrobeniusNorm();
1992  } else {
1993  const RCP<const Epetra_Operator> epOp = Thyra::get_Epetra_Operator(*op);
1994  const RCP<const Epetra_CrsMatrix> crsOp = rcp_dynamic_cast<const Epetra_CrsMatrix>(epOp,true);
1995  return crsOp->NormFrobenius();
1996  }
1997 }
1998 
1999 double oneNorm(const LinearOp & op)
2000 {
2001  if(Teko::TpetraHelpers::isTpetraLinearOp(op)){
2002  //const RCP<const Thyra::TpetraLinearOp<ST,LO,GO,NT> > tOp = rcp_dynamic_cast<const Thyra::TpetraLinearOp<ST,LO,GO,NT> >(op);
2003  //const RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > crsOp = rcp_dynamic_cast<const Tpetra::CrsMatrix<ST,LO,GO,NT> >(tOp->getConstTpetraOperator(),true);
2004  //return crsOp->getOneNorm();
2005  TEUCHOS_TEST_FOR_EXCEPTION(true,std::logic_error,"One norm not currently implemented for Tpetra matrices");
2006  } else {
2007  const RCP<const Epetra_Operator> epOp = Thyra::get_Epetra_Operator(*op);
2008  const RCP<const Epetra_CrsMatrix> crsOp = rcp_dynamic_cast<const Epetra_CrsMatrix>(epOp,true);
2009  return crsOp->NormOne();
2010  }
2011 }
2012 
2013 double infNorm(const LinearOp & op)
2014 {
2015  if(Teko::TpetraHelpers::isTpetraLinearOp(op)){
2016  //const RCP<const Thyra::TpetraLinearOp<ST,LO,GO,NT> > tOp = rcp_dynamic_cast<const Thyra::TpetraLinearOp<ST,LO,GO,NT> >(op);
2017  //const RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > crsOp = rcp_dynamic_cast<const Tpetra::CrsMatrix<ST,LO,GO,NT> >(tOp->getConstTpetraOperator(),true);
2018  //return crsOp->getInfNorm();
2019  TEUCHOS_TEST_FOR_EXCEPTION(true,std::logic_error,"Infinity norm not currently implemented for Tpetra matrices");
2020  } else {
2021  const RCP<const Epetra_Operator> epOp = Thyra::get_Epetra_Operator(*op);
2022  const RCP<const Epetra_CrsMatrix> crsOp = rcp_dynamic_cast<const Epetra_CrsMatrix>(epOp,true);
2023  return crsOp->NormInf();
2024  }
2025 }
2026 
2027 const LinearOp buildDiagonal(const MultiVector & src,const std::string & lbl)
2028 {
2029  RCP<Thyra::VectorBase<double> > dst = Thyra::createMember(src->range());
2030  Thyra::copy(*src->col(0),dst.ptr());
2031 
2032  return Thyra::diagonal<double>(dst,lbl);
2033 }
2034 
2035 const LinearOp buildInvDiagonal(const MultiVector & src,const std::string & lbl)
2036 {
2037  const RCP<const Thyra::VectorBase<double> > srcV = src->col(0);
2038  RCP<Thyra::VectorBase<double> > dst = Thyra::createMember(srcV->range());
2039  Thyra::reciprocal<double>(*srcV,dst.ptr());
2040 
2041  return Thyra::diagonal<double>(dst,lbl);
2042 }
2043 
2045 BlockedMultiVector buildBlockedMultiVector(const std::vector<MultiVector> & mvv)
2046 {
2047  Teuchos::Array<MultiVector> mvA;
2048  Teuchos::Array<VectorSpace> vsA;
2049 
2050  // build arrays of multi vectors and vector spaces
2051  std::vector<MultiVector>::const_iterator itr;
2052  for(itr=mvv.begin();itr!=mvv.end();++itr) {
2053  mvA.push_back(*itr);
2054  vsA.push_back((*itr)->range());
2055  }
2056 
2057  // construct the product vector space
2058  const RCP<const Thyra::DefaultProductVectorSpace<double> > vs
2059  = Thyra::productVectorSpace<double>(vsA);
2060 
2061  return Thyra::defaultProductMultiVector<double>(vs,mvA);
2062 }
2063 
2074 Teuchos::RCP<Thyra::VectorBase<double> > indicatorVector(
2075  const std::vector<int> & indices,
2076  const VectorSpace & vs,
2077  double onValue, double offValue)
2078 
2079 {
2080  using Teuchos::RCP;
2081 
2082  // create a new vector
2083  RCP<Thyra::VectorBase<double> > v = Thyra::createMember<double>(vs);
2084  Thyra::put_scalar<double>(offValue,v.ptr()); // fill it with "off" values
2085 
2086  // set on values
2087  for(std::size_t i=0;i<indices.size();i++)
2088  Thyra::set_ele<double>(indices[i],onValue,v.ptr());
2089 
2090  return v;
2091 }
2092 
2117 double computeSpectralRad(const RCP<const Thyra::LinearOpBase<double> > & A, double tol,
2118  bool isHermitian,int numBlocks,int restart,int verbosity)
2119 {
2120  typedef Thyra::LinearOpBase<double> OP;
2121  typedef Thyra::MultiVectorBase<double> MV;
2122 
2123  int startVectors = 1;
2124 
2125  // construct an initial guess
2126  const RCP<MV> ivec = Thyra::createMember(A->domain());
2127  Thyra::randomize(-1.0,1.0,ivec.ptr());
2128 
2129  RCP<Anasazi::BasicEigenproblem<double,MV,OP> > eigProb
2130  = rcp(new Anasazi::BasicEigenproblem<double,MV,OP>(A,ivec));
2131  eigProb->setNEV(1);
2132  eigProb->setHermitian(isHermitian);
2133 
2134  // set the problem up
2135  if(not eigProb->setProblem()) {
2136  // big time failure!
2137  return Teuchos::ScalarTraits<double>::nan();
2138  }
2139 
2140  // we want largert magnitude eigenvalue
2141  std::string which("LM"); // largest magnitude
2142 
2143  // Create the parameter list for the eigensolver
2144  // verbosity+=Anasazi::TimingDetails;
2145  Teuchos::ParameterList MyPL;
2146  MyPL.set( "Verbosity", verbosity );
2147  MyPL.set( "Which", which );
2148  MyPL.set( "Block Size", startVectors );
2149  MyPL.set( "Num Blocks", numBlocks );
2150  MyPL.set( "Maximum Restarts", restart );
2151  MyPL.set( "Convergence Tolerance", tol );
2152 
2153  // build status test manager
2154  // RCP<Anasazi::StatusTestMaxIters<double,MV,OP> > statTest
2155  // = rcp(new Anasazi::StatusTestMaxIters<double,MV,OP>(numBlocks*(restart+1)));
2156 
2157  // Create the Block Krylov Schur solver
2158  // This takes as inputs the eigenvalue problem and the solver parameters
2159  Anasazi::BlockKrylovSchurSolMgr<double,MV,OP> MyBlockKrylovSchur(eigProb, MyPL );
2160 
2161  // Solve the eigenvalue problem, and save the return code
2162  Anasazi::ReturnType solverreturn = MyBlockKrylovSchur.solve();
2163 
2164  if(solverreturn==Anasazi::Unconverged) {
2165  double real = MyBlockKrylovSchur.getRitzValues().begin()->realpart;
2166  double comp = MyBlockKrylovSchur.getRitzValues().begin()->imagpart;
2167 
2168  return -std::abs(std::sqrt(real*real+comp*comp));
2169 
2170  // cout << "Anasazi::BlockKrylovSchur::solve() did not converge!" << std::endl;
2171  // return -std::abs(MyBlockKrylovSchur.getRitzValues().begin()->realpart);
2172  }
2173  else { // solverreturn==Anasazi::Converged
2174  double real = eigProb->getSolution().Evals.begin()->realpart;
2175  double comp = eigProb->getSolution().Evals.begin()->imagpart;
2176 
2177  return std::abs(std::sqrt(real*real+comp*comp));
2178 
2179  // cout << "Anasazi::BlockKrylovSchur::solve() converged!" << endl;
2180  // return std::abs(eigProb->getSolution().Evals.begin()->realpart);
2181  }
2182 }
2183 
2207 double computeSmallestMagEig(const RCP<const Thyra::LinearOpBase<double> > & A, double tol,
2208  bool isHermitian,int numBlocks,int restart,int verbosity)
2209 {
2210  typedef Thyra::LinearOpBase<double> OP;
2211  typedef Thyra::MultiVectorBase<double> MV;
2212 
2213  int startVectors = 1;
2214 
2215  // construct an initial guess
2216  const RCP<MV> ivec = Thyra::createMember(A->domain());
2217  Thyra::randomize(-1.0,1.0,ivec.ptr());
2218 
2219  RCP<Anasazi::BasicEigenproblem<double,MV,OP> > eigProb
2220  = rcp(new Anasazi::BasicEigenproblem<double,MV,OP>(A,ivec));
2221  eigProb->setNEV(1);
2222  eigProb->setHermitian(isHermitian);
2223 
2224  // set the problem up
2225  if(not eigProb->setProblem()) {
2226  // big time failure!
2227  return Teuchos::ScalarTraits<double>::nan();
2228  }
2229 
2230  // we want largert magnitude eigenvalue
2231  std::string which("SM"); // smallest magnitude
2232 
2233  // Create the parameter list for the eigensolver
2234  Teuchos::ParameterList MyPL;
2235  MyPL.set( "Verbosity", verbosity );
2236  MyPL.set( "Which", which );
2237  MyPL.set( "Block Size", startVectors );
2238  MyPL.set( "Num Blocks", numBlocks );
2239  MyPL.set( "Maximum Restarts", restart );
2240  MyPL.set( "Convergence Tolerance", tol );
2241 
2242  // build status test manager
2243  // RCP<Anasazi::StatusTestMaxIters<double,MV,OP> > statTest
2244  // = rcp(new Anasazi::StatusTestMaxIters<double,MV,OP>(10));
2245 
2246  // Create the Block Krylov Schur solver
2247  // This takes as inputs the eigenvalue problem and the solver parameters
2248  Anasazi::BlockKrylovSchurSolMgr<double,MV,OP> MyBlockKrylovSchur(eigProb, MyPL );
2249 
2250  // Solve the eigenvalue problem, and save the return code
2251  Anasazi::ReturnType solverreturn = MyBlockKrylovSchur.solve();
2252 
2253  if(solverreturn==Anasazi::Unconverged) {
2254  // cout << "Anasazi::BlockKrylovSchur::solve() did not converge! " << std::endl;
2255  return -std::abs(MyBlockKrylovSchur.getRitzValues().begin()->realpart);
2256  }
2257  else { // solverreturn==Anasazi::Converged
2258  // cout << "Anasazi::BlockKrylovSchur::solve() converged!" << endl;
2259  return std::abs(eigProb->getSolution().Evals.begin()->realpart);
2260  }
2261 }
2262 
2271 ModifiableLinearOp getDiagonalOp(const Teko::LinearOp & A,const DiagonalType & dt)
2272 {
2273  switch(dt) {
2274  case Diagonal:
2275  return getDiagonalOp(A);
2276  case Lumped:
2277  return getLumpedMatrix(A);
2278  case AbsRowSum:
2279  return getAbsRowSumMatrix(A);
2280  case NotDiag:
2281  default:
2282  TEUCHOS_TEST_FOR_EXCEPT(true);
2283  };
2284 
2285  return Teuchos::null;
2286 }
2287 
2296 ModifiableLinearOp getInvDiagonalOp(const Teko::LinearOp & A,const Teko::DiagonalType & dt)
2297 {
2298  switch(dt) {
2299  case Diagonal:
2300  return getInvDiagonalOp(A);
2301  case Lumped:
2302  return getInvLumpedMatrix(A);
2303  case AbsRowSum:
2304  return getAbsRowSumInvMatrix(A);
2305  case NotDiag:
2306  default:
2307  TEUCHOS_TEST_FOR_EXCEPT(true);
2308  };
2309 
2310  return Teuchos::null;
2311 }
2312 
2319 std::string getDiagonalName(const DiagonalType & dt)
2320 {
2321  switch(dt) {
2322  case Diagonal:
2323  return "Diagonal";
2324  case Lumped:
2325  return "Lumped";
2326  case AbsRowSum:
2327  return "AbsRowSum";
2328  case NotDiag:
2329  return "NotDiag";
2330  case BlkDiag:
2331  return "BlkDiag";
2332  };
2333 
2334  return "<error>";
2335 }
2336 
2345 DiagonalType getDiagonalType(std::string name)
2346 {
2347  if(name=="Diagonal")
2348  return Diagonal;
2349  if(name=="Lumped")
2350  return Lumped;
2351  if(name=="AbsRowSum")
2352  return AbsRowSum;
2353  if(name=="BlkDiag")
2354  return BlkDiag;
2355 
2356  return NotDiag;
2357 }
2358 
2359 LinearOp probe(Teuchos::RCP<const Epetra_CrsGraph> &G,const LinearOp & Op){
2360 #ifdef Teko_ENABLE_Isorropia
2361  Teuchos::ParameterList probeList;
2362  Prober prober(G,probeList,true);
2363  Teuchos::RCP<Epetra_CrsMatrix> Mat=rcp(new Epetra_CrsMatrix(Copy,*G));
2365  prober.probe(Mwrap,*Mat);
2366  return Thyra::epetraLinearOp(Mat);
2367 #else
2368  TEUCHOS_TEST_FOR_EXCEPTION(true,std::runtime_error,"Probe requires Isorropia");
2369 #endif
2370 }
2371 
2372 double norm_1(const MultiVector & v,std::size_t col)
2373 {
2374  Teuchos::Array<double> n(v->domain()->dim());
2375  Thyra::norms_1<double>(*v,n);
2376 
2377  return n[col];
2378 }
2379 
2380 double norm_2(const MultiVector & v,std::size_t col)
2381 {
2382  Teuchos::Array<double> n(v->domain()->dim());
2383  Thyra::norms_2<double>(*v,n);
2384 
2385  return n[col];
2386 }
2387 
2388 void putScalar(const ModifiableLinearOp & op,double scalar)
2389 {
2390  try {
2391  // get Epetra_Operator
2392  RCP<Epetra_Operator> eOp = Thyra::get_Epetra_Operator(*op);
2393 
2394  // cast it to a CrsMatrix
2395  RCP<Epetra_CrsMatrix> eCrsOp = rcp_dynamic_cast<Epetra_CrsMatrix>(eOp,true);
2396 
2397  eCrsOp->PutScalar(scalar);
2398  }
2399  catch (std::exception & e) {
2400  RCP<Teuchos::FancyOStream> out = Teuchos::VerboseObjectBase::getDefaultOStream();
2401 
2402  *out << "Teko: putScalar requires an Epetra_CrsMatrix\n";
2403  *out << " Could not extract an Epetra_Operator from a \"" << op->description() << std::endl;
2404  *out << " OR\n";
2405  *out << " Could not cast an Epetra_Operator to a Epetra_CrsMatrix\n";
2406  *out << std::endl;
2407  *out << "*** THROWN EXCEPTION ***\n";
2408  *out << e.what() << std::endl;
2409  *out << "************************\n";
2410 
2411  throw e;
2412  }
2413 }
2414 
2415 void clipLower(MultiVector & v,double lowerBound)
2416 {
2417  using Teuchos::RCP;
2418  using Teuchos::rcp_dynamic_cast;
2419 
2420  // cast so entries are accessible
2421  // RCP<Thyra::SpmdMultiVectorBase<double> > spmdMVec
2422  // = rcp_dynamic_cast<Thyra::DefaultSpmdMultiVector<double> >(v);
2423 
2424  for(Thyra::Ordinal i=0;i<v->domain()->dim();i++) {
2425  RCP<Thyra::SpmdVectorBase<double> > spmdVec
2426  = rcp_dynamic_cast<Thyra::SpmdVectorBase<double> >(v->col(i),true);
2427 
2428  Teuchos::ArrayRCP<double> values;
2429  // spmdMVec->getNonconstLocalData(Teuchos::ptrFromRef(values),Teuchos::ptrFromRef(i));
2430  spmdVec->getNonconstLocalData(Teuchos::ptrFromRef(values));
2431  for(Teuchos::ArrayRCP<double>::size_type j=0;j<values.size();j++) {
2432  if(values[j]<lowerBound)
2433  values[j] = lowerBound;
2434  }
2435  }
2436 }
2437 
2438 void clipUpper(MultiVector & v,double upperBound)
2439 {
2440  using Teuchos::RCP;
2441  using Teuchos::rcp_dynamic_cast;
2442 
2443  // cast so entries are accessible
2444  // RCP<Thyra::SpmdMultiVectorBase<double> > spmdMVec
2445  // = rcp_dynamic_cast<Thyra::DefaultSpmdMultiVector<double> >(v);
2446  for(Thyra::Ordinal i=0;i<v->domain()->dim();i++) {
2447  RCP<Thyra::SpmdVectorBase<double> > spmdVec
2448  = rcp_dynamic_cast<Thyra::SpmdVectorBase<double> >(v->col(i),true);
2449 
2450  Teuchos::ArrayRCP<double> values;
2451  // spmdMVec->getNonconstLocalData(Teuchos::ptrFromRef(values),Teuchos::ptrFromRef(i));
2452  spmdVec->getNonconstLocalData(Teuchos::ptrFromRef(values));
2453  for(Teuchos::ArrayRCP<double>::size_type j=0;j<values.size();j++) {
2454  if(values[j]>upperBound)
2455  values[j] = upperBound;
2456  }
2457  }
2458 }
2459 
2460 void replaceValue(MultiVector & v,double currentValue,double newValue)
2461 {
2462  using Teuchos::RCP;
2463  using Teuchos::rcp_dynamic_cast;
2464 
2465  // cast so entries are accessible
2466  // RCP<Thyra::SpmdMultiVectorBase<double> > spmdMVec
2467  // = rcp_dynamic_cast<Thyra::SpmdMultiVectorBase<double> >(v,true);
2468  for(Thyra::Ordinal i=0;i<v->domain()->dim();i++) {
2469  RCP<Thyra::SpmdVectorBase<double> > spmdVec
2470  = rcp_dynamic_cast<Thyra::SpmdVectorBase<double> >(v->col(i),true);
2471 
2472  Teuchos::ArrayRCP<double> values;
2473  // spmdMVec->getNonconstLocalData(Teuchos::ptrFromRef(values),Teuchos::ptrFromRef(i));
2474  spmdVec->getNonconstLocalData(Teuchos::ptrFromRef(values));
2475  for(Teuchos::ArrayRCP<double>::size_type j=0;j<values.size();j++) {
2476  if(values[j]==currentValue)
2477  values[j] = newValue;
2478  }
2479  }
2480 }
2481 
2482 void columnAverages(const MultiVector & v,std::vector<double> & averages)
2483 {
2484  averages.resize(v->domain()->dim());
2485 
2486  // sum over each column
2487  Thyra::sums<double>(*v,averages);
2488 
2489  // build averages
2490  Thyra::Ordinal rows = v->range()->dim();
2491  for(std::size_t i=0;i<averages.size();i++)
2492  averages[i] = averages[i]/rows;
2493 }
2494 
2495 double average(const MultiVector & v)
2496 {
2497  Thyra::Ordinal rows = v->range()->dim();
2498  Thyra::Ordinal cols = v->domain()->dim();
2499 
2500  std::vector<double> averages;
2501  columnAverages(v,averages);
2502 
2503  double sum = 0.0;
2504  for(std::size_t i=0;i<averages.size();i++)
2505  sum += averages[i] * rows;
2506 
2507  return sum/(rows*cols);
2508 }
2509 
2510 bool isPhysicallyBlockedLinearOp(const LinearOp & op)
2511 {
2512  // See if the operator is a PBLO
2513  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > pblo = rcp_dynamic_cast<const Thyra::PhysicallyBlockedLinearOpBase<double> >(op);
2514  if (!pblo.is_null())
2515  return true;
2516 
2517  // See if the operator is a wrapped PBLO
2518  ST scalar = 0.0;
2519  Thyra::EOpTransp transp = Thyra::NOTRANS;
2520  RCP<const Thyra::LinearOpBase<ST> > wrapped_op;
2521  Thyra::unwrap(op, &scalar, &transp, &wrapped_op);
2522  pblo = rcp_dynamic_cast<const Thyra::PhysicallyBlockedLinearOpBase<double> >(wrapped_op);
2523  if (!pblo.is_null())
2524  return true;
2525 
2526  return false;
2527 }
2528 
2529 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > getPhysicallyBlockedLinearOp(const LinearOp & op, ST *scalar, bool *transp)
2530 {
2531  // If the operator is a TpetraLinearOp
2532  RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > pblo = rcp_dynamic_cast<const Thyra::PhysicallyBlockedLinearOpBase<double> >(op);
2533  if(!pblo.is_null()){
2534  *scalar = 1.0;
2535  *transp = false;
2536  return pblo;
2537  }
2538 
2539  // If the operator is a wrapped TpetraLinearOp
2540  RCP<const Thyra::LinearOpBase<ST> > wrapped_op;
2541  Thyra::EOpTransp eTransp = Thyra::NOTRANS;
2542  Thyra::unwrap(op, scalar, &eTransp, &wrapped_op);
2543  pblo = rcp_dynamic_cast<const Thyra::PhysicallyBlockedLinearOpBase<double> >(wrapped_op,true);
2544  if(!pblo.is_null()){
2545  *transp = true;
2546  if(eTransp == Thyra::NOTRANS)
2547  *transp = false;
2548  return pblo;
2549  }
2550 
2551  return Teuchos::null;
2552 }
2553 
2554 
2555 }
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 .