xref: /petsc/src/mat/tests/ex91.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests MatIncreaseOverlap(), MatCreateSubMatrices() for sequential MatSBAIJ format. Derived from ex51.c\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown int main(int argc,char **args)
7c4762a1bSJed Brown {
8c4762a1bSJed Brown   Mat            A,Atrans,sA,*submatA,*submatsA;
9c4762a1bSJed Brown   PetscInt       bs=1,m=43,ov=1,i,j,k,*rows,*cols,M,nd=5,*idx,mm,nn;
10c4762a1bSJed Brown   PetscMPIInt    size;
11c4762a1bSJed Brown   PetscScalar    *vals,rval,one=1.0;
12c4762a1bSJed Brown   IS             *is1,*is2;
13c4762a1bSJed Brown   PetscRandom    rand;
14c4762a1bSJed Brown   Vec            xx,s1,s2;
15c4762a1bSJed Brown   PetscReal      s1norm,s2norm,rnorm,tol = 10*PETSC_SMALL;
16c4762a1bSJed Brown   PetscBool      flg;
17c4762a1bSJed Brown 
18*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
195f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-mat_block_size",&bs,NULL));
205f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-mat_size",&m,NULL));
215f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-ov",&ov,NULL));
225f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-nd",&nd,NULL));
23c4762a1bSJed Brown 
24c4762a1bSJed Brown   /* create a SeqBAIJ matrix A */
25c4762a1bSJed Brown   M    = m*bs;
265f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSeqBAIJ(PETSC_COMM_SELF,bs,M,M,1,NULL,&A));
275f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE));
285f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PETSC_COMM_SELF,&rand));
295f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetFromOptions(rand));
30c4762a1bSJed Brown 
315f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(bs,&rows));
325f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(bs,&cols));
335f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(bs*bs,&vals));
345f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(M,&idx));
35c4762a1bSJed Brown 
36c4762a1bSJed Brown   /* Now set blocks of random values */
37c4762a1bSJed Brown   /* first, set diagonal blocks as zero */
38c4762a1bSJed Brown   for (j=0; j<bs*bs; j++) vals[j] = 0.0;
39c4762a1bSJed Brown   for (i=0; i<m; i++) {
40c4762a1bSJed Brown     cols[0] = i*bs; rows[0] = i*bs;
41c4762a1bSJed Brown     for (j=1; j<bs; j++) {
42c4762a1bSJed Brown       rows[j] = rows[j-1]+1;
43c4762a1bSJed Brown       cols[j] = cols[j-1]+1;
44c4762a1bSJed Brown     }
455f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(A,bs,rows,bs,cols,vals,ADD_VALUES));
46c4762a1bSJed Brown   }
47c4762a1bSJed Brown   /* second, add random blocks */
48c4762a1bSJed Brown   for (i=0; i<20*bs; i++) {
495f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomGetValue(rand,&rval));
50c4762a1bSJed Brown     cols[0] = bs*(int)(PetscRealPart(rval)*m);
515f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomGetValue(rand,&rval));
52c4762a1bSJed Brown     rows[0] = bs*(int)(PetscRealPart(rval)*m);
53c4762a1bSJed Brown     for (j=1; j<bs; j++) {
54c4762a1bSJed Brown       rows[j] = rows[j-1]+1;
55c4762a1bSJed Brown       cols[j] = cols[j-1]+1;
56c4762a1bSJed Brown     }
57c4762a1bSJed Brown 
58c4762a1bSJed Brown     for (j=0; j<bs*bs; j++) {
595f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscRandomGetValue(rand,&rval));
60c4762a1bSJed Brown       vals[j] = rval;
61c4762a1bSJed Brown     }
625f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(A,bs,rows,bs,cols,vals,ADD_VALUES));
63c4762a1bSJed Brown   }
64c4762a1bSJed Brown 
655f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
665f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
67c4762a1bSJed Brown 
68c4762a1bSJed Brown   /* make A a symmetric matrix: A <- A^T + A */
695f80ce2aSJacob Faibussowitsch   CHKERRQ(MatTranspose(A,MAT_INITIAL_MATRIX, &Atrans));
705f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAXPY(A,one,Atrans,DIFFERENT_NONZERO_PATTERN));
715f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&Atrans));
725f80ce2aSJacob Faibussowitsch   CHKERRQ(MatTranspose(A,MAT_INITIAL_MATRIX, &Atrans));
735f80ce2aSJacob Faibussowitsch   CHKERRQ(MatEqual(A, Atrans, &flg));
7428b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"A+A^T is non-symmetric");
755f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&Atrans));
76c4762a1bSJed Brown 
77c4762a1bSJed Brown   /* create a SeqSBAIJ matrix sA (= A) */
785f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE));
795f80ce2aSJacob Faibussowitsch   CHKERRQ(MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&sA));
80c4762a1bSJed Brown 
81c4762a1bSJed Brown   /* Test sA==A through MatMult() */
82c4762a1bSJed Brown   for (i=0; i<nd; i++) {
835f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetSize(A,&mm,&nn));
845f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,mm,&xx));
855f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(xx,&s1));
865f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(xx,&s2));
87c4762a1bSJed Brown     for (j=0; j<3; j++) {
885f80ce2aSJacob Faibussowitsch       CHKERRQ(VecSetRandom(xx,rand));
895f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMult(A,xx,s1));
905f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMult(sA,xx,s2));
915f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(s1,NORM_2,&s1norm));
925f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(s2,NORM_2,&s2norm));
93c4762a1bSJed Brown       rnorm = s2norm-s1norm;
94c4762a1bSJed Brown       if (rnorm<-tol || rnorm>tol) {
955f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n",(double)s1norm,(double)s2norm));
96c4762a1bSJed Brown       }
97c4762a1bSJed Brown     }
985f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&xx));
995f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&s1));
1005f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&s2));
101c4762a1bSJed Brown   }
102c4762a1bSJed Brown 
103c4762a1bSJed Brown   /* Test MatIncreaseOverlap() */
1045f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nd,&is1));
1055f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nd,&is2));
106c4762a1bSJed Brown 
107c4762a1bSJed Brown   for (i=0; i<nd; i++) {
1085f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomGetValue(rand,&rval));
109c4762a1bSJed Brown     size = (int)(PetscRealPart(rval)*m);
110c4762a1bSJed Brown     for (j=0; j<size; j++) {
1115f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscRandomGetValue(rand,&rval));
112c4762a1bSJed Brown       idx[j*bs] = bs*(int)(PetscRealPart(rval)*m);
113c4762a1bSJed Brown       for (k=1; k<bs; k++) idx[j*bs+k] = idx[j*bs]+k;
114c4762a1bSJed Brown     }
1155f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,size*bs,idx,PETSC_COPY_VALUES,is1+i));
1165f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,size*bs,idx,PETSC_COPY_VALUES,is2+i));
117c4762a1bSJed Brown   }
118c4762a1bSJed Brown   /* for debugging */
119c4762a1bSJed Brown   /*
1205f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_SELF));
1215f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(sA,PETSC_VIEWER_STDOUT_SELF));
122c4762a1bSJed Brown   */
123c4762a1bSJed Brown 
1245f80ce2aSJacob Faibussowitsch   CHKERRQ(MatIncreaseOverlap(A,nd,is1,ov));
1255f80ce2aSJacob Faibussowitsch   CHKERRQ(MatIncreaseOverlap(sA,nd,is2,ov));
126c4762a1bSJed Brown 
127c4762a1bSJed Brown   for (i=0; i<nd; ++i) {
1285f80ce2aSJacob Faibussowitsch     CHKERRQ(ISSort(is1[i]));
1295f80ce2aSJacob Faibussowitsch     CHKERRQ(ISSort(is2[i]));
130c4762a1bSJed Brown   }
131c4762a1bSJed Brown 
132c4762a1bSJed Brown   for (i=0; i<nd; ++i) {
1335f80ce2aSJacob Faibussowitsch     CHKERRQ(ISEqual(is1[i],is2[i],&flg));
13428b400f6SJacob Faibussowitsch     PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"i=%" PetscInt_FMT ", is1 != is2",i);
135c4762a1bSJed Brown   }
136c4762a1bSJed Brown 
1375f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSubMatrices(A,nd,is1,is1,MAT_INITIAL_MATRIX,&submatA));
1385f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSubMatrices(sA,nd,is2,is2,MAT_INITIAL_MATRIX,&submatsA));
139c4762a1bSJed Brown 
140c4762a1bSJed Brown   /* Test MatMult() */
141c4762a1bSJed Brown   for (i=0; i<nd; i++) {
1425f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetSize(submatA[i],&mm,&nn));
1435f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,mm,&xx));
1445f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(xx,&s1));
1455f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(xx,&s2));
146c4762a1bSJed Brown     for (j=0; j<3; j++) {
1475f80ce2aSJacob Faibussowitsch       CHKERRQ(VecSetRandom(xx,rand));
1485f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMult(submatA[i],xx,s1));
1495f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMult(submatsA[i],xx,s2));
1505f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(s1,NORM_2,&s1norm));
1515f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(s2,NORM_2,&s2norm));
152c4762a1bSJed Brown       rnorm = s2norm-s1norm;
153c4762a1bSJed Brown       if (rnorm<-tol || rnorm>tol) {
1545f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n",(double)s1norm,(double)s2norm));
155c4762a1bSJed Brown       }
156c4762a1bSJed Brown     }
1575f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&xx));
1585f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&s1));
1595f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&s2));
160c4762a1bSJed Brown   }
161c4762a1bSJed Brown 
162c4762a1bSJed Brown   /* Now test MatCreateSubmatrices with MAT_REUSE_MATRIX option */
1635f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSubMatrices(A,nd,is1,is1,MAT_REUSE_MATRIX,&submatA));
1645f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSubMatrices(sA,nd,is2,is2,MAT_REUSE_MATRIX,&submatsA));
165c4762a1bSJed Brown 
166c4762a1bSJed Brown   /* Test MatMult() */
167c4762a1bSJed Brown   for (i=0; i<nd; i++) {
1685f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetSize(submatA[i],&mm,&nn));
1695f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,mm,&xx));
1705f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(xx,&s1));
1715f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(xx,&s2));
172c4762a1bSJed Brown     for (j=0; j<3; j++) {
1735f80ce2aSJacob Faibussowitsch       CHKERRQ(VecSetRandom(xx,rand));
1745f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMult(submatA[i],xx,s1));
1755f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMult(submatsA[i],xx,s2));
1765f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(s1,NORM_2,&s1norm));
1775f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(s2,NORM_2,&s2norm));
178c4762a1bSJed Brown       rnorm = s2norm-s1norm;
179c4762a1bSJed Brown       if (rnorm<-tol || rnorm>tol) {
1805f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n",(double)s1norm,(double)s2norm));
181c4762a1bSJed Brown       }
182c4762a1bSJed Brown     }
1835f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&xx));
1845f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&s1));
1855f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&s2));
186c4762a1bSJed Brown   }
187c4762a1bSJed Brown 
188c4762a1bSJed Brown   /* Free allocated memory */
189c4762a1bSJed Brown   for (i=0; i<nd; ++i) {
1905f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&is1[i]));
1915f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&is2[i]));
192c4762a1bSJed Brown   }
1935f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroySubMatrices(nd,&submatA));
1945f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroySubMatrices(nd,&submatsA));
195c4762a1bSJed Brown 
1965f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(is1));
1975f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(is2));
1985f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(idx));
1995f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(rows));
2005f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(cols));
2015f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(vals));
2025f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
2035f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&sA));
2045f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&rand));
205*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
206*b122ec5aSJacob Faibussowitsch   return 0;
207c4762a1bSJed Brown }
208c4762a1bSJed Brown 
209c4762a1bSJed Brown /*TEST
210c4762a1bSJed Brown 
211c4762a1bSJed Brown    test:
212c4762a1bSJed Brown       args: -ov 2
213c4762a1bSJed Brown 
214c4762a1bSJed Brown TEST*/
215