xref: /petsc/src/mat/impls/aij/mpi/mmaij.c (revision a2f3521de2a9d7869700f4c17e26c23fcfeaa6f6)
1be1d678aSKris Buschelman 
28c79f6d3SBarry Smith /*
38c79f6d3SBarry Smith    Support for the parallel AIJ matrix vector multiply
48c79f6d3SBarry Smith */
5c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h>
68c79f6d3SBarry Smith 
74a2ae208SSatish Balay #undef __FUNCT__
84a2ae208SSatish Balay #define __FUNCT__ "MatSetUpMultiply_MPIAIJ"
9dfbe8321SBarry Smith PetscErrorCode MatSetUpMultiply_MPIAIJ(Mat mat)
108c79f6d3SBarry Smith {
1144a69424SLois Curfman McInnes   Mat_MPIAIJ         *aij = (Mat_MPIAIJ*)mat->data;
12ec8511deSBarry Smith   Mat_SeqAIJ         *B = (Mat_SeqAIJ*)(aij->B->data);
136849ba73SBarry Smith   PetscErrorCode     ierr;
14b1d57f15SBarry Smith   PetscInt           i,j,*aj = B->j,ec = 0,*garray;
151eb62cbbSBarry Smith   IS                 from,to;
161eb62cbbSBarry Smith   Vec                gvec;
17ace3abfcSBarry Smith   PetscBool          useblockis;
18aa482453SBarry Smith #if defined (PETSC_USE_CTABLE)
190f5bd95cSBarry Smith   PetscTable         gid1_lid1;
200f5bd95cSBarry Smith   PetscTablePosition tpos;
21b1d57f15SBarry Smith   PetscInt           gid,lid;
226f531f54SSatish Balay #else
23d0f46423SBarry Smith   PetscInt           N = mat->cmap->N,*indices;
242066d6f7SSatish Balay #endif
252066d6f7SSatish Balay 
263a40ed3dSBarry Smith   PetscFunctionBegin;
272066d6f7SSatish Balay 
28aa482453SBarry Smith #if defined (PETSC_USE_CTABLE)
29dc231df0SBarry Smith   /* use a table - Mark Adams */
30e23dfa41SBarry Smith   ierr = PetscTableCreate(aij->B->rmap->n,mat->cmap->N+1,&gid1_lid1);CHKERRQ(ierr);
31d0f46423SBarry Smith   for (i=0; i<aij->B->rmap->n; i++) {
322066d6f7SSatish Balay     for (j=0; j<B->ilen[i]; j++) {
33b1d57f15SBarry Smith       PetscInt data,gid1 = aj[B->i[i] + j] + 1;
340f5bd95cSBarry Smith       ierr = PetscTableFind(gid1_lid1,gid1,&data);CHKERRQ(ierr);
35fa46199cSSatish Balay       if (!data) {
362066d6f7SSatish Balay         /* one based table */
373861aac3SJed Brown         ierr = PetscTableAdd(gid1_lid1,gid1,++ec,INSERT_VALUES);CHKERRQ(ierr);
382066d6f7SSatish Balay       }
392066d6f7SSatish Balay     }
402066d6f7SSatish Balay   }
412066d6f7SSatish Balay   /* form array of columns we need */
42b1d57f15SBarry Smith   ierr = PetscMalloc((ec+1)*sizeof(PetscInt),&garray);CHKERRQ(ierr);
430f5bd95cSBarry Smith   ierr = PetscTableGetHeadPosition(gid1_lid1,&tpos);CHKERRQ(ierr);
442066d6f7SSatish Balay   while (tpos) {
450f5bd95cSBarry Smith     ierr = PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);CHKERRQ(ierr);
46b0a32e0cSBarry Smith     gid--;
47b0a32e0cSBarry Smith     lid--;
482066d6f7SSatish Balay     garray[lid] = gid;
492066d6f7SSatish Balay   }
500064e2bbSSatish Balay   ierr = PetscSortInt(ec,garray);CHKERRQ(ierr); /* sort, and rebuild */
510f5bd95cSBarry Smith   ierr = PetscTableRemoveAll(gid1_lid1);CHKERRQ(ierr);
522066d6f7SSatish Balay   for (i=0; i<ec; i++) {
533861aac3SJed Brown     ierr = PetscTableAdd(gid1_lid1,garray[i]+1,i+1,INSERT_VALUES);CHKERRQ(ierr);
542066d6f7SSatish Balay   }
552066d6f7SSatish Balay   /* compact out the extra columns in B */
56d0f46423SBarry Smith   for (i=0; i<aij->B->rmap->n; i++) {
572066d6f7SSatish Balay     for (j=0; j<B->ilen[i]; j++) {
58b1d57f15SBarry Smith       PetscInt gid1 = aj[B->i[i] + j] + 1;
590f5bd95cSBarry Smith       ierr = PetscTableFind(gid1_lid1,gid1,&lid);CHKERRQ(ierr);
60fa46199cSSatish Balay       lid --;
61b3fb0a6cSMatthew Knepley       aj[B->i[i] + j]  = lid;
622066d6f7SSatish Balay     }
632066d6f7SSatish Balay   }
64d0f46423SBarry Smith   aij->B->cmap->n = aij->B->cmap->N = ec;
6526283091SBarry Smith   ierr = PetscLayoutSetUp((aij->B->cmap));CHKERRQ(ierr);
666bc0bbbfSBarry Smith   ierr = PetscTableDestroy(&gid1_lid1);CHKERRQ(ierr);
672066d6f7SSatish Balay   /* Mark Adams */
682066d6f7SSatish Balay #else
6911285404SBarry Smith   /* Make an array as long as the number of columns */
701eb62cbbSBarry Smith   /* mark those columns that are in aij->B */
71b1d57f15SBarry Smith   ierr = PetscMalloc((N+1)*sizeof(PetscInt),&indices);CHKERRQ(ierr);
72b1d57f15SBarry Smith   ierr = PetscMemzero(indices,N*sizeof(PetscInt));CHKERRQ(ierr);
73d0f46423SBarry Smith   for (i=0; i<aij->B->rmap->n; i++) {
74d6dfbf8fSBarry Smith     for (j=0; j<B->ilen[i]; j++) {
75bfec09a0SHong Zhang       if (!indices[aj[B->i[i] + j] ]) ec++;
76bfec09a0SHong Zhang       indices[aj[B->i[i] + j] ] = 1;
77416022c9SBarry Smith     }
781eb62cbbSBarry Smith   }
798c79f6d3SBarry Smith 
801eb62cbbSBarry Smith   /* form array of columns we need */
81b1d57f15SBarry Smith   ierr = PetscMalloc((ec+1)*sizeof(PetscInt),&garray);CHKERRQ(ierr);
821eb62cbbSBarry Smith   ec = 0;
831eb62cbbSBarry Smith   for (i=0; i<N; i++) {
841eb62cbbSBarry Smith     if (indices[i]) garray[ec++] = i;
851eb62cbbSBarry Smith   }
861eb62cbbSBarry Smith 
871eb62cbbSBarry Smith   /* make indices now point into garray */
881eb62cbbSBarry Smith   for (i=0; i<ec; i++) {
89bfec09a0SHong Zhang     indices[garray[i]] = i;
901eb62cbbSBarry Smith   }
911eb62cbbSBarry Smith 
921eb62cbbSBarry Smith   /* compact out the extra columns in B */
93d0f46423SBarry Smith   for (i=0; i<aij->B->rmap->n; i++) {
94d6dfbf8fSBarry Smith     for (j=0; j<B->ilen[i]; j++) {
95bfec09a0SHong Zhang       aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
961eb62cbbSBarry Smith     }
97d6dfbf8fSBarry Smith   }
98d0f46423SBarry Smith   aij->B->cmap->n = aij->B->cmap->N = ec;
9926283091SBarry Smith   ierr = PetscLayoutSetUp((aij->B->cmap));CHKERRQ(ierr);
100606d414cSSatish Balay   ierr = PetscFree(indices);CHKERRQ(ierr);
1012066d6f7SSatish Balay #endif
1021eb62cbbSBarry Smith   /* create local vector that is used to scatter into */
103029af93fSBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,ec,&aij->lvec);CHKERRQ(ierr);
1041eb62cbbSBarry Smith 
105d6dfbf8fSBarry Smith   /* create two temporary Index sets for build scatter gather */
10603919abeSBarry Smith   /*  check for the special case where blocks are communicated for faster VecScatterXXX */
107*a2f3521dSMark F. Adams   useblockis = PETSC_FALSE;
108d0f46423SBarry Smith   if (mat->rmap->bs > 1) {
109d0f46423SBarry Smith     PetscInt bs = mat->rmap->bs,ibs,ga;
11003919abeSBarry Smith     if (!(ec % bs)) {
11103919abeSBarry Smith       for (i=0; i<ec/bs; i++) {
11203919abeSBarry Smith         if ((ga = garray[ibs = i*bs]) % bs) {
11303919abeSBarry Smith           useblockis = PETSC_FALSE;
11403919abeSBarry Smith           break;
11503919abeSBarry Smith         }
11603919abeSBarry Smith         for (j=1; j<bs; j++) {
11703919abeSBarry Smith           if (garray[ibs+j] != ga+j) {
11803919abeSBarry Smith             useblockis = PETSC_FALSE;
11903919abeSBarry Smith             break;
12003919abeSBarry Smith           }
12103919abeSBarry Smith         }
12203919abeSBarry Smith         if (!useblockis) break;
12303919abeSBarry Smith       }
12403919abeSBarry Smith     }
12503919abeSBarry Smith   }
12603919abeSBarry Smith   if (useblockis) {
127d0f46423SBarry Smith     PetscInt *ga,bs = mat->rmap->bs,iec = ec/bs;
128ae15b995SBarry Smith     ierr = PetscInfo(mat,"Using block index set to define scatter\n");
129d0f46423SBarry Smith     ierr = PetscMalloc((ec/mat->rmap->bs)*sizeof(PetscInt),&ga);CHKERRQ(ierr);
130e82e9f6bSBarry Smith     for (i=0; i<iec; i++) ga[i] = garray[i*bs]/bs;
131deff0451SBarry Smith     ierr = ISCreateBlock(((PetscObject)mat)->comm,bs,iec,ga,PETSC_OWN_POINTER,&from);CHKERRQ(ierr);
13203919abeSBarry Smith   } else {
13370b3c8c7SBarry Smith     ierr = ISCreateGeneral(((PetscObject)mat)->comm,ec,garray,PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
13403919abeSBarry Smith   }
135029af93fSBarry Smith   ierr = ISCreateStride(PETSC_COMM_SELF,ec,0,1,&to);CHKERRQ(ierr);
1361eb62cbbSBarry Smith 
1371eb62cbbSBarry Smith   /* create temporary global vector to generate scatter context */
138b5eb4454SBarry Smith   /* This does not allocate the array's memory so is efficient */
139778a2246SBarry Smith   ierr = VecCreateMPIWithArray(((PetscObject)mat)->comm,1,mat->cmap->n,mat->cmap->N,PETSC_NULL,&gvec);CHKERRQ(ierr);
1401eb62cbbSBarry Smith 
1412d336d48SLois Curfman McInnes   /* generate the scatter context */
14208480c60SBarry Smith   ierr = VecScatterCreate(gvec,from,aij->lvec,to,&aij->Mvctx);CHKERRQ(ierr);
14352e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,aij->Mvctx);CHKERRQ(ierr);
14452e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,aij->lvec);CHKERRQ(ierr);
14552e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,from);CHKERRQ(ierr);
14652e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,to);CHKERRQ(ierr);
1479e25ed09SBarry Smith   aij->garray = garray;
14852e6d16bSBarry Smith   ierr = PetscLogObjectMemory(mat,(ec+1)*sizeof(PetscInt));CHKERRQ(ierr);
1496bf464f9SBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
1506bf464f9SBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
1516bf464f9SBarry Smith   ierr = VecDestroy(&gvec);CHKERRQ(ierr);
1523a40ed3dSBarry Smith   PetscFunctionReturn(0);
1538c79f6d3SBarry Smith }
1549e25ed09SBarry Smith 
1559e25ed09SBarry Smith 
1564a2ae208SSatish Balay #undef __FUNCT__
157ab9863d7SBarry Smith #define __FUNCT__ "MatDisAssemble_MPIAIJ"
1582493cbb0SBarry Smith /*
1592493cbb0SBarry Smith      Takes the local part of an already assembled MPIAIJ matrix
1602493cbb0SBarry Smith    and disassembles it. This is to allow new nonzeros into the matrix
1612493cbb0SBarry Smith    that require more communication in the matrix vector multiply.
1622493cbb0SBarry Smith    Thus certain data-structures must be rebuilt.
1632493cbb0SBarry Smith 
1642493cbb0SBarry Smith    Kind of slow! But that's what application programmers get when
1652493cbb0SBarry Smith    they are sloppy.
1662493cbb0SBarry Smith */
167ab9863d7SBarry Smith PetscErrorCode MatDisAssemble_MPIAIJ(Mat A)
1682493cbb0SBarry Smith {
1692493cbb0SBarry Smith   Mat_MPIAIJ     *aij = (Mat_MPIAIJ*)A->data;
1702493cbb0SBarry Smith   Mat            B = aij->B,Bnew;
171ec8511deSBarry Smith   Mat_SeqAIJ     *Baij = (Mat_SeqAIJ*)B->data;
172dfbe8321SBarry Smith   PetscErrorCode ierr;
173d0f46423SBarry Smith   PetscInt       i,j,m = B->rmap->n,n = A->cmap->N,col,ct = 0,*garray = aij->garray,*nz,ec;
17487828ca2SBarry Smith   PetscScalar    v;
1752493cbb0SBarry Smith 
1763a40ed3dSBarry Smith   PetscFunctionBegin;
1772493cbb0SBarry Smith   /* free stuff related to matrix-vec multiply */
178b0a32e0cSBarry Smith   ierr = VecGetSize(aij->lvec,&ec);CHKERRQ(ierr); /* needed for PetscLogObjectMemory below */
1795e1f6667SBarry Smith   ierr = VecDestroy(&aij->lvec);CHKERRQ(ierr);
1805e1f6667SBarry Smith   ierr = VecScatterDestroy(&aij->Mvctx);CHKERRQ(ierr);
181464493b3SBarry Smith   if (aij->colmap) {
182aa482453SBarry Smith #if defined (PETSC_USE_CTABLE)
1836bc0bbbfSBarry Smith     ierr = PetscTableDestroy(&aij->colmap);CHKERRQ(ierr);
1842066d6f7SSatish Balay #else
185606d414cSSatish Balay     ierr = PetscFree(aij->colmap);CHKERRQ(ierr);
186d0f46423SBarry Smith     ierr = PetscLogObjectMemory(A,-aij->B->cmap->n*sizeof(PetscInt));CHKERRQ(ierr);
1872066d6f7SSatish Balay #endif
188464493b3SBarry Smith   }
1892493cbb0SBarry Smith 
1902493cbb0SBarry Smith   /* make sure that B is assembled so we can access its values */
1916d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
192fe2f2677SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1932493cbb0SBarry Smith 
1942493cbb0SBarry Smith   /* invent new B and copy stuff over */
195b1d57f15SBarry Smith   ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nz);CHKERRQ(ierr);
19648b35521SBarry Smith   for (i=0; i<m; i++) {
19748b35521SBarry Smith     nz[i] = Baij->i[i+1] - Baij->i[i];
19848b35521SBarry Smith   }
199f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,&Bnew);CHKERRQ(ierr);
200f69a0ea3SMatthew Knepley   ierr = MatSetSizes(Bnew,m,n,m,n);CHKERRQ(ierr);
201*a2f3521dSMark F. Adams   ierr = MatSetBlockSizes(Bnew,A->rmap->bs,A->cmap->bs);CHKERRQ(ierr);
2027adad957SLisandro Dalcin   ierr = MatSetType(Bnew,((PetscObject)B)->type_name);CHKERRQ(ierr);
203f204ca49SKris Buschelman   ierr = MatSeqAIJSetPreallocation(Bnew,0,nz);CHKERRQ(ierr);
2042576faa2SJed Brown   ((Mat_SeqAIJ*)Bnew->data)->nonew = Baij->nonew; /* Inherit insertion error options. */
205606d414cSSatish Balay   ierr = PetscFree(nz);CHKERRQ(ierr);
2062493cbb0SBarry Smith   for (i=0; i<m; i++) {
207bfec09a0SHong Zhang     for (j=Baij->i[i]; j<Baij->i[i+1]; j++) {
208bfec09a0SHong Zhang       col  = garray[Baij->j[ct]];
2092493cbb0SBarry Smith       v    = Baij->a[ct++];
21083271157SBarry Smith       ierr = MatSetValues(Bnew,1,&i,1,&col,&v,B->insertmode);CHKERRQ(ierr);
2112493cbb0SBarry Smith     }
2122493cbb0SBarry Smith   }
213606d414cSSatish Balay   ierr = PetscFree(aij->garray);CHKERRQ(ierr);
21452e6d16bSBarry Smith   ierr = PetscLogObjectMemory(A,-ec*sizeof(PetscInt));CHKERRQ(ierr);
2156bf464f9SBarry Smith   ierr = MatDestroy(&B);CHKERRQ(ierr);
21652e6d16bSBarry Smith   ierr = PetscLogObjectParent(A,Bnew);CHKERRQ(ierr);
2172493cbb0SBarry Smith   aij->B = Bnew;
218227d817aSBarry Smith   A->was_assembled = PETSC_FALSE;
2193a40ed3dSBarry Smith   PetscFunctionReturn(0);
2202493cbb0SBarry Smith }
2212493cbb0SBarry Smith 
2222cd6534aSBarry Smith /*      ugly stuff added for Glenn someday we should fix this up */
2232cd6534aSBarry Smith 
224b1d57f15SBarry Smith static PetscInt *auglyrmapd = 0,*auglyrmapo = 0;  /* mapping from the local ordering to the "diagonal" and "off-diagonal"
2252cd6534aSBarry Smith                                       parts of the local matrix */
2262cd6534aSBarry Smith static Vec auglydd = 0,auglyoo = 0;   /* work vectors used to scale the two parts of the local matrix */
2272cd6534aSBarry Smith 
2282cd6534aSBarry Smith 
2292cd6534aSBarry Smith #undef __FUNCT__
2302cd6534aSBarry Smith #define __FUNCT__ "MatMPIAIJDiagonalScaleLocalSetUp"
231dfbe8321SBarry Smith PetscErrorCode MatMPIAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale)
2322cd6534aSBarry Smith {
2332cd6534aSBarry Smith   Mat_MPIAIJ     *ina = (Mat_MPIAIJ*) inA->data; /*access private part of matrix */
234dfbe8321SBarry Smith   PetscErrorCode ierr;
235b1d57f15SBarry Smith   PetscInt       i,n,nt,cstart,cend,no,*garray = ina->garray,*lindices;
236b1d57f15SBarry Smith   PetscInt       *r_rmapd,*r_rmapo;
2372cd6534aSBarry Smith 
2382cd6534aSBarry Smith   PetscFunctionBegin;
2392cd6534aSBarry Smith   ierr = MatGetOwnershipRange(inA,&cstart,&cend);CHKERRQ(ierr);
2402cd6534aSBarry Smith   ierr = MatGetSize(ina->A,PETSC_NULL,&n);CHKERRQ(ierr);
241992144d0SBarry Smith   ierr = PetscMalloc((inA->rmap->mapping->n+1)*sizeof(PetscInt),&r_rmapd);CHKERRQ(ierr);
242992144d0SBarry Smith   ierr = PetscMemzero(r_rmapd,inA->rmap->mapping->n*sizeof(PetscInt));CHKERRQ(ierr);
2432cd6534aSBarry Smith   nt   = 0;
244992144d0SBarry Smith   for (i=0; i<inA->rmap->mapping->n; i++) {
245992144d0SBarry Smith     if (inA->rmap->mapping->indices[i] >= cstart && inA->rmap->mapping->indices[i] < cend) {
2462cd6534aSBarry Smith       nt++;
247992144d0SBarry Smith       r_rmapd[i] = inA->rmap->mapping->indices[i] + 1;
2482cd6534aSBarry Smith     }
2492cd6534aSBarry Smith   }
250e32f2f54SBarry Smith   if (nt != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D n %D",nt,n);
251b1d57f15SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&auglyrmapd);CHKERRQ(ierr);
252992144d0SBarry Smith   for (i=0; i<inA->rmap->mapping->n; i++) {
2532cd6534aSBarry Smith     if (r_rmapd[i]){
2542cd6534aSBarry Smith       auglyrmapd[(r_rmapd[i]-1)-cstart] = i;
2552cd6534aSBarry Smith     }
2562cd6534aSBarry Smith   }
2572cd6534aSBarry Smith   ierr = PetscFree(r_rmapd);CHKERRQ(ierr);
2582cd6534aSBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,n,&auglydd);CHKERRQ(ierr);
2592cd6534aSBarry Smith 
260d0f46423SBarry Smith   ierr = PetscMalloc((inA->cmap->N+1)*sizeof(PetscInt),&lindices);CHKERRQ(ierr);
261d0f46423SBarry Smith   ierr = PetscMemzero(lindices,inA->cmap->N*sizeof(PetscInt));CHKERRQ(ierr);
262d0f46423SBarry Smith   for (i=0; i<ina->B->cmap->n; i++) {
2632cd6534aSBarry Smith     lindices[garray[i]] = i+1;
2642cd6534aSBarry Smith   }
265992144d0SBarry Smith   no   = inA->rmap->mapping->n - nt;
266992144d0SBarry Smith   ierr = PetscMalloc((inA->rmap->mapping->n+1)*sizeof(PetscInt),&r_rmapo);CHKERRQ(ierr);
267992144d0SBarry Smith   ierr = PetscMemzero(r_rmapo,inA->rmap->mapping->n*sizeof(PetscInt));CHKERRQ(ierr);
2682cd6534aSBarry Smith   nt   = 0;
269992144d0SBarry Smith   for (i=0; i<inA->rmap->mapping->n; i++) {
270992144d0SBarry Smith     if (lindices[inA->rmap->mapping->indices[i]]) {
2712cd6534aSBarry Smith       nt++;
272992144d0SBarry Smith       r_rmapo[i] = lindices[inA->rmap->mapping->indices[i]];
2732cd6534aSBarry Smith     }
2742cd6534aSBarry Smith   }
275e32f2f54SBarry Smith   if (nt > no) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D no %D",nt,n);
2762cd6534aSBarry Smith   ierr = PetscFree(lindices);CHKERRQ(ierr);
277b1d57f15SBarry Smith   ierr = PetscMalloc((nt+1)*sizeof(PetscInt),&auglyrmapo);CHKERRQ(ierr);
278992144d0SBarry Smith   for (i=0; i<inA->rmap->mapping->n; i++) {
2792cd6534aSBarry Smith     if (r_rmapo[i]){
2802cd6534aSBarry Smith       auglyrmapo[(r_rmapo[i]-1)] = i;
2812cd6534aSBarry Smith     }
2822cd6534aSBarry Smith   }
2832cd6534aSBarry Smith   ierr = PetscFree(r_rmapo);CHKERRQ(ierr);
2842cd6534aSBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,nt,&auglyoo);CHKERRQ(ierr);
2852cd6534aSBarry Smith 
2862cd6534aSBarry Smith   PetscFunctionReturn(0);
2872cd6534aSBarry Smith }
2882cd6534aSBarry Smith 
2892cd6534aSBarry Smith #undef __FUNCT__
2902cd6534aSBarry Smith #define __FUNCT__ "MatMPIAIJDiagonalScaleLocal"
291dfbe8321SBarry Smith PetscErrorCode MatMPIAIJDiagonalScaleLocal(Mat A,Vec scale)
2922cd6534aSBarry Smith {
29392b32695SKris Buschelman   /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */
2944ac538c5SBarry Smith   PetscErrorCode ierr;
29592b32695SKris Buschelman 
29692b32695SKris Buschelman   PetscFunctionBegin;
2974ac538c5SBarry Smith   ierr = PetscTryMethod(A,"MatDiagonalScaleLocal_C",(Mat,Vec),(A,scale));CHKERRQ(ierr);
29892b32695SKris Buschelman   PetscFunctionReturn(0);
29992b32695SKris Buschelman }
30092b32695SKris Buschelman 
30192b32695SKris Buschelman EXTERN_C_BEGIN
30292b32695SKris Buschelman #undef __FUNCT__
30392b32695SKris Buschelman #define __FUNCT__ "MatDiagonalScaleLocal_MPIAIJ"
3047087cfbeSBarry Smith PetscErrorCode  MatDiagonalScaleLocal_MPIAIJ(Mat A,Vec scale)
30592b32695SKris Buschelman {
3062cd6534aSBarry Smith   Mat_MPIAIJ     *a = (Mat_MPIAIJ*) A->data; /*access private part of matrix */
307dfbe8321SBarry Smith   PetscErrorCode ierr;
308b1d57f15SBarry Smith   PetscInt       n,i;
3092cd6534aSBarry Smith   PetscScalar    *d,*o,*s;
3102cd6534aSBarry Smith 
3112cd6534aSBarry Smith   PetscFunctionBegin;
3122cd6534aSBarry Smith   if (!auglyrmapd) {
3132cd6534aSBarry Smith     ierr = MatMPIAIJDiagonalScaleLocalSetUp(A,scale);CHKERRQ(ierr);
3142cd6534aSBarry Smith   }
3152cd6534aSBarry Smith 
3161ebc52fbSHong Zhang   ierr = VecGetArray(scale,&s);CHKERRQ(ierr);
3172cd6534aSBarry Smith 
3182cd6534aSBarry Smith   ierr = VecGetLocalSize(auglydd,&n);CHKERRQ(ierr);
3191ebc52fbSHong Zhang   ierr = VecGetArray(auglydd,&d);CHKERRQ(ierr);
3202cd6534aSBarry Smith   for (i=0; i<n; i++) {
3212cd6534aSBarry Smith     d[i] = s[auglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */
3222cd6534aSBarry Smith   }
3231ebc52fbSHong Zhang   ierr = VecRestoreArray(auglydd,&d);CHKERRQ(ierr);
3242cd6534aSBarry Smith   /* column scale "diagonal" portion of local matrix */
3252cd6534aSBarry Smith   ierr = MatDiagonalScale(a->A,PETSC_NULL,auglydd);CHKERRQ(ierr);
3262cd6534aSBarry Smith 
3272cd6534aSBarry Smith   ierr = VecGetLocalSize(auglyoo,&n);CHKERRQ(ierr);
3281ebc52fbSHong Zhang   ierr = VecGetArray(auglyoo,&o);CHKERRQ(ierr);
3292cd6534aSBarry Smith   for (i=0; i<n; i++) {
3302cd6534aSBarry Smith     o[i] = s[auglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */
3312cd6534aSBarry Smith   }
3321ebc52fbSHong Zhang   ierr = VecRestoreArray(scale,&s);CHKERRQ(ierr);
3331ebc52fbSHong Zhang   ierr = VecRestoreArray(auglyoo,&o);CHKERRQ(ierr);
3342cd6534aSBarry Smith   /* column scale "off-diagonal" portion of local matrix */
3352cd6534aSBarry Smith   ierr = MatDiagonalScale(a->B,PETSC_NULL,auglyoo);CHKERRQ(ierr);
3362cd6534aSBarry Smith 
3372cd6534aSBarry Smith   PetscFunctionReturn(0);
3382cd6534aSBarry Smith }
33992b32695SKris Buschelman EXTERN_C_END
3402cd6534aSBarry Smith 
34148b35521SBarry Smith 
342