|
Teko Version of the Day
|
00001 /* 00002 // @HEADER 00003 // 00004 // *********************************************************************** 00005 // 00006 // Teko: A package for block and physics based preconditioning 00007 // Copyright 2010 Sandia Corporation 00008 // 00009 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation, 00010 // the U.S. Government retains certain rights in this software. 00011 // 00012 // Redistribution and use in source and binary forms, with or without 00013 // modification, are permitted provided that the following conditions are 00014 // met: 00015 // 00016 // 1. Redistributions of source code must retain the above copyright 00017 // notice, this list of conditions and the following disclaimer. 00018 // 00019 // 2. Redistributions in binary form must reproduce the above copyright 00020 // notice, this list of conditions and the following disclaimer in the 00021 // documentation and/or other materials provided with the distribution. 00022 // 00023 // 3. Neither the name of the Corporation nor the names of the 00024 // contributors may be used to endorse or promote products derived from 00025 // this software without specific prior written permission. 00026 // 00027 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY 00028 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE 00029 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR 00030 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE 00031 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, 00032 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, 00033 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR 00034 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF 00035 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING 00036 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS 00037 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE. 00038 // 00039 // Questions? Contact Eric C. Cyr (eccyr@sandia.gov) 00040 // 00041 // *********************************************************************** 00042 // 00043 // @HEADER 00044 00045 */ 00046 00047 #include "Epetra/Teko_StridedEpetraOperator.hpp" 00048 #include "Epetra/Teko_StridedMappingStrategy.hpp" 00049 #include "Epetra/Teko_ReorderedMappingStrategy.hpp" 00050 00051 #include "Teuchos_VerboseObject.hpp" 00052 00053 #include "Thyra_LinearOpBase.hpp" 00054 #include "Thyra_EpetraLinearOp.hpp" 00055 #include "Thyra_EpetraThyraWrappers.hpp" 00056 #include "Thyra_DefaultProductMultiVector.hpp" 00057 #include "Thyra_DefaultProductVectorSpace.hpp" 00058 #include "Thyra_DefaultBlockedLinearOp.hpp" 00059 #include "Thyra_get_Epetra_Operator.hpp" 00060 00061 #include "Epetra_Vector.h" 00062 #include "EpetraExt_MultiVectorOut.h" 00063 #include "EpetraExt_RowMatrixOut.h" 00064 00065 #include "Teko_Utilities.hpp" 00066 00067 namespace Teko { 00068 namespace Epetra { 00069 00070 using Teuchos::RCP; 00071 using Teuchos::rcp; 00072 using Teuchos::rcp_dynamic_cast; 00073 00074 StridedEpetraOperator::StridedEpetraOperator(int numVars,const Teuchos::RCP<Epetra_Operator> & content, 00075 const std::string & label) 00076 : Teko::Epetra::EpetraOperatorWrapper(), label_(label) 00077 { 00078 std::vector<int> vars; 00079 00080 // build vector describing the sub maps 00081 for(int i=0;i<numVars;i++) vars.push_back(1); 00082 00083 SetContent(vars,content); 00084 } 00085 00086 StridedEpetraOperator::StridedEpetraOperator(const std::vector<int> & vars,const Teuchos::RCP<Epetra_Operator> & content, 00087 const std::string & label) 00088 : Teko::Epetra::EpetraOperatorWrapper(), label_(label) 00089 { 00090 SetContent(vars,content); 00091 } 00092 00093 void StridedEpetraOperator::SetContent(const std::vector<int> & vars,const Teuchos::RCP<Epetra_Operator> & content) 00094 { 00095 fullContent_ = content; 00096 stridedMapping_ = rcp(new StridedMappingStrategy(vars,Teuchos::rcpFromRef(fullContent_->OperatorDomainMap()), 00097 fullContent_->Comm())); 00098 SetMapStrategy(stridedMapping_); 00099 00100 // build thyra operator 00101 BuildBlockedOperator(); 00102 } 00103 00104 void StridedEpetraOperator::BuildBlockedOperator() 00105 { 00106 TEUCHOS_ASSERT(stridedMapping_!=Teuchos::null); 00107 00108 // get a CRS matrix 00109 const RCP<const Epetra_CrsMatrix> crsContent = rcp_dynamic_cast<const Epetra_CrsMatrix>(fullContent_); 00110 00111 // ask the strategy to build the Thyra operator for you 00112 if(stridedOperator_==Teuchos::null) { 00113 stridedOperator_ = stridedMapping_->buildBlockedThyraOp(crsContent,label_); 00114 } 00115 else { 00116 const RCP<Thyra::BlockedLinearOpBase<double> > blkOp = rcp_dynamic_cast<Thyra::BlockedLinearOpBase<double> >(stridedOperator_,true); 00117 stridedMapping_->rebuildBlockedThyraOp(crsContent,blkOp); 00118 } 00119 00120 // set whatever is returned 00121 SetOperator(stridedOperator_,false); 00122 00123 // reorder if neccessary 00124 if(reorderManager_!=Teuchos::null) 00125 Reorder(*reorderManager_); 00126 } 00127 00128 const Teuchos::RCP<const Epetra_Operator> StridedEpetraOperator::GetBlock(int i,int j) const 00129 { 00130 const RCP<const Thyra::BlockedLinearOpBase<double> > blkOp 00131 = Teuchos::rcp_dynamic_cast<const Thyra::BlockedLinearOpBase<double> >(getThyraOp()); 00132 00133 return Thyra::get_Epetra_Operator(*blkOp->getBlock(i,j)); 00134 } 00135 00139 void StridedEpetraOperator::Reorder(const BlockReorderManager & brm) 00140 { 00141 reorderManager_ = rcp(new BlockReorderManager(brm)); 00142 00143 // build reordered objects 00144 RCP<const MappingStrategy> reorderMapping = rcp(new ReorderedMappingStrategy(*reorderManager_,stridedMapping_)); 00145 RCP<const Thyra::BlockedLinearOpBase<double> > blockOp 00146 = rcp_dynamic_cast<const Thyra::BlockedLinearOpBase<double> >(stridedOperator_); 00147 00148 RCP<const Thyra::LinearOpBase<double> > A = buildReorderedLinearOp(*reorderManager_,blockOp); 00149 00150 // set them as working values 00151 SetMapStrategy(reorderMapping); 00152 SetOperator(A,false); 00153 } 00154 00156 void StridedEpetraOperator::RemoveReording() 00157 { 00158 SetMapStrategy(stridedMapping_); 00159 SetOperator(stridedOperator_,false); 00160 reorderManager_ = Teuchos::null; 00161 } 00162 00165 void StridedEpetraOperator::WriteBlocks(const std::string & prefix) const 00166 { 00167 RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blockOp 00168 = rcp_dynamic_cast<Thyra::PhysicallyBlockedLinearOpBase<double> >(stridedOperator_); 00169 00170 // get size of strided block operator 00171 int rows = Teko::blockRowCount(blockOp); 00172 00173 for(int i=0;i<rows;i++) { 00174 for(int j=0;j<rows;j++) { 00175 // build the file name 00176 std::stringstream ss; 00177 ss << prefix << "_" << i << j << ".mm"; 00178 00179 // get the row matrix object (Note: can't use "GetBlock" method b/c matrix might be reordered) 00180 RCP<const Epetra_RowMatrix> mat 00181 = Teuchos::rcp_dynamic_cast<const Epetra_RowMatrix>(Thyra::get_Epetra_Operator(*blockOp->getBlock(i,j))); 00182 00183 // write to file 00184 EpetraExt::RowMatrixToMatrixMarketFile(ss.str().c_str(),*mat); 00185 } 00186 } 00187 } 00188 00196 std::string StridedEpetraOperator::PrintNorm(const eNormType & nrmType,const char newline) 00197 { 00198 BlockedLinearOp blockOp = toBlockedLinearOp(stridedOperator_); 00199 00200 // get size of strided block operator 00201 int rowCount = Teko::blockRowCount(blockOp); 00202 int colCount = Teko::blockRowCount(blockOp); 00203 00204 std::stringstream ss; 00205 ss.precision(4); 00206 ss << std::scientific; 00207 for(int row=0;row<rowCount;row++) { 00208 for(int col=0;col<colCount;col++) { 00209 // get the row matrix object (Note: can't use "GetBlock" method b/c matrix might be reordered) 00210 RCP<const Epetra_CrsMatrix> mat = Teuchos::rcp_dynamic_cast<const Epetra_CrsMatrix>( 00211 Thyra::get_Epetra_Operator(*blockOp->getBlock(row,col))); 00212 00213 // compute the norm 00214 double norm = 0.0; 00215 switch(nrmType) { 00216 case Inf: 00217 norm = mat->NormInf(); 00218 break; 00219 case One: 00220 norm = mat->NormOne(); 00221 break; 00222 case Frobenius: 00223 norm = mat->NormFrobenius(); 00224 break; 00225 default: 00226 TEST_FOR_EXCEPTION(true,std::runtime_error, 00227 "StridedEpetraOperator::eNormType incorrectly specified"); 00228 } 00229 00230 ss << norm << " "; 00231 } 00232 ss << newline; 00233 } 00234 00235 return ss.str(); 00236 } 00237 00238 #ifndef Teko_DEBUG_OFF 00239 bool StridedEpetraOperator::testAgainstFullOperator(int count,double tol) const 00240 { 00241 Epetra_Vector xf(OperatorRangeMap()); 00242 Epetra_Vector xs(OperatorRangeMap()); 00243 Epetra_Vector y(OperatorDomainMap()); 00244 00245 // test operator many times 00246 bool result = true; 00247 double diffNorm=0.0,trueNorm=0.0; 00248 for(int i=0;i<count;i++) { 00249 xf.PutScalar(0.0); 00250 xs.PutScalar(0.0); 00251 y.Random(); 00252 00253 // apply operator 00254 Apply(y,xs); // xs = A*y 00255 fullContent_->Apply(y,xf); // xf = A*y 00256 00257 // compute norms 00258 xs.Update(-1.0,xf,1.0); 00259 xs.Norm2(&diffNorm); 00260 xf.Norm2(&trueNorm); 00261 00262 // check result 00263 result &= (diffNorm/trueNorm < tol); 00264 } 00265 00266 return result; 00267 } 00268 #endif 00269 00270 } // end namespace Epetra 00271 } // end namespace Teko
1.7.4