IFPACK  Development
 All Classes Files Functions Variables Enumerations Friends
Ifpack_HIPS.cpp
00001 /*
00002 //@HEADER
00003 // ***********************************************************************
00004 //
00005 //       Ifpack: Object-Oriented Algebraic Preconditioner Package
00006 //                 Copyright (2002) Sandia Corporation
00007 //
00008 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
00009 // license for use of this work by or on behalf of the U.S. Government.
00010 //
00011 // Redistribution and use in source and binary forms, with or without
00012 // modification, are permitted provided that the following conditions are
00013 // met:
00014 //
00015 // 1. Redistributions of source code must retain the above copyright
00016 // notice, this list of conditions and the following disclaimer.
00017 //
00018 // 2. Redistributions in binary form must reproduce the above copyright
00019 // notice, this list of conditions and the following disclaimer in the
00020 // documentation and/or other materials provided with the distribution.
00021 //
00022 // 3. Neither the name of the Corporation nor the names of the
00023 // contributors may be used to endorse or promote products derived from
00024 // this software without specific prior written permission.
00025 //
00026 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
00027 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
00028 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
00029 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
00030 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
00031 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
00032 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
00033 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
00034 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
00035 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
00036 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
00037 //
00038 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
00039 //
00040 // ***********************************************************************
00041 //@HEADER
00042 */
00043 
00044 #include "Ifpack_HIPS.h"
00045 #if defined(HAVE_IFPACK_HIPS) && defined(HAVE_MPI)
00046 
00047 
00048 #include "Ifpack_Utils.h"
00049 extern "C" {
00050 #include "hips.h"
00051 }
00052 //#include "io.h"
00053 #include "Epetra_MpiComm.h"
00054 #include "Epetra_IntVector.h"
00055 #include "Epetra_Import.h"
00056 #include "Teuchos_ParameterList.hpp"
00057 #include "Teuchos_RefCountPtr.hpp"
00058 #ifdef IFPACK_NODE_AWARE_CODE
00059 extern int ML_NODE_ID;
00060 #endif
00061 
00062 using Teuchos::RefCountPtr;
00063 using Teuchos::rcp;
00064 
00065 
00066 
00067 Ifpack_HIPS::Ifpack_HIPS(Epetra_RowMatrix* A):
00068   A_(rcp(A,false)),
00069   HIPS_id(-1), // Assumes user will initialze HIPS outside 
00070   IsParallel_(false),
00071   IsInitialized_(false),
00072   IsComputed_(false),
00073   Label_(),
00074   Time_(A_->Comm())
00075 {
00076   
00077 }
00078 
00079 void Ifpack_HIPS::Destroy(){
00080   //NTS: Assume user will call HIPS_Finalize elsewhere - HIPS_Clean never needs
00081   //to be called if HIPS_Finalize is called at the end, unless you want to reuse
00082   //a slot.
00083 }
00084 
00085 
00086 
00087 int Ifpack_HIPS::Initialize(){
00088   if(Comm().NumProc() != 1) IsParallel_ = true;
00089   else IsParallel_ = false;
00090   IsInitialized_=true;
00091   return 0;
00092 }
00093 
00094 int Ifpack_HIPS::SetParameters(Teuchos::ParameterList& parameterlist){
00095   List_=parameterlist;
00096 
00097   // Grab the hips ID
00098   HIPS_id=List_.get("hips: id",-1);
00099   if(HIPS_id==-1) IFPACK_CHK_ERR(-1);
00100 
00101   // Set Defaults
00102   HIPS_SetDefaultOptions(HIPS_id,List_.get("hips: strategy",HIPS_ITERATIVE));  
00103   
00104   // Set the communicator
00105   const Epetra_MpiComm* MpiComm=dynamic_cast<const Epetra_MpiComm*>(&A_->Comm());
00106   if(!MpiComm) IFPACK_CHK_ERR(-2);  
00107   HIPS_SetCommunicator(HIPS_id,MpiComm->GetMpiComm());
00108   
00109   // Options
00110   HIPS_SetOptionINT(HIPS_id,HIPS_SYMMETRIC,List_.get("hips: symmetric",0));
00111   HIPS_SetOptionINT(HIPS_id,HIPS_VERBOSE,List_.get("hips: setup output",1));
00112   HIPS_SetOptionINT(HIPS_id,HIPS_LOCALLY,List_.get("hips: fill",0));
00113   // Sadly, this fill option doesn't work for HIPS_ITERATIVE mode, meaning the
00114   // only way to control fill-in is via the drop tolerance. 
00115 
00116   HIPS_SetOptionINT(HIPS_id,HIPS_REORDER,List_.get("hips: reorder",1));
00117   HIPS_SetOptionINT(HIPS_id,HIPS_GRAPH_SYM,List_.get("hips: graph symmetric",0));
00118   HIPS_SetOptionINT(HIPS_id,HIPS_FORTRAN_NUMBERING,List_.get("hips: fortran numbering",0));
00119   HIPS_SetOptionINT(HIPS_id,HIPS_DOF,List_.get("hips: dof per node",1));
00120 
00121   // This will disable the GMRES wrap around HIPS.
00122   //  HIPS_SetOptionINT(HIPS_id,HIPS_ITMAX,1);
00123   HIPS_SetOptionINT(HIPS_id,HIPS_ITMAX,-1);  
00124   //  HIPS_SetOptionINT(HIPS_id,HIPS_KRYLOV_METHOD,List_.get("hips: krylov",0));
00125   HIPS_SetOptionINT(HIPS_id,HIPS_KRYLOV_RESTART,-1);
00126   
00127   // Make sure the ILU always runs, by setting the internal tolerance really, really low.
00128   HIPS_SetOptionREAL(HIPS_id,HIPS_PREC,1e-100);
00129 
00130   // Options for Iterative only
00131   if(List_.get("hips: strategy",HIPS_ITERATIVE)==HIPS_ITERATIVE){  
00132     HIPS_SetOptionREAL(HIPS_id, HIPS_DROPTOL0, List_.get("hips: drop tolerance",1e-2));
00133     HIPS_SetOptionREAL(HIPS_id, HIPS_DROPTOL1, List_.get("hips: drop tolerance",1e-2));
00134     HIPS_SetOptionREAL(HIPS_id, HIPS_DROPTOLE, List_.get("hips: drop tolerance",1e-2));    
00135   }
00136   // NTS: This is only a subset of the actual HIPS options. 
00137   return 0;
00138 }
00139 
00140 
00141 int Ifpack_HIPS::Compute(){
00142   if(HIPS_id==-1) IFPACK_CHK_ERR(-1);
00143   int N=A_->NumMyRows(), nnz=A_->NumMyNonzeros();
00144   const Epetra_Comm &Comm=A_->Comm();
00145   int mypid=Comm.MyPID();
00146   
00147   // Pull the column indices, if possible
00148   int *rowptr,*colind,ierr,maxnr,Nr;
00149   double *values;  
00150   Epetra_CrsMatrix *Acrs=dynamic_cast<Epetra_CrsMatrix*>(&*A_);
00151   const Epetra_Map &RowMap=A_->RowMatrixRowMap();
00152   const Epetra_Map &ColMap=A_->RowMatrixColMap();
00153   if(Acrs) Acrs->ExtractCrsDataPointers(rowptr,colind,values);
00154   else{
00155     maxnr=A_->MaxNumEntries();
00156     colind=new int[maxnr];
00157     values=new double[maxnr];
00158   }
00159 
00160   // Create 0-to-N-1 consistent maps for rows and columns
00161   RowMap0_=rcp(new Epetra_Map(-1,N,0,Comm));
00162   Epetra_IntVector RowGIDs(View,RowMap,RowMap0_->MyGlobalElements());
00163   Epetra_IntVector ColGIDs(ColMap);
00164   Epetra_Import RowToCol(ColMap,RowMap);
00165   ColGIDs.Import(RowGIDs,RowToCol,Insert);
00166   ColMap0_=rcp(new Epetra_Map(-1,ColMap.NumMyElements(),ColGIDs.Values(),0,Comm));
00167 
00168   int *gcolind=0;
00169   if(Acrs){
00170     //  Global CIDs
00171     gcolind=new int[nnz];
00172     for(int j=0;j<nnz;j++) gcolind[j]=RowMap0_->GID(colind[j]);        
00173     ierr =  HIPS_GraphDistrCSR(HIPS_id,A_->NumGlobalRows(),A_->NumMyRows(),RowMap0_->MyGlobalElements(),
00174                                rowptr,gcolind);  
00175   }
00176   else{
00177     // Do things the hard way
00178     ierr=HIPS_GraphBegin(HIPS_id,A_->NumGlobalRows(),nnz);
00179     if(ierr!=HIPS_SUCCESS) IFPACK_CHK_ERR(-2);
00180 
00181     // Graph insert - RM mode
00182     for(int i=0;i<N;i++){
00183       A_->ExtractMyRowCopy(i,maxnr,Nr,values,colind);
00184       for(int j=0;j<Nr;j++){
00185         ierr=HIPS_GraphEdge(HIPS_id,RowMap0_->GID(i),ColMap0_->GID(colind[j]));
00186         if(ierr!=HIPS_SUCCESS) IFPACK_CHK_ERR(-3);        
00187       }
00188     }
00189     ierr=HIPS_GraphEnd(HIPS_id);
00190   }      
00191   if(ierr!=HIPS_SUCCESS) IFPACK_CHK_ERR(-4);  
00192   
00193   /*Have processor 0 send in the partition*/
00194   // NTS: This is really, really annoying.  Look at all this import/export
00195   // stuff.  This is mind-numbingly unnecessary.
00196   
00197   Epetra_Map OnePerProcMap(-1,1,0,Comm);
00198   Epetra_IntVector RowsPerProc(OnePerProcMap);
00199   Epetra_IntVector RowGID(View,*RowMap0_,RowMap0_->MyGlobalElements());
00200   
00201   // Get the RPP partial sums
00202   Comm.ScanSum(&N,&(RowsPerProc[0]),1); 
00203 
00204   // Build the maps for xfer to proc 0
00205   int OPP_els=0,RPP_els=0;
00206   if(!mypid){OPP_els=Comm.NumProc(); RPP_els=A_->NumGlobalRows();}
00207   Epetra_Map OPPMap_0(-1,OPP_els,0,Comm);
00208   Epetra_Map RPPMap_0(-1,RPP_els,0,Comm);
00209   Epetra_Import OPP_importer(OPPMap_0,OnePerProcMap);
00210   Epetra_Import RPP_importer(RPPMap_0,*RowMap0_);
00211   
00212   // Pull the vectors to proc 0
00213   Epetra_IntVector OPP_0(OPPMap_0);
00214   Epetra_IntVector RPP_0(RPPMap_0);
00215   OPP_0.Import(RowsPerProc,OPP_importer,Add);
00216   RPP_0.Import(RowGID,RPP_importer,Add);
00217  
00218   // Setup the partition
00219   if(!mypid){    
00220     int *mapptr=0;
00221     mapptr=new int[Comm.NumProc()+1];
00222     mapptr[0]=0;
00223     for(int i=0;i<Comm.NumProc();i++)
00224       mapptr[i+1]=OPP_0[i];
00225 
00226     // Call is only necessary on proc 0
00227     ierr=HIPS_SetPartition(HIPS_id,A_->Comm().NumProc(),mapptr,RPP_0.Values());
00228     HIPS_ExitOnError(ierr);
00229     delete [] mapptr;
00230   }      
00231 
00232   if(Acrs)
00233     ierr = HIPS_MatrixDistrCSR(HIPS_id,A_->NumMyRows(),RowMap0_->MyGlobalElements(),rowptr,gcolind,values,HIPS_ASSEMBLY_OVW,HIPS_ASSEMBLY_OVW,HIPS_ASSEMBLY_FOOL,0);  
00234   else{
00235     // Do things the hard way
00236      ierr = HIPS_AssemblyBegin(HIPS_id, nnz, HIPS_ASSEMBLY_OVW, HIPS_ASSEMBLY_OVW, HIPS_ASSEMBLY_FOOL,0);
00237      if(ierr!=HIPS_SUCCESS){HIPS_PrintError(ierr);IFPACK_CHK_ERR(-5);}
00238 
00239      // Matrix insert - RM Mode
00240      for(int i=0;i<N;i++){
00241        A_->ExtractMyRowCopy(i,maxnr,Nr,values,colind);
00242        for(int j=0;j<Nr;j++){      
00243          ierr = HIPS_AssemblySetValue(HIPS_id,RowMap0_->GID(i),ColMap0_->GID(colind[j]), values[j]);
00244          if(ierr!=HIPS_SUCCESS){HIPS_PrintError(ierr);IFPACK_CHK_ERR(-6);}
00245        }
00246      }
00247      ierr = HIPS_AssemblyEnd(HIPS_id);
00248   }
00249   if(ierr!=HIPS_SUCCESS){HIPS_PrintError(ierr);IFPACK_CHK_ERR(-7);}
00250 
00251   // Force factorization
00252   //NTS: This is odd.  There should be a better way to force this to happen.
00253   double *X=new double[A_->NumMyRows()];
00254   double *Y=new double[A_->NumMyRows()];
00255   for(int i=0;i<A_->NumMyRows();i++) X[i]=1.0;
00256 
00257   // Force HIPS to do it's own Import/Export
00258   ierr=HIPS_SetRHS(HIPS_id,A_->NumMyRows(),RowMap0_->MyGlobalElements(),X,HIPS_ASSEMBLY_OVW,HIPS_ASSEMBLY_OVW,HIPS_ASSEMBLY_FOOL);
00259   if(ierr!=HIPS_SUCCESS) IFPACK_CHK_ERR(-11);
00260   
00261   ierr=HIPS_GetSolution(HIPS_id,A_->NumMyRows(),RowMap0_->MyGlobalElements(),Y,HIPS_ASSEMBLY_FOOL);
00262   if(ierr!=HIPS_SUCCESS) {
00263     HIPS_PrintError(ierr);
00264     IFPACK_CHK_ERR(-12);
00265   }
00266   
00267   // Reset output for iteration
00268   HIPS_SetOptionINT(HIPS_id,HIPS_VERBOSE,List_.get("hips: iteration output",0));
00269 
00270   // Set Label
00271   int nnzP=0;
00272   HIPS_GetInfoINT(HIPS_id,HIPS_INFO_NNZ,&nnzP);
00273   if(nnzP>0) sprintf(Label_,"Ifpack_HIPS [dt=%4.1e fill=%4.2f]",List_.get("hips: drop tolerance",1e-2),(double)nnzP/(double)A_->NumGlobalNonzeros());
00274   else sprintf(Label_,"Ifpack_HIPS [dt=%4.1e]",List_.get("hips: drop tolerance",1e-2));
00275   // NTS: fill requires a HIPS debug level of at least 2
00276   
00277   IsComputed_=true;
00278 
00279   // Cleanup
00280   if(!Acrs){
00281     delete [] colind;
00282     delete [] values;
00283   }
00284   else
00285     delete [] gcolind;
00286   
00287   delete [] X;
00288   delete [] Y;
00289   return 0;
00290 }
00291 
00292 
00293 
00294 
00295 int Ifpack_HIPS::ApplyInverse(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const{
00296   int rv;
00297   if (!IsComputed())
00298     IFPACK_CHK_ERR(-3);
00299 
00300   if (X.NumVectors() != Y.NumVectors())
00301     IFPACK_CHK_ERR(-2);
00302 
00303   // HAQ: For now
00304   if(X.NumVectors()!=1) IFPACK_CHK_ERR(-42);
00305 
00306   // Wrapping for X==Y
00307   Teuchos::RefCountPtr< Epetra_MultiVector > X2;
00308   if (X.Pointers()[0] == Y.Pointers()[0])
00309     X2 = Teuchos::rcp( new Epetra_MultiVector(X) );    
00310   else
00311     X2 = Teuchos::rcp( (Epetra_MultiVector*)&X, false );
00312 
00313   Time_.ResetStartTime();
00314 
00315   // Force HIPS to do it's own Import/Export
00316   rv=HIPS_SetRHS(HIPS_id,A_->NumMyRows(),RowMap0_->MyGlobalElements(),(*X2)[0],HIPS_ASSEMBLY_OVW,HIPS_ASSEMBLY_OVW,HIPS_ASSEMBLY_FOOL);
00317   if(rv!=HIPS_SUCCESS) IFPACK_CHK_ERR(-11);
00318   
00319   rv=HIPS_GetSolution(HIPS_id,A_->NumMyRows(),RowMap0_->MyGlobalElements(),Y[0],HIPS_ASSEMBLY_FOOL);
00320   if(rv!=HIPS_SUCCESS) {
00321     HIPS_PrintError(rv);
00322     IFPACK_CHK_ERR(-12);
00323   }
00324   
00325   ++NumApplyInverse_;
00326   ApplyInverseTime_ += Time_.ElapsedTime();
00327 
00328   return(0);
00329 }
00330 
00331 
00332 ostream& Ifpack_HIPS::Print(ostream& os) const{
00333   os<<"Need to add meaningful output"<<endl;
00334   return os;
00335 }
00336 
00337 
00338 double Ifpack_HIPS::Condest(const Ifpack_CondestType CT, 
00339                              const int MaxIters,
00340                              const double Tol,
00341                              Epetra_RowMatrix* Matrix_in){
00342   return -1.0;
00343 }
00344 
00345 
00346 
00347 
00348 #endif
 All Classes Files Functions Variables Enumerations Friends