[petsc-users] DMRefine issue

Matthew Knepley knepley at gmail.com
Sun Feb 16 13:52:46 CST 2014


On Mon, Dec 30, 2013 at 7:23 PM, Yaakoub El Khamra
<yaakoub at tacc.utexas.edu>wrote:

>
> Greetings
> I am try to use DMRefine to refine a dm I create with DMPlexCreateBoxMesh
> with tonight's checkout. A quick glance at the DMView output from the DM's
> shows that the refinedMesh looks a lot like the original
> dm. DMGetRefineLevel returns 1 though. I must be missing something
> elementary here. Can someone please confirm if there are any glaring
> mistakes I made? Also if I want to refine and distribute, should I refine
> on a single processor then distribute or can I refine a couple of times,
> distribute, and refine some more?
>

1) Sorry this took a long time.

2) Refinement is somewhat complicated right now. If you use a mesh
generator to make the mesh (as BoxMesh does), then
    the default refinement is to use the generator itself. However, then
you must call

        ierr = DMPlexSetRefinementLimit(dm, refinementLimit);CHKERRQ(ierr);

   to set the maximum cell volume after refinement. If, however, the mesh
is constructed by hand the default is to use uniform
   refinement. You can also get uniform refinement using

      ierr = DMPlexSetRefinementUniform(dm,
refinementUniform);CHKERRQ(ierr);

3) Only uniform refinement works in parallel, so I tend to refine using the
generator in serial, distribute, and then uniformly refine.
     This is what SNES ex12 does.

4) Did you have another problem that I was supposed to look at on your
machine? I have lost my notes.

  Thanks,

     Matt


> The test code looks as follows:
>
>       program testRefine
>       implicit none
> #include "finclude/petsc.h90"
>       DM :: dm
>       DM :: refinedMesh
>       PetscInt :: r
>
>       PetscInt :: dim
>       PetscBool :: interpolate
>       PetscErrorCode :: ierr
>
>       call PetscInitialize(PETSC_NULL_CHARACTER, ierr)
>       CHKERRQ(ierr)
>       dim = 2
>       interpolate = PETSC_TRUE
> !     Create a mesh
>       call DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, interpolate, dm,
>      $     ierr)
>       CHKERRQ(ierr)
>
>       call DMView(dm, PETSC_VIEWER_STDOUT_WORLD, ierr)
>       CHKERRQ(ierr)
>
>
>       call DMRefine(dm, PETSC_COMM_WORLD, refinedMesh, ierr)
>       CHKERRQ(ierr)
>
>       call DMView(refinedMesh, PETSC_VIEWER_STDOUT_WORLD, ierr)
>       CHKERRQ(ierr)
>
>       call DMGetRefineLevel(refinedMesh, r, ierr)
>       write (*,*) "The value of r is:", r
>
>       call DMDestroy(dm, ierr)
>       call DMDestroy(refinedMesh, ierr)
>       CHKERRQ(ierr)
>
>       call PetscFinalize(ierr)
>       end program testRefine
>
>
>
> And the output:
> DM_0x1122bb0_0 in 2 dimensions:
>   0-cells: 9
>   1-cells: 16
>   2-cells: 8
> Labels:
>   marker: 2 strata of sizes (16, 5)
>   depth: 3 strata of sizes (9, 16, 8)
> DM_0x1122bb0_1 in 2 dimensions:
>   0-cells: 9
>   1-cells: 16
>   2-cells: 8
> Labels:
>   marker: 2 strata of sizes (16, 1)
>   depth: 3 strata of sizes (9, 16, 8)
>  The value of r is:           1
>
>
>
>
>
>
>
> Regards
> Yaakoub El Khamra
>



-- 
What most experimenters take for granted before they begin their
experiments is infinitely more interesting than any results to which their
experiments lead.
-- Norbert Wiener
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20140216/271a0bcf/attachment.html>


More information about the petsc-users mailing list