Main Page | Class Hierarchy | Class List | File List | Class Members

Prometheus_LinSysCore.h

00001 /* 00002 This is the Prometheus implementation of LinearSystemCore. 00003 $Id: Prometheus_LinSysCore.h,v 1.18 2004/05/07 20:47:41 adams Exp $ 00004 */ 00005 00006 #ifndef _Prometheus_LinSysCore_h_ 00007 #define _Prometheus_LinSysCore_h_ 00008 00009 #include <mpi.h> 00010 00011 #if !defined(PROM_NO_FEI) 00012 00013 /******************************** Sierra ********************************/ 00014 #include <fei_defs.h> 00015 #include <base/Data.h> 00016 #include <base/LinSysCore_flexible.h> 00017 00018 #else 00019 00020 /******************************* Olympus ********************************/ 00021 #include <stdio.h> 00022 class PromOptions; 00023 class PromPerfMonitor; 00024 class PromComm; 00025 typedef int GlobalID; 00026 typedef struct{ 00027 void setTypeName( char matrix_name[] ){} 00028 void setDataPtr( void * p ){ } 00029 void *getDataPtr()const{ return NULL; } 00030 char *getTypeName() const {return "K";} 00031 }Data; 00032 typedef int Lookup; 00033 class LinearSystemCore 00034 { 00035 public: 00036 LinearSystemCore(){} 00037 LinearSystemCore( MPI_Comm ){} 00038 virtual int parameters( int numParams, char** params ) = 0; 00039 virtual int setGlobalOffsets( int len, int* nodeOffsets, int* eqnOffsets, 00040 int* blkEqnOffsets) = 0; 00041 virtual int setConnectivities( GlobalID elemBlock, 00042 int numElements, 00043 int numNodesPerElem, 00044 const GlobalID* elemIDs, 00045 const int* const* connNodes) = 0; 00046 virtual int setStiffnessMatrices(GlobalID elemBlock, 00047 int numElems, 00048 const GlobalID* elemIDs, 00049 const double *const *const *stiff, 00050 int numEqnsPerElem, 00051 const int *const * eqnIndices) = 0; 00052 virtual int setLoadVectors( GlobalID elemBlock, 00053 int numElems, 00054 const GlobalID* elemIDs, 00055 const double *const * load, 00056 int numEqnsPerElem, 00057 const int *const * eqnIndices) = 0; 00058 virtual int setMultCREqns( int multCRSetID, 00059 int numCRs, int numNodesPerCR, 00060 int** nodeNumbers, int** eqnNumbers, 00061 int* fieldIDs, 00062 int* multiplierEqnNumbers) = 0; 00063 00064 virtual int setMatrixStructure( int** ptColIndices, int* ptRowLengths, 00065 int** blkColIndices, int* blkRowLengths, 00066 int* ptRowsPerBlkRow) = 0; 00067 /*int resetMatrixAndVector: 00068 don't destroy the structure of the matrix, but set the value 's' 00069 throughout the matrix and vectors. 00070 */ 00071 virtual int resetMatrixAndVector(double s) = 0; 00072 virtual int resetMatrix(double s) = 0; 00073 virtual int resetRHSVector(double s) = 0; 00074 00075 /*int sumIntoSystemMatrix: 00076 this is the primary assembly function. The coefficients 'values' 00077 are to be accumumlated into (added to any values already in place) 00078 global (0-based) equation 'row' of the matrix. 00079 */ 00080 virtual int sumIntoSystemMatrix( int numPtRows, const int* ptRows, 00081 int numPtCols, const int* ptCols, 00082 int numBlkRows, const int* blkRows, 00083 int numBlkCols, const int* blkCols, 00084 const double* const* values) = 0; 00085 00086 virtual int sumIntoSystemMatrix( int numPtRows, const int* ptRows, 00087 int numPtCols, const int* ptCols, 00088 const double* const* values) = 0; 00089 /*int sumIntoRHSVector: 00090 this is the rhs vector equivalent to sumIntoSystemMatrix above. 00091 */ 00092 virtual int sumIntoRHSVector(int num,const double*values,const int*indices)=0; 00093 virtual int putIntoRHSVector(int num,const double*values,const int*indices)=0; 00094 /* int matrixLoadComplete: 00095 do any internal synchronization/communication. 00096 */ 00097 virtual int matrixLoadComplete() = 0; 00098 00099 virtual int putNodalFieldData( int fieldID, int fieldSize, int* nodeNumbers, 00100 int numNodes, const double* data) = 0; 00101 /* functions for enforcing boundary conditions. */ 00102 virtual int enforceEssentialBC( int* globalEqn, 00103 double* alpha, 00104 double* gamma, int len) = 0; 00105 00106 virtual int enforceRemoteEssBCs( int numEqns, int* globalEqns, 00107 int** colIndices, int* colIndLen, 00108 double** coefs) = 0; 00109 00110 virtual int enforceOtherBC( int* globalEqn, double* alpha, 00111 double* beta, double* gamma, 00112 int len) = 0; 00113 00114 /* functions for managing multiple rhs vectors */ 00115 virtual int setNumRHSVectors(int numRHSs, const int* rhsIDs) = 0; 00116 00117 /* virtual int setRHSID: 00118 set the 'current' rhs context, assuming multiple rhs vectors. */ 00119 virtual int setRHSID(int rhsID) = 0; 00120 00121 /*int putInitialGuess: 00122 function for setting (a subset of) the initial-guess 00123 solution values (i.e., in the 'x' vector). 00124 */ 00125 virtual int putInitialGuess(const int* eqnNumbers, const double* values, int len) = 0; 00126 00127 /* function for getting all of the answers ('x' vector). */ 00128 virtual int getSolution(double* answers, int len) = 0; 00129 00130 /* function for getting the (single) entry at equation 00131 number 'eqnNumber'. */ 00132 virtual int getSolnEntry(int eqnNumber, double& answer) = 0; 00133 00134 /* function for launching the linear solver */ 00135 virtual int launchSolver( int &solveStatus, int &iterations) = 0; 00136 }; 00137 class LinSysCore_flexible : public LinearSystemCore 00138 { 00139 // LSC_flexible 00140 int resetConstraints( double s ); /* LSC_flexible */ 00141 int constraintsLoadComplete(){ return -1; } /* (old) LSC_flexible */ 00142 int setMultCRComplete(); /* LSC_flexible */ 00143 }; 00144 #endif 00145 00146 #define PROM_GET_VALUES -3 00147 #define PEN_INIT 1.e100 00148 00149 class PromCRVector; 00150 class PromCRMatrix; 00151 class Prometheus_EssBCData; 00152 class Prometheus_OthBCData; 00153 class Prometheus_ZeroEquations; 00154 class Prometheus_Constr; 00155 class Prometheus; 00156 class PromVector; 00157 class PromMatrix; 00158 class PromGrid; 00159 class PromMap; 00160 template <class T> class PromList; 00161 template <class T> class PromTable; 00162 00163 /************************************************************************* 00164 * class PromCRCache 00165 */ 00166 class PromCRCache 00167 { 00168 public: 00169 PromCRCache( int len, int geq, int* gids, double *vals ); 00170 ~PromCRCache(); 00171 int len_; 00172 int geq_; 00173 int *geqs_; 00174 double *vals_; 00175 PromCRCache *next_; 00176 }; 00177 00178 /************************************************************************* 00179 * class Prometheus_EssBCData 00180 */ 00181 class Prometheus_EssBCData 00182 { 00183 public: 00184 Prometheus_EssBCData(int* globalEqn, double* alpha, double* gamma, int len); 00185 ~Prometheus_EssBCData(); 00186 00187 int len_; 00188 int* globalEqn_; 00189 double* alpha_; 00190 double* gamma_; 00191 Prometheus_EssBCData *next_; 00192 }; 00193 00194 /************************************************************************* 00195 * class Prometheus_OthBCData 00196 */ 00197 class Prometheus_OthBCData 00198 { 00199 public: 00200 Prometheus_OthBCData( int* globalEqn, double* alpha, double* beta, 00201 double* gamma, int len); 00202 ~Prometheus_OthBCData(); 00203 00204 int len_; 00205 int* globalEqn_; 00206 double* alpha_; 00207 double* beta_; 00208 double* gamma_; 00209 Prometheus_OthBCData *next_; 00210 }; 00211 00212 /************************************************************************* 00213 * class Prometheus_ZeroEquations 00214 */ 00215 class Prometheus_ZeroEquations 00216 { 00217 public: 00218 Prometheus_ZeroEquations(int* localrows, int* sizes, int** idsarr, int len) : 00219 len_(len), zeroRowsLeq_(localrows), sizes_(sizes), idsarr_(idsarr), 00220 next_(NULL) {} 00221 ~Prometheus_ZeroEquations(){} 00222 int len_; 00223 int *sizes_; 00224 int *zeroRowsLeq_; 00225 int **idsarr_; 00226 Prometheus_ZeroEquations *next_; 00227 }; 00228 00229 class PromCR; 00230 class PromSolver; 00232 00235 class PromUzawaSolver 00236 { 00237 public: 00238 PromUzawaSolver() : pwork_(NULL), uwork_(NULL), uzawaC_(NULL),cr_tol_(0.0), 00239 Pen_stiff_(NULL), diag_(NULL), penalty_factor_(5.0), 00240 max_uzawa_its_(10) { } 00241 ~PromUzawaSolver(); 00242 // 00243 int Init( const PromMap &primmap, const PromMap &dualmap ); 00244 int Reset( const PromMap &dualmap ); 00245 int SetReg( PromMatrix *AA ); 00246 static int setRegularization_private( PromMatrix *const A, 00247 const PromMatrix *const C, 00248 PromVector *const diag, 00249 const double alpha, 00250 PromVector *const dd = NULL ); 00251 int SolveUzawa( PromSolver *const solver, const PromMatrix *const AA, 00252 PromCRVector *const RR, PromCRVector *const XX, 00253 const PromCRVector *const BB, PromVector *const x_0, 00254 const int verb, float normR_0, int &ret_out, 00255 int &its_out ) const; 00256 int setDummyReg( PromMatrix *const AA, 00257 const PromTable<PromCR*> *const CR_id ) const; 00258 // data 00259 PromVector *pwork_; 00260 PromVector *uwork_; 00261 PromMatrix *uzawaC_; 00262 PromVector *diag_; 00263 PromVector *Pen_stiff_; 00264 double penalty_factor_; 00265 double cr_tol_; 00266 short int max_uzawa_its_; 00267 }; 00268 00269 template <class T> class PromAList; 00270 class PromPC; 00271 class PromComm; 00272 class PromOptions; 00273 class PromPerfMonitor; 00275 00278 class Prometheus_LinSysCore: public LinSysCore_flexible, public PromUzawaSolver 00279 { 00280 public: 00281 Prometheus_LinSysCore( MPI_Comm comm ); 00282 int Init( const PromComm &comm, PromOptions &opts, 00283 PromPerfMonitor &perf ); 00284 virtual ~Prometheus_LinSysCore(); 00285 /*for creating another one, without knowing the run-time type 00286 of 'this' one. */ 00287 LinearSystemCore* clone(); 00288 00289 /* int parameters: 00290 for setting generic argc/argv style parameters. 00291 */ 00292 int parameters( int numParams, char** params ); 00293 00294 int setLookup( Lookup& lookup ); 00295 00296 int setGlobalOffsets( int len, int* nodeOffsets, int* eqnOffsets, 00297 int* blkEqnOffsets); 00298 00299 int setConnectivities( GlobalID elemBlock, 00300 int numElements, 00301 int numNodesPerElem, 00302 const GlobalID* elemIDs, 00303 const int* const* connNodes) ; 00304 00305 int setStiffnessMatrices(GlobalID elemBlock, 00306 int numElems, 00307 const GlobalID* elemIDs, 00308 const double *const *const *stiff, 00309 int numEqnsPerElem, 00310 const int *const * eqnIndices); 00311 00312 int setLoadVectors( GlobalID elemBlock, 00313 int numElems, 00314 const GlobalID* elemIDs, 00315 const double *const * load, 00316 int numEqnsPerElem, 00317 const int *const * eqnIndices); 00318 00319 int setMultCREqns( int multCRSetID, 00320 int numCRs, int numNodesPerCR, 00321 int** nodeNumbers, int** eqnNumbers, 00322 int* fieldIDs, 00323 int* multiplierEqnNumbers); 00324 00325 int setPenCREqns( int penCRSetID, 00326 int numCRs, int numNodesPerCR, 00327 int** nodeNumbers, int** eqnNumbers, 00328 int* fieldIDs); 00329 00330 int setMatrixStructure( int** ptColIndices, int* ptRowLengths, 00331 int** blkColIndices, int* blkRowLengths, 00332 int* ptRowsPerBlkRow) ; 00333 00334 /*int resetMatrixAndVector: 00335 don't destroy the structure of the matrix, but set the value 's' 00336 throughout the matrix and vectors. 00337 */ 00338 int resetMatrixAndVector(double s); 00339 int resetMatrix(double s); 00340 int resetRHSVector(double s); 00341 00342 /*int sumIntoSystemMatrix: 00343 this is the primary assembly function. The coefficients 'values' 00344 are to be accumumlated into (added to any values already in place) 00345 global (0-based) equation 'row' of the matrix. 00346 */ 00347 int sumIntoSystemMatrix( int numPtRows, const int* ptRows, 00348 int numPtCols, const int* ptCols, 00349 int numBlkRows, const int* blkRows, 00350 int numBlkCols, const int* blkCols, 00351 const double* const* values); 00352 00353 int sumIntoSystemMatrix( int numPtRows, const int* ptRows, 00354 int numPtCols, const int* ptCols, 00355 const double* const* values); 00356 00357 int putIntoSystemMatrix( int numPtRows, const int* ptRows, 00358 int numPtCols, const int* ptCols, 00359 const double* const* values); 00360 00361 int sumIntoSystemMatrix_private(int numPtRows, const int* ptRows, 00362 int numPtCols, const int* ptCols, 00363 const double* const* values, 00364 const int add_type, const bool image); 00365 int getMatrixRowLength( int row, int& length); 00366 00367 int getMatrixRow( int row, double* coefs, int* indices, 00368 int len, int& rowLength); 00369 00370 /*int sumIntoRHSVector: 00371 this is the rhs vector equivalent to sumIntoSystemMatrix above. 00372 */ 00373 int sumIntoRHSVector(int num, const double* values, const int* indices); 00374 int putIntoRHSVector(int num, const double* values, const int* indices); 00375 int sumIntoRHSVector_private( const int num, const double* values, 00376 const int* indices, const int add_type, 00377 const bool image ); 00378 int accessRHSVector_private( const int num, double *values, 00379 int *indices, const int add_type, 00380 const bool image ); 00381 00382 int getFromRHSVector( int num, double* values, const int* indices); 00383 /* int matrixLoadComplete: 00384 do any internal synchronization/communication. 00385 */ 00386 int matrixLoadComplete(); 00387 00388 int putNodalFieldData( int fieldID, int fieldSize, int* nodeNumbers, 00389 int numNodes, const double* data); 00390 private: 00391 int putCoords_private( int fieldID, int fieldSize, int* nodeNumbers, 00392 int numNodes, const double* data); 00393 public: 00394 /* functions for enforcing boundary conditions. */ 00395 int enforceEssentialBC( int* globalEqn, 00396 double* alpha, 00397 double* gamma, int len); 00398 00399 int enforceRemoteEssBCs( int numEqns, int* globalEqns, 00400 int** colIndices, int* colIndLen, 00401 double** coefs); 00402 00403 int enforceOtherBC( int* globalEqn, double* alpha, 00404 double* beta, double* gamma, 00405 int len); 00406 00407 /* functions for getting/setting matrix or vector pointers. */ 00408 00409 /*getMatrixPtr: 00410 obtain a pointer to the 'A' matrix. This should be considered a 00411 constant pointer -- i.e., this class remains responsible for the 00412 matrix (e.g., de-allocation upon destruction). 00413 */ 00414 int getMatrixPtr(Data& data); 00415 00416 /* copyInMatrix: 00417 replaces the internal matrix with a copy of the input argument, scaled 00418 by the coefficient 'scalar'. 00419 */ 00420 int copyInMatrix(double scalar, const Data& data); 00421 00422 /* copyOutMatrix: 00423 passes out a copy of the internal matrix, scaled by the coefficient 00424 'scalar'. 00425 */ 00426 int copyOutMatrix(double scalar, Data& data); 00427 00428 /* sumInMatrix: 00429 accumulate (sum) a copy of the input argument into the internal 00430 matrix, scaling the input by the coefficient 'scalar'. 00431 */ 00432 int sumInMatrix(double scalar, const Data& data); 00433 00434 /* get/setRHSVectorPtr: 00435 the same semantics apply here as for the matrixPtr functions above. 00436 */ 00437 int getRHSVectorPtr(Data& data); 00438 00439 /* copyInRHSVector/copyOutRHSVector/sumInRHSVector: 00440 the same semantics apply here as for the matrix functions above. 00441 */ 00442 int copyInRHSVector(double scalar, const Data& data); 00443 int copyOutRHSVector(double scalar, Data& data); 00444 int sumInRHSVector(double scalar, const Data& data); 00445 00446 /* destroyMatrixData/destroyVectorData: 00447 Utility function for destroying the matrix (or vector) in Data 00448 */ 00449 int destroyMatrixData(Data& data); 00450 int destroyVectorData(Data& data); 00451 00452 /* functions for managing multiple rhs vectors */ 00453 int setNumRHSVectors(int numRHSs, const int* rhsIDs); 00454 00455 /* int setRHSID: 00456 set the 'current' rhs context, assuming multiple rhs vectors. */ 00457 int setRHSID(int rhsID); 00458 00459 /*int putInitialGuess: 00460 function for setting (a subset of) the initial-guess 00461 solution values (i.e., in the 'x' vector). 00462 */ 00463 int putInitialGuess(const int* eqnNumbers, const double* values, int len); 00464 00465 /* function for getting all of the answers ('x' vector). */ 00466 int getSolution(double* answers, int len); 00467 00468 /* function for getting the (single) entry at equation 00469 number 'eqnNumber'. */ 00470 int getSolnEntry(int eqnNumber, double& answer); 00471 00472 /* function for obtaining the residual vector (using current solution or 00473 initial guess) */ 00474 int formResidual(double* values, int len); 00475 00476 /* function for launching the linear solver */ 00477 int launchSolver(int& solveStatus, int& iterations); 00478 00479 int writeSystem(const char* name); 00480 00481 // LSC_flexible 00482 int resetConstraints( double s ); /* LSC_flexible */ 00483 int constraintsLoadComplete(){ return -1; } /* (old) LSC_flexible */ 00484 int setMultCRComplete(); /* LSC_flexible */ 00485 00486 /* all-at-once MG support */ 00487 int CreateKKTObs( PromGrid *const grid ); 00488 private: 00489 int AggregateCRs( PromList<PromAList<int>*> &aggs, 00490 const PromMatrix*const Ti, const double tol2 = 0.2, 00491 const bool do_glob = (1==1), const bool set_with_local = (1==0) ) const; 00492 int MakeCRProlMap( PromList<PromAList<int>*> &aggs, 00493 const PromGrid * const grid, 00494 PromMap **clm_out, PromMatrix **prol_D_out) const; 00495 int MoveAggs(const PromMap &map, 00496 PromList<PromAList<int>*> &aggs, 00497 const int nextFactor, const PromMap*const CP_colMap, 00498 const PromMatrix*const CP )const; 00499 int SetBCforCRTable(const PromMatrix*const CC, const PromMatrix*const projC); 00500 int GetBCSends(); 00501 static int Mult( const PromMatrix * const AA, const PromMatrix *const BB, 00502 PromMatrix **C_out ); 00503 static int TransMult(const PromMatrix*const AA, const PromMatrix*const BB, 00504 const PromMap *const mp, PromMatrix**, 00505 const PromMap *const mpt=NULL, PromMatrix**t=NULL); 00506 int AggregateCRs2( PromList<PromAList<int>*> &aggs, 00507 const PromMatrix *const Ti, const double tol1, 00508 const double tol2, const bool do_glob, 00509 const bool set_with_local = (1==0) ) const; 00510 int MIS_CRs(int &ndone, const int nadj[], 00511 PromList<PromAList<int>*> &aggs, const PromMatrix *const Ti, 00512 const double tol2, const bool do_glob, 00513 PromAList<int> **lid_aggList, 00514 PromTable<PromAList<int>*> &ghost_list, 00515 PromTable<PromAList<int>*> &listID_list, 00516 int &locAggID, PromAList<int> *const dummy) const; 00517 int GlobalizeAggs(const PromMap &map, PromList<PromAList<int>*> &aggs, 00518 PromAList<int> **lid_aggList, 00519 PromTable<PromAList<int>*> &ghost_list, 00520 PromTable<PromAList<int>*> &listID_list, 00521 PromAList<int> *const dummy) const; 00522 public: 00523 int RemoveCRfromAllGrids(); 00524 PromVector *getRHS(); 00525 PromVector *getSolution(); 00526 int SetRTol( float tol ); 00527 int setMatrixStructure( const PromTable<int> &gid_ghostLid, 00528 const int ndf ); // Olympus interface 00529 private: 00530 int setMatrixStructure_private( const PromTable<int> &gid_ghostLid, 00531 const int ndf ); 00532 /* mapping for Langrange Multiplier problems */ 00533 inline int eq_map_private(const int geq,int *myeq,bool *prim,int *proc)const; 00534 inline int id_map_private(const int gid,int *myid,bool *prim,int *proc)const; 00535 int primGEq_lid_private( const int geq, int *lid_out ) const; 00536 int id_proc_private( const int gid, int &proc) const; 00537 int collectBCs(); 00538 int assembleCRs( PromTable<int> *projCREqs ); 00539 int applyBCs( PromMatrix *projC ); 00540 int enforceEssentialBC_private( int* globalEqn, double* alpha, double* gamma, 00541 int len, PromVector *bmod ); 00542 int enforceOtherBC_private( int* globalEqn, double* alpha, 00543 double* beta, double* gamma, int len); 00544 int constructPrometheus( const PromTable<int> &gid_ghostLid, 00545 const int ndf ); 00546 int ExtendElemList(); 00547 int addextraGhostsParents( PromGrid *lastg, PromGrid *nextg ); 00548 int ReFactorSolver(); 00549 int SymbolicSetUp(); 00550 int MakeSolver( const PromMap *const p_map ); 00551 int FactorKKT( PromCRMatrix *KK, PromGrid *grid, 00552 PromPC *pc, const PromCRVector *const xx); 00553 int LSC_Residual( PromMatrix *KK, const PromVector *const BB, 00554 PromVector *const XX, PromVector *const YY ); 00555 public: 00556 int SetGlobalResTol( double tol ){ normRTol_ = tol; return 0; } 00557 private: 00558 /* data */ 00559 PromMatrix *C_; // cache for KKTMG 'C' or projecton 'C' 00560 int *Ctnnz_; // cache for KKTMG 'Ct' 00561 PromMap *lmap_; // cache for uzawa and KKTMG constraints 00562 PromMap *proj_lmap_;// cache for projection constraints 00563 PromVector **bVecArr_; 00564 PromVector *xVec_; 00565 PromVector *work_; 00566 PromVector *init_guess_; 00567 PromVector *gamma_alpha_; 00568 void *zeroIS_; 00569 int *rhsIDs_; 00570 int *proc_gnode_; 00571 int *proc_globalEq_; 00572 int *proc_lmEq_; 00573 int *proc_primEq_; 00574 int *lid_primLEq_; 00575 PromList<int*> *elemList_; 00576 PromTable<double> *geq_fixity_; 00577 PromTable<PromTable<int>*> *proc_ghostCRIDs_; 00578 PromTable<PromCR*> *CR_id_; 00579 PromTable<PromCR*> *projCR_id_; 00580 00581 Prometheus_EssBCData *essbcdata_; 00582 Prometheus_OthBCData *othbcdata_; 00583 Prometheus_ZeroEquations *zeroEqs_; 00584 PromCRCache *CRCache_; 00585 PromCRCache *projCRCache_; 00586 public: 00587 Prometheus *prom_; 00588 private: 00589 PromComm *Comm_; 00590 PromOptions *options_; 00591 PromPerfMonitor *perf_mon_; 00592 float normRTol_; 00593 char root_[128]; 00594 int numGlobalRows_; 00595 int numGlobalRowsLm_; 00596 int numGlobalRowsPrim_; 00597 short int currentRHS_; 00598 short int dirty_guess_; /* keep track of initial solution guess in "work_"*/ 00599 short int mg_levels_; 00600 short int numRHSs_; 00601 short int verbose_; 00602 short int mype_; 00603 short int npe_; 00604 short int num_matrices_; 00605 short int projNdLimit_; 00606 short int nSolves_; 00607 bool use_uzawa_; 00608 bool matrix_setup_reuse_factor_; 00609 bool rhs_res_scale_; 00610 public: 00611 bool use_complex_; 00612 static const char version_no_[32]; 00613 MPI_Comm MPIComm_; 00614 // complex methods 00615 int setUseComplex( bool b = true ); 00616 int sumIntoSystemMatrix_i( int numPtRows, const int* ptRows, 00617 int numPtCols, const int* ptCols, 00618 int numBlkRows, const int* blkRows, 00619 int numBlkCols, const int* blkCols, 00620 const double* const* values); 00621 int sumIntoSystemMatrix_i( int numPtRows, const int* ptRows, 00622 int numPtCols, const int* ptCols, 00623 const double* const* values); 00624 int sumIntoRHSVector_i( int num, const double* values, const int* indices); 00625 int getFromRHSVector_i( int num, double* values, const int* indices); 00626 int putIntoRHSVector_i(int num, const double* values, const int* indices); 00627 int getSolution_i( double* answers, int len ); 00628 int getSolnEntry_i( int eqnNumber, double& answer ); 00629 int getSolution_private( double* answers, int len, const bool image ); 00630 int sumIntoSystemMatrix_private(int numPtRows, const int* ptRows, 00631 int numPtCols, const int* ptCols, 00632 int numBlkRows, const int* blkRows, 00633 int numBlkCols, const int* blkCols, 00634 const double* const* values, 00635 const bool image ); 00636 int formResidual_i(double* values, int len); 00637 int formResidual_private(double* values, int len, const bool image); 00638 }; 00639 00640 #endif

Generated on Fri May 21 14:17:55 2004 by doxygen 1.3.7