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