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
112namespace Teko {
113
114using Teuchos::rcp;
115using Teuchos::rcp_dynamic_cast;
116using Teuchos::RCP;
117#ifdef Teko_ENABLE_Isorropia
118using Isorropia::Epetra::Prober;
119#endif
120
121const 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
132inline 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
143inline 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
172RCP<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
230RCP<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 auto rowInd = typename Tpetra::CrsMatrix<ST,LO,GO,NT>::nonconst_global_inds_host_view_type(Kokkos::ViewAllocateWithoutInitializing("rowIndices"),stencil.getGlobalMaxNumRowEntries()+1);
238 auto rowData = typename Tpetra::CrsMatrix<ST,LO,GO,NT>::nonconst_values_host_view_type(Kokkos::ViewAllocateWithoutInitializing("rowIndices"),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,rowInd,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,rowInd,rowData);
281 }
282
283 gl->fillComplete();
284
285 return gl;
286}
287
310RCP<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
368RCP<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 auto rowInd = typename Tpetra::CrsMatrix<ST,LO,GO,NT>::nonconst_global_inds_host_view_type(Kokkos::ViewAllocateWithoutInitializing("rowIndices"),stencil.getGlobalMaxNumRowEntries()+1);
376 auto rowData = typename Tpetra::CrsMatrix<ST,LO,GO,NT>::nonconst_values_host_view_type(Kokkos::ViewAllocateWithoutInitializing("rowIndices"),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,rowInd,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,rowInd,rowData);
419 }
420
421 gl->fillComplete();
422
423 return gl;
424}
425
441void 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
461void 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
468void 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
482BlockedLinearOp 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
521BlockedLinearOp 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
578BlockedLinearOp 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
606bool isZeroOp(const LinearOp op)
607{
608 // if operator is null...then its zero!
609 if(op==Teuchos::null) return true;
610
611 // try to cast it to a zero linear operator
612 LinearOp test = rcp_dynamic_cast<const Thyra::ZeroLinearOpBase<double> >(op);
613
614 // if it works...then its zero...otherwise its null
615 if(test!=Teuchos::null) return true;
616
617 // See if the operator is a wrapped zero op
618 ST scalar = 0.0;
619 Thyra::EOpTransp transp = Thyra::NOTRANS;
620 RCP<const Thyra::LinearOpBase<ST> > wrapped_op;
621 Thyra::unwrap(op, &scalar, &transp, &wrapped_op);
622 test = rcp_dynamic_cast<const Thyra::ZeroLinearOpBase<double> >(wrapped_op);
623 return test!=Teuchos::null;
624}
625
634ModifiableLinearOp getAbsRowSumMatrix(const LinearOp & op)
635{
636 bool isTpetra = false;
637 RCP<const Epetra_CrsMatrix> eCrsOp;
638 RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOp;
639
640 try {
641 // get Epetra or Tpetra Operator
642 RCP<const Thyra::EpetraLinearOp > eOp = rcp_dynamic_cast<const Thyra::EpetraLinearOp >(op);
643 RCP<const Thyra::TpetraLinearOp<ST,LO,GO,NT> > tOp = rcp_dynamic_cast<const Thyra::TpetraLinearOp<ST,LO,GO,NT> >(op);
644
645 // cast it to a CrsMatrix
646 RCP<Teuchos::FancyOStream> out = Teuchos::VerboseObjectBase::getDefaultOStream();
647 if (!eOp.is_null()){
648 eCrsOp = rcp_dynamic_cast<const Epetra_CrsMatrix>(eOp->epetra_op(),true);
649 }
650 else if (!tOp.is_null()){
651 tCrsOp = rcp_dynamic_cast<const Tpetra::CrsMatrix<ST,LO,GO,NT> >(tOp->getConstTpetraOperator(),true);
652 isTpetra = true;
653 }
654 else
655 throw std::logic_error("Neither Epetra nor Tpetra");
656 }
657 catch (std::exception & e) {
658 RCP<Teuchos::FancyOStream> out = Teuchos::VerboseObjectBase::getDefaultOStream();
659
660 *out << "Teko: getAbsRowSumMatrix requires an Epetra_CrsMatrix or a Tpetra::CrsMatrix\n";
661 *out << " Could not extract an Epetra_Operator or a Tpetra_Operator from a \"" << op->description() << std::endl;
662 *out << " OR\n";
663 *out << " Could not cast an Epetra_Operator to a Epetra_CrsMatrix or a Tpetra_Operator to a Tpetra::CrsMatrix\n";
664 *out << std::endl;
665 *out << "*** THROWN EXCEPTION ***\n";
666 *out << e.what() << std::endl;
667 *out << "************************\n";
668
669 throw e;
670 }
671
672 if(!isTpetra){
673 // extract diagonal
674 const RCP<Epetra_Vector> ptrDiag = rcp(new Epetra_Vector(eCrsOp->RowMap()));
675 Epetra_Vector & diag = *ptrDiag;
676
677 // compute absolute value row sum
678 diag.PutScalar(0.0);
679 for(int i=0;i<eCrsOp->NumMyRows();i++) {
680 double * values = 0;
681 int numEntries;
682 eCrsOp->ExtractMyRowView(i,numEntries,values);
683
684 // build abs value row sum
685 for(int j=0;j<numEntries;j++)
686 diag[i] += std::abs(values[j]);
687 }
688
689 // build Thyra diagonal operator
690 return Teko::Epetra::thyraDiagOp(ptrDiag,eCrsOp->RowMap(),"absRowSum( " + op->getObjectLabel() + " )");
691 }
692 else {
693 // extract diagonal
694 const RCP<Tpetra::Vector<ST,LO,GO,NT> > ptrDiag = Tpetra::createVector<ST,LO,GO,NT>(tCrsOp->getRowMap());
695 Tpetra::Vector<ST,LO,GO,NT> & diag = *ptrDiag;
696
697 // compute absolute value row sum
698 diag.putScalar(0.0);
699 for(LO i=0;i<(LO) tCrsOp->getNodeNumRows();i++) {
700 LO numEntries = tCrsOp->getNumEntriesInLocalRow (i);
701 typename Tpetra::CrsMatrix<ST,LO,GO,NT>::local_inds_host_view_type indices;
702 typename Tpetra::CrsMatrix<ST,LO,GO,NT>::values_host_view_type values;
703 tCrsOp->getLocalRowView(i,indices,values);
704
705 // build abs value row sum
706 for(LO j=0;j<numEntries;j++)
707 diag.sumIntoLocalValue(i,std::abs(values(j)));
708 }
709
710 // build Thyra diagonal operator
711 return Teko::TpetraHelpers::thyraDiagOp(ptrDiag,*tCrsOp->getRowMap(),"absRowSum( " + op->getObjectLabel() + " ))");
712
713 }
714
715}
716
725ModifiableLinearOp getAbsRowSumInvMatrix(const LinearOp & op)
726{
727 // if this is a blocked operator, extract diagonals block by block
728 // FIXME: this does not add in values from off-diagonal blocks
729 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_op = rcp_dynamic_cast<const Thyra::PhysicallyBlockedLinearOpBase<double> >(op);
730 if(blocked_op != Teuchos::null){
731 int numRows = blocked_op->productRange()->numBlocks();
732 TEUCHOS_ASSERT(blocked_op->productDomain()->numBlocks() == numRows);
733 RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_diag = Thyra::defaultBlockedLinearOp<double>();
734 blocked_diag->beginBlockFill(numRows,numRows);
735 for(int r = 0; r < numRows; ++r){
736 for(int c = 0; c < numRows; ++c){
737 if(r==c)
738 blocked_diag->setNonconstBlock(r,c,getAbsRowSumInvMatrix(blocked_op->getBlock(r,c)));
739 else
740 blocked_diag->setBlock(r,c,Thyra::zero<double>(blocked_op->getBlock(r,c)->range(),blocked_op->getBlock(r,c)->domain()));
741 }
742 }
743 blocked_diag->endBlockFill();
744 return blocked_diag;
745 }
746
747 if(Teko::TpetraHelpers::isTpetraLinearOp(op)) {
748 ST scalar = 0.0;
749 bool transp = false;
750 RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOp = Teko::TpetraHelpers::getTpetraCrsMatrix(op, &scalar, &transp);
751
752 // extract diagonal
753 const RCP<Tpetra::Vector<ST,LO,GO,NT> > ptrDiag = Tpetra::createVector<ST,LO,GO,NT>(tCrsOp->getRowMap());
754 Tpetra::Vector<ST,LO,GO,NT> & diag = *ptrDiag;
755
756 // compute absolute value row sum
757 diag.putScalar(0.0);
758 for(LO i=0;i<(LO) tCrsOp->getNodeNumRows();i++) {
759 LO numEntries = tCrsOp->getNumEntriesInLocalRow (i);
760 typename Tpetra::CrsMatrix<ST,LO,GO,NT>::local_inds_host_view_type indices;
761 typename Tpetra::CrsMatrix<ST,LO,GO,NT>::values_host_view_type values;
762 tCrsOp->getLocalRowView(i,indices,values);
763
764 // build abs value row sum
765 for(LO j=0;j<numEntries;j++)
766 diag.sumIntoLocalValue(i,std::abs(values(j)));
767 }
768 diag.scale(scalar);
769 diag.reciprocal(diag); // invert entries
770
771 // build Thyra diagonal operator
772 return Teko::TpetraHelpers::thyraDiagOp(ptrDiag,*tCrsOp->getRowMap(),"absRowSum( " + op->getObjectLabel() + " ))");
773
774 }
775 else{
776 RCP<const Thyra::EpetraLinearOp > eOp = rcp_dynamic_cast<const Thyra::EpetraLinearOp >(op,true);
777 RCP<const Epetra_CrsMatrix> eCrsOp = rcp_dynamic_cast<const Epetra_CrsMatrix>(eOp->epetra_op(),true);
778
779 // extract diagonal
780 const RCP<Epetra_Vector> ptrDiag = rcp(new Epetra_Vector(eCrsOp->RowMap()));
781 Epetra_Vector & diag = *ptrDiag;
782
783 // compute absolute value row sum
784 diag.PutScalar(0.0);
785 for(int i=0;i<eCrsOp->NumMyRows();i++) {
786 double * values = 0;
787 int numEntries;
788 eCrsOp->ExtractMyRowView(i,numEntries,values);
789
790 // build abs value row sum
791 for(int j=0;j<numEntries;j++)
792 diag[i] += std::abs(values[j]);
793 }
794 diag.Reciprocal(diag); // invert entries
795
796 // build Thyra diagonal operator
797 return Teko::Epetra::thyraDiagOp(ptrDiag,eCrsOp->RowMap(),"absRowSum( " + op->getObjectLabel() + " )");
798 }
799
800}
801
809ModifiableLinearOp getLumpedMatrix(const LinearOp & op)
810{
811 RCP<Thyra::VectorBase<ST> > ones = Thyra::createMember(op->domain());
812 RCP<Thyra::VectorBase<ST> > diag = Thyra::createMember(op->range());
813
814 // set to all ones
815 Thyra::assign(ones.ptr(),1.0);
816
817 // compute lumped diagonal
818 // Thyra::apply(*op,Thyra::NONCONJ_ELE,*ones,&*diag);
819 Thyra::apply(*op,Thyra::NOTRANS,*ones,diag.ptr());
820
821 return rcp(new Thyra::DefaultDiagonalLinearOp<ST>(diag));
822}
823
832ModifiableLinearOp getInvLumpedMatrix(const LinearOp & op)
833{
834 RCP<Thyra::VectorBase<ST> > ones = Thyra::createMember(op->domain());
835 RCP<Thyra::VectorBase<ST> > diag = Thyra::createMember(op->range());
836
837 // set to all ones
838 Thyra::assign(ones.ptr(),1.0);
839
840 // compute lumped diagonal
841 Thyra::apply(*op,Thyra::NOTRANS,*ones,diag.ptr());
842 Thyra::reciprocal(*diag,diag.ptr());
843
844 return rcp(new Thyra::DefaultDiagonalLinearOp<ST>(diag));
845}
846
858const ModifiableLinearOp getDiagonalOp(const LinearOp & op)
859{
860 bool isTpetra = false;
861 RCP<const Epetra_CrsMatrix> eCrsOp;
862 RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOp;
863
864 try {
865 // get Epetra or Tpetra Operator
866 RCP<const Thyra::EpetraLinearOp > eOp = rcp_dynamic_cast<const Thyra::EpetraLinearOp >(op);
867 RCP<const Thyra::TpetraLinearOp<ST,LO,GO,NT> > tOp = rcp_dynamic_cast<const Thyra::TpetraLinearOp<ST,LO,GO,NT> >(op);
868
869 // cast it to a CrsMatrix
870 RCP<Teuchos::FancyOStream> out = Teuchos::VerboseObjectBase::getDefaultOStream();
871 if (!eOp.is_null()){
872 eCrsOp = rcp_dynamic_cast<const Epetra_CrsMatrix>(eOp->epetra_op(),true);
873 }
874 else if (!tOp.is_null()){
875 tCrsOp = rcp_dynamic_cast<const Tpetra::CrsMatrix<ST,LO,GO,NT> >(tOp->getConstTpetraOperator(),true);
876 isTpetra = true;
877 }
878 else
879 throw std::logic_error("Neither Epetra nor Tpetra");
880 }
881 catch (std::exception & e) {
882 RCP<Teuchos::FancyOStream> out = Teuchos::VerboseObjectBase::getDefaultOStream();
883
884 *out << "Teko: getDiagonalOp requires an Epetra_CrsMatrix or a Tpetra::CrsMatrix\n";
885 *out << " Could not extract an Epetra_Operator or a Tpetra_Operator from a \"" << op->description() << std::endl;
886 *out << " OR\n";
887 *out << " Could not cast an Epetra_Operator to a Epetra_CrsMatrix or a Tpetra_Operator to a Tpetra::CrsMatrix\n";
888 *out << std::endl;
889 *out << "*** THROWN EXCEPTION ***\n";
890 *out << e.what() << std::endl;
891 *out << "************************\n";
892
893 throw e;
894 }
895
896 if(!isTpetra){
897 // extract diagonal
898 const RCP<Epetra_Vector> diag = rcp(new Epetra_Vector(eCrsOp->RowMap()));
899 TEUCHOS_TEST_FOR_EXCEPT(eCrsOp->ExtractDiagonalCopy(*diag));
900
901 // build Thyra diagonal operator
902 return Teko::Epetra::thyraDiagOp(diag,eCrsOp->RowMap(),"inv(diag( " + op->getObjectLabel() + " ))");
903 }
904 else {
905 // extract diagonal
906 const RCP<Tpetra::Vector<ST,LO,GO,NT> > diag = Tpetra::createVector<ST,LO,GO,NT>(tCrsOp->getRowMap());
907 tCrsOp->getLocalDiagCopy(*diag);
908
909 // build Thyra diagonal operator
910 return Teko::TpetraHelpers::thyraDiagOp(diag,*tCrsOp->getRowMap(),"inv(diag( " + op->getObjectLabel() + " ))");
911
912 }
913}
914
915const MultiVector getDiagonal(const LinearOp & op)
916{
917 bool isTpetra = false;
918 RCP<const Epetra_CrsMatrix> eCrsOp;
919 RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOp;
920
921 try {
922 // get Epetra or Tpetra Operator
923 RCP<const Thyra::EpetraLinearOp > eOp = rcp_dynamic_cast<const Thyra::EpetraLinearOp >(op);
924 RCP<const Thyra::TpetraLinearOp<ST,LO,GO,NT> > tOp = rcp_dynamic_cast<const Thyra::TpetraLinearOp<ST,LO,GO,NT> >(op);
925
926 // cast it to a CrsMatrix
927 RCP<Teuchos::FancyOStream> out = Teuchos::VerboseObjectBase::getDefaultOStream();
928 if (!eOp.is_null()){
929 eCrsOp = rcp_dynamic_cast<const Epetra_CrsMatrix>(eOp->epetra_op(),true);
930 }
931 else if (!tOp.is_null()){
932 tCrsOp = rcp_dynamic_cast<const Tpetra::CrsMatrix<ST,LO,GO,NT> >(tOp->getConstTpetraOperator(),true);
933 isTpetra = true;
934 }
935 else
936 throw std::logic_error("Neither Epetra nor Tpetra");
937 }
938 catch (std::exception & e) {
939 RCP<Teuchos::FancyOStream> out = Teuchos::VerboseObjectBase::getDefaultOStream();
940
941 RCP<const Thyra::EpetraLinearOp > eOp = rcp_dynamic_cast<const Thyra::EpetraLinearOp >(op);
942 RCP<const Thyra::TpetraLinearOp<ST,LO,GO,NT> > tOp = rcp_dynamic_cast<const Thyra::TpetraLinearOp<ST,LO,GO,NT> >(op);
943 *out << eOp;
944 *out << tOp;
945
946 *out << "Teko: getDiagonal requires an Epetra_CrsMatrix or a Tpetra::CrsMatrix\n";
947 *out << " Could not extract an Epetra_Operator or a Tpetra_Operator from a \"" << op->description() << std::endl;
948 *out << " OR\n";
949 *out << " Could not cast an Epetra_Operator to a Epetra_CrsMatrix or a Tpetra_Operator to a Tpetra::CrsMatrix\n";
950 *out << std::endl;
951 *out << "*** THROWN EXCEPTION ***\n";
952 *out << e.what() << std::endl;
953 *out << "************************\n";
954
955 throw e;
956 }
957
958 if(!isTpetra){
959 // extract diagonal
960 const RCP<Epetra_Vector> diag = rcp(new Epetra_Vector(eCrsOp->RowMap()));
961 TEUCHOS_TEST_FOR_EXCEPT(eCrsOp->ExtractDiagonalCopy(*diag));
962
963 return Thyra::create_Vector(diag,Thyra::create_VectorSpace(Teuchos::rcpFromRef(eCrsOp->RowMap())));
964 }
965 else {
966 // extract diagonal
967 const RCP<Tpetra::Vector<ST,LO,GO,NT> > diag = Tpetra::createVector<ST,LO,GO,NT>(tCrsOp->getRowMap());
968 tCrsOp->getLocalDiagCopy(*diag);
969
970 // build Thyra diagonal operator
971 return Thyra::createVector<ST,LO,GO,NT>(diag,Thyra::createVectorSpace<ST,LO,GO,NT>(tCrsOp->getRowMap()));
972
973 }
974}
975
976const MultiVector getDiagonal(const Teko::LinearOp & A,const DiagonalType & dt)
977{
978 LinearOp diagOp = Teko::getDiagonalOp(A,dt);
979
980 Teuchos::RCP<const Thyra::MultiVectorBase<double> > v =
981 Teuchos::rcp_dynamic_cast<const Thyra::DiagonalLinearOpBase<double> >(diagOp)->getDiag();
982 return Teuchos::rcp_const_cast<Thyra::MultiVectorBase<double> >(v);
983}
984
996const ModifiableLinearOp getInvDiagonalOp(const LinearOp & op)
997{
998 // if this is a diagonal linear op already, just take the reciprocal
999 auto diagonal_op = rcp_dynamic_cast<const Thyra::DiagonalLinearOpBase<double>>(op);
1000 if(diagonal_op != Teuchos::null){
1001 auto diag = diagonal_op->getDiag();
1002 auto inv_diag = diag->clone_v();
1003 Thyra::reciprocal(*diag,inv_diag.ptr());
1004 return rcp(new Thyra::DefaultDiagonalLinearOp<double>(inv_diag));
1005 }
1006
1007 // if this is a blocked operator, extract diagonals block by block
1008 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_op = rcp_dynamic_cast<const Thyra::PhysicallyBlockedLinearOpBase<double> >(op);
1009 if(blocked_op != Teuchos::null){
1010 int numRows = blocked_op->productRange()->numBlocks();
1011 TEUCHOS_ASSERT(blocked_op->productDomain()->numBlocks() == numRows);
1012 RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_diag = Thyra::defaultBlockedLinearOp<double>();
1013 blocked_diag->beginBlockFill(numRows,numRows);
1014 for(int r = 0; r < numRows; ++r){
1015 for(int c = 0; c < numRows; ++c){
1016 if(r==c)
1017 blocked_diag->setNonconstBlock(r,c,getInvDiagonalOp(blocked_op->getBlock(r,c)));
1018 else
1019 blocked_diag->setBlock(r,c,Thyra::zero<double>(blocked_op->getBlock(r,c)->range(),blocked_op->getBlock(r,c)->domain()));
1020 }
1021 }
1022 blocked_diag->endBlockFill();
1023 return blocked_diag;
1024 }
1025
1026 if (Teko::TpetraHelpers::isTpetraLinearOp(op)){
1027 ST scalar = 0.0;
1028 bool transp = false;
1029 RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOp = Teko::TpetraHelpers::getTpetraCrsMatrix(op, &scalar, &transp);
1030
1031 // extract diagonal
1032 const RCP<Tpetra::Vector<ST,LO,GO,NT> > diag = Tpetra::createVector<ST,LO,GO,NT>(tCrsOp->getRowMap());
1033 diag->scale(scalar);
1034 tCrsOp->getLocalDiagCopy(*diag);
1035 diag->reciprocal(*diag);
1036
1037 // build Thyra diagonal operator
1038 return Teko::TpetraHelpers::thyraDiagOp(diag,*tCrsOp->getRowMap(),"inv(diag( " + op->getObjectLabel() + " ))");
1039
1040 }
1041 else {
1042 RCP<const Thyra::EpetraLinearOp > eOp = rcp_dynamic_cast<const Thyra::EpetraLinearOp >(op,true);
1043 RCP<const Epetra_CrsMatrix> eCrsOp = rcp_dynamic_cast<const Epetra_CrsMatrix>(eOp->epetra_op(),true);
1044
1045 // extract diagonal
1046 const RCP<Epetra_Vector> diag = rcp(new Epetra_Vector(eCrsOp->RowMap()));
1047 TEUCHOS_TEST_FOR_EXCEPT(eCrsOp->ExtractDiagonalCopy(*diag));
1048 diag->Reciprocal(*diag);
1049
1050 // build Thyra diagonal operator
1051 return Teko::Epetra::thyraDiagOp(diag,eCrsOp->RowMap(),"inv(diag( " + op->getObjectLabel() + " ))");
1052 }
1053}
1054
1067const LinearOp explicitMultiply(const LinearOp & opl,const LinearOp & opm,const LinearOp & opr)
1068{
1069 // if this is a blocked operator, multiply block by block
1070 // it is possible that not every factor in the product is blocked and these situations are handled separately
1071
1072 bool isBlockedL = isPhysicallyBlockedLinearOp(opl);
1073 bool isBlockedM = isPhysicallyBlockedLinearOp(opm);
1074 bool isBlockedR = isPhysicallyBlockedLinearOp(opr);
1075
1076 // all factors blocked
1077 if((isBlockedL && isBlockedM && isBlockedR)){
1078
1079 double scalarl = 0.0;
1080 bool transpl = false;
1081 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opl = getPhysicallyBlockedLinearOp(opl,&scalarl,&transpl);
1082 double scalarm = 0.0;
1083 bool transpm = false;
1084 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opm = getPhysicallyBlockedLinearOp(opm,&scalarm,&transpm);
1085 double scalarr = 0.0;
1086 bool transpr = false;
1087 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opr = getPhysicallyBlockedLinearOp(opr,&scalarr,&transpr);
1088 double scalar = scalarl*scalarm*scalarr;
1089
1090 int numRows = blocked_opl->productRange()->numBlocks();
1091 int numCols = blocked_opr->productDomain()->numBlocks();
1092 int numMiddle = blocked_opm->productRange()->numBlocks();
1093
1094 // Assume that the middle block is block nxn and that it's diagonal. Otherwise use the two argument explicitMultiply twice
1095 TEUCHOS_ASSERT(blocked_opm->productDomain()->numBlocks() == numMiddle);
1096 TEUCHOS_ASSERT(blocked_opl->productDomain()->numBlocks() == numMiddle);
1097 TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == numMiddle);
1098
1099 RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_product = Thyra::defaultBlockedLinearOp<double>();
1100 blocked_product->beginBlockFill(numRows,numCols);
1101 for(int r = 0; r < numRows; ++r){
1102 for(int c = 0; c < numCols; ++c){
1103 LinearOp product_rc = explicitMultiply(blocked_opl->getBlock(r,0),blocked_opm->getBlock(0,0),blocked_opr->getBlock(0,c));
1104 for(int m = 1; m < numMiddle; ++m){
1105 LinearOp product_m = explicitMultiply(blocked_opl->getBlock(r,m),blocked_opm->getBlock(m,m),blocked_opr->getBlock(m,c));
1106 product_rc = explicitAdd(product_rc,product_m);
1107 }
1108 blocked_product->setBlock(r,c,product_rc);
1109 }
1110 }
1111 blocked_product->endBlockFill();
1112 return Thyra::scale<double>(scalar,blocked_product.getConst());
1113 }
1114
1115 // left and right factors blocked
1116 if(isBlockedL && !isBlockedM && isBlockedR){
1117 double scalarl = 0.0;
1118 bool transpl = false;
1119 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opl = getPhysicallyBlockedLinearOp(opl,&scalarl,&transpl);
1120 double scalarr = 0.0;
1121 bool transpr = false;
1122 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opr = getPhysicallyBlockedLinearOp(opr,&scalarr,&transpr);
1123 double scalar = scalarl*scalarr;
1124
1125 int numRows = blocked_opl->productRange()->numBlocks();
1126 int numCols = blocked_opr->productDomain()->numBlocks();
1127 int numMiddle = 1;
1128
1129 // Assume that the middle block is 1x1 diagonal. Left must be rx1, right 1xc
1130 TEUCHOS_ASSERT(blocked_opl->productDomain()->numBlocks() == numMiddle);
1131 TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == numMiddle);
1132
1133 RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_product = Thyra::defaultBlockedLinearOp<double>();
1134 blocked_product->beginBlockFill(numRows,numCols);
1135 for(int r = 0; r < numRows; ++r){
1136 for(int c = 0; c < numCols; ++c){
1137 LinearOp product_rc = explicitMultiply(blocked_opl->getBlock(r,0),opm,blocked_opr->getBlock(0,c));
1138 blocked_product->setBlock(r,c,product_rc);
1139 }
1140 }
1141 blocked_product->endBlockFill();
1142 return Thyra::scale<double>(scalar,blocked_product.getConst());
1143 }
1144
1145 // only right factor blocked
1146 if(!isBlockedL && !isBlockedM && isBlockedR){
1147 double scalarr = 0.0;
1148 bool transpr = false;
1149 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opr = getPhysicallyBlockedLinearOp(opr,&scalarr,&transpr);
1150 double scalar = scalarr;
1151
1152 int numRows = 1;
1153 int numCols = blocked_opr->productDomain()->numBlocks();
1154 int numMiddle = 1;
1155
1156 // Assume that the middle block is 1x1 diagonal, left is 1x1. Right must be 1xc
1157 TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == numMiddle);
1158
1159 RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_product = Thyra::defaultBlockedLinearOp<double>();
1160 blocked_product->beginBlockFill(numRows,numCols);
1161 for(int c = 0; c < numCols; ++c){
1162 LinearOp product_c = explicitMultiply(opl,opm,blocked_opr->getBlock(0,c));
1163 blocked_product->setBlock(0,c,product_c);
1164 }
1165 blocked_product->endBlockFill();
1166 return Thyra::scale<double>(scalar,blocked_product.getConst());
1167 }
1168
1169 //TODO: three more cases (only non-blocked - blocked - non-blocked not possible)
1170
1171 bool isTpetral = Teko::TpetraHelpers::isTpetraLinearOp(opl);
1172 bool isTpetram = Teko::TpetraHelpers::isTpetraLinearOp(opm);
1173 bool isTpetrar = Teko::TpetraHelpers::isTpetraLinearOp(opr);
1174
1175 if(isTpetral && isTpetram && isTpetrar){ // Both operators are Tpetra matrices so explicitly multiply them
1176
1177 // Get left and right Tpetra crs operators
1178 ST scalarl = 0.0;
1179 bool transpl = false;
1180 RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpl = Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1181 ST scalarm = 0.0;
1182 bool transpm = false;
1183 RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpm = Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarm, &transpm);
1184 ST scalarr = 0.0;
1185 bool transpr = false;
1186 RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpr = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1187
1188 // Build output operator
1189 RCP<Thyra::LinearOpBase<ST> > explicitOp = rcp(new Thyra::TpetraLinearOp<ST,LO,GO,NT>());
1190 RCP<Thyra::TpetraLinearOp<ST,LO,GO,NT> > tExplicitOp = rcp_dynamic_cast<Thyra::TpetraLinearOp<ST,LO,GO,NT> >(explicitOp);
1191
1192 // Do explicit matrix-matrix multiply
1193 RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOplm = Tpetra::createCrsMatrix<ST,LO,GO,NT>(tCrsOpl->getRowMap());
1194 RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > explicitCrsOp = Tpetra::createCrsMatrix<ST,LO,GO,NT>(tCrsOpl->getRowMap());
1195 Tpetra::MatrixMatrix::Multiply<ST,LO,GO,NT>(*tCrsOpl,transpl,*tCrsOpm,transpm,*tCrsOplm);
1196 Tpetra::MatrixMatrix::Multiply<ST,LO,GO,NT>(*tCrsOplm,false,*tCrsOpr,transpr,*explicitCrsOp);
1197 explicitCrsOp->resumeFill();
1198 explicitCrsOp->scale(scalarl*scalarm*scalarr);
1199 explicitCrsOp->fillComplete(tCrsOpr->getDomainMap(),tCrsOpl->getRangeMap());
1200 tExplicitOp->initialize(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1201 return tExplicitOp;
1202
1203 } else if (isTpetral && !isTpetram && isTpetrar){ // Assume that the middle operator is diagonal
1204
1205 // Get left and right Tpetra crs operators
1206 ST scalarl = 0.0;
1207 bool transpl = false;
1208 RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpl = Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1209 ST scalarr = 0.0;
1210 bool transpr = false;
1211 RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpr = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1212
1213 RCP<const Tpetra::Vector<ST,LO,GO,NT> > diagPtr;
1214
1215 // Cast middle operator as DiagonalLinearOp and extract diagonal as Vector
1216 RCP<const Thyra::DiagonalLinearOpBase<ST> > dOpm = rcp_dynamic_cast<const Thyra::DiagonalLinearOpBase<ST> >(opm);
1217 if(dOpm != Teuchos::null){
1218 RCP<const Thyra::TpetraVector<ST,LO,GO,NT> > tPtr = rcp_dynamic_cast<const Thyra::TpetraVector<ST,LO,GO,NT> >(dOpm->getDiag(),true);
1219 diagPtr = rcp_dynamic_cast<const Tpetra::Vector<ST,LO,GO,NT> >(tPtr->getConstTpetraVector(),true);
1220 }
1221 // If it's not diagonal, maybe it's zero
1222 else if(rcp_dynamic_cast<const Thyra::ZeroLinearOpBase<ST> >(opm) != Teuchos::null){
1223 diagPtr = rcp(new Tpetra::Vector<ST,LO,GO,NT>(tCrsOpl->getDomainMap()));
1224 }
1225 else
1226 TEUCHOS_ASSERT(false);
1227
1228 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()));
1229
1230 // Do the diagonal scaling
1231 tCrsOplm->rightScale(*diagPtr);
1232
1233 // Build output operator
1234 RCP<Thyra::LinearOpBase<ST> > explicitOp = rcp(new Thyra::TpetraLinearOp<ST,LO,GO,NT>());
1235 RCP<Thyra::TpetraLinearOp<ST,LO,GO,NT> > tExplicitOp = rcp_dynamic_cast<Thyra::TpetraLinearOp<ST,LO,GO,NT> >(explicitOp);
1236
1237 // Do explicit matrix-matrix multiply
1238 RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > explicitCrsOp = Tpetra::createCrsMatrix<ST,LO,GO,NT>(tCrsOpl->getRowMap());
1239 Tpetra::MatrixMatrix::Multiply<ST,LO,GO,NT>(*tCrsOplm,false,*tCrsOpr,transpr,*explicitCrsOp);
1240 explicitCrsOp->resumeFill();
1241 explicitCrsOp->scale(scalarl*scalarr);
1242 explicitCrsOp->fillComplete(tCrsOpr->getDomainMap(),tCrsOpl->getRangeMap());
1243 tExplicitOp->initialize(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1244 return tExplicitOp;
1245
1246 } else { // Assume Epetra and we can use transformers
1247
1248 // build implicit multiply
1249 const LinearOp implicitOp = Thyra::multiply(opl,opm,opr);
1250
1251 // build transformer
1252 const RCP<Thyra::LinearOpTransformerBase<double> > prodTrans =
1253 Thyra::epetraExtDiagScaledMatProdTransformer();
1254
1255 // build operator and multiply
1256 const RCP<Thyra::LinearOpBase<double> > explicitOp = prodTrans->createOutputOp();
1257 prodTrans->transform(*implicitOp,explicitOp.ptr());
1258 explicitOp->setObjectLabel("explicit( " + opl->getObjectLabel() +
1259 " * " + opm->getObjectLabel() +
1260 " * " + opr->getObjectLabel() + " )");
1261
1262 return explicitOp;
1263
1264 }
1265}
1266
1281const ModifiableLinearOp explicitMultiply(const LinearOp & opl,const LinearOp & opm,const LinearOp & opr,
1282 const ModifiableLinearOp & destOp)
1283{
1284 bool isTpetral = Teko::TpetraHelpers::isTpetraLinearOp(opl);
1285 bool isTpetram = Teko::TpetraHelpers::isTpetraLinearOp(opm);
1286 bool isTpetrar = Teko::TpetraHelpers::isTpetraLinearOp(opr);
1287
1288 if(isTpetral && isTpetram && isTpetrar){ // Both operators are Tpetra matrices so explicitly multiply them
1289
1290 // Get left and right Tpetra crs operators
1291 ST scalarl = 0.0;
1292 bool transpl = false;
1293 RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpl = Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1294 ST scalarm = 0.0;
1295 bool transpm = false;
1296 RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpm = Teko::TpetraHelpers::getTpetraCrsMatrix(opm, &scalarm, &transpm);
1297 ST scalarr = 0.0;
1298 bool transpr = false;
1299 RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpr = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1300
1301 // Build output operator
1302 auto tExplicitOp = rcp_dynamic_cast<Thyra::TpetraLinearOp<ST,LO,GO,NT> >(destOp);
1303 if(tExplicitOp.is_null())
1304 tExplicitOp = rcp(new Thyra::TpetraLinearOp<ST,LO,GO,NT>());
1305
1306 // Do explicit matrix-matrix multiply
1307 RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOplm = Tpetra::createCrsMatrix<ST,LO,GO,NT>(tCrsOpl->getRowMap());
1308 RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > explicitCrsOp = Tpetra::createCrsMatrix<ST,LO,GO,NT>(tCrsOpl->getRowMap());
1309 Tpetra::MatrixMatrix::Multiply<ST,LO,GO,NT>(*tCrsOpl,transpl,*tCrsOpm,transpm,*tCrsOplm);
1310 Tpetra::MatrixMatrix::Multiply<ST,LO,GO,NT>(*tCrsOplm,false,*tCrsOpr,transpr,*explicitCrsOp);
1311 explicitCrsOp->resumeFill();
1312 explicitCrsOp->scale(scalarl*scalarm*scalarr);
1313 explicitCrsOp->fillComplete(tCrsOpr->getDomainMap(),tCrsOpl->getRangeMap());
1314 tExplicitOp->initialize(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1315 return tExplicitOp;
1316
1317 } else if (isTpetral && !isTpetram && isTpetrar){ // Assume that the middle operator is diagonal
1318
1319 // Get left and right Tpetra crs operators
1320 ST scalarl = 0.0;
1321 bool transpl = false;
1322 RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpl = Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1323 ST scalarr = 0.0;
1324 bool transpr = false;
1325 RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpr = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1326
1327 // Cast middle operator as DiagonalLinearOp and extract diagonal as Vector
1328 RCP<const Thyra::DiagonalLinearOpBase<ST> > dOpm = rcp_dynamic_cast<const Thyra::DiagonalLinearOpBase<ST> >(opm,true);
1329 RCP<const Thyra::TpetraVector<ST,LO,GO,NT> > tPtr = rcp_dynamic_cast<const Thyra::TpetraVector<ST,LO,GO,NT> >(dOpm->getDiag(),true);
1330 RCP<const Tpetra::Vector<ST,LO,GO,NT> > diagPtr = rcp_dynamic_cast<const Tpetra::Vector<ST,LO,GO,NT> >(tPtr->getConstTpetraVector(),true);
1331 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()));
1332
1333 // Do the diagonal scaling
1334 tCrsOplm->rightScale(*diagPtr);
1335
1336 // Build output operator
1337 RCP<Thyra::LinearOpBase<ST> > explicitOp = rcp(new Thyra::TpetraLinearOp<ST,LO,GO,NT>());
1338 RCP<Thyra::TpetraLinearOp<ST,LO,GO,NT> > tExplicitOp = rcp_dynamic_cast<Thyra::TpetraLinearOp<ST,LO,GO,NT> >(explicitOp);
1339
1340 // Do explicit matrix-matrix multiply
1341 RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > explicitCrsOp = Tpetra::createCrsMatrix<ST,LO,GO,NT>(tCrsOpl->getRowMap());
1342 Tpetra::MatrixMatrix::Multiply<ST,LO,GO,NT>(*tCrsOplm,false,*tCrsOpr,transpr,*explicitCrsOp);
1343 explicitCrsOp->resumeFill();
1344 explicitCrsOp->scale(scalarl*scalarr);
1345 explicitCrsOp->fillComplete(tCrsOpr->getDomainMap(),tCrsOpl->getRangeMap());
1346 tExplicitOp->initialize(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1347 return tExplicitOp;
1348
1349 } else { // Assume Epetra and we can use transformers
1350
1351 // build implicit multiply
1352 const LinearOp implicitOp = Thyra::multiply(opl,opm,opr);
1353
1354 // build transformer
1355 const RCP<Thyra::LinearOpTransformerBase<double> > prodTrans =
1356 Thyra::epetraExtDiagScaledMatProdTransformer();
1357
1358 // build operator destination operator
1359 ModifiableLinearOp explicitOp;
1360
1361 // if neccessary build a operator to put the explicit multiply into
1362 if(destOp==Teuchos::null)
1363 explicitOp = prodTrans->createOutputOp();
1364 else
1365 explicitOp = destOp;
1366
1367 // perform multiplication
1368 prodTrans->transform(*implicitOp,explicitOp.ptr());
1369
1370 // label it
1371 explicitOp->setObjectLabel("explicit( " + opl->getObjectLabel() +
1372 " * " + opm->getObjectLabel() +
1373 " * " + opr->getObjectLabel() + " )");
1374
1375 return explicitOp;
1376
1377 }
1378}
1379
1390const LinearOp explicitMultiply(const LinearOp & opl,const LinearOp & opr)
1391{
1392 // if this is a blocked operator, multiply block by block
1393 // it is possible that not every factor in the product is blocked and these situations are handled separately
1394
1395 bool isBlockedL = isPhysicallyBlockedLinearOp(opl);
1396 bool isBlockedR = isPhysicallyBlockedLinearOp(opr);
1397
1398 // both factors blocked
1399 if((isBlockedL && isBlockedR)){
1400 double scalarl = 0.0;
1401 bool transpl = false;
1402 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opl = getPhysicallyBlockedLinearOp(opl,&scalarl,&transpl);
1403 double scalarr = 0.0;
1404 bool transpr = false;
1405 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opr = getPhysicallyBlockedLinearOp(opr,&scalarr,&transpr);
1406 double scalar = scalarl*scalarr;
1407
1408 int numRows = blocked_opl->productRange()->numBlocks();
1409 int numCols = blocked_opr->productDomain()->numBlocks();
1410 int numMiddle = blocked_opl->productDomain()->numBlocks();
1411
1412 TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == numMiddle);
1413
1414 RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_product = Thyra::defaultBlockedLinearOp<double>();
1415 blocked_product->beginBlockFill(numRows,numCols);
1416 for(int r = 0; r < numRows; ++r){
1417 for(int c = 0; c < numCols; ++c){
1418 LinearOp product_rc = explicitMultiply(blocked_opl->getBlock(r,0),blocked_opr->getBlock(0,c));
1419 for(int m = 1; m < numMiddle; ++m){
1420 LinearOp product_m = explicitMultiply(blocked_opl->getBlock(r,m),blocked_opr->getBlock(m,c));
1421 product_rc = explicitAdd(product_rc,product_m);
1422 }
1423 blocked_product->setBlock(r,c,Thyra::scale(scalar,product_rc));
1424 }
1425 }
1426 blocked_product->endBlockFill();
1427 return blocked_product;
1428 }
1429
1430 // only left factor blocked
1431 if((isBlockedL && !isBlockedR)){
1432 double scalarl = 0.0;
1433 bool transpl = false;
1434 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opl = getPhysicallyBlockedLinearOp(opl,&scalarl,&transpl);
1435 double scalar = scalarl;
1436
1437 int numRows = blocked_opl->productRange()->numBlocks();
1438 int numCols = 1;
1439 int numMiddle = 1;
1440
1441 TEUCHOS_ASSERT(blocked_opl->productDomain()->numBlocks() == numMiddle);
1442
1443 RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_product = Thyra::defaultBlockedLinearOp<double>();
1444 blocked_product->beginBlockFill(numRows,numCols);
1445 for(int r = 0; r < numRows; ++r){
1446 LinearOp product_r = explicitMultiply(blocked_opl->getBlock(r,0),opr);
1447 blocked_product->setBlock(r,0,Thyra::scale(scalar,product_r));
1448 }
1449 blocked_product->endBlockFill();
1450 return blocked_product;
1451 }
1452
1453 // only right factor blocked
1454 if((!isBlockedL && isBlockedR)){
1455 double scalarr = 0.0;
1456 bool transpr = false;
1457 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opr = getPhysicallyBlockedLinearOp(opr,&scalarr,&transpr);
1458 double scalar = scalarr;
1459
1460 int numRows = 1;
1461 int numCols = blocked_opr->productDomain()->numBlocks();
1462 int numMiddle = 1;
1463
1464 TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == numMiddle);
1465
1466 RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_product = Thyra::defaultBlockedLinearOp<double>();
1467 blocked_product->beginBlockFill(numRows,numCols);
1468 for(int c = 0; c < numCols; ++c){
1469 LinearOp product_c = explicitMultiply(opl,blocked_opr->getBlock(0,c));
1470 blocked_product->setBlock(0,c,Thyra::scale(scalar,product_c));
1471 }
1472 blocked_product->endBlockFill();
1473 return blocked_product;
1474 }
1475
1476 bool isTpetral = Teko::TpetraHelpers::isTpetraLinearOp(opl);
1477 bool isTpetrar = Teko::TpetraHelpers::isTpetraLinearOp(opr);
1478
1479 if(isTpetral && isTpetrar){ // Both operators are Tpetra matrices so explicitly multiply them
1480 // Get left and right Tpetra crs operators
1481 ST scalarl = 0.0;
1482 bool transpl = false;
1483 RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpl = Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1484 ST scalarr = 0.0;
1485 bool transpr = false;
1486 RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpr = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1487
1488 // Build output operator
1489 RCP<Thyra::LinearOpBase<ST> > explicitOp = rcp(new Thyra::TpetraLinearOp<ST,LO,GO,NT>());
1490 RCP<Thyra::TpetraLinearOp<ST,LO,GO,NT> > tExplicitOp = rcp_dynamic_cast<Thyra::TpetraLinearOp<ST,LO,GO,NT> >(explicitOp);
1491
1492 // Do explicit matrix-matrix multiply
1493 RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > explicitCrsOp = Tpetra::createCrsMatrix<ST,LO,GO,NT>(tCrsOpl->getRowMap());
1494 Tpetra::MatrixMatrix::Multiply<ST,LO,GO,NT>(*tCrsOpl,transpl,*tCrsOpr,transpr,*explicitCrsOp);
1495 explicitCrsOp->resumeFill();
1496 explicitCrsOp->scale(scalarl*scalarr);
1497 explicitCrsOp->fillComplete(tCrsOpr->getDomainMap(),tCrsOpl->getRangeMap());
1498 tExplicitOp->initialize(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1499 return tExplicitOp;
1500
1501 } else if (isTpetral && !isTpetrar){ // Assume that the right operator is diagonal
1502
1503 // Get left Tpetra crs operator
1504 ST scalarl = 0.0;
1505 bool transpl = false;
1506 RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpl = Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1507
1508 // Cast right operator as DiagonalLinearOp and extract diagonal as Vector
1509 RCP<const Thyra::DiagonalLinearOpBase<ST> > dOpr = rcp_dynamic_cast<const Thyra::DiagonalLinearOpBase<ST> >(opr,true);
1510 RCP<const Thyra::TpetraVector<ST,LO,GO,NT> > tPtr = rcp_dynamic_cast<const Thyra::TpetraVector<ST,LO,GO,NT> >(dOpr->getDiag(),true);
1511 RCP<const Tpetra::Vector<ST,LO,GO,NT> > diagPtr = rcp_dynamic_cast<const Tpetra::Vector<ST,LO,GO,NT> >(tPtr->getConstTpetraVector(),true);
1512 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()));
1513
1514 explicitCrsOp->rightScale(*diagPtr);
1515 explicitCrsOp->resumeFill();
1516 explicitCrsOp->scale(scalarl);
1517 explicitCrsOp->fillComplete(tCrsOpl->getDomainMap(),tCrsOpl->getRangeMap());
1518
1519 return Thyra::constTpetraLinearOp<ST,LO,GO,NT>(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1520
1521 } else if (!isTpetral && isTpetrar){ // Assume that the left operator is diagonal
1522
1523 // Get right Tpetra crs operator
1524 ST scalarr = 0.0;
1525 bool transpr = false;
1526 RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpr = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1527
1528 RCP<const Tpetra::Vector<ST,LO,GO,NT> > diagPtr;
1529
1530 // Cast left operator as DiagonalLinearOp and extract diagonal as Vector
1531 RCP<const Thyra::DiagonalLinearOpBase<ST> > dOpl = rcp_dynamic_cast<const Thyra::DiagonalLinearOpBase<ST> >(opl);
1532 if(dOpl != Teuchos::null){
1533 RCP<const Thyra::TpetraVector<ST,LO,GO,NT> > tPtr = rcp_dynamic_cast<const Thyra::TpetraVector<ST,LO,GO,NT> >(dOpl->getDiag(),true);
1534 diagPtr = rcp_dynamic_cast<const Tpetra::Vector<ST,LO,GO,NT> >(tPtr->getConstTpetraVector(),true);
1535 }
1536 // If it's not diagonal, maybe it's zero
1537 else if(rcp_dynamic_cast<const Thyra::ZeroLinearOpBase<ST> >(opl) != Teuchos::null){
1538 diagPtr = rcp(new Tpetra::Vector<ST,LO,GO,NT>(tCrsOpr->getRangeMap()));
1539 }
1540 else
1541 TEUCHOS_ASSERT(false);
1542
1543 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()));
1544
1545 explicitCrsOp->leftScale(*diagPtr);
1546 explicitCrsOp->resumeFill();
1547 explicitCrsOp->scale(scalarr);
1548 explicitCrsOp->fillComplete(tCrsOpr->getDomainMap(),tCrsOpr->getRangeMap());
1549
1550 return Thyra::constTpetraLinearOp<ST,LO,GO,NT>(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1551
1552 } else { // Assume Epetra and we can use transformers
1553
1554 // build implicit multiply
1555 const LinearOp implicitOp = Thyra::multiply(opl,opr);
1556
1557 // build a scaling transformer
1558 RCP<Thyra::LinearOpTransformerBase<double> > prodTrans
1559 = Thyra::epetraExtDiagScalingTransformer();
1560
1561 // check to see if a scaling transformer works: if not use the
1562 // DiagScaledMatrixProduct transformer
1563 if(not prodTrans->isCompatible(*implicitOp))
1564 prodTrans = Thyra::epetraExtDiagScaledMatProdTransformer();
1565
1566 // build operator and multiply
1567 const RCP<Thyra::LinearOpBase<double> > explicitOp = prodTrans->createOutputOp();
1568 prodTrans->transform(*implicitOp,explicitOp.ptr());
1569 explicitOp->setObjectLabel("explicit( " + opl->getObjectLabel() +
1570 " * " + opr->getObjectLabel() + " )");
1571
1572 return explicitOp;
1573 }
1574}
1575
1589const ModifiableLinearOp explicitMultiply(const LinearOp & opl,const LinearOp & opr,
1590 const ModifiableLinearOp & destOp)
1591{
1592 // if this is a blocked operator, multiply block by block
1593 // it is possible that not every factor in the product is blocked and these situations are handled separately
1594
1595 bool isBlockedL = isPhysicallyBlockedLinearOp(opl);
1596 bool isBlockedR = isPhysicallyBlockedLinearOp(opr);
1597
1598 // both factors blocked
1599 if((isBlockedL && isBlockedR)){
1600 double scalarl = 0.0;
1601 bool transpl = false;
1602 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opl = getPhysicallyBlockedLinearOp(opl,&scalarl,&transpl);
1603 double scalarr = 0.0;
1604 bool transpr = false;
1605 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opr = getPhysicallyBlockedLinearOp(opr,&scalarr,&transpr);
1606 double scalar = scalarl*scalarr;
1607
1608 int numRows = blocked_opl->productRange()->numBlocks();
1609 int numCols = blocked_opr->productDomain()->numBlocks();
1610 int numMiddle = blocked_opl->productDomain()->numBlocks();
1611
1612 TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == numMiddle);
1613
1614 RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_product = Thyra::defaultBlockedLinearOp<double>();
1615 blocked_product->beginBlockFill(numRows,numCols);
1616 for(int r = 0; r < numRows; ++r){
1617 for(int c = 0; c < numCols; ++c){
1618
1619 LinearOp product_rc = explicitMultiply(blocked_opl->getBlock(r,0),blocked_opr->getBlock(0,c));
1620 for(int m = 1; m < numMiddle; ++m){
1621 LinearOp product_m = explicitMultiply(blocked_opl->getBlock(r,m),blocked_opr->getBlock(m,c));
1622 product_rc = explicitAdd(product_rc,product_m);
1623 }
1624 blocked_product->setBlock(r,c,Thyra::scale(scalar,product_rc));
1625 }
1626 }
1627 blocked_product->endBlockFill();
1628 return blocked_product;
1629 }
1630
1631 // only left factor blocked
1632 if((isBlockedL && !isBlockedR)){
1633 double scalarl = 0.0;
1634 bool transpl = false;
1635 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opl = getPhysicallyBlockedLinearOp(opl,&scalarl,&transpl);
1636 double scalar = scalarl;
1637
1638 int numRows = blocked_opl->productRange()->numBlocks();
1639 int numCols = 1;
1640 int numMiddle = 1;
1641
1642 TEUCHOS_ASSERT(blocked_opl->productDomain()->numBlocks() == numMiddle);
1643
1644 RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_product = Thyra::defaultBlockedLinearOp<double>();
1645 blocked_product->beginBlockFill(numRows,numCols);
1646 for(int r = 0; r < numRows; ++r){
1647 LinearOp product_r = explicitMultiply(blocked_opl->getBlock(r,0),opr);
1648 blocked_product->setBlock(r,0,Thyra::scale(scalar,product_r));
1649 }
1650 blocked_product->endBlockFill();
1651 return blocked_product;
1652 }
1653
1654 // only right factor blocked
1655 if((!isBlockedL && isBlockedR)){
1656 double scalarr = 0.0;
1657 bool transpr = false;
1658 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opr = getPhysicallyBlockedLinearOp(opr,&scalarr,&transpr);
1659 double scalar = scalarr;
1660
1661 int numRows = 1;
1662 int numCols = blocked_opr->productDomain()->numBlocks();
1663 int numMiddle = 1;
1664
1665 TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == numMiddle);
1666
1667 RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_product = Thyra::defaultBlockedLinearOp<double>();
1668 blocked_product->beginBlockFill(numRows,numCols);
1669 for(int c = 0; c < numCols; ++c){
1670 LinearOp product_c = explicitMultiply(opl,blocked_opr->getBlock(0,c));
1671 blocked_product->setBlock(0,c,Thyra::scale(scalar,product_c));
1672 }
1673 blocked_product->endBlockFill();
1674 return blocked_product;
1675 }
1676
1677 bool isTpetral = Teko::TpetraHelpers::isTpetraLinearOp(opl);
1678 bool isTpetrar = Teko::TpetraHelpers::isTpetraLinearOp(opr);
1679
1680 if(isTpetral && isTpetrar){ // Both operators are Tpetra matrices so use the explicit Tpetra matrix-matrix multiply
1681
1682 // Get left and right Tpetra crs operators
1683 ST scalarl = 0.0;
1684 bool transpl = false;
1685 RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpl = Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1686 ST scalarr = 0.0;
1687 bool transpr = false;
1688 RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpr = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1689
1690 // Build output operator
1691 RCP<Thyra::LinearOpBase<ST> > explicitOp;
1692 if(destOp!=Teuchos::null)
1693 explicitOp = destOp;
1694 else
1695 explicitOp = rcp(new Thyra::TpetraLinearOp<ST,LO,GO,NT>());
1696 RCP<Thyra::TpetraLinearOp<ST,LO,GO,NT> > tExplicitOp = rcp_dynamic_cast<Thyra::TpetraLinearOp<ST,LO,GO,NT> >(explicitOp);
1697
1698 // Do explicit matrix-matrix multiply
1699 RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > explicitCrsOp = Tpetra::createCrsMatrix<ST,LO,GO,NT>(tCrsOpl->getRowMap());
1700 Tpetra::MatrixMatrix::Multiply<ST,LO,GO,NT>(*tCrsOpl,transpl,*tCrsOpr,transpr,*explicitCrsOp);
1701 explicitCrsOp->resumeFill();
1702 explicitCrsOp->scale(scalarl*scalarr);
1703 explicitCrsOp->fillComplete(tCrsOpr->getDomainMap(),tCrsOpl->getRangeMap());
1704 tExplicitOp->initialize(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1705 return tExplicitOp;
1706
1707 } else if (isTpetral && !isTpetrar){ // Assume that the right operator is diagonal
1708
1709 // Get left Tpetra crs operator
1710 ST scalarl = 0.0;
1711 bool transpl = false;
1712 RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpl = Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1713
1714 // Cast right operator as DiagonalLinearOp and extract diagonal as Vector
1715 RCP<const Thyra::DiagonalLinearOpBase<ST> > dOpr = rcp_dynamic_cast<const Thyra::DiagonalLinearOpBase<ST> >(opr);
1716 RCP<const Thyra::TpetraVector<ST,LO,GO,NT> > tPtr = rcp_dynamic_cast<const Thyra::TpetraVector<ST,LO,GO,NT> >(dOpr->getDiag(),true);
1717 RCP<const Tpetra::Vector<ST,LO,GO,NT> > diagPtr = rcp_dynamic_cast<const Tpetra::Vector<ST,LO,GO,NT> >(tPtr->getConstTpetraVector(),true);
1718
1719 // Scale by the diagonal operator
1720 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()));
1721 explicitCrsOp->rightScale(*diagPtr);
1722 explicitCrsOp->resumeFill();
1723 explicitCrsOp->scale(scalarl);
1724 explicitCrsOp->fillComplete(tCrsOpl->getDomainMap(),tCrsOpl->getRangeMap());
1725 return Thyra::tpetraLinearOp<ST,LO,GO,NT>(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1726
1727 } else if (!isTpetral && isTpetrar){ // Assume that the left operator is diagonal
1728
1729 // Get right Tpetra crs operator
1730 ST scalarr = 0.0;
1731 bool transpr = false;
1732 RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpr = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1733
1734 // Cast leftt operator as DiagonalLinearOp and extract diagonal as Vector
1735 RCP<const Thyra::DiagonalLinearOpBase<ST> > dOpl = rcp_dynamic_cast<const Thyra::DiagonalLinearOpBase<ST> >(opl,true);
1736 RCP<const Thyra::TpetraVector<ST,LO,GO,NT> > tPtr = rcp_dynamic_cast<const Thyra::TpetraVector<ST,LO,GO,NT> >(dOpl->getDiag(),true);
1737 RCP<const Tpetra::Vector<ST,LO,GO,NT> > diagPtr = rcp_dynamic_cast<const Tpetra::Vector<ST,LO,GO,NT> >(tPtr->getConstTpetraVector(),true);
1738
1739 // Scale by the diagonal operator
1740 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()));
1741 explicitCrsOp->leftScale(*diagPtr);
1742 explicitCrsOp->resumeFill();
1743 explicitCrsOp->scale(scalarr);
1744 explicitCrsOp->fillComplete(tCrsOpr->getDomainMap(),tCrsOpr->getRangeMap());
1745 return Thyra::tpetraLinearOp<ST,LO,GO,NT>(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1746
1747 } else { // Assume Epetra and we can use transformers
1748
1749 // build implicit multiply
1750 const LinearOp implicitOp = Thyra::multiply(opl,opr);
1751
1752 // build a scaling transformer
1753
1754 RCP<Thyra::LinearOpTransformerBase<double> > prodTrans
1755 = Thyra::epetraExtDiagScalingTransformer();
1756
1757 // check to see if a scaling transformer works: if not use the
1758 // DiagScaledMatrixProduct transformer
1759 if(not prodTrans->isCompatible(*implicitOp))
1760 prodTrans = Thyra::epetraExtDiagScaledMatProdTransformer();
1761
1762 // build operator destination operator
1763 ModifiableLinearOp explicitOp;
1764
1765 // if neccessary build a operator to put the explicit multiply into
1766 if(destOp==Teuchos::null)
1767 explicitOp = prodTrans->createOutputOp();
1768 else
1769 explicitOp = destOp;
1770
1771 // perform multiplication
1772 prodTrans->transform(*implicitOp,explicitOp.ptr());
1773
1774 // label it
1775 explicitOp->setObjectLabel("explicit( " + opl->getObjectLabel() +
1776 " * " + opr->getObjectLabel() + " )");
1777
1778 return explicitOp;
1779 }
1780}
1781
1792const LinearOp explicitAdd(const LinearOp & opl_in,const LinearOp & opr_in)
1793{
1794 // if both blocked, add block by block
1795 if(isPhysicallyBlockedLinearOp(opl_in) && isPhysicallyBlockedLinearOp(opr_in)){
1796 double scalarl = 0.0;
1797 bool transpl = false;
1798 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opl = getPhysicallyBlockedLinearOp(opl_in, &scalarl, &transpl);
1799
1800 double scalarr = 0.0;
1801 bool transpr = false;
1802 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opr = getPhysicallyBlockedLinearOp(opr_in, &scalarr, &transpr);
1803
1804 int numRows = blocked_opl->productRange()->numBlocks();
1805 int numCols = blocked_opl->productDomain()->numBlocks();
1806 TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == numRows);
1807 TEUCHOS_ASSERT(blocked_opr->productDomain()->numBlocks() == numCols);
1808
1809 RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_sum = Thyra::defaultBlockedLinearOp<double>();
1810 blocked_sum->beginBlockFill(numRows,numCols);
1811 for(int r = 0; r < numRows; ++r)
1812 for(int c = 0; c < numCols; ++c)
1813 blocked_sum->setBlock(r,c,explicitAdd(Thyra::scale(scalarl,blocked_opl->getBlock(r,c)),Thyra::scale(scalarr,blocked_opr->getBlock(r,c))));
1814 blocked_sum->endBlockFill();
1815 return blocked_sum;
1816 }
1817
1818 // if only one is blocked, it must be 1x1
1819 LinearOp opl = opl_in;
1820 LinearOp opr = opr_in;
1821 if(isPhysicallyBlockedLinearOp(opl_in)){
1822 double scalarl = 0.0;
1823 bool transpl = false;
1824 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opl = getPhysicallyBlockedLinearOp(opl_in, &scalarl, &transpl);
1825 TEUCHOS_ASSERT(blocked_opl->productRange()->numBlocks() == 1);
1826 TEUCHOS_ASSERT(blocked_opl->productDomain()->numBlocks() == 1);
1827 opl = Thyra::scale(scalarl,blocked_opl->getBlock(0,0));
1828 }
1829 if(isPhysicallyBlockedLinearOp(opr_in)){
1830 double scalarr = 0.0;
1831 bool transpr = false;
1832 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opr = getPhysicallyBlockedLinearOp(opr_in, &scalarr, &transpr);
1833 TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == 1);
1834 TEUCHOS_ASSERT(blocked_opr->productDomain()->numBlocks() == 1);
1835 opr = Thyra::scale(scalarr,blocked_opr->getBlock(0,0));
1836 }
1837
1838 bool isTpetral = Teko::TpetraHelpers::isTpetraLinearOp(opl);
1839 bool isTpetrar = Teko::TpetraHelpers::isTpetraLinearOp(opr);
1840
1841 // if one of the operators in the sum is a thyra zero op
1842 if(isZeroOp(opl)){
1843 if(isZeroOp(opr))
1844 return opr; // return a zero op if both are zero
1845 if(isTpetrar){ // if other op is tpetra, replace this with a zero crs matrix
1846 ST scalar = 0.0;
1847 bool transp = false;
1848 auto crs_op = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalar, &transp);
1849 auto zero_crs = Tpetra::createCrsMatrix<ST,LO,GO,NT>(crs_op->getRowMap());
1850 zero_crs->fillComplete();
1851 opl = Thyra::constTpetraLinearOp<ST,LO,GO,NT>(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(crs_op->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(crs_op->getDomainMap()),zero_crs);
1852 isTpetral = true;
1853 } else
1854 return opr->clone();
1855 }
1856 if(isZeroOp(opr)){
1857 if(isTpetral){ // if other op is tpetra, replace this with a zero crs matrix
1858 ST scalar = 0.0;
1859 bool transp = false;
1860 auto crs_op = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalar, &transp);
1861 auto zero_crs = Tpetra::createCrsMatrix<ST,LO,GO,NT>(crs_op->getRowMap());
1862 zero_crs->fillComplete();
1863 opr = Thyra::constTpetraLinearOp<ST,LO,GO,NT>(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(crs_op->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(crs_op->getDomainMap()),zero_crs);
1864 isTpetrar = true;
1865 } else
1866 return opl->clone();
1867 }
1868
1869 if(isTpetral && isTpetrar){ // Both operators are Tpetra matrices so use the explicit Tpetra matrix-matrix add
1870
1871 // Get left and right Tpetra crs operators
1872 ST scalarl = 0.0;
1873 bool transpl = false;
1874 RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpl = Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1875 ST scalarr = 0.0;
1876 bool transpr = false;
1877 RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpr = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1878
1879 // Build output operator
1880 RCP<Thyra::LinearOpBase<ST> > explicitOp = rcp(new Thyra::TpetraLinearOp<ST,LO,GO,NT>());
1881 RCP<Thyra::TpetraLinearOp<ST,LO,GO,NT> > tExplicitOp = rcp_dynamic_cast<Thyra::TpetraLinearOp<ST,LO,GO,NT> >(explicitOp);
1882
1883 // Do explicit matrix-matrix add
1884 RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > explicitCrsOp = Tpetra::MatrixMatrix::add<ST,LO,GO,NT>(scalarl,transpl,*tCrsOpl,scalarr,transpr,*tCrsOpr);
1885 tExplicitOp->initialize(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1886 return tExplicitOp;
1887
1888 }else{//Assume Epetra
1889 // build implicit add
1890 const LinearOp implicitOp = Thyra::add(opl,opr);
1891
1892 // build transformer
1893 const RCP<Thyra::LinearOpTransformerBase<double> > prodTrans =
1894 Thyra::epetraExtAddTransformer();
1895
1896 // build operator and add
1897 const RCP<Thyra::LinearOpBase<double> > explicitOp = prodTrans->createOutputOp();
1898 prodTrans->transform(*implicitOp,explicitOp.ptr());
1899 explicitOp->setObjectLabel("explicit( " + opl->getObjectLabel() +
1900 " + " + opr->getObjectLabel() + " )");
1901
1902 return explicitOp;
1903 }
1904}
1905
1918const ModifiableLinearOp explicitAdd(const LinearOp & opl,const LinearOp & opr,
1919 const ModifiableLinearOp & destOp)
1920{
1921 // if blocked, add block by block
1922 if(isPhysicallyBlockedLinearOp(opl)){
1923 TEUCHOS_ASSERT(isPhysicallyBlockedLinearOp(opr));
1924
1925 double scalarl = 0.0;
1926 bool transpl = false;
1927 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opl = getPhysicallyBlockedLinearOp(opl, &scalarl, &transpl);
1928
1929 double scalarr = 0.0;
1930 bool transpr = false;
1931 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_opr = getPhysicallyBlockedLinearOp(opr, &scalarr, &transpr);
1932
1933 int numRows = blocked_opl->productRange()->numBlocks();
1934 int numCols = blocked_opl->productDomain()->numBlocks();
1935 TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == numRows);
1936 TEUCHOS_ASSERT(blocked_opr->productDomain()->numBlocks() == numCols);
1937
1938 RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_sum = Teuchos::rcp_dynamic_cast<Thyra::DefaultBlockedLinearOp<double>>(destOp);
1939 if(blocked_sum.is_null()) {
1940 // take care of the null case, this means we must alllocate memory
1941 blocked_sum = Thyra::defaultBlockedLinearOp<double>();
1942
1943 blocked_sum->beginBlockFill(numRows,numCols);
1944 for(int r = 0; r < numRows; ++r) {
1945 for(int c = 0; c < numCols; ++c) {
1946 auto block = explicitAdd(Thyra::scale(scalarl,blocked_opl->getBlock(r,c)),Thyra::scale(scalarr,blocked_opr->getBlock(r,c)),Teuchos::null);
1947 blocked_sum->setNonconstBlock(r,c,block);
1948 }
1949 }
1950 blocked_sum->endBlockFill();
1951
1952 }
1953 else {
1954 // in this case memory can be reused
1955 for(int r = 0; r < numRows; ++r)
1956 for(int c = 0; c < numCols; ++c)
1957 explicitAdd(Thyra::scale(scalarl,blocked_opl->getBlock(r,c)),Thyra::scale(scalarr,blocked_opr->getBlock(r,c)),blocked_sum->getNonconstBlock(r,c));
1958 }
1959
1960 return blocked_sum;
1961 }
1962
1963 bool isTpetral = Teko::TpetraHelpers::isTpetraLinearOp(opl);
1964 bool isTpetrar = Teko::TpetraHelpers::isTpetraLinearOp(opr);
1965
1966 if(isTpetral && isTpetrar){ // Both operators are Tpetra matrices so use the explicit Tpetra matrix-matrix add
1967
1968 // Get left and right Tpetra crs operators
1969 ST scalarl = 0.0;
1970 bool transpl = false;
1971 RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpl = Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1972 ST scalarr = 0.0;
1973 bool transpr = false;
1974 RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOpr = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1975
1976 // Build output operator
1977 RCP<Thyra::LinearOpBase<ST> > explicitOp;
1978 if(destOp!=Teuchos::null)
1979 explicitOp = destOp;
1980 else
1981 explicitOp = rcp(new Thyra::TpetraLinearOp<ST,LO,GO,NT>());
1982 RCP<Thyra::TpetraLinearOp<ST,LO,GO,NT> > tExplicitOp = rcp_dynamic_cast<Thyra::TpetraLinearOp<ST,LO,GO,NT> >(explicitOp);
1983
1984 // Do explicit matrix-matrix add
1985 RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > explicitCrsOp = Tpetra::MatrixMatrix::add<ST,LO,GO,NT>(scalarl,transpl,*tCrsOpl,scalarr,transpr,*tCrsOpr);
1986 tExplicitOp->initialize(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(explicitCrsOp->getDomainMap()),explicitCrsOp);
1987 return tExplicitOp;
1988
1989 }else{ // Assume Epetra
1990
1991 // build implicit add
1992 const LinearOp implicitOp = Thyra::add(opl,opr);
1993
1994 // build transformer
1995 const RCP<Thyra::LinearOpTransformerBase<double> > prodTrans =
1996 Thyra::epetraExtAddTransformer();
1997
1998 // build or reuse destination operator
1999 RCP<Thyra::LinearOpBase<double> > explicitOp;
2000 if(destOp!=Teuchos::null)
2001 explicitOp = destOp;
2002 else
2003 explicitOp = prodTrans->createOutputOp();
2004
2005 // perform add
2006 prodTrans->transform(*implicitOp,explicitOp.ptr());
2007 explicitOp->setObjectLabel("explicit( " + opl->getObjectLabel() +
2008 " + " + opr->getObjectLabel() + " )");
2009
2010 return explicitOp;
2011 }
2012}
2013
2018const ModifiableLinearOp explicitSum(const LinearOp & op,
2019 const ModifiableLinearOp & destOp)
2020{
2021 // convert operators to Epetra_CrsMatrix
2022 const RCP<const Epetra_CrsMatrix> epetraOp =
2023 rcp_dynamic_cast<const Epetra_CrsMatrix>(get_Epetra_Operator(*op), true);
2024
2025 if(destOp==Teuchos::null) {
2026 Teuchos::RCP<Epetra_Operator> epetraDest = Teuchos::rcp(new Epetra_CrsMatrix(*epetraOp));
2027
2028 return Thyra::nonconstEpetraLinearOp(epetraDest);
2029 }
2030
2031 const RCP<Epetra_CrsMatrix> epetraDest =
2032 rcp_dynamic_cast<Epetra_CrsMatrix>(get_Epetra_Operator(*destOp), true);
2033
2034 EpetraExt::MatrixMatrix::Add(*epetraOp,false,1.0,*epetraDest,1.0);
2035
2036 return destOp;
2037}
2038
2039const LinearOp explicitTranspose(const LinearOp & op)
2040{
2041 if(Teko::TpetraHelpers::isTpetraLinearOp(op)){
2042
2043 RCP<const Thyra::TpetraLinearOp<ST,LO,GO,NT> > tOp = rcp_dynamic_cast<const Thyra::TpetraLinearOp<ST,LO,GO,NT> >(op,true);
2044 RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOp = rcp_dynamic_cast<const Tpetra::CrsMatrix<ST,LO,GO,NT> >(tOp->getConstTpetraOperator(),true);
2045
2046 Tpetra::RowMatrixTransposer<ST,LO,GO,NT> transposer(tCrsOp);
2047 RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > transOp = transposer.createTranspose();
2048
2049 return Thyra::constTpetraLinearOp<ST,LO,GO,NT>(Thyra::tpetraVectorSpace<ST,LO,GO,NT>(transOp->getRangeMap()),Thyra::tpetraVectorSpace<ST,LO,GO,NT>(transOp->getDomainMap()),transOp);
2050
2051 } else {
2052
2053 RCP<const Epetra_Operator> eOp = Thyra::get_Epetra_Operator(*op);
2054 TEUCHOS_TEST_FOR_EXCEPTION(eOp==Teuchos::null,std::logic_error,
2055 "Teko::explicitTranspose Not an Epetra_Operator");
2056 RCP<const Epetra_RowMatrix> eRowMatrixOp
2057 = Teuchos::rcp_dynamic_cast<const Epetra_RowMatrix>(eOp);
2058 TEUCHOS_TEST_FOR_EXCEPTION(eRowMatrixOp==Teuchos::null,std::logic_error,
2059 "Teko::explicitTranspose Not an Epetra_RowMatrix");
2060
2061 // we now have a delete transpose operator
2062 EpetraExt::RowMatrix_Transpose tranposeOp;
2063 Epetra_RowMatrix & eMat = tranposeOp(const_cast<Epetra_RowMatrix &>(*eRowMatrixOp));
2064
2065 // this copy is because of a poor implementation of the EpetraExt::Transform
2066 // implementation
2067 Teuchos::RCP<Epetra_CrsMatrix> crsMat
2068 = Teuchos::rcp(new Epetra_CrsMatrix(dynamic_cast<Epetra_CrsMatrix &>(eMat)));
2069
2070 return Thyra::epetraLinearOp(crsMat);
2071 }
2072}
2073
2074double frobeniusNorm(const LinearOp & op_in)
2075{
2076 LinearOp op;
2077 double scalar = 1.0;
2078
2079 // if blocked, must be 1x1
2080 if(isPhysicallyBlockedLinearOp(op_in)){
2081 bool transp = false;
2082 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > blocked_op = getPhysicallyBlockedLinearOp(op_in,&scalar,&transp);
2083 TEUCHOS_ASSERT(blocked_op->productRange()->numBlocks() == 1);
2084 TEUCHOS_ASSERT(blocked_op->productDomain()->numBlocks() == 1);
2085 op = blocked_op->getBlock(0,0);
2086 } else
2087 op = op_in;
2088
2089 if(Teko::TpetraHelpers::isTpetraLinearOp(op)){
2090 const RCP<const Thyra::TpetraLinearOp<ST,LO,GO,NT> > tOp = rcp_dynamic_cast<const Thyra::TpetraLinearOp<ST,LO,GO,NT> >(op);
2091 const RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > crsOp = rcp_dynamic_cast<const Tpetra::CrsMatrix<ST,LO,GO,NT> >(tOp->getConstTpetraOperator(),true);
2092 return crsOp->getFrobeniusNorm();
2093 } else {
2094 const RCP<const Epetra_Operator> epOp = Thyra::get_Epetra_Operator(*op);
2095 const RCP<const Epetra_CrsMatrix> crsOp = rcp_dynamic_cast<const Epetra_CrsMatrix>(epOp,true);
2096 return crsOp->NormFrobenius();
2097 }
2098}
2099
2100double oneNorm(const LinearOp & op)
2101{
2102 if(Teko::TpetraHelpers::isTpetraLinearOp(op)){
2103 TEUCHOS_TEST_FOR_EXCEPTION(true,std::logic_error,"One norm not currently implemented for Tpetra matrices");
2104
2105 } else {
2106 const RCP<const Epetra_Operator> epOp = Thyra::get_Epetra_Operator(*op);
2107 const RCP<const Epetra_CrsMatrix> crsOp = rcp_dynamic_cast<const Epetra_CrsMatrix>(epOp,true);
2108 return crsOp->NormOne();
2109 }
2110}
2111
2112double infNorm(const LinearOp & op)
2113{
2114 if(Teko::TpetraHelpers::isTpetraLinearOp(op)){
2115 ST scalar = 0.0;
2116 bool transp = false;
2117 RCP<const Tpetra::CrsMatrix<ST,LO,GO,NT> > tCrsOp = Teko::TpetraHelpers::getTpetraCrsMatrix(op, &scalar, &transp);
2118
2119 // extract diagonal
2120 const RCP<Tpetra::Vector<ST,LO,GO,NT> > ptrDiag = Tpetra::createVector<ST,LO,GO,NT>(tCrsOp->getRowMap());
2121 Tpetra::Vector<ST,LO,GO,NT> & diag = *ptrDiag;
2122
2123 // compute absolute value row sum
2124 diag.putScalar(0.0);
2125 for(LO i=0;i<(LO) tCrsOp->getNodeNumRows();i++) {
2126 LO numEntries = tCrsOp->getNumEntriesInLocalRow (i);
2127 typename Tpetra::CrsMatrix<ST,LO,GO,NT>::local_inds_host_view_type indices;
2128 typename Tpetra::CrsMatrix<ST,LO,GO,NT>::values_host_view_type values;
2129 tCrsOp->getLocalRowView(i,indices,values);
2130
2131 // build abs value row sum
2132 for(LO j=0;j<numEntries;j++)
2133 diag.sumIntoLocalValue(i,std::abs(values(j)));
2134 }
2135 return diag.normInf()*scalar;
2136
2137 } else {
2138 const RCP<const Epetra_Operator> epOp = Thyra::get_Epetra_Operator(*op);
2139 const RCP<const Epetra_CrsMatrix> crsOp = rcp_dynamic_cast<const Epetra_CrsMatrix>(epOp,true);
2140 return crsOp->NormInf();
2141 }
2142}
2143
2144const LinearOp buildDiagonal(const MultiVector & src,const std::string & lbl)
2145{
2146 RCP<Thyra::VectorBase<double> > dst = Thyra::createMember(src->range());
2147 Thyra::copy(*src->col(0),dst.ptr());
2148
2149 return Thyra::diagonal<double>(dst,lbl);
2150}
2151
2152const LinearOp buildInvDiagonal(const MultiVector & src,const std::string & lbl)
2153{
2154 const RCP<const Thyra::VectorBase<double> > srcV = src->col(0);
2155 RCP<Thyra::VectorBase<double> > dst = Thyra::createMember(srcV->range());
2156 Thyra::reciprocal<double>(*srcV,dst.ptr());
2157
2158 return Thyra::diagonal<double>(dst,lbl);
2159}
2160
2162BlockedMultiVector buildBlockedMultiVector(const std::vector<MultiVector> & mvv)
2163{
2164 Teuchos::Array<MultiVector> mvA;
2165 Teuchos::Array<VectorSpace> vsA;
2166
2167 // build arrays of multi vectors and vector spaces
2168 std::vector<MultiVector>::const_iterator itr;
2169 for(itr=mvv.begin();itr!=mvv.end();++itr) {
2170 mvA.push_back(*itr);
2171 vsA.push_back((*itr)->range());
2172 }
2173
2174 // construct the product vector space
2175 const RCP<const Thyra::DefaultProductVectorSpace<double> > vs
2176 = Thyra::productVectorSpace<double>(vsA);
2177
2178 return Thyra::defaultProductMultiVector<double>(vs,mvA);
2179}
2180
2191Teuchos::RCP<Thyra::VectorBase<double> > indicatorVector(
2192 const std::vector<int> & indices,
2193 const VectorSpace & vs,
2194 double onValue, double offValue)
2195
2196{
2197 using Teuchos::RCP;
2198
2199 // create a new vector
2200 RCP<Thyra::VectorBase<double> > v = Thyra::createMember<double>(vs);
2201 Thyra::put_scalar<double>(offValue,v.ptr()); // fill it with "off" values
2202
2203 // set on values
2204 for(std::size_t i=0;i<indices.size();i++)
2205 Thyra::set_ele<double>(indices[i],onValue,v.ptr());
2206
2207 return v;
2208}
2209
2234double computeSpectralRad(const RCP<const Thyra::LinearOpBase<double> > & A, double tol,
2235 bool isHermitian,int numBlocks,int restart,int verbosity)
2236{
2237 typedef Thyra::LinearOpBase<double> OP;
2238 typedef Thyra::MultiVectorBase<double> MV;
2239
2240 int startVectors = 1;
2241
2242 // construct an initial guess
2243 const RCP<MV> ivec = Thyra::createMember(A->domain());
2244 Thyra::randomize(-1.0,1.0,ivec.ptr());
2245
2246 RCP<Anasazi::BasicEigenproblem<double,MV,OP> > eigProb
2247 = rcp(new Anasazi::BasicEigenproblem<double,MV,OP>(A,ivec));
2248 eigProb->setNEV(1);
2249 eigProb->setHermitian(isHermitian);
2250
2251 // set the problem up
2252 if(not eigProb->setProblem()) {
2253 // big time failure!
2254 return Teuchos::ScalarTraits<double>::nan();
2255 }
2256
2257 // we want largert magnitude eigenvalue
2258 std::string which("LM"); // largest magnitude
2259
2260 // Create the parameter list for the eigensolver
2261 // verbosity+=Anasazi::TimingDetails;
2262 Teuchos::ParameterList MyPL;
2263 MyPL.set( "Verbosity", verbosity );
2264 MyPL.set( "Which", which );
2265 MyPL.set( "Block Size", startVectors );
2266 MyPL.set( "Num Blocks", numBlocks );
2267 MyPL.set( "Maximum Restarts", restart );
2268 MyPL.set( "Convergence Tolerance", tol );
2269
2270 // build status test manager
2271 // RCP<Anasazi::StatusTestMaxIters<double,MV,OP> > statTest
2272 // = rcp(new Anasazi::StatusTestMaxIters<double,MV,OP>(numBlocks*(restart+1)));
2273
2274 // Create the Block Krylov Schur solver
2275 // This takes as inputs the eigenvalue problem and the solver parameters
2276 Anasazi::BlockKrylovSchurSolMgr<double,MV,OP> MyBlockKrylovSchur(eigProb, MyPL );
2277
2278 // Solve the eigenvalue problem, and save the return code
2279 Anasazi::ReturnType solverreturn = MyBlockKrylovSchur.solve();
2280
2281 if(solverreturn==Anasazi::Unconverged) {
2282 double real = MyBlockKrylovSchur.getRitzValues().begin()->realpart;
2283 double comp = MyBlockKrylovSchur.getRitzValues().begin()->imagpart;
2284
2285 return -std::abs(std::sqrt(real*real+comp*comp));
2286
2287 // cout << "Anasazi::BlockKrylovSchur::solve() did not converge!" << std::endl;
2288 // return -std::abs(MyBlockKrylovSchur.getRitzValues().begin()->realpart);
2289 }
2290 else { // solverreturn==Anasazi::Converged
2291 double real = eigProb->getSolution().Evals.begin()->realpart;
2292 double comp = eigProb->getSolution().Evals.begin()->imagpart;
2293
2294 return std::abs(std::sqrt(real*real+comp*comp));
2295
2296 // cout << "Anasazi::BlockKrylovSchur::solve() converged!" << endl;
2297 // return std::abs(eigProb->getSolution().Evals.begin()->realpart);
2298 }
2299}
2300
2324double computeSmallestMagEig(const RCP<const Thyra::LinearOpBase<double> > & A, double tol,
2325 bool isHermitian,int numBlocks,int restart,int verbosity)
2326{
2327 typedef Thyra::LinearOpBase<double> OP;
2328 typedef Thyra::MultiVectorBase<double> MV;
2329
2330 int startVectors = 1;
2331
2332 // construct an initial guess
2333 const RCP<MV> ivec = Thyra::createMember(A->domain());
2334 Thyra::randomize(-1.0,1.0,ivec.ptr());
2335
2336 RCP<Anasazi::BasicEigenproblem<double,MV,OP> > eigProb
2337 = rcp(new Anasazi::BasicEigenproblem<double,MV,OP>(A,ivec));
2338 eigProb->setNEV(1);
2339 eigProb->setHermitian(isHermitian);
2340
2341 // set the problem up
2342 if(not eigProb->setProblem()) {
2343 // big time failure!
2344 return Teuchos::ScalarTraits<double>::nan();
2345 }
2346
2347 // we want largert magnitude eigenvalue
2348 std::string which("SM"); // smallest magnitude
2349
2350 // Create the parameter list for the eigensolver
2351 Teuchos::ParameterList MyPL;
2352 MyPL.set( "Verbosity", verbosity );
2353 MyPL.set( "Which", which );
2354 MyPL.set( "Block Size", startVectors );
2355 MyPL.set( "Num Blocks", numBlocks );
2356 MyPL.set( "Maximum Restarts", restart );
2357 MyPL.set( "Convergence Tolerance", tol );
2358
2359 // build status test manager
2360 // RCP<Anasazi::StatusTestMaxIters<double,MV,OP> > statTest
2361 // = rcp(new Anasazi::StatusTestMaxIters<double,MV,OP>(10));
2362
2363 // Create the Block Krylov Schur solver
2364 // This takes as inputs the eigenvalue problem and the solver parameters
2365 Anasazi::BlockKrylovSchurSolMgr<double,MV,OP> MyBlockKrylovSchur(eigProb, MyPL );
2366
2367 // Solve the eigenvalue problem, and save the return code
2368 Anasazi::ReturnType solverreturn = MyBlockKrylovSchur.solve();
2369
2370 if(solverreturn==Anasazi::Unconverged) {
2371 // cout << "Anasazi::BlockKrylovSchur::solve() did not converge! " << std::endl;
2372 return -std::abs(MyBlockKrylovSchur.getRitzValues().begin()->realpart);
2373 }
2374 else { // solverreturn==Anasazi::Converged
2375 // cout << "Anasazi::BlockKrylovSchur::solve() converged!" << endl;
2376 return std::abs(eigProb->getSolution().Evals.begin()->realpart);
2377 }
2378}
2379
2388ModifiableLinearOp getDiagonalOp(const Teko::LinearOp & A,const DiagonalType & dt)
2389{
2390 switch(dt) {
2391 case Diagonal:
2392 return getDiagonalOp(A);
2393 case Lumped:
2394 return getLumpedMatrix(A);
2395 case AbsRowSum:
2396 return getAbsRowSumMatrix(A);
2397 case NotDiag:
2398 default:
2399 TEUCHOS_TEST_FOR_EXCEPT(true);
2400 };
2401
2402 return Teuchos::null;
2403}
2404
2413ModifiableLinearOp getInvDiagonalOp(const Teko::LinearOp & A,const Teko::DiagonalType & dt)
2414{
2415 switch(dt) {
2416 case Diagonal:
2417 return getInvDiagonalOp(A);
2418 case Lumped:
2419 return getInvLumpedMatrix(A);
2420 case AbsRowSum:
2421 return getAbsRowSumInvMatrix(A);
2422 case NotDiag:
2423 default:
2424 TEUCHOS_TEST_FOR_EXCEPT(true);
2425 };
2426
2427 return Teuchos::null;
2428}
2429
2436std::string getDiagonalName(const DiagonalType & dt)
2437{
2438 switch(dt) {
2439 case Diagonal:
2440 return "Diagonal";
2441 case Lumped:
2442 return "Lumped";
2443 case AbsRowSum:
2444 return "AbsRowSum";
2445 case NotDiag:
2446 return "NotDiag";
2447 case BlkDiag:
2448 return "BlkDiag";
2449 };
2450
2451 return "<error>";
2452}
2453
2462DiagonalType getDiagonalType(std::string name)
2463{
2464 if(name=="Diagonal")
2465 return Diagonal;
2466 if(name=="Lumped")
2467 return Lumped;
2468 if(name=="AbsRowSum")
2469 return AbsRowSum;
2470 if(name=="BlkDiag")
2471 return BlkDiag;
2472
2473 return NotDiag;
2474}
2475
2476LinearOp probe(Teuchos::RCP<const Epetra_CrsGraph> &G,const LinearOp & Op){
2477#ifdef Teko_ENABLE_Isorropia
2478 Teuchos::ParameterList probeList;
2479 Prober prober(G,probeList,true);
2480 Teuchos::RCP<Epetra_CrsMatrix> Mat=rcp(new Epetra_CrsMatrix(Copy,*G));
2482 prober.probe(Mwrap,*Mat);
2483 return Thyra::epetraLinearOp(Mat);
2484#else
2485 (void)G;
2486 (void)Op;
2487 TEUCHOS_TEST_FOR_EXCEPTION(true,std::runtime_error,"Probe requires Isorropia");
2488#endif
2489}
2490
2491double norm_1(const MultiVector & v,std::size_t col)
2492{
2493 Teuchos::Array<double> n(v->domain()->dim());
2494 Thyra::norms_1<double>(*v,n);
2495
2496 return n[col];
2497}
2498
2499double norm_2(const MultiVector & v,std::size_t col)
2500{
2501 Teuchos::Array<double> n(v->domain()->dim());
2502 Thyra::norms_2<double>(*v,n);
2503
2504 return n[col];
2505}
2506
2507void putScalar(const ModifiableLinearOp & op,double scalar)
2508{
2509 try {
2510 // get Epetra_Operator
2511 RCP<Epetra_Operator> eOp = Thyra::get_Epetra_Operator(*op);
2512
2513 // cast it to a CrsMatrix
2514 RCP<Epetra_CrsMatrix> eCrsOp = rcp_dynamic_cast<Epetra_CrsMatrix>(eOp,true);
2515
2516 eCrsOp->PutScalar(scalar);
2517 }
2518 catch (std::exception & e) {
2519 RCP<Teuchos::FancyOStream> out = Teuchos::VerboseObjectBase::getDefaultOStream();
2520
2521 *out << "Teko: putScalar requires an Epetra_CrsMatrix\n";
2522 *out << " Could not extract an Epetra_Operator from a \"" << op->description() << std::endl;
2523 *out << " OR\n";
2524 *out << " Could not cast an Epetra_Operator to a Epetra_CrsMatrix\n";
2525 *out << std::endl;
2526 *out << "*** THROWN EXCEPTION ***\n";
2527 *out << e.what() << std::endl;
2528 *out << "************************\n";
2529
2530 throw e;
2531 }
2532}
2533
2534void clipLower(MultiVector & v,double lowerBound)
2535{
2536 using Teuchos::RCP;
2537 using Teuchos::rcp_dynamic_cast;
2538
2539 // cast so entries are accessible
2540 // RCP<Thyra::SpmdMultiVectorBase<double> > spmdMVec
2541 // = rcp_dynamic_cast<Thyra::DefaultSpmdMultiVector<double> >(v);
2542
2543 for(Thyra::Ordinal i=0;i<v->domain()->dim();i++) {
2544 RCP<Thyra::SpmdVectorBase<double> > spmdVec
2545 = rcp_dynamic_cast<Thyra::SpmdVectorBase<double> >(v->col(i),true);
2546
2547 Teuchos::ArrayRCP<double> values;
2548 // spmdMVec->getNonconstLocalData(Teuchos::ptrFromRef(values),Teuchos::ptrFromRef(i));
2549 spmdVec->getNonconstLocalData(Teuchos::ptrFromRef(values));
2550 for(Teuchos::ArrayRCP<double>::size_type j=0;j<values.size();j++) {
2551 if(values[j]<lowerBound)
2552 values[j] = lowerBound;
2553 }
2554 }
2555}
2556
2557void clipUpper(MultiVector & v,double upperBound)
2558{
2559 using Teuchos::RCP;
2560 using Teuchos::rcp_dynamic_cast;
2561
2562 // cast so entries are accessible
2563 // RCP<Thyra::SpmdMultiVectorBase<double> > spmdMVec
2564 // = rcp_dynamic_cast<Thyra::DefaultSpmdMultiVector<double> >(v);
2565 for(Thyra::Ordinal i=0;i<v->domain()->dim();i++) {
2566 RCP<Thyra::SpmdVectorBase<double> > spmdVec
2567 = rcp_dynamic_cast<Thyra::SpmdVectorBase<double> >(v->col(i),true);
2568
2569 Teuchos::ArrayRCP<double> values;
2570 // spmdMVec->getNonconstLocalData(Teuchos::ptrFromRef(values),Teuchos::ptrFromRef(i));
2571 spmdVec->getNonconstLocalData(Teuchos::ptrFromRef(values));
2572 for(Teuchos::ArrayRCP<double>::size_type j=0;j<values.size();j++) {
2573 if(values[j]>upperBound)
2574 values[j] = upperBound;
2575 }
2576 }
2577}
2578
2579void replaceValue(MultiVector & v,double currentValue,double newValue)
2580{
2581 using Teuchos::RCP;
2582 using Teuchos::rcp_dynamic_cast;
2583
2584 // cast so entries are accessible
2585 // RCP<Thyra::SpmdMultiVectorBase<double> > spmdMVec
2586 // = rcp_dynamic_cast<Thyra::SpmdMultiVectorBase<double> >(v,true);
2587 for(Thyra::Ordinal i=0;i<v->domain()->dim();i++) {
2588 RCP<Thyra::SpmdVectorBase<double> > spmdVec
2589 = rcp_dynamic_cast<Thyra::SpmdVectorBase<double> >(v->col(i),true);
2590
2591 Teuchos::ArrayRCP<double> values;
2592 // spmdMVec->getNonconstLocalData(Teuchos::ptrFromRef(values),Teuchos::ptrFromRef(i));
2593 spmdVec->getNonconstLocalData(Teuchos::ptrFromRef(values));
2594 for(Teuchos::ArrayRCP<double>::size_type j=0;j<values.size();j++) {
2595 if(values[j]==currentValue)
2596 values[j] = newValue;
2597 }
2598 }
2599}
2600
2601void columnAverages(const MultiVector & v,std::vector<double> & averages)
2602{
2603 averages.resize(v->domain()->dim());
2604
2605 // sum over each column
2606 Thyra::sums<double>(*v,averages);
2607
2608 // build averages
2609 Thyra::Ordinal rows = v->range()->dim();
2610 for(std::size_t i=0;i<averages.size();i++)
2611 averages[i] = averages[i]/rows;
2612}
2613
2614double average(const MultiVector & v)
2615{
2616 Thyra::Ordinal rows = v->range()->dim();
2617 Thyra::Ordinal cols = v->domain()->dim();
2618
2619 std::vector<double> averages;
2620 columnAverages(v,averages);
2621
2622 double sum = 0.0;
2623 for(std::size_t i=0;i<averages.size();i++)
2624 sum += averages[i] * rows;
2625
2626 return sum/(rows*cols);
2627}
2628
2629bool isPhysicallyBlockedLinearOp(const LinearOp & op)
2630{
2631 // See if the operator is a PBLO
2632 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > pblo = rcp_dynamic_cast<const Thyra::PhysicallyBlockedLinearOpBase<double> >(op);
2633 if (!pblo.is_null())
2634 return true;
2635
2636 // See if the operator is a wrapped PBLO
2637 ST scalar = 0.0;
2638 Thyra::EOpTransp transp = Thyra::NOTRANS;
2639 RCP<const Thyra::LinearOpBase<ST> > wrapped_op;
2640 Thyra::unwrap(op, &scalar, &transp, &wrapped_op);
2641 pblo = rcp_dynamic_cast<const Thyra::PhysicallyBlockedLinearOpBase<double> >(wrapped_op);
2642 if (!pblo.is_null())
2643 return true;
2644
2645 return false;
2646}
2647
2648RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > getPhysicallyBlockedLinearOp(const LinearOp & op, ST *scalar, bool *transp)
2649{
2650 // If the operator is a TpetraLinearOp
2651 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > pblo = rcp_dynamic_cast<const Thyra::PhysicallyBlockedLinearOpBase<double> >(op);
2652 if(!pblo.is_null()){
2653 *scalar = 1.0;
2654 *transp = false;
2655 return pblo;
2656 }
2657
2658 // If the operator is a wrapped TpetraLinearOp
2659 RCP<const Thyra::LinearOpBase<ST> > wrapped_op;
2660 Thyra::EOpTransp eTransp = Thyra::NOTRANS;
2661 Thyra::unwrap(op, scalar, &eTransp, &wrapped_op);
2662 pblo = rcp_dynamic_cast<const Thyra::PhysicallyBlockedLinearOpBase<double> >(wrapped_op,true);
2663 if(!pblo.is_null()){
2664 *transp = true;
2665 if(eTransp == Thyra::NOTRANS)
2666 *transp = false;
2667 return pblo;
2668 }
2669
2670 return Teuchos::null;
2671}
2672
2673
2674}
void scale(const double alpha, MultiVector &x)
Scale a multivector by a constant.
DiagonalType
Type describing the type of diagonal to construct.
@ BlkDiag
Specifies that a block diagonal approximation is to be used.
@ NotDiag
For user convenience, if Teko recieves this value, exceptions will be thrown
@ AbsRowSum
Specifies that the diagonal entry is .
@ Diagonal
Specifies that just the diagonal is used.
@ Lumped
Specifies that row sum is used to form a diagonal.
int blockRowCount(const BlockedLinearOp &blo)
Get the row count in a block linear operator.
int blockColCount(const BlockedLinearOp &blo)
Get the column count in a block linear operator.
BlockedLinearOp createBlockedOp()
Build a new blocked linear operator.
Implements the Epetra_Operator interface with a Thyra LinearOperator. This enables the use of absrtac...