test_FiniteElementMesh.cpp 21.1 KB
Newer Older
1
2
3
4
5
6
#include "test_Physics.hpp"

inline constexpr size_t operator "" _zu(unsigned long long int i) {
    return (size_t) i;
}

7
8
9
class MasterTriangle : public ::testing::Test {
public:
    virtual void SetUp() {
JasonPries's avatar
JasonPries committed
10
11
12
        node.push_back(XY{0.0, 0.0});
        node.push_back(XY{1.0, 0.0});
        node.push_back(XY{0.0, 1.0});
13
14
15

        tri.push_back(Triangle<1>(0, 1_zu, 2_zu, 0_zu));

JasonPries's avatar
JasonPries committed
16
17
        reg.push_back(Region<2>(std::vector<size_t>{0}));

JasonPries's avatar
JasonPries committed
18
19
20
        bnd.push_back(Boundary<2>(std::vector<size_t>{0, 1}));
        bnd.push_back(Boundary<2>(std::vector<size_t>{0, 2}));
        bnd.push_back(Boundary<2>(std::vector<size_t>{1, 2}));
JasonPries's avatar
JasonPries committed
21

JasonPries's avatar
JasonPries committed
22
        femesh = FiniteElementMesh<2, 1>(node, tri, reg, bnd);
23
24
25
26
    }

    std::vector<XY> node;
    std::vector<Triangle<1>> tri;
JasonPries's avatar
JasonPries committed
27
    std::vector<Region<2>> reg;
JasonPries's avatar
JasonPries committed
28
    std::vector<Boundary<2>> bnd;
29

JasonPries's avatar
JasonPries committed
30
    FiniteElementMesh<2, 1> femesh;
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
};

TEST_F(MasterTriangle, mass_matrix) {
    auto mat = femesh.basis<2>();
    auto det = femesh.determinant<2>();

    Eigen::SparseMatrix<double> M = (mat[0] * det(0) * mat[0].transpose()) * TriangleQuadratureRule<2>::w[0]
                                    + (mat[1] * det(1) * mat[1].transpose()) * TriangleQuadratureRule<2>::w[1]
                                    + (mat[2] * det(2) * mat[2].transpose()) * TriangleQuadratureRule<2>::w[2];

    EXPECT_NEAR(M.coeffRef(0, 0), 1.0 / 12.0, FLT_EPSILON);
    EXPECT_NEAR(M.coeffRef(1, 0), 1.0 / 24.0, FLT_EPSILON);
    EXPECT_NEAR(M.coeffRef(2, 0), 1.0 / 24.0, FLT_EPSILON);
    EXPECT_NEAR(M.coeffRef(0, 1), 1.0 / 24.0, FLT_EPSILON);
    EXPECT_NEAR(M.coeffRef(1, 1), 1.0 / 12.0, FLT_EPSILON);
    EXPECT_NEAR(M.coeffRef(2, 1), 1.0 / 24.0, FLT_EPSILON);
    EXPECT_NEAR(M.coeffRef(0, 2), 1.0 / 24.0, FLT_EPSILON);
    EXPECT_NEAR(M.coeffRef(1, 2), 1.0 / 24.0, FLT_EPSILON);
    EXPECT_NEAR(M.coeffRef(2, 2), 1.0 / 12.0, FLT_EPSILON);
}

TEST_F(MasterTriangle, stiffness_matrix) {
    auto df = femesh.derivative<1>();
    auto det = femesh.determinant<1>();

    Eigen::SparseMatrix<double> K = (df.dx(0) * det(0) * df.dx(0).transpose() + df.dy(0) * det(0) * df.dy(0).transpose()) * TriangleQuadratureRule<1>::w[0];

    EXPECT_NEAR(K.coeffRef(0, 0), 1.0, FLT_EPSILON);
    EXPECT_NEAR(K.coeffRef(1, 0), -1.0 / 2.0, FLT_EPSILON);
    EXPECT_NEAR(K.coeffRef(2, 0), -1.0 / 2.0, FLT_EPSILON);
    EXPECT_NEAR(K.coeffRef(0, 1), -1.0 / 2.0, FLT_EPSILON);
    EXPECT_NEAR(K.coeffRef(1, 1), 1.0 / 2.0, FLT_EPSILON);
    EXPECT_NEAR(K.coeffRef(2, 1), 0.0, FLT_EPSILON);
    EXPECT_NEAR(K.coeffRef(0, 2), -1.0 / 2.0, FLT_EPSILON);
    EXPECT_NEAR(K.coeffRef(1, 2), 0.0, FLT_EPSILON);
    EXPECT_NEAR(K.coeffRef(2, 2), 1.0 / 2.0, FLT_EPSILON);
}

