Main Page   Class Hierarchy   Alphabetical List   Compound List   File List   Compound Members   File Members  

BandGenLinLapackSolver.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.1.1.1 $
00022 // $Date: 2000/09/15 08:23:28 $
00023 // $Source: /usr/local/cvs/OpenSees/SRC/system_of_eqn/linearSOE/bandGEN/BandGenLinLapackSolver.cpp,v $
00024                                                                         
00025                                                                         
00026 // File: ~/system_of_eqn/linearSOE/bandGEN/BandGenLinLapackSolver.h
00027 //
00028 // Written: fmk 
00029 // Created: Tue Sep 26 16:27:47: 1996
00030 // Revision: A
00031 //
00032 // Description: This file contains the class definition for 
00033 // BandGenLinLapackSolver. It solves the BandGenLinSOE object by calling
00034 // Lapack routines.
00035 //
00036 // What: "@(#) BandGenLinLapackSolver.h, revA"
00037 
00038 #include <BandGenLinLapackSolver.h>
00039 #include <BandGenLinSOE.h>
00040 #include <f2c.h>
00041 #include <math.h>
00042 
00043 
00044 BandGenLinLapackSolver::BandGenLinLapackSolver()
00045 :BandGenLinSolver(SOLVER_TAGS_BandGenLinLapackSolver),
00046  iPiv(0), iPivSize(0)
00047 {
00048     
00049 }
00050 
00051 BandGenLinLapackSolver::~BandGenLinLapackSolver()
00052 {
00053     if (iPiv != 0)
00054  delete [] iPiv;
00055 }
00056 
00057 #ifdef _WIN32
00058 
00059 extern "C" int _stdcall DGBSV(int *N, int *KL, int *KU, int *NRHS, double *A, 
00060          int *LDA, int *iPiv, double *B, int *LDB, 
00061          int *INFO);
00062 
00063 extern "C" int _stdcall DGBTRS(char *TRANS, unsigned int *sizeC, 
00064           int *N, int *KL, int *KU, int *NRHS,
00065           double *A, int *LDA, int *iPiv, 
00066           double *B, int *LDB, int *INFO);
00067 
00068 #else
00069 
00070 extern "C" int dgbsv_(int *N, int *KL, int *KU, int *NRHS, double *A, 
00071         int *LDA, int *iPiv, double *B, int *LDB, int *INFO);
00072         
00073 
00074 extern "C" int dgbtrs_(char *TRANS, int *N, int *KL, int *KU, int *NRHS, 
00075          double *A, int *LDA, int *iPiv, double *B, int *LDB, 
00076          int *INFO);
00077 #endif
00078 int
00079 BandGenLinLapackSolver::solve(void)
00080 {
00081     if (theSOE == 0) {
00082  cerr << "WARNING BandGenLinLapackSolver::solve(void)- ";
00083  cerr << " No LinearSOE object has been set\n";
00084  return -1;
00085     }
00086 
00087     int n = theSOE->size;    
00088     // check iPiv is large enough
00089     if (iPivSize < n) {
00090  cerr << "WARNING BandGenLinLapackSolver::solve(void)- ";
00091  cerr << " iPiv not large enough - has setSize() been called?\n";
00092  return -1;
00093     }     
00094 
00095     int kl = theSOE->numSubD;
00096     int ku = theSOE->numSuperD;
00097     int ldA = 2*kl + ku +1;
00098     int nrhs = 1;
00099     int ldB = n;
00100     int info;
00101     double *Aptr = theSOE->A;
00102     double *Xptr = theSOE->X;
00103     double *Bptr = theSOE->B;
00104     int    *iPIV = iPiv;
00105     
00106     // first copy B into X
00107     for (int i=0; i<n; i++) {
00108  *(Xptr++) = *(Bptr++);
00109     }
00110     Xptr = theSOE->X;
00111 
00112     // now solve AX = B
00113 
00114 #ifdef _WIN32
00115     {if (theSOE->factored == false)  
00116  // factor and solve 
00117  DGBSV(&n,&kl,&ku,&nrhs,Aptr,&ldA,iPIV,Xptr,&ldB,&info); 
00118     else  {
00119  // solve only using factored matrix
00120  unsigned int sizeC = 1;
00121  DGBTRS("N", &sizeC, &n,&kl,&ku,&nrhs,Aptr,&ldA,iPIV,Xptr,&ldB,&info);
00122     }}
00123 #else
00124     {if (theSOE->factored == false)      
00125  // factor and solve  
00126  dgbsv_(&n,&kl,&ku,&nrhs,Aptr,&ldA,iPIV,Xptr,&ldB,&info);
00127     else
00128  // solve only using factored matrix 
00129  dgbtrs_("N",&n,&kl,&ku,&nrhs,Aptr,&ldA,iPIV,Xptr,&ldB,&info);
00130     }
00131 #endif
00132     // check if successfull
00133     if (info != 0) {
00134  cerr << "WARNING BandGenLinLapackSolver::solve() -";
00135  cerr << "LAPACK routine returned " << info << endl;
00136  return -info;
00137     }
00138 
00139     theSOE->factored = true;
00140     return 0;
00141 }
00142     
00143 
00144 
00145 int
00146 BandGenLinLapackSolver::setSize()
00147 {
00148     // if iPiv not big enough, free it and get one large enough
00149     if (iPivSize < theSOE->size) {
00150  if (iPiv != 0)
00151      delete [] iPiv;
00152  
00153  iPiv = new int[theSOE->size];
00154  if (iPiv == 0) {
00155      cerr << "WARNING BandGenLinLapackSolver::setSize() ";
00156      cerr << " - ran out of memory for iPiv of size ";
00157      cerr << theSOE->size << endl;
00158      return -1;
00159  } else
00160      iPivSize = theSOE->size;
00161     }
00162  
00163     return 0;
00164 }
00165 
00166 int    
00167 BandGenLinLapackSolver::sendSelf(int commitTag, Channel &theChannel)
00168 {
00169     return 0;
00170 }
00171 
00172 int
00173 BandGenLinLapackSolver::recvSelf(int commitTag,
00174      Channel &theChannel, 
00175      FEM_ObjectBroker &theBroker)
00176 {
00177     // nothing to do
00178     return 0;
00179 }
Copyright Contact Us