-
Notifications
You must be signed in to change notification settings - Fork 37
/
Copy pathfault.cc
758 lines (633 loc) · 42.9 KB
/
fault.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
/*
Copyright (C) 2018 - 2020 by the authors of the World Builder code.
This file is part of the World Builder.
This program 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 of the License, or
(at your option) any later version.
This program 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 program. If not, see <https://www.gnu.org/licenses/>.
*/
#include <world_builder/features/fault.h>
#include <world_builder/utilities.h>
#include <world_builder/assert.h>
#include <world_builder/nan.h>
#include <world_builder/parameters.h>
#include <world_builder/types/array.h>
#include <world_builder/types/double.h>
#include <world_builder/types/segment.h>
#include <world_builder/types/point.h>
#include <world_builder/types/string.h>
#include <world_builder/types/unsigned_int.h>
#include <world_builder/types/plugin_system.h>
#include "rapidjson/prettywriter.h"
#include "rapidjson/stringbuffer.h"
#include "glm/glm.h"
using namespace std;
namespace WorldBuilder
{
using namespace Utilities;
namespace Features
{
Fault::Fault(WorldBuilder::World *world_)
:
reference_point(0,0,cartesian)
{
this->world = world_;
this->name = "fault";
}
Fault::~Fault()
{ }
void
Fault::declare_entries(Parameters &prm,
const std::string &parent_name,
const std::vector<std::string> & /*required_entries*/)
{
// This statment is needed because of the recursion associated with
// the sections entry.
if (parent_name == "items")
prm.enter_subsection("properties");
prm.declare_entry("min depth", Types::Double(0),
"The depth to which this feature is present");
prm.declare_entry("max depth", Types::Double(std::numeric_limits<double>::max()),
"The depth to which this feature is present");
prm.declare_entry("dip point", Types::Point<2>(),
"The depth to which this feature is present");
prm.declare_entry("segments", Types::Array(Types::Segment(0,Point<2>(0,0,invalid),Point<2>(0,0,invalid),Point<2>(0,0,invalid),
Types::PluginSystem("", Features::FaultModels::Temperature::Interface::declare_entries, {"model"}),
Types::PluginSystem("", Features::FaultModels::Composition::Interface::declare_entries, {"model"}),
Types::PluginSystem("", Features::FaultModels::Grains::Interface::declare_entries, {"model"}))),
"The depth to which this feature is present");
prm.declare_entry("temperature models",
Types::PluginSystem("", Features::FaultModels::Temperature::Interface::declare_entries, {"model"}),
"A list of temperature models.");
prm.declare_entry("composition models",
Types::PluginSystem("", Features::FaultModels::Composition::Interface::declare_entries, {"model"}),
"A list of composition models.");
prm.declare_entry("grains models",
Types::PluginSystem("", Features::FaultModels::Grains::Interface::declare_entries, {"model"}),
"A list of grains models.");
if (parent_name != "items")
{
// This only happens if we are not in sections
prm.declare_entry("sections", Types::Array(Types::PluginSystem("",Features::Fault::declare_entries, {"coordinate"}, false)),"A list of feature properties for a coordinate.");
}
else
{
// this only happens in sections
prm.declare_entry("coordinate", Types::UnsignedInt(0),
"The coordinate which should be overwritten");
prm.leave_subsection();
}
}
void
Fault::parse_entries(Parameters &prm)
{
const CoordinateSystem coordinate_system = prm.coordinate_system->natural_coordinate_system();
this->name = prm.get<std::string>("name");
this->get_coordinates("coordinates", prm, coordinate_system);
starting_depth = prm.get<double>("min depth");
maximum_depth = prm.get<double>("max depth");
const size_t n_sections = this->original_number_of_coordinates;
reference_point = prm.get<Point<2> >("dip point");
default_temperature_models.resize(0);
default_composition_models.resize(0);
default_grains_models.resize(0);
prm.get_shared_pointers<Features::FaultModels::Temperature::Interface>("temperature models", default_temperature_models);
prm.get_shared_pointers<Features::FaultModels::Composition::Interface>("composition models", default_composition_models);
prm.get_shared_pointers<Features::FaultModels::Grains::Interface>("grains models", default_grains_models);
// get the default segments.
default_segment_vector = prm.get_vector<Objects::Segment<Features::FaultModels::Temperature::Interface,
Features::FaultModels::Composition::Interface,
Features::FaultModels::Grains::Interface> >("segments", default_temperature_models, default_composition_models,default_grains_models);
// This vector stores segments to this coordiante/section.
// First used (raw) pointers to the segment relevant to this coordinate/section,
// but I do not trust it won't fail when memory is moved. So storing the all the data now.
segment_vector.resize(0);
segment_vector.resize(n_sections, default_segment_vector);
// now search whether a section is present, if so, replace the default segments.
std::vector<std::unique_ptr<Features::Fault> > sections_vector;
prm.get_unique_pointers("sections", sections_vector);
prm.enter_subsection("sections");
for (unsigned int i_section = 0; i_section < n_sections; ++i_section)
{
// first check whether this section/coordinate has a a special overwrite
for (unsigned int i_sector = 0; i_sector < sections_vector.size(); ++i_sector)
{
prm.enter_subsection(std::to_string(i_sector));
{
const unsigned int change_coord_number = prm.get<unsigned int>("coordinate");
WBAssertThrow(segment_vector.size() > change_coord_number, "Error: for subducting plate with name: '" << this->name
<< "', trying to change the section of coordinate " << change_coord_number
<< " while only " << segment_vector.size() << " coordinates are defined.");
std::vector<std::shared_ptr<Features::FaultModels::Temperature::Interface> > local_default_temperature_models;
std::vector<std::shared_ptr<Features::FaultModels::Composition::Interface> > local_default_composition_models;
std::vector<std::shared_ptr<Features::FaultModels::Grains::Interface> > local_default_grains_models;
if (prm.get_shared_pointers<Features::FaultModels::Temperature::Interface>("temperature models", local_default_temperature_models) == false)
{
// no local temperature model, use global default
local_default_temperature_models = default_temperature_models;
}
if (prm.get_shared_pointers<Features::FaultModels::Composition::Interface>("composition models", local_default_composition_models) == false)
{
// no local composition model, use global default
local_default_composition_models = default_composition_models;
}
if (prm.get_shared_pointers<Features::FaultModels::Grains::Interface>("grains models", local_default_grains_models) == false)
{
// no local grains model, use global default
local_default_grains_models = default_grains_models;
}
segment_vector[change_coord_number] = prm.get_vector<Objects::Segment<Features::FaultModels::Temperature::Interface,
Features::FaultModels::Composition::Interface,
Features::FaultModels::Grains::Interface> >("segments", local_default_temperature_models, local_default_composition_models, local_default_grains_models);
WBAssertThrow(segment_vector[change_coord_number].size() == default_segment_vector.size(),
"Error: There are not the same amount of segments in section with coordinate " << change_coord_number
<< " (" << segment_vector[change_coord_number].size() << " segments) as in the default segment ("
<< default_segment_vector.size() << " segments). This is not allowed.");
prm.enter_subsection("segments");
{
for (unsigned int i = 0; i < segment_vector[change_coord_number].size(); ++i)
{
prm.enter_subsection(std::to_string(i));
{
prm.enter_subsection("temperature models");
{
for (unsigned int j = 0; j < segment_vector[change_coord_number][i].temperature_systems.size(); ++j)
{
prm.enter_subsection(std::to_string(j));
{
segment_vector[change_coord_number][i].temperature_systems[j]->parse_entries(prm);
}
prm.leave_subsection();
}
}
prm.leave_subsection();
prm.enter_subsection("composition models");
{
for (unsigned int j = 0; j < segment_vector[change_coord_number][i].composition_systems.size(); ++j)
{
prm.enter_subsection(std::to_string(j));
{
segment_vector[change_coord_number][i].composition_systems[j]->parse_entries(prm);
}
prm.leave_subsection();
}
}
prm.leave_subsection();
prm.enter_subsection("grains models");
{
for (unsigned int j = 0; j < segment_vector[change_coord_number][i].grains_systems.size(); ++j)
{
prm.enter_subsection(std::to_string(j));
{
segment_vector[change_coord_number][i].grains_systems[j]->parse_entries(prm);
}
prm.leave_subsection();
}
}
prm.leave_subsection();
}
prm.leave_subsection();
}
}
prm.leave_subsection();
}
prm.leave_subsection();
}
}
prm.leave_subsection();
prm.enter_subsection("segments");
{
for (unsigned int i = 0; i < default_segment_vector.size(); ++i)
{
prm.enter_subsection(std::to_string(i));
{
prm.enter_subsection("temperature models");
{
for (unsigned int j = 0; j < default_segment_vector[i].temperature_systems.size(); ++j)
{
prm.enter_subsection(std::to_string(j));
{
default_segment_vector[i].temperature_systems[j]->parse_entries(prm);
}
prm.leave_subsection();
}
}
prm.leave_subsection();
prm.enter_subsection("composition models");
{
for (unsigned int j = 0; j < default_segment_vector[i].composition_systems.size(); ++j)
{
prm.enter_subsection(std::to_string(j));
{
default_segment_vector[i].composition_systems[j]->parse_entries(prm);
}
prm.leave_subsection();
}
}
prm.leave_subsection();
prm.enter_subsection("grains models");
{
for (unsigned int j = 0; j < default_segment_vector[i].grains_systems.size(); ++j)
{
prm.enter_subsection(std::to_string(j));
{
default_segment_vector[i].grains_systems[j]->parse_entries(prm);
}
prm.leave_subsection();
}
}
prm.leave_subsection();
}
prm.leave_subsection();
}
}
prm.leave_subsection();
maximum_slab_thickness = 0;
maximum_total_slab_length = 0;
total_slab_length.resize(original_number_of_coordinates);
slab_segment_lengths.resize(original_number_of_coordinates);
slab_segment_thickness.resize(original_number_of_coordinates);
slab_segment_top_truncation.resize(original_number_of_coordinates);
slab_segment_angles.resize(original_number_of_coordinates);
for (unsigned int i = 0; i < segment_vector.size(); ++i)
{
double local_total_slab_length = 0;
slab_segment_lengths[i].resize(segment_vector[i].size());
slab_segment_thickness[i].resize(segment_vector[i].size(), Point<2>(invalid));
slab_segment_top_truncation[i].resize(segment_vector[i].size(), Point<2>(invalid));
slab_segment_angles[i].resize(segment_vector[i].size(), Point<2>(invalid));
for (unsigned int j = 0; j < segment_vector[i].size(); ++j)
{
slab_segment_lengths[i][j] = segment_vector[i][j].value_length;
local_total_slab_length += segment_vector[i][j].value_length;
slab_segment_thickness[i][j] = segment_vector[i][j].value_thickness;
slab_segment_top_truncation[i][j] = segment_vector[i][j].value_top_truncation;
slab_segment_angles[i][j] = segment_vector[i][j].value_angle * (const_pi/180);
}
maximum_slab_thickness = std::max(maximum_slab_thickness, local_total_slab_length);
total_slab_length[i] = local_total_slab_length;
maximum_total_slab_length = std::max(maximum_total_slab_length, local_total_slab_length);
}
}
double
Fault::temperature(const Point<3> &position,
const double depth,
const double gravity_norm,
double temperature) const
{
WorldBuilder::Utilities::NaturalCoordinate natural_coordinate = WorldBuilder::Utilities::NaturalCoordinate(position,
*(world->parameters.coordinate_system));
// The depth variable is the distance from the surface to the position, the depth
// coordinate is the distance from the bottom of the model to the position and
// the starting radius is the distance from the bottom of the model to the surface.
const double starting_radius = natural_coordinate.get_depth_coordinate() + depth - starting_depth;
WBAssert(std::abs(starting_radius) > std::numeric_limits<double>::epsilon(), "World Builder error: starting_radius can not be zero. "
<< "Position = " << position[0] << ":" << position[1] << ":" << position[2]
<< ", natural_coordinate.get_depth_coordinate() = " << natural_coordinate.get_depth_coordinate()
<< ", depth = " << depth
<< ", starting_depth " << starting_depth
);
// todo: explain and check -starting_depth
if (depth <= maximum_depth && depth >= starting_depth && depth <= maximum_total_slab_length + maximum_slab_thickness)
{
// todo: explain
// This function only returns positive values, because we want
// the fault to be centered around the line provided by the user.
std::map<std::string,double> distance_from_planes =
WorldBuilder::Utilities::distance_point_from_curved_planes(position,
reference_point,
coordinates,
slab_segment_lengths,
slab_segment_angles,
starting_radius,
this->world->parameters.coordinate_system,
true,
one_dimensional_coordinates);
const double distance_from_plane = distance_from_planes["distanceFromPlane"];
const double distance_along_plane = distance_from_planes["distanceAlongPlane"];
const double section_fraction = distance_from_planes["sectionFraction"];
const size_t current_section = static_cast<size_t>(std::floor(one_dimensional_coordinates[static_cast<size_t>(distance_from_planes["section"])]));
const size_t next_section = current_section + 1;
const size_t current_segment = static_cast<size_t>(distance_from_planes["segment"]); // the original value was a unsigned in, converting it back.
//const size_t next_segment = current_segment + 1;
const double segment_fraction = distance_from_planes["segmentFraction"];
if (abs(distance_from_plane) < INFINITY || (distance_along_plane) < INFINITY)
{
// We want to do both section (horizontal) and segment (vertical) interpolation.
const double thickness_up = slab_segment_thickness[current_section][current_segment][0]
+ section_fraction
* (slab_segment_thickness[next_section][current_segment][0]
- slab_segment_thickness[current_section][current_segment][0]);
const double thickness_down = slab_segment_thickness[current_section][current_segment][1]
+ section_fraction
* (slab_segment_thickness[next_section][current_segment][1]
- slab_segment_thickness[current_section][current_segment][1]);
const double thickness_local = thickness_up + segment_fraction * (thickness_down - thickness_up);
const double max_slab_length = total_slab_length[current_section] +
section_fraction *
(total_slab_length[next_section] - total_slab_length[current_section]);
// secondly for top truncation
const double top_truncation_up = slab_segment_top_truncation[current_section][current_segment][0]
+ section_fraction
* (slab_segment_top_truncation[next_section][current_segment][0]
- slab_segment_top_truncation[current_section][current_segment][0]);
const double top_truncation_down = slab_segment_top_truncation[current_section][current_segment][1]
+ section_fraction
* (slab_segment_top_truncation[next_section][current_segment][1]
- slab_segment_top_truncation[current_section][current_segment][1]);
const double top_truncation_local = top_truncation_up + segment_fraction * (top_truncation_down - top_truncation_up);
// if the thickness is zero, we don't need to compute anything, so return.
if (std::fabs(thickness_local) < 2.0 * std::numeric_limits<double>::epsilon())
return temperature;
// if the thickness is smaller than what is truncated off at the top, we don't need to compute anything, so return.
if (thickness_local < top_truncation_local)
return temperature;
// Because both sides return positve values, we have to
// devide the thickness_local by two
if (std::fabs(distance_from_plane) > 0 &&
std::fabs(distance_from_plane) <= thickness_local * 0.5 &&
distance_along_plane > 0 &&
distance_along_plane <= max_slab_length)
{
// Inside the fault!
double temperature_current_section = temperature;
double temperature_next_section = temperature;
for (auto &temperature_model: segment_vector[current_section][current_segment].temperature_systems)
{
temperature_current_section = temperature_model->get_temperature(position,
depth,
gravity_norm,
temperature_current_section,
starting_depth,
maximum_depth,
distance_from_planes);
WBAssert(!std::isnan(temperature_current_section), "Temparture is not a number: " << temperature_current_section
<< ", based on a temperature model with the name " << temperature_model->get_name());
WBAssert(std::isfinite(temperature_current_section), "Temparture is not a finite: " << temperature_current_section
<< ", based on a temperature model with the name " << temperature_model->get_name());
}
for (auto &temperature_model: segment_vector[next_section][current_segment].temperature_systems)
{
temperature_next_section = temperature_model->get_temperature(position,
depth,
gravity_norm,
temperature_next_section,
starting_depth,
maximum_depth,
distance_from_planes);
WBAssert(!std::isnan(temperature_next_section), "Temparture is not a number: " << temperature_next_section
<< ", based on a temperature model with the name " << temperature_model->get_name());
WBAssert(std::isfinite(temperature_next_section), "Temparture is not a finite: " << temperature_next_section
<< ", based on a temperature model with the name " << temperature_model->get_name());
}
// linear interpolation between current and next section temperatures
temperature = temperature_current_section + section_fraction * (temperature_next_section - temperature_current_section);
}
}
}
return temperature;
}
double
Fault::composition(const Point<3> &position,
const double depth,
const unsigned int composition_number,
double composition) const
{
WorldBuilder::Utilities::NaturalCoordinate natural_coordinate = WorldBuilder::Utilities::NaturalCoordinate(position,
*(world->parameters.coordinate_system));
// todo: explain
const double starting_radius = natural_coordinate.get_depth_coordinate() + depth - starting_depth;
// todo: explain and check -starting_depth
if (depth <= maximum_depth && depth >= starting_depth && depth <= maximum_total_slab_length + maximum_slab_thickness)
{
// todo: explain
// This function only returns positive values, because we want
// the fault to be centered around the line provided by the user.
std::map<std::string,double> distance_from_planes =
WorldBuilder::Utilities::distance_point_from_curved_planes(position,
reference_point,
coordinates,
slab_segment_lengths,
slab_segment_angles,
starting_radius,
this->world->parameters.coordinate_system,
true,
one_dimensional_coordinates);
const double distance_from_plane = distance_from_planes["distanceFromPlane"];
const double distance_along_plane = distance_from_planes["distanceAlongPlane"];
const double section_fraction = distance_from_planes["sectionFraction"];
const size_t current_section = static_cast<size_t>(std::floor(one_dimensional_coordinates[static_cast<size_t>(distance_from_planes["section"])]));
const size_t next_section = current_section + 1;
const size_t current_segment = static_cast<size_t>(distance_from_planes["segment"]); // the original value was a unsigned int, turning it back.
//const size_t next_segment = current_segment + 1;
const double segment_fraction = distance_from_planes["segmentFraction"];
if (abs(distance_from_plane) < INFINITY || (distance_along_plane) < INFINITY)
{
// We want to do both section (horizontal) and segment (vertical) interpolation.
const double thickness_up = slab_segment_thickness[current_section][current_segment][0]
+ section_fraction
* (slab_segment_thickness[next_section][current_segment][0]
- slab_segment_thickness[current_section][current_segment][0]);
const double thickness_down = slab_segment_thickness[current_section][current_segment][1]
+ section_fraction
* (slab_segment_thickness[next_section][current_segment][1]
- slab_segment_thickness[current_section][current_segment][1]);
const double thickness_local = thickness_up + segment_fraction * (thickness_down - thickness_up);
// secondly for top truncation
const double top_truncation_up = slab_segment_top_truncation[current_section][current_segment][0]
+ section_fraction
* (slab_segment_top_truncation[next_section][current_segment][0]
- slab_segment_top_truncation[current_section][current_segment][0]);
const double top_truncation_down = slab_segment_top_truncation[current_section][current_segment][1]
+ section_fraction
* (slab_segment_top_truncation[next_section][current_segment][1]
- slab_segment_top_truncation[current_section][current_segment][1]);
const double top_truncation_local = top_truncation_up + segment_fraction * (top_truncation_down - top_truncation_up);
// if the thickness is zero, we don't need to compute anything, so return.
if (std::fabs(thickness_local) < 2.0 * std::numeric_limits<double>::epsilon())
return composition;
// if the thickness is smaller than what is truncated off at the top, we don't need to compute anything, so return.
if (thickness_local < top_truncation_local)
return composition;
const double max_slab_length = total_slab_length[current_section] +
section_fraction *
(total_slab_length[next_section] - total_slab_length[current_section]);
// Because both sides return positve values, we have to
// devide the thickness_local by two
if (std::fabs(distance_from_plane) > 0 &&
std::fabs(distance_from_plane) <= thickness_local * 0.5 &&
distance_along_plane > 0 &&
distance_along_plane <= max_slab_length)
{
// Inside the fault!
double composition_current_section = composition;
double composition_next_section = composition;
for (auto &composition_model: segment_vector[current_section][current_segment].composition_systems)
{
composition_current_section = composition_model->get_composition(position,
depth,
composition_number,
composition_current_section,
starting_depth,
maximum_depth,
distance_from_planes);
WBAssert(!std::isnan(composition_current_section), "Composition_current_section is not a number: " << composition_current_section
<< ", based on a temperature model with the name " << composition_model->get_name());
WBAssert(std::isfinite(composition_current_section), "Composition_current_section is not a finite: " << composition_current_section
<< ", based on a temperature model with the name " << composition_model->get_name());
}
for (auto &composition_model: segment_vector[next_section][current_segment].composition_systems)
{
composition_next_section = composition_model->get_composition(position,
depth,
composition_number,
composition_next_section,
starting_depth,
maximum_depth,
distance_from_planes);
WBAssert(!std::isnan(composition_next_section), "Composition_next_section is not a number: " << composition_next_section
<< ", based on a temperature model with the name " << composition_model->get_name());
WBAssert(std::isfinite(composition_next_section), "Composition_next_section is not a finite: " << composition_next_section
<< ", based on a temperature model with the name " << composition_model->get_name());
}
// linear interpolation between current and next section temperatures
composition = composition_current_section + section_fraction * (composition_next_section - composition_current_section);
}
}
}
return composition;
}
WorldBuilder::grains
Fault::grains(const Point<3> &position,
const double depth,
const unsigned int composition_number,
WorldBuilder::grains grains) const
{
WorldBuilder::Utilities::NaturalCoordinate natural_coordinate = WorldBuilder::Utilities::NaturalCoordinate(position,
*(world->parameters.coordinate_system));
// todo: explain
const double starting_radius = natural_coordinate.get_depth_coordinate() + depth - starting_depth;
// todo: explain and check -starting_depth
if (depth <= maximum_depth && depth >= starting_depth && depth <= maximum_total_slab_length + maximum_slab_thickness)
{
// todo: explain
// This function only returns positive values, because we want
// the fault to be centered around the line provided by the user.
std::map<std::string,double> distance_from_planes =
WorldBuilder::Utilities::distance_point_from_curved_planes(position,
reference_point,
coordinates,
slab_segment_lengths,
slab_segment_angles,
starting_radius,
this->world->parameters.coordinate_system,
true,
one_dimensional_coordinates);
const double distance_from_plane = distance_from_planes["distanceFromPlane"];
const double distance_along_plane = distance_from_planes["distanceAlongPlane"];
const double section_fraction = distance_from_planes["sectionFraction"];
const size_t current_section = static_cast<size_t>(std::floor(one_dimensional_coordinates[static_cast<size_t>(distance_from_planes["section"])]));
const size_t next_section = current_section + 1;
const size_t current_segment = static_cast<size_t>(distance_from_planes["segment"]); // the original value was a unsigned int, turning it back.
//const size_t next_segment = current_segment + 1;
const double segment_fraction = distance_from_planes["segmentFraction"];
if (abs(distance_from_plane) < INFINITY || (distance_along_plane) < INFINITY)
{
// We want to do both section (horizontal) and segment (vertical) interpolation.
const double thickness_up = slab_segment_thickness[current_section][current_segment][0]
+ section_fraction
* (slab_segment_thickness[next_section][current_segment][0]
- slab_segment_thickness[current_section][current_segment][0]);
const double thickness_down = slab_segment_thickness[current_section][current_segment][1]
+ section_fraction
* (slab_segment_thickness[next_section][current_segment][1]
- slab_segment_thickness[current_section][current_segment][1]);
const double thickness_local = thickness_up + segment_fraction * (thickness_down - thickness_up);
// secondly for top truncation
const double top_truncation_up = slab_segment_top_truncation[current_section][current_segment][0]
+ section_fraction
* (slab_segment_top_truncation[next_section][current_segment][0]
- slab_segment_top_truncation[current_section][current_segment][0]);
const double top_truncation_down = slab_segment_top_truncation[current_section][current_segment][1]
+ section_fraction
* (slab_segment_top_truncation[next_section][current_segment][1]
- slab_segment_top_truncation[current_section][current_segment][1]);
const double top_truncation_local = top_truncation_up + segment_fraction * (top_truncation_down - top_truncation_up);
// if the thickness is zero, we don't need to compute anything, so return.
if (std::fabs(thickness_local) < 2.0 * std::numeric_limits<double>::epsilon())
return grains;
// if the thickness is smaller than what is truncated off at the top, we don't need to compute anything, so return.
if (thickness_local < top_truncation_local)
return grains;
const double max_slab_length = total_slab_length[current_section] +
section_fraction *
(total_slab_length[next_section] - total_slab_length[current_section]);
// Because both sides return positve values, we have to
// devide the thickness_local by two
if (std::fabs(distance_from_plane) > 0 &&
std::fabs(distance_from_plane) <= thickness_local * 0.5 &&
distance_along_plane > 0 &&
distance_along_plane <= max_slab_length)
{
// Inside the fault!
WorldBuilder::grains grains_current_section = grains;
WorldBuilder::grains grains_next_section = grains;
for (auto &grains_model: segment_vector[current_section][current_segment].grains_systems)
{
grains_current_section = grains_model->get_grains(position,
depth,
composition_number,
grains_current_section,
starting_depth,
maximum_depth,
distance_from_planes);
/*WBAssert(!std::isnan(composition_current_section), "Composition_current_section is not a number: " << composition_current_section
<< ", based on a temperature model with the name " << composition_model->get_name());
WBAssert(std::isfinite(composition_current_section), "Composition_current_section is not a finite: " << composition_current_section
<< ", based on a temperature model with the name " << composition_model->get_name());*/
}
for (auto &grains_model: segment_vector[next_section][current_segment].grains_systems)
{
grains_next_section = grains_model->get_grains(position,
depth,
composition_number,
grains_next_section,
starting_depth,
maximum_depth,
distance_from_planes);
/*WBAssert(!std::isnan(composition_next_section), "Composition_next_section is not a number: " << composition_next_section
<< ", based on a temperature model with the name " << composition_model->get_name());
WBAssert(std::isfinite(composition_next_section), "Composition_next_section is not a finite: " << composition_next_section
<< ", based on a temperature model with the name " << composition_model->get_name());*/
}
// linear interpolation between current and next section temperatures
for (size_t i = 0; i < grains.sizes.size(); i++)
{
grains.sizes[i] = grains_current_section.sizes[i] + section_fraction * (grains_next_section.sizes[i] - grains_current_section.sizes[i]);
}
// average two rotations matrices throu quaternions.
for (size_t i = 0; i < grains_current_section.rotation_matrices.size(); i++)
{
glm::quaternion::quat quat_current = glm::quaternion::quat_cast(grains_current_section.rotation_matrices[i]);
glm::quaternion::quat quat_next = glm::quaternion::quat_cast(grains_next_section.rotation_matrices[i]);
glm::quaternion::quat quat_average = glm::quaternion::slerp(quat_current,quat_next,section_fraction);
grains.rotation_matrices[i] = glm::quaternion::mat3_cast(quat_average);
}
}
}
}
// todo: fix to an averged value
return grains;
}
/**
* Register plugin
*/
WB_REGISTER_FEATURE(Fault, fault)
}
}