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 = 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*swarm->vec_field_bs,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->vec_field_bs,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),bs,n*bs,array,&x);CHKERRQ(ierr); 127 } else { 128 ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)dm),bs,n*bs,PETSC_DETERMINE,array,&x);CHKERRQ(ierr); 129 } 130 PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"DMSwarmSharedField_%s",fieldname); 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,blocksize*size,NULL);CHKERRQ(ierr); 308 { 309 DataField gfield; 310 311 ierr = DataBucketGetDataFieldByName(swarm->db,fieldname,&gfield);CHKERRQ(ierr); 312 ierr = DataFieldSetBlockSize(gfield,blocksize);CHKERRQ(ierr); 313 } 314 swarm->db->field[swarm->db->nfields-1]->petsc_type = type; 315 316 PetscFunctionReturn(0); 317 } 318 319 #undef __FUNCT__ 320 #define __FUNCT__ "DMSwarmRegisterUserStructField" 321 PETSC_EXTERN PetscErrorCode DMSwarmRegisterUserStructField(DM dm,const char fieldname[],size_t size) 322 { 323 PetscErrorCode ierr; 324 DM_Swarm *swarm = (DM_Swarm*)dm->data; 325 326 ierr = DataBucketRegisterField(swarm->db,"DMSwarmRegisterUserStructField",fieldname,size,NULL);CHKERRQ(ierr); 327 swarm->db->field[swarm->db->nfields-1]->petsc_type = PETSC_STRUCT ; 328 329 PetscFunctionReturn(0); 330 } 331 332 #undef __FUNCT__ 333 #define __FUNCT__ "DMSwarmRegisterUserDatatypeField" 334 PETSC_EXTERN PetscErrorCode DMSwarmRegisterUserDatatypeField(DM dm,const char fieldname[],size_t size,PetscInt blocksize) 335 { 336 DM_Swarm *swarm = (DM_Swarm*)dm->data; 337 PetscErrorCode ierr; 338 339 ierr = DataBucketRegisterField(swarm->db,"DMSwarmRegisterUserDatatypeField",fieldname,blocksize*size,NULL);CHKERRQ(ierr); 340 { 341 DataField gfield; 342 343 ierr = DataBucketGetDataFieldByName(swarm->db,fieldname,&gfield);CHKERRQ(ierr); 344 ierr = DataFieldSetBlockSize(gfield,blocksize);CHKERRQ(ierr); 345 } 346 swarm->db->field[swarm->db->nfields-1]->petsc_type = PETSC_DATATYPE_UNKNOWN; 347 348 PetscFunctionReturn(0); 349 } 350 351 #undef __FUNCT__ 352 #define __FUNCT__ "DMSwarmGetField" 353 PETSC_EXTERN PetscErrorCode DMSwarmGetField(DM dm,const char fieldname[],PetscInt *blocksize,PetscDataType *type,void **data) 354 { 355 DM_Swarm *swarm = (DM_Swarm*)dm->data; 356 DataField gfield; 357 PetscErrorCode ierr; 358 359 if (!swarm->issetup) { ierr = DMSetUp(dm);CHKERRQ(ierr); } 360 361 ierr = DataBucketGetDataFieldByName(swarm->db,fieldname,&gfield);CHKERRQ(ierr); 362 ierr = DataFieldGetAccess(gfield);CHKERRQ(ierr); 363 ierr = DataFieldGetEntries(gfield,data);CHKERRQ(ierr); 364 if (blocksize) {*blocksize = gfield->bs; } 365 if (type) { *type = gfield->petsc_type; } 366 367 PetscFunctionReturn(0); 368 } 369 370 #undef __FUNCT__ 371 #define __FUNCT__ "DMSwarmRestoreField" 372 PETSC_EXTERN PetscErrorCode DMSwarmRestoreField(DM dm,const char fieldname[],PetscInt *blocksize,PetscDataType *type,void **data) 373 { 374 DM_Swarm *swarm = (DM_Swarm*)dm->data; 375 DataField gfield; 376 PetscErrorCode ierr; 377 378 ierr = DataBucketGetDataFieldByName(swarm->db,fieldname,&gfield);CHKERRQ(ierr); 379 ierr = DataFieldRestoreAccess(gfield);CHKERRQ(ierr); 380 if (data) *data = NULL; 381 382 PetscFunctionReturn(0); 383 } 384 385 #undef __FUNCT__ 386 #define __FUNCT__ "DMSwarmAddPoint" 387 PETSC_EXTERN PetscErrorCode DMSwarmAddPoint(DM dm) 388 { 389 DM_Swarm *swarm = (DM_Swarm*)dm->data; 390 PetscErrorCode ierr; 391 392 if (!swarm->issetup) { ierr = DMSetUp(dm);CHKERRQ(ierr); } 393 ierr = DataBucketAddPoint(swarm->db);CHKERRQ(ierr); 394 PetscFunctionReturn(0); 395 } 396 397 #undef __FUNCT__ 398 #define __FUNCT__ "DMSwarmAddNPoints" 399 PETSC_EXTERN PetscErrorCode DMSwarmAddNPoints(DM dm,PetscInt npoints) 400 { 401 DM_Swarm *swarm = (DM_Swarm*)dm->data; 402 PetscErrorCode ierr; 403 PetscInt nlocal; 404 405 ierr = DataBucketGetSizes(swarm->db,&nlocal,NULL,NULL);CHKERRQ(ierr); 406 nlocal = nlocal + npoints; 407 ierr = DataBucketSetSizes(swarm->db,nlocal,-1);CHKERRQ(ierr); 408 PetscFunctionReturn(0); 409 } 410 411 #undef __FUNCT__ 412 #define __FUNCT__ "DMSwarmRemovePoint" 413 PETSC_EXTERN PetscErrorCode DMSwarmRemovePoint(DM dm) 414 { 415 DM_Swarm *swarm = (DM_Swarm*)dm->data; 416 PetscErrorCode ierr; 417 418 ierr = DataBucketRemovePoint(swarm->db);CHKERRQ(ierr); 419 PetscFunctionReturn(0); 420 } 421 422 #undef __FUNCT__ 423 #define __FUNCT__ "DMSwarmRemovePointAtIndex" 424 PETSC_EXTERN PetscErrorCode DMSwarmRemovePointAtIndex(DM dm,PetscInt idx) 425 { 426 DM_Swarm *swarm = (DM_Swarm*)dm->data; 427 PetscErrorCode ierr; 428 429 ierr = DataBucketRemovePointAtIndex(swarm->db,idx);CHKERRQ(ierr); 430 PetscFunctionReturn(0); 431 } 432 433 #undef __FUNCT__ 434 #define __FUNCT__ "DMSwarmMigrate_Basic" 435 PetscErrorCode DMSwarmMigrate_Basic(DM dm,PetscBool remove_sent_points) 436 { 437 PetscErrorCode ierr; 438 ierr = DMSwarmMigrate_Push_Basic(dm,remove_sent_points);CHKERRQ(ierr); 439 PetscFunctionReturn(0); 440 } 441 442 #undef __FUNCT__ 443 #define __FUNCT__ "DMSwarmMigrate" 444 PETSC_EXTERN PetscErrorCode DMSwarmMigrate(DM dm,PetscBool remove_sent_points) 445 { 446 PetscErrorCode ierr; 447 ierr = DMSwarmMigrate_Basic(dm,remove_sent_points);CHKERRQ(ierr); 448 PetscFunctionReturn(0); 449 } 450 451 #undef __FUNCT__ 452 #define __FUNCT__ "DMSwarmCollectViewCreate" 453 PETSC_EXTERN PetscErrorCode DMSwarmCollectViewCreate(DM dm) 454 { 455 PetscErrorCode ierr; 456 DM_Swarm *swarm = (DM_Swarm*)dm->data; 457 PetscInt ng; 458 459 if (swarm->collect_view_active) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"CollectView currently active"); 460 461 ierr = DMSwarmGetLocalSize(dm,&ng);CHKERRQ(ierr); 462 switch (swarm->collect_type) { 463 case DMSWARM_COLLECT_BASIC: 464 ierr = DMSwarmMigrate_GlobalToLocal_Basic(dm,&ng);CHKERRQ(ierr); 465 break; 466 case DMSWARM_COLLECT_DMDABOUNDINGBOX: 467 //ierr = DMSwarmCollect_DMDABoundingBox(dm,&ng);CHKERRQ(ierr); 468 break; 469 case DMSWARM_COLLECT_GENERAL: 470 //ierr = DMSwarmCollect_General(dm,..,,..,&ng);CHKERRQ(ierr); 471 break; 472 473 default: 474 break; 475 } 476 477 swarm->collect_view_active = PETSC_TRUE; 478 swarm->collect_view_reset_nlocal = ng; 479 480 PetscFunctionReturn(0); 481 } 482 483 #undef __FUNCT__ 484 #define __FUNCT__ "DMSwarmCollectViewDestroy" 485 PETSC_EXTERN PetscErrorCode DMSwarmCollectViewDestroy(DM dm) 486 { 487 PetscErrorCode ierr; 488 DM_Swarm *swarm = (DM_Swarm*)dm->data; 489 490 if (!swarm->collect_view_active) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"CollectView is currently not active"); 491 ierr = DMSwarmSetLocalSizes(dm,swarm->collect_view_reset_nlocal,-1);CHKERRQ(ierr); 492 swarm->collect_view_active = PETSC_FALSE; 493 494 PetscFunctionReturn(0); 495 } 496 497 #undef __FUNCT__ 498 #define __FUNCT__ "DMSetup_Swarm" 499 PetscErrorCode DMSetup_Swarm(DM dm) 500 { 501 DM_Swarm *swarm = (DM_Swarm*)dm->data; 502 PetscErrorCode ierr; 503 PetscMPIInt rank; 504 PetscInt p,npoints,*rankval; 505 506 if (swarm->issetup) PetscFunctionReturn(0); 507 508 swarm->issetup = PETSC_TRUE; 509 510 /* check some fields were registered */ 511 if (swarm->db->nfields <= 2) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"At least one field user must be registered via DMSwarmRegisterXXX()"); 512 513 /* check local sizes were set */ 514 if (swarm->db->L == -1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Local sizes must be set via DMSwarmSetLocalSizes()"); 515 516 /* initialize values in pid and rank placeholders */ 517 /* TODO: [pid - use MPI_Scan] */ 518 519 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr); 520 ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 521 ierr = DMSwarmGetField(dm,"DMSwarm_rank",NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 522 for (p=0; p<npoints; p++) { 523 rankval[p] = (PetscInt)rank; 524 } 525 ierr = DMSwarmRestoreField(dm,"DMSwarm_rank",NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 526 527 PetscFunctionReturn(0); 528 } 529 530 #undef __FUNCT__ 531 #define __FUNCT__ "DMDestroy_Swarm" 532 PetscErrorCode DMDestroy_Swarm(DM dm) 533 { 534 DM_Swarm *swarm = (DM_Swarm*)dm->data; 535 PetscErrorCode ierr; 536 537 PetscFunctionBegin; 538 ierr = DataBucketDestroy(&swarm->db);CHKERRQ(ierr); 539 ierr = PetscFree(swarm);CHKERRQ(ierr); 540 PetscFunctionReturn(0); 541 } 542 543 #undef __FUNCT__ 544 #define __FUNCT__ "DMView_Swarm" 545 PetscErrorCode DMView_Swarm(DM dm, PetscViewer viewer) 546 { 547 DM_Swarm *swarm = (DM_Swarm*)dm->data; 548 PetscBool iascii,ibinary,ishdf5,isvtk; 549 PetscErrorCode ierr; 550 551 PetscFunctionBegin; 552 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 553 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 554 ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 555 ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY,&ibinary);CHKERRQ(ierr); 556 ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERVTK, &isvtk);CHKERRQ(ierr); 557 ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5);CHKERRQ(ierr); 558 if (iascii) { 559 ierr = DataBucketView(PetscObjectComm((PetscObject)dm),swarm->db,NULL,DATABUCKET_VIEW_STDOUT);CHKERRQ(ierr); 560 } else if (ibinary) { 561 SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"NO VTK support"); 562 } else if (ishdf5) { 563 #if defined(PETSC_HAVE_HDF5) 564 SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"NO HDF5 support"); 565 #else 566 SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"HDF5 not supported. Please reconfigure using --download-hdf5"); 567 #endif 568 } else if (isvtk) { 569 SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"NO VTK support"); 570 } 571 PetscFunctionReturn(0); 572 } 573 574 #undef __FUNCT__ 575 #define __FUNCT__ "DMCreate_Swarm" 576 PETSC_EXTERN PetscErrorCode DMCreate_Swarm(DM dm) 577 { 578 DM_Swarm *swarm; 579 PetscErrorCode ierr; 580 581 PetscFunctionBegin; 582 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 583 ierr = PetscNewLog(dm,&swarm);CHKERRQ(ierr); 584 585 ierr = DataBucketCreate(&swarm->db);CHKERRQ(ierr); 586 swarm->vec_field_set = PETSC_FALSE; 587 swarm->issetup = PETSC_FALSE; 588 swarm->swarm_type = DMSWARM_BASIC; 589 swarm->migrate_type = DMSWARM_MIGRATE_BASIC; 590 swarm->collect_type = DMSWARM_COLLECT_BASIC; 591 swarm->migrate_error_on_missing_point = PETSC_FALSE; 592 593 dm->dim = 0; 594 dm->data = swarm; 595 596 dm->ops->view = DMView_Swarm; 597 dm->ops->load = NULL; 598 dm->ops->setfromoptions = NULL; 599 dm->ops->clone = NULL; 600 dm->ops->setup = DMSetup_Swarm; 601 dm->ops->createdefaultsection = NULL; 602 dm->ops->createdefaultconstraints = NULL; 603 dm->ops->createglobalvector = DMCreateGlobalVector_Swarm; 604 dm->ops->createlocalvector = DMCreateLocalVector_Swarm; 605 dm->ops->getlocaltoglobalmapping = NULL; 606 dm->ops->createfieldis = NULL; 607 dm->ops->createcoordinatedm = NULL; 608 dm->ops->getcoloring = NULL; 609 dm->ops->creatematrix = NULL; 610 dm->ops->createinterpolation = NULL; 611 dm->ops->getaggregates = NULL; 612 dm->ops->getinjection = NULL; 613 dm->ops->refine = NULL; 614 dm->ops->coarsen = NULL; 615 dm->ops->refinehierarchy = NULL; 616 dm->ops->coarsenhierarchy = NULL; 617 dm->ops->globaltolocalbegin = NULL; 618 dm->ops->globaltolocalend = NULL; 619 dm->ops->localtoglobalbegin = NULL; 620 dm->ops->localtoglobalend = NULL; 621 dm->ops->destroy = DMDestroy_Swarm; 622 dm->ops->createsubdm = NULL; 623 dm->ops->getdimpoints = NULL; 624 dm->ops->locatepoints = NULL; 625 626 PetscFunctionReturn(0); 627 }