[petsc-users] How to access solution in /src/ksp/ksp/example/tutorial/ex45.c

Zhenglun (Alan) Wei zhenglun.wei at gmail.com
Mon Apr 9 17:37:45 CDT 2012


Dear All,
     I hope you're having a nice day.
     I'm working on ex45. It does not have any output subroutine inside. 
Therefore, I was trying to write out mine.
     I was trying to use DMDAVecGetArray; but I found that the DM object 
has been destroyed at very beginning of the main code. Therefore, I 
swithed to VecGetArray. However, it gives me warning like:
ex45.c:185: warning: passing argument 2 of 'VecGetArray' from 
incompatible pointer type
     Then, it shows a memory violation error during execution, which is 
attached as 'out'. I checked 'VecGetArray' examples in PETSc manual, it 
always shows how to use it to extract one-dimensional array. I wonder if 
it can do 3-dimensional?

     Also, there are other two minor questions:
1) in ex45.c, 'ComputeRHS', there is a if(x) {} statement, I write probe 
statements inside this if statement to see if the code ever go through 
it; the result shows it does not. I don't know why. ;(
2) probe statements are added in 'ComputeRHS' for both ex45.c and 
ex29.c. The size parameter used for DMDACreate3d in ex45 is -7, -7, -7; 
while for DMDACreate2d in ex29 is -7,-7. The probe statements show that 
the real 2d domain generated in ex29 is 25*25, which, as I understand, 
is because of multigrid. However, with multigrid scheme, the ex45 create 
a 7*7*7 3d domain. Any suggestions on this?

     I attached both ex45 and ex29 with this E-mail.

thanks in advance,
Alan

-------------- next part --------------

/*
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

   The same as ex22.c except it does not use DMMG, it uses its replacement.
   See src/snes/examples/tutorials/ex50.c

   Can also be run with -pc_type exotic -ksp_type fgmres

*/

static char help[] = "Solves 3D Laplacian using multigrid.\n\n";

#include <petscksp.h>
#include <petscdmda.h>
PetscInt i,j,k,mx,my,mz,xm,ym,zm,xs,ys,zs;

extern PetscErrorCode ComputeMatrix(DM,Vec,Mat,Mat,MatStructure*);
extern PetscErrorCode ComputeRHS(DM,Vec,Vec);
extern PetscErrorCode ComputeInitialGuess(DM, Vec);
extern PetscErrorCode OutputFunction(Vec);
#undef __FUNCT__
#define __FUNCT__ "main"
int main(int argc,char **argv)
{
  PetscErrorCode ierr;
  KSP            ksp;
  PetscReal      norm;
  DM             da;
  Vec            x,b,r;
  Mat            A;

  PetscInitialize(&argc,&argv,(char *)0,help);

  ierr = KSPCreate(PETSC_COMM_WORLD,&ksp);CHKERRQ(ierr);
  ierr = DMDACreate3d(PETSC_COMM_WORLD,DMDA_BOUNDARY_NONE,DMDA_BOUNDARY_NONE,DMDA_BOUNDARY_NONE,DMDA_STENCIL_STAR,-7,-7,-7,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,1,1,0,0,0,&da);CHKERRQ(ierr);  
  ierr = DMSetInitialGuess(da,ComputeInitialGuess);CHKERRQ(ierr);
  ierr = DMSetFunction(da,ComputeRHS);CHKERRQ(ierr);
  ierr = DMSetJacobian(da,ComputeMatrix);CHKERRQ(ierr);
  ierr = KSPSetDM(ksp,da);CHKERRQ(ierr);
  /*  ierr = KSPSetDMActive(ksp,PETSC_FALSE);CHKERRQ(ierr);*/
  ierr = DMDestroy(&da);CHKERRQ(ierr);

  ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr);
  ierr = KSPSetUp(ksp);CHKERRQ(ierr);
  ierr = KSPSolve(ksp,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
  ierr = KSPGetSolution(ksp,&x);CHKERRQ(ierr);
  ierr = KSPGetRhs(ksp,&b);CHKERRQ(ierr);
  ierr = VecDuplicate(b,&r);CHKERRQ(ierr);
  ierr = KSPGetOperators(ksp,&A,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);

  ierr = MatMult(A,x,r);CHKERRQ(ierr);
  ierr = VecAXPY(r,-1.0,b);CHKERRQ(ierr);
  ierr = VecNorm(r,NORM_2,&norm);CHKERRQ(ierr);
  ierr = PetscPrintf(PETSC_COMM_WORLD,"Residual norm %G\n",norm);CHKERRQ(ierr); 
  ierr = OutputFunction(x);

  ierr = VecDestroy(&r);CHKERRQ(ierr);
  ierr = KSPDestroy(&ksp);CHKERRQ(ierr);
  ierr = PetscFinalize();

  return 0;
}