TEST_F(MasterTriangle, magnetostatic) {
JasonPries's avatar
JasonPries committed
70
    //TODO: This isn't really testing anything...
JasonPries's avatar
JasonPries committed
71
72
73
    static constexpr size_t Q = 1;

    Magnetostatic<2, 1, Q, Variable::A> ms(femesh);
JasonPries's avatar
JasonPries committed
74
    ms.time(1.0);
JasonPries's avatar
JasonPries committed
75

JasonPries's avatar
JasonPries committed
76
    auto ff = [](double t) { return 1.0 * t; };
JasonPries's avatar
JasonPries committed
77

JasonPries's avatar
JasonPries committed
78
    ms.add_current_density(ff, std::vector<size_t>{0});
79
80
81

    ms.build_matrices();

JasonPries's avatar
JasonPries committed
82
83
84
85
86
87
88
89
90
    auto J = ms.init_unknown_matrix();
    auto v = ms.init_unknown_vector();
    auto r = ms.init_unknown_vector();
    auto f = ms.init_unknown_vector();
    auto Fx = ms.init_element_vector();
    auto Fy = ms.init_element_vector();
    auto dFxdx = ms.init_element_matrix();
    auto dFydy = ms.init_element_matrix();
    auto dFxdy = ms.init_element_matrix();
91

JasonPries's avatar
JasonPries committed
92
    ms.calculate_forcing(f);
JasonPries's avatar
JasonPries committed
93
    ms.linearize(J, v, r, f, Fx, Fy, dFxdx, dFydy, dFxdy);
94

JasonPries's avatar
JasonPries committed
95
    //auto ff = [](double t) { return 1.0 * t; };
96

JasonPries's avatar
JasonPries committed
97
    //HomogeneousForcing<2,1,Q> hf{ff, std::vector<size_t>{0}, ms.basis(), ms.weights(), ms.domain()};
98

JasonPries's avatar
JasonPries committed
99
    std::cout << f << std::endl;
100
101
102
}

class MasterTriangleRotationReflection : public ::testing::Test {
103
104
105
106
107
108
109
110
public:
    virtual void SetUp() {
        node.push_back(XY{0.0, 0.0});
        node.push_back(XY{1.0, 0.0});
        node.push_back(XY{0.0, 1.0});
        node.push_back(XY{-1.0, 0.0});
        node.push_back(XY{0.0, -1.0});

111
112
113
114
115
116
117
        tri.push_back(Triangle<1>(0, 1_zu, 2_zu, 0_zu)); // Master
        tri.push_back(Triangle<1>(1, 2_zu, 3_zu, 0_zu)); // 90 degree
        tri.push_back(Triangle<1>(2, 3_zu, 4_zu, 0_zu)); // 180 degree
        tri.push_back(Triangle<1>(3, 4_zu, 1_zu, 0_zu)); // 270 degree
        tri.push_back(Triangle<1>(4, 3_zu, 2_zu, 0_zu)); // y-axis reflection
        tri.push_back(Triangle<1>(5, 1_zu, 4_zu, 0_zu)); // x-axis reflection
        tri.push_back(Triangle<1>(6, 4_zu, 3_zu, 0_zu)); // x=y reflection
JasonPries's avatar
JasonPries committed
118

JasonPries's avatar
JasonPries committed
119
        reg.push_back(Region<2>(std::vector<size_t>{0, 1, 2, 3})); // Triangles 4,5,6 are duplicates
JasonPries's avatar
JasonPries committed
120

JasonPries's avatar
JasonPries committed
121
        femesh = FiniteElementMesh<2, 1>(node, tri, reg, bnd);
122
123
124
125
    }

    std::vector<XY> node;
    std::vector<Triangle<1>> tri;
JasonPries's avatar
JasonPries committed
126
    std::vector<Region<2>> reg;
JasonPries's avatar
JasonPries committed
127
128
    std::vector<Boundary<2>> bnd;

JasonPries's avatar
JasonPries committed
129
130

    FiniteElementMesh<2, 1> femesh;
131
132
};

