xref: /petsc/src/mat/tests/ex51.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests MatIncreaseOverlap(), MatCreateSubMatrices() for MatBAIJ format.\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown int main(int argc,char **args)
7c4762a1bSJed Brown {
8c4762a1bSJed Brown   Mat            A,B,*submatA,*submatB;
9c4762a1bSJed Brown   PetscInt       bs=1,m=43,ov=1,i,j,k,*rows,*cols,M,nd=5,*idx,mm,nn,lsize;
10c4762a1bSJed Brown   PetscScalar    *vals,rval;
11c4762a1bSJed Brown   IS             *is1,*is2;
12c4762a1bSJed Brown   PetscRandom    rdm;
13c4762a1bSJed Brown   Vec            xx,s1,s2;
14c4762a1bSJed Brown   PetscReal      s1norm,s2norm,rnorm,tol = PETSC_SQRT_MACHINE_EPSILON;
15c4762a1bSJed Brown   PetscBool      flg;
16c4762a1bSJed Brown 
17*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
185f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-mat_block_size",&bs,NULL));
195f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-mat_size",&m,NULL));
205f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-ov",&ov,NULL));
215f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-nd",&nd,NULL));
22c4762a1bSJed Brown   M    = m*bs;
23c4762a1bSJed Brown 
245f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSeqBAIJ(PETSC_COMM_SELF,bs,M,M,1,NULL,&A));
255f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE));
265f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSeqAIJ(PETSC_COMM_SELF,M,M,15,NULL,&B));
275f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE));
285f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PETSC_COMM_SELF,&rdm));
295f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetFromOptions(rdm));
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 values */
37c4762a1bSJed Brown   for (i=0; i<20*bs; i++) {
385f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomGetValue(rdm,&rval));
39c4762a1bSJed Brown     cols[0] = bs*(int)(PetscRealPart(rval)*m);
405f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomGetValue(rdm,&rval));
41c4762a1bSJed Brown     rows[0] = bs*(int)(PetscRealPart(rval)*m);
42c4762a1bSJed Brown     for (j=1; j<bs; j++) {
43c4762a1bSJed Brown       rows[j] = rows[j-1]+1;
44c4762a1bSJed Brown       cols[j] = cols[j-1]+1;
45c4762a1bSJed Brown     }
46c4762a1bSJed Brown 
47c4762a1bSJed Brown     for (j=0; j<bs*bs; j++) {
485f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscRandomGetValue(rdm,&rval));
49c4762a1bSJed Brown       vals[j] = rval;
50c4762a1bSJed Brown     }
515f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(A,bs,rows,bs,cols,vals,ADD_VALUES));
525f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(B,bs,rows,bs,cols,vals,ADD_VALUES));
53c4762a1bSJed Brown   }
54c4762a1bSJed Brown 
555f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
565f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
575f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
585f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
59c4762a1bSJed Brown 
60c4762a1bSJed Brown   /* Test MatIncreaseOverlap() */
615f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nd,&is1));
625f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nd,&is2));
63c4762a1bSJed Brown 
64c4762a1bSJed Brown   for (i=0; i<nd; i++) {
655f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomGetValue(rdm,&rval));
66c4762a1bSJed Brown     lsize = (int)(PetscRealPart(rval)*m);
67c4762a1bSJed Brown     for (j=0; j<lsize; j++) {
685f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscRandomGetValue(rdm,&rval));
69c4762a1bSJed Brown       idx[j*bs] = bs*(int)(PetscRealPart(rval)*m);
70c4762a1bSJed Brown       for (k=1; k<bs; k++) idx[j*bs+k] = idx[j*bs]+k;
71c4762a1bSJed Brown     }
725f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,lsize*bs,idx,PETSC_COPY_VALUES,is1+i));
735f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,lsize*bs,idx,PETSC_COPY_VALUES,is2+i));
74c4762a1bSJed Brown   }
755f80ce2aSJacob Faibussowitsch   CHKERRQ(MatIncreaseOverlap(A,nd,is1,ov));
765f80ce2aSJacob Faibussowitsch   CHKERRQ(MatIncreaseOverlap(B,nd,is2,ov));
77c4762a1bSJed Brown 
78c4762a1bSJed Brown   for (i=0; i<nd; ++i) {
795f80ce2aSJacob Faibussowitsch     CHKERRQ(ISEqual(is1[i],is2[i],&flg));
8028b400f6SJacob Faibussowitsch     PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"i=%" PetscInt_FMT ", flg =%d",i,(int)flg);
81c4762a1bSJed Brown   }
82c4762a1bSJed Brown 
83c4762a1bSJed Brown   for (i=0; i<nd; ++i) {
845f80ce2aSJacob Faibussowitsch     CHKERRQ(ISSort(is1[i]));
855f80ce2aSJacob Faibussowitsch     CHKERRQ(ISSort(is2[i]));
86c4762a1bSJed Brown   }
87c4762a1bSJed Brown 
885f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSubMatrices(A,nd,is1,is1,MAT_INITIAL_MATRIX,&submatA));
895f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSubMatrices(B,nd,is2,is2,MAT_INITIAL_MATRIX,&submatB));
90c4762a1bSJed Brown 
91c4762a1bSJed Brown   /* Test MatMult() */
92c4762a1bSJed Brown   for (i=0; i<nd; i++) {
935f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetSize(submatA[i],&mm,&nn));
945f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,mm,&xx));
955f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(xx,&s1));
965f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(xx,&s2));
97c4762a1bSJed Brown     for (j=0; j<3; j++) {
985f80ce2aSJacob Faibussowitsch       CHKERRQ(VecSetRandom(xx,rdm));
995f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMult(submatA[i],xx,s1));
1005f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMult(submatB[i],xx,s2));
1015f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(s1,NORM_2,&s1norm));
1025f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(s2,NORM_2,&s2norm));
103c4762a1bSJed Brown       rnorm = s2norm-s1norm;
104c4762a1bSJed Brown       if (rnorm<-tol || rnorm>tol) {
1055f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n",(double)s1norm,(double)s2norm));
106c4762a1bSJed Brown       }
107c4762a1bSJed Brown     }
1085f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&xx));
1095f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&s1));
1105f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&s2));
111c4762a1bSJed Brown   }
112c4762a1bSJed Brown   /* Now test MatCreateSubmatrices with MAT_REUSE_MATRIX option */
1135f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSubMatrices(A,nd,is1,is1,MAT_REUSE_MATRIX,&submatA));
1145f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSubMatrices(B,nd,is2,is2,MAT_REUSE_MATRIX,&submatB));
115c4762a1bSJed Brown 
116c4762a1bSJed Brown   /* Test MatMult() */
117c4762a1bSJed Brown   for (i=0; i<nd; i++) {
1185f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetSize(submatA[i],&mm,&nn));
1195f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,mm,&xx));
1205f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(xx,&s1));
1215f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(xx,&s2));
122c4762a1bSJed Brown     for (j=0; j<3; j++) {
1235f80ce2aSJacob Faibussowitsch       CHKERRQ(VecSetRandom(xx,rdm));
1245f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMult(submatA[i],xx,s1));
1255f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMult(submatB[i],xx,s2));
1265f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(s1,NORM_2,&s1norm));
1275f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(s2,NORM_2,&s2norm));
128c4762a1bSJed Brown       rnorm = s2norm-s1norm;
129c4762a1bSJed Brown       if (rnorm<-tol || rnorm>tol) {
1305f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n",(double)s1norm,(double)s2norm));
131c4762a1bSJed Brown       }
132c4762a1bSJed Brown     }
1335f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&xx));
1345f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&s1));
1355f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&s2));
136c4762a1bSJed Brown   }
137c4762a1bSJed Brown 
138c4762a1bSJed Brown   /* Free allocated memory */
139c4762a1bSJed Brown   for (i=0; i<nd; ++i) {
1405f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&is1[i]));
1415f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&is2[i]));
142c4762a1bSJed Brown   }
1435f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroySubMatrices(nd,&submatA));
1445f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroySubMatrices(nd,&submatB));
1455f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(is1));
1465f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(is2));
1475f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(idx));
1485f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(rows));
1495f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(cols));
1505f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(vals));
1515f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
1525f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
1535f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&rdm));
154*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
155*b122ec5aSJacob Faibussowitsch   return 0;
156c4762a1bSJed Brown }
157c4762a1bSJed Brown 
158c4762a1bSJed Brown /*TEST
159c4762a1bSJed Brown 
160c4762a1bSJed Brown    test:
161c4762a1bSJed Brown       args: -mat_block_size {{1 2  5 7 8}} -ov {{1 3}} -mat_size {{11 13}} -nd {{7}}
162c4762a1bSJed Brown 
163c4762a1bSJed Brown TEST*/
164