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