1 2 #include <private/characteristicimpl.h> /*I "petsccharacteristic.h" I*/ 3 4 PetscClassId CHARACTERISTIC_CLASSID; 5 PetscLogEvent CHARACTERISTIC_SetUp, CHARACTERISTIC_Solve, CHARACTERISTIC_QueueSetup, CHARACTERISTIC_DAUpdate; 6 PetscLogEvent CHARACTERISTIC_HalfTimeLocal, CHARACTERISTIC_HalfTimeRemote, CHARACTERISTIC_HalfTimeExchange; 7 PetscLogEvent CHARACTERISTIC_FullTimeLocal, CHARACTERISTIC_FullTimeRemote, CHARACTERISTIC_FullTimeExchange; 8 PetscBool CharacteristicRegisterAllCalled = PETSC_FALSE; 9 /* 10 Contains the list of registered characteristic routines 11 */ 12 PetscFList CharacteristicList = PETSC_NULL; 13 14 PetscErrorCode DMDAGetNeighborsRank(DM, PetscMPIInt []); 15 PetscInt DMDAGetNeighborRelative(DM, PassiveReal, PassiveReal); 16 PetscErrorCode DMDAMapToPeriodicDomain(DM, PetscScalar [] ); 17 18 PetscErrorCode HeapSort(Characteristic, Queue, PetscInt); 19 PetscErrorCode SiftDown(Characteristic, Queue, PetscInt, PetscInt); 20 21 #undef __FUNCT__ 22 #define __FUNCT__ "CharacteristicView" 23 PetscErrorCode CharacteristicView(Characteristic c, PetscViewer viewer) 24 { 25 PetscBool iascii; 26 PetscErrorCode ierr; 27 28 PetscFunctionBegin; 29 PetscValidHeaderSpecific(c, CHARACTERISTIC_CLASSID, 1); 30 if (!viewer) { 31 ierr = PetscViewerASCIIGetStdout(((PetscObject)c)->comm,&viewer);CHKERRQ(ierr); 32 } 33 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 34 PetscCheckSameComm(c, 1, viewer, 2); 35 36 ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 37 if (!iascii) { 38 if (c->ops->view) { 39 ierr = (*c->ops->view)(c, viewer);CHKERRQ(ierr); 40 } 41 } 42 PetscFunctionReturn(0); 43 } 44 45 #undef __FUNCT__ 46 #define __FUNCT__ "CharacteristicDestroy" 47 PetscErrorCode CharacteristicDestroy(Characteristic *c) 48 { 49 PetscErrorCode ierr; 50 51 PetscFunctionBegin; 52 if (!*c) PetscFunctionReturn(0); 53 PetscValidHeaderSpecific(*c, CHARACTERISTIC_CLASSID, 1); 54 if (--((PetscObject)(*c))->refct > 0) PetscFunctionReturn(0); 55 56 if ((*c)->ops->destroy) { 57 ierr = (*(*c)->ops->destroy)((*c));CHKERRQ(ierr); 58 } 59 ierr = MPI_Type_free(&(*c)->itemType);CHKERRQ(ierr); 60 ierr = PetscFree((*c)->queue);CHKERRQ(ierr); 61 ierr = PetscFree((*c)->queueLocal);CHKERRQ(ierr); 62 ierr = PetscFree((*c)->queueRemote);CHKERRQ(ierr); 63 ierr = PetscFree((*c)->neighbors);CHKERRQ(ierr); 64 ierr = PetscFree((*c)->needCount);CHKERRQ(ierr); 65 ierr = PetscFree((*c)->localOffsets);CHKERRQ(ierr); 66 ierr = PetscFree((*c)->fillCount);CHKERRQ(ierr); 67 ierr = PetscFree((*c)->remoteOffsets);CHKERRQ(ierr); 68 ierr = PetscFree((*c)->request);CHKERRQ(ierr); 69 ierr = PetscFree((*c)->status);CHKERRQ(ierr); 70 ierr = PetscHeaderDestroy(c);CHKERRQ(ierr); 71 PetscFunctionReturn(0); 72 } 73 74 #undef __FUNCT__ 75 #define __FUNCT__ "CharacteristicCreate" 76 PetscErrorCode CharacteristicCreate(MPI_Comm comm, Characteristic *c) 77 { 78 Characteristic newC; 79 PetscErrorCode ierr; 80 81 PetscFunctionBegin; 82 PetscValidPointer(c, 2); 83 *c = PETSC_NULL; 84 #ifndef PETSC_USE_DYNAMIC_LIBRARIES 85 ierr = CharacteristicInitializePackage(PETSC_NULL);CHKERRQ(ierr); 86 #endif 87 88 ierr = PetscHeaderCreate(newC, _p_Characteristic, struct _CharacteristicOps, CHARACTERISTIC_CLASSID, -1, "Characteristic", comm, CharacteristicDestroy, CharacteristicView);CHKERRQ(ierr); 89 ierr = PetscLogObjectCreate(newC);CHKERRQ(ierr); 90 *c = newC; 91 92 newC->structured = PETSC_TRUE; 93 newC->numIds = 0; 94 newC->velocityDA = PETSC_NULL; 95 newC->velocity = PETSC_NULL; 96 newC->velocityOld = PETSC_NULL; 97 newC->numVelocityComp = 0; 98 newC->velocityComp = PETSC_NULL; 99 newC->velocityInterp = PETSC_NULL; 100 newC->velocityInterpLocal = PETSC_NULL; 101 newC->velocityCtx = PETSC_NULL; 102 newC->fieldDA = PETSC_NULL; 103 newC->field = PETSC_NULL; 104 newC->numFieldComp = 0; 105 newC->fieldComp = PETSC_NULL; 106 newC->fieldInterp = PETSC_NULL; 107 newC->fieldInterpLocal = PETSC_NULL; 108 newC->fieldCtx = PETSC_NULL; 109 newC->itemType = PETSC_NULL; 110 newC->queue = PETSC_NULL; 111 newC->queueSize = 0; 112 newC->queueMax = 0; 113 newC->queueLocal = PETSC_NULL; 114 newC->queueLocalSize = 0; 115 newC->queueLocalMax = 0; 116 newC->queueRemote = PETSC_NULL; 117 newC->queueRemoteSize = 0; 118 newC->queueRemoteMax = 0; 119 newC->numNeighbors = 0; 120 newC->neighbors = PETSC_NULL; 121 newC->needCount = PETSC_NULL; 122 newC->localOffsets = PETSC_NULL; 123 newC->fillCount = PETSC_NULL; 124 newC->remoteOffsets = PETSC_NULL; 125 newC->request = PETSC_NULL; 126 newC->status = PETSC_NULL; 127 PetscFunctionReturn(0); 128 } 129 130 #undef __FUNCT__ 131 #define __FUNCT__ "CharacteristicSetType" 132 /*@C 133 CharacteristicSetType - Builds Characteristic for a particular solver. 134 135 Logically Collective on Characteristic 136 137 Input Parameters: 138 + c - the method of characteristics context 139 - type - a known method 140 141 Options Database Key: 142 . -characteristic_type <method> - Sets the method; use -help for a list 143 of available methods 144 145 Notes: 146 See "include/petsccharacteristic.h" for available methods 147 148 Normally, it is best to use the CharacteristicSetFromOptions() command and 149 then set the Characteristic type from the options database rather than by using 150 this routine. Using the options database provides the user with 151 maximum flexibility in evaluating the many different Krylov methods. 152 The CharacteristicSetType() routine is provided for those situations where it 153 is necessary to set the iterative solver independently of the command 154 line or options database. This might be the case, for example, when 155 the choice of iterative solver changes during the execution of the 156 program, and the user's application is taking responsibility for 157 choosing the appropriate method. In other words, this routine is 158 not for beginners. 159 160 Level: intermediate 161 162 .keywords: Characteristic, set, method 163 164 .seealso: CharacteristicType 165 166 @*/ 167 PetscErrorCode CharacteristicSetType(Characteristic c, const CharacteristicType type) 168 { 169 PetscErrorCode ierr, (*r)(Characteristic); 170 PetscBool match; 171 172 PetscFunctionBegin; 173 PetscValidHeaderSpecific(c, CHARACTERISTIC_CLASSID, 1); 174 PetscValidCharPointer(type, 2); 175 176 ierr = PetscTypeCompare((PetscObject) c, type, &match);CHKERRQ(ierr); 177 if (match) PetscFunctionReturn(0); 178 179 if (c->data) { 180 /* destroy the old private Characteristic context */ 181 ierr = (*c->ops->destroy)(c);CHKERRQ(ierr); 182 c->data = 0; 183 } 184 185 ierr = PetscFListFind(CharacteristicList, ((PetscObject)c)->comm,type,PETSC_TRUE, (void (**)(void)) &r);CHKERRQ(ierr); 186 if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown Characteristic type given: %s", type); 187 c->setupcalled = 0; 188 ierr = (*r)(c);CHKERRQ(ierr); 189 ierr = PetscObjectChangeTypeName((PetscObject) c, type);CHKERRQ(ierr); 190 PetscFunctionReturn(0); 191 } 192 193 #undef __FUNCT__ 194 #define __FUNCT__ "CharacteristicSetUp" 195 /*@ 196 CharacteristicSetUp - Sets up the internal data structures for the 197 later use of an iterative solver. 198 199 Collective on Characteristic 200 201 Input Parameter: 202 . ksp - iterative context obtained from CharacteristicCreate() 203 204 Level: developer 205 206 .keywords: Characteristic, setup 207 208 .seealso: CharacteristicCreate(), CharacteristicSolve(), CharacteristicDestroy() 209 @*/ 210 PetscErrorCode CharacteristicSetUp(Characteristic c) 211 { 212 PetscErrorCode ierr; 213 214 PetscFunctionBegin; 215 PetscValidHeaderSpecific(c, CHARACTERISTIC_CLASSID, 1); 216 217 if (!((PetscObject)c)->type_name){ 218 ierr = CharacteristicSetType(c, CHARACTERISTICDA);CHKERRQ(ierr); 219 } 220 221 if (c->setupcalled == 2) PetscFunctionReturn(0); 222 223 ierr = PetscLogEventBegin(CHARACTERISTIC_SetUp,c,0,0,0);CHKERRQ(ierr); 224 if (!c->setupcalled) { 225 ierr = (*c->ops->setup)(c);CHKERRQ(ierr); 226 } 227 ierr = PetscLogEventEnd(CHARACTERISTIC_SetUp,c,0,0,0);CHKERRQ(ierr); 228 c->setupcalled = 2; 229 PetscFunctionReturn(0); 230 } 231 232 #undef __FUNCT__ 233 #define __FUNCT__ "CharacteristicRegister" 234 /*@C 235 CharacteristicRegister - See CharacteristicRegisterDynamic() 236 237 Level: advanced 238 @*/ 239 PetscErrorCode CharacteristicRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(Characteristic)) 240 { 241 PetscErrorCode ierr; 242 char fullname[PETSC_MAX_PATH_LEN]; 243 244 PetscFunctionBegin; 245 ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr); 246 ierr = PetscFListAdd(&CharacteristicList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr); 247 PetscFunctionReturn(0); 248 } 249 250 #undef __FUNCT__ 251 #define __FUNCT__ "CharacteristicSetVelocityInterpolation" 252 PetscErrorCode CharacteristicSetVelocityInterpolation(Characteristic c, DM da, Vec v, Vec vOld, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(Vec, PetscReal [], PetscInt, PetscInt [], PetscScalar [], void *), void *ctx) 253 { 254 PetscFunctionBegin; 255 c->velocityDA = da; 256 c->velocity = v; 257 c->velocityOld = vOld; 258 c->numVelocityComp = numComponents; 259 c->velocityComp = components; 260 c->velocityInterp = interp; 261 c->velocityCtx = ctx; 262 PetscFunctionReturn(0); 263 } 264 265 #undef __FUNCT__ 266 #define __FUNCT__ "CharacteristicSetVelocityInterpolationLocal" 267 PetscErrorCode CharacteristicSetVelocityInterpolationLocal(Characteristic c, DM da, Vec v, Vec vOld, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(void *, PetscReal [], PetscInt, PetscInt [], PetscScalar [], void *), void *ctx) 268 { 269 PetscFunctionBegin; 270 c->velocityDA = da; 271 c->velocity = v; 272 c->velocityOld = vOld; 273 c->numVelocityComp = numComponents; 274 c->velocityComp = components; 275 c->velocityInterpLocal = interp; 276 c->velocityCtx = ctx; 277 PetscFunctionReturn(0); 278 } 279 280 #undef __FUNCT__ 281 #define __FUNCT__ "CharacteristicSetFieldInterpolation" 282 PetscErrorCode CharacteristicSetFieldInterpolation(Characteristic c, DM da, Vec v, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(Vec, PetscReal [], PetscInt, PetscInt [], PetscScalar [], void *), void *ctx) 283 { 284 PetscFunctionBegin; 285 #if 0 286 if (numComponents > 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "Fields with more than 2 components are not supported. Send mail to petsc-maint@mcs.anl.gov."); 287 #endif 288 c->fieldDA = da; 289 c->field = v; 290 c->numFieldComp = numComponents; 291 c->fieldComp = components; 292 c->fieldInterp = interp; 293 c->fieldCtx = ctx; 294 PetscFunctionReturn(0); 295 } 296 297 #undef __FUNCT__ 298 #define __FUNCT__ "CharacteristicSetFieldInterpolationLocal" 299 PetscErrorCode CharacteristicSetFieldInterpolationLocal(Characteristic c, DM da, Vec v, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(void *, PetscReal [], PetscInt, PetscInt [], PetscScalar [], void *), void *ctx) 300 { 301 PetscFunctionBegin; 302 #if 0 303 if (numComponents > 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "Fields with more than 2 components are not supported. Send mail to petsc-maint@mcs.anl.gov."); 304 #endif 305 c->fieldDA = da; 306 c->field = v; 307 c->numFieldComp = numComponents; 308 c->fieldComp = components; 309 c->fieldInterpLocal = interp; 310 c->fieldCtx = ctx; 311 PetscFunctionReturn(0); 312 } 313 314 #undef __FUNCT__ 315 #define __FUNCT__ "CharacteristicSolve" 316 PetscErrorCode CharacteristicSolve(Characteristic c, PetscReal dt, Vec solution) 317 { 318 CharacteristicPointDA2D Qi; 319 DM da = c->velocityDA; 320 Vec velocityLocal, velocityLocalOld; 321 Vec fieldLocal; 322 DMDALocalInfo info; 323 PetscScalar **solArray; 324 void *velocityArray; 325 void *velocityArrayOld; 326 void *fieldArray; 327 PassiveScalar *interpIndices; 328 PassiveScalar *velocityValues, *velocityValuesOld; 329 PassiveScalar *fieldValues; 330 PetscMPIInt rank; 331 PetscInt dim; 332 PetscMPIInt neighbors[9]; 333 PetscInt dof; 334 PetscInt gx, gy; 335 PetscInt n, is, ie, js, je, comp; 336 PetscErrorCode ierr; 337 338 PetscFunctionBegin; 339 c->queueSize = 0; 340 ierr = MPI_Comm_rank(((PetscObject)c)->comm, &rank);CHKERRQ(ierr); 341 ierr = DMDAGetNeighborsRank(da, neighbors);CHKERRQ(ierr); 342 ierr = CharacteristicSetNeighbors(c, 9, neighbors);CHKERRQ(ierr); 343 ierr = CharacteristicSetUp(c);CHKERRQ(ierr); 344 /* global and local grid info */ 345 ierr = DMDAGetInfo(da, &dim, &gx, &gy, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0);CHKERRQ(ierr); 346 ierr = DMDAGetLocalInfo(da, &info);CHKERRQ(ierr); 347 is = info.xs; ie = info.xs+info.xm; 348 js = info.ys; je = info.ys+info.ym; 349 /* Allocation */ 350 ierr = PetscMalloc(dim*sizeof(PetscScalar), &interpIndices);CHKERRQ(ierr); 351 ierr = PetscMalloc(c->numVelocityComp*sizeof(PetscScalar), &velocityValues);CHKERRQ(ierr); 352 ierr = PetscMalloc(c->numVelocityComp*sizeof(PetscScalar), &velocityValuesOld);CHKERRQ(ierr); 353 ierr = PetscMalloc(c->numFieldComp*sizeof(PetscScalar), &fieldValues);CHKERRQ(ierr); 354 ierr = PetscLogEventBegin(CHARACTERISTIC_Solve,0,0,0,0);CHKERRQ(ierr); 355 356 /* ----------------------------------------------------------------------- 357 PART 1, AT t-dt/2 358 -----------------------------------------------------------------------*/ 359 ierr = PetscLogEventBegin(CHARACTERISTIC_QueueSetup,0,0,0,0);CHKERRQ(ierr); 360 /* GET POSITION AT HALF TIME IN THE PAST */ 361 if (c->velocityInterpLocal) { 362 ierr = DMGetLocalVector(c->velocityDA, &velocityLocal);CHKERRQ(ierr); 363 ierr = DMGetLocalVector(c->velocityDA, &velocityLocalOld);CHKERRQ(ierr); 364 ierr = DMGlobalToLocalBegin(c->velocityDA, c->velocity, INSERT_VALUES, velocityLocal);CHKERRQ(ierr); 365 ierr = DMGlobalToLocalEnd(c->velocityDA, c->velocity, INSERT_VALUES, velocityLocal);CHKERRQ(ierr); 366 ierr = DMGlobalToLocalBegin(c->velocityDA, c->velocityOld, INSERT_VALUES, velocityLocalOld);CHKERRQ(ierr); 367 ierr = DMGlobalToLocalEnd(c->velocityDA, c->velocityOld, INSERT_VALUES, velocityLocalOld);CHKERRQ(ierr); 368 ierr = DMDAVecGetArray(c->velocityDA, velocityLocal, &velocityArray); CHKERRQ(ierr); 369 ierr = DMDAVecGetArray(c->velocityDA, velocityLocalOld, &velocityArrayOld);CHKERRQ(ierr); 370 } 371 ierr = PetscInfo(PETSC_NULL, "Calculating position at t_{n - 1/2}\n");CHKERRQ(ierr); 372 for(Qi.j = js; Qi.j < je; Qi.j++) { 373 for(Qi.i = is; Qi.i < ie; Qi.i++) { 374 interpIndices[0] = Qi.i; 375 interpIndices[1] = Qi.j; 376 if (c->velocityInterpLocal) { 377 c->velocityInterpLocal(velocityArray, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx); 378 } else { 379 c->velocityInterp(c->velocity, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx); 380 } 381 Qi.x = Qi.i - velocityValues[0]*dt/2.0; 382 Qi.y = Qi.j - velocityValues[1]*dt/2.0; 383 384 /* Determine whether the position at t - dt/2 is local */ 385 Qi.proc = DMDAGetNeighborRelative(da, Qi.x, Qi.y); 386 387 /* Check for Periodic boundaries and move all periodic points back onto the domain */ 388 ierr = DMDAMapCoordsToPeriodicDomain(da,&(Qi.x),&(Qi.y));CHKERRQ(ierr); 389 ierr = CharacteristicAddPoint(c, &Qi);CHKERRQ(ierr); 390 } 391 } 392 ierr = PetscLogEventEnd(CHARACTERISTIC_QueueSetup,0,0,0,0);CHKERRQ(ierr); 393 394 ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr); 395 ierr = CharacteristicSendCoordinatesBegin(c);CHKERRQ(ierr); 396 ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr); 397 398 ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeLocal,0,0,0,0);CHKERRQ(ierr); 399 /* Calculate velocity at t_n+1/2 (local values) */ 400 ierr = PetscInfo(PETSC_NULL, "Calculating local velocities at t_{n - 1/2}\n");CHKERRQ(ierr); 401 for(n = 0; n < c->queueSize; n++) { 402 Qi = c->queue[n]; 403 if (c->neighbors[Qi.proc] == rank) { 404 interpIndices[0] = Qi.x; 405 interpIndices[1] = Qi.y; 406 if (c->velocityInterpLocal) { 407 c->velocityInterpLocal(velocityArray, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx); 408 c->velocityInterpLocal(velocityArrayOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx); 409 } else { 410 c->velocityInterp(c->velocity, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx); 411 c->velocityInterp(c->velocityOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx); 412 } 413 Qi.x = 0.5*(velocityValues[0] + velocityValuesOld[0]); 414 Qi.y = 0.5*(velocityValues[1] + velocityValuesOld[1]); 415 } 416 c->queue[n] = Qi; 417 } 418 ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeLocal,0,0,0,0);CHKERRQ(ierr); 419 420 ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr); 421 ierr = CharacteristicSendCoordinatesEnd(c);CHKERRQ(ierr); 422 ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr); 423 424 425 /* Calculate velocity at t_n+1/2 (fill remote requests) */ 426 ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeRemote,0,0,0,0);CHKERRQ(ierr); 427 ierr = PetscInfo1(PETSC_NULL, "Calculating %d remote velocities at t_{n - 1/2}\n", c->queueRemoteSize);CHKERRQ(ierr); 428 for(n = 0; n < c->queueRemoteSize; n++) { 429 Qi = c->queueRemote[n]; 430 interpIndices[0] = Qi.x; 431 interpIndices[1] = Qi.y; 432 if (c->velocityInterpLocal) { 433 c->velocityInterpLocal(velocityArray, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx); 434 c->velocityInterpLocal(velocityArrayOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx); 435 } else { 436 c->velocityInterp(c->velocity, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx); 437 c->velocityInterp(c->velocityOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx); 438 } 439 Qi.x = 0.5*(velocityValues[0] + velocityValuesOld[0]); 440 Qi.y = 0.5*(velocityValues[1] + velocityValuesOld[1]); 441 c->queueRemote[n] = Qi; 442 } 443 ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeRemote,0,0,0,0);CHKERRQ(ierr); 444 ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr); 445 ierr = CharacteristicGetValuesBegin(c);CHKERRQ(ierr); 446 ierr = CharacteristicGetValuesEnd(c);CHKERRQ(ierr); 447 if (c->velocityInterpLocal) { 448 ierr = DMDAVecRestoreArray(c->velocityDA, velocityLocal, &velocityArray); CHKERRQ(ierr); 449 ierr = DMDAVecRestoreArray(c->velocityDA, velocityLocalOld, &velocityArrayOld);CHKERRQ(ierr); 450 ierr = DMRestoreLocalVector(c->velocityDA, &velocityLocal);CHKERRQ(ierr); 451 ierr = DMRestoreLocalVector(c->velocityDA, &velocityLocalOld);CHKERRQ(ierr); 452 } 453 ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr); 454 455 /* ----------------------------------------------------------------------- 456 PART 2, AT t-dt 457 -----------------------------------------------------------------------*/ 458 459 /* GET POSITION AT t_n (local values) */ 460 ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeLocal,0,0,0,0);CHKERRQ(ierr); 461 ierr = PetscInfo(PETSC_NULL, "Calculating position at t_{n}\n");CHKERRQ(ierr); 462 for(n = 0; n < c->queueSize; n++) { 463 Qi = c->queue[n]; 464 Qi.x = Qi.i - Qi.x*dt; 465 Qi.y = Qi.j - Qi.y*dt; 466 467 /* Determine whether the position at t-dt is local */ 468 Qi.proc = DMDAGetNeighborRelative(da, Qi.x, Qi.y); 469 470 /* Check for Periodic boundaries and move all periodic points back onto the domain */ 471 ierr = DMDAMapCoordsToPeriodicDomain(da,&(Qi.x),&(Qi.y));CHKERRQ(ierr); 472 473 c->queue[n] = Qi; 474 } 475 ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeLocal,0,0,0,0);CHKERRQ(ierr); 476 477 ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr); 478 ierr = CharacteristicSendCoordinatesBegin(c);CHKERRQ(ierr); 479 ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr); 480 481 /* GET VALUE AT FULL TIME IN THE PAST (LOCAL REQUESTS) */ 482 ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeLocal,0,0,0,0);CHKERRQ(ierr); 483 if (c->fieldInterpLocal) { 484 ierr = DMGetLocalVector(c->fieldDA, &fieldLocal);CHKERRQ(ierr); 485 ierr = DMGlobalToLocalBegin(c->fieldDA, c->field, INSERT_VALUES, fieldLocal);CHKERRQ(ierr); 486 ierr = DMGlobalToLocalEnd(c->fieldDA, c->field, INSERT_VALUES, fieldLocal);CHKERRQ(ierr); 487 ierr = DMDAVecGetArray(c->fieldDA, fieldLocal, &fieldArray);CHKERRQ(ierr); 488 } 489 ierr = PetscInfo(PETSC_NULL, "Calculating local field at t_{n}\n");CHKERRQ(ierr); 490 for(n = 0; n < c->queueSize; n++) { 491 if (c->neighbors[c->queue[n].proc] == rank) { 492 interpIndices[0] = c->queue[n].x; 493 interpIndices[1] = c->queue[n].y; 494 if (c->fieldInterpLocal) { 495 c->fieldInterpLocal(fieldArray, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx); 496 } else { 497 c->fieldInterp(c->field, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx); 498 } 499 for(comp = 0; comp < c->numFieldComp; comp++) { 500 c->queue[n].field[comp] = fieldValues[comp]; 501 } 502 } 503 } 504 ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeLocal,0,0,0,0);CHKERRQ(ierr); 505 506 ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr); 507 ierr = CharacteristicSendCoordinatesEnd(c);CHKERRQ(ierr); 508 ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr); 509 510 /* GET VALUE AT FULL TIME IN THE PAST (REMOTE REQUESTS) */ 511 ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeRemote,0,0,0,0);CHKERRQ(ierr); 512 ierr = PetscInfo1(PETSC_NULL, "Calculating %d remote field points at t_{n}\n", c->queueRemoteSize);CHKERRQ(ierr); 513 for(n = 0; n < c->queueRemoteSize; n++) { 514 interpIndices[0] = c->queueRemote[n].x; 515 interpIndices[1] = c->queueRemote[n].y; 516 517 /* for debugging purposes */ 518 if (1) { /* hacked bounds test...let's do better */ 519 PetscScalar im = interpIndices[0]; PetscScalar jm = interpIndices[1]; 520 521 if (( im < (PetscScalar) is - 1.) || (im > (PetscScalar) ie) || (jm < (PetscScalar) js - 1.) || (jm > (PetscScalar) je)) { 522 SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB, "Nonlocal point: (%g,%g)", im, jm); 523 } 524 } 525 526 if (c->fieldInterpLocal) { 527 c->fieldInterpLocal(fieldArray, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx); 528 } else { 529 c->fieldInterp(c->field, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx); 530 } 531 for(comp = 0; comp < c->numFieldComp; comp++) { 532 c->queueRemote[n].field[comp] = fieldValues[comp]; 533 } 534 } 535 ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeRemote,0,0,0,0);CHKERRQ(ierr); 536 537 ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr); 538 ierr = CharacteristicGetValuesBegin(c);CHKERRQ(ierr); 539 ierr = CharacteristicGetValuesEnd(c);CHKERRQ(ierr); 540 if (c->fieldInterpLocal) { 541 ierr = DMDAVecRestoreArray(c->fieldDA, fieldLocal, &fieldArray);CHKERRQ(ierr); 542 ierr = DMRestoreLocalVector(c->fieldDA, &fieldLocal);CHKERRQ(ierr); 543 } 544 ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr); 545 546 /* Return field of characteristics at t_n-1 */ 547 ierr = PetscLogEventBegin(CHARACTERISTIC_DAUpdate,0,0,0,0);CHKERRQ(ierr); 548 ierr = DMDAGetInfo(c->fieldDA,0,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); 549 ierr = DMDAVecGetArray(c->fieldDA, solution, &solArray);CHKERRQ(ierr); 550 for(n = 0; n < c->queueSize; n++) { 551 Qi = c->queue[n]; 552 for(comp = 0; comp < c->numFieldComp; comp++) { 553 solArray[Qi.j][Qi.i*dof+c->fieldComp[comp]] = Qi.field[comp]; 554 } 555 } 556 ierr = DMDAVecRestoreArray(c->fieldDA, solution, &solArray);CHKERRQ(ierr); 557 ierr = PetscLogEventEnd(CHARACTERISTIC_DAUpdate,0,0,0,0);CHKERRQ(ierr); 558 ierr = PetscLogEventEnd(CHARACTERISTIC_Solve,0,0,0,0);CHKERRQ(ierr); 559 560 /* Cleanup */ 561 ierr = PetscFree(interpIndices);CHKERRQ(ierr); 562 ierr = PetscFree(velocityValues);CHKERRQ(ierr); 563 ierr = PetscFree(velocityValuesOld);CHKERRQ(ierr); 564 ierr = PetscFree(fieldValues);CHKERRQ(ierr); 565 PetscFunctionReturn(0); 566 } 567 568 #undef __FUNCT__ 569 #define __FUNCT__ "CharacteristicSetNeighbors" 570 PetscErrorCode CharacteristicSetNeighbors(Characteristic c, PetscInt numNeighbors, PetscMPIInt neighbors[]) 571 { 572 PetscErrorCode ierr; 573 574 PetscFunctionBegin; 575 c->numNeighbors = numNeighbors; 576 ierr = PetscMalloc(numNeighbors * sizeof(PetscMPIInt), &c->neighbors);CHKERRQ(ierr); 577 ierr = PetscMemcpy(c->neighbors, neighbors, numNeighbors * sizeof(PetscMPIInt));CHKERRQ(ierr); 578 PetscFunctionReturn(0); 579 } 580 581 #undef __FUNCT__ 582 #define __FUNCT__ "CharacteristicAddPoint" 583 PetscErrorCode CharacteristicAddPoint(Characteristic c, CharacteristicPointDA2D *point) 584 { 585 PetscFunctionBegin; 586 if (c->queueSize >= c->queueMax) { 587 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Exceeeded maximum queue size %d", c->queueMax); 588 } 589 c->queue[c->queueSize++] = *point; 590 PetscFunctionReturn(0); 591 } 592 593 #undef __FUNCT__ 594 #define __FUNCT__ "CharacteristicSendCoordinatesBegin" 595 int CharacteristicSendCoordinatesBegin(Characteristic c) 596 { 597 PetscMPIInt rank, tag = 121; 598 PetscInt i, n; 599 PetscErrorCode ierr; 600 601 PetscFunctionBegin; 602 ierr = MPI_Comm_rank(((PetscObject)c)->comm, &rank);CHKERRQ(ierr); 603 ierr = HeapSort(c, c->queue, c->queueSize);CHKERRQ(ierr); 604 ierr = PetscMemzero(c->needCount, c->numNeighbors * sizeof(PetscInt));CHKERRQ(ierr); 605 for(i = 0; i < c->queueSize; i++) { 606 c->needCount[c->queue[i].proc]++; 607 } 608 c->fillCount[0] = 0; 609 for(n = 1; n < c->numNeighbors; n++) { 610 ierr = MPI_Irecv(&(c->fillCount[n]), 1, MPIU_INT, c->neighbors[n], tag, ((PetscObject)c)->comm, &(c->request[n-1]));CHKERRQ(ierr); 611 } 612 for(n = 1; n < c->numNeighbors; n++) { 613 ierr = MPI_Send(&(c->needCount[n]), 1, MPIU_INT, c->neighbors[n], tag, ((PetscObject)c)->comm);CHKERRQ(ierr); 614 } 615 ierr = MPI_Waitall(c->numNeighbors-1, c->request, c->status);CHKERRQ(ierr); 616 /* Initialize the remote queue */ 617 c->queueLocalMax = c->localOffsets[0] = 0; 618 c->queueRemoteMax = c->remoteOffsets[0] = 0; 619 for(n = 1; n < c->numNeighbors; n++) { 620 c->remoteOffsets[n] = c->queueRemoteMax; 621 c->queueRemoteMax += c->fillCount[n]; 622 c->localOffsets[n] = c->queueLocalMax; 623 c->queueLocalMax += c->needCount[n]; 624 } 625 /* HACK BEGIN */ 626 for(n = 1; n < c->numNeighbors; n++) { 627 c->localOffsets[n] += c->needCount[0]; 628 } 629 c->needCount[0] = 0; 630 /* HACK END */ 631 if (c->queueRemoteMax) { 632 ierr = PetscMalloc(sizeof(CharacteristicPointDA2D) * c->queueRemoteMax, &c->queueRemote);CHKERRQ(ierr); 633 } else { 634 c->queueRemote = PETSC_NULL; 635 } 636 c->queueRemoteSize = c->queueRemoteMax; 637 638 /* Send and Receive requests for values at t_n+1/2, giving the coordinates for interpolation */ 639 for(n = 1; n < c->numNeighbors; n++) { 640 ierr = PetscInfo2(PETSC_NULL, "Receiving %d requests for values from proc %d\n", c->fillCount[n], c->neighbors[n]);CHKERRQ(ierr); 641 ierr = MPI_Irecv(&(c->queueRemote[c->remoteOffsets[n]]), c->fillCount[n], c->itemType, c->neighbors[n], tag, ((PetscObject)c)->comm, &(c->request[n-1]));CHKERRQ(ierr); 642 } 643 for(n = 1; n < c->numNeighbors; n++) { 644 ierr = PetscInfo2(PETSC_NULL, "Sending %d requests for values from proc %d\n", c->needCount[n], c->neighbors[n]);CHKERRQ(ierr); 645 ierr = MPI_Send(&(c->queue[c->localOffsets[n]]), c->needCount[n], c->itemType, c->neighbors[n], tag, ((PetscObject)c)->comm);CHKERRQ(ierr); 646 } 647 PetscFunctionReturn(0); 648 } 649 650 #undef __FUNCT__ 651 #define __FUNCT__ "CharacteristicSendCoordinatesEnd" 652 PetscErrorCode CharacteristicSendCoordinatesEnd(Characteristic c) 653 { 654 #if 0 655 PetscMPIInt rank; 656 PetscInt n; 657 #endif 658 PetscErrorCode ierr; 659 660 PetscFunctionBegin; 661 ierr = MPI_Waitall(c->numNeighbors-1, c->request, c->status);CHKERRQ(ierr); 662 #if 0 663 ierr = MPI_Comm_rank(((PetscObject)c)->comm, &rank);CHKERRQ(ierr); 664 for(n = 0; n < c->queueRemoteSize; n++) { 665 if (c->neighbors[c->queueRemote[n].proc] == rank) { 666 SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB, "This is fucked up, n = %d proc = %d", n, c->queueRemote[n].proc); 667 } 668 } 669 #endif 670 PetscFunctionReturn(0); 671 } 672 673 #undef __FUNCT__ 674 #define __FUNCT__ "CharacteristicGetValuesBegin" 675 PetscErrorCode CharacteristicGetValuesBegin(Characteristic c) 676 { 677 PetscMPIInt tag = 121; 678 PetscInt n; 679 PetscErrorCode ierr; 680 681 PetscFunctionBegin; 682 /* SEND AND RECIEVE FILLED REQUESTS for velocities at t_n+1/2 */ 683 for(n = 1; n < c->numNeighbors; n++) { 684 ierr = MPI_Irecv(&(c->queue[c->localOffsets[n]]), c->needCount[n], c->itemType, c->neighbors[n], tag, ((PetscObject)c)->comm, &(c->request[n-1]));CHKERRQ(ierr); 685 } 686 for(n = 1; n < c->numNeighbors; n++) { 687 ierr = MPI_Send(&(c->queueRemote[c->remoteOffsets[n]]), c->fillCount[n], c->itemType, c->neighbors[n], tag, ((PetscObject)c)->comm);CHKERRQ(ierr); 688 } 689 PetscFunctionReturn(0); 690 } 691 692 #undef __FUNCT__ 693 #define __FUNCT__ "CharacteristicGetValuesEnd" 694 PetscErrorCode CharacteristicGetValuesEnd(Characteristic c) 695 { 696 PetscErrorCode ierr; 697 698 PetscFunctionBegin; 699 ierr = MPI_Waitall(c->numNeighbors-1, c->request, c->status);CHKERRQ(ierr); 700 /* Free queue of requests from other procs */ 701 ierr = PetscFree(c->queueRemote);CHKERRQ(ierr); 702 PetscFunctionReturn(0); 703 } 704 705 /*---------------------------------------------------------------------*/ 706 #undef __FUNCT__ 707 #define __FUNCT__ "HeapSort" 708 /* 709 Based on code from http://linux.wku.edu/~lamonml/algor/sort/heap.html 710 */ 711 int HeapSort(Characteristic c, Queue queue, PetscInt size) 712 /*---------------------------------------------------------------------*/ 713 { 714 CharacteristicPointDA2D temp; 715 int n; 716 717 if (0) { /* Check the order of the queue before sorting */ 718 PetscInfo(PETSC_NULL, "Before Heap sort\n"); 719 for (n=0; n<size; n++) { 720 PetscInfo2(PETSC_NULL,"%d %d\n",n,queue[n].proc); 721 } 722 } 723 724 /* SORTING PHASE */ 725 for (n = (size / 2)-1; n >= 0; n--) 726 SiftDown(c, queue, n, size-1); /* Rich had size-1 here, Matt had size*/ 727 for (n = size-1; n >= 1; n--) { 728 temp = queue[0]; 729 queue[0] = queue[n]; 730 queue[n] = temp; 731 SiftDown(c, queue, 0, n-1); 732 } 733 if (0) { /* Check the order of the queue after sorting */ 734 PetscInfo(PETSC_NULL, "Avter Heap sort\n"); 735 for (n=0; n<size; n++) { 736 PetscInfo2(PETSC_NULL,"%d %d\n",n,queue[n].proc); 737 } 738 } 739 return 0; 740 } 741 742 /*---------------------------------------------------------------------*/ 743 #undef __FUNCT__ 744 #define __FUNCT__ "SiftDown" 745 /* 746 Based on code from http://linux.wku.edu/~lamonml/algor/sort/heap.html 747 */ 748 PetscErrorCode SiftDown(Characteristic c, Queue queue, PetscInt root, PetscInt bottom) 749 /*---------------------------------------------------------------------*/ 750 { 751 PetscBool done = PETSC_FALSE; 752 PetscInt maxChild; 753 CharacteristicPointDA2D temp; 754 755 PetscFunctionBegin; 756 while ((root*2 <= bottom) && (!done)) { 757 if (root*2 == bottom) maxChild = root * 2; 758 else if (queue[root*2].proc > queue[root*2+1].proc) maxChild = root * 2; 759 else maxChild = root * 2 + 1; 760 761 if (queue[root].proc < queue[maxChild].proc) { 762 temp = queue[root]; 763 queue[root] = queue[maxChild]; 764 queue[maxChild] = temp; 765 root = maxChild; 766 } else 767 done = PETSC_TRUE; 768 } 769 PetscFunctionReturn(0); 770 } 771 772 #undef __FUNCT__ 773 #define __FUNCT__ "DMDAGetNeighborsRank" 774 /* [center, left, top-left, top, top-right, right, bottom-right, bottom, bottom-left] */ 775 PetscErrorCode DMDAGetNeighborsRank(DM da, PetscMPIInt neighbors[]) 776 { 777 DMDABoundaryType bx, by; 778 PetscBool IPeriodic = PETSC_FALSE, JPeriodic = PETSC_FALSE; 779 MPI_Comm comm; 780 PetscMPIInt rank; 781 PetscInt **procs,pi,pj,pim,pip,pjm,pjp,PI,PJ; 782 PetscErrorCode ierr; 783 784 PetscFunctionBegin; 785 ierr = PetscObjectGetComm((PetscObject) da, &comm);CHKERRQ(ierr); 786 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 787 ierr = DMDAGetInfo(da, 0, 0, 0, 0, &PI,&PJ, 0, 0, 0, &bx, &by,0, 0); 788 789 if (bx == DMDA_BOUNDARY_PERIODIC) { 790 IPeriodic = PETSC_TRUE; 791 } 792 if (by == DMDA_BOUNDARY_PERIODIC) { 793 JPeriodic = PETSC_TRUE; 794 } 795 796 neighbors[0] = rank; 797 rank = 0; 798 ierr = PetscMalloc(sizeof(PetscInt*)*PJ,&procs);CHKERRQ(ierr); 799 for (pj=0;pj<PJ;pj++) { 800 ierr = PetscMalloc(sizeof(PetscInt)*PI,&(procs[pj]));CHKERRQ(ierr); 801 for (pi=0;pi<PI;pi++) { 802 procs[pj][pi] = rank; 803 rank++; 804 } 805 } 806 807 pi = neighbors[0] % PI; 808 pj = neighbors[0] / PI; 809 pim = pi-1; if (pim<0) pim=PI-1; 810 pip = (pi+1)%PI; 811 pjm = pj-1; if (pjm<0) pjm=PJ-1; 812 pjp = (pj+1)%PJ; 813 814 neighbors[1] = procs[pj] [pim]; 815 neighbors[2] = procs[pjp][pim]; 816 neighbors[3] = procs[pjp][pi]; 817 neighbors[4] = procs[pjp][pip]; 818 neighbors[5] = procs[pj] [pip]; 819 neighbors[6] = procs[pjm][pip]; 820 neighbors[7] = procs[pjm][pi]; 821 neighbors[8] = procs[pjm][pim]; 822 823 if (!IPeriodic) { 824 if (pi==0) neighbors[1]=neighbors[2]=neighbors[8]=neighbors[0]; 825 if (pi==PI-1) neighbors[4]=neighbors[5]=neighbors[6]=neighbors[0]; 826 } 827 828 if (!JPeriodic) { 829 if (pj==0) neighbors[6]=neighbors[7]=neighbors[8]=neighbors[0]; 830 if (pj==PJ-1) neighbors[2]=neighbors[3]=neighbors[4]=neighbors[0]; 831 } 832 833 for(pj = 0; pj < PJ; pj++) { 834 ierr = PetscFree(procs[pj]);CHKERRQ(ierr); 835 } 836 ierr = PetscFree(procs);CHKERRQ(ierr); 837 PetscFunctionReturn(0); 838 } 839 840 #undef __FUNCT__ 841 #define __FUNCT__ "DMDAGetNeighborRelative" 842 /* 843 SUBDOMAIN NEIGHBORHOOD PROCESS MAP: 844 2 | 3 | 4 845 __|___|__ 846 1 | 0 | 5 847 __|___|__ 848 8 | 7 | 6 849 | | 850 */ 851 PetscInt DMDAGetNeighborRelative(DM da, PassiveReal ir, PassiveReal jr) 852 { 853 DMDALocalInfo info; 854 PassiveReal is,ie,js,je; 855 PetscErrorCode ierr; 856 857 ierr = DMDAGetLocalInfo(da, &info);CHKERRQ(ierr); 858 is = (PassiveReal) info.xs - 0.5; ie = (PassiveReal) info.xs + info.xm - 0.5; 859 js = (PassiveReal) info.ys - 0.5; je = (PassiveReal) info.ys + info.ym - 0.5; 860 861 if (ir >= is && ir <= ie) { /* center column */ 862 if (jr >= js && jr <= je) { 863 return 0; 864 } else if (jr < js) { 865 return 7; 866 } else { 867 return 3; 868 } 869 } else if (ir < is) { /* left column */ 870 if (jr >= js && jr <= je) { 871 return 1; 872 } else if (jr < js) { 873 return 8; 874 } else { 875 return 2; 876 } 877 } else { /* right column */ 878 if (jr >= js && jr <= je) { 879 return 5; 880 } else if (jr < js) { 881 return 6; 882 } else { 883 return 4; 884 } 885 } 886 } 887