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 {RUN_FULL, RUN_TEST} RunType; 26 27 typedef struct { 28 RunType runType; /* Whether to run tests, or solve the full problem */ 29 PetscBool useDualPenalty; /* Penalize deviation from both goals */ 30 } AppCtx; 31 32 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 33 { 34 const char *runTypes[2] = {"full", "test"}; 35 PetscInt run; 36 37 PetscFunctionBeginUser; 38 options->runType = RUN_FULL; 39 options->useDualPenalty = PETSC_FALSE; 40 PetscOptionsBegin(comm, "", "Inverse Problem Options", "DMPLEX"); 41 run = options->runType; 42 PetscCall(PetscOptionsEList("-run_type", "The run type", "ex2.c", runTypes, 2, runTypes[options->runType], &run, NULL)); 43 options->runType = (RunType) run; 44 PetscCall(PetscOptionsBool("-use_dual_penalty", "Penalize deviation from both goals", "ex2.c", options->useDualPenalty, &options->useDualPenalty, NULL)); 45 PetscOptionsEnd(); 46 PetscFunctionReturn(0); 47 } 48 49 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 50 { 51 PetscFunctionBeginUser; 52 PetscCall(DMCreate(comm, dm)); 53 PetscCall(DMSetType(*dm, DMPLEX)); 54 PetscCall(DMSetFromOptions(*dm)); 55 PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 56 PetscFunctionReturn(0); 57 } 58 59 void f0_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 60 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 61 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 62 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 63 { 64 f0[0] = (u[0] - (x[0]*x[0] + x[1]*x[1])); 65 } 66 void f0_u_full(PetscInt dim, PetscInt Nf, PetscInt NfAux, 67 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 68 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 69 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 70 { 71 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]))) + 72 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]))); 73 } 74 void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 75 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 76 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 77 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 78 { 79 PetscInt d; 80 for (d = 0; d < dim; ++d) f1[d] = u_x[dim*2+d]; 81 } 82 void g0_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 83 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 84 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 85 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 86 { 87 g0[0] = 1.0; 88 } 89 void g0_uu_full(PetscInt dim, PetscInt Nf, PetscInt NfAux, 90 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 91 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 92 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 93 { 94 g0[0] = PetscSqr(u[0] - sin(2.0*PETSC_PI*x[0]) * sin(2.0*PETSC_PI*x[1])) 95 + PetscSqr(u[0] - (x[0]*x[0] + x[1]*x[1])) 96 - 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] 97 + 4.0*(x[0]*x[0] + x[1]*x[1])*(sin(2.0*PETSC_PI*x[0]) * sin(2.0*PETSC_PI*x[1])); 98 } 99 void g3_ul(PetscInt dim, PetscInt Nf, PetscInt NfAux, 100 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 101 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 102 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 103 { 104 PetscInt d; 105 for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0; 106 } 107 108 void f0_a(PetscInt dim, PetscInt Nf, PetscInt NfAux, 109 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 110 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 111 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 112 { 113 f0[0] = u[1] - 4.0 /* 0.0 */ + u[2]; 114 } 115 void g0_aa(PetscInt dim, PetscInt Nf, PetscInt NfAux, 116 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 117 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 118 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 119 { 120 g0[0] = 1.0; 121 } 122 void g0_al(PetscInt dim, PetscInt Nf, PetscInt NfAux, 123 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 124 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 125 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 126 { 127 g0[0] = 1.0; 128 } 129 130 void f0_l(PetscInt dim, PetscInt Nf, PetscInt NfAux, 131 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 132 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 133 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 134 { 135 f0[0] = u[1]; 136 } 137 void f1_l(PetscInt dim, PetscInt Nf, PetscInt NfAux, 138 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 139 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 140 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 141 { 142 PetscInt d; 143 for (d = 0; d < dim; ++d) f1[d] = u_x[d]; 144 } 145 void g0_la(PetscInt dim, PetscInt Nf, PetscInt NfAux, 146 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 147 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 148 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 149 { 150 g0[0] = 1.0; 151 } 152 void g3_lu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 153 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 154 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 155 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 156 { 157 PetscInt d; 158 for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0; 159 } 160 161 /* 162 In 2D for Dirichlet conditions with a variable coefficient, we use exact solution: 163 164 u = x^2 + y^2 165 a = 4 166 d_A = 4 167 d_B = sin(2*pi*x[0]) * sin(2*pi*x[1]) 168 169 so that 170 171 -\Delta u + a = -4 + 4 = 0 172 */ 173 PetscErrorCode quadratic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 174 { 175 *u = x[0]*x[0] + x[1]*x[1]; 176 return 0; 177 } 178 PetscErrorCode constant_a_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *a, void *ctx) 179 { 180 *a = 4; 181 return 0; 182 } 183 PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *l, void *ctx) 184 { 185 *l = 0.0; 186 return 0; 187 } 188 189 PetscErrorCode SetupProblem(DM dm, AppCtx *user) 190 { 191 PetscDS ds; 192 DMLabel label; 193 const PetscInt id = 1; 194 195 PetscFunctionBeginUser; 196 PetscCall(DMGetDS(dm, &ds)); 197 PetscCall(PetscDSSetResidual(ds, 0, user->useDualPenalty == PETSC_TRUE ? f0_u_full : f0_u, f1_u)); 198 PetscCall(PetscDSSetResidual(ds, 1, f0_a, NULL)); 199 PetscCall(PetscDSSetResidual(ds, 2, f0_l, f1_l)); 200 PetscCall(PetscDSSetJacobian(ds, 0, 0, user->useDualPenalty == PETSC_TRUE ? g0_uu_full : g0_uu, NULL, NULL, NULL)); 201 PetscCall(PetscDSSetJacobian(ds, 0, 2, NULL, NULL, NULL, g3_ul)); 202 PetscCall(PetscDSSetJacobian(ds, 1, 1, g0_aa, NULL, NULL, NULL)); 203 PetscCall(PetscDSSetJacobian(ds, 1, 2, g0_al, NULL, NULL, NULL)); 204 PetscCall(PetscDSSetJacobian(ds, 2, 1, g0_la, NULL, NULL, NULL)); 205 PetscCall(PetscDSSetJacobian(ds, 2, 0, NULL, NULL, NULL, g3_lu)); 206 207 PetscCall(PetscDSSetExactSolution(ds, 0, quadratic_u_2d, NULL)); 208 PetscCall(PetscDSSetExactSolution(ds, 1, constant_a_2d, NULL)); 209 PetscCall(PetscDSSetExactSolution(ds, 2, zero, NULL)); 210 PetscCall(DMGetLabel(dm, "marker", &label)); 211 PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)()) quadratic_u_2d, NULL, user, NULL)); 212 PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 1, 0, NULL, (void (*)()) constant_a_2d, NULL, user, NULL)); 213 PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 2, 0, NULL, (void (*)()) zero, NULL, user, NULL)); 214 PetscFunctionReturn(0); 215 } 216 217 PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) 218 { 219 DM cdm = dm; 220 const PetscInt dim = 2; 221 PetscFE fe[3]; 222 PetscInt f; 223 MPI_Comm comm; 224 225 PetscFunctionBeginUser; 226 /* Create finite element */ 227 PetscCall(PetscObjectGetComm((PetscObject) dm, &comm)); 228 PetscCall(PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "potential_", -1, &fe[0])); 229 PetscCall(PetscObjectSetName((PetscObject) fe[0], "potential")); 230 PetscCall(PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "charge_", -1, &fe[1])); 231 PetscCall(PetscObjectSetName((PetscObject) fe[1], "charge")); 232 PetscCall(PetscFECopyQuadrature(fe[0], fe[1])); 233 PetscCall(PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "multiplier_", -1, &fe[2])); 234 PetscCall(PetscObjectSetName((PetscObject) fe[2], "multiplier")); 235 PetscCall(PetscFECopyQuadrature(fe[0], fe[2])); 236 /* Set discretization and boundary conditions for each mesh */ 237 for (f = 0; f < 3; ++f) PetscCall(DMSetField(dm, f, NULL, (PetscObject) fe[f])); 238 PetscCall(DMCreateDS(cdm)); 239 PetscCall(SetupProblem(dm, user)); 240 while (cdm) { 241 PetscCall(DMCopyDisc(dm, cdm)); 242 PetscCall(DMGetCoarseDM(cdm, &cdm)); 243 } 244 for (f = 0; f < 3; ++f) PetscCall(PetscFEDestroy(&fe[f])); 245 PetscFunctionReturn(0); 246 } 247 248 int main(int argc, char **argv) 249 { 250 DM dm; 251 SNES snes; 252 Vec u, r; 253 AppCtx user; 254 255 PetscCall(PetscInitialize(&argc, &argv, NULL,help)); 256 PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 257 PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 258 PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 259 PetscCall(SNESSetDM(snes, dm)); 260 PetscCall(SetupDiscretization(dm, &user)); 261 262 PetscCall(DMCreateGlobalVector(dm, &u)); 263 PetscCall(PetscObjectSetName((PetscObject) u, "solution")); 264 PetscCall(VecDuplicate(u, &r)); 265 PetscCall(DMPlexSetSNESLocalFEM(dm,&user,&user,&user)); 266 PetscCall(SNESSetFromOptions(snes)); 267 268 PetscCall(DMSNESCheckFromOptions(snes, u)); 269 if (user.runType == RUN_FULL) { 270 PetscDS ds; 271 PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 272 PetscErrorCode (*initialGuess[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx); 273 PetscReal error; 274 275 PetscCall(DMGetDS(dm, &ds)); 276 PetscCall(PetscDSGetExactSolution(ds, 0, &exactFuncs[0], NULL)); 277 PetscCall(PetscDSGetExactSolution(ds, 1, &exactFuncs[1], NULL)); 278 PetscCall(PetscDSGetExactSolution(ds, 2, &exactFuncs[2], NULL)); 279 initialGuess[0] = zero; 280 initialGuess[1] = zero; 281 initialGuess[2] = zero; 282 PetscCall(DMProjectFunction(dm, 0.0, initialGuess, NULL, INSERT_VALUES, u)); 283 PetscCall(VecViewFromOptions(u, NULL, "-initial_vec_view")); 284 PetscCall(DMComputeL2Diff(dm, 0.0, exactFuncs, NULL, u, &error)); 285 if (error < 1.0e-11) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Initial L_2 Error: < 1.0e-11\n")); 286 else PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Initial L_2 Error: %g\n", (double)error)); 287 PetscCall(SNESSolve(snes, NULL, u)); 288 PetscCall(DMComputeL2Diff(dm, 0.0, exactFuncs, NULL, u, &error)); 289 if (error < 1.0e-11) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Final L_2 Error: < 1.0e-11\n")); 290 else PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Final L_2 Error: %g\n", (double)error)); 291 } 292 PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view")); 293 294 PetscCall(VecDestroy(&u)); 295 PetscCall(VecDestroy(&r)); 296 PetscCall(SNESDestroy(&snes)); 297 PetscCall(DMDestroy(&dm)); 298 PetscCall(PetscFinalize()); 299 return 0; 300 } 301 302 /*TEST 303 304 build: 305 requires: !complex triangle 306 307 test: 308 suffix: 0 309 args: -run_type test -dmsnes_check -potential_petscspace_degree 2 -charge_petscspace_degree 1 -multiplier_petscspace_degree 1 310 311 test: 312 suffix: 1 313 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 314 315 test: 316 suffix: 2 317 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 318 319 TEST*/ 320