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 $y$ and $u$, given by 7c4762a1bSJed Brown \begin{align} 8c4762a1bSJed Brown 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) 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 null vector for the reference control $a_r$. Right now $\beta = 1$. 14c4762a1bSJed Brown 15c4762a1bSJed Brown The PDE will be the Laplace equation with homogeneous boundary conditions 16c4762a1bSJed Brown \begin{align} 17c4762a1bSJed Brown -Delta u = a 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 PetscBool useDualPenalty; /* Penalize deviation from both goals */ 30c4762a1bSJed Brown } AppCtx; 31c4762a1bSJed Brown 32c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 33c4762a1bSJed Brown { 34c4762a1bSJed Brown const char *runTypes[2] = {"full", "test"}; 35c4762a1bSJed Brown PetscInt run; 36c4762a1bSJed Brown 37c4762a1bSJed Brown PetscFunctionBeginUser; 38c4762a1bSJed Brown options->runType = RUN_FULL; 39c4762a1bSJed Brown options->useDualPenalty = PETSC_FALSE; 40d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Inverse Problem Options", "DMPLEX"); 41c4762a1bSJed Brown run = options->runType; 429566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-run_type", "The run type", "ex2.c", runTypes, 2, runTypes[options->runType], &run, NULL)); 43c4762a1bSJed Brown options->runType = (RunType) run; 449566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-use_dual_penalty", "Penalize deviation from both goals", "ex2.c", options->useDualPenalty, &options->useDualPenalty, NULL)); 45d0609cedSBarry Smith PetscOptionsEnd(); 46c4762a1bSJed Brown PetscFunctionReturn(0); 47c4762a1bSJed Brown } 48c4762a1bSJed Brown 49c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 50c4762a1bSJed Brown { 51c4762a1bSJed Brown PetscFunctionBeginUser; 529566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 539566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 549566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 559566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 56c4762a1bSJed Brown PetscFunctionReturn(0); 57c4762a1bSJed Brown } 58c4762a1bSJed Brown 59c4762a1bSJed Brown void f0_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 60c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 61c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 62c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 63c4762a1bSJed Brown { 64c4762a1bSJed Brown f0[0] = (u[0] - (x[0]*x[0] + x[1]*x[1])); 65c4762a1bSJed Brown } 66c4762a1bSJed Brown void f0_u_full(PetscInt dim, PetscInt Nf, PetscInt NfAux, 67c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 68c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 69c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 70c4762a1bSJed Brown { 71c4762a1bSJed Brown 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]))) + 72c4762a1bSJed Brown 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]))); 73c4762a1bSJed Brown } 74c4762a1bSJed Brown void f1_u(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, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 78c4762a1bSJed Brown { 79c4762a1bSJed Brown PetscInt d; 80c4762a1bSJed Brown for (d = 0; d < dim; ++d) f1[d] = u_x[dim*2+d]; 81c4762a1bSJed Brown } 82c4762a1bSJed Brown void g0_uu(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 g0[]) 86c4762a1bSJed Brown { 87c4762a1bSJed Brown g0[0] = 1.0; 88c4762a1bSJed Brown } 89c4762a1bSJed Brown void g0_uu_full(PetscInt dim, PetscInt Nf, PetscInt NfAux, 90c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 91c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 92c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 93c4762a1bSJed Brown { 94c4762a1bSJed Brown g0[0] = PetscSqr(u[0] - sin(2.0*PETSC_PI*x[0]) * sin(2.0*PETSC_PI*x[1])) 95c4762a1bSJed Brown + PetscSqr(u[0] - (x[0]*x[0] + x[1]*x[1])) 96c4762a1bSJed Brown - 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] 97c4762a1bSJed Brown + 4.0*(x[0]*x[0] + x[1]*x[1])*(sin(2.0*PETSC_PI*x[0]) * sin(2.0*PETSC_PI*x[1])); 98c4762a1bSJed Brown } 99c4762a1bSJed Brown void g3_ul(PetscInt dim, PetscInt Nf, PetscInt NfAux, 100c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 101c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 102c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 103c4762a1bSJed Brown { 104c4762a1bSJed Brown PetscInt d; 105c4762a1bSJed Brown for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0; 106c4762a1bSJed Brown } 107c4762a1bSJed Brown 108c4762a1bSJed Brown void f0_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 f0[]) 112c4762a1bSJed Brown { 113c4762a1bSJed Brown f0[0] = u[1] - 4.0 /* 0.0 */ + u[2]; 114c4762a1bSJed Brown } 115c4762a1bSJed Brown void g0_aa(PetscInt dim, PetscInt Nf, PetscInt NfAux, 116c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 117c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 118c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 119c4762a1bSJed Brown { 120c4762a1bSJed Brown g0[0] = 1.0; 121c4762a1bSJed Brown } 122c4762a1bSJed Brown void g0_al(PetscInt dim, PetscInt Nf, PetscInt NfAux, 123c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 124c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 125c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 126c4762a1bSJed Brown { 127c4762a1bSJed Brown g0[0] = 1.0; 128c4762a1bSJed Brown } 129c4762a1bSJed Brown 130c4762a1bSJed Brown void f0_l(PetscInt dim, PetscInt Nf, PetscInt NfAux, 131c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 132c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 133c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 134c4762a1bSJed Brown { 135c4762a1bSJed Brown f0[0] = u[1]; 136c4762a1bSJed Brown } 137c4762a1bSJed Brown void f1_l(PetscInt dim, PetscInt Nf, PetscInt NfAux, 138c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 139c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 140c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 141c4762a1bSJed Brown { 142c4762a1bSJed Brown PetscInt d; 143c4762a1bSJed Brown for (d = 0; d < dim; ++d) f1[d] = u_x[d]; 144c4762a1bSJed Brown } 145c4762a1bSJed Brown void g0_la(PetscInt dim, PetscInt Nf, PetscInt NfAux, 146c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 147c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 148c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 149c4762a1bSJed Brown { 150c4762a1bSJed Brown g0[0] = 1.0; 151c4762a1bSJed Brown } 152c4762a1bSJed Brown void g3_lu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 153c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 154c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 155c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 156c4762a1bSJed Brown { 157c4762a1bSJed Brown PetscInt d; 158c4762a1bSJed Brown for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0; 159c4762a1bSJed Brown } 160c4762a1bSJed Brown 161c4762a1bSJed Brown /* 162c4762a1bSJed Brown In 2D for Dirichlet conditions with a variable coefficient, we use exact solution: 163c4762a1bSJed Brown 164c4762a1bSJed Brown u = x^2 + y^2 165c4762a1bSJed Brown a = 4 166c4762a1bSJed Brown d_A = 4 167c4762a1bSJed Brown d_B = sin(2*pi*x[0]) * sin(2*pi*x[1]) 168c4762a1bSJed Brown 169c4762a1bSJed Brown so that 170c4762a1bSJed Brown 171c4762a1bSJed Brown -\Delta u + a = -4 + 4 = 0 172c4762a1bSJed Brown */ 173c4762a1bSJed Brown PetscErrorCode quadratic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 174c4762a1bSJed Brown { 175c4762a1bSJed Brown *u = x[0]*x[0] + x[1]*x[1]; 176c4762a1bSJed Brown return 0; 177c4762a1bSJed Brown } 178c4762a1bSJed Brown PetscErrorCode constant_a_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *a, void *ctx) 179c4762a1bSJed Brown { 180c4762a1bSJed Brown *a = 4; 181c4762a1bSJed Brown return 0; 182c4762a1bSJed Brown } 183c4762a1bSJed Brown PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *l, void *ctx) 184c4762a1bSJed Brown { 185c4762a1bSJed Brown *l = 0.0; 186c4762a1bSJed Brown return 0; 187c4762a1bSJed Brown } 188c4762a1bSJed Brown 189c4762a1bSJed Brown PetscErrorCode SetupProblem(DM dm, AppCtx *user) 190c4762a1bSJed Brown { 19145480ffeSMatthew G. Knepley PetscDS ds; 19245480ffeSMatthew G. Knepley DMLabel label; 193c4762a1bSJed Brown const PetscInt id = 1; 194c4762a1bSJed Brown 195c4762a1bSJed Brown PetscFunctionBeginUser; 1969566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 1979566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, user->useDualPenalty == PETSC_TRUE ? f0_u_full : f0_u, f1_u)); 1989566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 1, f0_a, NULL)); 1999566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 2, f0_l, f1_l)); 2009566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, user->useDualPenalty == PETSC_TRUE ? g0_uu_full : g0_uu, NULL, NULL, NULL)); 2019566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 2, NULL, NULL, NULL, g3_ul)); 2029566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 1, 1, g0_aa, NULL, NULL, NULL)); 2039566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 1, 2, g0_al, NULL, NULL, NULL)); 2049566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 2, 1, g0_la, NULL, NULL, NULL)); 2059566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 2, 0, NULL, NULL, NULL, g3_lu)); 206c4762a1bSJed Brown 2079566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, quadratic_u_2d, NULL)); 2089566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 1, constant_a_2d, NULL)); 2099566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 2, zero, NULL)); 2109566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 2119566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)()) quadratic_u_2d, NULL, user, NULL)); 2129566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 1, 0, NULL, (void (*)()) constant_a_2d, NULL, user, NULL)); 2139566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 2, 0, NULL, (void (*)()) zero, NULL, user, NULL)); 214c4762a1bSJed Brown PetscFunctionReturn(0); 215c4762a1bSJed Brown } 216c4762a1bSJed Brown 217c4762a1bSJed Brown PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) 218c4762a1bSJed Brown { 219c4762a1bSJed Brown DM cdm = dm; 220c4762a1bSJed Brown const PetscInt dim = 2; 221c4762a1bSJed Brown PetscFE fe[3]; 222c4762a1bSJed Brown PetscInt f; 223c4762a1bSJed Brown MPI_Comm comm; 224c4762a1bSJed Brown 225c4762a1bSJed Brown PetscFunctionBeginUser; 226c4762a1bSJed Brown /* Create finite element */ 2279566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject) dm, &comm)); 2289566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "potential_", -1, &fe[0])); 2299566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) fe[0], "potential")); 2309566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "charge_", -1, &fe[1])); 2319566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) fe[1], "charge")); 2329566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fe[0], fe[1])); 2339566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "multiplier_", -1, &fe[2])); 2349566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) fe[2], "multiplier")); 2359566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fe[0], fe[2])); 236c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 2379566063dSJacob Faibussowitsch for (f = 0; f < 3; ++f) PetscCall(DMSetField(dm, f, NULL, (PetscObject) fe[f])); 2389566063dSJacob Faibussowitsch PetscCall(DMCreateDS(cdm)); 2399566063dSJacob Faibussowitsch PetscCall(SetupProblem(dm, user)); 240c4762a1bSJed Brown while (cdm) { 2419566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm)); 2429566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm)); 243c4762a1bSJed Brown } 2449566063dSJacob Faibussowitsch for (f = 0; f < 3; ++f) PetscCall(PetscFEDestroy(&fe[f])); 245c4762a1bSJed Brown PetscFunctionReturn(0); 246c4762a1bSJed Brown } 247c4762a1bSJed Brown 248c4762a1bSJed Brown int main(int argc, char **argv) 249c4762a1bSJed Brown { 250c4762a1bSJed Brown DM dm; 251c4762a1bSJed Brown SNES snes; 252c4762a1bSJed Brown Vec u, r; 253c4762a1bSJed Brown AppCtx user; 254c4762a1bSJed Brown 255*327415f7SBarry Smith PetscFunctionBeginUser; 2569566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL,help)); 2579566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 2589566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 2599566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 2609566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, dm)); 2619566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dm, &user)); 262c4762a1bSJed Brown 2639566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u)); 2649566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) u, "solution")); 2659566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &r)); 2669566063dSJacob Faibussowitsch PetscCall(DMPlexSetSNESLocalFEM(dm,&user,&user,&user)); 2679566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 268c4762a1bSJed Brown 2699566063dSJacob Faibussowitsch PetscCall(DMSNESCheckFromOptions(snes, u)); 270c4762a1bSJed Brown if (user.runType == RUN_FULL) { 271348a1646SMatthew G. Knepley PetscDS ds; 272348a1646SMatthew G. Knepley PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 273c4762a1bSJed Brown PetscErrorCode (*initialGuess[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx); 274c4762a1bSJed Brown PetscReal error; 275c4762a1bSJed Brown 2769566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 2779566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 0, &exactFuncs[0], NULL)); 2789566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 1, &exactFuncs[1], NULL)); 2799566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 2, &exactFuncs[2], NULL)); 280c4762a1bSJed Brown initialGuess[0] = zero; 281c4762a1bSJed Brown initialGuess[1] = zero; 282c4762a1bSJed Brown initialGuess[2] = zero; 2839566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dm, 0.0, initialGuess, NULL, INSERT_VALUES, u)); 2849566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-initial_vec_view")); 2859566063dSJacob Faibussowitsch PetscCall(DMComputeL2Diff(dm, 0.0, exactFuncs, NULL, u, &error)); 2869566063dSJacob Faibussowitsch if (error < 1.0e-11) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Initial L_2 Error: < 1.0e-11\n")); 28763a3b9bcSJacob Faibussowitsch else PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Initial L_2 Error: %g\n", (double)error)); 2889566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, u)); 2899566063dSJacob Faibussowitsch PetscCall(DMComputeL2Diff(dm, 0.0, exactFuncs, NULL, u, &error)); 2909566063dSJacob Faibussowitsch if (error < 1.0e-11) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Final L_2 Error: < 1.0e-11\n")); 29163a3b9bcSJacob Faibussowitsch else PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Final L_2 Error: %g\n", (double)error)); 292c4762a1bSJed Brown } 2939566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view")); 294c4762a1bSJed Brown 2959566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 2969566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 2979566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 2989566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 2999566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 300b122ec5aSJacob Faibussowitsch return 0; 301c4762a1bSJed Brown } 302c4762a1bSJed Brown 303c4762a1bSJed Brown /*TEST 304c4762a1bSJed Brown 305c4762a1bSJed Brown build: 306c4762a1bSJed Brown requires: !complex triangle 307c4762a1bSJed Brown 308c4762a1bSJed Brown test: 309c4762a1bSJed Brown suffix: 0 310c4762a1bSJed Brown args: -run_type test -dmsnes_check -potential_petscspace_degree 2 -charge_petscspace_degree 1 -multiplier_petscspace_degree 1 311c4762a1bSJed Brown 312c4762a1bSJed Brown test: 313c4762a1bSJed Brown suffix: 1 314c4762a1bSJed Brown 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 315c4762a1bSJed Brown 316c4762a1bSJed Brown test: 317c4762a1bSJed Brown suffix: 2 318c4762a1bSJed Brown 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 319c4762a1bSJed Brown 320c4762a1bSJed Brown TEST*/ 321