-
Notifications
You must be signed in to change notification settings - Fork 0
/
allinone.cpp
214 lines (197 loc) · 9.32 KB
/
allinone.cpp
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
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
#include <iostream>
#include <string>
#include <vector>
#include <cmath>
#include <typeinfo>
#include <functional>
#include "include/elements.hpp"
#include "include/readfile.hpp"
#include "include/structure.hpp"
#include "include/quadrature.hpp"
#include "include/cgsolver.hpp"
#include "include/csrmatrixcpu.hpp"
#include "include/vectorcpu.hpp"
#include <omp.h>
#include <stdio.h>
using namespace std;
int main(int argc, char* argv[])
{
float walltime[2];
std::cout << ">- CPU: MESH FORMAT, Q2 -<" << std::endl;
std::vector<Vertex<float>> nodes;
std::vector<Element<float>*> elements;
file_to_mesh_all("../data/square_quadrilateral_q1_b.msh", nodes, elements);
std::cout << "num nodes: " << nodes.size() << std::endl;
std::cout << "num elements: " << elements.size() << std::endl;
CsrMatrixCpu<float> mat(nodes.size());
std::cout << "structure" << std::flush;
walltime[0] = omp_get_wtime();
structure(mat, elements);
walltime[0] -= omp_get_wtime();
std::cout << " - done (" << -walltime[0] * 1000.0 << ")" << std::endl;
std::cout << "assemble" << std::flush;
walltime[1] = omp_get_wtime();
float qp{std::sqrt(0.6)};
float weight[9] = {25.0, 40.0, 25.0,
40.0, 64.0, 40.0,
25.0, 40.0, 25.0};
float quadpoint[9][2] = {{-qp, -qp},
{-qp, 0.0},
{-qp, qp},
{0.0, -qp},
{0.0, 0.0},
{0.0, qp},
{ qp, -qp},
{ qp, 0.0},
{ qp, qp}};
// function begin: assemble
const size_t num_elem{elements.size()};
#pragma omp parallel for// num_threads(16)
for (size_t e=0; e < num_elem; e++)
{
const std::vector<size_t> vertexids{elements[e]->vertexids()}; // only needed for vertexids.size() -> always known (=8)
std::array<float, 8> coords = static_cast<QuadrilateralQ1<float>*>(elements[e])->get_pointcoords();
for (size_t i{0}; i < 4; ++i)
{
for (size_t j{0}; j < 4; ++j)
{
float val{0.0};
for (size_t p{0}; p < 9; ++p)
{
float xi = quadpoint[p][0];
float eta = quadpoint[p][1];
float B[2][2] =
{ { ( -(1.0f-eta)*coords[0] + (1.0f-eta)*coords[1] + (1.0f+eta)*coords[2] - (1.0f+eta)*coords[3] ) * 0.25f ,
( -(1.0f-xi )*coords[0] - (1.0f+xi )*coords[1] + (1.0f+xi )*coords[2] + (1.0f-xi )*coords[3] ) * 0.25f },
{ ( -(1.0f-eta)*coords[4] + (1.0f-eta)*coords[5] + (1.0f+eta)*coords[6] - (1.0f+eta)*coords[7] ) * 0.25f ,
( -(1.0f-xi )*coords[4] - (1.0f+xi )*coords[5] + (1.0f+xi )*coords[6] + (1.0f-xi )*coords[7] ) * 0.25f } };
// { { ( (eta-1.0)*coords[0] + (1.0-eta)*coords[1] + (1.0+eta)*coords[2] + (eta-1.0)*coords[3] ) * 0.25 ,
// ( (xi -1.0)*coords[0] + (xi -1.0)*coords[1] + (1.0+xi )*coords[2] + (1.0-xi )*coords[3] ) * 0.25 },
// { ( (eta-1.0)*coords[4] + (1.0-eta)*coords[5] + (1.0+eta)*coords[6] + (eta-1.0)*coords[7] ) * 0.25 ,
// ( (xi -1.0)*coords[4] + (xi -1.0)*coords[5] + (1.0+xi )*coords[6] + (1.0-xi )*coords[7] ) * 0.25 } };
// help vars
std::array<float, 2> grad1;
std::array<float, 2> grad2;
if (i == 0)
grad1 = {(1.0f - eta) * (-0.25f) ,
(1.0f - xi ) * (-0.25f) };
else if (i == 1)
grad1 = {(1.0f - eta) * 0.25f ,
(1.0f + xi ) * (-0.25f) };
else if (i == 2)
grad1 = {(1.0f + eta) * 0.25f ,
(1.0f + xi ) * 0.25f };
else //if (i == 3)
grad1 = {(1.0f + eta) * (-0.25f) ,
(1.0f - xi ) * 0.25f };
if (j == 0)
grad2 = {(1.0f - eta) * (-0.25f) ,
(1.0f - xi ) * (-0.25f) };
else if (j == 1)
grad2 = {(1.0f - eta) * 0.25f ,
(1.0f + xi ) * (-0.25f) };
else if (j == 2)
grad2 = {(1.0f + eta) * 0.25f ,
(1.0f + xi ) * 0.25f };
else //if (j == 3)
grad2 = {(1.0f + eta) * (-0.25f) ,
(1.0f - xi ) * 0.25f };
val += weight[p]
* ( ( B[1][1] * grad1[0] - B[1][0] * grad1[1]) * ( B[1][1] * grad2[0] - B[1][0] * grad2[1])
+ (-B[0][1] * grad1[0] + B[0][0] * grad1[1]) * (-B[0][1] * grad2[0] + B[0][0] * grad2[1]) )
/ std::abs(B[0][0] * B[1][1] - B[0][1] * B[1][0]);
} // end for p (quadrature point)
val /= 81.0f; // all weights are .../81
//val *= 0.0123456790123456790123;
mat.add(vertexids[i], vertexids[j], val);
} // end for j
} // end for i
} // end for elements
// function end: assemble
walltime[1] -= omp_get_wtime();
std::cout << " - done (" << -walltime[1] * 1000.0 << ")" << std::endl;
for (int k(0); k < 10; ++k)
std::cout << "(0, " << mat._colind[k] << ") = " << mat._values[k] << std::endl;
// assemble rhs
std::function<float(float, float)> f = [](float x, float y)
{ return static_cast<float>(2.0) * (x - x*x + y - y*y); };
size_t numvertices{nodes.size()};
VectorCpu rhs(numvertices, 0.0);
for (const auto& e : elements)
{
const std::vector<size_t> nodeids = e->vertexids();
Quadrature<Element, float> quad(e);
for (size_t i{0}; i < nodeids.size(); ++i)
// rhs.add(nodeids[i], f(nodes[nodeids[i]].x, nodes[nodeids[i]].y) * quad.integrate_basisfunction(2, i));
rhs.add(nodeids[i], f(nodes[nodeids[i]].x, nodes[nodeids[i]].y) * quad.integrate_basisfunction(3, i));
}
// dirichlet boundary
for (const auto& n : nodes)
{
if (n.x == 0.0 || n.y == 0.0 || n.x == 1.0 || n.y == 1.0)
{
for (size_t i{mat._rowptr[n.id]}; i < mat._rowptr[n.id + 1]; ++i)
mat._values[i] = 0.0f;
mat.set(n.id, n.id, 1.0f);
rhs.set(n.id, 0.0f);
}
}
// solve LGS
std::cout << "solve" << std::flush;
CgSolver<CsrMatrixCpu<float>, VectorCpu> solver(mat, rhs);
VectorCpu res(numvertices, 0.1);
solver.solve(res);
std::cout << " - done" << std::endl;
bool triangles_used = false; // false -> quadrilaterals
// write vtk-file
ofstream output("../data/allinone.vtk");
output << "# vtk DataFile Version 3.0" << std::endl;
output << "square q2" << std::endl;
output << "ASCII" << std::endl;
output << "DATASET UNSTRUCTURED_GRID" << std::endl;
output << std::endl;
output << "POINTS " << numvertices << (typeid(float) == typeid(float) ? " float" : " float") << std::endl;
for (const auto& n : nodes)
output << n.x << " " << n.y << " 0" << std::endl;
output << std::endl;
if (triangles_used) // bad style
output << "CELLS " << elements.size() << " " << 4*elements.size() << std::endl;
else
output << "CELLS " << elements.size() << " " << 5*elements.size() << std::endl;
for (const auto& e : elements)
{
//for (const auto id : e->vertexids())
//TODO
if (typeid(*e) == typeid(TriangleQ1<float>)
|| typeid(*e) == typeid(TriangleQ2<float>) )
output << "3 " << e->vertexids()[0] << " " << e->vertexids()[1] << " " << e->vertexids()[2] << std::endl;
else if (typeid(*e) == typeid(QuadrilateralQ1<float>)
|| typeid(*e) == typeid(QuadrilateralQ2<float>) )
output << "4 " << e->vertexids()[0] << " " << e->vertexids()[1] << " " << e->vertexids()[2] << " " << e->vertexids()[3] << std::endl;
else
assert(false);
}
output << std::endl;
output << "CELL_TYPES " << elements.size() << std::endl;
for (size_t i{0}; i < elements.size(); ++i)
{
if (typeid(*(elements[i])) == typeid(TriangleQ1<float>)
|| typeid(*(elements[i])) == typeid(TriangleQ2<float>) )
output << "5" << std::endl; // TriangleQ1
//output << "22" << std::endl; // TriangleQ2
else if (typeid(*(elements[i])) == typeid(QuadrilateralQ1<float>)
|| typeid(*(elements[i])) == typeid(QuadrilateralQ2<float>) )
output << "9" << std::endl; // QuadrilateralQ1
//output << "23" << std::endl; // QuadrilateralQ2
else
assert(false);
}
output << std::endl;
output << "POINT_DATA " << numvertices << std::endl;
output << "SCALARS u " << (typeid(float) == typeid(float) ? "float" : "float") << std::endl;
output << "LOOKUP_TABLE default" << std::endl;
for (size_t i{0}; i < numvertices; ++i)
output << (std::abs(res._values[i]) < 0.0001 ? 0 : res._values[i]) << std::endl;
output.close();
return 0;
}