xref: /petsc/src/mat/impls/baij/mpi/mmbaij.c (revision 6a719282e8b9d12a5dfd16d140b52bac3a4cf0eb)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2be1d678aSKris Buschelman 
38016bdd1SSatish Balay /*
4d9653453SSatish Balay    Support for the parallel BAIJ matrix vector multiply
58016bdd1SSatish Balay */
67c4f633dSBarry Smith #include "../src/mat/impls/baij/mpi/mpibaij.h"
7bba1ac68SSatish Balay 
8b24ad042SBarry Smith EXTERN PetscErrorCode MatSetValuesBlocked_SeqBAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
98016bdd1SSatish Balay 
104a2ae208SSatish Balay #undef __FUNCT__
114a2ae208SSatish Balay #define __FUNCT__ "MatSetUpMultiply_MPIBAIJ"
12dfbe8321SBarry Smith PetscErrorCode MatSetUpMultiply_MPIBAIJ(Mat mat)
138016bdd1SSatish Balay {
14d9653453SSatish Balay   Mat_MPIBAIJ    *baij = (Mat_MPIBAIJ*)mat->data;
15d9653453SSatish Balay   Mat_SeqBAIJ    *B = (Mat_SeqBAIJ*)(baij->B->data);
166849ba73SBarry Smith   PetscErrorCode ierr;
17b24ad042SBarry Smith   PetscInt       i,j,*aj = B->j,ec = 0,*garray;
18d0f46423SBarry Smith   PetscInt       bs = mat->rmap->bs,*stmp;
198016bdd1SSatish Balay   IS             from,to;
208016bdd1SSatish Balay   Vec            gvec;
21aa482453SBarry Smith #if defined (PETSC_USE_CTABLE)
220f5bd95cSBarry Smith   PetscTable     gid1_lid1;
230f5bd95cSBarry Smith   PetscTablePosition tpos;
24b24ad042SBarry Smith   PetscInt       gid,lid;
256f531f54SSatish Balay #else
26b24ad042SBarry Smith   PetscInt       Nbs = baij->Nbs,*indices;
2773a2e727SSatish Balay #endif
288016bdd1SSatish Balay 
293a40ed3dSBarry Smith   PetscFunctionBegin;
3073a2e727SSatish Balay 
31aa482453SBarry Smith #if defined (PETSC_USE_CTABLE)
3273a2e727SSatish Balay   /* use a table - Mark Adams */
33273d9f13SBarry Smith   ierr = PetscTableCreate(B->mbs,&gid1_lid1);CHKERRQ(ierr);
3473a2e727SSatish Balay   for (i=0; i<B->mbs; i++) {
3573a2e727SSatish Balay     for (j=0; j<B->ilen[i]; j++) {
36b24ad042SBarry Smith       PetscInt data,gid1 = aj[B->i[i]+j] + 1;
370f5bd95cSBarry Smith       ierr = PetscTableFind(gid1_lid1,gid1,&data);CHKERRQ(ierr);
38fa46199cSSatish Balay       if (!data) {
3973a2e727SSatish Balay         /* one based table */
400f5bd95cSBarry Smith         ierr = PetscTableAdd(gid1_lid1,gid1,++ec);CHKERRQ(ierr);
4173a2e727SSatish Balay       }
4273a2e727SSatish Balay     }
4373a2e727SSatish Balay   }
4473a2e727SSatish Balay   /* form array of columns we need */
45b24ad042SBarry Smith   ierr = PetscMalloc((ec+1)*sizeof(PetscInt),&garray);CHKERRQ(ierr);
460f5bd95cSBarry Smith   ierr = PetscTableGetHeadPosition(gid1_lid1,&tpos);CHKERRQ(ierr);
4773a2e727SSatish Balay   while (tpos) {
480f5bd95cSBarry Smith     ierr = PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);CHKERRQ(ierr);
4973a2e727SSatish Balay     gid--; lid--;
5073a2e727SSatish Balay     garray[lid] = gid;
5173a2e727SSatish Balay   }
520064e2bbSSatish Balay   ierr = PetscSortInt(ec,garray);CHKERRQ(ierr);
530f5bd95cSBarry Smith   ierr = PetscTableRemoveAll(gid1_lid1);CHKERRQ(ierr);
5473a2e727SSatish Balay   for (i=0; i<ec; i++) {
550f5bd95cSBarry Smith     ierr = PetscTableAdd(gid1_lid1,garray[i]+1,i+1);CHKERRQ(ierr);
5673a2e727SSatish Balay   }
5773a2e727SSatish Balay   /* compact out the extra columns in B */
5873a2e727SSatish Balay   for (i=0; i<B->mbs; i++) {
5973a2e727SSatish Balay     for (j=0; j<B->ilen[i]; j++) {
60b24ad042SBarry Smith       PetscInt gid1 = aj[B->i[i] + j] + 1;
610f5bd95cSBarry Smith       ierr = PetscTableFind(gid1_lid1,gid1,&lid);CHKERRQ(ierr);
62fa46199cSSatish Balay       lid --;
6373a2e727SSatish Balay       aj[B->i[i]+j] = lid;
6473a2e727SSatish Balay     }
6573a2e727SSatish Balay   }
6673a2e727SSatish Balay   B->nbs     = ec;
67d0f46423SBarry Smith   baij->B->cmap->n = baij->B->cmap->N = ec*mat->rmap->bs;
6826283091SBarry Smith   ierr = PetscLayoutSetUp((baij->B->cmap));CHKERRQ(ierr);
699c666560SBarry Smith   ierr = PetscTableDestroy(gid1_lid1);CHKERRQ(ierr);
7073a2e727SSatish Balay   /* Mark Adams */
7173a2e727SSatish Balay #else
72435da068SBarry Smith   /* Make an array as long as the number of columns */
73d9653453SSatish Balay   /* mark those columns that are in baij->B */
74b24ad042SBarry Smith   ierr = PetscMalloc((Nbs+1)*sizeof(PetscInt),&indices);CHKERRQ(ierr);
75b24ad042SBarry Smith   ierr = PetscMemzero(indices,Nbs*sizeof(PetscInt));CHKERRQ(ierr);
76d9653453SSatish Balay   for (i=0; i<B->mbs; i++) {
778016bdd1SSatish Balay     for (j=0; j<B->ilen[i]; j++) {
78d9653453SSatish Balay       if (!indices[aj[B->i[i] + j]]) ec++;
79d9653453SSatish Balay       indices[aj[B->i[i] + j]] = 1;
808016bdd1SSatish Balay     }
818016bdd1SSatish Balay   }
828016bdd1SSatish Balay 
838016bdd1SSatish Balay   /* form array of columns we need */
84b24ad042SBarry Smith   ierr = PetscMalloc((ec+1)*sizeof(PetscInt),&garray);CHKERRQ(ierr);
858016bdd1SSatish Balay   ec = 0;
86d9653453SSatish Balay   for (i=0; i<Nbs; i++) {
870bdbc534SSatish Balay     if (indices[i]) {
880bdbc534SSatish Balay       garray[ec++] = i;
890bdbc534SSatish Balay     }
908016bdd1SSatish Balay   }
918016bdd1SSatish Balay 
928016bdd1SSatish Balay   /* make indices now point into garray */
938016bdd1SSatish Balay   for (i=0; i<ec; i++) {
94d9653453SSatish Balay     indices[garray[i]] = i;
958016bdd1SSatish Balay   }
968016bdd1SSatish Balay 
978016bdd1SSatish Balay   /* compact out the extra columns in B */
98d9653453SSatish Balay   for (i=0; i<B->mbs; i++) {
998016bdd1SSatish Balay     for (j=0; j<B->ilen[i]; j++) {
100d9653453SSatish Balay       aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
1018016bdd1SSatish Balay     }
1028016bdd1SSatish Balay   }
103d9653453SSatish Balay   B->nbs       = ec;
104d0f46423SBarry Smith   baij->B->cmap->n =baij->B->cmap->N  = ec*mat->rmap->bs;
10526283091SBarry Smith   ierr = PetscLayoutSetUp((baij->B->cmap));CHKERRQ(ierr);
106606d414cSSatish Balay   ierr = PetscFree(indices);CHKERRQ(ierr);
10773a2e727SSatish Balay #endif
1088016bdd1SSatish Balay 
1098016bdd1SSatish Balay   /* create local vector that is used to scatter into */
110029af93fSBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,ec*bs,&baij->lvec);CHKERRQ(ierr);
1118016bdd1SSatish Balay 
112c16cb8f2SBarry Smith   /* create two temporary index sets for building scatter-gather */
113deff0451SBarry Smith   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
114c16cb8f2SBarry Smith 
115b24ad042SBarry Smith   ierr = PetscMalloc((ec+1)*sizeof(PetscInt),&stmp);CHKERRQ(ierr);
116e82e9f6bSBarry Smith   for (i=0; i<ec; i++) { stmp[i] = i; }
117deff0451SBarry Smith   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
1188016bdd1SSatish Balay 
1198016bdd1SSatish Balay   /* create temporary global vector to generate scatter context */
120d0f46423SBarry Smith   ierr = VecCreateMPIWithArray(((PetscObject)mat)->comm,mat->cmap->n,mat->cmap->N,PETSC_NULL,&gvec);CHKERRQ(ierr);
1218016bdd1SSatish Balay 
122d9653453SSatish Balay   ierr = VecScatterCreate(gvec,from,baij->lvec,to,&baij->Mvctx);CHKERRQ(ierr);
12390f02eecSBarry Smith 
12452e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,baij->Mvctx);CHKERRQ(ierr);
12552e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,baij->lvec);CHKERRQ(ierr);
12652e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,from);CHKERRQ(ierr);
12752e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,to);CHKERRQ(ierr);
128d9653453SSatish Balay   baij->garray = garray;
12952e6d16bSBarry Smith   ierr = PetscLogObjectMemory(mat,(ec+1)*sizeof(PetscInt));CHKERRQ(ierr);
1308016bdd1SSatish Balay   ierr = ISDestroy(from);CHKERRQ(ierr);
1318016bdd1SSatish Balay   ierr = ISDestroy(to);CHKERRQ(ierr);
132888f2ed8SSatish Balay   ierr = VecDestroy(gvec);CHKERRQ(ierr);
1333a40ed3dSBarry Smith   PetscFunctionReturn(0);
1348016bdd1SSatish Balay }
1358016bdd1SSatish Balay 
1368016bdd1SSatish Balay /*
137d9653453SSatish Balay      Takes the local part of an already assembled MPIBAIJ matrix
1388016bdd1SSatish Balay    and disassembles it. This is to allow new nonzeros into the matrix
1398016bdd1SSatish Balay    that require more communication in the matrix vector multiply.
1408016bdd1SSatish Balay    Thus certain data-structures must be rebuilt.
1418016bdd1SSatish Balay 
1428016bdd1SSatish Balay    Kind of slow! But that's what application programmers get when
1438016bdd1SSatish Balay    they are sloppy.
1448016bdd1SSatish Balay */
1454a2ae208SSatish Balay #undef __FUNCT__
1464a2ae208SSatish Balay #define __FUNCT__ "DisAssemble_MPIBAIJ"
147dfbe8321SBarry Smith PetscErrorCode DisAssemble_MPIBAIJ(Mat A)
1488016bdd1SSatish Balay {
149d9653453SSatish Balay   Mat_MPIBAIJ    *baij = (Mat_MPIBAIJ*)A->data;
150d9653453SSatish Balay   Mat            B = baij->B,Bnew;
151d9653453SSatish Balay   Mat_SeqBAIJ    *Bbaij = (Mat_SeqBAIJ*)B->data;
152dfbe8321SBarry Smith   PetscErrorCode ierr;
153d0f46423SBarry Smith   PetscInt       i,j,mbs=Bbaij->mbs,n = A->cmap->N,col,*garray=baij->garray;
154d0f46423SBarry Smith   PetscInt       bs2 = baij->bs2,*nz,ec,m = A->rmap->n;
1553eda8832SBarry Smith   MatScalar      *a = Bbaij->a;
156dd6ea824SBarry Smith   MatScalar      *atmp;
15797e5c40aSBarry Smith 
1588016bdd1SSatish Balay 
1593a40ed3dSBarry Smith   PetscFunctionBegin;
1608016bdd1SSatish Balay   /* free stuff related to matrix-vec multiply */
161b0a32e0cSBarry Smith   ierr = VecGetSize(baij->lvec,&ec);CHKERRQ(ierr); /* needed for PetscLogObjectMemory below */
162d9653453SSatish Balay   ierr = VecDestroy(baij->lvec);CHKERRQ(ierr); baij->lvec = 0;
163d9653453SSatish Balay   ierr = VecScatterDestroy(baij->Mvctx);CHKERRQ(ierr); baij->Mvctx = 0;
164d9653453SSatish Balay   if (baij->colmap) {
165aa482453SBarry Smith #if defined (PETSC_USE_CTABLE)
1669c666560SBarry Smith     ierr = PetscTableDestroy(baij->colmap); baij->colmap = 0;CHKERRQ(ierr);
16748e59246SSatish Balay #else
168606d414cSSatish Balay     ierr = PetscFree(baij->colmap);CHKERRQ(ierr);
169606d414cSSatish Balay     baij->colmap = 0;
17052e6d16bSBarry Smith     ierr = PetscLogObjectMemory(A,-Bbaij->nbs*sizeof(PetscInt));CHKERRQ(ierr);
17148e59246SSatish Balay #endif
1728016bdd1SSatish Balay   }
1738016bdd1SSatish Balay 
1748016bdd1SSatish Balay   /* make sure that B is assembled so we can access its values */
1756d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1763eda8832SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1778016bdd1SSatish Balay 
1788016bdd1SSatish Balay   /* invent new B and copy stuff over */
179b24ad042SBarry Smith   ierr = PetscMalloc(mbs*sizeof(PetscInt),&nz);CHKERRQ(ierr);
180d9653453SSatish Balay   for (i=0; i<mbs; i++) {
181d9653453SSatish Balay     nz[i] = Bbaij->i[i+1]-Bbaij->i[i];
1828016bdd1SSatish Balay   }
1837adad957SLisandro Dalcin   ierr = MatCreate(((PetscObject)B)->comm,&Bnew);CHKERRQ(ierr);
184f69a0ea3SMatthew Knepley   ierr = MatSetSizes(Bnew,m,n,m,n);CHKERRQ(ierr);
1857adad957SLisandro Dalcin   ierr = MatSetType(Bnew,((PetscObject)B)->type_name);CHKERRQ(ierr);
186d0f46423SBarry Smith   ierr = MatSeqBAIJSetPreallocation(Bnew,B->rmap->bs,0,nz);CHKERRQ(ierr);
1874e0d8c25SBarry Smith   ierr = MatSetOption(Bnew,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr);
188d9653453SSatish Balay 
189bba1ac68SSatish Balay   for (i=0; i<mbs; i++) {
190bba1ac68SSatish Balay     for (j=Bbaij->i[i]; j<Bbaij->i[i+1]; j++) {
191bba1ac68SSatish Balay       col  = garray[Bbaij->j[j]];
1923eda8832SBarry Smith       atmp = a + j*bs2;
193bba1ac68SSatish Balay       ierr = MatSetValuesBlocked_SeqBAIJ(Bnew,1,&i,1,&col,atmp,B->insertmode);CHKERRQ(ierr);
1948016bdd1SSatish Balay     }
1958016bdd1SSatish Balay   }
1964e0d8c25SBarry Smith   ierr = MatSetOption(Bnew,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr);
197bba1ac68SSatish Balay 
198bba1ac68SSatish Balay   ierr = PetscFree(nz);CHKERRQ(ierr);
199606d414cSSatish Balay   ierr = PetscFree(baij->garray);CHKERRQ(ierr);
200606d414cSSatish Balay   baij->garray = 0;
20152e6d16bSBarry Smith   ierr = PetscLogObjectMemory(A,-ec*sizeof(PetscInt));CHKERRQ(ierr);
2028016bdd1SSatish Balay   ierr = MatDestroy(B);CHKERRQ(ierr);
20352e6d16bSBarry Smith   ierr = PetscLogObjectParent(A,Bnew);CHKERRQ(ierr);
204d9653453SSatish Balay   baij->B = Bnew;
2058016bdd1SSatish Balay   A->was_assembled = PETSC_FALSE;
206*6a719282SBarry Smith   A->assembled     = PETSC_FALSE;
2073a40ed3dSBarry Smith   PetscFunctionReturn(0);
2088016bdd1SSatish Balay }
2098016bdd1SSatish Balay 
21004f1ad80SBarry Smith /*      ugly stuff added for Glenn someday we should fix this up */
21104f1ad80SBarry Smith 
212b24ad042SBarry Smith static PetscInt *uglyrmapd = 0,*uglyrmapo = 0;  /* mapping from the local ordering to the "diagonal" and "off-diagonal"
21304f1ad80SBarry Smith                                       parts of the local matrix */
2142cd6534aSBarry Smith static Vec uglydd = 0,uglyoo = 0;   /* work vectors used to scale the two parts of the local matrix */
21504f1ad80SBarry Smith 
21604f1ad80SBarry Smith 
217348d1a9dSSatish Balay #undef __FUNCT__
218348d1a9dSSatish Balay #define __FUNCT__ "MatMPIBAIJDiagonalScaleLocalSetUp"
219dfbe8321SBarry Smith PetscErrorCode MatMPIBAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale)
22004f1ad80SBarry Smith {
22104f1ad80SBarry Smith   Mat_MPIBAIJ    *ina = (Mat_MPIBAIJ*) inA->data; /*access private part of matrix */
22204f1ad80SBarry Smith   Mat_SeqBAIJ    *B = (Mat_SeqBAIJ*)ina->B->data;
223dfbe8321SBarry Smith   PetscErrorCode ierr;
224d0f46423SBarry Smith   PetscInt       bs = inA->rmap->bs,i,n,nt,j,cstart,cend,no,*garray = ina->garray,*lindices;
225b24ad042SBarry Smith   PetscInt       *r_rmapd,*r_rmapo;
22604f1ad80SBarry Smith 
22704f1ad80SBarry Smith   PetscFunctionBegin;
22804f1ad80SBarry Smith   ierr = MatGetOwnershipRange(inA,&cstart,&cend);CHKERRQ(ierr);
22904f1ad80SBarry Smith   ierr = MatGetSize(ina->A,PETSC_NULL,&n);CHKERRQ(ierr);
230b24ad042SBarry Smith   ierr = PetscMalloc((inA->bmapping->n+1)*sizeof(PetscInt),&r_rmapd);CHKERRQ(ierr);
231b24ad042SBarry Smith   ierr = PetscMemzero(r_rmapd,inA->bmapping->n*sizeof(PetscInt));CHKERRQ(ierr);
23204f1ad80SBarry Smith   nt   = 0;
23304f1ad80SBarry Smith   for (i=0; i<inA->bmapping->n; i++) {
23404f1ad80SBarry Smith     if (inA->bmapping->indices[i]*bs >= cstart && inA->bmapping->indices[i]*bs < cend) {
23504f1ad80SBarry Smith       nt++;
23604f1ad80SBarry Smith       r_rmapd[i] = inA->bmapping->indices[i] + 1;
23704f1ad80SBarry Smith     }
23804f1ad80SBarry Smith   }
239e32f2f54SBarry Smith   if (nt*bs != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt*bs %D n %D",nt*bs,n);
240b24ad042SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&uglyrmapd);CHKERRQ(ierr);
24104f1ad80SBarry Smith   for (i=0; i<inA->bmapping->n; i++) {
24204f1ad80SBarry Smith     if (r_rmapd[i]){
24304f1ad80SBarry Smith       for (j=0; j<bs; j++) {
24404f1ad80SBarry Smith         uglyrmapd[(r_rmapd[i]-1)*bs+j-cstart] = i*bs + j;
24504f1ad80SBarry Smith       }
24604f1ad80SBarry Smith     }
24704f1ad80SBarry Smith   }
24804f1ad80SBarry Smith   ierr = PetscFree(r_rmapd);CHKERRQ(ierr);
24904f1ad80SBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,n,&uglydd);CHKERRQ(ierr);
25004f1ad80SBarry Smith 
251b24ad042SBarry Smith   ierr = PetscMalloc((ina->Nbs+1)*sizeof(PetscInt),&lindices);CHKERRQ(ierr);
252b24ad042SBarry Smith   ierr = PetscMemzero(lindices,ina->Nbs*sizeof(PetscInt));CHKERRQ(ierr);
25304f1ad80SBarry Smith   for (i=0; i<B->nbs; i++) {
25404f1ad80SBarry Smith     lindices[garray[i]] = i+1;
25504f1ad80SBarry Smith   }
25604f1ad80SBarry Smith   no   = inA->bmapping->n - nt;
257b24ad042SBarry Smith   ierr = PetscMalloc((inA->bmapping->n+1)*sizeof(PetscInt),&r_rmapo);CHKERRQ(ierr);
258b24ad042SBarry Smith   ierr = PetscMemzero(r_rmapo,inA->bmapping->n*sizeof(PetscInt));CHKERRQ(ierr);
25904f1ad80SBarry Smith   nt   = 0;
26004f1ad80SBarry Smith   for (i=0; i<inA->bmapping->n; i++) {
26104f1ad80SBarry Smith     if (lindices[inA->bmapping->indices[i]]) {
26204f1ad80SBarry Smith       nt++;
26304f1ad80SBarry Smith       r_rmapo[i] = lindices[inA->bmapping->indices[i]];
26404f1ad80SBarry Smith     }
26504f1ad80SBarry Smith   }
266e32f2f54SBarry Smith   if (nt > no) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D no %D",nt,n);
26704f1ad80SBarry Smith   ierr = PetscFree(lindices);CHKERRQ(ierr);
268b24ad042SBarry Smith   ierr = PetscMalloc((nt*bs+1)*sizeof(PetscInt),&uglyrmapo);CHKERRQ(ierr);
26904f1ad80SBarry Smith   for (i=0; i<inA->bmapping->n; i++) {
27004f1ad80SBarry Smith     if (r_rmapo[i]){
27104f1ad80SBarry Smith       for (j=0; j<bs; j++) {
27204f1ad80SBarry Smith         uglyrmapo[(r_rmapo[i]-1)*bs+j] = i*bs + j;
27304f1ad80SBarry Smith       }
27404f1ad80SBarry Smith     }
27504f1ad80SBarry Smith   }
27604f1ad80SBarry Smith   ierr = PetscFree(r_rmapo);CHKERRQ(ierr);
27704f1ad80SBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,nt*bs,&uglyoo);CHKERRQ(ierr);
27804f1ad80SBarry Smith 
27904f1ad80SBarry Smith   PetscFunctionReturn(0);
28004f1ad80SBarry Smith }
28104f1ad80SBarry Smith 
282348d1a9dSSatish Balay #undef __FUNCT__
283348d1a9dSSatish Balay #define __FUNCT__ "MatMPIBAIJDiagonalScaleLocal"
284be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJDiagonalScaleLocal(Mat A,Vec scale)
28504f1ad80SBarry Smith {
28692b32695SKris Buschelman   /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */
2874ac538c5SBarry Smith   PetscErrorCode ierr;
28892b32695SKris Buschelman 
28992b32695SKris Buschelman   PetscFunctionBegin;
2904ac538c5SBarry Smith   ierr = PetscTryMethod(A,"MatDiagonalScaleLocal_C",(Mat,Vec),(A,scale));CHKERRQ(ierr);
29192b32695SKris Buschelman   PetscFunctionReturn(0);
29292b32695SKris Buschelman }
29392b32695SKris Buschelman 
29492b32695SKris Buschelman EXTERN_C_BEGIN
29592b32695SKris Buschelman #undef __FUNCT__
29692b32695SKris Buschelman #define __FUNCT__ "MatDiagonalScaleLocal_MPIBAIJ"
297be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScaleLocal_MPIBAIJ(Mat A,Vec scale)
29892b32695SKris Buschelman {
29904f1ad80SBarry Smith   Mat_MPIBAIJ    *a = (Mat_MPIBAIJ*) A->data; /*access private part of matrix */
300dfbe8321SBarry Smith   PetscErrorCode ierr;
301b24ad042SBarry Smith   PetscInt       n,i;
30204f1ad80SBarry Smith   PetscScalar    *d,*o,*s;
30304f1ad80SBarry Smith 
30404f1ad80SBarry Smith   PetscFunctionBegin;
3052cd6534aSBarry Smith   if (!uglyrmapd) {
3062cd6534aSBarry Smith     ierr = MatMPIBAIJDiagonalScaleLocalSetUp(A,scale);CHKERRQ(ierr);
3072cd6534aSBarry Smith   }
3082cd6534aSBarry Smith 
3091ebc52fbSHong Zhang   ierr = VecGetArray(scale,&s);CHKERRQ(ierr);
31004f1ad80SBarry Smith 
31104f1ad80SBarry Smith   ierr = VecGetLocalSize(uglydd,&n);CHKERRQ(ierr);
3121ebc52fbSHong Zhang   ierr = VecGetArray(uglydd,&d);CHKERRQ(ierr);
31304f1ad80SBarry Smith   for (i=0; i<n; i++) {
31404f1ad80SBarry Smith     d[i] = s[uglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */
31504f1ad80SBarry Smith   }
3161ebc52fbSHong Zhang   ierr = VecRestoreArray(uglydd,&d);CHKERRQ(ierr);
31704f1ad80SBarry Smith   /* column scale "diagonal" portion of local matrix */
31804f1ad80SBarry Smith   ierr = MatDiagonalScale(a->A,PETSC_NULL,uglydd);CHKERRQ(ierr);
31904f1ad80SBarry Smith 
32004f1ad80SBarry Smith   ierr = VecGetLocalSize(uglyoo,&n);CHKERRQ(ierr);
3211ebc52fbSHong Zhang   ierr = VecGetArray(uglyoo,&o);CHKERRQ(ierr);
32204f1ad80SBarry Smith   for (i=0; i<n; i++) {
32304f1ad80SBarry Smith     o[i] = s[uglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */
32404f1ad80SBarry Smith   }
3251ebc52fbSHong Zhang   ierr = VecRestoreArray(scale,&s);CHKERRQ(ierr);
3261ebc52fbSHong Zhang   ierr = VecRestoreArray(uglyoo,&o);CHKERRQ(ierr);
32704f1ad80SBarry Smith   /* column scale "off-diagonal" portion of local matrix */
32804f1ad80SBarry Smith   ierr = MatDiagonalScale(a->B,PETSC_NULL,uglyoo);CHKERRQ(ierr);
32904f1ad80SBarry Smith 
33004f1ad80SBarry Smith   PetscFunctionReturn(0);
33104f1ad80SBarry Smith }
33292b32695SKris Buschelman EXTERN_C_END
3338016bdd1SSatish Balay 
3343eda8832SBarry Smith 
335