xref: /petsc/src/dm/impls/swarm/swarmpic.c (revision 966bd95a39c2334d2e2ce17ad22128f3c1861eeb)
10e2ec84fSDave May #define PETSCDM_DLL
20e2ec84fSDave May #include <petsc/private/dmswarmimpl.h> /*I   "petscdmswarm.h"   I*/
30e2ec84fSDave May #include <petscsf.h>
4b799feefSDave May #include <petscdmda.h>
5b799feefSDave May #include <petscdmplex.h>
635b38c60SMatthew G. Knepley #include <petscdt.h>
7279f676cSBarry Smith #include "../src/dm/impls/swarm/data_bucket.h"
80e2ec84fSDave May 
935b38c60SMatthew G. Knepley #include <petsc/private/petscfeimpl.h> /* For CoordinatesRefToReal() */
1035b38c60SMatthew G. Knepley 
1119307e5cSMatthew G. Knepley PetscClassId DMSWARMCELLDM_CLASSID;
1219307e5cSMatthew G. Knepley 
1319307e5cSMatthew G. Knepley /*@
1419307e5cSMatthew G. Knepley   DMSwarmCellDMDestroy - destroy a `DMSwarmCellDM`
1519307e5cSMatthew G. Knepley 
1619307e5cSMatthew G. Knepley   Collective
1719307e5cSMatthew G. Knepley 
1819307e5cSMatthew G. Knepley   Input Parameter:
1919307e5cSMatthew G. Knepley . celldm - address of `DMSwarmCellDM`
2019307e5cSMatthew G. Knepley 
2119307e5cSMatthew G. Knepley   Level: advanced
2219307e5cSMatthew G. Knepley 
2319307e5cSMatthew G. Knepley .seealso: `DMSwarmCellDM`, `DMSwarmCellDMCreate()`
2419307e5cSMatthew G. Knepley @*/
2519307e5cSMatthew G. Knepley PetscErrorCode DMSwarmCellDMDestroy(DMSwarmCellDM *celldm)
2619307e5cSMatthew G. Knepley {
2719307e5cSMatthew G. Knepley   PetscFunctionBegin;
2819307e5cSMatthew G. Knepley   if (!*celldm) PetscFunctionReturn(PETSC_SUCCESS);
2919307e5cSMatthew G. Knepley   PetscValidHeaderSpecific(*celldm, DMSWARMCELLDM_CLASSID, 1);
3019307e5cSMatthew G. Knepley   if (--((PetscObject)*celldm)->refct > 0) {
3119307e5cSMatthew G. Knepley     *celldm = NULL;
3219307e5cSMatthew G. Knepley     PetscFunctionReturn(PETSC_SUCCESS);
3319307e5cSMatthew G. Knepley   }
3419307e5cSMatthew G. Knepley   PetscTryTypeMethod(*celldm, destroy);
3519307e5cSMatthew G. Knepley   for (PetscInt f = 0; f < (*celldm)->Nf; ++f) PetscCall(PetscFree((*celldm)->dmFields[f]));
3619307e5cSMatthew G. Knepley   PetscCall(PetscFree((*celldm)->dmFields));
3719307e5cSMatthew G. Knepley   for (PetscInt f = 0; f < (*celldm)->Nfc; ++f) PetscCall(PetscFree((*celldm)->coordFields[f]));
3819307e5cSMatthew G. Knepley   PetscCall(PetscFree((*celldm)->coordFields));
3919307e5cSMatthew G. Knepley   PetscCall(PetscFree((*celldm)->cellid));
4019307e5cSMatthew G. Knepley   PetscCall(DMSwarmSortDestroy(&(*celldm)->sort));
4119307e5cSMatthew G. Knepley   PetscCall(DMDestroy(&(*celldm)->dm));
4219307e5cSMatthew G. Knepley   PetscCall(PetscHeaderDestroy(celldm));
4319307e5cSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
4419307e5cSMatthew G. Knepley }
4519307e5cSMatthew G. Knepley 
4619307e5cSMatthew G. Knepley /*@
4719307e5cSMatthew G. Knepley   DMSwarmCellDMView - view a `DMSwarmCellDM`
4819307e5cSMatthew G. Knepley 
4919307e5cSMatthew G. Knepley   Collective
5019307e5cSMatthew G. Knepley 
5119307e5cSMatthew G. Knepley   Input Parameters:
5219307e5cSMatthew G. Knepley + celldm - `DMSwarmCellDM`
5319307e5cSMatthew G. Knepley - viewer - viewer to display field, for example `PETSC_VIEWER_STDOUT_WORLD`
5419307e5cSMatthew G. Knepley 
5519307e5cSMatthew G. Knepley   Level: advanced
5619307e5cSMatthew G. Knepley 
5719307e5cSMatthew G. Knepley .seealso: `DMSwarmCellDM`, `DMSwarmCellDMCreate()`
5819307e5cSMatthew G. Knepley @*/
5919307e5cSMatthew G. Knepley PetscErrorCode DMSwarmCellDMView(DMSwarmCellDM celldm, PetscViewer viewer)
6019307e5cSMatthew G. Knepley {
619f196a02SMartin Diehl   PetscBool isascii;
6219307e5cSMatthew G. Knepley 
6319307e5cSMatthew G. Knepley   PetscFunctionBegin;
6419307e5cSMatthew G. Knepley   PetscValidHeaderSpecific(celldm, DMSWARMCELLDM_CLASSID, 1);
6519307e5cSMatthew G. Knepley   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)celldm), &viewer));
6619307e5cSMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
6719307e5cSMatthew G. Knepley   PetscCheckSameComm(celldm, 1, viewer, 2);
689f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
699f196a02SMartin Diehl   if (isascii) {
7019307e5cSMatthew G. Knepley     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)celldm, viewer));
7119307e5cSMatthew G. Knepley     PetscCall(PetscViewerASCIIPushTab(viewer));
7219307e5cSMatthew G. Knepley     PetscCall(PetscViewerASCIIPrintf(viewer, "solution field%s:", celldm->Nf > 1 ? "s" : ""));
7319307e5cSMatthew G. Knepley     for (PetscInt f = 0; f < celldm->Nf; ++f) PetscCall(PetscViewerASCIIPrintf(viewer, " %s", celldm->dmFields[f]));
7419307e5cSMatthew G. Knepley     PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
7519307e5cSMatthew G. Knepley     PetscCall(PetscViewerASCIIPrintf(viewer, "coordinate field%s:", celldm->Nfc > 1 ? "s" : ""));
7619307e5cSMatthew G. Knepley     for (PetscInt f = 0; f < celldm->Nfc; ++f) PetscCall(PetscViewerASCIIPrintf(viewer, " %s", celldm->coordFields[f]));
7719307e5cSMatthew G. Knepley     PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
7819307e5cSMatthew G. Knepley     PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_DEFAULT));
7919307e5cSMatthew G. Knepley     PetscCall(DMView(celldm->dm, viewer));
8019307e5cSMatthew G. Knepley     PetscCall(PetscViewerPopFormat(viewer));
8119307e5cSMatthew G. Knepley   }
8219307e5cSMatthew G. Knepley   PetscTryTypeMethod(celldm, view, viewer);
839f196a02SMartin Diehl   if (isascii) PetscCall(PetscViewerASCIIPopTab(viewer));
8419307e5cSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
8519307e5cSMatthew G. Knepley }
8619307e5cSMatthew G. Knepley 
8719307e5cSMatthew G. Knepley /*@
8819307e5cSMatthew G. Knepley   DMSwarmCellDMGetDM - Returns the background `DM` for the `DMSwarm`
8919307e5cSMatthew G. Knepley 
9019307e5cSMatthew G. Knepley   Not Collective
9119307e5cSMatthew G. Knepley 
9219307e5cSMatthew G. Knepley   Input Parameter:
9319307e5cSMatthew G. Knepley . celldm - The `DMSwarmCellDM` object
9419307e5cSMatthew G. Knepley 
9519307e5cSMatthew G. Knepley   Output Parameter:
9619307e5cSMatthew G. Knepley . dm - The `DM` object
9719307e5cSMatthew G. Knepley 
9819307e5cSMatthew G. Knepley   Level: intermediate
9919307e5cSMatthew G. Knepley 
10019307e5cSMatthew G. Knepley .seealso: `DMSwarmCellDM`, `DM`, `DMSwarmSetCellDM()`
10119307e5cSMatthew G. Knepley @*/
10219307e5cSMatthew G. Knepley PetscErrorCode DMSwarmCellDMGetDM(DMSwarmCellDM celldm, DM *dm)
10319307e5cSMatthew G. Knepley {
10419307e5cSMatthew G. Knepley   PetscFunctionBegin;
10519307e5cSMatthew G. Knepley   PetscValidHeaderSpecific(celldm, DMSWARMCELLDM_CLASSID, 1);
10619307e5cSMatthew G. Knepley   PetscAssertPointer(dm, 2);
10719307e5cSMatthew G. Knepley   *dm = celldm->dm;
10819307e5cSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
10919307e5cSMatthew G. Knepley }
11019307e5cSMatthew G. Knepley 
11119307e5cSMatthew G. Knepley /*@C
11219307e5cSMatthew G. Knepley   DMSwarmCellDMGetFields - Returns the `DM` fields for the `DMSwarm`
11319307e5cSMatthew G. Knepley 
11419307e5cSMatthew G. Knepley   Not Collective
11519307e5cSMatthew G. Knepley 
11619307e5cSMatthew G. Knepley   Input Parameter:
11719307e5cSMatthew G. Knepley . celldm - The `DMSwarmCellDM` object
11819307e5cSMatthew G. Knepley 
11919307e5cSMatthew G. Knepley   Output Parameters:
12019307e5cSMatthew G. Knepley + Nf    - The number of fields
12119307e5cSMatthew G. Knepley - names - The array of field names in the `DMSWARM`
12219307e5cSMatthew G. Knepley 
12319307e5cSMatthew G. Knepley   Level: intermediate
12419307e5cSMatthew G. Knepley 
12519307e5cSMatthew G. Knepley .seealso: `DMSwarmCellDM`, `DM`, `DMSwarmSetCellDM()`
12619307e5cSMatthew G. Knepley @*/
12719307e5cSMatthew G. Knepley PetscErrorCode DMSwarmCellDMGetFields(DMSwarmCellDM celldm, PetscInt *Nf, const char **names[])
12819307e5cSMatthew G. Knepley {
12919307e5cSMatthew G. Knepley   PetscFunctionBegin;
13019307e5cSMatthew G. Knepley   PetscValidHeaderSpecific(celldm, DMSWARMCELLDM_CLASSID, 1);
13119307e5cSMatthew G. Knepley   if (Nf) {
13219307e5cSMatthew G. Knepley     PetscAssertPointer(Nf, 2);
13319307e5cSMatthew G. Knepley     *Nf = celldm->Nf;
13419307e5cSMatthew G. Knepley   }
13519307e5cSMatthew G. Knepley   if (names) {
13619307e5cSMatthew G. Knepley     PetscAssertPointer(names, 3);
13719307e5cSMatthew G. Knepley     *names = (const char **)celldm->dmFields;
13819307e5cSMatthew G. Knepley   }
13919307e5cSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
14019307e5cSMatthew G. Knepley }
14119307e5cSMatthew G. Knepley 
14219307e5cSMatthew G. Knepley /*@C
14319307e5cSMatthew G. Knepley   DMSwarmCellDMGetCoordinateFields - Returns the `DM` coordinate fields for the `DMSwarm`
14419307e5cSMatthew G. Knepley 
14519307e5cSMatthew G. Knepley   Not Collective
14619307e5cSMatthew G. Knepley 
14719307e5cSMatthew G. Knepley   Input Parameter:
14819307e5cSMatthew G. Knepley . celldm - The `DMSwarmCellDM` object
14919307e5cSMatthew G. Knepley 
15019307e5cSMatthew G. Knepley   Output Parameters:
15119307e5cSMatthew G. Knepley + Nfc   - The number of coordinate fields
15219307e5cSMatthew G. Knepley - names - The array of coordinate field names in the `DMSWARM`
15319307e5cSMatthew G. Knepley 
15419307e5cSMatthew G. Knepley   Level: intermediate
15519307e5cSMatthew G. Knepley 
15619307e5cSMatthew G. Knepley .seealso: `DMSwarmCellDM`, `DM`, `DMSwarmSetCellDM()`
15719307e5cSMatthew G. Knepley @*/
15819307e5cSMatthew G. Knepley PetscErrorCode DMSwarmCellDMGetCoordinateFields(DMSwarmCellDM celldm, PetscInt *Nfc, const char **names[])
15919307e5cSMatthew G. Knepley {
16019307e5cSMatthew G. Knepley   PetscFunctionBegin;
16119307e5cSMatthew G. Knepley   PetscValidHeaderSpecific(celldm, DMSWARMCELLDM_CLASSID, 1);
16219307e5cSMatthew G. Knepley   if (Nfc) {
16319307e5cSMatthew G. Knepley     PetscAssertPointer(Nfc, 2);
16419307e5cSMatthew G. Knepley     *Nfc = celldm->Nfc;
16519307e5cSMatthew G. Knepley   }
16619307e5cSMatthew G. Knepley   if (names) {
16719307e5cSMatthew G. Knepley     PetscAssertPointer(names, 3);
16819307e5cSMatthew G. Knepley     *names = (const char **)celldm->coordFields;
16919307e5cSMatthew G. Knepley   }
17019307e5cSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
17119307e5cSMatthew G. Knepley }
17219307e5cSMatthew G. Knepley 
17319307e5cSMatthew G. Knepley /*@C
17419307e5cSMatthew G. Knepley   DMSwarmCellDMGetCellID - Returns the cell id field name for the `DMSwarm`
17519307e5cSMatthew G. Knepley 
17619307e5cSMatthew G. Knepley   Not Collective
17719307e5cSMatthew G. Knepley 
17819307e5cSMatthew G. Knepley   Input Parameter:
17919307e5cSMatthew G. Knepley . celldm - The `DMSwarmCellDM` object
18019307e5cSMatthew G. Knepley 
18119307e5cSMatthew G. Knepley   Output Parameters:
18219307e5cSMatthew G. Knepley . cellid - The cell id field name in the `DMSWARM`
18319307e5cSMatthew G. Knepley 
18419307e5cSMatthew G. Knepley   Level: intermediate
18519307e5cSMatthew G. Knepley 
18619307e5cSMatthew G. Knepley .seealso: `DMSwarmCellDM`, `DM`, `DMSwarmSetCellDM()`
18719307e5cSMatthew G. Knepley @*/
18819307e5cSMatthew G. Knepley PetscErrorCode DMSwarmCellDMGetCellID(DMSwarmCellDM celldm, const char *cellid[])
18919307e5cSMatthew G. Knepley {
19019307e5cSMatthew G. Knepley   PetscFunctionBegin;
19119307e5cSMatthew G. Knepley   PetscValidHeaderSpecific(celldm, DMSWARMCELLDM_CLASSID, 1);
19219307e5cSMatthew G. Knepley   PetscAssertPointer(cellid, 2);
19319307e5cSMatthew G. Knepley   *cellid = celldm->cellid;
19419307e5cSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
19519307e5cSMatthew G. Knepley }
19619307e5cSMatthew G. Knepley 
19719307e5cSMatthew G. Knepley /*@C
19819307e5cSMatthew G. Knepley   DMSwarmCellDMGetSort - Returns the sort context over the active `DMSwarmCellDM` for the `DMSwarm`
19919307e5cSMatthew G. Knepley 
20019307e5cSMatthew G. Knepley   Not Collective
20119307e5cSMatthew G. Knepley 
20219307e5cSMatthew G. Knepley   Input Parameter:
20319307e5cSMatthew G. Knepley . celldm - The `DMSwarmCellDM` object
20419307e5cSMatthew G. Knepley 
20519307e5cSMatthew G. Knepley   Output Parameter:
20619307e5cSMatthew G. Knepley . sort - The `DMSwarmSort` object
20719307e5cSMatthew G. Knepley 
20819307e5cSMatthew G. Knepley   Level: intermediate
20919307e5cSMatthew G. Knepley 
21019307e5cSMatthew G. Knepley .seealso: `DMSwarmCellDM`, `DM`, `DMSwarmSetCellDM()`
21119307e5cSMatthew G. Knepley @*/
21219307e5cSMatthew G. Knepley PetscErrorCode DMSwarmCellDMGetSort(DMSwarmCellDM celldm, DMSwarmSort *sort)
21319307e5cSMatthew G. Knepley {
21419307e5cSMatthew G. Knepley   PetscFunctionBegin;
21519307e5cSMatthew G. Knepley   PetscValidHeaderSpecific(celldm, DMSWARMCELLDM_CLASSID, 1);
21619307e5cSMatthew G. Knepley   PetscAssertPointer(sort, 2);
21719307e5cSMatthew G. Knepley   *sort = celldm->sort;
21819307e5cSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
21919307e5cSMatthew G. Knepley }
22019307e5cSMatthew G. Knepley 
22119307e5cSMatthew G. Knepley /*@C
22219307e5cSMatthew G. Knepley   DMSwarmCellDMSetSort - Sets the sort context over the active `DMSwarmCellDM` for the `DMSwarm`
22319307e5cSMatthew G. Knepley 
22419307e5cSMatthew G. Knepley   Not Collective
22519307e5cSMatthew G. Knepley 
22619307e5cSMatthew G. Knepley   Input Parameters:
22719307e5cSMatthew G. Knepley + celldm - The `DMSwarmCellDM` object
22819307e5cSMatthew G. Knepley - sort   - The `DMSwarmSort` object
22919307e5cSMatthew G. Knepley 
23019307e5cSMatthew G. Knepley   Level: intermediate
23119307e5cSMatthew G. Knepley 
23219307e5cSMatthew G. Knepley .seealso: `DMSwarmCellDM`, `DM`, `DMSwarmSetCellDM()`
23319307e5cSMatthew G. Knepley @*/
23419307e5cSMatthew G. Knepley PetscErrorCode DMSwarmCellDMSetSort(DMSwarmCellDM celldm, DMSwarmSort sort)
23519307e5cSMatthew G. Knepley {
23619307e5cSMatthew G. Knepley   PetscFunctionBegin;
23719307e5cSMatthew G. Knepley   PetscValidHeaderSpecific(celldm, DMSWARMCELLDM_CLASSID, 1);
23819307e5cSMatthew G. Knepley   PetscAssertPointer(sort, 2);
23919307e5cSMatthew G. Knepley   celldm->sort = sort;
24019307e5cSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
24119307e5cSMatthew G. Knepley }
24219307e5cSMatthew G. Knepley 
24319307e5cSMatthew G. Knepley /*@
24419307e5cSMatthew G. Knepley   DMSwarmCellDMGetBlockSize - Returns the total blocksize for the `DM` fields
24519307e5cSMatthew G. Knepley 
24619307e5cSMatthew G. Knepley   Not Collective
24719307e5cSMatthew G. Knepley 
24819307e5cSMatthew G. Knepley   Input Parameters:
24919307e5cSMatthew G. Knepley + celldm - The `DMSwarmCellDM` object
25019307e5cSMatthew G. Knepley - sw     - The `DMSwarm` object
25119307e5cSMatthew G. Knepley 
25219307e5cSMatthew G. Knepley   Output Parameter:
25319307e5cSMatthew G. Knepley . bs - The total block size
25419307e5cSMatthew G. Knepley 
25519307e5cSMatthew G. Knepley   Level: intermediate
25619307e5cSMatthew G. Knepley 
25719307e5cSMatthew G. Knepley .seealso: `DMSwarmCellDM`, `DM`, `DMSwarmSetCellDM()`
25819307e5cSMatthew G. Knepley @*/
25919307e5cSMatthew G. Knepley PetscErrorCode DMSwarmCellDMGetBlockSize(DMSwarmCellDM celldm, DM sw, PetscInt *bs)
26019307e5cSMatthew G. Knepley {
26119307e5cSMatthew G. Knepley   PetscFunctionBegin;
26219307e5cSMatthew G. Knepley   PetscValidHeaderSpecific(celldm, DMSWARMCELLDM_CLASSID, 1);
26319307e5cSMatthew G. Knepley   PetscValidHeaderSpecific(sw, DM_CLASSID, 2);
26419307e5cSMatthew G. Knepley   PetscAssertPointer(bs, 3);
26519307e5cSMatthew G. Knepley   *bs = 0;
26619307e5cSMatthew G. Knepley   for (PetscInt f = 0; f < celldm->Nf; ++f) {
26719307e5cSMatthew G. Knepley     PetscInt fbs;
26819307e5cSMatthew G. Knepley 
26919307e5cSMatthew G. Knepley     PetscCall(DMSwarmGetFieldInfo(sw, celldm->dmFields[f], &fbs, NULL));
27019307e5cSMatthew G. Knepley     *bs += fbs;
27119307e5cSMatthew G. Knepley   }
27219307e5cSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
27319307e5cSMatthew G. Knepley }
27419307e5cSMatthew G. Knepley 
27519307e5cSMatthew G. Knepley /*@
27619307e5cSMatthew G. Knepley   DMSwarmCellDMCreate - create a `DMSwarmCellDM`
27719307e5cSMatthew G. Knepley 
27819307e5cSMatthew G. Knepley   Collective
27919307e5cSMatthew G. Knepley 
280d6ae5217SJose E. Roman   Input Parameters:
28119307e5cSMatthew G. Knepley + dm          - The background `DM` for the `DMSwarm`
28219307e5cSMatthew G. Knepley . Nf          - The number of swarm fields defined over `dm`
28319307e5cSMatthew G. Knepley . dmFields    - The swarm field names for the `dm` fields
28419307e5cSMatthew G. Knepley . Nfc         - The number of swarm fields to use for coordinates over `dm`
28519307e5cSMatthew G. Knepley - coordFields - The swarm field names for the `dm` coordinate fields
28619307e5cSMatthew G. Knepley 
28719307e5cSMatthew G. Knepley   Output Parameter:
28819307e5cSMatthew G. Knepley . celldm - The new `DMSwarmCellDM`
28919307e5cSMatthew G. Knepley 
29019307e5cSMatthew G. Knepley   Level: advanced
29119307e5cSMatthew G. Knepley 
29219307e5cSMatthew G. Knepley .seealso: `DMSwarmCellDM`, `DMSWARM`, `DMSetType()`
29319307e5cSMatthew G. Knepley @*/
29419307e5cSMatthew G. Knepley PetscErrorCode DMSwarmCellDMCreate(DM dm, PetscInt Nf, const char *dmFields[], PetscInt Nfc, const char *coordFields[], DMSwarmCellDM *celldm)
29519307e5cSMatthew G. Knepley {
29619307e5cSMatthew G. Knepley   DMSwarmCellDM b;
29719307e5cSMatthew G. Knepley   const char   *name;
29819307e5cSMatthew G. Knepley   char          cellid[PETSC_MAX_PATH_LEN];
29919307e5cSMatthew G. Knepley 
30019307e5cSMatthew G. Knepley   PetscFunctionBegin;
30119307e5cSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
30219307e5cSMatthew G. Knepley   if (Nf) PetscAssertPointer(dmFields, 3);
30319307e5cSMatthew G. Knepley   if (Nfc) PetscAssertPointer(coordFields, 5);
30419307e5cSMatthew G. Knepley   PetscCall(DMInitializePackage());
30519307e5cSMatthew G. Knepley 
30619307e5cSMatthew G. Knepley   PetscCall(PetscHeaderCreate(b, DMSWARMCELLDM_CLASSID, "DMSwarmCellDM", "Background DM for a Swarm", "DM", PetscObjectComm((PetscObject)dm), DMSwarmCellDMDestroy, DMSwarmCellDMView));
30719307e5cSMatthew G. Knepley   PetscCall(PetscObjectGetName((PetscObject)dm, &name));
30819307e5cSMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)b, name));
30919307e5cSMatthew G. Knepley   PetscCall(PetscObjectReference((PetscObject)dm));
31019307e5cSMatthew G. Knepley   b->dm  = dm;
31119307e5cSMatthew G. Knepley   b->Nf  = Nf;
31219307e5cSMatthew G. Knepley   b->Nfc = Nfc;
31319307e5cSMatthew G. Knepley   PetscCall(PetscMalloc1(b->Nf, &b->dmFields));
31419307e5cSMatthew G. Knepley   for (PetscInt f = 0; f < b->Nf; ++f) PetscCall(PetscStrallocpy(dmFields[f], &b->dmFields[f]));
31519307e5cSMatthew G. Knepley   PetscCall(PetscMalloc1(b->Nfc, &b->coordFields));
31619307e5cSMatthew G. Knepley   for (PetscInt f = 0; f < b->Nfc; ++f) PetscCall(PetscStrallocpy(coordFields[f], &b->coordFields[f]));
31719307e5cSMatthew G. Knepley   PetscCall(PetscSNPrintf(cellid, PETSC_MAX_PATH_LEN, "%s_cellid", name));
31819307e5cSMatthew G. Knepley   PetscCall(PetscStrallocpy(cellid, &b->cellid));
31919307e5cSMatthew G. Knepley   *celldm = b;
32019307e5cSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
32119307e5cSMatthew G. Knepley }
32219307e5cSMatthew G. Knepley 
3230e2ec84fSDave May /* Coordinate insertition/addition API */
324cc4c1da9SBarry Smith /*@
32520f4b53cSBarry Smith   DMSwarmSetPointsUniformCoordinates - Set point coordinates in a `DMSWARM` on a regular (ijk) grid
3260e2ec84fSDave May 
32720f4b53cSBarry Smith   Collective
3280e2ec84fSDave May 
32960225df5SJacob Faibussowitsch   Input Parameters:
33019307e5cSMatthew G. Knepley + sw      - the `DMSWARM`
3310e2ec84fSDave May . min     - minimum coordinate values in the x, y, z directions (array of length dim)
3320e2ec84fSDave May . max     - maximum coordinate values in the x, y, z directions (array of length dim)
3330e2ec84fSDave May . npoints - number of points in each spatial direction (array of length dim)
33420f4b53cSBarry Smith - mode    - indicates whether to append points to the swarm (`ADD_VALUES`), or over-ride existing points (`INSERT_VALUES`)
3350e2ec84fSDave May 
3360e2ec84fSDave May   Level: beginner
3370e2ec84fSDave May 
3380e2ec84fSDave May   Notes:
33920f4b53cSBarry Smith   When using mode = `INSERT_VALUES`, this method will reset the number of particles in the `DMSWARM`
340a3b724e8SBarry Smith   to be `npoints[0]` x `npoints[1]` (2D) or `npoints[0]` x `npoints[1]` x `npoints[2]` (3D). When using mode = `ADD_VALUES`,
34120f4b53cSBarry Smith   new points will be appended to any already existing in the `DMSWARM`
3420e2ec84fSDave May 
34320f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmSetType()`, `DMSwarmSetCellDM()`, `DMSwarmType`
3440e2ec84fSDave May @*/
34519307e5cSMatthew G. Knepley PetscErrorCode DMSwarmSetPointsUniformCoordinates(DM sw, PetscReal min[], PetscReal max[], PetscInt npoints[], InsertMode mode)
346d71ae5a4SJacob Faibussowitsch {
347c448b993SMatthew G. Knepley   PetscReal          lmin[] = {PETSC_MAX_REAL, PETSC_MAX_REAL, PETSC_MAX_REAL};
348c448b993SMatthew G. Knepley   PetscReal          lmax[] = {PETSC_MIN_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL};
34919307e5cSMatthew G. Knepley   PetscInt           i, j, k, bs, b, n_estimate, n_curr, n_new_est, p, n_found, Nfc;
3500e2ec84fSDave May   const PetscScalar *_coor;
35119307e5cSMatthew G. Knepley   DMSwarmCellDM      celldm;
35219307e5cSMatthew G. Knepley   DM                 dm;
3530e2ec84fSDave May   PetscReal          dx[3];
3542e3d3761SDave May   PetscInt           _npoints[] = {0, 0, 1};
3550e2ec84fSDave May   Vec                pos;
3560e2ec84fSDave May   PetscScalar       *_pos;
3570e2ec84fSDave May   PetscReal         *swarm_coor;
3580e2ec84fSDave May   PetscInt          *swarm_cellid;
3590e2ec84fSDave May   PetscSF            sfcell = NULL;
3600e2ec84fSDave May   const PetscSFNode *LA_sfcell;
36119307e5cSMatthew G. Knepley   const char       **coordFields, *cellid;
3620e2ec84fSDave May 
3630e2ec84fSDave May   PetscFunctionBegin;
36419307e5cSMatthew G. Knepley   DMSWARMPICVALID(sw);
36519307e5cSMatthew G. Knepley   PetscCall(DMSwarmGetCellDMActive(sw, &celldm));
36619307e5cSMatthew G. Knepley   PetscCall(DMSwarmCellDMGetCoordinateFields(celldm, &Nfc, &coordFields));
36719307e5cSMatthew G. Knepley   PetscCheck(Nfc == 1, PetscObjectComm((PetscObject)sw), PETSC_ERR_SUP, "We only support a single coordinate field right now, not %" PetscInt_FMT, Nfc);
36819307e5cSMatthew G. Knepley 
36919307e5cSMatthew G. Knepley   PetscCall(DMSwarmCellDMGetDM(celldm, &dm));
37019307e5cSMatthew G. Knepley   PetscCall(DMGetLocalBoundingBox(dm, lmin, lmax));
37119307e5cSMatthew G. Knepley   PetscCall(DMGetCoordinateDim(dm, &bs));
3720e2ec84fSDave May 
3730e2ec84fSDave May   for (b = 0; b < bs; b++) {
37471844bb8SDave May     if (npoints[b] > 1) {
3750e2ec84fSDave May       dx[b] = (max[b] - min[b]) / ((PetscReal)(npoints[b] - 1));
37671844bb8SDave May     } else {
37771844bb8SDave May       dx[b] = 0.0;
37871844bb8SDave May     }
3792e3d3761SDave May     _npoints[b] = npoints[b];
3800e2ec84fSDave May   }
3810e2ec84fSDave May 
3820e2ec84fSDave May   /* determine number of points living in the bounding box */
3830e2ec84fSDave May   n_estimate = 0;
3842e3d3761SDave May   for (k = 0; k < _npoints[2]; k++) {
3852e3d3761SDave May     for (j = 0; j < _npoints[1]; j++) {
3862e3d3761SDave May       for (i = 0; i < _npoints[0]; i++) {
3870e2ec84fSDave May         PetscReal xp[] = {0.0, 0.0, 0.0};
3880e2ec84fSDave May         PetscInt  ijk[3];
3890e2ec84fSDave May         PetscBool point_inside = PETSC_TRUE;
3900e2ec84fSDave May 
3910e2ec84fSDave May         ijk[0] = i;
3920e2ec84fSDave May         ijk[1] = j;
3930e2ec84fSDave May         ijk[2] = k;
394ad540459SPierre Jolivet         for (b = 0; b < bs; b++) xp[b] = min[b] + ijk[b] * dx[b];
3950e2ec84fSDave May         for (b = 0; b < bs; b++) {
396c448b993SMatthew G. Knepley           if (xp[b] < lmin[b]) point_inside = PETSC_FALSE;
397c448b993SMatthew G. Knepley           if (xp[b] > lmax[b]) point_inside = PETSC_FALSE;
3980e2ec84fSDave May         }
399ad540459SPierre Jolivet         if (point_inside) n_estimate++;
4000e2ec84fSDave May       }
4010e2ec84fSDave May     }
4020e2ec84fSDave May   }
4030e2ec84fSDave May 
4040e2ec84fSDave May   /* create candidate list */
4059566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_SELF, &pos));
4069566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(pos, bs * n_estimate, PETSC_DECIDE));
4079566063dSJacob Faibussowitsch   PetscCall(VecSetBlockSize(pos, bs));
4089566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(pos));
4099566063dSJacob Faibussowitsch   PetscCall(VecGetArray(pos, &_pos));
4100e2ec84fSDave May 
4110e2ec84fSDave May   n_estimate = 0;
4122e3d3761SDave May   for (k = 0; k < _npoints[2]; k++) {
4132e3d3761SDave May     for (j = 0; j < _npoints[1]; j++) {
4142e3d3761SDave May       for (i = 0; i < _npoints[0]; i++) {
4150e2ec84fSDave May         PetscReal xp[] = {0.0, 0.0, 0.0};
4160e2ec84fSDave May         PetscInt  ijk[3];
4170e2ec84fSDave May         PetscBool point_inside = PETSC_TRUE;
4180e2ec84fSDave May 
4190e2ec84fSDave May         ijk[0] = i;
4200e2ec84fSDave May         ijk[1] = j;
4210e2ec84fSDave May         ijk[2] = k;
422ad540459SPierre Jolivet         for (b = 0; b < bs; b++) xp[b] = min[b] + ijk[b] * dx[b];
4230e2ec84fSDave May         for (b = 0; b < bs; b++) {
424c448b993SMatthew G. Knepley           if (xp[b] < lmin[b]) point_inside = PETSC_FALSE;
425c448b993SMatthew G. Knepley           if (xp[b] > lmax[b]) point_inside = PETSC_FALSE;
4260e2ec84fSDave May         }
4270e2ec84fSDave May         if (point_inside) {
428ad540459SPierre Jolivet           for (b = 0; b < bs; b++) _pos[bs * n_estimate + b] = xp[b];
4290e2ec84fSDave May           n_estimate++;
4300e2ec84fSDave May         }
4310e2ec84fSDave May       }
4320e2ec84fSDave May     }
4330e2ec84fSDave May   }
4349566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(pos, &_pos));
4350e2ec84fSDave May 
4360e2ec84fSDave May   /* locate points */
43719307e5cSMatthew G. Knepley   PetscCall(DMLocatePoints(dm, pos, DM_POINTLOCATION_NONE, &sfcell));
4389566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell));
4390e2ec84fSDave May   n_found = 0;
4400e2ec84fSDave May   for (p = 0; p < n_estimate; p++) {
441ad540459SPierre Jolivet     if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) n_found++;
4420e2ec84fSDave May   }
4430e2ec84fSDave May 
4440e2ec84fSDave May   /* adjust size */
4450e2ec84fSDave May   if (mode == ADD_VALUES) {
44619307e5cSMatthew G. Knepley     PetscCall(DMSwarmGetLocalSize(sw, &n_curr));
4470e2ec84fSDave May     n_new_est = n_curr + n_found;
44819307e5cSMatthew G. Knepley     PetscCall(DMSwarmSetLocalSizes(sw, n_new_est, -1));
4490e2ec84fSDave May   }
4500e2ec84fSDave May   if (mode == INSERT_VALUES) {
4510e2ec84fSDave May     n_curr    = 0;
4520e2ec84fSDave May     n_new_est = n_found;
45319307e5cSMatthew G. Knepley     PetscCall(DMSwarmSetLocalSizes(sw, n_new_est, -1));
4540e2ec84fSDave May   }
4550e2ec84fSDave May 
4560e2ec84fSDave May   /* initialize new coords, cell owners, pid */
45719307e5cSMatthew G. Knepley   PetscCall(DMSwarmCellDMGetCellID(celldm, &cellid));
4589566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(pos, &_coor));
45919307e5cSMatthew G. Knepley   PetscCall(DMSwarmGetField(sw, coordFields[0], NULL, NULL, (void **)&swarm_coor));
46019307e5cSMatthew G. Knepley   PetscCall(DMSwarmGetField(sw, cellid, NULL, NULL, (void **)&swarm_cellid));
4610e2ec84fSDave May   n_found = 0;
4620e2ec84fSDave May   for (p = 0; p < n_estimate; p++) {
4630e2ec84fSDave May     if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) {
464ad540459SPierre Jolivet       for (b = 0; b < bs; b++) swarm_coor[bs * (n_curr + n_found) + b] = PetscRealPart(_coor[bs * p + b]);
4650e2ec84fSDave May       swarm_cellid[n_curr + n_found] = LA_sfcell[p].index;
4660e2ec84fSDave May       n_found++;
4670e2ec84fSDave May     }
4680e2ec84fSDave May   }
46919307e5cSMatthew G. Knepley   PetscCall(DMSwarmRestoreField(sw, cellid, NULL, NULL, (void **)&swarm_cellid));
47019307e5cSMatthew G. Knepley   PetscCall(DMSwarmRestoreField(sw, coordFields[0], NULL, NULL, (void **)&swarm_coor));
4719566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(pos, &_coor));
4720e2ec84fSDave May 
4739566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sfcell));
4749566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pos));
4753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4760e2ec84fSDave May }
4770e2ec84fSDave May 
478cc4c1da9SBarry Smith /*@
47920f4b53cSBarry Smith   DMSwarmSetPointCoordinates - Set point coordinates in a `DMSWARM` from a user defined list
4800e2ec84fSDave May 
48120f4b53cSBarry Smith   Collective
4820e2ec84fSDave May 
48360225df5SJacob Faibussowitsch   Input Parameters:
48419307e5cSMatthew G. Knepley + sw        - the `DMSWARM`
4850e2ec84fSDave May . npoints   - the number of points to insert
4860e2ec84fSDave May . coor      - the coordinate values
48720f4b53cSBarry Smith . redundant - if set to `PETSC_TRUE`, it is assumed that `npoints` and `coor` are only valid on rank 0 and should be broadcast to other ranks
48820f4b53cSBarry Smith - mode      - indicates whether to append points to the swarm (`ADD_VALUES`), or over-ride existing points (`INSERT_VALUES`)
4890e2ec84fSDave May 
4900e2ec84fSDave May   Level: beginner
4910e2ec84fSDave May 
4920e2ec84fSDave May   Notes:
49320f4b53cSBarry Smith   If the user has specified `redundant` as `PETSC_FALSE`, the cell `DM` will attempt to locate the coordinates provided by `coor` within
49420f4b53cSBarry Smith   its sub-domain. If they any values within `coor` are not located in the sub-domain, they will be ignored and will not get
49520f4b53cSBarry Smith   added to the `DMSWARM`.
4960e2ec84fSDave May 
49720f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmSetType()`, `DMSwarmSetCellDM()`, `DMSwarmType`, `DMSwarmSetPointsUniformCoordinates()`
4980e2ec84fSDave May @*/
49919307e5cSMatthew G. Knepley PetscErrorCode DMSwarmSetPointCoordinates(DM sw, PetscInt npoints, PetscReal coor[], PetscBool redundant, InsertMode mode)
500d71ae5a4SJacob Faibussowitsch {
5010e2ec84fSDave May   PetscReal          gmin[] = {PETSC_MAX_REAL, PETSC_MAX_REAL, PETSC_MAX_REAL};
5020e2ec84fSDave May   PetscReal          gmax[] = {PETSC_MIN_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL};
5030e2ec84fSDave May   PetscInt           i, N, bs, b, n_estimate, n_curr, n_new_est, p, n_found;
5040e2ec84fSDave May   Vec                coorlocal;
5050e2ec84fSDave May   const PetscScalar *_coor;
50619307e5cSMatthew G. Knepley   DMSwarmCellDM      celldm;
50719307e5cSMatthew G. Knepley   DM                 dm;
5080e2ec84fSDave May   Vec                pos;
5090e2ec84fSDave May   PetscScalar       *_pos;
5100e2ec84fSDave May   PetscReal         *swarm_coor;
5110e2ec84fSDave May   PetscInt          *swarm_cellid;
5120e2ec84fSDave May   PetscSF            sfcell = NULL;
5130e2ec84fSDave May   const PetscSFNode *LA_sfcell;
5140e2ec84fSDave May   PetscReal         *my_coor;
51519307e5cSMatthew G. Knepley   PetscInt           my_npoints, Nfc;
5160e2ec84fSDave May   PetscMPIInt        rank;
5170e2ec84fSDave May   MPI_Comm           comm;
51819307e5cSMatthew G. Knepley   const char       **coordFields, *cellid;
5190e2ec84fSDave May 
5200e2ec84fSDave May   PetscFunctionBegin;
52119307e5cSMatthew G. Knepley   DMSWARMPICVALID(sw);
52219307e5cSMatthew G. Knepley   PetscCall(PetscObjectGetComm((PetscObject)sw, &comm));
5239566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
5240e2ec84fSDave May 
52519307e5cSMatthew G. Knepley   PetscCall(DMSwarmGetCellDMActive(sw, &celldm));
52619307e5cSMatthew G. Knepley   PetscCall(DMSwarmCellDMGetCoordinateFields(celldm, &Nfc, &coordFields));
52719307e5cSMatthew G. Knepley   PetscCheck(Nfc == 1, PetscObjectComm((PetscObject)sw), PETSC_ERR_SUP, "We only support a single coordinate field right now, not %" PetscInt_FMT, Nfc);
52819307e5cSMatthew G. Knepley 
52919307e5cSMatthew G. Knepley   PetscCall(DMSwarmCellDMGetDM(celldm, &dm));
53019307e5cSMatthew G. Knepley   PetscCall(DMGetCoordinatesLocal(dm, &coorlocal));
5319566063dSJacob Faibussowitsch   PetscCall(VecGetSize(coorlocal, &N));
5329566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(coorlocal, &bs));
5330e2ec84fSDave May   N = N / bs;
5349566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(coorlocal, &_coor));
5350e2ec84fSDave May   for (i = 0; i < N; i++) {
5360e2ec84fSDave May     for (b = 0; b < bs; b++) {
537a5f152d1SDave May       gmin[b] = PetscMin(gmin[b], PetscRealPart(_coor[bs * i + b]));
538a5f152d1SDave May       gmax[b] = PetscMax(gmax[b], PetscRealPart(_coor[bs * i + b]));
5390e2ec84fSDave May     }
5400e2ec84fSDave May   }
5419566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(coorlocal, &_coor));
5420e2ec84fSDave May 
5430e2ec84fSDave May   /* broadcast points from rank 0 if requested */
5440e2ec84fSDave May   if (redundant) {
5456497c311SBarry Smith     PetscMPIInt imy;
5466497c311SBarry Smith 
5470e2ec84fSDave May     my_npoints = npoints;
5489566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(&my_npoints, 1, MPIU_INT, 0, comm));
5490e2ec84fSDave May 
5500e2ec84fSDave May     if (rank > 0) { /* allocate space */
5519566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(bs * my_npoints, &my_coor));
5520e2ec84fSDave May     } else {
5530e2ec84fSDave May       my_coor = coor;
5540e2ec84fSDave May     }
5556497c311SBarry Smith     PetscCall(PetscMPIIntCast(bs * my_npoints, &imy));
5566497c311SBarry Smith     PetscCallMPI(MPI_Bcast(my_coor, imy, MPIU_REAL, 0, comm));
5570e2ec84fSDave May   } else {
5580e2ec84fSDave May     my_npoints = npoints;
5590e2ec84fSDave May     my_coor    = coor;
5600e2ec84fSDave May   }
5610e2ec84fSDave May 
5620e2ec84fSDave May   /* determine the number of points living in the bounding box */
5630e2ec84fSDave May   n_estimate = 0;
5640e2ec84fSDave May   for (i = 0; i < my_npoints; i++) {
5650e2ec84fSDave May     PetscBool point_inside = PETSC_TRUE;
5660e2ec84fSDave May 
5670e2ec84fSDave May     for (b = 0; b < bs; b++) {
568ad540459SPierre Jolivet       if (my_coor[bs * i + b] < gmin[b]) point_inside = PETSC_FALSE;
569ad540459SPierre Jolivet       if (my_coor[bs * i + b] > gmax[b]) point_inside = PETSC_FALSE;
5700e2ec84fSDave May     }
571ad540459SPierre Jolivet     if (point_inside) n_estimate++;
5720e2ec84fSDave May   }
5730e2ec84fSDave May 
5740e2ec84fSDave May   /* create candidate list */
5759566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_SELF, &pos));
5769566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(pos, bs * n_estimate, PETSC_DECIDE));
5779566063dSJacob Faibussowitsch   PetscCall(VecSetBlockSize(pos, bs));
5789566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(pos));
5799566063dSJacob Faibussowitsch   PetscCall(VecGetArray(pos, &_pos));
5800e2ec84fSDave May 
5810e2ec84fSDave May   n_estimate = 0;
5820e2ec84fSDave May   for (i = 0; i < my_npoints; i++) {
5830e2ec84fSDave May     PetscBool point_inside = PETSC_TRUE;
5840e2ec84fSDave May 
5850e2ec84fSDave May     for (b = 0; b < bs; b++) {
586ad540459SPierre Jolivet       if (my_coor[bs * i + b] < gmin[b]) point_inside = PETSC_FALSE;
587ad540459SPierre Jolivet       if (my_coor[bs * i + b] > gmax[b]) point_inside = PETSC_FALSE;
5880e2ec84fSDave May     }
5890e2ec84fSDave May     if (point_inside) {
590ad540459SPierre Jolivet       for (b = 0; b < bs; b++) _pos[bs * n_estimate + b] = my_coor[bs * i + b];
5910e2ec84fSDave May       n_estimate++;
5920e2ec84fSDave May     }
5930e2ec84fSDave May   }
5949566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(pos, &_pos));
5950e2ec84fSDave May 
5960e2ec84fSDave May   /* locate points */
59719307e5cSMatthew G. Knepley   PetscCall(DMLocatePoints(dm, pos, DM_POINTLOCATION_NONE, &sfcell));
5980e2ec84fSDave May 
5999566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell));
6000e2ec84fSDave May   n_found = 0;
6010e2ec84fSDave May   for (p = 0; p < n_estimate; p++) {
602ad540459SPierre Jolivet     if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) n_found++;
6030e2ec84fSDave May   }
6040e2ec84fSDave May 
6050e2ec84fSDave May   /* adjust size */
6060e2ec84fSDave May   if (mode == ADD_VALUES) {
60719307e5cSMatthew G. Knepley     PetscCall(DMSwarmGetLocalSize(sw, &n_curr));
6080e2ec84fSDave May     n_new_est = n_curr + n_found;
60919307e5cSMatthew G. Knepley     PetscCall(DMSwarmSetLocalSizes(sw, n_new_est, -1));
6100e2ec84fSDave May   }
6110e2ec84fSDave May   if (mode == INSERT_VALUES) {
6120e2ec84fSDave May     n_curr    = 0;
6130e2ec84fSDave May     n_new_est = n_found;
61419307e5cSMatthew G. Knepley     PetscCall(DMSwarmSetLocalSizes(sw, n_new_est, -1));
6150e2ec84fSDave May   }
6160e2ec84fSDave May 
6170e2ec84fSDave May   /* initialize new coords, cell owners, pid */
61819307e5cSMatthew G. Knepley   PetscCall(DMSwarmCellDMGetCellID(celldm, &cellid));
6199566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(pos, &_coor));
62019307e5cSMatthew G. Knepley   PetscCall(DMSwarmGetField(sw, coordFields[0], NULL, NULL, (void **)&swarm_coor));
62119307e5cSMatthew G. Knepley   PetscCall(DMSwarmGetField(sw, cellid, NULL, NULL, (void **)&swarm_cellid));
6220e2ec84fSDave May   n_found = 0;
6230e2ec84fSDave May   for (p = 0; p < n_estimate; p++) {
6240e2ec84fSDave May     if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) {
625ad540459SPierre Jolivet       for (b = 0; b < bs; b++) swarm_coor[bs * (n_curr + n_found) + b] = PetscRealPart(_coor[bs * p + b]);
6260e2ec84fSDave May       swarm_cellid[n_curr + n_found] = LA_sfcell[p].index;
6270e2ec84fSDave May       n_found++;
6280e2ec84fSDave May     }
6290e2ec84fSDave May   }
63019307e5cSMatthew G. Knepley   PetscCall(DMSwarmRestoreField(sw, cellid, NULL, NULL, (void **)&swarm_cellid));
63119307e5cSMatthew G. Knepley   PetscCall(DMSwarmRestoreField(sw, coordFields[0], NULL, NULL, (void **)&swarm_coor));
6329566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(pos, &_coor));
6330e2ec84fSDave May 
6340e2ec84fSDave May   if (redundant) {
63548a46eb9SPierre Jolivet     if (rank > 0) PetscCall(PetscFree(my_coor));
6360e2ec84fSDave May   }
6379566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sfcell));
6389566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pos));
6393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6400e2ec84fSDave May }
6410e2ec84fSDave May 
6420e2ec84fSDave May extern PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_DA(DM, DM, DMSwarmPICLayoutType, PetscInt);
6430e2ec84fSDave May extern PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX(DM, DM, DMSwarmPICLayoutType, PetscInt);
6440e2ec84fSDave May 
645cc4c1da9SBarry Smith /*@
6460e2ec84fSDave May   DMSwarmInsertPointsUsingCellDM - Insert point coordinates within each cell
6470e2ec84fSDave May 
64820f4b53cSBarry Smith   Not Collective
6490e2ec84fSDave May 
65060225df5SJacob Faibussowitsch   Input Parameters:
65120f4b53cSBarry Smith + dm          - the `DMSWARM`
65220f4b53cSBarry Smith . layout_type - method used to fill each cell with the cell `DM`
6530e2ec84fSDave May - fill_param  - parameter controlling how many points per cell are added (the meaning of this parameter is dependent on the layout type)
6540e2ec84fSDave May 
6550e2ec84fSDave May   Level: beginner
6560e2ec84fSDave May 
6570e2ec84fSDave May   Notes:
65820f4b53cSBarry Smith   The insert method will reset any previous defined points within the `DMSWARM`.
659e7af74c8SDave May 
66020f4b53cSBarry Smith   When using a `DMDA` both 2D and 3D are supported for all layout types provided you are using `DMDA_ELEMENT_Q1`.
661e7af74c8SDave May 
662a4e35b19SJacob Faibussowitsch   When using a `DMPLEX` the following case are supported\:
66320f4b53cSBarry Smith .vb
664ea3d7275SDave May    (i) DMSWARMPIC_LAYOUT_REGULAR: 2D (triangle),
665ea3d7275SDave May    (ii) DMSWARMPIC_LAYOUT_GAUSS: 2D and 3D provided the cell is a tri/tet or a quad/hex,
666ea3d7275SDave May    (iii) DMSWARMPIC_LAYOUT_SUBDIVISION: 2D and 3D for quad/hex and 2D tri.
66720f4b53cSBarry Smith .ve
6680e2ec84fSDave May 
66920f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmPICLayoutType`, `DMSwarmSetType()`, `DMSwarmSetCellDM()`, `DMSwarmType`
6700e2ec84fSDave May @*/
671cc4c1da9SBarry Smith PetscErrorCode DMSwarmInsertPointsUsingCellDM(DM dm, DMSwarmPICLayoutType layout_type, PetscInt fill_param)
672d71ae5a4SJacob Faibussowitsch {
6730e2ec84fSDave May   DM        celldm;
6740e2ec84fSDave May   PetscBool isDA, isPLEX;
6750e2ec84fSDave May 
6760e2ec84fSDave May   PetscFunctionBegin;
6770e2ec84fSDave May   DMSWARMPICVALID(dm);
6789566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetCellDM(dm, &celldm));
6799566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMDA, &isDA));
6809566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMPLEX, &isPLEX));
6810e2ec84fSDave May   if (isDA) {
6829566063dSJacob Faibussowitsch     PetscCall(private_DMSwarmInsertPointsUsingCellDM_DA(dm, celldm, layout_type, fill_param));
6830e2ec84fSDave May   } else if (isPLEX) {
6849566063dSJacob Faibussowitsch     PetscCall(private_DMSwarmInsertPointsUsingCellDM_PLEX(dm, celldm, layout_type, fill_param));
6850e2ec84fSDave May   } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only supported for cell DMs of type DMDA and DMPLEX");
6863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6870e2ec84fSDave May }
6880e2ec84fSDave May 
689431382f2SDave May extern PetscErrorCode private_DMSwarmSetPointCoordinatesCellwise_PLEX(DM, DM, PetscInt, PetscReal *);
690431382f2SDave May 
691431382f2SDave May /*@C
692431382f2SDave May   DMSwarmSetPointCoordinatesCellwise - Insert point coordinates (defined over the reference cell) within each cell
693431382f2SDave May 
69420f4b53cSBarry Smith   Not Collective
695431382f2SDave May 
69660225df5SJacob Faibussowitsch   Input Parameters:
69720f4b53cSBarry Smith + dm      - the `DMSWARM`
698431382f2SDave May . npoints - the number of points to insert in each cell
699431382f2SDave May - xi      - the coordinates (defined in the local coordinate system for each cell) to insert
700431382f2SDave May 
701431382f2SDave May   Level: beginner
702431382f2SDave May 
703431382f2SDave May   Notes:
70420f4b53cSBarry Smith   The method will reset any previous defined points within the `DMSWARM`.
70520f4b53cSBarry Smith   Only supported for `DMPLEX`. If you are using a `DMDA` it is recommended to either use
70620f4b53cSBarry Smith   `DMSwarmInsertPointsUsingCellDM()`, or extract and set the coordinates yourself the following code
70720f4b53cSBarry Smith .vb
70820f4b53cSBarry Smith     PetscReal  *coor;
709d52c2f21SMatthew G. Knepley     const char *coordname;
710d52c2f21SMatthew G. Knepley     DMSwarmGetCoordinateField(dm, &coordname);
711d52c2f21SMatthew G. Knepley     DMSwarmGetField(dm,coordname,NULL,NULL,(void**)&coor);
71220f4b53cSBarry Smith     // user code to define the coordinates here
713d52c2f21SMatthew G. Knepley     DMSwarmRestoreField(dm,coordname,NULL,NULL,(void**)&coor);
71420f4b53cSBarry Smith .ve
715e7af74c8SDave May 
71620f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmSetCellDM()`, `DMSwarmInsertPointsUsingCellDM()`
717431382f2SDave May @*/
718cc4c1da9SBarry Smith PetscErrorCode DMSwarmSetPointCoordinatesCellwise(DM dm, PetscInt npoints, PetscReal xi[])
719d71ae5a4SJacob Faibussowitsch {
720431382f2SDave May   DM        celldm;
721431382f2SDave May   PetscBool isDA, isPLEX;
722431382f2SDave May 
7230e2ec84fSDave May   PetscFunctionBegin;
724431382f2SDave May   DMSWARMPICVALID(dm);
7259566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetCellDM(dm, &celldm));
7269566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMDA, &isDA));
7279566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMPLEX, &isPLEX));
72828b400f6SJacob Faibussowitsch   PetscCheck(!isDA, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only supported for cell DMs of type DMPLEX. Recommended you use DMSwarmInsertPointsUsingCellDM()");
729*966bd95aSPierre Jolivet   PetscCheck(isPLEX, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only supported for cell DMs of type DMDA and DMPLEX");
7309566063dSJacob Faibussowitsch   PetscCall(private_DMSwarmSetPointCoordinatesCellwise_PLEX(dm, celldm, npoints, xi));
7313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7320e2ec84fSDave May }
733431382f2SDave May 
734cc4c1da9SBarry Smith /*@
735b799feefSDave May   DMSwarmCreatePointPerCellCount - Count the number of points within all cells in the cell DM
7360e2ec84fSDave May 
73720f4b53cSBarry Smith   Not Collective
7380e2ec84fSDave May 
73960225df5SJacob Faibussowitsch   Input Parameter:
74019307e5cSMatthew G. Knepley . sw - the `DMSWARM`
7410e2ec84fSDave May 
74260225df5SJacob Faibussowitsch   Output Parameters:
74320f4b53cSBarry Smith + ncells - the number of cells in the cell `DM` (optional argument, pass `NULL` to ignore)
744b799feefSDave May - count  - array of length ncells containing the number of points per cell
7450e2ec84fSDave May 
7460e2ec84fSDave May   Level: beginner
7470e2ec84fSDave May 
7480e2ec84fSDave May   Notes:
7490e2ec84fSDave May   The array count is allocated internally and must be free'd by the user.
7500e2ec84fSDave May 
75120f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmSetType()`, `DMSwarmSetCellDM()`, `DMSwarmType`
7520e2ec84fSDave May @*/
75319307e5cSMatthew G. Knepley PetscErrorCode DMSwarmCreatePointPerCellCount(DM sw, PetscInt *ncells, PetscInt **count)
754d71ae5a4SJacob Faibussowitsch {
75519307e5cSMatthew G. Knepley   DMSwarmCellDM celldm;
756b799feefSDave May   PetscBool     isvalid;
757b799feefSDave May   PetscInt      nel;
758b799feefSDave May   PetscInt     *sum;
75919307e5cSMatthew G. Knepley   const char   *cellid;
760b799feefSDave May 
7610e2ec84fSDave May   PetscFunctionBegin;
76219307e5cSMatthew G. Knepley   PetscCall(DMSwarmSortGetIsValid(sw, &isvalid));
763b799feefSDave May   nel = 0;
764b799feefSDave May   if (isvalid) {
765b799feefSDave May     PetscInt e;
766b799feefSDave May 
76719307e5cSMatthew G. Knepley     PetscCall(DMSwarmSortGetSizes(sw, &nel, NULL));
768b799feefSDave May 
7699566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nel, &sum));
77019307e5cSMatthew G. Knepley     for (e = 0; e < nel; e++) PetscCall(DMSwarmSortGetNumberOfPointsPerCell(sw, e, &sum[e]));
771b799feefSDave May   } else {
77219307e5cSMatthew G. Knepley     DM        dm;
773b799feefSDave May     PetscBool isda, isplex, isshell;
774b799feefSDave May     PetscInt  p, npoints;
775b799feefSDave May     PetscInt *swarm_cellid;
776b799feefSDave May 
777b799feefSDave May     /* get the number of cells */
77819307e5cSMatthew G. Knepley     PetscCall(DMSwarmGetCellDMActive(sw, &celldm));
77919307e5cSMatthew G. Knepley     PetscCall(DMSwarmCellDMGetDM(celldm, &dm));
78019307e5cSMatthew G. Knepley     PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMDA, &isda));
78119307e5cSMatthew G. Knepley     PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMPLEX, &isplex));
78219307e5cSMatthew G. Knepley     PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
783b799feefSDave May     if (isda) {
784b799feefSDave May       PetscInt        _nel, _npe;
785b799feefSDave May       const PetscInt *_element;
786b799feefSDave May 
78719307e5cSMatthew G. Knepley       PetscCall(DMDAGetElements(dm, &_nel, &_npe, &_element));
788b799feefSDave May       nel = _nel;
78919307e5cSMatthew G. Knepley       PetscCall(DMDARestoreElements(dm, &_nel, &_npe, &_element));
790b799feefSDave May     } else if (isplex) {
791b799feefSDave May       PetscInt ps, pe;
792b799feefSDave May 
79319307e5cSMatthew G. Knepley       PetscCall(DMPlexGetHeightStratum(dm, 0, &ps, &pe));
794b799feefSDave May       nel = pe - ps;
795b799feefSDave May     } else if (isshell) {
796b799feefSDave May       PetscErrorCode (*method_DMShellGetNumberOfCells)(DM, PetscInt *);
797b799feefSDave May 
79819307e5cSMatthew G. Knepley       PetscCall(PetscObjectQueryFunction((PetscObject)dm, "DMGetNumberOfCells_C", &method_DMShellGetNumberOfCells));
799b799feefSDave May       if (method_DMShellGetNumberOfCells) {
80019307e5cSMatthew G. Knepley         PetscCall(method_DMShellGetNumberOfCells(dm, &nel));
8019371c9d4SSatish Balay       } else
80219307e5cSMatthew G. Knepley         SETERRQ(PetscObjectComm((PetscObject)sw), PETSC_ERR_SUP, "Cannot determine the number of cells for the DMSHELL object. User must provide a method via PetscObjectComposeFunction( (PetscObject)shelldm, \"DMGetNumberOfCells_C\", your_function_to_compute_number_of_cells);");
80319307e5cSMatthew G. Knepley     } else SETERRQ(PetscObjectComm((PetscObject)sw), PETSC_ERR_SUP, "Cannot determine the number of cells for a DM not of type DA, PLEX or SHELL");
804b799feefSDave May 
8059566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nel, &sum));
8069566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(sum, nel));
80719307e5cSMatthew G. Knepley     PetscCall(DMSwarmGetLocalSize(sw, &npoints));
80819307e5cSMatthew G. Knepley     PetscCall(DMSwarmCellDMGetCellID(celldm, &cellid));
80919307e5cSMatthew G. Knepley     PetscCall(DMSwarmGetField(sw, cellid, NULL, NULL, (void **)&swarm_cellid));
810b799feefSDave May     for (p = 0; p < npoints; p++) {
811ad540459SPierre Jolivet       if (swarm_cellid[p] != DMLOCATEPOINT_POINT_NOT_FOUND) sum[swarm_cellid[p]]++;
812b799feefSDave May     }
81319307e5cSMatthew G. Knepley     PetscCall(DMSwarmRestoreField(sw, cellid, NULL, NULL, (void **)&swarm_cellid));
814b799feefSDave May   }
815ad540459SPierre Jolivet   if (ncells) *ncells = nel;
816b799feefSDave May   *count = sum;
8173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8180e2ec84fSDave May }
81935b38c60SMatthew G. Knepley 
82035b38c60SMatthew G. Knepley /*@
82135b38c60SMatthew G. Knepley   DMSwarmGetNumSpecies - Get the number of particle species
82235b38c60SMatthew G. Knepley 
82320f4b53cSBarry Smith   Not Collective
82435b38c60SMatthew G. Knepley 
82560225df5SJacob Faibussowitsch   Input Parameter:
82660225df5SJacob Faibussowitsch . sw - the `DMSWARM`
82735b38c60SMatthew G. Knepley 
82860225df5SJacob Faibussowitsch   Output Parameters:
82935b38c60SMatthew G. Knepley . Ns - the number of species
83035b38c60SMatthew G. Knepley 
83135b38c60SMatthew G. Knepley   Level: intermediate
83235b38c60SMatthew G. Knepley 
83320f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmSetNumSpecies()`, `DMSwarmSetType()`, `DMSwarmType`
83435b38c60SMatthew G. Knepley @*/
835d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmGetNumSpecies(DM sw, PetscInt *Ns)
836d71ae5a4SJacob Faibussowitsch {
83735b38c60SMatthew G. Knepley   DM_Swarm *swarm = (DM_Swarm *)sw->data;
83835b38c60SMatthew G. Knepley 
83935b38c60SMatthew G. Knepley   PetscFunctionBegin;
84035b38c60SMatthew G. Knepley   *Ns = swarm->Ns;
8413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
84235b38c60SMatthew G. Knepley }
84335b38c60SMatthew G. Knepley 
84435b38c60SMatthew G. Knepley /*@
84535b38c60SMatthew G. Knepley   DMSwarmSetNumSpecies - Set the number of particle species
84635b38c60SMatthew G. Knepley 
84720f4b53cSBarry Smith   Not Collective
84835b38c60SMatthew G. Knepley 
84960225df5SJacob Faibussowitsch   Input Parameters:
85060225df5SJacob Faibussowitsch + sw - the `DMSWARM`
85135b38c60SMatthew G. Knepley - Ns - the number of species
85235b38c60SMatthew G. Knepley 
85335b38c60SMatthew G. Knepley   Level: intermediate
85435b38c60SMatthew G. Knepley 
85520f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmGetNumSpecies()`, `DMSwarmSetType()`, `DMSwarmType`
85635b38c60SMatthew G. Knepley @*/
857d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmSetNumSpecies(DM sw, PetscInt Ns)
858d71ae5a4SJacob Faibussowitsch {
85935b38c60SMatthew G. Knepley   DM_Swarm *swarm = (DM_Swarm *)sw->data;
86035b38c60SMatthew G. Knepley 
86135b38c60SMatthew G. Knepley   PetscFunctionBegin;
86235b38c60SMatthew G. Knepley   swarm->Ns = Ns;
8633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
86435b38c60SMatthew G. Knepley }
86535b38c60SMatthew G. Knepley 
86635b38c60SMatthew G. Knepley /*@C
8676c5a79ebSMatthew G. Knepley   DMSwarmGetCoordinateFunction - Get the function setting initial particle positions, if it exists
8686c5a79ebSMatthew G. Knepley 
86920f4b53cSBarry Smith   Not Collective
8706c5a79ebSMatthew G. Knepley 
87160225df5SJacob Faibussowitsch   Input Parameter:
87260225df5SJacob Faibussowitsch . sw - the `DMSWARM`
8736c5a79ebSMatthew G. Knepley 
8746c5a79ebSMatthew G. Knepley   Output Parameter:
8758434afd1SBarry Smith . coordFunc - the function setting initial particle positions, or `NULL`, see `PetscSimplePointFn` for the calling sequence
8766c5a79ebSMatthew G. Knepley 
8776c5a79ebSMatthew G. Knepley   Level: intermediate
8786c5a79ebSMatthew G. Knepley 
8798434afd1SBarry Smith .seealso: `DMSWARM`, `DMSwarmSetCoordinateFunction()`, `DMSwarmGetVelocityFunction()`, `DMSwarmInitializeCoordinates()`, `PetscSimplePointFn`
8806c5a79ebSMatthew G. Knepley @*/
8818434afd1SBarry Smith PetscErrorCode DMSwarmGetCoordinateFunction(DM sw, PetscSimplePointFn **coordFunc)
882d71ae5a4SJacob Faibussowitsch {
8836c5a79ebSMatthew G. Knepley   DM_Swarm *swarm = (DM_Swarm *)sw->data;
8846c5a79ebSMatthew G. Knepley 
8856c5a79ebSMatthew G. Knepley   PetscFunctionBegin;
8866c5a79ebSMatthew G. Knepley   PetscValidHeaderSpecific(sw, DM_CLASSID, 1);
8874f572ea9SToby Isaac   PetscAssertPointer(coordFunc, 2);
8886c5a79ebSMatthew G. Knepley   *coordFunc = swarm->coordFunc;
8893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8906c5a79ebSMatthew G. Knepley }
8916c5a79ebSMatthew G. Knepley 
8926c5a79ebSMatthew G. Knepley /*@C
8936c5a79ebSMatthew G. Knepley   DMSwarmSetCoordinateFunction - Set the function setting initial particle positions
8946c5a79ebSMatthew G. Knepley 
89520f4b53cSBarry Smith   Not Collective
8966c5a79ebSMatthew G. Knepley 
89760225df5SJacob Faibussowitsch   Input Parameters:
89860225df5SJacob Faibussowitsch + sw        - the `DMSWARM`
8998434afd1SBarry Smith - coordFunc - the function setting initial particle positions, see `PetscSimplePointFn` for the calling sequence
9006c5a79ebSMatthew G. Knepley 
9016c5a79ebSMatthew G. Knepley   Level: intermediate
9026c5a79ebSMatthew G. Knepley 
9038434afd1SBarry Smith .seealso: `DMSWARM`, `DMSwarmGetCoordinateFunction()`, `DMSwarmSetVelocityFunction()`, `DMSwarmInitializeCoordinates()`, `PetscSimplePointFn`
9046c5a79ebSMatthew G. Knepley @*/
9058434afd1SBarry Smith PetscErrorCode DMSwarmSetCoordinateFunction(DM sw, PetscSimplePointFn *coordFunc)
906d71ae5a4SJacob Faibussowitsch {
9076c5a79ebSMatthew G. Knepley   DM_Swarm *swarm = (DM_Swarm *)sw->data;
9086c5a79ebSMatthew G. Knepley 
9096c5a79ebSMatthew G. Knepley   PetscFunctionBegin;
9106c5a79ebSMatthew G. Knepley   PetscValidHeaderSpecific(sw, DM_CLASSID, 1);
9116c5a79ebSMatthew G. Knepley   PetscValidFunction(coordFunc, 2);
9126c5a79ebSMatthew G. Knepley   swarm->coordFunc = coordFunc;
9133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9146c5a79ebSMatthew G. Knepley }
9156c5a79ebSMatthew G. Knepley 
9166c5a79ebSMatthew G. Knepley /*@C
91760225df5SJacob Faibussowitsch   DMSwarmGetVelocityFunction - Get the function setting initial particle velocities, if it exists
9186c5a79ebSMatthew G. Knepley 
91920f4b53cSBarry Smith   Not Collective
9206c5a79ebSMatthew G. Knepley 
92160225df5SJacob Faibussowitsch   Input Parameter:
92260225df5SJacob Faibussowitsch . sw - the `DMSWARM`
9236c5a79ebSMatthew G. Knepley 
9246c5a79ebSMatthew G. Knepley   Output Parameter:
9258434afd1SBarry Smith . velFunc - the function setting initial particle velocities, or `NULL`, see `PetscSimplePointFn` for the calling sequence
9266c5a79ebSMatthew G. Knepley 
9276c5a79ebSMatthew G. Knepley   Level: intermediate
9286c5a79ebSMatthew G. Knepley 
9298434afd1SBarry Smith .seealso: `DMSWARM`, `DMSwarmSetVelocityFunction()`, `DMSwarmGetCoordinateFunction()`, `DMSwarmInitializeVelocities()`, `PetscSimplePointFn`
9306c5a79ebSMatthew G. Knepley @*/
9318434afd1SBarry Smith PetscErrorCode DMSwarmGetVelocityFunction(DM sw, PetscSimplePointFn **velFunc)
932d71ae5a4SJacob Faibussowitsch {
9336c5a79ebSMatthew G. Knepley   DM_Swarm *swarm = (DM_Swarm *)sw->data;
9346c5a79ebSMatthew G. Knepley 
9356c5a79ebSMatthew G. Knepley   PetscFunctionBegin;
9366c5a79ebSMatthew G. Knepley   PetscValidHeaderSpecific(sw, DM_CLASSID, 1);
9374f572ea9SToby Isaac   PetscAssertPointer(velFunc, 2);
9386c5a79ebSMatthew G. Knepley   *velFunc = swarm->velFunc;
9393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9406c5a79ebSMatthew G. Knepley }
9416c5a79ebSMatthew G. Knepley 
9426c5a79ebSMatthew G. Knepley /*@C
9436c5a79ebSMatthew G. Knepley   DMSwarmSetVelocityFunction - Set the function setting initial particle velocities
9446c5a79ebSMatthew G. Knepley 
94520f4b53cSBarry Smith   Not Collective
9466c5a79ebSMatthew G. Knepley 
94760225df5SJacob Faibussowitsch   Input Parameters:
948a4e35b19SJacob Faibussowitsch + sw      - the `DMSWARM`
9498434afd1SBarry Smith - velFunc - the function setting initial particle velocities, see `PetscSimplePointFn` for the calling sequence
9506c5a79ebSMatthew G. Knepley 
9516c5a79ebSMatthew G. Knepley   Level: intermediate
9526c5a79ebSMatthew G. Knepley 
9538434afd1SBarry Smith .seealso: `DMSWARM`, `DMSwarmGetVelocityFunction()`, `DMSwarmSetCoordinateFunction()`, `DMSwarmInitializeVelocities()`, `PetscSimplePointFn`
9546c5a79ebSMatthew G. Knepley @*/
9558434afd1SBarry Smith PetscErrorCode DMSwarmSetVelocityFunction(DM sw, PetscSimplePointFn *velFunc)
956d71ae5a4SJacob Faibussowitsch {
9576c5a79ebSMatthew G. Knepley   DM_Swarm *swarm = (DM_Swarm *)sw->data;
9586c5a79ebSMatthew G. Knepley 
9596c5a79ebSMatthew G. Knepley   PetscFunctionBegin;
9606c5a79ebSMatthew G. Knepley   PetscValidHeaderSpecific(sw, DM_CLASSID, 1);
9616c5a79ebSMatthew G. Knepley   PetscValidFunction(velFunc, 2);
9626c5a79ebSMatthew G. Knepley   swarm->velFunc = velFunc;
9633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9646c5a79ebSMatthew G. Knepley }
9656c5a79ebSMatthew G. Knepley 
9666c5a79ebSMatthew G. Knepley /*@C
96735b38c60SMatthew G. Knepley   DMSwarmComputeLocalSize - Compute the local number and distribution of particles based upon a density function
96835b38c60SMatthew G. Knepley 
96920f4b53cSBarry Smith   Not Collective
97035b38c60SMatthew G. Knepley 
97135b38c60SMatthew G. Knepley   Input Parameters:
97220f4b53cSBarry Smith + sw      - The `DMSWARM`
97335b38c60SMatthew G. Knepley . N       - The target number of particles
97435b38c60SMatthew G. Knepley - density - The density field for the particle layout, normalized to unity
97535b38c60SMatthew G. Knepley 
97635b38c60SMatthew G. Knepley   Level: advanced
97735b38c60SMatthew G. Knepley 
97820f4b53cSBarry Smith   Note:
97920f4b53cSBarry Smith   One particle will be created for each species.
98020f4b53cSBarry Smith 
98120f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmComputeLocalSizeFromOptions()`
98235b38c60SMatthew G. Knepley @*/
983f8662bd6SBarry Smith PetscErrorCode DMSwarmComputeLocalSize(DM sw, PetscInt N, PetscProbFn *density)
984d71ae5a4SJacob Faibussowitsch {
98535b38c60SMatthew G. Knepley   DM               dm;
98619307e5cSMatthew G. Knepley   DMSwarmCellDM    celldm;
98735b38c60SMatthew G. Knepley   PetscQuadrature  quad;
98835b38c60SMatthew G. Knepley   const PetscReal *xq, *wq;
989ea7032a0SMatthew G. Knepley   PetscReal       *n_int;
99019307e5cSMatthew G. Knepley   PetscInt        *npc_s, *swarm_cellid, Ni;
991ea7032a0SMatthew G. Knepley   PetscReal        gmin[3], gmax[3], xi0[3];
992ea7032a0SMatthew G. Knepley   PetscInt         Ns, cStart, cEnd, c, dim, d, Nq, q, Np = 0, p, s;
99335b38c60SMatthew G. Knepley   PetscBool        simplex;
99419307e5cSMatthew G. Knepley   const char      *cellid;
99535b38c60SMatthew G. Knepley 
99635b38c60SMatthew G. Knepley   PetscFunctionBegin;
9979566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetNumSpecies(sw, &Ns));
9989566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetCellDM(sw, &dm));
9999566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1000ea7032a0SMatthew G. Knepley   PetscCall(DMGetBoundingBox(dm, gmin, gmax));
10019566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
10029566063dSJacob Faibussowitsch   PetscCall(DMPlexIsSimplex(dm, &simplex));
10036858538eSMatthew G. Knepley   PetscCall(DMGetCoordinatesLocalSetUp(dm));
10049566063dSJacob Faibussowitsch   if (simplex) PetscCall(PetscDTStroudConicalQuadrature(dim, 1, 5, -1.0, 1.0, &quad));
10059566063dSJacob Faibussowitsch   else PetscCall(PetscDTGaussTensorQuadrature(dim, 1, 5, -1.0, 1.0, &quad));
10069566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, &xq, &wq));
1007ea7032a0SMatthew G. Knepley   PetscCall(PetscCalloc2(Ns, &n_int, (cEnd - cStart) * Ns, &npc_s));
100835b38c60SMatthew G. Knepley   /* Integrate the density function to get the number of particles in each cell */
100935b38c60SMatthew G. Knepley   for (d = 0; d < dim; ++d) xi0[d] = -1.0;
101035b38c60SMatthew G. Knepley   for (c = 0; c < cEnd - cStart; ++c) {
101135b38c60SMatthew G. Knepley     const PetscInt cell = c + cStart;
1012ea7032a0SMatthew G. Knepley     PetscReal      v0[3], J[9], invJ[9], detJ, detJp = 2. / (gmax[0] - gmin[0]), xr[3], den;
101335b38c60SMatthew G. Knepley 
1014ea7032a0SMatthew G. Knepley     /* Have to transform quadrature points/weights to cell domain */
10159566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeCellGeometryFEM(dm, cell, NULL, v0, J, invJ, &detJ));
1016ea7032a0SMatthew G. Knepley     PetscCall(PetscArrayzero(n_int, Ns));
101735b38c60SMatthew G. Knepley     for (q = 0; q < Nq; ++q) {
101835b38c60SMatthew G. Knepley       CoordinatesRefToReal(dim, dim, xi0, v0, J, &xq[q * dim], xr);
1019ea7032a0SMatthew G. Knepley       /* Have to transform mesh to domain of definition of PDF, [-1, 1], and weight PDF by |J|/2 */
1020ea7032a0SMatthew G. Knepley       xr[0] = detJp * (xr[0] - gmin[0]) - 1.;
1021ea7032a0SMatthew G. Knepley 
1022ea7032a0SMatthew G. Knepley       for (s = 0; s < Ns; ++s) {
1023ea7032a0SMatthew G. Knepley         PetscCall(density(xr, NULL, &den));
1024ea7032a0SMatthew G. Knepley         n_int[s] += (detJp * den) * (detJ * wq[q]) / (PetscReal)Ns;
102535b38c60SMatthew G. Knepley       }
1026ea7032a0SMatthew G. Knepley     }
1027ea7032a0SMatthew G. Knepley     for (s = 0; s < Ns; ++s) {
1028ea7032a0SMatthew G. Knepley       Ni = N;
1029f940b0e3Sdanofinn       npc_s[c * Ns + s] += (PetscInt)(Ni * n_int[s] + 0.5); // TODO Wish we wrapped round()
1030ea7032a0SMatthew G. Knepley       Np += npc_s[c * Ns + s];
1031ea7032a0SMatthew G. Knepley     }
103235b38c60SMatthew G. Knepley   }
10339566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureDestroy(&quad));
10349566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(sw, Np, 0));
103519307e5cSMatthew G. Knepley   PetscCall(DMSwarmGetCellDMActive(sw, &celldm));
103619307e5cSMatthew G. Knepley   PetscCall(DMSwarmCellDMGetCellID(celldm, &cellid));
103719307e5cSMatthew G. Knepley   PetscCall(DMSwarmGetField(sw, cellid, NULL, NULL, (void **)&swarm_cellid));
103835b38c60SMatthew G. Knepley   for (c = 0, p = 0; c < cEnd - cStart; ++c) {
1039ea7032a0SMatthew G. Knepley     for (s = 0; s < Ns; ++s) {
104019307e5cSMatthew G. Knepley       for (q = 0; q < npc_s[c * Ns + s]; ++q, ++p) swarm_cellid[p] = c;
1041ea7032a0SMatthew G. Knepley     }
104235b38c60SMatthew G. Knepley   }
104319307e5cSMatthew G. Knepley   PetscCall(DMSwarmRestoreField(sw, cellid, NULL, NULL, (void **)&swarm_cellid));
1044ea7032a0SMatthew G. Knepley   PetscCall(PetscFree2(n_int, npc_s));
10453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
104635b38c60SMatthew G. Knepley }
104735b38c60SMatthew G. Knepley 
104835b38c60SMatthew G. Knepley /*@
104935b38c60SMatthew G. Knepley   DMSwarmComputeLocalSizeFromOptions - Compute the local number and distribution of particles based upon a density function determined by options
105035b38c60SMatthew G. Knepley 
105120f4b53cSBarry Smith   Not Collective
105235b38c60SMatthew G. Knepley 
10532fe279fdSBarry Smith   Input Parameter:
10542fe279fdSBarry Smith . sw - The `DMSWARM`
105535b38c60SMatthew G. Knepley 
105635b38c60SMatthew G. Knepley   Level: advanced
105735b38c60SMatthew G. Knepley 
105820f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmComputeLocalSize()`
105935b38c60SMatthew G. Knepley @*/
1060d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmComputeLocalSizeFromOptions(DM sw)
1061d71ae5a4SJacob Faibussowitsch {
1062f8662bd6SBarry Smith   PetscProbFn *pdf;
106335b38c60SMatthew G. Knepley   const char  *prefix;
10646c5a79ebSMatthew G. Knepley   char         funcname[PETSC_MAX_PATH_LEN];
10656c5a79ebSMatthew G. Knepley   PetscInt    *N, Ns, dim, n;
10666c5a79ebSMatthew G. Knepley   PetscBool    flg;
10676c5a79ebSMatthew G. Knepley   PetscMPIInt  size, rank;
106835b38c60SMatthew G. Knepley 
106935b38c60SMatthew G. Knepley   PetscFunctionBegin;
10706c5a79ebSMatthew G. Knepley   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)sw), &size));
10716c5a79ebSMatthew G. Knepley   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sw), &rank));
10726c5a79ebSMatthew G. Knepley   PetscCall(PetscCalloc1(size, &N));
1073d0609cedSBarry Smith   PetscOptionsBegin(PetscObjectComm((PetscObject)sw), "", "DMSwarm Options", "DMSWARM");
10746c5a79ebSMatthew G. Knepley   n = size;
10756c5a79ebSMatthew G. Knepley   PetscCall(PetscOptionsIntArray("-dm_swarm_num_particles", "The target number of particles", "", N, &n, NULL));
10766c5a79ebSMatthew G. Knepley   PetscCall(DMSwarmGetNumSpecies(sw, &Ns));
10779566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-dm_swarm_num_species", "The number of species", "DMSwarmSetNumSpecies", Ns, &Ns, &flg));
10789566063dSJacob Faibussowitsch   if (flg) PetscCall(DMSwarmSetNumSpecies(sw, Ns));
10796c5a79ebSMatthew G. Knepley   PetscCall(PetscOptionsString("-dm_swarm_coordinate_function", "Function to determine particle coordinates", "DMSwarmSetCoordinateFunction", funcname, funcname, sizeof(funcname), &flg));
1080d0609cedSBarry Smith   PetscOptionsEnd();
10816c5a79ebSMatthew G. Knepley   if (flg) {
10828434afd1SBarry Smith     PetscSimplePointFn *coordFunc;
108335b38c60SMatthew G. Knepley 
10846c5a79ebSMatthew G. Knepley     PetscCall(DMSwarmGetNumSpecies(sw, &Ns));
10856c5a79ebSMatthew G. Knepley     PetscCall(PetscDLSym(NULL, funcname, (void **)&coordFunc));
10866c5a79ebSMatthew G. Knepley     PetscCheck(coordFunc, PetscObjectComm((PetscObject)sw), PETSC_ERR_ARG_WRONG, "Could not locate function %s", funcname);
10876c5a79ebSMatthew G. Knepley     PetscCall(DMSwarmGetNumSpecies(sw, &Ns));
10886c5a79ebSMatthew G. Knepley     PetscCall(DMSwarmSetLocalSizes(sw, N[rank] * Ns, 0));
10896c5a79ebSMatthew G. Knepley     PetscCall(DMSwarmSetCoordinateFunction(sw, coordFunc));
10906c5a79ebSMatthew G. Knepley   } else {
10919566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(sw, &dim));
10929566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetOptionsPrefix((PetscObject)sw, &prefix));
10939566063dSJacob Faibussowitsch     PetscCall(PetscProbCreateFromOptions(dim, prefix, "-dm_swarm_coordinate_density", &pdf, NULL, NULL));
10946c5a79ebSMatthew G. Knepley     PetscCall(DMSwarmComputeLocalSize(sw, N[rank], pdf));
10956c5a79ebSMatthew G. Knepley   }
10966c5a79ebSMatthew G. Knepley   PetscCall(PetscFree(N));
10973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
109835b38c60SMatthew G. Knepley }
109935b38c60SMatthew G. Knepley 
110035b38c60SMatthew G. Knepley /*@
110135b38c60SMatthew G. Knepley   DMSwarmInitializeCoordinates - Determine the initial coordinates of particles for a PIC method
110235b38c60SMatthew G. Knepley 
110320f4b53cSBarry Smith   Not Collective
110435b38c60SMatthew G. Knepley 
110520f4b53cSBarry Smith   Input Parameter:
110620f4b53cSBarry Smith . sw - The `DMSWARM`
110735b38c60SMatthew G. Knepley 
110835b38c60SMatthew G. Knepley   Level: advanced
110935b38c60SMatthew G. Knepley 
111020f4b53cSBarry Smith   Note:
111120f4b53cSBarry Smith   Currently, we randomly place particles in their assigned cell
111220f4b53cSBarry Smith 
111320f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmComputeLocalSize()`, `DMSwarmInitializeVelocities()`
111435b38c60SMatthew G. Knepley @*/
1115d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmInitializeCoordinates(DM sw)
1116d71ae5a4SJacob Faibussowitsch {
111719307e5cSMatthew G. Knepley   DMSwarmCellDM       celldm;
11188434afd1SBarry Smith   PetscSimplePointFn *coordFunc;
111935b38c60SMatthew G. Knepley   PetscScalar        *weight;
11206c5a79ebSMatthew G. Knepley   PetscReal          *x;
112135b38c60SMatthew G. Knepley   PetscInt           *species;
11226c5a79ebSMatthew G. Knepley   void               *ctx;
112335b38c60SMatthew G. Knepley   PetscBool           removePoints = PETSC_TRUE;
112435b38c60SMatthew G. Knepley   PetscDataType       dtype;
112519307e5cSMatthew G. Knepley   PetscInt            Nfc, Np, p, Ns, dim, d, bs;
112619307e5cSMatthew G. Knepley   const char        **coordFields;
112735b38c60SMatthew G. Knepley 
112835b38c60SMatthew G. Knepley   PetscFunctionBeginUser;
11296c5a79ebSMatthew G. Knepley   PetscCall(DMGetDimension(sw, &dim));
11306c5a79ebSMatthew G. Knepley   PetscCall(DMSwarmGetLocalSize(sw, &Np));
11319566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetNumSpecies(sw, &Ns));
11326c5a79ebSMatthew G. Knepley   PetscCall(DMSwarmGetCoordinateFunction(sw, &coordFunc));
11336c5a79ebSMatthew G. Knepley 
113419307e5cSMatthew G. Knepley   PetscCall(DMSwarmGetCellDMActive(sw, &celldm));
113519307e5cSMatthew G. Knepley   PetscCall(DMSwarmCellDMGetCoordinateFields(celldm, &Nfc, &coordFields));
113619307e5cSMatthew G. Knepley   PetscCheck(Nfc == 1, PetscObjectComm((PetscObject)sw), PETSC_ERR_SUP, "We only support a single coordinate field right now, not %" PetscInt_FMT, Nfc);
113719307e5cSMatthew G. Knepley 
113819307e5cSMatthew G. Knepley   PetscCall(DMSwarmGetField(sw, coordFields[0], &bs, &dtype, (void **)&x));
11396c5a79ebSMatthew G. Knepley   PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&weight));
11406c5a79ebSMatthew G. Knepley   PetscCall(DMSwarmGetField(sw, "species", NULL, NULL, (void **)&species));
11416c5a79ebSMatthew G. Knepley   if (coordFunc) {
11426c5a79ebSMatthew G. Knepley     PetscCall(DMGetApplicationContext(sw, &ctx));
11436c5a79ebSMatthew G. Knepley     for (p = 0; p < Np; ++p) {
11446c5a79ebSMatthew G. Knepley       PetscScalar X[3];
11456c5a79ebSMatthew G. Knepley 
11466c5a79ebSMatthew G. Knepley       PetscCall((*coordFunc)(dim, 0., NULL, p, X, ctx));
11476c5a79ebSMatthew G. Knepley       for (d = 0; d < dim; ++d) x[p * dim + d] = PetscRealPart(X[d]);
11486c5a79ebSMatthew G. Knepley       weight[p]  = 1.0;
11496c5a79ebSMatthew G. Knepley       species[p] = p % Ns;
11506c5a79ebSMatthew G. Knepley     }
11516c5a79ebSMatthew G. Knepley   } else {
11526c5a79ebSMatthew G. Knepley     DM          dm;
11536c5a79ebSMatthew G. Knepley     PetscRandom rnd;
11546c5a79ebSMatthew G. Knepley     PetscReal   xi0[3];
11556c5a79ebSMatthew G. Knepley     PetscInt    cStart, cEnd, c;
11566c5a79ebSMatthew G. Knepley 
11579566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetCellDM(sw, &dm));
11589566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
1159ea7032a0SMatthew G. Knepley     PetscCall(DMGetApplicationContext(sw, &ctx));
116035b38c60SMatthew G. Knepley 
116135b38c60SMatthew G. Knepley     /* Set particle position randomly in cell, set weights to 1 */
11629566063dSJacob Faibussowitsch     PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)dm), &rnd));
11639566063dSJacob Faibussowitsch     PetscCall(PetscRandomSetInterval(rnd, -1.0, 1.0));
11649566063dSJacob Faibussowitsch     PetscCall(PetscRandomSetFromOptions(rnd));
11659566063dSJacob Faibussowitsch     PetscCall(DMSwarmSortGetAccess(sw));
116635b38c60SMatthew G. Knepley     for (d = 0; d < dim; ++d) xi0[d] = -1.0;
116735b38c60SMatthew G. Knepley     for (c = cStart; c < cEnd; ++c) {
116835b38c60SMatthew G. Knepley       PetscReal v0[3], J[9], invJ[9], detJ;
116935b38c60SMatthew G. Knepley       PetscInt *pidx, Npc, q;
117035b38c60SMatthew G. Knepley 
11719566063dSJacob Faibussowitsch       PetscCall(DMSwarmSortGetPointsPerCell(sw, c, &Npc, &pidx));
11729566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ));
117335b38c60SMatthew G. Knepley       for (q = 0; q < Npc; ++q) {
117435b38c60SMatthew G. Knepley         const PetscInt p = pidx[q];
117535b38c60SMatthew G. Knepley         PetscReal      xref[3];
117635b38c60SMatthew G. Knepley 
11779566063dSJacob Faibussowitsch         for (d = 0; d < dim; ++d) PetscCall(PetscRandomGetValueReal(rnd, &xref[d]));
117835b38c60SMatthew G. Knepley         CoordinatesRefToReal(dim, dim, xi0, v0, J, xref, &x[p * dim]);
117935b38c60SMatthew G. Knepley 
1180ea7032a0SMatthew G. Knepley         weight[p]  = 1.0 / Np;
11816c5a79ebSMatthew G. Knepley         species[p] = p % Ns;
118235b38c60SMatthew G. Knepley       }
1183fe11354eSMatthew G. Knepley       PetscCall(DMSwarmSortRestorePointsPerCell(sw, c, &Npc, &pidx));
118435b38c60SMatthew G. Knepley     }
11859566063dSJacob Faibussowitsch     PetscCall(PetscRandomDestroy(&rnd));
11869566063dSJacob Faibussowitsch     PetscCall(DMSwarmSortRestoreAccess(sw));
11876c5a79ebSMatthew G. Knepley   }
118819307e5cSMatthew G. Knepley   PetscCall(DMSwarmRestoreField(sw, coordFields[0], NULL, NULL, (void **)&x));
11896c5a79ebSMatthew G. Knepley   PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&weight));
11909566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **)&species));
11916c5a79ebSMatthew G. Knepley 
11929566063dSJacob Faibussowitsch   PetscCall(DMSwarmMigrate(sw, removePoints));
11939566063dSJacob Faibussowitsch   PetscCall(DMLocalizeCoordinates(sw));
11943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
119535b38c60SMatthew G. Knepley }
119635b38c60SMatthew G. Knepley 
119735b38c60SMatthew G. Knepley /*@C
119835b38c60SMatthew G. Knepley   DMSwarmInitializeVelocities - Set the initial velocities of particles using a distribution.
119935b38c60SMatthew G. Knepley 
120020f4b53cSBarry Smith   Collective
120135b38c60SMatthew G. Knepley 
120235b38c60SMatthew G. Knepley   Input Parameters:
120320f4b53cSBarry Smith + sw      - The `DMSWARM` object
120435b38c60SMatthew G. Knepley . sampler - A function which uniformly samples the velocity PDF
120535b38c60SMatthew G. Knepley - v0      - The velocity scale for nondimensionalization for each species
120635b38c60SMatthew G. Knepley 
120735b38c60SMatthew G. Knepley   Level: advanced
120835b38c60SMatthew G. Knepley 
120920f4b53cSBarry Smith   Note:
121020f4b53cSBarry Smith   If `v0` is zero for the first species, all velocities are set to zero. If it is zero for any other species, the effect will be to give that species zero velocity.
121120f4b53cSBarry Smith 
121220f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmComputeLocalSize()`, `DMSwarmInitializeCoordinates()`, `DMSwarmInitializeVelocitiesFromOptions()`
121335b38c60SMatthew G. Knepley @*/
1214f8662bd6SBarry Smith PetscErrorCode DMSwarmInitializeVelocities(DM sw, PetscProbFn *sampler, const PetscReal v0[])
1215d71ae5a4SJacob Faibussowitsch {
12168434afd1SBarry Smith   PetscSimplePointFn *velFunc;
121735b38c60SMatthew G. Knepley   PetscReal          *v;
121835b38c60SMatthew G. Knepley   PetscInt           *species;
12196c5a79ebSMatthew G. Knepley   void               *ctx;
122035b38c60SMatthew G. Knepley   PetscInt            dim, Np, p;
122135b38c60SMatthew G. Knepley 
122235b38c60SMatthew G. Knepley   PetscFunctionBegin;
12236c5a79ebSMatthew G. Knepley   PetscCall(DMSwarmGetVelocityFunction(sw, &velFunc));
122435b38c60SMatthew G. Knepley 
12259566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(sw, &dim));
12269566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetLocalSize(sw, &Np));
12279566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(sw, "velocity", NULL, NULL, (void **)&v));
12289566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(sw, "species", NULL, NULL, (void **)&species));
12296c5a79ebSMatthew G. Knepley   if (v0[0] == 0.) {
12306c5a79ebSMatthew G. Knepley     PetscCall(PetscArrayzero(v, Np * dim));
12316c5a79ebSMatthew G. Knepley   } else if (velFunc) {
12326c5a79ebSMatthew G. Knepley     PetscCall(DMGetApplicationContext(sw, &ctx));
12336c5a79ebSMatthew G. Knepley     for (p = 0; p < Np; ++p) {
12346c5a79ebSMatthew G. Knepley       PetscInt    s = species[p], d;
12356c5a79ebSMatthew G. Knepley       PetscScalar vel[3];
12366c5a79ebSMatthew G. Knepley 
12376c5a79ebSMatthew G. Knepley       PetscCall((*velFunc)(dim, 0., NULL, p, vel, ctx));
12386c5a79ebSMatthew G. Knepley       for (d = 0; d < dim; ++d) v[p * dim + d] = (v0[s] / v0[0]) * PetscRealPart(vel[d]);
12396c5a79ebSMatthew G. Knepley     }
12406c5a79ebSMatthew G. Knepley   } else {
12416c5a79ebSMatthew G. Knepley     PetscRandom rnd;
12426c5a79ebSMatthew G. Knepley 
12436c5a79ebSMatthew G. Knepley     PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)sw), &rnd));
12446c5a79ebSMatthew G. Knepley     PetscCall(PetscRandomSetInterval(rnd, 0, 1.));
12456c5a79ebSMatthew G. Knepley     PetscCall(PetscRandomSetFromOptions(rnd));
12466c5a79ebSMatthew G. Knepley 
124735b38c60SMatthew G. Knepley     for (p = 0; p < Np; ++p) {
124835b38c60SMatthew G. Knepley       PetscInt  s = species[p], d;
124935b38c60SMatthew G. Knepley       PetscReal a[3], vel[3];
125035b38c60SMatthew G. Knepley 
12519566063dSJacob Faibussowitsch       for (d = 0; d < dim; ++d) PetscCall(PetscRandomGetValueReal(rnd, &a[d]));
12529566063dSJacob Faibussowitsch       PetscCall(sampler(a, NULL, vel));
1253ad540459SPierre Jolivet       for (d = 0; d < dim; ++d) v[p * dim + d] = (v0[s] / v0[0]) * vel[d];
125435b38c60SMatthew G. Knepley     }
12556c5a79ebSMatthew G. Knepley     PetscCall(PetscRandomDestroy(&rnd));
12566c5a79ebSMatthew G. Knepley   }
12579566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(sw, "velocity", NULL, NULL, (void **)&v));
12589566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **)&species));
12593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
126035b38c60SMatthew G. Knepley }
126135b38c60SMatthew G. Knepley 
126235b38c60SMatthew G. Knepley /*@
126335b38c60SMatthew G. Knepley   DMSwarmInitializeVelocitiesFromOptions - Set the initial velocities of particles using a distribution determined from options.
126435b38c60SMatthew G. Knepley 
126520f4b53cSBarry Smith   Collective
126635b38c60SMatthew G. Knepley 
126735b38c60SMatthew G. Knepley   Input Parameters:
126820f4b53cSBarry Smith + sw - The `DMSWARM` object
126935b38c60SMatthew G. Knepley - v0 - The velocity scale for nondimensionalization for each species
127035b38c60SMatthew G. Knepley 
127135b38c60SMatthew G. Knepley   Level: advanced
127235b38c60SMatthew G. Knepley 
127320f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmComputeLocalSize()`, `DMSwarmInitializeCoordinates()`, `DMSwarmInitializeVelocities()`
127435b38c60SMatthew G. Knepley @*/
1275d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmInitializeVelocitiesFromOptions(DM sw, const PetscReal v0[])
1276d71ae5a4SJacob Faibussowitsch {
1277f8662bd6SBarry Smith   PetscProbFn *sampler;
127835b38c60SMatthew G. Knepley   PetscInt     dim;
127935b38c60SMatthew G. Knepley   const char  *prefix;
12806c5a79ebSMatthew G. Knepley   char         funcname[PETSC_MAX_PATH_LEN];
12816c5a79ebSMatthew G. Knepley   PetscBool    flg;
128235b38c60SMatthew G. Knepley 
128335b38c60SMatthew G. Knepley   PetscFunctionBegin;
1284d0609cedSBarry Smith   PetscOptionsBegin(PetscObjectComm((PetscObject)sw), "", "DMSwarm Options", "DMSWARM");
12856c5a79ebSMatthew G. Knepley   PetscCall(PetscOptionsString("-dm_swarm_velocity_function", "Function to determine particle velocities", "DMSwarmSetVelocityFunction", funcname, funcname, sizeof(funcname), &flg));
1286d0609cedSBarry Smith   PetscOptionsEnd();
12876c5a79ebSMatthew G. Knepley   if (flg) {
12888434afd1SBarry Smith     PetscSimplePointFn *velFunc;
12896c5a79ebSMatthew G. Knepley 
12906c5a79ebSMatthew G. Knepley     PetscCall(PetscDLSym(NULL, funcname, (void **)&velFunc));
12916c5a79ebSMatthew G. Knepley     PetscCheck(velFunc, PetscObjectComm((PetscObject)sw), PETSC_ERR_ARG_WRONG, "Could not locate function %s", funcname);
12926c5a79ebSMatthew G. Knepley     PetscCall(DMSwarmSetVelocityFunction(sw, velFunc));
12936c5a79ebSMatthew G. Knepley   }
12949566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(sw, &dim));
12959566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)sw, &prefix));
12969566063dSJacob Faibussowitsch   PetscCall(PetscProbCreateFromOptions(dim, prefix, "-dm_swarm_velocity_density", NULL, NULL, &sampler));
12979566063dSJacob Faibussowitsch   PetscCall(DMSwarmInitializeVelocities(sw, sampler, v0));
12983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
129935b38c60SMatthew G. Knepley }
13009d50c5ebSMatthew G. Knepley 
13019d50c5ebSMatthew G. Knepley // The input vector U is assumed to be from a PetscFE. The Swarm fields are input as auxiliary values.
13029d50c5ebSMatthew G. Knepley PetscErrorCode DMProjectFieldLocal_Swarm(DM dm, PetscReal time, Vec U, PetscPointFn **funcs, InsertMode mode, Vec X)
13039d50c5ebSMatthew G. Knepley {
13049d50c5ebSMatthew G. Knepley   MPI_Comm         comm;
13059d50c5ebSMatthew G. Knepley   DM               dmIn;
13069d50c5ebSMatthew G. Knepley   PetscDS          ds;
13079d50c5ebSMatthew G. Knepley   PetscTabulation *T;
13089d50c5ebSMatthew G. Knepley   DMSwarmCellDM    celldm;
13099d50c5ebSMatthew G. Knepley   PetscScalar     *a, *val, *u, *u_x;
13109d50c5ebSMatthew G. Knepley   PetscFEGeom      fegeom;
13119d50c5ebSMatthew G. Knepley   PetscReal       *xi, *v0, *J, *invJ, detJ = 1.0, v0ref[3] = {-1.0, -1.0, -1.0};
13129d50c5ebSMatthew G. Knepley   PetscInt         dim, dE, Np, n, Nf, Nfc, Nfu, cStart, cEnd, maxC = 0, totbs = 0;
13139d50c5ebSMatthew G. Knepley   const char     **coordFields, **fields;
13149d50c5ebSMatthew G. Knepley   PetscReal      **coordVals, **vals;
13159d50c5ebSMatthew G. Knepley   PetscInt        *cbs, *bs, *uOff, *uOff_x;
13169d50c5ebSMatthew G. Knepley 
13179d50c5ebSMatthew G. Knepley   PetscFunctionBegin;
13189d50c5ebSMatthew G. Knepley   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
13199d50c5ebSMatthew G. Knepley   PetscCall(VecGetDM(U, &dmIn));
13209d50c5ebSMatthew G. Knepley   PetscCall(DMGetDimension(dmIn, &dim));
13219d50c5ebSMatthew G. Knepley   PetscCall(DMGetCoordinateDim(dmIn, &dE));
13229d50c5ebSMatthew G. Knepley   PetscCall(DMGetDS(dmIn, &ds));
13239d50c5ebSMatthew G. Knepley   PetscCall(PetscDSGetNumFields(ds, &Nfu));
13249d50c5ebSMatthew G. Knepley   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
13259d50c5ebSMatthew G. Knepley   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
13269d50c5ebSMatthew G. Knepley   PetscCall(PetscDSGetTabulation(ds, &T));
13279d50c5ebSMatthew G. Knepley   PetscCall(PetscDSGetEvaluationArrays(ds, &u, NULL, &u_x));
13289d50c5ebSMatthew G. Knepley   PetscCall(PetscMalloc3(dim, &v0, dim * dim, &J, dim * dim, &invJ));
13299d50c5ebSMatthew G. Knepley   PetscCall(DMPlexGetHeightStratum(dmIn, 0, &cStart, &cEnd));
13309d50c5ebSMatthew G. Knepley 
13319d50c5ebSMatthew G. Knepley   fegeom.dim      = dim;
13329d50c5ebSMatthew G. Knepley   fegeom.dimEmbed = dE;
13339d50c5ebSMatthew G. Knepley   fegeom.v        = v0;
13349d50c5ebSMatthew G. Knepley   fegeom.xi       = v0ref;
13359d50c5ebSMatthew G. Knepley   fegeom.J        = J;
13369d50c5ebSMatthew G. Knepley   fegeom.invJ     = invJ;
13379d50c5ebSMatthew G. Knepley   fegeom.detJ     = &detJ;
13389d50c5ebSMatthew G. Knepley 
13399d50c5ebSMatthew G. Knepley   PetscCall(DMSwarmGetLocalSize(dm, &Np));
13409d50c5ebSMatthew G. Knepley   PetscCall(VecGetLocalSize(X, &n));
13419d50c5ebSMatthew G. Knepley   PetscCheck(n == Np, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Output vector local size %" PetscInt_FMT " != %" PetscInt_FMT " number of local particles", n, Np);
13429d50c5ebSMatthew G. Knepley   PetscCall(DMSwarmGetCellDMActive(dm, &celldm));
13439d50c5ebSMatthew G. Knepley   PetscCall(DMSwarmCellDMGetCoordinateFields(celldm, &Nfc, &coordFields));
13449d50c5ebSMatthew G. Knepley   PetscCall(DMSwarmCellDMGetFields(celldm, &Nf, &fields));
13459d50c5ebSMatthew G. Knepley 
13469d50c5ebSMatthew G. Knepley   PetscCall(PetscMalloc2(Nfc, &coordVals, Nfc, &cbs));
13479d50c5ebSMatthew G. Knepley   for (PetscInt i = 0; i < Nfc; ++i) PetscCall(DMSwarmGetField(dm, coordFields[i], &cbs[i], NULL, (void **)&coordVals[i]));
13489d50c5ebSMatthew G. Knepley   PetscCall(PetscMalloc2(Nf, &vals, Nfc, &bs));
13499d50c5ebSMatthew G. Knepley   for (PetscInt i = 0; i < Nf; ++i) {
13509d50c5ebSMatthew G. Knepley     PetscCall(DMSwarmGetField(dm, fields[i], &bs[i], NULL, (void **)&vals[i]));
13519d50c5ebSMatthew G. Knepley     totbs += bs[i];
13529d50c5ebSMatthew G. Knepley   }
13539d50c5ebSMatthew G. Knepley 
13549d50c5ebSMatthew G. Knepley   PetscCall(DMSwarmSortGetAccess(dm));
13559d50c5ebSMatthew G. Knepley   for (PetscInt cell = cStart; cell < cEnd; ++cell) {
13569d50c5ebSMatthew G. Knepley     PetscInt *pindices, Npc;
13579d50c5ebSMatthew G. Knepley 
13589d50c5ebSMatthew G. Knepley     PetscCall(DMSwarmSortGetPointsPerCell(dm, cell, &Npc, &pindices));
13599d50c5ebSMatthew G. Knepley     maxC = PetscMax(maxC, Npc);
13609d50c5ebSMatthew G. Knepley     PetscCall(DMSwarmSortRestorePointsPerCell(dm, cell, &Npc, &pindices));
13619d50c5ebSMatthew G. Knepley   }
13629d50c5ebSMatthew G. Knepley   PetscCall(PetscMalloc3(maxC * dim, &xi, maxC * totbs, &val, Nfu, &T));
13639d50c5ebSMatthew G. Knepley   PetscCall(VecGetArray(X, &a));
13649d50c5ebSMatthew G. Knepley   {
13659d50c5ebSMatthew G. Knepley     for (PetscInt cell = cStart; cell < cEnd; ++cell) {
13669d50c5ebSMatthew G. Knepley       PetscInt *pindices, Npc;
13679d50c5ebSMatthew G. Knepley 
13689d50c5ebSMatthew G. Knepley       // TODO: Use DMField instead of assuming affine
13699d50c5ebSMatthew G. Knepley       PetscCall(DMPlexComputeCellGeometryFEM(dmIn, cell, NULL, v0, J, invJ, &detJ));
13709d50c5ebSMatthew G. Knepley       PetscCall(DMSwarmSortGetPointsPerCell(dm, cell, &Npc, &pindices));
13719d50c5ebSMatthew G. Knepley 
13729d50c5ebSMatthew G. Knepley       PetscScalar *closure = NULL;
13739d50c5ebSMatthew G. Knepley       PetscInt     Ncl;
13749d50c5ebSMatthew G. Knepley 
13759d50c5ebSMatthew G. Knepley       // Get fields from input vector and auxiliary fields from swarm
13769d50c5ebSMatthew G. Knepley       for (PetscInt p = 0; p < Npc; ++p) {
13779d50c5ebSMatthew G. Knepley         PetscReal xr[8];
13789d50c5ebSMatthew G. Knepley         PetscInt  off;
13799d50c5ebSMatthew G. Knepley 
13809d50c5ebSMatthew G. Knepley         off = 0;
13819d50c5ebSMatthew G. Knepley         for (PetscInt i = 0; i < Nfc; ++i) {
13829d50c5ebSMatthew G. Knepley           for (PetscInt b = 0; b < cbs[i]; ++b, ++off) xr[off] = coordVals[i][pindices[p] * cbs[i] + b];
13839d50c5ebSMatthew G. Knepley         }
13849d50c5ebSMatthew G. Knepley         PetscCheck(off == dim, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The total block size of coordinates is %" PetscInt_FMT " != %" PetscInt_FMT " the DM coordinate dimension", off, dim);
13859d50c5ebSMatthew G. Knepley         CoordinatesRealToRef(dE, dim, fegeom.xi, fegeom.v, fegeom.invJ, xr, &xi[p * dim]);
13869d50c5ebSMatthew G. Knepley         off = 0;
13879d50c5ebSMatthew G. Knepley         for (PetscInt i = 0; i < Nf; ++i) {
13889d50c5ebSMatthew G. Knepley           for (PetscInt b = 0; b < bs[i]; ++b, ++off) val[p * totbs + off] = vals[i][pindices[p] * bs[i] + b];
13899d50c5ebSMatthew G. Knepley         }
13909d50c5ebSMatthew G. Knepley         PetscCheck(off == totbs, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The total block size of swarm fields is %" PetscInt_FMT " != %" PetscInt_FMT " the computed total block size", off, totbs);
13919d50c5ebSMatthew G. Knepley       }
13929d50c5ebSMatthew G. Knepley       PetscCall(DMPlexVecGetClosure(dmIn, NULL, U, cell, &Ncl, &closure));
13939d50c5ebSMatthew G. Knepley       for (PetscInt field = 0; field < Nfu; ++field) {
13949d50c5ebSMatthew G. Knepley         PetscFE fe;
13959d50c5ebSMatthew G. Knepley 
13969d50c5ebSMatthew G. Knepley         PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe));
13979d50c5ebSMatthew G. Knepley         PetscCall(PetscFECreateTabulation(fe, 1, Npc, xi, 1, &T[field]));
13989d50c5ebSMatthew G. Knepley       }
13999d50c5ebSMatthew G. Knepley       for (PetscInt p = 0; p < Npc; ++p) {
14009d50c5ebSMatthew G. Knepley         // Get fields from input vector
14019d50c5ebSMatthew G. Knepley         PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nfu, 0, p, T, &fegeom, closure, NULL, u, u_x, NULL));
14029d50c5ebSMatthew G. Knepley         (*funcs[0])(dim, 1, 1, uOff, uOff_x, u, NULL, u_x, bs, NULL, &val[p * totbs], NULL, NULL, time, &xi[p * dim], 0, NULL, &a[pindices[p]]);
14039d50c5ebSMatthew G. Knepley       }
14049d50c5ebSMatthew G. Knepley       PetscCall(DMSwarmSortRestorePointsPerCell(dm, cell, &Npc, &pindices));
14059d50c5ebSMatthew G. Knepley       PetscCall(DMPlexVecRestoreClosure(dmIn, NULL, U, cell, &Ncl, &closure));
14069d50c5ebSMatthew G. Knepley       for (PetscInt field = 0; field < Nfu; ++field) PetscCall(PetscTabulationDestroy(&T[field]));
14079d50c5ebSMatthew G. Knepley     }
14089d50c5ebSMatthew G. Knepley   }
14099d50c5ebSMatthew G. Knepley   for (PetscInt i = 0; i < Nfc; ++i) PetscCall(DMSwarmRestoreField(dm, coordFields[i], &cbs[i], NULL, (void **)&coordVals[i]));
14109d50c5ebSMatthew G. Knepley   for (PetscInt i = 0; i < Nf; ++i) PetscCall(DMSwarmRestoreField(dm, fields[i], &bs[i], NULL, (void **)&vals[i]));
14119d50c5ebSMatthew G. Knepley   PetscCall(VecRestoreArray(X, &a));
14129d50c5ebSMatthew G. Knepley   PetscCall(DMSwarmSortRestoreAccess(dm));
14139d50c5ebSMatthew G. Knepley   PetscCall(PetscFree3(xi, val, T));
14149d50c5ebSMatthew G. Knepley   PetscCall(PetscFree3(v0, J, invJ));
14159d50c5ebSMatthew G. Knepley   PetscCall(PetscFree2(coordVals, cbs));
14169d50c5ebSMatthew G. Knepley   PetscCall(PetscFree2(vals, bs));
14179d50c5ebSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
14189d50c5ebSMatthew G. Knepley }
1419