1*c4762a1bSJed Brown 2*c4762a1bSJed Brown static const char help[] = "Solves PDE optimization problem using full-space method, treats state and adjoint variables separately.\n\n"; 3*c4762a1bSJed Brown 4*c4762a1bSJed Brown #include <petscdm.h> 5*c4762a1bSJed Brown #include <petscdmda.h> 6*c4762a1bSJed Brown #include <petscdmredundant.h> 7*c4762a1bSJed Brown #include <petscdmcomposite.h> 8*c4762a1bSJed Brown #include <petscpf.h> 9*c4762a1bSJed Brown #include <petscsnes.h> 10*c4762a1bSJed Brown 11*c4762a1bSJed Brown /* 12*c4762a1bSJed Brown 13*c4762a1bSJed Brown w - design variables (what we change to get an optimal solution) 14*c4762a1bSJed Brown u - state variables (i.e. the PDE solution) 15*c4762a1bSJed Brown lambda - the Lagrange multipliers 16*c4762a1bSJed Brown 17*c4762a1bSJed Brown U = (w u lambda) 18*c4762a1bSJed Brown 19*c4762a1bSJed Brown fu, fw, flambda contain the gradient of L(w,u,lambda) 20*c4762a1bSJed Brown 21*c4762a1bSJed Brown FU = (fw fu flambda) 22*c4762a1bSJed Brown 23*c4762a1bSJed Brown In this example the PDE is 24*c4762a1bSJed Brown Uxx = 2, 25*c4762a1bSJed Brown u(0) = w(0), thus this is the free parameter 26*c4762a1bSJed Brown u(1) = 0 27*c4762a1bSJed Brown the function we wish to minimize is 28*c4762a1bSJed Brown \integral u^{2} 29*c4762a1bSJed Brown 30*c4762a1bSJed Brown The exact solution for u is given by u(x) = x*x - 1.25*x + .25 31*c4762a1bSJed Brown 32*c4762a1bSJed Brown Use the usual centered finite differences. 33*c4762a1bSJed Brown 34*c4762a1bSJed Brown Note we treat the problem as non-linear though it happens to be linear 35*c4762a1bSJed Brown 36*c4762a1bSJed Brown See ex22.c for the same code, but that interlaces the u and the lambda 37*c4762a1bSJed Brown 38*c4762a1bSJed Brown */ 39*c4762a1bSJed Brown 40*c4762a1bSJed Brown typedef struct { 41*c4762a1bSJed Brown DM red1,da1,da2; 42*c4762a1bSJed Brown DM packer; 43*c4762a1bSJed Brown PetscViewer u_viewer,lambda_viewer; 44*c4762a1bSJed Brown PetscViewer fu_viewer,flambda_viewer; 45*c4762a1bSJed Brown } UserCtx; 46*c4762a1bSJed Brown 47*c4762a1bSJed Brown extern PetscErrorCode FormFunction(SNES,Vec,Vec,void*); 48*c4762a1bSJed Brown extern PetscErrorCode Monitor(SNES,PetscInt,PetscReal,void*); 49*c4762a1bSJed Brown 50*c4762a1bSJed Brown 51*c4762a1bSJed Brown int main(int argc,char **argv) 52*c4762a1bSJed Brown { 53*c4762a1bSJed Brown PetscErrorCode ierr; 54*c4762a1bSJed Brown PetscInt its; 55*c4762a1bSJed Brown Vec U,FU; 56*c4762a1bSJed Brown SNES snes; 57*c4762a1bSJed Brown UserCtx user; 58*c4762a1bSJed Brown 59*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 60*c4762a1bSJed Brown 61*c4762a1bSJed Brown /* Create a global vector that includes a single redundant array and two da arrays */ 62*c4762a1bSJed Brown ierr = DMCompositeCreate(PETSC_COMM_WORLD,&user.packer);CHKERRQ(ierr); 63*c4762a1bSJed Brown ierr = DMRedundantCreate(PETSC_COMM_WORLD,0,1,&user.red1);CHKERRQ(ierr); 64*c4762a1bSJed Brown ierr = DMCompositeAddDM(user.packer,user.red1);CHKERRQ(ierr); 65*c4762a1bSJed Brown ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,5,1,1,NULL,&user.da1);CHKERRQ(ierr); 66*c4762a1bSJed Brown ierr = DMSetFromOptions(user.da1);CHKERRQ(ierr); 67*c4762a1bSJed Brown ierr = DMSetUp(user.da1);CHKERRQ(ierr); 68*c4762a1bSJed Brown ierr = DMCompositeAddDM(user.packer,user.da1);CHKERRQ(ierr); 69*c4762a1bSJed Brown ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,5,1,1,NULL,&user.da2);CHKERRQ(ierr); 70*c4762a1bSJed Brown ierr = DMSetFromOptions(user.da2);CHKERRQ(ierr); 71*c4762a1bSJed Brown ierr = DMSetUp(user.da2);CHKERRQ(ierr); 72*c4762a1bSJed Brown ierr = DMDASetFieldName(user.da1,0,"u");CHKERRQ(ierr); 73*c4762a1bSJed Brown ierr = DMDASetFieldName(user.da2,0,"lambda");CHKERRQ(ierr); 74*c4762a1bSJed Brown ierr = DMCompositeAddDM(user.packer,user.da2);CHKERRQ(ierr); 75*c4762a1bSJed Brown ierr = DMCreateGlobalVector(user.packer,&U);CHKERRQ(ierr); 76*c4762a1bSJed Brown ierr = VecDuplicate(U,&FU);CHKERRQ(ierr); 77*c4762a1bSJed Brown 78*c4762a1bSJed Brown /* create graphics windows */ 79*c4762a1bSJed Brown ierr = PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"u - state variables",-1,-1,-1,-1,&user.u_viewer);CHKERRQ(ierr); 80*c4762a1bSJed Brown ierr = PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"lambda - Lagrange multipliers",-1,-1,-1,-1,&user.lambda_viewer);CHKERRQ(ierr); 81*c4762a1bSJed Brown ierr = PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"fu - derivate w.r.t. state variables",-1,-1,-1,-1,&user.fu_viewer);CHKERRQ(ierr); 82*c4762a1bSJed Brown ierr = PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"flambda - derivate w.r.t. Lagrange multipliers",-1,-1,-1,-1,&user.flambda_viewer);CHKERRQ(ierr); 83*c4762a1bSJed Brown 84*c4762a1bSJed Brown 85*c4762a1bSJed Brown /* create nonlinear solver */ 86*c4762a1bSJed Brown ierr = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(ierr); 87*c4762a1bSJed Brown ierr = SNESSetFunction(snes,FU,FormFunction,&user);CHKERRQ(ierr); 88*c4762a1bSJed Brown ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 89*c4762a1bSJed Brown ierr = SNESMonitorSet(snes,Monitor,&user,0);CHKERRQ(ierr); 90*c4762a1bSJed Brown ierr = SNESSolve(snes,NULL,U);CHKERRQ(ierr); 91*c4762a1bSJed Brown ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr); 92*c4762a1bSJed Brown ierr = SNESDestroy(&snes);CHKERRQ(ierr); 93*c4762a1bSJed Brown 94*c4762a1bSJed Brown ierr = DMDestroy(&user.red1);CHKERRQ(ierr); 95*c4762a1bSJed Brown ierr = DMDestroy(&user.da1);CHKERRQ(ierr); 96*c4762a1bSJed Brown ierr = DMDestroy(&user.da2);CHKERRQ(ierr); 97*c4762a1bSJed Brown ierr = DMDestroy(&user.packer);CHKERRQ(ierr); 98*c4762a1bSJed Brown ierr = VecDestroy(&U);CHKERRQ(ierr); 99*c4762a1bSJed Brown ierr = VecDestroy(&FU);CHKERRQ(ierr); 100*c4762a1bSJed Brown ierr = PetscViewerDestroy(&user.u_viewer);CHKERRQ(ierr); 101*c4762a1bSJed Brown ierr = PetscViewerDestroy(&user.lambda_viewer);CHKERRQ(ierr); 102*c4762a1bSJed Brown ierr = PetscViewerDestroy(&user.fu_viewer);CHKERRQ(ierr); 103*c4762a1bSJed Brown ierr = PetscViewerDestroy(&user.flambda_viewer);CHKERRQ(ierr); 104*c4762a1bSJed Brown ierr = PetscFinalize(); 105*c4762a1bSJed Brown return ierr; 106*c4762a1bSJed Brown } 107*c4762a1bSJed Brown 108*c4762a1bSJed Brown /* 109*c4762a1bSJed Brown Evaluates FU = Gradiant(L(w,u,lambda)) 110*c4762a1bSJed Brown 111*c4762a1bSJed Brown */ 112*c4762a1bSJed Brown PetscErrorCode FormFunction(SNES snes,Vec U,Vec FU,void *dummy) 113*c4762a1bSJed Brown { 114*c4762a1bSJed Brown UserCtx *user = (UserCtx*)dummy; 115*c4762a1bSJed Brown PetscErrorCode ierr; 116*c4762a1bSJed Brown PetscInt xs,xm,i,N; 117*c4762a1bSJed Brown PetscScalar *u,*lambda,*w,*fu,*fw,*flambda,d,h; 118*c4762a1bSJed Brown Vec vw,vu,vlambda,vfw,vfu,vflambda; 119*c4762a1bSJed Brown 120*c4762a1bSJed Brown PetscFunctionBeginUser; 121*c4762a1bSJed Brown ierr = DMCompositeGetLocalVectors(user->packer,&vw,&vu,&vlambda);CHKERRQ(ierr); 122*c4762a1bSJed Brown ierr = DMCompositeGetLocalVectors(user->packer,&vfw,&vfu,&vflambda);CHKERRQ(ierr); 123*c4762a1bSJed Brown ierr = DMCompositeScatter(user->packer,U,vw,vu,vlambda);CHKERRQ(ierr); 124*c4762a1bSJed Brown 125*c4762a1bSJed Brown ierr = DMDAGetCorners(user->da1,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr); 126*c4762a1bSJed Brown ierr = DMDAGetInfo(user->da1,0,&N,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr); 127*c4762a1bSJed Brown ierr = VecGetArray(vw,&w);CHKERRQ(ierr); 128*c4762a1bSJed Brown ierr = VecGetArray(vfw,&fw);CHKERRQ(ierr); 129*c4762a1bSJed Brown ierr = DMDAVecGetArray(user->da1,vu,&u);CHKERRQ(ierr); 130*c4762a1bSJed Brown ierr = DMDAVecGetArray(user->da1,vfu,&fu);CHKERRQ(ierr); 131*c4762a1bSJed Brown ierr = DMDAVecGetArray(user->da1,vlambda,&lambda);CHKERRQ(ierr); 132*c4762a1bSJed Brown ierr = DMDAVecGetArray(user->da1,vflambda,&flambda);CHKERRQ(ierr); 133*c4762a1bSJed Brown d = (N-1.0); 134*c4762a1bSJed Brown h = 1.0/d; 135*c4762a1bSJed Brown 136*c4762a1bSJed Brown /* derivative of L() w.r.t. w */ 137*c4762a1bSJed Brown if (xs == 0) { /* only first processor computes this */ 138*c4762a1bSJed Brown fw[0] = -2.*d*lambda[0]; 139*c4762a1bSJed Brown } 140*c4762a1bSJed Brown 141*c4762a1bSJed Brown /* derivative of L() w.r.t. u */ 142*c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 143*c4762a1bSJed Brown if (i == 0) flambda[0] = h*u[0] + 2.*d*lambda[0] - d*lambda[1]; 144*c4762a1bSJed Brown else if (i == 1) flambda[1] = 2.*h*u[1] + 2.*d*lambda[1] - d*lambda[2]; 145*c4762a1bSJed Brown else if (i == N-1) flambda[N-1] = h*u[N-1] + 2.*d*lambda[N-1] - d*lambda[N-2]; 146*c4762a1bSJed Brown else if (i == N-2) flambda[N-2] = 2.*h*u[N-2] + 2.*d*lambda[N-2] - d*lambda[N-3]; 147*c4762a1bSJed Brown else flambda[i] = 2.*h*u[i] - d*(lambda[i+1] - 2.0*lambda[i] + lambda[i-1]); 148*c4762a1bSJed Brown } 149*c4762a1bSJed Brown 150*c4762a1bSJed Brown /* derivative of L() w.r.t. lambda */ 151*c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 152*c4762a1bSJed Brown if (i == 0) fu[0] = 2.0*d*(u[0] - w[0]); 153*c4762a1bSJed Brown else if (i == N-1) fu[N-1] = 2.0*d*u[N-1]; 154*c4762a1bSJed Brown else fu[i] = -(d*(u[i+1] - 2.0*u[i] + u[i-1]) - 2.0*h); 155*c4762a1bSJed Brown } 156*c4762a1bSJed Brown 157*c4762a1bSJed Brown ierr = VecRestoreArray(vw,&w);CHKERRQ(ierr); 158*c4762a1bSJed Brown ierr = VecRestoreArray(vfw,&fw);CHKERRQ(ierr); 159*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(user->da1,vu,&u);CHKERRQ(ierr); 160*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(user->da1,vfu,&fu);CHKERRQ(ierr); 161*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(user->da1,vlambda,&lambda);CHKERRQ(ierr); 162*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(user->da1,vflambda,&flambda);CHKERRQ(ierr); 163*c4762a1bSJed Brown 164*c4762a1bSJed Brown ierr = DMCompositeGather(user->packer,INSERT_VALUES,FU,vfw,vfu,vflambda);CHKERRQ(ierr); 165*c4762a1bSJed Brown ierr = DMCompositeRestoreLocalVectors(user->packer,&vw,&vu,&vlambda);CHKERRQ(ierr); 166*c4762a1bSJed Brown ierr = DMCompositeRestoreLocalVectors(user->packer,&vfw,&vfu,&vflambda);CHKERRQ(ierr); 167*c4762a1bSJed Brown PetscFunctionReturn(0); 168*c4762a1bSJed Brown } 169*c4762a1bSJed Brown 170*c4762a1bSJed Brown PetscErrorCode Monitor(SNES snes,PetscInt its,PetscReal rnorm,void *dummy) 171*c4762a1bSJed Brown { 172*c4762a1bSJed Brown UserCtx *user = (UserCtx*)dummy; 173*c4762a1bSJed Brown PetscErrorCode ierr; 174*c4762a1bSJed Brown Vec w,u,lambda,U,F; 175*c4762a1bSJed Brown 176*c4762a1bSJed Brown PetscFunctionBeginUser; 177*c4762a1bSJed Brown ierr = SNESGetSolution(snes,&U);CHKERRQ(ierr); 178*c4762a1bSJed Brown ierr = DMCompositeGetAccess(user->packer,U,&w,&u,&lambda);CHKERRQ(ierr); 179*c4762a1bSJed Brown ierr = VecView(u,user->u_viewer);CHKERRQ(ierr); 180*c4762a1bSJed Brown ierr = VecView(lambda,user->lambda_viewer);CHKERRQ(ierr); 181*c4762a1bSJed Brown ierr = DMCompositeRestoreAccess(user->packer,U,&w,&u,&lambda);CHKERRQ(ierr); 182*c4762a1bSJed Brown 183*c4762a1bSJed Brown ierr = SNESGetFunction(snes,&F,0,0);CHKERRQ(ierr); 184*c4762a1bSJed Brown ierr = DMCompositeGetAccess(user->packer,F,&w,&u,&lambda);CHKERRQ(ierr); 185*c4762a1bSJed Brown ierr = VecView(u,user->fu_viewer);CHKERRQ(ierr); 186*c4762a1bSJed Brown ierr = VecView(lambda,user->flambda_viewer);CHKERRQ(ierr); 187*c4762a1bSJed Brown ierr = DMCompositeRestoreAccess(user->packer,F,&w,&u,&lambda);CHKERRQ(ierr); 188*c4762a1bSJed Brown PetscFunctionReturn(0); 189*c4762a1bSJed Brown } 190*c4762a1bSJed Brown 191*c4762a1bSJed Brown /*TEST 192*c4762a1bSJed Brown 193*c4762a1bSJed Brown test: 194*c4762a1bSJed Brown nsize: 4 195*c4762a1bSJed Brown args: -snes_linesearch_monitor -snes_mf -pc_type none -snes_monitor_short -nox -ksp_monitor_short -snes_converged_reason 196*c4762a1bSJed Brown requires: !single 197*c4762a1bSJed Brown 198*c4762a1bSJed Brown TEST*/ 199