Actual source code: aijbaij.c

  1: /*$Id: aijbaij.c,v 1.9 2001/08/07 03:02:55 balay Exp $*/

 3:  #include src/mat/impls/baij/seq/baij.h

  5: EXTERN_C_BEGIN
  8: int MatConvert_SeqBAIJ_SeqAIJ(Mat A,MatType newtype,Mat *newmat) {
  9:   Mat          B;
 10:   Mat_SeqBAIJ  *a = (Mat_SeqBAIJ*)A->data;
 11:   int          ierr,bs = a->bs,*ai = a->i,*aj = a->j,n = A->M/bs,i,j,k;
 12:   int          *rowlengths,*rows,*cols,maxlen = 0,ncols;
 13:   PetscScalar  *aa = a->a;

 16:   PetscMalloc(n*bs*sizeof(int),&rowlengths);
 17:   for (i=0; i<n; i++) {
 18:     maxlen = PetscMax(maxlen,(ai[i+1] - ai[i]));
 19:     for (j=0; j<bs; j++) {
 20:       rowlengths[i*bs+j] = bs*(ai[i+1] - ai[i]);
 21:     }
 22:   }
 23:   MatCreateSeqAIJ(PETSC_COMM_SELF,A->m,A->n,0,rowlengths,&B);
 24:   MatSetOption(B,MAT_COLUMN_ORIENTED);
 25:   MatSetOption(B,MAT_ROWS_SORTED);
 26:   MatSetOption(B,MAT_COLUMNS_SORTED);
 27:   PetscFree(rowlengths);

 29:   PetscMalloc(bs*sizeof(int),&rows);
 30:   PetscMalloc(bs*maxlen*sizeof(int),&cols);
 31:   for (i=0; i<n; i++) {
 32:     for (j=0; j<bs; j++) {
 33:       rows[j] = i*bs+j;
 34:     }
 35:     ncols = ai[i+1] - ai[i];
 36:     for (k=0; k<ncols; k++) {
 37:       for (j=0; j<bs; j++) {
 38:         cols[k*bs+j] = bs*(*aj) + j;
 39:       }
 40:       aj++;
 41:     }
 42:     MatSetValues(B,bs,rows,bs*ncols,cols,aa,INSERT_VALUES);
 43:     aa   += ncols*bs*bs;
 44:   }
 45:   PetscFree(cols);
 46:   PetscFree(rows);
 47:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 48:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);

 50:   /* Fake support for "inplace" convert. */
 51:   if (*newmat == A) {
 52:     MatDestroy(A);
 53:   }
 54:   *newmat = B;
 55:   return(0);
 56: }
 57: EXTERN_C_END

 59:  #include src/mat/impls/aij/seq/aij.h

 61: EXTERN_C_BEGIN
 64: int MatConvert_SeqAIJ_SeqBAIJ(Mat A,MatType newtype,Mat *newmat) {
 65:   Mat         B;
 66:   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)A->data;
 67:   Mat_SeqBAIJ *b;
 68:   int         ierr,*ai=a->i,m=A->M,n=A->N,i,*rowlengths;

 71:   if (n != m) SETERRQ(PETSC_ERR_ARG_WRONG,"Matrix must be square");

 73:   PetscMalloc(m*sizeof(int),&rowlengths);
 74:   for (i=0; i<m; i++) {
 75:     rowlengths[i] = ai[i+1] - ai[i];
 76:   }
 77:   MatCreateSeqBAIJ(PETSC_COMM_SELF,1,m,n,0,rowlengths,&B);
 78:   PetscFree(rowlengths);

 80:   MatSetOption(B,MAT_ROW_ORIENTED);
 81:   MatSetOption(B,MAT_ROWS_SORTED);
 82:   MatSetOption(B,MAT_COLUMNS_SORTED);
 83: 
 84:   b  = (Mat_SeqBAIJ*)(B->data);

 86:   PetscMemcpy(b->i,a->i,(m+1)*sizeof(int));
 87:   PetscMemcpy(b->ilen,a->ilen,m*sizeof(int));
 88:   PetscMemcpy(b->j,a->j,a->nz*sizeof(int));
 89:   PetscMemcpy(b->a,a->a,a->nz*sizeof(MatScalar));
 90: 
 91:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 92:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);

 94:   /* Fake support for "inplace" convert. */
 95:   if (*newmat == A) {
 96:     MatDestroy(A);
 97:   }
 98:   *newmat = B;
 99:   return(0);
100: }
101: EXTERN_C_END