#undef __FUNCT__
#define __FUNCT__ "ComputeRHS"
PetscErrorCode ComputeRHS(DM dm,Vec x,Vec b)
{
  PetscErrorCode ierr;
  PetscInt       mx,my,mz;
  PetscScalar    h;
  PetscScalar    ***xx,***bb;
  PetscScalar    Hx,Hy,Hz,HxHydHz,HyHzdHx,HxHzdHy;

  PetscFunctionBegin;
  ierr = DMDAGetInfo(dm,0,&mx,&my,&mz,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
  printf("mx = %d, my = %d, mz = %d\n", mx, my, mz);
  printf("xs = %d, ys = %d, zs = %d\n", xs, ys, zs);
  printf("xm = %d, ym = %d, zm = %d\n", xm, ym, zm);

  h    = 1.0/((mx-1)*(my-1)*(mz-1));
  ierr = VecSet(b,h);CHKERRQ(ierr);

  if (x) {
    ierr = DMDAVecGetArray(dm,x,&xx);CHKERRQ(ierr);
    ierr = DMDAVecGetArray(dm,b,&bb);CHKERRQ(ierr);

  ierr = DMDAGetInfo(dm,0,&mx,&my,&mz,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr); 
  printf("mx = %d, my = %d, mz = %d\n", mx, my, mz);
  printf("xs = %d, ys = %d, zs = %d\n", xs, ys, zs);
  printf("xm = %d, ym = %d, zm = %d\n", xm, ym, zm);
 
  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 = DMDAGetCorners(dm,&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++){
	if (!(i==0 || j==0 || k==0 || i==mx-1 || j==my-1 || k==mz-1)){
          bb[k][j][i] -= +HxHydHz*xx[k-1][j][i] + HxHzdHy*xx[k][j-1][i] + HyHzdHx*xx[k][j][i-1] + HyHzdHx*xx[k][j][i+1] + HxHzdHy*xx[k][j+1][i] + HxHydHz*xx[k+1][j][i];
        }
        bb[k][j][i] += 2.0*(HxHydHz + HxHzdHy + HyHzdHx)*xx[k][j][i];
      }
    }
  }
    ierr = DMDAVecRestoreArray(dm,x,&xx);CHKERRQ(ierr);
    ierr = DMDAVecRestoreArray(dm,b,&bb);CHKERRQ(ierr);
  }
  PetscFunctionReturn(0);
}
    
#undef __FUNCT__
#define __FUNCT__ "ComputeInitialGuess"
PetscErrorCode ComputeInitialGuess(DM dm,Vec b)
{
  PetscErrorCode ierr;

  PetscFunctionBegin;
  ierr = VecSet(b,0);CHKERRQ(ierr);
  PetscFunctionReturn(0);
}

