[petsc-users] questions about dmdacreaterestriction
Xiangdong
epscodes at gmail.com
Mon May 7 20:11:42 CDT 2018
Hi Barry,
Thanks for your message. When I use DMCreateInterpolation, I got the error
message like "[0]PETSC ERROR: Invalid argument
[0]PETSC ERROR: Coarsening factor in x must be 2". The error will be gone
if I change refine_x = 3 to refine_x=2. Why must the coarsening factor be a
factor of 2?
Thank you.
Xiangdong
Below is my test code:
int main(int argc, char **argv)
{
PetscErrorCode ierr;
PetscInitialize(&argc, &argv, PETSC_NULL, PETSC_NULL);
int Nx=6, Ny=4, Nz=1;
DM da;
ierr = DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,
DMDA_STENCIL_STAR,Nx,Ny,PETSC_DECIDE,
PETSC_DECIDE,1,1,NULL,NULL,&da);CHKERRQ(ierr);
int refine_x=3, refine_y=2, refine_z=1;
ierr = DMDASetRefinementFactor(da,refine_x,refine_y,refine_z);
ierr = DMDASetInterpolationType(da,DMDA_Q0);
ierr = DMSetFromOptions(da); CHKERRQ(ierr);
ierr = DMSetUp(da); CHKERRQ(ierr);
DM dac;
ierr = DMCoarsen(da,NULL,&dac);
Mat M;
Vec V;
ierr = DMCreateInterpolation(dac,da,&M,&V);
MatView(M,PETSC_VIEWER_STDOUT_WORLD);
ierr = PetscFinalize(); CHKERRQ(ierr);
PetscFunctionReturn(0);
}
On Mon, May 7, 2018 at 6:16 PM, Smith, Barry F. <bsmith at mcs.anl.gov> wrote:
>
>
>
> > On May 7, 2018, at 3:25 PM, Xiangdong <epscodes at gmail.com> wrote:
> >
> > Hello everyone,
> >
> > I have a 2D DMDA with Nx=6 and Ny=4. If I set refine_x=3 and refine_y=2,
> I can get a coarse DMDA by DMCoarsen(). However, when I tried to call the
> DMCreateRestriction between these two DMDAs, I got error about " No support
> for this operation for this object type. DMCreateRestriction not
> implemented for this type."
> >
> > Can you help me get around this?
>
> Us DMCreateInterpolation() and then use MatMultTranspose() for example
> to apply the restriction.
>
> >
> > Another simple question, where can I find the explanation/meaning of
> DMDA_Q0 and DMDA_Q1?
>
> This is finite element terminology, Q0 means piecewise constant basis
> functions (and hence piecewise constant interpolation) while Q1 means
> linear interpolation.
> For Q0 you should think of the unknowns as living on the centers of each
> element while for Q1 the unknowns live on the vertices of the elements.
>
> Barry
>
> >
> > Below is my test example.
> >
> > Thank you.
> >
> > Best,
> > Xiangdong
> >
> > #include <petsc.h>
> >
> > int main(int argc, char **argv)
> > {
> > PetscErrorCode ierr;
> >
> > PetscInitialize(&argc, &argv, PETSC_NULL, PETSC_NULL);
> >
> > int Nx=6, Ny=4, Nz=1;
> >
> > DM da;
> >
> > ierr = DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,
> DMDA_STENCIL_STAR,Nx,Ny,PETSC_DECIDE, PETSC_DECIDE,1,1,NULL,NULL,&
> da);CHKERRQ(ierr);
> >
> > int refine_x=3, refine_y=2, refine_z=1;
> > ierr = DMDASetRefinementFactor(da,refine_x,refine_y,refine_z);
> >
> > ierr = DMDASetInterpolationType(da,DMDA_Q0);
> >
> > ierr = DMSetFromOptions(da); CHKERRQ(ierr);
> > ierr = DMSetUp(da); CHKERRQ(ierr);
> >
> > DM dac;
> > ierr = DMCoarsen(da,NULL,&dac);
> >
> > Mat M;
> > Vec V;
> >
> > ierr = DMCreateRestriction(dac,da,&M);
> >
> > MatView(M,PETSC_VIEWER_STDOUT_WORLD);
> >
> >
> > ierr = PetscFinalize(); CHKERRQ(ierr);
> > PetscFunctionReturn(0);
> > }
> >
>
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20180507/95f8803a/attachment-0001.html>
More information about the petsc-users
mailing list