1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests for bugs in A->offloadmask consistency for GPU matrices\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6*d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 7*d71ae5a4SJacob Faibussowitsch { 8c4762a1bSJed Brown Mat A; 9c4762a1bSJed Brown PetscInt i, j, rstart, rend, m = 3; 10c4762a1bSJed Brown PetscScalar one = 1.0, zero = 0.0, negativeone = -1.0; 11c4762a1bSJed Brown PetscReal norm; 12c4762a1bSJed Brown Vec x, y; 13c4762a1bSJed Brown 14327415f7SBarry Smith PetscFunctionBeginUser; 159566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 169566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL)); 17c4762a1bSJed Brown 18c4762a1bSJed Brown for (i = 0; i < 2; i++) { 19c4762a1bSJed Brown /* Create the matrix and set it to contain explicit zero entries on the diagonal. */ 209566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 219566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, m * m, m * m)); 229566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 239566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 249566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &rstart, &rend)); 259566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &x, &y)); 269566063dSJacob Faibussowitsch PetscCall(VecSet(x, one)); 279566063dSJacob Faibussowitsch PetscCall(VecSet(y, zero)); 289566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(A, y, INSERT_VALUES)); 29c4762a1bSJed Brown 30c4762a1bSJed Brown /* Now set A to be the identity using various approaches. 31c4762a1bSJed Brown * Note that there may be other approaches that should be added here. */ 32c4762a1bSJed Brown switch (i) { 33*d71ae5a4SJacob Faibussowitsch case 0: 34*d71ae5a4SJacob Faibussowitsch PetscCall(MatDiagonalSet(A, x, INSERT_VALUES)); 35*d71ae5a4SJacob Faibussowitsch break; 36c4762a1bSJed Brown case 1: 3748a46eb9SPierre Jolivet for (j = rstart; j < rend; j++) PetscCall(MatSetValue(A, j, j, one, INSERT_VALUES)); 389566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 399566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 40c4762a1bSJed Brown break; 41c4762a1bSJed Brown case 2: 4248a46eb9SPierre Jolivet for (j = rstart; j < rend; j++) PetscCall(MatSetValuesRow(A, j, &one)); 439566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 449566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 45*d71ae5a4SJacob Faibussowitsch default: 46*d71ae5a4SJacob Faibussowitsch break; 47c4762a1bSJed Brown } 48c4762a1bSJed Brown 49c4762a1bSJed Brown /* Compute y <- A*x and verify that the difference between y and x is negligible, as it should be since A is the identity. */ 509566063dSJacob Faibussowitsch PetscCall(MatMult(A, x, y)); 519566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, negativeone, x)); 529566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &norm)); 5348a46eb9SPierre Jolivet if (norm > PETSC_SQRT_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test %" PetscInt_FMT ": Norm of error is %g, but should be near 0.\n", i, (double)norm)); 54c4762a1bSJed Brown 559566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 569566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 579566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 58c4762a1bSJed Brown } 59c4762a1bSJed Brown 609566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 61b122ec5aSJacob Faibussowitsch return 0; 62c4762a1bSJed Brown } 63c4762a1bSJed Brown 64c4762a1bSJed Brown /*TEST 65c4762a1bSJed Brown 66c4762a1bSJed Brown test: 67c4762a1bSJed Brown suffix: aijviennacl_1 68c4762a1bSJed Brown nsize: 1 69c4762a1bSJed Brown args: -mat_type aijviennacl 70c4762a1bSJed Brown requires: viennacl 71c4762a1bSJed Brown 72c4762a1bSJed Brown test: 73c4762a1bSJed Brown suffix: aijviennacl_2 74c4762a1bSJed Brown nsize: 2 75c4762a1bSJed Brown args: -mat_type aijviennacl 76c4762a1bSJed Brown requires: viennacl 77c4762a1bSJed Brown 78c4762a1bSJed Brown test: 79c4762a1bSJed Brown suffix: aijcusparse_1 80c4762a1bSJed Brown nsize: 1 81c4762a1bSJed Brown args: -mat_type aijcusparse 82c4762a1bSJed Brown requires: cuda 83c4762a1bSJed Brown 84c4762a1bSJed Brown test: 85c4762a1bSJed Brown suffix: aijcusparse_2 86c4762a1bSJed Brown nsize: 2 87c4762a1bSJed Brown args: -mat_type aijcusparse 88c4762a1bSJed Brown requires: cuda 89c4762a1bSJed Brown TEST*/ 90