#undef __FUNCT__
#define __FUNCT__ "ComputeMatrix"
PetscErrorCode ComputeMatrix(DM dm,Vec x,Mat jac,Mat B,MatStructure *stflg)
{
  DM             da = dm;
  PetscErrorCode ierr;
  PetscScalar    v[7],Hx,Hy,Hz,HxHydHz,HyHzdHx,HxHzdHy;
  MatStencil     row,col[7];

  ierr = DMDAGetInfo(da,0,&mx,&my,&mz,0,0,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 = DMDAGetCorners(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);
  *stflg = SAME_NONZERO_PATTERN;
  return 0;
}

#undef __FUNCT__
#define __FUNCT__ "OutputFunction"
PetscErrorCode OutputFunction(Vec x) 
{
  PetscErrorCode ierr;
  PetscScalar ***SolArray;
  int rank, countt, countn;
  char fname[20];
  FILE *fgrid, *fpplt3d;

  PetscFunctionBegin;

  if(x) { 
    ierr = VecGetArray(x, &SolArray);CHKERRQ(ierr);
  }

  printf("mx = %d, my = %d, mz = %d\n", mx, my, mz);
  printf("xs = %d, ys = %d, zs = %d\n", xs, ys, zs);
  printf("xm = %d, ym = %d, zm = %d\n", xm, ym, zm);

  MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
  sprintf(fname, "%s%d%s", "Rank", rank,"plt3d.q");
  fpplt3d = fopen(fname, "w");

  fprintf(fpplt3d, "%d  %d  %d\n", mx, my, mz);
  fprintf(fpplt3d, "%f  %f  %f  %f  %f  \n", 1.0, 1.0, 1.0, 1.0, 1.0);

  for(countn = 1; countn <= 4; countn++) {
    for(k = zs; k < zs+zm; k++) {
      for(j = ys; j < ys+ym; j++) {
	for(i = xs; i < xs+xm; i++) {
	  if(countn == 1) fprintf(fpplt3d, "%f  ", 1.0);
	  if(countn == 2) fprintf(fpplt3d, "%f  ", SolArray[k][j][i]);
	  if(countn == 3) fprintf(fpplt3d, "%f  ", SolArray[k][j][i]);
	  if(countn == 4) fprintf(fpplt3d, "%f  ", SolArray[k][j][i]);
	  if(countn == 5) fprintf(fpplt3d, "%f  ", SolArray[k][j][i]);
	  if(!(countt%6)) fprintf(fpplt3d, "\n");
	  countt++;
	}
      }
    }
  }
  fprintf(fpplt3d, "\n");
  fclose(fpplt3d);


  PetscFunctionReturn(0);
}
-------------- next part --------------
mx = 7, my = 7, mz = 7
xs = 0, ys = 0, zs = 0
xm = 0, ym = 0, zm = 0
Residual norm 7.06789e-07
mx = 7, my = 7, mz = 7
xs = 0, ys = 0, zs = 0
xm = 7, ym = 7, zm = 7
[0]PETSC ERROR: ------------------------------------------------------------------------
[0]PETSC ERROR: Caught signal number 11 SEGV: Segmentation Violation, probably memory access out of range
[0]PETSC ERROR: Try option -start_in_debugger or -on_error_attach_debugger
[0]PETSC ERROR: or see http://www.mcs.anl.gov/petsc/petsc-as/documentation/faq.html#valgrind[0]PETSC ERROR: or try http://valgrind.org on GNU/linux and Apple Mac OS X to find memory corruption errors
[0]PETSC ERROR: likely location of problem given in stack below
[0]PETSC ERROR: ---------------------  Stack Frames ------------------------------------
[0]PETSC ERROR: Note: The EXACT line numbers in the stack are not available,
[0]PETSC ERROR:       INSTEAD the line number of the start of the function
[0]PETSC ERROR:       is given.
[0]PETSC ERROR: [0] OutputFunction line 182 src/ksp/ksp/examples/tutorials/ex45.c
[0]PETSC ERROR: --------------------- Error Message ------------------------------------
[0]PETSC ERROR: Signal received!
[0]PETSC ERROR: ------------------------------------------------------------------------
[0]PETSC ERROR: Petsc Development HG revision: 85e6055943e0711fccdec2d08caeba48971d3d55  HG Date: Fri Sep 23 14:46:14 2011 -0700
[0]PETSC ERROR: See docs/changes/index.html for recent updates.
[0]PETSC ERROR: See docs/faq.html for hints about trouble shooting.
[0]PETSC ERROR: See docs/index.html for manual pages.
[0]PETSC ERROR: ------------------------------------------------------------------------
[0]PETSC ERROR: ./ex45 on a arch-linu named l2118a-linux.soecs.ku.edu by zlwei Mon Apr  9 18:19:47 2012
[0]PETSC ERROR: Libraries linked from /home/zlwei/soft/mercurial/petsc-dev/arch-linux2-c-debug/lib
[0]PETSC ERROR: Configure run at Fri Sep 23 17:13:32 2011
[0]PETSC ERROR: Configure options --download-f-blas-lapack=1 --download-mpich=1 --with-cc=gcc --with-fc=gfortran PETSC_ARCH=arch-linux2-c-debug
[0]PETSC ERROR: ------------------------------------------------------------------------
[0]PETSC ERROR: User provided function() line 0 in unknown directory unknown file
application called MPI_Abort(MPI_COMM_WORLD, 59) - process 0
[unset]: aborting job:
application called MPI_Abort(MPI_COMM_WORLD, 59) - process 0
-------------- next part --------------
/*T
   Concepts: KSP^solving a system of linear equations
   Concepts: KSP^Laplacian, 2d
   Processors: n
T*/

/*
Added at the request of Marc Garbey.

Inhomogeneous Laplacian in 2D. Modeled by the partial differential equation

   -div \rho grad u = f,  0 < x,y < 1,

with forcing function

   f = e^{-x^2/\nu} e^{-y^2/\nu}

with Dirichlet boundary conditions

   u = f(x,y) for x = 0, x = 1, y = 0, y = 1

or pure Neumman boundary conditions

This uses multigrid to solve the linear system
*/

static char help[] = "Solves 2D inhomogeneous Laplacian using multigrid.\n\n";

#include <petscdmda.h>
#include <petscksp.h>
#include <petscpcmg.h>
#include <petscdmmg.h>

extern PetscErrorCode ComputeMatrix(DMMG,Mat,Mat);
extern PetscErrorCode ComputeRHS(DMMG,Vec);
extern PetscErrorCode VecView_VTK(Vec, const char [], const char []);

typedef enum {DIRICHLET, NEUMANN} BCType;

typedef struct {
  PetscScalar   rho;
  PetscScalar   nu;
  BCType        bcType;
} UserContext;

#undef __FUNCT__
#define __FUNCT__ "main"
int main(int argc,char **argv)
{
  DMMG           *dmmg;
  DM             da;
  UserContext    user;
  PetscReal      norm;
  const char     *bcTypes[2] = {"dirichlet","neumann"};
  PetscErrorCode ierr;
  PetscInt       l,bc;

  PetscInitialize(&argc,&argv,(char *)0,help);

  ierr = DMMGCreate(PETSC_COMM_WORLD,3,PETSC_NULL,&dmmg);CHKERRQ(ierr);
  ierr = DMDACreate2d(PETSC_COMM_WORLD, DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_NONE,DMDA_STENCIL_STAR,-7,-7,PETSC_DECIDE,PETSC_DECIDE,1,1,0,0,&da);CHKERRQ(ierr);  
  ierr = DMMGSetDM(dmmg,(DM)da);CHKERRQ(ierr);
  ierr = DMDestroy(&da);CHKERRQ(ierr);
  for (l = 0; l < DMMGGetLevels(dmmg); l++) {
    ierr = DMMGSetUser(dmmg,l,&user);CHKERRQ(ierr);
  }

  ierr = PetscOptionsBegin(PETSC_COMM_WORLD, "", "Options for the inhomogeneous Poisson equation", "DMMG");
    user.rho    = 1.0;
    ierr        = PetscOptionsScalar("-rho", "The conductivity", "ex29.c", user.rho, &user.rho, PETSC_NULL);CHKERRQ(ierr);
    user.nu     = 0.1;
    ierr        = PetscOptionsScalar("-nu", "The width of the Gaussian source", "ex29.c", user.nu, &user.nu, PETSC_NULL);CHKERRQ(ierr);
    bc          = (PetscInt)DIRICHLET;
    ierr        = PetscOptionsEList("-bc_type","Type of boundary condition","ex29.c",bcTypes,2,bcTypes[0],&bc,PETSC_NULL);CHKERRQ(ierr);
    user.bcType = (BCType)bc;
  ierr = PetscOptionsEnd();

  ierr = DMMGSetKSP(dmmg,ComputeRHS,ComputeMatrix);CHKERRQ(ierr);
  if (user.bcType == NEUMANN) {
    ierr = DMMGSetNullSpace(dmmg,PETSC_TRUE,0,PETSC_NULL);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 = VecAssemblyBegin(DMMGGetx(dmmg));CHKERRQ(ierr);
  ierr = VecAssemblyEnd(DMMGGetx(dmmg));CHKERRQ(ierr);
//  ierr = VecView_VTK(DMMGGetRHS(dmmg), "rhs.vtk", bcTypes[user.bcType]);CHKERRQ(ierr);
//  ierr = VecView_VTK(DMMGGetx(dmmg), "solution.vtk", bcTypes[user.bcType]);CHKERRQ(ierr);

  ierr = DMMGDestroy(dmmg);CHKERRQ(ierr);
  ierr = PetscFinalize();

  return 0;
}

#undef __FUNCT__
#define __FUNCT__ "ComputeRHS"
PetscErrorCode ComputeRHS(DMMG dmmg, Vec b)
{
  DM             da = dmmg->dm;
  UserContext    *user = (UserContext *) dmmg->user;
  PetscErrorCode ierr;
  PetscInt       i,j,mx,my,xm,ym,xs,ys;
  PetscScalar    Hx,Hy;
  PetscScalar    **array;

  PetscFunctionBegin;
  ierr = DMDAGetInfo(da, 0, &mx, &my, 0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
  Hx   = 1.0 / (PetscReal)(mx-1);
  Hy   = 1.0 / (PetscReal)(my-1);
  ierr = DMDAGetCorners(da,&xs,&ys,0,&xm,&ym,0);CHKERRQ(ierr);
  printf("%d  %d  %d  %d  %d  %d  \n", mx, my, xs, ys, xm, ym);
  ierr = DMDAVecGetArray(da, b, &array);CHKERRQ(ierr);
  for (j=ys; j<ys+ym; j++){
    for(i=xs; i<xs+xm; i++){
      array[j][i] = PetscExpScalar(-((PetscReal)i*Hx)*((PetscReal)i*Hx)/user->nu)*PetscExpScalar(-((PetscReal)j*Hy)*((PetscReal)j*Hy)/user->nu)*Hx*Hy;
    }
  }
  ierr = DMDAVecRestoreArray(da, b, &array);CHKERRQ(ierr);
  ierr = VecAssemblyBegin(b);CHKERRQ(ierr);
  ierr = VecAssemblyEnd(b);CHKERRQ(ierr);

  /* force right hand side to be consistent for singular matrix */
  /* note this is really a hack, normally the model would provide you with a consistent right handside */
  if (user->bcType == NEUMANN) {
    MatNullSpace nullspace;

    ierr = KSPGetNullSpace(dmmg->ksp,&nullspace);CHKERRQ(ierr);
    ierr = MatNullSpaceRemove(nullspace,b,PETSC_NULL);CHKERRQ(ierr);
  }
  PetscFunctionReturn(0);
}

    
#undef __FUNCT__
#define __FUNCT__ "ComputeRho"
PetscErrorCode ComputeRho(PetscInt i, PetscInt j, PetscInt mx, PetscInt my, PetscScalar centerRho, PetscScalar *rho)
{
  PetscFunctionBegin;
  if ((i > mx/3.0) && (i < 2.0*mx/3.0) && (j > my/3.0) && (j < 2.0*my/3.0)) {
    *rho = centerRho;
  } else {
    *rho = 1.0;
  }
  PetscFunctionReturn(0);
}

#undef __FUNCT__
#define __FUNCT__ "ComputeMatrix"
PetscErrorCode ComputeMatrix(DMMG dmmg, Mat J,Mat jac)
{
  DM             da =  dmmg->dm;
  UserContext    *user = (UserContext *) dmmg->user;
  PetscScalar    centerRho = user->rho;
  PetscErrorCode ierr;
  PetscInt       i,j,mx,my,xm,ym,xs,ys,num;
  PetscScalar    v[5],Hx,Hy,HydHx,HxdHy,rho;
  MatStencil     row, col[5];

  PetscFunctionBegin;
  ierr = DMDAGetInfo(da,0,&mx,&my,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);  
  Hx    = 1.0 / (PetscReal)(mx-1);
  Hy    = 1.0 / (PetscReal)(my-1);
  HxdHy = Hx/Hy;
  HydHx = Hy/Hx;
  ierr = DMDAGetCorners(da,&xs,&ys,0,&xm,&ym,0);CHKERRQ(ierr);
  for (j=ys; j<ys+ym; j++){
    for(i=xs; i<xs+xm; i++){
      row.i = i; row.j = j;
      ierr = ComputeRho(i, j, mx, my, centerRho, &rho);CHKERRQ(ierr);
      if (i==0 || j==0 || i==mx-1 || j==my-1) {
        if (user->bcType == DIRICHLET) {
           v[0] = 2.0*rho*(HxdHy + HydHx);
          ierr = MatSetValuesStencil(jac,1,&row,1,&row,v,INSERT_VALUES);CHKERRQ(ierr);
        } else if (user->bcType == NEUMANN) {
          num = 0;
          if (j!=0) {
            v[num] = -rho*HxdHy;              col[num].i = i;   col[num].j = j-1;
            num++;
          }
          if (i!=0) {
            v[num] = -rho*HydHx;              col[num].i = i-1; col[num].j = j;
            num++;
          }
          if (i!=mx-1) {
            v[num] = -rho*HydHx;              col[num].i = i+1; col[num].j = j;
            num++;
          }
          if (j!=my-1) {
            v[num] = -rho*HxdHy;              col[num].i = i;   col[num].j = j+1;
            num++;
          }
          v[num]   = (num/2.0)*rho*(HxdHy + HydHx); col[num].i = i;   col[num].j = j;
          num++;
          ierr = MatSetValuesStencil(jac,1,&row,num,col,v,INSERT_VALUES);CHKERRQ(ierr);
        }
      } else {
        v[0] = -rho*HxdHy;              col[0].i = i;   col[0].j = j-1;
        v[1] = -rho*HydHx;              col[1].i = i-1; col[1].j = j;
        v[2] = 2.0*rho*(HxdHy + HydHx); col[2].i = i;   col[2].j = j;
        v[3] = -rho*HydHx;              col[3].i = i+1; col[3].j = j;
        v[4] = -rho*HxdHy;              col[4].i = i;   col[4].j = j+1;
        ierr = MatSetValuesStencil(jac,1,&row,5,col,v,INSERT_VALUES);CHKERRQ(ierr);
      }
    }
  }
  ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
  ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
  PetscFunctionReturn(0);
}
/*
#undef __FUNCT__
#define __FUNCT__ "VecView_VTK"
PetscErrorCode VecView_VTK(Vec x, const char filename[], const char bcName[])
{
  MPI_Comm           comm;
  DM                 da;
  Vec                coords;
  PetscViewer        viewer;
  PetscScalar       *array, *values;
  PetscInt           n, N, maxn, mx, my, dof;
  PetscInt           i, p;
  MPI_Status         status;
  PetscMPIInt        rank, size, tag, nn;
  PetscErrorCode     ierr;
  PetscBool          flg;

  PetscFunctionBegin;
  ierr = PetscObjectGetComm((PetscObject) x, &comm);CHKERRQ(ierr);
  ierr = PetscViewerASCIIOpen(comm, filename, &viewer);CHKERRQ(ierr);

  ierr = VecGetSize(x, &N);CHKERRQ(ierr);
  ierr = VecGetLocalSize(x, &n);CHKERRQ(ierr);
  ierr = PetscObjectQuery((PetscObject) x, "DM", (PetscObject *) &da);CHKERRQ(ierr);
  if (!da) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
  ierr = PetscTypeCompare((PetscObject)da,DMDA,&flg);CHKERRQ(ierr);
  if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");

  ierr = DMDAGetInfo(da, 0, &mx, &my, 0,0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr);

  ierr = PetscViewerASCIIPrintf(viewer, "# vtk DataFile Version 2.0\n");CHKERRQ(ierr);
  ierr = PetscViewerASCIIPrintf(viewer, "Inhomogeneous Poisson Equation with %s boundary conditions\n", bcName);CHKERRQ(ierr);
  ierr = PetscViewerASCIIPrintf(viewer, "ASCII\n");CHKERRQ(ierr);*/
  /* get coordinates of nodes */
/*  ierr = DMDAGetCoordinates(da, &coords);CHKERRQ(ierr);
  if (!coords) {
    ierr = DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 0.0);CHKERRQ(ierr);
    ierr = DMDAGetCoordinates(da, &coords);CHKERRQ(ierr);
  }
  ierr = PetscViewerASCIIPrintf(viewer, "DATASET RECTILINEAR_GRID\n");CHKERRQ(ierr);
  ierr = PetscViewerASCIIPrintf(viewer, "DIMENSIONS %d %d %d\n", mx, my, 1);CHKERRQ(ierr);
  ierr = VecGetArray(coords, &array);CHKERRQ(ierr);
  ierr = PetscViewerASCIIPrintf(viewer, "X_COORDINATES %d double\n", mx);CHKERRQ(ierr);
  for(i = 0; i < mx; i++) {
    ierr = PetscViewerASCIIPrintf(viewer, "%G ", PetscRealPart(array[i*2]));CHKERRQ(ierr);
  }
  ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr);
  ierr = PetscViewerASCIIPrintf(viewer, "Y_COORDINATES %d double\n", my);CHKERRQ(ierr);
  for(i = 0; i < my; i++) {
    ierr = PetscViewerASCIIPrintf(viewer, "%G ", PetscRealPart(array[i*mx*2+1]));CHKERRQ(ierr);
  }
  ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr);
  ierr = PetscViewerASCIIPrintf(viewer, "Z_COORDINATES %d double\n", 1);CHKERRQ(ierr);
  ierr = PetscViewerASCIIPrintf(viewer, "%G\n", 0.0);CHKERRQ(ierr);
  ierr = VecRestoreArray(coords, &array);CHKERRQ(ierr);
  ierr = PetscViewerASCIIPrintf(viewer, "POINT_DATA %d\n", N);CHKERRQ(ierr);
  ierr = PetscViewerASCIIPrintf(viewer, "SCALARS scalars double %d\n", dof);CHKERRQ(ierr);
  ierr = PetscViewerASCIIPrintf(viewer, "LOOKUP_TABLE default\n");CHKERRQ(ierr);
  ierr = VecGetArray(x, &array);CHKERRQ(ierr);
  /* Determine maximum message to arrive */
/*  ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
  ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
  ierr = MPI_Reduce(&n, &maxn, 1, MPIU_INT, MPI_MAX, 0, comm);CHKERRQ(ierr);
  tag  = ((PetscObject) viewer)->tag;
  if (!rank) {
    ierr = PetscMalloc((maxn+1) * sizeof(PetscScalar), &values);CHKERRQ(ierr);
    for(i = 0; i < n; i++) {
      ierr = PetscViewerASCIIPrintf(viewer, "%G\n", PetscRealPart(array[i]));CHKERRQ(ierr);
    }
    for(p = 1; p < size; p++) {
      ierr = MPI_Recv(values, (PetscMPIInt) n, MPIU_SCALAR, p, tag, comm, &status);CHKERRQ(ierr);
      ierr = MPI_Get_count(&status, MPIU_SCALAR, &nn);CHKERRQ(ierr);        
      for(i = 0; i < nn; i++) {
        ierr = PetscViewerASCIIPrintf(viewer, "%G\n", PetscRealPart(array[i]));CHKERRQ(ierr);
      }
    }
    ierr = PetscFree(values);CHKERRQ(ierr);
  } else {
    ierr = MPI_Send(array, n, MPIU_SCALAR, 0, tag, comm);CHKERRQ(ierr);
  }
  ierr = VecRestoreArray(x, &array);CHKERRQ(ierr);
  ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
  ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
  PetscFunctionReturn(0);
}*/


More information about the petsc-users mailing list