xref: /petsc/src/dm/impls/swarm/swarm.c (revision eb0c6899d8cb20e657fd4f83d7a531900bc59e0a)
157795646SDave May #define PETSCDM_DLL
257795646SDave May #include <petsc/private/dmswarmimpl.h> /*I   "petscdmswarm.h"   I*/
3e8f14785SLisandro Dalcin #include <petsc/private/hashsetij.h>
40643ed39SMatthew G. Knepley #include <petsc/private/petscfeimpl.h>
55917a6f0SStefano Zampini #include <petscviewer.h>
65917a6f0SStefano Zampini #include <petscdraw.h>
783c47955SMatthew G. Knepley #include <petscdmplex.h>
84cc7f7b2SMatthew G. Knepley #include <petscblaslapack.h>
9279f676cSBarry Smith #include "../src/dm/impls/swarm/data_bucket.h"
10d0c080abSJoseph Pusztay #include <petscdmlabel.h>
11d0c080abSJoseph Pusztay #include <petscsection.h>
1257795646SDave May 
13f2b2bee7SDave May PetscLogEvent DMSWARM_Migrate, DMSWARM_SetSizes, DMSWARM_AddPoints, DMSWARM_RemovePoints, DMSWARM_Sort;
14ed923d71SDave May PetscLogEvent DMSWARM_DataExchangerTopologySetup, DMSWARM_DataExchangerBegin, DMSWARM_DataExchangerEnd;
15ed923d71SDave May PetscLogEvent DMSWARM_DataExchangerSendCount, DMSWARM_DataExchangerPack;
16ed923d71SDave May 
17ea78f98cSLisandro Dalcin const char *DMSwarmTypeNames[]          = {"basic", "pic", NULL};
18ea78f98cSLisandro Dalcin const char *DMSwarmMigrateTypeNames[]   = {"basic", "dmcellnscatter", "dmcellexact", "user", NULL};
19ea78f98cSLisandro Dalcin const char *DMSwarmCollectTypeNames[]   = {"basic", "boundingbox", "general", "user", NULL};
20ea78f98cSLisandro Dalcin const char *DMSwarmPICLayoutTypeNames[] = {"regular", "gauss", "subdivision", NULL};
21f0cdbbbaSDave May 
22f0cdbbbaSDave May const char DMSwarmField_pid[]       = "DMSwarm_pid";
23f0cdbbbaSDave May const char DMSwarmField_rank[]      = "DMSwarm_rank";
24f0cdbbbaSDave May const char DMSwarmPICField_coor[]   = "DMSwarmPIC_coor";
25e2d107dbSDave May const char DMSwarmPICField_cellid[] = "DMSwarm_cellid";
26f0cdbbbaSDave May 
272e956fe4SStefano Zampini PetscInt SwarmDataFieldId = -1;
282e956fe4SStefano Zampini 
2974d0cae8SMatthew G. Knepley #if defined(PETSC_HAVE_HDF5)
3074d0cae8SMatthew G. Knepley   #include <petscviewerhdf5.h>
3174d0cae8SMatthew G. Knepley 
3266976f2fSJacob Faibussowitsch static PetscErrorCode VecView_Swarm_HDF5_Internal(Vec v, PetscViewer viewer)
33d71ae5a4SJacob Faibussowitsch {
3474d0cae8SMatthew G. Knepley   DM        dm;
3574d0cae8SMatthew G. Knepley   PetscReal seqval;
3674d0cae8SMatthew G. Knepley   PetscInt  seqnum, bs;
37*eb0c6899SMatthew G. Knepley   PetscBool isseq, ists;
3874d0cae8SMatthew G. Knepley 
3974d0cae8SMatthew G. Knepley   PetscFunctionBegin;
409566063dSJacob Faibussowitsch   PetscCall(VecGetDM(v, &dm));
419566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(v, &bs));
429566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5PushGroup(viewer, "/particle_fields"));
439566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)v, VECSEQ, &isseq));
449566063dSJacob Faibussowitsch   PetscCall(DMGetOutputSequenceNumber(dm, &seqnum, &seqval));
45*eb0c6899SMatthew G. Knepley   PetscCall(PetscViewerHDF5IsTimestepping(viewer, &ists));
46*eb0c6899SMatthew G. Knepley   if (ists) PetscCall(PetscViewerHDF5SetTimestep(viewer, seqnum));
479566063dSJacob Faibussowitsch   /* PetscCall(DMSequenceView_HDF5(dm, "time", seqnum, (PetscScalar) seqval, viewer)); */
489566063dSJacob Faibussowitsch   PetscCall(VecViewNative(v, viewer));
499566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject)v, "Nc", PETSC_INT, (void *)&bs));
509566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5PopGroup(viewer));
513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5274d0cae8SMatthew G. Knepley }
5374d0cae8SMatthew G. Knepley 
5466976f2fSJacob Faibussowitsch static PetscErrorCode DMSwarmView_HDF5(DM dm, PetscViewer viewer)
55d71ae5a4SJacob Faibussowitsch {
5674d0cae8SMatthew G. Knepley   Vec       coordinates;
5774d0cae8SMatthew G. Knepley   PetscInt  Np;
5874d0cae8SMatthew G. Knepley   PetscBool isseq;
5974d0cae8SMatthew G. Knepley 
6074d0cae8SMatthew G. Knepley   PetscFunctionBegin;
619566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetSize(dm, &Np));
629566063dSJacob Faibussowitsch   PetscCall(DMSwarmCreateGlobalVectorFromField(dm, DMSwarmPICField_coor, &coordinates));
639566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates"));
649566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5PushGroup(viewer, "/particles"));
659566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)coordinates, VECSEQ, &isseq));
669566063dSJacob Faibussowitsch   PetscCall(VecViewNative(coordinates, viewer));
679566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject)coordinates, "Np", PETSC_INT, (void *)&Np));
689566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5PopGroup(viewer));
699566063dSJacob Faibussowitsch   PetscCall(DMSwarmDestroyGlobalVectorFromField(dm, DMSwarmPICField_coor, &coordinates));
703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7174d0cae8SMatthew G. Knepley }
7274d0cae8SMatthew G. Knepley #endif
7374d0cae8SMatthew G. Knepley 
7466976f2fSJacob Faibussowitsch static PetscErrorCode VecView_Swarm(Vec v, PetscViewer viewer)
75d71ae5a4SJacob Faibussowitsch {
7674d0cae8SMatthew G. Knepley   DM dm;
77f9558f5cSBarry Smith #if defined(PETSC_HAVE_HDF5)
7874d0cae8SMatthew G. Knepley   PetscBool ishdf5;
79f9558f5cSBarry Smith #endif
8074d0cae8SMatthew G. Knepley 
8174d0cae8SMatthew G. Knepley   PetscFunctionBegin;
829566063dSJacob Faibussowitsch   PetscCall(VecGetDM(v, &dm));
8328b400f6SJacob Faibussowitsch   PetscCheck(dm, PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_WRONG, "Vector not generated from a DM");
84f9558f5cSBarry Smith #if defined(PETSC_HAVE_HDF5)
859566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
8674d0cae8SMatthew G. Knepley   if (ishdf5) {
879566063dSJacob Faibussowitsch     PetscCall(VecView_Swarm_HDF5_Internal(v, viewer));
883ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
8974d0cae8SMatthew G. Knepley   }
90f9558f5cSBarry Smith #endif
919566063dSJacob Faibussowitsch   PetscCall(VecViewNative(v, viewer));
923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9374d0cae8SMatthew G. Knepley }
9474d0cae8SMatthew G. Knepley 
95d3a51819SDave May /*@C
960bf7c1c5SMatthew G. Knepley   DMSwarmVectorGetField - Gets the field from which to define a `Vec` object
970bf7c1c5SMatthew G. Knepley   when `DMCreateLocalVector()`, or `DMCreateGlobalVector()` is called
980bf7c1c5SMatthew G. Knepley 
990bf7c1c5SMatthew G. Knepley   Not collective
1000bf7c1c5SMatthew G. Knepley 
1010bf7c1c5SMatthew G. Knepley   Input Parameter:
1020bf7c1c5SMatthew G. Knepley . dm - a `DMSWARM`
1030bf7c1c5SMatthew G. Knepley 
1040bf7c1c5SMatthew G. Knepley   Output Parameter:
1050bf7c1c5SMatthew G. Knepley . fieldname - the textual name given to a registered field, or the empty string if it has not been set
1060bf7c1c5SMatthew G. Knepley 
1070bf7c1c5SMatthew G. Knepley   Level: beginner
1080bf7c1c5SMatthew G. Knepley 
1090bf7c1c5SMatthew G. Knepley .seealso: `DM`, `DMSWARM`, `DMSwarmVectorDefineField()` `DMSwarmRegisterPetscDatatypeField()`, `DMCreateGlobalVector()`, `DMCreateLocalVector()`
1100bf7c1c5SMatthew G. Knepley @*/
1110bf7c1c5SMatthew G. Knepley PetscErrorCode DMSwarmVectorGetField(DM dm, const char *fieldname[])
1120bf7c1c5SMatthew G. Knepley {
1130bf7c1c5SMatthew G. Knepley   PetscFunctionBegin;
1140bf7c1c5SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
1150bf7c1c5SMatthew G. Knepley   PetscAssertPointer(fieldname, 2);
1160bf7c1c5SMatthew G. Knepley   *fieldname = ((DM_Swarm *)dm->data)->vec_field_name;
1170bf7c1c5SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
1180bf7c1c5SMatthew G. Knepley }
1190bf7c1c5SMatthew G. Knepley 
1200bf7c1c5SMatthew G. Knepley /*@C
12120f4b53cSBarry Smith   DMSwarmVectorDefineField - Sets the field from which to define a `Vec` object
12220f4b53cSBarry Smith   when `DMCreateLocalVector()`, or `DMCreateGlobalVector()` is called
12357795646SDave May 
12420f4b53cSBarry Smith   Collective
12557795646SDave May 
12660225df5SJacob Faibussowitsch   Input Parameters:
12720f4b53cSBarry Smith + dm        - a `DMSWARM`
12862741f57SDave May - fieldname - the textual name given to a registered field
12957795646SDave May 
130d3a51819SDave May   Level: beginner
13157795646SDave May 
132d3a51819SDave May   Notes:
13320f4b53cSBarry Smith   The field with name `fieldname` must be defined as having a data type of `PetscScalar`.
134e7af74c8SDave May 
13520f4b53cSBarry Smith   This function must be called prior to calling `DMCreateLocalVector()`, `DMCreateGlobalVector()`.
13620f4b53cSBarry Smith   Multiple calls to `DMSwarmVectorDefineField()` are permitted.
137e7af74c8SDave May 
1380bf7c1c5SMatthew G. Knepley .seealso: `DM`, `DMSWARM`, `DMSwarmVectorGetField()`, `DMSwarmRegisterPetscDatatypeField()`, `DMCreateGlobalVector()`, `DMCreateLocalVector()`
139d3a51819SDave May @*/
140d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmVectorDefineField(DM dm, const char fieldname[])
141d71ae5a4SJacob Faibussowitsch {
142b5bcf523SDave May   DM_Swarm     *swarm = (DM_Swarm *)dm->data;
143b5bcf523SDave May   PetscInt      bs, n;
144b5bcf523SDave May   PetscScalar  *array;
145b5bcf523SDave May   PetscDataType type;
146b5bcf523SDave May 
147a9cbaee5SMatthew G. Knepley   PetscFunctionBegin;
1480bf7c1c5SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
1490bf7c1c5SMatthew G. Knepley   if (fieldname) PetscAssertPointer(fieldname, 2);
1509566063dSJacob Faibussowitsch   if (!swarm->issetup) PetscCall(DMSetUp(dm));
1519566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &n, NULL, NULL));
1529566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, fieldname, &bs, &type, (void **)&array));
153b5bcf523SDave May 
154b5bcf523SDave May   /* Check all fields are of type PETSC_REAL or PETSC_SCALAR */
15508401ef6SPierre Jolivet   PetscCheck(type == PETSC_REAL, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only valid for PETSC_REAL");
1569566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(swarm->vec_field_name, PETSC_MAX_PATH_LEN - 1, "%s", fieldname));
157b5bcf523SDave May   swarm->vec_field_set    = PETSC_TRUE;
1581b1ea282SDave May   swarm->vec_field_bs     = bs;
159b5bcf523SDave May   swarm->vec_field_nlocal = n;
1609566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm, fieldname, &bs, &type, (void **)&array));
1613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
162b5bcf523SDave May }
163b5bcf523SDave May 
164cc651181SDave May /* requires DMSwarmDefineFieldVector has been called */
16566976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateGlobalVector_Swarm(DM dm, Vec *vec)
166d71ae5a4SJacob Faibussowitsch {
167b5bcf523SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
168b5bcf523SDave May   Vec       x;
169b5bcf523SDave May   char      name[PETSC_MAX_PATH_LEN];
170b5bcf523SDave May 
171a9cbaee5SMatthew G. Knepley   PetscFunctionBegin;
1729566063dSJacob Faibussowitsch   if (!swarm->issetup) PetscCall(DMSetUp(dm));
17328b400f6SJacob Faibussowitsch   PetscCheck(swarm->vec_field_set, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Must call DMSwarmVectorDefineField first");
17408401ef6SPierre Jolivet   PetscCheck(swarm->vec_field_nlocal == swarm->db->L, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "DMSwarm sizes have changed since last call to VectorDefineField first"); /* Stale data */
175cc651181SDave May 
1769566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN - 1, "DMSwarmField_%s", swarm->vec_field_name));
1779566063dSJacob Faibussowitsch   PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), &x));
1789566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)x, name));
1799566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(x, swarm->db->L * swarm->vec_field_bs, PETSC_DETERMINE));
1809566063dSJacob Faibussowitsch   PetscCall(VecSetBlockSize(x, swarm->vec_field_bs));
1819566063dSJacob Faibussowitsch   PetscCall(VecSetDM(x, dm));
1829566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(x));
1839566063dSJacob Faibussowitsch   PetscCall(VecSetDM(x, dm));
1849566063dSJacob Faibussowitsch   PetscCall(VecSetOperation(x, VECOP_VIEW, (void (*)(void))VecView_Swarm));
185b5bcf523SDave May   *vec = x;
1863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
187b5bcf523SDave May }
188b5bcf523SDave May 
189b5bcf523SDave May /* requires DMSwarmDefineFieldVector has been called */
19066976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateLocalVector_Swarm(DM dm, Vec *vec)
191d71ae5a4SJacob Faibussowitsch {
192b5bcf523SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
193b5bcf523SDave May   Vec       x;
194b5bcf523SDave May   char      name[PETSC_MAX_PATH_LEN];
195b5bcf523SDave May 
196a9cbaee5SMatthew G. Knepley   PetscFunctionBegin;
1979566063dSJacob Faibussowitsch   if (!swarm->issetup) PetscCall(DMSetUp(dm));
19828b400f6SJacob Faibussowitsch   PetscCheck(swarm->vec_field_set, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Must call DMSwarmVectorDefineField first");
19908401ef6SPierre Jolivet   PetscCheck(swarm->vec_field_nlocal == swarm->db->L, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "DMSwarm sizes have changed since last call to VectorDefineField first");
200cc651181SDave May 
2019566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN - 1, "DMSwarmField_%s", swarm->vec_field_name));
2029566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_SELF, &x));
2039566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)x, name));
2049566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(x, swarm->db->L * swarm->vec_field_bs, PETSC_DETERMINE));
2059566063dSJacob Faibussowitsch   PetscCall(VecSetBlockSize(x, swarm->vec_field_bs));
2069566063dSJacob Faibussowitsch   PetscCall(VecSetDM(x, dm));
2079566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(x));
208b5bcf523SDave May   *vec = x;
2093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
210b5bcf523SDave May }
211b5bcf523SDave May 
212d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSwarmDestroyVectorFromField_Private(DM dm, const char fieldname[], Vec *vec)
213d71ae5a4SJacob Faibussowitsch {
214fb1bcc12SMatthew G. Knepley   DM_Swarm        *swarm = (DM_Swarm *)dm->data;
21577048351SPatrick Sanan   DMSwarmDataField gfield;
2162e956fe4SStefano Zampini   PetscInt         bs, nlocal, fid = -1, cfid = -2;
2172e956fe4SStefano Zampini   PetscBool        flg;
218d3a51819SDave May 
219fb1bcc12SMatthew G. Knepley   PetscFunctionBegin;
2202e956fe4SStefano Zampini   /* check vector is an inplace array */
2212e956fe4SStefano Zampini   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldIdByName(swarm->db, fieldname, &fid));
2222e956fe4SStefano Zampini   PetscCall(PetscObjectComposedDataGetInt((PetscObject)*vec, SwarmDataFieldId, cfid, flg));
2232e956fe4SStefano Zampini   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);
2249566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(*vec, &nlocal));
2259566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(*vec, &bs));
22608401ef6SPierre 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");
2279566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield));
2289566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataFieldRestoreAccess(gfield));
2298df78e51SMark Adams   PetscCall(VecResetArray(*vec));
2309566063dSJacob Faibussowitsch   PetscCall(VecDestroy(vec));
2313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
232fb1bcc12SMatthew G. Knepley }
233fb1bcc12SMatthew G. Knepley 
234d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSwarmCreateVectorFromField_Private(DM dm, const char fieldname[], MPI_Comm comm, Vec *vec)
235d71ae5a4SJacob Faibussowitsch {
236fb1bcc12SMatthew G. Knepley   DM_Swarm     *swarm = (DM_Swarm *)dm->data;
237fb1bcc12SMatthew G. Knepley   PetscDataType type;
238fb1bcc12SMatthew G. Knepley   PetscScalar  *array;
2392e956fe4SStefano Zampini   PetscInt      bs, n, fid;
240fb1bcc12SMatthew G. Knepley   char          name[PETSC_MAX_PATH_LEN];
241e4fbd051SBarry Smith   PetscMPIInt   size;
2428df78e51SMark Adams   PetscBool     iscuda, iskokkos;
243fb1bcc12SMatthew G. Knepley 
244fb1bcc12SMatthew G. Knepley   PetscFunctionBegin;
2459566063dSJacob Faibussowitsch   if (!swarm->issetup) PetscCall(DMSetUp(dm));
2469566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &n, NULL, NULL));
2479566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, fieldname, &bs, &type, (void **)&array));
24808401ef6SPierre Jolivet   PetscCheck(type == PETSC_REAL, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only valid for PETSC_REAL");
249fb1bcc12SMatthew G. Knepley 
2509566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
2518df78e51SMark Adams   PetscCall(PetscStrcmp(dm->vectype, VECKOKKOS, &iskokkos));
2528df78e51SMark Adams   PetscCall(PetscStrcmp(dm->vectype, VECCUDA, &iscuda));
2538df78e51SMark Adams   PetscCall(VecCreate(comm, vec));
2548df78e51SMark Adams   PetscCall(VecSetSizes(*vec, n * bs, PETSC_DETERMINE));
2558df78e51SMark Adams   PetscCall(VecSetBlockSize(*vec, bs));
2568df78e51SMark Adams   if (iskokkos) PetscCall(VecSetType(*vec, VECKOKKOS));
2578df78e51SMark Adams   else if (iscuda) PetscCall(VecSetType(*vec, VECCUDA));
2588df78e51SMark Adams   else PetscCall(VecSetType(*vec, VECSTANDARD));
2598df78e51SMark Adams   PetscCall(VecPlaceArray(*vec, array));
2608df78e51SMark Adams 
2619566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN - 1, "DMSwarmSharedField_%s", fieldname));
2629566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)*vec, name));
263fb1bcc12SMatthew G. Knepley 
264fb1bcc12SMatthew G. Knepley   /* Set guard */
2652e956fe4SStefano Zampini   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldIdByName(swarm->db, fieldname, &fid));
2662e956fe4SStefano Zampini   PetscCall(PetscObjectComposedDataSetInt((PetscObject)*vec, SwarmDataFieldId, fid));
26774d0cae8SMatthew G. Knepley 
2689566063dSJacob Faibussowitsch   PetscCall(VecSetDM(*vec, dm));
2699566063dSJacob Faibussowitsch   PetscCall(VecSetOperation(*vec, VECOP_VIEW, (void (*)(void))VecView_Swarm));
2703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
271fb1bcc12SMatthew G. Knepley }
272fb1bcc12SMatthew G. Knepley 
2730643ed39SMatthew G. Knepley /* This creates a "mass matrix" between a finite element and particle space. If a finite element interpolant is given by
2740643ed39SMatthew G. Knepley 
2750643ed39SMatthew G. Knepley      \hat f = \sum_i f_i \phi_i
2760643ed39SMatthew G. Knepley 
2770643ed39SMatthew G. Knepley    and a particle function is given by
2780643ed39SMatthew G. Knepley 
2790643ed39SMatthew G. Knepley      f = \sum_i w_i \delta(x - x_i)
2800643ed39SMatthew G. Knepley 
2810643ed39SMatthew G. Knepley    then we want to require that
2820643ed39SMatthew G. Knepley 
2830643ed39SMatthew G. Knepley      M \hat f = M_p f
2840643ed39SMatthew G. Knepley 
2850643ed39SMatthew G. Knepley    where the particle mass matrix is given by
2860643ed39SMatthew G. Knepley 
2870643ed39SMatthew G. Knepley      (M_p)_{ij} = \int \phi_i \delta(x - x_j)
2880643ed39SMatthew G. Knepley 
2890643ed39SMatthew G. Knepley    The way Dave May does particles, they amount to quadratue weights rather than delta functions, so he has |J| is in
2900643ed39SMatthew G. Knepley    his integral. We allow this with the boolean flag.
2910643ed39SMatthew G. Knepley */
292d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSwarmComputeMassMatrix_Private(DM dmc, DM dmf, Mat mass, PetscBool useDeltaFunction, void *ctx)
293d71ae5a4SJacob Faibussowitsch {
29483c47955SMatthew G. Knepley   const char  *name = "Mass Matrix";
2950643ed39SMatthew G. Knepley   MPI_Comm     comm;
29683c47955SMatthew G. Knepley   PetscDS      prob;
29783c47955SMatthew G. Knepley   PetscSection fsection, globalFSection;
298e8f14785SLisandro Dalcin   PetscHSetIJ  ht;
2990643ed39SMatthew G. Knepley   PetscLayout  rLayout, colLayout;
30083c47955SMatthew G. Knepley   PetscInt    *dnz, *onz;
301adb2528bSMark Adams   PetscInt     locRows, locCols, rStart, colStart, colEnd, *rowIDXs;
3020643ed39SMatthew G. Knepley   PetscReal   *xi, *v0, *J, *invJ, detJ = 1.0, v0ref[3] = {-1.0, -1.0, -1.0};
30383c47955SMatthew G. Knepley   PetscScalar *elemMat;
3040bf7c1c5SMatthew G. Knepley   PetscInt     dim, Nf, field, cStart, cEnd, cell, totDim, maxC = 0, totNc = 0;
30583c47955SMatthew G. Knepley 
30683c47955SMatthew G. Knepley   PetscFunctionBegin;
3079566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)mass, &comm));
3089566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dmf, &dim));
3099566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dmf, &prob));
3109566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(prob, &Nf));
3119566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(prob, &totDim));
3129566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(dim, &v0, dim * dim, &J, dim * dim, &invJ));
3139566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dmf, &fsection));
3149566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalSection(dmf, &globalFSection));
3159566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd));
3169566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(mass, &locRows, &locCols));
31783c47955SMatthew G. Knepley 
3189566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreate(comm, &colLayout));
3199566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetLocalSize(colLayout, locCols));
3209566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(colLayout, 1));
3219566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(colLayout));
3229566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetRange(colLayout, &colStart, &colEnd));
3239566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&colLayout));
3240643ed39SMatthew G. Knepley 
3259566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreate(comm, &rLayout));
3269566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetLocalSize(rLayout, locRows));
3279566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(rLayout, 1));
3289566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(rLayout));
3299566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetRange(rLayout, &rStart, NULL));
3309566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&rLayout));
3310643ed39SMatthew G. Knepley 
3329566063dSJacob Faibussowitsch   PetscCall(PetscCalloc2(locRows, &dnz, locRows, &onz));
3339566063dSJacob Faibussowitsch   PetscCall(PetscHSetIJCreate(&ht));
33453e60ab4SJoseph Pusztay 
3359566063dSJacob Faibussowitsch   PetscCall(PetscSynchronizedFlush(comm, NULL));
3360bf7c1c5SMatthew G. Knepley   for (field = 0; field < Nf; ++field) {
3370bf7c1c5SMatthew G. Knepley     PetscObject  obj;
3380bf7c1c5SMatthew G. Knepley     PetscClassId id;
3390bf7c1c5SMatthew G. Knepley     PetscInt     Nc;
3400bf7c1c5SMatthew G. Knepley 
3410bf7c1c5SMatthew G. Knepley     PetscCall(PetscDSGetDiscretization(prob, field, &obj));
3420bf7c1c5SMatthew G. Knepley     PetscCall(PetscObjectGetClassId(obj, &id));
3430bf7c1c5SMatthew G. Knepley     if (id == PETSCFE_CLASSID) PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc));
3440bf7c1c5SMatthew G. Knepley     else PetscCall(PetscFVGetNumComponents((PetscFV)obj, &Nc));
3450bf7c1c5SMatthew G. Knepley     totNc += Nc;
3460bf7c1c5SMatthew G. Knepley   }
3470643ed39SMatthew G. Knepley   /* count non-zeros */
3489566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortGetAccess(dmc));
34983c47955SMatthew G. Knepley   for (field = 0; field < Nf; ++field) {
3500bf7c1c5SMatthew G. Knepley     PetscObject  obj;
3510bf7c1c5SMatthew G. Knepley     PetscClassId id;
3520bf7c1c5SMatthew G. Knepley     PetscInt     Nc;
3530bf7c1c5SMatthew G. Knepley 
3540bf7c1c5SMatthew G. Knepley     PetscCall(PetscDSGetDiscretization(prob, field, &obj));
3550bf7c1c5SMatthew G. Knepley     PetscCall(PetscObjectGetClassId(obj, &id));
3560bf7c1c5SMatthew G. Knepley     if (id == PETSCFE_CLASSID) PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc));
3570bf7c1c5SMatthew G. Knepley     else PetscCall(PetscFVGetNumComponents((PetscFV)obj, &Nc));
3580bf7c1c5SMatthew G. Knepley 
35983c47955SMatthew G. Knepley     for (cell = cStart; cell < cEnd; ++cell) {
3600643ed39SMatthew G. Knepley       PetscInt *findices, *cindices; /* fine is vertices, coarse is particles */
36183c47955SMatthew G. Knepley       PetscInt  numFIndices, numCIndices;
36283c47955SMatthew G. Knepley 
3639566063dSJacob Faibussowitsch       PetscCall(DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
3649566063dSJacob Faibussowitsch       PetscCall(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices));
365fc7c92abSMatthew G. Knepley       maxC = PetscMax(maxC, numCIndices);
36683c47955SMatthew G. Knepley       {
367e8f14785SLisandro Dalcin         PetscHashIJKey key;
368e8f14785SLisandro Dalcin         PetscBool      missing;
3690bf7c1c5SMatthew G. Knepley         for (PetscInt i = 0; i < numFIndices; ++i) {
370adb2528bSMark Adams           key.j = findices[i]; /* global column (from Plex) */
371adb2528bSMark Adams           if (key.j >= 0) {
37283c47955SMatthew G. Knepley             /* Get indices for coarse elements */
3730bf7c1c5SMatthew G. Knepley             for (PetscInt j = 0; j < numCIndices; ++j) {
3740bf7c1c5SMatthew G. Knepley               for (PetscInt c = 0; c < Nc; ++c) {
3750bf7c1c5SMatthew G. Knepley                 // TODO Need field offset on particle here
3760bf7c1c5SMatthew G. Knepley                 key.i = cindices[j] * totNc + c + rStart; /* global cols (from Swarm) */
377adb2528bSMark Adams                 if (key.i < 0) continue;
3789566063dSJacob Faibussowitsch                 PetscCall(PetscHSetIJQueryAdd(ht, key, &missing));
37983c47955SMatthew G. Knepley                 if (missing) {
3800643ed39SMatthew G. Knepley                   if ((key.j >= colStart) && (key.j < colEnd)) ++dnz[key.i - rStart];
381e8f14785SLisandro Dalcin                   else ++onz[key.i - rStart];
38263a3b9bcSJacob Faibussowitsch                 } else SETERRQ(PetscObjectComm((PetscObject)dmf), PETSC_ERR_SUP, "Set new value at %" PetscInt_FMT ",%" PetscInt_FMT, key.i, key.j);
38383c47955SMatthew G. Knepley               }
384fc7c92abSMatthew G. Knepley             }
385fc7c92abSMatthew G. Knepley           }
3860bf7c1c5SMatthew G. Knepley         }
3879566063dSJacob Faibussowitsch         PetscCall(PetscFree(cindices));
38883c47955SMatthew G. Knepley       }
3899566063dSJacob Faibussowitsch       PetscCall(DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
39083c47955SMatthew G. Knepley     }
39183c47955SMatthew G. Knepley   }
3929566063dSJacob Faibussowitsch   PetscCall(PetscHSetIJDestroy(&ht));
3939566063dSJacob Faibussowitsch   PetscCall(MatXAIJSetPreallocation(mass, 1, dnz, onz, NULL, NULL));
3949566063dSJacob Faibussowitsch   PetscCall(MatSetOption(mass, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
3959566063dSJacob Faibussowitsch   PetscCall(PetscFree2(dnz, onz));
3960bf7c1c5SMatthew G. Knepley   PetscCall(PetscMalloc3(maxC * totNc * totDim, &elemMat, maxC * totNc, &rowIDXs, maxC * dim, &xi));
39783c47955SMatthew G. Knepley   for (field = 0; field < Nf; ++field) {
398ef0bb6c7SMatthew G. Knepley     PetscTabulation Tcoarse;
39983c47955SMatthew G. Knepley     PetscObject     obj;
400ad9634fcSMatthew G. Knepley     PetscClassId    id;
401ea7032a0SMatthew G. Knepley     PetscReal      *fieldVals;
4020bf7c1c5SMatthew G. Knepley     PetscInt        Nc;
40383c47955SMatthew G. Knepley 
4049566063dSJacob Faibussowitsch     PetscCall(PetscDSGetDiscretization(prob, field, &obj));
405ad9634fcSMatthew G. Knepley     PetscCall(PetscObjectGetClassId(obj, &id));
406ad9634fcSMatthew G. Knepley     if (id == PETSCFE_CLASSID) PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc));
407ad9634fcSMatthew G. Knepley     else PetscCall(PetscFVGetNumComponents((PetscFV)obj, &Nc));
408ea7032a0SMatthew G. Knepley 
409ea7032a0SMatthew G. Knepley     PetscCall(DMSwarmGetField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **)&fieldVals));
41083c47955SMatthew G. Knepley     for (cell = cStart; cell < cEnd; ++cell) {
41183c47955SMatthew G. Knepley       PetscInt *findices, *cindices;
41283c47955SMatthew G. Knepley       PetscInt  numFIndices, numCIndices;
41383c47955SMatthew G. Knepley 
4140643ed39SMatthew G. Knepley       /* TODO: Use DMField instead of assuming affine */
4159566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ));
4169566063dSJacob Faibussowitsch       PetscCall(DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
4179566063dSJacob Faibussowitsch       PetscCall(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices));
4180bf7c1c5SMatthew G. Knepley       for (PetscInt j = 0; j < numCIndices; ++j) CoordinatesRealToRef(dim, dim, v0ref, v0, invJ, &fieldVals[cindices[j] * dim], &xi[j * dim]);
419ad9634fcSMatthew G. Knepley       if (id == PETSCFE_CLASSID) PetscCall(PetscFECreateTabulation((PetscFE)obj, 1, numCIndices, xi, 0, &Tcoarse));
420ad9634fcSMatthew G. Knepley       else PetscCall(PetscFVCreateTabulation((PetscFV)obj, 1, numCIndices, xi, 0, &Tcoarse));
42183c47955SMatthew G. Knepley       /* Get elemMat entries by multiplying by weight */
4220bf7c1c5SMatthew G. Knepley       PetscCall(PetscArrayzero(elemMat, numCIndices * Nc * totDim));
4230bf7c1c5SMatthew G. Knepley       for (PetscInt i = 0; i < numFIndices / Nc; ++i) {
4240bf7c1c5SMatthew G. Knepley         for (PetscInt j = 0; j < numCIndices; ++j) {
4250bf7c1c5SMatthew G. Knepley           for (PetscInt c = 0; c < Nc; ++c) {
4260bf7c1c5SMatthew G. Knepley             // TODO Need field offset on particle and field here
4270643ed39SMatthew G. Knepley             /* B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c */
4280bf7c1c5SMatthew G. Knepley             elemMat[(j * totNc + c) * numFIndices + i * Nc + c] += Tcoarse->T[0][(j * numFIndices + i * Nc + c) * Nc + c] * (useDeltaFunction ? 1.0 : detJ);
42983c47955SMatthew G. Knepley           }
4300643ed39SMatthew G. Knepley         }
4310643ed39SMatthew G. Knepley       }
4320bf7c1c5SMatthew G. Knepley       for (PetscInt j = 0; j < numCIndices; ++j)
4330bf7c1c5SMatthew G. Knepley         // TODO Need field offset on particle here
4340bf7c1c5SMatthew G. Knepley         for (PetscInt c = 0; c < Nc; ++c) rowIDXs[j * Nc + c] = cindices[j] * totNc + c + rStart;
4350bf7c1c5SMatthew G. Knepley       if (0) PetscCall(DMPrintCellMatrix(cell, name, numCIndices * Nc, numFIndices, elemMat));
4360bf7c1c5SMatthew G. Knepley       PetscCall(MatSetValues(mass, numCIndices * Nc, rowIDXs, numFIndices, findices, elemMat, ADD_VALUES));
4379566063dSJacob Faibussowitsch       PetscCall(PetscFree(cindices));
4389566063dSJacob Faibussowitsch       PetscCall(DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
4399566063dSJacob Faibussowitsch       PetscCall(PetscTabulationDestroy(&Tcoarse));
44083c47955SMatthew G. Knepley     }
441ea7032a0SMatthew G. Knepley     PetscCall(DMSwarmRestoreField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **)&fieldVals));
44283c47955SMatthew G. Knepley   }
4439566063dSJacob Faibussowitsch   PetscCall(PetscFree3(elemMat, rowIDXs, xi));
4449566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortRestoreAccess(dmc));
4459566063dSJacob Faibussowitsch   PetscCall(PetscFree3(v0, J, invJ));
4469566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mass, MAT_FINAL_ASSEMBLY));
4479566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mass, MAT_FINAL_ASSEMBLY));
4483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
44983c47955SMatthew G. Knepley }
45083c47955SMatthew G. Knepley 
451d0c080abSJoseph Pusztay /* Returns empty matrix for use with SNES FD */
452d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCreateMatrix_Swarm(DM sw, Mat *m)
453d71ae5a4SJacob Faibussowitsch {
454d0c080abSJoseph Pusztay   Vec      field;
455d0c080abSJoseph Pusztay   PetscInt size;
456d0c080abSJoseph Pusztay 
457d0c080abSJoseph Pusztay   PetscFunctionBegin;
4589566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(sw, &field));
4599566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(field, &size));
4609566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(sw, &field));
4619566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, m));
4629566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(*m));
4639566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*m, PETSC_DECIDE, PETSC_DECIDE, size, size));
4649566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(*m, 1, NULL));
4659566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(*m));
4669566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(*m, MAT_FINAL_ASSEMBLY));
4679566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(*m, MAT_FINAL_ASSEMBLY));
4689566063dSJacob Faibussowitsch   PetscCall(MatShift(*m, 1.0));
4699566063dSJacob Faibussowitsch   PetscCall(MatSetDM(*m, sw));
4703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
471d0c080abSJoseph Pusztay }
472d0c080abSJoseph Pusztay 
473adb2528bSMark Adams /* FEM cols, Particle rows */
474d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCreateMassMatrix_Swarm(DM dmCoarse, DM dmFine, Mat *mass)
475d71ae5a4SJacob Faibussowitsch {
476895a1698SMatthew G. Knepley   PetscSection gsf;
4770bf7c1c5SMatthew G. Knepley   PetscInt     m, n, Np, bs = 1;
47883c47955SMatthew G. Knepley   void        *ctx;
4790bf7c1c5SMatthew G. Knepley   PetscBool    set  = ((DM_Swarm *)dmCoarse->data)->vec_field_set;
4800bf7c1c5SMatthew G. Knepley   char        *name = ((DM_Swarm *)dmCoarse->data)->vec_field_name;
48183c47955SMatthew G. Knepley 
48283c47955SMatthew G. Knepley   PetscFunctionBegin;
4839566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalSection(dmFine, &gsf));
4849566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetConstrainedStorageSize(gsf, &m));
4850bf7c1c5SMatthew G. Knepley   PetscCall(DMSwarmGetLocalSize(dmCoarse, &Np));
4860bf7c1c5SMatthew G. Knepley   // TODO Include all fields
4870bf7c1c5SMatthew G. Knepley   if (set) PetscCall(DMSwarmGetFieldInfo(dmCoarse, name, &bs, NULL));
4880bf7c1c5SMatthew G. Knepley   n = Np * bs;
4899566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)dmCoarse), mass));
4909566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*mass, n, m, PETSC_DETERMINE, PETSC_DETERMINE));
4919566063dSJacob Faibussowitsch   PetscCall(MatSetType(*mass, dmCoarse->mattype));
4929566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(dmFine, &ctx));
49383c47955SMatthew G. Knepley 
4949566063dSJacob Faibussowitsch   PetscCall(DMSwarmComputeMassMatrix_Private(dmCoarse, dmFine, *mass, PETSC_TRUE, ctx));
4959566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(*mass, NULL, "-mass_mat_view"));
4963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
49783c47955SMatthew G. Knepley }
49883c47955SMatthew G. Knepley 
499d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSwarmComputeMassMatrixSquare_Private(DM dmc, DM dmf, Mat mass, PetscBool useDeltaFunction, void *ctx)
500d71ae5a4SJacob Faibussowitsch {
5014cc7f7b2SMatthew G. Knepley   const char  *name = "Mass Matrix Square";
5024cc7f7b2SMatthew G. Knepley   MPI_Comm     comm;
5034cc7f7b2SMatthew G. Knepley   PetscDS      prob;
5044cc7f7b2SMatthew G. Knepley   PetscSection fsection, globalFSection;
5054cc7f7b2SMatthew G. Knepley   PetscHSetIJ  ht;
5064cc7f7b2SMatthew G. Knepley   PetscLayout  rLayout, colLayout;
5074cc7f7b2SMatthew G. Knepley   PetscInt    *dnz, *onz, *adj, depth, maxConeSize, maxSupportSize, maxAdjSize;
5084cc7f7b2SMatthew G. Knepley   PetscInt     locRows, locCols, rStart, colStart, colEnd, *rowIDXs;
5094cc7f7b2SMatthew G. Knepley   PetscReal   *xi, *v0, *J, *invJ, detJ = 1.0, v0ref[3] = {-1.0, -1.0, -1.0};
5104cc7f7b2SMatthew G. Knepley   PetscScalar *elemMat, *elemMatSq;
5114cc7f7b2SMatthew G. Knepley   PetscInt     cdim, Nf, field, cStart, cEnd, cell, totDim, maxC = 0;
5124cc7f7b2SMatthew G. Knepley 
5134cc7f7b2SMatthew G. Knepley   PetscFunctionBegin;
5149566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)mass, &comm));
5159566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dmf, &cdim));
5169566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dmf, &prob));
5179566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(prob, &Nf));
5189566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(prob, &totDim));
5199566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(cdim, &v0, cdim * cdim, &J, cdim * cdim, &invJ));
5209566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dmf, &fsection));
5219566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalSection(dmf, &globalFSection));
5229566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd));
5239566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(mass, &locRows, &locCols));
5244cc7f7b2SMatthew G. Knepley 
5259566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreate(comm, &colLayout));
5269566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetLocalSize(colLayout, locCols));
5279566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(colLayout, 1));
5289566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(colLayout));
5299566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetRange(colLayout, &colStart, &colEnd));
5309566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&colLayout));
5314cc7f7b2SMatthew G. Knepley 
5329566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreate(comm, &rLayout));
5339566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetLocalSize(rLayout, locRows));
5349566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(rLayout, 1));
5359566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(rLayout));
5369566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetRange(rLayout, &rStart, NULL));
5379566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&rLayout));
5384cc7f7b2SMatthew G. Knepley 
5399566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(dmf, &depth));
5409566063dSJacob Faibussowitsch   PetscCall(DMPlexGetMaxSizes(dmf, &maxConeSize, &maxSupportSize));
5414cc7f7b2SMatthew G. Knepley   maxAdjSize = PetscPowInt(maxConeSize * maxSupportSize, depth);
5429566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(maxAdjSize, &adj));
5434cc7f7b2SMatthew G. Knepley 
5449566063dSJacob Faibussowitsch   PetscCall(PetscCalloc2(locRows, &dnz, locRows, &onz));
5459566063dSJacob Faibussowitsch   PetscCall(PetscHSetIJCreate(&ht));
5464cc7f7b2SMatthew G. Knepley   /* Count nonzeros
5474cc7f7b2SMatthew G. Knepley        This is just FVM++, but we cannot use the Plex P0 allocation since unknowns in a cell will not be contiguous
5484cc7f7b2SMatthew G. Knepley   */
5499566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortGetAccess(dmc));
5504cc7f7b2SMatthew G. Knepley   for (cell = cStart; cell < cEnd; ++cell) {
5514cc7f7b2SMatthew G. Knepley     PetscInt  i;
5524cc7f7b2SMatthew G. Knepley     PetscInt *cindices;
5534cc7f7b2SMatthew G. Knepley     PetscInt  numCIndices;
5544cc7f7b2SMatthew G. Knepley #if 0
5554cc7f7b2SMatthew G. Knepley     PetscInt  adjSize = maxAdjSize, a, j;
5564cc7f7b2SMatthew G. Knepley #endif
5574cc7f7b2SMatthew G. Knepley 
5589566063dSJacob Faibussowitsch     PetscCall(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices));
5594cc7f7b2SMatthew G. Knepley     maxC = PetscMax(maxC, numCIndices);
5604cc7f7b2SMatthew G. Knepley     /* Diagonal block */
561ad540459SPierre Jolivet     for (i = 0; i < numCIndices; ++i) dnz[cindices[i]] += numCIndices;
5624cc7f7b2SMatthew G. Knepley #if 0
5634cc7f7b2SMatthew G. Knepley     /* Off-diagonal blocks */
5649566063dSJacob Faibussowitsch     PetscCall(DMPlexGetAdjacency(dmf, cell, &adjSize, &adj));
5654cc7f7b2SMatthew G. Knepley     for (a = 0; a < adjSize; ++a) {
5664cc7f7b2SMatthew G. Knepley       if (adj[a] >= cStart && adj[a] < cEnd && adj[a] != cell) {
5674cc7f7b2SMatthew G. Knepley         const PetscInt ncell = adj[a];
5684cc7f7b2SMatthew G. Knepley         PetscInt      *ncindices;
5694cc7f7b2SMatthew G. Knepley         PetscInt       numNCIndices;
5704cc7f7b2SMatthew G. Knepley 
5719566063dSJacob Faibussowitsch         PetscCall(DMSwarmSortGetPointsPerCell(dmc, ncell, &numNCIndices, &ncindices));
5724cc7f7b2SMatthew G. Knepley         {
5734cc7f7b2SMatthew G. Knepley           PetscHashIJKey key;
5744cc7f7b2SMatthew G. Knepley           PetscBool      missing;
5754cc7f7b2SMatthew G. Knepley 
5764cc7f7b2SMatthew G. Knepley           for (i = 0; i < numCIndices; ++i) {
5774cc7f7b2SMatthew G. Knepley             key.i = cindices[i] + rStart; /* global rows (from Swarm) */
5784cc7f7b2SMatthew G. Knepley             if (key.i < 0) continue;
5794cc7f7b2SMatthew G. Knepley             for (j = 0; j < numNCIndices; ++j) {
5804cc7f7b2SMatthew G. Knepley               key.j = ncindices[j] + rStart; /* global column (from Swarm) */
5814cc7f7b2SMatthew G. Knepley               if (key.j < 0) continue;
5829566063dSJacob Faibussowitsch               PetscCall(PetscHSetIJQueryAdd(ht, key, &missing));
5834cc7f7b2SMatthew G. Knepley               if (missing) {
5844cc7f7b2SMatthew G. Knepley                 if ((key.j >= colStart) && (key.j < colEnd)) ++dnz[key.i - rStart];
5854cc7f7b2SMatthew G. Knepley                 else                                         ++onz[key.i - rStart];
5864cc7f7b2SMatthew G. Knepley               }
5874cc7f7b2SMatthew G. Knepley             }
5884cc7f7b2SMatthew G. Knepley           }
5894cc7f7b2SMatthew G. Knepley         }
5909566063dSJacob Faibussowitsch         PetscCall(PetscFree(ncindices));
5914cc7f7b2SMatthew G. Knepley       }
5924cc7f7b2SMatthew G. Knepley     }
5934cc7f7b2SMatthew G. Knepley #endif
5949566063dSJacob Faibussowitsch     PetscCall(PetscFree(cindices));
5954cc7f7b2SMatthew G. Knepley   }
5969566063dSJacob Faibussowitsch   PetscCall(PetscHSetIJDestroy(&ht));
5979566063dSJacob Faibussowitsch   PetscCall(MatXAIJSetPreallocation(mass, 1, dnz, onz, NULL, NULL));
5989566063dSJacob Faibussowitsch   PetscCall(MatSetOption(mass, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
5999566063dSJacob Faibussowitsch   PetscCall(PetscFree2(dnz, onz));
6009566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(maxC * totDim, &elemMat, maxC * maxC, &elemMatSq, maxC, &rowIDXs, maxC * cdim, &xi));
6014cc7f7b2SMatthew G. Knepley   /* Fill in values
6024cc7f7b2SMatthew G. Knepley        Each entry is a sum of terms \phi_i(x_p) \phi_i(x_q)
6034cc7f7b2SMatthew G. Knepley        Start just by producing block diagonal
6044cc7f7b2SMatthew G. Knepley        Could loop over adjacent cells
6054cc7f7b2SMatthew G. Knepley          Produce neighboring element matrix
6064cc7f7b2SMatthew G. Knepley          TODO Determine which columns and rows correspond to shared dual vector
6074cc7f7b2SMatthew G. Knepley          Do MatMatMult with rectangular matrices
6084cc7f7b2SMatthew G. Knepley          Insert block
6094cc7f7b2SMatthew G. Knepley   */
6104cc7f7b2SMatthew G. Knepley   for (field = 0; field < Nf; ++field) {
6114cc7f7b2SMatthew G. Knepley     PetscTabulation Tcoarse;
6124cc7f7b2SMatthew G. Knepley     PetscObject     obj;
6134cc7f7b2SMatthew G. Knepley     PetscReal      *coords;
6144cc7f7b2SMatthew G. Knepley     PetscInt        Nc, i;
6154cc7f7b2SMatthew G. Knepley 
6169566063dSJacob Faibussowitsch     PetscCall(PetscDSGetDiscretization(prob, field, &obj));
6179566063dSJacob Faibussowitsch     PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc));
61863a3b9bcSJacob Faibussowitsch     PetscCheck(Nc == 1, PetscObjectComm((PetscObject)dmf), PETSC_ERR_SUP, "Can only interpolate a scalar field from particles, Nc = %" PetscInt_FMT, Nc);
6199566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
6204cc7f7b2SMatthew G. Knepley     for (cell = cStart; cell < cEnd; ++cell) {
6214cc7f7b2SMatthew G. Knepley       PetscInt *findices, *cindices;
6224cc7f7b2SMatthew G. Knepley       PetscInt  numFIndices, numCIndices;
6234cc7f7b2SMatthew G. Knepley       PetscInt  p, c;
6244cc7f7b2SMatthew G. Knepley 
6254cc7f7b2SMatthew G. Knepley       /* TODO: Use DMField instead of assuming affine */
6269566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ));
6279566063dSJacob Faibussowitsch       PetscCall(DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
6289566063dSJacob Faibussowitsch       PetscCall(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices));
629ad540459SPierre Jolivet       for (p = 0; p < numCIndices; ++p) CoordinatesRealToRef(cdim, cdim, v0ref, v0, invJ, &coords[cindices[p] * cdim], &xi[p * cdim]);
6309566063dSJacob Faibussowitsch       PetscCall(PetscFECreateTabulation((PetscFE)obj, 1, numCIndices, xi, 0, &Tcoarse));
6314cc7f7b2SMatthew G. Knepley       /* Get elemMat entries by multiplying by weight */
6329566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(elemMat, numCIndices * totDim));
6334cc7f7b2SMatthew G. Knepley       for (i = 0; i < numFIndices; ++i) {
6344cc7f7b2SMatthew G. Knepley         for (p = 0; p < numCIndices; ++p) {
6354cc7f7b2SMatthew G. Knepley           for (c = 0; c < Nc; ++c) {
6364cc7f7b2SMatthew G. Knepley             /* B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c */
6374cc7f7b2SMatthew G. Knepley             elemMat[p * numFIndices + i] += Tcoarse->T[0][(p * numFIndices + i) * Nc + c] * (useDeltaFunction ? 1.0 : detJ);
6384cc7f7b2SMatthew G. Knepley           }
6394cc7f7b2SMatthew G. Knepley         }
6404cc7f7b2SMatthew G. Knepley       }
6419566063dSJacob Faibussowitsch       PetscCall(PetscTabulationDestroy(&Tcoarse));
6424cc7f7b2SMatthew G. Knepley       for (p = 0; p < numCIndices; ++p) rowIDXs[p] = cindices[p] + rStart;
6439566063dSJacob Faibussowitsch       if (0) PetscCall(DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat));
6444cc7f7b2SMatthew G. Knepley       /* Block diagonal */
64578884ca7SMatthew G. Knepley       if (numCIndices) {
6464cc7f7b2SMatthew G. Knepley         PetscBLASInt blasn, blask;
6474cc7f7b2SMatthew G. Knepley         PetscScalar  one = 1.0, zero = 0.0;
6484cc7f7b2SMatthew G. Knepley 
6499566063dSJacob Faibussowitsch         PetscCall(PetscBLASIntCast(numCIndices, &blasn));
6509566063dSJacob Faibussowitsch         PetscCall(PetscBLASIntCast(numFIndices, &blask));
651792fecdfSBarry Smith         PetscCallBLAS("BLASgemm", BLASgemm_("T", "N", &blasn, &blasn, &blask, &one, elemMat, &blask, elemMat, &blask, &zero, elemMatSq, &blasn));
6524cc7f7b2SMatthew G. Knepley       }
6539566063dSJacob Faibussowitsch       PetscCall(MatSetValues(mass, numCIndices, rowIDXs, numCIndices, rowIDXs, elemMatSq, ADD_VALUES));
6544cf0e950SBarry Smith       /* TODO off-diagonal */
6559566063dSJacob Faibussowitsch       PetscCall(PetscFree(cindices));
6569566063dSJacob Faibussowitsch       PetscCall(DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
6574cc7f7b2SMatthew G. Knepley     }
6589566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
6594cc7f7b2SMatthew G. Knepley   }
6609566063dSJacob Faibussowitsch   PetscCall(PetscFree4(elemMat, elemMatSq, rowIDXs, xi));
6619566063dSJacob Faibussowitsch   PetscCall(PetscFree(adj));
6629566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortRestoreAccess(dmc));
6639566063dSJacob Faibussowitsch   PetscCall(PetscFree3(v0, J, invJ));
6649566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mass, MAT_FINAL_ASSEMBLY));
6659566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mass, MAT_FINAL_ASSEMBLY));
6663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6674cc7f7b2SMatthew G. Knepley }
6684cc7f7b2SMatthew G. Knepley 
6694cc7f7b2SMatthew G. Knepley /*@C
6704cc7f7b2SMatthew G. Knepley   DMSwarmCreateMassMatrixSquare - Creates the block-diagonal of the square, M^T_p M_p, of the particle mass matrix M_p
6714cc7f7b2SMatthew G. Knepley 
67220f4b53cSBarry Smith   Collective
6734cc7f7b2SMatthew G. Knepley 
67460225df5SJacob Faibussowitsch   Input Parameters:
67520f4b53cSBarry Smith + dmCoarse - a `DMSWARM`
67620f4b53cSBarry Smith - dmFine   - a `DMPLEX`
6774cc7f7b2SMatthew G. Knepley 
67860225df5SJacob Faibussowitsch   Output Parameter:
6794cc7f7b2SMatthew G. Knepley . mass - the square of the particle mass matrix
6804cc7f7b2SMatthew G. Knepley 
6814cc7f7b2SMatthew G. Knepley   Level: advanced
6824cc7f7b2SMatthew G. Knepley 
68320f4b53cSBarry Smith   Note:
6844cc7f7b2SMatthew G. Knepley   We only compute the block diagonal since this provides a good preconditioner and is completely local. It would be possible in the
6854cc7f7b2SMatthew G. Knepley   future to compute the full normal equations.
6864cc7f7b2SMatthew G. Knepley 
68720f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMCreateMassMatrix()`
6884cc7f7b2SMatthew G. Knepley @*/
689d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmCreateMassMatrixSquare(DM dmCoarse, DM dmFine, Mat *mass)
690d71ae5a4SJacob Faibussowitsch {
6914cc7f7b2SMatthew G. Knepley   PetscInt n;
6924cc7f7b2SMatthew G. Knepley   void    *ctx;
6934cc7f7b2SMatthew G. Knepley 
6944cc7f7b2SMatthew G. Knepley   PetscFunctionBegin;
6959566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetLocalSize(dmCoarse, &n));
6969566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)dmCoarse), mass));
6979566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*mass, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
6989566063dSJacob Faibussowitsch   PetscCall(MatSetType(*mass, dmCoarse->mattype));
6999566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(dmFine, &ctx));
7004cc7f7b2SMatthew G. Knepley 
7019566063dSJacob Faibussowitsch   PetscCall(DMSwarmComputeMassMatrixSquare_Private(dmCoarse, dmFine, *mass, PETSC_TRUE, ctx));
7029566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(*mass, NULL, "-mass_sq_mat_view"));
7033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7044cc7f7b2SMatthew G. Knepley }
7054cc7f7b2SMatthew G. Knepley 
706fb1bcc12SMatthew G. Knepley /*@C
70720f4b53cSBarry Smith   DMSwarmCreateGlobalVectorFromField - Creates a `Vec` object sharing the array associated with a given field
708d3a51819SDave May 
70920f4b53cSBarry Smith   Collective
710d3a51819SDave May 
71160225df5SJacob Faibussowitsch   Input Parameters:
71220f4b53cSBarry Smith + dm        - a `DMSWARM`
71362741f57SDave May - fieldname - the textual name given to a registered field
714d3a51819SDave May 
71560225df5SJacob Faibussowitsch   Output Parameter:
716d3a51819SDave May . vec - the vector
717d3a51819SDave May 
718d3a51819SDave May   Level: beginner
719d3a51819SDave May 
7208b8a3813SDave May   Notes:
72120f4b53cSBarry Smith   The vector must be returned using a matching call to `DMSwarmDestroyGlobalVectorFromField()`.
7228b8a3813SDave May 
72320f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmDestroyGlobalVectorFromField()`
724d3a51819SDave May @*/
725d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmCreateGlobalVectorFromField(DM dm, const char fieldname[], Vec *vec)
726d71ae5a4SJacob Faibussowitsch {
727fb1bcc12SMatthew G. Knepley   MPI_Comm comm = PetscObjectComm((PetscObject)dm);
728b5bcf523SDave May 
729fb1bcc12SMatthew G. Knepley   PetscFunctionBegin;
730ea7032a0SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
7319566063dSJacob Faibussowitsch   PetscCall(DMSwarmCreateVectorFromField_Private(dm, fieldname, comm, vec));
7323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
733b5bcf523SDave May }
734b5bcf523SDave May 
735d3a51819SDave May /*@C
73620f4b53cSBarry Smith   DMSwarmDestroyGlobalVectorFromField - Destroys the `Vec` object which share the array associated with a given field
737d3a51819SDave May 
73820f4b53cSBarry Smith   Collective
739d3a51819SDave May 
74060225df5SJacob Faibussowitsch   Input Parameters:
74120f4b53cSBarry Smith + dm        - a `DMSWARM`
74262741f57SDave May - fieldname - the textual name given to a registered field
743d3a51819SDave May 
74460225df5SJacob Faibussowitsch   Output Parameter:
745d3a51819SDave May . vec - the vector
746d3a51819SDave May 
747d3a51819SDave May   Level: beginner
748d3a51819SDave May 
74920f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmCreateGlobalVectorFromField()`
750d3a51819SDave May @*/
751d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDestroyGlobalVectorFromField(DM dm, const char fieldname[], Vec *vec)
752d71ae5a4SJacob Faibussowitsch {
753fb1bcc12SMatthew G. Knepley   PetscFunctionBegin;
754ea7032a0SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
7559566063dSJacob Faibussowitsch   PetscCall(DMSwarmDestroyVectorFromField_Private(dm, fieldname, vec));
7563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
757b5bcf523SDave May }
758b5bcf523SDave May 
759fb1bcc12SMatthew G. Knepley /*@C
76020f4b53cSBarry Smith   DMSwarmCreateLocalVectorFromField - Creates a `Vec` object sharing the array associated with a given field
761fb1bcc12SMatthew G. Knepley 
76220f4b53cSBarry Smith   Collective
763fb1bcc12SMatthew G. Knepley 
76460225df5SJacob Faibussowitsch   Input Parameters:
76520f4b53cSBarry Smith + dm        - a `DMSWARM`
76662741f57SDave May - fieldname - the textual name given to a registered field
767fb1bcc12SMatthew G. Knepley 
76860225df5SJacob Faibussowitsch   Output Parameter:
769fb1bcc12SMatthew G. Knepley . vec - the vector
770fb1bcc12SMatthew G. Knepley 
771fb1bcc12SMatthew G. Knepley   Level: beginner
772fb1bcc12SMatthew G. Knepley 
77320f4b53cSBarry Smith   Note:
7748b8a3813SDave May   The vector must be returned using a matching call to DMSwarmDestroyLocalVectorFromField().
7758b8a3813SDave May 
77620f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmDestroyLocalVectorFromField()`
777fb1bcc12SMatthew G. Knepley @*/
778d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmCreateLocalVectorFromField(DM dm, const char fieldname[], Vec *vec)
779d71ae5a4SJacob Faibussowitsch {
780fb1bcc12SMatthew G. Knepley   MPI_Comm comm = PETSC_COMM_SELF;
781bbe8250bSMatthew G. Knepley 
782fb1bcc12SMatthew G. Knepley   PetscFunctionBegin;
7839566063dSJacob Faibussowitsch   PetscCall(DMSwarmCreateVectorFromField_Private(dm, fieldname, comm, vec));
7843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
785bbe8250bSMatthew G. Knepley }
786fb1bcc12SMatthew G. Knepley 
787fb1bcc12SMatthew G. Knepley /*@C
78820f4b53cSBarry Smith   DMSwarmDestroyLocalVectorFromField - Destroys the `Vec` object which share the array associated with a given field
789fb1bcc12SMatthew G. Knepley 
79020f4b53cSBarry Smith   Collective
791fb1bcc12SMatthew G. Knepley 
79260225df5SJacob Faibussowitsch   Input Parameters:
79320f4b53cSBarry Smith + dm        - a `DMSWARM`
79462741f57SDave May - fieldname - the textual name given to a registered field
795fb1bcc12SMatthew G. Knepley 
79660225df5SJacob Faibussowitsch   Output Parameter:
797fb1bcc12SMatthew G. Knepley . vec - the vector
798fb1bcc12SMatthew G. Knepley 
799fb1bcc12SMatthew G. Knepley   Level: beginner
800fb1bcc12SMatthew G. Knepley 
80120f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmCreateLocalVectorFromField()`
802fb1bcc12SMatthew G. Knepley @*/
803d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDestroyLocalVectorFromField(DM dm, const char fieldname[], Vec *vec)
804d71ae5a4SJacob Faibussowitsch {
805fb1bcc12SMatthew G. Knepley   PetscFunctionBegin;
806ea7032a0SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
8079566063dSJacob Faibussowitsch   PetscCall(DMSwarmDestroyVectorFromField_Private(dm, fieldname, vec));
8083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
809bbe8250bSMatthew G. Knepley }
810bbe8250bSMatthew G. Knepley 
811d3a51819SDave May /*@C
81220f4b53cSBarry Smith   DMSwarmInitializeFieldRegister - Initiates the registration of fields to a `DMSWARM`
813d3a51819SDave May 
81420f4b53cSBarry Smith   Collective
815d3a51819SDave May 
81660225df5SJacob Faibussowitsch   Input Parameter:
81720f4b53cSBarry Smith . dm - a `DMSWARM`
818d3a51819SDave May 
819d3a51819SDave May   Level: beginner
820d3a51819SDave May 
82120f4b53cSBarry Smith   Note:
82220f4b53cSBarry Smith   After all fields have been registered, you must call `DMSwarmFinalizeFieldRegister()`.
823d3a51819SDave May 
82420f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmFinalizeFieldRegister()`, `DMSwarmRegisterPetscDatatypeField()`,
825db781477SPatrick Sanan           `DMSwarmRegisterUserStructField()`, `DMSwarmRegisterUserDatatypeField()`
826d3a51819SDave May @*/
827d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmInitializeFieldRegister(DM dm)
828d71ae5a4SJacob Faibussowitsch {
8295f50eb2eSDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
8303454631fSDave May 
831521f74f9SMatthew G. Knepley   PetscFunctionBegin;
832cc651181SDave May   if (!swarm->field_registration_initialized) {
8335f50eb2eSDave May     swarm->field_registration_initialized = PETSC_TRUE;
834da81f932SPierre Jolivet     PetscCall(DMSwarmRegisterPetscDatatypeField(dm, DMSwarmField_pid, 1, PETSC_INT64)); /* unique identifier */
8359566063dSJacob Faibussowitsch     PetscCall(DMSwarmRegisterPetscDatatypeField(dm, DMSwarmField_rank, 1, PETSC_INT));  /* used for communication */
836cc651181SDave May   }
8373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8385f50eb2eSDave May }
8395f50eb2eSDave May 
84074d0cae8SMatthew G. Knepley /*@
84120f4b53cSBarry Smith   DMSwarmFinalizeFieldRegister - Finalizes the registration of fields to a `DMSWARM`
842d3a51819SDave May 
84320f4b53cSBarry Smith   Collective
844d3a51819SDave May 
84560225df5SJacob Faibussowitsch   Input Parameter:
84620f4b53cSBarry Smith . dm - a `DMSWARM`
847d3a51819SDave May 
848d3a51819SDave May   Level: beginner
849d3a51819SDave May 
85020f4b53cSBarry Smith   Note:
85120f4b53cSBarry Smith   After `DMSwarmFinalizeFieldRegister()` has been called, no new fields can be defined on the `DMSWARM`.
852d3a51819SDave May 
85320f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmInitializeFieldRegister()`, `DMSwarmRegisterPetscDatatypeField()`,
854db781477SPatrick Sanan           `DMSwarmRegisterUserStructField()`, `DMSwarmRegisterUserDatatypeField()`
855d3a51819SDave May @*/
856d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmFinalizeFieldRegister(DM dm)
857d71ae5a4SJacob Faibussowitsch {
8585f50eb2eSDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
8596845f8f5SDave May 
860521f74f9SMatthew G. Knepley   PetscFunctionBegin;
86148a46eb9SPierre Jolivet   if (!swarm->field_registration_finalized) PetscCall(DMSwarmDataBucketFinalize(swarm->db));
862f0cdbbbaSDave May   swarm->field_registration_finalized = PETSC_TRUE;
8633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8645f50eb2eSDave May }
8655f50eb2eSDave May 
86674d0cae8SMatthew G. Knepley /*@
86720f4b53cSBarry Smith   DMSwarmSetLocalSizes - Sets the length of all registered fields on the `DMSWARM`
868d3a51819SDave May 
86920f4b53cSBarry Smith   Not Collective
870d3a51819SDave May 
87160225df5SJacob Faibussowitsch   Input Parameters:
87220f4b53cSBarry Smith + dm     - a `DMSWARM`
873d3a51819SDave May . nlocal - the length of each registered field
87462741f57SDave May - buffer - the length of the buffer used to efficient dynamic re-sizing
875d3a51819SDave May 
876d3a51819SDave May   Level: beginner
877d3a51819SDave May 
87820f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmGetLocalSize()`
879d3a51819SDave May @*/
880d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmSetLocalSizes(DM dm, PetscInt nlocal, PetscInt buffer)
881d71ae5a4SJacob Faibussowitsch {
8825f50eb2eSDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
8835f50eb2eSDave May 
884521f74f9SMatthew G. Knepley   PetscFunctionBegin;
8859566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMSWARM_SetSizes, 0, 0, 0, 0));
8869566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketSetSizes(swarm->db, nlocal, buffer));
8879566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMSWARM_SetSizes, 0, 0, 0, 0));
8883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8895f50eb2eSDave May }
8905f50eb2eSDave May 
89174d0cae8SMatthew G. Knepley /*@
89220f4b53cSBarry Smith   DMSwarmSetCellDM - Attaches a `DM` to a `DMSWARM`
893d3a51819SDave May 
89420f4b53cSBarry Smith   Collective
895d3a51819SDave May 
89660225df5SJacob Faibussowitsch   Input Parameters:
89720f4b53cSBarry Smith + dm     - a `DMSWARM`
89820f4b53cSBarry Smith - dmcell - the `DM` to attach to the `DMSWARM`
899d3a51819SDave May 
900d3a51819SDave May   Level: beginner
901d3a51819SDave May 
90220f4b53cSBarry Smith   Note:
90320f4b53cSBarry Smith   The attached `DM` (dmcell) will be queried for point location and
90420f4b53cSBarry Smith   neighbor MPI-rank information if `DMSwarmMigrate()` is called.
905d3a51819SDave May 
90620f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmSetType()`, `DMSwarmGetCellDM()`, `DMSwarmMigrate()`
907d3a51819SDave May @*/
908d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmSetCellDM(DM dm, DM dmcell)
909d71ae5a4SJacob Faibussowitsch {
910b16650c8SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
911521f74f9SMatthew G. Knepley 
912521f74f9SMatthew G. Knepley   PetscFunctionBegin;
913b16650c8SDave May   swarm->dmcell = dmcell;
9143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
915b16650c8SDave May }
916b16650c8SDave May 
91774d0cae8SMatthew G. Knepley /*@
91820f4b53cSBarry Smith   DMSwarmGetCellDM - Fetches the attached cell `DM`
919d3a51819SDave May 
92020f4b53cSBarry Smith   Collective
921d3a51819SDave May 
92260225df5SJacob Faibussowitsch   Input Parameter:
92320f4b53cSBarry Smith . dm - a `DMSWARM`
924d3a51819SDave May 
92560225df5SJacob Faibussowitsch   Output Parameter:
92620f4b53cSBarry Smith . dmcell - the `DM` which was attached to the `DMSWARM`
927d3a51819SDave May 
928d3a51819SDave May   Level: beginner
929d3a51819SDave May 
93020f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmSetCellDM()`
931d3a51819SDave May @*/
932d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmGetCellDM(DM dm, DM *dmcell)
933d71ae5a4SJacob Faibussowitsch {
934fe39f135SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
935521f74f9SMatthew G. Knepley 
936521f74f9SMatthew G. Knepley   PetscFunctionBegin;
937fe39f135SDave May   *dmcell = swarm->dmcell;
9383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
939fe39f135SDave May }
940fe39f135SDave May 
94174d0cae8SMatthew G. Knepley /*@
942d5b43468SJose E. Roman   DMSwarmGetLocalSize - Retrieves the local length of fields registered
943d3a51819SDave May 
94420f4b53cSBarry Smith   Not Collective
945d3a51819SDave May 
94660225df5SJacob Faibussowitsch   Input Parameter:
94720f4b53cSBarry Smith . dm - a `DMSWARM`
948d3a51819SDave May 
94960225df5SJacob Faibussowitsch   Output Parameter:
950d3a51819SDave May . nlocal - the length of each registered field
951d3a51819SDave May 
952d3a51819SDave May   Level: beginner
953d3a51819SDave May 
95420f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmGetSize()`, `DMSwarmSetLocalSizes()`
955d3a51819SDave May @*/
956d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmGetLocalSize(DM dm, PetscInt *nlocal)
957d71ae5a4SJacob Faibussowitsch {
958dcf43ee8SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
959dcf43ee8SDave May 
960521f74f9SMatthew G. Knepley   PetscFunctionBegin;
9619566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, nlocal, NULL, NULL));
9623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
963dcf43ee8SDave May }
964dcf43ee8SDave May 
96574d0cae8SMatthew G. Knepley /*@
966d5b43468SJose E. Roman   DMSwarmGetSize - Retrieves the total length of fields registered
967d3a51819SDave May 
96820f4b53cSBarry Smith   Collective
969d3a51819SDave May 
97060225df5SJacob Faibussowitsch   Input Parameter:
97120f4b53cSBarry Smith . dm - a `DMSWARM`
972d3a51819SDave May 
97360225df5SJacob Faibussowitsch   Output Parameter:
974d3a51819SDave May . n - the total length of each registered field
975d3a51819SDave May 
976d3a51819SDave May   Level: beginner
977d3a51819SDave May 
978d3a51819SDave May   Note:
97920f4b53cSBarry Smith   This calls `MPI_Allreduce()` upon each call (inefficient but safe)
980d3a51819SDave May 
98120f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmGetLocalSize()`, `DMSwarmSetLocalSizes()`
982d3a51819SDave May @*/
983d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmGetSize(DM dm, PetscInt *n)
984d71ae5a4SJacob Faibussowitsch {
985dcf43ee8SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
9865627991aSBarry Smith   PetscInt  nlocal;
987dcf43ee8SDave May 
988521f74f9SMatthew G. Knepley   PetscFunctionBegin;
9899566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &nlocal, NULL, NULL));
990712fec58SPierre Jolivet   PetscCall(MPIU_Allreduce(&nlocal, n, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm)));
9913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
992dcf43ee8SDave May }
993dcf43ee8SDave May 
994d3a51819SDave May /*@C
99520f4b53cSBarry Smith   DMSwarmRegisterPetscDatatypeField - Register a field to a `DMSWARM` with a native PETSc data type
996d3a51819SDave May 
99720f4b53cSBarry Smith   Collective
998d3a51819SDave May 
99960225df5SJacob Faibussowitsch   Input Parameters:
100020f4b53cSBarry Smith + dm        - a `DMSWARM`
1001d3a51819SDave May . fieldname - the textual name to identify this field
1002d3a51819SDave May . blocksize - the number of each data type
100320f4b53cSBarry Smith - type      - a valid PETSc data type (`PETSC_CHAR`, `PETSC_SHORT`, `PETSC_INT`, `PETSC_FLOAT`, `PETSC_REAL`, `PETSC_LONG`)
1004d3a51819SDave May 
1005d3a51819SDave May   Level: beginner
1006d3a51819SDave May 
1007d3a51819SDave May   Notes:
10088b8a3813SDave May   The textual name for each registered field must be unique.
1009d3a51819SDave May 
101020f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterUserStructField()`, `DMSwarmRegisterUserDatatypeField()`
1011d3a51819SDave May @*/
1012d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmRegisterPetscDatatypeField(DM dm, const char fieldname[], PetscInt blocksize, PetscDataType type)
1013d71ae5a4SJacob Faibussowitsch {
1014b62e03f8SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1015b62e03f8SDave May   size_t    size;
1016b62e03f8SDave May 
1017521f74f9SMatthew G. Knepley   PetscFunctionBegin;
101828b400f6SJacob Faibussowitsch   PetscCheck(swarm->field_registration_initialized, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Must call DMSwarmInitializeFieldRegister() first");
101928b400f6SJacob Faibussowitsch   PetscCheck(!swarm->field_registration_finalized, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Cannot register additional fields after calling DMSwarmFinalizeFieldRegister() first");
10205f50eb2eSDave May 
102108401ef6SPierre Jolivet   PetscCheck(type != PETSC_OBJECT, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}");
102208401ef6SPierre Jolivet   PetscCheck(type != PETSC_FUNCTION, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}");
102308401ef6SPierre Jolivet   PetscCheck(type != PETSC_STRING, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}");
102408401ef6SPierre Jolivet   PetscCheck(type != PETSC_STRUCT, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}");
102508401ef6SPierre Jolivet   PetscCheck(type != PETSC_DATATYPE_UNKNOWN, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}");
1026b62e03f8SDave May 
10279566063dSJacob Faibussowitsch   PetscCall(PetscDataTypeGetSize(type, &size));
1028b62e03f8SDave May   /* Load a specific data type into data bucket, specifying textual name and its size in bytes */
10299566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketRegisterField(swarm->db, "DMSwarmRegisterPetscDatatypeField", fieldname, blocksize * size, NULL));
103052c3ed93SDave May   {
103177048351SPatrick Sanan     DMSwarmDataField gfield;
103252c3ed93SDave May 
10339566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield));
10349566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataFieldSetBlockSize(gfield, blocksize));
103552c3ed93SDave May   }
1036b62e03f8SDave May   swarm->db->field[swarm->db->nfields - 1]->petsc_type = type;
10373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1038b62e03f8SDave May }
1039b62e03f8SDave May 
1040d3a51819SDave May /*@C
104120f4b53cSBarry Smith   DMSwarmRegisterUserStructField - Register a user defined struct to a `DMSWARM`
1042d3a51819SDave May 
104320f4b53cSBarry Smith   Collective
1044d3a51819SDave May 
104560225df5SJacob Faibussowitsch   Input Parameters:
104620f4b53cSBarry Smith + dm        - a `DMSWARM`
1047d3a51819SDave May . fieldname - the textual name to identify this field
104862741f57SDave May - size      - the size in bytes of the user struct of each data type
1049d3a51819SDave May 
1050d3a51819SDave May   Level: beginner
1051d3a51819SDave May 
105220f4b53cSBarry Smith   Note:
10538b8a3813SDave May   The textual name for each registered field must be unique.
1054d3a51819SDave May 
105520f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmRegisterUserDatatypeField()`
1056d3a51819SDave May @*/
1057d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmRegisterUserStructField(DM dm, const char fieldname[], size_t size)
1058d71ae5a4SJacob Faibussowitsch {
1059b62e03f8SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1060b62e03f8SDave May 
1061521f74f9SMatthew G. Knepley   PetscFunctionBegin;
10629566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketRegisterField(swarm->db, "DMSwarmRegisterUserStructField", fieldname, size, NULL));
1063b62e03f8SDave May   swarm->db->field[swarm->db->nfields - 1]->petsc_type = PETSC_STRUCT;
10643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1065b62e03f8SDave May }
1066b62e03f8SDave May 
1067d3a51819SDave May /*@C
106820f4b53cSBarry Smith   DMSwarmRegisterUserDatatypeField - Register a user defined data type to a `DMSWARM`
1069d3a51819SDave May 
107020f4b53cSBarry Smith   Collective
1071d3a51819SDave May 
107260225df5SJacob Faibussowitsch   Input Parameters:
107320f4b53cSBarry Smith + dm        - a `DMSWARM`
1074d3a51819SDave May . fieldname - the textual name to identify this field
1075d3a51819SDave May . size      - the size in bytes of the user data type
107662741f57SDave May - blocksize - the number of each data type
1077d3a51819SDave May 
1078d3a51819SDave May   Level: beginner
1079d3a51819SDave May 
108020f4b53cSBarry Smith   Note:
10818b8a3813SDave May   The textual name for each registered field must be unique.
1082d3a51819SDave May 
108342747ad1SJacob Faibussowitsch .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmRegisterUserStructField()`
1084d3a51819SDave May @*/
1085d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmRegisterUserDatatypeField(DM dm, const char fieldname[], size_t size, PetscInt blocksize)
1086d71ae5a4SJacob Faibussowitsch {
1087b62e03f8SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1088b62e03f8SDave May 
1089521f74f9SMatthew G. Knepley   PetscFunctionBegin;
10909566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketRegisterField(swarm->db, "DMSwarmRegisterUserDatatypeField", fieldname, blocksize * size, NULL));
1091320740a0SDave May   {
109277048351SPatrick Sanan     DMSwarmDataField gfield;
1093320740a0SDave May 
10949566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield));
10959566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataFieldSetBlockSize(gfield, blocksize));
1096320740a0SDave May   }
1097b62e03f8SDave May   swarm->db->field[swarm->db->nfields - 1]->petsc_type = PETSC_DATATYPE_UNKNOWN;
10983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1099b62e03f8SDave May }
1100b62e03f8SDave May 
1101d3a51819SDave May /*@C
1102d3a51819SDave May   DMSwarmGetField - Get access to the underlying array storing all entries associated with a registered field
1103d3a51819SDave May 
110420f4b53cSBarry Smith   Not Collective
1105d3a51819SDave May 
110660225df5SJacob Faibussowitsch   Input Parameters:
110720f4b53cSBarry Smith + dm        - a `DMSWARM`
110862741f57SDave May - fieldname - the textual name to identify this field
1109d3a51819SDave May 
111060225df5SJacob Faibussowitsch   Output Parameters:
111162741f57SDave May + blocksize - the number of each data type
1112d3a51819SDave May . type      - the data type
111362741f57SDave May - data      - pointer to raw array
1114d3a51819SDave May 
1115d3a51819SDave May   Level: beginner
1116d3a51819SDave May 
1117d3a51819SDave May   Notes:
111820f4b53cSBarry Smith   The array must be returned using a matching call to `DMSwarmRestoreField()`.
1119d3a51819SDave May 
112020f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRestoreField()`
1121d3a51819SDave May @*/
1122d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmGetField(DM dm, const char fieldname[], PetscInt *blocksize, PetscDataType *type, void **data)
1123d71ae5a4SJacob Faibussowitsch {
1124b62e03f8SDave May   DM_Swarm        *swarm = (DM_Swarm *)dm->data;
112577048351SPatrick Sanan   DMSwarmDataField gfield;
1126b62e03f8SDave May 
1127521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1128ea7032a0SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
11299566063dSJacob Faibussowitsch   if (!swarm->issetup) PetscCall(DMSetUp(dm));
11309566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield));
11319566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataFieldGetAccess(gfield));
11329566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataFieldGetEntries(gfield, data));
1133ad540459SPierre Jolivet   if (blocksize) *blocksize = gfield->bs;
1134ad540459SPierre Jolivet   if (type) *type = gfield->petsc_type;
11353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1136b62e03f8SDave May }
1137b62e03f8SDave May 
1138d3a51819SDave May /*@C
1139d3a51819SDave May   DMSwarmRestoreField - Restore access to the underlying array storing all entries associated with a registered field
1140d3a51819SDave May 
114120f4b53cSBarry Smith   Not Collective
1142d3a51819SDave May 
114360225df5SJacob Faibussowitsch   Input Parameters:
114420f4b53cSBarry Smith + dm        - a `DMSWARM`
114562741f57SDave May - fieldname - the textual name to identify this field
1146d3a51819SDave May 
114760225df5SJacob Faibussowitsch   Output Parameters:
114862741f57SDave May + blocksize - the number of each data type
1149d3a51819SDave May . type      - the data type
115062741f57SDave May - data      - pointer to raw array
1151d3a51819SDave May 
1152d3a51819SDave May   Level: beginner
1153d3a51819SDave May 
1154d3a51819SDave May   Notes:
115520f4b53cSBarry Smith   The user must call `DMSwarmGetField()` prior to calling `DMSwarmRestoreField()`.
1156d3a51819SDave May 
115720f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmGetField()`
1158d3a51819SDave May @*/
1159d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmRestoreField(DM dm, const char fieldname[], PetscInt *blocksize, PetscDataType *type, void **data)
1160d71ae5a4SJacob Faibussowitsch {
1161b62e03f8SDave May   DM_Swarm        *swarm = (DM_Swarm *)dm->data;
116277048351SPatrick Sanan   DMSwarmDataField gfield;
1163b62e03f8SDave May 
1164521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1165ea7032a0SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
11669566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield));
11679566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataFieldRestoreAccess(gfield));
1168b62e03f8SDave May   if (data) *data = NULL;
11693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1170b62e03f8SDave May }
1171b62e03f8SDave May 
11720bf7c1c5SMatthew G. Knepley PetscErrorCode DMSwarmGetFieldInfo(DM dm, const char fieldname[], PetscInt *blocksize, PetscDataType *type)
11730bf7c1c5SMatthew G. Knepley {
11740bf7c1c5SMatthew G. Knepley   DM_Swarm        *swarm = (DM_Swarm *)dm->data;
11750bf7c1c5SMatthew G. Knepley   DMSwarmDataField gfield;
11760bf7c1c5SMatthew G. Knepley 
11770bf7c1c5SMatthew G. Knepley   PetscFunctionBegin;
11780bf7c1c5SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
11790bf7c1c5SMatthew G. Knepley   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield));
11800bf7c1c5SMatthew G. Knepley   if (blocksize) *blocksize = gfield->bs;
11810bf7c1c5SMatthew G. Knepley   if (type) *type = gfield->petsc_type;
11820bf7c1c5SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
11830bf7c1c5SMatthew G. Knepley }
11840bf7c1c5SMatthew G. Knepley 
118574d0cae8SMatthew G. Knepley /*@
118620f4b53cSBarry Smith   DMSwarmAddPoint - Add space for one new point in the `DMSWARM`
1187d3a51819SDave May 
118820f4b53cSBarry Smith   Not Collective
1189d3a51819SDave May 
119060225df5SJacob Faibussowitsch   Input Parameter:
119120f4b53cSBarry Smith . dm - a `DMSWARM`
1192d3a51819SDave May 
1193d3a51819SDave May   Level: beginner
1194d3a51819SDave May 
1195d3a51819SDave May   Notes:
11968b8a3813SDave May   The new point will have all fields initialized to zero.
1197d3a51819SDave May 
119820f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmAddNPoints()`
1199d3a51819SDave May @*/
1200d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmAddPoint(DM dm)
1201d71ae5a4SJacob Faibussowitsch {
1202cb1d1399SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1203cb1d1399SDave May 
1204521f74f9SMatthew G. Knepley   PetscFunctionBegin;
12059566063dSJacob Faibussowitsch   if (!swarm->issetup) PetscCall(DMSetUp(dm));
12069566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMSWARM_AddPoints, 0, 0, 0, 0));
12079566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketAddPoint(swarm->db));
12089566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMSWARM_AddPoints, 0, 0, 0, 0));
12093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1210cb1d1399SDave May }
1211cb1d1399SDave May 
121274d0cae8SMatthew G. Knepley /*@
121320f4b53cSBarry Smith   DMSwarmAddNPoints - Add space for a number of new points in the `DMSWARM`
1214d3a51819SDave May 
121520f4b53cSBarry Smith   Not Collective
1216d3a51819SDave May 
121760225df5SJacob Faibussowitsch   Input Parameters:
121820f4b53cSBarry Smith + dm      - a `DMSWARM`
121962741f57SDave May - npoints - the number of new points to add
1220d3a51819SDave May 
1221d3a51819SDave May   Level: beginner
1222d3a51819SDave May 
1223d3a51819SDave May   Notes:
12248b8a3813SDave May   The new point will have all fields initialized to zero.
1225d3a51819SDave May 
122620f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmAddPoint()`
1227d3a51819SDave May @*/
1228d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmAddNPoints(DM dm, PetscInt npoints)
1229d71ae5a4SJacob Faibussowitsch {
1230cb1d1399SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1231cb1d1399SDave May   PetscInt  nlocal;
1232cb1d1399SDave May 
1233521f74f9SMatthew G. Knepley   PetscFunctionBegin;
12349566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMSWARM_AddPoints, 0, 0, 0, 0));
12359566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &nlocal, NULL, NULL));
1236cb1d1399SDave May   nlocal = nlocal + npoints;
12379566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketSetSizes(swarm->db, nlocal, DMSWARM_DATA_BUCKET_BUFFER_DEFAULT));
12389566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMSWARM_AddPoints, 0, 0, 0, 0));
12393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1240cb1d1399SDave May }
1241cb1d1399SDave May 
124274d0cae8SMatthew G. Knepley /*@
124320f4b53cSBarry Smith   DMSwarmRemovePoint - Remove the last point from the `DMSWARM`
1244d3a51819SDave May 
124520f4b53cSBarry Smith   Not Collective
1246d3a51819SDave May 
124760225df5SJacob Faibussowitsch   Input Parameter:
124820f4b53cSBarry Smith . dm - a `DMSWARM`
1249d3a51819SDave May 
1250d3a51819SDave May   Level: beginner
1251d3a51819SDave May 
125220f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRemovePointAtIndex()`
1253d3a51819SDave May @*/
1254d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmRemovePoint(DM dm)
1255d71ae5a4SJacob Faibussowitsch {
1256cb1d1399SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1257cb1d1399SDave May 
1258521f74f9SMatthew G. Knepley   PetscFunctionBegin;
12599566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMSWARM_RemovePoints, 0, 0, 0, 0));
12609566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketRemovePoint(swarm->db));
12619566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMSWARM_RemovePoints, 0, 0, 0, 0));
12623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1263cb1d1399SDave May }
1264cb1d1399SDave May 
126574d0cae8SMatthew G. Knepley /*@
126620f4b53cSBarry Smith   DMSwarmRemovePointAtIndex - Removes a specific point from the `DMSWARM`
1267d3a51819SDave May 
126820f4b53cSBarry Smith   Not Collective
1269d3a51819SDave May 
127060225df5SJacob Faibussowitsch   Input Parameters:
127120f4b53cSBarry Smith + dm  - a `DMSWARM`
127262741f57SDave May - idx - index of point to remove
1273d3a51819SDave May 
1274d3a51819SDave May   Level: beginner
1275d3a51819SDave May 
127620f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRemovePoint()`
1277d3a51819SDave May @*/
1278d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmRemovePointAtIndex(DM dm, PetscInt idx)
1279d71ae5a4SJacob Faibussowitsch {
1280cb1d1399SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1281cb1d1399SDave May 
1282521f74f9SMatthew G. Knepley   PetscFunctionBegin;
12839566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMSWARM_RemovePoints, 0, 0, 0, 0));
12849566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketRemovePointAtIndex(swarm->db, idx));
12859566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMSWARM_RemovePoints, 0, 0, 0, 0));
12863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1287cb1d1399SDave May }
1288b62e03f8SDave May 
128974d0cae8SMatthew G. Knepley /*@
129020f4b53cSBarry Smith   DMSwarmCopyPoint - Copy point pj to point pi in the `DMSWARM`
1291ba4fc9c6SDave May 
129220f4b53cSBarry Smith   Not Collective
1293ba4fc9c6SDave May 
129460225df5SJacob Faibussowitsch   Input Parameters:
129520f4b53cSBarry Smith + dm - a `DMSWARM`
1296ba4fc9c6SDave May . pi - the index of the point to copy
1297ba4fc9c6SDave May - pj - the point index where the copy should be located
1298ba4fc9c6SDave May 
1299ba4fc9c6SDave May   Level: beginner
1300ba4fc9c6SDave May 
130120f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRemovePoint()`
1302ba4fc9c6SDave May @*/
1303d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmCopyPoint(DM dm, PetscInt pi, PetscInt pj)
1304d71ae5a4SJacob Faibussowitsch {
1305ba4fc9c6SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1306ba4fc9c6SDave May 
1307ba4fc9c6SDave May   PetscFunctionBegin;
13089566063dSJacob Faibussowitsch   if (!swarm->issetup) PetscCall(DMSetUp(dm));
13099566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketCopyPoint(swarm->db, pi, swarm->db, pj));
13103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1311ba4fc9c6SDave May }
1312ba4fc9c6SDave May 
131366976f2fSJacob Faibussowitsch static PetscErrorCode DMSwarmMigrate_Basic(DM dm, PetscBool remove_sent_points)
1314d71ae5a4SJacob Faibussowitsch {
1315521f74f9SMatthew G. Knepley   PetscFunctionBegin;
13169566063dSJacob Faibussowitsch   PetscCall(DMSwarmMigrate_Push_Basic(dm, remove_sent_points));
13173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13183454631fSDave May }
13193454631fSDave May 
132074d0cae8SMatthew G. Knepley /*@
132120f4b53cSBarry Smith   DMSwarmMigrate - Relocates points defined in the `DMSWARM` to other MPI-ranks
1322d3a51819SDave May 
132320f4b53cSBarry Smith   Collective
1324d3a51819SDave May 
132560225df5SJacob Faibussowitsch   Input Parameters:
132620f4b53cSBarry Smith + dm                 - the `DMSWARM`
132762741f57SDave May - remove_sent_points - flag indicating if sent points should be removed from the current MPI-rank
1328d3a51819SDave May 
1329d3a51819SDave May   Level: advanced
1330d3a51819SDave May 
133120f4b53cSBarry Smith   Notes:
133220f4b53cSBarry Smith   The `DM` will be modified to accommodate received points.
133320f4b53cSBarry Smith   If `remove_sent_points` is `PETSC_TRUE`, any points that were sent will be removed from the `DM`.
133420f4b53cSBarry Smith   Different styles of migration are supported. See `DMSwarmSetMigrateType()`.
133520f4b53cSBarry Smith 
133620f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmSetMigrateType()`
1337d3a51819SDave May @*/
1338d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmMigrate(DM dm, PetscBool remove_sent_points)
1339d71ae5a4SJacob Faibussowitsch {
1340f0cdbbbaSDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1341f0cdbbbaSDave May 
1342521f74f9SMatthew G. Knepley   PetscFunctionBegin;
13439566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMSWARM_Migrate, 0, 0, 0, 0));
1344f0cdbbbaSDave May   switch (swarm->migrate_type) {
1345d71ae5a4SJacob Faibussowitsch   case DMSWARM_MIGRATE_BASIC:
1346d71ae5a4SJacob Faibussowitsch     PetscCall(DMSwarmMigrate_Basic(dm, remove_sent_points));
1347d71ae5a4SJacob Faibussowitsch     break;
1348d71ae5a4SJacob Faibussowitsch   case DMSWARM_MIGRATE_DMCELLNSCATTER:
1349d71ae5a4SJacob Faibussowitsch     PetscCall(DMSwarmMigrate_CellDMScatter(dm, remove_sent_points));
1350d71ae5a4SJacob Faibussowitsch     break;
1351d71ae5a4SJacob Faibussowitsch   case DMSWARM_MIGRATE_DMCELLEXACT:
1352d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_MIGRATE_DMCELLEXACT not implemented");
1353d71ae5a4SJacob Faibussowitsch   case DMSWARM_MIGRATE_USER:
1354d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_MIGRATE_USER not implemented");
1355d71ae5a4SJacob Faibussowitsch   default:
1356d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_MIGRATE type unknown");
1357f0cdbbbaSDave May   }
13589566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMSWARM_Migrate, 0, 0, 0, 0));
13599566063dSJacob Faibussowitsch   PetscCall(DMClearGlobalVectors(dm));
13603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13613454631fSDave May }
13623454631fSDave May 
1363f0cdbbbaSDave May PetscErrorCode DMSwarmMigrate_GlobalToLocal_Basic(DM dm, PetscInt *globalsize);
1364f0cdbbbaSDave May 
1365d3a51819SDave May /*
1366d3a51819SDave May  DMSwarmCollectViewCreate
1367d3a51819SDave May 
1368d3a51819SDave May  * Applies a collection method and gathers point neighbour points into dm
1369d3a51819SDave May 
1370d3a51819SDave May  Notes:
13718b8a3813SDave May  Users should call DMSwarmCollectViewDestroy() after
1372d3a51819SDave May  they have finished computations associated with the collected points
1373d3a51819SDave May */
1374d3a51819SDave May 
137574d0cae8SMatthew G. Knepley /*@
1376d3a51819SDave May   DMSwarmCollectViewCreate - Applies a collection method and gathers points
137720f4b53cSBarry Smith   in neighbour ranks into the `DMSWARM`
1378d3a51819SDave May 
137920f4b53cSBarry Smith   Collective
1380d3a51819SDave May 
138160225df5SJacob Faibussowitsch   Input Parameter:
138220f4b53cSBarry Smith . dm - the `DMSWARM`
1383d3a51819SDave May 
1384d3a51819SDave May   Level: advanced
1385d3a51819SDave May 
138620f4b53cSBarry Smith   Notes:
138720f4b53cSBarry Smith   Users should call `DMSwarmCollectViewDestroy()` after
138820f4b53cSBarry Smith   they have finished computations associated with the collected points
138920f4b53cSBarry Smith   Different collect methods are supported. See `DMSwarmSetCollectType()`.
139020f4b53cSBarry Smith 
139120f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmCollectViewDestroy()`, `DMSwarmSetCollectType()`
1392d3a51819SDave May @*/
1393d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmCollectViewCreate(DM dm)
1394d71ae5a4SJacob Faibussowitsch {
13952712d1f2SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
13962712d1f2SDave May   PetscInt  ng;
13972712d1f2SDave May 
1398521f74f9SMatthew G. Knepley   PetscFunctionBegin;
139928b400f6SJacob Faibussowitsch   PetscCheck(!swarm->collect_view_active, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "CollectView currently active");
14009566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetLocalSize(dm, &ng));
1401480eef7bSDave May   switch (swarm->collect_type) {
1402d71ae5a4SJacob Faibussowitsch   case DMSWARM_COLLECT_BASIC:
1403d71ae5a4SJacob Faibussowitsch     PetscCall(DMSwarmMigrate_GlobalToLocal_Basic(dm, &ng));
1404d71ae5a4SJacob Faibussowitsch     break;
1405d71ae5a4SJacob Faibussowitsch   case DMSWARM_COLLECT_DMDABOUNDINGBOX:
1406d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_COLLECT_DMDABOUNDINGBOX not implemented");
1407d71ae5a4SJacob Faibussowitsch   case DMSWARM_COLLECT_GENERAL:
1408d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_COLLECT_GENERAL not implemented");
1409d71ae5a4SJacob Faibussowitsch   default:
1410d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_COLLECT type unknown");
1411480eef7bSDave May   }
1412480eef7bSDave May   swarm->collect_view_active       = PETSC_TRUE;
1413480eef7bSDave May   swarm->collect_view_reset_nlocal = ng;
14143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14152712d1f2SDave May }
14162712d1f2SDave May 
141774d0cae8SMatthew G. Knepley /*@
141820f4b53cSBarry Smith   DMSwarmCollectViewDestroy - Resets the `DMSWARM` to the size prior to calling `DMSwarmCollectViewCreate()`
1419d3a51819SDave May 
142020f4b53cSBarry Smith   Collective
1421d3a51819SDave May 
142260225df5SJacob Faibussowitsch   Input Parameters:
142320f4b53cSBarry Smith . dm - the `DMSWARM`
1424d3a51819SDave May 
1425d3a51819SDave May   Notes:
142620f4b53cSBarry Smith   Users should call `DMSwarmCollectViewCreate()` before this function is called.
1427d3a51819SDave May 
1428d3a51819SDave May   Level: advanced
1429d3a51819SDave May 
143020f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmCollectViewCreate()`, `DMSwarmSetCollectType()`
1431d3a51819SDave May @*/
1432d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmCollectViewDestroy(DM dm)
1433d71ae5a4SJacob Faibussowitsch {
14342712d1f2SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
14352712d1f2SDave May 
1436521f74f9SMatthew G. Knepley   PetscFunctionBegin;
143728b400f6SJacob Faibussowitsch   PetscCheck(swarm->collect_view_active, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "CollectView is currently not active");
14389566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(dm, swarm->collect_view_reset_nlocal, -1));
1439480eef7bSDave May   swarm->collect_view_active = PETSC_FALSE;
14403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14412712d1f2SDave May }
14423454631fSDave May 
144366976f2fSJacob Faibussowitsch static PetscErrorCode DMSwarmSetUpPIC(DM dm)
1444d71ae5a4SJacob Faibussowitsch {
1445f0cdbbbaSDave May   PetscInt dim;
1446f0cdbbbaSDave May 
1447521f74f9SMatthew G. Knepley   PetscFunctionBegin;
14489566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetNumSpecies(dm, 1));
14499566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
145063a3b9bcSJacob Faibussowitsch   PetscCheck(dim >= 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Dimension must be 1,2,3 - found %" PetscInt_FMT, dim);
145163a3b9bcSJacob Faibussowitsch   PetscCheck(dim <= 3, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Dimension must be 1,2,3 - found %" PetscInt_FMT, dim);
14529566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dm, DMSwarmPICField_coor, dim, PETSC_DOUBLE));
14539566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dm, DMSwarmPICField_cellid, 1, PETSC_INT));
14543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1455f0cdbbbaSDave May }
1456f0cdbbbaSDave May 
145774d0cae8SMatthew G. Knepley /*@
145831403186SMatthew G. Knepley   DMSwarmSetPointCoordinatesRandom - Sets initial coordinates for particles in each cell
145931403186SMatthew G. Knepley 
146020f4b53cSBarry Smith   Collective
146131403186SMatthew G. Knepley 
146260225df5SJacob Faibussowitsch   Input Parameters:
146320f4b53cSBarry Smith + dm  - the `DMSWARM`
146420f4b53cSBarry Smith - Npc - The number of particles per cell in the cell `DM`
146531403186SMatthew G. Knepley 
146631403186SMatthew G. Knepley   Level: intermediate
146731403186SMatthew G. Knepley 
146820f4b53cSBarry Smith   Notes:
146920f4b53cSBarry Smith   The user must use `DMSwarmSetCellDM()` to set the cell `DM` first. The particles are placed randomly inside each cell. If only
147020f4b53cSBarry Smith   one particle is in each cell, it is placed at the centroid.
147120f4b53cSBarry Smith 
147220f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmSetCellDM()`
147331403186SMatthew G. Knepley @*/
1474d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmSetPointCoordinatesRandom(DM dm, PetscInt Npc)
1475d71ae5a4SJacob Faibussowitsch {
147631403186SMatthew G. Knepley   DM             cdm;
147731403186SMatthew G. Knepley   PetscRandom    rnd;
147831403186SMatthew G. Knepley   DMPolytopeType ct;
147931403186SMatthew G. Knepley   PetscBool      simplex;
148031403186SMatthew G. Knepley   PetscReal     *centroid, *coords, *xi0, *v0, *J, *invJ, detJ;
148131403186SMatthew G. Knepley   PetscInt       dim, d, cStart, cEnd, c, p;
148231403186SMatthew G. Knepley 
148331403186SMatthew G. Knepley   PetscFunctionBeginUser;
14849566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)dm), &rnd));
14859566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetInterval(rnd, -1.0, 1.0));
14869566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetType(rnd, PETSCRAND48));
148731403186SMatthew G. Knepley 
14889566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetCellDM(dm, &cdm));
14899566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(cdm, &dim));
14909566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(cdm, 0, &cStart, &cEnd));
14919566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCellType(cdm, cStart, &ct));
149231403186SMatthew G. Knepley   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE;
149331403186SMatthew G. Knepley 
14949566063dSJacob Faibussowitsch   PetscCall(PetscMalloc5(dim, &centroid, dim, &xi0, dim, &v0, dim * dim, &J, dim * dim, &invJ));
149531403186SMatthew G. Knepley   for (d = 0; d < dim; ++d) xi0[d] = -1.0;
14969566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
149731403186SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
149831403186SMatthew G. Knepley     if (Npc == 1) {
14999566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFVM(cdm, c, NULL, centroid, NULL));
150031403186SMatthew G. Knepley       for (d = 0; d < dim; ++d) coords[c * dim + d] = centroid[d];
150131403186SMatthew G. Knepley     } else {
15029566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFEM(cdm, c, NULL, v0, J, invJ, &detJ)); /* affine */
150331403186SMatthew G. Knepley       for (p = 0; p < Npc; ++p) {
150431403186SMatthew G. Knepley         const PetscInt n   = c * Npc + p;
150531403186SMatthew G. Knepley         PetscReal      sum = 0.0, refcoords[3];
150631403186SMatthew G. Knepley 
150731403186SMatthew G. Knepley         for (d = 0; d < dim; ++d) {
15089566063dSJacob Faibussowitsch           PetscCall(PetscRandomGetValueReal(rnd, &refcoords[d]));
150931403186SMatthew G. Knepley           sum += refcoords[d];
151031403186SMatthew G. Knepley         }
15119371c9d4SSatish Balay         if (simplex && sum > 0.0)
15129371c9d4SSatish Balay           for (d = 0; d < dim; ++d) refcoords[d] -= PetscSqrtReal(dim) * sum;
151331403186SMatthew G. Knepley         CoordinatesRefToReal(dim, dim, xi0, v0, J, refcoords, &coords[n * dim]);
151431403186SMatthew G. Knepley       }
151531403186SMatthew G. Knepley     }
151631403186SMatthew G. Knepley   }
15179566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
15189566063dSJacob Faibussowitsch   PetscCall(PetscFree5(centroid, xi0, v0, J, invJ));
15199566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rnd));
15203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
152131403186SMatthew G. Knepley }
152231403186SMatthew G. Knepley 
152331403186SMatthew G. Knepley /*@
152420f4b53cSBarry Smith   DMSwarmSetType - Set particular flavor of `DMSWARM`
1525d3a51819SDave May 
152620f4b53cSBarry Smith   Collective
1527d3a51819SDave May 
152860225df5SJacob Faibussowitsch   Input Parameters:
152920f4b53cSBarry Smith + dm    - the `DMSWARM`
153020f4b53cSBarry Smith - stype - the `DMSWARM` type (e.g. `DMSWARM_PIC`)
1531d3a51819SDave May 
1532d3a51819SDave May   Level: advanced
1533d3a51819SDave May 
153420f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmSetMigrateType()`, `DMSwarmSetCollectType()`, `DMSwarmType`, `DMSWARM_PIC`, `DMSWARM_BASIC`
1535d3a51819SDave May @*/
1536d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmSetType(DM dm, DMSwarmType stype)
1537d71ae5a4SJacob Faibussowitsch {
1538f0cdbbbaSDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1539f0cdbbbaSDave May 
1540521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1541f0cdbbbaSDave May   swarm->swarm_type = stype;
154248a46eb9SPierre Jolivet   if (swarm->swarm_type == DMSWARM_PIC) PetscCall(DMSwarmSetUpPIC(dm));
15433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1544f0cdbbbaSDave May }
1545f0cdbbbaSDave May 
154666976f2fSJacob Faibussowitsch static PetscErrorCode DMSetup_Swarm(DM dm)
1547d71ae5a4SJacob Faibussowitsch {
15483454631fSDave May   DM_Swarm   *swarm = (DM_Swarm *)dm->data;
15493454631fSDave May   PetscMPIInt rank;
15503454631fSDave May   PetscInt    p, npoints, *rankval;
15513454631fSDave May 
1552521f74f9SMatthew G. Knepley   PetscFunctionBegin;
15533ba16761SJacob Faibussowitsch   if (swarm->issetup) PetscFunctionReturn(PETSC_SUCCESS);
15543454631fSDave May   swarm->issetup = PETSC_TRUE;
15553454631fSDave May 
1556f0cdbbbaSDave May   if (swarm->swarm_type == DMSWARM_PIC) {
1557f0cdbbbaSDave May     /* check dmcell exists */
155828b400f6SJacob Faibussowitsch     PetscCheck(swarm->dmcell, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "DMSWARM_PIC requires you call DMSwarmSetCellDM");
1559f0cdbbbaSDave May 
1560f0cdbbbaSDave May     if (swarm->dmcell->ops->locatepointssubdomain) {
1561f0cdbbbaSDave May       /* check methods exists for exact ownership identificiation */
15629566063dSJacob Faibussowitsch       PetscCall(PetscInfo(dm, "DMSWARM_PIC: Using method CellDM->ops->LocatePointsSubdomain\n"));
1563f0cdbbbaSDave May       swarm->migrate_type = DMSWARM_MIGRATE_DMCELLEXACT;
1564f0cdbbbaSDave May     } else {
1565f0cdbbbaSDave May       /* check methods exist for point location AND rank neighbor identification */
1566f0cdbbbaSDave May       if (swarm->dmcell->ops->locatepoints) {
15679566063dSJacob Faibussowitsch         PetscCall(PetscInfo(dm, "DMSWARM_PIC: Using method CellDM->LocatePoints\n"));
1568f0cdbbbaSDave May       } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "DMSWARM_PIC requires the method CellDM->ops->locatepoints be defined");
1569f0cdbbbaSDave May 
1570f0cdbbbaSDave May       if (swarm->dmcell->ops->getneighbors) {
15719566063dSJacob Faibussowitsch         PetscCall(PetscInfo(dm, "DMSWARM_PIC: Using method CellDM->GetNeigbors\n"));
1572f0cdbbbaSDave May       } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "DMSWARM_PIC requires the method CellDM->ops->getneighbors be defined");
1573f0cdbbbaSDave May 
1574f0cdbbbaSDave May       swarm->migrate_type = DMSWARM_MIGRATE_DMCELLNSCATTER;
1575f0cdbbbaSDave May     }
1576f0cdbbbaSDave May   }
1577f0cdbbbaSDave May 
15789566063dSJacob Faibussowitsch   PetscCall(DMSwarmFinalizeFieldRegister(dm));
1579f0cdbbbaSDave May 
15803454631fSDave May   /* check some fields were registered */
158108401ef6SPierre Jolivet   PetscCheck(swarm->db->nfields > 2, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "At least one field user must be registered via DMSwarmRegisterXXX()");
15823454631fSDave May 
15833454631fSDave May   /* check local sizes were set */
158408401ef6SPierre Jolivet   PetscCheck(swarm->db->L != -1, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Local sizes must be set via DMSwarmSetLocalSizes()");
15853454631fSDave May 
15863454631fSDave May   /* initialize values in pid and rank placeholders */
15873454631fSDave May   /* TODO: [pid - use MPI_Scan] */
15889566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
15899566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
15909566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
1591ad540459SPierre Jolivet   for (p = 0; p < npoints; p++) rankval[p] = (PetscInt)rank;
15929566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
15933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15943454631fSDave May }
15953454631fSDave May 
1596dc5f5ce9SDave May extern PetscErrorCode DMSwarmSortDestroy(DMSwarmSort *_ctx);
1597dc5f5ce9SDave May 
159866976f2fSJacob Faibussowitsch static PetscErrorCode DMDestroy_Swarm(DM dm)
1599d71ae5a4SJacob Faibussowitsch {
160057795646SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
160157795646SDave May 
160257795646SDave May   PetscFunctionBegin;
16033ba16761SJacob Faibussowitsch   if (--swarm->refct > 0) PetscFunctionReturn(PETSC_SUCCESS);
16049566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketDestroy(&swarm->db));
160548a46eb9SPierre Jolivet   if (swarm->sort_context) PetscCall(DMSwarmSortDestroy(&swarm->sort_context));
16069566063dSJacob Faibussowitsch   PetscCall(PetscFree(swarm));
16073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
160857795646SDave May }
160957795646SDave May 
161066976f2fSJacob Faibussowitsch static PetscErrorCode DMSwarmView_Draw(DM dm, PetscViewer viewer)
1611d71ae5a4SJacob Faibussowitsch {
1612a9ee3421SMatthew G. Knepley   DM         cdm;
1613a9ee3421SMatthew G. Knepley   PetscDraw  draw;
161431403186SMatthew G. Knepley   PetscReal *coords, oldPause, radius = 0.01;
1615a9ee3421SMatthew G. Knepley   PetscInt   Np, p, bs;
1616a9ee3421SMatthew G. Knepley 
1617a9ee3421SMatthew G. Knepley   PetscFunctionBegin;
16189566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, ((PetscObject)dm)->prefix, "-dm_view_swarm_radius", &radius, NULL));
16199566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
16209566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetCellDM(dm, &cdm));
16219566063dSJacob Faibussowitsch   PetscCall(PetscDrawGetPause(draw, &oldPause));
16229566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetPause(draw, 0.0));
16239566063dSJacob Faibussowitsch   PetscCall(DMView(cdm, viewer));
16249566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetPause(draw, oldPause));
1625a9ee3421SMatthew G. Knepley 
16269566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetLocalSize(dm, &Np));
16279566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, &bs, NULL, (void **)&coords));
1628a9ee3421SMatthew G. Knepley   for (p = 0; p < Np; ++p) {
1629a9ee3421SMatthew G. Knepley     const PetscInt i = p * bs;
1630a9ee3421SMatthew G. Knepley 
16319566063dSJacob Faibussowitsch     PetscCall(PetscDrawEllipse(draw, coords[i], coords[i + 1], radius, radius, PETSC_DRAW_BLUE));
1632a9ee3421SMatthew G. Knepley   }
16339566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, &bs, NULL, (void **)&coords));
16349566063dSJacob Faibussowitsch   PetscCall(PetscDrawFlush(draw));
16359566063dSJacob Faibussowitsch   PetscCall(PetscDrawPause(draw));
16369566063dSJacob Faibussowitsch   PetscCall(PetscDrawSave(draw));
16373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1638a9ee3421SMatthew G. Knepley }
1639a9ee3421SMatthew G. Knepley 
1640d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMView_Swarm_Ascii(DM dm, PetscViewer viewer)
1641d71ae5a4SJacob Faibussowitsch {
16426a5217c0SMatthew G. Knepley   PetscViewerFormat format;
16436a5217c0SMatthew G. Knepley   PetscInt         *sizes;
16446a5217c0SMatthew G. Knepley   PetscInt          dim, Np, maxSize = 17;
16456a5217c0SMatthew G. Knepley   MPI_Comm          comm;
16466a5217c0SMatthew G. Knepley   PetscMPIInt       rank, size, p;
16476a5217c0SMatthew G. Knepley   const char       *name;
16486a5217c0SMatthew G. Knepley 
16496a5217c0SMatthew G. Knepley   PetscFunctionBegin;
16506a5217c0SMatthew G. Knepley   PetscCall(PetscViewerGetFormat(viewer, &format));
16516a5217c0SMatthew G. Knepley   PetscCall(DMGetDimension(dm, &dim));
16526a5217c0SMatthew G. Knepley   PetscCall(DMSwarmGetLocalSize(dm, &Np));
16536a5217c0SMatthew G. Knepley   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
16546a5217c0SMatthew G. Knepley   PetscCallMPI(MPI_Comm_rank(comm, &rank));
16556a5217c0SMatthew G. Knepley   PetscCallMPI(MPI_Comm_size(comm, &size));
16566a5217c0SMatthew G. Knepley   PetscCall(PetscObjectGetName((PetscObject)dm, &name));
165763a3b9bcSJacob Faibussowitsch   if (name) PetscCall(PetscViewerASCIIPrintf(viewer, "%s in %" PetscInt_FMT " dimension%s:\n", name, dim, dim == 1 ? "" : "s"));
165863a3b9bcSJacob Faibussowitsch   else PetscCall(PetscViewerASCIIPrintf(viewer, "Swarm in %" PetscInt_FMT " dimension%s:\n", dim, dim == 1 ? "" : "s"));
16596a5217c0SMatthew G. Knepley   if (size < maxSize) PetscCall(PetscCalloc1(size, &sizes));
16606a5217c0SMatthew G. Knepley   else PetscCall(PetscCalloc1(3, &sizes));
16616a5217c0SMatthew G. Knepley   if (size < maxSize) {
16626a5217c0SMatthew G. Knepley     PetscCallMPI(MPI_Gather(&Np, 1, MPIU_INT, sizes, 1, MPIU_INT, 0, comm));
16636a5217c0SMatthew G. Knepley     PetscCall(PetscViewerASCIIPrintf(viewer, "  Number of particles per rank:"));
16646a5217c0SMatthew G. Knepley     for (p = 0; p < size; ++p) {
16656a5217c0SMatthew G. Knepley       if (rank == 0) PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT, sizes[p]));
16666a5217c0SMatthew G. Knepley     }
16676a5217c0SMatthew G. Knepley   } else {
16686a5217c0SMatthew G. Knepley     PetscInt locMinMax[2] = {Np, Np};
16696a5217c0SMatthew G. Knepley 
16706a5217c0SMatthew G. Knepley     PetscCall(PetscGlobalMinMaxInt(comm, locMinMax, sizes));
16716a5217c0SMatthew G. Knepley     PetscCall(PetscViewerASCIIPrintf(viewer, "  Min/Max of particles per rank: %" PetscInt_FMT "/%" PetscInt_FMT, sizes[0], sizes[1]));
16726a5217c0SMatthew G. Knepley   }
16736a5217c0SMatthew G. Knepley   PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
16746a5217c0SMatthew G. Knepley   PetscCall(PetscFree(sizes));
16756a5217c0SMatthew G. Knepley   if (format == PETSC_VIEWER_ASCII_INFO) {
16766a5217c0SMatthew G. Knepley     PetscInt *cell;
16776a5217c0SMatthew G. Knepley 
16786a5217c0SMatthew G. Knepley     PetscCall(PetscViewerASCIIPrintf(viewer, "  Cells containing each particle:\n"));
16796a5217c0SMatthew G. Knepley     PetscCall(PetscViewerASCIIPushSynchronized(viewer));
16806a5217c0SMatthew G. Knepley     PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&cell));
168163a3b9bcSJacob Faibussowitsch     for (p = 0; p < Np; ++p) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "  p%d: %" PetscInt_FMT "\n", p, cell[p]));
16826a5217c0SMatthew G. Knepley     PetscCall(PetscViewerFlush(viewer));
16836a5217c0SMatthew G. Knepley     PetscCall(PetscViewerASCIIPopSynchronized(viewer));
16846a5217c0SMatthew G. Knepley     PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&cell));
16856a5217c0SMatthew G. Knepley   }
16863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16876a5217c0SMatthew G. Knepley }
16886a5217c0SMatthew G. Knepley 
168966976f2fSJacob Faibussowitsch static PetscErrorCode DMView_Swarm(DM dm, PetscViewer viewer)
1690d71ae5a4SJacob Faibussowitsch {
16915f50eb2eSDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
16925627991aSBarry Smith   PetscBool iascii, ibinary, isvtk, isdraw;
16935627991aSBarry Smith #if defined(PETSC_HAVE_HDF5)
16945627991aSBarry Smith   PetscBool ishdf5;
16955627991aSBarry Smith #endif
16965f50eb2eSDave May 
16975f50eb2eSDave May   PetscFunctionBegin;
16985f50eb2eSDave May   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
16995f50eb2eSDave May   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
17009566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
17019566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &ibinary));
17029566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERVTK, &isvtk));
17035627991aSBarry Smith #if defined(PETSC_HAVE_HDF5)
17049566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
17055627991aSBarry Smith #endif
17069566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
17075f50eb2eSDave May   if (iascii) {
17086a5217c0SMatthew G. Knepley     PetscViewerFormat format;
17096a5217c0SMatthew G. Knepley 
17106a5217c0SMatthew G. Knepley     PetscCall(PetscViewerGetFormat(viewer, &format));
17116a5217c0SMatthew G. Knepley     switch (format) {
1712d71ae5a4SJacob Faibussowitsch     case PETSC_VIEWER_ASCII_INFO_DETAIL:
1713d71ae5a4SJacob Faibussowitsch       PetscCall(DMSwarmDataBucketView(PetscObjectComm((PetscObject)dm), swarm->db, NULL, DATABUCKET_VIEW_STDOUT));
1714d71ae5a4SJacob Faibussowitsch       break;
1715d71ae5a4SJacob Faibussowitsch     default:
1716d71ae5a4SJacob Faibussowitsch       PetscCall(DMView_Swarm_Ascii(dm, viewer));
17176a5217c0SMatthew G. Knepley     }
1718f7d195e4SLawrence Mitchell   } else {
17195f50eb2eSDave May #if defined(PETSC_HAVE_HDF5)
172048a46eb9SPierre Jolivet     if (ishdf5) PetscCall(DMSwarmView_HDF5(dm, viewer));
17215f50eb2eSDave May #endif
172248a46eb9SPierre Jolivet     if (isdraw) PetscCall(DMSwarmView_Draw(dm, viewer));
1723f7d195e4SLawrence Mitchell   }
17243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17255f50eb2eSDave May }
17265f50eb2eSDave May 
1727d0c080abSJoseph Pusztay /*@C
172820f4b53cSBarry Smith   DMSwarmGetCellSwarm - Extracts a single cell from the `DMSWARM` object, returns it as a single cell `DMSWARM`.
172920f4b53cSBarry 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.
1730d0c080abSJoseph Pusztay 
1731d0c080abSJoseph Pusztay   Noncollective
1732d0c080abSJoseph Pusztay 
173360225df5SJacob Faibussowitsch   Input Parameters:
173420f4b53cSBarry Smith + sw        - the `DMSWARM`
17355627991aSBarry Smith . cellID    - the integer id of the cell to be extracted and filtered
173620f4b53cSBarry Smith - cellswarm - The `DMSWARM` to receive the cell
1737d0c080abSJoseph Pusztay 
1738d0c080abSJoseph Pusztay   Level: beginner
1739d0c080abSJoseph Pusztay 
17405627991aSBarry Smith   Notes:
174120f4b53cSBarry Smith   This presently only supports `DMSWARM_PIC` type
17425627991aSBarry Smith 
174320f4b53cSBarry Smith   Should be restored with `DMSwarmRestoreCellSwarm()`
1744d0c080abSJoseph Pusztay 
174520f4b53cSBarry Smith   Changes to this cell of the swarm will be lost if they are made prior to restoring this cell.
174620f4b53cSBarry Smith 
174720f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRestoreCellSwarm()`
1748d0c080abSJoseph Pusztay @*/
1749d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMSwarmGetCellSwarm(DM sw, PetscInt cellID, DM cellswarm)
1750d71ae5a4SJacob Faibussowitsch {
1751d0c080abSJoseph Pusztay   DM_Swarm *original = (DM_Swarm *)sw->data;
1752d0c080abSJoseph Pusztay   DMLabel   label;
1753d0c080abSJoseph Pusztay   DM        dmc, subdmc;
1754d0c080abSJoseph Pusztay   PetscInt *pids, particles, dim;
1755d0c080abSJoseph Pusztay 
1756d0c080abSJoseph Pusztay   PetscFunctionBegin;
1757d0c080abSJoseph Pusztay   /* Configure new swarm */
17589566063dSJacob Faibussowitsch   PetscCall(DMSetType(cellswarm, DMSWARM));
17599566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(sw, &dim));
17609566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(cellswarm, dim));
17619566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetType(cellswarm, DMSWARM_PIC));
1762d0c080abSJoseph Pusztay   /* Destroy the unused, unconfigured data bucket to prevent stragglers in memory */
17639566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketDestroy(&((DM_Swarm *)cellswarm->data)->db));
17649566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortGetAccess(sw));
17659566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortGetNumberOfPointsPerCell(sw, cellID, &particles));
17669566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortGetPointsPerCell(sw, cellID, &particles, &pids));
17679566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketCreateFromSubset(original->db, particles, pids, &((DM_Swarm *)cellswarm->data)->db));
17689566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortRestoreAccess(sw));
17699566063dSJacob Faibussowitsch   PetscCall(PetscFree(pids));
17709566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetCellDM(sw, &dmc));
17719566063dSJacob Faibussowitsch   PetscCall(DMLabelCreate(PetscObjectComm((PetscObject)sw), "singlecell", &label));
17729566063dSJacob Faibussowitsch   PetscCall(DMAddLabel(dmc, label));
17739566063dSJacob Faibussowitsch   PetscCall(DMLabelSetValue(label, cellID, 1));
17749566063dSJacob Faibussowitsch   PetscCall(DMPlexFilter(dmc, label, 1, &subdmc));
17759566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetCellDM(cellswarm, subdmc));
17769566063dSJacob Faibussowitsch   PetscCall(DMLabelDestroy(&label));
17773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1778d0c080abSJoseph Pusztay }
1779d0c080abSJoseph Pusztay 
1780d0c080abSJoseph Pusztay /*@C
178120f4b53cSBarry Smith   DMSwarmRestoreCellSwarm - Restores a `DMSWARM` object obtained with `DMSwarmGetCellSwarm()`. All fields are copied back into the parent swarm.
1782d0c080abSJoseph Pusztay 
1783d0c080abSJoseph Pusztay   Noncollective
1784d0c080abSJoseph Pusztay 
178560225df5SJacob Faibussowitsch   Input Parameters:
178620f4b53cSBarry Smith + sw        - the parent `DMSWARM`
1787d0c080abSJoseph Pusztay . cellID    - the integer id of the cell to be copied back into the parent swarm
1788d0c080abSJoseph Pusztay - cellswarm - the cell swarm object
1789d0c080abSJoseph Pusztay 
1790d0c080abSJoseph Pusztay   Level: beginner
1791d0c080abSJoseph Pusztay 
17925627991aSBarry Smith   Note:
179320f4b53cSBarry Smith   This only supports `DMSWARM_PIC` types of `DMSWARM`s
1794d0c080abSJoseph Pusztay 
179520f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmGetCellSwarm()`
1796d0c080abSJoseph Pusztay @*/
1797d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMSwarmRestoreCellSwarm(DM sw, PetscInt cellID, DM cellswarm)
1798d71ae5a4SJacob Faibussowitsch {
1799d0c080abSJoseph Pusztay   DM        dmc;
1800d0c080abSJoseph Pusztay   PetscInt *pids, particles, p;
1801d0c080abSJoseph Pusztay 
1802d0c080abSJoseph Pusztay   PetscFunctionBegin;
18039566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortGetAccess(sw));
18049566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortGetPointsPerCell(sw, cellID, &particles, &pids));
18059566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortRestoreAccess(sw));
1806d0c080abSJoseph Pusztay   /* Pointwise copy of each particle based on pid. The parent swarm may not be altered during this process. */
180748a46eb9SPierre Jolivet   for (p = 0; p < particles; ++p) PetscCall(DMSwarmDataBucketCopyPoint(((DM_Swarm *)cellswarm->data)->db, pids[p], ((DM_Swarm *)sw->data)->db, pids[p]));
1808d0c080abSJoseph Pusztay   /* Free memory, destroy cell dm */
18099566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetCellDM(cellswarm, &dmc));
18109566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmc));
18119566063dSJacob Faibussowitsch   PetscCall(PetscFree(pids));
18123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1813d0c080abSJoseph Pusztay }
1814d0c080abSJoseph Pusztay 
1815d0c080abSJoseph Pusztay PETSC_INTERN PetscErrorCode DMClone_Swarm(DM, DM *);
1816d0c080abSJoseph Pusztay 
1817d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMInitialize_Swarm(DM sw)
1818d71ae5a4SJacob Faibussowitsch {
1819d0c080abSJoseph Pusztay   PetscFunctionBegin;
1820d0c080abSJoseph Pusztay   sw->dim                           = 0;
1821d0c080abSJoseph Pusztay   sw->ops->view                     = DMView_Swarm;
1822d0c080abSJoseph Pusztay   sw->ops->load                     = NULL;
1823d0c080abSJoseph Pusztay   sw->ops->setfromoptions           = NULL;
1824d0c080abSJoseph Pusztay   sw->ops->clone                    = DMClone_Swarm;
1825d0c080abSJoseph Pusztay   sw->ops->setup                    = DMSetup_Swarm;
1826d0c080abSJoseph Pusztay   sw->ops->createlocalsection       = NULL;
1827adc21957SMatthew G. Knepley   sw->ops->createsectionpermutation = NULL;
1828d0c080abSJoseph Pusztay   sw->ops->createdefaultconstraints = NULL;
1829d0c080abSJoseph Pusztay   sw->ops->createglobalvector       = DMCreateGlobalVector_Swarm;
1830d0c080abSJoseph Pusztay   sw->ops->createlocalvector        = DMCreateLocalVector_Swarm;
1831d0c080abSJoseph Pusztay   sw->ops->getlocaltoglobalmapping  = NULL;
1832d0c080abSJoseph Pusztay   sw->ops->createfieldis            = NULL;
1833d0c080abSJoseph Pusztay   sw->ops->createcoordinatedm       = NULL;
1834d0c080abSJoseph Pusztay   sw->ops->getcoloring              = NULL;
1835d0c080abSJoseph Pusztay   sw->ops->creatematrix             = DMCreateMatrix_Swarm;
1836d0c080abSJoseph Pusztay   sw->ops->createinterpolation      = NULL;
1837d0c080abSJoseph Pusztay   sw->ops->createinjection          = NULL;
1838d0c080abSJoseph Pusztay   sw->ops->createmassmatrix         = DMCreateMassMatrix_Swarm;
1839d0c080abSJoseph Pusztay   sw->ops->refine                   = NULL;
1840d0c080abSJoseph Pusztay   sw->ops->coarsen                  = NULL;
1841d0c080abSJoseph Pusztay   sw->ops->refinehierarchy          = NULL;
1842d0c080abSJoseph Pusztay   sw->ops->coarsenhierarchy         = NULL;
1843d0c080abSJoseph Pusztay   sw->ops->globaltolocalbegin       = NULL;
1844d0c080abSJoseph Pusztay   sw->ops->globaltolocalend         = NULL;
1845d0c080abSJoseph Pusztay   sw->ops->localtoglobalbegin       = NULL;
1846d0c080abSJoseph Pusztay   sw->ops->localtoglobalend         = NULL;
1847d0c080abSJoseph Pusztay   sw->ops->destroy                  = DMDestroy_Swarm;
1848d0c080abSJoseph Pusztay   sw->ops->createsubdm              = NULL;
1849d0c080abSJoseph Pusztay   sw->ops->getdimpoints             = NULL;
1850d0c080abSJoseph Pusztay   sw->ops->locatepoints             = NULL;
18513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1852d0c080abSJoseph Pusztay }
1853d0c080abSJoseph Pusztay 
1854d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode DMClone_Swarm(DM dm, DM *newdm)
1855d71ae5a4SJacob Faibussowitsch {
1856d0c080abSJoseph Pusztay   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1857d0c080abSJoseph Pusztay 
1858d0c080abSJoseph Pusztay   PetscFunctionBegin;
1859d0c080abSJoseph Pusztay   swarm->refct++;
1860d0c080abSJoseph Pusztay   (*newdm)->data = swarm;
18619566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)*newdm, DMSWARM));
18629566063dSJacob Faibussowitsch   PetscCall(DMInitialize_Swarm(*newdm));
18632e56dffeSJoe Pusztay   (*newdm)->dim = dm->dim;
18643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1865d0c080abSJoseph Pusztay }
1866d0c080abSJoseph Pusztay 
1867d3a51819SDave May /*MC
1868d3a51819SDave May 
186920f4b53cSBarry Smith  DMSWARM = "swarm" - A `DM` object used to represent arrays of data (fields) of arbitrary data type.
187062741f57SDave May  This implementation was designed for particle methods in which the underlying
1871d3a51819SDave May  data required to be represented is both (i) dynamic in length, (ii) and of arbitrary data type.
1872d3a51819SDave May 
187320f4b53cSBarry Smith  Level: intermediate
187420f4b53cSBarry Smith 
187520f4b53cSBarry Smith   Notes:
187620f4b53cSBarry Smith  User data can be represented by `DMSWARM` through a registering "fields".
187762741f57SDave May  To register a field, the user must provide;
187862741f57SDave May  (a) a unique name;
187962741f57SDave May  (b) the data type (or size in bytes);
188062741f57SDave May  (c) the block size of the data.
1881d3a51819SDave May 
1882d3a51819SDave May  For example, suppose the application requires a unique id, energy, momentum and density to be stored
1883c215a47eSMatthew Knepley  on a set of particles. Then the following code could be used
188420f4b53cSBarry Smith .vb
188520f4b53cSBarry Smith     DMSwarmInitializeFieldRegister(dm)
188620f4b53cSBarry Smith     DMSwarmRegisterPetscDatatypeField(dm,"uid",1,PETSC_LONG);
188720f4b53cSBarry Smith     DMSwarmRegisterPetscDatatypeField(dm,"energy",1,PETSC_REAL);
188820f4b53cSBarry Smith     DMSwarmRegisterPetscDatatypeField(dm,"momentum",3,PETSC_REAL);
188920f4b53cSBarry Smith     DMSwarmRegisterPetscDatatypeField(dm,"density",1,PETSC_FLOAT);
189020f4b53cSBarry Smith     DMSwarmFinalizeFieldRegister(dm)
189120f4b53cSBarry Smith .ve
1892d3a51819SDave May 
189320f4b53cSBarry Smith  The fields represented by `DMSWARM` are dynamic and can be re-sized at any time.
189420f4b53cSBarry Smith  The only restriction imposed by `DMSWARM` is that all fields contain the same number of points.
1895d3a51819SDave May 
1896d3a51819SDave May  To support particle methods, "migration" techniques are provided. These methods migrate data
18975627991aSBarry Smith  between ranks.
1898d3a51819SDave May 
189920f4b53cSBarry Smith  `DMSWARM` supports the methods `DMCreateGlobalVector()` and `DMCreateLocalVector()`.
190020f4b53cSBarry Smith  As a `DMSWARM` may internally define and store values of different data types,
190120f4b53cSBarry Smith  before calling `DMCreateGlobalVector()` or `DMCreateLocalVector()`, the user must inform `DMSWARM` which
190220f4b53cSBarry Smith  fields should be used to define a `Vec` object via
190320f4b53cSBarry Smith    `DMSwarmVectorDefineField()`
1904c215a47eSMatthew Knepley  The specified field can be changed at any time - thereby permitting vectors
1905c215a47eSMatthew Knepley  compatible with different fields to be created.
1906d3a51819SDave May 
190720f4b53cSBarry Smith  A dual representation of fields in the `DMSWARM` and a Vec object is permitted via
190820f4b53cSBarry Smith    `DMSwarmCreateGlobalVectorFromField()`
190920f4b53cSBarry Smith  Here the data defining the field in the `DMSWARM` is shared with a Vec.
1910d3a51819SDave May  This is inherently unsafe if you alter the size of the field at any time between
191120f4b53cSBarry Smith  calls to `DMSwarmCreateGlobalVectorFromField()` and `DMSwarmDestroyGlobalVectorFromField()`.
191220f4b53cSBarry Smith  If the local size of the `DMSWARM` does not match the local size of the global vector
191320f4b53cSBarry Smith  when `DMSwarmDestroyGlobalVectorFromField()` is called, an error is thrown.
1914d3a51819SDave May 
191562741f57SDave May  Additional high-level support is provided for Particle-In-Cell methods.
191620f4b53cSBarry Smith  Please refer to `DMSwarmSetType()`.
191762741f57SDave May 
191820f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMType`, `DMCreate()`, `DMSetType()`
1919d3a51819SDave May M*/
1920d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMCreate_Swarm(DM dm)
1921d71ae5a4SJacob Faibussowitsch {
192257795646SDave May   DM_Swarm *swarm;
192357795646SDave May 
192457795646SDave May   PetscFunctionBegin;
192557795646SDave May   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
19264dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&swarm));
1927f0cdbbbaSDave May   dm->data = swarm;
19289566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketCreate(&swarm->db));
19299566063dSJacob Faibussowitsch   PetscCall(DMSwarmInitializeFieldRegister(dm));
1930d0c080abSJoseph Pusztay   swarm->refct                          = 1;
1931b5bcf523SDave May   swarm->vec_field_set                  = PETSC_FALSE;
19323454631fSDave May   swarm->issetup                        = PETSC_FALSE;
1933480eef7bSDave May   swarm->swarm_type                     = DMSWARM_BASIC;
1934480eef7bSDave May   swarm->migrate_type                   = DMSWARM_MIGRATE_BASIC;
1935480eef7bSDave May   swarm->collect_type                   = DMSWARM_COLLECT_BASIC;
193640c453e9SDave May   swarm->migrate_error_on_missing_point = PETSC_FALSE;
1937f0cdbbbaSDave May   swarm->dmcell                         = NULL;
1938f0cdbbbaSDave May   swarm->collect_view_active            = PETSC_FALSE;
1939f0cdbbbaSDave May   swarm->collect_view_reset_nlocal      = -1;
19409566063dSJacob Faibussowitsch   PetscCall(DMInitialize_Swarm(dm));
19412e956fe4SStefano Zampini   if (SwarmDataFieldId == -1) PetscCall(PetscObjectComposedDataRegister(&SwarmDataFieldId));
19423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
194357795646SDave May }
1944