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 #include <DistributedSparseGenColLinSOE.h>
00032 #include <SparseGenColLinSolver.h>
00033 #include <Matrix.h>
00034 #include <Graph.h>
00035 #include <Vertex.h>
00036 #include <VertexIter.h>
00037 #include <f2c.h>
00038 #include <math.h>
00039 #include <Channel.h>
00040 #include <FEM_ObjectBroker.h>
00041
00042 DistributedSparseGenColLinSOE::DistributedSparseGenColLinSOE(SparseGenColLinSolver &theSolvr)
00043 :SparseGenColLinSOE(theSolvr, LinSOE_TAGS_DistributedSparseGenColLinSOE),
00044 processID(0), numChannels(0), theChannels(0), localCol(0), workArea(0), sizeWork(0), myB(0),
00045 myVectB(0)
00046
00047 {
00048 theSolvr.setLinearSOE(*this);
00049 }
00050
00051
00052 DistributedSparseGenColLinSOE::~DistributedSparseGenColLinSOE()
00053 {
00054 if (theChannels != 0)
00055 delete [] theChannels;
00056
00057 if (localCol != 0)
00058 for (int i=0; i<numChannels; i++)
00059 if (localCol[i] != 0)
00060 delete localCol[i];
00061 delete [] localCol;
00062
00063 if (myB != 0)
00064 delete [] myB;
00065
00066 if (myVectB != 0)
00067 delete myVectB;
00068 }
00069
00070
00071
00072 int
00073 DistributedSparseGenColLinSOE::setSize(Graph &theGraph)
00074 {
00075 int result = 0;
00076 int oldSize = size;
00077 int maxNumSubVertex = 0;
00078
00079
00080
00081 if (processID != 0) {
00082 Channel *theChannel = theChannels[0];
00083 theGraph.sendSelf(0, *theChannel);
00084
00085 static ID data(2);
00086 theChannel->recvID(0, 0, data);
00087 size = data(0);
00088 nnz = data(1);
00089
00090 ID *subMap = new ID(theGraph.getNumVertex());
00091 localCol[0] = subMap;
00092 Vertex *vertex;
00093 VertexIter &theSubVertices = theGraph.getVertices();
00094 int cnt = 0;
00095 while((vertex = theSubVertices()) != 0)
00096 (*subMap)(cnt++) = vertex->getTag();
00097
00098 theChannel->sendID(0, 0, *subMap);
00099
00100 if (nnz > Asize) {
00101 if (rowA != 0)
00102 delete [] rowA;
00103
00104 rowA = new int[nnz];
00105
00106 if (rowA == 0) {
00107 opserr << "WARNING SparseGenColLinSOE::SparseGenColLinSOE :";
00108 opserr << " ran out of memory for A and rowA with nnz = ";
00109 opserr << nnz << " \n";
00110 size = 0; Asize = 0; nnz = 0;
00111 result = -1;
00112 }
00113 }
00114
00115 if (size > Bsize) {
00116
00117 if (colStartA != 0)
00118 delete [] colStartA;
00119 colStartA = new int[size+1];
00120
00121 if (colStartA == 0) {
00122 opserr << "WARNING SparseGenColLinSOE::SparseGenColLinSOE :";
00123 opserr << " ran out of memory for vectors (size) (";
00124 opserr << size << ") \n";
00125 size = 0; Bsize = 0;
00126 result = -1;
00127 }
00128 }
00129
00130 ID rowAdata(rowA, nnz);
00131 ID colStartAdata(colStartA, size+1);
00132 theChannel->recvID(0, 0, rowAdata);
00133 theChannel->recvID(0, 0, colStartAdata);
00134 }
00135
00136
00137
00138
00139 else {
00140
00141
00142
00143 FEM_ObjectBroker theBroker;
00144 for (int j=0; j<numChannels; j++) {
00145 Channel *theChannel = theChannels[j];
00146 Graph theSubGraph;
00147 theSubGraph.recvSelf(0, *theChannel, theBroker);
00148 theGraph.merge(theSubGraph);
00149
00150 int numSubVertex = theSubGraph.getNumVertex();
00151 ID *subMap = new ID(numSubVertex);
00152 localCol[j] = subMap;
00153 if (numSubVertex > maxNumSubVertex)
00154 maxNumSubVertex = numSubVertex;
00155 }
00156
00157 size = theGraph.getNumVertex();
00158
00159
00160
00161
00162
00163 Vertex *theVertex;
00164 VertexIter &theVertices = theGraph.getVertices();
00165 nnz = 0;
00166 while ((theVertex = theVertices()) != 0) {
00167 const ID &theAdjacency = theVertex->getAdjacency();
00168 nnz += theAdjacency.Size() +1;
00169 }
00170
00171 static ID data(2);
00172 data(0) = size;
00173 data(1) = nnz;
00174
00175
00176
00177
00178 for (int j=0; j<numChannels; j++) {
00179 Channel *theChannel = theChannels[j];
00180 theChannel->sendID(0, 0, data);
00181
00182 ID *subMap = localCol[j];
00183 theChannel->recvID(0, 0, *subMap);
00184 }
00185
00186 if (nnz > Asize) {
00187 if (rowA != 0)
00188 delete [] rowA;
00189
00190 if (workArea != 0)
00191 delete [] workArea;
00192
00193 rowA = new int[nnz];
00194 workArea = new double[nnz];
00195
00196 sizeWork = nnz;
00197
00198 if (rowA == 0 || workArea == 0) {
00199 opserr << "WARNING SparseGenColLinSOE::SparseGenColLinSOE :";
00200 opserr << " ran out of memory for A and rowA with nnz = ";
00201 opserr << nnz << " \n";
00202 size = 0; Asize = 0; nnz = 0;
00203 result = -1;
00204 }
00205 }
00206
00207 if (size > Bsize) {
00208
00209 if (colStartA != 0)
00210 delete [] colStartA;
00211 colStartA = new int[size+1];
00212
00213 if (colStartA == 0) {
00214 opserr << "WARNING SparseGenColLinSOE::SparseGenColLinSOE :";
00215 opserr << " ran out of memory for vectors (size) (";
00216 opserr << size << ") \n";
00217 size = 0; Bsize = 0;
00218 result = -1;
00219 }
00220 }
00221
00222
00223 if (size != 0) {
00224 colStartA[0] = 0;
00225 int startLoc = 0;
00226 int lastLoc = 0;
00227 for (int a=0; a<size; a++) {
00228
00229 theVertex = theGraph.getVertexPtr(a);
00230 if (theVertex == 0) {
00231 opserr << "WARNING:SparseGenColLinSOE::setSize :";
00232 opserr << " vertex " << a << " not in graph! - size set to 0\n";
00233 size = 0;
00234 return -1;
00235 }
00236
00237 rowA[lastLoc++] = theVertex->getTag();
00238 const ID &theAdjacency = theVertex->getAdjacency();
00239 int idSize = theAdjacency.Size();
00240
00241
00242 for (int i=0; i<idSize; i++) {
00243
00244 int row = theAdjacency(i);
00245 bool foundPlace = false;
00246
00247 for (int j=startLoc; j<lastLoc; j++)
00248 if (rowA[j] > row) {
00249
00250
00251 for (int k=lastLoc; k>j; k--)
00252
00253 rowA[k] = rowA[k-1];
00254 rowA[j] = row;
00255 foundPlace = true;
00256 j = lastLoc;
00257 }
00258 if (foundPlace == false)
00259 rowA[lastLoc] = row;
00260
00261 lastLoc++;
00262 }
00263 colStartA[a+1] = lastLoc;;
00264 startLoc = lastLoc;
00265 }
00266 }
00267
00268
00269 ID rowAdata(rowA, nnz);
00270 ID colStartAdata(colStartA, size+1);
00271
00272 for (int j=0; j<numChannels; j++) {
00273 Channel *theChannel = theChannels[j];
00274 theChannel->sendID(0, 0, rowAdata);
00275 theChannel->sendID(0, 0, colStartAdata);
00276 }
00277 }
00278
00279 if (nnz > Asize) {
00280 if (A != 0)
00281 delete [] A;
00282
00283 A = new double[nnz];
00284
00285 if (A == 0 || rowA == 0) {
00286 opserr << "WARNING SparseGenColLinSOE::SparseGenColLinSOE :";
00287 opserr << " ran out of memory for A and rowA with nnz = ";
00288 opserr << nnz << " \n";
00289 size = 0; Asize = 0; nnz = 0;
00290 result = -1;
00291 }
00292
00293 Asize = nnz;
00294 }
00295
00296
00297 for (int i=0; i<Asize; i++)
00298 A[i] = 0;
00299
00300 factored = false;
00301
00302 if (size > Bsize) {
00303
00304
00305 if (B != 0) delete [] B;
00306 if (X != 0) delete [] X;
00307 if (myB != 0) delete [] myB;
00308
00309
00310 B = new double[size];
00311 X = new double[size];
00312 myB = new double[size];
00313
00314 if (B == 0 || X == 0 || colStartA == 0 || myB == 0) {
00315 opserr << "WARNING SparseGenColLinSOE::SparseGenColLinSOE :";
00316 opserr << " ran out of memory for vectors (size) (";
00317 opserr << size << ") \n";
00318 size = 0; Bsize = 0;
00319 result = -1;
00320 }
00321 else
00322 Bsize = size;
00323 }
00324
00325
00326 for (int j=0; j<size; j++) {
00327 B[j] = 0;
00328 X[j] = 0;
00329 myB[j] = 0;
00330 }
00331
00332
00333 if (size != oldSize) {
00334 if (vectX != 0)
00335 delete vectX;
00336
00337 if (vectB != 0)
00338 delete vectB;
00339
00340 vectX = new Vector(X,size);
00341 vectB = new Vector(B,size);
00342 myVectB = new Vector(myB, size);
00343 }
00344
00345 LinearSOESolver *theSolvr = this->getSolver();
00346 int solverOK = theSolvr->setSize();
00347 if (solverOK < 0) {
00348 opserr << "WARNING:DistributedSparseGenColLinSOE::setSize :";
00349 opserr << " solver failed setSize()\n";
00350 return solverOK;
00351 }
00352
00353 return result;
00354 }
00355
00356
00357 int
00358 DistributedSparseGenColLinSOE::addA(const Matrix &m, const ID &id, double fact)
00359 {
00360
00361 if (fact == 0.0)
00362 return 0;
00363
00364 int idSize = id.Size();
00365
00366
00367 if (idSize != m.noRows() && idSize != m.noCols()) {
00368 opserr << "SparseGenColLinSOE::addA() ";
00369 opserr << " - Matrix and ID not of similar sizes\n";
00370 return -1;
00371 }
00372
00373 if (fact == 1.0) {
00374 for (int i=0; i<idSize; i++) {
00375 int col = id(i);
00376 if (col < size && col >= 0) {
00377 int startColLoc = colStartA[col];
00378 int endColLoc = colStartA[col+1];
00379 for (int j=0; j<idSize; j++) {
00380 int row = id(j);
00381 if (row <size && row >= 0) {
00382
00383 for (int k=startColLoc; k<endColLoc; k++)
00384 if (rowA[k] == row) {
00385 A[k] += m(j,i);
00386 k = endColLoc;
00387 }
00388 }
00389 }
00390 }
00391 }
00392 } else {
00393 for (int i=0; i<idSize; i++) {
00394 int col = id(i);
00395 if (col < size && col >= 0) {
00396 int startColLoc = colStartA[col];
00397 int endColLoc = colStartA[col+1];
00398 for (int j=0; j<idSize; j++) {
00399 int row = id(j);
00400 if (row <size && row >= 0) {
00401
00402 for (int k=startColLoc; k<endColLoc; k++)
00403 if (rowA[k] == row) {
00404 A[k] += fact * m(j,i);
00405 k = endColLoc;
00406 }
00407 }
00408 }
00409 }
00410 }
00411 }
00412 return 0;
00413 }
00414
00415 int
00416 DistributedSparseGenColLinSOE::addB(const Vector &v, const ID &id, double fact)
00417 {
00418
00419 if (fact == 0.0) return 0;
00420
00421 int idSize = id.Size();
00422
00423 if (idSize != v.Size() ) {
00424 opserr << "SparseGenColLinSOE::addB() ";
00425 opserr << " - Vector and ID not of similar sizes\n";
00426 return -1;
00427 }
00428
00429 if (fact == 1.0) {
00430 for (int i=0; i<idSize; i++) {
00431 int pos = id(i);
00432 if (pos <size && pos >= 0)
00433 myB[pos] += v(i);
00434 }
00435 } else if (fact == -1.0) {
00436 for (int i=0; i<idSize; i++) {
00437 int pos = id(i);
00438 if (pos <size && pos >= 0)
00439 myB[pos] -= v(i);
00440 }
00441 } else {
00442 for (int i=0; i<idSize; i++) {
00443 int pos = id(i);
00444 if (pos <size && pos >= 0)
00445 myB[pos] += v(i) * fact;
00446 }
00447 }
00448 return 0;
00449 }
00450
00451
00452 int
00453 DistributedSparseGenColLinSOE::setB(const Vector &v, double fact)
00454 {
00455
00456 if (fact == 0.0) return 0;
00457
00458
00459 if (v.Size() != size) {
00460 opserr << "WARNING DistributedBandGenLinSOE::setB() -";
00461 opserr << " incomptable sizes " << size << " and " << v.Size() << endln;
00462 return -1;
00463 }
00464
00465 if (fact == 1.0) {
00466 for (int i=0; i<size; i++) {
00467 myB[i] = v(i);
00468 }
00469 } else if (fact == -1.0) {
00470 for (int i=0; i<size; i++) {
00471 myB[i] = -v(i);
00472 }
00473 } else {
00474 for (int i=0; i<size; i++) {
00475 myB[i] = v(i) * fact;
00476 }
00477 }
00478 return 0;
00479 }
00480
00481 void
00482 DistributedSparseGenColLinSOE::zeroB(void)
00483 {
00484 double *Bptr = myB;
00485 for (int i=0; i<size; i++)
00486 *Bptr++ = 0;
00487 }
00488
00489 int
00490 DistributedSparseGenColLinSOE::solve(void)
00491 {
00492 static ID result(1);
00493
00494
00495
00496
00497
00498 if (processID != 0) {
00499 Channel *theChannel = theChannels[0];
00500
00501
00502 theChannel->sendVector(0, 0, *myVectB);
00503
00504
00505 if (factored == false) {
00506 Vector vectA(A, nnz);
00507 theChannel->sendVector(0, 0, vectA);
00508 }
00509
00510 LinearSOESolver *theSoeSolver = this->getSolver();
00511 if (theSoeSolver->getClassTag() == SOLVER_TAGS_DistributedSuperLU)
00512 this->LinearSOE::solve();
00513
00514
00515 theChannel->recvVector(0, 0, *vectX);
00516 theChannel->recvVector(0, 0, *vectB);
00517 theChannel->recvID(0, 0, result);
00518 factored = true;
00519 }
00520
00521
00522
00523
00524
00525 else {
00526
00527
00528 *vectB = *myVectB;
00529
00530
00531 for (int j=0; j<numChannels; j++) {
00532
00533
00534 Channel *theChannel = theChannels[j];
00535 theChannel->recvVector(0, 0, *vectX);
00536 *vectB += *vectX;
00537
00538 if (factored == false) {
00539 Vector vectA(workArea, nnz);
00540 theChannel->recvVector(0, 0, vectA);
00541 for (int i=0; i<nnz; i++)
00542 A[i] += workArea[i];
00543 }
00544
00545
00546
00547
00548
00549
00550
00551
00552
00553
00554
00555
00556
00557
00558
00559 }
00560
00561
00562
00563 result(0) = this->LinearSOE::solve();
00564
00565
00566 for (int j=0; j<numChannels; j++) {
00567 Channel *theChannel = theChannels[j];
00568 theChannel->sendVector(0, 0, *vectX);
00569 theChannel->sendVector(0, 0, *vectB);
00570
00571 theChannel->sendID(0, 0, result);
00572 }
00573 }
00574
00575 return result(0);
00576 }
00577
00578
00579 const Vector &
00580 DistributedSparseGenColLinSOE::getB(void)
00581 {
00582
00583 if (processID != 0) {
00584 Channel *theChannel = theChannels[0];
00585
00586
00587 theChannel->sendVector(0, 0, *myVectB);
00588 theChannel->recvVector(0, 0, *vectB);
00589 }
00590
00591
00592
00593
00594
00595 else {
00596
00597 *vectB = *myVectB;
00598
00599 Vector remoteB(workArea, size);
00600
00601
00602 for (int j=0; j<numChannels; j++) {
00603
00604 Channel *theChannel = theChannels[j];
00605 theChannel->recvVector(0, 0, remoteB);
00606 *vectB += remoteB;
00607 }
00608
00609
00610 for (int j=0; j<numChannels; j++) {
00611 Channel *theChannel = theChannels[j];
00612 theChannel->sendVector(0, 0, *vectB);
00613 }
00614 }
00615
00616 return *vectB;
00617 }
00618
00619
00620 int
00621 DistributedSparseGenColLinSOE::sendSelf(int commitTag, Channel &theChannel)
00622 {
00623 int sendID =0;
00624
00625
00626
00627
00628
00629
00630 if (processID == 0) {
00631
00632
00633 bool found = false;
00634 for (int i=0; i<numChannels; i++)
00635 if (theChannels[i] == &theChannel) {
00636 sendID = i+1;
00637 found = true;
00638 }
00639
00640
00641 if (found == false) {
00642 int nextNumChannels = numChannels + 1;
00643 Channel **nextChannels = new Channel *[nextNumChannels];
00644 if (nextNumChannels == 0) {
00645 opserr << "DistributedSparseGenColLinSOE::sendSelf() - failed to allocate channel array of size: " <<
00646 nextNumChannels << endln;
00647 return -1;
00648 }
00649 for (int i=0; i<numChannels; i++)
00650 nextChannels[i] = theChannels[i];
00651 nextChannels[numChannels] = &theChannel;
00652
00653 numChannels = nextNumChannels;
00654
00655 if (theChannels != 0)
00656 delete [] theChannels;
00657
00658 theChannels = nextChannels;
00659
00660 if (localCol != 0)
00661 delete [] localCol;
00662 localCol = new ID *[numChannels];
00663 if (localCol == 0) {
00664 opserr << "DistributedSparseGenColLinSOE::sendSelf() - failed to allocate id array of size: " <<
00665 nextNumChannels << endln;
00666 return -1;
00667 }
00668 for (int i=0; i<numChannels; i++)
00669 localCol[i] = 0;
00670
00671
00672 sendID = numChannels;
00673 }
00674
00675 } else
00676 sendID = processID;
00677
00678
00679
00680 ID idData(1);
00681 idData(0) = sendID;
00682
00683 int res = theChannel.sendID(0, commitTag, idData);
00684 if (res < 0) {
00685 opserr <<"WARNING DistributedSparseGenColLinSOE::sendSelf() - failed to send data\n";
00686 return -1;
00687 }
00688
00689 return 0;
00690 }
00691
00692
00693 int
00694 DistributedSparseGenColLinSOE::recvSelf(int commitTag, Channel &theChannel, FEM_ObjectBroker &theBroker)
00695 {
00696 ID idData(1);
00697 int res = theChannel.recvID(0, commitTag, idData);
00698 if (res < 0) {
00699 opserr <<"WARNING DistributedSparseGenColLinSOE::recvSelf() - failed to send data\n";
00700 return -1;
00701 }
00702 processID = idData(0);
00703
00704 numChannels = 1;
00705 theChannels = new Channel *[1];
00706 theChannels[0] = &theChannel;
00707
00708
00709 localCol = new ID *[numChannels];
00710 for (int i=0; i<numChannels; i++)
00711 localCol[i] = 0;
00712
00713 return 0;
00714 }
00715
00716
00717 int
00718 DistributedSparseGenColLinSOE::setProcessID(int dTag)
00719 {
00720 processID = dTag;
00721 return 0;
00722 }
00723
00724 int
00725 DistributedSparseGenColLinSOE::setChannels(int nChannels, Channel **theC)
00726 {
00727 numChannels = nChannels;
00728
00729 if (theChannels != 0)
00730 delete [] theChannels;
00731
00732 theChannels = new Channel *[numChannels];
00733 for (int i=0; i<numChannels; i++)
00734 theChannels[i] = theC[i];
00735
00736
00737 localCol = new ID *[nChannels];
00738 for (int i=0; i<numChannels; i++)
00739 localCol[i] = 0;
00740
00741 return 0;
00742 }