-
Notifications
You must be signed in to change notification settings - Fork 292
/
Copy pathcell_prism6.C
648 lines (534 loc) · 16 KB
/
cell_prism6.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
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
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
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
315
316
317
318
319
320
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
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
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
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
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
642
643
644
645
646
647
648
// The libMesh Finite Element Library.
// Copyright (C) 2002-2025 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
// This library is free software; you can redistribute it and/or
// modify it under the terms of the GNU Lesser General Public
// License as published by the Free Software Foundation; either
// version 2.1 of the License, or (at your option) any later version.
// This library is distributed in the hope that it will be useful,
// but WITHOUT ANY WARRANTY; without even the implied warranty of
// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
// Lesser General Public License for more details.
// You should have received a copy of the GNU Lesser General Public
// License along with this library; if not, write to the Free Software
// Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
// Local includes
#include "libmesh/side.h"
#include "libmesh/cell_prism6.h"
#include "libmesh/edge_edge2.h"
#include "libmesh/face_quad4.h"
#include "libmesh/face_tri3.h"
#include "libmesh/enum_io_package.h"
#include "libmesh/enum_order.h"
#include "libmesh/fe_lagrange_shape_1D.h"
// C++ includes
#include <array>
// Helper functions for computing volume() and true_centroid()
namespace {
using namespace libMesh;
// Templated numerical quadrature implementation function used by both
// Prism6::volume() and Prism6::true_centroid() routines.
template <typename T>
void cell_integral(const Elem * elem, T & t);
/**
* Object passed to cell_integral() routine for computing the cell volume.
*/
struct VolumeIntegral
{
// API
void accumulate_qp(Real /*xi*/, Real /*eta*/, Real /*zeta*/, Real JxW)
{
vol += JxW;
};
Real vol = 0.;
};
/**
* Object passed to cell_integral() routine for computing the nodal
* volumes (one value per node)
*/
struct NodalVolumeIntegral
{
// API
void accumulate_qp(Real xi, Real eta, Real zeta, Real JxW)
{
// Copied from fe_lagrange_shape_3D.C
static const unsigned int i0[] = {0, 0, 0, 1, 1, 1};
static const unsigned int i1[] = {0, 1, 2, 0, 1, 2};
// Copied from fe_lagrange_shape_2D.C
auto tri3_phi = [](int i, Real x, Real y)
{
switch(i)
{
case 0:
return 1 - x - y;
case 1:
return x;
case 2:
return y;
default:
libmesh_error_msg("Invalid shape function index i = " << i);
}
};
for (int i=0; i<6; ++i)
{
Real phi =
tri3_phi(i1[i], xi, eta) * fe_lagrange_1D_linear_shape(i0[i], zeta);
V[i] += phi * JxW;
}
};
// nodal volumes
std::array<Real, 6> V {};
};
}
namespace libMesh
{
// ------------------------------------------------------------
// Prism6 class static member initializations
const int Prism6::num_nodes;
const int Prism6::nodes_per_side;
const int Prism6::nodes_per_edge;
const unsigned int Prism6::side_nodes_map[Prism6::num_sides][Prism6::nodes_per_side] =
{
{0, 2, 1, 99}, // Side 0
{0, 1, 4, 3}, // Side 1
{1, 2, 5, 4}, // Side 2
{2, 0, 3, 5}, // Side 3
{3, 4, 5, 99} // Side 4
};
const unsigned int Prism6::side_elems_map[Prism6::num_sides][Prism6::nodes_per_side] =
{
{0, 1, 2, 3}, // Side 0
{0, 1, 4, 5}, // Side 1
{1, 2, 5, 6}, // Side 2
{0, 2, 4, 6}, // Side 3
{4, 5, 6, 7} // Side 4
};
const unsigned int Prism6::edge_nodes_map[Prism6::num_edges][Prism6::nodes_per_edge] =
{
{0, 1}, // Edge 0
{1, 2}, // Edge 1
{0, 2}, // Edge 2
{0, 3}, // Edge 3
{1, 4}, // Edge 4
{2, 5}, // Edge 5
{3, 4}, // Edge 6
{4, 5}, // Edge 7
{3, 5} // Edge 8
};
// ------------------------------------------------------------
// Prism6 class member functions
bool Prism6::is_vertex(const unsigned int libmesh_dbg_var(n)) const
{
libmesh_assert_not_equal_to (n, invalid_uint);
return true;
}
bool Prism6::is_edge(const unsigned int) const
{
return false;
}
bool Prism6::is_face(const unsigned int) const
{
return false;
}
bool Prism6::is_node_on_side(const unsigned int n,
const unsigned int s) const
{
libmesh_assert_less (s, n_sides());
return std::find(std::begin(side_nodes_map[s]),
std::end(side_nodes_map[s]),
n) != std::end(side_nodes_map[s]);
}
std::vector<unsigned>
Prism6::nodes_on_side(const unsigned int s) const
{
libmesh_assert_less(s, n_sides());
auto trim = (s > 0 && s < 4) ? 0 : 1;
return {std::begin(side_nodes_map[s]), std::end(side_nodes_map[s]) - trim};
}
std::vector<unsigned>
Prism6::nodes_on_edge(const unsigned int e) const
{
libmesh_assert_less(e, n_edges());
return {std::begin(edge_nodes_map[e]), std::end(edge_nodes_map[e])};
}
bool Prism6::is_node_on_edge(const unsigned int n,
const unsigned int e) const
{
libmesh_assert_less (e, n_edges());
return std::find(std::begin(edge_nodes_map[e]),
std::end(edge_nodes_map[e]),
n) != std::end(edge_nodes_map[e]);
}
bool Prism6::has_affine_map() const
{
// Make sure z edges are affine
Point v = this->point(3) - this->point(0);
if (!v.relative_fuzzy_equals(this->point(4) - this->point(1), affine_tol) ||
!v.relative_fuzzy_equals(this->point(5) - this->point(2), affine_tol))
return false;
return true;
}
Order Prism6::default_order() const
{
return FIRST;
}
std::unique_ptr<Elem> Prism6::build_side_ptr (const unsigned int i,
bool proxy)
{
libmesh_assert_less (i, this->n_sides());
std::unique_ptr<Elem> face;
if (proxy)
{
#ifdef LIBMESH_ENABLE_DEPRECATED
libmesh_deprecated();
switch(i)
{
case 0:
case 4:
{
face = std::make_unique<Side<Tri3,Prism6>>(this,i);
break;
}
case 1:
case 2:
case 3:
{
face = std::make_unique<Side<Quad4,Prism6>>(this,i);
break;
}
default:
libmesh_error_msg("Invalid side i = " << i);
}
#else
libmesh_error();
#endif // LIBMESH_ENABLE_DEPRECATED
}
else
{
switch (i)
{
case 0: // the triangular face at z=-1
case 4: // the triangular face at z=1
{
face = std::make_unique<Tri3>();
break;
}
case 1: // the quad face at y=0
case 2: // the other quad face
case 3: // the quad face at x=0
{
face = std::make_unique<Quad4>();
break;
}
default:
libmesh_error_msg("Invalid side i = " << i);
}
// Set the nodes
for (auto n : face->node_index_range())
face->set_node(n) = this->node_ptr(Prism6::side_nodes_map[i][n]);
}
#ifdef LIBMESH_ENABLE_DEPRECATED
if (!proxy) // proxy sides used to leave parent() set
#endif
face->set_parent(nullptr);
face->set_interior_parent(this);
face->subdomain_id() = this->subdomain_id();
face->set_mapping_type(this->mapping_type());
#ifdef LIBMESH_ENABLE_AMR
face->set_p_level(this->p_level());
#endif
return face;
}
void Prism6::build_side_ptr (std::unique_ptr<Elem> & side,
const unsigned int i)
{
this->side_ptr(side, i);
side->set_interior_parent(this);
side->subdomain_id() = this->subdomain_id();
side->set_mapping_type(this->mapping_type());
#ifdef LIBMESH_ENABLE_AMR
side->set_p_level(this->p_level());
#endif
}
std::unique_ptr<Elem> Prism6::build_edge_ptr (const unsigned int i)
{
return this->simple_build_edge_ptr<Edge2,Prism6>(i);
}
void Prism6::build_edge_ptr (std::unique_ptr<Elem> & edge, const unsigned int i)
{
this->simple_build_edge_ptr<Prism6>(edge, i, EDGE2);
}
void Prism6::connectivity(const unsigned int libmesh_dbg_var(sc),
const IOPackage iop,
std::vector<dof_id_type> & conn) const
{
libmesh_assert(_nodes);
libmesh_assert_less (sc, this->n_sub_elem());
libmesh_assert_not_equal_to (iop, INVALID_IO_PACKAGE);
switch (iop)
{
case TECPLOT:
{
conn.resize(8);
conn[0] = this->node_id(0)+1;
conn[1] = this->node_id(1)+1;
conn[2] = this->node_id(2)+1;
conn[3] = this->node_id(2)+1;
conn[4] = this->node_id(3)+1;
conn[5] = this->node_id(4)+1;
conn[6] = this->node_id(5)+1;
conn[7] = this->node_id(5)+1;
return;
}
case VTK:
{
conn.resize(6);
conn[0] = this->node_id(0);
conn[1] = this->node_id(2);
conn[2] = this->node_id(1);
conn[3] = this->node_id(3);
conn[4] = this->node_id(5);
conn[5] = this->node_id(4);
return;
}
default:
libmesh_error_msg("Unsupported IO package " << iop);
}
}
#ifdef LIBMESH_ENABLE_AMR
const Real Prism6::_embedding_matrix[Prism6::num_children][Prism6::num_nodes][Prism6::num_nodes] =
{
// embedding matrix for child 0
{
// 0 1 2 3 4 5
{ 1.0, 0.0, 0.0, 0.0, 0.0, 0.0}, // 0
{ 0.5, 0.5, 0.0, 0.0, 0.0, 0.0}, // 1
{ 0.5, 0.0, 0.5, 0.0, 0.0, 0.0}, // 2
{ 0.5, 0.0, 0.0, 0.5, 0.0, 0.0}, // 3
{ .25, .25, 0.0, .25, .25, 0.0}, // 4
{ .25, 0.0, .25, .25, 0.0, .25} // 5
},
// embedding matrix for child 1
{
// 0 1 2 3 4 5
{ 0.5, 0.5, 0.0, 0.0, 0.0, 0.0}, // 0
{ 0.0, 1.0, 0.0, 0.0, 0.0, 0.0}, // 1
{ 0.0, 0.5, 0.5, 0.0, 0.0, 0.0}, // 2
{ .25, .25, 0.0, .25, .25, 0.0}, // 3
{ 0.0, 0.5, 0.0, 0.0, 0.5, 0.0}, // 4
{ 0.0, .25, .25, 0.0, .25, .25} // 5
},
// embedding matrix for child 2
{
// 0 1 2 3 4 5
{ 0.5, 0.0, 0.5, 0.0, 0.0, 0.0}, // 0
{ 0.0, 0.5, 0.5, 0.0, 0.0, 0.0}, // 1
{ 0.0, 0.0, 1.0, 0.0, 0.0, 0.0}, // 2
{ .25, 0.0, .25, .25, 0.0, .25}, // 3
{ 0.0, .25, .25, 0.0, .25, .25}, // 4
{ 0.0, 0.0, 0.5, 0.0, 0.0, 0.5} // 5
},
// embedding matrix for child 3
{
// 0 1 2 3 4 5
{ 0.5, 0.5, 0.0, 0.0, 0.0, 0.0}, // 0
{ 0.0, 0.5, 0.5, 0.0, 0.0, 0.0}, // 1
{ 0.5, 0.0, 0.5, 0.0, 0.0, 0.0}, // 2
{ .25, .25, 0.0, .25, .25, 0.0}, // 3
{ 0.0, .25, .25, 0.0, .25, .25}, // 4
{ .25, 0.0, .25, .25, 0.0, .25} // 5
},
// embedding matrix for child 4
{
// 0 1 2 3 4 5
{ 0.5, 0.0, 0.0, 0.5, 0.0, 0.0}, // 0
{ .25, .25, 0.0, .25, .25, 0.0}, // 1
{ .25, 0.0, .25, .25, 0.0, .25}, // 2
{ 0.0, 0.0, 0.0, 1.0, 0.0, 0.0}, // 3
{ 0.0, 0.0, 0.0, 0.5, 0.5, 0.0}, // 4
{ 0.0, 0.0, 0.0, 0.5, 0.0, 0.5} // 5
},
// embedding matrix for child 5
{
// 0 1 2 3 4 5
{ .25, .25, 0.0, .25, .25, 0.0}, // 0
{ 0.0, 0.5, 0.0, 0.0, 0.5, 0.0}, // 1
{ 0.0, .25, .25, 0.0, .25, .25}, // 2
{ 0.0, 0.0, 0.0, 0.5, 0.5, 0.0}, // 3
{ 0.0, 0.0, 0.0, 0.0, 1.0, 0.0}, // 4
{ 0.0, 0.0, 0.0, 0.0, 0.5, 0.5} // 5
},
// embedding matrix for child 6
{
// 0 1 2 3 4 5
{ .25, 0.0, .25, .25, 0.0, .25}, // 0
{ 0.0, .25, .25, 0.0, .25, .25}, // 1
{ 0.0, 0.0, 0.5, 0.0, 0.0, 0.5}, // 2
{ 0.0, 0.0, 0.0, 0.5, 0.0, 0.5}, // 3
{ 0.0, 0.0, 0.0, 0.0, 0.5, 0.5}, // 4
{ 0.0, 0.0, 0.0, 0.0, 0.0, 1.0} // 5
},
// embedding matrix for child 7
{
// 0 1 2 3 4 5
{ .25, .25, 0.0, .25, .25, 0.0}, // 0
{ 0.0, .25, .25, 0.0, .25, .25}, // 1
{ .25, 0.0, .25, .25, 0.0, .25}, // 2
{ 0.0, 0.0, 0.0, 0.5, 0.5, 0.0}, // 3
{ 0.0, 0.0, 0.0, 0.0, 0.5, 0.5}, // 4
{ 0.0, 0.0, 0.0, 0.5, 0.0, 0.5} // 5
}
};
#endif
Point Prism6::true_centroid () const
{
NodalVolumeIntegral nvi;
cell_integral(this, nvi);
// Compute centroid
Point cp;
Real vol = 0.;
for (int i=0; i<6; ++i)
{
cp += this->point(i) * nvi.V[i];
vol += nvi.V[i];
}
return cp / vol;
}
Real Prism6::volume () const
{
VolumeIntegral vi;
cell_integral(this, vi);
return vi.vol;
}
BoundingBox
Prism6::loose_bounding_box () const
{
return Elem::loose_bounding_box();
}
void
Prism6::permute(unsigned int perm_num)
{
libmesh_assert_less (perm_num, 6);
const unsigned int side = perm_num % 2;
const unsigned int rotate = perm_num / 2;
for (unsigned int i = 0; i != rotate; ++i)
{
swap3nodes(0,1,2);
swap3nodes(3,4,5);
swap3neighbors(1,2,3);
}
switch (side) {
case 0:
break;
case 1:
swap2nodes(1,3);
swap2nodes(0,4);
swap2nodes(2,5);
swap2neighbors(0,4);
swap2neighbors(2,3);
break;
default:
libmesh_error();
}
}
void
Prism6::flip(BoundaryInfo * boundary_info)
{
libmesh_assert(boundary_info);
swap2nodes(0,1);
swap2nodes(3,4);
swap2neighbors(2,3);
swap2boundarysides(2,3,boundary_info);
swap2boundaryedges(0,1,boundary_info);
swap2boundaryedges(3,4,boundary_info);
swap2boundaryedges(7,8,boundary_info);
}
ElemType
Prism6::side_type (const unsigned int s) const
{
libmesh_assert_less (s, 5);
if (s == 0 || s == 4)
return TRI3;
return QUAD4;
}
} // namespace libMesh
namespace // anonymous
{
template <typename T>
void cell_integral(const Elem * elem, T & t)
{
// Conveient references to our points.
const Point
&x0 = elem->point(0), &x1 = elem->point(1), &x2 = elem->point(2),
&x3 = elem->point(3), &x4 = elem->point(4), &x5 = elem->point(5);
// constant and zeta terms only. These are copied directly from a
// Python script.
Point dx_dxi[2] =
{
-x0/2 + x1/2 - x3/2 + x4/2, // constant
x0/2 - x1/2 - x3/2 + x4/2, // zeta
};
// constant and zeta terms only. These are copied directly from a
// Python script.
Point dx_deta[2] =
{
-x0/2 + x2/2 - x3/2 + x5/2, // constant
x0/2 - x2/2 - x3/2 + x5/2, // zeta
};
// Constant, xi, and eta terms
Point dx_dzeta[3] =
{
-x0/2 + x3/2, // constant
x0/2 - x2/2 - x3/2 + x5/2, // eta
x0/2 - x1/2 - x3/2 + x4/2 // xi
};
// The quadrature rule for the Prism6 is a tensor product between a
// four-point TRI3 rule (in xi, eta) and a two-point EDGE2 rule (in
// zeta) which is capable of integrating cubics exactly.
// Number of points in the 2D quadrature rule.
const int N2D = 4;
static const Real w2D[N2D] =
{
Real(1.5902069087198858469718450103758e-01L),
Real(9.0979309128011415302815498962418e-02L),
Real(1.5902069087198858469718450103758e-01L),
Real(9.0979309128011415302815498962418e-02L)
};
static const Real xi[N2D] =
{
Real(1.5505102572168219018027159252941e-01L),
Real(6.4494897427831780981972840747059e-01L),
Real(1.5505102572168219018027159252941e-01L),
Real(6.4494897427831780981972840747059e-01L)
};
static const Real eta[N2D] =
{
Real(1.7855872826361642311703513337422e-01L),
Real(7.5031110222608118177475598324603e-02L),
Real(6.6639024601470138670269327409637e-01L),
Real(2.8001991549907407200279599420481e-01L)
};
// Number of points in the 1D quadrature rule. The weights of the
// 1D quadrature rule are equal to 1.
const int N1D = 2;
// Points of the 1D quadrature rule
static const Real zeta[N1D] =
{
-std::sqrt(3.)/3,
std::sqrt(3.)/3.
};
for (int i=0; i<N2D; ++i)
{
// dx_dzeta depends only on the 2D quadrature rule points.
Point dx_dzeta_q = dx_dzeta[0] + eta[i]*dx_dzeta[1] + xi[i]*dx_dzeta[2];
for (int j=0; j<N1D; ++j)
{
// dx_dxi and dx_deta only depend on the 1D quadrature rule points.
Point
dx_dxi_q = dx_dxi[0] + zeta[j]*dx_dxi[1],
dx_deta_q = dx_deta[0] + zeta[j]*dx_deta[1];
// Compute t-specific contributions at current qp
t.accumulate_qp(xi[i], eta[i], zeta[j], w2D[i] * triple_product(dx_dxi_q, dx_deta_q, dx_dzeta_q));
}
}
}
}