-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathMesh.C
182 lines (160 loc) · 6.03 KB
/
Mesh.C
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
// Local header files
//#include "Mesh.h"
#include <vector>
#include "libmesh/elem.h"
#include "libmesh/face_tri3.h"
#include "libmesh/mesh.h"
#include "libmesh/mesh_generation.h"
#include "libmesh/mesh_tetgen_interface.h"
#include "libmesh/mesh_triangle_interface.h"
#include "libmesh/node.h"
#include "libmesh/serial_mesh.h"
// Bring in everything from the libMesh namespace using namespace libMesh;
using namespace libMesh;
// Helper routine for tetrahedralize_domain(). Adds the points and elements
// of a convex hull generated by TetGen to the input mesh
void add_sphere_convex_hull_to_mesh(MeshBase& mesh, libMesh::Real r_max, unsigned int points_on_sphere, std::vector<Node> geometry, std::string creator, const Real L, const unsigned int N);
unsigned int closest(std::vector<libMesh::Node> geom, libMesh::Point pt);
void tetrahedralise_sphere(UnstructuredMesh& mesh, std::vector<Node> geometry, std::string creator, Real r, int NrBall, Real VolConst, Real L, unsigned int N){
#ifdef LIBMESH_HAVE_TETGEN
//libMesh::Real r=10.;
add_sphere_convex_hull_to_mesh(mesh, r, NrBall, geometry, creator, L, N);
// 3.) Update neighbor information so that TetGen can verify there is a convex hull.
mesh.find_neighbors();
// 5.) Set parameters and tetrahedralize the domain
// 0 means "use TetGen default value"
//Real quality_constraint = 42*10; // practically no constraint
Real quality_constraint = 2;
// The volume constraint determines the max-allowed tetrahedral
// volume in the Mesh. TetGen will split cells which are larger than
// this size
Real volume_constraint = VolConst;
// Construct the Delaunay tetrahedralization
TetGenMeshInterface t(mesh);
t.pointset_convexhull();
mesh.find_neighbors();
t.triangulate_conformingDelaunayMesh(quality_constraint, volume_constraint);
// Find neighbors, etc in preparation for writing out the Mesh
mesh.prepare_for_use();
// Finally, write out the result --> is done in main already.
//mesh.write("sphere_3D.e");
#else
// Avoid compiler warnings
libmesh_ignore(mesh.comm);
#endif // LIBMESH_HAVE_TETGEN
}
std::vector<Point> fibonacci(unsigned int points_on_sphere){
std::vector<Point> points(points_on_sphere);
double x,y,z,r,phi;
double pi=3.1415926;
for (unsigned int i=0; i<points_on_sphere; i++){
y=2.*i/points_on_sphere-1.;
r=sqrt(1.-y*y);
phi=i*pi*(3.-sqrt(5.));
z=r*sin(phi);
x=r*cos(phi);
points[i]=Node(x,y,z);
//out<<points[i]<<std::endl;
}
//out<<std::endl<<std::endl;
return points;
}
std::vector<Point> archimedes(unsigned int points_on_sphere){
int n=floor(sqrt(points_on_sphere));
std::vector<Point> points(n*n);
double theta, phi;
for(int i=0; i<n; i++){
theta=i*6.2831852/n;
for(int j=1; j<=n; j++){
phi=acos(2.*j/n -1.);
points[i*n+j-1]=Node(cos(theta)*sin(phi),
sin(theta)*sin(phi),
2.*j/(n+1.)-1.);
//out<<points[i*n+j-1]<<std::endl;
}
}
return points;
}
std::vector<Point> spiral(unsigned int points_on_sphere){
std::vector<Point> points(points_on_sphere);
double pi=3.1415926;
double dl=pi*(3.-sqrt(5.));
double dz=2./(points_on_sphere +1.);
double z, l, r;
for(unsigned int i=1; i<=points_on_sphere; i++){
z=i*dz-1.;
l=i*dl;
r=sqrt(1.-z*z);
points[i-1]=Node(r*cos(l),r*sin(l),z);
}
return points;
}
void add_sphere_convex_hull_to_mesh(MeshBase& mesh, libMesh::Real r_max, unsigned int points_on_sphere, std::vector<Node> geometry, std::string creator, const Real L, const unsigned int N){
#ifdef LIBMESH_HAVE_TETGEN
std::vector<Point> point;
double x, scale;
// For each point in the map, insert it into the input mesh and copy it to all nuclear sites.
// keep track of the ID assigned.
unsigned int molsize=geometry.size();
unsigned int pts_circle;
// add one node each at the nuclear position.
for(unsigned int i=0; i<molsize; i++){
// add one point each at the nuclear position:
mesh.add_point( geometry[i]);
}
//now, add further points to the mesh: The nodes are located on spheres (circle)
// around the molecules; each sphere has a different radius (scale) in (0,r_max].
//outermost loop: over different spheres
for(unsigned int circle=0; circle<N; circle++){
x = (double)(2.*circle)/N - 1.;
scale = L*(1.-x)/(1.+x+2.*L/r_max);
if (scale == 0.)
continue;
// create a unit sphere of respective type with slowly incr. number of points
pts_circle=(int)(r_max*scale*0.3+1)*points_on_sphere;
//pts_circle=points_on_sphere;
if (creator=="fibonacci")
point=fibonacci(pts_circle);
else if (creator=="archimedes")
point=archimedes(pts_circle);
else if (creator=="spiral")
point=spiral(pts_circle);
// loop over atoms: each of them has a shpere around it
for(unsigned int i=0; i<molsize; i++){
// loop over the points on the sphere and add it respectively.
for (unsigned int it=0; it<point.size(); it++){
// shift and scale the coordinate as needed
point[it]*=scale;
point[it]+=geometry[i];
if (i==closest(geometry, point[it])){
mesh.add_point( point[it] );
}
// shift and scale back for use in the next round.
point[it]-= geometry[i];
point[it]/= scale;
}
}
}
// now, the point set is in mesh and we can call triangulate_pointset().
#else
// Avoid compiler warnings
libmesh_ignore(mesh);
libmesh_ignore(lower_limit);
libmesh_ignore(upper_limit);
#endif // LIBMESH_HAVE_TETGEN
}
unsigned int closest(std::vector<Node> geom, Point pt){
unsigned int molsize=geom.size();
//most outer loop: nuclear sites
unsigned int closest=0;
Point cl= geom[0]-pt;
Point ac;
for (unsigned int i=1; i<molsize; i++ ){
ac= geom[i]-pt;
if ( ac.size()<cl.size() ){
closest = i;
cl = ac;
}
}
return closest;
}