xref: /petsc/src/dm/impls/swarm/swarm.c (revision 9f196a0264fbaf0568fead3a30c861c7ae4cf663)
119307e5cSMatthew G. Knepley #include "petscdmswarm.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};
21fca3fa51SMatthew G. Knepley const char *DMSwarmRemapTypeNames[]     = {"none", "pfak", "colella", "DMSwarmRemapType", "DMSWARM_REMAP_", NULL};
22ea78f98cSLisandro Dalcin const char *DMSwarmPICLayoutTypeNames[] = {"regular", "gauss", "subdivision", NULL};
23f0cdbbbaSDave May 
24f0cdbbbaSDave May const char DMSwarmField_pid[]     = "DMSwarm_pid";
25f0cdbbbaSDave May const char DMSwarmField_rank[]    = "DMSwarm_rank";
26f0cdbbbaSDave May const char DMSwarmPICField_coor[] = "DMSwarmPIC_coor";
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 
3366976f2fSJacob Faibussowitsch static PetscErrorCode VecView_Swarm_HDF5_Internal(Vec v, PetscViewer viewer)
34d71ae5a4SJacob Faibussowitsch {
3574d0cae8SMatthew G. Knepley   DM        dm;
3674d0cae8SMatthew G. Knepley   PetscReal seqval;
3774d0cae8SMatthew G. Knepley   PetscInt  seqnum, bs;
38eb0c6899SMatthew G. Knepley   PetscBool isseq, ists;
3974d0cae8SMatthew G. Knepley 
4074d0cae8SMatthew G. Knepley   PetscFunctionBegin;
419566063dSJacob Faibussowitsch   PetscCall(VecGetDM(v, &dm));
429566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(v, &bs));
439566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5PushGroup(viewer, "/particle_fields"));
449566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)v, VECSEQ, &isseq));
459566063dSJacob Faibussowitsch   PetscCall(DMGetOutputSequenceNumber(dm, &seqnum, &seqval));
46eb0c6899SMatthew G. Knepley   PetscCall(PetscViewerHDF5IsTimestepping(viewer, &ists));
47eb0c6899SMatthew G. Knepley   if (ists) PetscCall(PetscViewerHDF5SetTimestep(viewer, seqnum));
489566063dSJacob Faibussowitsch   /* PetscCall(DMSequenceView_HDF5(dm, "time", seqnum, (PetscScalar) seqval, viewer)); */
499566063dSJacob Faibussowitsch   PetscCall(VecViewNative(v, viewer));
509566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject)v, "Nc", PETSC_INT, (void *)&bs));
519566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5PopGroup(viewer));
523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5374d0cae8SMatthew G. Knepley }
5474d0cae8SMatthew G. Knepley 
5566976f2fSJacob Faibussowitsch static PetscErrorCode DMSwarmView_HDF5(DM dm, PetscViewer viewer)
56d71ae5a4SJacob Faibussowitsch {
5719307e5cSMatthew G. Knepley   DMSwarmCellDM celldm;
5874d0cae8SMatthew G. Knepley   Vec           coordinates;
5919307e5cSMatthew G. Knepley   PetscInt      Np, Nfc;
6074d0cae8SMatthew G. Knepley   PetscBool     isseq;
6119307e5cSMatthew G. Knepley   const char  **coordFields;
6274d0cae8SMatthew G. Knepley 
6374d0cae8SMatthew G. Knepley   PetscFunctionBegin;
6419307e5cSMatthew G. Knepley   PetscCall(DMSwarmGetCellDMActive(dm, &celldm));
6519307e5cSMatthew G. Knepley   PetscCall(DMSwarmCellDMGetCoordinateFields(celldm, &Nfc, &coordFields));
6619307e5cSMatthew G. Knepley   PetscCheck(Nfc == 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "We only support a single coordinate field right now, not %" PetscInt_FMT, Nfc);
679566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetSize(dm, &Np));
6819307e5cSMatthew G. Knepley   PetscCall(DMSwarmCreateGlobalVectorFromField(dm, coordFields[0], &coordinates));
699566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates"));
709566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5PushGroup(viewer, "/particles"));
719566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)coordinates, VECSEQ, &isseq));
729566063dSJacob Faibussowitsch   PetscCall(VecViewNative(coordinates, viewer));
739566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject)coordinates, "Np", PETSC_INT, (void *)&Np));
749566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5PopGroup(viewer));
7519307e5cSMatthew G. Knepley   PetscCall(DMSwarmDestroyGlobalVectorFromField(dm, coordFields[0], &coordinates));
763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7774d0cae8SMatthew G. Knepley }
7874d0cae8SMatthew G. Knepley #endif
7974d0cae8SMatthew G. Knepley 
8066976f2fSJacob Faibussowitsch static PetscErrorCode VecView_Swarm(Vec v, PetscViewer viewer)
81d71ae5a4SJacob Faibussowitsch {
8274d0cae8SMatthew G. Knepley   DM dm;
83f9558f5cSBarry Smith #if defined(PETSC_HAVE_HDF5)
8474d0cae8SMatthew G. Knepley   PetscBool ishdf5;
85f9558f5cSBarry Smith #endif
8674d0cae8SMatthew G. Knepley 
8774d0cae8SMatthew G. Knepley   PetscFunctionBegin;
889566063dSJacob Faibussowitsch   PetscCall(VecGetDM(v, &dm));
8928b400f6SJacob Faibussowitsch   PetscCheck(dm, PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_WRONG, "Vector not generated from a DM");
90f9558f5cSBarry Smith #if defined(PETSC_HAVE_HDF5)
919566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
9274d0cae8SMatthew G. Knepley   if (ishdf5) {
939566063dSJacob Faibussowitsch     PetscCall(VecView_Swarm_HDF5_Internal(v, viewer));
943ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
9574d0cae8SMatthew G. Knepley   }
96f9558f5cSBarry Smith #endif
979566063dSJacob Faibussowitsch   PetscCall(VecViewNative(v, viewer));
983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9974d0cae8SMatthew G. Knepley }
10074d0cae8SMatthew G. Knepley 
101d52c2f21SMatthew G. Knepley /*@C
102d52c2f21SMatthew G. Knepley   DMSwarmVectorGetField - Gets the fields from which to define a `Vec` object
1030bf7c1c5SMatthew G. Knepley   when `DMCreateLocalVector()`, or `DMCreateGlobalVector()` is called
1040bf7c1c5SMatthew G. Knepley 
1050bf7c1c5SMatthew G. Knepley   Not collective
1060bf7c1c5SMatthew G. Knepley 
1070bf7c1c5SMatthew G. Knepley   Input Parameter:
10819307e5cSMatthew G. Knepley . sw - a `DMSWARM`
1090bf7c1c5SMatthew G. Knepley 
110d52c2f21SMatthew G. Knepley   Output Parameters:
111d52c2f21SMatthew G. Knepley + Nf         - the number of fields
112d52c2f21SMatthew G. Knepley - fieldnames - the textual name given to each registered field, or NULL if it has not been set
1130bf7c1c5SMatthew G. Knepley 
1140bf7c1c5SMatthew G. Knepley   Level: beginner
1150bf7c1c5SMatthew G. Knepley 
116d52c2f21SMatthew G. Knepley .seealso: `DM`, `DMSWARM`, `DMSwarmVectorDefineField()`, `DMSwarmRegisterPetscDatatypeField()`, `DMCreateGlobalVector()`, `DMCreateLocalVector()`
1170bf7c1c5SMatthew G. Knepley @*/
11819307e5cSMatthew G. Knepley PetscErrorCode DMSwarmVectorGetField(DM sw, PetscInt *Nf, const char **fieldnames[])
1190bf7c1c5SMatthew G. Knepley {
12019307e5cSMatthew G. Knepley   DMSwarmCellDM celldm;
12119307e5cSMatthew G. Knepley 
1220bf7c1c5SMatthew G. Knepley   PetscFunctionBegin;
12319307e5cSMatthew G. Knepley   PetscValidHeaderSpecificType(sw, DM_CLASSID, 1, DMSWARM);
12419307e5cSMatthew G. Knepley   PetscCall(DMSwarmGetCellDMActive(sw, &celldm));
12519307e5cSMatthew G. Knepley   PetscCall(DMSwarmCellDMGetFields(celldm, Nf, fieldnames));
1260bf7c1c5SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
1270bf7c1c5SMatthew G. Knepley }
1280bf7c1c5SMatthew G. Knepley 
129cc4c1da9SBarry Smith /*@
13020f4b53cSBarry Smith   DMSwarmVectorDefineField - Sets the field from which to define a `Vec` object
13120f4b53cSBarry Smith   when `DMCreateLocalVector()`, or `DMCreateGlobalVector()` is called
13257795646SDave May 
13320f4b53cSBarry Smith   Collective
13457795646SDave May 
13560225df5SJacob Faibussowitsch   Input Parameters:
13620f4b53cSBarry Smith + dm        - a `DMSWARM`
137d52c2f21SMatthew G. Knepley - fieldname - the textual name given to each registered field
13857795646SDave May 
139d3a51819SDave May   Level: beginner
14057795646SDave May 
141d3a51819SDave May   Notes:
14220f4b53cSBarry Smith   The field with name `fieldname` must be defined as having a data type of `PetscScalar`.
143e7af74c8SDave May 
14420f4b53cSBarry Smith   This function must be called prior to calling `DMCreateLocalVector()`, `DMCreateGlobalVector()`.
14520f4b53cSBarry Smith   Multiple calls to `DMSwarmVectorDefineField()` are permitted.
146e7af74c8SDave May 
147d52c2f21SMatthew G. Knepley .seealso: `DM`, `DMSWARM`, `DMSwarmVectorDefineFields()`, `DMSwarmVectorGetField()`, `DMSwarmRegisterPetscDatatypeField()`, `DMCreateGlobalVector()`, `DMCreateLocalVector()`
148d3a51819SDave May @*/
149d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmVectorDefineField(DM dm, const char fieldname[])
150d71ae5a4SJacob Faibussowitsch {
151d52c2f21SMatthew G. Knepley   PetscFunctionBegin;
152d52c2f21SMatthew G. Knepley   PetscCall(DMSwarmVectorDefineFields(dm, 1, &fieldname));
153d52c2f21SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
154d52c2f21SMatthew G. Knepley }
155d52c2f21SMatthew G. Knepley 
156d52c2f21SMatthew G. Knepley /*@C
157d52c2f21SMatthew G. Knepley   DMSwarmVectorDefineFields - Sets the fields from which to define a `Vec` object
158d52c2f21SMatthew G. Knepley   when `DMCreateLocalVector()`, or `DMCreateGlobalVector()` is called
159d52c2f21SMatthew G. Knepley 
160d52c2f21SMatthew G. Knepley   Collective, No Fortran support
161d52c2f21SMatthew G. Knepley 
162d52c2f21SMatthew G. Knepley   Input Parameters:
16319307e5cSMatthew G. Knepley + sw         - a `DMSWARM`
164d52c2f21SMatthew G. Knepley . Nf         - the number of fields
165d52c2f21SMatthew G. Knepley - fieldnames - the textual name given to each registered field
166d52c2f21SMatthew G. Knepley 
167d52c2f21SMatthew G. Knepley   Level: beginner
168d52c2f21SMatthew G. Knepley 
169d52c2f21SMatthew G. Knepley   Notes:
170d52c2f21SMatthew G. Knepley   Each field with name in `fieldnames` must be defined as having a data type of `PetscScalar`.
171d52c2f21SMatthew G. Knepley 
172d52c2f21SMatthew G. Knepley   This function must be called prior to calling `DMCreateLocalVector()`, `DMCreateGlobalVector()`.
173d52c2f21SMatthew G. Knepley   Multiple calls to `DMSwarmVectorDefineField()` are permitted.
174d52c2f21SMatthew G. Knepley 
175d52c2f21SMatthew G. Knepley .seealso: `DM`, `DMSWARM`, `DMSwarmVectorDefineField()`, `DMSwarmVectorGetField()`, `DMSwarmRegisterPetscDatatypeField()`, `DMCreateGlobalVector()`, `DMCreateLocalVector()`
176d52c2f21SMatthew G. Knepley @*/
17719307e5cSMatthew G. Knepley PetscErrorCode DMSwarmVectorDefineFields(DM sw, PetscInt Nf, const char *fieldnames[])
178d52c2f21SMatthew G. Knepley {
17919307e5cSMatthew G. Knepley   DM_Swarm     *swarm = (DM_Swarm *)sw->data;
18019307e5cSMatthew G. Knepley   DMSwarmCellDM celldm;
181b5bcf523SDave May 
182a9cbaee5SMatthew G. Knepley   PetscFunctionBegin;
18319307e5cSMatthew G. Knepley   PetscValidHeaderSpecificType(sw, DM_CLASSID, 1, DMSWARM);
184d52c2f21SMatthew G. Knepley   if (fieldnames) PetscAssertPointer(fieldnames, 3);
18519307e5cSMatthew G. Knepley   if (!swarm->issetup) PetscCall(DMSetUp(sw));
18619307e5cSMatthew G. Knepley   PetscCheck(Nf >= 0, PetscObjectComm((PetscObject)sw), PETSC_ERR_ARG_OUTOFRANGE, "Number of fields must be non-negative, not %" PetscInt_FMT, Nf);
18719307e5cSMatthew G. Knepley   // Create a dummy cell DM if none has been specified (I think we should not support this mode)
18819307e5cSMatthew G. Knepley   if (!swarm->activeCellDM) {
18919307e5cSMatthew G. Knepley     DM            dm;
19019307e5cSMatthew G. Knepley     DMSwarmCellDM celldm;
191b5bcf523SDave May 
19219307e5cSMatthew G. Knepley     PetscCall(DMCreate(PetscObjectComm((PetscObject)sw), &dm));
19319307e5cSMatthew G. Knepley     PetscCall(DMSetType(dm, DMSHELL));
19419307e5cSMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)dm, "dummy"));
19519307e5cSMatthew G. Knepley     PetscCall(DMSwarmCellDMCreate(dm, 0, NULL, 0, NULL, &celldm));
19619307e5cSMatthew G. Knepley     PetscCall(DMDestroy(&dm));
19719307e5cSMatthew G. Knepley     PetscCall(DMSwarmAddCellDM(sw, celldm));
19819307e5cSMatthew G. Knepley     PetscCall(DMSwarmCellDMDestroy(&celldm));
19919307e5cSMatthew G. Knepley     PetscCall(DMSwarmSetCellDMActive(sw, "dummy"));
20019307e5cSMatthew G. Knepley   }
20119307e5cSMatthew G. Knepley   PetscCall(DMSwarmGetCellDMActive(sw, &celldm));
20219307e5cSMatthew G. Knepley   for (PetscInt f = 0; f < celldm->Nf; ++f) PetscCall(PetscFree(celldm->dmFields[f]));
20319307e5cSMatthew G. Knepley   PetscCall(PetscFree(celldm->dmFields));
20419307e5cSMatthew G. Knepley 
20519307e5cSMatthew G. Knepley   celldm->Nf = Nf;
20619307e5cSMatthew G. Knepley   PetscCall(PetscMalloc1(Nf, &celldm->dmFields));
207d52c2f21SMatthew G. Knepley   for (PetscInt f = 0; f < Nf; ++f) {
208d52c2f21SMatthew G. Knepley     PetscDataType type;
209d52c2f21SMatthew G. Knepley 
210d52c2f21SMatthew G. Knepley     // Check all fields are of type PETSC_REAL or PETSC_SCALAR
21119307e5cSMatthew G. Knepley     PetscCall(DMSwarmGetFieldInfo(sw, fieldnames[f], NULL, &type));
21219307e5cSMatthew G. Knepley     PetscCheck(type == PETSC_REAL, PetscObjectComm((PetscObject)sw), PETSC_ERR_SUP, "Only valid for PETSC_REAL");
21319307e5cSMatthew G. Knepley     PetscCall(PetscStrallocpy(fieldnames[f], (char **)&celldm->dmFields[f]));
214d52c2f21SMatthew G. Knepley   }
2153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
216b5bcf523SDave May }
217b5bcf523SDave May 
218cc651181SDave May /* requires DMSwarmDefineFieldVector has been called */
21919307e5cSMatthew G. Knepley static PetscErrorCode DMCreateGlobalVector_Swarm(DM sw, Vec *vec)
220d71ae5a4SJacob Faibussowitsch {
22119307e5cSMatthew G. Knepley   DM_Swarm     *swarm = (DM_Swarm *)sw->data;
22219307e5cSMatthew G. Knepley   DMSwarmCellDM celldm;
223b5bcf523SDave May   Vec           x;
224b5bcf523SDave May   char          name[PETSC_MAX_PATH_LEN];
22519307e5cSMatthew G. Knepley   PetscInt      bs = 0, n;
226b5bcf523SDave May 
227a9cbaee5SMatthew G. Knepley   PetscFunctionBegin;
22819307e5cSMatthew G. Knepley   if (!swarm->issetup) PetscCall(DMSetUp(sw));
22919307e5cSMatthew G. Knepley   PetscCall(DMSwarmGetCellDMActive(sw, &celldm));
23019307e5cSMatthew G. Knepley   PetscCheck(celldm->Nf, PetscObjectComm((PetscObject)sw), PETSC_ERR_USER, "Active cell DM does not define any fields");
23119307e5cSMatthew G. Knepley   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &n, NULL, NULL));
232cc651181SDave May 
233d52c2f21SMatthew G. Knepley   PetscCall(PetscStrncpy(name, "DMSwarmField", PETSC_MAX_PATH_LEN));
23419307e5cSMatthew G. Knepley   for (PetscInt f = 0; f < celldm->Nf; ++f) {
23519307e5cSMatthew G. Knepley     PetscInt fbs;
236d52c2f21SMatthew G. Knepley     PetscCall(PetscStrlcat(name, "_", PETSC_MAX_PATH_LEN));
23719307e5cSMatthew G. Knepley     PetscCall(PetscStrlcat(name, celldm->dmFields[f], PETSC_MAX_PATH_LEN));
23819307e5cSMatthew G. Knepley     PetscCall(DMSwarmGetFieldInfo(sw, celldm->dmFields[f], &fbs, NULL));
23919307e5cSMatthew G. Knepley     bs += fbs;
240d52c2f21SMatthew G. Knepley   }
24119307e5cSMatthew G. Knepley   PetscCall(VecCreate(PetscObjectComm((PetscObject)sw), &x));
2429566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)x, name));
24319307e5cSMatthew G. Knepley   PetscCall(VecSetSizes(x, n * bs, PETSC_DETERMINE));
24419307e5cSMatthew G. Knepley   PetscCall(VecSetBlockSize(x, bs));
24519307e5cSMatthew G. Knepley   PetscCall(VecSetDM(x, sw));
2469566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(x));
2479566063dSJacob Faibussowitsch   PetscCall(VecSetOperation(x, VECOP_VIEW, (void (*)(void))VecView_Swarm));
248b5bcf523SDave May   *vec = x;
2493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
250b5bcf523SDave May }
251b5bcf523SDave May 
252b5bcf523SDave May /* requires DMSwarmDefineFieldVector has been called */
25319307e5cSMatthew G. Knepley static PetscErrorCode DMCreateLocalVector_Swarm(DM sw, Vec *vec)
254d71ae5a4SJacob Faibussowitsch {
25519307e5cSMatthew G. Knepley   DM_Swarm     *swarm = (DM_Swarm *)sw->data;
25619307e5cSMatthew G. Knepley   DMSwarmCellDM celldm;
257b5bcf523SDave May   Vec           x;
258b5bcf523SDave May   char          name[PETSC_MAX_PATH_LEN];
25919307e5cSMatthew G. Knepley   PetscInt      bs = 0, n;
260b5bcf523SDave May 
261a9cbaee5SMatthew G. Knepley   PetscFunctionBegin;
26219307e5cSMatthew G. Knepley   if (!swarm->issetup) PetscCall(DMSetUp(sw));
26319307e5cSMatthew G. Knepley   PetscCall(DMSwarmGetCellDMActive(sw, &celldm));
26419307e5cSMatthew G. Knepley   PetscCheck(celldm->Nf, PetscObjectComm((PetscObject)sw), PETSC_ERR_USER, "Active cell DM does not define any fields");
26519307e5cSMatthew G. Knepley   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &n, NULL, NULL));
266cc651181SDave May 
267d52c2f21SMatthew G. Knepley   PetscCall(PetscStrncpy(name, "DMSwarmField", PETSC_MAX_PATH_LEN));
26819307e5cSMatthew G. Knepley   for (PetscInt f = 0; f < celldm->Nf; ++f) {
26919307e5cSMatthew G. Knepley     PetscInt fbs;
270d52c2f21SMatthew G. Knepley     PetscCall(PetscStrlcat(name, "_", PETSC_MAX_PATH_LEN));
27119307e5cSMatthew G. Knepley     PetscCall(PetscStrlcat(name, celldm->dmFields[f], PETSC_MAX_PATH_LEN));
27219307e5cSMatthew G. Knepley     PetscCall(DMSwarmGetFieldInfo(sw, celldm->dmFields[f], &fbs, NULL));
27319307e5cSMatthew G. Knepley     bs += fbs;
274d52c2f21SMatthew G. Knepley   }
2759566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_SELF, &x));
2769566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)x, name));
27719307e5cSMatthew G. Knepley   PetscCall(VecSetSizes(x, n * bs, PETSC_DETERMINE));
27819307e5cSMatthew G. Knepley   PetscCall(VecSetBlockSize(x, bs));
27919307e5cSMatthew G. Knepley   PetscCall(VecSetDM(x, sw));
2809566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(x));
281b5bcf523SDave May   *vec = x;
2823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
283b5bcf523SDave May }
284b5bcf523SDave May 
285d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSwarmDestroyVectorFromField_Private(DM dm, const char fieldname[], Vec *vec)
286d71ae5a4SJacob Faibussowitsch {
287fb1bcc12SMatthew G. Knepley   DM_Swarm        *swarm = (DM_Swarm *)dm->data;
28877048351SPatrick Sanan   DMSwarmDataField gfield;
2892e956fe4SStefano Zampini   PetscInt         bs, nlocal, fid = -1, cfid = -2;
2902e956fe4SStefano Zampini   PetscBool        flg;
291d3a51819SDave May 
292fb1bcc12SMatthew G. Knepley   PetscFunctionBegin;
2932e956fe4SStefano Zampini   /* check vector is an inplace array */
2942e956fe4SStefano Zampini   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldIdByName(swarm->db, fieldname, &fid));
2952e956fe4SStefano Zampini   PetscCall(PetscObjectComposedDataGetInt((PetscObject)*vec, SwarmDataFieldId, cfid, flg));
296ea17275aSJose E. Roman   (void)flg; /* avoid compiler warning */
297e978a55eSPierre 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);
2989566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(*vec, &nlocal));
2999566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(*vec, &bs));
30008401ef6SPierre 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");
3019566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield));
3029566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataFieldRestoreAccess(gfield));
3038df78e51SMark Adams   PetscCall(VecResetArray(*vec));
3049566063dSJacob Faibussowitsch   PetscCall(VecDestroy(vec));
3053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
306fb1bcc12SMatthew G. Knepley }
307fb1bcc12SMatthew G. Knepley 
308d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSwarmCreateVectorFromField_Private(DM dm, const char fieldname[], MPI_Comm comm, Vec *vec)
309d71ae5a4SJacob Faibussowitsch {
310fb1bcc12SMatthew G. Knepley   DM_Swarm     *swarm = (DM_Swarm *)dm->data;
311fb1bcc12SMatthew G. Knepley   PetscDataType type;
312fb1bcc12SMatthew G. Knepley   PetscScalar  *array;
3132e956fe4SStefano Zampini   PetscInt      bs, n, fid;
314fb1bcc12SMatthew G. Knepley   char          name[PETSC_MAX_PATH_LEN];
315e4fbd051SBarry Smith   PetscMPIInt   size;
3167f92dde0SRylanor   PetscBool     iscuda, iskokkos, iship;
317fb1bcc12SMatthew G. Knepley 
318fb1bcc12SMatthew G. Knepley   PetscFunctionBegin;
3199566063dSJacob Faibussowitsch   if (!swarm->issetup) PetscCall(DMSetUp(dm));
3209566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &n, NULL, NULL));
3219566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, fieldname, &bs, &type, (void **)&array));
32208401ef6SPierre Jolivet   PetscCheck(type == PETSC_REAL, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only valid for PETSC_REAL");
323fb1bcc12SMatthew G. Knepley 
3249566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
3258df78e51SMark Adams   PetscCall(PetscStrcmp(dm->vectype, VECKOKKOS, &iskokkos));
3268df78e51SMark Adams   PetscCall(PetscStrcmp(dm->vectype, VECCUDA, &iscuda));
3277f92dde0SRylanor   PetscCall(PetscStrcmp(dm->vectype, VECHIP, &iship));
3288df78e51SMark Adams   PetscCall(VecCreate(comm, vec));
3298df78e51SMark Adams   PetscCall(VecSetSizes(*vec, n * bs, PETSC_DETERMINE));
3308df78e51SMark Adams   PetscCall(VecSetBlockSize(*vec, bs));
3318df78e51SMark Adams   if (iskokkos) PetscCall(VecSetType(*vec, VECKOKKOS));
3328df78e51SMark Adams   else if (iscuda) PetscCall(VecSetType(*vec, VECCUDA));
3337f92dde0SRylanor   else if (iship) PetscCall(VecSetType(*vec, VECHIP));
3348df78e51SMark Adams   else PetscCall(VecSetType(*vec, VECSTANDARD));
3358df78e51SMark Adams   PetscCall(VecPlaceArray(*vec, array));
3368df78e51SMark Adams 
3379566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN - 1, "DMSwarmSharedField_%s", fieldname));
3389566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)*vec, name));
339fb1bcc12SMatthew G. Knepley 
340fb1bcc12SMatthew G. Knepley   /* Set guard */
3412e956fe4SStefano Zampini   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldIdByName(swarm->db, fieldname, &fid));
3422e956fe4SStefano Zampini   PetscCall(PetscObjectComposedDataSetInt((PetscObject)*vec, SwarmDataFieldId, fid));
34374d0cae8SMatthew G. Knepley 
3449566063dSJacob Faibussowitsch   PetscCall(VecSetDM(*vec, dm));
3459566063dSJacob Faibussowitsch   PetscCall(VecSetOperation(*vec, VECOP_VIEW, (void (*)(void))VecView_Swarm));
3463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
347fb1bcc12SMatthew G. Knepley }
348fb1bcc12SMatthew G. Knepley 
34919307e5cSMatthew G. Knepley static PetscErrorCode DMSwarmDestroyVectorFromFields_Private(DM sw, PetscInt Nf, const char *fieldnames[], Vec *vec)
35019307e5cSMatthew G. Knepley {
35119307e5cSMatthew G. Knepley   DM_Swarm          *swarm = (DM_Swarm *)sw->data;
35219307e5cSMatthew G. Knepley   const PetscScalar *array;
35319307e5cSMatthew G. Knepley   PetscInt           bs, n, id = 0, cid = -2;
35419307e5cSMatthew G. Knepley   PetscBool          flg;
35519307e5cSMatthew G. Knepley 
35619307e5cSMatthew G. Knepley   PetscFunctionBegin;
35719307e5cSMatthew G. Knepley   for (PetscInt f = 0; f < Nf; ++f) {
35819307e5cSMatthew G. Knepley     PetscInt fid;
35919307e5cSMatthew G. Knepley 
36019307e5cSMatthew G. Knepley     PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldIdByName(swarm->db, fieldnames[f], &fid));
36119307e5cSMatthew G. Knepley     id += fid;
36219307e5cSMatthew G. Knepley   }
36319307e5cSMatthew G. Knepley   PetscCall(PetscObjectComposedDataGetInt((PetscObject)*vec, SwarmDataFieldId, cid, flg));
36419307e5cSMatthew G. Knepley   (void)flg; /* avoid compiler warning */
36519307e5cSMatthew G. Knepley   PetscCheck(cid == id, PetscObjectComm((PetscObject)sw), PETSC_ERR_USER, "Vector being destroyed was not created from DMSwarm field(%s)! %" PetscInt_FMT " != %" PetscInt_FMT, fieldnames[0], cid, id);
36619307e5cSMatthew G. Knepley   PetscCall(VecGetLocalSize(*vec, &n));
36719307e5cSMatthew G. Knepley   PetscCall(VecGetBlockSize(*vec, &bs));
36819307e5cSMatthew G. Knepley   n /= bs;
36919307e5cSMatthew G. Knepley   PetscCheck(n == swarm->db->L, PetscObjectComm((PetscObject)sw), PETSC_ERR_USER, "DMSwarm sizes have changed since vector was created - cannot ensure pointers are valid");
37019307e5cSMatthew G. Knepley   PetscCall(VecGetArrayRead(*vec, &array));
37119307e5cSMatthew G. Knepley   for (PetscInt f = 0, off = 0; f < Nf; ++f) {
37219307e5cSMatthew G. Knepley     PetscScalar  *farray;
37319307e5cSMatthew G. Knepley     PetscDataType ftype;
37419307e5cSMatthew G. Knepley     PetscInt      fbs;
37519307e5cSMatthew G. Knepley 
37619307e5cSMatthew G. Knepley     PetscCall(DMSwarmGetField(sw, fieldnames[f], &fbs, &ftype, (void **)&farray));
37719307e5cSMatthew G. Knepley     PetscCheck(off + fbs <= bs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid blocksize %" PetscInt_FMT " < %" PetscInt_FMT, bs, off + fbs);
37819307e5cSMatthew G. Knepley     for (PetscInt i = 0; i < n; ++i) {
37919307e5cSMatthew G. Knepley       for (PetscInt b = 0; b < fbs; ++b) farray[i * fbs + b] = array[i * bs + off + b];
38019307e5cSMatthew G. Knepley     }
38119307e5cSMatthew G. Knepley     off += fbs;
38219307e5cSMatthew G. Knepley     PetscCall(DMSwarmRestoreField(sw, fieldnames[f], &fbs, &ftype, (void **)&farray));
38319307e5cSMatthew G. Knepley   }
38419307e5cSMatthew G. Knepley   PetscCall(VecRestoreArrayRead(*vec, &array));
38519307e5cSMatthew G. Knepley   PetscCall(VecDestroy(vec));
38619307e5cSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
38719307e5cSMatthew G. Knepley }
38819307e5cSMatthew G. Knepley 
38919307e5cSMatthew G. Knepley static PetscErrorCode DMSwarmCreateVectorFromFields_Private(DM sw, PetscInt Nf, const char *fieldnames[], MPI_Comm comm, Vec *vec)
39019307e5cSMatthew G. Knepley {
39119307e5cSMatthew G. Knepley   DM_Swarm    *swarm = (DM_Swarm *)sw->data;
39219307e5cSMatthew G. Knepley   PetscScalar *array;
39319307e5cSMatthew G. Knepley   PetscInt     n, bs = 0, id = 0;
39419307e5cSMatthew G. Knepley   char         name[PETSC_MAX_PATH_LEN];
39519307e5cSMatthew G. Knepley 
39619307e5cSMatthew G. Knepley   PetscFunctionBegin;
39719307e5cSMatthew G. Knepley   if (!swarm->issetup) PetscCall(DMSetUp(sw));
39819307e5cSMatthew G. Knepley   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &n, NULL, NULL));
39919307e5cSMatthew G. Knepley   for (PetscInt f = 0; f < Nf; ++f) {
40019307e5cSMatthew G. Knepley     PetscDataType ftype;
40119307e5cSMatthew G. Knepley     PetscInt      fbs;
40219307e5cSMatthew G. Knepley 
40319307e5cSMatthew G. Knepley     PetscCall(DMSwarmGetFieldInfo(sw, fieldnames[f], &fbs, &ftype));
40419307e5cSMatthew G. Knepley     PetscCheck(ftype == PETSC_REAL, PetscObjectComm((PetscObject)sw), PETSC_ERR_SUP, "Only valid for PETSC_REAL");
40519307e5cSMatthew G. Knepley     bs += fbs;
40619307e5cSMatthew G. Knepley   }
40719307e5cSMatthew G. Knepley 
40819307e5cSMatthew G. Knepley   PetscCall(VecCreate(comm, vec));
40919307e5cSMatthew G. Knepley   PetscCall(VecSetSizes(*vec, n * bs, PETSC_DETERMINE));
41019307e5cSMatthew G. Knepley   PetscCall(VecSetBlockSize(*vec, bs));
41119307e5cSMatthew G. Knepley   PetscCall(VecSetType(*vec, sw->vectype));
41219307e5cSMatthew G. Knepley 
41319307e5cSMatthew G. Knepley   PetscCall(VecGetArrayWrite(*vec, &array));
41419307e5cSMatthew G. Knepley   for (PetscInt f = 0, off = 0; f < Nf; ++f) {
41519307e5cSMatthew G. Knepley     PetscScalar  *farray;
41619307e5cSMatthew G. Knepley     PetscDataType ftype;
41719307e5cSMatthew G. Knepley     PetscInt      fbs;
41819307e5cSMatthew G. Knepley 
41919307e5cSMatthew G. Knepley     PetscCall(DMSwarmGetField(sw, fieldnames[f], &fbs, &ftype, (void **)&farray));
42019307e5cSMatthew G. Knepley     for (PetscInt i = 0; i < n; ++i) {
42119307e5cSMatthew G. Knepley       for (PetscInt b = 0; b < fbs; ++b) array[i * bs + off + b] = farray[i * fbs + b];
42219307e5cSMatthew G. Knepley     }
42319307e5cSMatthew G. Knepley     off += fbs;
42419307e5cSMatthew G. Knepley     PetscCall(DMSwarmRestoreField(sw, fieldnames[f], &fbs, &ftype, (void **)&farray));
42519307e5cSMatthew G. Knepley   }
42619307e5cSMatthew G. Knepley   PetscCall(VecRestoreArrayWrite(*vec, &array));
42719307e5cSMatthew G. Knepley 
42819307e5cSMatthew G. Knepley   PetscCall(PetscStrncpy(name, "DMSwarmField", PETSC_MAX_PATH_LEN));
42919307e5cSMatthew G. Knepley   for (PetscInt f = 0; f < Nf; ++f) {
43019307e5cSMatthew G. Knepley     PetscCall(PetscStrlcat(name, "_", PETSC_MAX_PATH_LEN));
43119307e5cSMatthew G. Knepley     PetscCall(PetscStrlcat(name, fieldnames[f], PETSC_MAX_PATH_LEN));
43219307e5cSMatthew G. Knepley   }
43319307e5cSMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)*vec, name));
43419307e5cSMatthew G. Knepley 
43519307e5cSMatthew G. Knepley   for (PetscInt f = 0; f < Nf; ++f) {
43619307e5cSMatthew G. Knepley     PetscInt fid;
43719307e5cSMatthew G. Knepley 
43819307e5cSMatthew G. Knepley     PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldIdByName(swarm->db, fieldnames[f], &fid));
43919307e5cSMatthew G. Knepley     id += fid;
44019307e5cSMatthew G. Knepley   }
44119307e5cSMatthew G. Knepley   PetscCall(PetscObjectComposedDataSetInt((PetscObject)*vec, SwarmDataFieldId, id));
44219307e5cSMatthew G. Knepley 
44319307e5cSMatthew G. Knepley   PetscCall(VecSetDM(*vec, sw));
44419307e5cSMatthew G. Knepley   PetscCall(VecSetOperation(*vec, VECOP_VIEW, (void (*)(void))VecView_Swarm));
44519307e5cSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
44619307e5cSMatthew G. Knepley }
44719307e5cSMatthew G. Knepley 
4480643ed39SMatthew G. Knepley /* This creates a "mass matrix" between a finite element and particle space. If a finite element interpolant is given by
4490643ed39SMatthew G. Knepley 
4500643ed39SMatthew G. Knepley      \hat f = \sum_i f_i \phi_i
4510643ed39SMatthew G. Knepley 
4520643ed39SMatthew G. Knepley    and a particle function is given by
4530643ed39SMatthew G. Knepley 
4540643ed39SMatthew G. Knepley      f = \sum_i w_i \delta(x - x_i)
4550643ed39SMatthew G. Knepley 
4560643ed39SMatthew G. Knepley    then we want to require that
4570643ed39SMatthew G. Knepley 
4580643ed39SMatthew G. Knepley      M \hat f = M_p f
4590643ed39SMatthew G. Knepley 
4600643ed39SMatthew G. Knepley    where the particle mass matrix is given by
4610643ed39SMatthew G. Knepley 
4620643ed39SMatthew G. Knepley      (M_p)_{ij} = \int \phi_i \delta(x - x_j)
4630643ed39SMatthew G. Knepley 
4640643ed39SMatthew G. Knepley    The way Dave May does particles, they amount to quadratue weights rather than delta functions, so he has |J| is in
4650643ed39SMatthew G. Knepley    his integral. We allow this with the boolean flag.
4660643ed39SMatthew G. Knepley */
467d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSwarmComputeMassMatrix_Private(DM dmc, DM dmf, Mat mass, PetscBool useDeltaFunction, void *ctx)
468d71ae5a4SJacob Faibussowitsch {
46983c47955SMatthew G. Knepley   const char   *name = "Mass Matrix";
4700643ed39SMatthew G. Knepley   MPI_Comm      comm;
47119307e5cSMatthew G. Knepley   DMSwarmCellDM celldm;
47283c47955SMatthew G. Knepley   PetscDS       prob;
47383c47955SMatthew G. Knepley   PetscSection  fsection, globalFSection;
474e8f14785SLisandro Dalcin   PetscHSetIJ   ht;
4750643ed39SMatthew G. Knepley   PetscLayout   rLayout, colLayout;
47683c47955SMatthew G. Knepley   PetscInt     *dnz, *onz;
477adb2528bSMark Adams   PetscInt      locRows, locCols, rStart, colStart, colEnd, *rowIDXs;
4780643ed39SMatthew G. Knepley   PetscReal    *xi, *v0, *J, *invJ, detJ = 1.0, v0ref[3] = {-1.0, -1.0, -1.0};
47983c47955SMatthew G. Knepley   PetscScalar  *elemMat;
48019307e5cSMatthew G. Knepley   PetscInt      dim, Nf, Nfc, cStart, cEnd, totDim, maxC = 0, totNc = 0;
48119307e5cSMatthew G. Knepley   const char  **coordFields;
48219307e5cSMatthew G. Knepley   PetscReal   **coordVals;
48319307e5cSMatthew G. Knepley   PetscInt     *bs;
48483c47955SMatthew G. Knepley 
48583c47955SMatthew G. Knepley   PetscFunctionBegin;
4869566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)mass, &comm));
4879566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dmf, &dim));
4889566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dmf, &prob));
4899566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(prob, &Nf));
4909566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(prob, &totDim));
4919566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(dim, &v0, dim * dim, &J, dim * dim, &invJ));
4929566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dmf, &fsection));
4939566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalSection(dmf, &globalFSection));
4949566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd));
4959566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(mass, &locRows, &locCols));
49683c47955SMatthew G. Knepley 
49719307e5cSMatthew G. Knepley   PetscCall(DMSwarmGetCellDMActive(dmc, &celldm));
49819307e5cSMatthew G. Knepley   PetscCall(DMSwarmCellDMGetCoordinateFields(celldm, &Nfc, &coordFields));
49919307e5cSMatthew G. Knepley   PetscCall(PetscMalloc2(Nfc, &coordVals, Nfc, &bs));
500d52c2f21SMatthew G. Knepley 
5019566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreate(comm, &colLayout));
5029566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetLocalSize(colLayout, locCols));
5039566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(colLayout, 1));
5049566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(colLayout));
5059566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetRange(colLayout, &colStart, &colEnd));
5069566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&colLayout));
5070643ed39SMatthew G. Knepley 
5089566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreate(comm, &rLayout));
5099566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetLocalSize(rLayout, locRows));
5109566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(rLayout, 1));
5119566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(rLayout));
5129566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetRange(rLayout, &rStart, NULL));
5139566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&rLayout));
5140643ed39SMatthew G. Knepley 
5159566063dSJacob Faibussowitsch   PetscCall(PetscCalloc2(locRows, &dnz, locRows, &onz));
5169566063dSJacob Faibussowitsch   PetscCall(PetscHSetIJCreate(&ht));
51753e60ab4SJoseph Pusztay 
5189566063dSJacob Faibussowitsch   PetscCall(PetscSynchronizedFlush(comm, NULL));
51919307e5cSMatthew G. Knepley   for (PetscInt field = 0; field < Nf; ++field) {
5200bf7c1c5SMatthew G. Knepley     PetscObject  obj;
5210bf7c1c5SMatthew G. Knepley     PetscClassId id;
5220bf7c1c5SMatthew G. Knepley     PetscInt     Nc;
5230bf7c1c5SMatthew G. Knepley 
5240bf7c1c5SMatthew G. Knepley     PetscCall(PetscDSGetDiscretization(prob, field, &obj));
5250bf7c1c5SMatthew G. Knepley     PetscCall(PetscObjectGetClassId(obj, &id));
5260bf7c1c5SMatthew G. Knepley     if (id == PETSCFE_CLASSID) PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc));
5270bf7c1c5SMatthew G. Knepley     else PetscCall(PetscFVGetNumComponents((PetscFV)obj, &Nc));
5280bf7c1c5SMatthew G. Knepley     totNc += Nc;
5290bf7c1c5SMatthew G. Knepley   }
5300643ed39SMatthew G. Knepley   /* count non-zeros */
5319566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortGetAccess(dmc));
53219307e5cSMatthew G. Knepley   for (PetscInt field = 0; field < Nf; ++field) {
5330bf7c1c5SMatthew G. Knepley     PetscObject  obj;
5340bf7c1c5SMatthew G. Knepley     PetscClassId id;
5350bf7c1c5SMatthew G. Knepley     PetscInt     Nc;
5360bf7c1c5SMatthew G. Knepley 
5370bf7c1c5SMatthew G. Knepley     PetscCall(PetscDSGetDiscretization(prob, field, &obj));
5380bf7c1c5SMatthew G. Knepley     PetscCall(PetscObjectGetClassId(obj, &id));
5390bf7c1c5SMatthew G. Knepley     if (id == PETSCFE_CLASSID) PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc));
5400bf7c1c5SMatthew G. Knepley     else PetscCall(PetscFVGetNumComponents((PetscFV)obj, &Nc));
5410bf7c1c5SMatthew G. Knepley 
54219307e5cSMatthew G. Knepley     for (PetscInt cell = cStart; cell < cEnd; ++cell) {
5430643ed39SMatthew G. Knepley       PetscInt *findices, *cindices; /* fine is vertices, coarse is particles */
54483c47955SMatthew G. Knepley       PetscInt  numFIndices, numCIndices;
54583c47955SMatthew G. Knepley 
5469566063dSJacob Faibussowitsch       PetscCall(DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
5479566063dSJacob Faibussowitsch       PetscCall(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices));
548fc7c92abSMatthew G. Knepley       maxC = PetscMax(maxC, numCIndices);
54983c47955SMatthew G. Knepley       {
550e8f14785SLisandro Dalcin         PetscHashIJKey key;
551e8f14785SLisandro Dalcin         PetscBool      missing;
5520bf7c1c5SMatthew G. Knepley         for (PetscInt i = 0; i < numFIndices; ++i) {
553adb2528bSMark Adams           key.j = findices[i]; /* global column (from Plex) */
554adb2528bSMark Adams           if (key.j >= 0) {
55583c47955SMatthew G. Knepley             /* Get indices for coarse elements */
5560bf7c1c5SMatthew G. Knepley             for (PetscInt j = 0; j < numCIndices; ++j) {
5570bf7c1c5SMatthew G. Knepley               for (PetscInt c = 0; c < Nc; ++c) {
5580bf7c1c5SMatthew G. Knepley                 // TODO Need field offset on particle here
5590bf7c1c5SMatthew G. Knepley                 key.i = cindices[j] * totNc + c + rStart; /* global cols (from Swarm) */
560adb2528bSMark Adams                 if (key.i < 0) continue;
5619566063dSJacob Faibussowitsch                 PetscCall(PetscHSetIJQueryAdd(ht, key, &missing));
56283c47955SMatthew G. Knepley                 if (missing) {
5630643ed39SMatthew G. Knepley                   if ((key.j >= colStart) && (key.j < colEnd)) ++dnz[key.i - rStart];
564e8f14785SLisandro Dalcin                   else ++onz[key.i - rStart];
56563a3b9bcSJacob Faibussowitsch                 } else SETERRQ(PetscObjectComm((PetscObject)dmf), PETSC_ERR_SUP, "Set new value at %" PetscInt_FMT ",%" PetscInt_FMT, key.i, key.j);
56683c47955SMatthew G. Knepley               }
567fc7c92abSMatthew G. Knepley             }
568fc7c92abSMatthew G. Knepley           }
5690bf7c1c5SMatthew G. Knepley         }
570fe11354eSMatthew G. Knepley         PetscCall(DMSwarmSortRestorePointsPerCell(dmc, cell, &numCIndices, &cindices));
57183c47955SMatthew G. Knepley       }
5729566063dSJacob Faibussowitsch       PetscCall(DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
57383c47955SMatthew G. Knepley     }
57483c47955SMatthew G. Knepley   }
5759566063dSJacob Faibussowitsch   PetscCall(PetscHSetIJDestroy(&ht));
5769566063dSJacob Faibussowitsch   PetscCall(MatXAIJSetPreallocation(mass, 1, dnz, onz, NULL, NULL));
5779566063dSJacob Faibussowitsch   PetscCall(MatSetOption(mass, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
5789566063dSJacob Faibussowitsch   PetscCall(PetscFree2(dnz, onz));
5790bf7c1c5SMatthew G. Knepley   PetscCall(PetscMalloc3(maxC * totNc * totDim, &elemMat, maxC * totNc, &rowIDXs, maxC * dim, &xi));
58019307e5cSMatthew G. Knepley   for (PetscInt field = 0; field < Nf; ++field) {
581ef0bb6c7SMatthew G. Knepley     PetscTabulation Tcoarse;
58283c47955SMatthew G. Knepley     PetscObject     obj;
583ad9634fcSMatthew G. Knepley     PetscClassId    id;
58419307e5cSMatthew G. Knepley     PetscInt        Nc;
58583c47955SMatthew G. Knepley 
5869566063dSJacob Faibussowitsch     PetscCall(PetscDSGetDiscretization(prob, field, &obj));
587ad9634fcSMatthew G. Knepley     PetscCall(PetscObjectGetClassId(obj, &id));
588ad9634fcSMatthew G. Knepley     if (id == PETSCFE_CLASSID) PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc));
589ad9634fcSMatthew G. Knepley     else PetscCall(PetscFVGetNumComponents((PetscFV)obj, &Nc));
590ea7032a0SMatthew G. Knepley 
59119307e5cSMatthew G. Knepley     for (PetscInt i = 0; i < Nfc; ++i) PetscCall(DMSwarmGetField(dmc, coordFields[i], &bs[i], NULL, (void **)&coordVals[i]));
59219307e5cSMatthew G. Knepley     for (PetscInt cell = cStart; cell < cEnd; ++cell) {
59383c47955SMatthew G. Knepley       PetscInt *findices, *cindices;
59483c47955SMatthew G. Knepley       PetscInt  numFIndices, numCIndices;
59583c47955SMatthew G. Knepley 
5960643ed39SMatthew G. Knepley       /* TODO: Use DMField instead of assuming affine */
5979566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ));
5989566063dSJacob Faibussowitsch       PetscCall(DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
5999566063dSJacob Faibussowitsch       PetscCall(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices));
60019307e5cSMatthew G. Knepley       for (PetscInt j = 0; j < numCIndices; ++j) {
60119307e5cSMatthew G. Knepley         PetscReal xr[8];
60219307e5cSMatthew G. Knepley         PetscInt  off = 0;
60319307e5cSMatthew G. Knepley 
60419307e5cSMatthew G. Knepley         for (PetscInt i = 0; i < Nfc; ++i) {
60519307e5cSMatthew G. Knepley           for (PetscInt b = 0; b < bs[i]; ++b, ++off) xr[off] = coordVals[i][cindices[j] * bs[i] + b];
60619307e5cSMatthew G. Knepley         }
60719307e5cSMatthew G. Knepley         PetscCheck(off == dim, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The total block size of coordinates is %" PetscInt_FMT " != %" PetscInt_FMT " the DM coordinate dimension", off, dim);
60819307e5cSMatthew G. Knepley         CoordinatesRealToRef(dim, dim, v0ref, v0, invJ, xr, &xi[j * dim]);
60919307e5cSMatthew G. Knepley       }
610ad9634fcSMatthew G. Knepley       if (id == PETSCFE_CLASSID) PetscCall(PetscFECreateTabulation((PetscFE)obj, 1, numCIndices, xi, 0, &Tcoarse));
611ad9634fcSMatthew G. Knepley       else PetscCall(PetscFVCreateTabulation((PetscFV)obj, 1, numCIndices, xi, 0, &Tcoarse));
61283c47955SMatthew G. Knepley       /* Get elemMat entries by multiplying by weight */
6130bf7c1c5SMatthew G. Knepley       PetscCall(PetscArrayzero(elemMat, numCIndices * Nc * totDim));
6140bf7c1c5SMatthew G. Knepley       for (PetscInt i = 0; i < numFIndices / Nc; ++i) {
6150bf7c1c5SMatthew G. Knepley         for (PetscInt j = 0; j < numCIndices; ++j) {
6160bf7c1c5SMatthew G. Knepley           for (PetscInt c = 0; c < Nc; ++c) {
6170bf7c1c5SMatthew G. Knepley             // TODO Need field offset on particle and field here
6180643ed39SMatthew G. Knepley             /* B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c */
6190bf7c1c5SMatthew G. Knepley             elemMat[(j * totNc + c) * numFIndices + i * Nc + c] += Tcoarse->T[0][(j * numFIndices + i * Nc + c) * Nc + c] * (useDeltaFunction ? 1.0 : detJ);
62083c47955SMatthew G. Knepley           }
6210643ed39SMatthew G. Knepley         }
6220643ed39SMatthew G. Knepley       }
6230bf7c1c5SMatthew G. Knepley       for (PetscInt j = 0; j < numCIndices; ++j)
6240bf7c1c5SMatthew G. Knepley         // TODO Need field offset on particle here
6250bf7c1c5SMatthew G. Knepley         for (PetscInt c = 0; c < Nc; ++c) rowIDXs[j * Nc + c] = cindices[j] * totNc + c + rStart;
6260bf7c1c5SMatthew G. Knepley       if (0) PetscCall(DMPrintCellMatrix(cell, name, numCIndices * Nc, numFIndices, elemMat));
6270bf7c1c5SMatthew G. Knepley       PetscCall(MatSetValues(mass, numCIndices * Nc, rowIDXs, numFIndices, findices, elemMat, ADD_VALUES));
628fe11354eSMatthew G. Knepley       PetscCall(DMSwarmSortRestorePointsPerCell(dmc, cell, &numCIndices, &cindices));
6299566063dSJacob Faibussowitsch       PetscCall(DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
6309566063dSJacob Faibussowitsch       PetscCall(PetscTabulationDestroy(&Tcoarse));
63183c47955SMatthew G. Knepley     }
63219307e5cSMatthew G. Knepley     for (PetscInt i = 0; i < Nfc; ++i) PetscCall(DMSwarmRestoreField(dmc, coordFields[i], &bs[i], NULL, (void **)&coordVals[i]));
63383c47955SMatthew G. Knepley   }
6349566063dSJacob Faibussowitsch   PetscCall(PetscFree3(elemMat, rowIDXs, xi));
6359566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortRestoreAccess(dmc));
6369566063dSJacob Faibussowitsch   PetscCall(PetscFree3(v0, J, invJ));
63719307e5cSMatthew G. Knepley   PetscCall(PetscFree2(coordVals, bs));
6389566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mass, MAT_FINAL_ASSEMBLY));
6399566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mass, MAT_FINAL_ASSEMBLY));
6403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
64183c47955SMatthew G. Knepley }
64283c47955SMatthew G. Knepley 
643d0c080abSJoseph Pusztay /* Returns empty matrix for use with SNES FD */
644d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCreateMatrix_Swarm(DM sw, Mat *m)
645d71ae5a4SJacob Faibussowitsch {
646d0c080abSJoseph Pusztay   Vec      field;
647d0c080abSJoseph Pusztay   PetscInt size;
648d0c080abSJoseph Pusztay 
649d0c080abSJoseph Pusztay   PetscFunctionBegin;
6509566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(sw, &field));
6519566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(field, &size));
6529566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(sw, &field));
6539566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, m));
6549566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(*m));
6559566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*m, PETSC_DECIDE, PETSC_DECIDE, size, size));
6569566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(*m, 1, NULL));
6579566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(*m));
6589566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(*m, MAT_FINAL_ASSEMBLY));
6599566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(*m, MAT_FINAL_ASSEMBLY));
6609566063dSJacob Faibussowitsch   PetscCall(MatShift(*m, 1.0));
6619566063dSJacob Faibussowitsch   PetscCall(MatSetDM(*m, sw));
6623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
663d0c080abSJoseph Pusztay }
664d0c080abSJoseph Pusztay 
665adb2528bSMark Adams /* FEM cols, Particle rows */
666d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCreateMassMatrix_Swarm(DM dmCoarse, DM dmFine, Mat *mass)
667d71ae5a4SJacob Faibussowitsch {
66819307e5cSMatthew G. Knepley   DMSwarmCellDM celldm;
669895a1698SMatthew G. Knepley   PetscSection  gsf;
67019307e5cSMatthew G. Knepley   PetscInt      m, n, Np, bs;
67183c47955SMatthew G. Knepley   void         *ctx;
67283c47955SMatthew G. Knepley 
67383c47955SMatthew G. Knepley   PetscFunctionBegin;
67419307e5cSMatthew G. Knepley   PetscCall(DMSwarmGetCellDMActive(dmCoarse, &celldm));
67519307e5cSMatthew G. Knepley   PetscCheck(celldm->Nf, PetscObjectComm((PetscObject)dmCoarse), PETSC_ERR_USER, "Active cell DM does not define any fields");
6769566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalSection(dmFine, &gsf));
6779566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetConstrainedStorageSize(gsf, &m));
6780bf7c1c5SMatthew G. Knepley   PetscCall(DMSwarmGetLocalSize(dmCoarse, &Np));
67919307e5cSMatthew G. Knepley   PetscCall(DMSwarmCellDMGetBlockSize(celldm, dmCoarse, &bs));
68019307e5cSMatthew G. Knepley   n = Np * bs;
6819566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)dmCoarse), mass));
6829566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*mass, n, m, PETSC_DETERMINE, PETSC_DETERMINE));
6839566063dSJacob Faibussowitsch   PetscCall(MatSetType(*mass, dmCoarse->mattype));
6849566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(dmFine, &ctx));
68583c47955SMatthew G. Knepley 
6869566063dSJacob Faibussowitsch   PetscCall(DMSwarmComputeMassMatrix_Private(dmCoarse, dmFine, *mass, PETSC_TRUE, ctx));
6879566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(*mass, NULL, "-mass_mat_view"));
6883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
68983c47955SMatthew G. Knepley }
69083c47955SMatthew G. Knepley 
691d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSwarmComputeMassMatrixSquare_Private(DM dmc, DM dmf, Mat mass, PetscBool useDeltaFunction, void *ctx)
692d71ae5a4SJacob Faibussowitsch {
6934cc7f7b2SMatthew G. Knepley   const char   *name = "Mass Matrix Square";
6944cc7f7b2SMatthew G. Knepley   MPI_Comm      comm;
69519307e5cSMatthew G. Knepley   DMSwarmCellDM celldm;
6964cc7f7b2SMatthew G. Knepley   PetscDS       prob;
6974cc7f7b2SMatthew G. Knepley   PetscSection  fsection, globalFSection;
6984cc7f7b2SMatthew G. Knepley   PetscHSetIJ   ht;
6994cc7f7b2SMatthew G. Knepley   PetscLayout   rLayout, colLayout;
7004cc7f7b2SMatthew G. Knepley   PetscInt     *dnz, *onz, *adj, depth, maxConeSize, maxSupportSize, maxAdjSize;
7014cc7f7b2SMatthew G. Knepley   PetscInt      locRows, locCols, rStart, colStart, colEnd, *rowIDXs;
7024cc7f7b2SMatthew G. Knepley   PetscReal    *xi, *v0, *J, *invJ, detJ = 1.0, v0ref[3] = {-1.0, -1.0, -1.0};
7034cc7f7b2SMatthew G. Knepley   PetscScalar  *elemMat, *elemMatSq;
70419307e5cSMatthew G. Knepley   PetscInt      cdim, Nf, Nfc, cStart, cEnd, totDim, maxC = 0;
70519307e5cSMatthew G. Knepley   const char  **coordFields;
70619307e5cSMatthew G. Knepley   PetscReal   **coordVals;
70719307e5cSMatthew G. Knepley   PetscInt     *bs;
7084cc7f7b2SMatthew G. Knepley 
7094cc7f7b2SMatthew G. Knepley   PetscFunctionBegin;
7109566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)mass, &comm));
7119566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dmf, &cdim));
7129566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dmf, &prob));
7139566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(prob, &Nf));
7149566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(prob, &totDim));
7159566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(cdim, &v0, cdim * cdim, &J, cdim * cdim, &invJ));
7169566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dmf, &fsection));
7179566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalSection(dmf, &globalFSection));
7189566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd));
7199566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(mass, &locRows, &locCols));
7204cc7f7b2SMatthew G. Knepley 
72119307e5cSMatthew G. Knepley   PetscCall(DMSwarmGetCellDMActive(dmc, &celldm));
72219307e5cSMatthew G. Knepley   PetscCall(DMSwarmCellDMGetCoordinateFields(celldm, &Nfc, &coordFields));
72319307e5cSMatthew G. Knepley   PetscCall(PetscMalloc2(Nfc, &coordVals, Nfc, &bs));
724d52c2f21SMatthew G. Knepley 
7259566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreate(comm, &colLayout));
7269566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetLocalSize(colLayout, locCols));
7279566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(colLayout, 1));
7289566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(colLayout));
7299566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetRange(colLayout, &colStart, &colEnd));
7309566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&colLayout));
7314cc7f7b2SMatthew G. Knepley 
7329566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreate(comm, &rLayout));
7339566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetLocalSize(rLayout, locRows));
7349566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(rLayout, 1));
7359566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(rLayout));
7369566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetRange(rLayout, &rStart, NULL));
7379566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&rLayout));
7384cc7f7b2SMatthew G. Knepley 
7399566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(dmf, &depth));
7409566063dSJacob Faibussowitsch   PetscCall(DMPlexGetMaxSizes(dmf, &maxConeSize, &maxSupportSize));
7414cc7f7b2SMatthew G. Knepley   maxAdjSize = PetscPowInt(maxConeSize * maxSupportSize, depth);
7429566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(maxAdjSize, &adj));
7434cc7f7b2SMatthew G. Knepley 
7449566063dSJacob Faibussowitsch   PetscCall(PetscCalloc2(locRows, &dnz, locRows, &onz));
7459566063dSJacob Faibussowitsch   PetscCall(PetscHSetIJCreate(&ht));
7464cc7f7b2SMatthew G. Knepley   /* Count nonzeros
7474cc7f7b2SMatthew G. Knepley        This is just FVM++, but we cannot use the Plex P0 allocation since unknowns in a cell will not be contiguous
7484cc7f7b2SMatthew G. Knepley   */
7499566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortGetAccess(dmc));
75019307e5cSMatthew G. Knepley   for (PetscInt cell = cStart; cell < cEnd; ++cell) {
7514cc7f7b2SMatthew G. Knepley     PetscInt *cindices;
7524cc7f7b2SMatthew G. Knepley     PetscInt  numCIndices;
7534cc7f7b2SMatthew G. Knepley #if 0
7544cc7f7b2SMatthew G. Knepley     PetscInt  adjSize = maxAdjSize, a, j;
7554cc7f7b2SMatthew G. Knepley #endif
7564cc7f7b2SMatthew G. Knepley 
7579566063dSJacob Faibussowitsch     PetscCall(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices));
7584cc7f7b2SMatthew G. Knepley     maxC = PetscMax(maxC, numCIndices);
7594cc7f7b2SMatthew G. Knepley     /* Diagonal block */
76019307e5cSMatthew G. Knepley     for (PetscInt i = 0; i < numCIndices; ++i) dnz[cindices[i]] += numCIndices;
7614cc7f7b2SMatthew G. Knepley #if 0
7624cc7f7b2SMatthew G. Knepley     /* Off-diagonal blocks */
7639566063dSJacob Faibussowitsch     PetscCall(DMPlexGetAdjacency(dmf, cell, &adjSize, &adj));
7644cc7f7b2SMatthew G. Knepley     for (a = 0; a < adjSize; ++a) {
7654cc7f7b2SMatthew G. Knepley       if (adj[a] >= cStart && adj[a] < cEnd && adj[a] != cell) {
7664cc7f7b2SMatthew G. Knepley         const PetscInt ncell = adj[a];
7674cc7f7b2SMatthew G. Knepley         PetscInt      *ncindices;
7684cc7f7b2SMatthew G. Knepley         PetscInt       numNCIndices;
7694cc7f7b2SMatthew G. Knepley 
7709566063dSJacob Faibussowitsch         PetscCall(DMSwarmSortGetPointsPerCell(dmc, ncell, &numNCIndices, &ncindices));
7714cc7f7b2SMatthew G. Knepley         {
7724cc7f7b2SMatthew G. Knepley           PetscHashIJKey key;
7734cc7f7b2SMatthew G. Knepley           PetscBool      missing;
7744cc7f7b2SMatthew G. Knepley 
7754cc7f7b2SMatthew G. Knepley           for (i = 0; i < numCIndices; ++i) {
7764cc7f7b2SMatthew G. Knepley             key.i = cindices[i] + rStart; /* global rows (from Swarm) */
7774cc7f7b2SMatthew G. Knepley             if (key.i < 0) continue;
7784cc7f7b2SMatthew G. Knepley             for (j = 0; j < numNCIndices; ++j) {
7794cc7f7b2SMatthew G. Knepley               key.j = ncindices[j] + rStart; /* global column (from Swarm) */
7804cc7f7b2SMatthew G. Knepley               if (key.j < 0) continue;
7819566063dSJacob Faibussowitsch               PetscCall(PetscHSetIJQueryAdd(ht, key, &missing));
7824cc7f7b2SMatthew G. Knepley               if (missing) {
7834cc7f7b2SMatthew G. Knepley                 if ((key.j >= colStart) && (key.j < colEnd)) ++dnz[key.i - rStart];
7844cc7f7b2SMatthew G. Knepley                 else                                         ++onz[key.i - rStart];
7854cc7f7b2SMatthew G. Knepley               }
7864cc7f7b2SMatthew G. Knepley             }
7874cc7f7b2SMatthew G. Knepley           }
7884cc7f7b2SMatthew G. Knepley         }
789fe11354eSMatthew G. Knepley         PetscCall(DMSwarmSortRestorePointsPerCell(dmc, ncell, &numNCIndices, &ncindices));
7904cc7f7b2SMatthew G. Knepley       }
7914cc7f7b2SMatthew G. Knepley     }
7924cc7f7b2SMatthew G. Knepley #endif
793fe11354eSMatthew G. Knepley     PetscCall(DMSwarmSortRestorePointsPerCell(dmc, cell, &numCIndices, &cindices));
7944cc7f7b2SMatthew G. Knepley   }
7959566063dSJacob Faibussowitsch   PetscCall(PetscHSetIJDestroy(&ht));
7969566063dSJacob Faibussowitsch   PetscCall(MatXAIJSetPreallocation(mass, 1, dnz, onz, NULL, NULL));
7979566063dSJacob Faibussowitsch   PetscCall(MatSetOption(mass, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
7989566063dSJacob Faibussowitsch   PetscCall(PetscFree2(dnz, onz));
7999566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(maxC * totDim, &elemMat, maxC * maxC, &elemMatSq, maxC, &rowIDXs, maxC * cdim, &xi));
8004cc7f7b2SMatthew G. Knepley   /* Fill in values
8014cc7f7b2SMatthew G. Knepley        Each entry is a sum of terms \phi_i(x_p) \phi_i(x_q)
8024cc7f7b2SMatthew G. Knepley        Start just by producing block diagonal
8034cc7f7b2SMatthew G. Knepley        Could loop over adjacent cells
8044cc7f7b2SMatthew G. Knepley          Produce neighboring element matrix
8054cc7f7b2SMatthew G. Knepley          TODO Determine which columns and rows correspond to shared dual vector
8064cc7f7b2SMatthew G. Knepley          Do MatMatMult with rectangular matrices
8074cc7f7b2SMatthew G. Knepley          Insert block
8084cc7f7b2SMatthew G. Knepley   */
80919307e5cSMatthew G. Knepley   for (PetscInt field = 0; field < Nf; ++field) {
8104cc7f7b2SMatthew G. Knepley     PetscTabulation Tcoarse;
8114cc7f7b2SMatthew G. Knepley     PetscObject     obj;
81219307e5cSMatthew G. Knepley     PetscInt        Nc;
8134cc7f7b2SMatthew G. Knepley 
8149566063dSJacob Faibussowitsch     PetscCall(PetscDSGetDiscretization(prob, field, &obj));
8159566063dSJacob Faibussowitsch     PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc));
81663a3b9bcSJacob Faibussowitsch     PetscCheck(Nc == 1, PetscObjectComm((PetscObject)dmf), PETSC_ERR_SUP, "Can only interpolate a scalar field from particles, Nc = %" PetscInt_FMT, Nc);
81719307e5cSMatthew G. Knepley     for (PetscInt i = 0; i < Nfc; ++i) PetscCall(DMSwarmGetField(dmc, coordFields[i], &bs[i], NULL, (void **)&coordVals[i]));
81819307e5cSMatthew G. Knepley     for (PetscInt cell = cStart; cell < cEnd; ++cell) {
8194cc7f7b2SMatthew G. Knepley       PetscInt *findices, *cindices;
8204cc7f7b2SMatthew G. Knepley       PetscInt  numFIndices, numCIndices;
8214cc7f7b2SMatthew G. Knepley 
8224cc7f7b2SMatthew G. Knepley       /* TODO: Use DMField instead of assuming affine */
8239566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ));
8249566063dSJacob Faibussowitsch       PetscCall(DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
8259566063dSJacob Faibussowitsch       PetscCall(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices));
82619307e5cSMatthew G. Knepley       for (PetscInt p = 0; p < numCIndices; ++p) {
82719307e5cSMatthew G. Knepley         PetscReal xr[8];
82819307e5cSMatthew G. Knepley         PetscInt  off = 0;
82919307e5cSMatthew G. Knepley 
83019307e5cSMatthew G. Knepley         for (PetscInt i = 0; i < Nfc; ++i) {
83119307e5cSMatthew G. Knepley           for (PetscInt b = 0; b < bs[i]; ++b, ++off) xr[off] = coordVals[i][cindices[p] * bs[i] + b];
83219307e5cSMatthew G. Knepley         }
83319307e5cSMatthew G. Knepley         PetscCheck(off == cdim, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The total block size of coordinates is %" PetscInt_FMT " != %" PetscInt_FMT " the DM coordinate dimension", off, cdim);
83419307e5cSMatthew G. Knepley         CoordinatesRealToRef(cdim, cdim, v0ref, v0, invJ, xr, &xi[p * cdim]);
83519307e5cSMatthew G. Knepley       }
8369566063dSJacob Faibussowitsch       PetscCall(PetscFECreateTabulation((PetscFE)obj, 1, numCIndices, xi, 0, &Tcoarse));
8374cc7f7b2SMatthew G. Knepley       /* Get elemMat entries by multiplying by weight */
8389566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(elemMat, numCIndices * totDim));
83919307e5cSMatthew G. Knepley       for (PetscInt i = 0; i < numFIndices; ++i) {
84019307e5cSMatthew G. Knepley         for (PetscInt p = 0; p < numCIndices; ++p) {
84119307e5cSMatthew G. Knepley           for (PetscInt c = 0; c < Nc; ++c) {
8424cc7f7b2SMatthew G. Knepley             /* B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c */
8434cc7f7b2SMatthew G. Knepley             elemMat[p * numFIndices + i] += Tcoarse->T[0][(p * numFIndices + i) * Nc + c] * (useDeltaFunction ? 1.0 : detJ);
8444cc7f7b2SMatthew G. Knepley           }
8454cc7f7b2SMatthew G. Knepley         }
8464cc7f7b2SMatthew G. Knepley       }
8479566063dSJacob Faibussowitsch       PetscCall(PetscTabulationDestroy(&Tcoarse));
84819307e5cSMatthew G. Knepley       for (PetscInt p = 0; p < numCIndices; ++p) rowIDXs[p] = cindices[p] + rStart;
8499566063dSJacob Faibussowitsch       if (0) PetscCall(DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat));
8504cc7f7b2SMatthew G. Knepley       /* Block diagonal */
85178884ca7SMatthew G. Knepley       if (numCIndices) {
8524cc7f7b2SMatthew G. Knepley         PetscBLASInt blasn, blask;
8534cc7f7b2SMatthew G. Knepley         PetscScalar  one = 1.0, zero = 0.0;
8544cc7f7b2SMatthew G. Knepley 
8559566063dSJacob Faibussowitsch         PetscCall(PetscBLASIntCast(numCIndices, &blasn));
8569566063dSJacob Faibussowitsch         PetscCall(PetscBLASIntCast(numFIndices, &blask));
857792fecdfSBarry Smith         PetscCallBLAS("BLASgemm", BLASgemm_("T", "N", &blasn, &blasn, &blask, &one, elemMat, &blask, elemMat, &blask, &zero, elemMatSq, &blasn));
8584cc7f7b2SMatthew G. Knepley       }
8599566063dSJacob Faibussowitsch       PetscCall(MatSetValues(mass, numCIndices, rowIDXs, numCIndices, rowIDXs, elemMatSq, ADD_VALUES));
8604cf0e950SBarry Smith       /* TODO off-diagonal */
861fe11354eSMatthew G. Knepley       PetscCall(DMSwarmSortRestorePointsPerCell(dmc, cell, &numCIndices, &cindices));
8629566063dSJacob Faibussowitsch       PetscCall(DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
8634cc7f7b2SMatthew G. Knepley     }
86419307e5cSMatthew G. Knepley     for (PetscInt i = 0; i < Nfc; ++i) PetscCall(DMSwarmRestoreField(dmc, coordFields[i], &bs[i], NULL, (void **)&coordVals[i]));
8654cc7f7b2SMatthew G. Knepley   }
8669566063dSJacob Faibussowitsch   PetscCall(PetscFree4(elemMat, elemMatSq, rowIDXs, xi));
8679566063dSJacob Faibussowitsch   PetscCall(PetscFree(adj));
8689566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortRestoreAccess(dmc));
8699566063dSJacob Faibussowitsch   PetscCall(PetscFree3(v0, J, invJ));
87019307e5cSMatthew G. Knepley   PetscCall(PetscFree2(coordVals, bs));
8719566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mass, MAT_FINAL_ASSEMBLY));
8729566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mass, MAT_FINAL_ASSEMBLY));
8733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8744cc7f7b2SMatthew G. Knepley }
8754cc7f7b2SMatthew G. Knepley 
876cc4c1da9SBarry Smith /*@
8774cc7f7b2SMatthew G. Knepley   DMSwarmCreateMassMatrixSquare - Creates the block-diagonal of the square, M^T_p M_p, of the particle mass matrix M_p
8784cc7f7b2SMatthew G. Knepley 
87920f4b53cSBarry Smith   Collective
8804cc7f7b2SMatthew G. Knepley 
88160225df5SJacob Faibussowitsch   Input Parameters:
88220f4b53cSBarry Smith + dmCoarse - a `DMSWARM`
88320f4b53cSBarry Smith - dmFine   - a `DMPLEX`
8844cc7f7b2SMatthew G. Knepley 
88560225df5SJacob Faibussowitsch   Output Parameter:
8864cc7f7b2SMatthew G. Knepley . mass - the square of the particle mass matrix
8874cc7f7b2SMatthew G. Knepley 
8884cc7f7b2SMatthew G. Knepley   Level: advanced
8894cc7f7b2SMatthew G. Knepley 
89020f4b53cSBarry Smith   Note:
8914cc7f7b2SMatthew G. Knepley   We only compute the block diagonal since this provides a good preconditioner and is completely local. It would be possible in the
8924cc7f7b2SMatthew G. Knepley   future to compute the full normal equations.
8934cc7f7b2SMatthew G. Knepley 
89420f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMCreateMassMatrix()`
8954cc7f7b2SMatthew G. Knepley @*/
896d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmCreateMassMatrixSquare(DM dmCoarse, DM dmFine, Mat *mass)
897d71ae5a4SJacob Faibussowitsch {
8984cc7f7b2SMatthew G. Knepley   PetscInt n;
8994cc7f7b2SMatthew G. Knepley   void    *ctx;
9004cc7f7b2SMatthew G. Knepley 
9014cc7f7b2SMatthew G. Knepley   PetscFunctionBegin;
9029566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetLocalSize(dmCoarse, &n));
9039566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)dmCoarse), mass));
9049566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*mass, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
9059566063dSJacob Faibussowitsch   PetscCall(MatSetType(*mass, dmCoarse->mattype));
9069566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(dmFine, &ctx));
9074cc7f7b2SMatthew G. Knepley 
9089566063dSJacob Faibussowitsch   PetscCall(DMSwarmComputeMassMatrixSquare_Private(dmCoarse, dmFine, *mass, PETSC_TRUE, ctx));
9099566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(*mass, NULL, "-mass_sq_mat_view"));
9103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9114cc7f7b2SMatthew G. Knepley }
9124cc7f7b2SMatthew G. Knepley 
9131898fd5cSMatthew G. Knepley /* This creates a "gradient matrix" between a finite element and particle space, which is meant to enforce a weak divergence condition on the particle space. We are looking for a finite element field that has the same divergence as our particle field, so that
9141898fd5cSMatthew G. Knepley 
9151898fd5cSMatthew G. Knepley      \int_X \psi_i \nabla \cdot \hat f = \int_X \psi_i \nabla \cdot f
9161898fd5cSMatthew G. Knepley 
9171898fd5cSMatthew G. Knepley    and then integrate by parts
9181898fd5cSMatthew G. Knepley 
9191898fd5cSMatthew G. Knepley      \int_X \nabla \psi_i \cdot \hat f = \int_X \nabla \psi_i \cdot f
9201898fd5cSMatthew G. Knepley 
9211898fd5cSMatthew G. Knepley    where \psi is from a scalar FE space. If a finite element interpolant is given by
9221898fd5cSMatthew G. Knepley 
9231898fd5cSMatthew G. Knepley      \hat f^c = \sum_i f_i \phi^c_i
9241898fd5cSMatthew G. Knepley 
9251898fd5cSMatthew G. Knepley    and a particle function is given by
9261898fd5cSMatthew G. Knepley 
9271898fd5cSMatthew G. Knepley      f^c = \sum_p f^c_p \delta(x - x_p)
9281898fd5cSMatthew G. Knepley 
9291898fd5cSMatthew G. Knepley    then we want to require that
9301898fd5cSMatthew G. Knepley 
9311898fd5cSMatthew G. Knepley      D_f \hat f = D_p f
9321898fd5cSMatthew G. Knepley 
9331898fd5cSMatthew G. Knepley    where the gradient matrices are given by
9341898fd5cSMatthew G. Knepley 
9351898fd5cSMatthew G. Knepley      (D_f)_{i(jc)} = \int \partial_c \psi_i \phi_j
9361898fd5cSMatthew G. Knepley      (D_p)_{i(jc)} = \int \partial_c \psi_i \delta(x - x_j)
9371898fd5cSMatthew G. Knepley 
9381898fd5cSMatthew G. Knepley    Thus we need two finite element spaces, a scalar and a vector. The vector space holds the representer for the
9391898fd5cSMatthew G. Knepley    vector particle field. The scalar space holds the output of D_p or D_f, which is the weak divergence of the field.
9401898fd5cSMatthew G. Knepley 
9411898fd5cSMatthew G. Knepley    The way Dave May does particles, they amount to quadratue weights rather than delta functions, so he has |J| is in
9421898fd5cSMatthew G. Knepley    his integral. We allow this with the boolean flag.
9431898fd5cSMatthew G. Knepley */
9441898fd5cSMatthew G. Knepley static PetscErrorCode DMSwarmComputeGradientMatrix_Private(DM sw, DM dm, Mat derv, PetscBool useDeltaFunction, void *ctx)
9451898fd5cSMatthew G. Knepley {
9461898fd5cSMatthew G. Knepley   const char   *name = "Derivative Matrix";
9471898fd5cSMatthew G. Knepley   MPI_Comm      comm;
9481898fd5cSMatthew G. Knepley   DMSwarmCellDM celldm;
9491898fd5cSMatthew G. Knepley   PetscDS       ds;
9501898fd5cSMatthew G. Knepley   PetscSection  fsection, globalFSection;
9511898fd5cSMatthew G. Knepley   PetscLayout   rLayout;
9521898fd5cSMatthew G. Knepley   PetscInt      locRows, rStart, *rowIDXs;
9531898fd5cSMatthew G. Knepley   PetscReal    *xi, *v0, *J, *invJ, detJ = 1.0, v0ref[3] = {-1.0, -1.0, -1.0};
9541898fd5cSMatthew G. Knepley   PetscScalar  *elemMat;
9551898fd5cSMatthew G. Knepley   PetscInt      cdim, Nf, Nfc, cStart, cEnd, totDim, maxNpc = 0, totNc = 0;
9561898fd5cSMatthew G. Knepley   const char  **coordFields;
9571898fd5cSMatthew G. Knepley   PetscReal   **coordVals;
9581898fd5cSMatthew G. Knepley   PetscInt     *bs;
9591898fd5cSMatthew G. Knepley 
9601898fd5cSMatthew G. Knepley   PetscFunctionBegin;
9611898fd5cSMatthew G. Knepley   PetscCall(PetscObjectGetComm((PetscObject)derv, &comm));
9621898fd5cSMatthew G. Knepley   PetscCall(DMGetCoordinateDim(dm, &cdim));
9631898fd5cSMatthew G. Knepley   PetscCall(DMGetDS(dm, &ds));
9641898fd5cSMatthew G. Knepley   PetscCall(PetscDSGetNumFields(ds, &Nf));
9651898fd5cSMatthew G. Knepley   PetscCheck(Nf == 1, comm, PETSC_ERR_SUP, "Currently, we only support a single field");
9661898fd5cSMatthew G. Knepley   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
9671898fd5cSMatthew G. Knepley   PetscCall(PetscMalloc3(cdim, &v0, cdim * cdim, &J, cdim * cdim, &invJ));
9681898fd5cSMatthew G. Knepley   PetscCall(DMGetLocalSection(dm, &fsection));
9691898fd5cSMatthew G. Knepley   PetscCall(DMGetGlobalSection(dm, &globalFSection));
9701898fd5cSMatthew G. Knepley   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
9711898fd5cSMatthew G. Knepley   PetscCall(MatGetLocalSize(derv, &locRows, NULL));
9721898fd5cSMatthew G. Knepley 
9731898fd5cSMatthew G. Knepley   PetscCall(DMSwarmGetCellDMActive(sw, &celldm));
9741898fd5cSMatthew G. Knepley   PetscCall(DMSwarmCellDMGetCoordinateFields(celldm, &Nfc, &coordFields));
9751898fd5cSMatthew G. Knepley   PetscCheck(Nfc == 1, comm, PETSC_ERR_SUP, "Currently, we only support a single field");
9761898fd5cSMatthew G. Knepley   PetscCall(PetscMalloc2(Nfc, &coordVals, Nfc, &bs));
9771898fd5cSMatthew G. Knepley 
9781898fd5cSMatthew G. Knepley   PetscCall(PetscLayoutCreate(comm, &rLayout));
9791898fd5cSMatthew G. Knepley   PetscCall(PetscLayoutSetLocalSize(rLayout, locRows));
9801898fd5cSMatthew G. Knepley   PetscCall(PetscLayoutSetBlockSize(rLayout, cdim));
9811898fd5cSMatthew G. Knepley   PetscCall(PetscLayoutSetUp(rLayout));
9821898fd5cSMatthew G. Knepley   PetscCall(PetscLayoutGetRange(rLayout, &rStart, NULL));
9831898fd5cSMatthew G. Knepley   PetscCall(PetscLayoutDestroy(&rLayout));
9841898fd5cSMatthew G. Knepley 
9851898fd5cSMatthew G. Knepley   for (PetscInt field = 0; field < Nf; ++field) {
9861898fd5cSMatthew G. Knepley     PetscObject obj;
9871898fd5cSMatthew G. Knepley     PetscInt    Nc;
9881898fd5cSMatthew G. Knepley 
9891898fd5cSMatthew G. Knepley     PetscCall(PetscDSGetDiscretization(ds, field, &obj));
9901898fd5cSMatthew G. Knepley     PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc));
9911898fd5cSMatthew G. Knepley     totNc += Nc;
9921898fd5cSMatthew G. Knepley   }
9931898fd5cSMatthew G. Knepley   PetscCheck(totNc == 1, comm, PETSC_ERR_ARG_WRONG, "The number of field components %" PetscInt_FMT " != 1", totNc);
9941898fd5cSMatthew G. Knepley   /* count non-zeros */
9951898fd5cSMatthew G. Knepley   PetscCall(DMSwarmSortGetAccess(sw));
9961898fd5cSMatthew G. Knepley   for (PetscInt field = 0; field < Nf; ++field) {
9971898fd5cSMatthew G. Knepley     PetscObject obj;
9981898fd5cSMatthew G. Knepley     PetscInt    Nc;
9991898fd5cSMatthew G. Knepley 
10001898fd5cSMatthew G. Knepley     PetscCall(PetscDSGetDiscretization(ds, field, &obj));
10011898fd5cSMatthew G. Knepley     PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc));
10021898fd5cSMatthew G. Knepley     for (PetscInt cell = cStart; cell < cEnd; ++cell) {
10031898fd5cSMatthew G. Knepley       PetscInt *pind;
10041898fd5cSMatthew G. Knepley       PetscInt  Npc;
10051898fd5cSMatthew G. Knepley 
10061898fd5cSMatthew G. Knepley       PetscCall(DMSwarmSortGetPointsPerCell(sw, cell, &Npc, &pind));
10071898fd5cSMatthew G. Knepley       maxNpc = PetscMax(maxNpc, Npc);
10081898fd5cSMatthew G. Knepley       PetscCall(DMSwarmSortRestorePointsPerCell(sw, cell, &Npc, &pind));
10091898fd5cSMatthew G. Knepley     }
10101898fd5cSMatthew G. Knepley   }
10111898fd5cSMatthew G. Knepley   PetscCall(PetscMalloc3(maxNpc * cdim * totDim, &elemMat, maxNpc * cdim, &rowIDXs, maxNpc * cdim, &xi));
10121898fd5cSMatthew G. Knepley   for (PetscInt field = 0; field < Nf; ++field) {
10131898fd5cSMatthew G. Knepley     PetscTabulation Tcoarse;
10141898fd5cSMatthew G. Knepley     PetscFE         fe;
10151898fd5cSMatthew G. Knepley 
10161898fd5cSMatthew G. Knepley     PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe));
10171898fd5cSMatthew G. Knepley     for (PetscInt i = 0; i < Nfc; ++i) PetscCall(DMSwarmGetField(sw, coordFields[i], &bs[i], NULL, (void **)&coordVals[i]));
10181898fd5cSMatthew G. Knepley     for (PetscInt cell = cStart; cell < cEnd; ++cell) {
10191898fd5cSMatthew G. Knepley       PetscInt *findices, *pind;
10201898fd5cSMatthew G. Knepley       PetscInt  numFIndices, Npc;
10211898fd5cSMatthew G. Knepley 
10221898fd5cSMatthew G. Knepley       /* TODO: Use DMField instead of assuming affine */
10231898fd5cSMatthew G. Knepley       PetscCall(DMPlexComputeCellGeometryFEM(dm, cell, NULL, v0, J, invJ, &detJ));
10241898fd5cSMatthew G. Knepley       PetscCall(DMPlexGetClosureIndices(dm, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
10251898fd5cSMatthew G. Knepley       PetscCall(DMSwarmSortGetPointsPerCell(sw, cell, &Npc, &pind));
10261898fd5cSMatthew G. Knepley       for (PetscInt j = 0; j < Npc; ++j) {
10271898fd5cSMatthew G. Knepley         PetscReal xr[8];
10281898fd5cSMatthew G. Knepley         PetscInt  off = 0;
10291898fd5cSMatthew G. Knepley 
10301898fd5cSMatthew G. Knepley         for (PetscInt i = 0; i < Nfc; ++i) {
10311898fd5cSMatthew G. Knepley           for (PetscInt b = 0; b < bs[i]; ++b, ++off) xr[off] = coordVals[i][pind[j] * bs[i] + b];
10321898fd5cSMatthew G. Knepley         }
10331898fd5cSMatthew G. Knepley         PetscCheck(off == cdim, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The total block size of coordinates is %" PetscInt_FMT " != %" PetscInt_FMT " the DM coordinate dimension", off, cdim);
10341898fd5cSMatthew G. Knepley         CoordinatesRealToRef(cdim, cdim, v0ref, v0, invJ, xr, &xi[j * cdim]);
10351898fd5cSMatthew G. Knepley       }
10361898fd5cSMatthew G. Knepley       PetscCall(PetscFECreateTabulation(fe, 1, Npc, xi, 1, &Tcoarse));
10371898fd5cSMatthew G. Knepley       /* Get elemMat entries by multiplying by weight */
10381898fd5cSMatthew G. Knepley       PetscCall(PetscArrayzero(elemMat, Npc * cdim * totDim));
10391898fd5cSMatthew G. Knepley       for (PetscInt i = 0; i < numFIndices; ++i) {
10401898fd5cSMatthew G. Knepley         for (PetscInt j = 0; j < Npc; ++j) {
10411898fd5cSMatthew G. Knepley           /* D[((p*pdim + i)*Nc + c)*cdim + d] is the value at point p for basis function i, component c, derviative d */
10421898fd5cSMatthew G. Knepley           for (PetscInt d = 0; d < cdim; ++d) {
10431898fd5cSMatthew G. Knepley             xi[d] = 0.;
10441898fd5cSMatthew G. Knepley             for (PetscInt e = 0; e < cdim; ++e) xi[d] += invJ[e * cdim + d] * Tcoarse->T[1][(j * numFIndices + i) * cdim + e];
10451898fd5cSMatthew G. Knepley             elemMat[(j * cdim + d) * numFIndices + i] += xi[d] * (useDeltaFunction ? 1.0 : detJ);
10461898fd5cSMatthew G. Knepley           }
10471898fd5cSMatthew G. Knepley         }
10481898fd5cSMatthew G. Knepley       }
10491898fd5cSMatthew G. Knepley       for (PetscInt j = 0; j < Npc; ++j)
10501898fd5cSMatthew G. Knepley         for (PetscInt d = 0; d < cdim; ++d) rowIDXs[j * cdim + d] = pind[j] * cdim + d + rStart;
10511898fd5cSMatthew G. Knepley       if (0) PetscCall(DMPrintCellMatrix(cell, name, Npc * cdim, numFIndices, elemMat));
10521898fd5cSMatthew G. Knepley       PetscCall(MatSetValues(derv, Npc * cdim, rowIDXs, numFIndices, findices, elemMat, ADD_VALUES));
10531898fd5cSMatthew G. Knepley       PetscCall(DMSwarmSortRestorePointsPerCell(sw, cell, &Npc, &pind));
10541898fd5cSMatthew G. Knepley       PetscCall(DMPlexRestoreClosureIndices(dm, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
10551898fd5cSMatthew G. Knepley       PetscCall(PetscTabulationDestroy(&Tcoarse));
10561898fd5cSMatthew G. Knepley     }
10571898fd5cSMatthew G. Knepley     for (PetscInt i = 0; i < Nfc; ++i) PetscCall(DMSwarmRestoreField(sw, coordFields[i], &bs[i], NULL, (void **)&coordVals[i]));
10581898fd5cSMatthew G. Knepley   }
10591898fd5cSMatthew G. Knepley   PetscCall(PetscFree3(elemMat, rowIDXs, xi));
10601898fd5cSMatthew G. Knepley   PetscCall(DMSwarmSortRestoreAccess(sw));
10611898fd5cSMatthew G. Knepley   PetscCall(PetscFree3(v0, J, invJ));
10621898fd5cSMatthew G. Knepley   PetscCall(PetscFree2(coordVals, bs));
10631898fd5cSMatthew G. Knepley   PetscCall(MatAssemblyBegin(derv, MAT_FINAL_ASSEMBLY));
10641898fd5cSMatthew G. Knepley   PetscCall(MatAssemblyEnd(derv, MAT_FINAL_ASSEMBLY));
10651898fd5cSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
10661898fd5cSMatthew G. Knepley }
10671898fd5cSMatthew G. Knepley 
10681898fd5cSMatthew G. Knepley /* FEM cols:      this is a scalar space
10691898fd5cSMatthew G. Knepley    Particle rows: this is a vector space that contracts with the derivative
10701898fd5cSMatthew G. Knepley */
10711898fd5cSMatthew G. Knepley static PetscErrorCode DMCreateGradientMatrix_Swarm(DM sw, DM dm, Mat *derv)
10721898fd5cSMatthew G. Knepley {
10731898fd5cSMatthew G. Knepley   DMSwarmCellDM celldm;
10741898fd5cSMatthew G. Knepley   PetscSection  gs;
10751898fd5cSMatthew G. Knepley   PetscInt      cdim, m, n, Np, bs;
10761898fd5cSMatthew G. Knepley   void         *ctx;
10771898fd5cSMatthew G. Knepley   MPI_Comm      comm;
10781898fd5cSMatthew G. Knepley 
10791898fd5cSMatthew G. Knepley   PetscFunctionBegin;
10801898fd5cSMatthew G. Knepley   PetscCall(PetscObjectGetComm((PetscObject)sw, &comm));
10811898fd5cSMatthew G. Knepley   PetscCall(DMGetCoordinateDim(dm, &cdim));
10821898fd5cSMatthew G. Knepley   PetscCall(DMSwarmGetCellDMActive(sw, &celldm));
10831898fd5cSMatthew G. Knepley   PetscCheck(celldm->Nf, comm, PETSC_ERR_USER, "Active cell DM does not define any fields");
10841898fd5cSMatthew G. Knepley   PetscCall(DMGetGlobalSection(dm, &gs));
10851898fd5cSMatthew G. Knepley   PetscCall(PetscSectionGetConstrainedStorageSize(gs, &n));
10861898fd5cSMatthew G. Knepley   PetscCall(DMSwarmGetLocalSize(sw, &Np));
10871898fd5cSMatthew G. Knepley   PetscCall(DMSwarmCellDMGetBlockSize(celldm, sw, &bs));
10881898fd5cSMatthew G. Knepley   PetscCheck(cdim == bs, comm, PETSC_ERR_ARG_WRONG, "Coordinate dimension %" PetscInt_FMT " != %" PetscInt_FMT " swarm field block size", cdim, bs);
10891898fd5cSMatthew G. Knepley   m = Np * bs;
10901898fd5cSMatthew G. Knepley   PetscCall(MatCreate(PetscObjectComm((PetscObject)sw), derv));
10911898fd5cSMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)*derv, "Swarm Derivative Matrix"));
10921898fd5cSMatthew G. Knepley   PetscCall(MatSetSizes(*derv, m, n, PETSC_DETERMINE, PETSC_DETERMINE));
10931898fd5cSMatthew G. Knepley   PetscCall(MatSetType(*derv, sw->mattype));
10941898fd5cSMatthew G. Knepley   PetscCall(DMGetApplicationContext(dm, &ctx));
10951898fd5cSMatthew G. Knepley 
10961898fd5cSMatthew G. Knepley   PetscCall(DMSwarmComputeGradientMatrix_Private(sw, dm, *derv, PETSC_TRUE, ctx));
10971898fd5cSMatthew G. Knepley   PetscCall(MatViewFromOptions(*derv, NULL, "-gradient_mat_view"));
10981898fd5cSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
10991898fd5cSMatthew G. Knepley }
11001898fd5cSMatthew G. Knepley 
1101cc4c1da9SBarry Smith /*@
110220f4b53cSBarry Smith   DMSwarmCreateGlobalVectorFromField - Creates a `Vec` object sharing the array associated with a given field
1103d3a51819SDave May 
110420f4b53cSBarry Smith   Collective
1105d3a51819SDave May 
110660225df5SJacob Faibussowitsch   Input Parameters:
110720f4b53cSBarry Smith + dm        - a `DMSWARM`
110862741f57SDave May - fieldname - the textual name given to a registered field
1109d3a51819SDave May 
111060225df5SJacob Faibussowitsch   Output Parameter:
1111d3a51819SDave May . vec - the vector
1112d3a51819SDave May 
1113d3a51819SDave May   Level: beginner
1114d3a51819SDave May 
1115cc4c1da9SBarry Smith   Note:
111620f4b53cSBarry Smith   The vector must be returned using a matching call to `DMSwarmDestroyGlobalVectorFromField()`.
11178b8a3813SDave May 
111820f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmDestroyGlobalVectorFromField()`
1119d3a51819SDave May @*/
1120d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmCreateGlobalVectorFromField(DM dm, const char fieldname[], Vec *vec)
1121d71ae5a4SJacob Faibussowitsch {
1122fb1bcc12SMatthew G. Knepley   MPI_Comm comm = PetscObjectComm((PetscObject)dm);
1123b5bcf523SDave May 
1124fb1bcc12SMatthew G. Knepley   PetscFunctionBegin;
1125ea7032a0SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
11269566063dSJacob Faibussowitsch   PetscCall(DMSwarmCreateVectorFromField_Private(dm, fieldname, comm, vec));
11273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1128b5bcf523SDave May }
1129b5bcf523SDave May 
1130cc4c1da9SBarry Smith /*@
113120f4b53cSBarry Smith   DMSwarmDestroyGlobalVectorFromField - Destroys the `Vec` object which share the array associated with a given field
1132d3a51819SDave May 
113320f4b53cSBarry Smith   Collective
1134d3a51819SDave May 
113560225df5SJacob Faibussowitsch   Input Parameters:
113620f4b53cSBarry Smith + dm        - a `DMSWARM`
113762741f57SDave May - fieldname - the textual name given to a registered field
1138d3a51819SDave May 
113960225df5SJacob Faibussowitsch   Output Parameter:
1140d3a51819SDave May . vec - the vector
1141d3a51819SDave May 
1142d3a51819SDave May   Level: beginner
1143d3a51819SDave May 
114420f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmCreateGlobalVectorFromField()`
1145d3a51819SDave May @*/
1146d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDestroyGlobalVectorFromField(DM dm, const char fieldname[], Vec *vec)
1147d71ae5a4SJacob Faibussowitsch {
1148fb1bcc12SMatthew G. Knepley   PetscFunctionBegin;
1149ea7032a0SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
11509566063dSJacob Faibussowitsch   PetscCall(DMSwarmDestroyVectorFromField_Private(dm, fieldname, vec));
11513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1152b5bcf523SDave May }
1153b5bcf523SDave May 
1154cc4c1da9SBarry Smith /*@
115520f4b53cSBarry Smith   DMSwarmCreateLocalVectorFromField - Creates a `Vec` object sharing the array associated with a given field
1156fb1bcc12SMatthew G. Knepley 
115720f4b53cSBarry Smith   Collective
1158fb1bcc12SMatthew G. Knepley 
115960225df5SJacob Faibussowitsch   Input Parameters:
116020f4b53cSBarry Smith + dm        - a `DMSWARM`
116162741f57SDave May - fieldname - the textual name given to a registered field
1162fb1bcc12SMatthew G. Knepley 
116360225df5SJacob Faibussowitsch   Output Parameter:
1164fb1bcc12SMatthew G. Knepley . vec - the vector
1165fb1bcc12SMatthew G. Knepley 
1166fb1bcc12SMatthew G. Knepley   Level: beginner
1167fb1bcc12SMatthew G. Knepley 
116820f4b53cSBarry Smith   Note:
11698b8a3813SDave May   The vector must be returned using a matching call to DMSwarmDestroyLocalVectorFromField().
11708b8a3813SDave May 
117120f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmDestroyLocalVectorFromField()`
1172fb1bcc12SMatthew G. Knepley @*/
1173d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmCreateLocalVectorFromField(DM dm, const char fieldname[], Vec *vec)
1174d71ae5a4SJacob Faibussowitsch {
1175fb1bcc12SMatthew G. Knepley   MPI_Comm comm = PETSC_COMM_SELF;
1176bbe8250bSMatthew G. Knepley 
1177fb1bcc12SMatthew G. Knepley   PetscFunctionBegin;
11789566063dSJacob Faibussowitsch   PetscCall(DMSwarmCreateVectorFromField_Private(dm, fieldname, comm, vec));
11793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1180bbe8250bSMatthew G. Knepley }
1181fb1bcc12SMatthew G. Knepley 
1182cc4c1da9SBarry Smith /*@
118320f4b53cSBarry Smith   DMSwarmDestroyLocalVectorFromField - Destroys the `Vec` object which share the array associated with a given field
1184fb1bcc12SMatthew G. Knepley 
118520f4b53cSBarry Smith   Collective
1186fb1bcc12SMatthew G. Knepley 
118760225df5SJacob Faibussowitsch   Input Parameters:
118820f4b53cSBarry Smith + dm        - a `DMSWARM`
118962741f57SDave May - fieldname - the textual name given to a registered field
1190fb1bcc12SMatthew G. Knepley 
119160225df5SJacob Faibussowitsch   Output Parameter:
1192fb1bcc12SMatthew G. Knepley . vec - the vector
1193fb1bcc12SMatthew G. Knepley 
1194fb1bcc12SMatthew G. Knepley   Level: beginner
1195fb1bcc12SMatthew G. Knepley 
119620f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmCreateLocalVectorFromField()`
1197fb1bcc12SMatthew G. Knepley @*/
1198d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDestroyLocalVectorFromField(DM dm, const char fieldname[], Vec *vec)
1199d71ae5a4SJacob Faibussowitsch {
1200fb1bcc12SMatthew G. Knepley   PetscFunctionBegin;
1201ea7032a0SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
12029566063dSJacob Faibussowitsch   PetscCall(DMSwarmDestroyVectorFromField_Private(dm, fieldname, vec));
12033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1204bbe8250bSMatthew G. Knepley }
1205bbe8250bSMatthew G. Knepley 
1206cc4c1da9SBarry Smith /*@
120719307e5cSMatthew G. Knepley   DMSwarmCreateGlobalVectorFromFields - Creates a `Vec` object sharing the array associated with a given field set
120819307e5cSMatthew G. Knepley 
120919307e5cSMatthew G. Knepley   Collective
121019307e5cSMatthew G. Knepley 
121119307e5cSMatthew G. Knepley   Input Parameters:
121219307e5cSMatthew G. Knepley + dm         - a `DMSWARM`
121319307e5cSMatthew G. Knepley . Nf         - the number of fields
121419307e5cSMatthew G. Knepley - fieldnames - the textual names given to the registered fields
121519307e5cSMatthew G. Knepley 
121619307e5cSMatthew G. Knepley   Output Parameter:
121719307e5cSMatthew G. Knepley . vec - the vector
121819307e5cSMatthew G. Knepley 
121919307e5cSMatthew G. Knepley   Level: beginner
122019307e5cSMatthew G. Knepley 
122119307e5cSMatthew G. Knepley   Notes:
122219307e5cSMatthew G. Knepley   The vector must be returned using a matching call to `DMSwarmDestroyGlobalVectorFromFields()`.
122319307e5cSMatthew G. Knepley 
122419307e5cSMatthew G. Knepley   This vector is copyin-copyout, rather than a direct pointer like `DMSwarmCreateGlobalVectorFromField()`
122519307e5cSMatthew G. Knepley 
122619307e5cSMatthew G. Knepley .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmDestroyGlobalVectorFromFields()`
122719307e5cSMatthew G. Knepley @*/
122819307e5cSMatthew G. Knepley PetscErrorCode DMSwarmCreateGlobalVectorFromFields(DM dm, PetscInt Nf, const char *fieldnames[], Vec *vec)
122919307e5cSMatthew G. Knepley {
123019307e5cSMatthew G. Knepley   MPI_Comm comm = PetscObjectComm((PetscObject)dm);
123119307e5cSMatthew G. Knepley 
123219307e5cSMatthew G. Knepley   PetscFunctionBegin;
123319307e5cSMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
123419307e5cSMatthew G. Knepley   PetscCall(DMSwarmCreateVectorFromFields_Private(dm, Nf, fieldnames, comm, vec));
123519307e5cSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
123619307e5cSMatthew G. Knepley }
123719307e5cSMatthew G. Knepley 
123819307e5cSMatthew G. Knepley /*@
123919307e5cSMatthew G. Knepley   DMSwarmDestroyGlobalVectorFromFields - Destroys the `Vec` object which share the array associated with a given field set
124019307e5cSMatthew G. Knepley 
124119307e5cSMatthew G. Knepley   Collective
124219307e5cSMatthew G. Knepley 
124319307e5cSMatthew G. Knepley   Input Parameters:
124419307e5cSMatthew G. Knepley + dm         - a `DMSWARM`
124519307e5cSMatthew G. Knepley . Nf         - the number of fields
124619307e5cSMatthew G. Knepley - fieldnames - the textual names given to the registered fields
124719307e5cSMatthew G. Knepley 
124819307e5cSMatthew G. Knepley   Output Parameter:
124919307e5cSMatthew G. Knepley . vec - the vector
125019307e5cSMatthew G. Knepley 
125119307e5cSMatthew G. Knepley   Level: beginner
125219307e5cSMatthew G. Knepley 
125319307e5cSMatthew G. Knepley .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmCreateGlobalVectorFromField()`
125419307e5cSMatthew G. Knepley @*/
125519307e5cSMatthew G. Knepley PetscErrorCode DMSwarmDestroyGlobalVectorFromFields(DM dm, PetscInt Nf, const char *fieldnames[], Vec *vec)
125619307e5cSMatthew G. Knepley {
125719307e5cSMatthew G. Knepley   PetscFunctionBegin;
125819307e5cSMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
125919307e5cSMatthew G. Knepley   PetscCall(DMSwarmDestroyVectorFromFields_Private(dm, Nf, fieldnames, vec));
126019307e5cSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
126119307e5cSMatthew G. Knepley }
126219307e5cSMatthew G. Knepley 
126319307e5cSMatthew G. Knepley /*@
126419307e5cSMatthew G. Knepley   DMSwarmCreateLocalVectorFromFields - Creates a `Vec` object sharing the array associated with a given field set
126519307e5cSMatthew G. Knepley 
126619307e5cSMatthew G. Knepley   Collective
126719307e5cSMatthew G. Knepley 
126819307e5cSMatthew G. Knepley   Input Parameters:
126919307e5cSMatthew G. Knepley + dm         - a `DMSWARM`
127019307e5cSMatthew G. Knepley . Nf         - the number of fields
127119307e5cSMatthew G. Knepley - fieldnames - the textual names given to the registered fields
127219307e5cSMatthew G. Knepley 
127319307e5cSMatthew G. Knepley   Output Parameter:
127419307e5cSMatthew G. Knepley . vec - the vector
127519307e5cSMatthew G. Knepley 
127619307e5cSMatthew G. Knepley   Level: beginner
127719307e5cSMatthew G. Knepley 
127819307e5cSMatthew G. Knepley   Notes:
127919307e5cSMatthew G. Knepley   The vector must be returned using a matching call to DMSwarmDestroyLocalVectorFromField().
128019307e5cSMatthew G. Knepley 
128119307e5cSMatthew G. Knepley   This vector is copyin-copyout, rather than a direct pointer like `DMSwarmCreateLocalVectorFromField()`
128219307e5cSMatthew G. Knepley 
128319307e5cSMatthew G. Knepley .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmDestroyLocalVectorFromField()`
128419307e5cSMatthew G. Knepley @*/
128519307e5cSMatthew G. Knepley PetscErrorCode DMSwarmCreateLocalVectorFromFields(DM dm, PetscInt Nf, const char *fieldnames[], Vec *vec)
128619307e5cSMatthew G. Knepley {
128719307e5cSMatthew G. Knepley   MPI_Comm comm = PETSC_COMM_SELF;
128819307e5cSMatthew G. Knepley 
128919307e5cSMatthew G. Knepley   PetscFunctionBegin;
129019307e5cSMatthew G. Knepley   PetscCall(DMSwarmCreateVectorFromFields_Private(dm, Nf, fieldnames, comm, vec));
129119307e5cSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
129219307e5cSMatthew G. Knepley }
129319307e5cSMatthew G. Knepley 
129419307e5cSMatthew G. Knepley /*@
129519307e5cSMatthew G. Knepley   DMSwarmDestroyLocalVectorFromFields - Destroys the `Vec` object which share the array associated with a given field set
129619307e5cSMatthew G. Knepley 
129719307e5cSMatthew G. Knepley   Collective
129819307e5cSMatthew G. Knepley 
129919307e5cSMatthew G. Knepley   Input Parameters:
130019307e5cSMatthew G. Knepley + dm         - a `DMSWARM`
130119307e5cSMatthew G. Knepley . Nf         - the number of fields
130219307e5cSMatthew G. Knepley - fieldnames - the textual names given to the registered fields
130319307e5cSMatthew G. Knepley 
130419307e5cSMatthew G. Knepley   Output Parameter:
130519307e5cSMatthew G. Knepley . vec - the vector
130619307e5cSMatthew G. Knepley 
130719307e5cSMatthew G. Knepley   Level: beginner
130819307e5cSMatthew G. Knepley 
130919307e5cSMatthew G. Knepley .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmCreateLocalVectorFromFields()`
131019307e5cSMatthew G. Knepley @*/
131119307e5cSMatthew G. Knepley PetscErrorCode DMSwarmDestroyLocalVectorFromFields(DM dm, PetscInt Nf, const char *fieldnames[], Vec *vec)
131219307e5cSMatthew G. Knepley {
131319307e5cSMatthew G. Knepley   PetscFunctionBegin;
131419307e5cSMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
131519307e5cSMatthew G. Knepley   PetscCall(DMSwarmDestroyVectorFromFields_Private(dm, Nf, fieldnames, vec));
131619307e5cSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
131719307e5cSMatthew G. Knepley }
131819307e5cSMatthew G. Knepley 
131919307e5cSMatthew G. Knepley /*@
132020f4b53cSBarry Smith   DMSwarmInitializeFieldRegister - Initiates the registration of fields to a `DMSWARM`
1321d3a51819SDave May 
132220f4b53cSBarry Smith   Collective
1323d3a51819SDave May 
132460225df5SJacob Faibussowitsch   Input Parameter:
132520f4b53cSBarry Smith . dm - a `DMSWARM`
1326d3a51819SDave May 
1327d3a51819SDave May   Level: beginner
1328d3a51819SDave May 
132920f4b53cSBarry Smith   Note:
133020f4b53cSBarry Smith   After all fields have been registered, you must call `DMSwarmFinalizeFieldRegister()`.
1331d3a51819SDave May 
133220f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmFinalizeFieldRegister()`, `DMSwarmRegisterPetscDatatypeField()`,
1333db781477SPatrick Sanan           `DMSwarmRegisterUserStructField()`, `DMSwarmRegisterUserDatatypeField()`
1334d3a51819SDave May @*/
1335d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmInitializeFieldRegister(DM dm)
1336d71ae5a4SJacob Faibussowitsch {
13375f50eb2eSDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
13383454631fSDave May 
1339521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1340cc651181SDave May   if (!swarm->field_registration_initialized) {
13415f50eb2eSDave May     swarm->field_registration_initialized = PETSC_TRUE;
1342da81f932SPierre Jolivet     PetscCall(DMSwarmRegisterPetscDatatypeField(dm, DMSwarmField_pid, 1, PETSC_INT64)); /* unique identifier */
13439566063dSJacob Faibussowitsch     PetscCall(DMSwarmRegisterPetscDatatypeField(dm, DMSwarmField_rank, 1, PETSC_INT));  /* used for communication */
1344cc651181SDave May   }
13453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13465f50eb2eSDave May }
13475f50eb2eSDave May 
134874d0cae8SMatthew G. Knepley /*@
134920f4b53cSBarry Smith   DMSwarmFinalizeFieldRegister - Finalizes the registration of fields to a `DMSWARM`
1350d3a51819SDave May 
135120f4b53cSBarry Smith   Collective
1352d3a51819SDave May 
135360225df5SJacob Faibussowitsch   Input Parameter:
135420f4b53cSBarry Smith . dm - a `DMSWARM`
1355d3a51819SDave May 
1356d3a51819SDave May   Level: beginner
1357d3a51819SDave May 
135820f4b53cSBarry Smith   Note:
135920f4b53cSBarry Smith   After `DMSwarmFinalizeFieldRegister()` has been called, no new fields can be defined on the `DMSWARM`.
1360d3a51819SDave May 
136120f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmInitializeFieldRegister()`, `DMSwarmRegisterPetscDatatypeField()`,
1362db781477SPatrick Sanan           `DMSwarmRegisterUserStructField()`, `DMSwarmRegisterUserDatatypeField()`
1363d3a51819SDave May @*/
1364d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmFinalizeFieldRegister(DM dm)
1365d71ae5a4SJacob Faibussowitsch {
13665f50eb2eSDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
13676845f8f5SDave May 
1368521f74f9SMatthew G. Knepley   PetscFunctionBegin;
136948a46eb9SPierre Jolivet   if (!swarm->field_registration_finalized) PetscCall(DMSwarmDataBucketFinalize(swarm->db));
1370f0cdbbbaSDave May   swarm->field_registration_finalized = PETSC_TRUE;
13713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13725f50eb2eSDave May }
13735f50eb2eSDave May 
137474d0cae8SMatthew G. Knepley /*@
137520f4b53cSBarry Smith   DMSwarmSetLocalSizes - Sets the length of all registered fields on the `DMSWARM`
1376d3a51819SDave May 
137720f4b53cSBarry Smith   Not Collective
1378d3a51819SDave May 
137960225df5SJacob Faibussowitsch   Input Parameters:
1380fca3fa51SMatthew G. Knepley + sw     - a `DMSWARM`
1381d3a51819SDave May . nlocal - the length of each registered field
138262741f57SDave May - buffer - the length of the buffer used to efficient dynamic re-sizing
1383d3a51819SDave May 
1384d3a51819SDave May   Level: beginner
1385d3a51819SDave May 
138620f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmGetLocalSize()`
1387d3a51819SDave May @*/
1388fca3fa51SMatthew G. Knepley PetscErrorCode DMSwarmSetLocalSizes(DM sw, PetscInt nlocal, PetscInt buffer)
1389d71ae5a4SJacob Faibussowitsch {
1390fca3fa51SMatthew G. Knepley   DM_Swarm   *swarm = (DM_Swarm *)sw->data;
1391fca3fa51SMatthew G. Knepley   PetscMPIInt rank;
1392fca3fa51SMatthew G. Knepley   PetscInt   *rankval;
13935f50eb2eSDave May 
1394521f74f9SMatthew G. Knepley   PetscFunctionBegin;
13959566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMSWARM_SetSizes, 0, 0, 0, 0));
13969566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketSetSizes(swarm->db, nlocal, buffer));
13979566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMSWARM_SetSizes, 0, 0, 0, 0));
1398fca3fa51SMatthew G. Knepley 
1399fca3fa51SMatthew G. Knepley   // Initialize values in pid and rank placeholders
1400fca3fa51SMatthew G. Knepley   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sw), &rank));
1401fca3fa51SMatthew G. Knepley   PetscCall(DMSwarmGetField(sw, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
1402fca3fa51SMatthew G. Knepley   for (PetscInt p = 0; p < nlocal; p++) rankval[p] = rank;
1403fca3fa51SMatthew G. Knepley   PetscCall(DMSwarmRestoreField(sw, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
1404fca3fa51SMatthew G. Knepley   /* TODO: [pid - use MPI_Scan] */
14053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14065f50eb2eSDave May }
14075f50eb2eSDave May 
140874d0cae8SMatthew G. Knepley /*@
140920f4b53cSBarry Smith   DMSwarmSetCellDM - Attaches a `DM` to a `DMSWARM`
1410d3a51819SDave May 
141120f4b53cSBarry Smith   Collective
1412d3a51819SDave May 
141360225df5SJacob Faibussowitsch   Input Parameters:
141419307e5cSMatthew G. Knepley + sw - a `DMSWARM`
141519307e5cSMatthew G. Knepley - dm - the `DM` to attach to the `DMSWARM`
1416d3a51819SDave May 
1417d3a51819SDave May   Level: beginner
1418d3a51819SDave May 
141920f4b53cSBarry Smith   Note:
142019307e5cSMatthew G. Knepley   The attached `DM` (dm) will be queried for point location and
142120f4b53cSBarry Smith   neighbor MPI-rank information if `DMSwarmMigrate()` is called.
1422d3a51819SDave May 
142320f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmSetType()`, `DMSwarmGetCellDM()`, `DMSwarmMigrate()`
1424d3a51819SDave May @*/
142519307e5cSMatthew G. Knepley PetscErrorCode DMSwarmSetCellDM(DM sw, DM dm)
1426d71ae5a4SJacob Faibussowitsch {
142719307e5cSMatthew G. Knepley   DMSwarmCellDM celldm;
142819307e5cSMatthew G. Knepley   const char   *name;
142919307e5cSMatthew G. Knepley   char         *coordName;
1430d52c2f21SMatthew G. Knepley 
1431d52c2f21SMatthew G. Knepley   PetscFunctionBegin;
1432d52c2f21SMatthew G. Knepley   PetscValidHeaderSpecific(sw, DM_CLASSID, 1);
143319307e5cSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
143419307e5cSMatthew G. Knepley   PetscCall(PetscStrallocpy(DMSwarmPICField_coor, &coordName));
143519307e5cSMatthew G. Knepley   PetscCall(DMSwarmCellDMCreate(dm, 0, NULL, 1, (const char **)&coordName, &celldm));
143619307e5cSMatthew G. Knepley   PetscCall(PetscFree(coordName));
143719307e5cSMatthew G. Knepley   PetscCall(PetscObjectGetName((PetscObject)celldm, &name));
143819307e5cSMatthew G. Knepley   PetscCall(DMSwarmAddCellDM(sw, celldm));
143919307e5cSMatthew G. Knepley   PetscCall(DMSwarmCellDMDestroy(&celldm));
144019307e5cSMatthew G. Knepley   PetscCall(DMSwarmSetCellDMActive(sw, name));
1441d52c2f21SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
1442d52c2f21SMatthew G. Knepley }
1443d52c2f21SMatthew G. Knepley 
1444d52c2f21SMatthew G. Knepley /*@
144519307e5cSMatthew G. Knepley   DMSwarmGetCellDM - Fetches the active cell `DM`
1446d52c2f21SMatthew G. Knepley 
1447d52c2f21SMatthew G. Knepley   Collective
1448d52c2f21SMatthew G. Knepley 
1449d52c2f21SMatthew G. Knepley   Input Parameter:
1450d52c2f21SMatthew G. Knepley . sw - a `DMSWARM`
1451d52c2f21SMatthew G. Knepley 
145219307e5cSMatthew G. Knepley   Output Parameter:
145319307e5cSMatthew G. Knepley . dm - the active `DM` for the `DMSWARM`
145419307e5cSMatthew G. Knepley 
1455d52c2f21SMatthew G. Knepley   Level: beginner
1456d52c2f21SMatthew G. Knepley 
145719307e5cSMatthew G. Knepley .seealso: `DM`, `DMSWARM`, `DMSwarmSetCellDM()`
1458d52c2f21SMatthew G. Knepley @*/
145919307e5cSMatthew G. Knepley PetscErrorCode DMSwarmGetCellDM(DM sw, DM *dm)
1460d52c2f21SMatthew G. Knepley {
1461d52c2f21SMatthew G. Knepley   DM_Swarm     *swarm = (DM_Swarm *)sw->data;
146219307e5cSMatthew G. Knepley   DMSwarmCellDM celldm;
146319307e5cSMatthew G. Knepley 
146419307e5cSMatthew G. Knepley   PetscFunctionBegin;
1465fca3fa51SMatthew G. Knepley   PetscValidHeaderSpecific(sw, DM_CLASSID, 1);
146619307e5cSMatthew G. Knepley   PetscCall(PetscObjectListFind(swarm->cellDMs, swarm->activeCellDM, (PetscObject *)&celldm));
146719307e5cSMatthew G. Knepley   PetscCheck(celldm, PetscObjectComm((PetscObject)sw), PETSC_ERR_ARG_WRONG, "There is no cell DM named %s in this Swarm", swarm->activeCellDM);
146819307e5cSMatthew G. Knepley   PetscCall(DMSwarmCellDMGetDM(celldm, dm));
146919307e5cSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
147019307e5cSMatthew G. Knepley }
147119307e5cSMatthew G. Knepley 
1472fca3fa51SMatthew G. Knepley /*@C
1473fca3fa51SMatthew G. Knepley   DMSwarmGetCellDMNames - Get the list of cell `DM` names
1474fca3fa51SMatthew G. Knepley 
1475fca3fa51SMatthew G. Knepley   Not collective
1476fca3fa51SMatthew G. Knepley 
1477fca3fa51SMatthew G. Knepley   Input Parameter:
1478fca3fa51SMatthew G. Knepley . sw - a `DMSWARM`
1479fca3fa51SMatthew G. Knepley 
1480fca3fa51SMatthew G. Knepley   Output Parameters:
1481fca3fa51SMatthew G. Knepley + Ndm     - the number of `DMSwarmCellDM` in the `DMSWARM`
1482fca3fa51SMatthew G. Knepley - celldms - the name of each `DMSwarmCellDM`
1483fca3fa51SMatthew G. Knepley 
1484fca3fa51SMatthew G. Knepley   Level: beginner
1485fca3fa51SMatthew G. Knepley 
1486fca3fa51SMatthew G. Knepley .seealso: `DM`, `DMSWARM`, `DMSwarmSetCellDM()`, `DMSwarmGetCellDMByName()`
1487fca3fa51SMatthew G. Knepley @*/
1488fca3fa51SMatthew G. Knepley PetscErrorCode DMSwarmGetCellDMNames(DM sw, PetscInt *Ndm, const char **celldms[])
1489fca3fa51SMatthew G. Knepley {
1490fca3fa51SMatthew G. Knepley   DM_Swarm       *swarm = (DM_Swarm *)sw->data;
1491fca3fa51SMatthew G. Knepley   PetscObjectList next  = swarm->cellDMs;
1492fca3fa51SMatthew G. Knepley   PetscInt        n     = 0;
1493fca3fa51SMatthew G. Knepley 
1494fca3fa51SMatthew G. Knepley   PetscFunctionBegin;
1495fca3fa51SMatthew G. Knepley   PetscValidHeaderSpecific(sw, DM_CLASSID, 1);
1496fca3fa51SMatthew G. Knepley   PetscAssertPointer(Ndm, 2);
1497fca3fa51SMatthew G. Knepley   PetscAssertPointer(celldms, 3);
1498fca3fa51SMatthew G. Knepley   while (next) {
1499fca3fa51SMatthew G. Knepley     next = next->next;
1500fca3fa51SMatthew G. Knepley     ++n;
1501fca3fa51SMatthew G. Knepley   }
1502fca3fa51SMatthew G. Knepley   PetscCall(PetscMalloc1(n, celldms));
1503fca3fa51SMatthew G. Knepley   next = swarm->cellDMs;
1504fca3fa51SMatthew G. Knepley   n    = 0;
1505fca3fa51SMatthew G. Knepley   while (next) {
1506fca3fa51SMatthew G. Knepley     (*celldms)[n] = (const char *)next->obj->name;
1507fca3fa51SMatthew G. Knepley     next          = next->next;
1508fca3fa51SMatthew G. Knepley     ++n;
1509fca3fa51SMatthew G. Knepley   }
1510fca3fa51SMatthew G. Knepley   *Ndm = n;
1511fca3fa51SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
1512fca3fa51SMatthew G. Knepley }
1513fca3fa51SMatthew G. Knepley 
151419307e5cSMatthew G. Knepley /*@
151519307e5cSMatthew G. Knepley   DMSwarmSetCellDMActive - Activates a cell `DM` for a `DMSWARM`
151619307e5cSMatthew G. Knepley 
151719307e5cSMatthew G. Knepley   Collective
151819307e5cSMatthew G. Knepley 
151919307e5cSMatthew G. Knepley   Input Parameters:
152019307e5cSMatthew G. Knepley + sw   - a `DMSWARM`
152119307e5cSMatthew G. Knepley - name - name of the cell `DM` to active for the `DMSWARM`
152219307e5cSMatthew G. Knepley 
152319307e5cSMatthew G. Knepley   Level: beginner
152419307e5cSMatthew G. Knepley 
152519307e5cSMatthew G. Knepley   Note:
152619307e5cSMatthew G. Knepley   The attached `DM` (dmcell) will be queried for point location and
152719307e5cSMatthew G. Knepley   neighbor MPI-rank information if `DMSwarmMigrate()` is called.
152819307e5cSMatthew G. Knepley 
152919307e5cSMatthew G. Knepley .seealso: `DM`, `DMSWARM`, `DMSwarmCellDM`, `DMSwarmSetType()`, `DMSwarmAddCellDM()`, `DMSwarmSetCellDM()`, `DMSwarmMigrate()`
153019307e5cSMatthew G. Knepley @*/
153119307e5cSMatthew G. Knepley PetscErrorCode DMSwarmSetCellDMActive(DM sw, const char name[])
153219307e5cSMatthew G. Knepley {
153319307e5cSMatthew G. Knepley   DM_Swarm     *swarm = (DM_Swarm *)sw->data;
153419307e5cSMatthew G. Knepley   DMSwarmCellDM celldm;
1535d52c2f21SMatthew G. Knepley 
1536d52c2f21SMatthew G. Knepley   PetscFunctionBegin;
1537d52c2f21SMatthew G. Knepley   PetscValidHeaderSpecific(sw, DM_CLASSID, 1);
153819307e5cSMatthew G. Knepley   PetscCall(PetscInfo(sw, "Setting cell DM to %s\n", name));
153919307e5cSMatthew G. Knepley   PetscCall(PetscFree(swarm->activeCellDM));
154019307e5cSMatthew G. Knepley   PetscCall(PetscStrallocpy(name, (char **)&swarm->activeCellDM));
154119307e5cSMatthew G. Knepley   PetscCall(DMSwarmGetCellDMActive(sw, &celldm));
154219307e5cSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
1543d52c2f21SMatthew G. Knepley }
154419307e5cSMatthew G. Knepley 
154519307e5cSMatthew G. Knepley /*@
154619307e5cSMatthew G. Knepley   DMSwarmGetCellDMActive - Returns the active cell `DM` for a `DMSWARM`
154719307e5cSMatthew G. Knepley 
154819307e5cSMatthew G. Knepley   Collective
154919307e5cSMatthew G. Knepley 
155019307e5cSMatthew G. Knepley   Input Parameter:
155119307e5cSMatthew G. Knepley . sw - a `DMSWARM`
155219307e5cSMatthew G. Knepley 
155319307e5cSMatthew G. Knepley   Output Parameter:
155419307e5cSMatthew G. Knepley . celldm - the active `DMSwarmCellDM`
155519307e5cSMatthew G. Knepley 
155619307e5cSMatthew G. Knepley   Level: beginner
155719307e5cSMatthew G. Knepley 
155819307e5cSMatthew G. Knepley .seealso: `DM`, `DMSWARM`, `DMSwarmCellDM`, `DMSwarmSetType()`, `DMSwarmAddCellDM()`, `DMSwarmSetCellDM()`, `DMSwarmMigrate()`
155919307e5cSMatthew G. Knepley @*/
156019307e5cSMatthew G. Knepley PetscErrorCode DMSwarmGetCellDMActive(DM sw, DMSwarmCellDM *celldm)
156119307e5cSMatthew G. Knepley {
156219307e5cSMatthew G. Knepley   DM_Swarm *swarm = (DM_Swarm *)sw->data;
156319307e5cSMatthew G. Knepley 
156419307e5cSMatthew G. Knepley   PetscFunctionBegin;
156519307e5cSMatthew G. Knepley   PetscValidHeaderSpecific(sw, DM_CLASSID, 1);
1566fca3fa51SMatthew G. Knepley   PetscAssertPointer(celldm, 2);
156719307e5cSMatthew G. Knepley   PetscCheck(swarm->activeCellDM, PetscObjectComm((PetscObject)sw), PETSC_ERR_ARG_WRONGSTATE, "Swarm has no active cell DM");
156819307e5cSMatthew G. Knepley   PetscCall(PetscObjectListFind(swarm->cellDMs, swarm->activeCellDM, (PetscObject *)celldm));
1569fca3fa51SMatthew G. Knepley   PetscCheck(*celldm, PetscObjectComm((PetscObject)sw), PETSC_ERR_ARG_WRONGSTATE, "Swarm has no valid cell DM for %s", swarm->activeCellDM);
1570fca3fa51SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
1571fca3fa51SMatthew G. Knepley }
1572fca3fa51SMatthew G. Knepley 
1573fca3fa51SMatthew G. Knepley /*@C
1574fca3fa51SMatthew G. Knepley   DMSwarmGetCellDMByName - Get a `DMSwarmCellDM` from its name
1575fca3fa51SMatthew G. Knepley 
1576fca3fa51SMatthew G. Knepley   Not collective
1577fca3fa51SMatthew G. Knepley 
1578fca3fa51SMatthew G. Knepley   Input Parameters:
1579fca3fa51SMatthew G. Knepley + sw   - a `DMSWARM`
1580fca3fa51SMatthew G. Knepley - name - the name
1581fca3fa51SMatthew G. Knepley 
1582fca3fa51SMatthew G. Knepley   Output Parameter:
1583fca3fa51SMatthew G. Knepley . celldm - the `DMSwarmCellDM`
1584fca3fa51SMatthew G. Knepley 
1585fca3fa51SMatthew G. Knepley   Level: beginner
1586fca3fa51SMatthew G. Knepley 
1587fca3fa51SMatthew G. Knepley .seealso: `DM`, `DMSWARM`, `DMSwarmSetCellDM()`, `DMSwarmGetCellDMNames()`
1588fca3fa51SMatthew G. Knepley @*/
1589fca3fa51SMatthew G. Knepley PetscErrorCode DMSwarmGetCellDMByName(DM sw, const char name[], DMSwarmCellDM *celldm)
1590fca3fa51SMatthew G. Knepley {
1591fca3fa51SMatthew G. Knepley   DM_Swarm *swarm = (DM_Swarm *)sw->data;
1592fca3fa51SMatthew G. Knepley 
1593fca3fa51SMatthew G. Knepley   PetscFunctionBegin;
1594fca3fa51SMatthew G. Knepley   PetscValidHeaderSpecific(sw, DM_CLASSID, 1);
1595fca3fa51SMatthew G. Knepley   PetscAssertPointer(name, 2);
1596fca3fa51SMatthew G. Knepley   PetscAssertPointer(celldm, 3);
1597fca3fa51SMatthew G. Knepley   PetscCall(PetscObjectListFind(swarm->cellDMs, name, (PetscObject *)celldm));
1598fca3fa51SMatthew G. Knepley   PetscCheck(*celldm, PetscObjectComm((PetscObject)sw), PETSC_ERR_ARG_WRONGSTATE, "Swarm has no valid cell DM for %s", name);
159919307e5cSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
160019307e5cSMatthew G. Knepley }
160119307e5cSMatthew G. Knepley 
160219307e5cSMatthew G. Knepley /*@
160319307e5cSMatthew G. Knepley   DMSwarmAddCellDM - Adds a cell `DM` to the `DMSWARM`
160419307e5cSMatthew G. Knepley 
160519307e5cSMatthew G. Knepley   Collective
160619307e5cSMatthew G. Knepley 
160719307e5cSMatthew G. Knepley   Input Parameters:
160819307e5cSMatthew G. Knepley + sw     - a `DMSWARM`
160919307e5cSMatthew G. Knepley - celldm - the `DMSwarmCellDM`
161019307e5cSMatthew G. Knepley 
161119307e5cSMatthew G. Knepley   Level: beginner
161219307e5cSMatthew G. Knepley 
161319307e5cSMatthew G. Knepley   Note:
161419307e5cSMatthew G. Knepley   Cell DMs with the same name will share the cellid field
161519307e5cSMatthew G. Knepley 
161619307e5cSMatthew G. Knepley .seealso: `DM`, `DMSWARM`, `DMSwarmSetType()`, `DMSwarmPushCellDM()`, `DMSwarmSetCellDM()`, `DMSwarmMigrate()`
161719307e5cSMatthew G. Knepley @*/
161819307e5cSMatthew G. Knepley PetscErrorCode DMSwarmAddCellDM(DM sw, DMSwarmCellDM celldm)
161919307e5cSMatthew G. Knepley {
162019307e5cSMatthew G. Knepley   DM_Swarm   *swarm = (DM_Swarm *)sw->data;
162119307e5cSMatthew G. Knepley   const char *name;
162219307e5cSMatthew G. Knepley   PetscInt    dim;
162319307e5cSMatthew G. Knepley   PetscBool   flg;
162419307e5cSMatthew G. Knepley   MPI_Comm    comm;
162519307e5cSMatthew G. Knepley 
162619307e5cSMatthew G. Knepley   PetscFunctionBegin;
162719307e5cSMatthew G. Knepley   PetscValidHeaderSpecific(sw, DM_CLASSID, 1);
162819307e5cSMatthew G. Knepley   PetscCall(PetscObjectGetComm((PetscObject)sw, &comm));
162919307e5cSMatthew G. Knepley   PetscValidHeaderSpecific(celldm, DMSWARMCELLDM_CLASSID, 2);
163019307e5cSMatthew G. Knepley   PetscCall(PetscObjectGetName((PetscObject)celldm, &name));
163119307e5cSMatthew G. Knepley   PetscCall(PetscObjectListAdd(&swarm->cellDMs, name, (PetscObject)celldm));
163219307e5cSMatthew G. Knepley   PetscCall(DMGetDimension(sw, &dim));
163319307e5cSMatthew G. Knepley   for (PetscInt f = 0; f < celldm->Nfc; ++f) {
163419307e5cSMatthew G. Knepley     PetscCall(DMSwarmDataFieldStringInList(celldm->coordFields[f], swarm->db->nfields, (const DMSwarmDataField *)swarm->db->field, &flg));
163519307e5cSMatthew G. Knepley     if (!flg) {
163619307e5cSMatthew G. Knepley       PetscCall(DMSwarmRegisterPetscDatatypeField(sw, celldm->coordFields[f], dim, PETSC_DOUBLE));
163719307e5cSMatthew G. Knepley     } else {
163819307e5cSMatthew G. Knepley       PetscDataType dt;
163919307e5cSMatthew G. Knepley       PetscInt      bs;
164019307e5cSMatthew G. Knepley 
164119307e5cSMatthew G. Knepley       PetscCall(DMSwarmGetFieldInfo(sw, celldm->coordFields[f], &bs, &dt));
164219307e5cSMatthew G. Knepley       PetscCheck(bs == dim, comm, PETSC_ERR_ARG_WRONG, "Coordinate field %s has blocksize %" PetscInt_FMT " != %" PetscInt_FMT " spatial dimension", celldm->coordFields[f], bs, dim);
164319307e5cSMatthew G. Knepley       PetscCheck(dt == PETSC_DOUBLE, comm, PETSC_ERR_ARG_WRONG, "Coordinate field %s has datatype %s != PETSC_DOUBLE", celldm->coordFields[f], PetscDataTypes[dt]);
164419307e5cSMatthew G. Knepley     }
164519307e5cSMatthew G. Knepley   }
164619307e5cSMatthew G. Knepley   // Assume that DMs with the same name share the cellid field
164719307e5cSMatthew G. Knepley   PetscCall(DMSwarmDataFieldStringInList(celldm->cellid, swarm->db->nfields, (const DMSwarmDataField *)swarm->db->field, &flg));
164819307e5cSMatthew G. Knepley   if (!flg) {
164919307e5cSMatthew G. Knepley     PetscBool   isShell, isDummy;
165019307e5cSMatthew G. Knepley     const char *name;
165119307e5cSMatthew G. Knepley 
165219307e5cSMatthew G. Knepley     // Allow dummy DMSHELL (I don't think we should support this mode)
165319307e5cSMatthew G. Knepley     PetscCall(PetscObjectTypeCompare((PetscObject)celldm->dm, DMSHELL, &isShell));
165419307e5cSMatthew G. Knepley     PetscCall(PetscObjectGetName((PetscObject)celldm->dm, &name));
165519307e5cSMatthew G. Knepley     PetscCall(PetscStrcmp(name, "dummy", &isDummy));
165619307e5cSMatthew G. Knepley     if (!isShell || !isDummy) PetscCall(DMSwarmRegisterPetscDatatypeField(sw, celldm->cellid, 1, PETSC_INT));
165719307e5cSMatthew G. Knepley   }
165819307e5cSMatthew G. Knepley   PetscCall(DMSwarmSetCellDMActive(sw, name));
16593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1660fe39f135SDave May }
1661fe39f135SDave May 
166274d0cae8SMatthew G. Knepley /*@
1663d5b43468SJose E. Roman   DMSwarmGetLocalSize - Retrieves the local length of fields registered
1664d3a51819SDave May 
166520f4b53cSBarry Smith   Not Collective
1666d3a51819SDave May 
166760225df5SJacob Faibussowitsch   Input Parameter:
166820f4b53cSBarry Smith . dm - a `DMSWARM`
1669d3a51819SDave May 
167060225df5SJacob Faibussowitsch   Output Parameter:
1671d3a51819SDave May . nlocal - the length of each registered field
1672d3a51819SDave May 
1673d3a51819SDave May   Level: beginner
1674d3a51819SDave May 
167520f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmGetSize()`, `DMSwarmSetLocalSizes()`
1676d3a51819SDave May @*/
1677d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmGetLocalSize(DM dm, PetscInt *nlocal)
1678d71ae5a4SJacob Faibussowitsch {
1679dcf43ee8SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1680dcf43ee8SDave May 
1681521f74f9SMatthew G. Knepley   PetscFunctionBegin;
16829566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, nlocal, NULL, NULL));
16833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1684dcf43ee8SDave May }
1685dcf43ee8SDave May 
168674d0cae8SMatthew G. Knepley /*@
1687d5b43468SJose E. Roman   DMSwarmGetSize - Retrieves the total length of fields registered
1688d3a51819SDave May 
168920f4b53cSBarry Smith   Collective
1690d3a51819SDave May 
169160225df5SJacob Faibussowitsch   Input Parameter:
169220f4b53cSBarry Smith . dm - a `DMSWARM`
1693d3a51819SDave May 
169460225df5SJacob Faibussowitsch   Output Parameter:
1695d3a51819SDave May . n - the total length of each registered field
1696d3a51819SDave May 
1697d3a51819SDave May   Level: beginner
1698d3a51819SDave May 
1699d3a51819SDave May   Note:
170020f4b53cSBarry Smith   This calls `MPI_Allreduce()` upon each call (inefficient but safe)
1701d3a51819SDave May 
170220f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmGetLocalSize()`, `DMSwarmSetLocalSizes()`
1703d3a51819SDave May @*/
1704d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmGetSize(DM dm, PetscInt *n)
1705d71ae5a4SJacob Faibussowitsch {
1706dcf43ee8SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
17075627991aSBarry Smith   PetscInt  nlocal;
1708dcf43ee8SDave May 
1709521f74f9SMatthew G. Knepley   PetscFunctionBegin;
17109566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &nlocal, NULL, NULL));
1711462c564dSBarry Smith   PetscCallMPI(MPIU_Allreduce(&nlocal, n, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm)));
17123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1713dcf43ee8SDave May }
1714dcf43ee8SDave May 
1715ce78bad3SBarry Smith /*@C
171620f4b53cSBarry Smith   DMSwarmRegisterPetscDatatypeField - Register a field to a `DMSWARM` with a native PETSc data type
1717d3a51819SDave May 
171820f4b53cSBarry Smith   Collective
1719d3a51819SDave May 
172060225df5SJacob Faibussowitsch   Input Parameters:
172120f4b53cSBarry Smith + dm        - a `DMSWARM`
1722d3a51819SDave May . fieldname - the textual name to identify this field
1723d3a51819SDave May . blocksize - the number of each data type
172420f4b53cSBarry Smith - type      - a valid PETSc data type (`PETSC_CHAR`, `PETSC_SHORT`, `PETSC_INT`, `PETSC_FLOAT`, `PETSC_REAL`, `PETSC_LONG`)
1725d3a51819SDave May 
1726d3a51819SDave May   Level: beginner
1727d3a51819SDave May 
1728d3a51819SDave May   Notes:
17298b8a3813SDave May   The textual name for each registered field must be unique.
1730d3a51819SDave May 
173120f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterUserStructField()`, `DMSwarmRegisterUserDatatypeField()`
1732d3a51819SDave May @*/
1733d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmRegisterPetscDatatypeField(DM dm, const char fieldname[], PetscInt blocksize, PetscDataType type)
1734d71ae5a4SJacob Faibussowitsch {
1735b62e03f8SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1736b62e03f8SDave May   size_t    size;
1737b62e03f8SDave May 
1738521f74f9SMatthew G. Knepley   PetscFunctionBegin;
173928b400f6SJacob Faibussowitsch   PetscCheck(swarm->field_registration_initialized, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Must call DMSwarmInitializeFieldRegister() first");
174028b400f6SJacob Faibussowitsch   PetscCheck(!swarm->field_registration_finalized, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Cannot register additional fields after calling DMSwarmFinalizeFieldRegister() first");
17415f50eb2eSDave May 
174208401ef6SPierre Jolivet   PetscCheck(type != PETSC_OBJECT, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}");
174308401ef6SPierre Jolivet   PetscCheck(type != PETSC_FUNCTION, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}");
174408401ef6SPierre Jolivet   PetscCheck(type != PETSC_STRING, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}");
174508401ef6SPierre Jolivet   PetscCheck(type != PETSC_STRUCT, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}");
174608401ef6SPierre Jolivet   PetscCheck(type != PETSC_DATATYPE_UNKNOWN, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}");
1747b62e03f8SDave May 
17489566063dSJacob Faibussowitsch   PetscCall(PetscDataTypeGetSize(type, &size));
1749b62e03f8SDave May   /* Load a specific data type into data bucket, specifying textual name and its size in bytes */
17509566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketRegisterField(swarm->db, "DMSwarmRegisterPetscDatatypeField", fieldname, blocksize * size, NULL));
175152c3ed93SDave May   {
175277048351SPatrick Sanan     DMSwarmDataField gfield;
175352c3ed93SDave May 
17549566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield));
17559566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataFieldSetBlockSize(gfield, blocksize));
175652c3ed93SDave May   }
1757b62e03f8SDave May   swarm->db->field[swarm->db->nfields - 1]->petsc_type = type;
17583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1759b62e03f8SDave May }
1760b62e03f8SDave May 
1761d3a51819SDave May /*@C
176220f4b53cSBarry Smith   DMSwarmRegisterUserStructField - Register a user defined struct to a `DMSWARM`
1763d3a51819SDave May 
176420f4b53cSBarry Smith   Collective
1765d3a51819SDave May 
176660225df5SJacob Faibussowitsch   Input Parameters:
176720f4b53cSBarry Smith + dm        - a `DMSWARM`
1768d3a51819SDave May . fieldname - the textual name to identify this field
176962741f57SDave May - size      - the size in bytes of the user struct of each data type
1770d3a51819SDave May 
1771d3a51819SDave May   Level: beginner
1772d3a51819SDave May 
177320f4b53cSBarry Smith   Note:
17748b8a3813SDave May   The textual name for each registered field must be unique.
1775d3a51819SDave May 
177620f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmRegisterUserDatatypeField()`
1777d3a51819SDave May @*/
1778d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmRegisterUserStructField(DM dm, const char fieldname[], size_t size)
1779d71ae5a4SJacob Faibussowitsch {
1780b62e03f8SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1781b62e03f8SDave May 
1782521f74f9SMatthew G. Knepley   PetscFunctionBegin;
17839566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketRegisterField(swarm->db, "DMSwarmRegisterUserStructField", fieldname, size, NULL));
1784b62e03f8SDave May   swarm->db->field[swarm->db->nfields - 1]->petsc_type = PETSC_STRUCT;
17853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1786b62e03f8SDave May }
1787b62e03f8SDave May 
1788d3a51819SDave May /*@C
178920f4b53cSBarry Smith   DMSwarmRegisterUserDatatypeField - Register a user defined data type to a `DMSWARM`
1790d3a51819SDave May 
179120f4b53cSBarry Smith   Collective
1792d3a51819SDave May 
179360225df5SJacob Faibussowitsch   Input Parameters:
179420f4b53cSBarry Smith + dm        - a `DMSWARM`
1795d3a51819SDave May . fieldname - the textual name to identify this field
1796d3a51819SDave May . size      - the size in bytes of the user data type
179762741f57SDave May - blocksize - the number of each data type
1798d3a51819SDave May 
1799d3a51819SDave May   Level: beginner
1800d3a51819SDave May 
180120f4b53cSBarry Smith   Note:
18028b8a3813SDave May   The textual name for each registered field must be unique.
1803d3a51819SDave May 
180442747ad1SJacob Faibussowitsch .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmRegisterUserStructField()`
1805d3a51819SDave May @*/
1806d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmRegisterUserDatatypeField(DM dm, const char fieldname[], size_t size, PetscInt blocksize)
1807d71ae5a4SJacob Faibussowitsch {
1808b62e03f8SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1809b62e03f8SDave May 
1810521f74f9SMatthew G. Knepley   PetscFunctionBegin;
18119566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketRegisterField(swarm->db, "DMSwarmRegisterUserDatatypeField", fieldname, blocksize * size, NULL));
1812320740a0SDave May   {
181377048351SPatrick Sanan     DMSwarmDataField gfield;
1814320740a0SDave May 
18159566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield));
18169566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataFieldSetBlockSize(gfield, blocksize));
1817320740a0SDave May   }
1818b62e03f8SDave May   swarm->db->field[swarm->db->nfields - 1]->petsc_type = PETSC_DATATYPE_UNKNOWN;
18193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1820b62e03f8SDave May }
1821b62e03f8SDave May 
1822d3a51819SDave May /*@C
1823d3a51819SDave May   DMSwarmGetField - Get access to the underlying array storing all entries associated with a registered field
1824d3a51819SDave May 
1825cc4c1da9SBarry Smith   Not Collective, No Fortran Support
1826d3a51819SDave May 
182760225df5SJacob Faibussowitsch   Input Parameters:
182820f4b53cSBarry Smith + dm        - a `DMSWARM`
182962741f57SDave May - fieldname - the textual name to identify this field
1830d3a51819SDave May 
183160225df5SJacob Faibussowitsch   Output Parameters:
183262741f57SDave May + blocksize - the number of each data type
1833d3a51819SDave May . type      - the data type
183462741f57SDave May - data      - pointer to raw array
1835d3a51819SDave May 
1836d3a51819SDave May   Level: beginner
1837d3a51819SDave May 
1838d3a51819SDave May   Notes:
183920f4b53cSBarry Smith   The array must be returned using a matching call to `DMSwarmRestoreField()`.
1840d3a51819SDave May 
1841ce78bad3SBarry Smith   Fortran Note:
1842ce78bad3SBarry Smith   Only works for `type` of `PETSC_SCALAR`
1843ce78bad3SBarry Smith 
184420f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRestoreField()`
1845d3a51819SDave May @*/
1846ce78bad3SBarry Smith PetscErrorCode DMSwarmGetField(DM dm, const char fieldname[], PetscInt *blocksize, PetscDataType *type, void **data) PeNS
1847d71ae5a4SJacob Faibussowitsch {
1848b62e03f8SDave May   DM_Swarm        *swarm = (DM_Swarm *)dm->data;
184977048351SPatrick Sanan   DMSwarmDataField gfield;
1850b62e03f8SDave May 
1851521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1852ea7032a0SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
18539566063dSJacob Faibussowitsch   if (!swarm->issetup) PetscCall(DMSetUp(dm));
18549566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield));
18559566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataFieldGetAccess(gfield));
18569566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataFieldGetEntries(gfield, data));
1857ad540459SPierre Jolivet   if (blocksize) *blocksize = gfield->bs;
1858ad540459SPierre Jolivet   if (type) *type = gfield->petsc_type;
18593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1860b62e03f8SDave May }
1861b62e03f8SDave May 
1862d3a51819SDave May /*@C
1863d3a51819SDave May   DMSwarmRestoreField - Restore access to the underlying array storing all entries associated with a registered field
1864d3a51819SDave May 
1865ce78bad3SBarry Smith   Not Collective
1866d3a51819SDave May 
186760225df5SJacob Faibussowitsch   Input Parameters:
186820f4b53cSBarry Smith + dm        - a `DMSWARM`
186962741f57SDave May - fieldname - the textual name to identify this field
1870d3a51819SDave May 
187160225df5SJacob Faibussowitsch   Output Parameters:
187262741f57SDave May + blocksize - the number of each data type
1873d3a51819SDave May . type      - the data type
187462741f57SDave May - data      - pointer to raw array
1875d3a51819SDave May 
1876d3a51819SDave May   Level: beginner
1877d3a51819SDave May 
1878d3a51819SDave May   Notes:
187920f4b53cSBarry Smith   The user must call `DMSwarmGetField()` prior to calling `DMSwarmRestoreField()`.
1880d3a51819SDave May 
1881ce78bad3SBarry Smith   Fortran Note:
1882ce78bad3SBarry Smith   Only works for `type` of `PETSC_SCALAR`
1883ce78bad3SBarry Smith 
188420f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmGetField()`
1885d3a51819SDave May @*/
1886ce78bad3SBarry Smith PetscErrorCode DMSwarmRestoreField(DM dm, const char fieldname[], PetscInt *blocksize, PetscDataType *type, void **data) PeNS
1887d71ae5a4SJacob Faibussowitsch {
1888b62e03f8SDave May   DM_Swarm        *swarm = (DM_Swarm *)dm->data;
188977048351SPatrick Sanan   DMSwarmDataField gfield;
1890b62e03f8SDave May 
1891521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1892ea7032a0SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
18939566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield));
18949566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataFieldRestoreAccess(gfield));
1895b62e03f8SDave May   if (data) *data = NULL;
18963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1897b62e03f8SDave May }
1898b62e03f8SDave May 
18990bf7c1c5SMatthew G. Knepley PetscErrorCode DMSwarmGetFieldInfo(DM dm, const char fieldname[], PetscInt *blocksize, PetscDataType *type)
19000bf7c1c5SMatthew G. Knepley {
19010bf7c1c5SMatthew G. Knepley   DM_Swarm        *swarm = (DM_Swarm *)dm->data;
19020bf7c1c5SMatthew G. Knepley   DMSwarmDataField gfield;
19030bf7c1c5SMatthew G. Knepley 
19040bf7c1c5SMatthew G. Knepley   PetscFunctionBegin;
19050bf7c1c5SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
19060bf7c1c5SMatthew G. Knepley   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield));
19070bf7c1c5SMatthew G. Knepley   if (blocksize) *blocksize = gfield->bs;
19080bf7c1c5SMatthew G. Knepley   if (type) *type = gfield->petsc_type;
19090bf7c1c5SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
19100bf7c1c5SMatthew G. Knepley }
19110bf7c1c5SMatthew G. Knepley 
191274d0cae8SMatthew G. Knepley /*@
191320f4b53cSBarry Smith   DMSwarmAddPoint - Add space for one new point in the `DMSWARM`
1914d3a51819SDave May 
191520f4b53cSBarry Smith   Not Collective
1916d3a51819SDave May 
191760225df5SJacob Faibussowitsch   Input Parameter:
191820f4b53cSBarry Smith . dm - a `DMSWARM`
1919d3a51819SDave May 
1920d3a51819SDave May   Level: beginner
1921d3a51819SDave May 
1922d3a51819SDave May   Notes:
19238b8a3813SDave May   The new point will have all fields initialized to zero.
1924d3a51819SDave May 
192520f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmAddNPoints()`
1926d3a51819SDave May @*/
1927d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmAddPoint(DM dm)
1928d71ae5a4SJacob Faibussowitsch {
1929cb1d1399SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1930cb1d1399SDave May 
1931521f74f9SMatthew G. Knepley   PetscFunctionBegin;
19329566063dSJacob Faibussowitsch   if (!swarm->issetup) PetscCall(DMSetUp(dm));
19339566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMSWARM_AddPoints, 0, 0, 0, 0));
19349566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketAddPoint(swarm->db));
19359566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMSWARM_AddPoints, 0, 0, 0, 0));
19363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1937cb1d1399SDave May }
1938cb1d1399SDave May 
193974d0cae8SMatthew G. Knepley /*@
194020f4b53cSBarry Smith   DMSwarmAddNPoints - Add space for a number of new points in the `DMSWARM`
1941d3a51819SDave May 
194220f4b53cSBarry Smith   Not Collective
1943d3a51819SDave May 
194460225df5SJacob Faibussowitsch   Input Parameters:
194520f4b53cSBarry Smith + dm      - a `DMSWARM`
194662741f57SDave May - npoints - the number of new points to add
1947d3a51819SDave May 
1948d3a51819SDave May   Level: beginner
1949d3a51819SDave May 
1950d3a51819SDave May   Notes:
19518b8a3813SDave May   The new point will have all fields initialized to zero.
1952d3a51819SDave May 
195320f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmAddPoint()`
1954d3a51819SDave May @*/
1955d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmAddNPoints(DM dm, PetscInt npoints)
1956d71ae5a4SJacob Faibussowitsch {
1957cb1d1399SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1958cb1d1399SDave May   PetscInt  nlocal;
1959cb1d1399SDave May 
1960521f74f9SMatthew G. Knepley   PetscFunctionBegin;
19619566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMSWARM_AddPoints, 0, 0, 0, 0));
19629566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &nlocal, NULL, NULL));
1963670292f5SMatthew Knepley   nlocal = PetscMax(nlocal, 0) + npoints;
19649566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketSetSizes(swarm->db, nlocal, DMSWARM_DATA_BUCKET_BUFFER_DEFAULT));
19659566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMSWARM_AddPoints, 0, 0, 0, 0));
19663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1967cb1d1399SDave May }
1968cb1d1399SDave May 
196974d0cae8SMatthew G. Knepley /*@
197020f4b53cSBarry Smith   DMSwarmRemovePoint - Remove the last point from the `DMSWARM`
1971d3a51819SDave May 
197220f4b53cSBarry Smith   Not Collective
1973d3a51819SDave May 
197460225df5SJacob Faibussowitsch   Input Parameter:
197520f4b53cSBarry Smith . dm - a `DMSWARM`
1976d3a51819SDave May 
1977d3a51819SDave May   Level: beginner
1978d3a51819SDave May 
197920f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRemovePointAtIndex()`
1980d3a51819SDave May @*/
1981d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmRemovePoint(DM dm)
1982d71ae5a4SJacob Faibussowitsch {
1983cb1d1399SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1984cb1d1399SDave May 
1985521f74f9SMatthew G. Knepley   PetscFunctionBegin;
19869566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMSWARM_RemovePoints, 0, 0, 0, 0));
19879566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketRemovePoint(swarm->db));
19889566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMSWARM_RemovePoints, 0, 0, 0, 0));
19893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1990cb1d1399SDave May }
1991cb1d1399SDave May 
199274d0cae8SMatthew G. Knepley /*@
199320f4b53cSBarry Smith   DMSwarmRemovePointAtIndex - Removes a specific point from the `DMSWARM`
1994d3a51819SDave May 
199520f4b53cSBarry Smith   Not Collective
1996d3a51819SDave May 
199760225df5SJacob Faibussowitsch   Input Parameters:
199820f4b53cSBarry Smith + dm  - a `DMSWARM`
199962741f57SDave May - idx - index of point to remove
2000d3a51819SDave May 
2001d3a51819SDave May   Level: beginner
2002d3a51819SDave May 
200320f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRemovePoint()`
2004d3a51819SDave May @*/
2005d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmRemovePointAtIndex(DM dm, PetscInt idx)
2006d71ae5a4SJacob Faibussowitsch {
2007cb1d1399SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
2008cb1d1399SDave May 
2009521f74f9SMatthew G. Knepley   PetscFunctionBegin;
20109566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMSWARM_RemovePoints, 0, 0, 0, 0));
20119566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketRemovePointAtIndex(swarm->db, idx));
20129566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMSWARM_RemovePoints, 0, 0, 0, 0));
20133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2014cb1d1399SDave May }
2015b62e03f8SDave May 
201674d0cae8SMatthew G. Knepley /*@
201720f4b53cSBarry Smith   DMSwarmCopyPoint - Copy point pj to point pi in the `DMSWARM`
2018ba4fc9c6SDave May 
201920f4b53cSBarry Smith   Not Collective
2020ba4fc9c6SDave May 
202160225df5SJacob Faibussowitsch   Input Parameters:
202220f4b53cSBarry Smith + dm - a `DMSWARM`
2023ba4fc9c6SDave May . pi - the index of the point to copy
2024ba4fc9c6SDave May - pj - the point index where the copy should be located
2025ba4fc9c6SDave May 
2026ba4fc9c6SDave May   Level: beginner
2027ba4fc9c6SDave May 
202820f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRemovePoint()`
2029ba4fc9c6SDave May @*/
2030d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmCopyPoint(DM dm, PetscInt pi, PetscInt pj)
2031d71ae5a4SJacob Faibussowitsch {
2032ba4fc9c6SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
2033ba4fc9c6SDave May 
2034ba4fc9c6SDave May   PetscFunctionBegin;
20359566063dSJacob Faibussowitsch   if (!swarm->issetup) PetscCall(DMSetUp(dm));
20369566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketCopyPoint(swarm->db, pi, swarm->db, pj));
20373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2038ba4fc9c6SDave May }
2039ba4fc9c6SDave May 
204066976f2fSJacob Faibussowitsch static PetscErrorCode DMSwarmMigrate_Basic(DM dm, PetscBool remove_sent_points)
2041d71ae5a4SJacob Faibussowitsch {
2042521f74f9SMatthew G. Knepley   PetscFunctionBegin;
20439566063dSJacob Faibussowitsch   PetscCall(DMSwarmMigrate_Push_Basic(dm, remove_sent_points));
20443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
20453454631fSDave May }
20463454631fSDave May 
204774d0cae8SMatthew G. Knepley /*@
204820f4b53cSBarry Smith   DMSwarmMigrate - Relocates points defined in the `DMSWARM` to other MPI-ranks
2049d3a51819SDave May 
205020f4b53cSBarry Smith   Collective
2051d3a51819SDave May 
205260225df5SJacob Faibussowitsch   Input Parameters:
205320f4b53cSBarry Smith + dm                 - the `DMSWARM`
205462741f57SDave May - remove_sent_points - flag indicating if sent points should be removed from the current MPI-rank
2055d3a51819SDave May 
2056d3a51819SDave May   Level: advanced
2057d3a51819SDave May 
205820f4b53cSBarry Smith   Notes:
205920f4b53cSBarry Smith   The `DM` will be modified to accommodate received points.
206020f4b53cSBarry Smith   If `remove_sent_points` is `PETSC_TRUE`, any points that were sent will be removed from the `DM`.
206120f4b53cSBarry Smith   Different styles of migration are supported. See `DMSwarmSetMigrateType()`.
206220f4b53cSBarry Smith 
206320f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmSetMigrateType()`
2064d3a51819SDave May @*/
2065d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmMigrate(DM dm, PetscBool remove_sent_points)
2066d71ae5a4SJacob Faibussowitsch {
2067f0cdbbbaSDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
2068f0cdbbbaSDave May 
2069521f74f9SMatthew G. Knepley   PetscFunctionBegin;
20709566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMSWARM_Migrate, 0, 0, 0, 0));
2071f0cdbbbaSDave May   switch (swarm->migrate_type) {
2072d71ae5a4SJacob Faibussowitsch   case DMSWARM_MIGRATE_BASIC:
2073d71ae5a4SJacob Faibussowitsch     PetscCall(DMSwarmMigrate_Basic(dm, remove_sent_points));
2074d71ae5a4SJacob Faibussowitsch     break;
2075d71ae5a4SJacob Faibussowitsch   case DMSWARM_MIGRATE_DMCELLNSCATTER:
2076d71ae5a4SJacob Faibussowitsch     PetscCall(DMSwarmMigrate_CellDMScatter(dm, remove_sent_points));
2077d71ae5a4SJacob Faibussowitsch     break;
2078d71ae5a4SJacob Faibussowitsch   case DMSWARM_MIGRATE_DMCELLEXACT:
2079d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_MIGRATE_DMCELLEXACT not implemented");
2080d71ae5a4SJacob Faibussowitsch   case DMSWARM_MIGRATE_USER:
2081d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_MIGRATE_USER not implemented");
2082d71ae5a4SJacob Faibussowitsch   default:
2083d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_MIGRATE type unknown");
2084f0cdbbbaSDave May   }
20859566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMSWARM_Migrate, 0, 0, 0, 0));
20869566063dSJacob Faibussowitsch   PetscCall(DMClearGlobalVectors(dm));
20873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
20883454631fSDave May }
20893454631fSDave May 
2090f0cdbbbaSDave May PetscErrorCode DMSwarmMigrate_GlobalToLocal_Basic(DM dm, PetscInt *globalsize);
2091f0cdbbbaSDave May 
2092d3a51819SDave May /*
2093d3a51819SDave May  DMSwarmCollectViewCreate
2094d3a51819SDave May 
2095d3a51819SDave May  * Applies a collection method and gathers point neighbour points into dm
2096d3a51819SDave May 
2097d3a51819SDave May  Notes:
20988b8a3813SDave May  Users should call DMSwarmCollectViewDestroy() after
2099d3a51819SDave May  they have finished computations associated with the collected points
2100d3a51819SDave May */
2101d3a51819SDave May 
210274d0cae8SMatthew G. Knepley /*@
2103d3a51819SDave May   DMSwarmCollectViewCreate - Applies a collection method and gathers points
210420f4b53cSBarry Smith   in neighbour ranks into the `DMSWARM`
2105d3a51819SDave May 
210620f4b53cSBarry Smith   Collective
2107d3a51819SDave May 
210860225df5SJacob Faibussowitsch   Input Parameter:
210920f4b53cSBarry Smith . dm - the `DMSWARM`
2110d3a51819SDave May 
2111d3a51819SDave May   Level: advanced
2112d3a51819SDave May 
211320f4b53cSBarry Smith   Notes:
211420f4b53cSBarry Smith   Users should call `DMSwarmCollectViewDestroy()` after
211520f4b53cSBarry Smith   they have finished computations associated with the collected points
21160764c050SBarry Smith 
211720f4b53cSBarry Smith   Different collect methods are supported. See `DMSwarmSetCollectType()`.
211820f4b53cSBarry Smith 
21190764c050SBarry Smith   Developer Note:
21200764c050SBarry Smith   Create and Destroy routines create new objects that can get destroyed, they do not change the state
21210764c050SBarry Smith   of the current object.
21220764c050SBarry Smith 
212320f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmCollectViewDestroy()`, `DMSwarmSetCollectType()`
2124d3a51819SDave May @*/
2125d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmCollectViewCreate(DM dm)
2126d71ae5a4SJacob Faibussowitsch {
21272712d1f2SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
21282712d1f2SDave May   PetscInt  ng;
21292712d1f2SDave May 
2130521f74f9SMatthew G. Knepley   PetscFunctionBegin;
213128b400f6SJacob Faibussowitsch   PetscCheck(!swarm->collect_view_active, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "CollectView currently active");
21329566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetLocalSize(dm, &ng));
2133480eef7bSDave May   switch (swarm->collect_type) {
2134d71ae5a4SJacob Faibussowitsch   case DMSWARM_COLLECT_BASIC:
2135d71ae5a4SJacob Faibussowitsch     PetscCall(DMSwarmMigrate_GlobalToLocal_Basic(dm, &ng));
2136d71ae5a4SJacob Faibussowitsch     break;
2137d71ae5a4SJacob Faibussowitsch   case DMSWARM_COLLECT_DMDABOUNDINGBOX:
2138d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_COLLECT_DMDABOUNDINGBOX not implemented");
2139d71ae5a4SJacob Faibussowitsch   case DMSWARM_COLLECT_GENERAL:
2140d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_COLLECT_GENERAL not implemented");
2141d71ae5a4SJacob Faibussowitsch   default:
2142d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_COLLECT type unknown");
2143480eef7bSDave May   }
2144480eef7bSDave May   swarm->collect_view_active       = PETSC_TRUE;
2145480eef7bSDave May   swarm->collect_view_reset_nlocal = ng;
21463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
21472712d1f2SDave May }
21482712d1f2SDave May 
214974d0cae8SMatthew G. Knepley /*@
215020f4b53cSBarry Smith   DMSwarmCollectViewDestroy - Resets the `DMSWARM` to the size prior to calling `DMSwarmCollectViewCreate()`
2151d3a51819SDave May 
215220f4b53cSBarry Smith   Collective
2153d3a51819SDave May 
215460225df5SJacob Faibussowitsch   Input Parameters:
215520f4b53cSBarry Smith . dm - the `DMSWARM`
2156d3a51819SDave May 
2157d3a51819SDave May   Notes:
215820f4b53cSBarry Smith   Users should call `DMSwarmCollectViewCreate()` before this function is called.
2159d3a51819SDave May 
2160d3a51819SDave May   Level: advanced
2161d3a51819SDave May 
21620764c050SBarry Smith   Developer Note:
21630764c050SBarry Smith   Create and Destroy routines create new objects that can get destroyed, they do not change the state
21640764c050SBarry Smith   of the current object.
21650764c050SBarry Smith 
216620f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmCollectViewCreate()`, `DMSwarmSetCollectType()`
2167d3a51819SDave May @*/
2168d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmCollectViewDestroy(DM dm)
2169d71ae5a4SJacob Faibussowitsch {
21702712d1f2SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
21712712d1f2SDave May 
2172521f74f9SMatthew G. Knepley   PetscFunctionBegin;
217328b400f6SJacob Faibussowitsch   PetscCheck(swarm->collect_view_active, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "CollectView is currently not active");
21749566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(dm, swarm->collect_view_reset_nlocal, -1));
2175480eef7bSDave May   swarm->collect_view_active = PETSC_FALSE;
21763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
21772712d1f2SDave May }
21783454631fSDave May 
217966976f2fSJacob Faibussowitsch static PetscErrorCode DMSwarmSetUpPIC(DM dm)
2180d71ae5a4SJacob Faibussowitsch {
2181f0cdbbbaSDave May   PetscInt dim;
2182f0cdbbbaSDave May 
2183521f74f9SMatthew G. Knepley   PetscFunctionBegin;
21849566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetNumSpecies(dm, 1));
21859566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
218663a3b9bcSJacob Faibussowitsch   PetscCheck(dim >= 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Dimension must be 1,2,3 - found %" PetscInt_FMT, dim);
218763a3b9bcSJacob Faibussowitsch   PetscCheck(dim <= 3, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Dimension must be 1,2,3 - found %" PetscInt_FMT, dim);
21883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2189f0cdbbbaSDave May }
2190f0cdbbbaSDave May 
219174d0cae8SMatthew G. Knepley /*@
219231403186SMatthew G. Knepley   DMSwarmSetPointCoordinatesRandom - Sets initial coordinates for particles in each cell
219331403186SMatthew G. Knepley 
219420f4b53cSBarry Smith   Collective
219531403186SMatthew G. Knepley 
219660225df5SJacob Faibussowitsch   Input Parameters:
219720f4b53cSBarry Smith + dm  - the `DMSWARM`
219820f4b53cSBarry Smith - Npc - The number of particles per cell in the cell `DM`
219931403186SMatthew G. Knepley 
220031403186SMatthew G. Knepley   Level: intermediate
220131403186SMatthew G. Knepley 
220220f4b53cSBarry Smith   Notes:
220320f4b53cSBarry Smith   The user must use `DMSwarmSetCellDM()` to set the cell `DM` first. The particles are placed randomly inside each cell. If only
220420f4b53cSBarry Smith   one particle is in each cell, it is placed at the centroid.
220520f4b53cSBarry Smith 
220620f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmSetCellDM()`
220731403186SMatthew G. Knepley @*/
2208d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmSetPointCoordinatesRandom(DM dm, PetscInt Npc)
2209d71ae5a4SJacob Faibussowitsch {
221031403186SMatthew G. Knepley   DM             cdm;
221119307e5cSMatthew G. Knepley   DMSwarmCellDM  celldm;
221231403186SMatthew G. Knepley   PetscRandom    rnd;
221331403186SMatthew G. Knepley   DMPolytopeType ct;
221431403186SMatthew G. Knepley   PetscBool      simplex;
221531403186SMatthew G. Knepley   PetscReal     *centroid, *coords, *xi0, *v0, *J, *invJ, detJ;
221619307e5cSMatthew G. Knepley   PetscInt       dim, d, cStart, cEnd, c, p, Nfc;
221719307e5cSMatthew G. Knepley   const char   **coordFields;
221831403186SMatthew G. Knepley 
221931403186SMatthew G. Knepley   PetscFunctionBeginUser;
22209566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)dm), &rnd));
22219566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetInterval(rnd, -1.0, 1.0));
22229566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetType(rnd, PETSCRAND48));
222331403186SMatthew G. Knepley 
222419307e5cSMatthew G. Knepley   PetscCall(DMSwarmGetCellDMActive(dm, &celldm));
222519307e5cSMatthew G. Knepley   PetscCall(DMSwarmCellDMGetCoordinateFields(celldm, &Nfc, &coordFields));
222619307e5cSMatthew G. Knepley   PetscCheck(Nfc == 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "We only support a single coordinate field right now, not %" PetscInt_FMT, Nfc);
22279566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetCellDM(dm, &cdm));
22289566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(cdm, &dim));
22299566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(cdm, 0, &cStart, &cEnd));
22309566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCellType(cdm, cStart, &ct));
223131403186SMatthew G. Knepley   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE;
223231403186SMatthew G. Knepley 
22339566063dSJacob Faibussowitsch   PetscCall(PetscMalloc5(dim, &centroid, dim, &xi0, dim, &v0, dim * dim, &J, dim * dim, &invJ));
223431403186SMatthew G. Knepley   for (d = 0; d < dim; ++d) xi0[d] = -1.0;
223519307e5cSMatthew G. Knepley   PetscCall(DMSwarmGetField(dm, coordFields[0], NULL, NULL, (void **)&coords));
223631403186SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
223731403186SMatthew G. Knepley     if (Npc == 1) {
22389566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFVM(cdm, c, NULL, centroid, NULL));
223931403186SMatthew G. Knepley       for (d = 0; d < dim; ++d) coords[c * dim + d] = centroid[d];
224031403186SMatthew G. Knepley     } else {
22419566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFEM(cdm, c, NULL, v0, J, invJ, &detJ)); /* affine */
224231403186SMatthew G. Knepley       for (p = 0; p < Npc; ++p) {
224331403186SMatthew G. Knepley         const PetscInt n   = c * Npc + p;
224431403186SMatthew G. Knepley         PetscReal      sum = 0.0, refcoords[3];
224531403186SMatthew G. Knepley 
224631403186SMatthew G. Knepley         for (d = 0; d < dim; ++d) {
22479566063dSJacob Faibussowitsch           PetscCall(PetscRandomGetValueReal(rnd, &refcoords[d]));
224831403186SMatthew G. Knepley           sum += refcoords[d];
224931403186SMatthew G. Knepley         }
22509371c9d4SSatish Balay         if (simplex && sum > 0.0)
22519371c9d4SSatish Balay           for (d = 0; d < dim; ++d) refcoords[d] -= PetscSqrtReal(dim) * sum;
225231403186SMatthew G. Knepley         CoordinatesRefToReal(dim, dim, xi0, v0, J, refcoords, &coords[n * dim]);
225331403186SMatthew G. Knepley       }
225431403186SMatthew G. Knepley     }
225531403186SMatthew G. Knepley   }
225619307e5cSMatthew G. Knepley   PetscCall(DMSwarmRestoreField(dm, coordFields[0], NULL, NULL, (void **)&coords));
22579566063dSJacob Faibussowitsch   PetscCall(PetscFree5(centroid, xi0, v0, J, invJ));
22589566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rnd));
22593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
226031403186SMatthew G. Knepley }
226131403186SMatthew G. Knepley 
226231403186SMatthew G. Knepley /*@
2263fca3fa51SMatthew G. Knepley   DMSwarmGetType - Get particular flavor of `DMSWARM`
2264fca3fa51SMatthew G. Knepley 
2265fca3fa51SMatthew G. Knepley   Collective
2266fca3fa51SMatthew G. Knepley 
2267fca3fa51SMatthew G. Knepley   Input Parameter:
2268fca3fa51SMatthew G. Knepley . sw - the `DMSWARM`
2269fca3fa51SMatthew G. Knepley 
2270fca3fa51SMatthew G. Knepley   Output Parameter:
2271fca3fa51SMatthew G. Knepley . stype - the `DMSWARM` type (e.g. `DMSWARM_PIC`)
2272fca3fa51SMatthew G. Knepley 
2273fca3fa51SMatthew G. Knepley   Level: advanced
2274fca3fa51SMatthew G. Knepley 
2275fca3fa51SMatthew G. Knepley .seealso: `DM`, `DMSWARM`, `DMSwarmSetMigrateType()`, `DMSwarmSetCollectType()`, `DMSwarmType`, `DMSWARM_PIC`, `DMSWARM_BASIC`
2276fca3fa51SMatthew G. Knepley @*/
2277fca3fa51SMatthew G. Knepley PetscErrorCode DMSwarmGetType(DM sw, DMSwarmType *stype)
2278fca3fa51SMatthew G. Knepley {
2279fca3fa51SMatthew G. Knepley   DM_Swarm *swarm = (DM_Swarm *)sw->data;
2280fca3fa51SMatthew G. Knepley 
2281fca3fa51SMatthew G. Knepley   PetscFunctionBegin;
2282fca3fa51SMatthew G. Knepley   PetscValidHeaderSpecific(sw, DM_CLASSID, 1);
2283fca3fa51SMatthew G. Knepley   PetscAssertPointer(stype, 2);
2284fca3fa51SMatthew G. Knepley   *stype = swarm->swarm_type;
2285fca3fa51SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
2286fca3fa51SMatthew G. Knepley }
2287fca3fa51SMatthew G. Knepley 
2288fca3fa51SMatthew G. Knepley /*@
228920f4b53cSBarry Smith   DMSwarmSetType - Set particular flavor of `DMSWARM`
2290d3a51819SDave May 
229120f4b53cSBarry Smith   Collective
2292d3a51819SDave May 
229360225df5SJacob Faibussowitsch   Input Parameters:
2294fca3fa51SMatthew G. Knepley + sw    - the `DMSWARM`
229520f4b53cSBarry Smith - stype - the `DMSWARM` type (e.g. `DMSWARM_PIC`)
2296d3a51819SDave May 
2297d3a51819SDave May   Level: advanced
2298d3a51819SDave May 
229920f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmSetMigrateType()`, `DMSwarmSetCollectType()`, `DMSwarmType`, `DMSWARM_PIC`, `DMSWARM_BASIC`
2300d3a51819SDave May @*/
2301fca3fa51SMatthew G. Knepley PetscErrorCode DMSwarmSetType(DM sw, DMSwarmType stype)
2302d71ae5a4SJacob Faibussowitsch {
2303fca3fa51SMatthew G. Knepley   DM_Swarm *swarm = (DM_Swarm *)sw->data;
2304f0cdbbbaSDave May 
2305521f74f9SMatthew G. Knepley   PetscFunctionBegin;
2306fca3fa51SMatthew G. Knepley   PetscValidHeaderSpecific(sw, DM_CLASSID, 1);
2307f0cdbbbaSDave May   swarm->swarm_type = stype;
2308fca3fa51SMatthew G. Knepley   if (swarm->swarm_type == DMSWARM_PIC) PetscCall(DMSwarmSetUpPIC(sw));
23093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2310f0cdbbbaSDave May }
2311f0cdbbbaSDave May 
2312fca3fa51SMatthew G. Knepley static PetscErrorCode DMSwarmCreateRemapDM_Private(DM sw, DM *rdm)
2313d71ae5a4SJacob Faibussowitsch {
2314fca3fa51SMatthew G. Knepley   PetscFE        fe;
2315fca3fa51SMatthew G. Knepley   DMPolytopeType ct;
2316fca3fa51SMatthew G. Knepley   PetscInt       dim, cStart;
2317fca3fa51SMatthew G. Knepley   const char    *prefix = "remap_";
2318fca3fa51SMatthew G. Knepley 
2319fca3fa51SMatthew G. Knepley   PetscFunctionBegin;
2320fca3fa51SMatthew G. Knepley   PetscCall(DMCreate(PetscObjectComm((PetscObject)sw), rdm));
2321fca3fa51SMatthew G. Knepley   PetscCall(DMSetType(*rdm, DMPLEX));
2322fca3fa51SMatthew G. Knepley   PetscCall(DMPlexSetOptionsPrefix(*rdm, prefix));
2323fca3fa51SMatthew G. Knepley   PetscCall(DMSetFromOptions(*rdm));
2324fca3fa51SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)*rdm, "remap"));
2325fca3fa51SMatthew G. Knepley   PetscCall(DMViewFromOptions(*rdm, NULL, "-dm_view"));
2326fca3fa51SMatthew G. Knepley 
2327fca3fa51SMatthew G. Knepley   PetscCall(DMGetDimension(*rdm, &dim));
2328fca3fa51SMatthew G. Knepley   PetscCall(DMPlexGetHeightStratum(*rdm, 0, &cStart, NULL));
2329fca3fa51SMatthew G. Knepley   PetscCall(DMPlexGetCellType(*rdm, cStart, &ct));
2330fca3fa51SMatthew G. Knepley   PetscCall(PetscFECreateByCell(PETSC_COMM_SELF, dim, 1, ct, prefix, PETSC_DETERMINE, &fe));
2331fca3fa51SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)fe, "distribution"));
2332fca3fa51SMatthew G. Knepley   PetscCall(DMSetField(*rdm, 0, NULL, (PetscObject)fe));
2333fca3fa51SMatthew G. Knepley   PetscCall(DMCreateDS(*rdm));
2334fca3fa51SMatthew G. Knepley   PetscCall(PetscFEDestroy(&fe));
2335fca3fa51SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
2336fca3fa51SMatthew G. Knepley }
2337fca3fa51SMatthew G. Knepley 
2338fca3fa51SMatthew G. Knepley static PetscErrorCode DMSetup_Swarm(DM sw)
2339fca3fa51SMatthew G. Knepley {
2340fca3fa51SMatthew G. Knepley   DM_Swarm *swarm = (DM_Swarm *)sw->data;
23413454631fSDave May 
2342521f74f9SMatthew G. Knepley   PetscFunctionBegin;
23433ba16761SJacob Faibussowitsch   if (swarm->issetup) PetscFunctionReturn(PETSC_SUCCESS);
23443454631fSDave May   swarm->issetup = PETSC_TRUE;
23453454631fSDave May 
2346fca3fa51SMatthew G. Knepley   if (swarm->remap_type != DMSWARM_REMAP_NONE) {
2347fca3fa51SMatthew G. Knepley     DMSwarmCellDM celldm;
2348fca3fa51SMatthew G. Knepley     DM            rdm;
2349fca3fa51SMatthew G. Knepley     const char   *fieldnames[2]  = {DMSwarmPICField_coor, "velocity"};
2350fca3fa51SMatthew G. Knepley     const char   *vfieldnames[1] = {"w_q"};
2351fca3fa51SMatthew G. Knepley 
2352fca3fa51SMatthew G. Knepley     PetscCall(DMSwarmCreateRemapDM_Private(sw, &rdm));
2353fca3fa51SMatthew G. Knepley     PetscCall(DMSwarmCellDMCreate(rdm, 1, vfieldnames, 2, fieldnames, &celldm));
2354fca3fa51SMatthew G. Knepley     PetscCall(DMSwarmAddCellDM(sw, celldm));
2355fca3fa51SMatthew G. Knepley     PetscCall(DMSwarmCellDMDestroy(&celldm));
2356fca3fa51SMatthew G. Knepley     PetscCall(DMDestroy(&rdm));
2357fca3fa51SMatthew G. Knepley   }
2358fca3fa51SMatthew G. Knepley 
2359f0cdbbbaSDave May   if (swarm->swarm_type == DMSWARM_PIC) {
236019307e5cSMatthew G. Knepley     DMSwarmCellDM celldm;
2361f0cdbbbaSDave May 
2362fca3fa51SMatthew G. Knepley     PetscCall(DMSwarmGetCellDMActive(sw, &celldm));
2363fca3fa51SMatthew G. Knepley     PetscCheck(celldm, PetscObjectComm((PetscObject)sw), PETSC_ERR_USER, "No active cell DM. DMSWARM_PIC requires you call DMSwarmSetCellDM() or DMSwarmAddCellDM()");
236419307e5cSMatthew G. Knepley     if (celldm->dm->ops->locatepointssubdomain) {
2365f0cdbbbaSDave May       /* check methods exists for exact ownership identificiation */
2366fca3fa51SMatthew G. Knepley       PetscCall(PetscInfo(sw, "DMSWARM_PIC: Using method CellDM->ops->LocatePointsSubdomain\n"));
2367f0cdbbbaSDave May       swarm->migrate_type = DMSWARM_MIGRATE_DMCELLEXACT;
2368f0cdbbbaSDave May     } else {
2369f0cdbbbaSDave May       /* check methods exist for point location AND rank neighbor identification */
237019307e5cSMatthew G. Knepley       if (celldm->dm->ops->locatepoints) {
2371fca3fa51SMatthew G. Knepley         PetscCall(PetscInfo(sw, "DMSWARM_PIC: Using method CellDM->LocatePoints\n"));
2372fca3fa51SMatthew G. Knepley       } else SETERRQ(PetscObjectComm((PetscObject)sw), PETSC_ERR_USER, "DMSWARM_PIC requires the method CellDM->ops->locatepoints be defined");
2373f0cdbbbaSDave May 
237419307e5cSMatthew G. Knepley       if (celldm->dm->ops->getneighbors) {
2375fca3fa51SMatthew G. Knepley         PetscCall(PetscInfo(sw, "DMSWARM_PIC: Using method CellDM->GetNeigbors\n"));
2376fca3fa51SMatthew G. Knepley       } else SETERRQ(PetscObjectComm((PetscObject)sw), PETSC_ERR_USER, "DMSWARM_PIC requires the method CellDM->ops->getneighbors be defined");
2377f0cdbbbaSDave May 
2378f0cdbbbaSDave May       swarm->migrate_type = DMSWARM_MIGRATE_DMCELLNSCATTER;
2379f0cdbbbaSDave May     }
2380f0cdbbbaSDave May   }
2381f0cdbbbaSDave May 
2382fca3fa51SMatthew G. Knepley   PetscCall(DMSwarmFinalizeFieldRegister(sw));
2383f0cdbbbaSDave May 
23843454631fSDave May   /* check some fields were registered */
2385fca3fa51SMatthew G. Knepley   PetscCheck(swarm->db->nfields > 2, PetscObjectComm((PetscObject)sw), PETSC_ERR_USER, "At least one field user must be registered via DMSwarmRegisterXXX()");
23863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
23873454631fSDave May }
23883454631fSDave May 
238966976f2fSJacob Faibussowitsch static PetscErrorCode DMDestroy_Swarm(DM dm)
2390d71ae5a4SJacob Faibussowitsch {
239157795646SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
239257795646SDave May 
239357795646SDave May   PetscFunctionBegin;
23943ba16761SJacob Faibussowitsch   if (--swarm->refct > 0) PetscFunctionReturn(PETSC_SUCCESS);
239519307e5cSMatthew G. Knepley   PetscCall(PetscObjectListDestroy(&swarm->cellDMs));
239619307e5cSMatthew G. Knepley   PetscCall(PetscFree(swarm->activeCellDM));
23979566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketDestroy(&swarm->db));
23989566063dSJacob Faibussowitsch   PetscCall(PetscFree(swarm));
23993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
240057795646SDave May }
240157795646SDave May 
240266976f2fSJacob Faibussowitsch static PetscErrorCode DMSwarmView_Draw(DM dm, PetscViewer viewer)
2403d71ae5a4SJacob Faibussowitsch {
2404a9ee3421SMatthew G. Knepley   DM            cdm;
240519307e5cSMatthew G. Knepley   DMSwarmCellDM celldm;
2406a9ee3421SMatthew G. Knepley   PetscDraw     draw;
240731403186SMatthew G. Knepley   PetscReal    *coords, oldPause, radius = 0.01;
240819307e5cSMatthew G. Knepley   PetscInt      Np, p, bs, Nfc;
240919307e5cSMatthew G. Knepley   const char  **coordFields;
2410a9ee3421SMatthew G. Knepley 
2411a9ee3421SMatthew G. Knepley   PetscFunctionBegin;
24129566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, ((PetscObject)dm)->prefix, "-dm_view_swarm_radius", &radius, NULL));
24139566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
24149566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetCellDM(dm, &cdm));
24159566063dSJacob Faibussowitsch   PetscCall(PetscDrawGetPause(draw, &oldPause));
24169566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetPause(draw, 0.0));
24179566063dSJacob Faibussowitsch   PetscCall(DMView(cdm, viewer));
24189566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetPause(draw, oldPause));
2419a9ee3421SMatthew G. Knepley 
242019307e5cSMatthew G. Knepley   PetscCall(DMSwarmGetCellDMActive(dm, &celldm));
242119307e5cSMatthew G. Knepley   PetscCall(DMSwarmCellDMGetCoordinateFields(celldm, &Nfc, &coordFields));
242219307e5cSMatthew G. Knepley   PetscCheck(Nfc == 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "We only support a single coordinate field right now, not %" PetscInt_FMT, Nfc);
24239566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetLocalSize(dm, &Np));
242419307e5cSMatthew G. Knepley   PetscCall(DMSwarmGetField(dm, coordFields[0], &bs, NULL, (void **)&coords));
2425a9ee3421SMatthew G. Knepley   for (p = 0; p < Np; ++p) {
2426a9ee3421SMatthew G. Knepley     const PetscInt i = p * bs;
2427a9ee3421SMatthew G. Knepley 
24289566063dSJacob Faibussowitsch     PetscCall(PetscDrawEllipse(draw, coords[i], coords[i + 1], radius, radius, PETSC_DRAW_BLUE));
2429a9ee3421SMatthew G. Knepley   }
243019307e5cSMatthew G. Knepley   PetscCall(DMSwarmRestoreField(dm, coordFields[0], &bs, NULL, (void **)&coords));
24319566063dSJacob Faibussowitsch   PetscCall(PetscDrawFlush(draw));
24329566063dSJacob Faibussowitsch   PetscCall(PetscDrawPause(draw));
24339566063dSJacob Faibussowitsch   PetscCall(PetscDrawSave(draw));
24343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2435a9ee3421SMatthew G. Knepley }
2436a9ee3421SMatthew G. Knepley 
2437d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMView_Swarm_Ascii(DM dm, PetscViewer viewer)
2438d71ae5a4SJacob Faibussowitsch {
24396a5217c0SMatthew G. Knepley   PetscViewerFormat format;
24406a5217c0SMatthew G. Knepley   PetscInt         *sizes;
24416a5217c0SMatthew G. Knepley   PetscInt          dim, Np, maxSize = 17;
24426a5217c0SMatthew G. Knepley   MPI_Comm          comm;
24436a5217c0SMatthew G. Knepley   PetscMPIInt       rank, size, p;
244419307e5cSMatthew G. Knepley   const char       *name, *cellid;
24456a5217c0SMatthew G. Knepley 
24466a5217c0SMatthew G. Knepley   PetscFunctionBegin;
24476a5217c0SMatthew G. Knepley   PetscCall(PetscViewerGetFormat(viewer, &format));
24486a5217c0SMatthew G. Knepley   PetscCall(DMGetDimension(dm, &dim));
24496a5217c0SMatthew G. Knepley   PetscCall(DMSwarmGetLocalSize(dm, &Np));
24506a5217c0SMatthew G. Knepley   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
24516a5217c0SMatthew G. Knepley   PetscCallMPI(MPI_Comm_rank(comm, &rank));
24526a5217c0SMatthew G. Knepley   PetscCallMPI(MPI_Comm_size(comm, &size));
24536a5217c0SMatthew G. Knepley   PetscCall(PetscObjectGetName((PetscObject)dm, &name));
245463a3b9bcSJacob Faibussowitsch   if (name) PetscCall(PetscViewerASCIIPrintf(viewer, "%s in %" PetscInt_FMT " dimension%s:\n", name, dim, dim == 1 ? "" : "s"));
245563a3b9bcSJacob Faibussowitsch   else PetscCall(PetscViewerASCIIPrintf(viewer, "Swarm in %" PetscInt_FMT " dimension%s:\n", dim, dim == 1 ? "" : "s"));
24566a5217c0SMatthew G. Knepley   if (size < maxSize) PetscCall(PetscCalloc1(size, &sizes));
24576a5217c0SMatthew G. Knepley   else PetscCall(PetscCalloc1(3, &sizes));
24586a5217c0SMatthew G. Knepley   if (size < maxSize) {
24596a5217c0SMatthew G. Knepley     PetscCallMPI(MPI_Gather(&Np, 1, MPIU_INT, sizes, 1, MPIU_INT, 0, comm));
24606a5217c0SMatthew G. Knepley     PetscCall(PetscViewerASCIIPrintf(viewer, "  Number of particles per rank:"));
24616a5217c0SMatthew G. Knepley     for (p = 0; p < size; ++p) {
24626a5217c0SMatthew G. Knepley       if (rank == 0) PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT, sizes[p]));
24636a5217c0SMatthew G. Knepley     }
24646a5217c0SMatthew G. Knepley   } else {
24656a5217c0SMatthew G. Knepley     PetscInt locMinMax[2] = {Np, Np};
24666a5217c0SMatthew G. Knepley 
24676a5217c0SMatthew G. Knepley     PetscCall(PetscGlobalMinMaxInt(comm, locMinMax, sizes));
24686a5217c0SMatthew G. Knepley     PetscCall(PetscViewerASCIIPrintf(viewer, "  Min/Max of particles per rank: %" PetscInt_FMT "/%" PetscInt_FMT, sizes[0], sizes[1]));
24696a5217c0SMatthew G. Knepley   }
24706a5217c0SMatthew G. Knepley   PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
24716a5217c0SMatthew G. Knepley   PetscCall(PetscFree(sizes));
24726a5217c0SMatthew G. Knepley   if (format == PETSC_VIEWER_ASCII_INFO) {
247319307e5cSMatthew G. Knepley     DMSwarmCellDM celldm;
24746a5217c0SMatthew G. Knepley     PetscInt     *cell;
24756a5217c0SMatthew G. Knepley 
24766a5217c0SMatthew G. Knepley     PetscCall(PetscViewerASCIIPrintf(viewer, "  Cells containing each particle:\n"));
24776a5217c0SMatthew G. Knepley     PetscCall(PetscViewerASCIIPushSynchronized(viewer));
247819307e5cSMatthew G. Knepley     PetscCall(DMSwarmGetCellDMActive(dm, &celldm));
247919307e5cSMatthew G. Knepley     PetscCall(DMSwarmCellDMGetCellID(celldm, &cellid));
248019307e5cSMatthew G. Knepley     PetscCall(DMSwarmGetField(dm, cellid, NULL, NULL, (void **)&cell));
248163a3b9bcSJacob Faibussowitsch     for (p = 0; p < Np; ++p) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "  p%d: %" PetscInt_FMT "\n", p, cell[p]));
24826a5217c0SMatthew G. Knepley     PetscCall(PetscViewerFlush(viewer));
24836a5217c0SMatthew G. Knepley     PetscCall(PetscViewerASCIIPopSynchronized(viewer));
248419307e5cSMatthew G. Knepley     PetscCall(DMSwarmRestoreField(dm, cellid, NULL, NULL, (void **)&cell));
24856a5217c0SMatthew G. Knepley   }
24863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
24876a5217c0SMatthew G. Knepley }
24886a5217c0SMatthew G. Knepley 
248966976f2fSJacob Faibussowitsch static PetscErrorCode DMView_Swarm(DM dm, PetscViewer viewer)
2490d71ae5a4SJacob Faibussowitsch {
24915f50eb2eSDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
2492*9f196a02SMartin Diehl   PetscBool isascii, ibinary, isvtk, isdraw, ispython;
24935627991aSBarry Smith #if defined(PETSC_HAVE_HDF5)
24945627991aSBarry Smith   PetscBool ishdf5;
24955627991aSBarry Smith #endif
24965f50eb2eSDave May 
24975f50eb2eSDave May   PetscFunctionBegin;
24985f50eb2eSDave May   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
24995f50eb2eSDave May   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
2500*9f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
25019566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &ibinary));
25029566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERVTK, &isvtk));
25035627991aSBarry Smith #if defined(PETSC_HAVE_HDF5)
25049566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
25055627991aSBarry Smith #endif
25069566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
25074fc69c0aSMatthew G. Knepley   PetscCall(PetscObjectHasFunction((PetscObject)viewer, "PetscViewerPythonViewObject_C", &ispython));
2508*9f196a02SMartin Diehl   if (isascii) {
25096a5217c0SMatthew G. Knepley     PetscViewerFormat format;
25106a5217c0SMatthew G. Knepley 
25116a5217c0SMatthew G. Knepley     PetscCall(PetscViewerGetFormat(viewer, &format));
25126a5217c0SMatthew G. Knepley     switch (format) {
2513d71ae5a4SJacob Faibussowitsch     case PETSC_VIEWER_ASCII_INFO_DETAIL:
2514d71ae5a4SJacob Faibussowitsch       PetscCall(DMSwarmDataBucketView(PetscObjectComm((PetscObject)dm), swarm->db, NULL, DATABUCKET_VIEW_STDOUT));
2515d71ae5a4SJacob Faibussowitsch       break;
2516d71ae5a4SJacob Faibussowitsch     default:
2517d71ae5a4SJacob Faibussowitsch       PetscCall(DMView_Swarm_Ascii(dm, viewer));
25186a5217c0SMatthew G. Knepley     }
2519f7d195e4SLawrence Mitchell   } else {
25205f50eb2eSDave May #if defined(PETSC_HAVE_HDF5)
252148a46eb9SPierre Jolivet     if (ishdf5) PetscCall(DMSwarmView_HDF5(dm, viewer));
25225f50eb2eSDave May #endif
252348a46eb9SPierre Jolivet     if (isdraw) PetscCall(DMSwarmView_Draw(dm, viewer));
25244fc69c0aSMatthew G. Knepley     if (ispython) PetscCall(PetscViewerPythonViewObject(viewer, (PetscObject)dm));
2525f7d195e4SLawrence Mitchell   }
25263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
25275f50eb2eSDave May }
25285f50eb2eSDave May 
2529cc4c1da9SBarry Smith /*@
253020f4b53cSBarry Smith   DMSwarmGetCellSwarm - Extracts a single cell from the `DMSWARM` object, returns it as a single cell `DMSWARM`.
253120f4b53cSBarry 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.
2532d0c080abSJoseph Pusztay 
2533d0c080abSJoseph Pusztay   Noncollective
2534d0c080abSJoseph Pusztay 
253560225df5SJacob Faibussowitsch   Input Parameters:
253620f4b53cSBarry Smith + sw        - the `DMSWARM`
25375627991aSBarry Smith . cellID    - the integer id of the cell to be extracted and filtered
253820f4b53cSBarry Smith - cellswarm - The `DMSWARM` to receive the cell
2539d0c080abSJoseph Pusztay 
2540d0c080abSJoseph Pusztay   Level: beginner
2541d0c080abSJoseph Pusztay 
25425627991aSBarry Smith   Notes:
254320f4b53cSBarry Smith   This presently only supports `DMSWARM_PIC` type
25445627991aSBarry Smith 
254520f4b53cSBarry Smith   Should be restored with `DMSwarmRestoreCellSwarm()`
2546d0c080abSJoseph Pusztay 
254720f4b53cSBarry Smith   Changes to this cell of the swarm will be lost if they are made prior to restoring this cell.
254820f4b53cSBarry Smith 
254920f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRestoreCellSwarm()`
2550d0c080abSJoseph Pusztay @*/
2551cc4c1da9SBarry Smith PetscErrorCode DMSwarmGetCellSwarm(DM sw, PetscInt cellID, DM cellswarm)
2552d71ae5a4SJacob Faibussowitsch {
2553d0c080abSJoseph Pusztay   DM_Swarm   *original = (DM_Swarm *)sw->data;
2554d0c080abSJoseph Pusztay   DMLabel     label;
2555d0c080abSJoseph Pusztay   DM          dmc, subdmc;
2556d0c080abSJoseph Pusztay   PetscInt   *pids, particles, dim;
255719307e5cSMatthew G. Knepley   const char *name;
2558d0c080abSJoseph Pusztay 
2559d0c080abSJoseph Pusztay   PetscFunctionBegin;
2560d0c080abSJoseph Pusztay   /* Configure new swarm */
25619566063dSJacob Faibussowitsch   PetscCall(DMSetType(cellswarm, DMSWARM));
25629566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(sw, &dim));
25639566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(cellswarm, dim));
25649566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetType(cellswarm, DMSWARM_PIC));
2565d0c080abSJoseph Pusztay   /* Destroy the unused, unconfigured data bucket to prevent stragglers in memory */
25669566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketDestroy(&((DM_Swarm *)cellswarm->data)->db));
25679566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortGetAccess(sw));
25689566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortGetNumberOfPointsPerCell(sw, cellID, &particles));
25699566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortGetPointsPerCell(sw, cellID, &particles, &pids));
25709566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketCreateFromSubset(original->db, particles, pids, &((DM_Swarm *)cellswarm->data)->db));
25719566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortRestoreAccess(sw));
2572fe11354eSMatthew G. Knepley   PetscCall(DMSwarmSortRestorePointsPerCell(sw, cellID, &particles, &pids));
25739566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetCellDM(sw, &dmc));
25749566063dSJacob Faibussowitsch   PetscCall(DMLabelCreate(PetscObjectComm((PetscObject)sw), "singlecell", &label));
25759566063dSJacob Faibussowitsch   PetscCall(DMAddLabel(dmc, label));
25769566063dSJacob Faibussowitsch   PetscCall(DMLabelSetValue(label, cellID, 1));
257730cbcd5dSksagiyam   PetscCall(DMPlexFilter(dmc, label, 1, PETSC_FALSE, PETSC_FALSE, NULL, &subdmc));
257819307e5cSMatthew G. Knepley   PetscCall(PetscObjectGetName((PetscObject)dmc, &name));
257919307e5cSMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)subdmc, name));
25809566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetCellDM(cellswarm, subdmc));
25819566063dSJacob Faibussowitsch   PetscCall(DMLabelDestroy(&label));
25823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2583d0c080abSJoseph Pusztay }
2584d0c080abSJoseph Pusztay 
2585cc4c1da9SBarry Smith /*@
258620f4b53cSBarry Smith   DMSwarmRestoreCellSwarm - Restores a `DMSWARM` object obtained with `DMSwarmGetCellSwarm()`. All fields are copied back into the parent swarm.
2587d0c080abSJoseph Pusztay 
2588d0c080abSJoseph Pusztay   Noncollective
2589d0c080abSJoseph Pusztay 
259060225df5SJacob Faibussowitsch   Input Parameters:
259120f4b53cSBarry Smith + sw        - the parent `DMSWARM`
2592d0c080abSJoseph Pusztay . cellID    - the integer id of the cell to be copied back into the parent swarm
2593d0c080abSJoseph Pusztay - cellswarm - the cell swarm object
2594d0c080abSJoseph Pusztay 
2595d0c080abSJoseph Pusztay   Level: beginner
2596d0c080abSJoseph Pusztay 
25975627991aSBarry Smith   Note:
259820f4b53cSBarry Smith   This only supports `DMSWARM_PIC` types of `DMSWARM`s
2599d0c080abSJoseph Pusztay 
260020f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmGetCellSwarm()`
2601d0c080abSJoseph Pusztay @*/
2602cc4c1da9SBarry Smith PetscErrorCode DMSwarmRestoreCellSwarm(DM sw, PetscInt cellID, DM cellswarm)
2603d71ae5a4SJacob Faibussowitsch {
2604d0c080abSJoseph Pusztay   DM        dmc;
2605d0c080abSJoseph Pusztay   PetscInt *pids, particles, p;
2606d0c080abSJoseph Pusztay 
2607d0c080abSJoseph Pusztay   PetscFunctionBegin;
26089566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortGetAccess(sw));
26099566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortGetPointsPerCell(sw, cellID, &particles, &pids));
26109566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortRestoreAccess(sw));
2611d0c080abSJoseph Pusztay   /* Pointwise copy of each particle based on pid. The parent swarm may not be altered during this process. */
261248a46eb9SPierre Jolivet   for (p = 0; p < particles; ++p) PetscCall(DMSwarmDataBucketCopyPoint(((DM_Swarm *)cellswarm->data)->db, pids[p], ((DM_Swarm *)sw->data)->db, pids[p]));
2613d0c080abSJoseph Pusztay   /* Free memory, destroy cell dm */
26149566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetCellDM(cellswarm, &dmc));
26159566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmc));
2616fe11354eSMatthew G. Knepley   PetscCall(DMSwarmSortRestorePointsPerCell(sw, cellID, &particles, &pids));
26173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2618d0c080abSJoseph Pusztay }
2619d0c080abSJoseph Pusztay 
2620d52c2f21SMatthew G. Knepley /*@
2621d52c2f21SMatthew G. Knepley   DMSwarmComputeMoments - Compute the first three particle moments for a given field
2622d52c2f21SMatthew G. Knepley 
2623d52c2f21SMatthew G. Knepley   Noncollective
2624d52c2f21SMatthew G. Knepley 
2625d52c2f21SMatthew G. Knepley   Input Parameters:
2626d52c2f21SMatthew G. Knepley + sw         - the `DMSWARM`
2627d52c2f21SMatthew G. Knepley . coordinate - the coordinate field name
2628d52c2f21SMatthew G. Knepley - weight     - the weight field name
2629d52c2f21SMatthew G. Knepley 
2630d52c2f21SMatthew G. Knepley   Output Parameter:
2631d52c2f21SMatthew G. Knepley . moments - the field moments
2632d52c2f21SMatthew G. Knepley 
2633d52c2f21SMatthew G. Knepley   Level: intermediate
2634d52c2f21SMatthew G. Knepley 
2635d52c2f21SMatthew G. Knepley   Notes:
2636d52c2f21SMatthew G. Knepley   The `moments` array should be of length bs + 2, where bs is the block size of the coordinate field.
2637d52c2f21SMatthew G. Knepley 
2638d52c2f21SMatthew G. Knepley   The weight field must be a scalar, having blocksize 1.
2639d52c2f21SMatthew G. Knepley 
2640d52c2f21SMatthew G. Knepley .seealso: `DM`, `DMSWARM`, `DMPlexComputeMoments()`
2641d52c2f21SMatthew G. Knepley @*/
2642d52c2f21SMatthew G. Knepley PetscErrorCode DMSwarmComputeMoments(DM sw, const char coordinate[], const char weight[], PetscReal moments[])
2643d52c2f21SMatthew G. Knepley {
2644d52c2f21SMatthew G. Knepley   const PetscReal *coords;
2645d52c2f21SMatthew G. Knepley   const PetscReal *w;
2646d52c2f21SMatthew G. Knepley   PetscReal       *mom;
2647d52c2f21SMatthew G. Knepley   PetscDataType    dtc, dtw;
2648d52c2f21SMatthew G. Knepley   PetscInt         bsc, bsw, Np;
2649d52c2f21SMatthew G. Knepley   MPI_Comm         comm;
2650d52c2f21SMatthew G. Knepley 
2651d52c2f21SMatthew G. Knepley   PetscFunctionBegin;
2652d52c2f21SMatthew G. Knepley   PetscValidHeaderSpecific(sw, DM_CLASSID, 1);
2653d52c2f21SMatthew G. Knepley   PetscAssertPointer(coordinate, 2);
2654d52c2f21SMatthew G. Knepley   PetscAssertPointer(weight, 3);
2655d52c2f21SMatthew G. Knepley   PetscAssertPointer(moments, 4);
2656d52c2f21SMatthew G. Knepley   PetscCall(PetscObjectGetComm((PetscObject)sw, &comm));
2657d52c2f21SMatthew G. Knepley   PetscCall(DMSwarmGetField(sw, coordinate, &bsc, &dtc, (void **)&coords));
2658d52c2f21SMatthew G. Knepley   PetscCall(DMSwarmGetField(sw, weight, &bsw, &dtw, (void **)&w));
2659d52c2f21SMatthew G. Knepley   PetscCheck(dtc == PETSC_REAL, comm, PETSC_ERR_ARG_WRONG, "Coordinate field %s must be real, not %s", coordinate, PetscDataTypes[dtc]);
2660d52c2f21SMatthew G. Knepley   PetscCheck(dtw == PETSC_REAL, comm, PETSC_ERR_ARG_WRONG, "Weight field %s must be real, not %s", weight, PetscDataTypes[dtw]);
2661d52c2f21SMatthew G. Knepley   PetscCheck(bsw == 1, comm, PETSC_ERR_ARG_WRONG, "Weight field %s must be a scalar, not blocksize %" PetscInt_FMT, weight, bsw);
2662d52c2f21SMatthew G. Knepley   PetscCall(DMSwarmGetLocalSize(sw, &Np));
2663d52c2f21SMatthew G. Knepley   PetscCall(DMGetWorkArray(sw, bsc + 2, MPIU_REAL, &mom));
2664d52c2f21SMatthew G. Knepley   PetscCall(PetscArrayzero(mom, bsc + 2));
2665d52c2f21SMatthew G. Knepley   for (PetscInt p = 0; p < Np; ++p) {
2666d52c2f21SMatthew G. Knepley     const PetscReal *c  = &coords[p * bsc];
2667d52c2f21SMatthew G. Knepley     const PetscReal  wp = w[p];
2668d52c2f21SMatthew G. Knepley 
2669d52c2f21SMatthew G. Knepley     mom[0] += wp;
2670d52c2f21SMatthew G. Knepley     for (PetscInt d = 0; d < bsc; ++d) {
2671d52c2f21SMatthew G. Knepley       mom[d + 1] += wp * c[d];
2672d52c2f21SMatthew G. Knepley       mom[d + bsc + 1] += wp * PetscSqr(c[d]);
2673d52c2f21SMatthew G. Knepley     }
2674d52c2f21SMatthew G. Knepley   }
2675d52c2f21SMatthew G. Knepley   PetscCall(DMSwarmRestoreField(sw, "velocity", NULL, NULL, (void **)&coords));
2676d52c2f21SMatthew G. Knepley   PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&w));
2677d52c2f21SMatthew G. Knepley   PetscCallMPI(MPIU_Allreduce(mom, moments, bsc + 2, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)sw)));
2678d52c2f21SMatthew G. Knepley   PetscCall(DMRestoreWorkArray(sw, bsc + 2, MPIU_REAL, &mom));
2679d0c080abSJoseph Pusztay   PetscFunctionReturn(PETSC_SUCCESS);
2680d0c080abSJoseph Pusztay }
2681d0c080abSJoseph Pusztay 
2682ce78bad3SBarry Smith static PetscErrorCode DMSetFromOptions_Swarm(DM dm, PetscOptionItems PetscOptionsObject)
2683fca3fa51SMatthew G. Knepley {
2684fca3fa51SMatthew G. Knepley   DM_Swarm *swarm = (DM_Swarm *)dm->data;
2685fca3fa51SMatthew G. Knepley 
2686fca3fa51SMatthew G. Knepley   PetscFunctionBegin;
2687fca3fa51SMatthew G. Knepley   PetscOptionsHeadBegin(PetscOptionsObject, "DMSwarm Options");
2688fca3fa51SMatthew G. Knepley   PetscCall(PetscOptionsEnum("-dm_swarm_remap_type", "Remap algorithm", "DMSwarmSetRemapType", DMSwarmRemapTypeNames, (PetscEnum)swarm->remap_type, (PetscEnum *)&swarm->remap_type, NULL));
2689fca3fa51SMatthew G. Knepley   PetscOptionsHeadEnd();
2690fca3fa51SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
2691fca3fa51SMatthew G. Knepley }
2692fca3fa51SMatthew G. Knepley 
2693d0c080abSJoseph Pusztay PETSC_INTERN PetscErrorCode DMClone_Swarm(DM, DM *);
2694d0c080abSJoseph Pusztay 
2695d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMInitialize_Swarm(DM sw)
2696d71ae5a4SJacob Faibussowitsch {
2697d0c080abSJoseph Pusztay   PetscFunctionBegin;
2698d0c080abSJoseph Pusztay   sw->ops->view                     = DMView_Swarm;
2699d0c080abSJoseph Pusztay   sw->ops->load                     = NULL;
2700fca3fa51SMatthew G. Knepley   sw->ops->setfromoptions           = DMSetFromOptions_Swarm;
2701d0c080abSJoseph Pusztay   sw->ops->clone                    = DMClone_Swarm;
2702d0c080abSJoseph Pusztay   sw->ops->setup                    = DMSetup_Swarm;
2703d0c080abSJoseph Pusztay   sw->ops->createlocalsection       = NULL;
2704adc21957SMatthew G. Knepley   sw->ops->createsectionpermutation = NULL;
2705d0c080abSJoseph Pusztay   sw->ops->createdefaultconstraints = NULL;
2706d0c080abSJoseph Pusztay   sw->ops->createglobalvector       = DMCreateGlobalVector_Swarm;
2707d0c080abSJoseph Pusztay   sw->ops->createlocalvector        = DMCreateLocalVector_Swarm;
2708d0c080abSJoseph Pusztay   sw->ops->getlocaltoglobalmapping  = NULL;
2709d0c080abSJoseph Pusztay   sw->ops->createfieldis            = NULL;
2710d0c080abSJoseph Pusztay   sw->ops->createcoordinatedm       = NULL;
271199acd26cSksagiyam   sw->ops->createcellcoordinatedm   = NULL;
2712d0c080abSJoseph Pusztay   sw->ops->getcoloring              = NULL;
2713d0c080abSJoseph Pusztay   sw->ops->creatematrix             = DMCreateMatrix_Swarm;
2714d0c080abSJoseph Pusztay   sw->ops->createinterpolation      = NULL;
2715d0c080abSJoseph Pusztay   sw->ops->createinjection          = NULL;
2716d0c080abSJoseph Pusztay   sw->ops->createmassmatrix         = DMCreateMassMatrix_Swarm;
27171898fd5cSMatthew G. Knepley   sw->ops->creategradientmatrix     = DMCreateGradientMatrix_Swarm;
2718d0c080abSJoseph Pusztay   sw->ops->refine                   = NULL;
2719d0c080abSJoseph Pusztay   sw->ops->coarsen                  = NULL;
2720d0c080abSJoseph Pusztay   sw->ops->refinehierarchy          = NULL;
2721d0c080abSJoseph Pusztay   sw->ops->coarsenhierarchy         = NULL;
27229d50c5ebSMatthew G. Knepley   sw->ops->globaltolocalbegin       = DMGlobalToLocalBegin_Swarm;
27239d50c5ebSMatthew G. Knepley   sw->ops->globaltolocalend         = DMGlobalToLocalEnd_Swarm;
27249d50c5ebSMatthew G. Knepley   sw->ops->localtoglobalbegin       = DMLocalToGlobalBegin_Swarm;
27259d50c5ebSMatthew G. Knepley   sw->ops->localtoglobalend         = DMLocalToGlobalEnd_Swarm;
2726d0c080abSJoseph Pusztay   sw->ops->destroy                  = DMDestroy_Swarm;
2727d0c080abSJoseph Pusztay   sw->ops->createsubdm              = NULL;
2728d0c080abSJoseph Pusztay   sw->ops->getdimpoints             = NULL;
2729d0c080abSJoseph Pusztay   sw->ops->locatepoints             = NULL;
27309d50c5ebSMatthew G. Knepley   sw->ops->projectfieldlocal        = DMProjectFieldLocal_Swarm;
27313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2732d0c080abSJoseph Pusztay }
2733d0c080abSJoseph Pusztay 
2734d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode DMClone_Swarm(DM dm, DM *newdm)
2735d71ae5a4SJacob Faibussowitsch {
2736d0c080abSJoseph Pusztay   DM_Swarm *swarm = (DM_Swarm *)dm->data;
2737d0c080abSJoseph Pusztay 
2738d0c080abSJoseph Pusztay   PetscFunctionBegin;
2739d0c080abSJoseph Pusztay   swarm->refct++;
2740d0c080abSJoseph Pusztay   (*newdm)->data = swarm;
27419566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)*newdm, DMSWARM));
27429566063dSJacob Faibussowitsch   PetscCall(DMInitialize_Swarm(*newdm));
27432e56dffeSJoe Pusztay   (*newdm)->dim = dm->dim;
27443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2745d0c080abSJoseph Pusztay }
2746d0c080abSJoseph Pusztay 
2747d3a51819SDave May /*MC
27480b4b7b1cSBarry Smith  DMSWARM = "swarm" - A `DM` object for particle methods, such as particle-in-cell (PIC), in which the underlying
27490b4b7b1cSBarry Smith            data is both (i) dynamic in length, (ii) and of arbitrary data type.
2750d3a51819SDave May 
275120f4b53cSBarry Smith  Level: intermediate
275220f4b53cSBarry Smith 
275320f4b53cSBarry Smith  Notes:
27540b4b7b1cSBarry Smith  User data can be represented by `DMSWARM` through a registering "fields" which are to be stored on particles.
275562741f57SDave May  To register a field, the user must provide;
275662741f57SDave May  (a) a unique name;
275762741f57SDave May  (b) the data type (or size in bytes);
275862741f57SDave May  (c) the block size of the data.
2759d3a51819SDave May 
2760d3a51819SDave May  For example, suppose the application requires a unique id, energy, momentum and density to be stored
2761c215a47eSMatthew Knepley  on a set of particles. Then the following code could be used
276220f4b53cSBarry Smith .vb
276320f4b53cSBarry Smith     DMSwarmInitializeFieldRegister(dm)
276420f4b53cSBarry Smith     DMSwarmRegisterPetscDatatypeField(dm,"uid",1,PETSC_LONG);
276520f4b53cSBarry Smith     DMSwarmRegisterPetscDatatypeField(dm,"energy",1,PETSC_REAL);
276620f4b53cSBarry Smith     DMSwarmRegisterPetscDatatypeField(dm,"momentum",3,PETSC_REAL);
276720f4b53cSBarry Smith     DMSwarmRegisterPetscDatatypeField(dm,"density",1,PETSC_FLOAT);
276820f4b53cSBarry Smith     DMSwarmFinalizeFieldRegister(dm)
276920f4b53cSBarry Smith .ve
2770d3a51819SDave May 
277120f4b53cSBarry Smith  The fields represented by `DMSWARM` are dynamic and can be re-sized at any time.
27720b4b7b1cSBarry Smith  The only restriction imposed by `DMSWARM` is that all fields contain the same number of particles.
2773d3a51819SDave May 
2774d3a51819SDave May  To support particle methods, "migration" techniques are provided. These methods migrate data
27755627991aSBarry Smith  between ranks.
2776d3a51819SDave May 
277720f4b53cSBarry Smith  `DMSWARM` supports the methods `DMCreateGlobalVector()` and `DMCreateLocalVector()`.
277820f4b53cSBarry Smith  As a `DMSWARM` may internally define and store values of different data types,
277920f4b53cSBarry Smith  before calling `DMCreateGlobalVector()` or `DMCreateLocalVector()`, the user must inform `DMSWARM` which
27800b4b7b1cSBarry Smith  fields should be used to define a `Vec` object via `DMSwarmVectorDefineField()`
2781c215a47eSMatthew Knepley  The specified field can be changed at any time - thereby permitting vectors
2782c215a47eSMatthew Knepley  compatible with different fields to be created.
2783d3a51819SDave May 
27840b4b7b1cSBarry Smith  A dual representation of fields in the `DMSWARM` and a Vec object is permitted via `DMSwarmCreateGlobalVectorFromField()`
27850b4b7b1cSBarry Smith  Here the data defining the field in the `DMSWARM` is shared with a `Vec`.
2786d3a51819SDave May  This is inherently unsafe if you alter the size of the field at any time between
278720f4b53cSBarry Smith  calls to `DMSwarmCreateGlobalVectorFromField()` and `DMSwarmDestroyGlobalVectorFromField()`.
278820f4b53cSBarry Smith  If the local size of the `DMSWARM` does not match the local size of the global vector
278920f4b53cSBarry Smith  when `DMSwarmDestroyGlobalVectorFromField()` is called, an error is thrown.
2790d3a51819SDave May 
27910b4b7b1cSBarry Smith  Additional high-level support is provided for Particle-In-Cell methods. Refer to `DMSwarmSetType()`.
279262741f57SDave May 
27930b4b7b1cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMType`, `DMCreate()`, `DMSetType()`, `DMSwarmSetType()`, `DMSwarmType`, `DMSwarmCreateGlobalVectorFromField()`,
27940b4b7b1cSBarry Smith          `DMCreateGlobalVector()`, `DMCreateLocalVector()`
2795d3a51819SDave May M*/
2796cc4c1da9SBarry Smith 
2797d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMCreate_Swarm(DM dm)
2798d71ae5a4SJacob Faibussowitsch {
279957795646SDave May   DM_Swarm *swarm;
280057795646SDave May 
280157795646SDave May   PetscFunctionBegin;
280257795646SDave May   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
28034dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&swarm));
2804f0cdbbbaSDave May   dm->data = swarm;
28059566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketCreate(&swarm->db));
28069566063dSJacob Faibussowitsch   PetscCall(DMSwarmInitializeFieldRegister(dm));
2807d52c2f21SMatthew G. Knepley   dm->dim                               = 0;
2808d0c080abSJoseph Pusztay   swarm->refct                          = 1;
28093454631fSDave May   swarm->issetup                        = PETSC_FALSE;
2810480eef7bSDave May   swarm->swarm_type                     = DMSWARM_BASIC;
2811480eef7bSDave May   swarm->migrate_type                   = DMSWARM_MIGRATE_BASIC;
2812480eef7bSDave May   swarm->collect_type                   = DMSWARM_COLLECT_BASIC;
281340c453e9SDave May   swarm->migrate_error_on_missing_point = PETSC_FALSE;
2814f0cdbbbaSDave May   swarm->collect_view_active            = PETSC_FALSE;
2815f0cdbbbaSDave May   swarm->collect_view_reset_nlocal      = -1;
28169566063dSJacob Faibussowitsch   PetscCall(DMInitialize_Swarm(dm));
28172e956fe4SStefano Zampini   if (SwarmDataFieldId == -1) PetscCall(PetscObjectComposedDataRegister(&SwarmDataFieldId));
28183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
281957795646SDave May }
282019307e5cSMatthew G. Knepley 
282119307e5cSMatthew G. Knepley /* Replace dm with the contents of ndm, and then destroy ndm
282219307e5cSMatthew G. Knepley    - Share the DM_Swarm structure
282319307e5cSMatthew G. Knepley */
282419307e5cSMatthew G. Knepley PetscErrorCode DMSwarmReplace(DM dm, DM *ndm)
282519307e5cSMatthew G. Knepley {
282619307e5cSMatthew G. Knepley   DM               dmNew = *ndm;
282719307e5cSMatthew G. Knepley   const PetscReal *maxCell, *Lstart, *L;
282819307e5cSMatthew G. Knepley   PetscInt         dim;
282919307e5cSMatthew G. Knepley 
283019307e5cSMatthew G. Knepley   PetscFunctionBegin;
283119307e5cSMatthew G. Knepley   if (dm == dmNew) {
283219307e5cSMatthew G. Knepley     PetscCall(DMDestroy(ndm));
283319307e5cSMatthew G. Knepley     PetscFunctionReturn(PETSC_SUCCESS);
283419307e5cSMatthew G. Knepley   }
283519307e5cSMatthew G. Knepley   dm->setupcalled = dmNew->setupcalled;
283619307e5cSMatthew G. Knepley   if (!dm->hdr.name) {
283719307e5cSMatthew G. Knepley     const char *name;
283819307e5cSMatthew G. Knepley 
283919307e5cSMatthew G. Knepley     PetscCall(PetscObjectGetName((PetscObject)*ndm, &name));
284019307e5cSMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)dm, name));
284119307e5cSMatthew G. Knepley   }
284219307e5cSMatthew G. Knepley   PetscCall(DMGetDimension(dmNew, &dim));
284319307e5cSMatthew G. Knepley   PetscCall(DMSetDimension(dm, dim));
284419307e5cSMatthew G. Knepley   PetscCall(DMGetPeriodicity(dmNew, &maxCell, &Lstart, &L));
284519307e5cSMatthew G. Knepley   PetscCall(DMSetPeriodicity(dm, maxCell, Lstart, L));
284619307e5cSMatthew G. Knepley   PetscCall(DMDestroy_Swarm(dm));
284719307e5cSMatthew G. Knepley   PetscCall(DMInitialize_Swarm(dm));
284819307e5cSMatthew G. Knepley   dm->data = dmNew->data;
284919307e5cSMatthew G. Knepley   ((DM_Swarm *)dmNew->data)->refct++;
285019307e5cSMatthew G. Knepley   PetscCall(DMDestroy(ndm));
285119307e5cSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
285219307e5cSMatthew G. Knepley }
2853fca3fa51SMatthew G. Knepley 
2854fca3fa51SMatthew G. Knepley /*@
2855fca3fa51SMatthew G. Knepley   DMSwarmDuplicate - Creates a new `DMSWARM` with the same fields and cell `DM`s but no particles
2856fca3fa51SMatthew G. Knepley 
2857fca3fa51SMatthew G. Knepley   Collective
2858fca3fa51SMatthew G. Knepley 
2859fca3fa51SMatthew G. Knepley   Input Parameter:
2860fca3fa51SMatthew G. Knepley . sw - the `DMSWARM`
2861fca3fa51SMatthew G. Knepley 
2862fca3fa51SMatthew G. Knepley   Output Parameter:
2863fca3fa51SMatthew G. Knepley . nsw - the new `DMSWARM`
2864fca3fa51SMatthew G. Knepley 
2865fca3fa51SMatthew G. Knepley   Level: beginner
2866fca3fa51SMatthew G. Knepley 
2867fca3fa51SMatthew G. Knepley .seealso: `DM`, `DMSWARM`, `DMSwarmCreate()`, `DMClone()`
2868fca3fa51SMatthew G. Knepley @*/
2869fca3fa51SMatthew G. Knepley PetscErrorCode DMSwarmDuplicate(DM sw, DM *nsw)
2870fca3fa51SMatthew G. Knepley {
2871fca3fa51SMatthew G. Knepley   DM_Swarm         *swarm = (DM_Swarm *)sw->data;
2872fca3fa51SMatthew G. Knepley   DMSwarmDataField *fields;
2873fca3fa51SMatthew G. Knepley   DMSwarmCellDM     celldm, ncelldm;
2874fca3fa51SMatthew G. Knepley   DMSwarmType       stype;
2875fca3fa51SMatthew G. Knepley   const char       *name, **celldmnames;
2876fca3fa51SMatthew G. Knepley   void             *ctx;
2877fca3fa51SMatthew G. Knepley   PetscInt          dim, Nf, Ndm;
2878fca3fa51SMatthew G. Knepley   PetscBool         flg;
2879fca3fa51SMatthew G. Knepley 
2880fca3fa51SMatthew G. Knepley   PetscFunctionBegin;
2881fca3fa51SMatthew G. Knepley   PetscCall(DMCreate(PetscObjectComm((PetscObject)sw), nsw));
2882fca3fa51SMatthew G. Knepley   PetscCall(DMSetType(*nsw, DMSWARM));
2883fca3fa51SMatthew G. Knepley   PetscCall(PetscObjectGetName((PetscObject)sw, &name));
2884fca3fa51SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)*nsw, name));
2885fca3fa51SMatthew G. Knepley   PetscCall(DMGetDimension(sw, &dim));
2886fca3fa51SMatthew G. Knepley   PetscCall(DMSetDimension(*nsw, dim));
2887fca3fa51SMatthew G. Knepley   PetscCall(DMSwarmGetType(sw, &stype));
2888fca3fa51SMatthew G. Knepley   PetscCall(DMSwarmSetType(*nsw, stype));
2889fca3fa51SMatthew G. Knepley   PetscCall(DMGetApplicationContext(sw, &ctx));
2890fca3fa51SMatthew G. Knepley   PetscCall(DMSetApplicationContext(*nsw, ctx));
2891fca3fa51SMatthew G. Knepley 
2892fca3fa51SMatthew G. Knepley   PetscCall(DMSwarmDataBucketGetDMSwarmDataFields(swarm->db, &Nf, &fields));
2893fca3fa51SMatthew G. Knepley   for (PetscInt f = 0; f < Nf; ++f) {
2894fca3fa51SMatthew G. Knepley     PetscCall(DMSwarmDataFieldStringInList(fields[f]->name, ((DM_Swarm *)(*nsw)->data)->db->nfields, (const DMSwarmDataField *)((DM_Swarm *)(*nsw)->data)->db->field, &flg));
2895fca3fa51SMatthew G. Knepley     if (!flg) PetscCall(DMSwarmRegisterPetscDatatypeField(*nsw, fields[f]->name, fields[f]->bs, fields[f]->petsc_type));
2896fca3fa51SMatthew G. Knepley   }
2897fca3fa51SMatthew G. Knepley 
2898fca3fa51SMatthew G. Knepley   PetscCall(DMSwarmGetCellDMNames(sw, &Ndm, &celldmnames));
2899fca3fa51SMatthew G. Knepley   for (PetscInt c = 0; c < Ndm; ++c) {
2900fca3fa51SMatthew G. Knepley     DM           dm;
2901fca3fa51SMatthew G. Knepley     PetscInt     Ncf;
2902fca3fa51SMatthew G. Knepley     const char **coordfields, **fields;
2903fca3fa51SMatthew G. Knepley 
2904fca3fa51SMatthew G. Knepley     PetscCall(DMSwarmGetCellDMByName(sw, celldmnames[c], &celldm));
2905fca3fa51SMatthew G. Knepley     PetscCall(DMSwarmCellDMGetDM(celldm, &dm));
2906fca3fa51SMatthew G. Knepley     PetscCall(DMSwarmCellDMGetCoordinateFields(celldm, &Ncf, &coordfields));
2907fca3fa51SMatthew G. Knepley     PetscCall(DMSwarmCellDMGetFields(celldm, &Nf, &fields));
2908fca3fa51SMatthew G. Knepley     PetscCall(DMSwarmCellDMCreate(dm, Nf, fields, Ncf, coordfields, &ncelldm));
2909fca3fa51SMatthew G. Knepley     PetscCall(DMSwarmAddCellDM(*nsw, ncelldm));
2910fca3fa51SMatthew G. Knepley     PetscCall(DMSwarmCellDMDestroy(&ncelldm));
2911fca3fa51SMatthew G. Knepley   }
2912fca3fa51SMatthew G. Knepley   PetscCall(PetscFree(celldmnames));
2913fca3fa51SMatthew G. Knepley 
2914fca3fa51SMatthew G. Knepley   PetscCall(DMSetFromOptions(*nsw));
2915fca3fa51SMatthew G. Knepley   PetscCall(DMSetUp(*nsw));
2916fca3fa51SMatthew G. Knepley   PetscCall(DMSwarmGetCellDMActive(sw, &celldm));
2917fca3fa51SMatthew G. Knepley   PetscCall(PetscObjectGetName((PetscObject)celldm, &name));
2918fca3fa51SMatthew G. Knepley   PetscCall(DMSwarmSetCellDMActive(*nsw, name));
2919fca3fa51SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
2920fca3fa51SMatthew G. Knepley }
29219d50c5ebSMatthew G. Knepley 
29229d50c5ebSMatthew G. Knepley PetscErrorCode DMLocalToGlobalBegin_Swarm(DM dm, Vec l, InsertMode mode, Vec g)
29239d50c5ebSMatthew G. Knepley {
29249d50c5ebSMatthew G. Knepley   PetscFunctionBegin;
29259d50c5ebSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
29269d50c5ebSMatthew G. Knepley }
29279d50c5ebSMatthew G. Knepley 
29289d50c5ebSMatthew G. Knepley PetscErrorCode DMLocalToGlobalEnd_Swarm(DM dm, Vec l, InsertMode mode, Vec g)
29299d50c5ebSMatthew G. Knepley {
29309d50c5ebSMatthew G. Knepley   PetscFunctionBegin;
29319d50c5ebSMatthew G. Knepley   switch (mode) {
29329d50c5ebSMatthew G. Knepley   case INSERT_VALUES:
29339d50c5ebSMatthew G. Knepley     PetscCall(VecCopy(l, g));
29349d50c5ebSMatthew G. Knepley     break;
29359d50c5ebSMatthew G. Knepley   case ADD_VALUES:
29369d50c5ebSMatthew G. Knepley     PetscCall(VecAXPY(g, 1., l));
29379d50c5ebSMatthew G. Knepley     break;
29389d50c5ebSMatthew G. Knepley   default:
29399d50c5ebSMatthew G. Knepley     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mode not supported: %d", mode);
29409d50c5ebSMatthew G. Knepley   }
29419d50c5ebSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
29429d50c5ebSMatthew G. Knepley }
29439d50c5ebSMatthew G. Knepley 
29449d50c5ebSMatthew G. Knepley PetscErrorCode DMGlobalToLocalBegin_Swarm(DM dm, Vec g, InsertMode mode, Vec l)
29459d50c5ebSMatthew G. Knepley {
29469d50c5ebSMatthew G. Knepley   PetscFunctionBegin;
29479d50c5ebSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
29489d50c5ebSMatthew G. Knepley }
29499d50c5ebSMatthew G. Knepley 
29509d50c5ebSMatthew G. Knepley PetscErrorCode DMGlobalToLocalEnd_Swarm(DM dm, Vec g, InsertMode mode, Vec l)
29519d50c5ebSMatthew G. Knepley {
29529d50c5ebSMatthew G. Knepley   PetscFunctionBegin;
29539d50c5ebSMatthew G. Knepley   PetscCall(DMLocalToGlobalEnd_Swarm(dm, g, mode, l));
29549d50c5ebSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
29559d50c5ebSMatthew G. Knepley }
2956