-
Notifications
You must be signed in to change notification settings - Fork 227
/
chunk.cc
1055 lines (884 loc) · 36.5 KB
/
chunk.cc
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) 2011 - 2018 by the authors of the ASPECT code.
This file is part of ASPECT.
ASPECT is free software; you can redistribute it and/or modify
it under the terms of the GNU General Public License as published by
the Free Software Foundation; either version 2, or (at your option)
any later version.
ASPECT 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 General Public License for more details.
You should have received a copy of the GNU General Public License
along with ASPECT; see the file LICENSE. If not see
<http://www.gnu.org/licenses/>.
*/
#include <aspect/geometry_model/chunk.h>
#include <aspect/geometry_model/initial_topography_model/zero_topography.h>
#include <aspect/geometry_model/initial_topography_model/ascii_data.h>
#include <aspect/simulator_signals.h>
#include <deal.II/grid/grid_generator.h>
#include <deal.II/grid/tria_iterator.h>
#include <deal.II/grid/tria_accessor.h>
#include <deal.II/grid/grid_tools.h>
#if !DEAL_II_VERSION_GTE(9,0,0)
#include <deal.II/grid/tria_boundary_lib.h>
#endif
#if DEAL_II_VERSION_GTE(9,0,0)
#include <deal.II/grid/tria.h>
#endif
namespace aspect
{
namespace GeometryModel
{
template <int dim>
Chunk<dim>::ChunkGeometry::ChunkGeometry()
:
point1_lon(0.0),
inner_radius(3471e3),
max_depth(2900e3)
{}
template <int dim>
Chunk<dim>::ChunkGeometry::ChunkGeometry(const ChunkGeometry &other)
:
ChartManifold<dim,dim>(other),
point1_lon(other.point1_lon),
inner_radius(other.inner_radius),
max_depth(other.max_depth)
{
this->initialize(other.topo);
}
template <int dim>
void
Chunk<dim>::ChunkGeometry::initialize(const InitialTopographyModel::Interface<dim> *topo_pointer)
{
topo = topo_pointer;
}
template <int dim>
DerivativeForm<1,dim,dim>
Chunk<dim>::ChunkGeometry::
push_forward_gradient(const Point<dim> &chart_point) const
{
const double R = chart_point[0]; // Radius
Assert (R > 0.0, ExcMessage("Negative radius for given point."));
Tensor<2,dim> DX;
// prior to deal.II 9, we do not apply initial topography
#if !DEAL_II_VERSION_GTE(9,0,0)
const double phi = chart_point[1]; // Longitude
switch (dim)
{
case 2:
{
DX[0][0] = std::cos(phi);
DX[0][1] = -R * std::sin(phi);
DX[1][0] = std::sin(phi);
DX[1][1] = R * std::cos(phi);
break;
}
case 3:
{
const double theta = chart_point[2]; // Latitude (not colatitude)
DX[0][0] = std::cos(theta) * std::cos(phi);
DX[0][1] = -R * std::cos(theta) * std::sin(phi);
DX[0][2] = -R * std::sin(theta) * std::cos(phi);
DX[1][0] = std::cos(theta) * std::sin(phi);
DX[1][1] = R * std::cos(theta) * std::cos(phi);
DX[1][2] = -R * std::sin(theta) * std::sin(phi);
DX[2][0] = std::sin(theta);
DX[2][1] = 0;
DX[2][2] = R * std::cos(theta);
break;
}
default:
Assert (false, ExcNotImplemented ());
}
return DerivativeForm<1,dim,dim>(DX);
#else
// The initial topography derivatives (dtopo/dphi and dtopo/dtheta)
// They are zero for the ZeroTopography model
Tensor<1,dim-1> topo_derivatives;
if (const InitialTopographyModel::AsciiData<dim> *itm = dynamic_cast<const InitialTopographyModel::AsciiData<dim> *> (topo))
topo_derivatives = itm->vector_gradient(push_forward_sphere(chart_point));
// Construct surface point in lon(,lat) coordinates
Point<dim-1> surface_point;
for (unsigned int d=0; d<dim-1; d++)
surface_point[d] = chart_point[d+1];
// Convert latitude to colatitude
if (dim == 3)
surface_point[1] = 0.5*numbers::PI - surface_point[1];
// get the maximum topography at the surface point
const double d_topo = topo->value(surface_point);
// get the spherical point including topography
const Point<dim> topo_point = push_forward_topo(chart_point);
const double R_topo = topo_point[0];
const double phi_topo = topo_point[1];
// The derivatives of topo_point to chart_point
Tensor<2, dim> Dtopo;
// The derivatives of the cartesian point to chart_point
DerivativeForm<1, dim, dim> Dtotal;
switch (dim)
{
case 2:
{
// R_topo = R + topo(phi) * ((R-R_0)/(R_1-R_0)) = R + topo(phi) * ((R-R_0)/max_depth)
// phi_topo = phi
//dR_topo/dR
Dtopo[0][0] = (d_topo / max_depth) + 1.;
//dR_topo/dphi = dR_topo/dtopo * dtopo/dphi
Dtopo[0][1] = (R-inner_radius)/max_depth * topo_derivatives[0];
//dphi_topo/dR
Dtopo[1][0] = 0.;
//dphi_topo/dphi
Dtopo[1][1] = 1.;
//dx/dR_topo
DX[0][0] = std::cos(phi_topo);
//dx/dphi_topo
DX[0][1] = -R_topo * std::sin(phi_topo);
//dy/dR_topo
DX[1][0] = std::sin(phi_topo);
//dy/dphi_topo
DX[1][1] = R_topo * std::cos(phi_topo);
break;
}
case 3:
{
// R_topo = R + topo(phi,theta) * ((R-R_0)/(R_1-R_0))
// phi_topo = phi
// theta_topo = theta
//dR_topo/dR
Dtopo[0][0] = (d_topo / max_depth) + 1.;
//dR_topo/dphi
Dtopo[0][1] = (R - inner_radius) / max_depth * topo_derivatives[0];
//dR_topo/dtheta
Dtopo[0][2] = (R - inner_radius) / max_depth * topo_derivatives[1];
//dphi_topo/dR
Dtopo[1][0] = 0.;
//dphi_topo/dphi
Dtopo[1][1] = 1.;
//dphi_topo/dtheta
Dtopo[1][2] = 0.;
//dtheta_topo/dR
Dtopo[2][0] = 0.;
//dtheta_topo/dphi
Dtopo[2][1] = 0.;
//dtheta_topo/dtheta
Dtopo[2][2] = 1.;
const double theta_topo = topo_point[2];
// The derivatives of the cartesian points to topo_point
DX[0][0] = std::cos(theta_topo) * std::cos(phi_topo);
DX[0][1] = -R_topo * std::cos(theta_topo) * std::sin(phi_topo);
DX[0][2] = -R_topo * std::sin(theta_topo) * std::cos(phi_topo); //reorder
DX[1][0] = std::cos(theta_topo) * std::sin(phi_topo);
DX[1][1] = R_topo * std::cos(theta_topo) * std::cos(phi_topo);
DX[1][2] = -R_topo * std::sin(theta_topo) * std::sin(phi_topo);
DX[2][0] = std::sin(theta_topo);
DX[2][1] = 0;
DX[2][2] = R_topo * std::cos(theta_topo);
break;
}
default:
Assert (false, ExcNotImplemented ());
}
Dtotal = DerivativeForm<1,dim,dim>(DX * Dtopo);
return Dtotal;
#endif
}
template <int dim>
Point<dim>
Chunk<dim>::ChunkGeometry::
push_forward(const Point<dim> &r_phi_theta) const
{
#if !DEAL_II_VERSION_GTE(9,0,0)
return push_forward_sphere(r_phi_theta);
#else
// Only take into account topography when we're not using the ZeroTopography plugin
if (dynamic_cast<const InitialTopographyModel::ZeroTopography<dim>*>(topo) != nullptr)
return push_forward_sphere(r_phi_theta);
else
return push_forward_sphere(push_forward_topo(r_phi_theta));
#endif
}
template <int dim>
Point<dim>
Chunk<dim>::ChunkGeometry::
pull_back(const Point<dim> &x_y_z) const
{
#if !DEAL_II_VERSION_GTE(9,0,0)
return pull_back_sphere(x_y_z);
#else
// Only take into account topography when we're not using the ZeroTopography plugin
if (dynamic_cast<const InitialTopographyModel::ZeroTopography<dim>*>(topo) != nullptr)
return pull_back_sphere(x_y_z);
else
return pull_back_topo(pull_back_sphere(x_y_z));
#endif
}
template <int dim>
Point<dim>
Chunk<dim>::ChunkGeometry::
push_forward_sphere(const Point<dim> &input_vertex) const
{
Point<dim> output_vertex;
switch (dim)
{
case 2:
{
output_vertex[0] = input_vertex[0]*std::cos(input_vertex[1]); // x=rcosphi
output_vertex[1] = input_vertex[0]*std::sin(input_vertex[1]); // z=rsinphi
break;
}
case 3:
{
output_vertex[0] = input_vertex[0]*std::cos(input_vertex[2])*std::cos(input_vertex[1]); // x=rsinthetacosphi
output_vertex[1] = input_vertex[0]*std::cos(input_vertex[2])*std::sin(input_vertex[1]); // y=rsinthetasinphi
output_vertex[2] = input_vertex[0]*std::sin(input_vertex[2]); // z=rcostheta
break;
}
default:
Assert (false, ExcNotImplemented ());
}
return output_vertex;
}
template <int dim>
Point<dim>
Chunk<dim>::ChunkGeometry::
pull_back_sphere(const Point<dim> &v) const
{
Point<dim> output_vertex;
switch (dim)
{
case 2:
{
output_vertex[1] = std::atan2(v[1], v[0]);
output_vertex[0] = v.norm();
// We must guarantee that all returned points have a longitude coordinate
// value that is larger than (or equal to) the longitude of point1.
// For example:
// If the domain runs from longitude -10 to 200 degrees,
// atan2 will also return a negative value (-180 to -160) for the points
// with longitude 180 to 200. These values must be corrected
// so that they are larger than the minimum longitude value of -10,
// by adding 360 degrees.
// A 100*epsilon ensures we catch all cases.
if (output_vertex[1] < 0.0)
if (output_vertex[1] < point1_lon - 100 * std::abs(point1_lon)*std::numeric_limits<double>::epsilon())
output_vertex[1] += 2.0 * numbers::PI;
break;
}
case 3:
{
const double radius=v.norm();
output_vertex[0] = radius;
output_vertex[1] = std::atan2(v[1], v[0]);
// See 2D case
if (output_vertex[1] < 0.0)
if (output_vertex[1] < point1_lon - 100 * std::abs(point1_lon)*std::numeric_limits<double>::epsilon())
output_vertex[1] += 2.0 * numbers::PI;
output_vertex[2] = std::asin(v[2]/radius);
break;
}
default:
Assert (false, ExcNotImplemented ());
}
return output_vertex;
}
template <int dim>
void
Chunk<dim>::ChunkGeometry::
set_min_longitude(const double p1_lon)
{
point1_lon = p1_lon;
}
#if DEAL_II_VERSION_GTE(9,0,0)
template <int dim>
std::unique_ptr<Manifold<dim,dim> >
Chunk<dim>::ChunkGeometry::
clone() const
{
return std_cxx14::make_unique<ChunkGeometry>(*this);
}
#endif
template <int dim>
Point<dim>
Chunk<dim>::ChunkGeometry::
push_forward_topo(const Point<dim> &r_phi_theta) const
{
// the radius of the current point without topography
const double radius = r_phi_theta[0];
// Grab lon,lat coordinates
Point<dim-1> surface_point;
for (unsigned int d=0; d<dim-1; d++)
surface_point[d] = r_phi_theta[d+1];
// Convert latitude to colatitude
if (dim == 3)
surface_point[1] = 0.5*numbers::PI - surface_point[1];
const double topography = topo->value(surface_point);
// adjust the radius based on the radius of the point
// through a linear interpolation between 0 at max depth and
// "topography" at the surface
const double topo_radius = std::max(inner_radius,radius + (radius-inner_radius)/max_depth*topography);
// return the point with adjusted radius
Point<dim> topor_phi_theta = r_phi_theta;
topor_phi_theta[0] = topo_radius;
return topor_phi_theta;
}
template <int dim>
Point<dim>
Chunk<dim>::ChunkGeometry::
pull_back_topo(const Point<dim> &topor_phi_theta) const
{
// the radius of the point with topography
const double topo_radius = topor_phi_theta[0];
// Grab lon,lat coordinates
Point<dim-1> surface_point;
for (unsigned int d=0; d<dim-1; d++)
surface_point[d] = topor_phi_theta[d+1];
// Convert latitude to colatitude
if (dim == 3)
surface_point[1] = 0.5*numbers::PI - surface_point[1];
const double topography = topo->value(surface_point);
// remove the topography (which scales with radius)
const double radius = std::max(inner_radius,(topo_radius*max_depth+inner_radius*topography)/(max_depth+topography));
// return the point without topography
Point<dim> r_phi_theta = topor_phi_theta;
r_phi_theta[0] = radius;
return r_phi_theta;
}
template <int dim>
double
Chunk<dim>::ChunkGeometry::
get_radius(const Point<dim> &x_y_z) const
{
const Point<dim> r_phi_theta = pull_back(x_y_z);
Point<dim-1> surface_point;
for (unsigned int d=0; d<dim-1; d++)
surface_point[d] = r_phi_theta[d+1];
// Convert latitude to colatitude
if (dim == 3)
surface_point[1] = 0.5*numbers::PI - surface_point[1];
const double topography = topo->value(surface_point);
// return the outer radius at this phi, theta point including topography
return topography + inner_radius + max_depth;
}
template <int dim>
void
Chunk<dim>::ChunkGeometry::
set_min_radius(const double p1_rad)
{
inner_radius = p1_rad;
}
template <int dim>
void
Chunk<dim>::ChunkGeometry::
set_max_depth(const double p2_p1_rad)
{
max_depth = p2_p1_rad;
}
template <int dim>
void
Chunk<dim>::initialize ()
{
#if !DEAL_II_VERSION_GTE(9,0,0)
// Call function to connect the set/clear manifold id functions
// to the right signal
connect_to_signal(this->get_signals());
AssertThrow(dynamic_cast<const InitialTopographyModel::ZeroTopography<dim>*>(&this->get_initial_topography_model()) != nullptr,
ExcMessage("Only with deal.II 9 or higher, an initial topography model can be used."));
#else
AssertThrow(dynamic_cast<const InitialTopographyModel::AsciiData<dim>*>(&this->get_initial_topography_model()) != nullptr ||
dynamic_cast<const InitialTopographyModel::ZeroTopography<dim>*>(&this->get_initial_topography_model()) != nullptr,
ExcMessage("At the moment, only the Zero or AsciiData initial topography model can be used."));
#endif
manifold.initialize(&(this->get_initial_topography_model()));
}
template <int dim>
void
Chunk<dim>::
create_coarse_mesh (parallel::distributed::Triangulation<dim> &coarse_grid) const
{
std::vector<unsigned int> rep_vec(repetitions, repetitions+dim);
GridGenerator::subdivided_hyper_rectangle (coarse_grid,
rep_vec,
point1,
point2,
true);
#if !DEAL_II_VERSION_GTE(9,0,0)
// At this point, all boundary faces have their correct boundary
// indicators, but the edges do not. We want the edges of curved
// faces to be curved as well, so we set the edge boundary indicators
// to the same boundary indicators as their faces.
for (typename Triangulation<dim>::active_cell_iterator
cell = coarse_grid.begin_active();
cell != coarse_grid.end(); ++cell)
for (unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; ++f)
if (cell->face(f)->at_boundary())
// Set edges on the radial faces; both adjacent faces
// should agree on where new points along the boundary lie
// for these edges, so the order of the boundaries does not matter
if ((cell->face(f)->boundary_id() == 2)
||
(cell->face(f)->boundary_id() == 3))
cell->face(f)->set_all_boundary_ids(cell->face(f)->boundary_id());
for (typename Triangulation<dim>::active_cell_iterator
cell = coarse_grid.begin_active();
cell != coarse_grid.end(); ++cell)
for (unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; ++f)
if (cell->face(f)->at_boundary())
// Set edges on the radial faces; both adjacent faces
// should agree on where new points along the boundary lie
// for these edges, so the order of the boundaries does not matter
if ((cell->face(f)->boundary_id() == 4)
||
(cell->face(f)->boundary_id() == 5))
cell->face(f)->set_all_boundary_ids(cell->face(f)->boundary_id());
for (typename Triangulation<dim>::active_cell_iterator
cell = coarse_grid.begin_active();
cell != coarse_grid.end(); ++cell)
for (unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; ++f)
if (cell->face(f)->at_boundary())
// (Re-)Set edges on the spherical shells to ensure that
// they are all curved as expected
if ((cell->face(f)->boundary_id() == 0)
||
(cell->face(f)->boundary_id() == 1))
cell->face(f)->set_all_boundary_ids(cell->face(f)->boundary_id());
#endif
// Transform box into spherical chunk
GridTools::transform (
[&](const Point<dim> &p) -> Point<dim>
{
return manifold.push_forward(p);
},
coarse_grid);
// Deal with a curved mesh
// Attach the real manifold to slot 15.
coarse_grid.set_manifold (15, manifold);
for (typename Triangulation<dim>::active_cell_iterator cell =
coarse_grid.begin_active(); cell != coarse_grid.end(); ++cell)
cell->set_all_manifold_ids (15);
#if !DEAL_II_VERSION_GTE(9,0,0)
// On the boundary faces, set boundary objects.
// The east and west boundaries are straight,
// the inner and outer boundary are part of
// a spherical shell. In 3D, the north and south boundaries
// are part of a cone with the tip at the origin.
static const StraightBoundary<dim> boundary_straight;
// Attach boundary objects to the straight east and west boundaries
coarse_grid.set_boundary(2, boundary_straight);
coarse_grid.set_boundary(3, boundary_straight);
if (dim == 3)
{
// Define the center point of the greater radius end of the
// north and south boundary cones.
// These lie along the z-axis.
Point<dim> center;
Point<dim> north, south;
const double outer_radius = point2[0];
north[dim-1] = outer_radius * std::sin(point2[2]);
south[dim-1] = outer_radius * std::sin(point1[2]);
// Define the radius of the cones
const double north_radius = std::sqrt(outer_radius*outer_radius-north[dim-1]*north[dim-1]);
const double south_radius = std::sqrt(outer_radius*outer_radius-south[dim-1]*south[dim-1]);
static const ConeBoundary<dim> boundary_cone_north(0.0,north_radius,center,north);
static const ConeBoundary<dim> boundary_cone_south(0.0,south_radius,center,south);
// Attach boundary objects to the conical north and south boundaries
// If one of the boundaries lies at the equator,
// just use the straight boundary.
if (point2[2] != 0.0)
coarse_grid.set_boundary (5, boundary_cone_north);
else
coarse_grid.set_boundary (5, boundary_straight);
if (point1[2] != 0.0)
coarse_grid.set_boundary (4, boundary_cone_south);
else
coarse_grid.set_boundary (4, boundary_straight);
}
// Attach shell boundary objects to the curved inner and outer boundaries
static const HyperShellBoundary<dim> boundary_shell;
coarse_grid.set_boundary (0, boundary_shell);
coarse_grid.set_boundary (1, boundary_shell);
#endif
}
#if !DEAL_II_VERSION_GTE(9,0,0)
template <int dim>
void
Chunk<dim>::set_manifold_ids (typename parallel::distributed::Triangulation<dim> &triangulation)
{
// Set all cells, faces and edges to manifold_id 15
for (typename Triangulation<dim>::active_cell_iterator cell =
triangulation.begin_active(); cell != triangulation.end(); ++cell)
cell->set_all_manifold_ids (15);
}
template <int dim>
void
Chunk<dim>::clear_manifold_ids (typename parallel::distributed::Triangulation<dim> &triangulation)
{
// Clear the manifold_id from the faces and edges at the boundary
// so that the boundary objects can be used
for (typename Triangulation<dim>::active_cell_iterator cell =
triangulation.begin_active(); cell != triangulation.end(); ++cell)
for (unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; ++f)
if (cell->face(f)->at_boundary())
cell->face(f)->set_all_manifold_ids (numbers::flat_manifold_id);
}
template <int dim>
void
Chunk<dim>::connect_to_signal (SimulatorSignals<dim> &signals)
{
// Connect the topography function to the signal
signals.pre_compute_no_normal_flux_constraints.connect (
[&](typename parallel::distributed::Triangulation<dim> &tria)
{
this->clear_manifold_ids(tria);
});
signals.post_compute_no_normal_flux_constraints.connect (
[&](typename parallel::distributed::Triangulation<dim> &tria)
{
this->set_manifold_ids(tria);
});
}
#endif
template <int dim>
std::set<types::boundary_id>
Chunk<dim>::
get_used_boundary_indicators () const
{
// boundary indicators are zero through 2*dim-1
std::set<types::boundary_id> s;
for (unsigned int i=0; i<2*dim; ++i)
s.insert (i);
return s;
}
template <int dim>
std::map<std::string,types::boundary_id>
Chunk<dim>::
get_symbolic_boundary_names_map () const
{
switch (dim)
{
case 2:
{
static const std::pair<std::string,types::boundary_id> mapping[]
= { std::pair<std::string,types::boundary_id>("bottom", 0),
std::pair<std::string,types::boundary_id>("top", 1),
std::pair<std::string,types::boundary_id>("west", 2),
std::pair<std::string,types::boundary_id>("east", 3)
};
return std::map<std::string,types::boundary_id> (&mapping[0],
&mapping[sizeof(mapping)/sizeof(mapping[0])]);
}
case 3:
{
static const std::pair<std::string,types::boundary_id> mapping[]
= { std::pair<std::string,types::boundary_id>("bottom", 0),
std::pair<std::string,types::boundary_id>("top", 1),
std::pair<std::string,types::boundary_id>("west", 2),
std::pair<std::string,types::boundary_id>("east", 3),
std::pair<std::string,types::boundary_id>("south", 4),
std::pair<std::string,types::boundary_id>("north", 5)
};
return std::map<std::string,types::boundary_id> (&mapping[0],
&mapping[sizeof(mapping)/sizeof(mapping[0])]);
}
}
Assert (false, ExcNotImplemented());
return std::map<std::string,types::boundary_id>();
}
template <int dim>
double
Chunk<dim>::
length_scale () const
{
// As described in the first ASPECT paper, a length scale of
// 10km = 1e4m works well for the pressure scaling for earth
// sized spherical shells. use a length scale that
// yields this value for the R0,R1 corresponding to earth
// but otherwise scales like (R1-R0)
return 1e4 * maximal_depth() / (6336000.-3481000.);
}
template <int dim>
double
Chunk<dim>::depth(const Point<dim> &position) const
{
// depth is defined wrt the reference surface point2[0]
// negative depth is not allowed
return std::max (0., std::min (point2[0]-position.norm(), maximal_depth()));
}
template <int dim>
double
Chunk<dim>::depth_wrt_topo(const Point<dim> &position) const
{
// depth is defined wrt the reference surface point2[0] + the topography
// depth is therefore always positive
const double outer_radius = manifold.get_radius(position);
const Point<dim> rtopo_phi_theta = manifold.pull_back_sphere(position);
Assert (rtopo_phi_theta[0] <= outer_radius, ExcMessage("The radius is bigger than the maximum radius."));
return std::max(0.0, outer_radius - rtopo_phi_theta[0]);
}
template <int dim>
double
Chunk<dim>::height_above_reference_surface(const Point<dim> &position) const
{
return position.norm()-point2[0];
}
template <int dim>
Point<dim>
Chunk<dim>::representative_point(const double depth) const
{
Assert (depth >= 0,
ExcMessage ("Given depth must be positive or zero."));
Assert (depth <= maximal_depth(),
ExcMessage ("Given depth must be less than or equal to the maximal depth of this geometry."));
// Choose a point at the mean longitude (and latitude)
Point<dim> p = 0.5*(point2+point1);
// at a depth beneath the top surface
p[0] = point2[0]-depth;
// Now convert to Cartesian coordinates
return manifold.push_forward_sphere(p);
}
template <int dim>
double
Chunk<dim>::west_longitude () const
{
return point1[1];
}
template <int dim>
double
Chunk<dim>::east_longitude () const
{
return point2[1];
}
template <int dim>
double
Chunk<dim>::longitude_range () const
{
return point2[1] - point1[1];
}
template <int dim>
double
Chunk<dim>::south_latitude () const
{
if (dim == 3)
return point1[2];
else
return 0;
}
template <int dim>
double
Chunk<dim>::north_latitude () const
{
if (dim==3)
return point2[2];
else
return 0;
}
template <int dim>
double
Chunk<dim>::latitude_range () const
{
if (dim==3)
return point2[2] - point1[2];
else
return 0;
}
template <int dim>
double
Chunk<dim>::maximal_depth() const
{
return point2[0]-point1[0];
}
template <int dim>
double
Chunk<dim>::inner_radius () const
{
return point1[0];
}
template <int dim>
double
Chunk<dim>::outer_radius () const
{
return point2[0];
}
template <int dim>
bool
Chunk<dim>::has_curved_elements() const
{
return true;
}
template <int dim>
bool
Chunk<dim>::point_is_in_domain(const Point<dim> &point) const
{
AssertThrow(this->get_free_surface_boundary_indicators().size() == 0 ||
// we are still before the first time step has started
this->get_timestep_number() == numbers::invalid_unsigned_int,
ExcMessage("After displacement of the free surface, this function can no longer be used to determine whether a point lies in the domain or not."));
AssertThrow(dynamic_cast<const InitialTopographyModel::ZeroTopography<dim>*>(&this->get_initial_topography_model()) != nullptr,
ExcMessage("After adding topography, this function can no longer be used to determine whether a point lies in the domain or not."));
const Point<dim> spherical_point = manifold.pull_back(point);
for (unsigned int d = 0; d < dim; d++)
if (spherical_point[d] > point2[d]+std::numeric_limits<double>::epsilon()*std::abs(point2[d]) ||
spherical_point[d] < point1[d]-std::numeric_limits<double>::epsilon()*std::abs(point2[d]))
return false;
return true;
}
template <int dim>
std::array<double,dim>
Chunk<dim>::cartesian_to_natural_coordinates(const Point<dim> &position_point) const
{
// the chunk manifold has a order of radius, longitude, latitude.
// This is exactly what we need.
const Point<dim> transformed_point = manifold.pull_back(position_point);
std::array<double,dim> position_array;
for (unsigned int i = 0; i < dim; i++)
position_array[i] = transformed_point(i);
return position_array;
}
template <int dim>
aspect::Utilities::Coordinates::CoordinateSystem
Chunk<dim>::natural_coordinate_system() const
{
return aspect::Utilities::Coordinates::CoordinateSystem::spherical;
}
template <int dim>
Point<dim>
Chunk<dim>::natural_to_cartesian_coordinates(const std::array<double,dim> &position_tensor) const
{
Point<dim> position_point;
for (unsigned int i = 0; i < dim; i++)
position_point[i] = position_tensor[i];
const Point<dim> transformed_point = manifold.push_forward(position_point);
return transformed_point;
}
template <int dim>
void
Chunk<dim>::
declare_parameters (ParameterHandler &prm)
{
prm.enter_subsection("Geometry model");
{
prm.enter_subsection("Chunk");
{
prm.declare_entry ("Chunk inner radius", "0",
Patterns::Double (0),
"Radius at the bottom surface of the chunk. Units: m.");
prm.declare_entry ("Chunk outer radius", "1",
Patterns::Double (0),
"Radius at the top surface of the chunk. Units: m.");
prm.declare_entry ("Chunk minimum longitude", "0",
Patterns::Double (-180, 360), // enables crossing of either hemisphere
"Minimum longitude of the chunk. Units: degrees.");
prm.declare_entry ("Chunk maximum longitude", "1",
Patterns::Double (-180, 360), // enables crossing of either hemisphere
"Maximum longitude of the chunk. Units: degrees.");
prm.declare_entry ("Chunk minimum latitude", "0",
Patterns::Double (-90, 90),
"Minimum latitude of the chunk. This value is ignored "
"if the simulation is in 2d. Units: degrees.");
prm.declare_entry ("Chunk maximum latitude", "1",
Patterns::Double (-90, 90),
"Maximum latitude of the chunk. This value is ignored "
"if the simulation is in 2d. Units: degrees.");
prm.declare_entry ("Radius repetitions", "1",
Patterns::Integer (1),
"Number of cells in radius.");
prm.declare_entry ("Longitude repetitions", "1",
Patterns::Integer (1),
"Number of cells in longitude.");
prm.declare_entry ("Latitude repetitions", "1",
Patterns::Integer (1),
"Number of cells in latitude. This value is ignored "
"if the simulation is in 2d");
}
prm.leave_subsection();
}
prm.leave_subsection();
}
template <int dim>
void
Chunk<dim>::parse_parameters (ParameterHandler &prm)
{
prm.enter_subsection("Geometry model");
{
prm.enter_subsection("Chunk");
{
const double degtorad = dealii::numbers::PI/180;
Assert (dim >= 2, ExcInternalError());
Assert (dim <= 3, ExcInternalError());
if (dim >= 2)
{
point1[0] = prm.get_double ("Chunk inner radius");
point2[0] = prm.get_double ("Chunk outer radius");
repetitions[0] = prm.get_integer ("Radius repetitions");
point1[1] = prm.get_double ("Chunk minimum longitude") * degtorad;
point2[1] = prm.get_double ("Chunk maximum longitude") * degtorad;
repetitions[1] = prm.get_integer ("Longitude repetitions");
AssertThrow (point1[0] < point2[0],
ExcMessage ("Inner radius must be less than outer radius."));
AssertThrow (point1[1] < point2[1],
ExcMessage ("Minimum longitude must be less than maximum longitude."));
AssertThrow (point2[1] - point1[1] < 2.*numbers::PI,
ExcMessage ("Maximum - minimum longitude should be less than 360 degrees."));
}
// Inform the manifold about the minimum longitude
manifold.set_min_longitude(point1[1]);
// Inform the manifold about the minimum radius
manifold.set_min_radius(point1[0]);
// Inform the manifold about the maximum depth (without topo)
manifold.set_max_depth(point2[0]-point1[0]);
if (dim == 3)
{
point1[2] = prm.get_double ("Chunk minimum latitude") * degtorad;
point2[2] = prm.get_double ("Chunk maximum latitude") * degtorad;
repetitions[2] = prm.get_integer ("Latitude repetitions");