1c4762a1bSJed Brown 2c4762a1bSJed Brown static const char help[] = "Solves PDE optimization problem using full-space method, treats state and adjoint variables separately.\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscdm.h> 5c4762a1bSJed Brown #include <petscdmda.h> 6c4762a1bSJed Brown #include <petscdmredundant.h> 7c4762a1bSJed Brown #include <petscdmcomposite.h> 8c4762a1bSJed Brown #include <petscpf.h> 9c4762a1bSJed Brown #include <petscsnes.h> 10c4762a1bSJed Brown 11c4762a1bSJed Brown /* 12c4762a1bSJed Brown 13c4762a1bSJed Brown w - design variables (what we change to get an optimal solution) 14c4762a1bSJed Brown u - state variables (i.e. the PDE solution) 15c4762a1bSJed Brown lambda - the Lagrange multipliers 16c4762a1bSJed Brown 17c4762a1bSJed Brown U = (w u lambda) 18c4762a1bSJed Brown 19c4762a1bSJed Brown fu, fw, flambda contain the gradient of L(w,u,lambda) 20c4762a1bSJed Brown 21c4762a1bSJed Brown FU = (fw fu flambda) 22c4762a1bSJed Brown 23c4762a1bSJed Brown In this example the PDE is 24c4762a1bSJed Brown Uxx = 2, 25c4762a1bSJed Brown u(0) = w(0), thus this is the free parameter 26c4762a1bSJed Brown u(1) = 0 27c4762a1bSJed Brown the function we wish to minimize is 28c4762a1bSJed Brown \integral u^{2} 29c4762a1bSJed Brown 30c4762a1bSJed Brown The exact solution for u is given by u(x) = x*x - 1.25*x + .25 31c4762a1bSJed Brown 32c4762a1bSJed Brown Use the usual centered finite differences. 33c4762a1bSJed Brown 34c4762a1bSJed Brown Note we treat the problem as non-linear though it happens to be linear 35c4762a1bSJed Brown 36c4762a1bSJed Brown See ex22.c for the same code, but that interlaces the u and the lambda 37c4762a1bSJed Brown 38c4762a1bSJed Brown */ 39c4762a1bSJed Brown 40c4762a1bSJed Brown typedef struct { 41c4762a1bSJed Brown DM red1, da1, da2; 42c4762a1bSJed Brown DM packer; 43c4762a1bSJed Brown PetscViewer u_viewer, lambda_viewer; 44c4762a1bSJed Brown PetscViewer fu_viewer, flambda_viewer; 45c4762a1bSJed Brown } UserCtx; 46c4762a1bSJed Brown 47c4762a1bSJed Brown extern PetscErrorCode FormFunction(SNES, Vec, Vec, void *); 48c4762a1bSJed Brown extern PetscErrorCode Monitor(SNES, PetscInt, PetscReal, void *); 49c4762a1bSJed Brown 50d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 51d71ae5a4SJacob Faibussowitsch { 52c4762a1bSJed Brown PetscInt its; 53c4762a1bSJed Brown Vec U, FU; 54c4762a1bSJed Brown SNES snes; 55c4762a1bSJed Brown UserCtx user; 56c4762a1bSJed Brown 57327415f7SBarry Smith PetscFunctionBeginUser; 589566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 59c4762a1bSJed Brown 60c4762a1bSJed Brown /* Create a global vector that includes a single redundant array and two da arrays */ 619566063dSJacob Faibussowitsch PetscCall(DMCompositeCreate(PETSC_COMM_WORLD, &user.packer)); 629566063dSJacob Faibussowitsch PetscCall(DMRedundantCreate(PETSC_COMM_WORLD, 0, 1, &user.red1)); 639566063dSJacob Faibussowitsch PetscCall(DMCompositeAddDM(user.packer, user.red1)); 649566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, 5, 1, 1, NULL, &user.da1)); 659566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(user.da1)); 669566063dSJacob Faibussowitsch PetscCall(DMSetUp(user.da1)); 679566063dSJacob Faibussowitsch PetscCall(DMCompositeAddDM(user.packer, user.da1)); 689566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, 5, 1, 1, NULL, &user.da2)); 699566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(user.da2)); 709566063dSJacob Faibussowitsch PetscCall(DMSetUp(user.da2)); 719566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(user.da1, 0, "u")); 729566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(user.da2, 0, "lambda")); 739566063dSJacob Faibussowitsch PetscCall(DMCompositeAddDM(user.packer, user.da2)); 749566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(user.packer, &U)); 759566063dSJacob Faibussowitsch PetscCall(VecDuplicate(U, &FU)); 76c4762a1bSJed Brown 77c4762a1bSJed Brown /* create graphics windows */ 789566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawOpen(PETSC_COMM_WORLD, 0, "u - state variables", -1, -1, -1, -1, &user.u_viewer)); 799566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawOpen(PETSC_COMM_WORLD, 0, "lambda - Lagrange multipliers", -1, -1, -1, -1, &user.lambda_viewer)); 809566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawOpen(PETSC_COMM_WORLD, 0, "fu - derivate w.r.t. state variables", -1, -1, -1, -1, &user.fu_viewer)); 819566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawOpen(PETSC_COMM_WORLD, 0, "flambda - derivate w.r.t. Lagrange multipliers", -1, -1, -1, -1, &user.flambda_viewer)); 82c4762a1bSJed Brown 83c4762a1bSJed Brown /* create nonlinear solver */ 849566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 859566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes, FU, FormFunction, &user)); 869566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 879566063dSJacob Faibussowitsch PetscCall(SNESMonitorSet(snes, Monitor, &user, 0)); 889566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, U)); 899566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &its)); 909566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 91c4762a1bSJed Brown 929566063dSJacob Faibussowitsch PetscCall(DMDestroy(&user.red1)); 939566063dSJacob Faibussowitsch PetscCall(DMDestroy(&user.da1)); 949566063dSJacob Faibussowitsch PetscCall(DMDestroy(&user.da2)); 959566063dSJacob Faibussowitsch PetscCall(DMDestroy(&user.packer)); 969566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U)); 979566063dSJacob Faibussowitsch PetscCall(VecDestroy(&FU)); 989566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&user.u_viewer)); 999566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&user.lambda_viewer)); 1009566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&user.fu_viewer)); 1019566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&user.flambda_viewer)); 1029566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 103b122ec5aSJacob Faibussowitsch return 0; 104c4762a1bSJed Brown } 105c4762a1bSJed Brown 106c4762a1bSJed Brown /* 107*d5b43468SJose E. Roman Evaluates FU = Gradient(L(w,u,lambda)) 108c4762a1bSJed Brown 109c4762a1bSJed Brown */ 110d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunction(SNES snes, Vec U, Vec FU, void *dummy) 111d71ae5a4SJacob Faibussowitsch { 112c4762a1bSJed Brown UserCtx *user = (UserCtx *)dummy; 113c4762a1bSJed Brown PetscInt xs, xm, i, N; 114c4762a1bSJed Brown PetscScalar *u, *lambda, *w, *fu, *fw, *flambda, d, h; 115c4762a1bSJed Brown Vec vw, vu, vlambda, vfw, vfu, vflambda; 116c4762a1bSJed Brown 117c4762a1bSJed Brown PetscFunctionBeginUser; 1189566063dSJacob Faibussowitsch PetscCall(DMCompositeGetLocalVectors(user->packer, &vw, &vu, &vlambda)); 1199566063dSJacob Faibussowitsch PetscCall(DMCompositeGetLocalVectors(user->packer, &vfw, &vfu, &vflambda)); 1209566063dSJacob Faibussowitsch PetscCall(DMCompositeScatter(user->packer, U, vw, vu, vlambda)); 121c4762a1bSJed Brown 1229566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(user->da1, &xs, NULL, NULL, &xm, NULL, NULL)); 1239566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(user->da1, 0, &N, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0)); 1249566063dSJacob Faibussowitsch PetscCall(VecGetArray(vw, &w)); 1259566063dSJacob Faibussowitsch PetscCall(VecGetArray(vfw, &fw)); 1269566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(user->da1, vu, &u)); 1279566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(user->da1, vfu, &fu)); 1289566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(user->da1, vlambda, &lambda)); 1299566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(user->da1, vflambda, &flambda)); 130c4762a1bSJed Brown d = (N - 1.0); 131c4762a1bSJed Brown h = 1.0 / d; 132c4762a1bSJed Brown 133c4762a1bSJed Brown /* derivative of L() w.r.t. w */ 134c4762a1bSJed Brown if (xs == 0) { /* only first processor computes this */ 135c4762a1bSJed Brown fw[0] = -2. * d * lambda[0]; 136c4762a1bSJed Brown } 137c4762a1bSJed Brown 138c4762a1bSJed Brown /* derivative of L() w.r.t. u */ 139c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 140c4762a1bSJed Brown if (i == 0) flambda[0] = h * u[0] + 2. * d * lambda[0] - d * lambda[1]; 141c4762a1bSJed Brown else if (i == 1) flambda[1] = 2. * h * u[1] + 2. * d * lambda[1] - d * lambda[2]; 142c4762a1bSJed Brown else if (i == N - 1) flambda[N - 1] = h * u[N - 1] + 2. * d * lambda[N - 1] - d * lambda[N - 2]; 143c4762a1bSJed Brown else if (i == N - 2) flambda[N - 2] = 2. * h * u[N - 2] + 2. * d * lambda[N - 2] - d * lambda[N - 3]; 144c4762a1bSJed Brown else flambda[i] = 2. * h * u[i] - d * (lambda[i + 1] - 2.0 * lambda[i] + lambda[i - 1]); 145c4762a1bSJed Brown } 146c4762a1bSJed Brown 147c4762a1bSJed Brown /* derivative of L() w.r.t. lambda */ 148c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 149c4762a1bSJed Brown if (i == 0) fu[0] = 2.0 * d * (u[0] - w[0]); 150c4762a1bSJed Brown else if (i == N - 1) fu[N - 1] = 2.0 * d * u[N - 1]; 151c4762a1bSJed Brown else fu[i] = -(d * (u[i + 1] - 2.0 * u[i] + u[i - 1]) - 2.0 * h); 152c4762a1bSJed Brown } 153c4762a1bSJed Brown 1549566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(vw, &w)); 1559566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(vfw, &fw)); 1569566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(user->da1, vu, &u)); 1579566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(user->da1, vfu, &fu)); 1589566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(user->da1, vlambda, &lambda)); 1599566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(user->da1, vflambda, &flambda)); 160c4762a1bSJed Brown 1619566063dSJacob Faibussowitsch PetscCall(DMCompositeGather(user->packer, INSERT_VALUES, FU, vfw, vfu, vflambda)); 1629566063dSJacob Faibussowitsch PetscCall(DMCompositeRestoreLocalVectors(user->packer, &vw, &vu, &vlambda)); 1639566063dSJacob Faibussowitsch PetscCall(DMCompositeRestoreLocalVectors(user->packer, &vfw, &vfu, &vflambda)); 164c4762a1bSJed Brown PetscFunctionReturn(0); 165c4762a1bSJed Brown } 166c4762a1bSJed Brown 167d71ae5a4SJacob Faibussowitsch PetscErrorCode Monitor(SNES snes, PetscInt its, PetscReal rnorm, void *dummy) 168d71ae5a4SJacob Faibussowitsch { 169c4762a1bSJed Brown UserCtx *user = (UserCtx *)dummy; 170c4762a1bSJed Brown Vec w, u, lambda, U, F; 171c4762a1bSJed Brown 172c4762a1bSJed Brown PetscFunctionBeginUser; 1739566063dSJacob Faibussowitsch PetscCall(SNESGetSolution(snes, &U)); 1749566063dSJacob Faibussowitsch PetscCall(DMCompositeGetAccess(user->packer, U, &w, &u, &lambda)); 1759566063dSJacob Faibussowitsch PetscCall(VecView(u, user->u_viewer)); 1769566063dSJacob Faibussowitsch PetscCall(VecView(lambda, user->lambda_viewer)); 1779566063dSJacob Faibussowitsch PetscCall(DMCompositeRestoreAccess(user->packer, U, &w, &u, &lambda)); 178c4762a1bSJed Brown 1799566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(snes, &F, 0, 0)); 1809566063dSJacob Faibussowitsch PetscCall(DMCompositeGetAccess(user->packer, F, &w, &u, &lambda)); 1819566063dSJacob Faibussowitsch PetscCall(VecView(u, user->fu_viewer)); 1829566063dSJacob Faibussowitsch PetscCall(VecView(lambda, user->flambda_viewer)); 1839566063dSJacob Faibussowitsch PetscCall(DMCompositeRestoreAccess(user->packer, F, &w, &u, &lambda)); 184c4762a1bSJed Brown PetscFunctionReturn(0); 185c4762a1bSJed Brown } 186c4762a1bSJed Brown 187c4762a1bSJed Brown /*TEST 188c4762a1bSJed Brown 189c4762a1bSJed Brown test: 190c4762a1bSJed Brown nsize: 4 191c4762a1bSJed Brown args: -snes_linesearch_monitor -snes_mf -pc_type none -snes_monitor_short -nox -ksp_monitor_short -snes_converged_reason 192c4762a1bSJed Brown requires: !single 193c4762a1bSJed Brown 194c4762a1bSJed Brown TEST*/ 195