xref: /petsc/src/snes/tutorials/ex25.c (revision d0e6bf2ad94dcc89b258ce16c7987200a4714786)
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 
main(int argc,char ** argv)28d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
29d71ae5a4SJacob Faibussowitsch {
30c4762a1bSJed Brown   SNES      snes;
31c4762a1bSJed Brown   PetscInt  its, lits;
32c4762a1bSJed Brown   PetscReal litspit;
33c4762a1bSJed Brown   DM        da;
34c4762a1bSJed Brown 
35327415f7SBarry Smith   PetscFunctionBeginUser;
369566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
37c4762a1bSJed Brown   /*
38c4762a1bSJed Brown       Set the DMDA (grid structure) for the grids.
39c4762a1bSJed Brown   */
409566063dSJacob 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));
419566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
429566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
439566063dSJacob Faibussowitsch   PetscCall(DMDASNESSetFunctionLocal(da, INSERT_VALUES, (PetscErrorCode (*)(DMDALocalInfo *, void *, void *, void *))FormFunctionLocal, NULL));
449566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
459566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(snes, da));
469566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
47c4762a1bSJed Brown 
489566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
49c4762a1bSJed Brown 
509566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, 0, 0));
519566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes, &its));
529566063dSJacob Faibussowitsch   PetscCall(SNESGetLinearSolveIterations(snes, &lits));
53c4762a1bSJed Brown   litspit = ((PetscReal)lits) / ((PetscReal)its);
5463a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Number of SNES iterations = %" PetscInt_FMT "\n", its));
5563a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Number of Linear iterations = %" PetscInt_FMT "\n", lits));
569566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Average Linear its / SNES = %e\n", (double)litspit));
57c4762a1bSJed Brown 
589566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
599566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
60b122ec5aSJacob Faibussowitsch   return 0;
61c4762a1bSJed Brown }
62c4762a1bSJed Brown 
FormFunctionLocal(DMDALocalInfo * info,PetscScalar ** t,PetscScalar ** f,void * ptr)63d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunctionLocal(DMDALocalInfo *info, PetscScalar **t, PetscScalar **f, void *ptr)
64d71ae5a4SJacob Faibussowitsch {
65c4762a1bSJed Brown   PetscInt    i, j;
66c4762a1bSJed Brown   PetscScalar hx, hy;
67c4762a1bSJed Brown   PetscScalar gradup, graddown, gradleft, gradright, gradx, grady;
68c4762a1bSJed Brown   PetscScalar coeffup, coeffdown, coeffleft, coeffright;
69c4762a1bSJed Brown 
70c4762a1bSJed Brown   PetscFunctionBeginUser;
719371c9d4SSatish Balay   hx = 1.0 / (PetscReal)(info->mx - 1);
729371c9d4SSatish Balay   hy = 1.0 / (PetscReal)(info->my - 1);
73c4762a1bSJed Brown 
74c4762a1bSJed Brown   /* Evaluate function */
75c4762a1bSJed Brown   for (j = info->ys; j < info->ys + info->ym; j++) {
76c4762a1bSJed Brown     for (i = info->xs; i < info->xs + info->xm; i++) {
77c4762a1bSJed Brown       if (i == 0 || i == info->mx - 1 || j == 0 || j == info->my - 1) {
78c4762a1bSJed Brown         f[j][i] = t[j][i] - (1.0 - (2.0 * hx * (PetscReal)i - 1.0) * (2.0 * hx * (PetscReal)i - 1.0));
79c4762a1bSJed Brown       } else {
80c4762a1bSJed Brown         gradup    = (t[j + 1][i] - t[j][i]) / hy;
81c4762a1bSJed Brown         graddown  = (t[j][i] - t[j - 1][i]) / hy;
82c4762a1bSJed Brown         gradright = (t[j][i + 1] - t[j][i]) / hx;
83c4762a1bSJed Brown         gradleft  = (t[j][i] - t[j][i - 1]) / hx;
84c4762a1bSJed Brown 
85c4762a1bSJed Brown         gradx = .5 * (t[j][i + 1] - t[j][i - 1]) / hx;
86c4762a1bSJed Brown         grady = .5 * (t[j + 1][i] - t[j - 1][i]) / hy;
87c4762a1bSJed Brown 
88c4762a1bSJed Brown         coeffup   = 1.0 / PetscSqrtScalar(1.0 + gradup * gradup + gradx * gradx);
89c4762a1bSJed Brown         coeffdown = 1.0 / PetscSqrtScalar(1.0 + graddown * graddown + gradx * gradx);
90c4762a1bSJed Brown 
91c4762a1bSJed Brown         coeffleft  = 1.0 / PetscSqrtScalar(1.0 + gradleft * gradleft + grady * grady);
92c4762a1bSJed Brown         coeffright = 1.0 / PetscSqrtScalar(1.0 + gradright * gradright + grady * grady);
93c4762a1bSJed Brown 
94c4762a1bSJed Brown         f[j][i] = (coeffup * gradup - coeffdown * graddown) * hx + (coeffright * gradright - coeffleft * gradleft) * hy;
95c4762a1bSJed Brown       }
96c4762a1bSJed Brown     }
97c4762a1bSJed Brown   }
983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
99c4762a1bSJed Brown }
100c4762a1bSJed Brown 
101c4762a1bSJed Brown /*TEST
102c4762a1bSJed Brown 
103c4762a1bSJed Brown    test:
104c4762a1bSJed Brown       args: -pc_type mg -da_refine 1 -ksp_type fgmres
105c4762a1bSJed Brown 
106c4762a1bSJed Brown    test:
107c4762a1bSJed Brown       suffix: 2
108c4762a1bSJed Brown       nsize: 2
109c4762a1bSJed Brown       args: -pc_type mg -da_refine 1 -ksp_type fgmres
110c4762a1bSJed Brown 
11141ba4c6cSHeeho Park    test:
11241ba4c6cSHeeho Park       suffix: 3
11341ba4c6cSHeeho Park       nsize: 2
11441ba4c6cSHeeho Park       args: -pc_type mg -da_refine 1 -ksp_type fgmres -snes_type newtontrdc -snes_trdc_use_cauchy false
11541ba4c6cSHeeho Park 
11641ba4c6cSHeeho Park    test:
11741ba4c6cSHeeho Park       suffix: 4
11841ba4c6cSHeeho Park       nsize: 2
11941ba4c6cSHeeho Park       args: -pc_type mg -da_refine 1 -ksp_type fgmres -snes_type newtontrdc
120*c89294f2SStefano Zampini       filter: sed -e "s/SNES iterations = 1[0-3]/SNES iterations = 13/g" |sed -e "s/Linear iterations = 2[7-9]/Linear iterations = 29/g" |sed -e "s/Linear iterations = 3[0-1]/Linear iterations = 29/g"
12141ba4c6cSHeeho Park 
122c4762a1bSJed Brown TEST*/
123