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 $a$ and $u$, given by 7 \begin{align} 8 L(u, a, \lambda) = \frac{1}{2} || Qu - d ||^2 + \frac{1}{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 exact control for the reference $a_r$. 14 15 The PDE will be the Laplace equation with homogeneous boundary conditions 16 \begin{align} 17 -nabla \cdot a \nabla u = f 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 } AppCtx; 30 31 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 32 { 33 const char *runTypes[2] = {"full", "test"}; 34 PetscInt run; 35 PetscErrorCode ierr; 36 37 PetscFunctionBeginUser; 38 options->runType = RUN_FULL; 39 40 ierr = PetscOptionsBegin(comm, "", "Inverse Problem Options", "DMPLEX");PetscCall(ierr); 41 run = options->runType; 42 PetscCall(PetscOptionsEList("-run_type", "The run type", "ex1.c", runTypes, 2, runTypes[options->runType], &run, NULL)); 43 options->runType = (RunType) run; 44 ierr = PetscOptionsEnd();PetscCall(ierr); 45 PetscFunctionReturn(0); 46 } 47 48 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 49 { 50 PetscFunctionBeginUser; 51 PetscCall(DMCreate(comm, dm)); 52 PetscCall(DMSetType(*dm, DMPLEX)); 53 PetscCall(DMSetFromOptions(*dm)); 54 PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 55 PetscFunctionReturn(0); 56 } 57 58 /* u - (x^2 + y^2) */ 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 /* a \nabla\lambda */ 67 void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 68 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 69 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 70 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 71 { 72 PetscInt d; 73 for (d = 0; d < dim; ++d) f1[d] = u[1]*u_x[dim*2+d]; 74 } 75 /* I */ 76 void g0_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 77 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 78 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 79 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 80 { 81 g0[0] = 1.0; 82 } 83 /* \nabla */ 84 void g2_ua(PetscInt dim, PetscInt Nf, PetscInt NfAux, 85 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 86 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 87 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]) 88 { 89 PetscInt d; 90 for (d = 0; d < dim; ++d) g2[d] = u_x[dim*2+d]; 91 } 92 /* a */ 93 void g3_ul(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 g3[]) 97 { 98 PetscInt d; 99 for (d = 0; d < dim; ++d) g3[d*dim+d] = u[1]; 100 } 101 /* a - (x + y) */ 102 void f0_a(PetscInt dim, PetscInt Nf, PetscInt NfAux, 103 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 104 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 105 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 106 { 107 f0[0] = u[1] - (x[0] + x[1]); 108 } 109 /* \lambda \nabla u */ 110 void f1_a(PetscInt dim, PetscInt Nf, PetscInt NfAux, 111 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 112 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 113 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 114 { 115 PetscInt d; 116 for (d = 0; d < dim; ++d) f1[d] = u[2]*u_x[d]; 117 } 118 /* I */ 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 /* 6 (x + y) */ 127 void f0_l(PetscInt dim, PetscInt Nf, PetscInt NfAux, 128 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 129 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 130 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 131 { 132 f0[0] = 6.0*(x[0] + x[1]); 133 } 134 /* a \nabla u */ 135 void f1_l(PetscInt dim, PetscInt Nf, PetscInt NfAux, 136 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 137 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 138 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 139 { 140 PetscInt d; 141 for (d = 0; d < dim; ++d) f1[d] = u[1]*u_x[d]; 142 } 143 /* \nabla u */ 144 void g2_la(PetscInt dim, PetscInt Nf, PetscInt NfAux, 145 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 146 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 147 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]) 148 { 149 PetscInt d; 150 for (d = 0; d < dim; ++d) g2[d] = u_x[d]; 151 } 152 /* a */ 153 void g3_lu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 154 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 155 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 156 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 157 { 158 PetscInt d; 159 for (d = 0; d < dim; ++d) g3[d*dim+d] = u[1]; 160 } 161 162 /* 163 In 2D for Dirichlet conditions with a variable coefficient, we use exact solution: 164 165 u = x^2 + y^2 166 f = 6 (x + y) 167 kappa(a) = a = (x + y) 168 169 so that 170 171 -\div \kappa(a) \grad u + f = -6 (x + y) + 6 (x + y) = 0 172 */ 173 PetscErrorCode quadratic_u_2d(PetscInt dim, PetscReal t, 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 linear_a_2d(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar *a, void *ctx) 179 { 180 *a = x[0] + x[1]; 181 return 0; 182 } 183 PetscErrorCode zero(PetscInt dim, PetscReal t, 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, f0_u, f1_u)); 198 PetscCall(PetscDSSetResidual(ds, 1, f0_a, f1_a)); 199 PetscCall(PetscDSSetResidual(ds, 2, f0_l, f1_l)); 200 PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_uu, NULL, NULL, NULL)); 201 PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_ua, NULL)); 202 PetscCall(PetscDSSetJacobian(ds, 0, 2, NULL, NULL, NULL, g3_ul)); 203 PetscCall(PetscDSSetJacobian(ds, 1, 1, g0_aa, NULL, NULL, NULL)); 204 PetscCall(PetscDSSetJacobian(ds, 2, 1, NULL, NULL, g2_la, 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, linear_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 (*)(void)) quadratic_u_2d, NULL, user, NULL)); 212 PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 1, 0, NULL, (void (*)(void)) linear_a_2d, NULL, user, NULL)); 213 PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 2, 0, NULL, (void (*)(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, "conductivity_", -1, &fe[1])); 231 PetscCall(PetscObjectSetName((PetscObject) fe[1], "conductivity")); 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(dm)); 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 t, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 272 PetscErrorCode (*initialGuess[3])(PetscInt dim, PetscReal t, 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", 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", 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 306 307 test: 308 suffix: 0 309 requires: triangle 310 args: -run_type test -dmsnes_check -potential_petscspace_degree 2 -conductivity_petscspace_degree 1 -multiplier_petscspace_degree 2 311 312 test: 313 suffix: 1 314 requires: triangle 315 args: -potential_petscspace_degree 2 -conductivity_petscspace_degree 1 -multiplier_petscspace_degree 2 -snes_monitor -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 -fieldsplit_multiplier_ksp_rtol 1.0e-10 -fieldsplit_multiplier_pc_type lu -sol_vec_view 316 317 TEST*/ 318