133
TEST_F(MasterTriangleRotationReflection, jacobian_1) {
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
    auto J = tri[0].jacobian<1>(node);

    EXPECT_NEAR(J(0, 0), 1.0, FLT_EPSILON);
    EXPECT_NEAR(J(0, 1), 0.0, FLT_EPSILON);
    EXPECT_NEAR(J(1, 0), 0.0, FLT_EPSILON);
    EXPECT_NEAR(J(1, 1), 1.0, FLT_EPSILON);

    J = tri[1].jacobian<1>(node);

    EXPECT_NEAR(J(0, 0), 0.0, FLT_EPSILON);
    EXPECT_NEAR(J(0, 1), -1.0, FLT_EPSILON);
    EXPECT_NEAR(J(1, 0), 1.0, FLT_EPSILON);
    EXPECT_NEAR(J(1, 1), 0.0, FLT_EPSILON);

    J = tri[2].jacobian<1>(node);

    EXPECT_NEAR(J(0, 0), -1.0, FLT_EPSILON);
    EXPECT_NEAR(J(0, 1), 0.0, FLT_EPSILON);
    EXPECT_NEAR(J(1, 0), 0.0, FLT_EPSILON);
    EXPECT_NEAR(J(1, 1), -1.0, FLT_EPSILON);

    J = tri[3].jacobian<1>(node);

    EXPECT_NEAR(J(0, 0), 0.0, FLT_EPSILON);
    EXPECT_NEAR(J(0, 1), 1.0, FLT_EPSILON);
    EXPECT_NEAR(J(1, 0), -1.0, FLT_EPSILON);
    EXPECT_NEAR(J(1, 1), 0.0, FLT_EPSILON);

    J = tri[4].jacobian<1>(node);

    EXPECT_NEAR(J(0, 0), -1.0, FLT_EPSILON);
    EXPECT_NEAR(J(0, 1), 0.0, FLT_EPSILON);
    EXPECT_NEAR(J(1, 0), 0.0, FLT_EPSILON);
    EXPECT_NEAR(J(1, 1), 1.0, FLT_EPSILON);

    J = tri[5].jacobian<1>(node);

    EXPECT_NEAR(J(0, 0), 1.0, FLT_EPSILON);
    EXPECT_NEAR(J(0, 1), 0.0, FLT_EPSILON);
    EXPECT_NEAR(J(1, 0), 0.0, FLT_EPSILON);
    EXPECT_NEAR(J(1, 1), -1.0, FLT_EPSILON);

    J = tri[6].jacobian<1>(node);

    EXPECT_NEAR(J(0, 0), 0.0, FLT_EPSILON);
    EXPECT_NEAR(J(0, 1), -1.0, FLT_EPSILON);
    EXPECT_NEAR(J(1, 0), -1.0, FLT_EPSILON);
    EXPECT_NEAR(J(1, 1), 0.0, FLT_EPSILON);
}

184
TEST_F(MasterTriangleRotationReflection, basis_1) {
JasonPries's avatar
JasonPries committed
185
    auto mat = femesh.basis<1>();
186
187

    for (size_t i = 0; i != tri.size(); ++i) {
JasonPries's avatar
JasonPries committed
188
189
190
        EXPECT_NEAR(mat(0, tri[i].node(0), i), 1.0 / 3.0, FLT_EPSILON);
        EXPECT_NEAR(mat(0, tri[i].node(1), i), 1.0 / 3.0, FLT_EPSILON);
        EXPECT_NEAR(mat(0, tri[i].node(2), i), 1.0 / 3.0, FLT_EPSILON);
191
192
193
    }
}

194
TEST_F(MasterTriangleRotationReflection, derivative_1) {
JasonPries's avatar
JasonPries committed
195
    auto df = femesh.derivative<1>();
196

JasonPries's avatar
JasonPries committed
197
198
199
200
201
202
    EXPECT_NEAR(df.dx(0, tri[0].node(0), 0), 1.0, FLT_EPSILON);
    EXPECT_NEAR(df.dy(0, tri[0].node(0), 0), 0.0, FLT_EPSILON);
    EXPECT_NEAR(df.dx(0, tri[0].node(1), 0), 0.0, FLT_EPSILON);
    EXPECT_NEAR(df.dy(0, tri[0].node(1), 0), 1.0, FLT_EPSILON);
    EXPECT_NEAR(df.dx(0, tri[0].node(2), 0), -1.0, FLT_EPSILON);
    EXPECT_NEAR(df.dy(0, tri[0].node(2), 0), -1.0, FLT_EPSILON);
203
204
}

205
TEST_F(MasterTriangleRotationReflection, basis_2) {
JasonPries's avatar
JasonPries committed
206
    auto mat = femesh.basis<2>();
207
208

    for (size_t i = 0; i != tri.size(); ++i) {
JasonPries's avatar
JasonPries committed
209
210
211
212
213
214
215
        EXPECT_NEAR(mat(0, tri[i].node(0), i), 2.0 / 3.0, FLT_EPSILON);
        EXPECT_NEAR(mat(0, tri[i].node(1), i), 1.0 / 6.0, FLT_EPSILON);
        EXPECT_NEAR(mat(0, tri[i].node(1), i), 1.0 / 6.0, FLT_EPSILON);

        EXPECT_NEAR(mat(1, tri[i].node(0), i), 1.0 / 6.0, FLT_EPSILON);
        EXPECT_NEAR(mat(1, tri[i].node(1), i), 2.0 / 3.0, FLT_EPSILON);
        EXPECT_NEAR(mat(1, tri[i].node(2), i), 1.0 / 6.0, FLT_EPSILON);
216

JasonPries's avatar
JasonPries committed
217
218
219
220
221
        EXPECT_NEAR(mat(2, tri[i].node(0), i), 1.0 / 6.0, FLT_EPSILON);
        EXPECT_NEAR(mat(2, tri[i].node(1), i), 1.0 / 6.0, FLT_EPSILON);
        EXPECT_NEAR(mat(2, tri[i].node(2), i), 2.0 / 3.0, FLT_EPSILON);
    }
}
222

