[petsc-users] DMMG with PBC
khalid ashraf
khalid_eee at yahoo.com
Tue Apr 19 18:26:00 CDT 2011
>>How much of a difference?
With the applied XYZPeriodic, if I keep the follwoing lines
if (i==0 || j==0 || k==0 || i==mx-1 || j==my-1 || k==mz-1){
v[0] = 2.0*(HxHydHz + HxHzdHy + HyHzdHx);
ierr =
MatSetValuesStencil(B,1,&row,1,&row,v,INSERT_VALUES);CHKERRQ(ierr);
} else
Then the error between 1proc and 4 procs is only after 5th decimal point.
However, if I comment out the
above lines then the results are completely different for 1 and 4 procs.
I am attaching the output of some last data points of a 10X10X8 grid.
1 proc output:
-4.40214
-4.39202
-4.38693
-4.38547
-4.38687
-4.39047
4 proc output:
0.000188031
0.000169784
0.000157229
0.000178713
0.000179637
0.000188031
0.000169784
0.000157229
0.000178713
0.000179637
I am attaching the faulty code here for your review.
Thanks.
Khalid
static char help[] = "Solves 3D Laplacian using multigrid.\n\n";
#include "petscda.h"
#include "petscksp.h"
#include "petscdmmg.h"
#include "myHeaderfile.h"
extern PetscErrorCode ComputeMatrix(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,1,PETSC_NULL,&dmmg);CHKERRQ(ierr);
ierr =
DACreate3d(PETSC_COMM_WORLD,DA_XYZPERIODIC,DA_STENCIL_STAR,10,10,8,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,ComputeMatrix);CHKERRQ(ierr);
ierr =DMMGSetNullSpace(dmmg,PETSC_TRUE,0,PETSC_NULL);
ierr = DMMGSetUp(dmmg);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=VecView_VTK(DMMGGetx(dmmg),"X",&appctx);
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 = 10.0/((mx-1)*(my-1)*(mz-1));
ierr = VecSet(b,h);CHKERRQ(ierr);
PetscFunctionReturn(0);
}
#undef __FUNCT__
#define __FUNCT__ "ComputeMatrix"
PetscErrorCode ComputeMatrix(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<zs+zm; k++){
for (j=ys; j<ys+ym; j++){
for(i=xs; i<xs+xm; i++){
row.i = i; row.j = j; row.k = k;
/* if (i==0 || j==0 || k==0 || i==mx-1 || j==my-1 || k==mz-1){
v[0] = 2.0*(HxHydHz + HxHzdHy + HyHzdHx);
ierr =
MatSetValuesStencil(B,1,&row,1,&row,v,INSERT_VALUES);CHKERRQ(ierr);
} else */
{
v[0] = -HxHydHz;col[0].i = i; col[0].j = j; col[0].k = k-1;
v[1] = -HxHzdHy;col[1].i = i; col[1].j = j-1; col[1].k = k;
v[2] = -HyHzdHx;col[2].i = i-1; col[2].j = j; col[2].k = k;
v[3] = 2.0*(HxHydHz + HxHzdHy + HyHzdHx);col[3].i = row.i; col[3].j =
row.j; col[3].k = row.k;
v[4] = -HyHzdHx;col[4].i = i+1; col[4].j = j; col[4].k = k;
v[5] = -HxHzdHy;col[5].i = i; col[5].j = j+1; col[5].k = k;
v[6] = -HxHydHz;col[6].i = i; col[6].j = j; col[6].k = k+1;
ierr = MatSetValuesStencil(B,1,&row,7,col,v,INSERT_VALUES);CHKERRQ(ierr);
}
}
}
}
ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
return 0;
}
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20110419/d2e304ad/attachment.htm>
More information about the petsc-users
mailing list