157795646SDave May #define PETSCDM_DLL 257795646SDave May #include <petsc/private/dmswarmimpl.h> /*I "petscdmswarm.h" I*/ 3e8f14785SLisandro Dalcin #include <petsc/private/hashsetij.h> 40643ed39SMatthew G. Knepley #include <petsc/private/petscfeimpl.h> 55917a6f0SStefano Zampini #include <petscviewer.h> 65917a6f0SStefano Zampini #include <petscdraw.h> 783c47955SMatthew G. Knepley #include <petscdmplex.h> 84cc7f7b2SMatthew G. Knepley #include <petscblaslapack.h> 9279f676cSBarry Smith #include "../src/dm/impls/swarm/data_bucket.h" 10d0c080abSJoseph Pusztay #include <petscdmlabel.h> 11d0c080abSJoseph Pusztay #include <petscsection.h> 1257795646SDave May 13f2b2bee7SDave May PetscLogEvent DMSWARM_Migrate, DMSWARM_SetSizes, DMSWARM_AddPoints, DMSWARM_RemovePoints, DMSWARM_Sort; 14ed923d71SDave May PetscLogEvent DMSWARM_DataExchangerTopologySetup, DMSWARM_DataExchangerBegin, DMSWARM_DataExchangerEnd; 15ed923d71SDave May PetscLogEvent DMSWARM_DataExchangerSendCount, DMSWARM_DataExchangerPack; 16ed923d71SDave May 17ea78f98cSLisandro Dalcin const char *DMSwarmTypeNames[] = {"basic", "pic", NULL}; 18ea78f98cSLisandro Dalcin const char *DMSwarmMigrateTypeNames[] = {"basic", "dmcellnscatter", "dmcellexact", "user", NULL}; 19ea78f98cSLisandro Dalcin const char *DMSwarmCollectTypeNames[] = {"basic", "boundingbox", "general", "user", NULL}; 20ea78f98cSLisandro Dalcin const char *DMSwarmPICLayoutTypeNames[] = {"regular", "gauss", "subdivision", NULL}; 21f0cdbbbaSDave May 22f0cdbbbaSDave May const char DMSwarmField_pid[] = "DMSwarm_pid"; 23f0cdbbbaSDave May const char DMSwarmField_rank[] = "DMSwarm_rank"; 24f0cdbbbaSDave May const char DMSwarmPICField_coor[] = "DMSwarmPIC_coor"; 25e2d107dbSDave May const char DMSwarmPICField_cellid[] = "DMSwarm_cellid"; 26f0cdbbbaSDave May 272e956fe4SStefano Zampini PetscInt SwarmDataFieldId = -1; 282e956fe4SStefano Zampini 2974d0cae8SMatthew G. Knepley #if defined(PETSC_HAVE_HDF5) 3074d0cae8SMatthew G. Knepley #include <petscviewerhdf5.h> 3174d0cae8SMatthew G. Knepley 3266976f2fSJacob Faibussowitsch static PetscErrorCode VecView_Swarm_HDF5_Internal(Vec v, PetscViewer viewer) 33d71ae5a4SJacob Faibussowitsch { 3474d0cae8SMatthew G. Knepley DM dm; 3574d0cae8SMatthew G. Knepley PetscReal seqval; 3674d0cae8SMatthew G. Knepley PetscInt seqnum, bs; 37eb0c6899SMatthew G. Knepley PetscBool isseq, ists; 3874d0cae8SMatthew G. Knepley 3974d0cae8SMatthew G. Knepley PetscFunctionBegin; 409566063dSJacob Faibussowitsch PetscCall(VecGetDM(v, &dm)); 419566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(v, &bs)); 429566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5PushGroup(viewer, "/particle_fields")); 439566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)v, VECSEQ, &isseq)); 449566063dSJacob Faibussowitsch PetscCall(DMGetOutputSequenceNumber(dm, &seqnum, &seqval)); 45eb0c6899SMatthew G. Knepley PetscCall(PetscViewerHDF5IsTimestepping(viewer, &ists)); 46eb0c6899SMatthew G. Knepley if (ists) PetscCall(PetscViewerHDF5SetTimestep(viewer, seqnum)); 479566063dSJacob Faibussowitsch /* PetscCall(DMSequenceView_HDF5(dm, "time", seqnum, (PetscScalar) seqval, viewer)); */ 489566063dSJacob Faibussowitsch PetscCall(VecViewNative(v, viewer)); 499566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject)v, "Nc", PETSC_INT, (void *)&bs)); 509566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5PopGroup(viewer)); 513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5274d0cae8SMatthew G. Knepley } 5374d0cae8SMatthew G. Knepley 5466976f2fSJacob Faibussowitsch static PetscErrorCode DMSwarmView_HDF5(DM dm, PetscViewer viewer) 55d71ae5a4SJacob Faibussowitsch { 5674d0cae8SMatthew G. Knepley Vec coordinates; 5774d0cae8SMatthew G. Knepley PetscInt Np; 5874d0cae8SMatthew G. Knepley PetscBool isseq; 5974d0cae8SMatthew G. Knepley 6074d0cae8SMatthew G. Knepley PetscFunctionBegin; 619566063dSJacob Faibussowitsch PetscCall(DMSwarmGetSize(dm, &Np)); 629566063dSJacob Faibussowitsch PetscCall(DMSwarmCreateGlobalVectorFromField(dm, DMSwarmPICField_coor, &coordinates)); 639566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates")); 649566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5PushGroup(viewer, "/particles")); 659566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)coordinates, VECSEQ, &isseq)); 669566063dSJacob Faibussowitsch PetscCall(VecViewNative(coordinates, viewer)); 679566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject)coordinates, "Np", PETSC_INT, (void *)&Np)); 689566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5PopGroup(viewer)); 699566063dSJacob Faibussowitsch PetscCall(DMSwarmDestroyGlobalVectorFromField(dm, DMSwarmPICField_coor, &coordinates)); 703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7174d0cae8SMatthew G. Knepley } 7274d0cae8SMatthew G. Knepley #endif 7374d0cae8SMatthew G. Knepley 7466976f2fSJacob Faibussowitsch static PetscErrorCode VecView_Swarm(Vec v, PetscViewer viewer) 75d71ae5a4SJacob Faibussowitsch { 7674d0cae8SMatthew G. Knepley DM dm; 77f9558f5cSBarry Smith #if defined(PETSC_HAVE_HDF5) 7874d0cae8SMatthew G. Knepley PetscBool ishdf5; 79f9558f5cSBarry Smith #endif 8074d0cae8SMatthew G. Knepley 8174d0cae8SMatthew G. Knepley PetscFunctionBegin; 829566063dSJacob Faibussowitsch PetscCall(VecGetDM(v, &dm)); 8328b400f6SJacob Faibussowitsch PetscCheck(dm, PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_WRONG, "Vector not generated from a DM"); 84f9558f5cSBarry Smith #if defined(PETSC_HAVE_HDF5) 859566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5)); 8674d0cae8SMatthew G. Knepley if (ishdf5) { 879566063dSJacob Faibussowitsch PetscCall(VecView_Swarm_HDF5_Internal(v, viewer)); 883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8974d0cae8SMatthew G. Knepley } 90f9558f5cSBarry Smith #endif 919566063dSJacob Faibussowitsch PetscCall(VecViewNative(v, viewer)); 923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9374d0cae8SMatthew G. Knepley } 9474d0cae8SMatthew G. Knepley 95cc4c1da9SBarry Smith /*@ 960bf7c1c5SMatthew G. Knepley DMSwarmVectorGetField - Gets the field from which to define a `Vec` object 970bf7c1c5SMatthew G. Knepley when `DMCreateLocalVector()`, or `DMCreateGlobalVector()` is called 980bf7c1c5SMatthew G. Knepley 990bf7c1c5SMatthew G. Knepley Not collective 1000bf7c1c5SMatthew G. Knepley 1010bf7c1c5SMatthew G. Knepley Input Parameter: 1020bf7c1c5SMatthew G. Knepley . dm - a `DMSWARM` 1030bf7c1c5SMatthew G. Knepley 1040bf7c1c5SMatthew G. Knepley Output Parameter: 1050bf7c1c5SMatthew G. Knepley . fieldname - the textual name given to a registered field, or the empty string if it has not been set 1060bf7c1c5SMatthew G. Knepley 1070bf7c1c5SMatthew G. Knepley Level: beginner 1080bf7c1c5SMatthew G. Knepley 1090bf7c1c5SMatthew G. Knepley .seealso: `DM`, `DMSWARM`, `DMSwarmVectorDefineField()` `DMSwarmRegisterPetscDatatypeField()`, `DMCreateGlobalVector()`, `DMCreateLocalVector()` 1100bf7c1c5SMatthew G. Knepley @*/ 1110bf7c1c5SMatthew G. Knepley PetscErrorCode DMSwarmVectorGetField(DM dm, const char *fieldname[]) 1120bf7c1c5SMatthew G. Knepley { 1130bf7c1c5SMatthew G. Knepley PetscFunctionBegin; 1140bf7c1c5SMatthew G. Knepley PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM); 1150bf7c1c5SMatthew G. Knepley PetscAssertPointer(fieldname, 2); 1160bf7c1c5SMatthew G. Knepley *fieldname = ((DM_Swarm *)dm->data)->vec_field_name; 1170bf7c1c5SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 1180bf7c1c5SMatthew G. Knepley } 1190bf7c1c5SMatthew G. Knepley 120cc4c1da9SBarry Smith /*@ 12120f4b53cSBarry Smith DMSwarmVectorDefineField - Sets the field from which to define a `Vec` object 12220f4b53cSBarry Smith when `DMCreateLocalVector()`, or `DMCreateGlobalVector()` is called 12357795646SDave May 12420f4b53cSBarry Smith Collective 12557795646SDave May 12660225df5SJacob Faibussowitsch Input Parameters: 12720f4b53cSBarry Smith + dm - a `DMSWARM` 12862741f57SDave May - fieldname - the textual name given to a registered field 12957795646SDave May 130d3a51819SDave May Level: beginner 13157795646SDave May 132d3a51819SDave May Notes: 13320f4b53cSBarry Smith The field with name `fieldname` must be defined as having a data type of `PetscScalar`. 134e7af74c8SDave May 13520f4b53cSBarry Smith This function must be called prior to calling `DMCreateLocalVector()`, `DMCreateGlobalVector()`. 13620f4b53cSBarry Smith Multiple calls to `DMSwarmVectorDefineField()` are permitted. 137e7af74c8SDave May 1380bf7c1c5SMatthew G. Knepley .seealso: `DM`, `DMSWARM`, `DMSwarmVectorGetField()`, `DMSwarmRegisterPetscDatatypeField()`, `DMCreateGlobalVector()`, `DMCreateLocalVector()` 139d3a51819SDave May @*/ 140d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmVectorDefineField(DM dm, const char fieldname[]) 141d71ae5a4SJacob Faibussowitsch { 142b5bcf523SDave May DM_Swarm *swarm = (DM_Swarm *)dm->data; 143b5bcf523SDave May PetscInt bs, n; 144b5bcf523SDave May PetscScalar *array; 145b5bcf523SDave May PetscDataType type; 146b5bcf523SDave May 147a9cbaee5SMatthew G. Knepley PetscFunctionBegin; 1480bf7c1c5SMatthew G. Knepley PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM); 1490bf7c1c5SMatthew G. Knepley if (fieldname) PetscAssertPointer(fieldname, 2); 1509566063dSJacob Faibussowitsch if (!swarm->issetup) PetscCall(DMSetUp(dm)); 1519566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &n, NULL, NULL)); 1529566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, fieldname, &bs, &type, (void **)&array)); 153b5bcf523SDave May 154b5bcf523SDave May /* Check all fields are of type PETSC_REAL or PETSC_SCALAR */ 15508401ef6SPierre Jolivet PetscCheck(type == PETSC_REAL, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only valid for PETSC_REAL"); 1569566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(swarm->vec_field_name, PETSC_MAX_PATH_LEN - 1, "%s", fieldname)); 157b5bcf523SDave May swarm->vec_field_set = PETSC_TRUE; 1581b1ea282SDave May swarm->vec_field_bs = bs; 159b5bcf523SDave May swarm->vec_field_nlocal = n; 1609566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, fieldname, &bs, &type, (void **)&array)); 1613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 162b5bcf523SDave May } 163b5bcf523SDave May 164cc651181SDave May /* requires DMSwarmDefineFieldVector has been called */ 16566976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateGlobalVector_Swarm(DM dm, Vec *vec) 166d71ae5a4SJacob Faibussowitsch { 167b5bcf523SDave May DM_Swarm *swarm = (DM_Swarm *)dm->data; 168b5bcf523SDave May Vec x; 169b5bcf523SDave May char name[PETSC_MAX_PATH_LEN]; 170b5bcf523SDave May 171a9cbaee5SMatthew G. Knepley PetscFunctionBegin; 1729566063dSJacob Faibussowitsch if (!swarm->issetup) PetscCall(DMSetUp(dm)); 17328b400f6SJacob Faibussowitsch PetscCheck(swarm->vec_field_set, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Must call DMSwarmVectorDefineField first"); 17408401ef6SPierre Jolivet PetscCheck(swarm->vec_field_nlocal == swarm->db->L, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "DMSwarm sizes have changed since last call to VectorDefineField first"); /* Stale data */ 175cc651181SDave May 1769566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN - 1, "DMSwarmField_%s", swarm->vec_field_name)); 1779566063dSJacob Faibussowitsch PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), &x)); 1789566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)x, name)); 1799566063dSJacob Faibussowitsch PetscCall(VecSetSizes(x, swarm->db->L * swarm->vec_field_bs, PETSC_DETERMINE)); 1809566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(x, swarm->vec_field_bs)); 1819566063dSJacob Faibussowitsch PetscCall(VecSetDM(x, dm)); 1829566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(x)); 1839566063dSJacob Faibussowitsch PetscCall(VecSetDM(x, dm)); 1849566063dSJacob Faibussowitsch PetscCall(VecSetOperation(x, VECOP_VIEW, (void (*)(void))VecView_Swarm)); 185b5bcf523SDave May *vec = x; 1863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 187b5bcf523SDave May } 188b5bcf523SDave May 189b5bcf523SDave May /* requires DMSwarmDefineFieldVector has been called */ 19066976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateLocalVector_Swarm(DM dm, Vec *vec) 191d71ae5a4SJacob Faibussowitsch { 192b5bcf523SDave May DM_Swarm *swarm = (DM_Swarm *)dm->data; 193b5bcf523SDave May Vec x; 194b5bcf523SDave May char name[PETSC_MAX_PATH_LEN]; 195b5bcf523SDave May 196a9cbaee5SMatthew G. Knepley PetscFunctionBegin; 1979566063dSJacob Faibussowitsch if (!swarm->issetup) PetscCall(DMSetUp(dm)); 19828b400f6SJacob Faibussowitsch PetscCheck(swarm->vec_field_set, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Must call DMSwarmVectorDefineField first"); 19908401ef6SPierre Jolivet PetscCheck(swarm->vec_field_nlocal == swarm->db->L, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "DMSwarm sizes have changed since last call to VectorDefineField first"); 200cc651181SDave May 2019566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN - 1, "DMSwarmField_%s", swarm->vec_field_name)); 2029566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, &x)); 2039566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)x, name)); 2049566063dSJacob Faibussowitsch PetscCall(VecSetSizes(x, swarm->db->L * swarm->vec_field_bs, PETSC_DETERMINE)); 2059566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(x, swarm->vec_field_bs)); 2069566063dSJacob Faibussowitsch PetscCall(VecSetDM(x, dm)); 2079566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(x)); 208b5bcf523SDave May *vec = x; 2093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 210b5bcf523SDave May } 211b5bcf523SDave May 212d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSwarmDestroyVectorFromField_Private(DM dm, const char fieldname[], Vec *vec) 213d71ae5a4SJacob Faibussowitsch { 214fb1bcc12SMatthew G. Knepley DM_Swarm *swarm = (DM_Swarm *)dm->data; 21577048351SPatrick Sanan DMSwarmDataField gfield; 2162e956fe4SStefano Zampini PetscInt bs, nlocal, fid = -1, cfid = -2; 2172e956fe4SStefano Zampini PetscBool flg; 218d3a51819SDave May 219fb1bcc12SMatthew G. Knepley PetscFunctionBegin; 2202e956fe4SStefano Zampini /* check vector is an inplace array */ 2212e956fe4SStefano Zampini PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldIdByName(swarm->db, fieldname, &fid)); 2222e956fe4SStefano Zampini PetscCall(PetscObjectComposedDataGetInt((PetscObject)*vec, SwarmDataFieldId, cfid, flg)); 223ea17275aSJose E. Roman (void)flg; /* avoid compiler warning */ 224e978a55eSPierre Jolivet PetscCheck(cfid == fid, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Vector being destroyed was not created from DMSwarm field(%s)! %" PetscInt_FMT " != %" PetscInt_FMT, fieldname, cfid, fid); 2259566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(*vec, &nlocal)); 2269566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(*vec, &bs)); 22708401ef6SPierre Jolivet PetscCheck(nlocal / bs == swarm->db->L, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "DMSwarm sizes have changed since vector was created - cannot ensure pointers are valid"); 2289566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield)); 2299566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldRestoreAccess(gfield)); 2308df78e51SMark Adams PetscCall(VecResetArray(*vec)); 2319566063dSJacob Faibussowitsch PetscCall(VecDestroy(vec)); 2323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 233fb1bcc12SMatthew G. Knepley } 234fb1bcc12SMatthew G. Knepley 235d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSwarmCreateVectorFromField_Private(DM dm, const char fieldname[], MPI_Comm comm, Vec *vec) 236d71ae5a4SJacob Faibussowitsch { 237fb1bcc12SMatthew G. Knepley DM_Swarm *swarm = (DM_Swarm *)dm->data; 238fb1bcc12SMatthew G. Knepley PetscDataType type; 239fb1bcc12SMatthew G. Knepley PetscScalar *array; 2402e956fe4SStefano Zampini PetscInt bs, n, fid; 241fb1bcc12SMatthew G. Knepley char name[PETSC_MAX_PATH_LEN]; 242e4fbd051SBarry Smith PetscMPIInt size; 2438df78e51SMark Adams PetscBool iscuda, iskokkos; 244fb1bcc12SMatthew G. Knepley 245fb1bcc12SMatthew G. Knepley PetscFunctionBegin; 2469566063dSJacob Faibussowitsch if (!swarm->issetup) PetscCall(DMSetUp(dm)); 2479566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &n, NULL, NULL)); 2489566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, fieldname, &bs, &type, (void **)&array)); 24908401ef6SPierre Jolivet PetscCheck(type == PETSC_REAL, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only valid for PETSC_REAL"); 250fb1bcc12SMatthew G. Knepley 2519566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 2528df78e51SMark Adams PetscCall(PetscStrcmp(dm->vectype, VECKOKKOS, &iskokkos)); 2538df78e51SMark Adams PetscCall(PetscStrcmp(dm->vectype, VECCUDA, &iscuda)); 2548df78e51SMark Adams PetscCall(VecCreate(comm, vec)); 2558df78e51SMark Adams PetscCall(VecSetSizes(*vec, n * bs, PETSC_DETERMINE)); 2568df78e51SMark Adams PetscCall(VecSetBlockSize(*vec, bs)); 2578df78e51SMark Adams if (iskokkos) PetscCall(VecSetType(*vec, VECKOKKOS)); 2588df78e51SMark Adams else if (iscuda) PetscCall(VecSetType(*vec, VECCUDA)); 2598df78e51SMark Adams else PetscCall(VecSetType(*vec, VECSTANDARD)); 2608df78e51SMark Adams PetscCall(VecPlaceArray(*vec, array)); 2618df78e51SMark Adams 2629566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN - 1, "DMSwarmSharedField_%s", fieldname)); 2639566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*vec, name)); 264fb1bcc12SMatthew G. Knepley 265fb1bcc12SMatthew G. Knepley /* Set guard */ 2662e956fe4SStefano Zampini PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldIdByName(swarm->db, fieldname, &fid)); 2672e956fe4SStefano Zampini PetscCall(PetscObjectComposedDataSetInt((PetscObject)*vec, SwarmDataFieldId, fid)); 26874d0cae8SMatthew G. Knepley 2699566063dSJacob Faibussowitsch PetscCall(VecSetDM(*vec, dm)); 2709566063dSJacob Faibussowitsch PetscCall(VecSetOperation(*vec, VECOP_VIEW, (void (*)(void))VecView_Swarm)); 2713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 272fb1bcc12SMatthew G. Knepley } 273fb1bcc12SMatthew G. Knepley 2740643ed39SMatthew G. Knepley /* This creates a "mass matrix" between a finite element and particle space. If a finite element interpolant is given by 2750643ed39SMatthew G. Knepley 2760643ed39SMatthew G. Knepley \hat f = \sum_i f_i \phi_i 2770643ed39SMatthew G. Knepley 2780643ed39SMatthew G. Knepley and a particle function is given by 2790643ed39SMatthew G. Knepley 2800643ed39SMatthew G. Knepley f = \sum_i w_i \delta(x - x_i) 2810643ed39SMatthew G. Knepley 2820643ed39SMatthew G. Knepley then we want to require that 2830643ed39SMatthew G. Knepley 2840643ed39SMatthew G. Knepley M \hat f = M_p f 2850643ed39SMatthew G. Knepley 2860643ed39SMatthew G. Knepley where the particle mass matrix is given by 2870643ed39SMatthew G. Knepley 2880643ed39SMatthew G. Knepley (M_p)_{ij} = \int \phi_i \delta(x - x_j) 2890643ed39SMatthew G. Knepley 2900643ed39SMatthew G. Knepley The way Dave May does particles, they amount to quadratue weights rather than delta functions, so he has |J| is in 2910643ed39SMatthew G. Knepley his integral. We allow this with the boolean flag. 2920643ed39SMatthew G. Knepley */ 293d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSwarmComputeMassMatrix_Private(DM dmc, DM dmf, Mat mass, PetscBool useDeltaFunction, void *ctx) 294d71ae5a4SJacob Faibussowitsch { 29583c47955SMatthew G. Knepley const char *name = "Mass Matrix"; 2960643ed39SMatthew G. Knepley MPI_Comm comm; 29783c47955SMatthew G. Knepley PetscDS prob; 29883c47955SMatthew G. Knepley PetscSection fsection, globalFSection; 299e8f14785SLisandro Dalcin PetscHSetIJ ht; 3000643ed39SMatthew G. Knepley PetscLayout rLayout, colLayout; 30183c47955SMatthew G. Knepley PetscInt *dnz, *onz; 302adb2528bSMark Adams PetscInt locRows, locCols, rStart, colStart, colEnd, *rowIDXs; 3030643ed39SMatthew G. Knepley PetscReal *xi, *v0, *J, *invJ, detJ = 1.0, v0ref[3] = {-1.0, -1.0, -1.0}; 30483c47955SMatthew G. Knepley PetscScalar *elemMat; 3050bf7c1c5SMatthew G. Knepley PetscInt dim, Nf, field, cStart, cEnd, cell, totDim, maxC = 0, totNc = 0; 30683c47955SMatthew G. Knepley 30783c47955SMatthew G. Knepley PetscFunctionBegin; 3089566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)mass, &comm)); 3099566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dmf, &dim)); 3109566063dSJacob Faibussowitsch PetscCall(DMGetDS(dmf, &prob)); 3119566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(prob, &Nf)); 3129566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(prob, &totDim)); 3139566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(dim, &v0, dim * dim, &J, dim * dim, &invJ)); 3149566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dmf, &fsection)); 3159566063dSJacob Faibussowitsch PetscCall(DMGetGlobalSection(dmf, &globalFSection)); 3169566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd)); 3179566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(mass, &locRows, &locCols)); 31883c47955SMatthew G. Knepley 3199566063dSJacob Faibussowitsch PetscCall(PetscLayoutCreate(comm, &colLayout)); 3209566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetLocalSize(colLayout, locCols)); 3219566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetBlockSize(colLayout, 1)); 3229566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(colLayout)); 3239566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRange(colLayout, &colStart, &colEnd)); 3249566063dSJacob Faibussowitsch PetscCall(PetscLayoutDestroy(&colLayout)); 3250643ed39SMatthew G. Knepley 3269566063dSJacob Faibussowitsch PetscCall(PetscLayoutCreate(comm, &rLayout)); 3279566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetLocalSize(rLayout, locRows)); 3289566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetBlockSize(rLayout, 1)); 3299566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(rLayout)); 3309566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRange(rLayout, &rStart, NULL)); 3319566063dSJacob Faibussowitsch PetscCall(PetscLayoutDestroy(&rLayout)); 3320643ed39SMatthew G. Knepley 3339566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(locRows, &dnz, locRows, &onz)); 3349566063dSJacob Faibussowitsch PetscCall(PetscHSetIJCreate(&ht)); 33553e60ab4SJoseph Pusztay 3369566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(comm, NULL)); 3370bf7c1c5SMatthew G. Knepley for (field = 0; field < Nf; ++field) { 3380bf7c1c5SMatthew G. Knepley PetscObject obj; 3390bf7c1c5SMatthew G. Knepley PetscClassId id; 3400bf7c1c5SMatthew G. Knepley PetscInt Nc; 3410bf7c1c5SMatthew G. Knepley 3420bf7c1c5SMatthew G. Knepley PetscCall(PetscDSGetDiscretization(prob, field, &obj)); 3430bf7c1c5SMatthew G. Knepley PetscCall(PetscObjectGetClassId(obj, &id)); 3440bf7c1c5SMatthew G. Knepley if (id == PETSCFE_CLASSID) PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc)); 3450bf7c1c5SMatthew G. Knepley else PetscCall(PetscFVGetNumComponents((PetscFV)obj, &Nc)); 3460bf7c1c5SMatthew G. Knepley totNc += Nc; 3470bf7c1c5SMatthew G. Knepley } 3480643ed39SMatthew G. Knepley /* count non-zeros */ 3499566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetAccess(dmc)); 35083c47955SMatthew G. Knepley for (field = 0; field < Nf; ++field) { 3510bf7c1c5SMatthew G. Knepley PetscObject obj; 3520bf7c1c5SMatthew G. Knepley PetscClassId id; 3530bf7c1c5SMatthew G. Knepley PetscInt Nc; 3540bf7c1c5SMatthew G. Knepley 3550bf7c1c5SMatthew G. Knepley PetscCall(PetscDSGetDiscretization(prob, field, &obj)); 3560bf7c1c5SMatthew G. Knepley PetscCall(PetscObjectGetClassId(obj, &id)); 3570bf7c1c5SMatthew G. Knepley if (id == PETSCFE_CLASSID) PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc)); 3580bf7c1c5SMatthew G. Knepley else PetscCall(PetscFVGetNumComponents((PetscFV)obj, &Nc)); 3590bf7c1c5SMatthew G. Knepley 36083c47955SMatthew G. Knepley for (cell = cStart; cell < cEnd; ++cell) { 3610643ed39SMatthew G. Knepley PetscInt *findices, *cindices; /* fine is vertices, coarse is particles */ 36283c47955SMatthew G. Knepley PetscInt numFIndices, numCIndices; 36383c47955SMatthew G. Knepley 3649566063dSJacob Faibussowitsch PetscCall(DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL)); 3659566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices)); 366fc7c92abSMatthew G. Knepley maxC = PetscMax(maxC, numCIndices); 36783c47955SMatthew G. Knepley { 368e8f14785SLisandro Dalcin PetscHashIJKey key; 369e8f14785SLisandro Dalcin PetscBool missing; 3700bf7c1c5SMatthew G. Knepley for (PetscInt i = 0; i < numFIndices; ++i) { 371adb2528bSMark Adams key.j = findices[i]; /* global column (from Plex) */ 372adb2528bSMark Adams if (key.j >= 0) { 37383c47955SMatthew G. Knepley /* Get indices for coarse elements */ 3740bf7c1c5SMatthew G. Knepley for (PetscInt j = 0; j < numCIndices; ++j) { 3750bf7c1c5SMatthew G. Knepley for (PetscInt c = 0; c < Nc; ++c) { 3760bf7c1c5SMatthew G. Knepley // TODO Need field offset on particle here 3770bf7c1c5SMatthew G. Knepley key.i = cindices[j] * totNc + c + rStart; /* global cols (from Swarm) */ 378adb2528bSMark Adams if (key.i < 0) continue; 3799566063dSJacob Faibussowitsch PetscCall(PetscHSetIJQueryAdd(ht, key, &missing)); 38083c47955SMatthew G. Knepley if (missing) { 3810643ed39SMatthew G. Knepley if ((key.j >= colStart) && (key.j < colEnd)) ++dnz[key.i - rStart]; 382e8f14785SLisandro Dalcin else ++onz[key.i - rStart]; 38363a3b9bcSJacob Faibussowitsch } else SETERRQ(PetscObjectComm((PetscObject)dmf), PETSC_ERR_SUP, "Set new value at %" PetscInt_FMT ",%" PetscInt_FMT, key.i, key.j); 38483c47955SMatthew G. Knepley } 385fc7c92abSMatthew G. Knepley } 386fc7c92abSMatthew G. Knepley } 3870bf7c1c5SMatthew G. Knepley } 3889566063dSJacob Faibussowitsch PetscCall(PetscFree(cindices)); 38983c47955SMatthew G. Knepley } 3909566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL)); 39183c47955SMatthew G. Knepley } 39283c47955SMatthew G. Knepley } 3939566063dSJacob Faibussowitsch PetscCall(PetscHSetIJDestroy(&ht)); 3949566063dSJacob Faibussowitsch PetscCall(MatXAIJSetPreallocation(mass, 1, dnz, onz, NULL, NULL)); 3959566063dSJacob Faibussowitsch PetscCall(MatSetOption(mass, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE)); 3969566063dSJacob Faibussowitsch PetscCall(PetscFree2(dnz, onz)); 3970bf7c1c5SMatthew G. Knepley PetscCall(PetscMalloc3(maxC * totNc * totDim, &elemMat, maxC * totNc, &rowIDXs, maxC * dim, &xi)); 39883c47955SMatthew G. Knepley for (field = 0; field < Nf; ++field) { 399ef0bb6c7SMatthew G. Knepley PetscTabulation Tcoarse; 40083c47955SMatthew G. Knepley PetscObject obj; 401ad9634fcSMatthew G. Knepley PetscClassId id; 402ea7032a0SMatthew G. Knepley PetscReal *fieldVals; 4030bf7c1c5SMatthew G. Knepley PetscInt Nc; 40483c47955SMatthew G. Knepley 4059566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(prob, field, &obj)); 406ad9634fcSMatthew G. Knepley PetscCall(PetscObjectGetClassId(obj, &id)); 407ad9634fcSMatthew G. Knepley if (id == PETSCFE_CLASSID) PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc)); 408ad9634fcSMatthew G. Knepley else PetscCall(PetscFVGetNumComponents((PetscFV)obj, &Nc)); 409ea7032a0SMatthew G. Knepley 410ea7032a0SMatthew G. Knepley PetscCall(DMSwarmGetField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **)&fieldVals)); 41183c47955SMatthew G. Knepley for (cell = cStart; cell < cEnd; ++cell) { 41283c47955SMatthew G. Knepley PetscInt *findices, *cindices; 41383c47955SMatthew G. Knepley PetscInt numFIndices, numCIndices; 41483c47955SMatthew G. Knepley 4150643ed39SMatthew G. Knepley /* TODO: Use DMField instead of assuming affine */ 4169566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ)); 4179566063dSJacob Faibussowitsch PetscCall(DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL)); 4189566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices)); 4190bf7c1c5SMatthew G. Knepley for (PetscInt j = 0; j < numCIndices; ++j) CoordinatesRealToRef(dim, dim, v0ref, v0, invJ, &fieldVals[cindices[j] * dim], &xi[j * dim]); 420ad9634fcSMatthew G. Knepley if (id == PETSCFE_CLASSID) PetscCall(PetscFECreateTabulation((PetscFE)obj, 1, numCIndices, xi, 0, &Tcoarse)); 421ad9634fcSMatthew G. Knepley else PetscCall(PetscFVCreateTabulation((PetscFV)obj, 1, numCIndices, xi, 0, &Tcoarse)); 42283c47955SMatthew G. Knepley /* Get elemMat entries by multiplying by weight */ 4230bf7c1c5SMatthew G. Knepley PetscCall(PetscArrayzero(elemMat, numCIndices * Nc * totDim)); 4240bf7c1c5SMatthew G. Knepley for (PetscInt i = 0; i < numFIndices / Nc; ++i) { 4250bf7c1c5SMatthew G. Knepley for (PetscInt j = 0; j < numCIndices; ++j) { 4260bf7c1c5SMatthew G. Knepley for (PetscInt c = 0; c < Nc; ++c) { 4270bf7c1c5SMatthew G. Knepley // TODO Need field offset on particle and field here 4280643ed39SMatthew G. Knepley /* B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c */ 4290bf7c1c5SMatthew G. Knepley elemMat[(j * totNc + c) * numFIndices + i * Nc + c] += Tcoarse->T[0][(j * numFIndices + i * Nc + c) * Nc + c] * (useDeltaFunction ? 1.0 : detJ); 43083c47955SMatthew G. Knepley } 4310643ed39SMatthew G. Knepley } 4320643ed39SMatthew G. Knepley } 4330bf7c1c5SMatthew G. Knepley for (PetscInt j = 0; j < numCIndices; ++j) 4340bf7c1c5SMatthew G. Knepley // TODO Need field offset on particle here 4350bf7c1c5SMatthew G. Knepley for (PetscInt c = 0; c < Nc; ++c) rowIDXs[j * Nc + c] = cindices[j] * totNc + c + rStart; 4360bf7c1c5SMatthew G. Knepley if (0) PetscCall(DMPrintCellMatrix(cell, name, numCIndices * Nc, numFIndices, elemMat)); 4370bf7c1c5SMatthew G. Knepley PetscCall(MatSetValues(mass, numCIndices * Nc, rowIDXs, numFIndices, findices, elemMat, ADD_VALUES)); 4389566063dSJacob Faibussowitsch PetscCall(PetscFree(cindices)); 4399566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL)); 4409566063dSJacob Faibussowitsch PetscCall(PetscTabulationDestroy(&Tcoarse)); 44183c47955SMatthew G. Knepley } 442ea7032a0SMatthew G. Knepley PetscCall(DMSwarmRestoreField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **)&fieldVals)); 44383c47955SMatthew G. Knepley } 4449566063dSJacob Faibussowitsch PetscCall(PetscFree3(elemMat, rowIDXs, xi)); 4459566063dSJacob Faibussowitsch PetscCall(DMSwarmSortRestoreAccess(dmc)); 4469566063dSJacob Faibussowitsch PetscCall(PetscFree3(v0, J, invJ)); 4479566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(mass, MAT_FINAL_ASSEMBLY)); 4489566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(mass, MAT_FINAL_ASSEMBLY)); 4493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 45083c47955SMatthew G. Knepley } 45183c47955SMatthew G. Knepley 452d0c080abSJoseph Pusztay /* Returns empty matrix for use with SNES FD */ 453d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCreateMatrix_Swarm(DM sw, Mat *m) 454d71ae5a4SJacob Faibussowitsch { 455d0c080abSJoseph Pusztay Vec field; 456d0c080abSJoseph Pusztay PetscInt size; 457d0c080abSJoseph Pusztay 458d0c080abSJoseph Pusztay PetscFunctionBegin; 4599566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(sw, &field)); 4609566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(field, &size)); 4619566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(sw, &field)); 4629566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, m)); 4639566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(*m)); 4649566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*m, PETSC_DECIDE, PETSC_DECIDE, size, size)); 4659566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(*m, 1, NULL)); 4669566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(*m)); 4679566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(*m, MAT_FINAL_ASSEMBLY)); 4689566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(*m, MAT_FINAL_ASSEMBLY)); 4699566063dSJacob Faibussowitsch PetscCall(MatShift(*m, 1.0)); 4709566063dSJacob Faibussowitsch PetscCall(MatSetDM(*m, sw)); 4713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 472d0c080abSJoseph Pusztay } 473d0c080abSJoseph Pusztay 474adb2528bSMark Adams /* FEM cols, Particle rows */ 475d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCreateMassMatrix_Swarm(DM dmCoarse, DM dmFine, Mat *mass) 476d71ae5a4SJacob Faibussowitsch { 477895a1698SMatthew G. Knepley PetscSection gsf; 4780bf7c1c5SMatthew G. Knepley PetscInt m, n, Np, bs = 1; 47983c47955SMatthew G. Knepley void *ctx; 4800bf7c1c5SMatthew G. Knepley PetscBool set = ((DM_Swarm *)dmCoarse->data)->vec_field_set; 4810bf7c1c5SMatthew G. Knepley char *name = ((DM_Swarm *)dmCoarse->data)->vec_field_name; 48283c47955SMatthew G. Knepley 48383c47955SMatthew G. Knepley PetscFunctionBegin; 4849566063dSJacob Faibussowitsch PetscCall(DMGetGlobalSection(dmFine, &gsf)); 4859566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstrainedStorageSize(gsf, &m)); 4860bf7c1c5SMatthew G. Knepley PetscCall(DMSwarmGetLocalSize(dmCoarse, &Np)); 4870bf7c1c5SMatthew G. Knepley // TODO Include all fields 4880bf7c1c5SMatthew G. Knepley if (set) PetscCall(DMSwarmGetFieldInfo(dmCoarse, name, &bs, NULL)); 4890bf7c1c5SMatthew G. Knepley n = Np * bs; 4909566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)dmCoarse), mass)); 4919566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*mass, n, m, PETSC_DETERMINE, PETSC_DETERMINE)); 4929566063dSJacob Faibussowitsch PetscCall(MatSetType(*mass, dmCoarse->mattype)); 4939566063dSJacob Faibussowitsch PetscCall(DMGetApplicationContext(dmFine, &ctx)); 49483c47955SMatthew G. Knepley 4959566063dSJacob Faibussowitsch PetscCall(DMSwarmComputeMassMatrix_Private(dmCoarse, dmFine, *mass, PETSC_TRUE, ctx)); 4969566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(*mass, NULL, "-mass_mat_view")); 4973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 49883c47955SMatthew G. Knepley } 49983c47955SMatthew G. Knepley 500d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSwarmComputeMassMatrixSquare_Private(DM dmc, DM dmf, Mat mass, PetscBool useDeltaFunction, void *ctx) 501d71ae5a4SJacob Faibussowitsch { 5024cc7f7b2SMatthew G. Knepley const char *name = "Mass Matrix Square"; 5034cc7f7b2SMatthew G. Knepley MPI_Comm comm; 5044cc7f7b2SMatthew G. Knepley PetscDS prob; 5054cc7f7b2SMatthew G. Knepley PetscSection fsection, globalFSection; 5064cc7f7b2SMatthew G. Knepley PetscHSetIJ ht; 5074cc7f7b2SMatthew G. Knepley PetscLayout rLayout, colLayout; 5084cc7f7b2SMatthew G. Knepley PetscInt *dnz, *onz, *adj, depth, maxConeSize, maxSupportSize, maxAdjSize; 5094cc7f7b2SMatthew G. Knepley PetscInt locRows, locCols, rStart, colStart, colEnd, *rowIDXs; 5104cc7f7b2SMatthew G. Knepley PetscReal *xi, *v0, *J, *invJ, detJ = 1.0, v0ref[3] = {-1.0, -1.0, -1.0}; 5114cc7f7b2SMatthew G. Knepley PetscScalar *elemMat, *elemMatSq; 5124cc7f7b2SMatthew G. Knepley PetscInt cdim, Nf, field, cStart, cEnd, cell, totDim, maxC = 0; 5134cc7f7b2SMatthew G. Knepley 5144cc7f7b2SMatthew G. Knepley PetscFunctionBegin; 5159566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)mass, &comm)); 5169566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dmf, &cdim)); 5179566063dSJacob Faibussowitsch PetscCall(DMGetDS(dmf, &prob)); 5189566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(prob, &Nf)); 5199566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(prob, &totDim)); 5209566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(cdim, &v0, cdim * cdim, &J, cdim * cdim, &invJ)); 5219566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dmf, &fsection)); 5229566063dSJacob Faibussowitsch PetscCall(DMGetGlobalSection(dmf, &globalFSection)); 5239566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd)); 5249566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(mass, &locRows, &locCols)); 5254cc7f7b2SMatthew G. Knepley 5269566063dSJacob Faibussowitsch PetscCall(PetscLayoutCreate(comm, &colLayout)); 5279566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetLocalSize(colLayout, locCols)); 5289566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetBlockSize(colLayout, 1)); 5299566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(colLayout)); 5309566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRange(colLayout, &colStart, &colEnd)); 5319566063dSJacob Faibussowitsch PetscCall(PetscLayoutDestroy(&colLayout)); 5324cc7f7b2SMatthew G. Knepley 5339566063dSJacob Faibussowitsch PetscCall(PetscLayoutCreate(comm, &rLayout)); 5349566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetLocalSize(rLayout, locRows)); 5359566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetBlockSize(rLayout, 1)); 5369566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(rLayout)); 5379566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRange(rLayout, &rStart, NULL)); 5389566063dSJacob Faibussowitsch PetscCall(PetscLayoutDestroy(&rLayout)); 5394cc7f7b2SMatthew G. Knepley 5409566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepth(dmf, &depth)); 5419566063dSJacob Faibussowitsch PetscCall(DMPlexGetMaxSizes(dmf, &maxConeSize, &maxSupportSize)); 5424cc7f7b2SMatthew G. Knepley maxAdjSize = PetscPowInt(maxConeSize * maxSupportSize, depth); 5439566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(maxAdjSize, &adj)); 5444cc7f7b2SMatthew G. Knepley 5459566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(locRows, &dnz, locRows, &onz)); 5469566063dSJacob Faibussowitsch PetscCall(PetscHSetIJCreate(&ht)); 5474cc7f7b2SMatthew G. Knepley /* Count nonzeros 5484cc7f7b2SMatthew G. Knepley This is just FVM++, but we cannot use the Plex P0 allocation since unknowns in a cell will not be contiguous 5494cc7f7b2SMatthew G. Knepley */ 5509566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetAccess(dmc)); 5514cc7f7b2SMatthew G. Knepley for (cell = cStart; cell < cEnd; ++cell) { 5524cc7f7b2SMatthew G. Knepley PetscInt i; 5534cc7f7b2SMatthew G. Knepley PetscInt *cindices; 5544cc7f7b2SMatthew G. Knepley PetscInt numCIndices; 5554cc7f7b2SMatthew G. Knepley #if 0 5564cc7f7b2SMatthew G. Knepley PetscInt adjSize = maxAdjSize, a, j; 5574cc7f7b2SMatthew G. Knepley #endif 5584cc7f7b2SMatthew G. Knepley 5599566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices)); 5604cc7f7b2SMatthew G. Knepley maxC = PetscMax(maxC, numCIndices); 5614cc7f7b2SMatthew G. Knepley /* Diagonal block */ 562ad540459SPierre Jolivet for (i = 0; i < numCIndices; ++i) dnz[cindices[i]] += numCIndices; 5634cc7f7b2SMatthew G. Knepley #if 0 5644cc7f7b2SMatthew G. Knepley /* Off-diagonal blocks */ 5659566063dSJacob Faibussowitsch PetscCall(DMPlexGetAdjacency(dmf, cell, &adjSize, &adj)); 5664cc7f7b2SMatthew G. Knepley for (a = 0; a < adjSize; ++a) { 5674cc7f7b2SMatthew G. Knepley if (adj[a] >= cStart && adj[a] < cEnd && adj[a] != cell) { 5684cc7f7b2SMatthew G. Knepley const PetscInt ncell = adj[a]; 5694cc7f7b2SMatthew G. Knepley PetscInt *ncindices; 5704cc7f7b2SMatthew G. Knepley PetscInt numNCIndices; 5714cc7f7b2SMatthew G. Knepley 5729566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetPointsPerCell(dmc, ncell, &numNCIndices, &ncindices)); 5734cc7f7b2SMatthew G. Knepley { 5744cc7f7b2SMatthew G. Knepley PetscHashIJKey key; 5754cc7f7b2SMatthew G. Knepley PetscBool missing; 5764cc7f7b2SMatthew G. Knepley 5774cc7f7b2SMatthew G. Knepley for (i = 0; i < numCIndices; ++i) { 5784cc7f7b2SMatthew G. Knepley key.i = cindices[i] + rStart; /* global rows (from Swarm) */ 5794cc7f7b2SMatthew G. Knepley if (key.i < 0) continue; 5804cc7f7b2SMatthew G. Knepley for (j = 0; j < numNCIndices; ++j) { 5814cc7f7b2SMatthew G. Knepley key.j = ncindices[j] + rStart; /* global column (from Swarm) */ 5824cc7f7b2SMatthew G. Knepley if (key.j < 0) continue; 5839566063dSJacob Faibussowitsch PetscCall(PetscHSetIJQueryAdd(ht, key, &missing)); 5844cc7f7b2SMatthew G. Knepley if (missing) { 5854cc7f7b2SMatthew G. Knepley if ((key.j >= colStart) && (key.j < colEnd)) ++dnz[key.i - rStart]; 5864cc7f7b2SMatthew G. Knepley else ++onz[key.i - rStart]; 5874cc7f7b2SMatthew G. Knepley } 5884cc7f7b2SMatthew G. Knepley } 5894cc7f7b2SMatthew G. Knepley } 5904cc7f7b2SMatthew G. Knepley } 5919566063dSJacob Faibussowitsch PetscCall(PetscFree(ncindices)); 5924cc7f7b2SMatthew G. Knepley } 5934cc7f7b2SMatthew G. Knepley } 5944cc7f7b2SMatthew G. Knepley #endif 5959566063dSJacob Faibussowitsch PetscCall(PetscFree(cindices)); 5964cc7f7b2SMatthew G. Knepley } 5979566063dSJacob Faibussowitsch PetscCall(PetscHSetIJDestroy(&ht)); 5989566063dSJacob Faibussowitsch PetscCall(MatXAIJSetPreallocation(mass, 1, dnz, onz, NULL, NULL)); 5999566063dSJacob Faibussowitsch PetscCall(MatSetOption(mass, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE)); 6009566063dSJacob Faibussowitsch PetscCall(PetscFree2(dnz, onz)); 6019566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(maxC * totDim, &elemMat, maxC * maxC, &elemMatSq, maxC, &rowIDXs, maxC * cdim, &xi)); 6024cc7f7b2SMatthew G. Knepley /* Fill in values 6034cc7f7b2SMatthew G. Knepley Each entry is a sum of terms \phi_i(x_p) \phi_i(x_q) 6044cc7f7b2SMatthew G. Knepley Start just by producing block diagonal 6054cc7f7b2SMatthew G. Knepley Could loop over adjacent cells 6064cc7f7b2SMatthew G. Knepley Produce neighboring element matrix 6074cc7f7b2SMatthew G. Knepley TODO Determine which columns and rows correspond to shared dual vector 6084cc7f7b2SMatthew G. Knepley Do MatMatMult with rectangular matrices 6094cc7f7b2SMatthew G. Knepley Insert block 6104cc7f7b2SMatthew G. Knepley */ 6114cc7f7b2SMatthew G. Knepley for (field = 0; field < Nf; ++field) { 6124cc7f7b2SMatthew G. Knepley PetscTabulation Tcoarse; 6134cc7f7b2SMatthew G. Knepley PetscObject obj; 6144cc7f7b2SMatthew G. Knepley PetscReal *coords; 6154cc7f7b2SMatthew G. Knepley PetscInt Nc, i; 6164cc7f7b2SMatthew G. Knepley 6179566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(prob, field, &obj)); 6189566063dSJacob Faibussowitsch PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc)); 61963a3b9bcSJacob Faibussowitsch PetscCheck(Nc == 1, PetscObjectComm((PetscObject)dmf), PETSC_ERR_SUP, "Can only interpolate a scalar field from particles, Nc = %" PetscInt_FMT, Nc); 6209566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 6214cc7f7b2SMatthew G. Knepley for (cell = cStart; cell < cEnd; ++cell) { 6224cc7f7b2SMatthew G. Knepley PetscInt *findices, *cindices; 6234cc7f7b2SMatthew G. Knepley PetscInt numFIndices, numCIndices; 6244cc7f7b2SMatthew G. Knepley PetscInt p, c; 6254cc7f7b2SMatthew G. Knepley 6264cc7f7b2SMatthew G. Knepley /* TODO: Use DMField instead of assuming affine */ 6279566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ)); 6289566063dSJacob Faibussowitsch PetscCall(DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL)); 6299566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices)); 630ad540459SPierre Jolivet for (p = 0; p < numCIndices; ++p) CoordinatesRealToRef(cdim, cdim, v0ref, v0, invJ, &coords[cindices[p] * cdim], &xi[p * cdim]); 6319566063dSJacob Faibussowitsch PetscCall(PetscFECreateTabulation((PetscFE)obj, 1, numCIndices, xi, 0, &Tcoarse)); 6324cc7f7b2SMatthew G. Knepley /* Get elemMat entries by multiplying by weight */ 6339566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(elemMat, numCIndices * totDim)); 6344cc7f7b2SMatthew G. Knepley for (i = 0; i < numFIndices; ++i) { 6354cc7f7b2SMatthew G. Knepley for (p = 0; p < numCIndices; ++p) { 6364cc7f7b2SMatthew G. Knepley for (c = 0; c < Nc; ++c) { 6374cc7f7b2SMatthew G. Knepley /* B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c */ 6384cc7f7b2SMatthew G. Knepley elemMat[p * numFIndices + i] += Tcoarse->T[0][(p * numFIndices + i) * Nc + c] * (useDeltaFunction ? 1.0 : detJ); 6394cc7f7b2SMatthew G. Knepley } 6404cc7f7b2SMatthew G. Knepley } 6414cc7f7b2SMatthew G. Knepley } 6429566063dSJacob Faibussowitsch PetscCall(PetscTabulationDestroy(&Tcoarse)); 6434cc7f7b2SMatthew G. Knepley for (p = 0; p < numCIndices; ++p) rowIDXs[p] = cindices[p] + rStart; 6449566063dSJacob Faibussowitsch if (0) PetscCall(DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat)); 6454cc7f7b2SMatthew G. Knepley /* Block diagonal */ 64678884ca7SMatthew G. Knepley if (numCIndices) { 6474cc7f7b2SMatthew G. Knepley PetscBLASInt blasn, blask; 6484cc7f7b2SMatthew G. Knepley PetscScalar one = 1.0, zero = 0.0; 6494cc7f7b2SMatthew G. Knepley 6509566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(numCIndices, &blasn)); 6519566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(numFIndices, &blask)); 652792fecdfSBarry Smith PetscCallBLAS("BLASgemm", BLASgemm_("T", "N", &blasn, &blasn, &blask, &one, elemMat, &blask, elemMat, &blask, &zero, elemMatSq, &blasn)); 6534cc7f7b2SMatthew G. Knepley } 6549566063dSJacob Faibussowitsch PetscCall(MatSetValues(mass, numCIndices, rowIDXs, numCIndices, rowIDXs, elemMatSq, ADD_VALUES)); 6554cf0e950SBarry Smith /* TODO off-diagonal */ 6569566063dSJacob Faibussowitsch PetscCall(PetscFree(cindices)); 6579566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL)); 6584cc7f7b2SMatthew G. Knepley } 6599566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 6604cc7f7b2SMatthew G. Knepley } 6619566063dSJacob Faibussowitsch PetscCall(PetscFree4(elemMat, elemMatSq, rowIDXs, xi)); 6629566063dSJacob Faibussowitsch PetscCall(PetscFree(adj)); 6639566063dSJacob Faibussowitsch PetscCall(DMSwarmSortRestoreAccess(dmc)); 6649566063dSJacob Faibussowitsch PetscCall(PetscFree3(v0, J, invJ)); 6659566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(mass, MAT_FINAL_ASSEMBLY)); 6669566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(mass, MAT_FINAL_ASSEMBLY)); 6673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6684cc7f7b2SMatthew G. Knepley } 6694cc7f7b2SMatthew G. Knepley 670cc4c1da9SBarry Smith /*@ 6714cc7f7b2SMatthew G. Knepley DMSwarmCreateMassMatrixSquare - Creates the block-diagonal of the square, M^T_p M_p, of the particle mass matrix M_p 6724cc7f7b2SMatthew G. Knepley 67320f4b53cSBarry Smith Collective 6744cc7f7b2SMatthew G. Knepley 67560225df5SJacob Faibussowitsch Input Parameters: 67620f4b53cSBarry Smith + dmCoarse - a `DMSWARM` 67720f4b53cSBarry Smith - dmFine - a `DMPLEX` 6784cc7f7b2SMatthew G. Knepley 67960225df5SJacob Faibussowitsch Output Parameter: 6804cc7f7b2SMatthew G. Knepley . mass - the square of the particle mass matrix 6814cc7f7b2SMatthew G. Knepley 6824cc7f7b2SMatthew G. Knepley Level: advanced 6834cc7f7b2SMatthew G. Knepley 68420f4b53cSBarry Smith Note: 6854cc7f7b2SMatthew G. Knepley We only compute the block diagonal since this provides a good preconditioner and is completely local. It would be possible in the 6864cc7f7b2SMatthew G. Knepley future to compute the full normal equations. 6874cc7f7b2SMatthew G. Knepley 68820f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMCreateMassMatrix()` 6894cc7f7b2SMatthew G. Knepley @*/ 690d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmCreateMassMatrixSquare(DM dmCoarse, DM dmFine, Mat *mass) 691d71ae5a4SJacob Faibussowitsch { 6924cc7f7b2SMatthew G. Knepley PetscInt n; 6934cc7f7b2SMatthew G. Knepley void *ctx; 6944cc7f7b2SMatthew G. Knepley 6954cc7f7b2SMatthew G. Knepley PetscFunctionBegin; 6969566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(dmCoarse, &n)); 6979566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)dmCoarse), mass)); 6989566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*mass, n, n, PETSC_DETERMINE, PETSC_DETERMINE)); 6999566063dSJacob Faibussowitsch PetscCall(MatSetType(*mass, dmCoarse->mattype)); 7009566063dSJacob Faibussowitsch PetscCall(DMGetApplicationContext(dmFine, &ctx)); 7014cc7f7b2SMatthew G. Knepley 7029566063dSJacob Faibussowitsch PetscCall(DMSwarmComputeMassMatrixSquare_Private(dmCoarse, dmFine, *mass, PETSC_TRUE, ctx)); 7039566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(*mass, NULL, "-mass_sq_mat_view")); 7043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7054cc7f7b2SMatthew G. Knepley } 7064cc7f7b2SMatthew G. Knepley 707cc4c1da9SBarry Smith /*@ 70820f4b53cSBarry Smith DMSwarmCreateGlobalVectorFromField - Creates a `Vec` object sharing the array associated with a given field 709d3a51819SDave May 71020f4b53cSBarry Smith Collective 711d3a51819SDave May 71260225df5SJacob Faibussowitsch Input Parameters: 71320f4b53cSBarry Smith + dm - a `DMSWARM` 71462741f57SDave May - fieldname - the textual name given to a registered field 715d3a51819SDave May 71660225df5SJacob Faibussowitsch Output Parameter: 717d3a51819SDave May . vec - the vector 718d3a51819SDave May 719d3a51819SDave May Level: beginner 720d3a51819SDave May 721cc4c1da9SBarry Smith Note: 72220f4b53cSBarry Smith The vector must be returned using a matching call to `DMSwarmDestroyGlobalVectorFromField()`. 7238b8a3813SDave May 72420f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmDestroyGlobalVectorFromField()` 725d3a51819SDave May @*/ 726d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmCreateGlobalVectorFromField(DM dm, const char fieldname[], Vec *vec) 727d71ae5a4SJacob Faibussowitsch { 728fb1bcc12SMatthew G. Knepley MPI_Comm comm = PetscObjectComm((PetscObject)dm); 729b5bcf523SDave May 730fb1bcc12SMatthew G. Knepley PetscFunctionBegin; 731ea7032a0SMatthew G. Knepley PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM); 7329566063dSJacob Faibussowitsch PetscCall(DMSwarmCreateVectorFromField_Private(dm, fieldname, comm, vec)); 7333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 734b5bcf523SDave May } 735b5bcf523SDave May 736cc4c1da9SBarry Smith /*@ 73720f4b53cSBarry Smith DMSwarmDestroyGlobalVectorFromField - Destroys the `Vec` object which share the array associated with a given field 738d3a51819SDave May 73920f4b53cSBarry Smith Collective 740d3a51819SDave May 74160225df5SJacob Faibussowitsch Input Parameters: 74220f4b53cSBarry Smith + dm - a `DMSWARM` 74362741f57SDave May - fieldname - the textual name given to a registered field 744d3a51819SDave May 74560225df5SJacob Faibussowitsch Output Parameter: 746d3a51819SDave May . vec - the vector 747d3a51819SDave May 748d3a51819SDave May Level: beginner 749d3a51819SDave May 75020f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmCreateGlobalVectorFromField()` 751d3a51819SDave May @*/ 752d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDestroyGlobalVectorFromField(DM dm, const char fieldname[], Vec *vec) 753d71ae5a4SJacob Faibussowitsch { 754fb1bcc12SMatthew G. Knepley PetscFunctionBegin; 755ea7032a0SMatthew G. Knepley PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM); 7569566063dSJacob Faibussowitsch PetscCall(DMSwarmDestroyVectorFromField_Private(dm, fieldname, vec)); 7573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 758b5bcf523SDave May } 759b5bcf523SDave May 760cc4c1da9SBarry Smith /*@ 76120f4b53cSBarry Smith DMSwarmCreateLocalVectorFromField - Creates a `Vec` object sharing the array associated with a given field 762fb1bcc12SMatthew G. Knepley 76320f4b53cSBarry Smith Collective 764fb1bcc12SMatthew G. Knepley 76560225df5SJacob Faibussowitsch Input Parameters: 76620f4b53cSBarry Smith + dm - a `DMSWARM` 76762741f57SDave May - fieldname - the textual name given to a registered field 768fb1bcc12SMatthew G. Knepley 76960225df5SJacob Faibussowitsch Output Parameter: 770fb1bcc12SMatthew G. Knepley . vec - the vector 771fb1bcc12SMatthew G. Knepley 772fb1bcc12SMatthew G. Knepley Level: beginner 773fb1bcc12SMatthew G. Knepley 77420f4b53cSBarry Smith Note: 7758b8a3813SDave May The vector must be returned using a matching call to DMSwarmDestroyLocalVectorFromField(). 7768b8a3813SDave May 77720f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmDestroyLocalVectorFromField()` 778fb1bcc12SMatthew G. Knepley @*/ 779d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmCreateLocalVectorFromField(DM dm, const char fieldname[], Vec *vec) 780d71ae5a4SJacob Faibussowitsch { 781fb1bcc12SMatthew G. Knepley MPI_Comm comm = PETSC_COMM_SELF; 782bbe8250bSMatthew G. Knepley 783fb1bcc12SMatthew G. Knepley PetscFunctionBegin; 7849566063dSJacob Faibussowitsch PetscCall(DMSwarmCreateVectorFromField_Private(dm, fieldname, comm, vec)); 7853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 786bbe8250bSMatthew G. Knepley } 787fb1bcc12SMatthew G. Knepley 788cc4c1da9SBarry Smith /*@ 78920f4b53cSBarry Smith DMSwarmDestroyLocalVectorFromField - Destroys the `Vec` object which share the array associated with a given field 790fb1bcc12SMatthew G. Knepley 79120f4b53cSBarry Smith Collective 792fb1bcc12SMatthew G. Knepley 79360225df5SJacob Faibussowitsch Input Parameters: 79420f4b53cSBarry Smith + dm - a `DMSWARM` 79562741f57SDave May - fieldname - the textual name given to a registered field 796fb1bcc12SMatthew G. Knepley 79760225df5SJacob Faibussowitsch Output Parameter: 798fb1bcc12SMatthew G. Knepley . vec - the vector 799fb1bcc12SMatthew G. Knepley 800fb1bcc12SMatthew G. Knepley Level: beginner 801fb1bcc12SMatthew G. Knepley 80220f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmCreateLocalVectorFromField()` 803fb1bcc12SMatthew G. Knepley @*/ 804d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDestroyLocalVectorFromField(DM dm, const char fieldname[], Vec *vec) 805d71ae5a4SJacob Faibussowitsch { 806fb1bcc12SMatthew G. Knepley PetscFunctionBegin; 807ea7032a0SMatthew G. Knepley PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM); 8089566063dSJacob Faibussowitsch PetscCall(DMSwarmDestroyVectorFromField_Private(dm, fieldname, vec)); 8093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 810bbe8250bSMatthew G. Knepley } 811bbe8250bSMatthew G. Knepley 812cc4c1da9SBarry Smith /*@ 81320f4b53cSBarry Smith DMSwarmInitializeFieldRegister - Initiates the registration of fields to a `DMSWARM` 814d3a51819SDave May 81520f4b53cSBarry Smith Collective 816d3a51819SDave May 81760225df5SJacob Faibussowitsch Input Parameter: 81820f4b53cSBarry Smith . dm - a `DMSWARM` 819d3a51819SDave May 820d3a51819SDave May Level: beginner 821d3a51819SDave May 82220f4b53cSBarry Smith Note: 82320f4b53cSBarry Smith After all fields have been registered, you must call `DMSwarmFinalizeFieldRegister()`. 824d3a51819SDave May 82520f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmFinalizeFieldRegister()`, `DMSwarmRegisterPetscDatatypeField()`, 826db781477SPatrick Sanan `DMSwarmRegisterUserStructField()`, `DMSwarmRegisterUserDatatypeField()` 827d3a51819SDave May @*/ 828d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmInitializeFieldRegister(DM dm) 829d71ae5a4SJacob Faibussowitsch { 8305f50eb2eSDave May DM_Swarm *swarm = (DM_Swarm *)dm->data; 8313454631fSDave May 832521f74f9SMatthew G. Knepley PetscFunctionBegin; 833cc651181SDave May if (!swarm->field_registration_initialized) { 8345f50eb2eSDave May swarm->field_registration_initialized = PETSC_TRUE; 835da81f932SPierre Jolivet PetscCall(DMSwarmRegisterPetscDatatypeField(dm, DMSwarmField_pid, 1, PETSC_INT64)); /* unique identifier */ 8369566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(dm, DMSwarmField_rank, 1, PETSC_INT)); /* used for communication */ 837cc651181SDave May } 8383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8395f50eb2eSDave May } 8405f50eb2eSDave May 84174d0cae8SMatthew G. Knepley /*@ 84220f4b53cSBarry Smith DMSwarmFinalizeFieldRegister - Finalizes the registration of fields to a `DMSWARM` 843d3a51819SDave May 84420f4b53cSBarry Smith Collective 845d3a51819SDave May 84660225df5SJacob Faibussowitsch Input Parameter: 84720f4b53cSBarry Smith . dm - a `DMSWARM` 848d3a51819SDave May 849d3a51819SDave May Level: beginner 850d3a51819SDave May 85120f4b53cSBarry Smith Note: 85220f4b53cSBarry Smith After `DMSwarmFinalizeFieldRegister()` has been called, no new fields can be defined on the `DMSWARM`. 853d3a51819SDave May 85420f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmInitializeFieldRegister()`, `DMSwarmRegisterPetscDatatypeField()`, 855db781477SPatrick Sanan `DMSwarmRegisterUserStructField()`, `DMSwarmRegisterUserDatatypeField()` 856d3a51819SDave May @*/ 857d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmFinalizeFieldRegister(DM dm) 858d71ae5a4SJacob Faibussowitsch { 8595f50eb2eSDave May DM_Swarm *swarm = (DM_Swarm *)dm->data; 8606845f8f5SDave May 861521f74f9SMatthew G. Knepley PetscFunctionBegin; 86248a46eb9SPierre Jolivet if (!swarm->field_registration_finalized) PetscCall(DMSwarmDataBucketFinalize(swarm->db)); 863f0cdbbbaSDave May swarm->field_registration_finalized = PETSC_TRUE; 8643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8655f50eb2eSDave May } 8665f50eb2eSDave May 86774d0cae8SMatthew G. Knepley /*@ 86820f4b53cSBarry Smith DMSwarmSetLocalSizes - Sets the length of all registered fields on the `DMSWARM` 869d3a51819SDave May 87020f4b53cSBarry Smith Not Collective 871d3a51819SDave May 87260225df5SJacob Faibussowitsch Input Parameters: 87320f4b53cSBarry Smith + dm - a `DMSWARM` 874d3a51819SDave May . nlocal - the length of each registered field 87562741f57SDave May - buffer - the length of the buffer used to efficient dynamic re-sizing 876d3a51819SDave May 877d3a51819SDave May Level: beginner 878d3a51819SDave May 87920f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmGetLocalSize()` 880d3a51819SDave May @*/ 881d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmSetLocalSizes(DM dm, PetscInt nlocal, PetscInt buffer) 882d71ae5a4SJacob Faibussowitsch { 8835f50eb2eSDave May DM_Swarm *swarm = (DM_Swarm *)dm->data; 8845f50eb2eSDave May 885521f74f9SMatthew G. Knepley PetscFunctionBegin; 8869566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(DMSWARM_SetSizes, 0, 0, 0, 0)); 8879566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketSetSizes(swarm->db, nlocal, buffer)); 8889566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(DMSWARM_SetSizes, 0, 0, 0, 0)); 8893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8905f50eb2eSDave May } 8915f50eb2eSDave May 89274d0cae8SMatthew G. Knepley /*@ 89320f4b53cSBarry Smith DMSwarmSetCellDM - Attaches a `DM` to a `DMSWARM` 894d3a51819SDave May 89520f4b53cSBarry Smith Collective 896d3a51819SDave May 89760225df5SJacob Faibussowitsch Input Parameters: 89820f4b53cSBarry Smith + dm - a `DMSWARM` 89920f4b53cSBarry Smith - dmcell - the `DM` to attach to the `DMSWARM` 900d3a51819SDave May 901d3a51819SDave May Level: beginner 902d3a51819SDave May 90320f4b53cSBarry Smith Note: 90420f4b53cSBarry Smith The attached `DM` (dmcell) will be queried for point location and 90520f4b53cSBarry Smith neighbor MPI-rank information if `DMSwarmMigrate()` is called. 906d3a51819SDave May 90720f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmSetType()`, `DMSwarmGetCellDM()`, `DMSwarmMigrate()` 908d3a51819SDave May @*/ 909d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmSetCellDM(DM dm, DM dmcell) 910d71ae5a4SJacob Faibussowitsch { 911b16650c8SDave May DM_Swarm *swarm = (DM_Swarm *)dm->data; 912521f74f9SMatthew G. Knepley 913521f74f9SMatthew G. Knepley PetscFunctionBegin; 914b16650c8SDave May swarm->dmcell = dmcell; 9153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 916b16650c8SDave May } 917b16650c8SDave May 91874d0cae8SMatthew G. Knepley /*@ 91920f4b53cSBarry Smith DMSwarmGetCellDM - Fetches the attached cell `DM` 920d3a51819SDave May 92120f4b53cSBarry Smith Collective 922d3a51819SDave May 92360225df5SJacob Faibussowitsch Input Parameter: 92420f4b53cSBarry Smith . dm - a `DMSWARM` 925d3a51819SDave May 92660225df5SJacob Faibussowitsch Output Parameter: 92720f4b53cSBarry Smith . dmcell - the `DM` which was attached to the `DMSWARM` 928d3a51819SDave May 929d3a51819SDave May Level: beginner 930d3a51819SDave May 93120f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmSetCellDM()` 932d3a51819SDave May @*/ 933d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmGetCellDM(DM dm, DM *dmcell) 934d71ae5a4SJacob Faibussowitsch { 935fe39f135SDave May DM_Swarm *swarm = (DM_Swarm *)dm->data; 936521f74f9SMatthew G. Knepley 937521f74f9SMatthew G. Knepley PetscFunctionBegin; 938fe39f135SDave May *dmcell = swarm->dmcell; 9393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 940fe39f135SDave May } 941fe39f135SDave May 94274d0cae8SMatthew G. Knepley /*@ 943d5b43468SJose E. Roman DMSwarmGetLocalSize - Retrieves the local length of fields registered 944d3a51819SDave May 94520f4b53cSBarry Smith Not Collective 946d3a51819SDave May 94760225df5SJacob Faibussowitsch Input Parameter: 94820f4b53cSBarry Smith . dm - a `DMSWARM` 949d3a51819SDave May 95060225df5SJacob Faibussowitsch Output Parameter: 951d3a51819SDave May . nlocal - the length of each registered field 952d3a51819SDave May 953d3a51819SDave May Level: beginner 954d3a51819SDave May 95520f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmGetSize()`, `DMSwarmSetLocalSizes()` 956d3a51819SDave May @*/ 957d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmGetLocalSize(DM dm, PetscInt *nlocal) 958d71ae5a4SJacob Faibussowitsch { 959dcf43ee8SDave May DM_Swarm *swarm = (DM_Swarm *)dm->data; 960dcf43ee8SDave May 961521f74f9SMatthew G. Knepley PetscFunctionBegin; 9629566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(swarm->db, nlocal, NULL, NULL)); 9633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 964dcf43ee8SDave May } 965dcf43ee8SDave May 96674d0cae8SMatthew G. Knepley /*@ 967d5b43468SJose E. Roman DMSwarmGetSize - Retrieves the total length of fields registered 968d3a51819SDave May 96920f4b53cSBarry Smith Collective 970d3a51819SDave May 97160225df5SJacob Faibussowitsch Input Parameter: 97220f4b53cSBarry Smith . dm - a `DMSWARM` 973d3a51819SDave May 97460225df5SJacob Faibussowitsch Output Parameter: 975d3a51819SDave May . n - the total length of each registered field 976d3a51819SDave May 977d3a51819SDave May Level: beginner 978d3a51819SDave May 979d3a51819SDave May Note: 98020f4b53cSBarry Smith This calls `MPI_Allreduce()` upon each call (inefficient but safe) 981d3a51819SDave May 98220f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmGetLocalSize()`, `DMSwarmSetLocalSizes()` 983d3a51819SDave May @*/ 984d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmGetSize(DM dm, PetscInt *n) 985d71ae5a4SJacob Faibussowitsch { 986dcf43ee8SDave May DM_Swarm *swarm = (DM_Swarm *)dm->data; 9875627991aSBarry Smith PetscInt nlocal; 988dcf43ee8SDave May 989521f74f9SMatthew G. Knepley PetscFunctionBegin; 9909566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &nlocal, NULL, NULL)); 991462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(&nlocal, n, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm))); 9923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 993dcf43ee8SDave May } 994dcf43ee8SDave May 995cc4c1da9SBarry Smith /*@ 99620f4b53cSBarry Smith DMSwarmRegisterPetscDatatypeField - Register a field to a `DMSWARM` with a native PETSc data type 997d3a51819SDave May 99820f4b53cSBarry Smith Collective 999d3a51819SDave May 100060225df5SJacob Faibussowitsch Input Parameters: 100120f4b53cSBarry Smith + dm - a `DMSWARM` 1002d3a51819SDave May . fieldname - the textual name to identify this field 1003d3a51819SDave May . blocksize - the number of each data type 100420f4b53cSBarry Smith - type - a valid PETSc data type (`PETSC_CHAR`, `PETSC_SHORT`, `PETSC_INT`, `PETSC_FLOAT`, `PETSC_REAL`, `PETSC_LONG`) 1005d3a51819SDave May 1006d3a51819SDave May Level: beginner 1007d3a51819SDave May 1008d3a51819SDave May Notes: 10098b8a3813SDave May The textual name for each registered field must be unique. 1010d3a51819SDave May 101120f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterUserStructField()`, `DMSwarmRegisterUserDatatypeField()` 1012d3a51819SDave May @*/ 1013d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmRegisterPetscDatatypeField(DM dm, const char fieldname[], PetscInt blocksize, PetscDataType type) 1014d71ae5a4SJacob Faibussowitsch { 1015b62e03f8SDave May DM_Swarm *swarm = (DM_Swarm *)dm->data; 1016b62e03f8SDave May size_t size; 1017b62e03f8SDave May 1018521f74f9SMatthew G. Knepley PetscFunctionBegin; 101928b400f6SJacob Faibussowitsch PetscCheck(swarm->field_registration_initialized, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Must call DMSwarmInitializeFieldRegister() first"); 102028b400f6SJacob Faibussowitsch PetscCheck(!swarm->field_registration_finalized, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Cannot register additional fields after calling DMSwarmFinalizeFieldRegister() first"); 10215f50eb2eSDave May 102208401ef6SPierre Jolivet PetscCheck(type != PETSC_OBJECT, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}"); 102308401ef6SPierre Jolivet PetscCheck(type != PETSC_FUNCTION, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}"); 102408401ef6SPierre Jolivet PetscCheck(type != PETSC_STRING, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}"); 102508401ef6SPierre Jolivet PetscCheck(type != PETSC_STRUCT, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}"); 102608401ef6SPierre Jolivet PetscCheck(type != PETSC_DATATYPE_UNKNOWN, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}"); 1027b62e03f8SDave May 10289566063dSJacob Faibussowitsch PetscCall(PetscDataTypeGetSize(type, &size)); 1029b62e03f8SDave May /* Load a specific data type into data bucket, specifying textual name and its size in bytes */ 10309566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketRegisterField(swarm->db, "DMSwarmRegisterPetscDatatypeField", fieldname, blocksize * size, NULL)); 103152c3ed93SDave May { 103277048351SPatrick Sanan DMSwarmDataField gfield; 103352c3ed93SDave May 10349566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield)); 10359566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldSetBlockSize(gfield, blocksize)); 103652c3ed93SDave May } 1037b62e03f8SDave May swarm->db->field[swarm->db->nfields - 1]->petsc_type = type; 10383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1039b62e03f8SDave May } 1040b62e03f8SDave May 1041d3a51819SDave May /*@C 104220f4b53cSBarry Smith DMSwarmRegisterUserStructField - Register a user defined struct to a `DMSWARM` 1043d3a51819SDave May 104420f4b53cSBarry Smith Collective 1045d3a51819SDave May 104660225df5SJacob Faibussowitsch Input Parameters: 104720f4b53cSBarry Smith + dm - a `DMSWARM` 1048d3a51819SDave May . fieldname - the textual name to identify this field 104962741f57SDave May - size - the size in bytes of the user struct of each data type 1050d3a51819SDave May 1051d3a51819SDave May Level: beginner 1052d3a51819SDave May 105320f4b53cSBarry Smith Note: 10548b8a3813SDave May The textual name for each registered field must be unique. 1055d3a51819SDave May 105620f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmRegisterUserDatatypeField()` 1057d3a51819SDave May @*/ 1058d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmRegisterUserStructField(DM dm, const char fieldname[], size_t size) 1059d71ae5a4SJacob Faibussowitsch { 1060b62e03f8SDave May DM_Swarm *swarm = (DM_Swarm *)dm->data; 1061b62e03f8SDave May 1062521f74f9SMatthew G. Knepley PetscFunctionBegin; 10639566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketRegisterField(swarm->db, "DMSwarmRegisterUserStructField", fieldname, size, NULL)); 1064b62e03f8SDave May swarm->db->field[swarm->db->nfields - 1]->petsc_type = PETSC_STRUCT; 10653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1066b62e03f8SDave May } 1067b62e03f8SDave May 1068d3a51819SDave May /*@C 106920f4b53cSBarry Smith DMSwarmRegisterUserDatatypeField - Register a user defined data type to a `DMSWARM` 1070d3a51819SDave May 107120f4b53cSBarry Smith Collective 1072d3a51819SDave May 107360225df5SJacob Faibussowitsch Input Parameters: 107420f4b53cSBarry Smith + dm - a `DMSWARM` 1075d3a51819SDave May . fieldname - the textual name to identify this field 1076d3a51819SDave May . size - the size in bytes of the user data type 107762741f57SDave May - blocksize - the number of each data type 1078d3a51819SDave May 1079d3a51819SDave May Level: beginner 1080d3a51819SDave May 108120f4b53cSBarry Smith Note: 10828b8a3813SDave May The textual name for each registered field must be unique. 1083d3a51819SDave May 108442747ad1SJacob Faibussowitsch .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmRegisterUserStructField()` 1085d3a51819SDave May @*/ 1086d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmRegisterUserDatatypeField(DM dm, const char fieldname[], size_t size, PetscInt blocksize) 1087d71ae5a4SJacob Faibussowitsch { 1088b62e03f8SDave May DM_Swarm *swarm = (DM_Swarm *)dm->data; 1089b62e03f8SDave May 1090521f74f9SMatthew G. Knepley PetscFunctionBegin; 10919566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketRegisterField(swarm->db, "DMSwarmRegisterUserDatatypeField", fieldname, blocksize * size, NULL)); 1092320740a0SDave May { 109377048351SPatrick Sanan DMSwarmDataField gfield; 1094320740a0SDave May 10959566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield)); 10969566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldSetBlockSize(gfield, blocksize)); 1097320740a0SDave May } 1098b62e03f8SDave May swarm->db->field[swarm->db->nfields - 1]->petsc_type = PETSC_DATATYPE_UNKNOWN; 10993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1100b62e03f8SDave May } 1101b62e03f8SDave May 1102d3a51819SDave May /*@C 1103d3a51819SDave May DMSwarmGetField - Get access to the underlying array storing all entries associated with a registered field 1104d3a51819SDave May 1105cc4c1da9SBarry Smith Not Collective, No Fortran Support 1106d3a51819SDave May 110760225df5SJacob Faibussowitsch Input Parameters: 110820f4b53cSBarry Smith + dm - a `DMSWARM` 110962741f57SDave May - fieldname - the textual name to identify this field 1110d3a51819SDave May 111160225df5SJacob Faibussowitsch Output Parameters: 111262741f57SDave May + blocksize - the number of each data type 1113d3a51819SDave May . type - the data type 111462741f57SDave May - data - pointer to raw array 1115d3a51819SDave May 1116d3a51819SDave May Level: beginner 1117d3a51819SDave May 1118d3a51819SDave May Notes: 111920f4b53cSBarry Smith The array must be returned using a matching call to `DMSwarmRestoreField()`. 1120d3a51819SDave May 112120f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRestoreField()` 1122d3a51819SDave May @*/ 1123d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmGetField(DM dm, const char fieldname[], PetscInt *blocksize, PetscDataType *type, void **data) 1124d71ae5a4SJacob Faibussowitsch { 1125b62e03f8SDave May DM_Swarm *swarm = (DM_Swarm *)dm->data; 112677048351SPatrick Sanan DMSwarmDataField gfield; 1127b62e03f8SDave May 1128521f74f9SMatthew G. Knepley PetscFunctionBegin; 1129ea7032a0SMatthew G. Knepley PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM); 11309566063dSJacob Faibussowitsch if (!swarm->issetup) PetscCall(DMSetUp(dm)); 11319566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield)); 11329566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldGetAccess(gfield)); 11339566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldGetEntries(gfield, data)); 1134ad540459SPierre Jolivet if (blocksize) *blocksize = gfield->bs; 1135ad540459SPierre Jolivet if (type) *type = gfield->petsc_type; 11363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1137b62e03f8SDave May } 1138b62e03f8SDave May 1139d3a51819SDave May /*@C 1140d3a51819SDave May DMSwarmRestoreField - Restore access to the underlying array storing all entries associated with a registered field 1141d3a51819SDave May 1142cc4c1da9SBarry Smith Not Collective, No Fortran Support 1143d3a51819SDave May 114460225df5SJacob Faibussowitsch Input Parameters: 114520f4b53cSBarry Smith + dm - a `DMSWARM` 114662741f57SDave May - fieldname - the textual name to identify this field 1147d3a51819SDave May 114860225df5SJacob Faibussowitsch Output Parameters: 114962741f57SDave May + blocksize - the number of each data type 1150d3a51819SDave May . type - the data type 115162741f57SDave May - data - pointer to raw array 1152d3a51819SDave May 1153d3a51819SDave May Level: beginner 1154d3a51819SDave May 1155d3a51819SDave May Notes: 115620f4b53cSBarry Smith The user must call `DMSwarmGetField()` prior to calling `DMSwarmRestoreField()`. 1157d3a51819SDave May 115820f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmGetField()` 1159d3a51819SDave May @*/ 1160d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmRestoreField(DM dm, const char fieldname[], PetscInt *blocksize, PetscDataType *type, void **data) 1161d71ae5a4SJacob Faibussowitsch { 1162b62e03f8SDave May DM_Swarm *swarm = (DM_Swarm *)dm->data; 116377048351SPatrick Sanan DMSwarmDataField gfield; 1164b62e03f8SDave May 1165521f74f9SMatthew G. Knepley PetscFunctionBegin; 1166ea7032a0SMatthew G. Knepley PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM); 11679566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield)); 11689566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldRestoreAccess(gfield)); 1169b62e03f8SDave May if (data) *data = NULL; 11703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1171b62e03f8SDave May } 1172b62e03f8SDave May 11730bf7c1c5SMatthew G. Knepley PetscErrorCode DMSwarmGetFieldInfo(DM dm, const char fieldname[], PetscInt *blocksize, PetscDataType *type) 11740bf7c1c5SMatthew G. Knepley { 11750bf7c1c5SMatthew G. Knepley DM_Swarm *swarm = (DM_Swarm *)dm->data; 11760bf7c1c5SMatthew G. Knepley DMSwarmDataField gfield; 11770bf7c1c5SMatthew G. Knepley 11780bf7c1c5SMatthew G. Knepley PetscFunctionBegin; 11790bf7c1c5SMatthew G. Knepley PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM); 11800bf7c1c5SMatthew G. Knepley PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield)); 11810bf7c1c5SMatthew G. Knepley if (blocksize) *blocksize = gfield->bs; 11820bf7c1c5SMatthew G. Knepley if (type) *type = gfield->petsc_type; 11830bf7c1c5SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 11840bf7c1c5SMatthew G. Knepley } 11850bf7c1c5SMatthew G. Knepley 118674d0cae8SMatthew G. Knepley /*@ 118720f4b53cSBarry Smith DMSwarmAddPoint - Add space for one new point in the `DMSWARM` 1188d3a51819SDave May 118920f4b53cSBarry Smith Not Collective 1190d3a51819SDave May 119160225df5SJacob Faibussowitsch Input Parameter: 119220f4b53cSBarry Smith . dm - a `DMSWARM` 1193d3a51819SDave May 1194d3a51819SDave May Level: beginner 1195d3a51819SDave May 1196d3a51819SDave May Notes: 11978b8a3813SDave May The new point will have all fields initialized to zero. 1198d3a51819SDave May 119920f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmAddNPoints()` 1200d3a51819SDave May @*/ 1201d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmAddPoint(DM dm) 1202d71ae5a4SJacob Faibussowitsch { 1203cb1d1399SDave May DM_Swarm *swarm = (DM_Swarm *)dm->data; 1204cb1d1399SDave May 1205521f74f9SMatthew G. Knepley PetscFunctionBegin; 12069566063dSJacob Faibussowitsch if (!swarm->issetup) PetscCall(DMSetUp(dm)); 12079566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(DMSWARM_AddPoints, 0, 0, 0, 0)); 12089566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketAddPoint(swarm->db)); 12099566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(DMSWARM_AddPoints, 0, 0, 0, 0)); 12103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1211cb1d1399SDave May } 1212cb1d1399SDave May 121374d0cae8SMatthew G. Knepley /*@ 121420f4b53cSBarry Smith DMSwarmAddNPoints - Add space for a number of new points in the `DMSWARM` 1215d3a51819SDave May 121620f4b53cSBarry Smith Not Collective 1217d3a51819SDave May 121860225df5SJacob Faibussowitsch Input Parameters: 121920f4b53cSBarry Smith + dm - a `DMSWARM` 122062741f57SDave May - npoints - the number of new points to add 1221d3a51819SDave May 1222d3a51819SDave May Level: beginner 1223d3a51819SDave May 1224d3a51819SDave May Notes: 12258b8a3813SDave May The new point will have all fields initialized to zero. 1226d3a51819SDave May 122720f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmAddPoint()` 1228d3a51819SDave May @*/ 1229d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmAddNPoints(DM dm, PetscInt npoints) 1230d71ae5a4SJacob Faibussowitsch { 1231cb1d1399SDave May DM_Swarm *swarm = (DM_Swarm *)dm->data; 1232cb1d1399SDave May PetscInt nlocal; 1233cb1d1399SDave May 1234521f74f9SMatthew G. Knepley PetscFunctionBegin; 12359566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(DMSWARM_AddPoints, 0, 0, 0, 0)); 12369566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &nlocal, NULL, NULL)); 1237cb1d1399SDave May nlocal = nlocal + npoints; 12389566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketSetSizes(swarm->db, nlocal, DMSWARM_DATA_BUCKET_BUFFER_DEFAULT)); 12399566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(DMSWARM_AddPoints, 0, 0, 0, 0)); 12403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1241cb1d1399SDave May } 1242cb1d1399SDave May 124374d0cae8SMatthew G. Knepley /*@ 124420f4b53cSBarry Smith DMSwarmRemovePoint - Remove the last point from the `DMSWARM` 1245d3a51819SDave May 124620f4b53cSBarry Smith Not Collective 1247d3a51819SDave May 124860225df5SJacob Faibussowitsch Input Parameter: 124920f4b53cSBarry Smith . dm - a `DMSWARM` 1250d3a51819SDave May 1251d3a51819SDave May Level: beginner 1252d3a51819SDave May 125320f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRemovePointAtIndex()` 1254d3a51819SDave May @*/ 1255d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmRemovePoint(DM dm) 1256d71ae5a4SJacob Faibussowitsch { 1257cb1d1399SDave May DM_Swarm *swarm = (DM_Swarm *)dm->data; 1258cb1d1399SDave May 1259521f74f9SMatthew G. Knepley PetscFunctionBegin; 12609566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(DMSWARM_RemovePoints, 0, 0, 0, 0)); 12619566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketRemovePoint(swarm->db)); 12629566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(DMSWARM_RemovePoints, 0, 0, 0, 0)); 12633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1264cb1d1399SDave May } 1265cb1d1399SDave May 126674d0cae8SMatthew G. Knepley /*@ 126720f4b53cSBarry Smith DMSwarmRemovePointAtIndex - Removes a specific point from the `DMSWARM` 1268d3a51819SDave May 126920f4b53cSBarry Smith Not Collective 1270d3a51819SDave May 127160225df5SJacob Faibussowitsch Input Parameters: 127220f4b53cSBarry Smith + dm - a `DMSWARM` 127362741f57SDave May - idx - index of point to remove 1274d3a51819SDave May 1275d3a51819SDave May Level: beginner 1276d3a51819SDave May 127720f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRemovePoint()` 1278d3a51819SDave May @*/ 1279d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmRemovePointAtIndex(DM dm, PetscInt idx) 1280d71ae5a4SJacob Faibussowitsch { 1281cb1d1399SDave May DM_Swarm *swarm = (DM_Swarm *)dm->data; 1282cb1d1399SDave May 1283521f74f9SMatthew G. Knepley PetscFunctionBegin; 12849566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(DMSWARM_RemovePoints, 0, 0, 0, 0)); 12859566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketRemovePointAtIndex(swarm->db, idx)); 12869566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(DMSWARM_RemovePoints, 0, 0, 0, 0)); 12873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1288cb1d1399SDave May } 1289b62e03f8SDave May 129074d0cae8SMatthew G. Knepley /*@ 129120f4b53cSBarry Smith DMSwarmCopyPoint - Copy point pj to point pi in the `DMSWARM` 1292ba4fc9c6SDave May 129320f4b53cSBarry Smith Not Collective 1294ba4fc9c6SDave May 129560225df5SJacob Faibussowitsch Input Parameters: 129620f4b53cSBarry Smith + dm - a `DMSWARM` 1297ba4fc9c6SDave May . pi - the index of the point to copy 1298ba4fc9c6SDave May - pj - the point index where the copy should be located 1299ba4fc9c6SDave May 1300ba4fc9c6SDave May Level: beginner 1301ba4fc9c6SDave May 130220f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRemovePoint()` 1303ba4fc9c6SDave May @*/ 1304d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmCopyPoint(DM dm, PetscInt pi, PetscInt pj) 1305d71ae5a4SJacob Faibussowitsch { 1306ba4fc9c6SDave May DM_Swarm *swarm = (DM_Swarm *)dm->data; 1307ba4fc9c6SDave May 1308ba4fc9c6SDave May PetscFunctionBegin; 13099566063dSJacob Faibussowitsch if (!swarm->issetup) PetscCall(DMSetUp(dm)); 13109566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketCopyPoint(swarm->db, pi, swarm->db, pj)); 13113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1312ba4fc9c6SDave May } 1313ba4fc9c6SDave May 131466976f2fSJacob Faibussowitsch static PetscErrorCode DMSwarmMigrate_Basic(DM dm, PetscBool remove_sent_points) 1315d71ae5a4SJacob Faibussowitsch { 1316521f74f9SMatthew G. Knepley PetscFunctionBegin; 13179566063dSJacob Faibussowitsch PetscCall(DMSwarmMigrate_Push_Basic(dm, remove_sent_points)); 13183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13193454631fSDave May } 13203454631fSDave May 132174d0cae8SMatthew G. Knepley /*@ 132220f4b53cSBarry Smith DMSwarmMigrate - Relocates points defined in the `DMSWARM` to other MPI-ranks 1323d3a51819SDave May 132420f4b53cSBarry Smith Collective 1325d3a51819SDave May 132660225df5SJacob Faibussowitsch Input Parameters: 132720f4b53cSBarry Smith + dm - the `DMSWARM` 132862741f57SDave May - remove_sent_points - flag indicating if sent points should be removed from the current MPI-rank 1329d3a51819SDave May 1330d3a51819SDave May Level: advanced 1331d3a51819SDave May 133220f4b53cSBarry Smith Notes: 133320f4b53cSBarry Smith The `DM` will be modified to accommodate received points. 133420f4b53cSBarry Smith If `remove_sent_points` is `PETSC_TRUE`, any points that were sent will be removed from the `DM`. 133520f4b53cSBarry Smith Different styles of migration are supported. See `DMSwarmSetMigrateType()`. 133620f4b53cSBarry Smith 133720f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmSetMigrateType()` 1338d3a51819SDave May @*/ 1339d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmMigrate(DM dm, PetscBool remove_sent_points) 1340d71ae5a4SJacob Faibussowitsch { 1341f0cdbbbaSDave May DM_Swarm *swarm = (DM_Swarm *)dm->data; 1342f0cdbbbaSDave May 1343521f74f9SMatthew G. Knepley PetscFunctionBegin; 13449566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(DMSWARM_Migrate, 0, 0, 0, 0)); 1345f0cdbbbaSDave May switch (swarm->migrate_type) { 1346d71ae5a4SJacob Faibussowitsch case DMSWARM_MIGRATE_BASIC: 1347d71ae5a4SJacob Faibussowitsch PetscCall(DMSwarmMigrate_Basic(dm, remove_sent_points)); 1348d71ae5a4SJacob Faibussowitsch break; 1349d71ae5a4SJacob Faibussowitsch case DMSWARM_MIGRATE_DMCELLNSCATTER: 1350d71ae5a4SJacob Faibussowitsch PetscCall(DMSwarmMigrate_CellDMScatter(dm, remove_sent_points)); 1351d71ae5a4SJacob Faibussowitsch break; 1352d71ae5a4SJacob Faibussowitsch case DMSWARM_MIGRATE_DMCELLEXACT: 1353d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_MIGRATE_DMCELLEXACT not implemented"); 1354d71ae5a4SJacob Faibussowitsch case DMSWARM_MIGRATE_USER: 1355d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_MIGRATE_USER not implemented"); 1356d71ae5a4SJacob Faibussowitsch default: 1357d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_MIGRATE type unknown"); 1358f0cdbbbaSDave May } 13599566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(DMSWARM_Migrate, 0, 0, 0, 0)); 13609566063dSJacob Faibussowitsch PetscCall(DMClearGlobalVectors(dm)); 13613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13623454631fSDave May } 13633454631fSDave May 1364f0cdbbbaSDave May PetscErrorCode DMSwarmMigrate_GlobalToLocal_Basic(DM dm, PetscInt *globalsize); 1365f0cdbbbaSDave May 1366d3a51819SDave May /* 1367d3a51819SDave May DMSwarmCollectViewCreate 1368d3a51819SDave May 1369d3a51819SDave May * Applies a collection method and gathers point neighbour points into dm 1370d3a51819SDave May 1371d3a51819SDave May Notes: 13728b8a3813SDave May Users should call DMSwarmCollectViewDestroy() after 1373d3a51819SDave May they have finished computations associated with the collected points 1374d3a51819SDave May */ 1375d3a51819SDave May 137674d0cae8SMatthew G. Knepley /*@ 1377d3a51819SDave May DMSwarmCollectViewCreate - Applies a collection method and gathers points 137820f4b53cSBarry Smith in neighbour ranks into the `DMSWARM` 1379d3a51819SDave May 138020f4b53cSBarry Smith Collective 1381d3a51819SDave May 138260225df5SJacob Faibussowitsch Input Parameter: 138320f4b53cSBarry Smith . dm - the `DMSWARM` 1384d3a51819SDave May 1385d3a51819SDave May Level: advanced 1386d3a51819SDave May 138720f4b53cSBarry Smith Notes: 138820f4b53cSBarry Smith Users should call `DMSwarmCollectViewDestroy()` after 138920f4b53cSBarry Smith they have finished computations associated with the collected points 13900764c050SBarry Smith 139120f4b53cSBarry Smith Different collect methods are supported. See `DMSwarmSetCollectType()`. 139220f4b53cSBarry Smith 13930764c050SBarry Smith Developer Note: 13940764c050SBarry Smith Create and Destroy routines create new objects that can get destroyed, they do not change the state 13950764c050SBarry Smith of the current object. 13960764c050SBarry Smith 139720f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmCollectViewDestroy()`, `DMSwarmSetCollectType()` 1398d3a51819SDave May @*/ 1399d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmCollectViewCreate(DM dm) 1400d71ae5a4SJacob Faibussowitsch { 14012712d1f2SDave May DM_Swarm *swarm = (DM_Swarm *)dm->data; 14022712d1f2SDave May PetscInt ng; 14032712d1f2SDave May 1404521f74f9SMatthew G. Knepley PetscFunctionBegin; 140528b400f6SJacob Faibussowitsch PetscCheck(!swarm->collect_view_active, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "CollectView currently active"); 14069566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(dm, &ng)); 1407480eef7bSDave May switch (swarm->collect_type) { 1408d71ae5a4SJacob Faibussowitsch case DMSWARM_COLLECT_BASIC: 1409d71ae5a4SJacob Faibussowitsch PetscCall(DMSwarmMigrate_GlobalToLocal_Basic(dm, &ng)); 1410d71ae5a4SJacob Faibussowitsch break; 1411d71ae5a4SJacob Faibussowitsch case DMSWARM_COLLECT_DMDABOUNDINGBOX: 1412d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_COLLECT_DMDABOUNDINGBOX not implemented"); 1413d71ae5a4SJacob Faibussowitsch case DMSWARM_COLLECT_GENERAL: 1414d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_COLLECT_GENERAL not implemented"); 1415d71ae5a4SJacob Faibussowitsch default: 1416d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_COLLECT type unknown"); 1417480eef7bSDave May } 1418480eef7bSDave May swarm->collect_view_active = PETSC_TRUE; 1419480eef7bSDave May swarm->collect_view_reset_nlocal = ng; 14203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14212712d1f2SDave May } 14222712d1f2SDave May 142374d0cae8SMatthew G. Knepley /*@ 142420f4b53cSBarry Smith DMSwarmCollectViewDestroy - Resets the `DMSWARM` to the size prior to calling `DMSwarmCollectViewCreate()` 1425d3a51819SDave May 142620f4b53cSBarry Smith Collective 1427d3a51819SDave May 142860225df5SJacob Faibussowitsch Input Parameters: 142920f4b53cSBarry Smith . dm - the `DMSWARM` 1430d3a51819SDave May 1431d3a51819SDave May Notes: 143220f4b53cSBarry Smith Users should call `DMSwarmCollectViewCreate()` before this function is called. 1433d3a51819SDave May 1434d3a51819SDave May Level: advanced 1435d3a51819SDave May 14360764c050SBarry Smith Developer Note: 14370764c050SBarry Smith Create and Destroy routines create new objects that can get destroyed, they do not change the state 14380764c050SBarry Smith of the current object. 14390764c050SBarry Smith 144020f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmCollectViewCreate()`, `DMSwarmSetCollectType()` 1441d3a51819SDave May @*/ 1442d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmCollectViewDestroy(DM dm) 1443d71ae5a4SJacob Faibussowitsch { 14442712d1f2SDave May DM_Swarm *swarm = (DM_Swarm *)dm->data; 14452712d1f2SDave May 1446521f74f9SMatthew G. Knepley PetscFunctionBegin; 144728b400f6SJacob Faibussowitsch PetscCheck(swarm->collect_view_active, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "CollectView is currently not active"); 14489566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(dm, swarm->collect_view_reset_nlocal, -1)); 1449480eef7bSDave May swarm->collect_view_active = PETSC_FALSE; 14503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14512712d1f2SDave May } 14523454631fSDave May 145366976f2fSJacob Faibussowitsch static PetscErrorCode DMSwarmSetUpPIC(DM dm) 1454d71ae5a4SJacob Faibussowitsch { 1455f0cdbbbaSDave May PetscInt dim; 1456f0cdbbbaSDave May 1457521f74f9SMatthew G. Knepley PetscFunctionBegin; 14589566063dSJacob Faibussowitsch PetscCall(DMSwarmSetNumSpecies(dm, 1)); 14599566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 146063a3b9bcSJacob Faibussowitsch PetscCheck(dim >= 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Dimension must be 1,2,3 - found %" PetscInt_FMT, dim); 146163a3b9bcSJacob Faibussowitsch PetscCheck(dim <= 3, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Dimension must be 1,2,3 - found %" PetscInt_FMT, dim); 14629566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(dm, DMSwarmPICField_coor, dim, PETSC_DOUBLE)); 14639566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(dm, DMSwarmPICField_cellid, 1, PETSC_INT)); 14643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1465f0cdbbbaSDave May } 1466f0cdbbbaSDave May 146774d0cae8SMatthew G. Knepley /*@ 146831403186SMatthew G. Knepley DMSwarmSetPointCoordinatesRandom - Sets initial coordinates for particles in each cell 146931403186SMatthew G. Knepley 147020f4b53cSBarry Smith Collective 147131403186SMatthew G. Knepley 147260225df5SJacob Faibussowitsch Input Parameters: 147320f4b53cSBarry Smith + dm - the `DMSWARM` 147420f4b53cSBarry Smith - Npc - The number of particles per cell in the cell `DM` 147531403186SMatthew G. Knepley 147631403186SMatthew G. Knepley Level: intermediate 147731403186SMatthew G. Knepley 147820f4b53cSBarry Smith Notes: 147920f4b53cSBarry Smith The user must use `DMSwarmSetCellDM()` to set the cell `DM` first. The particles are placed randomly inside each cell. If only 148020f4b53cSBarry Smith one particle is in each cell, it is placed at the centroid. 148120f4b53cSBarry Smith 148220f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmSetCellDM()` 148331403186SMatthew G. Knepley @*/ 1484d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmSetPointCoordinatesRandom(DM dm, PetscInt Npc) 1485d71ae5a4SJacob Faibussowitsch { 148631403186SMatthew G. Knepley DM cdm; 148731403186SMatthew G. Knepley PetscRandom rnd; 148831403186SMatthew G. Knepley DMPolytopeType ct; 148931403186SMatthew G. Knepley PetscBool simplex; 149031403186SMatthew G. Knepley PetscReal *centroid, *coords, *xi0, *v0, *J, *invJ, detJ; 149131403186SMatthew G. Knepley PetscInt dim, d, cStart, cEnd, c, p; 149231403186SMatthew G. Knepley 149331403186SMatthew G. Knepley PetscFunctionBeginUser; 14949566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)dm), &rnd)); 14959566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rnd, -1.0, 1.0)); 14969566063dSJacob Faibussowitsch PetscCall(PetscRandomSetType(rnd, PETSCRAND48)); 149731403186SMatthew G. Knepley 14989566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(dm, &cdm)); 14999566063dSJacob Faibussowitsch PetscCall(DMGetDimension(cdm, &dim)); 15009566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(cdm, 0, &cStart, &cEnd)); 15019566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(cdm, cStart, &ct)); 150231403186SMatthew G. Knepley simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE; 150331403186SMatthew G. Knepley 15049566063dSJacob Faibussowitsch PetscCall(PetscMalloc5(dim, ¢roid, dim, &xi0, dim, &v0, dim * dim, &J, dim * dim, &invJ)); 150531403186SMatthew G. Knepley for (d = 0; d < dim; ++d) xi0[d] = -1.0; 15069566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 150731403186SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 150831403186SMatthew G. Knepley if (Npc == 1) { 15099566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFVM(cdm, c, NULL, centroid, NULL)); 151031403186SMatthew G. Knepley for (d = 0; d < dim; ++d) coords[c * dim + d] = centroid[d]; 151131403186SMatthew G. Knepley } else { 15129566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFEM(cdm, c, NULL, v0, J, invJ, &detJ)); /* affine */ 151331403186SMatthew G. Knepley for (p = 0; p < Npc; ++p) { 151431403186SMatthew G. Knepley const PetscInt n = c * Npc + p; 151531403186SMatthew G. Knepley PetscReal sum = 0.0, refcoords[3]; 151631403186SMatthew G. Knepley 151731403186SMatthew G. Knepley for (d = 0; d < dim; ++d) { 15189566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValueReal(rnd, &refcoords[d])); 151931403186SMatthew G. Knepley sum += refcoords[d]; 152031403186SMatthew G. Knepley } 15219371c9d4SSatish Balay if (simplex && sum > 0.0) 15229371c9d4SSatish Balay for (d = 0; d < dim; ++d) refcoords[d] -= PetscSqrtReal(dim) * sum; 152331403186SMatthew G. Knepley CoordinatesRefToReal(dim, dim, xi0, v0, J, refcoords, &coords[n * dim]); 152431403186SMatthew G. Knepley } 152531403186SMatthew G. Knepley } 152631403186SMatthew G. Knepley } 15279566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 15289566063dSJacob Faibussowitsch PetscCall(PetscFree5(centroid, xi0, v0, J, invJ)); 15299566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rnd)); 15303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 153131403186SMatthew G. Knepley } 153231403186SMatthew G. Knepley 153331403186SMatthew G. Knepley /*@ 153420f4b53cSBarry Smith DMSwarmSetType - Set particular flavor of `DMSWARM` 1535d3a51819SDave May 153620f4b53cSBarry Smith Collective 1537d3a51819SDave May 153860225df5SJacob Faibussowitsch Input Parameters: 153920f4b53cSBarry Smith + dm - the `DMSWARM` 154020f4b53cSBarry Smith - stype - the `DMSWARM` type (e.g. `DMSWARM_PIC`) 1541d3a51819SDave May 1542d3a51819SDave May Level: advanced 1543d3a51819SDave May 154420f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmSetMigrateType()`, `DMSwarmSetCollectType()`, `DMSwarmType`, `DMSWARM_PIC`, `DMSWARM_BASIC` 1545d3a51819SDave May @*/ 1546d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmSetType(DM dm, DMSwarmType stype) 1547d71ae5a4SJacob Faibussowitsch { 1548f0cdbbbaSDave May DM_Swarm *swarm = (DM_Swarm *)dm->data; 1549f0cdbbbaSDave May 1550521f74f9SMatthew G. Knepley PetscFunctionBegin; 1551f0cdbbbaSDave May swarm->swarm_type = stype; 155248a46eb9SPierre Jolivet if (swarm->swarm_type == DMSWARM_PIC) PetscCall(DMSwarmSetUpPIC(dm)); 15533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1554f0cdbbbaSDave May } 1555f0cdbbbaSDave May 155666976f2fSJacob Faibussowitsch static PetscErrorCode DMSetup_Swarm(DM dm) 1557d71ae5a4SJacob Faibussowitsch { 15583454631fSDave May DM_Swarm *swarm = (DM_Swarm *)dm->data; 15593454631fSDave May PetscMPIInt rank; 15603454631fSDave May PetscInt p, npoints, *rankval; 15613454631fSDave May 1562521f74f9SMatthew G. Knepley PetscFunctionBegin; 15633ba16761SJacob Faibussowitsch if (swarm->issetup) PetscFunctionReturn(PETSC_SUCCESS); 15643454631fSDave May swarm->issetup = PETSC_TRUE; 15653454631fSDave May 1566f0cdbbbaSDave May if (swarm->swarm_type == DMSWARM_PIC) { 1567f0cdbbbaSDave May /* check dmcell exists */ 156828b400f6SJacob Faibussowitsch PetscCheck(swarm->dmcell, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "DMSWARM_PIC requires you call DMSwarmSetCellDM"); 1569f0cdbbbaSDave May 1570f0cdbbbaSDave May if (swarm->dmcell->ops->locatepointssubdomain) { 1571f0cdbbbaSDave May /* check methods exists for exact ownership identificiation */ 15729566063dSJacob Faibussowitsch PetscCall(PetscInfo(dm, "DMSWARM_PIC: Using method CellDM->ops->LocatePointsSubdomain\n")); 1573f0cdbbbaSDave May swarm->migrate_type = DMSWARM_MIGRATE_DMCELLEXACT; 1574f0cdbbbaSDave May } else { 1575f0cdbbbaSDave May /* check methods exist for point location AND rank neighbor identification */ 1576f0cdbbbaSDave May if (swarm->dmcell->ops->locatepoints) { 15779566063dSJacob Faibussowitsch PetscCall(PetscInfo(dm, "DMSWARM_PIC: Using method CellDM->LocatePoints\n")); 1578f0cdbbbaSDave May } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "DMSWARM_PIC requires the method CellDM->ops->locatepoints be defined"); 1579f0cdbbbaSDave May 1580f0cdbbbaSDave May if (swarm->dmcell->ops->getneighbors) { 15819566063dSJacob Faibussowitsch PetscCall(PetscInfo(dm, "DMSWARM_PIC: Using method CellDM->GetNeigbors\n")); 1582f0cdbbbaSDave May } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "DMSWARM_PIC requires the method CellDM->ops->getneighbors be defined"); 1583f0cdbbbaSDave May 1584f0cdbbbaSDave May swarm->migrate_type = DMSWARM_MIGRATE_DMCELLNSCATTER; 1585f0cdbbbaSDave May } 1586f0cdbbbaSDave May } 1587f0cdbbbaSDave May 15889566063dSJacob Faibussowitsch PetscCall(DMSwarmFinalizeFieldRegister(dm)); 1589f0cdbbbaSDave May 15903454631fSDave May /* check some fields were registered */ 159108401ef6SPierre Jolivet PetscCheck(swarm->db->nfields > 2, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "At least one field user must be registered via DMSwarmRegisterXXX()"); 15923454631fSDave May 15933454631fSDave May /* check local sizes were set */ 159408401ef6SPierre Jolivet PetscCheck(swarm->db->L != -1, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Local sizes must be set via DMSwarmSetLocalSizes()"); 15953454631fSDave May 15963454631fSDave May /* initialize values in pid and rank placeholders */ 15973454631fSDave May /* TODO: [pid - use MPI_Scan] */ 15989566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 15999566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL)); 16009566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval)); 1601ad540459SPierre Jolivet for (p = 0; p < npoints; p++) rankval[p] = (PetscInt)rank; 16029566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval)); 16033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16043454631fSDave May } 16053454631fSDave May 1606dc5f5ce9SDave May extern PetscErrorCode DMSwarmSortDestroy(DMSwarmSort *_ctx); 1607dc5f5ce9SDave May 160866976f2fSJacob Faibussowitsch static PetscErrorCode DMDestroy_Swarm(DM dm) 1609d71ae5a4SJacob Faibussowitsch { 161057795646SDave May DM_Swarm *swarm = (DM_Swarm *)dm->data; 161157795646SDave May 161257795646SDave May PetscFunctionBegin; 16133ba16761SJacob Faibussowitsch if (--swarm->refct > 0) PetscFunctionReturn(PETSC_SUCCESS); 16149566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketDestroy(&swarm->db)); 161548a46eb9SPierre Jolivet if (swarm->sort_context) PetscCall(DMSwarmSortDestroy(&swarm->sort_context)); 16169566063dSJacob Faibussowitsch PetscCall(PetscFree(swarm)); 16173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 161857795646SDave May } 161957795646SDave May 162066976f2fSJacob Faibussowitsch static PetscErrorCode DMSwarmView_Draw(DM dm, PetscViewer viewer) 1621d71ae5a4SJacob Faibussowitsch { 1622a9ee3421SMatthew G. Knepley DM cdm; 1623a9ee3421SMatthew G. Knepley PetscDraw draw; 162431403186SMatthew G. Knepley PetscReal *coords, oldPause, radius = 0.01; 1625a9ee3421SMatthew G. Knepley PetscInt Np, p, bs; 1626a9ee3421SMatthew G. Knepley 1627a9ee3421SMatthew G. Knepley PetscFunctionBegin; 16289566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, ((PetscObject)dm)->prefix, "-dm_view_swarm_radius", &radius, NULL)); 16299566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw)); 16309566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(dm, &cdm)); 16319566063dSJacob Faibussowitsch PetscCall(PetscDrawGetPause(draw, &oldPause)); 16329566063dSJacob Faibussowitsch PetscCall(PetscDrawSetPause(draw, 0.0)); 16339566063dSJacob Faibussowitsch PetscCall(DMView(cdm, viewer)); 16349566063dSJacob Faibussowitsch PetscCall(PetscDrawSetPause(draw, oldPause)); 1635a9ee3421SMatthew G. Knepley 16369566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(dm, &Np)); 16379566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, &bs, NULL, (void **)&coords)); 1638a9ee3421SMatthew G. Knepley for (p = 0; p < Np; ++p) { 1639a9ee3421SMatthew G. Knepley const PetscInt i = p * bs; 1640a9ee3421SMatthew G. Knepley 16419566063dSJacob Faibussowitsch PetscCall(PetscDrawEllipse(draw, coords[i], coords[i + 1], radius, radius, PETSC_DRAW_BLUE)); 1642a9ee3421SMatthew G. Knepley } 16439566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, &bs, NULL, (void **)&coords)); 16449566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 16459566063dSJacob Faibussowitsch PetscCall(PetscDrawPause(draw)); 16469566063dSJacob Faibussowitsch PetscCall(PetscDrawSave(draw)); 16473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1648a9ee3421SMatthew G. Knepley } 1649a9ee3421SMatthew G. Knepley 1650d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMView_Swarm_Ascii(DM dm, PetscViewer viewer) 1651d71ae5a4SJacob Faibussowitsch { 16526a5217c0SMatthew G. Knepley PetscViewerFormat format; 16536a5217c0SMatthew G. Knepley PetscInt *sizes; 16546a5217c0SMatthew G. Knepley PetscInt dim, Np, maxSize = 17; 16556a5217c0SMatthew G. Knepley MPI_Comm comm; 16566a5217c0SMatthew G. Knepley PetscMPIInt rank, size, p; 16576a5217c0SMatthew G. Knepley const char *name; 16586a5217c0SMatthew G. Knepley 16596a5217c0SMatthew G. Knepley PetscFunctionBegin; 16606a5217c0SMatthew G. Knepley PetscCall(PetscViewerGetFormat(viewer, &format)); 16616a5217c0SMatthew G. Knepley PetscCall(DMGetDimension(dm, &dim)); 16626a5217c0SMatthew G. Knepley PetscCall(DMSwarmGetLocalSize(dm, &Np)); 16636a5217c0SMatthew G. Knepley PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 16646a5217c0SMatthew G. Knepley PetscCallMPI(MPI_Comm_rank(comm, &rank)); 16656a5217c0SMatthew G. Knepley PetscCallMPI(MPI_Comm_size(comm, &size)); 16666a5217c0SMatthew G. Knepley PetscCall(PetscObjectGetName((PetscObject)dm, &name)); 166763a3b9bcSJacob Faibussowitsch if (name) PetscCall(PetscViewerASCIIPrintf(viewer, "%s in %" PetscInt_FMT " dimension%s:\n", name, dim, dim == 1 ? "" : "s")); 166863a3b9bcSJacob Faibussowitsch else PetscCall(PetscViewerASCIIPrintf(viewer, "Swarm in %" PetscInt_FMT " dimension%s:\n", dim, dim == 1 ? "" : "s")); 16696a5217c0SMatthew G. Knepley if (size < maxSize) PetscCall(PetscCalloc1(size, &sizes)); 16706a5217c0SMatthew G. Knepley else PetscCall(PetscCalloc1(3, &sizes)); 16716a5217c0SMatthew G. Knepley if (size < maxSize) { 16726a5217c0SMatthew G. Knepley PetscCallMPI(MPI_Gather(&Np, 1, MPIU_INT, sizes, 1, MPIU_INT, 0, comm)); 16736a5217c0SMatthew G. Knepley PetscCall(PetscViewerASCIIPrintf(viewer, " Number of particles per rank:")); 16746a5217c0SMatthew G. Knepley for (p = 0; p < size; ++p) { 16756a5217c0SMatthew G. Knepley if (rank == 0) PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT, sizes[p])); 16766a5217c0SMatthew G. Knepley } 16776a5217c0SMatthew G. Knepley } else { 16786a5217c0SMatthew G. Knepley PetscInt locMinMax[2] = {Np, Np}; 16796a5217c0SMatthew G. Knepley 16806a5217c0SMatthew G. Knepley PetscCall(PetscGlobalMinMaxInt(comm, locMinMax, sizes)); 16816a5217c0SMatthew G. Knepley PetscCall(PetscViewerASCIIPrintf(viewer, " Min/Max of particles per rank: %" PetscInt_FMT "/%" PetscInt_FMT, sizes[0], sizes[1])); 16826a5217c0SMatthew G. Knepley } 16836a5217c0SMatthew G. Knepley PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 16846a5217c0SMatthew G. Knepley PetscCall(PetscFree(sizes)); 16856a5217c0SMatthew G. Knepley if (format == PETSC_VIEWER_ASCII_INFO) { 16866a5217c0SMatthew G. Knepley PetscInt *cell; 16876a5217c0SMatthew G. Knepley 16886a5217c0SMatthew G. Knepley PetscCall(PetscViewerASCIIPrintf(viewer, " Cells containing each particle:\n")); 16896a5217c0SMatthew G. Knepley PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 16906a5217c0SMatthew G. Knepley PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&cell)); 169163a3b9bcSJacob Faibussowitsch for (p = 0; p < Np; ++p) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " p%d: %" PetscInt_FMT "\n", p, cell[p])); 16926a5217c0SMatthew G. Knepley PetscCall(PetscViewerFlush(viewer)); 16936a5217c0SMatthew G. Knepley PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 16946a5217c0SMatthew G. Knepley PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&cell)); 16956a5217c0SMatthew G. Knepley } 16963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16976a5217c0SMatthew G. Knepley } 16986a5217c0SMatthew G. Knepley 169966976f2fSJacob Faibussowitsch static PetscErrorCode DMView_Swarm(DM dm, PetscViewer viewer) 1700d71ae5a4SJacob Faibussowitsch { 17015f50eb2eSDave May DM_Swarm *swarm = (DM_Swarm *)dm->data; 17025627991aSBarry Smith PetscBool iascii, ibinary, isvtk, isdraw; 17035627991aSBarry Smith #if defined(PETSC_HAVE_HDF5) 17045627991aSBarry Smith PetscBool ishdf5; 17055627991aSBarry Smith #endif 17065f50eb2eSDave May 17075f50eb2eSDave May PetscFunctionBegin; 17085f50eb2eSDave May PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 17095f50eb2eSDave May PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 17109566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 17119566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &ibinary)); 17129566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERVTK, &isvtk)); 17135627991aSBarry Smith #if defined(PETSC_HAVE_HDF5) 17149566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5)); 17155627991aSBarry Smith #endif 17169566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw)); 17175f50eb2eSDave May if (iascii) { 17186a5217c0SMatthew G. Knepley PetscViewerFormat format; 17196a5217c0SMatthew G. Knepley 17206a5217c0SMatthew G. Knepley PetscCall(PetscViewerGetFormat(viewer, &format)); 17216a5217c0SMatthew G. Knepley switch (format) { 1722d71ae5a4SJacob Faibussowitsch case PETSC_VIEWER_ASCII_INFO_DETAIL: 1723d71ae5a4SJacob Faibussowitsch PetscCall(DMSwarmDataBucketView(PetscObjectComm((PetscObject)dm), swarm->db, NULL, DATABUCKET_VIEW_STDOUT)); 1724d71ae5a4SJacob Faibussowitsch break; 1725d71ae5a4SJacob Faibussowitsch default: 1726d71ae5a4SJacob Faibussowitsch PetscCall(DMView_Swarm_Ascii(dm, viewer)); 17276a5217c0SMatthew G. Knepley } 1728f7d195e4SLawrence Mitchell } else { 17295f50eb2eSDave May #if defined(PETSC_HAVE_HDF5) 173048a46eb9SPierre Jolivet if (ishdf5) PetscCall(DMSwarmView_HDF5(dm, viewer)); 17315f50eb2eSDave May #endif 173248a46eb9SPierre Jolivet if (isdraw) PetscCall(DMSwarmView_Draw(dm, viewer)); 1733f7d195e4SLawrence Mitchell } 17343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 17355f50eb2eSDave May } 17365f50eb2eSDave May 1737cc4c1da9SBarry Smith /*@ 173820f4b53cSBarry Smith DMSwarmGetCellSwarm - Extracts a single cell from the `DMSWARM` object, returns it as a single cell `DMSWARM`. 173920f4b53cSBarry Smith The cell `DM` is filtered for fields of that cell, and the filtered `DM` is used as the cell `DM` of the new swarm object. 1740d0c080abSJoseph Pusztay 1741d0c080abSJoseph Pusztay Noncollective 1742d0c080abSJoseph Pusztay 174360225df5SJacob Faibussowitsch Input Parameters: 174420f4b53cSBarry Smith + sw - the `DMSWARM` 17455627991aSBarry Smith . cellID - the integer id of the cell to be extracted and filtered 174620f4b53cSBarry Smith - cellswarm - The `DMSWARM` to receive the cell 1747d0c080abSJoseph Pusztay 1748d0c080abSJoseph Pusztay Level: beginner 1749d0c080abSJoseph Pusztay 17505627991aSBarry Smith Notes: 175120f4b53cSBarry Smith This presently only supports `DMSWARM_PIC` type 17525627991aSBarry Smith 175320f4b53cSBarry Smith Should be restored with `DMSwarmRestoreCellSwarm()` 1754d0c080abSJoseph Pusztay 175520f4b53cSBarry Smith Changes to this cell of the swarm will be lost if they are made prior to restoring this cell. 175620f4b53cSBarry Smith 175720f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmRestoreCellSwarm()` 1758d0c080abSJoseph Pusztay @*/ 1759cc4c1da9SBarry Smith PetscErrorCode DMSwarmGetCellSwarm(DM sw, PetscInt cellID, DM cellswarm) 1760d71ae5a4SJacob Faibussowitsch { 1761d0c080abSJoseph Pusztay DM_Swarm *original = (DM_Swarm *)sw->data; 1762d0c080abSJoseph Pusztay DMLabel label; 1763d0c080abSJoseph Pusztay DM dmc, subdmc; 1764d0c080abSJoseph Pusztay PetscInt *pids, particles, dim; 1765d0c080abSJoseph Pusztay 1766d0c080abSJoseph Pusztay PetscFunctionBegin; 1767d0c080abSJoseph Pusztay /* Configure new swarm */ 17689566063dSJacob Faibussowitsch PetscCall(DMSetType(cellswarm, DMSWARM)); 17699566063dSJacob Faibussowitsch PetscCall(DMGetDimension(sw, &dim)); 17709566063dSJacob Faibussowitsch PetscCall(DMSetDimension(cellswarm, dim)); 17719566063dSJacob Faibussowitsch PetscCall(DMSwarmSetType(cellswarm, DMSWARM_PIC)); 1772d0c080abSJoseph Pusztay /* Destroy the unused, unconfigured data bucket to prevent stragglers in memory */ 17739566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketDestroy(&((DM_Swarm *)cellswarm->data)->db)); 17749566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetAccess(sw)); 17759566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetNumberOfPointsPerCell(sw, cellID, &particles)); 17769566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetPointsPerCell(sw, cellID, &particles, &pids)); 17779566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketCreateFromSubset(original->db, particles, pids, &((DM_Swarm *)cellswarm->data)->db)); 17789566063dSJacob Faibussowitsch PetscCall(DMSwarmSortRestoreAccess(sw)); 17799566063dSJacob Faibussowitsch PetscCall(PetscFree(pids)); 17809566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(sw, &dmc)); 17819566063dSJacob Faibussowitsch PetscCall(DMLabelCreate(PetscObjectComm((PetscObject)sw), "singlecell", &label)); 17829566063dSJacob Faibussowitsch PetscCall(DMAddLabel(dmc, label)); 17839566063dSJacob Faibussowitsch PetscCall(DMLabelSetValue(label, cellID, 1)); 178430cbcd5dSksagiyam PetscCall(DMPlexFilter(dmc, label, 1, PETSC_FALSE, PETSC_FALSE, NULL, &subdmc)); 17859566063dSJacob Faibussowitsch PetscCall(DMSwarmSetCellDM(cellswarm, subdmc)); 17869566063dSJacob Faibussowitsch PetscCall(DMLabelDestroy(&label)); 17873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1788d0c080abSJoseph Pusztay } 1789d0c080abSJoseph Pusztay 1790cc4c1da9SBarry Smith /*@ 179120f4b53cSBarry Smith DMSwarmRestoreCellSwarm - Restores a `DMSWARM` object obtained with `DMSwarmGetCellSwarm()`. All fields are copied back into the parent swarm. 1792d0c080abSJoseph Pusztay 1793d0c080abSJoseph Pusztay Noncollective 1794d0c080abSJoseph Pusztay 179560225df5SJacob Faibussowitsch Input Parameters: 179620f4b53cSBarry Smith + sw - the parent `DMSWARM` 1797d0c080abSJoseph Pusztay . cellID - the integer id of the cell to be copied back into the parent swarm 1798d0c080abSJoseph Pusztay - cellswarm - the cell swarm object 1799d0c080abSJoseph Pusztay 1800d0c080abSJoseph Pusztay Level: beginner 1801d0c080abSJoseph Pusztay 18025627991aSBarry Smith Note: 180320f4b53cSBarry Smith This only supports `DMSWARM_PIC` types of `DMSWARM`s 1804d0c080abSJoseph Pusztay 180520f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmGetCellSwarm()` 1806d0c080abSJoseph Pusztay @*/ 1807cc4c1da9SBarry Smith PetscErrorCode DMSwarmRestoreCellSwarm(DM sw, PetscInt cellID, DM cellswarm) 1808d71ae5a4SJacob Faibussowitsch { 1809d0c080abSJoseph Pusztay DM dmc; 1810d0c080abSJoseph Pusztay PetscInt *pids, particles, p; 1811d0c080abSJoseph Pusztay 1812d0c080abSJoseph Pusztay PetscFunctionBegin; 18139566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetAccess(sw)); 18149566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetPointsPerCell(sw, cellID, &particles, &pids)); 18159566063dSJacob Faibussowitsch PetscCall(DMSwarmSortRestoreAccess(sw)); 1816d0c080abSJoseph Pusztay /* Pointwise copy of each particle based on pid. The parent swarm may not be altered during this process. */ 181748a46eb9SPierre Jolivet for (p = 0; p < particles; ++p) PetscCall(DMSwarmDataBucketCopyPoint(((DM_Swarm *)cellswarm->data)->db, pids[p], ((DM_Swarm *)sw->data)->db, pids[p])); 1818d0c080abSJoseph Pusztay /* Free memory, destroy cell dm */ 18199566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(cellswarm, &dmc)); 18209566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmc)); 18219566063dSJacob Faibussowitsch PetscCall(PetscFree(pids)); 18223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1823d0c080abSJoseph Pusztay } 1824d0c080abSJoseph Pusztay 1825d0c080abSJoseph Pusztay PETSC_INTERN PetscErrorCode DMClone_Swarm(DM, DM *); 1826d0c080abSJoseph Pusztay 1827d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMInitialize_Swarm(DM sw) 1828d71ae5a4SJacob Faibussowitsch { 1829d0c080abSJoseph Pusztay PetscFunctionBegin; 1830d0c080abSJoseph Pusztay sw->dim = 0; 1831d0c080abSJoseph Pusztay sw->ops->view = DMView_Swarm; 1832d0c080abSJoseph Pusztay sw->ops->load = NULL; 1833d0c080abSJoseph Pusztay sw->ops->setfromoptions = NULL; 1834d0c080abSJoseph Pusztay sw->ops->clone = DMClone_Swarm; 1835d0c080abSJoseph Pusztay sw->ops->setup = DMSetup_Swarm; 1836d0c080abSJoseph Pusztay sw->ops->createlocalsection = NULL; 1837adc21957SMatthew G. Knepley sw->ops->createsectionpermutation = NULL; 1838d0c080abSJoseph Pusztay sw->ops->createdefaultconstraints = NULL; 1839d0c080abSJoseph Pusztay sw->ops->createglobalvector = DMCreateGlobalVector_Swarm; 1840d0c080abSJoseph Pusztay sw->ops->createlocalvector = DMCreateLocalVector_Swarm; 1841d0c080abSJoseph Pusztay sw->ops->getlocaltoglobalmapping = NULL; 1842d0c080abSJoseph Pusztay sw->ops->createfieldis = NULL; 1843d0c080abSJoseph Pusztay sw->ops->createcoordinatedm = NULL; 1844d0c080abSJoseph Pusztay sw->ops->getcoloring = NULL; 1845d0c080abSJoseph Pusztay sw->ops->creatematrix = DMCreateMatrix_Swarm; 1846d0c080abSJoseph Pusztay sw->ops->createinterpolation = NULL; 1847d0c080abSJoseph Pusztay sw->ops->createinjection = NULL; 1848d0c080abSJoseph Pusztay sw->ops->createmassmatrix = DMCreateMassMatrix_Swarm; 1849d0c080abSJoseph Pusztay sw->ops->refine = NULL; 1850d0c080abSJoseph Pusztay sw->ops->coarsen = NULL; 1851d0c080abSJoseph Pusztay sw->ops->refinehierarchy = NULL; 1852d0c080abSJoseph Pusztay sw->ops->coarsenhierarchy = NULL; 1853d0c080abSJoseph Pusztay sw->ops->globaltolocalbegin = NULL; 1854d0c080abSJoseph Pusztay sw->ops->globaltolocalend = NULL; 1855d0c080abSJoseph Pusztay sw->ops->localtoglobalbegin = NULL; 1856d0c080abSJoseph Pusztay sw->ops->localtoglobalend = NULL; 1857d0c080abSJoseph Pusztay sw->ops->destroy = DMDestroy_Swarm; 1858d0c080abSJoseph Pusztay sw->ops->createsubdm = NULL; 1859d0c080abSJoseph Pusztay sw->ops->getdimpoints = NULL; 1860d0c080abSJoseph Pusztay sw->ops->locatepoints = NULL; 18613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1862d0c080abSJoseph Pusztay } 1863d0c080abSJoseph Pusztay 1864d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode DMClone_Swarm(DM dm, DM *newdm) 1865d71ae5a4SJacob Faibussowitsch { 1866d0c080abSJoseph Pusztay DM_Swarm *swarm = (DM_Swarm *)dm->data; 1867d0c080abSJoseph Pusztay 1868d0c080abSJoseph Pusztay PetscFunctionBegin; 1869d0c080abSJoseph Pusztay swarm->refct++; 1870d0c080abSJoseph Pusztay (*newdm)->data = swarm; 18719566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)*newdm, DMSWARM)); 18729566063dSJacob Faibussowitsch PetscCall(DMInitialize_Swarm(*newdm)); 18732e56dffeSJoe Pusztay (*newdm)->dim = dm->dim; 18743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1875d0c080abSJoseph Pusztay } 1876d0c080abSJoseph Pusztay 1877d3a51819SDave May /*MC 1878*0b4b7b1cSBarry Smith DMSWARM = "swarm" - A `DM` object for particle methods, such as particle-in-cell (PIC), in which the underlying 1879*0b4b7b1cSBarry Smith data is both (i) dynamic in length, (ii) and of arbitrary data type. 1880d3a51819SDave May 188120f4b53cSBarry Smith Level: intermediate 188220f4b53cSBarry Smith 188320f4b53cSBarry Smith Notes: 1884*0b4b7b1cSBarry Smith User data can be represented by `DMSWARM` through a registering "fields" which are to be stored on particles. 188562741f57SDave May To register a field, the user must provide; 188662741f57SDave May (a) a unique name; 188762741f57SDave May (b) the data type (or size in bytes); 188862741f57SDave May (c) the block size of the data. 1889d3a51819SDave May 1890d3a51819SDave May For example, suppose the application requires a unique id, energy, momentum and density to be stored 1891c215a47eSMatthew Knepley on a set of particles. Then the following code could be used 189220f4b53cSBarry Smith .vb 189320f4b53cSBarry Smith DMSwarmInitializeFieldRegister(dm) 189420f4b53cSBarry Smith DMSwarmRegisterPetscDatatypeField(dm,"uid",1,PETSC_LONG); 189520f4b53cSBarry Smith DMSwarmRegisterPetscDatatypeField(dm,"energy",1,PETSC_REAL); 189620f4b53cSBarry Smith DMSwarmRegisterPetscDatatypeField(dm,"momentum",3,PETSC_REAL); 189720f4b53cSBarry Smith DMSwarmRegisterPetscDatatypeField(dm,"density",1,PETSC_FLOAT); 189820f4b53cSBarry Smith DMSwarmFinalizeFieldRegister(dm) 189920f4b53cSBarry Smith .ve 1900d3a51819SDave May 190120f4b53cSBarry Smith The fields represented by `DMSWARM` are dynamic and can be re-sized at any time. 1902*0b4b7b1cSBarry Smith The only restriction imposed by `DMSWARM` is that all fields contain the same number of particles. 1903d3a51819SDave May 1904d3a51819SDave May To support particle methods, "migration" techniques are provided. These methods migrate data 19055627991aSBarry Smith between ranks. 1906d3a51819SDave May 190720f4b53cSBarry Smith `DMSWARM` supports the methods `DMCreateGlobalVector()` and `DMCreateLocalVector()`. 190820f4b53cSBarry Smith As a `DMSWARM` may internally define and store values of different data types, 190920f4b53cSBarry Smith before calling `DMCreateGlobalVector()` or `DMCreateLocalVector()`, the user must inform `DMSWARM` which 1910*0b4b7b1cSBarry Smith fields should be used to define a `Vec` object via `DMSwarmVectorDefineField()` 1911c215a47eSMatthew Knepley The specified field can be changed at any time - thereby permitting vectors 1912c215a47eSMatthew Knepley compatible with different fields to be created. 1913d3a51819SDave May 1914*0b4b7b1cSBarry Smith A dual representation of fields in the `DMSWARM` and a Vec object is permitted via `DMSwarmCreateGlobalVectorFromField()` 1915*0b4b7b1cSBarry Smith Here the data defining the field in the `DMSWARM` is shared with a `Vec`. 1916d3a51819SDave May This is inherently unsafe if you alter the size of the field at any time between 191720f4b53cSBarry Smith calls to `DMSwarmCreateGlobalVectorFromField()` and `DMSwarmDestroyGlobalVectorFromField()`. 191820f4b53cSBarry Smith If the local size of the `DMSWARM` does not match the local size of the global vector 191920f4b53cSBarry Smith when `DMSwarmDestroyGlobalVectorFromField()` is called, an error is thrown. 1920d3a51819SDave May 1921*0b4b7b1cSBarry Smith Additional high-level support is provided for Particle-In-Cell methods. Refer to `DMSwarmSetType()`. 192262741f57SDave May 1923*0b4b7b1cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMType`, `DMCreate()`, `DMSetType()`, `DMSwarmSetType()`, `DMSwarmType`, `DMSwarmCreateGlobalVectorFromField()`, 1924*0b4b7b1cSBarry Smith `DMCreateGlobalVector()`, `DMCreateLocalVector()` 1925d3a51819SDave May M*/ 1926cc4c1da9SBarry Smith 1927d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMCreate_Swarm(DM dm) 1928d71ae5a4SJacob Faibussowitsch { 192957795646SDave May DM_Swarm *swarm; 193057795646SDave May 193157795646SDave May PetscFunctionBegin; 193257795646SDave May PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 19334dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&swarm)); 1934f0cdbbbaSDave May dm->data = swarm; 19359566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketCreate(&swarm->db)); 19369566063dSJacob Faibussowitsch PetscCall(DMSwarmInitializeFieldRegister(dm)); 1937d0c080abSJoseph Pusztay swarm->refct = 1; 1938b5bcf523SDave May swarm->vec_field_set = PETSC_FALSE; 19393454631fSDave May swarm->issetup = PETSC_FALSE; 1940480eef7bSDave May swarm->swarm_type = DMSWARM_BASIC; 1941480eef7bSDave May swarm->migrate_type = DMSWARM_MIGRATE_BASIC; 1942480eef7bSDave May swarm->collect_type = DMSWARM_COLLECT_BASIC; 194340c453e9SDave May swarm->migrate_error_on_missing_point = PETSC_FALSE; 1944f0cdbbbaSDave May swarm->dmcell = NULL; 1945f0cdbbbaSDave May swarm->collect_view_active = PETSC_FALSE; 1946f0cdbbbaSDave May swarm->collect_view_reset_nlocal = -1; 19479566063dSJacob Faibussowitsch PetscCall(DMInitialize_Swarm(dm)); 19482e956fe4SStefano Zampini if (SwarmDataFieldId == -1) PetscCall(PetscObjectComposedDataRegister(&SwarmDataFieldId)); 19493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 195057795646SDave May } 1951