00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024 #include "map.h"
00025
00026
00027
00028 static void checkEquality(const char *what,
00029 int v1,const char *src1,
00030 int v2,const char *src2)
00031 {
00032 if (v1!=v2) {
00033 CkPrintf("FEM> %s value %d, from %d, doesn't match %d, from %d!\n",
00034 what, v1,src1, v2,src2);
00035 CkAbort("FEM Equality assertation failed");
00036 }
00037 }
00038 static void checkRange(const char *what,int v,int max)
00039 {
00040 if ((v<0)||(v>=max)) {
00041 CkPrintf("FEM> %s value %d should be between 0 and %d!\n",
00042 what,v,max);
00043 CkAbort("FEM Range assertation failed");
00044 }
00045 }
00046
00047 static void checkArrayEntries(const int *arr,int nArr,int max,const char *what)
00048 {
00049 #if CMK_ERROR_CHECKING
00050
00051 for (int e=0;e<nArr;e++) checkRange(what,arr[e],max);
00052 #endif
00053 }
00054
00055
00056 static void check(const FEM_Mesh *mesh) {
00057 #if 0 //FIXME
00058 for (int t=0;t<mesh->elem.size();t++)
00059 if (mesh->elem.has(t)) {
00060
00061 checkArrayEntries(mesh->elem[t].getConn().getData(),
00062 mesh->elem[t].size()*mesh->elem[t].getNodesPer(),
00063 mesh->node.size(), "element connectivity, from FEM_Set_Elem_Conn,");
00064 }
00065 for (int s=0;s<mesh->nSparse();s++) {
00066
00067 const FEM_Sparse &src=mesh->getSparse(s);
00068 checkArrayEntries(src.getNodes(0),src.size()*src.getNodesPer(),
00069 mesh->node.size(), "sparse data nodes, from FEM_Set_Sparse,");
00070 if (src.getElem())
00071 checkArrayEntries(src.getElem(),src.size()*2,
00072 mesh->nElems(),"sparse data elements, from FEM_Set_Sparse_Elem,");
00073 }
00074 #endif
00075 }
00076 static void check(const FEM_Partition &partition,const FEM_Mesh *mesh) {
00077 const int *canon=partition.getCanon();
00078 if (canon!=NULL)
00079 checkArrayEntries(canon,mesh->node.size(),mesh->node.size(),
00080 "node canonicalization array, from FEM_Set_Symmetries");
00081 }
00082
00084 void FEM_Ghost_Stencil::check(const FEM_Mesh &mesh) const {
00085 const FEM_Elem &elem=mesh.elem[mesh.chkET(elType)];
00086 checkEquality("number of elements",
00087 n,"FEM_Ghost_stencil",
00088 elem.size(),"mesh");
00089 int i,prevEnd=0;
00090 for (i=0;i<n;i++) {
00091 int nElts=ends[i]-prevEnd;
00092 checkRange("FEM_Ghost_stencil index array",nElts,n);
00093 prevEnd=ends[i];
00094 }
00095 int m=ends[n-1];
00096 for (i=0;i<m;i++) {
00097 int t=adj[2*i+0];
00098 mesh.chkET(t);
00099 checkRange("FEM_Ghost_stencil neighbor array",
00100 adj[2*i+1],mesh.elem[t].size());
00101 }
00102 }
00103
00104
00105
00117 class dynChunk {
00118 public:
00119 typedef CkVec<int> dynList;
00120 NumberedVec<dynList> elem;
00121 dynList node;
00122
00123
00124 int addRealElement(int type,int globalNo) {
00125 elem[type].push_back(globalNo);
00126 return elem[type].size()-1;
00127 }
00128 int addRealNode(int globalNo) {
00129 node.push_back(globalNo);
00130 return node.size()-1;
00131 }
00132 };
00133
00138 class splitter {
00139
00140 FEM_Mesh *mesh;
00141 const int *elem2chunk;
00142 int nChunks;
00143
00144 FEM_Mesh **chunks;
00145
00146
00147 chunkList *gNode;
00148 CkVec<chunkList *> gElem;
00149 dynChunk *dyn;
00150
00155 void renumberNodesLocal(int row, BasicTable2d<int> &table,
00156 int chunk, FEM_Symmetries_t sym)
00157 {
00158 int *nodes=table.getRow(row);
00159 int nNodes=table.width();
00160 for (int i=0;i<nNodes;i++)
00161 nodes[i]=gNode[nodes[i]].localOnChunk(chunk,sym);
00162 }
00163
00164
00165 FEM_Sparse **sparseDest;
00166
00168 void copySparse(const FEM_Sparse &src,int s,int chunk,FEM_Symmetries_t sym) {
00169 FEM_Sparse *dest=sparseDest[chunk];
00170 int d=dest->push_back(src,s);
00171 renumberNodesLocal(d, dest->setConn(), chunk, sym);
00172
00173 if (dest->hasElements()) {
00174
00175
00176
00177 }
00178 }
00179
00181 void copySparseChunks(const FEM_Sparse &src,int s,bool forGhost);
00182
00183
00184
00185 unsigned char *ghostNode;
00186 const int *canon;
00187 const FEM_Symmetries_t *sym;
00188 int curGhostLayerNo;
00189 int totGhostElem,totGhostNode;
00190
00192 void addStencil(const FEM_Ghost_Stencil &s,const FEM_Partition &partition);
00193
00195 void addLayer(const FEM_Ghost_Layer &g,const FEM_Partition &partition);
00196
00198 bool hasGhostNodes(const int *conn,int nodesPer)
00199 {
00200 for (int i=0;i<nodesPer;i++)
00201 if (ghostNode[conn[i]])
00202 return true;
00203 return false;
00204 }
00205
00207 bool addTuple(int *dest,FEM_Symmetries_t *destSym,
00208 const int *elem2tuple,int nodesPerTuple,const int *conn) const;
00209
00211 void addSymmetryGhost(const elemList &a);
00213 void addGlobalGhost(int srcType,int srcNum, int destType,int destNum, bool addNodes);
00214
00217 void addGhostPair(const elemList &src,const elemList &dest,bool addNodes);
00218
00220 int addGhostElement(int t,int gNo,int srcChunk,int destChunk,FEM_Symmetries_t sym);
00221
00223 int addGhostNode(int gnNo,int srcChunk,int destChunk,FEM_Symmetries_t sym);
00224
00226 int addGhostInner(const FEM_Entity &gEnt,int gNo, chunkList &gDest,
00227 int srcChunk,FEM_Entity &srcEnt, int destChunk,FEM_Entity &destEnt,
00228 FEM_Symmetries_t sym, int isNode, int t);
00229
00231 void addElemElems(const FEM_Partition &partition);
00232 void buildElemElemData(const FEM_ElemAdj_Layer &g,const FEM_Partition &partition);
00233
00234
00235
00236
00237 void bad(const char *why) {
00238 CkAbort(why);
00239 }
00240 void equal(int is,int should,const char *what) {
00241 if (is!=should) {
00242 CkPrintf("ERROR! Expected %s to be %d, got %d!\n",
00243 what,should,is);
00244 bad("Internal FEM data structure corruption! (unequal)\n");
00245 }
00246 }
00247 void range(int value,int lo,int hi,const char *what) {
00248 if (value<lo || value>=hi) {
00249 CkPrintf("ERROR! %s: %d is out of range (%d,%d)!\n",what,value,lo,hi);
00250 bad("Internal FEM data structure corruption! (out of range)\n");
00251 }
00252 }
00253 void nonnegative(int value,const char *what) {
00254 if (value<0) {
00255 CkPrintf("ERROR! Expected %s to be non-negative, got %d!\n",
00256 what,value);
00257 bad("Internal FEM data structure corruption! (negative)\n");
00258 }
00259 }
00260
00261 public:
00262 splitter(FEM_Mesh *mesh_,const int *elem2chunk_,int nChunks_);
00263 ~splitter();
00264
00265
00266 void buildCommLists(void);
00267
00268
00269 void addGhosts(const FEM_Partition &partition);
00270
00271
00272 void separateSparse(bool forGhost);
00273
00274
00275 void aboutToCreate(void);
00276
00277
00278 FEM_Mesh *createMesh(int c);
00279
00280 void consistencyCheck(void);
00281
00282 };
00283
00284
00285 splitter::splitter(FEM_Mesh *mesh_,const int *elem2chunk_,int nChunks_)
00286 :mesh(mesh_), elem2chunk(elem2chunk_),nChunks(nChunks_)
00287 {
00288 chunks=new FEM_Mesh* [nChunks];
00289 dyn=new dynChunk[nChunks];
00290 int c;
00291 for (c=0;c<nChunks;c++) {
00292 chunks[c]=new FEM_Mesh;
00293 chunks[c]->copyShape(*mesh);
00294 dyn[c].elem.makeLonger(mesh->elem.size()-1);
00295 }
00296
00297 gNode=new chunkList[mesh->node.size()];
00298 gElem.resize(mesh->elem.size());
00299 for (int t=0;t<mesh->elem.size();t++) {
00300 if (mesh->elem.has(t))
00301 gElem[t]=new chunkList[mesh->elem[t].size()];
00302 else
00303 gElem[t]=NULL;
00304 }
00305
00306 sparseDest=new FEM_Sparse* [nChunks];
00307 ghostNode=NULL;
00308 canon=NULL;
00309 }
00310
00311 splitter::~splitter() {
00312 delete[] gNode;
00313 delete[] dyn;
00314 for (int t=0;t<mesh->elem.size();t++)
00315 delete[] gElem[t];
00316 for (int c=0;c<nChunks;c++)
00317 delete chunks[c];
00318 delete[] chunks;
00319 delete[] sparseDest;
00320 }
00321
00322 void splitter::buildCommLists(void)
00323 {
00324
00325
00326 int t;
00327 int e;
00328 int n;
00329 for (t=0;t<mesh->elem.size();t++)
00330 if (mesh->elem.has(t)) {
00331 int typeStart=mesh->nElems(t);
00332 const FEM_Elem &src=mesh->elem[t];
00333 for (e=0;e<src.size();e++) {
00334 int c=elem2chunk[typeStart+e];
00335 const int *srcConn=src.getConn().getRow(e);
00336 gElem[t][e].set(c,dyn[c].addRealElement(t,e),noSymmetries,-1);
00337 int lNo=dyn[c].node.size();
00338 for (n=0;n<src.getNodesPer();n++) {
00339 int gNo=srcConn[n];
00340 if (gNode[gNo].addchunk(c,lNo,noSymmetries,-1)) {
00341
00342 dyn[c].addRealNode(gNo);
00343 lNo++;
00344 }
00345 }
00346 }
00347 }
00348
00349
00350 for (n=0;n<mesh->node.size();n++) {
00351 if (gNode[n].isShared())
00352 {
00353
00354 int len=gNode[n].length();
00355 for (int bi=0;bi<len;bi++)
00356 for (int ai=0;ai<bi;ai++)
00357 {
00358 chunkList &a=gNode[n][ai];
00359 chunkList &b=gNode[n][bi];
00360 chunks[a.chunk]->node.shared.add(
00361 a.chunk,a.localNo,
00362 b.chunk,b.localNo,
00363 chunks[b.chunk]->node.shared);
00364 }
00365 }
00366 }
00367 }
00368
00369
00370
00371 void splitter::aboutToCreate()
00372 {
00373 for (int t=0;t<mesh->elem.size();t++)
00374 {delete[] gElem[t];gElem[t]=NULL;}
00375 }
00376
00377
00378
00379
00380 FEM_Mesh *splitter::createMesh(int c)
00381 {
00382 int t;
00383 FEM_Mesh &dest=*chunks[c];
00384
00385 dest.udata=mesh->udata;
00386
00387
00388 FEM_Node *destNode=&dest.node;
00389 const dynChunk::dynList &srcIdx=dyn[c].node;
00390 int nNode=srcIdx.size();
00391 destNode->setLength(nNode);
00392 for (int lNo=0;lNo<nNode;lNo++) {
00393 int gNo=srcIdx[lNo];
00394
00395 destNode->copyEntity(lNo, mesh->node, gNo);
00396 destNode->setSymmetries(lNo,noSymmetries);
00397
00398 chunkList *head=&gNode[gNo];
00399 const chunkList *cur=head->onChunk(c,noSymmetries);
00400 destNode->setPrimary(lNo,head==cur);
00401 }
00402
00403
00404 for (t=0;t<dest.elem.size();t++)
00405 if (dest.elem.has(t))
00406 {
00407 FEM_Elem *destElem=&dest.elem[t];
00408 const dynChunk::dynList &srcIdx=dyn[c].elem[t];
00409 int nEl=srcIdx.size();
00410 destElem->setLength(nEl);
00411 for (int lNo=0;lNo<nEl;lNo++) {
00412 int gNo=srcIdx[lNo];
00413
00414 destElem->copyEntity(lNo, mesh->elem[t], gNo);
00415
00416 renumberNodesLocal(lNo, destElem->setConn(), c,noSymmetries);
00417 }
00418 }
00419
00420 chunks[c]=NULL;
00421 dest.becomeGetting();
00422 return &dest;
00423 }
00424
00425
00426
00427
00428
00429 void FEM_Mesh_split(FEM_Mesh *mesh,int nChunks,
00430 const FEM_Partition &partition,FEM_Mesh_Output *out)
00431 {
00432 const int *elem2chunk=partition.getPartition(mesh,nChunks);
00433
00434 CkThresholdTimer time("FEM Split> Setting up",1.0);
00435 checkArrayEntries(elem2chunk,mesh->nElems(),nChunks,
00436 "elem2chunk, from FEM_Set_Partition or metis,");
00437 check(mesh);
00438 check(partition,mesh);
00439
00440 mesh->setSymList(partition.getSymList());
00441 splitter s(mesh,elem2chunk,nChunks);
00442
00443 time.start("FEM Split> Finding comm lists");
00444 s.buildCommLists();
00445
00446 time.start("FEM Split> Finding ghosts");
00447 s.separateSparse(false);
00448 s.addGhosts(partition);
00449 s.separateSparse(true);
00450
00451 time.start("FEM Split> Copying mesh data");
00452
00453 s.aboutToCreate();
00454 for (int c=0;c<nChunks;c++)
00455 out->accept(c,s.createMesh(c));
00456 }
00457
00458
00459
00461 void splitter::copySparseChunks(const FEM_Sparse &src,int s,bool forGhost)
00462 {
00463 if (src.hasElements())
00464 {
00465 const int *elemRec=src.getElem().getRow(s);
00466 int t=elemRec[0];
00467 int e=elemRec[1];
00468 if (!forGhost)
00469 {
00470 copySparse(src,s,elem2chunk[mesh->getGlobalElem(t,e)],noSymmetries);
00471 } else {
00472 chunkList *cur=&gElem[t][e];
00473 while (cur!=NULL) {
00474 if (FEM_Is_ghost_index(cur->localNo))
00475 copySparse(src,s,cur->chunk,cur->sym);
00476 cur=cur->next;
00477 }
00478 }
00479 }
00480 else
00481 {
00482 int nNodes=src.getConn().width();
00483 if (nNodes==0) CkAbort("Registered an FEM sparse without nodes or elements");
00484 FEM_Symmetries_t sym=noSymmetries;
00485 const int *nodes=src.getConn().getRow(s);
00486 int nCopied=0;
00487
00488
00489 for (const chunkList *cur=&gNode[nodes[0]];cur!=NULL;cur=cur->next)
00490 {
00491 int testchunk=cur->chunk;
00492 bool hasGhost=false;
00493 for (int i=0;i<nNodes;i++) {
00494 const chunkList *li=gNode[nodes[i]].onChunk(testchunk,sym);
00495 if (li==NULL)
00496 { testchunk=-1; break; }
00497 if (FEM_Is_ghost_index(li->localNo))
00498 hasGhost=true;
00499 }
00500 if (testchunk==-1) continue;
00501 if (forGhost && !hasGhost) continue;
00502 copySparse(src,s,testchunk,sym);
00503 nCopied++;
00504 }
00505 if (nCopied==0 && !forGhost) FEM_Abort("copySparseChunks",
00506 "Sparse record %d does not lie on any chunk!",s);
00507 }
00508 }
00509
00510
00511
00512 void splitter::separateSparse(bool forGhost)
00513 {
00514 for (int t=0;t<mesh->sparse.size();t++)
00515 if (mesh->sparse.has(t))
00516 {
00517
00518 const FEM_Sparse &src=mesh->sparse.get(t);
00519 for (int c=0;c<nChunks;c++)
00520 {
00521 FEM_Sparse *d=&chunks[c]->sparse.set(t);
00522 if (forGhost) d=(FEM_Sparse *)d->getGhost();
00523 sparseDest[c]=d;
00524 }
00525
00526 for (int s=0;s<src.size();s++) copySparseChunks(src,s,forGhost);
00527 }
00528 }
00529
00530
00531
00532
00533
00534
00535
00536
00537
00538
00539
00540
00541
00542
00543
00544
00545
00546
00547
00548
00549
00550
00551
00552
00553
00554
00555
00556 CkHashCode CkHashFunction_ints(const void *keyData,size_t keyLen)
00557 {
00558 const int *d=(const int *)keyData;
00559 int l=keyLen/sizeof(int);
00560 CkHashCode ret=d[0];
00561 for (int i=1;i<l;i++)
00562 ret=ret^circleShift(d[i],i*23);
00563 return ret;
00564 }
00565
00566 int CkHashCompare_ints(const void *k1,const void *k2,size_t keyLen)
00567 {
00568 const int *d1=(const int *)k1;
00569 const int *d2=(const int *)k2;
00570 int l=keyLen/sizeof(int);
00571 for (int i=0;i<l;i++) if (d1[i]!=d2[i]) return 0;
00572 return 1;
00573 }
00574
00575
00576 extern "C" int ck_fem_map_compare_int(const void *a, const void *b)
00577 {
00578 return (*(const int *)a)-(*(const int *)b);
00579 }
00580
00581
00582 void splitter::addGhosts(const FEM_Partition &partition)
00583 {
00584 if (partition.getRegions()==0) return;
00585
00586
00587 ghostNode=new unsigned char[mesh->node.size()];
00588 int n,nNode=mesh->node.size();
00589 for (n=0;n<nNode;n++) {
00590 ghostNode[n]=(gNode[n].isShared());
00591 }
00592
00593
00594 canon=partition.getCanon();
00595 sym=partition.getSymmetries();
00596 if (sym!=NULL)
00597 for (n=0;n<nNode;n++)
00598 if (sym[n]!=(FEM_Symmetries_t)0)
00599 ghostNode[n]=1;
00600
00601
00602 curGhostLayerNo=0;
00603 for (int regNo=0;regNo<partition.getRegions();regNo++)
00604 {
00605 const FEM_Ghost_Region &r=partition.getRegion(regNo);
00606 totGhostElem=0,totGhostNode=0;
00607 int thisLayer=curGhostLayerNo;
00608
00609 if (r.layer!=0) {
00610 addLayer(*r.layer,partition);
00611 curGhostLayerNo++;
00612 }
00613 else if (r.stencil!=0)
00614 {
00615 addStencil(*r.stencil,partition);
00616 }
00617 else {
00618 curGhostLayerNo++;
00619 }
00620
00621 if (totGhostElem>0)
00622 CkPrintf("FEM Ghost %s %d> %d new ghost elements, %d new ghost nodes\n",
00623 r.layer?"layer":"stencil",
00624 1+thisLayer,
00625 totGhostElem,totGhostNode);
00626 }
00627
00628 delete[] ghostNode; ghostNode=NULL;
00629 }
00630
00631
00632 void splitter::addStencil(const FEM_Ghost_Stencil &s,const FEM_Partition &partition)
00633 {
00634 int t=s.getType();
00635 s.check(*mesh);
00636 int i,j,n=mesh->elem[t].size();
00637 for (i=0;i<n;i++) {
00638 const int *nbor;
00639 j=0;
00640 while (NULL!= (nbor=s.getNeighbor(i,j++))) {
00641
00642
00643 addGlobalGhost(nbor[0],nbor[1],t,i, s.wantNodes());
00644 }
00645 }
00646 }
00647
00650 void splitter::addGlobalGhost(
00651 int srcType,int srcNum, int destType,int destNum, bool doAddNodes)
00652 {
00653
00654 chunkList &srcC=gElem[srcType][srcNum];
00655
00656 for (chunkList *destC=&gElem[destType][destNum];
00657 destC!=NULL;
00658 destC=destC->next)
00659 if (srcC.chunk!=destC->chunk
00660 && destC->layerNo<curGhostLayerNo
00661 )
00662 {
00663 elemList src(srcC.chunk, srcC.localNo, srcType, (FEM_Symmetries_t)0);
00664 elemList dest(destC->chunk, destC->localNo, destType, (FEM_Symmetries_t)0);
00665 addGhostPair(src,dest,doAddNodes);
00666 }
00667 }
00668
00669
00670 void splitter::addLayer(const FEM_Ghost_Layer &g,const FEM_Partition &partition)
00671 {
00672 tupleTable table(g.nodesPerTuple);
00673
00674
00675 for (int t=0;t<mesh->elem.size();t++)
00676 if (mesh->elem.has(t)) {
00677 if (!g.elem[t].add) continue;
00678
00679 int gElemCount=mesh->elem[t].size();
00680 for (int gNo=0;gNo<gElemCount;gNo++)
00681 {
00682 const int *conn=mesh->elem[t].connFor(gNo);
00683 if (hasGhostNodes(conn,mesh->elem[t].getNodesPer()))
00684 {
00685 for (chunkList *cur=&gElem[t][gNo];cur!=NULL;cur=cur->next)
00686 for (int u=0;u<g.elem[t].tuplesPerElem;u++)
00687 {
00688 int tuple[tupleTable::MAX_TUPLE];
00689 FEM_Symmetries_t allSym;
00690 if (addTuple(tuple,&allSym,
00691 &g.elem[t].elem2tuple[u*g.nodesPerTuple],
00692 g.nodesPerTuple,conn)) {
00693 table.addTuple(tuple,new elemList(cur->chunk,cur->localNo,t,allSym));
00694 }
00695 }
00696 }
00697 }
00698 }
00699
00700
00701 table.beginLookup();
00702 elemList *l;
00703 while (NULL!=(l=table.lookupNext())) {
00704 if (l->next==NULL)
00705 addSymmetryGhost(*l);
00706 else {
00707
00708 for (const elemList *a=l;a!=NULL;a=a->next)
00709 for (const elemList *b=l;b!=NULL;b=b->next)
00710 if (a!=b && a->localNo>=0)
00711 addGhostPair(*a,*b,g.addNodes);
00712 }
00713 }
00714 }
00715
00716
00717
00718
00719 bool splitter::addTuple(int *dest,FEM_Symmetries_t *destSym,const int *elem2tuple,
00720 int nodesPerTuple,const int *conn) const
00721 {
00722 FEM_Symmetries_t allSym=(FEM_Symmetries_t)(~0);
00723 for (int i=0;i<nodesPerTuple;i++) {
00724 int eidx=elem2tuple[i];
00725 if (eidx==-1) {
00726 dest[i]=-1;
00727 } else {
00728 int n=conn[eidx];
00729 if (!ghostNode[n])
00730 return false;
00731 if (sym!=NULL) {
00732 allSym=allSym & sym[n];
00733 n=canon[n];
00734 }
00735 dest[i]=n;
00736 }
00737 }
00738
00739 if (sym!=NULL) *destSym=allSym; else *destSym=0;
00740 return true;
00741 }
00742
00743
00744 void splitter::addSymmetryGhost(const elemList &a)
00745 {
00746 if (a.sym==0) return;
00747 CkAbort("FEM map> Mirror symmetry ghosts not yet supported");
00748 }
00749
00750
00759 void splitter::addGhostPair(const elemList &src,const elemList &dest,bool addNodes)
00760 {
00761
00762 int srcChunk=src.chunk;
00763 int destChunk=dest.chunk;
00764 int elemSym=dest.sym;
00765 bool isSymmetry=(dest.sym!=0);
00766
00767 if ((!isSymmetry) && srcChunk==destChunk)
00768 return;
00769
00770 int t=src.type;
00771 if (src.localNo<0) FEM_Abort("addGhostPair","Cannot add a ghost of a ghost (src num=%d)",src.localNo);
00772 int gNo=dyn[srcChunk].elem[t][src.localNo];
00773
00774 int newNo=addGhostElement(t,gNo,srcChunk,destChunk,elemSym);
00775 if (newNo==-1) return;
00776
00777
00778
00779 int srcConnN=mesh->elem[t].getNodesPer();
00780 const int *srcConn=mesh->elem[t].connFor(gNo);
00781
00782 int dt=dest.type;
00783 int destConnN=mesh->elem[dt].getNodesPer();
00784 const int *destConn=NULL;
00785 FEM_Elem *destElemGhosts=(FEM_Elem *)chunks[destChunk]->elem[t].getGhost();
00786 int *newConn=destElemGhosts->connFor(newNo);
00787 if (isSymmetry) {
00788 destConn=mesh->elem[dt].connFor(dyn[destChunk].elem[dt][dest.localNo]);
00789 }
00790
00791
00792 for (int sn=0;sn<srcConnN;sn++)
00793 {
00794 FEM_Symmetries_t nodeSym=noSymmetries;
00795 int gnNo=srcConn[sn];
00796 ghostNode[gnNo]=1;
00797 if (isSymmetry)
00798 {
00799
00800
00801 nodeSym=dest.sym;
00802 for (int dn=0;dn<destConnN;dn++)
00803 if (canon[destConn[dn]]==canon[gnNo])
00804 {
00805 gnNo=destConn[dn];
00806
00807 nodeSym=noSymmetries;
00808 break;
00809 }
00810 }
00811 if (addNodes)
00812 addGhostNode(gnNo,srcChunk,destChunk,nodeSym);
00813
00814
00815 newConn[sn]=gNode[gnNo].localOnChunk(destChunk,nodeSym);
00816 }
00817 }
00818
00819
00821 int splitter::addGhostElement(int t,int gNo,int srcChunk,int destChunk,FEM_Symmetries_t sym)
00822 {
00823 FEM_Elem &destEnt=*(FEM_Elem *)chunks[destChunk]->elem[t].getGhost();
00824 int destNo=addGhostInner(mesh->elem[t],gNo,gElem[t][gNo],
00825 srcChunk,chunks[srcChunk]->elem[t],
00826 destChunk, destEnt,
00827 sym, 0, t);
00828 if (destNo!=-1)
00829 {
00830 totGhostElem++;
00831 return destNo;
00832 }
00833 return -1;
00834 }
00836 int splitter::addGhostNode(int gnNo,int srcChunk,int destChunk,FEM_Symmetries_t sym)
00837 {
00838 int destNo=addGhostInner(mesh->node,gnNo,gNode[gnNo],
00839 srcChunk,chunks[srcChunk]->node,
00840 destChunk,*chunks[destChunk]->node.getGhost(),
00841 sym, 1, 0);
00842 if (destNo!=-1)
00843 {
00844 totGhostNode++;
00845 return destNo;
00846 }
00847 return -1;
00848 }
00849
00860 int splitter::addGhostInner(const FEM_Entity &gEnt,int gNo, chunkList &gDest,
00861 int srcChunk,FEM_Entity &srcEnt, int destChunk,FEM_Entity &destEnt,
00862 FEM_Symmetries_t sym, int isNode, int t)
00863 {
00864 if (gDest.onChunk(destChunk,sym))
00865 return -1;
00866 int srcNo=gDest.localOnChunk(srcChunk,noSymmetries);
00867 if (srcNo<0)
00868 return -1;
00869
00870
00871
00872 int destNo=destEnt.push_back(gEnt,gNo);
00873 destEnt.setSymmetries(destNo,sym);
00874 gDest.addchunk(destChunk,FEM_To_ghost_index(destNo),sym,curGhostLayerNo);
00875 chunkList *l = &gDest;
00876 while(l->next != NULL) {
00877 if(l->chunk != srcChunk && l->chunk != destChunk && l->localNo>=0) {
00878 if(isNode)
00879 chunks[l->chunk]->node.setGhostSend().add(l->chunk,l->localNo,destChunk,destNo,destEnt.setGhostRecv());
00880 else
00881 chunks[l->chunk]->elem[t].setGhostSend().add(l->chunk,l->localNo,destChunk,destNo,destEnt.setGhostRecv());
00882 }
00883 l = l->next;
00884 }
00885 srcEnt.setGhostSend().add(srcChunk,srcNo,destChunk,destNo,
00886 destEnt.setGhostRecv());
00887 return destNo;
00888 }
00889
00890
00891
00892 void splitter::consistencyCheck(void)
00893 {
00894 #if CMK_ERROR_CHECKING
00895 bool skipCheck=false;
00896 #else
00897 bool skipCheck=true;
00898 #endif
00899 if (skipCheck) return;
00900
00901 CkThresholdTimer time("FEM Split> Performing consistency check",1.0);
00902
00903 int t,c;
00904 FEM_Symmetries_t sym=noSymmetries;
00905
00906 for (c=0;c<nChunks;c++) {
00907 for (t=0;t<mesh->elem.size();t++)
00908 if (mesh->elem.has(t)) {
00909 const dynChunk::dynList &srcIdx=dyn[c].elem[t];
00910 int nEl=srcIdx.size();
00911 for (int lNo=0;lNo<nEl;lNo++) {
00912 int gNo=srcIdx[lNo];
00913 range(gNo,0,mesh->elem[t].size(),"global element number");
00914 equal(gElem[t][gNo].localOnChunk(c,sym),lNo,"gElem[t] local number");
00915 }
00916 }
00917 const dynChunk::dynList &srcIdx=dyn[c].node;
00918 int nNo=srcIdx.size();
00919 for (int lNo=0;lNo<nNo;lNo++) {
00920 int gNo=srcIdx[lNo];
00921 range(gNo,0,mesh->node.size(),"global node number");
00922 equal(gNode[gNo].localOnChunk(c,sym),lNo,"gNode[] local number");
00923 }
00924 }
00925
00926
00927 for (int gNo=0;gNo<mesh->node.size();gNo++) {
00928 for (chunkList *l=&gNode[gNo];l!=NULL;l=l->next) {
00929 if (l->chunk<0)
00930 FEM_Abort("partitioning","Global node %d (0-based) is not connected to any elements!",gNo);
00931 range(l->chunk,0,nChunks,"gNode chunk");
00932 if (!FEM_Is_ghost_index(l->localNo))
00933 equal(dyn[l->chunk].node[l->localNo],gNo,"chunk local node");
00934 }
00935 }
00936
00937 for (t=0;t<mesh->elem.size();t++)
00938 if (mesh->elem.has(t)) {
00939 for (int gNo=0;gNo<mesh->elem[t].size();gNo++) {
00940 for (chunkList *l=&gElem[t][gNo];l!=NULL;l=l->next) {
00941 range(l->chunk,0,nChunks,"gElem chunk");
00942 if (!FEM_Is_ghost_index(l->localNo))
00943 equal(dyn[l->chunk].elem[t][l->localNo],gNo,"chunk local element");
00944 }
00945 }
00946 }
00947
00948
00949 }
00950
00951
00952
00953
00954
00955
00956
00957
00958
00959
00960
00961
00962
00964 class FEM_Entity_numberer {
00965 CkVec<unsigned char> occupied;
00966 int occupiedBefore;
00967
00968
00969 int nextUnoccupied(void) {
00970 while (occupiedBefore<occupied.size()) {
00971 if (occupied[occupiedBefore]==0)
00972 return occupiedBefore;
00973 else occupiedBefore++;
00974 }
00975
00976 occupied.push_back(1);
00977 return occupiedBefore;
00978 }
00979
00980 public:
00981 FEM_Entity_numberer() { occupiedBefore=0; }
00982
00984 void mark(FEM_Entity &src) {
00985 int l,len=src.size();
00986 for (l=0;l<len;l++) {
00987 int g=src.getGlobalno(l);
00988 if (g!=-1) {
00989 while (occupied.size()<=g) {
00990 occupied.push_back(0);
00991 }
00992
00993
00994 occupied[g]=1;
00995 }
00996 }
00997 }
00998
01000 void add(FEM_Entity &l_src,FEM_Entity &g_dest) {
01001 int l,len=l_src.size();
01002 for (l=0;l<len;l++) add(l_src,l,g_dest);
01003 }
01004
01006 int add(FEM_Entity &l_src,int l,FEM_Entity &g_dest) {
01007 int g=l_src.getGlobalno(l);
01008 if (g==-1) {
01009 g=nextUnoccupied();
01010 occupied[g]=1;
01011
01012 l_src.setGlobalno(l,g);
01013 }
01014 if (g_dest.size()<=g) g_dest.setLength(g+1);
01015 g_dest.copyEntity(g,l_src,l);
01016 g_dest.setGlobalno(g,g);
01017 return g;
01018 }
01019 };
01020
01021
01022 static void renumberConn(const FEM_Elem &src_e,int l,FEM_Elem &dest_e,int g,
01023 const FEM_Mesh &mesh)
01024 {
01025 const int *srcConn=src_e.connFor(l);
01026 int *dstConn=dest_e.setConn().getRow(g);
01027 for (int n=0;n<src_e.getNodesPer();n++)
01028 dstConn[n]=mesh.node.getGlobalno(srcConn[n]);
01029 }
01030
01031 FEM_Mesh *FEM_Mesh_assemble(int nChunks,FEM_Mesh **chunks)
01032 {
01033 int t,c;
01034 FEM_Mesh *m=new FEM_Mesh;
01035 for(c=0; c<nChunks;c++)
01036 m->copyShape(*chunks[c]);
01037
01038 m->udata=chunks[0]->udata;
01039
01040
01041 FEM_Entity_numberer nodeNum;
01042 for(c=0; c<nChunks;c++) nodeNum.mark(chunks[c]->node);
01043 for(c=0; c<nChunks;c++) nodeNum.add(chunks[c]->node,m->node);
01044
01045
01046 int nElemTypes=m->elem.size();
01047 for (t=0;t<nElemTypes;t++) {
01048 FEM_Entity_numberer elemNum;
01049 for(c=0; c<nChunks;c++)
01050 if (chunks[c]->elem.has(t))
01051 elemNum.mark(chunks[c]->elem[t]);
01052
01053 for(c=0; c<nChunks;c++)
01054 if (chunks[c]->elem.has(t)) {
01055 FEM_Elem &src_e=chunks[c]->elem[t];
01056 if (!m->elem.has(t)) m->elem.set(t).copyShape(src_e);
01057 FEM_Elem &dest_e=m->elem[t];
01058 for (int l=0;l<src_e.size();l++) {
01059 int g=elemNum.add(src_e,l,dest_e);
01060 renumberConn(src_e,l,dest_e,g,*chunks[c]);
01061 }
01062 }
01063 }
01064
01065
01066 int nSparseTypes=m->sparse.size();
01067 for (t=0;t<nSparseTypes;t++) {
01068 FEM_Entity_numberer sparseNum;
01069 for(c=0; c<nChunks;c++)
01070 if (chunks[c]->sparse.has(t))
01071 sparseNum.mark(chunks[c]->sparse[t]);
01072
01073 for(c=0; c<nChunks;c++)
01074 if (chunks[c]->sparse.has(t)) {
01075 FEM_Sparse &src_e=chunks[c]->sparse[t];
01076 if (!m->sparse.has(t)) m->sparse.set(t).copyShape(src_e);
01077 FEM_Sparse &dest_e=m->sparse[t];
01078 for (int l=0;l<src_e.size();l++) {
01079 int g=sparseNum.add(src_e,l,dest_e);
01080 renumberConn(src_e,l,dest_e,g,*chunks[c]);
01081
01082 }
01083 }
01084 }
01085
01086 m->becomeGetting();
01087 return m;
01088 }