xref: /petsc/src/dm/impls/swarm/swarm.c (revision 6a5217c03994f2d95bb2e6dbd8bed42381aeb015)
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 
2774d0cae8SMatthew G. Knepley #if defined(PETSC_HAVE_HDF5)
2874d0cae8SMatthew G. Knepley #include <petscviewerhdf5.h>
2974d0cae8SMatthew G. Knepley 
3074d0cae8SMatthew G. Knepley PetscErrorCode VecView_Swarm_HDF5_Internal(Vec v, PetscViewer viewer)
3174d0cae8SMatthew G. Knepley {
3274d0cae8SMatthew G. Knepley   DM             dm;
3374d0cae8SMatthew G. Knepley   PetscReal      seqval;
3474d0cae8SMatthew G. Knepley   PetscInt       seqnum, bs;
3574d0cae8SMatthew G. Knepley   PetscBool      isseq;
3674d0cae8SMatthew G. Knepley 
3774d0cae8SMatthew G. Knepley   PetscFunctionBegin;
389566063dSJacob Faibussowitsch   PetscCall(VecGetDM(v, &dm));
399566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(v, &bs));
409566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5PushGroup(viewer, "/particle_fields"));
419566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject) v, VECSEQ, &isseq));
429566063dSJacob Faibussowitsch   PetscCall(DMGetOutputSequenceNumber(dm, &seqnum, &seqval));
439566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5SetTimestep(viewer, seqnum));
449566063dSJacob Faibussowitsch   /* PetscCall(DMSequenceView_HDF5(dm, "time", seqnum, (PetscScalar) seqval, viewer)); */
459566063dSJacob Faibussowitsch   PetscCall(VecViewNative(v, viewer));
469566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject) v, "Nc", PETSC_INT, (void *) &bs));
479566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5PopGroup(viewer));
4874d0cae8SMatthew G. Knepley   PetscFunctionReturn(0);
4974d0cae8SMatthew G. Knepley }
5074d0cae8SMatthew G. Knepley 
5174d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmView_HDF5(DM dm, PetscViewer viewer)
5274d0cae8SMatthew G. Knepley {
5374d0cae8SMatthew G. Knepley   Vec            coordinates;
5474d0cae8SMatthew G. Knepley   PetscInt       Np;
5574d0cae8SMatthew G. Knepley   PetscBool      isseq;
5674d0cae8SMatthew G. Knepley 
5774d0cae8SMatthew G. Knepley   PetscFunctionBegin;
589566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetSize(dm, &Np));
599566063dSJacob Faibussowitsch   PetscCall(DMSwarmCreateGlobalVectorFromField(dm, DMSwarmPICField_coor, &coordinates));
609566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) coordinates, "coordinates"));
619566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5PushGroup(viewer, "/particles"));
629566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject) coordinates, VECSEQ, &isseq));
639566063dSJacob Faibussowitsch   PetscCall(VecViewNative(coordinates, viewer));
649566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject) coordinates, "Np", PETSC_INT, (void *) &Np));
659566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5PopGroup(viewer));
669566063dSJacob Faibussowitsch   PetscCall(DMSwarmDestroyGlobalVectorFromField(dm, DMSwarmPICField_coor, &coordinates));
6774d0cae8SMatthew G. Knepley   PetscFunctionReturn(0);
6874d0cae8SMatthew G. Knepley }
6974d0cae8SMatthew G. Knepley #endif
7074d0cae8SMatthew G. Knepley 
7174d0cae8SMatthew G. Knepley PetscErrorCode VecView_Swarm(Vec v, PetscViewer viewer)
7274d0cae8SMatthew G. Knepley {
7374d0cae8SMatthew G. Knepley   DM             dm;
74f9558f5cSBarry Smith #if defined(PETSC_HAVE_HDF5)
7574d0cae8SMatthew G. Knepley   PetscBool      ishdf5;
76f9558f5cSBarry Smith #endif
7774d0cae8SMatthew G. Knepley 
7874d0cae8SMatthew G. Knepley   PetscFunctionBegin;
799566063dSJacob Faibussowitsch   PetscCall(VecGetDM(v, &dm));
8028b400f6SJacob Faibussowitsch   PetscCheck(dm,PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_WRONG, "Vector not generated from a DM");
81f9558f5cSBarry Smith #if defined(PETSC_HAVE_HDF5)
829566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERHDF5,  &ishdf5));
8374d0cae8SMatthew G. Knepley   if (ishdf5) {
849566063dSJacob Faibussowitsch       PetscCall(VecView_Swarm_HDF5_Internal(v, viewer));
85f9558f5cSBarry Smith       PetscFunctionReturn(0);
8674d0cae8SMatthew G. Knepley   }
87f9558f5cSBarry Smith #endif
889566063dSJacob Faibussowitsch   PetscCall(VecViewNative(v, viewer));
8974d0cae8SMatthew G. Knepley   PetscFunctionReturn(0);
9074d0cae8SMatthew G. Knepley }
9174d0cae8SMatthew G. Knepley 
92d3a51819SDave May /*@C
9362741f57SDave May    DMSwarmVectorDefineField - Sets the field from which to define a Vec object
9462741f57SDave May                              when DMCreateLocalVector(), or DMCreateGlobalVector() is called
9557795646SDave May 
96d083f849SBarry Smith    Collective on dm
9757795646SDave May 
98d3a51819SDave May    Input parameters:
9962741f57SDave May +  dm - a DMSwarm
10062741f57SDave May -  fieldname - the textual name given to a registered field
10157795646SDave May 
102d3a51819SDave May    Level: beginner
10357795646SDave May 
104d3a51819SDave May    Notes:
105e7af74c8SDave May 
10662741f57SDave May    The field with name fieldname must be defined as having a data type of PetscScalar.
107e7af74c8SDave May 
108d3a51819SDave May    This function must be called prior to calling DMCreateLocalVector(), DMCreateGlobalVector().
109d3a51819SDave May    Mutiple calls to DMSwarmVectorDefineField() are permitted.
11057795646SDave May 
1118b8a3813SDave May .seealso: DMSwarmRegisterPetscDatatypeField(), DMCreateGlobalVector(), DMCreateLocalVector()
112d3a51819SDave May @*/
11374d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmVectorDefineField(DM dm,const char fieldname[])
114b5bcf523SDave May {
115b5bcf523SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
116b5bcf523SDave May   PetscInt       bs,n;
117b5bcf523SDave May   PetscScalar    *array;
118b5bcf523SDave May   PetscDataType  type;
119b5bcf523SDave May 
120a9cbaee5SMatthew G. Knepley   PetscFunctionBegin;
1219566063dSJacob Faibussowitsch   if (!swarm->issetup) PetscCall(DMSetUp(dm));
1229566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db,&n,NULL,NULL));
1239566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm,fieldname,&bs,&type,(void**)&array));
124b5bcf523SDave May 
125b5bcf523SDave May   /* Check all fields are of type PETSC_REAL or PETSC_SCALAR */
12608401ef6SPierre Jolivet   PetscCheck(type == PETSC_REAL,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only valid for PETSC_REAL");
1279566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(swarm->vec_field_name,PETSC_MAX_PATH_LEN-1,"%s",fieldname));
128b5bcf523SDave May   swarm->vec_field_set = PETSC_TRUE;
1291b1ea282SDave May   swarm->vec_field_bs = bs;
130b5bcf523SDave May   swarm->vec_field_nlocal = n;
1319566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm,fieldname,&bs,&type,(void**)&array));
132b5bcf523SDave May   PetscFunctionReturn(0);
133b5bcf523SDave May }
134b5bcf523SDave May 
135cc651181SDave May /* requires DMSwarmDefineFieldVector has been called */
136b5bcf523SDave May PetscErrorCode DMCreateGlobalVector_Swarm(DM dm,Vec *vec)
137b5bcf523SDave May {
138b5bcf523SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
139b5bcf523SDave May   Vec            x;
140b5bcf523SDave May   char           name[PETSC_MAX_PATH_LEN];
141b5bcf523SDave May 
142a9cbaee5SMatthew G. Knepley   PetscFunctionBegin;
1439566063dSJacob Faibussowitsch   if (!swarm->issetup) PetscCall(DMSetUp(dm));
14428b400f6SJacob Faibussowitsch   PetscCheck(swarm->vec_field_set,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMSwarmVectorDefineField first");
14508401ef6SPierre 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 */
146cc651181SDave May 
1479566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"DMSwarmField_%s",swarm->vec_field_name));
1489566063dSJacob Faibussowitsch   PetscCall(VecCreate(PetscObjectComm((PetscObject)dm),&x));
1499566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)x,name));
1509566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(x,swarm->db->L*swarm->vec_field_bs,PETSC_DETERMINE));
1519566063dSJacob Faibussowitsch   PetscCall(VecSetBlockSize(x,swarm->vec_field_bs));
1529566063dSJacob Faibussowitsch   PetscCall(VecSetDM(x,dm));
1539566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(x));
1549566063dSJacob Faibussowitsch   PetscCall(VecSetDM(x, dm));
1559566063dSJacob Faibussowitsch   PetscCall(VecSetOperation(x, VECOP_VIEW, (void (*)(void)) VecView_Swarm));
156b5bcf523SDave May   *vec = x;
157b5bcf523SDave May   PetscFunctionReturn(0);
158b5bcf523SDave May }
159b5bcf523SDave May 
160b5bcf523SDave May /* requires DMSwarmDefineFieldVector has been called */
161b5bcf523SDave May PetscErrorCode DMCreateLocalVector_Swarm(DM dm,Vec *vec)
162b5bcf523SDave May {
163b5bcf523SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
164b5bcf523SDave May   Vec            x;
165b5bcf523SDave May   char           name[PETSC_MAX_PATH_LEN];
166b5bcf523SDave May 
167a9cbaee5SMatthew G. Knepley   PetscFunctionBegin;
1689566063dSJacob Faibussowitsch   if (!swarm->issetup) PetscCall(DMSetUp(dm));
16928b400f6SJacob Faibussowitsch   PetscCheck(swarm->vec_field_set,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMSwarmVectorDefineField first");
17008401ef6SPierre 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");
171cc651181SDave May 
1729566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"DMSwarmField_%s",swarm->vec_field_name));
1739566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_SELF,&x));
1749566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)x,name));
1759566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(x,swarm->db->L*swarm->vec_field_bs,PETSC_DETERMINE));
1769566063dSJacob Faibussowitsch   PetscCall(VecSetBlockSize(x,swarm->vec_field_bs));
1779566063dSJacob Faibussowitsch   PetscCall(VecSetDM(x,dm));
1789566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(x));
179b5bcf523SDave May   *vec = x;
180b5bcf523SDave May   PetscFunctionReturn(0);
181b5bcf523SDave May }
182b5bcf523SDave May 
183fb1bcc12SMatthew G. Knepley static PetscErrorCode DMSwarmDestroyVectorFromField_Private(DM dm, const char fieldname[], Vec *vec)
184fb1bcc12SMatthew G. Knepley {
185fb1bcc12SMatthew G. Knepley   DM_Swarm         *swarm = (DM_Swarm *) dm->data;
18677048351SPatrick Sanan   DMSwarmDataField gfield;
187fb1bcc12SMatthew G. Knepley   void             (*fptr)(void);
188fb1bcc12SMatthew G. Knepley   PetscInt         bs, nlocal;
189fb1bcc12SMatthew G. Knepley   char             name[PETSC_MAX_PATH_LEN];
190d3a51819SDave May 
191fb1bcc12SMatthew G. Knepley   PetscFunctionBegin;
1929566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(*vec, &nlocal));
1939566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(*vec, &bs));
19408401ef6SPierre 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");
1959566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield));
196fb1bcc12SMatthew G. Knepley   /* check vector is an inplace array */
1979566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN-1, "DMSwarm_VecFieldInPlace_%s", fieldname));
1989566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject) *vec, name, &fptr));
19928b400f6SJacob Faibussowitsch   PetscCheck(fptr,PetscObjectComm((PetscObject) dm), PETSC_ERR_USER, "Vector being destroyed was not created from DMSwarm field(%s)", fieldname);
2009566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataFieldRestoreAccess(gfield));
2018df78e51SMark Adams   PetscCall(VecResetArray(*vec));
2029566063dSJacob Faibussowitsch   PetscCall(VecDestroy(vec));
203fb1bcc12SMatthew G. Knepley   PetscFunctionReturn(0);
204fb1bcc12SMatthew G. Knepley }
205fb1bcc12SMatthew G. Knepley 
206fb1bcc12SMatthew G. Knepley static PetscErrorCode DMSwarmCreateVectorFromField_Private(DM dm, const char fieldname[], MPI_Comm comm, Vec *vec)
207fb1bcc12SMatthew G. Knepley {
208fb1bcc12SMatthew G. Knepley   DM_Swarm      *swarm = (DM_Swarm *) dm->data;
209fb1bcc12SMatthew G. Knepley   PetscDataType  type;
210fb1bcc12SMatthew G. Knepley   PetscScalar   *array;
211fb1bcc12SMatthew G. Knepley   PetscInt       bs, n;
212fb1bcc12SMatthew G. Knepley   char           name[PETSC_MAX_PATH_LEN];
213e4fbd051SBarry Smith   PetscMPIInt    size;
2148df78e51SMark Adams   PetscBool      iscuda,iskokkos;
215fb1bcc12SMatthew G. Knepley 
216fb1bcc12SMatthew G. Knepley   PetscFunctionBegin;
2179566063dSJacob Faibussowitsch   if (!swarm->issetup) PetscCall(DMSetUp(dm));
2189566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &n, NULL, NULL));
2199566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, fieldname, &bs, &type, (void **) &array));
22008401ef6SPierre Jolivet   PetscCheck(type == PETSC_REAL,PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Only valid for PETSC_REAL");
221fb1bcc12SMatthew G. Knepley 
2229566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
2238df78e51SMark Adams   PetscCall(PetscStrcmp(dm->vectype,VECKOKKOS,&iskokkos));
2248df78e51SMark Adams   PetscCall(PetscStrcmp(dm->vectype,VECCUDA,&iscuda));
2258df78e51SMark Adams   PetscCall(VecCreate(comm,vec));
2268df78e51SMark Adams   PetscCall(VecSetSizes(*vec,n*bs,PETSC_DETERMINE));
2278df78e51SMark Adams   PetscCall(VecSetBlockSize(*vec,bs));
2288df78e51SMark Adams   if (iskokkos) PetscCall(VecSetType(*vec,VECKOKKOS));
2298df78e51SMark Adams   else if (iscuda) PetscCall(VecSetType(*vec,VECCUDA));
2308df78e51SMark Adams   else PetscCall(VecSetType(*vec,VECSTANDARD));
2318df78e51SMark Adams   PetscCall(VecPlaceArray(*vec,array));
2328df78e51SMark Adams 
2339566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN-1, "DMSwarmSharedField_%s", fieldname));
2349566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) *vec, name));
235fb1bcc12SMatthew G. Knepley 
236fb1bcc12SMatthew G. Knepley   /* Set guard */
2379566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN-1, "DMSwarm_VecFieldInPlace_%s", fieldname));
2389566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject) *vec, name, DMSwarmDestroyVectorFromField_Private));
23974d0cae8SMatthew G. Knepley 
2409566063dSJacob Faibussowitsch   PetscCall(VecSetDM(*vec, dm));
2419566063dSJacob Faibussowitsch   PetscCall(VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Swarm));
242fb1bcc12SMatthew G. Knepley   PetscFunctionReturn(0);
243fb1bcc12SMatthew G. Knepley }
244fb1bcc12SMatthew G. Knepley 
2450643ed39SMatthew G. Knepley /* This creates a "mass matrix" between a finite element and particle space. If a finite element interpolant is given by
2460643ed39SMatthew G. Knepley 
2470643ed39SMatthew G. Knepley      \hat f = \sum_i f_i \phi_i
2480643ed39SMatthew G. Knepley 
2490643ed39SMatthew G. Knepley    and a particle function is given by
2500643ed39SMatthew G. Knepley 
2510643ed39SMatthew G. Knepley      f = \sum_i w_i \delta(x - x_i)
2520643ed39SMatthew G. Knepley 
2530643ed39SMatthew G. Knepley    then we want to require that
2540643ed39SMatthew G. Knepley 
2550643ed39SMatthew G. Knepley      M \hat f = M_p f
2560643ed39SMatthew G. Knepley 
2570643ed39SMatthew G. Knepley    where the particle mass matrix is given by
2580643ed39SMatthew G. Knepley 
2590643ed39SMatthew G. Knepley      (M_p)_{ij} = \int \phi_i \delta(x - x_j)
2600643ed39SMatthew G. Knepley 
2610643ed39SMatthew G. Knepley    The way Dave May does particles, they amount to quadratue weights rather than delta functions, so he has |J| is in
2620643ed39SMatthew G. Knepley    his integral. We allow this with the boolean flag.
2630643ed39SMatthew G. Knepley */
2640643ed39SMatthew G. Knepley static PetscErrorCode DMSwarmComputeMassMatrix_Private(DM dmc, DM dmf, Mat mass, PetscBool useDeltaFunction, void *ctx)
26583c47955SMatthew G. Knepley {
26683c47955SMatthew G. Knepley   const char    *name = "Mass Matrix";
2670643ed39SMatthew G. Knepley   MPI_Comm       comm;
26883c47955SMatthew G. Knepley   PetscDS        prob;
26983c47955SMatthew G. Knepley   PetscSection   fsection, globalFSection;
270e8f14785SLisandro Dalcin   PetscHSetIJ    ht;
2710643ed39SMatthew G. Knepley   PetscLayout    rLayout, colLayout;
27283c47955SMatthew G. Knepley   PetscInt      *dnz, *onz;
273adb2528bSMark Adams   PetscInt       locRows, locCols, rStart, colStart, colEnd, *rowIDXs;
2740643ed39SMatthew G. Knepley   PetscReal     *xi, *v0, *J, *invJ, detJ = 1.0, v0ref[3] = {-1.0, -1.0, -1.0};
27583c47955SMatthew G. Knepley   PetscScalar   *elemMat;
2760643ed39SMatthew G. Knepley   PetscInt       dim, Nf, field, cStart, cEnd, cell, totDim, maxC = 0;
27783c47955SMatthew G. Knepley 
27883c47955SMatthew G. Knepley   PetscFunctionBegin;
2799566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject) mass, &comm));
2809566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dmf, &dim));
2819566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dmf, &prob));
2829566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(prob, &Nf));
2839566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(prob, &totDim));
2849566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(dim, &v0, dim*dim, &J, dim*dim,&invJ));
2859566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dmf, &fsection));
2869566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalSection(dmf, &globalFSection));
2879566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd));
2889566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(mass, &locRows, &locCols));
28983c47955SMatthew G. Knepley 
2909566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreate(comm, &colLayout));
2919566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetLocalSize(colLayout, locCols));
2929566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(colLayout, 1));
2939566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(colLayout));
2949566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetRange(colLayout, &colStart, &colEnd));
2959566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&colLayout));
2960643ed39SMatthew G. Knepley 
2979566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreate(comm, &rLayout));
2989566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetLocalSize(rLayout, locRows));
2999566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(rLayout, 1));
3009566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(rLayout));
3019566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetRange(rLayout, &rStart, NULL));
3029566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&rLayout));
3030643ed39SMatthew G. Knepley 
3049566063dSJacob Faibussowitsch   PetscCall(PetscCalloc2(locRows, &dnz, locRows, &onz));
3059566063dSJacob Faibussowitsch   PetscCall(PetscHSetIJCreate(&ht));
30653e60ab4SJoseph Pusztay 
3079566063dSJacob Faibussowitsch   PetscCall(PetscSynchronizedFlush(comm, NULL));
3080643ed39SMatthew G. Knepley   /* count non-zeros */
3099566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortGetAccess(dmc));
31083c47955SMatthew G. Knepley   for (field = 0; field < Nf; ++field) {
31183c47955SMatthew G. Knepley     for (cell = cStart; cell < cEnd; ++cell) {
3120643ed39SMatthew G. Knepley       PetscInt  c, i;
3130643ed39SMatthew G. Knepley       PetscInt *findices,   *cindices; /* fine is vertices, coarse is particles */
31483c47955SMatthew G. Knepley       PetscInt  numFIndices, numCIndices;
31583c47955SMatthew G. Knepley 
3169566063dSJacob Faibussowitsch       PetscCall(DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
3179566063dSJacob Faibussowitsch       PetscCall(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices));
318fc7c92abSMatthew G. Knepley       maxC = PetscMax(maxC, numCIndices);
31983c47955SMatthew G. Knepley       {
320e8f14785SLisandro Dalcin         PetscHashIJKey key;
321e8f14785SLisandro Dalcin         PetscBool      missing;
32283c47955SMatthew G. Knepley         for (i = 0; i < numFIndices; ++i) {
323adb2528bSMark Adams           key.j = findices[i]; /* global column (from Plex) */
324adb2528bSMark Adams           if (key.j >= 0) {
32583c47955SMatthew G. Knepley             /* Get indices for coarse elements */
32683c47955SMatthew G. Knepley             for (c = 0; c < numCIndices; ++c) {
327adb2528bSMark Adams               key.i = cindices[c] + rStart; /* global cols (from Swarm) */
328adb2528bSMark Adams               if (key.i < 0) continue;
3299566063dSJacob Faibussowitsch               PetscCall(PetscHSetIJQueryAdd(ht, key, &missing));
33083c47955SMatthew G. Knepley               if (missing) {
3310643ed39SMatthew G. Knepley                 if ((key.j >= colStart) && (key.j < colEnd)) ++dnz[key.i - rStart];
332e8f14785SLisandro Dalcin                 else                                         ++onz[key.i - rStart];
33398921bdaSJacob Faibussowitsch               } else SETERRQ(PetscObjectComm((PetscObject) dmf), PETSC_ERR_SUP, "Set new value at %D,%D", key.i, key.j);
33483c47955SMatthew G. Knepley             }
335fc7c92abSMatthew G. Knepley           }
336fc7c92abSMatthew G. Knepley         }
3379566063dSJacob Faibussowitsch         PetscCall(PetscFree(cindices));
33883c47955SMatthew G. Knepley       }
3399566063dSJacob Faibussowitsch       PetscCall(DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
34083c47955SMatthew G. Knepley     }
34183c47955SMatthew G. Knepley   }
3429566063dSJacob Faibussowitsch   PetscCall(PetscHSetIJDestroy(&ht));
3439566063dSJacob Faibussowitsch   PetscCall(MatXAIJSetPreallocation(mass, 1, dnz, onz, NULL, NULL));
3449566063dSJacob Faibussowitsch   PetscCall(MatSetOption(mass, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
3459566063dSJacob Faibussowitsch   PetscCall(PetscFree2(dnz, onz));
3469566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(maxC*totDim, &elemMat, maxC, &rowIDXs, maxC*dim, &xi));
34783c47955SMatthew G. Knepley   for (field = 0; field < Nf; ++field) {
348ef0bb6c7SMatthew G. Knepley     PetscTabulation Tcoarse;
34983c47955SMatthew G. Knepley     PetscObject     obj;
350ef0bb6c7SMatthew G. Knepley     PetscReal       *coords;
3510643ed39SMatthew G. Knepley     PetscInt        Nc, i;
35283c47955SMatthew G. Knepley 
3539566063dSJacob Faibussowitsch     PetscCall(PetscDSGetDiscretization(prob, field, &obj));
3549566063dSJacob Faibussowitsch     PetscCall(PetscFEGetNumComponents((PetscFE) obj, &Nc));
35508401ef6SPierre Jolivet     PetscCheck(Nc == 1,PetscObjectComm((PetscObject) dmf), PETSC_ERR_SUP, "Can only interpolate a scalar field from particles, Nc = %D", Nc);
3569566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **) &coords));
35783c47955SMatthew G. Knepley     for (cell = cStart; cell < cEnd; ++cell) {
35883c47955SMatthew G. Knepley       PetscInt *findices  , *cindices;
35983c47955SMatthew G. Knepley       PetscInt  numFIndices, numCIndices;
3600643ed39SMatthew G. Knepley       PetscInt  p, c;
36183c47955SMatthew G. Knepley 
3620643ed39SMatthew G. Knepley       /* TODO: Use DMField instead of assuming affine */
3639566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ));
3649566063dSJacob Faibussowitsch       PetscCall(DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
3659566063dSJacob Faibussowitsch       PetscCall(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices));
3660643ed39SMatthew G. Knepley       for (p = 0; p < numCIndices; ++p) {
3670643ed39SMatthew G. Knepley         CoordinatesRealToRef(dim, dim, v0ref, v0, invJ, &coords[cindices[p]*dim], &xi[p*dim]);
3680643ed39SMatthew G. Knepley       }
3699566063dSJacob Faibussowitsch       PetscCall(PetscFECreateTabulation((PetscFE) obj, 1, numCIndices, xi, 0, &Tcoarse));
37083c47955SMatthew G. Knepley       /* Get elemMat entries by multiplying by weight */
3719566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(elemMat, numCIndices*totDim));
37283c47955SMatthew G. Knepley       for (i = 0; i < numFIndices; ++i) {
3730643ed39SMatthew G. Knepley         for (p = 0; p < numCIndices; ++p) {
3740643ed39SMatthew G. Knepley           for (c = 0; c < Nc; ++c) {
3750643ed39SMatthew G. Knepley             /* B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c */
376ef0bb6c7SMatthew G. Knepley             elemMat[p*numFIndices+i] += Tcoarse->T[0][(p*numFIndices + i)*Nc + c]*(useDeltaFunction ? 1.0 : detJ);
37783c47955SMatthew G. Knepley           }
3780643ed39SMatthew G. Knepley         }
3790643ed39SMatthew G. Knepley       }
380adb2528bSMark Adams       for (p = 0; p < numCIndices; ++p) rowIDXs[p] = cindices[p] + rStart;
3819566063dSJacob Faibussowitsch       if (0) PetscCall(DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat));
3829566063dSJacob Faibussowitsch       PetscCall(MatSetValues(mass, numCIndices, rowIDXs, numFIndices, findices, elemMat, ADD_VALUES));
3839566063dSJacob Faibussowitsch       PetscCall(PetscFree(cindices));
3849566063dSJacob Faibussowitsch       PetscCall(DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
3859566063dSJacob Faibussowitsch       PetscCall(PetscTabulationDestroy(&Tcoarse));
38683c47955SMatthew G. Knepley     }
3879566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **) &coords));
38883c47955SMatthew G. Knepley   }
3899566063dSJacob Faibussowitsch   PetscCall(PetscFree3(elemMat, rowIDXs, xi));
3909566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortRestoreAccess(dmc));
3919566063dSJacob Faibussowitsch   PetscCall(PetscFree3(v0, J, invJ));
3929566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mass, MAT_FINAL_ASSEMBLY));
3939566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mass, MAT_FINAL_ASSEMBLY));
39483c47955SMatthew G. Knepley   PetscFunctionReturn(0);
39583c47955SMatthew G. Knepley }
39683c47955SMatthew G. Knepley 
397d0c080abSJoseph Pusztay /* Returns empty matrix for use with SNES FD */
39870a7d78aSStefano Zampini static PetscErrorCode DMCreateMatrix_Swarm(DM sw, Mat* m)
39970a7d78aSStefano Zampini {
400d0c080abSJoseph Pusztay   Vec            field;
401d0c080abSJoseph Pusztay   PetscInt       size;
402d0c080abSJoseph Pusztay 
403d0c080abSJoseph Pusztay   PetscFunctionBegin;
4049566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(sw, &field));
4059566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(field, &size));
4069566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(sw, &field));
4079566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, m));
4089566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(*m));
4099566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*m, PETSC_DECIDE, PETSC_DECIDE, size, size));
4109566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(*m, 1, NULL));
4119566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(*m));
4129566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(*m, MAT_FINAL_ASSEMBLY));
4139566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(*m, MAT_FINAL_ASSEMBLY));
4149566063dSJacob Faibussowitsch   PetscCall(MatShift(*m, 1.0));
4159566063dSJacob Faibussowitsch   PetscCall(MatSetDM(*m, sw));
416d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
417d0c080abSJoseph Pusztay }
418d0c080abSJoseph Pusztay 
419adb2528bSMark Adams /* FEM cols, Particle rows */
42083c47955SMatthew G. Knepley static PetscErrorCode DMCreateMassMatrix_Swarm(DM dmCoarse, DM dmFine, Mat *mass)
42183c47955SMatthew G. Knepley {
422895a1698SMatthew G. Knepley   PetscSection   gsf;
42383c47955SMatthew G. Knepley   PetscInt       m, n;
42483c47955SMatthew G. Knepley   void          *ctx;
42583c47955SMatthew G. Knepley 
42683c47955SMatthew G. Knepley   PetscFunctionBegin;
4279566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalSection(dmFine, &gsf));
4289566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetConstrainedStorageSize(gsf, &m));
4299566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetLocalSize(dmCoarse, &n));
4309566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject) dmCoarse), mass));
4319566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*mass, n, m, PETSC_DETERMINE, PETSC_DETERMINE));
4329566063dSJacob Faibussowitsch   PetscCall(MatSetType(*mass, dmCoarse->mattype));
4339566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(dmFine, &ctx));
43483c47955SMatthew G. Knepley 
4359566063dSJacob Faibussowitsch   PetscCall(DMSwarmComputeMassMatrix_Private(dmCoarse, dmFine, *mass, PETSC_TRUE, ctx));
4369566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(*mass, NULL, "-mass_mat_view"));
43783c47955SMatthew G. Knepley   PetscFunctionReturn(0);
43883c47955SMatthew G. Knepley }
43983c47955SMatthew G. Knepley 
4404cc7f7b2SMatthew G. Knepley static PetscErrorCode DMSwarmComputeMassMatrixSquare_Private(DM dmc, DM dmf, Mat mass, PetscBool useDeltaFunction, void *ctx)
4414cc7f7b2SMatthew G. Knepley {
4424cc7f7b2SMatthew G. Knepley   const char    *name = "Mass Matrix Square";
4434cc7f7b2SMatthew G. Knepley   MPI_Comm       comm;
4444cc7f7b2SMatthew G. Knepley   PetscDS        prob;
4454cc7f7b2SMatthew G. Knepley   PetscSection   fsection, globalFSection;
4464cc7f7b2SMatthew G. Knepley   PetscHSetIJ    ht;
4474cc7f7b2SMatthew G. Knepley   PetscLayout    rLayout, colLayout;
4484cc7f7b2SMatthew G. Knepley   PetscInt      *dnz, *onz, *adj, depth, maxConeSize, maxSupportSize, maxAdjSize;
4494cc7f7b2SMatthew G. Knepley   PetscInt       locRows, locCols, rStart, colStart, colEnd, *rowIDXs;
4504cc7f7b2SMatthew G. Knepley   PetscReal     *xi, *v0, *J, *invJ, detJ = 1.0, v0ref[3] = {-1.0, -1.0, -1.0};
4514cc7f7b2SMatthew G. Knepley   PetscScalar   *elemMat, *elemMatSq;
4524cc7f7b2SMatthew G. Knepley   PetscInt       cdim, Nf, field, cStart, cEnd, cell, totDim, maxC = 0;
4534cc7f7b2SMatthew G. Knepley 
4544cc7f7b2SMatthew G. Knepley   PetscFunctionBegin;
4559566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject) mass, &comm));
4569566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dmf, &cdim));
4579566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dmf, &prob));
4589566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(prob, &Nf));
4599566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(prob, &totDim));
4609566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(cdim, &v0, cdim*cdim, &J, cdim*cdim,&invJ));
4619566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dmf, &fsection));
4629566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalSection(dmf, &globalFSection));
4639566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd));
4649566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(mass, &locRows, &locCols));
4654cc7f7b2SMatthew G. Knepley 
4669566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreate(comm, &colLayout));
4679566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetLocalSize(colLayout, locCols));
4689566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(colLayout, 1));
4699566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(colLayout));
4709566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetRange(colLayout, &colStart, &colEnd));
4719566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&colLayout));
4724cc7f7b2SMatthew G. Knepley 
4739566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreate(comm, &rLayout));
4749566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetLocalSize(rLayout, locRows));
4759566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(rLayout, 1));
4769566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(rLayout));
4779566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetRange(rLayout, &rStart, NULL));
4789566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&rLayout));
4794cc7f7b2SMatthew G. Knepley 
4809566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(dmf, &depth));
4819566063dSJacob Faibussowitsch   PetscCall(DMPlexGetMaxSizes(dmf, &maxConeSize, &maxSupportSize));
4824cc7f7b2SMatthew G. Knepley   maxAdjSize = PetscPowInt(maxConeSize*maxSupportSize, depth);
4839566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(maxAdjSize, &adj));
4844cc7f7b2SMatthew G. Knepley 
4859566063dSJacob Faibussowitsch   PetscCall(PetscCalloc2(locRows, &dnz, locRows, &onz));
4869566063dSJacob Faibussowitsch   PetscCall(PetscHSetIJCreate(&ht));
4874cc7f7b2SMatthew G. Knepley   /* Count nonzeros
4884cc7f7b2SMatthew G. Knepley        This is just FVM++, but we cannot use the Plex P0 allocation since unknowns in a cell will not be contiguous
4894cc7f7b2SMatthew G. Knepley   */
4909566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortGetAccess(dmc));
4914cc7f7b2SMatthew G. Knepley   for (cell = cStart; cell < cEnd; ++cell) {
4924cc7f7b2SMatthew G. Knepley     PetscInt  i;
4934cc7f7b2SMatthew G. Knepley     PetscInt *cindices;
4944cc7f7b2SMatthew G. Knepley     PetscInt  numCIndices;
4954cc7f7b2SMatthew G. Knepley   #if 0
4964cc7f7b2SMatthew G. Knepley     PetscInt  adjSize = maxAdjSize, a, j;
4974cc7f7b2SMatthew G. Knepley   #endif
4984cc7f7b2SMatthew G. Knepley 
4999566063dSJacob Faibussowitsch     PetscCall(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices));
5004cc7f7b2SMatthew G. Knepley     maxC = PetscMax(maxC, numCIndices);
5014cc7f7b2SMatthew G. Knepley     /* Diagonal block */
5024cc7f7b2SMatthew G. Knepley     for (i = 0; i < numCIndices; ++i) {dnz[cindices[i]] += numCIndices;}
5034cc7f7b2SMatthew G. Knepley #if 0
5044cc7f7b2SMatthew G. Knepley     /* Off-diagonal blocks */
5059566063dSJacob Faibussowitsch     PetscCall(DMPlexGetAdjacency(dmf, cell, &adjSize, &adj));
5064cc7f7b2SMatthew G. Knepley     for (a = 0; a < adjSize; ++a) {
5074cc7f7b2SMatthew G. Knepley       if (adj[a] >= cStart && adj[a] < cEnd && adj[a] != cell) {
5084cc7f7b2SMatthew G. Knepley         const PetscInt ncell = adj[a];
5094cc7f7b2SMatthew G. Knepley         PetscInt      *ncindices;
5104cc7f7b2SMatthew G. Knepley         PetscInt       numNCIndices;
5114cc7f7b2SMatthew G. Knepley 
5129566063dSJacob Faibussowitsch         PetscCall(DMSwarmSortGetPointsPerCell(dmc, ncell, &numNCIndices, &ncindices));
5134cc7f7b2SMatthew G. Knepley         {
5144cc7f7b2SMatthew G. Knepley           PetscHashIJKey key;
5154cc7f7b2SMatthew G. Knepley           PetscBool      missing;
5164cc7f7b2SMatthew G. Knepley 
5174cc7f7b2SMatthew G. Knepley           for (i = 0; i < numCIndices; ++i) {
5184cc7f7b2SMatthew G. Knepley             key.i = cindices[i] + rStart; /* global rows (from Swarm) */
5194cc7f7b2SMatthew G. Knepley             if (key.i < 0) continue;
5204cc7f7b2SMatthew G. Knepley             for (j = 0; j < numNCIndices; ++j) {
5214cc7f7b2SMatthew G. Knepley               key.j = ncindices[j] + rStart; /* global column (from Swarm) */
5224cc7f7b2SMatthew G. Knepley               if (key.j < 0) continue;
5239566063dSJacob Faibussowitsch               PetscCall(PetscHSetIJQueryAdd(ht, key, &missing));
5244cc7f7b2SMatthew G. Knepley               if (missing) {
5254cc7f7b2SMatthew G. Knepley                 if ((key.j >= colStart) && (key.j < colEnd)) ++dnz[key.i - rStart];
5264cc7f7b2SMatthew G. Knepley                 else                                         ++onz[key.i - rStart];
5274cc7f7b2SMatthew G. Knepley               }
5284cc7f7b2SMatthew G. Knepley             }
5294cc7f7b2SMatthew G. Knepley           }
5304cc7f7b2SMatthew G. Knepley         }
5319566063dSJacob Faibussowitsch         PetscCall(PetscFree(ncindices));
5324cc7f7b2SMatthew G. Knepley       }
5334cc7f7b2SMatthew G. Knepley     }
5344cc7f7b2SMatthew G. Knepley #endif
5359566063dSJacob Faibussowitsch     PetscCall(PetscFree(cindices));
5364cc7f7b2SMatthew G. Knepley   }
5379566063dSJacob Faibussowitsch   PetscCall(PetscHSetIJDestroy(&ht));
5389566063dSJacob Faibussowitsch   PetscCall(MatXAIJSetPreallocation(mass, 1, dnz, onz, NULL, NULL));
5399566063dSJacob Faibussowitsch   PetscCall(MatSetOption(mass, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
5409566063dSJacob Faibussowitsch   PetscCall(PetscFree2(dnz, onz));
5419566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(maxC*totDim, &elemMat, maxC*maxC, &elemMatSq, maxC, &rowIDXs, maxC*cdim, &xi));
5424cc7f7b2SMatthew G. Knepley   /* Fill in values
5434cc7f7b2SMatthew G. Knepley        Each entry is a sum of terms \phi_i(x_p) \phi_i(x_q)
5444cc7f7b2SMatthew G. Knepley        Start just by producing block diagonal
5454cc7f7b2SMatthew G. Knepley        Could loop over adjacent cells
5464cc7f7b2SMatthew G. Knepley          Produce neighboring element matrix
5474cc7f7b2SMatthew G. Knepley          TODO Determine which columns and rows correspond to shared dual vector
5484cc7f7b2SMatthew G. Knepley          Do MatMatMult with rectangular matrices
5494cc7f7b2SMatthew G. Knepley          Insert block
5504cc7f7b2SMatthew G. Knepley   */
5514cc7f7b2SMatthew G. Knepley   for (field = 0; field < Nf; ++field) {
5524cc7f7b2SMatthew G. Knepley     PetscTabulation Tcoarse;
5534cc7f7b2SMatthew G. Knepley     PetscObject     obj;
5544cc7f7b2SMatthew G. Knepley     PetscReal       *coords;
5554cc7f7b2SMatthew G. Knepley     PetscInt        Nc, i;
5564cc7f7b2SMatthew G. Knepley 
5579566063dSJacob Faibussowitsch     PetscCall(PetscDSGetDiscretization(prob, field, &obj));
5589566063dSJacob Faibussowitsch     PetscCall(PetscFEGetNumComponents((PetscFE) obj, &Nc));
55908401ef6SPierre Jolivet     PetscCheck(Nc == 1,PetscObjectComm((PetscObject) dmf), PETSC_ERR_SUP, "Can only interpolate a scalar field from particles, Nc = %D", Nc);
5609566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **) &coords));
5614cc7f7b2SMatthew G. Knepley     for (cell = cStart; cell < cEnd; ++cell) {
5624cc7f7b2SMatthew G. Knepley       PetscInt *findices  , *cindices;
5634cc7f7b2SMatthew G. Knepley       PetscInt  numFIndices, numCIndices;
5644cc7f7b2SMatthew G. Knepley       PetscInt  p, c;
5654cc7f7b2SMatthew G. Knepley 
5664cc7f7b2SMatthew G. Knepley       /* TODO: Use DMField instead of assuming affine */
5679566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ));
5689566063dSJacob Faibussowitsch       PetscCall(DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
5699566063dSJacob Faibussowitsch       PetscCall(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices));
5704cc7f7b2SMatthew G. Knepley       for (p = 0; p < numCIndices; ++p) {
5714cc7f7b2SMatthew G. Knepley         CoordinatesRealToRef(cdim, cdim, v0ref, v0, invJ, &coords[cindices[p]*cdim], &xi[p*cdim]);
5724cc7f7b2SMatthew G. Knepley       }
5739566063dSJacob Faibussowitsch       PetscCall(PetscFECreateTabulation((PetscFE) obj, 1, numCIndices, xi, 0, &Tcoarse));
5744cc7f7b2SMatthew G. Knepley       /* Get elemMat entries by multiplying by weight */
5759566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(elemMat, numCIndices*totDim));
5764cc7f7b2SMatthew G. Knepley       for (i = 0; i < numFIndices; ++i) {
5774cc7f7b2SMatthew G. Knepley         for (p = 0; p < numCIndices; ++p) {
5784cc7f7b2SMatthew G. Knepley           for (c = 0; c < Nc; ++c) {
5794cc7f7b2SMatthew G. Knepley             /* B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c */
5804cc7f7b2SMatthew G. Knepley             elemMat[p*numFIndices+i] += Tcoarse->T[0][(p*numFIndices + i)*Nc + c]*(useDeltaFunction ? 1.0 : detJ);
5814cc7f7b2SMatthew G. Knepley           }
5824cc7f7b2SMatthew G. Knepley         }
5834cc7f7b2SMatthew G. Knepley       }
5849566063dSJacob Faibussowitsch       PetscCall(PetscTabulationDestroy(&Tcoarse));
5854cc7f7b2SMatthew G. Knepley       for (p = 0; p < numCIndices; ++p) rowIDXs[p] = cindices[p] + rStart;
5869566063dSJacob Faibussowitsch       if (0) PetscCall(DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat));
5874cc7f7b2SMatthew G. Knepley       /* Block diagonal */
58878884ca7SMatthew G. Knepley       if (numCIndices) {
5894cc7f7b2SMatthew G. Knepley         PetscBLASInt blasn, blask;
5904cc7f7b2SMatthew G. Knepley         PetscScalar  one = 1.0, zero = 0.0;
5914cc7f7b2SMatthew G. Knepley 
5929566063dSJacob Faibussowitsch         PetscCall(PetscBLASIntCast(numCIndices, &blasn));
5939566063dSJacob Faibussowitsch         PetscCall(PetscBLASIntCast(numFIndices, &blask));
5944cc7f7b2SMatthew G. Knepley         PetscStackCallBLAS("BLASgemm",BLASgemm_("T","N",&blasn,&blasn,&blask,&one,elemMat,&blask,elemMat,&blask,&zero,elemMatSq,&blasn));
5954cc7f7b2SMatthew G. Knepley       }
5969566063dSJacob Faibussowitsch       PetscCall(MatSetValues(mass, numCIndices, rowIDXs, numCIndices, rowIDXs, elemMatSq, ADD_VALUES));
5974cc7f7b2SMatthew G. Knepley       /* TODO Off-diagonal */
5989566063dSJacob Faibussowitsch       PetscCall(PetscFree(cindices));
5999566063dSJacob Faibussowitsch       PetscCall(DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
6004cc7f7b2SMatthew G. Knepley     }
6019566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **) &coords));
6024cc7f7b2SMatthew G. Knepley   }
6039566063dSJacob Faibussowitsch   PetscCall(PetscFree4(elemMat, elemMatSq, rowIDXs, xi));
6049566063dSJacob Faibussowitsch   PetscCall(PetscFree(adj));
6059566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortRestoreAccess(dmc));
6069566063dSJacob Faibussowitsch   PetscCall(PetscFree3(v0, J, invJ));
6079566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mass, MAT_FINAL_ASSEMBLY));
6089566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mass, MAT_FINAL_ASSEMBLY));
6094cc7f7b2SMatthew G. Knepley   PetscFunctionReturn(0);
6104cc7f7b2SMatthew G. Knepley }
6114cc7f7b2SMatthew G. Knepley 
6124cc7f7b2SMatthew G. Knepley /*@C
6134cc7f7b2SMatthew G. Knepley   DMSwarmCreateMassMatrixSquare - Creates the block-diagonal of the square, M^T_p M_p, of the particle mass matrix M_p
6144cc7f7b2SMatthew G. Knepley 
6154cc7f7b2SMatthew G. Knepley   Collective on dmCoarse
6164cc7f7b2SMatthew G. Knepley 
6174cc7f7b2SMatthew G. Knepley   Input parameters:
6184cc7f7b2SMatthew G. Knepley + dmCoarse - a DMSwarm
6194cc7f7b2SMatthew G. Knepley - dmFine   - a DMPlex
6204cc7f7b2SMatthew G. Knepley 
6214cc7f7b2SMatthew G. Knepley   Output parameter:
6224cc7f7b2SMatthew G. Knepley . mass     - the square of the particle mass matrix
6234cc7f7b2SMatthew G. Knepley 
6244cc7f7b2SMatthew G. Knepley   Level: advanced
6254cc7f7b2SMatthew G. Knepley 
6264cc7f7b2SMatthew G. Knepley   Notes:
6274cc7f7b2SMatthew G. Knepley   We only compute the block diagonal since this provides a good preconditioner and is completely local. It would be possible in the
6284cc7f7b2SMatthew G. Knepley   future to compute the full normal equations.
6294cc7f7b2SMatthew G. Knepley 
6304cc7f7b2SMatthew G. Knepley .seealso: DMCreateMassMatrix()
6314cc7f7b2SMatthew G. Knepley @*/
6324cc7f7b2SMatthew G. Knepley PetscErrorCode DMSwarmCreateMassMatrixSquare(DM dmCoarse, DM dmFine, Mat *mass)
6334cc7f7b2SMatthew G. Knepley {
6344cc7f7b2SMatthew G. Knepley   PetscInt       n;
6354cc7f7b2SMatthew G. Knepley   void          *ctx;
6364cc7f7b2SMatthew G. Knepley 
6374cc7f7b2SMatthew G. Knepley   PetscFunctionBegin;
6389566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetLocalSize(dmCoarse, &n));
6399566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject) dmCoarse), mass));
6409566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*mass, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
6419566063dSJacob Faibussowitsch   PetscCall(MatSetType(*mass, dmCoarse->mattype));
6429566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(dmFine, &ctx));
6434cc7f7b2SMatthew G. Knepley 
6449566063dSJacob Faibussowitsch   PetscCall(DMSwarmComputeMassMatrixSquare_Private(dmCoarse, dmFine, *mass, PETSC_TRUE, ctx));
6459566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(*mass, NULL, "-mass_sq_mat_view"));
6464cc7f7b2SMatthew G. Knepley   PetscFunctionReturn(0);
6474cc7f7b2SMatthew G. Knepley }
6484cc7f7b2SMatthew G. Knepley 
649fb1bcc12SMatthew G. Knepley /*@C
650d3a51819SDave May    DMSwarmCreateGlobalVectorFromField - Creates a Vec object sharing the array associated with a given field
651d3a51819SDave May 
652d083f849SBarry Smith    Collective on dm
653d3a51819SDave May 
654d3a51819SDave May    Input parameters:
65562741f57SDave May +  dm - a DMSwarm
65662741f57SDave May -  fieldname - the textual name given to a registered field
657d3a51819SDave May 
6588b8a3813SDave May    Output parameter:
659d3a51819SDave May .  vec - the vector
660d3a51819SDave May 
661d3a51819SDave May    Level: beginner
662d3a51819SDave May 
6638b8a3813SDave May    Notes:
6648b8a3813SDave May    The vector must be returned using a matching call to DMSwarmDestroyGlobalVectorFromField().
6658b8a3813SDave May 
6668b8a3813SDave May .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmDestroyGlobalVectorFromField()
667d3a51819SDave May @*/
66874d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmCreateGlobalVectorFromField(DM dm,const char fieldname[],Vec *vec)
669b5bcf523SDave May {
670fb1bcc12SMatthew G. Knepley   MPI_Comm       comm = PetscObjectComm((PetscObject) dm);
671b5bcf523SDave May 
672fb1bcc12SMatthew G. Knepley   PetscFunctionBegin;
6739566063dSJacob Faibussowitsch   PetscCall(DMSwarmCreateVectorFromField_Private(dm, fieldname, comm, vec));
674b5bcf523SDave May   PetscFunctionReturn(0);
675b5bcf523SDave May }
676b5bcf523SDave May 
677d3a51819SDave May /*@C
678d3a51819SDave May    DMSwarmDestroyGlobalVectorFromField - Destroys the Vec object which share the array associated with a given field
679d3a51819SDave May 
680d083f849SBarry Smith    Collective on dm
681d3a51819SDave May 
682d3a51819SDave May    Input parameters:
68362741f57SDave May +  dm - a DMSwarm
68462741f57SDave May -  fieldname - the textual name given to a registered field
685d3a51819SDave May 
6868b8a3813SDave May    Output parameter:
687d3a51819SDave May .  vec - the vector
688d3a51819SDave May 
689d3a51819SDave May    Level: beginner
690d3a51819SDave May 
6918b8a3813SDave May .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmCreateGlobalVectorFromField()
692d3a51819SDave May @*/
69374d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmDestroyGlobalVectorFromField(DM dm,const char fieldname[],Vec *vec)
694b5bcf523SDave May {
695fb1bcc12SMatthew G. Knepley   PetscFunctionBegin;
6969566063dSJacob Faibussowitsch   PetscCall(DMSwarmDestroyVectorFromField_Private(dm, fieldname, vec));
697b5bcf523SDave May   PetscFunctionReturn(0);
698b5bcf523SDave May }
699b5bcf523SDave May 
700fb1bcc12SMatthew G. Knepley /*@C
701fb1bcc12SMatthew G. Knepley    DMSwarmCreateLocalVectorFromField - Creates a Vec object sharing the array associated with a given field
702fb1bcc12SMatthew G. Knepley 
703d083f849SBarry Smith    Collective on dm
704fb1bcc12SMatthew G. Knepley 
705fb1bcc12SMatthew G. Knepley    Input parameters:
70662741f57SDave May +  dm - a DMSwarm
70762741f57SDave May -  fieldname - the textual name given to a registered field
708fb1bcc12SMatthew G. Knepley 
7098b8a3813SDave May    Output parameter:
710fb1bcc12SMatthew G. Knepley .  vec - the vector
711fb1bcc12SMatthew G. Knepley 
712fb1bcc12SMatthew G. Knepley    Level: beginner
713fb1bcc12SMatthew G. Knepley 
7148b8a3813SDave May    Notes:
7158b8a3813SDave May    The vector must be returned using a matching call to DMSwarmDestroyLocalVectorFromField().
7168b8a3813SDave May 
7178b8a3813SDave May .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmDestroyLocalVectorFromField()
718fb1bcc12SMatthew G. Knepley @*/
71974d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmCreateLocalVectorFromField(DM dm,const char fieldname[],Vec *vec)
720bbe8250bSMatthew G. Knepley {
721fb1bcc12SMatthew G. Knepley   MPI_Comm       comm = PETSC_COMM_SELF;
722bbe8250bSMatthew G. Knepley 
723fb1bcc12SMatthew G. Knepley   PetscFunctionBegin;
7249566063dSJacob Faibussowitsch   PetscCall(DMSwarmCreateVectorFromField_Private(dm, fieldname, comm, vec));
725fb1bcc12SMatthew G. Knepley   PetscFunctionReturn(0);
726bbe8250bSMatthew G. Knepley }
727fb1bcc12SMatthew G. Knepley 
728fb1bcc12SMatthew G. Knepley /*@C
729fb1bcc12SMatthew G. Knepley    DMSwarmDestroyLocalVectorFromField - Destroys the Vec object which share the array associated with a given field
730fb1bcc12SMatthew G. Knepley 
731d083f849SBarry Smith    Collective on dm
732fb1bcc12SMatthew G. Knepley 
733fb1bcc12SMatthew G. Knepley    Input parameters:
73462741f57SDave May +  dm - a DMSwarm
73562741f57SDave May -  fieldname - the textual name given to a registered field
736fb1bcc12SMatthew G. Knepley 
7378b8a3813SDave May    Output parameter:
738fb1bcc12SMatthew G. Knepley .  vec - the vector
739fb1bcc12SMatthew G. Knepley 
740fb1bcc12SMatthew G. Knepley    Level: beginner
741fb1bcc12SMatthew G. Knepley 
7428b8a3813SDave May .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmCreateLocalVectorFromField()
743fb1bcc12SMatthew G. Knepley @*/
74474d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmDestroyLocalVectorFromField(DM dm,const char fieldname[],Vec *vec)
745fb1bcc12SMatthew G. Knepley {
746fb1bcc12SMatthew G. Knepley   PetscFunctionBegin;
7479566063dSJacob Faibussowitsch   PetscCall(DMSwarmDestroyVectorFromField_Private(dm, fieldname, vec));
748bbe8250bSMatthew G. Knepley   PetscFunctionReturn(0);
749bbe8250bSMatthew G. Knepley }
750bbe8250bSMatthew G. Knepley 
751d3a51819SDave May /*@C
752d3a51819SDave May    DMSwarmInitializeFieldRegister - Initiates the registration of fields to a DMSwarm
753d3a51819SDave May 
754d083f849SBarry Smith    Collective on dm
755d3a51819SDave May 
756d3a51819SDave May    Input parameter:
757d3a51819SDave May .  dm - a DMSwarm
758d3a51819SDave May 
759d3a51819SDave May    Level: beginner
760d3a51819SDave May 
761d3a51819SDave May    Notes:
7628b8a3813SDave May    After all fields have been registered, you must call DMSwarmFinalizeFieldRegister().
763d3a51819SDave May 
764d3a51819SDave May .seealso: DMSwarmFinalizeFieldRegister(), DMSwarmRegisterPetscDatatypeField(),
765d3a51819SDave May           DMSwarmRegisterUserStructField(), DMSwarmRegisterUserDatatypeField()
766d3a51819SDave May @*/
76774d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmInitializeFieldRegister(DM dm)
7685f50eb2eSDave May {
7695f50eb2eSDave May   DM_Swarm      *swarm = (DM_Swarm *) dm->data;
7703454631fSDave May 
771521f74f9SMatthew G. Knepley   PetscFunctionBegin;
772cc651181SDave May   if (!swarm->field_registration_initialized) {
7735f50eb2eSDave May     swarm->field_registration_initialized = PETSC_TRUE;
7749566063dSJacob Faibussowitsch     PetscCall(DMSwarmRegisterPetscDatatypeField(dm,DMSwarmField_pid,1,PETSC_INT64)); /* unique identifer */
7759566063dSJacob Faibussowitsch     PetscCall(DMSwarmRegisterPetscDatatypeField(dm,DMSwarmField_rank,1,PETSC_INT)); /* used for communication */
776cc651181SDave May   }
7775f50eb2eSDave May   PetscFunctionReturn(0);
7785f50eb2eSDave May }
7795f50eb2eSDave May 
78074d0cae8SMatthew G. Knepley /*@
781d3a51819SDave May    DMSwarmFinalizeFieldRegister - Finalizes the registration of fields to a DMSwarm
782d3a51819SDave May 
783d083f849SBarry Smith    Collective on dm
784d3a51819SDave May 
785d3a51819SDave May    Input parameter:
786d3a51819SDave May .  dm - a DMSwarm
787d3a51819SDave May 
788d3a51819SDave May    Level: beginner
789d3a51819SDave May 
790d3a51819SDave May    Notes:
79162741f57SDave May    After DMSwarmFinalizeFieldRegister() has been called, no new fields can be defined on the DMSwarm.
792d3a51819SDave May 
793d3a51819SDave May .seealso: DMSwarmInitializeFieldRegister(), DMSwarmRegisterPetscDatatypeField(),
794d3a51819SDave May           DMSwarmRegisterUserStructField(), DMSwarmRegisterUserDatatypeField()
795d3a51819SDave May @*/
79674d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmFinalizeFieldRegister(DM dm)
7975f50eb2eSDave May {
7985f50eb2eSDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
7996845f8f5SDave May 
800521f74f9SMatthew G. Knepley   PetscFunctionBegin;
801f0cdbbbaSDave May   if (!swarm->field_registration_finalized) {
8029566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketFinalize(swarm->db));
803f0cdbbbaSDave May   }
804f0cdbbbaSDave May   swarm->field_registration_finalized = PETSC_TRUE;
8055f50eb2eSDave May   PetscFunctionReturn(0);
8065f50eb2eSDave May }
8075f50eb2eSDave May 
80874d0cae8SMatthew G. Knepley /*@
809d3a51819SDave May    DMSwarmSetLocalSizes - Sets the length of all registered fields on the DMSwarm
810d3a51819SDave May 
811d3a51819SDave May    Not collective
812d3a51819SDave May 
813d3a51819SDave May    Input parameters:
81462741f57SDave May +  dm - a DMSwarm
815d3a51819SDave May .  nlocal - the length of each registered field
81662741f57SDave May -  buffer - the length of the buffer used to efficient dynamic re-sizing
817d3a51819SDave May 
818d3a51819SDave May    Level: beginner
819d3a51819SDave May 
820d3a51819SDave May .seealso: DMSwarmGetLocalSize()
821d3a51819SDave May @*/
82274d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmSetLocalSizes(DM dm,PetscInt nlocal,PetscInt buffer)
8235f50eb2eSDave May {
8245f50eb2eSDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
8255f50eb2eSDave May 
826521f74f9SMatthew G. Knepley   PetscFunctionBegin;
8279566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMSWARM_SetSizes,0,0,0,0));
8289566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketSetSizes(swarm->db,nlocal,buffer));
8299566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMSWARM_SetSizes,0,0,0,0));
8305f50eb2eSDave May   PetscFunctionReturn(0);
8315f50eb2eSDave May }
8325f50eb2eSDave May 
83374d0cae8SMatthew G. Knepley /*@
834d3a51819SDave May    DMSwarmSetCellDM - Attachs a DM to a DMSwarm
835d3a51819SDave May 
836d083f849SBarry Smith    Collective on dm
837d3a51819SDave May 
838d3a51819SDave May    Input parameters:
83962741f57SDave May +  dm - a DMSwarm
84062741f57SDave May -  dmcell - the DM to attach to the DMSwarm
841d3a51819SDave May 
842d3a51819SDave May    Level: beginner
843d3a51819SDave May 
844d3a51819SDave May    Notes:
845d3a51819SDave May    The attached DM (dmcell) will be queried for point location and
8468b8a3813SDave May    neighbor MPI-rank information if DMSwarmMigrate() is called.
847d3a51819SDave May 
8488b8a3813SDave May .seealso: DMSwarmSetType(), DMSwarmGetCellDM(), DMSwarmMigrate()
849d3a51819SDave May @*/
85074d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmSetCellDM(DM dm,DM dmcell)
851b16650c8SDave May {
852b16650c8SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
853521f74f9SMatthew G. Knepley 
854521f74f9SMatthew G. Knepley   PetscFunctionBegin;
855b16650c8SDave May   swarm->dmcell = dmcell;
856b16650c8SDave May   PetscFunctionReturn(0);
857b16650c8SDave May }
858b16650c8SDave May 
85974d0cae8SMatthew G. Knepley /*@
860d3a51819SDave May    DMSwarmGetCellDM - Fetches the attached cell DM
861d3a51819SDave May 
862d083f849SBarry Smith    Collective on dm
863d3a51819SDave May 
864d3a51819SDave May    Input parameter:
865d3a51819SDave May .  dm - a DMSwarm
866d3a51819SDave May 
867d3a51819SDave May    Output parameter:
868d3a51819SDave May .  dmcell - the DM which was attached to the DMSwarm
869d3a51819SDave May 
870d3a51819SDave May    Level: beginner
871d3a51819SDave May 
872d3a51819SDave May .seealso: DMSwarmSetCellDM()
873d3a51819SDave May @*/
87474d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmGetCellDM(DM dm,DM *dmcell)
875fe39f135SDave May {
876fe39f135SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
877521f74f9SMatthew G. Knepley 
878521f74f9SMatthew G. Knepley   PetscFunctionBegin;
879fe39f135SDave May   *dmcell = swarm->dmcell;
880fe39f135SDave May   PetscFunctionReturn(0);
881fe39f135SDave May }
882fe39f135SDave May 
88374d0cae8SMatthew G. Knepley /*@
884d3a51819SDave May    DMSwarmGetLocalSize - Retrives the local length of fields registered
885d3a51819SDave May 
886d3a51819SDave May    Not collective
887d3a51819SDave May 
888d3a51819SDave May    Input parameter:
889d3a51819SDave May .  dm - a DMSwarm
890d3a51819SDave May 
891d3a51819SDave May    Output parameter:
892d3a51819SDave May .  nlocal - the length of each registered field
893d3a51819SDave May 
894d3a51819SDave May    Level: beginner
895d3a51819SDave May 
8968b8a3813SDave May .seealso: DMSwarmGetSize(), DMSwarmSetLocalSizes()
897d3a51819SDave May @*/
89874d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmGetLocalSize(DM dm,PetscInt *nlocal)
899dcf43ee8SDave May {
900dcf43ee8SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
901dcf43ee8SDave May 
902521f74f9SMatthew G. Knepley   PetscFunctionBegin;
9039566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db,nlocal,NULL,NULL));
904dcf43ee8SDave May   PetscFunctionReturn(0);
905dcf43ee8SDave May }
906dcf43ee8SDave May 
90774d0cae8SMatthew G. Knepley /*@
908d3a51819SDave May    DMSwarmGetSize - Retrives the total length of fields registered
909d3a51819SDave May 
910d083f849SBarry Smith    Collective on dm
911d3a51819SDave May 
912d3a51819SDave May    Input parameter:
913d3a51819SDave May .  dm - a DMSwarm
914d3a51819SDave May 
915d3a51819SDave May    Output parameter:
916d3a51819SDave May .  n - the total length of each registered field
917d3a51819SDave May 
918d3a51819SDave May    Level: beginner
919d3a51819SDave May 
920d3a51819SDave May    Note:
921d3a51819SDave May    This calls MPI_Allreduce upon each call (inefficient but safe)
922d3a51819SDave May 
9238b8a3813SDave May .seealso: DMSwarmGetLocalSize(), DMSwarmSetLocalSizes()
924d3a51819SDave May @*/
92574d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmGetSize(DM dm,PetscInt *n)
926dcf43ee8SDave May {
927dcf43ee8SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
9285627991aSBarry Smith   PetscInt       nlocal;
929dcf43ee8SDave May 
930521f74f9SMatthew G. Knepley   PetscFunctionBegin;
9319566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db,&nlocal,NULL,NULL));
9329566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allreduce(&nlocal,n,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm)));
933dcf43ee8SDave May   PetscFunctionReturn(0);
934dcf43ee8SDave May }
935dcf43ee8SDave May 
936d3a51819SDave May /*@C
9378b8a3813SDave May    DMSwarmRegisterPetscDatatypeField - Register a field to a DMSwarm with a native PETSc data type
938d3a51819SDave May 
939d083f849SBarry Smith    Collective on dm
940d3a51819SDave May 
941d3a51819SDave May    Input parameters:
94262741f57SDave May +  dm - a DMSwarm
943d3a51819SDave May .  fieldname - the textual name to identify this field
944d3a51819SDave May .  blocksize - the number of each data type
94562741f57SDave May -  type - a valid PETSc data type (PETSC_CHAR, PETSC_SHORT, PETSC_INT, PETSC_FLOAT, PETSC_REAL, PETSC_LONG)
946d3a51819SDave May 
947d3a51819SDave May    Level: beginner
948d3a51819SDave May 
949d3a51819SDave May    Notes:
9508b8a3813SDave May    The textual name for each registered field must be unique.
951d3a51819SDave May 
952d3a51819SDave May .seealso: DMSwarmRegisterUserStructField(), DMSwarmRegisterUserDatatypeField()
953d3a51819SDave May @*/
95474d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmRegisterPetscDatatypeField(DM dm,const char fieldname[],PetscInt blocksize,PetscDataType type)
955b62e03f8SDave May {
956b62e03f8SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
957b62e03f8SDave May   size_t         size;
958b62e03f8SDave May 
959521f74f9SMatthew G. Knepley   PetscFunctionBegin;
96028b400f6SJacob Faibussowitsch   PetscCheck(swarm->field_registration_initialized,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMSwarmInitializeFieldRegister() first");
96128b400f6SJacob Faibussowitsch   PetscCheck(!swarm->field_registration_finalized,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Cannot register additional fields after calling DMSwarmFinalizeFieldRegister() first");
9625f50eb2eSDave May 
96308401ef6SPierre Jolivet   PetscCheck(type != PETSC_OBJECT,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
96408401ef6SPierre Jolivet   PetscCheck(type != PETSC_FUNCTION,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
96508401ef6SPierre Jolivet   PetscCheck(type != PETSC_STRING,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
96608401ef6SPierre Jolivet   PetscCheck(type != PETSC_STRUCT,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
96708401ef6SPierre Jolivet   PetscCheck(type != PETSC_DATATYPE_UNKNOWN,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
968b62e03f8SDave May 
9699566063dSJacob Faibussowitsch   PetscCall(PetscDataTypeGetSize(type, &size));
970b62e03f8SDave May   /* Load a specific data type into data bucket, specifying textual name and its size in bytes */
9719566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketRegisterField(swarm->db,"DMSwarmRegisterPetscDatatypeField",fieldname,blocksize*size,NULL));
97252c3ed93SDave May   {
97377048351SPatrick Sanan     DMSwarmDataField gfield;
97452c3ed93SDave May 
9759566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,fieldname,&gfield));
9769566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataFieldSetBlockSize(gfield,blocksize));
97752c3ed93SDave May   }
978b62e03f8SDave May   swarm->db->field[swarm->db->nfields-1]->petsc_type = type;
979b62e03f8SDave May   PetscFunctionReturn(0);
980b62e03f8SDave May }
981b62e03f8SDave May 
982d3a51819SDave May /*@C
983d3a51819SDave May    DMSwarmRegisterUserStructField - Register a user defined struct to a DMSwarm
984d3a51819SDave May 
985d083f849SBarry Smith    Collective on dm
986d3a51819SDave May 
987d3a51819SDave May    Input parameters:
98862741f57SDave May +  dm - a DMSwarm
989d3a51819SDave May .  fieldname - the textual name to identify this field
99062741f57SDave May -  size - the size in bytes of the user struct of each data type
991d3a51819SDave May 
992d3a51819SDave May    Level: beginner
993d3a51819SDave May 
994d3a51819SDave May    Notes:
9958b8a3813SDave May    The textual name for each registered field must be unique.
996d3a51819SDave May 
997d3a51819SDave May .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmRegisterUserDatatypeField()
998d3a51819SDave May @*/
99974d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmRegisterUserStructField(DM dm,const char fieldname[],size_t size)
1000b62e03f8SDave May {
1001b62e03f8SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
1002b62e03f8SDave May 
1003521f74f9SMatthew G. Knepley   PetscFunctionBegin;
10049566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketRegisterField(swarm->db,"DMSwarmRegisterUserStructField",fieldname,size,NULL));
1005b62e03f8SDave May   swarm->db->field[swarm->db->nfields-1]->petsc_type = PETSC_STRUCT ;
1006b62e03f8SDave May   PetscFunctionReturn(0);
1007b62e03f8SDave May }
1008b62e03f8SDave May 
1009d3a51819SDave May /*@C
1010d3a51819SDave May    DMSwarmRegisterUserDatatypeField - Register a user defined data type to a DMSwarm
1011d3a51819SDave May 
1012d083f849SBarry Smith    Collective on dm
1013d3a51819SDave May 
1014d3a51819SDave May    Input parameters:
101562741f57SDave May +  dm - a DMSwarm
1016d3a51819SDave May .  fieldname - the textual name to identify this field
1017d3a51819SDave May .  size - the size in bytes of the user data type
101862741f57SDave May -  blocksize - the number of each data type
1019d3a51819SDave May 
1020d3a51819SDave May    Level: beginner
1021d3a51819SDave May 
1022d3a51819SDave May    Notes:
10238b8a3813SDave May    The textual name for each registered field must be unique.
1024d3a51819SDave May 
1025d3a51819SDave May .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmRegisterUserStructField(), DMSwarmRegisterUserDatatypeField()
1026d3a51819SDave May @*/
102774d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmRegisterUserDatatypeField(DM dm,const char fieldname[],size_t size,PetscInt blocksize)
1028b62e03f8SDave May {
1029b62e03f8SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
1030b62e03f8SDave May 
1031521f74f9SMatthew G. Knepley   PetscFunctionBegin;
10329566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketRegisterField(swarm->db,"DMSwarmRegisterUserDatatypeField",fieldname,blocksize*size,NULL));
1033320740a0SDave May   {
103477048351SPatrick Sanan     DMSwarmDataField gfield;
1035320740a0SDave May 
10369566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,fieldname,&gfield));
10379566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataFieldSetBlockSize(gfield,blocksize));
1038320740a0SDave May   }
1039b62e03f8SDave May   swarm->db->field[swarm->db->nfields-1]->petsc_type = PETSC_DATATYPE_UNKNOWN;
1040b62e03f8SDave May   PetscFunctionReturn(0);
1041b62e03f8SDave May }
1042b62e03f8SDave May 
1043d3a51819SDave May /*@C
1044d3a51819SDave May    DMSwarmGetField - Get access to the underlying array storing all entries associated with a registered field
1045d3a51819SDave May 
1046d3a51819SDave May    Not collective
1047d3a51819SDave May 
1048d3a51819SDave May    Input parameters:
104962741f57SDave May +  dm - a DMSwarm
105062741f57SDave May -  fieldname - the textual name to identify this field
1051d3a51819SDave May 
1052d3a51819SDave May    Output parameters:
105362741f57SDave May +  blocksize - the number of each data type
1054d3a51819SDave May .  type - the data type
105562741f57SDave May -  data - pointer to raw array
1056d3a51819SDave May 
1057d3a51819SDave May    Level: beginner
1058d3a51819SDave May 
1059d3a51819SDave May    Notes:
10608b8a3813SDave May    The array must be returned using a matching call to DMSwarmRestoreField().
1061d3a51819SDave May 
1062d3a51819SDave May .seealso: DMSwarmRestoreField()
1063d3a51819SDave May @*/
106474d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmGetField(DM dm,const char fieldname[],PetscInt *blocksize,PetscDataType *type,void **data)
1065b62e03f8SDave May {
1066b62e03f8SDave May   DM_Swarm         *swarm = (DM_Swarm*)dm->data;
106777048351SPatrick Sanan   DMSwarmDataField gfield;
1068b62e03f8SDave May 
1069521f74f9SMatthew G. Knepley   PetscFunctionBegin;
10709566063dSJacob Faibussowitsch   if (!swarm->issetup) PetscCall(DMSetUp(dm));
10719566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,fieldname,&gfield));
10729566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataFieldGetAccess(gfield));
10739566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataFieldGetEntries(gfield,data));
10741b1ea282SDave May   if (blocksize) {*blocksize = gfield->bs; }
1075b5bcf523SDave May   if (type) { *type = gfield->petsc_type; }
1076b62e03f8SDave May   PetscFunctionReturn(0);
1077b62e03f8SDave May }
1078b62e03f8SDave May 
1079d3a51819SDave May /*@C
1080d3a51819SDave May    DMSwarmRestoreField - Restore access to the underlying array storing all entries associated with a registered field
1081d3a51819SDave May 
1082d3a51819SDave May    Not collective
1083d3a51819SDave May 
1084d3a51819SDave May    Input parameters:
108562741f57SDave May +  dm - a DMSwarm
108662741f57SDave May -  fieldname - the textual name to identify this field
1087d3a51819SDave May 
1088d3a51819SDave May    Output parameters:
108962741f57SDave May +  blocksize - the number of each data type
1090d3a51819SDave May .  type - the data type
109162741f57SDave May -  data - pointer to raw array
1092d3a51819SDave May 
1093d3a51819SDave May    Level: beginner
1094d3a51819SDave May 
1095d3a51819SDave May    Notes:
10968b8a3813SDave May    The user must call DMSwarmGetField() prior to calling DMSwarmRestoreField().
1097d3a51819SDave May 
1098d3a51819SDave May .seealso: DMSwarmGetField()
1099d3a51819SDave May @*/
110074d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmRestoreField(DM dm,const char fieldname[],PetscInt *blocksize,PetscDataType *type,void **data)
1101b62e03f8SDave May {
1102b62e03f8SDave May   DM_Swarm         *swarm = (DM_Swarm*)dm->data;
110377048351SPatrick Sanan   DMSwarmDataField gfield;
1104b62e03f8SDave May 
1105521f74f9SMatthew G. Knepley   PetscFunctionBegin;
11069566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,fieldname,&gfield));
11079566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataFieldRestoreAccess(gfield));
1108b62e03f8SDave May   if (data) *data = NULL;
1109b62e03f8SDave May   PetscFunctionReturn(0);
1110b62e03f8SDave May }
1111b62e03f8SDave May 
111274d0cae8SMatthew G. Knepley /*@
1113d3a51819SDave May    DMSwarmAddPoint - Add space for one new point in the DMSwarm
1114d3a51819SDave May 
1115d3a51819SDave May    Not collective
1116d3a51819SDave May 
1117d3a51819SDave May    Input parameter:
1118d3a51819SDave May .  dm - a DMSwarm
1119d3a51819SDave May 
1120d3a51819SDave May    Level: beginner
1121d3a51819SDave May 
1122d3a51819SDave May    Notes:
11238b8a3813SDave May    The new point will have all fields initialized to zero.
1124d3a51819SDave May 
1125d3a51819SDave May .seealso: DMSwarmAddNPoints()
1126d3a51819SDave May @*/
112774d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmAddPoint(DM dm)
1128cb1d1399SDave May {
1129cb1d1399SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
1130cb1d1399SDave May 
1131521f74f9SMatthew G. Knepley   PetscFunctionBegin;
11329566063dSJacob Faibussowitsch   if (!swarm->issetup) PetscCall(DMSetUp(dm));
11339566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMSWARM_AddPoints,0,0,0,0));
11349566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketAddPoint(swarm->db));
11359566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMSWARM_AddPoints,0,0,0,0));
1136cb1d1399SDave May   PetscFunctionReturn(0);
1137cb1d1399SDave May }
1138cb1d1399SDave May 
113974d0cae8SMatthew G. Knepley /*@
1140d3a51819SDave May    DMSwarmAddNPoints - Add space for a number of new points in the DMSwarm
1141d3a51819SDave May 
1142d3a51819SDave May    Not collective
1143d3a51819SDave May 
1144d3a51819SDave May    Input parameters:
114562741f57SDave May +  dm - a DMSwarm
114662741f57SDave May -  npoints - the number of new points to add
1147d3a51819SDave May 
1148d3a51819SDave May    Level: beginner
1149d3a51819SDave May 
1150d3a51819SDave May    Notes:
11518b8a3813SDave May    The new point will have all fields initialized to zero.
1152d3a51819SDave May 
1153d3a51819SDave May .seealso: DMSwarmAddPoint()
1154d3a51819SDave May @*/
115574d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmAddNPoints(DM dm,PetscInt npoints)
1156cb1d1399SDave May {
1157cb1d1399SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
1158cb1d1399SDave May   PetscInt       nlocal;
1159cb1d1399SDave May 
1160521f74f9SMatthew G. Knepley   PetscFunctionBegin;
11619566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMSWARM_AddPoints,0,0,0,0));
11629566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db,&nlocal,NULL,NULL));
1163cb1d1399SDave May   nlocal = nlocal + npoints;
11649566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketSetSizes(swarm->db,nlocal,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT));
11659566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMSWARM_AddPoints,0,0,0,0));
1166cb1d1399SDave May   PetscFunctionReturn(0);
1167cb1d1399SDave May }
1168cb1d1399SDave May 
116974d0cae8SMatthew G. Knepley /*@
1170d3a51819SDave May    DMSwarmRemovePoint - Remove the last point from the DMSwarm
1171d3a51819SDave May 
1172d3a51819SDave May    Not collective
1173d3a51819SDave May 
1174d3a51819SDave May    Input parameter:
1175d3a51819SDave May .  dm - a DMSwarm
1176d3a51819SDave May 
1177d3a51819SDave May    Level: beginner
1178d3a51819SDave May 
1179d3a51819SDave May .seealso: DMSwarmRemovePointAtIndex()
1180d3a51819SDave May @*/
118174d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmRemovePoint(DM dm)
1182cb1d1399SDave May {
1183cb1d1399SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
1184cb1d1399SDave May 
1185521f74f9SMatthew G. Knepley   PetscFunctionBegin;
11869566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMSWARM_RemovePoints,0,0,0,0));
11879566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketRemovePoint(swarm->db));
11889566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMSWARM_RemovePoints,0,0,0,0));
1189cb1d1399SDave May   PetscFunctionReturn(0);
1190cb1d1399SDave May }
1191cb1d1399SDave May 
119274d0cae8SMatthew G. Knepley /*@
1193d3a51819SDave May    DMSwarmRemovePointAtIndex - Removes a specific point from the DMSwarm
1194d3a51819SDave May 
1195d3a51819SDave May    Not collective
1196d3a51819SDave May 
1197d3a51819SDave May    Input parameters:
119862741f57SDave May +  dm - a DMSwarm
119962741f57SDave May -  idx - index of point to remove
1200d3a51819SDave May 
1201d3a51819SDave May    Level: beginner
1202d3a51819SDave May 
1203d3a51819SDave May .seealso: DMSwarmRemovePoint()
1204d3a51819SDave May @*/
120574d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmRemovePointAtIndex(DM dm,PetscInt idx)
1206cb1d1399SDave May {
1207cb1d1399SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
1208cb1d1399SDave May 
1209521f74f9SMatthew G. Knepley   PetscFunctionBegin;
12109566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMSWARM_RemovePoints,0,0,0,0));
12119566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketRemovePointAtIndex(swarm->db,idx));
12129566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMSWARM_RemovePoints,0,0,0,0));
1213cb1d1399SDave May   PetscFunctionReturn(0);
1214cb1d1399SDave May }
1215b62e03f8SDave May 
121674d0cae8SMatthew G. Knepley /*@
1217ba4fc9c6SDave May    DMSwarmCopyPoint - Copy point pj to point pi in the DMSwarm
1218ba4fc9c6SDave May 
1219ba4fc9c6SDave May    Not collective
1220ba4fc9c6SDave May 
1221ba4fc9c6SDave May    Input parameters:
1222ba4fc9c6SDave May +  dm - a DMSwarm
1223ba4fc9c6SDave May .  pi - the index of the point to copy
1224ba4fc9c6SDave May -  pj - the point index where the copy should be located
1225ba4fc9c6SDave May 
1226ba4fc9c6SDave May  Level: beginner
1227ba4fc9c6SDave May 
1228ba4fc9c6SDave May .seealso: DMSwarmRemovePoint()
1229ba4fc9c6SDave May @*/
123074d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmCopyPoint(DM dm,PetscInt pi,PetscInt pj)
1231ba4fc9c6SDave May {
1232ba4fc9c6SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
1233ba4fc9c6SDave May 
1234ba4fc9c6SDave May   PetscFunctionBegin;
12359566063dSJacob Faibussowitsch   if (!swarm->issetup) PetscCall(DMSetUp(dm));
12369566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketCopyPoint(swarm->db,pi,swarm->db,pj));
1237ba4fc9c6SDave May   PetscFunctionReturn(0);
1238ba4fc9c6SDave May }
1239ba4fc9c6SDave May 
1240095059a4SDave May PetscErrorCode DMSwarmMigrate_Basic(DM dm,PetscBool remove_sent_points)
12413454631fSDave May {
1242521f74f9SMatthew G. Knepley   PetscFunctionBegin;
12439566063dSJacob Faibussowitsch   PetscCall(DMSwarmMigrate_Push_Basic(dm,remove_sent_points));
12443454631fSDave May   PetscFunctionReturn(0);
12453454631fSDave May }
12463454631fSDave May 
124774d0cae8SMatthew G. Knepley /*@
1248d3a51819SDave May    DMSwarmMigrate - Relocates points defined in the DMSwarm to other MPI-ranks
1249d3a51819SDave May 
1250d083f849SBarry Smith    Collective on dm
1251d3a51819SDave May 
1252d3a51819SDave May    Input parameters:
125362741f57SDave May +  dm - the DMSwarm
125462741f57SDave May -  remove_sent_points - flag indicating if sent points should be removed from the current MPI-rank
1255d3a51819SDave May 
1256d3a51819SDave May    Notes:
1257a5b23f4aSJose E. Roman    The DM will be modified to accommodate received points.
12588b8a3813SDave May    If remove_sent_points = PETSC_TRUE, any points that were sent will be removed from the DM.
12598b8a3813SDave May    Different styles of migration are supported. See DMSwarmSetMigrateType().
1260d3a51819SDave May 
1261d3a51819SDave May    Level: advanced
1262d3a51819SDave May 
1263d3a51819SDave May .seealso: DMSwarmSetMigrateType()
1264d3a51819SDave May @*/
126574d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmMigrate(DM dm,PetscBool remove_sent_points)
12663454631fSDave May {
1267f0cdbbbaSDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
1268f0cdbbbaSDave May 
1269521f74f9SMatthew G. Knepley   PetscFunctionBegin;
12709566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMSWARM_Migrate,0,0,0,0));
1271f0cdbbbaSDave May   switch (swarm->migrate_type) {
1272f0cdbbbaSDave May     case DMSWARM_MIGRATE_BASIC:
12739566063dSJacob Faibussowitsch       PetscCall(DMSwarmMigrate_Basic(dm,remove_sent_points));
1274f0cdbbbaSDave May       break;
1275f0cdbbbaSDave May     case DMSWARM_MIGRATE_DMCELLNSCATTER:
12769566063dSJacob Faibussowitsch       PetscCall(DMSwarmMigrate_CellDMScatter(dm,remove_sent_points));
1277f0cdbbbaSDave May       break;
1278f0cdbbbaSDave May     case DMSWARM_MIGRATE_DMCELLEXACT:
1279f0cdbbbaSDave May       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_MIGRATE_DMCELLEXACT not implemented");
1280f0cdbbbaSDave May     case DMSWARM_MIGRATE_USER:
1281f0cdbbbaSDave May       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_MIGRATE_USER not implemented");
1282f0cdbbbaSDave May     default:
1283f0cdbbbaSDave May       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_MIGRATE type unknown");
1284f0cdbbbaSDave May   }
12859566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMSWARM_Migrate,0,0,0,0));
12869566063dSJacob Faibussowitsch   PetscCall(DMClearGlobalVectors(dm));
12873454631fSDave May   PetscFunctionReturn(0);
12883454631fSDave May }
12893454631fSDave May 
1290f0cdbbbaSDave May PetscErrorCode DMSwarmMigrate_GlobalToLocal_Basic(DM dm,PetscInt *globalsize);
1291f0cdbbbaSDave May 
1292d3a51819SDave May /*
1293d3a51819SDave May  DMSwarmCollectViewCreate
1294d3a51819SDave May 
1295d3a51819SDave May  * Applies a collection method and gathers point neighbour points into dm
1296d3a51819SDave May 
1297d3a51819SDave May  Notes:
12988b8a3813SDave May  Users should call DMSwarmCollectViewDestroy() after
1299d3a51819SDave May  they have finished computations associated with the collected points
1300d3a51819SDave May */
1301d3a51819SDave May 
130274d0cae8SMatthew G. Knepley /*@
1303d3a51819SDave May    DMSwarmCollectViewCreate - Applies a collection method and gathers points
13045627991aSBarry Smith                               in neighbour ranks into the DMSwarm
1305d3a51819SDave May 
1306d083f849SBarry Smith    Collective on dm
1307d3a51819SDave May 
1308d3a51819SDave May    Input parameter:
1309d3a51819SDave May .  dm - the DMSwarm
1310d3a51819SDave May 
1311d3a51819SDave May    Notes:
1312d3a51819SDave May    Users should call DMSwarmCollectViewDestroy() after
1313d3a51819SDave May    they have finished computations associated with the collected points
13148b8a3813SDave May    Different collect methods are supported. See DMSwarmSetCollectType().
1315d3a51819SDave May 
1316d3a51819SDave May    Level: advanced
1317d3a51819SDave May 
1318d3a51819SDave May .seealso: DMSwarmCollectViewDestroy(), DMSwarmSetCollectType()
1319d3a51819SDave May @*/
132074d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmCollectViewCreate(DM dm)
13212712d1f2SDave May {
13222712d1f2SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
13232712d1f2SDave May   PetscInt       ng;
13242712d1f2SDave May 
1325521f74f9SMatthew G. Knepley   PetscFunctionBegin;
132628b400f6SJacob Faibussowitsch   PetscCheck(!swarm->collect_view_active,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"CollectView currently active");
13279566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetLocalSize(dm,&ng));
1328480eef7bSDave May   switch (swarm->collect_type) {
1329f0cdbbbaSDave May 
1330480eef7bSDave May     case DMSWARM_COLLECT_BASIC:
13319566063dSJacob Faibussowitsch       PetscCall(DMSwarmMigrate_GlobalToLocal_Basic(dm,&ng));
1332480eef7bSDave May       break;
1333480eef7bSDave May     case DMSWARM_COLLECT_DMDABOUNDINGBOX:
1334f0cdbbbaSDave May       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_COLLECT_DMDABOUNDINGBOX not implemented");
1335480eef7bSDave May     case DMSWARM_COLLECT_GENERAL:
1336f0cdbbbaSDave May       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_COLLECT_GENERAL not implemented");
1337480eef7bSDave May     default:
1338f0cdbbbaSDave May       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_COLLECT type unknown");
1339480eef7bSDave May   }
1340480eef7bSDave May   swarm->collect_view_active = PETSC_TRUE;
1341480eef7bSDave May   swarm->collect_view_reset_nlocal = ng;
13422712d1f2SDave May   PetscFunctionReturn(0);
13432712d1f2SDave May }
13442712d1f2SDave May 
134574d0cae8SMatthew G. Knepley /*@
1346d3a51819SDave May    DMSwarmCollectViewDestroy - Resets the DMSwarm to the size prior to calling DMSwarmCollectViewCreate()
1347d3a51819SDave May 
1348d083f849SBarry Smith    Collective on dm
1349d3a51819SDave May 
1350d3a51819SDave May    Input parameters:
1351d3a51819SDave May .  dm - the DMSwarm
1352d3a51819SDave May 
1353d3a51819SDave May    Notes:
1354d3a51819SDave May    Users should call DMSwarmCollectViewCreate() before this function is called.
1355d3a51819SDave May 
1356d3a51819SDave May    Level: advanced
1357d3a51819SDave May 
1358d3a51819SDave May .seealso: DMSwarmCollectViewCreate(), DMSwarmSetCollectType()
1359d3a51819SDave May @*/
136074d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmCollectViewDestroy(DM dm)
13612712d1f2SDave May {
13622712d1f2SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
13632712d1f2SDave May 
1364521f74f9SMatthew G. Knepley   PetscFunctionBegin;
136528b400f6SJacob Faibussowitsch   PetscCheck(swarm->collect_view_active,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"CollectView is currently not active");
13669566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(dm,swarm->collect_view_reset_nlocal,-1));
1367480eef7bSDave May   swarm->collect_view_active = PETSC_FALSE;
13682712d1f2SDave May   PetscFunctionReturn(0);
13692712d1f2SDave May }
13703454631fSDave May 
1371f0cdbbbaSDave May PetscErrorCode DMSwarmSetUpPIC(DM dm)
1372f0cdbbbaSDave May {
1373f0cdbbbaSDave May   PetscInt       dim;
1374f0cdbbbaSDave May 
1375521f74f9SMatthew G. Knepley   PetscFunctionBegin;
13769566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetNumSpecies(dm, 1));
13779566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm,&dim));
137808401ef6SPierre Jolivet   PetscCheck(dim >= 1,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Dimension must be 1,2,3 - found %D",dim);
137908401ef6SPierre Jolivet   PetscCheck(dim <= 3,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Dimension must be 1,2,3 - found %D",dim);
13809566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dm,DMSwarmPICField_coor,dim,PETSC_DOUBLE));
13819566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dm,DMSwarmPICField_cellid,1,PETSC_INT));
1382f0cdbbbaSDave May   PetscFunctionReturn(0);
1383f0cdbbbaSDave May }
1384f0cdbbbaSDave May 
138574d0cae8SMatthew G. Knepley /*@
138631403186SMatthew G. Knepley   DMSwarmSetPointCoordinatesRandom - Sets initial coordinates for particles in each cell
138731403186SMatthew G. Knepley 
138831403186SMatthew G. Knepley   Collective on dm
138931403186SMatthew G. Knepley 
139031403186SMatthew G. Knepley   Input parameters:
139131403186SMatthew G. Knepley + dm  - the DMSwarm
139231403186SMatthew G. Knepley - Npc - The number of particles per cell in the cell DM
139331403186SMatthew G. Knepley 
139431403186SMatthew G. Knepley   Notes:
139531403186SMatthew G. Knepley   The user must use DMSwarmSetCellDM() to set the cell DM first. The particles are placed randomly inside each cell. If only
139631403186SMatthew G. Knepley   one particle is in each cell, it is placed at the centroid.
139731403186SMatthew G. Knepley 
139831403186SMatthew G. Knepley   Level: intermediate
139931403186SMatthew G. Knepley 
140031403186SMatthew G. Knepley .seealso: DMSwarmSetCellDM()
140131403186SMatthew G. Knepley @*/
140231403186SMatthew G. Knepley PetscErrorCode DMSwarmSetPointCoordinatesRandom(DM dm, PetscInt Npc)
140331403186SMatthew G. Knepley {
140431403186SMatthew G. Knepley   DM             cdm;
140531403186SMatthew G. Knepley   PetscRandom    rnd;
140631403186SMatthew G. Knepley   DMPolytopeType ct;
140731403186SMatthew G. Knepley   PetscBool      simplex;
140831403186SMatthew G. Knepley   PetscReal     *centroid, *coords, *xi0, *v0, *J, *invJ, detJ;
140931403186SMatthew G. Knepley   PetscInt       dim, d, cStart, cEnd, c, p;
141031403186SMatthew G. Knepley 
141131403186SMatthew G. Knepley   PetscFunctionBeginUser;
14129566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject) dm), &rnd));
14139566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetInterval(rnd, -1.0, 1.0));
14149566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetType(rnd, PETSCRAND48));
141531403186SMatthew G. Knepley 
14169566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetCellDM(dm, &cdm));
14179566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(cdm, &dim));
14189566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(cdm, 0, &cStart, &cEnd));
14199566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCellType(cdm, cStart, &ct));
142031403186SMatthew G. Knepley   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
142131403186SMatthew G. Knepley 
14229566063dSJacob Faibussowitsch   PetscCall(PetscMalloc5(dim, &centroid, dim, &xi0, dim, &v0, dim*dim, &J, dim*dim, &invJ));
142331403186SMatthew G. Knepley   for (d = 0; d < dim; ++d) xi0[d] = -1.0;
14249566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, NULL, NULL, (void **) &coords));
142531403186SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
142631403186SMatthew G. Knepley     if (Npc == 1) {
14279566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFVM(cdm, c, NULL, centroid, NULL));
142831403186SMatthew G. Knepley       for (d = 0; d < dim; ++d) coords[c*dim+d] = centroid[d];
142931403186SMatthew G. Knepley     } else {
14309566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFEM(cdm, c, NULL, v0, J, invJ, &detJ)); /* affine */
143131403186SMatthew G. Knepley       for (p = 0; p < Npc; ++p) {
143231403186SMatthew G. Knepley         const PetscInt n   = c*Npc + p;
143331403186SMatthew G. Knepley         PetscReal      sum = 0.0, refcoords[3];
143431403186SMatthew G. Knepley 
143531403186SMatthew G. Knepley         for (d = 0; d < dim; ++d) {
14369566063dSJacob Faibussowitsch           PetscCall(PetscRandomGetValueReal(rnd, &refcoords[d]));
143731403186SMatthew G. Knepley           sum += refcoords[d];
143831403186SMatthew G. Knepley         }
143931403186SMatthew G. Knepley         if (simplex && sum > 0.0) for (d = 0; d < dim; ++d) refcoords[d] -= PetscSqrtReal(dim)*sum;
144031403186SMatthew G. Knepley         CoordinatesRefToReal(dim, dim, xi0, v0, J, refcoords, &coords[n*dim]);
144131403186SMatthew G. Knepley       }
144231403186SMatthew G. Knepley     }
144331403186SMatthew G. Knepley   }
14449566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, NULL, NULL, (void **) &coords));
14459566063dSJacob Faibussowitsch   PetscCall(PetscFree5(centroid, xi0, v0, J, invJ));
14469566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rnd));
144731403186SMatthew G. Knepley   PetscFunctionReturn(0);
144831403186SMatthew G. Knepley }
144931403186SMatthew G. Knepley 
145031403186SMatthew G. Knepley /*@
1451d3a51819SDave May    DMSwarmSetType - Set particular flavor of DMSwarm
1452d3a51819SDave May 
1453d083f849SBarry Smith    Collective on dm
1454d3a51819SDave May 
1455d3a51819SDave May    Input parameters:
145662741f57SDave May +  dm - the DMSwarm
145762741f57SDave May -  stype - the DMSwarm type (e.g. DMSWARM_PIC)
1458d3a51819SDave May 
1459d3a51819SDave May    Level: advanced
1460d3a51819SDave May 
14615627991aSBarry Smith .seealso: DMSwarmSetMigrateType(), DMSwarmSetCollectType(), DMSwarmType, DMSWARM_PIC, DMSWARM_BASIC
1462d3a51819SDave May @*/
146374d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmSetType(DM dm,DMSwarmType stype)
1464f0cdbbbaSDave May {
1465f0cdbbbaSDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
1466f0cdbbbaSDave May 
1467521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1468f0cdbbbaSDave May   swarm->swarm_type = stype;
1469f0cdbbbaSDave May   if (swarm->swarm_type == DMSWARM_PIC) {
14709566063dSJacob Faibussowitsch     PetscCall(DMSwarmSetUpPIC(dm));
1471f0cdbbbaSDave May   }
1472f0cdbbbaSDave May   PetscFunctionReturn(0);
1473f0cdbbbaSDave May }
1474f0cdbbbaSDave May 
14753454631fSDave May PetscErrorCode DMSetup_Swarm(DM dm)
14763454631fSDave May {
14773454631fSDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
14783454631fSDave May   PetscMPIInt    rank;
14793454631fSDave May   PetscInt       p,npoints,*rankval;
14803454631fSDave May 
1481521f74f9SMatthew G. Knepley   PetscFunctionBegin;
14823454631fSDave May   if (swarm->issetup) PetscFunctionReturn(0);
14833454631fSDave May   swarm->issetup = PETSC_TRUE;
14843454631fSDave May 
1485f0cdbbbaSDave May   if (swarm->swarm_type == DMSWARM_PIC) {
1486f0cdbbbaSDave May     /* check dmcell exists */
148728b400f6SJacob Faibussowitsch     PetscCheck(swarm->dmcell,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"DMSWARM_PIC requires you call DMSwarmSetCellDM");
1488f0cdbbbaSDave May 
1489f0cdbbbaSDave May     if (swarm->dmcell->ops->locatepointssubdomain) {
1490f0cdbbbaSDave May       /* check methods exists for exact ownership identificiation */
14919566063dSJacob Faibussowitsch       PetscCall(PetscInfo(dm, "DMSWARM_PIC: Using method CellDM->ops->LocatePointsSubdomain\n"));
1492f0cdbbbaSDave May       swarm->migrate_type = DMSWARM_MIGRATE_DMCELLEXACT;
1493f0cdbbbaSDave May     } else {
1494f0cdbbbaSDave May       /* check methods exist for point location AND rank neighbor identification */
1495f0cdbbbaSDave May       if (swarm->dmcell->ops->locatepoints) {
14969566063dSJacob Faibussowitsch         PetscCall(PetscInfo(dm, "DMSWARM_PIC: Using method CellDM->LocatePoints\n"));
1497f0cdbbbaSDave May       } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"DMSWARM_PIC requires the method CellDM->ops->locatepoints be defined");
1498f0cdbbbaSDave May 
1499f0cdbbbaSDave May       if (swarm->dmcell->ops->getneighbors) {
15009566063dSJacob Faibussowitsch         PetscCall(PetscInfo(dm, "DMSWARM_PIC: Using method CellDM->GetNeigbors\n"));
1501f0cdbbbaSDave May       } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"DMSWARM_PIC requires the method CellDM->ops->getneighbors be defined");
1502f0cdbbbaSDave May 
1503f0cdbbbaSDave May       swarm->migrate_type = DMSWARM_MIGRATE_DMCELLNSCATTER;
1504f0cdbbbaSDave May     }
1505f0cdbbbaSDave May   }
1506f0cdbbbaSDave May 
15079566063dSJacob Faibussowitsch   PetscCall(DMSwarmFinalizeFieldRegister(dm));
1508f0cdbbbaSDave May 
15093454631fSDave May   /* check some fields were registered */
151008401ef6SPierre Jolivet   PetscCheck(swarm->db->nfields > 2,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"At least one field user must be registered via DMSwarmRegisterXXX()");
15113454631fSDave May 
15123454631fSDave May   /* check local sizes were set */
151308401ef6SPierre Jolivet   PetscCheck(swarm->db->L != -1,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Local sizes must be set via DMSwarmSetLocalSizes()");
15143454631fSDave May 
15153454631fSDave May   /* initialize values in pid and rank placeholders */
15163454631fSDave May   /* TODO: [pid - use MPI_Scan] */
15179566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank));
15189566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL));
15199566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval));
15203454631fSDave May   for (p=0; p<npoints; p++) {
15213454631fSDave May     rankval[p] = (PetscInt)rank;
15223454631fSDave May   }
15239566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval));
15243454631fSDave May   PetscFunctionReturn(0);
15253454631fSDave May }
15263454631fSDave May 
1527dc5f5ce9SDave May extern PetscErrorCode DMSwarmSortDestroy(DMSwarmSort *_ctx);
1528dc5f5ce9SDave May 
152957795646SDave May PetscErrorCode DMDestroy_Swarm(DM dm)
153057795646SDave May {
153157795646SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
153257795646SDave May 
153357795646SDave May   PetscFunctionBegin;
1534d0c080abSJoseph Pusztay   if (--swarm->refct > 0) PetscFunctionReturn(0);
15359566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketDestroy(&swarm->db));
1536dc5f5ce9SDave May   if (swarm->sort_context) {
15379566063dSJacob Faibussowitsch     PetscCall(DMSwarmSortDestroy(&swarm->sort_context));
1538dc5f5ce9SDave May   }
15399566063dSJacob Faibussowitsch   PetscCall(PetscFree(swarm));
154057795646SDave May   PetscFunctionReturn(0);
154157795646SDave May }
154257795646SDave May 
1543a9ee3421SMatthew G. Knepley PetscErrorCode DMSwarmView_Draw(DM dm, PetscViewer viewer)
1544a9ee3421SMatthew G. Knepley {
1545a9ee3421SMatthew G. Knepley   DM             cdm;
1546a9ee3421SMatthew G. Knepley   PetscDraw      draw;
154731403186SMatthew G. Knepley   PetscReal     *coords, oldPause, radius = 0.01;
1548a9ee3421SMatthew G. Knepley   PetscInt       Np, p, bs;
1549a9ee3421SMatthew G. Knepley 
1550a9ee3421SMatthew G. Knepley   PetscFunctionBegin;
15519566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, ((PetscObject) dm)->prefix, "-dm_view_swarm_radius", &radius, NULL));
15529566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
15539566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetCellDM(dm, &cdm));
15549566063dSJacob Faibussowitsch   PetscCall(PetscDrawGetPause(draw, &oldPause));
15559566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetPause(draw, 0.0));
15569566063dSJacob Faibussowitsch   PetscCall(DMView(cdm, viewer));
15579566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetPause(draw, oldPause));
1558a9ee3421SMatthew G. Knepley 
15599566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetLocalSize(dm, &Np));
15609566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, &bs, NULL, (void **) &coords));
1561a9ee3421SMatthew G. Knepley   for (p = 0; p < Np; ++p) {
1562a9ee3421SMatthew G. Knepley     const PetscInt i = p*bs;
1563a9ee3421SMatthew G. Knepley 
15649566063dSJacob Faibussowitsch     PetscCall(PetscDrawEllipse(draw, coords[i], coords[i+1], radius, radius, PETSC_DRAW_BLUE));
1565a9ee3421SMatthew G. Knepley   }
15669566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, &bs, NULL, (void **) &coords));
15679566063dSJacob Faibussowitsch   PetscCall(PetscDrawFlush(draw));
15689566063dSJacob Faibussowitsch   PetscCall(PetscDrawPause(draw));
15699566063dSJacob Faibussowitsch   PetscCall(PetscDrawSave(draw));
1570a9ee3421SMatthew G. Knepley   PetscFunctionReturn(0);
1571a9ee3421SMatthew G. Knepley }
1572a9ee3421SMatthew G. Knepley 
1573*6a5217c0SMatthew G. Knepley static PetscErrorCode DMView_Swarm_Ascii(DM dm, PetscViewer viewer)
1574*6a5217c0SMatthew G. Knepley {
1575*6a5217c0SMatthew G. Knepley   PetscViewerFormat format;
1576*6a5217c0SMatthew G. Knepley   PetscInt         *sizes;
1577*6a5217c0SMatthew G. Knepley   PetscInt          dim, Np, maxSize = 17;
1578*6a5217c0SMatthew G. Knepley   MPI_Comm          comm;
1579*6a5217c0SMatthew G. Knepley   PetscMPIInt       rank, size, p;
1580*6a5217c0SMatthew G. Knepley   const char       *name;
1581*6a5217c0SMatthew G. Knepley 
1582*6a5217c0SMatthew G. Knepley   PetscFunctionBegin;
1583*6a5217c0SMatthew G. Knepley   PetscCall(PetscViewerGetFormat(viewer, &format));
1584*6a5217c0SMatthew G. Knepley   PetscCall(DMGetDimension(dm, &dim));
1585*6a5217c0SMatthew G. Knepley   PetscCall(DMSwarmGetLocalSize(dm, &Np));
1586*6a5217c0SMatthew G. Knepley   PetscCall(PetscObjectGetComm((PetscObject) dm, &comm));
1587*6a5217c0SMatthew G. Knepley   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1588*6a5217c0SMatthew G. Knepley   PetscCallMPI(MPI_Comm_size(comm, &size));
1589*6a5217c0SMatthew G. Knepley   PetscCall(PetscObjectGetName((PetscObject) dm, &name));
1590*6a5217c0SMatthew G. Knepley   if (name) PetscCall(PetscViewerASCIIPrintf(viewer, "%s in %D dimension%s:\n", name, dim, dim == 1 ? "" : "s"));
1591*6a5217c0SMatthew G. Knepley   else      PetscCall(PetscViewerASCIIPrintf(viewer, "Swarm in %D dimension%s:\n", dim, dim == 1 ? "" : "s"));
1592*6a5217c0SMatthew G. Knepley   if (size < maxSize) PetscCall(PetscCalloc1(size, &sizes));
1593*6a5217c0SMatthew G. Knepley   else                PetscCall(PetscCalloc1(3,    &sizes));
1594*6a5217c0SMatthew G. Knepley   if (size < maxSize) {
1595*6a5217c0SMatthew G. Knepley     PetscCallMPI(MPI_Gather(&Np, 1, MPIU_INT, sizes, 1, MPIU_INT, 0, comm));
1596*6a5217c0SMatthew G. Knepley     PetscCall(PetscViewerASCIIPrintf(viewer, "  Number of particles per rank:"));
1597*6a5217c0SMatthew G. Knepley     for (p = 0; p < size; ++p) {
1598*6a5217c0SMatthew G. Knepley       if (rank == 0) PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT, sizes[p]));
1599*6a5217c0SMatthew G. Knepley     }
1600*6a5217c0SMatthew G. Knepley   } else {
1601*6a5217c0SMatthew G. Knepley     PetscInt locMinMax[2] = {Np, Np};
1602*6a5217c0SMatthew G. Knepley 
1603*6a5217c0SMatthew G. Knepley     PetscCall(PetscGlobalMinMaxInt(comm, locMinMax, sizes));
1604*6a5217c0SMatthew G. Knepley     PetscCall(PetscViewerASCIIPrintf(viewer, "  Min/Max of particles per rank: %" PetscInt_FMT "/%" PetscInt_FMT, sizes[0], sizes[1]));
1605*6a5217c0SMatthew G. Knepley   }
1606*6a5217c0SMatthew G. Knepley   PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1607*6a5217c0SMatthew G. Knepley   PetscCall(PetscFree(sizes));
1608*6a5217c0SMatthew G. Knepley   if (format == PETSC_VIEWER_ASCII_INFO) {
1609*6a5217c0SMatthew G. Knepley     PetscInt *cell;
1610*6a5217c0SMatthew G. Knepley 
1611*6a5217c0SMatthew G. Knepley     PetscCall(PetscViewerASCIIPrintf(viewer, "  Cells containing each particle:\n"));
1612*6a5217c0SMatthew G. Knepley     PetscCall(PetscViewerASCIIPushSynchronized(viewer));
1613*6a5217c0SMatthew G. Knepley     PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void**) &cell));
1614*6a5217c0SMatthew G. Knepley     for (p = 0; p < Np; ++p) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "  p%" PetscInt_FMT ": %" PetscInt_FMT "\n", p, cell[p]));
1615*6a5217c0SMatthew G. Knepley     PetscCall(PetscViewerFlush(viewer));
1616*6a5217c0SMatthew G. Knepley     PetscCall(PetscViewerASCIIPopSynchronized(viewer));
1617*6a5217c0SMatthew G. Knepley     PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void**) &cell));
1618*6a5217c0SMatthew G. Knepley   }
1619*6a5217c0SMatthew G. Knepley   PetscFunctionReturn(0);
1620*6a5217c0SMatthew G. Knepley }
1621*6a5217c0SMatthew G. Knepley 
16225f50eb2eSDave May PetscErrorCode DMView_Swarm(DM dm, PetscViewer viewer)
16235f50eb2eSDave May {
16245f50eb2eSDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
16255627991aSBarry Smith   PetscBool      iascii,ibinary,isvtk,isdraw;
16265627991aSBarry Smith #if defined(PETSC_HAVE_HDF5)
16275627991aSBarry Smith   PetscBool      ishdf5;
16285627991aSBarry Smith #endif
16295f50eb2eSDave May 
16305f50eb2eSDave May   PetscFunctionBegin;
16315f50eb2eSDave May   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
16325f50eb2eSDave May   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
16339566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
16349566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY,&ibinary));
16359566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERVTK,   &isvtk));
16365627991aSBarry Smith #if defined(PETSC_HAVE_HDF5)
16379566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5,  &ishdf5));
16385627991aSBarry Smith #endif
16399566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW,  &isdraw));
16405f50eb2eSDave May   if (iascii) {
1641*6a5217c0SMatthew G. Knepley     PetscViewerFormat format;
1642*6a5217c0SMatthew G. Knepley 
1643*6a5217c0SMatthew G. Knepley     PetscCall(PetscViewerGetFormat(viewer, &format));
1644*6a5217c0SMatthew G. Knepley     switch (format) {
1645*6a5217c0SMatthew G. Knepley       case PETSC_VIEWER_ASCII_INFO_DETAIL:
1646*6a5217c0SMatthew G. Knepley         PetscCall(DMSwarmDataBucketView(PetscObjectComm((PetscObject) dm), swarm->db, NULL, DATABUCKET_VIEW_STDOUT));break;
1647*6a5217c0SMatthew G. Knepley       default: PetscCall(DMView_Swarm_Ascii(dm, viewer));
1648*6a5217c0SMatthew G. Knepley     }
164928b400f6SJacob Faibussowitsch   } else PetscCheck(!ibinary,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"NO Binary support");
16505f50eb2eSDave May #if defined(PETSC_HAVE_HDF5)
165174d0cae8SMatthew G. Knepley   else if (ishdf5) {
16529566063dSJacob Faibussowitsch     PetscCall(DMSwarmView_HDF5(dm, viewer));
165374d0cae8SMatthew G. Knepley   }
16545f50eb2eSDave May #endif
1655298827fbSBarry Smith   else if (isdraw) {
16569566063dSJacob Faibussowitsch     PetscCall(DMSwarmView_Draw(dm, viewer));
16575f50eb2eSDave May   }
16585f50eb2eSDave May   PetscFunctionReturn(0);
16595f50eb2eSDave May }
16605f50eb2eSDave May 
1661d0c080abSJoseph Pusztay /*@C
1662d0c080abSJoseph Pusztay    DMSwarmGetCellSwarm - Extracts a single cell from the DMSwarm object, returns it as a single cell DMSwarm.
1663d0c080abSJoseph Pusztay    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.
1664d0c080abSJoseph Pusztay 
1665d0c080abSJoseph Pusztay    Important: Changes to this cell of the swarm will be lost if they are made prior to restoring this cell.
1666d0c080abSJoseph Pusztay 
1667d0c080abSJoseph Pusztay    Noncollective
1668d0c080abSJoseph Pusztay 
1669d0c080abSJoseph Pusztay    Input parameters:
1670d0c080abSJoseph Pusztay +  sw - the DMSwarm
16715627991aSBarry Smith .  cellID - the integer id of the cell to be extracted and filtered
16725627991aSBarry Smith -  cellswarm - The DMSwarm to receive the cell
1673d0c080abSJoseph Pusztay 
1674d0c080abSJoseph Pusztay    Level: beginner
1675d0c080abSJoseph Pusztay 
16765627991aSBarry Smith    Notes:
16775627991aSBarry Smith       This presently only supports DMSWARM_PIC type
16785627991aSBarry Smith 
16795627991aSBarry Smith       Should be restored with DMSwarmRestoreCellSwarm()
1680d0c080abSJoseph Pusztay 
1681d0c080abSJoseph Pusztay .seealso: DMSwarmRestoreCellSwarm()
1682d0c080abSJoseph Pusztay @*/
1683d0c080abSJoseph Pusztay PETSC_EXTERN PetscErrorCode DMSwarmGetCellSwarm(DM sw, PetscInt cellID, DM cellswarm)
1684d0c080abSJoseph Pusztay {
1685d0c080abSJoseph Pusztay   DM_Swarm      *original = (DM_Swarm*) sw->data;
1686d0c080abSJoseph Pusztay   DMLabel        label;
1687d0c080abSJoseph Pusztay   DM             dmc, subdmc;
1688d0c080abSJoseph Pusztay   PetscInt      *pids, particles, dim;
1689d0c080abSJoseph Pusztay 
1690d0c080abSJoseph Pusztay   PetscFunctionBegin;
1691d0c080abSJoseph Pusztay   /* Configure new swarm */
16929566063dSJacob Faibussowitsch   PetscCall(DMSetType(cellswarm, DMSWARM));
16939566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(sw, &dim));
16949566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(cellswarm, dim));
16959566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetType(cellswarm, DMSWARM_PIC));
1696d0c080abSJoseph Pusztay   /* Destroy the unused, unconfigured data bucket to prevent stragglers in memory */
16979566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketDestroy(&((DM_Swarm*)cellswarm->data)->db));
16989566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortGetAccess(sw));
16999566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortGetNumberOfPointsPerCell(sw, cellID, &particles));
17009566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortGetPointsPerCell(sw, cellID, &particles, &pids));
17019566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketCreateFromSubset(original->db, particles, pids, &((DM_Swarm*)cellswarm->data)->db));
17029566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortRestoreAccess(sw));
17039566063dSJacob Faibussowitsch   PetscCall(PetscFree(pids));
17049566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetCellDM(sw, &dmc));
17059566063dSJacob Faibussowitsch   PetscCall(DMLabelCreate(PetscObjectComm((PetscObject)sw), "singlecell", &label));
17069566063dSJacob Faibussowitsch   PetscCall(DMAddLabel(dmc, label));
17079566063dSJacob Faibussowitsch   PetscCall(DMLabelSetValue(label, cellID, 1));
17089566063dSJacob Faibussowitsch   PetscCall(DMPlexFilter(dmc, label, 1, &subdmc));
17099566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetCellDM(cellswarm, subdmc));
17109566063dSJacob Faibussowitsch   PetscCall(DMLabelDestroy(&label));
1711d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
1712d0c080abSJoseph Pusztay }
1713d0c080abSJoseph Pusztay 
1714d0c080abSJoseph Pusztay /*@C
17155627991aSBarry Smith    DMSwarmRestoreCellSwarm - Restores a DMSwarm object obtained with DMSwarmGetCellSwarm(). All fields are copied back into the parent swarm.
1716d0c080abSJoseph Pusztay 
1717d0c080abSJoseph Pusztay    Noncollective
1718d0c080abSJoseph Pusztay 
1719d0c080abSJoseph Pusztay    Input parameters:
1720d0c080abSJoseph Pusztay +  sw - the parent DMSwarm
1721d0c080abSJoseph Pusztay .  cellID - the integer id of the cell to be copied back into the parent swarm
1722d0c080abSJoseph Pusztay -  cellswarm - the cell swarm object
1723d0c080abSJoseph Pusztay 
1724d0c080abSJoseph Pusztay    Level: beginner
1725d0c080abSJoseph Pusztay 
17265627991aSBarry Smith    Note:
17275627991aSBarry Smith     This only supports DMSWARM_PIC types of DMSwarms
1728d0c080abSJoseph Pusztay 
1729d0c080abSJoseph Pusztay .seealso: DMSwarmGetCellSwarm()
1730d0c080abSJoseph Pusztay @*/
1731d0c080abSJoseph Pusztay PETSC_EXTERN PetscErrorCode DMSwarmRestoreCellSwarm(DM sw, PetscInt cellID, DM cellswarm)
1732d0c080abSJoseph Pusztay {
1733d0c080abSJoseph Pusztay   DM             dmc;
1734d0c080abSJoseph Pusztay   PetscInt       *pids, particles, p;
1735d0c080abSJoseph Pusztay 
1736d0c080abSJoseph Pusztay   PetscFunctionBegin;
17379566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortGetAccess(sw));
17389566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortGetPointsPerCell(sw, cellID, &particles, &pids));
17399566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortRestoreAccess(sw));
1740d0c080abSJoseph Pusztay   /* Pointwise copy of each particle based on pid. The parent swarm may not be altered during this process. */
1741d0c080abSJoseph Pusztay   for (p=0; p<particles; ++p) {
17429566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketCopyPoint(((DM_Swarm*)cellswarm->data)->db,pids[p],((DM_Swarm*)sw->data)->db,pids[p]));
1743d0c080abSJoseph Pusztay   }
1744d0c080abSJoseph Pusztay   /* Free memory, destroy cell dm */
17459566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetCellDM(cellswarm, &dmc));
17469566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmc));
17479566063dSJacob Faibussowitsch   PetscCall(PetscFree(pids));
1748d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
1749d0c080abSJoseph Pusztay }
1750d0c080abSJoseph Pusztay 
1751d0c080abSJoseph Pusztay PETSC_INTERN PetscErrorCode DMClone_Swarm(DM, DM *);
1752d0c080abSJoseph Pusztay 
1753d0c080abSJoseph Pusztay static PetscErrorCode DMInitialize_Swarm(DM sw)
1754d0c080abSJoseph Pusztay {
1755d0c080abSJoseph Pusztay   PetscFunctionBegin;
1756d0c080abSJoseph Pusztay   sw->dim  = 0;
1757d0c080abSJoseph Pusztay   sw->ops->view                            = DMView_Swarm;
1758d0c080abSJoseph Pusztay   sw->ops->load                            = NULL;
1759d0c080abSJoseph Pusztay   sw->ops->setfromoptions                  = NULL;
1760d0c080abSJoseph Pusztay   sw->ops->clone                           = DMClone_Swarm;
1761d0c080abSJoseph Pusztay   sw->ops->setup                           = DMSetup_Swarm;
1762d0c080abSJoseph Pusztay   sw->ops->createlocalsection              = NULL;
1763d0c080abSJoseph Pusztay   sw->ops->createdefaultconstraints        = NULL;
1764d0c080abSJoseph Pusztay   sw->ops->createglobalvector              = DMCreateGlobalVector_Swarm;
1765d0c080abSJoseph Pusztay   sw->ops->createlocalvector               = DMCreateLocalVector_Swarm;
1766d0c080abSJoseph Pusztay   sw->ops->getlocaltoglobalmapping         = NULL;
1767d0c080abSJoseph Pusztay   sw->ops->createfieldis                   = NULL;
1768d0c080abSJoseph Pusztay   sw->ops->createcoordinatedm              = NULL;
1769d0c080abSJoseph Pusztay   sw->ops->getcoloring                     = NULL;
1770d0c080abSJoseph Pusztay   sw->ops->creatematrix                    = DMCreateMatrix_Swarm;
1771d0c080abSJoseph Pusztay   sw->ops->createinterpolation             = NULL;
1772d0c080abSJoseph Pusztay   sw->ops->createinjection                 = NULL;
1773d0c080abSJoseph Pusztay   sw->ops->createmassmatrix                = DMCreateMassMatrix_Swarm;
1774d0c080abSJoseph Pusztay   sw->ops->refine                          = NULL;
1775d0c080abSJoseph Pusztay   sw->ops->coarsen                         = NULL;
1776d0c080abSJoseph Pusztay   sw->ops->refinehierarchy                 = NULL;
1777d0c080abSJoseph Pusztay   sw->ops->coarsenhierarchy                = NULL;
1778d0c080abSJoseph Pusztay   sw->ops->globaltolocalbegin              = NULL;
1779d0c080abSJoseph Pusztay   sw->ops->globaltolocalend                = NULL;
1780d0c080abSJoseph Pusztay   sw->ops->localtoglobalbegin              = NULL;
1781d0c080abSJoseph Pusztay   sw->ops->localtoglobalend                = NULL;
1782d0c080abSJoseph Pusztay   sw->ops->destroy                         = DMDestroy_Swarm;
1783d0c080abSJoseph Pusztay   sw->ops->createsubdm                     = NULL;
1784d0c080abSJoseph Pusztay   sw->ops->getdimpoints                    = NULL;
1785d0c080abSJoseph Pusztay   sw->ops->locatepoints                    = NULL;
1786d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
1787d0c080abSJoseph Pusztay }
1788d0c080abSJoseph Pusztay 
1789d0c080abSJoseph Pusztay PETSC_INTERN PetscErrorCode DMClone_Swarm(DM dm, DM *newdm)
1790d0c080abSJoseph Pusztay {
1791d0c080abSJoseph Pusztay   DM_Swarm       *swarm = (DM_Swarm *) dm->data;
1792d0c080abSJoseph Pusztay 
1793d0c080abSJoseph Pusztay   PetscFunctionBegin;
1794d0c080abSJoseph Pusztay   swarm->refct++;
1795d0c080abSJoseph Pusztay   (*newdm)->data = swarm;
17969566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject) *newdm, DMSWARM));
17979566063dSJacob Faibussowitsch   PetscCall(DMInitialize_Swarm(*newdm));
17982e56dffeSJoe Pusztay   (*newdm)->dim = dm->dim;
1799d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
1800d0c080abSJoseph Pusztay }
1801d0c080abSJoseph Pusztay 
1802d3a51819SDave May /*MC
1803d3a51819SDave May 
1804d3a51819SDave May  DMSWARM = "swarm" - A DM object used to represent arrays of data (fields) of arbitrary data type.
180562741f57SDave May  This implementation was designed for particle methods in which the underlying
1806d3a51819SDave May  data required to be represented is both (i) dynamic in length, (ii) and of arbitrary data type.
1807d3a51819SDave May 
180862741f57SDave May  User data can be represented by DMSwarm through a registering "fields".
180962741f57SDave May  To register a field, the user must provide;
181062741f57SDave May  (a) a unique name;
181162741f57SDave May  (b) the data type (or size in bytes);
181262741f57SDave May  (c) the block size of the data.
1813d3a51819SDave May 
1814d3a51819SDave May  For example, suppose the application requires a unique id, energy, momentum and density to be stored
1815c215a47eSMatthew Knepley  on a set of particles. Then the following code could be used
1816d3a51819SDave May 
181762741f57SDave May $    DMSwarmInitializeFieldRegister(dm)
181862741f57SDave May $    DMSwarmRegisterPetscDatatypeField(dm,"uid",1,PETSC_LONG);
181962741f57SDave May $    DMSwarmRegisterPetscDatatypeField(dm,"energy",1,PETSC_REAL);
182062741f57SDave May $    DMSwarmRegisterPetscDatatypeField(dm,"momentum",3,PETSC_REAL);
182162741f57SDave May $    DMSwarmRegisterPetscDatatypeField(dm,"density",1,PETSC_FLOAT);
182262741f57SDave May $    DMSwarmFinalizeFieldRegister(dm)
1823d3a51819SDave May 
1824d3a51819SDave May  The fields represented by DMSwarm are dynamic and can be re-sized at any time.
182562741f57SDave May  The only restriction imposed by DMSwarm is that all fields contain the same number of points.
1826d3a51819SDave May 
1827d3a51819SDave May  To support particle methods, "migration" techniques are provided. These methods migrate data
18285627991aSBarry Smith  between ranks.
1829d3a51819SDave May 
1830d3a51819SDave May  DMSwarm supports the methods DMCreateGlobalVector() and DMCreateLocalVector().
1831d3a51819SDave May  As a DMSwarm may internally define and store values of different data types,
183262741f57SDave May  before calling DMCreateGlobalVector() or DMCreateLocalVector(), the user must inform DMSwarm which
1833d3a51819SDave May  fields should be used to define a Vec object via
1834d3a51819SDave May    DMSwarmVectorDefineField()
1835c215a47eSMatthew Knepley  The specified field can be changed at any time - thereby permitting vectors
1836c215a47eSMatthew Knepley  compatible with different fields to be created.
1837d3a51819SDave May 
183862741f57SDave May  A dual representation of fields in the DMSwarm and a Vec object is permitted via
1839d3a51819SDave May    DMSwarmCreateGlobalVectorFromField()
1840d3a51819SDave May  Here the data defining the field in the DMSwarm is shared with a Vec.
1841d3a51819SDave May  This is inherently unsafe if you alter the size of the field at any time between
1842d3a51819SDave May  calls to DMSwarmCreateGlobalVectorFromField() and DMSwarmDestroyGlobalVectorFromField().
1843cc651181SDave May  If the local size of the DMSwarm does not match the local size of the global vector
1844cc651181SDave May  when DMSwarmDestroyGlobalVectorFromField() is called, an error is thrown.
1845d3a51819SDave May 
184662741f57SDave May  Additional high-level support is provided for Particle-In-Cell methods.
184762741f57SDave May  Please refer to the man page for DMSwarmSetType().
184862741f57SDave May 
1849d3a51819SDave May  Level: beginner
1850d3a51819SDave May 
1851d3a51819SDave May .seealso: DMType, DMCreate(), DMSetType()
1852d3a51819SDave May M*/
185357795646SDave May PETSC_EXTERN PetscErrorCode DMCreate_Swarm(DM dm)
185457795646SDave May {
185557795646SDave May   DM_Swarm      *swarm;
185657795646SDave May 
185757795646SDave May   PetscFunctionBegin;
185857795646SDave May   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
18599566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(dm,&swarm));
1860f0cdbbbaSDave May   dm->data = swarm;
18619566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketCreate(&swarm->db));
18629566063dSJacob Faibussowitsch   PetscCall(DMSwarmInitializeFieldRegister(dm));
1863d0c080abSJoseph Pusztay   swarm->refct = 1;
1864b5bcf523SDave May   swarm->vec_field_set                  = PETSC_FALSE;
18653454631fSDave May   swarm->issetup                        = PETSC_FALSE;
1866480eef7bSDave May   swarm->swarm_type                     = DMSWARM_BASIC;
1867480eef7bSDave May   swarm->migrate_type                   = DMSWARM_MIGRATE_BASIC;
1868480eef7bSDave May   swarm->collect_type                   = DMSWARM_COLLECT_BASIC;
186940c453e9SDave May   swarm->migrate_error_on_missing_point = PETSC_FALSE;
1870f0cdbbbaSDave May   swarm->dmcell                         = NULL;
1871f0cdbbbaSDave May   swarm->collect_view_active            = PETSC_FALSE;
1872f0cdbbbaSDave May   swarm->collect_view_reset_nlocal      = -1;
18739566063dSJacob Faibussowitsch   PetscCall(DMInitialize_Swarm(dm));
187457795646SDave May   PetscFunctionReturn(0);
187557795646SDave May }
1876