/** \brief Function that is called to build the preconditioner
  *        for the linear operator that is passed in.
  */
LinearOp IterativePreconditionerFactory::buildPreconditionerOperator(LinearOp & lo,PreconditionerState & state) const
{
   TEUCHOS_TEST_FOR_EXCEPTION(precFactory_==Teuchos::null,std::runtime_error,
                      "ERROR: Teko::IterativePreconditionerFactory::buildPreconditionerOperator requires that a "
                   << "preconditioner factory has been set. Currently it is null!");

   // build user specified preconditioner
   ModifiableLinearOp & invP = state.getModifiableOp("prec");
   if(invP==Teuchos::null)
      invP = Teko::buildInverse(*precFactory_,lo);
   else
      Teko::rebuildInverse(*precFactory_,lo,invP);

   // // no repititions are required
   // if(correctionNum_==0) return invP;

   // now build correction operator
   LinearOp I = Thyra::identity(lo->range(),"I");
   LinearOp AiP = multiply(lo,invP.getConst(),"AiP");
   LinearOp correction = add(I,scale(-1.0,AiP)); // I - A * iPA

   LinearOp resMap = I; // will map r_n to r_{n+1}
   for(unsigned int i=0;i<correctionNum_;i++) 
      resMap = add(I,multiply(resMap,correction)); // resMap = I + resMap*(I-A*iP)
   
   // iP = (I-A*iP)^{correctionNum}
   return multiply(invP.getConst(),resMap);
}
Exemplo n.º 2
0
//! Get the domain space of a linear operator
inline VectorSpace domainSpace(const LinearOp & lo)
{ return lo->domain(); }
Exemplo n.º 3
0
//! Get the range space of a linear operator
inline VectorSpace rangeSpace(const LinearOp & lo)
{ return lo->range(); }
  bool tSIMPLEPreconditionerFactory::test_diagonal(int verbosity,std::ostream & os,int use_blocking)
{
   // make sure the preconditioner is working by testing against the identity matrix
   typedef RCP<const Thyra::LinearOpBase<double> > LinearOp;

   bool status = false;
   bool allPassed = true;
   double vec[2];
   double diff = 0.0;

   // build 4x4 matrix with block 2x2 diagonal subblocks
   //
   //            [ 1 0 7 0 ]
   // [ F G ] =  [ 0 2 0 8 ]
   // [ D C ]    [ 5 0 3 0 ]
   //            [ 0 6 0 4 ]
   //

   vec[0] = 1.0; vec[1] = 2.0;
   LinearOp F = Teko::Test::DiagMatrix(2,vec,"F");

   vec[0] = 7.0; vec[1] = 8.0;
   LinearOp G = Teko::Test::DiagMatrix(2,vec,"G");

   vec[0] = 5.0; vec[1] = 6.0;
   LinearOp D = Teko::Test::DiagMatrix(2,vec,"D");

   vec[0] = 3.0; vec[1] = 4.0;
   LinearOp C = Teko::Test::DiagMatrix(2,vec,"C");

   vec[0] = 1.0; vec[1] = 0.5;
   LinearOp iF = Teko::Test::DiagMatrix(2,vec,"inv(F)");

   vec[0] = -0.03125; vec[1] = -0.05;
   LinearOp iS = Teko::Test::DiagMatrix(2,vec,"inv(S)");

   RCP<Teko::InverseFactory> invF = rcp(new Teko::StaticOpInverseFactory(iF));
   RCP<Teko::InverseFactory> invS = rcp(new Teko::StaticOpInverseFactory(iS));

   LinearOp A = Thyra::block2x2(F,G,D,C);
   RCP<SIMPLEPreconditionerFactory> sFactory = rcp(new SIMPLEPreconditionerFactory(invF,invS,0.9));
   const RCP<const Thyra::PreconditionerFactoryBase<double> > precFactory =sFactory;
   RCP<Thyra::PreconditionerBase<double> > prec = Thyra::prec<double>(*precFactory,A);

   if(use_blocking==1){
     // parameter list for (1,1) block
     Teuchos::ParameterList List,BlkList;   
     BlkList.set("number of local blocks",1);
     BlkList.set("block start index",&*block_starts_);
     BlkList.set("block entry gids",&*block_gids_);
     List.set("H options",BlkList);
     List.set("Explicit Velocity Inverse Type","BlkDiag");
     List.set("Inverse Pressure Type","Amesos");
     sFactory->initializeFromParameterList(List);
   }
   else if(use_blocking==2){
     Teuchos::ParameterList List,BlkList;   
     BlkList.set("contiguous block size",2);
     List.set("H options",BlkList);
     List.set("Explicit Velocity Inverse Type","BlkDiag");
     List.set("Inverse Pressure Type","Amesos");
     sFactory->initializeFromParameterList(List);
   }

   // build linear operator
   RCP<const Thyra::LinearOpBase<double> > precOp = prec->getUnspecifiedPrecOp();

   const RCP<Epetra_Map> map = rcp(new Epetra_Map(2,0,*comm));
   // construct a couple of vectors
   Epetra_Vector ea(*map),eb(*map);
   Epetra_Vector ef(*map),eg(*map);
   const RCP<const Thyra::MultiVectorBase<double> > x = BlockVector(ea,eb,A->domain());
   const RCP<const Thyra::MultiVectorBase<double> > z = BlockVector(ef,eg,A->domain());
   const RCP<Thyra::MultiVectorBase<double> > y = Thyra::createMembers(A->range(),1); 

   // now checks of the preconditioner (should be exact!)
   /////////////////////////////////////////////////////////////////////////

   // test vector [0 1 1 3]
   ea[0] = 0.0; ea[1] = 1.0; eb[0] = 1.0; eb[1] = 3.0;
   ef[0] =  0.21875; ef[1]  =  0.5;  
   eg[0] = -0.028125; eg[1] =  0.0;
   Thyra::apply(*precOp,Thyra::NOTRANS,*x,y.ptr());
   status = ((diff = Teko::Test::Difference(y,z)/Thyra::norm_2(*z->col(0)))<tolerance_);
   if(not status || verbosity>=10 ) { 
      os << std::endl << "   tSIMPLEPreconditionerFactory::test_diagonal " << toString(status) << ":  (y=inv(A)*x) != z (|y-z|_2/|z|_2 = " 
                      << diff << ")" << std::endl;
      os << "      "; Print(os,"x",x);
      os << "      "; Print(os,"y",y);
      os << "      "; Print(os,"z",z);
   }
   allPassed &= status;

   // test vector [-2 4 7 9]
   ea[0] =-2.0; ea[1] = 4.0; eb[0] = 7.0; eb[1] = 9.0;
   ef[0] = 1.71875; ef[1] =  1.4;
   eg[0] = -0.478125; eg[1] = 0.135;
   Thyra::apply(*precOp,Thyra::NOTRANS,*x,y.ptr());
   status = ((diff = Teko::Test::Difference(y,z)/Thyra::norm_2(*z->col(0)))<tolerance_);
   if(not status || verbosity>=10 ) { 
      os << std::endl << "   tSIMPLEPreconditionerFactory::test_diagonal " << toString(status) << ":  (y=inv(A)*x) != z (|y-z|_2/|z|_2 = " 
                      << diff << ")" << std::endl;
      os << "      "; Print(os,"x",x);
      os << "      "; Print(os,"y",y);
      os << "      "; Print(os,"z",z);
   }
   allPassed &= status;

   // test vector [1 0 0 -5]
   ea[0] = 1.0; ea[1] = 0.0; eb[0] = 0.0; eb[1] =-5.0;
   ef[0] = -0.09375; ef[1] = -1.0;
   eg[0] =  0.140625; eg[1] =  0.225;
   Thyra::apply(*precOp,Thyra::NOTRANS,*x,y.ptr());
   status = ((diff = Teko::Test::Difference(y,z)/Thyra::norm_2(*z->col(0)))<tolerance_);
   if(not status || verbosity>=10 ) { 
      os << std::endl << "   tSIMPLEPreconditionerFactory::test_diagonal " << toString(status) << ":  (y=inv(A)*x) != z (|y-z|_2/|z|_2 = " 
                      << diff << ")" << std::endl;
      os << "      "; Print(os,"x",x);
      os << "      "; Print(os,"y",y);
      os << "      "; Print(os,"z",z);
   }
   allPassed &= status;

   // test vector [4 -4 6 12]
   ea[0] = 4.0; ea[1] =-4.0; eb[0] = 6.0; eb[1] =12.0;
   ef[0] = 0.9375; ef[1] =  2.800000000000001;
   eg[0] = 0.39375; eg[1] = -1.08;
   Thyra::apply(*precOp,Thyra::NOTRANS,*x,y.ptr());
   status = ((diff = Teko::Test::Difference(y,z)/Thyra::norm_2(*z->col(0)))<tolerance_);
   if(not status || verbosity>=10 ) { 
      os << std::endl << "   tSIMPLEPreconditionerFactory::test_diagonal " << toString(status) << ":  (y=inv(A)*x) != z (|y-z|_2/|z|_2 = " 
                      << diff << ")" << std::endl;
      os << "      "; Print(os,"x",x);
      os << "      "; Print(os,"y",y);
      os << "      "; Print(os,"z",z);
   }
   allPassed &= status;

   return allPassed;
}
bool tLSCStablePreconditionerFactory::test_diagonal(int verbosity,std::ostream & os)
{
   // make sure the preconditioner is working by testing against the identity matrix
   typedef RCP<const Thyra::LinearOpBase<double> > LinearOp;

   bool status = false;
   bool allPassed = true;
   double vec[2];
   double diff = 0.0;

   // build 4x4 matrix with block 2x2 diagonal subblocks
   //
   //            [ 1 0 7 0 ]
   // [ F G ] =  [ 0 2 0 8 ]
   // [ D C ]    [ 5 0 0 0 ]
   //            [ 0 6 0 0 ]
   //

   vec[0] = 1.0; vec[1] = 2.0;
   LinearOp F = Teko::Test::DiagMatrix(2,vec);

   vec[0] = 7.0; vec[1] = 8.0;
   LinearOp G = Teko::Test::DiagMatrix(2,vec);

   vec[0] = 5.0; vec[1] = 6.0;
   LinearOp D = Teko::Test::DiagMatrix(2,vec);

   vec[0] = 0.0; vec[1] = 0.0;
   LinearOp C = Teko::Test::DiagMatrix(2,vec);

   vec[0] = 1.0; vec[1] = 0.5;
   LinearOp iF = Teko::Test::DiagMatrix(2,vec);

   // S = -C+D*iF*G
   vec[0] = 0.028571428571429; vec[1] = 0.020833333333333; 
   LinearOp iBBt = Teko::Test::DiagMatrix(2,vec);

   LinearOp A = Thyra::block2x2(F,G,D,C);
   const RCP<const Thyra::PreconditionerFactoryBase<double> > precFactory 
         = rcp(new LSCPreconditionerFactory(iF,iBBt,Teuchos::null));
   RCP<Thyra::PreconditionerBase<double> > prec = Thyra::prec<double>(*precFactory,A);

   // build linear operator
   RCP<const Thyra::LinearOpBase<double> > precOp = prec->getUnspecifiedPrecOp();

   const RCP<Epetra_Map> map = rcp(new Epetra_Map(2,0,*comm));
   // construct a couple of vectors
   Epetra_Vector ea(*map),eb(*map);
   Epetra_Vector ef(*map),eg(*map);
   const RCP<const Thyra::MultiVectorBase<double> > x = BlockVector(ea,eb,A->domain());
   const RCP<const Thyra::MultiVectorBase<double> > z = BlockVector(ef,eg,A->domain());
   const RCP<Thyra::MultiVectorBase<double> > y = Thyra::createMembers(A->range(),1); 

   // now checks of the preconditioner (should be exact!)
   /////////////////////////////////////////////////////////////////////////

   // test vector [0 1 1 3]
   ea[0] = 0.0; ea[1] = 1.0; eb[0] = 1.0; eb[1] = 3.0;
   // ef[0] =  0.200000000000000; ef[1] =  0.500000000000000;  
   // eg[0] = -0.028571428571429; eg[1] =  0;
   ef[0] =  0.200000000000000; ef[1] =  1.000000000000000;  
   eg[0] = -0.028571428571429; eg[1] =  -0.125;
   Thyra::apply(*precOp,Thyra::NOTRANS,*x,y.ptr());
   status = ((diff = Teko::Test::Difference(y,z)/Thyra::norm_2(*z->col(0)))<tolerance_);
   if(not status || verbosity>=10 ) { 
      os << std::endl << "   tLSCStablePreconditionerFactory::test_diagonal " << toString(status) << ":  (y=inv(A)*x) != z (|y-z|_2/|z|_2 = " 
                      << diff << ")" << std::endl;
      os << "      "; Print(os,"x",x);
      os << "      "; Print(os,"y",y);
      os << "      "; Print(os,"z",z);
   }
   allPassed &= status;

   // test vector [-2 4 7 9]
   ea[0] =-2.0; ea[1] = 4.0; eb[0] = 7.0; eb[1] = 9.0;
   // ef[0] =  1.400000000000000; ef[1] =  1.500000000000000;
   // eg[0] = -0.485714285714286; eg[1] =  0.125000000000000;
   ef[0] = -0.600000000000000; ef[1] =  3.500000000000000;
   eg[0] = -0.200000000000000; eg[1] = -0.375000000000000;
   Thyra::apply(*precOp,Thyra::NOTRANS,*x,y.ptr());
   status = ((diff = Teko::Test::Difference(y,z)/Thyra::norm_2(*z->col(0)))<tolerance_);
   if(not status || verbosity>=10 ) { 
      os << std::endl << "   tLSCStablePreconditionerFactory::test_diagonal " << toString(status) << ":  (y=inv(A)*x) != z (|y-z|_2/|z|_2 = " 
                      << diff << ")" << std::endl;
      os << "      "; Print(os,"x",x);
      os << "      "; Print(os,"y",y);
      os << "      "; Print(os,"z",z);
   }
   allPassed &= status;

   // test vector [1 0 0 -5]
   ea[0] = 1.0; ea[1] = 0.0; eb[0] = 0.0; eb[1] =-5.0;
   // ef[0] =  0.000000000000000; ef[1] = -0.833333333333333;
   // eg[0] =  0.142857142857143; eg[1] =  0.208333333333333;
   ef[0] =  1.000000000000000; ef[1] = -0.833333333333333;
   eg[0] =  0.000000000000000; eg[1] =  0.208333333333333;
   Thyra::apply(*precOp,Thyra::NOTRANS,*x,y.ptr());
   status = ((diff = Teko::Test::Difference(y,z)/Thyra::norm_2(*z->col(0)))<tolerance_);
   if(not status || verbosity>=10 ) { 
      os << std::endl << "   tLSCStablePreconditionerFactory::test_diagonal " << toString(status) << ":  (y=inv(A)*x) != z (|y-z|_2/|z|_2 = " 
                      << diff << ")" << std::endl;
      os << "      "; Print(os,"x",x);
      os << "      "; Print(os,"y",y);
      os << "      "; Print(os,"z",z);
   }
   allPassed &= status;

   // test vector [4 -4 6 12]
   ea[0] = 4.0; ea[1] =-4.0; eb[0] = 6.0; eb[1] =12.0;
   // ef[0] =  1.200000000000000; ef[1] =  2.000000000000000;
   // eg[0] =  0.400000000000000; eg[1] = -1.000000000000000;
   ef[0] = 05.200000000000000; ef[1] =  0.000000000000000;
   eg[0] = -0.171428571428571; eg[1] = -0.500000000000000;
   Thyra::apply(*precOp,Thyra::NOTRANS,*x,y.ptr());
   status = ((diff = Teko::Test::Difference(y,z)/Thyra::norm_2(*z->col(0)))<tolerance_);
   if(not status || verbosity>=10 ) { 
      os << std::endl << "   tLSCStablePreconditionerFactory::test_diagonal " << toString(status) << ":  (y=inv(A)*x) != z (|y-z|_2/|z|_2 = " 
                      << diff << ")" << std::endl;
      os << "      "; Print(os,"x",x);
      os << "      "; Print(os,"y",y);
      os << "      "; Print(os,"z",z);
   }
   allPassed &= status;

   return allPassed;
}
bool tLSCStablePreconditionerFactory::test_identity(int verbosity,std::ostream & os)
{
   // make sure the preconditioner is working by testing against the identity matrix
   typedef RCP<const Thyra::LinearOpBase<double> > LinearOp;

   bool status = false;
   bool allPassed = true;

   LinearOp Iu = Thyra::identity<double>(invF_->range());
   LinearOp Ip = Thyra::identity<double>(invBQBt_->range());
   LinearOp Zp = Thyra::zero<double>(invBQBt_->range(),invF_->domain());
   LinearOp invBQBt = Ip;

   LinearOp A = Thyra::block2x2(Iu,Ip,Iu,Zp);
   const RCP<const Thyra::PreconditionerFactoryBase<double> > precFactory 
         = rcp(new LSCPreconditionerFactory(Iu,invBQBt,Teuchos::null));
   RCP<Thyra::PreconditionerBase<double> > prec = Thyra::prec<double>(*precFactory,A);

   // build linear operator
   RCP<const Thyra::LinearOpBase<double> > precOp = prec->getUnspecifiedPrecOp();

   const RCP<Epetra_Map> map = rcp(new Epetra_Map(2,0,*comm));
   // construct a couple of vectors
   Epetra_Vector ea(*map),eb(*map);
   Epetra_Vector ef(*map),eg(*map);
   const RCP<const Thyra::MultiVectorBase<double> > x = BlockVector(ea,eb,A->domain());
   const RCP<const Thyra::MultiVectorBase<double> > z = BlockVector(ef,eg,A->domain());
   const RCP<Thyra::MultiVectorBase<double> > y = Thyra::createMembers(A->range(),1); 

   // test vector [0 1 1 3]
   ea[0] = 0.0; ea[1] = 1.0; eb[0] = 1.0; eb[1] = 3.0;
   // ef[0] = eb[0]; ef[1] = eb[1]; eg[0] = ea[0]-eb[0]; eg[1] = ea[1]-eb[1];
   ef[0] = ea[0]+eb[0]; ef[1] = ea[1]+eb[1]; eg[0] = -eb[0]; eg[1] = -eb[1];
   Thyra::apply(*precOp,Thyra::NOTRANS,*x,y.ptr());
   status = (Teko::Test::Difference(y,z)<tolerance_);
   if(not status || verbosity>=10) { 
      os << std::endl << "   tLSCStablePreconditionerFactory::test_Identity " << toString(status) << ": (y=inv(A)*x) != z" << std::endl;
      os << "      "; Print(os,"x",x);
      os << "      "; Print(os,"y",y);
      os << "      "; Print(os,"z",z);
   }
   allPassed &= status;

   // test vector [-2 4 7 9]
   ea[0] =-2.0; ea[1] = 4.0; eb[0] = 7.0; eb[1] = 9.0;
   // ef[0] = eb[0]; ef[1] = eb[1]; eg[0] = ea[0]-eb[0]; eg[1] = ea[1]-eb[1];
   ef[0] = ea[0]+eb[0]; ef[1] = ea[1]+eb[1]; eg[0] = -eb[0]; eg[1] = -eb[1];
   Thyra::apply(*precOp,Thyra::NOTRANS,*x,y.ptr());
   status = (Teko::Test::Difference(y,z)<tolerance_);
   if(not status || verbosity>=10) { 
      os << std::endl << "   tLSCStablePreconditionerFactory::test_Identity " << toString(status) << ": (y=inv(A)*x) != z" << std::endl;
      os << "      "; Print(os,"x",x);
      os << "      "; Print(os,"y",y);
      os << "      "; Print(os,"z",z);
   }
   allPassed &= status;

   // test vector [1 0 0 -5]
   ea[0] = 1.0; ea[1] = 0.0; eb[0] = 0.0; eb[1] =-5.0;
   // ef[0] = eb[0]; ef[1] = eb[1]; eg[0] = ea[0]-eb[0]; eg[1] = ea[1]-eb[1];
   ef[0] = ea[0]+eb[0]; ef[1] = ea[1]+eb[1]; eg[0] = -eb[0]; eg[1] = -eb[1];
   Thyra::apply(*precOp,Thyra::NOTRANS,*x,y.ptr());
   status = (Teko::Test::Difference(y,z)<tolerance_);
   if(not status || verbosity>=10) { 
      os << std::endl << "   tLSCStablePreconditionerFactory::test_Identity " << toString(status) << ": (y=inv(A)*x) != z" << std::endl;
      os << "      "; Print(os,"x",x);
      os << "      "; Print(os,"y",y);
      os << "      "; Print(os,"z",z);
   }
   allPassed &= status;

   // test vector [4 -4 6 12]
   ea[0] = 4.0; ea[1] =-4.0; eb[0] = 6.0; eb[1] =12.0;
   // ef[0] = eb[0]; ef[1] = eb[1]; eg[0] = ea[0]-eb[0]; eg[1] = ea[1]-eb[1];
   ef[0] = ea[0]+eb[0]; ef[1] = ea[1]+eb[1]; eg[0] = -eb[0]; eg[1] = -eb[1];
   Thyra::apply(*precOp,Thyra::NOTRANS,*x,y.ptr());
   status = (Teko::Test::Difference(y,z)<tolerance_);
   if(not status || verbosity>=10) { 
      os << std::endl << "   tLSCStablePreconditionerFactory::test_Identity " << toString(status) << ": (y=inv(A)*x) != z" << std::endl;
      os << "      "; Print(os,"x",x);
      os << "      "; Print(os,"y",y);
      os << "      "; Print(os,"z",z);
   }
   allPassed &= status;

   return allPassed;
}
bool tJacobi2x2PreconditionerFactory::test_identity(int verbosity,std::ostream & os)
{
   // make sure the preconditioner is working by testing against the identity matrix
   typedef RCP<const Thyra::VectorBase<double> > Vector;
   typedef RCP<const Thyra::VectorSpaceBase<double> > VectorSpace;
   typedef RCP<const Thyra::LinearOpBase<double> > LinearOp;

   bool status = false;
   bool allPassed = true;
   double diff = 0.0;

   LinearOp F = Thyra::identity<double>(invF_->range());
   LinearOp G = Thyra::identity<double>(invF_->range());
   LinearOp D = Thyra::identity<double>(invC_->range());
   LinearOp C = Thyra::identity<double>(invC_->range());

   LinearOp A = Thyra::block2x2(F,G,D,C);
   RCP<Thyra::PreconditionerFactoryBase<double> > precFactory 
         = rcp(new JacobiPreconditionerFactory(F,C));
   RCP<Thyra::PreconditionerBase<double> > prec = Thyra::prec<double>(*precFactory,A);

   // build linear operator
   RCP<const Thyra::LinearOpBase<double> > precOp = prec->getUnspecifiedPrecOp();

   const RCP<Epetra_Map> map = rcp(new Epetra_Map(2,0,*comm));
   // construct a couple of vectors
   Epetra_Vector ea(*map),eb(*map);
   const RCP<const Thyra::MultiVectorBase<double> > x = BlockVector(ea,eb,A->domain());
   const RCP<Thyra::MultiVectorBase<double> > y = Thyra::createMembers(A->range(),1); 

   // test vector [0 1 1 3]
   ea[0] = 0.0; ea[1] = 1.0; eb[0] = 1.0; eb[1] = 3.0;
   Thyra::apply(*precOp,Thyra::NOTRANS,*x,y.ptr());
   status = ((diff = Teko::Test::Difference(y,x)/Thyra::norm_2(*x->col(0)))<tolerance_);
   if(not status || verbosity>=10) { 
      os << std::endl << "   tJacobi2x2PreconditionerFactory::test_identity " << toString(status) << ":  (y=inv(A)*x) != z (|y-z|_2/|z|_2 = " 
                      << diff << ")" << std::endl;
      os << "      "; Print(os,"x",x);
      os << "      "; Print(os,"y",y);
   }
   allPassed &= status;

   // test vector [-2 4 7 9]
   ea[0] =-2.0; ea[1] = 4.0; eb[0] = 7.0; eb[1] = 9.0;
   Thyra::apply(*precOp,Thyra::NOTRANS,*x,y.ptr());
   status = ((diff = Teko::Test::Difference(y,x)/Thyra::norm_2(*x->col(0)))<tolerance_);
   if(not status || verbosity>=10) { 
      os << std::endl << "   tJacobi2x2PreconditionerFactory::test_identity " << toString(status) << ":  (y=inv(A)*x) != z (|y-z|_2/|z|_2 = " 
                      << diff << ")" << std::endl;
      os << "      "; Print(os,"x",x);
      os << "      "; Print(os,"y",y);
   }
   allPassed &= status;

   // test vector [1 0 0 -5]
   ea[0] = 1.0; ea[1] = 0.0; eb[0] = 0.0; eb[1] =-5.0;
   Thyra::apply(*precOp,Thyra::NOTRANS,*x,y.ptr());
   status = ((diff = Teko::Test::Difference(y,x)/Thyra::norm_2(*x->col(0)))<tolerance_);
   if(not status || verbosity>=10) { 
      os << std::endl << "   tJacobi2x2PreconditionerFactory::test_identity " << toString(status) << ":  (y=inv(A)*x) != z (|y-z|_2/|z|_2 = " 
                      << diff << ")" << std::endl;
      os << "      "; Print(os,"x",x);
      os << "      "; Print(os,"y",y);
   }
   allPassed &= status;

   // test vector [4 -4 6 12]
   ea[0] = 4.0; ea[1] =-4.0; eb[0] = 6.0; eb[1] =12.0;
   Thyra::apply(*precOp,Thyra::NOTRANS,*x,y.ptr());
   status = ((diff = Teko::Test::Difference(y,x)/Thyra::norm_2(*x->col(0)))<tolerance_);
   if(not status || verbosity>=10) { 
      os << std::endl << "   tJacobi2x2PreconditionerFactory::test_identity " << toString(status) << ":  (y=inv(A)*x) != z (|y-z|_2/|z|_2 = " 
                      << diff << ")" << std::endl;
      os << "      "; Print(os,"x",x);
      os << "      "; Print(os,"y",y);
   }
   allPassed &= status;

   return allPassed;
}