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