1c4762a1bSJed Brown 2c4762a1bSJed Brown /* 3c4762a1bSJed Brown Laplacian in 3D. Use for testing MatSolve routines. 4c4762a1bSJed Brown Modeled by the partial differential equation 5c4762a1bSJed Brown 6c4762a1bSJed Brown - Laplacian u = 1,0 < x,y,z < 1, 7c4762a1bSJed Brown 8c4762a1bSJed Brown with boundary conditions 9c4762a1bSJed Brown u = 1 for x = 0, x = 1, y = 0, y = 1, z = 0, z = 1. 10c4762a1bSJed Brown */ 11c4762a1bSJed Brown 12c4762a1bSJed Brown static char help[] = "This example is for testing different MatSolve routines :MatSolve(), MatSolveAdd(), MatSolveTranspose(), MatSolveTransposeAdd(), and MatMatSolve().\n\ 13c4762a1bSJed Brown Example usage: ./ex129 -mat_type aij -dof 2\n\n"; 14c4762a1bSJed Brown 15c4762a1bSJed Brown #include <petscdm.h> 16c4762a1bSJed Brown #include <petscdmda.h> 17c4762a1bSJed Brown 18c4762a1bSJed Brown extern PetscErrorCode ComputeMatrix(DM, Mat); 19c4762a1bSJed Brown extern PetscErrorCode ComputeRHS(DM, Vec); 20c4762a1bSJed Brown extern PetscErrorCode ComputeRHSMatrix(PetscInt, PetscInt, Mat *); 21c4762a1bSJed Brown 22d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 23d71ae5a4SJacob Faibussowitsch { 24c4762a1bSJed Brown PetscMPIInt size; 25c4762a1bSJed Brown Vec x, b, y, b1; 26c4762a1bSJed Brown DM da; 27c4762a1bSJed Brown Mat A, F, RHS, X, C1; 28c4762a1bSJed Brown MatFactorInfo info; 29c4762a1bSJed Brown IS perm, iperm; 30c4762a1bSJed Brown PetscInt dof = 1, M = 8, m, n, nrhs; 31c4762a1bSJed Brown PetscScalar one = 1.0; 32c4762a1bSJed Brown PetscReal norm, tol = 1000 * PETSC_MACHINE_EPSILON; 33c4762a1bSJed Brown PetscBool InplaceLU = PETSC_FALSE; 34c4762a1bSJed Brown 35327415f7SBarry Smith PetscFunctionBeginUser; 369566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 379566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 38be096a46SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only"); 399566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-dof", &dof, NULL)); 409566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-M", &M, NULL)); 41c4762a1bSJed Brown 429566063dSJacob Faibussowitsch PetscCall(DMDACreate(PETSC_COMM_WORLD, &da)); 439566063dSJacob Faibussowitsch PetscCall(DMSetDimension(da, 3)); 449566063dSJacob Faibussowitsch PetscCall(DMDASetBoundaryType(da, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE)); 459566063dSJacob Faibussowitsch PetscCall(DMDASetStencilType(da, DMDA_STENCIL_STAR)); 469566063dSJacob Faibussowitsch PetscCall(DMDASetSizes(da, M, M, M)); 479566063dSJacob Faibussowitsch PetscCall(DMDASetNumProcs(da, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE)); 489566063dSJacob Faibussowitsch PetscCall(DMDASetDof(da, dof)); 499566063dSJacob Faibussowitsch PetscCall(DMDASetStencilWidth(da, 1)); 509566063dSJacob Faibussowitsch PetscCall(DMDASetOwnershipRanges(da, NULL, NULL, NULL)); 519566063dSJacob Faibussowitsch PetscCall(DMSetMatType(da, MATBAIJ)); 529566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 539566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 54c4762a1bSJed Brown 559566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &x)); 569566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &b)); 579566063dSJacob Faibussowitsch PetscCall(VecDuplicate(b, &y)); 589566063dSJacob Faibussowitsch PetscCall(ComputeRHS(da, b)); 599566063dSJacob Faibussowitsch PetscCall(VecSet(y, one)); 609566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(da, &A)); 619566063dSJacob Faibussowitsch PetscCall(ComputeMatrix(da, A)); 629566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &m, &n)); 63c4762a1bSJed Brown nrhs = 2; 649566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-nrhs", &nrhs, NULL)); 659566063dSJacob Faibussowitsch PetscCall(ComputeRHSMatrix(m, nrhs, &RHS)); 669566063dSJacob Faibussowitsch PetscCall(MatDuplicate(RHS, MAT_DO_NOT_COPY_VALUES, &X)); 67c4762a1bSJed Brown 689566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(A, MATORDERINGND, &perm, &iperm)); 69c4762a1bSJed Brown 709566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-inplacelu", &InplaceLU, NULL)); 719566063dSJacob Faibussowitsch PetscCall(MatFactorInfoInitialize(&info)); 72c4762a1bSJed Brown if (!InplaceLU) { 739566063dSJacob Faibussowitsch PetscCall(MatGetFactor(A, MATSOLVERPETSC, MAT_FACTOR_LU, &F)); 74c4762a1bSJed Brown info.fill = 5.0; 759566063dSJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(F, A, perm, iperm, &info)); 769566063dSJacob Faibussowitsch PetscCall(MatLUFactorNumeric(F, A, &info)); 77c4762a1bSJed Brown } else { /* Test inplace factorization */ 789566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &F)); 799566063dSJacob Faibussowitsch PetscCall(MatLUFactor(F, perm, iperm, &info)); 80c4762a1bSJed Brown } 81c4762a1bSJed Brown 829566063dSJacob Faibussowitsch PetscCall(VecDuplicate(y, &b1)); 83c4762a1bSJed Brown 84c4762a1bSJed Brown /* MatSolve */ 859566063dSJacob Faibussowitsch PetscCall(MatSolve(F, b, x)); 869566063dSJacob Faibussowitsch PetscCall(MatMult(A, x, b1)); 879566063dSJacob Faibussowitsch PetscCall(VecAXPY(b1, -1.0, b)); 889566063dSJacob Faibussowitsch PetscCall(VecNorm(b1, NORM_2, &norm)); 8948a46eb9SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "MatSolve : Error of norm %g\n", (double)norm)); 90c4762a1bSJed Brown 91c4762a1bSJed Brown /* MatSolveTranspose */ 929566063dSJacob Faibussowitsch PetscCall(MatSolveTranspose(F, b, x)); 939566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, x, b1)); 949566063dSJacob Faibussowitsch PetscCall(VecAXPY(b1, -1.0, b)); 959566063dSJacob Faibussowitsch PetscCall(VecNorm(b1, NORM_2, &norm)); 9648a46eb9SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "MatSolveTranspose : Error of norm %g\n", (double)norm)); 97c4762a1bSJed Brown 98c4762a1bSJed Brown /* MatSolveAdd */ 999566063dSJacob Faibussowitsch PetscCall(MatSolveAdd(F, b, y, x)); 1009566063dSJacob Faibussowitsch PetscCall(MatMult(A, y, b1)); 1019566063dSJacob Faibussowitsch PetscCall(VecScale(b1, -1.0)); 1029566063dSJacob Faibussowitsch PetscCall(MatMultAdd(A, x, b1, b1)); 1039566063dSJacob Faibussowitsch PetscCall(VecAXPY(b1, -1.0, b)); 1049566063dSJacob Faibussowitsch PetscCall(VecNorm(b1, NORM_2, &norm)); 10548a46eb9SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "MatSolveAdd : Error of norm %g\n", (double)norm)); 106c4762a1bSJed Brown 107c4762a1bSJed Brown /* MatSolveTransposeAdd */ 1089566063dSJacob Faibussowitsch PetscCall(MatSolveTransposeAdd(F, b, y, x)); 1099566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, y, b1)); 1109566063dSJacob Faibussowitsch PetscCall(VecScale(b1, -1.0)); 1119566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(A, x, b1, b1)); 1129566063dSJacob Faibussowitsch PetscCall(VecAXPY(b1, -1.0, b)); 1139566063dSJacob Faibussowitsch PetscCall(VecNorm(b1, NORM_2, &norm)); 11448a46eb9SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "MatSolveTransposeAdd : Error of norm %g\n", (double)norm)); 115c4762a1bSJed Brown 116c4762a1bSJed Brown /* MatMatSolve */ 1179566063dSJacob Faibussowitsch PetscCall(MatMatSolve(F, RHS, X)); 1189566063dSJacob Faibussowitsch PetscCall(MatMatMult(A, X, MAT_INITIAL_MATRIX, 2.0, &C1)); 1199566063dSJacob Faibussowitsch PetscCall(MatAXPY(C1, -1.0, RHS, SAME_NONZERO_PATTERN)); 1209566063dSJacob Faibussowitsch PetscCall(MatNorm(C1, NORM_FROBENIUS, &norm)); 12148a46eb9SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "MatMatSolve : Error of norm %g\n", (double)norm)); 122c4762a1bSJed Brown 1239566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1249566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 1259566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b1)); 1269566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 1279566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1289566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F)); 1299566063dSJacob Faibussowitsch PetscCall(MatDestroy(&RHS)); 1309566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C1)); 1319566063dSJacob Faibussowitsch PetscCall(MatDestroy(&X)); 1329566063dSJacob Faibussowitsch PetscCall(ISDestroy(&perm)); 1339566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iperm)); 1349566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 1359566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 136b122ec5aSJacob Faibussowitsch return 0; 137c4762a1bSJed Brown } 138c4762a1bSJed Brown 139d71ae5a4SJacob Faibussowitsch PetscErrorCode ComputeRHS(DM da, Vec b) 140d71ae5a4SJacob Faibussowitsch { 141c4762a1bSJed Brown PetscInt mx, my, mz; 142c4762a1bSJed Brown PetscScalar h; 143c4762a1bSJed Brown 144c4762a1bSJed Brown PetscFunctionBegin; 1459566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &mx, &my, &mz, 0, 0, 0, 0, 0, 0, 0, 0, 0)); 146c4762a1bSJed Brown h = 1.0 / ((mx - 1) * (my - 1) * (mz - 1)); 1479566063dSJacob Faibussowitsch PetscCall(VecSet(b, h)); 148*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 149c4762a1bSJed Brown } 150c4762a1bSJed Brown 151d71ae5a4SJacob Faibussowitsch PetscErrorCode ComputeRHSMatrix(PetscInt m, PetscInt nrhs, Mat *C) 152d71ae5a4SJacob Faibussowitsch { 153c4762a1bSJed Brown PetscRandom rand; 154c4762a1bSJed Brown Mat RHS; 155c4762a1bSJed Brown PetscScalar *array, rval; 156c4762a1bSJed Brown PetscInt i, k; 157c4762a1bSJed Brown 158c4762a1bSJed Brown PetscFunctionBegin; 1599566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &RHS)); 1609566063dSJacob Faibussowitsch PetscCall(MatSetSizes(RHS, m, PETSC_DECIDE, PETSC_DECIDE, nrhs)); 1619566063dSJacob Faibussowitsch PetscCall(MatSetType(RHS, MATSEQDENSE)); 1629566063dSJacob Faibussowitsch PetscCall(MatSetUp(RHS)); 163c4762a1bSJed Brown 1649566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rand)); 1659566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rand)); 1669566063dSJacob Faibussowitsch PetscCall(MatDenseGetArray(RHS, &array)); 167c4762a1bSJed Brown for (i = 0; i < m; i++) { 1689566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand, &rval)); 169c4762a1bSJed Brown array[i] = rval; 170c4762a1bSJed Brown } 171c4762a1bSJed Brown if (nrhs > 1) { 172c4762a1bSJed Brown for (k = 1; k < nrhs; k++) { 173ad540459SPierre Jolivet for (i = 0; i < m; i++) array[m * k + i] = array[i]; 174c4762a1bSJed Brown } 175c4762a1bSJed Brown } 1769566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArray(RHS, &array)); 1779566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(RHS, MAT_FINAL_ASSEMBLY)); 1789566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(RHS, MAT_FINAL_ASSEMBLY)); 179c4762a1bSJed Brown *C = RHS; 1809566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rand)); 181*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 182c4762a1bSJed Brown } 183c4762a1bSJed Brown 184d71ae5a4SJacob Faibussowitsch PetscErrorCode ComputeMatrix(DM da, Mat B) 185d71ae5a4SJacob Faibussowitsch { 186c4762a1bSJed Brown PetscInt i, j, k, mx, my, mz, xm, ym, zm, xs, ys, zs, dof, k1, k2, k3; 187c4762a1bSJed Brown PetscScalar *v, *v_neighbor, Hx, Hy, Hz, HxHydHz, HyHzdHx, HxHzdHy, r1, r2; 188c4762a1bSJed Brown MatStencil row, col; 189c4762a1bSJed Brown PetscRandom rand; 190c4762a1bSJed Brown 191c4762a1bSJed Brown PetscFunctionBegin; 1929566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rand)); 1939566063dSJacob Faibussowitsch PetscCall(PetscRandomSetSeed(rand, 1)); 1949566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rand, -.001, .001)); 1959566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rand)); 196c4762a1bSJed Brown 1979566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &mx, &my, &mz, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); 198c4762a1bSJed Brown /* For simplicity, this example only works on mx=my=mz */ 199e00437b9SBarry Smith PetscCheck(mx == my && mx == mz, PETSC_COMM_SELF, PETSC_ERR_SUP, "This example only works with mx %" PetscInt_FMT " = my %" PetscInt_FMT " = mz %" PetscInt_FMT, mx, my, mz); 200c4762a1bSJed Brown 2019371c9d4SSatish Balay Hx = 1.0 / (PetscReal)(mx - 1); 2029371c9d4SSatish Balay Hy = 1.0 / (PetscReal)(my - 1); 2039371c9d4SSatish Balay Hz = 1.0 / (PetscReal)(mz - 1); 2049371c9d4SSatish Balay HxHydHz = Hx * Hy / Hz; 2059371c9d4SSatish Balay HxHzdHy = Hx * Hz / Hy; 2069371c9d4SSatish Balay HyHzdHx = Hy * Hz / Hx; 207c4762a1bSJed Brown 2089566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2 * dof * dof + 1, &v)); 209c4762a1bSJed Brown v_neighbor = v + dof * dof; 2109566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(v, 2 * dof * dof + 1)); 211c4762a1bSJed Brown k3 = 0; 212c4762a1bSJed Brown for (k1 = 0; k1 < dof; k1++) { 213c4762a1bSJed Brown for (k2 = 0; k2 < dof; k2++) { 214c4762a1bSJed Brown if (k1 == k2) { 215c4762a1bSJed Brown v[k3] = 2.0 * (HxHydHz + HxHzdHy + HyHzdHx); 216c4762a1bSJed Brown v_neighbor[k3] = -HxHydHz; 217c4762a1bSJed Brown } else { 2189566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand, &r1)); 2199566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand, &r2)); 220c4762a1bSJed Brown 221c4762a1bSJed Brown v[k3] = r1; 222c4762a1bSJed Brown v_neighbor[k3] = r2; 223c4762a1bSJed Brown } 224c4762a1bSJed Brown k3++; 225c4762a1bSJed Brown } 226c4762a1bSJed Brown } 2279566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm)); 228c4762a1bSJed Brown 229c4762a1bSJed Brown for (k = zs; k < zs + zm; k++) { 230c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 231c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 2329371c9d4SSatish Balay row.i = i; 2339371c9d4SSatish Balay row.j = j; 2349371c9d4SSatish Balay row.k = k; 235a5b23f4aSJose E. Roman if (i == 0 || j == 0 || k == 0 || i == mx - 1 || j == my - 1 || k == mz - 1) { /* boundary points */ 2369566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlockedStencil(B, 1, &row, 1, &row, v, INSERT_VALUES)); 237c4762a1bSJed Brown } else { /* interior points */ 238c4762a1bSJed Brown /* center */ 2399371c9d4SSatish Balay col.i = i; 2409371c9d4SSatish Balay col.j = j; 2419371c9d4SSatish Balay col.k = k; 2429566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlockedStencil(B, 1, &row, 1, &col, v, INSERT_VALUES)); 243c4762a1bSJed Brown 244c4762a1bSJed Brown /* x neighbors */ 2459371c9d4SSatish Balay col.i = i - 1; 2469371c9d4SSatish Balay col.j = j; 2479371c9d4SSatish Balay col.k = k; 2489566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlockedStencil(B, 1, &row, 1, &col, v_neighbor, INSERT_VALUES)); 2499371c9d4SSatish Balay col.i = i + 1; 2509371c9d4SSatish Balay col.j = j; 2519371c9d4SSatish Balay col.k = k; 2529566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlockedStencil(B, 1, &row, 1, &col, v_neighbor, INSERT_VALUES)); 253c4762a1bSJed Brown 254c4762a1bSJed Brown /* y neighbors */ 2559371c9d4SSatish Balay col.i = i; 2569371c9d4SSatish Balay col.j = j - 1; 2579371c9d4SSatish Balay col.k = k; 2589566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlockedStencil(B, 1, &row, 1, &col, v_neighbor, INSERT_VALUES)); 2599371c9d4SSatish Balay col.i = i; 2609371c9d4SSatish Balay col.j = j + 1; 2619371c9d4SSatish Balay col.k = k; 2629566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlockedStencil(B, 1, &row, 1, &col, v_neighbor, INSERT_VALUES)); 263c4762a1bSJed Brown 264c4762a1bSJed Brown /* z neighbors */ 2659371c9d4SSatish Balay col.i = i; 2669371c9d4SSatish Balay col.j = j; 2679371c9d4SSatish Balay col.k = k - 1; 2689566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlockedStencil(B, 1, &row, 1, &col, v_neighbor, INSERT_VALUES)); 2699371c9d4SSatish Balay col.i = i; 2709371c9d4SSatish Balay col.j = j; 2719371c9d4SSatish Balay col.k = k + 1; 2729566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlockedStencil(B, 1, &row, 1, &col, v_neighbor, INSERT_VALUES)); 273c4762a1bSJed Brown } 274c4762a1bSJed Brown } 275c4762a1bSJed Brown } 276c4762a1bSJed Brown } 2779566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 2789566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 2799566063dSJacob Faibussowitsch PetscCall(PetscFree(v)); 2809566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rand)); 281*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 282c4762a1bSJed Brown } 283c4762a1bSJed Brown 284c4762a1bSJed Brown /*TEST 285c4762a1bSJed Brown 286c4762a1bSJed Brown test: 287c4762a1bSJed Brown args: -dm_mat_type aij -dof 1 288c4762a1bSJed Brown output_file: output/ex129.out 289c4762a1bSJed Brown 290c4762a1bSJed Brown test: 291c4762a1bSJed Brown suffix: 2 292c4762a1bSJed Brown args: -dm_mat_type aij -dof 1 -inplacelu 293c4762a1bSJed Brown output_file: output/ex129.out 294c4762a1bSJed Brown 295c4762a1bSJed Brown TEST*/ 296