1*736140d5SStefano Zampini static const char help[] = "\ 2*736140d5SStefano Zampini Minimum surface area problem in 2D using DMPLEX.\n\ 3*736140d5SStefano Zampini It solves the unconstrained minimization problem \n\ 4*736140d5SStefano Zampini \n\ 5*736140d5SStefano Zampini argmin \\int_\\Omega (1 + ||u||^2), u = g on \\partial\\Omega.\n\ 6*736140d5SStefano Zampini \n\ 7*736140d5SStefano Zampini This example shows how to setup an optimization problem using DMPLEX FEM routines.\n\ 8*736140d5SStefano Zampini It supports nonlinear domain decomposition and multilevel solvers.\n\ 9*736140d5SStefano Zampini \n\n"; 10*736140d5SStefano Zampini 11*736140d5SStefano Zampini #include <petscdmplex.h> 12*736140d5SStefano Zampini #include <petscsnes.h> 13*736140d5SStefano Zampini #include <petscds.h> 14*736140d5SStefano Zampini 15*736140d5SStefano Zampini /* The pointwise evaluation routine for the objective function */ 16*736140d5SStefano Zampini static void objective_2d(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar obj[]) 17*736140d5SStefano Zampini { 18*736140d5SStefano Zampini const PetscScalar ux2 = PetscSqr(u_x[0]); 19*736140d5SStefano Zampini const PetscScalar uy2 = PetscSqr(u_x[1]); 20*736140d5SStefano Zampini 21*736140d5SStefano Zampini obj[0] = PetscSqrtReal(PetscAbsScalar(1.0 + ux2 + uy2)); 22*736140d5SStefano Zampini } 23*736140d5SStefano Zampini 24*736140d5SStefano Zampini /* The pointwise evaluation routine for the gradient wrt the gradients of the FEM basis */ 25*736140d5SStefano Zampini static void gradient_1_2d(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f[]) 26*736140d5SStefano Zampini { 27*736140d5SStefano Zampini const PetscScalar ux2 = PetscSqr(u_x[0]); 28*736140d5SStefano Zampini const PetscScalar uy2 = PetscSqr(u_x[1]); 29*736140d5SStefano Zampini const PetscScalar v = 1. / PetscSqrtReal(PetscAbsScalar(1.0 + ux2 + uy2)); 30*736140d5SStefano Zampini 31*736140d5SStefano Zampini f[0] = v * u_x[0]; 32*736140d5SStefano Zampini f[1] = v * u_x[1]; 33*736140d5SStefano Zampini } 34*736140d5SStefano Zampini 35*736140d5SStefano Zampini /* The pointwise evaluation routine for the hessian wrt the gradients of the FEM basis */ 36*736140d5SStefano Zampini static void hessian_11_2d(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar jac[]) 37*736140d5SStefano Zampini { 38*736140d5SStefano Zampini const PetscScalar ux2 = PetscSqr(u_x[0]); 39*736140d5SStefano Zampini const PetscScalar uy2 = PetscSqr(u_x[1]); 40*736140d5SStefano Zampini const PetscScalar uxy = u_x[0] * u_x[1]; 41*736140d5SStefano Zampini const PetscScalar v1 = 1. / PetscSqrtReal(PetscAbsScalar(1.0 + ux2 + uy2)); 42*736140d5SStefano Zampini const PetscScalar v2 = v1 / (1.0 + ux2 + uy2); 43*736140d5SStefano Zampini 44*736140d5SStefano Zampini jac[0] = v1 - v2 * ux2; 45*736140d5SStefano Zampini jac[1] = -v2 * uxy; 46*736140d5SStefano Zampini jac[2] = -v2 * uxy; 47*736140d5SStefano Zampini jac[3] = v1 - v2 * uy2; 48*736140d5SStefano Zampini } 49*736140d5SStefano Zampini 50*736140d5SStefano Zampini /* The boundary conditions we impose */ 51*736140d5SStefano Zampini static PetscErrorCode sins_2d(PetscInt dim, PetscReal time, const PetscReal xx[], PetscInt Nc, PetscScalar *u, void *ctx) 52*736140d5SStefano Zampini { 53*736140d5SStefano Zampini PetscFunctionBeginUser; 54*736140d5SStefano Zampini const PetscReal pi2 = PETSC_PI * 2; 55*736140d5SStefano Zampini const PetscReal x = xx[0]; 56*736140d5SStefano Zampini const PetscReal y = xx[1]; 57*736140d5SStefano Zampini 58*736140d5SStefano Zampini *u = (y - 0.5) * PetscSinReal(pi2 * x) + (x - 0.5) * PetscSinReal(pi2 * y); 59*736140d5SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 60*736140d5SStefano Zampini } 61*736140d5SStefano Zampini 62*736140d5SStefano Zampini static PetscErrorCode CreateBCLabel(DM dm, const char name[]) 63*736140d5SStefano Zampini { 64*736140d5SStefano Zampini DM plex; 65*736140d5SStefano Zampini DMLabel label; 66*736140d5SStefano Zampini 67*736140d5SStefano Zampini PetscFunctionBeginUser; 68*736140d5SStefano Zampini PetscCall(DMCreateLabel(dm, name)); 69*736140d5SStefano Zampini PetscCall(DMGetLabel(dm, name, &label)); 70*736140d5SStefano Zampini PetscCall(DMConvert(dm, DMPLEX, &plex)); 71*736140d5SStefano Zampini PetscCall(DMPlexMarkBoundaryFaces(plex, 1, label)); 72*736140d5SStefano Zampini PetscCall(DMDestroy(&plex)); 73*736140d5SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 74*736140d5SStefano Zampini } 75*736140d5SStefano Zampini 76*736140d5SStefano Zampini static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm) 77*736140d5SStefano Zampini { 78*736140d5SStefano Zampini PetscFunctionBeginUser; 79*736140d5SStefano Zampini PetscCall(DMCreate(comm, dm)); 80*736140d5SStefano Zampini PetscCall(DMSetType(*dm, DMPLEX)); 81*736140d5SStefano Zampini PetscCall(DMSetFromOptions(*dm)); 82*736140d5SStefano Zampini PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 83*736140d5SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 84*736140d5SStefano Zampini } 85*736140d5SStefano Zampini 86*736140d5SStefano Zampini static PetscErrorCode SetupProblem(DM dm) 87*736140d5SStefano Zampini { 88*736140d5SStefano Zampini PetscDS ds; 89*736140d5SStefano Zampini DMLabel label; 90*736140d5SStefano Zampini PetscInt dim; 91*736140d5SStefano Zampini const PetscInt id = 1; 92*736140d5SStefano Zampini 93*736140d5SStefano Zampini PetscFunctionBeginUser; 94*736140d5SStefano Zampini PetscCall(DMGetDS(dm, &ds)); 95*736140d5SStefano Zampini PetscCall(DMGetDimension(dm, &dim)); 96*736140d5SStefano Zampini PetscCheck(dim == 2, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only for 2-D\n"); 97*736140d5SStefano Zampini PetscCall(PetscDSSetObjective(ds, 0, objective_2d)); 98*736140d5SStefano Zampini PetscCall(PetscDSSetResidual(ds, 0, NULL, gradient_1_2d)); 99*736140d5SStefano Zampini PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, hessian_11_2d)); 100*736140d5SStefano Zampini PetscCall(DMGetLabel(dm, "marker", &label)); 101*736140d5SStefano Zampini PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "data", label, 1, &id, 0, 0, NULL, (void (*)(void))sins_2d, NULL, NULL, NULL)); 102*736140d5SStefano Zampini PetscCall(DMPlexSetSNESLocalFEM(dm, PETSC_TRUE, NULL)); 103*736140d5SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 104*736140d5SStefano Zampini } 105*736140d5SStefano Zampini 106*736140d5SStefano Zampini static PetscErrorCode SetupDiscretization(DM dm) 107*736140d5SStefano Zampini { 108*736140d5SStefano Zampini DM plex, cdm = dm; 109*736140d5SStefano Zampini PetscFE fe; 110*736140d5SStefano Zampini PetscBool simplex; 111*736140d5SStefano Zampini PetscInt dim; 112*736140d5SStefano Zampini 113*736140d5SStefano Zampini PetscFunctionBeginUser; 114*736140d5SStefano Zampini PetscCall(DMGetDimension(dm, &dim)); 115*736140d5SStefano Zampini PetscCall(DMConvert(dm, DMPLEX, &plex)); 116*736140d5SStefano Zampini PetscCall(DMPlexIsSimplex(plex, &simplex)); 117*736140d5SStefano Zampini PetscCall(DMDestroy(&plex)); 118*736140d5SStefano Zampini PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, NULL, -1, &fe)); 119*736140d5SStefano Zampini PetscCall(PetscObjectSetName((PetscObject)fe, "potential")); 120*736140d5SStefano Zampini PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe)); 121*736140d5SStefano Zampini PetscCall(DMCreateDS(dm)); 122*736140d5SStefano Zampini PetscCall(SetupProblem(dm)); 123*736140d5SStefano Zampini while (cdm) { 124*736140d5SStefano Zampini PetscBool hasLabel; 125*736140d5SStefano Zampini 126*736140d5SStefano Zampini PetscCall(DMHasLabel(cdm, "marker", &hasLabel)); 127*736140d5SStefano Zampini if (!hasLabel) PetscCall(CreateBCLabel(cdm, "marker")); 128*736140d5SStefano Zampini PetscCall(DMCopyDisc(dm, cdm)); 129*736140d5SStefano Zampini PetscCall(DMGetCoarseDM(cdm, &cdm)); 130*736140d5SStefano Zampini } 131*736140d5SStefano Zampini PetscCall(PetscFEDestroy(&fe)); 132*736140d5SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 133*736140d5SStefano Zampini } 134*736140d5SStefano Zampini 135*736140d5SStefano Zampini int main(int argc, char **argv) 136*736140d5SStefano Zampini { 137*736140d5SStefano Zampini DM dm; /* Problem specification */ 138*736140d5SStefano Zampini SNES snes; /* nonlinear solver */ 139*736140d5SStefano Zampini Vec u; /* solution vector */ 140*736140d5SStefano Zampini 141*736140d5SStefano Zampini PetscFunctionBeginUser; 142*736140d5SStefano Zampini PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 143*736140d5SStefano Zampini PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 144*736140d5SStefano Zampini PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm)); 145*736140d5SStefano Zampini PetscCall(SNESSetDM(snes, dm)); 146*736140d5SStefano Zampini 147*736140d5SStefano Zampini PetscCall(SetupDiscretization(dm)); 148*736140d5SStefano Zampini 149*736140d5SStefano Zampini PetscCall(SNESSetFromOptions(snes)); 150*736140d5SStefano Zampini 151*736140d5SStefano Zampini PetscCall(DMCreateGlobalVector(dm, &u)); 152*736140d5SStefano Zampini PetscCall(SNESSolve(snes, NULL, u)); 153*736140d5SStefano Zampini 154*736140d5SStefano Zampini PetscCall(VecDestroy(&u)); 155*736140d5SStefano Zampini PetscCall(SNESDestroy(&snes)); 156*736140d5SStefano Zampini PetscCall(DMDestroy(&dm)); 157*736140d5SStefano Zampini PetscCall(PetscFinalize()); 158*736140d5SStefano Zampini return 0; 159*736140d5SStefano Zampini } 160*736140d5SStefano Zampini 161*736140d5SStefano Zampini /*TEST 162*736140d5SStefano Zampini 163*736140d5SStefano Zampini test: 164*736140d5SStefano Zampini requires: !single 165*736140d5SStefano Zampini suffix: qn_nasm 166*736140d5SStefano Zampini filter: sed -e "s/CONVERGED_FNORM_ABS/CONVERGED_FNORM_RELATIVE/g" 167*736140d5SStefano Zampini nsize: 4 168*736140d5SStefano Zampini args: -petscspace_degree 1 -dm_refine 2 -snes_type qn -snes_npc_side {{left right}separate output} -npc_snes_type nasm -npc_snes_nasm_type restrict -npc_sub_snes_linesearch_order 1 -npc_sub_snes_linesearch_type bt -dm_plex_dd_overlap 1 -snes_linesearch_type bt -snes_linesearch_order 1 -npc_sub_pc_type lu -npc_sub_ksp_type preonly -snes_converged_reason -snes_monitor_short -petscpartitioner_type simple -npc_sub_snes_max_it 1 -dm_plex_simplex 0 -snes_rtol 1.e-6 169*736140d5SStefano Zampini 170*736140d5SStefano Zampini TEST*/ 171