-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathdynamic-muscle.cc
4932 lines (4239 loc) · 195 KB
/
dynamic-muscle.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
/* ---------------------------------------------------------------------
*
* Flexodeal Lite
* Copyright (C) 2024 Neuromuscular Mechanics Laboratory
*
* This file is part of the Flexodeal library
*
* The Flexodeal library is free software; you can use it, redistribute
* it, and/or modify it under the terms of the GNU Lesser General
* Public License as published by the Free Software Foundation; either
* version 2.1 of the License, or (at your option) any later version.
* The full text of the license can be found in the file LICENSE.md at
* the top level directory of Flexodeal.
*
* ---------------------------------------------------------------------
*
* Author: Javier Almonacid
* PhD Candidate, Applied and Computational Mathematics
* Neuromuscular Mechanics Laboratory (NML)
* Simon Fraser University
* Spring 2024
*
* This software has been created based on the "muscle code" developed by
* members of the NML since 2014:
*
* Harshil Pathak
* Ryan N. Konno
* Cassidy Tam
* Sebastian A. Dominguez-Rivera
* Stephanie A. Ross
* David Ryan
* Hadi Rahemi
* Prof. James M. Wakeling (SFU Biomedical Physiology and Kinesiology)
* Prof. Nilima Nigam (SFU Mathematics)
*
* Furthermore, the structure of the code is based on deal.II's step-44 tutorial:
*
* Pelteret, J.-P., & McBride, A. (2012). The deal.II tutorial step-44:
* Three-field formulation for non-linear solid mechanics.
* Zenodo. https://doi.org/10.5281/zenodo.439772
*
* Some comments have been preserved from step-44 to increase readability of the
* code. The mathematical details of the muscle model in use are available here:
*
* Almonacid, J. A., Domínguez-Rivera, S. A., Konno, R. N., Nigam, N.,
* Ross, S. A., Tam, C., & Wakeling, J. M. (2024).
* A three-dimensional model of skeletal muscle tissues.
* SIAM Journal on Applied Mathematics, S538-S566.
*
*
*/
// First, we include all the headers necessary for this code. All the headers
// (with the exception of fe_dgq.h) have already been discussed in step-44.
#include <deal.II/base/function.h>
#include <deal.II/base/parameter_handler.h>
#include <deal.II/base/point.h>
#include <deal.II/base/quadrature_lib.h>
#include <deal.II/base/symmetric_tensor.h>
#include <deal.II/base/tensor.h>
#include <deal.II/base/timer.h>
#include <deal.II/base/work_stream.h>
#include <deal.II/dofs/dof_renumbering.h>
#include <deal.II/dofs/dof_tools.h>
#include <deal.II/base/quadrature_point_data.h>
#include <deal.II/grid/grid_generator.h>
#include <deal.II/grid/grid_tools.h>
#include <deal.II/grid/grid_in.h>
#include <deal.II/grid/grid_out.h>
#include <deal.II/grid/tria.h>
#include <deal.II/fe/fe_dgp_monomial.h>
#include <deal.II/fe/fe_dgq.h>
#include <deal.II/fe/fe_q.h>
#include <deal.II/fe/fe_system.h>
#include <deal.II/fe/fe_tools.h>
#include <deal.II/fe/fe_values.h>
#include <deal.II/fe/mapping_q_eulerian.h>
#include <deal.II/lac/block_sparse_matrix.h>
#include <deal.II/lac/block_vector.h>
#include <deal.II/lac/dynamic_sparsity_pattern.h>
#include <deal.II/lac/full_matrix.h>
#include <deal.II/lac/precondition_selector.h>
#include <deal.II/lac/solver_cg.h>
#include <deal.II/lac/solver_gmres.h>
#include <deal.II/lac/solver_selector.h>
#include <deal.II/lac/sparse_direct.h>
#include <deal.II/lac/affine_constraints.h>
#include <deal.II/lac/linear_operator.h>
#include <deal.II/lac/packaged_operation.h>
#include <deal.II/numerics/data_out.h>
#include <deal.II/numerics/vector_tools.h>
#include <deal.II/physics/elasticity/kinematics.h>
#include <deal.II/physics/elasticity/standard_tensors.h>
#include <sys/stat.h>
#include <sys/types.h>
#include <fstream>
// Then, we place all functions and classes inside a namespace of its own.
namespace Flexodeal
{
using namespace dealii;
// @sect3{Run-time parameters}
//
// There are several parameters that can be set in the code so we set up a
// ParameterHandler object to read in the choices at run-time.
namespace Parameters
{
// @sect4{Finite Element system}
// As mentioned in the introduction, a different order interpolation should
// be used for the displacement $\mathbf{u}$ than for the pressure
// $\widetilde{p}$ and the dilatation $\widetilde{J}$. Choosing
// $\widetilde{p}$ and $\widetilde{J}$ as discontinuous (constant) functions
// at the element level leads to the mean-dilatation method. The
// discontinuous approximation allows $\widetilde{p}$ and $\widetilde{J}$ to
// be condensed out and a classical displacement based method is recovered.
// Here we specify the polynomial order used to approximate the solution.
// The quadrature order should be adjusted accordingly.
struct FESystem
{
unsigned int poly_degree;
unsigned int quad_order;
std::string type_of_simulation;
static void declare_parameters(ParameterHandler &prm);
void parse_parameters(ParameterHandler &prm);
};
void FESystem::declare_parameters(ParameterHandler &prm)
{
prm.enter_subsection("Finite element system");
{
prm.declare_entry("Polynomial degree",
"2",
Patterns::Integer(0),
"Displacement system polynomial order");
prm.declare_entry("Quadrature order",
"3",
Patterns::Integer(0),
"Gauss quadrature order");
prm.declare_entry("Type of simulation",
"quasi-static",
Patterns::Selection("quasi-static|dynamic"),
"Type of simulation (quasi-static or dynamic)");
}
prm.leave_subsection();
}
void FESystem::parse_parameters(ParameterHandler &prm)
{
prm.enter_subsection("Finite element system");
{
poly_degree = prm.get_integer("Polynomial degree");
quad_order = prm.get_integer("Quadrature order");
type_of_simulation = prm.get("Type of simulation");
}
prm.leave_subsection();
}
// @sect4{Geometry}
// Next, we process the parameters related to the geometry of the problem.
// To keep things simple, we consider a block of muscle tissue of given
// length, width, and height. Note that this structure will probably have
// to be adjusted whenever a different geometry is considered.
struct Geometry
{
unsigned int global_refinement;
double length;
double width;
double height;
double scale;
double p_p0;
static void declare_parameters(ParameterHandler &prm);
void parse_parameters(ParameterHandler &prm);
};
void Geometry::declare_parameters(ParameterHandler &prm)
{
prm.enter_subsection("Geometry");
{
prm.declare_entry("Global refinement",
"2",
Patterns::Integer(0),
"Global refinement level");
prm.declare_entry("Length", "1.0",
Patterns::Double(0.0),
"Length in the x direction");
prm.declare_entry("Width", "1.0",
Patterns::Double(0.0),
"Width in the y direction");
prm.declare_entry("Height", "1.0",
Patterns::Double(0.0),
"Height in the z direction");
prm.declare_entry("Grid scale",
"1e-3",
Patterns::Double(0.0),
"Global grid scaling factor");
prm.declare_entry("Pressure ratio p/p0",
"100",
Patterns::Selection("0|20|40|60|80|100"),
"Ratio of applied pressure to reference pressure");
}
prm.leave_subsection();
}
void Geometry::parse_parameters(ParameterHandler &prm)
{
prm.enter_subsection("Geometry");
{
global_refinement = prm.get_integer("Global refinement");
length = prm.get_double("Length");
width = prm.get_double("Width");
height = prm.get_double("Height");
scale = prm.get_double("Grid scale");
p_p0 = prm.get_double("Pressure ratio p/p0");
}
prm.leave_subsection();
}
// @sect4{Muscle properties}
// Then, we process all the intrinsic properties of muscle tissue.
// In this model, muscle is viewed as a bundle of muscle fibres
// surrounded by a base material.
struct MuscleProperties
{
double muscle_density;
// Fibre properties
double max_iso_stress_muscle;
double kappa_muscle;
double max_strain_rate;
double muscle_fibre_orientation_x;
double muscle_fibre_orientation_y;
double muscle_fibre_orientation_z;
// Base material properties
double max_iso_stress_basematerial;
double muscle_basematerial_factor;
double muscle_basemat_c1;
double muscle_basemat_c2;
double muscle_basemat_c3;
static void
declare_parameters(ParameterHandler &prm);
void
parse_parameters(ParameterHandler &prm);
};
void MuscleProperties::declare_parameters(ParameterHandler &prm)
{
prm.enter_subsection("Materials");
{
prm.declare_entry("Muscle density", "1060",
Patterns::Double(),
"Muscle tissue density");
// Fibre properties
prm.declare_entry("Sigma naught muscle", "2.0e5",
Patterns::Double(),
"Muscle maximum isometric stress");
prm.declare_entry("Bulk modulus muscle", "1.0e6",
Patterns::Double(),
"Muscle kappa value");
prm.declare_entry("Max strain rate", "0.0",
Patterns::Double(),
"Maximum muscle fibre strain rate");
prm.declare_entry("Muscle x component", "1.0",
Patterns::Double(0.0),
"Muscle fibre orientation x direction");
prm.declare_entry("Muscle y component", "0.0",
Patterns::Double(0.0),
"Muscle fibre orientation y direction");
prm.declare_entry("Muscle z component", "0.0",
Patterns::Double(0.0),
"Muscle fibre orientation z direction");
// Base material properties
prm.declare_entry("Sigma naught base material", "2.0e5",
Patterns::Double(),
"Base material maximum isometric stress");
prm.declare_entry("Muscle base material factor", "1.0e1",
Patterns::Double(),
"Fictitious muscle base material multiplier");
prm.declare_entry("Muscle base material constant 1", "0.0",
Patterns::Double(),
"Muscle base material constant 1");
prm.declare_entry("Muscle base material constant 2", "0.0",
Patterns::Double(),
"Muscle base material constant 2");
prm.declare_entry("Muscle base material constant 3", "0.0",
Patterns::Double(),
"Muscle base material constant 3");
}
prm.leave_subsection();
}
void MuscleProperties::parse_parameters(ParameterHandler &prm)
{
prm.enter_subsection("Materials");
{
muscle_density = prm.get_double("Muscle density");
max_iso_stress_muscle = prm.get_double("Sigma naught muscle");
kappa_muscle = prm.get_double("Bulk modulus muscle");
max_strain_rate = prm.get_double("Max strain rate");
muscle_fibre_orientation_x = prm.get_double("Muscle x component");
muscle_fibre_orientation_y = prm.get_double("Muscle y component");
muscle_fibre_orientation_z = prm.get_double("Muscle z component");
max_iso_stress_basematerial = prm.get_double("Sigma naught base material");
muscle_basematerial_factor = prm.get_double("Muscle base material factor");
muscle_basemat_c1 = prm.get_double("Muscle base material constant 1");
muscle_basemat_c2 = prm.get_double("Muscle base material constant 2");
muscle_basemat_c3 = prm.get_double("Muscle base material constant 3");
}
prm.leave_subsection();
}
// @sect4{Linear solver}
// Next, we choose both solver and preconditioner settings. The use of an
// effective preconditioner is critical to ensure convergence when a large
// nonlinear motion occurs within a Newton increment.
struct LinearSolver
{
std::string type_lin;
double tol_lin;
double max_iterations_lin;
bool use_static_condensation;
std::string preconditioner_type;
double preconditioner_relaxation;
static void declare_parameters(ParameterHandler &prm);
void parse_parameters(ParameterHandler &prm);
};
void LinearSolver::declare_parameters(ParameterHandler &prm)
{
prm.enter_subsection("Linear solver");
{
prm.declare_entry("Solver type",
"CG",
Patterns::Selection("GMRES|CG|Direct"),
"Type of solver used to solve the linear system");
prm.declare_entry("Residual",
"1e-6",
Patterns::Double(0.0),
"Linear solver residual (scaled by residual norm)");
prm.declare_entry(
"Max iteration multiplier",
"1",
Patterns::Double(0.0),
"Linear solver iterations (multiples of the system matrix size)");
prm.declare_entry("Use static condensation",
"true",
Patterns::Bool(),
"Solve the full block system or a reduced problem");
prm.declare_entry("Preconditioner type",
"ssor",
Patterns::Selection("jacobi|ssor"),
"Type of preconditioner");
prm.declare_entry("Preconditioner relaxation",
"0.65",
Patterns::Double(0.0),
"Preconditioner relaxation value");
}
prm.leave_subsection();
}
void LinearSolver::parse_parameters(ParameterHandler &prm)
{
prm.enter_subsection("Linear solver");
{
type_lin = prm.get("Solver type");
tol_lin = prm.get_double("Residual");
max_iterations_lin = prm.get_double("Max iteration multiplier");
use_static_condensation = prm.get_bool("Use static condensation");
preconditioner_type = prm.get("Preconditioner type");
preconditioner_relaxation = prm.get_double("Preconditioner relaxation");
}
prm.leave_subsection();
}
// @sect4{Nonlinear solver}
// A Newton scheme is used to solve the nonlinear system of governing
// equations. Below we process the stopping criteria for the solver.
struct NonlinearSolver
{
std::string type_nonlinear_solver;
unsigned int max_iterations_NR;
double tol_f;
double tol_u;
static void declare_parameters(ParameterHandler &prm);
void parse_parameters(ParameterHandler &prm);
};
void NonlinearSolver::declare_parameters(ParameterHandler &prm)
{
prm.enter_subsection("Nonlinear solver");
{
prm.declare_entry("Nonlinear solver type", "classicNewton",
Patterns::Selection("classicNewton|acceleratedNewton"),
"Type of nonlinear iteration");
prm.declare_entry("Max iterations Newton-Raphson",
"10",
Patterns::Integer(0),
"Number of Newton-Raphson iterations allowed");
prm.declare_entry("Tolerance force",
"1.0e-9",
Patterns::Double(0.0),
"Force residual tolerance");
prm.declare_entry("Tolerance displacement",
"1.0e-6",
Patterns::Double(0.0),
"Displacement error tolerance");
}
prm.leave_subsection();
}
void NonlinearSolver::parse_parameters(ParameterHandler &prm)
{
prm.enter_subsection("Nonlinear solver");
{
type_nonlinear_solver = prm.get("Nonlinear solver type");
max_iterations_NR = prm.get_integer("Max iterations Newton-Raphson");
tol_f = prm.get_double("Tolerance force");
tol_u = prm.get_double("Tolerance displacement");
}
prm.leave_subsection();
}
// @sect4{Time}
// Set the timestep size $ \varDelta t $ and the simulation end-time.
struct Time
{
double delta_t;
double end_time;
static void declare_parameters(ParameterHandler &prm);
void parse_parameters(ParameterHandler &prm);
};
void Time::declare_parameters(ParameterHandler &prm)
{
prm.enter_subsection("Time");
{
prm.declare_entry("End time", "1", Patterns::Double(), "End time");
prm.declare_entry("Time step size",
"0.1",
Patterns::Double(),
"Time step size");
}
prm.leave_subsection();
}
void Time::parse_parameters(ParameterHandler &prm)
{
prm.enter_subsection("Time");
{
end_time = prm.get_double("End time");
delta_t = prm.get_double("Time step size");
}
prm.leave_subsection();
}
// @sect4{PrescribedDisplacement}
// Set the parameters for the prescribed displacement. Because the profile
// itself is given from control_points_strain.dat file, the only thing we
// keep track here is the ID of the face we are pulling/pushing from. This
// is critical to correctly identify the face for which forces will be
// reported as a time series.
struct PrescribedDisplacement
{
unsigned int pulling_face_id;
static void
declare_parameters(ParameterHandler &prm);
void
parse_parameters(ParameterHandler &prm);
};
void PrescribedDisplacement::declare_parameters(ParameterHandler &prm)
{
prm.enter_subsection("Prescribed displacement");
{
prm.declare_entry("Pulling face ID", "1",
Patterns::Integer(0),
"Boundary ID of face being pulled/pushed");
}
prm.leave_subsection();
}
void PrescribedDisplacement::parse_parameters(ParameterHandler &prm)
{
prm.enter_subsection("Prescribed displacement");
{
pulling_face_id = prm.get_integer("Pulling face ID");
}
prm.leave_subsection();
}
// @sect4{Measuring locations}
// We use this section to indicate the file which contains a list
// of markers. We call a "marker" a mesh vertex that contains
// displacement degrees of freedom, so that we can output displacement
// information at particular points. Older versions of this code had only
// support for three points (left, mid, and right), but now all we need
// is a list containing a label and the 3D point. For example, a file
// called markers.dat could contain the following information:
// p1 0.0075 0.005 0.005
// p2 0.015 0.005 0.005
// p3 0.0225 0.005 0.005
// Note that the markers MUST corresponds to places with availability of
// displacement degrees of freedom (for Q2 and Q1 elements, these are mesh
// vertices).
struct MeasuringLocations
{
std::string markers_list_file;
static void
declare_parameters(ParameterHandler &prm);
void
parse_parameters(ParameterHandler &prm);
};
void Parameters::MeasuringLocations::declare_parameters(ParameterHandler &prm)
{
prm.enter_subsection("Measuring locations");
{
prm.declare_entry("Markers list file", "",
Patterns::FileName(),
"File including list of markers");
}
prm.leave_subsection();
}
void Parameters::MeasuringLocations::parse_parameters(ParameterHandler &prm)
{
prm.enter_subsection("Measuring locations");
{
markers_list_file = prm.get("Markers list file");
}
prm.leave_subsection();
}
// @sect4{All parameters}
// Finally we consolidate all of the above structures into a single
// container that holds all of our run-time selections.
struct AllParameters : public FESystem,
public Geometry,
public MuscleProperties,
public LinearSolver,
public NonlinearSolver,
public Time,
public PrescribedDisplacement,
public MeasuringLocations
{
AllParameters(const std::string &input_file);
static void declare_parameters(ParameterHandler &prm);
void parse_parameters(ParameterHandler &prm);
// We decide to keep track of the parameter handler
// object since we want to call prm.print_parameters()
// after the constructor has been called and we need
// to first set up the directories with the current
// timestamp.
ParameterHandler prm;
};
AllParameters::AllParameters(const std::string &input_file)
{
declare_parameters(prm);
prm.parse_input(input_file);
parse_parameters(prm);
}
void AllParameters::declare_parameters(ParameterHandler &prm)
{
FESystem::declare_parameters(prm);
Geometry::declare_parameters(prm);
MuscleProperties::declare_parameters(prm);
LinearSolver::declare_parameters(prm);
NonlinearSolver::declare_parameters(prm);
Time::declare_parameters(prm);
PrescribedDisplacement::declare_parameters(prm);
MeasuringLocations::declare_parameters(prm);
}
void AllParameters::parse_parameters(ParameterHandler &prm)
{
FESystem::parse_parameters(prm);
Geometry::parse_parameters(prm);
MuscleProperties::parse_parameters(prm);
LinearSolver::parse_parameters(prm);
NonlinearSolver::parse_parameters(prm);
Time::parse_parameters(prm);
PrescribedDisplacement::parse_parameters(prm);
MeasuringLocations::parse_parameters(prm);
}
} // namespace Parameters
// @sect3{Time class}
// A simple class to store time data. Its functioning is transparent so no
// discussion is necessary. For simplicity we assume a constant time step
// size.
class Time
{
public:
Time(const double time_end, const double delta_t)
: timestep(0)
, time_previous(0.0)
, time_current(0.0)
, time_end(time_end)
, delta_t(delta_t)
{}
virtual ~Time() = default;
double current() const
{
return time_current;
}
double previous() const
{
return time_previous;
}
double end() const
{
return time_end;
}
double get_delta_t() const
{
return delta_t;
}
unsigned int get_timestep() const
{
return timestep;
}
void increment()
{
time_previous = time_current;
time_current += delta_t;
++timestep;
}
private:
unsigned int timestep;
double time_previous;
double time_current;
const double time_end;
const double delta_t;
};
// @sect3{Function from tabular data and interpolation tools}
// This class can be used to read data from a table of
// x,y values (separated by a space or tab) and construct
// a continuous function from them using linear interpolation.
// This will be particularly useful to read activations and
// muscle length profiles from data. Although it may seem like
// an overkill, this class can also be used to implement simple
// linear ramps. For instance, a linear ramp between (0,0) and
// (0.5, 1.0) only requires a .dat file with 2 lines:
//
// 0.0 0.0
// 0.5 1.0
//
class TabularFunction
{
public:
TabularFunction(const std::string filename)
{
initialize_map(filename);
}
double operator()(const double t) const;
private:
std::map<double,double> table_values;
void initialize_map(const std::string filename);
};
// Read data (control points) and store them in the table_values map
void TabularFunction::initialize_map(const std::string filename)
{
std::ifstream infile(filename);
// Raise an exception if file cannot be open (perhaps it does not
// even exist!).
if (infile.fail())
throw std::invalid_argument("Cannot open file: " + filename +
". Make sure the file exists and it has read permissions.");
double x, y;
while (infile >> x >> y)
table_values.insert({x, y});
}
// Evaluate the data:
// - Interpolate if not found and t (independent variable)
// is in data range
// - Retrieve the original y coordinate if found
// - Output a constant value if t exceeds data range
double TabularFunction::operator()(const double t) const
{
double out = 1000000; /*Bogus value*/
auto iter_t = table_values.find(t);
if (iter_t == table_values.end()) /* Value not found: interpolate! */
{
if (t <= table_values.rbegin()->first)
{
auto t1 = table_values.upper_bound(t);
auto t0 = std::prev(t1);
out = ((t1->second - t0->second)/(t1->first - t0->first)) * (t - t0->first) + t0->second;
}
else
out = table_values.rbegin()->second; /* Constant after last point */
}
else /* Value found! Return this value */
out = iter_t->second;
return out;
}
// @sect3{Marker class (measuring locations)}
// This is an extension of the Point<dim> class. Additional attributes are
// a label, a boolean used to track whether this point was found in the
// as a vertex in the mesh or not, and a global_dof_index which will be
// set during Solid<dim>::system_setup.
template <int dim>
struct Marker : public Point<dim>
{
Marker(const std::string label,
const Tensor<1, dim> location)
: Point<dim>(location)
, label(label)
, found(false)
{}
const std::string label;
bool found;
std::vector<dealii::types::global_dof_index> global_dof_index;
};
// @sect3{Muscle tissue within a three-field formulation}
// Muscle can be described as a quasi-incompressible fibre-reinforced
// material. Similar to the <code>Material_Compressible_Neo_Hook_Three_Field</code>
// in step-44, the material here can be described by a strain-energy
// function $ \Psi = \Psi_{\text{iso}}(\overline{\mathbf{b}}) +
// \Psi_{\text{vol}}(\widetilde{J})$, where $\Psi_{iso}$ is composed of a
// fibre and a base material component $\Psi_{fibre}$ and $\Psi_{base}$,
// respectively.
//
// In the current configuration, the fibre component takes the form
// $\boldsymbol{\sigma} = \sigma_{Hill} \mathbf{a} \otimes \mathbf{a}$,
// where $\sigma_{Hill} = \sigma_0 \left\{ a(t) \sigma_L(\lambda)
// \sigma_V(\epsilon) + \sigma_P(\lambda) \right\}$.
//
// In turn, the base material component is simply given by a Yeoh SEF.
//
// A note on quasi-incompressibility of muscle:
//
// Technically, muscle is compressible (Baskin & Paolini 1966). However, the
// magnitude of this volume change is so small that most physiologists would
// consider muscle as "incompressible". Therefore, the volume changes that we
// expect here will (and have to) be small and the dilation J will hover around
// the value of 1. From a numerical perspective, though these changes are small,
// they are extremely important to prevent locking effects.
template <int dim>
class Muscle_Tissues_Three_Field
{
public:
Muscle_Tissues_Three_Field(const std::string type_contraction,
const double max_iso_stress_muscle,
const double kappa_muscle,
const double max_strain_rate,
const double initial_fibre_orientation_x,
const double initial_fibre_orientation_y,
const double initial_fibre_orientation_z,
const double max_iso_stress_basematerial,
const double muscle_basematerial_factor,
const double muscle_basemat_c1,
const double muscle_basemat_c2,
const double muscle_basemat_c3)
:
type_of_contraction(type_contraction),
sigma_naught_muscle(max_iso_stress_muscle),
kappa_muscle(kappa_muscle),
strain_rate_naught(max_strain_rate),
initial_fibre_orientation({initial_fibre_orientation_x,
initial_fibre_orientation_y,
initial_fibre_orientation_z}),
sigma_naught_basematerial(max_iso_stress_basematerial),
s_base_muscle(muscle_basematerial_factor),
c1_basematerial_muscle(muscle_basemat_c1),
c2_basematerial_muscle(muscle_basemat_c2),
c3_basematerial_muscle(muscle_basemat_c3),
/* Physiological variables */
stretch_bar(1.0),
strain_rate_bar(0.0),
fibre_time_activation(0.0),
orientation(initial_fibre_orientation),
/* Mechanical variables */
det_F(1.0),
p_tilde(0.0),
J_tilde(1.0),
b_bar(Physics::Elasticity::StandardTensors<dim>::I),
trace_b_bar(3.0),
trace_d(0.0),
delta_t(0.0)
{
Assert(kappa_muscle > 0,
ExcMessage("Bulk modulus must be positive!"));
Assert(initial_fibre_orientation.norm() != 0,
ExcMessage("Initial fibre orientation must be a nonzero vector!"));
}
// We update the material model with various deformation dependent data
// based on $F$ and the pressure $\widetilde{p}$ and dilatation
// $\widetilde{J}$, and at the end of the function include a physical
// check for internal consistency:
void update_material_data(const Tensor<2, dim> &F,
const double p_tilde_in,
const double J_tilde_in,
const double fibre_time_activation_in,
const Tensor<2, dim> &grad_velocity,
const double delta_t_in)
{
// First compute the determinant of the deformation tensor
// and stop the program immediately if it is negative.
det_F = determinant(F);
AssertThrow(det_F > 0, ExcInternalError());
// Then, update the rest of the variables.
p_tilde = p_tilde_in;
J_tilde = J_tilde_in;
fibre_time_activation = fibre_time_activation_in;
delta_t = delta_t_in;
const Tensor<2, dim> F_bar = Physics::Elasticity::Kinematics::F_iso(F);
b_bar = Physics::Elasticity::Kinematics::b(F_bar);
trace_b_bar = first_invariant(b_bar);
// Update stretch_bar and strain_rate_bar using the current Newton iterate.
const SymmetricTensor<2,dim> symm_grad_velocity = Physics::Elasticity::Kinematics::d(F,grad_velocity);
const Tensor<2,dim> dev_symm_grad_velocity = Physics::Elasticity::StandardTensors<dim>::dev_P * symm_grad_velocity;
trace_d = first_invariant(symm_grad_velocity);
orientation = F_bar * initial_fibre_orientation;
stretch_bar = std::sqrt(orientation * orientation);
strain_rate_bar = (1.0 / strain_rate_naught) *
orientation * (dev_symm_grad_velocity * orientation) / stretch_bar;
}
// The second function determines the Kirchhoff stress $\boldsymbol{\tau}
// = \boldsymbol{\tau}_{\textrm{iso}} + \boldsymbol{\tau}_{\textrm{vol}}$
SymmetricTensor<2, dim> get_tau()
{
return get_tau_iso() + get_tau_vol();
}
// The following set of functions determine each contribution to the
// Kirchhoff stress in detail. We make this functions public as they
// need to be accessed by PointHistory when computing forces and
// energies. The first one determines the volumetric Kirchhoff stress
// $\boldsymbol{\tau}_{\textrm{vol}}$:
SymmetricTensor<2, dim> get_tau_vol() const
{
return p_tilde * det_F * Physics::Elasticity::StandardTensors<dim>::I;
}
// Next, determine the isochoric Kirchhoff stress
// $\boldsymbol{\tau}_{\textrm{iso}} =
// \mathcal{P}:\overline{\boldsymbol{\tau}}$:
SymmetricTensor<2, dim> get_tau_iso() const
{
return Physics::Elasticity::StandardTensors<dim>::dev_P * get_tau_bar();
}
// Just like the SEF of the system, the isochoric
// Kirchhoff stress is made of two parts:
// a fibre component and a base material component.
// In particular, we decide to subdivide the
// fibre component into its active and passive
// components.
SymmetricTensor<2,dim> get_tau_iso_muscle_active()
{
return Physics::Elasticity::StandardTensors<dim>::dev_P * get_tau_muscle_active_bar();
}
SymmetricTensor<2,dim> get_tau_iso_muscle_passive()
{
return Physics::Elasticity::StandardTensors<dim>::dev_P * get_tau_muscle_passive_bar();
}
SymmetricTensor<2,dim> get_tau_iso_muscle_basematerial()
{
return Physics::Elasticity::StandardTensors<dim>::dev_P * get_tau_muscle_basematerial_bar();
}
// The fourth-order elasticity tensor in the spatial setting
// $\mathfrak{c}$ is calculated from the SEF $\Psi$ as $ J
// \mathfrak{c}_{ijkl} = F_{iA} F_{jB} \mathfrak{C}_{ABCD} F_{kC} F_{lD}$
// where $ \mathfrak{C} = 4 \frac{\partial^2 \Psi(\mathbf{C})}{\partial
// \mathbf{C} \partial \mathbf{C}}$
SymmetricTensor<4, dim> get_Jc() const
{
return get_Jc_vol() + get_Jc_iso();
}
// Derivative of the volumetric free energy with respect to
// $\widetilde{J}$ return $\frac{\partial
// \Psi_{\text{vol}}(\widetilde{J})}{\partial \widetilde{J}}$
double get_dPsi_vol_dJ() const
{
return (kappa_muscle / 2.0) * (J_tilde - 1.0 / J_tilde);
}
// Second derivative of the volumetric free energy wrt $\widetilde{J}$. We
// need the following computation explicitly in the tangent so we make it
// public. We calculate $\frac{\partial^2
// \Psi_{\textrm{vol}}(\widetilde{J})}{\partial \widetilde{J} \partial
// \widetilde{J}}$
double get_d2Psi_vol_dJ2() const
{
return ((kappa_muscle / 2.0) * (1.0 + 1.0 / (J_tilde * J_tilde)));
}
// The next few functions return various data that we choose to store with
// the material:
double get_det_F() const
{
return det_F;
}
double get_p_tilde() const
{
return p_tilde;
}
double get_J_tilde() const