223
TEST_F(MasterTriangleRotationReflection, derivative_2) {
JasonPries's avatar
JasonPries committed
224
    auto df = femesh.derivative<2>();
225

JasonPries's avatar
JasonPries committed
226
227
228
229
230
231
232
    for (size_t q = 0; q != 3; ++q) {
        EXPECT_NEAR(df.dx(q, tri[0].node(0), 0), 1.0, FLT_EPSILON);
        EXPECT_NEAR(df.dy(q, tri[0].node(0), 0), 0.0, FLT_EPSILON);
        EXPECT_NEAR(df.dx(q, tri[0].node(1), 0), 0.0, FLT_EPSILON);
        EXPECT_NEAR(df.dy(q, tri[0].node(1), 0), 1.0, FLT_EPSILON);
        EXPECT_NEAR(df.dx(q, tri[0].node(2), 0), -1.0, FLT_EPSILON);
        EXPECT_NEAR(df.dy(q, tri[0].node(2), 0), -1.0, FLT_EPSILON);
233
234
235
    }
}

236
class SimpleSquare : public ::testing::Test {
JasonPries's avatar
JasonPries committed
237
238
239
240
241
242
243
public:
    virtual void SetUp() {
        node.push_back(XY{0.0, 0.0});
        node.push_back(XY{1.0, 0.0});
        node.push_back(XY{0.0, 1.0});
        node.push_back(XY{1.0, 1.0});

244
245
        tri.push_back(Triangle<1>(0, 0_zu, 1_zu, 2_zu));
        tri.push_back(Triangle<1>(1, 1_zu, 3_zu, 2_zu));
246

JasonPries's avatar
JasonPries committed
247
248
249
        reg.push_back(Region<2>(std::vector<size_t>{0}));
        reg.push_back(Region<2>(std::vector<size_t>{1}));

JasonPries's avatar
JasonPries committed
250
251
252
253
254
        bnd.push_back(Boundary<2>(std::vector<size_t>{0, 1}));
        bnd.push_back(Boundary<2>(std::vector<size_t>{0, 2}));
        bnd.push_back(Boundary<2>(std::vector<size_t>{1, 2}));
        bnd.push_back(Boundary<2>(std::vector<size_t>{1, 3}));
        bnd.push_back(Boundary<2>(std::vector<size_t>{2, 3}));
JasonPries's avatar
JasonPries committed
255
256

        femesh = FiniteElementMesh<2, 1>(node, tri, reg, bnd);
JasonPries's avatar
JasonPries committed
257
258
259
260
261

        vx_m << 1.0, 0.0, 1.0, 0.0;
        vx_p << 0.0, 1.0, 0.0, 1.0;
        vy_m << 1.0, 1.0, 0.0, 0.0;
        vy_p << 0.0, 0.0, 1.0, 1.0;
262
263
    }

JasonPries's avatar
JasonPries committed
264
265
    std::vector<XY> node;
    std::vector<Triangle<1>> tri;
JasonPries's avatar
JasonPries committed
266
    std::vector<Region<2>> reg;
JasonPries's avatar
JasonPries committed
267
    std::vector<Boundary<2>> bnd;
JasonPries's avatar
JasonPries committed
268
269
270
271
272
273
274
275
276

    Eigen::Vector4d vx_m;
    Eigen::Vector4d vx_p;
    Eigen::Vector4d vy_m;
    Eigen::Vector4d vy_p;

    FiniteElementMesh<2, 1> femesh;
};

277
TEST_F(SimpleSquare, derivative_1_temp) {
JasonPries's avatar
JasonPries committed
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
    auto df = femesh.derivative<1>();

    Eigen::Vector2d tmp;

    tmp = df.dx(0).transpose() * vx_p;
    EXPECT_NEAR(tmp(0), 1.0, FLT_EPSILON);
    EXPECT_NEAR(tmp(0), 1.0, FLT_EPSILON);

    tmp = df.dx(0).transpose() * vx_m;
    EXPECT_NEAR(tmp(0), -1.0, FLT_EPSILON);
    EXPECT_NEAR(tmp(0), -1.0, FLT_EPSILON);

    tmp = df.dx(0).transpose() * vy_m;
    EXPECT_NEAR(tmp(0), 0.0, FLT_EPSILON);
    EXPECT_NEAR(tmp(0), 0.0, FLT_EPSILON);

    tmp = df.dx(0).transpose() * vy_m;
    EXPECT_NEAR(tmp(0), 0.0, FLT_EPSILON);
    EXPECT_NEAR(tmp(0), 0.0, FLT_EPSILON);

    tmp = df.dy(0).transpose() * vy_p;
    EXPECT_NEAR(tmp(0), 1.0, FLT_EPSILON);
    EXPECT_NEAR(tmp(0), 1.0, FLT_EPSILON);

    tmp = df.dy(0).transpose() * vy_m;
    EXPECT_NEAR(tmp(0), -1.0, FLT_EPSILON);
    EXPECT_NEAR(tmp(0), -1.0, FLT_EPSILON);

    tmp = df.dy(0).transpose() * vx_m;
    EXPECT_NEAR(tmp(0), 0.0, FLT_EPSILON);
    EXPECT_NEAR(tmp(0), 0.0, FLT_EPSILON);

    tmp = df.dy(0).transpose() * vx_m;
    EXPECT_NEAR(tmp(0), 0.0, FLT_EPSILON);
    EXPECT_NEAR(tmp(0), 0.0, FLT_EPSILON);
}

