[petsc-users] questions about dmdacreaterestriction

Smith, Barry F. bsmith at mcs.anl.gov
Mon May 7 17:16:33 CDT 2018


   

> 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);
> }
> 



More information about the petsc-users mailing list