1a30f8f8cSSatish Balay 2a30f8f8cSSatish Balay /* 3a30f8f8cSSatish Balay Support for the parallel SBAIJ matrix vector multiply 4a30f8f8cSSatish Balay */ 5c6db04a5SJed Brown #include <../src/mat/impls/sbaij/mpi/mpisbaij.h> 67cff1425SSatish Balay 7dfbe8321SBarry Smith PetscErrorCode MatSetUpMultiply_MPISBAIJ(Mat mat) 8a30f8f8cSSatish Balay { 940781036SHong Zhang Mat_MPISBAIJ *sbaij = (Mat_MPISBAIJ*)mat->data; 1040781036SHong Zhang Mat_SeqBAIJ *B = (Mat_SeqBAIJ*)(sbaij->B->data); 116849ba73SBarry Smith PetscErrorCode ierr; 1213f74950SBarry Smith PetscInt Nbs = sbaij->Nbs,i,j,*indices,*aj = B->j,ec = 0,*garray,*sgarray; 13d0f46423SBarry Smith PetscInt bs = mat->rmap->bs,*stmp,mbs=sbaij->mbs, vec_size,nt; 1440781036SHong Zhang IS from,to; 1540781036SHong Zhang Vec gvec; 1613f74950SBarry Smith PetscMPIInt rank =sbaij->rank,lsize,size=sbaij->size; 17077aedafSJed Brown PetscInt *owners=sbaij->rangebs,*ec_owner,k; 18077aedafSJed Brown const PetscInt *sowners; 1940781036SHong Zhang PetscScalar *ptr; 2040781036SHong Zhang 2140781036SHong Zhang PetscFunctionBegin; 226bf464f9SBarry Smith ierr = VecScatterDestroy(&sbaij->sMvctx);CHKERRQ(ierr); 2340781036SHong Zhang 2440781036SHong Zhang /* For the first stab we make an array as long as the number of columns */ 2540781036SHong Zhang /* mark those columns that are in sbaij->B */ 261795a4d1SJed Brown ierr = PetscCalloc1(Nbs,&indices);CHKERRQ(ierr); 2740781036SHong Zhang for (i=0; i<mbs; i++) { 2840781036SHong Zhang for (j=0; j<B->ilen[i]; j++) { 2940781036SHong Zhang if (!indices[aj[B->i[i] + j]]) ec++; 3040781036SHong Zhang indices[aj[B->i[i] + j]] = 1; 3140781036SHong Zhang } 3240781036SHong Zhang } 3340781036SHong Zhang 3440781036SHong Zhang /* form arrays of columns we need */ 35785e854fSJed Brown ierr = PetscMalloc1(ec,&garray);CHKERRQ(ierr); 36dcca6d9dSJed Brown ierr = PetscMalloc2(2*ec,&sgarray,ec,&ec_owner);CHKERRQ(ierr); 3740781036SHong Zhang 3840781036SHong Zhang ec = 0; 3940781036SHong Zhang for (j=0; j<size; j++) { 4040781036SHong Zhang for (i=owners[j]; i<owners[j+1]; i++) { 4140781036SHong Zhang if (indices[i]) { 4240781036SHong Zhang garray[ec] = i; 4340781036SHong Zhang ec_owner[ec] = j; 4440781036SHong Zhang ec++; 4540781036SHong Zhang } 4640781036SHong Zhang } 4740781036SHong Zhang } 4840781036SHong Zhang 4940781036SHong Zhang /* make indices now point into garray */ 50b424e231SHong Zhang for (i=0; i<ec; i++) indices[garray[i]] = i; 5140781036SHong Zhang 5240781036SHong Zhang /* compact out the extra columns in B */ 5340781036SHong Zhang for (i=0; i<mbs; i++) { 5440781036SHong Zhang for (j=0; j<B->ilen[i]; j++) aj[B->i[i] + j] = indices[aj[B->i[i] + j]]; 5540781036SHong Zhang } 5640781036SHong Zhang B->nbs = ec; 57*ca5434daSLawrence Mitchell ierr = PetscLayoutDestroy(&sbaij->B->cmap);CHKERRQ(ierr); 58*ca5434daSLawrence Mitchell ierr = PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)sbaij->B),ec*mat->rmap->bs,ec*mat->rmap->bs,mat->rmap->bs,&sbaij->B->cmap);CHKERRQ(ierr); 5940781036SHong Zhang ierr = PetscFree(indices);CHKERRQ(ierr); 6040781036SHong Zhang 6140781036SHong Zhang /* create local vector that is used to scatter into */ 6240781036SHong Zhang ierr = VecCreateSeq(PETSC_COMM_SELF,ec*bs,&sbaij->lvec);CHKERRQ(ierr); 6340781036SHong Zhang 6440781036SHong Zhang /* create two temporary index sets for building scatter-gather */ 65785e854fSJed Brown ierr = PetscMalloc1(2*ec,&stmp);CHKERRQ(ierr); 66deff0451SBarry Smith ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,PETSC_COPY_VALUES,&from);CHKERRQ(ierr); 6726fbe8dcSKarl Rupp for (i=0; i<ec; i++) stmp[i] = i; 68deff0451SBarry Smith ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,PETSC_COPY_VALUES,&to);CHKERRQ(ierr); 6940781036SHong Zhang 7093d1592dSHong Zhang /* generate the scatter context 714c500f23SPierre Jolivet -- Mvctx and lvec are not used by MatMult_MPISBAIJ(), but usefull for some applications */ 72ce94432eSBarry Smith ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat),1,mat->cmap->n,mat->cmap->N,NULL,&gvec);CHKERRQ(ierr); 739448b7f1SJunchao Zhang ierr = VecScatterCreate(gvec,from,sbaij->lvec,to,&sbaij->Mvctx);CHKERRQ(ierr); 746bf464f9SBarry Smith ierr = VecDestroy(&gvec);CHKERRQ(ierr); 7540781036SHong Zhang 7640781036SHong Zhang sbaij->garray = garray; 7726fbe8dcSKarl Rupp 783bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->Mvctx);CHKERRQ(ierr); 793bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->lvec);CHKERRQ(ierr); 803bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)from);CHKERRQ(ierr); 813bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)to);CHKERRQ(ierr); 8240781036SHong Zhang 836bf464f9SBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 846bf464f9SBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 8540781036SHong Zhang 8640781036SHong Zhang /* create parallel vector that is used by SBAIJ matrix to scatter from/into */ 8740781036SHong Zhang lsize = (mbs + ec)*bs; 88ce94432eSBarry Smith ierr = VecCreateMPI(PetscObjectComm((PetscObject)mat),lsize,PETSC_DETERMINE,&sbaij->slvec0);CHKERRQ(ierr); 8940781036SHong Zhang ierr = VecDuplicate(sbaij->slvec0,&sbaij->slvec1);CHKERRQ(ierr); 9040781036SHong Zhang ierr = VecGetSize(sbaij->slvec0,&vec_size);CHKERRQ(ierr); 9140781036SHong Zhang 92077aedafSJed Brown ierr = VecGetOwnershipRanges(sbaij->slvec0,&sowners);CHKERRQ(ierr); 9340781036SHong Zhang 9440781036SHong Zhang /* x index in the IS sfrom */ 9540781036SHong Zhang for (i=0; i<ec; i++) { 9640781036SHong Zhang j = ec_owner[i]; 9740781036SHong Zhang sgarray[i] = garray[i] + (sowners[j]/bs - owners[j]); 9840781036SHong Zhang } 9940781036SHong Zhang /* b index in the IS sfrom */ 10040781036SHong Zhang k = sowners[rank]/bs + mbs; 10140781036SHong Zhang for (i=ec,j=0; i< 2*ec; i++,j++) sgarray[i] = k + j; 102deff0451SBarry Smith ierr = ISCreateBlock(PETSC_COMM_SELF,bs,2*ec,sgarray,PETSC_COPY_VALUES,&from);CHKERRQ(ierr); 10340781036SHong Zhang 10440781036SHong Zhang /* x index in the IS sto */ 10540781036SHong Zhang k = sowners[rank]/bs + mbs; 106e82e9f6bSBarry Smith for (i=0; i<ec; i++) stmp[i] = (k + i); 10740781036SHong Zhang /* b index in the IS sto */ 108e82e9f6bSBarry Smith for (i=ec; i<2*ec; i++) stmp[i] = sgarray[i-ec]; 10940781036SHong Zhang 110deff0451SBarry Smith ierr = ISCreateBlock(PETSC_COMM_SELF,bs,2*ec,stmp,PETSC_COPY_VALUES,&to);CHKERRQ(ierr); 11140781036SHong Zhang 1129448b7f1SJunchao Zhang ierr = VecScatterCreate(sbaij->slvec0,from,sbaij->slvec1,to,&sbaij->sMvctx);CHKERRQ(ierr); 11340781036SHong Zhang 11440781036SHong Zhang ierr = VecGetLocalSize(sbaij->slvec1,&nt);CHKERRQ(ierr); 1151ebc52fbSHong Zhang ierr = VecGetArray(sbaij->slvec1,&ptr);CHKERRQ(ierr); 116778a2246SBarry Smith ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,bs*mbs,ptr,&sbaij->slvec1a);CHKERRQ(ierr); 117778a2246SBarry Smith ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,nt-bs*mbs,ptr+bs*mbs,&sbaij->slvec1b);CHKERRQ(ierr); 1181ebc52fbSHong Zhang ierr = VecRestoreArray(sbaij->slvec1,&ptr);CHKERRQ(ierr); 11940781036SHong Zhang 1201ebc52fbSHong Zhang ierr = VecGetArray(sbaij->slvec0,&ptr);CHKERRQ(ierr); 121778a2246SBarry Smith ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,nt-bs*mbs,ptr+bs*mbs,&sbaij->slvec0b);CHKERRQ(ierr); 1221ebc52fbSHong Zhang ierr = VecRestoreArray(sbaij->slvec0,&ptr);CHKERRQ(ierr); 12340781036SHong Zhang 12440781036SHong Zhang ierr = PetscFree(stmp);CHKERRQ(ierr); 125ce94432eSBarry Smith ierr = MPI_Barrier(PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 12640781036SHong Zhang 1273bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->sMvctx);CHKERRQ(ierr); 1283bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec0);CHKERRQ(ierr); 1293bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec1);CHKERRQ(ierr); 1303bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec0b);CHKERRQ(ierr); 1313bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec1a);CHKERRQ(ierr); 1323bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec1b);CHKERRQ(ierr); 1333bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)from);CHKERRQ(ierr); 1343bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)to);CHKERRQ(ierr); 13540781036SHong Zhang 1363bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)mat,(ec+1)*sizeof(PetscInt));CHKERRQ(ierr); 1376bf464f9SBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 1386bf464f9SBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 13974ed9c26SBarry Smith ierr = PetscFree2(sgarray,ec_owner);CHKERRQ(ierr); 14040781036SHong Zhang PetscFunctionReturn(0); 14140781036SHong Zhang } 14240781036SHong Zhang 143dfbe8321SBarry Smith PetscErrorCode MatSetUpMultiply_MPISBAIJ_2comm(Mat mat) 14440781036SHong Zhang { 1452896befeSSatish Balay Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data; 146a30f8f8cSSatish Balay Mat_SeqBAIJ *B = (Mat_SeqBAIJ*)(baij->B->data); 1476849ba73SBarry Smith PetscErrorCode ierr; 14813f74950SBarry Smith PetscInt i,j,*aj = B->j,ec = 0,*garray; 149d0f46423SBarry Smith PetscInt bs = mat->rmap->bs,*stmp; 150a30f8f8cSSatish Balay IS from,to; 151a30f8f8cSSatish Balay Vec gvec; 152a30f8f8cSSatish Balay #if defined(PETSC_USE_CTABLE) 153a30f8f8cSSatish Balay PetscTable gid1_lid1; 154a30f8f8cSSatish Balay PetscTablePosition tpos; 15513f74950SBarry Smith PetscInt gid,lid; 1566f531f54SSatish Balay #else 15713f74950SBarry Smith PetscInt Nbs = baij->Nbs,*indices; 158a30f8f8cSSatish Balay #endif 159a30f8f8cSSatish Balay 160a30f8f8cSSatish Balay PetscFunctionBegin; 161a30f8f8cSSatish Balay #if defined(PETSC_USE_CTABLE) 162a30f8f8cSSatish Balay /* use a table - Mark Adams */ 163e23dfa41SBarry Smith PetscTableCreate(B->mbs,baij->Nbs+1,&gid1_lid1); 164a30f8f8cSSatish Balay for (i=0; i<B->mbs; i++) { 165a30f8f8cSSatish Balay for (j=0; j<B->ilen[i]; j++) { 16613f74950SBarry Smith PetscInt data,gid1 = aj[B->i[i]+j] + 1; 167a30f8f8cSSatish Balay ierr = PetscTableFind(gid1_lid1,gid1,&data);CHKERRQ(ierr); 168a30f8f8cSSatish Balay if (!data) { 169a30f8f8cSSatish Balay /* one based table */ 1703861aac3SJed Brown ierr = PetscTableAdd(gid1_lid1,gid1,++ec,INSERT_VALUES);CHKERRQ(ierr); 171a30f8f8cSSatish Balay } 172a30f8f8cSSatish Balay } 173a30f8f8cSSatish Balay } 174a30f8f8cSSatish Balay /* form array of columns we need */ 175785e854fSJed Brown ierr = PetscMalloc1(ec,&garray);CHKERRQ(ierr); 176a30f8f8cSSatish Balay ierr = PetscTableGetHeadPosition(gid1_lid1,&tpos);CHKERRQ(ierr); 177a30f8f8cSSatish Balay while (tpos) { 178a30f8f8cSSatish Balay ierr = PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);CHKERRQ(ierr); 179a30f8f8cSSatish Balay gid--; lid--; 180a30f8f8cSSatish Balay garray[lid] = gid; 181a30f8f8cSSatish Balay } 182a30f8f8cSSatish Balay ierr = PetscSortInt(ec,garray);CHKERRQ(ierr); 183a30f8f8cSSatish Balay ierr = PetscTableRemoveAll(gid1_lid1);CHKERRQ(ierr); 184a30f8f8cSSatish Balay for (i=0; i<ec; i++) { 1853861aac3SJed Brown ierr = PetscTableAdd(gid1_lid1,garray[i]+1,i+1,INSERT_VALUES);CHKERRQ(ierr); 186a30f8f8cSSatish Balay } 187a30f8f8cSSatish Balay /* compact out the extra columns in B */ 188a30f8f8cSSatish Balay for (i=0; i<B->mbs; i++) { 189a30f8f8cSSatish Balay for (j=0; j<B->ilen[i]; j++) { 19013f74950SBarry Smith PetscInt gid1 = aj[B->i[i] + j] + 1; 191a30f8f8cSSatish Balay ierr = PetscTableFind(gid1_lid1,gid1,&lid);CHKERRQ(ierr); 192a30f8f8cSSatish Balay lid--; 193a30f8f8cSSatish Balay aj[B->i[i]+j] = lid; 194a30f8f8cSSatish Balay } 195a30f8f8cSSatish Balay } 196a30f8f8cSSatish Balay B->nbs = ec; 197*ca5434daSLawrence Mitchell ierr = PetscLayoutDestroy(&baij->B->cmap);CHKERRQ(ierr); 198*ca5434daSLawrence Mitchell ierr = PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)baij->B),ec*mat->rmap->bs,ec*mat->rmap->bs,mat->rmap->bs,&baij->B->cmap);CHKERRQ(ierr); 1996bc0bbbfSBarry Smith ierr = PetscTableDestroy(&gid1_lid1);CHKERRQ(ierr); 200a30f8f8cSSatish Balay #else 201a30f8f8cSSatish Balay /* For the first stab we make an array as long as the number of columns */ 202a30f8f8cSSatish Balay /* mark those columns that are in baij->B */ 2031795a4d1SJed Brown ierr = PetscCalloc1(Nbs,&indices);CHKERRQ(ierr); 204a30f8f8cSSatish Balay for (i=0; i<B->mbs; i++) { 205a30f8f8cSSatish Balay for (j=0; j<B->ilen[i]; j++) { 206a30f8f8cSSatish Balay if (!indices[aj[B->i[i] + j]]) ec++; 207a30f8f8cSSatish Balay indices[aj[B->i[i] + j]] = 1; 208a30f8f8cSSatish Balay } 209a30f8f8cSSatish Balay } 210a30f8f8cSSatish Balay 211a30f8f8cSSatish Balay /* form array of columns we need */ 212785e854fSJed Brown ierr = PetscMalloc1(ec,&garray);CHKERRQ(ierr); 213f3ef73ceSBarry Smith ec = 0; 214f3ef73ceSBarry Smith for (i=0; i<Nbs; i++) { 215f3ef73ceSBarry Smith if (indices[i]) { 216f3ef73ceSBarry Smith garray[ec++] = i; 217f3ef73ceSBarry Smith } 218f3ef73ceSBarry Smith } 219a30f8f8cSSatish Balay 220a30f8f8cSSatish Balay /* make indices now point into garray */ 22126fbe8dcSKarl Rupp for (i=0; i<ec; i++) indices[garray[i]] = i; 222a30f8f8cSSatish Balay 223a30f8f8cSSatish Balay /* compact out the extra columns in B */ 224a30f8f8cSSatish Balay for (i=0; i<B->mbs; i++) { 225a30f8f8cSSatish Balay for (j=0; j<B->ilen[i]; j++) { 226a30f8f8cSSatish Balay aj[B->i[i] + j] = indices[aj[B->i[i] + j]]; 227a30f8f8cSSatish Balay } 228a30f8f8cSSatish Balay } 229a30f8f8cSSatish Balay B->nbs = ec; 23026fbe8dcSKarl Rupp 231d0f46423SBarry Smith baij->B->cmap->n = ec*mat->rmap->bs; 23226fbe8dcSKarl Rupp 233a30f8f8cSSatish Balay ierr = PetscFree(indices);CHKERRQ(ierr); 234a30f8f8cSSatish Balay #endif 235633e10c7SHong Zhang 236a30f8f8cSSatish Balay /* create local vector that is used to scatter into */ 237a30f8f8cSSatish Balay ierr = VecCreateSeq(PETSC_COMM_SELF,ec*bs,&baij->lvec);CHKERRQ(ierr); 238a30f8f8cSSatish Balay 239a30f8f8cSSatish Balay /* create two temporary index sets for building scatter-gather */ 240deff0451SBarry Smith ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,PETSC_COPY_VALUES,&from);CHKERRQ(ierr); 241a30f8f8cSSatish Balay 242785e854fSJed Brown ierr = PetscMalloc1(ec,&stmp);CHKERRQ(ierr); 24326fbe8dcSKarl Rupp for (i=0; i<ec; i++) stmp[i] = i; 244deff0451SBarry Smith ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,PETSC_OWN_POINTER,&to);CHKERRQ(ierr); 245a30f8f8cSSatish Balay 246a30f8f8cSSatish Balay /* create temporary global vector to generate scatter context */ 247ce94432eSBarry Smith ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat),1,mat->cmap->n,mat->cmap->N,NULL,&gvec);CHKERRQ(ierr); 2489448b7f1SJunchao Zhang ierr = VecScatterCreate(gvec,from,baij->lvec,to,&baij->Mvctx);CHKERRQ(ierr); 2496bf464f9SBarry Smith ierr = VecDestroy(&gvec);CHKERRQ(ierr); 250a30f8f8cSSatish Balay 2513bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)baij->Mvctx);CHKERRQ(ierr); 2523bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)baij->lvec);CHKERRQ(ierr); 2533bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)from);CHKERRQ(ierr); 2543bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)to);CHKERRQ(ierr); 25526fbe8dcSKarl Rupp 256a30f8f8cSSatish Balay baij->garray = garray; 25726fbe8dcSKarl Rupp 2583bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)mat,(ec+1)*sizeof(PetscInt));CHKERRQ(ierr); 2596bf464f9SBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 2606bf464f9SBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 261a30f8f8cSSatish Balay PetscFunctionReturn(0); 262a30f8f8cSSatish Balay } 263a30f8f8cSSatish Balay 264a30f8f8cSSatish Balay /* 26501b2bd88SHong Zhang Takes the local part of an already assembled MPISBAIJ matrix 266a30f8f8cSSatish Balay and disassembles it. This is to allow new nonzeros into the matrix 267a30f8f8cSSatish Balay that require more communication in the matrix vector multiply. 268a30f8f8cSSatish Balay Thus certain data-structures must be rebuilt. 269a30f8f8cSSatish Balay 270a30f8f8cSSatish Balay Kind of slow! But that's what application programmers get when 271a30f8f8cSSatish Balay they are sloppy. 272a30f8f8cSSatish Balay */ 273ab9863d7SBarry Smith PetscErrorCode MatDisAssemble_MPISBAIJ(Mat A) 274a30f8f8cSSatish Balay { 2752896befeSSatish Balay Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)A->data; 276a30f8f8cSSatish Balay Mat B = baij->B,Bnew; 277a30f8f8cSSatish Balay Mat_SeqBAIJ *Bbaij = (Mat_SeqBAIJ*)B->data; 278dfbe8321SBarry Smith PetscErrorCode ierr; 279d0f46423SBarry Smith PetscInt i,j,mbs=Bbaij->mbs,n = A->cmap->N,col,*garray=baij->garray; 280d0f46423SBarry Smith PetscInt k,bs=A->rmap->bs,bs2=baij->bs2,*rvals,*nz,ec,m=A->rmap->n; 281a30f8f8cSSatish Balay MatScalar *a = Bbaij->a; 28287828ca2SBarry Smith PetscScalar *atmp; 283ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_MAT_SINGLE) 28413f74950SBarry Smith PetscInt l; 285a30f8f8cSSatish Balay #endif 286a30f8f8cSSatish Balay 287a30f8f8cSSatish Balay PetscFunctionBegin; 288ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_MAT_SINGLE) 289785e854fSJed Brown ierr = PetscMalloc1(A->rmap->bs,&atmp); 29082502324SSatish Balay #endif 291a30f8f8cSSatish Balay /* free stuff related to matrix-vec multiply */ 292b0a32e0cSBarry Smith ierr = VecGetSize(baij->lvec,&ec);CHKERRQ(ierr); /* needed for PetscLogObjectMemory below */ 2936bf464f9SBarry Smith ierr = VecDestroy(&baij->lvec);CHKERRQ(ierr); 2946bf464f9SBarry Smith ierr = VecScatterDestroy(&baij->Mvctx);CHKERRQ(ierr); 29501b2bd88SHong Zhang 2966bf464f9SBarry Smith ierr = VecDestroy(&baij->slvec0);CHKERRQ(ierr); 2976bf464f9SBarry Smith ierr = VecDestroy(&baij->slvec0b);CHKERRQ(ierr); 2986bf464f9SBarry Smith ierr = VecDestroy(&baij->slvec1);CHKERRQ(ierr); 2996bf464f9SBarry Smith ierr = VecDestroy(&baij->slvec1a);CHKERRQ(ierr); 3006bf464f9SBarry Smith ierr = VecDestroy(&baij->slvec1b);CHKERRQ(ierr); 30101b2bd88SHong Zhang 302a30f8f8cSSatish Balay if (baij->colmap) { 303a30f8f8cSSatish Balay #if defined(PETSC_USE_CTABLE) 3046bc0bbbfSBarry Smith ierr = PetscTableDestroy(&baij->colmap);CHKERRQ(ierr); 305a30f8f8cSSatish Balay #else 306a30f8f8cSSatish Balay ierr = PetscFree(baij->colmap);CHKERRQ(ierr); 3073bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)A,-Bbaij->nbs*sizeof(PetscInt));CHKERRQ(ierr); 308a30f8f8cSSatish Balay #endif 309a30f8f8cSSatish Balay } 310a30f8f8cSSatish Balay 311a30f8f8cSSatish Balay /* make sure that B is assembled so we can access its values */ 312a30f8f8cSSatish Balay ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 313a30f8f8cSSatish Balay ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 314a30f8f8cSSatish Balay 315a30f8f8cSSatish Balay /* invent new B and copy stuff over */ 316785e854fSJed Brown ierr = PetscMalloc1(mbs,&nz);CHKERRQ(ierr); 317a30f8f8cSSatish Balay for (i=0; i<mbs; i++) { 318a30f8f8cSSatish Balay nz[i] = Bbaij->i[i+1]-Bbaij->i[i]; 319a30f8f8cSSatish Balay } 320f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,&Bnew);CHKERRQ(ierr); 321f69a0ea3SMatthew Knepley ierr = MatSetSizes(Bnew,m,n,m,n);CHKERRQ(ierr); 3227adad957SLisandro Dalcin ierr = MatSetType(Bnew,((PetscObject)B)->type_name);CHKERRQ(ierr); 323d0f46423SBarry Smith ierr = MatSeqBAIJSetPreallocation(Bnew,B->rmap->bs,0,nz);CHKERRQ(ierr); 324a30f8f8cSSatish Balay ierr = PetscFree(nz);CHKERRQ(ierr); 325a30f8f8cSSatish Balay 326b38c15b3SStefano Zampini if (Bbaij->nonew >= 0) { /* Inherit insertion error options (if positive). */ 327b38c15b3SStefano Zampini ((Mat_SeqSBAIJ*)Bnew->data)->nonew = Bbaij->nonew; 328b38c15b3SStefano Zampini } 329b38c15b3SStefano Zampini 33077341eacSDmitry Karpeev /* 33177341eacSDmitry Karpeev Ensure that B's nonzerostate is monotonically increasing. 33277341eacSDmitry Karpeev Or should this follow the MatSetValues() loop to preserve B's nonzerstate across a MatDisAssemble() call? 33377341eacSDmitry Karpeev */ 33477341eacSDmitry Karpeev Bnew->nonzerostate = B->nonzerostate; 335785e854fSJed Brown ierr = PetscMalloc1(bs,&rvals);CHKERRQ(ierr); 336a30f8f8cSSatish Balay for (i=0; i<mbs; i++) { 337a30f8f8cSSatish Balay rvals[0] = bs*i; 33826fbe8dcSKarl Rupp for (j=1; j<bs; j++) rvals[j] = rvals[j-1] + 1; 339a30f8f8cSSatish Balay for (j=Bbaij->i[i]; j<Bbaij->i[i+1]; j++) { 340a30f8f8cSSatish Balay col = garray[Bbaij->j[j]]*bs; 341a30f8f8cSSatish Balay for (k=0; k<bs; k++) { 342ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_MAT_SINGLE) 343a30f8f8cSSatish Balay for (l=0; l<bs; l++) atmp[l] = a[j*bs2+l]; 344a30f8f8cSSatish Balay #else 34571730473SSatish Balay atmp = a+j*bs2 + k*bs; 346a30f8f8cSSatish Balay #endif 347c8407628SSatish Balay ierr = MatSetValues_SeqSBAIJ(Bnew,bs,rvals,1,&col,atmp,B->insertmode);CHKERRQ(ierr); 348a30f8f8cSSatish Balay col++; 349a30f8f8cSSatish Balay } 350a30f8f8cSSatish Balay } 351a30f8f8cSSatish Balay } 352ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_MAT_SINGLE) 353a30f8f8cSSatish Balay ierr = PetscFree(atmp);CHKERRQ(ierr); 354a30f8f8cSSatish Balay #endif 355a30f8f8cSSatish Balay ierr = PetscFree(baij->garray);CHKERRQ(ierr); 35626fbe8dcSKarl Rupp 357a30f8f8cSSatish Balay baij->garray = 0; 35826fbe8dcSKarl Rupp 359a30f8f8cSSatish Balay ierr = PetscFree(rvals);CHKERRQ(ierr); 3603bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)A,-ec*sizeof(PetscInt));CHKERRQ(ierr); 3616bf464f9SBarry Smith ierr = MatDestroy(&B);CHKERRQ(ierr); 3623bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)A,(PetscObject)Bnew);CHKERRQ(ierr); 36326fbe8dcSKarl Rupp 364a30f8f8cSSatish Balay baij->B = Bnew; 36526fbe8dcSKarl Rupp 366a30f8f8cSSatish Balay A->was_assembled = PETSC_FALSE; 367a30f8f8cSSatish Balay PetscFunctionReturn(0); 368a30f8f8cSSatish Balay } 369