1af33a6ddSJed Brown 230a76a96SBarry Smith #include <private/characteristicimpl.h> /*I "petsccharacteristic.h" I*/ 3af33a6ddSJed Brown 4af33a6ddSJed Brown PetscClassId CHARACTERISTIC_CLASSID; 5af33a6ddSJed Brown PetscLogEvent CHARACTERISTIC_SetUp, CHARACTERISTIC_Solve, CHARACTERISTIC_QueueSetup, CHARACTERISTIC_DAUpdate; 6af33a6ddSJed Brown PetscLogEvent CHARACTERISTIC_HalfTimeLocal, CHARACTERISTIC_HalfTimeRemote, CHARACTERISTIC_HalfTimeExchange; 7af33a6ddSJed Brown PetscLogEvent CHARACTERISTIC_FullTimeLocal, CHARACTERISTIC_FullTimeRemote, CHARACTERISTIC_FullTimeExchange; 8af33a6ddSJed Brown PetscBool CharacteristicRegisterAllCalled = PETSC_FALSE; 9af33a6ddSJed Brown /* 10af33a6ddSJed Brown Contains the list of registered characteristic routines 11af33a6ddSJed Brown */ 12af33a6ddSJed Brown PetscFList CharacteristicList = PETSC_NULL; 13af33a6ddSJed Brown 14af33a6ddSJed Brown PetscErrorCode DMDAGetNeighborsRank(DM, PetscMPIInt []); 15af33a6ddSJed Brown PetscInt DMDAGetNeighborRelative(DM, PassiveReal, PassiveReal); 16af33a6ddSJed Brown PetscErrorCode DMDAMapToPeriodicDomain(DM, PetscScalar [] ); 17af33a6ddSJed Brown 18af33a6ddSJed Brown PetscErrorCode HeapSort(Characteristic, Queue, PetscInt); 19af33a6ddSJed Brown PetscErrorCode SiftDown(Characteristic, Queue, PetscInt, PetscInt); 20af33a6ddSJed Brown 21af33a6ddSJed Brown #undef __FUNCT__ 22af33a6ddSJed Brown #define __FUNCT__ "CharacteristicView" 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) { 31af33a6ddSJed Brown ierr = PetscViewerASCIIGetStdout(((PetscObject)c)->comm,&viewer);CHKERRQ(ierr); 32af33a6ddSJed Brown } 33af33a6ddSJed Brown PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 34af33a6ddSJed Brown PetscCheckSameComm(c, 1, viewer, 2); 35af33a6ddSJed Brown 36af33a6ddSJed Brown ierr = PetscTypeCompare((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 45af33a6ddSJed Brown #undef __FUNCT__ 46af33a6ddSJed Brown #define __FUNCT__ "CharacteristicDestroy" 476bf464f9SBarry Smith PetscErrorCode CharacteristicDestroy(Characteristic *c) 48af33a6ddSJed Brown { 49af33a6ddSJed Brown PetscErrorCode ierr; 50af33a6ddSJed Brown 51af33a6ddSJed Brown PetscFunctionBegin; 526bf464f9SBarry Smith if (!*c) PetscFunctionReturn(0); 536bf464f9SBarry Smith PetscValidHeaderSpecific(*c, CHARACTERISTIC_CLASSID, 1); 546bf464f9SBarry Smith if (--((PetscObject)(*c))->refct > 0) PetscFunctionReturn(0); 55af33a6ddSJed Brown 566bf464f9SBarry Smith if ((*c)->ops->destroy) { 576bf464f9SBarry Smith ierr = (*(*c)->ops->destroy)((*c));CHKERRQ(ierr); 58af33a6ddSJed Brown } 596bf464f9SBarry Smith ierr = MPI_Type_free(&(*c)->itemType);CHKERRQ(ierr); 606bf464f9SBarry Smith ierr = PetscFree((*c)->queue);CHKERRQ(ierr); 616bf464f9SBarry Smith ierr = PetscFree((*c)->queueLocal);CHKERRQ(ierr); 626bf464f9SBarry Smith ierr = PetscFree((*c)->queueRemote);CHKERRQ(ierr); 636bf464f9SBarry Smith ierr = PetscFree((*c)->neighbors);CHKERRQ(ierr); 646bf464f9SBarry Smith ierr = PetscFree((*c)->needCount);CHKERRQ(ierr); 656bf464f9SBarry Smith ierr = PetscFree((*c)->localOffsets);CHKERRQ(ierr); 666bf464f9SBarry Smith ierr = PetscFree((*c)->fillCount);CHKERRQ(ierr); 676bf464f9SBarry Smith ierr = PetscFree((*c)->remoteOffsets);CHKERRQ(ierr); 686bf464f9SBarry Smith ierr = PetscFree((*c)->request);CHKERRQ(ierr); 696bf464f9SBarry Smith ierr = PetscFree((*c)->status);CHKERRQ(ierr); 70af33a6ddSJed Brown ierr = PetscHeaderDestroy(c);CHKERRQ(ierr); 71af33a6ddSJed Brown PetscFunctionReturn(0); 72af33a6ddSJed Brown } 73af33a6ddSJed Brown 74af33a6ddSJed Brown #undef __FUNCT__ 75af33a6ddSJed Brown #define __FUNCT__ "CharacteristicCreate" 76af33a6ddSJed Brown PetscErrorCode CharacteristicCreate(MPI_Comm comm, Characteristic *c) 77af33a6ddSJed Brown { 78af33a6ddSJed Brown Characteristic newC; 79af33a6ddSJed Brown PetscErrorCode ierr; 80af33a6ddSJed Brown 81af33a6ddSJed Brown PetscFunctionBegin; 82af33a6ddSJed Brown PetscValidPointer(c, 2); 83af33a6ddSJed Brown *c = PETSC_NULL; 84af33a6ddSJed Brown #ifndef PETSC_USE_DYNAMIC_LIBRARIES 85af33a6ddSJed Brown ierr = CharacteristicInitializePackage(PETSC_NULL);CHKERRQ(ierr); 86af33a6ddSJed Brown #endif 87af33a6ddSJed Brown 88af33a6ddSJed Brown ierr = PetscHeaderCreate(newC, _p_Characteristic, struct _CharacteristicOps, CHARACTERISTIC_CLASSID, -1, "Characteristic", comm, CharacteristicDestroy, CharacteristicView);CHKERRQ(ierr); 89af33a6ddSJed Brown ierr = PetscLogObjectCreate(newC);CHKERRQ(ierr); 90af33a6ddSJed Brown *c = newC; 91af33a6ddSJed Brown 92af33a6ddSJed Brown newC->structured = PETSC_TRUE; 93af33a6ddSJed Brown newC->numIds = 0; 94af33a6ddSJed Brown newC->velocityDA = PETSC_NULL; 95af33a6ddSJed Brown newC->velocity = PETSC_NULL; 96af33a6ddSJed Brown newC->velocityOld = PETSC_NULL; 97af33a6ddSJed Brown newC->numVelocityComp = 0; 98af33a6ddSJed Brown newC->velocityComp = PETSC_NULL; 99af33a6ddSJed Brown newC->velocityInterp = PETSC_NULL; 100af33a6ddSJed Brown newC->velocityInterpLocal = PETSC_NULL; 101af33a6ddSJed Brown newC->velocityCtx = PETSC_NULL; 102af33a6ddSJed Brown newC->fieldDA = PETSC_NULL; 103af33a6ddSJed Brown newC->field = PETSC_NULL; 104af33a6ddSJed Brown newC->numFieldComp = 0; 105af33a6ddSJed Brown newC->fieldComp = PETSC_NULL; 106af33a6ddSJed Brown newC->fieldInterp = PETSC_NULL; 107af33a6ddSJed Brown newC->fieldInterpLocal = PETSC_NULL; 108af33a6ddSJed Brown newC->fieldCtx = PETSC_NULL; 109af33a6ddSJed Brown newC->itemType = PETSC_NULL; 110af33a6ddSJed Brown newC->queue = PETSC_NULL; 111af33a6ddSJed Brown newC->queueSize = 0; 112af33a6ddSJed Brown newC->queueMax = 0; 113af33a6ddSJed Brown newC->queueLocal = PETSC_NULL; 114af33a6ddSJed Brown newC->queueLocalSize = 0; 115af33a6ddSJed Brown newC->queueLocalMax = 0; 116af33a6ddSJed Brown newC->queueRemote = PETSC_NULL; 117af33a6ddSJed Brown newC->queueRemoteSize = 0; 118af33a6ddSJed Brown newC->queueRemoteMax = 0; 119af33a6ddSJed Brown newC->numNeighbors = 0; 120af33a6ddSJed Brown newC->neighbors = PETSC_NULL; 121af33a6ddSJed Brown newC->needCount = PETSC_NULL; 122af33a6ddSJed Brown newC->localOffsets = PETSC_NULL; 123af33a6ddSJed Brown newC->fillCount = PETSC_NULL; 124af33a6ddSJed Brown newC->remoteOffsets = PETSC_NULL; 125af33a6ddSJed Brown newC->request = PETSC_NULL; 126af33a6ddSJed Brown newC->status = PETSC_NULL; 127af33a6ddSJed Brown PetscFunctionReturn(0); 128af33a6ddSJed Brown } 129af33a6ddSJed Brown 130af33a6ddSJed Brown #undef __FUNCT__ 131af33a6ddSJed Brown #define __FUNCT__ "CharacteristicSetType" 132af33a6ddSJed Brown /*@C 133af33a6ddSJed Brown CharacteristicSetType - Builds Characteristic for a particular solver. 134af33a6ddSJed Brown 135af33a6ddSJed Brown Logically Collective on Characteristic 136af33a6ddSJed Brown 137af33a6ddSJed Brown Input Parameters: 138af33a6ddSJed Brown + c - the method of characteristics context 139af33a6ddSJed Brown - type - a known method 140af33a6ddSJed Brown 141af33a6ddSJed Brown Options Database Key: 142af33a6ddSJed Brown . -characteristic_type <method> - Sets the method; use -help for a list 143af33a6ddSJed Brown of available methods 144af33a6ddSJed Brown 145af33a6ddSJed Brown Notes: 14630a76a96SBarry Smith See "include/petsccharacteristic.h" for available methods 147af33a6ddSJed Brown 148af33a6ddSJed Brown Normally, it is best to use the CharacteristicSetFromOptions() command and 149af33a6ddSJed Brown then set the Characteristic type from the options database rather than by using 150af33a6ddSJed Brown this routine. Using the options database provides the user with 151af33a6ddSJed Brown maximum flexibility in evaluating the many different Krylov methods. 152af33a6ddSJed Brown The CharacteristicSetType() routine is provided for those situations where it 153af33a6ddSJed Brown is necessary to set the iterative solver independently of the command 154af33a6ddSJed Brown line or options database. This might be the case, for example, when 155af33a6ddSJed Brown the choice of iterative solver changes during the execution of the 156af33a6ddSJed Brown program, and the user's application is taking responsibility for 157af33a6ddSJed Brown choosing the appropriate method. In other words, this routine is 158af33a6ddSJed Brown not for beginners. 159af33a6ddSJed Brown 160af33a6ddSJed Brown Level: intermediate 161af33a6ddSJed Brown 162af33a6ddSJed Brown .keywords: Characteristic, set, method 163af33a6ddSJed Brown 164af33a6ddSJed Brown .seealso: CharacteristicType 165af33a6ddSJed Brown 166af33a6ddSJed Brown @*/ 167af33a6ddSJed Brown PetscErrorCode CharacteristicSetType(Characteristic c, const CharacteristicType type) 168af33a6ddSJed Brown { 169af33a6ddSJed Brown PetscErrorCode ierr, (*r)(Characteristic); 170af33a6ddSJed Brown PetscBool match; 171af33a6ddSJed Brown 172af33a6ddSJed Brown PetscFunctionBegin; 173af33a6ddSJed Brown PetscValidHeaderSpecific(c, CHARACTERISTIC_CLASSID, 1); 174af33a6ddSJed Brown PetscValidCharPointer(type, 2); 175af33a6ddSJed Brown 176af33a6ddSJed Brown ierr = PetscTypeCompare((PetscObject) c, type, &match);CHKERRQ(ierr); 177af33a6ddSJed Brown if (match) PetscFunctionReturn(0); 178af33a6ddSJed Brown 179af33a6ddSJed Brown if (c->data) { 180af33a6ddSJed Brown /* destroy the old private Characteristic context */ 181af33a6ddSJed Brown ierr = (*c->ops->destroy)(c);CHKERRQ(ierr); 182af33a6ddSJed Brown c->data = 0; 183af33a6ddSJed Brown } 184af33a6ddSJed Brown 1854b91b6eaSBarry Smith ierr = PetscFListFind(CharacteristicList, ((PetscObject)c)->comm,type,PETSC_TRUE, (void (**)(void)) &r);CHKERRQ(ierr); 186af33a6ddSJed Brown if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown Characteristic type given: %s", type); 187af33a6ddSJed Brown c->setupcalled = 0; 188af33a6ddSJed Brown ierr = (*r)(c);CHKERRQ(ierr); 189af33a6ddSJed Brown ierr = PetscObjectChangeTypeName((PetscObject) c, type);CHKERRQ(ierr); 190af33a6ddSJed Brown PetscFunctionReturn(0); 191af33a6ddSJed Brown } 192af33a6ddSJed Brown 193af33a6ddSJed Brown #undef __FUNCT__ 194af33a6ddSJed Brown #define __FUNCT__ "CharacteristicSetUp" 195af33a6ddSJed Brown /*@ 196af33a6ddSJed Brown CharacteristicSetUp - Sets up the internal data structures for the 197af33a6ddSJed Brown later use of an iterative solver. 198af33a6ddSJed Brown 199af33a6ddSJed Brown Collective on Characteristic 200af33a6ddSJed Brown 201af33a6ddSJed Brown Input Parameter: 202af33a6ddSJed Brown . ksp - iterative context obtained from CharacteristicCreate() 203af33a6ddSJed Brown 204af33a6ddSJed Brown Level: developer 205af33a6ddSJed Brown 206af33a6ddSJed Brown .keywords: Characteristic, setup 207af33a6ddSJed Brown 208af33a6ddSJed Brown .seealso: CharacteristicCreate(), CharacteristicSolve(), CharacteristicDestroy() 209af33a6ddSJed Brown @*/ 210af33a6ddSJed Brown PetscErrorCode CharacteristicSetUp(Characteristic c) 211af33a6ddSJed Brown { 212af33a6ddSJed Brown PetscErrorCode ierr; 213af33a6ddSJed Brown 214af33a6ddSJed Brown PetscFunctionBegin; 215af33a6ddSJed Brown PetscValidHeaderSpecific(c, CHARACTERISTIC_CLASSID, 1); 216af33a6ddSJed Brown 217af33a6ddSJed Brown if (!((PetscObject)c)->type_name){ 218af33a6ddSJed Brown ierr = CharacteristicSetType(c, CHARACTERISTICDA);CHKERRQ(ierr); 219af33a6ddSJed Brown } 220af33a6ddSJed Brown 221af33a6ddSJed Brown if (c->setupcalled == 2) PetscFunctionReturn(0); 222af33a6ddSJed Brown 223af33a6ddSJed Brown ierr = PetscLogEventBegin(CHARACTERISTIC_SetUp,c,0,0,0);CHKERRQ(ierr); 224af33a6ddSJed Brown if (!c->setupcalled) { 225af33a6ddSJed Brown ierr = (*c->ops->setup)(c);CHKERRQ(ierr); 226af33a6ddSJed Brown } 227af33a6ddSJed Brown ierr = PetscLogEventEnd(CHARACTERISTIC_SetUp,c,0,0,0);CHKERRQ(ierr); 228af33a6ddSJed Brown c->setupcalled = 2; 229af33a6ddSJed Brown PetscFunctionReturn(0); 230af33a6ddSJed Brown } 231af33a6ddSJed Brown 232af33a6ddSJed Brown #undef __FUNCT__ 233af33a6ddSJed Brown #define __FUNCT__ "CharacteristicRegister" 234af33a6ddSJed Brown /*@C 235af33a6ddSJed Brown CharacteristicRegister - See CharacteristicRegisterDynamic() 236af33a6ddSJed Brown 237af33a6ddSJed Brown Level: advanced 238af33a6ddSJed Brown @*/ 239af33a6ddSJed Brown PetscErrorCode CharacteristicRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(Characteristic)) 240af33a6ddSJed Brown { 241af33a6ddSJed Brown PetscErrorCode ierr; 242af33a6ddSJed Brown char fullname[PETSC_MAX_PATH_LEN]; 243af33a6ddSJed Brown 244af33a6ddSJed Brown PetscFunctionBegin; 245af33a6ddSJed Brown ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr); 246af33a6ddSJed Brown ierr = PetscFListAdd(&CharacteristicList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr); 247af33a6ddSJed Brown PetscFunctionReturn(0); 248af33a6ddSJed Brown } 249af33a6ddSJed Brown 250af33a6ddSJed Brown #undef __FUNCT__ 251af33a6ddSJed Brown #define __FUNCT__ "CharacteristicSetVelocityInterpolation" 252af33a6ddSJed 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) 253af33a6ddSJed Brown { 254af33a6ddSJed Brown PetscFunctionBegin; 255af33a6ddSJed Brown c->velocityDA = da; 256af33a6ddSJed Brown c->velocity = v; 257af33a6ddSJed Brown c->velocityOld = vOld; 258af33a6ddSJed Brown c->numVelocityComp = numComponents; 259af33a6ddSJed Brown c->velocityComp = components; 260af33a6ddSJed Brown c->velocityInterp = interp; 261af33a6ddSJed Brown c->velocityCtx = ctx; 262af33a6ddSJed Brown PetscFunctionReturn(0); 263af33a6ddSJed Brown } 264af33a6ddSJed Brown 265af33a6ddSJed Brown #undef __FUNCT__ 266af33a6ddSJed Brown #define __FUNCT__ "CharacteristicSetVelocityInterpolationLocal" 267af33a6ddSJed 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) 268af33a6ddSJed Brown { 269af33a6ddSJed Brown PetscFunctionBegin; 270af33a6ddSJed Brown c->velocityDA = da; 271af33a6ddSJed Brown c->velocity = v; 272af33a6ddSJed Brown c->velocityOld = vOld; 273af33a6ddSJed Brown c->numVelocityComp = numComponents; 274af33a6ddSJed Brown c->velocityComp = components; 275af33a6ddSJed Brown c->velocityInterpLocal = interp; 276af33a6ddSJed Brown c->velocityCtx = ctx; 277af33a6ddSJed Brown PetscFunctionReturn(0); 278af33a6ddSJed Brown } 279af33a6ddSJed Brown 280af33a6ddSJed Brown #undef __FUNCT__ 281af33a6ddSJed Brown #define __FUNCT__ "CharacteristicSetFieldInterpolation" 282af33a6ddSJed Brown PetscErrorCode CharacteristicSetFieldInterpolation(Characteristic c, DM da, Vec v, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(Vec, PetscReal [], PetscInt, PetscInt [], PetscScalar [], void *), void *ctx) 283af33a6ddSJed Brown { 284af33a6ddSJed Brown PetscFunctionBegin; 285af33a6ddSJed Brown #if 0 286af33a6ddSJed 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."); 287af33a6ddSJed Brown #endif 288af33a6ddSJed Brown c->fieldDA = da; 289af33a6ddSJed Brown c->field = v; 290af33a6ddSJed Brown c->numFieldComp = numComponents; 291af33a6ddSJed Brown c->fieldComp = components; 292af33a6ddSJed Brown c->fieldInterp = interp; 293af33a6ddSJed Brown c->fieldCtx = ctx; 294af33a6ddSJed Brown PetscFunctionReturn(0); 295af33a6ddSJed Brown } 296af33a6ddSJed Brown 297af33a6ddSJed Brown #undef __FUNCT__ 298af33a6ddSJed Brown #define __FUNCT__ "CharacteristicSetFieldInterpolationLocal" 299af33a6ddSJed Brown PetscErrorCode CharacteristicSetFieldInterpolationLocal(Characteristic c, DM da, Vec v, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(void *, PetscReal [], PetscInt, PetscInt [], PetscScalar [], void *), void *ctx) 300af33a6ddSJed Brown { 301af33a6ddSJed Brown PetscFunctionBegin; 302af33a6ddSJed Brown #if 0 303af33a6ddSJed 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."); 304af33a6ddSJed Brown #endif 305af33a6ddSJed Brown c->fieldDA = da; 306af33a6ddSJed Brown c->field = v; 307af33a6ddSJed Brown c->numFieldComp = numComponents; 308af33a6ddSJed Brown c->fieldComp = components; 309af33a6ddSJed Brown c->fieldInterpLocal = interp; 310af33a6ddSJed Brown c->fieldCtx = ctx; 311af33a6ddSJed Brown PetscFunctionReturn(0); 312af33a6ddSJed Brown } 313af33a6ddSJed Brown 314af33a6ddSJed Brown #undef __FUNCT__ 315af33a6ddSJed Brown #define __FUNCT__ "CharacteristicSolve" 316af33a6ddSJed Brown PetscErrorCode CharacteristicSolve(Characteristic c, PetscReal dt, Vec solution) 317af33a6ddSJed Brown { 318af33a6ddSJed Brown CharacteristicPointDA2D Qi; 319af33a6ddSJed Brown DM da = c->velocityDA; 320af33a6ddSJed Brown Vec velocityLocal, velocityLocalOld; 321af33a6ddSJed Brown Vec fieldLocal; 322af33a6ddSJed Brown DMDALocalInfo info; 323af33a6ddSJed Brown PetscScalar **solArray; 324af33a6ddSJed Brown void *velocityArray; 325af33a6ddSJed Brown void *velocityArrayOld; 326af33a6ddSJed Brown void *fieldArray; 327af33a6ddSJed Brown PassiveScalar *interpIndices; 328af33a6ddSJed Brown PassiveScalar *velocityValues, *velocityValuesOld; 329af33a6ddSJed Brown PassiveScalar *fieldValues; 330af33a6ddSJed Brown PetscMPIInt rank; 331af33a6ddSJed Brown PetscInt dim; 332af33a6ddSJed Brown PetscMPIInt neighbors[9]; 333af33a6ddSJed Brown PetscInt dof; 334af33a6ddSJed Brown PetscInt gx, gy; 335af33a6ddSJed Brown PetscInt n, is, ie, js, je, comp; 336af33a6ddSJed Brown PetscErrorCode ierr; 337af33a6ddSJed Brown 338af33a6ddSJed Brown PetscFunctionBegin; 339af33a6ddSJed Brown c->queueSize = 0; 340af33a6ddSJed Brown ierr = MPI_Comm_rank(((PetscObject)c)->comm, &rank);CHKERRQ(ierr); 341af33a6ddSJed Brown ierr = DMDAGetNeighborsRank(da, neighbors);CHKERRQ(ierr); 342af33a6ddSJed Brown ierr = CharacteristicSetNeighbors(c, 9, neighbors);CHKERRQ(ierr); 343af33a6ddSJed Brown ierr = CharacteristicSetUp(c);CHKERRQ(ierr); 344af33a6ddSJed Brown /* global and local grid info */ 345af33a6ddSJed Brown ierr = DMDAGetInfo(da, &dim, &gx, &gy, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0);CHKERRQ(ierr); 346af33a6ddSJed Brown ierr = DMDAGetLocalInfo(da, &info);CHKERRQ(ierr); 347af33a6ddSJed Brown is = info.xs; ie = info.xs+info.xm; 348af33a6ddSJed Brown js = info.ys; je = info.ys+info.ym; 349af33a6ddSJed Brown /* Allocation */ 350af33a6ddSJed Brown ierr = PetscMalloc(dim*sizeof(PetscScalar), &interpIndices);CHKERRQ(ierr); 351af33a6ddSJed Brown ierr = PetscMalloc(c->numVelocityComp*sizeof(PetscScalar), &velocityValues);CHKERRQ(ierr); 352af33a6ddSJed Brown ierr = PetscMalloc(c->numVelocityComp*sizeof(PetscScalar), &velocityValuesOld);CHKERRQ(ierr); 353af33a6ddSJed Brown ierr = PetscMalloc(c->numFieldComp*sizeof(PetscScalar), &fieldValues);CHKERRQ(ierr); 354af33a6ddSJed Brown ierr = PetscLogEventBegin(CHARACTERISTIC_Solve,0,0,0,0);CHKERRQ(ierr); 355af33a6ddSJed Brown 356af33a6ddSJed Brown /* ----------------------------------------------------------------------- 357af33a6ddSJed Brown PART 1, AT t-dt/2 358af33a6ddSJed Brown -----------------------------------------------------------------------*/ 359af33a6ddSJed Brown ierr = PetscLogEventBegin(CHARACTERISTIC_QueueSetup,0,0,0,0);CHKERRQ(ierr); 360af33a6ddSJed Brown /* GET POSITION AT HALF TIME IN THE PAST */ 361af33a6ddSJed Brown if (c->velocityInterpLocal) { 362af33a6ddSJed Brown ierr = DMGetLocalVector(c->velocityDA, &velocityLocal);CHKERRQ(ierr); 363af33a6ddSJed Brown ierr = DMGetLocalVector(c->velocityDA, &velocityLocalOld);CHKERRQ(ierr); 364af33a6ddSJed Brown ierr = DMGlobalToLocalBegin(c->velocityDA, c->velocity, INSERT_VALUES, velocityLocal);CHKERRQ(ierr); 365af33a6ddSJed Brown ierr = DMGlobalToLocalEnd(c->velocityDA, c->velocity, INSERT_VALUES, velocityLocal);CHKERRQ(ierr); 366af33a6ddSJed Brown ierr = DMGlobalToLocalBegin(c->velocityDA, c->velocityOld, INSERT_VALUES, velocityLocalOld);CHKERRQ(ierr); 367af33a6ddSJed Brown ierr = DMGlobalToLocalEnd(c->velocityDA, c->velocityOld, INSERT_VALUES, velocityLocalOld);CHKERRQ(ierr); 368af33a6ddSJed Brown ierr = DMDAVecGetArray(c->velocityDA, velocityLocal, &velocityArray); CHKERRQ(ierr); 369af33a6ddSJed Brown ierr = DMDAVecGetArray(c->velocityDA, velocityLocalOld, &velocityArrayOld);CHKERRQ(ierr); 370af33a6ddSJed Brown } 371af33a6ddSJed Brown ierr = PetscInfo(PETSC_NULL, "Calculating position at t_{n - 1/2}\n");CHKERRQ(ierr); 372af33a6ddSJed Brown for(Qi.j = js; Qi.j < je; Qi.j++) { 373af33a6ddSJed Brown for(Qi.i = is; Qi.i < ie; Qi.i++) { 374af33a6ddSJed Brown interpIndices[0] = Qi.i; 375af33a6ddSJed Brown interpIndices[1] = Qi.j; 376af33a6ddSJed Brown if (c->velocityInterpLocal) { 377af33a6ddSJed Brown c->velocityInterpLocal(velocityArray, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx); 378af33a6ddSJed Brown } else { 379af33a6ddSJed Brown c->velocityInterp(c->velocity, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx); 380af33a6ddSJed Brown } 381af33a6ddSJed Brown Qi.x = Qi.i - velocityValues[0]*dt/2.0; 382af33a6ddSJed Brown Qi.y = Qi.j - velocityValues[1]*dt/2.0; 383af33a6ddSJed Brown 384af33a6ddSJed Brown /* Determine whether the position at t - dt/2 is local */ 385af33a6ddSJed Brown Qi.proc = DMDAGetNeighborRelative(da, Qi.x, Qi.y); 386af33a6ddSJed Brown 387af33a6ddSJed Brown /* Check for Periodic boundaries and move all periodic points back onto the domain */ 388af33a6ddSJed Brown ierr = DMDAMapCoordsToPeriodicDomain(da,&(Qi.x),&(Qi.y));CHKERRQ(ierr); 389af33a6ddSJed Brown ierr = CharacteristicAddPoint(c, &Qi);CHKERRQ(ierr); 390af33a6ddSJed Brown } 391af33a6ddSJed Brown } 392af33a6ddSJed Brown ierr = PetscLogEventEnd(CHARACTERISTIC_QueueSetup,0,0,0,0);CHKERRQ(ierr); 393af33a6ddSJed Brown 394af33a6ddSJed Brown ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr); 395af33a6ddSJed Brown ierr = CharacteristicSendCoordinatesBegin(c);CHKERRQ(ierr); 396af33a6ddSJed Brown ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr); 397af33a6ddSJed Brown 398af33a6ddSJed Brown ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeLocal,0,0,0,0);CHKERRQ(ierr); 399af33a6ddSJed Brown /* Calculate velocity at t_n+1/2 (local values) */ 400af33a6ddSJed Brown ierr = PetscInfo(PETSC_NULL, "Calculating local velocities at t_{n - 1/2}\n");CHKERRQ(ierr); 401af33a6ddSJed Brown for(n = 0; n < c->queueSize; n++) { 402af33a6ddSJed Brown Qi = c->queue[n]; 403af33a6ddSJed Brown if (c->neighbors[Qi.proc] == rank) { 404af33a6ddSJed Brown interpIndices[0] = Qi.x; 405af33a6ddSJed Brown interpIndices[1] = Qi.y; 406af33a6ddSJed Brown if (c->velocityInterpLocal) { 407af33a6ddSJed Brown c->velocityInterpLocal(velocityArray, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx); 408af33a6ddSJed Brown c->velocityInterpLocal(velocityArrayOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx); 409af33a6ddSJed Brown } else { 410af33a6ddSJed Brown c->velocityInterp(c->velocity, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx); 411af33a6ddSJed Brown c->velocityInterp(c->velocityOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx); 412af33a6ddSJed Brown } 413af33a6ddSJed Brown Qi.x = 0.5*(velocityValues[0] + velocityValuesOld[0]); 414af33a6ddSJed Brown Qi.y = 0.5*(velocityValues[1] + velocityValuesOld[1]); 415af33a6ddSJed Brown } 416af33a6ddSJed Brown c->queue[n] = Qi; 417af33a6ddSJed Brown } 418af33a6ddSJed Brown ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeLocal,0,0,0,0);CHKERRQ(ierr); 419af33a6ddSJed Brown 420af33a6ddSJed Brown ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr); 421af33a6ddSJed Brown ierr = CharacteristicSendCoordinatesEnd(c);CHKERRQ(ierr); 422af33a6ddSJed Brown ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr); 423af33a6ddSJed Brown 424af33a6ddSJed Brown 425af33a6ddSJed Brown /* Calculate velocity at t_n+1/2 (fill remote requests) */ 426af33a6ddSJed Brown ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeRemote,0,0,0,0);CHKERRQ(ierr); 427af33a6ddSJed Brown ierr = PetscInfo1(PETSC_NULL, "Calculating %d remote velocities at t_{n - 1/2}\n", c->queueRemoteSize);CHKERRQ(ierr); 428af33a6ddSJed Brown for(n = 0; n < c->queueRemoteSize; n++) { 429af33a6ddSJed Brown Qi = c->queueRemote[n]; 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 c->queueRemote[n] = Qi; 442af33a6ddSJed Brown } 443af33a6ddSJed Brown ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeRemote,0,0,0,0);CHKERRQ(ierr); 444af33a6ddSJed Brown ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr); 445af33a6ddSJed Brown ierr = CharacteristicGetValuesBegin(c);CHKERRQ(ierr); 446af33a6ddSJed Brown ierr = CharacteristicGetValuesEnd(c);CHKERRQ(ierr); 447af33a6ddSJed Brown if (c->velocityInterpLocal) { 448af33a6ddSJed Brown ierr = DMDAVecRestoreArray(c->velocityDA, velocityLocal, &velocityArray); CHKERRQ(ierr); 449af33a6ddSJed Brown ierr = DMDAVecRestoreArray(c->velocityDA, velocityLocalOld, &velocityArrayOld);CHKERRQ(ierr); 450af33a6ddSJed Brown ierr = DMRestoreLocalVector(c->velocityDA, &velocityLocal);CHKERRQ(ierr); 451af33a6ddSJed Brown ierr = DMRestoreLocalVector(c->velocityDA, &velocityLocalOld);CHKERRQ(ierr); 452af33a6ddSJed Brown } 453af33a6ddSJed Brown ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr); 454af33a6ddSJed Brown 455af33a6ddSJed Brown /* ----------------------------------------------------------------------- 456af33a6ddSJed Brown PART 2, AT t-dt 457af33a6ddSJed Brown -----------------------------------------------------------------------*/ 458af33a6ddSJed Brown 459af33a6ddSJed Brown /* GET POSITION AT t_n (local values) */ 460af33a6ddSJed Brown ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeLocal,0,0,0,0);CHKERRQ(ierr); 461af33a6ddSJed Brown ierr = PetscInfo(PETSC_NULL, "Calculating position at t_{n}\n");CHKERRQ(ierr); 462af33a6ddSJed Brown for(n = 0; n < c->queueSize; n++) { 463af33a6ddSJed Brown Qi = c->queue[n]; 464af33a6ddSJed Brown Qi.x = Qi.i - Qi.x*dt; 465af33a6ddSJed Brown Qi.y = Qi.j - Qi.y*dt; 466af33a6ddSJed Brown 467af33a6ddSJed Brown /* Determine whether the position at t-dt is local */ 468af33a6ddSJed Brown Qi.proc = DMDAGetNeighborRelative(da, Qi.x, Qi.y); 469af33a6ddSJed Brown 470af33a6ddSJed Brown /* Check for Periodic boundaries and move all periodic points back onto the domain */ 471af33a6ddSJed Brown ierr = DMDAMapCoordsToPeriodicDomain(da,&(Qi.x),&(Qi.y));CHKERRQ(ierr); 472af33a6ddSJed Brown 473af33a6ddSJed Brown c->queue[n] = Qi; 474af33a6ddSJed Brown } 475af33a6ddSJed Brown ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeLocal,0,0,0,0);CHKERRQ(ierr); 476af33a6ddSJed Brown 477af33a6ddSJed Brown ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr); 478af33a6ddSJed Brown ierr = CharacteristicSendCoordinatesBegin(c);CHKERRQ(ierr); 479af33a6ddSJed Brown ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr); 480af33a6ddSJed Brown 481af33a6ddSJed Brown /* GET VALUE AT FULL TIME IN THE PAST (LOCAL REQUESTS) */ 482af33a6ddSJed Brown ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeLocal,0,0,0,0);CHKERRQ(ierr); 483af33a6ddSJed Brown if (c->fieldInterpLocal) { 484af33a6ddSJed Brown ierr = DMGetLocalVector(c->fieldDA, &fieldLocal);CHKERRQ(ierr); 485af33a6ddSJed Brown ierr = DMGlobalToLocalBegin(c->fieldDA, c->field, INSERT_VALUES, fieldLocal);CHKERRQ(ierr); 486af33a6ddSJed Brown ierr = DMGlobalToLocalEnd(c->fieldDA, c->field, INSERT_VALUES, fieldLocal);CHKERRQ(ierr); 487af33a6ddSJed Brown ierr = DMDAVecGetArray(c->fieldDA, fieldLocal, &fieldArray);CHKERRQ(ierr); 488af33a6ddSJed Brown } 489af33a6ddSJed Brown ierr = PetscInfo(PETSC_NULL, "Calculating local field at t_{n}\n");CHKERRQ(ierr); 490af33a6ddSJed Brown for(n = 0; n < c->queueSize; n++) { 491af33a6ddSJed Brown if (c->neighbors[c->queue[n].proc] == rank) { 492af33a6ddSJed Brown interpIndices[0] = c->queue[n].x; 493af33a6ddSJed Brown interpIndices[1] = c->queue[n].y; 494af33a6ddSJed Brown if (c->fieldInterpLocal) { 495af33a6ddSJed Brown c->fieldInterpLocal(fieldArray, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx); 496af33a6ddSJed Brown } else { 497af33a6ddSJed Brown c->fieldInterp(c->field, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx); 498af33a6ddSJed Brown } 499af33a6ddSJed Brown for(comp = 0; comp < c->numFieldComp; comp++) { 500af33a6ddSJed Brown c->queue[n].field[comp] = fieldValues[comp]; 501af33a6ddSJed Brown } 502af33a6ddSJed Brown } 503af33a6ddSJed Brown } 504af33a6ddSJed Brown ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeLocal,0,0,0,0);CHKERRQ(ierr); 505af33a6ddSJed Brown 506af33a6ddSJed Brown ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr); 507af33a6ddSJed Brown ierr = CharacteristicSendCoordinatesEnd(c);CHKERRQ(ierr); 508af33a6ddSJed Brown ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr); 509af33a6ddSJed Brown 510af33a6ddSJed Brown /* GET VALUE AT FULL TIME IN THE PAST (REMOTE REQUESTS) */ 511af33a6ddSJed Brown ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeRemote,0,0,0,0);CHKERRQ(ierr); 512af33a6ddSJed Brown ierr = PetscInfo1(PETSC_NULL, "Calculating %d remote field points at t_{n}\n", c->queueRemoteSize);CHKERRQ(ierr); 513af33a6ddSJed Brown for(n = 0; n < c->queueRemoteSize; n++) { 514af33a6ddSJed Brown interpIndices[0] = c->queueRemote[n].x; 515af33a6ddSJed Brown interpIndices[1] = c->queueRemote[n].y; 516af33a6ddSJed Brown 517af33a6ddSJed Brown /* for debugging purposes */ 518af33a6ddSJed Brown if (1) { /* hacked bounds test...let's do better */ 519af33a6ddSJed Brown PetscScalar im = interpIndices[0]; PetscScalar jm = interpIndices[1]; 520af33a6ddSJed Brown 521af33a6ddSJed Brown if (( im < (PetscScalar) is - 1.) || (im > (PetscScalar) ie) || (jm < (PetscScalar) js - 1.) || (jm > (PetscScalar) je)) { 522af33a6ddSJed Brown SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB, "Nonlocal point: (%g,%g)", im, jm); 523af33a6ddSJed Brown } 524af33a6ddSJed Brown } 525af33a6ddSJed Brown 526af33a6ddSJed Brown if (c->fieldInterpLocal) { 527af33a6ddSJed Brown c->fieldInterpLocal(fieldArray, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx); 528af33a6ddSJed Brown } else { 529af33a6ddSJed Brown c->fieldInterp(c->field, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx); 530af33a6ddSJed Brown } 531af33a6ddSJed Brown for(comp = 0; comp < c->numFieldComp; comp++) { 532af33a6ddSJed Brown c->queueRemote[n].field[comp] = fieldValues[comp]; 533af33a6ddSJed Brown } 534af33a6ddSJed Brown } 535af33a6ddSJed Brown ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeRemote,0,0,0,0);CHKERRQ(ierr); 536af33a6ddSJed Brown 537af33a6ddSJed Brown ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr); 538af33a6ddSJed Brown ierr = CharacteristicGetValuesBegin(c);CHKERRQ(ierr); 539af33a6ddSJed Brown ierr = CharacteristicGetValuesEnd(c);CHKERRQ(ierr); 540af33a6ddSJed Brown if (c->fieldInterpLocal) { 541af33a6ddSJed Brown ierr = DMDAVecRestoreArray(c->fieldDA, fieldLocal, &fieldArray);CHKERRQ(ierr); 542af33a6ddSJed Brown ierr = DMRestoreLocalVector(c->fieldDA, &fieldLocal);CHKERRQ(ierr); 543af33a6ddSJed Brown } 544af33a6ddSJed Brown ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr); 545af33a6ddSJed Brown 546af33a6ddSJed Brown /* Return field of characteristics at t_n-1 */ 547af33a6ddSJed Brown ierr = PetscLogEventBegin(CHARACTERISTIC_DAUpdate,0,0,0,0);CHKERRQ(ierr); 548af33a6ddSJed Brown ierr = DMDAGetInfo(c->fieldDA,0,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); 549af33a6ddSJed Brown ierr = DMDAVecGetArray(c->fieldDA, solution, &solArray);CHKERRQ(ierr); 550af33a6ddSJed Brown for(n = 0; n < c->queueSize; n++) { 551af33a6ddSJed Brown Qi = c->queue[n]; 552af33a6ddSJed Brown for(comp = 0; comp < c->numFieldComp; comp++) { 553af33a6ddSJed Brown solArray[Qi.j][Qi.i*dof+c->fieldComp[comp]] = Qi.field[comp]; 554af33a6ddSJed Brown } 555af33a6ddSJed Brown } 556af33a6ddSJed Brown ierr = DMDAVecRestoreArray(c->fieldDA, solution, &solArray);CHKERRQ(ierr); 557af33a6ddSJed Brown ierr = PetscLogEventEnd(CHARACTERISTIC_DAUpdate,0,0,0,0);CHKERRQ(ierr); 558af33a6ddSJed Brown ierr = PetscLogEventEnd(CHARACTERISTIC_Solve,0,0,0,0);CHKERRQ(ierr); 559af33a6ddSJed Brown 560af33a6ddSJed Brown /* Cleanup */ 561af33a6ddSJed Brown ierr = PetscFree(interpIndices);CHKERRQ(ierr); 562af33a6ddSJed Brown ierr = PetscFree(velocityValues);CHKERRQ(ierr); 563af33a6ddSJed Brown ierr = PetscFree(velocityValuesOld);CHKERRQ(ierr); 564af33a6ddSJed Brown ierr = PetscFree(fieldValues);CHKERRQ(ierr); 565af33a6ddSJed Brown PetscFunctionReturn(0); 566af33a6ddSJed Brown } 567af33a6ddSJed Brown 568af33a6ddSJed Brown #undef __FUNCT__ 569af33a6ddSJed Brown #define __FUNCT__ "CharacteristicSetNeighbors" 570af33a6ddSJed Brown PetscErrorCode CharacteristicSetNeighbors(Characteristic c, PetscInt numNeighbors, PetscMPIInt neighbors[]) 571af33a6ddSJed Brown { 572af33a6ddSJed Brown PetscErrorCode ierr; 573af33a6ddSJed Brown 574af33a6ddSJed Brown PetscFunctionBegin; 575af33a6ddSJed Brown c->numNeighbors = numNeighbors; 576*6eaa7dc8SBarry Smith ierr = PetscFree(c->neighbors);CHKERRQ(ierr); 577af33a6ddSJed Brown ierr = PetscMalloc(numNeighbors * sizeof(PetscMPIInt), &c->neighbors);CHKERRQ(ierr); 578af33a6ddSJed Brown ierr = PetscMemcpy(c->neighbors, neighbors, numNeighbors * sizeof(PetscMPIInt));CHKERRQ(ierr); 579af33a6ddSJed Brown PetscFunctionReturn(0); 580af33a6ddSJed Brown } 581af33a6ddSJed Brown 582af33a6ddSJed Brown #undef __FUNCT__ 583af33a6ddSJed Brown #define __FUNCT__ "CharacteristicAddPoint" 584af33a6ddSJed Brown PetscErrorCode CharacteristicAddPoint(Characteristic c, CharacteristicPointDA2D *point) 585af33a6ddSJed Brown { 586af33a6ddSJed Brown PetscFunctionBegin; 587af33a6ddSJed Brown if (c->queueSize >= c->queueMax) { 588af33a6ddSJed Brown SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Exceeeded maximum queue size %d", c->queueMax); 589af33a6ddSJed Brown } 590af33a6ddSJed Brown c->queue[c->queueSize++] = *point; 591af33a6ddSJed Brown PetscFunctionReturn(0); 592af33a6ddSJed Brown } 593af33a6ddSJed Brown 594af33a6ddSJed Brown #undef __FUNCT__ 595af33a6ddSJed Brown #define __FUNCT__ "CharacteristicSendCoordinatesBegin" 596af33a6ddSJed Brown int CharacteristicSendCoordinatesBegin(Characteristic c) 597af33a6ddSJed Brown { 598af33a6ddSJed Brown PetscMPIInt rank, tag = 121; 599af33a6ddSJed Brown PetscInt i, n; 600af33a6ddSJed Brown PetscErrorCode ierr; 601af33a6ddSJed Brown 602af33a6ddSJed Brown PetscFunctionBegin; 603af33a6ddSJed Brown ierr = MPI_Comm_rank(((PetscObject)c)->comm, &rank);CHKERRQ(ierr); 604af33a6ddSJed Brown ierr = HeapSort(c, c->queue, c->queueSize);CHKERRQ(ierr); 605af33a6ddSJed Brown ierr = PetscMemzero(c->needCount, c->numNeighbors * sizeof(PetscInt));CHKERRQ(ierr); 606af33a6ddSJed Brown for(i = 0; i < c->queueSize; i++) { 607af33a6ddSJed Brown c->needCount[c->queue[i].proc]++; 608af33a6ddSJed Brown } 609af33a6ddSJed Brown c->fillCount[0] = 0; 610af33a6ddSJed Brown for(n = 1; n < c->numNeighbors; n++) { 611af33a6ddSJed Brown ierr = MPI_Irecv(&(c->fillCount[n]), 1, MPIU_INT, c->neighbors[n], tag, ((PetscObject)c)->comm, &(c->request[n-1]));CHKERRQ(ierr); 612af33a6ddSJed Brown } 613af33a6ddSJed Brown for(n = 1; n < c->numNeighbors; n++) { 614af33a6ddSJed Brown ierr = MPI_Send(&(c->needCount[n]), 1, MPIU_INT, c->neighbors[n], tag, ((PetscObject)c)->comm);CHKERRQ(ierr); 615af33a6ddSJed Brown } 616af33a6ddSJed Brown ierr = MPI_Waitall(c->numNeighbors-1, c->request, c->status);CHKERRQ(ierr); 617af33a6ddSJed Brown /* Initialize the remote queue */ 618af33a6ddSJed Brown c->queueLocalMax = c->localOffsets[0] = 0; 619af33a6ddSJed Brown c->queueRemoteMax = c->remoteOffsets[0] = 0; 620af33a6ddSJed Brown for(n = 1; n < c->numNeighbors; n++) { 621af33a6ddSJed Brown c->remoteOffsets[n] = c->queueRemoteMax; 622af33a6ddSJed Brown c->queueRemoteMax += c->fillCount[n]; 623af33a6ddSJed Brown c->localOffsets[n] = c->queueLocalMax; 624af33a6ddSJed Brown c->queueLocalMax += c->needCount[n]; 625af33a6ddSJed Brown } 626af33a6ddSJed Brown /* HACK BEGIN */ 627af33a6ddSJed Brown for(n = 1; n < c->numNeighbors; n++) { 628af33a6ddSJed Brown c->localOffsets[n] += c->needCount[0]; 629af33a6ddSJed Brown } 630af33a6ddSJed Brown c->needCount[0] = 0; 631af33a6ddSJed Brown /* HACK END */ 632af33a6ddSJed Brown if (c->queueRemoteMax) { 633af33a6ddSJed Brown ierr = PetscMalloc(sizeof(CharacteristicPointDA2D) * c->queueRemoteMax, &c->queueRemote);CHKERRQ(ierr); 634af33a6ddSJed Brown } else { 635af33a6ddSJed Brown c->queueRemote = PETSC_NULL; 636af33a6ddSJed Brown } 637af33a6ddSJed Brown c->queueRemoteSize = c->queueRemoteMax; 638af33a6ddSJed Brown 639af33a6ddSJed Brown /* Send and Receive requests for values at t_n+1/2, giving the coordinates for interpolation */ 640af33a6ddSJed Brown for(n = 1; n < c->numNeighbors; n++) { 641af33a6ddSJed Brown ierr = PetscInfo2(PETSC_NULL, "Receiving %d requests for values from proc %d\n", c->fillCount[n], c->neighbors[n]);CHKERRQ(ierr); 642af33a6ddSJed Brown ierr = MPI_Irecv(&(c->queueRemote[c->remoteOffsets[n]]), c->fillCount[n], c->itemType, c->neighbors[n], tag, ((PetscObject)c)->comm, &(c->request[n-1]));CHKERRQ(ierr); 643af33a6ddSJed Brown } 644af33a6ddSJed Brown for(n = 1; n < c->numNeighbors; n++) { 645af33a6ddSJed Brown ierr = PetscInfo2(PETSC_NULL, "Sending %d requests for values from proc %d\n", c->needCount[n], c->neighbors[n]);CHKERRQ(ierr); 646af33a6ddSJed Brown ierr = MPI_Send(&(c->queue[c->localOffsets[n]]), c->needCount[n], c->itemType, c->neighbors[n], tag, ((PetscObject)c)->comm);CHKERRQ(ierr); 647af33a6ddSJed Brown } 648af33a6ddSJed Brown PetscFunctionReturn(0); 649af33a6ddSJed Brown } 650af33a6ddSJed Brown 651af33a6ddSJed Brown #undef __FUNCT__ 652af33a6ddSJed Brown #define __FUNCT__ "CharacteristicSendCoordinatesEnd" 653af33a6ddSJed Brown PetscErrorCode CharacteristicSendCoordinatesEnd(Characteristic c) 654af33a6ddSJed Brown { 655af33a6ddSJed Brown #if 0 656af33a6ddSJed Brown PetscMPIInt rank; 657af33a6ddSJed Brown PetscInt n; 658af33a6ddSJed Brown #endif 659af33a6ddSJed Brown PetscErrorCode ierr; 660af33a6ddSJed Brown 661af33a6ddSJed Brown PetscFunctionBegin; 662af33a6ddSJed Brown ierr = MPI_Waitall(c->numNeighbors-1, c->request, c->status);CHKERRQ(ierr); 663af33a6ddSJed Brown #if 0 664af33a6ddSJed Brown ierr = MPI_Comm_rank(((PetscObject)c)->comm, &rank);CHKERRQ(ierr); 665af33a6ddSJed Brown for(n = 0; n < c->queueRemoteSize; n++) { 666af33a6ddSJed Brown if (c->neighbors[c->queueRemote[n].proc] == rank) { 667af33a6ddSJed Brown SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB, "This is fucked up, n = %d proc = %d", n, c->queueRemote[n].proc); 668af33a6ddSJed Brown } 669af33a6ddSJed Brown } 670af33a6ddSJed Brown #endif 671af33a6ddSJed Brown PetscFunctionReturn(0); 672af33a6ddSJed Brown } 673af33a6ddSJed Brown 674af33a6ddSJed Brown #undef __FUNCT__ 675af33a6ddSJed Brown #define __FUNCT__ "CharacteristicGetValuesBegin" 676af33a6ddSJed Brown PetscErrorCode CharacteristicGetValuesBegin(Characteristic c) 677af33a6ddSJed Brown { 678af33a6ddSJed Brown PetscMPIInt tag = 121; 679af33a6ddSJed Brown PetscInt n; 680af33a6ddSJed Brown PetscErrorCode ierr; 681af33a6ddSJed Brown 682af33a6ddSJed Brown PetscFunctionBegin; 683af33a6ddSJed Brown /* SEND AND RECIEVE FILLED REQUESTS for velocities at t_n+1/2 */ 684af33a6ddSJed Brown for(n = 1; n < c->numNeighbors; n++) { 685af33a6ddSJed Brown ierr = MPI_Irecv(&(c->queue[c->localOffsets[n]]), c->needCount[n], c->itemType, c->neighbors[n], tag, ((PetscObject)c)->comm, &(c->request[n-1]));CHKERRQ(ierr); 686af33a6ddSJed Brown } 687af33a6ddSJed Brown for(n = 1; n < c->numNeighbors; n++) { 688af33a6ddSJed Brown ierr = MPI_Send(&(c->queueRemote[c->remoteOffsets[n]]), c->fillCount[n], c->itemType, c->neighbors[n], tag, ((PetscObject)c)->comm);CHKERRQ(ierr); 689af33a6ddSJed Brown } 690af33a6ddSJed Brown PetscFunctionReturn(0); 691af33a6ddSJed Brown } 692af33a6ddSJed Brown 693af33a6ddSJed Brown #undef __FUNCT__ 694af33a6ddSJed Brown #define __FUNCT__ "CharacteristicGetValuesEnd" 695af33a6ddSJed Brown PetscErrorCode CharacteristicGetValuesEnd(Characteristic c) 696af33a6ddSJed Brown { 697af33a6ddSJed Brown PetscErrorCode ierr; 698af33a6ddSJed Brown 699af33a6ddSJed Brown PetscFunctionBegin; 700af33a6ddSJed Brown ierr = MPI_Waitall(c->numNeighbors-1, c->request, c->status);CHKERRQ(ierr); 701af33a6ddSJed Brown /* Free queue of requests from other procs */ 702af33a6ddSJed Brown ierr = PetscFree(c->queueRemote);CHKERRQ(ierr); 703af33a6ddSJed Brown PetscFunctionReturn(0); 704af33a6ddSJed Brown } 705af33a6ddSJed Brown 706af33a6ddSJed Brown /*---------------------------------------------------------------------*/ 707af33a6ddSJed Brown #undef __FUNCT__ 708af33a6ddSJed Brown #define __FUNCT__ "HeapSort" 709af33a6ddSJed Brown /* 710af33a6ddSJed Brown Based on code from http://linux.wku.edu/~lamonml/algor/sort/heap.html 711af33a6ddSJed Brown */ 712af33a6ddSJed Brown int HeapSort(Characteristic c, Queue queue, PetscInt size) 713af33a6ddSJed Brown /*---------------------------------------------------------------------*/ 714af33a6ddSJed Brown { 715af33a6ddSJed Brown CharacteristicPointDA2D temp; 716af33a6ddSJed Brown int n; 717af33a6ddSJed Brown 718af33a6ddSJed Brown if (0) { /* Check the order of the queue before sorting */ 719af33a6ddSJed Brown PetscInfo(PETSC_NULL, "Before Heap sort\n"); 720af33a6ddSJed Brown for (n=0; n<size; n++) { 721af33a6ddSJed Brown PetscInfo2(PETSC_NULL,"%d %d\n",n,queue[n].proc); 722af33a6ddSJed Brown } 723af33a6ddSJed Brown } 724af33a6ddSJed Brown 725af33a6ddSJed Brown /* SORTING PHASE */ 726af33a6ddSJed Brown for (n = (size / 2)-1; n >= 0; n--) 727af33a6ddSJed Brown SiftDown(c, queue, n, size-1); /* Rich had size-1 here, Matt had size*/ 728af33a6ddSJed Brown for (n = size-1; n >= 1; n--) { 729af33a6ddSJed Brown temp = queue[0]; 730af33a6ddSJed Brown queue[0] = queue[n]; 731af33a6ddSJed Brown queue[n] = temp; 732af33a6ddSJed Brown SiftDown(c, queue, 0, n-1); 733af33a6ddSJed Brown } 734af33a6ddSJed Brown if (0) { /* Check the order of the queue after sorting */ 735af33a6ddSJed Brown PetscInfo(PETSC_NULL, "Avter Heap sort\n"); 736af33a6ddSJed Brown for (n=0; n<size; n++) { 737af33a6ddSJed Brown PetscInfo2(PETSC_NULL,"%d %d\n",n,queue[n].proc); 738af33a6ddSJed Brown } 739af33a6ddSJed Brown } 740af33a6ddSJed Brown return 0; 741af33a6ddSJed Brown } 742af33a6ddSJed Brown 743af33a6ddSJed Brown /*---------------------------------------------------------------------*/ 744af33a6ddSJed Brown #undef __FUNCT__ 745af33a6ddSJed Brown #define __FUNCT__ "SiftDown" 746af33a6ddSJed Brown /* 747af33a6ddSJed Brown Based on code from http://linux.wku.edu/~lamonml/algor/sort/heap.html 748af33a6ddSJed Brown */ 749af33a6ddSJed Brown PetscErrorCode SiftDown(Characteristic c, Queue queue, PetscInt root, PetscInt bottom) 750af33a6ddSJed Brown /*---------------------------------------------------------------------*/ 751af33a6ddSJed Brown { 752af33a6ddSJed Brown PetscBool done = PETSC_FALSE; 753af33a6ddSJed Brown PetscInt maxChild; 754af33a6ddSJed Brown CharacteristicPointDA2D temp; 755af33a6ddSJed Brown 756af33a6ddSJed Brown PetscFunctionBegin; 757af33a6ddSJed Brown while ((root*2 <= bottom) && (!done)) { 758af33a6ddSJed Brown if (root*2 == bottom) maxChild = root * 2; 759af33a6ddSJed Brown else if (queue[root*2].proc > queue[root*2+1].proc) maxChild = root * 2; 760af33a6ddSJed Brown else maxChild = root * 2 + 1; 761af33a6ddSJed Brown 762af33a6ddSJed Brown if (queue[root].proc < queue[maxChild].proc) { 763af33a6ddSJed Brown temp = queue[root]; 764af33a6ddSJed Brown queue[root] = queue[maxChild]; 765af33a6ddSJed Brown queue[maxChild] = temp; 766af33a6ddSJed Brown root = maxChild; 767af33a6ddSJed Brown } else 768af33a6ddSJed Brown done = PETSC_TRUE; 769af33a6ddSJed Brown } 770af33a6ddSJed Brown PetscFunctionReturn(0); 771af33a6ddSJed Brown } 772af33a6ddSJed Brown 773af33a6ddSJed Brown #undef __FUNCT__ 774af33a6ddSJed Brown #define __FUNCT__ "DMDAGetNeighborsRank" 775af33a6ddSJed Brown /* [center, left, top-left, top, top-right, right, bottom-right, bottom, bottom-left] */ 776af33a6ddSJed Brown PetscErrorCode DMDAGetNeighborsRank(DM da, PetscMPIInt neighbors[]) 777af33a6ddSJed Brown { 778af33a6ddSJed Brown DMDABoundaryType bx, by; 779af33a6ddSJed Brown PetscBool IPeriodic = PETSC_FALSE, JPeriodic = PETSC_FALSE; 780af33a6ddSJed Brown MPI_Comm comm; 781af33a6ddSJed Brown PetscMPIInt rank; 782af33a6ddSJed Brown PetscInt **procs,pi,pj,pim,pip,pjm,pjp,PI,PJ; 783af33a6ddSJed Brown PetscErrorCode ierr; 784af33a6ddSJed Brown 785af33a6ddSJed Brown PetscFunctionBegin; 786af33a6ddSJed Brown ierr = PetscObjectGetComm((PetscObject) da, &comm);CHKERRQ(ierr); 787af33a6ddSJed Brown ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 788af33a6ddSJed Brown ierr = DMDAGetInfo(da, 0, 0, 0, 0, &PI,&PJ, 0, 0, 0, &bx, &by,0, 0); 789af33a6ddSJed Brown 790af33a6ddSJed Brown if (bx == DMDA_BOUNDARY_PERIODIC) { 791af33a6ddSJed Brown IPeriodic = PETSC_TRUE; 792af33a6ddSJed Brown } 793af33a6ddSJed Brown if (by == DMDA_BOUNDARY_PERIODIC) { 794af33a6ddSJed Brown JPeriodic = PETSC_TRUE; 795af33a6ddSJed Brown } 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 */ 863af33a6ddSJed Brown if (jr >= js && jr <= je) { 864af33a6ddSJed Brown return 0; 865af33a6ddSJed Brown } else if (jr < js) { 866af33a6ddSJed Brown return 7; 867af33a6ddSJed Brown } else { 868af33a6ddSJed Brown return 3; 869af33a6ddSJed Brown } 870af33a6ddSJed Brown } else if (ir < is) { /* left column */ 871af33a6ddSJed Brown if (jr >= js && jr <= je) { 872af33a6ddSJed Brown return 1; 873af33a6ddSJed Brown } else if (jr < js) { 874af33a6ddSJed Brown return 8; 875af33a6ddSJed Brown } else { 876af33a6ddSJed Brown return 2; 877af33a6ddSJed Brown } 878af33a6ddSJed Brown } else { /* right column */ 879af33a6ddSJed Brown if (jr >= js && jr <= je) { 880af33a6ddSJed Brown return 5; 881af33a6ddSJed Brown } else if (jr < js) { 882af33a6ddSJed Brown return 6; 883af33a6ddSJed Brown } else { 884af33a6ddSJed Brown return 4; 885af33a6ddSJed Brown } 886af33a6ddSJed Brown } 887af33a6ddSJed Brown } 888