-
Notifications
You must be signed in to change notification settings - Fork 708
/
mapping_q_internal.h
2490 lines (2214 loc) · 105 KB
/
mapping_q_internal.h
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
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
// ---------------------------------------------------------------------
//
// Copyright (C) 2020 - 2022 by the deal.II authors
//
// This file is part of the deal.II library.
//
// The deal.II library is free software; you can use it, 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.
// The full text of the license can be found in the file LICENSE.md at
// the top level directory of deal.II.
//
// ---------------------------------------------------------------------
#ifndef dealii_mapping_q_internal_h
#define dealii_mapping_q_internal_h
#include <deal.II/base/config.h>
#include <deal.II/base/array_view.h>
#include <deal.II/base/derivative_form.h>
#include <deal.II/base/geometry_info.h>
#include <deal.II/base/point.h>
#include <deal.II/base/qprojector.h>
#include <deal.II/base/table.h>
#include <deal.II/base/tensor.h>
#include <deal.II/fe/fe_tools.h>
#include <deal.II/fe/fe_update_flags.h>
#include <deal.II/fe/fe_values.h>
#include <deal.II/fe/mapping_q.h>
#include <deal.II/grid/grid_tools.h>
#include <deal.II/matrix_free/evaluation_flags.h>
#include <deal.II/matrix_free/evaluation_template_factory.h>
#include <deal.II/matrix_free/fe_evaluation_data.h>
#include <deal.II/matrix_free/mapping_info_storage.h>
#include <deal.II/matrix_free/shape_info.h>
#include <deal.II/matrix_free/tensor_product_kernels.h>
#include <array>
#include <limits>
DEAL_II_NAMESPACE_OPEN
namespace internal
{
/**
* Internal namespace to implement methods specific to MappingQ1, in
* particular an explicit formula for the transformation from the real to
* the unit cell in 2d.
*/
namespace MappingQ1
{
// These are left as templates on the spatial dimension (even though dim
// == spacedim must be true for them to make sense) because templates are
// expanded before the compiler eliminates code due to the 'if (dim ==
// spacedim)' statement (see the body of the general
// transform_real_to_unit_cell).
template <int spacedim>
inline Point<1>
transform_real_to_unit_cell(
const std::array<Point<spacedim>, GeometryInfo<1>::vertices_per_cell>
& vertices,
const Point<spacedim> &p)
{
Assert(spacedim == 1, ExcInternalError());
return Point<1>((p[0] - vertices[0](0)) /
(vertices[1](0) - vertices[0](0)));
}
template <int spacedim>
inline Point<2>
transform_real_to_unit_cell(
const std::array<Point<spacedim>, GeometryInfo<2>::vertices_per_cell>
& vertices,
const Point<spacedim> &p)
{
Assert(spacedim == 2, ExcInternalError());
// For accuracy reasons, we do all arithmetic in extended precision
// (long double). This has a noticeable effect on the hit rate for
// borderline cases and thus makes the algorithm more robust.
const long double x = p(0);
const long double y = p(1);
const long double x0 = vertices[0](0);
const long double x1 = vertices[1](0);
const long double x2 = vertices[2](0);
const long double x3 = vertices[3](0);
const long double y0 = vertices[0](1);
const long double y1 = vertices[1](1);
const long double y2 = vertices[2](1);
const long double y3 = vertices[3](1);
const long double a = (x1 - x3) * (y0 - y2) - (x0 - x2) * (y1 - y3);
const long double b = -(x0 - x1 - x2 + x3) * y + (x - 2 * x1 + x3) * y0 -
(x - 2 * x0 + x2) * y1 - (x - x1) * y2 +
(x - x0) * y3;
const long double c = (x0 - x1) * y - (x - x1) * y0 + (x - x0) * y1;
const long double discriminant = b * b - 4 * a * c;
// exit if the point is not in the cell (this is the only case where the
// discriminant is negative)
AssertThrow(
discriminant > 0.0,
(typename Mapping<spacedim, spacedim>::ExcTransformationFailed()));
long double eta1;
long double eta2;
const long double sqrt_discriminant = std::sqrt(discriminant);
// special case #1: if a is near-zero to make the discriminant exactly
// equal b, then use the linear formula
if (b != 0.0 && std::abs(b) == sqrt_discriminant)
{
eta1 = -c / b;
eta2 = -c / b;
}
// special case #2: a is zero for parallelograms and very small for
// near-parallelograms:
else if (std::abs(a) < 1e-8 * std::abs(b))
{
// if both a and c are very small then the root should be near
// zero: this first case will capture that
eta1 = 2 * c / (-b - sqrt_discriminant);
eta2 = 2 * c / (-b + sqrt_discriminant);
}
// finally, use the plain version:
else
{
eta1 = (-b - sqrt_discriminant) / (2 * a);
eta2 = (-b + sqrt_discriminant) / (2 * a);
}
// pick the one closer to the center of the cell.
const long double eta =
(std::abs(eta1 - 0.5) < std::abs(eta2 - 0.5)) ? eta1 : eta2;
/*
* There are two ways to compute xi from eta, but either one may have a
* zero denominator.
*/
const long double subexpr0 = -eta * x2 + x0 * (eta - 1);
const long double xi_denominator0 = eta * x3 - x1 * (eta - 1) + subexpr0;
const long double max_x = std::max(std::max(std::abs(x0), std::abs(x1)),
std::max(std::abs(x2), std::abs(x3)));
if (std::abs(xi_denominator0) > 1e-10 * max_x)
{
const double xi = (x + subexpr0) / xi_denominator0;
return {xi, static_cast<double>(eta)};
}
else
{
const long double max_y =
std::max(std::max(std::abs(y0), std::abs(y1)),
std::max(std::abs(y2), std::abs(y3)));
const long double subexpr1 = -eta * y2 + y0 * (eta - 1);
const long double xi_denominator1 =
eta * y3 - y1 * (eta - 1) + subexpr1;
if (std::abs(xi_denominator1) > 1e-10 * max_y)
{
const double xi = (subexpr1 + y) / xi_denominator1;
return {xi, static_cast<double>(eta)};
}
else // give up and try Newton iteration
{
AssertThrow(
false,
(typename Mapping<spacedim,
spacedim>::ExcTransformationFailed()));
}
}
// bogus return to placate compiler. It should not be possible to get
// here.
Assert(false, ExcInternalError());
return {std::numeric_limits<double>::quiet_NaN(),
std::numeric_limits<double>::quiet_NaN()};
}
template <int spacedim>
inline Point<3>
transform_real_to_unit_cell(
const std::array<Point<spacedim>, GeometryInfo<3>::vertices_per_cell>
& /*vertices*/,
const Point<spacedim> & /*p*/)
{
// It should not be possible to get here
Assert(false, ExcInternalError());
return {std::numeric_limits<double>::quiet_NaN(),
std::numeric_limits<double>::quiet_NaN(),
std::numeric_limits<double>::quiet_NaN()};
}
} // namespace MappingQ1
/**
* Internal namespace to implement methods of MappingQ, such as the
* evaluation of the mapping and the transformation between real and unit
* cell.
*/
namespace MappingQImplementation
{
/**
* This function generates the reference cell support points from the 1d
* support points by expanding the tensor product.
*/
template <int dim>
std::vector<Point<dim>>
unit_support_points(const std::vector<Point<1>> & line_support_points,
const std::vector<unsigned int> &renumbering)
{
AssertDimension(Utilities::pow(line_support_points.size(), dim),
renumbering.size());
std::vector<Point<dim>> points(renumbering.size());
const unsigned int n1 = line_support_points.size();
for (unsigned int q2 = 0, q = 0; q2 < (dim > 2 ? n1 : 1); ++q2)
for (unsigned int q1 = 0; q1 < (dim > 1 ? n1 : 1); ++q1)
for (unsigned int q0 = 0; q0 < n1; ++q0, ++q)
{
points[renumbering[q]][0] = line_support_points[q0][0];
if (dim > 1)
points[renumbering[q]][1] = line_support_points[q1][0];
if (dim > 2)
points[renumbering[q]][2] = line_support_points[q2][0];
}
return points;
}
/**
* This function is needed by the constructor of
* <tt>MappingQ<dim,spacedim></tt> for <tt>dim=</tt> 2 and 3.
*
* For the definition of the @p support_point_weights_on_quad please
* refer to the description of TransfiniteInterpolationManifold.
*/
inline dealii::Table<2, double>
compute_support_point_weights_on_quad(const unsigned int polynomial_degree)
{
dealii::Table<2, double> loqvs;
// we are asked to compute weights for interior support points, but
// there are no interior points if degree==1
if (polynomial_degree == 1)
return loqvs;
const unsigned int M = polynomial_degree - 1;
const unsigned int n_inner_2d = M * M;
const unsigned int n_outer_2d = 4 + 4 * M;
// set the weights of transfinite interpolation
loqvs.reinit(n_inner_2d, n_outer_2d);
QGaussLobatto<2> gl(polynomial_degree + 1);
for (unsigned int i = 0; i < M; ++i)
for (unsigned int j = 0; j < M; ++j)
{
const Point<2> &p =
gl.point((i + 1) * (polynomial_degree + 1) + (j + 1));
const unsigned int index_table = i * M + j;
for (unsigned int v = 0; v < 4; ++v)
loqvs(index_table, v) =
-GeometryInfo<2>::d_linear_shape_function(p, v);
loqvs(index_table, 4 + i) = 1. - p[0];
loqvs(index_table, 4 + i + M) = p[0];
loqvs(index_table, 4 + j + 2 * M) = 1. - p[1];
loqvs(index_table, 4 + j + 3 * M) = p[1];
}
// the sum of weights of the points at the outer rim should be one.
// check this
for (unsigned int unit_point = 0; unit_point < n_inner_2d; ++unit_point)
Assert(std::fabs(std::accumulate(loqvs[unit_point].begin(),
loqvs[unit_point].end(),
0.) -
1) < 1e-13 * polynomial_degree,
ExcInternalError());
return loqvs;
}
/**
* This function is needed by the constructor of <tt>MappingQ<3></tt>.
*
* For the definition of the @p support_point_weights_on_quad please
* refer to the description of TransfiniteInterpolationManifold.
*/
inline dealii::Table<2, double>
compute_support_point_weights_on_hex(const unsigned int polynomial_degree)
{
dealii::Table<2, double> lohvs;
// we are asked to compute weights for interior support points, but
// there are no interior points if degree==1
if (polynomial_degree == 1)
return lohvs;
const unsigned int M = polynomial_degree - 1;
const unsigned int n_inner = Utilities::fixed_power<3>(M);
const unsigned int n_outer = 8 + 12 * M + 6 * M * M;
// set the weights of transfinite interpolation
lohvs.reinit(n_inner, n_outer);
QGaussLobatto<3> gl(polynomial_degree + 1);
for (unsigned int i = 0; i < M; ++i)
for (unsigned int j = 0; j < M; ++j)
for (unsigned int k = 0; k < M; ++k)
{
const Point<3> & p = gl.point((i + 1) * (M + 2) * (M + 2) +
(j + 1) * (M + 2) + (k + 1));
const unsigned int index_table = i * M * M + j * M + k;
// vertices
for (unsigned int v = 0; v < 8; ++v)
lohvs(index_table, v) =
GeometryInfo<3>::d_linear_shape_function(p, v);
// lines
{
constexpr std::array<unsigned int, 4> line_coordinates_y(
{{0, 1, 4, 5}});
const Point<2> py(p[0], p[2]);
for (unsigned int l = 0; l < 4; ++l)
lohvs(index_table, 8 + line_coordinates_y[l] * M + j) =
-GeometryInfo<2>::d_linear_shape_function(py, l);
}
{
constexpr std::array<unsigned int, 4> line_coordinates_x(
{{2, 3, 6, 7}});
const Point<2> px(p[1], p[2]);
for (unsigned int l = 0; l < 4; ++l)
lohvs(index_table, 8 + line_coordinates_x[l] * M + k) =
-GeometryInfo<2>::d_linear_shape_function(px, l);
}
{
constexpr std::array<unsigned int, 4> line_coordinates_z(
{{8, 9, 10, 11}});
const Point<2> pz(p[0], p[1]);
for (unsigned int l = 0; l < 4; ++l)
lohvs(index_table, 8 + line_coordinates_z[l] * M + i) =
-GeometryInfo<2>::d_linear_shape_function(pz, l);
}
// quads
lohvs(index_table, 8 + 12 * M + 0 * M * M + i * M + j) =
1. - p[0];
lohvs(index_table, 8 + 12 * M + 1 * M * M + i * M + j) = p[0];
lohvs(index_table, 8 + 12 * M + 2 * M * M + k * M + i) =
1. - p[1];
lohvs(index_table, 8 + 12 * M + 3 * M * M + k * M + i) = p[1];
lohvs(index_table, 8 + 12 * M + 4 * M * M + j * M + k) =
1. - p[2];
lohvs(index_table, 8 + 12 * M + 5 * M * M + j * M + k) = p[2];
}
// the sum of weights of the points at the outer rim should be one.
// check this
for (unsigned int unit_point = 0; unit_point < n_inner; ++unit_point)
Assert(std::fabs(std::accumulate(lohvs[unit_point].begin(),
lohvs[unit_point].end(),
0.) -
1) < 1e-13 * polynomial_degree,
ExcInternalError());
return lohvs;
}
/**
* This function collects the output of
* compute_support_point_weights_on_{quad,hex} in a single data structure.
*/
inline std::vector<dealii::Table<2, double>>
compute_support_point_weights_perimeter_to_interior(
const unsigned int polynomial_degree,
const unsigned int dim)
{
Assert(dim > 0 && dim <= 3, ExcImpossibleInDim(dim));
std::vector<dealii::Table<2, double>> output(dim);
if (polynomial_degree <= 1)
return output;
// fill the 1d interior weights
QGaussLobatto<1> quadrature(polynomial_degree + 1);
output[0].reinit(polynomial_degree - 1,
GeometryInfo<1>::vertices_per_cell);
for (unsigned int q = 0; q < polynomial_degree - 1; ++q)
for (const unsigned int i : GeometryInfo<1>::vertex_indices())
output[0](q, i) =
GeometryInfo<1>::d_linear_shape_function(quadrature.point(q + 1),
i);
if (dim > 1)
output[1] = compute_support_point_weights_on_quad(polynomial_degree);
if (dim > 2)
output[2] = compute_support_point_weights_on_hex(polynomial_degree);
return output;
}
/**
* Collects all interior points for the various dimensions.
*/
template <int dim>
inline dealii::Table<2, double>
compute_support_point_weights_cell(const unsigned int polynomial_degree)
{
Assert(dim > 0 && dim <= 3, ExcImpossibleInDim(dim));
if (polynomial_degree <= 1)
return dealii::Table<2, double>();
QGaussLobatto<dim> quadrature(polynomial_degree + 1);
const std::vector<unsigned int> h2l =
FETools::hierarchic_to_lexicographic_numbering<dim>(polynomial_degree);
dealii::Table<2, double> output(quadrature.size() -
GeometryInfo<dim>::vertices_per_cell,
GeometryInfo<dim>::vertices_per_cell);
for (unsigned int q = 0; q < output.size(0); ++q)
for (const unsigned int i : GeometryInfo<dim>::vertex_indices())
output(q, i) = GeometryInfo<dim>::d_linear_shape_function(
quadrature.point(h2l[q + GeometryInfo<dim>::vertices_per_cell]), i);
return output;
}
/**
* Using the relative weights of the shape functions evaluated at
* one point on the reference cell (and stored in data.shape_values
* and accessed via data.shape(0,i)) and the locations of mapping
* support points (stored in data.mapping_support_points), compute
* the mapped location of that point in real space.
*/
template <int dim, int spacedim>
inline Point<spacedim>
compute_mapped_location_of_point(
const typename dealii::MappingQ<dim, spacedim>::InternalData &data)
{
AssertDimension(data.shape_values.size(),
data.mapping_support_points.size());
// use now the InternalData to compute the point in real space.
Point<spacedim> p_real;
for (unsigned int i = 0; i < data.mapping_support_points.size(); ++i)
p_real += data.mapping_support_points[i] * data.shape(0, i);
return p_real;
}
/**
* Implementation of transform_real_to_unit_cell for either type double
* or VectorizedArray<double>
*/
template <int dim, int spacedim, typename Number>
inline Point<dim, Number>
do_transform_real_to_unit_cell_internal(
const Point<spacedim, Number> & p,
const Point<dim, Number> & initial_p_unit,
const std::vector<Point<spacedim>> & points,
const std::vector<Polynomials::Polynomial<double>> &polynomials_1d,
const std::vector<unsigned int> & renumber,
const bool print_iterations_to_deallog = false)
{
if (print_iterations_to_deallog)
deallog << "Start MappingQ::do_transform_real_to_unit_cell for real "
<< "point [ " << p << " ] " << std::endl;
AssertDimension(points.size(),
Utilities::pow(polynomials_1d.size(), dim));
// Newton iteration to solve
// f(x)=p(x)-p=0
// where we are looking for 'x' and p(x) is the forward transformation
// from unit to real cell. We solve this using a Newton iteration
// x_{n+1}=x_n-[f'(x)]^{-1}f(x)
// The start value is set to be the linear approximation to the cell
// The shape values and derivatives of the mapping at this point are
// previously computed.
Point<dim, Number> p_unit = initial_p_unit;
auto p_real = internal::evaluate_tensor_product_value_and_gradient(
polynomials_1d, points, p_unit, polynomials_1d.size() == 2, renumber);
Tensor<1, spacedim, Number> f = p_real.first - p;
// early out if we already have our point in all SIMD lanes, i.e.,
// f.norm_square() < 1e-24 * p_real.second[0].norm_square(). To enable
// this step also for VectorizedArray where we do not have operator <,
// compare the result to zero which is defined for SIMD types
if (std::max(Number(0.),
f.norm_square() - 1e-24 * p_real.second[0].norm_square()) ==
Number(0.))
return p_unit;
// we need to compare the position of the computed p(x) against the
// given point 'p'. We will terminate the iteration and return 'x' if
// they are less than eps apart. The question is how to choose eps --
// or, put maybe more generally: in which norm we want these 'p' and
// 'p(x)' to be eps apart.
//
// the question is difficult since we may have to deal with very
// elongated cells where we may achieve 1e-12*h for the distance of
// these two points in the 'long' direction, but achieving this
// tolerance in the 'short' direction of the cell may not be possible
//
// what we do instead is then to terminate iterations if
// \| p(x) - p \|_A < eps
// where the A-norm is somehow induced by the transformation of the
// cell. in particular, we want to measure distances relative to the
// sizes of the cell in its principal directions.
//
// to define what exactly A should be, note that to first order we have
// the following (assuming that x* is the solution of the problem, i.e.,
// p(x*)=p):
// p(x) - p = p(x) - p(x*)
// = -grad p(x) * (x*-x) + higher order terms
// This suggest to measure with a norm that corresponds to
// A = {[grad p(x)]^T [grad p(x)]}^{-1}
// because then
// \| p(x) - p \|_A \approx \| x - x* \|
// Consequently, we will try to enforce that
// \| p(x) - p \|_A = \| f \| <= eps
//
// Note that using this norm is a bit dangerous since the norm changes
// in every iteration (A isn't fixed by depending on xk). However, if
// the cell is not too deformed (it may be stretched, but not twisted)
// then the mapping is almost linear and A is indeed constant or
// nearly so.
const double eps = 1.e-11;
const unsigned int newton_iteration_limit = 20;
Point<dim, Number> invalid_point;
invalid_point[0] = std::numeric_limits<double>::infinity();
bool tried_project_to_unit_cell = false;
unsigned int newton_iteration = 0;
Number f_weighted_norm_square = 1.;
Number last_f_weighted_norm_square = 1.;
do
{
if (print_iterations_to_deallog)
deallog << "Newton iteration " << newton_iteration
<< " for unit point guess " << p_unit << std::endl;
// f'(x)
Tensor<2, spacedim, Number> df;
for (unsigned int d = 0; d < spacedim; ++d)
for (unsigned int e = 0; e < dim; ++e)
df[d][e] = p_real.second[e][d];
// check determinand(df) > 0 on all SIMD lanes
if (!(std::min(determinant(df),
Number(std::numeric_limits<double>::min())) ==
Number(std::numeric_limits<double>::min())))
{
// We allow to enter this function with an initial guess
// outside the unit cell. We might have run into invalid
// Jacobians due to that, so we should at least try once to go
// back to the unit cell and go on with the Newton iteration
// from there. Since the outside case is unlikely, we can
// afford spending the extra effort at this place.
if (tried_project_to_unit_cell == false)
{
p_unit = GeometryInfo<dim>::project_to_unit_cell(p_unit);
p_real = internal::evaluate_tensor_product_value_and_gradient(
polynomials_1d,
points,
p_unit,
polynomials_1d.size() == 2,
renumber);
f = p_real.first - p;
f_weighted_norm_square = 1.;
last_f_weighted_norm_square = 1;
tried_project_to_unit_cell = true;
continue;
}
else
return invalid_point;
}
// Solve [f'(x)]d=f(x)
const Tensor<2, spacedim, Number> df_inverse = invert(df);
const Tensor<1, spacedim, Number> delta = df_inverse * f;
last_f_weighted_norm_square = delta.norm_square();
if (print_iterations_to_deallog)
deallog << " delta=" << delta << std::endl;
// do a line search
double step_length = 1.0;
do
{
// update of p_unit. The spacedim-th component of transformed
// point is simply ignored in codimension one case. When this
// component is not zero, then we are projecting the point to
// the surface or curve identified by the cell.
Point<dim, Number> p_unit_trial = p_unit;
for (unsigned int i = 0; i < dim; ++i)
p_unit_trial[i] -= step_length * delta[i];
// shape values and derivatives at new p_unit point
const auto p_real_trial =
internal::evaluate_tensor_product_value_and_gradient(
polynomials_1d,
points,
p_unit_trial,
polynomials_1d.size() == 2,
renumber);
const Tensor<1, spacedim, Number> f_trial =
p_real_trial.first - p;
f_weighted_norm_square = (df_inverse * f_trial).norm_square();
if (print_iterations_to_deallog)
{
deallog << " step_length=" << step_length << std::endl;
if (step_length == 1.0)
deallog << " ||f || =" << f.norm() << std::endl;
deallog << " ||f*|| =" << f_trial.norm() << std::endl
<< " ||f*||_A ="
<< std::sqrt(f_weighted_norm_square) << std::endl;
}
// See if we are making progress with the current step length
// and if not, reduce it by a factor of two and try again.
//
// Strictly speaking, we should probably use the same norm as we
// use for the outer algorithm. In practice, line search is just
// a crutch to find a "reasonable" step length, and so using the
// l2 norm is probably just fine.
//
// check f_trial.norm() < f.norm() in SIMD form. This is a bit
// more complicated because some SIMD lanes might not be doing
// any progress any more as they have already reached roundoff
// accuracy. We define that as the case of updates less than
// 1e-6. The tolerance might seem coarse but since we are
// dealing with a Newton iteration of a polynomial function we
// either converge quadratically or not any more. Thus, our
// selection is to terminate if either the norm of f is
// decreasing or that threshold of 1e-6 is reached.
if (std::max(f_weighted_norm_square - 1e-6 * 1e-6, Number(0.)) *
std::max(f_trial.norm_square() - f.norm_square(),
Number(0.)) ==
Number(0.))
{
p_real = p_real_trial;
p_unit = p_unit_trial;
f = f_trial;
break;
}
else if (step_length > 0.05)
step_length *= 0.5;
else
break;
}
while (true);
// In case we terminated the line search due to the step becoming
// too small, we give the iteration another try with the
// projection of the initial guess to the unit cell before we give
// up, just like for the negative determinant case.
if (step_length <= 0.05 && tried_project_to_unit_cell == false)
{
p_unit = GeometryInfo<dim>::project_to_unit_cell(p_unit);
p_real = internal::evaluate_tensor_product_value_and_gradient(
polynomials_1d,
points,
p_unit,
polynomials_1d.size() == 2,
renumber);
f = p_real.first - p;
f_weighted_norm_square = 1.;
last_f_weighted_norm_square = 1;
tried_project_to_unit_cell = true;
continue;
}
else if (step_length <= 0.05)
return invalid_point;
++newton_iteration;
if (newton_iteration > newton_iteration_limit)
return invalid_point;
}
// Stop if f_weighted_norm_square <= eps^2 on all SIMD lanes or if the
// weighted norm is less than 1e-6 and the improvement against the
// previous step was less than a factor of 10 (in that regime, we
// either have quadratic convergence or roundoff errors due to a bad
// mapping)
while (
!(std::max(f_weighted_norm_square - eps * eps, Number(0.)) *
std::max(last_f_weighted_norm_square -
std::min(f_weighted_norm_square, Number(1e-6 * 1e-6)) *
100.,
Number(0.)) ==
Number(0.)));
if (print_iterations_to_deallog)
deallog << "Iteration converged for p_unit = [ " << p_unit
<< " ] and iteration error "
<< std::sqrt(f_weighted_norm_square) << std::endl;
return p_unit;
}
/**
* Implementation of transform_real_to_unit_cell for dim==spacedim-1
*/
template <int dim>
inline Point<dim>
do_transform_real_to_unit_cell_internal_codim1(
const typename dealii::Triangulation<dim, dim + 1>::cell_iterator &cell,
const Point<dim + 1> & p,
const Point<dim> & initial_p_unit,
typename dealii::MappingQ<dim, dim + 1>::InternalData &mdata)
{
const unsigned int spacedim = dim + 1;
const unsigned int n_shapes = mdata.shape_values.size();
(void)n_shapes;
Assert(n_shapes != 0, ExcInternalError());
Assert(mdata.shape_derivatives.size() == n_shapes, ExcInternalError());
Assert(mdata.shape_second_derivatives.size() == n_shapes,
ExcInternalError());
std::vector<Point<spacedim>> &points = mdata.mapping_support_points;
Assert(points.size() == n_shapes, ExcInternalError());
Point<spacedim> p_minus_F;
Tensor<1, spacedim> DF[dim];
Tensor<1, spacedim> D2F[dim][dim];
Point<dim> p_unit = initial_p_unit;
Point<dim> f;
Tensor<2, dim> df;
// Evaluate first and second derivatives
mdata.compute_shape_function_values(std::vector<Point<dim>>(1, p_unit));
for (unsigned int k = 0; k < mdata.n_shape_functions; ++k)
{
const Tensor<1, dim> & grad_phi_k = mdata.derivative(0, k);
const Tensor<2, dim> & hessian_k = mdata.second_derivative(0, k);
const Point<spacedim> &point_k = points[k];
for (unsigned int j = 0; j < dim; ++j)
{
DF[j] += grad_phi_k[j] * point_k;
for (unsigned int l = 0; l < dim; ++l)
D2F[j][l] += hessian_k[j][l] * point_k;
}
}
p_minus_F = p;
p_minus_F -= compute_mapped_location_of_point<dim, spacedim>(mdata);
for (unsigned int j = 0; j < dim; ++j)
f[j] = DF[j] * p_minus_F;
for (unsigned int j = 0; j < dim; ++j)
{
f[j] = DF[j] * p_minus_F;
for (unsigned int l = 0; l < dim; ++l)
df[j][l] = -DF[j] * DF[l] + D2F[j][l] * p_minus_F;
}
const double eps = 1.e-12 * cell->diameter();
const unsigned int loop_limit = 10;
unsigned int loop = 0;
while (f.norm() > eps && loop++ < loop_limit)
{
// Solve [df(x)]d=f(x)
const Tensor<1, dim> d =
invert(df) * static_cast<const Tensor<1, dim> &>(f);
p_unit -= d;
for (unsigned int j = 0; j < dim; ++j)
{
DF[j].clear();
for (unsigned int l = 0; l < dim; ++l)
D2F[j][l].clear();
}
mdata.compute_shape_function_values(
std::vector<Point<dim>>(1, p_unit));
for (unsigned int k = 0; k < mdata.n_shape_functions; ++k)
{
const Tensor<1, dim> & grad_phi_k = mdata.derivative(0, k);
const Tensor<2, dim> & hessian_k = mdata.second_derivative(0, k);
const Point<spacedim> &point_k = points[k];
for (unsigned int j = 0; j < dim; ++j)
{
DF[j] += grad_phi_k[j] * point_k;
for (unsigned int l = 0; l < dim; ++l)
D2F[j][l] += hessian_k[j][l] * point_k;
}
}
// TODO: implement a line search here in much the same way as for
// the corresponding function above that does so for dim==spacedim
p_minus_F = p;
p_minus_F -= compute_mapped_location_of_point<dim, spacedim>(mdata);
for (unsigned int j = 0; j < dim; ++j)
{
f[j] = DF[j] * p_minus_F;
for (unsigned int l = 0; l < dim; ++l)
df[j][l] = -DF[j] * DF[l] + D2F[j][l] * p_minus_F;
}
}
// Here we check that in the last execution of while the first
// condition was already wrong, meaning the residual was below
// eps. Only if the first condition failed, loop will have been
// increased and tested, and thus have reached the limit.
AssertThrow(loop < loop_limit,
(typename Mapping<dim, spacedim>::ExcTransformationFailed()));
return p_unit;
}
/**
* A class to compute a quadratic approximation to the inverse map from
* real to unit points by a least-squares fit along the mapping support
* points. The least squares fit is special in the sense that the
* approximation is constructed for the inverse function of a
* MappingQ, which is generally a rational function. This allows
* for a very cheap evaluation of the inverse map by a simple polynomial
* interpolation, which can be used as a better initial guess for
* transforming points from real to unit coordinates than an affine
* approximation.
*
* Far away outside the unit cell, this approximation can become
* inaccurate for non-affine cell shapes. This must be expected from a
* fit of a polynomial to a rational function, and due to the fact that
* the region of the least squares fit, the unit cell, is left. Hence,
* use this function with care in those situations.
*/
template <int dim, int spacedim>
class InverseQuadraticApproximation
{
public:
/**
* Number of basis functions in the quadratic approximation.
*/
static constexpr unsigned int n_functions =
(spacedim == 1 ? 3 : (spacedim == 2 ? 6 : 10));
/**
* Constructor.
*
* @param real_support_points The position of the mapping support points
* in real space, queried by
* MappingQ::compute_mapping_support_points().
*
* @param unit_support_points The location of the support points in
* reference coordinates $[0, 1]^d$ that map to the mapping support
* points in real space by a polynomial map.
*/
InverseQuadraticApproximation(
const std::vector<Point<spacedim>> &real_support_points,
const std::vector<Point<dim>> & unit_support_points)
: normalization_shift(real_support_points[0])
, normalization_length(
1. / real_support_points[0].distance(real_support_points[1]))
, is_affine(true)
{
AssertDimension(real_support_points.size(), unit_support_points.size());
// For the bi-/trilinear approximation, we cannot build a quadratic
// polynomial due to a lack of points (interpolation matrix would
// get singular). Similarly, it is not entirely clear how to gather
// enough information for the case dim < spacedim.
//
// In both cases we require the vector real_support_points to
// contain the vertex positions and fall back to an affine
// approximation:
Assert(dim == spacedim || real_support_points.size() ==
GeometryInfo<dim>::vertices_per_cell,
ExcInternalError());
if (real_support_points.size() == GeometryInfo<dim>::vertices_per_cell)
{
const auto affine = GridTools::affine_cell_approximation<dim>(
make_array_view(real_support_points));
DerivativeForm<1, spacedim, dim> A_inv =
affine.first.covariant_form().transpose();
// The code for evaluation assumes an additional transformation of
// the form (x - normalization_shift) * normalization_length --
// account for this in the definition of the coefficients.
coefficients[0] =
apply_transformation(A_inv, normalization_shift - affine.second);
for (unsigned int d = 0; d < spacedim; ++d)
for (unsigned int e = 0; e < dim; ++e)
coefficients[1 + d][e] =
A_inv[e][d] * (1.0 / normalization_length);
is_affine = true;
return;
}
Tensor<2, n_functions> matrix;
std::array<double, n_functions> shape_values;
for (unsigned int q = 0; q < unit_support_points.size(); ++q)
{
// Evaluate quadratic shape functions in point, with the
// normalization applied in order to avoid roundoff issues with
// scaling far away from 1.
shape_values[0] = 1.;
const Tensor<1, spacedim> p_scaled =
normalization_length *
(real_support_points[q] - normalization_shift);
for (unsigned int d = 0; d < spacedim; ++d)
shape_values[1 + d] = p_scaled[d];
for (unsigned int d = 0, c = 0; d < spacedim; ++d)
for (unsigned int e = 0; e <= d; ++e, ++c)
shape_values[1 + spacedim + c] = p_scaled[d] * p_scaled[e];
// Build lower diagonal of least squares matrix and rhs, the
// essential part being that we construct the matrix with the
// real points and the right hand side by comparing to the
// reference point positions which sets up an inverse
// interpolation.
for (unsigned int i = 0; i < n_functions; ++i)
for (unsigned int j = 0; j < n_functions; ++j)
matrix[i][j] += shape_values[i] * shape_values[j];
for (unsigned int i = 0; i < n_functions; ++i)
coefficients[i] += shape_values[i] * unit_support_points[q];
}
// Factorize the matrix A = L * L^T in-place with the
// Cholesky-Banachiewicz algorithm. The implementation is similar to
// FullMatrix::cholesky() but re-implemented to avoid memory
// allocations and some unnecessary divisions which we can do here as
// we only need to solve with dim right hand sides.
for (unsigned int i = 0; i < n_functions; ++i)
{
double Lij_sum = 0;
for (unsigned int j = 0; j < i; ++j)
{
double Lik_Ljk_sum = 0;
for (unsigned int k = 0; k < j; ++k)
Lik_Ljk_sum += matrix[i][k] * matrix[j][k];
matrix[i][j] = matrix[j][j] * (matrix[i][j] - Lik_Ljk_sum);
Lij_sum += matrix[i][j] * matrix[i][j];
}
AssertThrow(matrix[i][i] - Lij_sum >= 0,
ExcMessage("Matrix of normal equations not positive "
"definite"));
// Store the inverse in the diagonal since that is the quantity
// needed later in the factorization as well as the forward and
// backward substitution, minimizing the number of divisions.
matrix[i][i] = 1. / std::sqrt(matrix[i][i] - Lij_sum);
}
// Solve lower triangular part, L * y = rhs.
for (unsigned int i = 0; i < n_functions; ++i)
{
Point<dim> sum = coefficients[i];
for (unsigned int j = 0; j < i; ++j)
sum -= matrix[i][j] * coefficients[j];
coefficients[i] = sum * matrix[i][i];
}
// Solve upper triangular part, L^T * x = y (i.e., x = A^{-1} * rhs)
for (unsigned int i = n_functions; i > 0;)