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 #include <MumpsParallelSOE.h>
00031 #include <MumpsParallelSolver.h>
00032 #include <Matrix.h>
00033 #include <Graph.h>
00034 #include <Vertex.h>
00035 #include <VertexIter.h>
00036 #include <f2c.h>
00037 #include <math.h>
00038 #include <Channel.h>
00039 #include <FEM_ObjectBroker.h>
00040
00041 MumpsParallelSOE::MumpsParallelSOE(MumpsParallelSolver &theSolvr, int matType)
00042 :MumpsSOE(theSolvr, LinSOE_TAGS_MumpsParallelSOE, matType),
00043 processID(0), numChannels(0), theChannels(0), localCol(0), workArea(0),
00044 sizeWork(0), myB(0), myVectB(0)
00045 {
00046 theSolvr.setLinearSOE(*this);
00047 }
00048
00049
00050 MumpsParallelSOE::~MumpsParallelSOE()
00051 {
00052 if (theChannels != 0)
00053 delete [] theChannels;
00054
00055 if (localCol != 0)
00056 for (int i=0; i<numChannels; i++)
00057 if (localCol[i] != 0)
00058 delete localCol[i];
00059 delete [] localCol;
00060
00061 if (myB != 0)
00062 delete [] myB;
00063
00064 if (myVectB != 0)
00065 delete myVectB;
00066 }
00067
00068 int
00069 MumpsParallelSOE::setSize(Graph &theGraph)
00070 {
00071 int result = 0;
00072 int oldSize = size;
00073 int maxNumSubVertex = 0;
00074
00075
00076 int maxVertexTag = -1;
00077 Vertex *theVertex;
00078 int newNNZ = 0;
00079 size = theGraph.getNumVertex();
00080 int mySize = size;
00081
00082 VertexIter &theVertices = theGraph.getVertices();
00083 while ((theVertex = theVertices()) != 0) {
00084 int vertexTag = theVertex->getTag();
00085 if (vertexTag > maxVertexTag)
00086 maxVertexTag = vertexTag;
00087 const ID &theAdjacency = theVertex->getAdjacency();
00088 newNNZ += theAdjacency.Size() +1;
00089 }
00090
00091 if (matType != 0) {
00092
00093
00094 newNNZ -= size;
00095 newNNZ /= 2;
00096 newNNZ += size;
00097 }
00098
00099 nnz = newNNZ;
00100
00101 if (processID != 0) {
00102
00103
00104
00105
00106
00107 static ID data(1);
00108
00109 data(0) = maxVertexTag;
00110 Channel *theChannel = theChannels[0];
00111 theChannel->sendID(0, 0, data);
00112 theChannel->recvID(0, 0, data);
00113
00114 size = data(0);
00115
00116 } else {
00117
00118
00119
00120
00121
00122 static ID data(1);
00123 FEM_ObjectBroker theBroker;
00124 for (int j=0; j<numChannels; j++) {
00125 Channel *theChannel = theChannels[j];
00126 theChannel->recvID(0, 0, data);
00127 if (data(0) > maxVertexTag)
00128 maxVertexTag = data(0);
00129 }
00130
00131 data(0) = maxVertexTag;
00132
00133 for (int j=0; j<numChannels; j++) {
00134 Channel *theChannel = theChannels[j];
00135 theChannel->sendID(0, 0, data);
00136 }
00137 size = maxVertexTag;
00138 }
00139
00140 size+=1;
00141
00142 if (nnz > Asize) {
00143 if (A != 0) delete [] A;
00144 if (rowA != 0) delete [] rowA;
00145 if (colA != 0) delete [] colA;
00146
00147 A = new double[nnz];
00148 rowA = new int[nnz];
00149 colA = new int[nnz];
00150
00151 if (rowA == 0 || A == 0 || colA == 0) {
00152 opserr << "WARNING SparseGenColLinSOE::SparseGenColLinSOE :";
00153 opserr << " ran out of memory for A and rowA with nnz = ";
00154 opserr << nnz << " \n";
00155 size = 0; Asize = 0; nnz = 0;
00156 result = -1;
00157 }
00158 Asize = nnz;
00159 }
00160
00161
00162 if (size > Bsize) {
00163
00164
00165 if (B != 0) delete [] B;
00166 if (X != 0) delete [] X;
00167 if (myB != 0) delete [] myB;
00168 if (workArea != 0) delete [] workArea;
00169 if (colStartA != 0) delete [] colStartA;
00170
00171
00172 B = new double[size];
00173 X = new double[size];
00174 myB = new double[size];
00175 workArea = new double[size];
00176 colStartA = new int[size+1];
00177
00178 if (B == 0 || X == 0 || colStartA == 0 || workArea == 0 || myB == 0) {
00179 opserr << "WARNING MumpsSOE::MumpsSOE :";
00180 opserr << " ran out of memory for vectors (size) (";
00181 opserr << size << ") \n";
00182 size = 0; Bsize = 0;
00183 result = -1;
00184 }
00185 else
00186 Bsize = size;
00187 }
00188
00189
00190 for (int j=0; j<size; j++) {
00191 B[j] = 0;
00192 X[j] = 0;
00193 myB[j] = 0;
00194 }
00195
00196
00197 if (size != oldSize) {
00198 if (vectX != 0) delete vectX;
00199 if (vectB != 0) delete vectB;
00200 if (myVectB != 0) delete myVectB;
00201
00202 vectX = new Vector(X,size);
00203 vectB = new Vector(B,size);
00204 myVectB = new Vector(myB, size);
00205 }
00206
00207
00208 if (size != 0) {
00209 colStartA[0] = 0;
00210 int startLoc = 0;
00211 int lastLoc = 0;
00212 for (int a=0; a<size; a++) {
00213
00214 theVertex = theGraph.getVertexPtr(a);
00215 if (theVertex != 0) {
00216
00217 int vertexTag = theVertex->getTag();
00218 rowA[lastLoc++] = vertexTag;
00219 const ID &theAdjacency = theVertex->getAdjacency();
00220 int idSize = theAdjacency.Size();
00221
00222
00223
00224 if (matType != 0) {
00225
00226
00227 for (int i=0; i<idSize; i++) {
00228 int row = theAdjacency(i);
00229 if (row > vertexTag) {
00230 bool foundPlace = false;
00231
00232 for (int j=startLoc; j<lastLoc; j++)
00233 if (rowA[j] > row) {
00234
00235
00236 for (int k=lastLoc; k>j; k--)
00237 rowA[k] = rowA[k-1];
00238 rowA[j] = row;
00239 foundPlace = true;
00240 j = lastLoc;
00241 }
00242
00243 if (foundPlace == false)
00244 rowA[lastLoc] = row;
00245 lastLoc++;
00246 }
00247 }
00248
00249 } else {
00250
00251
00252 for (int i=0; i<idSize; i++) {
00253 int row = theAdjacency(i);
00254 bool foundPlace = false;
00255
00256 for (int j=startLoc; j<lastLoc; j++)
00257 if (rowA[j] > row) {
00258
00259
00260 for (int k=lastLoc; k>j; k--)
00261 rowA[k] = rowA[k-1];
00262 rowA[j] = row;
00263 foundPlace = true;
00264 j = lastLoc;
00265 }
00266 if (foundPlace == false)
00267 rowA[lastLoc] = row;
00268
00269 lastLoc++;
00270 }
00271 }
00272 }
00273 colStartA[a+1] = lastLoc;
00274 startLoc = lastLoc;
00275 }
00276 }
00277
00278
00279 int count = 0;
00280 for (int i=0; i<size; i++) {
00281 for (int k=colStartA[i]; k<colStartA[i+1]; k++)
00282 colA[count++] = i;
00283 }
00284
00285 LinearSOESolver *theSolvr = this->getSolver();
00286 int solverOK = theSolvr->setSize();
00287 if (solverOK < 0) {
00288 opserr << "WARNING:MumpsParallelSOE::setSize :";
00289 opserr << " solver failed setSize()\n";
00290 return solverOK;
00291 }
00292
00293 return result;
00294 }
00295
00296
00297 int
00298 MumpsParallelSOE::solve(void)
00299 {
00300 int resSolver = 0;
00301
00302
00303
00304
00305
00306 if (processID != 0) {
00307
00308
00309 Channel *theChannel = theChannels[0];
00310 theChannel->sendVector(0, 0, *myVectB);
00311 LinearSOESolver *theSoeSolver = this->getSolver();
00312 resSolver = this->LinearSOE::solve();
00313
00314 if (resSolver == 0) {
00315
00316 theChannel->recvVector(0, 0, *vectX);
00317 theChannel->recvVector(0, 0, *vectB);
00318 factored = true;
00319 }
00320 }
00321
00322
00323
00324
00325
00326 else {
00327
00328
00329 *vectB = *myVectB;
00330
00331
00332 for (int j=0; j<numChannels; j++) {
00333
00334 Channel *theChannel = theChannels[j];
00335 theChannel->recvVector(0, 0, *vectX);
00336 *vectB += *vectX;
00337 }
00338
00339
00340 resSolver = this->LinearSOE::solve();
00341
00342
00343 if (resSolver == 0) {
00344 for (int j=0; j<numChannels; j++) {
00345 Channel *theChannel = theChannels[j];
00346 theChannel->sendVector(0, 0, *vectX);
00347 theChannel->sendVector(0, 0, *vectB);
00348 }
00349 }
00350 }
00351
00352 return resSolver;
00353 }
00354
00355
00356
00357 int
00358 MumpsParallelSOE::addB(const Vector &v, const ID &id, double fact)
00359 {
00360
00361 if (fact == 0.0) return 0;
00362
00363 int idSize = id.Size();
00364
00365 if (idSize != v.Size() ) {
00366 opserr << "SparseGenColLinSOE::addB() ";
00367 opserr << " - Vector and ID not of similar sizes\n";
00368 return -1;
00369 }
00370
00371 if (fact == 1.0) {
00372 for (int i=0; i<idSize; i++) {
00373 int pos = id(i);
00374 if (pos <size && pos >= 0)
00375 myB[pos] += v(i);
00376 }
00377 } else if (fact == -1.0) {
00378 for (int i=0; i<idSize; i++) {
00379 int pos = id(i);
00380 if (pos <size && pos >= 0)
00381 myB[pos] -= v(i);
00382 }
00383 } else {
00384 for (int i=0; i<idSize; i++) {
00385 int pos = id(i);
00386 if (pos <size && pos >= 0)
00387 myB[pos] += v(i) * fact;
00388 }
00389 }
00390 return 0;
00391 }
00392
00393
00394 int
00395 MumpsParallelSOE::setB(const Vector &v, double fact)
00396 {
00397
00398 if (fact == 0.0) return 0;
00399
00400
00401 if (v.Size() != size) {
00402 opserr << "WARNING DistributedBandGenLinSOE::setB() -";
00403 opserr << " incomptable sizes " << size << " and " << v.Size() << endln;
00404 return -1;
00405 }
00406
00407 if (fact == 1.0) {
00408 for (int i=0; i<size; i++) {
00409 myB[i] = v(i);
00410 }
00411 } else if (fact == -1.0) {
00412 for (int i=0; i<size; i++) {
00413 myB[i] = -v(i);
00414 }
00415 } else {
00416 for (int i=0; i<size; i++) {
00417 myB[i] = v(i) * fact;
00418 }
00419 }
00420 return 0;
00421 }
00422
00423 void
00424 MumpsParallelSOE::zeroB(void)
00425 {
00426 double *Bptr = myB;
00427 for (int i=0; i<size; i++)
00428 *Bptr++ = 0;
00429 }
00430
00431
00432 const Vector &
00433 MumpsParallelSOE::getB(void)
00434 {
00435
00436 if (processID != 0) {
00437 Channel *theChannel = theChannels[0];
00438
00439
00440 theChannel->sendVector(0, 0, *myVectB);
00441 theChannel->recvVector(0, 0, *vectB);
00442 }
00443
00444
00445
00446
00447
00448 else {
00449
00450 *vectB = *myVectB;
00451
00452 Vector remoteB(workArea, size);
00453
00454
00455 for (int j=0; j<numChannels; j++) {
00456
00457 Channel *theChannel = theChannels[j];
00458 theChannel->recvVector(0, 0, remoteB);
00459 *vectB += remoteB;
00460 }
00461
00462
00463 for (int j=0; j<numChannels; j++) {
00464 Channel *theChannel = theChannels[j];
00465 theChannel->sendVector(0, 0, *vectB);
00466 }
00467 }
00468
00469 return *vectB;
00470 }
00471
00472
00473 int
00474 MumpsParallelSOE::sendSelf(int commitTag, Channel &theChannel)
00475 {
00476 int sendID =0;
00477
00478
00479
00480
00481
00482
00483 if (processID == 0) {
00484
00485
00486 bool found = false;
00487 for (int i=0; i<numChannels; i++)
00488 if (theChannels[i] == &theChannel) {
00489 sendID = i+1;
00490 found = true;
00491 }
00492
00493
00494 if (found == false) {
00495 int nextNumChannels = numChannels + 1;
00496 Channel **nextChannels = new Channel *[nextNumChannels];
00497 if (nextNumChannels == 0) {
00498 opserr << "MumpsParallelSOE::sendSelf() - failed to allocate channel array of size: " <<
00499 nextNumChannels << endln;
00500 return -1;
00501 }
00502 for (int i=0; i<numChannels; i++)
00503 nextChannels[i] = theChannels[i];
00504 nextChannels[numChannels] = &theChannel;
00505
00506 numChannels = nextNumChannels;
00507
00508 if (theChannels != 0)
00509 delete [] theChannels;
00510
00511 theChannels = nextChannels;
00512
00513 if (localCol != 0)
00514 delete [] localCol;
00515 localCol = new ID *[numChannels];
00516 if (localCol == 0) {
00517 opserr << "MumpsParallelSOE::sendSelf() - failed to allocate id array of size: " <<
00518 nextNumChannels << endln;
00519 return -1;
00520 }
00521 for (int i=0; i<numChannels; i++)
00522 localCol[i] = 0;
00523
00524
00525 sendID = numChannels;
00526 }
00527
00528 } else
00529 sendID = processID;
00530
00531
00532
00533 ID idData(2);
00534 idData(0) = sendID;
00535 idData(1) = matType;
00536
00537 int res = theChannel.sendID(0, commitTag, idData);
00538 if (res < 0) {
00539 opserr <<"WARNING MumpsParallelSOE::sendSelf() - failed to send data\n";
00540 return -1;
00541 }
00542
00543 return 0;
00544 }
00545
00546
00547 int
00548 MumpsParallelSOE::recvSelf(int commitTag, Channel &theChannel, FEM_ObjectBroker &theBroker)
00549 {
00550 ID idData(2);
00551 int res = theChannel.recvID(0, commitTag, idData);
00552 if (res < 0) {
00553 opserr <<"WARNING MumpsParallelSOE::recvSelf() - failed to send data\n";
00554 return -1;
00555 }
00556 processID = idData(0);
00557 matType = idData(1);
00558
00559 numChannels = 1;
00560 theChannels = new Channel *[1];
00561 theChannels[0] = &theChannel;
00562
00563 localCol = new ID *[numChannels];
00564 for (int i=0; i<numChannels; i++)
00565 localCol[i] = 0;
00566
00567 return 0;
00568 }
00569
00570
00571 int
00572 MumpsParallelSOE::setProcessID(int dTag)
00573 {
00574 processID = dTag;
00575 return 0;
00576 }
00577
00578 int
00579 MumpsParallelSOE::setChannels(int nChannels, Channel **theC)
00580 {
00581 numChannels = nChannels;
00582
00583 if (theChannels != 0)
00584 delete [] theChannels;
00585
00586 theChannels = new Channel *[numChannels];
00587 for (int i=0; i<numChannels; i++)
00588 theChannels[i] = theC[i];
00589
00590
00591 localCol = new ID *[nChannels];
00592 for (int i=0; i<numChannels; i++)
00593 localCol[i] = 0;
00594
00595 return 0;
00596 }