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 50c4762a1bSJed Brown int main(int argc,char **argv) 51c4762a1bSJed Brown { 52c4762a1bSJed Brown PetscErrorCode ierr; 53c4762a1bSJed Brown PetscInt its; 54c4762a1bSJed Brown Vec U,FU; 55c4762a1bSJed Brown SNES snes; 56c4762a1bSJed Brown UserCtx user; 57c4762a1bSJed Brown 58c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 59c4762a1bSJed Brown 60c4762a1bSJed Brown /* Create a global vector that includes a single redundant array and two da arrays */ 61*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCompositeCreate(PETSC_COMM_WORLD,&user.packer)); 62*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMRedundantCreate(PETSC_COMM_WORLD,0,1,&user.red1)); 63*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCompositeAddDM(user.packer,user.red1)); 64*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,5,1,1,NULL,&user.da1)); 65*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(user.da1)); 66*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(user.da1)); 67*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCompositeAddDM(user.packer,user.da1)); 68*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,5,1,1,NULL,&user.da2)); 69*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(user.da2)); 70*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(user.da2)); 71*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetFieldName(user.da1,0,"u")); 72*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetFieldName(user.da2,0,"lambda")); 73*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCompositeAddDM(user.packer,user.da2)); 74*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(user.packer,&U)); 75*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(U,&FU)); 76c4762a1bSJed Brown 77c4762a1bSJed Brown /* create graphics windows */ 78*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"u - state variables",-1,-1,-1,-1,&user.u_viewer)); 79*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"lambda - Lagrange multipliers",-1,-1,-1,-1,&user.lambda_viewer)); 80*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"fu - derivate w.r.t. state variables",-1,-1,-1,-1,&user.fu_viewer)); 81*5f80ce2aSJacob Faibussowitsch CHKERRQ(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 */ 84*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESCreate(PETSC_COMM_WORLD,&snes)); 85*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetFunction(snes,FU,FormFunction,&user)); 86*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetFromOptions(snes)); 87*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESMonitorSet(snes,Monitor,&user,0)); 88*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSolve(snes,NULL,U)); 89*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetIterationNumber(snes,&its)); 90*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESDestroy(&snes)); 91c4762a1bSJed Brown 92*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&user.red1)); 93*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&user.da1)); 94*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&user.da2)); 95*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&user.packer)); 96*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&U)); 97*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&FU)); 98*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&user.u_viewer)); 99*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&user.lambda_viewer)); 100*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&user.fu_viewer)); 101*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&user.flambda_viewer)); 102c4762a1bSJed Brown ierr = PetscFinalize(); 103c4762a1bSJed Brown return ierr; 104c4762a1bSJed Brown } 105c4762a1bSJed Brown 106c4762a1bSJed Brown /* 107c4762a1bSJed Brown Evaluates FU = Gradiant(L(w,u,lambda)) 108c4762a1bSJed Brown 109c4762a1bSJed Brown */ 110c4762a1bSJed Brown PetscErrorCode FormFunction(SNES snes,Vec U,Vec FU,void *dummy) 111c4762a1bSJed Brown { 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; 118*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCompositeGetLocalVectors(user->packer,&vw,&vu,&vlambda)); 119*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCompositeGetLocalVectors(user->packer,&vfw,&vfu,&vflambda)); 120*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCompositeScatter(user->packer,U,vw,vu,vlambda)); 121c4762a1bSJed Brown 122*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(user->da1,&xs,NULL,NULL,&xm,NULL,NULL)); 123*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetInfo(user->da1,0,&N,0,0,0,0,0,0,0,0,0,0,0)); 124*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(vw,&w)); 125*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(vfw,&fw)); 126*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(user->da1,vu,&u)); 127*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(user->da1,vfu,&fu)); 128*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(user->da1,vlambda,&lambda)); 129*5f80ce2aSJacob Faibussowitsch CHKERRQ(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 154*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(vw,&w)); 155*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(vfw,&fw)); 156*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(user->da1,vu,&u)); 157*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(user->da1,vfu,&fu)); 158*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(user->da1,vlambda,&lambda)); 159*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(user->da1,vflambda,&flambda)); 160c4762a1bSJed Brown 161*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCompositeGather(user->packer,INSERT_VALUES,FU,vfw,vfu,vflambda)); 162*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCompositeRestoreLocalVectors(user->packer,&vw,&vu,&vlambda)); 163*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCompositeRestoreLocalVectors(user->packer,&vfw,&vfu,&vflambda)); 164c4762a1bSJed Brown PetscFunctionReturn(0); 165c4762a1bSJed Brown } 166c4762a1bSJed Brown 167c4762a1bSJed Brown PetscErrorCode Monitor(SNES snes,PetscInt its,PetscReal rnorm,void *dummy) 168c4762a1bSJed Brown { 169c4762a1bSJed Brown UserCtx *user = (UserCtx*)dummy; 170c4762a1bSJed Brown Vec w,u,lambda,U,F; 171c4762a1bSJed Brown 172c4762a1bSJed Brown PetscFunctionBeginUser; 173*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetSolution(snes,&U)); 174*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCompositeGetAccess(user->packer,U,&w,&u,&lambda)); 175*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(u,user->u_viewer)); 176*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(lambda,user->lambda_viewer)); 177*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCompositeRestoreAccess(user->packer,U,&w,&u,&lambda)); 178c4762a1bSJed Brown 179*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetFunction(snes,&F,0,0)); 180*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCompositeGetAccess(user->packer,F,&w,&u,&lambda)); 181*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(u,user->fu_viewer)); 182*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(lambda,user->flambda_viewer)); 183*5f80ce2aSJacob Faibussowitsch CHKERRQ(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