Actual source code: ex37.c

  2: static char help[] = "Tests MatCopy() and MatStore/RetrieveValues().\n\n";

 4:  #include petscmat.h

  8: int main(int argc,char **args)
  9: {
 10:   Mat            C,A;
 11:   PetscInt       i, n = 10,midx[3];
 13:   PetscScalar    v[3];
 14:   PetscTruth     flg;

 16:   PetscInitialize(&argc,&args,(char *)0,help);
 17:   PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);

 19:   MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,&C);
 20:   MatSetFromOptions(C);

 22:   v[0] = -1.; v[1] = 2.; v[2] = -1.;
 23:   for (i=1; i<n-1; i++){
 24:     midx[2] = i-1; midx[1] = i; midx[0] = i+1;
 25:     MatSetValues(C,1,&i,3,midx,v,INSERT_VALUES);
 26:   }
 27:   i = 0; midx[0] = 0; midx[1] = 1;
 28:   v[0] = 2.0; v[1] = -1.;
 29:   MatSetValues(C,1,&i,2,midx,v,INSERT_VALUES);
 30:   i = n-1; midx[0] = n-2; midx[1] = n-1;
 31:   v[0] = -1.0; v[1] = 2.;
 32:   MatSetValues(C,1,&i,2,midx,v,INSERT_VALUES);

 34:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 35:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

 37:   MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,&A);
 38:   MatSetFromOptions(A);
 39:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 40:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 42:   /* test matrices with different nonzero patterns */
 43:   MatCopy(C,A,DIFFERENT_NONZERO_PATTERN);

 45:   /* Now C and A have the same nonzero pattern */
 46:   MatSetOption(C,MAT_NO_NEW_NONZERO_LOCATIONS);
 47:   MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);
 48:   MatCopy(C,A,SAME_NONZERO_PATTERN);

 50:   MatView(C,PETSC_VIEWER_STDOUT_WORLD);
 51:   MatView(A,PETSC_VIEWER_STDOUT_WORLD);

 53:   MatEqual(A,C,&flg);
 54:   if (flg) {
 55:     PetscPrintf(PETSC_COMM_WORLD,"Matrices are equal\n");
 56:   } else {
 57:     SETERRQ(1,"Matrices are NOT equal");
 58:   }

 60:   MatStoreValues(A);
 61:   MatZeroEntries(A);
 62:   MatRetrieveValues(A);
 63:   MatEqual(A,C,&flg);
 64:   if (flg) {
 65:     PetscPrintf(PETSC_COMM_WORLD,"Matrices are equal\n");
 66:   } else {
 67:     SETERRQ(1,"Matrices are NOT equal");
 68:   }

 70:   MatDestroy(C);
 71:   MatDestroy(A);

 73:   PetscFinalize();
 74:   return 0;
 75: }

 77: