00001
00002
00003
00004
00005
00006
00007
00008 #ifndef __ParFUM_Adapt_Algs_H
00009 #define __ParFUM_Adapt_Algs_H
00010
00011 #define SLIVERAREA 1.0e-18
00012 #define REFINE_TOL 1.1 // Refine elements with average edge length >
00013
00014 #define COARSEN_TOL 0.8 // Coarsen element with average edge length <
00015
00016 #define QUALITY_MIN 0.7
00017
00018 class FEM_Adapt_Algs;
00019 CtvExtern(FEM_Adapt_Algs *, _adaptAlgs);
00020
00021 class femMeshModify;
00022 class FEM_Adapt;
00023 class FEM_AdaptL;
00024
00026
00029 class FEM_Adapt_Algs {
00030 friend class FEM_AdaptL;
00031 friend class FEM_Adapt;
00032 friend class femMeshModify;
00033 friend class FEM_Interpolate;
00034 friend class FEM_MUtil;
00035
00036 public:
00038 int coord_attr;
00040 int bc_attr;
00041
00042 protected:
00044 FEM_Mesh *theMesh;
00046 femMeshModify *theMod;
00048 FEM_AdaptL *theAdaptor;
00050 int numNodes;
00052 int numElements;
00054 int dim;
00056 typedef struct {
00057 int elID;
00058 double len;
00059 } elemHeap;
00061 elemHeap *coarsenElements;
00063 elemHeap *refineElements;
00065 elemHeap *refineStack;
00067 int coarsenHeapSize;
00069 int refineHeapSize;
00071 int refineTop;
00073 int *elemConn;
00075 double *coordsn1, *coordsn2, *coordsn3;
00076
00077 public:
00079 FEM_Adapt_Algs() {
00080 theMesh = NULL; theMod = NULL; theAdaptor = NULL;
00081 }
00083 FEM_Adapt_Algs(FEM_Mesh *m, femMeshModify *fm);
00085 FEM_Adapt_Algs(femMeshModify *fm);
00087 void FEM_AdaptAlgsSetMesh(FEM_Mesh *m) {theMesh = m;}
00089 void FEM_Adapt_Algs_Init(int coord_at, int bc_at, int dimension) {
00090 coord_attr = coord_at;
00091 bc_attr = bc_at;
00092 dim = dimension;
00093 }
00095 ~FEM_Adapt_Algs();
00097 void pup(PUP::er &p) {
00098 p|coord_attr;
00099 p|bc_attr;
00100 p|dim;
00101 }
00102
00104 void FEM_Refine(int qm, int method, double factor, double *sizes);
00106 void FEM_Coarsen(int qm, int method, double factor, double *sizes);
00108 void FEM_AdaptMesh(int qm, int method, double factor, double *sizes);
00110 void FEM_Smooth(int qm, int method);
00112 void FEM_mesh_smooth(FEM_Mesh *meshP, int *nodes, int nNodes, int attrNo);
00114 void FEM_Repair(int qm);
00116 void FEM_Remesh(int qm, int method, double factor, double *sizes);
00117
00119 void SetReferenceMesh();
00121 void GradateMesh(double smoothness);
00122
00123 private:
00124
00126 int Refine(int qm, int method, double factor, double *sizes);
00128 int Coarsen(int qm, int method, double factor, double *sizes);
00129
00131 void SetMeshSize(int method, double factor, double *sizes);
00133 void Insert(int elID, double len, int cFlag);
00135 int Delete_Min(int cFlag);
00136
00137 public:
00139 int refine_element_leb(int e);
00141 void refine_flip_element_leb(int e, int p, int n1, int n2,
00142 double le);
00143
00145 int simple_refine(double targetA, double xmin=0.0, double ymin=0.0, double xmax=1.0, double ymax=1.0);
00147 int simple_coarsen(double targetA, double xmin=0.0, double ymin=0.0, double xmax=1.0, double ymax=1.0);
00148
00150 double length(int n1, int n2);
00152 double length(double *n1_coord, double *n2_coord);
00154 double getArea(int n1, int n2, int n3);
00156 double getArea(double *n1_coord, double *n2_coord, double *n3_coord);
00158 double getSignedArea(int n1, int n2, int n3);
00160 double getSignedArea(double *n1_coord, double *n2_coord, double *n3_coord);
00162 int getCoord(int n1, double *crds);
00164 int getShortestEdge(int n1, int n2, int n3, int* shortestEdge);
00165
00167 double getAreaQuality(int elem);
00169 void ensureQuality(int n1, int n2, int n3);
00171 bool controlQualityF(int n1, int n2, int n3, int n4);
00173 bool controlQualityR(int n1, int n2, int n3, int n4);
00175 bool controlQualityR(double *n1_coord, double *n2_coord, double *n3_coord);
00177 bool controlQualityC(int n1, int n2, int n3, double *n4_coord);
00179 bool controlQualityC(double *n1_coord, double *n2_coord, double *n3_coord, double *n4_coord);
00180
00182 bool flipOrBisect(int elId, int n1, int n2, int maxEdgeIdx, double maxlen);
00184 void tests(bool b);
00185 };
00186
00187
00188
00189 #endif