157795646SDave May 257795646SDave May #define PETSCDM_DLL 357795646SDave May #include <petsc/private/dmswarmimpl.h> /*I "petscdmswarm.h" I*/ 4b62e03f8SDave May #include "data_bucket.h" 557795646SDave May 657795646SDave May //typedef PetscErrorCode (*swarm_project)(DM,DM,Vec) DMSwarmProjectMethod; /* swarm, geometry, result */ 757795646SDave May 857795646SDave May //typedef enum { PROJECT_DMDA_AQ1=0, PROJECT_DMDA_P0 } DMSwarmDMDAProjectionType; 957795646SDave May 1057795646SDave May #if 0 1157795646SDave May 1257795646SDave May /* Defines what the local space will be */ 13b62e03f8SDave May PetscErrorCode DMSwarmSetOverlap(void) 1457795646SDave May { 1557795646SDave May 16b62e03f8SDave May PetscFunctionReturn(0); 1757795646SDave May } 1857795646SDave May 1957795646SDave May 2057795646SDave May /* coordinates */ 2157795646SDave May /* 2257795646SDave May DMGetCoordinateDM returns self 2357795646SDave May DMGetCoordinates and DMGetCoordinatesLocal return same thing 2457795646SDave May Local view could be used to define overlapping information 2557795646SDave May */ 2657795646SDave May 2757795646SDave May #endif 2857795646SDave May 2957795646SDave May #undef __FUNCT__ 30*b5bcf523SDave May #define __FUNCT__ "DMSwarmVectorDefineField" 31*b5bcf523SDave May PETSC_EXTERN PetscErrorCode DMSwarmVectorDefineField(DM dm,const char fieldname[]) 32*b5bcf523SDave May { 33*b5bcf523SDave May DM_Swarm *swarm = (DM_Swarm*)dm->data; 34*b5bcf523SDave May PetscErrorCode ierr; 35*b5bcf523SDave May PetscInt bs,n; 36*b5bcf523SDave May PetscScalar *array; 37*b5bcf523SDave May PetscDataType type; 38*b5bcf523SDave May 39*b5bcf523SDave May DataBucketGetSizes(swarm->db,&n,NULL,NULL); 40*b5bcf523SDave May ierr = DMSwarmGetField(dm,fieldname,&bs,&type,(void**)&array);CHKERRQ(ierr); 41*b5bcf523SDave May 42*b5bcf523SDave May /* Check all fields are of type PETSC_REAL or PETSC_SCALAR */ 43*b5bcf523SDave May if (type != PETSC_REAL) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only valid for PETSC_REAL"); 44*b5bcf523SDave May 45*b5bcf523SDave May PetscSNPrintf(swarm->vec_field_name,PETSC_MAX_PATH_LEN-1,"%s",fieldname); 46*b5bcf523SDave May swarm->vec_field_set = PETSC_TRUE; 47*b5bcf523SDave May swarm->vec_field_bs = 1;//bs; 48*b5bcf523SDave May swarm->vec_field_nlocal = n; 49*b5bcf523SDave May 50*b5bcf523SDave May PetscFunctionReturn(0); 51*b5bcf523SDave May } 52*b5bcf523SDave May 53*b5bcf523SDave May /* requires DMSwarmDefineFieldVector has been called */ 54*b5bcf523SDave May #undef __FUNCT__ 55*b5bcf523SDave May #define __FUNCT__ "DMCreateGlobalVector_Swarm" 56*b5bcf523SDave May PetscErrorCode DMCreateGlobalVector_Swarm(DM dm,Vec *vec) 57*b5bcf523SDave May { 58*b5bcf523SDave May DM_Swarm *swarm = (DM_Swarm*)dm->data; 59*b5bcf523SDave May PetscErrorCode ierr; 60*b5bcf523SDave May Vec x; 61*b5bcf523SDave May char name[PETSC_MAX_PATH_LEN]; 62*b5bcf523SDave May 63*b5bcf523SDave May if (!swarm->vec_field_set) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMSwarmVectorDefineField first"); 64*b5bcf523SDave May PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"DMSwarmField_%s",swarm->vec_field_name); 65*b5bcf523SDave May ierr = VecCreate(PetscObjectComm((PetscObject)dm),&x);CHKERRQ(ierr); 66*b5bcf523SDave May ierr = PetscObjectSetName((PetscObject)x,name);CHKERRQ(ierr); 67*b5bcf523SDave May ierr = VecSetSizes(x,swarm->vec_field_nlocal,PETSC_DETERMINE);CHKERRQ(ierr); 68*b5bcf523SDave May ierr = VecSetBlockSize(x,swarm->vec_field_bs);CHKERRQ(ierr); 69*b5bcf523SDave May ierr = VecSetFromOptions(x);CHKERRQ(ierr); 70*b5bcf523SDave May *vec = x; 71*b5bcf523SDave May 72*b5bcf523SDave May PetscFunctionReturn(0); 73*b5bcf523SDave May } 74*b5bcf523SDave May 75*b5bcf523SDave May /* requires DMSwarmDefineFieldVector has been called */ 76*b5bcf523SDave May #undef __FUNCT__ 77*b5bcf523SDave May #define __FUNCT__ "DMCreateLocalVector_Swarm" 78*b5bcf523SDave May PetscErrorCode DMCreateLocalVector_Swarm(DM dm,Vec *vec) 79*b5bcf523SDave May { 80*b5bcf523SDave May DM_Swarm *swarm = (DM_Swarm*)dm->data; 81*b5bcf523SDave May PetscErrorCode ierr; 82*b5bcf523SDave May Vec x; 83*b5bcf523SDave May 84*b5bcf523SDave May char name[PETSC_MAX_PATH_LEN]; 85*b5bcf523SDave May 86*b5bcf523SDave May if (!swarm->vec_field_set) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMSwarmVectorDefineField first"); 87*b5bcf523SDave May PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"DMSwarmField_%s",swarm->vec_field_name); 88*b5bcf523SDave May ierr = VecCreate(PETSC_COMM_SELF,&x);CHKERRQ(ierr); 89*b5bcf523SDave May ierr = PetscObjectSetName((PetscObject)x,name);CHKERRQ(ierr); 90*b5bcf523SDave May ierr = VecSetSizes(x,swarm->vec_field_nlocal,swarm->vec_field_nlocal);CHKERRQ(ierr); 91*b5bcf523SDave May ierr = VecSetBlockSize(x,swarm->vec_field_bs);CHKERRQ(ierr); 92*b5bcf523SDave May ierr = VecSetFromOptions(x);CHKERRQ(ierr); 93*b5bcf523SDave May *vec = x; 94*b5bcf523SDave May 95*b5bcf523SDave May PetscFunctionReturn(0); 96*b5bcf523SDave May } 97*b5bcf523SDave May 98*b5bcf523SDave May #undef __FUNCT__ 99*b5bcf523SDave May #define __FUNCT__ "DMSwarmCreateGlobalVectorFromField" 100*b5bcf523SDave May PETSC_EXTERN PetscErrorCode DMSwarmCreateGlobalVectorFromField(DM dm,const char fieldname[],Vec *vec) 101*b5bcf523SDave May { 102*b5bcf523SDave May DM_Swarm *swarm = (DM_Swarm*)dm->data; 103*b5bcf523SDave May PetscErrorCode ierr; 104*b5bcf523SDave May PetscInt bs,n; 105*b5bcf523SDave May PetscScalar *array; 106*b5bcf523SDave May Vec x; 107*b5bcf523SDave May PetscDataType type; 108*b5bcf523SDave May char name[PETSC_MAX_PATH_LEN]; 109*b5bcf523SDave May 110*b5bcf523SDave May DataBucketGetSizes(swarm->db,&n,NULL,NULL); 111*b5bcf523SDave May ierr = DMSwarmGetField(dm,fieldname,&bs,&type,(void**)&array);CHKERRQ(ierr); 112*b5bcf523SDave May 113*b5bcf523SDave May /* Check all fields are of type PETSC_REAL or PETSC_SCALAR */ 114*b5bcf523SDave May if (type != PETSC_REAL) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only valid for PETSC_REAL"); 115*b5bcf523SDave May 116*b5bcf523SDave May ierr = VecCreateSeqWithArray(PetscObjectComm((PetscObject)dm),1,n,array,&x);CHKERRQ(ierr); 117*b5bcf523SDave May 118*b5bcf523SDave May PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"DMSwarm_VecFieldInPlace_%s",fieldname); 119*b5bcf523SDave May ierr = PetscObjectComposeFunction((PetscObject)x,name,DMSwarmDestroyGlobalVectorFromField);CHKERRQ(ierr); 120*b5bcf523SDave May 121*b5bcf523SDave May /* Set guard */ 122*b5bcf523SDave May 123*b5bcf523SDave May *vec = x; 124*b5bcf523SDave May PetscFunctionReturn(0); 125*b5bcf523SDave May } 126*b5bcf523SDave May 127*b5bcf523SDave May #undef __FUNCT__ 128*b5bcf523SDave May #define __FUNCT__ "DMSwarmDestroyGlobalVectorFromField" 129*b5bcf523SDave May PETSC_EXTERN PetscErrorCode DMSwarmDestroyGlobalVectorFromField(DM dm,const char fieldname[],Vec *vec) 130*b5bcf523SDave May { 131*b5bcf523SDave May DM_Swarm *swarm = (DM_Swarm*)dm->data; 132*b5bcf523SDave May PetscErrorCode ierr; 133*b5bcf523SDave May DataField gfield; 134*b5bcf523SDave May char name[PETSC_MAX_PATH_LEN]; 135*b5bcf523SDave May void (*fptr)(void); 136*b5bcf523SDave May 137*b5bcf523SDave May /* get data field */ 138*b5bcf523SDave May DataBucketGetDataFieldByName(swarm->db,fieldname,&gfield); 139*b5bcf523SDave May 140*b5bcf523SDave May /* check vector is an inplace array */ 141*b5bcf523SDave May PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"DMSwarm_VecFieldInPlace_%s",fieldname); 142*b5bcf523SDave May ierr = PetscObjectQueryFunction((PetscObject)(*vec),name,&fptr);CHKERRQ(ierr); 143*b5bcf523SDave May if (!fptr) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Vector being destroyed was not created from DMSwarm field(%s)",fieldname); 144*b5bcf523SDave May 145*b5bcf523SDave May /* restore data field */ 146*b5bcf523SDave May DataFieldRestoreAccess(gfield); 147*b5bcf523SDave May 148*b5bcf523SDave May ierr = VecDestroy(vec);CHKERRQ(ierr); 149*b5bcf523SDave May 150*b5bcf523SDave May PetscFunctionReturn(0); 151*b5bcf523SDave May } 152*b5bcf523SDave May 153*b5bcf523SDave May /* 154*b5bcf523SDave May PETSC_EXTERN PetscErrorCode DMSwarmCreateGlobalVectorFromFields(DM dm,const PetscInt nf,const char *fieldnames[],Vec *vec) 155*b5bcf523SDave May { 156*b5bcf523SDave May PetscFunctionReturn(0); 157*b5bcf523SDave May } 158*b5bcf523SDave May 159*b5bcf523SDave May PETSC_EXTERN PetscErrorCode DMSwarmRestoreGlobalVectorFromFields(DM dm,Vec *vec) 160*b5bcf523SDave May { 161*b5bcf523SDave May PetscFunctionReturn(0); 162*b5bcf523SDave May } 163*b5bcf523SDave May */ 164*b5bcf523SDave May 165*b5bcf523SDave May #undef __FUNCT__ 1665f50eb2eSDave May #define __FUNCT__ "DMSwarmInitializeFieldRegister" 1675f50eb2eSDave May PETSC_EXTERN PetscErrorCode DMSwarmInitializeFieldRegister(DM dm) 1685f50eb2eSDave May { 1695f50eb2eSDave May DM_Swarm *swarm = (DM_Swarm*)dm->data; 1705f50eb2eSDave May swarm->field_registration_initialized = PETSC_TRUE; 1715f50eb2eSDave May PetscFunctionReturn(0); 1725f50eb2eSDave May } 1735f50eb2eSDave May 1745f50eb2eSDave May #undef __FUNCT__ 1755f50eb2eSDave May #define __FUNCT__ "DMSwarmFinalizeFieldRegister" 1765f50eb2eSDave May PETSC_EXTERN PetscErrorCode DMSwarmFinalizeFieldRegister(DM dm) 1775f50eb2eSDave May { 1785f50eb2eSDave May DM_Swarm *swarm = (DM_Swarm*)dm->data; 1795f50eb2eSDave May swarm->field_registration_finalized = PETSC_TRUE; 1805f50eb2eSDave May DataBucketFinalize(swarm->db); 1815f50eb2eSDave May PetscFunctionReturn(0); 1825f50eb2eSDave May } 1835f50eb2eSDave May 1845f50eb2eSDave May #undef __FUNCT__ 1855f50eb2eSDave May #define __FUNCT__ "DMSwarmSetLocalSizes" 1865f50eb2eSDave May PETSC_EXTERN PetscErrorCode DMSwarmSetLocalSizes(DM dm,PetscInt nlocal,PetscInt buffer) 1875f50eb2eSDave May { 1885f50eb2eSDave May DM_Swarm *swarm = (DM_Swarm*)dm->data; 1895f50eb2eSDave May 1905f50eb2eSDave May DataBucketSetSizes(swarm->db,nlocal,buffer); 1915f50eb2eSDave May 1925f50eb2eSDave May PetscFunctionReturn(0); 1935f50eb2eSDave May } 1945f50eb2eSDave May 1955f50eb2eSDave May #undef __FUNCT__ 196b62e03f8SDave May #define __FUNCT__ "DMSwarmRegisterPetscDatatypeField" 1975f50eb2eSDave May PETSC_EXTERN PetscErrorCode DMSwarmRegisterPetscDatatypeField(DM dm,const char fieldname[],PetscInt blocksize,PetscDataType type) 198b62e03f8SDave May { 199b62e03f8SDave May DM_Swarm *swarm = (DM_Swarm*)dm->data; 200b62e03f8SDave May size_t size; 201b62e03f8SDave May 2025f50eb2eSDave May if (!swarm->field_registration_initialized) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMSwarmInitializeFieldRegister() first"); 2035f50eb2eSDave May if (swarm->field_registration_finalized) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Cannot register additional fields after calling DMSwarmFinalizeFieldRegister() first"); 2045f50eb2eSDave May 2055f50eb2eSDave May if (type == PETSC_OBJECT) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}"); 2065f50eb2eSDave May if (type == PETSC_FUNCTION) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}"); 2075f50eb2eSDave May if (type == PETSC_STRING) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}"); 2085f50eb2eSDave May if (type == PETSC_STRUCT) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}"); 2095f50eb2eSDave May if (type == PETSC_DATATYPE_UNKNOWN) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}"); 210b62e03f8SDave May 211b62e03f8SDave May switch (type) { 212b62e03f8SDave May case PETSC_CHAR: 213b62e03f8SDave May size = sizeof(PetscChar); 214b62e03f8SDave May break; 215b62e03f8SDave May case PETSC_SHORT: 216b62e03f8SDave May size = sizeof(PetscShort); 217b62e03f8SDave May break; 218b62e03f8SDave May case PETSC_INT: 219b62e03f8SDave May size = sizeof(PetscInt); 220b62e03f8SDave May break; 221b62e03f8SDave May case PETSC_LONG: 222b62e03f8SDave May size = sizeof(Petsc64bitInt); 223b62e03f8SDave May break; 224b62e03f8SDave May case PETSC_FLOAT: 225b62e03f8SDave May size = sizeof(PetscFloat); 226b62e03f8SDave May break; 227b62e03f8SDave May case PETSC_DOUBLE: 228b62e03f8SDave May size = sizeof(PetscReal); 229b62e03f8SDave May break; 230b62e03f8SDave May 231b62e03f8SDave May default: 2325f50eb2eSDave May SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}"); 233b62e03f8SDave May break; 234b62e03f8SDave May } 235b62e03f8SDave May 236b62e03f8SDave May /* Load a specific data type into data bucket, specifying textual name and its size in bytes */ 237b62e03f8SDave May DataBucketRegisterField(swarm->db,fieldname,size,NULL); 238b62e03f8SDave May swarm->db->field[swarm->db->nfields-1]->petsc_type = type; 239b62e03f8SDave May 240b62e03f8SDave May PetscFunctionReturn(0); 241b62e03f8SDave May } 242b62e03f8SDave May 243b62e03f8SDave May #undef __FUNCT__ 244b62e03f8SDave May #define __FUNCT__ "DMSwarmRegisterUserStructField" 2455f50eb2eSDave May PETSC_EXTERN PetscErrorCode DMSwarmRegisterUserStructField(DM dm,const char fieldname[],size_t size) 246b62e03f8SDave May { 247b62e03f8SDave May DM_Swarm *swarm = (DM_Swarm*)dm->data; 248b62e03f8SDave May 249b62e03f8SDave May DataBucketRegisterField(swarm->db,fieldname,size,NULL); 250b62e03f8SDave May swarm->db->field[swarm->db->nfields-1]->petsc_type = PETSC_STRUCT ; 251b62e03f8SDave May 252b62e03f8SDave May PetscFunctionReturn(0); 253b62e03f8SDave May } 254b62e03f8SDave May 255b62e03f8SDave May #undef __FUNCT__ 256b62e03f8SDave May #define __FUNCT__ "DMSwarmRegisterUserDatatypeField" 2575f50eb2eSDave May PETSC_EXTERN PetscErrorCode DMSwarmRegisterUserDatatypeField(DM dm,const char fieldname[],size_t size) 258b62e03f8SDave May { 259b62e03f8SDave May DM_Swarm *swarm = (DM_Swarm*)dm->data; 260b62e03f8SDave May 261b62e03f8SDave May DataBucketRegisterField(swarm->db,fieldname,size,NULL); 262b62e03f8SDave May swarm->db->field[swarm->db->nfields-1]->petsc_type = PETSC_DATATYPE_UNKNOWN; 263b62e03f8SDave May 264b62e03f8SDave May PetscFunctionReturn(0); 265b62e03f8SDave May } 266b62e03f8SDave May 267b62e03f8SDave May #undef __FUNCT__ 268b62e03f8SDave May #define __FUNCT__ "DMSwarmGetField" 2695f50eb2eSDave May PETSC_EXTERN PetscErrorCode DMSwarmGetField(DM dm,const char fieldname[],PetscInt *blocksize,PetscDataType *type,void **data) 270b62e03f8SDave May { 271b62e03f8SDave May DM_Swarm *swarm = (DM_Swarm*)dm->data; 272b62e03f8SDave May DataField gfield; 273b62e03f8SDave May 274b62e03f8SDave May DataBucketGetDataFieldByName(swarm->db,fieldname,&gfield); 275b62e03f8SDave May DataFieldGetAccess(gfield); 276b62e03f8SDave May DataFieldGetEntries(gfield,data); 277*b5bcf523SDave May if (blocksize) {} 278*b5bcf523SDave May if (type) { *type = gfield->petsc_type; } 279b62e03f8SDave May 280b62e03f8SDave May PetscFunctionReturn(0); 281b62e03f8SDave May } 282b62e03f8SDave May 283b62e03f8SDave May #undef __FUNCT__ 284b62e03f8SDave May #define __FUNCT__ "DMSwarmRestoreField" 2855f50eb2eSDave May PETSC_EXTERN PetscErrorCode DMSwarmRestoreField(DM dm,const char fieldname[],PetscInt *blocksize,PetscDataType *type,void **data) 286b62e03f8SDave May { 287b62e03f8SDave May DM_Swarm *swarm = (DM_Swarm*)dm->data; 288b62e03f8SDave May DataField gfield; 289b62e03f8SDave May 290b62e03f8SDave May DataBucketGetDataFieldByName(swarm->db,fieldname,&gfield); 291b62e03f8SDave May DataFieldRestoreAccess(gfield); 292b62e03f8SDave May if (data) *data = NULL; 293b62e03f8SDave May 294b62e03f8SDave May PetscFunctionReturn(0); 295b62e03f8SDave May } 296b62e03f8SDave May 297b62e03f8SDave May 298b62e03f8SDave May #undef __FUNCT__ 29957795646SDave May #define __FUNCT__ "DMDestroy_Swarm" 30057795646SDave May PetscErrorCode DMDestroy_Swarm(DM dm) 30157795646SDave May { 30257795646SDave May DM_Swarm *swarm = (DM_Swarm*)dm->data; 30357795646SDave May PetscErrorCode ierr; 30457795646SDave May 30557795646SDave May PetscFunctionBegin; 306b62e03f8SDave May DataBucketDestroy(&swarm->db); 30757795646SDave May ierr = PetscFree(swarm);CHKERRQ(ierr); 30857795646SDave May PetscFunctionReturn(0); 30957795646SDave May } 31057795646SDave May 31157795646SDave May #undef __FUNCT__ 3125f50eb2eSDave May #define __FUNCT__ "DMView_Swarm" 3135f50eb2eSDave May PetscErrorCode DMView_Swarm(DM dm, PetscViewer viewer) 3145f50eb2eSDave May { 3155f50eb2eSDave May DM_Swarm *swarm = (DM_Swarm*)dm->data; 3165f50eb2eSDave May PetscBool iascii,ibinary,ishdf5,isvtk; 3175f50eb2eSDave May PetscErrorCode ierr; 3185f50eb2eSDave May 3195f50eb2eSDave May PetscFunctionBegin; 3205f50eb2eSDave May PetscValidHeaderSpecific(dm,DM_CLASSID,1); 3215f50eb2eSDave May PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 3225f50eb2eSDave May ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 3235f50eb2eSDave May ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY,&ibinary);CHKERRQ(ierr); 3245f50eb2eSDave May ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERVTK, &isvtk);CHKERRQ(ierr); 3255f50eb2eSDave May ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5);CHKERRQ(ierr); 3265f50eb2eSDave May if (iascii) { 3275f50eb2eSDave May DataBucketView(PetscObjectComm((PetscObject)dm),swarm->db,NULL,DATABUCKET_VIEW_STDOUT); 3285f50eb2eSDave May } else if (ibinary) { 3295f50eb2eSDave May SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"NO VTK support"); 3305f50eb2eSDave May } else if (ishdf5) { 3315f50eb2eSDave May #if defined(PETSC_HAVE_HDF5) 3325f50eb2eSDave May SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"NO HDF5 support"); 3335f50eb2eSDave May #else 3345f50eb2eSDave May SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"HDF5 not supported. Please reconfigure using --download-hdf5"); 3355f50eb2eSDave May #endif 3365f50eb2eSDave May } else if (isvtk) { 3375f50eb2eSDave May SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"NO VTK support"); 3385f50eb2eSDave May } 3395f50eb2eSDave May PetscFunctionReturn(0); 3405f50eb2eSDave May } 3415f50eb2eSDave May 3425f50eb2eSDave May #undef __FUNCT__ 34357795646SDave May #define __FUNCT__ "DMCreate_Swarm" 34457795646SDave May PETSC_EXTERN PetscErrorCode DMCreate_Swarm(DM dm) 34557795646SDave May { 34657795646SDave May DM_Swarm *swarm; 34757795646SDave May PetscErrorCode ierr; 34857795646SDave May 34957795646SDave May PetscFunctionBegin; 35057795646SDave May PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 35157795646SDave May ierr = PetscNewLog(dm,&swarm);CHKERRQ(ierr); 35257795646SDave May 353b62e03f8SDave May DataBucketCreate(&swarm->db); 354*b5bcf523SDave May swarm->vec_field_set = PETSC_FALSE; 355b62e03f8SDave May 35657795646SDave May dm->dim = 0; 35757795646SDave May dm->data = swarm; 35857795646SDave May 3595f50eb2eSDave May dm->ops->view = DMView_Swarm; 36057795646SDave May dm->ops->load = NULL; 36157795646SDave May dm->ops->setfromoptions = NULL; 36257795646SDave May dm->ops->clone = NULL; 36357795646SDave May dm->ops->setup = NULL; 36457795646SDave May dm->ops->createdefaultsection = NULL; 36557795646SDave May dm->ops->createdefaultconstraints = NULL; 366*b5bcf523SDave May dm->ops->createglobalvector = DMCreateGlobalVector_Swarm; 367*b5bcf523SDave May dm->ops->createlocalvector = DMCreateLocalVector_Swarm; 36857795646SDave May dm->ops->getlocaltoglobalmapping = NULL; 36957795646SDave May dm->ops->createfieldis = NULL; 37057795646SDave May dm->ops->createcoordinatedm = NULL; 37157795646SDave May dm->ops->getcoloring = NULL; 37257795646SDave May dm->ops->creatematrix = NULL; 37357795646SDave May dm->ops->createinterpolation = NULL; 37457795646SDave May dm->ops->getaggregates = NULL; 37557795646SDave May dm->ops->getinjection = NULL; 37657795646SDave May dm->ops->refine = NULL; 37757795646SDave May dm->ops->coarsen = NULL; 37857795646SDave May dm->ops->refinehierarchy = NULL; 37957795646SDave May dm->ops->coarsenhierarchy = NULL; 38057795646SDave May dm->ops->globaltolocalbegin = NULL; 38157795646SDave May dm->ops->globaltolocalend = NULL; 38257795646SDave May dm->ops->localtoglobalbegin = NULL; 38357795646SDave May dm->ops->localtoglobalend = NULL; 38457795646SDave May dm->ops->destroy = DMDestroy_Swarm; 38557795646SDave May dm->ops->createsubdm = NULL; 38657795646SDave May dm->ops->getdimpoints = NULL; 38757795646SDave May dm->ops->locatepoints = NULL; 38857795646SDave May 38957795646SDave May PetscFunctionReturn(0); 39057795646SDave May }