xref: /petsc/src/dm/impls/swarm/swarm.c (revision d52c2f216e65d48d17427a896822b982b0e2da6e)
1*d52c2f21SMatthew G. Knepley #include "petscsys.h"
257795646SDave May #define PETSCDM_DLL
357795646SDave May #include <petsc/private/dmswarmimpl.h> /*I   "petscdmswarm.h"   I*/
4e8f14785SLisandro Dalcin #include <petsc/private/hashsetij.h>
50643ed39SMatthew G. Knepley #include <petsc/private/petscfeimpl.h>
65917a6f0SStefano Zampini #include <petscviewer.h>
75917a6f0SStefano Zampini #include <petscdraw.h>
883c47955SMatthew G. Knepley #include <petscdmplex.h>
94cc7f7b2SMatthew G. Knepley #include <petscblaslapack.h>
10279f676cSBarry Smith #include "../src/dm/impls/swarm/data_bucket.h"
11d0c080abSJoseph Pusztay #include <petscdmlabel.h>
12d0c080abSJoseph Pusztay #include <petscsection.h>
1357795646SDave May 
14f2b2bee7SDave May PetscLogEvent DMSWARM_Migrate, DMSWARM_SetSizes, DMSWARM_AddPoints, DMSWARM_RemovePoints, DMSWARM_Sort;
15ed923d71SDave May PetscLogEvent DMSWARM_DataExchangerTopologySetup, DMSWARM_DataExchangerBegin, DMSWARM_DataExchangerEnd;
16ed923d71SDave May PetscLogEvent DMSWARM_DataExchangerSendCount, DMSWARM_DataExchangerPack;
17ed923d71SDave May 
18ea78f98cSLisandro Dalcin const char *DMSwarmTypeNames[]          = {"basic", "pic", NULL};
19ea78f98cSLisandro Dalcin const char *DMSwarmMigrateTypeNames[]   = {"basic", "dmcellnscatter", "dmcellexact", "user", NULL};
20ea78f98cSLisandro Dalcin const char *DMSwarmCollectTypeNames[]   = {"basic", "boundingbox", "general", "user", NULL};
21ea78f98cSLisandro Dalcin const char *DMSwarmPICLayoutTypeNames[] = {"regular", "gauss", "subdivision", NULL};
22f0cdbbbaSDave May 
23f0cdbbbaSDave May const char DMSwarmField_pid[]       = "DMSwarm_pid";
24f0cdbbbaSDave May const char DMSwarmField_rank[]      = "DMSwarm_rank";
25f0cdbbbaSDave May const char DMSwarmPICField_coor[]   = "DMSwarmPIC_coor";
26e2d107dbSDave May const char DMSwarmPICField_cellid[] = "DMSwarm_cellid";
27f0cdbbbaSDave May 
282e956fe4SStefano Zampini PetscInt SwarmDataFieldId = -1;
292e956fe4SStefano Zampini 
3074d0cae8SMatthew G. Knepley #if defined(PETSC_HAVE_HDF5)
3174d0cae8SMatthew G. Knepley   #include <petscviewerhdf5.h>
3274d0cae8SMatthew G. Knepley 
33*d52c2f21SMatthew G. Knepley static PetscErrorCode DMInitialize_Swarm(DM);
34*d52c2f21SMatthew G. Knepley static PetscErrorCode DMDestroy_Swarm(DM);
35*d52c2f21SMatthew G. Knepley 
36*d52c2f21SMatthew G. Knepley /* Replace dm with the contents of ndm, and then destroy ndm
37*d52c2f21SMatthew G. Knepley    - Share the DM_Swarm structure
38*d52c2f21SMatthew G. Knepley */
39*d52c2f21SMatthew G. Knepley PetscErrorCode DMSwarmReplace_Internal(DM dm, DM *ndm)
40*d52c2f21SMatthew G. Knepley {
41*d52c2f21SMatthew G. Knepley   DM               dmNew = *ndm;
42*d52c2f21SMatthew G. Knepley   const PetscReal *maxCell, *Lstart, *L;
43*d52c2f21SMatthew G. Knepley   PetscInt         dim;
44*d52c2f21SMatthew G. Knepley 
45*d52c2f21SMatthew G. Knepley   PetscFunctionBegin;
46*d52c2f21SMatthew G. Knepley   if (dm == dmNew) {
47*d52c2f21SMatthew G. Knepley     PetscCall(DMDestroy(ndm));
48*d52c2f21SMatthew G. Knepley     PetscFunctionReturn(PETSC_SUCCESS);
49*d52c2f21SMatthew G. Knepley   }
50*d52c2f21SMatthew G. Knepley   dm->setupcalled = dmNew->setupcalled;
51*d52c2f21SMatthew G. Knepley   if (!dm->hdr.name) {
52*d52c2f21SMatthew G. Knepley     const char *name;
53*d52c2f21SMatthew G. Knepley 
54*d52c2f21SMatthew G. Knepley     PetscCall(PetscObjectGetName((PetscObject)*ndm, &name));
55*d52c2f21SMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)dm, name));
56*d52c2f21SMatthew G. Knepley   }
57*d52c2f21SMatthew G. Knepley   PetscCall(DMGetDimension(dmNew, &dim));
58*d52c2f21SMatthew G. Knepley   PetscCall(DMSetDimension(dm, dim));
59*d52c2f21SMatthew G. Knepley   PetscCall(DMGetPeriodicity(dmNew, &maxCell, &Lstart, &L));
60*d52c2f21SMatthew G. Knepley   PetscCall(DMSetPeriodicity(dm, maxCell, Lstart, L));
61*d52c2f21SMatthew G. Knepley   PetscCall(DMDestroy_Swarm(dm));
62*d52c2f21SMatthew G. Knepley   PetscCall(DMInitialize_Swarm(dm));
63*d52c2f21SMatthew G. Knepley   dm->data = dmNew->data;
64*d52c2f21SMatthew G. Knepley   ((DM_Swarm *)dmNew->data)->refct++;
65*d52c2f21SMatthew G. Knepley   PetscCall(DMDestroy(ndm));
66*d52c2f21SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
67*d52c2f21SMatthew G. Knepley }
68*d52c2f21SMatthew G. Knepley 
6966976f2fSJacob Faibussowitsch static PetscErrorCode VecView_Swarm_HDF5_Internal(Vec v, PetscViewer viewer)
70d71ae5a4SJacob Faibussowitsch {
7174d0cae8SMatthew G. Knepley   DM        dm;
7274d0cae8SMatthew G. Knepley   PetscReal seqval;
7374d0cae8SMatthew G. Knepley   PetscInt  seqnum, bs;
74eb0c6899SMatthew G. Knepley   PetscBool isseq, ists;
7574d0cae8SMatthew G. Knepley 
7674d0cae8SMatthew G. Knepley   PetscFunctionBegin;
779566063dSJacob Faibussowitsch   PetscCall(VecGetDM(v, &dm));
789566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(v, &bs));
799566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5PushGroup(viewer, "/particle_fields"));
809566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)v, VECSEQ, &isseq));
819566063dSJacob Faibussowitsch   PetscCall(DMGetOutputSequenceNumber(dm, &seqnum, &seqval));
82eb0c6899SMatthew G. Knepley   PetscCall(PetscViewerHDF5IsTimestepping(viewer, &ists));
83eb0c6899SMatthew G. Knepley   if (ists) PetscCall(PetscViewerHDF5SetTimestep(viewer, seqnum));
849566063dSJacob Faibussowitsch   /* PetscCall(DMSequenceView_HDF5(dm, "time", seqnum, (PetscScalar) seqval, viewer)); */
859566063dSJacob Faibussowitsch   PetscCall(VecViewNative(v, viewer));
869566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject)v, "Nc", PETSC_INT, (void *)&bs));
879566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5PopGroup(viewer));
883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8974d0cae8SMatthew G. Knepley }
9074d0cae8SMatthew G. Knepley 
9166976f2fSJacob Faibussowitsch static PetscErrorCode DMSwarmView_HDF5(DM dm, PetscViewer viewer)
92d71ae5a4SJacob Faibussowitsch {
9374d0cae8SMatthew G. Knepley   Vec         coordinates;
9474d0cae8SMatthew G. Knepley   PetscInt    Np;
9574d0cae8SMatthew G. Knepley   PetscBool   isseq;
96*d52c2f21SMatthew G. Knepley   const char *coordname;
9774d0cae8SMatthew G. Knepley 
9874d0cae8SMatthew G. Knepley   PetscFunctionBegin;
99*d52c2f21SMatthew G. Knepley   PetscCall(DMSwarmGetCoordinateField(dm, &coordname));
1009566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetSize(dm, &Np));
101*d52c2f21SMatthew G. Knepley   PetscCall(DMSwarmCreateGlobalVectorFromField(dm, coordname, &coordinates));
1029566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates"));
1039566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5PushGroup(viewer, "/particles"));
1049566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)coordinates, VECSEQ, &isseq));
1059566063dSJacob Faibussowitsch   PetscCall(VecViewNative(coordinates, viewer));
1069566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject)coordinates, "Np", PETSC_INT, (void *)&Np));
1079566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5PopGroup(viewer));
108*d52c2f21SMatthew G. Knepley   PetscCall(DMSwarmDestroyGlobalVectorFromField(dm, coordname, &coordinates));
1093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11074d0cae8SMatthew G. Knepley }
11174d0cae8SMatthew G. Knepley #endif
11274d0cae8SMatthew G. Knepley 
11366976f2fSJacob Faibussowitsch static PetscErrorCode VecView_Swarm(Vec v, PetscViewer viewer)
114d71ae5a4SJacob Faibussowitsch {
11574d0cae8SMatthew G. Knepley   DM dm;
116f9558f5cSBarry Smith #if defined(PETSC_HAVE_HDF5)
11774d0cae8SMatthew G. Knepley   PetscBool ishdf5;
118f9558f5cSBarry Smith #endif
11974d0cae8SMatthew G. Knepley 
12074d0cae8SMatthew G. Knepley   PetscFunctionBegin;
1219566063dSJacob Faibussowitsch   PetscCall(VecGetDM(v, &dm));
12228b400f6SJacob Faibussowitsch   PetscCheck(dm, PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_WRONG, "Vector not generated from a DM");
123f9558f5cSBarry Smith #if defined(PETSC_HAVE_HDF5)
1249566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
12574d0cae8SMatthew G. Knepley   if (ishdf5) {
1269566063dSJacob Faibussowitsch     PetscCall(VecView_Swarm_HDF5_Internal(v, viewer));
1273ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
12874d0cae8SMatthew G. Knepley   }
129f9558f5cSBarry Smith #endif
1309566063dSJacob Faibussowitsch   PetscCall(VecViewNative(v, viewer));
1313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13274d0cae8SMatthew G. Knepley }
13374d0cae8SMatthew G. Knepley 
134*d52c2f21SMatthew G. Knepley /*@C
135*d52c2f21SMatthew G. Knepley   DMSwarmVectorGetField - Gets the fields from which to define a `Vec` object
1360bf7c1c5SMatthew G. Knepley   when `DMCreateLocalVector()`, or `DMCreateGlobalVector()` is called
1370bf7c1c5SMatthew G. Knepley 
1380bf7c1c5SMatthew G. Knepley   Not collective
1390bf7c1c5SMatthew G. Knepley 
1400bf7c1c5SMatthew G. Knepley   Input Parameter:
1410bf7c1c5SMatthew G. Knepley . dm - a `DMSWARM`
1420bf7c1c5SMatthew G. Knepley 
143*d52c2f21SMatthew G. Knepley   Output Parameters:
144*d52c2f21SMatthew G. Knepley + Nf         - the number of fields
145*d52c2f21SMatthew G. Knepley - fieldnames - the textual name given to each registered field, or NULL if it has not been set
1460bf7c1c5SMatthew G. Knepley 
1470bf7c1c5SMatthew G. Knepley   Level: beginner
1480bf7c1c5SMatthew G. Knepley 
149*d52c2f21SMatthew G. Knepley .seealso: `DM`, `DMSWARM`, `DMSwarmVectorDefineField()`, `DMSwarmRegisterPetscDatatypeField()`, `DMCreateGlobalVector()`, `DMCreateLocalVector()`
1500bf7c1c5SMatthew G. Knepley @*/
151*d52c2f21SMatthew G. Knepley PetscErrorCode DMSwarmVectorGetField(DM dm, PetscInt *Nf, const char **fieldnames[])
1520bf7c1c5SMatthew G. Knepley {
1530bf7c1c5SMatthew G. Knepley   PetscFunctionBegin;
1540bf7c1c5SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
155*d52c2f21SMatthew G. Knepley   if (Nf) {
156*d52c2f21SMatthew G. Knepley     PetscAssertPointer(Nf, 2);
157*d52c2f21SMatthew G. Knepley     *Nf = ((DM_Swarm *)dm->data)->vec_field_num;
158*d52c2f21SMatthew G. Knepley   }
159*d52c2f21SMatthew G. Knepley   if (fieldnames) {
160*d52c2f21SMatthew G. Knepley     PetscAssertPointer(fieldnames, 3);
161*d52c2f21SMatthew G. Knepley     *fieldnames = ((DM_Swarm *)dm->data)->vec_field_names;
162*d52c2f21SMatthew G. Knepley   }
1630bf7c1c5SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
1640bf7c1c5SMatthew G. Knepley }
1650bf7c1c5SMatthew G. Knepley 
166cc4c1da9SBarry Smith /*@
16720f4b53cSBarry Smith   DMSwarmVectorDefineField - Sets the field from which to define a `Vec` object
16820f4b53cSBarry Smith   when `DMCreateLocalVector()`, or `DMCreateGlobalVector()` is called
16957795646SDave May 
17020f4b53cSBarry Smith   Collective
17157795646SDave May 
17260225df5SJacob Faibussowitsch   Input Parameters:
17320f4b53cSBarry Smith + dm        - a `DMSWARM`
174*d52c2f21SMatthew G. Knepley - fieldname - the textual name given to each registered field
17557795646SDave May 
176d3a51819SDave May   Level: beginner
17757795646SDave May 
178d3a51819SDave May   Notes:
17920f4b53cSBarry Smith   The field with name `fieldname` must be defined as having a data type of `PetscScalar`.
180e7af74c8SDave May 
18120f4b53cSBarry Smith   This function must be called prior to calling `DMCreateLocalVector()`, `DMCreateGlobalVector()`.
18220f4b53cSBarry Smith   Multiple calls to `DMSwarmVectorDefineField()` are permitted.
183e7af74c8SDave May 
184*d52c2f21SMatthew G. Knepley .seealso: `DM`, `DMSWARM`, `DMSwarmVectorDefineFields()`, `DMSwarmVectorGetField()`, `DMSwarmRegisterPetscDatatypeField()`, `DMCreateGlobalVector()`, `DMCreateLocalVector()`
185d3a51819SDave May @*/
186d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmVectorDefineField(DM dm, const char fieldname[])
187d71ae5a4SJacob Faibussowitsch {
188*d52c2f21SMatthew G. Knepley   PetscFunctionBegin;
189*d52c2f21SMatthew G. Knepley   PetscCall(DMSwarmVectorDefineFields(dm, 1, &fieldname));
190*d52c2f21SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
191*d52c2f21SMatthew G. Knepley }
192*d52c2f21SMatthew G. Knepley 
193*d52c2f21SMatthew G. Knepley /*@C
194*d52c2f21SMatthew G. Knepley   DMSwarmVectorDefineFields - Sets the fields from which to define a `Vec` object
195*d52c2f21SMatthew G. Knepley   when `DMCreateLocalVector()`, or `DMCreateGlobalVector()` is called
196*d52c2f21SMatthew G. Knepley 
197*d52c2f21SMatthew G. Knepley   Collective, No Fortran support
198*d52c2f21SMatthew G. Knepley 
199*d52c2f21SMatthew G. Knepley   Input Parameters:
200*d52c2f21SMatthew G. Knepley + dm         - a `DMSWARM`
201*d52c2f21SMatthew G. Knepley . Nf         - the number of fields
202*d52c2f21SMatthew G. Knepley - fieldnames - the textual name given to each registered field
203*d52c2f21SMatthew G. Knepley 
204*d52c2f21SMatthew G. Knepley   Level: beginner
205*d52c2f21SMatthew G. Knepley 
206*d52c2f21SMatthew G. Knepley   Notes:
207*d52c2f21SMatthew G. Knepley   Each field with name in `fieldnames` must be defined as having a data type of `PetscScalar`.
208*d52c2f21SMatthew G. Knepley 
209*d52c2f21SMatthew G. Knepley   This function must be called prior to calling `DMCreateLocalVector()`, `DMCreateGlobalVector()`.
210*d52c2f21SMatthew G. Knepley   Multiple calls to `DMSwarmVectorDefineField()` are permitted.
211*d52c2f21SMatthew G. Knepley 
212*d52c2f21SMatthew G. Knepley .seealso: `DM`, `DMSWARM`, `DMSwarmVectorDefineField()`, `DMSwarmVectorGetField()`, `DMSwarmRegisterPetscDatatypeField()`, `DMCreateGlobalVector()`, `DMCreateLocalVector()`
213*d52c2f21SMatthew G. Knepley @*/
214*d52c2f21SMatthew G. Knepley PetscErrorCode DMSwarmVectorDefineFields(DM dm, PetscInt Nf, const char *fieldnames[])
215*d52c2f21SMatthew G. Knepley {
216*d52c2f21SMatthew G. Knepley   DM_Swarm *sw = (DM_Swarm *)dm->data;
217b5bcf523SDave May 
218a9cbaee5SMatthew G. Knepley   PetscFunctionBegin;
2190bf7c1c5SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
220*d52c2f21SMatthew G. Knepley   if (fieldnames) PetscAssertPointer(fieldnames, 3);
221*d52c2f21SMatthew G. Knepley   if (!sw->issetup) PetscCall(DMSetUp(dm));
222*d52c2f21SMatthew G. Knepley   PetscCheck(Nf >= 0, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Number of fields must be non-negative, not %" PetscInt_FMT, Nf);
223*d52c2f21SMatthew G. Knepley   for (PetscInt f = 0; f < sw->vec_field_num; ++f) PetscCall(PetscFree(sw->vec_field_names[f]));
224*d52c2f21SMatthew G. Knepley   PetscCall(PetscFree(sw->vec_field_names));
225b5bcf523SDave May 
226*d52c2f21SMatthew G. Knepley   sw->vec_field_num = Nf;
227*d52c2f21SMatthew G. Knepley   sw->vec_field_bs  = 0;
228*d52c2f21SMatthew G. Knepley   PetscCall(DMSwarmDataBucketGetSizes(sw->db, &sw->vec_field_nlocal, NULL, NULL));
229*d52c2f21SMatthew G. Knepley   PetscCall(PetscMalloc1(Nf, &sw->vec_field_names));
230*d52c2f21SMatthew G. Knepley   for (PetscInt f = 0; f < Nf; ++f) {
231*d52c2f21SMatthew G. Knepley     PetscInt      bs;
232*d52c2f21SMatthew G. Knepley     PetscScalar  *array;
233*d52c2f21SMatthew G. Knepley     PetscDataType type;
234*d52c2f21SMatthew G. Knepley 
235*d52c2f21SMatthew G. Knepley     PetscCall(DMSwarmGetField(dm, fieldnames[f], &bs, &type, (void **)&array));
236*d52c2f21SMatthew G. Knepley     // Check all fields are of type PETSC_REAL or PETSC_SCALAR
23708401ef6SPierre Jolivet     PetscCheck(type == PETSC_REAL, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only valid for PETSC_REAL");
238*d52c2f21SMatthew G. Knepley     sw->vec_field_bs += bs;
239*d52c2f21SMatthew G. Knepley     PetscCall(DMSwarmRestoreField(dm, fieldnames[f], &bs, &type, (void **)&array));
240*d52c2f21SMatthew G. Knepley     PetscCall(PetscStrallocpy(fieldnames[f], (char **)&sw->vec_field_names[f]));
241*d52c2f21SMatthew G. Knepley   }
2423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
243b5bcf523SDave May }
244b5bcf523SDave May 
245cc651181SDave May /* requires DMSwarmDefineFieldVector has been called */
24666976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateGlobalVector_Swarm(DM dm, Vec *vec)
247d71ae5a4SJacob Faibussowitsch {
248b5bcf523SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
249b5bcf523SDave May   Vec       x;
250b5bcf523SDave May   char      name[PETSC_MAX_PATH_LEN];
251b5bcf523SDave May 
252a9cbaee5SMatthew G. Knepley   PetscFunctionBegin;
2539566063dSJacob Faibussowitsch   if (!swarm->issetup) PetscCall(DMSetUp(dm));
254*d52c2f21SMatthew G. Knepley   PetscCheck(swarm->vec_field_num, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Must call DMSwarmVectorDefineField(s) first");
255*d52c2f21SMatthew G. Knepley   PetscCheck(swarm->vec_field_nlocal == swarm->db->L, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "DMSwarm sizes have changed since last call to VectorDefineField(s) first"); /* Stale data */
256cc651181SDave May 
257*d52c2f21SMatthew G. Knepley   PetscCall(PetscStrncpy(name, "DMSwarmField", PETSC_MAX_PATH_LEN));
258*d52c2f21SMatthew G. Knepley   for (PetscInt f = 0; f < swarm->vec_field_num; ++f) {
259*d52c2f21SMatthew G. Knepley     PetscCall(PetscStrlcat(name, "_", PETSC_MAX_PATH_LEN));
260*d52c2f21SMatthew G. Knepley     PetscCall(PetscStrlcat(name, swarm->vec_field_names[f], PETSC_MAX_PATH_LEN));
261*d52c2f21SMatthew G. Knepley   }
2629566063dSJacob Faibussowitsch   PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), &x));
2639566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)x, name));
2649566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(x, swarm->db->L * swarm->vec_field_bs, PETSC_DETERMINE));
2659566063dSJacob Faibussowitsch   PetscCall(VecSetBlockSize(x, swarm->vec_field_bs));
2669566063dSJacob Faibussowitsch   PetscCall(VecSetDM(x, dm));
2679566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(x));
2689566063dSJacob Faibussowitsch   PetscCall(VecSetDM(x, dm));
2699566063dSJacob Faibussowitsch   PetscCall(VecSetOperation(x, VECOP_VIEW, (void (*)(void))VecView_Swarm));
270b5bcf523SDave May   *vec = x;
2713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
272b5bcf523SDave May }
273b5bcf523SDave May 
274b5bcf523SDave May /* requires DMSwarmDefineFieldVector has been called */
27566976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateLocalVector_Swarm(DM dm, Vec *vec)
276d71ae5a4SJacob Faibussowitsch {
277b5bcf523SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
278b5bcf523SDave May   Vec       x;
279b5bcf523SDave May   char      name[PETSC_MAX_PATH_LEN];
280b5bcf523SDave May 
281a9cbaee5SMatthew G. Knepley   PetscFunctionBegin;
2829566063dSJacob Faibussowitsch   if (!swarm->issetup) PetscCall(DMSetUp(dm));
283*d52c2f21SMatthew G. Knepley   PetscCheck(swarm->vec_field_num, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Must call DMSwarmVectorDefineField(s) first");
284*d52c2f21SMatthew G. Knepley   PetscCheck(swarm->vec_field_nlocal == swarm->db->L, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "DMSwarm sizes have changed since last call to VectorDefineField(s) first");
285cc651181SDave May 
286*d52c2f21SMatthew G. Knepley   PetscCall(PetscStrncpy(name, "DMSwarmField", PETSC_MAX_PATH_LEN));
287*d52c2f21SMatthew G. Knepley   for (PetscInt f = 0; f < swarm->vec_field_num; ++f) {
288*d52c2f21SMatthew G. Knepley     PetscCall(PetscStrlcat(name, "_", PETSC_MAX_PATH_LEN));
289*d52c2f21SMatthew G. Knepley     PetscCall(PetscStrlcat(name, swarm->vec_field_names[f], PETSC_MAX_PATH_LEN));
290*d52c2f21SMatthew G. Knepley   }
2919566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_SELF, &x));
2929566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)x, name));
2939566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(x, swarm->db->L * swarm->vec_field_bs, PETSC_DETERMINE));
2949566063dSJacob Faibussowitsch   PetscCall(VecSetBlockSize(x, swarm->vec_field_bs));
2959566063dSJacob Faibussowitsch   PetscCall(VecSetDM(x, dm));
2969566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(x));
297b5bcf523SDave May   *vec = x;
2983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
299b5bcf523SDave May }
300b5bcf523SDave May 
301d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSwarmDestroyVectorFromField_Private(DM dm, const char fieldname[], Vec *vec)
302d71ae5a4SJacob Faibussowitsch {
303fb1bcc12SMatthew G. Knepley   DM_Swarm        *swarm = (DM_Swarm *)dm->data;
30477048351SPatrick Sanan   DMSwarmDataField gfield;
3052e956fe4SStefano Zampini   PetscInt         bs, nlocal, fid = -1, cfid = -2;
3062e956fe4SStefano Zampini   PetscBool        flg;
307d3a51819SDave May 
308fb1bcc12SMatthew G. Knepley   PetscFunctionBegin;
3092e956fe4SStefano Zampini   /* check vector is an inplace array */
3102e956fe4SStefano Zampini   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldIdByName(swarm->db, fieldname, &fid));
3112e956fe4SStefano Zampini   PetscCall(PetscObjectComposedDataGetInt((PetscObject)*vec, SwarmDataFieldId, cfid, flg));
312ea17275aSJose E. Roman   (void)flg; /* avoid compiler warning */
313e978a55eSPierre Jolivet   PetscCheck(cfid == fid, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Vector being destroyed was not created from DMSwarm field(%s)! %" PetscInt_FMT " != %" PetscInt_FMT, fieldname, cfid, fid);
3149566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(*vec, &nlocal));
3159566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(*vec, &bs));
31608401ef6SPierre Jolivet   PetscCheck(nlocal / bs == swarm->db->L, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "DMSwarm sizes have changed since vector was created - cannot ensure pointers are valid");
3179566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield));
3189566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataFieldRestoreAccess(gfield));
3198df78e51SMark Adams   PetscCall(VecResetArray(*vec));
3209566063dSJacob Faibussowitsch   PetscCall(VecDestroy(vec));
3213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
322fb1bcc12SMatthew G. Knepley }
323fb1bcc12SMatthew G. Knepley 
324d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSwarmCreateVectorFromField_Private(DM dm, const char fieldname[], MPI_Comm comm, Vec *vec)
325d71ae5a4SJacob Faibussowitsch {
326fb1bcc12SMatthew G. Knepley   DM_Swarm     *swarm = (DM_Swarm *)dm->data;
327fb1bcc12SMatthew G. Knepley   PetscDataType type;
328fb1bcc12SMatthew G. Knepley   PetscScalar  *array;
3292e956fe4SStefano Zampini   PetscInt      bs, n, fid;
330fb1bcc12SMatthew G. Knepley   char          name[PETSC_MAX_PATH_LEN];
331e4fbd051SBarry Smith   PetscMPIInt   size;
3328df78e51SMark Adams   PetscBool     iscuda, iskokkos;
333fb1bcc12SMatthew G. Knepley 
334fb1bcc12SMatthew G. Knepley   PetscFunctionBegin;
3359566063dSJacob Faibussowitsch   if (!swarm->issetup) PetscCall(DMSetUp(dm));
3369566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &n, NULL, NULL));
3379566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, fieldname, &bs, &type, (void **)&array));
33808401ef6SPierre Jolivet   PetscCheck(type == PETSC_REAL, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only valid for PETSC_REAL");
339fb1bcc12SMatthew G. Knepley 
3409566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
3418df78e51SMark Adams   PetscCall(PetscStrcmp(dm->vectype, VECKOKKOS, &iskokkos));
3428df78e51SMark Adams   PetscCall(PetscStrcmp(dm->vectype, VECCUDA, &iscuda));
3438df78e51SMark Adams   PetscCall(VecCreate(comm, vec));
3448df78e51SMark Adams   PetscCall(VecSetSizes(*vec, n * bs, PETSC_DETERMINE));
3458df78e51SMark Adams   PetscCall(VecSetBlockSize(*vec, bs));
3468df78e51SMark Adams   if (iskokkos) PetscCall(VecSetType(*vec, VECKOKKOS));
3478df78e51SMark Adams   else if (iscuda) PetscCall(VecSetType(*vec, VECCUDA));
3488df78e51SMark Adams   else PetscCall(VecSetType(*vec, VECSTANDARD));
3498df78e51SMark Adams   PetscCall(VecPlaceArray(*vec, array));
3508df78e51SMark Adams 
3519566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN - 1, "DMSwarmSharedField_%s", fieldname));
3529566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)*vec, name));
353fb1bcc12SMatthew G. Knepley 
354fb1bcc12SMatthew G. Knepley   /* Set guard */
3552e956fe4SStefano Zampini   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldIdByName(swarm->db, fieldname, &fid));
3562e956fe4SStefano Zampini   PetscCall(PetscObjectComposedDataSetInt((PetscObject)*vec, SwarmDataFieldId, fid));
35774d0cae8SMatthew G. Knepley 
3589566063dSJacob Faibussowitsch   PetscCall(VecSetDM(*vec, dm));
3599566063dSJacob Faibussowitsch   PetscCall(VecSetOperation(*vec, VECOP_VIEW, (void (*)(void))VecView_Swarm));
3603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
361fb1bcc12SMatthew G. Knepley }
362fb1bcc12SMatthew G. Knepley 
3630643ed39SMatthew G. Knepley /* This creates a "mass matrix" between a finite element and particle space. If a finite element interpolant is given by
3640643ed39SMatthew G. Knepley 
3650643ed39SMatthew G. Knepley      \hat f = \sum_i f_i \phi_i
3660643ed39SMatthew G. Knepley 
3670643ed39SMatthew G. Knepley    and a particle function is given by
3680643ed39SMatthew G. Knepley 
3690643ed39SMatthew G. Knepley      f = \sum_i w_i \delta(x - x_i)
3700643ed39SMatthew G. Knepley 
3710643ed39SMatthew G. Knepley    then we want to require that
3720643ed39SMatthew G. Knepley 
3730643ed39SMatthew G. Knepley      M \hat f = M_p f
3740643ed39SMatthew G. Knepley 
3750643ed39SMatthew G. Knepley    where the particle mass matrix is given by
3760643ed39SMatthew G. Knepley 
3770643ed39SMatthew G. Knepley      (M_p)_{ij} = \int \phi_i \delta(x - x_j)
3780643ed39SMatthew G. Knepley 
3790643ed39SMatthew G. Knepley    The way Dave May does particles, they amount to quadratue weights rather than delta functions, so he has |J| is in
3800643ed39SMatthew G. Knepley    his integral. We allow this with the boolean flag.
3810643ed39SMatthew G. Knepley */
382d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSwarmComputeMassMatrix_Private(DM dmc, DM dmf, Mat mass, PetscBool useDeltaFunction, void *ctx)
383d71ae5a4SJacob Faibussowitsch {
38483c47955SMatthew G. Knepley   const char  *name = "Mass Matrix";
3850643ed39SMatthew G. Knepley   MPI_Comm     comm;
38683c47955SMatthew G. Knepley   PetscDS      prob;
38783c47955SMatthew G. Knepley   PetscSection fsection, globalFSection;
388e8f14785SLisandro Dalcin   PetscHSetIJ  ht;
3890643ed39SMatthew G. Knepley   PetscLayout  rLayout, colLayout;
39083c47955SMatthew G. Knepley   PetscInt    *dnz, *onz;
391adb2528bSMark Adams   PetscInt     locRows, locCols, rStart, colStart, colEnd, *rowIDXs;
3920643ed39SMatthew G. Knepley   PetscReal   *xi, *v0, *J, *invJ, detJ = 1.0, v0ref[3] = {-1.0, -1.0, -1.0};
39383c47955SMatthew G. Knepley   PetscScalar *elemMat;
3940bf7c1c5SMatthew G. Knepley   PetscInt     dim, Nf, field, cStart, cEnd, cell, totDim, maxC = 0, totNc = 0;
395*d52c2f21SMatthew G. Knepley   const char  *coordname;
39683c47955SMatthew G. Knepley 
39783c47955SMatthew G. Knepley   PetscFunctionBegin;
3989566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)mass, &comm));
3999566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dmf, &dim));
4009566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dmf, &prob));
4019566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(prob, &Nf));
4029566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(prob, &totDim));
4039566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(dim, &v0, dim * dim, &J, dim * dim, &invJ));
4049566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dmf, &fsection));
4059566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalSection(dmf, &globalFSection));
4069566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd));
4079566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(mass, &locRows, &locCols));
40883c47955SMatthew G. Knepley 
409*d52c2f21SMatthew G. Knepley   PetscCall(DMSwarmGetCoordinateField(dmc, &coordname));
410*d52c2f21SMatthew G. Knepley 
4119566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreate(comm, &colLayout));
4129566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetLocalSize(colLayout, locCols));
4139566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(colLayout, 1));
4149566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(colLayout));
4159566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetRange(colLayout, &colStart, &colEnd));
4169566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&colLayout));
4170643ed39SMatthew G. Knepley 
4189566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreate(comm, &rLayout));
4199566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetLocalSize(rLayout, locRows));
4209566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(rLayout, 1));
4219566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(rLayout));
4229566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetRange(rLayout, &rStart, NULL));
4239566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&rLayout));
4240643ed39SMatthew G. Knepley 
4259566063dSJacob Faibussowitsch   PetscCall(PetscCalloc2(locRows, &dnz, locRows, &onz));
4269566063dSJacob Faibussowitsch   PetscCall(PetscHSetIJCreate(&ht));
42753e60ab4SJoseph Pusztay 
4289566063dSJacob Faibussowitsch   PetscCall(PetscSynchronizedFlush(comm, NULL));
4290bf7c1c5SMatthew G. Knepley   for (field = 0; field < Nf; ++field) {
4300bf7c1c5SMatthew G. Knepley     PetscObject  obj;
4310bf7c1c5SMatthew G. Knepley     PetscClassId id;
4320bf7c1c5SMatthew G. Knepley     PetscInt     Nc;
4330bf7c1c5SMatthew G. Knepley 
4340bf7c1c5SMatthew G. Knepley     PetscCall(PetscDSGetDiscretization(prob, field, &obj));
4350bf7c1c5SMatthew G. Knepley     PetscCall(PetscObjectGetClassId(obj, &id));
4360bf7c1c5SMatthew G. Knepley     if (id == PETSCFE_CLASSID) PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc));
4370bf7c1c5SMatthew G. Knepley     else PetscCall(PetscFVGetNumComponents((PetscFV)obj, &Nc));
4380bf7c1c5SMatthew G. Knepley     totNc += Nc;
4390bf7c1c5SMatthew G. Knepley   }
4400643ed39SMatthew G. Knepley   /* count non-zeros */
4419566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortGetAccess(dmc));
44283c47955SMatthew G. Knepley   for (field = 0; field < Nf; ++field) {
4430bf7c1c5SMatthew G. Knepley     PetscObject  obj;
4440bf7c1c5SMatthew G. Knepley     PetscClassId id;
4450bf7c1c5SMatthew G. Knepley     PetscInt     Nc;
4460bf7c1c5SMatthew G. Knepley 
4470bf7c1c5SMatthew G. Knepley     PetscCall(PetscDSGetDiscretization(prob, field, &obj));
4480bf7c1c5SMatthew G. Knepley     PetscCall(PetscObjectGetClassId(obj, &id));
4490bf7c1c5SMatthew G. Knepley     if (id == PETSCFE_CLASSID) PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc));
4500bf7c1c5SMatthew G. Knepley     else PetscCall(PetscFVGetNumComponents((PetscFV)obj, &Nc));
4510bf7c1c5SMatthew G. Knepley 
45283c47955SMatthew G. Knepley     for (cell = cStart; cell < cEnd; ++cell) {
4530643ed39SMatthew G. Knepley       PetscInt *findices, *cindices; /* fine is vertices, coarse is particles */
45483c47955SMatthew G. Knepley       PetscInt  numFIndices, numCIndices;
45583c47955SMatthew G. Knepley 
4569566063dSJacob Faibussowitsch       PetscCall(DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
4579566063dSJacob Faibussowitsch       PetscCall(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices));
458fc7c92abSMatthew G. Knepley       maxC = PetscMax(maxC, numCIndices);
45983c47955SMatthew G. Knepley       {
460e8f14785SLisandro Dalcin         PetscHashIJKey key;
461e8f14785SLisandro Dalcin         PetscBool      missing;
4620bf7c1c5SMatthew G. Knepley         for (PetscInt i = 0; i < numFIndices; ++i) {
463adb2528bSMark Adams           key.j = findices[i]; /* global column (from Plex) */
464adb2528bSMark Adams           if (key.j >= 0) {
46583c47955SMatthew G. Knepley             /* Get indices for coarse elements */
4660bf7c1c5SMatthew G. Knepley             for (PetscInt j = 0; j < numCIndices; ++j) {
4670bf7c1c5SMatthew G. Knepley               for (PetscInt c = 0; c < Nc; ++c) {
4680bf7c1c5SMatthew G. Knepley                 // TODO Need field offset on particle here
4690bf7c1c5SMatthew G. Knepley                 key.i = cindices[j] * totNc + c + rStart; /* global cols (from Swarm) */
470adb2528bSMark Adams                 if (key.i < 0) continue;
4719566063dSJacob Faibussowitsch                 PetscCall(PetscHSetIJQueryAdd(ht, key, &missing));
47283c47955SMatthew G. Knepley                 if (missing) {
4730643ed39SMatthew G. Knepley                   if ((key.j >= colStart) && (key.j < colEnd)) ++dnz[key.i - rStart];
474e8f14785SLisandro Dalcin                   else ++onz[key.i - rStart];
47563a3b9bcSJacob Faibussowitsch                 } else SETERRQ(PetscObjectComm((PetscObject)dmf), PETSC_ERR_SUP, "Set new value at %" PetscInt_FMT ",%" PetscInt_FMT, key.i, key.j);
47683c47955SMatthew G. Knepley               }
477fc7c92abSMatthew G. Knepley             }
478fc7c92abSMatthew G. Knepley           }
4790bf7c1c5SMatthew G. Knepley         }
480fe11354eSMatthew G. Knepley         PetscCall(DMSwarmSortRestorePointsPerCell(dmc, cell, &numCIndices, &cindices));
48183c47955SMatthew G. Knepley       }
4829566063dSJacob Faibussowitsch       PetscCall(DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
48383c47955SMatthew G. Knepley     }
48483c47955SMatthew G. Knepley   }
4859566063dSJacob Faibussowitsch   PetscCall(PetscHSetIJDestroy(&ht));
4869566063dSJacob Faibussowitsch   PetscCall(MatXAIJSetPreallocation(mass, 1, dnz, onz, NULL, NULL));
4879566063dSJacob Faibussowitsch   PetscCall(MatSetOption(mass, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
4889566063dSJacob Faibussowitsch   PetscCall(PetscFree2(dnz, onz));
4890bf7c1c5SMatthew G. Knepley   PetscCall(PetscMalloc3(maxC * totNc * totDim, &elemMat, maxC * totNc, &rowIDXs, maxC * dim, &xi));
49083c47955SMatthew G. Knepley   for (field = 0; field < Nf; ++field) {
491ef0bb6c7SMatthew G. Knepley     PetscTabulation Tcoarse;
49283c47955SMatthew G. Knepley     PetscObject     obj;
493ad9634fcSMatthew G. Knepley     PetscClassId    id;
494ea7032a0SMatthew G. Knepley     PetscReal      *fieldVals;
495*d52c2f21SMatthew G. Knepley     PetscInt        Nc, bs;
49683c47955SMatthew G. Knepley 
4979566063dSJacob Faibussowitsch     PetscCall(PetscDSGetDiscretization(prob, field, &obj));
498ad9634fcSMatthew G. Knepley     PetscCall(PetscObjectGetClassId(obj, &id));
499ad9634fcSMatthew G. Knepley     if (id == PETSCFE_CLASSID) PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc));
500ad9634fcSMatthew G. Knepley     else PetscCall(PetscFVGetNumComponents((PetscFV)obj, &Nc));
501ea7032a0SMatthew G. Knepley 
502*d52c2f21SMatthew G. Knepley     PetscCall(DMSwarmGetField(dmc, coordname, &bs, NULL, (void **)&fieldVals));
50383c47955SMatthew G. Knepley     for (cell = cStart; cell < cEnd; ++cell) {
50483c47955SMatthew G. Knepley       PetscInt *findices, *cindices;
50583c47955SMatthew G. Knepley       PetscInt  numFIndices, numCIndices;
50683c47955SMatthew G. Knepley 
5070643ed39SMatthew G. Knepley       /* TODO: Use DMField instead of assuming affine */
5089566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ));
5099566063dSJacob Faibussowitsch       PetscCall(DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
5109566063dSJacob Faibussowitsch       PetscCall(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices));
511*d52c2f21SMatthew G. Knepley       for (PetscInt j = 0; j < numCIndices; ++j) CoordinatesRealToRef(dim, dim, v0ref, v0, invJ, &fieldVals[cindices[j] * bs], &xi[j * dim]);
512ad9634fcSMatthew G. Knepley       if (id == PETSCFE_CLASSID) PetscCall(PetscFECreateTabulation((PetscFE)obj, 1, numCIndices, xi, 0, &Tcoarse));
513ad9634fcSMatthew G. Knepley       else PetscCall(PetscFVCreateTabulation((PetscFV)obj, 1, numCIndices, xi, 0, &Tcoarse));
51483c47955SMatthew G. Knepley       /* Get elemMat entries by multiplying by weight */
5150bf7c1c5SMatthew G. Knepley       PetscCall(PetscArrayzero(elemMat, numCIndices * Nc * totDim));
5160bf7c1c5SMatthew G. Knepley       for (PetscInt i = 0; i < numFIndices / Nc; ++i) {
5170bf7c1c5SMatthew G. Knepley         for (PetscInt j = 0; j < numCIndices; ++j) {
5180bf7c1c5SMatthew G. Knepley           for (PetscInt c = 0; c < Nc; ++c) {
5190bf7c1c5SMatthew G. Knepley             // TODO Need field offset on particle and field here
5200643ed39SMatthew G. Knepley             /* B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c */
5210bf7c1c5SMatthew G. Knepley             elemMat[(j * totNc + c) * numFIndices + i * Nc + c] += Tcoarse->T[0][(j * numFIndices + i * Nc + c) * Nc + c] * (useDeltaFunction ? 1.0 : detJ);
52283c47955SMatthew G. Knepley           }
5230643ed39SMatthew G. Knepley         }
5240643ed39SMatthew G. Knepley       }
5250bf7c1c5SMatthew G. Knepley       for (PetscInt j = 0; j < numCIndices; ++j)
5260bf7c1c5SMatthew G. Knepley         // TODO Need field offset on particle here
5270bf7c1c5SMatthew G. Knepley         for (PetscInt c = 0; c < Nc; ++c) rowIDXs[j * Nc + c] = cindices[j] * totNc + c + rStart;
5280bf7c1c5SMatthew G. Knepley       if (0) PetscCall(DMPrintCellMatrix(cell, name, numCIndices * Nc, numFIndices, elemMat));
5290bf7c1c5SMatthew G. Knepley       PetscCall(MatSetValues(mass, numCIndices * Nc, rowIDXs, numFIndices, findices, elemMat, ADD_VALUES));
530fe11354eSMatthew G. Knepley       PetscCall(DMSwarmSortRestorePointsPerCell(dmc, cell, &numCIndices, &cindices));
5319566063dSJacob Faibussowitsch       PetscCall(DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
5329566063dSJacob Faibussowitsch       PetscCall(PetscTabulationDestroy(&Tcoarse));
53383c47955SMatthew G. Knepley     }
534*d52c2f21SMatthew G. Knepley     PetscCall(DMSwarmRestoreField(dmc, coordname, &bs, NULL, (void **)&fieldVals));
53583c47955SMatthew G. Knepley   }
5369566063dSJacob Faibussowitsch   PetscCall(PetscFree3(elemMat, rowIDXs, xi));
5379566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortRestoreAccess(dmc));
5389566063dSJacob Faibussowitsch   PetscCall(PetscFree3(v0, J, invJ));
5399566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mass, MAT_FINAL_ASSEMBLY));
5409566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mass, MAT_FINAL_ASSEMBLY));
5413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
54283c47955SMatthew G. Knepley }
54383c47955SMatthew G. Knepley 
544d0c080abSJoseph Pusztay /* Returns empty matrix for use with SNES FD */
545d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCreateMatrix_Swarm(DM sw, Mat *m)
546d71ae5a4SJacob Faibussowitsch {
547d0c080abSJoseph Pusztay   Vec      field;
548d0c080abSJoseph Pusztay   PetscInt size;
549d0c080abSJoseph Pusztay 
550d0c080abSJoseph Pusztay   PetscFunctionBegin;
5519566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(sw, &field));
5529566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(field, &size));
5539566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(sw, &field));
5549566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, m));
5559566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(*m));
5569566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*m, PETSC_DECIDE, PETSC_DECIDE, size, size));
5579566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(*m, 1, NULL));
5589566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(*m));
5599566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(*m, MAT_FINAL_ASSEMBLY));
5609566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(*m, MAT_FINAL_ASSEMBLY));
5619566063dSJacob Faibussowitsch   PetscCall(MatShift(*m, 1.0));
5629566063dSJacob Faibussowitsch   PetscCall(MatSetDM(*m, sw));
5633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
564d0c080abSJoseph Pusztay }
565d0c080abSJoseph Pusztay 
566adb2528bSMark Adams /* FEM cols, Particle rows */
567d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCreateMassMatrix_Swarm(DM dmCoarse, DM dmFine, Mat *mass)
568d71ae5a4SJacob Faibussowitsch {
569*d52c2f21SMatthew G. Knepley   DM_Swarm    *swarm = (DM_Swarm *)dmCoarse->data;
570895a1698SMatthew G. Knepley   PetscSection gsf;
571*d52c2f21SMatthew G. Knepley   PetscInt     m, n, Np;
57283c47955SMatthew G. Knepley   void        *ctx;
57383c47955SMatthew G. Knepley 
57483c47955SMatthew G. Knepley   PetscFunctionBegin;
575*d52c2f21SMatthew G. Knepley   PetscCheck(swarm->vec_field_num, PetscObjectComm((PetscObject)dmCoarse), PETSC_ERR_USER, "Must call DMSwarmVectorDefineField(s) first");
5769566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalSection(dmFine, &gsf));
5779566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetConstrainedStorageSize(gsf, &m));
5780bf7c1c5SMatthew G. Knepley   PetscCall(DMSwarmGetLocalSize(dmCoarse, &Np));
579*d52c2f21SMatthew G. Knepley   n = Np * swarm->vec_field_bs;
5809566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)dmCoarse), mass));
5819566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*mass, n, m, PETSC_DETERMINE, PETSC_DETERMINE));
5829566063dSJacob Faibussowitsch   PetscCall(MatSetType(*mass, dmCoarse->mattype));
5839566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(dmFine, &ctx));
58483c47955SMatthew G. Knepley 
5859566063dSJacob Faibussowitsch   PetscCall(DMSwarmComputeMassMatrix_Private(dmCoarse, dmFine, *mass, PETSC_TRUE, ctx));
5869566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(*mass, NULL, "-mass_mat_view"));
5873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
58883c47955SMatthew G. Knepley }
58983c47955SMatthew G. Knepley 
590d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSwarmComputeMassMatrixSquare_Private(DM dmc, DM dmf, Mat mass, PetscBool useDeltaFunction, void *ctx)
591d71ae5a4SJacob Faibussowitsch {
5924cc7f7b2SMatthew G. Knepley   const char  *name = "Mass Matrix Square";
5934cc7f7b2SMatthew G. Knepley   MPI_Comm     comm;
5944cc7f7b2SMatthew G. Knepley   PetscDS      prob;
5954cc7f7b2SMatthew G. Knepley   PetscSection fsection, globalFSection;
5964cc7f7b2SMatthew G. Knepley   PetscHSetIJ  ht;
5974cc7f7b2SMatthew G. Knepley   PetscLayout  rLayout, colLayout;
5984cc7f7b2SMatthew G. Knepley   PetscInt    *dnz, *onz, *adj, depth, maxConeSize, maxSupportSize, maxAdjSize;
5994cc7f7b2SMatthew G. Knepley   PetscInt     locRows, locCols, rStart, colStart, colEnd, *rowIDXs;
6004cc7f7b2SMatthew G. Knepley   PetscReal   *xi, *v0, *J, *invJ, detJ = 1.0, v0ref[3] = {-1.0, -1.0, -1.0};
6014cc7f7b2SMatthew G. Knepley   PetscScalar *elemMat, *elemMatSq;
6024cc7f7b2SMatthew G. Knepley   PetscInt     cdim, Nf, field, cStart, cEnd, cell, totDim, maxC = 0;
603*d52c2f21SMatthew G. Knepley   const char  *coordname;
6044cc7f7b2SMatthew G. Knepley 
6054cc7f7b2SMatthew G. Knepley   PetscFunctionBegin;
6069566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)mass, &comm));
6079566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dmf, &cdim));
6089566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dmf, &prob));
6099566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(prob, &Nf));
6109566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(prob, &totDim));
6119566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(cdim, &v0, cdim * cdim, &J, cdim * cdim, &invJ));
6129566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dmf, &fsection));
6139566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalSection(dmf, &globalFSection));
6149566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd));
6159566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(mass, &locRows, &locCols));
6164cc7f7b2SMatthew G. Knepley 
617*d52c2f21SMatthew G. Knepley   PetscCall(DMSwarmGetCoordinateField(dmc, &coordname));
618*d52c2f21SMatthew G. Knepley 
6199566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreate(comm, &colLayout));
6209566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetLocalSize(colLayout, locCols));
6219566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(colLayout, 1));
6229566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(colLayout));
6239566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetRange(colLayout, &colStart, &colEnd));
6249566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&colLayout));
6254cc7f7b2SMatthew G. Knepley 
6269566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreate(comm, &rLayout));
6279566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetLocalSize(rLayout, locRows));
6289566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(rLayout, 1));
6299566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(rLayout));
6309566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetRange(rLayout, &rStart, NULL));
6319566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&rLayout));
6324cc7f7b2SMatthew G. Knepley 
6339566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(dmf, &depth));
6349566063dSJacob Faibussowitsch   PetscCall(DMPlexGetMaxSizes(dmf, &maxConeSize, &maxSupportSize));
6354cc7f7b2SMatthew G. Knepley   maxAdjSize = PetscPowInt(maxConeSize * maxSupportSize, depth);
6369566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(maxAdjSize, &adj));
6374cc7f7b2SMatthew G. Knepley 
6389566063dSJacob Faibussowitsch   PetscCall(PetscCalloc2(locRows, &dnz, locRows, &onz));
6399566063dSJacob Faibussowitsch   PetscCall(PetscHSetIJCreate(&ht));
6404cc7f7b2SMatthew G. Knepley   /* Count nonzeros
6414cc7f7b2SMatthew G. Knepley        This is just FVM++, but we cannot use the Plex P0 allocation since unknowns in a cell will not be contiguous
6424cc7f7b2SMatthew G. Knepley   */
6439566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortGetAccess(dmc));
6444cc7f7b2SMatthew G. Knepley   for (cell = cStart; cell < cEnd; ++cell) {
6454cc7f7b2SMatthew G. Knepley     PetscInt  i;
6464cc7f7b2SMatthew G. Knepley     PetscInt *cindices;
6474cc7f7b2SMatthew G. Knepley     PetscInt  numCIndices;
6484cc7f7b2SMatthew G. Knepley #if 0
6494cc7f7b2SMatthew G. Knepley     PetscInt  adjSize = maxAdjSize, a, j;
6504cc7f7b2SMatthew G. Knepley #endif
6514cc7f7b2SMatthew G. Knepley 
6529566063dSJacob Faibussowitsch     PetscCall(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices));
6534cc7f7b2SMatthew G. Knepley     maxC = PetscMax(maxC, numCIndices);
6544cc7f7b2SMatthew G. Knepley     /* Diagonal block */
655ad540459SPierre Jolivet     for (i = 0; i < numCIndices; ++i) dnz[cindices[i]] += numCIndices;
6564cc7f7b2SMatthew G. Knepley #if 0
6574cc7f7b2SMatthew G. Knepley     /* Off-diagonal blocks */
6589566063dSJacob Faibussowitsch     PetscCall(DMPlexGetAdjacency(dmf, cell, &adjSize, &adj));
6594cc7f7b2SMatthew G. Knepley     for (a = 0; a < adjSize; ++a) {
6604cc7f7b2SMatthew G. Knepley       if (adj[a] >= cStart && adj[a] < cEnd && adj[a] != cell) {
6614cc7f7b2SMatthew G. Knepley         const PetscInt ncell = adj[a];
6624cc7f7b2SMatthew G. Knepley         PetscInt      *ncindices;
6634cc7f7b2SMatthew G. Knepley         PetscInt       numNCIndices;
6644cc7f7b2SMatthew G. Knepley 
6659566063dSJacob Faibussowitsch         PetscCall(DMSwarmSortGetPointsPerCell(dmc, ncell, &numNCIndices, &ncindices));
6664cc7f7b2SMatthew G. Knepley         {
6674cc7f7b2SMatthew G. Knepley           PetscHashIJKey key;
6684cc7f7b2SMatthew G. Knepley           PetscBool      missing;
6694cc7f7b2SMatthew G. Knepley 
6704cc7f7b2SMatthew G. Knepley           for (i = 0; i < numCIndices; ++i) {
6714cc7f7b2SMatthew G. Knepley             key.i = cindices[i] + rStart; /* global rows (from Swarm) */
6724cc7f7b2SMatthew G. Knepley             if (key.i < 0) continue;
6734cc7f7b2SMatthew G. Knepley             for (j = 0; j < numNCIndices; ++j) {
6744cc7f7b2SMatthew G. Knepley               key.j = ncindices[j] + rStart; /* global column (from Swarm) */
6754cc7f7b2SMatthew G. Knepley               if (key.j < 0) continue;
6769566063dSJacob Faibussowitsch               PetscCall(PetscHSetIJQueryAdd(ht, key, &missing));
6774cc7f7b2SMatthew G. Knepley               if (missing) {
6784cc7f7b2SMatthew G. Knepley                 if ((key.j >= colStart) && (key.j < colEnd)) ++dnz[key.i - rStart];
6794cc7f7b2SMatthew G. Knepley                 else                                         ++onz[key.i - rStart];
6804cc7f7b2SMatthew G. Knepley               }
6814cc7f7b2SMatthew G. Knepley             }
6824cc7f7b2SMatthew G. Knepley           }
6834cc7f7b2SMatthew G. Knepley         }
684fe11354eSMatthew G. Knepley         PetscCall(DMSwarmSortRestorePointsPerCell(dmc, ncell, &numNCIndices, &ncindices));
6854cc7f7b2SMatthew G. Knepley       }
6864cc7f7b2SMatthew G. Knepley     }
6874cc7f7b2SMatthew G. Knepley #endif
688fe11354eSMatthew G. Knepley     PetscCall(DMSwarmSortRestorePointsPerCell(dmc, cell, &numCIndices, &cindices));
6894cc7f7b2SMatthew G. Knepley   }
6909566063dSJacob Faibussowitsch   PetscCall(PetscHSetIJDestroy(&ht));
6919566063dSJacob Faibussowitsch   PetscCall(MatXAIJSetPreallocation(mass, 1, dnz, onz, NULL, NULL));
6929566063dSJacob Faibussowitsch   PetscCall(MatSetOption(mass, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
6939566063dSJacob Faibussowitsch   PetscCall(PetscFree2(dnz, onz));
6949566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(maxC * totDim, &elemMat, maxC * maxC, &elemMatSq, maxC, &rowIDXs, maxC * cdim, &xi));
6954cc7f7b2SMatthew G. Knepley   /* Fill in values
6964cc7f7b2SMatthew G. Knepley        Each entry is a sum of terms \phi_i(x_p) \phi_i(x_q)
6974cc7f7b2SMatthew G. Knepley        Start just by producing block diagonal
6984cc7f7b2SMatthew G. Knepley        Could loop over adjacent cells
6994cc7f7b2SMatthew G. Knepley          Produce neighboring element matrix
7004cc7f7b2SMatthew G. Knepley          TODO Determine which columns and rows correspond to shared dual vector
7014cc7f7b2SMatthew G. Knepley          Do MatMatMult with rectangular matrices
7024cc7f7b2SMatthew G. Knepley          Insert block
7034cc7f7b2SMatthew G. Knepley   */
7044cc7f7b2SMatthew G. Knepley   for (field = 0; field < Nf; ++field) {
7054cc7f7b2SMatthew G. Knepley     PetscTabulation Tcoarse;
7064cc7f7b2SMatthew G. Knepley     PetscObject     obj;
7074cc7f7b2SMatthew G. Knepley     PetscReal      *coords;
7084cc7f7b2SMatthew G. Knepley     PetscInt        Nc, i;
7094cc7f7b2SMatthew G. Knepley 
7109566063dSJacob Faibussowitsch     PetscCall(PetscDSGetDiscretization(prob, field, &obj));
7119566063dSJacob Faibussowitsch     PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc));
71263a3b9bcSJacob Faibussowitsch     PetscCheck(Nc == 1, PetscObjectComm((PetscObject)dmf), PETSC_ERR_SUP, "Can only interpolate a scalar field from particles, Nc = %" PetscInt_FMT, Nc);
713*d52c2f21SMatthew G. Knepley     PetscCall(DMSwarmGetField(dmc, coordname, NULL, NULL, (void **)&coords));
7144cc7f7b2SMatthew G. Knepley     for (cell = cStart; cell < cEnd; ++cell) {
7154cc7f7b2SMatthew G. Knepley       PetscInt *findices, *cindices;
7164cc7f7b2SMatthew G. Knepley       PetscInt  numFIndices, numCIndices;
7174cc7f7b2SMatthew G. Knepley       PetscInt  p, c;
7184cc7f7b2SMatthew G. Knepley 
7194cc7f7b2SMatthew G. Knepley       /* TODO: Use DMField instead of assuming affine */
7209566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ));
7219566063dSJacob Faibussowitsch       PetscCall(DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
7229566063dSJacob Faibussowitsch       PetscCall(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices));
723ad540459SPierre Jolivet       for (p = 0; p < numCIndices; ++p) CoordinatesRealToRef(cdim, cdim, v0ref, v0, invJ, &coords[cindices[p] * cdim], &xi[p * cdim]);
7249566063dSJacob Faibussowitsch       PetscCall(PetscFECreateTabulation((PetscFE)obj, 1, numCIndices, xi, 0, &Tcoarse));
7254cc7f7b2SMatthew G. Knepley       /* Get elemMat entries by multiplying by weight */
7269566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(elemMat, numCIndices * totDim));
7274cc7f7b2SMatthew G. Knepley       for (i = 0; i < numFIndices; ++i) {
7284cc7f7b2SMatthew G. Knepley         for (p = 0; p < numCIndices; ++p) {
7294cc7f7b2SMatthew G. Knepley           for (c = 0; c < Nc; ++c) {
7304cc7f7b2SMatthew G. Knepley             /* B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c */
7314cc7f7b2SMatthew G. Knepley             elemMat[p * numFIndices + i] += Tcoarse->T[0][(p * numFIndices + i) * Nc + c] * (useDeltaFunction ? 1.0 : detJ);
7324cc7f7b2SMatthew G. Knepley           }
7334cc7f7b2SMatthew G. Knepley         }
7344cc7f7b2SMatthew G. Knepley       }
7359566063dSJacob Faibussowitsch       PetscCall(PetscTabulationDestroy(&Tcoarse));
7364cc7f7b2SMatthew G. Knepley       for (p = 0; p < numCIndices; ++p) rowIDXs[p] = cindices[p] + rStart;
7379566063dSJacob Faibussowitsch       if (0) PetscCall(DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat));
7384cc7f7b2SMatthew G. Knepley       /* Block diagonal */
73978884ca7SMatthew G. Knepley       if (numCIndices) {
7404cc7f7b2SMatthew G. Knepley         PetscBLASInt blasn, blask;
7414cc7f7b2SMatthew G. Knepley         PetscScalar  one = 1.0, zero = 0.0;
7424cc7f7b2SMatthew G. Knepley 
7439566063dSJacob Faibussowitsch         PetscCall(PetscBLASIntCast(numCIndices, &blasn));
7449566063dSJacob Faibussowitsch         PetscCall(PetscBLASIntCast(numFIndices, &blask));
745792fecdfSBarry Smith         PetscCallBLAS("BLASgemm", BLASgemm_("T", "N", &blasn, &blasn, &blask, &one, elemMat, &blask, elemMat, &blask, &zero, elemMatSq, &blasn));
7464cc7f7b2SMatthew G. Knepley       }
7479566063dSJacob Faibussowitsch       PetscCall(MatSetValues(mass, numCIndices, rowIDXs, numCIndices, rowIDXs, elemMatSq, ADD_VALUES));
7484cf0e950SBarry Smith       /* TODO off-diagonal */
749fe11354eSMatthew G. Knepley       PetscCall(DMSwarmSortRestorePointsPerCell(dmc, cell, &numCIndices, &cindices));
7509566063dSJacob Faibussowitsch       PetscCall(DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
7514cc7f7b2SMatthew G. Knepley     }
752*d52c2f21SMatthew G. Knepley     PetscCall(DMSwarmRestoreField(dmc, coordname, NULL, NULL, (void **)&coords));
7534cc7f7b2SMatthew G. Knepley   }
7549566063dSJacob Faibussowitsch   PetscCall(PetscFree4(elemMat, elemMatSq, rowIDXs, xi));
7559566063dSJacob Faibussowitsch   PetscCall(PetscFree(adj));
7569566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortRestoreAccess(dmc));
7579566063dSJacob Faibussowitsch   PetscCall(PetscFree3(v0, J, invJ));
7589566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mass, MAT_FINAL_ASSEMBLY));
7599566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mass, MAT_FINAL_ASSEMBLY));
7603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7614cc7f7b2SMatthew G. Knepley }
7624cc7f7b2SMatthew G. Knepley 
763cc4c1da9SBarry Smith /*@
7644cc7f7b2SMatthew G. Knepley   DMSwarmCreateMassMatrixSquare - Creates the block-diagonal of the square, M^T_p M_p, of the particle mass matrix M_p
7654cc7f7b2SMatthew G. Knepley 
76620f4b53cSBarry Smith   Collective
7674cc7f7b2SMatthew G. Knepley 
76860225df5SJacob Faibussowitsch   Input Parameters:
76920f4b53cSBarry Smith + dmCoarse - a `DMSWARM`
77020f4b53cSBarry Smith - dmFine   - a `DMPLEX`
7714cc7f7b2SMatthew G. Knepley 
77260225df5SJacob Faibussowitsch   Output Parameter:
7734cc7f7b2SMatthew G. Knepley . mass - the square of the particle mass matrix
7744cc7f7b2SMatthew G. Knepley 
7754cc7f7b2SMatthew G. Knepley   Level: advanced
7764cc7f7b2SMatthew G. Knepley 
77720f4b53cSBarry Smith   Note:
7784cc7f7b2SMatthew G. Knepley   We only compute the block diagonal since this provides a good preconditioner and is completely local. It would be possible in the
7794cc7f7b2SMatthew G. Knepley   future to compute the full normal equations.
7804cc7f7b2SMatthew G. Knepley 
78120f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMCreateMassMatrix()`
7824cc7f7b2SMatthew G. Knepley @*/
783d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmCreateMassMatrixSquare(DM dmCoarse, DM dmFine, Mat *mass)
784d71ae5a4SJacob Faibussowitsch {
7854cc7f7b2SMatthew G. Knepley   PetscInt n;
7864cc7f7b2SMatthew G. Knepley   void    *ctx;
7874cc7f7b2SMatthew G. Knepley 
7884cc7f7b2SMatthew G. Knepley   PetscFunctionBegin;
7899566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetLocalSize(dmCoarse, &n));
7909566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)dmCoarse), mass));
7919566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*mass, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
7929566063dSJacob Faibussowitsch   PetscCall(MatSetType(*mass, dmCoarse->mattype));
7939566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(dmFine, &ctx));
7944cc7f7b2SMatthew G. Knepley 
7959566063dSJacob Faibussowitsch   PetscCall(DMSwarmComputeMassMatrixSquare_Private(dmCoarse, dmFine, *mass, PETSC_TRUE, ctx));
7969566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(*mass, NULL, "-mass_sq_mat_view"));
7973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7984cc7f7b2SMatthew G. Knepley }
7994cc7f7b2SMatthew G. Knepley 
800cc4c1da9SBarry Smith /*@
80120f4b53cSBarry Smith   DMSwarmCreateGlobalVectorFromField - Creates a `Vec` object sharing the array associated with a given field
802d3a51819SDave May 
80320f4b53cSBarry Smith   Collective
804d3a51819SDave May 
80560225df5SJacob Faibussowitsch   Input Parameters:
80620f4b53cSBarry Smith + dm        - a `DMSWARM`
80762741f57SDave May - fieldname - the textual name given to a registered field
808d3a51819SDave May 
80960225df5SJacob Faibussowitsch   Output Parameter:
810d3a51819SDave May . vec - the vector
811d3a51819SDave May 
812d3a51819SDave May   Level: beginner
813d3a51819SDave May 
814cc4c1da9SBarry Smith   Note:
81520f4b53cSBarry Smith   The vector must be returned using a matching call to `DMSwarmDestroyGlobalVectorFromField()`.
8168b8a3813SDave May 
81720f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmDestroyGlobalVectorFromField()`
818d3a51819SDave May @*/
819d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmCreateGlobalVectorFromField(DM dm, const char fieldname[], Vec *vec)
820d71ae5a4SJacob Faibussowitsch {
821fb1bcc12SMatthew G. Knepley   MPI_Comm comm = PetscObjectComm((PetscObject)dm);
822b5bcf523SDave May 
823fb1bcc12SMatthew G. Knepley   PetscFunctionBegin;
824ea7032a0SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
8259566063dSJacob Faibussowitsch   PetscCall(DMSwarmCreateVectorFromField_Private(dm, fieldname, comm, vec));
8263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
827b5bcf523SDave May }
828b5bcf523SDave May 
829cc4c1da9SBarry Smith /*@
83020f4b53cSBarry Smith   DMSwarmDestroyGlobalVectorFromField - Destroys the `Vec` object which share the array associated with a given field
831d3a51819SDave May 
83220f4b53cSBarry Smith   Collective
833d3a51819SDave May 
83460225df5SJacob Faibussowitsch   Input Parameters:
83520f4b53cSBarry Smith + dm        - a `DMSWARM`
83662741f57SDave May - fieldname - the textual name given to a registered field
837d3a51819SDave May 
83860225df5SJacob Faibussowitsch   Output Parameter:
839d3a51819SDave May . vec - the vector
840d3a51819SDave May 
841d3a51819SDave May   Level: beginner
842d3a51819SDave May 
84320f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmCreateGlobalVectorFromField()`
844d3a51819SDave May @*/
845d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDestroyGlobalVectorFromField(DM dm, const char fieldname[], Vec *vec)
846d71ae5a4SJacob Faibussowitsch {
847fb1bcc12SMatthew G. Knepley   PetscFunctionBegin;
848ea7032a0SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
8499566063dSJacob Faibussowitsch   PetscCall(DMSwarmDestroyVectorFromField_Private(dm, fieldname, vec));
8503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
851b5bcf523SDave May }
852b5bcf523SDave May 
853cc4c1da9SBarry Smith /*@
85420f4b53cSBarry Smith   DMSwarmCreateLocalVectorFromField - Creates a `Vec` object sharing the array associated with a given field
855fb1bcc12SMatthew G. Knepley 
85620f4b53cSBarry Smith   Collective
857fb1bcc12SMatthew G. Knepley 
85860225df5SJacob Faibussowitsch   Input Parameters:
85920f4b53cSBarry Smith + dm        - a `DMSWARM`
86062741f57SDave May - fieldname - the textual name given to a registered field
861fb1bcc12SMatthew G. Knepley 
86260225df5SJacob Faibussowitsch   Output Parameter:
863fb1bcc12SMatthew G. Knepley . vec - the vector
864fb1bcc12SMatthew G. Knepley 
865fb1bcc12SMatthew G. Knepley   Level: beginner
866fb1bcc12SMatthew G. Knepley 
86720f4b53cSBarry Smith   Note:
8688b8a3813SDave May   The vector must be returned using a matching call to DMSwarmDestroyLocalVectorFromField().
8698b8a3813SDave May 
87020f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmDestroyLocalVectorFromField()`
871fb1bcc12SMatthew G. Knepley @*/
872d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmCreateLocalVectorFromField(DM dm, const char fieldname[], Vec *vec)
873d71ae5a4SJacob Faibussowitsch {
874fb1bcc12SMatthew G. Knepley   MPI_Comm comm = PETSC_COMM_SELF;
875bbe8250bSMatthew G. Knepley 
876fb1bcc12SMatthew G. Knepley   PetscFunctionBegin;
8779566063dSJacob Faibussowitsch   PetscCall(DMSwarmCreateVectorFromField_Private(dm, fieldname, comm, vec));
8783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
879bbe8250bSMatthew G. Knepley }
880fb1bcc12SMatthew G. Knepley 
881cc4c1da9SBarry Smith /*@
88220f4b53cSBarry Smith   DMSwarmDestroyLocalVectorFromField - Destroys the `Vec` object which share the array associated with a given field
883fb1bcc12SMatthew G. Knepley 
88420f4b53cSBarry Smith   Collective
885fb1bcc12SMatthew G. Knepley 
88660225df5SJacob Faibussowitsch   Input Parameters:
88720f4b53cSBarry Smith + dm        - a `DMSWARM`
88862741f57SDave May - fieldname - the textual name given to a registered field
889fb1bcc12SMatthew G. Knepley 
89060225df5SJacob Faibussowitsch   Output Parameter:
891fb1bcc12SMatthew G. Knepley . vec - the vector
892fb1bcc12SMatthew G. Knepley 
893fb1bcc12SMatthew G. Knepley   Level: beginner
894fb1bcc12SMatthew G. Knepley 
89520f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmCreateLocalVectorFromField()`
896fb1bcc12SMatthew G. Knepley @*/
897d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDestroyLocalVectorFromField(DM dm, const char fieldname[], Vec *vec)
898d71ae5a4SJacob Faibussowitsch {
899fb1bcc12SMatthew G. Knepley   PetscFunctionBegin;
900ea7032a0SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
9019566063dSJacob Faibussowitsch   PetscCall(DMSwarmDestroyVectorFromField_Private(dm, fieldname, vec));
9023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
903bbe8250bSMatthew G. Knepley }
904bbe8250bSMatthew G. Knepley 
905cc4c1da9SBarry Smith /*@
90620f4b53cSBarry Smith   DMSwarmInitializeFieldRegister - Initiates the registration of fields to a `DMSWARM`
907d3a51819SDave May 
90820f4b53cSBarry Smith   Collective
909d3a51819SDave May 
91060225df5SJacob Faibussowitsch   Input Parameter:
91120f4b53cSBarry Smith . dm - a `DMSWARM`
912d3a51819SDave May 
913d3a51819SDave May   Level: beginner
914d3a51819SDave May 
91520f4b53cSBarry Smith   Note:
91620f4b53cSBarry Smith   After all fields have been registered, you must call `DMSwarmFinalizeFieldRegister()`.
917d3a51819SDave May 
91820f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmFinalizeFieldRegister()`, `DMSwarmRegisterPetscDatatypeField()`,
919db781477SPatrick Sanan           `DMSwarmRegisterUserStructField()`, `DMSwarmRegisterUserDatatypeField()`
920d3a51819SDave May @*/
921d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmInitializeFieldRegister(DM dm)
922d71ae5a4SJacob Faibussowitsch {
9235f50eb2eSDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
9243454631fSDave May 
925521f74f9SMatthew G. Knepley   PetscFunctionBegin;
926cc651181SDave May   if (!swarm->field_registration_initialized) {
9275f50eb2eSDave May     swarm->field_registration_initialized = PETSC_TRUE;
928da81f932SPierre Jolivet     PetscCall(DMSwarmRegisterPetscDatatypeField(dm, DMSwarmField_pid, 1, PETSC_INT64)); /* unique identifier */
9299566063dSJacob Faibussowitsch     PetscCall(DMSwarmRegisterPetscDatatypeField(dm, DMSwarmField_rank, 1, PETSC_INT));  /* used for communication */
930cc651181SDave May   }
9313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9325f50eb2eSDave May }
9335f50eb2eSDave May 
93474d0cae8SMatthew G. Knepley /*@
93520f4b53cSBarry Smith   DMSwarmFinalizeFieldRegister - Finalizes the registration of fields to a `DMSWARM`
936d3a51819SDave May 
93720f4b53cSBarry Smith   Collective
938d3a51819SDave May 
93960225df5SJacob Faibussowitsch   Input Parameter:
94020f4b53cSBarry Smith . dm - a `DMSWARM`
941d3a51819SDave May 
942d3a51819SDave May   Level: beginner
943d3a51819SDave May 
94420f4b53cSBarry Smith   Note:
94520f4b53cSBarry Smith   After `DMSwarmFinalizeFieldRegister()` has been called, no new fields can be defined on the `DMSWARM`.
946d3a51819SDave May 
94720f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmInitializeFieldRegister()`, `DMSwarmRegisterPetscDatatypeField()`,
948db781477SPatrick Sanan           `DMSwarmRegisterUserStructField()`, `DMSwarmRegisterUserDatatypeField()`
949d3a51819SDave May @*/
950d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmFinalizeFieldRegister(DM dm)
951d71ae5a4SJacob Faibussowitsch {
9525f50eb2eSDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
9536845f8f5SDave May 
954521f74f9SMatthew G. Knepley   PetscFunctionBegin;
95548a46eb9SPierre Jolivet   if (!swarm->field_registration_finalized) PetscCall(DMSwarmDataBucketFinalize(swarm->db));
956f0cdbbbaSDave May   swarm->field_registration_finalized = PETSC_TRUE;
9573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9585f50eb2eSDave May }
9595f50eb2eSDave May 
96074d0cae8SMatthew G. Knepley /*@
96120f4b53cSBarry Smith   DMSwarmSetLocalSizes - Sets the length of all registered fields on the `DMSWARM`
962d3a51819SDave May 
96320f4b53cSBarry Smith   Not Collective
964d3a51819SDave May 
96560225df5SJacob Faibussowitsch   Input Parameters:
96620f4b53cSBarry Smith + dm     - a `DMSWARM`
967d3a51819SDave May . nlocal - the length of each registered field
96862741f57SDave May - buffer - the length of the buffer used to efficient dynamic re-sizing
969d3a51819SDave May 
970d3a51819SDave May   Level: beginner
971d3a51819SDave May 
97220f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmGetLocalSize()`
973d3a51819SDave May @*/
974d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmSetLocalSizes(DM dm, PetscInt nlocal, PetscInt buffer)
975d71ae5a4SJacob Faibussowitsch {
9765f50eb2eSDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
9775f50eb2eSDave May 
978521f74f9SMatthew G. Knepley   PetscFunctionBegin;
9799566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMSWARM_SetSizes, 0, 0, 0, 0));
9809566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketSetSizes(swarm->db, nlocal, buffer));
9819566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMSWARM_SetSizes, 0, 0, 0, 0));
9823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9835f50eb2eSDave May }
9845f50eb2eSDave May 
98574d0cae8SMatthew G. Knepley /*@
98620f4b53cSBarry Smith   DMSwarmSetCellDM - Attaches a `DM` to a `DMSWARM`
987d3a51819SDave May 
98820f4b53cSBarry Smith   Collective
989d3a51819SDave May 
99060225df5SJacob Faibussowitsch   Input Parameters:
99120f4b53cSBarry Smith + dm     - a `DMSWARM`
99220f4b53cSBarry Smith - dmcell - the `DM` to attach to the `DMSWARM`
993d3a51819SDave May 
994d3a51819SDave May   Level: beginner
995d3a51819SDave May 
99620f4b53cSBarry Smith   Note:
99720f4b53cSBarry Smith   The attached `DM` (dmcell) will be queried for point location and
99820f4b53cSBarry Smith   neighbor MPI-rank information if `DMSwarmMigrate()` is called.
999d3a51819SDave May 
100020f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmSetType()`, `DMSwarmGetCellDM()`, `DMSwarmMigrate()`
1001d3a51819SDave May @*/
1002d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmSetCellDM(DM dm, DM dmcell)
1003d71ae5a4SJacob Faibussowitsch {
1004521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1005*d52c2f21SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1006*d52c2f21SMatthew G. Knepley   PetscValidHeaderSpecific(dmcell, DM_CLASSID, 2);
1007*d52c2f21SMatthew G. Knepley   PetscCall(DMSwarmPushCellDM(dm, dmcell, 0, NULL, DMSwarmPICField_coor));
10083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1009b16650c8SDave May }
1010b16650c8SDave May 
101174d0cae8SMatthew G. Knepley /*@
101220f4b53cSBarry Smith   DMSwarmGetCellDM - Fetches the attached cell `DM`
1013d3a51819SDave May 
101420f4b53cSBarry Smith   Collective
1015d3a51819SDave May 
101660225df5SJacob Faibussowitsch   Input Parameter:
101720f4b53cSBarry Smith . dm - a `DMSWARM`
1018d3a51819SDave May 
101960225df5SJacob Faibussowitsch   Output Parameter:
102020f4b53cSBarry Smith . dmcell - the `DM` which was attached to the `DMSWARM`
1021d3a51819SDave May 
1022d3a51819SDave May   Level: beginner
1023d3a51819SDave May 
102420f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmSetCellDM()`
1025d3a51819SDave May @*/
1026d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmGetCellDM(DM dm, DM *dmcell)
1027d71ae5a4SJacob Faibussowitsch {
1028fe39f135SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1029521f74f9SMatthew G. Knepley 
1030521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1031*d52c2f21SMatthew G. Knepley   *dmcell = swarm->cellinfo ? swarm->cellinfo[0].dm : NULL;
1032*d52c2f21SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
1033*d52c2f21SMatthew G. Knepley }
1034*d52c2f21SMatthew G. Knepley 
1035*d52c2f21SMatthew G. Knepley static PetscErrorCode CellDMInfoDestroy(CellDMInfo *info)
1036*d52c2f21SMatthew G. Knepley {
1037*d52c2f21SMatthew G. Knepley   PetscFunctionBegin;
1038*d52c2f21SMatthew G. Knepley   for (PetscInt f = 0; f < (*info)->Nf; ++f) PetscCall(PetscFree((*info)->dmFields[f]));
1039*d52c2f21SMatthew G. Knepley   PetscCall(PetscFree((*info)->dmFields));
1040*d52c2f21SMatthew G. Knepley   PetscCall(PetscFree((*info)->coordField));
1041*d52c2f21SMatthew G. Knepley   PetscCall(DMDestroy(&(*info)->dm));
1042*d52c2f21SMatthew G. Knepley   PetscCall(PetscFree(*info));
1043*d52c2f21SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
1044*d52c2f21SMatthew G. Knepley }
1045*d52c2f21SMatthew G. Knepley 
1046*d52c2f21SMatthew G. Knepley /*@
1047*d52c2f21SMatthew G. Knepley   DMSwarmPushCellDM - Attaches a `DM` to a `DMSWARM`
1048*d52c2f21SMatthew G. Knepley 
1049*d52c2f21SMatthew G. Knepley   Collective
1050*d52c2f21SMatthew G. Knepley 
1051*d52c2f21SMatthew G. Knepley   Input Parameters:
1052*d52c2f21SMatthew G. Knepley + sw         - a `DMSWARM`
1053*d52c2f21SMatthew G. Knepley . dmcell     - the `DM` to attach to the `DMSWARM`
1054*d52c2f21SMatthew G. Knepley . Nf         - the number of swarm fields defining the `DM`
1055*d52c2f21SMatthew G. Knepley . dmFields   - an array of field names for the fields defining the `DM`
1056*d52c2f21SMatthew G. Knepley - coordField - the name for the swarm field to use for particle coordinates on the cell `DM`
1057*d52c2f21SMatthew G. Knepley 
1058*d52c2f21SMatthew G. Knepley   Level: beginner
1059*d52c2f21SMatthew G. Knepley 
1060*d52c2f21SMatthew G. Knepley   Note:
1061*d52c2f21SMatthew G. Knepley   The attached `DM` (dmcell) will be queried for point location and
1062*d52c2f21SMatthew G. Knepley   neighbor MPI-rank information if `DMSwarmMigrate()` is called.
1063*d52c2f21SMatthew G. Knepley 
1064*d52c2f21SMatthew G. Knepley .seealso: `DM`, `DMSWARM`, `DMSwarmSetType()`, `DMSwarmPopCellDM()`, `DMSwarmSetCellDM()`, `DMSwarmMigrate()`
1065*d52c2f21SMatthew G. Knepley @*/
1066*d52c2f21SMatthew G. Knepley PetscErrorCode DMSwarmPushCellDM(DM sw, DM dmcell, PetscInt Nf, const char *dmFields[], const char coordField[])
1067*d52c2f21SMatthew G. Knepley {
1068*d52c2f21SMatthew G. Knepley   DM_Swarm  *swarm = (DM_Swarm *)sw->data;
1069*d52c2f21SMatthew G. Knepley   CellDMInfo info;
1070*d52c2f21SMatthew G. Knepley   PetscBool  rebin = swarm->cellinfo ? PETSC_TRUE : PETSC_FALSE;
1071*d52c2f21SMatthew G. Knepley 
1072*d52c2f21SMatthew G. Knepley   PetscFunctionBegin;
1073*d52c2f21SMatthew G. Knepley   PetscValidHeaderSpecific(sw, DM_CLASSID, 1);
1074*d52c2f21SMatthew G. Knepley   PetscValidHeaderSpecific(dmcell, DM_CLASSID, 2);
1075*d52c2f21SMatthew G. Knepley   if (Nf) PetscAssertPointer(dmFields, 4);
1076*d52c2f21SMatthew G. Knepley   PetscAssertPointer(coordField, 5);
1077*d52c2f21SMatthew G. Knepley   PetscCall(PetscNew(&info));
1078*d52c2f21SMatthew G. Knepley   PetscCall(PetscObjectReference((PetscObject)dmcell));
1079*d52c2f21SMatthew G. Knepley   info->dm        = dmcell;
1080*d52c2f21SMatthew G. Knepley   info->Nf        = Nf;
1081*d52c2f21SMatthew G. Knepley   info->next      = swarm->cellinfo;
1082*d52c2f21SMatthew G. Knepley   swarm->cellinfo = info;
1083*d52c2f21SMatthew G. Knepley   // Define the DM fields
1084*d52c2f21SMatthew G. Knepley   PetscCall(PetscMalloc1(info->Nf, &info->dmFields));
1085*d52c2f21SMatthew G. Knepley   for (PetscInt f = 0; f < info->Nf; ++f) PetscCall(PetscStrallocpy(dmFields[f], &info->dmFields[f]));
1086*d52c2f21SMatthew G. Knepley   if (info->Nf) PetscCall(DMSwarmVectorDefineFields(sw, info->Nf, (const char **)info->dmFields));
1087*d52c2f21SMatthew G. Knepley   // Set the coordinate field
1088*d52c2f21SMatthew G. Knepley   PetscCall(PetscStrallocpy(coordField, &info->coordField));
1089*d52c2f21SMatthew G. Knepley   if (info->coordField) PetscCall(DMSwarmSetCoordinateField(sw, info->coordField));
1090*d52c2f21SMatthew G. Knepley   // Rebin the cells and set cell_id field
1091*d52c2f21SMatthew G. Knepley   if (rebin) PetscCall(DMSwarmMigrate(sw, PETSC_FALSE));
1092*d52c2f21SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
1093*d52c2f21SMatthew G. Knepley }
1094*d52c2f21SMatthew G. Knepley 
1095*d52c2f21SMatthew G. Knepley /*@
1096*d52c2f21SMatthew G. Knepley   DMSwarmPopCellDM - Discard the current cell `DM` and restore the previous one, if it exists
1097*d52c2f21SMatthew G. Knepley 
1098*d52c2f21SMatthew G. Knepley   Collective
1099*d52c2f21SMatthew G. Knepley 
1100*d52c2f21SMatthew G. Knepley   Input Parameter:
1101*d52c2f21SMatthew G. Knepley . sw - a `DMSWARM`
1102*d52c2f21SMatthew G. Knepley 
1103*d52c2f21SMatthew G. Knepley   Level: beginner
1104*d52c2f21SMatthew G. Knepley 
1105*d52c2f21SMatthew G. Knepley .seealso: `DM`, `DMSWARM`, `DMSwarmSetType()`, `DMSwarmPushCellDM()`, `DMSwarmSetCellDM()`, `DMSwarmMigrate()`
1106*d52c2f21SMatthew G. Knepley @*/
1107*d52c2f21SMatthew G. Knepley PetscErrorCode DMSwarmPopCellDM(DM sw)
1108*d52c2f21SMatthew G. Knepley {
1109*d52c2f21SMatthew G. Knepley   DM_Swarm  *swarm = (DM_Swarm *)sw->data;
1110*d52c2f21SMatthew G. Knepley   CellDMInfo info  = swarm->cellinfo;
1111*d52c2f21SMatthew G. Knepley 
1112*d52c2f21SMatthew G. Knepley   PetscFunctionBegin;
1113*d52c2f21SMatthew G. Knepley   PetscValidHeaderSpecific(sw, DM_CLASSID, 1);
1114*d52c2f21SMatthew G. Knepley   if (!swarm->cellinfo) PetscFunctionReturn(PETSC_SUCCESS);
1115*d52c2f21SMatthew G. Knepley   if (info->next) {
1116*d52c2f21SMatthew G. Knepley     CellDMInfo newinfo = info->next;
1117*d52c2f21SMatthew G. Knepley 
1118*d52c2f21SMatthew G. Knepley     swarm->cellinfo = info->next;
1119*d52c2f21SMatthew G. Knepley     // Define the DM fields
1120*d52c2f21SMatthew G. Knepley     PetscCall(DMSwarmVectorDefineFields(sw, newinfo->Nf, (const char **)newinfo->dmFields));
1121*d52c2f21SMatthew G. Knepley     // Set the coordinate field
1122*d52c2f21SMatthew G. Knepley     PetscCall(DMSwarmSetCoordinateField(sw, newinfo->coordField));
1123*d52c2f21SMatthew G. Knepley     // Rebin the cells and set cell_id field
1124*d52c2f21SMatthew G. Knepley     PetscCall(DMSwarmMigrate(sw, PETSC_FALSE));
1125*d52c2f21SMatthew G. Knepley   }
1126*d52c2f21SMatthew G. Knepley   PetscCall(CellDMInfoDestroy(&info));
11273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1128fe39f135SDave May }
1129fe39f135SDave May 
113074d0cae8SMatthew G. Knepley /*@
1131d5b43468SJose E. Roman   DMSwarmGetLocalSize - Retrieves the local length of fields registered
1132d3a51819SDave May 
113320f4b53cSBarry Smith   Not Collective
1134d3a51819SDave May 
113560225df5SJacob Faibussowitsch   Input Parameter:
113620f4b53cSBarry Smith . dm - a `DMSWARM`
1137d3a51819SDave May 
113860225df5SJacob Faibussowitsch   Output Parameter:
1139d3a51819SDave May . nlocal - the length of each registered field
1140d3a51819SDave May 
1141d3a51819SDave May   Level: beginner
1142d3a51819SDave May 
114320f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmGetSize()`, `DMSwarmSetLocalSizes()`
1144d3a51819SDave May @*/
1145d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmGetLocalSize(DM dm, PetscInt *nlocal)
1146d71ae5a4SJacob Faibussowitsch {
1147dcf43ee8SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1148dcf43ee8SDave May 
1149521f74f9SMatthew G. Knepley   PetscFunctionBegin;
11509566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, nlocal, NULL, NULL));
11513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1152dcf43ee8SDave May }
1153dcf43ee8SDave May 
115474d0cae8SMatthew G. Knepley /*@
1155d5b43468SJose E. Roman   DMSwarmGetSize - Retrieves the total length of fields registered
1156d3a51819SDave May 
115720f4b53cSBarry Smith   Collective
1158d3a51819SDave May 
115960225df5SJacob Faibussowitsch   Input Parameter:
116020f4b53cSBarry Smith . dm - a `DMSWARM`
1161d3a51819SDave May 
116260225df5SJacob Faibussowitsch   Output Parameter:
1163d3a51819SDave May . n - the total length of each registered field
1164d3a51819SDave May 
1165d3a51819SDave May   Level: beginner
1166d3a51819SDave May 
1167d3a51819SDave May   Note:
116820f4b53cSBarry Smith   This calls `MPI_Allreduce()` upon each call (inefficient but safe)
1169d3a51819SDave May 
117020f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmGetLocalSize()`, `DMSwarmSetLocalSizes()`
1171d3a51819SDave May @*/
1172d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmGetSize(DM dm, PetscInt *n)
1173d71ae5a4SJacob Faibussowitsch {
1174dcf43ee8SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
11755627991aSBarry Smith   PetscInt  nlocal;
1176dcf43ee8SDave May 
1177521f74f9SMatthew G. Knepley   PetscFunctionBegin;
11789566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &nlocal, NULL, NULL));
1179462c564dSBarry Smith   PetscCallMPI(MPIU_Allreduce(&nlocal, n, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm)));
11803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1181dcf43ee8SDave May }
1182dcf43ee8SDave May 
1183cc4c1da9SBarry Smith /*@
118420f4b53cSBarry Smith   DMSwarmRegisterPetscDatatypeField - Register a field to a `DMSWARM` with a native PETSc data type
1185d3a51819SDave May 
118620f4b53cSBarry Smith   Collective
1187d3a51819SDave May 
118860225df5SJacob Faibussowitsch   Input Parameters:
118920f4b53cSBarry Smith + dm        - a `DMSWARM`
1190d3a51819SDave May . fieldname - the textual name to identify this field
1191d3a51819SDave May . blocksize - the number of each data type
119220f4b53cSBarry Smith - type      - a valid PETSc data type (`PETSC_CHAR`, `PETSC_SHORT`, `PETSC_INT`, `PETSC_FLOAT`, `PETSC_REAL`, `PETSC_LONG`)
1193d3a51819SDave May 
1194d3a51819SDave May   Level: beginner
1195d3a51819SDave May 
1196d3a51819SDave May   Notes:
11978b8a3813SDave May   The textual name for each registered field must be unique.
1198d3a51819SDave May 
119920f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterUserStructField()`, `DMSwarmRegisterUserDatatypeField()`
1200d3a51819SDave May @*/
1201d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmRegisterPetscDatatypeField(DM dm, const char fieldname[], PetscInt blocksize, PetscDataType type)
1202d71ae5a4SJacob Faibussowitsch {
1203b62e03f8SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1204b62e03f8SDave May   size_t    size;
1205b62e03f8SDave May 
1206521f74f9SMatthew G. Knepley   PetscFunctionBegin;
120728b400f6SJacob Faibussowitsch   PetscCheck(swarm->field_registration_initialized, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Must call DMSwarmInitializeFieldRegister() first");
120828b400f6SJacob Faibussowitsch   PetscCheck(!swarm->field_registration_finalized, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Cannot register additional fields after calling DMSwarmFinalizeFieldRegister() first");
12095f50eb2eSDave May 
121008401ef6SPierre Jolivet   PetscCheck(type != PETSC_OBJECT, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}");
121108401ef6SPierre Jolivet   PetscCheck(type != PETSC_FUNCTION, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}");
121208401ef6SPierre Jolivet   PetscCheck(type != PETSC_STRING, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}");
121308401ef6SPierre Jolivet   PetscCheck(type != PETSC_STRUCT, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}");
121408401ef6SPierre Jolivet   PetscCheck(type != PETSC_DATATYPE_UNKNOWN, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}");
1215b62e03f8SDave May 
12169566063dSJacob Faibussowitsch   PetscCall(PetscDataTypeGetSize(type, &size));
1217b62e03f8SDave May   /* Load a specific data type into data bucket, specifying textual name and its size in bytes */
12189566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketRegisterField(swarm->db, "DMSwarmRegisterPetscDatatypeField", fieldname, blocksize * size, NULL));
121952c3ed93SDave May   {
122077048351SPatrick Sanan     DMSwarmDataField gfield;
122152c3ed93SDave May 
12229566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield));
12239566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataFieldSetBlockSize(gfield, blocksize));
122452c3ed93SDave May   }
1225b62e03f8SDave May   swarm->db->field[swarm->db->nfields - 1]->petsc_type = type;
12263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1227b62e03f8SDave May }
1228b62e03f8SDave May 
1229d3a51819SDave May /*@C
123020f4b53cSBarry Smith   DMSwarmRegisterUserStructField - Register a user defined struct to a `DMSWARM`
1231d3a51819SDave May 
123220f4b53cSBarry Smith   Collective
1233d3a51819SDave May 
123460225df5SJacob Faibussowitsch   Input Parameters:
123520f4b53cSBarry Smith + dm        - a `DMSWARM`
1236d3a51819SDave May . fieldname - the textual name to identify this field
123762741f57SDave May - size      - the size in bytes of the user struct of each data type
1238d3a51819SDave May 
1239d3a51819SDave May   Level: beginner
1240d3a51819SDave May 
124120f4b53cSBarry Smith   Note:
12428b8a3813SDave May   The textual name for each registered field must be unique.
1243d3a51819SDave May 
124420f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmRegisterUserDatatypeField()`
1245d3a51819SDave May @*/
1246d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmRegisterUserStructField(DM dm, const char fieldname[], size_t size)
1247d71ae5a4SJacob Faibussowitsch {
1248b62e03f8SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1249b62e03f8SDave May 
1250521f74f9SMatthew G. Knepley   PetscFunctionBegin;
12519566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketRegisterField(swarm->db, "DMSwarmRegisterUserStructField", fieldname, size, NULL));
1252b62e03f8SDave May   swarm->db->field[swarm->db->nfields - 1]->petsc_type = PETSC_STRUCT;
12533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1254b62e03f8SDave May }
1255b62e03f8SDave May 
1256d3a51819SDave May /*@C
125720f4b53cSBarry Smith   DMSwarmRegisterUserDatatypeField - Register a user defined data type to a `DMSWARM`
1258d3a51819SDave May 
125920f4b53cSBarry Smith   Collective
1260d3a51819SDave May 
126160225df5SJacob Faibussowitsch   Input Parameters:
126220f4b53cSBarry Smith + dm        - a `DMSWARM`
1263d3a51819SDave May . fieldname - the textual name to identify this field
1264d3a51819SDave May . size      - the size in bytes of the user data type
126562741f57SDave May - blocksize - the number of each data type
1266d3a51819SDave May 
1267d3a51819SDave May   Level: beginner
1268d3a51819SDave May 
126920f4b53cSBarry Smith   Note:
12708b8a3813SDave May   The textual name for each registered field must be unique.
1271d3a51819SDave May 
127242747ad1SJacob Faibussowitsch .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmRegisterUserStructField()`
1273d3a51819SDave May @*/
1274d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmRegisterUserDatatypeField(DM dm, const char fieldname[], size_t size, PetscInt blocksize)
1275d71ae5a4SJacob Faibussowitsch {
1276b62e03f8SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1277b62e03f8SDave May 
1278521f74f9SMatthew G. Knepley   PetscFunctionBegin;
12799566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketRegisterField(swarm->db, "DMSwarmRegisterUserDatatypeField", fieldname, blocksize * size, NULL));
1280320740a0SDave May   {
128177048351SPatrick Sanan     DMSwarmDataField gfield;
1282320740a0SDave May 
12839566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield));
12849566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataFieldSetBlockSize(gfield, blocksize));
1285320740a0SDave May   }
1286b62e03f8SDave May   swarm->db->field[swarm->db->nfields - 1]->petsc_type = PETSC_DATATYPE_UNKNOWN;
12873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1288b62e03f8SDave May }
1289b62e03f8SDave May 
1290d3a51819SDave May /*@C
1291d3a51819SDave May   DMSwarmGetField - Get access to the underlying array storing all entries associated with a registered field
1292d3a51819SDave May 
1293cc4c1da9SBarry Smith   Not Collective, No Fortran Support
1294d3a51819SDave May 
129560225df5SJacob Faibussowitsch   Input Parameters:
129620f4b53cSBarry Smith + dm        - a `DMSWARM`
129762741f57SDave May - fieldname - the textual name to identify this field
1298d3a51819SDave May 
129960225df5SJacob Faibussowitsch   Output Parameters:
130062741f57SDave May + blocksize - the number of each data type
1301d3a51819SDave May . type      - the data type
130262741f57SDave May - data      - pointer to raw array
1303d3a51819SDave May 
1304d3a51819SDave May   Level: beginner
1305d3a51819SDave May 
1306d3a51819SDave May   Notes:
130720f4b53cSBarry Smith   The array must be returned using a matching call to `DMSwarmRestoreField()`.
1308d3a51819SDave May 
130920f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRestoreField()`
1310d3a51819SDave May @*/
1311d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmGetField(DM dm, const char fieldname[], PetscInt *blocksize, PetscDataType *type, void **data)
1312d71ae5a4SJacob Faibussowitsch {
1313b62e03f8SDave May   DM_Swarm        *swarm = (DM_Swarm *)dm->data;
131477048351SPatrick Sanan   DMSwarmDataField gfield;
1315b62e03f8SDave May 
1316521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1317ea7032a0SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
13189566063dSJacob Faibussowitsch   if (!swarm->issetup) PetscCall(DMSetUp(dm));
13199566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield));
13209566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataFieldGetAccess(gfield));
13219566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataFieldGetEntries(gfield, data));
1322ad540459SPierre Jolivet   if (blocksize) *blocksize = gfield->bs;
1323ad540459SPierre Jolivet   if (type) *type = gfield->petsc_type;
13243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1325b62e03f8SDave May }
1326b62e03f8SDave May 
1327d3a51819SDave May /*@C
1328d3a51819SDave May   DMSwarmRestoreField - Restore access to the underlying array storing all entries associated with a registered field
1329d3a51819SDave May 
1330cc4c1da9SBarry Smith   Not Collective, No Fortran Support
1331d3a51819SDave May 
133260225df5SJacob Faibussowitsch   Input Parameters:
133320f4b53cSBarry Smith + dm        - a `DMSWARM`
133462741f57SDave May - fieldname - the textual name to identify this field
1335d3a51819SDave May 
133660225df5SJacob Faibussowitsch   Output Parameters:
133762741f57SDave May + blocksize - the number of each data type
1338d3a51819SDave May . type      - the data type
133962741f57SDave May - data      - pointer to raw array
1340d3a51819SDave May 
1341d3a51819SDave May   Level: beginner
1342d3a51819SDave May 
1343d3a51819SDave May   Notes:
134420f4b53cSBarry Smith   The user must call `DMSwarmGetField()` prior to calling `DMSwarmRestoreField()`.
1345d3a51819SDave May 
134620f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmGetField()`
1347d3a51819SDave May @*/
1348d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmRestoreField(DM dm, const char fieldname[], PetscInt *blocksize, PetscDataType *type, void **data)
1349d71ae5a4SJacob Faibussowitsch {
1350b62e03f8SDave May   DM_Swarm        *swarm = (DM_Swarm *)dm->data;
135177048351SPatrick Sanan   DMSwarmDataField gfield;
1352b62e03f8SDave May 
1353521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1354ea7032a0SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
13559566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield));
13569566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataFieldRestoreAccess(gfield));
1357b62e03f8SDave May   if (data) *data = NULL;
13583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1359b62e03f8SDave May }
1360b62e03f8SDave May 
1361*d52c2f21SMatthew G. Knepley /*@
1362*d52c2f21SMatthew G. Knepley   DMSwarmGetCoordinateField - Get the name of the field holding particle coordinates
1363*d52c2f21SMatthew G. Knepley 
1364*d52c2f21SMatthew G. Knepley   Not Collective
1365*d52c2f21SMatthew G. Knepley 
1366*d52c2f21SMatthew G. Knepley   Input Parameter:
1367*d52c2f21SMatthew G. Knepley . sw - a `DMSWARM`
1368*d52c2f21SMatthew G. Knepley 
1369*d52c2f21SMatthew G. Knepley   Output Parameter:
1370*d52c2f21SMatthew G. Knepley . fieldname - the name of  the coordinate field
1371*d52c2f21SMatthew G. Knepley 
1372*d52c2f21SMatthew G. Knepley   Level: intermediate
1373*d52c2f21SMatthew G. Knepley 
1374*d52c2f21SMatthew G. Knepley .seealso: `DM`, `DMSWARM`, `DMSwarmSetCoordinateField()`, `DMSwarmGetField()`, `DMSwarmSetCellDM()`, `DMSwarmMigrate()`
1375*d52c2f21SMatthew G. Knepley @*/
1376*d52c2f21SMatthew G. Knepley PetscErrorCode DMSwarmGetCoordinateField(DM sw, const char *fieldname[])
1377*d52c2f21SMatthew G. Knepley {
1378*d52c2f21SMatthew G. Knepley   DM_Swarm *swarm = (DM_Swarm *)sw->data;
1379*d52c2f21SMatthew G. Knepley 
1380*d52c2f21SMatthew G. Knepley   PetscFunctionBegin;
1381*d52c2f21SMatthew G. Knepley   PetscValidHeaderSpecificType(sw, DM_CLASSID, 1, DMSWARM);
1382*d52c2f21SMatthew G. Knepley   PetscAssertPointer(fieldname, 2);
1383*d52c2f21SMatthew G. Knepley   *fieldname = swarm->coord_name;
1384*d52c2f21SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
1385*d52c2f21SMatthew G. Knepley }
1386*d52c2f21SMatthew G. Knepley 
1387*d52c2f21SMatthew G. Knepley /*@
1388*d52c2f21SMatthew G. Knepley   DMSwarmSetCoordinateField - Set the name of the field holding particle coordinates
1389*d52c2f21SMatthew G. Knepley 
1390*d52c2f21SMatthew G. Knepley   Not Collective
1391*d52c2f21SMatthew G. Knepley 
1392*d52c2f21SMatthew G. Knepley   Input Parameters:
1393*d52c2f21SMatthew G. Knepley + sw        - a `DMSWARM`
1394*d52c2f21SMatthew G. Knepley - fieldname - the name of  the coordinate field
1395*d52c2f21SMatthew G. Knepley 
1396*d52c2f21SMatthew G. Knepley   Level: intermediate
1397*d52c2f21SMatthew G. Knepley 
1398*d52c2f21SMatthew G. Knepley .seealso: `DM`, `DMSWARM`, `DMSwarmGetCoordinateField()`, `DMSwarmGetField()`, `DMSwarmSetCellDM()`, `DMSwarmMigrate()`
1399*d52c2f21SMatthew G. Knepley @*/
1400*d52c2f21SMatthew G. Knepley PetscErrorCode DMSwarmSetCoordinateField(DM sw, const char fieldname[])
1401*d52c2f21SMatthew G. Knepley {
1402*d52c2f21SMatthew G. Knepley   DM_Swarm *swarm = (DM_Swarm *)sw->data;
1403*d52c2f21SMatthew G. Knepley 
1404*d52c2f21SMatthew G. Knepley   PetscFunctionBegin;
1405*d52c2f21SMatthew G. Knepley   PetscValidHeaderSpecificType(sw, DM_CLASSID, 1, DMSWARM);
1406*d52c2f21SMatthew G. Knepley   PetscAssertPointer(fieldname, 2);
1407*d52c2f21SMatthew G. Knepley   PetscCall(PetscFree(swarm->coord_name));
1408*d52c2f21SMatthew G. Knepley   PetscCall(PetscStrallocpy(fieldname, (char **)&swarm->coord_name));
1409*d52c2f21SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
1410*d52c2f21SMatthew G. Knepley }
1411*d52c2f21SMatthew G. Knepley 
14120bf7c1c5SMatthew G. Knepley PetscErrorCode DMSwarmGetFieldInfo(DM dm, const char fieldname[], PetscInt *blocksize, PetscDataType *type)
14130bf7c1c5SMatthew G. Knepley {
14140bf7c1c5SMatthew G. Knepley   DM_Swarm        *swarm = (DM_Swarm *)dm->data;
14150bf7c1c5SMatthew G. Knepley   DMSwarmDataField gfield;
14160bf7c1c5SMatthew G. Knepley 
14170bf7c1c5SMatthew G. Knepley   PetscFunctionBegin;
14180bf7c1c5SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
14190bf7c1c5SMatthew G. Knepley   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield));
14200bf7c1c5SMatthew G. Knepley   if (blocksize) *blocksize = gfield->bs;
14210bf7c1c5SMatthew G. Knepley   if (type) *type = gfield->petsc_type;
14220bf7c1c5SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
14230bf7c1c5SMatthew G. Knepley }
14240bf7c1c5SMatthew G. Knepley 
142574d0cae8SMatthew G. Knepley /*@
142620f4b53cSBarry Smith   DMSwarmAddPoint - Add space for one new point in the `DMSWARM`
1427d3a51819SDave May 
142820f4b53cSBarry Smith   Not Collective
1429d3a51819SDave May 
143060225df5SJacob Faibussowitsch   Input Parameter:
143120f4b53cSBarry Smith . dm - a `DMSWARM`
1432d3a51819SDave May 
1433d3a51819SDave May   Level: beginner
1434d3a51819SDave May 
1435d3a51819SDave May   Notes:
14368b8a3813SDave May   The new point will have all fields initialized to zero.
1437d3a51819SDave May 
143820f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmAddNPoints()`
1439d3a51819SDave May @*/
1440d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmAddPoint(DM dm)
1441d71ae5a4SJacob Faibussowitsch {
1442cb1d1399SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1443cb1d1399SDave May 
1444521f74f9SMatthew G. Knepley   PetscFunctionBegin;
14459566063dSJacob Faibussowitsch   if (!swarm->issetup) PetscCall(DMSetUp(dm));
14469566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMSWARM_AddPoints, 0, 0, 0, 0));
14479566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketAddPoint(swarm->db));
14489566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMSWARM_AddPoints, 0, 0, 0, 0));
14493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1450cb1d1399SDave May }
1451cb1d1399SDave May 
145274d0cae8SMatthew G. Knepley /*@
145320f4b53cSBarry Smith   DMSwarmAddNPoints - Add space for a number of new points in the `DMSWARM`
1454d3a51819SDave May 
145520f4b53cSBarry Smith   Not Collective
1456d3a51819SDave May 
145760225df5SJacob Faibussowitsch   Input Parameters:
145820f4b53cSBarry Smith + dm      - a `DMSWARM`
145962741f57SDave May - npoints - the number of new points to add
1460d3a51819SDave May 
1461d3a51819SDave May   Level: beginner
1462d3a51819SDave May 
1463d3a51819SDave May   Notes:
14648b8a3813SDave May   The new point will have all fields initialized to zero.
1465d3a51819SDave May 
146620f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmAddPoint()`
1467d3a51819SDave May @*/
1468d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmAddNPoints(DM dm, PetscInt npoints)
1469d71ae5a4SJacob Faibussowitsch {
1470cb1d1399SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1471cb1d1399SDave May   PetscInt  nlocal;
1472cb1d1399SDave May 
1473521f74f9SMatthew G. Knepley   PetscFunctionBegin;
14749566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMSWARM_AddPoints, 0, 0, 0, 0));
14759566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &nlocal, NULL, NULL));
1476cb1d1399SDave May   nlocal = nlocal + npoints;
14779566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketSetSizes(swarm->db, nlocal, DMSWARM_DATA_BUCKET_BUFFER_DEFAULT));
14789566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMSWARM_AddPoints, 0, 0, 0, 0));
14793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1480cb1d1399SDave May }
1481cb1d1399SDave May 
148274d0cae8SMatthew G. Knepley /*@
148320f4b53cSBarry Smith   DMSwarmRemovePoint - Remove the last point from the `DMSWARM`
1484d3a51819SDave May 
148520f4b53cSBarry Smith   Not Collective
1486d3a51819SDave May 
148760225df5SJacob Faibussowitsch   Input Parameter:
148820f4b53cSBarry Smith . dm - a `DMSWARM`
1489d3a51819SDave May 
1490d3a51819SDave May   Level: beginner
1491d3a51819SDave May 
149220f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRemovePointAtIndex()`
1493d3a51819SDave May @*/
1494d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmRemovePoint(DM dm)
1495d71ae5a4SJacob Faibussowitsch {
1496cb1d1399SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1497cb1d1399SDave May 
1498521f74f9SMatthew G. Knepley   PetscFunctionBegin;
14999566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMSWARM_RemovePoints, 0, 0, 0, 0));
15009566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketRemovePoint(swarm->db));
15019566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMSWARM_RemovePoints, 0, 0, 0, 0));
15023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1503cb1d1399SDave May }
1504cb1d1399SDave May 
150574d0cae8SMatthew G. Knepley /*@
150620f4b53cSBarry Smith   DMSwarmRemovePointAtIndex - Removes a specific point from the `DMSWARM`
1507d3a51819SDave May 
150820f4b53cSBarry Smith   Not Collective
1509d3a51819SDave May 
151060225df5SJacob Faibussowitsch   Input Parameters:
151120f4b53cSBarry Smith + dm  - a `DMSWARM`
151262741f57SDave May - idx - index of point to remove
1513d3a51819SDave May 
1514d3a51819SDave May   Level: beginner
1515d3a51819SDave May 
151620f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRemovePoint()`
1517d3a51819SDave May @*/
1518d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmRemovePointAtIndex(DM dm, PetscInt idx)
1519d71ae5a4SJacob Faibussowitsch {
1520cb1d1399SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1521cb1d1399SDave May 
1522521f74f9SMatthew G. Knepley   PetscFunctionBegin;
15239566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMSWARM_RemovePoints, 0, 0, 0, 0));
15249566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketRemovePointAtIndex(swarm->db, idx));
15259566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMSWARM_RemovePoints, 0, 0, 0, 0));
15263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1527cb1d1399SDave May }
1528b62e03f8SDave May 
152974d0cae8SMatthew G. Knepley /*@
153020f4b53cSBarry Smith   DMSwarmCopyPoint - Copy point pj to point pi in the `DMSWARM`
1531ba4fc9c6SDave May 
153220f4b53cSBarry Smith   Not Collective
1533ba4fc9c6SDave May 
153460225df5SJacob Faibussowitsch   Input Parameters:
153520f4b53cSBarry Smith + dm - a `DMSWARM`
1536ba4fc9c6SDave May . pi - the index of the point to copy
1537ba4fc9c6SDave May - pj - the point index where the copy should be located
1538ba4fc9c6SDave May 
1539ba4fc9c6SDave May   Level: beginner
1540ba4fc9c6SDave May 
154120f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRemovePoint()`
1542ba4fc9c6SDave May @*/
1543d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmCopyPoint(DM dm, PetscInt pi, PetscInt pj)
1544d71ae5a4SJacob Faibussowitsch {
1545ba4fc9c6SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1546ba4fc9c6SDave May 
1547ba4fc9c6SDave May   PetscFunctionBegin;
15489566063dSJacob Faibussowitsch   if (!swarm->issetup) PetscCall(DMSetUp(dm));
15499566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketCopyPoint(swarm->db, pi, swarm->db, pj));
15503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1551ba4fc9c6SDave May }
1552ba4fc9c6SDave May 
155366976f2fSJacob Faibussowitsch static PetscErrorCode DMSwarmMigrate_Basic(DM dm, PetscBool remove_sent_points)
1554d71ae5a4SJacob Faibussowitsch {
1555521f74f9SMatthew G. Knepley   PetscFunctionBegin;
15569566063dSJacob Faibussowitsch   PetscCall(DMSwarmMigrate_Push_Basic(dm, remove_sent_points));
15573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15583454631fSDave May }
15593454631fSDave May 
156074d0cae8SMatthew G. Knepley /*@
156120f4b53cSBarry Smith   DMSwarmMigrate - Relocates points defined in the `DMSWARM` to other MPI-ranks
1562d3a51819SDave May 
156320f4b53cSBarry Smith   Collective
1564d3a51819SDave May 
156560225df5SJacob Faibussowitsch   Input Parameters:
156620f4b53cSBarry Smith + dm                 - the `DMSWARM`
156762741f57SDave May - remove_sent_points - flag indicating if sent points should be removed from the current MPI-rank
1568d3a51819SDave May 
1569d3a51819SDave May   Level: advanced
1570d3a51819SDave May 
157120f4b53cSBarry Smith   Notes:
157220f4b53cSBarry Smith   The `DM` will be modified to accommodate received points.
157320f4b53cSBarry Smith   If `remove_sent_points` is `PETSC_TRUE`, any points that were sent will be removed from the `DM`.
157420f4b53cSBarry Smith   Different styles of migration are supported. See `DMSwarmSetMigrateType()`.
157520f4b53cSBarry Smith 
157620f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmSetMigrateType()`
1577d3a51819SDave May @*/
1578d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmMigrate(DM dm, PetscBool remove_sent_points)
1579d71ae5a4SJacob Faibussowitsch {
1580f0cdbbbaSDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1581f0cdbbbaSDave May 
1582521f74f9SMatthew G. Knepley   PetscFunctionBegin;
15839566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMSWARM_Migrate, 0, 0, 0, 0));
1584f0cdbbbaSDave May   switch (swarm->migrate_type) {
1585d71ae5a4SJacob Faibussowitsch   case DMSWARM_MIGRATE_BASIC:
1586d71ae5a4SJacob Faibussowitsch     PetscCall(DMSwarmMigrate_Basic(dm, remove_sent_points));
1587d71ae5a4SJacob Faibussowitsch     break;
1588d71ae5a4SJacob Faibussowitsch   case DMSWARM_MIGRATE_DMCELLNSCATTER:
1589d71ae5a4SJacob Faibussowitsch     PetscCall(DMSwarmMigrate_CellDMScatter(dm, remove_sent_points));
1590d71ae5a4SJacob Faibussowitsch     break;
1591d71ae5a4SJacob Faibussowitsch   case DMSWARM_MIGRATE_DMCELLEXACT:
1592d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_MIGRATE_DMCELLEXACT not implemented");
1593d71ae5a4SJacob Faibussowitsch   case DMSWARM_MIGRATE_USER:
1594d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_MIGRATE_USER not implemented");
1595d71ae5a4SJacob Faibussowitsch   default:
1596d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_MIGRATE type unknown");
1597f0cdbbbaSDave May   }
15989566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMSWARM_Migrate, 0, 0, 0, 0));
15999566063dSJacob Faibussowitsch   PetscCall(DMClearGlobalVectors(dm));
16003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16013454631fSDave May }
16023454631fSDave May 
1603f0cdbbbaSDave May PetscErrorCode DMSwarmMigrate_GlobalToLocal_Basic(DM dm, PetscInt *globalsize);
1604f0cdbbbaSDave May 
1605d3a51819SDave May /*
1606d3a51819SDave May  DMSwarmCollectViewCreate
1607d3a51819SDave May 
1608d3a51819SDave May  * Applies a collection method and gathers point neighbour points into dm
1609d3a51819SDave May 
1610d3a51819SDave May  Notes:
16118b8a3813SDave May  Users should call DMSwarmCollectViewDestroy() after
1612d3a51819SDave May  they have finished computations associated with the collected points
1613d3a51819SDave May */
1614d3a51819SDave May 
161574d0cae8SMatthew G. Knepley /*@
1616d3a51819SDave May   DMSwarmCollectViewCreate - Applies a collection method and gathers points
161720f4b53cSBarry Smith   in neighbour ranks into the `DMSWARM`
1618d3a51819SDave May 
161920f4b53cSBarry Smith   Collective
1620d3a51819SDave May 
162160225df5SJacob Faibussowitsch   Input Parameter:
162220f4b53cSBarry Smith . dm - the `DMSWARM`
1623d3a51819SDave May 
1624d3a51819SDave May   Level: advanced
1625d3a51819SDave May 
162620f4b53cSBarry Smith   Notes:
162720f4b53cSBarry Smith   Users should call `DMSwarmCollectViewDestroy()` after
162820f4b53cSBarry Smith   they have finished computations associated with the collected points
16290764c050SBarry Smith 
163020f4b53cSBarry Smith   Different collect methods are supported. See `DMSwarmSetCollectType()`.
163120f4b53cSBarry Smith 
16320764c050SBarry Smith   Developer Note:
16330764c050SBarry Smith   Create and Destroy routines create new objects that can get destroyed, they do not change the state
16340764c050SBarry Smith   of the current object.
16350764c050SBarry Smith 
163620f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmCollectViewDestroy()`, `DMSwarmSetCollectType()`
1637d3a51819SDave May @*/
1638d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmCollectViewCreate(DM dm)
1639d71ae5a4SJacob Faibussowitsch {
16402712d1f2SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
16412712d1f2SDave May   PetscInt  ng;
16422712d1f2SDave May 
1643521f74f9SMatthew G. Knepley   PetscFunctionBegin;
164428b400f6SJacob Faibussowitsch   PetscCheck(!swarm->collect_view_active, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "CollectView currently active");
16459566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetLocalSize(dm, &ng));
1646480eef7bSDave May   switch (swarm->collect_type) {
1647d71ae5a4SJacob Faibussowitsch   case DMSWARM_COLLECT_BASIC:
1648d71ae5a4SJacob Faibussowitsch     PetscCall(DMSwarmMigrate_GlobalToLocal_Basic(dm, &ng));
1649d71ae5a4SJacob Faibussowitsch     break;
1650d71ae5a4SJacob Faibussowitsch   case DMSWARM_COLLECT_DMDABOUNDINGBOX:
1651d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_COLLECT_DMDABOUNDINGBOX not implemented");
1652d71ae5a4SJacob Faibussowitsch   case DMSWARM_COLLECT_GENERAL:
1653d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_COLLECT_GENERAL not implemented");
1654d71ae5a4SJacob Faibussowitsch   default:
1655d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_COLLECT type unknown");
1656480eef7bSDave May   }
1657480eef7bSDave May   swarm->collect_view_active       = PETSC_TRUE;
1658480eef7bSDave May   swarm->collect_view_reset_nlocal = ng;
16593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16602712d1f2SDave May }
16612712d1f2SDave May 
166274d0cae8SMatthew G. Knepley /*@
166320f4b53cSBarry Smith   DMSwarmCollectViewDestroy - Resets the `DMSWARM` to the size prior to calling `DMSwarmCollectViewCreate()`
1664d3a51819SDave May 
166520f4b53cSBarry Smith   Collective
1666d3a51819SDave May 
166760225df5SJacob Faibussowitsch   Input Parameters:
166820f4b53cSBarry Smith . dm - the `DMSWARM`
1669d3a51819SDave May 
1670d3a51819SDave May   Notes:
167120f4b53cSBarry Smith   Users should call `DMSwarmCollectViewCreate()` before this function is called.
1672d3a51819SDave May 
1673d3a51819SDave May   Level: advanced
1674d3a51819SDave May 
16750764c050SBarry Smith   Developer Note:
16760764c050SBarry Smith   Create and Destroy routines create new objects that can get destroyed, they do not change the state
16770764c050SBarry Smith   of the current object.
16780764c050SBarry Smith 
167920f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmCollectViewCreate()`, `DMSwarmSetCollectType()`
1680d3a51819SDave May @*/
1681d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmCollectViewDestroy(DM dm)
1682d71ae5a4SJacob Faibussowitsch {
16832712d1f2SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
16842712d1f2SDave May 
1685521f74f9SMatthew G. Knepley   PetscFunctionBegin;
168628b400f6SJacob Faibussowitsch   PetscCheck(swarm->collect_view_active, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "CollectView is currently not active");
16879566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(dm, swarm->collect_view_reset_nlocal, -1));
1688480eef7bSDave May   swarm->collect_view_active = PETSC_FALSE;
16893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16902712d1f2SDave May }
16913454631fSDave May 
169266976f2fSJacob Faibussowitsch static PetscErrorCode DMSwarmSetUpPIC(DM dm)
1693d71ae5a4SJacob Faibussowitsch {
1694f0cdbbbaSDave May   PetscInt dim;
1695f0cdbbbaSDave May 
1696521f74f9SMatthew G. Knepley   PetscFunctionBegin;
16979566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetNumSpecies(dm, 1));
16989566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
169963a3b9bcSJacob Faibussowitsch   PetscCheck(dim >= 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Dimension must be 1,2,3 - found %" PetscInt_FMT, dim);
170063a3b9bcSJacob Faibussowitsch   PetscCheck(dim <= 3, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Dimension must be 1,2,3 - found %" PetscInt_FMT, dim);
17019566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dm, DMSwarmPICField_coor, dim, PETSC_DOUBLE));
17029566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dm, DMSwarmPICField_cellid, 1, PETSC_INT));
1703*d52c2f21SMatthew G. Knepley   PetscCall(DMSwarmSetCoordinateField(dm, DMSwarmPICField_coor));
17043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1705f0cdbbbaSDave May }
1706f0cdbbbaSDave May 
170774d0cae8SMatthew G. Knepley /*@
170831403186SMatthew G. Knepley   DMSwarmSetPointCoordinatesRandom - Sets initial coordinates for particles in each cell
170931403186SMatthew G. Knepley 
171020f4b53cSBarry Smith   Collective
171131403186SMatthew G. Knepley 
171260225df5SJacob Faibussowitsch   Input Parameters:
171320f4b53cSBarry Smith + dm  - the `DMSWARM`
171420f4b53cSBarry Smith - Npc - The number of particles per cell in the cell `DM`
171531403186SMatthew G. Knepley 
171631403186SMatthew G. Knepley   Level: intermediate
171731403186SMatthew G. Knepley 
171820f4b53cSBarry Smith   Notes:
171920f4b53cSBarry Smith   The user must use `DMSwarmSetCellDM()` to set the cell `DM` first. The particles are placed randomly inside each cell. If only
172020f4b53cSBarry Smith   one particle is in each cell, it is placed at the centroid.
172120f4b53cSBarry Smith 
172220f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmSetCellDM()`
172331403186SMatthew G. Knepley @*/
1724d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmSetPointCoordinatesRandom(DM dm, PetscInt Npc)
1725d71ae5a4SJacob Faibussowitsch {
172631403186SMatthew G. Knepley   DM             cdm;
172731403186SMatthew G. Knepley   PetscRandom    rnd;
172831403186SMatthew G. Knepley   DMPolytopeType ct;
172931403186SMatthew G. Knepley   PetscBool      simplex;
173031403186SMatthew G. Knepley   PetscReal     *centroid, *coords, *xi0, *v0, *J, *invJ, detJ;
173131403186SMatthew G. Knepley   PetscInt       dim, d, cStart, cEnd, c, p;
1732*d52c2f21SMatthew G. Knepley   const char    *coordname;
173331403186SMatthew G. Knepley 
173431403186SMatthew G. Knepley   PetscFunctionBeginUser;
17359566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)dm), &rnd));
17369566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetInterval(rnd, -1.0, 1.0));
17379566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetType(rnd, PETSCRAND48));
173831403186SMatthew G. Knepley 
1739*d52c2f21SMatthew G. Knepley   PetscCall(DMSwarmGetCoordinateField(dm, &coordname));
17409566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetCellDM(dm, &cdm));
17419566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(cdm, &dim));
17429566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(cdm, 0, &cStart, &cEnd));
17439566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCellType(cdm, cStart, &ct));
174431403186SMatthew G. Knepley   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE;
174531403186SMatthew G. Knepley 
17469566063dSJacob Faibussowitsch   PetscCall(PetscMalloc5(dim, &centroid, dim, &xi0, dim, &v0, dim * dim, &J, dim * dim, &invJ));
174731403186SMatthew G. Knepley   for (d = 0; d < dim; ++d) xi0[d] = -1.0;
1748*d52c2f21SMatthew G. Knepley   PetscCall(DMSwarmGetField(dm, coordname, NULL, NULL, (void **)&coords));
174931403186SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
175031403186SMatthew G. Knepley     if (Npc == 1) {
17519566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFVM(cdm, c, NULL, centroid, NULL));
175231403186SMatthew G. Knepley       for (d = 0; d < dim; ++d) coords[c * dim + d] = centroid[d];
175331403186SMatthew G. Knepley     } else {
17549566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFEM(cdm, c, NULL, v0, J, invJ, &detJ)); /* affine */
175531403186SMatthew G. Knepley       for (p = 0; p < Npc; ++p) {
175631403186SMatthew G. Knepley         const PetscInt n   = c * Npc + p;
175731403186SMatthew G. Knepley         PetscReal      sum = 0.0, refcoords[3];
175831403186SMatthew G. Knepley 
175931403186SMatthew G. Knepley         for (d = 0; d < dim; ++d) {
17609566063dSJacob Faibussowitsch           PetscCall(PetscRandomGetValueReal(rnd, &refcoords[d]));
176131403186SMatthew G. Knepley           sum += refcoords[d];
176231403186SMatthew G. Knepley         }
17639371c9d4SSatish Balay         if (simplex && sum > 0.0)
17649371c9d4SSatish Balay           for (d = 0; d < dim; ++d) refcoords[d] -= PetscSqrtReal(dim) * sum;
176531403186SMatthew G. Knepley         CoordinatesRefToReal(dim, dim, xi0, v0, J, refcoords, &coords[n * dim]);
176631403186SMatthew G. Knepley       }
176731403186SMatthew G. Knepley     }
176831403186SMatthew G. Knepley   }
1769*d52c2f21SMatthew G. Knepley   PetscCall(DMSwarmRestoreField(dm, coordname, NULL, NULL, (void **)&coords));
17709566063dSJacob Faibussowitsch   PetscCall(PetscFree5(centroid, xi0, v0, J, invJ));
17719566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rnd));
17723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
177331403186SMatthew G. Knepley }
177431403186SMatthew G. Knepley 
177531403186SMatthew G. Knepley /*@
177620f4b53cSBarry Smith   DMSwarmSetType - Set particular flavor of `DMSWARM`
1777d3a51819SDave May 
177820f4b53cSBarry Smith   Collective
1779d3a51819SDave May 
178060225df5SJacob Faibussowitsch   Input Parameters:
178120f4b53cSBarry Smith + dm    - the `DMSWARM`
178220f4b53cSBarry Smith - stype - the `DMSWARM` type (e.g. `DMSWARM_PIC`)
1783d3a51819SDave May 
1784d3a51819SDave May   Level: advanced
1785d3a51819SDave May 
178620f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmSetMigrateType()`, `DMSwarmSetCollectType()`, `DMSwarmType`, `DMSWARM_PIC`, `DMSWARM_BASIC`
1787d3a51819SDave May @*/
1788d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmSetType(DM dm, DMSwarmType stype)
1789d71ae5a4SJacob Faibussowitsch {
1790f0cdbbbaSDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1791f0cdbbbaSDave May 
1792521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1793f0cdbbbaSDave May   swarm->swarm_type = stype;
179448a46eb9SPierre Jolivet   if (swarm->swarm_type == DMSWARM_PIC) PetscCall(DMSwarmSetUpPIC(dm));
17953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1796f0cdbbbaSDave May }
1797f0cdbbbaSDave May 
179866976f2fSJacob Faibussowitsch static PetscErrorCode DMSetup_Swarm(DM dm)
1799d71ae5a4SJacob Faibussowitsch {
18003454631fSDave May   DM_Swarm   *swarm = (DM_Swarm *)dm->data;
18013454631fSDave May   PetscMPIInt rank;
18023454631fSDave May   PetscInt    p, npoints, *rankval;
18033454631fSDave May 
1804521f74f9SMatthew G. Knepley   PetscFunctionBegin;
18053ba16761SJacob Faibussowitsch   if (swarm->issetup) PetscFunctionReturn(PETSC_SUCCESS);
18063454631fSDave May   swarm->issetup = PETSC_TRUE;
18073454631fSDave May 
1808f0cdbbbaSDave May   if (swarm->swarm_type == DMSWARM_PIC) {
1809f0cdbbbaSDave May     /* check dmcell exists */
1810*d52c2f21SMatthew G. Knepley     PetscCheck(swarm->cellinfo && swarm->cellinfo[0].dm, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "DMSWARM_PIC requires you call DMSwarmSetCellDM() or DMSwarmPushCellDM()");
1811f0cdbbbaSDave May 
1812*d52c2f21SMatthew G. Knepley     if (swarm->cellinfo[0].dm->ops->locatepointssubdomain) {
1813f0cdbbbaSDave May       /* check methods exists for exact ownership identificiation */
18149566063dSJacob Faibussowitsch       PetscCall(PetscInfo(dm, "DMSWARM_PIC: Using method CellDM->ops->LocatePointsSubdomain\n"));
1815f0cdbbbaSDave May       swarm->migrate_type = DMSWARM_MIGRATE_DMCELLEXACT;
1816f0cdbbbaSDave May     } else {
1817f0cdbbbaSDave May       /* check methods exist for point location AND rank neighbor identification */
1818*d52c2f21SMatthew G. Knepley       if (swarm->cellinfo[0].dm->ops->locatepoints) {
18199566063dSJacob Faibussowitsch         PetscCall(PetscInfo(dm, "DMSWARM_PIC: Using method CellDM->LocatePoints\n"));
1820f0cdbbbaSDave May       } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "DMSWARM_PIC requires the method CellDM->ops->locatepoints be defined");
1821f0cdbbbaSDave May 
1822*d52c2f21SMatthew G. Knepley       if (swarm->cellinfo[0].dm->ops->getneighbors) {
18239566063dSJacob Faibussowitsch         PetscCall(PetscInfo(dm, "DMSWARM_PIC: Using method CellDM->GetNeigbors\n"));
1824f0cdbbbaSDave May       } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "DMSWARM_PIC requires the method CellDM->ops->getneighbors be defined");
1825f0cdbbbaSDave May 
1826f0cdbbbaSDave May       swarm->migrate_type = DMSWARM_MIGRATE_DMCELLNSCATTER;
1827f0cdbbbaSDave May     }
1828f0cdbbbaSDave May   }
1829f0cdbbbaSDave May 
18309566063dSJacob Faibussowitsch   PetscCall(DMSwarmFinalizeFieldRegister(dm));
1831f0cdbbbaSDave May 
18323454631fSDave May   /* check some fields were registered */
183308401ef6SPierre Jolivet   PetscCheck(swarm->db->nfields > 2, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "At least one field user must be registered via DMSwarmRegisterXXX()");
18343454631fSDave May 
18353454631fSDave May   /* check local sizes were set */
183608401ef6SPierre Jolivet   PetscCheck(swarm->db->L != -1, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Local sizes must be set via DMSwarmSetLocalSizes()");
18373454631fSDave May 
18383454631fSDave May   /* initialize values in pid and rank placeholders */
18393454631fSDave May   /* TODO: [pid - use MPI_Scan] */
18409566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
18419566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
18429566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
1843ad540459SPierre Jolivet   for (p = 0; p < npoints; p++) rankval[p] = (PetscInt)rank;
18449566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
18453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
18463454631fSDave May }
18473454631fSDave May 
1848dc5f5ce9SDave May extern PetscErrorCode DMSwarmSortDestroy(DMSwarmSort *_ctx);
1849dc5f5ce9SDave May 
185066976f2fSJacob Faibussowitsch static PetscErrorCode DMDestroy_Swarm(DM dm)
1851d71ae5a4SJacob Faibussowitsch {
185257795646SDave May   DM_Swarm  *swarm = (DM_Swarm *)dm->data;
1853*d52c2f21SMatthew G. Knepley   CellDMInfo info  = swarm->cellinfo;
185457795646SDave May 
185557795646SDave May   PetscFunctionBegin;
18563ba16761SJacob Faibussowitsch   if (--swarm->refct > 0) PetscFunctionReturn(PETSC_SUCCESS);
18579566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketDestroy(&swarm->db));
1858*d52c2f21SMatthew G. Knepley   for (PetscInt f = 0; f < swarm->vec_field_num; ++f) PetscCall(PetscFree(swarm->vec_field_names[f]));
1859*d52c2f21SMatthew G. Knepley   PetscCall(PetscFree(swarm->vec_field_names));
1860*d52c2f21SMatthew G. Knepley   PetscCall(PetscFree(swarm->coord_name));
186148a46eb9SPierre Jolivet   if (swarm->sort_context) PetscCall(DMSwarmSortDestroy(&swarm->sort_context));
1862*d52c2f21SMatthew G. Knepley   while (info) {
1863*d52c2f21SMatthew G. Knepley     CellDMInfo tmp = info;
1864*d52c2f21SMatthew G. Knepley 
1865*d52c2f21SMatthew G. Knepley     info = info->next;
1866*d52c2f21SMatthew G. Knepley     PetscCall(CellDMInfoDestroy(&tmp));
1867*d52c2f21SMatthew G. Knepley   }
18689566063dSJacob Faibussowitsch   PetscCall(PetscFree(swarm));
18693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
187057795646SDave May }
187157795646SDave May 
187266976f2fSJacob Faibussowitsch static PetscErrorCode DMSwarmView_Draw(DM dm, PetscViewer viewer)
1873d71ae5a4SJacob Faibussowitsch {
1874a9ee3421SMatthew G. Knepley   DM          cdm;
1875a9ee3421SMatthew G. Knepley   PetscDraw   draw;
187631403186SMatthew G. Knepley   PetscReal  *coords, oldPause, radius = 0.01;
1877a9ee3421SMatthew G. Knepley   PetscInt    Np, p, bs;
1878*d52c2f21SMatthew G. Knepley   const char *coordname;
1879a9ee3421SMatthew G. Knepley 
1880a9ee3421SMatthew G. Knepley   PetscFunctionBegin;
18819566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, ((PetscObject)dm)->prefix, "-dm_view_swarm_radius", &radius, NULL));
18829566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
18839566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetCellDM(dm, &cdm));
18849566063dSJacob Faibussowitsch   PetscCall(PetscDrawGetPause(draw, &oldPause));
18859566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetPause(draw, 0.0));
18869566063dSJacob Faibussowitsch   PetscCall(DMView(cdm, viewer));
18879566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetPause(draw, oldPause));
1888a9ee3421SMatthew G. Knepley 
1889*d52c2f21SMatthew G. Knepley   PetscCall(DMSwarmGetCoordinateField(dm, &coordname));
18909566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetLocalSize(dm, &Np));
1891*d52c2f21SMatthew G. Knepley   PetscCall(DMSwarmGetField(dm, coordname, &bs, NULL, (void **)&coords));
1892a9ee3421SMatthew G. Knepley   for (p = 0; p < Np; ++p) {
1893a9ee3421SMatthew G. Knepley     const PetscInt i = p * bs;
1894a9ee3421SMatthew G. Knepley 
18959566063dSJacob Faibussowitsch     PetscCall(PetscDrawEllipse(draw, coords[i], coords[i + 1], radius, radius, PETSC_DRAW_BLUE));
1896a9ee3421SMatthew G. Knepley   }
1897*d52c2f21SMatthew G. Knepley   PetscCall(DMSwarmRestoreField(dm, coordname, &bs, NULL, (void **)&coords));
18989566063dSJacob Faibussowitsch   PetscCall(PetscDrawFlush(draw));
18999566063dSJacob Faibussowitsch   PetscCall(PetscDrawPause(draw));
19009566063dSJacob Faibussowitsch   PetscCall(PetscDrawSave(draw));
19013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1902a9ee3421SMatthew G. Knepley }
1903a9ee3421SMatthew G. Knepley 
1904d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMView_Swarm_Ascii(DM dm, PetscViewer viewer)
1905d71ae5a4SJacob Faibussowitsch {
19066a5217c0SMatthew G. Knepley   PetscViewerFormat format;
19076a5217c0SMatthew G. Knepley   PetscInt         *sizes;
19086a5217c0SMatthew G. Knepley   PetscInt          dim, Np, maxSize = 17;
19096a5217c0SMatthew G. Knepley   MPI_Comm          comm;
19106a5217c0SMatthew G. Knepley   PetscMPIInt       rank, size, p;
19116a5217c0SMatthew G. Knepley   const char       *name;
19126a5217c0SMatthew G. Knepley 
19136a5217c0SMatthew G. Knepley   PetscFunctionBegin;
19146a5217c0SMatthew G. Knepley   PetscCall(PetscViewerGetFormat(viewer, &format));
19156a5217c0SMatthew G. Knepley   PetscCall(DMGetDimension(dm, &dim));
19166a5217c0SMatthew G. Knepley   PetscCall(DMSwarmGetLocalSize(dm, &Np));
19176a5217c0SMatthew G. Knepley   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
19186a5217c0SMatthew G. Knepley   PetscCallMPI(MPI_Comm_rank(comm, &rank));
19196a5217c0SMatthew G. Knepley   PetscCallMPI(MPI_Comm_size(comm, &size));
19206a5217c0SMatthew G. Knepley   PetscCall(PetscObjectGetName((PetscObject)dm, &name));
192163a3b9bcSJacob Faibussowitsch   if (name) PetscCall(PetscViewerASCIIPrintf(viewer, "%s in %" PetscInt_FMT " dimension%s:\n", name, dim, dim == 1 ? "" : "s"));
192263a3b9bcSJacob Faibussowitsch   else PetscCall(PetscViewerASCIIPrintf(viewer, "Swarm in %" PetscInt_FMT " dimension%s:\n", dim, dim == 1 ? "" : "s"));
19236a5217c0SMatthew G. Knepley   if (size < maxSize) PetscCall(PetscCalloc1(size, &sizes));
19246a5217c0SMatthew G. Knepley   else PetscCall(PetscCalloc1(3, &sizes));
19256a5217c0SMatthew G. Knepley   if (size < maxSize) {
19266a5217c0SMatthew G. Knepley     PetscCallMPI(MPI_Gather(&Np, 1, MPIU_INT, sizes, 1, MPIU_INT, 0, comm));
19276a5217c0SMatthew G. Knepley     PetscCall(PetscViewerASCIIPrintf(viewer, "  Number of particles per rank:"));
19286a5217c0SMatthew G. Knepley     for (p = 0; p < size; ++p) {
19296a5217c0SMatthew G. Knepley       if (rank == 0) PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT, sizes[p]));
19306a5217c0SMatthew G. Knepley     }
19316a5217c0SMatthew G. Knepley   } else {
19326a5217c0SMatthew G. Knepley     PetscInt locMinMax[2] = {Np, Np};
19336a5217c0SMatthew G. Knepley 
19346a5217c0SMatthew G. Knepley     PetscCall(PetscGlobalMinMaxInt(comm, locMinMax, sizes));
19356a5217c0SMatthew G. Knepley     PetscCall(PetscViewerASCIIPrintf(viewer, "  Min/Max of particles per rank: %" PetscInt_FMT "/%" PetscInt_FMT, sizes[0], sizes[1]));
19366a5217c0SMatthew G. Knepley   }
19376a5217c0SMatthew G. Knepley   PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
19386a5217c0SMatthew G. Knepley   PetscCall(PetscFree(sizes));
19396a5217c0SMatthew G. Knepley   if (format == PETSC_VIEWER_ASCII_INFO) {
19406a5217c0SMatthew G. Knepley     PetscInt *cell;
19416a5217c0SMatthew G. Knepley 
19426a5217c0SMatthew G. Knepley     PetscCall(PetscViewerASCIIPrintf(viewer, "  Cells containing each particle:\n"));
19436a5217c0SMatthew G. Knepley     PetscCall(PetscViewerASCIIPushSynchronized(viewer));
19446a5217c0SMatthew G. Knepley     PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&cell));
194563a3b9bcSJacob Faibussowitsch     for (p = 0; p < Np; ++p) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "  p%d: %" PetscInt_FMT "\n", p, cell[p]));
19466a5217c0SMatthew G. Knepley     PetscCall(PetscViewerFlush(viewer));
19476a5217c0SMatthew G. Knepley     PetscCall(PetscViewerASCIIPopSynchronized(viewer));
19486a5217c0SMatthew G. Knepley     PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&cell));
19496a5217c0SMatthew G. Knepley   }
19503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19516a5217c0SMatthew G. Knepley }
19526a5217c0SMatthew G. Knepley 
195366976f2fSJacob Faibussowitsch static PetscErrorCode DMView_Swarm(DM dm, PetscViewer viewer)
1954d71ae5a4SJacob Faibussowitsch {
19555f50eb2eSDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
19565627991aSBarry Smith   PetscBool iascii, ibinary, isvtk, isdraw;
19575627991aSBarry Smith #if defined(PETSC_HAVE_HDF5)
19585627991aSBarry Smith   PetscBool ishdf5;
19595627991aSBarry Smith #endif
19605f50eb2eSDave May 
19615f50eb2eSDave May   PetscFunctionBegin;
19625f50eb2eSDave May   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
19635f50eb2eSDave May   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
19649566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
19659566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &ibinary));
19669566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERVTK, &isvtk));
19675627991aSBarry Smith #if defined(PETSC_HAVE_HDF5)
19689566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
19695627991aSBarry Smith #endif
19709566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
19715f50eb2eSDave May   if (iascii) {
19726a5217c0SMatthew G. Knepley     PetscViewerFormat format;
19736a5217c0SMatthew G. Knepley 
19746a5217c0SMatthew G. Knepley     PetscCall(PetscViewerGetFormat(viewer, &format));
19756a5217c0SMatthew G. Knepley     switch (format) {
1976d71ae5a4SJacob Faibussowitsch     case PETSC_VIEWER_ASCII_INFO_DETAIL:
1977d71ae5a4SJacob Faibussowitsch       PetscCall(DMSwarmDataBucketView(PetscObjectComm((PetscObject)dm), swarm->db, NULL, DATABUCKET_VIEW_STDOUT));
1978d71ae5a4SJacob Faibussowitsch       break;
1979d71ae5a4SJacob Faibussowitsch     default:
1980d71ae5a4SJacob Faibussowitsch       PetscCall(DMView_Swarm_Ascii(dm, viewer));
19816a5217c0SMatthew G. Knepley     }
1982f7d195e4SLawrence Mitchell   } else {
19835f50eb2eSDave May #if defined(PETSC_HAVE_HDF5)
198448a46eb9SPierre Jolivet     if (ishdf5) PetscCall(DMSwarmView_HDF5(dm, viewer));
19855f50eb2eSDave May #endif
198648a46eb9SPierre Jolivet     if (isdraw) PetscCall(DMSwarmView_Draw(dm, viewer));
1987f7d195e4SLawrence Mitchell   }
19883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19895f50eb2eSDave May }
19905f50eb2eSDave May 
1991cc4c1da9SBarry Smith /*@
199220f4b53cSBarry Smith   DMSwarmGetCellSwarm - Extracts a single cell from the `DMSWARM` object, returns it as a single cell `DMSWARM`.
199320f4b53cSBarry Smith   The cell `DM` is filtered for fields of that cell, and the filtered `DM` is used as the cell `DM` of the new swarm object.
1994d0c080abSJoseph Pusztay 
1995d0c080abSJoseph Pusztay   Noncollective
1996d0c080abSJoseph Pusztay 
199760225df5SJacob Faibussowitsch   Input Parameters:
199820f4b53cSBarry Smith + sw        - the `DMSWARM`
19995627991aSBarry Smith . cellID    - the integer id of the cell to be extracted and filtered
200020f4b53cSBarry Smith - cellswarm - The `DMSWARM` to receive the cell
2001d0c080abSJoseph Pusztay 
2002d0c080abSJoseph Pusztay   Level: beginner
2003d0c080abSJoseph Pusztay 
20045627991aSBarry Smith   Notes:
200520f4b53cSBarry Smith   This presently only supports `DMSWARM_PIC` type
20065627991aSBarry Smith 
200720f4b53cSBarry Smith   Should be restored with `DMSwarmRestoreCellSwarm()`
2008d0c080abSJoseph Pusztay 
200920f4b53cSBarry Smith   Changes to this cell of the swarm will be lost if they are made prior to restoring this cell.
201020f4b53cSBarry Smith 
201120f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRestoreCellSwarm()`
2012d0c080abSJoseph Pusztay @*/
2013cc4c1da9SBarry Smith PetscErrorCode DMSwarmGetCellSwarm(DM sw, PetscInt cellID, DM cellswarm)
2014d71ae5a4SJacob Faibussowitsch {
2015d0c080abSJoseph Pusztay   DM_Swarm *original = (DM_Swarm *)sw->data;
2016d0c080abSJoseph Pusztay   DMLabel   label;
2017d0c080abSJoseph Pusztay   DM        dmc, subdmc;
2018d0c080abSJoseph Pusztay   PetscInt *pids, particles, dim;
2019d0c080abSJoseph Pusztay 
2020d0c080abSJoseph Pusztay   PetscFunctionBegin;
2021d0c080abSJoseph Pusztay   /* Configure new swarm */
20229566063dSJacob Faibussowitsch   PetscCall(DMSetType(cellswarm, DMSWARM));
20239566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(sw, &dim));
20249566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(cellswarm, dim));
20259566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetType(cellswarm, DMSWARM_PIC));
2026d0c080abSJoseph Pusztay   /* Destroy the unused, unconfigured data bucket to prevent stragglers in memory */
20279566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketDestroy(&((DM_Swarm *)cellswarm->data)->db));
20289566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortGetAccess(sw));
20299566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortGetNumberOfPointsPerCell(sw, cellID, &particles));
20309566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortGetPointsPerCell(sw, cellID, &particles, &pids));
20319566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketCreateFromSubset(original->db, particles, pids, &((DM_Swarm *)cellswarm->data)->db));
20329566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortRestoreAccess(sw));
2033fe11354eSMatthew G. Knepley   PetscCall(DMSwarmSortRestorePointsPerCell(sw, cellID, &particles, &pids));
20349566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetCellDM(sw, &dmc));
20359566063dSJacob Faibussowitsch   PetscCall(DMLabelCreate(PetscObjectComm((PetscObject)sw), "singlecell", &label));
20369566063dSJacob Faibussowitsch   PetscCall(DMAddLabel(dmc, label));
20379566063dSJacob Faibussowitsch   PetscCall(DMLabelSetValue(label, cellID, 1));
203830cbcd5dSksagiyam   PetscCall(DMPlexFilter(dmc, label, 1, PETSC_FALSE, PETSC_FALSE, NULL, &subdmc));
20399566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetCellDM(cellswarm, subdmc));
20409566063dSJacob Faibussowitsch   PetscCall(DMLabelDestroy(&label));
20413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2042d0c080abSJoseph Pusztay }
2043d0c080abSJoseph Pusztay 
2044cc4c1da9SBarry Smith /*@
204520f4b53cSBarry Smith   DMSwarmRestoreCellSwarm - Restores a `DMSWARM` object obtained with `DMSwarmGetCellSwarm()`. All fields are copied back into the parent swarm.
2046d0c080abSJoseph Pusztay 
2047d0c080abSJoseph Pusztay   Noncollective
2048d0c080abSJoseph Pusztay 
204960225df5SJacob Faibussowitsch   Input Parameters:
205020f4b53cSBarry Smith + sw        - the parent `DMSWARM`
2051d0c080abSJoseph Pusztay . cellID    - the integer id of the cell to be copied back into the parent swarm
2052d0c080abSJoseph Pusztay - cellswarm - the cell swarm object
2053d0c080abSJoseph Pusztay 
2054d0c080abSJoseph Pusztay   Level: beginner
2055d0c080abSJoseph Pusztay 
20565627991aSBarry Smith   Note:
205720f4b53cSBarry Smith   This only supports `DMSWARM_PIC` types of `DMSWARM`s
2058d0c080abSJoseph Pusztay 
205920f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmGetCellSwarm()`
2060d0c080abSJoseph Pusztay @*/
2061cc4c1da9SBarry Smith PetscErrorCode DMSwarmRestoreCellSwarm(DM sw, PetscInt cellID, DM cellswarm)
2062d71ae5a4SJacob Faibussowitsch {
2063d0c080abSJoseph Pusztay   DM        dmc;
2064d0c080abSJoseph Pusztay   PetscInt *pids, particles, p;
2065d0c080abSJoseph Pusztay 
2066d0c080abSJoseph Pusztay   PetscFunctionBegin;
20679566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortGetAccess(sw));
20689566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortGetPointsPerCell(sw, cellID, &particles, &pids));
20699566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortRestoreAccess(sw));
2070d0c080abSJoseph Pusztay   /* Pointwise copy of each particle based on pid. The parent swarm may not be altered during this process. */
207148a46eb9SPierre Jolivet   for (p = 0; p < particles; ++p) PetscCall(DMSwarmDataBucketCopyPoint(((DM_Swarm *)cellswarm->data)->db, pids[p], ((DM_Swarm *)sw->data)->db, pids[p]));
2072d0c080abSJoseph Pusztay   /* Free memory, destroy cell dm */
20739566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetCellDM(cellswarm, &dmc));
20749566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmc));
2075fe11354eSMatthew G. Knepley   PetscCall(DMSwarmSortRestorePointsPerCell(sw, cellID, &particles, &pids));
20763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2077d0c080abSJoseph Pusztay }
2078d0c080abSJoseph Pusztay 
2079*d52c2f21SMatthew G. Knepley /*@
2080*d52c2f21SMatthew G. Knepley   DMSwarmComputeMoments - Compute the first three particle moments for a given field
2081*d52c2f21SMatthew G. Knepley 
2082*d52c2f21SMatthew G. Knepley   Noncollective
2083*d52c2f21SMatthew G. Knepley 
2084*d52c2f21SMatthew G. Knepley   Input Parameters:
2085*d52c2f21SMatthew G. Knepley + sw         - the `DMSWARM`
2086*d52c2f21SMatthew G. Knepley . coordinate - the coordinate field name
2087*d52c2f21SMatthew G. Knepley - weight     - the weight field name
2088*d52c2f21SMatthew G. Knepley 
2089*d52c2f21SMatthew G. Knepley   Output Parameter:
2090*d52c2f21SMatthew G. Knepley . moments - the field moments
2091*d52c2f21SMatthew G. Knepley 
2092*d52c2f21SMatthew G. Knepley   Level: intermediate
2093*d52c2f21SMatthew G. Knepley 
2094*d52c2f21SMatthew G. Knepley   Notes:
2095*d52c2f21SMatthew G. Knepley   The `moments` array should be of length bs + 2, where bs is the block size of the coordinate field.
2096*d52c2f21SMatthew G. Knepley 
2097*d52c2f21SMatthew G. Knepley   The weight field must be a scalar, having blocksize 1.
2098*d52c2f21SMatthew G. Knepley 
2099*d52c2f21SMatthew G. Knepley .seealso: `DM`, `DMSWARM`, `DMPlexComputeMoments()`
2100*d52c2f21SMatthew G. Knepley @*/
2101*d52c2f21SMatthew G. Knepley PetscErrorCode DMSwarmComputeMoments(DM sw, const char coordinate[], const char weight[], PetscReal moments[])
2102*d52c2f21SMatthew G. Knepley {
2103*d52c2f21SMatthew G. Knepley   const PetscReal *coords;
2104*d52c2f21SMatthew G. Knepley   const PetscReal *w;
2105*d52c2f21SMatthew G. Knepley   PetscReal       *mom;
2106*d52c2f21SMatthew G. Knepley   PetscDataType    dtc, dtw;
2107*d52c2f21SMatthew G. Knepley   PetscInt         bsc, bsw, Np;
2108*d52c2f21SMatthew G. Knepley   MPI_Comm         comm;
2109*d52c2f21SMatthew G. Knepley 
2110*d52c2f21SMatthew G. Knepley   PetscFunctionBegin;
2111*d52c2f21SMatthew G. Knepley   PetscValidHeaderSpecific(sw, DM_CLASSID, 1);
2112*d52c2f21SMatthew G. Knepley   PetscAssertPointer(coordinate, 2);
2113*d52c2f21SMatthew G. Knepley   PetscAssertPointer(weight, 3);
2114*d52c2f21SMatthew G. Knepley   PetscAssertPointer(moments, 4);
2115*d52c2f21SMatthew G. Knepley   PetscCall(PetscObjectGetComm((PetscObject)sw, &comm));
2116*d52c2f21SMatthew G. Knepley   PetscCall(DMSwarmGetField(sw, coordinate, &bsc, &dtc, (void **)&coords));
2117*d52c2f21SMatthew G. Knepley   PetscCall(DMSwarmGetField(sw, weight, &bsw, &dtw, (void **)&w));
2118*d52c2f21SMatthew G. Knepley   PetscCheck(dtc == PETSC_REAL, comm, PETSC_ERR_ARG_WRONG, "Coordinate field %s must be real, not %s", coordinate, PetscDataTypes[dtc]);
2119*d52c2f21SMatthew G. Knepley   PetscCheck(dtw == PETSC_REAL, comm, PETSC_ERR_ARG_WRONG, "Weight field %s must be real, not %s", weight, PetscDataTypes[dtw]);
2120*d52c2f21SMatthew G. Knepley   PetscCheck(bsw == 1, comm, PETSC_ERR_ARG_WRONG, "Weight field %s must be a scalar, not blocksize %" PetscInt_FMT, weight, bsw);
2121*d52c2f21SMatthew G. Knepley   PetscCall(DMSwarmGetLocalSize(sw, &Np));
2122*d52c2f21SMatthew G. Knepley   PetscCall(DMGetWorkArray(sw, bsc + 2, MPIU_REAL, &mom));
2123*d52c2f21SMatthew G. Knepley   PetscCall(PetscArrayzero(mom, bsc + 2));
2124*d52c2f21SMatthew G. Knepley   for (PetscInt p = 0; p < Np; ++p) {
2125*d52c2f21SMatthew G. Knepley     const PetscReal *c  = &coords[p * bsc];
2126*d52c2f21SMatthew G. Knepley     const PetscReal  wp = w[p];
2127*d52c2f21SMatthew G. Knepley 
2128*d52c2f21SMatthew G. Knepley     mom[0] += wp;
2129*d52c2f21SMatthew G. Knepley     for (PetscInt d = 0; d < bsc; ++d) {
2130*d52c2f21SMatthew G. Knepley       mom[d + 1] += wp * c[d];
2131*d52c2f21SMatthew G. Knepley       mom[d + bsc + 1] += wp * PetscSqr(c[d]);
2132*d52c2f21SMatthew G. Knepley     }
2133*d52c2f21SMatthew G. Knepley   }
2134*d52c2f21SMatthew G. Knepley   PetscCall(DMSwarmRestoreField(sw, "velocity", NULL, NULL, (void **)&coords));
2135*d52c2f21SMatthew G. Knepley   PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&w));
2136*d52c2f21SMatthew G. Knepley   PetscCallMPI(MPIU_Allreduce(mom, moments, bsc + 2, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)sw)));
2137*d52c2f21SMatthew G. Knepley   PetscCall(DMRestoreWorkArray(sw, bsc + 2, MPIU_REAL, &mom));
2138*d52c2f21SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
2139*d52c2f21SMatthew G. Knepley }
2140*d52c2f21SMatthew G. Knepley 
2141d0c080abSJoseph Pusztay PETSC_INTERN PetscErrorCode DMClone_Swarm(DM, DM *);
2142d0c080abSJoseph Pusztay 
2143d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMInitialize_Swarm(DM sw)
2144d71ae5a4SJacob Faibussowitsch {
2145d0c080abSJoseph Pusztay   PetscFunctionBegin;
2146d0c080abSJoseph Pusztay   sw->ops->view                     = DMView_Swarm;
2147d0c080abSJoseph Pusztay   sw->ops->load                     = NULL;
2148d0c080abSJoseph Pusztay   sw->ops->setfromoptions           = NULL;
2149d0c080abSJoseph Pusztay   sw->ops->clone                    = DMClone_Swarm;
2150d0c080abSJoseph Pusztay   sw->ops->setup                    = DMSetup_Swarm;
2151d0c080abSJoseph Pusztay   sw->ops->createlocalsection       = NULL;
2152adc21957SMatthew G. Knepley   sw->ops->createsectionpermutation = NULL;
2153d0c080abSJoseph Pusztay   sw->ops->createdefaultconstraints = NULL;
2154d0c080abSJoseph Pusztay   sw->ops->createglobalvector       = DMCreateGlobalVector_Swarm;
2155d0c080abSJoseph Pusztay   sw->ops->createlocalvector        = DMCreateLocalVector_Swarm;
2156d0c080abSJoseph Pusztay   sw->ops->getlocaltoglobalmapping  = NULL;
2157d0c080abSJoseph Pusztay   sw->ops->createfieldis            = NULL;
2158d0c080abSJoseph Pusztay   sw->ops->createcoordinatedm       = NULL;
2159d0c080abSJoseph Pusztay   sw->ops->getcoloring              = NULL;
2160d0c080abSJoseph Pusztay   sw->ops->creatematrix             = DMCreateMatrix_Swarm;
2161d0c080abSJoseph Pusztay   sw->ops->createinterpolation      = NULL;
2162d0c080abSJoseph Pusztay   sw->ops->createinjection          = NULL;
2163d0c080abSJoseph Pusztay   sw->ops->createmassmatrix         = DMCreateMassMatrix_Swarm;
2164d0c080abSJoseph Pusztay   sw->ops->refine                   = NULL;
2165d0c080abSJoseph Pusztay   sw->ops->coarsen                  = NULL;
2166d0c080abSJoseph Pusztay   sw->ops->refinehierarchy          = NULL;
2167d0c080abSJoseph Pusztay   sw->ops->coarsenhierarchy         = NULL;
2168d0c080abSJoseph Pusztay   sw->ops->globaltolocalbegin       = NULL;
2169d0c080abSJoseph Pusztay   sw->ops->globaltolocalend         = NULL;
2170d0c080abSJoseph Pusztay   sw->ops->localtoglobalbegin       = NULL;
2171d0c080abSJoseph Pusztay   sw->ops->localtoglobalend         = NULL;
2172d0c080abSJoseph Pusztay   sw->ops->destroy                  = DMDestroy_Swarm;
2173d0c080abSJoseph Pusztay   sw->ops->createsubdm              = NULL;
2174d0c080abSJoseph Pusztay   sw->ops->getdimpoints             = NULL;
2175d0c080abSJoseph Pusztay   sw->ops->locatepoints             = NULL;
21763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2177d0c080abSJoseph Pusztay }
2178d0c080abSJoseph Pusztay 
2179d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode DMClone_Swarm(DM dm, DM *newdm)
2180d71ae5a4SJacob Faibussowitsch {
2181d0c080abSJoseph Pusztay   DM_Swarm *swarm = (DM_Swarm *)dm->data;
2182d0c080abSJoseph Pusztay 
2183d0c080abSJoseph Pusztay   PetscFunctionBegin;
2184d0c080abSJoseph Pusztay   swarm->refct++;
2185d0c080abSJoseph Pusztay   (*newdm)->data = swarm;
21869566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)*newdm, DMSWARM));
21879566063dSJacob Faibussowitsch   PetscCall(DMInitialize_Swarm(*newdm));
21882e56dffeSJoe Pusztay   (*newdm)->dim = dm->dim;
21893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2190d0c080abSJoseph Pusztay }
2191d0c080abSJoseph Pusztay 
2192d3a51819SDave May /*MC
219320f4b53cSBarry Smith  DMSWARM = "swarm" - A `DM` object used to represent arrays of data (fields) of arbitrary data type.
219462741f57SDave May  This implementation was designed for particle methods in which the underlying
2195d3a51819SDave May  data required to be represented is both (i) dynamic in length, (ii) and of arbitrary data type.
2196d3a51819SDave May 
219720f4b53cSBarry Smith  Level: intermediate
219820f4b53cSBarry Smith 
219920f4b53cSBarry Smith   Notes:
220020f4b53cSBarry Smith  User data can be represented by `DMSWARM` through a registering "fields".
220162741f57SDave May  To register a field, the user must provide;
220262741f57SDave May  (a) a unique name;
220362741f57SDave May  (b) the data type (or size in bytes);
220462741f57SDave May  (c) the block size of the data.
2205d3a51819SDave May 
2206d3a51819SDave May  For example, suppose the application requires a unique id, energy, momentum and density to be stored
2207c215a47eSMatthew Knepley  on a set of particles. Then the following code could be used
220820f4b53cSBarry Smith .vb
220920f4b53cSBarry Smith     DMSwarmInitializeFieldRegister(dm)
221020f4b53cSBarry Smith     DMSwarmRegisterPetscDatatypeField(dm,"uid",1,PETSC_LONG);
221120f4b53cSBarry Smith     DMSwarmRegisterPetscDatatypeField(dm,"energy",1,PETSC_REAL);
221220f4b53cSBarry Smith     DMSwarmRegisterPetscDatatypeField(dm,"momentum",3,PETSC_REAL);
221320f4b53cSBarry Smith     DMSwarmRegisterPetscDatatypeField(dm,"density",1,PETSC_FLOAT);
221420f4b53cSBarry Smith     DMSwarmFinalizeFieldRegister(dm)
221520f4b53cSBarry Smith .ve
2216d3a51819SDave May 
221720f4b53cSBarry Smith  The fields represented by `DMSWARM` are dynamic and can be re-sized at any time.
221820f4b53cSBarry Smith  The only restriction imposed by `DMSWARM` is that all fields contain the same number of points.
2219d3a51819SDave May 
2220d3a51819SDave May  To support particle methods, "migration" techniques are provided. These methods migrate data
22215627991aSBarry Smith  between ranks.
2222d3a51819SDave May 
222320f4b53cSBarry Smith  `DMSWARM` supports the methods `DMCreateGlobalVector()` and `DMCreateLocalVector()`.
222420f4b53cSBarry Smith  As a `DMSWARM` may internally define and store values of different data types,
222520f4b53cSBarry Smith  before calling `DMCreateGlobalVector()` or `DMCreateLocalVector()`, the user must inform `DMSWARM` which
222620f4b53cSBarry Smith  fields should be used to define a `Vec` object via
222720f4b53cSBarry Smith    `DMSwarmVectorDefineField()`
2228c215a47eSMatthew Knepley  The specified field can be changed at any time - thereby permitting vectors
2229c215a47eSMatthew Knepley  compatible with different fields to be created.
2230d3a51819SDave May 
223120f4b53cSBarry Smith  A dual representation of fields in the `DMSWARM` and a Vec object is permitted via
223220f4b53cSBarry Smith    `DMSwarmCreateGlobalVectorFromField()`
223320f4b53cSBarry Smith  Here the data defining the field in the `DMSWARM` is shared with a Vec.
2234d3a51819SDave May  This is inherently unsafe if you alter the size of the field at any time between
223520f4b53cSBarry Smith  calls to `DMSwarmCreateGlobalVectorFromField()` and `DMSwarmDestroyGlobalVectorFromField()`.
223620f4b53cSBarry Smith  If the local size of the `DMSWARM` does not match the local size of the global vector
223720f4b53cSBarry Smith  when `DMSwarmDestroyGlobalVectorFromField()` is called, an error is thrown.
2238d3a51819SDave May 
223962741f57SDave May  Additional high-level support is provided for Particle-In-Cell methods.
224020f4b53cSBarry Smith  Please refer to `DMSwarmSetType()`.
224162741f57SDave May 
224220f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMType`, `DMCreate()`, `DMSetType()`
2243d3a51819SDave May M*/
2244cc4c1da9SBarry Smith 
2245d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMCreate_Swarm(DM dm)
2246d71ae5a4SJacob Faibussowitsch {
224757795646SDave May   DM_Swarm *swarm;
224857795646SDave May 
224957795646SDave May   PetscFunctionBegin;
225057795646SDave May   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
22514dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&swarm));
2252f0cdbbbaSDave May   dm->data = swarm;
22539566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketCreate(&swarm->db));
22549566063dSJacob Faibussowitsch   PetscCall(DMSwarmInitializeFieldRegister(dm));
2255*d52c2f21SMatthew G. Knepley   dm->dim                               = 0;
2256d0c080abSJoseph Pusztay   swarm->refct                          = 1;
22573454631fSDave May   swarm->issetup                        = PETSC_FALSE;
2258480eef7bSDave May   swarm->swarm_type                     = DMSWARM_BASIC;
2259480eef7bSDave May   swarm->migrate_type                   = DMSWARM_MIGRATE_BASIC;
2260480eef7bSDave May   swarm->collect_type                   = DMSWARM_COLLECT_BASIC;
226140c453e9SDave May   swarm->migrate_error_on_missing_point = PETSC_FALSE;
2262*d52c2f21SMatthew G. Knepley   swarm->cellinfo                       = NULL;
2263f0cdbbbaSDave May   swarm->collect_view_active            = PETSC_FALSE;
2264f0cdbbbaSDave May   swarm->collect_view_reset_nlocal      = -1;
22659566063dSJacob Faibussowitsch   PetscCall(DMInitialize_Swarm(dm));
22662e956fe4SStefano Zampini   if (SwarmDataFieldId == -1) PetscCall(PetscObjectComposedDataRegister(&SwarmDataFieldId));
22673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
226857795646SDave May }
2269