1af33a6ddSJed Brown 2af0996ceSBarry Smith #include <petsc/private/characteristicimpl.h> /*I "petsccharacteristic.h" I*/ 31e25c274SJed Brown #include <petscdmda.h> 4665c2dedSJed Brown #include <petscviewer.h> 5af33a6ddSJed Brown 6af33a6ddSJed Brown PetscClassId CHARACTERISTIC_CLASSID; 7af33a6ddSJed Brown PetscLogEvent CHARACTERISTIC_SetUp, CHARACTERISTIC_Solve, CHARACTERISTIC_QueueSetup, CHARACTERISTIC_DAUpdate; 8af33a6ddSJed Brown PetscLogEvent CHARACTERISTIC_HalfTimeLocal, CHARACTERISTIC_HalfTimeRemote, CHARACTERISTIC_HalfTimeExchange; 9af33a6ddSJed Brown PetscLogEvent CHARACTERISTIC_FullTimeLocal, CHARACTERISTIC_FullTimeRemote, CHARACTERISTIC_FullTimeExchange; 10af33a6ddSJed Brown /* 11af33a6ddSJed Brown Contains the list of registered characteristic routines 12af33a6ddSJed Brown */ 130298fd71SBarry Smith PetscFunctionList CharacteristicList = NULL; 14140e18c1SBarry Smith PetscBool CharacteristicRegisterAllCalled = PETSC_FALSE; 15af33a6ddSJed Brown 16af33a6ddSJed Brown PetscErrorCode DMDAGetNeighborsRank(DM, PetscMPIInt []); 171cc8b206SBarry Smith PetscInt DMDAGetNeighborRelative(DM, PetscReal, PetscReal); 18af33a6ddSJed Brown PetscErrorCode DMDAMapToPeriodicDomain(DM, PetscScalar []); 19af33a6ddSJed Brown 200b8f38c8SBarry Smith PetscErrorCode CharacteristicHeapSort(Characteristic, Queue, PetscInt); 210b8f38c8SBarry Smith PetscErrorCode CharacteristicSiftDown(Characteristic, Queue, PetscInt, PetscInt); 22af33a6ddSJed Brown 23af33a6ddSJed Brown PetscErrorCode CharacteristicView(Characteristic c, PetscViewer viewer) 24af33a6ddSJed Brown { 25af33a6ddSJed Brown PetscBool iascii; 26af33a6ddSJed Brown PetscErrorCode ierr; 27af33a6ddSJed Brown 28af33a6ddSJed Brown PetscFunctionBegin; 29af33a6ddSJed Brown PetscValidHeaderSpecific(c, CHARACTERISTIC_CLASSID, 1); 30af33a6ddSJed Brown if (!viewer) { 31ce94432eSBarry Smith ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)c),&viewer);CHKERRQ(ierr); 32af33a6ddSJed Brown } 33af33a6ddSJed Brown PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 34af33a6ddSJed Brown PetscCheckSameComm(c, 1, viewer, 2); 35af33a6ddSJed Brown 36251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 37af33a6ddSJed Brown if (!iascii) { 38af33a6ddSJed Brown if (c->ops->view) { 39af33a6ddSJed Brown ierr = (*c->ops->view)(c, viewer);CHKERRQ(ierr); 40af33a6ddSJed Brown } 41af33a6ddSJed Brown } 42af33a6ddSJed Brown PetscFunctionReturn(0); 43af33a6ddSJed Brown } 44af33a6ddSJed Brown 456bf464f9SBarry Smith PetscErrorCode CharacteristicDestroy(Characteristic *c) 46af33a6ddSJed Brown { 47af33a6ddSJed Brown PetscErrorCode ierr; 48af33a6ddSJed Brown 49af33a6ddSJed Brown PetscFunctionBegin; 506bf464f9SBarry Smith if (!*c) PetscFunctionReturn(0); 516bf464f9SBarry Smith PetscValidHeaderSpecific(*c, CHARACTERISTIC_CLASSID, 1); 526bf464f9SBarry Smith if (--((PetscObject)(*c))->refct > 0) PetscFunctionReturn(0); 53af33a6ddSJed Brown 546bf464f9SBarry Smith if ((*c)->ops->destroy) { 556bf464f9SBarry Smith ierr = (*(*c)->ops->destroy)((*c));CHKERRQ(ierr); 56af33a6ddSJed Brown } 576bf464f9SBarry Smith ierr = MPI_Type_free(&(*c)->itemType);CHKERRQ(ierr); 586bf464f9SBarry Smith ierr = PetscFree((*c)->queue);CHKERRQ(ierr); 596bf464f9SBarry Smith ierr = PetscFree((*c)->queueLocal);CHKERRQ(ierr); 606bf464f9SBarry Smith ierr = PetscFree((*c)->queueRemote);CHKERRQ(ierr); 616bf464f9SBarry Smith ierr = PetscFree((*c)->neighbors);CHKERRQ(ierr); 626bf464f9SBarry Smith ierr = PetscFree((*c)->needCount);CHKERRQ(ierr); 636bf464f9SBarry Smith ierr = PetscFree((*c)->localOffsets);CHKERRQ(ierr); 646bf464f9SBarry Smith ierr = PetscFree((*c)->fillCount);CHKERRQ(ierr); 656bf464f9SBarry Smith ierr = PetscFree((*c)->remoteOffsets);CHKERRQ(ierr); 666bf464f9SBarry Smith ierr = PetscFree((*c)->request);CHKERRQ(ierr); 676bf464f9SBarry Smith ierr = PetscFree((*c)->status);CHKERRQ(ierr); 68af33a6ddSJed Brown ierr = PetscHeaderDestroy(c);CHKERRQ(ierr); 69af33a6ddSJed Brown PetscFunctionReturn(0); 70af33a6ddSJed Brown } 71af33a6ddSJed Brown 72af33a6ddSJed Brown PetscErrorCode CharacteristicCreate(MPI_Comm comm, Characteristic *c) 73af33a6ddSJed Brown { 74af33a6ddSJed Brown Characteristic newC; 75af33a6ddSJed Brown PetscErrorCode ierr; 76af33a6ddSJed Brown 77af33a6ddSJed Brown PetscFunctionBegin; 78af33a6ddSJed Brown PetscValidPointer(c, 2); 790298fd71SBarry Smith *c = NULL; 80607a6623SBarry Smith ierr = CharacteristicInitializePackage();CHKERRQ(ierr); 81af33a6ddSJed Brown 821b266c99SBarry Smith ierr = PetscHeaderCreate(newC, CHARACTERISTIC_CLASSID, "Characteristic", "Characteristic", "Characteristic", comm, CharacteristicDestroy, CharacteristicView);CHKERRQ(ierr); 83af33a6ddSJed Brown *c = newC; 84af33a6ddSJed Brown 85af33a6ddSJed Brown newC->structured = PETSC_TRUE; 86af33a6ddSJed Brown newC->numIds = 0; 870298fd71SBarry Smith newC->velocityDA = NULL; 880298fd71SBarry Smith newC->velocity = NULL; 890298fd71SBarry Smith newC->velocityOld = NULL; 90af33a6ddSJed Brown newC->numVelocityComp = 0; 910298fd71SBarry Smith newC->velocityComp = NULL; 920298fd71SBarry Smith newC->velocityInterp = NULL; 930298fd71SBarry Smith newC->velocityInterpLocal = NULL; 940298fd71SBarry Smith newC->velocityCtx = NULL; 950298fd71SBarry Smith newC->fieldDA = NULL; 960298fd71SBarry Smith newC->field = NULL; 97af33a6ddSJed Brown newC->numFieldComp = 0; 980298fd71SBarry Smith newC->fieldComp = NULL; 990298fd71SBarry Smith newC->fieldInterp = NULL; 1000298fd71SBarry Smith newC->fieldInterpLocal = NULL; 1010298fd71SBarry Smith newC->fieldCtx = NULL; 1020298fd71SBarry Smith newC->itemType = 0; 1030298fd71SBarry Smith newC->queue = NULL; 104af33a6ddSJed Brown newC->queueSize = 0; 105af33a6ddSJed Brown newC->queueMax = 0; 1060298fd71SBarry Smith newC->queueLocal = NULL; 107af33a6ddSJed Brown newC->queueLocalSize = 0; 108af33a6ddSJed Brown newC->queueLocalMax = 0; 1090298fd71SBarry Smith newC->queueRemote = NULL; 110af33a6ddSJed Brown newC->queueRemoteSize = 0; 111af33a6ddSJed Brown newC->queueRemoteMax = 0; 112af33a6ddSJed Brown newC->numNeighbors = 0; 1130298fd71SBarry Smith newC->neighbors = NULL; 1140298fd71SBarry Smith newC->needCount = NULL; 1150298fd71SBarry Smith newC->localOffsets = NULL; 1160298fd71SBarry Smith newC->fillCount = NULL; 1170298fd71SBarry Smith newC->remoteOffsets = NULL; 1180298fd71SBarry Smith newC->request = NULL; 1190298fd71SBarry Smith newC->status = NULL; 120af33a6ddSJed Brown PetscFunctionReturn(0); 121af33a6ddSJed Brown } 122af33a6ddSJed Brown 123af33a6ddSJed Brown /*@C 124af33a6ddSJed Brown CharacteristicSetType - Builds Characteristic for a particular solver. 125af33a6ddSJed Brown 126af33a6ddSJed Brown Logically Collective on Characteristic 127af33a6ddSJed Brown 128af33a6ddSJed Brown Input Parameters: 129af33a6ddSJed Brown + c - the method of characteristics context 130af33a6ddSJed Brown - type - a known method 131af33a6ddSJed Brown 132af33a6ddSJed Brown Options Database Key: 133af33a6ddSJed Brown . -characteristic_type <method> - Sets the method; use -help for a list 134af33a6ddSJed Brown of available methods 135af33a6ddSJed Brown 136af33a6ddSJed Brown Notes: 13730a76a96SBarry Smith See "include/petsccharacteristic.h" for available methods 138af33a6ddSJed Brown 139af33a6ddSJed Brown Normally, it is best to use the CharacteristicSetFromOptions() command and 140af33a6ddSJed Brown then set the Characteristic type from the options database rather than by using 141af33a6ddSJed Brown this routine. Using the options database provides the user with 142af33a6ddSJed Brown maximum flexibility in evaluating the many different Krylov methods. 143af33a6ddSJed Brown The CharacteristicSetType() routine is provided for those situations where it 144af33a6ddSJed Brown is necessary to set the iterative solver independently of the command 145af33a6ddSJed Brown line or options database. This might be the case, for example, when 146af33a6ddSJed Brown the choice of iterative solver changes during the execution of the 147af33a6ddSJed Brown program, and the user's application is taking responsibility for 148af33a6ddSJed Brown choosing the appropriate method. In other words, this routine is 149af33a6ddSJed Brown not for beginners. 150af33a6ddSJed Brown 151af33a6ddSJed Brown Level: intermediate 152af33a6ddSJed Brown 153af33a6ddSJed Brown .keywords: Characteristic, set, method 154af33a6ddSJed Brown 155af33a6ddSJed Brown .seealso: CharacteristicType 156af33a6ddSJed Brown 157af33a6ddSJed Brown @*/ 15819fd82e9SBarry Smith PetscErrorCode CharacteristicSetType(Characteristic c, CharacteristicType type) 159af33a6ddSJed Brown { 160af33a6ddSJed Brown PetscErrorCode ierr, (*r)(Characteristic); 161af33a6ddSJed Brown PetscBool match; 162af33a6ddSJed Brown 163af33a6ddSJed Brown PetscFunctionBegin; 164af33a6ddSJed Brown PetscValidHeaderSpecific(c, CHARACTERISTIC_CLASSID, 1); 165af33a6ddSJed Brown PetscValidCharPointer(type, 2); 166af33a6ddSJed Brown 167251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject) c, type, &match);CHKERRQ(ierr); 168af33a6ddSJed Brown if (match) PetscFunctionReturn(0); 169af33a6ddSJed Brown 170af33a6ddSJed Brown if (c->data) { 171af33a6ddSJed Brown /* destroy the old private Characteristic context */ 172af33a6ddSJed Brown ierr = (*c->ops->destroy)(c);CHKERRQ(ierr); 1730298fd71SBarry Smith c->ops->destroy = NULL; 174af33a6ddSJed Brown c->data = 0; 175af33a6ddSJed Brown } 176af33a6ddSJed Brown 1771c9cd337SJed Brown ierr = PetscFunctionListFind(CharacteristicList,type,&r);CHKERRQ(ierr); 178af33a6ddSJed Brown if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown Characteristic type given: %s", type); 179af33a6ddSJed Brown c->setupcalled = 0; 180af33a6ddSJed Brown ierr = (*r)(c);CHKERRQ(ierr); 181af33a6ddSJed Brown ierr = PetscObjectChangeTypeName((PetscObject) c, type);CHKERRQ(ierr); 182af33a6ddSJed Brown PetscFunctionReturn(0); 183af33a6ddSJed Brown } 184af33a6ddSJed Brown 185af33a6ddSJed Brown /*@ 186af33a6ddSJed Brown CharacteristicSetUp - Sets up the internal data structures for the 187af33a6ddSJed Brown later use of an iterative solver. 188af33a6ddSJed Brown 189af33a6ddSJed Brown Collective on Characteristic 190af33a6ddSJed Brown 191af33a6ddSJed Brown Input Parameter: 192af33a6ddSJed Brown . ksp - iterative context obtained from CharacteristicCreate() 193af33a6ddSJed Brown 194af33a6ddSJed Brown Level: developer 195af33a6ddSJed Brown 196af33a6ddSJed Brown .keywords: Characteristic, setup 197af33a6ddSJed Brown 198af33a6ddSJed Brown .seealso: CharacteristicCreate(), CharacteristicSolve(), CharacteristicDestroy() 199af33a6ddSJed Brown @*/ 200af33a6ddSJed Brown PetscErrorCode CharacteristicSetUp(Characteristic c) 201af33a6ddSJed Brown { 202af33a6ddSJed Brown PetscErrorCode ierr; 203af33a6ddSJed Brown 204af33a6ddSJed Brown PetscFunctionBegin; 205af33a6ddSJed Brown PetscValidHeaderSpecific(c, CHARACTERISTIC_CLASSID, 1); 206af33a6ddSJed Brown 207af33a6ddSJed Brown if (!((PetscObject)c)->type_name) { 208af33a6ddSJed Brown ierr = CharacteristicSetType(c, CHARACTERISTICDA);CHKERRQ(ierr); 209af33a6ddSJed Brown } 210af33a6ddSJed Brown 211af33a6ddSJed Brown if (c->setupcalled == 2) PetscFunctionReturn(0); 212af33a6ddSJed Brown 213af33a6ddSJed Brown ierr = PetscLogEventBegin(CHARACTERISTIC_SetUp,c,0,0,0);CHKERRQ(ierr); 214af33a6ddSJed Brown if (!c->setupcalled) { 215af33a6ddSJed Brown ierr = (*c->ops->setup)(c);CHKERRQ(ierr); 216af33a6ddSJed Brown } 217af33a6ddSJed Brown ierr = PetscLogEventEnd(CHARACTERISTIC_SetUp,c,0,0,0);CHKERRQ(ierr); 218af33a6ddSJed Brown c->setupcalled = 2; 219af33a6ddSJed Brown PetscFunctionReturn(0); 220af33a6ddSJed Brown } 221af33a6ddSJed Brown 222af33a6ddSJed Brown /*@C 2231c84c290SBarry Smith CharacteristicRegister - Adds a solver to the method of characteristics package. 2241c84c290SBarry Smith 2251c84c290SBarry Smith Not Collective 2261c84c290SBarry Smith 2271c84c290SBarry Smith Input Parameters: 2281c84c290SBarry Smith + name_solver - name of a new user-defined solver 2291c84c290SBarry Smith - routine_create - routine to create method context 2301c84c290SBarry Smith 2311c84c290SBarry Smith Sample usage: 2321c84c290SBarry Smith .vb 233bdf89e91SBarry Smith CharacteristicRegister("my_char", MyCharCreate); 2341c84c290SBarry Smith .ve 2351c84c290SBarry Smith 2361c84c290SBarry Smith Then, your Characteristic type can be chosen with the procedural interface via 2371c84c290SBarry Smith .vb 2381c84c290SBarry Smith CharacteristicCreate(MPI_Comm, Characteristic* &char); 2391c84c290SBarry Smith CharacteristicSetType(char,"my_char"); 2401c84c290SBarry Smith .ve 2411c84c290SBarry Smith or at runtime via the option 2421c84c290SBarry Smith .vb 2431c84c290SBarry Smith -characteristic_type my_char 2441c84c290SBarry Smith .ve 2451c84c290SBarry Smith 2461c84c290SBarry Smith Notes: 2471c84c290SBarry Smith CharacteristicRegister() may be called multiple times to add several user-defined solvers. 2481c84c290SBarry Smith 2491c84c290SBarry Smith .keywords: Characteristic, register 2501c84c290SBarry Smith 2511c84c290SBarry Smith .seealso: CharacteristicRegisterAll(), CharacteristicRegisterDestroy() 252af33a6ddSJed Brown 253af33a6ddSJed Brown Level: advanced 254af33a6ddSJed Brown @*/ 255bdf89e91SBarry Smith PetscErrorCode CharacteristicRegister(const char sname[],PetscErrorCode (*function)(Characteristic)) 256af33a6ddSJed Brown { 257af33a6ddSJed Brown PetscErrorCode ierr; 258af33a6ddSJed Brown 259af33a6ddSJed Brown PetscFunctionBegin; 260a240a19fSJed Brown ierr = PetscFunctionListAdd(&CharacteristicList,sname,function);CHKERRQ(ierr); 261af33a6ddSJed Brown PetscFunctionReturn(0); 262af33a6ddSJed Brown } 263af33a6ddSJed Brown 264af33a6ddSJed Brown PetscErrorCode CharacteristicSetVelocityInterpolation(Characteristic c, DM da, Vec v, Vec vOld, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(Vec, PetscReal[], PetscInt, PetscInt[], PetscScalar[], void*), void *ctx) 265af33a6ddSJed Brown { 266af33a6ddSJed Brown PetscFunctionBegin; 267af33a6ddSJed Brown c->velocityDA = da; 268af33a6ddSJed Brown c->velocity = v; 269af33a6ddSJed Brown c->velocityOld = vOld; 270af33a6ddSJed Brown c->numVelocityComp = numComponents; 271af33a6ddSJed Brown c->velocityComp = components; 272af33a6ddSJed Brown c->velocityInterp = interp; 273af33a6ddSJed Brown c->velocityCtx = ctx; 274af33a6ddSJed Brown PetscFunctionReturn(0); 275af33a6ddSJed Brown } 276af33a6ddSJed Brown 277af33a6ddSJed Brown PetscErrorCode CharacteristicSetVelocityInterpolationLocal(Characteristic c, DM da, Vec v, Vec vOld, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(void*, PetscReal [], PetscInt, PetscInt[], PetscScalar[], void*), void *ctx) 278af33a6ddSJed Brown { 279af33a6ddSJed Brown PetscFunctionBegin; 280af33a6ddSJed Brown c->velocityDA = da; 281af33a6ddSJed Brown c->velocity = v; 282af33a6ddSJed Brown c->velocityOld = vOld; 283af33a6ddSJed Brown c->numVelocityComp = numComponents; 284af33a6ddSJed Brown c->velocityComp = components; 285af33a6ddSJed Brown c->velocityInterpLocal = interp; 286af33a6ddSJed Brown c->velocityCtx = ctx; 287af33a6ddSJed Brown PetscFunctionReturn(0); 288af33a6ddSJed Brown } 289af33a6ddSJed Brown 290af33a6ddSJed Brown PetscErrorCode CharacteristicSetFieldInterpolation(Characteristic c, DM da, Vec v, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(Vec, PetscReal[], PetscInt, PetscInt[], PetscScalar[], void*), void *ctx) 291af33a6ddSJed Brown { 292af33a6ddSJed Brown PetscFunctionBegin; 293af33a6ddSJed Brown #if 0 294af33a6ddSJed Brown 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."); 295af33a6ddSJed Brown #endif 296af33a6ddSJed Brown c->fieldDA = da; 297af33a6ddSJed Brown c->field = v; 298af33a6ddSJed Brown c->numFieldComp = numComponents; 299af33a6ddSJed Brown c->fieldComp = components; 300af33a6ddSJed Brown c->fieldInterp = interp; 301af33a6ddSJed Brown c->fieldCtx = ctx; 302af33a6ddSJed Brown PetscFunctionReturn(0); 303af33a6ddSJed Brown } 304af33a6ddSJed Brown 305af33a6ddSJed Brown PetscErrorCode CharacteristicSetFieldInterpolationLocal(Characteristic c, DM da, Vec v, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(void*, PetscReal[], PetscInt, PetscInt[], PetscScalar [], void*), void *ctx) 306af33a6ddSJed Brown { 307af33a6ddSJed Brown PetscFunctionBegin; 308af33a6ddSJed Brown #if 0 309af33a6ddSJed Brown 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."); 310af33a6ddSJed Brown #endif 311af33a6ddSJed Brown c->fieldDA = da; 312af33a6ddSJed Brown c->field = v; 313af33a6ddSJed Brown c->numFieldComp = numComponents; 314af33a6ddSJed Brown c->fieldComp = components; 315af33a6ddSJed Brown c->fieldInterpLocal = interp; 316af33a6ddSJed Brown c->fieldCtx = ctx; 317af33a6ddSJed Brown PetscFunctionReturn(0); 318af33a6ddSJed Brown } 319af33a6ddSJed Brown 320af33a6ddSJed Brown PetscErrorCode CharacteristicSolve(Characteristic c, PetscReal dt, Vec solution) 321af33a6ddSJed Brown { 322af33a6ddSJed Brown CharacteristicPointDA2D Qi; 323af33a6ddSJed Brown DM da = c->velocityDA; 324af33a6ddSJed Brown Vec velocityLocal, velocityLocalOld; 325af33a6ddSJed Brown Vec fieldLocal; 326af33a6ddSJed Brown DMDALocalInfo info; 327af33a6ddSJed Brown PetscScalar **solArray; 328af33a6ddSJed Brown void *velocityArray; 329af33a6ddSJed Brown void *velocityArrayOld; 330af33a6ddSJed Brown void *fieldArray; 3311cc8b206SBarry Smith PetscScalar *interpIndices; 3321cc8b206SBarry Smith PetscScalar *velocityValues, *velocityValuesOld; 3331cc8b206SBarry Smith PetscScalar *fieldValues; 334af33a6ddSJed Brown PetscMPIInt rank; 335af33a6ddSJed Brown PetscInt dim; 336af33a6ddSJed Brown PetscMPIInt neighbors[9]; 337af33a6ddSJed Brown PetscInt dof; 338af33a6ddSJed Brown PetscInt gx, gy; 339af33a6ddSJed Brown PetscInt n, is, ie, js, je, comp; 340af33a6ddSJed Brown PetscErrorCode ierr; 341af33a6ddSJed Brown 342af33a6ddSJed Brown PetscFunctionBegin; 343af33a6ddSJed Brown c->queueSize = 0; 344ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank);CHKERRQ(ierr); 345af33a6ddSJed Brown ierr = DMDAGetNeighborsRank(da, neighbors);CHKERRQ(ierr); 346af33a6ddSJed Brown ierr = CharacteristicSetNeighbors(c, 9, neighbors);CHKERRQ(ierr); 347af33a6ddSJed Brown ierr = CharacteristicSetUp(c);CHKERRQ(ierr); 348af33a6ddSJed Brown /* global and local grid info */ 349af33a6ddSJed Brown ierr = DMDAGetInfo(da, &dim, &gx, &gy, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0);CHKERRQ(ierr); 350af33a6ddSJed Brown ierr = DMDAGetLocalInfo(da, &info);CHKERRQ(ierr); 351af33a6ddSJed Brown is = info.xs; ie = info.xs+info.xm; 352af33a6ddSJed Brown js = info.ys; je = info.ys+info.ym; 353af33a6ddSJed Brown /* Allocation */ 354785e854fSJed Brown ierr = PetscMalloc1(dim, &interpIndices);CHKERRQ(ierr); 355785e854fSJed Brown ierr = PetscMalloc1(c->numVelocityComp, &velocityValues);CHKERRQ(ierr); 356785e854fSJed Brown ierr = PetscMalloc1(c->numVelocityComp, &velocityValuesOld);CHKERRQ(ierr); 357785e854fSJed Brown ierr = PetscMalloc1(c->numFieldComp, &fieldValues);CHKERRQ(ierr); 358af33a6ddSJed Brown ierr = PetscLogEventBegin(CHARACTERISTIC_Solve,0,0,0,0);CHKERRQ(ierr); 359af33a6ddSJed Brown 360af33a6ddSJed Brown /* ----------------------------------------------------------------------- 361af33a6ddSJed Brown PART 1, AT t-dt/2 362af33a6ddSJed Brown -----------------------------------------------------------------------*/ 363af33a6ddSJed Brown ierr = PetscLogEventBegin(CHARACTERISTIC_QueueSetup,0,0,0,0);CHKERRQ(ierr); 364af33a6ddSJed Brown /* GET POSITION AT HALF TIME IN THE PAST */ 365af33a6ddSJed Brown if (c->velocityInterpLocal) { 366af33a6ddSJed Brown ierr = DMGetLocalVector(c->velocityDA, &velocityLocal);CHKERRQ(ierr); 367af33a6ddSJed Brown ierr = DMGetLocalVector(c->velocityDA, &velocityLocalOld);CHKERRQ(ierr); 368af33a6ddSJed Brown ierr = DMGlobalToLocalBegin(c->velocityDA, c->velocity, INSERT_VALUES, velocityLocal);CHKERRQ(ierr); 369af33a6ddSJed Brown ierr = DMGlobalToLocalEnd(c->velocityDA, c->velocity, INSERT_VALUES, velocityLocal);CHKERRQ(ierr); 370af33a6ddSJed Brown ierr = DMGlobalToLocalBegin(c->velocityDA, c->velocityOld, INSERT_VALUES, velocityLocalOld);CHKERRQ(ierr); 371af33a6ddSJed Brown ierr = DMGlobalToLocalEnd(c->velocityDA, c->velocityOld, INSERT_VALUES, velocityLocalOld);CHKERRQ(ierr); 372af33a6ddSJed Brown ierr = DMDAVecGetArray(c->velocityDA, velocityLocal, &velocityArray);CHKERRQ(ierr); 373af33a6ddSJed Brown ierr = DMDAVecGetArray(c->velocityDA, velocityLocalOld, &velocityArrayOld);CHKERRQ(ierr); 374af33a6ddSJed Brown } 3750298fd71SBarry Smith ierr = PetscInfo(NULL, "Calculating position at t_{n - 1/2}\n");CHKERRQ(ierr); 376af33a6ddSJed Brown for (Qi.j = js; Qi.j < je; Qi.j++) { 377af33a6ddSJed Brown for (Qi.i = is; Qi.i < ie; Qi.i++) { 378af33a6ddSJed Brown interpIndices[0] = Qi.i; 379af33a6ddSJed Brown interpIndices[1] = Qi.j; 38061ba4676SBarry Smith if (c->velocityInterpLocal) {ierr = c->velocityInterpLocal(velocityArray, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx);CHKERRQ(ierr);} 38161ba4676SBarry Smith else {ierr = c->velocityInterp(c->velocity, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx);CHKERRQ(ierr);} 382af33a6ddSJed Brown Qi.x = Qi.i - velocityValues[0]*dt/2.0; 383af33a6ddSJed Brown Qi.y = Qi.j - velocityValues[1]*dt/2.0; 384af33a6ddSJed Brown 385af33a6ddSJed Brown /* Determine whether the position at t - dt/2 is local */ 386af33a6ddSJed Brown Qi.proc = DMDAGetNeighborRelative(da, Qi.x, Qi.y); 387af33a6ddSJed Brown 388af33a6ddSJed Brown /* Check for Periodic boundaries and move all periodic points back onto the domain */ 389af33a6ddSJed Brown ierr = DMDAMapCoordsToPeriodicDomain(da,&(Qi.x),&(Qi.y));CHKERRQ(ierr); 390af33a6ddSJed Brown ierr = CharacteristicAddPoint(c, &Qi);CHKERRQ(ierr); 391af33a6ddSJed Brown } 392af33a6ddSJed Brown } 393af33a6ddSJed Brown ierr = PetscLogEventEnd(CHARACTERISTIC_QueueSetup,0,0,0,0);CHKERRQ(ierr); 394af33a6ddSJed Brown 395af33a6ddSJed Brown ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr); 396af33a6ddSJed Brown ierr = CharacteristicSendCoordinatesBegin(c);CHKERRQ(ierr); 397af33a6ddSJed Brown ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr); 398af33a6ddSJed Brown 399af33a6ddSJed Brown ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeLocal,0,0,0,0);CHKERRQ(ierr); 400af33a6ddSJed Brown /* Calculate velocity at t_n+1/2 (local values) */ 4010298fd71SBarry Smith ierr = PetscInfo(NULL, "Calculating local velocities at t_{n - 1/2}\n");CHKERRQ(ierr); 402af33a6ddSJed Brown for (n = 0; n < c->queueSize; n++) { 403af33a6ddSJed Brown Qi = c->queue[n]; 404af33a6ddSJed Brown if (c->neighbors[Qi.proc] == rank) { 405af33a6ddSJed Brown interpIndices[0] = Qi.x; 406af33a6ddSJed Brown interpIndices[1] = Qi.y; 407af33a6ddSJed Brown if (c->velocityInterpLocal) { 40861ba4676SBarry Smith ierr = c->velocityInterpLocal(velocityArray, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx);CHKERRQ(ierr); 40961ba4676SBarry Smith ierr = c->velocityInterpLocal(velocityArrayOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx);CHKERRQ(ierr); 410af33a6ddSJed Brown } else { 41161ba4676SBarry Smith ierr = c->velocityInterp(c->velocity, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx);CHKERRQ(ierr); 41261528463SBarry Smith ierr = c->velocityInterp(c->velocityOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx);CHKERRQ(ierr); 413af33a6ddSJed Brown } 414af33a6ddSJed Brown Qi.x = 0.5*(velocityValues[0] + velocityValuesOld[0]); 415af33a6ddSJed Brown Qi.y = 0.5*(velocityValues[1] + velocityValuesOld[1]); 416af33a6ddSJed Brown } 417af33a6ddSJed Brown c->queue[n] = Qi; 418af33a6ddSJed Brown } 419af33a6ddSJed Brown ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeLocal,0,0,0,0);CHKERRQ(ierr); 420af33a6ddSJed Brown 421af33a6ddSJed Brown ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr); 422af33a6ddSJed Brown ierr = CharacteristicSendCoordinatesEnd(c);CHKERRQ(ierr); 423af33a6ddSJed Brown ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr); 424af33a6ddSJed Brown 425af33a6ddSJed Brown 426af33a6ddSJed Brown /* Calculate velocity at t_n+1/2 (fill remote requests) */ 427af33a6ddSJed Brown ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeRemote,0,0,0,0);CHKERRQ(ierr); 4280298fd71SBarry Smith ierr = PetscInfo1(NULL, "Calculating %d remote velocities at t_{n - 1/2}\n", c->queueRemoteSize);CHKERRQ(ierr); 429af33a6ddSJed Brown for (n = 0; n < c->queueRemoteSize; n++) { 430af33a6ddSJed Brown Qi = c->queueRemote[n]; 431af33a6ddSJed Brown interpIndices[0] = Qi.x; 432af33a6ddSJed Brown interpIndices[1] = Qi.y; 433af33a6ddSJed Brown if (c->velocityInterpLocal) { 43461ba4676SBarry Smith ierr = c->velocityInterpLocal(velocityArray, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx);CHKERRQ(ierr); 43561ba4676SBarry Smith ierr = c->velocityInterpLocal(velocityArrayOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx);CHKERRQ(ierr); 436af33a6ddSJed Brown } else { 43761ba4676SBarry Smith ierr = c->velocityInterp(c->velocity, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx);CHKERRQ(ierr); 43861ba4676SBarry Smith ierr = c->velocityInterp(c->velocityOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx);CHKERRQ(ierr); 439af33a6ddSJed Brown } 440af33a6ddSJed Brown Qi.x = 0.5*(velocityValues[0] + velocityValuesOld[0]); 441af33a6ddSJed Brown Qi.y = 0.5*(velocityValues[1] + velocityValuesOld[1]); 442af33a6ddSJed Brown c->queueRemote[n] = Qi; 443af33a6ddSJed Brown } 444af33a6ddSJed Brown ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeRemote,0,0,0,0);CHKERRQ(ierr); 445af33a6ddSJed Brown ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr); 446af33a6ddSJed Brown ierr = CharacteristicGetValuesBegin(c);CHKERRQ(ierr); 447af33a6ddSJed Brown ierr = CharacteristicGetValuesEnd(c);CHKERRQ(ierr); 448af33a6ddSJed Brown if (c->velocityInterpLocal) { 449af33a6ddSJed Brown ierr = DMDAVecRestoreArray(c->velocityDA, velocityLocal, &velocityArray);CHKERRQ(ierr); 450af33a6ddSJed Brown ierr = DMDAVecRestoreArray(c->velocityDA, velocityLocalOld, &velocityArrayOld);CHKERRQ(ierr); 451af33a6ddSJed Brown ierr = DMRestoreLocalVector(c->velocityDA, &velocityLocal);CHKERRQ(ierr); 452af33a6ddSJed Brown ierr = DMRestoreLocalVector(c->velocityDA, &velocityLocalOld);CHKERRQ(ierr); 453af33a6ddSJed Brown } 454af33a6ddSJed Brown ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr); 455af33a6ddSJed Brown 456af33a6ddSJed Brown /* ----------------------------------------------------------------------- 457af33a6ddSJed Brown PART 2, AT t-dt 458af33a6ddSJed Brown -----------------------------------------------------------------------*/ 459af33a6ddSJed Brown 460af33a6ddSJed Brown /* GET POSITION AT t_n (local values) */ 461af33a6ddSJed Brown ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeLocal,0,0,0,0);CHKERRQ(ierr); 4620298fd71SBarry Smith ierr = PetscInfo(NULL, "Calculating position at t_{n}\n");CHKERRQ(ierr); 463af33a6ddSJed Brown for (n = 0; n < c->queueSize; n++) { 464af33a6ddSJed Brown Qi = c->queue[n]; 465af33a6ddSJed Brown Qi.x = Qi.i - Qi.x*dt; 466af33a6ddSJed Brown Qi.y = Qi.j - Qi.y*dt; 467af33a6ddSJed Brown 468af33a6ddSJed Brown /* Determine whether the position at t-dt is local */ 469af33a6ddSJed Brown Qi.proc = DMDAGetNeighborRelative(da, Qi.x, Qi.y); 470af33a6ddSJed Brown 471af33a6ddSJed Brown /* Check for Periodic boundaries and move all periodic points back onto the domain */ 472af33a6ddSJed Brown ierr = DMDAMapCoordsToPeriodicDomain(da,&(Qi.x),&(Qi.y));CHKERRQ(ierr); 473af33a6ddSJed Brown 474af33a6ddSJed Brown c->queue[n] = Qi; 475af33a6ddSJed Brown } 476af33a6ddSJed Brown ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeLocal,0,0,0,0);CHKERRQ(ierr); 477af33a6ddSJed Brown 478af33a6ddSJed Brown ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr); 479af33a6ddSJed Brown ierr = CharacteristicSendCoordinatesBegin(c);CHKERRQ(ierr); 480af33a6ddSJed Brown ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr); 481af33a6ddSJed Brown 482af33a6ddSJed Brown /* GET VALUE AT FULL TIME IN THE PAST (LOCAL REQUESTS) */ 483af33a6ddSJed Brown ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeLocal,0,0,0,0);CHKERRQ(ierr); 484af33a6ddSJed Brown if (c->fieldInterpLocal) { 485af33a6ddSJed Brown ierr = DMGetLocalVector(c->fieldDA, &fieldLocal);CHKERRQ(ierr); 486af33a6ddSJed Brown ierr = DMGlobalToLocalBegin(c->fieldDA, c->field, INSERT_VALUES, fieldLocal);CHKERRQ(ierr); 487af33a6ddSJed Brown ierr = DMGlobalToLocalEnd(c->fieldDA, c->field, INSERT_VALUES, fieldLocal);CHKERRQ(ierr); 488af33a6ddSJed Brown ierr = DMDAVecGetArray(c->fieldDA, fieldLocal, &fieldArray);CHKERRQ(ierr); 489af33a6ddSJed Brown } 4900298fd71SBarry Smith ierr = PetscInfo(NULL, "Calculating local field at t_{n}\n");CHKERRQ(ierr); 491af33a6ddSJed Brown for (n = 0; n < c->queueSize; n++) { 492af33a6ddSJed Brown if (c->neighbors[c->queue[n].proc] == rank) { 493af33a6ddSJed Brown interpIndices[0] = c->queue[n].x; 494af33a6ddSJed Brown interpIndices[1] = c->queue[n].y; 49561ba4676SBarry Smith if (c->fieldInterpLocal) {ierr = c->fieldInterpLocal(fieldArray, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx);CHKERRQ(ierr);} 49661ba4676SBarry Smith else {ierr = c->fieldInterp(c->field, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx);CHKERRQ(ierr);} 497bbd56ea5SKarl Rupp for (comp = 0; comp < c->numFieldComp; comp++) c->queue[n].field[comp] = fieldValues[comp]; 498af33a6ddSJed Brown } 499af33a6ddSJed Brown } 500af33a6ddSJed Brown ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeLocal,0,0,0,0);CHKERRQ(ierr); 501af33a6ddSJed Brown 502af33a6ddSJed Brown ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr); 503af33a6ddSJed Brown ierr = CharacteristicSendCoordinatesEnd(c);CHKERRQ(ierr); 504af33a6ddSJed Brown ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr); 505af33a6ddSJed Brown 506af33a6ddSJed Brown /* GET VALUE AT FULL TIME IN THE PAST (REMOTE REQUESTS) */ 507af33a6ddSJed Brown ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeRemote,0,0,0,0);CHKERRQ(ierr); 5080298fd71SBarry Smith ierr = PetscInfo1(NULL, "Calculating %d remote field points at t_{n}\n", c->queueRemoteSize);CHKERRQ(ierr); 509af33a6ddSJed Brown for (n = 0; n < c->queueRemoteSize; n++) { 510af33a6ddSJed Brown interpIndices[0] = c->queueRemote[n].x; 511af33a6ddSJed Brown interpIndices[1] = c->queueRemote[n].y; 512af33a6ddSJed Brown 513af33a6ddSJed Brown /* for debugging purposes */ 514af33a6ddSJed Brown if (1) { /* hacked bounds test...let's do better */ 515af33a6ddSJed Brown PetscScalar im = interpIndices[0]; PetscScalar jm = interpIndices[1]; 516af33a6ddSJed Brown 517bbd56ea5SKarl Rupp 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); 518af33a6ddSJed Brown } 519af33a6ddSJed Brown 52061ba4676SBarry Smith if (c->fieldInterpLocal) {ierr = c->fieldInterpLocal(fieldArray, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx);CHKERRQ(ierr);} 52161528463SBarry Smith else {ierr = c->fieldInterp(c->field, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx);CHKERRQ(ierr);} 522bbd56ea5SKarl Rupp for (comp = 0; comp < c->numFieldComp; comp++) c->queueRemote[n].field[comp] = fieldValues[comp]; 523af33a6ddSJed Brown } 524af33a6ddSJed Brown ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeRemote,0,0,0,0);CHKERRQ(ierr); 525af33a6ddSJed Brown 526af33a6ddSJed Brown ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr); 527af33a6ddSJed Brown ierr = CharacteristicGetValuesBegin(c);CHKERRQ(ierr); 528af33a6ddSJed Brown ierr = CharacteristicGetValuesEnd(c);CHKERRQ(ierr); 529af33a6ddSJed Brown if (c->fieldInterpLocal) { 530af33a6ddSJed Brown ierr = DMDAVecRestoreArray(c->fieldDA, fieldLocal, &fieldArray);CHKERRQ(ierr); 531af33a6ddSJed Brown ierr = DMRestoreLocalVector(c->fieldDA, &fieldLocal);CHKERRQ(ierr); 532af33a6ddSJed Brown } 533af33a6ddSJed Brown ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr); 534af33a6ddSJed Brown 535af33a6ddSJed Brown /* Return field of characteristics at t_n-1 */ 536af33a6ddSJed Brown ierr = PetscLogEventBegin(CHARACTERISTIC_DAUpdate,0,0,0,0);CHKERRQ(ierr); 537af33a6ddSJed Brown ierr = DMDAGetInfo(c->fieldDA,0,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); 538af33a6ddSJed Brown ierr = DMDAVecGetArray(c->fieldDA, solution, &solArray);CHKERRQ(ierr); 539af33a6ddSJed Brown for (n = 0; n < c->queueSize; n++) { 540af33a6ddSJed Brown Qi = c->queue[n]; 541bbd56ea5SKarl Rupp for (comp = 0; comp < c->numFieldComp; comp++) solArray[Qi.j][Qi.i*dof+c->fieldComp[comp]] = Qi.field[comp]; 542af33a6ddSJed Brown } 543af33a6ddSJed Brown ierr = DMDAVecRestoreArray(c->fieldDA, solution, &solArray);CHKERRQ(ierr); 544af33a6ddSJed Brown ierr = PetscLogEventEnd(CHARACTERISTIC_DAUpdate,0,0,0,0);CHKERRQ(ierr); 545af33a6ddSJed Brown ierr = PetscLogEventEnd(CHARACTERISTIC_Solve,0,0,0,0);CHKERRQ(ierr); 546af33a6ddSJed Brown 547af33a6ddSJed Brown /* Cleanup */ 548af33a6ddSJed Brown ierr = PetscFree(interpIndices);CHKERRQ(ierr); 549af33a6ddSJed Brown ierr = PetscFree(velocityValues);CHKERRQ(ierr); 550af33a6ddSJed Brown ierr = PetscFree(velocityValuesOld);CHKERRQ(ierr); 551af33a6ddSJed Brown ierr = PetscFree(fieldValues);CHKERRQ(ierr); 552af33a6ddSJed Brown PetscFunctionReturn(0); 553af33a6ddSJed Brown } 554af33a6ddSJed Brown 555af33a6ddSJed Brown PetscErrorCode CharacteristicSetNeighbors(Characteristic c, PetscInt numNeighbors, PetscMPIInt neighbors[]) 556af33a6ddSJed Brown { 557af33a6ddSJed Brown PetscErrorCode ierr; 558af33a6ddSJed Brown 559af33a6ddSJed Brown PetscFunctionBegin; 560af33a6ddSJed Brown c->numNeighbors = numNeighbors; 5616eaa7dc8SBarry Smith ierr = PetscFree(c->neighbors);CHKERRQ(ierr); 562785e854fSJed Brown ierr = PetscMalloc1(numNeighbors, &c->neighbors);CHKERRQ(ierr); 563af33a6ddSJed Brown ierr = PetscMemcpy(c->neighbors, neighbors, numNeighbors * sizeof(PetscMPIInt));CHKERRQ(ierr); 564af33a6ddSJed Brown PetscFunctionReturn(0); 565af33a6ddSJed Brown } 566af33a6ddSJed Brown 567af33a6ddSJed Brown PetscErrorCode CharacteristicAddPoint(Characteristic c, CharacteristicPointDA2D *point) 568af33a6ddSJed Brown { 569af33a6ddSJed Brown PetscFunctionBegin; 570f23aa3ddSBarry Smith if (c->queueSize >= c->queueMax) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Exceeeded maximum queue size %d", c->queueMax); 571af33a6ddSJed Brown c->queue[c->queueSize++] = *point; 572af33a6ddSJed Brown PetscFunctionReturn(0); 573af33a6ddSJed Brown } 574af33a6ddSJed Brown 575af33a6ddSJed Brown int CharacteristicSendCoordinatesBegin(Characteristic c) 576af33a6ddSJed Brown { 577af33a6ddSJed Brown PetscMPIInt rank, tag = 121; 578af33a6ddSJed Brown PetscInt i, n; 579af33a6ddSJed Brown PetscErrorCode ierr; 580af33a6ddSJed Brown 581af33a6ddSJed Brown PetscFunctionBegin; 582ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank);CHKERRQ(ierr); 5830b8f38c8SBarry Smith ierr = CharacteristicHeapSort(c, c->queue, c->queueSize);CHKERRQ(ierr); 584af33a6ddSJed Brown ierr = PetscMemzero(c->needCount, c->numNeighbors * sizeof(PetscInt));CHKERRQ(ierr); 585bbd56ea5SKarl Rupp for (i = 0; i < c->queueSize; i++) c->needCount[c->queue[i].proc]++; 586af33a6ddSJed Brown c->fillCount[0] = 0; 587af33a6ddSJed Brown for (n = 1; n < c->numNeighbors; n++) { 588ce94432eSBarry Smith ierr = MPI_Irecv(&(c->fillCount[n]), 1, MPIU_INT, c->neighbors[n], tag, PetscObjectComm((PetscObject)c), &(c->request[n-1]));CHKERRQ(ierr); 589af33a6ddSJed Brown } 590af33a6ddSJed Brown for (n = 1; n < c->numNeighbors; n++) { 591ce94432eSBarry Smith ierr = MPI_Send(&(c->needCount[n]), 1, MPIU_INT, c->neighbors[n], tag, PetscObjectComm((PetscObject)c));CHKERRQ(ierr); 592af33a6ddSJed Brown } 593af33a6ddSJed Brown ierr = MPI_Waitall(c->numNeighbors-1, c->request, c->status);CHKERRQ(ierr); 594af33a6ddSJed Brown /* Initialize the remote queue */ 595af33a6ddSJed Brown c->queueLocalMax = c->localOffsets[0] = 0; 596af33a6ddSJed Brown c->queueRemoteMax = c->remoteOffsets[0] = 0; 597af33a6ddSJed Brown for (n = 1; n < c->numNeighbors; n++) { 598af33a6ddSJed Brown c->remoteOffsets[n] = c->queueRemoteMax; 599af33a6ddSJed Brown c->queueRemoteMax += c->fillCount[n]; 600af33a6ddSJed Brown c->localOffsets[n] = c->queueLocalMax; 601af33a6ddSJed Brown c->queueLocalMax += c->needCount[n]; 602af33a6ddSJed Brown } 603af33a6ddSJed Brown /* HACK BEGIN */ 604bbd56ea5SKarl Rupp for (n = 1; n < c->numNeighbors; n++) c->localOffsets[n] += c->needCount[0]; 605af33a6ddSJed Brown c->needCount[0] = 0; 606af33a6ddSJed Brown /* HACK END */ 607af33a6ddSJed Brown if (c->queueRemoteMax) { 60895dccacaSBarry Smith ierr = PetscMalloc1(c->queueRemoteMax, &c->queueRemote);CHKERRQ(ierr); 6090298fd71SBarry Smith } else c->queueRemote = NULL; 610af33a6ddSJed Brown c->queueRemoteSize = c->queueRemoteMax; 611af33a6ddSJed Brown 612af33a6ddSJed Brown /* Send and Receive requests for values at t_n+1/2, giving the coordinates for interpolation */ 613af33a6ddSJed Brown for (n = 1; n < c->numNeighbors; n++) { 6140298fd71SBarry Smith ierr = PetscInfo2(NULL, "Receiving %d requests for values from proc %d\n", c->fillCount[n], c->neighbors[n]);CHKERRQ(ierr); 615ce94432eSBarry Smith 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); 616af33a6ddSJed Brown } 617af33a6ddSJed Brown for (n = 1; n < c->numNeighbors; n++) { 6180298fd71SBarry Smith ierr = PetscInfo2(NULL, "Sending %d requests for values from proc %d\n", c->needCount[n], c->neighbors[n]);CHKERRQ(ierr); 619ce94432eSBarry Smith ierr = MPI_Send(&(c->queue[c->localOffsets[n]]), c->needCount[n], c->itemType, c->neighbors[n], tag, PetscObjectComm((PetscObject)c));CHKERRQ(ierr); 620af33a6ddSJed Brown } 621af33a6ddSJed Brown PetscFunctionReturn(0); 622af33a6ddSJed Brown } 623af33a6ddSJed Brown 624af33a6ddSJed Brown PetscErrorCode CharacteristicSendCoordinatesEnd(Characteristic c) 625af33a6ddSJed Brown { 626af33a6ddSJed Brown #if 0 627af33a6ddSJed Brown PetscMPIInt rank; 628af33a6ddSJed Brown PetscInt n; 629af33a6ddSJed Brown #endif 630af33a6ddSJed Brown PetscErrorCode ierr; 631af33a6ddSJed Brown 632af33a6ddSJed Brown PetscFunctionBegin; 633af33a6ddSJed Brown ierr = MPI_Waitall(c->numNeighbors-1, c->request, c->status);CHKERRQ(ierr); 634af33a6ddSJed Brown #if 0 635ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank);CHKERRQ(ierr); 636af33a6ddSJed Brown for (n = 0; n < c->queueRemoteSize; n++) { 637f23aa3ddSBarry Smith 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); 638af33a6ddSJed Brown } 639af33a6ddSJed Brown #endif 640af33a6ddSJed Brown PetscFunctionReturn(0); 641af33a6ddSJed Brown } 642af33a6ddSJed Brown 643af33a6ddSJed Brown PetscErrorCode CharacteristicGetValuesBegin(Characteristic c) 644af33a6ddSJed Brown { 645af33a6ddSJed Brown PetscMPIInt tag = 121; 646af33a6ddSJed Brown PetscInt n; 647af33a6ddSJed Brown PetscErrorCode ierr; 648af33a6ddSJed Brown 649af33a6ddSJed Brown PetscFunctionBegin; 650af33a6ddSJed Brown /* SEND AND RECIEVE FILLED REQUESTS for velocities at t_n+1/2 */ 651af33a6ddSJed Brown for (n = 1; n < c->numNeighbors; n++) { 652ce94432eSBarry Smith 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); 653af33a6ddSJed Brown } 654af33a6ddSJed Brown for (n = 1; n < c->numNeighbors; n++) { 655ce94432eSBarry Smith ierr = MPI_Send(&(c->queueRemote[c->remoteOffsets[n]]), c->fillCount[n], c->itemType, c->neighbors[n], tag, PetscObjectComm((PetscObject)c));CHKERRQ(ierr); 656af33a6ddSJed Brown } 657af33a6ddSJed Brown PetscFunctionReturn(0); 658af33a6ddSJed Brown } 659af33a6ddSJed Brown 660af33a6ddSJed Brown PetscErrorCode CharacteristicGetValuesEnd(Characteristic c) 661af33a6ddSJed Brown { 662af33a6ddSJed Brown PetscErrorCode ierr; 663af33a6ddSJed Brown 664af33a6ddSJed Brown PetscFunctionBegin; 665af33a6ddSJed Brown ierr = MPI_Waitall(c->numNeighbors-1, c->request, c->status);CHKERRQ(ierr); 666af33a6ddSJed Brown /* Free queue of requests from other procs */ 667af33a6ddSJed Brown ierr = PetscFree(c->queueRemote);CHKERRQ(ierr); 668af33a6ddSJed Brown PetscFunctionReturn(0); 669af33a6ddSJed Brown } 670af33a6ddSJed Brown 671af33a6ddSJed Brown /*---------------------------------------------------------------------*/ 672af33a6ddSJed Brown /* 673af33a6ddSJed Brown Based on code from http://linux.wku.edu/~lamonml/algor/sort/heap.html 674af33a6ddSJed Brown */ 6750b8f38c8SBarry Smith PetscErrorCode CharacteristicHeapSort(Characteristic c, Queue queue, PetscInt size) 676af33a6ddSJed Brown /*---------------------------------------------------------------------*/ 677af33a6ddSJed Brown { 6780b8f38c8SBarry Smith PetscErrorCode ierr; 679af33a6ddSJed Brown CharacteristicPointDA2D temp; 6800b8f38c8SBarry Smith PetscInt n; 681af33a6ddSJed Brown 6820b8f38c8SBarry Smith PetscFunctionBegin; 683af33a6ddSJed Brown if (0) { /* Check the order of the queue before sorting */ 684*994fe344SLisandro Dalcin ierr = PetscInfo(NULL, "Before Heap sort\n");CHKERRQ(ierr); 685af33a6ddSJed Brown for (n=0; n<size; n++) { 6860298fd71SBarry Smith ierr = PetscInfo2(NULL,"%d %d\n",n,queue[n].proc);CHKERRQ(ierr); 687af33a6ddSJed Brown } 688af33a6ddSJed Brown } 689af33a6ddSJed Brown 690af33a6ddSJed Brown /* SORTING PHASE */ 6910b8f38c8SBarry Smith for (n = (size / 2)-1; n >= 0; n--) { 6920b8f38c8SBarry Smith ierr = CharacteristicSiftDown(c, queue, n, size-1);CHKERRQ(ierr); /* Rich had size-1 here, Matt had size*/ 6930b8f38c8SBarry Smith } 694af33a6ddSJed Brown for (n = size-1; n >= 1; n--) { 695af33a6ddSJed Brown temp = queue[0]; 696af33a6ddSJed Brown queue[0] = queue[n]; 697af33a6ddSJed Brown queue[n] = temp; 6980b8f38c8SBarry Smith ierr = CharacteristicSiftDown(c, queue, 0, n-1);CHKERRQ(ierr); 699af33a6ddSJed Brown } 700af33a6ddSJed Brown if (0) { /* Check the order of the queue after sorting */ 7010298fd71SBarry Smith ierr = PetscInfo(NULL, "Avter Heap sort\n");CHKERRQ(ierr); 702af33a6ddSJed Brown for (n=0; n<size; n++) { 7030298fd71SBarry Smith ierr = PetscInfo2(NULL,"%d %d\n",n,queue[n].proc);CHKERRQ(ierr); 704af33a6ddSJed Brown } 705af33a6ddSJed Brown } 7060b8f38c8SBarry Smith PetscFunctionReturn(0); 707af33a6ddSJed Brown } 708af33a6ddSJed Brown 709af33a6ddSJed Brown /*---------------------------------------------------------------------*/ 710af33a6ddSJed Brown /* 711af33a6ddSJed Brown Based on code from http://linux.wku.edu/~lamonml/algor/sort/heap.html 712af33a6ddSJed Brown */ 7130b8f38c8SBarry Smith PetscErrorCode CharacteristicSiftDown(Characteristic c, Queue queue, PetscInt root, PetscInt bottom) 714af33a6ddSJed Brown /*---------------------------------------------------------------------*/ 715af33a6ddSJed Brown { 716af33a6ddSJed Brown PetscBool done = PETSC_FALSE; 717af33a6ddSJed Brown PetscInt maxChild; 718af33a6ddSJed Brown CharacteristicPointDA2D temp; 719af33a6ddSJed Brown 720af33a6ddSJed Brown PetscFunctionBegin; 721af33a6ddSJed Brown while ((root*2 <= bottom) && (!done)) { 722af33a6ddSJed Brown if (root*2 == bottom) maxChild = root * 2; 723af33a6ddSJed Brown else if (queue[root*2].proc > queue[root*2+1].proc) maxChild = root * 2; 724af33a6ddSJed Brown else maxChild = root * 2 + 1; 725af33a6ddSJed Brown 726af33a6ddSJed Brown if (queue[root].proc < queue[maxChild].proc) { 727af33a6ddSJed Brown temp = queue[root]; 728af33a6ddSJed Brown queue[root] = queue[maxChild]; 729af33a6ddSJed Brown queue[maxChild] = temp; 730af33a6ddSJed Brown root = maxChild; 731db4deed7SKarl Rupp } else done = PETSC_TRUE; 732af33a6ddSJed Brown } 733af33a6ddSJed Brown PetscFunctionReturn(0); 734af33a6ddSJed Brown } 735af33a6ddSJed Brown 736af33a6ddSJed Brown /* [center, left, top-left, top, top-right, right, bottom-right, bottom, bottom-left] */ 737af33a6ddSJed Brown PetscErrorCode DMDAGetNeighborsRank(DM da, PetscMPIInt neighbors[]) 738af33a6ddSJed Brown { 739bff4a2f0SMatthew G. Knepley DMBoundaryType bx, by; 740af33a6ddSJed Brown PetscBool IPeriodic = PETSC_FALSE, JPeriodic = PETSC_FALSE; 741af33a6ddSJed Brown MPI_Comm comm; 742af33a6ddSJed Brown PetscMPIInt rank; 743af33a6ddSJed Brown PetscInt **procs,pi,pj,pim,pip,pjm,pjp,PI,PJ; 744af33a6ddSJed Brown PetscErrorCode ierr; 745af33a6ddSJed Brown 746af33a6ddSJed Brown PetscFunctionBegin; 747af33a6ddSJed Brown ierr = PetscObjectGetComm((PetscObject) da, &comm);CHKERRQ(ierr); 748af33a6ddSJed Brown ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 74922d28d08SBarry Smith ierr = DMDAGetInfo(da, 0, 0, 0, 0, &PI,&PJ, 0, 0, 0, &bx, &by,0, 0);CHKERRQ(ierr); 750af33a6ddSJed Brown 751bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC) IPeriodic = PETSC_TRUE; 752bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_PERIODIC) JPeriodic = PETSC_TRUE; 753af33a6ddSJed Brown 754af33a6ddSJed Brown neighbors[0] = rank; 755af33a6ddSJed Brown rank = 0; 756854ce69bSBarry Smith ierr = PetscMalloc1(PJ,&procs);CHKERRQ(ierr); 757af33a6ddSJed Brown for (pj=0; pj<PJ; pj++) { 758854ce69bSBarry Smith ierr = PetscMalloc1(PI,&(procs[pj]));CHKERRQ(ierr); 759af33a6ddSJed Brown for (pi=0; pi<PI; pi++) { 760af33a6ddSJed Brown procs[pj][pi] = rank; 761af33a6ddSJed Brown rank++; 762af33a6ddSJed Brown } 763af33a6ddSJed Brown } 764af33a6ddSJed Brown 765af33a6ddSJed Brown pi = neighbors[0] % PI; 766af33a6ddSJed Brown pj = neighbors[0] / PI; 767af33a6ddSJed Brown pim = pi-1; if (pim<0) pim=PI-1; 768af33a6ddSJed Brown pip = (pi+1)%PI; 769af33a6ddSJed Brown pjm = pj-1; if (pjm<0) pjm=PJ-1; 770af33a6ddSJed Brown pjp = (pj+1)%PJ; 771af33a6ddSJed Brown 772af33a6ddSJed Brown neighbors[1] = procs[pj] [pim]; 773af33a6ddSJed Brown neighbors[2] = procs[pjp][pim]; 774af33a6ddSJed Brown neighbors[3] = procs[pjp][pi]; 775af33a6ddSJed Brown neighbors[4] = procs[pjp][pip]; 776af33a6ddSJed Brown neighbors[5] = procs[pj] [pip]; 777af33a6ddSJed Brown neighbors[6] = procs[pjm][pip]; 778af33a6ddSJed Brown neighbors[7] = procs[pjm][pi]; 779af33a6ddSJed Brown neighbors[8] = procs[pjm][pim]; 780af33a6ddSJed Brown 781af33a6ddSJed Brown if (!IPeriodic) { 782af33a6ddSJed Brown if (pi==0) neighbors[1]=neighbors[2]=neighbors[8]=neighbors[0]; 783af33a6ddSJed Brown if (pi==PI-1) neighbors[4]=neighbors[5]=neighbors[6]=neighbors[0]; 784af33a6ddSJed Brown } 785af33a6ddSJed Brown 786af33a6ddSJed Brown if (!JPeriodic) { 787af33a6ddSJed Brown if (pj==0) neighbors[6]=neighbors[7]=neighbors[8]=neighbors[0]; 788af33a6ddSJed Brown if (pj==PJ-1) neighbors[2]=neighbors[3]=neighbors[4]=neighbors[0]; 789af33a6ddSJed Brown } 790af33a6ddSJed Brown 791af33a6ddSJed Brown for (pj = 0; pj < PJ; pj++) { 792af33a6ddSJed Brown ierr = PetscFree(procs[pj]);CHKERRQ(ierr); 793af33a6ddSJed Brown } 794af33a6ddSJed Brown ierr = PetscFree(procs);CHKERRQ(ierr); 795af33a6ddSJed Brown PetscFunctionReturn(0); 796af33a6ddSJed Brown } 797af33a6ddSJed Brown 798af33a6ddSJed Brown /* 799af33a6ddSJed Brown SUBDOMAIN NEIGHBORHOOD PROCESS MAP: 800af33a6ddSJed Brown 2 | 3 | 4 801af33a6ddSJed Brown __|___|__ 802af33a6ddSJed Brown 1 | 0 | 5 803af33a6ddSJed Brown __|___|__ 804af33a6ddSJed Brown 8 | 7 | 6 805af33a6ddSJed Brown | | 806af33a6ddSJed Brown */ 8071cc8b206SBarry Smith PetscInt DMDAGetNeighborRelative(DM da, PetscReal ir, PetscReal jr) 808af33a6ddSJed Brown { 809af33a6ddSJed Brown DMDALocalInfo info; 8101cc8b206SBarry Smith PetscReal is,ie,js,je; 811af33a6ddSJed Brown PetscErrorCode ierr; 812af33a6ddSJed Brown 813af33a6ddSJed Brown ierr = DMDAGetLocalInfo(da, &info);CHKERRQ(ierr); 8141cc8b206SBarry Smith is = (PetscReal) info.xs - 0.5; ie = (PetscReal) info.xs + info.xm - 0.5; 8151cc8b206SBarry Smith js = (PetscReal) info.ys - 0.5; je = (PetscReal) info.ys + info.ym - 0.5; 816af33a6ddSJed Brown 817af33a6ddSJed Brown if (ir >= is && ir <= ie) { /* center column */ 818bbd56ea5SKarl Rupp if (jr >= js && jr <= je) return 0; 819bbd56ea5SKarl Rupp else if (jr < js) return 7; 820bbd56ea5SKarl Rupp else return 3; 821af33a6ddSJed Brown } else if (ir < is) { /* left column */ 822bbd56ea5SKarl Rupp if (jr >= js && jr <= je) return 1; 823bbd56ea5SKarl Rupp else if (jr < js) return 8; 824bbd56ea5SKarl Rupp else return 2; 825af33a6ddSJed Brown } else { /* right column */ 826bbd56ea5SKarl Rupp if (jr >= js && jr <= je) return 5; 827bbd56ea5SKarl Rupp else if (jr < js) return 6; 828bbd56ea5SKarl Rupp else return 4; 829af33a6ddSJed Brown } 830af33a6ddSJed Brown } 831