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(CharacteristicList,type,&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 - routine_create - routine to create method context 245 246 Sample usage: 247 .vb 248 CharacteristicRegister("my_char", MyCharCreate); 249 .ve 250 251 Then, your Characteristic type can be chosen with the procedural interface via 252 .vb 253 CharacteristicCreate(MPI_Comm, Characteristic* &char); 254 CharacteristicSetType(char,"my_char"); 255 .ve 256 or at runtime via the option 257 .vb 258 -characteristic_type my_char 259 .ve 260 261 Notes: 262 CharacteristicRegister() may be called multiple times to add several user-defined solvers. 263 264 .keywords: Characteristic, register 265 266 .seealso: CharacteristicRegisterAll(), CharacteristicRegisterDestroy() 267 268 Level: advanced 269 @*/ 270 PetscErrorCode CharacteristicRegister(const char sname[],PetscErrorCode (*function)(Characteristic)) 271 { 272 PetscErrorCode ierr; 273 274 PetscFunctionBegin; 275 ierr = PetscFunctionListAdd(&CharacteristicList,sname,function);CHKERRQ(ierr); 276 PetscFunctionReturn(0); 277 } 278 279 #undef __FUNCT__ 280 #define __FUNCT__ "CharacteristicSetVelocityInterpolation" 281 PetscErrorCode CharacteristicSetVelocityInterpolation(Characteristic c, DM da, Vec v, Vec vOld, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(Vec, PetscReal[], PetscInt, PetscInt[], PetscScalar[], void*), void *ctx) 282 { 283 PetscFunctionBegin; 284 c->velocityDA = da; 285 c->velocity = v; 286 c->velocityOld = vOld; 287 c->numVelocityComp = numComponents; 288 c->velocityComp = components; 289 c->velocityInterp = interp; 290 c->velocityCtx = ctx; 291 PetscFunctionReturn(0); 292 } 293 294 #undef __FUNCT__ 295 #define __FUNCT__ "CharacteristicSetVelocityInterpolationLocal" 296 PetscErrorCode CharacteristicSetVelocityInterpolationLocal(Characteristic c, DM da, Vec v, Vec vOld, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(void*, PetscReal [], PetscInt, PetscInt[], PetscScalar[], void*), void *ctx) 297 { 298 PetscFunctionBegin; 299 c->velocityDA = da; 300 c->velocity = v; 301 c->velocityOld = vOld; 302 c->numVelocityComp = numComponents; 303 c->velocityComp = components; 304 c->velocityInterpLocal = interp; 305 c->velocityCtx = ctx; 306 PetscFunctionReturn(0); 307 } 308 309 #undef __FUNCT__ 310 #define __FUNCT__ "CharacteristicSetFieldInterpolation" 311 PetscErrorCode CharacteristicSetFieldInterpolation(Characteristic c, DM da, Vec v, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(Vec, PetscReal[], PetscInt, PetscInt[], PetscScalar[], void*), void *ctx) 312 { 313 PetscFunctionBegin; 314 #if 0 315 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."); 316 #endif 317 c->fieldDA = da; 318 c->field = v; 319 c->numFieldComp = numComponents; 320 c->fieldComp = components; 321 c->fieldInterp = interp; 322 c->fieldCtx = ctx; 323 PetscFunctionReturn(0); 324 } 325 326 #undef __FUNCT__ 327 #define __FUNCT__ "CharacteristicSetFieldInterpolationLocal" 328 PetscErrorCode CharacteristicSetFieldInterpolationLocal(Characteristic c, DM da, Vec v, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(void*, PetscReal[], PetscInt, PetscInt[], PetscScalar [], void*), void *ctx) 329 { 330 PetscFunctionBegin; 331 #if 0 332 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."); 333 #endif 334 c->fieldDA = da; 335 c->field = v; 336 c->numFieldComp = numComponents; 337 c->fieldComp = components; 338 c->fieldInterpLocal = interp; 339 c->fieldCtx = ctx; 340 PetscFunctionReturn(0); 341 } 342 343 #undef __FUNCT__ 344 #define __FUNCT__ "CharacteristicSolve" 345 PetscErrorCode CharacteristicSolve(Characteristic c, PetscReal dt, Vec solution) 346 { 347 CharacteristicPointDA2D Qi; 348 DM da = c->velocityDA; 349 Vec velocityLocal, velocityLocalOld; 350 Vec fieldLocal; 351 DMDALocalInfo info; 352 PetscScalar **solArray; 353 void *velocityArray; 354 void *velocityArrayOld; 355 void *fieldArray; 356 PassiveScalar *interpIndices; 357 PassiveScalar *velocityValues, *velocityValuesOld; 358 PassiveScalar *fieldValues; 359 PetscMPIInt rank; 360 PetscInt dim; 361 PetscMPIInt neighbors[9]; 362 PetscInt dof; 363 PetscInt gx, gy; 364 PetscInt n, is, ie, js, je, comp; 365 PetscErrorCode ierr; 366 367 PetscFunctionBegin; 368 c->queueSize = 0; 369 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank);CHKERRQ(ierr); 370 ierr = DMDAGetNeighborsRank(da, neighbors);CHKERRQ(ierr); 371 ierr = CharacteristicSetNeighbors(c, 9, neighbors);CHKERRQ(ierr); 372 ierr = CharacteristicSetUp(c);CHKERRQ(ierr); 373 /* global and local grid info */ 374 ierr = DMDAGetInfo(da, &dim, &gx, &gy, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0);CHKERRQ(ierr); 375 ierr = DMDAGetLocalInfo(da, &info);CHKERRQ(ierr); 376 is = info.xs; ie = info.xs+info.xm; 377 js = info.ys; je = info.ys+info.ym; 378 /* Allocation */ 379 ierr = PetscMalloc(dim*sizeof(PetscScalar), &interpIndices);CHKERRQ(ierr); 380 ierr = PetscMalloc(c->numVelocityComp*sizeof(PetscScalar), &velocityValues);CHKERRQ(ierr); 381 ierr = PetscMalloc(c->numVelocityComp*sizeof(PetscScalar), &velocityValuesOld);CHKERRQ(ierr); 382 ierr = PetscMalloc(c->numFieldComp*sizeof(PetscScalar), &fieldValues);CHKERRQ(ierr); 383 ierr = PetscLogEventBegin(CHARACTERISTIC_Solve,0,0,0,0);CHKERRQ(ierr); 384 385 /* ----------------------------------------------------------------------- 386 PART 1, AT t-dt/2 387 -----------------------------------------------------------------------*/ 388 ierr = PetscLogEventBegin(CHARACTERISTIC_QueueSetup,0,0,0,0);CHKERRQ(ierr); 389 /* GET POSITION AT HALF TIME IN THE PAST */ 390 if (c->velocityInterpLocal) { 391 ierr = DMGetLocalVector(c->velocityDA, &velocityLocal);CHKERRQ(ierr); 392 ierr = DMGetLocalVector(c->velocityDA, &velocityLocalOld);CHKERRQ(ierr); 393 ierr = DMGlobalToLocalBegin(c->velocityDA, c->velocity, INSERT_VALUES, velocityLocal);CHKERRQ(ierr); 394 ierr = DMGlobalToLocalEnd(c->velocityDA, c->velocity, INSERT_VALUES, velocityLocal);CHKERRQ(ierr); 395 ierr = DMGlobalToLocalBegin(c->velocityDA, c->velocityOld, INSERT_VALUES, velocityLocalOld);CHKERRQ(ierr); 396 ierr = DMGlobalToLocalEnd(c->velocityDA, c->velocityOld, INSERT_VALUES, velocityLocalOld);CHKERRQ(ierr); 397 ierr = DMDAVecGetArray(c->velocityDA, velocityLocal, &velocityArray);CHKERRQ(ierr); 398 ierr = DMDAVecGetArray(c->velocityDA, velocityLocalOld, &velocityArrayOld);CHKERRQ(ierr); 399 } 400 ierr = PetscInfo(NULL, "Calculating position at t_{n - 1/2}\n");CHKERRQ(ierr); 401 for (Qi.j = js; Qi.j < je; Qi.j++) { 402 for (Qi.i = is; Qi.i < ie; Qi.i++) { 403 interpIndices[0] = Qi.i; 404 interpIndices[1] = Qi.j; 405 if (c->velocityInterpLocal) c->velocityInterpLocal(velocityArray, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx); 406 else c->velocityInterp(c->velocity, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx); 407 Qi.x = Qi.i - velocityValues[0]*dt/2.0; 408 Qi.y = Qi.j - velocityValues[1]*dt/2.0; 409 410 /* Determine whether the position at t - dt/2 is local */ 411 Qi.proc = DMDAGetNeighborRelative(da, Qi.x, Qi.y); 412 413 /* Check for Periodic boundaries and move all periodic points back onto the domain */ 414 ierr = DMDAMapCoordsToPeriodicDomain(da,&(Qi.x),&(Qi.y));CHKERRQ(ierr); 415 ierr = CharacteristicAddPoint(c, &Qi);CHKERRQ(ierr); 416 } 417 } 418 ierr = PetscLogEventEnd(CHARACTERISTIC_QueueSetup,0,0,0,0);CHKERRQ(ierr); 419 420 ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr); 421 ierr = CharacteristicSendCoordinatesBegin(c);CHKERRQ(ierr); 422 ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr); 423 424 ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeLocal,0,0,0,0);CHKERRQ(ierr); 425 /* Calculate velocity at t_n+1/2 (local values) */ 426 ierr = PetscInfo(NULL, "Calculating local velocities at t_{n - 1/2}\n");CHKERRQ(ierr); 427 for (n = 0; n < c->queueSize; n++) { 428 Qi = c->queue[n]; 429 if (c->neighbors[Qi.proc] == rank) { 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 } 442 c->queue[n] = Qi; 443 } 444 ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeLocal,0,0,0,0);CHKERRQ(ierr); 445 446 ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr); 447 ierr = CharacteristicSendCoordinatesEnd(c);CHKERRQ(ierr); 448 ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr); 449 450 451 /* Calculate velocity at t_n+1/2 (fill remote requests) */ 452 ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeRemote,0,0,0,0);CHKERRQ(ierr); 453 ierr = PetscInfo1(NULL, "Calculating %d remote velocities at t_{n - 1/2}\n", c->queueRemoteSize);CHKERRQ(ierr); 454 for (n = 0; n < c->queueRemoteSize; n++) { 455 Qi = c->queueRemote[n]; 456 interpIndices[0] = Qi.x; 457 interpIndices[1] = Qi.y; 458 if (c->velocityInterpLocal) { 459 c->velocityInterpLocal(velocityArray, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx); 460 c->velocityInterpLocal(velocityArrayOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx); 461 } else { 462 c->velocityInterp(c->velocity, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx); 463 c->velocityInterp(c->velocityOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx); 464 } 465 Qi.x = 0.5*(velocityValues[0] + velocityValuesOld[0]); 466 Qi.y = 0.5*(velocityValues[1] + velocityValuesOld[1]); 467 c->queueRemote[n] = Qi; 468 } 469 ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeRemote,0,0,0,0);CHKERRQ(ierr); 470 ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr); 471 ierr = CharacteristicGetValuesBegin(c);CHKERRQ(ierr); 472 ierr = CharacteristicGetValuesEnd(c);CHKERRQ(ierr); 473 if (c->velocityInterpLocal) { 474 ierr = DMDAVecRestoreArray(c->velocityDA, velocityLocal, &velocityArray);CHKERRQ(ierr); 475 ierr = DMDAVecRestoreArray(c->velocityDA, velocityLocalOld, &velocityArrayOld);CHKERRQ(ierr); 476 ierr = DMRestoreLocalVector(c->velocityDA, &velocityLocal);CHKERRQ(ierr); 477 ierr = DMRestoreLocalVector(c->velocityDA, &velocityLocalOld);CHKERRQ(ierr); 478 } 479 ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr); 480 481 /* ----------------------------------------------------------------------- 482 PART 2, AT t-dt 483 -----------------------------------------------------------------------*/ 484 485 /* GET POSITION AT t_n (local values) */ 486 ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeLocal,0,0,0,0);CHKERRQ(ierr); 487 ierr = PetscInfo(NULL, "Calculating position at t_{n}\n");CHKERRQ(ierr); 488 for (n = 0; n < c->queueSize; n++) { 489 Qi = c->queue[n]; 490 Qi.x = Qi.i - Qi.x*dt; 491 Qi.y = Qi.j - Qi.y*dt; 492 493 /* Determine whether the position at t-dt is local */ 494 Qi.proc = DMDAGetNeighborRelative(da, Qi.x, Qi.y); 495 496 /* Check for Periodic boundaries and move all periodic points back onto the domain */ 497 ierr = DMDAMapCoordsToPeriodicDomain(da,&(Qi.x),&(Qi.y));CHKERRQ(ierr); 498 499 c->queue[n] = Qi; 500 } 501 ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeLocal,0,0,0,0);CHKERRQ(ierr); 502 503 ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr); 504 ierr = CharacteristicSendCoordinatesBegin(c);CHKERRQ(ierr); 505 ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr); 506 507 /* GET VALUE AT FULL TIME IN THE PAST (LOCAL REQUESTS) */ 508 ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeLocal,0,0,0,0);CHKERRQ(ierr); 509 if (c->fieldInterpLocal) { 510 ierr = DMGetLocalVector(c->fieldDA, &fieldLocal);CHKERRQ(ierr); 511 ierr = DMGlobalToLocalBegin(c->fieldDA, c->field, INSERT_VALUES, fieldLocal);CHKERRQ(ierr); 512 ierr = DMGlobalToLocalEnd(c->fieldDA, c->field, INSERT_VALUES, fieldLocal);CHKERRQ(ierr); 513 ierr = DMDAVecGetArray(c->fieldDA, fieldLocal, &fieldArray);CHKERRQ(ierr); 514 } 515 ierr = PetscInfo(NULL, "Calculating local field at t_{n}\n");CHKERRQ(ierr); 516 for (n = 0; n < c->queueSize; n++) { 517 if (c->neighbors[c->queue[n].proc] == rank) { 518 interpIndices[0] = c->queue[n].x; 519 interpIndices[1] = c->queue[n].y; 520 if (c->fieldInterpLocal) c->fieldInterpLocal(fieldArray, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx); 521 else c->fieldInterp(c->field, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx); 522 for (comp = 0; comp < c->numFieldComp; comp++) c->queue[n].field[comp] = fieldValues[comp]; 523 } 524 } 525 ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeLocal,0,0,0,0);CHKERRQ(ierr); 526 527 ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr); 528 ierr = CharacteristicSendCoordinatesEnd(c);CHKERRQ(ierr); 529 ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr); 530 531 /* GET VALUE AT FULL TIME IN THE PAST (REMOTE REQUESTS) */ 532 ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeRemote,0,0,0,0);CHKERRQ(ierr); 533 ierr = PetscInfo1(NULL, "Calculating %d remote field points at t_{n}\n", c->queueRemoteSize);CHKERRQ(ierr); 534 for (n = 0; n < c->queueRemoteSize; n++) { 535 interpIndices[0] = c->queueRemote[n].x; 536 interpIndices[1] = c->queueRemote[n].y; 537 538 /* for debugging purposes */ 539 if (1) { /* hacked bounds test...let's do better */ 540 PetscScalar im = interpIndices[0]; PetscScalar jm = interpIndices[1]; 541 542 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); 543 } 544 545 if (c->fieldInterpLocal) c->fieldInterpLocal(fieldArray, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx); 546 else c->fieldInterp(c->field, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx); 547 for (comp = 0; comp < c->numFieldComp; comp++) c->queueRemote[n].field[comp] = fieldValues[comp]; 548 } 549 ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeRemote,0,0,0,0);CHKERRQ(ierr); 550 551 ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr); 552 ierr = CharacteristicGetValuesBegin(c);CHKERRQ(ierr); 553 ierr = CharacteristicGetValuesEnd(c);CHKERRQ(ierr); 554 if (c->fieldInterpLocal) { 555 ierr = DMDAVecRestoreArray(c->fieldDA, fieldLocal, &fieldArray);CHKERRQ(ierr); 556 ierr = DMRestoreLocalVector(c->fieldDA, &fieldLocal);CHKERRQ(ierr); 557 } 558 ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr); 559 560 /* Return field of characteristics at t_n-1 */ 561 ierr = PetscLogEventBegin(CHARACTERISTIC_DAUpdate,0,0,0,0);CHKERRQ(ierr); 562 ierr = DMDAGetInfo(c->fieldDA,0,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); 563 ierr = DMDAVecGetArray(c->fieldDA, solution, &solArray);CHKERRQ(ierr); 564 for (n = 0; n < c->queueSize; n++) { 565 Qi = c->queue[n]; 566 for (comp = 0; comp < c->numFieldComp; comp++) solArray[Qi.j][Qi.i*dof+c->fieldComp[comp]] = Qi.field[comp]; 567 } 568 ierr = DMDAVecRestoreArray(c->fieldDA, solution, &solArray);CHKERRQ(ierr); 569 ierr = PetscLogEventEnd(CHARACTERISTIC_DAUpdate,0,0,0,0);CHKERRQ(ierr); 570 ierr = PetscLogEventEnd(CHARACTERISTIC_Solve,0,0,0,0);CHKERRQ(ierr); 571 572 /* Cleanup */ 573 ierr = PetscFree(interpIndices);CHKERRQ(ierr); 574 ierr = PetscFree(velocityValues);CHKERRQ(ierr); 575 ierr = PetscFree(velocityValuesOld);CHKERRQ(ierr); 576 ierr = PetscFree(fieldValues);CHKERRQ(ierr); 577 PetscFunctionReturn(0); 578 } 579 580 #undef __FUNCT__ 581 #define __FUNCT__ "CharacteristicSetNeighbors" 582 PetscErrorCode CharacteristicSetNeighbors(Characteristic c, PetscInt numNeighbors, PetscMPIInt neighbors[]) 583 { 584 PetscErrorCode ierr; 585 586 PetscFunctionBegin; 587 c->numNeighbors = numNeighbors; 588 ierr = PetscFree(c->neighbors);CHKERRQ(ierr); 589 ierr = PetscMalloc(numNeighbors * sizeof(PetscMPIInt), &c->neighbors);CHKERRQ(ierr); 590 ierr = PetscMemcpy(c->neighbors, neighbors, numNeighbors * sizeof(PetscMPIInt));CHKERRQ(ierr); 591 PetscFunctionReturn(0); 592 } 593 594 #undef __FUNCT__ 595 #define __FUNCT__ "CharacteristicAddPoint" 596 PetscErrorCode CharacteristicAddPoint(Characteristic c, CharacteristicPointDA2D *point) 597 { 598 PetscFunctionBegin; 599 if (c->queueSize >= c->queueMax) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Exceeeded maximum queue size %d", c->queueMax); 600 c->queue[c->queueSize++] = *point; 601 PetscFunctionReturn(0); 602 } 603 604 #undef __FUNCT__ 605 #define __FUNCT__ "CharacteristicSendCoordinatesBegin" 606 int CharacteristicSendCoordinatesBegin(Characteristic c) 607 { 608 PetscMPIInt rank, tag = 121; 609 PetscInt i, n; 610 PetscErrorCode ierr; 611 612 PetscFunctionBegin; 613 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank);CHKERRQ(ierr); 614 ierr = CharacteristicHeapSort(c, c->queue, c->queueSize);CHKERRQ(ierr); 615 ierr = PetscMemzero(c->needCount, c->numNeighbors * sizeof(PetscInt));CHKERRQ(ierr); 616 for (i = 0; i < c->queueSize; i++) c->needCount[c->queue[i].proc]++; 617 c->fillCount[0] = 0; 618 for (n = 1; n < c->numNeighbors; n++) { 619 ierr = MPI_Irecv(&(c->fillCount[n]), 1, MPIU_INT, c->neighbors[n], tag, PetscObjectComm((PetscObject)c), &(c->request[n-1]));CHKERRQ(ierr); 620 } 621 for (n = 1; n < c->numNeighbors; n++) { 622 ierr = MPI_Send(&(c->needCount[n]), 1, MPIU_INT, c->neighbors[n], tag, PetscObjectComm((PetscObject)c));CHKERRQ(ierr); 623 } 624 ierr = MPI_Waitall(c->numNeighbors-1, c->request, c->status);CHKERRQ(ierr); 625 /* Initialize the remote queue */ 626 c->queueLocalMax = c->localOffsets[0] = 0; 627 c->queueRemoteMax = c->remoteOffsets[0] = 0; 628 for (n = 1; n < c->numNeighbors; n++) { 629 c->remoteOffsets[n] = c->queueRemoteMax; 630 c->queueRemoteMax += c->fillCount[n]; 631 c->localOffsets[n] = c->queueLocalMax; 632 c->queueLocalMax += c->needCount[n]; 633 } 634 /* HACK BEGIN */ 635 for (n = 1; n < c->numNeighbors; n++) c->localOffsets[n] += c->needCount[0]; 636 c->needCount[0] = 0; 637 /* HACK END */ 638 if (c->queueRemoteMax) { 639 ierr = PetscMalloc(sizeof(CharacteristicPointDA2D) * c->queueRemoteMax, &c->queueRemote);CHKERRQ(ierr); 640 } else c->queueRemote = NULL; 641 c->queueRemoteSize = c->queueRemoteMax; 642 643 /* Send and Receive requests for values at t_n+1/2, giving the coordinates for interpolation */ 644 for (n = 1; n < c->numNeighbors; n++) { 645 ierr = PetscInfo2(NULL, "Receiving %d requests for values from proc %d\n", c->fillCount[n], c->neighbors[n]);CHKERRQ(ierr); 646 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); 647 } 648 for (n = 1; n < c->numNeighbors; n++) { 649 ierr = PetscInfo2(NULL, "Sending %d requests for values from proc %d\n", c->needCount[n], c->neighbors[n]);CHKERRQ(ierr); 650 ierr = MPI_Send(&(c->queue[c->localOffsets[n]]), c->needCount[n], c->itemType, c->neighbors[n], tag, PetscObjectComm((PetscObject)c));CHKERRQ(ierr); 651 } 652 PetscFunctionReturn(0); 653 } 654 655 #undef __FUNCT__ 656 #define __FUNCT__ "CharacteristicSendCoordinatesEnd" 657 PetscErrorCode CharacteristicSendCoordinatesEnd(Characteristic c) 658 { 659 #if 0 660 PetscMPIInt rank; 661 PetscInt n; 662 #endif 663 PetscErrorCode ierr; 664 665 PetscFunctionBegin; 666 ierr = MPI_Waitall(c->numNeighbors-1, c->request, c->status);CHKERRQ(ierr); 667 #if 0 668 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank);CHKERRQ(ierr); 669 for (n = 0; n < c->queueRemoteSize; n++) { 670 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); 671 } 672 #endif 673 PetscFunctionReturn(0); 674 } 675 676 #undef __FUNCT__ 677 #define __FUNCT__ "CharacteristicGetValuesBegin" 678 PetscErrorCode CharacteristicGetValuesBegin(Characteristic c) 679 { 680 PetscMPIInt tag = 121; 681 PetscInt n; 682 PetscErrorCode ierr; 683 684 PetscFunctionBegin; 685 /* SEND AND RECIEVE FILLED REQUESTS for velocities at t_n+1/2 */ 686 for (n = 1; n < c->numNeighbors; n++) { 687 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); 688 } 689 for (n = 1; n < c->numNeighbors; n++) { 690 ierr = MPI_Send(&(c->queueRemote[c->remoteOffsets[n]]), c->fillCount[n], c->itemType, c->neighbors[n], tag, PetscObjectComm((PetscObject)c));CHKERRQ(ierr); 691 } 692 PetscFunctionReturn(0); 693 } 694 695 #undef __FUNCT__ 696 #define __FUNCT__ "CharacteristicGetValuesEnd" 697 PetscErrorCode CharacteristicGetValuesEnd(Characteristic c) 698 { 699 PetscErrorCode ierr; 700 701 PetscFunctionBegin; 702 ierr = MPI_Waitall(c->numNeighbors-1, c->request, c->status);CHKERRQ(ierr); 703 /* Free queue of requests from other procs */ 704 ierr = PetscFree(c->queueRemote);CHKERRQ(ierr); 705 PetscFunctionReturn(0); 706 } 707 708 /*---------------------------------------------------------------------*/ 709 #undef __FUNCT__ 710 #define __FUNCT__ "CharacteristicHeapSort" 711 /* 712 Based on code from http://linux.wku.edu/~lamonml/algor/sort/heap.html 713 */ 714 PetscErrorCode CharacteristicHeapSort(Characteristic c, Queue queue, PetscInt size) 715 /*---------------------------------------------------------------------*/ 716 { 717 PetscErrorCode ierr; 718 CharacteristicPointDA2D temp; 719 PetscInt n; 720 721 PetscFunctionBegin; 722 if (0) { /* Check the order of the queue before sorting */ 723 PetscInfo(NULL, "Before Heap sort\n"); 724 for (n=0; n<size; n++) { 725 ierr = PetscInfo2(NULL,"%d %d\n",n,queue[n].proc);CHKERRQ(ierr); 726 } 727 } 728 729 /* SORTING PHASE */ 730 for (n = (size / 2)-1; n >= 0; n--) { 731 ierr = CharacteristicSiftDown(c, queue, n, size-1);CHKERRQ(ierr); /* Rich had size-1 here, Matt had size*/ 732 } 733 for (n = size-1; n >= 1; n--) { 734 temp = queue[0]; 735 queue[0] = queue[n]; 736 queue[n] = temp; 737 ierr = CharacteristicSiftDown(c, queue, 0, n-1);CHKERRQ(ierr); 738 } 739 if (0) { /* Check the order of the queue after sorting */ 740 ierr = PetscInfo(NULL, "Avter Heap sort\n");CHKERRQ(ierr); 741 for (n=0; n<size; n++) { 742 ierr = PetscInfo2(NULL,"%d %d\n",n,queue[n].proc);CHKERRQ(ierr); 743 } 744 } 745 PetscFunctionReturn(0); 746 } 747 748 /*---------------------------------------------------------------------*/ 749 #undef __FUNCT__ 750 #define __FUNCT__ "CharacteristicSiftDown" 751 /* 752 Based on code from http://linux.wku.edu/~lamonml/algor/sort/heap.html 753 */ 754 PetscErrorCode CharacteristicSiftDown(Characteristic c, Queue queue, PetscInt root, PetscInt bottom) 755 /*---------------------------------------------------------------------*/ 756 { 757 PetscBool done = PETSC_FALSE; 758 PetscInt maxChild; 759 CharacteristicPointDA2D temp; 760 761 PetscFunctionBegin; 762 while ((root*2 <= bottom) && (!done)) { 763 if (root*2 == bottom) maxChild = root * 2; 764 else if (queue[root*2].proc > queue[root*2+1].proc) maxChild = root * 2; 765 else maxChild = root * 2 + 1; 766 767 if (queue[root].proc < queue[maxChild].proc) { 768 temp = queue[root]; 769 queue[root] = queue[maxChild]; 770 queue[maxChild] = temp; 771 root = maxChild; 772 } else done = PETSC_TRUE; 773 } 774 PetscFunctionReturn(0); 775 } 776 777 #undef __FUNCT__ 778 #define __FUNCT__ "DMDAGetNeighborsRank" 779 /* [center, left, top-left, top, top-right, right, bottom-right, bottom, bottom-left] */ 780 PetscErrorCode DMDAGetNeighborsRank(DM da, PetscMPIInt neighbors[]) 781 { 782 DMDABoundaryType bx, by; 783 PetscBool IPeriodic = PETSC_FALSE, JPeriodic = PETSC_FALSE; 784 MPI_Comm comm; 785 PetscMPIInt rank; 786 PetscInt **procs,pi,pj,pim,pip,pjm,pjp,PI,PJ; 787 PetscErrorCode ierr; 788 789 PetscFunctionBegin; 790 ierr = PetscObjectGetComm((PetscObject) da, &comm);CHKERRQ(ierr); 791 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 792 ierr = DMDAGetInfo(da, 0, 0, 0, 0, &PI,&PJ, 0, 0, 0, &bx, &by,0, 0);CHKERRQ(ierr); 793 794 if (bx == DMDA_BOUNDARY_PERIODIC) IPeriodic = PETSC_TRUE; 795 if (by == DMDA_BOUNDARY_PERIODIC) JPeriodic = PETSC_TRUE; 796 797 neighbors[0] = rank; 798 rank = 0; 799 ierr = PetscMalloc(sizeof(PetscInt*)*PJ,&procs);CHKERRQ(ierr); 800 for (pj=0; pj<PJ; pj++) { 801 ierr = PetscMalloc(sizeof(PetscInt)*PI,&(procs[pj]));CHKERRQ(ierr); 802 for (pi=0; pi<PI; pi++) { 803 procs[pj][pi] = rank; 804 rank++; 805 } 806 } 807 808 pi = neighbors[0] % PI; 809 pj = neighbors[0] / PI; 810 pim = pi-1; if (pim<0) pim=PI-1; 811 pip = (pi+1)%PI; 812 pjm = pj-1; if (pjm<0) pjm=PJ-1; 813 pjp = (pj+1)%PJ; 814 815 neighbors[1] = procs[pj] [pim]; 816 neighbors[2] = procs[pjp][pim]; 817 neighbors[3] = procs[pjp][pi]; 818 neighbors[4] = procs[pjp][pip]; 819 neighbors[5] = procs[pj] [pip]; 820 neighbors[6] = procs[pjm][pip]; 821 neighbors[7] = procs[pjm][pi]; 822 neighbors[8] = procs[pjm][pim]; 823 824 if (!IPeriodic) { 825 if (pi==0) neighbors[1]=neighbors[2]=neighbors[8]=neighbors[0]; 826 if (pi==PI-1) neighbors[4]=neighbors[5]=neighbors[6]=neighbors[0]; 827 } 828 829 if (!JPeriodic) { 830 if (pj==0) neighbors[6]=neighbors[7]=neighbors[8]=neighbors[0]; 831 if (pj==PJ-1) neighbors[2]=neighbors[3]=neighbors[4]=neighbors[0]; 832 } 833 834 for (pj = 0; pj < PJ; pj++) { 835 ierr = PetscFree(procs[pj]);CHKERRQ(ierr); 836 } 837 ierr = PetscFree(procs);CHKERRQ(ierr); 838 PetscFunctionReturn(0); 839 } 840 841 #undef __FUNCT__ 842 #define __FUNCT__ "DMDAGetNeighborRelative" 843 /* 844 SUBDOMAIN NEIGHBORHOOD PROCESS MAP: 845 2 | 3 | 4 846 __|___|__ 847 1 | 0 | 5 848 __|___|__ 849 8 | 7 | 6 850 | | 851 */ 852 PetscInt DMDAGetNeighborRelative(DM da, PassiveReal ir, PassiveReal jr) 853 { 854 DMDALocalInfo info; 855 PassiveReal is,ie,js,je; 856 PetscErrorCode ierr; 857 858 ierr = DMDAGetLocalInfo(da, &info);CHKERRQ(ierr); 859 is = (PassiveReal) info.xs - 0.5; ie = (PassiveReal) info.xs + info.xm - 0.5; 860 js = (PassiveReal) info.ys - 0.5; je = (PassiveReal) info.ys + info.ym - 0.5; 861 862 if (ir >= is && ir <= ie) { /* center column */ 863 if (jr >= js && jr <= je) return 0; 864 else if (jr < js) return 7; 865 else return 3; 866 } else if (ir < is) { /* left column */ 867 if (jr >= js && jr <= je) return 1; 868 else if (jr < js) return 8; 869 else return 2; 870 } else { /* right column */ 871 if (jr >= js && jr <= je) return 5; 872 else if (jr < js) return 6; 873 else return 4; 874 } 875 } 876