1c4762a1bSJed Brown static char help[] = "Tests DMSLICED operations\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscdmsliced.h> 4c4762a1bSJed Brown 5*9371c9d4SSatish Balay int main(int argc, char *argv[]) { 6c4762a1bSJed Brown char mat_type[256] = MATAIJ; /* default matrix type */ 7c4762a1bSJed Brown MPI_Comm comm; 8c4762a1bSJed Brown PetscMPIInt rank, size; 9c4762a1bSJed Brown DM slice; 10c4762a1bSJed Brown PetscInt i, bs = 1, N = 5, n, m, rstart, ghosts[2], *d_nnz, *o_nnz, dfill[4] = {1, 0, 0, 1}, ofill[4] = {1, 1, 1, 1}; 11c4762a1bSJed Brown PetscReal alpha = 1, K = 1, rho0 = 1, u0 = 0, sigma = 0.2; 12c4762a1bSJed Brown PetscBool useblock = PETSC_TRUE; 13c4762a1bSJed Brown PetscScalar *xx; 14c4762a1bSJed Brown Mat A; 15c4762a1bSJed Brown Vec x, b, lf; 16c4762a1bSJed Brown 17327415f7SBarry Smith PetscFunctionBeginUser; 189566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, 0, help)); 19c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 209566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 22c4762a1bSJed Brown 23d0609cedSBarry Smith PetscOptionsBegin(comm, 0, "Options for DMSliced test", 0); 24c4762a1bSJed Brown { 259566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-n", "Global number of nodes", "", N, &N, NULL)); 269566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-bs", "Block size (1 or 2)", "", bs, &bs, NULL)); 27c4762a1bSJed Brown if (bs != 1) { 2808401ef6SPierre Jolivet PetscCheck(bs == 2, PETSC_COMM_WORLD, PETSC_ERR_SUP, "Block size must be 1 or 2"); 299566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-alpha", "Inverse time step for wave operator", "", alpha, &alpha, NULL)); 309566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-K", "Bulk modulus of compressibility", "", K, &K, NULL)); 319566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-rho0", "Reference density", "", rho0, &rho0, NULL)); 329566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-u0", "Reference velocity", "", u0, &u0, NULL)); 339566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-sigma", "Width of Gaussian density perturbation", "", sigma, &sigma, NULL)); 349566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-block", "Use block matrix assembly", "", useblock, &useblock, NULL)); 35c4762a1bSJed Brown } 369566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-sliced_mat_type", "Matrix type to use (aij or baij)", "", mat_type, mat_type, sizeof(mat_type), NULL)); 37c4762a1bSJed Brown } 38d0609cedSBarry Smith PetscOptionsEnd(); 39c4762a1bSJed Brown 40c4762a1bSJed Brown /* Split ownership, set up periodic grid in 1D */ 41c4762a1bSJed Brown n = PETSC_DECIDE; 429566063dSJacob Faibussowitsch PetscCall(PetscSplitOwnership(comm, &n, &N)); 43c4762a1bSJed Brown rstart = 0; 449566063dSJacob Faibussowitsch PetscCallMPI(MPI_Scan(&n, &rstart, 1, MPIU_INT, MPI_SUM, comm)); 45c4762a1bSJed Brown rstart -= n; 46c4762a1bSJed Brown ghosts[0] = (N + rstart - 1) % N; 47c4762a1bSJed Brown ghosts[1] = (rstart + n) % N; 48c4762a1bSJed Brown 499566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(n, &d_nnz, n, &o_nnz)); 50c4762a1bSJed Brown for (i = 0; i < n; i++) { 51c4762a1bSJed Brown if (size > 1 && (i == 0 || i == n - 1)) { 52c4762a1bSJed Brown d_nnz[i] = 2; 53c4762a1bSJed Brown o_nnz[i] = 1; 54c4762a1bSJed Brown } else { 55c4762a1bSJed Brown d_nnz[i] = 3; 56c4762a1bSJed Brown o_nnz[i] = 0; 57c4762a1bSJed Brown } 58c4762a1bSJed Brown } 599566063dSJacob Faibussowitsch PetscCall(DMSlicedCreate(comm, bs, n, 2, ghosts, d_nnz, o_nnz, &slice)); /* Currently does not copy X_nnz so we can't free them until after DMSlicedGetMatrix */ 60c4762a1bSJed Brown 619566063dSJacob Faibussowitsch if (!useblock) PetscCall(DMSlicedSetBlockFills(slice, dfill, ofill)); /* Irrelevant for baij formats */ 629566063dSJacob Faibussowitsch PetscCall(DMSetMatType(slice, mat_type)); 639566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(slice, &A)); 649566063dSJacob Faibussowitsch PetscCall(PetscFree2(d_nnz, o_nnz)); 659566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE)); 66c4762a1bSJed Brown 679566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(slice, &x)); 689566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &b)); 69c4762a1bSJed Brown 709566063dSJacob Faibussowitsch PetscCall(VecGhostGetLocalForm(x, &lf)); 719566063dSJacob Faibussowitsch PetscCall(VecGetSize(lf, &m)); 7263a3b9bcSJacob Faibussowitsch PetscCheck(m == (n + 2) * bs, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "size of local form %" PetscInt_FMT ", expected %" PetscInt_FMT, m, (n + 2) * bs); 739566063dSJacob Faibussowitsch PetscCall(VecGetArray(lf, &xx)); 74c4762a1bSJed Brown for (i = 0; i < n; i++) { 75c4762a1bSJed Brown PetscInt row[2], col[9], im, ip; 76c4762a1bSJed Brown PetscScalar v[12]; 77c4762a1bSJed Brown const PetscReal xref = 2.0 * (rstart + i) / N - 1; /* [-1,1] */ 78c4762a1bSJed Brown const PetscReal h = 1.0 / N; /* grid spacing */ 79c4762a1bSJed Brown im = (i == 0) ? n : i - 1; 80c4762a1bSJed Brown ip = (i == n - 1) ? n + 1 : i + 1; 81c4762a1bSJed Brown switch (bs) { 82c4762a1bSJed Brown case 1: /* Laplacian with periodic boundaries */ 83*9371c9d4SSatish Balay col[0] = im; 84*9371c9d4SSatish Balay col[1] = i; 85*9371c9d4SSatish Balay col[2] = ip; 86*9371c9d4SSatish Balay v[0] = -h; 87*9371c9d4SSatish Balay v[1] = 2 * h; 88*9371c9d4SSatish Balay v[2] = -h; 899566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(A, 1, &i, 3, col, v, INSERT_VALUES)); 90c4762a1bSJed Brown xx[i] = PetscSinReal(xref * PETSC_PI); 91c4762a1bSJed Brown break; 92c4762a1bSJed Brown case 2: /* Linear acoustic wave operator in variables [rho, u], central differences, periodic, timestep 1/alpha */ 93*9371c9d4SSatish Balay v[0] = -0.5 * u0; 94*9371c9d4SSatish Balay v[1] = -0.5 * K; 95*9371c9d4SSatish Balay v[2] = alpha; 96*9371c9d4SSatish Balay v[3] = 0; 97*9371c9d4SSatish Balay v[4] = 0.5 * u0; 98*9371c9d4SSatish Balay v[5] = 0.5 * K; 99*9371c9d4SSatish Balay v[6] = -0.5 / rho0; 100*9371c9d4SSatish Balay v[7] = -0.5 * u0; 101*9371c9d4SSatish Balay v[8] = 0; 102*9371c9d4SSatish Balay v[9] = alpha; 103*9371c9d4SSatish Balay v[10] = 0.5 / rho0; 104*9371c9d4SSatish Balay v[11] = 0.5 * u0; 105c4762a1bSJed Brown if (useblock) { 106*9371c9d4SSatish Balay row[0] = i; 107*9371c9d4SSatish Balay col[0] = im; 108*9371c9d4SSatish Balay col[1] = i; 109*9371c9d4SSatish Balay col[2] = ip; 1109566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlockedLocal(A, 1, row, 3, col, v, INSERT_VALUES)); 111c4762a1bSJed Brown } else { 112*9371c9d4SSatish Balay row[0] = 2 * i; 113*9371c9d4SSatish Balay row[1] = 2 * i + 1; 114*9371c9d4SSatish Balay col[0] = 2 * im; 115*9371c9d4SSatish Balay col[1] = 2 * im + 1; 116*9371c9d4SSatish Balay col[2] = 2 * i; 117*9371c9d4SSatish Balay col[3] = 2 * ip; 118*9371c9d4SSatish Balay col[4] = 2 * ip + 1; 119*9371c9d4SSatish Balay v[3] = v[4]; 120*9371c9d4SSatish Balay v[4] = v[5]; /* pack values in first row */ 1219566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(A, 1, row, 5, col, v, INSERT_VALUES)); 122c4762a1bSJed Brown col[2] = 2 * i + 1; 123*9371c9d4SSatish Balay v[8] = v[9]; 124*9371c9d4SSatish Balay v[9] = v[10]; 125*9371c9d4SSatish Balay v[10] = v[11]; /* pack values in second row */ 1269566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(A, 1, row + 1, 5, col, v + 6, INSERT_VALUES)); 127c4762a1bSJed Brown } 128c4762a1bSJed Brown /* Set current state (gaussian density perturbation) */ 129c4762a1bSJed Brown xx[2 * i] = 0.2 * PetscExpReal(-PetscSqr(xref) / (2 * PetscSqr(sigma))); 130c4762a1bSJed Brown xx[2 * i + 1] = 0; 131c4762a1bSJed Brown break; 13263a3b9bcSJacob Faibussowitsch default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "not implemented for block size %" PetscInt_FMT, bs); 133c4762a1bSJed Brown } 134c4762a1bSJed Brown } 1359566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(lf, &xx)); 1369566063dSJacob Faibussowitsch PetscCall(VecGhostRestoreLocalForm(x, &lf)); 1379566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1389566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 139c4762a1bSJed Brown 1409566063dSJacob Faibussowitsch PetscCall(MatMult(A, x, b)); 1419566063dSJacob Faibussowitsch PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); 1429566063dSJacob Faibussowitsch PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD)); 1439566063dSJacob Faibussowitsch PetscCall(VecView(b, PETSC_VIEWER_STDOUT_WORLD)); 144c4762a1bSJed Brown 145c4762a1bSJed Brown /* Update the ghosted values, view the result on rank 0. */ 1469566063dSJacob Faibussowitsch PetscCall(VecGhostUpdateBegin(b, INSERT_VALUES, SCATTER_FORWARD)); 1479566063dSJacob Faibussowitsch PetscCall(VecGhostUpdateEnd(b, INSERT_VALUES, SCATTER_FORWARD)); 148dd400576SPatrick Sanan if (rank == 0) { 1499566063dSJacob Faibussowitsch PetscCall(VecGhostGetLocalForm(b, &lf)); 1509566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF, "Local form of b on rank 0, last two nodes are ghost nodes\n")); 1519566063dSJacob Faibussowitsch PetscCall(VecView(lf, PETSC_VIEWER_STDOUT_SELF)); 1529566063dSJacob Faibussowitsch PetscCall(VecGhostRestoreLocalForm(b, &lf)); 153c4762a1bSJed Brown } 154c4762a1bSJed Brown 1559566063dSJacob Faibussowitsch PetscCall(DMDestroy(&slice)); 1569566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1579566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 1589566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1599566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 160b122ec5aSJacob Faibussowitsch return 0; 161c4762a1bSJed Brown } 162c4762a1bSJed Brown 163c4762a1bSJed Brown /*TEST 164c4762a1bSJed Brown 165c4762a1bSJed Brown test: 166c4762a1bSJed Brown nsize: 2 167c4762a1bSJed Brown args: -bs 2 -block 0 -sliced_mat_type baij -alpha 10 -u0 0.1 168c4762a1bSJed Brown 169c4762a1bSJed Brown test: 170c4762a1bSJed Brown suffix: 2 171c4762a1bSJed Brown nsize: 2 172c4762a1bSJed Brown args: -bs 2 -block 1 -sliced_mat_type aij -alpha 10 -u0 0.1 173c4762a1bSJed Brown 174c4762a1bSJed Brown test: 175c4762a1bSJed Brown suffix: 3 176c4762a1bSJed Brown nsize: 2 177c4762a1bSJed Brown args: -bs 2 -block 0 -sliced_mat_type aij -alpha 10 -u0 0.1 178c4762a1bSJed Brown 179c4762a1bSJed Brown TEST*/ 180