1 2 #define PETSCDM_DLL 3 #include <petsc/private/dmswarmimpl.h> /*I "petscdmswarm.h" I*/ 4 #include "data_bucket.h" 5 6 PetscErrorCode DMSwarmMigrate_Push_Basic(DM dm,PetscBool remove_sent_points); 7 8 9 //typedef PetscErrorCode (*swarm_project)(DM,DM,Vec) DMSwarmProjectMethod; /* swarm, geometry, result */ 10 11 //typedef enum { PROJECT_DMDA_AQ1=0, PROJECT_DMDA_P0 } DMSwarmDMDAProjectionType; 12 13 #if 0 14 15 /* Defines what the local space will be */ 16 PetscErrorCode DMSwarmSetOverlap(void) 17 { 18 19 PetscFunctionReturn(0); 20 } 21 22 23 /* coordinates */ 24 /* 25 DMGetCoordinateDM returns self 26 DMGetCoordinates and DMGetCoordinatesLocal return same thing 27 Local view could be used to define overlapping information 28 */ 29 30 #endif 31 32 #undef __FUNCT__ 33 #define __FUNCT__ "DMSwarmVectorDefineField" 34 PETSC_EXTERN PetscErrorCode DMSwarmVectorDefineField(DM dm,const char fieldname[]) 35 { 36 DM_Swarm *swarm = (DM_Swarm*)dm->data; 37 PetscErrorCode ierr; 38 PetscInt bs,n; 39 PetscScalar *array; 40 PetscDataType type; 41 42 if (!swarm->issetup) { ierr = DMSetUp(dm);CHKERRQ(ierr); } 43 ierr = DataBucketGetSizes(swarm->db,&n,NULL,NULL);CHKERRQ(ierr); 44 ierr = DMSwarmGetField(dm,fieldname,&bs,&type,(void**)&array);CHKERRQ(ierr); 45 46 /* Check all fields are of type PETSC_REAL or PETSC_SCALAR */ 47 if (type != PETSC_REAL) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only valid for PETSC_REAL"); 48 49 PetscSNPrintf(swarm->vec_field_name,PETSC_MAX_PATH_LEN-1,"%s",fieldname); 50 swarm->vec_field_set = PETSC_TRUE; 51 swarm->vec_field_bs = 1;//bs; 52 swarm->vec_field_nlocal = n; 53 ierr = DMSwarmRestoreField(dm,fieldname,&bs,&type,(void**)&array);CHKERRQ(ierr); 54 55 PetscFunctionReturn(0); 56 } 57 58 /* requires DMSwarmDefineFieldVector has been called */ 59 #undef __FUNCT__ 60 #define __FUNCT__ "DMCreateGlobalVector_Swarm" 61 PetscErrorCode DMCreateGlobalVector_Swarm(DM dm,Vec *vec) 62 { 63 DM_Swarm *swarm = (DM_Swarm*)dm->data; 64 PetscErrorCode ierr; 65 Vec x; 66 char name[PETSC_MAX_PATH_LEN]; 67 68 if (!swarm->issetup) { ierr = DMSetUp(dm);CHKERRQ(ierr); } 69 if (!swarm->vec_field_set) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMSwarmVectorDefineField first"); 70 PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"DMSwarmField_%s",swarm->vec_field_name); 71 ierr = VecCreate(PetscObjectComm((PetscObject)dm),&x);CHKERRQ(ierr); 72 ierr = PetscObjectSetName((PetscObject)x,name);CHKERRQ(ierr); 73 ierr = VecSetSizes(x,swarm->db->L,PETSC_DETERMINE);CHKERRQ(ierr); 74 ierr = VecSetBlockSize(x,swarm->vec_field_bs);CHKERRQ(ierr); 75 ierr = VecSetFromOptions(x);CHKERRQ(ierr); 76 *vec = x; 77 78 PetscFunctionReturn(0); 79 } 80 81 /* requires DMSwarmDefineFieldVector has been called */ 82 #undef __FUNCT__ 83 #define __FUNCT__ "DMCreateLocalVector_Swarm" 84 PetscErrorCode DMCreateLocalVector_Swarm(DM dm,Vec *vec) 85 { 86 DM_Swarm *swarm = (DM_Swarm*)dm->data; 87 PetscErrorCode ierr; 88 Vec x; 89 char name[PETSC_MAX_PATH_LEN]; 90 91 if (!swarm->issetup) { ierr = DMSetUp(dm);CHKERRQ(ierr); } 92 if (!swarm->vec_field_set) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMSwarmVectorDefineField first"); 93 PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"DMSwarmField_%s",swarm->vec_field_name); 94 ierr = VecCreate(PETSC_COMM_SELF,&x);CHKERRQ(ierr); 95 ierr = PetscObjectSetName((PetscObject)x,name);CHKERRQ(ierr); 96 ierr = VecSetSizes(x,swarm->db->L,swarm->db->L);CHKERRQ(ierr); 97 ierr = VecSetBlockSize(x,swarm->vec_field_bs);CHKERRQ(ierr); 98 ierr = VecSetFromOptions(x);CHKERRQ(ierr); 99 *vec = x; 100 101 PetscFunctionReturn(0); 102 } 103 104 #undef __FUNCT__ 105 #define __FUNCT__ "DMSwarmCreateGlobalVectorFromField" 106 PETSC_EXTERN PetscErrorCode DMSwarmCreateGlobalVectorFromField(DM dm,const char fieldname[],Vec *vec) 107 { 108 DM_Swarm *swarm = (DM_Swarm*)dm->data; 109 PetscErrorCode ierr; 110 PetscInt bs,n; 111 PetscScalar *array; 112 Vec x; 113 PetscDataType type; 114 char name[PETSC_MAX_PATH_LEN]; 115 PetscMPIInt commsize; 116 117 if (!swarm->issetup) { ierr = DMSetUp(dm);CHKERRQ(ierr); } 118 ierr = DataBucketGetSizes(swarm->db,&n,NULL,NULL);CHKERRQ(ierr); 119 ierr = DMSwarmGetField(dm,fieldname,&bs,&type,(void**)&array);CHKERRQ(ierr); 120 121 /* Check all fields are of type PETSC_REAL or PETSC_SCALAR */ 122 if (type != PETSC_REAL) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only valid for PETSC_REAL"); 123 124 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&commsize);CHKERRQ(ierr); 125 if (commsize == 1) { 126 ierr = VecCreateSeqWithArray(PetscObjectComm((PetscObject)dm),1,n,array,&x);CHKERRQ(ierr); 127 } else { 128 ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)dm),1,n,PETSC_DETERMINE,array,&x);CHKERRQ(ierr); 129 } 130 PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"DMSwarmSharedField_%s",swarm->vec_field_name); 131 ierr = PetscObjectSetName((PetscObject)x,name);CHKERRQ(ierr); 132 133 /* Set guard */ 134 PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"DMSwarm_VecFieldInPlace_%s",fieldname); 135 ierr = PetscObjectComposeFunction((PetscObject)x,name,DMSwarmDestroyGlobalVectorFromField);CHKERRQ(ierr); 136 137 *vec = x; 138 PetscFunctionReturn(0); 139 } 140 141 #undef __FUNCT__ 142 #define __FUNCT__ "DMSwarmDestroyGlobalVectorFromField" 143 PETSC_EXTERN PetscErrorCode DMSwarmDestroyGlobalVectorFromField(DM dm,const char fieldname[],Vec *vec) 144 { 145 DM_Swarm *swarm = (DM_Swarm*)dm->data; 146 PetscErrorCode ierr; 147 DataField gfield; 148 char name[PETSC_MAX_PATH_LEN]; 149 void (*fptr)(void); 150 151 /* get data field */ 152 ierr = DataBucketGetDataFieldByName(swarm->db,fieldname,&gfield);CHKERRQ(ierr); 153 154 /* check vector is an inplace array */ 155 PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"DMSwarm_VecFieldInPlace_%s",fieldname); 156 ierr = PetscObjectQueryFunction((PetscObject)(*vec),name,&fptr);CHKERRQ(ierr); 157 if (!fptr) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Vector being destroyed was not created from DMSwarm field(%s)",fieldname); 158 159 /* restore data field */ 160 ierr = DataFieldRestoreAccess(gfield);CHKERRQ(ierr); 161 162 ierr = VecDestroy(vec);CHKERRQ(ierr); 163 164 PetscFunctionReturn(0); 165 } 166 167 /* 168 PETSC_EXTERN PetscErrorCode DMSwarmCreateGlobalVectorFromFields(DM dm,const PetscInt nf,const char *fieldnames[],Vec *vec) 169 { 170 PetscFunctionReturn(0); 171 } 172 173 PETSC_EXTERN PetscErrorCode DMSwarmRestoreGlobalVectorFromFields(DM dm,Vec *vec) 174 { 175 PetscFunctionReturn(0); 176 } 177 */ 178 179 #undef __FUNCT__ 180 #define __FUNCT__ "DMSwarmInitializeFieldRegister" 181 PETSC_EXTERN PetscErrorCode DMSwarmInitializeFieldRegister(DM dm) 182 { 183 DM_Swarm *swarm = (DM_Swarm*)dm->data; 184 PetscErrorCode ierr; 185 186 swarm->field_registration_initialized = PETSC_TRUE; 187 188 ierr = DMSwarmRegisterPetscDatatypeField(dm,"DMSwarm_pid",1,PETSC_LONG);CHKERRQ(ierr); /* unique identifer */ 189 ierr = DMSwarmRegisterPetscDatatypeField(dm,"DMSwarm_rank",1,PETSC_INT);CHKERRQ(ierr); /* used for communication */ 190 191 PetscFunctionReturn(0); 192 } 193 194 #undef __FUNCT__ 195 #define __FUNCT__ "DMSwarmFinalizeFieldRegister" 196 PETSC_EXTERN PetscErrorCode DMSwarmFinalizeFieldRegister(DM dm) 197 { 198 DM_Swarm *swarm = (DM_Swarm*)dm->data; 199 PetscErrorCode ierr; 200 201 swarm->field_registration_finalized = PETSC_TRUE; 202 ierr = DataBucketFinalize(swarm->db);CHKERRQ(ierr); 203 PetscFunctionReturn(0); 204 } 205 206 #undef __FUNCT__ 207 #define __FUNCT__ "DMSwarmSetLocalSizes" 208 PETSC_EXTERN PetscErrorCode DMSwarmSetLocalSizes(DM dm,PetscInt nlocal,PetscInt buffer) 209 { 210 DM_Swarm *swarm = (DM_Swarm*)dm->data; 211 PetscErrorCode ierr; 212 213 ierr = DataBucketSetSizes(swarm->db,nlocal,buffer);CHKERRQ(ierr); 214 215 PetscFunctionReturn(0); 216 } 217 218 #undef __FUNCT__ 219 #define __FUNCT__ "DMSwarmSetCellDM" 220 PETSC_EXTERN PetscErrorCode DMSwarmSetCellDM(DM dm,DM dmcell) 221 { 222 DM_Swarm *swarm = (DM_Swarm*)dm->data; 223 swarm->dmcell = dmcell; 224 PetscFunctionReturn(0); 225 } 226 227 #undef __FUNCT__ 228 #define __FUNCT__ "DMSwarmGetLocalSize" 229 PETSC_EXTERN PetscErrorCode DMSwarmGetLocalSize(DM dm,PetscInt *nlocal) 230 { 231 DM_Swarm *swarm = (DM_Swarm*)dm->data; 232 PetscErrorCode ierr; 233 234 if (nlocal) { 235 ierr = DataBucketGetSizes(swarm->db,nlocal,NULL,NULL);CHKERRQ(ierr); 236 } 237 238 PetscFunctionReturn(0); 239 } 240 241 #undef __FUNCT__ 242 #define __FUNCT__ "DMSwarmGetSize" 243 PETSC_EXTERN PetscErrorCode DMSwarmGetSize(DM dm,PetscInt *n) 244 { 245 DM_Swarm *swarm = (DM_Swarm*)dm->data; 246 PetscErrorCode ierr; 247 PetscInt nlocal,ng; 248 249 ierr = DataBucketGetSizes(swarm->db,&nlocal,NULL,NULL);CHKERRQ(ierr); 250 ierr = MPI_Allreduce(&nlocal,&ng,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 251 if (n) { *n = ng; } 252 PetscFunctionReturn(0); 253 } 254 255 #undef __FUNCT__ 256 #define __FUNCT__ "DMSwarmRegisterPetscDatatypeField" 257 PETSC_EXTERN PetscErrorCode DMSwarmRegisterPetscDatatypeField(DM dm,const char fieldname[],PetscInt blocksize,PetscDataType type) 258 { 259 PetscErrorCode ierr; 260 DM_Swarm *swarm = (DM_Swarm*)dm->data; 261 size_t size; 262 263 if (!swarm->field_registration_initialized) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMSwarmInitializeFieldRegister() first"); 264 if (swarm->field_registration_finalized) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Cannot register additional fields after calling DMSwarmFinalizeFieldRegister() first"); 265 266 if (type == PETSC_OBJECT) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}"); 267 if (type == PETSC_FUNCTION) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}"); 268 if (type == PETSC_STRING) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}"); 269 if (type == PETSC_STRUCT) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}"); 270 if (type == PETSC_DATATYPE_UNKNOWN) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}"); 271 272 switch (type) { 273 case PETSC_CHAR: 274 size = sizeof(PetscChar); 275 break; 276 case PETSC_SHORT: 277 size = sizeof(PetscShort); 278 break; 279 case PETSC_INT: 280 size = sizeof(PetscInt); 281 break; 282 case PETSC_LONG: 283 size = sizeof(Petsc64bitInt); 284 break; 285 case PETSC_FLOAT: 286 size = sizeof(PetscFloat); 287 break; 288 case PETSC_DOUBLE: 289 size = sizeof(PetscReal); 290 break; 291 292 default: 293 SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}"); 294 break; 295 } 296 297 /* Load a specific data type into data bucket, specifying textual name and its size in bytes */ 298 ierr = DataBucketRegisterField(swarm->db,"DMSwarmRegisterPetscDatatypeField",fieldname,size,NULL);CHKERRQ(ierr); 299 swarm->db->field[swarm->db->nfields-1]->petsc_type = type; 300 301 PetscFunctionReturn(0); 302 } 303 304 #undef __FUNCT__ 305 #define __FUNCT__ "DMSwarmRegisterUserStructField" 306 PETSC_EXTERN PetscErrorCode DMSwarmRegisterUserStructField(DM dm,const char fieldname[],size_t size) 307 { 308 PetscErrorCode ierr; 309 DM_Swarm *swarm = (DM_Swarm*)dm->data; 310 311 ierr = DataBucketRegisterField(swarm->db,"DMSwarmRegisterUserStructField",fieldname,size,NULL);CHKERRQ(ierr); 312 swarm->db->field[swarm->db->nfields-1]->petsc_type = PETSC_STRUCT ; 313 314 PetscFunctionReturn(0); 315 } 316 317 #undef __FUNCT__ 318 #define __FUNCT__ "DMSwarmRegisterUserDatatypeField" 319 PETSC_EXTERN PetscErrorCode DMSwarmRegisterUserDatatypeField(DM dm,const char fieldname[],size_t size) 320 { 321 DM_Swarm *swarm = (DM_Swarm*)dm->data; 322 PetscErrorCode ierr; 323 324 ierr = DataBucketRegisterField(swarm->db,"DMSwarmRegisterUserDatatypeField",fieldname,size,NULL);CHKERRQ(ierr); 325 swarm->db->field[swarm->db->nfields-1]->petsc_type = PETSC_DATATYPE_UNKNOWN; 326 327 PetscFunctionReturn(0); 328 } 329 330 #undef __FUNCT__ 331 #define __FUNCT__ "DMSwarmGetField" 332 PETSC_EXTERN PetscErrorCode DMSwarmGetField(DM dm,const char fieldname[],PetscInt *blocksize,PetscDataType *type,void **data) 333 { 334 DM_Swarm *swarm = (DM_Swarm*)dm->data; 335 DataField gfield; 336 PetscErrorCode ierr; 337 338 if (!swarm->issetup) { ierr = DMSetUp(dm);CHKERRQ(ierr); } 339 340 ierr = DataBucketGetDataFieldByName(swarm->db,fieldname,&gfield);CHKERRQ(ierr); 341 ierr = DataFieldGetAccess(gfield);CHKERRQ(ierr); 342 ierr = DataFieldGetEntries(gfield,data);CHKERRQ(ierr); 343 if (blocksize) {} 344 if (type) { *type = gfield->petsc_type; } 345 346 PetscFunctionReturn(0); 347 } 348 349 #undef __FUNCT__ 350 #define __FUNCT__ "DMSwarmRestoreField" 351 PETSC_EXTERN PetscErrorCode DMSwarmRestoreField(DM dm,const char fieldname[],PetscInt *blocksize,PetscDataType *type,void **data) 352 { 353 DM_Swarm *swarm = (DM_Swarm*)dm->data; 354 DataField gfield; 355 PetscErrorCode ierr; 356 357 ierr = DataBucketGetDataFieldByName(swarm->db,fieldname,&gfield);CHKERRQ(ierr); 358 ierr = DataFieldRestoreAccess(gfield);CHKERRQ(ierr); 359 if (data) *data = NULL; 360 361 PetscFunctionReturn(0); 362 } 363 364 #undef __FUNCT__ 365 #define __FUNCT__ "DMSwarmAddPoint" 366 PETSC_EXTERN PetscErrorCode DMSwarmAddPoint(DM dm) 367 { 368 DM_Swarm *swarm = (DM_Swarm*)dm->data; 369 PetscErrorCode ierr; 370 371 if (!swarm->issetup) { ierr = DMSetUp(dm);CHKERRQ(ierr); } 372 ierr = DataBucketAddPoint(swarm->db);CHKERRQ(ierr); 373 PetscFunctionReturn(0); 374 } 375 376 #undef __FUNCT__ 377 #define __FUNCT__ "DMSwarmAddNPoints" 378 PETSC_EXTERN PetscErrorCode DMSwarmAddNPoints(DM dm,PetscInt npoints) 379 { 380 DM_Swarm *swarm = (DM_Swarm*)dm->data; 381 PetscErrorCode ierr; 382 PetscInt nlocal; 383 384 ierr = DataBucketGetSizes(swarm->db,&nlocal,NULL,NULL);CHKERRQ(ierr); 385 nlocal = nlocal + npoints; 386 ierr = DataBucketSetSizes(swarm->db,nlocal,-1);CHKERRQ(ierr); 387 PetscFunctionReturn(0); 388 } 389 390 #undef __FUNCT__ 391 #define __FUNCT__ "DMSwarmRemovePoint" 392 PETSC_EXTERN PetscErrorCode DMSwarmRemovePoint(DM dm) 393 { 394 DM_Swarm *swarm = (DM_Swarm*)dm->data; 395 PetscErrorCode ierr; 396 397 ierr = DataBucketRemovePoint(swarm->db);CHKERRQ(ierr); 398 PetscFunctionReturn(0); 399 } 400 401 #undef __FUNCT__ 402 #define __FUNCT__ "DMSwarmRemovePointAtIndex" 403 PETSC_EXTERN PetscErrorCode DMSwarmRemovePointAtIndex(DM dm,PetscInt idx) 404 { 405 DM_Swarm *swarm = (DM_Swarm*)dm->data; 406 PetscErrorCode ierr; 407 408 ierr = DataBucketRemovePointAtIndex(swarm->db,idx);CHKERRQ(ierr); 409 PetscFunctionReturn(0); 410 } 411 412 #undef __FUNCT__ 413 #define __FUNCT__ "DMSwarmMigrate_Basic" 414 PetscErrorCode DMSwarmMigrate_Basic(DM dm,PetscBool remove_sent_points) 415 { 416 PetscErrorCode ierr; 417 ierr = DMSwarmMigrate_Push_Basic(dm,remove_sent_points);CHKERRQ(ierr); 418 PetscFunctionReturn(0); 419 } 420 421 #undef __FUNCT__ 422 #define __FUNCT__ "DMSwarmMigrate" 423 PETSC_EXTERN PetscErrorCode DMSwarmMigrate(DM dm,PetscBool remove_sent_points) 424 { 425 PetscErrorCode ierr; 426 ierr = DMSwarmMigrate_Basic(dm,remove_sent_points);CHKERRQ(ierr); 427 PetscFunctionReturn(0); 428 } 429 430 #undef __FUNCT__ 431 #define __FUNCT__ "DMSwarmGlobalToLocalViewCreate" 432 PETSC_EXTERN PetscErrorCode DMSwarmGlobalToLocalViewCreate(DM dm,InsertMode mode) 433 { 434 PetscErrorCode ierr; 435 DM_Swarm *swarm = (DM_Swarm*)dm->data; 436 PetscInt ng; 437 438 if (mode != INSERT_VALUES) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only mode INSERT_VALUES is supported"); 439 440 ierr = DMSwarmMigrate_GlobalToLocal_Basic(dm,&ng);CHKERRQ(ierr); 441 swarm->view_ng = ng; 442 443 PetscFunctionReturn(0); 444 } 445 446 #undef __FUNCT__ 447 #define __FUNCT__ "DMSwarmGlobalToLocalViewDestroy" 448 PETSC_EXTERN PetscErrorCode DMSwarmGlobalToLocalViewDestroy(DM dm,InsertMode mode) 449 { 450 PetscErrorCode ierr; 451 DM_Swarm *swarm = (DM_Swarm*)dm->data; 452 453 if (mode != INSERT_VALUES) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only mode INSERT_VALUES is supported"); 454 455 ierr = DMSwarmSetLocalSizes(dm,swarm->view_ng,-1);CHKERRQ(ierr); 456 457 PetscFunctionReturn(0); 458 } 459 460 #undef __FUNCT__ 461 #define __FUNCT__ "DMSetup_Swarm" 462 PetscErrorCode DMSetup_Swarm(DM dm) 463 { 464 DM_Swarm *swarm = (DM_Swarm*)dm->data; 465 PetscErrorCode ierr; 466 PetscMPIInt rank; 467 PetscInt p,npoints,*rankval; 468 469 if (swarm->issetup) PetscFunctionReturn(0); 470 471 PetscPrintf(PETSC_COMM_SELF,"Swarm setup \n"); 472 swarm->issetup = PETSC_TRUE; 473 474 /* check some fields were registered */ 475 if (swarm->db->nfields <= 2) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"At least one field user must be registered via DMSwarmRegisterXXX()"); 476 477 /* check local sizes were set */ 478 if (swarm->db->L == -1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Local sizes must be set via DMSwarmSetLocalSizes()"); 479 480 /* initialize values in pid and rank placeholders */ 481 /* TODO: [pid - use MPI_Scan] */ 482 483 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr); 484 ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 485 ierr = DMSwarmGetField(dm,"DMSwarm_rank",NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 486 for (p=0; p<npoints; p++) { 487 rankval[p] = (PetscInt)rank; 488 } 489 ierr = DMSwarmRestoreField(dm,"DMSwarm_rank",NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 490 491 PetscFunctionReturn(0); 492 } 493 494 #undef __FUNCT__ 495 #define __FUNCT__ "DMDestroy_Swarm" 496 PetscErrorCode DMDestroy_Swarm(DM dm) 497 { 498 DM_Swarm *swarm = (DM_Swarm*)dm->data; 499 PetscErrorCode ierr; 500 501 PetscFunctionBegin; 502 ierr = DataBucketDestroy(&swarm->db);CHKERRQ(ierr); 503 ierr = PetscFree(swarm);CHKERRQ(ierr); 504 PetscFunctionReturn(0); 505 } 506 507 #undef __FUNCT__ 508 #define __FUNCT__ "DMView_Swarm" 509 PetscErrorCode DMView_Swarm(DM dm, PetscViewer viewer) 510 { 511 DM_Swarm *swarm = (DM_Swarm*)dm->data; 512 PetscBool iascii,ibinary,ishdf5,isvtk; 513 PetscErrorCode ierr; 514 515 PetscFunctionBegin; 516 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 517 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 518 ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 519 ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY,&ibinary);CHKERRQ(ierr); 520 ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERVTK, &isvtk);CHKERRQ(ierr); 521 ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5);CHKERRQ(ierr); 522 if (iascii) { 523 ierr = DataBucketView(PetscObjectComm((PetscObject)dm),swarm->db,NULL,DATABUCKET_VIEW_STDOUT);CHKERRQ(ierr); 524 } else if (ibinary) { 525 SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"NO VTK support"); 526 } else if (ishdf5) { 527 #if defined(PETSC_HAVE_HDF5) 528 SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"NO HDF5 support"); 529 #else 530 SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"HDF5 not supported. Please reconfigure using --download-hdf5"); 531 #endif 532 } else if (isvtk) { 533 SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"NO VTK support"); 534 } 535 PetscFunctionReturn(0); 536 } 537 538 #undef __FUNCT__ 539 #define __FUNCT__ "DMCreate_Swarm" 540 PETSC_EXTERN PetscErrorCode DMCreate_Swarm(DM dm) 541 { 542 DM_Swarm *swarm; 543 PetscErrorCode ierr; 544 545 PetscFunctionBegin; 546 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 547 ierr = PetscNewLog(dm,&swarm);CHKERRQ(ierr); 548 549 ierr = DataBucketCreate(&swarm->db);CHKERRQ(ierr); 550 swarm->vec_field_set = PETSC_FALSE; 551 swarm->issetup = PETSC_FALSE; 552 553 dm->dim = 0; 554 dm->data = swarm; 555 556 dm->ops->view = DMView_Swarm; 557 dm->ops->load = NULL; 558 dm->ops->setfromoptions = NULL; 559 dm->ops->clone = NULL; 560 dm->ops->setup = DMSetup_Swarm; 561 dm->ops->createdefaultsection = NULL; 562 dm->ops->createdefaultconstraints = NULL; 563 dm->ops->createglobalvector = DMCreateGlobalVector_Swarm; 564 dm->ops->createlocalvector = DMCreateLocalVector_Swarm; 565 dm->ops->getlocaltoglobalmapping = NULL; 566 dm->ops->createfieldis = NULL; 567 dm->ops->createcoordinatedm = NULL; 568 dm->ops->getcoloring = NULL; 569 dm->ops->creatematrix = NULL; 570 dm->ops->createinterpolation = NULL; 571 dm->ops->getaggregates = NULL; 572 dm->ops->getinjection = NULL; 573 dm->ops->refine = NULL; 574 dm->ops->coarsen = NULL; 575 dm->ops->refinehierarchy = NULL; 576 dm->ops->coarsenhierarchy = NULL; 577 dm->ops->globaltolocalbegin = NULL; 578 dm->ops->globaltolocalend = NULL; 579 dm->ops->localtoglobalbegin = NULL; 580 dm->ops->localtoglobalend = NULL; 581 dm->ops->destroy = DMDestroy_Swarm; 582 dm->ops->createsubdm = NULL; 583 dm->ops->getdimpoints = NULL; 584 dm->ops->locatepoints = NULL; 585 586 PetscFunctionReturn(0); 587 }