1af0996ceSBarry Smith #include <petsc/private/characteristicimpl.h> /*I "petsccharacteristic.h" I*/ 21e25c274SJed Brown #include <petscdmda.h> 3665c2dedSJed Brown #include <petscviewer.h> 4af33a6ddSJed Brown 5af33a6ddSJed Brown PetscClassId CHARACTERISTIC_CLASSID; 6af33a6ddSJed Brown PetscLogEvent CHARACTERISTIC_SetUp, CHARACTERISTIC_Solve, CHARACTERISTIC_QueueSetup, CHARACTERISTIC_DAUpdate; 7af33a6ddSJed Brown PetscLogEvent CHARACTERISTIC_HalfTimeLocal, CHARACTERISTIC_HalfTimeRemote, CHARACTERISTIC_HalfTimeExchange; 8af33a6ddSJed Brown PetscLogEvent CHARACTERISTIC_FullTimeLocal, CHARACTERISTIC_FullTimeRemote, CHARACTERISTIC_FullTimeExchange; 9af33a6ddSJed Brown /* 10af33a6ddSJed Brown Contains the list of registered characteristic routines 11af33a6ddSJed Brown */ 120298fd71SBarry Smith PetscFunctionList CharacteristicList = NULL; 13140e18c1SBarry Smith PetscBool CharacteristicRegisterAllCalled = PETSC_FALSE; 14af33a6ddSJed Brown 1566976f2fSJacob Faibussowitsch static PetscErrorCode DMDAGetNeighborsRank(DM, PetscMPIInt[]); 1666976f2fSJacob Faibussowitsch static PetscInt DMDAGetNeighborRelative(DM, PetscReal, PetscReal); 17af33a6ddSJed Brown 1866976f2fSJacob Faibussowitsch static PetscErrorCode CharacteristicHeapSort(Characteristic, Queue, PetscInt); 1966976f2fSJacob Faibussowitsch static PetscErrorCode CharacteristicSiftDown(Characteristic, Queue, PetscInt, PetscInt); 20af33a6ddSJed Brown 2166976f2fSJacob Faibussowitsch static PetscErrorCode CharacteristicView(Characteristic c, PetscViewer viewer) 22d71ae5a4SJacob Faibussowitsch { 23af33a6ddSJed Brown PetscBool iascii; 24af33a6ddSJed Brown 25af33a6ddSJed Brown PetscFunctionBegin; 26af33a6ddSJed Brown PetscValidHeaderSpecific(c, CHARACTERISTIC_CLASSID, 1); 2748a46eb9SPierre Jolivet if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)c), &viewer)); 28af33a6ddSJed Brown PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 29af33a6ddSJed Brown PetscCheckSameComm(c, 1, viewer, 2); 30af33a6ddSJed Brown 319566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 32ad540459SPierre Jolivet if (!iascii) PetscTryTypeMethod(c, view, viewer); 333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 34af33a6ddSJed Brown } 35af33a6ddSJed Brown 36d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicDestroy(Characteristic *c) 37d71ae5a4SJacob Faibussowitsch { 38af33a6ddSJed Brown PetscFunctionBegin; 393ba16761SJacob Faibussowitsch if (!*c) PetscFunctionReturn(PETSC_SUCCESS); 406bf464f9SBarry Smith PetscValidHeaderSpecific(*c, CHARACTERISTIC_CLASSID, 1); 41*f4f49eeaSPierre Jolivet if (--((PetscObject)*c)->refct > 0) PetscFunctionReturn(PETSC_SUCCESS); 42af33a6ddSJed Brown 43213acdd3SPierre Jolivet PetscTryTypeMethod(*c, destroy); 449566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&(*c)->itemType)); 459566063dSJacob Faibussowitsch PetscCall(PetscFree((*c)->queue)); 469566063dSJacob Faibussowitsch PetscCall(PetscFree((*c)->queueLocal)); 479566063dSJacob Faibussowitsch PetscCall(PetscFree((*c)->queueRemote)); 489566063dSJacob Faibussowitsch PetscCall(PetscFree((*c)->neighbors)); 499566063dSJacob Faibussowitsch PetscCall(PetscFree((*c)->needCount)); 509566063dSJacob Faibussowitsch PetscCall(PetscFree((*c)->localOffsets)); 519566063dSJacob Faibussowitsch PetscCall(PetscFree((*c)->fillCount)); 529566063dSJacob Faibussowitsch PetscCall(PetscFree((*c)->remoteOffsets)); 539566063dSJacob Faibussowitsch PetscCall(PetscFree((*c)->request)); 549566063dSJacob Faibussowitsch PetscCall(PetscFree((*c)->status)); 559566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(c)); 563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 57af33a6ddSJed Brown } 58af33a6ddSJed Brown 59d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicCreate(MPI_Comm comm, Characteristic *c) 60d71ae5a4SJacob Faibussowitsch { 61af33a6ddSJed Brown Characteristic newC; 62af33a6ddSJed Brown 63af33a6ddSJed Brown PetscFunctionBegin; 644f572ea9SToby Isaac PetscAssertPointer(c, 2); 650298fd71SBarry Smith *c = NULL; 669566063dSJacob Faibussowitsch PetscCall(CharacteristicInitializePackage()); 67af33a6ddSJed Brown 689566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(newC, CHARACTERISTIC_CLASSID, "Characteristic", "Characteristic", "Characteristic", comm, CharacteristicDestroy, CharacteristicView)); 69af33a6ddSJed Brown *c = newC; 70af33a6ddSJed Brown 71af33a6ddSJed Brown newC->structured = PETSC_TRUE; 72af33a6ddSJed Brown newC->numIds = 0; 730298fd71SBarry Smith newC->velocityDA = NULL; 740298fd71SBarry Smith newC->velocity = NULL; 750298fd71SBarry Smith newC->velocityOld = NULL; 76af33a6ddSJed Brown newC->numVelocityComp = 0; 770298fd71SBarry Smith newC->velocityComp = NULL; 780298fd71SBarry Smith newC->velocityInterp = NULL; 790298fd71SBarry Smith newC->velocityInterpLocal = NULL; 800298fd71SBarry Smith newC->velocityCtx = NULL; 810298fd71SBarry Smith newC->fieldDA = NULL; 820298fd71SBarry Smith newC->field = NULL; 83af33a6ddSJed Brown newC->numFieldComp = 0; 840298fd71SBarry Smith newC->fieldComp = NULL; 850298fd71SBarry Smith newC->fieldInterp = NULL; 860298fd71SBarry Smith newC->fieldInterpLocal = NULL; 870298fd71SBarry Smith newC->fieldCtx = NULL; 880298fd71SBarry Smith newC->itemType = 0; 890298fd71SBarry Smith newC->queue = NULL; 90af33a6ddSJed Brown newC->queueSize = 0; 91af33a6ddSJed Brown newC->queueMax = 0; 920298fd71SBarry Smith newC->queueLocal = NULL; 93af33a6ddSJed Brown newC->queueLocalSize = 0; 94af33a6ddSJed Brown newC->queueLocalMax = 0; 950298fd71SBarry Smith newC->queueRemote = NULL; 96af33a6ddSJed Brown newC->queueRemoteSize = 0; 97af33a6ddSJed Brown newC->queueRemoteMax = 0; 98af33a6ddSJed Brown newC->numNeighbors = 0; 990298fd71SBarry Smith newC->neighbors = NULL; 1000298fd71SBarry Smith newC->needCount = NULL; 1010298fd71SBarry Smith newC->localOffsets = NULL; 1020298fd71SBarry Smith newC->fillCount = NULL; 1030298fd71SBarry Smith newC->remoteOffsets = NULL; 1040298fd71SBarry Smith newC->request = NULL; 1050298fd71SBarry Smith newC->status = NULL; 1063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 107af33a6ddSJed Brown } 108af33a6ddSJed Brown 109af33a6ddSJed Brown /*@C 110af33a6ddSJed Brown CharacteristicSetType - Builds Characteristic for a particular solver. 111af33a6ddSJed Brown 112c3339decSBarry Smith Logically Collective 113af33a6ddSJed Brown 114af33a6ddSJed Brown Input Parameters: 115af33a6ddSJed Brown + c - the method of characteristics context 116af33a6ddSJed Brown - type - a known method 117af33a6ddSJed Brown 118af33a6ddSJed Brown Options Database Key: 119af33a6ddSJed Brown . -characteristic_type <method> - Sets the method; use -help for a list 120af33a6ddSJed Brown of available methods 121af33a6ddSJed Brown 122bcf0153eSBarry Smith Level: intermediate 123bcf0153eSBarry Smith 124af33a6ddSJed Brown Notes: 12530a76a96SBarry Smith See "include/petsccharacteristic.h" for available methods 126af33a6ddSJed Brown 127af33a6ddSJed Brown Normally, it is best to use the CharacteristicSetFromOptions() command and 128af33a6ddSJed Brown then set the Characteristic type from the options database rather than by using 129af33a6ddSJed Brown this routine. Using the options database provides the user with 130af33a6ddSJed Brown maximum flexibility in evaluating the many different Krylov methods. 131af33a6ddSJed Brown The CharacteristicSetType() routine is provided for those situations where it 132af33a6ddSJed Brown is necessary to set the iterative solver independently of the command 133af33a6ddSJed Brown line or options database. This might be the case, for example, when 134af33a6ddSJed Brown the choice of iterative solver changes during the execution of the 135af33a6ddSJed Brown program, and the user's application is taking responsibility for 136af33a6ddSJed Brown choosing the appropriate method. In other words, this routine is 137af33a6ddSJed Brown not for beginners. 138af33a6ddSJed Brown 1391cc06b55SBarry Smith .seealso: [](ch_ts), `CharacteristicType` 140af33a6ddSJed Brown @*/ 141d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicSetType(Characteristic c, CharacteristicType type) 142d71ae5a4SJacob Faibussowitsch { 143af33a6ddSJed Brown PetscBool match; 1445f80ce2aSJacob Faibussowitsch PetscErrorCode (*r)(Characteristic); 145af33a6ddSJed Brown 146af33a6ddSJed Brown PetscFunctionBegin; 147af33a6ddSJed Brown PetscValidHeaderSpecific(c, CHARACTERISTIC_CLASSID, 1); 1484f572ea9SToby Isaac PetscAssertPointer(type, 2); 149af33a6ddSJed Brown 1509566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)c, type, &match)); 1513ba16761SJacob Faibussowitsch if (match) PetscFunctionReturn(PETSC_SUCCESS); 152af33a6ddSJed Brown 153af33a6ddSJed Brown if (c->data) { 154af33a6ddSJed Brown /* destroy the old private Characteristic context */ 155dbbe0bcdSBarry Smith PetscUseTypeMethod(c, destroy); 1560298fd71SBarry Smith c->ops->destroy = NULL; 1572120983fSLisandro Dalcin c->data = NULL; 158af33a6ddSJed Brown } 159af33a6ddSJed Brown 1609566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(CharacteristicList, type, &r)); 1616adde796SStefano Zampini PetscCheck(r, PetscObjectComm((PetscObject)c), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown Characteristic type given: %s", type); 162af33a6ddSJed Brown c->setupcalled = 0; 1639566063dSJacob Faibussowitsch PetscCall((*r)(c)); 1649566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)c, type)); 1653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 166af33a6ddSJed Brown } 167af33a6ddSJed Brown 168af33a6ddSJed Brown /*@ 169af33a6ddSJed Brown CharacteristicSetUp - Sets up the internal data structures for the 170af33a6ddSJed Brown later use of an iterative solver. 171af33a6ddSJed Brown 172c3339decSBarry Smith Collective 173af33a6ddSJed Brown 174af33a6ddSJed Brown Input Parameter: 175b43aa488SJacob Faibussowitsch . c - iterative context obtained from CharacteristicCreate() 176af33a6ddSJed Brown 177af33a6ddSJed Brown Level: developer 178af33a6ddSJed Brown 1791cc06b55SBarry Smith .seealso: [](ch_ts), `CharacteristicCreate()`, `CharacteristicSolve()`, `CharacteristicDestroy()` 180af33a6ddSJed Brown @*/ 181d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicSetUp(Characteristic c) 182d71ae5a4SJacob Faibussowitsch { 183af33a6ddSJed Brown PetscFunctionBegin; 184af33a6ddSJed Brown PetscValidHeaderSpecific(c, CHARACTERISTIC_CLASSID, 1); 185af33a6ddSJed Brown 18648a46eb9SPierre Jolivet if (!((PetscObject)c)->type_name) PetscCall(CharacteristicSetType(c, CHARACTERISTICDA)); 187af33a6ddSJed Brown 1883ba16761SJacob Faibussowitsch if (c->setupcalled == 2) PetscFunctionReturn(PETSC_SUCCESS); 189af33a6ddSJed Brown 1909566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(CHARACTERISTIC_SetUp, c, NULL, NULL, NULL)); 191ad540459SPierre Jolivet if (!c->setupcalled) PetscUseTypeMethod(c, setup); 1929566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(CHARACTERISTIC_SetUp, c, NULL, NULL, NULL)); 193af33a6ddSJed Brown c->setupcalled = 2; 1943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 195af33a6ddSJed Brown } 196af33a6ddSJed Brown 197af33a6ddSJed Brown /*@C 1981c84c290SBarry Smith CharacteristicRegister - Adds a solver to the method of characteristics package. 1991c84c290SBarry Smith 2001c84c290SBarry Smith Not Collective 2011c84c290SBarry Smith 2021c84c290SBarry Smith Input Parameters: 2032fe279fdSBarry Smith + sname - name of a new user-defined solver 2042fe279fdSBarry Smith - function - routine to create method context 2051c84c290SBarry Smith 206bcf0153eSBarry Smith Level: advanced 207bcf0153eSBarry Smith 208b43aa488SJacob Faibussowitsch Example Usage: 2091c84c290SBarry Smith .vb 210bdf89e91SBarry Smith CharacteristicRegister("my_char", MyCharCreate); 2111c84c290SBarry Smith .ve 2121c84c290SBarry Smith 2131c84c290SBarry Smith Then, your Characteristic type can be chosen with the procedural interface via 2141c84c290SBarry Smith .vb 2151c84c290SBarry Smith CharacteristicCreate(MPI_Comm, Characteristic* &char); 2161c84c290SBarry Smith CharacteristicSetType(char,"my_char"); 2171c84c290SBarry Smith .ve 2181c84c290SBarry Smith or at runtime via the option 2191c84c290SBarry Smith .vb 2201c84c290SBarry Smith -characteristic_type my_char 2211c84c290SBarry Smith .ve 2221c84c290SBarry Smith 2231c84c290SBarry Smith Notes: 2241c84c290SBarry Smith CharacteristicRegister() may be called multiple times to add several user-defined solvers. 2251c84c290SBarry Smith 2261cc06b55SBarry Smith .seealso: [](ch_ts), `CharacteristicRegisterAll()`, `CharacteristicRegisterDestroy()` 227af33a6ddSJed Brown @*/ 228d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicRegister(const char sname[], PetscErrorCode (*function)(Characteristic)) 229d71ae5a4SJacob Faibussowitsch { 230af33a6ddSJed Brown PetscFunctionBegin; 2319566063dSJacob Faibussowitsch PetscCall(CharacteristicInitializePackage()); 2329566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&CharacteristicList, sname, function)); 2333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 234af33a6ddSJed Brown } 235af33a6ddSJed Brown 236d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicSetVelocityInterpolation(Characteristic c, DM da, Vec v, Vec vOld, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(Vec, PetscReal[], PetscInt, PetscInt[], PetscScalar[], void *), void *ctx) 237d71ae5a4SJacob Faibussowitsch { 238af33a6ddSJed Brown PetscFunctionBegin; 239af33a6ddSJed Brown c->velocityDA = da; 240af33a6ddSJed Brown c->velocity = v; 241af33a6ddSJed Brown c->velocityOld = vOld; 242af33a6ddSJed Brown c->numVelocityComp = numComponents; 243af33a6ddSJed Brown c->velocityComp = components; 244af33a6ddSJed Brown c->velocityInterp = interp; 245af33a6ddSJed Brown c->velocityCtx = ctx; 2463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 247af33a6ddSJed Brown } 248af33a6ddSJed Brown 249d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicSetVelocityInterpolationLocal(Characteristic c, DM da, Vec v, Vec vOld, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(void *, PetscReal[], PetscInt, PetscInt[], PetscScalar[], void *), void *ctx) 250d71ae5a4SJacob Faibussowitsch { 251af33a6ddSJed Brown PetscFunctionBegin; 252af33a6ddSJed Brown c->velocityDA = da; 253af33a6ddSJed Brown c->velocity = v; 254af33a6ddSJed Brown c->velocityOld = vOld; 255af33a6ddSJed Brown c->numVelocityComp = numComponents; 256af33a6ddSJed Brown c->velocityComp = components; 257af33a6ddSJed Brown c->velocityInterpLocal = interp; 258af33a6ddSJed Brown c->velocityCtx = ctx; 2593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 260af33a6ddSJed Brown } 261af33a6ddSJed Brown 262d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicSetFieldInterpolation(Characteristic c, DM da, Vec v, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(Vec, PetscReal[], PetscInt, PetscInt[], PetscScalar[], void *), void *ctx) 263d71ae5a4SJacob Faibussowitsch { 264af33a6ddSJed Brown PetscFunctionBegin; 265af33a6ddSJed Brown #if 0 2663c633725SBarry Smith PetscCheck(numComponents <= 2,PETSC_COMM_SELF,PETSC_ERR_SUP, "Fields with more than 2 components are not supported. Send mail to petsc-maint@mcs.anl.gov."); 267af33a6ddSJed Brown #endif 268af33a6ddSJed Brown c->fieldDA = da; 269af33a6ddSJed Brown c->field = v; 270af33a6ddSJed Brown c->numFieldComp = numComponents; 271af33a6ddSJed Brown c->fieldComp = components; 272af33a6ddSJed Brown c->fieldInterp = interp; 273af33a6ddSJed Brown c->fieldCtx = ctx; 2743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 275af33a6ddSJed Brown } 276af33a6ddSJed Brown 277d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicSetFieldInterpolationLocal(Characteristic c, DM da, Vec v, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(void *, PetscReal[], PetscInt, PetscInt[], PetscScalar[], void *), void *ctx) 278d71ae5a4SJacob Faibussowitsch { 279af33a6ddSJed Brown PetscFunctionBegin; 280af33a6ddSJed Brown #if 0 2813c633725SBarry Smith PetscCheck(numComponents <= 2,PETSC_COMM_SELF,PETSC_ERR_SUP, "Fields with more than 2 components are not supported. Send mail to petsc-maint@mcs.anl.gov."); 282af33a6ddSJed Brown #endif 283af33a6ddSJed Brown c->fieldDA = da; 284af33a6ddSJed Brown c->field = v; 285af33a6ddSJed Brown c->numFieldComp = numComponents; 286af33a6ddSJed Brown c->fieldComp = components; 287af33a6ddSJed Brown c->fieldInterpLocal = interp; 288af33a6ddSJed Brown c->fieldCtx = ctx; 2893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 290af33a6ddSJed Brown } 291af33a6ddSJed Brown 292d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicSolve(Characteristic c, PetscReal dt, Vec solution) 293d71ae5a4SJacob Faibussowitsch { 294af33a6ddSJed Brown CharacteristicPointDA2D Qi; 295af33a6ddSJed Brown DM da = c->velocityDA; 296af33a6ddSJed Brown Vec velocityLocal, velocityLocalOld; 297af33a6ddSJed Brown Vec fieldLocal; 298af33a6ddSJed Brown DMDALocalInfo info; 299af33a6ddSJed Brown PetscScalar **solArray; 300af33a6ddSJed Brown void *velocityArray; 301af33a6ddSJed Brown void *velocityArrayOld; 302af33a6ddSJed Brown void *fieldArray; 3031cc8b206SBarry Smith PetscScalar *interpIndices; 3041cc8b206SBarry Smith PetscScalar *velocityValues, *velocityValuesOld; 3051cc8b206SBarry Smith PetscScalar *fieldValues; 306af33a6ddSJed Brown PetscMPIInt rank; 307af33a6ddSJed Brown PetscInt dim; 308af33a6ddSJed Brown PetscMPIInt neighbors[9]; 309af33a6ddSJed Brown PetscInt dof; 310af33a6ddSJed Brown PetscInt gx, gy; 311af33a6ddSJed Brown PetscInt n, is, ie, js, je, comp; 312af33a6ddSJed Brown 313af33a6ddSJed Brown PetscFunctionBegin; 314af33a6ddSJed Brown c->queueSize = 0; 3159566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank)); 3169566063dSJacob Faibussowitsch PetscCall(DMDAGetNeighborsRank(da, neighbors)); 3179566063dSJacob Faibussowitsch PetscCall(CharacteristicSetNeighbors(c, 9, neighbors)); 3189566063dSJacob Faibussowitsch PetscCall(CharacteristicSetUp(c)); 319af33a6ddSJed Brown /* global and local grid info */ 3209566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, &gx, &gy, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 3219566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info)); 3229371c9d4SSatish Balay is = info.xs; 3239371c9d4SSatish Balay ie = info.xs + info.xm; 3249371c9d4SSatish Balay js = info.ys; 3259371c9d4SSatish Balay je = info.ys + info.ym; 326af33a6ddSJed Brown /* Allocation */ 3279566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dim, &interpIndices)); 3289566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(c->numVelocityComp, &velocityValues)); 3299566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(c->numVelocityComp, &velocityValuesOld)); 3309566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(c->numFieldComp, &fieldValues)); 3319566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(CHARACTERISTIC_Solve, NULL, NULL, NULL, NULL)); 332af33a6ddSJed Brown 3331d27aa22SBarry Smith /* 334af33a6ddSJed Brown PART 1, AT t-dt/2 3351d27aa22SBarry Smith */ 3369566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(CHARACTERISTIC_QueueSetup, NULL, NULL, NULL, NULL)); 337af33a6ddSJed Brown /* GET POSITION AT HALF TIME IN THE PAST */ 338af33a6ddSJed Brown if (c->velocityInterpLocal) { 3399566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(c->velocityDA, &velocityLocal)); 3409566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(c->velocityDA, &velocityLocalOld)); 3419566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(c->velocityDA, c->velocity, INSERT_VALUES, velocityLocal)); 3429566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(c->velocityDA, c->velocity, INSERT_VALUES, velocityLocal)); 3439566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(c->velocityDA, c->velocityOld, INSERT_VALUES, velocityLocalOld)); 3449566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(c->velocityDA, c->velocityOld, INSERT_VALUES, velocityLocalOld)); 3459566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(c->velocityDA, velocityLocal, &velocityArray)); 3469566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(c->velocityDA, velocityLocalOld, &velocityArrayOld)); 347af33a6ddSJed Brown } 3489566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Calculating position at t_{n - 1/2}\n")); 349af33a6ddSJed Brown for (Qi.j = js; Qi.j < je; Qi.j++) { 350af33a6ddSJed Brown for (Qi.i = is; Qi.i < ie; Qi.i++) { 351af33a6ddSJed Brown interpIndices[0] = Qi.i; 352af33a6ddSJed Brown interpIndices[1] = Qi.j; 3539566063dSJacob Faibussowitsch if (c->velocityInterpLocal) PetscCall(c->velocityInterpLocal(velocityArray, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx)); 3549566063dSJacob Faibussowitsch else PetscCall(c->velocityInterp(c->velocity, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx)); 355af33a6ddSJed Brown Qi.x = Qi.i - velocityValues[0] * dt / 2.0; 356af33a6ddSJed Brown Qi.y = Qi.j - velocityValues[1] * dt / 2.0; 357af33a6ddSJed Brown 358af33a6ddSJed Brown /* Determine whether the position at t - dt/2 is local */ 359af33a6ddSJed Brown Qi.proc = DMDAGetNeighborRelative(da, Qi.x, Qi.y); 360af33a6ddSJed Brown 361af33a6ddSJed Brown /* Check for Periodic boundaries and move all periodic points back onto the domain */ 362*f4f49eeaSPierre Jolivet PetscCall(DMDAMapCoordsToPeriodicDomain(da, &Qi.x, &Qi.y)); 3639566063dSJacob Faibussowitsch PetscCall(CharacteristicAddPoint(c, &Qi)); 364af33a6ddSJed Brown } 365af33a6ddSJed Brown } 3669566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(CHARACTERISTIC_QueueSetup, NULL, NULL, NULL, NULL)); 367af33a6ddSJed Brown 3689566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange, NULL, NULL, NULL, NULL)); 3699566063dSJacob Faibussowitsch PetscCall(CharacteristicSendCoordinatesBegin(c)); 3709566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange, NULL, NULL, NULL, NULL)); 371af33a6ddSJed Brown 3729566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(CHARACTERISTIC_HalfTimeLocal, NULL, NULL, NULL, NULL)); 373af33a6ddSJed Brown /* Calculate velocity at t_n+1/2 (local values) */ 3749566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Calculating local velocities at t_{n - 1/2}\n")); 375af33a6ddSJed Brown for (n = 0; n < c->queueSize; n++) { 376af33a6ddSJed Brown Qi = c->queue[n]; 377af33a6ddSJed Brown if (c->neighbors[Qi.proc] == rank) { 378af33a6ddSJed Brown interpIndices[0] = Qi.x; 379af33a6ddSJed Brown interpIndices[1] = Qi.y; 380af33a6ddSJed Brown if (c->velocityInterpLocal) { 3819566063dSJacob Faibussowitsch PetscCall(c->velocityInterpLocal(velocityArray, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx)); 3829566063dSJacob Faibussowitsch PetscCall(c->velocityInterpLocal(velocityArrayOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx)); 383af33a6ddSJed Brown } else { 3849566063dSJacob Faibussowitsch PetscCall(c->velocityInterp(c->velocity, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx)); 3859566063dSJacob Faibussowitsch PetscCall(c->velocityInterp(c->velocityOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx)); 386af33a6ddSJed Brown } 387af33a6ddSJed Brown Qi.x = 0.5 * (velocityValues[0] + velocityValuesOld[0]); 388af33a6ddSJed Brown Qi.y = 0.5 * (velocityValues[1] + velocityValuesOld[1]); 389af33a6ddSJed Brown } 390af33a6ddSJed Brown c->queue[n] = Qi; 391af33a6ddSJed Brown } 3929566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(CHARACTERISTIC_HalfTimeLocal, NULL, NULL, NULL, NULL)); 393af33a6ddSJed Brown 3949566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange, NULL, NULL, NULL, NULL)); 3959566063dSJacob Faibussowitsch PetscCall(CharacteristicSendCoordinatesEnd(c)); 3969566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange, NULL, NULL, NULL, NULL)); 397af33a6ddSJed Brown 398af33a6ddSJed Brown /* Calculate velocity at t_n+1/2 (fill remote requests) */ 3999566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(CHARACTERISTIC_HalfTimeRemote, NULL, NULL, NULL, NULL)); 40063a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Calculating %" PetscInt_FMT " remote velocities at t_{n - 1/2}\n", c->queueRemoteSize)); 401af33a6ddSJed Brown for (n = 0; n < c->queueRemoteSize; n++) { 402af33a6ddSJed Brown Qi = c->queueRemote[n]; 403af33a6ddSJed Brown interpIndices[0] = Qi.x; 404af33a6ddSJed Brown interpIndices[1] = Qi.y; 405af33a6ddSJed Brown if (c->velocityInterpLocal) { 4069566063dSJacob Faibussowitsch PetscCall(c->velocityInterpLocal(velocityArray, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx)); 4079566063dSJacob Faibussowitsch PetscCall(c->velocityInterpLocal(velocityArrayOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx)); 408af33a6ddSJed Brown } else { 4099566063dSJacob Faibussowitsch PetscCall(c->velocityInterp(c->velocity, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx)); 4109566063dSJacob Faibussowitsch PetscCall(c->velocityInterp(c->velocityOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx)); 411af33a6ddSJed Brown } 412af33a6ddSJed Brown Qi.x = 0.5 * (velocityValues[0] + velocityValuesOld[0]); 413af33a6ddSJed Brown Qi.y = 0.5 * (velocityValues[1] + velocityValuesOld[1]); 414af33a6ddSJed Brown c->queueRemote[n] = Qi; 415af33a6ddSJed Brown } 4169566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(CHARACTERISTIC_HalfTimeRemote, NULL, NULL, NULL, NULL)); 4179566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange, NULL, NULL, NULL, NULL)); 4189566063dSJacob Faibussowitsch PetscCall(CharacteristicGetValuesBegin(c)); 4199566063dSJacob Faibussowitsch PetscCall(CharacteristicGetValuesEnd(c)); 420af33a6ddSJed Brown if (c->velocityInterpLocal) { 4219566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(c->velocityDA, velocityLocal, &velocityArray)); 4229566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(c->velocityDA, velocityLocalOld, &velocityArrayOld)); 4239566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(c->velocityDA, &velocityLocal)); 4249566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(c->velocityDA, &velocityLocalOld)); 425af33a6ddSJed Brown } 4269566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange, NULL, NULL, NULL, NULL)); 427af33a6ddSJed Brown 4281d27aa22SBarry Smith /* 429af33a6ddSJed Brown PART 2, AT t-dt 4301d27aa22SBarry Smith */ 431af33a6ddSJed Brown 432af33a6ddSJed Brown /* GET POSITION AT t_n (local values) */ 4339566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(CHARACTERISTIC_FullTimeLocal, NULL, NULL, NULL, NULL)); 4349566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Calculating position at t_{n}\n")); 435af33a6ddSJed Brown for (n = 0; n < c->queueSize; n++) { 436af33a6ddSJed Brown Qi = c->queue[n]; 437af33a6ddSJed Brown Qi.x = Qi.i - Qi.x * dt; 438af33a6ddSJed Brown Qi.y = Qi.j - Qi.y * dt; 439af33a6ddSJed Brown 440af33a6ddSJed Brown /* Determine whether the position at t-dt is local */ 441af33a6ddSJed Brown Qi.proc = DMDAGetNeighborRelative(da, Qi.x, Qi.y); 442af33a6ddSJed Brown 443af33a6ddSJed Brown /* Check for Periodic boundaries and move all periodic points back onto the domain */ 444*f4f49eeaSPierre Jolivet PetscCall(DMDAMapCoordsToPeriodicDomain(da, &Qi.x, &Qi.y)); 445af33a6ddSJed Brown 446af33a6ddSJed Brown c->queue[n] = Qi; 447af33a6ddSJed Brown } 4489566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(CHARACTERISTIC_FullTimeLocal, NULL, NULL, NULL, NULL)); 449af33a6ddSJed Brown 4509566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange, NULL, NULL, NULL, NULL)); 4519566063dSJacob Faibussowitsch PetscCall(CharacteristicSendCoordinatesBegin(c)); 4529566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange, NULL, NULL, NULL, NULL)); 453af33a6ddSJed Brown 454af33a6ddSJed Brown /* GET VALUE AT FULL TIME IN THE PAST (LOCAL REQUESTS) */ 4559566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(CHARACTERISTIC_FullTimeLocal, NULL, NULL, NULL, NULL)); 456af33a6ddSJed Brown if (c->fieldInterpLocal) { 4579566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(c->fieldDA, &fieldLocal)); 4589566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(c->fieldDA, c->field, INSERT_VALUES, fieldLocal)); 4599566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(c->fieldDA, c->field, INSERT_VALUES, fieldLocal)); 4609566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(c->fieldDA, fieldLocal, &fieldArray)); 461af33a6ddSJed Brown } 4629566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Calculating local field at t_{n}\n")); 463af33a6ddSJed Brown for (n = 0; n < c->queueSize; n++) { 464af33a6ddSJed Brown if (c->neighbors[c->queue[n].proc] == rank) { 465af33a6ddSJed Brown interpIndices[0] = c->queue[n].x; 466af33a6ddSJed Brown interpIndices[1] = c->queue[n].y; 4679566063dSJacob Faibussowitsch if (c->fieldInterpLocal) PetscCall(c->fieldInterpLocal(fieldArray, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx)); 4689566063dSJacob Faibussowitsch else PetscCall(c->fieldInterp(c->field, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx)); 469bbd56ea5SKarl Rupp for (comp = 0; comp < c->numFieldComp; comp++) c->queue[n].field[comp] = fieldValues[comp]; 470af33a6ddSJed Brown } 471af33a6ddSJed Brown } 4729566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(CHARACTERISTIC_FullTimeLocal, NULL, NULL, NULL, NULL)); 473af33a6ddSJed Brown 4749566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange, NULL, NULL, NULL, NULL)); 4759566063dSJacob Faibussowitsch PetscCall(CharacteristicSendCoordinatesEnd(c)); 4769566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange, NULL, NULL, NULL, NULL)); 477af33a6ddSJed Brown 478af33a6ddSJed Brown /* GET VALUE AT FULL TIME IN THE PAST (REMOTE REQUESTS) */ 4799566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(CHARACTERISTIC_FullTimeRemote, NULL, NULL, NULL, NULL)); 48063a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Calculating %" PetscInt_FMT " remote field points at t_{n}\n", c->queueRemoteSize)); 481af33a6ddSJed Brown for (n = 0; n < c->queueRemoteSize; n++) { 482af33a6ddSJed Brown interpIndices[0] = c->queueRemote[n].x; 483af33a6ddSJed Brown interpIndices[1] = c->queueRemote[n].y; 484af33a6ddSJed Brown 485af33a6ddSJed Brown /* for debugging purposes */ 486af33a6ddSJed Brown if (1) { /* hacked bounds test...let's do better */ 4879371c9d4SSatish Balay PetscScalar im = interpIndices[0]; 4889371c9d4SSatish Balay PetscScalar jm = interpIndices[1]; 489af33a6ddSJed Brown 49063a3b9bcSJacob Faibussowitsch PetscCheck((im >= (PetscScalar)is - 1.) && (im <= (PetscScalar)ie) && (jm >= (PetscScalar)js - 1.) && (jm <= (PetscScalar)je), PETSC_COMM_SELF, PETSC_ERR_LIB, "Nonlocal point: (%g,%g)", (double)PetscAbsScalar(im), (double)PetscAbsScalar(jm)); 491af33a6ddSJed Brown } 492af33a6ddSJed Brown 4939566063dSJacob Faibussowitsch if (c->fieldInterpLocal) PetscCall(c->fieldInterpLocal(fieldArray, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx)); 4949566063dSJacob Faibussowitsch else PetscCall(c->fieldInterp(c->field, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx)); 495bbd56ea5SKarl Rupp for (comp = 0; comp < c->numFieldComp; comp++) c->queueRemote[n].field[comp] = fieldValues[comp]; 496af33a6ddSJed Brown } 4979566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(CHARACTERISTIC_FullTimeRemote, NULL, NULL, NULL, NULL)); 498af33a6ddSJed Brown 4999566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange, NULL, NULL, NULL, NULL)); 5009566063dSJacob Faibussowitsch PetscCall(CharacteristicGetValuesBegin(c)); 5019566063dSJacob Faibussowitsch PetscCall(CharacteristicGetValuesEnd(c)); 502af33a6ddSJed Brown if (c->fieldInterpLocal) { 5039566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(c->fieldDA, fieldLocal, &fieldArray)); 5049566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(c->fieldDA, &fieldLocal)); 505af33a6ddSJed Brown } 5069566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange, NULL, NULL, NULL, NULL)); 507af33a6ddSJed Brown 508af33a6ddSJed Brown /* Return field of characteristics at t_n-1 */ 5099566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(CHARACTERISTIC_DAUpdate, NULL, NULL, NULL, NULL)); 5109566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(c->fieldDA, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL)); 5119566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(c->fieldDA, solution, &solArray)); 512af33a6ddSJed Brown for (n = 0; n < c->queueSize; n++) { 513af33a6ddSJed Brown Qi = c->queue[n]; 514bbd56ea5SKarl Rupp for (comp = 0; comp < c->numFieldComp; comp++) solArray[Qi.j][Qi.i * dof + c->fieldComp[comp]] = Qi.field[comp]; 515af33a6ddSJed Brown } 5169566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(c->fieldDA, solution, &solArray)); 5179566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(CHARACTERISTIC_DAUpdate, NULL, NULL, NULL, NULL)); 5189566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(CHARACTERISTIC_Solve, NULL, NULL, NULL, NULL)); 519af33a6ddSJed Brown 520af33a6ddSJed Brown /* Cleanup */ 5219566063dSJacob Faibussowitsch PetscCall(PetscFree(interpIndices)); 5229566063dSJacob Faibussowitsch PetscCall(PetscFree(velocityValues)); 5239566063dSJacob Faibussowitsch PetscCall(PetscFree(velocityValuesOld)); 5249566063dSJacob Faibussowitsch PetscCall(PetscFree(fieldValues)); 5253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 526af33a6ddSJed Brown } 527af33a6ddSJed Brown 528d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicSetNeighbors(Characteristic c, PetscInt numNeighbors, PetscMPIInt neighbors[]) 529d71ae5a4SJacob Faibussowitsch { 530af33a6ddSJed Brown PetscFunctionBegin; 531af33a6ddSJed Brown c->numNeighbors = numNeighbors; 5329566063dSJacob Faibussowitsch PetscCall(PetscFree(c->neighbors)); 5339566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numNeighbors, &c->neighbors)); 5349566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(c->neighbors, neighbors, numNeighbors)); 5353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 536af33a6ddSJed Brown } 537af33a6ddSJed Brown 538d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicAddPoint(Characteristic c, CharacteristicPointDA2D *point) 539d71ae5a4SJacob Faibussowitsch { 540af33a6ddSJed Brown PetscFunctionBegin; 54163a3b9bcSJacob Faibussowitsch PetscCheck(c->queueSize < c->queueMax, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Exceeded maximum queue size %" PetscInt_FMT, c->queueMax); 542af33a6ddSJed Brown c->queue[c->queueSize++] = *point; 5433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 544af33a6ddSJed Brown } 545af33a6ddSJed Brown 5463ba16761SJacob Faibussowitsch PetscErrorCode CharacteristicSendCoordinatesBegin(Characteristic c) 547d71ae5a4SJacob Faibussowitsch { 548af33a6ddSJed Brown PetscMPIInt rank, tag = 121; 549af33a6ddSJed Brown PetscInt i, n; 550af33a6ddSJed Brown 551af33a6ddSJed Brown PetscFunctionBegin; 5529566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank)); 5539566063dSJacob Faibussowitsch PetscCall(CharacteristicHeapSort(c, c->queue, c->queueSize)); 5549566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(c->needCount, c->numNeighbors)); 555bbd56ea5SKarl Rupp for (i = 0; i < c->queueSize; i++) c->needCount[c->queue[i].proc]++; 556af33a6ddSJed Brown c->fillCount[0] = 0; 557*f4f49eeaSPierre Jolivet for (n = 1; n < c->numNeighbors; n++) PetscCallMPI(MPI_Irecv(&c->fillCount[n], 1, MPIU_INT, c->neighbors[n], tag, PetscObjectComm((PetscObject)c), &c->request[n - 1])); 558*f4f49eeaSPierre Jolivet for (n = 1; n < c->numNeighbors; n++) PetscCallMPI(MPI_Send(&c->needCount[n], 1, MPIU_INT, c->neighbors[n], tag, PetscObjectComm((PetscObject)c))); 5599566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(c->numNeighbors - 1, c->request, c->status)); 560af33a6ddSJed Brown /* Initialize the remote queue */ 561af33a6ddSJed Brown c->queueLocalMax = c->localOffsets[0] = 0; 562af33a6ddSJed Brown c->queueRemoteMax = c->remoteOffsets[0] = 0; 563af33a6ddSJed Brown for (n = 1; n < c->numNeighbors; n++) { 564af33a6ddSJed Brown c->remoteOffsets[n] = c->queueRemoteMax; 565af33a6ddSJed Brown c->queueRemoteMax += c->fillCount[n]; 566af33a6ddSJed Brown c->localOffsets[n] = c->queueLocalMax; 567af33a6ddSJed Brown c->queueLocalMax += c->needCount[n]; 568af33a6ddSJed Brown } 569af33a6ddSJed Brown /* HACK BEGIN */ 570bbd56ea5SKarl Rupp for (n = 1; n < c->numNeighbors; n++) c->localOffsets[n] += c->needCount[0]; 571af33a6ddSJed Brown c->needCount[0] = 0; 572af33a6ddSJed Brown /* HACK END */ 573af33a6ddSJed Brown if (c->queueRemoteMax) { 5749566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(c->queueRemoteMax, &c->queueRemote)); 5750298fd71SBarry Smith } else c->queueRemote = NULL; 576af33a6ddSJed Brown c->queueRemoteSize = c->queueRemoteMax; 577af33a6ddSJed Brown 578af33a6ddSJed Brown /* Send and Receive requests for values at t_n+1/2, giving the coordinates for interpolation */ 579af33a6ddSJed Brown for (n = 1; n < c->numNeighbors; n++) { 58063a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Receiving %" PetscInt_FMT " requests for values from proc %d\n", c->fillCount[n], c->neighbors[n])); 581*f4f49eeaSPierre Jolivet PetscCallMPI(MPI_Irecv(&(c->queueRemote[c->remoteOffsets[n]]), c->fillCount[n], c->itemType, c->neighbors[n], tag, PetscObjectComm((PetscObject)c), &c->request[n - 1])); 582af33a6ddSJed Brown } 583af33a6ddSJed Brown for (n = 1; n < c->numNeighbors; n++) { 58463a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Sending %" PetscInt_FMT " requests for values from proc %d\n", c->needCount[n], c->neighbors[n])); 5859566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send(&(c->queue[c->localOffsets[n]]), c->needCount[n], c->itemType, c->neighbors[n], tag, PetscObjectComm((PetscObject)c))); 586af33a6ddSJed Brown } 5873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 588af33a6ddSJed Brown } 589af33a6ddSJed Brown 590d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicSendCoordinatesEnd(Characteristic c) 591d71ae5a4SJacob Faibussowitsch { 592af33a6ddSJed Brown #if 0 593af33a6ddSJed Brown PetscMPIInt rank; 594af33a6ddSJed Brown PetscInt n; 595af33a6ddSJed Brown #endif 596af33a6ddSJed Brown 597af33a6ddSJed Brown PetscFunctionBegin; 5989566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(c->numNeighbors - 1, c->request, c->status)); 599af33a6ddSJed Brown #if 0 6009566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank)); 601af33a6ddSJed Brown for (n = 0; n < c->queueRemoteSize; n++) { 6023c633725SBarry Smith PetscCheck(c->neighbors[c->queueRemote[n].proc] != rank,PETSC_COMM_SELF,PETSC_ERR_PLIB, "This is messed up, n = %d proc = %d", n, c->queueRemote[n].proc); 603af33a6ddSJed Brown } 604af33a6ddSJed Brown #endif 6053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 606af33a6ddSJed Brown } 607af33a6ddSJed Brown 608d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicGetValuesBegin(Characteristic c) 609d71ae5a4SJacob Faibussowitsch { 610af33a6ddSJed Brown PetscMPIInt tag = 121; 611af33a6ddSJed Brown PetscInt n; 612af33a6ddSJed Brown 613af33a6ddSJed Brown PetscFunctionBegin; 614d5b43468SJose E. Roman /* SEND AND RECEIVE FILLED REQUESTS for velocities at t_n+1/2 */ 615*f4f49eeaSPierre Jolivet for (n = 1; n < c->numNeighbors; n++) PetscCallMPI(MPI_Irecv(&(c->queue[c->localOffsets[n]]), c->needCount[n], c->itemType, c->neighbors[n], tag, PetscObjectComm((PetscObject)c), &c->request[n - 1])); 61648a46eb9SPierre Jolivet for (n = 1; n < c->numNeighbors; n++) PetscCallMPI(MPI_Send(&(c->queueRemote[c->remoteOffsets[n]]), c->fillCount[n], c->itemType, c->neighbors[n], tag, PetscObjectComm((PetscObject)c))); 6173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 618af33a6ddSJed Brown } 619af33a6ddSJed Brown 620d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicGetValuesEnd(Characteristic c) 621d71ae5a4SJacob Faibussowitsch { 622af33a6ddSJed Brown PetscFunctionBegin; 6239566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(c->numNeighbors - 1, c->request, c->status)); 624af33a6ddSJed Brown /* Free queue of requests from other procs */ 6259566063dSJacob Faibussowitsch PetscCall(PetscFree(c->queueRemote)); 6263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 627af33a6ddSJed Brown } 628af33a6ddSJed Brown 629af33a6ddSJed Brown /* 630af33a6ddSJed Brown Based on code from http://linux.wku.edu/~lamonml/algor/sort/heap.html 631af33a6ddSJed Brown */ 63266976f2fSJacob Faibussowitsch static PetscErrorCode CharacteristicHeapSort(Characteristic c, Queue queue, PetscInt size) 633af33a6ddSJed Brown { 634af33a6ddSJed Brown CharacteristicPointDA2D temp; 6350b8f38c8SBarry Smith PetscInt n; 636af33a6ddSJed Brown 6370b8f38c8SBarry Smith PetscFunctionBegin; 638af33a6ddSJed Brown if (0) { /* Check the order of the queue before sorting */ 6399566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Before Heap sort\n")); 64048a46eb9SPierre Jolivet for (n = 0; n < size; n++) PetscCall(PetscInfo(NULL, "%" PetscInt_FMT " %d\n", n, queue[n].proc)); 641af33a6ddSJed Brown } 642af33a6ddSJed Brown 643af33a6ddSJed Brown /* SORTING PHASE */ 6449371c9d4SSatish Balay for (n = (size / 2) - 1; n >= 0; n--) { PetscCall(CharacteristicSiftDown(c, queue, n, size - 1)); /* Rich had size-1 here, Matt had size*/ } 645af33a6ddSJed Brown for (n = size - 1; n >= 1; n--) { 646af33a6ddSJed Brown temp = queue[0]; 647af33a6ddSJed Brown queue[0] = queue[n]; 648af33a6ddSJed Brown queue[n] = temp; 6499566063dSJacob Faibussowitsch PetscCall(CharacteristicSiftDown(c, queue, 0, n - 1)); 650af33a6ddSJed Brown } 651af33a6ddSJed Brown if (0) { /* Check the order of the queue after sorting */ 6529566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Avter Heap sort\n")); 65348a46eb9SPierre Jolivet for (n = 0; n < size; n++) PetscCall(PetscInfo(NULL, "%" PetscInt_FMT " %d\n", n, queue[n].proc)); 654af33a6ddSJed Brown } 6553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 656af33a6ddSJed Brown } 657af33a6ddSJed Brown 658af33a6ddSJed Brown /* 659af33a6ddSJed Brown Based on code from http://linux.wku.edu/~lamonml/algor/sort/heap.html 660af33a6ddSJed Brown */ 66166976f2fSJacob Faibussowitsch static PetscErrorCode CharacteristicSiftDown(Characteristic c, Queue queue, PetscInt root, PetscInt bottom) 662af33a6ddSJed Brown { 663af33a6ddSJed Brown PetscBool done = PETSC_FALSE; 664af33a6ddSJed Brown PetscInt maxChild; 665af33a6ddSJed Brown CharacteristicPointDA2D temp; 666af33a6ddSJed Brown 667af33a6ddSJed Brown PetscFunctionBegin; 668af33a6ddSJed Brown while ((root * 2 <= bottom) && (!done)) { 669af33a6ddSJed Brown if (root * 2 == bottom) maxChild = root * 2; 670af33a6ddSJed Brown else if (queue[root * 2].proc > queue[root * 2 + 1].proc) maxChild = root * 2; 671af33a6ddSJed Brown else maxChild = root * 2 + 1; 672af33a6ddSJed Brown 673af33a6ddSJed Brown if (queue[root].proc < queue[maxChild].proc) { 674af33a6ddSJed Brown temp = queue[root]; 675af33a6ddSJed Brown queue[root] = queue[maxChild]; 676af33a6ddSJed Brown queue[maxChild] = temp; 677af33a6ddSJed Brown root = maxChild; 678db4deed7SKarl Rupp } else done = PETSC_TRUE; 679af33a6ddSJed Brown } 6803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 681af33a6ddSJed Brown } 682af33a6ddSJed Brown 683af33a6ddSJed Brown /* [center, left, top-left, top, top-right, right, bottom-right, bottom, bottom-left] */ 68466976f2fSJacob Faibussowitsch static PetscErrorCode DMDAGetNeighborsRank(DM da, PetscMPIInt neighbors[]) 685d71ae5a4SJacob Faibussowitsch { 686bff4a2f0SMatthew G. Knepley DMBoundaryType bx, by; 687af33a6ddSJed Brown PetscBool IPeriodic = PETSC_FALSE, JPeriodic = PETSC_FALSE; 688af33a6ddSJed Brown MPI_Comm comm; 689af33a6ddSJed Brown PetscMPIInt rank; 690af33a6ddSJed Brown PetscInt **procs, pi, pj, pim, pip, pjm, pjp, PI, PJ; 691af33a6ddSJed Brown 692af33a6ddSJed Brown PetscFunctionBegin; 6939566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)da, &comm)); 6949566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 6959566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, NULL, NULL, NULL, NULL, &PI, &PJ, NULL, NULL, NULL, &bx, &by, NULL, NULL)); 696af33a6ddSJed Brown 697bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC) IPeriodic = PETSC_TRUE; 698bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_PERIODIC) JPeriodic = PETSC_TRUE; 699af33a6ddSJed Brown 700af33a6ddSJed Brown neighbors[0] = rank; 701af33a6ddSJed Brown rank = 0; 7029566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(PJ, &procs)); 703af33a6ddSJed Brown for (pj = 0; pj < PJ; pj++) { 704*f4f49eeaSPierre Jolivet PetscCall(PetscMalloc1(PI, &procs[pj])); 705af33a6ddSJed Brown for (pi = 0; pi < PI; pi++) { 706af33a6ddSJed Brown procs[pj][pi] = rank; 707af33a6ddSJed Brown rank++; 708af33a6ddSJed Brown } 709af33a6ddSJed Brown } 710af33a6ddSJed Brown 711af33a6ddSJed Brown pi = neighbors[0] % PI; 712af33a6ddSJed Brown pj = neighbors[0] / PI; 7139371c9d4SSatish Balay pim = pi - 1; 7149371c9d4SSatish Balay if (pim < 0) pim = PI - 1; 715af33a6ddSJed Brown pip = (pi + 1) % PI; 7169371c9d4SSatish Balay pjm = pj - 1; 7179371c9d4SSatish Balay if (pjm < 0) pjm = PJ - 1; 718af33a6ddSJed Brown pjp = (pj + 1) % PJ; 719af33a6ddSJed Brown 720af33a6ddSJed Brown neighbors[1] = procs[pj][pim]; 721af33a6ddSJed Brown neighbors[2] = procs[pjp][pim]; 722af33a6ddSJed Brown neighbors[3] = procs[pjp][pi]; 723af33a6ddSJed Brown neighbors[4] = procs[pjp][pip]; 724af33a6ddSJed Brown neighbors[5] = procs[pj][pip]; 725af33a6ddSJed Brown neighbors[6] = procs[pjm][pip]; 726af33a6ddSJed Brown neighbors[7] = procs[pjm][pi]; 727af33a6ddSJed Brown neighbors[8] = procs[pjm][pim]; 728af33a6ddSJed Brown 729af33a6ddSJed Brown if (!IPeriodic) { 730af33a6ddSJed Brown if (pi == 0) neighbors[1] = neighbors[2] = neighbors[8] = neighbors[0]; 731af33a6ddSJed Brown if (pi == PI - 1) neighbors[4] = neighbors[5] = neighbors[6] = neighbors[0]; 732af33a6ddSJed Brown } 733af33a6ddSJed Brown 734af33a6ddSJed Brown if (!JPeriodic) { 735af33a6ddSJed Brown if (pj == 0) neighbors[6] = neighbors[7] = neighbors[8] = neighbors[0]; 736af33a6ddSJed Brown if (pj == PJ - 1) neighbors[2] = neighbors[3] = neighbors[4] = neighbors[0]; 737af33a6ddSJed Brown } 738af33a6ddSJed Brown 73948a46eb9SPierre Jolivet for (pj = 0; pj < PJ; pj++) PetscCall(PetscFree(procs[pj])); 7409566063dSJacob Faibussowitsch PetscCall(PetscFree(procs)); 7413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 742af33a6ddSJed Brown } 743af33a6ddSJed Brown 744af33a6ddSJed Brown /* 745af33a6ddSJed Brown SUBDOMAIN NEIGHBORHOOD PROCESS MAP: 746af33a6ddSJed Brown 2 | 3 | 4 747af33a6ddSJed Brown __|___|__ 748af33a6ddSJed Brown 1 | 0 | 5 749af33a6ddSJed Brown __|___|__ 750af33a6ddSJed Brown 8 | 7 | 6 751af33a6ddSJed Brown | | 752af33a6ddSJed Brown */ 75366976f2fSJacob Faibussowitsch static PetscInt DMDAGetNeighborRelative(DM da, PetscReal ir, PetscReal jr) 754d71ae5a4SJacob Faibussowitsch { 755af33a6ddSJed Brown DMDALocalInfo info; 7561cc8b206SBarry Smith PetscReal is, ie, js, je; 757af33a6ddSJed Brown 7589566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info)); 7599371c9d4SSatish Balay is = (PetscReal)info.xs - 0.5; 7609371c9d4SSatish Balay ie = (PetscReal)info.xs + info.xm - 0.5; 7619371c9d4SSatish Balay js = (PetscReal)info.ys - 0.5; 7629371c9d4SSatish Balay je = (PetscReal)info.ys + info.ym - 0.5; 763af33a6ddSJed Brown 764af33a6ddSJed Brown if (ir >= is && ir <= ie) { /* center column */ 765bbd56ea5SKarl Rupp if (jr >= js && jr <= je) return 0; 766bbd56ea5SKarl Rupp else if (jr < js) return 7; 767bbd56ea5SKarl Rupp else return 3; 768af33a6ddSJed Brown } else if (ir < is) { /* left column */ 769bbd56ea5SKarl Rupp if (jr >= js && jr <= je) return 1; 770bbd56ea5SKarl Rupp else if (jr < js) return 8; 771bbd56ea5SKarl Rupp else return 2; 772af33a6ddSJed Brown } else { /* right column */ 773bbd56ea5SKarl Rupp if (jr >= js && jr <= je) return 5; 774bbd56ea5SKarl Rupp else if (jr < js) return 6; 775bbd56ea5SKarl Rupp else return 4; 776af33a6ddSJed Brown } 777af33a6ddSJed Brown } 778