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