315
TEST_F(SimpleSquare, derivative_2_temp) {
JasonPries's avatar
JasonPries committed
316
317
318
    auto df = femesh.derivative<2>();

    Eigen::Vector2d tmp;
319
320

    for (size_t i = 0; i != 3; ++i) {
JasonPries's avatar
JasonPries committed
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
        tmp = df.dx(i).transpose() * vx_p;
        EXPECT_NEAR(tmp(0), 1.0, FLT_EPSILON);
        EXPECT_NEAR(tmp(0), 1.0, FLT_EPSILON);

        tmp = df.dx(i).transpose() * vx_m;
        EXPECT_NEAR(tmp(0), -1.0, FLT_EPSILON);
        EXPECT_NEAR(tmp(0), -1.0, FLT_EPSILON);

        tmp = df.dx(i).transpose() * vy_m;
        EXPECT_NEAR(tmp(0), 0.0, FLT_EPSILON);
        EXPECT_NEAR(tmp(0), 0.0, FLT_EPSILON);

        tmp = df.dx(i).transpose() * vy_m;
        EXPECT_NEAR(tmp(0), 0.0, FLT_EPSILON);
        EXPECT_NEAR(tmp(0), 0.0, FLT_EPSILON);

        tmp = df.dy(i).transpose() * vy_p;
        EXPECT_NEAR(tmp(0), 1.0, FLT_EPSILON);
        EXPECT_NEAR(tmp(0), 1.0, FLT_EPSILON);

        tmp = df.dy(i).transpose() * vy_m;
        EXPECT_NEAR(tmp(0), -1.0, FLT_EPSILON);
        EXPECT_NEAR(tmp(0), -1.0, FLT_EPSILON);

        tmp = df.dy(i).transpose() * vx_m;
        EXPECT_NEAR(tmp(0), 0.0, FLT_EPSILON);
        EXPECT_NEAR(tmp(0), 0.0, FLT_EPSILON);

        tmp = df.dy(i).transpose() * vx_m;
        EXPECT_NEAR(tmp(0), 0.0, FLT_EPSILON);
        EXPECT_NEAR(tmp(0), 0.0, FLT_EPSILON);
352
    }
JasonPries's avatar
JasonPries committed
353
354
355
356
357
358
359
360
361
362
363
}

