1c4762a1bSJed Brown static const char help[] = "Minimum surface problem in 2D.\n\ 2c4762a1bSJed Brown Uses 2-dimensional distributed arrays.\n\ 3c4762a1bSJed Brown \n\ 4c4762a1bSJed Brown Solves the linear systems via multilevel methods \n\ 5c4762a1bSJed Brown \n\n"; 6c4762a1bSJed Brown 7c4762a1bSJed Brown /* 8c4762a1bSJed Brown 9c4762a1bSJed Brown This example models the partial differential equation 10c4762a1bSJed Brown 11c4762a1bSJed Brown - Div((1 + ||GRAD T||^2)^(1/2) (GRAD T)) = 0. 12c4762a1bSJed Brown 13c4762a1bSJed Brown in the unit square, which is uniformly discretized in each of x and 14c4762a1bSJed Brown y in this simple encoding. The degrees of freedom are vertex centered 15c4762a1bSJed Brown 16c4762a1bSJed Brown A finite difference approximation with the usual 5-point stencil 17c4762a1bSJed Brown is used to discretize the boundary value problem to obtain a 18c4762a1bSJed Brown nonlinear system of equations. 19c4762a1bSJed Brown 20c4762a1bSJed Brown */ 21c4762a1bSJed Brown 22c4762a1bSJed Brown #include <petscsnes.h> 23c4762a1bSJed Brown #include <petscdm.h> 24c4762a1bSJed Brown #include <petscdmda.h> 25c4762a1bSJed Brown 26c4762a1bSJed Brown extern PetscErrorCode FormFunctionLocal(DMDALocalInfo *, PetscScalar **, PetscScalar **, void *); 27c4762a1bSJed Brown 28*9371c9d4SSatish Balay int main(int argc, char **argv) { 29c4762a1bSJed Brown SNES snes; 30c4762a1bSJed Brown PetscInt its, lits; 31c4762a1bSJed Brown PetscReal litspit; 32c4762a1bSJed Brown DM da; 33c4762a1bSJed Brown 34327415f7SBarry Smith PetscFunctionBeginUser; 359566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 36c4762a1bSJed Brown /* 37c4762a1bSJed Brown Set the DMDA (grid structure) for the grids. 38c4762a1bSJed Brown */ 399566063dSJacob Faibussowitsch PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, 5, 5, PETSC_DECIDE, PETSC_DECIDE, 1, 1, 0, 0, &da)); 409566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 419566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 429566063dSJacob Faibussowitsch PetscCall(DMDASNESSetFunctionLocal(da, INSERT_VALUES, (PetscErrorCode(*)(DMDALocalInfo *, void *, void *, void *))FormFunctionLocal, NULL)); 439566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 449566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, da)); 459566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 46c4762a1bSJed Brown 479566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 48c4762a1bSJed Brown 499566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, 0, 0)); 509566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &its)); 519566063dSJacob Faibussowitsch PetscCall(SNESGetLinearSolveIterations(snes, &lits)); 52c4762a1bSJed Brown litspit = ((PetscReal)lits) / ((PetscReal)its); 5363a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Number of SNES iterations = %" PetscInt_FMT "\n", its)); 5463a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Number of Linear iterations = %" PetscInt_FMT "\n", lits)); 559566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Average Linear its / SNES = %e\n", (double)litspit)); 56c4762a1bSJed Brown 579566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 589566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 59b122ec5aSJacob Faibussowitsch return 0; 60c4762a1bSJed Brown } 61c4762a1bSJed Brown 62*9371c9d4SSatish Balay PetscErrorCode FormFunctionLocal(DMDALocalInfo *info, PetscScalar **t, PetscScalar **f, void *ptr) { 63c4762a1bSJed Brown PetscInt i, j; 64c4762a1bSJed Brown PetscScalar hx, hy; 65c4762a1bSJed Brown PetscScalar gradup, graddown, gradleft, gradright, gradx, grady; 66c4762a1bSJed Brown PetscScalar coeffup, coeffdown, coeffleft, coeffright; 67c4762a1bSJed Brown 68c4762a1bSJed Brown PetscFunctionBeginUser; 69*9371c9d4SSatish Balay hx = 1.0 / (PetscReal)(info->mx - 1); 70*9371c9d4SSatish Balay hy = 1.0 / (PetscReal)(info->my - 1); 71c4762a1bSJed Brown 72c4762a1bSJed Brown /* Evaluate function */ 73c4762a1bSJed Brown for (j = info->ys; j < info->ys + info->ym; j++) { 74c4762a1bSJed Brown for (i = info->xs; i < info->xs + info->xm; i++) { 75c4762a1bSJed Brown if (i == 0 || i == info->mx - 1 || j == 0 || j == info->my - 1) { 76c4762a1bSJed Brown f[j][i] = t[j][i] - (1.0 - (2.0 * hx * (PetscReal)i - 1.0) * (2.0 * hx * (PetscReal)i - 1.0)); 77c4762a1bSJed Brown } else { 78c4762a1bSJed Brown gradup = (t[j + 1][i] - t[j][i]) / hy; 79c4762a1bSJed Brown graddown = (t[j][i] - t[j - 1][i]) / hy; 80c4762a1bSJed Brown gradright = (t[j][i + 1] - t[j][i]) / hx; 81c4762a1bSJed Brown gradleft = (t[j][i] - t[j][i - 1]) / hx; 82c4762a1bSJed Brown 83c4762a1bSJed Brown gradx = .5 * (t[j][i + 1] - t[j][i - 1]) / hx; 84c4762a1bSJed Brown grady = .5 * (t[j + 1][i] - t[j - 1][i]) / hy; 85c4762a1bSJed Brown 86c4762a1bSJed Brown coeffup = 1.0 / PetscSqrtScalar(1.0 + gradup * gradup + gradx * gradx); 87c4762a1bSJed Brown coeffdown = 1.0 / PetscSqrtScalar(1.0 + graddown * graddown + gradx * gradx); 88c4762a1bSJed Brown 89c4762a1bSJed Brown coeffleft = 1.0 / PetscSqrtScalar(1.0 + gradleft * gradleft + grady * grady); 90c4762a1bSJed Brown coeffright = 1.0 / PetscSqrtScalar(1.0 + gradright * gradright + grady * grady); 91c4762a1bSJed Brown 92c4762a1bSJed Brown f[j][i] = (coeffup * gradup - coeffdown * graddown) * hx + (coeffright * gradright - coeffleft * gradleft) * hy; 93c4762a1bSJed Brown } 94c4762a1bSJed Brown } 95c4762a1bSJed Brown } 96c4762a1bSJed Brown PetscFunctionReturn(0); 97c4762a1bSJed Brown } 98c4762a1bSJed Brown 99c4762a1bSJed Brown /*TEST 100c4762a1bSJed Brown 101c4762a1bSJed Brown test: 102c4762a1bSJed Brown args: -pc_type mg -da_refine 1 -ksp_type fgmres 103c4762a1bSJed Brown 104c4762a1bSJed Brown test: 105c4762a1bSJed Brown suffix: 2 106c4762a1bSJed Brown nsize: 2 107c4762a1bSJed Brown args: -pc_type mg -da_refine 1 -ksp_type fgmres 108c4762a1bSJed Brown 10941ba4c6cSHeeho Park test: 11041ba4c6cSHeeho Park suffix: 3 11141ba4c6cSHeeho Park nsize: 2 11241ba4c6cSHeeho Park args: -pc_type mg -da_refine 1 -ksp_type fgmres -snes_type newtontrdc -snes_trdc_use_cauchy false 11341ba4c6cSHeeho Park 11441ba4c6cSHeeho Park test: 11541ba4c6cSHeeho Park suffix: 4 11641ba4c6cSHeeho Park nsize: 2 11741ba4c6cSHeeho Park args: -pc_type mg -da_refine 1 -ksp_type fgmres -snes_type newtontrdc 11841ba4c6cSHeeho Park filter: sed -e "s/SNES iterations = 1[1-3]/SNES iterations = 13/g" |sed -e "s/Linear iterations = 2[8-9]/Linear iterations = 29/g" |sed -e "s/Linear iterations = 3[0-1]/Linear iterations = 29/g" 11941ba4c6cSHeeho Park 120c4762a1bSJed Brown TEST*/ 121