1 static char help[] = "One-Shot Multigrid for Parameter Estimation Problem for the Poisson Equation.\n\ 2 Using the Interior Point Method.\n\n\n"; 3 4 /*F 5 We are solving the parameter estimation problem for the Laplacian. We will ask to minimize a Lagrangian 6 function over $y$ and $u$, given by 7 \begin{align} 8 L(u, a, \lambda) = \frac{1}{2} || Qu - d_A ||^2 || Qu - d_B ||^2 + \frac{\beta}{2} || L (a - a_r) ||^2 + \lambda F(u; a) 9 \end{align} 10 where $Q$ is a sampling operator, $L$ is a regularization operator, $F$ defines the PDE. 11 12 Currently, we have perfect information, meaning $Q = I$, and then we need no regularization, $L = I$. We 13 also give the null vector for the reference control $a_r$. Right now $\beta = 1$. 14 15 The PDE will be the Laplace equation with homogeneous boundary conditions 16 \begin{align} 17 -Delta u = a 18 \end{align} 19 20 F*/ 21 22 #include <petsc.h> 23 #include <petscfe.h> 24 25 typedef enum { 26 RUN_FULL, 27 RUN_TEST 28 } RunType; 29 30 typedef struct { 31 RunType runType; /* Whether to run tests, or solve the full problem */ 32 PetscBool useDualPenalty; /* Penalize deviation from both goals */ 33 } AppCtx; 34 35 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) { 36 const char *runTypes[2] = {"full", "test"}; 37 PetscInt run; 38 39 PetscFunctionBeginUser; 40 options->runType = RUN_FULL; 41 options->useDualPenalty = PETSC_FALSE; 42 PetscOptionsBegin(comm, "", "Inverse Problem Options", "DMPLEX"); 43 run = options->runType; 44 PetscCall(PetscOptionsEList("-run_type", "The run type", "ex2.c", runTypes, 2, runTypes[options->runType], &run, NULL)); 45 options->runType = (RunType)run; 46 PetscCall(PetscOptionsBool("-use_dual_penalty", "Penalize deviation from both goals", "ex2.c", options->useDualPenalty, &options->useDualPenalty, NULL)); 47 PetscOptionsEnd(); 48 PetscFunctionReturn(0); 49 } 50 51 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) { 52 PetscFunctionBeginUser; 53 PetscCall(DMCreate(comm, dm)); 54 PetscCall(DMSetType(*dm, DMPLEX)); 55 PetscCall(DMSetFromOptions(*dm)); 56 PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 57 PetscFunctionReturn(0); 58 } 59 60 void f0_u(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 f0[]) { 61 f0[0] = (u[0] - (x[0] * x[0] + x[1] * x[1])); 62 } 63 void f0_u_full(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 f0[]) { 64 f0[0] = (u[0] - (x[0] * x[0] + x[1] * x[1])) * PetscSqr(u[0] - (sin(2.0 * PETSC_PI * x[0]) * sin(2.0 * PETSC_PI * x[1]))) + PetscSqr(u[0] - (x[0] * x[0] + x[1] * x[1])) * (u[0] - (sin(2.0 * PETSC_PI * x[0]) * sin(2.0 * PETSC_PI * x[1]))); 65 } 66 void f1_u(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 f1[]) { 67 PetscInt d; 68 for (d = 0; d < dim; ++d) f1[d] = u_x[dim * 2 + d]; 69 } 70 void g0_uu(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 g0[]) { 71 g0[0] = 1.0; 72 } 73 void g0_uu_full(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 g0[]) { 74 g0[0] = PetscSqr(u[0] - sin(2.0 * PETSC_PI * x[0]) * sin(2.0 * PETSC_PI * x[1])) + PetscSqr(u[0] - (x[0] * x[0] + x[1] * x[1])) - 2.0 * ((x[0] * x[0] + x[1] * x[1]) + (sin(2.0 * PETSC_PI * x[0]) * sin(2.0 * PETSC_PI * x[1]))) * u[0] + 4.0 * (x[0] * x[0] + x[1] * x[1]) * (sin(2.0 * PETSC_PI * x[0]) * sin(2.0 * PETSC_PI * x[1])); 75 } 76 void g3_ul(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 g3[]) { 77 PetscInt d; 78 for (d = 0; d < dim; ++d) g3[d * dim + d] = 1.0; 79 } 80 81 void f0_a(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 f0[]) { 82 f0[0] = u[1] - 4.0 /* 0.0 */ + u[2]; 83 } 84 void g0_aa(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 g0[]) { 85 g0[0] = 1.0; 86 } 87 void g0_al(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 g0[]) { 88 g0[0] = 1.0; 89 } 90 91 void f0_l(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 f0[]) { 92 f0[0] = u[1]; 93 } 94 void f1_l(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 f1[]) { 95 PetscInt d; 96 for (d = 0; d < dim; ++d) f1[d] = u_x[d]; 97 } 98 void g0_la(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 g0[]) { 99 g0[0] = 1.0; 100 } 101 void g3_lu(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 g3[]) { 102 PetscInt d; 103 for (d = 0; d < dim; ++d) g3[d * dim + d] = 1.0; 104 } 105 106 /* 107 In 2D for Dirichlet conditions with a variable coefficient, we use exact solution: 108 109 u = x^2 + y^2 110 a = 4 111 d_A = 4 112 d_B = sin(2*pi*x[0]) * sin(2*pi*x[1]) 113 114 so that 115 116 -\Delta u + a = -4 + 4 = 0 117 */ 118 PetscErrorCode quadratic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) { 119 *u = x[0] * x[0] + x[1] * x[1]; 120 return 0; 121 } 122 PetscErrorCode constant_a_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *a, void *ctx) { 123 *a = 4; 124 return 0; 125 } 126 PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *l, void *ctx) { 127 *l = 0.0; 128 return 0; 129 } 130 131 PetscErrorCode SetupProblem(DM dm, AppCtx *user) { 132 PetscDS ds; 133 DMLabel label; 134 const PetscInt id = 1; 135 136 PetscFunctionBeginUser; 137 PetscCall(DMGetDS(dm, &ds)); 138 PetscCall(PetscDSSetResidual(ds, 0, user->useDualPenalty == PETSC_TRUE ? f0_u_full : f0_u, f1_u)); 139 PetscCall(PetscDSSetResidual(ds, 1, f0_a, NULL)); 140 PetscCall(PetscDSSetResidual(ds, 2, f0_l, f1_l)); 141 PetscCall(PetscDSSetJacobian(ds, 0, 0, user->useDualPenalty == PETSC_TRUE ? g0_uu_full : g0_uu, NULL, NULL, NULL)); 142 PetscCall(PetscDSSetJacobian(ds, 0, 2, NULL, NULL, NULL, g3_ul)); 143 PetscCall(PetscDSSetJacobian(ds, 1, 1, g0_aa, NULL, NULL, NULL)); 144 PetscCall(PetscDSSetJacobian(ds, 1, 2, g0_al, NULL, NULL, NULL)); 145 PetscCall(PetscDSSetJacobian(ds, 2, 1, g0_la, NULL, NULL, NULL)); 146 PetscCall(PetscDSSetJacobian(ds, 2, 0, NULL, NULL, NULL, g3_lu)); 147 148 PetscCall(PetscDSSetExactSolution(ds, 0, quadratic_u_2d, NULL)); 149 PetscCall(PetscDSSetExactSolution(ds, 1, constant_a_2d, NULL)); 150 PetscCall(PetscDSSetExactSolution(ds, 2, zero, NULL)); 151 PetscCall(DMGetLabel(dm, "marker", &label)); 152 PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)())quadratic_u_2d, NULL, user, NULL)); 153 PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 1, 0, NULL, (void (*)())constant_a_2d, NULL, user, NULL)); 154 PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 2, 0, NULL, (void (*)())zero, NULL, user, NULL)); 155 PetscFunctionReturn(0); 156 } 157 158 PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) { 159 DM cdm = dm; 160 const PetscInt dim = 2; 161 PetscFE fe[3]; 162 PetscInt f; 163 MPI_Comm comm; 164 165 PetscFunctionBeginUser; 166 /* Create finite element */ 167 PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 168 PetscCall(PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "potential_", -1, &fe[0])); 169 PetscCall(PetscObjectSetName((PetscObject)fe[0], "potential")); 170 PetscCall(PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "charge_", -1, &fe[1])); 171 PetscCall(PetscObjectSetName((PetscObject)fe[1], "charge")); 172 PetscCall(PetscFECopyQuadrature(fe[0], fe[1])); 173 PetscCall(PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "multiplier_", -1, &fe[2])); 174 PetscCall(PetscObjectSetName((PetscObject)fe[2], "multiplier")); 175 PetscCall(PetscFECopyQuadrature(fe[0], fe[2])); 176 /* Set discretization and boundary conditions for each mesh */ 177 for (f = 0; f < 3; ++f) PetscCall(DMSetField(dm, f, NULL, (PetscObject)fe[f])); 178 PetscCall(DMCreateDS(cdm)); 179 PetscCall(SetupProblem(dm, user)); 180 while (cdm) { 181 PetscCall(DMCopyDisc(dm, cdm)); 182 PetscCall(DMGetCoarseDM(cdm, &cdm)); 183 } 184 for (f = 0; f < 3; ++f) PetscCall(PetscFEDestroy(&fe[f])); 185 PetscFunctionReturn(0); 186 } 187 188 int main(int argc, char **argv) { 189 DM dm; 190 SNES snes; 191 Vec u, r; 192 AppCtx user; 193 194 PetscFunctionBeginUser; 195 PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 196 PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 197 PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 198 PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 199 PetscCall(SNESSetDM(snes, dm)); 200 PetscCall(SetupDiscretization(dm, &user)); 201 202 PetscCall(DMCreateGlobalVector(dm, &u)); 203 PetscCall(PetscObjectSetName((PetscObject)u, "solution")); 204 PetscCall(VecDuplicate(u, &r)); 205 PetscCall(DMPlexSetSNESLocalFEM(dm, &user, &user, &user)); 206 PetscCall(SNESSetFromOptions(snes)); 207 208 PetscCall(DMSNESCheckFromOptions(snes, u)); 209 if (user.runType == RUN_FULL) { 210 PetscDS ds; 211 PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 212 PetscErrorCode (*initialGuess[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx); 213 PetscReal error; 214 215 PetscCall(DMGetDS(dm, &ds)); 216 PetscCall(PetscDSGetExactSolution(ds, 0, &exactFuncs[0], NULL)); 217 PetscCall(PetscDSGetExactSolution(ds, 1, &exactFuncs[1], NULL)); 218 PetscCall(PetscDSGetExactSolution(ds, 2, &exactFuncs[2], NULL)); 219 initialGuess[0] = zero; 220 initialGuess[1] = zero; 221 initialGuess[2] = zero; 222 PetscCall(DMProjectFunction(dm, 0.0, initialGuess, NULL, INSERT_VALUES, u)); 223 PetscCall(VecViewFromOptions(u, NULL, "-initial_vec_view")); 224 PetscCall(DMComputeL2Diff(dm, 0.0, exactFuncs, NULL, u, &error)); 225 if (error < 1.0e-11) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Initial L_2 Error: < 1.0e-11\n")); 226 else PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Initial L_2 Error: %g\n", (double)error)); 227 PetscCall(SNESSolve(snes, NULL, u)); 228 PetscCall(DMComputeL2Diff(dm, 0.0, exactFuncs, NULL, u, &error)); 229 if (error < 1.0e-11) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Final L_2 Error: < 1.0e-11\n")); 230 else PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Final L_2 Error: %g\n", (double)error)); 231 } 232 PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view")); 233 234 PetscCall(VecDestroy(&u)); 235 PetscCall(VecDestroy(&r)); 236 PetscCall(SNESDestroy(&snes)); 237 PetscCall(DMDestroy(&dm)); 238 PetscCall(PetscFinalize()); 239 return 0; 240 } 241 242 /*TEST 243 244 build: 245 requires: !complex triangle 246 247 test: 248 suffix: 0 249 args: -run_type test -dmsnes_check -potential_petscspace_degree 2 -charge_petscspace_degree 1 -multiplier_petscspace_degree 1 250 251 test: 252 suffix: 1 253 args: -potential_petscspace_degree 2 -charge_petscspace_degree 1 -multiplier_petscspace_degree 1 -snes_monitor -snes_converged_reason -pc_type fieldsplit -pc_fieldsplit_0_fields 0,1 -pc_fieldsplit_1_fields 2 -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition selfp -fieldsplit_0_pc_type lu -sol_vec_view 254 255 test: 256 suffix: 2 257 args: -potential_petscspace_degree 2 -charge_petscspace_degree 1 -multiplier_petscspace_degree 1 -snes_monitor -snes_converged_reason -snes_fd -pc_type fieldsplit -pc_fieldsplit_0_fields 0,1 -pc_fieldsplit_1_fields 2 -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition selfp -fieldsplit_0_pc_type lu -sol_vec_view 258 259 TEST*/ 260