10e2ec84fSDave May #include <petsc/private/dmswarmimpl.h> /*I "petscdmswarm.h" I*/ 20e2ec84fSDave May #include <petscsf.h> 3b799feefSDave May #include <petscdmda.h> 4b799feefSDave May #include <petscdmplex.h> 535b38c60SMatthew G. Knepley #include <petscdt.h> 6279f676cSBarry Smith #include "../src/dm/impls/swarm/data_bucket.h" 70e2ec84fSDave May 835b38c60SMatthew G. Knepley #include <petsc/private/petscfeimpl.h> /* For CoordinatesRefToReal() */ 935b38c60SMatthew G. Knepley 1019307e5cSMatthew G. Knepley PetscClassId DMSWARMCELLDM_CLASSID; 1119307e5cSMatthew G. Knepley 1219307e5cSMatthew G. Knepley /*@ 1319307e5cSMatthew G. Knepley DMSwarmCellDMDestroy - destroy a `DMSwarmCellDM` 1419307e5cSMatthew G. Knepley 1519307e5cSMatthew G. Knepley Collective 1619307e5cSMatthew G. Knepley 1719307e5cSMatthew G. Knepley Input Parameter: 1819307e5cSMatthew G. Knepley . celldm - address of `DMSwarmCellDM` 1919307e5cSMatthew G. Knepley 2019307e5cSMatthew G. Knepley Level: advanced 2119307e5cSMatthew G. Knepley 2219307e5cSMatthew G. Knepley .seealso: `DMSwarmCellDM`, `DMSwarmCellDMCreate()` 2319307e5cSMatthew G. Knepley @*/ 2419307e5cSMatthew G. Knepley PetscErrorCode DMSwarmCellDMDestroy(DMSwarmCellDM *celldm) 2519307e5cSMatthew G. Knepley { 2619307e5cSMatthew G. Knepley PetscFunctionBegin; 2719307e5cSMatthew G. Knepley if (!*celldm) PetscFunctionReturn(PETSC_SUCCESS); 2819307e5cSMatthew G. Knepley PetscValidHeaderSpecific(*celldm, DMSWARMCELLDM_CLASSID, 1); 2919307e5cSMatthew G. Knepley if (--((PetscObject)*celldm)->refct > 0) { 3019307e5cSMatthew G. Knepley *celldm = NULL; 3119307e5cSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 3219307e5cSMatthew G. Knepley } 3319307e5cSMatthew G. Knepley PetscTryTypeMethod(*celldm, destroy); 3419307e5cSMatthew G. Knepley for (PetscInt f = 0; f < (*celldm)->Nf; ++f) PetscCall(PetscFree((*celldm)->dmFields[f])); 3519307e5cSMatthew G. Knepley PetscCall(PetscFree((*celldm)->dmFields)); 3619307e5cSMatthew G. Knepley for (PetscInt f = 0; f < (*celldm)->Nfc; ++f) PetscCall(PetscFree((*celldm)->coordFields[f])); 3719307e5cSMatthew G. Knepley PetscCall(PetscFree((*celldm)->coordFields)); 3819307e5cSMatthew G. Knepley PetscCall(PetscFree((*celldm)->cellid)); 3919307e5cSMatthew G. Knepley PetscCall(DMSwarmSortDestroy(&(*celldm)->sort)); 4019307e5cSMatthew G. Knepley PetscCall(DMDestroy(&(*celldm)->dm)); 4119307e5cSMatthew G. Knepley PetscCall(PetscHeaderDestroy(celldm)); 4219307e5cSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 4319307e5cSMatthew G. Knepley } 4419307e5cSMatthew G. Knepley 4519307e5cSMatthew G. Knepley /*@ 4619307e5cSMatthew G. Knepley DMSwarmCellDMView - view a `DMSwarmCellDM` 4719307e5cSMatthew G. Knepley 4819307e5cSMatthew G. Knepley Collective 4919307e5cSMatthew G. Knepley 5019307e5cSMatthew G. Knepley Input Parameters: 5119307e5cSMatthew G. Knepley + celldm - `DMSwarmCellDM` 5219307e5cSMatthew G. Knepley - viewer - viewer to display field, for example `PETSC_VIEWER_STDOUT_WORLD` 5319307e5cSMatthew G. Knepley 5419307e5cSMatthew G. Knepley Level: advanced 5519307e5cSMatthew G. Knepley 5619307e5cSMatthew G. Knepley .seealso: `DMSwarmCellDM`, `DMSwarmCellDMCreate()` 5719307e5cSMatthew G. Knepley @*/ 5819307e5cSMatthew G. Knepley PetscErrorCode DMSwarmCellDMView(DMSwarmCellDM celldm, PetscViewer viewer) 5919307e5cSMatthew G. Knepley { 609f196a02SMartin Diehl PetscBool isascii; 6119307e5cSMatthew G. Knepley 6219307e5cSMatthew G. Knepley PetscFunctionBegin; 6319307e5cSMatthew G. Knepley PetscValidHeaderSpecific(celldm, DMSWARMCELLDM_CLASSID, 1); 6419307e5cSMatthew G. Knepley if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)celldm), &viewer)); 6519307e5cSMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 6619307e5cSMatthew G. Knepley PetscCheckSameComm(celldm, 1, viewer, 2); 679f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 689f196a02SMartin Diehl if (isascii) { 6919307e5cSMatthew G. Knepley PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)celldm, viewer)); 7019307e5cSMatthew G. Knepley PetscCall(PetscViewerASCIIPushTab(viewer)); 7119307e5cSMatthew G. Knepley PetscCall(PetscViewerASCIIPrintf(viewer, "solution field%s:", celldm->Nf > 1 ? "s" : "")); 7219307e5cSMatthew G. Knepley for (PetscInt f = 0; f < celldm->Nf; ++f) PetscCall(PetscViewerASCIIPrintf(viewer, " %s", celldm->dmFields[f])); 7319307e5cSMatthew G. Knepley PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 7419307e5cSMatthew G. Knepley PetscCall(PetscViewerASCIIPrintf(viewer, "coordinate field%s:", celldm->Nfc > 1 ? "s" : "")); 7519307e5cSMatthew G. Knepley for (PetscInt f = 0; f < celldm->Nfc; ++f) PetscCall(PetscViewerASCIIPrintf(viewer, " %s", celldm->coordFields[f])); 7619307e5cSMatthew G. Knepley PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 7719307e5cSMatthew G. Knepley PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_DEFAULT)); 7819307e5cSMatthew G. Knepley PetscCall(DMView(celldm->dm, viewer)); 7919307e5cSMatthew G. Knepley PetscCall(PetscViewerPopFormat(viewer)); 8019307e5cSMatthew G. Knepley } 8119307e5cSMatthew G. Knepley PetscTryTypeMethod(celldm, view, viewer); 829f196a02SMartin Diehl if (isascii) PetscCall(PetscViewerASCIIPopTab(viewer)); 8319307e5cSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 8419307e5cSMatthew G. Knepley } 8519307e5cSMatthew G. Knepley 8619307e5cSMatthew G. Knepley /*@ 8719307e5cSMatthew G. Knepley DMSwarmCellDMGetDM - Returns the background `DM` for the `DMSwarm` 8819307e5cSMatthew G. Knepley 8919307e5cSMatthew G. Knepley Not Collective 9019307e5cSMatthew G. Knepley 9119307e5cSMatthew G. Knepley Input Parameter: 9219307e5cSMatthew G. Knepley . celldm - The `DMSwarmCellDM` object 9319307e5cSMatthew G. Knepley 9419307e5cSMatthew G. Knepley Output Parameter: 9519307e5cSMatthew G. Knepley . dm - The `DM` object 9619307e5cSMatthew G. Knepley 9719307e5cSMatthew G. Knepley Level: intermediate 9819307e5cSMatthew G. Knepley 9919307e5cSMatthew G. Knepley .seealso: `DMSwarmCellDM`, `DM`, `DMSwarmSetCellDM()` 10019307e5cSMatthew G. Knepley @*/ 10119307e5cSMatthew G. Knepley PetscErrorCode DMSwarmCellDMGetDM(DMSwarmCellDM celldm, DM *dm) 10219307e5cSMatthew G. Knepley { 10319307e5cSMatthew G. Knepley PetscFunctionBegin; 10419307e5cSMatthew G. Knepley PetscValidHeaderSpecific(celldm, DMSWARMCELLDM_CLASSID, 1); 10519307e5cSMatthew G. Knepley PetscAssertPointer(dm, 2); 10619307e5cSMatthew G. Knepley *dm = celldm->dm; 10719307e5cSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 10819307e5cSMatthew G. Knepley } 10919307e5cSMatthew G. Knepley 11019307e5cSMatthew G. Knepley /*@C 11119307e5cSMatthew G. Knepley DMSwarmCellDMGetFields - Returns the `DM` fields for the `DMSwarm` 11219307e5cSMatthew G. Knepley 11319307e5cSMatthew G. Knepley Not Collective 11419307e5cSMatthew G. Knepley 11519307e5cSMatthew G. Knepley Input Parameter: 11619307e5cSMatthew G. Knepley . celldm - The `DMSwarmCellDM` object 11719307e5cSMatthew G. Knepley 11819307e5cSMatthew G. Knepley Output Parameters: 11919307e5cSMatthew G. Knepley + Nf - The number of fields 12019307e5cSMatthew G. Knepley - names - The array of field names in the `DMSWARM` 12119307e5cSMatthew G. Knepley 12219307e5cSMatthew G. Knepley Level: intermediate 12319307e5cSMatthew G. Knepley 12419307e5cSMatthew G. Knepley .seealso: `DMSwarmCellDM`, `DM`, `DMSwarmSetCellDM()` 12519307e5cSMatthew G. Knepley @*/ 12619307e5cSMatthew G. Knepley PetscErrorCode DMSwarmCellDMGetFields(DMSwarmCellDM celldm, PetscInt *Nf, const char **names[]) 12719307e5cSMatthew G. Knepley { 12819307e5cSMatthew G. Knepley PetscFunctionBegin; 12919307e5cSMatthew G. Knepley PetscValidHeaderSpecific(celldm, DMSWARMCELLDM_CLASSID, 1); 13019307e5cSMatthew G. Knepley if (Nf) { 13119307e5cSMatthew G. Knepley PetscAssertPointer(Nf, 2); 13219307e5cSMatthew G. Knepley *Nf = celldm->Nf; 13319307e5cSMatthew G. Knepley } 13419307e5cSMatthew G. Knepley if (names) { 13519307e5cSMatthew G. Knepley PetscAssertPointer(names, 3); 13619307e5cSMatthew G. Knepley *names = (const char **)celldm->dmFields; 13719307e5cSMatthew G. Knepley } 13819307e5cSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 13919307e5cSMatthew G. Knepley } 14019307e5cSMatthew G. Knepley 14119307e5cSMatthew G. Knepley /*@C 14219307e5cSMatthew G. Knepley DMSwarmCellDMGetCoordinateFields - Returns the `DM` coordinate fields for the `DMSwarm` 14319307e5cSMatthew G. Knepley 14419307e5cSMatthew G. Knepley Not Collective 14519307e5cSMatthew G. Knepley 14619307e5cSMatthew G. Knepley Input Parameter: 14719307e5cSMatthew G. Knepley . celldm - The `DMSwarmCellDM` object 14819307e5cSMatthew G. Knepley 14919307e5cSMatthew G. Knepley Output Parameters: 15019307e5cSMatthew G. Knepley + Nfc - The number of coordinate fields 15119307e5cSMatthew G. Knepley - names - The array of coordinate field names in the `DMSWARM` 15219307e5cSMatthew G. Knepley 15319307e5cSMatthew G. Knepley Level: intermediate 15419307e5cSMatthew G. Knepley 15519307e5cSMatthew G. Knepley .seealso: `DMSwarmCellDM`, `DM`, `DMSwarmSetCellDM()` 15619307e5cSMatthew G. Knepley @*/ 15719307e5cSMatthew G. Knepley PetscErrorCode DMSwarmCellDMGetCoordinateFields(DMSwarmCellDM celldm, PetscInt *Nfc, const char **names[]) 15819307e5cSMatthew G. Knepley { 15919307e5cSMatthew G. Knepley PetscFunctionBegin; 16019307e5cSMatthew G. Knepley PetscValidHeaderSpecific(celldm, DMSWARMCELLDM_CLASSID, 1); 16119307e5cSMatthew G. Knepley if (Nfc) { 16219307e5cSMatthew G. Knepley PetscAssertPointer(Nfc, 2); 16319307e5cSMatthew G. Knepley *Nfc = celldm->Nfc; 16419307e5cSMatthew G. Knepley } 16519307e5cSMatthew G. Knepley if (names) { 16619307e5cSMatthew G. Knepley PetscAssertPointer(names, 3); 16719307e5cSMatthew G. Knepley *names = (const char **)celldm->coordFields; 16819307e5cSMatthew G. Knepley } 16919307e5cSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 17019307e5cSMatthew G. Knepley } 17119307e5cSMatthew G. Knepley 17219307e5cSMatthew G. Knepley /*@C 17319307e5cSMatthew G. Knepley DMSwarmCellDMGetCellID - Returns the cell id field name for the `DMSwarm` 17419307e5cSMatthew G. Knepley 17519307e5cSMatthew G. Knepley Not Collective 17619307e5cSMatthew G. Knepley 17719307e5cSMatthew G. Knepley Input Parameter: 17819307e5cSMatthew G. Knepley . celldm - The `DMSwarmCellDM` object 17919307e5cSMatthew G. Knepley 18019307e5cSMatthew G. Knepley Output Parameters: 18119307e5cSMatthew G. Knepley . cellid - The cell id field name in the `DMSWARM` 18219307e5cSMatthew G. Knepley 18319307e5cSMatthew G. Knepley Level: intermediate 18419307e5cSMatthew G. Knepley 18519307e5cSMatthew G. Knepley .seealso: `DMSwarmCellDM`, `DM`, `DMSwarmSetCellDM()` 18619307e5cSMatthew G. Knepley @*/ 18719307e5cSMatthew G. Knepley PetscErrorCode DMSwarmCellDMGetCellID(DMSwarmCellDM celldm, const char *cellid[]) 18819307e5cSMatthew G. Knepley { 18919307e5cSMatthew G. Knepley PetscFunctionBegin; 19019307e5cSMatthew G. Knepley PetscValidHeaderSpecific(celldm, DMSWARMCELLDM_CLASSID, 1); 19119307e5cSMatthew G. Knepley PetscAssertPointer(cellid, 2); 19219307e5cSMatthew G. Knepley *cellid = celldm->cellid; 19319307e5cSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 19419307e5cSMatthew G. Knepley } 19519307e5cSMatthew G. Knepley 19619307e5cSMatthew G. Knepley /*@C 19719307e5cSMatthew G. Knepley DMSwarmCellDMGetSort - Returns the sort context over the active `DMSwarmCellDM` for the `DMSwarm` 19819307e5cSMatthew G. Knepley 19919307e5cSMatthew G. Knepley Not Collective 20019307e5cSMatthew G. Knepley 20119307e5cSMatthew G. Knepley Input Parameter: 20219307e5cSMatthew G. Knepley . celldm - The `DMSwarmCellDM` object 20319307e5cSMatthew G. Knepley 20419307e5cSMatthew G. Knepley Output Parameter: 20519307e5cSMatthew G. Knepley . sort - The `DMSwarmSort` object 20619307e5cSMatthew G. Knepley 20719307e5cSMatthew G. Knepley Level: intermediate 20819307e5cSMatthew G. Knepley 20919307e5cSMatthew G. Knepley .seealso: `DMSwarmCellDM`, `DM`, `DMSwarmSetCellDM()` 21019307e5cSMatthew G. Knepley @*/ 21119307e5cSMatthew G. Knepley PetscErrorCode DMSwarmCellDMGetSort(DMSwarmCellDM celldm, DMSwarmSort *sort) 21219307e5cSMatthew G. Knepley { 21319307e5cSMatthew G. Knepley PetscFunctionBegin; 21419307e5cSMatthew G. Knepley PetscValidHeaderSpecific(celldm, DMSWARMCELLDM_CLASSID, 1); 21519307e5cSMatthew G. Knepley PetscAssertPointer(sort, 2); 21619307e5cSMatthew G. Knepley *sort = celldm->sort; 21719307e5cSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 21819307e5cSMatthew G. Knepley } 21919307e5cSMatthew G. Knepley 22019307e5cSMatthew G. Knepley /*@C 22119307e5cSMatthew G. Knepley DMSwarmCellDMSetSort - Sets the sort context over the active `DMSwarmCellDM` for the `DMSwarm` 22219307e5cSMatthew G. Knepley 22319307e5cSMatthew G. Knepley Not Collective 22419307e5cSMatthew G. Knepley 22519307e5cSMatthew G. Knepley Input Parameters: 22619307e5cSMatthew G. Knepley + celldm - The `DMSwarmCellDM` object 22719307e5cSMatthew G. Knepley - sort - The `DMSwarmSort` object 22819307e5cSMatthew G. Knepley 22919307e5cSMatthew G. Knepley Level: intermediate 23019307e5cSMatthew G. Knepley 23119307e5cSMatthew G. Knepley .seealso: `DMSwarmCellDM`, `DM`, `DMSwarmSetCellDM()` 23219307e5cSMatthew G. Knepley @*/ 23319307e5cSMatthew G. Knepley PetscErrorCode DMSwarmCellDMSetSort(DMSwarmCellDM celldm, DMSwarmSort sort) 23419307e5cSMatthew G. Knepley { 23519307e5cSMatthew G. Knepley PetscFunctionBegin; 23619307e5cSMatthew G. Knepley PetscValidHeaderSpecific(celldm, DMSWARMCELLDM_CLASSID, 1); 23719307e5cSMatthew G. Knepley PetscAssertPointer(sort, 2); 23819307e5cSMatthew G. Knepley celldm->sort = sort; 23919307e5cSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 24019307e5cSMatthew G. Knepley } 24119307e5cSMatthew G. Knepley 24219307e5cSMatthew G. Knepley /*@ 24319307e5cSMatthew G. Knepley DMSwarmCellDMGetBlockSize - Returns the total blocksize for the `DM` fields 24419307e5cSMatthew G. Knepley 24519307e5cSMatthew G. Knepley Not Collective 24619307e5cSMatthew G. Knepley 24719307e5cSMatthew G. Knepley Input Parameters: 24819307e5cSMatthew G. Knepley + celldm - The `DMSwarmCellDM` object 24919307e5cSMatthew G. Knepley - sw - The `DMSwarm` object 25019307e5cSMatthew G. Knepley 25119307e5cSMatthew G. Knepley Output Parameter: 25219307e5cSMatthew G. Knepley . bs - The total block size 25319307e5cSMatthew G. Knepley 25419307e5cSMatthew G. Knepley Level: intermediate 25519307e5cSMatthew G. Knepley 25619307e5cSMatthew G. Knepley .seealso: `DMSwarmCellDM`, `DM`, `DMSwarmSetCellDM()` 25719307e5cSMatthew G. Knepley @*/ 25819307e5cSMatthew G. Knepley PetscErrorCode DMSwarmCellDMGetBlockSize(DMSwarmCellDM celldm, DM sw, PetscInt *bs) 25919307e5cSMatthew G. Knepley { 26019307e5cSMatthew G. Knepley PetscFunctionBegin; 26119307e5cSMatthew G. Knepley PetscValidHeaderSpecific(celldm, DMSWARMCELLDM_CLASSID, 1); 26219307e5cSMatthew G. Knepley PetscValidHeaderSpecific(sw, DM_CLASSID, 2); 26319307e5cSMatthew G. Knepley PetscAssertPointer(bs, 3); 26419307e5cSMatthew G. Knepley *bs = 0; 26519307e5cSMatthew G. Knepley for (PetscInt f = 0; f < celldm->Nf; ++f) { 26619307e5cSMatthew G. Knepley PetscInt fbs; 26719307e5cSMatthew G. Knepley 26819307e5cSMatthew G. Knepley PetscCall(DMSwarmGetFieldInfo(sw, celldm->dmFields[f], &fbs, NULL)); 26919307e5cSMatthew G. Knepley *bs += fbs; 27019307e5cSMatthew G. Knepley } 27119307e5cSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 27219307e5cSMatthew G. Knepley } 27319307e5cSMatthew G. Knepley 27419307e5cSMatthew G. Knepley /*@ 27519307e5cSMatthew G. Knepley DMSwarmCellDMCreate - create a `DMSwarmCellDM` 27619307e5cSMatthew G. Knepley 27719307e5cSMatthew G. Knepley Collective 27819307e5cSMatthew G. Knepley 279d6ae5217SJose E. Roman Input Parameters: 28019307e5cSMatthew G. Knepley + dm - The background `DM` for the `DMSwarm` 28119307e5cSMatthew G. Knepley . Nf - The number of swarm fields defined over `dm` 28219307e5cSMatthew G. Knepley . dmFields - The swarm field names for the `dm` fields 28319307e5cSMatthew G. Knepley . Nfc - The number of swarm fields to use for coordinates over `dm` 28419307e5cSMatthew G. Knepley - coordFields - The swarm field names for the `dm` coordinate fields 28519307e5cSMatthew G. Knepley 28619307e5cSMatthew G. Knepley Output Parameter: 28719307e5cSMatthew G. Knepley . celldm - The new `DMSwarmCellDM` 28819307e5cSMatthew G. Knepley 28919307e5cSMatthew G. Knepley Level: advanced 29019307e5cSMatthew G. Knepley 29119307e5cSMatthew G. Knepley .seealso: `DMSwarmCellDM`, `DMSWARM`, `DMSetType()` 29219307e5cSMatthew G. Knepley @*/ 29319307e5cSMatthew G. Knepley PetscErrorCode DMSwarmCellDMCreate(DM dm, PetscInt Nf, const char *dmFields[], PetscInt Nfc, const char *coordFields[], DMSwarmCellDM *celldm) 29419307e5cSMatthew G. Knepley { 29519307e5cSMatthew G. Knepley DMSwarmCellDM b; 29619307e5cSMatthew G. Knepley const char *name; 29719307e5cSMatthew G. Knepley char cellid[PETSC_MAX_PATH_LEN]; 29819307e5cSMatthew G. Knepley 29919307e5cSMatthew G. Knepley PetscFunctionBegin; 30019307e5cSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 30119307e5cSMatthew G. Knepley if (Nf) PetscAssertPointer(dmFields, 3); 30219307e5cSMatthew G. Knepley if (Nfc) PetscAssertPointer(coordFields, 5); 30319307e5cSMatthew G. Knepley PetscCall(DMInitializePackage()); 30419307e5cSMatthew G. Knepley 30519307e5cSMatthew G. Knepley PetscCall(PetscHeaderCreate(b, DMSWARMCELLDM_CLASSID, "DMSwarmCellDM", "Background DM for a Swarm", "DM", PetscObjectComm((PetscObject)dm), DMSwarmCellDMDestroy, DMSwarmCellDMView)); 30619307e5cSMatthew G. Knepley PetscCall(PetscObjectGetName((PetscObject)dm, &name)); 30719307e5cSMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)b, name)); 30819307e5cSMatthew G. Knepley PetscCall(PetscObjectReference((PetscObject)dm)); 30919307e5cSMatthew G. Knepley b->dm = dm; 31019307e5cSMatthew G. Knepley b->Nf = Nf; 31119307e5cSMatthew G. Knepley b->Nfc = Nfc; 31219307e5cSMatthew G. Knepley PetscCall(PetscMalloc1(b->Nf, &b->dmFields)); 31319307e5cSMatthew G. Knepley for (PetscInt f = 0; f < b->Nf; ++f) PetscCall(PetscStrallocpy(dmFields[f], &b->dmFields[f])); 31419307e5cSMatthew G. Knepley PetscCall(PetscMalloc1(b->Nfc, &b->coordFields)); 31519307e5cSMatthew G. Knepley for (PetscInt f = 0; f < b->Nfc; ++f) PetscCall(PetscStrallocpy(coordFields[f], &b->coordFields[f])); 31619307e5cSMatthew G. Knepley PetscCall(PetscSNPrintf(cellid, PETSC_MAX_PATH_LEN, "%s_cellid", name)); 31719307e5cSMatthew G. Knepley PetscCall(PetscStrallocpy(cellid, &b->cellid)); 31819307e5cSMatthew G. Knepley *celldm = b; 31919307e5cSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 32019307e5cSMatthew G. Knepley } 32119307e5cSMatthew G. Knepley 3220e2ec84fSDave May /* Coordinate insertition/addition API */ 323cc4c1da9SBarry Smith /*@ 32420f4b53cSBarry Smith DMSwarmSetPointsUniformCoordinates - Set point coordinates in a `DMSWARM` on a regular (ijk) grid 3250e2ec84fSDave May 32620f4b53cSBarry Smith Collective 3270e2ec84fSDave May 32860225df5SJacob Faibussowitsch Input Parameters: 32919307e5cSMatthew G. Knepley + sw - the `DMSWARM` 3300e2ec84fSDave May . min - minimum coordinate values in the x, y, z directions (array of length dim) 3310e2ec84fSDave May . max - maximum coordinate values in the x, y, z directions (array of length dim) 3320e2ec84fSDave May . npoints - number of points in each spatial direction (array of length dim) 33320f4b53cSBarry Smith - mode - indicates whether to append points to the swarm (`ADD_VALUES`), or over-ride existing points (`INSERT_VALUES`) 3340e2ec84fSDave May 3350e2ec84fSDave May Level: beginner 3360e2ec84fSDave May 3370e2ec84fSDave May Notes: 33820f4b53cSBarry Smith When using mode = `INSERT_VALUES`, this method will reset the number of particles in the `DMSWARM` 339a3b724e8SBarry Smith to be `npoints[0]` x `npoints[1]` (2D) or `npoints[0]` x `npoints[1]` x `npoints[2]` (3D). When using mode = `ADD_VALUES`, 34020f4b53cSBarry Smith new points will be appended to any already existing in the `DMSWARM` 3410e2ec84fSDave May 34220f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmSetType()`, `DMSwarmSetCellDM()`, `DMSwarmType` 3430e2ec84fSDave May @*/ 34419307e5cSMatthew G. Knepley PetscErrorCode DMSwarmSetPointsUniformCoordinates(DM sw, PetscReal min[], PetscReal max[], PetscInt npoints[], InsertMode mode) 345d71ae5a4SJacob Faibussowitsch { 346c448b993SMatthew G. Knepley PetscReal lmin[] = {PETSC_MAX_REAL, PETSC_MAX_REAL, PETSC_MAX_REAL}; 347c448b993SMatthew G. Knepley PetscReal lmax[] = {PETSC_MIN_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL}; 34819307e5cSMatthew G. Knepley PetscInt i, j, k, bs, b, n_estimate, n_curr, n_new_est, p, n_found, Nfc; 3490e2ec84fSDave May const PetscScalar *_coor; 35019307e5cSMatthew G. Knepley DMSwarmCellDM celldm; 35119307e5cSMatthew G. Knepley DM dm; 3520e2ec84fSDave May PetscReal dx[3]; 3532e3d3761SDave May PetscInt _npoints[] = {0, 0, 1}; 3540e2ec84fSDave May Vec pos; 3550e2ec84fSDave May PetscScalar *_pos; 3560e2ec84fSDave May PetscReal *swarm_coor; 3570e2ec84fSDave May PetscInt *swarm_cellid; 3580e2ec84fSDave May PetscSF sfcell = NULL; 3590e2ec84fSDave May const PetscSFNode *LA_sfcell; 36019307e5cSMatthew G. Knepley const char **coordFields, *cellid; 3610e2ec84fSDave May 3620e2ec84fSDave May PetscFunctionBegin; 36319307e5cSMatthew G. Knepley DMSWARMPICVALID(sw); 36419307e5cSMatthew G. Knepley PetscCall(DMSwarmGetCellDMActive(sw, &celldm)); 36519307e5cSMatthew G. Knepley PetscCall(DMSwarmCellDMGetCoordinateFields(celldm, &Nfc, &coordFields)); 36619307e5cSMatthew G. Knepley PetscCheck(Nfc == 1, PetscObjectComm((PetscObject)sw), PETSC_ERR_SUP, "We only support a single coordinate field right now, not %" PetscInt_FMT, Nfc); 36719307e5cSMatthew G. Knepley 36819307e5cSMatthew G. Knepley PetscCall(DMSwarmCellDMGetDM(celldm, &dm)); 36919307e5cSMatthew G. Knepley PetscCall(DMGetLocalBoundingBox(dm, lmin, lmax)); 37019307e5cSMatthew G. Knepley PetscCall(DMGetCoordinateDim(dm, &bs)); 3710e2ec84fSDave May 3720e2ec84fSDave May for (b = 0; b < bs; b++) { 37371844bb8SDave May if (npoints[b] > 1) { 3740e2ec84fSDave May dx[b] = (max[b] - min[b]) / ((PetscReal)(npoints[b] - 1)); 37571844bb8SDave May } else { 37671844bb8SDave May dx[b] = 0.0; 37771844bb8SDave May } 3782e3d3761SDave May _npoints[b] = npoints[b]; 3790e2ec84fSDave May } 3800e2ec84fSDave May 3810e2ec84fSDave May /* determine number of points living in the bounding box */ 3820e2ec84fSDave May n_estimate = 0; 3832e3d3761SDave May for (k = 0; k < _npoints[2]; k++) { 3842e3d3761SDave May for (j = 0; j < _npoints[1]; j++) { 3852e3d3761SDave May for (i = 0; i < _npoints[0]; i++) { 3860e2ec84fSDave May PetscReal xp[] = {0.0, 0.0, 0.0}; 3870e2ec84fSDave May PetscInt ijk[3]; 3880e2ec84fSDave May PetscBool point_inside = PETSC_TRUE; 3890e2ec84fSDave May 3900e2ec84fSDave May ijk[0] = i; 3910e2ec84fSDave May ijk[1] = j; 3920e2ec84fSDave May ijk[2] = k; 393ad540459SPierre Jolivet for (b = 0; b < bs; b++) xp[b] = min[b] + ijk[b] * dx[b]; 3940e2ec84fSDave May for (b = 0; b < bs; b++) { 395c448b993SMatthew G. Knepley if (xp[b] < lmin[b]) point_inside = PETSC_FALSE; 396c448b993SMatthew G. Knepley if (xp[b] > lmax[b]) point_inside = PETSC_FALSE; 3970e2ec84fSDave May } 398ad540459SPierre Jolivet if (point_inside) n_estimate++; 3990e2ec84fSDave May } 4000e2ec84fSDave May } 4010e2ec84fSDave May } 4020e2ec84fSDave May 4030e2ec84fSDave May /* create candidate list */ 4049566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, &pos)); 4059566063dSJacob Faibussowitsch PetscCall(VecSetSizes(pos, bs * n_estimate, PETSC_DECIDE)); 4069566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(pos, bs)); 4079566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(pos)); 4089566063dSJacob Faibussowitsch PetscCall(VecGetArray(pos, &_pos)); 4090e2ec84fSDave May 4100e2ec84fSDave May n_estimate = 0; 4112e3d3761SDave May for (k = 0; k < _npoints[2]; k++) { 4122e3d3761SDave May for (j = 0; j < _npoints[1]; j++) { 4132e3d3761SDave May for (i = 0; i < _npoints[0]; i++) { 4140e2ec84fSDave May PetscReal xp[] = {0.0, 0.0, 0.0}; 4150e2ec84fSDave May PetscInt ijk[3]; 4160e2ec84fSDave May PetscBool point_inside = PETSC_TRUE; 4170e2ec84fSDave May 4180e2ec84fSDave May ijk[0] = i; 4190e2ec84fSDave May ijk[1] = j; 4200e2ec84fSDave May ijk[2] = k; 421ad540459SPierre Jolivet for (b = 0; b < bs; b++) xp[b] = min[b] + ijk[b] * dx[b]; 4220e2ec84fSDave May for (b = 0; b < bs; b++) { 423c448b993SMatthew G. Knepley if (xp[b] < lmin[b]) point_inside = PETSC_FALSE; 424c448b993SMatthew G. Knepley if (xp[b] > lmax[b]) point_inside = PETSC_FALSE; 4250e2ec84fSDave May } 4260e2ec84fSDave May if (point_inside) { 427ad540459SPierre Jolivet for (b = 0; b < bs; b++) _pos[bs * n_estimate + b] = xp[b]; 4280e2ec84fSDave May n_estimate++; 4290e2ec84fSDave May } 4300e2ec84fSDave May } 4310e2ec84fSDave May } 4320e2ec84fSDave May } 4339566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(pos, &_pos)); 4340e2ec84fSDave May 4350e2ec84fSDave May /* locate points */ 43619307e5cSMatthew G. Knepley PetscCall(DMLocatePoints(dm, pos, DM_POINTLOCATION_NONE, &sfcell)); 4379566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell)); 4380e2ec84fSDave May n_found = 0; 4390e2ec84fSDave May for (p = 0; p < n_estimate; p++) { 440ad540459SPierre Jolivet if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) n_found++; 4410e2ec84fSDave May } 4420e2ec84fSDave May 4430e2ec84fSDave May /* adjust size */ 4440e2ec84fSDave May if (mode == ADD_VALUES) { 44519307e5cSMatthew G. Knepley PetscCall(DMSwarmGetLocalSize(sw, &n_curr)); 4460e2ec84fSDave May n_new_est = n_curr + n_found; 44719307e5cSMatthew G. Knepley PetscCall(DMSwarmSetLocalSizes(sw, n_new_est, -1)); 4480e2ec84fSDave May } 4490e2ec84fSDave May if (mode == INSERT_VALUES) { 4500e2ec84fSDave May n_curr = 0; 4510e2ec84fSDave May n_new_est = n_found; 45219307e5cSMatthew G. Knepley PetscCall(DMSwarmSetLocalSizes(sw, n_new_est, -1)); 4530e2ec84fSDave May } 4540e2ec84fSDave May 4550e2ec84fSDave May /* initialize new coords, cell owners, pid */ 45619307e5cSMatthew G. Knepley PetscCall(DMSwarmCellDMGetCellID(celldm, &cellid)); 4579566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pos, &_coor)); 45819307e5cSMatthew G. Knepley PetscCall(DMSwarmGetField(sw, coordFields[0], NULL, NULL, (void **)&swarm_coor)); 45919307e5cSMatthew G. Knepley PetscCall(DMSwarmGetField(sw, cellid, NULL, NULL, (void **)&swarm_cellid)); 4600e2ec84fSDave May n_found = 0; 4610e2ec84fSDave May for (p = 0; p < n_estimate; p++) { 4620e2ec84fSDave May if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) { 463ad540459SPierre Jolivet for (b = 0; b < bs; b++) swarm_coor[bs * (n_curr + n_found) + b] = PetscRealPart(_coor[bs * p + b]); 4640e2ec84fSDave May swarm_cellid[n_curr + n_found] = LA_sfcell[p].index; 4650e2ec84fSDave May n_found++; 4660e2ec84fSDave May } 4670e2ec84fSDave May } 46819307e5cSMatthew G. Knepley PetscCall(DMSwarmRestoreField(sw, cellid, NULL, NULL, (void **)&swarm_cellid)); 46919307e5cSMatthew G. Knepley PetscCall(DMSwarmRestoreField(sw, coordFields[0], NULL, NULL, (void **)&swarm_coor)); 4709566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pos, &_coor)); 4710e2ec84fSDave May 4729566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfcell)); 4739566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pos)); 4743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4750e2ec84fSDave May } 4760e2ec84fSDave May 477cc4c1da9SBarry Smith /*@ 47820f4b53cSBarry Smith DMSwarmSetPointCoordinates - Set point coordinates in a `DMSWARM` from a user defined list 4790e2ec84fSDave May 48020f4b53cSBarry Smith Collective 4810e2ec84fSDave May 48260225df5SJacob Faibussowitsch Input Parameters: 48319307e5cSMatthew G. Knepley + sw - the `DMSWARM` 4840e2ec84fSDave May . npoints - the number of points to insert 4850e2ec84fSDave May . coor - the coordinate values 48620f4b53cSBarry Smith . redundant - if set to `PETSC_TRUE`, it is assumed that `npoints` and `coor` are only valid on rank 0 and should be broadcast to other ranks 48720f4b53cSBarry Smith - mode - indicates whether to append points to the swarm (`ADD_VALUES`), or over-ride existing points (`INSERT_VALUES`) 4880e2ec84fSDave May 4890e2ec84fSDave May Level: beginner 4900e2ec84fSDave May 4910e2ec84fSDave May Notes: 49220f4b53cSBarry Smith If the user has specified `redundant` as `PETSC_FALSE`, the cell `DM` will attempt to locate the coordinates provided by `coor` within 49320f4b53cSBarry Smith its sub-domain. If they any values within `coor` are not located in the sub-domain, they will be ignored and will not get 49420f4b53cSBarry Smith added to the `DMSWARM`. 4950e2ec84fSDave May 49620f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmSetType()`, `DMSwarmSetCellDM()`, `DMSwarmType`, `DMSwarmSetPointsUniformCoordinates()` 4970e2ec84fSDave May @*/ 49819307e5cSMatthew G. Knepley PetscErrorCode DMSwarmSetPointCoordinates(DM sw, PetscInt npoints, PetscReal coor[], PetscBool redundant, InsertMode mode) 499d71ae5a4SJacob Faibussowitsch { 5000e2ec84fSDave May PetscReal gmin[] = {PETSC_MAX_REAL, PETSC_MAX_REAL, PETSC_MAX_REAL}; 5010e2ec84fSDave May PetscReal gmax[] = {PETSC_MIN_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL}; 5020e2ec84fSDave May PetscInt i, N, bs, b, n_estimate, n_curr, n_new_est, p, n_found; 5030e2ec84fSDave May Vec coorlocal; 5040e2ec84fSDave May const PetscScalar *_coor; 50519307e5cSMatthew G. Knepley DMSwarmCellDM celldm; 50619307e5cSMatthew G. Knepley DM dm; 5070e2ec84fSDave May Vec pos; 5080e2ec84fSDave May PetscScalar *_pos; 5090e2ec84fSDave May PetscReal *swarm_coor; 5100e2ec84fSDave May PetscInt *swarm_cellid; 5110e2ec84fSDave May PetscSF sfcell = NULL; 5120e2ec84fSDave May const PetscSFNode *LA_sfcell; 5130e2ec84fSDave May PetscReal *my_coor; 51419307e5cSMatthew G. Knepley PetscInt my_npoints, Nfc; 5150e2ec84fSDave May PetscMPIInt rank; 5160e2ec84fSDave May MPI_Comm comm; 51719307e5cSMatthew G. Knepley const char **coordFields, *cellid; 5180e2ec84fSDave May 5190e2ec84fSDave May PetscFunctionBegin; 52019307e5cSMatthew G. Knepley DMSWARMPICVALID(sw); 52119307e5cSMatthew G. Knepley PetscCall(PetscObjectGetComm((PetscObject)sw, &comm)); 5229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 5230e2ec84fSDave May 52419307e5cSMatthew G. Knepley PetscCall(DMSwarmGetCellDMActive(sw, &celldm)); 52519307e5cSMatthew G. Knepley PetscCall(DMSwarmCellDMGetCoordinateFields(celldm, &Nfc, &coordFields)); 52619307e5cSMatthew G. Knepley PetscCheck(Nfc == 1, PetscObjectComm((PetscObject)sw), PETSC_ERR_SUP, "We only support a single coordinate field right now, not %" PetscInt_FMT, Nfc); 52719307e5cSMatthew G. Knepley 52819307e5cSMatthew G. Knepley PetscCall(DMSwarmCellDMGetDM(celldm, &dm)); 52919307e5cSMatthew G. Knepley PetscCall(DMGetCoordinatesLocal(dm, &coorlocal)); 5309566063dSJacob Faibussowitsch PetscCall(VecGetSize(coorlocal, &N)); 5319566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(coorlocal, &bs)); 5320e2ec84fSDave May N = N / bs; 5339566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(coorlocal, &_coor)); 5340e2ec84fSDave May for (i = 0; i < N; i++) { 5350e2ec84fSDave May for (b = 0; b < bs; b++) { 536a5f152d1SDave May gmin[b] = PetscMin(gmin[b], PetscRealPart(_coor[bs * i + b])); 537a5f152d1SDave May gmax[b] = PetscMax(gmax[b], PetscRealPart(_coor[bs * i + b])); 5380e2ec84fSDave May } 5390e2ec84fSDave May } 5409566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(coorlocal, &_coor)); 5410e2ec84fSDave May 5420e2ec84fSDave May /* broadcast points from rank 0 if requested */ 5430e2ec84fSDave May if (redundant) { 5446497c311SBarry Smith PetscMPIInt imy; 5456497c311SBarry Smith 5460e2ec84fSDave May my_npoints = npoints; 5479566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&my_npoints, 1, MPIU_INT, 0, comm)); 5480e2ec84fSDave May 5490e2ec84fSDave May if (rank > 0) { /* allocate space */ 5509566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs * my_npoints, &my_coor)); 5510e2ec84fSDave May } else { 5520e2ec84fSDave May my_coor = coor; 5530e2ec84fSDave May } 5546497c311SBarry Smith PetscCall(PetscMPIIntCast(bs * my_npoints, &imy)); 5556497c311SBarry Smith PetscCallMPI(MPI_Bcast(my_coor, imy, MPIU_REAL, 0, comm)); 5560e2ec84fSDave May } else { 5570e2ec84fSDave May my_npoints = npoints; 5580e2ec84fSDave May my_coor = coor; 5590e2ec84fSDave May } 5600e2ec84fSDave May 5610e2ec84fSDave May /* determine the number of points living in the bounding box */ 5620e2ec84fSDave May n_estimate = 0; 5630e2ec84fSDave May for (i = 0; i < my_npoints; i++) { 5640e2ec84fSDave May PetscBool point_inside = PETSC_TRUE; 5650e2ec84fSDave May 5660e2ec84fSDave May for (b = 0; b < bs; b++) { 567ad540459SPierre Jolivet if (my_coor[bs * i + b] < gmin[b]) point_inside = PETSC_FALSE; 568ad540459SPierre Jolivet if (my_coor[bs * i + b] > gmax[b]) point_inside = PETSC_FALSE; 5690e2ec84fSDave May } 570ad540459SPierre Jolivet if (point_inside) n_estimate++; 5710e2ec84fSDave May } 5720e2ec84fSDave May 5730e2ec84fSDave May /* create candidate list */ 5749566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, &pos)); 5759566063dSJacob Faibussowitsch PetscCall(VecSetSizes(pos, bs * n_estimate, PETSC_DECIDE)); 5769566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(pos, bs)); 5779566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(pos)); 5789566063dSJacob Faibussowitsch PetscCall(VecGetArray(pos, &_pos)); 5790e2ec84fSDave May 5800e2ec84fSDave May n_estimate = 0; 5810e2ec84fSDave May for (i = 0; i < my_npoints; i++) { 5820e2ec84fSDave May PetscBool point_inside = PETSC_TRUE; 5830e2ec84fSDave May 5840e2ec84fSDave May for (b = 0; b < bs; b++) { 585ad540459SPierre Jolivet if (my_coor[bs * i + b] < gmin[b]) point_inside = PETSC_FALSE; 586ad540459SPierre Jolivet if (my_coor[bs * i + b] > gmax[b]) point_inside = PETSC_FALSE; 5870e2ec84fSDave May } 5880e2ec84fSDave May if (point_inside) { 589ad540459SPierre Jolivet for (b = 0; b < bs; b++) _pos[bs * n_estimate + b] = my_coor[bs * i + b]; 5900e2ec84fSDave May n_estimate++; 5910e2ec84fSDave May } 5920e2ec84fSDave May } 5939566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(pos, &_pos)); 5940e2ec84fSDave May 5950e2ec84fSDave May /* locate points */ 59619307e5cSMatthew G. Knepley PetscCall(DMLocatePoints(dm, pos, DM_POINTLOCATION_NONE, &sfcell)); 5970e2ec84fSDave May 5989566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell)); 5990e2ec84fSDave May n_found = 0; 6000e2ec84fSDave May for (p = 0; p < n_estimate; p++) { 601ad540459SPierre Jolivet if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) n_found++; 6020e2ec84fSDave May } 6030e2ec84fSDave May 6040e2ec84fSDave May /* adjust size */ 6050e2ec84fSDave May if (mode == ADD_VALUES) { 60619307e5cSMatthew G. Knepley PetscCall(DMSwarmGetLocalSize(sw, &n_curr)); 6070e2ec84fSDave May n_new_est = n_curr + n_found; 60819307e5cSMatthew G. Knepley PetscCall(DMSwarmSetLocalSizes(sw, n_new_est, -1)); 6090e2ec84fSDave May } 6100e2ec84fSDave May if (mode == INSERT_VALUES) { 6110e2ec84fSDave May n_curr = 0; 6120e2ec84fSDave May n_new_est = n_found; 61319307e5cSMatthew G. Knepley PetscCall(DMSwarmSetLocalSizes(sw, n_new_est, -1)); 6140e2ec84fSDave May } 6150e2ec84fSDave May 6160e2ec84fSDave May /* initialize new coords, cell owners, pid */ 61719307e5cSMatthew G. Knepley PetscCall(DMSwarmCellDMGetCellID(celldm, &cellid)); 6189566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pos, &_coor)); 61919307e5cSMatthew G. Knepley PetscCall(DMSwarmGetField(sw, coordFields[0], NULL, NULL, (void **)&swarm_coor)); 62019307e5cSMatthew G. Knepley PetscCall(DMSwarmGetField(sw, cellid, NULL, NULL, (void **)&swarm_cellid)); 6210e2ec84fSDave May n_found = 0; 6220e2ec84fSDave May for (p = 0; p < n_estimate; p++) { 6230e2ec84fSDave May if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) { 624ad540459SPierre Jolivet for (b = 0; b < bs; b++) swarm_coor[bs * (n_curr + n_found) + b] = PetscRealPart(_coor[bs * p + b]); 6250e2ec84fSDave May swarm_cellid[n_curr + n_found] = LA_sfcell[p].index; 6260e2ec84fSDave May n_found++; 6270e2ec84fSDave May } 6280e2ec84fSDave May } 62919307e5cSMatthew G. Knepley PetscCall(DMSwarmRestoreField(sw, cellid, NULL, NULL, (void **)&swarm_cellid)); 63019307e5cSMatthew G. Knepley PetscCall(DMSwarmRestoreField(sw, coordFields[0], NULL, NULL, (void **)&swarm_coor)); 6319566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pos, &_coor)); 6320e2ec84fSDave May 6330e2ec84fSDave May if (redundant) { 63448a46eb9SPierre Jolivet if (rank > 0) PetscCall(PetscFree(my_coor)); 6350e2ec84fSDave May } 6369566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfcell)); 6379566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pos)); 6383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6390e2ec84fSDave May } 6400e2ec84fSDave May 6410e2ec84fSDave May extern PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_DA(DM, DM, DMSwarmPICLayoutType, PetscInt); 6420e2ec84fSDave May extern PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX(DM, DM, DMSwarmPICLayoutType, PetscInt); 6430e2ec84fSDave May 644cc4c1da9SBarry Smith /*@ 6450e2ec84fSDave May DMSwarmInsertPointsUsingCellDM - Insert point coordinates within each cell 6460e2ec84fSDave May 64720f4b53cSBarry Smith Not Collective 6480e2ec84fSDave May 64960225df5SJacob Faibussowitsch Input Parameters: 65020f4b53cSBarry Smith + dm - the `DMSWARM` 65120f4b53cSBarry Smith . layout_type - method used to fill each cell with the cell `DM` 6520e2ec84fSDave May - fill_param - parameter controlling how many points per cell are added (the meaning of this parameter is dependent on the layout type) 6530e2ec84fSDave May 6540e2ec84fSDave May Level: beginner 6550e2ec84fSDave May 6560e2ec84fSDave May Notes: 65720f4b53cSBarry Smith The insert method will reset any previous defined points within the `DMSWARM`. 658e7af74c8SDave May 65920f4b53cSBarry Smith When using a `DMDA` both 2D and 3D are supported for all layout types provided you are using `DMDA_ELEMENT_Q1`. 660e7af74c8SDave May 661a4e35b19SJacob Faibussowitsch When using a `DMPLEX` the following case are supported\: 66220f4b53cSBarry Smith .vb 663ea3d7275SDave May (i) DMSWARMPIC_LAYOUT_REGULAR: 2D (triangle), 664ea3d7275SDave May (ii) DMSWARMPIC_LAYOUT_GAUSS: 2D and 3D provided the cell is a tri/tet or a quad/hex, 665ea3d7275SDave May (iii) DMSWARMPIC_LAYOUT_SUBDIVISION: 2D and 3D for quad/hex and 2D tri. 66620f4b53cSBarry Smith .ve 6670e2ec84fSDave May 66820f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmPICLayoutType`, `DMSwarmSetType()`, `DMSwarmSetCellDM()`, `DMSwarmType` 6690e2ec84fSDave May @*/ 670cc4c1da9SBarry Smith PetscErrorCode DMSwarmInsertPointsUsingCellDM(DM dm, DMSwarmPICLayoutType layout_type, PetscInt fill_param) 671d71ae5a4SJacob Faibussowitsch { 6720e2ec84fSDave May DM celldm; 6730e2ec84fSDave May PetscBool isDA, isPLEX; 6740e2ec84fSDave May 6750e2ec84fSDave May PetscFunctionBegin; 6760e2ec84fSDave May DMSWARMPICVALID(dm); 6779566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(dm, &celldm)); 6789566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMDA, &isDA)); 6799566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMPLEX, &isPLEX)); 6800e2ec84fSDave May if (isDA) { 6819566063dSJacob Faibussowitsch PetscCall(private_DMSwarmInsertPointsUsingCellDM_DA(dm, celldm, layout_type, fill_param)); 6820e2ec84fSDave May } else if (isPLEX) { 6839566063dSJacob Faibussowitsch PetscCall(private_DMSwarmInsertPointsUsingCellDM_PLEX(dm, celldm, layout_type, fill_param)); 6840e2ec84fSDave May } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only supported for cell DMs of type DMDA and DMPLEX"); 6853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6860e2ec84fSDave May } 6870e2ec84fSDave May 688431382f2SDave May extern PetscErrorCode private_DMSwarmSetPointCoordinatesCellwise_PLEX(DM, DM, PetscInt, PetscReal *); 689431382f2SDave May 690431382f2SDave May /*@C 691431382f2SDave May DMSwarmSetPointCoordinatesCellwise - Insert point coordinates (defined over the reference cell) within each cell 692431382f2SDave May 69320f4b53cSBarry Smith Not Collective 694431382f2SDave May 69560225df5SJacob Faibussowitsch Input Parameters: 69620f4b53cSBarry Smith + dm - the `DMSWARM` 697431382f2SDave May . npoints - the number of points to insert in each cell 698431382f2SDave May - xi - the coordinates (defined in the local coordinate system for each cell) to insert 699431382f2SDave May 700431382f2SDave May Level: beginner 701431382f2SDave May 702431382f2SDave May Notes: 70320f4b53cSBarry Smith The method will reset any previous defined points within the `DMSWARM`. 70420f4b53cSBarry Smith Only supported for `DMPLEX`. If you are using a `DMDA` it is recommended to either use 70520f4b53cSBarry Smith `DMSwarmInsertPointsUsingCellDM()`, or extract and set the coordinates yourself the following code 70620f4b53cSBarry Smith .vb 70720f4b53cSBarry Smith PetscReal *coor; 708d52c2f21SMatthew G. Knepley const char *coordname; 709d52c2f21SMatthew G. Knepley DMSwarmGetCoordinateField(dm, &coordname); 710d52c2f21SMatthew G. Knepley DMSwarmGetField(dm,coordname,NULL,NULL,(void**)&coor); 71120f4b53cSBarry Smith // user code to define the coordinates here 712d52c2f21SMatthew G. Knepley DMSwarmRestoreField(dm,coordname,NULL,NULL,(void**)&coor); 71320f4b53cSBarry Smith .ve 714e7af74c8SDave May 71520f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmSetCellDM()`, `DMSwarmInsertPointsUsingCellDM()` 716431382f2SDave May @*/ 717cc4c1da9SBarry Smith PetscErrorCode DMSwarmSetPointCoordinatesCellwise(DM dm, PetscInt npoints, PetscReal xi[]) 718d71ae5a4SJacob Faibussowitsch { 719431382f2SDave May DM celldm; 720431382f2SDave May PetscBool isDA, isPLEX; 721431382f2SDave May 7220e2ec84fSDave May PetscFunctionBegin; 723431382f2SDave May DMSWARMPICVALID(dm); 7249566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(dm, &celldm)); 7259566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMDA, &isDA)); 7269566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMPLEX, &isPLEX)); 72728b400f6SJacob Faibussowitsch PetscCheck(!isDA, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only supported for cell DMs of type DMPLEX. Recommended you use DMSwarmInsertPointsUsingCellDM()"); 728*966bd95aSPierre Jolivet PetscCheck(isPLEX, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only supported for cell DMs of type DMDA and DMPLEX"); 7299566063dSJacob Faibussowitsch PetscCall(private_DMSwarmSetPointCoordinatesCellwise_PLEX(dm, celldm, npoints, xi)); 7303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7310e2ec84fSDave May } 732431382f2SDave May 733cc4c1da9SBarry Smith /*@ 734b799feefSDave May DMSwarmCreatePointPerCellCount - Count the number of points within all cells in the cell DM 7350e2ec84fSDave May 73620f4b53cSBarry Smith Not Collective 7370e2ec84fSDave May 73860225df5SJacob Faibussowitsch Input Parameter: 73919307e5cSMatthew G. Knepley . sw - the `DMSWARM` 7400e2ec84fSDave May 74160225df5SJacob Faibussowitsch Output Parameters: 74220f4b53cSBarry Smith + ncells - the number of cells in the cell `DM` (optional argument, pass `NULL` to ignore) 743b799feefSDave May - count - array of length ncells containing the number of points per cell 7440e2ec84fSDave May 7450e2ec84fSDave May Level: beginner 7460e2ec84fSDave May 7470e2ec84fSDave May Notes: 7480e2ec84fSDave May The array count is allocated internally and must be free'd by the user. 7490e2ec84fSDave May 75020f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmSetType()`, `DMSwarmSetCellDM()`, `DMSwarmType` 7510e2ec84fSDave May @*/ 75219307e5cSMatthew G. Knepley PetscErrorCode DMSwarmCreatePointPerCellCount(DM sw, PetscInt *ncells, PetscInt **count) 753d71ae5a4SJacob Faibussowitsch { 75419307e5cSMatthew G. Knepley DMSwarmCellDM celldm; 755b799feefSDave May PetscBool isvalid; 756b799feefSDave May PetscInt nel; 757b799feefSDave May PetscInt *sum; 75819307e5cSMatthew G. Knepley const char *cellid; 759b799feefSDave May 7600e2ec84fSDave May PetscFunctionBegin; 76119307e5cSMatthew G. Knepley PetscCall(DMSwarmSortGetIsValid(sw, &isvalid)); 762b799feefSDave May nel = 0; 763b799feefSDave May if (isvalid) { 764b799feefSDave May PetscInt e; 765b799feefSDave May 76619307e5cSMatthew G. Knepley PetscCall(DMSwarmSortGetSizes(sw, &nel, NULL)); 767b799feefSDave May 7689566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nel, &sum)); 76919307e5cSMatthew G. Knepley for (e = 0; e < nel; e++) PetscCall(DMSwarmSortGetNumberOfPointsPerCell(sw, e, &sum[e])); 770b799feefSDave May } else { 77119307e5cSMatthew G. Knepley DM dm; 772b799feefSDave May PetscBool isda, isplex, isshell; 773b799feefSDave May PetscInt p, npoints; 774b799feefSDave May PetscInt *swarm_cellid; 775b799feefSDave May 776b799feefSDave May /* get the number of cells */ 77719307e5cSMatthew G. Knepley PetscCall(DMSwarmGetCellDMActive(sw, &celldm)); 77819307e5cSMatthew G. Knepley PetscCall(DMSwarmCellDMGetDM(celldm, &dm)); 77919307e5cSMatthew G. Knepley PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMDA, &isda)); 78019307e5cSMatthew G. Knepley PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMPLEX, &isplex)); 78119307e5cSMatthew G. Knepley PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell)); 782b799feefSDave May if (isda) { 783b799feefSDave May PetscInt _nel, _npe; 784b799feefSDave May const PetscInt *_element; 785b799feefSDave May 78619307e5cSMatthew G. Knepley PetscCall(DMDAGetElements(dm, &_nel, &_npe, &_element)); 787b799feefSDave May nel = _nel; 78819307e5cSMatthew G. Knepley PetscCall(DMDARestoreElements(dm, &_nel, &_npe, &_element)); 789b799feefSDave May } else if (isplex) { 790b799feefSDave May PetscInt ps, pe; 791b799feefSDave May 79219307e5cSMatthew G. Knepley PetscCall(DMPlexGetHeightStratum(dm, 0, &ps, &pe)); 793b799feefSDave May nel = pe - ps; 794b799feefSDave May } else if (isshell) { 795b799feefSDave May PetscErrorCode (*method_DMShellGetNumberOfCells)(DM, PetscInt *); 796b799feefSDave May 79719307e5cSMatthew G. Knepley PetscCall(PetscObjectQueryFunction((PetscObject)dm, "DMGetNumberOfCells_C", &method_DMShellGetNumberOfCells)); 798b799feefSDave May if (method_DMShellGetNumberOfCells) { 79919307e5cSMatthew G. Knepley PetscCall(method_DMShellGetNumberOfCells(dm, &nel)); 8009371c9d4SSatish Balay } else 80119307e5cSMatthew G. Knepley SETERRQ(PetscObjectComm((PetscObject)sw), PETSC_ERR_SUP, "Cannot determine the number of cells for the DMSHELL object. User must provide a method via PetscObjectComposeFunction( (PetscObject)shelldm, \"DMGetNumberOfCells_C\", your_function_to_compute_number_of_cells);"); 80219307e5cSMatthew G. Knepley } else SETERRQ(PetscObjectComm((PetscObject)sw), PETSC_ERR_SUP, "Cannot determine the number of cells for a DM not of type DA, PLEX or SHELL"); 803b799feefSDave May 8049566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nel, &sum)); 8059566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(sum, nel)); 80619307e5cSMatthew G. Knepley PetscCall(DMSwarmGetLocalSize(sw, &npoints)); 80719307e5cSMatthew G. Knepley PetscCall(DMSwarmCellDMGetCellID(celldm, &cellid)); 80819307e5cSMatthew G. Knepley PetscCall(DMSwarmGetField(sw, cellid, NULL, NULL, (void **)&swarm_cellid)); 809b799feefSDave May for (p = 0; p < npoints; p++) { 810ad540459SPierre Jolivet if (swarm_cellid[p] != DMLOCATEPOINT_POINT_NOT_FOUND) sum[swarm_cellid[p]]++; 811b799feefSDave May } 81219307e5cSMatthew G. Knepley PetscCall(DMSwarmRestoreField(sw, cellid, NULL, NULL, (void **)&swarm_cellid)); 813b799feefSDave May } 814ad540459SPierre Jolivet if (ncells) *ncells = nel; 815b799feefSDave May *count = sum; 8163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8170e2ec84fSDave May } 81835b38c60SMatthew G. Knepley 81935b38c60SMatthew G. Knepley /*@ 82035b38c60SMatthew G. Knepley DMSwarmGetNumSpecies - Get the number of particle species 82135b38c60SMatthew G. Knepley 82220f4b53cSBarry Smith Not Collective 82335b38c60SMatthew G. Knepley 82460225df5SJacob Faibussowitsch Input Parameter: 82560225df5SJacob Faibussowitsch . sw - the `DMSWARM` 82635b38c60SMatthew G. Knepley 82760225df5SJacob Faibussowitsch Output Parameters: 82835b38c60SMatthew G. Knepley . Ns - the number of species 82935b38c60SMatthew G. Knepley 83035b38c60SMatthew G. Knepley Level: intermediate 83135b38c60SMatthew G. Knepley 83220f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmSetNumSpecies()`, `DMSwarmSetType()`, `DMSwarmType` 83335b38c60SMatthew G. Knepley @*/ 834d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmGetNumSpecies(DM sw, PetscInt *Ns) 835d71ae5a4SJacob Faibussowitsch { 83635b38c60SMatthew G. Knepley DM_Swarm *swarm = (DM_Swarm *)sw->data; 83735b38c60SMatthew G. Knepley 83835b38c60SMatthew G. Knepley PetscFunctionBegin; 83935b38c60SMatthew G. Knepley *Ns = swarm->Ns; 8403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 84135b38c60SMatthew G. Knepley } 84235b38c60SMatthew G. Knepley 84335b38c60SMatthew G. Knepley /*@ 84435b38c60SMatthew G. Knepley DMSwarmSetNumSpecies - Set the number of particle species 84535b38c60SMatthew G. Knepley 84620f4b53cSBarry Smith Not Collective 84735b38c60SMatthew G. Knepley 84860225df5SJacob Faibussowitsch Input Parameters: 84960225df5SJacob Faibussowitsch + sw - the `DMSWARM` 85035b38c60SMatthew G. Knepley - Ns - the number of species 85135b38c60SMatthew G. Knepley 85235b38c60SMatthew G. Knepley Level: intermediate 85335b38c60SMatthew G. Knepley 85420f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmGetNumSpecies()`, `DMSwarmSetType()`, `DMSwarmType` 85535b38c60SMatthew G. Knepley @*/ 856d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmSetNumSpecies(DM sw, PetscInt Ns) 857d71ae5a4SJacob Faibussowitsch { 85835b38c60SMatthew G. Knepley DM_Swarm *swarm = (DM_Swarm *)sw->data; 85935b38c60SMatthew G. Knepley 86035b38c60SMatthew G. Knepley PetscFunctionBegin; 86135b38c60SMatthew G. Knepley swarm->Ns = Ns; 8623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 86335b38c60SMatthew G. Knepley } 86435b38c60SMatthew G. Knepley 86535b38c60SMatthew G. Knepley /*@C 8666c5a79ebSMatthew G. Knepley DMSwarmGetCoordinateFunction - Get the function setting initial particle positions, if it exists 8676c5a79ebSMatthew G. Knepley 86820f4b53cSBarry Smith Not Collective 8696c5a79ebSMatthew G. Knepley 87060225df5SJacob Faibussowitsch Input Parameter: 87160225df5SJacob Faibussowitsch . sw - the `DMSWARM` 8726c5a79ebSMatthew G. Knepley 8736c5a79ebSMatthew G. Knepley Output Parameter: 8748434afd1SBarry Smith . coordFunc - the function setting initial particle positions, or `NULL`, see `PetscSimplePointFn` for the calling sequence 8756c5a79ebSMatthew G. Knepley 8766c5a79ebSMatthew G. Knepley Level: intermediate 8776c5a79ebSMatthew G. Knepley 8788434afd1SBarry Smith .seealso: `DMSWARM`, `DMSwarmSetCoordinateFunction()`, `DMSwarmGetVelocityFunction()`, `DMSwarmInitializeCoordinates()`, `PetscSimplePointFn` 8796c5a79ebSMatthew G. Knepley @*/ 8808434afd1SBarry Smith PetscErrorCode DMSwarmGetCoordinateFunction(DM sw, PetscSimplePointFn **coordFunc) 881d71ae5a4SJacob Faibussowitsch { 8826c5a79ebSMatthew G. Knepley DM_Swarm *swarm = (DM_Swarm *)sw->data; 8836c5a79ebSMatthew G. Knepley 8846c5a79ebSMatthew G. Knepley PetscFunctionBegin; 8856c5a79ebSMatthew G. Knepley PetscValidHeaderSpecific(sw, DM_CLASSID, 1); 8864f572ea9SToby Isaac PetscAssertPointer(coordFunc, 2); 8876c5a79ebSMatthew G. Knepley *coordFunc = swarm->coordFunc; 8883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8896c5a79ebSMatthew G. Knepley } 8906c5a79ebSMatthew G. Knepley 8916c5a79ebSMatthew G. Knepley /*@C 8926c5a79ebSMatthew G. Knepley DMSwarmSetCoordinateFunction - Set the function setting initial particle positions 8936c5a79ebSMatthew G. Knepley 89420f4b53cSBarry Smith Not Collective 8956c5a79ebSMatthew G. Knepley 89660225df5SJacob Faibussowitsch Input Parameters: 89760225df5SJacob Faibussowitsch + sw - the `DMSWARM` 8988434afd1SBarry Smith - coordFunc - the function setting initial particle positions, see `PetscSimplePointFn` for the calling sequence 8996c5a79ebSMatthew G. Knepley 9006c5a79ebSMatthew G. Knepley Level: intermediate 9016c5a79ebSMatthew G. Knepley 9028434afd1SBarry Smith .seealso: `DMSWARM`, `DMSwarmGetCoordinateFunction()`, `DMSwarmSetVelocityFunction()`, `DMSwarmInitializeCoordinates()`, `PetscSimplePointFn` 9036c5a79ebSMatthew G. Knepley @*/ 9048434afd1SBarry Smith PetscErrorCode DMSwarmSetCoordinateFunction(DM sw, PetscSimplePointFn *coordFunc) 905d71ae5a4SJacob Faibussowitsch { 9066c5a79ebSMatthew G. Knepley DM_Swarm *swarm = (DM_Swarm *)sw->data; 9076c5a79ebSMatthew G. Knepley 9086c5a79ebSMatthew G. Knepley PetscFunctionBegin; 9096c5a79ebSMatthew G. Knepley PetscValidHeaderSpecific(sw, DM_CLASSID, 1); 9106c5a79ebSMatthew G. Knepley PetscValidFunction(coordFunc, 2); 9116c5a79ebSMatthew G. Knepley swarm->coordFunc = coordFunc; 9123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9136c5a79ebSMatthew G. Knepley } 9146c5a79ebSMatthew G. Knepley 9156c5a79ebSMatthew G. Knepley /*@C 91660225df5SJacob Faibussowitsch DMSwarmGetVelocityFunction - Get the function setting initial particle velocities, if it exists 9176c5a79ebSMatthew G. Knepley 91820f4b53cSBarry Smith Not Collective 9196c5a79ebSMatthew G. Knepley 92060225df5SJacob Faibussowitsch Input Parameter: 92160225df5SJacob Faibussowitsch . sw - the `DMSWARM` 9226c5a79ebSMatthew G. Knepley 9236c5a79ebSMatthew G. Knepley Output Parameter: 9248434afd1SBarry Smith . velFunc - the function setting initial particle velocities, or `NULL`, see `PetscSimplePointFn` for the calling sequence 9256c5a79ebSMatthew G. Knepley 9266c5a79ebSMatthew G. Knepley Level: intermediate 9276c5a79ebSMatthew G. Knepley 9288434afd1SBarry Smith .seealso: `DMSWARM`, `DMSwarmSetVelocityFunction()`, `DMSwarmGetCoordinateFunction()`, `DMSwarmInitializeVelocities()`, `PetscSimplePointFn` 9296c5a79ebSMatthew G. Knepley @*/ 9308434afd1SBarry Smith PetscErrorCode DMSwarmGetVelocityFunction(DM sw, PetscSimplePointFn **velFunc) 931d71ae5a4SJacob Faibussowitsch { 9326c5a79ebSMatthew G. Knepley DM_Swarm *swarm = (DM_Swarm *)sw->data; 9336c5a79ebSMatthew G. Knepley 9346c5a79ebSMatthew G. Knepley PetscFunctionBegin; 9356c5a79ebSMatthew G. Knepley PetscValidHeaderSpecific(sw, DM_CLASSID, 1); 9364f572ea9SToby Isaac PetscAssertPointer(velFunc, 2); 9376c5a79ebSMatthew G. Knepley *velFunc = swarm->velFunc; 9383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9396c5a79ebSMatthew G. Knepley } 9406c5a79ebSMatthew G. Knepley 9416c5a79ebSMatthew G. Knepley /*@C 9426c5a79ebSMatthew G. Knepley DMSwarmSetVelocityFunction - Set the function setting initial particle velocities 9436c5a79ebSMatthew G. Knepley 94420f4b53cSBarry Smith Not Collective 9456c5a79ebSMatthew G. Knepley 94660225df5SJacob Faibussowitsch Input Parameters: 947a4e35b19SJacob Faibussowitsch + sw - the `DMSWARM` 9488434afd1SBarry Smith - velFunc - the function setting initial particle velocities, see `PetscSimplePointFn` for the calling sequence 9496c5a79ebSMatthew G. Knepley 9506c5a79ebSMatthew G. Knepley Level: intermediate 9516c5a79ebSMatthew G. Knepley 9528434afd1SBarry Smith .seealso: `DMSWARM`, `DMSwarmGetVelocityFunction()`, `DMSwarmSetCoordinateFunction()`, `DMSwarmInitializeVelocities()`, `PetscSimplePointFn` 9536c5a79ebSMatthew G. Knepley @*/ 9548434afd1SBarry Smith PetscErrorCode DMSwarmSetVelocityFunction(DM sw, PetscSimplePointFn *velFunc) 955d71ae5a4SJacob Faibussowitsch { 9566c5a79ebSMatthew G. Knepley DM_Swarm *swarm = (DM_Swarm *)sw->data; 9576c5a79ebSMatthew G. Knepley 9586c5a79ebSMatthew G. Knepley PetscFunctionBegin; 9596c5a79ebSMatthew G. Knepley PetscValidHeaderSpecific(sw, DM_CLASSID, 1); 9606c5a79ebSMatthew G. Knepley PetscValidFunction(velFunc, 2); 9616c5a79ebSMatthew G. Knepley swarm->velFunc = velFunc; 9623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9636c5a79ebSMatthew G. Knepley } 9646c5a79ebSMatthew G. Knepley 9656c5a79ebSMatthew G. Knepley /*@C 96635b38c60SMatthew G. Knepley DMSwarmComputeLocalSize - Compute the local number and distribution of particles based upon a density function 96735b38c60SMatthew G. Knepley 96820f4b53cSBarry Smith Not Collective 96935b38c60SMatthew G. Knepley 97035b38c60SMatthew G. Knepley Input Parameters: 97120f4b53cSBarry Smith + sw - The `DMSWARM` 97235b38c60SMatthew G. Knepley . N - The target number of particles 97335b38c60SMatthew G. Knepley - density - The density field for the particle layout, normalized to unity 97435b38c60SMatthew G. Knepley 97535b38c60SMatthew G. Knepley Level: advanced 97635b38c60SMatthew G. Knepley 97720f4b53cSBarry Smith Note: 97820f4b53cSBarry Smith One particle will be created for each species. 97920f4b53cSBarry Smith 98020f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmComputeLocalSizeFromOptions()` 98135b38c60SMatthew G. Knepley @*/ 982f8662bd6SBarry Smith PetscErrorCode DMSwarmComputeLocalSize(DM sw, PetscInt N, PetscProbFn *density) 983d71ae5a4SJacob Faibussowitsch { 98435b38c60SMatthew G. Knepley DM dm; 98519307e5cSMatthew G. Knepley DMSwarmCellDM celldm; 98635b38c60SMatthew G. Knepley PetscQuadrature quad; 98735b38c60SMatthew G. Knepley const PetscReal *xq, *wq; 988ea7032a0SMatthew G. Knepley PetscReal *n_int; 98919307e5cSMatthew G. Knepley PetscInt *npc_s, *swarm_cellid, Ni; 990ea7032a0SMatthew G. Knepley PetscReal gmin[3], gmax[3], xi0[3]; 991ea7032a0SMatthew G. Knepley PetscInt Ns, cStart, cEnd, c, dim, d, Nq, q, Np = 0, p, s; 99235b38c60SMatthew G. Knepley PetscBool simplex; 99319307e5cSMatthew G. Knepley const char *cellid; 99435b38c60SMatthew G. Knepley 99535b38c60SMatthew G. Knepley PetscFunctionBegin; 9969566063dSJacob Faibussowitsch PetscCall(DMSwarmGetNumSpecies(sw, &Ns)); 9979566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(sw, &dm)); 9989566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 999ea7032a0SMatthew G. Knepley PetscCall(DMGetBoundingBox(dm, gmin, gmax)); 10009566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 10019566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(dm, &simplex)); 10026858538eSMatthew G. Knepley PetscCall(DMGetCoordinatesLocalSetUp(dm)); 10039566063dSJacob Faibussowitsch if (simplex) PetscCall(PetscDTStroudConicalQuadrature(dim, 1, 5, -1.0, 1.0, &quad)); 10049566063dSJacob Faibussowitsch else PetscCall(PetscDTGaussTensorQuadrature(dim, 1, 5, -1.0, 1.0, &quad)); 10059566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, &xq, &wq)); 1006ea7032a0SMatthew G. Knepley PetscCall(PetscCalloc2(Ns, &n_int, (cEnd - cStart) * Ns, &npc_s)); 100735b38c60SMatthew G. Knepley /* Integrate the density function to get the number of particles in each cell */ 100835b38c60SMatthew G. Knepley for (d = 0; d < dim; ++d) xi0[d] = -1.0; 100935b38c60SMatthew G. Knepley for (c = 0; c < cEnd - cStart; ++c) { 101035b38c60SMatthew G. Knepley const PetscInt cell = c + cStart; 1011ea7032a0SMatthew G. Knepley PetscReal v0[3], J[9], invJ[9], detJ, detJp = 2. / (gmax[0] - gmin[0]), xr[3], den; 101235b38c60SMatthew G. Knepley 1013ea7032a0SMatthew G. Knepley /* Have to transform quadrature points/weights to cell domain */ 10149566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFEM(dm, cell, NULL, v0, J, invJ, &detJ)); 1015ea7032a0SMatthew G. Knepley PetscCall(PetscArrayzero(n_int, Ns)); 101635b38c60SMatthew G. Knepley for (q = 0; q < Nq; ++q) { 101735b38c60SMatthew G. Knepley CoordinatesRefToReal(dim, dim, xi0, v0, J, &xq[q * dim], xr); 1018ea7032a0SMatthew G. Knepley /* Have to transform mesh to domain of definition of PDF, [-1, 1], and weight PDF by |J|/2 */ 1019ea7032a0SMatthew G. Knepley xr[0] = detJp * (xr[0] - gmin[0]) - 1.; 1020ea7032a0SMatthew G. Knepley 1021ea7032a0SMatthew G. Knepley for (s = 0; s < Ns; ++s) { 1022ea7032a0SMatthew G. Knepley PetscCall(density(xr, NULL, &den)); 1023ea7032a0SMatthew G. Knepley n_int[s] += (detJp * den) * (detJ * wq[q]) / (PetscReal)Ns; 102435b38c60SMatthew G. Knepley } 1025ea7032a0SMatthew G. Knepley } 1026ea7032a0SMatthew G. Knepley for (s = 0; s < Ns; ++s) { 1027ea7032a0SMatthew G. Knepley Ni = N; 1028f940b0e3Sdanofinn npc_s[c * Ns + s] += (PetscInt)(Ni * n_int[s] + 0.5); // TODO Wish we wrapped round() 1029ea7032a0SMatthew G. Knepley Np += npc_s[c * Ns + s]; 1030ea7032a0SMatthew G. Knepley } 103135b38c60SMatthew G. Knepley } 10329566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&quad)); 10339566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(sw, Np, 0)); 103419307e5cSMatthew G. Knepley PetscCall(DMSwarmGetCellDMActive(sw, &celldm)); 103519307e5cSMatthew G. Knepley PetscCall(DMSwarmCellDMGetCellID(celldm, &cellid)); 103619307e5cSMatthew G. Knepley PetscCall(DMSwarmGetField(sw, cellid, NULL, NULL, (void **)&swarm_cellid)); 103735b38c60SMatthew G. Knepley for (c = 0, p = 0; c < cEnd - cStart; ++c) { 1038ea7032a0SMatthew G. Knepley for (s = 0; s < Ns; ++s) { 103919307e5cSMatthew G. Knepley for (q = 0; q < npc_s[c * Ns + s]; ++q, ++p) swarm_cellid[p] = c; 1040ea7032a0SMatthew G. Knepley } 104135b38c60SMatthew G. Knepley } 104219307e5cSMatthew G. Knepley PetscCall(DMSwarmRestoreField(sw, cellid, NULL, NULL, (void **)&swarm_cellid)); 1043ea7032a0SMatthew G. Knepley PetscCall(PetscFree2(n_int, npc_s)); 10443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 104535b38c60SMatthew G. Knepley } 104635b38c60SMatthew G. Knepley 104735b38c60SMatthew G. Knepley /*@ 104835b38c60SMatthew G. Knepley DMSwarmComputeLocalSizeFromOptions - Compute the local number and distribution of particles based upon a density function determined by options 104935b38c60SMatthew G. Knepley 105020f4b53cSBarry Smith Not Collective 105135b38c60SMatthew G. Knepley 10522fe279fdSBarry Smith Input Parameter: 10532fe279fdSBarry Smith . sw - The `DMSWARM` 105435b38c60SMatthew G. Knepley 105535b38c60SMatthew G. Knepley Level: advanced 105635b38c60SMatthew G. Knepley 105720f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmComputeLocalSize()` 105835b38c60SMatthew G. Knepley @*/ 1059d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmComputeLocalSizeFromOptions(DM sw) 1060d71ae5a4SJacob Faibussowitsch { 1061f8662bd6SBarry Smith PetscProbFn *pdf; 106235b38c60SMatthew G. Knepley const char *prefix; 10636c5a79ebSMatthew G. Knepley char funcname[PETSC_MAX_PATH_LEN]; 10646c5a79ebSMatthew G. Knepley PetscInt *N, Ns, dim, n; 10656c5a79ebSMatthew G. Knepley PetscBool flg; 10666c5a79ebSMatthew G. Knepley PetscMPIInt size, rank; 106735b38c60SMatthew G. Knepley 106835b38c60SMatthew G. Knepley PetscFunctionBegin; 10696c5a79ebSMatthew G. Knepley PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)sw), &size)); 10706c5a79ebSMatthew G. Knepley PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sw), &rank)); 10716c5a79ebSMatthew G. Knepley PetscCall(PetscCalloc1(size, &N)); 1072d0609cedSBarry Smith PetscOptionsBegin(PetscObjectComm((PetscObject)sw), "", "DMSwarm Options", "DMSWARM"); 10736c5a79ebSMatthew G. Knepley n = size; 10746c5a79ebSMatthew G. Knepley PetscCall(PetscOptionsIntArray("-dm_swarm_num_particles", "The target number of particles", "", N, &n, NULL)); 10756c5a79ebSMatthew G. Knepley PetscCall(DMSwarmGetNumSpecies(sw, &Ns)); 10769566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-dm_swarm_num_species", "The number of species", "DMSwarmSetNumSpecies", Ns, &Ns, &flg)); 10779566063dSJacob Faibussowitsch if (flg) PetscCall(DMSwarmSetNumSpecies(sw, Ns)); 10786c5a79ebSMatthew G. Knepley PetscCall(PetscOptionsString("-dm_swarm_coordinate_function", "Function to determine particle coordinates", "DMSwarmSetCoordinateFunction", funcname, funcname, sizeof(funcname), &flg)); 1079d0609cedSBarry Smith PetscOptionsEnd(); 10806c5a79ebSMatthew G. Knepley if (flg) { 10818434afd1SBarry Smith PetscSimplePointFn *coordFunc; 108235b38c60SMatthew G. Knepley 10836c5a79ebSMatthew G. Knepley PetscCall(DMSwarmGetNumSpecies(sw, &Ns)); 10846c5a79ebSMatthew G. Knepley PetscCall(PetscDLSym(NULL, funcname, (void **)&coordFunc)); 10856c5a79ebSMatthew G. Knepley PetscCheck(coordFunc, PetscObjectComm((PetscObject)sw), PETSC_ERR_ARG_WRONG, "Could not locate function %s", funcname); 10866c5a79ebSMatthew G. Knepley PetscCall(DMSwarmGetNumSpecies(sw, &Ns)); 10876c5a79ebSMatthew G. Knepley PetscCall(DMSwarmSetLocalSizes(sw, N[rank] * Ns, 0)); 10886c5a79ebSMatthew G. Knepley PetscCall(DMSwarmSetCoordinateFunction(sw, coordFunc)); 10896c5a79ebSMatthew G. Knepley } else { 10909566063dSJacob Faibussowitsch PetscCall(DMGetDimension(sw, &dim)); 10919566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)sw, &prefix)); 10929566063dSJacob Faibussowitsch PetscCall(PetscProbCreateFromOptions(dim, prefix, "-dm_swarm_coordinate_density", &pdf, NULL, NULL)); 10936c5a79ebSMatthew G. Knepley PetscCall(DMSwarmComputeLocalSize(sw, N[rank], pdf)); 10946c5a79ebSMatthew G. Knepley } 10956c5a79ebSMatthew G. Knepley PetscCall(PetscFree(N)); 10963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 109735b38c60SMatthew G. Knepley } 109835b38c60SMatthew G. Knepley 109935b38c60SMatthew G. Knepley /*@ 110035b38c60SMatthew G. Knepley DMSwarmInitializeCoordinates - Determine the initial coordinates of particles for a PIC method 110135b38c60SMatthew G. Knepley 110220f4b53cSBarry Smith Not Collective 110335b38c60SMatthew G. Knepley 110420f4b53cSBarry Smith Input Parameter: 110520f4b53cSBarry Smith . sw - The `DMSWARM` 110635b38c60SMatthew G. Knepley 110735b38c60SMatthew G. Knepley Level: advanced 110835b38c60SMatthew G. Knepley 110920f4b53cSBarry Smith Note: 111020f4b53cSBarry Smith Currently, we randomly place particles in their assigned cell 111120f4b53cSBarry Smith 111220f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmComputeLocalSize()`, `DMSwarmInitializeVelocities()` 111335b38c60SMatthew G. Knepley @*/ 1114d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmInitializeCoordinates(DM sw) 1115d71ae5a4SJacob Faibussowitsch { 111619307e5cSMatthew G. Knepley DMSwarmCellDM celldm; 11178434afd1SBarry Smith PetscSimplePointFn *coordFunc; 111835b38c60SMatthew G. Knepley PetscScalar *weight; 11196c5a79ebSMatthew G. Knepley PetscReal *x; 112035b38c60SMatthew G. Knepley PetscInt *species; 11216c5a79ebSMatthew G. Knepley void *ctx; 112235b38c60SMatthew G. Knepley PetscBool removePoints = PETSC_TRUE; 112335b38c60SMatthew G. Knepley PetscDataType dtype; 112419307e5cSMatthew G. Knepley PetscInt Nfc, Np, p, Ns, dim, d, bs; 112519307e5cSMatthew G. Knepley const char **coordFields; 112635b38c60SMatthew G. Knepley 112735b38c60SMatthew G. Knepley PetscFunctionBeginUser; 11286c5a79ebSMatthew G. Knepley PetscCall(DMGetDimension(sw, &dim)); 11296c5a79ebSMatthew G. Knepley PetscCall(DMSwarmGetLocalSize(sw, &Np)); 11309566063dSJacob Faibussowitsch PetscCall(DMSwarmGetNumSpecies(sw, &Ns)); 11316c5a79ebSMatthew G. Knepley PetscCall(DMSwarmGetCoordinateFunction(sw, &coordFunc)); 11326c5a79ebSMatthew G. Knepley 113319307e5cSMatthew G. Knepley PetscCall(DMSwarmGetCellDMActive(sw, &celldm)); 113419307e5cSMatthew G. Knepley PetscCall(DMSwarmCellDMGetCoordinateFields(celldm, &Nfc, &coordFields)); 113519307e5cSMatthew G. Knepley PetscCheck(Nfc == 1, PetscObjectComm((PetscObject)sw), PETSC_ERR_SUP, "We only support a single coordinate field right now, not %" PetscInt_FMT, Nfc); 113619307e5cSMatthew G. Knepley 113719307e5cSMatthew G. Knepley PetscCall(DMSwarmGetField(sw, coordFields[0], &bs, &dtype, (void **)&x)); 11386c5a79ebSMatthew G. Knepley PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&weight)); 11396c5a79ebSMatthew G. Knepley PetscCall(DMSwarmGetField(sw, "species", NULL, NULL, (void **)&species)); 11406c5a79ebSMatthew G. Knepley if (coordFunc) { 11416c5a79ebSMatthew G. Knepley PetscCall(DMGetApplicationContext(sw, &ctx)); 11426c5a79ebSMatthew G. Knepley for (p = 0; p < Np; ++p) { 11436c5a79ebSMatthew G. Knepley PetscScalar X[3]; 11446c5a79ebSMatthew G. Knepley 11456c5a79ebSMatthew G. Knepley PetscCall((*coordFunc)(dim, 0., NULL, p, X, ctx)); 11466c5a79ebSMatthew G. Knepley for (d = 0; d < dim; ++d) x[p * dim + d] = PetscRealPart(X[d]); 11476c5a79ebSMatthew G. Knepley weight[p] = 1.0; 11486c5a79ebSMatthew G. Knepley species[p] = p % Ns; 11496c5a79ebSMatthew G. Knepley } 11506c5a79ebSMatthew G. Knepley } else { 11516c5a79ebSMatthew G. Knepley DM dm; 11526c5a79ebSMatthew G. Knepley PetscRandom rnd; 11536c5a79ebSMatthew G. Knepley PetscReal xi0[3]; 11546c5a79ebSMatthew G. Knepley PetscInt cStart, cEnd, c; 11556c5a79ebSMatthew G. Knepley 11569566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(sw, &dm)); 11579566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 1158ea7032a0SMatthew G. Knepley PetscCall(DMGetApplicationContext(sw, &ctx)); 115935b38c60SMatthew G. Knepley 116035b38c60SMatthew G. Knepley /* Set particle position randomly in cell, set weights to 1 */ 11619566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)dm), &rnd)); 11629566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rnd, -1.0, 1.0)); 11639566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rnd)); 11649566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetAccess(sw)); 116535b38c60SMatthew G. Knepley for (d = 0; d < dim; ++d) xi0[d] = -1.0; 116635b38c60SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 116735b38c60SMatthew G. Knepley PetscReal v0[3], J[9], invJ[9], detJ; 116835b38c60SMatthew G. Knepley PetscInt *pidx, Npc, q; 116935b38c60SMatthew G. Knepley 11709566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetPointsPerCell(sw, c, &Npc, &pidx)); 11719566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ)); 117235b38c60SMatthew G. Knepley for (q = 0; q < Npc; ++q) { 117335b38c60SMatthew G. Knepley const PetscInt p = pidx[q]; 117435b38c60SMatthew G. Knepley PetscReal xref[3]; 117535b38c60SMatthew G. Knepley 11769566063dSJacob Faibussowitsch for (d = 0; d < dim; ++d) PetscCall(PetscRandomGetValueReal(rnd, &xref[d])); 117735b38c60SMatthew G. Knepley CoordinatesRefToReal(dim, dim, xi0, v0, J, xref, &x[p * dim]); 117835b38c60SMatthew G. Knepley 1179ea7032a0SMatthew G. Knepley weight[p] = 1.0 / Np; 11806c5a79ebSMatthew G. Knepley species[p] = p % Ns; 118135b38c60SMatthew G. Knepley } 1182fe11354eSMatthew G. Knepley PetscCall(DMSwarmSortRestorePointsPerCell(sw, c, &Npc, &pidx)); 118335b38c60SMatthew G. Knepley } 11849566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rnd)); 11859566063dSJacob Faibussowitsch PetscCall(DMSwarmSortRestoreAccess(sw)); 11866c5a79ebSMatthew G. Knepley } 118719307e5cSMatthew G. Knepley PetscCall(DMSwarmRestoreField(sw, coordFields[0], NULL, NULL, (void **)&x)); 11886c5a79ebSMatthew G. Knepley PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&weight)); 11899566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **)&species)); 11906c5a79ebSMatthew G. Knepley 11919566063dSJacob Faibussowitsch PetscCall(DMSwarmMigrate(sw, removePoints)); 11929566063dSJacob Faibussowitsch PetscCall(DMLocalizeCoordinates(sw)); 11933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 119435b38c60SMatthew G. Knepley } 119535b38c60SMatthew G. Knepley 119635b38c60SMatthew G. Knepley /*@C 119735b38c60SMatthew G. Knepley DMSwarmInitializeVelocities - Set the initial velocities of particles using a distribution. 119835b38c60SMatthew G. Knepley 119920f4b53cSBarry Smith Collective 120035b38c60SMatthew G. Knepley 120135b38c60SMatthew G. Knepley Input Parameters: 120220f4b53cSBarry Smith + sw - The `DMSWARM` object 120335b38c60SMatthew G. Knepley . sampler - A function which uniformly samples the velocity PDF 120435b38c60SMatthew G. Knepley - v0 - The velocity scale for nondimensionalization for each species 120535b38c60SMatthew G. Knepley 120635b38c60SMatthew G. Knepley Level: advanced 120735b38c60SMatthew G. Knepley 120820f4b53cSBarry Smith Note: 120920f4b53cSBarry Smith If `v0` is zero for the first species, all velocities are set to zero. If it is zero for any other species, the effect will be to give that species zero velocity. 121020f4b53cSBarry Smith 121120f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmComputeLocalSize()`, `DMSwarmInitializeCoordinates()`, `DMSwarmInitializeVelocitiesFromOptions()` 121235b38c60SMatthew G. Knepley @*/ 1213f8662bd6SBarry Smith PetscErrorCode DMSwarmInitializeVelocities(DM sw, PetscProbFn *sampler, const PetscReal v0[]) 1214d71ae5a4SJacob Faibussowitsch { 12158434afd1SBarry Smith PetscSimplePointFn *velFunc; 121635b38c60SMatthew G. Knepley PetscReal *v; 121735b38c60SMatthew G. Knepley PetscInt *species; 12186c5a79ebSMatthew G. Knepley void *ctx; 121935b38c60SMatthew G. Knepley PetscInt dim, Np, p; 122035b38c60SMatthew G. Knepley 122135b38c60SMatthew G. Knepley PetscFunctionBegin; 12226c5a79ebSMatthew G. Knepley PetscCall(DMSwarmGetVelocityFunction(sw, &velFunc)); 122335b38c60SMatthew G. Knepley 12249566063dSJacob Faibussowitsch PetscCall(DMGetDimension(sw, &dim)); 12259566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(sw, &Np)); 12269566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(sw, "velocity", NULL, NULL, (void **)&v)); 12279566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(sw, "species", NULL, NULL, (void **)&species)); 12286c5a79ebSMatthew G. Knepley if (v0[0] == 0.) { 12296c5a79ebSMatthew G. Knepley PetscCall(PetscArrayzero(v, Np * dim)); 12306c5a79ebSMatthew G. Knepley } else if (velFunc) { 12316c5a79ebSMatthew G. Knepley PetscCall(DMGetApplicationContext(sw, &ctx)); 12326c5a79ebSMatthew G. Knepley for (p = 0; p < Np; ++p) { 12336c5a79ebSMatthew G. Knepley PetscInt s = species[p], d; 12346c5a79ebSMatthew G. Knepley PetscScalar vel[3]; 12356c5a79ebSMatthew G. Knepley 12366c5a79ebSMatthew G. Knepley PetscCall((*velFunc)(dim, 0., NULL, p, vel, ctx)); 12376c5a79ebSMatthew G. Knepley for (d = 0; d < dim; ++d) v[p * dim + d] = (v0[s] / v0[0]) * PetscRealPart(vel[d]); 12386c5a79ebSMatthew G. Knepley } 12396c5a79ebSMatthew G. Knepley } else { 12406c5a79ebSMatthew G. Knepley PetscRandom rnd; 12416c5a79ebSMatthew G. Knepley 12426c5a79ebSMatthew G. Knepley PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)sw), &rnd)); 12436c5a79ebSMatthew G. Knepley PetscCall(PetscRandomSetInterval(rnd, 0, 1.)); 12446c5a79ebSMatthew G. Knepley PetscCall(PetscRandomSetFromOptions(rnd)); 12456c5a79ebSMatthew G. Knepley 124635b38c60SMatthew G. Knepley for (p = 0; p < Np; ++p) { 124735b38c60SMatthew G. Knepley PetscInt s = species[p], d; 124835b38c60SMatthew G. Knepley PetscReal a[3], vel[3]; 124935b38c60SMatthew G. Knepley 12509566063dSJacob Faibussowitsch for (d = 0; d < dim; ++d) PetscCall(PetscRandomGetValueReal(rnd, &a[d])); 12519566063dSJacob Faibussowitsch PetscCall(sampler(a, NULL, vel)); 1252ad540459SPierre Jolivet for (d = 0; d < dim; ++d) v[p * dim + d] = (v0[s] / v0[0]) * vel[d]; 125335b38c60SMatthew G. Knepley } 12546c5a79ebSMatthew G. Knepley PetscCall(PetscRandomDestroy(&rnd)); 12556c5a79ebSMatthew G. Knepley } 12569566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sw, "velocity", NULL, NULL, (void **)&v)); 12579566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **)&species)); 12583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 125935b38c60SMatthew G. Knepley } 126035b38c60SMatthew G. Knepley 126135b38c60SMatthew G. Knepley /*@ 126235b38c60SMatthew G. Knepley DMSwarmInitializeVelocitiesFromOptions - Set the initial velocities of particles using a distribution determined from options. 126335b38c60SMatthew G. Knepley 126420f4b53cSBarry Smith Collective 126535b38c60SMatthew G. Knepley 126635b38c60SMatthew G. Knepley Input Parameters: 126720f4b53cSBarry Smith + sw - The `DMSWARM` object 126835b38c60SMatthew G. Knepley - v0 - The velocity scale for nondimensionalization for each species 126935b38c60SMatthew G. Knepley 127035b38c60SMatthew G. Knepley Level: advanced 127135b38c60SMatthew G. Knepley 127220f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmComputeLocalSize()`, `DMSwarmInitializeCoordinates()`, `DMSwarmInitializeVelocities()` 127335b38c60SMatthew G. Knepley @*/ 1274d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmInitializeVelocitiesFromOptions(DM sw, const PetscReal v0[]) 1275d71ae5a4SJacob Faibussowitsch { 1276f8662bd6SBarry Smith PetscProbFn *sampler; 127735b38c60SMatthew G. Knepley PetscInt dim; 127835b38c60SMatthew G. Knepley const char *prefix; 12796c5a79ebSMatthew G. Knepley char funcname[PETSC_MAX_PATH_LEN]; 12806c5a79ebSMatthew G. Knepley PetscBool flg; 128135b38c60SMatthew G. Knepley 128235b38c60SMatthew G. Knepley PetscFunctionBegin; 1283d0609cedSBarry Smith PetscOptionsBegin(PetscObjectComm((PetscObject)sw), "", "DMSwarm Options", "DMSWARM"); 12846c5a79ebSMatthew G. Knepley PetscCall(PetscOptionsString("-dm_swarm_velocity_function", "Function to determine particle velocities", "DMSwarmSetVelocityFunction", funcname, funcname, sizeof(funcname), &flg)); 1285d0609cedSBarry Smith PetscOptionsEnd(); 12866c5a79ebSMatthew G. Knepley if (flg) { 12878434afd1SBarry Smith PetscSimplePointFn *velFunc; 12886c5a79ebSMatthew G. Knepley 12896c5a79ebSMatthew G. Knepley PetscCall(PetscDLSym(NULL, funcname, (void **)&velFunc)); 12906c5a79ebSMatthew G. Knepley PetscCheck(velFunc, PetscObjectComm((PetscObject)sw), PETSC_ERR_ARG_WRONG, "Could not locate function %s", funcname); 12916c5a79ebSMatthew G. Knepley PetscCall(DMSwarmSetVelocityFunction(sw, velFunc)); 12926c5a79ebSMatthew G. Knepley } 12939566063dSJacob Faibussowitsch PetscCall(DMGetDimension(sw, &dim)); 12949566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)sw, &prefix)); 12959566063dSJacob Faibussowitsch PetscCall(PetscProbCreateFromOptions(dim, prefix, "-dm_swarm_velocity_density", NULL, NULL, &sampler)); 12969566063dSJacob Faibussowitsch PetscCall(DMSwarmInitializeVelocities(sw, sampler, v0)); 12973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 129835b38c60SMatthew G. Knepley } 12999d50c5ebSMatthew G. Knepley 13009d50c5ebSMatthew G. Knepley // The input vector U is assumed to be from a PetscFE. The Swarm fields are input as auxiliary values. 13019d50c5ebSMatthew G. Knepley PetscErrorCode DMProjectFieldLocal_Swarm(DM dm, PetscReal time, Vec U, PetscPointFn **funcs, InsertMode mode, Vec X) 13029d50c5ebSMatthew G. Knepley { 13039d50c5ebSMatthew G. Knepley MPI_Comm comm; 13049d50c5ebSMatthew G. Knepley DM dmIn; 13059d50c5ebSMatthew G. Knepley PetscDS ds; 13069d50c5ebSMatthew G. Knepley PetscTabulation *T; 13079d50c5ebSMatthew G. Knepley DMSwarmCellDM celldm; 13089d50c5ebSMatthew G. Knepley PetscScalar *a, *val, *u, *u_x; 13099d50c5ebSMatthew G. Knepley PetscFEGeom fegeom; 13109d50c5ebSMatthew G. Knepley PetscReal *xi, *v0, *J, *invJ, detJ = 1.0, v0ref[3] = {-1.0, -1.0, -1.0}; 13119d50c5ebSMatthew G. Knepley PetscInt dim, dE, Np, n, Nf, Nfc, Nfu, cStart, cEnd, maxC = 0, totbs = 0; 13129d50c5ebSMatthew G. Knepley const char **coordFields, **fields; 13139d50c5ebSMatthew G. Knepley PetscReal **coordVals, **vals; 13149d50c5ebSMatthew G. Knepley PetscInt *cbs, *bs, *uOff, *uOff_x; 13159d50c5ebSMatthew G. Knepley 13169d50c5ebSMatthew G. Knepley PetscFunctionBegin; 13179d50c5ebSMatthew G. Knepley PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 13189d50c5ebSMatthew G. Knepley PetscCall(VecGetDM(U, &dmIn)); 13199d50c5ebSMatthew G. Knepley PetscCall(DMGetDimension(dmIn, &dim)); 13209d50c5ebSMatthew G. Knepley PetscCall(DMGetCoordinateDim(dmIn, &dE)); 13219d50c5ebSMatthew G. Knepley PetscCall(DMGetDS(dmIn, &ds)); 13229d50c5ebSMatthew G. Knepley PetscCall(PetscDSGetNumFields(ds, &Nfu)); 13239d50c5ebSMatthew G. Knepley PetscCall(PetscDSGetComponentOffsets(ds, &uOff)); 13249d50c5ebSMatthew G. Knepley PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x)); 13259d50c5ebSMatthew G. Knepley PetscCall(PetscDSGetTabulation(ds, &T)); 13269d50c5ebSMatthew G. Knepley PetscCall(PetscDSGetEvaluationArrays(ds, &u, NULL, &u_x)); 13279d50c5ebSMatthew G. Knepley PetscCall(PetscMalloc3(dim, &v0, dim * dim, &J, dim * dim, &invJ)); 13289d50c5ebSMatthew G. Knepley PetscCall(DMPlexGetHeightStratum(dmIn, 0, &cStart, &cEnd)); 13299d50c5ebSMatthew G. Knepley 13309d50c5ebSMatthew G. Knepley fegeom.dim = dim; 13319d50c5ebSMatthew G. Knepley fegeom.dimEmbed = dE; 13329d50c5ebSMatthew G. Knepley fegeom.v = v0; 13339d50c5ebSMatthew G. Knepley fegeom.xi = v0ref; 13349d50c5ebSMatthew G. Knepley fegeom.J = J; 13359d50c5ebSMatthew G. Knepley fegeom.invJ = invJ; 13369d50c5ebSMatthew G. Knepley fegeom.detJ = &detJ; 13379d50c5ebSMatthew G. Knepley 13389d50c5ebSMatthew G. Knepley PetscCall(DMSwarmGetLocalSize(dm, &Np)); 13399d50c5ebSMatthew G. Knepley PetscCall(VecGetLocalSize(X, &n)); 13409d50c5ebSMatthew G. Knepley PetscCheck(n == Np, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Output vector local size %" PetscInt_FMT " != %" PetscInt_FMT " number of local particles", n, Np); 13419d50c5ebSMatthew G. Knepley PetscCall(DMSwarmGetCellDMActive(dm, &celldm)); 13429d50c5ebSMatthew G. Knepley PetscCall(DMSwarmCellDMGetCoordinateFields(celldm, &Nfc, &coordFields)); 13439d50c5ebSMatthew G. Knepley PetscCall(DMSwarmCellDMGetFields(celldm, &Nf, &fields)); 13449d50c5ebSMatthew G. Knepley 13459d50c5ebSMatthew G. Knepley PetscCall(PetscMalloc2(Nfc, &coordVals, Nfc, &cbs)); 13469d50c5ebSMatthew G. Knepley for (PetscInt i = 0; i < Nfc; ++i) PetscCall(DMSwarmGetField(dm, coordFields[i], &cbs[i], NULL, (void **)&coordVals[i])); 13479d50c5ebSMatthew G. Knepley PetscCall(PetscMalloc2(Nf, &vals, Nfc, &bs)); 13489d50c5ebSMatthew G. Knepley for (PetscInt i = 0; i < Nf; ++i) { 13499d50c5ebSMatthew G. Knepley PetscCall(DMSwarmGetField(dm, fields[i], &bs[i], NULL, (void **)&vals[i])); 13509d50c5ebSMatthew G. Knepley totbs += bs[i]; 13519d50c5ebSMatthew G. Knepley } 13529d50c5ebSMatthew G. Knepley 13539d50c5ebSMatthew G. Knepley PetscCall(DMSwarmSortGetAccess(dm)); 13549d50c5ebSMatthew G. Knepley for (PetscInt cell = cStart; cell < cEnd; ++cell) { 13559d50c5ebSMatthew G. Knepley PetscInt *pindices, Npc; 13569d50c5ebSMatthew G. Knepley 13579d50c5ebSMatthew G. Knepley PetscCall(DMSwarmSortGetPointsPerCell(dm, cell, &Npc, &pindices)); 13589d50c5ebSMatthew G. Knepley maxC = PetscMax(maxC, Npc); 13599d50c5ebSMatthew G. Knepley PetscCall(DMSwarmSortRestorePointsPerCell(dm, cell, &Npc, &pindices)); 13609d50c5ebSMatthew G. Knepley } 13619d50c5ebSMatthew G. Knepley PetscCall(PetscMalloc3(maxC * dim, &xi, maxC * totbs, &val, Nfu, &T)); 13629d50c5ebSMatthew G. Knepley PetscCall(VecGetArray(X, &a)); 13639d50c5ebSMatthew G. Knepley { 13649d50c5ebSMatthew G. Knepley for (PetscInt cell = cStart; cell < cEnd; ++cell) { 13659d50c5ebSMatthew G. Knepley PetscInt *pindices, Npc; 13669d50c5ebSMatthew G. Knepley 13679d50c5ebSMatthew G. Knepley // TODO: Use DMField instead of assuming affine 13689d50c5ebSMatthew G. Knepley PetscCall(DMPlexComputeCellGeometryFEM(dmIn, cell, NULL, v0, J, invJ, &detJ)); 13699d50c5ebSMatthew G. Knepley PetscCall(DMSwarmSortGetPointsPerCell(dm, cell, &Npc, &pindices)); 13709d50c5ebSMatthew G. Knepley 13719d50c5ebSMatthew G. Knepley PetscScalar *closure = NULL; 13729d50c5ebSMatthew G. Knepley PetscInt Ncl; 13739d50c5ebSMatthew G. Knepley 13749d50c5ebSMatthew G. Knepley // Get fields from input vector and auxiliary fields from swarm 13759d50c5ebSMatthew G. Knepley for (PetscInt p = 0; p < Npc; ++p) { 13769d50c5ebSMatthew G. Knepley PetscReal xr[8]; 13779d50c5ebSMatthew G. Knepley PetscInt off; 13789d50c5ebSMatthew G. Knepley 13799d50c5ebSMatthew G. Knepley off = 0; 13809d50c5ebSMatthew G. Knepley for (PetscInt i = 0; i < Nfc; ++i) { 13819d50c5ebSMatthew G. Knepley for (PetscInt b = 0; b < cbs[i]; ++b, ++off) xr[off] = coordVals[i][pindices[p] * cbs[i] + b]; 13829d50c5ebSMatthew G. Knepley } 13839d50c5ebSMatthew G. Knepley PetscCheck(off == dim, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The total block size of coordinates is %" PetscInt_FMT " != %" PetscInt_FMT " the DM coordinate dimension", off, dim); 13849d50c5ebSMatthew G. Knepley CoordinatesRealToRef(dE, dim, fegeom.xi, fegeom.v, fegeom.invJ, xr, &xi[p * dim]); 13859d50c5ebSMatthew G. Knepley off = 0; 13869d50c5ebSMatthew G. Knepley for (PetscInt i = 0; i < Nf; ++i) { 13879d50c5ebSMatthew G. Knepley for (PetscInt b = 0; b < bs[i]; ++b, ++off) val[p * totbs + off] = vals[i][pindices[p] * bs[i] + b]; 13889d50c5ebSMatthew G. Knepley } 13899d50c5ebSMatthew G. Knepley PetscCheck(off == totbs, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The total block size of swarm fields is %" PetscInt_FMT " != %" PetscInt_FMT " the computed total block size", off, totbs); 13909d50c5ebSMatthew G. Knepley } 13919d50c5ebSMatthew G. Knepley PetscCall(DMPlexVecGetClosure(dmIn, NULL, U, cell, &Ncl, &closure)); 13929d50c5ebSMatthew G. Knepley for (PetscInt field = 0; field < Nfu; ++field) { 13939d50c5ebSMatthew G. Knepley PetscFE fe; 13949d50c5ebSMatthew G. Knepley 13959d50c5ebSMatthew G. Knepley PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe)); 13969d50c5ebSMatthew G. Knepley PetscCall(PetscFECreateTabulation(fe, 1, Npc, xi, 1, &T[field])); 13979d50c5ebSMatthew G. Knepley } 13989d50c5ebSMatthew G. Knepley for (PetscInt p = 0; p < Npc; ++p) { 13999d50c5ebSMatthew G. Knepley // Get fields from input vector 14009d50c5ebSMatthew G. Knepley PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nfu, 0, p, T, &fegeom, closure, NULL, u, u_x, NULL)); 14019d50c5ebSMatthew G. Knepley (*funcs[0])(dim, 1, 1, uOff, uOff_x, u, NULL, u_x, bs, NULL, &val[p * totbs], NULL, NULL, time, &xi[p * dim], 0, NULL, &a[pindices[p]]); 14029d50c5ebSMatthew G. Knepley } 14039d50c5ebSMatthew G. Knepley PetscCall(DMSwarmSortRestorePointsPerCell(dm, cell, &Npc, &pindices)); 14049d50c5ebSMatthew G. Knepley PetscCall(DMPlexVecRestoreClosure(dmIn, NULL, U, cell, &Ncl, &closure)); 14059d50c5ebSMatthew G. Knepley for (PetscInt field = 0; field < Nfu; ++field) PetscCall(PetscTabulationDestroy(&T[field])); 14069d50c5ebSMatthew G. Knepley } 14079d50c5ebSMatthew G. Knepley } 14089d50c5ebSMatthew G. Knepley for (PetscInt i = 0; i < Nfc; ++i) PetscCall(DMSwarmRestoreField(dm, coordFields[i], &cbs[i], NULL, (void **)&coordVals[i])); 14099d50c5ebSMatthew G. Knepley for (PetscInt i = 0; i < Nf; ++i) PetscCall(DMSwarmRestoreField(dm, fields[i], &bs[i], NULL, (void **)&vals[i])); 14109d50c5ebSMatthew G. Knepley PetscCall(VecRestoreArray(X, &a)); 14119d50c5ebSMatthew G. Knepley PetscCall(DMSwarmSortRestoreAccess(dm)); 14129d50c5ebSMatthew G. Knepley PetscCall(PetscFree3(xi, val, T)); 14139d50c5ebSMatthew G. Knepley PetscCall(PetscFree3(v0, J, invJ)); 14149d50c5ebSMatthew G. Knepley PetscCall(PetscFree2(coordVals, cbs)); 14159d50c5ebSMatthew G. Knepley PetscCall(PetscFree2(vals, bs)); 14169d50c5ebSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 14179d50c5ebSMatthew G. Knepley } 1418