TEST_F(SimpleSquare, magnetostatic_forcing) {
    // TODO: Write Physics::add_condition() (templated?)
    // TODO: BoundaryCondition, ForcingCondition, MotionCondition?
    // TODO: ?then template parameters of conditions can be hidden?
    // TODO: Hierarchy of conditions: Those altering RHS only (forcing, some bcs), those altering matrices in a fixed way (some bcs), those altering the mesh (motion, some bcs)
    //      TODO: Applied in reverse order, Motion->BoundaryCondition->Forcing,
    //      TODO: Or, TimeVarying->TimeInvariant->Forcing?

    Magnetostatic<2, 1, 2, Variable::A> ms{femesh};
JasonPries's avatar
JasonPries committed
364
    ms.time(1.0);
JasonPries's avatar
JasonPries committed
365
366
367
    ms.build_matrices();

    auto ff = [](double t) { return 1.0 * t; };
JasonPries's avatar
JasonPries committed
368
    auto f = ms.init_unknown_vector();
JasonPries's avatar
JasonPries committed
369
370
371

    ms.add_current_density(ff, std::vector<size_t>{0});
    ms.calculate_forcing(f);
JasonPries's avatar
JasonPries committed
372
373
374
375
    EXPECT_NEAR(f(0), 1.0 / 6.0, FLT_EPSILON);
    EXPECT_NEAR(f(1), 1.0 / 6.0, FLT_EPSILON);
    EXPECT_NEAR(f(2), 1.0 / 6.0, FLT_EPSILON);
    EXPECT_NEAR(f(3), 0.0 / 6.0, FLT_EPSILON);
JasonPries's avatar
JasonPries committed
376
377
378
379
    ms.remove_current_density(0);

    ms.add_current_density(ff, std::vector<size_t>{1});
    ms.calculate_forcing(f);
JasonPries's avatar
JasonPries committed
380
381
382
383
    EXPECT_NEAR(f(0), 0.0 / 6.0, FLT_EPSILON);
    EXPECT_NEAR(f(1), 1.0 / 6.0, FLT_EPSILON);
    EXPECT_NEAR(f(2), 1.0 / 6.0, FLT_EPSILON);
    EXPECT_NEAR(f(3), 1.0 / 6.0, FLT_EPSILON);
JasonPries's avatar
JasonPries committed
384
385
    ms.remove_current_density(0);

JasonPries's avatar
JasonPries committed
386
    ms.add_current_density(ff, std::vector<size_t>{0, 1});
JasonPries's avatar
JasonPries committed
387
    ms.calculate_forcing(f);
JasonPries's avatar
JasonPries committed
388
389
390
391
    EXPECT_NEAR(f(0), 1.0 / 6.0, FLT_EPSILON);
    EXPECT_NEAR(f(1), 2.0 / 6.0, FLT_EPSILON);
    EXPECT_NEAR(f(2), 2.0 / 6.0, FLT_EPSILON);
    EXPECT_NEAR(f(3), 1.0 / 6.0, FLT_EPSILON);
JasonPries's avatar
JasonPries committed
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
    ms.remove_current_density(0);

    /*
    auto J = ms.init_matrix();
    auto v = ms.init_vector();
    auto r = ms.init_vector();
    v.setZero();

    ms.linearize(J,v,r,f);

    //Eigen::SimplicialLDLT<Eigen::SparseMatrix<double>> solver;
    //solver.compute(J);
    //if(solver.info()!=Eigen::Success) {
    //    std::cout << "Failed" << std::endl;
    //} else {
    //    std::cout << "Success" << std::endl;
    //}

    Eigen::SparseMatrix<double> J3 = J.block(1,1,3,3);
    Eigen::SparseMatrix<double> J2 = J.block(1,0,3,1);
    Eigen::SparseMatrix<double> J1 = J.block(0,1,1,3);
    Eigen::SparseMatrix<double> J0 = J.block(0,0,1,1);

    Eigen::SimplicialLDLT<Eigen::SparseMatrix<double>> solver;
    Eigen::SparseMatrix<double> Jb = J3 - (J2 * (J0 \ J1));

    std::cout << J3 << std::endl;
    std::cout << J2 << std::endl;
    std::cout << J1 << std::endl;
    std::cout << J0 << std::endl;
    std::cout << Jb << std::endl;

    Eigen::ConjugateGradient<Eigen::SparseMatrix<double>, Eigen::Lower|Eigen::Upper> cg;
    cg.compute(J);
    if(cg.info()!=Eigen::Success) {
        std::cout << "Failed" << std::endl;
    } else {
        std::cout << "Success" << std::endl;
    }
    v -= cg.solve(r);

    std::cout << J << std::endl;
    std::cout << v << std::endl;
    std::cout << r << std::endl;
    std::cout << f << std::endl;
    std::cout << J * v - r << std::endl;


    ms.linearize(J,v,r,f);
    std::cout << J << std::endl;
    std::cout << v << std::endl;
    std::cout << r << std::endl;
    std::cout << f << std::endl;

    std::cout << ms.derivatives().dx(0).transpose() * v << std::endl;
    std::cout << ms.derivatives().dy(0).transpose() * v << std::endl;
    */
}

