Actual source code: dmswarmimpl.h
1: #pragma once
3: #include <petscvec.h>
4: #include <petscmat.h>
5: #include <petscdmswarm.h>
6: #include <petsc/private/dmimpl.h>
8: PETSC_EXTERN PetscBool SwarmProjcite;
9: PETSC_EXTERN const char SwarmProjCitation[];
11: PETSC_EXTERN PetscLogEvent DMSWARM_Migrate;
12: PETSC_EXTERN PetscLogEvent DMSWARM_SetSizes;
13: PETSC_EXTERN PetscLogEvent DMSWARM_AddPoints;
14: PETSC_EXTERN PetscLogEvent DMSWARM_RemovePoints;
15: PETSC_EXTERN PetscLogEvent DMSWARM_Sort;
16: PETSC_EXTERN PetscLogEvent DMSWARM_DataExchangerTopologySetup;
17: PETSC_EXTERN PetscLogEvent DMSWARM_DataExchangerBegin;
18: PETSC_EXTERN PetscLogEvent DMSWARM_DataExchangerEnd;
19: PETSC_EXTERN PetscLogEvent DMSWARM_DataExchangerSendCount;
20: PETSC_EXTERN PetscLogEvent DMSWARM_DataExchangerPack;
22: typedef struct _DMSwarmCellDMOps *DMSwarmCellDMOps;
23: struct _DMSwarmCellDMOps {
24: PetscErrorCode (*create)(DMSwarmCellDM);
25: PetscErrorCode (*destroy)(DMSwarmCellDM);
26: PetscErrorCode (*setfromoptions)(PetscOptionItems *, DMSwarmCellDM);
27: PetscErrorCode (*setup)(DMSwarmCellDM);
28: PetscErrorCode (*view)(DMSwarmCellDM, PetscViewer);
29: };
31: struct _p_DMSwarmCellDM {
32: PETSCHEADER(struct _DMSwarmCellDMOps);
33: DM dm; // The cell DM
34: PetscInt Nf; // The number of DM fields
35: char **dmFields; // Swarm fields defining this DM
36: PetscInt Nfc; // The number of DM coordinate fields
37: char **coordFields; // Swarm field for coordinates on this DM
38: char *cellid; // Swarm field with cell ids for particles
39: DMSwarmSort sort; // Ordering of particles by enclosing cell
40: };
42: /*
43: Error checking to ensure the swarm type is correct and that a cell DM has been set
44: */
45: #define DMSWARMPICVALID(obj) \
46: do { \
47: DM_Swarm *_swarm = (DM_Swarm *)(obj)->data; \
48: PetscCheck(_swarm->swarm_type == DMSWARM_PIC, PetscObjectComm((PetscObject)(obj)), PETSC_ERR_SUP, "Valid only for DMSwarm-PIC. You must call DMSwarmSetType(dm,DMSWARM_PIC)"); \
49: PetscCheck(_swarm->activeCellDM, PetscObjectComm((PetscObject)(obj)), PETSC_ERR_SUP, "Valid only for DMSwarmPIC if the cell DM is set. You must call DMSwarmSetCellDM() or DMSwarmAddCellDM()"); \
50: } while (0)
52: typedef struct {
53: DMSwarmDataBucket db;
54: PetscInt refct;
55: PetscBool field_registration_initialized;
56: PetscBool field_registration_finalized;
57: /* DMSwarmProjectMethod *swarm_project;*/ /* swarm, geometry, result */
59: /* PetscInt overlap; */
60: /* PetscErrorCode (*update_overlap)(void); */
62: PetscObjectList cellDMs;
63: const char *activeCellDM;
65: PetscBool issetup;
66: DMSwarmType swarm_type;
67: DMSwarmMigrateType migrate_type;
68: DMSwarmCollectType collect_type;
69: DMSwarmRemapType remap_type;
71: PetscBool migrate_error_on_missing_point;
73: PetscBool collect_view_active;
74: PetscInt collect_view_reset_nlocal;
76: /* Support for PIC */
77: PetscInt Ns; /* The number of particle species */
79: PetscSimplePointFn *coordFunc; /* Function to set particle coordinates */
80: PetscSimplePointFn *velFunc; /* Function to set particle velocities */
82: /* Debugging */
83: PetscInt printCoords;
84: PetscInt printWeights;
85: } DM_Swarm;
87: typedef struct {
88: PetscInt point_index;
89: PetscInt cell_index;
90: } SwarmPoint;
92: struct _p_DMSwarmSort {
93: PetscBool isvalid;
94: PetscInt ncells, npoints;
95: PetscInt *pcell_offsets;
96: SwarmPoint *list;
97: };
99: PETSC_INTERN PetscErrorCode DMSwarmMigrate_Push_Basic(DM, PetscBool);
100: PETSC_INTERN PetscErrorCode DMSwarmMigrate_CellDMScatter(DM, PetscBool);
101: PETSC_INTERN PetscErrorCode DMSwarmMigrate_CellDMExact(DM, PetscBool);