xref: /petsc/src/mat/tests/ex81.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
18c78258cSHong Zhang 
28c78258cSHong Zhang static char help[] = "Tests MatOption MAT_FORCE_DIAGONAL_ENTRIES.\n\n";
38c78258cSHong Zhang 
48c78258cSHong Zhang #include <petscmat.h>
58c78258cSHong Zhang 
68c78258cSHong Zhang int main(int argc,char **args)
78c78258cSHong Zhang {
88c78258cSHong Zhang   Mat            A,B;
98c78258cSHong Zhang   Vec            diag;
108c78258cSHong Zhang   PetscInt       i,n = 10,col[3],test;
118c78258cSHong Zhang   PetscScalar    v[3];
128c78258cSHong Zhang 
13*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
145f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
158c78258cSHong Zhang 
168c78258cSHong Zhang   /* Create A which has empty 0-th row and column */
175f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
185f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n));
195f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(A,MATAIJ));
205f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(A));
215f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(A));
228c78258cSHong Zhang 
238c78258cSHong Zhang   v[0] = -1.; v[1] = 2.; v[2] = -1.;
248c78258cSHong Zhang   for (i=2; i<n-1; i++) {
258c78258cSHong Zhang     col[0] = i-1; col[1] = i; col[2] = i+1;
265f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(A,1,&i,3,col,v,INSERT_VALUES));
278c78258cSHong Zhang   }
288c78258cSHong Zhang   i    = 1; col[0] = 1; col[1] = 2;
295f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(A,1,&i,2,col,v+1,INSERT_VALUES));
308c78258cSHong Zhang   i    = n - 1; col[0] = n - 2; col[1] = n - 1;
315f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(A,1,&i,2,col,v,INSERT_VALUES));
325f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
335f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
348c78258cSHong Zhang 
358c78258cSHong Zhang   for (test = 0; test < 2; test++) {
365f80ce2aSJacob Faibussowitsch     CHKERRQ(MatProductCreate(A,A,NULL,&B));
378c78258cSHong Zhang 
388c78258cSHong Zhang     if (test == 0) {
398c78258cSHong Zhang       /* Compute B = A*A; B misses 0-th diagonal */
405f80ce2aSJacob Faibussowitsch       CHKERRQ(MatProductSetType(B,MATPRODUCT_AB));
415f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetOptionsPrefix(B,"AB_"));
428c78258cSHong Zhang     } else {
438c78258cSHong Zhang       /* Compute B = A^t*A; B misses 0-th diagonal */
445f80ce2aSJacob Faibussowitsch       CHKERRQ(MatProductSetType(B,MATPRODUCT_AtB));
455f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetOptionsPrefix(B,"AtB_"));
468c78258cSHong Zhang     }
478c78258cSHong Zhang 
488c78258cSHong Zhang     /* Force allocate missing diagonal entries of B */
495f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetOption(B,MAT_FORCE_DIAGONAL_ENTRIES,PETSC_TRUE));
505f80ce2aSJacob Faibussowitsch     CHKERRQ(MatProductSetFromOptions(B));
518c78258cSHong Zhang 
525f80ce2aSJacob Faibussowitsch     CHKERRQ(MatProductSymbolic(B));
535f80ce2aSJacob Faibussowitsch     CHKERRQ(MatProductNumeric(B));
548c78258cSHong Zhang 
555f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE));
568c78258cSHong Zhang 
578c78258cSHong Zhang     /* Insert entries to diagonal of B */
585f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateVecs(B,NULL,&diag));
595f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetDiagonal(B,diag));
605f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetValue(diag,0,100.0,INSERT_VALUES));
615f80ce2aSJacob Faibussowitsch     CHKERRQ(VecAssemblyBegin(diag));
625f80ce2aSJacob Faibussowitsch     CHKERRQ(VecAssemblyEnd(diag));
638c78258cSHong Zhang 
645f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDiagonalSet(B,diag,INSERT_VALUES));
658c78258cSHong Zhang     if (test == 1) {
665f80ce2aSJacob Faibussowitsch       CHKERRQ(MatView(B,PETSC_VIEWER_STDOUT_WORLD));
678c78258cSHong Zhang     }
685f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&B));
695f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&diag));
708c78258cSHong Zhang   }
718c78258cSHong Zhang 
725f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
73*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
74*b122ec5aSJacob Faibussowitsch   return 0;
758c78258cSHong Zhang }
768c78258cSHong Zhang 
778c78258cSHong Zhang /*TEST
788c78258cSHong Zhang 
798c78258cSHong Zhang    test:
808c78258cSHong Zhang      output_file: output/ex81_1.out
818c78258cSHong Zhang 
828c78258cSHong Zhang    test:
838c78258cSHong Zhang      suffix: 2
843e662e0bSHong Zhang      args: -AtB_mat_product_algorithm at*b
858c78258cSHong Zhang      output_file: output/ex81_1.out
868c78258cSHong Zhang 
878c78258cSHong Zhang    test:
888c78258cSHong Zhang      suffix: 3
893e662e0bSHong Zhang      args: -AtB_mat_product_algorithm outerproduct
908c78258cSHong Zhang      output_file: output/ex81_1.out
918c78258cSHong Zhang 
928c78258cSHong Zhang    test:
938c78258cSHong Zhang      suffix: 4
948c78258cSHong Zhang      nsize: 3
953e662e0bSHong Zhang      args: -AtB_mat_product_algorithm nonscalable
968c78258cSHong Zhang      output_file: output/ex81_3.out
978c78258cSHong Zhang 
988c78258cSHong Zhang    test:
998c78258cSHong Zhang      suffix: 5
1008c78258cSHong Zhang      nsize: 3
1013e662e0bSHong Zhang      args: -AtB_mat_product_algorithm scalable
1028c78258cSHong Zhang      output_file: output/ex81_3.out
1038c78258cSHong Zhang 
1048c78258cSHong Zhang    test:
1058c78258cSHong Zhang      suffix: 6
1068c78258cSHong Zhang      nsize: 3
1073e662e0bSHong Zhang      args: -AtB_mat_product_algorithm at*b
1088c78258cSHong Zhang      output_file: output/ex81_3.out
1098c78258cSHong Zhang 
1108c78258cSHong Zhang TEST*/
111