Line | Branch | Exec | Source |
---|---|---|---|
1 | // This file is part of Bembel, the higher order C++ boundary element library. | ||
2 | // | ||
3 | // Copyright (C) 2024 see <http://www.bembel.eu> | ||
4 | // | ||
5 | // It was written as part of a cooperation of J. Doelz, H. Harbrecht, S. Kurz, | ||
6 | // M. Multerer, S. Schoeps, and F. Wolf at Technische Universitaet Darmstadt, | ||
7 | // Universitaet Basel, and Universita della Svizzera italiana, Lugano. This | ||
8 | // source code is subject to the GNU General Public License version 3 and | ||
9 | // provided WITHOUT ANY WARRANTY, see <http://www.bembel.eu> for further | ||
10 | // information. | ||
11 | |||
12 | #ifndef BEMBEL_SRC_GEOMETRY_GEOMETRYIO_HPP_ | ||
13 | #define BEMBEL_SRC_GEOMETRY_GEOMETRYIO_HPP_ | ||
14 | |||
15 | namespace Bembel { | ||
16 | |||
17 | /** | ||
18 | * \ingroup Geometry | ||
19 | * \brief loads geometry from file with GEOPDE-format. Note that the direction | ||
20 | * of the normals must be consistent | ||
21 | * | ||
22 | * \param file_name path/filename pointing to the geometry file | ||
23 | * \return std::vector of NURBS::Patch describing geometry | ||
24 | */ | ||
25 | 26 | inline std::vector<Patch> LoadGeometryFileDAT( | |
26 | const std::string &file_name) noexcept { | ||
27 | 26 | std::vector<Bembel::Patch> out; | |
28 | 26 | std::stringstream iss; | |
29 | 26 | std::string word; | |
30 | 26 | std::string row; | |
31 | int infoInt[5]; | ||
32 | 26 | std::ifstream file; | |
33 | 26 | file.open(file_name); | |
34 |
1/2✗ Branch 1 not taken.
✓ Branch 2 taken 26 times.
|
26 | if (!file) { |
35 | ✗ | std::cerr << "File " << file_name << " doesn't exist!"; | |
36 | ✗ | exit(1); | |
37 | } | ||
38 | // first 4 rows are irrelevant | ||
39 |
2/2✓ Branch 0 taken 130 times.
✓ Branch 1 taken 26 times.
|
156 | for (int i = 0; i < 5; i++) { |
40 | 130 | getline(file, row); | |
41 | } | ||
42 | // main information | ||
43 | 26 | iss.str(row); | |
44 |
2/2✓ Branch 0 taken 130 times.
✓ Branch 1 taken 26 times.
|
156 | for (int i = 0; i < 5; i++) { |
45 | 130 | iss >> word; | |
46 | |||
47 | 130 | infoInt[i] = stoi(word); | |
48 | } | ||
49 | // main loop - patches | ||
50 |
2/2✓ Branch 0 taken 88 times.
✓ Branch 1 taken 26 times.
|
114 | for (int patchNr = 1; patchNr <= infoInt[2]; patchNr++) { |
51 | 88 | Bembel::Patch tempPatch; | |
52 | 88 | std::vector<double> tempknt1; | |
53 | 88 | std::vector<double> tempknt2; | |
54 | 88 | std::vector<int> info; // p and ncp / 0,1-p 2,3-ncp | |
55 | 88 | std::vector<Eigen::MatrixXd> tmp; | |
56 | |||
57 | 88 | getline(file, row); // file_name | |
58 | |||
59 | 88 | getline(file, row); | |
60 | 88 | iss.str(row); | |
61 |
2/2✓ Branch 2 taken 176 times.
✓ Branch 3 taken 88 times.
|
264 | while (iss >> word) { |
62 | 176 | info.push_back(stoi(word)); | |
63 | } | ||
64 | 88 | getline(file, row); | |
65 | |||
66 | 88 | word = ""; | |
67 | 88 | iss.clear(); | |
68 | 88 | iss.str(row); | |
69 |
2/2✓ Branch 2 taken 176 times.
✓ Branch 3 taken 88 times.
|
264 | while (iss >> word) { |
70 | 176 | info.push_back(stoi(word)); | |
71 | } | ||
72 | |||
73 | 88 | word = ""; | |
74 | 88 | iss.clear(); | |
75 | |||
76 | // In textfiles are only 2 knotVectors and 4 Matrices | ||
77 | |||
78 | 88 | getline(file, row); | |
79 | 88 | iss.str(row); | |
80 | |||
81 | // first knotVector | ||
82 | |||
83 |
2/2✓ Branch 2 taken 629 times.
✓ Branch 3 taken 88 times.
|
717 | for (int k = 0; k < info[0] + info[2] + 1; k++) { |
84 | 629 | iss >> word; | |
85 | 629 | tempknt1.push_back(atof(word.c_str())); | |
86 | } | ||
87 | 88 | iss.clear(); | |
88 | 88 | getline(file, row); | |
89 | 88 | iss.str(row); | |
90 | 88 | word = ""; | |
91 | // second knotVector | ||
92 |
2/2✓ Branch 2 taken 631 times.
✓ Branch 3 taken 88 times.
|
719 | for (int k = 0; k < info[1] + info[3] + 1; k++) { |
93 | 631 | iss >> word; | |
94 | 631 | tempknt2.push_back(atof(word.c_str())); | |
95 | } | ||
96 | 88 | iss.clear(); | |
97 | // 4 Matrices | ||
98 | 88 | int N = info[2]; | |
99 | 88 | int M = info[3]; | |
100 |
2/2✓ Branch 0 taken 352 times.
✓ Branch 1 taken 88 times.
|
440 | for (int k = 0; k < 4; k++) { |
101 | Eigen::Matrix<double, -1, -1> tempMatrix( | ||
102 | 352 | M, N); // == Eigen::MatrixXd tempMatrix(M,N); | |
103 | 352 | getline(file, row); | |
104 | 352 | iss.str(row); | |
105 |
2/2✓ Branch 0 taken 1264 times.
✓ Branch 1 taken 352 times.
|
1616 | for (int i = 0; i < M; i++) |
106 |
2/2✓ Branch 0 taken 5312 times.
✓ Branch 1 taken 1264 times.
|
6576 | for (int j = 0; j < N; j++) { |
107 | 5312 | iss >> word; | |
108 | |||
109 | 5312 | tempMatrix(i, j) = atof(word.c_str()); | |
110 | } | ||
111 | |||
112 | 352 | tmp.push_back(tempMatrix); | |
113 | 352 | iss.clear(); | |
114 | 352 | } | |
115 | // Important | ||
116 | 88 | tempPatch.init_Patch(tmp, tempknt1, tempknt2); | |
117 | 88 | out.push_back(tempPatch); | |
118 | 88 | } | |
119 | |||
120 | 26 | file.close(); | |
121 | 26 | return out; | |
122 | 26 | } | |
123 | |||
124 | /** | ||
125 | * \ingroup Gemetry | ||
126 | * \brief method to generate textfile for the geometry | ||
127 | * \param file_name name of new textfile | ||
128 | * \param number_of_patches overall number of patches | ||
129 | **/ | ||
130 | 1 | inline void MakeFile(const std::string &file_name, | |
131 | int number_of_patches) noexcept { | ||
132 | 1 | std::ofstream file; | |
133 | 1 | file.open(file_name); | |
134 | file << "# nurbs mesh v.2.1" | ||
135 | 1 | << "\r\n"; | |
136 | 1 | file << "# " << file_name << "\r\n"; | |
137 | file << "# Generated by BEMBEL, see www.bembel.eu" | ||
138 | 1 | << "\r\n"; | |
139 | file << "#" | ||
140 | 1 | << "\r\n"; | |
141 | 1 | file << "2 3 " << number_of_patches << " 0 0 " | |
142 | 1 | << "\r\n"; | |
143 | 1 | file.close(); | |
144 | 1 | } | |
145 | /** | ||
146 | * \ingroup Gemetry | ||
147 | * \brief method to write Patch information into textfile | ||
148 | * \param knt1 knotVector1 | ||
149 | * \param knt2 knotVector2 | ||
150 | * \param xyzw Vector with x,y,z,w Matices | ||
151 | * \param file_name filename | ||
152 | * \param current_patch_number current patch number | ||
153 | **/ | ||
154 | 2 | void WritePatch(const std::string &file_name, int current_patch_number, | |
155 | const std::vector<Eigen::MatrixXd> &xyzw, | ||
156 | const std::vector<double> &knt1, | ||
157 | const std::vector<double> &knt2) noexcept { | ||
158 | 2 | std::ofstream file(file_name, std::ios::app); | |
159 | 2 | int N = xyzw[0].cols(); // ncp | |
160 | 2 | int M = xyzw[0].rows(); // ncp | |
161 | 2 | int p1 = knt1.size() - N - 1; | |
162 | 2 | int p2 = knt2.size() - M - 1; | |
163 | 2 | file << "PATCH " << current_patch_number << " \r\n"; | |
164 | 2 | file << p1 << " " << p2 << " \r\n"; | |
165 | 2 | file << N << " " << M << " \r\n"; | |
166 | // knotVectors | ||
167 |
2/2✓ Branch 1 taken 18 times.
✓ Branch 2 taken 2 times.
|
20 | for (unsigned int i = 0; i < knt1.size() - 1; i++) { |
168 | 18 | file << std::fixed << std::setprecision(15) << knt1[i] << " "; | |
169 | } | ||
170 | 2 | file << std::fixed << std::setprecision(15) << knt1[knt1.size() - 1] | |
171 | 2 | << " \r\n"; | |
172 |
2/2✓ Branch 1 taken 20 times.
✓ Branch 2 taken 2 times.
|
22 | for (unsigned int i = 0; i < knt2.size() - 1; i++) { |
173 | 20 | file << std::fixed << std::setprecision(15) << knt2[i] << " "; | |
174 | } | ||
175 | 2 | file << std::fixed << std::setprecision(15) << knt2[knt2.size() - 1] | |
176 | 2 | << " \r\n"; | |
177 | |||
178 | // Matrices | ||
179 |
2/2✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
|
10 | for (unsigned int n = 0; n < xyzw.size(); n++) { |
180 |
2/2✓ Branch 0 taken 44 times.
✓ Branch 1 taken 8 times.
|
52 | for (int i = 0; i < M; i++) |
181 |
2/2✓ Branch 0 taken 220 times.
✓ Branch 1 taken 44 times.
|
264 | for (int j = 0; j < N; j++) { |
182 | 220 | file << std::fixed << std::setprecision(15) << xyzw[n](i, j); | |
183 |
2/2✓ Branch 0 taken 8 times.
✓ Branch 1 taken 212 times.
|
220 | if (N * i + j == M * N - 1) { |
184 | 8 | file << " \r\n"; | |
185 | 8 | continue; | |
186 | } | ||
187 | 212 | file << " "; | |
188 | } | ||
189 | } | ||
190 | |||
191 | 2 | file.close(); | |
192 | 2 | } | |
193 | |||
194 | /** | ||
195 | * \ingroup Geometry | ||
196 | * \brief exports a geometry from Bembel to a .dat file. | ||
197 | * | ||
198 | * Limitation: This functions assumes a p-open knot vector without internal | ||
199 | * knots. | ||
200 | * | ||
201 | * \param Geometry std::vector of NURBS::Patch describing geometry | ||
202 | * \param file_name path/filename to be written to | ||
203 | */ | ||
204 | 1 | void WriteDATFile(const std::vector<Patch> &Geometry, | |
205 | const std::string &file_name) { | ||
206 | 1 | const int number_of_patches = Geometry.size(); | |
207 | 1 | MakeFile(file_name, number_of_patches); | |
208 | |||
209 | 1 | int patch_count = 0; | |
210 |
2/2✓ Branch 5 taken 2 times.
✓ Branch 6 taken 1 times.
|
3 | for (auto &patch : Geometry) { |
211 |
1/2✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
|
2 | assert(patch.unique_knots_x_.size() == 2 && |
212 | "I assume 0 and 1 as unique knots!"); | ||
213 |
1/2✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
|
2 | assert(patch.unique_knots_y_.size() == 2 && |
214 | "I assume 0 and 1 as unique knots!"); | ||
215 | |||
216 | 2 | const int polynomial_degree_x = patch.polynomial_degree_x_; | |
217 | 2 | const int polynomial_degree_y = patch.polynomial_degree_y_; | |
218 | |||
219 |
1/2✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
|
2 | std::vector<double> knt_x(2 * polynomial_degree_x, 1.0); |
220 |
1/2✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
|
2 | std::vector<double> knt_y(2 * polynomial_degree_y, 1.0); |
221 |
2/2✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
|
12 | for (auto i = 0; i < polynomial_degree_x; ++i) knt_x[i] = 0.0; |
222 |
2/2✓ Branch 1 taken 11 times.
✓ Branch 2 taken 2 times.
|
13 | for (auto i = 0; i < polynomial_degree_y; ++i) knt_y[i] = 0.0; |
223 | |||
224 | 2 | const int number_of_points_x = knt_x.size() - polynomial_degree_x; | |
225 | 2 | const int number_of_points_y = knt_y.size() - polynomial_degree_y; | |
226 | |||
227 | std::vector<Eigen::MatrixXd> data( | ||
228 |
3/6✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
|
4 | 4, Eigen::MatrixXd::Zero(number_of_points_y, number_of_points_x)); |
229 | |||
230 | 2 | const int matrix_size = number_of_points_x * number_of_points_y; | |
231 |
2/2✓ Branch 0 taken 55 times.
✓ Branch 1 taken 2 times.
|
57 | for (auto i = 0; i < matrix_size; ++i) { |
232 | 55 | const int rowIndex = i % number_of_points_y; | |
233 | 55 | const int colIndex = i / number_of_points_y; | |
234 |
1/2✓ Branch 3 taken 55 times.
✗ Branch 4 not taken.
|
55 | data[0](rowIndex, colIndex) = patch.data_[4 * i]; |
235 |
1/2✓ Branch 3 taken 55 times.
✗ Branch 4 not taken.
|
55 | data[1](rowIndex, colIndex) = patch.data_[4 * i + 1]; |
236 |
1/2✓ Branch 3 taken 55 times.
✗ Branch 4 not taken.
|
55 | data[2](rowIndex, colIndex) = patch.data_[4 * i + 2]; |
237 |
1/2✓ Branch 3 taken 55 times.
✗ Branch 4 not taken.
|
55 | data[3](rowIndex, colIndex) = patch.data_[4 * i + 3]; |
238 | } | ||
239 | |||
240 | 2 | WritePatch(file_name, patch_count, data, knt_x, knt_y); | |
241 | 2 | ++patch_count; | |
242 | 2 | } | |
243 | |||
244 | 1 | return; | |
245 | } | ||
246 | |||
247 | } // namespace Bembel | ||
248 | #endif // BEMBEL_SRC_GEOMETRY_GEOMETRYIO_HPP_ | ||
249 |