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 PetscErrorCode (**exactFuncs)(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 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 41 ierr = PetscOptionsBegin(comm, "", "Inverse Problem Options", "DMPLEX");CHKERRQ(ierr); 42 run = options->runType; 43 ierr = PetscOptionsEList("-run_type", "The run type", "ex1.c", runTypes, 2, runTypes[options->runType], &run, NULL);CHKERRQ(ierr); 44 options->runType = (RunType) run; 45 ierr = PetscOptionsEnd();CHKERRQ(ierr); 46 PetscFunctionReturn(0); 47 } 48 49 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 50 { 51 DM distributedMesh = NULL; 52 PetscErrorCode ierr; 53 54 PetscFunctionBeginUser; 55 ierr = DMPlexCreateBoxMesh(comm, 2, PETSC_TRUE, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr); 56 ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr); 57 ierr = DMPlexDistribute(*dm, 0, NULL, &distributedMesh);CHKERRQ(ierr); 58 if (distributedMesh) { 59 ierr = DMDestroy(dm);CHKERRQ(ierr); 60 *dm = distributedMesh; 61 } 62 ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 63 ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 64 PetscFunctionReturn(0); 65 } 66 67 /* u - (x^2 + y^2) */ 68 void f0_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 69 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 70 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 71 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 72 { 73 f0[0] = u[0] - (x[0]*x[0] + x[1]*x[1]); 74 } 75 /* a \nabla\lambda */ 76 void f1_u(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, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 80 { 81 PetscInt d; 82 for (d = 0; d < dim; ++d) f1[d] = u[1]*u_x[dim*2+d]; 83 } 84 /* I */ 85 void g0_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 86 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 87 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 88 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 89 { 90 g0[0] = 1.0; 91 } 92 /* \nabla */ 93 void g2_ua(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 g2[]) 97 { 98 PetscInt d; 99 for (d = 0; d < dim; ++d) g2[d] = u_x[dim*2+d]; 100 } 101 /* a */ 102 void g3_ul(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, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 106 { 107 PetscInt d; 108 for (d = 0; d < dim; ++d) g3[d*dim+d] = u[1]; 109 } 110 /* a - (x + y) */ 111 void f0_a(PetscInt dim, PetscInt Nf, PetscInt NfAux, 112 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 113 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 114 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 115 { 116 f0[0] = u[1] - (x[0] + x[1]); 117 } 118 /* \lambda \nabla u */ 119 void f1_a(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, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 123 { 124 PetscInt d; 125 for (d = 0; d < dim; ++d) f1[d] = u[2]*u_x[d]; 126 } 127 /* I */ 128 void g0_aa(PetscInt dim, PetscInt Nf, PetscInt NfAux, 129 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 130 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 131 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 132 { 133 g0[0] = 1.0; 134 } 135 /* 6 (x + y) */ 136 void f0_l(PetscInt dim, PetscInt Nf, PetscInt NfAux, 137 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 138 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 139 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 140 { 141 f0[0] = 6.0*(x[0] + x[1]); 142 } 143 /* a \nabla u */ 144 void f1_l(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, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 148 { 149 PetscInt d; 150 for (d = 0; d < dim; ++d) f1[d] = u[1]*u_x[d]; 151 } 152 /* \nabla u */ 153 void g2_la(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 g2[]) 157 { 158 PetscInt d; 159 for (d = 0; d < dim; ++d) g2[d] = u_x[d]; 160 } 161 /* a */ 162 void g3_lu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 163 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 164 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 165 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 166 { 167 PetscInt d; 168 for (d = 0; d < dim; ++d) g3[d*dim+d] = u[1]; 169 } 170 171 /* 172 In 2D for Dirichlet conditions with a variable coefficient, we use exact solution: 173 174 u = x^2 + y^2 175 f = 6 (x + y) 176 kappa(a) = a = (x + y) 177 178 so that 179 180 -\div \kappa(a) \grad u + f = -6 (x + y) + 6 (x + y) = 0 181 */ 182 PetscErrorCode quadratic_u_2d(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 183 { 184 *u = x[0]*x[0] + x[1]*x[1]; 185 return 0; 186 } 187 PetscErrorCode linear_a_2d(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar *a, void *ctx) 188 { 189 *a = x[0] + x[1]; 190 return 0; 191 } 192 PetscErrorCode zero(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar *l, void *ctx) 193 { 194 *l = 0.0; 195 return 0; 196 } 197 198 PetscErrorCode SetupProblem(DM dm, AppCtx *user) 199 { 200 PetscDS prob; 201 const PetscInt id = 1; 202 PetscErrorCode ierr; 203 204 PetscFunctionBeginUser; 205 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 206 ierr = PetscDSSetResidual(prob, 0, f0_u, f1_u);CHKERRQ(ierr); 207 ierr = PetscDSSetResidual(prob, 1, f0_a, f1_a);CHKERRQ(ierr); 208 ierr = PetscDSSetResidual(prob, 2, f0_l, f1_l);CHKERRQ(ierr); 209 ierr = PetscDSSetJacobian(prob, 0, 0, g0_uu, NULL, NULL, NULL);CHKERRQ(ierr); 210 ierr = PetscDSSetJacobian(prob, 0, 1, NULL, NULL, g2_ua, NULL);CHKERRQ(ierr); 211 ierr = PetscDSSetJacobian(prob, 0, 2, NULL, NULL, NULL, g3_ul);CHKERRQ(ierr); 212 ierr = PetscDSSetJacobian(prob, 1, 1, g0_aa, NULL, NULL, NULL);CHKERRQ(ierr); 213 ierr = PetscDSSetJacobian(prob, 2, 1, NULL, NULL, g2_la, NULL);CHKERRQ(ierr); 214 ierr = PetscDSSetJacobian(prob, 2, 0, NULL, NULL, NULL, g3_lu);CHKERRQ(ierr); 215 216 user->exactFuncs[0] = quadratic_u_2d; 217 user->exactFuncs[1] = linear_a_2d; 218 user->exactFuncs[2] = zero; 219 ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "marker", 0, 0, NULL, (void (*)(void)) user->exactFuncs[0], 1, &id, user);CHKERRQ(ierr); 220 ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "marker", 1, 0, NULL, (void (*)(void)) user->exactFuncs[1], 1, &id, user);CHKERRQ(ierr); 221 ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "marker", 2, 0, NULL, (void (*)(void)) user->exactFuncs[2], 1, &id, user);CHKERRQ(ierr); 222 PetscFunctionReturn(0); 223 } 224 225 PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) 226 { 227 DM cdm = dm; 228 const PetscInt dim = 2; 229 PetscFE fe[3]; 230 PetscInt f; 231 MPI_Comm comm; 232 PetscErrorCode ierr; 233 234 PetscFunctionBeginUser; 235 /* Create finite element */ 236 ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 237 ierr = PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "potential_", -1, &fe[0]);CHKERRQ(ierr); 238 ierr = PetscObjectSetName((PetscObject) fe[0], "potential");CHKERRQ(ierr); 239 ierr = PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "conductivity_", -1, &fe[1]);CHKERRQ(ierr); 240 ierr = PetscObjectSetName((PetscObject) fe[1], "conductivity");CHKERRQ(ierr); 241 ierr = PetscFECopyQuadrature(fe[0], fe[1]);CHKERRQ(ierr); 242 ierr = PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "multiplier_", -1, &fe[2]);CHKERRQ(ierr); 243 ierr = PetscObjectSetName((PetscObject) fe[2], "multiplier");CHKERRQ(ierr); 244 ierr = PetscFECopyQuadrature(fe[0], fe[2]);CHKERRQ(ierr); 245 /* Set discretization and boundary conditions for each mesh */ 246 for (f = 0; f < 3; ++f) {ierr = DMSetField(dm, f, NULL, (PetscObject) fe[f]);CHKERRQ(ierr);} 247 ierr = DMCreateDS(dm);CHKERRQ(ierr); 248 ierr = SetupProblem(dm, user);CHKERRQ(ierr); 249 while (cdm) { 250 ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr); 251 ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr); 252 } 253 for (f = 0; f < 3; ++f) {ierr = PetscFEDestroy(&fe[f]);CHKERRQ(ierr);} 254 PetscFunctionReturn(0); 255 } 256 257 int main(int argc, char **argv) 258 { 259 DM dm; 260 SNES snes; 261 Vec u, r; 262 AppCtx user; 263 PetscErrorCode ierr; 264 265 ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 266 ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr); 267 ierr = SNESCreate(PETSC_COMM_WORLD, &snes);CHKERRQ(ierr); 268 ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr); 269 ierr = SNESSetDM(snes, dm);CHKERRQ(ierr); 270 271 ierr = PetscMalloc(3 * sizeof(void (*)()), &user.exactFuncs);CHKERRQ(ierr); 272 ierr = SetupDiscretization(dm, &user);CHKERRQ(ierr); 273 274 ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr); 275 ierr = PetscObjectSetName((PetscObject) u, "solution");CHKERRQ(ierr); 276 ierr = VecDuplicate(u, &r);CHKERRQ(ierr); 277 ierr = DMPlexSetSNESLocalFEM(dm,&user,&user,&user);CHKERRQ(ierr); 278 ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 279 280 ierr = DMProjectFunction(dm, 0.0, user.exactFuncs, NULL, INSERT_ALL_VALUES, u);CHKERRQ(ierr); 281 ierr = DMSNESCheckFromOptions(snes, u, user.exactFuncs, NULL);CHKERRQ(ierr); 282 if (user.runType == RUN_FULL) { 283 PetscErrorCode (*initialGuess[3])(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx); 284 PetscReal error; 285 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, user.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, user.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 = PetscFree(user.exactFuncs);CHKERRQ(ierr); 306 ierr = PetscFinalize(); 307 return ierr; 308 } 309 310 /*TEST 311 312 build: 313 requires: !complex 314 315 test: 316 suffix: 0 317 requires: triangle 318 args: -run_type test -dmsnes_check -potential_petscspace_degree 2 -conductivity_petscspace_degree 1 -multiplier_petscspace_degree 2 319 320 test: 321 suffix: 1 322 requires: triangle 323 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 324 325 TEST*/ 326