/* Laplacian in 3D. Modeled by the partial differential equation - Laplacian u = 1,0 < x,y,z < 1, with boundary conditions u = 1 for x = 0, x = 1, y = 0, y = 1, z = 0, z = 1. This uses multigrid to solve the linear system */ static char help[] = "Solves 3D Laplacian using multigrid.\n\n"; #include "petscda.h" #include "petscksp.h" #include "petscdmmg.h" extern PetscErrorCode ComputeJacobian(DMMG,Mat,Mat); extern PetscErrorCode ComputeRHS(DMMG,Vec); #undef __FUNCT__ #define __FUNCT__ "main" int main(int argc,char **argv) { PetscErrorCode ierr; DMMG *dmmg; PetscReal norm; DA da; PetscInitialize(&argc,&argv,(char *)0,help); ierr = DMMGCreate(PETSC_COMM_WORLD,3,PETSC_NULL,&dmmg);CHKERRQ(ierr); ierr = DACreate3d(PETSC_COMM_WORLD,DA_NONPERIODIC,DA_STENCIL_STAR,-3,-3,-3,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,1,1,0,0,0,&da);CHKERRQ(ierr); ierr = DMMGSetDM(dmmg,(DM)da);CHKERRQ(ierr); ierr = DADestroy(da);CHKERRQ(ierr); ierr = DMMGSetKSP(dmmg,ComputeRHS,ComputeJacobian);CHKERRQ(ierr); ierr = DMMGSolve(dmmg);CHKERRQ(ierr); ierr = MatMult(DMMGGetJ(dmmg),DMMGGetx(dmmg),DMMGGetr(dmmg));CHKERRQ(ierr); ierr = VecAXPY(DMMGGetr(dmmg),-1.0,DMMGGetRHS(dmmg));CHKERRQ(ierr); ierr = VecNorm(DMMGGetr(dmmg),NORM_2,&norm);CHKERRQ(ierr); /* ierr = PetscPrintf(PETSC_COMM_WORLD,"Residual norm %G\n",norm);CHKERRQ(ierr); */ ierr = DMMGDestroy(dmmg);CHKERRQ(ierr); ierr = PetscFinalize();CHKERRQ(ierr); return 0; } #undef __FUNCT__ #define __FUNCT__ "ComputeRHS" PetscErrorCode ComputeRHS(DMMG dmmg,Vec b) { PetscErrorCode ierr; PetscInt mx,my,mz; PetscScalar h; PetscFunctionBegin; ierr = DAGetInfo((DA)dmmg->dm,0,&mx,&my,&mz,0,0,0,0,0,0,0);CHKERRQ(ierr); h = 1.0/((mx-1)*(my-1)*(mz-1)); ierr = VecSet(b,h);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "ComputeJacobian" PetscErrorCode ComputeJacobian(DMMG dmmg,Mat jac,Mat B) { DA da = (DA)dmmg->dm; PetscErrorCode ierr; PetscInt i,j,k,mx,my,mz,xm,ym,zm,xs,ys,zs; PetscScalar v[7],Hx,Hy,Hz,HxHydHz,HyHzdHx,HxHzdHy; MatStencil row,col[7]; ierr = DAGetInfo(da,0,&mx,&my,&mz,0,0,0,0,0,0,0);CHKERRQ(ierr); Hx = 1.0 / (PetscReal)(mx-1); Hy = 1.0 / (PetscReal)(my-1); Hz = 1.0 / (PetscReal)(mz-1); HxHydHz = Hx*Hy/Hz; HxHzdHy = Hx*Hz/Hy; HyHzdHx = Hy*Hz/Hx; ierr = DAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr); for (k=zs; k