// Read this mesh from the FEM framework's mesh m and include ghosts void readGhostFEM(int m,TetMesh &t) { int nNode=FEM_Mesh_get_length(m,FEM_NODE); int ngNode=FEM_Mesh_get_length(m,FEM_NODE+FEM_GHOST); int nTet=FEM_Mesh_get_length(m,FEM_ELEM+0); int ngTet=FEM_Mesh_get_length(m,FEM_ELEM+0+FEM_GHOST); t.allocate(nTet+ngTet, nNode+ngNode); FEM_Mesh_data(m,FEM_NODE,FEM_COORD,t.getPointArray(),0,nNode,FEM_DOUBLE,3); FEM_Mesh_data(m,FEM_ELEM+0,FEM_CONN,t.getTetConn(),0,nTet,FEM_INDEX_0,4); // add ghost data if it exists if (ngTet > 0) { FEM_Mesh_data(m,FEM_ELEM+0+FEM_GHOST,FEM_CONN,t.getTet(nTet),0,ngTet, FEM_INDEX_0,4); if (ngNode > 0) FEM_Mesh_data(m,FEM_NODE+FEM_GHOST,FEM_COORD,t.getPoint(nNode),0,ngNode, FEM_DOUBLE,3); } t.nonGhostPt = nNode; t.nonGhostTet = nTet; // convert ghost node indices from negative to positive //printf("nNode=%d ngNode=%d nTet=%d ngTet=%d\n", nNode, ngNode, nTet, ngTet); if (ngTet > 0) { int *tconn; for (int i=nTet; i<ngTet+nTet; i++) { tconn = t.getTet(i); for (int j=0; j<4; j++) { //printf("Tet=%d conn[%d]=%d before conversion.\n", i, j, tconn[j]); if (tconn[j] < -1) tconn[j] = (tconn[j]*-1)-2+nNode; //printf("Tet=%d conn[%d]=%d after conversion.\n", i, j, tconn[j]); assert(tconn[j] < nNode+ngNode); assert(!(tconn[j] < -1)); } } } }
// Read this mesh from the FEM framework's mesh m void readFEM(int m,TetMesh &t) { int nNode=FEM_Mesh_get_length(m,FEM_NODE); int nTet=FEM_Mesh_get_length(m,FEM_ELEM+0); t.allocate(nTet, nNode); FEM_Mesh_data(m,FEM_NODE,FEM_COORD,t.getPointArray(),0,nNode,FEM_DOUBLE,3); FEM_Mesh_data(m,FEM_ELEM+0,FEM_CONN,t.getTetConn(),0,nTet,FEM_INDEX_0,4); }
void FEM_mesh_smooth(int mesh, int *nodes, int nNodes, int attrNo) { vector2d newPos, *coords, *ghostCoords,theCoord; int nNod, nGn, *boundVals, nodesInChunk; int neighbors[3], *adjnodes; int gIdxN; FEM_Mesh *meshP = FEM_Mesh_lookup(mesh, "driver"); nodesInChunk = FEM_Mesh_get_length(mesh,FEM_NODE); boundVals = new int[nodesInChunk]; nGn = FEM_Mesh_get_length(mesh, FEM_GHOST + FEM_NODE); coords = new vector2d[nodesInChunk+nGn]; FEM_Mesh_data(mesh, FEM_NODE, FEM_BOUNDARY, (int*) boundVals, 0, nodesInChunk, FEM_INT, 1); FEM_Mesh_data(mesh, FEM_NODE, attrNo, (double*)coords, 0, nodesInChunk, FEM_DOUBLE, 2); IDXL_Layout_t coord_layout = IDXL_Layout_create(IDXL_DOUBLE, 2); FEM_Update_ghost_field(coord_layout,-1, coords); ghostCoords = &(coords[nodesInChunk]); // FEM_Mesh_data(FEM_Mesh_default_write(), FEM_GHOST+FEM_NODE, attrNo, (double*)ghostCoords, 0, nGn, FEM_DOUBLE, 2); for (int i=0; i<nNodes; i++) { newPos.x=0; newPos.y=0; CkAssert(nodes[i]<nodesInChunk); if (FEM_is_valid(mesh, FEM_NODE, i) && boundVals[i]>-1) //node must be internal { meshP->n2n_getAll(i, &adjnodes, &nNod); // for (int j=0; j<nNod; j++) CkPrintf("node[%d]: %d\n", i,adjnodes[j]); for (int j=0; j<nNod; j++) { //for all adjacent nodes, find coords if (adjnodes[j]<-1) { gIdxN = FEM_From_ghost_index(adjnodes[j]); newPos.x += theCoord.x=ghostCoords[gIdxN].x; newPos.y += theCoord.y=ghostCoords[gIdxN].y; } else { newPos.x += theCoord.x=coords[adjnodes[j]].x; newPos.y += theCoord.y=coords[adjnodes[j]].y; } int rank=FEM_My_partition(); if (rank==0 && i==3 || rank==1 && i==17) { CkPrintf("node[%d]: %d\n", i,adjnodes[j]); CkPrintf("%d: (%f, %f)\n", j, theCoord.x, theCoord.y); } } newPos.x/=nNod; newPos.y/=nNod; FEM_set_entity_coord2(mesh, FEM_NODE, nodes[i], newPos.x, newPos.y); delete [] adjnodes; } } // FEM_Update_field(coord_layout, coords); // FEM_Mesh_data(FEM_Mesh_default_write(), FEM_NODE, attrNo, (double*)coords, 0, nodesInChunk, FEM_DOUBLE, 2); if (coords) delete [] coords; delete [] boundVals; }
void repeat_after_split(void *data){ myGlobals *g = (myGlobals *)data; g->nelems = FEM_Mesh_get_length(FEM_Mesh_default_read(),FEM_ELEM); g->nnodes = FEM_Mesh_get_length(FEM_Mesh_default_read(),FEM_NODE); for(int k=0;k<g->nnodes;k++){ printf(" node %d ( %.6f %.6f )\n",k,g->coord[k].x,g->coord[k].y); } calcMasses(*g); };
void FEM_REFINE2D_Newmesh(int meshID,int nodeID,int elemID,int nodeBoundary){ int nelems = FEM_Mesh_get_length(meshID,elemID); int nghost = FEM_Mesh_get_length(meshID,elemID+FEM_GHOST); int total = nghost + nelems; int *tempMesh = new int[3*total]; int nnodes = FEM_Mesh_get_length(meshID,nodeID); int *tempBoundaries=NULL; if(nodeBoundary){ tempBoundaries = new int[nnodes]; } FEM_Mesh_data(meshID,elemID,FEM_CONN,&tempMesh[0],0,nelems,FEM_INDEX_0,3); FEM_Mesh_data(meshID,elemID+FEM_GHOST,FEM_CONN,&tempMesh[3*nelems],0,nghost,FEM_INDEX_0,3); for(int t=nelems;t<total;t++){ for(int j=0;j<3;j++){ if(FEM_Is_ghost_index(tempMesh[3*t+j])){ tempMesh[3*t+j] += nelems; } } } /*Set up the global ID's, for refinement*/ int myID = FEM_My_partition(); int *gid=new int[2*total]; for (int i=0;i<nelems;i++) { gid[2*i+0]=myID; //Local element-- my chunk gid[2*i+1]=i; //Local number } int gid_fid=FEM_Create_field(FEM_INT,2,0,2*sizeof(int)); FEM_Update_ghost_field(gid_fid,0,gid); if(nodeBoundary){ FEM_Mesh_data(meshID,nodeID,FEM_BOUNDARY,tempBoundaries,0,nnodes,FEM_INT,1); printf("NODE BOUNDARIES-------------\n"); for(int i=0;i<nnodes;i++){ printf("%d %d \n",i,tempBoundaries[i]); } printf("------------------------\n"); } /*Set up refinement framework*/ /*FIX ME! PASS IN EDGE BOUNDARIES! */ const int *edgeConn = NULL; const int *edgeBounds = NULL; int nEdges = 0; REFINE2D_NewMesh(meshID,nelems,total,nnodes,(int *)tempMesh,gid,tempBoundaries,edgeBounds, edgeConn, nEdges); if(tempBoundaries){ delete [] tempBoundaries; } delete [] gid; delete [] tempMesh; }
void BulkAdapt::dumpConn() { FEM_Elem &elem = meshPtr->elem[0]; // elem is local elements int nelems=FEM_Mesh_get_length(meshID, FEM_ELEM+0); // Get number of elements for (int i=0; i<nelems; i++) { int *conn = elem.connFor(i); CkPrintf("[%d] %d: %d %d %d %d\n", partitionID, i, conn[0], conn[1], conn[2], conn[3]); } }
/** Extract all FEM communication information into Roccom format. */ static CkVec<int> getRoccomPconn(int fem_mesh,int *ghost_len,const int *paneFmChunk) { CkVec<int> pconn; // Shared nodes come first: getRoccomPconn(IDXL_Get_send(FEM_Comm_shared(fem_mesh,FEM_NODE)),0,pconn,paneFmChunk); int realLen=pconn.size(); // Sent ghost nodes: getRoccomPconn(IDXL_Get_send(FEM_Comm_ghost(fem_mesh,FEM_NODE)),0,pconn,paneFmChunk); // Received ghost nodes (use bias to switch to Roccom ghost node numbering) getRoccomPconn(IDXL_Get_recv(FEM_Comm_ghost(fem_mesh,FEM_NODE)), FEM_Mesh_get_length(fem_mesh,FEM_NODE),pconn,paneFmChunk); // Handle elements (much tougher!) // Find list of element types int elems[1024]; int e, ne=FEM_Mesh_get_entities(fem_mesh,elems); for (e=0;e<ne;e++) if (elems[e]<FEM_ELEM || elems[e]>=FEM_SPARSE) elems[e--]=elems[--ne]; // swap out bad entity with the end // Make one output IDXL that combines all element types: IDXL_t elghost=IDXL_Create(); int out_r=0, out_g=0; // output indices for real; ghost for (e=0;e<ne;e++) { IDXL_Combine(elghost,FEM_Comm_ghost(fem_mesh,elems[e]), out_r,out_g); out_r+=FEM_Mesh_get_length(fem_mesh,elems[e]); out_g+=FEM_Mesh_get_length(fem_mesh,elems[e]+FEM_GHOST); } // Sent ghost elements: getRoccomPconn(IDXL_Get_send(elghost),0,pconn,paneFmChunk); // Received ghost elements (shift all ghosts to start after real elements) getRoccomPconn(IDXL_Get_recv(elghost),out_r,pconn,paneFmChunk); IDXL_Destroy(elghost); if (ghost_len) *ghost_len=pconn.size()-realLen; return pconn; }
CDECL void FEM_Update_ghost_field(int fid, int elType, void *v_data) { int mesh=FEM_Mesh_default_read(); int entity; if (elType==-1) entity=FEM_NODE; else entity=FEM_ELEM+elType; IDXL_t src=FEM_Comm_ghost(mesh,entity); int nReal=FEM_Mesh_get_length(mesh,entity); int bytesPerRec=IDXL_Get_layout_distance(fid); char *data=(char *)v_data; IDXL_Comm_t comm=IDXL_Comm_begin(1,0); IDXL_Comm_send(comm,src,fid,data); IDXL_Comm_recv(comm,src,fid,&data[nReal*bytesPerRec]); //Begin recv'ing ghosts after all reals IDXL_Comm_wait(comm); }
void rebuildArrays (int mesh, myGlobals &g) { CkPrintf("Rebuilding arrays. \n"); delete [] g.conn; delete [] g.coord; delete [] g.vCoord; delete [] g.vConn; g.nelems=FEM_Mesh_get_length(mesh, FEM_ELEM); g.nnodes=FEM_Mesh_get_length(mesh, FEM_NODE); g.nVnodes = FEM_count_valid(mesh, FEM_NODE); g.nVelems = FEM_count_valid(mesh, FEM_ELEM); g.coord=new vector2d[g.nnodes]; g.conn=new connRec[g.nelems]; g.vConn = new connRec[g.nVelems]; g.vCoord = new vector2d[g.nVnodes]; FEM_Mesh_data(mesh, FEM_NODE, FEM_DATA, (double *)g.coord, 0, g.nnodes, FEM_DOUBLE, 2); int j=0; for (int i=0; i<g.nnodes;i++) { if (FEM_is_valid(mesh, FEM_NODE, i)) { // CkPrintf("g.coord[%d]: (%f, %f) \n", i, g.coord[i].x, g.coord[i].y); g.vCoord[j]=g.coord[i]; j++; } } j=0; for (int i=0; i<g.nelems;i++) { FEM_Mesh_lookup(mesh, "driver")->e2n_getAll(i, (int *)g.conn[i]); if (FEM_is_valid(mesh, FEM_ELEM, i)) { // CkPrintf("g.conn[%d]: (%d, %d, %d) \n", i, g.conn[i][0], g.conn[i][1], g.conn[i][2]); for (int k=0; k<3; k++) g.vConn[j][k]=g.conn[i][k]; j++; } } }
/// Recreate the shared nodes. An alternate incorrect version can be enabled by undefining CORRECT_COORD_COMPARISON void ParFUM_recreateSharedNodes(int meshid, int dim, MPI_Comm newComm) { #define CORRECT_COORD_COMPARISON MPI_Comm comm = newComm; int rank, nParts; int send_count=0; // sanity check int recv_count=0; // sanity check MPI_Comm_rank(comm, &rank); MPI_Comm_size(comm, &nParts); #if SUPER_FAST_SPECIFIC_TORUS #define TORUSY 15 #define TORUSZ 15 CkPrintf("rank %d is manually configuring the IDXL lists to make the shared node generation fast\n"); FEM_Mesh *mesh = (FEM_chunk::get("ParFUM_recreateSharedNodes"))->lookup(meshid,"ParFUM_recreateSharedNodes"); IDXL_Side &shared = mesh->node.shared; int low = (rank-1+nParts) % nParts; int high = (rank+1) % nParts; IDXL_List &list1 = shared.addList(low); IDXL_List &list2 = shared.addList(high); int nodesInPlane = TORUSY * TORUSZ; int numNodes = FEM_Mesh_get_length(meshid,FEM_NODE); // vp - 1 for(int j=0;j<nodesInPlane;j++){ list1.push_back(j); } // vp + 1 for(int j=0;j<nodesInPlane;j++){ list2.push_back(numNodes - nodesInPlane +j); } return; #else // Shared data will be temporarily stored in the following structure int *sharedNodeCounts; // sharedCounts[i] = number of nodes shared with rank i int **sharedNodeLists; // sharedNodes[i] is the list of nodes shared with rank i // Initialize shared data sharedNodeCounts = (int *)malloc(nParts*sizeof(int)); sharedNodeLists = (int **)malloc(nParts*sizeof(int *)); for (int i=0; i<nParts; i++) { sharedNodeLists[i] = NULL; sharedNodeCounts[i] = 0; } // Get local node count and coordinates int numNodes; int coord_msg_tag=42, sharedlist_msg_tag=43; double *nodeCoords; numNodes = FEM_Mesh_get_length(meshid,FEM_NODE); nodeCoords = (double *)malloc(dim*numNodes*sizeof(double)); FEM_Mesh_become_get(meshid); FEM_Mesh_data(meshid,FEM_NODE,FEM_COORD, nodeCoords, 0, numNodes,FEM_DOUBLE, dim); //MPI_Barrier(MPI_COMM_WORLD); if (rank==0) CkPrintf("Extracted node data...\n"); // Begin exchange of node coordinates to determine shared nodes // FIX ME: compute bounding box, only exchange when bounding boxes collide /// The highest partition # to which I send my coordinates(wraps around) int sendUpperBound; if(nParts %2==0){ sendUpperBound = rank + (nParts/2) - (rank%2); } else { sendUpperBound = rank + (nParts/2) ; } /// The lowest partition # to which I send my coordinates(wraps around) int sendLowerBound; if(nParts %2==0){ sendLowerBound = rank - (nParts/2) + ((rank+1)%2); } else { sendLowerBound = rank - (nParts/2); } // Special case optimization for when the mesh is generated in such a way that only neighboring partitions share nodes // look for command line argument #ifdef SHARED_NODES_ONLY_NEIGHBOR //#warning "ParFUM_recreateSharedNodes only allows adjacent partitions(rank +/- 1) to have shared nodes" sendUpperBound = rank + 1; sendLowerBound = rank - 1; #endif for (int i=rank+1; i<=sendUpperBound; i++) { //send nodeCoords to rank i MPI_Send(nodeCoords, dim*numNodes, MPI_DOUBLE, i%nParts, coord_msg_tag, comm); send_count ++; // printf("[%d] Sending %d doubles to rank %d \n",rank,dim*numNodes,i%nParts); } // Receive coordinates from the appropriate number of other partitions // These can be received in any order for (int i=sendLowerBound; i<rank; i++) { std::vector<int> remoteSharedNodes, localSharedNodes; double *recvNodeCoords; MPI_Status status; int source, length; // Probe for a coordinate message from any source; extract source and msg length MPI_Probe(MPI_ANY_SOURCE, coord_msg_tag, comm, &status); source = status.MPI_SOURCE; length = status.MPI_LENGTH/sizeof(double); // printf("[%d] Receiving %d doubles from rank %d \n",rank,length,source); recv_count ++; // Receive whatever data was available according to probe recvNodeCoords = (double *)malloc(length*sizeof(double)); MPI_Recv((void*)recvNodeCoords, length, MPI_DOUBLE, source, coord_msg_tag, comm, &status); // Match coords between local nodes and received coords int recvNodeCount = length/dim; // PERFORM THE NODE COMPARISONS #ifdef SHARED_NODES_ONLY_NEIGHBOR int borderNodes = BORDERNODES; //#warning "Only the first and last BORDERNODES nodes on each partition are candidates for being shared nodes" // indices are inclusive int myBottomLow = 0; int myBottomHigh = borderNodes; int myTopLow = numNodes - borderNodes; int myTopHigh = numNodes-1; int recvBottomLow = 0; int recvBottomHigh = borderNodes; int recvTopLow = recvNodeCount - borderNodes; int recvTopHigh = recvNodeCount-1; CkPrintf("[%d] rank=%d myBottomLow=%d myBottomHigh=%d myTopLow=%d myTopHigh=%d recvBottomLow=%d recvBottomHigh=%d recvTopLow=%d recvTopHigh=%d\n", CkMyPe(), rank, myBottomLow, myBottomHigh, myTopLow, myTopHigh, recvBottomLow, recvBottomHigh, recvTopLow, recvTopHigh); // make sure the top region is non-negative if(myTopLow < 0) myTopLow = 0; if(recvTopLow < 0) recvTopLow = 0; // make the two regions be non-overlapping if(myBottomHigh >= myTopLow) myTopLow = myTopLow-1; if(recvBottomHigh >= recvTopLow) recvTopLow = recvTopLow-1; for (int j=myBottomLow; j<=myBottomHigh; j++) { for (int k=recvBottomLow; k<=recvBottomHigh; k++) { if (coordEqual(&nodeCoords[j*dim], &recvNodeCoords[k*dim], dim)) { localSharedNodes.push_back(j); remoteSharedNodes.push_back(k); break; } } } for (int j=myTopLow; j<=myBottomHigh; j++) { for (int k=recvTopLow; k<=recvTopHigh; k++) { if (coordEqual(&nodeCoords[j*dim], &recvNodeCoords[k*dim], dim)) { localSharedNodes.push_back(j); remoteSharedNodes.push_back(k); break; } } } for (int j=myTopLow; j<=myTopHigh; j++) { for (int k=recvBottomLow; k<=recvBottomHigh; k++) { if (coordEqual(&nodeCoords[j*dim], &recvNodeCoords[k*dim], dim)) { localSharedNodes.push_back(j); remoteSharedNodes.push_back(k); break; } } } for (int j=myBottomLow; j<=myTopHigh; j++) { for (int k=recvTopLow; k<=recvTopHigh; k++) { if (coordEqual(&nodeCoords[j*dim], &recvNodeCoords[k*dim], dim)) { localSharedNodes.push_back(j); remoteSharedNodes.push_back(k); break; } } } #else // CkPrintf("Comparing %d nodes with %d received nodes\n", numNodes, recvNodeCount); for (int j=0; j<numNodes; j++) { for (int k=0; k<recvNodeCount; k++) { if (coordEqual(&nodeCoords[j*dim], &recvNodeCoords[k*dim], dim)) { localSharedNodes.push_back(j); remoteSharedNodes.push_back(k); //printf("[%d] found local node %d to match with remote node %d \n",rank,j,k); break; } } } #endif // Copy local nodes that were shared with source into the data structure int *localSharedNodeList = (int *)malloc(localSharedNodes.size()*sizeof(int)); for (int m=0; m<localSharedNodes.size(); m++) { localSharedNodeList[m] = localSharedNodes[m]; } sharedNodeCounts[source] = localSharedNodes.size(); sharedNodeLists[source] = localSharedNodeList; // do not delete localSharedNodeList as a pointer to it is stored // Send remote nodes that were shared with this partition to remote partition MPI_Send((int *)&remoteSharedNodes[0], remoteSharedNodes.size(), MPI_INT, source, sharedlist_msg_tag, comm); free(recvNodeCoords); } for (int i=rank+1; i<=sendUpperBound; i++) { // recv shared node lists (from the partitions in any order) int *sharedNodes; MPI_Status status; int source, length; // Probe for a shared node list from any source; extract source and msg length MPI_Probe(MPI_ANY_SOURCE, sharedlist_msg_tag, comm, &status); source = status.MPI_SOURCE; length = status.MPI_LENGTH/sizeof(int); // Recv the shared node list the probe revealed was available sharedNodes = (int *)malloc(length*sizeof(int)); MPI_Recv((void*)sharedNodes, length, MPI_INT, source, sharedlist_msg_tag, comm, &status); // Store the shared node list in the data structure sharedNodeCounts[source] = length; sharedNodeLists[source] = sharedNodes; // don't delete sharedNodes! we kept a pointer to it! } if (rank==0) CkPrintf("Received new shared node lists...\n"); // IMPLEMENT ME: use sharedNodeLists and sharedNodeCounts to move shared node data // to IDXL FEM_Mesh *mesh = (FEM_chunk::get("ParFUM_recreateSharedNodes"))->lookup(meshid,"ParFUM_recreateSharedNodes"); IDXL_Side &shared = mesh->node.shared; for(int i=0;i<nParts;i++){ if(i == rank) continue; if(sharedNodeCounts[i] != 0){ IDXL_List &list = shared.addList(i); for(int j=0;j<sharedNodeCounts[i];j++){ list.push_back(sharedNodeLists[i][j]); } } } MPI_Barrier(MPI_COMM_WORLD); if (rank==0) CkPrintf("Recreation of shared nodes complete...\n"); //printf("After recreating shared nodes %d \n",rank); //shared.print(); #ifdef SHARED_NODES_ONLY_NEIGHBOR CkAssert(send_count + recv_count == 2); #else CkAssert(send_count + recv_count == nParts-1); #endif // Clean up free(nodeCoords); free(sharedNodeCounts); for (int i=0; i<nParts; i++) { if (sharedNodeLists[i]) free(sharedNodeLists[i]); } free(sharedNodeLists); #endif // normal mode, not super fast mesh specific one }
extern "C" void driver(void) { int ignored; int i, count; int myChunk=FEM_My_partition(); /*Add a refinement object to FEM array*/ CkPrintf("[%d] begin init\n",myChunk); FEM_REFINE2D_Init(); CkPrintf("[%d] end init\n",myChunk); myGlobals g; FEM_Register(&g,(FEM_PupFn)pup_myGlobals); init_myGlobal(&g); g.nnodes = FEM_Mesh_get_length(FEM_Mesh_default_read(),FEM_NODE); int maxNodes = g.nnodes; g.maxnodes=2*maxNodes; g.m_i_fid=FEM_Create_field(FEM_DOUBLE,1,0,sizeof(double)); resize_nodes((void *)&g,&g.nnodes,&maxNodes); int nghost=0; g.nelems=FEM_Mesh_get_length(FEM_Mesh_default_read(),FEM_ELEM); g.maxelems=g.nelems; resize_elems((void *)&g,&g.nelems,&g.maxelems); FEM_REFINE2D_Newmesh(FEM_Mesh_default_read(),FEM_NODE,FEM_ELEM); //Initialize associated data for (i=0;i<g.maxnodes;i++) { g.R_net[i]=g.d[i]=g.v[i]=g.a[i]=vector2d(0.0); } //Apply a small initial perturbation to positions for (i=0;i<g.nnodes;i++) { const double max=1.0e-15/15.0; //Tiny perturbation g.d[i].x+=max*(i&15); g.d[i].y+=max*((i+5)&15); } int fid=FEM_Create_field(FEM_DOUBLE,2,0,sizeof(vector2d)); for (i=0;i<g.nelems;i++){ checkTriangle(g,i); } sleep(5); //Timeloop if (CkMyPe()==0){ CkPrintf("Entering timeloop\n"); } // int tSteps=0x70FF00FF; int tSteps=4; int z=13; calcMasses(g); double startTime=CkWallTimer(); double curArea=2.5e-5/1024; int t = 0; // THIS IS THE INITIAL MESH SENT TO NetFEM if (1) { //Publish data to the net publishMeshToNetFEM(g,myChunk,t); } double desiredArea; /* //should not be necessary as it would have been set in the init for (i=0; i<g.nnodes; i++) { g.validNode[i] = 1; } for (i=0; i<g.nelems; i++) { g.validElem[i] = 1; }*/ double avgArea = 0.0; for (i=0;i<g.nelems;i++) { avgArea += calcArea(g, i); } avgArea /= g.nelems; for (t=1;t<=tSteps;t++) { /* if (1) { //Structural mechanics //Compute forces on nodes exerted by elements CST_NL(g.coord,g.conn,g.R_net,g.d,matConst,g.nnodes,g.nelems,g.S11,g.S22,g.S12); //Communicate net force on shared nodes FEM_Update_field(fid,g.R_net); //Advance node positions advanceNodes(dt,g.nnodes,g.coord,g.R_net,g.a,g.v,g.d,g.m_i,(t%4)==0); }*/ //Debugging/perf. output double curTime=CkWallTimer(); double total=curTime-startTime; startTime=curTime; vector2d *loc; double *areas; // prepare to coarsen loc=new vector2d[2*g.nnodes]; for (i=0;i<g.nnodes;i++) { loc[i]=g.coord[i];//+g.d[i]; } areas=new double[g.nelems]; for (i=0;i<g.nelems;i++) { areas[i] = avgArea; } //coarsen one element at a time //int coarseIdx = (23 + 4*t)%g.nnodes; //areas[coarseIdx] = calcArea(g,coarseIdx)*2.5; CkPrintf("[%d] Starting coarsening step: %d nodes, %d elements\n", myChunk,countValidEntities(g.validNode,g.nnodes),countValidEntities(g.validElem,g.nelems)); FEM_REFINE2D_Coarsen(FEM_Mesh_default_read(),FEM_NODE,(double *)g.coord,FEM_ELEM,areas,FEM_SPARSE); repeat_after_split((void *)&g); g.nelems = FEM_Mesh_get_length(FEM_Mesh_default_read(),FEM_ELEM); g.nnodes = FEM_Mesh_get_length(FEM_Mesh_default_read(),FEM_NODE); CkPrintf("[%d] Done with coarsening step: %d nodes, %d elements\n", myChunk,countValidEntities(g.validNode,g.nnodes),countValidEntities(g.validElem,g.nelems)); delete [] loc; delete[] areas; // THIS IS THE COARSENED MESH SENT TO NetFEM if (1) { //Publish data to the net publishMeshToNetFEM(g,myChunk,2*t-1); } //prepare to refine loc=new vector2d[2*g.nnodes]; for (i=0;i<g.nnodes;i++) { loc[i]=g.coord[i];//+g.d[i]; } areas=new double[g.nelems]; for (i=0;i<g.nelems;i++) { areas[i] = avgArea; } //refine one element at a time //int refIdx = (13 + 3*t)%g.nnodes; //areas[refIdx] = calcArea(g,refIdx)/1.5; CkPrintf("[%d] Starting refinement step: %d nodes, %d elements\n", myChunk,countValidEntities(g.validNode,g.nnodes),countValidEntities(g.validElem,g.nelems)); FEM_REFINE2D_Split(FEM_Mesh_default_read(),FEM_NODE,(double *)loc,FEM_ELEM,areas,FEM_SPARSE); repeat_after_split((void *)&g); g.nelems = FEM_Mesh_get_length(FEM_Mesh_default_read(),FEM_ELEM); g.nnodes = FEM_Mesh_get_length(FEM_Mesh_default_read(),FEM_NODE); CkPrintf("[%d] Done with refinement step: %d nodes, %d elements\n", myChunk,countValidEntities(g.validNode,g.nnodes),countValidEntities(g.validElem,g.nelems)); delete [] loc; delete[] areas; // THIS IS THE REFINED MESH SENT TO NetFEM if (1) { //Publish data to the net publishMeshToNetFEM(g,myChunk,2*t); } } if (CkMyPe()==0) CkPrintf("Driver finished\n"); }
void publish_data_netfem(int i, myGlobals g) { MPI_Barrier(MPI_COMM_WORLD); if (1) { //Publish data to the net int mesh=FEM_Mesh_default_read(); // Tell framework we are reading data from the mesh int rank = FEM_My_partition(); g.nnodes=FEM_Mesh_get_length(mesh,FEM_NODE); // Get number of nodes g.coord=new vector3d[g.nnodes]; int count = 0; vector3d *coord = new vector3d[g.nnodes]; int *maptovalid = new int[g.nnodes]; double *nodeid = new double[g.nnodes]; // Get node positions FEM_Mesh_data(mesh, FEM_NODE, FEM_DATA+0, (double*)g.coord, 0, g.nnodes, FEM_DOUBLE, 3); for(int j=0; j<g.nnodes; j++) { if(FEM_is_valid(mesh,FEM_NODE+0,j)) { coord[count].x = g.coord[j].x; coord[count].y = g.coord[j].y; coord[count].z = g.coord[j].z; maptovalid[j] = count; nodeid[count] = j; count++; } } NetFEM n=NetFEM_Begin(rank,i,3,NetFEM_WRITE); NetFEM_Nodes(n,count,(double *)coord,"Position (m)"); //NetFEM_Nodes(n,g.nnodes,(double *)g.coord,"Position (m)"); NetFEM_Scalar(n,nodeid,1,"Node ID"); // Get element data g.nelems=FEM_Mesh_get_length(mesh,FEM_ELEM+0); // Get number of elements g.conn=new connRec[g.nelems]; connRec *conn = new connRec[g.nelems]; double *elid = new double[g.nelems]; count = 0; // Get connectivity for elements FEM_Mesh_data(mesh, FEM_ELEM+0, FEM_CONN, (int *)g.conn, 0, g.nelems, FEM_INDEX_0, 4); double totalVolume = 0.0; for(int j=0; j<g.nelems; j++) { if(FEM_is_valid(mesh,FEM_ELEM+0,j)) { conn[count][0] = maptovalid[g.conn[j][0]]; conn[count][1] = maptovalid[g.conn[j][1]]; conn[count][2] = maptovalid[g.conn[j][2]]; conn[count][3] = maptovalid[g.conn[j][3]]; elid[count] = j; totalVolume += getVolume(coord[conn[count][0]],coord[conn[count][1]],coord[conn[count][2]],coord[conn[count][3]]); if(totalVolume != totalVolume) { CkPrintf("NAN\n"); } count++; } } NetFEM_Elements(n,count,4,(int *)conn,"Tets"); //NetFEM_Elements(n,g.nelems,4,(int *)g.conn,"Tets"); NetFEM_Scalar(n,elid,1,"Element ID"); NetFEM_End(n); double finalVolume; CkPrintf("Chunk[%d]: local volume: %.12f\n",rank,totalVolume); MPI_Reduce((void*)&totalVolume,(void*)&finalVolume,1,MPI_DOUBLE,MPI_SUM,0,MPI_COMM_WORLD); if(rank == 0) CkPrintf("Chunk[%d]: total volume: %.12f\n",rank,finalVolume); delete [] g.coord; delete [] g.conn; delete [] coord; delete [] conn; delete [] maptovalid; delete [] nodeid; delete [] elid; } }
// A driver() function // driver() is required in all FEM programs extern "C" void driver(void) { int nnodes,nelems,nelems2,ignored; int i, myId=FEM_My_partition(); myGlobals g; FEM_Register(&g,(FEM_PupFn)pup_myGlobals); _registerFEMMeshModify(); printf("partition %d is in driver\n", myId); int mesh=FEM_Mesh_default_read(); // Tell framework we are reading data from the mesh // Get node data nnodes=FEM_Mesh_get_length(mesh,FEM_NODE); // Get number of nodes g.nnodes=nnodes; g.coord=new vector3d[nnodes]; // Get node positions FEM_Mesh_data(mesh, FEM_NODE, FEM_DATA+0, (double*)g.coord, 0, nnodes, FEM_DOUBLE, 3); // Get element data nelems=FEM_Mesh_get_length(mesh,FEM_ELEM+0); // Get number of elements g.nelems=nelems; g.conn=new connRec[nelems]; // Get connectivity for elements FEM_Mesh_data(mesh, FEM_ELEM+0, FEM_CONN, (int *)g.conn, 0, nelems, FEM_INDEX_0, 4); double* tridata =new double[nelems]; FEM_Mesh_data(mesh, FEM_ELEM+0, FEM_DATA, (int *)tridata, 0, nelems, FEM_DOUBLE, 1); int nelemsghost=FEM_Mesh_get_length(mesh,FEM_ELEM+0+FEM_GHOST); double* trighostdata =new double[nelemsghost]; FEM_Mesh_data(mesh, FEM_ELEM+0+FEM_GHOST, FEM_DATA, (int *)trighostdata, 0, nelemsghost, FEM_DOUBLE, 1); int nnodesghost=FEM_Mesh_get_length(mesh,FEM_NODE+0+FEM_GHOST); double* nodeghostdata =new double[3*nnodesghost]; FEM_Mesh_data(mesh, FEM_NODE+0+FEM_GHOST, FEM_DATA, (int *)nodeghostdata, 0, nnodesghost, FEM_DOUBLE, 3); { const int triangleFaces[12] = {0,1,2,1,2,3,2,3,0,0,3,1}; FEM_Add_elem2face_tuples(mesh, 0, 3, 4, triangleFaces); FEM_Mesh_create_elem_elem_adjacency(mesh); FEM_Mesh_allocate_valid_attr(mesh, FEM_ELEM+0); FEM_Mesh_allocate_valid_attr(mesh, FEM_NODE); FEM_Mesh_create_node_elem_adjacency(mesh); FEM_Mesh_create_node_node_adjacency(mesh); int netIndex = 0; publish_data_netfem(netIndex,g); netIndex++; #ifdef SUMMARY_ON FEM_Print_Mesh_Summary(mesh); #endif int rank = 0; MPI_Comm comm = MPI_COMM_WORLD; MPI_Comm_rank(comm,&rank); MPI_Barrier(comm); FEM_REF_INIT(mesh); FEM_Mesh *meshP = FEM_Mesh_lookup(FEM_Mesh_default_read(),"driver"); FEM_AdaptL *ada = meshP->getfmMM()->getfmAdaptL(); int ret_op = -1; FEM_Adapt_Algs *adaptAlgs= meshP->getfmMM()->getfmAdaptAlgs(); adaptAlgs->FEM_Adapt_Algs_Init(FEM_DATA+0,FEM_DATA+4); FEM_Interpolate *interp = meshP->getfmMM()->getfmInp(); //interp->FEM_SetInterpolateNodeEdgeFnPtr(interpolate); MPI_Barrier(comm); //CkPrintf("Shadow arrays have been bound\n"); /* #ifdef SUMMARY_ON FEM_Print_Mesh_Summary(mesh); #endif CkPrintf("Marking 5 nodes and one element as invalid\n"); FEM_set_entity_invalid(mesh, FEM_NODE, 5); FEM_set_entity_invalid(mesh, FEM_NODE, 6); FEM_set_entity_invalid(mesh, FEM_NODE, 7); FEM_set_entity_invalid(mesh, FEM_NODE, 8); FEM_set_entity_invalid(mesh, FEM_NODE, 9); FEM_set_entity_invalid(mesh, FEM_ELEM, 9); #ifdef SUMMARY_ON FEM_Print_Mesh_Summary(mesh); #endif CkPrintf("Marking 5 nodes and one element as valid again\n"); FEM_set_entity_valid(mesh, FEM_NODE, 5); FEM_set_entity_valid(mesh, FEM_NODE, 6); FEM_set_entity_valid(mesh, FEM_NODE, 7); FEM_set_entity_valid(mesh, FEM_NODE, 8); FEM_set_entity_valid(mesh, FEM_NODE, 9); FEM_set_entity_valid(mesh, FEM_ELEM, 9); */ #ifdef SUMMARY_ON FEM_Print_Mesh_Summary(mesh); #endif int adjs[4]; int elemid; if(rank == 0) { adjs[0] = 0; adjs[1] = 1; adjs[2] = 2; adjs[3] = 3; elemid = 0; } else if(rank == 1) { adjs[0] = 19; adjs[1] = 5; adjs[2] = 7; adjs[3] = 21; elemid = 21; } else if(rank == 2) { adjs[0] = 8; adjs[1] = 11; adjs[2] = 6; adjs[3] = 21; elemid = 7; } else { adjs[0] = 0; adjs[1] = 1; adjs[2] = 2; adjs[3] = 3; elemid = 0; } int newel1 = 0; #ifdef SUMMARY_ON FEM_Print_Mesh_Summary(mesh); FEM_Print_n2e(mesh,adjs[0]); FEM_Print_n2e(mesh,adjs[1]); FEM_Print_n2e(mesh,adjs[2]); FEM_Print_n2e(mesh,adjs[3]); FEM_Print_n2n(mesh,adjs[0]); FEM_Print_n2n(mesh,adjs[1]); FEM_Print_n2n(mesh,adjs[2]); FEM_Print_n2n(mesh,adjs[3]); FEM_Print_e2n(mesh,elemid); FEM_Print_e2e(mesh,elemid); #endif //FEM_Modify_Lock(mesh, adjs, 3, adjs, 0); if(rank == 0) { FEM_remove_element(mesh, elemid, 0, 1); } //FEM_Modify_Unlock(mesh); MPI_Barrier(comm); #ifdef SUMMARY_ON FEM_Print_Mesh_Summary(mesh); #endif //FEM_Modify_Lock(mesh, adjs, 3, adjs, 0); if(rank == 0) { newel1 = FEM_add_element(mesh,adjs,4,0,0); CkPrintf("New Element\n"); } //FEM_Modify_Unlock(mesh); publish_data_netfem(netIndex,g); netIndex++; #ifdef SUMMARY_ON FEM_Print_Mesh_Summary(mesh); FEM_Print_n2e(mesh,adjs[0]); FEM_Print_n2e(mesh,adjs[1]); FEM_Print_n2e(mesh,adjs[2]); FEM_Print_n2e(mesh,adjs[3]); FEM_Print_n2n(mesh,adjs[0]); FEM_Print_n2n(mesh,adjs[1]); FEM_Print_n2n(mesh,adjs[2]); FEM_Print_n2n(mesh,adjs[3]); FEM_Print_e2n(mesh,newel1); FEM_Print_e2e(mesh,newel1); #endif /* if(rank==0){ FEM_Print_Mesh_Summary(mesh); CkPrintf("%d: Removing element \n", rank); int nelemsghost =FEM_Mesh_get_length(mesh,FEM_ELEM+0+FEM_GHOST); int numvalidghost =FEM_count_valid(mesh,FEM_ELEM+0+FEM_GHOST); CkPrintf("nelemsghost=%d numvalidghost=%d\n", nelemsghost, numvalidghost); for(int i=1;i<20;i++){ if(FEM_is_valid(mesh, FEM_ELEM+FEM_GHOST, i)){ double data[1]; FEM_Mesh_data(mesh, FEM_ELEM+FEM_GHOST, FEM_DATA, (int *)data, i, 1, FEM_DOUBLE, 1); CkPrintf("%d: Eating ghost element %d with value %f\n", rank, i, data[1]); int conn[3]; FEM_Mesh_data(mesh, FEM_ELEM+FEM_GHOST, FEM_CONN, (int *)conn, i, 1, FEM_INDEX_0, 3); CkPrintf("conn for element is: %d %d %d\n", conn[0], conn[1], conn[2]); FEM_Modify_Lock(mesh, conn, 3, conn, 0); FEM_remove_element(mesh, FEM_From_ghost_index(i), 0, 1); FEM_Modify_Unlock(mesh); MPI_Barrier(comm); FEM_Print_Mesh_Summary(mesh); FEM_Modify_Lock(mesh, conn, 3, conn, 0); FEM_add_element(mesh, conn, 3, 0, rank); // add locally FEM_Modify_Unlock(mesh); CkPrintf("New conn for element is: %d %d %d\n", conn[0], conn[1], conn[2]); publish_data_netfem(netIndex,g); netIndex++; FEM_Print_Mesh_Summary(mesh); } else{ // CkPrintf("invalid element %d\n", i); } } } else { CkPrintf("Rank %d\n", rank); for(int i=1;i<20;i++){ MPI_Barrier(comm); FEM_Print_Mesh_Summary(mesh); publish_data_netfem(netIndex,g); netIndex++; FEM_Print_Mesh_Summary(mesh); } } publish_data_netfem(netIndex,g); netIndex++; */ /* CkPrintf("Starting Local edge flips on individual chunks\n"); int flip[4]; if(rank == 0) { flip[0] = 1; flip[1] = 2; flip[2] = 0; flip[3] = 3; } else if(rank == 1) { flip[0] = 1; flip[1] = 2; flip[2] = 0; flip[3] = 4; } else if(rank == 2) { flip[0] = 13; flip[1] = 14; flip[2] = 15; flip[3] = 7; } else { flip[0] = 0; flip[1] = 1; flip[2] = 2; flip[3] = 3; } #ifdef SUMMARY_ON FEM_Print_Mesh_Summary(mesh); #endif ret_op = ada->edge_flip(flip[0],flip[1]); publish_data_netfem(netIndex,g); netIndex++; ret_op = ada->edge_flip(flip[2],flip[3]); publish_data_netfem(netIndex,g); netIndex++; #ifdef SUMMARY_ON FEM_Print_Mesh_Summary(mesh); #endif CkPrintf("Starting shared edge flips on individual chunks\n"); int sflip[4]; if(rank == 0) { sflip[0] = 19; sflip[1] = 18; sflip[2] = 1; sflip[3] = -4; } else if(rank == 1) { sflip[0] = 5; sflip[1] = 6; sflip[2] = 7; sflip[3] = -5; } else if(rank == 2) { sflip[0] = 11; sflip[1] = 2; sflip[2] = 0; sflip[3] = -2; } else { sflip[0] = 0; sflip[1] = 1; sflip[2] = 2; sflip[3] = 3; } #ifdef SUMMARY_ON FEM_Print_Mesh_Summary(mesh); #endif ret_op = ada->edge_flip(sflip[0],sflip[1]); publish_data_netfem(netIndex,g); netIndex++; if(ret_op > 0) { if(sflip[2]<0) sflip[2] = ret_op; else if(sflip[3]<0) sflip[3] = ret_op; } ret_op = ada->edge_flip(sflip[2],sflip[3]); publish_data_netfem(netIndex,g); netIndex++; #ifdef SUMMARY_ON FEM_Print_Mesh_Summary(mesh); #endif CkPrintf("Starting Local edge bisect on individual chunks\n"); int bisect[2]; if(rank == 0) { bisect[0] = 16; bisect[1] = 21; } else if(rank == 1) { bisect[0] = 5; bisect[1] = 6; } else if(rank == 2) { bisect[0] = 8; bisect[1] = 11; } else { bisect[0] = 0; bisect[1] = 1; } if(rank==2) ret_op = ada->edge_bisect(bisect[0],bisect[1]); publish_data_netfem(netIndex,g); netIndex++; adaptAlgs->tests(); MPI_Barrier(comm); #ifdef SUMMARY_ON FEM_Print_Mesh_Summary(mesh); #endif CkPrintf("Starting Local vertex remove on individual chunks\n"); int vr[2]; if(rank == 0) { vr[0] = ret_op; vr[1] = 6; } else if(rank == 1) { vr[0] = ret_op; vr[1] = 13; } else if(rank == 2) { vr[0] = ret_op; vr[1] = 21; } else { vr[0] = ret_op; vr[1] = 1; } #ifdef SUMMARY_ON FEM_Print_Mesh_Summary(mesh); #endif ret_op = ada->vertex_remove(vr[0],vr[1]); publish_data_netfem(netIndex,g); netIndex++; #ifdef SUMMARY_ON FEM_Print_Mesh_Summary(mesh); #endif CkPrintf("Starting shared edge bisect on individual chunks\n"); int sbisect[2]; if(rank == 0) { sbisect[0] = 1;//4;//21; sbisect[1] = 19;//20; } else if(rank == 1) { sbisect[0] = 0; sbisect[1] = 21; } else if(rank == 2) { sbisect[0] = 20;//9;//3; sbisect[1] = 8;//2;//19; } else { sbisect[0] = 0; sbisect[1] = 1; } #ifdef SUMMARY_ON FEM_Print_Mesh_Summary(mesh); #endif ret_op = ada->edge_bisect(sbisect[0],sbisect[1]); publish_data_netfem(netIndex,g); netIndex++; CkPrintf("Starting shared vertex remove on individual chunks\n"); int svr[2]; if(rank == 0) { svr[0] = ret_op; svr[1] = 19; } else if(rank == 1) { svr[0] = ret_op; svr[1] = 21; } else if(rank == 2) { svr[0] = ret_op; svr[1] = 20; } else { svr[0] = ret_op; svr[1] = 1; } #ifdef SUMMARY_ON FEM_Print_Mesh_Summary(mesh); #endif ret_op = ada->vertex_remove(svr[0],svr[1]); publish_data_netfem(netIndex,g); netIndex++; #ifdef SUMMARY_ON FEM_Print_Mesh_Summary(mesh); #endif CkPrintf("Starting Local edge contract on individual chunks\n"); int contract[2]; if(rank == 0) { contract[1] = 16; contract[0] = 21; } else if(rank == 1) { contract[0] = 5; contract[1] = 6; } else if(rank == 2) { contract[0] = 19; contract[1] = 17; } else { contract[0] = 0; contract[1] = 1; } #ifdef SUMMARY_ON FEM_Print_Mesh_Summary(mesh); #endif if(rank==0) ret_op = ada->edge_contraction(contract[0],contract[1]); publish_data_netfem(netIndex,g); netIndex++; //CkPrintf("Starting Local vertex split on individual chunks\n"); /* int vs[3]; if(rank == 0) { vs[0] = ret_op; vs[1] = 9; vs[2] = 13; } else if(rank == 1) { vs[0] = ret_op; vs[1] = 8; vs[2] = 7; } else if(rank == 2) { vs[0] = ret_op; vs[1] = 14; vs[2] = 23; } else { vs[0] = ret_op; vs[1] = 2; vs[2] = 3; } #ifdef SUMMARY_ON //FEM_Print_Mesh_Summary(mesh); #endif //ret_op = ada->vertex_split(vs[0],vs[1],vs[2]); //publish_data_netfem(netIndex,g); netIndex++; #ifdef SUMMARY_ON FEM_Print_Mesh_Summary(mesh); #endif */ /* CkPrintf("Starting shared edge contract on individual chunks\n"); int scontract[2]; if(rank == 0) { scontract[0] = 9; scontract[1] = 10; } else if(rank == 1) { scontract[0] = 5; scontract[1] = 6; } else if(rank == 2) { scontract[0] = 11; scontract[1] = 2; } else { scontract[0] = 0; scontract[1] = 1; } #ifdef SUMMARY_ON FEM_Print_Mesh_Summary(mesh); #endif ret_op = ada->edge_contraction(scontract[0],scontract[1]); publish_data_netfem(netIndex,g); netIndex++; /* //CkPrintf("Starting shared vertex split on individual chunks\n"); int svs[3]; if(rank == 0) { svs[0] = ret_op; svs[1] = 1; svs[2] = -6; } else if(rank == 1) { svs[0] = ret_op; svs[1] = 7; svs[2] = 7; } else if(rank == 2) { svs[0] = ret_op; svs[1] = 0; svs[2] = -2; } else { svs[0] = ret_op; svs[1] = 2; svs[2] = 3; } #ifdef SUMMARY_ON //FEM_Print_Mesh_Summary(mesh); #endif //ret_op = ada->vertex_split(svs[0],svs[1],svs[2]); //publish_data_netfem(netIndex,g); netIndex++; #ifdef SUMMARY_ON FEM_Print_Mesh_Summary(mesh); #endif */ /* CkPrintf("Starting LEB on individual chunks\n"); int *leb_elem = (int*)malloc(1*sizeof(int)); if(rank==0) { leb_elem[0] = 2; } else if(rank==1) { leb_elem[0] = 13; //4; } else if(rank==2) { leb_elem[0] = 20; //26; } else if (rank == 3){ leb_elem[0] = 14; } else { leb_elem[0] = 0; } adaptAlgs->refine_element_leb(leb_elem[0]); publish_data_netfem(netIndex,g); netIndex++; */ /* int nEle; //for(int tstep = 0; tstep < 2; tstep++) { nEle = FEM_Mesh_get_length(mesh, FEM_ELEM); for (int i=0; i<nEle; i++) if (FEM_is_valid(mesh, FEM_ELEM, i)) adaptAlgs->refine_element_leb(i); publish_data_netfem(netIndex,g); netIndex++; FEM_Print_Mesh_Summary(mesh); //} */ double targetArea = 0.00004; for(int tstep = 0; tstep < 0; tstep++) { adaptAlgs->simple_refine(targetArea); publish_data_netfem(netIndex,g); netIndex++; adaptAlgs->tests(); MPI_Barrier(comm); //int *nodes = new int[g.nnodes]; //for (int i=0; i<g.nnodes; i++) nodes[i]=i; //FEM_mesh_smooth(mesh, nodes, g.nnodes, FEM_DATA+0); //publish_data_netfem(netIndex,g); netIndex++; //delete [] nodes; #ifdef SUMMARY_ON FEM_Print_Mesh_Summary(mesh); #endif adaptAlgs->simple_coarsen(targetArea); publish_data_netfem(netIndex,g); netIndex++; adaptAlgs->tests(); MPI_Barrier(comm); //int *nodes = new int[g.nnodes]; //for (int i=0; i<g.nnodes; i++) nodes[i]=i; //FEM_mesh_smooth(mesh, nodes, g.nnodes, FEM_DATA+0); //publish_data_netfem(netIndex,g); netIndex++; //delete [] nodes; #ifdef SUMMARY_ON FEM_Print_Mesh_Summary(mesh); #endif } targetArea = 0.00004; for(int tstep = 0; tstep < 0; tstep++) { adaptAlgs->simple_refine(targetArea); publish_data_netfem(netIndex,g); netIndex++; adaptAlgs->tests(); MPI_Barrier(comm); //int *nodes = new int[g.nnodes]; //for (int i=0; i<g.nnodes; i++) nodes[i]=i; //FEM_mesh_smooth(mesh, nodes, g.nnodes, FEM_DATA+0); //publish_data_netfem(netIndex,g); netIndex++; //delete [] nodes; targetArea *= 0.5; #ifdef SUMMARY_ON FEM_Print_Mesh_Summary(mesh); #endif } targetArea /= 0.4; for(int tstep = 0; tstep < 0; tstep++) { adaptAlgs->simple_coarsen(targetArea); publish_data_netfem(netIndex,g); netIndex++; adaptAlgs->tests(); MPI_Barrier(comm); //int *nodes = new int[g.nnodes]; //for (int i=0; i<g.nnodes; i++) nodes[i]=i; //FEM_mesh_smooth(mesh, nodes, g.nnodes, FEM_DATA+0); //publish_data_netfem(netIndex,g); netIndex++; //delete [] nodes; targetArea *= 1.5; #ifdef SUMMARY_ON FEM_Print_Mesh_Summary(mesh); #endif } //wave propagation on a bar targetArea = 0.00004; double xmin = 0.00; double xmax = 0.1; double ymin = 0.00; double ymax = 0.01; for(int tstep = 0; tstep < 0; tstep++) { targetArea = 0.000002; adaptAlgs->simple_refine(targetArea, xmin, ymin, xmax, ymax); publish_data_netfem(netIndex,g); netIndex++; #ifdef SUMMARY_ON FEM_Print_Mesh_Summary(mesh); #endif targetArea = 0.0000014; adaptAlgs->simple_coarsen(targetArea, xmin, ymin, xmax, ymax); publish_data_netfem(netIndex,g); netIndex++; #ifdef SUMMARY_ON FEM_Print_Mesh_Summary(mesh); #endif ymin += 0.01; ymax += 0.01; } //crack propagation on a block targetArea = 0.00004; xmin = 0.00; xmax = 0.2; double xcrackmin = 0.09; double xcrackmax = 0.10; ymin = 0.00; ymax = 0.02; for(int tstep = 0; tstep < 0; tstep++) { targetArea = 0.000025; adaptAlgs->simple_refine(targetArea, xmin, ymin, xmax, ymax); publish_data_netfem(netIndex,g); netIndex++; #ifdef SUMMARY_ON FEM_Print_Mesh_Summary(mesh); #endif targetArea = 0.00005; adaptAlgs->simple_coarsen(targetArea, xmin, ymin, xmax, ymax); //publish_data_netfem(netIndex,g); netIndex++; #ifdef SUMMARY_ON FEM_Print_Mesh_Summary(mesh); #endif /*if(tstep > 2) { targetArea = 0.000025; adaptAlgs->simple_refine(targetArea, xcrackmin, ymin, xcrackmax, ymax); //publish_data_netfem(netIndex,g); netIndex++; #ifdef SUMMARY_ON FEM_Print_Mesh_Summary(mesh); #endif xcrackmin -= 0.004; xcrackmax += 0.004; } */ ymin += 0.02; ymax += 0.02; } CkPrintf("chunk %d Waiting for Synchronization\n",rank); MPI_Barrier(comm); CkPrintf("Synchronized\n"); #ifdef SUMMARY_ON FEM_Print_Mesh_Summary(mesh); #endif publish_data_netfem(netIndex,g); netIndex++; CkExit(); } }
void FEM_mesh_smooth(int mesh, int *nodes, int nNodes, int attrNo) { vector3d *centroids, newPos, *coords, *ghostCoords, *vGcoords; int nEle, nGn, *boundVals, nodesInChunk, nVg; int neighbors[3], *adjelems; int gIdxN; int j=0; double x[3], y[3]; FEM_Mesh *meshP = FEM_Mesh_lookup(mesh, "driver"); nodesInChunk = FEM_Mesh_get_length(mesh,FEM_NODE); boundVals = new int[nodesInChunk]; nGn = FEM_Mesh_get_length(mesh, FEM_GHOST + FEM_NODE); coords = new vector3d[nodesInChunk+nGn]; FEM_Mesh_data(mesh, FEM_NODE, FEM_BOUNDARY, (int*) boundVals, 0, nodesInChunk, FEM_INT, 1); FEM_Mesh_data(mesh, FEM_NODE, attrNo, (double*)coords, 0, nodesInChunk, FEM_DOUBLE, 2); for (int i=0; i<(nodesInChunk); i++) { //CkPrintf(" coords[%d]: (%f, %f)\n", i, coords[i].x, coords[i].y); } IDXL_Layout_t coord_layout = IDXL_Layout_create(IDXL_DOUBLE, 2); FEM_Update_ghost_field(coord_layout,-1, coords); ghostCoords = &(coords[nodesInChunk]); /* for (int i=0; i<nGn;i++) { if (FEM_is_valid(mesh, FEM_GHOST+FEM_NODE, i)) { CkPrintf("ghost %d is valid \n", i); // vGcoords[j]=ghostCoords[i]; //j++; } else CkPrintf("ghost %d is invalid \n", i); } */ for (int i=0; i<(nodesInChunk+nGn); i++) { //CkPrintf(" coords[%d]: (%f, %f)\n", i, coords[i].x, coords[i].y); } // FEM_Mesh_data(FEM_Mesh_default_write(), FEM_GHOST+FEM_NODE, attrNo, (double*)ghostCoords, 0, nGn, FEM_DOUBLE, 2); for (int i=0; i<nNodes; i++) { newPos.x=0; newPos.y=0; CkAssert(nodes[i]<nodesInChunk); if (FEM_is_valid(mesh, FEM_NODE, i) && boundVals[i]>-1) //node must be internal { meshP->n2e_getAll(i, &adjelems, &nEle); centroids = new vector3d[nEle]; for (int j=0; j<nEle; j++) { //for all adjacent elements, find centroids meshP->e2n_getAll(adjelems[j], neighbors); for (int k=0; k<3; k++) { if (neighbors[k]<-1) { gIdxN = FEM_From_ghost_index(neighbors[k]); x[k] = ghostCoords[gIdxN].x; y[k] = ghostCoords[gIdxN].y; } else { x[k] = coords[neighbors[k]].x; y[k] = coords[neighbors[k]].y; } } centroids[j].x=(x[0]+x[1]+x[2])/3.0; centroids[j].y=(y[0]+y[1]+y[2])/3.0; newPos.x += centroids[j].x; newPos.y += centroids[j].y; } newPos.x/=nEle; newPos.y/=nEle; FEM_set_entity_coord2(mesh, FEM_NODE, nodes[i], newPos.x, newPos.y); delete [] centroids; delete [] adjelems; } } delete [] coords; delete [] boundVals; }
extern "C" void driver(void) { int ignored; int i; int myChunk=FEM_My_partition(); /*Add a refinement object to FEM array*/ CkPrintf("[%d] begin init\n",myChunk); FEM_REFINE2D_Init(); CkPrintf("[%d] end init\n",myChunk); myGlobals g; FEM_Register(&g,(FEM_PupFn)pup_myGlobals); init_myGlobal(&g); g.nnodes = FEM_Mesh_get_length(FEM_Mesh_default_read(),FEM_NODE); int maxNodes = g.nnodes; g.maxnodes=2*maxNodes; g.m_i_fid=FEM_Create_field(FEM_DOUBLE,1,0,sizeof(double)); resize_nodes((void *)&g,&g.nnodes,&maxNodes); int nghost=0; g.nelems=FEM_Mesh_get_length(FEM_Mesh_default_read(),FEM_ELEM); g.maxelems=g.nelems; resize_elems((void *)&g,&g.nelems,&g.maxelems); g.nedges = FEM_Mesh_get_length(FEM_Mesh_default_read(),FEM_SPARSE); g.maxedges = g.nedges; resize_edges((void *)&g,&g.nedges,&g.maxedges); FEM_REFINE2D_Newmesh(FEM_Mesh_default_read(),FEM_NODE,FEM_ELEM); //Initialize associated data for (i=0;i<g.maxnodes;i++){ g.R_net[i]=g.d[i]=g.v[i]=g.a[i]=vector2d(0.0); } //Apply a small initial perturbation to positions for (i=0;i<g.nnodes;i++) { const double max=1.0e-15/15.0; //Tiny perturbation g.d[i].x+=max*(i&15); g.d[i].y+=max*((i+5)&15); } int fid=FEM_Create_field(FEM_DOUBLE,2,0,sizeof(vector2d)); for (i=0;i<g.nelems;i++){ checkTriangle(g,i); } //Timeloop if (CkMyPe()==0){ CkPrintf("Entering timeloop\n"); } int tSteps=0x70FF00FF; calcMasses(g); double startTime=CkWallTimer(); double curArea=2.0e-5; for (int t=0;t<tSteps;t++) { if (1) { //Structural mechanics //Compute forces on nodes exerted by elements CST_NL(g.coord,g.conn,g.R_net,g.d,matConst,g.nnodes,g.nelems,g.S11,g.S22,g.S12); //Communicate net force on shared nodes FEM_Update_field(fid,g.R_net); //Advance node positions advanceNodes(dt,g.nnodes,g.coord,g.R_net,g.a,g.v,g.d,g.m_i,(t%4)==0); } //Debugging/perf. output double curTime=CkWallTimer(); double total=curTime-startTime; startTime=curTime; if (CkMyPe()==0 && (t%64==0)) CkPrintf("%d %.6f sec for loop %d \n",CkNumPes(),total,t); /* if (0 && t%16==0) { CkPrintf(" Triangle 0:\n"); for (int j=0;j<3;j++) { int n=g.conn[0][j]; CkPrintf(" Node %d: coord=(%.4f,%.4f) d=(%.4g,%.4g)\n", n,g.coord[n].x,g.coord[n].y,g.d[n].x,g.d[n].y); } }*/ // if (t%512==0) // FEM_Migrate(); if (t%128==0) { //Refinement: vector2d *loc=new vector2d[2*g.nnodes]; for (i=0;i<g.nnodes;i++) { loc[i]=g.coord[i];//+g.d[i]; } double *areas=new double[g.nelems]; curArea=curArea*0.98; for (i=0;i<g.nelems;i++) { #if 0 double origArea=8e-8; //Typical triangle size if (fabs(g.S12[i])>1.0e8) areas[i]=origArea*0.9; //Refine stuff that's stressed else areas[i]=origArea; //Leave everything else big #endif areas[i]=curArea; } CkPrintf("[%d] Starting refinement step: %d nodes, %d elements to %.3g\n", myChunk,g.nnodes,g.nelems,curArea); FEM_REFINE2D_Split(FEM_Mesh_default_read(),FEM_NODE,(double *)loc,FEM_ELEM,areas,FEM_SPARSE); repeat_after_split((void *)&g); g.nelems = FEM_Mesh_get_length(FEM_Mesh_default_read(),FEM_ELEM); g.nnodes = FEM_Mesh_get_length(FEM_Mesh_default_read(),FEM_NODE); CkPrintf("[%d] Done with refinement step: %d nodes, %d elements\n", myChunk,g.nnodes,g.nelems); } if (1) { //Publish data to the net NetFEM n=NetFEM_Begin(myChunk,t,2,NetFEM_POINTAT); NetFEM_Nodes(n,g.nnodes,(double *)g.coord,"Position (m)"); NetFEM_Vector(n,(double *)g.d,"Displacement (m)"); NetFEM_Vector(n,(double *)g.v,"Velocity (m/s)"); NetFEM_Elements(n,g.nelems,3,(int *)g.conn,"Triangles"); NetFEM_Scalar(n,g.S11,1,"X Stress (pure)"); NetFEM_Scalar(n,g.S22,1,"Y Stress (pure)"); NetFEM_Scalar(n,g.S12,1,"Shear Stress (pure)"); NetFEM_End(n); } } if (CkMyPe()==0) CkPrintf("Driver finished\n"); }
void FEM_REFINE2D_Split(int meshID,int nodeID,double *coord,int elemID,double *desiredAreas,int sparseID){ int nnodes = FEM_Mesh_get_length(meshID,nodeID); int nelems = FEM_Mesh_get_length(meshID,elemID); int actual_nodes = nnodes, actual_elems = nelems; FEM_Refine_Operation_Data refine_data; refine_data.meshID = meshID; refine_data.nodeID = nodeID; refine_data.sparseID = sparseID; refine_data.elemID = elemID; refine_data.cur_nodes = FEM_Mesh_get_length(meshID,nodeID); /*Copy the cordinates of the nodes into a vector, the cordinates of the new nodes will be inserted into this vector and will be used to sort all the nodes on the basis of the distance from origin */ CkVec<double> coordVec; for(int i=0;i<nnodes*2;i++){ coordVec.push_back(coord[i]); } refine_data.coordVec = &coordVec; refine_data.coord = coord; /*find out the attributes of the node */ FEM_Entity *e=refine_data.node = FEM_Entity_lookup(meshID,nodeID,"REFINE2D_Mesh"); CkVec<FEM_Attribute *> *attrs = refine_data.attrs = e->getAttrVec(); /* FEM_DataAttribute *boundaryAttr = (FEM_DataAttribute *)e->lookup(FEM_BOUNDARY,"split"); if(boundaryAttr != NULL){ AllocTable2d<int> &boundaryTable = boundaryAttr->getInt(); printf(" Node Boundary flags \n"); for(int i=0;i<nnodes;i++){ printf("Node %d flag %d \n",i,boundaryTable[i][0]); } } */ FEM_Entity *elem = refine_data.elem = FEM_Entity_lookup(meshID,elemID,"REFIN2D_Mesh_elem"); CkVec<FEM_Attribute *> *elemattrs = refine_data.elemattrs = elem->getAttrVec(); FEM_Attribute *connAttr = elem->lookup(FEM_CONN,"REFINE2D_Mesh"); if(connAttr == NULL){ CkAbort("Grrrr element without connectivity \n"); } AllocTable2d<int> &connTable = ((FEM_IndexAttribute *)connAttr)->get(); refine_data.connTable = &connTable; //hashtable to store the new node number as a function of the two old numbers CkHashtableT<intdual,int> newnodes(nnodes); refine_data.newnodes = &newnodes; /* Get the FEM_BOUNDARY data of sparse elements and load it into a hashtable indexed by the 2 node ids that make up the edge. The data in the hashtable is the index number of the sparse element */ FEM_Entity *sparse; CkVec<FEM_Attribute *> *sparseattrs; FEM_Attribute *sparseConnAttr, *sparseBoundaryAttr; AllocTable2d<int> *sparseConnTable, *sparseBoundaryTable; CkHashtableT<intdual,int> nodes2sparse(nelems); refine_data.nodes2sparse = &nodes2sparse; if(sparseID != -1){ sparse = refine_data.sparse = FEM_Entity_lookup(meshID,sparseID,"REFINE2D_Mesh_sparse"); refine_data.sparseattrs = sparseattrs = sparse->getAttrVec(); refine_data.sparseConnAttr = sparseConnAttr = sparse->lookup(FEM_CONN,"REFINE2D_Mesh_sparse"); sparseConnTable = &(((FEM_IndexAttribute *)sparseConnAttr)->get()); refine_data.sparseBoundaryAttr = sparseBoundaryAttr = sparse->lookup(FEM_BOUNDARY,"REFINE2D_Mesh_sparse"); if(sparseBoundaryAttr == NULL){ CkAbort("Specified sparse elements without boundary conditions"); } FEM_DataAttribute *validEdgeAttribute = (FEM_DataAttribute *)sparse->lookup(FEM_VALID,"REFINE2D_Mesh_sparse"); if(validEdgeAttribute){ refine_data.validEdge = &(validEdgeAttribute->getInt()); }else{ refine_data.validEdge = NULL; } /* since the default value in the hashtable is 0, to distinguish between uninserted keys and the sparse element with index 0, the index of the sparse elements is incremented by 1 while inserting. */ // printf("[%d] Sparse elements\n",FEM_My_partition()); for(int j=0;j<sparse->size();j++){ if(refine_data.validEdge == NULL || (*(refine_data.validEdge))[j][0]){ int *cdata = (*sparseConnTable)[j]; // printf("%d < %d,%d > \n",j,cdata[0],cdata[1]); nodes2sparse.put(intdual(cdata[0],cdata[1])) = j+1; } } }else{ printf("Edge boundary conditions not passed into FEM_REFINE2D_Split \n"); } //count the actual number of nodes and elements if(refine_data.node->lookup(FEM_VALID,"refine2D_splilt") != NULL){ AllocTable2d<int> &validNodeTable = ((FEM_DataAttribute *)(refine_data.node->lookup(FEM_VALID,"refine2D_splilt")))->getInt(); actual_nodes = countValidEntities(validNodeTable.getData(),nnodes); } if(refine_data.elem->lookup(FEM_VALID,"refine2D_splilt") != NULL){ AllocTable2d<int> &validElemTable = ((FEM_DataAttribute *)(refine_data.elem->lookup(FEM_VALID,"refine2D_splilt")))->getInt(); actual_elems = countValidEntities(validElemTable.getData(),nelems); } DEBUGINT(printf("%d %d \n",nnodes,nelems)); REFINE2D_Split(actual_nodes,coord,actual_elems,desiredAreas,&refine_data); int nSplits= refine_data.nSplits = REFINE2D_Get_Split_Length(); DEBUGINT(printf("called REFINE2D_Split nSplits = %d \n",nSplits)); if(nSplits == 0){ return; } for(int split = 0;split < nSplits;split++){ refineData op; REFINE2D_Get_Split(split,&op); FEM_Refine_Operation(&refine_data,op); } DEBUGINT(printf("Cordinate list length %d according to FEM %d\n",coordVec.size()/2,FEM_Mesh_get_length(meshID,nodeID))); IDXL_Sort_2d(FEM_Comm_shared(meshID,nodeID),coordVec.getVec()); int read = FEM_Mesh_is_get(meshID) ; assert(read); }
void FEM_REFINE2D_Coarsen(int meshID,int nodeID,double *coord,int elemID,double *desiredAreas,int sparseID){ int nnodes = FEM_Mesh_get_length(meshID,nodeID); int nelems = FEM_Mesh_get_length(meshID,elemID); int nodeCount=0,elemCount=0; /* The attributes of the different entities */ FEM_Entity *node=FEM_Entity_lookup(meshID,nodeID,"REFINE2D_Mesh"); CkVec<FEM_Attribute *> *attrs = node->getAttrVec(); FEM_Attribute *validNodeAttr = node->lookup(FEM_VALID,"FEM_COARSEN"); FEM_Attribute *nodeBoundaryAttr = node->lookup(FEM_BOUNDARY,"FEM_COARSEN"); if(!nodeBoundaryAttr){ printf("Warning:- no boundary flags for nodes \n"); } FEM_Entity *elem = FEM_Entity_lookup(meshID,elemID,"REFIN2D_Mesh_elem"); CkVec<FEM_Attribute *> *elemattrs = elem->getAttrVec(); FEM_Attribute *validElemAttr = elem->lookup(FEM_VALID,"FEM_COARSEN"); FEM_Attribute *connAttr = elem->lookup(FEM_CONN,"REFINE2D_Mesh"); if(connAttr == NULL){ CkAbort("Grrrr element without connectivity \n"); } AllocTable2d<int> &connTable = ((FEM_IndexAttribute *)connAttr)->get(); AllocTable2d<int>&validNodeTable = ((FEM_DataAttribute *)validNodeAttr)->getInt(); AllocTable2d<int>&validElemTable = ((FEM_DataAttribute *)validElemAttr)->getInt(); FEM_Operation_Data *coarsen_data = new FEM_Operation_Data; coarsen_data->node = (FEM_Node *)node; coarsen_data->coord = coord; int *connData = coarsen_data->connData= connTable.getData(); int *validNodeData = coarsen_data->validNodeData = validNodeTable.getData(); int *validElemData = coarsen_data->validElemData = validElemTable.getData(); /* Extract the data for node boundaries */ int *nodeBoundaryData= coarsen_data->nodeBoundaryData =NULL; if(nodeBoundaryAttr){ AllocTable2d<int> &nodeBoundaryTable = ((FEM_DataAttribute *)nodeBoundaryAttr)->getInt(); nodeBoundaryData = coarsen_data->nodeBoundaryData = nodeBoundaryTable.getData(); } for(int k=0;k<nnodes;k++){ int valid = validNodeData[k]; if(validNodeData[k]){ nodeCount++; } } for(int k=0;k<nelems;k++){ if(validElemData[k]){ elemCount++; } } /* populate the hashtable from nodes2edges with the valid edges */ CkHashtableT<intdual,int> nodes2sparse(nelems); coarsen_data->nodes2sparse = &nodes2sparse; FEM_Attribute *sparseBoundaryAttr; if(sparseID != -1){ FEM_Entity *sparse = FEM_Entity_lookup(meshID,sparseID,"Coarsen_sparse"); FEM_DataAttribute *validEdgeAttribute = (FEM_DataAttribute *)sparse->lookup(FEM_VALID,"Coarsen_sparse"); FEM_IndexAttribute *sparseConnAttribute = (FEM_IndexAttribute *)sparse->lookup(FEM_CONN,"Coarsen_sparse"); AllocTable2d<int> *sparseConnTable = coarsen_data->sparseConnTable = &(sparseConnAttribute->get()); coarsen_data->sparseBoundaryAttr = sparseBoundaryAttr = sparse->lookup(FEM_BOUNDARY,"REFINE2D_Mesh_sparse"); if(sparseBoundaryAttr == NULL){ CkAbort("Specified sparse elements without boundary conditions"); } if(validEdgeAttribute){ coarsen_data->validEdge = &(validEdgeAttribute->getInt()); }else{ coarsen_data->validEdge = NULL; } /* since the default value in the hashtable is 0, to distinguish between uninserted keys and the sparse element with index 0, the index of the sparse elements is incremented by 1 while inserting. */ printf("[%d] Sparse elements\n",FEM_My_partition()); for(int j=0;j<sparse->size();j++){ if(coarsen_data->validEdge == NULL || (*(coarsen_data->validEdge))[j][0]){ int *cdata = (*sparseConnTable)[j]; printf("%d < %d,%d > \n",j,cdata[0],cdata[1]); nodes2sparse.put(intdual(cdata[0],cdata[1])) = j+1; } } }else{ coarsen_data->validEdge = NULL; } DEBUGINT(printf("coarsen %d %d \n",nodeCount,elemCount)); REFINE2D_Coarsen(nodeCount,coord,elemCount,desiredAreas,coarsen_data); int nCollapses = REFINE2D_Get_Collapse_Length(); for(int i=0;i<nnodes;i++){ int valid = validNodeData[i]; } delete coarsen_data; }
void FEM_Refine_Operation(FEM_Refine_Operation_Data *data,refineData &op){ int meshID = data->meshID; int nodeID = data->nodeID; int sparseID = data->sparseID; CkVec<FEM_Attribute *> *attrs = data->attrs; CkVec<FEM_Attribute *> *elemattrs = data->elemattrs; CkHashtableT<intdual,int> *newnodes=data->newnodes; CkHashtableT<intdual,int> *nodes2sparse = data->nodes2sparse; FEM_Attribute *sparseConnAttr = data->sparseConnAttr, *sparseBoundaryAttr = data->sparseBoundaryAttr; AllocTable2d<int> *connTable = data->connTable; double *coord = data->coord; AllocTable2d<int> *sparseConnTable, *sparseBoundaryTable; CkVec<FEM_Attribute *> *sparseattrs = data->sparseattrs; /* FEM_DataAttribute *boundaryAttr = (FEM_DataAttribute *)data->node->lookup(FEM_BOUNDARY,"split"); if(boundaryAttr != NULL){ AllocTable2d<int> &boundaryTable = boundaryAttr->getInt(); printf(" Node Boundary flags \n"); for(int i=0;i<data->node->size();i++){ printf("Node %d flag %d \n",i,boundaryTable[i][0]); } } */ int tri=op.tri,A=op.A,B=op.B,C=op.C,D=op.D; double frac=op.frac; // current number of nodes in the mesh int *connData = connTable->getData(); int flags=op.flag; if((flags & 0x1) || (flags & 0x2)){ //new node DEBUGINT(CkPrintf("---- Adding node %d\n",D)); /* lastA=A; lastB=B; lastD=D;*/ if (A>=data->cur_nodes) CkAbort("Calculated A is invalid!"); if (B>=data->cur_nodes) CkAbort("Calculated B is invalid!"); if(D >= data->cur_nodes){ data->node->setLength(D+1); data->cur_nodes = D+1; } for(int i=0;i<attrs->size();i++){ FEM_Attribute *a = (FEM_Attribute *)(*attrs)[i]; if(a->getAttr()<FEM_ATTRIB_TAG_MAX){ FEM_DataAttribute *d = (FEM_DataAttribute *)a; d->interpolate(A,B,D,frac); }else{ /* Boundary value of new node D should be boundary value of the edge (sparse element) that contains the two nodes A & B. */ if(a->getAttr() == FEM_BOUNDARY){ if(sparseID != -1){ int sidx = nodes2sparse->get(intdual(A,B))-1; if(sidx == -1){ CkAbort("no sparse element between these 2 nodes, are they really connected ??"); } sparseBoundaryTable = &(((FEM_DataAttribute *)sparseBoundaryAttr)->getInt()); int boundaryVal = ((*sparseBoundaryTable)[sidx])[0]; (((FEM_DataAttribute *)a)->getInt()[D])[0] = boundaryVal; }else{ /* if sparse elements don't exist, just do simple interpolation */ FEM_DataAttribute *d = (FEM_DataAttribute *)a; d->interpolate(A,B,D,frac); } } } } int AandB[2]; AandB[0]=A; AandB[1]=B; /* Add a new node D between A and B */ IDXL_Add_entity(FEM_Comm_shared(meshID,nodeID),D,2,AandB); double Dx = coord[2*A]*(1-frac)+frac*coord[2*B]; double Dy = coord[2*A+1]*(1-frac)+frac*coord[2*B+1]; data->coordVec->push_back(Dx); data->coordVec->push_back(Dy); newnodes->put(intdual(A,B))=D; /* add the new sparse element <D,B> and modify the connectivity of the old one from <A,B> to <A,D> and change the hashtable to reflect that change */ if(sparseID != -1){ int oldsidx = nodes2sparse->get(intdual(A,B))-1; int newsidx = data->sparse->size(); data->sparse->setLength(newsidx+1); for(int satt = 0;satt<sparseattrs->size();satt++){ if((*sparseattrs)[satt]->getAttr() == FEM_CONN){ /* change the conn of the old sparse to A,D and new one to B,D */ sparseConnTable = &(((FEM_IndexAttribute *)sparseConnAttr)->get()); int *oldconn = (*sparseConnTable)[oldsidx]; int *newconn = (*sparseConnTable)[newsidx]; oldconn[0] = A; oldconn[1] = D; newconn[0] = D; newconn[1] = B; //printf("<%d,%d> edge being split into <%d,%d> <%d,%d> \n",A,B,A,D,D,B); }else{ /* apart from conn copy everything else */ FEM_Attribute *attr = (FEM_Attribute *)(*sparseattrs)[satt]; attr->copyEntity(newsidx,*attr,oldsidx); } } /* modify the hashtable - delete the old edge and the new ones */ nodes2sparse->remove(intdual(A,B)); nodes2sparse->put(intdual(A,D)) = oldsidx+1; nodes2sparse->put(intdual(D,B)) = newsidx+1; } } //add a new triangle /*TODO: replace FEM_ELEM with parameter*/ int newTri = op._new; int cur_elements = FEM_Mesh_get_length(data->meshID,data->elemID); DEBUGINT(CkPrintf("---- Adding triangle %d after splitting %d \n",newTri,tri)); if(newTri >= cur_elements){ data->elem->setLength(newTri+1); } D = newnodes->get(intdual(A,B)); for(int j=0;j<elemattrs->size();j++){ if((*elemattrs)[j]->getAttr() == FEM_CONN){ DEBUGINT(CkPrintf("elem attr conn code %d \n",(*elemattrs)[j]->getAttr())); //it is a connectivity attribute.. get the connectivity right FEM_IndexAttribute *connAttr = (FEM_IndexAttribute *)(*elemattrs)[j]; AllocTable2d<int> &table = connAttr->get(); DEBUGINT(CkPrintf("Table of connectivity attribute starts at %p width %d \n",table[0],connAttr->getWidth())); int *oldRow = table[tri]; int *newRow = table[newTri]; for (int i=0;i<3;i++){ if (oldRow[i]==A){ oldRow[i]=D; DEBUGINT(CkPrintf("In triangle %d %d replaced by %d \n",tri,A,D)); } } for (int i=0; i<3; i++) { if (oldRow[i] == B){ newRow[i] = D; } else if (oldRow[i] == C){ newRow[i] = C; } else if (oldRow[i] == D){ newRow[i] = A; } } DEBUGINT(CkPrintf("New Triangle %d (%d %d %d) conn %p\n",newTri,newRow[0],newRow[1],newRow[2],newRow)); }else{ FEM_Attribute *elattr = (FEM_Attribute *)(*elemattrs)[j]; if(elattr->getAttr() < FEM_ATTRIB_FIRST){ elattr->copyEntity(newTri,*elattr,tri); } } } if(sparseID != -1){ /* add sparse element (edge between C and D) */ int cdidx = data->sparse->size(); data->sparse->setLength(cdidx+1); for(int satt = 0; satt < sparseattrs->size();satt++){ if((*sparseattrs)[satt]->getAttr() == FEM_CONN){ sparseConnTable = &(((FEM_IndexAttribute *)sparseConnAttr)->get()); int *cdconn = (*sparseConnTable)[cdidx]; cdconn[0]=C; cdconn[1]=D; } if((*sparseattrs)[satt]->getAttr() == FEM_BOUNDARY){ /* An edge connecting C and D has to be an internal edge */ sparseBoundaryTable = &(((FEM_DataAttribute *)sparseBoundaryAttr)->getInt()); ((*sparseBoundaryTable)[cdidx])[0] = 0; } } if(data->validEdge){ (*(data->validEdge))[cdidx][0] = 1; } nodes2sparse->put(intdual(C,D)) = cdidx+1; } }