class TwoRegionHex : public ::testing::Test {
public:
    virtual void SetUp() {
        node.push_back(XY{0.0, 0.0});

        double r{1.0};
        for (size_t i = 0; i != 6; ++i) {
            double a = (i * M_PI) / 3.0;
            node.push_back(XY{r * std::cos(a), r * std::sin(a)});
        }

        r = 2.0;
        for (size_t i = 0; i != 6; ++i) {
            double a = (i * M_PI) / 3.0 + M_PI / 6.0;
            node.push_back(XY{r * std::cos(a), r * std::sin(a)});
        }

        tri.push_back(Triangle<1>(0, 0_zu, 1_zu, 2_zu));
        tri.push_back(Triangle<1>(1, 0_zu, 2_zu, 3_zu));
        tri.push_back(Triangle<1>(2, 0_zu, 3_zu, 4_zu));
        tri.push_back(Triangle<1>(3, 0_zu, 4_zu, 5_zu));
        tri.push_back(Triangle<1>(4, 0_zu, 5_zu, 6_zu));
        tri.push_back(Triangle<1>(5, 0_zu, 6_zu, 1_zu));

JasonPries's avatar
JasonPries committed
475
        reg.push_back(Region<2>(std::vector<size_t>{0, 1, 2, 3, 4, 5}));
JasonPries's avatar
JasonPries committed
476

JasonPries's avatar
JasonPries committed
477
        bnd.push_back(Boundary<2>(std::vector<size_t>{1, 2, 3, 4, 5, 6}));
JasonPries's avatar
JasonPries committed
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496

        tri.push_back(Triangle<1>(6, 1_zu, 7_zu, 2_zu));
        tri.push_back(Triangle<1>(7, 2_zu, 7_zu, 8_zu));

        tri.push_back(Triangle<1>(8, 2_zu, 8_zu, 3_zu));
        tri.push_back(Triangle<1>(9, 3_zu, 8_zu, 9_zu));

        tri.push_back(Triangle<1>(10, 3_zu, 9_zu, 4_zu));
        tri.push_back(Triangle<1>(11, 4_zu, 9_zu, 10_zu));

        tri.push_back(Triangle<1>(12, 4_zu, 10_zu, 5_zu));
        tri.push_back(Triangle<1>(13, 5_zu, 10_zu, 11_zu));

        tri.push_back(Triangle<1>(14, 5_zu, 11_zu, 6_zu));
        tri.push_back(Triangle<1>(15, 6_zu, 11_zu, 12_zu));

        tri.push_back(Triangle<1>(16, 6_zu, 12_zu, 1_zu));
        tri.push_back(Triangle<1>(17, 1_zu, 12_zu, 7_zu));

JasonPries's avatar
JasonPries committed
497
        reg.push_back(Region<2>(std::vector<size_t>{6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17}));
JasonPries's avatar
JasonPries committed
498

JasonPries's avatar
JasonPries committed
499
        bnd.push_back(Boundary<2>(std::vector<size_t>{7, 8, 9, 10, 11, 12}));
JasonPries's avatar
JasonPries committed
500
501
502
503
504
505
506
507
508
509
510
511
512

        femesh = FiniteElementMesh<2, 1>(node, tri, reg, bnd);
    }

    std::vector<XY> node;
    std::vector<Triangle<1>> tri;
    std::vector<Region<2>> reg;
    std::vector<Boundary<2>> bnd;

    FiniteElementMesh<2, 1> femesh;
};

