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