PetscSparseSeqSolver.cpp

Go to the documentation of this file.
00001 /* ****************************************************************** **
00002 **    OpenSees - Open System for Earthquake Engineering Simulation    **
00003 **          Pacific Earthquake Engineering Research Center            **
00004 **                                                                    **
00005 **                                                                    **
00006 ** (C) Copyright 1999, The Regents of the University of California    **
00007 ** All Rights Reserved.                                               **
00008 **                                                                    **
00009 ** Commercial use of this program without express permission of the   **
00010 ** University of California, Berkeley, is strictly prohibited.  See   **
00011 ** file 'COPYRIGHT'  in main directory for information on usage and   **
00012 ** redistribution,  and for a DISCLAIMER OF ALL WARRANTIES.           **
00013 **                                                                    **
00014 ** Developed by:                                                      **
00015 **   Frank McKenna (fmckenna@ce.berkeley.edu)                         **
00016 **   Gregory L. Fenves (fenves@ce.berkeley.edu)                       **
00017 **   Filip C. Filippou (filippou@ce.berkeley.edu)                     **
00018 **                                                                    **
00019 ** ****************************************************************** */
00020                                                                         
00021 // $Revision: 1.2 $
00022 // $Date: 2005/05/18 19:26:59 $
00023 // $Source: /usr/local/cvs/OpenSees/SRC/system_of_eqn/linearSOE/petsc/PetscSparseSeqSolver.cpp,v $
00024                                                                         
00025 // Written: fmk 
00026 // Created: 04/05
00027                                                                         
00028 // What: "@(#) PetscSparseSeqSolver.h, revA"
00029 
00030 #include <PetscSparseSeqSolver.h>
00031 #include <SparseGenRowLinSOE.h>
00032 #include <f2c.h>
00033 #include <math.h>
00034 #include <Channel.h>
00035 #include <FEM_ObjectBroker.h>
00036 #include <Timer.h>
00037 
00038 PetscSparseSeqSolver::PetscSparseSeqSolver(KSPType meth, PCType pre)
00039   :SparseGenRowLinSolver(SOLVER_TAGS_PetscSparseSeqSolver), 
00040    rTol(PETSC_DEFAULT), aTol(PETSC_DEFAULT), dTol(PETSC_DEFAULT), maxIts(PETSC_DEFAULT)
00041 {
00042   PetscInitialize(0, PETSC_NULL, (char *)0, PETSC_NULL);
00043 
00044   method = meth;
00045   preconditioner = pre;
00046 }
00047 
00048 PetscSparseSeqSolver::PetscSparseSeqSolver(KSPType meth, PCType pre, double relTol, double absTol, double divTol, int maxIterations)
00049   :SparseGenRowLinSolver(SOLVER_TAGS_PetscSparseSeqSolver), 
00050    rTol(relTol), aTol(absTol), dTol(divTol), maxIts(maxIterations)
00051 {
00052   PetscInitialize(0, PETSC_NULL, (char *)0, PETSC_NULL);
00053 
00054   method = meth;
00055   preconditioner = pre;
00056 }
00057 
00058 
00059 PetscSparseSeqSolver::~PetscSparseSeqSolver()
00060 {
00061   MatDestroy(A);
00062   VecDestroy(x);
00063   VecDestroy(b);
00064   KSPDestroy(ksp);
00065 }
00066 
00067 
00068 
00069 int
00070 PetscSparseSeqSolver::solve(void)
00071 {
00072   PetscErrorCode ierr = KSPSolve(ksp, b, x); CHKERRQ(ierr); 
00073 
00074   return ierr;
00075 }
00076 
00077 
00078 int
00079 PetscSparseSeqSolver::getNumIterations(void)
00080 {
00081   PetscErrorCode ierr = KSPGetIterationNumber(ksp, &its);
00082     
00083   return its;
00084 }
00085 
00086 
00087 int
00088 PetscSparseSeqSolver::setSize()
00089 {
00090   int n = theSOE->size;
00091   int nnz = theSOE->nnz;
00092   double *Adata = theSOE->A;
00093   double *Xdata = theSOE->X;
00094   double *Bdata = theSOE->B;
00095   int *colA = theSOE->colA;
00096   int *rowStartA = theSOE->rowStartA;
00097 
00098 
00099   /* 
00100    * Call Petsc VecCreate & MatCreate; NOTE: using previously allocated storage
00101    * 
00102    */
00103 
00104   PetscTruth flg;
00105   PetscErrorCode ierr = PetscOptionsGetInt(PETSC_NULL,"-n", &n, &flg); CHKERRQ(ierr);
00106   ierr = VecCreateSeqWithArray(PETSC_COMM_WORLD, n, Xdata, &x); CHKERRQ(ierr); 
00107   ierr = VecCreateSeqWithArray(PETSC_COMM_WORLD, n, Bdata, &b); CHKERRQ(ierr); 
00108   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_WORLD, n, n, rowStartA, colA,  Adata , &A);
00109 
00110   /* 
00111    * Create linear solver context
00112    */
00113 
00114    KSPCreate(PETSC_COMM_WORLD, &ksp);
00115 
00116    /* 
00117     *  Set operators. NOTE: matrix that defines the linear system
00118     *  also serves as the preconditioning matrix.
00119     */
00120 
00121    KSPSetOperators(ksp, A, A, DIFFERENT_NONZERO_PATTERN);
00122 
00123    /* 
00124     *  Set solution scheme & tolerances
00125     */
00126 
00127    ierr = KSPSetType(ksp, method); CHKERRQ(ierr); 
00128    ierr = KSPSetTolerances(ksp, rTol, aTol, dTol, maxIts); 
00129 
00130    /* 
00131     *  Set preconditioning scheme
00132     */
00133 
00134    KSPGetPC(ksp, &pc);
00135    ierr = PCSetType(pc,  preconditioner); CHKERRQ(ierr); 
00136 
00137    /* 
00138     *  Finally mark so that uses last solution as initial guess in each solution
00139     *    NOTE: maybe change this as another user supplied option
00140     */
00141 
00142    ierr = KSPSetInitialGuessNonzero(ksp, PETSC_TRUE); CHKERRQ(ierr); 
00143 
00144    return ierr;
00145 }
00146 
00147 
00148 int 
00149 PetscSparseSeqSolver::setLinearSOE(SparseGenRowLinSOE &theSystem)
00150 {
00151   theSOE = &theSystem;
00152   return 0;
00153 }
00154 
00155 
00156 int
00157 PetscSparseSeqSolver::sendSelf(int cTag, Channel &theChannel)
00158 {
00159     // nothing to do
00160     return 0;
00161 }
00162 
00163 int
00164 PetscSparseSeqSolver::recvSelf(int cTag, Channel &theChannel, 
00165                       FEM_ObjectBroker &theBroker)
00166 {
00167     // nothing to do
00168     return 0;
00169 }
00170 
00171 
00172 

Generated on Mon Oct 23 15:05:29 2006 for OpenSees by doxygen 1.5.0