PPL Logo

libs/ck-libs/TMRC2D/refine.h File Reference

Go to the source code of this file.

Data Structures

struct  refineData
struct  collapseData
struct  updateData
struct  replacedeleteData
struct  coarsenData

Enumerations

enum  { INVALID = 0, COLLAPSE, UPDATE, REPLACE }

Functions

void REFINE2D_Init (void)
 Create a new refinement object for this virtual processor.
void REFINE2D_NewMesh (int meshID, int nEl, int nGhost, int nnodes, const int *conn, const int *gid, const int *boundaries, const int *edgeBounds, const int *edgeConn, int nEdges)
 Push a new mesh into the refinement system.
void REFINE2D_Split (int nNode, double *coord, int nEl, double *desiredArea, FEM_Refine_Operation_Data *refine_data)
int REFINE2D_Get_Split_Length (void)
 Get the number of split triangles.
void REFINE2D_Get_Split (int splitNo, refineData *d)
 Return one split triangle.
void REFINE2D_Coarsen (int nNode, double *coord, int nEl, double *desiredArea, FEM_Operation_Data *data)
void REFINE2D_Get_Collapse (int i, coarsenData *output)
int REFINE2D_Get_Collapse_Length ()
void REFINE2D_Check (int nEle, const int *conn, int nNodes)
 Check to make sure our connectivity and the refine connectivity agree.


Enumeration Type Documentation

anonymous enum

Enumerator:
INVALID 
COLLAPSE 
UPDATE 
REPLACE 

Definition at line 142 of file refine.h.


Function Documentation

void REFINE2D_Init ( void   ) 

Create a new refinement object for this virtual processor.

Must be called exactly once at startup.

Definition at line 18 of file refine.C.

References CkArrayOptions::bindTo(), comm, TCharm::get(), TCharm::getNumElements(), TCharm::getProxy(), mesh, MPI_Bcast, MPI_Comm_rank, chunkMsg::myThreads, chunkMsg::nChunks, rank, and TCharm::suspend().

Referenced by FEM_REFINE2D_Init(), and FTN_NAME().

Here is the call graph for this function:

Here is the caller graph for this function:

void REFINE2D_NewMesh ( int  meshID,
int  nEl,
int  nGhost,
int  nnodes,
const int conn,
const int gid,
const int boundaries,
const int edgeBounds,
const int edgeConn,
int  nEdges 
)

Push a new mesh into the refinement system.

This is the first call user programs make into the refinement system. This call need *not* be repeated when the refinement system changes the mesh; only when the user changes the mesh (e.g., to coarsen it).

Conn is row-major, and maps an element number to three node numbers. Hence conn[i*3+j] gives the local number of node j of element i. Because of shared nodes, every node of every local element will be local. Ghost elements may not have the complete set of nodes-- some of their nodes may have the invalid number -1.

Elements with numbers between 0 and nEl-1 (inclusive) are local. Elements with numbers between nEl and nGhost-1 (inclusive) are "ghosts"-- elements that are actually local on another processor. There are guaranteed to be enough ghosts that every local element's non-boundary edge will face a ghost element.

gid maps an element number to a chunk number and local number on that chunk. These are stored at gid[i*2+0] (chunk number) and gid[i*2+1] (local number).

boundaries is the boundary flags for nodes. boundaries may be null, if the user doesnt specify boundary flags. nnodes specifies the number of nodes

Definition at line 46 of file refine.C.

References CkWaitQD(), and MPI_Barrier.

Here is the call graph for this function:

void REFINE2D_Split ( int  nNode,
double *  coord,
int  nEl,
double *  desiredArea,
FEM_Refine_Operation_Data refine_data 
)

Definition at line 244 of file refine.C.

References C, chunk::multipleRefine(), chunk::refineResultsStorage, and chunk::updateNodeCoords().

Here is the call graph for this function:

int REFINE2D_Get_Split_Length ( void   ) 

Get the number of split triangles.

Definition at line 176 of file refine.C.

References refineResults::countResults(), and getResults().

Referenced by FEM_REFINE2D_Split(), and FTN_NAME().

Here is the call graph for this function:

Here is the caller graph for this function:

void REFINE2D_Get_Split ( int  splitNo,
refineData d 
)

Return one split triangle.

A and B are the nodes along the splitting edge:

C C / \ /|\ / \ / | \ / \ => / | \ / tri \ / | \ / \ /tri | new\ B --------- A B --- D --- A

The original triangle's node A should be replaced by D; while a new triangle should be inserted with nodes CAD.

The new node D's location should equal A*(1-frac)+B*frac. For a simple splitter, frac will always be 0.5.

If nodes A and B are shared with some other processor, that processor will also receive a "split" call for the same edge. If nodes A and B are shared by some other local triangle, that triangle will immediately receive a "split" call for the same edge.

Parameters: -splitNo is the number of this split. Pass splitNo in increasing order from 0 to REFINE2D_Get_Split_Length()-1. -conn is the triangle connectivity array, used to look up the node numbers A, B, and C. This array is *not* modified. -tri returns the number of the old triangle being split; as labelled above. -A,B,C returns the numbers of the nodes in the above diagram. -frac returns the weighting for D between A and B; for now, this is always 0.5.

Client's responsibilities: -Add the new node D. Since both sides of a shared local edge will receive a "split" call, you must ensure the node is not added twice; you can do this by checking this split's nodes A and B against the previous split. -Update connectivity for source triangle -Add new triangle.

Definition at line 299 of file refine.C.

References refineResults::extract(), getResults(), and r.

Here is the call graph for this function:

void REFINE2D_Coarsen ( int  nNode,
double *  coord,
int  nEl,
double *  desiredArea,
FEM_Operation_Data data 
)

Definition at line 257 of file refine.C.

References C, chunk::coarsenResultsStorage, chunk::multipleCoarsen(), and chunk::updateNodeCoords().

Referenced by FEM_REFINE2D_Coarsen().

Here is the call graph for this function:

Here is the caller graph for this function:

void REFINE2D_Get_Collapse ( int  i,
coarsenData output 
)

Definition at line 332 of file refine.C.

References coarsenResults::extract(), and getCoarsenResults().

Here is the call graph for this function:

int REFINE2D_Get_Collapse_Length (  ) 

Definition at line 324 of file refine.C.

References coarsenResults::countResults(), and getCoarsenResults().

Referenced by FEM_REFINE2D_Coarsen().

Here is the call graph for this function:

Here is the caller graph for this function:

void REFINE2D_Check ( int  nEle,
const int conn,
int  nNodes 
)

Check to make sure our connectivity and the refine connectivity agree.

Definition at line 256 of file refine.C.

References checkConn().

Here is the call graph for this function:


Generated on Mon Sep 21 08:11:26 2020 for Charm++ by  doxygen 1.5.5