xref: /petsc/src/mat/tests/ex81.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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   PetscErrorCode ierr;
128c78258cSHong Zhang   PetscScalar    v[3];
138c78258cSHong Zhang 
148c78258cSHong Zhang   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
15*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
168c78258cSHong Zhang 
178c78258cSHong Zhang   /* Create A which has empty 0-th row and column */
18*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
19*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n));
20*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(A,MATAIJ));
21*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(A));
22*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(A));
238c78258cSHong Zhang 
248c78258cSHong Zhang   v[0] = -1.; v[1] = 2.; v[2] = -1.;
258c78258cSHong Zhang   for (i=2; i<n-1; i++) {
268c78258cSHong Zhang     col[0] = i-1; col[1] = i; col[2] = i+1;
27*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(A,1,&i,3,col,v,INSERT_VALUES));
288c78258cSHong Zhang   }
298c78258cSHong Zhang   i    = 1; col[0] = 1; col[1] = 2;
30*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(A,1,&i,2,col,v+1,INSERT_VALUES));
318c78258cSHong Zhang   i    = n - 1; col[0] = n - 2; col[1] = n - 1;
32*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(A,1,&i,2,col,v,INSERT_VALUES));
33*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
34*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
358c78258cSHong Zhang 
368c78258cSHong Zhang   for (test = 0; test < 2; test++) {
37*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatProductCreate(A,A,NULL,&B));
388c78258cSHong Zhang 
398c78258cSHong Zhang     if (test == 0) {
408c78258cSHong Zhang       /* Compute B = A*A; B misses 0-th diagonal */
41*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatProductSetType(B,MATPRODUCT_AB));
42*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetOptionsPrefix(B,"AB_"));
438c78258cSHong Zhang     } else {
448c78258cSHong Zhang       /* Compute B = A^t*A; B misses 0-th diagonal */
45*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatProductSetType(B,MATPRODUCT_AtB));
46*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetOptionsPrefix(B,"AtB_"));
478c78258cSHong Zhang     }
488c78258cSHong Zhang 
498c78258cSHong Zhang     /* Force allocate missing diagonal entries of B */
50*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetOption(B,MAT_FORCE_DIAGONAL_ENTRIES,PETSC_TRUE));
51*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatProductSetFromOptions(B));
528c78258cSHong Zhang 
53*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatProductSymbolic(B));
54*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatProductNumeric(B));
558c78258cSHong Zhang 
56*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE));
578c78258cSHong Zhang 
588c78258cSHong Zhang     /* Insert entries to diagonal of B */
59*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateVecs(B,NULL,&diag));
60*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetDiagonal(B,diag));
61*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetValue(diag,0,100.0,INSERT_VALUES));
62*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecAssemblyBegin(diag));
63*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecAssemblyEnd(diag));
648c78258cSHong Zhang 
65*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDiagonalSet(B,diag,INSERT_VALUES));
668c78258cSHong Zhang     if (test == 1) {
67*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatView(B,PETSC_VIEWER_STDOUT_WORLD));
688c78258cSHong Zhang     }
69*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&B));
70*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&diag));
718c78258cSHong Zhang   }
728c78258cSHong Zhang 
73*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
748c78258cSHong Zhang   ierr = PetscFinalize();
758c78258cSHong Zhang   return ierr;
768c78258cSHong Zhang }
778c78258cSHong Zhang 
788c78258cSHong Zhang /*TEST
798c78258cSHong Zhang 
808c78258cSHong Zhang    test:
818c78258cSHong Zhang      output_file: output/ex81_1.out
828c78258cSHong Zhang 
838c78258cSHong Zhang    test:
848c78258cSHong Zhang      suffix: 2
853e662e0bSHong Zhang      args: -AtB_mat_product_algorithm at*b
868c78258cSHong Zhang      output_file: output/ex81_1.out
878c78258cSHong Zhang 
888c78258cSHong Zhang    test:
898c78258cSHong Zhang      suffix: 3
903e662e0bSHong Zhang      args: -AtB_mat_product_algorithm outerproduct
918c78258cSHong Zhang      output_file: output/ex81_1.out
928c78258cSHong Zhang 
938c78258cSHong Zhang    test:
948c78258cSHong Zhang      suffix: 4
958c78258cSHong Zhang      nsize: 3
963e662e0bSHong Zhang      args: -AtB_mat_product_algorithm nonscalable
978c78258cSHong Zhang      output_file: output/ex81_3.out
988c78258cSHong Zhang 
998c78258cSHong Zhang    test:
1008c78258cSHong Zhang      suffix: 5
1018c78258cSHong Zhang      nsize: 3
1023e662e0bSHong Zhang      args: -AtB_mat_product_algorithm scalable
1038c78258cSHong Zhang      output_file: output/ex81_3.out
1048c78258cSHong Zhang 
1058c78258cSHong Zhang    test:
1068c78258cSHong Zhang      suffix: 6
1078c78258cSHong Zhang      nsize: 3
1083e662e0bSHong Zhang      args: -AtB_mat_product_algorithm at*b
1098c78258cSHong Zhang      output_file: output/ex81_3.out
1108c78258cSHong Zhang 
1118c78258cSHong Zhang TEST*/
112