-
Notifications
You must be signed in to change notification settings - Fork 265
Expand file tree
/
Copy pathhelper_functions.cc
More file actions
2885 lines (2493 loc) · 141 KB
/
helper_functions.cc
File metadata and controls
2885 lines (2493 loc) · 141 KB
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 - 2024 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/simulator.h>
#include <aspect/advection_field.h>
#include <aspect/melt.h>
#include <aspect/volume_of_fluid/handler.h>
#include <aspect/newton.h>
#include <aspect/global.h>
#include <aspect/geometry_model/interface.h>
#include <aspect/geometry_model/ellipsoidal_chunk.h>
#include <aspect/heating_model/interface.h>
#include <aspect/heating_model/adiabatic_heating.h>
#include <aspect/material_model/interface.h>
#include <aspect/mesh_deformation/interface.h>
#include <aspect/particle/generator/interface.h>
#include <aspect/particle/integrator/interface.h>
#include <aspect/particle/interpolator/interface.h>
#include <aspect/particle/property/interface.h>
#include <aspect/postprocess/visualization.h>
#include <aspect/prescribed_solution/interface.h>
#include <aspect/prescribed_stokes_solution/interface.h>
#include <deal.II/base/index_set.h>
#include <deal.II/base/conditional_ostream.h>
#include <deal.II/base/quadrature_lib.h>
#include <deal.II/base/signaling_nan.h>
#include <deal.II/lac/block_sparsity_pattern.h>
#include <deal.II/grid/grid_tools.h>
#include <deal.II/dofs/dof_renumbering.h>
#include <deal.II/dofs/dof_accessor.h>
#include <deal.II/dofs/dof_tools.h>
#include <deal.II/fe/fe_q.h>
#include <deal.II/fe/fe_dgq.h>
#include <deal.II/fe/fe_dgp.h>
#include <deal.II/fe/fe_values.h>
#include <deal.II/sundials/arkode.h>
#include <fstream>
#include <iostream>
#include <string>
namespace aspect
{
template <int dim>
void Simulator<dim>::write_plugin_graph (std::ostream &out) const
{
// write the preamble
out << "digraph Plugins\n"
"{\n"
" splines=line;\n"
" splines=true;\n"
" overlap=false;\n"
" edge [fontname=\"FreeSans\",\n"
" fontsize=\"10\",\n"
" labelfontname=\"FreeSans\",\n"
" labelfontsize=\"10\",\n"
" color=\"black\",\n"
" style=\"solid\"];\n"
" node [fontname=\"FreeSans\",\n"
" fontsize=\"10\",\n"
" shape=\"rectangle\",\n"
" height=0.2,\n"
" width=0.4,\n"
" color=\"black\",\n"
" fillcolor=\"white\",\n"
" style=\"filled\"];\n"
" layout=neato;\n"
"\n";
// then also write nodes for the Simulator and SimulatorAccess classes,
// and an arrow from the former to the latter to indicate flow of
// information
out << " Simulator [height=1.5,width=2,shape=\"octagon\",fillcolor=\"yellow\"];\n";
out << " SimulatorAccess [height=1.2,width=1.2,shape=\"rect\",fillcolor=\"yellow\"];\n";
out << " Simulator -> SimulatorAccess [len=1, weight=100];\n";
// then go through all plugin systems and output everything we have
AdiabaticConditions::write_plugin_graph<dim>(out);
BoundaryComposition::Manager<dim>::write_plugin_graph(out);
BoundaryFluidPressure::write_plugin_graph<dim>(out);
BoundaryTemperature::Manager<dim>::write_plugin_graph(out);
BoundaryTraction::Manager<dim>::write_plugin_graph(out);
BoundaryVelocity::Manager<dim>::write_plugin_graph(out);
InitialTopographyModel::write_plugin_graph<dim>(out);
GeometryModel::write_plugin_graph<dim>(out);
GravityModel::write_plugin_graph<dim>(out);
HeatingModel::Manager<dim>::write_plugin_graph(out);
InitialComposition::Manager<dim>::write_plugin_graph(out);
InitialTemperature::Manager<dim>::write_plugin_graph(out);
MaterialModel::write_plugin_graph<dim>(out);
MeshRefinement::Manager<dim>::write_plugin_graph(out);
Particle::Generator::write_plugin_graph<dim>(out);
Particle::Integrator::write_plugin_graph<dim>(out);
Particle::Interpolator::write_plugin_graph<dim>(out);
Particle::Property::Manager<dim>::write_plugin_graph(out);
Postprocess::Manager<dim>::write_plugin_graph(out);
Postprocess::Visualization<dim>::write_plugin_graph(out);
PrescribedSolution::Manager<dim>::write_plugin_graph(out);
PrescribedStokesSolution::write_plugin_graph<dim>(out);
TerminationCriteria::Manager<dim>::write_plugin_graph(out);
// end the graph
out << '}'
<< std::endl;
}
template <int dim>
void Simulator<dim>::output_statistics()
{
// Only write the statistics file from processor zero
if (Utilities::MPI::this_mpi_process(mpi_communicator)!=0)
return;
// Formatting the table we're about to output and writing the
// actual file may take some time, so do it on a separate
// thread. We do this using a lambda function that takes
// a copy of the statistics object to make sure that whatever
// we do to the 'real' statistics object at the time of
// writing data doesn't affect what we write.
//
// Before we can start working on a new thread, we need to
// make sure that the previous thread is done or they'll
// step on each other's feet.
if (output_statistics_thread.joinable())
output_statistics_thread.join();
// Write data in the background through a lambda function.
// This happening in the background means that we have
// to create a copy of the statistics table, since whatever is
// running in the foreground may continue to add entries to the
// statistics table at the same time.
output_statistics_thread = std::thread (
[statistics_copy_ptr = std::make_unique<TableHandler>(statistics),this]()
{
// First write everything into a string in memory
std::ostringstream stream;
statistics_copy_ptr->write_text (stream,
TableHandler::table_with_separate_column_description);
stream.flush();
const std::string statistics_contents = stream.str();
// Next find out whether we need to write everything into
// the statistics file, or whether it is enough to just write
// the last few bytes that were added since we wrote to that
// file again. The way we do that is by checking whether the
// first few bytes of the string we just created match what we
// had previously written. One might think that they always should,
// but the statistics object automatically sizes the column widths
// of its output to match what is being written, and so if a later
// entry requires more width, then even the first columns are
// changed -- in that case, we will have to write everything,
// not just append one line.
const bool write_everything
= ( // We may have never written anything. More precisely, this
// case happens if the statistics_last_write_size is at the
// value initialized by the Simulator::Simulator()
// constructor, and this can happen in two situations:
// (i) At the end of the first time step; and (ii) upon restart
// since the variable we query here is not serialized. It is clear
// that in both situations, we want to write the
// entire contents of the statistics object. For the second
// case, this is also appropriate since someone may have
// previously restarted from a checkpoint, run a couple of
// time steps that have added to the statistics file, but then
// aborted the run again; a later restart from the same
// checkpoint then requires overwriting the statistics file
// on disk with what we have when this function is called for
// the first time after the restart. The same situation
// happens if the simulation kept running for some time after
// a checkpoint, but is resumed from that checkpoint (i.e.,
// at an earlier time step than when the statistics file was
// written to last). In these situations, we effectively want
// to "truncate" the file to the state stored in the checkpoint,
// and we do that by just overwriting the entire file.
(statistics_last_write_size == 0)
||
// Or the size of the statistics file may have
// shrunk mysteriously -- this shouldn't happen
// but if it did we'd get into trouble with the
// .substr() call in the next check.
(statistics_last_write_size > statistics_contents.size())
||
// Or the hash of what we wrote last time doesn't match
// the hash of the first part of what we want to write
(statistics_last_hash
!=
std::hash<std::string>()(statistics_contents.substr(0, statistics_last_write_size))) );
const std::string stat_file_name = parameters.output_directory + "statistics";
if (write_everything)
{
// Write what we have into a tmp file, then move that into
// place
const std::string tmp_file_name = stat_file_name + ".tmp";
{
std::ofstream tmp_file (tmp_file_name);
tmp_file << statistics_contents;
}
std::rename(tmp_file_name.c_str(), stat_file_name.c_str());
}
else
{
// If we don't have to write everything, then the first part of what
// we want to write matches what's already on disk. In that case,
// we just have to append what's new.
std::ofstream stat_file (stat_file_name, std::ios::app);
stat_file << statistics_contents.substr(statistics_last_write_size, std::string::npos);
}
// Now update the size and hash of what we just wrote so that
// we can compare against it next time we get here. Note that we do
// not need to guard access to these variables with a mutex because
// this is the only function that touches the variables, and
// this function runs only once at a time (on a different
// thread, but it's not started a second time while the previous
// run hasn't finished).
statistics_last_write_size = statistics_contents.size();
statistics_last_hash = std::hash<std::string>()(statistics_contents);
});
}
template <int dim>
double
Simulator<dim>::
compute_pressure_scaling_factor() const
{
// Determine how to treat the pressure. We have to scale it for the solver
// to make velocities and pressures of roughly the same (numerical) size,
// and we may have to fix up the right hand side vector before solving for
// compressible models if there are no in-/outflow boundaries
//
// We do this by scaling the divergence equation by a constant factor that
// is equal to a reference viscosity divided by a length scale.
// We get the latter from the geometry model, and the former
// by looping over all cells and averaging the "order of magnitude"
// of the viscosity. The order of magnitude is the logarithm of
// the viscosity, so
//
// \eta_{ref} = exp( 1/N * (log(eta_1) + log(eta_2) + ... + log(eta_N))
//
// where the \eta_i are typical viscosities on the cells of the mesh.
// For this, we just take the viscosity at the cell center.
//
// The formula above computes the exponential of the average of the
// logarithms. It is easy to verify that this is equivalent to
// computing the *geometric* mean of the viscosities, but the
// formula above is numerically stable.
// Do return NaN if we do not solve the Stokes equation. Technically, there is
// no harm in computing the factor anyway, but the viscosity has to be positive
// for this function to work, and some models may set viscosity to 0 if not solving
// the Stokes equation. Returning NaN guarantees this value cannot be used.
if (parameters.nonlinear_solver == NonlinearSolver::no_Advection_no_Stokes ||
parameters.nonlinear_solver == NonlinearSolver::single_Advection_no_Stokes)
return numbers::signaling_nan<double>();
const QMidpoint<dim> quadrature_formula;
FEValues<dim> fe_values (*mapping,
finite_element,
quadrature_formula,
update_values |
update_gradients |
update_quadrature_points |
update_JxW_values);
const unsigned int n_q_points = quadrature_formula.size();
double local_integrated_viscosity_logarithm = 0.;
double local_volume = 0.;
MaterialModel::MaterialModelInputs<dim> in(n_q_points,
introspection.n_compositional_fields);
MaterialModel::MaterialModelOutputs<dim> out(n_q_points,
introspection.n_compositional_fields);
// We do not need to compute anything but the viscosity
in.requested_properties = MaterialModel::MaterialProperties::viscosity;
for (const auto &cell : dof_handler.active_cell_iterators())
if (cell->is_locally_owned())
{
fe_values.reinit (cell);
in.reinit(fe_values,
cell,
introspection,
solution);
// We do not call the cell-wise average function of the
// material model, because we average globally below
material_model->evaluate(in, out);
// Evaluate viscosity at the mid-point of each cell and
// calculate the volume weighted harmonic average of all cells
for (unsigned int q=0; q<n_q_points; ++q)
{
Assert(out.viscosities[q] > 0,
ExcMessage ("The viscosity needs to be a "
"positive quantity."));
const double JxW = fe_values.JxW(q);
local_integrated_viscosity_logarithm += std::log(out.viscosities[q]) * JxW;
local_volume += JxW;
}
}
// vector for packing local values before MPI summing them
double values[2] = {local_integrated_viscosity_logarithm, local_volume};
Utilities::MPI::sum(values, mpi_communicator, values);
const double reference_viscosity = std::exp(values[0]/values[1]);
const double length_scale = geometry_model->length_scale();
// Allow the user to inspect and/or overwrite our result:
const double result = reference_viscosity / length_scale;
if (boost::optional<double> value = signals.modify_pressure_scaling(result, reference_viscosity, length_scale))
return *value;
else
return result;
}
template <int dim>
double
Simulator<dim>::
get_maximal_velocity (const LinearAlgebra::BlockVector &solution) const
{
// use a quadrature formula that has one point at
// the location of each degree of freedom in the
// velocity element
const QIterated<dim> quadrature_formula (QTrapezoid<1>(),
parameters.stokes_velocity_degree);
const unsigned int n_q_points = quadrature_formula.size();
FEValues<dim> fe_values (*mapping, finite_element, quadrature_formula, update_values);
std::vector<Tensor<1,dim>> velocity_values(n_q_points);
double max_local_velocity = 0;
// loop over all locally owned cells and evaluate the velocities at each
// quadrature point (i.e. each node). keep a running tally of the largest
// such velocity
for (const auto &cell : dof_handler.active_cell_iterators())
if (cell->is_locally_owned())
{
fe_values.reinit (cell);
fe_values[introspection.extractors.velocities].get_function_values (solution,
velocity_values);
for (unsigned int q=0; q<n_q_points; ++q)
max_local_velocity = std::max (max_local_velocity,
velocity_values[q].norm());
}
// return the largest value over all processors
return Utilities::MPI::max (max_local_velocity, mpi_communicator);
}
template <int dim>
bool Simulator<dim>::maybe_do_initial_refinement (const unsigned int max_refinement_level)
{
if (pre_refinement_step < parameters.initial_adaptive_refinement)
{
if (parameters.timing_output_frequency == 0)
{
computing_timer.print_summary ();
pcout << "-- Total wallclock time elapsed including restarts: "
<< std::round(wall_timer.wall_time()+total_walltime_until_last_snapshot)
<< 's' << std::endl;
}
output_statistics();
// we only want to do the postprocessing here if it is not already done in
// the nonlinear iteration scheme, which is the case if we run postprocessors
// on all nonlinear iterations
if (parameters.run_postprocessors_on_initial_refinement && (!parameters.run_postprocessors_on_nonlinear_iterations))
postprocess ();
refine_mesh (max_refinement_level);
++pre_refinement_step;
return true;
}
else
{
// invalidate the value of pre_refinement_step since it will no longer be used from here on
pre_refinement_step = std::numeric_limits<unsigned int>::max();
return false;
}
}
template <int dim>
void Simulator<dim>::exchange_refinement_flags ()
{
// Communicate refinement flags on ghost cells from the owner of the
// cell. This is necessary to get consistent refinement, as mesh
// smoothing would undo some of the requested coarsening/refinement.
auto pack
= [] (const typename DoFHandler<dim>::active_cell_iterator &cell) -> std::uint8_t
{
if (cell->refine_flag_set())
return 1;
if (cell->coarsen_flag_set())
return 2;
return 0;
};
auto unpack
= [] (const typename DoFHandler<dim>::active_cell_iterator &cell, const std::uint8_t &flag) -> void
{
cell->clear_coarsen_flag();
cell->clear_refine_flag();
if (flag==1)
cell->set_refine_flag();
else if (flag==2)
cell->set_coarsen_flag();
};
GridTools::exchange_cell_data_to_ghosts<std::uint8_t, DoFHandler<dim>>
(dof_handler, pack, unpack);
}
template <int dim>
void Simulator<dim>::maybe_refine_mesh (const double new_time_step,
unsigned int &max_refinement_level)
{
/*
* see if this is an additional refinement cycle. An additional refinement
* cycle differs from a regular, because the maximal refinement level allowed
* is increased by one from this time on.
*/
if ((parameters.additional_refinement_times.size() > 0)
&&
(parameters.additional_refinement_times.front () < time+new_time_step))
{
// loop over as many times as this is necessary
while ((parameters.additional_refinement_times.size() > 0)
&&
(parameters.additional_refinement_times.front () < time+new_time_step))
{
++max_refinement_level;
refine_mesh (max_refinement_level);
parameters.additional_refinement_times
.erase (parameters.additional_refinement_times.begin());
}
}
// see if this is a time step where regular refinement is requested
else if ((timestep_number > 0
&&
parameters.adaptive_refinement_interval > 0
&&
timestep_number % parameters.adaptive_refinement_interval == 0)
||
(timestep_number == 0 && parameters.adaptive_refinement_interval == 1)
)
{
refine_mesh (max_refinement_level);
}
}
template <int dim>
void Simulator<dim>::maybe_write_timing_output () const
{
bool write_timing_output = false;
if (parameters.timing_output_frequency <= 1)
write_timing_output = true;
else if ((timestep_number > 0) &&
(timestep_number % parameters.timing_output_frequency == 0))
write_timing_output = true;
// if requested output a summary of the current timing information
if (write_timing_output)
{
computing_timer.print_summary ();
pcout << "-- Total wallclock time elapsed including restarts: "
<< std::round(wall_timer.wall_time()+total_walltime_until_last_snapshot)
<< 's' << std::endl;
}
}
template <int dim>
bool Simulator<dim>::maybe_write_checkpoint (const std::time_t last_checkpoint_time,
const bool force_writing_checkpoint)
{
// Do a checkpoint if this is the end of simulation,
// and the termination criteria say to checkpoint at the end.
bool write_checkpoint = force_writing_checkpoint;
// If we base checkpoint frequency on timing, measure the time at process 0
// This prevents race conditions where some processes will checkpoint and others won't
if (!write_checkpoint && parameters.checkpoint_time_secs > 0)
{
const bool global_do_checkpoint = Utilities::MPI::broadcast(mpi_communicator,
(std::time(nullptr)-last_checkpoint_time) >=
parameters.checkpoint_time_secs,
/* root= */0);
if (global_do_checkpoint)
write_checkpoint = true;
}
// If we base checkpoint frequency on steps, see if it's time for another checkpoint
if (!write_checkpoint &&
(parameters.checkpoint_time_secs == 0) &&
(parameters.checkpoint_steps > 0) &&
(timestep_number % parameters.checkpoint_steps == 0))
write_checkpoint = true;
// Do a checkpoint if indicated by checkpoint parameters
if (write_checkpoint)
{
create_snapshot();
// matrices will be regenerated after a resume, so do that here too
// to be consistent. otherwise we would get different results
// for a restarted computation than for one that ran straight
// through
rebuild_stokes_matrix =
rebuild_stokes_preconditioner = true;
}
return write_checkpoint;
}
template <int dim>
void Simulator<dim>::advance_time (const double step_size)
{
old_time_step = time_step;
time_step = step_size;
time += time_step;
++timestep_number;
// prepare for the next time step by shifting solution vectors
// by one time step. In timestep 0 (just increased in the
// line above) initialize both old_solution
// and old_old_solution with the currently computed solution.
if (timestep_number == 1)
{
old_old_solution = solution;
old_solution = solution;
}
else
{
old_old_solution = old_solution;
old_solution = solution;
}
}
template <int dim>
std::pair<double,double>
Simulator<dim>::
get_extrapolated_advection_field_range (const AdvectionField &advection_field) const
{
const QIterated<dim> quadrature_formula (QTrapezoid<1>(),
advection_field.polynomial_degree(introspection));
const unsigned int n_q_points = quadrature_formula.size();
const FEValuesExtractors::Scalar field = advection_field.scalar_extractor(introspection);
FEValues<dim> fe_values (*mapping, finite_element, quadrature_formula,
update_values);
std::vector<double> old_field_values(n_q_points);
std::vector<double> old_old_field_values(n_q_points);
// This presets the minimum with a bigger
// and the maximum with a smaller number
// than one that is going to appear. Will
// be overwritten in the cell loop or in
// the communication step at the
// latest.
double min_local_field = std::numeric_limits<double>::max(),
max_local_field = std::numeric_limits<double>::lowest();
if (timestep_number > 1)
{
for (const auto &cell : dof_handler.active_cell_iterators())
if (cell->is_locally_owned())
{
fe_values.reinit (cell);
fe_values[field].get_function_values (old_solution,
old_field_values);
fe_values[field].get_function_values (old_old_solution,
old_old_field_values);
for (unsigned int q=0; q<n_q_points; ++q)
{
const double extrapolated_field =
(1. + time_step/old_time_step) * old_field_values[q]-
time_step/old_time_step * old_old_field_values[q];
min_local_field = std::min (min_local_field,
extrapolated_field);
max_local_field = std::max (max_local_field,
extrapolated_field);
}
}
}
else
{
for (const auto &cell : dof_handler.active_cell_iterators())
if (cell->is_locally_owned())
{
fe_values.reinit (cell);
fe_values[field].get_function_values (old_solution,
old_field_values);
for (unsigned int q=0; q<n_q_points; ++q)
{
const double extrapolated_field = old_field_values[q];
min_local_field = std::min (min_local_field,
extrapolated_field);
max_local_field = std::max (max_local_field,
extrapolated_field);
}
}
}
return std::make_pair(Utilities::MPI::min (min_local_field,
mpi_communicator),
Utilities::MPI::max (max_local_field,
mpi_communicator));
}
template <int dim>
void Simulator<dim>::interpolate_onto_velocity_system(const TensorFunction<1,dim> &func,
LinearAlgebra::Vector &vec) const
{
Assert(introspection.block_indices.velocities == 0, ExcNotImplemented());
const std::vector<Point<dim>> mesh_support_points = finite_element.base_element(introspection.base_elements.velocities).get_unit_support_points();
const unsigned int n_velocity_dofs_per_cell = finite_element.base_element(introspection.base_elements.velocities).dofs_per_cell;
FEValues<dim> mesh_points (*mapping, finite_element, mesh_support_points, update_quadrature_points);
std::vector<types::global_dof_index> cell_dof_indices (finite_element.dofs_per_cell);
for (const auto &cell : dof_handler.active_cell_iterators())
if (cell->is_locally_owned())
{
mesh_points.reinit(cell);
cell->get_dof_indices (cell_dof_indices);
for (unsigned int j=0; j<n_velocity_dofs_per_cell; ++j)
for (unsigned int dir=0; dir<dim; ++dir)
{
const unsigned int support_point_index
= finite_element.component_to_system_index(/*velocity component=*/ introspection.component_indices.velocities[dir],
/*dof index within component=*/ j);
vec[cell_dof_indices[support_point_index]] = func.value(mesh_points.quadrature_point(j))[dir];
}
}
vec.compress(VectorOperation::insert);
#if DEAL_II_VERSION_GTE(9,7,0)
AffineConstraints<double> hanging_node_constraints(introspection.index_sets.system_relevant_set,
introspection.index_sets.system_relevant_set);
#else
AffineConstraints<double> hanging_node_constraints(introspection.index_sets.system_relevant_set);
#endif
DoFTools::make_hanging_node_constraints(dof_handler, hanging_node_constraints);
hanging_node_constraints.close();
// Create a view of all constraints that only pertains to the
// Stokes subset of degrees of freedom. We can then use this later
// to call constraints.distribute(), constraints.set_zero(), etc.,
// on those block vectors that only have the Stokes components in
// them.
//
// For the moment, assume that the Stokes degrees are first in the
// overall vector, so that they form a contiguous range starting
// at zero. The assertion checks this, but this could easily be
// generalized if the Stokes block were not starting at zero.
#if DEAL_II_VERSION_GTE(9,6,0)
Assert (introspection.block_indices.velocities == 0,
ExcNotImplemented());
if (parameters.use_direct_stokes_solver == false)
Assert (introspection.block_indices.pressure == 1,
ExcNotImplemented());
IndexSet stokes_dofs (dof_handler.n_dofs());
stokes_dofs.add_range (0, vec.size());
const AffineConstraints<double> stokes_hanging_node_constraints
= hanging_node_constraints.get_view (stokes_dofs);
#else
const AffineConstraints<double> &stokes_hanging_node_constraints = hanging_node_constraints;
#endif
stokes_hanging_node_constraints.distribute(vec);
}
template <int dim>
double Simulator<dim>::normalize_pressure (LinearAlgebra::BlockVector &vector) const
{
if (parameters.pressure_normalization == "no")
return 0;
const FEValuesExtractors::Scalar &extractor_pressure =
(parameters.include_melt_transport ?
introspection.variable("fluid pressure").extractor_scalar()
: introspection.extractors.pressure);
double my_pressure = 0.0;
double my_area = 0.0;
if (parameters.pressure_normalization == "surface")
{
const types::boundary_id top_boundary_id = geometry_model->translate_symbolic_boundary_name_to_id("top");
const Quadrature<dim-1> &quadrature = this->introspection.face_quadratures.pressure;
const unsigned int n_q_points = quadrature.size();
FEFaceValues<dim> fe_face_values (*mapping, finite_element, quadrature,
update_JxW_values | update_values);
std::vector<double> pressure_values(n_q_points);
for (const auto &cell : dof_handler.active_cell_iterators())
if (cell->is_locally_owned())
{
for (const unsigned int face_no : cell->face_indices())
{
const typename DoFHandler<dim>::face_iterator face = cell->face (face_no);
if (face->at_boundary() && face->boundary_id() == top_boundary_id)
{
fe_face_values.reinit (cell, face_no);
fe_face_values[extractor_pressure].get_function_values(vector,
pressure_values);
for (unsigned int q = 0; q < n_q_points; ++q)
{
const double JxW = fe_face_values.JxW(q);
my_pressure += pressure_values[q]
* JxW;
my_area += JxW;
}
}
}
}
}
else if (parameters.pressure_normalization == "volume")
{
const Quadrature<dim> &quadrature = this->introspection.quadratures.pressure;
const unsigned int n_q_points = quadrature.size();
FEValues<dim> fe_values (*mapping, finite_element, quadrature,
update_JxW_values | update_values);
std::vector<double> pressure_values(n_q_points);
for (const auto &cell : dof_handler.active_cell_iterators())
if (cell->is_locally_owned())
{
fe_values.reinit (cell);
fe_values[extractor_pressure].get_function_values(vector,
pressure_values);
for (unsigned int q = 0; q < n_q_points; ++q)
{
my_pressure += pressure_values[q]
* fe_values.JxW (q);
my_area += fe_values.JxW (q);
}
}
}
else
AssertThrow (false, ExcMessage("Invalid pressure normalization method: " +
parameters.pressure_normalization));
// sum up the integrals from each processor and compute the result we care about
double pressure_adjustment = numbers::signaling_nan<double>();
{
const double my_temp[2] = {my_pressure, my_area};
double temp[2];
Utilities::MPI::sum (my_temp, mpi_communicator, temp);
const double pressure = temp[0];
const double area = temp[1];
Assert (area > 0,
ExcMessage("While computing the average pressure, the area/volume "
"to integrate over was found to be zero or negative. This "
"indicates that no appropriate surface faces were found, "
"which is typically the case if the geometry model is not "
"set up correctly."));
if (parameters.pressure_normalization == "surface")
pressure_adjustment = -pressure/area + parameters.surface_pressure;
else if (parameters.pressure_normalization == "volume")
pressure_adjustment = -pressure/area;
else
AssertThrow(false, ExcNotImplemented());
}
// A complication is that we can't modify individual
// elements of the solution vector since that one has ghost element.
// rather, we first need to localize it and then distribute back
LinearAlgebra::BlockVector distributed_vector (introspection.index_sets.system_partitioning,
mpi_communicator);
distributed_vector = vector;
if (parameters.use_locally_conservative_discretization == false)
{
if (introspection.block_indices.velocities != introspection.block_indices.pressure
&& !parameters.include_melt_transport)
distributed_vector.block(introspection.block_indices.pressure).add(pressure_adjustment);
else
{
// pressure is not in a separate block, so we have to modify the values manually
const unsigned int pressure_component = (parameters.include_melt_transport ?
introspection.variable("fluid pressure").first_component_index
: introspection.component_indices.pressure);
const unsigned int n_local_pressure_dofs = (parameters.include_melt_transport ?
finite_element.base_element(introspection.variable("fluid pressure").base_index).dofs_per_cell
: finite_element.base_element(introspection.base_elements.pressure).dofs_per_cell);
std::vector<types::global_dof_index> local_dof_indices (finite_element.dofs_per_cell);
for (const auto &cell : dof_handler.active_cell_iterators())
if (cell->is_locally_owned())
{
cell->get_dof_indices (local_dof_indices);
for (unsigned int j=0; j<n_local_pressure_dofs; ++j)
{
unsigned int support_point_index
= finite_element.component_to_system_index(pressure_component,
/*dof index within component=*/ j);
// then adjust its value. Note that because we end up touching
// entries more than once, we are not simply incrementing
// distributed_vector but copy from the unchanged vector.
distributed_vector(local_dof_indices[support_point_index]) = vector(local_dof_indices[support_point_index]) + pressure_adjustment;
}
}
distributed_vector.compress(VectorOperation::insert);
}
}
else
{
// this case is a bit more complicated: if the condition above is false
// then we use the FE_DGP element for which the shape functions do not
// add up to one; consequently, adding a constant to all degrees of
// freedom does not alter the overall function by that constant, but
// by something different
//
// we can work around this by using the documented property of the
// FE_DGP element that the first shape function is constant.
// consequently, adding the adjustment to the global function is
// achieved by adding the adjustment to the first pressure degree
// of freedom on each cell.
Assert (dynamic_cast<const FE_DGP<dim>*>(&finite_element.base_element(introspection.base_elements.pressure)) != nullptr,
ExcInternalError());
const unsigned int pressure_component = (parameters.include_melt_transport ?
introspection.variable("fluid pressure").first_component_index
: introspection.component_indices.pressure);
std::vector<types::global_dof_index> local_dof_indices (finite_element.dofs_per_cell);
for (const auto &cell : dof_handler.active_cell_iterators())
if (cell->is_locally_owned())
{
// identify the first pressure dof
cell->get_dof_indices (local_dof_indices);
const unsigned int first_pressure_dof
= finite_element.component_to_system_index (pressure_component, 0);
// make sure that this DoF is really owned by the current processor
// and that it is in fact a pressure dof
Assert (dof_handler.locally_owned_dofs().is_element(local_dof_indices[first_pressure_dof]),
ExcInternalError());
// then adjust its value
distributed_vector(local_dof_indices[first_pressure_dof]) = vector(local_dof_indices[first_pressure_dof])
+ pressure_adjustment;
}
distributed_vector.compress(VectorOperation::insert);
}
// now get back to the original vector and return the adjustment used
// in the computations above
vector = distributed_vector;
return pressure_adjustment;
}
template <int dim>
void
Simulator<dim>::
denormalize_pressure (const double pressure_adjustment,
LinearAlgebra::BlockVector &vector) const
{
if (parameters.pressure_normalization == "no")
return;
if (parameters.use_locally_conservative_discretization == false)
{
if ((introspection.block_indices.velocities != introspection.block_indices.pressure)
&& !parameters.include_melt_transport)
vector.block(introspection.block_indices.pressure).add(-1.0 * pressure_adjustment);
else
{
// pressure is not in a separate block so we have to modify the values manually
const unsigned int pressure_component = (parameters.include_melt_transport ?
introspection.variable("fluid pressure").first_component_index
: introspection.component_indices.pressure);
const unsigned int n_local_pressure_dofs = (parameters.include_melt_transport ?
finite_element.base_element(introspection.variable("fluid pressure").base_index).dofs_per_cell
: finite_element.base_element(introspection.base_elements.pressure).dofs_per_cell);
// We may touch the same DoF multiple times, so we need to copy the
// vector before modifying it to have access to the original value.
LinearAlgebra::BlockVector vector_backup;
vector_backup.reinit(vector, /* omit_zeroing_entries = */ true);
const unsigned int pressure_block_index =
parameters.include_melt_transport ?
introspection.variable("fluid pressure").block_index
:
introspection.block_indices.pressure;
vector_backup.block(pressure_block_index) = vector.block(pressure_block_index);
std::vector<types::global_dof_index> local_dof_indices (finite_element.dofs_per_cell);
for (const auto &cell : dof_handler.active_cell_iterators())
if (cell->is_locally_owned())
{
cell->get_dof_indices (local_dof_indices);
for (unsigned int j=0; j<n_local_pressure_dofs; ++j)
{
const unsigned int local_dof_index
= finite_element.component_to_system_index(pressure_component,
/*dof index within component=*/ j);
// Then adjust its value. vector could be a vector with ghost elements
// or a fully distributed vector. In the latter case only access dofs that are
// locally owned. Note that because we end up touching
// entries more than once, we are not simply incrementing
// vector but copy from the vector_backup copy.
if (vector.has_ghost_elements() ||
dof_handler.locally_owned_dofs().is_element(local_dof_indices[local_dof_index]))
vector(local_dof_indices[local_dof_index])
= vector_backup(local_dof_indices[local_dof_index]) - pressure_adjustment;
}
}
vector.compress(VectorOperation::insert);
}
}
else
{
// this case is a bit more complicated: if the condition above is false
// then we use the FE_DGP element for which the shape functions do not
// add up to one; consequently, adding a constant to all degrees of
// freedom does not alter the overall function by that constant, but
// by something different
//
// we can work around this by using the documented property of the
// FE_DGP element that the first shape function is constant.