xref: /petsc/src/mat/tests/ex52.c (revision 1baa6e3354bfe224b33a0c158f482508792a8a6e)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests various routines in MatMPIBAIJ 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;
9c4762a1bSJed Brown   PetscInt       m=2,bs=1,M,row,col,start,end,i,j,k;
10c4762a1bSJed Brown   PetscMPIInt    rank,size;
11c4762a1bSJed Brown   PetscScalar    data=100;
12c4762a1bSJed Brown   PetscBool      flg;
13c4762a1bSJed Brown 
149566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
159566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
169566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
17c4762a1bSJed Brown 
18c4762a1bSJed Brown   /* Test MatSetValues() and MatGetValues() */
199566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-mat_block_size",&bs,NULL));
209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-mat_size",&m,NULL));
21c4762a1bSJed Brown 
22c4762a1bSJed Brown   M    = m*bs*size;
239566063dSJacob Faibussowitsch   PetscCall(MatCreateBAIJ(PETSC_COMM_WORLD,bs,PETSC_DECIDE,PETSC_DECIDE,M,M,PETSC_DECIDE,NULL,PETSC_DECIDE,NULL,&A));
24c4762a1bSJed Brown 
259566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A,&start,&end));
269566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL,"-column_oriented",&flg));
27*1baa6e33SBarry Smith   if (flg) PetscCall(MatSetOption(A,MAT_ROW_ORIENTED,PETSC_FALSE));
28c4762a1bSJed Brown 
29c4762a1bSJed Brown   /* inproc assembly */
30c4762a1bSJed Brown   for (row=start; row<end; row++) {
31c4762a1bSJed Brown     for (col=start; col<end; col++,data+=1) {
329566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&row,1,&col,&data,INSERT_VALUES));
33c4762a1bSJed Brown     }
34c4762a1bSJed Brown   }
359566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
369566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
37c4762a1bSJed Brown 
38c4762a1bSJed Brown   /* offproc assembly */
39c4762a1bSJed Brown   data = 5.0;
40c4762a1bSJed Brown   row  = (M+start-1)%M;
41c4762a1bSJed Brown   for (col=0; col<M; col++) {
429566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A,1,&row,1,&col,&data,ADD_VALUES));
43c4762a1bSJed Brown   }
449566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
459566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
46c4762a1bSJed Brown 
47c4762a1bSJed Brown   /* Test MatSetValuesBlocked() */
489566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL,"-test_setvaluesblocked",&flg));
49c4762a1bSJed Brown   if (flg) {
50c4762a1bSJed Brown     PetscScalar *bval;
51c4762a1bSJed Brown     row /= bs;
52c4762a1bSJed Brown     col  = start/bs;
539566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(bs*bs,&bval));
54c4762a1bSJed Brown     k = 1;
55a5b23f4aSJose E. Roman     /* row oriented - default */
56c4762a1bSJed Brown     for (i=0; i<bs; i++) {
57c4762a1bSJed Brown       for (j=0; j<bs; j++) {
58c4762a1bSJed Brown         bval[i*bs+j] = (PetscScalar)k; k++;
59c4762a1bSJed Brown       }
60c4762a1bSJed Brown     }
619566063dSJacob Faibussowitsch     PetscCall(MatSetValuesBlocked(A,1,&row,1,&col,bval,INSERT_VALUES));
629566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
639566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
649566063dSJacob Faibussowitsch     PetscCall(PetscFree(bval));
65c4762a1bSJed Brown   }
66c4762a1bSJed Brown 
679566063dSJacob Faibussowitsch   PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
689566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
699566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
70b122ec5aSJacob Faibussowitsch   return 0;
71c4762a1bSJed Brown }
72c4762a1bSJed Brown 
73c4762a1bSJed Brown /*TEST
74c4762a1bSJed Brown 
75c4762a1bSJed Brown    test:
76c4762a1bSJed Brown       suffix: 1
77c4762a1bSJed Brown       nsize: 3
78c4762a1bSJed Brown       args: -mat_block_size 2 -test_setvaluesblocked
79c4762a1bSJed Brown 
80c4762a1bSJed Brown    test:
81c4762a1bSJed Brown       suffix: 2
82c4762a1bSJed Brown       nsize: 3
83c4762a1bSJed Brown       args: -mat_block_size 2 -test_setvaluesblocked -column_oriented
84c4762a1bSJed Brown 
85c4762a1bSJed Brown    test:
86c4762a1bSJed Brown       suffix: 3
87c4762a1bSJed Brown       nsize: 3
88c4762a1bSJed Brown       args: -mat_block_size 1 -test_setvaluesblocked
89c4762a1bSJed Brown 
90c4762a1bSJed Brown    test:
91c4762a1bSJed Brown       suffix: 4
92c4762a1bSJed Brown       nsize: 3
93c4762a1bSJed Brown       args: -mat_block_size 1 -test_setvaluesblocked -column_oriented
94c4762a1bSJed Brown 
95c4762a1bSJed Brown TEST*/
96