On 04/11/2022 02:43, Matthew Knepley wrote:
On Thu, Nov 3, 2022 at 8:36 PM Matthew Knepley <knep...@gmail.com> wrote:

    On Thu, Oct 27, 2022 at 11:57 AM Semplice Matteo
    <matteo.sempl...@uninsubria.it> wrote:

        Dear Petsc developers,
        I am trying to use a DMSwarm to locate a cloud of points with
        respect to a background mesh. In the real application the
        points will be loaded from disk, but I have created a small
        demo in which

          * each processor creates Npart particles, all within the
            domain covered by the mesh, but not all in the local
            portion of the mesh
          * migrate the particles

        After migration most particles are not any more in the DMSwarm
        (how many and which ones seems to depend on the number of
        cpus, but it never happens that all particle survive the
        migration process).

        I am clearly missing some step, since I'd expect that a DMDA
        would be able to locate particles without the need to go
        through a DMShell as it is done in
        src/dm/tutorials/swarm_ex3.c.html
        
<https://eur01.safelinks.protection.outlook.com/?url=https%3A%2F%2Fpetsc.org%2Fmain%2Fsrc%2Fdm%2Ftutorials%2Fswarm_ex3.c.html&data=05%7C01%7Cmatteo.semplice%40uninsubria.it%7C98335c5548fd4df3ec9208dabe060043%7C9252ed8bdffc401c86ca6237da9991fa%7C0%7C0%7C638031230248558895%7CUnknown%7CTWFpbGZsb3d8eyJWIjoiMC4wLjAwMDAiLCJQIjoiV2luMzIiLCJBTiI6Ik1haWwiLCJXVCI6Mn0%3D%7C3000%7C%7C%7C&sdata=j23JDuWGg3BwgLdSNbO0biatt%2FfyslLISMn2mM6MxrI%3D&reserved=0>

        I attach my demo code.

        Could someone give me a hint?


    Thanks for sending this. I found the problem. Someone has some
    overly fancy code inside DMDA to figure out the local bounding box
    from the coordinates.
    It is broken for DM_BOUNDARY_GHOSTED, but we never tested with
    this. I will fix it.


Okay, I think this fix is correct

https://gitlab.com/petsc/petsc/-/merge_requests/5802 <https://eur01.safelinks.protection.outlook.com/?url=https%3A%2F%2Fgitlab.com%2Fpetsc%2Fpetsc%2F-%2Fmerge_requests%2F5802&data=05%7C01%7Cmatteo.semplice%40uninsubria.it%7C98335c5548fd4df3ec9208dabe060043%7C9252ed8bdffc401c86ca6237da9991fa%7C0%7C0%7C638031230248558895%7CUnknown%7CTWFpbGZsb3d8eyJWIjoiMC4wLjAwMDAiLCJQIjoiV2luMzIiLCJBTiI6Ik1haWwiLCJXVCI6Mn0%3D%7C3000%7C%7C%7C&sdata=qNjWHwveQoOIgvwLq1WBmMRlBr%2FLH%2FQQphOzoxNP7NQ%3D&reserved=0>

I incorporated your test as src/dm/impls/da/tests/ex1.c. Can you take a look and see if this fixes your issue?

Yes, we have tested 2d and 3d, with various combinations of DM_BOUNDARY_* along different directions and it works like a charm.

On a side note, neither DMSwarmViewXDMF nor DMSwarmMigrate seem to be implemented for 1d: I get

[0]PETSC ERROR: No support for this operation for this object type[0]PETSC ERROR: Support not provided for 1D

However, currently I have no need for this feature.

Finally, if the test is meant to stay in the source, you may remove the call to DMSwarmRegisterPetscDatatypeField as in the attached patch.

Thanks a lot!!

    Matteo and Silvia


  Thanks,

     Matt

      Thanks,

         Matt

        Best
            Matteo



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

    https://www.cse.buffalo.edu/~knepley/
    
<https://eur01.safelinks.protection.outlook.com/?url=http%3A%2F%2Fwww.cse.buffalo.edu%2F~knepley%2F&data=05%7C01%7Cmatteo.semplice%40uninsubria.it%7C98335c5548fd4df3ec9208dabe060043%7C9252ed8bdffc401c86ca6237da9991fa%7C0%7C0%7C638031230248558895%7CUnknown%7CTWFpbGZsb3d8eyJWIjoiMC4wLjAwMDAiLCJQIjoiV2luMzIiLCJBTiI6Ik1haWwiLCJXVCI6Mn0%3D%7C3000%7C%7C%7C&sdata=NUrLdRYvJVG9FK%2B66ku2yE6gNsX5xDMsccNsdhSQXHA%3D&reserved=0>



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

https://www.cse.buffalo.edu/~knepley/ <https://eur01.safelinks.protection.outlook.com/?url=http%3A%2F%2Fwww.cse.buffalo.edu%2F~knepley%2F&data=05%7C01%7Cmatteo.semplice%40uninsubria.it%7C98335c5548fd4df3ec9208dabe060043%7C9252ed8bdffc401c86ca6237da9991fa%7C0%7C0%7C638031230248558895%7CUnknown%7CTWFpbGZsb3d8eyJWIjoiMC4wLjAwMDAiLCJQIjoiV2luMzIiLCJBTiI6Ik1haWwiLCJXVCI6Mn0%3D%7C3000%7C%7C%7C&sdata=NUrLdRYvJVG9FK%2B66ku2yE6gNsX5xDMsccNsdhSQXHA%3D&reserved=0>

--
Prof. Matteo Semplice
Università degli Studi dell’Insubria
Dipartimento di Scienza e Alta Tecnologia – DiSAT
Professore Associato
Via Valleggio, 11 – 22100 Como (CO) – Italia
tel.: +39 031 2386316
diff --git a/src/dm/impls/da/tests/ex1.c b/src/dm/impls/da/tests/ex1.c
index d7d77ead3d5..3926aeaa52f 100644
--- a/src/dm/impls/da/tests/ex1.c
+++ b/src/dm/impls/da/tests/ex1.c
@@ -53,7 +53,6 @@ int main(int argc, char **argv)
   PetscCall(DMSetFromOptions(sw));
   PetscCall(DMSwarmSetCellDM(sw, dm));
   PetscCall(DMSwarmInitializeFieldRegister(sw));
-  PetscCall(DMSwarmRegisterPetscDatatypeField(sw, "u", 1, PETSC_SCALAR));
   PetscCall(DMSwarmFinalizeFieldRegister(sw));
   PetscCall(DMSwarmSetLocalSizes(sw, Np, 2));
   PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));

Reply via email to