1736140d5SStefano Zampini static const char help[] = "\
2736140d5SStefano Zampini Minimum surface area problem in 2D using DMPLEX.\n\
3736140d5SStefano Zampini It solves the unconstrained minimization problem \n\
4736140d5SStefano Zampini \n\
5736140d5SStefano Zampini argmin \\int_\\Omega (1 + ||u||^2), u = g on \\partial\\Omega.\n\
6736140d5SStefano Zampini \n\
7736140d5SStefano Zampini This example shows how to setup an optimization problem using DMPLEX FEM routines.\n\
8736140d5SStefano Zampini It supports nonlinear domain decomposition and multilevel solvers.\n\
9736140d5SStefano Zampini \n\n";
10736140d5SStefano Zampini
11736140d5SStefano Zampini #include <petscdmplex.h>
12736140d5SStefano Zampini #include <petscsnes.h>
13736140d5SStefano Zampini #include <petscds.h>
14736140d5SStefano Zampini
15736140d5SStefano Zampini /* The pointwise evaluation routine for the objective function */
objective_2d(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar obj[])16736140d5SStefano Zampini static void objective_2d(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar obj[])
17736140d5SStefano Zampini {
18736140d5SStefano Zampini const PetscScalar ux2 = PetscSqr(u_x[0]);
19736140d5SStefano Zampini const PetscScalar uy2 = PetscSqr(u_x[1]);
20736140d5SStefano Zampini
216497c311SBarry Smith obj[0] = PetscSqrtReal(PetscAbsScalar(1 + ux2 + uy2));
22736140d5SStefano Zampini }
23736140d5SStefano Zampini
24736140d5SStefano Zampini /* The pointwise evaluation routine for the gradient wrt the gradients of the FEM basis */
gradient_1_2d(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar f[])25736140d5SStefano Zampini static void gradient_1_2d(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f[])
26736140d5SStefano Zampini {
27736140d5SStefano Zampini const PetscScalar ux2 = PetscSqr(u_x[0]);
28736140d5SStefano Zampini const PetscScalar uy2 = PetscSqr(u_x[1]);
296497c311SBarry Smith const PetscScalar v = 1. / PetscSqrtReal(PetscAbsScalar(1 + ux2 + uy2));
30736140d5SStefano Zampini
31736140d5SStefano Zampini f[0] = v * u_x[0];
32736140d5SStefano Zampini f[1] = v * u_x[1];
33736140d5SStefano Zampini }
34736140d5SStefano Zampini
35736140d5SStefano Zampini /* The pointwise evaluation routine for the hessian wrt the gradients of the FEM basis */
hessian_11_2d(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,PetscReal u_tShift,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar jac[])36736140d5SStefano Zampini static void hessian_11_2d(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar jac[])
37736140d5SStefano Zampini {
38736140d5SStefano Zampini const PetscScalar ux2 = PetscSqr(u_x[0]);
39736140d5SStefano Zampini const PetscScalar uy2 = PetscSqr(u_x[1]);
40736140d5SStefano Zampini const PetscScalar uxy = u_x[0] * u_x[1];
416497c311SBarry Smith const PetscScalar v1 = 1. / PetscSqrtReal(PetscAbsScalar(1 + ux2 + uy2));
426497c311SBarry Smith const PetscScalar v2 = v1 / (1 + ux2 + uy2);
43736140d5SStefano Zampini
44736140d5SStefano Zampini jac[0] = v1 - v2 * ux2;
45736140d5SStefano Zampini jac[1] = -v2 * uxy;
46736140d5SStefano Zampini jac[2] = -v2 * uxy;
47736140d5SStefano Zampini jac[3] = v1 - v2 * uy2;
48736140d5SStefano Zampini }
49736140d5SStefano Zampini
50736140d5SStefano Zampini /* The boundary conditions we impose */
sins_2d(PetscInt dim,PetscReal time,const PetscReal xx[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)51*2a8381b2SBarry Smith static PetscErrorCode sins_2d(PetscInt dim, PetscReal time, const PetscReal xx[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
52736140d5SStefano Zampini {
53736140d5SStefano Zampini PetscFunctionBeginUser;
54736140d5SStefano Zampini const PetscReal pi2 = PETSC_PI * 2;
55736140d5SStefano Zampini const PetscReal x = xx[0];
56736140d5SStefano Zampini const PetscReal y = xx[1];
57736140d5SStefano Zampini
58736140d5SStefano Zampini *u = (y - 0.5) * PetscSinReal(pi2 * x) + (x - 0.5) * PetscSinReal(pi2 * y);
59736140d5SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS);
60736140d5SStefano Zampini }
61736140d5SStefano Zampini
CreateBCLabel(DM dm,const char name[])62736140d5SStefano Zampini static PetscErrorCode CreateBCLabel(DM dm, const char name[])
63736140d5SStefano Zampini {
64736140d5SStefano Zampini DM plex;
65736140d5SStefano Zampini DMLabel label;
66736140d5SStefano Zampini
67736140d5SStefano Zampini PetscFunctionBeginUser;
68736140d5SStefano Zampini PetscCall(DMCreateLabel(dm, name));
69736140d5SStefano Zampini PetscCall(DMGetLabel(dm, name, &label));
70736140d5SStefano Zampini PetscCall(DMConvert(dm, DMPLEX, &plex));
71736140d5SStefano Zampini PetscCall(DMPlexMarkBoundaryFaces(plex, 1, label));
72736140d5SStefano Zampini PetscCall(DMDestroy(&plex));
73736140d5SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS);
74736140d5SStefano Zampini }
75736140d5SStefano Zampini
CreateMesh(MPI_Comm comm,DM * dm)76736140d5SStefano Zampini static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm)
77736140d5SStefano Zampini {
78736140d5SStefano Zampini PetscFunctionBeginUser;
79736140d5SStefano Zampini PetscCall(DMCreate(comm, dm));
80736140d5SStefano Zampini PetscCall(DMSetType(*dm, DMPLEX));
81736140d5SStefano Zampini PetscCall(DMSetFromOptions(*dm));
82736140d5SStefano Zampini PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
83736140d5SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS);
84736140d5SStefano Zampini }
85736140d5SStefano Zampini
SetupProblem(DM dm)86736140d5SStefano Zampini static PetscErrorCode SetupProblem(DM dm)
87736140d5SStefano Zampini {
88736140d5SStefano Zampini PetscDS ds;
89736140d5SStefano Zampini DMLabel label;
90736140d5SStefano Zampini PetscInt dim;
91736140d5SStefano Zampini const PetscInt id = 1;
92736140d5SStefano Zampini
93736140d5SStefano Zampini PetscFunctionBeginUser;
94736140d5SStefano Zampini PetscCall(DMGetDS(dm, &ds));
95736140d5SStefano Zampini PetscCall(DMGetDimension(dm, &dim));
9600045ab3SPierre Jolivet PetscCheck(dim == 2, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only for 2-D");
97736140d5SStefano Zampini PetscCall(PetscDSSetObjective(ds, 0, objective_2d));
98736140d5SStefano Zampini PetscCall(PetscDSSetResidual(ds, 0, NULL, gradient_1_2d));
99736140d5SStefano Zampini PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, hessian_11_2d));
100736140d5SStefano Zampini PetscCall(DMGetLabel(dm, "marker", &label));
10157d50842SBarry Smith PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "data", label, 1, &id, 0, 0, NULL, (PetscVoidFn *)sins_2d, NULL, NULL, NULL));
102736140d5SStefano Zampini PetscCall(DMPlexSetSNESLocalFEM(dm, PETSC_TRUE, NULL));
103736140d5SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS);
104736140d5SStefano Zampini }
105736140d5SStefano Zampini
SetupDiscretization(DM dm)106736140d5SStefano Zampini static PetscErrorCode SetupDiscretization(DM dm)
107736140d5SStefano Zampini {
108736140d5SStefano Zampini DM plex, cdm = dm;
109736140d5SStefano Zampini PetscFE fe;
110736140d5SStefano Zampini PetscBool simplex;
111736140d5SStefano Zampini PetscInt dim;
112736140d5SStefano Zampini
113736140d5SStefano Zampini PetscFunctionBeginUser;
114736140d5SStefano Zampini PetscCall(DMGetDimension(dm, &dim));
115736140d5SStefano Zampini PetscCall(DMConvert(dm, DMPLEX, &plex));
116736140d5SStefano Zampini PetscCall(DMPlexIsSimplex(plex, &simplex));
117736140d5SStefano Zampini PetscCall(DMDestroy(&plex));
118736140d5SStefano Zampini PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, NULL, -1, &fe));
119736140d5SStefano Zampini PetscCall(PetscObjectSetName((PetscObject)fe, "potential"));
120736140d5SStefano Zampini PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe));
121736140d5SStefano Zampini PetscCall(DMCreateDS(dm));
122736140d5SStefano Zampini PetscCall(SetupProblem(dm));
123736140d5SStefano Zampini while (cdm) {
124736140d5SStefano Zampini PetscBool hasLabel;
125736140d5SStefano Zampini
126736140d5SStefano Zampini PetscCall(DMHasLabel(cdm, "marker", &hasLabel));
127736140d5SStefano Zampini if (!hasLabel) PetscCall(CreateBCLabel(cdm, "marker"));
128736140d5SStefano Zampini PetscCall(DMCopyDisc(dm, cdm));
129736140d5SStefano Zampini PetscCall(DMGetCoarseDM(cdm, &cdm));
130736140d5SStefano Zampini }
131736140d5SStefano Zampini PetscCall(PetscFEDestroy(&fe));
132736140d5SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS);
133736140d5SStefano Zampini }
134736140d5SStefano Zampini
main(int argc,char ** argv)135736140d5SStefano Zampini int main(int argc, char **argv)
136736140d5SStefano Zampini {
137736140d5SStefano Zampini DM dm; /* Problem specification */
138736140d5SStefano Zampini SNES snes; /* nonlinear solver */
139736140d5SStefano Zampini Vec u; /* solution vector */
140736140d5SStefano Zampini
141736140d5SStefano Zampini PetscFunctionBeginUser;
142736140d5SStefano Zampini PetscCall(PetscInitialize(&argc, &argv, NULL, help));
143736140d5SStefano Zampini PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
144736140d5SStefano Zampini PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm));
145736140d5SStefano Zampini PetscCall(SNESSetDM(snes, dm));
146736140d5SStefano Zampini
147736140d5SStefano Zampini PetscCall(SetupDiscretization(dm));
148736140d5SStefano Zampini
149736140d5SStefano Zampini PetscCall(SNESSetFromOptions(snes));
150736140d5SStefano Zampini
151736140d5SStefano Zampini PetscCall(DMCreateGlobalVector(dm, &u));
152736140d5SStefano Zampini PetscCall(SNESSolve(snes, NULL, u));
153736140d5SStefano Zampini
154736140d5SStefano Zampini PetscCall(VecDestroy(&u));
155736140d5SStefano Zampini PetscCall(SNESDestroy(&snes));
156736140d5SStefano Zampini PetscCall(DMDestroy(&dm));
157736140d5SStefano Zampini PetscCall(PetscFinalize());
158736140d5SStefano Zampini return 0;
159736140d5SStefano Zampini }
160736140d5SStefano Zampini
161736140d5SStefano Zampini /*TEST
162736140d5SStefano Zampini
163736140d5SStefano Zampini test:
164736140d5SStefano Zampini requires: !single
165736140d5SStefano Zampini suffix: qn_nasm
166736140d5SStefano Zampini filter: sed -e "s/CONVERGED_FNORM_ABS/CONVERGED_FNORM_RELATIVE/g"
167736140d5SStefano Zampini nsize: 4
168736140d5SStefano Zampini args: -petscspace_degree 1 -dm_refine 2 -snes_type qn -snes_npc_side {{left right}separate output} -npc_snes_type nasm -npc_snes_nasm_type restrict -npc_sub_snes_linesearch_order 1 -npc_sub_snes_linesearch_type bt -dm_plex_dd_overlap 1 -snes_linesearch_type bt -snes_linesearch_order 1 -npc_sub_pc_type lu -npc_sub_ksp_type preonly -snes_converged_reason -snes_monitor_short -petscpartitioner_type simple -npc_sub_snes_max_it 1 -dm_plex_simplex 0 -snes_rtol 1.e-6
169736140d5SStefano Zampini
170736140d5SStefano Zampini TEST*/
171