xref: /petsc/src/dm/impls/swarm/swarm.c (revision ea78f98c112368f404cd6d4fff6d4dfe73e5a1e7)
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>
8279f676cSBarry Smith #include "../src/dm/impls/swarm/data_bucket.h"
957795646SDave May 
10f2b2bee7SDave May PetscLogEvent DMSWARM_Migrate, DMSWARM_SetSizes, DMSWARM_AddPoints, DMSWARM_RemovePoints, DMSWARM_Sort;
11ed923d71SDave May PetscLogEvent DMSWARM_DataExchangerTopologySetup, DMSWARM_DataExchangerBegin, DMSWARM_DataExchangerEnd;
12ed923d71SDave May PetscLogEvent DMSWARM_DataExchangerSendCount, DMSWARM_DataExchangerPack;
13ed923d71SDave May 
14*ea78f98cSLisandro Dalcin const char* DMSwarmTypeNames[] = { "basic", "pic", NULL };
15*ea78f98cSLisandro Dalcin const char* DMSwarmMigrateTypeNames[] = { "basic", "dmcellnscatter", "dmcellexact", "user", NULL };
16*ea78f98cSLisandro Dalcin const char* DMSwarmCollectTypeNames[] = { "basic", "boundingbox", "general", "user", NULL  };
17*ea78f98cSLisandro Dalcin const char* DMSwarmPICLayoutTypeNames[] = { "regular", "gauss", "subdivision", NULL  };
18f0cdbbbaSDave May 
19f0cdbbbaSDave May const char DMSwarmField_pid[] = "DMSwarm_pid";
20f0cdbbbaSDave May const char DMSwarmField_rank[] = "DMSwarm_rank";
21f0cdbbbaSDave May const char DMSwarmPICField_coor[] = "DMSwarmPIC_coor";
22e2d107dbSDave May const char DMSwarmPICField_cellid[] = "DMSwarm_cellid";
23f0cdbbbaSDave May 
2474d0cae8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode VecView_Seq(Vec, PetscViewer);
2574d0cae8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode VecView_MPI(Vec, PetscViewer);
2674d0cae8SMatthew G. Knepley 
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   PetscErrorCode ierr;
3774d0cae8SMatthew G. Knepley 
3874d0cae8SMatthew G. Knepley   PetscFunctionBegin;
3974d0cae8SMatthew G. Knepley   ierr = VecGetDM(v, &dm);CHKERRQ(ierr);
4074d0cae8SMatthew G. Knepley   ierr = VecGetBlockSize(v, &bs);CHKERRQ(ierr);
4174d0cae8SMatthew G. Knepley   ierr = PetscViewerHDF5PushGroup(viewer, "/particle_fields");CHKERRQ(ierr);
4274d0cae8SMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) v, VECSEQ, &isseq);CHKERRQ(ierr);
4374d0cae8SMatthew G. Knepley   ierr = DMGetOutputSequenceNumber(dm, &seqnum, &seqval);CHKERRQ(ierr);
4474d0cae8SMatthew G. Knepley   ierr = PetscViewerHDF5SetTimestep(viewer, seqnum);CHKERRQ(ierr);
4574d0cae8SMatthew G. Knepley   /* ierr = DMSequenceView_HDF5(dm, "time", seqnum, (PetscScalar) seqval, viewer);CHKERRQ(ierr); */
4674d0cae8SMatthew G. Knepley   if (isseq) {ierr = VecView_Seq(v, viewer);CHKERRQ(ierr);}
4774d0cae8SMatthew G. Knepley   else       {ierr = VecView_MPI(v, viewer);CHKERRQ(ierr);}
4874d0cae8SMatthew G. Knepley   ierr = PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject) v, "Nc", PETSC_INT, (void *) &bs);CHKERRQ(ierr);
4974d0cae8SMatthew G. Knepley   ierr = PetscViewerHDF5PopGroup(viewer);CHKERRQ(ierr);
5074d0cae8SMatthew G. Knepley   PetscFunctionReturn(0);
5174d0cae8SMatthew G. Knepley }
5274d0cae8SMatthew G. Knepley 
5374d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmView_HDF5(DM dm, PetscViewer viewer)
5474d0cae8SMatthew G. Knepley {
5574d0cae8SMatthew G. Knepley   Vec            coordinates;
5674d0cae8SMatthew G. Knepley   PetscInt       Np;
5774d0cae8SMatthew G. Knepley   PetscBool      isseq;
5874d0cae8SMatthew G. Knepley   PetscErrorCode ierr;
5974d0cae8SMatthew G. Knepley 
6074d0cae8SMatthew G. Knepley   PetscFunctionBegin;
6174d0cae8SMatthew G. Knepley   ierr = DMSwarmGetSize(dm, &Np);CHKERRQ(ierr);
6274d0cae8SMatthew G. Knepley   ierr = DMSwarmCreateGlobalVectorFromField(dm, DMSwarmPICField_coor, &coordinates);CHKERRQ(ierr);
6374d0cae8SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
6474d0cae8SMatthew G. Knepley   ierr = PetscViewerHDF5PushGroup(viewer, "/particles");CHKERRQ(ierr);
6574d0cae8SMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) coordinates, VECSEQ, &isseq);CHKERRQ(ierr);
6674d0cae8SMatthew G. Knepley   if (isseq) {ierr = VecView_Seq(coordinates, viewer);CHKERRQ(ierr);}
6774d0cae8SMatthew G. Knepley   else       {ierr = VecView_MPI(coordinates, viewer);CHKERRQ(ierr);}
6874d0cae8SMatthew G. Knepley   ierr = PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject) coordinates, "Np", PETSC_INT, (void *) &Np);CHKERRQ(ierr);
6974d0cae8SMatthew G. Knepley   ierr = PetscViewerHDF5PopGroup(viewer);CHKERRQ(ierr);
7074d0cae8SMatthew G. Knepley   ierr = DMSwarmDestroyGlobalVectorFromField(dm, DMSwarmPICField_coor, &coordinates);CHKERRQ(ierr);
7174d0cae8SMatthew G. Knepley   PetscFunctionReturn(0);
7274d0cae8SMatthew G. Knepley }
7374d0cae8SMatthew G. Knepley #endif
7474d0cae8SMatthew G. Knepley 
7574d0cae8SMatthew G. Knepley PetscErrorCode VecView_Swarm(Vec v, PetscViewer viewer)
7674d0cae8SMatthew G. Knepley {
7774d0cae8SMatthew G. Knepley   DM             dm;
7874d0cae8SMatthew G. Knepley   PetscBool      ishdf5;
7974d0cae8SMatthew G. Knepley   PetscErrorCode ierr;
8074d0cae8SMatthew G. Knepley 
8174d0cae8SMatthew G. Knepley   PetscFunctionBegin;
8274d0cae8SMatthew G. Knepley   ierr = VecGetDM(v, &dm);CHKERRQ(ierr);
8374d0cae8SMatthew G. Knepley   if (!dm) SETERRQ(PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_WRONG, "Vector not generated from a DM");
8474d0cae8SMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERHDF5,  &ishdf5);CHKERRQ(ierr);
8574d0cae8SMatthew G. Knepley   if (ishdf5) {
8674d0cae8SMatthew G. Knepley #if defined(PETSC_HAVE_HDF5)
8774d0cae8SMatthew G. Knepley     ierr = VecView_Swarm_HDF5_Internal(v, viewer);CHKERRQ(ierr);
8874d0cae8SMatthew G. Knepley #else
8974d0cae8SMatthew G. Knepley     SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
9074d0cae8SMatthew G. Knepley #endif
9174d0cae8SMatthew G. Knepley   } else {
9274d0cae8SMatthew G. Knepley     PetscBool isseq;
9374d0cae8SMatthew G. Knepley 
9474d0cae8SMatthew G. Knepley     ierr = PetscObjectTypeCompare((PetscObject) v, VECSEQ, &isseq);CHKERRQ(ierr);
9574d0cae8SMatthew G. Knepley     if (isseq) {ierr = VecView_Seq(v, viewer);CHKERRQ(ierr);}
9674d0cae8SMatthew G. Knepley     else       {ierr = VecView_MPI(v, viewer);CHKERRQ(ierr);}
9774d0cae8SMatthew G. Knepley   }
9874d0cae8SMatthew G. Knepley   PetscFunctionReturn(0);
9974d0cae8SMatthew G. Knepley }
10074d0cae8SMatthew G. Knepley 
101d3a51819SDave May /*@C
10262741f57SDave May    DMSwarmVectorDefineField - Sets the field from which to define a Vec object
10362741f57SDave May                              when DMCreateLocalVector(), or DMCreateGlobalVector() is called
10457795646SDave May 
105d083f849SBarry Smith    Collective on dm
10657795646SDave May 
107d3a51819SDave May    Input parameters:
10862741f57SDave May +  dm - a DMSwarm
10962741f57SDave May -  fieldname - the textual name given to a registered field
11057795646SDave May 
111d3a51819SDave May    Level: beginner
11257795646SDave May 
113d3a51819SDave May    Notes:
114e7af74c8SDave May 
11562741f57SDave May    The field with name fieldname must be defined as having a data type of PetscScalar.
116e7af74c8SDave May 
117d3a51819SDave May    This function must be called prior to calling DMCreateLocalVector(), DMCreateGlobalVector().
118d3a51819SDave May    Mutiple calls to DMSwarmVectorDefineField() are permitted.
11957795646SDave May 
1208b8a3813SDave May .seealso: DMSwarmRegisterPetscDatatypeField(), DMCreateGlobalVector(), DMCreateLocalVector()
121d3a51819SDave May @*/
12274d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmVectorDefineField(DM dm,const char fieldname[])
123b5bcf523SDave May {
124b5bcf523SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
125b5bcf523SDave May   PetscErrorCode ierr;
126b5bcf523SDave May   PetscInt       bs,n;
127b5bcf523SDave May   PetscScalar    *array;
128b5bcf523SDave May   PetscDataType  type;
129b5bcf523SDave May 
130a9cbaee5SMatthew G. Knepley   PetscFunctionBegin;
1313454631fSDave May   if (!swarm->issetup) { ierr = DMSetUp(dm);CHKERRQ(ierr); }
13277048351SPatrick Sanan   ierr = DMSwarmDataBucketGetSizes(swarm->db,&n,NULL,NULL);CHKERRQ(ierr);
133b5bcf523SDave May   ierr = DMSwarmGetField(dm,fieldname,&bs,&type,(void**)&array);CHKERRQ(ierr);
134b5bcf523SDave May 
135b5bcf523SDave May   /* Check all fields are of type PETSC_REAL or PETSC_SCALAR */
136b5bcf523SDave May   if (type != PETSC_REAL) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only valid for PETSC_REAL");
137521f74f9SMatthew G. Knepley   ierr = PetscSNPrintf(swarm->vec_field_name,PETSC_MAX_PATH_LEN-1,"%s",fieldname);CHKERRQ(ierr);
138b5bcf523SDave May   swarm->vec_field_set = PETSC_TRUE;
1391b1ea282SDave May   swarm->vec_field_bs = bs;
140b5bcf523SDave May   swarm->vec_field_nlocal = n;
141dcf43ee8SDave May   ierr = DMSwarmRestoreField(dm,fieldname,&bs,&type,(void**)&array);CHKERRQ(ierr);
142b5bcf523SDave May   PetscFunctionReturn(0);
143b5bcf523SDave May }
144b5bcf523SDave May 
145cc651181SDave May /* requires DMSwarmDefineFieldVector has been called */
146b5bcf523SDave May PetscErrorCode DMCreateGlobalVector_Swarm(DM dm,Vec *vec)
147b5bcf523SDave May {
148b5bcf523SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
149b5bcf523SDave May   PetscErrorCode ierr;
150b5bcf523SDave May   Vec            x;
151b5bcf523SDave May   char           name[PETSC_MAX_PATH_LEN];
152b5bcf523SDave May 
153a9cbaee5SMatthew G. Knepley   PetscFunctionBegin;
1543454631fSDave May   if (!swarm->issetup) { ierr = DMSetUp(dm);CHKERRQ(ierr); }
155b5bcf523SDave May   if (!swarm->vec_field_set) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMSwarmVectorDefineField first");
156cc651181SDave May   if (swarm->vec_field_nlocal != swarm->db->L) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"DMSwarm sizes have changed since last call to VectorDefineField first"); /* Stale data */
157cc651181SDave May 
158521f74f9SMatthew G. Knepley   ierr = PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"DMSwarmField_%s",swarm->vec_field_name);CHKERRQ(ierr);
159b5bcf523SDave May   ierr = VecCreate(PetscObjectComm((PetscObject)dm),&x);CHKERRQ(ierr);
160b5bcf523SDave May   ierr = PetscObjectSetName((PetscObject)x,name);CHKERRQ(ierr);
1611b1ea282SDave May   ierr = VecSetSizes(x,swarm->db->L*swarm->vec_field_bs,PETSC_DETERMINE);CHKERRQ(ierr);
162b5bcf523SDave May   ierr = VecSetBlockSize(x,swarm->vec_field_bs);CHKERRQ(ierr);
163c6b011d8SStefano Zampini   ierr = VecSetDM(x,dm);CHKERRQ(ierr);
164b5bcf523SDave May   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
16574d0cae8SMatthew G. Knepley   ierr = VecSetDM(x, dm);CHKERRQ(ierr);
16674d0cae8SMatthew G. Knepley   ierr = VecSetOperation(x, VECOP_VIEW, (void (*)(void)) VecView_Swarm);CHKERRQ(ierr);
167b5bcf523SDave May   *vec = x;
168b5bcf523SDave May   PetscFunctionReturn(0);
169b5bcf523SDave May }
170b5bcf523SDave May 
171b5bcf523SDave May /* requires DMSwarmDefineFieldVector has been called */
172b5bcf523SDave May PetscErrorCode DMCreateLocalVector_Swarm(DM dm,Vec *vec)
173b5bcf523SDave May {
174b5bcf523SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
175b5bcf523SDave May   PetscErrorCode ierr;
176b5bcf523SDave May   Vec x;
177b5bcf523SDave May   char name[PETSC_MAX_PATH_LEN];
178b5bcf523SDave May 
179a9cbaee5SMatthew G. Knepley   PetscFunctionBegin;
1803454631fSDave May   if (!swarm->issetup) { ierr = DMSetUp(dm);CHKERRQ(ierr); }
181b5bcf523SDave May   if (!swarm->vec_field_set) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMSwarmVectorDefineField first");
182cc651181SDave May   if (swarm->vec_field_nlocal != swarm->db->L) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"DMSwarm sizes have changed since last call to VectorDefineField first"); /* Stale data */
183cc651181SDave May 
184521f74f9SMatthew G. Knepley   ierr = PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"DMSwarmField_%s",swarm->vec_field_name);CHKERRQ(ierr);
185b5bcf523SDave May   ierr = VecCreate(PETSC_COMM_SELF,&x);CHKERRQ(ierr);
186b5bcf523SDave May   ierr = PetscObjectSetName((PetscObject)x,name);CHKERRQ(ierr);
187071900c8SMatthew G. Knepley   ierr = VecSetSizes(x,swarm->db->L*swarm->vec_field_bs,PETSC_DETERMINE);CHKERRQ(ierr);
188b5bcf523SDave May   ierr = VecSetBlockSize(x,swarm->vec_field_bs);CHKERRQ(ierr);
189c6b011d8SStefano Zampini   ierr = VecSetDM(x,dm);CHKERRQ(ierr);
190b5bcf523SDave May   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
191b5bcf523SDave May   *vec = x;
192b5bcf523SDave May   PetscFunctionReturn(0);
193b5bcf523SDave May }
194b5bcf523SDave May 
195fb1bcc12SMatthew G. Knepley static PetscErrorCode DMSwarmDestroyVectorFromField_Private(DM dm, const char fieldname[], Vec *vec)
196fb1bcc12SMatthew G. Knepley {
197fb1bcc12SMatthew G. Knepley   DM_Swarm      *swarm = (DM_Swarm *) dm->data;
19877048351SPatrick Sanan   DMSwarmDataField      gfield;
199fb1bcc12SMatthew G. Knepley   void         (*fptr)(void);
200fb1bcc12SMatthew G. Knepley   PetscInt       bs, nlocal;
201fb1bcc12SMatthew G. Knepley   char           name[PETSC_MAX_PATH_LEN];
202fb1bcc12SMatthew G. Knepley   PetscErrorCode ierr;
203d3a51819SDave May 
204fb1bcc12SMatthew G. Knepley   PetscFunctionBegin;
205fb1bcc12SMatthew G. Knepley   ierr = VecGetLocalSize(*vec, &nlocal);CHKERRQ(ierr);
206fb1bcc12SMatthew G. Knepley   ierr = VecGetBlockSize(*vec, &bs);CHKERRQ(ierr);
207fb1bcc12SMatthew G. Knepley   if (nlocal/bs != swarm->db->L) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"DMSwarm sizes have changed since vector was created - cannot ensure pointers are valid"); /* Stale data */
20877048351SPatrick Sanan   ierr = DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield);CHKERRQ(ierr);
209fb1bcc12SMatthew G. Knepley   /* check vector is an inplace array */
210521f74f9SMatthew G. Knepley   ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN-1, "DMSwarm_VecFieldInPlace_%s", fieldname);CHKERRQ(ierr);
211fb1bcc12SMatthew G. Knepley   ierr = PetscObjectQueryFunction((PetscObject) *vec, name, &fptr);CHKERRQ(ierr);
212fb1bcc12SMatthew G. Knepley   if (!fptr) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_USER, "Vector being destroyed was not created from DMSwarm field(%s)", fieldname);
21377048351SPatrick Sanan   ierr = DMSwarmDataFieldRestoreAccess(gfield);CHKERRQ(ierr);
214fb1bcc12SMatthew G. Knepley   ierr = VecDestroy(vec);CHKERRQ(ierr);
215fb1bcc12SMatthew G. Knepley   PetscFunctionReturn(0);
216fb1bcc12SMatthew G. Knepley }
217fb1bcc12SMatthew G. Knepley 
218fb1bcc12SMatthew G. Knepley static PetscErrorCode DMSwarmCreateVectorFromField_Private(DM dm, const char fieldname[], MPI_Comm comm, Vec *vec)
219fb1bcc12SMatthew G. Knepley {
220fb1bcc12SMatthew G. Knepley   DM_Swarm      *swarm = (DM_Swarm *) dm->data;
221fb1bcc12SMatthew G. Knepley   PetscDataType  type;
222fb1bcc12SMatthew G. Knepley   PetscScalar   *array;
223fb1bcc12SMatthew G. Knepley   PetscInt       bs, n;
224fb1bcc12SMatthew G. Knepley   char           name[PETSC_MAX_PATH_LEN];
225e4fbd051SBarry Smith   PetscMPIInt    size;
226fb1bcc12SMatthew G. Knepley   PetscErrorCode ierr;
227fb1bcc12SMatthew G. Knepley 
228fb1bcc12SMatthew G. Knepley   PetscFunctionBegin;
229fb1bcc12SMatthew G. Knepley   if (!swarm->issetup) {ierr = DMSetUp(dm);CHKERRQ(ierr);}
23077048351SPatrick Sanan   ierr = DMSwarmDataBucketGetSizes(swarm->db, &n, NULL, NULL);CHKERRQ(ierr);
231fb1bcc12SMatthew G. Knepley   ierr = DMSwarmGetField(dm, fieldname, &bs, &type, (void **) &array);CHKERRQ(ierr);
232fb1bcc12SMatthew G. Knepley   if (type != PETSC_REAL) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Only valid for PETSC_REAL");
233fb1bcc12SMatthew G. Knepley 
234e4fbd051SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
235e4fbd051SBarry Smith   if (size == 1) {
236fb1bcc12SMatthew G. Knepley     ierr = VecCreateSeqWithArray(comm, bs, n*bs, array, vec);CHKERRQ(ierr);
237fb1bcc12SMatthew G. Knepley   } else {
238fb1bcc12SMatthew G. Knepley     ierr = VecCreateMPIWithArray(comm, bs, n*bs, PETSC_DETERMINE, array, vec);CHKERRQ(ierr);
239fb1bcc12SMatthew G. Knepley   }
240fb1bcc12SMatthew G. Knepley   ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN-1, "DMSwarmSharedField_%s", fieldname);CHKERRQ(ierr);
241fb1bcc12SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) *vec, name);CHKERRQ(ierr);
242fb1bcc12SMatthew G. Knepley 
243fb1bcc12SMatthew G. Knepley   /* Set guard */
244fb1bcc12SMatthew G. Knepley   ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN-1, "DMSwarm_VecFieldInPlace_%s", fieldname);CHKERRQ(ierr);
245fb1bcc12SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject) *vec, name, DMSwarmDestroyVectorFromField_Private);CHKERRQ(ierr);
24674d0cae8SMatthew G. Knepley 
24774d0cae8SMatthew G. Knepley   ierr = VecSetDM(*vec, dm);CHKERRQ(ierr);
24874d0cae8SMatthew G. Knepley   ierr = VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Swarm);CHKERRQ(ierr);
249fb1bcc12SMatthew G. Knepley   PetscFunctionReturn(0);
250fb1bcc12SMatthew G. Knepley }
251fb1bcc12SMatthew G. Knepley 
2520643ed39SMatthew G. Knepley /* This creates a "mass matrix" between a finite element and particle space. If a finite element interpolant is given by
2530643ed39SMatthew G. Knepley 
2540643ed39SMatthew G. Knepley      \hat f = \sum_i f_i \phi_i
2550643ed39SMatthew G. Knepley 
2560643ed39SMatthew G. Knepley    and a particle function is given by
2570643ed39SMatthew G. Knepley 
2580643ed39SMatthew G. Knepley      f = \sum_i w_i \delta(x - x_i)
2590643ed39SMatthew G. Knepley 
2600643ed39SMatthew G. Knepley    then we want to require that
2610643ed39SMatthew G. Knepley 
2620643ed39SMatthew G. Knepley      M \hat f = M_p f
2630643ed39SMatthew G. Knepley 
2640643ed39SMatthew G. Knepley    where the particle mass matrix is given by
2650643ed39SMatthew G. Knepley 
2660643ed39SMatthew G. Knepley      (M_p)_{ij} = \int \phi_i \delta(x - x_j)
2670643ed39SMatthew G. Knepley 
2680643ed39SMatthew G. Knepley    The way Dave May does particles, they amount to quadratue weights rather than delta functions, so he has |J| is in
2690643ed39SMatthew G. Knepley    his integral. We allow this with the boolean flag.
2700643ed39SMatthew G. Knepley */
2710643ed39SMatthew G. Knepley static PetscErrorCode DMSwarmComputeMassMatrix_Private(DM dmc, DM dmf, Mat mass, PetscBool useDeltaFunction, void *ctx)
27283c47955SMatthew G. Knepley {
27383c47955SMatthew G. Knepley   const char    *name = "Mass Matrix";
2740643ed39SMatthew G. Knepley   MPI_Comm       comm;
27583c47955SMatthew G. Knepley   PetscDS        prob;
27683c47955SMatthew G. Knepley   PetscSection   fsection, globalFSection;
277e8f14785SLisandro Dalcin   PetscHSetIJ    ht;
2780643ed39SMatthew G. Knepley   PetscLayout    rLayout, colLayout;
27983c47955SMatthew G. Knepley   PetscInt      *dnz, *onz;
280adb2528bSMark Adams   PetscInt       locRows, locCols, rStart, colStart, colEnd, *rowIDXs;
2810643ed39SMatthew G. Knepley   PetscReal     *xi, *v0, *J, *invJ, detJ = 1.0, v0ref[3] = {-1.0, -1.0, -1.0};
28283c47955SMatthew G. Knepley   PetscScalar   *elemMat;
2830643ed39SMatthew G. Knepley   PetscInt       dim, Nf, field, cStart, cEnd, cell, totDim, maxC = 0;
28483c47955SMatthew G. Knepley   PetscErrorCode ierr;
28583c47955SMatthew G. Knepley 
28683c47955SMatthew G. Knepley   PetscFunctionBegin;
2870643ed39SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) mass, &comm);CHKERRQ(ierr);
28883c47955SMatthew G. Knepley   ierr = DMGetCoordinateDim(dmf, &dim);CHKERRQ(ierr);
28983c47955SMatthew G. Knepley   ierr = DMGetDS(dmf, &prob);CHKERRQ(ierr);
29083c47955SMatthew G. Knepley   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
2910643ed39SMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
29283c47955SMatthew G. Knepley   ierr = PetscMalloc3(dim, &v0, dim*dim, &J, dim*dim,&invJ);CHKERRQ(ierr);
29392fd8e1eSJed Brown   ierr = DMGetLocalSection(dmf, &fsection);CHKERRQ(ierr);
294e87a4003SBarry Smith   ierr = DMGetGlobalSection(dmf, &globalFSection);CHKERRQ(ierr);
29583c47955SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd);CHKERRQ(ierr);
2960643ed39SMatthew G. Knepley   ierr = MatGetLocalSize(mass, &locRows, &locCols);CHKERRQ(ierr);
29783c47955SMatthew G. Knepley 
2980643ed39SMatthew G. Knepley   ierr = PetscLayoutCreate(comm, &colLayout);CHKERRQ(ierr);
2990643ed39SMatthew G. Knepley   ierr = PetscLayoutSetLocalSize(colLayout, locCols);CHKERRQ(ierr);
3000643ed39SMatthew G. Knepley   ierr = PetscLayoutSetBlockSize(colLayout, 1);CHKERRQ(ierr);
3010643ed39SMatthew G. Knepley   ierr = PetscLayoutSetUp(colLayout);CHKERRQ(ierr);
3020643ed39SMatthew G. Knepley   ierr = PetscLayoutGetRange(colLayout, &colStart, &colEnd);CHKERRQ(ierr);
3030643ed39SMatthew G. Knepley   ierr = PetscLayoutDestroy(&colLayout);CHKERRQ(ierr);
3040643ed39SMatthew G. Knepley 
3050643ed39SMatthew G. Knepley   ierr = PetscLayoutCreate(comm, &rLayout);CHKERRQ(ierr);
30683c47955SMatthew G. Knepley   ierr = PetscLayoutSetLocalSize(rLayout, locRows);CHKERRQ(ierr);
30783c47955SMatthew G. Knepley   ierr = PetscLayoutSetBlockSize(rLayout, 1);CHKERRQ(ierr);
30883c47955SMatthew G. Knepley   ierr = PetscLayoutSetUp(rLayout);CHKERRQ(ierr);
3090643ed39SMatthew G. Knepley   ierr = PetscLayoutGetRange(rLayout, &rStart, NULL);CHKERRQ(ierr);
31083c47955SMatthew G. Knepley   ierr = PetscLayoutDestroy(&rLayout);CHKERRQ(ierr);
3110643ed39SMatthew G. Knepley 
31283c47955SMatthew G. Knepley   ierr = PetscCalloc2(locRows, &dnz, locRows, &onz);CHKERRQ(ierr);
313e8f14785SLisandro Dalcin   ierr = PetscHSetIJCreate(&ht);CHKERRQ(ierr);
31453e60ab4SJoseph Pusztay 
3150643ed39SMatthew G. Knepley   ierr = PetscSynchronizedFlush(comm, NULL);CHKERRQ(ierr);
3160643ed39SMatthew G. Knepley   /* count non-zeros */
3170643ed39SMatthew G. Knepley   ierr = DMSwarmSortGetAccess(dmc);CHKERRQ(ierr);
31883c47955SMatthew G. Knepley   for (field = 0; field < Nf; ++field) {
31983c47955SMatthew G. Knepley     for (cell = cStart; cell < cEnd; ++cell) {
3200643ed39SMatthew G. Knepley       PetscInt  c, i;
3210643ed39SMatthew G. Knepley       PetscInt *findices,   *cindices; /* fine is vertices, coarse is particles */
32283c47955SMatthew G. Knepley       PetscInt  numFIndices, numCIndices;
32383c47955SMatthew G. Knepley 
32471f0bbf9SMatthew G. Knepley       ierr = DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL);CHKERRQ(ierr);
325fc7c92abSMatthew G. Knepley       ierr = DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices);CHKERRQ(ierr);
326fc7c92abSMatthew G. Knepley       maxC = PetscMax(maxC, numCIndices);
32783c47955SMatthew G. Knepley       {
328e8f14785SLisandro Dalcin         PetscHashIJKey key;
329e8f14785SLisandro Dalcin         PetscBool      missing;
33083c47955SMatthew G. Knepley         for (i = 0; i < numFIndices; ++i) {
331adb2528bSMark Adams           key.j = findices[i]; /* global column (from Plex) */
332adb2528bSMark Adams           if (key.j >= 0) {
33383c47955SMatthew G. Knepley             /* Get indices for coarse elements */
33483c47955SMatthew G. Knepley             for (c = 0; c < numCIndices; ++c) {
335adb2528bSMark Adams               key.i = cindices[c] + rStart; /* global cols (from Swarm) */
336adb2528bSMark Adams               if (key.i < 0) continue;
337e8f14785SLisandro Dalcin               ierr = PetscHSetIJQueryAdd(ht, key, &missing);CHKERRQ(ierr);
33883c47955SMatthew G. Knepley               if (missing) {
3390643ed39SMatthew G. Knepley                 if ((key.j >= colStart) && (key.j < colEnd)) ++dnz[key.i - rStart];
340e8f14785SLisandro Dalcin                 else                                         ++onz[key.i - rStart];
3410643ed39SMatthew G. Knepley               } else SETERRQ2(PetscObjectComm((PetscObject) dmf), PETSC_ERR_SUP, "Set new value at %D,%D", key.i, key.j);
34283c47955SMatthew G. Knepley             }
343fc7c92abSMatthew G. Knepley           }
344fc7c92abSMatthew G. Knepley         }
34583c47955SMatthew G. Knepley         ierr = PetscFree(cindices);CHKERRQ(ierr);
34683c47955SMatthew G. Knepley       }
34771f0bbf9SMatthew G. Knepley       ierr = DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL);CHKERRQ(ierr);
34883c47955SMatthew G. Knepley     }
34983c47955SMatthew G. Knepley   }
350e8f14785SLisandro Dalcin   ierr = PetscHSetIJDestroy(&ht);CHKERRQ(ierr);
35183c47955SMatthew G. Knepley   ierr = MatXAIJSetPreallocation(mass, 1, dnz, onz, NULL, NULL);CHKERRQ(ierr);
35283c47955SMatthew G. Knepley   ierr = MatSetOption(mass, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE);CHKERRQ(ierr);
35383c47955SMatthew G. Knepley   ierr = PetscFree2(dnz, onz);CHKERRQ(ierr);
354adb2528bSMark Adams   ierr = PetscMalloc3(maxC*totDim, &elemMat, maxC, &rowIDXs, maxC*dim, &xi);CHKERRQ(ierr);
35583c47955SMatthew G. Knepley   for (field = 0; field < Nf; ++field) {
356ef0bb6c7SMatthew G. Knepley     PetscTabulation Tcoarse;
35783c47955SMatthew G. Knepley     PetscObject       obj;
358ef0bb6c7SMatthew G. Knepley     PetscReal        *coords;
3590643ed39SMatthew G. Knepley     PetscInt          Nc, i;
36083c47955SMatthew G. Knepley 
36183c47955SMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, field, &obj);CHKERRQ(ierr);
3620643ed39SMatthew G. Knepley     ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);
3630643ed39SMatthew G. Knepley     if (Nc != 1) SETERRQ1(PetscObjectComm((PetscObject) dmf), PETSC_ERR_SUP, "Can only interpolate a scalar field from particles, Nc = %D", Nc);
3640643ed39SMatthew G. Knepley     ierr = DMSwarmGetField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr);
36583c47955SMatthew G. Knepley     for (cell = cStart; cell < cEnd; ++cell) {
36683c47955SMatthew G. Knepley       PetscInt *findices  , *cindices;
36783c47955SMatthew G. Knepley       PetscInt  numFIndices, numCIndices;
3680643ed39SMatthew G. Knepley       PetscInt  p, c;
36983c47955SMatthew G. Knepley 
3700643ed39SMatthew G. Knepley       /* TODO: Use DMField instead of assuming affine */
37183c47955SMatthew G. Knepley       ierr = DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
37271f0bbf9SMatthew G. Knepley       ierr = DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL);CHKERRQ(ierr);
37383c47955SMatthew G. Knepley       ierr = DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices);CHKERRQ(ierr);
3740643ed39SMatthew G. Knepley       for (p = 0; p < numCIndices; ++p) {
3750643ed39SMatthew G. Knepley         CoordinatesRealToRef(dim, dim, v0ref, v0, invJ, &coords[cindices[p]*dim], &xi[p*dim]);
3760643ed39SMatthew G. Knepley       }
377ef0bb6c7SMatthew G. Knepley       ierr = PetscFECreateTabulation((PetscFE) obj, 1, numCIndices, xi, 0, &Tcoarse);CHKERRQ(ierr);
37883c47955SMatthew G. Knepley       /* Get elemMat entries by multiplying by weight */
379580bdb30SBarry Smith       ierr = PetscArrayzero(elemMat, numCIndices*totDim);CHKERRQ(ierr);
38083c47955SMatthew G. Knepley       for (i = 0; i < numFIndices; ++i) {
3810643ed39SMatthew G. Knepley         for (p = 0; p < numCIndices; ++p) {
3820643ed39SMatthew G. Knepley           for (c = 0; c < Nc; ++c) {
3830643ed39SMatthew G. Knepley             /* B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c */
384ef0bb6c7SMatthew G. Knepley             elemMat[p*numFIndices+i] += Tcoarse->T[0][(p*numFIndices + i)*Nc + c]*(useDeltaFunction ? 1.0 : detJ);
38583c47955SMatthew G. Knepley           }
3860643ed39SMatthew G. Knepley         }
3870643ed39SMatthew G. Knepley       }
388adb2528bSMark Adams       for (p = 0; p < numCIndices; ++p) rowIDXs[p] = cindices[p] + rStart;
38983c47955SMatthew G. Knepley       if (0) {ierr = DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat);CHKERRQ(ierr);}
390adb2528bSMark Adams       ierr = MatSetValues(mass, numCIndices, rowIDXs, numFIndices, findices, elemMat, ADD_VALUES);CHKERRQ(ierr);
39183c47955SMatthew G. Knepley       ierr = PetscFree(cindices);CHKERRQ(ierr);
39271f0bbf9SMatthew G. Knepley       ierr = DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL);CHKERRQ(ierr);
393ef0bb6c7SMatthew G. Knepley       ierr = PetscTabulationDestroy(&Tcoarse);CHKERRQ(ierr);
39483c47955SMatthew G. Knepley     }
3950643ed39SMatthew G. Knepley     ierr = DMSwarmRestoreField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr);
39683c47955SMatthew G. Knepley   }
397adb2528bSMark Adams   ierr = PetscFree3(elemMat, rowIDXs, xi);CHKERRQ(ierr);
39883c47955SMatthew G. Knepley   ierr = DMSwarmSortRestoreAccess(dmc);CHKERRQ(ierr);
39983c47955SMatthew G. Knepley   ierr = PetscFree3(v0, J, invJ);CHKERRQ(ierr);
40083c47955SMatthew G. Knepley   ierr = MatAssemblyBegin(mass, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
40183c47955SMatthew G. Knepley   ierr = MatAssemblyEnd(mass, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
40283c47955SMatthew G. Knepley   PetscFunctionReturn(0);
40383c47955SMatthew G. Knepley }
40483c47955SMatthew G. Knepley 
405adb2528bSMark Adams /* FEM cols, Particle rows */
40683c47955SMatthew G. Knepley static PetscErrorCode DMCreateMassMatrix_Swarm(DM dmCoarse, DM dmFine, Mat *mass)
40783c47955SMatthew G. Knepley {
408895a1698SMatthew G. Knepley   PetscSection   gsf;
40983c47955SMatthew G. Knepley   PetscInt       m, n;
41083c47955SMatthew G. Knepley   void          *ctx;
41183c47955SMatthew G. Knepley   PetscErrorCode ierr;
41283c47955SMatthew G. Knepley 
41383c47955SMatthew G. Knepley   PetscFunctionBegin;
414e87a4003SBarry Smith   ierr = DMGetGlobalSection(dmFine, &gsf);CHKERRQ(ierr);
41583c47955SMatthew G. Knepley   ierr = PetscSectionGetConstrainedStorageSize(gsf, &m);CHKERRQ(ierr);
416895a1698SMatthew G. Knepley   ierr = DMSwarmGetLocalSize(dmCoarse, &n);CHKERRQ(ierr);
41783c47955SMatthew G. Knepley   ierr = MatCreate(PetscObjectComm((PetscObject) dmCoarse), mass);CHKERRQ(ierr);
418adb2528bSMark Adams   ierr = MatSetSizes(*mass, n, m, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
41983c47955SMatthew G. Knepley   ierr = MatSetType(*mass, dmCoarse->mattype);CHKERRQ(ierr);
42083c47955SMatthew G. Knepley   ierr = DMGetApplicationContext(dmFine, &ctx);CHKERRQ(ierr);
42183c47955SMatthew G. Knepley 
4220643ed39SMatthew G. Knepley   ierr = DMSwarmComputeMassMatrix_Private(dmCoarse, dmFine, *mass, PETSC_TRUE, ctx);CHKERRQ(ierr);
42383c47955SMatthew G. Knepley   ierr = MatViewFromOptions(*mass, NULL, "-mass_mat_view");CHKERRQ(ierr);
42483c47955SMatthew G. Knepley   PetscFunctionReturn(0);
42583c47955SMatthew G. Knepley }
42683c47955SMatthew G. Knepley 
427fb1bcc12SMatthew G. Knepley /*@C
428d3a51819SDave May    DMSwarmCreateGlobalVectorFromField - Creates a Vec object sharing the array associated with a given field
429d3a51819SDave May 
430d083f849SBarry Smith    Collective on dm
431d3a51819SDave May 
432d3a51819SDave May    Input parameters:
43362741f57SDave May +  dm - a DMSwarm
43462741f57SDave May -  fieldname - the textual name given to a registered field
435d3a51819SDave May 
4368b8a3813SDave May    Output parameter:
437d3a51819SDave May .  vec - the vector
438d3a51819SDave May 
439d3a51819SDave May    Level: beginner
440d3a51819SDave May 
4418b8a3813SDave May    Notes:
4428b8a3813SDave May    The vector must be returned using a matching call to DMSwarmDestroyGlobalVectorFromField().
4438b8a3813SDave May 
4448b8a3813SDave May .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmDestroyGlobalVectorFromField()
445d3a51819SDave May @*/
44674d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmCreateGlobalVectorFromField(DM dm,const char fieldname[],Vec *vec)
447b5bcf523SDave May {
448fb1bcc12SMatthew G. Knepley   MPI_Comm       comm = PetscObjectComm((PetscObject) dm);
449b5bcf523SDave May   PetscErrorCode ierr;
450b5bcf523SDave May 
451fb1bcc12SMatthew G. Knepley   PetscFunctionBegin;
452fb1bcc12SMatthew G. Knepley   ierr = DMSwarmCreateVectorFromField_Private(dm, fieldname, comm, vec);CHKERRQ(ierr);
453b5bcf523SDave May   PetscFunctionReturn(0);
454b5bcf523SDave May }
455b5bcf523SDave May 
456d3a51819SDave May /*@C
457d3a51819SDave May    DMSwarmDestroyGlobalVectorFromField - Destroys the Vec object which share the array associated with a given field
458d3a51819SDave May 
459d083f849SBarry Smith    Collective on dm
460d3a51819SDave May 
461d3a51819SDave May    Input parameters:
46262741f57SDave May +  dm - a DMSwarm
46362741f57SDave May -  fieldname - the textual name given to a registered field
464d3a51819SDave May 
4658b8a3813SDave May    Output parameter:
466d3a51819SDave May .  vec - the vector
467d3a51819SDave May 
468d3a51819SDave May    Level: beginner
469d3a51819SDave May 
4708b8a3813SDave May .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmCreateGlobalVectorFromField()
471d3a51819SDave May @*/
47274d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmDestroyGlobalVectorFromField(DM dm,const char fieldname[],Vec *vec)
473b5bcf523SDave May {
474b5bcf523SDave May   PetscErrorCode ierr;
475cc651181SDave May 
476fb1bcc12SMatthew G. Knepley   PetscFunctionBegin;
477fb1bcc12SMatthew G. Knepley   ierr = DMSwarmDestroyVectorFromField_Private(dm, fieldname, vec);CHKERRQ(ierr);
478b5bcf523SDave May   PetscFunctionReturn(0);
479b5bcf523SDave May }
480b5bcf523SDave May 
481fb1bcc12SMatthew G. Knepley /*@C
482fb1bcc12SMatthew G. Knepley    DMSwarmCreateLocalVectorFromField - Creates a Vec object sharing the array associated with a given field
483fb1bcc12SMatthew G. Knepley 
484d083f849SBarry Smith    Collective on dm
485fb1bcc12SMatthew G. Knepley 
486fb1bcc12SMatthew G. Knepley    Input parameters:
48762741f57SDave May +  dm - a DMSwarm
48862741f57SDave May -  fieldname - the textual name given to a registered field
489fb1bcc12SMatthew G. Knepley 
4908b8a3813SDave May    Output parameter:
491fb1bcc12SMatthew G. Knepley .  vec - the vector
492fb1bcc12SMatthew G. Knepley 
493fb1bcc12SMatthew G. Knepley    Level: beginner
494fb1bcc12SMatthew G. Knepley 
4958b8a3813SDave May    Notes:
4968b8a3813SDave May    The vector must be returned using a matching call to DMSwarmDestroyLocalVectorFromField().
4978b8a3813SDave May 
4988b8a3813SDave May .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmDestroyLocalVectorFromField()
499fb1bcc12SMatthew G. Knepley @*/
50074d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmCreateLocalVectorFromField(DM dm,const char fieldname[],Vec *vec)
501bbe8250bSMatthew G. Knepley {
502fb1bcc12SMatthew G. Knepley   MPI_Comm       comm = PETSC_COMM_SELF;
503bbe8250bSMatthew G. Knepley   PetscErrorCode ierr;
504bbe8250bSMatthew G. Knepley 
505fb1bcc12SMatthew G. Knepley   PetscFunctionBegin;
506fb1bcc12SMatthew G. Knepley   ierr = DMSwarmCreateVectorFromField_Private(dm, fieldname, comm, vec);CHKERRQ(ierr);
507fb1bcc12SMatthew G. Knepley   PetscFunctionReturn(0);
508bbe8250bSMatthew G. Knepley }
509fb1bcc12SMatthew G. Knepley 
510fb1bcc12SMatthew G. Knepley /*@C
511fb1bcc12SMatthew G. Knepley    DMSwarmDestroyLocalVectorFromField - Destroys the Vec object which share the array associated with a given field
512fb1bcc12SMatthew G. Knepley 
513d083f849SBarry Smith    Collective on dm
514fb1bcc12SMatthew G. Knepley 
515fb1bcc12SMatthew G. Knepley    Input parameters:
51662741f57SDave May +  dm - a DMSwarm
51762741f57SDave May -  fieldname - the textual name given to a registered field
518fb1bcc12SMatthew G. Knepley 
5198b8a3813SDave May    Output parameter:
520fb1bcc12SMatthew G. Knepley .  vec - the vector
521fb1bcc12SMatthew G. Knepley 
522fb1bcc12SMatthew G. Knepley    Level: beginner
523fb1bcc12SMatthew G. Knepley 
5248b8a3813SDave May .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmCreateLocalVectorFromField()
525fb1bcc12SMatthew G. Knepley @*/
52674d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmDestroyLocalVectorFromField(DM dm,const char fieldname[],Vec *vec)
527fb1bcc12SMatthew G. Knepley {
528fb1bcc12SMatthew G. Knepley   PetscErrorCode ierr;
529fb1bcc12SMatthew G. Knepley 
530fb1bcc12SMatthew G. Knepley   PetscFunctionBegin;
531fb1bcc12SMatthew G. Knepley   ierr = DMSwarmDestroyVectorFromField_Private(dm, fieldname, vec);CHKERRQ(ierr);
532bbe8250bSMatthew G. Knepley   PetscFunctionReturn(0);
533bbe8250bSMatthew G. Knepley }
534bbe8250bSMatthew G. Knepley 
535b5bcf523SDave May /*
53674d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmCreateGlobalVectorFromFields(DM dm,const PetscInt nf,const char *fieldnames[],Vec *vec)
537b5bcf523SDave May {
538b5bcf523SDave May   PetscFunctionReturn(0);
539b5bcf523SDave May }
540b5bcf523SDave May 
54174d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmRestoreGlobalVectorFromFields(DM dm,Vec *vec)
542b5bcf523SDave May {
543b5bcf523SDave May   PetscFunctionReturn(0);
544b5bcf523SDave May }
545b5bcf523SDave May */
546b5bcf523SDave May 
547d3a51819SDave May /*@C
548d3a51819SDave May    DMSwarmInitializeFieldRegister - Initiates the registration of fields to a DMSwarm
549d3a51819SDave May 
550d083f849SBarry Smith    Collective on dm
551d3a51819SDave May 
552d3a51819SDave May    Input parameter:
553d3a51819SDave May .  dm - a DMSwarm
554d3a51819SDave May 
555d3a51819SDave May    Level: beginner
556d3a51819SDave May 
557d3a51819SDave May    Notes:
5588b8a3813SDave May    After all fields have been registered, you must call DMSwarmFinalizeFieldRegister().
559d3a51819SDave May 
560d3a51819SDave May .seealso: DMSwarmFinalizeFieldRegister(), DMSwarmRegisterPetscDatatypeField(),
561d3a51819SDave May  DMSwarmRegisterUserStructField(), DMSwarmRegisterUserDatatypeField()
562d3a51819SDave May @*/
56374d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmInitializeFieldRegister(DM dm)
5645f50eb2eSDave May {
5655f50eb2eSDave May   DM_Swarm      *swarm = (DM_Swarm *) dm->data;
5663454631fSDave May   PetscErrorCode ierr;
5673454631fSDave May 
568521f74f9SMatthew G. Knepley   PetscFunctionBegin;
569cc651181SDave May   if (!swarm->field_registration_initialized) {
5705f50eb2eSDave May     swarm->field_registration_initialized = PETSC_TRUE;
57143f984edSMatthew G. Knepley     ierr = DMSwarmRegisterPetscDatatypeField(dm,DMSwarmField_pid,1,PETSC_INT64);CHKERRQ(ierr); /* unique identifer */
572f0cdbbbaSDave May     ierr = DMSwarmRegisterPetscDatatypeField(dm,DMSwarmField_rank,1,PETSC_INT);CHKERRQ(ierr); /* used for communication */
573cc651181SDave May   }
5745f50eb2eSDave May   PetscFunctionReturn(0);
5755f50eb2eSDave May }
5765f50eb2eSDave May 
57774d0cae8SMatthew G. Knepley /*@
578d3a51819SDave May    DMSwarmFinalizeFieldRegister - Finalizes the registration of fields to a DMSwarm
579d3a51819SDave May 
580d083f849SBarry Smith    Collective on dm
581d3a51819SDave May 
582d3a51819SDave May    Input parameter:
583d3a51819SDave May .  dm - a DMSwarm
584d3a51819SDave May 
585d3a51819SDave May    Level: beginner
586d3a51819SDave May 
587d3a51819SDave May    Notes:
58862741f57SDave May    After DMSwarmFinalizeFieldRegister() has been called, no new fields can be defined on the DMSwarm.
589d3a51819SDave May 
590d3a51819SDave May .seealso: DMSwarmInitializeFieldRegister(), DMSwarmRegisterPetscDatatypeField(),
591d3a51819SDave May  DMSwarmRegisterUserStructField(), DMSwarmRegisterUserDatatypeField()
592d3a51819SDave May @*/
59374d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmFinalizeFieldRegister(DM dm)
5945f50eb2eSDave May {
5955f50eb2eSDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
5966845f8f5SDave May   PetscErrorCode ierr;
5976845f8f5SDave May 
598521f74f9SMatthew G. Knepley   PetscFunctionBegin;
599f0cdbbbaSDave May   if (!swarm->field_registration_finalized) {
60077048351SPatrick Sanan     ierr = DMSwarmDataBucketFinalize(swarm->db);CHKERRQ(ierr);
601f0cdbbbaSDave May   }
602f0cdbbbaSDave May   swarm->field_registration_finalized = PETSC_TRUE;
6035f50eb2eSDave May   PetscFunctionReturn(0);
6045f50eb2eSDave May }
6055f50eb2eSDave May 
60674d0cae8SMatthew G. Knepley /*@
607d3a51819SDave May    DMSwarmSetLocalSizes - Sets the length of all registered fields on the DMSwarm
608d3a51819SDave May 
609d3a51819SDave May    Not collective
610d3a51819SDave May 
611d3a51819SDave May    Input parameters:
61262741f57SDave May +  dm - a DMSwarm
613d3a51819SDave May .  nlocal - the length of each registered field
61462741f57SDave May -  buffer - the length of the buffer used to efficient dynamic re-sizing
615d3a51819SDave May 
616d3a51819SDave May    Level: beginner
617d3a51819SDave May 
618d3a51819SDave May .seealso: DMSwarmGetLocalSize()
619d3a51819SDave May @*/
62074d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmSetLocalSizes(DM dm,PetscInt nlocal,PetscInt buffer)
6215f50eb2eSDave May {
6225f50eb2eSDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
6236845f8f5SDave May   PetscErrorCode ierr;
6245f50eb2eSDave May 
625521f74f9SMatthew G. Knepley   PetscFunctionBegin;
626f2b2bee7SDave May   ierr = PetscLogEventBegin(DMSWARM_SetSizes,0,0,0,0);CHKERRQ(ierr);
62777048351SPatrick Sanan   ierr = DMSwarmDataBucketSetSizes(swarm->db,nlocal,buffer);CHKERRQ(ierr);
628f2b2bee7SDave May   ierr = PetscLogEventEnd(DMSWARM_SetSizes,0,0,0,0);CHKERRQ(ierr);
6295f50eb2eSDave May   PetscFunctionReturn(0);
6305f50eb2eSDave May }
6315f50eb2eSDave May 
63274d0cae8SMatthew G. Knepley /*@
633d3a51819SDave May    DMSwarmSetCellDM - Attachs a DM to a DMSwarm
634d3a51819SDave May 
635d083f849SBarry Smith    Collective on dm
636d3a51819SDave May 
637d3a51819SDave May    Input parameters:
63862741f57SDave May +  dm - a DMSwarm
63962741f57SDave May -  dmcell - the DM to attach to the DMSwarm
640d3a51819SDave May 
641d3a51819SDave May    Level: beginner
642d3a51819SDave May 
643d3a51819SDave May    Notes:
644d3a51819SDave May    The attached DM (dmcell) will be queried for point location and
6458b8a3813SDave May    neighbor MPI-rank information if DMSwarmMigrate() is called.
646d3a51819SDave May 
6478b8a3813SDave May .seealso: DMSwarmSetType(), DMSwarmGetCellDM(), DMSwarmMigrate()
648d3a51819SDave May @*/
64974d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmSetCellDM(DM dm,DM dmcell)
650b16650c8SDave May {
651b16650c8SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
652521f74f9SMatthew G. Knepley 
653521f74f9SMatthew G. Knepley   PetscFunctionBegin;
654b16650c8SDave May   swarm->dmcell = dmcell;
655b16650c8SDave May   PetscFunctionReturn(0);
656b16650c8SDave May }
657b16650c8SDave May 
65874d0cae8SMatthew G. Knepley /*@
659d3a51819SDave May    DMSwarmGetCellDM - Fetches the attached cell DM
660d3a51819SDave May 
661d083f849SBarry Smith    Collective on dm
662d3a51819SDave May 
663d3a51819SDave May    Input parameter:
664d3a51819SDave May .  dm - a DMSwarm
665d3a51819SDave May 
666d3a51819SDave May    Output parameter:
667d3a51819SDave May .  dmcell - the DM which was attached to the DMSwarm
668d3a51819SDave May 
669d3a51819SDave May    Level: beginner
670d3a51819SDave May 
671d3a51819SDave May .seealso: DMSwarmSetCellDM()
672d3a51819SDave May @*/
67374d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmGetCellDM(DM dm,DM *dmcell)
674fe39f135SDave May {
675fe39f135SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
676521f74f9SMatthew G. Knepley 
677521f74f9SMatthew G. Knepley   PetscFunctionBegin;
678fe39f135SDave May   *dmcell = swarm->dmcell;
679fe39f135SDave May   PetscFunctionReturn(0);
680fe39f135SDave May }
681fe39f135SDave May 
68274d0cae8SMatthew G. Knepley /*@
683d3a51819SDave May    DMSwarmGetLocalSize - Retrives the local length of fields registered
684d3a51819SDave May 
685d3a51819SDave May    Not collective
686d3a51819SDave May 
687d3a51819SDave May    Input parameter:
688d3a51819SDave May .  dm - a DMSwarm
689d3a51819SDave May 
690d3a51819SDave May    Output parameter:
691d3a51819SDave May .  nlocal - the length of each registered field
692d3a51819SDave May 
693d3a51819SDave May    Level: beginner
694d3a51819SDave May 
6958b8a3813SDave May .seealso: DMSwarmGetSize(), DMSwarmSetLocalSizes()
696d3a51819SDave May @*/
69774d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmGetLocalSize(DM dm,PetscInt *nlocal)
698dcf43ee8SDave May {
699dcf43ee8SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
700dcf43ee8SDave May   PetscErrorCode ierr;
701dcf43ee8SDave May 
702521f74f9SMatthew G. Knepley   PetscFunctionBegin;
70377048351SPatrick Sanan   if (nlocal) {ierr = DMSwarmDataBucketGetSizes(swarm->db,nlocal,NULL,NULL);CHKERRQ(ierr);}
704dcf43ee8SDave May   PetscFunctionReturn(0);
705dcf43ee8SDave May }
706dcf43ee8SDave May 
70774d0cae8SMatthew G. Knepley /*@
708d3a51819SDave May    DMSwarmGetSize - Retrives the total length of fields registered
709d3a51819SDave May 
710d083f849SBarry Smith    Collective on dm
711d3a51819SDave May 
712d3a51819SDave May    Input parameter:
713d3a51819SDave May .  dm - a DMSwarm
714d3a51819SDave May 
715d3a51819SDave May    Output parameter:
716d3a51819SDave May .  n - the total length of each registered field
717d3a51819SDave May 
718d3a51819SDave May    Level: beginner
719d3a51819SDave May 
720d3a51819SDave May    Note:
721d3a51819SDave May    This calls MPI_Allreduce upon each call (inefficient but safe)
722d3a51819SDave May 
7238b8a3813SDave May .seealso: DMSwarmGetLocalSize(), DMSwarmSetLocalSizes()
724d3a51819SDave May @*/
72574d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmGetSize(DM dm,PetscInt *n)
726dcf43ee8SDave May {
727dcf43ee8SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
728dcf43ee8SDave May   PetscErrorCode ierr;
729dcf43ee8SDave May   PetscInt       nlocal,ng;
730dcf43ee8SDave May 
731521f74f9SMatthew G. Knepley   PetscFunctionBegin;
73277048351SPatrick Sanan   ierr = DMSwarmDataBucketGetSizes(swarm->db,&nlocal,NULL,NULL);CHKERRQ(ierr);
733dcf43ee8SDave May   ierr = MPI_Allreduce(&nlocal,&ng,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
734dcf43ee8SDave May   if (n) { *n = ng; }
735dcf43ee8SDave May   PetscFunctionReturn(0);
736dcf43ee8SDave May }
737dcf43ee8SDave May 
738d3a51819SDave May /*@C
7398b8a3813SDave May    DMSwarmRegisterPetscDatatypeField - Register a field to a DMSwarm with a native PETSc data type
740d3a51819SDave May 
741d083f849SBarry Smith    Collective on dm
742d3a51819SDave May 
743d3a51819SDave May    Input parameters:
74462741f57SDave May +  dm - a DMSwarm
745d3a51819SDave May .  fieldname - the textual name to identify this field
746d3a51819SDave May .  blocksize - the number of each data type
74762741f57SDave May -  type - a valid PETSc data type (PETSC_CHAR, PETSC_SHORT, PETSC_INT, PETSC_FLOAT, PETSC_REAL, PETSC_LONG)
748d3a51819SDave May 
749d3a51819SDave May    Level: beginner
750d3a51819SDave May 
751d3a51819SDave May    Notes:
7528b8a3813SDave May    The textual name for each registered field must be unique.
753d3a51819SDave May 
754d3a51819SDave May .seealso: DMSwarmRegisterUserStructField(), DMSwarmRegisterUserDatatypeField()
755d3a51819SDave May @*/
75674d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmRegisterPetscDatatypeField(DM dm,const char fieldname[],PetscInt blocksize,PetscDataType type)
757b62e03f8SDave May {
7582eac95f8SDave May   PetscErrorCode ierr;
759b62e03f8SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
760b62e03f8SDave May   size_t         size;
761b62e03f8SDave May 
762521f74f9SMatthew G. Knepley   PetscFunctionBegin;
7635f50eb2eSDave May   if (!swarm->field_registration_initialized) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMSwarmInitializeFieldRegister() first");
7645f50eb2eSDave May   if (swarm->field_registration_finalized) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Cannot register additional fields after calling DMSwarmFinalizeFieldRegister() first");
7655f50eb2eSDave May 
7665f50eb2eSDave May   if (type == PETSC_OBJECT) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
7675f50eb2eSDave May   if (type == PETSC_FUNCTION) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
7685f50eb2eSDave May   if (type == PETSC_STRING) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
7695f50eb2eSDave May   if (type == PETSC_STRUCT) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
7705f50eb2eSDave May   if (type == PETSC_DATATYPE_UNKNOWN) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
771b62e03f8SDave May 
7722ddcf43eSMatthew G. Knepley   ierr = PetscDataTypeGetSize(type, &size);CHKERRQ(ierr);
773b62e03f8SDave May   /* Load a specific data type into data bucket, specifying textual name and its size in bytes */
77477048351SPatrick Sanan   ierr = DMSwarmDataBucketRegisterField(swarm->db,"DMSwarmRegisterPetscDatatypeField",fieldname,blocksize*size,NULL);CHKERRQ(ierr);
77552c3ed93SDave May   {
77677048351SPatrick Sanan     DMSwarmDataField gfield;
77752c3ed93SDave May 
77877048351SPatrick Sanan     ierr = DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,fieldname,&gfield);CHKERRQ(ierr);
77977048351SPatrick Sanan     ierr = DMSwarmDataFieldSetBlockSize(gfield,blocksize);CHKERRQ(ierr);
78052c3ed93SDave May   }
781b62e03f8SDave May   swarm->db->field[swarm->db->nfields-1]->petsc_type = type;
782b62e03f8SDave May   PetscFunctionReturn(0);
783b62e03f8SDave May }
784b62e03f8SDave May 
785d3a51819SDave May /*@C
786d3a51819SDave May    DMSwarmRegisterUserStructField - Register a user defined struct to a DMSwarm
787d3a51819SDave May 
788d083f849SBarry Smith    Collective on dm
789d3a51819SDave May 
790d3a51819SDave May    Input parameters:
79162741f57SDave May +  dm - a DMSwarm
792d3a51819SDave May .  fieldname - the textual name to identify this field
79362741f57SDave May -  size - the size in bytes of the user struct of each data type
794d3a51819SDave May 
795d3a51819SDave May    Level: beginner
796d3a51819SDave May 
797d3a51819SDave May    Notes:
7988b8a3813SDave May    The textual name for each registered field must be unique.
799d3a51819SDave May 
800d3a51819SDave May .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmRegisterUserDatatypeField()
801d3a51819SDave May @*/
80274d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmRegisterUserStructField(DM dm,const char fieldname[],size_t size)
803b62e03f8SDave May {
8042eac95f8SDave May   PetscErrorCode ierr;
805b62e03f8SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
806b62e03f8SDave May 
807521f74f9SMatthew G. Knepley   PetscFunctionBegin;
80877048351SPatrick Sanan   ierr = DMSwarmDataBucketRegisterField(swarm->db,"DMSwarmRegisterUserStructField",fieldname,size,NULL);CHKERRQ(ierr);
809b62e03f8SDave May   swarm->db->field[swarm->db->nfields-1]->petsc_type = PETSC_STRUCT ;
810b62e03f8SDave May   PetscFunctionReturn(0);
811b62e03f8SDave May }
812b62e03f8SDave May 
813d3a51819SDave May /*@C
814d3a51819SDave May    DMSwarmRegisterUserDatatypeField - Register a user defined data type to a DMSwarm
815d3a51819SDave May 
816d083f849SBarry Smith    Collective on dm
817d3a51819SDave May 
818d3a51819SDave May    Input parameters:
81962741f57SDave May +  dm - a DMSwarm
820d3a51819SDave May .  fieldname - the textual name to identify this field
821d3a51819SDave May .  size - the size in bytes of the user data type
82262741f57SDave May -  blocksize - the number of each data type
823d3a51819SDave May 
824d3a51819SDave May    Level: beginner
825d3a51819SDave May 
826d3a51819SDave May    Notes:
8278b8a3813SDave May    The textual name for each registered field must be unique.
828d3a51819SDave May 
829d3a51819SDave May .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmRegisterUserStructField(), DMSwarmRegisterUserDatatypeField()
830d3a51819SDave May @*/
83174d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmRegisterUserDatatypeField(DM dm,const char fieldname[],size_t size,PetscInt blocksize)
832b62e03f8SDave May {
833b62e03f8SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
8346845f8f5SDave May   PetscErrorCode ierr;
835b62e03f8SDave May 
836521f74f9SMatthew G. Knepley   PetscFunctionBegin;
83777048351SPatrick Sanan   ierr = DMSwarmDataBucketRegisterField(swarm->db,"DMSwarmRegisterUserDatatypeField",fieldname,blocksize*size,NULL);CHKERRQ(ierr);
838320740a0SDave May   {
83977048351SPatrick Sanan     DMSwarmDataField gfield;
840320740a0SDave May 
84177048351SPatrick Sanan     ierr = DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,fieldname,&gfield);CHKERRQ(ierr);
84277048351SPatrick Sanan     ierr = DMSwarmDataFieldSetBlockSize(gfield,blocksize);CHKERRQ(ierr);
843320740a0SDave May   }
844b62e03f8SDave May   swarm->db->field[swarm->db->nfields-1]->petsc_type = PETSC_DATATYPE_UNKNOWN;
845b62e03f8SDave May   PetscFunctionReturn(0);
846b62e03f8SDave May }
847b62e03f8SDave May 
848d3a51819SDave May /*@C
849d3a51819SDave May    DMSwarmGetField - Get access to the underlying array storing all entries associated with a registered field
850d3a51819SDave May 
851d3a51819SDave May    Not collective
852d3a51819SDave May 
853d3a51819SDave May    Input parameters:
85462741f57SDave May +  dm - a DMSwarm
85562741f57SDave May -  fieldname - the textual name to identify this field
856d3a51819SDave May 
857d3a51819SDave May    Output parameters:
85862741f57SDave May +  blocksize - the number of each data type
859d3a51819SDave May .  type - the data type
86062741f57SDave May -  data - pointer to raw array
861d3a51819SDave May 
862d3a51819SDave May    Level: beginner
863d3a51819SDave May 
864d3a51819SDave May    Notes:
8658b8a3813SDave May    The array must be returned using a matching call to DMSwarmRestoreField().
866d3a51819SDave May 
867d3a51819SDave May .seealso: DMSwarmRestoreField()
868d3a51819SDave May @*/
86974d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmGetField(DM dm,const char fieldname[],PetscInt *blocksize,PetscDataType *type,void **data)
870b62e03f8SDave May {
871b62e03f8SDave May   DM_Swarm         *swarm = (DM_Swarm*)dm->data;
87277048351SPatrick Sanan   DMSwarmDataField gfield;
8732eac95f8SDave May   PetscErrorCode ierr;
874b62e03f8SDave May 
875521f74f9SMatthew G. Knepley   PetscFunctionBegin;
8763454631fSDave May   if (!swarm->issetup) { ierr = DMSetUp(dm);CHKERRQ(ierr); }
87777048351SPatrick Sanan   ierr = DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,fieldname,&gfield);CHKERRQ(ierr);
87877048351SPatrick Sanan   ierr = DMSwarmDataFieldGetAccess(gfield);CHKERRQ(ierr);
87977048351SPatrick Sanan   ierr = DMSwarmDataFieldGetEntries(gfield,data);CHKERRQ(ierr);
8801b1ea282SDave May   if (blocksize) {*blocksize = gfield->bs; }
881b5bcf523SDave May   if (type) { *type = gfield->petsc_type; }
882b62e03f8SDave May   PetscFunctionReturn(0);
883b62e03f8SDave May }
884b62e03f8SDave May 
885d3a51819SDave May /*@C
886d3a51819SDave May    DMSwarmRestoreField - Restore access to the underlying array storing all entries associated with a registered field
887d3a51819SDave May 
888d3a51819SDave May    Not collective
889d3a51819SDave May 
890d3a51819SDave May    Input parameters:
89162741f57SDave May +  dm - a DMSwarm
89262741f57SDave May -  fieldname - the textual name to identify this field
893d3a51819SDave May 
894d3a51819SDave May    Output parameters:
89562741f57SDave May +  blocksize - the number of each data type
896d3a51819SDave May .  type - the data type
89762741f57SDave May -  data - pointer to raw array
898d3a51819SDave May 
899d3a51819SDave May    Level: beginner
900d3a51819SDave May 
901d3a51819SDave May    Notes:
9028b8a3813SDave May    The user must call DMSwarmGetField() prior to calling DMSwarmRestoreField().
903d3a51819SDave May 
904d3a51819SDave May .seealso: DMSwarmGetField()
905d3a51819SDave May @*/
90674d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmRestoreField(DM dm,const char fieldname[],PetscInt *blocksize,PetscDataType *type,void **data)
907b62e03f8SDave May {
908b62e03f8SDave May   DM_Swarm         *swarm = (DM_Swarm*)dm->data;
90977048351SPatrick Sanan   DMSwarmDataField gfield;
9102eac95f8SDave May   PetscErrorCode ierr;
911b62e03f8SDave May 
912521f74f9SMatthew G. Knepley   PetscFunctionBegin;
91377048351SPatrick Sanan   ierr = DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,fieldname,&gfield);CHKERRQ(ierr);
91477048351SPatrick Sanan   ierr = DMSwarmDataFieldRestoreAccess(gfield);CHKERRQ(ierr);
915b62e03f8SDave May   if (data) *data = NULL;
916b62e03f8SDave May   PetscFunctionReturn(0);
917b62e03f8SDave May }
918b62e03f8SDave May 
91974d0cae8SMatthew G. Knepley /*@
920d3a51819SDave May    DMSwarmAddPoint - Add space for one new point in the DMSwarm
921d3a51819SDave May 
922d3a51819SDave May    Not collective
923d3a51819SDave May 
924d3a51819SDave May    Input parameter:
925d3a51819SDave May .  dm - a DMSwarm
926d3a51819SDave May 
927d3a51819SDave May    Level: beginner
928d3a51819SDave May 
929d3a51819SDave May    Notes:
9308b8a3813SDave May    The new point will have all fields initialized to zero.
931d3a51819SDave May 
932d3a51819SDave May .seealso: DMSwarmAddNPoints()
933d3a51819SDave May @*/
93474d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmAddPoint(DM dm)
935cb1d1399SDave May {
936cb1d1399SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
937cb1d1399SDave May   PetscErrorCode ierr;
938cb1d1399SDave May 
939521f74f9SMatthew G. Knepley   PetscFunctionBegin;
9403454631fSDave May   if (!swarm->issetup) {ierr = DMSetUp(dm);CHKERRQ(ierr);}
941f2b2bee7SDave May   ierr = PetscLogEventBegin(DMSWARM_AddPoints,0,0,0,0);CHKERRQ(ierr);
94277048351SPatrick Sanan   ierr = DMSwarmDataBucketAddPoint(swarm->db);CHKERRQ(ierr);
943f2b2bee7SDave May   ierr = PetscLogEventEnd(DMSWARM_AddPoints,0,0,0,0);CHKERRQ(ierr);
944cb1d1399SDave May   PetscFunctionReturn(0);
945cb1d1399SDave May }
946cb1d1399SDave May 
94774d0cae8SMatthew G. Knepley /*@
948d3a51819SDave May    DMSwarmAddNPoints - Add space for a number of new points in the DMSwarm
949d3a51819SDave May 
950d3a51819SDave May    Not collective
951d3a51819SDave May 
952d3a51819SDave May    Input parameters:
95362741f57SDave May +  dm - a DMSwarm
95462741f57SDave May -  npoints - the number of new points to add
955d3a51819SDave May 
956d3a51819SDave May    Level: beginner
957d3a51819SDave May 
958d3a51819SDave May    Notes:
9598b8a3813SDave May    The new point will have all fields initialized to zero.
960d3a51819SDave May 
961d3a51819SDave May .seealso: DMSwarmAddPoint()
962d3a51819SDave May @*/
96374d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmAddNPoints(DM dm,PetscInt npoints)
964cb1d1399SDave May {
965cb1d1399SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
966cb1d1399SDave May   PetscErrorCode ierr;
967cb1d1399SDave May   PetscInt       nlocal;
968cb1d1399SDave May 
969521f74f9SMatthew G. Knepley   PetscFunctionBegin;
970f2b2bee7SDave May   ierr = PetscLogEventBegin(DMSWARM_AddPoints,0,0,0,0);CHKERRQ(ierr);
97177048351SPatrick Sanan   ierr = DMSwarmDataBucketGetSizes(swarm->db,&nlocal,NULL,NULL);CHKERRQ(ierr);
972cb1d1399SDave May   nlocal = nlocal + npoints;
97377048351SPatrick Sanan   ierr = DMSwarmDataBucketSetSizes(swarm->db,nlocal,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT);CHKERRQ(ierr);
974f2b2bee7SDave May   ierr = PetscLogEventEnd(DMSWARM_AddPoints,0,0,0,0);CHKERRQ(ierr);
975cb1d1399SDave May   PetscFunctionReturn(0);
976cb1d1399SDave May }
977cb1d1399SDave May 
97874d0cae8SMatthew G. Knepley /*@
979d3a51819SDave May    DMSwarmRemovePoint - Remove the last point from the DMSwarm
980d3a51819SDave May 
981d3a51819SDave May    Not collective
982d3a51819SDave May 
983d3a51819SDave May    Input parameter:
984d3a51819SDave May .  dm - a DMSwarm
985d3a51819SDave May 
986d3a51819SDave May    Level: beginner
987d3a51819SDave May 
988d3a51819SDave May .seealso: DMSwarmRemovePointAtIndex()
989d3a51819SDave May @*/
99074d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmRemovePoint(DM dm)
991cb1d1399SDave May {
992cb1d1399SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
993cb1d1399SDave May   PetscErrorCode ierr;
994cb1d1399SDave May 
995521f74f9SMatthew G. Knepley   PetscFunctionBegin;
996f2b2bee7SDave May   ierr = PetscLogEventBegin(DMSWARM_RemovePoints,0,0,0,0);CHKERRQ(ierr);
99777048351SPatrick Sanan   ierr = DMSwarmDataBucketRemovePoint(swarm->db);CHKERRQ(ierr);
998f2b2bee7SDave May   ierr = PetscLogEventEnd(DMSWARM_RemovePoints,0,0,0,0);CHKERRQ(ierr);
999cb1d1399SDave May   PetscFunctionReturn(0);
1000cb1d1399SDave May }
1001cb1d1399SDave May 
100274d0cae8SMatthew G. Knepley /*@
1003d3a51819SDave May    DMSwarmRemovePointAtIndex - Removes a specific point from the DMSwarm
1004d3a51819SDave May 
1005d3a51819SDave May    Not collective
1006d3a51819SDave May 
1007d3a51819SDave May    Input parameters:
100862741f57SDave May +  dm - a DMSwarm
100962741f57SDave May -  idx - index of point to remove
1010d3a51819SDave May 
1011d3a51819SDave May    Level: beginner
1012d3a51819SDave May 
1013d3a51819SDave May .seealso: DMSwarmRemovePoint()
1014d3a51819SDave May @*/
101574d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmRemovePointAtIndex(DM dm,PetscInt idx)
1016cb1d1399SDave May {
1017cb1d1399SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
1018cb1d1399SDave May   PetscErrorCode ierr;
1019cb1d1399SDave May 
1020521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1021f2b2bee7SDave May   ierr = PetscLogEventBegin(DMSWARM_RemovePoints,0,0,0,0);CHKERRQ(ierr);
102277048351SPatrick Sanan   ierr = DMSwarmDataBucketRemovePointAtIndex(swarm->db,idx);CHKERRQ(ierr);
1023f2b2bee7SDave May   ierr = PetscLogEventEnd(DMSWARM_RemovePoints,0,0,0,0);CHKERRQ(ierr);
1024cb1d1399SDave May   PetscFunctionReturn(0);
1025cb1d1399SDave May }
1026b62e03f8SDave May 
102774d0cae8SMatthew G. Knepley /*@
1028ba4fc9c6SDave May    DMSwarmCopyPoint - Copy point pj to point pi in the DMSwarm
1029ba4fc9c6SDave May 
1030ba4fc9c6SDave May    Not collective
1031ba4fc9c6SDave May 
1032ba4fc9c6SDave May    Input parameters:
1033ba4fc9c6SDave May +  dm - a DMSwarm
1034ba4fc9c6SDave May .  pi - the index of the point to copy
1035ba4fc9c6SDave May -  pj - the point index where the copy should be located
1036ba4fc9c6SDave May 
1037ba4fc9c6SDave May  Level: beginner
1038ba4fc9c6SDave May 
1039ba4fc9c6SDave May .seealso: DMSwarmRemovePoint()
1040ba4fc9c6SDave May @*/
104174d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmCopyPoint(DM dm,PetscInt pi,PetscInt pj)
1042ba4fc9c6SDave May {
1043ba4fc9c6SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
1044ba4fc9c6SDave May   PetscErrorCode ierr;
1045ba4fc9c6SDave May 
1046ba4fc9c6SDave May   PetscFunctionBegin;
1047ba4fc9c6SDave May   if (!swarm->issetup) {ierr = DMSetUp(dm);CHKERRQ(ierr);}
104877048351SPatrick Sanan   ierr = DMSwarmDataBucketCopyPoint(swarm->db,pi,swarm->db,pj);CHKERRQ(ierr);
1049ba4fc9c6SDave May   PetscFunctionReturn(0);
1050ba4fc9c6SDave May }
1051ba4fc9c6SDave May 
1052095059a4SDave May PetscErrorCode DMSwarmMigrate_Basic(DM dm,PetscBool remove_sent_points)
10533454631fSDave May {
1054dcf43ee8SDave May   PetscErrorCode ierr;
1055521f74f9SMatthew G. Knepley 
1056521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1057dcf43ee8SDave May   ierr = DMSwarmMigrate_Push_Basic(dm,remove_sent_points);CHKERRQ(ierr);
10583454631fSDave May   PetscFunctionReturn(0);
10593454631fSDave May }
10603454631fSDave May 
106174d0cae8SMatthew G. Knepley /*@
1062d3a51819SDave May    DMSwarmMigrate - Relocates points defined in the DMSwarm to other MPI-ranks
1063d3a51819SDave May 
1064d083f849SBarry Smith    Collective on dm
1065d3a51819SDave May 
1066d3a51819SDave May    Input parameters:
106762741f57SDave May +  dm - the DMSwarm
106862741f57SDave May -  remove_sent_points - flag indicating if sent points should be removed from the current MPI-rank
1069d3a51819SDave May 
1070d3a51819SDave May    Notes:
10718b8a3813SDave May    The DM will be modified to accomodate received points.
10728b8a3813SDave May    If remove_sent_points = PETSC_TRUE, any points that were sent will be removed from the DM.
10738b8a3813SDave May    Different styles of migration are supported. See DMSwarmSetMigrateType().
1074d3a51819SDave May 
1075d3a51819SDave May    Level: advanced
1076d3a51819SDave May 
1077d3a51819SDave May .seealso: DMSwarmSetMigrateType()
1078d3a51819SDave May @*/
107974d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmMigrate(DM dm,PetscBool remove_sent_points)
10803454631fSDave May {
1081f0cdbbbaSDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
10823454631fSDave May   PetscErrorCode ierr;
1083f0cdbbbaSDave May 
1084521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1085ed923d71SDave May   ierr = PetscLogEventBegin(DMSWARM_Migrate,0,0,0,0);CHKERRQ(ierr);
1086f0cdbbbaSDave May   switch (swarm->migrate_type) {
1087f0cdbbbaSDave May     case DMSWARM_MIGRATE_BASIC:
1088095059a4SDave May       ierr = DMSwarmMigrate_Basic(dm,remove_sent_points);CHKERRQ(ierr);
1089f0cdbbbaSDave May       break;
1090f0cdbbbaSDave May     case DMSWARM_MIGRATE_DMCELLNSCATTER:
1091f0cdbbbaSDave May       ierr = DMSwarmMigrate_CellDMScatter(dm,remove_sent_points);CHKERRQ(ierr);
1092f0cdbbbaSDave May       break;
1093f0cdbbbaSDave May     case DMSWARM_MIGRATE_DMCELLEXACT:
1094f0cdbbbaSDave May       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_MIGRATE_DMCELLEXACT not implemented");
1095f0cdbbbaSDave May       break;
1096f0cdbbbaSDave May     case DMSWARM_MIGRATE_USER:
1097f0cdbbbaSDave May       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_MIGRATE_USER not implemented");
1098f0cdbbbaSDave May       break;
1099f0cdbbbaSDave May     default:
1100f0cdbbbaSDave May       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_MIGRATE type unknown");
1101f0cdbbbaSDave May       break;
1102f0cdbbbaSDave May   }
1103ed923d71SDave May   ierr = PetscLogEventEnd(DMSWARM_Migrate,0,0,0,0);CHKERRQ(ierr);
11043454631fSDave May   PetscFunctionReturn(0);
11053454631fSDave May }
11063454631fSDave May 
1107f0cdbbbaSDave May PetscErrorCode DMSwarmMigrate_GlobalToLocal_Basic(DM dm,PetscInt *globalsize);
1108f0cdbbbaSDave May 
1109d3a51819SDave May /*
1110d3a51819SDave May  DMSwarmCollectViewCreate
1111d3a51819SDave May 
1112d3a51819SDave May  * Applies a collection method and gathers point neighbour points into dm
1113d3a51819SDave May 
1114d3a51819SDave May  Notes:
11158b8a3813SDave May  Users should call DMSwarmCollectViewDestroy() after
1116d3a51819SDave May  they have finished computations associated with the collected points
1117d3a51819SDave May */
1118d3a51819SDave May 
111974d0cae8SMatthew G. Knepley /*@
1120d3a51819SDave May    DMSwarmCollectViewCreate - Applies a collection method and gathers points
1121d3a51819SDave May    in neighbour MPI-ranks into the DMSwarm
1122d3a51819SDave May 
1123d083f849SBarry Smith    Collective on dm
1124d3a51819SDave May 
1125d3a51819SDave May    Input parameter:
1126d3a51819SDave May .  dm - the DMSwarm
1127d3a51819SDave May 
1128d3a51819SDave May    Notes:
1129d3a51819SDave May    Users should call DMSwarmCollectViewDestroy() after
1130d3a51819SDave May    they have finished computations associated with the collected points
11318b8a3813SDave May    Different collect methods are supported. See DMSwarmSetCollectType().
1132d3a51819SDave May 
1133d3a51819SDave May    Level: advanced
1134d3a51819SDave May 
1135d3a51819SDave May .seealso: DMSwarmCollectViewDestroy(), DMSwarmSetCollectType()
1136d3a51819SDave May @*/
113774d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmCollectViewCreate(DM dm)
11382712d1f2SDave May {
11392712d1f2SDave May   PetscErrorCode ierr;
11402712d1f2SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
11412712d1f2SDave May   PetscInt ng;
11422712d1f2SDave May 
1143521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1144480eef7bSDave May   if (swarm->collect_view_active) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"CollectView currently active");
1145480eef7bSDave May   ierr = DMSwarmGetLocalSize(dm,&ng);CHKERRQ(ierr);
1146480eef7bSDave May   switch (swarm->collect_type) {
1147f0cdbbbaSDave May 
1148480eef7bSDave May     case DMSWARM_COLLECT_BASIC:
11492712d1f2SDave May       ierr = DMSwarmMigrate_GlobalToLocal_Basic(dm,&ng);CHKERRQ(ierr);
1150480eef7bSDave May       break;
1151480eef7bSDave May     case DMSWARM_COLLECT_DMDABOUNDINGBOX:
1152f0cdbbbaSDave May       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_COLLECT_DMDABOUNDINGBOX not implemented");
1153480eef7bSDave May       break;
1154480eef7bSDave May     case DMSWARM_COLLECT_GENERAL:
1155f0cdbbbaSDave May       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_COLLECT_GENERAL not implemented");
1156480eef7bSDave May       break;
1157480eef7bSDave May     default:
1158f0cdbbbaSDave May       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_COLLECT type unknown");
1159480eef7bSDave May       break;
1160480eef7bSDave May   }
1161480eef7bSDave May   swarm->collect_view_active = PETSC_TRUE;
1162480eef7bSDave May   swarm->collect_view_reset_nlocal = ng;
11632712d1f2SDave May   PetscFunctionReturn(0);
11642712d1f2SDave May }
11652712d1f2SDave May 
116674d0cae8SMatthew G. Knepley /*@
1167d3a51819SDave May    DMSwarmCollectViewDestroy - Resets the DMSwarm to the size prior to calling DMSwarmCollectViewCreate()
1168d3a51819SDave May 
1169d083f849SBarry Smith    Collective on dm
1170d3a51819SDave May 
1171d3a51819SDave May    Input parameters:
1172d3a51819SDave May .  dm - the DMSwarm
1173d3a51819SDave May 
1174d3a51819SDave May    Notes:
1175d3a51819SDave May    Users should call DMSwarmCollectViewCreate() before this function is called.
1176d3a51819SDave May 
1177d3a51819SDave May    Level: advanced
1178d3a51819SDave May 
1179d3a51819SDave May .seealso: DMSwarmCollectViewCreate(), DMSwarmSetCollectType()
1180d3a51819SDave May @*/
118174d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmCollectViewDestroy(DM dm)
11822712d1f2SDave May {
11832712d1f2SDave May   PetscErrorCode ierr;
11842712d1f2SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
11852712d1f2SDave May 
1186521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1187480eef7bSDave May   if (!swarm->collect_view_active) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"CollectView is currently not active");
1188480eef7bSDave May   ierr = DMSwarmSetLocalSizes(dm,swarm->collect_view_reset_nlocal,-1);CHKERRQ(ierr);
1189480eef7bSDave May   swarm->collect_view_active = PETSC_FALSE;
11902712d1f2SDave May   PetscFunctionReturn(0);
11912712d1f2SDave May }
11923454631fSDave May 
1193f0cdbbbaSDave May PetscErrorCode DMSwarmSetUpPIC(DM dm)
1194f0cdbbbaSDave May {
1195f0cdbbbaSDave May   PetscInt       dim;
1196f0cdbbbaSDave May   PetscErrorCode ierr;
1197f0cdbbbaSDave May 
1198521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1199f0cdbbbaSDave May   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
1200f0cdbbbaSDave May   if (dim < 1) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Dimension must be 1,2,3 - found %D",dim);
1201f0cdbbbaSDave May   if (dim > 3) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Dimension must be 1,2,3 - found %D",dim);
1202f0cdbbbaSDave May   ierr = DMSwarmRegisterPetscDatatypeField(dm,DMSwarmPICField_coor,dim,PETSC_DOUBLE);CHKERRQ(ierr);
1203e2d107dbSDave May   ierr = DMSwarmRegisterPetscDatatypeField(dm,DMSwarmPICField_cellid,1,PETSC_INT);CHKERRQ(ierr);
1204f0cdbbbaSDave May   PetscFunctionReturn(0);
1205f0cdbbbaSDave May }
1206f0cdbbbaSDave May 
120774d0cae8SMatthew G. Knepley /*@
1208d3a51819SDave May    DMSwarmSetType - Set particular flavor of DMSwarm
1209d3a51819SDave May 
1210d083f849SBarry Smith    Collective on dm
1211d3a51819SDave May 
1212d3a51819SDave May    Input parameters:
121362741f57SDave May +  dm - the DMSwarm
121462741f57SDave May -  stype - the DMSwarm type (e.g. DMSWARM_PIC)
1215d3a51819SDave May 
1216d3a51819SDave May    Level: advanced
1217d3a51819SDave May 
1218d3a51819SDave May .seealso: DMSwarmSetMigrateType(), DMSwarmSetCollectType()
1219d3a51819SDave May @*/
122074d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmSetType(DM dm,DMSwarmType stype)
1221f0cdbbbaSDave May {
1222f0cdbbbaSDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
1223f0cdbbbaSDave May   PetscErrorCode ierr;
1224f0cdbbbaSDave May 
1225521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1226f0cdbbbaSDave May   swarm->swarm_type = stype;
1227f0cdbbbaSDave May   if (swarm->swarm_type == DMSWARM_PIC) {
1228f0cdbbbaSDave May     ierr = DMSwarmSetUpPIC(dm);CHKERRQ(ierr);
1229f0cdbbbaSDave May   }
1230f0cdbbbaSDave May   PetscFunctionReturn(0);
1231f0cdbbbaSDave May }
1232f0cdbbbaSDave May 
12333454631fSDave May PetscErrorCode DMSetup_Swarm(DM dm)
12343454631fSDave May {
12353454631fSDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
12363454631fSDave May   PetscErrorCode ierr;
12373454631fSDave May   PetscMPIInt    rank;
12383454631fSDave May   PetscInt       p,npoints,*rankval;
12393454631fSDave May 
1240521f74f9SMatthew G. Knepley   PetscFunctionBegin;
12413454631fSDave May   if (swarm->issetup) PetscFunctionReturn(0);
12423454631fSDave May 
12433454631fSDave May   swarm->issetup = PETSC_TRUE;
12443454631fSDave May 
1245f0cdbbbaSDave May   if (swarm->swarm_type == DMSWARM_PIC) {
1246f0cdbbbaSDave May     /* check dmcell exists */
1247f0cdbbbaSDave May     if (!swarm->dmcell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"DMSWARM_PIC requires you call DMSwarmSetCellDM");
1248f0cdbbbaSDave May 
1249f0cdbbbaSDave May     if (swarm->dmcell->ops->locatepointssubdomain) {
1250f0cdbbbaSDave May       /* check methods exists for exact ownership identificiation */
125177b6ec44SMatthew G. Knepley       ierr = PetscInfo(dm, "DMSWARM_PIC: Using method CellDM->ops->LocatePointsSubdomain\n");CHKERRQ(ierr);
1252f0cdbbbaSDave May       swarm->migrate_type = DMSWARM_MIGRATE_DMCELLEXACT;
1253f0cdbbbaSDave May     } else {
1254f0cdbbbaSDave May       /* check methods exist for point location AND rank neighbor identification */
1255f0cdbbbaSDave May       if (swarm->dmcell->ops->locatepoints) {
125677b6ec44SMatthew G. Knepley         ierr = PetscInfo(dm, "DMSWARM_PIC: Using method CellDM->LocatePoints\n");CHKERRQ(ierr);
1257f0cdbbbaSDave May       } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"DMSWARM_PIC requires the method CellDM->ops->locatepoints be defined");
1258f0cdbbbaSDave May 
1259f0cdbbbaSDave May       if (swarm->dmcell->ops->getneighbors) {
126077b6ec44SMatthew G. Knepley         ierr = PetscInfo(dm, "DMSWARM_PIC: Using method CellDM->GetNeigbors\n");CHKERRQ(ierr);
1261f0cdbbbaSDave May       } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"DMSWARM_PIC requires the method CellDM->ops->getneighbors be defined");
1262f0cdbbbaSDave May 
1263f0cdbbbaSDave May       swarm->migrate_type = DMSWARM_MIGRATE_DMCELLNSCATTER;
1264f0cdbbbaSDave May     }
1265f0cdbbbaSDave May   }
1266f0cdbbbaSDave May 
1267f0cdbbbaSDave May   ierr = DMSwarmFinalizeFieldRegister(dm);CHKERRQ(ierr);
1268f0cdbbbaSDave May 
12693454631fSDave May   /* check some fields were registered */
12703454631fSDave May   if (swarm->db->nfields <= 2) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"At least one field user must be registered via DMSwarmRegisterXXX()");
12713454631fSDave May 
12723454631fSDave May   /* check local sizes were set */
12733454631fSDave May   if (swarm->db->L == -1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Local sizes must be set via DMSwarmSetLocalSizes()");
12743454631fSDave May 
12753454631fSDave May   /* initialize values in pid and rank placeholders */
12763454631fSDave May   /* TODO: [pid - use MPI_Scan] */
12773454631fSDave May   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
127877048351SPatrick Sanan   ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
1279f0cdbbbaSDave May   ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
12803454631fSDave May   for (p=0; p<npoints; p++) {
12813454631fSDave May     rankval[p] = (PetscInt)rank;
12823454631fSDave May   }
1283f0cdbbbaSDave May   ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
12843454631fSDave May   PetscFunctionReturn(0);
12853454631fSDave May }
12863454631fSDave May 
1287dc5f5ce9SDave May extern PetscErrorCode DMSwarmSortDestroy(DMSwarmSort *_ctx);
1288dc5f5ce9SDave May 
128957795646SDave May PetscErrorCode DMDestroy_Swarm(DM dm)
129057795646SDave May {
129157795646SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
129257795646SDave May   PetscErrorCode ierr;
129357795646SDave May 
129457795646SDave May   PetscFunctionBegin;
129577048351SPatrick Sanan   ierr = DMSwarmDataBucketDestroy(&swarm->db);CHKERRQ(ierr);
1296dc5f5ce9SDave May   if (swarm->sort_context) {
1297dc5f5ce9SDave May     ierr = DMSwarmSortDestroy(&swarm->sort_context);CHKERRQ(ierr);
1298dc5f5ce9SDave May   }
129957795646SDave May   ierr = PetscFree(swarm);CHKERRQ(ierr);
130057795646SDave May   PetscFunctionReturn(0);
130157795646SDave May }
130257795646SDave May 
1303a9ee3421SMatthew G. Knepley PetscErrorCode DMSwarmView_Draw(DM dm, PetscViewer viewer)
1304a9ee3421SMatthew G. Knepley {
1305a9ee3421SMatthew G. Knepley   DM             cdm;
1306a9ee3421SMatthew G. Knepley   PetscDraw      draw;
1307a9ee3421SMatthew G. Knepley   PetscReal     *coords, oldPause;
1308a9ee3421SMatthew G. Knepley   PetscInt       Np, p, bs;
1309a9ee3421SMatthew G. Knepley   PetscErrorCode ierr;
1310a9ee3421SMatthew G. Knepley 
1311a9ee3421SMatthew G. Knepley   PetscFunctionBegin;
1312a9ee3421SMatthew G. Knepley   ierr = PetscViewerDrawGetDraw(viewer, 0, &draw);CHKERRQ(ierr);
1313a9ee3421SMatthew G. Knepley   ierr = DMSwarmGetCellDM(dm, &cdm);CHKERRQ(ierr);
1314a9ee3421SMatthew G. Knepley   ierr = PetscDrawGetPause(draw, &oldPause);CHKERRQ(ierr);
1315a9ee3421SMatthew G. Knepley   ierr = PetscDrawSetPause(draw, 0.0);CHKERRQ(ierr);
1316a9ee3421SMatthew G. Knepley   ierr = DMView(cdm, viewer);CHKERRQ(ierr);
1317a9ee3421SMatthew G. Knepley   ierr = PetscDrawSetPause(draw, oldPause);CHKERRQ(ierr);
1318a9ee3421SMatthew G. Knepley 
1319a9ee3421SMatthew G. Knepley   ierr = DMSwarmGetLocalSize(dm, &Np);CHKERRQ(ierr);
1320a9ee3421SMatthew G. Knepley   ierr = DMSwarmGetField(dm, DMSwarmPICField_coor, &bs, NULL, (void **) &coords);CHKERRQ(ierr);
1321a9ee3421SMatthew G. Knepley   for (p = 0; p < Np; ++p) {
1322a9ee3421SMatthew G. Knepley     const PetscInt i = p*bs;
1323a9ee3421SMatthew G. Knepley 
1324a9ee3421SMatthew G. Knepley     ierr = PetscDrawEllipse(draw, coords[i], coords[i+1], 0.01, 0.01, PETSC_DRAW_BLUE);CHKERRQ(ierr);
1325a9ee3421SMatthew G. Knepley   }
1326a9ee3421SMatthew G. Knepley   ierr = DMSwarmRestoreField(dm, DMSwarmPICField_coor, &bs, NULL, (void **) &coords);CHKERRQ(ierr);
1327a9ee3421SMatthew G. Knepley   ierr = PetscDrawFlush(draw);CHKERRQ(ierr);
1328a9ee3421SMatthew G. Knepley   ierr = PetscDrawPause(draw);CHKERRQ(ierr);
1329a9ee3421SMatthew G. Knepley   ierr = PetscDrawSave(draw);CHKERRQ(ierr);
1330a9ee3421SMatthew G. Knepley   PetscFunctionReturn(0);
1331a9ee3421SMatthew G. Knepley }
1332a9ee3421SMatthew G. Knepley 
13335f50eb2eSDave May PetscErrorCode DMView_Swarm(DM dm, PetscViewer viewer)
13345f50eb2eSDave May {
13355f50eb2eSDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
1336a9ee3421SMatthew G. Knepley   PetscBool      iascii,ibinary,ishdf5,isvtk,isdraw;
13375f50eb2eSDave May   PetscErrorCode ierr;
13385f50eb2eSDave May 
13395f50eb2eSDave May   PetscFunctionBegin;
13405f50eb2eSDave May   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
13415f50eb2eSDave May   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
13425f50eb2eSDave May   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
13435f50eb2eSDave May   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY,&ibinary);CHKERRQ(ierr);
13445f50eb2eSDave May   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERVTK,   &isvtk);CHKERRQ(ierr);
13455f50eb2eSDave May   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5,  &ishdf5);CHKERRQ(ierr);
1346a9ee3421SMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW,  &isdraw);CHKERRQ(ierr);
13475f50eb2eSDave May   if (iascii) {
134877048351SPatrick Sanan     ierr = DMSwarmDataBucketView(PetscObjectComm((PetscObject)dm),swarm->db,NULL,DATABUCKET_VIEW_STDOUT);CHKERRQ(ierr);
1349298827fbSBarry Smith   } else if (ibinary) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"NO Binary support");
13505f50eb2eSDave May #if defined(PETSC_HAVE_HDF5)
135174d0cae8SMatthew G. Knepley   else if (ishdf5) {
135274d0cae8SMatthew G. Knepley     ierr = DMSwarmView_HDF5(dm, viewer);CHKERRQ(ierr);
135374d0cae8SMatthew G. Knepley   }
13545f50eb2eSDave May #else
1355298827fbSBarry Smith   else if (ishdf5) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"HDF5 not supported. Please reconfigure using --download-hdf5");
13565f50eb2eSDave May #endif
1357298827fbSBarry Smith   else if (isvtk) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"NO VTK support");
1358298827fbSBarry Smith   else if (isdraw) {
1359a9ee3421SMatthew G. Knepley     ierr = DMSwarmView_Draw(dm, viewer);CHKERRQ(ierr);
13605f50eb2eSDave May   }
13615f50eb2eSDave May   PetscFunctionReturn(0);
13625f50eb2eSDave May }
13635f50eb2eSDave May 
1364d3a51819SDave May /*MC
1365d3a51819SDave May 
1366d3a51819SDave May  DMSWARM = "swarm" - A DM object used to represent arrays of data (fields) of arbitrary data type.
136762741f57SDave May  This implementation was designed for particle methods in which the underlying
1368d3a51819SDave May  data required to be represented is both (i) dynamic in length, (ii) and of arbitrary data type.
1369d3a51819SDave May 
137062741f57SDave May  User data can be represented by DMSwarm through a registering "fields".
137162741f57SDave May  To register a field, the user must provide;
137262741f57SDave May  (a) a unique name;
137362741f57SDave May  (b) the data type (or size in bytes);
137462741f57SDave May  (c) the block size of the data.
1375d3a51819SDave May 
1376d3a51819SDave May  For example, suppose the application requires a unique id, energy, momentum and density to be stored
1377c215a47eSMatthew Knepley  on a set of particles. Then the following code could be used
1378d3a51819SDave May 
137962741f57SDave May $    DMSwarmInitializeFieldRegister(dm)
138062741f57SDave May $    DMSwarmRegisterPetscDatatypeField(dm,"uid",1,PETSC_LONG);
138162741f57SDave May $    DMSwarmRegisterPetscDatatypeField(dm,"energy",1,PETSC_REAL);
138262741f57SDave May $    DMSwarmRegisterPetscDatatypeField(dm,"momentum",3,PETSC_REAL);
138362741f57SDave May $    DMSwarmRegisterPetscDatatypeField(dm,"density",1,PETSC_FLOAT);
138462741f57SDave May $    DMSwarmFinalizeFieldRegister(dm)
1385d3a51819SDave May 
1386d3a51819SDave May  The fields represented by DMSwarm are dynamic and can be re-sized at any time.
138762741f57SDave May  The only restriction imposed by DMSwarm is that all fields contain the same number of points.
1388d3a51819SDave May 
1389d3a51819SDave May  To support particle methods, "migration" techniques are provided. These methods migrate data
1390d3a51819SDave May  between MPI-ranks.
1391d3a51819SDave May 
1392d3a51819SDave May  DMSwarm supports the methods DMCreateGlobalVector() and DMCreateLocalVector().
1393d3a51819SDave May  As a DMSwarm may internally define and store values of different data types,
139462741f57SDave May  before calling DMCreateGlobalVector() or DMCreateLocalVector(), the user must inform DMSwarm which
1395d3a51819SDave May  fields should be used to define a Vec object via
1396d3a51819SDave May    DMSwarmVectorDefineField()
1397c215a47eSMatthew Knepley  The specified field can be changed at any time - thereby permitting vectors
1398c215a47eSMatthew Knepley  compatible with different fields to be created.
1399d3a51819SDave May 
140062741f57SDave May  A dual representation of fields in the DMSwarm and a Vec object is permitted via
1401d3a51819SDave May    DMSwarmCreateGlobalVectorFromField()
1402d3a51819SDave May  Here the data defining the field in the DMSwarm is shared with a Vec.
1403d3a51819SDave May  This is inherently unsafe if you alter the size of the field at any time between
1404d3a51819SDave May  calls to DMSwarmCreateGlobalVectorFromField() and DMSwarmDestroyGlobalVectorFromField().
1405cc651181SDave May  If the local size of the DMSwarm does not match the local size of the global vector
1406cc651181SDave May  when DMSwarmDestroyGlobalVectorFromField() is called, an error is thrown.
1407d3a51819SDave May 
140862741f57SDave May  Additional high-level support is provided for Particle-In-Cell methods.
140962741f57SDave May  Please refer to the man page for DMSwarmSetType().
141062741f57SDave May 
1411d3a51819SDave May  Level: beginner
1412d3a51819SDave May 
1413d3a51819SDave May .seealso: DMType, DMCreate(), DMSetType()
1414d3a51819SDave May M*/
141557795646SDave May PETSC_EXTERN PetscErrorCode DMCreate_Swarm(DM dm)
141657795646SDave May {
141757795646SDave May   DM_Swarm      *swarm;
141857795646SDave May   PetscErrorCode ierr;
141957795646SDave May 
142057795646SDave May   PetscFunctionBegin;
142157795646SDave May   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
142257795646SDave May   ierr     = PetscNewLog(dm,&swarm);CHKERRQ(ierr);
1423f0cdbbbaSDave May   dm->data = swarm;
142457795646SDave May 
142577048351SPatrick Sanan   ierr = DMSwarmDataBucketCreate(&swarm->db);CHKERRQ(ierr);
1426f0cdbbbaSDave May   ierr = DMSwarmInitializeFieldRegister(dm);CHKERRQ(ierr);
1427f0cdbbbaSDave May 
1428b5bcf523SDave May   swarm->vec_field_set = PETSC_FALSE;
14293454631fSDave May   swarm->issetup = PETSC_FALSE;
1430480eef7bSDave May   swarm->swarm_type = DMSWARM_BASIC;
1431480eef7bSDave May   swarm->migrate_type = DMSWARM_MIGRATE_BASIC;
1432480eef7bSDave May   swarm->collect_type = DMSWARM_COLLECT_BASIC;
143340c453e9SDave May   swarm->migrate_error_on_missing_point = PETSC_FALSE;
1434b62e03f8SDave May 
1435f0cdbbbaSDave May   swarm->dmcell = NULL;
1436f0cdbbbaSDave May   swarm->collect_view_active = PETSC_FALSE;
1437f0cdbbbaSDave May   swarm->collect_view_reset_nlocal = -1;
143857795646SDave May 
1439f0cdbbbaSDave May   dm->dim  = 0;
14405f50eb2eSDave May   dm->ops->view                            = DMView_Swarm;
144157795646SDave May   dm->ops->load                            = NULL;
144257795646SDave May   dm->ops->setfromoptions                  = NULL;
144357795646SDave May   dm->ops->clone                           = NULL;
14443454631fSDave May   dm->ops->setup                           = DMSetup_Swarm;
14451bb6d2a8SBarry Smith   dm->ops->createlocalsection              = NULL;
144657795646SDave May   dm->ops->createdefaultconstraints        = NULL;
1447b5bcf523SDave May   dm->ops->createglobalvector              = DMCreateGlobalVector_Swarm;
1448b5bcf523SDave May   dm->ops->createlocalvector               = DMCreateLocalVector_Swarm;
144957795646SDave May   dm->ops->getlocaltoglobalmapping         = NULL;
145057795646SDave May   dm->ops->createfieldis                   = NULL;
145157795646SDave May   dm->ops->createcoordinatedm              = NULL;
145257795646SDave May   dm->ops->getcoloring                     = NULL;
145357795646SDave May   dm->ops->creatematrix                    = NULL;
145457795646SDave May   dm->ops->createinterpolation             = NULL;
14555a84ad33SLisandro Dalcin   dm->ops->createinjection                 = NULL;
145683c47955SMatthew G. Knepley   dm->ops->createmassmatrix                = DMCreateMassMatrix_Swarm;
145757795646SDave May   dm->ops->refine                          = NULL;
145857795646SDave May   dm->ops->coarsen                         = NULL;
145957795646SDave May   dm->ops->refinehierarchy                 = NULL;
146057795646SDave May   dm->ops->coarsenhierarchy                = NULL;
146157795646SDave May   dm->ops->globaltolocalbegin              = NULL;
146257795646SDave May   dm->ops->globaltolocalend                = NULL;
146357795646SDave May   dm->ops->localtoglobalbegin              = NULL;
146457795646SDave May   dm->ops->localtoglobalend                = NULL;
146557795646SDave May   dm->ops->destroy                         = DMDestroy_Swarm;
146657795646SDave May   dm->ops->createsubdm                     = NULL;
146757795646SDave May   dm->ops->getdimpoints                    = NULL;
146857795646SDave May   dm->ops->locatepoints                    = NULL;
146957795646SDave May   PetscFunctionReturn(0);
147057795646SDave May }
1471