xref: /petsc/src/dm/impls/swarm/swarm.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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;
38*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetDM(v, &dm));
39*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetBlockSize(v, &bs));
40*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerHDF5PushGroup(viewer, "/particle_fields"));
41*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject) v, VECSEQ, &isseq));
42*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetOutputSequenceNumber(dm, &seqnum, &seqval));
43*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerHDF5SetTimestep(viewer, seqnum));
44*5f80ce2aSJacob Faibussowitsch   /* CHKERRQ(DMSequenceView_HDF5(dm, "time", seqnum, (PetscScalar) seqval, viewer)); */
45*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecViewNative(v, viewer));
46*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject) v, "Nc", PETSC_INT, (void *) &bs));
47*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
58*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetSize(dm, &Np));
59*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmCreateGlobalVectorFromField(dm, DMSwarmPICField_coor, &coordinates));
60*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) coordinates, "coordinates"));
61*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerHDF5PushGroup(viewer, "/particles"));
62*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject) coordinates, VECSEQ, &isseq));
63*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecViewNative(coordinates, viewer));
64*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject) coordinates, "Np", PETSC_INT, (void *) &Np));
65*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerHDF5PopGroup(viewer));
66*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
79*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetDM(v, &dm));
802c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!dm,PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_WRONG, "Vector not generated from a DM");
81f9558f5cSBarry Smith #if defined(PETSC_HAVE_HDF5)
82*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERHDF5,  &ishdf5));
8374d0cae8SMatthew G. Knepley   if (ishdf5) {
84*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecView_Swarm_HDF5_Internal(v, viewer));
85f9558f5cSBarry Smith       PetscFunctionReturn(0);
8674d0cae8SMatthew G. Knepley   }
87f9558f5cSBarry Smith #endif
88*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
121*5f80ce2aSJacob Faibussowitsch   if (!swarm->issetup) CHKERRQ(DMSetUp(dm));
122*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketGetSizes(swarm->db,&n,NULL,NULL));
123*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetField(dm,fieldname,&bs,&type,(void**)&array));
124b5bcf523SDave May 
125b5bcf523SDave May   /* Check all fields are of type PETSC_REAL or PETSC_SCALAR */
1262c71b3e2SJacob Faibussowitsch   PetscCheckFalse(type != PETSC_REAL,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only valid for PETSC_REAL");
127*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
131*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
143*5f80ce2aSJacob Faibussowitsch   if (!swarm->issetup) CHKERRQ(DMSetUp(dm));
1442c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!swarm->vec_field_set,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMSwarmVectorDefineField first");
1452c71b3e2SJacob Faibussowitsch   PetscCheckFalse(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 
147*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"DMSwarmField_%s",swarm->vec_field_name));
148*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PetscObjectComm((PetscObject)dm),&x));
149*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)x,name));
150*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(x,swarm->db->L*swarm->vec_field_bs,PETSC_DETERMINE));
151*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetBlockSize(x,swarm->vec_field_bs));
152*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetDM(x,dm));
153*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(x));
154*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetDM(x, dm));
155*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
168*5f80ce2aSJacob Faibussowitsch   if (!swarm->issetup) CHKERRQ(DMSetUp(dm));
1692c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!swarm->vec_field_set,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMSwarmVectorDefineField first");
1702c71b3e2SJacob Faibussowitsch   PetscCheckFalse(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 
172*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"DMSwarmField_%s",swarm->vec_field_name));
173*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_SELF,&x));
174*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)x,name));
175*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(x,swarm->db->L*swarm->vec_field_bs,PETSC_DETERMINE));
176*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetBlockSize(x,swarm->vec_field_bs));
177*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetDM(x,dm));
178*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
192*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetLocalSize(*vec, &nlocal));
193*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetBlockSize(*vec, &bs));
1942c71b3e2SJacob Faibussowitsch   PetscCheckFalse(nlocal/bs != swarm->db->L,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"DMSwarm sizes have changed since vector was created - cannot ensure pointers are valid");
195*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield));
196fb1bcc12SMatthew G. Knepley   /* check vector is an inplace array */
197*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSNPrintf(name, PETSC_MAX_PATH_LEN-1, "DMSwarm_VecFieldInPlace_%s", fieldname));
198*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectQueryFunction((PetscObject) *vec, name, &fptr));
1992c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!fptr,PetscObjectComm((PetscObject) dm), PETSC_ERR_USER, "Vector being destroyed was not created from DMSwarm field(%s)", fieldname);
200*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataFieldRestoreAccess(gfield));
201*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(vec));
202fb1bcc12SMatthew G. Knepley   PetscFunctionReturn(0);
203fb1bcc12SMatthew G. Knepley }
204fb1bcc12SMatthew G. Knepley 
205fb1bcc12SMatthew G. Knepley static PetscErrorCode DMSwarmCreateVectorFromField_Private(DM dm, const char fieldname[], MPI_Comm comm, Vec *vec)
206fb1bcc12SMatthew G. Knepley {
207fb1bcc12SMatthew G. Knepley   DM_Swarm      *swarm = (DM_Swarm *) dm->data;
208fb1bcc12SMatthew G. Knepley   PetscDataType  type;
209fb1bcc12SMatthew G. Knepley   PetscScalar   *array;
210fb1bcc12SMatthew G. Knepley   PetscInt       bs, n;
211fb1bcc12SMatthew G. Knepley   char           name[PETSC_MAX_PATH_LEN];
212e4fbd051SBarry Smith   PetscMPIInt    size;
213fb1bcc12SMatthew G. Knepley 
214fb1bcc12SMatthew G. Knepley   PetscFunctionBegin;
215*5f80ce2aSJacob Faibussowitsch   if (!swarm->issetup) CHKERRQ(DMSetUp(dm));
216*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketGetSizes(swarm->db, &n, NULL, NULL));
217*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetField(dm, fieldname, &bs, &type, (void **) &array));
2182c71b3e2SJacob Faibussowitsch   PetscCheckFalse(type != PETSC_REAL,PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Only valid for PETSC_REAL");
219fb1bcc12SMatthew G. Knepley 
220*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(comm, &size));
221e4fbd051SBarry Smith   if (size == 1) {
222*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCreateSeqWithArray(comm, bs, n*bs, array, vec));
223fb1bcc12SMatthew G. Knepley   } else {
224*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCreateMPIWithArray(comm, bs, n*bs, PETSC_DETERMINE, array, vec));
225fb1bcc12SMatthew G. Knepley   }
226*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSNPrintf(name, PETSC_MAX_PATH_LEN-1, "DMSwarmSharedField_%s", fieldname));
227*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) *vec, name));
228fb1bcc12SMatthew G. Knepley 
229fb1bcc12SMatthew G. Knepley   /* Set guard */
230*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSNPrintf(name, PETSC_MAX_PATH_LEN-1, "DMSwarm_VecFieldInPlace_%s", fieldname));
231*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject) *vec, name, DMSwarmDestroyVectorFromField_Private));
23274d0cae8SMatthew G. Knepley 
233*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetDM(*vec, dm));
234*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Swarm));
235fb1bcc12SMatthew G. Knepley   PetscFunctionReturn(0);
236fb1bcc12SMatthew G. Knepley }
237fb1bcc12SMatthew G. Knepley 
2380643ed39SMatthew G. Knepley /* This creates a "mass matrix" between a finite element and particle space. If a finite element interpolant is given by
2390643ed39SMatthew G. Knepley 
2400643ed39SMatthew G. Knepley      \hat f = \sum_i f_i \phi_i
2410643ed39SMatthew G. Knepley 
2420643ed39SMatthew G. Knepley    and a particle function is given by
2430643ed39SMatthew G. Knepley 
2440643ed39SMatthew G. Knepley      f = \sum_i w_i \delta(x - x_i)
2450643ed39SMatthew G. Knepley 
2460643ed39SMatthew G. Knepley    then we want to require that
2470643ed39SMatthew G. Knepley 
2480643ed39SMatthew G. Knepley      M \hat f = M_p f
2490643ed39SMatthew G. Knepley 
2500643ed39SMatthew G. Knepley    where the particle mass matrix is given by
2510643ed39SMatthew G. Knepley 
2520643ed39SMatthew G. Knepley      (M_p)_{ij} = \int \phi_i \delta(x - x_j)
2530643ed39SMatthew G. Knepley 
2540643ed39SMatthew G. Knepley    The way Dave May does particles, they amount to quadratue weights rather than delta functions, so he has |J| is in
2550643ed39SMatthew G. Knepley    his integral. We allow this with the boolean flag.
2560643ed39SMatthew G. Knepley */
2570643ed39SMatthew G. Knepley static PetscErrorCode DMSwarmComputeMassMatrix_Private(DM dmc, DM dmf, Mat mass, PetscBool useDeltaFunction, void *ctx)
25883c47955SMatthew G. Knepley {
25983c47955SMatthew G. Knepley   const char    *name = "Mass Matrix";
2600643ed39SMatthew G. Knepley   MPI_Comm       comm;
26183c47955SMatthew G. Knepley   PetscDS        prob;
26283c47955SMatthew G. Knepley   PetscSection   fsection, globalFSection;
263e8f14785SLisandro Dalcin   PetscHSetIJ    ht;
2640643ed39SMatthew G. Knepley   PetscLayout    rLayout, colLayout;
26583c47955SMatthew G. Knepley   PetscInt      *dnz, *onz;
266adb2528bSMark Adams   PetscInt       locRows, locCols, rStart, colStart, colEnd, *rowIDXs;
2670643ed39SMatthew G. Knepley   PetscReal     *xi, *v0, *J, *invJ, detJ = 1.0, v0ref[3] = {-1.0, -1.0, -1.0};
26883c47955SMatthew G. Knepley   PetscScalar   *elemMat;
2690643ed39SMatthew G. Knepley   PetscInt       dim, Nf, field, cStart, cEnd, cell, totDim, maxC = 0;
27083c47955SMatthew G. Knepley 
27183c47955SMatthew G. Knepley   PetscFunctionBegin;
272*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetComm((PetscObject) mass, &comm));
273*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateDim(dmf, &dim));
274*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDS(dmf, &prob));
275*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSGetNumFields(prob, &Nf));
276*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSGetTotalDimension(prob, &totDim));
277*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc3(dim, &v0, dim*dim, &J, dim*dim,&invJ));
278*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalSection(dmf, &fsection));
279*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetGlobalSection(dmf, &globalFSection));
280*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd));
281*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetLocalSize(mass, &locRows, &locCols));
28283c47955SMatthew G. Knepley 
283*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutCreate(comm, &colLayout));
284*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutSetLocalSize(colLayout, locCols));
285*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutSetBlockSize(colLayout, 1));
286*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutSetUp(colLayout));
287*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutGetRange(colLayout, &colStart, &colEnd));
288*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutDestroy(&colLayout));
2890643ed39SMatthew G. Knepley 
290*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutCreate(comm, &rLayout));
291*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutSetLocalSize(rLayout, locRows));
292*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutSetBlockSize(rLayout, 1));
293*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutSetUp(rLayout));
294*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutGetRange(rLayout, &rStart, NULL));
295*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutDestroy(&rLayout));
2960643ed39SMatthew G. Knepley 
297*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCalloc2(locRows, &dnz, locRows, &onz));
298*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIJCreate(&ht));
29953e60ab4SJoseph Pusztay 
300*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSynchronizedFlush(comm, NULL));
3010643ed39SMatthew G. Knepley   /* count non-zeros */
302*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSortGetAccess(dmc));
30383c47955SMatthew G. Knepley   for (field = 0; field < Nf; ++field) {
30483c47955SMatthew G. Knepley     for (cell = cStart; cell < cEnd; ++cell) {
3050643ed39SMatthew G. Knepley       PetscInt  c, i;
3060643ed39SMatthew G. Knepley       PetscInt *findices,   *cindices; /* fine is vertices, coarse is particles */
30783c47955SMatthew G. Knepley       PetscInt  numFIndices, numCIndices;
30883c47955SMatthew G. Knepley 
309*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
310*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices));
311fc7c92abSMatthew G. Knepley       maxC = PetscMax(maxC, numCIndices);
31283c47955SMatthew G. Knepley       {
313e8f14785SLisandro Dalcin         PetscHashIJKey key;
314e8f14785SLisandro Dalcin         PetscBool      missing;
31583c47955SMatthew G. Knepley         for (i = 0; i < numFIndices; ++i) {
316adb2528bSMark Adams           key.j = findices[i]; /* global column (from Plex) */
317adb2528bSMark Adams           if (key.j >= 0) {
31883c47955SMatthew G. Knepley             /* Get indices for coarse elements */
31983c47955SMatthew G. Knepley             for (c = 0; c < numCIndices; ++c) {
320adb2528bSMark Adams               key.i = cindices[c] + rStart; /* global cols (from Swarm) */
321adb2528bSMark Adams               if (key.i < 0) continue;
322*5f80ce2aSJacob Faibussowitsch               CHKERRQ(PetscHSetIJQueryAdd(ht, key, &missing));
32383c47955SMatthew G. Knepley               if (missing) {
3240643ed39SMatthew G. Knepley                 if ((key.j >= colStart) && (key.j < colEnd)) ++dnz[key.i - rStart];
325e8f14785SLisandro Dalcin                 else                                         ++onz[key.i - rStart];
32698921bdaSJacob Faibussowitsch               } else SETERRQ(PetscObjectComm((PetscObject) dmf), PETSC_ERR_SUP, "Set new value at %D,%D", key.i, key.j);
32783c47955SMatthew G. Knepley             }
328fc7c92abSMatthew G. Knepley           }
329fc7c92abSMatthew G. Knepley         }
330*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscFree(cindices));
33183c47955SMatthew G. Knepley       }
332*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
33383c47955SMatthew G. Knepley     }
33483c47955SMatthew G. Knepley   }
335*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIJDestroy(&ht));
336*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatXAIJSetPreallocation(mass, 1, dnz, onz, NULL, NULL));
337*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOption(mass, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
338*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree2(dnz, onz));
339*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc3(maxC*totDim, &elemMat, maxC, &rowIDXs, maxC*dim, &xi));
34083c47955SMatthew G. Knepley   for (field = 0; field < Nf; ++field) {
341ef0bb6c7SMatthew G. Knepley     PetscTabulation Tcoarse;
34283c47955SMatthew G. Knepley     PetscObject     obj;
343ef0bb6c7SMatthew G. Knepley     PetscReal       *coords;
3440643ed39SMatthew G. Knepley     PetscInt        Nc, i;
34583c47955SMatthew G. Knepley 
346*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSGetDiscretization(prob, field, &obj));
347*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFEGetNumComponents((PetscFE) obj, &Nc));
3482c71b3e2SJacob Faibussowitsch     PetscCheckFalse(Nc != 1,PetscObjectComm((PetscObject) dmf), PETSC_ERR_SUP, "Can only interpolate a scalar field from particles, Nc = %D", Nc);
349*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmGetField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **) &coords));
35083c47955SMatthew G. Knepley     for (cell = cStart; cell < cEnd; ++cell) {
35183c47955SMatthew G. Knepley       PetscInt *findices  , *cindices;
35283c47955SMatthew G. Knepley       PetscInt  numFIndices, numCIndices;
3530643ed39SMatthew G. Knepley       PetscInt  p, c;
35483c47955SMatthew G. Knepley 
3550643ed39SMatthew G. Knepley       /* TODO: Use DMField instead of assuming affine */
356*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ));
357*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
358*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices));
3590643ed39SMatthew G. Knepley       for (p = 0; p < numCIndices; ++p) {
3600643ed39SMatthew G. Knepley         CoordinatesRealToRef(dim, dim, v0ref, v0, invJ, &coords[cindices[p]*dim], &xi[p*dim]);
3610643ed39SMatthew G. Knepley       }
362*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFECreateTabulation((PetscFE) obj, 1, numCIndices, xi, 0, &Tcoarse));
36383c47955SMatthew G. Knepley       /* Get elemMat entries by multiplying by weight */
364*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscArrayzero(elemMat, numCIndices*totDim));
36583c47955SMatthew G. Knepley       for (i = 0; i < numFIndices; ++i) {
3660643ed39SMatthew G. Knepley         for (p = 0; p < numCIndices; ++p) {
3670643ed39SMatthew G. Knepley           for (c = 0; c < Nc; ++c) {
3680643ed39SMatthew G. Knepley             /* B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c */
369ef0bb6c7SMatthew G. Knepley             elemMat[p*numFIndices+i] += Tcoarse->T[0][(p*numFIndices + i)*Nc + c]*(useDeltaFunction ? 1.0 : detJ);
37083c47955SMatthew G. Knepley           }
3710643ed39SMatthew G. Knepley         }
3720643ed39SMatthew G. Knepley       }
373adb2528bSMark Adams       for (p = 0; p < numCIndices; ++p) rowIDXs[p] = cindices[p] + rStart;
374*5f80ce2aSJacob Faibussowitsch       if (0) CHKERRQ(DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat));
375*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValues(mass, numCIndices, rowIDXs, numFIndices, findices, elemMat, ADD_VALUES));
376*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFree(cindices));
377*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
378*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscTabulationDestroy(&Tcoarse));
37983c47955SMatthew G. Knepley     }
380*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmRestoreField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **) &coords));
38183c47955SMatthew G. Knepley   }
382*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree3(elemMat, rowIDXs, xi));
383*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSortRestoreAccess(dmc));
384*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree3(v0, J, invJ));
385*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(mass, MAT_FINAL_ASSEMBLY));
386*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(mass, MAT_FINAL_ASSEMBLY));
38783c47955SMatthew G. Knepley   PetscFunctionReturn(0);
38883c47955SMatthew G. Knepley }
38983c47955SMatthew G. Knepley 
390d0c080abSJoseph Pusztay /* Returns empty matrix for use with SNES FD */
39170a7d78aSStefano Zampini static PetscErrorCode DMCreateMatrix_Swarm(DM sw, Mat* m)
39270a7d78aSStefano Zampini {
393d0c080abSJoseph Pusztay   Vec            field;
394d0c080abSJoseph Pusztay   PetscInt       size;
395d0c080abSJoseph Pusztay 
396d0c080abSJoseph Pusztay   PetscFunctionBegin;
397*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetGlobalVector(sw, &field));
398*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetLocalSize(field, &size));
399*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreGlobalVector(sw, &field));
400*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD, m));
401*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(*m));
402*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(*m, PETSC_DECIDE, PETSC_DECIDE, size, size));
403*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqAIJSetPreallocation(*m, 1, NULL));
404*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatZeroEntries(*m));
405*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(*m, MAT_FINAL_ASSEMBLY));
406*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(*m, MAT_FINAL_ASSEMBLY));
407*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShift(*m, 1.0));
408*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetDM(*m, sw));
409d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
410d0c080abSJoseph Pusztay }
411d0c080abSJoseph Pusztay 
412adb2528bSMark Adams /* FEM cols, Particle rows */
41383c47955SMatthew G. Knepley static PetscErrorCode DMCreateMassMatrix_Swarm(DM dmCoarse, DM dmFine, Mat *mass)
41483c47955SMatthew G. Knepley {
415895a1698SMatthew G. Knepley   PetscSection   gsf;
41683c47955SMatthew G. Knepley   PetscInt       m, n;
41783c47955SMatthew G. Knepley   void          *ctx;
41883c47955SMatthew G. Knepley 
41983c47955SMatthew G. Knepley   PetscFunctionBegin;
420*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetGlobalSection(dmFine, &gsf));
421*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionGetConstrainedStorageSize(gsf, &m));
422*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetLocalSize(dmCoarse, &n));
423*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PetscObjectComm((PetscObject) dmCoarse), mass));
424*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(*mass, n, m, PETSC_DETERMINE, PETSC_DETERMINE));
425*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(*mass, dmCoarse->mattype));
426*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetApplicationContext(dmFine, &ctx));
42783c47955SMatthew G. Knepley 
428*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmComputeMassMatrix_Private(dmCoarse, dmFine, *mass, PETSC_TRUE, ctx));
429*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatViewFromOptions(*mass, NULL, "-mass_mat_view"));
43083c47955SMatthew G. Knepley   PetscFunctionReturn(0);
43183c47955SMatthew G. Knepley }
43283c47955SMatthew G. Knepley 
4334cc7f7b2SMatthew G. Knepley static PetscErrorCode DMSwarmComputeMassMatrixSquare_Private(DM dmc, DM dmf, Mat mass, PetscBool useDeltaFunction, void *ctx)
4344cc7f7b2SMatthew G. Knepley {
4354cc7f7b2SMatthew G. Knepley   const char    *name = "Mass Matrix Square";
4364cc7f7b2SMatthew G. Knepley   MPI_Comm       comm;
4374cc7f7b2SMatthew G. Knepley   PetscDS        prob;
4384cc7f7b2SMatthew G. Knepley   PetscSection   fsection, globalFSection;
4394cc7f7b2SMatthew G. Knepley   PetscHSetIJ    ht;
4404cc7f7b2SMatthew G. Knepley   PetscLayout    rLayout, colLayout;
4414cc7f7b2SMatthew G. Knepley   PetscInt      *dnz, *onz, *adj, depth, maxConeSize, maxSupportSize, maxAdjSize;
4424cc7f7b2SMatthew G. Knepley   PetscInt       locRows, locCols, rStart, colStart, colEnd, *rowIDXs;
4434cc7f7b2SMatthew G. Knepley   PetscReal     *xi, *v0, *J, *invJ, detJ = 1.0, v0ref[3] = {-1.0, -1.0, -1.0};
4444cc7f7b2SMatthew G. Knepley   PetscScalar   *elemMat, *elemMatSq;
4454cc7f7b2SMatthew G. Knepley   PetscInt       cdim, Nf, field, cStart, cEnd, cell, totDim, maxC = 0;
4464cc7f7b2SMatthew G. Knepley 
4474cc7f7b2SMatthew G. Knepley   PetscFunctionBegin;
448*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetComm((PetscObject) mass, &comm));
449*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateDim(dmf, &cdim));
450*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDS(dmf, &prob));
451*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSGetNumFields(prob, &Nf));
452*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSGetTotalDimension(prob, &totDim));
453*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc3(cdim, &v0, cdim*cdim, &J, cdim*cdim,&invJ));
454*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalSection(dmf, &fsection));
455*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetGlobalSection(dmf, &globalFSection));
456*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd));
457*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetLocalSize(mass, &locRows, &locCols));
4584cc7f7b2SMatthew G. Knepley 
459*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutCreate(comm, &colLayout));
460*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutSetLocalSize(colLayout, locCols));
461*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutSetBlockSize(colLayout, 1));
462*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutSetUp(colLayout));
463*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutGetRange(colLayout, &colStart, &colEnd));
464*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutDestroy(&colLayout));
4654cc7f7b2SMatthew G. Knepley 
466*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutCreate(comm, &rLayout));
467*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutSetLocalSize(rLayout, locRows));
468*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutSetBlockSize(rLayout, 1));
469*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutSetUp(rLayout));
470*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutGetRange(rLayout, &rStart, NULL));
471*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutDestroy(&rLayout));
4724cc7f7b2SMatthew G. Knepley 
473*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetDepth(dmf, &depth));
474*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetMaxSizes(dmf, &maxConeSize, &maxSupportSize));
4754cc7f7b2SMatthew G. Knepley   maxAdjSize = PetscPowInt(maxConeSize*maxSupportSize, depth);
476*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(maxAdjSize, &adj));
4774cc7f7b2SMatthew G. Knepley 
478*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCalloc2(locRows, &dnz, locRows, &onz));
479*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIJCreate(&ht));
4804cc7f7b2SMatthew G. Knepley   /* Count nonzeros
4814cc7f7b2SMatthew G. Knepley        This is just FVM++, but we cannot use the Plex P0 allocation since unknowns in a cell will not be contiguous
4824cc7f7b2SMatthew G. Knepley   */
483*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSortGetAccess(dmc));
4844cc7f7b2SMatthew G. Knepley   for (cell = cStart; cell < cEnd; ++cell) {
4854cc7f7b2SMatthew G. Knepley     PetscInt  i;
4864cc7f7b2SMatthew G. Knepley     PetscInt *cindices;
4874cc7f7b2SMatthew G. Knepley     PetscInt  numCIndices;
4884cc7f7b2SMatthew G. Knepley   #if 0
4894cc7f7b2SMatthew G. Knepley     PetscInt  adjSize = maxAdjSize, a, j;
4904cc7f7b2SMatthew G. Knepley   #endif
4914cc7f7b2SMatthew G. Knepley 
492*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices));
4934cc7f7b2SMatthew G. Knepley     maxC = PetscMax(maxC, numCIndices);
4944cc7f7b2SMatthew G. Knepley     /* Diagonal block */
4954cc7f7b2SMatthew G. Knepley     for (i = 0; i < numCIndices; ++i) {dnz[cindices[i]] += numCIndices;}
4964cc7f7b2SMatthew G. Knepley #if 0
4974cc7f7b2SMatthew G. Knepley     /* Off-diagonal blocks */
498*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetAdjacency(dmf, cell, &adjSize, &adj));
4994cc7f7b2SMatthew G. Knepley     for (a = 0; a < adjSize; ++a) {
5004cc7f7b2SMatthew G. Knepley       if (adj[a] >= cStart && adj[a] < cEnd && adj[a] != cell) {
5014cc7f7b2SMatthew G. Knepley         const PetscInt ncell = adj[a];
5024cc7f7b2SMatthew G. Knepley         PetscInt      *ncindices;
5034cc7f7b2SMatthew G. Knepley         PetscInt       numNCIndices;
5044cc7f7b2SMatthew G. Knepley 
505*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMSwarmSortGetPointsPerCell(dmc, ncell, &numNCIndices, &ncindices));
5064cc7f7b2SMatthew G. Knepley         {
5074cc7f7b2SMatthew G. Knepley           PetscHashIJKey key;
5084cc7f7b2SMatthew G. Knepley           PetscBool      missing;
5094cc7f7b2SMatthew G. Knepley 
5104cc7f7b2SMatthew G. Knepley           for (i = 0; i < numCIndices; ++i) {
5114cc7f7b2SMatthew G. Knepley             key.i = cindices[i] + rStart; /* global rows (from Swarm) */
5124cc7f7b2SMatthew G. Knepley             if (key.i < 0) continue;
5134cc7f7b2SMatthew G. Knepley             for (j = 0; j < numNCIndices; ++j) {
5144cc7f7b2SMatthew G. Knepley               key.j = ncindices[j] + rStart; /* global column (from Swarm) */
5154cc7f7b2SMatthew G. Knepley               if (key.j < 0) continue;
516*5f80ce2aSJacob Faibussowitsch               CHKERRQ(PetscHSetIJQueryAdd(ht, key, &missing));
5174cc7f7b2SMatthew G. Knepley               if (missing) {
5184cc7f7b2SMatthew G. Knepley                 if ((key.j >= colStart) && (key.j < colEnd)) ++dnz[key.i - rStart];
5194cc7f7b2SMatthew G. Knepley                 else                                         ++onz[key.i - rStart];
5204cc7f7b2SMatthew G. Knepley               }
5214cc7f7b2SMatthew G. Knepley             }
5224cc7f7b2SMatthew G. Knepley           }
5234cc7f7b2SMatthew G. Knepley         }
524*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscFree(ncindices));
5254cc7f7b2SMatthew G. Knepley       }
5264cc7f7b2SMatthew G. Knepley     }
5274cc7f7b2SMatthew G. Knepley #endif
528*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(cindices));
5294cc7f7b2SMatthew G. Knepley   }
530*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIJDestroy(&ht));
531*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatXAIJSetPreallocation(mass, 1, dnz, onz, NULL, NULL));
532*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOption(mass, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
533*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree2(dnz, onz));
534*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc4(maxC*totDim, &elemMat, maxC*maxC, &elemMatSq, maxC, &rowIDXs, maxC*cdim, &xi));
5354cc7f7b2SMatthew G. Knepley   /* Fill in values
5364cc7f7b2SMatthew G. Knepley        Each entry is a sum of terms \phi_i(x_p) \phi_i(x_q)
5374cc7f7b2SMatthew G. Knepley        Start just by producing block diagonal
5384cc7f7b2SMatthew G. Knepley        Could loop over adjacent cells
5394cc7f7b2SMatthew G. Knepley          Produce neighboring element matrix
5404cc7f7b2SMatthew G. Knepley          TODO Determine which columns and rows correspond to shared dual vector
5414cc7f7b2SMatthew G. Knepley          Do MatMatMult with rectangular matrices
5424cc7f7b2SMatthew G. Knepley          Insert block
5434cc7f7b2SMatthew G. Knepley   */
5444cc7f7b2SMatthew G. Knepley   for (field = 0; field < Nf; ++field) {
5454cc7f7b2SMatthew G. Knepley     PetscTabulation Tcoarse;
5464cc7f7b2SMatthew G. Knepley     PetscObject     obj;
5474cc7f7b2SMatthew G. Knepley     PetscReal       *coords;
5484cc7f7b2SMatthew G. Knepley     PetscInt        Nc, i;
5494cc7f7b2SMatthew G. Knepley 
550*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSGetDiscretization(prob, field, &obj));
551*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFEGetNumComponents((PetscFE) obj, &Nc));
5522c71b3e2SJacob Faibussowitsch     PetscCheckFalse(Nc != 1,PetscObjectComm((PetscObject) dmf), PETSC_ERR_SUP, "Can only interpolate a scalar field from particles, Nc = %D", Nc);
553*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmGetField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **) &coords));
5544cc7f7b2SMatthew G. Knepley     for (cell = cStart; cell < cEnd; ++cell) {
5554cc7f7b2SMatthew G. Knepley       PetscInt *findices  , *cindices;
5564cc7f7b2SMatthew G. Knepley       PetscInt  numFIndices, numCIndices;
5574cc7f7b2SMatthew G. Knepley       PetscInt  p, c;
5584cc7f7b2SMatthew G. Knepley 
5594cc7f7b2SMatthew G. Knepley       /* TODO: Use DMField instead of assuming affine */
560*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ));
561*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
562*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices));
5634cc7f7b2SMatthew G. Knepley       for (p = 0; p < numCIndices; ++p) {
5644cc7f7b2SMatthew G. Knepley         CoordinatesRealToRef(cdim, cdim, v0ref, v0, invJ, &coords[cindices[p]*cdim], &xi[p*cdim]);
5654cc7f7b2SMatthew G. Knepley       }
566*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFECreateTabulation((PetscFE) obj, 1, numCIndices, xi, 0, &Tcoarse));
5674cc7f7b2SMatthew G. Knepley       /* Get elemMat entries by multiplying by weight */
568*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscArrayzero(elemMat, numCIndices*totDim));
5694cc7f7b2SMatthew G. Knepley       for (i = 0; i < numFIndices; ++i) {
5704cc7f7b2SMatthew G. Knepley         for (p = 0; p < numCIndices; ++p) {
5714cc7f7b2SMatthew G. Knepley           for (c = 0; c < Nc; ++c) {
5724cc7f7b2SMatthew G. Knepley             /* B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c */
5734cc7f7b2SMatthew G. Knepley             elemMat[p*numFIndices+i] += Tcoarse->T[0][(p*numFIndices + i)*Nc + c]*(useDeltaFunction ? 1.0 : detJ);
5744cc7f7b2SMatthew G. Knepley           }
5754cc7f7b2SMatthew G. Knepley         }
5764cc7f7b2SMatthew G. Knepley       }
577*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscTabulationDestroy(&Tcoarse));
5784cc7f7b2SMatthew G. Knepley       for (p = 0; p < numCIndices; ++p) rowIDXs[p] = cindices[p] + rStart;
579*5f80ce2aSJacob Faibussowitsch       if (0) CHKERRQ(DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat));
5804cc7f7b2SMatthew G. Knepley       /* Block diagonal */
58178884ca7SMatthew G. Knepley       if (numCIndices) {
5824cc7f7b2SMatthew G. Knepley         PetscBLASInt blasn, blask;
5834cc7f7b2SMatthew G. Knepley         PetscScalar  one = 1.0, zero = 0.0;
5844cc7f7b2SMatthew G. Knepley 
585*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscBLASIntCast(numCIndices, &blasn));
586*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscBLASIntCast(numFIndices, &blask));
5874cc7f7b2SMatthew G. Knepley         PetscStackCallBLAS("BLASgemm",BLASgemm_("T","N",&blasn,&blasn,&blask,&one,elemMat,&blask,elemMat,&blask,&zero,elemMatSq,&blasn));
5884cc7f7b2SMatthew G. Knepley       }
589*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValues(mass, numCIndices, rowIDXs, numCIndices, rowIDXs, elemMatSq, ADD_VALUES));
5904cc7f7b2SMatthew G. Knepley       /* TODO Off-diagonal */
591*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFree(cindices));
592*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
5934cc7f7b2SMatthew G. Knepley     }
594*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmRestoreField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **) &coords));
5954cc7f7b2SMatthew G. Knepley   }
596*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree4(elemMat, elemMatSq, rowIDXs, xi));
597*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(adj));
598*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSortRestoreAccess(dmc));
599*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree3(v0, J, invJ));
600*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(mass, MAT_FINAL_ASSEMBLY));
601*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(mass, MAT_FINAL_ASSEMBLY));
6024cc7f7b2SMatthew G. Knepley   PetscFunctionReturn(0);
6034cc7f7b2SMatthew G. Knepley }
6044cc7f7b2SMatthew G. Knepley 
6054cc7f7b2SMatthew G. Knepley /*@C
6064cc7f7b2SMatthew G. Knepley   DMSwarmCreateMassMatrixSquare - Creates the block-diagonal of the square, M^T_p M_p, of the particle mass matrix M_p
6074cc7f7b2SMatthew G. Knepley 
6084cc7f7b2SMatthew G. Knepley   Collective on dmCoarse
6094cc7f7b2SMatthew G. Knepley 
6104cc7f7b2SMatthew G. Knepley   Input parameters:
6114cc7f7b2SMatthew G. Knepley + dmCoarse - a DMSwarm
6124cc7f7b2SMatthew G. Knepley - dmFine   - a DMPlex
6134cc7f7b2SMatthew G. Knepley 
6144cc7f7b2SMatthew G. Knepley   Output parameter:
6154cc7f7b2SMatthew G. Knepley . mass     - the square of the particle mass matrix
6164cc7f7b2SMatthew G. Knepley 
6174cc7f7b2SMatthew G. Knepley   Level: advanced
6184cc7f7b2SMatthew G. Knepley 
6194cc7f7b2SMatthew G. Knepley   Notes:
6204cc7f7b2SMatthew G. Knepley   We only compute the block diagonal since this provides a good preconditioner and is completely local. It would be possible in the
6214cc7f7b2SMatthew G. Knepley   future to compute the full normal equations.
6224cc7f7b2SMatthew G. Knepley 
6234cc7f7b2SMatthew G. Knepley .seealso: DMCreateMassMatrix()
6244cc7f7b2SMatthew G. Knepley @*/
6254cc7f7b2SMatthew G. Knepley PetscErrorCode DMSwarmCreateMassMatrixSquare(DM dmCoarse, DM dmFine, Mat *mass)
6264cc7f7b2SMatthew G. Knepley {
6274cc7f7b2SMatthew G. Knepley   PetscInt       n;
6284cc7f7b2SMatthew G. Knepley   void          *ctx;
6294cc7f7b2SMatthew G. Knepley 
6304cc7f7b2SMatthew G. Knepley   PetscFunctionBegin;
631*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetLocalSize(dmCoarse, &n));
632*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PetscObjectComm((PetscObject) dmCoarse), mass));
633*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(*mass, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
634*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(*mass, dmCoarse->mattype));
635*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetApplicationContext(dmFine, &ctx));
6364cc7f7b2SMatthew G. Knepley 
637*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmComputeMassMatrixSquare_Private(dmCoarse, dmFine, *mass, PETSC_TRUE, ctx));
638*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatViewFromOptions(*mass, NULL, "-mass_sq_mat_view"));
6394cc7f7b2SMatthew G. Knepley   PetscFunctionReturn(0);
6404cc7f7b2SMatthew G. Knepley }
6414cc7f7b2SMatthew G. Knepley 
642fb1bcc12SMatthew G. Knepley /*@C
643d3a51819SDave May    DMSwarmCreateGlobalVectorFromField - Creates a Vec object sharing the array associated with a given field
644d3a51819SDave May 
645d083f849SBarry Smith    Collective on dm
646d3a51819SDave May 
647d3a51819SDave May    Input parameters:
64862741f57SDave May +  dm - a DMSwarm
64962741f57SDave May -  fieldname - the textual name given to a registered field
650d3a51819SDave May 
6518b8a3813SDave May    Output parameter:
652d3a51819SDave May .  vec - the vector
653d3a51819SDave May 
654d3a51819SDave May    Level: beginner
655d3a51819SDave May 
6568b8a3813SDave May    Notes:
6578b8a3813SDave May    The vector must be returned using a matching call to DMSwarmDestroyGlobalVectorFromField().
6588b8a3813SDave May 
6598b8a3813SDave May .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmDestroyGlobalVectorFromField()
660d3a51819SDave May @*/
66174d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmCreateGlobalVectorFromField(DM dm,const char fieldname[],Vec *vec)
662b5bcf523SDave May {
663fb1bcc12SMatthew G. Knepley   MPI_Comm       comm = PetscObjectComm((PetscObject) dm);
664b5bcf523SDave May 
665fb1bcc12SMatthew G. Knepley   PetscFunctionBegin;
666*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmCreateVectorFromField_Private(dm, fieldname, comm, vec));
667b5bcf523SDave May   PetscFunctionReturn(0);
668b5bcf523SDave May }
669b5bcf523SDave May 
670d3a51819SDave May /*@C
671d3a51819SDave May    DMSwarmDestroyGlobalVectorFromField - Destroys the Vec object which share the array associated with a given field
672d3a51819SDave May 
673d083f849SBarry Smith    Collective on dm
674d3a51819SDave May 
675d3a51819SDave May    Input parameters:
67662741f57SDave May +  dm - a DMSwarm
67762741f57SDave May -  fieldname - the textual name given to a registered field
678d3a51819SDave May 
6798b8a3813SDave May    Output parameter:
680d3a51819SDave May .  vec - the vector
681d3a51819SDave May 
682d3a51819SDave May    Level: beginner
683d3a51819SDave May 
6848b8a3813SDave May .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmCreateGlobalVectorFromField()
685d3a51819SDave May @*/
68674d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmDestroyGlobalVectorFromField(DM dm,const char fieldname[],Vec *vec)
687b5bcf523SDave May {
688fb1bcc12SMatthew G. Knepley   PetscFunctionBegin;
689*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDestroyVectorFromField_Private(dm, fieldname, vec));
690b5bcf523SDave May   PetscFunctionReturn(0);
691b5bcf523SDave May }
692b5bcf523SDave May 
693fb1bcc12SMatthew G. Knepley /*@C
694fb1bcc12SMatthew G. Knepley    DMSwarmCreateLocalVectorFromField - Creates a Vec object sharing the array associated with a given field
695fb1bcc12SMatthew G. Knepley 
696d083f849SBarry Smith    Collective on dm
697fb1bcc12SMatthew G. Knepley 
698fb1bcc12SMatthew G. Knepley    Input parameters:
69962741f57SDave May +  dm - a DMSwarm
70062741f57SDave May -  fieldname - the textual name given to a registered field
701fb1bcc12SMatthew G. Knepley 
7028b8a3813SDave May    Output parameter:
703fb1bcc12SMatthew G. Knepley .  vec - the vector
704fb1bcc12SMatthew G. Knepley 
705fb1bcc12SMatthew G. Knepley    Level: beginner
706fb1bcc12SMatthew G. Knepley 
7078b8a3813SDave May    Notes:
7088b8a3813SDave May    The vector must be returned using a matching call to DMSwarmDestroyLocalVectorFromField().
7098b8a3813SDave May 
7108b8a3813SDave May .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmDestroyLocalVectorFromField()
711fb1bcc12SMatthew G. Knepley @*/
71274d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmCreateLocalVectorFromField(DM dm,const char fieldname[],Vec *vec)
713bbe8250bSMatthew G. Knepley {
714fb1bcc12SMatthew G. Knepley   MPI_Comm       comm = PETSC_COMM_SELF;
715bbe8250bSMatthew G. Knepley 
716fb1bcc12SMatthew G. Knepley   PetscFunctionBegin;
717*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmCreateVectorFromField_Private(dm, fieldname, comm, vec));
718fb1bcc12SMatthew G. Knepley   PetscFunctionReturn(0);
719bbe8250bSMatthew G. Knepley }
720fb1bcc12SMatthew G. Knepley 
721fb1bcc12SMatthew G. Knepley /*@C
722fb1bcc12SMatthew G. Knepley    DMSwarmDestroyLocalVectorFromField - Destroys the Vec object which share the array associated with a given field
723fb1bcc12SMatthew G. Knepley 
724d083f849SBarry Smith    Collective on dm
725fb1bcc12SMatthew G. Knepley 
726fb1bcc12SMatthew G. Knepley    Input parameters:
72762741f57SDave May +  dm - a DMSwarm
72862741f57SDave May -  fieldname - the textual name given to a registered field
729fb1bcc12SMatthew G. Knepley 
7308b8a3813SDave May    Output parameter:
731fb1bcc12SMatthew G. Knepley .  vec - the vector
732fb1bcc12SMatthew G. Knepley 
733fb1bcc12SMatthew G. Knepley    Level: beginner
734fb1bcc12SMatthew G. Knepley 
7358b8a3813SDave May .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmCreateLocalVectorFromField()
736fb1bcc12SMatthew G. Knepley @*/
73774d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmDestroyLocalVectorFromField(DM dm,const char fieldname[],Vec *vec)
738fb1bcc12SMatthew G. Knepley {
739fb1bcc12SMatthew G. Knepley   PetscFunctionBegin;
740*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDestroyVectorFromField_Private(dm, fieldname, vec));
741bbe8250bSMatthew G. Knepley   PetscFunctionReturn(0);
742bbe8250bSMatthew G. Knepley }
743bbe8250bSMatthew G. Knepley 
744d3a51819SDave May /*@C
745d3a51819SDave May    DMSwarmInitializeFieldRegister - Initiates the registration of fields to a DMSwarm
746d3a51819SDave May 
747d083f849SBarry Smith    Collective on dm
748d3a51819SDave May 
749d3a51819SDave May    Input parameter:
750d3a51819SDave May .  dm - a DMSwarm
751d3a51819SDave May 
752d3a51819SDave May    Level: beginner
753d3a51819SDave May 
754d3a51819SDave May    Notes:
7558b8a3813SDave May    After all fields have been registered, you must call DMSwarmFinalizeFieldRegister().
756d3a51819SDave May 
757d3a51819SDave May .seealso: DMSwarmFinalizeFieldRegister(), DMSwarmRegisterPetscDatatypeField(),
758d3a51819SDave May           DMSwarmRegisterUserStructField(), DMSwarmRegisterUserDatatypeField()
759d3a51819SDave May @*/
76074d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmInitializeFieldRegister(DM dm)
7615f50eb2eSDave May {
7625f50eb2eSDave May   DM_Swarm      *swarm = (DM_Swarm *) dm->data;
7633454631fSDave May 
764521f74f9SMatthew G. Knepley   PetscFunctionBegin;
765cc651181SDave May   if (!swarm->field_registration_initialized) {
7665f50eb2eSDave May     swarm->field_registration_initialized = PETSC_TRUE;
767*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmRegisterPetscDatatypeField(dm,DMSwarmField_pid,1,PETSC_INT64)); /* unique identifer */
768*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmRegisterPetscDatatypeField(dm,DMSwarmField_rank,1,PETSC_INT)); /* used for communication */
769cc651181SDave May   }
7705f50eb2eSDave May   PetscFunctionReturn(0);
7715f50eb2eSDave May }
7725f50eb2eSDave May 
77374d0cae8SMatthew G. Knepley /*@
774d3a51819SDave May    DMSwarmFinalizeFieldRegister - Finalizes the registration of fields to a DMSwarm
775d3a51819SDave May 
776d083f849SBarry Smith    Collective on dm
777d3a51819SDave May 
778d3a51819SDave May    Input parameter:
779d3a51819SDave May .  dm - a DMSwarm
780d3a51819SDave May 
781d3a51819SDave May    Level: beginner
782d3a51819SDave May 
783d3a51819SDave May    Notes:
78462741f57SDave May    After DMSwarmFinalizeFieldRegister() has been called, no new fields can be defined on the DMSwarm.
785d3a51819SDave May 
786d3a51819SDave May .seealso: DMSwarmInitializeFieldRegister(), DMSwarmRegisterPetscDatatypeField(),
787d3a51819SDave May           DMSwarmRegisterUserStructField(), DMSwarmRegisterUserDatatypeField()
788d3a51819SDave May @*/
78974d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmFinalizeFieldRegister(DM dm)
7905f50eb2eSDave May {
7915f50eb2eSDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
7926845f8f5SDave May 
793521f74f9SMatthew G. Knepley   PetscFunctionBegin;
794f0cdbbbaSDave May   if (!swarm->field_registration_finalized) {
795*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataBucketFinalize(swarm->db));
796f0cdbbbaSDave May   }
797f0cdbbbaSDave May   swarm->field_registration_finalized = PETSC_TRUE;
7985f50eb2eSDave May   PetscFunctionReturn(0);
7995f50eb2eSDave May }
8005f50eb2eSDave May 
80174d0cae8SMatthew G. Knepley /*@
802d3a51819SDave May    DMSwarmSetLocalSizes - Sets the length of all registered fields on the DMSwarm
803d3a51819SDave May 
804d3a51819SDave May    Not collective
805d3a51819SDave May 
806d3a51819SDave May    Input parameters:
80762741f57SDave May +  dm - a DMSwarm
808d3a51819SDave May .  nlocal - the length of each registered field
80962741f57SDave May -  buffer - the length of the buffer used to efficient dynamic re-sizing
810d3a51819SDave May 
811d3a51819SDave May    Level: beginner
812d3a51819SDave May 
813d3a51819SDave May .seealso: DMSwarmGetLocalSize()
814d3a51819SDave May @*/
81574d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmSetLocalSizes(DM dm,PetscInt nlocal,PetscInt buffer)
8165f50eb2eSDave May {
8175f50eb2eSDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
8185f50eb2eSDave May 
819521f74f9SMatthew G. Knepley   PetscFunctionBegin;
820*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventBegin(DMSWARM_SetSizes,0,0,0,0));
821*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketSetSizes(swarm->db,nlocal,buffer));
822*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventEnd(DMSWARM_SetSizes,0,0,0,0));
8235f50eb2eSDave May   PetscFunctionReturn(0);
8245f50eb2eSDave May }
8255f50eb2eSDave May 
82674d0cae8SMatthew G. Knepley /*@
827d3a51819SDave May    DMSwarmSetCellDM - Attachs a DM to a DMSwarm
828d3a51819SDave May 
829d083f849SBarry Smith    Collective on dm
830d3a51819SDave May 
831d3a51819SDave May    Input parameters:
83262741f57SDave May +  dm - a DMSwarm
83362741f57SDave May -  dmcell - the DM to attach to the DMSwarm
834d3a51819SDave May 
835d3a51819SDave May    Level: beginner
836d3a51819SDave May 
837d3a51819SDave May    Notes:
838d3a51819SDave May    The attached DM (dmcell) will be queried for point location and
8398b8a3813SDave May    neighbor MPI-rank information if DMSwarmMigrate() is called.
840d3a51819SDave May 
8418b8a3813SDave May .seealso: DMSwarmSetType(), DMSwarmGetCellDM(), DMSwarmMigrate()
842d3a51819SDave May @*/
84374d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmSetCellDM(DM dm,DM dmcell)
844b16650c8SDave May {
845b16650c8SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
846521f74f9SMatthew G. Knepley 
847521f74f9SMatthew G. Knepley   PetscFunctionBegin;
848b16650c8SDave May   swarm->dmcell = dmcell;
849b16650c8SDave May   PetscFunctionReturn(0);
850b16650c8SDave May }
851b16650c8SDave May 
85274d0cae8SMatthew G. Knepley /*@
853d3a51819SDave May    DMSwarmGetCellDM - Fetches the attached cell DM
854d3a51819SDave May 
855d083f849SBarry Smith    Collective on dm
856d3a51819SDave May 
857d3a51819SDave May    Input parameter:
858d3a51819SDave May .  dm - a DMSwarm
859d3a51819SDave May 
860d3a51819SDave May    Output parameter:
861d3a51819SDave May .  dmcell - the DM which was attached to the DMSwarm
862d3a51819SDave May 
863d3a51819SDave May    Level: beginner
864d3a51819SDave May 
865d3a51819SDave May .seealso: DMSwarmSetCellDM()
866d3a51819SDave May @*/
86774d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmGetCellDM(DM dm,DM *dmcell)
868fe39f135SDave May {
869fe39f135SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
870521f74f9SMatthew G. Knepley 
871521f74f9SMatthew G. Knepley   PetscFunctionBegin;
872fe39f135SDave May   *dmcell = swarm->dmcell;
873fe39f135SDave May   PetscFunctionReturn(0);
874fe39f135SDave May }
875fe39f135SDave May 
87674d0cae8SMatthew G. Knepley /*@
877d3a51819SDave May    DMSwarmGetLocalSize - Retrives the local length of fields registered
878d3a51819SDave May 
879d3a51819SDave May    Not collective
880d3a51819SDave May 
881d3a51819SDave May    Input parameter:
882d3a51819SDave May .  dm - a DMSwarm
883d3a51819SDave May 
884d3a51819SDave May    Output parameter:
885d3a51819SDave May .  nlocal - the length of each registered field
886d3a51819SDave May 
887d3a51819SDave May    Level: beginner
888d3a51819SDave May 
8898b8a3813SDave May .seealso: DMSwarmGetSize(), DMSwarmSetLocalSizes()
890d3a51819SDave May @*/
89174d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmGetLocalSize(DM dm,PetscInt *nlocal)
892dcf43ee8SDave May {
893dcf43ee8SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
894dcf43ee8SDave May 
895521f74f9SMatthew G. Knepley   PetscFunctionBegin;
896*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketGetSizes(swarm->db,nlocal,NULL,NULL));
897dcf43ee8SDave May   PetscFunctionReturn(0);
898dcf43ee8SDave May }
899dcf43ee8SDave May 
90074d0cae8SMatthew G. Knepley /*@
901d3a51819SDave May    DMSwarmGetSize - Retrives the total length of fields registered
902d3a51819SDave May 
903d083f849SBarry Smith    Collective on dm
904d3a51819SDave May 
905d3a51819SDave May    Input parameter:
906d3a51819SDave May .  dm - a DMSwarm
907d3a51819SDave May 
908d3a51819SDave May    Output parameter:
909d3a51819SDave May .  n - the total length of each registered field
910d3a51819SDave May 
911d3a51819SDave May    Level: beginner
912d3a51819SDave May 
913d3a51819SDave May    Note:
914d3a51819SDave May    This calls MPI_Allreduce upon each call (inefficient but safe)
915d3a51819SDave May 
9168b8a3813SDave May .seealso: DMSwarmGetLocalSize(), DMSwarmSetLocalSizes()
917d3a51819SDave May @*/
91874d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmGetSize(DM dm,PetscInt *n)
919dcf43ee8SDave May {
920dcf43ee8SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
9215627991aSBarry Smith   PetscInt       nlocal;
922dcf43ee8SDave May 
923521f74f9SMatthew G. Knepley   PetscFunctionBegin;
924*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketGetSizes(swarm->db,&nlocal,NULL,NULL));
925*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Allreduce(&nlocal,n,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm)));
926dcf43ee8SDave May   PetscFunctionReturn(0);
927dcf43ee8SDave May }
928dcf43ee8SDave May 
929d3a51819SDave May /*@C
9308b8a3813SDave May    DMSwarmRegisterPetscDatatypeField - Register a field to a DMSwarm with a native PETSc data type
931d3a51819SDave May 
932d083f849SBarry Smith    Collective on dm
933d3a51819SDave May 
934d3a51819SDave May    Input parameters:
93562741f57SDave May +  dm - a DMSwarm
936d3a51819SDave May .  fieldname - the textual name to identify this field
937d3a51819SDave May .  blocksize - the number of each data type
93862741f57SDave May -  type - a valid PETSc data type (PETSC_CHAR, PETSC_SHORT, PETSC_INT, PETSC_FLOAT, PETSC_REAL, PETSC_LONG)
939d3a51819SDave May 
940d3a51819SDave May    Level: beginner
941d3a51819SDave May 
942d3a51819SDave May    Notes:
9438b8a3813SDave May    The textual name for each registered field must be unique.
944d3a51819SDave May 
945d3a51819SDave May .seealso: DMSwarmRegisterUserStructField(), DMSwarmRegisterUserDatatypeField()
946d3a51819SDave May @*/
94774d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmRegisterPetscDatatypeField(DM dm,const char fieldname[],PetscInt blocksize,PetscDataType type)
948b62e03f8SDave May {
949b62e03f8SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
950b62e03f8SDave May   size_t         size;
951b62e03f8SDave May 
952521f74f9SMatthew G. Knepley   PetscFunctionBegin;
9532c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!swarm->field_registration_initialized,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMSwarmInitializeFieldRegister() first");
9542c71b3e2SJacob Faibussowitsch   PetscCheckFalse(swarm->field_registration_finalized,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Cannot register additional fields after calling DMSwarmFinalizeFieldRegister() first");
9555f50eb2eSDave May 
9562c71b3e2SJacob Faibussowitsch   PetscCheckFalse(type == PETSC_OBJECT,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
9572c71b3e2SJacob Faibussowitsch   PetscCheckFalse(type == PETSC_FUNCTION,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
9582c71b3e2SJacob Faibussowitsch   PetscCheckFalse(type == PETSC_STRING,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
9592c71b3e2SJacob Faibussowitsch   PetscCheckFalse(type == PETSC_STRUCT,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
9602c71b3e2SJacob Faibussowitsch   PetscCheckFalse(type == PETSC_DATATYPE_UNKNOWN,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
961b62e03f8SDave May 
962*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDataTypeGetSize(type, &size));
963b62e03f8SDave May   /* Load a specific data type into data bucket, specifying textual name and its size in bytes */
964*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketRegisterField(swarm->db,"DMSwarmRegisterPetscDatatypeField",fieldname,blocksize*size,NULL));
96552c3ed93SDave May   {
96677048351SPatrick Sanan     DMSwarmDataField gfield;
96752c3ed93SDave May 
968*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,fieldname,&gfield));
969*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataFieldSetBlockSize(gfield,blocksize));
97052c3ed93SDave May   }
971b62e03f8SDave May   swarm->db->field[swarm->db->nfields-1]->petsc_type = type;
972b62e03f8SDave May   PetscFunctionReturn(0);
973b62e03f8SDave May }
974b62e03f8SDave May 
975d3a51819SDave May /*@C
976d3a51819SDave May    DMSwarmRegisterUserStructField - Register a user defined struct to a DMSwarm
977d3a51819SDave May 
978d083f849SBarry Smith    Collective on dm
979d3a51819SDave May 
980d3a51819SDave May    Input parameters:
98162741f57SDave May +  dm - a DMSwarm
982d3a51819SDave May .  fieldname - the textual name to identify this field
98362741f57SDave May -  size - the size in bytes of the user struct of each data type
984d3a51819SDave May 
985d3a51819SDave May    Level: beginner
986d3a51819SDave May 
987d3a51819SDave May    Notes:
9888b8a3813SDave May    The textual name for each registered field must be unique.
989d3a51819SDave May 
990d3a51819SDave May .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmRegisterUserDatatypeField()
991d3a51819SDave May @*/
99274d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmRegisterUserStructField(DM dm,const char fieldname[],size_t size)
993b62e03f8SDave May {
994b62e03f8SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
995b62e03f8SDave May 
996521f74f9SMatthew G. Knepley   PetscFunctionBegin;
997*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketRegisterField(swarm->db,"DMSwarmRegisterUserStructField",fieldname,size,NULL));
998b62e03f8SDave May   swarm->db->field[swarm->db->nfields-1]->petsc_type = PETSC_STRUCT ;
999b62e03f8SDave May   PetscFunctionReturn(0);
1000b62e03f8SDave May }
1001b62e03f8SDave May 
1002d3a51819SDave May /*@C
1003d3a51819SDave May    DMSwarmRegisterUserDatatypeField - Register a user defined data type to a DMSwarm
1004d3a51819SDave May 
1005d083f849SBarry Smith    Collective on dm
1006d3a51819SDave May 
1007d3a51819SDave May    Input parameters:
100862741f57SDave May +  dm - a DMSwarm
1009d3a51819SDave May .  fieldname - the textual name to identify this field
1010d3a51819SDave May .  size - the size in bytes of the user data type
101162741f57SDave May -  blocksize - the number of each data type
1012d3a51819SDave May 
1013d3a51819SDave May    Level: beginner
1014d3a51819SDave May 
1015d3a51819SDave May    Notes:
10168b8a3813SDave May    The textual name for each registered field must be unique.
1017d3a51819SDave May 
1018d3a51819SDave May .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmRegisterUserStructField(), DMSwarmRegisterUserDatatypeField()
1019d3a51819SDave May @*/
102074d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmRegisterUserDatatypeField(DM dm,const char fieldname[],size_t size,PetscInt blocksize)
1021b62e03f8SDave May {
1022b62e03f8SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
1023b62e03f8SDave May 
1024521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1025*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketRegisterField(swarm->db,"DMSwarmRegisterUserDatatypeField",fieldname,blocksize*size,NULL));
1026320740a0SDave May   {
102777048351SPatrick Sanan     DMSwarmDataField gfield;
1028320740a0SDave May 
1029*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,fieldname,&gfield));
1030*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataFieldSetBlockSize(gfield,blocksize));
1031320740a0SDave May   }
1032b62e03f8SDave May   swarm->db->field[swarm->db->nfields-1]->petsc_type = PETSC_DATATYPE_UNKNOWN;
1033b62e03f8SDave May   PetscFunctionReturn(0);
1034b62e03f8SDave May }
1035b62e03f8SDave May 
1036d3a51819SDave May /*@C
1037d3a51819SDave May    DMSwarmGetField - Get access to the underlying array storing all entries associated with a registered field
1038d3a51819SDave May 
1039d3a51819SDave May    Not collective
1040d3a51819SDave May 
1041d3a51819SDave May    Input parameters:
104262741f57SDave May +  dm - a DMSwarm
104362741f57SDave May -  fieldname - the textual name to identify this field
1044d3a51819SDave May 
1045d3a51819SDave May    Output parameters:
104662741f57SDave May +  blocksize - the number of each data type
1047d3a51819SDave May .  type - the data type
104862741f57SDave May -  data - pointer to raw array
1049d3a51819SDave May 
1050d3a51819SDave May    Level: beginner
1051d3a51819SDave May 
1052d3a51819SDave May    Notes:
10538b8a3813SDave May    The array must be returned using a matching call to DMSwarmRestoreField().
1054d3a51819SDave May 
1055d3a51819SDave May .seealso: DMSwarmRestoreField()
1056d3a51819SDave May @*/
105774d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmGetField(DM dm,const char fieldname[],PetscInt *blocksize,PetscDataType *type,void **data)
1058b62e03f8SDave May {
1059b62e03f8SDave May   DM_Swarm         *swarm = (DM_Swarm*)dm->data;
106077048351SPatrick Sanan   DMSwarmDataField gfield;
1061b62e03f8SDave May 
1062521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1063*5f80ce2aSJacob Faibussowitsch   if (!swarm->issetup) CHKERRQ(DMSetUp(dm));
1064*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,fieldname,&gfield));
1065*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataFieldGetAccess(gfield));
1066*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataFieldGetEntries(gfield,data));
10671b1ea282SDave May   if (blocksize) {*blocksize = gfield->bs; }
1068b5bcf523SDave May   if (type) { *type = gfield->petsc_type; }
1069b62e03f8SDave May   PetscFunctionReturn(0);
1070b62e03f8SDave May }
1071b62e03f8SDave May 
1072d3a51819SDave May /*@C
1073d3a51819SDave May    DMSwarmRestoreField - Restore access to the underlying array storing all entries associated with a registered field
1074d3a51819SDave May 
1075d3a51819SDave May    Not collective
1076d3a51819SDave May 
1077d3a51819SDave May    Input parameters:
107862741f57SDave May +  dm - a DMSwarm
107962741f57SDave May -  fieldname - the textual name to identify this field
1080d3a51819SDave May 
1081d3a51819SDave May    Output parameters:
108262741f57SDave May +  blocksize - the number of each data type
1083d3a51819SDave May .  type - the data type
108462741f57SDave May -  data - pointer to raw array
1085d3a51819SDave May 
1086d3a51819SDave May    Level: beginner
1087d3a51819SDave May 
1088d3a51819SDave May    Notes:
10898b8a3813SDave May    The user must call DMSwarmGetField() prior to calling DMSwarmRestoreField().
1090d3a51819SDave May 
1091d3a51819SDave May .seealso: DMSwarmGetField()
1092d3a51819SDave May @*/
109374d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmRestoreField(DM dm,const char fieldname[],PetscInt *blocksize,PetscDataType *type,void **data)
1094b62e03f8SDave May {
1095b62e03f8SDave May   DM_Swarm         *swarm = (DM_Swarm*)dm->data;
109677048351SPatrick Sanan   DMSwarmDataField gfield;
1097b62e03f8SDave May 
1098521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1099*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,fieldname,&gfield));
1100*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataFieldRestoreAccess(gfield));
1101b62e03f8SDave May   if (data) *data = NULL;
1102b62e03f8SDave May   PetscFunctionReturn(0);
1103b62e03f8SDave May }
1104b62e03f8SDave May 
110574d0cae8SMatthew G. Knepley /*@
1106d3a51819SDave May    DMSwarmAddPoint - Add space for one new point in the DMSwarm
1107d3a51819SDave May 
1108d3a51819SDave May    Not collective
1109d3a51819SDave May 
1110d3a51819SDave May    Input parameter:
1111d3a51819SDave May .  dm - a DMSwarm
1112d3a51819SDave May 
1113d3a51819SDave May    Level: beginner
1114d3a51819SDave May 
1115d3a51819SDave May    Notes:
11168b8a3813SDave May    The new point will have all fields initialized to zero.
1117d3a51819SDave May 
1118d3a51819SDave May .seealso: DMSwarmAddNPoints()
1119d3a51819SDave May @*/
112074d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmAddPoint(DM dm)
1121cb1d1399SDave May {
1122cb1d1399SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
1123cb1d1399SDave May 
1124521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1125*5f80ce2aSJacob Faibussowitsch   if (!swarm->issetup) CHKERRQ(DMSetUp(dm));
1126*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventBegin(DMSWARM_AddPoints,0,0,0,0));
1127*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketAddPoint(swarm->db));
1128*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventEnd(DMSWARM_AddPoints,0,0,0,0));
1129cb1d1399SDave May   PetscFunctionReturn(0);
1130cb1d1399SDave May }
1131cb1d1399SDave May 
113274d0cae8SMatthew G. Knepley /*@
1133d3a51819SDave May    DMSwarmAddNPoints - Add space for a number of new points in the DMSwarm
1134d3a51819SDave May 
1135d3a51819SDave May    Not collective
1136d3a51819SDave May 
1137d3a51819SDave May    Input parameters:
113862741f57SDave May +  dm - a DMSwarm
113962741f57SDave May -  npoints - the number of new points to add
1140d3a51819SDave May 
1141d3a51819SDave May    Level: beginner
1142d3a51819SDave May 
1143d3a51819SDave May    Notes:
11448b8a3813SDave May    The new point will have all fields initialized to zero.
1145d3a51819SDave May 
1146d3a51819SDave May .seealso: DMSwarmAddPoint()
1147d3a51819SDave May @*/
114874d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmAddNPoints(DM dm,PetscInt npoints)
1149cb1d1399SDave May {
1150cb1d1399SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
1151cb1d1399SDave May   PetscInt       nlocal;
1152cb1d1399SDave May 
1153521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1154*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventBegin(DMSWARM_AddPoints,0,0,0,0));
1155*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketGetSizes(swarm->db,&nlocal,NULL,NULL));
1156cb1d1399SDave May   nlocal = nlocal + npoints;
1157*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketSetSizes(swarm->db,nlocal,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT));
1158*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventEnd(DMSWARM_AddPoints,0,0,0,0));
1159cb1d1399SDave May   PetscFunctionReturn(0);
1160cb1d1399SDave May }
1161cb1d1399SDave May 
116274d0cae8SMatthew G. Knepley /*@
1163d3a51819SDave May    DMSwarmRemovePoint - Remove the last point from the DMSwarm
1164d3a51819SDave May 
1165d3a51819SDave May    Not collective
1166d3a51819SDave May 
1167d3a51819SDave May    Input parameter:
1168d3a51819SDave May .  dm - a DMSwarm
1169d3a51819SDave May 
1170d3a51819SDave May    Level: beginner
1171d3a51819SDave May 
1172d3a51819SDave May .seealso: DMSwarmRemovePointAtIndex()
1173d3a51819SDave May @*/
117474d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmRemovePoint(DM dm)
1175cb1d1399SDave May {
1176cb1d1399SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
1177cb1d1399SDave May 
1178521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1179*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventBegin(DMSWARM_RemovePoints,0,0,0,0));
1180*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketRemovePoint(swarm->db));
1181*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventEnd(DMSWARM_RemovePoints,0,0,0,0));
1182cb1d1399SDave May   PetscFunctionReturn(0);
1183cb1d1399SDave May }
1184cb1d1399SDave May 
118574d0cae8SMatthew G. Knepley /*@
1186d3a51819SDave May    DMSwarmRemovePointAtIndex - Removes a specific point from the DMSwarm
1187d3a51819SDave May 
1188d3a51819SDave May    Not collective
1189d3a51819SDave May 
1190d3a51819SDave May    Input parameters:
119162741f57SDave May +  dm - a DMSwarm
119262741f57SDave May -  idx - index of point to remove
1193d3a51819SDave May 
1194d3a51819SDave May    Level: beginner
1195d3a51819SDave May 
1196d3a51819SDave May .seealso: DMSwarmRemovePoint()
1197d3a51819SDave May @*/
119874d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmRemovePointAtIndex(DM dm,PetscInt idx)
1199cb1d1399SDave May {
1200cb1d1399SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
1201cb1d1399SDave May 
1202521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1203*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventBegin(DMSWARM_RemovePoints,0,0,0,0));
1204*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketRemovePointAtIndex(swarm->db,idx));
1205*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventEnd(DMSWARM_RemovePoints,0,0,0,0));
1206cb1d1399SDave May   PetscFunctionReturn(0);
1207cb1d1399SDave May }
1208b62e03f8SDave May 
120974d0cae8SMatthew G. Knepley /*@
1210ba4fc9c6SDave May    DMSwarmCopyPoint - Copy point pj to point pi in the DMSwarm
1211ba4fc9c6SDave May 
1212ba4fc9c6SDave May    Not collective
1213ba4fc9c6SDave May 
1214ba4fc9c6SDave May    Input parameters:
1215ba4fc9c6SDave May +  dm - a DMSwarm
1216ba4fc9c6SDave May .  pi - the index of the point to copy
1217ba4fc9c6SDave May -  pj - the point index where the copy should be located
1218ba4fc9c6SDave May 
1219ba4fc9c6SDave May  Level: beginner
1220ba4fc9c6SDave May 
1221ba4fc9c6SDave May .seealso: DMSwarmRemovePoint()
1222ba4fc9c6SDave May @*/
122374d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmCopyPoint(DM dm,PetscInt pi,PetscInt pj)
1224ba4fc9c6SDave May {
1225ba4fc9c6SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
1226ba4fc9c6SDave May 
1227ba4fc9c6SDave May   PetscFunctionBegin;
1228*5f80ce2aSJacob Faibussowitsch   if (!swarm->issetup) CHKERRQ(DMSetUp(dm));
1229*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketCopyPoint(swarm->db,pi,swarm->db,pj));
1230ba4fc9c6SDave May   PetscFunctionReturn(0);
1231ba4fc9c6SDave May }
1232ba4fc9c6SDave May 
1233095059a4SDave May PetscErrorCode DMSwarmMigrate_Basic(DM dm,PetscBool remove_sent_points)
12343454631fSDave May {
1235521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1236*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmMigrate_Push_Basic(dm,remove_sent_points));
12373454631fSDave May   PetscFunctionReturn(0);
12383454631fSDave May }
12393454631fSDave May 
124074d0cae8SMatthew G. Knepley /*@
1241d3a51819SDave May    DMSwarmMigrate - Relocates points defined in the DMSwarm to other MPI-ranks
1242d3a51819SDave May 
1243d083f849SBarry Smith    Collective on dm
1244d3a51819SDave May 
1245d3a51819SDave May    Input parameters:
124662741f57SDave May +  dm - the DMSwarm
124762741f57SDave May -  remove_sent_points - flag indicating if sent points should be removed from the current MPI-rank
1248d3a51819SDave May 
1249d3a51819SDave May    Notes:
1250a5b23f4aSJose E. Roman    The DM will be modified to accommodate received points.
12518b8a3813SDave May    If remove_sent_points = PETSC_TRUE, any points that were sent will be removed from the DM.
12528b8a3813SDave May    Different styles of migration are supported. See DMSwarmSetMigrateType().
1253d3a51819SDave May 
1254d3a51819SDave May    Level: advanced
1255d3a51819SDave May 
1256d3a51819SDave May .seealso: DMSwarmSetMigrateType()
1257d3a51819SDave May @*/
125874d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmMigrate(DM dm,PetscBool remove_sent_points)
12593454631fSDave May {
1260f0cdbbbaSDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
1261f0cdbbbaSDave May 
1262521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1263*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventBegin(DMSWARM_Migrate,0,0,0,0));
1264f0cdbbbaSDave May   switch (swarm->migrate_type) {
1265f0cdbbbaSDave May     case DMSWARM_MIGRATE_BASIC:
1266*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMSwarmMigrate_Basic(dm,remove_sent_points));
1267f0cdbbbaSDave May       break;
1268f0cdbbbaSDave May     case DMSWARM_MIGRATE_DMCELLNSCATTER:
1269*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMSwarmMigrate_CellDMScatter(dm,remove_sent_points));
1270f0cdbbbaSDave May       break;
1271f0cdbbbaSDave May     case DMSWARM_MIGRATE_DMCELLEXACT:
1272f0cdbbbaSDave May       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_MIGRATE_DMCELLEXACT not implemented");
1273f0cdbbbaSDave May     case DMSWARM_MIGRATE_USER:
1274f0cdbbbaSDave May       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_MIGRATE_USER not implemented");
1275f0cdbbbaSDave May     default:
1276f0cdbbbaSDave May       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_MIGRATE type unknown");
1277f0cdbbbaSDave May   }
1278*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventEnd(DMSWARM_Migrate,0,0,0,0));
1279*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMClearGlobalVectors(dm));
12803454631fSDave May   PetscFunctionReturn(0);
12813454631fSDave May }
12823454631fSDave May 
1283f0cdbbbaSDave May PetscErrorCode DMSwarmMigrate_GlobalToLocal_Basic(DM dm,PetscInt *globalsize);
1284f0cdbbbaSDave May 
1285d3a51819SDave May /*
1286d3a51819SDave May  DMSwarmCollectViewCreate
1287d3a51819SDave May 
1288d3a51819SDave May  * Applies a collection method and gathers point neighbour points into dm
1289d3a51819SDave May 
1290d3a51819SDave May  Notes:
12918b8a3813SDave May  Users should call DMSwarmCollectViewDestroy() after
1292d3a51819SDave May  they have finished computations associated with the collected points
1293d3a51819SDave May */
1294d3a51819SDave May 
129574d0cae8SMatthew G. Knepley /*@
1296d3a51819SDave May    DMSwarmCollectViewCreate - Applies a collection method and gathers points
12975627991aSBarry Smith                               in neighbour ranks into the DMSwarm
1298d3a51819SDave May 
1299d083f849SBarry Smith    Collective on dm
1300d3a51819SDave May 
1301d3a51819SDave May    Input parameter:
1302d3a51819SDave May .  dm - the DMSwarm
1303d3a51819SDave May 
1304d3a51819SDave May    Notes:
1305d3a51819SDave May    Users should call DMSwarmCollectViewDestroy() after
1306d3a51819SDave May    they have finished computations associated with the collected points
13078b8a3813SDave May    Different collect methods are supported. See DMSwarmSetCollectType().
1308d3a51819SDave May 
1309d3a51819SDave May    Level: advanced
1310d3a51819SDave May 
1311d3a51819SDave May .seealso: DMSwarmCollectViewDestroy(), DMSwarmSetCollectType()
1312d3a51819SDave May @*/
131374d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmCollectViewCreate(DM dm)
13142712d1f2SDave May {
13152712d1f2SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
13162712d1f2SDave May   PetscInt       ng;
13172712d1f2SDave May 
1318521f74f9SMatthew G. Knepley   PetscFunctionBegin;
13192c71b3e2SJacob Faibussowitsch   PetscCheckFalse(swarm->collect_view_active,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"CollectView currently active");
1320*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetLocalSize(dm,&ng));
1321480eef7bSDave May   switch (swarm->collect_type) {
1322f0cdbbbaSDave May 
1323480eef7bSDave May     case DMSWARM_COLLECT_BASIC:
1324*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMSwarmMigrate_GlobalToLocal_Basic(dm,&ng));
1325480eef7bSDave May       break;
1326480eef7bSDave May     case DMSWARM_COLLECT_DMDABOUNDINGBOX:
1327f0cdbbbaSDave May       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_COLLECT_DMDABOUNDINGBOX not implemented");
1328480eef7bSDave May     case DMSWARM_COLLECT_GENERAL:
1329f0cdbbbaSDave May       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_COLLECT_GENERAL not implemented");
1330480eef7bSDave May     default:
1331f0cdbbbaSDave May       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_COLLECT type unknown");
1332480eef7bSDave May   }
1333480eef7bSDave May   swarm->collect_view_active = PETSC_TRUE;
1334480eef7bSDave May   swarm->collect_view_reset_nlocal = ng;
13352712d1f2SDave May   PetscFunctionReturn(0);
13362712d1f2SDave May }
13372712d1f2SDave May 
133874d0cae8SMatthew G. Knepley /*@
1339d3a51819SDave May    DMSwarmCollectViewDestroy - Resets the DMSwarm to the size prior to calling DMSwarmCollectViewCreate()
1340d3a51819SDave May 
1341d083f849SBarry Smith    Collective on dm
1342d3a51819SDave May 
1343d3a51819SDave May    Input parameters:
1344d3a51819SDave May .  dm - the DMSwarm
1345d3a51819SDave May 
1346d3a51819SDave May    Notes:
1347d3a51819SDave May    Users should call DMSwarmCollectViewCreate() before this function is called.
1348d3a51819SDave May 
1349d3a51819SDave May    Level: advanced
1350d3a51819SDave May 
1351d3a51819SDave May .seealso: DMSwarmCollectViewCreate(), DMSwarmSetCollectType()
1352d3a51819SDave May @*/
135374d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmCollectViewDestroy(DM dm)
13542712d1f2SDave May {
13552712d1f2SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
13562712d1f2SDave May 
1357521f74f9SMatthew G. Knepley   PetscFunctionBegin;
13582c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!swarm->collect_view_active,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"CollectView is currently not active");
1359*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSetLocalSizes(dm,swarm->collect_view_reset_nlocal,-1));
1360480eef7bSDave May   swarm->collect_view_active = PETSC_FALSE;
13612712d1f2SDave May   PetscFunctionReturn(0);
13622712d1f2SDave May }
13633454631fSDave May 
1364f0cdbbbaSDave May PetscErrorCode DMSwarmSetUpPIC(DM dm)
1365f0cdbbbaSDave May {
1366f0cdbbbaSDave May   PetscInt       dim;
1367f0cdbbbaSDave May 
1368521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1369*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSetNumSpecies(dm, 1));
1370*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(dm,&dim));
13712c71b3e2SJacob Faibussowitsch   PetscCheckFalse(dim < 1,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Dimension must be 1,2,3 - found %D",dim);
13722c71b3e2SJacob Faibussowitsch   PetscCheckFalse(dim > 3,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Dimension must be 1,2,3 - found %D",dim);
1373*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRegisterPetscDatatypeField(dm,DMSwarmPICField_coor,dim,PETSC_DOUBLE));
1374*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRegisterPetscDatatypeField(dm,DMSwarmPICField_cellid,1,PETSC_INT));
1375f0cdbbbaSDave May   PetscFunctionReturn(0);
1376f0cdbbbaSDave May }
1377f0cdbbbaSDave May 
137874d0cae8SMatthew G. Knepley /*@
137931403186SMatthew G. Knepley   DMSwarmSetPointCoordinatesRandom - Sets initial coordinates for particles in each cell
138031403186SMatthew G. Knepley 
138131403186SMatthew G. Knepley   Collective on dm
138231403186SMatthew G. Knepley 
138331403186SMatthew G. Knepley   Input parameters:
138431403186SMatthew G. Knepley + dm  - the DMSwarm
138531403186SMatthew G. Knepley - Npc - The number of particles per cell in the cell DM
138631403186SMatthew G. Knepley 
138731403186SMatthew G. Knepley   Notes:
138831403186SMatthew G. Knepley   The user must use DMSwarmSetCellDM() to set the cell DM first. The particles are placed randomly inside each cell. If only
138931403186SMatthew G. Knepley   one particle is in each cell, it is placed at the centroid.
139031403186SMatthew G. Knepley 
139131403186SMatthew G. Knepley   Level: intermediate
139231403186SMatthew G. Knepley 
139331403186SMatthew G. Knepley .seealso: DMSwarmSetCellDM()
139431403186SMatthew G. Knepley @*/
139531403186SMatthew G. Knepley PetscErrorCode DMSwarmSetPointCoordinatesRandom(DM dm, PetscInt Npc)
139631403186SMatthew G. Knepley {
139731403186SMatthew G. Knepley   DM             cdm;
139831403186SMatthew G. Knepley   PetscRandom    rnd;
139931403186SMatthew G. Knepley   DMPolytopeType ct;
140031403186SMatthew G. Knepley   PetscBool      simplex;
140131403186SMatthew G. Knepley   PetscReal     *centroid, *coords, *xi0, *v0, *J, *invJ, detJ;
140231403186SMatthew G. Knepley   PetscInt       dim, d, cStart, cEnd, c, p;
140331403186SMatthew G. Knepley 
140431403186SMatthew G. Knepley   PetscFunctionBeginUser;
1405*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PetscObjectComm((PetscObject) dm), &rnd));
1406*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetInterval(rnd, -1.0, 1.0));
1407*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetType(rnd, PETSCRAND48));
140831403186SMatthew G. Knepley 
1409*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetCellDM(dm, &cdm));
1410*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(cdm, &dim));
1411*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetHeightStratum(cdm, 0, &cStart, &cEnd));
1412*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetCellType(cdm, cStart, &ct));
141331403186SMatthew G. Knepley   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
141431403186SMatthew G. Knepley 
1415*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc5(dim, &centroid, dim, &xi0, dim, &v0, dim*dim, &J, dim*dim, &invJ));
141631403186SMatthew G. Knepley   for (d = 0; d < dim; ++d) xi0[d] = -1.0;
1417*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetField(dm, DMSwarmPICField_coor, NULL, NULL, (void **) &coords));
141831403186SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
141931403186SMatthew G. Knepley     if (Npc == 1) {
1420*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexComputeCellGeometryFVM(cdm, c, NULL, centroid, NULL));
142131403186SMatthew G. Knepley       for (d = 0; d < dim; ++d) coords[c*dim+d] = centroid[d];
142231403186SMatthew G. Knepley     } else {
1423*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexComputeCellGeometryFEM(cdm, c, NULL, v0, J, invJ, &detJ)); /* affine */
142431403186SMatthew G. Knepley       for (p = 0; p < Npc; ++p) {
142531403186SMatthew G. Knepley         const PetscInt n   = c*Npc + p;
142631403186SMatthew G. Knepley         PetscReal      sum = 0.0, refcoords[3];
142731403186SMatthew G. Knepley 
142831403186SMatthew G. Knepley         for (d = 0; d < dim; ++d) {
1429*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscRandomGetValueReal(rnd, &refcoords[d]));
143031403186SMatthew G. Knepley           sum += refcoords[d];
143131403186SMatthew G. Knepley         }
143231403186SMatthew G. Knepley         if (simplex && sum > 0.0) for (d = 0; d < dim; ++d) refcoords[d] -= PetscSqrtReal(dim)*sum;
143331403186SMatthew G. Knepley         CoordinatesRefToReal(dim, dim, xi0, v0, J, refcoords, &coords[n*dim]);
143431403186SMatthew G. Knepley       }
143531403186SMatthew G. Knepley     }
143631403186SMatthew G. Knepley   }
1437*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRestoreField(dm, DMSwarmPICField_coor, NULL, NULL, (void **) &coords));
1438*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree5(centroid, xi0, v0, J, invJ));
1439*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&rnd));
144031403186SMatthew G. Knepley   PetscFunctionReturn(0);
144131403186SMatthew G. Knepley }
144231403186SMatthew G. Knepley 
144331403186SMatthew G. Knepley /*@
1444d3a51819SDave May    DMSwarmSetType - Set particular flavor of DMSwarm
1445d3a51819SDave May 
1446d083f849SBarry Smith    Collective on dm
1447d3a51819SDave May 
1448d3a51819SDave May    Input parameters:
144962741f57SDave May +  dm - the DMSwarm
145062741f57SDave May -  stype - the DMSwarm type (e.g. DMSWARM_PIC)
1451d3a51819SDave May 
1452d3a51819SDave May    Level: advanced
1453d3a51819SDave May 
14545627991aSBarry Smith .seealso: DMSwarmSetMigrateType(), DMSwarmSetCollectType(), DMSwarmType, DMSWARM_PIC, DMSWARM_BASIC
1455d3a51819SDave May @*/
145674d0cae8SMatthew G. Knepley PetscErrorCode DMSwarmSetType(DM dm,DMSwarmType stype)
1457f0cdbbbaSDave May {
1458f0cdbbbaSDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
1459f0cdbbbaSDave May 
1460521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1461f0cdbbbaSDave May   swarm->swarm_type = stype;
1462f0cdbbbaSDave May   if (swarm->swarm_type == DMSWARM_PIC) {
1463*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmSetUpPIC(dm));
1464f0cdbbbaSDave May   }
1465f0cdbbbaSDave May   PetscFunctionReturn(0);
1466f0cdbbbaSDave May }
1467f0cdbbbaSDave May 
14683454631fSDave May PetscErrorCode DMSetup_Swarm(DM dm)
14693454631fSDave May {
14703454631fSDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
14713454631fSDave May   PetscMPIInt    rank;
14723454631fSDave May   PetscInt       p,npoints,*rankval;
14733454631fSDave May 
1474521f74f9SMatthew G. Knepley   PetscFunctionBegin;
14753454631fSDave May   if (swarm->issetup) PetscFunctionReturn(0);
14763454631fSDave May   swarm->issetup = PETSC_TRUE;
14773454631fSDave May 
1478f0cdbbbaSDave May   if (swarm->swarm_type == DMSWARM_PIC) {
1479f0cdbbbaSDave May     /* check dmcell exists */
14802c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!swarm->dmcell,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"DMSWARM_PIC requires you call DMSwarmSetCellDM");
1481f0cdbbbaSDave May 
1482f0cdbbbaSDave May     if (swarm->dmcell->ops->locatepointssubdomain) {
1483f0cdbbbaSDave May       /* check methods exists for exact ownership identificiation */
1484*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscInfo(dm, "DMSWARM_PIC: Using method CellDM->ops->LocatePointsSubdomain\n"));
1485f0cdbbbaSDave May       swarm->migrate_type = DMSWARM_MIGRATE_DMCELLEXACT;
1486f0cdbbbaSDave May     } else {
1487f0cdbbbaSDave May       /* check methods exist for point location AND rank neighbor identification */
1488f0cdbbbaSDave May       if (swarm->dmcell->ops->locatepoints) {
1489*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscInfo(dm, "DMSWARM_PIC: Using method CellDM->LocatePoints\n"));
1490f0cdbbbaSDave May       } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"DMSWARM_PIC requires the method CellDM->ops->locatepoints be defined");
1491f0cdbbbaSDave May 
1492f0cdbbbaSDave May       if (swarm->dmcell->ops->getneighbors) {
1493*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscInfo(dm, "DMSWARM_PIC: Using method CellDM->GetNeigbors\n"));
1494f0cdbbbaSDave May       } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"DMSWARM_PIC requires the method CellDM->ops->getneighbors be defined");
1495f0cdbbbaSDave May 
1496f0cdbbbaSDave May       swarm->migrate_type = DMSWARM_MIGRATE_DMCELLNSCATTER;
1497f0cdbbbaSDave May     }
1498f0cdbbbaSDave May   }
1499f0cdbbbaSDave May 
1500*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmFinalizeFieldRegister(dm));
1501f0cdbbbaSDave May 
15023454631fSDave May   /* check some fields were registered */
15032c71b3e2SJacob Faibussowitsch   PetscCheckFalse(swarm->db->nfields <= 2,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"At least one field user must be registered via DMSwarmRegisterXXX()");
15043454631fSDave May 
15053454631fSDave May   /* check local sizes were set */
15062c71b3e2SJacob Faibussowitsch   PetscCheckFalse(swarm->db->L == -1,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Local sizes must be set via DMSwarmSetLocalSizes()");
15073454631fSDave May 
15083454631fSDave May   /* initialize values in pid and rank placeholders */
15093454631fSDave May   /* TODO: [pid - use MPI_Scan] */
1510*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank));
1511*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL));
1512*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval));
15133454631fSDave May   for (p=0; p<npoints; p++) {
15143454631fSDave May     rankval[p] = (PetscInt)rank;
15153454631fSDave May   }
1516*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval));
15173454631fSDave May   PetscFunctionReturn(0);
15183454631fSDave May }
15193454631fSDave May 
1520dc5f5ce9SDave May extern PetscErrorCode DMSwarmSortDestroy(DMSwarmSort *_ctx);
1521dc5f5ce9SDave May 
152257795646SDave May PetscErrorCode DMDestroy_Swarm(DM dm)
152357795646SDave May {
152457795646SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
152557795646SDave May 
152657795646SDave May   PetscFunctionBegin;
1527d0c080abSJoseph Pusztay   if (--swarm->refct > 0) PetscFunctionReturn(0);
1528*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketDestroy(&swarm->db));
1529dc5f5ce9SDave May   if (swarm->sort_context) {
1530*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmSortDestroy(&swarm->sort_context));
1531dc5f5ce9SDave May   }
1532*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(swarm));
153357795646SDave May   PetscFunctionReturn(0);
153457795646SDave May }
153557795646SDave May 
1536a9ee3421SMatthew G. Knepley PetscErrorCode DMSwarmView_Draw(DM dm, PetscViewer viewer)
1537a9ee3421SMatthew G. Knepley {
1538a9ee3421SMatthew G. Knepley   DM             cdm;
1539a9ee3421SMatthew G. Knepley   PetscDraw      draw;
154031403186SMatthew G. Knepley   PetscReal     *coords, oldPause, radius = 0.01;
1541a9ee3421SMatthew G. Knepley   PetscInt       Np, p, bs;
1542a9ee3421SMatthew G. Knepley 
1543a9ee3421SMatthew G. Knepley   PetscFunctionBegin;
1544*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL, ((PetscObject) dm)->prefix, "-dm_view_swarm_radius", &radius, NULL));
1545*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDrawGetDraw(viewer, 0, &draw));
1546*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetCellDM(dm, &cdm));
1547*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDrawGetPause(draw, &oldPause));
1548*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDrawSetPause(draw, 0.0));
1549*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMView(cdm, viewer));
1550*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDrawSetPause(draw, oldPause));
1551a9ee3421SMatthew G. Knepley 
1552*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetLocalSize(dm, &Np));
1553*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetField(dm, DMSwarmPICField_coor, &bs, NULL, (void **) &coords));
1554a9ee3421SMatthew G. Knepley   for (p = 0; p < Np; ++p) {
1555a9ee3421SMatthew G. Knepley     const PetscInt i = p*bs;
1556a9ee3421SMatthew G. Knepley 
1557*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDrawEllipse(draw, coords[i], coords[i+1], radius, radius, PETSC_DRAW_BLUE));
1558a9ee3421SMatthew G. Knepley   }
1559*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRestoreField(dm, DMSwarmPICField_coor, &bs, NULL, (void **) &coords));
1560*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDrawFlush(draw));
1561*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDrawPause(draw));
1562*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDrawSave(draw));
1563a9ee3421SMatthew G. Knepley   PetscFunctionReturn(0);
1564a9ee3421SMatthew G. Knepley }
1565a9ee3421SMatthew G. Knepley 
15665f50eb2eSDave May PetscErrorCode DMView_Swarm(DM dm, PetscViewer viewer)
15675f50eb2eSDave May {
15685f50eb2eSDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
15695627991aSBarry Smith   PetscBool      iascii,ibinary,isvtk,isdraw;
15705627991aSBarry Smith #if defined(PETSC_HAVE_HDF5)
15715627991aSBarry Smith   PetscBool      ishdf5;
15725627991aSBarry Smith #endif
15735f50eb2eSDave May 
15745f50eb2eSDave May   PetscFunctionBegin;
15755f50eb2eSDave May   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
15765f50eb2eSDave May   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1577*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1578*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY,&ibinary));
1579*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERVTK,   &isvtk));
15805627991aSBarry Smith #if defined(PETSC_HAVE_HDF5)
1581*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5,  &ishdf5));
15825627991aSBarry Smith #endif
1583*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW,  &isdraw));
15845f50eb2eSDave May   if (iascii) {
1585*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataBucketView(PetscObjectComm((PetscObject)dm),swarm->db,NULL,DATABUCKET_VIEW_STDOUT));
15862c71b3e2SJacob Faibussowitsch   } else PetscCheckFalse(ibinary,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"NO Binary support");
15875f50eb2eSDave May #if defined(PETSC_HAVE_HDF5)
158874d0cae8SMatthew G. Knepley   else if (ishdf5) {
1589*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmView_HDF5(dm, viewer));
159074d0cae8SMatthew G. Knepley   }
15915f50eb2eSDave May #endif
1592298827fbSBarry Smith   else if (isdraw) {
1593*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmView_Draw(dm, viewer));
15945f50eb2eSDave May   }
15955f50eb2eSDave May   PetscFunctionReturn(0);
15965f50eb2eSDave May }
15975f50eb2eSDave May 
1598d0c080abSJoseph Pusztay /*@C
1599d0c080abSJoseph Pusztay    DMSwarmGetCellSwarm - Extracts a single cell from the DMSwarm object, returns it as a single cell DMSwarm.
1600d0c080abSJoseph 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.
1601d0c080abSJoseph Pusztay 
1602d0c080abSJoseph Pusztay    Important: Changes to this cell of the swarm will be lost if they are made prior to restoring this cell.
1603d0c080abSJoseph Pusztay 
1604d0c080abSJoseph Pusztay    Noncollective
1605d0c080abSJoseph Pusztay 
1606d0c080abSJoseph Pusztay    Input parameters:
1607d0c080abSJoseph Pusztay +  sw - the DMSwarm
16085627991aSBarry Smith .  cellID - the integer id of the cell to be extracted and filtered
16095627991aSBarry Smith -  cellswarm - The DMSwarm to receive the cell
1610d0c080abSJoseph Pusztay 
1611d0c080abSJoseph Pusztay    Level: beginner
1612d0c080abSJoseph Pusztay 
16135627991aSBarry Smith    Notes:
16145627991aSBarry Smith       This presently only supports DMSWARM_PIC type
16155627991aSBarry Smith 
16165627991aSBarry Smith       Should be restored with DMSwarmRestoreCellSwarm()
1617d0c080abSJoseph Pusztay 
1618d0c080abSJoseph Pusztay .seealso: DMSwarmRestoreCellSwarm()
1619d0c080abSJoseph Pusztay @*/
1620d0c080abSJoseph Pusztay PETSC_EXTERN PetscErrorCode DMSwarmGetCellSwarm(DM sw, PetscInt cellID, DM cellswarm)
1621d0c080abSJoseph Pusztay {
1622d0c080abSJoseph Pusztay   DM_Swarm      *original = (DM_Swarm*) sw->data;
1623d0c080abSJoseph Pusztay   DMLabel        label;
1624d0c080abSJoseph Pusztay   DM             dmc, subdmc;
1625d0c080abSJoseph Pusztay   PetscInt      *pids, particles, dim;
1626d0c080abSJoseph Pusztay 
1627d0c080abSJoseph Pusztay   PetscFunctionBegin;
1628d0c080abSJoseph Pusztay   /* Configure new swarm */
1629*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetType(cellswarm, DMSWARM));
1630*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(sw, &dim));
1631*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetDimension(cellswarm, dim));
1632*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSetType(cellswarm, DMSWARM_PIC));
1633d0c080abSJoseph Pusztay   /* Destroy the unused, unconfigured data bucket to prevent stragglers in memory */
1634*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketDestroy(&((DM_Swarm*)cellswarm->data)->db));
1635*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSortGetAccess(sw));
1636*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSortGetNumberOfPointsPerCell(sw, cellID, &particles));
1637*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSortGetPointsPerCell(sw, cellID, &particles, &pids));
1638*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketCreateFromSubset(original->db, particles, pids, &((DM_Swarm*)cellswarm->data)->db));
1639*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSortRestoreAccess(sw));
1640*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(pids));
1641*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetCellDM(sw, &dmc));
1642*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLabelCreate(PetscObjectComm((PetscObject)sw), "singlecell", &label));
1643*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMAddLabel(dmc, label));
1644*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLabelSetValue(label, cellID, 1));
1645*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexFilter(dmc, label, 1, &subdmc));
1646*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSetCellDM(cellswarm, subdmc));
1647*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLabelDestroy(&label));
1648d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
1649d0c080abSJoseph Pusztay }
1650d0c080abSJoseph Pusztay 
1651d0c080abSJoseph Pusztay /*@C
16525627991aSBarry Smith    DMSwarmRestoreCellSwarm - Restores a DMSwarm object obtained with DMSwarmGetCellSwarm(). All fields are copied back into the parent swarm.
1653d0c080abSJoseph Pusztay 
1654d0c080abSJoseph Pusztay    Noncollective
1655d0c080abSJoseph Pusztay 
1656d0c080abSJoseph Pusztay    Input parameters:
1657d0c080abSJoseph Pusztay +  sw - the parent DMSwarm
1658d0c080abSJoseph Pusztay .  cellID - the integer id of the cell to be copied back into the parent swarm
1659d0c080abSJoseph Pusztay -  cellswarm - the cell swarm object
1660d0c080abSJoseph Pusztay 
1661d0c080abSJoseph Pusztay    Level: beginner
1662d0c080abSJoseph Pusztay 
16635627991aSBarry Smith    Note:
16645627991aSBarry Smith     This only supports DMSWARM_PIC types of DMSwarms
1665d0c080abSJoseph Pusztay 
1666d0c080abSJoseph Pusztay .seealso: DMSwarmGetCellSwarm()
1667d0c080abSJoseph Pusztay @*/
1668d0c080abSJoseph Pusztay PETSC_EXTERN PetscErrorCode DMSwarmRestoreCellSwarm(DM sw, PetscInt cellID, DM cellswarm)
1669d0c080abSJoseph Pusztay {
1670d0c080abSJoseph Pusztay   DM             dmc;
1671d0c080abSJoseph Pusztay   PetscInt       *pids, particles, p;
1672d0c080abSJoseph Pusztay 
1673d0c080abSJoseph Pusztay   PetscFunctionBegin;
1674*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSortGetAccess(sw));
1675*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSortGetPointsPerCell(sw, cellID, &particles, &pids));
1676*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSortRestoreAccess(sw));
1677d0c080abSJoseph Pusztay   /* Pointwise copy of each particle based on pid. The parent swarm may not be altered during this process. */
1678d0c080abSJoseph Pusztay   for (p=0; p<particles; ++p) {
1679*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataBucketCopyPoint(((DM_Swarm*)cellswarm->data)->db,pids[p],((DM_Swarm*)sw->data)->db,pids[p]));
1680d0c080abSJoseph Pusztay   }
1681d0c080abSJoseph Pusztay   /* Free memory, destroy cell dm */
1682*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetCellDM(cellswarm, &dmc));
1683*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dmc));
1684*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(pids));
1685d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
1686d0c080abSJoseph Pusztay }
1687d0c080abSJoseph Pusztay 
1688d0c080abSJoseph Pusztay PETSC_INTERN PetscErrorCode DMClone_Swarm(DM, DM *);
1689d0c080abSJoseph Pusztay 
1690d0c080abSJoseph Pusztay static PetscErrorCode DMInitialize_Swarm(DM sw)
1691d0c080abSJoseph Pusztay {
1692d0c080abSJoseph Pusztay   PetscFunctionBegin;
1693d0c080abSJoseph Pusztay   sw->dim  = 0;
1694d0c080abSJoseph Pusztay   sw->ops->view                            = DMView_Swarm;
1695d0c080abSJoseph Pusztay   sw->ops->load                            = NULL;
1696d0c080abSJoseph Pusztay   sw->ops->setfromoptions                  = NULL;
1697d0c080abSJoseph Pusztay   sw->ops->clone                           = DMClone_Swarm;
1698d0c080abSJoseph Pusztay   sw->ops->setup                           = DMSetup_Swarm;
1699d0c080abSJoseph Pusztay   sw->ops->createlocalsection              = NULL;
1700d0c080abSJoseph Pusztay   sw->ops->createdefaultconstraints        = NULL;
1701d0c080abSJoseph Pusztay   sw->ops->createglobalvector              = DMCreateGlobalVector_Swarm;
1702d0c080abSJoseph Pusztay   sw->ops->createlocalvector               = DMCreateLocalVector_Swarm;
1703d0c080abSJoseph Pusztay   sw->ops->getlocaltoglobalmapping         = NULL;
1704d0c080abSJoseph Pusztay   sw->ops->createfieldis                   = NULL;
1705d0c080abSJoseph Pusztay   sw->ops->createcoordinatedm              = NULL;
1706d0c080abSJoseph Pusztay   sw->ops->getcoloring                     = NULL;
1707d0c080abSJoseph Pusztay   sw->ops->creatematrix                    = DMCreateMatrix_Swarm;
1708d0c080abSJoseph Pusztay   sw->ops->createinterpolation             = NULL;
1709d0c080abSJoseph Pusztay   sw->ops->createinjection                 = NULL;
1710d0c080abSJoseph Pusztay   sw->ops->createmassmatrix                = DMCreateMassMatrix_Swarm;
1711d0c080abSJoseph Pusztay   sw->ops->refine                          = NULL;
1712d0c080abSJoseph Pusztay   sw->ops->coarsen                         = NULL;
1713d0c080abSJoseph Pusztay   sw->ops->refinehierarchy                 = NULL;
1714d0c080abSJoseph Pusztay   sw->ops->coarsenhierarchy                = NULL;
1715d0c080abSJoseph Pusztay   sw->ops->globaltolocalbegin              = NULL;
1716d0c080abSJoseph Pusztay   sw->ops->globaltolocalend                = NULL;
1717d0c080abSJoseph Pusztay   sw->ops->localtoglobalbegin              = NULL;
1718d0c080abSJoseph Pusztay   sw->ops->localtoglobalend                = NULL;
1719d0c080abSJoseph Pusztay   sw->ops->destroy                         = DMDestroy_Swarm;
1720d0c080abSJoseph Pusztay   sw->ops->createsubdm                     = NULL;
1721d0c080abSJoseph Pusztay   sw->ops->getdimpoints                    = NULL;
1722d0c080abSJoseph Pusztay   sw->ops->locatepoints                    = NULL;
1723d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
1724d0c080abSJoseph Pusztay }
1725d0c080abSJoseph Pusztay 
1726d0c080abSJoseph Pusztay PETSC_INTERN PetscErrorCode DMClone_Swarm(DM dm, DM *newdm)
1727d0c080abSJoseph Pusztay {
1728d0c080abSJoseph Pusztay   DM_Swarm       *swarm = (DM_Swarm *) dm->data;
1729d0c080abSJoseph Pusztay 
1730d0c080abSJoseph Pusztay   PetscFunctionBegin;
1731d0c080abSJoseph Pusztay   swarm->refct++;
1732d0c080abSJoseph Pusztay   (*newdm)->data = swarm;
1733*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectChangeTypeName((PetscObject) *newdm, DMSWARM));
1734*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMInitialize_Swarm(*newdm));
17352e56dffeSJoe Pusztay   (*newdm)->dim = dm->dim;
1736d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
1737d0c080abSJoseph Pusztay }
1738d0c080abSJoseph Pusztay 
1739d3a51819SDave May /*MC
1740d3a51819SDave May 
1741d3a51819SDave May  DMSWARM = "swarm" - A DM object used to represent arrays of data (fields) of arbitrary data type.
174262741f57SDave May  This implementation was designed for particle methods in which the underlying
1743d3a51819SDave May  data required to be represented is both (i) dynamic in length, (ii) and of arbitrary data type.
1744d3a51819SDave May 
174562741f57SDave May  User data can be represented by DMSwarm through a registering "fields".
174662741f57SDave May  To register a field, the user must provide;
174762741f57SDave May  (a) a unique name;
174862741f57SDave May  (b) the data type (or size in bytes);
174962741f57SDave May  (c) the block size of the data.
1750d3a51819SDave May 
1751d3a51819SDave May  For example, suppose the application requires a unique id, energy, momentum and density to be stored
1752c215a47eSMatthew Knepley  on a set of particles. Then the following code could be used
1753d3a51819SDave May 
175462741f57SDave May $    DMSwarmInitializeFieldRegister(dm)
175562741f57SDave May $    DMSwarmRegisterPetscDatatypeField(dm,"uid",1,PETSC_LONG);
175662741f57SDave May $    DMSwarmRegisterPetscDatatypeField(dm,"energy",1,PETSC_REAL);
175762741f57SDave May $    DMSwarmRegisterPetscDatatypeField(dm,"momentum",3,PETSC_REAL);
175862741f57SDave May $    DMSwarmRegisterPetscDatatypeField(dm,"density",1,PETSC_FLOAT);
175962741f57SDave May $    DMSwarmFinalizeFieldRegister(dm)
1760d3a51819SDave May 
1761d3a51819SDave May  The fields represented by DMSwarm are dynamic and can be re-sized at any time.
176262741f57SDave May  The only restriction imposed by DMSwarm is that all fields contain the same number of points.
1763d3a51819SDave May 
1764d3a51819SDave May  To support particle methods, "migration" techniques are provided. These methods migrate data
17655627991aSBarry Smith  between ranks.
1766d3a51819SDave May 
1767d3a51819SDave May  DMSwarm supports the methods DMCreateGlobalVector() and DMCreateLocalVector().
1768d3a51819SDave May  As a DMSwarm may internally define and store values of different data types,
176962741f57SDave May  before calling DMCreateGlobalVector() or DMCreateLocalVector(), the user must inform DMSwarm which
1770d3a51819SDave May  fields should be used to define a Vec object via
1771d3a51819SDave May    DMSwarmVectorDefineField()
1772c215a47eSMatthew Knepley  The specified field can be changed at any time - thereby permitting vectors
1773c215a47eSMatthew Knepley  compatible with different fields to be created.
1774d3a51819SDave May 
177562741f57SDave May  A dual representation of fields in the DMSwarm and a Vec object is permitted via
1776d3a51819SDave May    DMSwarmCreateGlobalVectorFromField()
1777d3a51819SDave May  Here the data defining the field in the DMSwarm is shared with a Vec.
1778d3a51819SDave May  This is inherently unsafe if you alter the size of the field at any time between
1779d3a51819SDave May  calls to DMSwarmCreateGlobalVectorFromField() and DMSwarmDestroyGlobalVectorFromField().
1780cc651181SDave May  If the local size of the DMSwarm does not match the local size of the global vector
1781cc651181SDave May  when DMSwarmDestroyGlobalVectorFromField() is called, an error is thrown.
1782d3a51819SDave May 
178362741f57SDave May  Additional high-level support is provided for Particle-In-Cell methods.
178462741f57SDave May  Please refer to the man page for DMSwarmSetType().
178562741f57SDave May 
1786d3a51819SDave May  Level: beginner
1787d3a51819SDave May 
1788d3a51819SDave May .seealso: DMType, DMCreate(), DMSetType()
1789d3a51819SDave May M*/
179057795646SDave May PETSC_EXTERN PetscErrorCode DMCreate_Swarm(DM dm)
179157795646SDave May {
179257795646SDave May   DM_Swarm      *swarm;
179357795646SDave May 
179457795646SDave May   PetscFunctionBegin;
179557795646SDave May   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1796*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscNewLog(dm,&swarm));
1797f0cdbbbaSDave May   dm->data = swarm;
1798*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketCreate(&swarm->db));
1799*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmInitializeFieldRegister(dm));
1800d0c080abSJoseph Pusztay   swarm->refct = 1;
1801b5bcf523SDave May   swarm->vec_field_set                  = PETSC_FALSE;
18023454631fSDave May   swarm->issetup                        = PETSC_FALSE;
1803480eef7bSDave May   swarm->swarm_type                     = DMSWARM_BASIC;
1804480eef7bSDave May   swarm->migrate_type                   = DMSWARM_MIGRATE_BASIC;
1805480eef7bSDave May   swarm->collect_type                   = DMSWARM_COLLECT_BASIC;
180640c453e9SDave May   swarm->migrate_error_on_missing_point = PETSC_FALSE;
1807f0cdbbbaSDave May   swarm->dmcell                         = NULL;
1808f0cdbbbaSDave May   swarm->collect_view_active            = PETSC_FALSE;
1809f0cdbbbaSDave May   swarm->collect_view_reset_nlocal      = -1;
1810*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMInitialize_Swarm(dm));
181157795646SDave May   PetscFunctionReturn(0);
181257795646SDave May }
1813