1c4762a1bSJed Brown static char help[] = "One-Shot Multigrid for Parameter Estimation Problem for the Poisson Equation.\n\ 2c4762a1bSJed Brown Using the Interior Point Method.\n\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown /*F 5c4762a1bSJed Brown We are solving the parameter estimation problem for the Laplacian. We will ask to minimize a Lagrangian 6c4762a1bSJed Brown function over $a$ and $u$, given by 7c4762a1bSJed Brown \begin{align} 8c4762a1bSJed Brown L(u, a, \lambda) = \frac{1}{2} || Qu - d ||^2 + \frac{1}{2} || L (a - a_r) ||^2 + \lambda F(u; a) 9c4762a1bSJed Brown \end{align} 10c4762a1bSJed Brown where $Q$ is a sampling operator, $L$ is a regularization operator, $F$ defines the PDE. 11c4762a1bSJed Brown 12c4762a1bSJed Brown Currently, we have perfect information, meaning $Q = I$, and then we need no regularization, $L = I$. We 13c4762a1bSJed Brown also give the exact control for the reference $a_r$. 14c4762a1bSJed Brown 15c4762a1bSJed Brown The PDE will be the Laplace equation with homogeneous boundary conditions 16c4762a1bSJed Brown \begin{align} 17c4762a1bSJed Brown -nabla \cdot a \nabla u = f 18c4762a1bSJed Brown \end{align} 19c4762a1bSJed Brown 20c4762a1bSJed Brown F*/ 21c4762a1bSJed Brown 22c4762a1bSJed Brown #include <petsc.h> 23c4762a1bSJed Brown #include <petscfe.h> 24c4762a1bSJed Brown 25c4762a1bSJed Brown typedef enum {RUN_FULL, RUN_TEST} RunType; 26c4762a1bSJed Brown 27c4762a1bSJed Brown typedef struct { 28c4762a1bSJed Brown RunType runType; /* Whether to run tests, or solve the full problem */ 29c4762a1bSJed Brown } AppCtx; 30c4762a1bSJed Brown 31c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 32c4762a1bSJed Brown { 33c4762a1bSJed Brown const char *runTypes[2] = {"full", "test"}; 34c4762a1bSJed Brown PetscInt run; 35c4762a1bSJed Brown 36c4762a1bSJed Brown PetscFunctionBeginUser; 37c4762a1bSJed Brown options->runType = RUN_FULL; 38d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Inverse Problem Options", "DMPLEX"); 39c4762a1bSJed Brown run = options->runType; 409566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-run_type", "The run type", "ex1.c", runTypes, 2, runTypes[options->runType], &run, NULL)); 41c4762a1bSJed Brown options->runType = (RunType) run; 42d0609cedSBarry Smith PetscOptionsEnd(); 43c4762a1bSJed Brown PetscFunctionReturn(0); 44c4762a1bSJed Brown } 45c4762a1bSJed Brown 46c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 47c4762a1bSJed Brown { 48c4762a1bSJed Brown PetscFunctionBeginUser; 499566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 509566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 519566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 529566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 53c4762a1bSJed Brown PetscFunctionReturn(0); 54c4762a1bSJed Brown } 55c4762a1bSJed Brown 56c4762a1bSJed Brown /* u - (x^2 + y^2) */ 57c4762a1bSJed Brown void f0_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 58c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 59c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 60c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 61c4762a1bSJed Brown { 62c4762a1bSJed Brown f0[0] = u[0] - (x[0]*x[0] + x[1]*x[1]); 63c4762a1bSJed Brown } 64c4762a1bSJed Brown /* a \nabla\lambda */ 65c4762a1bSJed Brown void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 66c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 67c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 68c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 69c4762a1bSJed Brown { 70c4762a1bSJed Brown PetscInt d; 71c4762a1bSJed Brown for (d = 0; d < dim; ++d) f1[d] = u[1]*u_x[dim*2+d]; 72c4762a1bSJed Brown } 73c4762a1bSJed Brown /* I */ 74c4762a1bSJed Brown void g0_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 75c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 76c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 77c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 78c4762a1bSJed Brown { 79c4762a1bSJed Brown g0[0] = 1.0; 80c4762a1bSJed Brown } 81c4762a1bSJed Brown /* \nabla */ 82c4762a1bSJed Brown void g2_ua(PetscInt dim, PetscInt Nf, PetscInt NfAux, 83c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 84c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 85c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]) 86c4762a1bSJed Brown { 87c4762a1bSJed Brown PetscInt d; 88c4762a1bSJed Brown for (d = 0; d < dim; ++d) g2[d] = u_x[dim*2+d]; 89c4762a1bSJed Brown } 90c4762a1bSJed Brown /* a */ 91c4762a1bSJed Brown void g3_ul(PetscInt dim, PetscInt Nf, PetscInt NfAux, 92c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 93c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 94c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 95c4762a1bSJed Brown { 96c4762a1bSJed Brown PetscInt d; 97c4762a1bSJed Brown for (d = 0; d < dim; ++d) g3[d*dim+d] = u[1]; 98c4762a1bSJed Brown } 99c4762a1bSJed Brown /* a - (x + y) */ 100c4762a1bSJed Brown void f0_a(PetscInt dim, PetscInt Nf, PetscInt NfAux, 101c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 102c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 103c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 104c4762a1bSJed Brown { 105c4762a1bSJed Brown f0[0] = u[1] - (x[0] + x[1]); 106c4762a1bSJed Brown } 107c4762a1bSJed Brown /* \lambda \nabla u */ 108c4762a1bSJed Brown void f1_a(PetscInt dim, PetscInt Nf, PetscInt NfAux, 109c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 110c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 111c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 112c4762a1bSJed Brown { 113c4762a1bSJed Brown PetscInt d; 114c4762a1bSJed Brown for (d = 0; d < dim; ++d) f1[d] = u[2]*u_x[d]; 115c4762a1bSJed Brown } 116c4762a1bSJed Brown /* I */ 117c4762a1bSJed Brown void g0_aa(PetscInt dim, PetscInt Nf, PetscInt NfAux, 118c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 119c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 120c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 121c4762a1bSJed Brown { 122c4762a1bSJed Brown g0[0] = 1.0; 123c4762a1bSJed Brown } 124c4762a1bSJed Brown /* 6 (x + y) */ 125c4762a1bSJed Brown void f0_l(PetscInt dim, PetscInt Nf, PetscInt NfAux, 126c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 127c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 128c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 129c4762a1bSJed Brown { 130c4762a1bSJed Brown f0[0] = 6.0*(x[0] + x[1]); 131c4762a1bSJed Brown } 132c4762a1bSJed Brown /* a \nabla u */ 133c4762a1bSJed Brown void f1_l(PetscInt dim, PetscInt Nf, PetscInt NfAux, 134c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 135c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 136c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 137c4762a1bSJed Brown { 138c4762a1bSJed Brown PetscInt d; 139c4762a1bSJed Brown for (d = 0; d < dim; ++d) f1[d] = u[1]*u_x[d]; 140c4762a1bSJed Brown } 141c4762a1bSJed Brown /* \nabla u */ 142c4762a1bSJed Brown void g2_la(PetscInt dim, PetscInt Nf, PetscInt NfAux, 143c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 144c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 145c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]) 146c4762a1bSJed Brown { 147c4762a1bSJed Brown PetscInt d; 148c4762a1bSJed Brown for (d = 0; d < dim; ++d) g2[d] = u_x[d]; 149c4762a1bSJed Brown } 150c4762a1bSJed Brown /* a */ 151c4762a1bSJed Brown void g3_lu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 152c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 153c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 154c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 155c4762a1bSJed Brown { 156c4762a1bSJed Brown PetscInt d; 157c4762a1bSJed Brown for (d = 0; d < dim; ++d) g3[d*dim+d] = u[1]; 158c4762a1bSJed Brown } 159c4762a1bSJed Brown 160c4762a1bSJed Brown /* 161c4762a1bSJed Brown In 2D for Dirichlet conditions with a variable coefficient, we use exact solution: 162c4762a1bSJed Brown 163c4762a1bSJed Brown u = x^2 + y^2 164c4762a1bSJed Brown f = 6 (x + y) 165c4762a1bSJed Brown kappa(a) = a = (x + y) 166c4762a1bSJed Brown 167c4762a1bSJed Brown so that 168c4762a1bSJed Brown 169c4762a1bSJed Brown -\div \kappa(a) \grad u + f = -6 (x + y) + 6 (x + y) = 0 170c4762a1bSJed Brown */ 171c4762a1bSJed Brown PetscErrorCode quadratic_u_2d(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 172c4762a1bSJed Brown { 173c4762a1bSJed Brown *u = x[0]*x[0] + x[1]*x[1]; 174c4762a1bSJed Brown return 0; 175c4762a1bSJed Brown } 176c4762a1bSJed Brown PetscErrorCode linear_a_2d(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar *a, void *ctx) 177c4762a1bSJed Brown { 178c4762a1bSJed Brown *a = x[0] + x[1]; 179c4762a1bSJed Brown return 0; 180c4762a1bSJed Brown } 181c4762a1bSJed Brown PetscErrorCode zero(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar *l, void *ctx) 182c4762a1bSJed Brown { 183c4762a1bSJed Brown *l = 0.0; 184c4762a1bSJed Brown return 0; 185c4762a1bSJed Brown } 186c4762a1bSJed Brown 187c4762a1bSJed Brown PetscErrorCode SetupProblem(DM dm, AppCtx *user) 188c4762a1bSJed Brown { 18945480ffeSMatthew G. Knepley PetscDS ds; 19045480ffeSMatthew G. Knepley DMLabel label; 191c4762a1bSJed Brown const PetscInt id = 1; 192c4762a1bSJed Brown 193c4762a1bSJed Brown PetscFunctionBeginUser; 1949566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 1959566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_u, f1_u)); 1969566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 1, f0_a, f1_a)); 1979566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 2, f0_l, f1_l)); 1989566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_uu, NULL, NULL, NULL)); 1999566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_ua, NULL)); 2009566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 2, NULL, NULL, NULL, g3_ul)); 2019566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 1, 1, g0_aa, NULL, NULL, NULL)); 2029566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 2, 1, NULL, NULL, g2_la, NULL)); 2039566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 2, 0, NULL, NULL, NULL, g3_lu)); 204c4762a1bSJed Brown 2059566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, quadratic_u_2d, NULL)); 2069566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 1, linear_a_2d, NULL)); 2079566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 2, zero, NULL)); 2089566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 2099566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) quadratic_u_2d, NULL, user, NULL)); 2109566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 1, 0, NULL, (void (*)(void)) linear_a_2d, NULL, user, NULL)); 2119566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 2, 0, NULL, (void (*)(void)) zero, NULL, user, NULL)); 212c4762a1bSJed Brown PetscFunctionReturn(0); 213c4762a1bSJed Brown } 214c4762a1bSJed Brown 215c4762a1bSJed Brown PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) 216c4762a1bSJed Brown { 217c4762a1bSJed Brown DM cdm = dm; 218c4762a1bSJed Brown const PetscInt dim = 2; 219c4762a1bSJed Brown PetscFE fe[3]; 220c4762a1bSJed Brown PetscInt f; 221c4762a1bSJed Brown MPI_Comm comm; 222c4762a1bSJed Brown 223c4762a1bSJed Brown PetscFunctionBeginUser; 224c4762a1bSJed Brown /* Create finite element */ 2259566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject) dm, &comm)); 2269566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "potential_", -1, &fe[0])); 2279566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) fe[0], "potential")); 2289566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "conductivity_", -1, &fe[1])); 2299566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) fe[1], "conductivity")); 2309566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fe[0], fe[1])); 2319566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "multiplier_", -1, &fe[2])); 2329566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) fe[2], "multiplier")); 2339566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fe[0], fe[2])); 234c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 2359566063dSJacob Faibussowitsch for (f = 0; f < 3; ++f) PetscCall(DMSetField(dm, f, NULL, (PetscObject) fe[f])); 2369566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 2379566063dSJacob Faibussowitsch PetscCall(SetupProblem(dm, user)); 238c4762a1bSJed Brown while (cdm) { 2399566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm)); 2409566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm)); 241c4762a1bSJed Brown } 2429566063dSJacob Faibussowitsch for (f = 0; f < 3; ++f) PetscCall(PetscFEDestroy(&fe[f])); 243c4762a1bSJed Brown PetscFunctionReturn(0); 244c4762a1bSJed Brown } 245c4762a1bSJed Brown 246c4762a1bSJed Brown int main(int argc, char **argv) 247c4762a1bSJed Brown { 248c4762a1bSJed Brown DM dm; 249c4762a1bSJed Brown SNES snes; 250c4762a1bSJed Brown Vec u, r; 251c4762a1bSJed Brown AppCtx user; 252c4762a1bSJed Brown 253*327415f7SBarry Smith PetscFunctionBeginUser; 2549566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL,help)); 2559566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 2569566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 2579566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 2589566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, dm)); 2599566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dm, &user)); 260c4762a1bSJed Brown 2619566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u)); 2629566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) u, "solution")); 2639566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &r)); 2649566063dSJacob Faibussowitsch PetscCall(DMPlexSetSNESLocalFEM(dm,&user,&user,&user)); 2659566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 266c4762a1bSJed Brown 2679566063dSJacob Faibussowitsch PetscCall(DMSNESCheckFromOptions(snes, u)); 268c4762a1bSJed Brown if (user.runType == RUN_FULL) { 269348a1646SMatthew G. Knepley PetscDS ds; 270348a1646SMatthew G. Knepley PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 271c4762a1bSJed Brown PetscErrorCode (*initialGuess[3])(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx); 272c4762a1bSJed Brown PetscReal error; 273c4762a1bSJed Brown 2749566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 2759566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 0, &exactFuncs[0], NULL)); 2769566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 1, &exactFuncs[1], NULL)); 2779566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 2, &exactFuncs[2], NULL)); 278c4762a1bSJed Brown initialGuess[0] = zero; 279c4762a1bSJed Brown initialGuess[1] = zero; 280c4762a1bSJed Brown initialGuess[2] = zero; 2819566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dm, 0.0, initialGuess, NULL, INSERT_VALUES, u)); 2829566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-initial_vec_view")); 2839566063dSJacob Faibussowitsch PetscCall(DMComputeL2Diff(dm, 0.0, exactFuncs, NULL, u, &error)); 2849566063dSJacob Faibussowitsch if (error < 1.0e-11) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Initial L_2 Error: < 1.0e-11\n")); 28563a3b9bcSJacob Faibussowitsch else PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Initial L_2 Error: %g\n", (double)error)); 2869566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, u)); 2879566063dSJacob Faibussowitsch PetscCall(DMComputeL2Diff(dm, 0.0, exactFuncs, NULL, u, &error)); 2889566063dSJacob Faibussowitsch if (error < 1.0e-11) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Final L_2 Error: < 1.0e-11\n")); 28963a3b9bcSJacob Faibussowitsch else PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Final L_2 Error: %g\n", (double)error)); 290c4762a1bSJed Brown } 2919566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view")); 292c4762a1bSJed Brown 2939566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 2949566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 2959566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 2969566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 2979566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 298b122ec5aSJacob Faibussowitsch return 0; 299c4762a1bSJed Brown } 300c4762a1bSJed Brown 301c4762a1bSJed Brown /*TEST 302c4762a1bSJed Brown 303c4762a1bSJed Brown build: 304c4762a1bSJed Brown requires: !complex 305c4762a1bSJed Brown 306c4762a1bSJed Brown test: 307c4762a1bSJed Brown suffix: 0 308c4762a1bSJed Brown requires: triangle 309c4762a1bSJed Brown args: -run_type test -dmsnes_check -potential_petscspace_degree 2 -conductivity_petscspace_degree 1 -multiplier_petscspace_degree 2 310c4762a1bSJed Brown 311c4762a1bSJed Brown test: 312c4762a1bSJed Brown suffix: 1 313c4762a1bSJed Brown requires: triangle 314c4762a1bSJed Brown 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 315c4762a1bSJed Brown 316c4762a1bSJed Brown TEST*/ 317