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__ "DMSwarmGetCellDM" 229 PETSC_EXTERN PetscErrorCode DMSwarmGetCellDM(DM dm,DM *dmcell) 230 { 231 DM_Swarm *swarm = (DM_Swarm*)dm->data; 232 *dmcell = swarm->dmcell; 233 PetscFunctionReturn(0); 234 } 235 236 #undef __FUNCT__ 237 #define __FUNCT__ "DMSwarmGetLocalSize" 238 PETSC_EXTERN PetscErrorCode DMSwarmGetLocalSize(DM dm,PetscInt *nlocal) 239 { 240 DM_Swarm *swarm = (DM_Swarm*)dm->data; 241 PetscErrorCode ierr; 242 243 if (nlocal) { 244 ierr = DataBucketGetSizes(swarm->db,nlocal,NULL,NULL);CHKERRQ(ierr); 245 } 246 247 PetscFunctionReturn(0); 248 } 249 250 #undef __FUNCT__ 251 #define __FUNCT__ "DMSwarmGetSize" 252 PETSC_EXTERN PetscErrorCode DMSwarmGetSize(DM dm,PetscInt *n) 253 { 254 DM_Swarm *swarm = (DM_Swarm*)dm->data; 255 PetscErrorCode ierr; 256 PetscInt nlocal,ng; 257 258 ierr = DataBucketGetSizes(swarm->db,&nlocal,NULL,NULL);CHKERRQ(ierr); 259 ierr = MPI_Allreduce(&nlocal,&ng,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 260 if (n) { *n = ng; } 261 PetscFunctionReturn(0); 262 } 263 264 #undef __FUNCT__ 265 #define __FUNCT__ "DMSwarmRegisterPetscDatatypeField" 266 PETSC_EXTERN PetscErrorCode DMSwarmRegisterPetscDatatypeField(DM dm,const char fieldname[],PetscInt blocksize,PetscDataType type) 267 { 268 PetscErrorCode ierr; 269 DM_Swarm *swarm = (DM_Swarm*)dm->data; 270 size_t size; 271 272 if (!swarm->field_registration_initialized) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMSwarmInitializeFieldRegister() first"); 273 if (swarm->field_registration_finalized) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Cannot register additional fields after calling DMSwarmFinalizeFieldRegister() first"); 274 275 if (type == PETSC_OBJECT) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}"); 276 if (type == PETSC_FUNCTION) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}"); 277 if (type == PETSC_STRING) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}"); 278 if (type == PETSC_STRUCT) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}"); 279 if (type == PETSC_DATATYPE_UNKNOWN) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}"); 280 281 switch (type) { 282 case PETSC_CHAR: 283 size = sizeof(PetscChar); 284 break; 285 case PETSC_SHORT: 286 size = sizeof(PetscShort); 287 break; 288 case PETSC_INT: 289 size = sizeof(PetscInt); 290 break; 291 case PETSC_LONG: 292 size = sizeof(Petsc64bitInt); 293 break; 294 case PETSC_FLOAT: 295 size = sizeof(PetscFloat); 296 break; 297 case PETSC_DOUBLE: 298 size = sizeof(PetscReal); 299 break; 300 301 default: 302 SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}"); 303 break; 304 } 305 306 /* Load a specific data type into data bucket, specifying textual name and its size in bytes */ 307 ierr = DataBucketRegisterField(swarm->db,"DMSwarmRegisterPetscDatatypeField",fieldname,size,NULL);CHKERRQ(ierr); 308 swarm->db->field[swarm->db->nfields-1]->petsc_type = type; 309 310 PetscFunctionReturn(0); 311 } 312 313 #undef __FUNCT__ 314 #define __FUNCT__ "DMSwarmRegisterUserStructField" 315 PETSC_EXTERN PetscErrorCode DMSwarmRegisterUserStructField(DM dm,const char fieldname[],size_t size) 316 { 317 PetscErrorCode ierr; 318 DM_Swarm *swarm = (DM_Swarm*)dm->data; 319 320 ierr = DataBucketRegisterField(swarm->db,"DMSwarmRegisterUserStructField",fieldname,size,NULL);CHKERRQ(ierr); 321 swarm->db->field[swarm->db->nfields-1]->petsc_type = PETSC_STRUCT ; 322 323 PetscFunctionReturn(0); 324 } 325 326 #undef __FUNCT__ 327 #define __FUNCT__ "DMSwarmRegisterUserDatatypeField" 328 PETSC_EXTERN PetscErrorCode DMSwarmRegisterUserDatatypeField(DM dm,const char fieldname[],size_t size) 329 { 330 DM_Swarm *swarm = (DM_Swarm*)dm->data; 331 PetscErrorCode ierr; 332 333 ierr = DataBucketRegisterField(swarm->db,"DMSwarmRegisterUserDatatypeField",fieldname,size,NULL);CHKERRQ(ierr); 334 swarm->db->field[swarm->db->nfields-1]->petsc_type = PETSC_DATATYPE_UNKNOWN; 335 336 PetscFunctionReturn(0); 337 } 338 339 #undef __FUNCT__ 340 #define __FUNCT__ "DMSwarmGetField" 341 PETSC_EXTERN PetscErrorCode DMSwarmGetField(DM dm,const char fieldname[],PetscInt *blocksize,PetscDataType *type,void **data) 342 { 343 DM_Swarm *swarm = (DM_Swarm*)dm->data; 344 DataField gfield; 345 PetscErrorCode ierr; 346 347 if (!swarm->issetup) { ierr = DMSetUp(dm);CHKERRQ(ierr); } 348 349 ierr = DataBucketGetDataFieldByName(swarm->db,fieldname,&gfield);CHKERRQ(ierr); 350 ierr = DataFieldGetAccess(gfield);CHKERRQ(ierr); 351 ierr = DataFieldGetEntries(gfield,data);CHKERRQ(ierr); 352 if (blocksize) {} 353 if (type) { *type = gfield->petsc_type; } 354 355 PetscFunctionReturn(0); 356 } 357 358 #undef __FUNCT__ 359 #define __FUNCT__ "DMSwarmRestoreField" 360 PETSC_EXTERN PetscErrorCode DMSwarmRestoreField(DM dm,const char fieldname[],PetscInt *blocksize,PetscDataType *type,void **data) 361 { 362 DM_Swarm *swarm = (DM_Swarm*)dm->data; 363 DataField gfield; 364 PetscErrorCode ierr; 365 366 ierr = DataBucketGetDataFieldByName(swarm->db,fieldname,&gfield);CHKERRQ(ierr); 367 ierr = DataFieldRestoreAccess(gfield);CHKERRQ(ierr); 368 if (data) *data = NULL; 369 370 PetscFunctionReturn(0); 371 } 372 373 #undef __FUNCT__ 374 #define __FUNCT__ "DMSwarmAddPoint" 375 PETSC_EXTERN PetscErrorCode DMSwarmAddPoint(DM dm) 376 { 377 DM_Swarm *swarm = (DM_Swarm*)dm->data; 378 PetscErrorCode ierr; 379 380 if (!swarm->issetup) { ierr = DMSetUp(dm);CHKERRQ(ierr); } 381 ierr = DataBucketAddPoint(swarm->db);CHKERRQ(ierr); 382 PetscFunctionReturn(0); 383 } 384 385 #undef __FUNCT__ 386 #define __FUNCT__ "DMSwarmAddNPoints" 387 PETSC_EXTERN PetscErrorCode DMSwarmAddNPoints(DM dm,PetscInt npoints) 388 { 389 DM_Swarm *swarm = (DM_Swarm*)dm->data; 390 PetscErrorCode ierr; 391 PetscInt nlocal; 392 393 ierr = DataBucketGetSizes(swarm->db,&nlocal,NULL,NULL);CHKERRQ(ierr); 394 nlocal = nlocal + npoints; 395 ierr = DataBucketSetSizes(swarm->db,nlocal,-1);CHKERRQ(ierr); 396 PetscFunctionReturn(0); 397 } 398 399 #undef __FUNCT__ 400 #define __FUNCT__ "DMSwarmRemovePoint" 401 PETSC_EXTERN PetscErrorCode DMSwarmRemovePoint(DM dm) 402 { 403 DM_Swarm *swarm = (DM_Swarm*)dm->data; 404 PetscErrorCode ierr; 405 406 ierr = DataBucketRemovePoint(swarm->db);CHKERRQ(ierr); 407 PetscFunctionReturn(0); 408 } 409 410 #undef __FUNCT__ 411 #define __FUNCT__ "DMSwarmRemovePointAtIndex" 412 PETSC_EXTERN PetscErrorCode DMSwarmRemovePointAtIndex(DM dm,PetscInt idx) 413 { 414 DM_Swarm *swarm = (DM_Swarm*)dm->data; 415 PetscErrorCode ierr; 416 417 ierr = DataBucketRemovePointAtIndex(swarm->db,idx);CHKERRQ(ierr); 418 PetscFunctionReturn(0); 419 } 420 421 #undef __FUNCT__ 422 #define __FUNCT__ "DMSwarmMigrate_Basic" 423 PetscErrorCode DMSwarmMigrate_Basic(DM dm,PetscBool remove_sent_points) 424 { 425 PetscErrorCode ierr; 426 ierr = DMSwarmMigrate_Push_Basic(dm,remove_sent_points);CHKERRQ(ierr); 427 PetscFunctionReturn(0); 428 } 429 430 #undef __FUNCT__ 431 #define __FUNCT__ "DMSwarmMigrate" 432 PETSC_EXTERN PetscErrorCode DMSwarmMigrate(DM dm,PetscBool remove_sent_points) 433 { 434 PetscErrorCode ierr; 435 ierr = DMSwarmMigrate_Basic(dm,remove_sent_points);CHKERRQ(ierr); 436 PetscFunctionReturn(0); 437 } 438 439 #undef __FUNCT__ 440 #define __FUNCT__ "DMSwarmCollectViewCreate" 441 PETSC_EXTERN PetscErrorCode DMSwarmCollectViewCreate(DM dm) 442 { 443 PetscErrorCode ierr; 444 DM_Swarm *swarm = (DM_Swarm*)dm->data; 445 PetscInt ng; 446 447 if (swarm->collect_view_active) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"CollectView currently active"); 448 449 ierr = DMSwarmGetLocalSize(dm,&ng);CHKERRQ(ierr); 450 switch (swarm->collect_type) { 451 case DMSWARM_COLLECT_BASIC: 452 ierr = DMSwarmMigrate_GlobalToLocal_Basic(dm,&ng);CHKERRQ(ierr); 453 break; 454 case DMSWARM_COLLECT_DMDABOUNDINGBOX: 455 //ierr = DMSwarmCollect_DMDABoundingBox(dm,&ng);CHKERRQ(ierr); 456 break; 457 case DMSWARM_COLLECT_GENERAL: 458 //ierr = DMSwarmCollect_General(dm,..,,..,&ng);CHKERRQ(ierr); 459 break; 460 461 default: 462 break; 463 } 464 465 swarm->collect_view_active = PETSC_TRUE; 466 swarm->collect_view_reset_nlocal = ng; 467 468 PetscFunctionReturn(0); 469 } 470 471 #undef __FUNCT__ 472 #define __FUNCT__ "DMSwarmCollectViewDestroy" 473 PETSC_EXTERN PetscErrorCode DMSwarmCollectViewDestroy(DM dm) 474 { 475 PetscErrorCode ierr; 476 DM_Swarm *swarm = (DM_Swarm*)dm->data; 477 478 if (!swarm->collect_view_active) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"CollectView is currently not active"); 479 ierr = DMSwarmSetLocalSizes(dm,swarm->collect_view_reset_nlocal,-1);CHKERRQ(ierr); 480 swarm->collect_view_active = PETSC_FALSE; 481 482 PetscFunctionReturn(0); 483 } 484 485 #undef __FUNCT__ 486 #define __FUNCT__ "DMSetup_Swarm" 487 PetscErrorCode DMSetup_Swarm(DM dm) 488 { 489 DM_Swarm *swarm = (DM_Swarm*)dm->data; 490 PetscErrorCode ierr; 491 PetscMPIInt rank; 492 PetscInt p,npoints,*rankval; 493 494 if (swarm->issetup) PetscFunctionReturn(0); 495 496 swarm->issetup = PETSC_TRUE; 497 498 /* check some fields were registered */ 499 if (swarm->db->nfields <= 2) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"At least one field user must be registered via DMSwarmRegisterXXX()"); 500 501 /* check local sizes were set */ 502 if (swarm->db->L == -1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Local sizes must be set via DMSwarmSetLocalSizes()"); 503 504 /* initialize values in pid and rank placeholders */ 505 /* TODO: [pid - use MPI_Scan] */ 506 507 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr); 508 ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 509 ierr = DMSwarmGetField(dm,"DMSwarm_rank",NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 510 for (p=0; p<npoints; p++) { 511 rankval[p] = (PetscInt)rank; 512 } 513 ierr = DMSwarmRestoreField(dm,"DMSwarm_rank",NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 514 515 PetscFunctionReturn(0); 516 } 517 518 #undef __FUNCT__ 519 #define __FUNCT__ "DMDestroy_Swarm" 520 PetscErrorCode DMDestroy_Swarm(DM dm) 521 { 522 DM_Swarm *swarm = (DM_Swarm*)dm->data; 523 PetscErrorCode ierr; 524 525 PetscFunctionBegin; 526 ierr = DataBucketDestroy(&swarm->db);CHKERRQ(ierr); 527 ierr = PetscFree(swarm);CHKERRQ(ierr); 528 PetscFunctionReturn(0); 529 } 530 531 #undef __FUNCT__ 532 #define __FUNCT__ "DMView_Swarm" 533 PetscErrorCode DMView_Swarm(DM dm, PetscViewer viewer) 534 { 535 DM_Swarm *swarm = (DM_Swarm*)dm->data; 536 PetscBool iascii,ibinary,ishdf5,isvtk; 537 PetscErrorCode ierr; 538 539 PetscFunctionBegin; 540 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 541 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 542 ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 543 ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY,&ibinary);CHKERRQ(ierr); 544 ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERVTK, &isvtk);CHKERRQ(ierr); 545 ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5);CHKERRQ(ierr); 546 if (iascii) { 547 ierr = DataBucketView(PetscObjectComm((PetscObject)dm),swarm->db,NULL,DATABUCKET_VIEW_STDOUT);CHKERRQ(ierr); 548 } else if (ibinary) { 549 SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"NO VTK support"); 550 } else if (ishdf5) { 551 #if defined(PETSC_HAVE_HDF5) 552 SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"NO HDF5 support"); 553 #else 554 SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"HDF5 not supported. Please reconfigure using --download-hdf5"); 555 #endif 556 } else if (isvtk) { 557 SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"NO VTK support"); 558 } 559 PetscFunctionReturn(0); 560 } 561 562 #undef __FUNCT__ 563 #define __FUNCT__ "DMCreate_Swarm" 564 PETSC_EXTERN PetscErrorCode DMCreate_Swarm(DM dm) 565 { 566 DM_Swarm *swarm; 567 PetscErrorCode ierr; 568 569 PetscFunctionBegin; 570 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 571 ierr = PetscNewLog(dm,&swarm);CHKERRQ(ierr); 572 573 ierr = DataBucketCreate(&swarm->db);CHKERRQ(ierr); 574 swarm->vec_field_set = PETSC_FALSE; 575 swarm->issetup = PETSC_FALSE; 576 swarm->swarm_type = DMSWARM_BASIC; 577 swarm->migrate_type = DMSWARM_MIGRATE_BASIC; 578 swarm->collect_type = DMSWARM_COLLECT_BASIC; 579 580 dm->dim = 0; 581 dm->data = swarm; 582 583 dm->ops->view = DMView_Swarm; 584 dm->ops->load = NULL; 585 dm->ops->setfromoptions = NULL; 586 dm->ops->clone = NULL; 587 dm->ops->setup = DMSetup_Swarm; 588 dm->ops->createdefaultsection = NULL; 589 dm->ops->createdefaultconstraints = NULL; 590 dm->ops->createglobalvector = DMCreateGlobalVector_Swarm; 591 dm->ops->createlocalvector = DMCreateLocalVector_Swarm; 592 dm->ops->getlocaltoglobalmapping = NULL; 593 dm->ops->createfieldis = NULL; 594 dm->ops->createcoordinatedm = NULL; 595 dm->ops->getcoloring = NULL; 596 dm->ops->creatematrix = NULL; 597 dm->ops->createinterpolation = NULL; 598 dm->ops->getaggregates = NULL; 599 dm->ops->getinjection = NULL; 600 dm->ops->refine = NULL; 601 dm->ops->coarsen = NULL; 602 dm->ops->refinehierarchy = NULL; 603 dm->ops->coarsenhierarchy = NULL; 604 dm->ops->globaltolocalbegin = NULL; 605 dm->ops->globaltolocalend = NULL; 606 dm->ops->localtoglobalbegin = NULL; 607 dm->ops->localtoglobalend = NULL; 608 dm->ops->destroy = DMDestroy_Swarm; 609 dm->ops->createsubdm = NULL; 610 dm->ops->getdimpoints = NULL; 611 dm->ops->locatepoints = NULL; 612 613 PetscFunctionReturn(0); 614 }