Exemplo n.º 1
0
int 
SymArpackSOE::setSize(Graph &theGraph)
{
    int result = 0;
    int oldSize = size;
    size = theGraph.getNumVertex();

    // fist itearte through the vertices of the graph to get nnz
    Vertex *theVertex;
    int newNNZ = 0;
    VertexIter &theVertices = theGraph.getVertices();
    while ((theVertex = theVertices()) != 0) {
        const ID &theAdjacency = theVertex->getAdjacency();
	newNNZ += theAdjacency.Size(); 
    }
    nnz = newNNZ;

    if (colA != 0)
      delete [] colA;

    colA = new int[newNNZ];	
    if (colA == 0) {
	opserr << "WARNING SymArpackSOE::SymArpackSOE :";
	opserr << " ran out of memory for colA with nnz = ";
      	opserr << newNNZ << " \n";
       	size = 0; nnz = 0;
       	result =  -1;
    } 
	
    factored = false;

    if (rowStartA != 0) 
      delete [] rowStartA;
    rowStartA = new int[size+1]; 
    if (rowStartA == 0) {
        opserr << "SymArpackSOE::ran out of memory for rowStartA." << endln;
        result = -1;
    }

    // fill in rowStartA and colA
    if (size != 0) {
        rowStartA[0] = 0;
        int startLoc = 0;
	int lastLoc = 0;

	for (int a=0; a<size; a++) {
	   theVertex = theGraph.getVertexPtr(a);
	   if (theVertex == 0) {
	        opserr << "WARNING:SymArpackSOE::setSize :";
	        opserr << " vertex " << a << " not in graph! - size set to 0\n";
	        size = 0;
	        return -1;
	   }

	   const ID &theAdjacency = theVertex->getAdjacency();
	   int idSize = theAdjacency.Size();
	
	// now we have to place the entries in the ID into order in colA
	   for (int i=0; i<idSize; i++) {
	      int row = theAdjacency(i);
	      bool foundPlace = false;
	 
	      for (int j=startLoc; j<lastLoc; j++)
	          if (colA[j] > row) { 
	      // move the entries already there one further on
	      // and place col in current location
	              for (int k=lastLoc; k>j; k--)
		          colA[k] = colA[k-1];
                      colA[j] = row;
		      foundPlace = true;
    	              j = lastLoc;
		  }
		  
	      if (foundPlace == false) // put in at the end
	      	   colA[lastLoc] = row;

	      lastLoc++;
	   }
	   rowStartA[a+1] = lastLoc;;	    
	   startLoc = lastLoc;
	}
    }

    // begin to choose different ordering schema.
    //   cout << "Enter DOF Numberer Type: \n";
    //   cout << "[1] Minimum Degree, [2] Nested Dissection, [3] RCM: ";
    int LSPARSE = 1;
    //   cin >> LSPARSE;

// call "C" function to form elimination tree and do the symbolic factorization.
//    nblks = symFactorization(rowStartA, colA, size, LSPARSE);
    nblks = symFactorization(rowStartA, colA, size, LSPARSE,
			     &xblk, &invp, &rowblks, &begblk, &first, &penv, &diag);

    // invoke setSize() on the Solver
    EigenSolver *theSolvr = this->getSolver();
    int solverOK = theSolvr->setSize();
    if (solverOK < 0) {
	opserr << "WARNING:BandArpackSOE::setSize :";
	opserr << " solver failed setSize()\n";
	return solverOK;
    } 


    return result;
}
int 
DistributedSparseGenColLinSOE::setSize(Graph &theGraph)
{
  int result = 0;
  int oldSize = size;
  int maxNumSubVertex = 0;

  // if subprocess, collect graph, send it off, 
  // vector back containing size of system, etc.
  if (processID != 0) {
    Channel *theChannel = theChannels[0];
    theGraph.sendSelf(0, *theChannel);
    
    static ID data(2);
    theChannel->recvID(0, 0, data);
    size = data(0);
    nnz = data(1);

    ID *subMap = new ID(theGraph.getNumVertex());
    localCol[0] = subMap;
    Vertex *vertex;
    VertexIter &theSubVertices = theGraph.getVertices();
    int cnt = 0;
    while((vertex = theSubVertices()) != 0) 
      (*subMap)(cnt++) = vertex->getTag();

    theChannel->sendID(0, 0, *subMap);

    if (nnz > Asize) { // we have to get more space for A and rowA
      if (rowA != 0)
	delete [] rowA;
      
      rowA = new int[nnz];
	
      if (rowA == 0) {
	opserr << "WARNING SparseGenColLinSOE::SparseGenColLinSOE :";
	opserr << " ran out of memory for A and rowA with nnz = ";
	opserr << nnz << " \n";
	size = 0; Asize = 0; nnz = 0;
	result =  -1;
      } 
    }

    if (size > Bsize) { // we have to get space for the vectors
	
      if (colStartA != 0) 
	delete [] colStartA;
      colStartA = new int[size+1]; 
      
      if (colStartA == 0) {
	opserr << "WARNING SparseGenColLinSOE::SparseGenColLinSOE :";
	opserr << " ran out of memory for vectors (size) (";
	opserr << size << ") \n";
	size = 0; Bsize = 0;
	result =  -1;
      }    
    }

    ID rowAdata(rowA, nnz);
    ID colStartAdata(colStartA, size+1);
    theChannel->recvID(0, 0, rowAdata);
    theChannel->recvID(0, 0, colStartAdata);
  } 
  
  // if main domain, collect graphs from all subdomains,
  // merge into 1, number this one, send to subdomains the
  // id containing dof tags & start id's.
  else {

    // from each distributed soe recv it's graph
    // and merge them into master graph
    FEM_ObjectBroker theBroker;
    for (int j=0; j<numChannels; j++) {
      Channel *theChannel = theChannels[j];
      Graph theSubGraph;
      theSubGraph.recvSelf(0, *theChannel, theBroker);
      theGraph.merge(theSubGraph);

      int numSubVertex = theSubGraph.getNumVertex();
      ID *subMap = new ID(numSubVertex);
      localCol[j] = subMap;
      if (numSubVertex > maxNumSubVertex)
	maxNumSubVertex = numSubVertex;
    }

    size = theGraph.getNumVertex();
  
    //
    // determine the number of non-zeros
    //

    Vertex *theVertex;
    VertexIter &theVertices = theGraph.getVertices();
    nnz = 0;
    while ((theVertex = theVertices()) != 0) {
	const ID &theAdjacency = theVertex->getAdjacency();
	nnz += theAdjacency.Size() +1; // the +1 is for the diag entry
    }

    static ID data(2);
    data(0) = size;
    data(1) = nnz;

    // to each distributed soe send the size data
    // and merge them into master graph

    for (int j=0; j<numChannels; j++) {
      Channel *theChannel = theChannels[j];
      theChannel->sendID(0, 0, data);

      ID *subMap = localCol[j];
      theChannel->recvID(0, 0, *subMap);
    }    

    if (nnz > Asize) { // we have to get more space for A and rowA
      if (rowA != 0)
	delete [] rowA;

      if (workArea != 0) 
	delete [] workArea;

      rowA = new int[nnz];
      workArea = new double[nnz];

      sizeWork = nnz;
	
      if (rowA == 0 || workArea == 0) {
	opserr << "WARNING SparseGenColLinSOE::SparseGenColLinSOE :";
	opserr << " ran out of memory for A and rowA with nnz = ";
	opserr << nnz << " \n";
	size = 0; Asize = 0; nnz = 0;
	result =  -1;
      } 
    }

    if (size > Bsize) { // we have to get space for the vectors
	
      if (colStartA != 0) 
	delete [] colStartA;
      colStartA = new int[size+1]; 
      
      if (colStartA == 0) {
	opserr << "WARNING SparseGenColLinSOE::SparseGenColLinSOE :";
	opserr << " ran out of memory for vectors (size) (";
	opserr << size << ") \n";
	size = 0; Bsize = 0;
	result =  -1;
      }    
    }

    // fill in colStartA and rowA
    if (size != 0) {
      colStartA[0] = 0;
      int startLoc = 0;
      int lastLoc = 0;
      for (int a=0; a<size; a++) {
	
	theVertex = theGraph.getVertexPtr(a);
	if (theVertex == 0) {
	  opserr << "WARNING:SparseGenColLinSOE::setSize :";
	  opserr << " vertex " << a << " not in graph! - size set to 0\n";
	  size = 0;
	  return -1;
	}
	
	rowA[lastLoc++] = theVertex->getTag(); // place diag in first
	const ID &theAdjacency = theVertex->getAdjacency();
	int idSize = theAdjacency.Size();
	
	// now we have to place the entries in the ID into order in rowA
	for (int i=0; i<idSize; i++) {
	  
	  int row = theAdjacency(i);
	  bool foundPlace = false;
	  // find a place in rowA for current col
	  for (int j=startLoc; j<lastLoc; j++)
	    if (rowA[j] > row) { 
	      // move the entries already there one further on
	      // and place col in current location
	      for (int k=lastLoc; k>j; k--)
		
		rowA[k] = rowA[k-1];
	      rowA[j] = row;
	      foundPlace = true;
	      j = lastLoc;
	    }
	  if (foundPlace == false) // put in at the end
	    rowA[lastLoc] = row;
	  
	  lastLoc++;
	}
	colStartA[a+1] = lastLoc;;	    
	startLoc = lastLoc;
      }
    }


    ID rowAdata(rowA, nnz);
    ID colStartAdata(colStartA, size+1);

    for (int j=0; j<numChannels; j++) {
      Channel *theChannel = theChannels[j];
      theChannel->sendID(0, 0, rowAdata);
      theChannel->sendID(0, 0, colStartAdata);
    }
  }

  if (nnz > Asize) { // we have to get more space for A and rowA
    if (A != 0) 
      delete [] A;
    
    A = new double[nnz];
    
    if (A == 0 || rowA == 0) {
      opserr << "WARNING SparseGenColLinSOE::SparseGenColLinSOE :";
      opserr << " ran out of memory for A and rowA with nnz = ";
      opserr << nnz << " \n";
      size = 0; Asize = 0; nnz = 0;
      result =  -1;
    } 
    
    Asize = nnz;
  }
  
  // zero the matrix
  for (int i=0; i<Asize; i++)
    A[i] = 0;

  factored = false;

  if (size > Bsize) { // we have to get space for the vectors
    
    // delete the old	
    if (B != 0) delete [] B;
    if (X != 0) delete [] X;
    if (myB != 0) delete [] myB;
    
    // create the new
    B = new double[size];
    X = new double[size];
    myB = new double[size];
    
    if (B == 0 || X == 0 || colStartA == 0 || myB == 0) {
      opserr << "WARNING SparseGenColLinSOE::SparseGenColLinSOE :";
      opserr << " ran out of memory for vectors (size) (";
      opserr << size << ") \n";
      size = 0; Bsize = 0;
      result =  -1;
    }
    else
      Bsize = size;
  }

  // zero the vectors
  for (int j=0; j<size; j++) {
    B[j] = 0;
    X[j] = 0;
    myB[j] = 0;
  }

  // create new Vectors objects
  if (size != oldSize) {
    if (vectX != 0)
      delete vectX;
    
    if (vectB != 0)
      delete vectB;
    
    vectX = new Vector(X,size);
    vectB = new Vector(B,size);	
    myVectB = new Vector(myB, size);
  }

  LinearSOESolver *theSolvr = this->getSolver();
  int solverOK = theSolvr->setSize();
  if (solverOK < 0) {
    opserr << "WARNING:DistributedSparseGenColLinSOE::setSize :";
    opserr << " solver failed setSize()\n";
    return solverOK;
  }    

  return result;    
}
Exemplo n.º 3
0
int 
SkylineSPD::setSize(Graph &theGraph)
{
    int oldSize = size;
    int result = 0;
    size = theGraph.getNumVertex();

    // check we have enough space in iDiagLoc and iLastCol
    // if not delete old and create new
    if (size > Bsize) { 
      if (iDiagLoc != 0) delete [] iDiagLoc;
      if (RowTop != 0) delete [] RowTop;
      if (topRowPtr != 0) delete [] topRowPtr;
      if (invD != 0) delete [] invD;
      // if (topRowPtr != 0) free((void *)topRowPtr);
      
      iDiagLoc = new int[size];
      RowTop = new int[size];
      invD = new double[size]; 
      //      topRowPtr = new double *[size] ;
      topRowPtr = (double **)malloc(size *sizeof(double *));


      if (iDiagLoc == 0 || RowTop == 0 || topRowPtr == 0 || invD == 0) {
	opserr << "WARNING SkylineSPD::setSize() : ";
	opserr << " - ran out of memory for iDiagLoc\n";
	size = 0; Asize = 0;
	result = -1;
      }
    }
    
    // zero out iDiagLoc 
    for (int i=0; i<size; i++) {
      iDiagLoc[i] = 0;
    }
    
    // now we go through the vertices to find the height of each col and
    // width of each row from the connectivity information.
    
    Vertex *vertexPtr;
    VertexIter &theVertices = theGraph.getVertices();
    
    while ((vertexPtr = theVertices()) != 0) {
      int vertexNum = vertexPtr->getTag();
      const ID &theAdjacency = vertexPtr->getAdjacency();
      int iiDiagLoc = iDiagLoc[vertexNum];
      int *iiDiagLocPtr = &(iDiagLoc[vertexNum]);
      
      for (int i=0; i<theAdjacency.Size(); i++) {
	int otherNum = theAdjacency(i);
	int diff = vertexNum-otherNum;
	if (diff > 0) {
	  if (iiDiagLoc < diff) {
	    iiDiagLoc = diff;
	    *iiDiagLocPtr = diff;
	  }
	} 
      }
    }
    
    
    // now go through iDiagLoc, adding 1 for the diagonal element
    // and then adding previous entry to give current location.
    if (iDiagLoc != 0)
      iDiagLoc[0] = 1; // NOTE FORTRAN ARRAY LOCATION
    
    for (int j=1; j<size; j++)
      iDiagLoc[j] = iDiagLoc[j] + 1 + iDiagLoc[j-1];
    
    if (iDiagLoc != 0)       
      profileSize = iDiagLoc[size-1];
    
    // check if we need more space to hold A
    // if so then go get it
    if (profileSize > Asize) { 
      
      // delete old space
      if (A != 0)
	delete [] A;
      
      // get new space
      A = new double[profileSize];
      
      if (A == 0) {
	opserr << "SkylineSPD::SkylineSPD :";
	opserr << " ran out of memory for A (size,Profile) (";
	opserr << size <<", " << profileSize << ") \n";
	size = 0;  Asize = 0;  profileSize = 0;
	result = -1;
      }
      else 
	Asize = profileSize;
    }
    
    // zero the matrix
    for (int k=0; k<profileSize; k++)
      A[k] = 0;
    
    isAfactored = false;
    isAcondensed = false;    
    
    if (size > Bsize) { // we have to get another space for A
      
      // delete the old
      if (B != 0) delete [] B;
      if (X != 0) delete [] X;
      
      // create the new	
      B = new double[size];
      X = new double[size];
      
      if (B == 0 || X == 0 ) {
	opserr << "SkylineSPD::SkylineSPD :";
	opserr << " ran out of memory for vectors (size) (";
	opserr << size << ") \n";
	size = 0; Bsize = 0;
	result = -1;
      }
    }
    
    // zero the vectors
    for (int l=0; l<size; l++) {
      B[l] = 0;
      X[l] = 0;
    }
    
    if (size != oldSize) {
      
      if (vectX != 0)
	delete vectX;
      if (vectB != 0)
	delete vectB;
      
      vectX = new Vector(X,size);
      vectB = new Vector(B,size);
      
      if (size > Bsize)
	Bsize = size;
    }
    
    RowTop[0] = 0;
    topRowPtr[0] = A;
    for (int j=1; j<size; j++) {
      int icolsz = iDiagLoc[j] - iDiagLoc[j-1];
      RowTop[j] = j - icolsz +  1;
      topRowPtr[j] = &A[iDiagLoc[j-1]]; // FORTRAN array indexing in iDiagLoc
    }
    
    return result;
}
Exemplo n.º 4
0
int 
ShadowPetscSOE::setSize(Graph &theGraph)
{
  int n = theGraph.getNumVertex();
  int size = n;
  int N = n;

  // fist itearte through the vertices of the graph to get nnz,
  // the number of non-zeros.
  Vertex *theVertex;
  int NNZ = 0;
  VertexIter &theVertices = theGraph.getVertices();
  while ((theVertex = theVertices()) != 0) {
    const ID &theAdjacency = theVertex->getAdjacency();
    NNZ += theAdjacency.Size() +1; // the +1 is for the diag entry
  }

  // we determine the number of non-zeros & number of nonzero's
  // in each row of A
  // create two integer arrays. One containing the column indices for each
  // entry in A stored in order by row and column number. 
  // The other a pointer into this array for each row.
  int *rowStartA = new int[size];  // start locations of rows in colA
  int *colA = new int[NNZ]; // column locations, stored row-wise

  // fill in rowStartA and colA
  if (size != 0) {
    rowStartA[0] = 0;
    int startLoc = 0;
    int lastLoc = 0;
    for (int a=0; a<size; a++) {

      theVertex = theGraph.getVertexPtr(a);
      if (theVertex == 0) {
	opserr << "WARNING:SparseGenLinSOE::setSize :";
	opserr << " vertex " << a << " not in graph! - size set to 0\n";
	size = 0;
	return -1;
      }

      colA[lastLoc++] = theVertex->getTag(); // place diag in first
      const ID &theAdjacency = theVertex->getAdjacency();
      int idSize = theAdjacency.Size();
	
      // now we have to place the entries in the ID into order in colA
      for (int i=0; i<idSize; i++) {

	int row = theAdjacency(i);
	bool foundPlace = false;
	// find a place in colA for current col
	for (int j=startLoc; j<lastLoc; j++)
	  if (colA[j] > row) { 
	    // move the entries already there one further on
	    // and place col in current location
	    for (int k=lastLoc; k>j; k--)
	      colA[k] = colA[k-1];
	    colA[j] = row;
	    foundPlace = true;
	      j = lastLoc;
	  }
	if (foundPlace == false) // put in at the end
	  colA[lastLoc] = row;

	lastLoc++;
      }
      rowStartA[a+1] = lastLoc;;	    
      startLoc = lastLoc;
    }
  }

  // now for each of the SOE's we determine how to invoke setSizeParallel()
  sendData[0] = 2;
  MPI_Bcast(sendBuffer, 3, MPI_INT, 0, PETSC_COMM_WORLD);

  int dnz = 0;
  int onz = 0;
  int numRowsTyp = N/blockSize/numProcessors;
  numRowsTyp *= blockSize;
  int numRowsLast = numRowsTyp + (N - numProcessors*numRowsTyp); // lastProc extra rows
  int *dnnz = new int[numRowsLast];
  int *onnz = new int[numRowsLast];

  // first determine start and end rows of diag block
  int numRows = numRowsLast;
  int endRow = N-1;
  int startRow = endRow-numRows +1; 
  
  int result = 0;
  // we have to go last to first
  for (int i=numProcessors-1; i>=0; i--) {

    // for each processor determine onnz and dnnz info in colA[]
    for (int j=0; j<numRows; j++) {
      dnnz[j] = 0; 
      onnz[j] = 0;
      int rowStart = rowStartA[startRow+j];
      int nextRowStart = rowStartA[startRow+j+1];
      for (int k=rowStart; k<nextRowStart; k++) {
	int col = colA[k];
	if (col < startRow || col > endRow)
	  onnz[j] += 1;
	else
	  dnnz[j] += 1;
      }
    }

    // now send the data
    if (i != 0) {  // remote object
      sendData[0] = 2;
      sendData[1] = numRows;
      sendData[2] = n;
      
      int tag = 100;
      MPI_Send(sendBuffer, 3, MPI_INT, i, tag, PETSC_COMM_WORLD);

      tag = 101;
      void *buffer = (void *)dnnz;
      MPI_Send(buffer, numRows, MPI_INT, i, tag, PETSC_COMM_WORLD);

      tag = 102;
      buffer = (void *)onnz;
      MPI_Send(buffer, numRows, MPI_INT, i, tag, PETSC_COMM_WORLD);
    }      

    // increment startRow, endRow and numRows if necessary
    endRow = startRow-1;
    numRows = numRowsTyp;
    startRow = endRow-numRows +1; 
  }

  // we broadcast again before we start setSizeParallel()
  // this is because Petsc all processes need to call setup at same
  // time .. if don't we hang
  MPI_Bcast(sendBuffer, 3, MPI_INT, 0, PETSC_COMM_WORLD);  
  result = theSOE.setSizeParallel(numRows, n, dnz, dnnz, onz, onnz);

  // free up trhe memory used
  /*
  delete [] colA;
  delete [] rowStartA;
  delete [] onnz;
  delete [] dnnz;
  */

  return result;
}
Exemplo n.º 5
0
int 
SparseGenRowLinSOE::setSize(Graph &theGraph)
{

    int result = 0;
    int oldSize = size;
    size = theGraph.getNumVertex();

    // fist itearte through the vertices of the graph to get nnz
    Vertex *theVertex;
    int newNNZ = 0;
    VertexIter &theVertices = theGraph.getVertices();
    while ((theVertex = theVertices()) != 0) {
	const ID &theAdjacency = theVertex->getAdjacency();
	newNNZ += theAdjacency.Size() +1; // the +1 is for the diag entry
    }
    nnz = newNNZ;

    if (newNNZ > Asize) { // we have to get more space for A and colA
	if (A != 0) 
	    delete [] A;
	if (colA != 0)
	    delete [] colA;
	
	A = new double[newNNZ];
	colA = new int[newNNZ];
	
        if (A == 0 || colA == 0) {
            opserr << "WARNING SparseGenRowLinSOE::SparseGenRowLinSOE :";
	    opserr << " ran out of memory for A and colA with nnz = ";
	    opserr << newNNZ << " \n";
	    size = 0; Asize = 0; nnz = 0;
	    result =  -1;
        } 
	
	Asize = newNNZ;
    }

    // zero the matrix
    for (int i=0; i<Asize; i++)
	A[i] = 0;
	
    factored = false;
    
    if (size > Bsize) { // we have to get space for the vectors
	
	// delete the old	
	if (B != 0) delete [] B;
	if (X != 0) delete [] X;
	if (rowStartA != 0) delete [] rowStartA;

	// create the new
	B = new double[size];
	X = new double[size];
	rowStartA = new int[size+1]; 
	
        if (B == 0 || X == 0 || rowStartA == 0) {
            opserr << "WARNING SparseGenRowLinSOE::SparseGenRowLinSOE :";
	    opserr << " ran out of memory for vectors (size) (";
	    opserr << size << ") \n";
	    size = 0; Bsize = 0;
	    result =  -1;
        }
	else
	    Bsize = size;
    }

    // zero the vectors
    for (int j=0; j<size; j++) {
	B[j] = 0;
	X[j] = 0;
    }
    
    // create new Vectors objects
    if (size != oldSize) {
	if (vectX != 0)
	    delete vectX;

	if (vectB != 0)
	    delete vectB;
	
	vectX = new Vector(X,size);
	vectB = new Vector(B,size);	
    }

    // fill in rowStartA and colA
    if (size != 0) {
      rowStartA[0] = 0;
      int startLoc = 0;
      int lastLoc = 0;
      for (int a=0; a<size; a++) {

	theVertex = theGraph.getVertexPtr(a);
	if (theVertex == 0) {
	  opserr << "WARNING:SparseGenRowLinSOE::setSize :";
	  opserr << " vertex " << a << " not in graph! - size set to 0\n";
	  size = 0;
	  return -1;
	}

	colA[lastLoc++] = theVertex->getTag(); // place diag in first
	const ID &theAdjacency = theVertex->getAdjacency();
	int idSize = theAdjacency.Size();
	
	// now we have to place the entries in the ID into order in colA
	for (int i=0; i<idSize; i++) {

	  int row = theAdjacency(i);
	  bool foundPlace = false;
	  // find a place in colA for current col
	  for (int j=startLoc; j<lastLoc; j++)
	    if (colA[j] > row) { 
	      // move the entries already there one further on
	      // and place col in current location
	      for (int k=lastLoc; k>j; k--)
		
		colA[k] = colA[k-1];
	      colA[j] = row;
	      foundPlace = true;
	      j = lastLoc;
	    }
	  if (foundPlace == false) // put in at the end
	    colA[lastLoc] = row;

	  lastLoc++;
	}
	rowStartA[a+1] = lastLoc;;	    
	startLoc = lastLoc;
      }
    }

    // invoke setSize() on the Solver   
     LinearSOESolver *the_Solver = this->getSolver();
    int solverOK = the_Solver->setSize();
    if (solverOK < 0) {
	opserr << "WARNING:SparseGenRowLinSOE::setSize :";
	opserr << " solver failed setSize()\n";
	return solverOK;
    }    
    return result;
}
Exemplo n.º 6
0
int 
BandGenLinSOE::setSize(Graph &theGraph)
{
    int result = 0;
    int oldSize = size;
    size = theGraph.getNumVertex();
    
    /*
     * determine the number of superdiagonals and subdiagonals
     */
    
    numSubD = 0;
    numSuperD = 0;

    Vertex *vertexPtr;
    VertexIter &theVertices = theGraph.getVertices();
    
    while ((vertexPtr = theVertices()) != 0) {
	int vertexNum = vertexPtr->getTag();
	const ID &theAdjacency = vertexPtr->getAdjacency();
	for (int i=0; i<theAdjacency.Size(); i++) {
	    int otherNum = theAdjacency(i);
	    int diff = vertexNum - otherNum;
	    if (diff > 0) {
		if (diff > numSuperD)
		    numSuperD = diff;
	    } else 
		if (diff < numSubD)
		    numSubD = diff;
	}
    }
    numSubD *= -1;

    int newSize = size * (2*numSubD + numSuperD +1);
    if (newSize > Asize) { // we have to get another space for A

	if (A != 0) 
	    delete [] A;

	A = new double[newSize];
	
        if (A == 0) {
            opserr << "WARNING BandGenLinSOE::BandGenLinSOE :";
	    opserr << " ran out of memory for A (size,super,sub) (";
	    opserr << size <<", " << numSuperD << ", " << numSubD << ") \n";
	    Asize = 0; size = 0; numSubD = 0; numSuperD = 0;
	    result= -1;
        }
	else  
	    Asize = newSize;
    }

    // zero the matrix
    for (int i=0; i<Asize; i++)
	A[i] = 0;
	
    factored = false;
    
    if (size > Bsize) { // we have to get space for the vectors
	
	// delete the old	
	if (B != 0) delete [] B;
	if (X != 0) delete [] X;

	// create the new
	B = new double[size];
	X = new double[size];
	
        if (B == 0 || X == 0) {
            opserr << "WARNING BandGenLinSOE::BandGenLinSOE :";
	    opserr << " ran out of memory for vectors (size) (";
	    opserr << size << ") \n";
	    Bsize = 0; size = 0; numSubD = 0; numSuperD = 0;
	    result = -1;
        }
	else 
	    Bsize = size;
    }

    // zero the vectors
    for (int j=0; j<size; j++) {
	B[j] = 0;
	X[j] = 0;
    }

    // get new Vector objects if size has changes
    if (oldSize != size) {
	if (vectX != 0) 
	    delete vectX;

	if (vectB != 0) 
	    delete vectB;
		
	vectX = new Vector(X,size);
	vectB = new Vector(B,size);
    }
    
    // invoke setSize() on the Solver
    LinearSOESolver *theSolvr = this->getSolver();
    int solverOK = theSolvr->setSize();
    if (solverOK < 0) {
	opserr << "WARNING:BandGenLinSOE::setSize :";
	opserr << " solver failed setSize()\n";
	return solverOK;
    }    

    return result;    
}
Exemplo n.º 7
0
int 
MumpsParallelSOE::setSize(Graph &theGraph)
{
  int result = 0;
  int oldSize = size;
  int maxNumSubVertex = 0;
  
  // fist itearte through the vertices of the graph to get nnzLoc and n
  int maxVertexTag = -1;
  Vertex *theVertex;
  int newNNZ = 0;
  size = theGraph.getNumVertex();
  int mySize = size;
  //opserr << "MumpsParallelSOE: size : " << size << endln;

  VertexIter &theVertices = theGraph.getVertices();
  while ((theVertex = theVertices()) != 0) {
    int vertexTag = theVertex->getTag();
    if (vertexTag > maxVertexTag)
      maxVertexTag = vertexTag;
    const ID &theAdjacency = theVertex->getAdjacency();
    newNNZ += theAdjacency.Size() +1; // the +1 is for the diag entry
  }

  if (matType !=  0) {

    // symmetric - allows us to reduce nnz by almost half
    newNNZ -= size;
    newNNZ /= 2;
    newNNZ += size;
  }

  nnz = newNNZ;

  if (processID != 0) {

    //
    // if subprocess, send local max vertexTag (n)
    // recv ax n from P0
    //
    static ID data(1);

    data(0) = maxVertexTag;
    Channel *theChannel = theChannels[0];
    theChannel->sendID(0, 0, data);
    theChannel->recvID(0, 0, data);
    
    size = data(0);

  } else {

    //
    // from each distributed soe recv it's max n and compare; return max n to all
    //

    static ID data(1);
    FEM_ObjectBroker theBroker;
    for (int j=0; j<numChannels; j++) {
      Channel *theChannel = theChannels[j];
      theChannel->recvID(0, 0, data);
      if (data(0) > maxVertexTag)
	maxVertexTag = data(0);
    }

    data(0) = maxVertexTag;

    for (int j=0; j<numChannels; j++) {
      Channel *theChannel = theChannels[j];
      theChannel->sendID(0, 0, data);
    }
    size = maxVertexTag;
  }

  size+=1; // vertices numbered 0 through n-1

  if (nnz > Asize) { // we have to get more space for A and rowA and colA

    if (A != 0) delete [] A;
    if (rowA != 0) delete [] rowA;
    if (colA != 0) delete [] colA;

    A = new double[nnz];    
    rowA = new int[nnz];
    colA = new int[nnz];
      
    for (int i=0; i<nnz; i++) {
      A[i]=0;
      rowA[i]=0;
      colA[i]=0;
    }

    if (rowA == 0 || A == 0 || colA == 0) {
      opserr << "WARNING SparseGenColLinSOE::SparseGenColLinSOE :";
      opserr << " ran out of memory for A and rowA with nnz = ";
      opserr << nnz << " \n";
      size = 0; Asize = 0; nnz = 0;
	result =  -1;
    } 
    Asize = nnz;
  }

  
  if (size > Bsize) { // we have to get space for the vectors

    if (B != 0) delete [] B;
    if (X != 0) delete [] X;    
    if (myB != 0) delete [] myB;
    if (workArea != 0) delete [] workArea;
    if (colStartA != 0)  delete [] colStartA;

    // create the new
    B = new double[size];
    X = new double[size];
    myB = new double[size];
    workArea = new double[size];
    colStartA = new int[size+1]; 
    
    if (B == 0 || X == 0 || colStartA == 0 || workArea == 0 || myB == 0) {
      opserr << "WARNING MumpsSOE::MumpsSOE :";
      opserr << " ran out of memory for vectors (size) (";
      opserr << size << ") \n";
      size = 0; Bsize = 0;
      result =  -1;
    }
    else
      Bsize = size;

  }
  
  // zero the vectors
  for (int j=0; j<size; j++) {
    B[j] = 0;
    X[j] = 0;
    myB[j] = 0;
  }

  // create new Vectors objects
  if (size != oldSize) {
    if (vectX != 0) delete vectX;
    if (vectB != 0) delete vectB;
    if (myVectB != 0) delete myVectB;
    
    vectX = new Vector(X,size);
    vectB = new Vector(B,size);	
    myVectB = new Vector(myB, size);
  }

  // fill in colStartA and rowA
  if (size != 0) {
    colStartA[0] = 0;
    int startLoc = 0;
    int lastLoc = 0;
    for (int a=0; a<size; a++) {
      
      theVertex = theGraph.getVertexPtr(a);
      if (theVertex != 0) {
	
	int vertexTag = theVertex->getTag();
	rowA[lastLoc++] = vertexTag; // place diag in first
	const ID &theAdjacency = theVertex->getAdjacency();
	int idSize = theAdjacency.Size();
	
	// now we have to place the entries in the ID into order in rowA
	
	if (matType != 0) {
	  
	  // symmetric
	  for (int i=0; i<idSize; i++) {
	    int row = theAdjacency(i);
	    if (row > vertexTag) {
	      bool foundPlace = false;
	      // find a place in rowA for current col
	      for (int j=startLoc; j<lastLoc; j++)
		if (rowA[j] > row) { 
		  // move the entries already there one further on
		  // and place col in current location
		  for (int k=lastLoc; k>j; k--)
		    rowA[k] = rowA[k-1];
		  rowA[j] = row;
		  foundPlace = true;
		  j = lastLoc;
		}
	      
	      if (foundPlace == false) // put in at the end
		rowA[lastLoc] = row;
	      lastLoc++;
	    }
	  }
	  
	} else {

	  // unsymmetric	  
	  for (int i=0; i<idSize; i++) {
	    int row = theAdjacency(i);
	    bool foundPlace = false;
	    // find a place in rowA for current col
	    for (int j=startLoc; j<lastLoc; j++)
	      if (rowA[j] > row) { 
		// move the entries already there one further on
		// and place col in current location
		for (int k=lastLoc; k>j; k--)
		  rowA[k] = rowA[k-1];
		rowA[j] = row;
		foundPlace = true;
		j = lastLoc;
	      }
	    if (foundPlace == false) // put in at the end
	      rowA[lastLoc] = row;
	    
	    lastLoc++;
	  }
	}
      }
      colStartA[a+1] = lastLoc;
      startLoc = lastLoc;
    }
  }

  // fill in colA
  int count = 0;
  for (int i=0; i<size; i++) {
    for (int k=colStartA[i]; k<colStartA[i+1]; k++)
      colA[count++] = i;
  }

  LinearSOESolver *theSolvr = this->getSolver();

  int solverOK = theSolvr->setSize();
  if (solverOK < 0) {
    opserr << "WARNING:MumpsParallelSOE::setSize :";
    opserr << " solver failed setSize()\n";
    return solverOK;
  }    

  return result;    
}