114277c92SJacob Faibussowitsch static const char help[] = "Tests MatGetDiagonal().\n\n"; 214277c92SJacob Faibussowitsch 314277c92SJacob Faibussowitsch #include <petscmat.h> 414277c92SJacob Faibussowitsch 514277c92SJacob Faibussowitsch static PetscErrorCode CheckDiagonal(Mat A, Vec diag, PetscScalar dval) 614277c92SJacob Faibussowitsch { 714277c92SJacob Faibussowitsch static PetscBool first_time = PETSC_TRUE; 814277c92SJacob Faibussowitsch const PetscReal rtol = 1e-10, atol = PETSC_SMALL; 914277c92SJacob Faibussowitsch PetscInt rstart, rend, n; 1014277c92SJacob Faibussowitsch const PetscScalar *arr; 1114277c92SJacob Faibussowitsch 1214277c92SJacob Faibussowitsch PetscFunctionBegin; 1314277c92SJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &rstart, &rend)); 14baca6076SPierre Jolivet // If matrix is AIJ, MatSetRandom() will have randomly chosen the locations of nonzeros, 1514277c92SJacob Faibussowitsch // which may not be on the diagonal. So a reallocation is not necessarily a bad thing here. 1614277c92SJacob Faibussowitsch if (first_time) PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 1714277c92SJacob Faibussowitsch for (PetscInt i = rstart; i < rend; ++i) PetscCall(MatSetValue(A, i, i, dval, INSERT_VALUES)); 1814277c92SJacob Faibussowitsch if (first_time) { 1914277c92SJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE)); 2014277c92SJacob Faibussowitsch first_time = PETSC_FALSE; 2114277c92SJacob Faibussowitsch } 2214277c92SJacob Faibussowitsch 2314277c92SJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 2414277c92SJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 2514277c92SJacob Faibussowitsch PetscCall(MatViewFromOptions(A, NULL, "-mat_view_assembled")); 2614277c92SJacob Faibussowitsch 2714277c92SJacob Faibussowitsch PetscCall(VecSetRandom(diag, NULL)); 2814277c92SJacob Faibussowitsch PetscCall(MatGetDiagonal(A, diag)); 2914277c92SJacob Faibussowitsch PetscCall(VecViewFromOptions(diag, NULL, "-diag_vec_view")); 3014277c92SJacob Faibussowitsch 3114277c92SJacob Faibussowitsch PetscCall(VecGetLocalSize(diag, &n)); 3214277c92SJacob Faibussowitsch PetscCall(VecGetArrayRead(diag, &arr)); 3314277c92SJacob Faibussowitsch for (PetscInt i = 0; i < n; ++i) { 3414277c92SJacob Faibussowitsch const PetscScalar lhs = arr[i]; 3514277c92SJacob Faibussowitsch 3614277c92SJacob Faibussowitsch if (!PetscIsCloseAtTolScalar(lhs, dval, rtol, atol)) { 3714277c92SJacob Faibussowitsch const PetscReal lhs_r = PetscRealPart(lhs); 3814277c92SJacob Faibussowitsch const PetscReal lhs_i = PetscImaginaryPart(lhs); 3914277c92SJacob Faibussowitsch const PetscReal dval_r = PetscRealPart(dval); 4014277c92SJacob Faibussowitsch const PetscReal dval_i = PetscImaginaryPart(dval); 4114277c92SJacob Faibussowitsch 4214277c92SJacob Faibussowitsch PetscCheck(PetscIsCloseAtTol(lhs_r, dval_r, rtol, atol), PETSC_COMM_SELF, PETSC_ERR_PLIB, "Real component actual[%" PetscInt_FMT "] %g != expected[%" PetscInt_FMT "] %g", i, (double)lhs_r, i, (double)dval_r); 4314277c92SJacob Faibussowitsch PetscCheck(PetscIsCloseAtTol(lhs_i, dval_i, rtol, atol), PETSC_COMM_SELF, PETSC_ERR_PLIB, "Imaginary component actual[%" PetscInt_FMT "] %g != expected[%" PetscInt_FMT "] %g", i, (double)lhs_i, i, (double)dval_i); 4414277c92SJacob Faibussowitsch } 4514277c92SJacob Faibussowitsch } 4614277c92SJacob Faibussowitsch PetscCall(VecRestoreArrayRead(diag, &arr)); 4714277c92SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4814277c92SJacob Faibussowitsch } 4914277c92SJacob Faibussowitsch 5014277c92SJacob Faibussowitsch static PetscErrorCode InitializeMatrix(Mat A) 5114277c92SJacob Faibussowitsch { 5214277c92SJacob Faibussowitsch MatType mtype; 5314277c92SJacob Faibussowitsch char *tmp; 5414277c92SJacob Faibussowitsch 5514277c92SJacob Faibussowitsch PetscFunctionBegin; 5614277c92SJacob Faibussowitsch PetscCall(MatSetUp(A)); 5714277c92SJacob Faibussowitsch PetscCall(MatGetType(A, &mtype)); 5814277c92SJacob Faibussowitsch PetscCall(PetscStrstr(mtype, "aij", &tmp)); 5914277c92SJacob Faibussowitsch if (tmp) { 6014277c92SJacob Faibussowitsch PetscInt rows, cols, diag_nnz, offdiag_nnz; 6114277c92SJacob Faibussowitsch PetscInt *dnnz, *onnz; 6214277c92SJacob Faibussowitsch 6314277c92SJacob Faibussowitsch // is some kind of AIJ 6414277c92SJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &rows, &cols)); 6514277c92SJacob Faibussowitsch // at least 3 nonzeros in diagonal block 6614277c92SJacob Faibussowitsch diag_nnz = PetscMin(cols, 3); 6714277c92SJacob Faibussowitsch // leave at least 3 *zeros* per row 6814277c92SJacob Faibussowitsch offdiag_nnz = PetscMax(cols - diag_nnz - 3, 0); 6914277c92SJacob Faibussowitsch PetscCall(PetscMalloc2(rows, &dnnz, rows, &onnz)); 7014277c92SJacob Faibussowitsch for (PetscInt i = 0; i < rows; ++i) { 7114277c92SJacob Faibussowitsch dnnz[i] = diag_nnz; 7214277c92SJacob Faibussowitsch onnz[i] = offdiag_nnz; 7314277c92SJacob Faibussowitsch } 7414277c92SJacob Faibussowitsch PetscCall(MatXAIJSetPreallocation(A, PETSC_DECIDE, dnnz, onnz, NULL, NULL)); 7514277c92SJacob Faibussowitsch PetscCall(PetscFree2(dnnz, onnz)); 7614277c92SJacob Faibussowitsch } 7714277c92SJacob Faibussowitsch 7814277c92SJacob Faibussowitsch PetscCall(MatSetRandom(A, NULL)); 7914277c92SJacob Faibussowitsch PetscCall(MatViewFromOptions(A, NULL, "-mat_view_setup")); 8014277c92SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8114277c92SJacob Faibussowitsch } 8214277c92SJacob Faibussowitsch 8314277c92SJacob Faibussowitsch int main(int argc, char **argv) 8414277c92SJacob Faibussowitsch { 8514277c92SJacob Faibussowitsch Mat A; 8614277c92SJacob Faibussowitsch Vec diag; 8714277c92SJacob Faibussowitsch 8814277c92SJacob Faibussowitsch PetscFunctionBeginUser; 8914277c92SJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 9014277c92SJacob Faibussowitsch 9114277c92SJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 9214277c92SJacob Faibussowitsch PetscCall(MatSetSizes(A, 10, 10, PETSC_DECIDE, PETSC_DECIDE)); 9314277c92SJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 9414277c92SJacob Faibussowitsch 9514277c92SJacob Faibussowitsch PetscCall(InitializeMatrix(A)); 9614277c92SJacob Faibussowitsch 9714277c92SJacob Faibussowitsch PetscCall(MatCreateVecs(A, &diag, NULL)); 9814277c92SJacob Faibussowitsch 9914277c92SJacob Faibussowitsch PetscCall(CheckDiagonal(A, diag, 0.0)); 10014277c92SJacob Faibussowitsch PetscCall(CheckDiagonal(A, diag, 1.0)); 10114277c92SJacob Faibussowitsch PetscCall(CheckDiagonal(A, diag, 2.0)); 10214277c92SJacob Faibussowitsch 10314277c92SJacob Faibussowitsch PetscCall(VecDestroy(&diag)); 10414277c92SJacob Faibussowitsch PetscCall(MatDestroy(&A)); 10514277c92SJacob Faibussowitsch PetscCall(PetscFinalize()); 10614277c92SJacob Faibussowitsch return 0; 10714277c92SJacob Faibussowitsch } 10814277c92SJacob Faibussowitsch 10914277c92SJacob Faibussowitsch /*TEST 11014277c92SJacob Faibussowitsch 11114277c92SJacob Faibussowitsch testset: 11214277c92SJacob Faibussowitsch output_file: ./output/empty.out 11314277c92SJacob Faibussowitsch nsize: {{1 2}} 11414277c92SJacob Faibussowitsch suffix: dense 11514277c92SJacob Faibussowitsch test: 11614277c92SJacob Faibussowitsch suffix: standard 11714277c92SJacob Faibussowitsch args: -mat_type dense 11814277c92SJacob Faibussowitsch test: 11914277c92SJacob Faibussowitsch suffix: cuda 12014277c92SJacob Faibussowitsch requires: cuda 12114277c92SJacob Faibussowitsch args: -mat_type densecuda 12214277c92SJacob Faibussowitsch test: 12314277c92SJacob Faibussowitsch suffix: hip 12414277c92SJacob Faibussowitsch requires: hip 12514277c92SJacob Faibussowitsch args: -mat_type densehip 12614277c92SJacob Faibussowitsch 12714277c92SJacob Faibussowitsch testset: 12814277c92SJacob Faibussowitsch output_file: ./output/empty.out 12914277c92SJacob Faibussowitsch nsize: {{1 2}} 13014277c92SJacob Faibussowitsch suffix: aij 13114277c92SJacob Faibussowitsch test: 13214277c92SJacob Faibussowitsch suffix: standard 13314277c92SJacob Faibussowitsch args: -mat_type aij 13414277c92SJacob Faibussowitsch test: 13514277c92SJacob Faibussowitsch suffix: cuda 13614277c92SJacob Faibussowitsch requires: cuda 13714277c92SJacob Faibussowitsch args: -mat_type aijcusparse 13814277c92SJacob Faibussowitsch test: 13914277c92SJacob Faibussowitsch suffix: hip 14014277c92SJacob Faibussowitsch requires: hip 14114277c92SJacob Faibussowitsch args: -mat_type aijhipsparse 14214277c92SJacob Faibussowitsch test: 14314277c92SJacob Faibussowitsch suffix: kokkos 144*4c073473SPierre Jolivet requires: kokkos kokkos_kernels 14514277c92SJacob Faibussowitsch args: -mat_type aijkokkos 14614277c92SJacob Faibussowitsch 14714277c92SJacob Faibussowitsch TEST*/ 148