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: /*
23: Error checking to ensure the swarm type is correct and that a cell DM has been set
24: */
25: #define DMSWARMPICVALID(dm) \
26: do { \
27: DM_Swarm *_swarm = (DM_Swarm *)(dm)->data; \
28: PetscCheck(_swarm->swarm_type == DMSWARM_PIC, PetscObjectComm((PetscObject)(dm)), PETSC_ERR_SUP, "Valid only for DMSwarm-PIC. You must call DMSwarmSetType(dm,DMSWARM_PIC)"); \
29: PetscCheck(_swarm->dmcell, PetscObjectComm((PetscObject)(dm)), PETSC_ERR_SUP, "Valid only for DMSwarmPIC if the cell DM is set. You must call DMSwarmSetCellDM(dm,celldm)"); \
30: } while (0)
32: typedef struct {
33: DMSwarmDataBucket db;
34: PetscInt refct;
35: PetscBool field_registration_initialized;
36: PetscBool field_registration_finalized;
37: /* DMSwarmProjectMethod *swarm_project;*/ /* swarm, geometry, result */
39: /* PetscInt overlap; */
40: /* PetscErrorCode (*update_overlap)(void); */
42: char vec_field_name[PETSC_MAX_PATH_LEN];
43: PetscBool vec_field_set;
44: PetscInt vec_field_bs, vec_field_nlocal;
46: PetscBool issetup;
47: DMSwarmType swarm_type;
48: DMSwarmMigrateType migrate_type;
49: DMSwarmCollectType collect_type;
51: DM dmcell;
53: PetscBool migrate_error_on_missing_point;
55: PetscBool collect_view_active;
56: PetscInt collect_view_reset_nlocal;
57: DMSwarmSort sort_context;
59: /* Support for PIC */
60: PetscInt Ns; /* The number of particle species */
62: PetscSimplePointFn *coordFunc; /* Function to set particle coordinates */
63: PetscSimplePointFn *velFunc; /* Function to set particle velocities */
64: } DM_Swarm;
66: typedef struct {
67: PetscInt point_index;
68: PetscInt cell_index;
69: } SwarmPoint;
71: struct _p_DMSwarmSort {
72: PetscBool isvalid;
73: PetscInt ncells, npoints;
74: PetscInt *pcell_offsets;
75: SwarmPoint *list;
76: };
78: PETSC_INTERN PetscErrorCode DMSwarmMigrate_Push_Basic(DM, PetscBool);
79: PETSC_INTERN PetscErrorCode DMSwarmMigrate_CellDMScatter(DM, PetscBool);
80: PETSC_INTERN PetscErrorCode DMSwarmMigrate_CellDMExact(DM, PetscBool);