00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038 #include <fstream.h>
00039
00040 #include <ItpackLinSOE.h>
00041 #include <ItpackLinSolver.h>
00042 #include <Matrix.h>
00043 #include <Graph.h>
00044 #include <Vertex.h>
00045 #include <VertexIter.h>
00046 #include <math.h>
00047
00048 #include <Channel.h>
00049 #include <FEM_ObjectBroker.h>
00050
00051 ItpackLinSOE::ItpackLinSOE(ItpackLinSolver &the_Solver)
00052 :LinearSOE(the_Solver, LinSOE_TAGS_ItpackLinSOE),
00053 size(0), nnz(0), A(0), B(0), X(0), colA(0), rowStartA(0),
00054 vectX(0), vectB(0), Asize(0), Bsize(0), Aformed(false)
00055 {
00056 the_Solver.setLinearSOE(*this);
00057 }
00058
00059 ItpackLinSOE::~ItpackLinSOE()
00060 {
00061 if (A != 0) delete [] A;
00062 if (B != 0) delete [] B;
00063 if (X != 0) delete [] X;
00064 if (rowStartA != 0) delete [] rowStartA;
00065 if (colA != 0) delete [] colA;
00066 if (vectX != 0) delete vectX;
00067 if (vectB != 0) delete vectB;
00068 }
00069
00070 int
00071 ItpackLinSOE::getNumEqn(void) const
00072 {
00073 return size;
00074 }
00075
00076 int
00077 ItpackLinSOE::setSize(Graph &theGraph)
00078 {
00079 int result = 0;
00080 int oldSize = size;
00081 size = theGraph.getNumVertex();
00082
00083
00084 Vertex *theVertex;
00085 int newNNZ = 0;
00086 VertexIter &theVertices = theGraph.getVertices();
00087 while ((theVertex = theVertices()) != 0) {
00088 const ID &theAdjacency = theVertex->getAdjacency();
00089 newNNZ += theAdjacency.Size() +1;
00090 }
00091 nnz = newNNZ;
00092
00093 opserr << "ItpackLinSOE::setSize - n " << size << " nnz " << nnz << endln;
00094
00095 if (nnz > Asize) {
00096
00097 if (A != 0)
00098 delete [] A;
00099 if (colA != 0)
00100 delete [] colA;
00101
00102 A = new double[nnz];
00103 colA = new int[nnz];
00104
00105 if (A == 0 || colA == 0) {
00106 opserr << "WARNING ItpackLinSOE::ItpackLinSOE :";
00107 opserr << " ran out of memory for A and colA with nnz = ";
00108 opserr << newNNZ << " \n";
00109 size = 0; Asize = 0; nnz = 0;
00110 result = -1;
00111 }
00112
00113 Asize = nnz;
00114 }
00115
00116
00117 for (int i=0; i<Asize; i++)
00118 A[i] = 0;
00119
00120 if (size > Bsize) {
00121
00122
00123 if (B != 0) delete [] B;
00124 if (X != 0) delete [] X;
00125 if (rowStartA != 0) delete [] rowStartA;
00126
00127
00128 B = new double[size];
00129 X = new double[size];
00130 rowStartA = new int[size+1];
00131
00132 if (B == 0 || X == 0 || rowStartA == 0) {
00133 opserr << "WARNING ItpackLinSOE::ItpackLinSOE :";
00134 opserr << " ran out of memory for vectors (size) (";
00135 opserr << size << ") \n";
00136 size = 0; Bsize = 0;
00137 result = -1;
00138 }
00139 else
00140 Bsize = size;
00141 }
00142
00143
00144 for (int j=0; j<size; j++) {
00145 B[j] = 0;
00146 X[j] = 0;
00147 }
00148
00149
00150 if (size != oldSize) {
00151 if (vectX != 0)
00152 delete vectX;
00153
00154 if (vectB != 0)
00155 delete vectB;
00156
00157 vectX = new Vector(X,size);
00158 vectB = new Vector(B,size);
00159 }
00160
00161
00162 if (size != 0) {
00163 rowStartA[0] = 0;
00164 int startLoc = 0;
00165 int lastLoc = 0;
00166 for (int a=0; a<size; a++) {
00167
00168 theVertex = theGraph.getVertexPtr(a);
00169 if (theVertex == 0) {
00170 opserr << "WARNING:ItpackLinSOE::setSize :";
00171 opserr << " vertex " << a << " not in graph! - size set to 0\n";
00172 size = 0;
00173 return -1;
00174 }
00175
00176 colA[lastLoc++] = theVertex->getTag();
00177 const ID &theAdjacency = theVertex->getAdjacency();
00178 int idSize = theAdjacency.Size();
00179
00180
00181 for (int i=0; i<idSize; i++) {
00182
00183 int row = theAdjacency(i);
00184 bool foundPlace = false;
00185
00186 for (int j=startLoc; j<lastLoc; j++)
00187 if (colA[j] > row) {
00188
00189
00190 for (int k=lastLoc; k>j; k--)
00191 colA[k] = colA[k-1];
00192 colA[j] = row;
00193 foundPlace = true;
00194 j = lastLoc;
00195 }
00196 if (foundPlace == false)
00197 colA[lastLoc] = row;
00198
00199 lastLoc++;
00200 }
00201 rowStartA[a+1] = lastLoc;
00202 startLoc = lastLoc;
00203 }
00204 }
00205
00206
00207
00208
00209
00210
00211
00212
00213
00214
00215
00216
00217
00218
00219
00220
00221
00222 Aformed = false;
00223
00224
00225 LinearSOESolver *the_Solver = this->getSolver();
00226 int solverOK = the_Solver->setSize();
00227 if (solverOK < 0) {
00228 opserr << "WARNING:ItpackLinSOE::setSize :";
00229 opserr << " solver failed setSize()\n";
00230 return solverOK;
00231 }
00232 return result;
00233 }
00234
00235 int
00236 ItpackLinSOE::addA(const Matrix &m, const ID &id, double fact)
00237 {
00238
00239 if (fact == 0.0)
00240 return 0;
00241
00242 int idSize = id.Size();
00243
00244
00245 if (idSize != m.noRows() && idSize != m.noCols()) {
00246 opserr << "ItpackLinSOE::addA() ";
00247 opserr << " - Matrix and ID not of similar sizes\n";
00248 return -1;
00249 }
00250
00251 if (fact == 1.0) {
00252 for (int i=0; i<idSize; i++) {
00253 int row = id(i);
00254 if (row < size && row >= 0) {
00255 int startRowLoc = rowStartA[row];
00256 int endRowLoc = rowStartA[row+1];
00257 for (int j=0; j<idSize; j++) {
00258 int col = id(j);
00259 if (col <size && col >= 0) {
00260
00261 for (int k=startRowLoc; k<endRowLoc; k++)
00262 if (colA[k] == col) {
00263 A[k] += m(i,j);
00264 k = endRowLoc;
00265 }
00266 }
00267 }
00268 }
00269 }
00270 }
00271 else {
00272 for (int i=0; i<idSize; i++) {
00273 int row = id(i);
00274 if (row < size && row >= 0) {
00275 int startRowLoc = rowStartA[row];
00276 int endRowLoc = rowStartA[row+1];
00277 for (int j=0; j<idSize; j++) {
00278 int col = id(j);
00279 if (col <size && col >= 0) {
00280
00281 for (int k=startRowLoc; k<endRowLoc; k++)
00282 if (colA[k] == col) {
00283 A[k] += fact * m(i,j);
00284 k = endRowLoc;
00285 }
00286 }
00287 }
00288 }
00289 }
00290 }
00291
00292 return 0;
00293 }
00294
00295 int
00296 ItpackLinSOE::addB(const Vector &v, const ID &id, double fact)
00297 {
00298
00299
00300 if (fact == 0.0)
00301 return 0;
00302
00303 int idSize = id.Size();
00304
00305 if (idSize != v.Size() ) {
00306 opserr << "ItpackLinSOE::addB() ";
00307 opserr << " - Vector and ID not of similar sizes\n";
00308 return -1;
00309 }
00310
00311 if (fact == 1.0) {
00312 for (int i=0; i<idSize; i++) {
00313 int pos = id(i);
00314 if (pos <size && pos >= 0)
00315 B[pos] += v(i);
00316 }
00317 }
00318 else if (fact == -1.0) {
00319 for (int i=0; i<idSize; i++) {
00320 int pos = id(i);
00321 if (pos <size && pos >= 0)
00322 B[pos] -= v(i);
00323 }
00324 }
00325 else {
00326 for (int i=0; i<idSize; i++) {
00327 int pos = id(i);
00328 if (pos <size && pos >= 0)
00329 B[pos] += v(i) * fact;
00330 }
00331 }
00332
00333 return 0;
00334 }
00335
00336 int
00337 ItpackLinSOE::setB(const Vector &v, double fact)
00338 {
00339 if (v.Size() != size) {
00340 opserr << "WARNING ItpackLinSOE::setB() -";
00341 opserr << " incomptable sizes " << size << " and " << v.Size() << endln;
00342 return -1;
00343 }
00344
00345 if (fact == 0.0) {
00346 for (int i=0; i<size; i++)
00347 B[i] = 0.0;
00348 }
00349 else if (fact == 1.0) {
00350 for (int i=0; i<size; i++)
00351 B[i] = v(i);
00352 }
00353 else if (fact == -1.0) {
00354 for (int i=0; i<size; i++)
00355 B[i] = -v(i);
00356 }
00357 else {
00358 for (int i=0; i<size; i++)
00359 B[i] = v(i) * fact;
00360 }
00361 return 0;
00362 }
00363
00364 void
00365 ItpackLinSOE::zeroA(void)
00366 {
00367 double *Aptr = A;
00368 for (int i=0; i<Asize; i++)
00369 *Aptr++ = 0;
00370 Aformed = false;
00371 }
00372
00373 void
00374 ItpackLinSOE::zeroB(void)
00375 {
00376 double *Bptr = B;
00377 for (int i=0; i<size; i++)
00378 *Bptr++ = 0;
00379 }
00380
00381 void
00382 ItpackLinSOE::setX(int loc, double value)
00383 {
00384 if (loc < size && loc >=0)
00385 X[loc] = value;
00386 }
00387
00388 void
00389 ItpackLinSOE::setX(const Vector &x)
00390 {
00391 if (x.Size() != size) {
00392 opserr << "WARNING ItpackLinSOE::setX() -";
00393 opserr << " incomptable sizes " << size << " and " << x.Size() << endln;
00394 }
00395 else {
00396 for (int i = 0; i < size; i++)
00397 X[i] = x(i);
00398 }
00399 }
00400
00401 const Vector &
00402 ItpackLinSOE::getX(void)
00403 {
00404 if (vectX == 0) {
00405 opserr << "FATAL ItpackLinSOE::getX - vectX == 0";
00406 exit(-1);
00407 }
00408 return *vectX;
00409 }
00410
00411 const Vector &
00412 ItpackLinSOE::getB(void)
00413 {
00414 if (vectB == 0) {
00415 opserr << "FATAL ItpackLinSOE::getB - vectB == 0";
00416 exit(-1);
00417 }
00418 return *vectB;
00419 }
00420
00421 double
00422 ItpackLinSOE::normRHS(void)
00423 {
00424 double norm =0.0;
00425 for (int i=0; i<size; i++) {
00426 double Yi = B[i];
00427 norm += Yi*Yi;
00428 }
00429 return sqrt(norm);
00430 }
00431
00432 int
00433 ItpackLinSOE::setItpackLinSolver(ItpackLinSolver &newSolver)
00434 {
00435 newSolver.setLinearSOE(*this);
00436
00437 if (size != 0) {
00438 int solverOK = newSolver.setSize();
00439 if (solverOK < 0) {
00440 opserr << "WARNING:ItpackLinSOE::setSolver :";
00441 opserr << "the new solver could not setSeize() - staying with old\n";
00442 return -1;
00443 }
00444 }
00445
00446 return this->LinearSOE::setSolver(newSolver);
00447 }
00448
00449 int
00450 ItpackLinSOE::sendSelf(int cTag, Channel &theChannel)
00451 {
00452 return -1;
00453 }
00454
00455 int
00456 ItpackLinSOE::recvSelf(int cTag, Channel &theChannel,
00457 FEM_ObjectBroker &theBroker)
00458 {
00459 return -1;
00460 }