cppcheck cleanup sweep for unowned modules
[charm.git] / src / libs / ck-libs / ParFUM / mesh_modify.C
blobbe59bac840bb1296cde359b4fbb0b9808af9d8ff
1 /**
2  * \addtogroup ParFUM
3 */
4 /*@{*/
7 /** File: fem_mesh_modify.C
8  * Authors: Nilesh Choudhury, Isaac Dooley
9  * 
10  * This file contains a set of functions, which allow primitive operations upon meshes in parallel.
11  * See the assumptions listed in fem_mesh_modify.h before using these functions.
12  * 
13  * It also contains a shadow array class that is attached to a chunk and
14  * handles all remote communication between chunks (needed for adaptivity)
15  */
17 #include "ParFUM.h"
18 #include "ParFUM_internals.h"
19 #include "adapt_adj.h"
21 CProxy_femMeshModify meshMod;
24 /** A wrapper to simplify the lookup to determine whether a node is shared
25  */
26 inline int is_shared(FEM_Mesh *m, int node){
27   return m->getfmMM()->getfmUtil()->isShared(node);
31 /* Some helper print functions that can be used for debugging 
32    helper functions that print the various node/element connectivities
33    calls the corresponding functions in the FEM_MUtil class to do the printing */
34 CDECL void FEM_Print_n2n(int mesh, int nodeid){
35   FEM_Mesh *m=FEM_Mesh_lookup(mesh,"FEM_Print_Mesh_Summary");
36   m->getfmMM()->getfmUtil()->FEM_Print_n2n(m, nodeid);
38 CDECL void FEM_Print_n2e(int mesh, int eid){
39   FEM_Mesh *m=FEM_Mesh_lookup(mesh,"FEM_Print_Mesh_Summary");
40   m->getfmMM()->getfmUtil()->FEM_Print_n2e(m, eid);
42 CDECL void FEM_Print_e2n(int mesh, int eid){
43   FEM_Mesh *m=FEM_Mesh_lookup(mesh,"FEM_Print_Mesh_Summary");
44   m->getfmMM()->getfmUtil()->FEM_Print_e2n(m, eid);
46 CDECL void FEM_Print_e2e(int mesh, int eid){
47   FEM_Mesh *m=FEM_Mesh_lookup(mesh,"FEM_Print_Mesh_Summary");
48   m->getfmMM()->getfmUtil()->FEM_Print_e2e(m, eid);
51 /** Prints the mesh summary, i.e. number of valid nodes/elements (local/ghost)
52  */
53 CDECL void FEM_Print_Mesh_Summary(int mesh){
54   CkPrintf("---------------FEM_Print_Mesh_Summary-------------\n");
55   FEM_Mesh *m=FEM_Mesh_lookup(mesh,"FEM_Print_Mesh_Summary");
56   // Print Node information
57   CkPrintf("     Nodes: %d/%d and ghost nodes: %d/%d\n", m->node.count_valid(), m->node.size(),m->node.getGhost()->count_valid(),m->node.getGhost()->size());
58   // Print Element information
59   CkPrintf("     Element Types: %d\n", m->elem.size());
60   for (int t=0;t<m->elem.size();t++) {// for each element type t
61     if (m->elem.has(t)) {
62       int numEl = m->elem[t].size();
63       int numElG = m->elem[t].getGhost()->size();
64       int numValidEl = m->elem[t].count_valid();
65       int numValidElG = m->elem[t].getGhost()->count_valid();
66       CkPrintf("     Element type %d contains %d/%d elements and %d/%d ghosts\n", t, numValidEl, numEl, numValidElG, numElG);
67     }
68   }
69   CkPrintf("\n");
74 //========================Basic Adaptivity operations=======================
75 /*These are the basic add/remove node/element operations */
76 //The following functions translate the calls from meshId to mesh pointer
77 CDECL int FEM_add_node(int mesh, int* adjacent_nodes, int num_adjacent_nodes, int *chunks, int numChunks, int forceShared){
78   return FEM_add_node(FEM_Mesh_lookup(mesh,"FEM_add_node"), adjacent_nodes, num_adjacent_nodes, chunks, numChunks, forceShared);
80 CDECL void FEM_remove_node(int mesh,int node){
81   return FEM_remove_node(FEM_Mesh_lookup(mesh,"FEM_remove_node"), node);
83 CDECL int FEM_add_element(int mesh, int* conn, int conn_size, int elem_type, int chunkNo){
84   return FEM_add_element(FEM_Mesh_lookup(mesh,"FEM_add_element"), conn, conn_size, elem_type, chunkNo);
86 CDECL int FEM_remove_element(int mesh, int element, int elem_type, int permanent){
87   return FEM_remove_element(FEM_Mesh_lookup(mesh,"FEM_remove_element"), element, elem_type, permanent);
89 CDECL int FEM_purge_element(int mesh, int element, int elem_type) {
90   return FEM_purge_element(FEM_Mesh_lookup(mesh,"FEM_remove_element"), element, elem_type);
94 /* Implementation of the add/remove/purge node/element functions 
95    (along with helpers)
96    Add node and helpers -- add_local/ add_remote/ add_shared, etc*/
98 /** Add a node to this mesh 'm', could be ghost/local decided by 'addGhost'
99     basically, this function finds an empty row in the nodelist
100     and clears its adjacencies
101     then, it creates a lock for this node, if this node is new (i.e. unused)
102     otherwise, the existing lock is reset
104 int FEM_add_node_local(FEM_Mesh *m, bool addGhost, bool doLocking, bool doAdjacencies){
105   int newNode;
106   if(addGhost){
107     // find a place for new node in tables, reusing old invalid nodes if possible
108     newNode = m->node.getGhost()->get_next_invalid(m,true,true); 
109         if(doAdjacencies){
110           m->n2e_removeAll(FEM_To_ghost_index(newNode));    // initialize element adjacencies
111           m->n2n_removeAll(FEM_To_ghost_index(newNode));    // initialize node adjacencies
112         }
113   }
114   else{
115     newNode = m->node.get_next_invalid(m,true,false);
116     if(doAdjacencies){
117           m->n2e_removeAll(newNode);    // initialize element adjacencies
118           m->n2n_removeAll(newNode);    // initialize node adjacencies
119         }
120     //add a lock
121     if(doLocking){
122           if(newNode >= m->getfmMM()->fmLockN.size()) {
123                 m->getfmMM()->fmLockN.push_back(FEM_lockN(newNode,m->getfmMM()));
124           }
125           else {
126                 m->getfmMM()->fmLockN[newNode].reset(newNode,m->getfmMM());
127           }
128         }
129   }
130 #ifdef DEBUG_2
131   CkPrintf("[%d]Adding Node %d\n",m->getfmMM()->getfmUtil()->getIdx(),(addGhost==0)?newNode:FEM_To_ghost_index(newNode));
132 #endif
133   return newNode;  // return a new index
137 /** 'adjacentNodes' specifies the set of nodes between which this new node
138     is to be added.
139     'forceshared' is used to convey extra information. Sometimes, we know that
140     the new node will not be shared, i.e. if the two nodes are shared, but the
141     two elements on two sides of this node belong to the same chunk, then this node
142     will not be shared. (in spite of the two nodes being shared by the same set of chunks
143     'numChunks' is the number of entries in 'chunks'.
144     'chunks' contains the chunk index where this node is to be added.
145     This is populated when the operation originates from edge_bisect and the 
146     operation knows which all chunks this node needs to be added to.
147     Returns the index of the new node that was added
148  */
149 int FEM_add_node(FEM_Mesh *m, int* adjacentNodes, int numAdjacentNodes, int*chunks, int numChunks, int forceShared){
150   // add local node
151   //should be used only when all the adjacent nodes are shared 
152   //but you know that the new node should not
153   //be shared, but should belong to only one of the chunks sharing the face/edge.
154   //usually it is safe to use -1, except in some weird cases.
155   int index = m->getfmMM()->getIdx();
156   if(numChunks==1 && chunks[0]!=index) {
157     //translate the indices.. the new node will be a ghost
158     int chunkNo = chunks[0];
159     addNodeMsg *am = new (numAdjacentNodes,numChunks,0) addNodeMsg;
160     am->chk = index;
161     am->nBetween = numAdjacentNodes;
162     for(int i=0; i<numAdjacentNodes; i++) {
163       am->between[i] = m->getfmMM()->getfmUtil()->exists_in_IDXL(m,adjacentNodes[i],chunkNo,0);
164       CkAssert(am->between[i]!=-1);
165     }
166     am->chunks[0] = chunkNo; //there should be only 1 chunk
167     am->numChunks = numChunks;
168     intMsg *imsg = new intMsg(-1);
169     m->getfmMM()->getfmUtil()->idxllock(m,chunkNo,0);
170     imsg = meshMod[chunkNo].addNodeRemote(am);
171     int newghost = FEM_add_node_local(m,1);
172     m->node.ghost->ghostRecv.addNode(newghost,chunkNo);
173     m->getfmMM()->getfmUtil()->idxlunlock(m,chunkNo,0);
174     //this is the ghostsend index on that chunk, translate it from ghostrecv idxl table
175     //return FEM_To_ghost_index(m->getfmMM()->getfmUtil()->lookup_in_IDXL(m,imsg->i,chunkNo,2));
176     //rather this is same as
177     delete imsg;
178     return FEM_To_ghost_index(newghost);
179   }
180   int newNode = FEM_add_node_local(m, 0);
181   int sharedCount = 0;
182   //interpolate the node data on the new node, based on the ones from which it is created
183   FEM_Interpolate *inp = m->getfmMM()->getfmInp();
184   FEM_Interpolate::NodalArgs nm;
185   nm.n = newNode;
186   for(int i=0; i<numAdjacentNodes; i++) {
187     nm.nodes[i] = adjacentNodes[i];
188   }
189   nm.frac = 0.5;
190   nm.addNode = true;
191   inp->FEM_InterpolateNodeOnEdge(nm);
192   //if any of the adjacent nodes are shared...
193   if(numChunks>1) {
194     // for each adjacent node, if the node is shared
195     for(int i=0;i<numAdjacentNodes;i++){ //a newly added node is shared only if all the 
196       //nodes between which it is added are shared
197       if(is_shared(m, adjacentNodes[i]))
198         {
199           sharedCount++;
200           // lookup adjacent_nodes[i] in IDXL, to find all remote chunks which share this node
201           // call_shared_node_remote() on all chunks for which the shared node exists
202           // we must make sure that we only call the remote entry method once for each remote chunk
203         }
204     }
205     //this is the entry in the IDXL.
206     //just a sanity check, not required really
207     if((sharedCount==numAdjacentNodes && numAdjacentNodes!=0) && forceShared!=-1 ) {
208       //m->getfmMM()->getfmUtil()->splitEntityAll(m, newNode, numAdjacentNodes, adjacentNodes, 0);
209       m->getfmMM()->getfmUtil()->splitEntitySharing(m, newNode, numAdjacentNodes, adjacentNodes, numChunks, chunks);
210     }
211   }
212   return newNode;
215 /** The function called by the entry method on the remote chunk
216     Adds a node on this chunk between the nodes specified in 'between'
217  */
218 void FEM_add_shared_node_remote(FEM_Mesh *m, int chk, int nBetween, int *between){
219   // create local node
220   int newnode = FEM_add_node_local(m, 0);
221   // must negotiate the common IDXL number for the new node, 
222   // and store it in appropriate IDXL tables
223   //note that these are the shared indices, local indices need to be calculated
224   m->getfmMM()->getfmUtil()->splitEntityRemote(m, chk, newnode, nBetween, between);
229 /*Remove node along with helpers -- remove_local/ remove_remote, etc */
230 /** This function invalidates a node entry and if there are any 
231     entries in the IDXL lists (ghost/local) then, these entries are cleared
232     both on the local as well as the remote chunk
233     and then the node is invalidated and any lock it might be holding is reset
234     The pre-condition which is asserted here that the node should not
235     have any adjacency entries, n2e or n2n.
237 void FEM_remove_node_local(FEM_Mesh *m, int node) {
238   // if node is local:
239   int numAdjNodes, numAdjElts;
240   int *adjNodes, *adjElts;
241   m->n2n_getAll(node, adjNodes, numAdjNodes);
242   m->n2e_getAll(node, adjElts, numAdjElts);
243   // mark node as deleted/invalid
244   if(FEM_Is_ghost_index(node)){
245     //CkAssert((numAdjNodes==0) && (numAdjElts==0));
246     //otherwise this ghost node is connected to some element in another chunk, 
247     //which the chunk that just informed us doesn't know abt
248     //look up the ghostrecv idxl list & clean up all instances
249     const IDXL_Rec *irec = m->node.ghost->ghostRecv.getRec(FEM_To_ghost_index(node));
250     //this ghost node might be sent from a chunk without any elems
251     int size = 0;
252     if(irec) size = irec->getShared();
253     int *chunks1, *inds1;
254     if(size>0) {
255       chunks1 = new int[size];
256       inds1 = new int[size];
257     }
258     for(int i=0; i<size; i++) {
259       chunks1[i] = irec->getChk(i);
260       inds1[i] = irec->getIdx(i);
261     }
262     int index = m->getfmMM()->getIdx();
263     for(int i=0; i<size; i++) {
264       int chk = chunks1[i];
265       int sharedIdx = inds1[i];
266       m->node.ghost->ghostRecv.removeNode(FEM_To_ghost_index(node),chk);
267       meshMod[chk].removeIDXLRemote(index,sharedIdx,1);
268     }
269     if(size>0) {
270       delete [] chunks1;
271       delete [] inds1;
272     }
273     m->node.ghost->set_invalid(FEM_To_ghost_index(node),true);
274   }
275   else {
276     //look it up on the idxl list and delete any instances in it
277     const IDXL_Rec *irec = m->node.ghostSend.getRec(node);
278     int size = 0;
279     if(irec) size = irec->getShared();
280     if(size > 0) {
281       int *chknos = (int*)malloc(size*sizeof(int));
282       int *sharedIndices = (int*)malloc(size*sizeof(int));
283       for(int i=0; i<size; i++) {
284         chknos[i] = irec->getChk(i);
285         sharedIndices[i] = irec->getIdx(i);
286       }
287       for(int chkno=0; chkno<size; chkno++) {
288         int remoteChunk = chknos[chkno];
289         int sharedIdx = sharedIndices[chkno];
290         //go to that remote chunk & delete this idxl list and this ghost node
291         meshMod[remoteChunk].removeGhostNode(m->getfmMM()->getIdx(), sharedIdx);
292         m->node.ghostSend.removeNode(node, remoteChunk);
293       }
294       free(chknos);
295       free(sharedIndices);
296     }
297     /*if(!((numAdjNodes==0) && (numAdjElts==0))) {
298       CkPrintf("Error: Node %d cannot be removed, it is connected to :\n",node);
299       m->getfmMM()->fmUtil->FEM_Print_n2e(m,node);
300       m->getfmMM()->fmUtil->FEM_Print_n2n(m,node);
301       }*/
302     // we shouldn't be removing a node away that is connected to something
303     CkAssert((numAdjNodes==0) && (numAdjElts==0));
304     //should be done for the locked node only, i.e. the node on the smallest chunk no
305     m->node.set_invalid(node,true);
306     m->getfmMM()->fmLockN[node].reset(node,m->getfmMM());
307   }
308   if(numAdjNodes != 0) delete[] adjNodes;
309   if(numAdjElts != 0) delete[] adjElts;
310 #ifdef DEBUG_2
311   CkPrintf("[%d]Removing Node %d\n",m->getfmMM()->getfmUtil()->getIdx(),node);
312 #endif
313   return;
316 /** remove a local or shared node, but NOT a ghost node
317     Should probably be able to handle ghosts someday, but here
318     we interpret it as deleting the entry in the ghost list.
320 void FEM_remove_node(FEM_Mesh *m, int node){
321   CkAssert(node >= 0);
322   //someone might actually want to remove a ghost node... 
323   //when there is a ghost edge with both end nodes shared
324   //we will have to notice such a situation and call it on the remotechunk
325   if(FEM_Is_ghost_index(node)) {
326     //this is interpreted as a chunk trying to delete this ghost node, 
327     //this is just to get rid of its version of the node entity
328     CkAssert(m->node.ghost->is_valid(FEM_To_ghost_index(node))); // make sure the node is still there
329     FEM_remove_node_local(m,node);
330   }
331   else {
332     CkAssert(m->node.is_valid(node)); // make sure the node is still there
333     // if node is shared:
334     if(is_shared(m, node)){
335       // verify it is not adjacent to any elements locally
336       int numAdjNodes, numAdjElts;
337       int *adjNodes, *adjElts;
338       m->n2n_getAll(node, adjNodes, numAdjNodes);
339       m->n2e_getAll(node, adjElts, numAdjElts);
340       // we shouldn't be removing a node away that is connected to anything
341       CkAssert((numAdjNodes==0) && (numAdjElts==0)); 
342       // verify it is not adjacent to any elements on any of the associated chunks
343       // delete it on remote chunks(shared and ghost), update IDXL tables
344       m->getfmMM()->getfmUtil()->removeNodeAll(m, node);
345       // mark node as deleted/invalid locally
346       //FEM_remove_node_local(m,node);
347       //m->node.set_invalid(node,true);
348       if(numAdjNodes!=0) delete[] adjNodes;
349       if(numAdjElts!=0) delete[] adjElts;
350     }
351     else {
352       FEM_remove_node_local(m,node);
353     }
354   }
359 /*Add an element, along with helpers -- add_remote, add_local, update_e2e etc */
360 /** A helper function for FEM_add_element_local below
361     Will only work with the same element type as the one given, may crash otherwise
363 void update_new_element_e2e(FEM_Mesh *m, int newEl, int elemType);
364 void FEM_update_new_element_e2e(FEM_Mesh *m, int newEl, int elemType){
365     update_new_element_e2e(m, newEl, elemType);
368 /** Regenerate the e2e connectivity for the given element.
369  *   Currently supports multiple element types.
370  */
371 void update_new_element_e2e(FEM_Mesh *m, int newEl, int elemType) {
373         // Create tuple table
374         FEM_ElemAdj_Layer *g = m->getElemAdjLayer();
375         CkAssert(g->initialized);
376         tupleTable table(g->nodesPerTuple);
378         FEM_Symmetries_t allSym;
379         // locate all elements adjacent to the nodes adjacent to the new 
380         // element, including ghosts, and the new element itself
381         const int nodesPerElem = m->elem[elemType].getNodesPer();
382         int *adjnodes = new int[nodesPerElem];
384         /// A list of all elements adjacent to one of the nodes
385         CkVec<ElemID> nodeAdjacentElements;
387         m->e2n_getAll(newEl, adjnodes, elemType);
388         for(int i=0;i<nodesPerElem;i++){
389                 CkVec<ElemID> adjElements = m->n2e_getAll(adjnodes[i]);
390                 for(int j=0;j<adjElements.size();j++){
391                         int found=0;
392                         // only insert if it is not already in the list
393                         // we use a slow linear scan of the CkVec
394                         for(int i=0;i<nodeAdjacentElements.length();i++) {
395                                 if(nodeAdjacentElements[i] == adjElements[j]) found=1;
396                         }
397                         if(!found){
398                                 nodeAdjacentElements.push_back(adjElements[j]);
399                         }
400                 }
401         }
403         delete[] adjnodes;
404         // Add all the potentially adjacent elements to the tuple table
405         for(int i=0;i<nodeAdjacentElements.length();i++){
406                 int nextElem = nodeAdjacentElements[i].getSignedId();
407                 int nextElemType = nodeAdjacentElements[i].getUnsignedType();
409                 int tuple[tupleTable::MAX_TUPLE];
410                 int *conn;
412                 if(FEM_Is_ghost_index(nextElem))
413                         conn=((FEM_Elem*)m->elem[nextElemType].getGhost())->connFor(FEM_To_ghost_index(nextElem));
414                 else
415                         conn=m->elem[nextElemType].connFor(nextElem);
417                 for (int u=0;u< g->elem[nextElemType].tuplesPerElem;u++) {
418                         for (int i=0;i< g->nodesPerTuple; i++) {
419                                 int eidx=g->elem[nextElemType].elem2tuple[i+u*g->nodesPerTuple];
420                                 if (eidx==-1)  //"not-there" node--
421                                         tuple[i]=-1; //Don't map via connectivity
422                                 else           //Ordinary node
423                                         tuple[i]=conn[eidx]; 
424                         }
425                         table.addTuple(tuple,new elemList(0,nextElem,nextElemType,allSym,u)); 
426                 }
427         }
429         // extract true adjacencies from table and update all e2e tables for both newEl and the others
430         table.beginLookup();
431         // look through each elemList that is returned by the tuple table
432         elemList *l;
433         FEM_IndexAttribute *elemAdjTypesAttr = (FEM_IndexAttribute *)m->elem[elemType].lookup(FEM_ELEM_ELEM_ADJ_TYPES,"update_new_element_e2e");
434         FEM_IndexAttribute *elemAdjAttr = (FEM_IndexAttribute *)m->elem[elemType].lookup(FEM_ELEM_ELEM_ADJACENCY,"update_new_element_e2e");
435         FEM_IndexAttribute *elemAdjTypesAttrGhost = (FEM_IndexAttribute *)m->elem[elemType].getGhost()->lookup(FEM_ELEM_ELEM_ADJ_TYPES,"update_new_element_e2e");
436         FEM_IndexAttribute *elemAdjAttrGhost = (FEM_IndexAttribute *)m->elem[elemType].getGhost()->lookup(FEM_ELEM_ELEM_ADJACENCY,"update_new_element_e2e");
437         //allocate some data structures
438         AllocTable2d<int> &adjTable = elemAdjAttr->get();
439         int *adjs = adjTable.getData();
440         AllocTable2d<int> &adjTypesTable = elemAdjTypesAttr->get();
441         int *adjTypes = adjTypesTable.getData();
442         AllocTable2d<int> &adjTableGhost = elemAdjAttrGhost->get();
443         int *adjsGhost = adjTableGhost.getData();
444         AllocTable2d<int> &adjTypesTableGhost = elemAdjTypesAttrGhost->get();
445         int *adjTypesGhost = adjTypesTableGhost.getData();
446         while (NULL!=(l=table.lookupNext())) {
447                 if (l->next==NULL) {} // One-entry list: must be a symmetry
448                 else { /* Several elements in list: normal case */
449                         // for each a,b pair of adjacent edges
450                         for (const elemList *a=l;a!=NULL;a=a->next){
451                                 for (const elemList *b=a->next;b!=NULL;b=b->next){
452                                         // if a and b are different elements
453                                         if((a->localNo != b->localNo) || (a->type != b->type)){
454                                                 //CkPrintf("%d:%d:%d adj to %d:%d:%d\n", a->type, a->localNo, a->tupleNo, b->type, b->localNo, b->tupleNo);
455                                                 // Put b in a's adjacency list
456                                                 if(FEM_Is_ghost_index(a->localNo)){
457                                                         const int j = FEM_To_ghost_index(a->localNo)* g->elem[a->type].tuplesPerElem + a->tupleNo;
458                                                         adjsGhost[j] = b->localNo;
459                                                         adjTypesGhost[j] = b->type;
460                                                 }
461                                                 else{
462                                                         const int j= a->localNo*g->elem[a->type].tuplesPerElem + a->tupleNo;
463                                                         adjs[j] = b->localNo;
464                                                         adjTypes[j] = b->type;
465                                                 }
466                                                 // Put a in b's adjacency list
467                                                 if(FEM_Is_ghost_index(b->localNo)){
468                                                         const int j = FEM_To_ghost_index(b->localNo)*g->elem[b->type].tuplesPerElem + b->tupleNo;
469                                                         adjsGhost[j] = a->localNo;
470                                                         adjTypesGhost[j] = a->type;
471                                                 }
472                                                 else{
473                                                         const int j= b->localNo*g->elem[b->type].tuplesPerElem + b->tupleNo;
474                                                         adjs[j] = a->localNo;
475                                                         adjTypes[j] = a->type;
476                                                 }
477                                         }
478                                 }
479                         }
480                 }
481         }
484 /** A helper function that adds the local element, and updates adjacencies
485  */
486 int FEM_add_element_local(FEM_Mesh *m, int *conn, int connSize, int elemType, bool addGhost, bool create_adjacencies){
487   // lengthen element attributes
488   int newEl;
489   if(addGhost){
490     // find a place in the array for the new el
491     newEl = m->elem[elemType].getGhost()->get_next_invalid(m,false,true); 
492     // update element's conn, i.e. e2n table
493     ((FEM_Elem*)m->elem[elemType].getGhost())->connIs(newEl,conn);
494     // get the signed ghost value
495     newEl = FEM_From_ghost_index(newEl);
496   }
497   else{
498     newEl = m->elem[elemType].get_next_invalid(m,false,false);
499     /*#ifdef FEM_ELEMSORDERED
500     //correct the orientation, which might have been lost if it moves across chunks
501     if(m->getfmMM()->fmAdaptAlgs->getSignedArea(conn[0],conn[1],conn[2])>0) {
502       int tmp = conn[1];
503       conn[1] = conn[2];
504       conn[2] = tmp;
505     }
506     #endif*/
507     // update element's conn, i.e. e2n table
508     m->elem[elemType].connIs(newEl,conn); 
509   }
511   if(create_adjacencies){
512         // add to n2e and n2n table
513         for(int i=0;i<connSize;i++){
514           m->n2e_add(conn[i],newEl);
515           for(int j=i+1;j<connSize;j++){
516                 if(! m->n2n_exists(conn[i],conn[j]))
517                   m->n2n_add(conn[i],conn[j]);
518                 if(! m->n2n_exists(conn[j],conn[i]))
519                   m->n2n_add(conn[j],conn[i]);
520           }
521         }
522         // update e2e table -- too complicated, so it gets is own function
523         m->e2e_removeAll(newEl);
524         update_new_element_e2e(m,newEl,elemType);
525         int *adjes = new int[connSize];
526         m->e2e_getAll(newEl, adjes, 0);
527         
528 //      CkAssert( (m->elem[elemType].size()<2) || !((adjes[0]==adjes[1] && adjes[0]!=-1) || (adjes[1]==adjes[2] && adjes[1]!=-1) || (adjes[2]==adjes[0] && adjes[2]!=-1)));
530 #ifdef DEBUG_2
531         CkPrintf("[%d]Adding element %d(%d,%d,%d)\n",m->getfmMM()->getfmUtil()->getIdx(),newEl,conn[0],conn[1],conn[2]);
532 #endif
533 #ifdef FEM_ELEMSORDERED
534         if(!FEM_Is_ghost_index(newEl)) {
535           m->getfmMM()->fmAdaptAlgs->ensureQuality(conn[0],conn[1],conn[2]);
536         }
537 #endif
538         delete[] adjes;
539   }
540   
541   return newEl;
544 /** This addds an element with the given 'conn' connecitivity and 'connSize' number of nodes in e2n
545     to the mesh 'm'. 'chunkNo' denotes the chunk the new element should be added to.
546     If it is -1, then this function decides which chunk this new element should belong to,
547     based on its node connectivity.
548     It returns the index of the newly created element
550 int FEM_add_element(FEM_Mesh *m, int* conn, int connSize, int elemType, int chunkNo){
551   //the first step is to figure out the type of each node (local/shared/ghost)
552   int newEl = -1;
553   int index = m->getfmMM()->getIdx();
554   int buildGhosts = 0; //this determines if we need to modify ghosts locally or remotely
555                        //if all nodes are local, there is no need to update ghosts anywhere
556   int sharedcount=0;
557   int ghostcount=0;
558   int localcount=0;
559   CkAssert(conn[0]!=conn[1] &&conn[1]!=conn[2] &&conn[2]!=conn[0]);
560   int *nodetype = (int *)malloc(connSize *sizeof(int)); //0 -- local, 1 -- shared, 2--ghost
562 #ifdef DEBUG_2
563   CkPrintf("addElement, line %d\n", __LINE__);
564 #endif
565   for(int i=0;i<connSize;i++){
566     //make sure that every connectivity is valid
567     CkAssert(conn[i]!=-1); 
568     if(is_shared(m,conn[i])) {
569       nodetype[i] = 1;
570       sharedcount++;
571     }
572     else if(FEM_Is_ghost_index(conn[i])) {
573       nodetype[i] = 2;
574       ghostcount++;
575     }
576     else {
577       nodetype[i] = 0;
578     }
579   }
580   //if it is not shared or ghost, it should be local
581   localcount = connSize - (sharedcount + ghostcount);
583 #ifdef DEBUG_2
584   CkPrintf("addElement, line %d\n", __LINE__);
585 #endif
586   if(sharedcount==0 && ghostcount==0){
587     // add a local elem with all local nodes
588     newEl = FEM_add_element_local(m,conn,connSize,elemType,0);
589     //no modifications required for ghostsend or ghostrecv of nodes or elements
590   }
591   else if(ghostcount==0 && sharedcount > 0 && localcount > 0){
592     // a local elem with not ALL shared nodes or ALL local nodes
593     newEl = FEM_add_element_local(m,conn,connSize,elemType,0);
594     buildGhosts = 1;
595   }
596   else if(ghostcount==0 && sharedcount > 0 && localcount == 0){
597     // a local elem with ALL shared nodes
598     //it is interesting to note that such a situation can only occur between two chunks,
599     //in any number of dimensions.
600     //So, the solution is to figure out, which chunk it belongs to...
601     if(!(chunkNo==index || chunkNo==-1)) {
602       //if chunkNo points to a remotechunk, then this element should be added to that chunk
603       addElemMsg *am = new (connSize, 0, 0) addElemMsg;
604       int chk = index;
605       am->chk = chk;
606       am->elemtype = elemType;
607       am->connSize = connSize;
608       am->numGhostIndex = 0;
609       for(int i=0; i<connSize; i++) {
610         //translate the node connectivity local to this chunk to shared idxl entries
611         CkAssert(nodetype[i] == 1);
612         am->conn[i] = m->getfmMM()->getfmUtil()->exists_in_IDXL(m,conn[i],chunkNo,0);
613         CkAssert(am->conn[i] >= 0);
614       }
615       intMsg* imsg = meshMod[chunkNo].addElementRemote(am);
616       int shidx = imsg->i;
617       delete imsg;
618       const IDXL_List ilist = m->elem[elemType].ghost->ghostRecv.addList(chunkNo);
619       newEl = ilist[shidx];
620       newEl = FEM_To_ghost_index(newEl);
621       free(nodetype);
622       return newEl;
623     }
624     //otherwise, it will be added on this chunk
625     newEl = FEM_add_element_local(m,conn,connSize,elemType,0);
626     buildGhosts = 1;
627   }
628   else if(ghostcount > 0 /*&& localcount == 0 && sharedcount > 0*/) { 
629     //I will assume that this is a correct call.. the caller never gives junk ghost elements
630     //it is a remote elem with some shared nodes
631     //this is the part when it eats a ghost element
632 #ifdef DEBUG_2
633   CkPrintf("addElement, line %d\n", __LINE__);
634 #endif
635     if((chunkNo!=-1) && (chunkNo==index)) { 
636       //this is because the chunk doing this operation is supposed to eat into someone else.
637       //change all ghost nodes to shared nodes
638       //this also involves going to allchunks that it was local/shared 
639       //on and make it shared and add this chunk to the list of chunks it is shared to.
640 #ifdef DEBUG_2
641   CkPrintf("addElement, line %d\n", __LINE__);
642 #endif
643       for(int i=0; i<connSize; i++) {
644         if(nodetype[i]==2) {
645           //build up the list of chunks it is shared/local to
646           int numchunks;
647           IDXL_Share **chunks1;
648           //make sure that the ghostRecv table is completely populated...
649           //if this is a ghost->local transition, then it will be local on only one
650           //chunk, i.e. on the chunk where it was local before this eat operation started.
651 #ifdef DEBUG_2
652   CkPrintf("addElement, line %d\n", __LINE__);
653 #endif
654           m->getfmMM()->getfmUtil()->getChunkNos(0,conn[i],&numchunks,&chunks1);
655           //add a new node with the same attributes as this ghost node,
656           //do not remove the ghost node yet
657 #ifdef DEBUG_2
658   CkPrintf("addElement, line %d\n", __LINE__);
659 #endif
660           int newN = m->getfmMM()->getfmUtil()->Replace_node_local(m, conn[i], -1);
661           //need to know if the newly added node will be local or shared
662           //if shared, add index to the shared list of this node on all the chunks
663           //if local, remove the oldIdx on the other chunk, if the oldIdx
664           //will be local only on this chunk after this operation
665           //(happens during a local->ghost transition)
666 #ifdef DEBUG_2
667   CkPrintf("addElement, line %d\n", __LINE__);
668 #endif
669           for(int j=0; j<numchunks; j++) {
670             int chk = chunks1[j]->chk;
671             if(chk==index) continue;
672             //coarser lock
673 #ifdef DEBUG_2
674   CkPrintf("addElement, line %d\n", __LINE__);
675 #endif
676             m->getfmMM()->getfmUtil()->idxllock(m,chk,0);
677 #ifdef DEBUG_2
678   CkPrintf("addElement, line %d\n", __LINE__);
679 #endif
680             //this new node, might be a shared node or a local node... 
681             m->node.shared.addNode(newN,chk);
682 #ifdef DEBUG_2
683   CkPrintf("addElement, line %d\n", __LINE__);
684 #endif
685             //find out what chk calls this node (from the ghostrecv idxl list)
686             int idx = m->getfmMM()->getfmUtil()->exists_in_IDXL(m, conn[i], chk, 2);
687             m->node.ghost->ghostRecv.removeNode(FEM_To_ghost_index(conn[i]), chk);
688 #ifdef DEBUG_2
689   CkPrintf("addElement, line %d\n", __LINE__);
690 #endif
691             //when doing this add, make all idxl lists consistent, and return the list 
692             //of n2e connections for this node, which are local to this chunk, 
693             //and also a list of the connectivity
694             meshMod[chk].addToSharedList(index, idx); 
695 #ifdef DEBUG_2
696   CkPrintf("addElement, line %d\n", __LINE__);
697 #endif
698             //when adding a new node, always verify if this new node is a corner on
699             //the other chunk, if so add it to the list of fixed nodes (corners)
700             int idx1 = m->getfmMM()->getfmUtil()->exists_in_IDXL(m, newN, chk, 0);
701 #ifdef DEBUG_2
702   CkPrintf("addElement, line %d\n", __LINE__);
703 #endif
704             boolMsg* bmsg1 = meshMod[chk].isFixedNodeRemote(index,idx1);
705             bool isfixed = bmsg1->b;
706             delete bmsg1;
707             if(isfixed) m->getfmMM()->fmfixedNodes.push_back(newN);
708 #ifdef DEBUG_2
709   CkPrintf("addElement, line %d\n", __LINE__);
710 #endif
711             //coarser lock
712             m->getfmMM()->getfmUtil()->idxlunlock(m, chk, 0);
713 #ifdef DEBUG_2
714   CkPrintf("addElement, line %d\n", __LINE__);
715 #endif
716           }
717           nodetype[i] = 1;
718           sharedcount++;
719           ghostcount--;
720           //remove the ghost node
721           FEM_remove_node_local(m,conn[i]);
722           conn[i] = newN;
723           //lock the newly formed node, if it needs to be
724 #ifdef DEBUG_2
725   CkPrintf("addElement, line %d\n", __LINE__);
726 #endif
727 #ifndef CPSD
728           FEM_Modify_LockUpdate(m,conn[i]);
729 #endif
730 #ifdef DEBUG_2
731   CkPrintf("addElement, line %d\n", __LINE__);
732 #endif
733           for(int j=0; j<numchunks; j++) {
734             delete chunks1[j];
735           }
736           if(numchunks != 0) free(chunks1);
737         }
738       }
739 #ifdef DEBUG_2
740   CkPrintf("addElement, line %d\n", __LINE__);
741 #endif
742       newEl = FEM_add_element_local(m,conn,connSize,elemType,0);
743 #ifdef DEBUG_2
744   CkPrintf("addElement, line %d\n", __LINE__);
745 #endif
746       buildGhosts = 1;
747     }
748     else {
749 #ifdef DEBUG_2
750   CkPrintf("addElement, line %d\n", __LINE__);
751 #endif
752       int numSharedChunks = 0;
753       int remoteChunk = -1;
754       if(chunkNo==-1) {
755         CkAssert(false); //shouldn't be here
756         //figure out which chunk it is local to
757         //among all chunks who share some of the nodes or from whom this chunk receives ghost nodes
758         //if there is any chunk which owns a ghost node which is not shared, then that is a local node 
759         //to that chunk and that chunk owns that element. However, only that chunk knows abt it.
760         //So, just go to the owner of every ghost node and figure out who all share that node.
761         //Build up this table of nodes owned by which all chunks.
762         //The chunk that is in the table corresponding to all nodes wins the element.
763         CkVec<int> **allShared;
764         CkVec<int> *allChunks = NULL;
765         int **sharedConn; 
766         m->getfmMM()->getfmUtil()->buildChunkToNodeTable(nodetype, sharedcount, ghostcount, localcount, conn, connSize, &allShared, &numSharedChunks, &allChunks, &sharedConn);   
767         //we are looking for a chunk which does not have a ghost node
768         for(int i=0; i<numSharedChunks; i++) {
769           remoteChunk = i;
770           for(int j=0; j<connSize; j++) {
771             if(sharedConn[i][j] == -1 || sharedConn[i][j] == 2) {
772               remoteChunk = -1;
773               break; //this chunk has a ghost node
774             }
775             //if(sharedConn[i][j] == 0) {
776             //  break; //this is a local node, hence it is the remotechunk
777             //}
778           }
779           if(remoteChunk == i) break;
780           else remoteChunk = -1;
781         }
782         if(remoteChunk==-1) {
783           //every chunk has a ghost node.
784           if(chunkNo != -1) {
785             //upgrade the ghost node to a shared node on chunkNo
786             remoteChunk = chunkNo;
787             for(int k=0; k<numSharedChunks; k++) {
788               if(chunkNo == (*allChunks)[k]) {
789                 for(int l=0; l<connSize; l++) {
790                   if(sharedConn[k][l]==-1 || sharedConn[k][l]==2) {
791                     //FIXME: upgrade this ghost node
792                   }
793                 }
794               }
795             }
796           }
797           else {
798             CkPrintf("[%d]ERROR: Can not derive where (%d,%d,%d) belongs to\n",index,conn[0],conn[1],conn[2]);
799             CkAssert(false);
800           }
801         }
802         remoteChunk = (*allChunks)[remoteChunk];
803         //convert all connections to the shared IDXL indices. We should also tell which are ghost indices
804         CkPrintf("[%d]Error: I derived it should go to chunk %d, which is not %d\n",index,remoteChunk,chunkNo);
805         for(int k=0; k<numSharedChunks; k++) {
806           free(sharedConn[k]);
807         }
808         if(numSharedChunks!=0) free(sharedConn);
809         for(int k=0; k<connSize; k++) {
810           delete allShared[k];
811         }
812         free(allShared);
813         allChunks->free();
814         if(allChunks!=NULL) free(allChunks);
815       }
816       else {
817         remoteChunk=chunkNo;
818       }
819       int numGhostNodes = 0;
820       for(int i=0; i<connSize; i++) {
821         if(nodetype[i] == 2) { //a ghost node
822           numGhostNodes++;
823         }
824       }
825       CkAssert(numGhostNodes > 0);
826       addElemMsg *am = new (connSize, numGhostNodes, 0) addElemMsg;
827       int chk = index;
828       am->chk = chk;
829       am->elemtype = elemType;
830       am->connSize = connSize;
831       am->numGhostIndex = numGhostNodes;
832       int j = 0;
833       for(int i=0; i<connSize; i++) {
834         if(nodetype[i] == 1) {
835           am->conn[i] = m->getfmMM()->getfmUtil()->exists_in_IDXL(m,conn[i],remoteChunk,0);
836           CkAssert(am->conn[i] >= 0);
837         }
838         else if(nodetype[i] == 2) {
839           am->conn[i] = m->getfmMM()->getfmUtil()->exists_in_IDXL(m,conn[i],remoteChunk,2);
840           CkAssert(am->conn[i] >= 0);
841           am->ghostIndices[j] = i;
842           j++;
843         }
844       }
845       intMsg* imsg = meshMod[remoteChunk].addElementRemote(am);
846       int shidx = imsg->i;
847       delete imsg;
848       const IDXL_List ilist = m->elem[elemType].ghost->ghostRecv.addList(remoteChunk);
849       newEl = ilist[shidx];
850       newEl = FEM_To_ghost_index(newEl);
851     }
852 #ifdef DEBUG_2
853   CkPrintf("addElement, line %d\n", __LINE__);
854 #endif
855   }
856   else if(ghostcount > 0 && localcount == 0 && sharedcount == 0) { 
857     // it is a remote elem with no shared nodes
858     //I guess such a situation will never occur
859     //figure out which chunk it is local to -- this would be really difficult to do in such a case.
860     //so, for now, we do not allow a chunk to add an element 
861     //for which it does not share even a single node.
862     // almost the code in the preceeding else case will work for this, but its better not to allow this
863   }
864   else if(ghostcount > 0 && localcount > 0 && sharedcount > 0){
865     // it is a flip operation
866     //actually this can be generalized as an operation which moves the boundary,
867     //to make the ghost node shared
868     //if one uses FEM_elem_acquire, then this condition gets distributed 
869     //across other conditions already done.
870     //   promote ghosts to shared on others, requesting new ghosts
871     //   grow local element and attribute tables if needed
872     //   add to local elem[elemType] table, and update IDXL if needed
873     //   update remote adjacencies
874     //   update local adjacencies
875   }
876   else if(ghostcount > 0 && localcount > 0 && sharedcount == 0) { 
877     //this is an impossible case
878   }
879   //start building ghosts if required
880   if(buildGhosts==1) {
881     //   make this element ghost on all others, updating all IDXL's
882     //   also in same remote entry method, update adjacencies on all others
883     //   grow local element and attribute tables if needed
884     //   add to local elem[elemType] table, and update IDXL if needed
885     //   update local adjacencies
886     //   return the new element id
887     //build a mapping of all shared chunks to all nodes in this element    
888 #ifdef DEBUG_2
889   CkPrintf("addElement, line %d\n", __LINE__);
890 #endif
891     CkVec<int> **allShared;
892     int numSharedChunks = 0;
893     CkVec<int> *allChunks = NULL;
894     int **sharedConn; 
895     m->getfmMM()->getfmUtil()->buildChunkToNodeTable(nodetype, sharedcount, ghostcount, localcount, conn, connSize, &allShared, &numSharedChunks, &allChunks, &sharedConn);   
896     //add all the local nodes in this element to the ghost list, if they did not exist already
897     for(int i=0; i<numSharedChunks; i++) {
898       int chk = (*allChunks)[i];
899       if(chk == index) continue; //it is this chunk
900       //it is a new element so it could not have existed as a ghost on that chunk. Just add it.
901       m->getfmMM()->getfmUtil()->idxllock(m, chk, 0);
902       m->elem[elemType].ghostSend.addNode(newEl,chk);
903       //ghost nodes should be added only if they were not already present as ghosts on that chunk.
904       int *sharedIndices = (int *)malloc((connSize)*sizeof(int));
905       int *typeOfIndex = (int *)malloc((connSize)*sizeof(int));
906       for(int j=0; j<connSize; j++) {
907         int sharedNode = m->getfmMM()->getfmUtil()->exists_in_IDXL(m,conn[j],chk,0);
908         if(sharedNode == -1) {
909           //node 'j' is a ghost on chunk 'i'
910           int sharedGhost = m->getfmMM()->getfmUtil()->exists_in_IDXL(m,conn[j],chk,1);
911           if( sharedGhost == -1) {
912             //this might be a new ghost, figure out if any of the chunks sharing
913             //this node has created this as a ghost on 'chk'
914             const IDXL_Rec *irec = m->node.shared.getRec(conn[j]);
915             if(irec) {
916               int noShared = irec->getShared();
917               for(int sharedck=0; sharedck<noShared; sharedck++) {
918                 int ckshared = irec->getChk(sharedck);
919                 int idxshared = irec->getIdx(sharedck);
920                 if(ckshared == chk) continue;
921                 CkAssert(chk!=index && chk!=ckshared && ckshared!=index);
922                 intMsg* imsg = meshMod[ckshared].getIdxGhostSend(index,idxshared,chk);
923                 int idxghostsend = imsg->i;
924                 delete imsg;
925                 if(idxghostsend != -1) {
926                   m->node.ghostSend.addNode(conn[j],chk);
927                   meshMod[chk].updateIdxlList(index,idxghostsend,ckshared);
928                   sharedGhost = m->getfmMM()->getfmUtil()->exists_in_IDXL(m,conn[j],chk,1);
929                   CkAssert(sharedGhost != -1);
930                   break; //found a chunk that sends it out, update my tables
931                 }
932                 //Chunk 'ckshared' does not send this to Chunk 'chk' as ghost
933               }
934             }
935             //else it is a new ghost
936           }
937           if(sharedGhost == -1) {
938             sharedIndices[j] = conn[j];
939             typeOfIndex[j] = -1;
940           }
941           else {
942             //it is a shared ghost
943             sharedIndices[j] = sharedGhost;
944             typeOfIndex[j] = 1;
945           }
946         }
947         else {
948           //it is a shared node
949           sharedIndices[j] = sharedNode;
950           typeOfIndex[j] = 0;
951         }
952       }
953       addGhostElemMsg *fm = new (connSize, connSize, 0)addGhostElemMsg;
954       fm->chk = index;
955       fm->elemType = elemType;
956 #ifdef DEBUG_2
957   CkPrintf("addElement, line %d\n", __LINE__);
958 #endif
959       for(int j=0; j<connSize; j++) {
960         fm->indices[j] = sharedIndices[j];
961         fm->typeOfIndex[j] = typeOfIndex[j];
962         if(typeOfIndex[j]==-1) {
963           m->node.ghostSend.addNode(sharedIndices[j],chk);
964         }
965       }
966       fm->connSize = connSize;
967       meshMod[chk].addGhostElem(fm); //newEl, m->fmMM->idx, elemType;
968       m->getfmMM()->getfmUtil()->idxlunlock(m, chk, 0);
969       free(sharedIndices);
970       free(typeOfIndex);
971     }
972     for(int k=0; k<numSharedChunks; k++) {
973       free(sharedConn[k]);
974     }
975     if(numSharedChunks!=0) free(sharedConn);
976     for(int k=0; k<connSize; k++) {
977       delete allShared[k];
978     }
979     free(allShared);
980     allChunks->free();
981     if(allChunks!=NULL) free(allChunks);
982 #ifdef DEBUG_2
983   CkPrintf("addElement, line %d\n", __LINE__);
984 #endif
985   }
986   free(nodetype);
987 #ifdef DEBUG_2
988   CkPrintf("addElement, line %d\n", __LINE__);
989 #endif
990   return newEl;
995 /*Remove an element with helpers --remove_local, etc */
996 /** Remove a local element from any of the adjacency tables that use it,
997     namely the n2e, e2e and e2n tables of the neighbors. 
998     e2n of the neighbors does not change
999     all tables for this element get deleted anyway when it is marked invalid, which is 
1000     now done in the purge step
1002 void FEM_remove_element_local(FEM_Mesh *m, int element, int etype){
1003   // replace this element with -1 in adjacent nodes' adjacencies
1004   // nodesPerEl be the number of nodes that can be adjacent to this element
1005   const int nodesPerEl = m->elem[etype].getConn().width(); 
1006   int *adjnodes = new int[nodesPerEl];
1007   m->e2n_getAll(element, adjnodes, etype);
1008 #ifdef DEBUG_2
1009   CkPrintf("[%d]Deleting element %d(%d,%d,%d)\n",m->getfmMM()->getfmUtil()->getIdx(),element,adjnodes[0],adjnodes[1],adjnodes[2]);
1010 #endif
1011   for(int i=0;i<nodesPerEl;i++) {
1012     //if(adjnodes[i] != -1) //if an element is local, then an adjacent node should not be -1
1013     m->n2e_remove(adjnodes[i],element);
1014   }
1015   // replace this element with -1 in adjacent elements' adjacencies
1016   const int numAdjElts = nodesPerEl;  
1017   int *adjelts = new int[numAdjElts]; 
1018   m->e2e_getAll(element, adjelts, etype);
1019   for(int i=0;i<numAdjElts;i++){
1020     m->e2e_replace(adjelts[i],element,-1);
1021   }
1022   // We must now remove any n2n adjacencies which existed because of the 
1023   // element that is now being removed. This is done by removing all 
1024   // n2n adjacencies for the nodes of this element, and recreating those
1025   // that still exist by using the neighboring elements.
1026   //FIXME: I think, if we just consider on which all faces (for a 2D, a face is an edge),
1027   //we have valid neighboring elements, then none of the edges (n2n conn) on that
1028   //face should go off, this would be more efficient.
1029   for(int i=0;i<nodesPerEl;i++) {
1030     for(int j=i+1;j<nodesPerEl;j++){
1031       m->n2n_remove(adjnodes[i],adjnodes[j]);
1032       m->n2n_remove(adjnodes[j],adjnodes[i]);
1033     }
1034   }
1035   for(int i=0;i<numAdjElts;i++){ // for each neighboring element
1036     if(adjelts[i] != -1){
1037       int *adjnodes2 = new int[nodesPerEl];
1038       m->e2n_getAll(adjelts[i], adjnodes2, etype);
1039       // for each j,k pair of nodes adjacent to the neighboring element
1040       for(int j=0;j<nodesPerEl;j++){ 
1041         for(int k=j+1;k<nodesPerEl;k++){   
1042           if(!m->n2n_exists(adjnodes2[j],adjnodes2[k]))
1043             m->n2n_add(adjnodes2[j],adjnodes2[k]);
1044           if(!m->n2n_exists(adjnodes2[k],adjnodes2[j]))
1045             m->n2n_add(adjnodes2[k],adjnodes2[j]);
1046         }
1047       }
1048       delete[] adjnodes2;
1049     }
1050   }
1051   //done in purge now
1052   /*if(FEM_Is_ghost_index(element)){
1053     m->elem[etype].getGhost()->set_invalid(FEM_To_ghost_index(element),false);
1054   }
1055   else {
1056     m->elem[etype].set_invalid(element,false);
1057     }*/
1058   
1059   //cleanup
1060   delete[] adjelts;
1061   delete[] adjnodes;
1062   return;
1065 /** Can be called on local or ghost elements
1066     removes the connectivity of the element on this chunk as well as all chunks on which it is a ghost
1067     'permanent' is -1 if the element is to be deleted, e.g. in a coarsen operation, one element
1068     is deleted permanently.. for a refine however the existing elements are eaten by the new ones
1069     'permanent' is the id of the chunk which wants to be the owner (i.e. wants to eat it)
1070     this is done just to save deleting some things when it would be immediately added again, 
1071     when the element is deleted on one and added on another chunk 
1072     The return value is the chunk on which this element was local
1073     Removing an element is a two step process, FEM_remove_element and FEM_purge_element
1074     This step will get rid of all adjacencies info, but the idxl information
1075     will not be deleted, just ghostsend on this junk and receive info on the remote chunk.
1076     The idxl info will be deleted in Purge, that way shared 
1077     chunks could still refer to it, needed for solution transfer
1079 int FEM_remove_element(FEM_Mesh *m, int elementid, int elemtype, int permanent, bool aggressive_node_removal){
1080   //CkAssert(elementid != -1);
1081 #ifdef DEBUG_2
1082   CkPrintf("removeElement, line %d\n", __LINE__);
1083 #endif
1084   if(elementid == -1) return -1;
1085   int index = m->getfmMM()->getfmUtil()->getIdx();
1086   if(FEM_Is_ghost_index(elementid)){
1087     //an element can come as a ghost from only one chunk, so just convert args and call it on that chunk
1088     int ghostid = FEM_To_ghost_index(elementid);
1089     const IDXL_Rec *irec = m->elem[elemtype].ghost->ghostRecv.getRec(ghostid);
1090     int size = irec->getShared();
1091     CkAssert(size == 1);
1092     int remoteChunk = irec->getChk(0);
1093     int sharedIdx = irec->getIdx(0);
1094     //the message that should go to the remote chunk, where this element lives
1095     removeElemMsg *rm = new removeElemMsg;
1096     rm->chk = index;
1097     rm->elementid = sharedIdx;
1098     rm->elemtype = elemtype;
1099     rm->permanent = permanent;
1100     rm->aggressive_node_removal = aggressive_node_removal;
1101     meshMod[remoteChunk].removeElementRemote(rm);
1102     // remove local ghost element, now done in purge
1103 #ifdef DEBUG_2
1104   CkPrintf("removeElement, line %d\n", __LINE__);
1105 #endif
1106     return remoteChunk;
1107   }
1108   else {
1109     //this is a local element on this chunk
1110     //if it is a ghost element on some other chunks, then for all ghost nodes in this element,
1111     //we should individually go to each ghost and figure out if it has any connectivity in
1112     //the ghost layer after the removal. If not, then it will be removed from their ghost layers.
1113     //This is something similar to an add_element, where we do a similar analaysis on each 
1114     //involved chunk neighbor and add ghost nodes.
1115     //get a list of chunks for which this element is a ghost
1116     //this has to be done as the idxl lists change in this loop, so we
1117     //if we call getRec multiple times, we will get a new irec every time.
1118     const IDXL_Rec *irec = m->elem[elemtype].ghostSend.getRec(elementid);
1119     if(irec){
1120       int numSharedChunks = irec->getShared();
1121       int connSize = m->elem[elemtype].getConn().width();
1122       if(numSharedChunks>0) {
1123         int *chknos = (int*)malloc(numSharedChunks*sizeof(int));
1124         int *inds = (int*)malloc(numSharedChunks*sizeof(int));
1125         for(int i=0; i<numSharedChunks; i++) {
1126           chknos[i] = irec->getChk(i);
1127           inds[i] = irec->getIdx(i);
1128         }
1129         int *newghost, *ghostidx, *losingThisNode, *nodes, *willBeGhost;
1130         //will someone acquire this element soon?
1131         if(permanent>=0) { //this element will be reassigned or change connectivity
1132           newghost = new int[connSize];
1133           ghostidx = new int[connSize]; //need to create new ghosts
1134           losingThisNode = new int[connSize];
1135           nodes = new int[connSize];
1136           willBeGhost = new int[connSize]; 
1137           //willBeGhost is needed because not necessarily a node that is lost
1138           //will be a ghost, think of what happens when an isolated element is deleted
1139           int *elems;
1140           int numElems;
1141           m->e2n_getAll(elementid,nodes,elemtype);
1142           for(int i=0; i<connSize; i++) {
1143             newghost[i] = -1;
1144             ghostidx[i] = -1;
1145             //Is this chunk losing this node?
1146             losingThisNode[i] = 1;
1147             //if this node is connected to any element which is local and is not being deleted now
1148             //then this chunk will not lose this node now
1149             m->n2e_getAll(nodes[i], elems, numElems);
1150             for(int k=0; k<numElems; k++) {
1151               if((!FEM_Is_ghost_index(elems[k])) && (elems[k]!=elementid)) {
1152                 losingThisNode[i] = 0;
1153                 break;
1154               }
1155             }
1156             if(losingThisNode[i]) {
1157               willBeGhost[i] = 1; 
1158               //current allotment, it might be modified, see below
1159             }
1160             free(elems);
1161           }
1162           //verify if it is losing all nodes, this means that this was an isolated element
1163           //and this element does not have to be a ghost
1164           if(losingThisNode[0]==1 && losingThisNode[1]==1 && losingThisNode[2]==1) {
1165             //only if it is connected to any other node of this chunk, it will be a ghost
1166             for(int i=0; i<connSize; i++) {
1167               int *ndnbrs;
1168               int numndnbrs;
1169               m->n2n_getAll(nodes[i], ndnbrs, numndnbrs);
1170               willBeGhost[i]=0;
1171               //if this node is connected to any node that is local other than the nodes of this
1172               //element, only then will this node be a ghost.
1173               for(int n1=0; n1<numndnbrs; n1++) {
1174                 if(ndnbrs[n1]>=0 && ndnbrs[n1]!=nodes[(i+1)%connSize] && ndnbrs[n1]!=nodes[(i+2)%connSize]) {
1175                   willBeGhost[i]=1;
1176                 }
1177               }
1178             }
1179           }
1180           //add a new ghost for every node that should be lost and a ghost added
1181           for(int i=0; i<connSize; i++) {
1182             if(losingThisNode[i]) {
1183               //lock it on the min chunk on which this node is local
1184 #ifndef CPSD
1185               FEM_Modify_LockAll(m,nodes[i],false);
1186 #endif
1187             }
1188             if(losingThisNode[i]==1 && willBeGhost[i]==1) {
1189               newghost[i] = FEM_add_node_local(m,1);
1190               ghostidx[i] = FEM_To_ghost_index(newghost[i]); //translate the index
1191             }
1192           }
1193         }
1194         //the above code just figured out what to do with all the nodes of this element
1195         //and corrected the locks as well
1196         /*//a correction to deal with physical corners
1197         //if losing a node which is a corner, then upgrade that ghostnode on the chunk 
1198         //to a shared node before performing this entire operation
1199         //I can only think of this happening for 1 chunk
1200         for(int i=0; i<connSize; i++) {
1201           if(losingThisNode[i]==1 && m->node.shared.getRec(nodes[i]==NULL)) {
1202           }
1203           }*/
1204         int deleteNodeLater = -1; //if a local node becomes ghost, then one
1205         //has to keep track of that and delete it only after the remote chunk
1206         //has copied the attributes, this interger remembers which node is to be
1207         //deleted later
1208         //Note that the chunk to which it loses is always 'permanent'
1209         for(int i=0; i<numSharedChunks; i++) {
1210           bool lockedRemoteIDXL = false;
1211           int chk = chknos[i];
1212           int sharedIdx = inds[i];
1213           int numGhostNodes = 0;
1214           int *ghostIndices = (int*)malloc(connSize*sizeof(int));
1215           int numGhostRN = 0;
1216           CkVec<int> ghostRNIndices;
1217           int numGhostRE = 0;
1218           CkVec<int> ghostREIndices;
1219           int numSharedNodes = 0;
1220           int *sharedIndices = (int*)malloc(connSize*sizeof(int));
1221           //purge will do this now
1222           //m->elem[elemtype].ghostSend.removeNode(elementid, chk);
1223           if(permanent>=0) {
1224             //get the list of n2e for all nodes of this element. 
1225             //If any node has only this element in its list.
1226             //it no longer should be a ghost on chk
1227             //the next step is to figure out all the ghost elements 
1228             //& nodes that it no longer needs to have
1229             CkVec<int> testelems;
1230             const IDXL_List ll = m->elem[elemtype].ghostSend.addList(chk);
1231             int size = ll.size();
1232             const IDXL_List ln = m->node.ghostSend.addList(chk);
1233             int sizeN = ln.size();
1234             const IDXL_List lre = m->elem[elemtype].ghost->ghostRecv.addList(chk);
1235             int lresize = lre.size();
1236             const IDXL_List lrn = m->node.ghost->ghostRecv.addList(chk);
1237             int lrnsize = lrn.size();
1238             //iterate over the nodes connected to this element
1239             for(int j=0; j<connSize; j++) {
1240               int *elems;
1241               int numElems;
1242               //we are tyring to figure out if this local node, which might be a ghost
1243               //on some chunk be deleted on that chunk
1244               m->n2e_getAll(nodes[j], elems, numElems);
1245               //do not delete ghost nodes on the eating chunk, and be sure 
1246               //that the node is indeed a ghost on the chunk being considered
1247               if((chk != permanent) && (m->getfmMM()->fmUtil->exists_in_IDXL(m,nodes[j],chk,1)!=-1)) {
1248                 //'chk' is not the chunk that is eating this element
1249                 //if any of these elems is a ghost on chk then do not delete this ghost node
1250                 bool shouldBeDeleted = true;
1251                 for(int k=0; k<numElems; k++) {
1252                   if(elems[k]==elementid) continue;
1253                   if(elems[k]>0) {
1254                     if(m->getfmMM()->fmUtil->exists_in_IDXL(m,elems[k],chk,3)!=-1) {
1255                       shouldBeDeleted = false;
1256                       break;
1257                     }
1258                   }
1259                 }
1260                 //find out what the other shared chunks think abt it
1261                 //if every shared chunk believes it should be deleted, get rid of it
1262                 if(shouldBeDeleted) {
1263                   const IDXL_Rec *irecsh = m->node.shared.getRec(nodes[j]);
1264                   if(irecsh!=NULL) {
1265                     for(int k=0; k<irecsh->getShared(); k++) {
1266                       if(shouldBeDeleted) {
1267                         boolMsg* bmsg = meshMod[irecsh->getChk(k)].shouldLoseGhost(index,irecsh->getIdx(k),chk);
1268                         shouldBeDeleted = bmsg->b;
1269                         delete bmsg;
1270                       }
1271                     }
1272                   }
1273                 }
1274                 //since we have a remote call in between, to make the program thread-safe
1275                 //we will need to verify if the connectivity of this node changed
1276                 int *elems1, numElems1;
1277                 m->n2e_getAll(nodes[j], elems1, numElems1);
1278                 if(numElems1!=numElems) {
1279                   CkAssert(false);
1280                 }
1281                 else {
1282                   if(numElems1!=0) delete[] elems1;
1283                 }
1285                 //add this to the list of ghost nodes to be deleted on the remote chunk
1286                 if(shouldBeDeleted) {
1287                   //convert this local index to a shared index
1288                   int shidx = m->getfmMM()->fmUtil->exists_in_IDXL(m,nodes[j],chk,1);
1289                   if(shidx!=-1) {
1290 #ifdef DEBUG_2
1291   CkPrintf("removeElement, line %d\n", __LINE__);
1292 #endif
1293                     m->node.ghostSend.removeNode(nodes[j], chk);
1294 #ifdef DEBUG_2
1295   CkPrintf("removeElement, line %d\n", __LINE__);
1296 #endif
1297                     ghostIndices[numGhostNodes] = shidx;
1298                     numGhostNodes++;
1299                   }
1300                 }
1301               }
1302               //try to find out which ghosts this chunk will be losing
1303               //if I am losing this node, then only will I lose some ghosts
1304               if(losingThisNode[j]) {
1305                 //was either a shared node, but now it should be a ghost node
1306                 //or it was a local node(corner), but now should be a ghost
1307                 for(int k=0; k<numElems; k++) {
1308                   int *nds = (int*)malloc(m->elem[elemtype].getConn().width()*sizeof(int));
1309                   if(FEM_Is_ghost_index(elems[k]) && m->getfmMM()->getfmUtil()->exists_in_IDXL(m,elems[k],chk,4,elemtype)!=-1) {
1310                     //this ghost element might be deleted, if it is not connected to any
1311                     //other node that is local or shared to this chunk
1312                     m->e2n_getAll(elems[k], nds, elemtype);
1313                     int geShouldBeDeleted = 1;
1314                     //look at the connectivity of this element
1315                     for(int l=0; l<connSize; l++) {
1316                       //if a node is not a ghost but is about to be lost anyway, 
1317                       //then it should be deleted
1318                       if(!FEM_Is_ghost_index(nds[l]) && (nodes[j]!=nds[l])) {
1319                         if(losingThisNode[(j+1)%connSize]==1 && nds[l]==nodes[(j+1)%connSize]) continue;
1320                         else if(losingThisNode[(j+2)%connSize]==1 && nds[l]==nodes[(j+2)%connSize]) continue;
1321                         geShouldBeDeleted = 0;
1322                       }
1323                     }
1324                     //the element will be deleted only if all the nodes are about to be lost
1325                     //here we delete this ghost element and clean the idxl lists as well
1326                     //and mark the element invalid in the fem_elem list
1327                     if(geShouldBeDeleted) {
1328                       int sge = m->getfmMM()->getfmUtil()->exists_in_IDXL(m,elems[k],chk,4,elemtype);
1329 #ifdef DEBUG_2
1330   CkPrintf("removeElement, line %d\n", __LINE__);
1331 #endif
1332                       m->elem[elemtype].ghost->ghostRecv.removeNode(FEM_To_ghost_index(elems[k]), chk);
1333 #ifdef DEBUG_2
1334   CkPrintf("removeElement, line %d\n", __LINE__);
1335 #endif
1336                       FEM_remove_element_local(m,elems[k],elemtype);
1337                       m->elem[elemtype].ghost->set_invalid(FEM_From_ghost_index(elems[k]),false);
1338                       ghostREIndices.push_back(sge);
1339                       numGhostRE++;
1340                       //find out if any nodes need to be removed
1341                       for(int l=0; l<connSize; l++) {
1342                         int *elts;
1343                         int numElts;
1344                         m->n2e_getAll(nds[l], elts, numElts);
1345                         //if this is no longer connected to this chunk through 
1346                         //any of its neighboring elements, it should no longer be a ghost
1347 // NOTE:
1348 // for edge flips, this will remove needed nodes on processor boundaries. The right
1349 // solution is to properly figure out when that situation applies, but for now just keep
1350 // the nodes always.
1351 //#ifndef CPSD
1352                         if(nds[l]<-1) { //it is a ghost
1353                           bool removeflag = true;
1354                           for(int lm=0; lm<numElts; lm++) {
1355                             //if it is connected to any element on this chunk, be it 
1356                             //local or ghost, other than the element that is about
1357                             //to be deleted, then it will still be a ghost
1358                             if(elts[lm]!=elems[k]) {
1359                               removeflag = false;
1360                             }
1361                           }
1362                           if(removeflag) {
1363                             //remove this ghost node on this chunk
1364                             int sgn = m->getfmMM()->getfmUtil()->exists_in_IDXL(m,nds[l],chk,2);
1365 #ifdef DEBUG_2
1366   CkPrintf("removeElement, line %d\n", __LINE__);
1367 #endif
1368                             m->node.ghost->ghostRecv.removeNode(FEM_To_ghost_index(nds[l]), chk);
1369 #ifdef DEBUG_2
1370   CkPrintf("removeElement, line %d\n", __LINE__);
1371 #endif
1372                             if(numElts==0) FEM_remove_node_local(m,nds[l]);
1373                             ghostRNIndices.push_back(sgn);
1374                             numGhostRN++;
1375                           }
1376                         }
1377 //#endif
1378                         if(numElts!=0) delete[] elts;
1379                       }
1380                     }
1381                   }
1382                   free(nds);
1383                 }
1384                 //if it is shared, tell that chunk, it no longer exists on me
1385                 int ssn = m->getfmMM()->getfmUtil()->exists_in_IDXL(m,nodes[j],chk,0);
1386                 //add a new coarse idxl lock, should lock a remote idxl only once
1387                 //even if it loses more than one nodes
1388                 if(!lockedRemoteIDXL) m->getfmMM()->getfmUtil()->idxllock(m,chk,0);
1389                 lockedRemoteIDXL = true;
1390                 bool losinglocal = false;
1391                 if(ssn==-1 && chk==permanent) {
1392                   //seems to me that if a local node becomes a ghost because of an element
1393                   //being eaten by a chunk, then it has to be a corner
1394                   //the other chunk must have this as a ghost
1395                   losinglocal = true;
1396                   ssn = -1000000000;
1397                   ssn += m->getfmMM()->getfmUtil()->exists_in_IDXL(m,nodes[j],chk,1);
1398                   if(willBeGhost[j]==1) {
1399                     m->node.ghost->ghostRecv.addNode(newghost[j],chk);
1400                   }
1401                   else { //happens when it was local but won't even be a ghost now
1402                     ssn -= 500000000; //to tell the other chunk, not to add it in ghostsend
1403                   }
1404                   //m->node.ghostSend.removeNode(nodes[j], chk);
1405                   sharedIndices[numSharedNodes] = ssn;
1406                   numSharedNodes++;
1407                 }
1408                 else 
1409                     if(ssn!=-1) {
1410                   if(willBeGhost[j]==1) {
1411                     m->node.ghost->ghostRecv.addNode(newghost[j],chk);
1412                   }
1413                   else { //happens when it was shared but won't even be a ghost now
1414                     ssn -= 500000000; //to tell the other chunk, not to add it in ghostsend
1415                   }
1416 #ifdef DEBUG_2
1417   CkPrintf("removeElement, line %d\n", __LINE__);
1418 #endif
1419                   m->node.shared.removeNode(nodes[j], chk);
1420 #ifdef DEBUG_2
1421   CkPrintf("removeElement, line %d\n", __LINE__);
1422 #endif
1423                   sharedIndices[numSharedNodes] = ssn;
1424                   numSharedNodes++;
1425                 }
1426                 if(i==numSharedChunks-1) { //this is the last time
1427 #ifdef DEBUG_2
1428   CkPrintf("removeElement, line %d\n", __LINE__);
1429 #endif
1430                   //so, update this element's adjacencies now
1431                   //since this node is about to be lost, so update the neighboring 
1432                   //elements and nodes that this node has a connectivity with,
1433                   //about the new node that will replace it
1434                   //add only to the elems which still exist
1435                   //if this node will not be a ghost, then these neighboring nodes and elems
1436                   //to this node do not really need updating, but I believe it will not
1437                   //hurt, because the connectivity in this case will be -1
1438                   int *elems1, numElems1;
1439                   m->n2e_getAll(nodes[j], elems1, numElems1);
1440                   for(int k=0; k<numElems1; k++) {
1441                     bool flagElem = false;
1442                     if(FEM_Is_ghost_index(elems1[k])) {
1443                       if(m->elem[0].ghost->is_valid(FEM_From_ghost_index(elems1[k]))) {
1444                         flagElem = true;
1445                       }
1446                     }
1447                     else {
1448                       if(m->elem[0].is_valid(elems1[k])) flagElem = true;
1449                     }
1450                     if(flagElem) {
1451                       m->e2n_replace(elems1[k],nodes[j],ghostidx[j],elemtype);
1452                       m->n2e_remove(nodes[j],elems1[k]);
1453                       m->n2e_add(ghostidx[j],elems1[k]);
1454                       testelems.push_back(elems1[k]);
1455                     }
1456                   }
1457                   if(numElems1!=0) delete[] elems1;
1458                   int *n2ns;
1459                   int numn2ns;
1460                   m->n2n_getAll(nodes[j],n2ns,numn2ns);
1461                   for(int k=0; k<numn2ns; k++) {
1462                     m->n2n_replace(n2ns[k],nodes[j],ghostidx[j]);
1463                     m->n2n_remove(nodes[j],n2ns[k]);
1464                     m->n2n_add(ghostidx[j],n2ns[k]);
1465                   }
1466                   if(numn2ns!=0) delete[] n2ns;
1467                   //if it is not shared, then all shared entries have been deleted
1468                   const IDXL_Rec *irec1 = m->node.shared.getRec(nodes[j]);
1469 #ifdef DEBUG_2
1470   CkPrintf("removeElement, line %d\n", __LINE__);
1471 #endif
1472 // This code removes nodes that we still need in cases involving acquisition/flip for
1473 // CPSD. I don't really follow the logic here, so I'm just cutting this bit out for now
1474 //#ifndef CPSD
1475           CkPrintf("removeElement, aggressive_node_removal: %d\n", aggressive_node_removal);
1476           if (aggressive_node_removal) {
1477               if(!irec1) {
1478                   if(!losinglocal) {
1479                       FEM_remove_node_local(m,nodes[j]);
1480                   }
1481                   else {
1482                       deleteNodeLater = nodes[j];
1483                   }
1484 #ifdef DEBUG_2
1485                   CkPrintf("removeElement, line %d\n", __LINE__);
1486 #endif
1487                   //if losing a local node, then can remove it only after its attributes 
1488                   //have been copied, since this is the only copy.. 
1489                   //(as no other chunk has this node as local/shared)
1490                   //so we'll need to delete the ghostsend idxl entry and the node later
1491               }
1492           }
1493 //#endif
1494                 }
1495               }
1496               if(numElems!=0) delete[] elems;
1497             }
1498             //now that all ghost nodes to be removed have been decided, we add the elem,
1499             // and call the entry method test if all the elems added are valid
1500             //sanity testing code START
1501 #ifdef DEBUG_2
1502   CkPrintf("removeElement, line %d\n", __LINE__);
1503 #endif
1504             for(int testelemscnt=0; testelemscnt <testelems.size(); testelemscnt++) {
1505               int el = testelems[testelemscnt];
1506               if(FEM_Is_ghost_index(el)) {
1507                 CkAssert(m->elem[0].ghost->is_valid(FEM_From_ghost_index(el))==1);
1508               }
1509               else {
1510                 CkAssert(m->elem[0].is_valid(el)==1);
1511               }
1512             }
1513 #ifdef DEBUG_2
1514   CkPrintf("removeElement, line %d\n", __LINE__);
1515 #endif
1516             //sanity testing code END
1517           }
1518 #ifdef DEBUG_2
1519   CkPrintf("removeElement, line %d\n", __LINE__);
1520 #endif
1521           removeGhostElemMsg *rm = new (numGhostNodes, numGhostRN, numGhostRE, numSharedNodes, 0) removeGhostElemMsg;
1522           rm->chk = index;
1523           rm->elemtype = elemtype;
1524           rm->elementid = sharedIdx;
1525           rm->numGhostIndex = numGhostNodes;
1526           for(int j=0; j<numGhostNodes; j++) {
1527             rm->ghostIndices[j] = ghostIndices[j];
1528           }
1529           rm->numGhostRNIndex = numGhostRN;
1530           for(int j=0; j<numGhostRN; j++) {
1531             rm->ghostRNIndices[j] = ghostRNIndices[j];
1532           }
1533           rm->numGhostREIndex = numGhostRE;
1534           for(int j=0; j<numGhostRE; j++) {
1535             rm->ghostREIndices[j] = ghostREIndices[j];
1536           }
1537           rm->numSharedIndex = numSharedNodes;
1538           for(int j=0; j<numSharedNodes; j++) {
1539             rm->sharedIndices[j] = sharedIndices[j];
1540           }
1541 #ifdef DEBUG_2
1542   CkPrintf("removeElement, line %d\n", __LINE__);
1543 #endif
1544           meshMod[chk].removeGhostElem(rm);  //update the ghosts on all shared chunks
1545 #ifdef DEBUG_2
1546   CkPrintf("removeElement, line %d\n", __LINE__);
1547 #endif
1548           //can delete the ghostSend and node now, if it was not deleted earlier
1549           if((i==numSharedChunks-1) && (deleteNodeLater!=-1)) {
1550 #ifdef DEBUG_2
1551   CkPrintf("removeElement, line %d\n", __LINE__);
1552 #endif
1553             m->node.ghostSend.removeNode(deleteNodeLater, permanent);
1554             FEM_remove_node_local(m,deleteNodeLater);
1555 #ifdef DEBUG_2
1556   CkPrintf("leaving removeElement, line %d\n", __LINE__);
1557 #endif
1558             //finally if it is on the fmfixednodes list, remove it
1559             //an O(n) search, but fmfixedNodes is a small list and called rarely
1560             for(int cnt=0; cnt<m->getfmMM()->fmfixedNodes.size(); cnt++) {
1561               if(m->getfmMM()->fmfixedNodes[cnt]==deleteNodeLater) {
1562                 m->getfmMM()->fmfixedNodes.remove(cnt);
1563                 break;
1564               }
1565             }
1566           }
1567           if(lockedRemoteIDXL) m->getfmMM()->getfmUtil()->idxlunlock(m,chk,0);
1568           free(ghostIndices);
1569           free(sharedIndices);
1570         }
1571         free(chknos);
1572         free(inds);
1573         if(permanent>=0) {
1574           delete [] newghost;
1575           delete [] ghostidx;
1576           delete [] losingThisNode;
1577           delete [] nodes;
1578           delete [] willBeGhost;
1579         }
1580       }
1581     }
1582     // remove local element
1583 #ifdef DEBUG_2
1584   CkPrintf("removeElement, line %d\n", __LINE__);
1585 #endif
1586     FEM_remove_element_local(m, elementid, elemtype);
1587 #ifdef DEBUG_2
1588   CkPrintf("removeElement, line %d\n", __LINE__);
1589 #endif
1590   }
1591 #ifdef DEBUG_2
1592   CkPrintf("leaving removeElement, line %d\n", __LINE__);
1593 #endif
1594   return index;
1599 /** Purge an element
1600     works for both local and ghost elements. Gets rid of any IDXL connectivity that might
1601     be left behind and cleans it on all chunks. Then it invalidates the element entry
1602     in the element table
1604 int FEM_purge_element(FEM_Mesh *m, int elementid, int elemtype) {
1605   if(elementid==-1) return 1;
1606   int index = m->getfmMM()->idx;
1607   if(FEM_Is_ghost_index(elementid)) {
1608     const IDXL_Rec *irec1 = m->elem[elemtype].ghost->ghostRecv.getRec(FEM_To_ghost_index(elementid));
1609     int remotechk = irec1->getChk(0);
1610     int sharedIdx = irec1->getIdx(0);
1611     meshMod[remotechk].purgeElement(index,sharedIdx);
1612   }
1613   else {
1614     //figure out all chunks where it is a ghost on
1615     const IDXL_Rec *irec = m->elem[elemtype].ghostSend.getRec(elementid);
1616     if(irec){
1617       //copy the (chunk,index) tuple
1618       int numSharedChunks = irec->getShared();
1619       int connSize = m->elem[elemtype].getConn().width();
1620       int *chknos, *inds;
1621       if(numSharedChunks>0) {
1622         chknos = (int*)malloc(numSharedChunks*sizeof(int));
1623         inds = (int*)malloc(numSharedChunks*sizeof(int));
1624         for(int i=0; i<numSharedChunks; i++) {
1625           chknos[i] = irec->getChk(i);
1626           inds[i] = irec->getIdx(i);
1627         }
1628       }
1629       //go and clean it up on each of these chunks
1630       for(int i=0; i<numSharedChunks; i++) {
1631         meshMod[chknos[i]].cleanupIDXL(index,inds[i]);
1632         m->elem[elemtype].ghostSend.removeNode(elementid, chknos[i]);
1633       }
1634       //free up allocated memory
1635       if(numSharedChunks>0) {
1636         free(chknos);
1637         free(inds);
1638       }
1639     }
1640   }
1641   // delete element by marking invalid
1642   if(!FEM_Is_ghost_index(elementid)){
1643     m->elem[elemtype].set_invalid(elementid,false);
1644   }
1645   return 1;
1647 //========================Basic Adaptivity operations=======================
1654 //========================Locking operations=======================
1655 /* These are the locking operations on nodes or chunks, depending on the locking scheme */
1656 //these are based on chunk locking: DEPRECATED
1657 CDECL int FEM_Modify_Lock(int mesh, int* affectedNodes, int numAffectedNodes, int* affectedElts, int numAffectedElts, int elemtype){
1658   return FEM_Modify_Lock(FEM_Mesh_lookup(mesh,"FEM_Modify_Lock"), affectedNodes, numAffectedNodes, affectedElts, numAffectedElts, elemtype);
1660 CDECL int FEM_Modify_Unlock(int mesh){
1661   return FEM_Modify_Unlock(FEM_Mesh_lookup(mesh,"FEM_Modify_Unlock"));
1663 //these are based on node locking
1664 CDECL int FEM_Modify_LockN(int mesh, int nodeId, int readLock) {
1665   return FEM_Modify_LockN(FEM_Mesh_lookup(mesh,"FEM_Modify_LockN"),nodeId, readLock);
1667 CDECL int FEM_Modify_UnlockN(int mesh, int nodeId, int readLock) {
1668   return FEM_Modify_UnlockN(FEM_Mesh_lookup(mesh,"FEM_Modify_UnlockN"),nodeId, readLock);
1672 /** coarse lock (lock chunks)
1673     find out which all chunks these nodes belong to and lock them all
1675 int FEM_Modify_Lock(FEM_Mesh *m, int* affectedNodes, int numAffectedNodes, int* affectedElts, int numAffectedElts, int elemtype) {
1676   return m->getfmMM()->getfmLock()->lock(numAffectedNodes, affectedNodes, numAffectedElts, affectedElts, elemtype);
1679 /** Unlock all chunks that was being held by a lock originating from here
1680  */
1681 int FEM_Modify_Unlock(FEM_Mesh *m) {
1682   return m->getfmMM()->getfmLock()->unlock();
1685 /** Lock a node
1686     get a 'read/write lock on 'nodeId'
1687     we find the chunk with the smallest id that owns/shares this node and lock this node
1688     on this chunk. This decreases the amount of communication and the number of 
1689     contention points for a single lock, as well as any chance of a deadlock
1690     since there is only a single resource to be acquired for a request
1692 int FEM_Modify_LockN(FEM_Mesh *m, int nodeId, int readlock) {
1693   int ret = -1;
1694   if(is_shared(m,nodeId)) {
1695     //find the index of the least chunk it is shared on and try to lock it
1696     int index = m->getfmMM()->getIdx();
1697     int numchunks;
1698     IDXL_Share **chunks1;
1699     m->getfmMM()->getfmUtil()->getChunkNos(0,nodeId,&numchunks,&chunks1);
1700     int minChunk = MAX_CHUNK;
1701     for(int j=0; j<numchunks; j++) {
1702       int pchk = chunks1[j]->chk;
1703       if(pchk < minChunk) minChunk = pchk;
1704     }
1705     for(int j=0; j<numchunks; j++) {
1706       delete chunks1[j];
1707     }
1708     if(numchunks!=0) free(chunks1);
1709     CkAssert(minChunk!=MAX_CHUNK);
1710     if(minChunk==index) {
1711       if(readlock) {
1712         ret = m->getfmMM()->fmLockN[nodeId].rlock();
1713       } else {
1714         ret = m->getfmMM()->fmLockN[nodeId].wlock(index);
1715       }
1716       if(ret==1) {
1717         m->getfmMM()->fmLockN[nodeId].verifyLock();
1718       }
1719       return ret;
1720     }
1721     else {
1722       CkAssert(minChunk!=MAX_CHUNK);
1723       int sharedIdx = m->getfmMM()->getfmUtil()->exists_in_IDXL(m,nodeId,minChunk,0);
1724       if(sharedIdx < 0) return -1;
1725       intMsg* imsg;
1726       if(readlock) {
1727         imsg = meshMod[minChunk].lockRemoteNode(sharedIdx, index, 0, 1);
1728       } else {
1729         imsg = meshMod[minChunk].lockRemoteNode(sharedIdx, index, 0, 0);
1730       }
1731       ret = imsg->i;
1732       delete imsg;
1733       if(ret==1) {
1734         boolMsg* bmsg = meshMod[minChunk].verifyLock(index,sharedIdx,0);
1735         delete bmsg;
1736       }
1737       return ret;
1738     }
1739   }
1740   else if(FEM_Is_ghost_index(nodeId)) {
1741     //find the index of the least chunk that owns it & try to lock it
1742     int index = m->getfmMM()->getIdx();
1743     int numchunks;
1744     IDXL_Share **chunks1;
1745     m->getfmMM()->getfmUtil()->getChunkNos(0,nodeId,&numchunks,&chunks1);
1746     int minChunk = MAX_CHUNK;
1747     for(int j=0; j<numchunks; j++) {
1748       int pchk = chunks1[j]->chk;
1749       if(pchk == index) continue;
1750       if(pchk < minChunk) minChunk = pchk;
1751     }
1752     for(int j=0; j<numchunks; j++) {
1753       delete chunks1[j];
1754     }
1755     if(numchunks!=0) free(chunks1);
1756     if(minChunk==MAX_CHUNK) return -1;
1757     int sharedIdx = m->getfmMM()->getfmUtil()->exists_in_IDXL(m,nodeId,minChunk,2);
1758     if(sharedIdx < 0) return -1;
1759     intMsg* imsg;
1760     if(readlock) {
1761       imsg = meshMod[minChunk].lockRemoteNode(sharedIdx, index, 1, 1);
1762     } else {
1763       imsg = meshMod[minChunk].lockRemoteNode(sharedIdx, index, 1, 0);
1764     }
1765     ret = imsg->i;
1766     delete imsg;
1767     /*if(ret==1) {
1768       if(nodeId==-8 && index==4) {
1769         CkPrintf("Locking node %d on chunk %d\n",nodeId, minChunk);
1770       }
1771       boolMsg* bmsg = meshMod[minChunk].verifyLock(index, sharedIdx, 1);
1772       delete bmsg;
1773       }*/
1774     return ret;
1775   }
1776   else {
1777     if(readlock) {
1778       ret = m->getfmMM()->fmLockN[nodeId].rlock();
1779     } else {
1780       int index = m->getfmMM()->getIdx();
1781       ret = m->getfmMM()->fmLockN[nodeId].wlock(index);
1782     }
1783     /*if(ret==1) {
1784       m->getfmMM()->fmLockN[nodeId].verifyLock();
1785       }*/
1786     return ret;
1787   }
1788   return -1; //should not reach here
1791 /** Unlock this node
1792     find the chunk that holds the lock and unlock it on that chunk
1793     the type of the lock is also considered while unlocking a lock
1795 int FEM_Modify_UnlockN(FEM_Mesh *m, int nodeId, int readlock) {
1796   if(is_shared(m,nodeId)) {
1797     //find the index of the least chunk it is shared on and try to unlock it
1798     int index = m->getfmMM()->getIdx();
1799     int numchunks;
1800     IDXL_Share **chunks1;
1801     m->getfmMM()->getfmUtil()->getChunkNos(0,nodeId,&numchunks,&chunks1);
1802     int minChunk = MAX_CHUNK;
1803     for(int j=0; j<numchunks; j++) {
1804       int pchk = chunks1[j]->chk;
1805       if(pchk < minChunk) minChunk = pchk;
1806     }
1807     for(int j=0; j<numchunks; j++) {
1808       delete chunks1[j];
1809     }
1810     if(numchunks!=0) free(chunks1);
1811     if(minChunk==index) {
1812       if(readlock) {
1813         return m->getfmMM()->fmLockN[nodeId].runlock();
1814       } else {
1815         return m->getfmMM()->fmLockN[nodeId].wunlock(index);
1816       }
1817     }
1818     else {
1819       CkAssert(minChunk!=MAX_CHUNK);
1820       int sharedIdx = m->getfmMM()->getfmUtil()->exists_in_IDXL(m,nodeId,minChunk,0);
1821       intMsg* imsg;
1822       if(readlock) {
1823         imsg = meshMod[minChunk].unlockRemoteNode(sharedIdx, index, 0, 1);
1824       } else {
1825         imsg = meshMod[minChunk].unlockRemoteNode(sharedIdx, index, 0, 0);
1826       }
1827       int ret = imsg->i;
1828       delete imsg;
1829       return ret;
1830     }
1831   }
1832   else if(FEM_Is_ghost_index(nodeId)) {
1833     //find the index of the least chunk that owns it & try to unlock it
1834     int index = m->getfmMM()->getIdx();
1835     int numchunks;
1836     IDXL_Share **chunks1;
1837     m->getfmMM()->getfmUtil()->getChunkNos(0,nodeId,&numchunks,&chunks1);
1838     int minChunk = MAX_CHUNK;
1839     for(int j=0; j<numchunks; j++) {
1840       int pchk = chunks1[j]->chk;
1841       if(pchk == index) continue;
1842       if(pchk < minChunk) minChunk = pchk;
1843     }
1844     for(int j=0; j<numchunks; j++) {
1845       delete chunks1[j];
1846     }
1847     if(numchunks!=0) free(chunks1);
1848     if(minChunk==MAX_CHUNK) return -1;
1849     int sharedIdx = m->getfmMM()->getfmUtil()->exists_in_IDXL(m,nodeId,minChunk,2);
1850     intMsg* imsg;
1851     if(readlock) {
1852       imsg = meshMod[minChunk].unlockRemoteNode(sharedIdx, index, 1, 1);
1853     } else {
1854       imsg = meshMod[minChunk].unlockRemoteNode(sharedIdx, index, 1, 0);
1855     }
1856     int ret = imsg->i;
1857     delete imsg;
1858     return ret;
1859   }
1860   else {
1861     if(readlock) {
1862       return m->getfmMM()->fmLockN[nodeId].runlock();
1863     } else {
1864       int index = m->getfmMM()->getIdx();
1865       return m->getfmMM()->fmLockN[nodeId].wunlock(index);
1866     }
1867   }
1868   return -1; //should not reach here
1873 /** When a chunk is losing a node, the lock might need to be reassigned
1874     lock it on the minchunk other than 'this chunk'
1875     should only be called when 'this chunk' is losing the node as local/shared
1877 void FEM_Modify_LockAll(FEM_Mesh*m, int nodeId, bool lockall) {
1878   int index = m->getfmMM()->getIdx();
1879   int numchunks;
1880   const IDXL_Rec *irec = m->node.shared.getRec(nodeId);
1881   //if it is a corner node, it might not have been a shared node, so can't lock it on anything else
1882   if(irec) {
1883     numchunks = irec->getShared();
1884     if(!lockall) {
1885       int minchunk=MAX_CHUNK;
1886       int sharedIdx = -1;
1887       for(int i=0; i<numchunks; i++) {
1888         int pchk = irec->getChk(i); 
1889         if(pchk<minchunk) {
1890           minchunk = pchk;
1891           sharedIdx = irec->getIdx(i);
1892         }
1893       }
1894       CkAssert(minchunk!=MAX_CHUNK && sharedIdx!=-1);
1895       //lock it on this chunk, if not already locked
1896       intMsg* imsg = meshMod[minchunk].lockRemoteNode(sharedIdx, index, 0, 0);
1897       int done = imsg->i;
1898       delete imsg;
1899       //if done=-1, then it is already locked, otherwise we just locked it
1900     }
1901     else {
1902       for(int i=0; i<numchunks; i++) {
1903         int pchk = irec->getChk(i); 
1904         int sharedIdx = irec->getIdx(i);
1905         intMsg* imsg = meshMod[pchk].lockRemoteNode(sharedIdx, index, 0, 0);
1906         int done = imsg->i;
1907         delete imsg;
1908       }
1909       m->getfmMM()->fmLockN[nodeId].wlock(index);
1910     }
1911   }
1912   return;
1915 /**Deprecated: Update the lock on this node (by locking the newly formed node
1916    The locks will be cleared later
1917    must be a local node, lock it & then unlock it if needed
1918  */
1919 void FEM_Modify_LockUpdate(FEM_Mesh*m, int nodeId, bool lockall) {
1920   int index = m->getfmMM()->getIdx();
1921   int numchunks;
1922   const IDXL_Rec *irec = m->node.shared.getRec(nodeId);
1923   //if it is a corner, the new node might not be shared
1924   //nothing was locked, hence nothing needs to be unlocked too
1925   if(irec) {
1926     numchunks = irec->getShared();
1927     int minchunk=MAX_CHUNK;
1928     int minI=-1;
1929     for(int i=0; i<numchunks; i++) {
1930       int pchk = irec->getChk(i); 
1931       if(pchk<minchunk) {
1932         minchunk = pchk;
1933         minI=i;
1934       }
1935     }
1936     if(!lockall) {
1937       if(minchunk>index) {
1938         int prevminchunk=minchunk;
1939         minchunk=index;
1940         int sharedIdx = irec->getIdx(minI);
1941         CkAssert(prevminchunk!=MAX_CHUNK && sharedIdx!=-1);
1942         intMsg* imsg = meshMod[prevminchunk].unlockRemoteNode(sharedIdx, index, 0, 0);
1943         delete imsg;
1944       }
1945       else if(minchunk < index) {
1946         //unlock the previously acquired lock
1947         int sharedIdx = irec->getIdx(minI);
1948         intMsg* imsg = meshMod[minchunk].lockRemoteNode(sharedIdx, index, 0, 0);
1949         delete imsg;
1950         m->getfmMM()->fmLockN[nodeId].wunlock(index);
1951       }
1952     }
1953     else {
1954       if(minchunk>index) minchunk=index;
1955       if(minchunk!=index) {
1956         m->getfmMM()->fmLockN[nodeId].wunlock(index);
1957       }
1958       for(int i=0; i<numchunks; i++) {
1959         int pchk = irec->getChk(i);
1960         if(pchk!=minchunk) {
1961           int sharedIdx = irec->getIdx(i);
1962           intMsg* imsg = meshMod[pchk].unlockRemoteNode(sharedIdx, index, 0, 0);
1963           delete imsg;
1964         }
1965       }
1966     }
1967   }
1968   return;
1971 /** Deprecated: For the newly acquired node, correct the lock by removing superfluous locks
1972     Should always be called on the new node to correct the locks
1973     should make sure that it is called before unlocking, i.e. before the entire operation completes.
1974     gets rid of the extra safety lock & makes the locks consistent
1976 void FEM_Modify_correctLockN(FEM_Mesh *m, int nodeId) {
1977   int index = m->getfmMM()->getIdx();
1978   int numchunks;
1979   IDXL_Share **chunks1;
1980   m->getfmMM()->getfmUtil()->getChunkNos(0,nodeId,&numchunks,&chunks1);
1981   int minChunk = MAX_CHUNK;
1982   int owner = -1;
1983   for(int j=0; j<numchunks; j++) {
1984     int pchk = chunks1[j]->chk;
1985     if(pchk < minChunk) minChunk = pchk;
1986   }
1987   if(is_shared(m,nodeId)) {
1988     for(int j=0; j<numchunks; j++) {
1989       int pchk = chunks1[j]->chk;
1990       if(pchk == index) owner = m->getfmMM()->fmLockN[nodeId].lockOwner();
1991       else {
1992         int sharedIdx = m->getfmMM()->getfmUtil()->exists_in_IDXL(m,nodeId,pchk,0);
1993         intMsg* imsg = meshMod[pchk].hasLockRemoteNode(sharedIdx, index, 0);
1994         owner = imsg->i;
1995         delete imsg;
1996       }
1997       if(owner != -1) { //this node is locked
1998         if(pchk == minChunk) {
1999           //the lock is good
2000           break;
2001         }
2002         else {
2003           //unlock the node on pchk & lock it on minChunk
2004           int locknodes = nodeId;
2005           int gotlocks = 1;
2006           int done = -1;
2007           if(pchk==index) m->getfmMM()->fmLockN[nodeId].wunlock(index);
2008           else {
2009             int sharedIdx = m->getfmMM()->getfmUtil()->exists_in_IDXL(m,nodeId,pchk,0);
2010             intMsg* imsg = meshMod[pchk].unlockRemoteNode(sharedIdx, index, 0, 0);
2011             delete imsg;
2012           }
2013           if(minChunk==index) done = m->getfmMM()->fmLockN[nodeId].wlock(index);
2014           else {
2015             CkAssert(minChunk!=MAX_CHUNK);
2016             int sharedIdx = m->getfmMM()->getfmUtil()->exists_in_IDXL(m,nodeId,minChunk,0);
2017             intMsg* imsg = meshMod[minChunk].lockRemoteNode(sharedIdx, index, 0, 0);
2018             done = imsg->i;
2019             delete imsg;
2020           }
2021           break;
2022         }
2023       }
2024     }
2025   }
2026   else if(FEM_Is_ghost_index(nodeId)) {
2027     //someone must have the lock
2028     for(int j=0; j<numchunks; j++) {
2029       int pchk = chunks1[j]->chk;
2030       int sharedIdx = m->getfmMM()->getfmUtil()->exists_in_IDXL(m,nodeId,pchk,2);
2031       intMsg* imsg = meshMod[pchk].hasLockRemoteNode(sharedIdx, index, 1);
2032       owner = imsg->i;
2033       delete imsg;
2034       if(owner != -1) { //this node is locked
2035         if(pchk == minChunk) {
2036           //the lock is good
2037           break;
2038         }
2039         else {
2040           //unlock the node on pchk & lock it on minChunk
2041           int locknodes = nodeId;
2042           int gotlocks = 1;
2043           int done = -1;
2044           int sharedIdx = m->getfmMM()->getfmUtil()->exists_in_IDXL(m,nodeId,pchk,2);
2045           intMsg* imsg = meshMod[pchk].unlockRemoteNode(sharedIdx, index, 1, 0);
2046           delete imsg;
2047           gotlocks=-1;
2048           while(done==-1) {
2049             done = m->getfmMM()->fmAdaptL->lockNodes(&gotlocks, &locknodes, 0, &locknodes, 1);
2050           }
2051           break;
2052         }
2053       }
2054     }
2055     if(owner==-1) {
2056       int locknodes = nodeId;
2057       int done = -1;
2058       int gotlocks=-1;
2059       while(done==-1) {
2060         done = m->getfmMM()->fmAdaptL->lockNodes(&gotlocks, &locknodes, 0, &locknodes, 1);
2061       }
2062     }
2063   }
2064   for(int j=0; j<numchunks; j++) {
2065     delete chunks1[j];
2066   }
2067   if(numchunks!=0) free(chunks1);
2068   return;
2070 //========================Locking operations=======================
2075 //=======================Read/Write operations=======================
2076 /* The following are functions to read/write data from/to the mesh */
2077 /** wrapper function that takes a 'mesh pointer' instead of meshId, and gets/sets the 
2078     data corresponding to the this 'attr' of this 'entity'
2079     'firstItem' is the starting point in the list, 'length' is the number of entries
2080     'datatype' is the type of data while 'width' is the size of each data entry
2082 void FEM_Mesh_dataP(FEM_Mesh *fem_mesh,int entity,int attr,void *data, int firstItem,int length, int datatype,int width) {
2083   IDXL_Layout lo(datatype,width);
2084   FEM_Mesh_data_layoutP(fem_mesh,entity,attr,data,firstItem,length,lo);
2087 /** similar to the above function, but takes a 'layout' instead of datatype and width
2088  */
2089 void FEM_Mesh_data_layoutP(FEM_Mesh *fem_mesh,int entity,int attr,void *data, int firstItem,int length, IDXL_Layout_t layout) {
2090   const char *caller="FEM_Mesh_data_layout";
2091   FEM_Mesh_data_layoutP(fem_mesh,entity,attr,data,firstItem,length,IDXL_Layout_List::get().get(layout,caller));
2094 /** The helper function that both the above functions call
2095     does the actual setting/getting of the data in/from the mesh
2097 void FEM_Mesh_data_layoutP(FEM_Mesh *m,int entity,int attr,void *data, int firstItem,int length, const IDXL_Layout &layout) {
2098   const char *caller="FEM_Mesh_data";
2099   //FEMAPI(caller);
2100   //FEM_Mesh *m=FEM_Mesh_lookup(fem_mesh,caller);
2101   FEM_Attribute *a=m->lookup(entity,caller)->lookup(attr,caller);
2102   if (m->isSetting()) 
2103     a->set(data,firstItem,length,layout,caller);
2104   else /* m->isGetting()*/
2105     a->get(data,firstItem,length,layout,caller);
2110 /** Copy the essential attributes for this ghost node from remote chunks
2111     This would be only done for attributes of a node, currently for coords & boundary
2112  */
2113 void FEM_Ghost_Essential_attributes(FEM_Mesh *m, int coord_attr, int bc_attr, int nodeid) {
2114   femMeshModify *theMod = m->getfmMM();
2115   FEM_MUtil *util = theMod->getfmUtil();
2116   int index = theMod->idx;
2117   if(FEM_Is_ghost_index(nodeid)) {
2118     //build up the list of chunks it is ghost to
2119     int numchunks;
2120     IDXL_Share **chunks1;
2121     util->getChunkNos(0,nodeid,&numchunks,&chunks1);
2122     for(int j=0; j<numchunks; j++) {
2123       int chk = chunks1[j]->chk;
2124       if(chk==index) continue;
2125       int ghostidx = util->exists_in_IDXL(m,nodeid,chk,2);
2126       double2Msg *d = meshMod[chk].getRemoteCoord(index,ghostidx);
2127       intMsg *im = meshMod[chk].getRemoteBound(index,ghostidx);
2128       double *coord = new double[2];
2129       coord[0] = d->i; coord[1] = d->j;
2130       int bound = im->i;
2131       CkVec<FEM_Attribute *>*attrs = m->node.ghost->getAttrVec();
2132       for (int i=0; i<attrs->size(); i++) {
2133         FEM_Attribute *a = (FEM_Attribute *)(*attrs)[i];
2134         if (a->getAttr() == theMod->fmAdaptAlgs->coord_attr) {
2135           FEM_DataAttribute *d = (FEM_DataAttribute *)a;
2136           d->getDouble().setRow(FEM_From_ghost_index(nodeid),coord,0);
2137         }
2138         else if(a->getAttr() == FEM_BOUNDARY) {
2139           FEM_DataAttribute *d = (FEM_DataAttribute *)a;
2140           d->getInt().setRow(FEM_From_ghost_index(nodeid),bound);
2141         }
2142       }
2143       delete [] coord;
2144       for(int j=0; j<numchunks; j++) {
2145         delete chunks1[j];
2146       }
2147       if(numchunks != 0) free(chunks1);
2148       delete im;
2149       delete d;
2150       break;
2151     }
2152   }
2153   return;
2155 //=======================Read/Write operations=======================
2163 //=======================femMeshModify: shadow array=======================
2164 femMeshModify::femMeshModify(femMeshModMsg *fm) {
2165   numChunks = fm->numChunks;
2166   idx = fm->myChunk;
2167   fmLock = new FEM_lock(idx, this);
2168   fmUtil = new FEM_MUtil(idx, this);
2169   fmAdapt = NULL;
2170   fmAdaptL = NULL;
2171   fmAdaptAlgs = NULL;
2172   fmInp = NULL;
2173   fmMesh = NULL;
2174   delete fm;
2177 femMeshModify::femMeshModify(CkMigrateMessage *m) {
2178   tc = NULL;
2179   fmMesh = NULL;
2180   fmLock = new FEM_lock(this);
2181   fmUtil = new FEM_MUtil(this);
2182   fmInp = new FEM_Interpolate(this);
2183   fmAdapt = new FEM_Adapt(this);
2184   fmAdaptL = new FEM_AdaptL(this);
2185   fmAdaptAlgs = new FEM_Adapt_Algs(this);
2188 femMeshModify::~femMeshModify() {
2189   if(fmLock != NULL) {
2190     delete fmLock;
2191   }
2192   if(fmUtil != NULL) {
2193     delete fmUtil;
2194   }
2199 void femMeshModify::pup(PUP::er &p) {
2200   p|numChunks;
2201   p|idx;
2202   p|tproxy;
2203   fmLock->pup(p);
2204   p|fmLockN;
2205   p|fmIdxlLock;
2206   p|fmfixedNodes;
2207   fmUtil->pup(p);
2208   fmInp->pup(p);
2209   fmAdapt->pup(p);
2210   fmAdaptL->pup(p);
2211   fmAdaptAlgs->pup(p);
2214 enum {FEM_globalID=33};
2215 void femMeshModify::ckJustMigrated(void) {
2216   ArrayElement1D::ckJustMigrated();
2217   //set the pointer to fmMM
2218   tc = tproxy[idx].ckLocal();
2219   CkVec<TCharm::UserData> &v=tc->sud;
2220   FEM_chunk *c = (FEM_chunk*)(v[FEM_globalID].getData());
2221   fmMesh = c->getMesh("ckJustMigrated");
2222   fmMesh->fmMM = this;
2223   setPointersAfterMigrate(fmMesh);
2226 void femMeshModify::setPointersAfterMigrate(FEM_Mesh *m) {
2227   fmMesh = m;
2228   fmInp->FEM_InterpolateSetMesh(fmMesh);
2229   fmAdapt->FEM_AdaptSetMesh(fmMesh);
2230   fmAdaptL->FEM_AdaptLSetMesh(fmMesh);
2231   fmAdaptAlgs->FEM_AdaptAlgsSetMesh(fmMesh);
2232   for(int i=0; i<fmLockN.size(); i++) fmLockN[i].setMeshModify(this);
2237 /** Part of the initialization phase, create all the new objects
2238     Populate all data structures, include all locks
2239     It also computes all the fixed nodes and populates a data structure with it
2240  */
2241 void femMeshModify::setFemMesh(FEMMeshMsg *fm) {
2242   fmMesh = fm->m;
2243   tc = fm->t;
2244   tproxy = tc->getProxy();
2245   fmMesh->setFemMeshModify(this);
2246   fmAdapt = new FEM_Adapt(fmMesh, this);
2247   fmAdaptL = new FEM_AdaptL(fmMesh, this);
2248   fmAdaptAlgs = new FEM_Adapt_Algs(fmMesh, this);
2249   fmInp = new FEM_Interpolate(fmMesh, this);
2250   //populate the node locks
2251   int nsize = fmMesh->node.size();
2252   for(int i=0; i<nsize; i++) {
2253     fmLockN.push_back(FEM_lockN(i,this));
2254   }
2255   /*int gsize = fmMesh->node.ghost->size();
2256     for(int i=0; i<gsize; i++) {
2257     fmgLockN.push_back(new FEM_lockN(FEM_To_ghost_index(i),this));
2258     }*/
2259   for(int i=0; i<numChunks; i++) {
2260     fmIdxlLock.push_back(false);
2261   }
2262   //compute all the fixed nodes
2263   for(int i=0; i<nsize; i++) {
2264     if(fmAdaptL->isCorner(i)) {
2265       fmfixedNodes.push_back(i);
2266     }
2267   }
2268   delete fm;
2269   return;
2274 /* Coarse chunk locks */
2275 intMsg *femMeshModify::lockRemoteChunk(int2Msg *msg) {
2276   CtvAccess(_curTCharm) = tc;
2277   intMsg *imsg = new intMsg(0);
2278   int ret = fmLock->lock(msg->i, msg->j);
2279   imsg->i = ret;
2280   delete msg;
2281   return imsg;
2284 intMsg *femMeshModify::unlockRemoteChunk(int2Msg *msg) {
2285   CtvAccess(_curTCharm) = tc;
2286   intMsg *imsg = new intMsg(0);
2287   int ret = fmLock->unlock(msg->i, msg->j);
2288   imsg->i = ret;
2289   delete msg;
2290   return imsg;
2293 /* node locks */
2294 intMsg *femMeshModify::lockRemoteNode(int sharedIdx, int fromChk, int isGhost, int readLock) {
2295   CtvAccess(_curTCharm) = tc;
2296   int localIdx;
2297   if(isGhost) {
2298     localIdx = fmUtil->lookup_in_IDXL(fmMesh, sharedIdx, fromChk, 1);
2299   } else {
2300     localIdx = fmUtil->lookup_in_IDXL(fmMesh, sharedIdx, fromChk, 0);
2301   }
2302   //CkAssert(localIdx != -1);
2303   intMsg *imsg = new intMsg(0);
2304   int ret;
2305   if(localIdx == -1) {
2306     ret = -1;
2307   }
2308   else {
2309     if(readLock) {
2310       ret = fmLockN[localIdx].rlock();
2311     } else {
2312       ret = fmLockN[localIdx].wlock(fromChk);
2313     }
2314   }
2315   imsg->i = ret;
2316   return imsg;
2319 intMsg *femMeshModify::unlockRemoteNode(int sharedIdx, int fromChk, int isGhost, int readLock) {
2320   CtvAccess(_curTCharm) = tc;
2321   int localIdx;
2322   if(isGhost) {
2323     localIdx = fmUtil->lookup_in_IDXL(fmMesh, sharedIdx, fromChk, 1);
2324   } else {
2325     localIdx = fmUtil->lookup_in_IDXL(fmMesh, sharedIdx, fromChk, 0);
2326   }
2327   CkAssert(localIdx != -1);
2328   intMsg *imsg = new intMsg(0);
2329   int ret;
2330   if(readLock) {
2331     ret = fmLockN[localIdx].runlock();
2332   } else {
2333     ret = fmLockN[localIdx].wunlock(fromChk);
2334   }
2335   imsg->i = ret;
2336   return imsg;
2341 chunkListMsg *femMeshModify::getChunksSharingGhostNode(int2Msg *i2m) {
2342   CtvAccess(_curTCharm) = tc;
2343   chunkListMsg *clm;
2344   clm = fmUtil->getChunksSharingGhostNodeRemote(fmMesh, i2m->i, i2m->j);
2345   delete i2m;
2346   return clm;
2351 intMsg *femMeshModify::addNodeRemote(addNodeMsg *msg) {
2352   CtvAccess(_curTCharm) = tc;
2353   intMsg *imsg = new intMsg(-1);
2354   //translate the indices
2355   int *localIndices = (int*)malloc(msg->nBetween*sizeof(int));
2356   int *chunks = (int*)malloc(msg->numChunks*sizeof(int));
2357   CkAssert(msg->numChunks==1);
2358   chunks[0] = msg->chunks[0];
2359   for(int i=0; i<msg->nBetween; i++) {
2360     localIndices[i] = fmUtil->lookup_in_IDXL(fmMesh, msg->between[i], chunks[0], 0);
2361     CkAssert(localIndices[i] != -1);
2362   }
2363   int ret = FEM_add_node(fmMesh, localIndices, msg->nBetween, chunks, msg->numChunks, msg->forceShared);
2364   //this is a ghost on that chunk,
2365   //add it to the idxl & update that guys idxl list
2366   fmMesh->node.ghostSend.addNode(ret,chunks[0]);
2367   imsg->i = fmUtil->exists_in_IDXL(fmMesh,ret,chunks[0],1);
2368   free(localIndices);
2369   delete msg;
2370   return imsg;
2373 void femMeshModify::addSharedNodeRemote(sharedNodeMsg *fm) {
2374   CtvAccess(_curTCharm) = tc;
2375   FEM_add_shared_node_remote(fmMesh, fm->chk, fm->nBetween, fm->between);
2376   delete fm;
2377   return;
2380 void femMeshModify::removeSharedNodeRemote(removeSharedNodeMsg *fm) {
2381   CtvAccess(_curTCharm) = tc;
2382   fmUtil->removeNodeRemote(fmMesh, fm->chk, fm->index);
2383   delete fm;
2384   return;
2387 void femMeshModify::removeGhostNode(int fromChk, int sharedIdx) {
2388   CtvAccess(_curTCharm) = tc;
2389   fmUtil->removeGhostNodeRemote(fmMesh, fromChk, sharedIdx);
2390   return;
2395 void femMeshModify::addGhostElem(addGhostElemMsg *fm) {
2396   CtvAccess(_curTCharm) = tc;
2397   fmUtil->addGhostElementRemote(fmMesh, fm->chk, fm->elemType, fm->indices, fm->typeOfIndex, fm->connSize);
2398   delete fm;
2399   return;
2402 intMsg *femMeshModify::addElementRemote(addElemMsg *fm) {
2403   CtvAccess(_curTCharm) = tc;
2404   int newEl = fmUtil->addElemRemote(fmMesh, fm->chk, fm->elemtype, fm->connSize, fm->conn, fm->numGhostIndex, fm->ghostIndices);
2405   //expecting that this element will always be local to this chunk
2406   CkAssert(!FEM_Is_ghost_index(newEl));
2407   int shidx = fmUtil->exists_in_IDXL(fmMesh,newEl,fm->chk,3);
2408   intMsg *imsg = new intMsg(shidx);
2409   delete fm;
2410   return imsg;
2413 void femMeshModify::removeGhostElem(removeGhostElemMsg *fm) {
2414   CtvAccess(_curTCharm) = tc;
2415   fmUtil->removeGhostElementRemote(fmMesh, fm->chk, fm->elementid, fm->elemtype, fm->numGhostIndex, fm->ghostIndices, fm->numGhostRNIndex, fm->ghostRNIndices, fm->numGhostREIndex, fm->ghostREIndices, fm->numSharedIndex, fm->sharedIndices);
2416   delete fm;
2417   return;
2420 void femMeshModify::removeElementRemote(removeElemMsg *fm) {
2421   CtvAccess(_curTCharm) = tc;
2422   fmUtil->removeElemRemote(fmMesh, fm->chk, fm->elementid, fm->elemtype, fm->permanent, fm->aggressive_node_removal);
2423   delete fm;
2424   return;
2429 intMsg *femMeshModify::eatIntoElement(int fromChk, int sharedIdx) {
2430   CtvAccess(_curTCharm) = tc;
2431   intMsg *imsg = new intMsg(-1);
2432   int localIdx = fmUtil->lookup_in_IDXL(fmMesh, sharedIdx, fromChk, 4);
2433   if(localIdx==-1) return imsg;
2434   int newEl = fmUtil->eatIntoElement(FEM_To_ghost_index(localIdx));
2435   int returnIdx = fmUtil->exists_in_IDXL(fmMesh, newEl, fromChk, 3);
2436   //returnIdx=-1 means that it was a discontinuous part
2437   if(returnIdx==-1) {
2438     //the other chunk has lost the entire region, unlock all the nodes involved
2439     //if there is a node 'fromChk' doesn't know about, 
2440     //but it is holding a lock on it, unlock it
2441     // nodesPerEl should be the number of nodes that can be adjacent to this element
2442     const int nodesPerEl = fmMesh->elem[0].getConn().width();
2443     int *adjnodes = new int[nodesPerEl];
2444     int numLocks = nodesPerEl + 1; //for 2D.. will be different for 3D
2445     int *lockednodes = new int[numLocks];
2446     fmMesh->e2n_getAll(newEl, adjnodes, 0);
2447     int newNode = -1;
2448     bool foundNewNode = false;
2449     //get the other node, which will be on an element that is an e2e of newEl
2450     //and which fromChk doesn't know abt, but has a lock on
2451     /*int *adjelems = new int[nodesPerEl];
2452     fmMesh->e2e_getAll(newEl, adjelems, 0);
2453     int *nnds = new int[nodesPerEl];
2454     for(int i=0; i<nodesPerEl; i++) {
2455       if(adjelems[i]!=-1) {
2456         fmMesh->e2n_getAll(adjelems[i], nnds, 0);
2457         for(int j=0; j<nodesPerEl; j++) {
2458           bool istarget = true;
2459           for(int k=0; k<nodesPerEl; k++) {
2460             if(nnds[j]==adjnodes[k]) {
2461               istarget = false;
2462             }
2463           }
2464           //get the owner of the lock
2465           if(istarget) {
2466             int lockowner = fmUtil->getLockOwner(nnds[j]);
2467             if(lockowner==fromChk) {
2468               bool knows = fmUtil->knowsAbtNode(fromChk,nnds[j]);
2469               //fromChk doesn't know abt this node
2470               if(!knows) {
2471                 newNode = nnds[j];
2472                 foundNewNode = true;
2473                 break;
2474               }
2475             }
2476           }
2477           if(foundNewNode) break;
2478         }
2479       }
2480       if(foundNewNode) break;
2481     }
2482     delete[] adjelems;
2483     delete[] nnds;
2484     */
2485     int *gotlocks = new int[numLocks];
2486     for(int i=0; i<nodesPerEl; i++) {
2487       lockednodes[i] = adjnodes[i];
2488       gotlocks[i] = 1;
2489     }
2490     if(foundNewNode) {
2491       lockednodes[nodesPerEl] = newNode;
2492       gotlocks[nodesPerEl] = 1;
2493     }
2494     else numLocks--;
2495     fmAdaptL->unlockNodes(gotlocks, lockednodes, 0, lockednodes, numLocks);
2496     delete[] lockednodes;
2497     delete[] gotlocks;
2498     delete[] adjnodes;
2499   }
2500   imsg->i = returnIdx;
2501   return imsg;
2506 intMsg *femMeshModify::getLockOwner(int fromChk, int sharedIdx) {
2507   CtvAccess(_curTCharm) = tc;
2508   int localIdx = fmUtil->lookup_in_IDXL(fmMesh, sharedIdx, fromChk, 0);
2509   int ret = fmUtil->getLockOwner(localIdx);
2510   intMsg *imsg = new intMsg(ret);
2511   return imsg;
2514 boolMsg *femMeshModify::knowsAbtNode(int fromChk, int toChk, int sharedIdx) {
2515   CtvAccess(_curTCharm) = tc;
2516   int localIdx = fmUtil->lookup_in_IDXL(fmMesh, sharedIdx, fromChk, 1);
2517   bool ret = fmUtil->knowsAbtNode(toChk,localIdx);
2518   boolMsg *bmsg = new boolMsg(ret);
2519   return bmsg;
2524 void femMeshModify::refine_flip_element_leb(int fromChk, int propElemT, int propNodeT,
2525                                             int newNodeT, int nbrOpNodeT,
2526                                             int nbrghost, double longEdgeLen)
2528   CtvAccess(_curTCharm) = tc;
2529   int propElem, propNode, newNode, nbrOpNode;
2530   propElem = getfmUtil()->lookup_in_IDXL(fmMesh, propElemT, fromChk, 3, 0);
2531   propNode = getfmUtil()->lookup_in_IDXL(fmMesh, propNodeT, fromChk, 0, -1);
2532   newNode = getfmUtil()->lookup_in_IDXL(fmMesh, newNodeT, fromChk, 0, -1);
2533   if(nbrghost) {
2534     nbrOpNode = getfmUtil()->lookup_in_IDXL(fmMesh, nbrOpNodeT, fromChk, 1,-1); 
2535   }
2536   else {
2537     nbrOpNode = getfmUtil()->lookup_in_IDXL(fmMesh, nbrOpNodeT, fromChk, 0,-1); 
2538   }
2539   fmAdaptAlgs->refine_flip_element_leb(propElem, propNode, newNode, nbrOpNode, longEdgeLen);
2540   return;
2545 void femMeshModify::addToSharedList(int fromChk, int sharedIdx) {
2546   CtvAccess(_curTCharm) = tc;
2547   fmUtil->addToSharedList(fmMesh, fromChk, sharedIdx);
2550 void femMeshModify::updateAttrs(updateAttrsMsg* umsg) {
2551   CtvAccess(_curTCharm) = tc;
2552   int localIdx = -1;
2553   if(!umsg->isnode) {
2554     localIdx = fmUtil->lookup_in_IDXL(fmMesh, umsg->sharedIdx, umsg->fromChk, 3);
2555   }
2556   else {
2557     if(!umsg->isGhost) {
2558       localIdx = fmUtil->lookup_in_IDXL(fmMesh, umsg->sharedIdx, umsg->fromChk, 0);
2559     }
2560     else localIdx = fmUtil->lookup_in_IDXL(fmMesh, umsg->sharedIdx, umsg->fromChk, 1);
2561   }
2562   CkAssert(localIdx!=-1);
2563   fmUtil->updateAttrs(umsg->data,umsg->datasize,localIdx,umsg->isnode,umsg->elemType);
2564   delete umsg;
2565   return;
2568 void femMeshModify::updateghostsend(verifyghostsendMsg *vmsg) {
2569   CtvAccess(_curTCharm) = tc;
2570   int localIdx = fmUtil->lookup_in_IDXL(fmMesh, vmsg->sharedIdx, vmsg->fromChk, 0);
2571   fmUtil->UpdateGhostSend(localIdx, vmsg->chunks, vmsg->numchks);
2572   delete vmsg;
2575 findgsMsg *femMeshModify::findghostsend(int fromChk, int sharedIdx) {
2576   CtvAccess(_curTCharm) = tc;
2577   int localIdx = fmUtil->lookup_in_IDXL(fmMesh, sharedIdx, fromChk, 0);
2578   int *chkl, numchkl=0;
2579   fmUtil->findGhostSend(localIdx, chkl, numchkl);
2580   findgsMsg *fmsg = new(numchkl)findgsMsg();
2581   fmsg->numchks = numchkl;
2582   for(int i=0; i<numchkl; i++) fmsg->chunks[i] = chkl[i];
2583   if(numchkl>0) delete[] chkl;
2584   return fmsg;
2587 intMsg *femMeshModify::getIdxGhostSend(int fromChk, int idxshared, int toChk) {
2588   CtvAccess(_curTCharm) = tc;
2589   int localIdx = fmUtil->lookup_in_IDXL(fmMesh, idxshared, fromChk, 0);
2590   int idxghostsend = -1;
2591   if(localIdx != -1) { 
2592     const IDXL_Rec *irec = fmMesh->node.ghostSend.getRec(localIdx);
2593     if(irec) {
2594       for(int i=0; i<irec->getShared(); i++) {
2595         if(irec->getChk(i) == toChk) {
2596           idxghostsend = fmUtil->exists_in_IDXL(fmMesh, localIdx, toChk, 1);
2597           break;
2598         }
2599       }
2600     }
2601   }
2602   intMsg *d = new intMsg(idxghostsend);
2603   return d;
2608 double2Msg *femMeshModify::getRemoteCoord(int fromChk, int ghostIdx) {
2609   CtvAccess(_curTCharm) = tc;
2610   //int localIdx = fmUtil->lookup_in_IDXL(fmMesh, ghostIdx, fromChk, 1);
2611   if(ghostIdx == fmMesh->node.ghostSend.addList(fromChk).size()) {
2612     double2Msg *d = new double2Msg(-2.0,-2.0);
2613     return d;
2614   }
2615   else {
2616     int localIdx = fmMesh->node.ghostSend.addList(fromChk)[ghostIdx];
2617     double coord[2];
2618     FEM_Mesh_dataP(fmMesh, FEM_NODE, fmAdaptAlgs->coord_attr, coord, localIdx, 1, FEM_DOUBLE, 2);
2619     double2Msg *d = new double2Msg(coord[0], coord[1]);
2620     return d;
2621   }
2624 intMsg *femMeshModify::getRemoteBound(int fromChk, int ghostIdx) {
2625   CtvAccess(_curTCharm) = tc;
2626   int localIdx = fmUtil->lookup_in_IDXL(fmMesh, ghostIdx, fromChk, 1);
2627   int bound;
2628   FEM_Mesh_dataP(fmMesh, FEM_NODE, FEM_BOUNDARY, &bound, localIdx, 1, FEM_INT, 1);
2629   intMsg *d = new intMsg(bound);
2630   return d;
2635 void femMeshModify::updateIdxlList(int fromChk, int idxTrans, int transChk) {
2636   CtvAccess(_curTCharm) = tc;
2637   int idxghostrecv = fmUtil->lookup_in_IDXL(fmMesh, idxTrans, transChk, 2);
2638   CkAssert(idxghostrecv != -1);
2639   fmMesh->node.ghost->ghostRecv.addNode(idxghostrecv,fromChk);
2640   return;
2643 void femMeshModify::removeIDXLRemote(int fromChk, int sharedIdx, int type) {
2644   CtvAccess(_curTCharm) = tc;
2645   int localIdx = fmUtil->lookup_in_IDXL(fmMesh, sharedIdx, fromChk, type);
2646   //CkAssert(localIdx != -1);
2647   if(localIdx!=-1) {
2648     fmMesh->node.ghostSend.removeNode(localIdx,fromChk);
2649   }
2650   return;
2653 void femMeshModify::addTransIDXLRemote(int fromChk, int sharedIdx, int transChk) {
2654   CtvAccess(_curTCharm) = tc;
2655   int localIdx = fmUtil->lookup_in_IDXL(fmMesh, sharedIdx, transChk, 0);
2656   CkAssert(localIdx != -1);
2657   fmMesh->node.ghostSend.addNode(localIdx,fromChk);
2658   return;
2661 void femMeshModify::verifyIdxlList(int fromChk, int size, int type) {
2662   CtvAccess(_curTCharm) = tc;
2663   fmUtil->verifyIdxlListRemote(fmMesh, fromChk, size, type);
2664   return;
2669 void femMeshModify::idxllockRemote(int fromChk, int type) {
2670   CtvAccess(_curTCharm) = tc;
2671   if(type==1) type = 2;
2672   else if(type ==2) type = 1;
2673   else if(type ==3) type = 4;
2674   else if(type ==4) type = 3;
2675   fmUtil->idxllockLocal(fmMesh, fromChk, type);
2676   return;
2679 void femMeshModify::idxlunlockRemote(int fromChk, int type) {
2680   CtvAccess(_curTCharm) = tc;
2681   if(type==1) type = 2;
2682   else if(type ==2) type = 1;
2683   else if(type ==3) type = 4;
2684   else if(type ==4) type = 3;
2685   fmUtil->idxlunlockLocal(fmMesh, fromChk, type);
2686   return;
2691 intMsg *femMeshModify::hasLockRemoteNode(int sharedIdx, int fromChk, int isGhost) {
2692   CtvAccess(_curTCharm) = tc;
2693   int localIdx;
2694   if(isGhost) {
2695     localIdx = fmUtil->lookup_in_IDXL(fmMesh, sharedIdx, fromChk, 1);
2696   } else {
2697     localIdx = fmUtil->lookup_in_IDXL(fmMesh, sharedIdx, fromChk, 0);
2698   }
2699   CkAssert(localIdx != -1);
2700   intMsg *imsg = new intMsg(0);
2701   int ret = fmLockN[localIdx].lockOwner();
2702   imsg->i = ret;
2703   return imsg;
2706 void femMeshModify::modifyLockAll(int fromChk, int sharedIdx) {
2707   CtvAccess(_curTCharm) = tc;
2708   int localIdx = fmUtil->lookup_in_IDXL(fmMesh, sharedIdx, fromChk, 1);
2709   FEM_Modify_LockAll(fmMesh, localIdx);
2710   return;
2713 boolMsg *femMeshModify::verifyLock(int fromChk, int sharedIdx, int isGhost) {
2714   CtvAccess(_curTCharm) = tc;
2715   int localIdx;
2716   if(isGhost) {
2717     localIdx = fmUtil->lookup_in_IDXL(fmMesh, sharedIdx, fromChk, 1);
2718   } else {
2719     localIdx = fmUtil->lookup_in_IDXL(fmMesh, sharedIdx, fromChk, 0);
2720   }
2721   //CkAssert(localIdx != -1);
2722   boolMsg *bmsg = new boolMsg(0);
2723   bool ret;
2724   if(localIdx == -1) {
2725     ret = false;
2726   }
2727   else {
2728     ret = fmLockN[localIdx].verifyLock();
2729   }
2730   bmsg->b = ret;
2731   return bmsg;
2734 void femMeshModify::verifyghostsend(verifyghostsendMsg *vmsg) {
2735   CtvAccess(_curTCharm) = tc;
2736   int localIdx = fmUtil->lookup_in_IDXL(fmMesh, vmsg->sharedIdx, vmsg->fromChk, 1);
2737   const IDXL_Rec *irec = fmMesh->node.shared.getRec(localIdx);
2738   if (irec!=NULL) {
2739     int numsh = irec->getShared();
2740     CkAssert(numsh==vmsg->numchks-1);
2741     for(int i=0; i<numsh; i++) {
2742       int ckl = irec->getChk(i);
2743       bool found = false;
2744       for(int j=0; j<numsh+1; j++) {
2745         if(vmsg->chunks[j]==ckl) {
2746           found = true; break;
2747         }
2748       }
2749       CkAssert(found);
2750     }
2751   }
2752   delete vmsg;
2755 boolMsg *femMeshModify::shouldLoseGhost(int fromChk, int sharedIdx, int toChk) {
2756   CtvAccess(_curTCharm) = tc;
2757   int localIdx = fmUtil->lookup_in_IDXL(fmMesh, sharedIdx, fromChk, 0);
2758   int *elems, numElems;
2759   fmMesh->n2e_getAll(localIdx, elems, numElems);
2760   bool shouldBeDeleted = true;
2761   for(int k=0; k<numElems; k++) {
2762     if(elems[k]>=0) {
2763       if(fmMesh->getfmMM()->fmUtil->exists_in_IDXL(fmMesh,elems[k],toChk,3)!=-1) {
2764         shouldBeDeleted = false;
2765         break;
2766       }
2767     }
2768   }
2769   if(numElems>0) delete[] elems;
2770   boolMsg *bmsg = new boolMsg(shouldBeDeleted);
2771   return bmsg;
2776 /** The node index is the sharedIdx with 'fromChk'
2777     This local node index is then to be send as a ghost to 'toChk'
2779 void femMeshModify::addghostsendl(int fromChk, int sharedIdx, int toChk, int transIdx) {
2780   CtvAccess(_curTCharm) = tc;
2781   int localIdx = fmUtil->lookup_in_IDXL(fmMesh, sharedIdx, fromChk, 0);
2782   int sharedghost = fmUtil->exists_in_IDXL(fmMesh,localIdx,toChk,1);
2783   if(sharedghost==-1) { //it needs to be added from this chunk
2784     //lock idxl
2785     fmMesh->node.ghostSend.addNode(localIdx,toChk);
2786     meshMod[toChk].addghostsendl1(idx,fromChk,transIdx);
2787     //unlock idxl
2788   }
2791 /** The node index is obtained as a the ghost received from 'transChk' at 'transIdx'
2792  */
2793 void femMeshModify::addghostsendl1(int fromChk, int transChk, int transIdx) {
2794   CtvAccess(_curTCharm) = tc;
2795   int localIdx = fmUtil->lookup_in_IDXL(fmMesh, transIdx, transChk, 2);
2796   fmMesh->node.ghost->ghostRecv.addNode(localIdx,fromChk);
2799 /** The node is the obtained on the ghost recv list from 'fromChk' at 'sharedIdx' index.
2800  */
2801 void femMeshModify::addghostsendr(int fromChk, int sharedIdx, int toChk, int transIdx) {
2802   CtvAccess(_curTCharm) = tc;
2803   int localIdx = fmUtil->lookup_in_IDXL(fmMesh, sharedIdx, fromChk, 2);
2804   int sharedghost = fmUtil->exists_in_IDXL(fmMesh,FEM_To_ghost_index(localIdx),toChk,2);
2805   if(sharedghost==-1) { //it needs to be added from this chunk
2806     //lock idxl
2807     fmUtil->idxllock(fmMesh,toChk,0);
2808     fmMesh->node.ghost->ghostRecv.addNode(localIdx,toChk);
2809     meshMod[toChk].addghostsendr1(idx,fromChk,transIdx);
2810     fmUtil->idxlunlock(fmMesh,toChk,0);
2811     //unlock idxl
2812   }
2815 /** The local node index is obtained on the shared list of 'transChk' at 'transIdx'
2816  */
2817 void femMeshModify::addghostsendr1(int fromChk, int transChk, int transIdx) {
2818   CtvAccess(_curTCharm) = tc;
2819   int localIdx = fmUtil->lookup_in_IDXL(fmMesh, transIdx, transChk, 0);
2820   fmMesh->node.ghostSend.addNode(localIdx,fromChk);
2823 boolMsg *femMeshModify::willItLose(int fromChk, int sharedIdx) {
2824   CtvAccess(_curTCharm) = tc;
2825   int localIdx = fmUtil->lookup_in_IDXL(fmMesh, sharedIdx, fromChk, 3);
2826   int nnbrs[3];
2827   fmMesh->e2n_getAll(localIdx,nnbrs,0);
2828   //if it loses any node if it loses this element, then it should let fromChk, acquire this elem
2829   bool willlose = true;
2830   for(int i=0; i<3; i++) {
2831     int *enbrs, numenbrs;
2832     fmMesh->n2e_getAll(nnbrs[i],enbrs,numenbrs);
2833     willlose = true;
2834     for(int j=0; j<numenbrs; j++) {
2835       if(enbrs[j]>=0 && enbrs[j]!=localIdx) {
2836         willlose = false;
2837         break;
2838       }
2839     }
2840     if(numenbrs>0) delete [] enbrs;
2841     if(willlose) break;
2842   }
2843   boolMsg *bmsg = new boolMsg(willlose);
2844   return bmsg;
2849 /** The two local elements indices are obtained from the ghost send list to 'fromChk'
2850     Data is copied from the first element to the second
2852 void femMeshModify::interpolateElemCopy(int fromChk, int sharedIdx1, int sharedIdx2) {
2853   CtvAccess(_curTCharm) = tc;
2854   int localIdx1 = fmUtil->lookup_in_IDXL(fmMesh, sharedIdx1, fromChk, 3);
2855   int localIdx2 = fmUtil->lookup_in_IDXL(fmMesh, sharedIdx2, fromChk, 3);
2856   CkAssert(localIdx1!=-1 && localIdx2!=-1);
2857   fmUtil->copyElemData(0,localIdx1,localIdx2);
2860 void femMeshModify::cleanupIDXL(int fromChk, int sharedIdx) {
2861   CtvAccess(_curTCharm) = tc;
2862   int localIdx = fmUtil->lookup_in_IDXL(fmMesh, sharedIdx, fromChk, 4);
2863   CkAssert(fmUtil->exists_in_IDXL(fmMesh,FEM_To_ghost_index(localIdx),fromChk,4)!=-1);
2864   fmMesh->elem[0].ghost->ghostRecv.removeNode(localIdx, fromChk);
2865   fmMesh->elem[0].getGhost()->set_invalid(localIdx,false);
2868 void femMeshModify::purgeElement(int fromChk, int sharedIdx) {
2869   CtvAccess(_curTCharm) = tc;
2870   int localIdx = fmUtil->lookup_in_IDXL(fmMesh, sharedIdx, fromChk, 3);
2871   CkAssert(localIdx!=-1);
2872   FEM_purge_element(fmMesh,localIdx,0);
2875 entDataMsg *femMeshModify::packEntData(int fromChk, int sharedIdx, bool isnode, int elemType) {
2876   CtvAccess(_curTCharm) = tc;
2877   int localIdx=-1;
2878   if(isnode) {
2879     localIdx = fmUtil->lookup_in_IDXL(fmMesh, sharedIdx, fromChk, 1);
2880   }
2881   else {
2882     localIdx = fmUtil->lookup_in_IDXL(fmMesh, sharedIdx, fromChk, 3);
2883   }
2884   CkAssert(localIdx!=-1);
2885   char *data; 
2886   int size, count;
2887   fmUtil->packEntData(&data, &size, &count, localIdx, isnode, elemType);
2888   entDataMsg *edm = new (size, 0) entDataMsg(count,size);
2889   memcpy(edm->data, data, size);
2890   free(data);
2891   return edm;
2894 boolMsg *femMeshModify::isFixedNodeRemote(int fromChk, int sharedIdx) {
2895   CtvAccess(_curTCharm) = tc;
2896   int localIdx = fmUtil->lookup_in_IDXL(fmMesh, sharedIdx, fromChk, 0);
2897   for(int i=0; i<fmfixedNodes.size(); i++) {
2898     if(fmfixedNodes[i]==localIdx) {
2899       return new boolMsg(true);
2900     }
2901   }
2902   return new boolMsg(false);
2907 //need these functions for debugging
2908 void femMeshModify::finish(void) {
2909   CkPrintf("[%d]finished this program!\n",idx);
2910   return;
2913 void femMeshModify::finish1(void) {
2914   meshMod[0].finish();
2915   return;
2918 #include "ParFUM_Adapt.def.h"
2919 /*@}*/