xref: /petsc/src/dm/impls/swarm/swarm.c (revision 0bf7c1c561ba610fc9f083972c61de6e58d2aa49)
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;
3774d0cae8SMatthew G. Knepley   PetscBool isseq;
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));
459566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5SetTimestep(viewer, seqnum));
469566063dSJacob Faibussowitsch   /* PetscCall(DMSequenceView_HDF5(dm, "time", seqnum, (PetscScalar) seqval, viewer)); */
479566063dSJacob Faibussowitsch   PetscCall(VecViewNative(v, viewer));
489566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject)v, "Nc", PETSC_INT, (void *)&bs));
499566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5PopGroup(viewer));
503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5174d0cae8SMatthew G. Knepley }
5274d0cae8SMatthew G. Knepley 
5366976f2fSJacob Faibussowitsch static PetscErrorCode DMSwarmView_HDF5(DM dm, PetscViewer viewer)
54d71ae5a4SJacob Faibussowitsch {
5574d0cae8SMatthew G. Knepley   Vec       coordinates;
5674d0cae8SMatthew G. Knepley   PetscInt  Np;
5774d0cae8SMatthew G. Knepley   PetscBool isseq;
5874d0cae8SMatthew G. Knepley 
5974d0cae8SMatthew G. Knepley   PetscFunctionBegin;
609566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetSize(dm, &Np));
619566063dSJacob Faibussowitsch   PetscCall(DMSwarmCreateGlobalVectorFromField(dm, DMSwarmPICField_coor, &coordinates));
629566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates"));
639566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5PushGroup(viewer, "/particles"));
649566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)coordinates, VECSEQ, &isseq));
659566063dSJacob Faibussowitsch   PetscCall(VecViewNative(coordinates, viewer));
669566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject)coordinates, "Np", PETSC_INT, (void *)&Np));
679566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5PopGroup(viewer));
689566063dSJacob Faibussowitsch   PetscCall(DMSwarmDestroyGlobalVectorFromField(dm, DMSwarmPICField_coor, &coordinates));
693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7074d0cae8SMatthew G. Knepley }
7174d0cae8SMatthew G. Knepley #endif
7274d0cae8SMatthew G. Knepley 
7366976f2fSJacob Faibussowitsch static PetscErrorCode VecView_Swarm(Vec v, PetscViewer viewer)
74d71ae5a4SJacob Faibussowitsch {
7574d0cae8SMatthew G. Knepley   DM dm;
76f9558f5cSBarry Smith #if defined(PETSC_HAVE_HDF5)
7774d0cae8SMatthew G. Knepley   PetscBool ishdf5;
78f9558f5cSBarry Smith #endif
7974d0cae8SMatthew G. Knepley 
8074d0cae8SMatthew G. Knepley   PetscFunctionBegin;
819566063dSJacob Faibussowitsch   PetscCall(VecGetDM(v, &dm));
8228b400f6SJacob Faibussowitsch   PetscCheck(dm, PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_WRONG, "Vector not generated from a DM");
83f9558f5cSBarry Smith #if defined(PETSC_HAVE_HDF5)
849566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
8574d0cae8SMatthew G. Knepley   if (ishdf5) {
869566063dSJacob Faibussowitsch     PetscCall(VecView_Swarm_HDF5_Internal(v, viewer));
873ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
8874d0cae8SMatthew G. Knepley   }
89f9558f5cSBarry Smith #endif
909566063dSJacob Faibussowitsch   PetscCall(VecViewNative(v, viewer));
913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9274d0cae8SMatthew G. Knepley }
9374d0cae8SMatthew G. Knepley 
94d3a51819SDave May /*@C
95*0bf7c1c5SMatthew G. Knepley   DMSwarmVectorGetField - Gets the field from which to define a `Vec` object
96*0bf7c1c5SMatthew G. Knepley   when `DMCreateLocalVector()`, or `DMCreateGlobalVector()` is called
97*0bf7c1c5SMatthew G. Knepley 
98*0bf7c1c5SMatthew G. Knepley   Not collective
99*0bf7c1c5SMatthew G. Knepley 
100*0bf7c1c5SMatthew G. Knepley   Input Parameter:
101*0bf7c1c5SMatthew G. Knepley . dm - a `DMSWARM`
102*0bf7c1c5SMatthew G. Knepley 
103*0bf7c1c5SMatthew G. Knepley   Output Parameter:
104*0bf7c1c5SMatthew G. Knepley . fieldname - the textual name given to a registered field, or the empty string if it has not been set
105*0bf7c1c5SMatthew G. Knepley 
106*0bf7c1c5SMatthew G. Knepley   Level: beginner
107*0bf7c1c5SMatthew G. Knepley 
108*0bf7c1c5SMatthew G. Knepley .seealso: `DM`, `DMSWARM`, `DMSwarmVectorDefineField()` `DMSwarmRegisterPetscDatatypeField()`, `DMCreateGlobalVector()`, `DMCreateLocalVector()`
109*0bf7c1c5SMatthew G. Knepley @*/
110*0bf7c1c5SMatthew G. Knepley PetscErrorCode DMSwarmVectorGetField(DM dm, const char *fieldname[])
111*0bf7c1c5SMatthew G. Knepley {
112*0bf7c1c5SMatthew G. Knepley   PetscFunctionBegin;
113*0bf7c1c5SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
114*0bf7c1c5SMatthew G. Knepley   PetscAssertPointer(fieldname, 2);
115*0bf7c1c5SMatthew G. Knepley   *fieldname = ((DM_Swarm *)dm->data)->vec_field_name;
116*0bf7c1c5SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
117*0bf7c1c5SMatthew G. Knepley }
118*0bf7c1c5SMatthew G. Knepley 
119*0bf7c1c5SMatthew G. Knepley /*@C
12020f4b53cSBarry Smith   DMSwarmVectorDefineField - Sets the field from which to define a `Vec` object
12120f4b53cSBarry Smith   when `DMCreateLocalVector()`, or `DMCreateGlobalVector()` is called
12257795646SDave May 
12320f4b53cSBarry Smith   Collective
12457795646SDave May 
12560225df5SJacob Faibussowitsch   Input Parameters:
12620f4b53cSBarry Smith + dm        - a `DMSWARM`
12762741f57SDave May - fieldname - the textual name given to a registered field
12857795646SDave May 
129d3a51819SDave May   Level: beginner
13057795646SDave May 
131d3a51819SDave May   Notes:
13220f4b53cSBarry Smith   The field with name `fieldname` must be defined as having a data type of `PetscScalar`.
133e7af74c8SDave May 
13420f4b53cSBarry Smith   This function must be called prior to calling `DMCreateLocalVector()`, `DMCreateGlobalVector()`.
13520f4b53cSBarry Smith   Multiple calls to `DMSwarmVectorDefineField()` are permitted.
136e7af74c8SDave May 
137*0bf7c1c5SMatthew G. Knepley .seealso: `DM`, `DMSWARM`, `DMSwarmVectorGetField()`, `DMSwarmRegisterPetscDatatypeField()`, `DMCreateGlobalVector()`, `DMCreateLocalVector()`
138d3a51819SDave May @*/
139d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmVectorDefineField(DM dm, const char fieldname[])
140d71ae5a4SJacob Faibussowitsch {
141b5bcf523SDave May   DM_Swarm     *swarm = (DM_Swarm *)dm->data;
142b5bcf523SDave May   PetscInt      bs, n;
143b5bcf523SDave May   PetscScalar  *array;
144b5bcf523SDave May   PetscDataType type;
145b5bcf523SDave May 
146a9cbaee5SMatthew G. Knepley   PetscFunctionBegin;
147*0bf7c1c5SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
148*0bf7c1c5SMatthew G. Knepley   if (fieldname) PetscAssertPointer(fieldname, 2);
1499566063dSJacob Faibussowitsch   if (!swarm->issetup) PetscCall(DMSetUp(dm));
1509566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &n, NULL, NULL));
1519566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, fieldname, &bs, &type, (void **)&array));
152b5bcf523SDave May 
153b5bcf523SDave May   /* Check all fields are of type PETSC_REAL or PETSC_SCALAR */
15408401ef6SPierre Jolivet   PetscCheck(type == PETSC_REAL, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only valid for PETSC_REAL");
1559566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(swarm->vec_field_name, PETSC_MAX_PATH_LEN - 1, "%s", fieldname));
156b5bcf523SDave May   swarm->vec_field_set    = PETSC_TRUE;
1571b1ea282SDave May   swarm->vec_field_bs     = bs;
158b5bcf523SDave May   swarm->vec_field_nlocal = n;
1599566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm, fieldname, &bs, &type, (void **)&array));
1603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
161b5bcf523SDave May }
162b5bcf523SDave May 
163cc651181SDave May /* requires DMSwarmDefineFieldVector has been called */
16466976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateGlobalVector_Swarm(DM dm, Vec *vec)
165d71ae5a4SJacob Faibussowitsch {
166b5bcf523SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
167b5bcf523SDave May   Vec       x;
168b5bcf523SDave May   char      name[PETSC_MAX_PATH_LEN];
169b5bcf523SDave May 
170a9cbaee5SMatthew G. Knepley   PetscFunctionBegin;
1719566063dSJacob Faibussowitsch   if (!swarm->issetup) PetscCall(DMSetUp(dm));
17228b400f6SJacob Faibussowitsch   PetscCheck(swarm->vec_field_set, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Must call DMSwarmVectorDefineField first");
17308401ef6SPierre 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 */
174cc651181SDave May 
1759566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN - 1, "DMSwarmField_%s", swarm->vec_field_name));
1769566063dSJacob Faibussowitsch   PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), &x));
1779566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)x, name));
1789566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(x, swarm->db->L * swarm->vec_field_bs, PETSC_DETERMINE));
1799566063dSJacob Faibussowitsch   PetscCall(VecSetBlockSize(x, swarm->vec_field_bs));
1809566063dSJacob Faibussowitsch   PetscCall(VecSetDM(x, dm));
1819566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(x));
1829566063dSJacob Faibussowitsch   PetscCall(VecSetDM(x, dm));
1839566063dSJacob Faibussowitsch   PetscCall(VecSetOperation(x, VECOP_VIEW, (void (*)(void))VecView_Swarm));
184b5bcf523SDave May   *vec = x;
1853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
186b5bcf523SDave May }
187b5bcf523SDave May 
188b5bcf523SDave May /* requires DMSwarmDefineFieldVector has been called */
18966976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateLocalVector_Swarm(DM dm, Vec *vec)
190d71ae5a4SJacob Faibussowitsch {
191b5bcf523SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
192b5bcf523SDave May   Vec       x;
193b5bcf523SDave May   char      name[PETSC_MAX_PATH_LEN];
194b5bcf523SDave May 
195a9cbaee5SMatthew G. Knepley   PetscFunctionBegin;
1969566063dSJacob Faibussowitsch   if (!swarm->issetup) PetscCall(DMSetUp(dm));
19728b400f6SJacob Faibussowitsch   PetscCheck(swarm->vec_field_set, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Must call DMSwarmVectorDefineField first");
19808401ef6SPierre 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");
199cc651181SDave May 
2009566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN - 1, "DMSwarmField_%s", swarm->vec_field_name));
2019566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_SELF, &x));
2029566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)x, name));
2039566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(x, swarm->db->L * swarm->vec_field_bs, PETSC_DETERMINE));
2049566063dSJacob Faibussowitsch   PetscCall(VecSetBlockSize(x, swarm->vec_field_bs));
2059566063dSJacob Faibussowitsch   PetscCall(VecSetDM(x, dm));
2069566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(x));
207b5bcf523SDave May   *vec = x;
2083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
209b5bcf523SDave May }
210b5bcf523SDave May 
211d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSwarmDestroyVectorFromField_Private(DM dm, const char fieldname[], Vec *vec)
212d71ae5a4SJacob Faibussowitsch {
213fb1bcc12SMatthew G. Knepley   DM_Swarm        *swarm = (DM_Swarm *)dm->data;
21477048351SPatrick Sanan   DMSwarmDataField gfield;
2152e956fe4SStefano Zampini   PetscInt         bs, nlocal, fid = -1, cfid = -2;
2162e956fe4SStefano Zampini   PetscBool        flg;
217d3a51819SDave May 
218fb1bcc12SMatthew G. Knepley   PetscFunctionBegin;
2192e956fe4SStefano Zampini   /* check vector is an inplace array */
2202e956fe4SStefano Zampini   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldIdByName(swarm->db, fieldname, &fid));
2212e956fe4SStefano Zampini   PetscCall(PetscObjectComposedDataGetInt((PetscObject)*vec, SwarmDataFieldId, cfid, flg));
2222e956fe4SStefano 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);
2239566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(*vec, &nlocal));
2249566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(*vec, &bs));
22508401ef6SPierre 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");
2269566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield));
2279566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataFieldRestoreAccess(gfield));
2288df78e51SMark Adams   PetscCall(VecResetArray(*vec));
2299566063dSJacob Faibussowitsch   PetscCall(VecDestroy(vec));
2303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
231fb1bcc12SMatthew G. Knepley }
232fb1bcc12SMatthew G. Knepley 
233d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSwarmCreateVectorFromField_Private(DM dm, const char fieldname[], MPI_Comm comm, Vec *vec)
234d71ae5a4SJacob Faibussowitsch {
235fb1bcc12SMatthew G. Knepley   DM_Swarm     *swarm = (DM_Swarm *)dm->data;
236fb1bcc12SMatthew G. Knepley   PetscDataType type;
237fb1bcc12SMatthew G. Knepley   PetscScalar  *array;
2382e956fe4SStefano Zampini   PetscInt      bs, n, fid;
239fb1bcc12SMatthew G. Knepley   char          name[PETSC_MAX_PATH_LEN];
240e4fbd051SBarry Smith   PetscMPIInt   size;
2418df78e51SMark Adams   PetscBool     iscuda, iskokkos;
242fb1bcc12SMatthew G. Knepley 
243fb1bcc12SMatthew G. Knepley   PetscFunctionBegin;
2449566063dSJacob Faibussowitsch   if (!swarm->issetup) PetscCall(DMSetUp(dm));
2459566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &n, NULL, NULL));
2469566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, fieldname, &bs, &type, (void **)&array));
24708401ef6SPierre Jolivet   PetscCheck(type == PETSC_REAL, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only valid for PETSC_REAL");
248fb1bcc12SMatthew G. Knepley 
2499566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
2508df78e51SMark Adams   PetscCall(PetscStrcmp(dm->vectype, VECKOKKOS, &iskokkos));
2518df78e51SMark Adams   PetscCall(PetscStrcmp(dm->vectype, VECCUDA, &iscuda));
2528df78e51SMark Adams   PetscCall(VecCreate(comm, vec));
2538df78e51SMark Adams   PetscCall(VecSetSizes(*vec, n * bs, PETSC_DETERMINE));
2548df78e51SMark Adams   PetscCall(VecSetBlockSize(*vec, bs));
2558df78e51SMark Adams   if (iskokkos) PetscCall(VecSetType(*vec, VECKOKKOS));
2568df78e51SMark Adams   else if (iscuda) PetscCall(VecSetType(*vec, VECCUDA));
2578df78e51SMark Adams   else PetscCall(VecSetType(*vec, VECSTANDARD));
2588df78e51SMark Adams   PetscCall(VecPlaceArray(*vec, array));
2598df78e51SMark Adams 
2609566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN - 1, "DMSwarmSharedField_%s", fieldname));
2619566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)*vec, name));
262fb1bcc12SMatthew G. Knepley 
263fb1bcc12SMatthew G. Knepley   /* Set guard */
2642e956fe4SStefano Zampini   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldIdByName(swarm->db, fieldname, &fid));
2652e956fe4SStefano Zampini   PetscCall(PetscObjectComposedDataSetInt((PetscObject)*vec, SwarmDataFieldId, fid));
26674d0cae8SMatthew G. Knepley 
2679566063dSJacob Faibussowitsch   PetscCall(VecSetDM(*vec, dm));
2689566063dSJacob Faibussowitsch   PetscCall(VecSetOperation(*vec, VECOP_VIEW, (void (*)(void))VecView_Swarm));
2693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
270fb1bcc12SMatthew G. Knepley }
271fb1bcc12SMatthew G. Knepley 
2720643ed39SMatthew G. Knepley /* This creates a "mass matrix" between a finite element and particle space. If a finite element interpolant is given by
2730643ed39SMatthew G. Knepley 
2740643ed39SMatthew G. Knepley      \hat f = \sum_i f_i \phi_i
2750643ed39SMatthew G. Knepley 
2760643ed39SMatthew G. Knepley    and a particle function is given by
2770643ed39SMatthew G. Knepley 
2780643ed39SMatthew G. Knepley      f = \sum_i w_i \delta(x - x_i)
2790643ed39SMatthew G. Knepley 
2800643ed39SMatthew G. Knepley    then we want to require that
2810643ed39SMatthew G. Knepley 
2820643ed39SMatthew G. Knepley      M \hat f = M_p f
2830643ed39SMatthew G. Knepley 
2840643ed39SMatthew G. Knepley    where the particle mass matrix is given by
2850643ed39SMatthew G. Knepley 
2860643ed39SMatthew G. Knepley      (M_p)_{ij} = \int \phi_i \delta(x - x_j)
2870643ed39SMatthew G. Knepley 
2880643ed39SMatthew G. Knepley    The way Dave May does particles, they amount to quadratue weights rather than delta functions, so he has |J| is in
2890643ed39SMatthew G. Knepley    his integral. We allow this with the boolean flag.
2900643ed39SMatthew G. Knepley */
291d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSwarmComputeMassMatrix_Private(DM dmc, DM dmf, Mat mass, PetscBool useDeltaFunction, void *ctx)
292d71ae5a4SJacob Faibussowitsch {
29383c47955SMatthew G. Knepley   const char  *name = "Mass Matrix";
2940643ed39SMatthew G. Knepley   MPI_Comm     comm;
29583c47955SMatthew G. Knepley   PetscDS      prob;
29683c47955SMatthew G. Knepley   PetscSection fsection, globalFSection;
297e8f14785SLisandro Dalcin   PetscHSetIJ  ht;
2980643ed39SMatthew G. Knepley   PetscLayout  rLayout, colLayout;
29983c47955SMatthew G. Knepley   PetscInt    *dnz, *onz;
300adb2528bSMark Adams   PetscInt     locRows, locCols, rStart, colStart, colEnd, *rowIDXs;
3010643ed39SMatthew G. Knepley   PetscReal   *xi, *v0, *J, *invJ, detJ = 1.0, v0ref[3] = {-1.0, -1.0, -1.0};
30283c47955SMatthew G. Knepley   PetscScalar *elemMat;
303*0bf7c1c5SMatthew G. Knepley   PetscInt     dim, Nf, field, cStart, cEnd, cell, totDim, maxC = 0, totNc = 0;
30483c47955SMatthew G. Knepley 
30583c47955SMatthew G. Knepley   PetscFunctionBegin;
3069566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)mass, &comm));
3079566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dmf, &dim));
3089566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dmf, &prob));
3099566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(prob, &Nf));
3109566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(prob, &totDim));
3119566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(dim, &v0, dim * dim, &J, dim * dim, &invJ));
3129566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dmf, &fsection));
3139566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalSection(dmf, &globalFSection));
3149566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd));
3159566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(mass, &locRows, &locCols));
31683c47955SMatthew G. Knepley 
3179566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreate(comm, &colLayout));
3189566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetLocalSize(colLayout, locCols));
3199566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(colLayout, 1));
3209566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(colLayout));
3219566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetRange(colLayout, &colStart, &colEnd));
3229566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&colLayout));
3230643ed39SMatthew G. Knepley 
3249566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreate(comm, &rLayout));
3259566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetLocalSize(rLayout, locRows));
3269566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(rLayout, 1));
3279566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(rLayout));
3289566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetRange(rLayout, &rStart, NULL));
3299566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&rLayout));
3300643ed39SMatthew G. Knepley 
3319566063dSJacob Faibussowitsch   PetscCall(PetscCalloc2(locRows, &dnz, locRows, &onz));
3329566063dSJacob Faibussowitsch   PetscCall(PetscHSetIJCreate(&ht));
33353e60ab4SJoseph Pusztay 
3349566063dSJacob Faibussowitsch   PetscCall(PetscSynchronizedFlush(comm, NULL));
335*0bf7c1c5SMatthew G. Knepley   for (field = 0; field < Nf; ++field) {
336*0bf7c1c5SMatthew G. Knepley     PetscObject  obj;
337*0bf7c1c5SMatthew G. Knepley     PetscClassId id;
338*0bf7c1c5SMatthew G. Knepley     PetscInt     Nc;
339*0bf7c1c5SMatthew G. Knepley 
340*0bf7c1c5SMatthew G. Knepley     PetscCall(PetscDSGetDiscretization(prob, field, &obj));
341*0bf7c1c5SMatthew G. Knepley     PetscCall(PetscObjectGetClassId(obj, &id));
342*0bf7c1c5SMatthew G. Knepley     if (id == PETSCFE_CLASSID) PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc));
343*0bf7c1c5SMatthew G. Knepley     else PetscCall(PetscFVGetNumComponents((PetscFV)obj, &Nc));
344*0bf7c1c5SMatthew G. Knepley     totNc += Nc;
345*0bf7c1c5SMatthew G. Knepley   }
3460643ed39SMatthew G. Knepley   /* count non-zeros */
3479566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortGetAccess(dmc));
34883c47955SMatthew G. Knepley   for (field = 0; field < Nf; ++field) {
349*0bf7c1c5SMatthew G. Knepley     PetscObject  obj;
350*0bf7c1c5SMatthew G. Knepley     PetscClassId id;
351*0bf7c1c5SMatthew G. Knepley     PetscInt     Nc;
352*0bf7c1c5SMatthew G. Knepley 
353*0bf7c1c5SMatthew G. Knepley     PetscCall(PetscDSGetDiscretization(prob, field, &obj));
354*0bf7c1c5SMatthew G. Knepley     PetscCall(PetscObjectGetClassId(obj, &id));
355*0bf7c1c5SMatthew G. Knepley     if (id == PETSCFE_CLASSID) PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc));
356*0bf7c1c5SMatthew G. Knepley     else PetscCall(PetscFVGetNumComponents((PetscFV)obj, &Nc));
357*0bf7c1c5SMatthew G. Knepley 
35883c47955SMatthew G. Knepley     for (cell = cStart; cell < cEnd; ++cell) {
3590643ed39SMatthew G. Knepley       PetscInt *findices, *cindices; /* fine is vertices, coarse is particles */
36083c47955SMatthew G. Knepley       PetscInt  numFIndices, numCIndices;
36183c47955SMatthew G. Knepley 
3629566063dSJacob Faibussowitsch       PetscCall(DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
3639566063dSJacob Faibussowitsch       PetscCall(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices));
364fc7c92abSMatthew G. Knepley       maxC = PetscMax(maxC, numCIndices);
36583c47955SMatthew G. Knepley       {
366e8f14785SLisandro Dalcin         PetscHashIJKey key;
367e8f14785SLisandro Dalcin         PetscBool      missing;
368*0bf7c1c5SMatthew G. Knepley         for (PetscInt i = 0; i < numFIndices; ++i) {
369adb2528bSMark Adams           key.j = findices[i]; /* global column (from Plex) */
370adb2528bSMark Adams           if (key.j >= 0) {
37183c47955SMatthew G. Knepley             /* Get indices for coarse elements */
372*0bf7c1c5SMatthew G. Knepley             for (PetscInt j = 0; j < numCIndices; ++j) {
373*0bf7c1c5SMatthew G. Knepley               for (PetscInt c = 0; c < Nc; ++c) {
374*0bf7c1c5SMatthew G. Knepley                 // TODO Need field offset on particle here
375*0bf7c1c5SMatthew G. Knepley                 key.i = cindices[j] * totNc + c + rStart; /* global cols (from Swarm) */
376adb2528bSMark Adams                 if (key.i < 0) continue;
3779566063dSJacob Faibussowitsch                 PetscCall(PetscHSetIJQueryAdd(ht, key, &missing));
37883c47955SMatthew G. Knepley                 if (missing) {
3790643ed39SMatthew G. Knepley                   if ((key.j >= colStart) && (key.j < colEnd)) ++dnz[key.i - rStart];
380e8f14785SLisandro Dalcin                   else ++onz[key.i - rStart];
38163a3b9bcSJacob Faibussowitsch                 } else SETERRQ(PetscObjectComm((PetscObject)dmf), PETSC_ERR_SUP, "Set new value at %" PetscInt_FMT ",%" PetscInt_FMT, key.i, key.j);
38283c47955SMatthew G. Knepley               }
383fc7c92abSMatthew G. Knepley             }
384fc7c92abSMatthew G. Knepley           }
385*0bf7c1c5SMatthew G. Knepley         }
3869566063dSJacob Faibussowitsch         PetscCall(PetscFree(cindices));
38783c47955SMatthew G. Knepley       }
3889566063dSJacob Faibussowitsch       PetscCall(DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
38983c47955SMatthew G. Knepley     }
39083c47955SMatthew G. Knepley   }
3919566063dSJacob Faibussowitsch   PetscCall(PetscHSetIJDestroy(&ht));
3929566063dSJacob Faibussowitsch   PetscCall(MatXAIJSetPreallocation(mass, 1, dnz, onz, NULL, NULL));
3939566063dSJacob Faibussowitsch   PetscCall(MatSetOption(mass, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
3949566063dSJacob Faibussowitsch   PetscCall(PetscFree2(dnz, onz));
395*0bf7c1c5SMatthew G. Knepley   PetscCall(PetscMalloc3(maxC * totNc * totDim, &elemMat, maxC * totNc, &rowIDXs, maxC * dim, &xi));
39683c47955SMatthew G. Knepley   for (field = 0; field < Nf; ++field) {
397ef0bb6c7SMatthew G. Knepley     PetscTabulation Tcoarse;
39883c47955SMatthew G. Knepley     PetscObject     obj;
399ad9634fcSMatthew G. Knepley     PetscClassId    id;
400ea7032a0SMatthew G. Knepley     PetscReal      *fieldVals;
401*0bf7c1c5SMatthew G. Knepley     PetscInt        Nc;
40283c47955SMatthew G. Knepley 
4039566063dSJacob Faibussowitsch     PetscCall(PetscDSGetDiscretization(prob, field, &obj));
404ad9634fcSMatthew G. Knepley     PetscCall(PetscObjectGetClassId(obj, &id));
405ad9634fcSMatthew G. Knepley     if (id == PETSCFE_CLASSID) PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc));
406ad9634fcSMatthew G. Knepley     else PetscCall(PetscFVGetNumComponents((PetscFV)obj, &Nc));
407ea7032a0SMatthew G. Knepley 
408ea7032a0SMatthew G. Knepley     PetscCall(DMSwarmGetField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **)&fieldVals));
40983c47955SMatthew G. Knepley     for (cell = cStart; cell < cEnd; ++cell) {
41083c47955SMatthew G. Knepley       PetscInt *findices, *cindices;
41183c47955SMatthew G. Knepley       PetscInt  numFIndices, numCIndices;
41283c47955SMatthew G. Knepley 
4130643ed39SMatthew G. Knepley       /* TODO: Use DMField instead of assuming affine */
4149566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ));
4159566063dSJacob Faibussowitsch       PetscCall(DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
4169566063dSJacob Faibussowitsch       PetscCall(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices));
417*0bf7c1c5SMatthew G. Knepley       for (PetscInt j = 0; j < numCIndices; ++j) CoordinatesRealToRef(dim, dim, v0ref, v0, invJ, &fieldVals[cindices[j] * dim], &xi[j * dim]);
418ad9634fcSMatthew G. Knepley       if (id == PETSCFE_CLASSID) PetscCall(PetscFECreateTabulation((PetscFE)obj, 1, numCIndices, xi, 0, &Tcoarse));
419ad9634fcSMatthew G. Knepley       else PetscCall(PetscFVCreateTabulation((PetscFV)obj, 1, numCIndices, xi, 0, &Tcoarse));
42083c47955SMatthew G. Knepley       /* Get elemMat entries by multiplying by weight */
421*0bf7c1c5SMatthew G. Knepley       PetscCall(PetscArrayzero(elemMat, numCIndices * Nc * totDim));
422*0bf7c1c5SMatthew G. Knepley       for (PetscInt i = 0; i < numFIndices / Nc; ++i) {
423*0bf7c1c5SMatthew G. Knepley         for (PetscInt j = 0; j < numCIndices; ++j) {
424*0bf7c1c5SMatthew G. Knepley           for (PetscInt c = 0; c < Nc; ++c) {
425*0bf7c1c5SMatthew G. Knepley             // TODO Need field offset on particle and field here
4260643ed39SMatthew G. Knepley             /* B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c */
427*0bf7c1c5SMatthew G. Knepley             elemMat[(j * totNc + c) * numFIndices + i * Nc + c] += Tcoarse->T[0][(j * numFIndices + i * Nc + c) * Nc + c] * (useDeltaFunction ? 1.0 : detJ);
42883c47955SMatthew G. Knepley           }
4290643ed39SMatthew G. Knepley         }
4300643ed39SMatthew G. Knepley       }
431*0bf7c1c5SMatthew G. Knepley       for (PetscInt j = 0; j < numCIndices; ++j)
432*0bf7c1c5SMatthew G. Knepley         // TODO Need field offset on particle here
433*0bf7c1c5SMatthew G. Knepley         for (PetscInt c = 0; c < Nc; ++c) rowIDXs[j * Nc + c] = cindices[j] * totNc + c + rStart;
434*0bf7c1c5SMatthew G. Knepley       if (0) PetscCall(DMPrintCellMatrix(cell, name, numCIndices * Nc, numFIndices, elemMat));
435*0bf7c1c5SMatthew G. Knepley       PetscCall(MatSetValues(mass, numCIndices * Nc, rowIDXs, numFIndices, findices, elemMat, ADD_VALUES));
4369566063dSJacob Faibussowitsch       PetscCall(PetscFree(cindices));
4379566063dSJacob Faibussowitsch       PetscCall(DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
4389566063dSJacob Faibussowitsch       PetscCall(PetscTabulationDestroy(&Tcoarse));
43983c47955SMatthew G. Knepley     }
440ea7032a0SMatthew G. Knepley     PetscCall(DMSwarmRestoreField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **)&fieldVals));
44183c47955SMatthew G. Knepley   }
4429566063dSJacob Faibussowitsch   PetscCall(PetscFree3(elemMat, rowIDXs, xi));
4439566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortRestoreAccess(dmc));
4449566063dSJacob Faibussowitsch   PetscCall(PetscFree3(v0, J, invJ));
4459566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mass, MAT_FINAL_ASSEMBLY));
4469566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mass, MAT_FINAL_ASSEMBLY));
4473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
44883c47955SMatthew G. Knepley }
44983c47955SMatthew G. Knepley 
450d0c080abSJoseph Pusztay /* Returns empty matrix for use with SNES FD */
451d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCreateMatrix_Swarm(DM sw, Mat *m)
452d71ae5a4SJacob Faibussowitsch {
453d0c080abSJoseph Pusztay   Vec      field;
454d0c080abSJoseph Pusztay   PetscInt size;
455d0c080abSJoseph Pusztay 
456d0c080abSJoseph Pusztay   PetscFunctionBegin;
4579566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(sw, &field));
4589566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(field, &size));
4599566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(sw, &field));
4609566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, m));
4619566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(*m));
4629566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*m, PETSC_DECIDE, PETSC_DECIDE, size, size));
4639566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(*m, 1, NULL));
4649566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(*m));
4659566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(*m, MAT_FINAL_ASSEMBLY));
4669566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(*m, MAT_FINAL_ASSEMBLY));
4679566063dSJacob Faibussowitsch   PetscCall(MatShift(*m, 1.0));
4689566063dSJacob Faibussowitsch   PetscCall(MatSetDM(*m, sw));
4693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
470d0c080abSJoseph Pusztay }
471d0c080abSJoseph Pusztay 
472adb2528bSMark Adams /* FEM cols, Particle rows */
473d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCreateMassMatrix_Swarm(DM dmCoarse, DM dmFine, Mat *mass)
474d71ae5a4SJacob Faibussowitsch {
475895a1698SMatthew G. Knepley   PetscSection gsf;
476*0bf7c1c5SMatthew G. Knepley   PetscInt     m, n, Np, bs = 1;
47783c47955SMatthew G. Knepley   void        *ctx;
478*0bf7c1c5SMatthew G. Knepley   PetscBool    set  = ((DM_Swarm *)dmCoarse->data)->vec_field_set;
479*0bf7c1c5SMatthew G. Knepley   char        *name = ((DM_Swarm *)dmCoarse->data)->vec_field_name;
48083c47955SMatthew G. Knepley 
48183c47955SMatthew G. Knepley   PetscFunctionBegin;
4829566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalSection(dmFine, &gsf));
4839566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetConstrainedStorageSize(gsf, &m));
484*0bf7c1c5SMatthew G. Knepley   PetscCall(DMSwarmGetLocalSize(dmCoarse, &Np));
485*0bf7c1c5SMatthew G. Knepley   // TODO Include all fields
486*0bf7c1c5SMatthew G. Knepley   if (set) PetscCall(DMSwarmGetFieldInfo(dmCoarse, name, &bs, NULL));
487*0bf7c1c5SMatthew G. Knepley   n = Np * bs;
4889566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)dmCoarse), mass));
4899566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*mass, n, m, PETSC_DETERMINE, PETSC_DETERMINE));
4909566063dSJacob Faibussowitsch   PetscCall(MatSetType(*mass, dmCoarse->mattype));
4919566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(dmFine, &ctx));
49283c47955SMatthew G. Knepley 
4939566063dSJacob Faibussowitsch   PetscCall(DMSwarmComputeMassMatrix_Private(dmCoarse, dmFine, *mass, PETSC_TRUE, ctx));
4949566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(*mass, NULL, "-mass_mat_view"));
4953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
49683c47955SMatthew G. Knepley }
49783c47955SMatthew G. Knepley 
498d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSwarmComputeMassMatrixSquare_Private(DM dmc, DM dmf, Mat mass, PetscBool useDeltaFunction, void *ctx)
499d71ae5a4SJacob Faibussowitsch {
5004cc7f7b2SMatthew G. Knepley   const char  *name = "Mass Matrix Square";
5014cc7f7b2SMatthew G. Knepley   MPI_Comm     comm;
5024cc7f7b2SMatthew G. Knepley   PetscDS      prob;
5034cc7f7b2SMatthew G. Knepley   PetscSection fsection, globalFSection;
5044cc7f7b2SMatthew G. Knepley   PetscHSetIJ  ht;
5054cc7f7b2SMatthew G. Knepley   PetscLayout  rLayout, colLayout;
5064cc7f7b2SMatthew G. Knepley   PetscInt    *dnz, *onz, *adj, depth, maxConeSize, maxSupportSize, maxAdjSize;
5074cc7f7b2SMatthew G. Knepley   PetscInt     locRows, locCols, rStart, colStart, colEnd, *rowIDXs;
5084cc7f7b2SMatthew G. Knepley   PetscReal   *xi, *v0, *J, *invJ, detJ = 1.0, v0ref[3] = {-1.0, -1.0, -1.0};
5094cc7f7b2SMatthew G. Knepley   PetscScalar *elemMat, *elemMatSq;
5104cc7f7b2SMatthew G. Knepley   PetscInt     cdim, Nf, field, cStart, cEnd, cell, totDim, maxC = 0;
5114cc7f7b2SMatthew G. Knepley 
5124cc7f7b2SMatthew G. Knepley   PetscFunctionBegin;
5139566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)mass, &comm));
5149566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dmf, &cdim));
5159566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dmf, &prob));
5169566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(prob, &Nf));
5179566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(prob, &totDim));
5189566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(cdim, &v0, cdim * cdim, &J, cdim * cdim, &invJ));
5199566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dmf, &fsection));
5209566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalSection(dmf, &globalFSection));
5219566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd));
5229566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(mass, &locRows, &locCols));
5234cc7f7b2SMatthew G. Knepley 
5249566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreate(comm, &colLayout));
5259566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetLocalSize(colLayout, locCols));
5269566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(colLayout, 1));
5279566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(colLayout));
5289566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetRange(colLayout, &colStart, &colEnd));
5299566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&colLayout));
5304cc7f7b2SMatthew G. Knepley 
5319566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreate(comm, &rLayout));
5329566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetLocalSize(rLayout, locRows));
5339566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(rLayout, 1));
5349566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(rLayout));
5359566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetRange(rLayout, &rStart, NULL));
5369566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&rLayout));
5374cc7f7b2SMatthew G. Knepley 
5389566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(dmf, &depth));
5399566063dSJacob Faibussowitsch   PetscCall(DMPlexGetMaxSizes(dmf, &maxConeSize, &maxSupportSize));
5404cc7f7b2SMatthew G. Knepley   maxAdjSize = PetscPowInt(maxConeSize * maxSupportSize, depth);
5419566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(maxAdjSize, &adj));
5424cc7f7b2SMatthew G. Knepley 
5439566063dSJacob Faibussowitsch   PetscCall(PetscCalloc2(locRows, &dnz, locRows, &onz));
5449566063dSJacob Faibussowitsch   PetscCall(PetscHSetIJCreate(&ht));
5454cc7f7b2SMatthew G. Knepley   /* Count nonzeros
5464cc7f7b2SMatthew G. Knepley        This is just FVM++, but we cannot use the Plex P0 allocation since unknowns in a cell will not be contiguous
5474cc7f7b2SMatthew G. Knepley   */
5489566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortGetAccess(dmc));
5494cc7f7b2SMatthew G. Knepley   for (cell = cStart; cell < cEnd; ++cell) {
5504cc7f7b2SMatthew G. Knepley     PetscInt  i;
5514cc7f7b2SMatthew G. Knepley     PetscInt *cindices;
5524cc7f7b2SMatthew G. Knepley     PetscInt  numCIndices;
5534cc7f7b2SMatthew G. Knepley #if 0
5544cc7f7b2SMatthew G. Knepley     PetscInt  adjSize = maxAdjSize, a, j;
5554cc7f7b2SMatthew G. Knepley #endif
5564cc7f7b2SMatthew G. Knepley 
5579566063dSJacob Faibussowitsch     PetscCall(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices));
5584cc7f7b2SMatthew G. Knepley     maxC = PetscMax(maxC, numCIndices);
5594cc7f7b2SMatthew G. Knepley     /* Diagonal block */
560ad540459SPierre Jolivet     for (i = 0; i < numCIndices; ++i) dnz[cindices[i]] += numCIndices;
5614cc7f7b2SMatthew G. Knepley #if 0
5624cc7f7b2SMatthew G. Knepley     /* Off-diagonal blocks */
5639566063dSJacob Faibussowitsch     PetscCall(DMPlexGetAdjacency(dmf, cell, &adjSize, &adj));
5644cc7f7b2SMatthew G. Knepley     for (a = 0; a < adjSize; ++a) {
5654cc7f7b2SMatthew G. Knepley       if (adj[a] >= cStart && adj[a] < cEnd && adj[a] != cell) {
5664cc7f7b2SMatthew G. Knepley         const PetscInt ncell = adj[a];
5674cc7f7b2SMatthew G. Knepley         PetscInt      *ncindices;
5684cc7f7b2SMatthew G. Knepley         PetscInt       numNCIndices;
5694cc7f7b2SMatthew G. Knepley 
5709566063dSJacob Faibussowitsch         PetscCall(DMSwarmSortGetPointsPerCell(dmc, ncell, &numNCIndices, &ncindices));
5714cc7f7b2SMatthew G. Knepley         {
5724cc7f7b2SMatthew G. Knepley           PetscHashIJKey key;
5734cc7f7b2SMatthew G. Knepley           PetscBool      missing;
5744cc7f7b2SMatthew G. Knepley 
5754cc7f7b2SMatthew G. Knepley           for (i = 0; i < numCIndices; ++i) {
5764cc7f7b2SMatthew G. Knepley             key.i = cindices[i] + rStart; /* global rows (from Swarm) */
5774cc7f7b2SMatthew G. Knepley             if (key.i < 0) continue;
5784cc7f7b2SMatthew G. Knepley             for (j = 0; j < numNCIndices; ++j) {
5794cc7f7b2SMatthew G. Knepley               key.j = ncindices[j] + rStart; /* global column (from Swarm) */
5804cc7f7b2SMatthew G. Knepley               if (key.j < 0) continue;
5819566063dSJacob Faibussowitsch               PetscCall(PetscHSetIJQueryAdd(ht, key, &missing));
5824cc7f7b2SMatthew G. Knepley               if (missing) {
5834cc7f7b2SMatthew G. Knepley                 if ((key.j >= colStart) && (key.j < colEnd)) ++dnz[key.i - rStart];
5844cc7f7b2SMatthew G. Knepley                 else                                         ++onz[key.i - rStart];
5854cc7f7b2SMatthew G. Knepley               }
5864cc7f7b2SMatthew G. Knepley             }
5874cc7f7b2SMatthew G. Knepley           }
5884cc7f7b2SMatthew G. Knepley         }
5899566063dSJacob Faibussowitsch         PetscCall(PetscFree(ncindices));
5904cc7f7b2SMatthew G. Knepley       }
5914cc7f7b2SMatthew G. Knepley     }
5924cc7f7b2SMatthew G. Knepley #endif
5939566063dSJacob Faibussowitsch     PetscCall(PetscFree(cindices));
5944cc7f7b2SMatthew G. Knepley   }
5959566063dSJacob Faibussowitsch   PetscCall(PetscHSetIJDestroy(&ht));
5969566063dSJacob Faibussowitsch   PetscCall(MatXAIJSetPreallocation(mass, 1, dnz, onz, NULL, NULL));
5979566063dSJacob Faibussowitsch   PetscCall(MatSetOption(mass, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
5989566063dSJacob Faibussowitsch   PetscCall(PetscFree2(dnz, onz));
5999566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(maxC * totDim, &elemMat, maxC * maxC, &elemMatSq, maxC, &rowIDXs, maxC * cdim, &xi));
6004cc7f7b2SMatthew G. Knepley   /* Fill in values
6014cc7f7b2SMatthew G. Knepley        Each entry is a sum of terms \phi_i(x_p) \phi_i(x_q)
6024cc7f7b2SMatthew G. Knepley        Start just by producing block diagonal
6034cc7f7b2SMatthew G. Knepley        Could loop over adjacent cells
6044cc7f7b2SMatthew G. Knepley          Produce neighboring element matrix
6054cc7f7b2SMatthew G. Knepley          TODO Determine which columns and rows correspond to shared dual vector
6064cc7f7b2SMatthew G. Knepley          Do MatMatMult with rectangular matrices
6074cc7f7b2SMatthew G. Knepley          Insert block
6084cc7f7b2SMatthew G. Knepley   */
6094cc7f7b2SMatthew G. Knepley   for (field = 0; field < Nf; ++field) {
6104cc7f7b2SMatthew G. Knepley     PetscTabulation Tcoarse;
6114cc7f7b2SMatthew G. Knepley     PetscObject     obj;
6124cc7f7b2SMatthew G. Knepley     PetscReal      *coords;
6134cc7f7b2SMatthew G. Knepley     PetscInt        Nc, i;
6144cc7f7b2SMatthew G. Knepley 
6159566063dSJacob Faibussowitsch     PetscCall(PetscDSGetDiscretization(prob, field, &obj));
6169566063dSJacob Faibussowitsch     PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc));
61763a3b9bcSJacob Faibussowitsch     PetscCheck(Nc == 1, PetscObjectComm((PetscObject)dmf), PETSC_ERR_SUP, "Can only interpolate a scalar field from particles, Nc = %" PetscInt_FMT, Nc);
6189566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
6194cc7f7b2SMatthew G. Knepley     for (cell = cStart; cell < cEnd; ++cell) {
6204cc7f7b2SMatthew G. Knepley       PetscInt *findices, *cindices;
6214cc7f7b2SMatthew G. Knepley       PetscInt  numFIndices, numCIndices;
6224cc7f7b2SMatthew G. Knepley       PetscInt  p, c;
6234cc7f7b2SMatthew G. Knepley 
6244cc7f7b2SMatthew G. Knepley       /* TODO: Use DMField instead of assuming affine */
6259566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ));
6269566063dSJacob Faibussowitsch       PetscCall(DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
6279566063dSJacob Faibussowitsch       PetscCall(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices));
628ad540459SPierre Jolivet       for (p = 0; p < numCIndices; ++p) CoordinatesRealToRef(cdim, cdim, v0ref, v0, invJ, &coords[cindices[p] * cdim], &xi[p * cdim]);
6299566063dSJacob Faibussowitsch       PetscCall(PetscFECreateTabulation((PetscFE)obj, 1, numCIndices, xi, 0, &Tcoarse));
6304cc7f7b2SMatthew G. Knepley       /* Get elemMat entries by multiplying by weight */
6319566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(elemMat, numCIndices * totDim));
6324cc7f7b2SMatthew G. Knepley       for (i = 0; i < numFIndices; ++i) {
6334cc7f7b2SMatthew G. Knepley         for (p = 0; p < numCIndices; ++p) {
6344cc7f7b2SMatthew G. Knepley           for (c = 0; c < Nc; ++c) {
6354cc7f7b2SMatthew G. Knepley             /* B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c */
6364cc7f7b2SMatthew G. Knepley             elemMat[p * numFIndices + i] += Tcoarse->T[0][(p * numFIndices + i) * Nc + c] * (useDeltaFunction ? 1.0 : detJ);
6374cc7f7b2SMatthew G. Knepley           }
6384cc7f7b2SMatthew G. Knepley         }
6394cc7f7b2SMatthew G. Knepley       }
6409566063dSJacob Faibussowitsch       PetscCall(PetscTabulationDestroy(&Tcoarse));
6414cc7f7b2SMatthew G. Knepley       for (p = 0; p < numCIndices; ++p) rowIDXs[p] = cindices[p] + rStart;
6429566063dSJacob Faibussowitsch       if (0) PetscCall(DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat));
6434cc7f7b2SMatthew G. Knepley       /* Block diagonal */
64478884ca7SMatthew G. Knepley       if (numCIndices) {
6454cc7f7b2SMatthew G. Knepley         PetscBLASInt blasn, blask;
6464cc7f7b2SMatthew G. Knepley         PetscScalar  one = 1.0, zero = 0.0;
6474cc7f7b2SMatthew G. Knepley 
6489566063dSJacob Faibussowitsch         PetscCall(PetscBLASIntCast(numCIndices, &blasn));
6499566063dSJacob Faibussowitsch         PetscCall(PetscBLASIntCast(numFIndices, &blask));
650792fecdfSBarry Smith         PetscCallBLAS("BLASgemm", BLASgemm_("T", "N", &blasn, &blasn, &blask, &one, elemMat, &blask, elemMat, &blask, &zero, elemMatSq, &blasn));
6514cc7f7b2SMatthew G. Knepley       }
6529566063dSJacob Faibussowitsch       PetscCall(MatSetValues(mass, numCIndices, rowIDXs, numCIndices, rowIDXs, elemMatSq, ADD_VALUES));
6534cf0e950SBarry Smith       /* TODO off-diagonal */
6549566063dSJacob Faibussowitsch       PetscCall(PetscFree(cindices));
6559566063dSJacob Faibussowitsch       PetscCall(DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
6564cc7f7b2SMatthew G. Knepley     }
6579566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
6584cc7f7b2SMatthew G. Knepley   }
6599566063dSJacob Faibussowitsch   PetscCall(PetscFree4(elemMat, elemMatSq, rowIDXs, xi));
6609566063dSJacob Faibussowitsch   PetscCall(PetscFree(adj));
6619566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortRestoreAccess(dmc));
6629566063dSJacob Faibussowitsch   PetscCall(PetscFree3(v0, J, invJ));
6639566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mass, MAT_FINAL_ASSEMBLY));
6649566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mass, MAT_FINAL_ASSEMBLY));
6653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6664cc7f7b2SMatthew G. Knepley }
6674cc7f7b2SMatthew G. Knepley 
6684cc7f7b2SMatthew G. Knepley /*@C
6694cc7f7b2SMatthew G. Knepley   DMSwarmCreateMassMatrixSquare - Creates the block-diagonal of the square, M^T_p M_p, of the particle mass matrix M_p
6704cc7f7b2SMatthew G. Knepley 
67120f4b53cSBarry Smith   Collective
6724cc7f7b2SMatthew G. Knepley 
67360225df5SJacob Faibussowitsch   Input Parameters:
67420f4b53cSBarry Smith + dmCoarse - a `DMSWARM`
67520f4b53cSBarry Smith - dmFine   - a `DMPLEX`
6764cc7f7b2SMatthew G. Knepley 
67760225df5SJacob Faibussowitsch   Output Parameter:
6784cc7f7b2SMatthew G. Knepley . mass - the square of the particle mass matrix
6794cc7f7b2SMatthew G. Knepley 
6804cc7f7b2SMatthew G. Knepley   Level: advanced
6814cc7f7b2SMatthew G. Knepley 
68220f4b53cSBarry Smith   Note:
6834cc7f7b2SMatthew G. Knepley   We only compute the block diagonal since this provides a good preconditioner and is completely local. It would be possible in the
6844cc7f7b2SMatthew G. Knepley   future to compute the full normal equations.
6854cc7f7b2SMatthew G. Knepley 
68620f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMCreateMassMatrix()`
6874cc7f7b2SMatthew G. Knepley @*/
688d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmCreateMassMatrixSquare(DM dmCoarse, DM dmFine, Mat *mass)
689d71ae5a4SJacob Faibussowitsch {
6904cc7f7b2SMatthew G. Knepley   PetscInt n;
6914cc7f7b2SMatthew G. Knepley   void    *ctx;
6924cc7f7b2SMatthew G. Knepley 
6934cc7f7b2SMatthew G. Knepley   PetscFunctionBegin;
6949566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetLocalSize(dmCoarse, &n));
6959566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)dmCoarse), mass));
6969566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*mass, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
6979566063dSJacob Faibussowitsch   PetscCall(MatSetType(*mass, dmCoarse->mattype));
6989566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(dmFine, &ctx));
6994cc7f7b2SMatthew G. Knepley 
7009566063dSJacob Faibussowitsch   PetscCall(DMSwarmComputeMassMatrixSquare_Private(dmCoarse, dmFine, *mass, PETSC_TRUE, ctx));
7019566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(*mass, NULL, "-mass_sq_mat_view"));
7023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7034cc7f7b2SMatthew G. Knepley }
7044cc7f7b2SMatthew G. Knepley 
705fb1bcc12SMatthew G. Knepley /*@C
70620f4b53cSBarry Smith   DMSwarmCreateGlobalVectorFromField - Creates a `Vec` object sharing the array associated with a given field
707d3a51819SDave May 
70820f4b53cSBarry Smith   Collective
709d3a51819SDave May 
71060225df5SJacob Faibussowitsch   Input Parameters:
71120f4b53cSBarry Smith + dm        - a `DMSWARM`
71262741f57SDave May - fieldname - the textual name given to a registered field
713d3a51819SDave May 
71460225df5SJacob Faibussowitsch   Output Parameter:
715d3a51819SDave May . vec - the vector
716d3a51819SDave May 
717d3a51819SDave May   Level: beginner
718d3a51819SDave May 
7198b8a3813SDave May   Notes:
72020f4b53cSBarry Smith   The vector must be returned using a matching call to `DMSwarmDestroyGlobalVectorFromField()`.
7218b8a3813SDave May 
72220f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmDestroyGlobalVectorFromField()`
723d3a51819SDave May @*/
724d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmCreateGlobalVectorFromField(DM dm, const char fieldname[], Vec *vec)
725d71ae5a4SJacob Faibussowitsch {
726fb1bcc12SMatthew G. Knepley   MPI_Comm comm = PetscObjectComm((PetscObject)dm);
727b5bcf523SDave May 
728fb1bcc12SMatthew G. Knepley   PetscFunctionBegin;
729ea7032a0SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
7309566063dSJacob Faibussowitsch   PetscCall(DMSwarmCreateVectorFromField_Private(dm, fieldname, comm, vec));
7313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
732b5bcf523SDave May }
733b5bcf523SDave May 
734d3a51819SDave May /*@C
73520f4b53cSBarry Smith   DMSwarmDestroyGlobalVectorFromField - Destroys the `Vec` object which share the array associated with a given field
736d3a51819SDave May 
73720f4b53cSBarry Smith   Collective
738d3a51819SDave May 
73960225df5SJacob Faibussowitsch   Input Parameters:
74020f4b53cSBarry Smith + dm        - a `DMSWARM`
74162741f57SDave May - fieldname - the textual name given to a registered field
742d3a51819SDave May 
74360225df5SJacob Faibussowitsch   Output Parameter:
744d3a51819SDave May . vec - the vector
745d3a51819SDave May 
746d3a51819SDave May   Level: beginner
747d3a51819SDave May 
74820f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmCreateGlobalVectorFromField()`
749d3a51819SDave May @*/
750d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDestroyGlobalVectorFromField(DM dm, const char fieldname[], Vec *vec)
751d71ae5a4SJacob Faibussowitsch {
752fb1bcc12SMatthew G. Knepley   PetscFunctionBegin;
753ea7032a0SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
7549566063dSJacob Faibussowitsch   PetscCall(DMSwarmDestroyVectorFromField_Private(dm, fieldname, vec));
7553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
756b5bcf523SDave May }
757b5bcf523SDave May 
758fb1bcc12SMatthew G. Knepley /*@C
75920f4b53cSBarry Smith   DMSwarmCreateLocalVectorFromField - Creates a `Vec` object sharing the array associated with a given field
760fb1bcc12SMatthew G. Knepley 
76120f4b53cSBarry Smith   Collective
762fb1bcc12SMatthew G. Knepley 
76360225df5SJacob Faibussowitsch   Input Parameters:
76420f4b53cSBarry Smith + dm        - a `DMSWARM`
76562741f57SDave May - fieldname - the textual name given to a registered field
766fb1bcc12SMatthew G. Knepley 
76760225df5SJacob Faibussowitsch   Output Parameter:
768fb1bcc12SMatthew G. Knepley . vec - the vector
769fb1bcc12SMatthew G. Knepley 
770fb1bcc12SMatthew G. Knepley   Level: beginner
771fb1bcc12SMatthew G. Knepley 
77220f4b53cSBarry Smith   Note:
7738b8a3813SDave May   The vector must be returned using a matching call to DMSwarmDestroyLocalVectorFromField().
7748b8a3813SDave May 
77520f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmDestroyLocalVectorFromField()`
776fb1bcc12SMatthew G. Knepley @*/
777d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmCreateLocalVectorFromField(DM dm, const char fieldname[], Vec *vec)
778d71ae5a4SJacob Faibussowitsch {
779fb1bcc12SMatthew G. Knepley   MPI_Comm comm = PETSC_COMM_SELF;
780bbe8250bSMatthew G. Knepley 
781fb1bcc12SMatthew G. Knepley   PetscFunctionBegin;
7829566063dSJacob Faibussowitsch   PetscCall(DMSwarmCreateVectorFromField_Private(dm, fieldname, comm, vec));
7833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
784bbe8250bSMatthew G. Knepley }
785fb1bcc12SMatthew G. Knepley 
786fb1bcc12SMatthew G. Knepley /*@C
78720f4b53cSBarry Smith   DMSwarmDestroyLocalVectorFromField - Destroys the `Vec` object which share the array associated with a given field
788fb1bcc12SMatthew G. Knepley 
78920f4b53cSBarry Smith   Collective
790fb1bcc12SMatthew G. Knepley 
79160225df5SJacob Faibussowitsch   Input Parameters:
79220f4b53cSBarry Smith + dm        - a `DMSWARM`
79362741f57SDave May - fieldname - the textual name given to a registered field
794fb1bcc12SMatthew G. Knepley 
79560225df5SJacob Faibussowitsch   Output Parameter:
796fb1bcc12SMatthew G. Knepley . vec - the vector
797fb1bcc12SMatthew G. Knepley 
798fb1bcc12SMatthew G. Knepley   Level: beginner
799fb1bcc12SMatthew G. Knepley 
80020f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmCreateLocalVectorFromField()`
801fb1bcc12SMatthew G. Knepley @*/
802d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDestroyLocalVectorFromField(DM dm, const char fieldname[], Vec *vec)
803d71ae5a4SJacob Faibussowitsch {
804fb1bcc12SMatthew G. Knepley   PetscFunctionBegin;
805ea7032a0SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
8069566063dSJacob Faibussowitsch   PetscCall(DMSwarmDestroyVectorFromField_Private(dm, fieldname, vec));
8073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
808bbe8250bSMatthew G. Knepley }
809bbe8250bSMatthew G. Knepley 
810d3a51819SDave May /*@C
81120f4b53cSBarry Smith   DMSwarmInitializeFieldRegister - Initiates the registration of fields to a `DMSWARM`
812d3a51819SDave May 
81320f4b53cSBarry Smith   Collective
814d3a51819SDave May 
81560225df5SJacob Faibussowitsch   Input Parameter:
81620f4b53cSBarry Smith . dm - a `DMSWARM`
817d3a51819SDave May 
818d3a51819SDave May   Level: beginner
819d3a51819SDave May 
82020f4b53cSBarry Smith   Note:
82120f4b53cSBarry Smith   After all fields have been registered, you must call `DMSwarmFinalizeFieldRegister()`.
822d3a51819SDave May 
82320f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmFinalizeFieldRegister()`, `DMSwarmRegisterPetscDatatypeField()`,
824db781477SPatrick Sanan           `DMSwarmRegisterUserStructField()`, `DMSwarmRegisterUserDatatypeField()`
825d3a51819SDave May @*/
826d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmInitializeFieldRegister(DM dm)
827d71ae5a4SJacob Faibussowitsch {
8285f50eb2eSDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
8293454631fSDave May 
830521f74f9SMatthew G. Knepley   PetscFunctionBegin;
831cc651181SDave May   if (!swarm->field_registration_initialized) {
8325f50eb2eSDave May     swarm->field_registration_initialized = PETSC_TRUE;
833da81f932SPierre Jolivet     PetscCall(DMSwarmRegisterPetscDatatypeField(dm, DMSwarmField_pid, 1, PETSC_INT64)); /* unique identifier */
8349566063dSJacob Faibussowitsch     PetscCall(DMSwarmRegisterPetscDatatypeField(dm, DMSwarmField_rank, 1, PETSC_INT));  /* used for communication */
835cc651181SDave May   }
8363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8375f50eb2eSDave May }
8385f50eb2eSDave May 
83974d0cae8SMatthew G. Knepley /*@
84020f4b53cSBarry Smith   DMSwarmFinalizeFieldRegister - Finalizes the registration of fields to a `DMSWARM`
841d3a51819SDave May 
84220f4b53cSBarry Smith   Collective
843d3a51819SDave May 
84460225df5SJacob Faibussowitsch   Input Parameter:
84520f4b53cSBarry Smith . dm - a `DMSWARM`
846d3a51819SDave May 
847d3a51819SDave May   Level: beginner
848d3a51819SDave May 
84920f4b53cSBarry Smith   Note:
85020f4b53cSBarry Smith   After `DMSwarmFinalizeFieldRegister()` has been called, no new fields can be defined on the `DMSWARM`.
851d3a51819SDave May 
85220f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmInitializeFieldRegister()`, `DMSwarmRegisterPetscDatatypeField()`,
853db781477SPatrick Sanan           `DMSwarmRegisterUserStructField()`, `DMSwarmRegisterUserDatatypeField()`
854d3a51819SDave May @*/
855d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmFinalizeFieldRegister(DM dm)
856d71ae5a4SJacob Faibussowitsch {
8575f50eb2eSDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
8586845f8f5SDave May 
859521f74f9SMatthew G. Knepley   PetscFunctionBegin;
86048a46eb9SPierre Jolivet   if (!swarm->field_registration_finalized) PetscCall(DMSwarmDataBucketFinalize(swarm->db));
861f0cdbbbaSDave May   swarm->field_registration_finalized = PETSC_TRUE;
8623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8635f50eb2eSDave May }
8645f50eb2eSDave May 
86574d0cae8SMatthew G. Knepley /*@
86620f4b53cSBarry Smith   DMSwarmSetLocalSizes - Sets the length of all registered fields on the `DMSWARM`
867d3a51819SDave May 
86820f4b53cSBarry Smith   Not Collective
869d3a51819SDave May 
87060225df5SJacob Faibussowitsch   Input Parameters:
87120f4b53cSBarry Smith + dm     - a `DMSWARM`
872d3a51819SDave May . nlocal - the length of each registered field
87362741f57SDave May - buffer - the length of the buffer used to efficient dynamic re-sizing
874d3a51819SDave May 
875d3a51819SDave May   Level: beginner
876d3a51819SDave May 
87720f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmGetLocalSize()`
878d3a51819SDave May @*/
879d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmSetLocalSizes(DM dm, PetscInt nlocal, PetscInt buffer)
880d71ae5a4SJacob Faibussowitsch {
8815f50eb2eSDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
8825f50eb2eSDave May 
883521f74f9SMatthew G. Knepley   PetscFunctionBegin;
8849566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMSWARM_SetSizes, 0, 0, 0, 0));
8859566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketSetSizes(swarm->db, nlocal, buffer));
8869566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMSWARM_SetSizes, 0, 0, 0, 0));
8873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8885f50eb2eSDave May }
8895f50eb2eSDave May 
89074d0cae8SMatthew G. Knepley /*@
89120f4b53cSBarry Smith   DMSwarmSetCellDM - Attaches a `DM` to a `DMSWARM`
892d3a51819SDave May 
89320f4b53cSBarry Smith   Collective
894d3a51819SDave May 
89560225df5SJacob Faibussowitsch   Input Parameters:
89620f4b53cSBarry Smith + dm     - a `DMSWARM`
89720f4b53cSBarry Smith - dmcell - the `DM` to attach to the `DMSWARM`
898d3a51819SDave May 
899d3a51819SDave May   Level: beginner
900d3a51819SDave May 
90120f4b53cSBarry Smith   Note:
90220f4b53cSBarry Smith   The attached `DM` (dmcell) will be queried for point location and
90320f4b53cSBarry Smith   neighbor MPI-rank information if `DMSwarmMigrate()` is called.
904d3a51819SDave May 
90520f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmSetType()`, `DMSwarmGetCellDM()`, `DMSwarmMigrate()`
906d3a51819SDave May @*/
907d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmSetCellDM(DM dm, DM dmcell)
908d71ae5a4SJacob Faibussowitsch {
909b16650c8SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
910521f74f9SMatthew G. Knepley 
911521f74f9SMatthew G. Knepley   PetscFunctionBegin;
912b16650c8SDave May   swarm->dmcell = dmcell;
9133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
914b16650c8SDave May }
915b16650c8SDave May 
91674d0cae8SMatthew G. Knepley /*@
91720f4b53cSBarry Smith   DMSwarmGetCellDM - Fetches the attached cell `DM`
918d3a51819SDave May 
91920f4b53cSBarry Smith   Collective
920d3a51819SDave May 
92160225df5SJacob Faibussowitsch   Input Parameter:
92220f4b53cSBarry Smith . dm - a `DMSWARM`
923d3a51819SDave May 
92460225df5SJacob Faibussowitsch   Output Parameter:
92520f4b53cSBarry Smith . dmcell - the `DM` which was attached to the `DMSWARM`
926d3a51819SDave May 
927d3a51819SDave May   Level: beginner
928d3a51819SDave May 
92920f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmSetCellDM()`
930d3a51819SDave May @*/
931d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmGetCellDM(DM dm, DM *dmcell)
932d71ae5a4SJacob Faibussowitsch {
933fe39f135SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
934521f74f9SMatthew G. Knepley 
935521f74f9SMatthew G. Knepley   PetscFunctionBegin;
936fe39f135SDave May   *dmcell = swarm->dmcell;
9373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
938fe39f135SDave May }
939fe39f135SDave May 
94074d0cae8SMatthew G. Knepley /*@
941d5b43468SJose E. Roman   DMSwarmGetLocalSize - Retrieves the local length of fields registered
942d3a51819SDave May 
94320f4b53cSBarry Smith   Not Collective
944d3a51819SDave May 
94560225df5SJacob Faibussowitsch   Input Parameter:
94620f4b53cSBarry Smith . dm - a `DMSWARM`
947d3a51819SDave May 
94860225df5SJacob Faibussowitsch   Output Parameter:
949d3a51819SDave May . nlocal - the length of each registered field
950d3a51819SDave May 
951d3a51819SDave May   Level: beginner
952d3a51819SDave May 
95320f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmGetSize()`, `DMSwarmSetLocalSizes()`
954d3a51819SDave May @*/
955d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmGetLocalSize(DM dm, PetscInt *nlocal)
956d71ae5a4SJacob Faibussowitsch {
957dcf43ee8SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
958dcf43ee8SDave May 
959521f74f9SMatthew G. Knepley   PetscFunctionBegin;
9609566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, nlocal, NULL, NULL));
9613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
962dcf43ee8SDave May }
963dcf43ee8SDave May 
96474d0cae8SMatthew G. Knepley /*@
965d5b43468SJose E. Roman   DMSwarmGetSize - Retrieves the total length of fields registered
966d3a51819SDave May 
96720f4b53cSBarry Smith   Collective
968d3a51819SDave May 
96960225df5SJacob Faibussowitsch   Input Parameter:
97020f4b53cSBarry Smith . dm - a `DMSWARM`
971d3a51819SDave May 
97260225df5SJacob Faibussowitsch   Output Parameter:
973d3a51819SDave May . n - the total length of each registered field
974d3a51819SDave May 
975d3a51819SDave May   Level: beginner
976d3a51819SDave May 
977d3a51819SDave May   Note:
97820f4b53cSBarry Smith   This calls `MPI_Allreduce()` upon each call (inefficient but safe)
979d3a51819SDave May 
98020f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmGetLocalSize()`, `DMSwarmSetLocalSizes()`
981d3a51819SDave May @*/
982d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmGetSize(DM dm, PetscInt *n)
983d71ae5a4SJacob Faibussowitsch {
984dcf43ee8SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
9855627991aSBarry Smith   PetscInt  nlocal;
986dcf43ee8SDave May 
987521f74f9SMatthew G. Knepley   PetscFunctionBegin;
9889566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &nlocal, NULL, NULL));
989712fec58SPierre Jolivet   PetscCall(MPIU_Allreduce(&nlocal, n, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm)));
9903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
991dcf43ee8SDave May }
992dcf43ee8SDave May 
993d3a51819SDave May /*@C
99420f4b53cSBarry Smith   DMSwarmRegisterPetscDatatypeField - Register a field to a `DMSWARM` with a native PETSc data type
995d3a51819SDave May 
99620f4b53cSBarry Smith   Collective
997d3a51819SDave May 
99860225df5SJacob Faibussowitsch   Input Parameters:
99920f4b53cSBarry Smith + dm        - a `DMSWARM`
1000d3a51819SDave May . fieldname - the textual name to identify this field
1001d3a51819SDave May . blocksize - the number of each data type
100220f4b53cSBarry Smith - type      - a valid PETSc data type (`PETSC_CHAR`, `PETSC_SHORT`, `PETSC_INT`, `PETSC_FLOAT`, `PETSC_REAL`, `PETSC_LONG`)
1003d3a51819SDave May 
1004d3a51819SDave May   Level: beginner
1005d3a51819SDave May 
1006d3a51819SDave May   Notes:
10078b8a3813SDave May   The textual name for each registered field must be unique.
1008d3a51819SDave May 
100920f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterUserStructField()`, `DMSwarmRegisterUserDatatypeField()`
1010d3a51819SDave May @*/
1011d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmRegisterPetscDatatypeField(DM dm, const char fieldname[], PetscInt blocksize, PetscDataType type)
1012d71ae5a4SJacob Faibussowitsch {
1013b62e03f8SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1014b62e03f8SDave May   size_t    size;
1015b62e03f8SDave May 
1016521f74f9SMatthew G. Knepley   PetscFunctionBegin;
101728b400f6SJacob Faibussowitsch   PetscCheck(swarm->field_registration_initialized, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Must call DMSwarmInitializeFieldRegister() first");
101828b400f6SJacob Faibussowitsch   PetscCheck(!swarm->field_registration_finalized, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Cannot register additional fields after calling DMSwarmFinalizeFieldRegister() first");
10195f50eb2eSDave May 
102008401ef6SPierre Jolivet   PetscCheck(type != PETSC_OBJECT, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}");
102108401ef6SPierre Jolivet   PetscCheck(type != PETSC_FUNCTION, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}");
102208401ef6SPierre Jolivet   PetscCheck(type != PETSC_STRING, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}");
102308401ef6SPierre Jolivet   PetscCheck(type != PETSC_STRUCT, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}");
102408401ef6SPierre Jolivet   PetscCheck(type != PETSC_DATATYPE_UNKNOWN, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}");
1025b62e03f8SDave May 
10269566063dSJacob Faibussowitsch   PetscCall(PetscDataTypeGetSize(type, &size));
1027b62e03f8SDave May   /* Load a specific data type into data bucket, specifying textual name and its size in bytes */
10289566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketRegisterField(swarm->db, "DMSwarmRegisterPetscDatatypeField", fieldname, blocksize * size, NULL));
102952c3ed93SDave May   {
103077048351SPatrick Sanan     DMSwarmDataField gfield;
103152c3ed93SDave May 
10329566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield));
10339566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataFieldSetBlockSize(gfield, blocksize));
103452c3ed93SDave May   }
1035b62e03f8SDave May   swarm->db->field[swarm->db->nfields - 1]->petsc_type = type;
10363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1037b62e03f8SDave May }
1038b62e03f8SDave May 
1039d3a51819SDave May /*@C
104020f4b53cSBarry Smith   DMSwarmRegisterUserStructField - Register a user defined struct to a `DMSWARM`
1041d3a51819SDave May 
104220f4b53cSBarry Smith   Collective
1043d3a51819SDave May 
104460225df5SJacob Faibussowitsch   Input Parameters:
104520f4b53cSBarry Smith + dm        - a `DMSWARM`
1046d3a51819SDave May . fieldname - the textual name to identify this field
104762741f57SDave May - size      - the size in bytes of the user struct of each data type
1048d3a51819SDave May 
1049d3a51819SDave May   Level: beginner
1050d3a51819SDave May 
105120f4b53cSBarry Smith   Note:
10528b8a3813SDave May   The textual name for each registered field must be unique.
1053d3a51819SDave May 
105420f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmRegisterUserDatatypeField()`
1055d3a51819SDave May @*/
1056d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmRegisterUserStructField(DM dm, const char fieldname[], size_t size)
1057d71ae5a4SJacob Faibussowitsch {
1058b62e03f8SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1059b62e03f8SDave May 
1060521f74f9SMatthew G. Knepley   PetscFunctionBegin;
10619566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketRegisterField(swarm->db, "DMSwarmRegisterUserStructField", fieldname, size, NULL));
1062b62e03f8SDave May   swarm->db->field[swarm->db->nfields - 1]->petsc_type = PETSC_STRUCT;
10633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1064b62e03f8SDave May }
1065b62e03f8SDave May 
1066d3a51819SDave May /*@C
106720f4b53cSBarry Smith   DMSwarmRegisterUserDatatypeField - Register a user defined data type to a `DMSWARM`
1068d3a51819SDave May 
106920f4b53cSBarry Smith   Collective
1070d3a51819SDave May 
107160225df5SJacob Faibussowitsch   Input Parameters:
107220f4b53cSBarry Smith + dm        - a `DMSWARM`
1073d3a51819SDave May . fieldname - the textual name to identify this field
1074d3a51819SDave May . size      - the size in bytes of the user data type
107562741f57SDave May - blocksize - the number of each data type
1076d3a51819SDave May 
1077d3a51819SDave May   Level: beginner
1078d3a51819SDave May 
107920f4b53cSBarry Smith   Note:
10808b8a3813SDave May   The textual name for each registered field must be unique.
1081d3a51819SDave May 
108242747ad1SJacob Faibussowitsch .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmRegisterUserStructField()`
1083d3a51819SDave May @*/
1084d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmRegisterUserDatatypeField(DM dm, const char fieldname[], size_t size, PetscInt blocksize)
1085d71ae5a4SJacob Faibussowitsch {
1086b62e03f8SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1087b62e03f8SDave May 
1088521f74f9SMatthew G. Knepley   PetscFunctionBegin;
10899566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketRegisterField(swarm->db, "DMSwarmRegisterUserDatatypeField", fieldname, blocksize * size, NULL));
1090320740a0SDave May   {
109177048351SPatrick Sanan     DMSwarmDataField gfield;
1092320740a0SDave May 
10939566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield));
10949566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataFieldSetBlockSize(gfield, blocksize));
1095320740a0SDave May   }
1096b62e03f8SDave May   swarm->db->field[swarm->db->nfields - 1]->petsc_type = PETSC_DATATYPE_UNKNOWN;
10973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1098b62e03f8SDave May }
1099b62e03f8SDave May 
1100d3a51819SDave May /*@C
1101d3a51819SDave May   DMSwarmGetField - Get access to the underlying array storing all entries associated with a registered field
1102d3a51819SDave May 
110320f4b53cSBarry Smith   Not Collective
1104d3a51819SDave May 
110560225df5SJacob Faibussowitsch   Input Parameters:
110620f4b53cSBarry Smith + dm        - a `DMSWARM`
110762741f57SDave May - fieldname - the textual name to identify this field
1108d3a51819SDave May 
110960225df5SJacob Faibussowitsch   Output Parameters:
111062741f57SDave May + blocksize - the number of each data type
1111d3a51819SDave May . type      - the data type
111262741f57SDave May - data      - pointer to raw array
1113d3a51819SDave May 
1114d3a51819SDave May   Level: beginner
1115d3a51819SDave May 
1116d3a51819SDave May   Notes:
111720f4b53cSBarry Smith   The array must be returned using a matching call to `DMSwarmRestoreField()`.
1118d3a51819SDave May 
111920f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRestoreField()`
1120d3a51819SDave May @*/
1121d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmGetField(DM dm, const char fieldname[], PetscInt *blocksize, PetscDataType *type, void **data)
1122d71ae5a4SJacob Faibussowitsch {
1123b62e03f8SDave May   DM_Swarm        *swarm = (DM_Swarm *)dm->data;
112477048351SPatrick Sanan   DMSwarmDataField gfield;
1125b62e03f8SDave May 
1126521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1127ea7032a0SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
11289566063dSJacob Faibussowitsch   if (!swarm->issetup) PetscCall(DMSetUp(dm));
11299566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield));
11309566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataFieldGetAccess(gfield));
11319566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataFieldGetEntries(gfield, data));
1132ad540459SPierre Jolivet   if (blocksize) *blocksize = gfield->bs;
1133ad540459SPierre Jolivet   if (type) *type = gfield->petsc_type;
11343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1135b62e03f8SDave May }
1136b62e03f8SDave May 
1137d3a51819SDave May /*@C
1138d3a51819SDave May   DMSwarmRestoreField - Restore access to the underlying array storing all entries associated with a registered field
1139d3a51819SDave May 
114020f4b53cSBarry Smith   Not Collective
1141d3a51819SDave May 
114260225df5SJacob Faibussowitsch   Input Parameters:
114320f4b53cSBarry Smith + dm        - a `DMSWARM`
114462741f57SDave May - fieldname - the textual name to identify this field
1145d3a51819SDave May 
114660225df5SJacob Faibussowitsch   Output Parameters:
114762741f57SDave May + blocksize - the number of each data type
1148d3a51819SDave May . type      - the data type
114962741f57SDave May - data      - pointer to raw array
1150d3a51819SDave May 
1151d3a51819SDave May   Level: beginner
1152d3a51819SDave May 
1153d3a51819SDave May   Notes:
115420f4b53cSBarry Smith   The user must call `DMSwarmGetField()` prior to calling `DMSwarmRestoreField()`.
1155d3a51819SDave May 
115620f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmGetField()`
1157d3a51819SDave May @*/
1158d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmRestoreField(DM dm, const char fieldname[], PetscInt *blocksize, PetscDataType *type, void **data)
1159d71ae5a4SJacob Faibussowitsch {
1160b62e03f8SDave May   DM_Swarm        *swarm = (DM_Swarm *)dm->data;
116177048351SPatrick Sanan   DMSwarmDataField gfield;
1162b62e03f8SDave May 
1163521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1164ea7032a0SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
11659566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield));
11669566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataFieldRestoreAccess(gfield));
1167b62e03f8SDave May   if (data) *data = NULL;
11683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1169b62e03f8SDave May }
1170b62e03f8SDave May 
1171*0bf7c1c5SMatthew G. Knepley PetscErrorCode DMSwarmGetFieldInfo(DM dm, const char fieldname[], PetscInt *blocksize, PetscDataType *type)
1172*0bf7c1c5SMatthew G. Knepley {
1173*0bf7c1c5SMatthew G. Knepley   DM_Swarm        *swarm = (DM_Swarm *)dm->data;
1174*0bf7c1c5SMatthew G. Knepley   DMSwarmDataField gfield;
1175*0bf7c1c5SMatthew G. Knepley 
1176*0bf7c1c5SMatthew G. Knepley   PetscFunctionBegin;
1177*0bf7c1c5SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
1178*0bf7c1c5SMatthew G. Knepley   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield));
1179*0bf7c1c5SMatthew G. Knepley   if (blocksize) *blocksize = gfield->bs;
1180*0bf7c1c5SMatthew G. Knepley   if (type) *type = gfield->petsc_type;
1181*0bf7c1c5SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
1182*0bf7c1c5SMatthew G. Knepley }
1183*0bf7c1c5SMatthew G. Knepley 
118474d0cae8SMatthew G. Knepley /*@
118520f4b53cSBarry Smith   DMSwarmAddPoint - Add space for one new point in the `DMSWARM`
1186d3a51819SDave May 
118720f4b53cSBarry Smith   Not Collective
1188d3a51819SDave May 
118960225df5SJacob Faibussowitsch   Input Parameter:
119020f4b53cSBarry Smith . dm - a `DMSWARM`
1191d3a51819SDave May 
1192d3a51819SDave May   Level: beginner
1193d3a51819SDave May 
1194d3a51819SDave May   Notes:
11958b8a3813SDave May   The new point will have all fields initialized to zero.
1196d3a51819SDave May 
119720f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmAddNPoints()`
1198d3a51819SDave May @*/
1199d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmAddPoint(DM dm)
1200d71ae5a4SJacob Faibussowitsch {
1201cb1d1399SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1202cb1d1399SDave May 
1203521f74f9SMatthew G. Knepley   PetscFunctionBegin;
12049566063dSJacob Faibussowitsch   if (!swarm->issetup) PetscCall(DMSetUp(dm));
12059566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMSWARM_AddPoints, 0, 0, 0, 0));
12069566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketAddPoint(swarm->db));
12079566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMSWARM_AddPoints, 0, 0, 0, 0));
12083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1209cb1d1399SDave May }
1210cb1d1399SDave May 
121174d0cae8SMatthew G. Knepley /*@
121220f4b53cSBarry Smith   DMSwarmAddNPoints - Add space for a number of new points in the `DMSWARM`
1213d3a51819SDave May 
121420f4b53cSBarry Smith   Not Collective
1215d3a51819SDave May 
121660225df5SJacob Faibussowitsch   Input Parameters:
121720f4b53cSBarry Smith + dm      - a `DMSWARM`
121862741f57SDave May - npoints - the number of new points to add
1219d3a51819SDave May 
1220d3a51819SDave May   Level: beginner
1221d3a51819SDave May 
1222d3a51819SDave May   Notes:
12238b8a3813SDave May   The new point will have all fields initialized to zero.
1224d3a51819SDave May 
122520f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmAddPoint()`
1226d3a51819SDave May @*/
1227d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmAddNPoints(DM dm, PetscInt npoints)
1228d71ae5a4SJacob Faibussowitsch {
1229cb1d1399SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1230cb1d1399SDave May   PetscInt  nlocal;
1231cb1d1399SDave May 
1232521f74f9SMatthew G. Knepley   PetscFunctionBegin;
12339566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMSWARM_AddPoints, 0, 0, 0, 0));
12349566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &nlocal, NULL, NULL));
1235cb1d1399SDave May   nlocal = nlocal + npoints;
12369566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketSetSizes(swarm->db, nlocal, DMSWARM_DATA_BUCKET_BUFFER_DEFAULT));
12379566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMSWARM_AddPoints, 0, 0, 0, 0));
12383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1239cb1d1399SDave May }
1240cb1d1399SDave May 
124174d0cae8SMatthew G. Knepley /*@
124220f4b53cSBarry Smith   DMSwarmRemovePoint - Remove the last point from the `DMSWARM`
1243d3a51819SDave May 
124420f4b53cSBarry Smith   Not Collective
1245d3a51819SDave May 
124660225df5SJacob Faibussowitsch   Input Parameter:
124720f4b53cSBarry Smith . dm - a `DMSWARM`
1248d3a51819SDave May 
1249d3a51819SDave May   Level: beginner
1250d3a51819SDave May 
125120f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRemovePointAtIndex()`
1252d3a51819SDave May @*/
1253d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmRemovePoint(DM dm)
1254d71ae5a4SJacob Faibussowitsch {
1255cb1d1399SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1256cb1d1399SDave May 
1257521f74f9SMatthew G. Knepley   PetscFunctionBegin;
12589566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMSWARM_RemovePoints, 0, 0, 0, 0));
12599566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketRemovePoint(swarm->db));
12609566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMSWARM_RemovePoints, 0, 0, 0, 0));
12613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1262cb1d1399SDave May }
1263cb1d1399SDave May 
126474d0cae8SMatthew G. Knepley /*@
126520f4b53cSBarry Smith   DMSwarmRemovePointAtIndex - Removes a specific point from the `DMSWARM`
1266d3a51819SDave May 
126720f4b53cSBarry Smith   Not Collective
1268d3a51819SDave May 
126960225df5SJacob Faibussowitsch   Input Parameters:
127020f4b53cSBarry Smith + dm  - a `DMSWARM`
127162741f57SDave May - idx - index of point to remove
1272d3a51819SDave May 
1273d3a51819SDave May   Level: beginner
1274d3a51819SDave May 
127520f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRemovePoint()`
1276d3a51819SDave May @*/
1277d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmRemovePointAtIndex(DM dm, PetscInt idx)
1278d71ae5a4SJacob Faibussowitsch {
1279cb1d1399SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1280cb1d1399SDave May 
1281521f74f9SMatthew G. Knepley   PetscFunctionBegin;
12829566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMSWARM_RemovePoints, 0, 0, 0, 0));
12839566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketRemovePointAtIndex(swarm->db, idx));
12849566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMSWARM_RemovePoints, 0, 0, 0, 0));
12853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1286cb1d1399SDave May }
1287b62e03f8SDave May 
128874d0cae8SMatthew G. Knepley /*@
128920f4b53cSBarry Smith   DMSwarmCopyPoint - Copy point pj to point pi in the `DMSWARM`
1290ba4fc9c6SDave May 
129120f4b53cSBarry Smith   Not Collective
1292ba4fc9c6SDave May 
129360225df5SJacob Faibussowitsch   Input Parameters:
129420f4b53cSBarry Smith + dm - a `DMSWARM`
1295ba4fc9c6SDave May . pi - the index of the point to copy
1296ba4fc9c6SDave May - pj - the point index where the copy should be located
1297ba4fc9c6SDave May 
1298ba4fc9c6SDave May   Level: beginner
1299ba4fc9c6SDave May 
130020f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRemovePoint()`
1301ba4fc9c6SDave May @*/
1302d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmCopyPoint(DM dm, PetscInt pi, PetscInt pj)
1303d71ae5a4SJacob Faibussowitsch {
1304ba4fc9c6SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1305ba4fc9c6SDave May 
1306ba4fc9c6SDave May   PetscFunctionBegin;
13079566063dSJacob Faibussowitsch   if (!swarm->issetup) PetscCall(DMSetUp(dm));
13089566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketCopyPoint(swarm->db, pi, swarm->db, pj));
13093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1310ba4fc9c6SDave May }
1311ba4fc9c6SDave May 
131266976f2fSJacob Faibussowitsch static PetscErrorCode DMSwarmMigrate_Basic(DM dm, PetscBool remove_sent_points)
1313d71ae5a4SJacob Faibussowitsch {
1314521f74f9SMatthew G. Knepley   PetscFunctionBegin;
13159566063dSJacob Faibussowitsch   PetscCall(DMSwarmMigrate_Push_Basic(dm, remove_sent_points));
13163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13173454631fSDave May }
13183454631fSDave May 
131974d0cae8SMatthew G. Knepley /*@
132020f4b53cSBarry Smith   DMSwarmMigrate - Relocates points defined in the `DMSWARM` to other MPI-ranks
1321d3a51819SDave May 
132220f4b53cSBarry Smith   Collective
1323d3a51819SDave May 
132460225df5SJacob Faibussowitsch   Input Parameters:
132520f4b53cSBarry Smith + dm                 - the `DMSWARM`
132662741f57SDave May - remove_sent_points - flag indicating if sent points should be removed from the current MPI-rank
1327d3a51819SDave May 
1328d3a51819SDave May   Level: advanced
1329d3a51819SDave May 
133020f4b53cSBarry Smith   Notes:
133120f4b53cSBarry Smith   The `DM` will be modified to accommodate received points.
133220f4b53cSBarry Smith   If `remove_sent_points` is `PETSC_TRUE`, any points that were sent will be removed from the `DM`.
133320f4b53cSBarry Smith   Different styles of migration are supported. See `DMSwarmSetMigrateType()`.
133420f4b53cSBarry Smith 
133520f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmSetMigrateType()`
1336d3a51819SDave May @*/
1337d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmMigrate(DM dm, PetscBool remove_sent_points)
1338d71ae5a4SJacob Faibussowitsch {
1339f0cdbbbaSDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1340f0cdbbbaSDave May 
1341521f74f9SMatthew G. Knepley   PetscFunctionBegin;
13429566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMSWARM_Migrate, 0, 0, 0, 0));
1343f0cdbbbaSDave May   switch (swarm->migrate_type) {
1344d71ae5a4SJacob Faibussowitsch   case DMSWARM_MIGRATE_BASIC:
1345d71ae5a4SJacob Faibussowitsch     PetscCall(DMSwarmMigrate_Basic(dm, remove_sent_points));
1346d71ae5a4SJacob Faibussowitsch     break;
1347d71ae5a4SJacob Faibussowitsch   case DMSWARM_MIGRATE_DMCELLNSCATTER:
1348d71ae5a4SJacob Faibussowitsch     PetscCall(DMSwarmMigrate_CellDMScatter(dm, remove_sent_points));
1349d71ae5a4SJacob Faibussowitsch     break;
1350d71ae5a4SJacob Faibussowitsch   case DMSWARM_MIGRATE_DMCELLEXACT:
1351d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_MIGRATE_DMCELLEXACT not implemented");
1352d71ae5a4SJacob Faibussowitsch   case DMSWARM_MIGRATE_USER:
1353d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_MIGRATE_USER not implemented");
1354d71ae5a4SJacob Faibussowitsch   default:
1355d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_MIGRATE type unknown");
1356f0cdbbbaSDave May   }
13579566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMSWARM_Migrate, 0, 0, 0, 0));
13589566063dSJacob Faibussowitsch   PetscCall(DMClearGlobalVectors(dm));
13593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13603454631fSDave May }
13613454631fSDave May 
1362f0cdbbbaSDave May PetscErrorCode DMSwarmMigrate_GlobalToLocal_Basic(DM dm, PetscInt *globalsize);
1363f0cdbbbaSDave May 
1364d3a51819SDave May /*
1365d3a51819SDave May  DMSwarmCollectViewCreate
1366d3a51819SDave May 
1367d3a51819SDave May  * Applies a collection method and gathers point neighbour points into dm
1368d3a51819SDave May 
1369d3a51819SDave May  Notes:
13708b8a3813SDave May  Users should call DMSwarmCollectViewDestroy() after
1371d3a51819SDave May  they have finished computations associated with the collected points
1372d3a51819SDave May */
1373d3a51819SDave May 
137474d0cae8SMatthew G. Knepley /*@
1375d3a51819SDave May   DMSwarmCollectViewCreate - Applies a collection method and gathers points
137620f4b53cSBarry Smith   in neighbour ranks into the `DMSWARM`
1377d3a51819SDave May 
137820f4b53cSBarry Smith   Collective
1379d3a51819SDave May 
138060225df5SJacob Faibussowitsch   Input Parameter:
138120f4b53cSBarry Smith . dm - the `DMSWARM`
1382d3a51819SDave May 
1383d3a51819SDave May   Level: advanced
1384d3a51819SDave May 
138520f4b53cSBarry Smith   Notes:
138620f4b53cSBarry Smith   Users should call `DMSwarmCollectViewDestroy()` after
138720f4b53cSBarry Smith   they have finished computations associated with the collected points
138820f4b53cSBarry Smith   Different collect methods are supported. See `DMSwarmSetCollectType()`.
138920f4b53cSBarry Smith 
139020f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmCollectViewDestroy()`, `DMSwarmSetCollectType()`
1391d3a51819SDave May @*/
1392d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmCollectViewCreate(DM dm)
1393d71ae5a4SJacob Faibussowitsch {
13942712d1f2SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
13952712d1f2SDave May   PetscInt  ng;
13962712d1f2SDave May 
1397521f74f9SMatthew G. Knepley   PetscFunctionBegin;
139828b400f6SJacob Faibussowitsch   PetscCheck(!swarm->collect_view_active, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "CollectView currently active");
13999566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetLocalSize(dm, &ng));
1400480eef7bSDave May   switch (swarm->collect_type) {
1401d71ae5a4SJacob Faibussowitsch   case DMSWARM_COLLECT_BASIC:
1402d71ae5a4SJacob Faibussowitsch     PetscCall(DMSwarmMigrate_GlobalToLocal_Basic(dm, &ng));
1403d71ae5a4SJacob Faibussowitsch     break;
1404d71ae5a4SJacob Faibussowitsch   case DMSWARM_COLLECT_DMDABOUNDINGBOX:
1405d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_COLLECT_DMDABOUNDINGBOX not implemented");
1406d71ae5a4SJacob Faibussowitsch   case DMSWARM_COLLECT_GENERAL:
1407d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_COLLECT_GENERAL not implemented");
1408d71ae5a4SJacob Faibussowitsch   default:
1409d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_COLLECT type unknown");
1410480eef7bSDave May   }
1411480eef7bSDave May   swarm->collect_view_active       = PETSC_TRUE;
1412480eef7bSDave May   swarm->collect_view_reset_nlocal = ng;
14133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14142712d1f2SDave May }
14152712d1f2SDave May 
141674d0cae8SMatthew G. Knepley /*@
141720f4b53cSBarry Smith   DMSwarmCollectViewDestroy - Resets the `DMSWARM` to the size prior to calling `DMSwarmCollectViewCreate()`
1418d3a51819SDave May 
141920f4b53cSBarry Smith   Collective
1420d3a51819SDave May 
142160225df5SJacob Faibussowitsch   Input Parameters:
142220f4b53cSBarry Smith . dm - the `DMSWARM`
1423d3a51819SDave May 
1424d3a51819SDave May   Notes:
142520f4b53cSBarry Smith   Users should call `DMSwarmCollectViewCreate()` before this function is called.
1426d3a51819SDave May 
1427d3a51819SDave May   Level: advanced
1428d3a51819SDave May 
142920f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmCollectViewCreate()`, `DMSwarmSetCollectType()`
1430d3a51819SDave May @*/
1431d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmCollectViewDestroy(DM dm)
1432d71ae5a4SJacob Faibussowitsch {
14332712d1f2SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
14342712d1f2SDave May 
1435521f74f9SMatthew G. Knepley   PetscFunctionBegin;
143628b400f6SJacob Faibussowitsch   PetscCheck(swarm->collect_view_active, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "CollectView is currently not active");
14379566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(dm, swarm->collect_view_reset_nlocal, -1));
1438480eef7bSDave May   swarm->collect_view_active = PETSC_FALSE;
14393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14402712d1f2SDave May }
14413454631fSDave May 
144266976f2fSJacob Faibussowitsch static PetscErrorCode DMSwarmSetUpPIC(DM dm)
1443d71ae5a4SJacob Faibussowitsch {
1444f0cdbbbaSDave May   PetscInt dim;
1445f0cdbbbaSDave May 
1446521f74f9SMatthew G. Knepley   PetscFunctionBegin;
14479566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetNumSpecies(dm, 1));
14489566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
144963a3b9bcSJacob Faibussowitsch   PetscCheck(dim >= 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Dimension must be 1,2,3 - found %" PetscInt_FMT, dim);
145063a3b9bcSJacob Faibussowitsch   PetscCheck(dim <= 3, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Dimension must be 1,2,3 - found %" PetscInt_FMT, dim);
14519566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dm, DMSwarmPICField_coor, dim, PETSC_DOUBLE));
14529566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dm, DMSwarmPICField_cellid, 1, PETSC_INT));
14533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1454f0cdbbbaSDave May }
1455f0cdbbbaSDave May 
145674d0cae8SMatthew G. Knepley /*@
145731403186SMatthew G. Knepley   DMSwarmSetPointCoordinatesRandom - Sets initial coordinates for particles in each cell
145831403186SMatthew G. Knepley 
145920f4b53cSBarry Smith   Collective
146031403186SMatthew G. Knepley 
146160225df5SJacob Faibussowitsch   Input Parameters:
146220f4b53cSBarry Smith + dm  - the `DMSWARM`
146320f4b53cSBarry Smith - Npc - The number of particles per cell in the cell `DM`
146431403186SMatthew G. Knepley 
146531403186SMatthew G. Knepley   Level: intermediate
146631403186SMatthew G. Knepley 
146720f4b53cSBarry Smith   Notes:
146820f4b53cSBarry Smith   The user must use `DMSwarmSetCellDM()` to set the cell `DM` first. The particles are placed randomly inside each cell. If only
146920f4b53cSBarry Smith   one particle is in each cell, it is placed at the centroid.
147020f4b53cSBarry Smith 
147120f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmSetCellDM()`
147231403186SMatthew G. Knepley @*/
1473d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmSetPointCoordinatesRandom(DM dm, PetscInt Npc)
1474d71ae5a4SJacob Faibussowitsch {
147531403186SMatthew G. Knepley   DM             cdm;
147631403186SMatthew G. Knepley   PetscRandom    rnd;
147731403186SMatthew G. Knepley   DMPolytopeType ct;
147831403186SMatthew G. Knepley   PetscBool      simplex;
147931403186SMatthew G. Knepley   PetscReal     *centroid, *coords, *xi0, *v0, *J, *invJ, detJ;
148031403186SMatthew G. Knepley   PetscInt       dim, d, cStart, cEnd, c, p;
148131403186SMatthew G. Knepley 
148231403186SMatthew G. Knepley   PetscFunctionBeginUser;
14839566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)dm), &rnd));
14849566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetInterval(rnd, -1.0, 1.0));
14859566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetType(rnd, PETSCRAND48));
148631403186SMatthew G. Knepley 
14879566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetCellDM(dm, &cdm));
14889566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(cdm, &dim));
14899566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(cdm, 0, &cStart, &cEnd));
14909566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCellType(cdm, cStart, &ct));
149131403186SMatthew G. Knepley   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE;
149231403186SMatthew G. Knepley 
14939566063dSJacob Faibussowitsch   PetscCall(PetscMalloc5(dim, &centroid, dim, &xi0, dim, &v0, dim * dim, &J, dim * dim, &invJ));
149431403186SMatthew G. Knepley   for (d = 0; d < dim; ++d) xi0[d] = -1.0;
14959566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
149631403186SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
149731403186SMatthew G. Knepley     if (Npc == 1) {
14989566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFVM(cdm, c, NULL, centroid, NULL));
149931403186SMatthew G. Knepley       for (d = 0; d < dim; ++d) coords[c * dim + d] = centroid[d];
150031403186SMatthew G. Knepley     } else {
15019566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFEM(cdm, c, NULL, v0, J, invJ, &detJ)); /* affine */
150231403186SMatthew G. Knepley       for (p = 0; p < Npc; ++p) {
150331403186SMatthew G. Knepley         const PetscInt n   = c * Npc + p;
150431403186SMatthew G. Knepley         PetscReal      sum = 0.0, refcoords[3];
150531403186SMatthew G. Knepley 
150631403186SMatthew G. Knepley         for (d = 0; d < dim; ++d) {
15079566063dSJacob Faibussowitsch           PetscCall(PetscRandomGetValueReal(rnd, &refcoords[d]));
150831403186SMatthew G. Knepley           sum += refcoords[d];
150931403186SMatthew G. Knepley         }
15109371c9d4SSatish Balay         if (simplex && sum > 0.0)
15119371c9d4SSatish Balay           for (d = 0; d < dim; ++d) refcoords[d] -= PetscSqrtReal(dim) * sum;
151231403186SMatthew G. Knepley         CoordinatesRefToReal(dim, dim, xi0, v0, J, refcoords, &coords[n * dim]);
151331403186SMatthew G. Knepley       }
151431403186SMatthew G. Knepley     }
151531403186SMatthew G. Knepley   }
15169566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
15179566063dSJacob Faibussowitsch   PetscCall(PetscFree5(centroid, xi0, v0, J, invJ));
15189566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rnd));
15193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
152031403186SMatthew G. Knepley }
152131403186SMatthew G. Knepley 
152231403186SMatthew G. Knepley /*@
152320f4b53cSBarry Smith   DMSwarmSetType - Set particular flavor of `DMSWARM`
1524d3a51819SDave May 
152520f4b53cSBarry Smith   Collective
1526d3a51819SDave May 
152760225df5SJacob Faibussowitsch   Input Parameters:
152820f4b53cSBarry Smith + dm    - the `DMSWARM`
152920f4b53cSBarry Smith - stype - the `DMSWARM` type (e.g. `DMSWARM_PIC`)
1530d3a51819SDave May 
1531d3a51819SDave May   Level: advanced
1532d3a51819SDave May 
153320f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmSetMigrateType()`, `DMSwarmSetCollectType()`, `DMSwarmType`, `DMSWARM_PIC`, `DMSWARM_BASIC`
1534d3a51819SDave May @*/
1535d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmSetType(DM dm, DMSwarmType stype)
1536d71ae5a4SJacob Faibussowitsch {
1537f0cdbbbaSDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1538f0cdbbbaSDave May 
1539521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1540f0cdbbbaSDave May   swarm->swarm_type = stype;
154148a46eb9SPierre Jolivet   if (swarm->swarm_type == DMSWARM_PIC) PetscCall(DMSwarmSetUpPIC(dm));
15423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1543f0cdbbbaSDave May }
1544f0cdbbbaSDave May 
154566976f2fSJacob Faibussowitsch static PetscErrorCode DMSetup_Swarm(DM dm)
1546d71ae5a4SJacob Faibussowitsch {
15473454631fSDave May   DM_Swarm   *swarm = (DM_Swarm *)dm->data;
15483454631fSDave May   PetscMPIInt rank;
15493454631fSDave May   PetscInt    p, npoints, *rankval;
15503454631fSDave May 
1551521f74f9SMatthew G. Knepley   PetscFunctionBegin;
15523ba16761SJacob Faibussowitsch   if (swarm->issetup) PetscFunctionReturn(PETSC_SUCCESS);
15533454631fSDave May   swarm->issetup = PETSC_TRUE;
15543454631fSDave May 
1555f0cdbbbaSDave May   if (swarm->swarm_type == DMSWARM_PIC) {
1556f0cdbbbaSDave May     /* check dmcell exists */
155728b400f6SJacob Faibussowitsch     PetscCheck(swarm->dmcell, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "DMSWARM_PIC requires you call DMSwarmSetCellDM");
1558f0cdbbbaSDave May 
1559f0cdbbbaSDave May     if (swarm->dmcell->ops->locatepointssubdomain) {
1560f0cdbbbaSDave May       /* check methods exists for exact ownership identificiation */
15619566063dSJacob Faibussowitsch       PetscCall(PetscInfo(dm, "DMSWARM_PIC: Using method CellDM->ops->LocatePointsSubdomain\n"));
1562f0cdbbbaSDave May       swarm->migrate_type = DMSWARM_MIGRATE_DMCELLEXACT;
1563f0cdbbbaSDave May     } else {
1564f0cdbbbaSDave May       /* check methods exist for point location AND rank neighbor identification */
1565f0cdbbbaSDave May       if (swarm->dmcell->ops->locatepoints) {
15669566063dSJacob Faibussowitsch         PetscCall(PetscInfo(dm, "DMSWARM_PIC: Using method CellDM->LocatePoints\n"));
1567f0cdbbbaSDave May       } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "DMSWARM_PIC requires the method CellDM->ops->locatepoints be defined");
1568f0cdbbbaSDave May 
1569f0cdbbbaSDave May       if (swarm->dmcell->ops->getneighbors) {
15709566063dSJacob Faibussowitsch         PetscCall(PetscInfo(dm, "DMSWARM_PIC: Using method CellDM->GetNeigbors\n"));
1571f0cdbbbaSDave May       } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "DMSWARM_PIC requires the method CellDM->ops->getneighbors be defined");
1572f0cdbbbaSDave May 
1573f0cdbbbaSDave May       swarm->migrate_type = DMSWARM_MIGRATE_DMCELLNSCATTER;
1574f0cdbbbaSDave May     }
1575f0cdbbbaSDave May   }
1576f0cdbbbaSDave May 
15779566063dSJacob Faibussowitsch   PetscCall(DMSwarmFinalizeFieldRegister(dm));
1578f0cdbbbaSDave May 
15793454631fSDave May   /* check some fields were registered */
158008401ef6SPierre Jolivet   PetscCheck(swarm->db->nfields > 2, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "At least one field user must be registered via DMSwarmRegisterXXX()");
15813454631fSDave May 
15823454631fSDave May   /* check local sizes were set */
158308401ef6SPierre Jolivet   PetscCheck(swarm->db->L != -1, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Local sizes must be set via DMSwarmSetLocalSizes()");
15843454631fSDave May 
15853454631fSDave May   /* initialize values in pid and rank placeholders */
15863454631fSDave May   /* TODO: [pid - use MPI_Scan] */
15879566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
15889566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
15899566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
1590ad540459SPierre Jolivet   for (p = 0; p < npoints; p++) rankval[p] = (PetscInt)rank;
15919566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
15923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15933454631fSDave May }
15943454631fSDave May 
1595dc5f5ce9SDave May extern PetscErrorCode DMSwarmSortDestroy(DMSwarmSort *_ctx);
1596dc5f5ce9SDave May 
159766976f2fSJacob Faibussowitsch static PetscErrorCode DMDestroy_Swarm(DM dm)
1598d71ae5a4SJacob Faibussowitsch {
159957795646SDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
160057795646SDave May 
160157795646SDave May   PetscFunctionBegin;
16023ba16761SJacob Faibussowitsch   if (--swarm->refct > 0) PetscFunctionReturn(PETSC_SUCCESS);
16039566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketDestroy(&swarm->db));
160448a46eb9SPierre Jolivet   if (swarm->sort_context) PetscCall(DMSwarmSortDestroy(&swarm->sort_context));
16059566063dSJacob Faibussowitsch   PetscCall(PetscFree(swarm));
16063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
160757795646SDave May }
160857795646SDave May 
160966976f2fSJacob Faibussowitsch static PetscErrorCode DMSwarmView_Draw(DM dm, PetscViewer viewer)
1610d71ae5a4SJacob Faibussowitsch {
1611a9ee3421SMatthew G. Knepley   DM         cdm;
1612a9ee3421SMatthew G. Knepley   PetscDraw  draw;
161331403186SMatthew G. Knepley   PetscReal *coords, oldPause, radius = 0.01;
1614a9ee3421SMatthew G. Knepley   PetscInt   Np, p, bs;
1615a9ee3421SMatthew G. Knepley 
1616a9ee3421SMatthew G. Knepley   PetscFunctionBegin;
16179566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, ((PetscObject)dm)->prefix, "-dm_view_swarm_radius", &radius, NULL));
16189566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
16199566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetCellDM(dm, &cdm));
16209566063dSJacob Faibussowitsch   PetscCall(PetscDrawGetPause(draw, &oldPause));
16219566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetPause(draw, 0.0));
16229566063dSJacob Faibussowitsch   PetscCall(DMView(cdm, viewer));
16239566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetPause(draw, oldPause));
1624a9ee3421SMatthew G. Knepley 
16259566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetLocalSize(dm, &Np));
16269566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, &bs, NULL, (void **)&coords));
1627a9ee3421SMatthew G. Knepley   for (p = 0; p < Np; ++p) {
1628a9ee3421SMatthew G. Knepley     const PetscInt i = p * bs;
1629a9ee3421SMatthew G. Knepley 
16309566063dSJacob Faibussowitsch     PetscCall(PetscDrawEllipse(draw, coords[i], coords[i + 1], radius, radius, PETSC_DRAW_BLUE));
1631a9ee3421SMatthew G. Knepley   }
16329566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, &bs, NULL, (void **)&coords));
16339566063dSJacob Faibussowitsch   PetscCall(PetscDrawFlush(draw));
16349566063dSJacob Faibussowitsch   PetscCall(PetscDrawPause(draw));
16359566063dSJacob Faibussowitsch   PetscCall(PetscDrawSave(draw));
16363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1637a9ee3421SMatthew G. Knepley }
1638a9ee3421SMatthew G. Knepley 
1639d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMView_Swarm_Ascii(DM dm, PetscViewer viewer)
1640d71ae5a4SJacob Faibussowitsch {
16416a5217c0SMatthew G. Knepley   PetscViewerFormat format;
16426a5217c0SMatthew G. Knepley   PetscInt         *sizes;
16436a5217c0SMatthew G. Knepley   PetscInt          dim, Np, maxSize = 17;
16446a5217c0SMatthew G. Knepley   MPI_Comm          comm;
16456a5217c0SMatthew G. Knepley   PetscMPIInt       rank, size, p;
16466a5217c0SMatthew G. Knepley   const char       *name;
16476a5217c0SMatthew G. Knepley 
16486a5217c0SMatthew G. Knepley   PetscFunctionBegin;
16496a5217c0SMatthew G. Knepley   PetscCall(PetscViewerGetFormat(viewer, &format));
16506a5217c0SMatthew G. Knepley   PetscCall(DMGetDimension(dm, &dim));
16516a5217c0SMatthew G. Knepley   PetscCall(DMSwarmGetLocalSize(dm, &Np));
16526a5217c0SMatthew G. Knepley   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
16536a5217c0SMatthew G. Knepley   PetscCallMPI(MPI_Comm_rank(comm, &rank));
16546a5217c0SMatthew G. Knepley   PetscCallMPI(MPI_Comm_size(comm, &size));
16556a5217c0SMatthew G. Knepley   PetscCall(PetscObjectGetName((PetscObject)dm, &name));
165663a3b9bcSJacob Faibussowitsch   if (name) PetscCall(PetscViewerASCIIPrintf(viewer, "%s in %" PetscInt_FMT " dimension%s:\n", name, dim, dim == 1 ? "" : "s"));
165763a3b9bcSJacob Faibussowitsch   else PetscCall(PetscViewerASCIIPrintf(viewer, "Swarm in %" PetscInt_FMT " dimension%s:\n", dim, dim == 1 ? "" : "s"));
16586a5217c0SMatthew G. Knepley   if (size < maxSize) PetscCall(PetscCalloc1(size, &sizes));
16596a5217c0SMatthew G. Knepley   else PetscCall(PetscCalloc1(3, &sizes));
16606a5217c0SMatthew G. Knepley   if (size < maxSize) {
16616a5217c0SMatthew G. Knepley     PetscCallMPI(MPI_Gather(&Np, 1, MPIU_INT, sizes, 1, MPIU_INT, 0, comm));
16626a5217c0SMatthew G. Knepley     PetscCall(PetscViewerASCIIPrintf(viewer, "  Number of particles per rank:"));
16636a5217c0SMatthew G. Knepley     for (p = 0; p < size; ++p) {
16646a5217c0SMatthew G. Knepley       if (rank == 0) PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT, sizes[p]));
16656a5217c0SMatthew G. Knepley     }
16666a5217c0SMatthew G. Knepley   } else {
16676a5217c0SMatthew G. Knepley     PetscInt locMinMax[2] = {Np, Np};
16686a5217c0SMatthew G. Knepley 
16696a5217c0SMatthew G. Knepley     PetscCall(PetscGlobalMinMaxInt(comm, locMinMax, sizes));
16706a5217c0SMatthew G. Knepley     PetscCall(PetscViewerASCIIPrintf(viewer, "  Min/Max of particles per rank: %" PetscInt_FMT "/%" PetscInt_FMT, sizes[0], sizes[1]));
16716a5217c0SMatthew G. Knepley   }
16726a5217c0SMatthew G. Knepley   PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
16736a5217c0SMatthew G. Knepley   PetscCall(PetscFree(sizes));
16746a5217c0SMatthew G. Knepley   if (format == PETSC_VIEWER_ASCII_INFO) {
16756a5217c0SMatthew G. Knepley     PetscInt *cell;
16766a5217c0SMatthew G. Knepley 
16776a5217c0SMatthew G. Knepley     PetscCall(PetscViewerASCIIPrintf(viewer, "  Cells containing each particle:\n"));
16786a5217c0SMatthew G. Knepley     PetscCall(PetscViewerASCIIPushSynchronized(viewer));
16796a5217c0SMatthew G. Knepley     PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&cell));
168063a3b9bcSJacob Faibussowitsch     for (p = 0; p < Np; ++p) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "  p%d: %" PetscInt_FMT "\n", p, cell[p]));
16816a5217c0SMatthew G. Knepley     PetscCall(PetscViewerFlush(viewer));
16826a5217c0SMatthew G. Knepley     PetscCall(PetscViewerASCIIPopSynchronized(viewer));
16836a5217c0SMatthew G. Knepley     PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&cell));
16846a5217c0SMatthew G. Knepley   }
16853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16866a5217c0SMatthew G. Knepley }
16876a5217c0SMatthew G. Knepley 
168866976f2fSJacob Faibussowitsch static PetscErrorCode DMView_Swarm(DM dm, PetscViewer viewer)
1689d71ae5a4SJacob Faibussowitsch {
16905f50eb2eSDave May   DM_Swarm *swarm = (DM_Swarm *)dm->data;
16915627991aSBarry Smith   PetscBool iascii, ibinary, isvtk, isdraw;
16925627991aSBarry Smith #if defined(PETSC_HAVE_HDF5)
16935627991aSBarry Smith   PetscBool ishdf5;
16945627991aSBarry Smith #endif
16955f50eb2eSDave May 
16965f50eb2eSDave May   PetscFunctionBegin;
16975f50eb2eSDave May   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
16985f50eb2eSDave May   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
16999566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
17009566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &ibinary));
17019566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERVTK, &isvtk));
17025627991aSBarry Smith #if defined(PETSC_HAVE_HDF5)
17039566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
17045627991aSBarry Smith #endif
17059566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
17065f50eb2eSDave May   if (iascii) {
17076a5217c0SMatthew G. Knepley     PetscViewerFormat format;
17086a5217c0SMatthew G. Knepley 
17096a5217c0SMatthew G. Knepley     PetscCall(PetscViewerGetFormat(viewer, &format));
17106a5217c0SMatthew G. Knepley     switch (format) {
1711d71ae5a4SJacob Faibussowitsch     case PETSC_VIEWER_ASCII_INFO_DETAIL:
1712d71ae5a4SJacob Faibussowitsch       PetscCall(DMSwarmDataBucketView(PetscObjectComm((PetscObject)dm), swarm->db, NULL, DATABUCKET_VIEW_STDOUT));
1713d71ae5a4SJacob Faibussowitsch       break;
1714d71ae5a4SJacob Faibussowitsch     default:
1715d71ae5a4SJacob Faibussowitsch       PetscCall(DMView_Swarm_Ascii(dm, viewer));
17166a5217c0SMatthew G. Knepley     }
1717f7d195e4SLawrence Mitchell   } else {
17185f50eb2eSDave May #if defined(PETSC_HAVE_HDF5)
171948a46eb9SPierre Jolivet     if (ishdf5) PetscCall(DMSwarmView_HDF5(dm, viewer));
17205f50eb2eSDave May #endif
172148a46eb9SPierre Jolivet     if (isdraw) PetscCall(DMSwarmView_Draw(dm, viewer));
1722f7d195e4SLawrence Mitchell   }
17233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17245f50eb2eSDave May }
17255f50eb2eSDave May 
1726d0c080abSJoseph Pusztay /*@C
172720f4b53cSBarry Smith   DMSwarmGetCellSwarm - Extracts a single cell from the `DMSWARM` object, returns it as a single cell `DMSWARM`.
172820f4b53cSBarry 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.
1729d0c080abSJoseph Pusztay 
1730d0c080abSJoseph Pusztay   Noncollective
1731d0c080abSJoseph Pusztay 
173260225df5SJacob Faibussowitsch   Input Parameters:
173320f4b53cSBarry Smith + sw        - the `DMSWARM`
17345627991aSBarry Smith . cellID    - the integer id of the cell to be extracted and filtered
173520f4b53cSBarry Smith - cellswarm - The `DMSWARM` to receive the cell
1736d0c080abSJoseph Pusztay 
1737d0c080abSJoseph Pusztay   Level: beginner
1738d0c080abSJoseph Pusztay 
17395627991aSBarry Smith   Notes:
174020f4b53cSBarry Smith   This presently only supports `DMSWARM_PIC` type
17415627991aSBarry Smith 
174220f4b53cSBarry Smith   Should be restored with `DMSwarmRestoreCellSwarm()`
1743d0c080abSJoseph Pusztay 
174420f4b53cSBarry Smith   Changes to this cell of the swarm will be lost if they are made prior to restoring this cell.
174520f4b53cSBarry Smith 
174620f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRestoreCellSwarm()`
1747d0c080abSJoseph Pusztay @*/
1748d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMSwarmGetCellSwarm(DM sw, PetscInt cellID, DM cellswarm)
1749d71ae5a4SJacob Faibussowitsch {
1750d0c080abSJoseph Pusztay   DM_Swarm *original = (DM_Swarm *)sw->data;
1751d0c080abSJoseph Pusztay   DMLabel   label;
1752d0c080abSJoseph Pusztay   DM        dmc, subdmc;
1753d0c080abSJoseph Pusztay   PetscInt *pids, particles, dim;
1754d0c080abSJoseph Pusztay 
1755d0c080abSJoseph Pusztay   PetscFunctionBegin;
1756d0c080abSJoseph Pusztay   /* Configure new swarm */
17579566063dSJacob Faibussowitsch   PetscCall(DMSetType(cellswarm, DMSWARM));
17589566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(sw, &dim));
17599566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(cellswarm, dim));
17609566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetType(cellswarm, DMSWARM_PIC));
1761d0c080abSJoseph Pusztay   /* Destroy the unused, unconfigured data bucket to prevent stragglers in memory */
17629566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketDestroy(&((DM_Swarm *)cellswarm->data)->db));
17639566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortGetAccess(sw));
17649566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortGetNumberOfPointsPerCell(sw, cellID, &particles));
17659566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortGetPointsPerCell(sw, cellID, &particles, &pids));
17669566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketCreateFromSubset(original->db, particles, pids, &((DM_Swarm *)cellswarm->data)->db));
17679566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortRestoreAccess(sw));
17689566063dSJacob Faibussowitsch   PetscCall(PetscFree(pids));
17699566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetCellDM(sw, &dmc));
17709566063dSJacob Faibussowitsch   PetscCall(DMLabelCreate(PetscObjectComm((PetscObject)sw), "singlecell", &label));
17719566063dSJacob Faibussowitsch   PetscCall(DMAddLabel(dmc, label));
17729566063dSJacob Faibussowitsch   PetscCall(DMLabelSetValue(label, cellID, 1));
17739566063dSJacob Faibussowitsch   PetscCall(DMPlexFilter(dmc, label, 1, &subdmc));
17749566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetCellDM(cellswarm, subdmc));
17759566063dSJacob Faibussowitsch   PetscCall(DMLabelDestroy(&label));
17763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1777d0c080abSJoseph Pusztay }
1778d0c080abSJoseph Pusztay 
1779d0c080abSJoseph Pusztay /*@C
178020f4b53cSBarry Smith   DMSwarmRestoreCellSwarm - Restores a `DMSWARM` object obtained with `DMSwarmGetCellSwarm()`. All fields are copied back into the parent swarm.
1781d0c080abSJoseph Pusztay 
1782d0c080abSJoseph Pusztay   Noncollective
1783d0c080abSJoseph Pusztay 
178460225df5SJacob Faibussowitsch   Input Parameters:
178520f4b53cSBarry Smith + sw        - the parent `DMSWARM`
1786d0c080abSJoseph Pusztay . cellID    - the integer id of the cell to be copied back into the parent swarm
1787d0c080abSJoseph Pusztay - cellswarm - the cell swarm object
1788d0c080abSJoseph Pusztay 
1789d0c080abSJoseph Pusztay   Level: beginner
1790d0c080abSJoseph Pusztay 
17915627991aSBarry Smith   Note:
179220f4b53cSBarry Smith   This only supports `DMSWARM_PIC` types of `DMSWARM`s
1793d0c080abSJoseph Pusztay 
179420f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmGetCellSwarm()`
1795d0c080abSJoseph Pusztay @*/
1796d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMSwarmRestoreCellSwarm(DM sw, PetscInt cellID, DM cellswarm)
1797d71ae5a4SJacob Faibussowitsch {
1798d0c080abSJoseph Pusztay   DM        dmc;
1799d0c080abSJoseph Pusztay   PetscInt *pids, particles, p;
1800d0c080abSJoseph Pusztay 
1801d0c080abSJoseph Pusztay   PetscFunctionBegin;
18029566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortGetAccess(sw));
18039566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortGetPointsPerCell(sw, cellID, &particles, &pids));
18049566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortRestoreAccess(sw));
1805d0c080abSJoseph Pusztay   /* Pointwise copy of each particle based on pid. The parent swarm may not be altered during this process. */
180648a46eb9SPierre Jolivet   for (p = 0; p < particles; ++p) PetscCall(DMSwarmDataBucketCopyPoint(((DM_Swarm *)cellswarm->data)->db, pids[p], ((DM_Swarm *)sw->data)->db, pids[p]));
1807d0c080abSJoseph Pusztay   /* Free memory, destroy cell dm */
18089566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetCellDM(cellswarm, &dmc));
18099566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmc));
18109566063dSJacob Faibussowitsch   PetscCall(PetscFree(pids));
18113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1812d0c080abSJoseph Pusztay }
1813d0c080abSJoseph Pusztay 
1814d0c080abSJoseph Pusztay PETSC_INTERN PetscErrorCode DMClone_Swarm(DM, DM *);
1815d0c080abSJoseph Pusztay 
1816d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMInitialize_Swarm(DM sw)
1817d71ae5a4SJacob Faibussowitsch {
1818d0c080abSJoseph Pusztay   PetscFunctionBegin;
1819d0c080abSJoseph Pusztay   sw->dim                           = 0;
1820d0c080abSJoseph Pusztay   sw->ops->view                     = DMView_Swarm;
1821d0c080abSJoseph Pusztay   sw->ops->load                     = NULL;
1822d0c080abSJoseph Pusztay   sw->ops->setfromoptions           = NULL;
1823d0c080abSJoseph Pusztay   sw->ops->clone                    = DMClone_Swarm;
1824d0c080abSJoseph Pusztay   sw->ops->setup                    = DMSetup_Swarm;
1825d0c080abSJoseph Pusztay   sw->ops->createlocalsection       = NULL;
1826adc21957SMatthew G. Knepley   sw->ops->createsectionpermutation = NULL;
1827d0c080abSJoseph Pusztay   sw->ops->createdefaultconstraints = NULL;
1828d0c080abSJoseph Pusztay   sw->ops->createglobalvector       = DMCreateGlobalVector_Swarm;
1829d0c080abSJoseph Pusztay   sw->ops->createlocalvector        = DMCreateLocalVector_Swarm;
1830d0c080abSJoseph Pusztay   sw->ops->getlocaltoglobalmapping  = NULL;
1831d0c080abSJoseph Pusztay   sw->ops->createfieldis            = NULL;
1832d0c080abSJoseph Pusztay   sw->ops->createcoordinatedm       = NULL;
1833d0c080abSJoseph Pusztay   sw->ops->getcoloring              = NULL;
1834d0c080abSJoseph Pusztay   sw->ops->creatematrix             = DMCreateMatrix_Swarm;
1835d0c080abSJoseph Pusztay   sw->ops->createinterpolation      = NULL;
1836d0c080abSJoseph Pusztay   sw->ops->createinjection          = NULL;
1837d0c080abSJoseph Pusztay   sw->ops->createmassmatrix         = DMCreateMassMatrix_Swarm;
1838d0c080abSJoseph Pusztay   sw->ops->refine                   = NULL;
1839d0c080abSJoseph Pusztay   sw->ops->coarsen                  = NULL;
1840d0c080abSJoseph Pusztay   sw->ops->refinehierarchy          = NULL;
1841d0c080abSJoseph Pusztay   sw->ops->coarsenhierarchy         = NULL;
1842d0c080abSJoseph Pusztay   sw->ops->globaltolocalbegin       = NULL;
1843d0c080abSJoseph Pusztay   sw->ops->globaltolocalend         = NULL;
1844d0c080abSJoseph Pusztay   sw->ops->localtoglobalbegin       = NULL;
1845d0c080abSJoseph Pusztay   sw->ops->localtoglobalend         = NULL;
1846d0c080abSJoseph Pusztay   sw->ops->destroy                  = DMDestroy_Swarm;
1847d0c080abSJoseph Pusztay   sw->ops->createsubdm              = NULL;
1848d0c080abSJoseph Pusztay   sw->ops->getdimpoints             = NULL;
1849d0c080abSJoseph Pusztay   sw->ops->locatepoints             = NULL;
18503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1851d0c080abSJoseph Pusztay }
1852d0c080abSJoseph Pusztay 
1853d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode DMClone_Swarm(DM dm, DM *newdm)
1854d71ae5a4SJacob Faibussowitsch {
1855d0c080abSJoseph Pusztay   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1856d0c080abSJoseph Pusztay 
1857d0c080abSJoseph Pusztay   PetscFunctionBegin;
1858d0c080abSJoseph Pusztay   swarm->refct++;
1859d0c080abSJoseph Pusztay   (*newdm)->data = swarm;
18609566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)*newdm, DMSWARM));
18619566063dSJacob Faibussowitsch   PetscCall(DMInitialize_Swarm(*newdm));
18622e56dffeSJoe Pusztay   (*newdm)->dim = dm->dim;
18633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1864d0c080abSJoseph Pusztay }
1865d0c080abSJoseph Pusztay 
1866d3a51819SDave May /*MC
1867d3a51819SDave May 
186820f4b53cSBarry Smith  DMSWARM = "swarm" - A `DM` object used to represent arrays of data (fields) of arbitrary data type.
186962741f57SDave May  This implementation was designed for particle methods in which the underlying
1870d3a51819SDave May  data required to be represented is both (i) dynamic in length, (ii) and of arbitrary data type.
1871d3a51819SDave May 
187220f4b53cSBarry Smith  Level: intermediate
187320f4b53cSBarry Smith 
187420f4b53cSBarry Smith   Notes:
187520f4b53cSBarry Smith  User data can be represented by `DMSWARM` through a registering "fields".
187662741f57SDave May  To register a field, the user must provide;
187762741f57SDave May  (a) a unique name;
187862741f57SDave May  (b) the data type (or size in bytes);
187962741f57SDave May  (c) the block size of the data.
1880d3a51819SDave May 
1881d3a51819SDave May  For example, suppose the application requires a unique id, energy, momentum and density to be stored
1882c215a47eSMatthew Knepley  on a set of particles. Then the following code could be used
188320f4b53cSBarry Smith .vb
188420f4b53cSBarry Smith     DMSwarmInitializeFieldRegister(dm)
188520f4b53cSBarry Smith     DMSwarmRegisterPetscDatatypeField(dm,"uid",1,PETSC_LONG);
188620f4b53cSBarry Smith     DMSwarmRegisterPetscDatatypeField(dm,"energy",1,PETSC_REAL);
188720f4b53cSBarry Smith     DMSwarmRegisterPetscDatatypeField(dm,"momentum",3,PETSC_REAL);
188820f4b53cSBarry Smith     DMSwarmRegisterPetscDatatypeField(dm,"density",1,PETSC_FLOAT);
188920f4b53cSBarry Smith     DMSwarmFinalizeFieldRegister(dm)
189020f4b53cSBarry Smith .ve
1891d3a51819SDave May 
189220f4b53cSBarry Smith  The fields represented by `DMSWARM` are dynamic and can be re-sized at any time.
189320f4b53cSBarry Smith  The only restriction imposed by `DMSWARM` is that all fields contain the same number of points.
1894d3a51819SDave May 
1895d3a51819SDave May  To support particle methods, "migration" techniques are provided. These methods migrate data
18965627991aSBarry Smith  between ranks.
1897d3a51819SDave May 
189820f4b53cSBarry Smith  `DMSWARM` supports the methods `DMCreateGlobalVector()` and `DMCreateLocalVector()`.
189920f4b53cSBarry Smith  As a `DMSWARM` may internally define and store values of different data types,
190020f4b53cSBarry Smith  before calling `DMCreateGlobalVector()` or `DMCreateLocalVector()`, the user must inform `DMSWARM` which
190120f4b53cSBarry Smith  fields should be used to define a `Vec` object via
190220f4b53cSBarry Smith    `DMSwarmVectorDefineField()`
1903c215a47eSMatthew Knepley  The specified field can be changed at any time - thereby permitting vectors
1904c215a47eSMatthew Knepley  compatible with different fields to be created.
1905d3a51819SDave May 
190620f4b53cSBarry Smith  A dual representation of fields in the `DMSWARM` and a Vec object is permitted via
190720f4b53cSBarry Smith    `DMSwarmCreateGlobalVectorFromField()`
190820f4b53cSBarry Smith  Here the data defining the field in the `DMSWARM` is shared with a Vec.
1909d3a51819SDave May  This is inherently unsafe if you alter the size of the field at any time between
191020f4b53cSBarry Smith  calls to `DMSwarmCreateGlobalVectorFromField()` and `DMSwarmDestroyGlobalVectorFromField()`.
191120f4b53cSBarry Smith  If the local size of the `DMSWARM` does not match the local size of the global vector
191220f4b53cSBarry Smith  when `DMSwarmDestroyGlobalVectorFromField()` is called, an error is thrown.
1913d3a51819SDave May 
191462741f57SDave May  Additional high-level support is provided for Particle-In-Cell methods.
191520f4b53cSBarry Smith  Please refer to `DMSwarmSetType()`.
191662741f57SDave May 
191720f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMType`, `DMCreate()`, `DMSetType()`
1918d3a51819SDave May M*/
1919d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMCreate_Swarm(DM dm)
1920d71ae5a4SJacob Faibussowitsch {
192157795646SDave May   DM_Swarm *swarm;
192257795646SDave May 
192357795646SDave May   PetscFunctionBegin;
192457795646SDave May   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
19254dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&swarm));
1926f0cdbbbaSDave May   dm->data = swarm;
19279566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketCreate(&swarm->db));
19289566063dSJacob Faibussowitsch   PetscCall(DMSwarmInitializeFieldRegister(dm));
1929d0c080abSJoseph Pusztay   swarm->refct                          = 1;
1930b5bcf523SDave May   swarm->vec_field_set                  = PETSC_FALSE;
19313454631fSDave May   swarm->issetup                        = PETSC_FALSE;
1932480eef7bSDave May   swarm->swarm_type                     = DMSWARM_BASIC;
1933480eef7bSDave May   swarm->migrate_type                   = DMSWARM_MIGRATE_BASIC;
1934480eef7bSDave May   swarm->collect_type                   = DMSWARM_COLLECT_BASIC;
193540c453e9SDave May   swarm->migrate_error_on_missing_point = PETSC_FALSE;
1936f0cdbbbaSDave May   swarm->dmcell                         = NULL;
1937f0cdbbbaSDave May   swarm->collect_view_active            = PETSC_FALSE;
1938f0cdbbbaSDave May   swarm->collect_view_reset_nlocal      = -1;
19399566063dSJacob Faibussowitsch   PetscCall(DMInitialize_Swarm(dm));
19402e956fe4SStefano Zampini   if (SwarmDataFieldId == -1) PetscCall(PetscObjectComposedDataRegister(&SwarmDataFieldId));
19413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
194257795646SDave May }
1943