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