-
Notifications
You must be signed in to change notification settings - Fork 62
/
bioelectric_finite_elasticity_routines.f90
executable file
·3386 lines (2992 loc) · 199 KB
/
bioelectric_finite_elasticity_routines.f90
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
!> \file
!> \authors Thomas Heidlauf
!> \brief This module handles all routines pertaining to bioelectrics coupled with finite elasticity.
!>
!> \section LICENSE
!>
!> Version: MPL 1.1/GPL 2.0/LGPL 2.1
!>
!> The contents of this file are subject to the Mozilla Public License
!> Version 1.1 (the "License"); you may not use this file except in
!> compliance with the License. You may obtain a copy of the License at
!> http://www.mozilla.org/MPL/
!>
!> Software distributed under the License is distributed on an "AS IS"
!> basis, WITHOUT WARRANTY OF ANY KIND, either express or implied. See the
!> License for the specific language governing rights and limitations
!> under the License.
!>
!> The Original Code is OpenCMISS
!>
!> The Initial Developer of the Original Code is University of Auckland,
!> Auckland, New Zealand, the University of Oxford, Oxford, United
!> Kingdom and King's College, London, United Kingdom. Portions created
!> by the University of Auckland, the University of Oxford and King's
!> College, London are Copyright (C) 2007-2010 by the University of
!> Auckland, the University of Oxford and King's College, London.
!> All Rights Reserved.
!>
!> Contributor(s): Chris Bradley, Thomas Klotz
!>
!> Alternatively, the contents of this file may be used under the terms of
!> either the GNU General Public License Version 2 or later (the "GPL"), or
!> the GNU Lesser General Public License Version 2.1 or later (the "LGPL"),
!> in which case the provisions of the GPL or the LGPL are applicable instead
!> of those above. If you wish to allow use of your version of this file only
!> under the terms of either the GPL or the LGPL, and not to allow others to
!> use your version of this file under the terms of the MPL, indicate your
!> decision by deleting the provisions above and replace them with the notice
!> and other provisions required by the GPL or the LGPL. If you do not delete
!> the provisions above, a recipient may use your version of this file under
!> the terms of any one of the MPL, the GPL or the LGPL.
!>
!>This module handles all routines pertaining to bioelectrics coupled with finite elasticity.
MODULE BIOELECTRIC_FINITE_ELASTICITY_ROUTINES
USE BaseRoutines
USE BasisRoutines
USE BasisAccessRoutines
USE BIOELECTRIC_ROUTINES
USE BIODOMAIN_EQUATION_ROUTINES
USE Constants
USE CONTROL_LOOP_ROUTINES
USE ControlLoopAccessRoutines
USE EquationsRoutines
USE EquationsAccessRoutines
USE EquationsMappingRoutines
USE EquationsMappingAccessRoutines
USE EquationsMatricesRoutines
USE EquationsSetConstants
USE FIELD_IO_ROUTINES
USE FIELD_ROUTINES
USE FieldAccessRoutines
USE FINITE_ELASTICITY_ROUTINES
USE INPUT_OUTPUT
USE ISO_VARYING_STRING
USE Kinds
USE Maths
USE PROBLEM_CONSTANTS
USE Strings
USE SOLVER_ROUTINES
USE SolverAccessRoutines
USE Types
#include "macros.h"
IMPLICIT NONE
PUBLIC BioelectricFiniteElasticity_EquationsSetSetup
PUBLIC BioelectricFiniteElasticity_EquationsSetSolutionMethodSet
PUBLIC BIOELECTRIC_FINITE_ELASTICITY_PROBLEM_SETUP
PUBLIC BioelectricFiniteElasticity_ProblemSpecificationSet
PUBLIC BioelectricFiniteElasticity_FiniteElementCalculate
PUBLIC BIOELECTRIC_FINITE_ELASTICITY_PRE_SOLVE
PUBLIC BIOELECTRIC_FINITE_ELASTICITY_POST_SOLVE
PUBLIC BioelectricFiniteElasticity_ControlLoopPreLoop
PUBLIC BioelectricFiniteElasticity_ControlLoopPostLoop
PUBLIC BioelectricFiniteElasticity_UpdateGeometricField
CONTAINS
!
!================================================================================================================================
!
!>Sets/changes the solution method for a bioelectrics finite elasticity equation type of a multi physics equations set class.
SUBROUTINE BioelectricFiniteElasticity_EquationsSetSolutionMethodSet(EQUATIONS_SET,SOLUTION_METHOD,ERR,ERROR,*)
!Argument variables
TYPE(EQUATIONS_SET_TYPE), POINTER :: EQUATIONS_SET !<A pointer to the equations set to set the solution method for
INTEGER(INTG), INTENT(IN) :: SOLUTION_METHOD !<The solution method to set
INTEGER(INTG), INTENT(OUT) :: ERR !<The error code
TYPE(VARYING_STRING), INTENT(OUT) :: ERROR !<The error string
!Local Variables
TYPE(VARYING_STRING) :: LOCAL_ERROR
ENTERS("BioelectricFiniteElasticity_EquationsSetSolutionMethodSet",ERR,ERROR,*999)
IF(ASSOCIATED(EQUATIONS_SET)) THEN
IF(.NOT.ALLOCATED(EQUATIONS_SET%SPECIFICATION)) THEN
CALL FlagError("Equations set specification is not allocated.",err,error,*999)
ELSE IF(SIZE(EQUATIONS_SET%SPECIFICATION,1)/=3) THEN
CALL FlagError("Equations set specification must have three entries for a "// &
& "Bioelectric-finite elasticity type equations set.",err,error,*999)
END IF
SELECT CASE(EQUATIONS_SET%SPECIFICATION(3))
CASE(EQUATIONS_SET_STANDARD_MONODOMAIN_ELASTICITY_SUBTYPE, &
& EQUATIONS_SET_1D3D_MONODOMAIN_ELASTICITY_SUBTYPE, &
& EQUATIONS_SET_MONODOMAIN_ELASTICITY_W_TITIN_SUBTYPE, &
& EQUATIONS_SET_MONODOMAIN_ELASTICITY_VELOCITY_SUBTYPE, &
& EQUATIONS_SET_1D3D_MONODOMAIN_ACTIVE_STRAIN_SUBTYPE, &
& EQUATIONS_SET_MONODOMAIN_ELASTICITY_MUSCLE_TENDON_SUBTYPE)
SELECT CASE(SOLUTION_METHOD)
CASE(EQUATIONS_SET_FEM_SOLUTION_METHOD)
EQUATIONS_SET%SOLUTION_METHOD=EQUATIONS_SET_FEM_SOLUTION_METHOD
CASE(EQUATIONS_SET_BEM_SOLUTION_METHOD)
CALL FlagError("Not implemented.",ERR,ERROR,*999)
CASE(EQUATIONS_SET_FD_SOLUTION_METHOD)
CALL FlagError("Not implemented.",ERR,ERROR,*999)
CASE(EQUATIONS_SET_FV_SOLUTION_METHOD)
CALL FlagError("Not implemented.",ERR,ERROR,*999)
CASE(EQUATIONS_SET_GFEM_SOLUTION_METHOD)
CALL FlagError("Not implemented.",ERR,ERROR,*999)
CASE(EQUATIONS_SET_GFV_SOLUTION_METHOD)
CALL FlagError("Not implemented.",ERR,ERROR,*999)
CASE DEFAULT
LOCAL_ERROR="The specified solution method of "//TRIM(NUMBER_TO_VSTRING(SOLUTION_METHOD,"*",ERR,ERROR))//" is invalid."
CALL FlagError(LOCAL_ERROR,ERR,ERROR,*999)
END SELECT
CASE DEFAULT
LOCAL_ERROR="Equations set subtype of "//TRIM(NUMBER_TO_VSTRING(EQUATIONS_SET%SPECIFICATION(3),"*",ERR,ERROR))// &
& " is not valid for a bioelectrics finite elasticity equation type of a multi physics equations set class."
CALL FlagError(LOCAL_ERROR,ERR,ERROR,*999)
END SELECT
ELSE
CALL FlagError("Equations set is not associated.",ERR,ERROR,*999)
ENDIF
EXITS("BioelectricFiniteElasticity_EquationsSetSolutionMethodSet")
RETURN
999 ERRORS("BioelectricFiniteElasticity_EquationsSetSolutionMethodSet",ERR,ERROR)
EXITS("BioelectricFiniteElasticity_EquationsSetSolutionMethodSet")
RETURN 1
END SUBROUTINE BioelectricFiniteElasticity_EquationsSetSolutionMethodSet
!
!================================================================================================================================
!
!>Sets up the bioelectrics finite elasticity equation.
SUBROUTINE BioelectricFiniteElasticity_EquationsSetSetup(EQUATIONS_SET,EQUATIONS_SET_SETUP,ERR,ERROR,*)
!Argument variables
TYPE(EQUATIONS_SET_TYPE), POINTER :: EQUATIONS_SET !<A pointer to the equations set to setup
TYPE(EQUATIONS_SET_SETUP_TYPE), INTENT(INOUT) :: EQUATIONS_SET_SETUP !<The equations set setup information
INTEGER(INTG), INTENT(OUT) :: ERR !<The error code
TYPE(VARYING_STRING), INTENT(OUT) :: ERROR !<The error string
ENTERS("BioelectricFiniteElasticity_EquationsSetSetup",ERR,ERROR,*999)
CALL FlagError("BioelectricFiniteElasticity_EquationsSetSetup is not implemented.",ERR,ERROR,*999)
EXITS("BioelectricFiniteElasticity_EquationsSetSetup")
RETURN
999 ERRORS("BioelectricFiniteElasticity_EquationsSetSetup",ERR,ERROR)
EXITS("BioelectricFiniteElasticity_EquationsSetSetup")
RETURN 1
END SUBROUTINE BioelectricFiniteElasticity_EquationsSetSetup
!
!================================================================================================================================
!
!>Calculates the element stiffness matrices and RHS for a bioelectrics finite elasticity equation finite element equations set.
SUBROUTINE BioelectricFiniteElasticity_FiniteElementCalculate(EQUATIONS_SET,ELEMENT_NUMBER,ERR,ERROR,*)
!Argument variables
TYPE(EQUATIONS_SET_TYPE), POINTER :: EQUATIONS_SET !<A pointer to the equations set to perform the finite element calculations on
INTEGER(INTG), INTENT(IN) :: ELEMENT_NUMBER !<The element number to calculate
INTEGER(INTG), INTENT(OUT) :: ERR !<The error code
TYPE(VARYING_STRING), INTENT(OUT) :: ERROR !<The error string
ENTERS("BioelectricFiniteElasticity_FiniteElementCalculate",ERR,ERROR,*999)
CALL FlagError("BioelectricFiniteElasticity_FiniteElementCalculate is not implemented.",ERR,ERROR,*999)
EXITS("BioelectricFiniteElasticity_FiniteElementCalculate")
RETURN
999 ERRORS("BioelectricFiniteElasticity_FiniteElementCalculate",ERR,ERROR)
EXITS("BioelectricFiniteElasticity_FiniteElementCalculate")
RETURN 1
END SUBROUTINE BioelectricFiniteElasticity_FiniteElementCalculate
!
!================================================================================================================================
!
!>Sets the problem specification for a bioelectric finite elasticity problem type .
SUBROUTINE BioelectricFiniteElasticity_ProblemSpecificationSet(problem,problemSpecification,err,error,*)
!Argument variables
TYPE(PROBLEM_TYPE), POINTER :: problem !<A pointer to the problem to set the problem specification for
INTEGER(INTG), INTENT(IN) :: problemSpecification(:) !<The problem specification to set
INTEGER(INTG), INTENT(OUT) :: err !<The error code
TYPE(VARYING_STRING), INTENT(OUT) :: error !<The error string
!Local Variables
TYPE(VARYING_STRING) :: localError
INTEGER(INTG) :: problemSubtype
ENTERS("BioelectricFiniteElasticity_ProblemSpecificationSet",err,error,*999)
IF(ASSOCIATED(PROBLEM)) THEN
IF(SIZE(problemSpecification,1)==3) THEN
problemSubtype=problemSpecification(3)
SELECT CASE(problemSubtype)
CASE(PROBLEM_GUDUNOV_MONODOMAIN_SIMPLE_ELASTICITY_SUBTYPE, &
& PROBLEM_GUDUNOV_MONODOMAIN_1D3D_ELASTICITY_SUBTYPE, &
& PROBLEM_MONODOMAIN_ELASTICITY_W_TITIN_SUBTYPE, &
& PROBLEM_MONODOMAIN_ELASTICITY_VELOCITY_SUBTYPE, &
& PROBLEM_MONODOMAIN_1D3D_ACTIVE_STRAIN_SUBTYPE, &
& PROBLEM_MONODOMAIN_ELASTICITY_MUSCLE_TENDON_SUBTYPE)
!ok
CASE DEFAULT
localError="The third problem specification of "//TRIM(NumberToVstring(problemSubtype,"*",err,error))// &
& " is not valid for a bioelectric finite elasticity type of a multi physics problem."
CALL FlagError(localError,err,error,*999)
END SELECT
IF(ALLOCATED(problem%specification)) THEN
CALL FlagError("Problem specification is already allocated.",err,error,*999)
ELSE
ALLOCATE(problem%specification(3),stat=err)
IF(err/=0) CALL FlagError("Could not allocate problem specification.",err,error,*999)
END IF
problem%specification(1:3)=[PROBLEM_MULTI_PHYSICS_CLASS,PROBLEM_BIOELECTRIC_FINITE_ELASTICITY_TYPE, &
& problemSubtype]
ELSE
CALL FlagError("Bioelectric finite elasticity problem specification must have three entries.",err,error,*999)
END IF
ELSE
CALL FlagError("Problem is not associated.",err,error,*999)
END IF
EXITS("BioelectricFiniteElasticity_ProblemSpecificationSet")
RETURN
999 ERRORS("BioelectricFiniteElasticity_ProblemSpecificationSet",err,error)
EXITS("BioelectricFiniteElasticity_ProblemSpecificationSet")
RETURN 1
END SUBROUTINE BioelectricFiniteElasticity_ProblemSpecificationSet
!
!================================================================================================================================
!
!>Sets up the bioelectric finite elasticity problem.
SUBROUTINE BIOELECTRIC_FINITE_ELASTICITY_PROBLEM_SETUP(PROBLEM,PROBLEM_SETUP,ERR,ERROR,*)
!Argument variables
TYPE(PROBLEM_TYPE), POINTER :: PROBLEM !<A pointer to the problem to setup
TYPE(PROBLEM_SETUP_TYPE), INTENT(INOUT) :: PROBLEM_SETUP !<The problem setup information
INTEGER(INTG), INTENT(OUT) :: ERR !<The error code
TYPE(VARYING_STRING), INTENT(OUT) :: ERROR !<The error string
!Local Variables
TYPE(CONTROL_LOOP_TYPE), POINTER :: CONTROL_LOOP,CONTROL_LOOP_ROOT
TYPE(CONTROL_LOOP_TYPE), POINTER :: MONODOMAIN_SUB_LOOP,ELASTICITY_SUB_LOOP
TYPE(SOLVER_TYPE), POINTER :: SOLVER
TYPE(SOLVER_EQUATIONS_TYPE), POINTER :: SOLVER_EQUATIONS
TYPE(CELLML_EQUATIONS_TYPE), POINTER :: CELLML_EQUATIONS
TYPE(SOLVERS_TYPE), POINTER :: SOLVERS,MONODOMAIN_SOLVERS,ELASTICITY_SOLVERS
TYPE(VARYING_STRING) :: LOCAL_ERROR
ENTERS("BIOELECTRIC_FINITE_ELASTICITY_PROBLEM_SETUP",ERR,ERROR,*999)
NULLIFY(CONTROL_LOOP)
NULLIFY(MONODOMAIN_SUB_LOOP)
NULLIFY(ELASTICITY_SUB_LOOP)
NULLIFY(SOLVER)
NULLIFY(SOLVERS)
NULLIFY(MONODOMAIN_SOLVERS)
NULLIFY(ELASTICITY_SOLVERS)
NULLIFY(SOLVER_EQUATIONS)
NULLIFY(CELLML_EQUATIONS)
IF(ASSOCIATED(PROBLEM)) THEN
IF(.NOT.ALLOCATED(problem%specification)) THEN
CALL FlagError("Problem specification is not allocated.",err,error,*999)
ELSE IF(SIZE(problem%specification,1)<3) THEN
CALL FlagError("Problem specification must have three entries for a bioelectric-finite elasticity problem.",err,error,*999)
END IF
SELECT CASE(PROBLEM%SPECIFICATION(3))
!--------------------------------------------------------------------
! Transient Gudunov monodomain, simple finite elasticity
!--------------------------------------------------------------------
CASE(PROBLEM_GUDUNOV_MONODOMAIN_SIMPLE_ELASTICITY_SUBTYPE,PROBLEM_GUDUNOV_MONODOMAIN_1D3D_ELASTICITY_SUBTYPE, &
& PROBLEM_MONODOMAIN_ELASTICITY_W_TITIN_SUBTYPE,PROBLEM_MONODOMAIN_ELASTICITY_VELOCITY_SUBTYPE, &
& PROBLEM_MONODOMAIN_1D3D_ACTIVE_STRAIN_SUBTYPE,PROBLEM_MONODOMAIN_ELASTICITY_MUSCLE_TENDON_SUBTYPE)
SELECT CASE(PROBLEM_SETUP%SETUP_TYPE)
CASE(PROBLEM_SETUP_INITIAL_TYPE)
SELECT CASE(PROBLEM_SETUP%ACTION_TYPE)
CASE(PROBLEM_SETUP_START_ACTION)
!Do nothing
CASE(PROBLEM_SETUP_FINISH_ACTION)
!Do nothing
CASE DEFAULT
LOCAL_ERROR="The action type of "//TRIM(NUMBER_TO_VSTRING(PROBLEM_SETUP%ACTION_TYPE,"*",ERR,ERROR))// &
& " for a setup type of "//TRIM(NUMBER_TO_VSTRING(PROBLEM_SETUP%SETUP_TYPE,"*",ERR,ERROR))// &
& " is invalid for a bioelectrics finite elasticity equation."
CALL FlagError(LOCAL_ERROR,ERR,ERROR,*999)
END SELECT
CASE(PROBLEM_SETUP_CONTROL_TYPE)
SELECT CASE(PROBLEM_SETUP%ACTION_TYPE)
CASE(PROBLEM_SETUP_START_ACTION)
!Set up a time control loop
CALL CONTROL_LOOP_CREATE_START(PROBLEM,CONTROL_LOOP,ERR,ERROR,*999)
CALL CONTROL_LOOP_TYPE_SET(CONTROL_LOOP,PROBLEM_CONTROL_TIME_LOOP_TYPE,ERR,ERROR,*999)
CALL CONTROL_LOOP_NUMBER_OF_SUB_LOOPS_SET(CONTROL_LOOP,2,ERR,ERROR,*999)
CALL CONTROL_LOOP_OUTPUT_TYPE_SET(CONTROL_LOOP,CONTROL_LOOP_PROGRESS_OUTPUT,ERR,ERROR,*999)
!Set up the control sub loop for monodomain
CALL CONTROL_LOOP_SUB_LOOP_GET(CONTROL_LOOP,1,MONODOMAIN_SUB_LOOP,ERR,ERROR,*999)
CALL CONTROL_LOOP_LABEL_SET(MONODOMAIN_SUB_LOOP,'MONODOMAIN_TIME_LOOP',ERR,ERROR,*999)
CALL CONTROL_LOOP_TYPE_SET(MONODOMAIN_SUB_LOOP,PROBLEM_CONTROL_TIME_LOOP_TYPE,ERR,ERROR,*999)
CALL CONTROL_LOOP_OUTPUT_TYPE_SET(MONODOMAIN_SUB_LOOP,CONTROL_LOOP_PROGRESS_OUTPUT,ERR,ERROR,*999)
IF(PROBLEM%specification(3)==PROBLEM_MONODOMAIN_ELASTICITY_VELOCITY_SUBTYPE) THEN
!Set up the control sub loop for finite elasicity
CALL CONTROL_LOOP_SUB_LOOP_GET(CONTROL_LOOP,2,ELASTICITY_SUB_LOOP,ERR,ERROR,*999)
CALL CONTROL_LOOP_LABEL_SET(ELASTICITY_SUB_LOOP,'ELASTICITY_WHILE_LOOP',ERR,ERROR,*999)
CALL CONTROL_LOOP_TYPE_SET(ELASTICITY_SUB_LOOP,PROBLEM_CONTROL_WHILE_LOOP_TYPE,ERR,ERROR,*999)
CALL CONTROL_LOOP_OUTPUT_TYPE_SET(ELASTICITY_SUB_LOOP,CONTROL_LOOP_PROGRESS_OUTPUT,ERR,ERROR,*999)
ELSE
!Set up the control sub loop for finite elasicity
CALL CONTROL_LOOP_SUB_LOOP_GET(CONTROL_LOOP,2,ELASTICITY_SUB_LOOP,ERR,ERROR,*999)
CALL CONTROL_LOOP_LABEL_SET(ELASTICITY_SUB_LOOP,'ELASTICITY_LOAD_INCREMENT_LOOP',ERR,ERROR,*999)
CALL CONTROL_LOOP_TYPE_SET(ELASTICITY_SUB_LOOP,PROBLEM_CONTROL_LOAD_INCREMENT_LOOP_TYPE,ERR,ERROR,*999)
CALL CONTROL_LOOP_OUTPUT_TYPE_SET(ELASTICITY_SUB_LOOP,CONTROL_LOOP_PROGRESS_OUTPUT,ERR,ERROR,*999)
ENDIF
CASE(PROBLEM_SETUP_FINISH_ACTION)
!Finish the control loops
CONTROL_LOOP_ROOT=>PROBLEM%CONTROL_LOOP
CALL CONTROL_LOOP_GET(CONTROL_LOOP_ROOT,CONTROL_LOOP_NODE,CONTROL_LOOP,ERR,ERROR,*999)
CALL CONTROL_LOOP_CREATE_FINISH(CONTROL_LOOP,ERR,ERROR,*999)
!Sub-loops are finished when parent is finished
CASE DEFAULT
LOCAL_ERROR="The action type of "//TRIM(NUMBER_TO_VSTRING(PROBLEM_SETUP%ACTION_TYPE,"*",ERR,ERROR))// &
& " for a setup type of "//TRIM(NUMBER_TO_VSTRING(PROBLEM_SETUP%SETUP_TYPE,"*",ERR,ERROR))// &
& " is invalid for a bioelectrics finite elasticity equation."
CALL FlagError(LOCAL_ERROR,ERR,ERROR,*999)
END SELECT
CASE(PROBLEM_SETUP_SOLVERS_TYPE)
!Get the control loop
CONTROL_LOOP_ROOT=>PROBLEM%CONTROL_LOOP
CALL CONTROL_LOOP_GET(CONTROL_LOOP_ROOT,CONTROL_LOOP_NODE,CONTROL_LOOP,ERR,ERROR,*999)
SELECT CASE(PROBLEM_SETUP%ACTION_TYPE)
CASE(PROBLEM_SETUP_START_ACTION)
!Get the monodomain sub loop
CALL CONTROL_LOOP_SUB_LOOP_GET(CONTROL_LOOP,1,MONODOMAIN_SUB_LOOP,ERR,ERROR,*999)
!Start the solvers creation
CALL SOLVERS_CREATE_START(MONODOMAIN_SUB_LOOP,MONODOMAIN_SOLVERS,ERR,ERROR,*999)
CALL SOLVERS_NUMBER_SET(MONODOMAIN_SOLVERS,2,ERR,ERROR,*999)
!Set the first solver to be a differential-algebraic equations solver
NULLIFY(SOLVER)
CALL SOLVERS_SOLVER_GET(MONODOMAIN_SOLVERS,1,SOLVER,ERR,ERROR,*999)
CALL SOLVER_TYPE_SET(SOLVER,SOLVER_DAE_TYPE,ERR,ERROR,*999)
CALL SOLVER_LABEL_SET(SOLVER,"ODE Solver",ERR,ERROR,*999)
CALL SOLVER_LIBRARY_TYPE_SET(SOLVER,SOLVER_CMISS_LIBRARY,ERR,ERROR,*999)
!Set the second solver to be a dynamic solver
NULLIFY(SOLVER)
CALL SOLVERS_SOLVER_GET(MONODOMAIN_SOLVERS,2,SOLVER,ERR,ERROR,*999)
CALL SOLVER_TYPE_SET(SOLVER,SOLVER_DYNAMIC_TYPE,ERR,ERROR,*999)
CALL SOLVER_DYNAMIC_ORDER_SET(SOLVER,SOLVER_DYNAMIC_FIRST_ORDER,ERR,ERROR,*999)
CALL SOLVER_LABEL_SET(SOLVER,"Parabolic solver",ERR,ERROR,*999)
CALL SOLVER_DYNAMIC_DEGREE_SET(SOLVER,SOLVER_DYNAMIC_FIRST_DEGREE,ERR,ERROR,*999)
CALL SOLVER_DYNAMIC_SCHEME_SET(SOLVER,SOLVER_DYNAMIC_CRANK_NICOLSON_SCHEME,ERR,ERROR,*999)
CALL SOLVER_DYNAMIC_RESTART_SET(SOLVER,.TRUE.,ERR,ERROR,*999)
CALL SOLVER_LIBRARY_TYPE_SET(SOLVER,SOLVER_CMISS_LIBRARY,ERR,ERROR,*999)
!Get the finite elasticity sub loop
CALL CONTROL_LOOP_SUB_LOOP_GET(CONTROL_LOOP,2,ELASTICITY_SUB_LOOP,ERR,ERROR,*999)
!Start the solvers creation
CALL SOLVERS_CREATE_START(ELASTICITY_SUB_LOOP,ELASTICITY_SOLVERS,ERR,ERROR,*999)
CALL SOLVERS_NUMBER_SET(ELASTICITY_SOLVERS,1,ERR,ERROR,*999)
!Set the finite elasticity solver to be a nonlinear solver
NULLIFY(SOLVER)
CALL SOLVERS_SOLVER_GET(ELASTICITY_SOLVERS,1,SOLVER,ERR,ERROR,*999)
CALL SOLVER_TYPE_SET(SOLVER,SOLVER_NONLINEAR_TYPE,ERR,ERROR,*999)
CALL SOLVER_LIBRARY_TYPE_SET(SOLVER,SOLVER_PETSC_LIBRARY,ERR,ERROR,*999)
CASE(PROBLEM_SETUP_FINISH_ACTION)
!Get the monodomain solvers
CALL CONTROL_LOOP_SUB_LOOP_GET(CONTROL_LOOP,1,MONODOMAIN_SUB_LOOP,ERR,ERROR,*999)
CALL CONTROL_LOOP_SOLVERS_GET(MONODOMAIN_SUB_LOOP,MONODOMAIN_SOLVERS,ERR,ERROR,*999)
!Finish the solvers creation
CALL SOLVERS_CREATE_FINISH(MONODOMAIN_SOLVERS,ERR,ERROR,*999)
!Get the finite elasticity solvers
CALL CONTROL_LOOP_SUB_LOOP_GET(CONTROL_LOOP,2,ELASTICITY_SUB_LOOP,ERR,ERROR,*999)
CALL CONTROL_LOOP_SOLVERS_GET(ELASTICITY_SUB_LOOP,ELASTICITY_SOLVERS,ERR,ERROR,*999)
!Finish the solvers creation
CALL SOLVERS_CREATE_FINISH(ELASTICITY_SOLVERS,ERR,ERROR,*999)
CASE DEFAULT
LOCAL_ERROR="The action type of "//TRIM(NUMBER_TO_VSTRING(PROBLEM_SETUP%ACTION_TYPE,"*",ERR,ERROR))// &
& " for a setup type of "//TRIM(NUMBER_TO_VSTRING(PROBLEM_SETUP%SETUP_TYPE,"*",ERR,ERROR))// &
& " is invalid for a bioelectrics finite elasticity equation."
CALL FlagError(LOCAL_ERROR,ERR,ERROR,*999)
END SELECT
CASE(PROBLEM_SETUP_SOLVER_EQUATIONS_TYPE)
SELECT CASE(PROBLEM_SETUP%ACTION_TYPE)
CASE(PROBLEM_SETUP_START_ACTION)
!Get the control loop and solvers
CONTROL_LOOP_ROOT=>PROBLEM%CONTROL_LOOP
CALL CONTROL_LOOP_GET(CONTROL_LOOP_ROOT,CONTROL_LOOP_NODE,CONTROL_LOOP,ERR,ERROR,*999)
!Get the monodomain sub loop and solvers
CALL CONTROL_LOOP_SUB_LOOP_GET(CONTROL_LOOP,1,MONODOMAIN_SUB_LOOP,ERR,ERROR,*999)
CALL CONTROL_LOOP_SOLVERS_GET(MONODOMAIN_SUB_LOOP,MONODOMAIN_SOLVERS,ERR,ERROR,*999)
!Create the solver equations for the second (parabolic) solver
NULLIFY(SOLVER)
CALL SOLVERS_SOLVER_GET(MONODOMAIN_SOLVERS,2,SOLVER,ERR,ERROR,*999)
CALL SOLVER_EQUATIONS_CREATE_START(SOLVER,SOLVER_EQUATIONS,ERR,ERROR,*999)
CALL SOLVER_EQUATIONS_LINEARITY_TYPE_SET(SOLVER_EQUATIONS,SOLVER_EQUATIONS_LINEAR,ERR,ERROR,*999)
CALL SOLVER_EQUATIONS_TIME_DEPENDENCE_TYPE_SET(SOLVER_EQUATIONS,SOLVER_EQUATIONS_FIRST_ORDER_DYNAMIC,ERR,ERROR,*999)
CALL SOLVER_EQUATIONS_SPARSITY_TYPE_SET(SOLVER_EQUATIONS,SOLVER_SPARSE_MATRICES,ERR,ERROR,*999)
!Get the finite elasticity sub loop and solvers
CALL CONTROL_LOOP_SUB_LOOP_GET(CONTROL_LOOP,2,ELASTICITY_SUB_LOOP,ERR,ERROR,*999)
CALL CONTROL_LOOP_SOLVERS_GET(ELASTICITY_SUB_LOOP,ELASTICITY_SOLVERS,ERR,ERROR,*999)
!Get the finite elasticity solver and create the finite elasticity solver equations
NULLIFY(SOLVER)
NULLIFY(SOLVER_EQUATIONS)
CALL SOLVERS_SOLVER_GET(ELASTICITY_SOLVERS,1,SOLVER,ERR,ERROR,*999)
CALL SOLVER_EQUATIONS_CREATE_START(SOLVER,SOLVER_EQUATIONS,ERR,ERROR,*999)
CALL SOLVER_EQUATIONS_LINEARITY_TYPE_SET(SOLVER_EQUATIONS,SOLVER_EQUATIONS_NONLINEAR,ERR,ERROR,*999)
CALL SOLVER_EQUATIONS_TIME_DEPENDENCE_TYPE_SET(SOLVER_EQUATIONS,SOLVER_EQUATIONS_STATIC,ERR,ERROR,*999)
CALL SOLVER_EQUATIONS_SPARSITY_TYPE_SET(SOLVER_EQUATIONS,SOLVER_SPARSE_MATRICES,ERR,ERROR,*999)
CASE(PROBLEM_SETUP_FINISH_ACTION)
!Get the control loop
CONTROL_LOOP_ROOT=>PROBLEM%CONTROL_LOOP
CALL CONTROL_LOOP_GET(CONTROL_LOOP_ROOT,CONTROL_LOOP_NODE,CONTROL_LOOP,ERR,ERROR,*999)
!Get the monodomain sub loop and solvers
CALL CONTROL_LOOP_SUB_LOOP_GET(CONTROL_LOOP,1,MONODOMAIN_SUB_LOOP,ERR,ERROR,*999)
CALL CONTROL_LOOP_SOLVERS_GET(MONODOMAIN_SUB_LOOP,MONODOMAIN_SOLVERS,ERR,ERROR,*999)
!Get the solver equations for the second (parabolic) solver
NULLIFY(SOLVER)
NULLIFY(SOLVER_EQUATIONS)
CALL SOLVERS_SOLVER_GET(MONODOMAIN_SOLVERS,2,SOLVER,ERR,ERROR,*999)
CALL SOLVER_SOLVER_EQUATIONS_GET(SOLVER,SOLVER_EQUATIONS,ERR,ERROR,*999)
!Finish the solver equations creation
CALL SOLVER_EQUATIONS_CREATE_FINISH(SOLVER_EQUATIONS,ERR,ERROR,*999)
!Get the finite elasticity sub loop and solvers
CALL CONTROL_LOOP_SUB_LOOP_GET(CONTROL_LOOP,2,ELASTICITY_SUB_LOOP,ERR,ERROR,*999)
CALL CONTROL_LOOP_SOLVERS_GET(ELASTICITY_SUB_LOOP,ELASTICITY_SOLVERS,ERR,ERROR,*999)
!Finish the creation of the finite elasticity solver equations
NULLIFY(SOLVER)
NULLIFY(SOLVER_EQUATIONS)
CALL SOLVERS_SOLVER_GET(ELASTICITY_SOLVERS,1,SOLVER,ERR,ERROR,*999)
CALL SOLVER_SOLVER_EQUATIONS_GET(SOLVER,SOLVER_EQUATIONS,ERR,ERROR,*999)
CALL SOLVER_EQUATIONS_CREATE_FINISH(SOLVER_EQUATIONS,ERR,ERROR,*999)
CASE DEFAULT
LOCAL_ERROR="The action type of "//TRIM(NUMBER_TO_VSTRING(PROBLEM_SETUP%ACTION_TYPE,"*",ERR,ERROR))// &
& " for a setup type of "//TRIM(NUMBER_TO_VSTRING(PROBLEM_SETUP%SETUP_TYPE,"*",ERR,ERROR))// &
& " is invalid for a bioelectrics finite elasticity equation."
CALL FlagError(LOCAL_ERROR,ERR,ERROR,*999)
END SELECT
CASE(PROBLEM_SETUP_CELLML_EQUATIONS_TYPE)
SELECT CASE(PROBLEM_SETUP%ACTION_TYPE)
CASE(PROBLEM_SETUP_START_ACTION)
!Get the control loop
CONTROL_LOOP_ROOT=>PROBLEM%CONTROL_LOOP
CALL CONTROL_LOOP_GET(CONTROL_LOOP_ROOT,CONTROL_LOOP_NODE,CONTROL_LOOP,ERR,ERROR,*999)
CALL CONTROL_LOOP_SUB_LOOP_GET(CONTROL_LOOP,1,MONODOMAIN_SUB_LOOP,ERR,ERROR,*999)
!Get the solvers
CALL CONTROL_LOOP_SOLVERS_GET(MONODOMAIN_SUB_LOOP,MONODOMAIN_SOLVERS,ERR,ERROR,*999)
!Create the CellML equations for the first DAE solver
CALL SOLVERS_SOLVER_GET(MONODOMAIN_SOLVERS,1,SOLVER,ERR,ERROR,*999)
CALL CELLML_EQUATIONS_CREATE_START(SOLVER,CELLML_EQUATIONS,ERR,ERROR,*999)
CASE(PROBLEM_SETUP_FINISH_ACTION)
!Get the control loop
CONTROL_LOOP_ROOT=>PROBLEM%CONTROL_LOOP
CALL CONTROL_LOOP_GET(CONTROL_LOOP_ROOT,CONTROL_LOOP_NODE,CONTROL_LOOP,ERR,ERROR,*999)
CALL CONTROL_LOOP_SUB_LOOP_GET(CONTROL_LOOP,1,MONODOMAIN_SUB_LOOP,ERR,ERROR,*999)
!Get the solvers
CALL CONTROL_LOOP_SOLVERS_GET(MONODOMAIN_SUB_LOOP,MONODOMAIN_SOLVERS,ERR,ERROR,*999)
!Get the CellML equations for the first DAE solver
CALL SOLVERS_SOLVER_GET(MONODOMAIN_SOLVERS,1,SOLVER,ERR,ERROR,*999)
CALL SOLVER_CELLML_EQUATIONS_GET(SOLVER,CELLML_EQUATIONS,ERR,ERROR,*999)
!Finish the CellML equations creation
CALL CELLML_EQUATIONS_CREATE_FINISH(CELLML_EQUATIONS,ERR,ERROR,*999)
CASE DEFAULT
LOCAL_ERROR="The action type of "//TRIM(NUMBER_TO_VSTRING(PROBLEM_SETUP%ACTION_TYPE,"*",ERR,ERROR))// &
& " for a setup type of "//TRIM(NUMBER_TO_VSTRING(PROBLEM_SETUP%SETUP_TYPE,"*",ERR,ERROR))// &
& " is invalid for a bioelectrics finite elasticity equation."
CALL FlagError(LOCAL_ERROR,ERR,ERROR,*999)
END SELECT
CASE DEFAULT
LOCAL_ERROR="The setup type of "//TRIM(NUMBER_TO_VSTRING(PROBLEM_SETUP%SETUP_TYPE,"*",ERR,ERROR))// &
& " is invalid for a bioelectrics finite elasticity equation."
CALL FlagError(LOCAL_ERROR,ERR,ERROR,*999)
END SELECT
CASE DEFAULT
LOCAL_ERROR="The problem subtype of "//TRIM(NUMBER_TO_VSTRING(PROBLEM%SPECIFICATION(3),"*",ERR,ERROR))// &
& " does not equal a transient monodomain quasistatic finite elasticity equation subtype."
CALL FlagError(LOCAL_ERROR,ERR,ERROR,*999)
END SELECT
ELSE
CALL FlagError("Problem is not associated.",ERR,ERROR,*999)
ENDIF
EXITS("BIOELECTRIC_FINITE_ELASTICITY_PROBLEM_SETUP")
RETURN
999 ERRORSEXITS("BIOELECTRIC_FINITE_ELASTICITY_PROBLEM_SETUP",ERR,ERROR)
RETURN 1
END SUBROUTINE BIOELECTRIC_FINITE_ELASTICITY_PROBLEM_SETUP
!
!================================================================================================================================
!
!>Sets up the bioelectrics finite elasticity problem pre-solve.
SUBROUTINE BIOELECTRIC_FINITE_ELASTICITY_PRE_SOLVE(CONTROL_LOOP,SOLVER,ERR,ERROR,*)
!Argument variables
TYPE(CONTROL_LOOP_TYPE), POINTER :: CONTROL_LOOP !<A pointer to the control loop to solve.
TYPE(SOLVER_TYPE), POINTER :: SOLVER !<A pointer to the solver
INTEGER(INTG), INTENT(OUT) :: ERR !<The error code
TYPE(VARYING_STRING), INTENT(OUT) :: ERROR !<The error string
!Local Variables
TYPE(VARYING_STRING) :: LOCAL_ERROR
ENTERS("BIOELECTRIC_FINITE_ELASTICITY_PRE_SOLVE",ERR,ERROR,*999)
IF(ASSOCIATED(CONTROL_LOOP)) THEN
IF(ASSOCIATED(SOLVER)) THEN
IF(ASSOCIATED(CONTROL_LOOP%PROBLEM)) THEN
IF(.NOT.ALLOCATED(control_loop%problem%specification)) THEN
CALL FlagError("Problem specification is not allocated.",err,error,*999)
ELSE IF(SIZE(control_loop%problem%specification,1)<3) THEN
CALL FlagError("Problem specification must have three entries for a bioelectric-finite elasticity problem.", &
& err,error,*999)
END IF
SELECT CASE(CONTROL_LOOP%PROBLEM%SPECIFICATION(3))
CASE(PROBLEM_GUDUNOV_MONODOMAIN_SIMPLE_ELASTICITY_SUBTYPE,PROBLEM_GUDUNOV_MONODOMAIN_1D3D_ELASTICITY_SUBTYPE, &
& PROBLEM_MONODOMAIN_ELASTICITY_W_TITIN_SUBTYPE,PROBLEM_MONODOMAIN_1D3D_ACTIVE_STRAIN_SUBTYPE, &
& PROBLEM_MONODOMAIN_ELASTICITY_MUSCLE_TENDON_SUBTYPE)
SELECT CASE(CONTROL_LOOP%LOOP_TYPE)
CASE(PROBLEM_CONTROL_TIME_LOOP_TYPE)
CALL BIODOMAIN_PRE_SOLVE(SOLVER,ERR,ERROR,*999)
CASE(PROBLEM_CONTROL_LOAD_INCREMENT_LOOP_TYPE)
CALL FiniteElasticity_PreSolve(solver,err,error,*999)
CASE DEFAULT
LOCAL_ERROR="Control loop loop type "//TRIM(NUMBER_TO_VSTRING(CONTROL_LOOP%LOOP_TYPE,"*",ERR,ERROR))// &
& " is not valid for a bioelectrics finite elasticity type of a multi physics problem class."
CALL FlagError(LOCAL_ERROR,ERR,ERROR,*999)
END SELECT
CASE(PROBLEM_MONODOMAIN_ELASTICITY_VELOCITY_SUBTYPE)
SELECT CASE(CONTROL_LOOP%LOOP_TYPE)
CASE(PROBLEM_CONTROL_TIME_LOOP_TYPE)
CALL BIODOMAIN_PRE_SOLVE(SOLVER,ERR,ERROR,*999)
CASE(PROBLEM_CONTROL_WHILE_LOOP_TYPE)
CALL FiniteElasticity_PreSolve(solver,err,error,*999)
CASE DEFAULT
LOCAL_ERROR="Control loop loop type "//TRIM(NUMBER_TO_VSTRING(CONTROL_LOOP%LOOP_TYPE,"*",ERR,ERROR))// &
& " is not valid for a bioelectrics finite elasticity type of a multi physics problem class."
CALL FlagError(LOCAL_ERROR,ERR,ERROR,*999)
END SELECT
CASE DEFAULT
LOCAL_ERROR="Problem subtype "//TRIM(NUMBER_TO_VSTRING(CONTROL_LOOP%PROBLEM%SPECIFICATION(3),"*",ERR,ERROR))// &
& " is not valid for a bioelectrics finite elasticity type of a multi physics problem class."
CALL FlagError(LOCAL_ERROR,ERR,ERROR,*999)
END SELECT
ELSE
CALL FlagError("Problem is not associated.",ERR,ERROR,*999)
ENDIF
ELSE
CALL FlagError("Solver is not associated.",ERR,ERROR,*999)
ENDIF
ELSE
CALL FlagError("Control loop is not associated.",ERR,ERROR,*999)
ENDIF
EXITS("BIOELECTRIC_FINITE_ELASTICITY_PRE_SOLVE")
RETURN
999 ERRORSEXITS("BIOELECTRIC_FINITE_ELASTICITY_PRE_SOLVE",ERR,ERROR)
RETURN 1
END SUBROUTINE BIOELECTRIC_FINITE_ELASTICITY_PRE_SOLVE
!
!================================================================================================================================
!
!>Sets up the bioelectrics finite elasticity problem post solve.
SUBROUTINE BIOELECTRIC_FINITE_ELASTICITY_POST_SOLVE(CONTROL_LOOP,SOLVER,ERR,ERROR,*)
!Argument variables
TYPE(CONTROL_LOOP_TYPE), POINTER :: CONTROL_LOOP !<A pointer to the control loop to solve.
TYPE(SOLVER_TYPE), POINTER :: SOLVER!<A pointer to the solver
INTEGER(INTG), INTENT(OUT) :: ERR !<The error code
TYPE(VARYING_STRING), INTENT(OUT) :: ERROR !<The error string
!Local Variables
TYPE(VARYING_STRING) :: LOCAL_ERROR
ENTERS("BIOELECTRIC_FINITE_ELASTICITY_POST_SOLVE",ERR,ERROR,*999)
IF(ASSOCIATED(CONTROL_LOOP)) THEN
IF(ASSOCIATED(SOLVER)) THEN
IF(ASSOCIATED(CONTROL_LOOP%PROBLEM)) THEN
IF(.NOT.ALLOCATED(control_loop%problem%specification)) THEN
CALL FlagError("Problem specification is not allocated.",err,error,*999)
ELSE IF(SIZE(control_loop%problem%specification,1)<3) THEN
CALL FlagError("Problem specification must have three entries for a bioelectric-finite elasticity problem.", &
& err,error,*999)
END IF
SELECT CASE(CONTROL_LOOP%PROBLEM%SPECIFICATION(3))
CASE(PROBLEM_GUDUNOV_MONODOMAIN_SIMPLE_ELASTICITY_SUBTYPE,PROBLEM_GUDUNOV_MONODOMAIN_1D3D_ELASTICITY_SUBTYPE, &
& PROBLEM_MONODOMAIN_ELASTICITY_W_TITIN_SUBTYPE,PROBLEM_MONODOMAIN_ELASTICITY_VELOCITY_SUBTYPE, &
& PROBLEM_MONODOMAIN_1D3D_ACTIVE_STRAIN_SUBTYPE,PROBLEM_MONODOMAIN_ELASTICITY_MUSCLE_TENDON_SUBTYPE)
SELECT CASE(SOLVER%SOLVE_TYPE)
CASE(SOLVER_DAE_TYPE)
CALL BIOELECTRIC_POST_SOLVE(SOLVER,ERR,ERROR,*999)
CASE(SOLVER_DYNAMIC_TYPE)
CALL BIOELECTRIC_POST_SOLVE(SOLVER,ERR,ERROR,*999)
CASE(SOLVER_NONLINEAR_TYPE)
CALL FiniteElasticity_PostSolve(solver,err,error,*999)
CASE DEFAULT
LOCAL_ERROR="Solver solve type "//TRIM(NUMBER_TO_VSTRING(SOLVER%SOLVE_TYPE,"*",ERR,ERROR))// &
& " is not valid for a bioelectrics finite elasticity type of a multi physics problem class."
CALL FlagError(LOCAL_ERROR,ERR,ERROR,*999)
END SELECT
CASE DEFAULT
LOCAL_ERROR="Problem subtype "//TRIM(NUMBER_TO_VSTRING(CONTROL_LOOP%PROBLEM%SPECIFICATION(3),"*",ERR,ERROR))// &
& " is not valid for a bioelectrics finite elasticity type of a multi physics problem class."
CALL FlagError(LOCAL_ERROR,ERR,ERROR,*999)
END SELECT
ELSE
CALL FlagError("Problem is not associated.",ERR,ERROR,*999)
ENDIF
ELSE
CALL FlagError("Solver is not associated.",ERR,ERROR,*999)
ENDIF
ELSE
CALL FlagError("Control loop is not associated.",ERR,ERROR,*999)
ENDIF
EXITS("BIOELECTRIC_FINITE_ELASTICITY_POST_SOLVE")
RETURN
999 ERRORSEXITS("BIOELECTRIC_FINITE_ELASTICITY_POST_SOLVE",ERR,ERROR)
RETURN 1
END SUBROUTINE BIOELECTRIC_FINITE_ELASTICITY_POST_SOLVE
!
!================================================================================================================================
!
!>Sets up the bioelectrics finite elasticity problem pre-control loop.
SUBROUTINE BioelectricFiniteElasticity_ControlLoopPreLoop(CONTROL_LOOP,ERR,ERROR,*)
!Argument variables
TYPE(CONTROL_LOOP_TYPE), POINTER :: CONTROL_LOOP !<A pointer to the control loop to solve.
INTEGER(INTG), INTENT(OUT) :: ERR !<The error code
TYPE(VARYING_STRING), INTENT(OUT) :: ERROR !<The error string
!Local Variables
TYPE(PROBLEM_TYPE), POINTER :: PROBLEM
TYPE(VARYING_STRING) :: LOCAL_ERROR
ENTERS("BioelectricFiniteElasticity_ControlLoopPreLoop",ERR,ERROR,*999)
IF(ASSOCIATED(CONTROL_LOOP)) THEN
IF(CONTROL_LOOP%NUMBER_OF_SUB_LOOPS==0) THEN
PROBLEM=>CONTROL_LOOP%PROBLEM
IF(ASSOCIATED(PROBLEM)) THEN
IF(.NOT.ALLOCATED(control_loop%problem%specification)) THEN
CALL FlagError("Problem specification is not allocated.",err,error,*999)
ELSE IF(SIZE(control_loop%problem%specification,1)<3) THEN
CALL FlagError("Problem specification must have three entries for a bioelectric-finite elasticity problem.", &
& err,error,*999)
END IF
SELECT CASE(PROBLEM%SPECIFICATION(2))
CASE(PROBLEM_BIOELECTRIC_FINITE_ELASTICITY_TYPE)
SELECT CASE(CONTROL_LOOP%LOOP_TYPE)
CASE(PROBLEM_CONTROL_TIME_LOOP_TYPE)
!do nothing ???
CASE(PROBLEM_CONTROL_LOAD_INCREMENT_LOOP_TYPE)
SELECT CASE(PROBLEM%SPECIFICATION(3))
CASE(PROBLEM_GUDUNOV_MONODOMAIN_1D3D_ELASTICITY_SUBTYPE,PROBLEM_MONODOMAIN_1D3D_ACTIVE_STRAIN_SUBTYPE)
CALL BioelectricFiniteElasticity_IndependentFieldInterpolate(CONTROL_LOOP,ERR,ERROR,*999)
CASE(PROBLEM_MONODOMAIN_ELASTICITY_W_TITIN_SUBTYPE,PROBLEM_MONODOMAIN_ELASTICITY_MUSCLE_TENDON_SUBTYPE)
CALL BIOELECTRIC_FINITE_ELASTICITY_COMPUTE_TITIN(CONTROL_LOOP,ERR,ERROR,*999)
CALL BioelectricFiniteElasticity_IndependentFieldInterpolate(CONTROL_LOOP,ERR,ERROR,*999)
CASE DEFAULT
LOCAL_ERROR="The third problem specification of "// &
& TRIM(NUMBER_TO_VSTRING(PROBLEM%SPECIFICATION(3),"*",ERR,ERROR))// &
& " is not valid for bioelectric finite elasticity problem."
CALL FlagError(LOCAL_ERROR,ERR,ERROR,*999)
END SELECT
CALL FiniteElasticity_ControlTimeLoopPreLoop(CONTROL_LOOP,ERR,ERROR,*999)
CASE(PROBLEM_CONTROL_WHILE_LOOP_TYPE)
SELECT CASE(PROBLEM%SPECIFICATION(3))
CASE(PROBLEM_MONODOMAIN_ELASTICITY_VELOCITY_SUBTYPE)
IF(CONTROL_LOOP%WHILE_LOOP%ITERATION_NUMBER==1) THEN
CALL BioelectricFiniteElasticity_IndependentFieldInterpolate(CONTROL_LOOP,ERR,ERROR,*999)
ENDIF
CALL BioelectricFiniteElasticity_ComputeFibreStretch(CONTROL_LOOP,ERR,ERROR,*999)
CALL BioelectricFiniteElasticity_ForceLengthVelocityRelation(CONTROL_LOOP,ERR,ERROR,*999)
CASE DEFAULT
LOCAL_ERROR="The third problem specification of "// &
& TRIM(NUMBER_TO_VSTRING(PROBLEM%SPECIFICATION(3),"*",ERR,ERROR))// &
& " is not valid for bioelectric finite elasticity problem."
CALL FlagError(LOCAL_ERROR,ERR,ERROR,*999)
END SELECT
CALL FiniteElasticity_ControlTimeLoopPreLoop(CONTROL_LOOP,ERR,ERROR,*999)
CASE DEFAULT
LOCAL_ERROR="Control loop loop type "//TRIM(NUMBER_TO_VSTRING(CONTROL_LOOP%LOOP_TYPE,"*",ERR,ERROR))// &
& " is not valid for bioelectric finite elasticity problem type."
CALL FlagError(LOCAL_ERROR,ERR,ERROR,*999)
END SELECT
CASE DEFAULT
LOCAL_ERROR="The second problem specification of "// &
& TRIM(NUMBER_TO_VSTRING(PROBLEM%SPECIFICATION(2),"*",ERR,ERROR))// &
& " is not valid for a multi physics problem."
CALL FlagError(LOCAL_ERROR,ERR,ERROR,*999)
END SELECT
ELSE
CALL FlagError("Control loop problem is not associated.",ERR,ERROR,*999)
ENDIF
ELSE
!the main time loop - do nothing!
ENDIF
ELSE
CALL FlagError("Control loop is not associated.",ERR,ERROR,*999)
ENDIF
EXITS("BioelectricFiniteElasticity_ControlLoopPreLoop")
RETURN
999 ERRORS("BioelectricFiniteElasticity_ControlLoopPreLoop",ERR,ERROR)
EXITS("BioelectricFiniteElasticity_ControlLoopPreLoop")
RETURN 1
END SUBROUTINE BioelectricFiniteElasticity_ControlLoopPreLoop
!
!================================================================================================================================
!
!>Computes the fibre stretch for bioelectrics finite elasticity problems.
SUBROUTINE BioelectricFiniteElasticity_ComputeFibreStretch(CONTROL_LOOP,ERR,ERROR,*)
!Argument variables
TYPE(CONTROL_LOOP_TYPE), POINTER :: CONTROL_LOOP !<A pointer to the control loop to solve.
INTEGER(INTG), INTENT(OUT) :: ERR !<The error code
TYPE(VARYING_STRING), INTENT(OUT) :: ERROR !<The error string
!Local Variables
TYPE(EQUATIONS_SET_TYPE), POINTER :: EQUATIONS_SET
TYPE(FIELD_TYPE), POINTER :: dependentField,fibreField,geometricField,independentField
TYPE(SOLVER_TYPE), POINTER :: SOLVER
TYPE(SOLVER_EQUATIONS_TYPE), POINTER :: SOLVER_EQUATIONS
TYPE(SOLVER_MAPPING_TYPE), POINTER :: SOLVER_MAPPING
TYPE(SOLVERS_TYPE), POINTER :: SOLVERS
TYPE(VARYING_STRING) :: LOCAL_ERROR
TYPE(BASIS_TYPE), POINTER :: DEPENDENT_BASIS
TYPE(EquationsType), POINTER :: equations
TYPE(EquationsMappingVectorType), POINTER :: vectorMapping
TYPE(EquationsMappingNonlinearType), POINTER :: nonlinearMapping
TYPE(EquationsVectorType), POINTER :: vectorEquations
TYPE(QUADRATURE_SCHEME_TYPE), POINTER :: DEPENDENT_QUADRATURE_SCHEME
TYPE(FIELD_INTERPOLATION_PARAMETERS_TYPE), POINTER :: GEOMETRIC_INTERPOLATION_PARAMETERS, &
& FIBRE_INTERPOLATION_PARAMETERS,DEPENDENT_INTERPOLATION_PARAMETERS
TYPE(FIELD_INTERPOLATED_POINT_TYPE), POINTER :: GEOMETRIC_INTERPOLATED_POINT,FIBRE_INTERPOLATED_POINT, &
& DEPENDENT_INTERPOLATED_POINT
TYPE(FIELD_INTERPOLATED_POINT_METRICS_TYPE), POINTER :: GEOMETRIC_INTERPOLATED_POINT_METRICS, &
& DEPENDENT_INTERPOLATED_POINT_METRICS
TYPE(BASIS_TYPE), POINTER :: GEOMETRIC_BASIS
TYPE(DECOMPOSITION_TYPE), POINTER :: DECOMPOSITION
TYPE(FIELD_VARIABLE_TYPE), POINTER :: FIELD_VAR_U1
INTEGER(INTG) :: DEPENDENT_NUMBER_OF_GAUSS_POINTS
INTEGER(INTG) :: MESH_COMPONENT_NUMBER,NUMBER_OF_ELEMENTS,FIELD_VAR_TYPE
INTEGER(INTG) :: equations_set_idx,gauss_idx,dof_idx,element_idx,idx
REAL(DP) :: DZDNU(3,3),DZDNUT(3,3),AZL(3,3)
ENTERS("BioelectricFiniteElasticity_ComputeFibreStretch",ERR,ERROR,*999)
NULLIFY(SOLVERS)
NULLIFY(SOLVER)
NULLIFY(SOLVER_EQUATIONS)
NULLIFY(FIELD_VAR_U1)
NULLIFY(DEPENDENT_BASIS)
NULLIFY(EQUATIONS)
NULLIFY(dependentField,fibreField,geometricField,independentField)
NULLIFY(DEPENDENT_QUADRATURE_SCHEME)
NULLIFY(GEOMETRIC_INTERPOLATION_PARAMETERS,FIBRE_INTERPOLATION_PARAMETERS,DEPENDENT_INTERPOLATION_PARAMETERS)
NULLIFY(GEOMETRIC_INTERPOLATED_POINT,FIBRE_INTERPOLATED_POINT,DEPENDENT_INTERPOLATED_POINT)
NULLIFY(GEOMETRIC_INTERPOLATED_POINT_METRICS,DEPENDENT_INTERPOLATED_POINT_METRICS)
IF(ASSOCIATED(CONTROL_LOOP)) THEN
IF(CONTROL_LOOP%NUMBER_OF_SUB_LOOPS==0) THEN
SELECT CASE(CONTROL_LOOP%LOOP_TYPE)
CASE(PROBLEM_CONTROL_TIME_LOOP_TYPE)
!do nothing
CASE(PROBLEM_CONTROL_LOAD_INCREMENT_LOOP_TYPE)
!do nothing
CASE(PROBLEM_CONTROL_WHILE_LOOP_TYPE)
CALL CONTROL_LOOP_SOLVERS_GET(CONTROL_LOOP,SOLVERS,ERR,ERROR,*999)
CALL SOLVERS_SOLVER_GET(SOLVERS,1,SOLVER,ERR,ERROR,*999)
CALL SOLVER_SOLVER_EQUATIONS_GET(SOLVER,SOLVER_EQUATIONS,ERR,ERROR,*999)
!Loop over the equations sets associated with the solver
IF(ASSOCIATED(SOLVER_EQUATIONS)) THEN
SOLVER_MAPPING=>SOLVER_EQUATIONS%SOLVER_MAPPING
IF(ASSOCIATED(SOLVER_MAPPING)) THEN
DO equations_set_idx=1,SOLVER_MAPPING%NUMBER_OF_EQUATIONS_SETS
EQUATIONS_SET=>SOLVER_MAPPING%EQUATIONS_SETS(equations_set_idx)%ptr
IF(ASSOCIATED(EQUATIONS_SET)) THEN
geometricField=>EQUATIONS_SET%GEOMETRY%GEOMETRIC_FIELD
IF(.NOT.ASSOCIATED(geometricField)) THEN
LOCAL_ERROR="Geometric field is not associated for equations set index "// &
& TRIM(NUMBER_TO_VSTRING(equations_set_idx,"*",ERR,ERROR))//" in the solver mapping."
CALL FlagError(LOCAL_ERROR,ERR,ERROR,*999)
ENDIF
dependentField=>EQUATIONS_SET%DEPENDENT%DEPENDENT_FIELD
IF(.NOT.ASSOCIATED(dependentField)) THEN
LOCAL_ERROR="Dependent field is not associated for equations set index "// &
& TRIM(NUMBER_TO_VSTRING(equations_set_idx,"*",ERR,ERROR))//" in the solver mapping."
CALL FlagError(LOCAL_ERROR,ERR,ERROR,*999)
ENDIF
independentField=>EQUATIONS_SET%INDEPENDENT%independent_Field
IF(.NOT.ASSOCIATED(independentField)) THEN
LOCAL_ERROR="Independent field is not associated for equations set index "// &
& TRIM(NUMBER_TO_VSTRING(equations_set_idx,"*",ERR,ERROR))//" in the solver mapping."
CALL FlagError(LOCAL_ERROR,ERR,ERROR,*999)
ENDIF
EQUATIONS=>EQUATIONS_SET%EQUATIONS
IF(.NOT.ASSOCIATED(EQUATIONS)) THEN
LOCAL_ERROR="Equations is not associated for equations set index "// &
& TRIM(NUMBER_TO_VSTRING(equations_set_idx,"*",ERR,ERROR))//" in the solver mapping."
CALL FlagError(LOCAL_ERROR,ERR,ERROR,*999)
ENDIF
NULLIFY(vectorEquations)
CALL Equations_VectorEquationsGet(equations,vectorEquations,err,error,*999)
NULLIFY(vectorMapping)
CALL EquationsVector_VectorMappingGet(vectorEquations,vectorMapping,err,error,*999)
NULLIFY(nonlinearMapping)
CALL EquationsMappingVector_NonlinearMappingGet(vectorMapping,nonlinearMapping,err,error,*999)
fibreField=>equations%interpolation%fibreField
IF(.NOT.ASSOCIATED(fibreField)) THEN
LOCAL_ERROR="Fibre field is not associated for equations set index "// &
& TRIM(NUMBER_TO_VSTRING(equations_set_idx,"*",ERR,ERROR))//" in the solver mapping."
CALL FlagError(LOCAL_ERROR,ERR,ERROR,*999)
ENDIF
ELSE
LOCAL_ERROR="Equations set is not associated for equations set index "// &
& TRIM(NUMBER_TO_VSTRING(equations_set_idx,"*",ERR,ERROR))//" in the solver mapping."
CALL FlagError(LOCAL_ERROR,ERR,ERROR,*999)
ENDIF
ENDDO !equations_set_idx
ELSE
CALL FlagError("Solver equations solver mapping is not associated.",ERR,ERROR,*999)
ENDIF
ELSE
CALL FlagError("Solver solver equations are not associated.",ERR,ERROR,*999)
ENDIF
CALL Field_VariableGet(independentField,FIELD_U1_VARIABLE_TYPE,FIELD_VAR_U1,ERR,ERROR,*999)
DECOMPOSITION=>dependentField%DECOMPOSITION
MESH_COMPONENT_NUMBER=DECOMPOSITION%MESH_COMPONENT_NUMBER
IF(CONTROL_LOOP%WHILE_LOOP%ITERATION_NUMBER==1) THEN
!copy the old fibre stretch to the previous values parameter set
CALL FIELD_PARAMETER_SETS_COPY(independentField,FIELD_U1_VARIABLE_TYPE, &
& FIELD_VALUES_SET_TYPE,FIELD_PREVIOUS_VALUES_SET_TYPE,1.0_DP,ERR,ERROR,*999)
ENDIF
NUMBER_OF_ELEMENTS=dependentField%GEOMETRIC_FIELD%DECOMPOSITION%TOPOLOGY%ELEMENTS%NUMBER_OF_ELEMENTS
!Grab interpolation parameters
FIELD_VAR_TYPE=nonlinearMapping%residualVariables(1)%ptr%VARIABLE_TYPE
DEPENDENT_INTERPOLATION_PARAMETERS=>equations%interpolation%dependentInterpParameters(FIELD_VAR_TYPE)%ptr
GEOMETRIC_INTERPOLATION_PARAMETERS=>equations%interpolation%geometricInterpParameters(FIELD_U_VARIABLE_TYPE)%ptr
FIBRE_INTERPOLATION_PARAMETERS=>equations%interpolation%fibreInterpParameters(FIELD_U_VARIABLE_TYPE)%ptr
!Point interpolation pointer
GEOMETRIC_INTERPOLATED_POINT=>equations%interpolation%geometricInterpPoint(FIELD_U_VARIABLE_TYPE)%ptr
GEOMETRIC_INTERPOLATED_POINT_METRICS=>equations%interpolation%geometricInterpPointMetrics(FIELD_U_VARIABLE_TYPE)%ptr
FIBRE_INTERPOLATED_POINT=>equations%interpolation%fibreInterpPoint(FIELD_U_VARIABLE_TYPE)%ptr
DEPENDENT_INTERPOLATED_POINT=>equations%interpolation%dependentInterpPoint(FIELD_VAR_TYPE)%ptr
DEPENDENT_INTERPOLATED_POINT_METRICS=>equations%interpolation%dependentInterpPointMetrics(FIELD_VAR_TYPE)%ptr
!loop over the elements of the finite elasticity mesh (internal and boundary elements)
DO element_idx=1,NUMBER_OF_ELEMENTS
DEPENDENT_BASIS=>DECOMPOSITION%DOMAIN(MESH_COMPONENT_NUMBER)%ptr%TOPOLOGY%ELEMENTS%ELEMENTS(element_idx)%BASIS
DEPENDENT_QUADRATURE_SCHEME=>DEPENDENT_BASIS%QUADRATURE%QUADRATURE_SCHEME_MAP(BASIS_DEFAULT_QUADRATURE_SCHEME)%ptr
DEPENDENT_NUMBER_OF_GAUSS_POINTS=DEPENDENT_QUADRATURE_SCHEME%NUMBER_OF_GAUSS
GEOMETRIC_BASIS=>geometricField%DECOMPOSITION%DOMAIN(geometricField%DECOMPOSITION%MESH_COMPONENT_NUMBER)%ptr% &
& TOPOLOGY%ELEMENTS%ELEMENTS(element_idx)%BASIS
!Initialise tensors and matrices
DZDNU=0.0_DP
DO idx=1,3
DZDNU(idx,idx)=1.0_DP
ENDDO
CALL FIELD_INTERPOLATION_PARAMETERS_ELEMENT_GET(FIELD_VALUES_SET_TYPE,element_idx, &
& GEOMETRIC_INTERPOLATION_PARAMETERS,ERR,ERROR,*999)
CALL FIELD_INTERPOLATION_PARAMETERS_ELEMENT_GET(FIELD_VALUES_SET_TYPE,element_idx, &
& FIBRE_INTERPOLATION_PARAMETERS,ERR,ERROR,*999)
CALL FIELD_INTERPOLATION_PARAMETERS_ELEMENT_GET(FIELD_VALUES_SET_TYPE,element_idx, &
& DEPENDENT_INTERPOLATION_PARAMETERS,ERR,ERROR,*999)
!Loop over gauss points
DO gauss_idx=1,DEPENDENT_NUMBER_OF_GAUSS_POINTS
!Interpolate dependent, geometric, fibre and materials fields
CALL FIELD_INTERPOLATE_GAUSS(FIRST_PART_DERIV,BASIS_DEFAULT_QUADRATURE_SCHEME,gauss_idx, &
& DEPENDENT_INTERPOLATED_POINT,ERR,ERROR,*999)
CALL FIELD_INTERPOLATED_POINT_METRICS_CALCULATE(DEPENDENT_BASIS%NUMBER_OF_XI, &
& DEPENDENT_INTERPOLATED_POINT_METRICS,ERR,ERROR,*999)
CALL FIELD_INTERPOLATE_GAUSS(FIRST_PART_DERIV,BASIS_DEFAULT_QUADRATURE_SCHEME,gauss_idx, &
& GEOMETRIC_INTERPOLATED_POINT,ERR,ERROR,*999)
CALL FIELD_INTERPOLATED_POINT_METRICS_CALCULATE(GEOMETRIC_BASIS%NUMBER_OF_XI, &
& GEOMETRIC_INTERPOLATED_POINT_METRICS,ERR,ERROR,*999)
CALL FIELD_INTERPOLATE_GAUSS(FIRST_PART_DERIV,BASIS_DEFAULT_QUADRATURE_SCHEME,gauss_idx, &
& FIBRE_INTERPOLATED_POINT,ERR,ERROR,*999)
!Calculate F=dZ/dNU, the deformation gradient tensor at the gauss point
CALL FiniteElasticity_GaussDeformationGradientTensor(DEPENDENT_INTERPOLATED_POINT_METRICS, &
& GEOMETRIC_INTERPOLATED_POINT_METRICS,FIBRE_INTERPOLATED_POINT,DZDNU,ERR,ERROR,*999)
!compute C=F^T F
CALL MatrixTranspose(DZDNU,DZDNUT,ERR,ERROR,*999)
CALL MatrixProduct(DZDNUT,DZDNU,AZL,ERR,ERROR,*999)
!store the fibre stretch lambda_f, i.e., sqrt(C_11) or AZL(1,1)
dof_idx=FIELD_VAR_U1%COMPONENTS(1)%PARAM_TO_DOF_MAP%GAUSS_POINT_PARAM2DOF_MAP%GAUSS_POINTS(gauss_idx,element_idx)
CALL FIELD_PARAMETER_SET_UPDATE_LOCAL_DOF(independentField,FIELD_U1_VARIABLE_TYPE,FIELD_VALUES_SET_TYPE, &
& dof_idx,SQRT(AZL(1,1)),ERR,ERROR,*999)
ENDDO !gauss_idx
ENDDO !element_idx
!now the ghost elements -- get the relevant info from the other computational nodes
CALL FIELD_PARAMETER_SET_UPDATE_START(independentField,FIELD_U1_VARIABLE_TYPE,FIELD_VALUES_SET_TYPE,ERR,ERROR,*999)
CALL FIELD_PARAMETER_SET_UPDATE_FINISH(independentField,FIELD_U1_VARIABLE_TYPE,FIELD_VALUES_SET_TYPE,ERR,ERROR,*999)
CASE DEFAULT
LOCAL_ERROR="Control loop type "//TRIM(NUMBER_TO_VSTRING(CONTROL_LOOP%LOOP_TYPE,"*",ERR,ERROR))// &
& " is not valid for a bioelectrics finite elasticity type of a multi physics problem class."
CALL FlagError(LOCAL_ERROR,ERR,ERROR,*999)
END SELECT
ELSE !the control loop contains subloops
!do nothing
ENDIF
ELSE
CALL FlagError("Control loop is not associated.",ERR,ERROR,*999)
ENDIF
EXITS("BioelectricFiniteElasticity_ComputeFibreStretch")
RETURN
999 ERRORS("BioelectricFiniteElasticity_ComputeFibreStretch",ERR,ERROR)
EXITS("BioelectricFiniteElasticity_ComputeFibreStretch")
RETURN 1
END SUBROUTINE BioelectricFiniteElasticity_ComputeFibreStretch
!
!================================================================================================================================
!
!>Computes the bioelectrics finite elasticity force-length and force_velocity relations.
SUBROUTINE BioelectricFiniteElasticity_ForceLengthVelocityRelation(CONTROL_LOOP,ERR,ERROR,*)
!Argument variables
TYPE(CONTROL_LOOP_TYPE), POINTER :: CONTROL_LOOP !<A pointer to the control loop to solve.
INTEGER(INTG), INTENT(OUT) :: ERR !<The error code
TYPE(VARYING_STRING), INTENT(OUT) :: ERROR !<The error string
!Local Variables