00001
00002
00003
00004
00005
00006
00007 #include "ParFUM.h"
00008 #include "ParFUM_internals.h"
00009
00010
00011 CLINKAGE void FEM_REF_INIT(int mesh) {
00012 CkArrayID femRefId;
00013 int cid;
00014 int size;
00015 TCharm *tc=TCharm::get();
00016 MPI_Comm comm = MPI_COMM_WORLD;
00017 MPI_Comm_rank(comm,&cid);
00018 MPI_Comm_size(comm,&size);
00019 if(cid==0) {
00020 CkArrayOptions opts;
00021 opts.bindTo(tc->getProxy());
00022 femMeshModMsg *fm = new femMeshModMsg;
00023 femRefId = CProxy_femMeshModify::ckNew(fm, opts);
00024 }
00025 MPI_Bcast(&femRefId, sizeof(CkArrayID), MPI_BYTE, 0, comm);
00026 meshMod = femRefId;
00027 femMeshModMsg *fm = new femMeshModMsg(size,cid);
00028 meshMod[cid].insert(fm);
00029 FEM_Mesh *m=FEM_Mesh_lookup(mesh,"FEM_REF_INIT");
00030 FEMMeshMsg *msg = new FEMMeshMsg(m,tc);
00031 meshMod[cid].setFemMesh(msg);
00032 return;
00033 }
00034
00035 void FEM_ADAPT_Init(int meshID) {
00036 const int triangleFaces[6] = {0,1,1,2,2,0};
00037 FEM_Add_elem2face_tuples(meshID, 0, 2, 3, triangleFaces);
00038 FEM_Mesh_allocate_valid_attr(meshID, FEM_ELEM+0);
00039 FEM_Mesh_allocate_valid_attr(meshID, FEM_NODE);
00040 FEM_Mesh_create_elem_elem_adjacency(meshID);
00041 FEM_Mesh_create_node_elem_adjacency(meshID);
00042 FEM_Mesh_create_node_node_adjacency(meshID);
00043
00044 FEM_REF_INIT(meshID);
00045 FEM_Mesh *meshP = FEM_Mesh_lookup(meshID, "FEM_ADAPT_Init");
00046 CtvInitialize(FEM_Adapt_Algs *, _adaptAlgs);
00047 CtvAccess(_adaptAlgs) = meshP->getfmMM()->getfmAdaptAlgs();
00048 CtvAccess(_adaptAlgs)->FEM_Adapt_Algs_Init(FEM_DATA+0, FEM_BOUNDARY,2);
00049 }
00050 FLINKAGE void FTN_NAME(FEM_ADAPT_INIT,fem_adapt_init)(int *meshID)
00051 {
00052 FEM_ADAPT_Init(*meshID);
00053 }
00054
00055
00056
00057 void FEM_ADAPT_Refine(int meshID, int qm, int method, double factor,
00058 double *sizes) {
00059 FEM_Mesh* mesh = FEM_Mesh_lookup(meshID, "FEM_ADAPT_Refine");
00060 mesh->getfmMM()->getfmAdaptAlgs()->FEM_Refine(qm, method, factor, sizes);
00061 }
00062 FLINKAGE void FTN_NAME(FEM_ADAPT_REFINE,fem_adapt_refine)(int* meshID,
00063 int *qm, int *method, double *factor, double *sizes)
00064 {
00065 FEM_ADAPT_Refine(*meshID, *qm, *method, *factor, sizes);
00066 }
00067
00068
00069 void FEM_ADAPT_Coarsen(int meshID, int qm, int method, double factor,
00070 double *sizes) {
00071 FEM_Mesh* mesh = FEM_Mesh_lookup(meshID, "FEM_ADAPT_Coarsen");
00072 mesh->getfmMM()->getfmAdaptAlgs()->FEM_Coarsen(qm, method, factor, sizes);
00073 }
00074 FLINKAGE void FTN_NAME(FEM_ADAPT_COARSEN,fem_adapt_coarsen)(int* meshID,
00075 int *qm, int *method, double *factor, double *sizes)
00076 {
00077 FEM_ADAPT_Coarsen(*meshID, *qm, *method, *factor, sizes);
00078 }
00079
00080
00081 void FEM_ADAPT_AdaptMesh(int meshID, int qm, int method, double factor,
00082 double *sizes) {
00083 FEM_Mesh* mesh = FEM_Mesh_lookup(meshID, "FEM_ADAPT_AdaptMesh");
00084 mesh->getfmMM()->getfmAdaptAlgs()->FEM_AdaptMesh(qm, method, factor, sizes);
00085 }
00086 FLINKAGE void FTN_NAME(FEM_ADAPT_ADAPTMESH,fem_adapt_adaptmesh)(int* meshID,
00087 int *qm, int *method, double *factor, double *sizes)
00088 {
00089 FEM_ADAPT_AdaptMesh(*meshID, *qm, *method, *factor, sizes);
00090 }
00091
00092
00093 void FEM_ADAPT_SetElementSizeField(int meshID, int elem, double size) {
00094 FEM_Mesh *meshP = FEM_Mesh_lookup(meshID, "FEM_ADAPT_SetElementSizeField");
00095 meshP->elem[0].setMeshSizing(elem, size);
00096 }
00097 FLINKAGE void FTN_NAME(FEM_ADAPT_SETELEMENTSIZEFIELD,fem_adapt_setelementsizefield)(int *meshID, int *elem, double *size)
00098 {
00099 FEM_ADAPT_SetElementSizeField(*meshID, *elem, *size);
00100 }
00101
00102
00103 void FEM_ADAPT_SetElementsSizeField(int meshID, double *sizes) {
00104 FEM_Mesh *meshP = FEM_Mesh_lookup(meshID, "FEM_ADAPT_SetElementsSizeField");
00105 int numElements = meshP->elem[0].size();
00106 for (int i=0; i<numElements; i++) {
00107 meshP->elem[0].setMeshSizing(i, sizes[i]);
00108 }
00109 }
00110 FLINKAGE void FTN_NAME(FEM_ADAPT_SETELEMENTSSIZEFIELD,fem_adapt_setelementssizefield)(int *meshID, double *sizes)
00111 {
00112 FEM_ADAPT_SetElementsSizeField(*meshID, sizes);
00113 }
00114
00115
00116 void FEM_ADAPT_SetReferenceMesh(int meshID) {
00117 FEM_Mesh* mesh = FEM_Mesh_lookup(meshID, "FEM_ADAPT_SetReferenceMesh");
00118 mesh->getfmMM()->getfmAdaptAlgs()->SetReferenceMesh();
00119 }
00120 FLINKAGE void FTN_NAME(FEM_ADAPT_SETREFERENCEMESH, fem_adapt_setreferencemesh)(int* meshID)
00121 {
00122 FEM_ADAPT_SetReferenceMesh(*meshID);
00123 }
00124
00125
00126 void FEM_ADAPT_GradateMesh(int meshID, double smoothness)
00127 {
00128 FEM_Mesh* mesh = FEM_Mesh_lookup(meshID, "FEM_ADAPT_GradateMesh");
00129 mesh->getfmMM()->getfmAdaptAlgs()->GradateMesh(smoothness);
00130 }
00131 FLINKAGE void FTN_NAME(FEM_ADAPT_GRADATEMESH, fem_adapt_gradatemesh)(int* meshID, double* smoothness)
00132 {
00133 FEM_ADAPT_GradateMesh(*meshID, *smoothness);
00134 }
00135
00136 void FEM_ADAPT_TestMesh(int meshID) {
00137 FEM_Mesh* mesh = FEM_Mesh_lookup(meshID, "FEM_ADAPT_GradateMesh");
00138 mesh->getfmMM()->getfmAdaptAlgs()->tests(true);
00139 }
00140 FLINKAGE void FTN_NAME(FEM_ADAPT_TESTMESH, fem_adapt_testmesh)(int* meshID)
00141 {
00142 FEM_ADAPT_TestMesh(*meshID);
00143 }
00144
00145 int FEM_ADAPT_SimpleRefineMesh(int meshID, double targetA, double xmin, double ymin, double xmax, double ymax) {
00146 FEM_Mesh* mesh = FEM_Mesh_lookup(meshID, "FEM_ADAPT_GradateMesh");
00147 return mesh->getfmMM()->getfmAdaptAlgs()->simple_refine(targetA,xmin,ymin,xmax,ymax);
00148 }
00149 FLINKAGE void FTN_NAME(FEM_ADAPT_SIMPLEREFINEMESH, fem_adapt_simplerefinemesh)(int* meshID, double* targetA, double* xmin, double* ymin, double* xmax, double* ymax)
00150 {
00151 FEM_ADAPT_SimpleRefineMesh(*meshID,*targetA,*xmin,*ymin,*xmax,*ymax);
00152 }
00153
00154 int FEM_ADAPT_SimpleCoarsenMesh(int meshID, double targetA, double xmin, double ymin, double xmax, double ymax) {
00155 FEM_Mesh* mesh = FEM_Mesh_lookup(meshID, "FEM_ADAPT_GradateMesh");
00156 return mesh->getfmMM()->getfmAdaptAlgs()->simple_coarsen(targetA,xmin,ymin,xmax,ymax);
00157 }
00158 FLINKAGE void FTN_NAME(FEM_ADAPT_SIMPLECOARSENMESH, fem_adapt_simplecoarsenmesh)(int* meshID, double* targetA, double* xmin, double* ymin, double* xmax, double* ymax)
00159 {
00160 FEM_ADAPT_SimpleCoarsenMesh(*meshID,*targetA,*xmin,*ymin,*xmax,*ymax);
00161 }