TEST_F(TwoRegionHex, magnetostatic_forcing) {
JasonPries's avatar
JasonPries committed
513
    Magnetostatic<2, 1, 1, Variable::A> ms{femesh};
JasonPries's avatar
JasonPries committed
514
    ms.time(0.0);
JasonPries's avatar
JasonPries committed
515
516

    // Test without boundary conditions
JasonPries's avatar
JasonPries committed
517
    ms.build_matrices();
JasonPries's avatar
JasonPries committed
518

JasonPries's avatar
JasonPries committed
519
    auto ff = [](double t) { return 1.0; };
JasonPries's avatar
JasonPries committed
520
    auto f = ms.init_unknown_vector();
JasonPries's avatar
JasonPries committed
521
522
523
524
525
526
527

    ms.add_current_density(ff, std::vector<size_t>{0});
    ms.add_magnetic_insulation(std::vector<size_t>{1});
    ms.calculate_forcing(f);

    double nodal_current{sqrt(3.0) / 4.0 / 3.0};
    EXPECT_NEAR(f(0), nodal_current * 6.0, FLT_EPSILON);
JasonPries's avatar
JasonPries committed
528
529
530
    for (size_t i = 1; i != 7; ++i) {
        EXPECT_NEAR(f(i), nodal_current * 2.0, FLT_EPSILON);
    }
JasonPries's avatar
JasonPries committed
531

JasonPries's avatar
JasonPries committed
532
533
534
535
536
537
538
539
540
    auto J = ms.init_unknown_matrix();
    auto v = ms.init_unknown_vector();
    auto r = ms.init_unknown_vector();
    auto Fx = ms.init_element_vector();
    auto Fy = ms.init_element_vector();
    auto dFxdx = ms.init_element_matrix();
    auto dFydy = ms.init_element_matrix();
    auto dFxdy = ms.init_element_matrix();

JasonPries's avatar
JasonPries committed
541
542
    v.setZero();

JasonPries's avatar
JasonPries committed
543
    ms.linearize(J, v, r, f, Fx, Fy, dFxdx, dFydy, dFxdy);
JasonPries's avatar
JasonPries committed
544
545
546
547

    //std::cout << J << std::endl;
    //std::cout << r << std::endl;

JasonPries's avatar
JasonPries committed
548
    // Test with boundary conditions
JasonPries's avatar
JasonPries committed
549
550
    ms.apply_conditions();

JasonPries's avatar
JasonPries committed
551
552
553
554
555
556
557
558
559
560
    J = ms.init_unknown_matrix();
    v = ms.init_unknown_vector();
    r = ms.init_unknown_vector();
    f = ms.init_unknown_vector();
    Fx = ms.init_element_vector();
    Fy = ms.init_element_vector();
    dFxdx = ms.init_element_matrix();
    dFydy = ms.init_element_matrix();
    dFxdy = ms.init_element_matrix();

JasonPries's avatar
JasonPries committed
561
562
563
564
    v.setZero();

    ms.calculate_forcing(f);
    EXPECT_NEAR(f(0), nodal_current * 6.0, FLT_EPSILON);
JasonPries's avatar
JasonPries committed
565
566
567
    for (size_t i = 1; i != 7; ++i) {
        EXPECT_NEAR(f(i), nodal_current * 2.0, FLT_EPSILON);
    }
JasonPries's avatar
JasonPries committed
568

JasonPries's avatar
JasonPries committed
569
    ms.linearize(J, v, r, f, Fx, Fy, dFxdx, dFydy, dFxdy);
JasonPries's avatar
JasonPries committed
570
571
572
573
574
    Eigen::SimplicialLDLT<Eigen::SparseMatrix<double>> ldlt;
    ldlt.compute(J);
    ASSERT_EQ(ldlt.info(), Eigen::Success);
    v -= ldlt.solve(r);

JasonPries's avatar
JasonPries committed
575
    // Verify equation is solved
JasonPries's avatar
JasonPries committed
576
    ms.linearize(J, v, r, f, Fx, Fy, dFxdx, dFydy, dFxdy);
JasonPries's avatar
JasonPries committed
577
578
579
580
581
582
583
584
585
    for (size_t i = 0; i != r.size(); ++i) {
        EXPECT_NEAR(r(i), 0.0, FLT_EPSILON);
    }

    // Fundamental solution is 0.25*r^2 + C in forced region (0)
    EXPECT_NEAR(v(0), v(1) + 0.25, FLT_EPSILON);
    for (size_t i = 1; i != 7; ++i) {
        EXPECT_NEAR(v(1), v(i), FLT_EPSILON);
    }
JasonPries's avatar
JasonPries committed
586

JasonPries's avatar
JasonPries committed
587
588
589
590
591
    // Test flux-density values
    // Values are exact in forced region (0)
    // Values are approximated in unforced region (1)
    Eigen::VectorXd Bx = ms.derivatives().dy(0).transpose() * v;
    Eigen::VectorXd By = -ms.derivatives().dx(0).transpose() * v;
JasonPries's avatar
JasonPries committed
592

JasonPries's avatar
JasonPries committed
593
594
    Eigen::VectorXd Bmag(By.size());
    Eigen::VectorXd Bang(By.size());
JasonPries's avatar
JasonPries committed
595

JasonPries's avatar
JasonPries committed
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
    for (size_t i = 0; i != By.size(); ++i) {
        Bmag(i) = hypot(Bx(i), By(i));
        Bang(i) = atan2(By(i), Bx(i)) * 180.0 / M_PI;
    }

    double x_q = (1.0 + cos(M_PI / 3.0)) / 3.0;
    double y_q = sin(M_PI / 3.0) / 3.0;
    double r_q = hypot(x_q, y_q);
    double B_q = r_q / 2.0;

    for (size_t i = 0; i != 6; ++i) {
        EXPECT_NEAR(Bmag(i), B_q, FLT_EPSILON);

        double ang = 120.0 + 60.0 * i;
        ang -= 360.0 * (ang > 180.0);

        EXPECT_NEAR(Bang(i), ang, FLT_EPSILON);
    }

    x_q = (1.0 + cos(M_PI / 3.0)) / 2.0;
    y_q = sin(M_PI / 3.0) / 2.0;
    r_q = 2.0 - hypot(x_q, y_q);
    B_q = v(1) / r_q;

    for (size_t i = 6; i < 18; i += 2) {
        EXPECT_NEAR(Bmag(i), B_q, FLT_EPSILON);

        double ang = 120.0 + 60.0 * (i - 6) / 2.0;
        ang -= 360.0 * (ang > 180.0);

        EXPECT_NEAR(Bang(i), ang, FLT_EPSILON);
    }

    x_q = (2.0 * cos(M_PI / 6.0) + 2.0 * cos(M_PI / 6.0 * 3.0)) / 2.0;
    y_q = (2.0 * sin(M_PI / 6.0) + 2.0 * sin(M_PI / 6.0 * 3.0)) / 2.0;
    r_q = hypot(x_q, y_q) - 1.0;
    B_q = v(1) / r_q;

    for (size_t i = 7; i < 18; i += 2) {
        EXPECT_NEAR(Bmag(i), B_q, FLT_EPSILON);

        double ang = 120.0 + 60.0 * (i - 6) / 2.0;
        ang -= 360.0 * (ang > 180.0);

        EXPECT_NEAR(Bang(i), ang, FLT_EPSILON);
    }
642
}