-
Notifications
You must be signed in to change notification settings - Fork 1
/
nano_PF_singleVariant_calibratedInstability.cc
2328 lines (1937 loc) · 86.1 KB
/
nano_PF_singleVariant_calibratedInstability.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
/* Author: Hamed Babaei, 2016 */
#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/base/conditional_ostream.h>
#include <deal.II/base/std_cxx11/shared_ptr.h>
#include <deal.II/base/logstream.h>
#include <deal.II/base/utilities.h>
#include <deal.II/base/index_set.h>
#include <deal.II/lac/generic_linear_algebra.h>
#include <deal.II/lac/vector.h>
#include <deal.II/lac/full_matrix.h>
#include <deal.II/lac/sparse_matrix.h>
#include <deal.II/lac/dynamic_sparsity_pattern.h>
#include <deal.II/lac/solver_cg.h>
#include <deal.II/lac/precondition.h>
#include <deal.II/lac/constraint_matrix.h>
#include <deal.II/lac/precondition_selector.h>
#include <deal.II/lac/sparse_direct.h>
#include <deal.II/lac/sparsity_tools.h>
#include <deal.II/lac/constrained_linear_operator.h>
#include <deal.II/lac/trilinos_sparse_matrix.h>
#include <deal.II/lac/trilinos_vector.h>
#include <deal.II/lac/trilinos_precondition.h>
#include <deal.II/lac/trilinos_solver.h>
#include <deal.II/lac/petsc_parallel_sparse_matrix.h>
#include <deal.II/lac/petsc_parallel_vector.h>
#include <deal.II/lac/petsc_solver.h>
#include <deal.II/lac/petsc_precondition.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/tria.h>
#include <deal.II/grid/tria_boundary_lib.h>
#include <deal.II/grid/tria_accessor.h>
#include <deal.II/grid/tria_iterator.h>
#include <deal.II/dofs/dof_renumbering.h>
#include <deal.II/dofs/dof_tools.h>
#include <deal.II/dofs/dof_handler.h>
#include <deal.II/dofs/dof_accessor.h>
#include <deal.II/fe/fe_q.h>
#include <deal.II/fe/fe_dgq.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/fe/mapping_q.h>
#include <deal.II/numerics/data_out.h>
#include <deal.II/numerics/vector_tools.h>
#include <deal.II/numerics/matrix_tools.h>
#include <deal.II/distributed/tria.h>
#include <deal.II/distributed/grid_refinement.h>
#include <iostream>
#include <fstream>
//////////////////////////////////////////////////////////////////
namespace PhaseField
{
using namespace dealii;
// INPUT OF PARAMETERS
namespace Parameters
{
struct FESystem
{
unsigned int poly_degree;
unsigned int quad_order;
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.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");
}
prm.leave_subsection();
}
////////////////////////////////////////////////////
struct Geometry
{
unsigned int refinement;
double scale;
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", "4",
Patterns::Integer(0),
"Global refinement level");
prm.declare_entry("Grid scale", "1",
Patterns::Double(0.0),
"Global grid scaling factor");
}
prm.leave_subsection();
}
void Geometry::parse_parameters(ParameterHandler &prm)
{
prm.enter_subsection("Geometry");
{
refinement = prm.get_integer("Global refinement");
scale = prm.get_double("Grid scale");
}
prm.leave_subsection();
}
/////////////////////////////////////////////////
struct Materials
{
double lambda0; // austenite phase
double mu0; // austenite phase
double lambda1; // martensite phase
double mu1; // martensite phase
double L; // interface mobility
double beta; // gradient energy coefficient
double A0; // parameter for barrier height
double theta; // temperature
double thetac; // critical temperature
double thetae; // equilibrium temperature
double thetan; // Crank-Nicolson parameter
static void
declare_parameters(ParameterHandler &prm);
void
parse_parameters(ParameterHandler &prm);
};
void Materials::declare_parameters(ParameterHandler &prm)
{
prm.enter_subsection("Material properties");
{
prm.declare_entry("lambda austenite", "144.0", /* */
Patterns::Double(),
"lambda austenite");
prm.declare_entry("mu austenite", "74.0",
Patterns::Double(0.0),
"mu austenite");
prm.declare_entry("lambda martensite", "379.0",
Patterns::Double(),
"lambda martensite");
prm.declare_entry("mu martensite", "134.0",
Patterns::Double(0.0),
"mu martensite");
prm.declare_entry("kinetic coeff", "2.6",
Patterns::Double(0.0),
"kinetic coeff");
prm.declare_entry("energy coeff", "0.1",
Patterns::Double(0.0),
"energy coeff");
prm.declare_entry("barrier height", "0.028",
Patterns::Double(),
"barrier height");
prm.declare_entry("temperature", "50.0",
Patterns::Double(),
"temperature");
prm.declare_entry("temperature crit", "-183.0",
Patterns::Double(),
"temperature crit");
prm.declare_entry("temperature eql", "215.0",
Patterns::Double(),
"temperature eql");
prm.declare_entry("crank-nicolson parameter", "0.5",
Patterns::Double(),
"crank-nicolson parameter");
}
prm.leave_subsection();
}
void Materials::parse_parameters(ParameterHandler &prm)
{
prm.enter_subsection("Material properties");
{
lambda0 = prm.get_double("lambda austenite");
mu0 = prm.get_double("mu austenite");
lambda1 = prm.get_double("lambda martensite");
mu1 = prm.get_double("mu martensite");
L = prm.get_double("kinetic coeff");
beta = prm.get_double("energy coeff");
A0 = prm.get_double("barrier height");
theta = prm.get_double("temperature");
thetac = prm.get_double("temperature crit");
thetae = prm.get_double("temperature eql");
thetan = prm.get_double("crank-nicolson parameter");
}
prm.leave_subsection();
}
/////////////////////////////////////////////////
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", "5",
Patterns::Double(),
"End time");
prm.declare_entry("Time step size", "0.01",
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();
}
///////////////////////////////////////////////////////
struct AllParameters : public FESystem,
public Geometry,
public Materials,
public Time
{
AllParameters(const std::string &input_file);
static void
declare_parameters(ParameterHandler &prm);
void
parse_parameters(ParameterHandler &prm);
};
AllParameters::AllParameters(const std::string &input_file)
{
ParameterHandler prm;
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);
Materials::declare_parameters(prm);
Time::declare_parameters(prm);
}
void AllParameters::parse_parameters(ParameterHandler &prm)
{
FESystem::parse_parameters(prm);
Geometry::parse_parameters(prm);
Materials::parse_parameters(prm);
Time::parse_parameters(prm);
}
}
// DEFINE SECOND ORDER IDENTITY, AND TWO FOURTH ORDER IDENTITY TENSORS
template <int dim>
class StandardTensors
{
public:
static const SymmetricTensor<2, dim> I;
static const SymmetricTensor<4, dim> IxI;
static const SymmetricTensor<4, dim> II;
};
template <int dim>
const SymmetricTensor<2, dim>
StandardTensors<dim>::I = unit_symmetric_tensor<dim>();
template <int dim>
const SymmetricTensor<4, dim>
StandardTensors<dim>::IxI = outer_product(I, I);
template <int dim>
const SymmetricTensor<4, dim>
StandardTensors<dim>::II = identity_tensor<dim>();
// DEFINE TIME STEP, CURRENT TIME ETC.
class Time
{
public:
Time (const double time_end,
const double delta_t)
:
timestep(0),
time_current(0.0),
time_end(time_end),
delta_t(delta_t),
strainstep(0)
{}
virtual ~Time()
{}
double current() const
{
return time_current;
}
double end() const
{
return time_end;
}
double get_delta_t() const
{
return delta_t;
}
unsigned int get_timestep() const
{
return timestep;
}
void increment()
{
time_current += delta_t;
++timestep;
}
unsigned int get_strainstep() const
{
return strainstep;
}
void strainincrement()
{
++strainstep;
}
private:
unsigned int timestep;
double time_current;
const double time_end;
const double delta_t;
unsigned int strainstep;
};
//////////////////////////////////////////////////////////
template <int dim>
class Material_Constitutive
{
public:
Material_Constitutive(const double A0,
const double theta,
const double thetac,
const double thetae)
:
//kinetic tensors
det_F(1.0),
Fe(Tensor<2, dim>()),
ge(StandardTensors<dim>::I),
Be(SymmetricTensor<2, dim>()),
I1(0.0),
Ge(StandardTensors<dim>::I),
Ee(Tensor<2, dim>()),
FeEe(Tensor<2, dim>()),
EeEe(Tensor<2, dim>()),
//elastic constants
C0(Vector<double> (9)),
C1(Vector<double> (9)),
lambda0(Vector<double> (3)),
lambda1(Vector<double> (3)),
lambda (Vector<double> (3)),
mu0(Vector<double> (3)),
mu1(Vector<double> (3)),
mu (Vector<double> (3)),
nu0(Vector<double> (3)),
nu1(Vector<double> (3)),
nu (Vector<double> (3)),
dnu (Vector<double> (3)),
dmu (Vector<double> (3)),
dlambda (Vector<double> (3)),
ddnu (Vector<double> (3)),
ddmu (Vector<double> (3)),
ddlambda (Vector<double> (3)),
//Phase-field constants
A00(A0),
theta1(theta),
thetac1(thetac),
thetae1(thetae)
{}
~Material_Constitutive()
{}
void update_material_data (const Tensor<2, dim> &F, const Tensor<2, dim> &F_e,
const double phi, const double dphi,const double ddphi,
const double dphi_doublewell, const double ddphi_doublewell,
const double dphi_deltaG,const double ddphi_deltaG,
const double phi_elastic, const double dphi_elastic, const double ddphi_elastic,
const double A, const double deltaG)
{
// Kinetic tensors
Fe=F_e;
det_F = determinant(F); // determinant of total F
ge = symmetrize(Fe*transpose(Fe)); // ge = Fe.Fe^T
Be = 0.5*(ge-StandardTensors<dim>::I); // Be = 1/2(ge-I)
I1 = trace(Be); // first invariant of Be
Ge = symmetrize(transpose(Fe)*Fe); //Ge=Fe^T.Fe
Ee = 0.5*(Ge-StandardTensors<dim>::I); //Ee=1/2(Ge-I)
FeEe= Fe*Ee;
EeEe= Ee*Ee;
// Elastic constants for Orthotropic materials
C0[0]= 167.5;//C0_11
C0[1]= 167.5;//C0_22
C0[2]= 167.5;//C0_33
C0[3]= 80.1;//C0_44
C0[4]= 80.1;//C0_55
C0[5]= 80.1;//C0_66
C0[6]= 65.0;//C0_12
C0[7]= 65.0;//C0_13
C0[8]= 65.0;//C0_23
C1[0]= 174.76;//C1_11
C1[1]= 136.68;//C1_22
C1[2]= 174.76;//C1_33
C1[3]= 60.24;//C1_44
C1[4]= 42.22;//C1_55
C1[5]= 60.24;//C1_66
C1[6]= 68.00;//C1_12
C1[7]= 102.00;//C1_13
C1[8]= 68.00;//C1_23
lambda0[0]= C0[0]+C0[8]+2*C0[3]-(C0[6]+C0[7]+2*C0[4]+2*C0[5]);
lambda0[1]= C0[1]+C0[7]+2*C0[4]-(C0[6]+C0[8]+2*C0[3]+2*C0[5]);
lambda0[2]= C0[2]+C0[6]+2*C0[5]-(C0[7]+C0[8]+2*C0[3]+2*C0[4]);
lambda1[0]= C1[0]+C1[8]+2*C1[3]-(C1[6]+C1[7]+2*C1[4]+2*C1[5]);
lambda1[1]= C1[1]+C1[7]+2*C1[4]-(C1[6]+C1[8]+2*C1[3]+2*C1[5]);
lambda1[2]= C1[2]+C1[6]+2*C1[5]-(C1[7]+C1[8]+2*C1[3]+2*C1[4]);
mu0[0]= 0.5*(C0[6]+C0[7]-C0[8]);
mu0[1]= 0.5*(C0[6]+C0[8]-C0[7]);
mu0[2]= 0.5*(C0[7]+C0[8]-C0[6]);
mu1[0]= 0.5*(C1[6]+C1[7]-C1[8]);
mu1[1]= 0.5*(C1[6]+C1[8]-C1[7]);
mu1[2]= 0.5*(C1[7]+C1[8]-C1[6]);
nu0[0]= 0.5*(C0[4]+C0[5]-C0[3]);
nu0[1]= 0.5*(C0[3]+C0[5]-C0[4]);
nu0[2]= 0.5*(C0[3]+C0[4]-C0[5]);
nu1[0]= 0.5*(C1[4]+C1[5]-C1[3]);
nu1[1]= 0.5*(C1[3]+C1[5]-C1[4]);
nu1[2]= 0.5*(C1[3]+C1[4]-C1[5]);
for (unsigned int n=0; n<3; ++n)
{
lambda[n] = lambda0[n]+(lambda1[n]-lambda0[n])*phi_elastic;
mu[n] = mu0[n]+(mu1[n]-mu0[n])*phi_elastic;
nu[n] = nu0[n]+(nu1[n]-nu0[n])*phi_elastic;
dlambda[n] = (lambda1[n]-lambda0[n])*dphi_elastic;
dmu[n] = (mu1[n]-mu0[n])*dphi_elastic;
dnu[n] = (nu1[n]-nu0[n])*dphi_elastic;
ddlambda[n] = (lambda1[n]-lambda0[n])*ddphi_elastic;
ddmu[n] = (mu1[n]-mu0[n])*ddphi_elastic;
ddnu[n] = (nu1[n]-nu0[n])*ddphi_elastic;
}
// deltaG = A00*(theta1-thetae1)/3.0; // difference in thermal energy
// A = A00*(theta1-thetac1); // compute barrier height
deltaG1= deltaG;//Si parameters
A1= A;
// derivatives of interpolation functions
dphi1 = dphi;
ddphi1 = ddphi;
dphi_doublewell1 = dphi_doublewell;
ddphi_doublewell1 = ddphi_doublewell;
dphi_deltaG1=dphi_deltaG;
ddphi_deltaG1=ddphi_deltaG;
Assert(det_F > 0, ExcInternalError());
}
SymmetricTensor<2, dim> get_tau() // tau is the Kirchhoff stress := J*Cauchy stress
{
return get_tau_kirchhoff();
}
SymmetricTensor<4, dim> get_Jc() const // comupte the fourth order modulus tensor in the deformed config.
{
return get_Jc_modulus();
}
double get_det_F() const // determinant of deformation gradient tensor
{
return det_F;
}
double get_driving_force_noStress() const // driving force excluding the transformational work
{
return get_driv_forc_noStress ();
}
double get_deri_driving_force_noStress() const // driving force excluding the transformational work
{
return get_deri_driv_forc_noStress ();
}
protected:
double det_F;
double I1;
Tensor<2, dim> Fe;
SymmetricTensor<2, dim> ge;
SymmetricTensor<2, dim> Be;
SymmetricTensor<2, dim> Ge;
Tensor<2, dim> Ee;
Tensor<2, dim> FeEe;
Tensor<2, dim> EeEe;
Vector<double> C0;
Vector<double> C1;
Vector<double> lambda0;
Vector<double> lambda1;
Vector<double> lambda;
Vector<double> mu0;
Vector<double> mu1;
Vector<double> mu;
Vector<double> nu0;
Vector<double> nu1;
Vector<double> nu;
Vector<double> dnu;
Vector<double> dmu;
Vector<double> dlambda;
Vector<double> ddnu;
Vector<double> ddmu;
Vector<double> ddlambda;
double A00;
double theta1;
double thetac1;
double thetae1;
double deltaG1;
double A1;
double dphi1;
double ddphi1;
double dphi_doublewell1;
double ddphi_doublewell1;
double dphi_deltaG1;
double ddphi_deltaG1;
double phi_elastic1;
// compute the Kirchhoff stress
SymmetricTensor<2, dim> get_tau_kirchhoff() const
{
SymmetricTensor<2, dim> kirchhoff_stress;
for (unsigned int n=0; n<dim; ++n)
for (unsigned int i=0; i<dim; ++i)
for (unsigned int j=0; j<dim; ++j)
kirchhoff_stress[i][j] += lambda[n]*Ee[n][n]*Fe[i][n]*Fe[j][n]+
mu[n]*(I1*Fe[i][n]*Fe[j][n]+Ee[n][n]*ge[i][j])+
2*nu[n]*(Fe[i][n]*FeEe[j][n]+FeEe[i][n]*Fe[j][n]);
return kirchhoff_stress;
}
// compute the modulus J*d_ijkl
SymmetricTensor<4, dim> get_Jc_modulus() const
{
SymmetricTensor<4,dim> elasticityTensor;
for (unsigned int n=0; n<dim; ++n)
for (unsigned int i=0; i<dim; ++i)
for (unsigned int j=0; j<dim; ++j)
for (unsigned int k=0; k<dim; ++k)
for (unsigned int l=0; l<dim; ++l)
elasticityTensor[i][j][k][l] +=
lambda[n]*Fe[i][n]*Fe[j][n]*Fe[k][n]*Fe[l][n]+
mu[n]*(Fe[i][n]*Fe[j][n]*ge[k][l]+ ge[i][j]*Fe[k][n]*Fe[l][n])+
nu[n]*(Fe[i][n]*ge[j][k]*Fe[l][n]+ Fe[j][n]*ge[i][k]*Fe[l][n]+
Fe[i][n]*ge[j][l]*Fe[k][n]+ Fe[j][n]*ge[i][l]*Fe[k][n]);
return elasticityTensor;
}
// compute the driving force excluding the transformational work
double get_driv_forc_noStress () const
{
double d_elastic_energy (0.0);
for (unsigned int n=0; n<dim; ++n)
d_elastic_energy += dlambda[n]*Ee[n][n]*Ee[n][n] + 2*dmu[n]*Ee[n][n]*I1 + 4*dnu[n]*EeEe[n][n];
return det_F*(-d_elastic_energy -A1*dphi_doublewell1-deltaG1*dphi1);
}
double get_deri_driv_forc_noStress () const
{
double dd_elastic_energy (0.0);
for (unsigned int n=0; n<dim; ++n)
dd_elastic_energy += 2*ddlambda[n]*Ee[n][n]*Ee[n][n] + 4*ddmu[n]*Ee[n][n]*I1 + 8*ddnu[n]*EeEe[n][n];
return det_F*(-dd_elastic_energy -A1*ddphi_doublewell1-deltaG1*ddphi1);
}
};
//////////////////////////////////////////////////////////////////////////////////////////
// updates the quadrature point history
template <int dim>
class PointHistory
{
public:
PointHistory()
:
material(NULL),
F_inv(StandardTensors<dim>::I),
tau(SymmetricTensor<2, dim>()),
Jc(SymmetricTensor<4, dim>())
{}
virtual ~PointHistory()
{
delete material;
material = NULL;
}
void setup_lqp (const Parameters::AllParameters ¶meters)
{
material = new Material_Constitutive<dim>(parameters.A0,
parameters.theta, parameters.thetac, parameters.thetae);
update_values(Tensor<2, dim>(), double (), double (), double(), double(),double() );
}
void update_values (const Tensor<2, dim> &Grad_u_n,
const double new_eta,
const double old_eta,
const double eta_laplacians,
const double thetan,
const double beta)
{
// updating order parameter
double eta = thetan*new_eta +(1-thetan)*old_eta;
// Interpolation functions and their derivatives for traditional interpolation function
//(This function is not used in the current version of theory anymore)
const double phi = 3*pow (eta, 2.0)-2*pow (eta, 3.0); // interpolation function
const double dphi = 6*eta-6*pow (eta, 2.0); // derivative of interpolation fn
const double ddphi = 6-12*eta;
// Interpolation functions and their derivatives for interpolation of double-well barriar
const double dphi_doublewell = 2*eta-6*pow (eta, 2.0)+4*pow (eta, 3.0);
const double ddphi_doublewell =2-12*eta+12*pow (eta, 2.0);
// Interpolation functions and their derivatives for interpolation of jump in thermal energy delta_G
const double dphi_deltaG=6*eta-9.3*pow (eta, 2.0)+13.2*pow (eta, 3.0)-16.5*pow (eta, 4.0)+6.6*pow (eta, 5.0);
const double ddphi_deltaG=6-18.6*eta+39.6*pow (eta, 2.0)-66*pow (eta, 3.0)+33*pow (eta, 4.0);
// Interpolation functions and their derivatives for interpolation of elastic constants
const double phi_elastic= 10*pow(eta, 3.0)-15*pow(eta, 4.0)+6*pow(eta, 5.0);
const double dphi_elastic= 30*pow(eta, 2.0)-60*pow(eta, 3.0)+30*pow(eta, 4.0);
const double ddphi_elastic= 60*eta -180*pow(eta, 2.0)+120*pow(eta, 3.0);
// Interpolation functions and their derivatives for interpolation of transformation strain
// Follwoing constants are obtaned after calibration of instability criteria with atomistic simulations
double a3;
double w3;
double a1;
double w1;
double A;
double deltaG;
a3= 3.3528;
w3= -2.6471;
a1= 0.9211*a3;
w1= 1.0406*w3;
A= 12.4192;
deltaG= 2;
const double phi_1 = a1*pow (eta, 2.0)+(10-3*a1+w1)*pow(eta,3.0)+(3*a1-2*w1-15)*pow(eta,4.0)+(6-a1+w1)*pow(eta,5.0);
const double dphi_1 = 2*a1*eta+3*(10-3*a1+w1)*pow(eta,2.0)+4*(3*a1-2*w1-15)*pow(eta,3.0)+5*(6-a1+w1)*pow(eta,4.0);
const double ddphi_1 = 2*a1+6*(10-3*a1+w1)*eta+12*(3*a1-2*w1-15)*pow(eta,2.0)+20*(6-a1+w1)*pow(eta,3.0);
const double phi_3 = a3*pow (eta, 2.0)+(10-3*a3+w3)*pow(eta,3.0)+(3*a3-2*w3-15)*pow(eta,4.0)+(6-a3+w3)*pow(eta,5.0);
const double dphi_3 = 2*a3*eta+3*(10-3*a3+w3)*pow(eta,2.0)+4*(3*a3-2*w3-15)*pow(eta,3.0)+5*(6-a3+w3)*pow(eta,4.0);
const double ddphi_3 = 2*a3+6*(10-3*a3+w3)*eta+12*(3*a3-2*w3-15)*pow(eta,2.0)+20*(6-a3+w3)*pow(eta,3.0);
Tensor<2, dim> eps_t; // transformation strain tensor
eps_t[0][0] = 0.1753;
eps_t[1][1] = -0.447;
eps_t[2][2] = 0.1753;
Tensor<2, dim> eps_t_phi;
eps_t_phi[0][0] = 0.1753*phi_1;
eps_t_phi[1][1] = -0.447 *phi_3;
eps_t_phi[2][2] = 0.1753*phi_1;
Tensor<2, dim> eps_t_dphi;
eps_t_dphi[0][0] = 0.1753*dphi_1;
eps_t_dphi[1][1] = -0.447 *dphi_3;
eps_t_dphi[2][2] = 0.1753*dphi_1;
Tensor<2, dim> eps_t_ddphi;
eps_t_ddphi[0][0] = 0.1753*ddphi_1;
eps_t_ddphi[1][1] = -0.447 *ddphi_3;
eps_t_ddphi[2][2] = 0.1753*ddphi_1;
Tensor<2, dim> dinverse_Ft;
dinverse_Ft[0][0] = -0.1753*dphi_1/pow((1+0.1753*phi_1), 2.0);
dinverse_Ft[1][1] = 0.447 *dphi_3/pow((1- 0.447*phi_3), 2.0);
dinverse_Ft[2][2] = -0.1753*dphi_1/pow((1+0.1753*phi_1), 2.0);
F = (Tensor<2, dim>(StandardTensors<dim>::I) + Grad_u_n); // total F
Ft = (Tensor<2, dim>(StandardTensors<dim>::I) + Tensor<2, dim>(eps_t_phi)); //transformation deformation gradient
Fe = F * invert(Ft); // elastic part of deformation gradient
material->update_material_data(F, Fe,
phi, dphi, ddphi,
dphi_doublewell, ddphi_doublewell,
dphi_deltaG,ddphi_deltaG,
phi_elastic, dphi_elastic, ddphi_elastic,
A,deltaG);
F_inv = invert(F);
F_inv_tr = transpose(F_inv);
tau = material->get_tau(); // extracting kirchhoff stress
Jc = material->get_Jc(); // extracting J*d_ijkl
driving_force_noStress = material->get_driving_force_noStress(); // extracting driving force with no stress
deri_driving_force_noStress = material->get_deri_driving_force_noStress();
const Tensor<2, dim> temp_tensor = F_inv*Tensor<2, dim>(tau);
const Tensor<2, dim> temp_tensor1 = temp_tensor * Fe;
const Tensor<2, dim> temp_tensor2 = Tensor<2, dim>(tau)* dinverse_Ft;
get_driv_forc = scalar_product(temp_tensor1, eps_t_dphi) +
driving_force_noStress; // computing total driving force
get_deri_driv_forc = scalar_product(temp_tensor2, eps_t_dphi)+
scalar_product(temp_tensor1, eps_t_ddphi) +
deri_driving_force_noStress;
get_driv_forc_total= get_driv_forc+beta* eta_laplacians;
Assert(determinant(F_inv) > 0, ExcInternalError());
//////////////////////////////////////////////////////////////////////////////
}
const Tensor<2, dim> &get_F() const
{
return F;
}
double get_det_F() const
{
return material->get_det_F();
}
const Tensor<2, dim> &get_F_inv() const
{
return F_inv;
}
const Tensor<2, dim> &get_F_inv_tr() const
{
return F_inv_tr;
}
const SymmetricTensor<2, dim> &get_tau() const
{
return tau;
}
const SymmetricTensor<4, dim> &get_Jc() const
{
return Jc;
}
double get_driving_force() const
{
return get_driv_forc;
}
double get_deri_driving_force() const
{
return get_deri_driv_forc;
}
double get_driving_force_total() const
{
return get_driv_forc_total;
}
private:
Material_Constitutive<dim> *material;
Tensor<2, dim> F;
Tensor<2, dim> F_inv;
Tensor<2, dim> F_inv_tr;
Tensor<2, dim> Ft;
SymmetricTensor<2, dim> tau;
SymmetricTensor<4, dim> Jc;
Tensor<2, dim> Fe;
double driving_force_noStress;
double deri_driving_force_noStress;
double dphi;
double ddphi;
double get_driv_forc;
double get_deri_driv_forc;
double get_driv_forc_total;
};
///////////////////////////////////////////////////////////////
template <int dim>
class Solid
{
public:
Solid(const std::string &input_file);
virtual
~Solid();
void
run();
private:
void make_grid();
void system_setup();
void assemble_system();
void make_constraints(const int &it_nr);
void solve_nonlinear_timestep();
unsigned int solve();
void assemble_system_eta();
void solve_nonlinear_timestep_eta();
unsigned int solve_eta();
void setup_qph();
void update_qph_incremental();
void output_results() const;
void output_resultant_stress();
Parameters::AllParameters parameters;
double vol_reference; // volume of the reference config
double vol_current; // volume of the current config
Time time; // variable of type class 'Time
MPI_Comm mpi_communicator;
parallel::distributed::Triangulation<dim> triangulation;
ConditionalOStream pcout;
const unsigned int degree; // degree of polynomial of shape functions
const FESystem<dim> fe; // fe object
DoFHandler<dim> dof_handler; // we have used two dof_handler: one for mechanics another for order parameter
const unsigned int dofs_per_cell; // no of dofs per cell for the mechanics problem
const FEValuesExtractors::Vector u_fe;
const QGauss<dim> qf_cell; // quadrature points in the cell
const QGauss<dim - 1> qf_face; // quadrature points at the face
const unsigned int n_q_points; // no of quadrature points in the cell
const unsigned int n_q_points_f; // no of quadrature points at the face
ConstraintMatrix constraints; // constraint object
FE_DGQ<dim> history_fe;
DoFHandler<dim> history_dof_handler;
std::vector<PointHistory<dim> > quadrature_point_history;
IndexSet locally_owned_dofs;
IndexSet locally_relevant_dofs;
TrilinosWrappers::SparseMatrix tangent_matrix; // tangent stiffenss matrix
TrilinosWrappers::MPI::Vector system_rhs; // system right hand side or residual of mechanics problem
TrilinosWrappers::MPI::Vector solution; // solution vector for displacement
TrilinosWrappers::MPI::Vector solution_update; // another vector containing the displacement soln
const unsigned int degree_eta; // degree of polynomial for eta
FE_Q<dim> fe_eta; // fe object for eta
DoFHandler<dim> dof_handler_eta; //another dof_handler for eta
const unsigned int dofs_per_cell_eta; // dof per eta cell
const QGauss<dim> qf_cell_eta;
const unsigned int n_q_points_eta;
ConstraintMatrix constraints_eta;
IndexSet locally_owned_dofs_eta;
IndexSet locally_relevant_dofs_eta;
PETScWrappers::MPI::SparseMatrix mass_matrix; // mass matrix out of Ginxburg-Landau eqn
PETScWrappers::MPI::SparseMatrix laplace_matrix; // Laplace matrix out of Ginxburg-Landau eqn
PETScWrappers::MPI::SparseMatrix system_matrix_eta;
PETScWrappers::MPI::SparseMatrix nl_matrix;
PETScWrappers::MPI::SparseMatrix tmp_matrix;
PETScWrappers::MPI::Vector nl_term;
PETScWrappers::MPI::Vector system_rhs_eta;
PETScWrappers::MPI::Vector solution_eta; // solution vector for eta
PETScWrappers::MPI::Vector old_solution_eta;
PETScWrappers::MPI::Vector solution_update_eta;
PETScWrappers::MPI::Vector tmp_vector; // a vector used for solving eta
Vector<double> resultant_stress;
Vector<double> static_stress;
bool apply_strain;
unsigned int load_step;
};
/////////////////////////////////////////////////////////
// defines the initial condition for the order parameter
template <int dim>
class InitialValues : public Function<dim>
{
public:
InitialValues (const int ×tep )
:
Function<dim>(),
time_step (timestep)
{}
virtual double value(const Point<dim> &p,
const unsigned int /*component = 0*/) const;
private:
const int time_step ;
};
template <int dim>
double InitialValues<dim>::value (const Point<dim> &p,
const unsigned int /*component*/) const
{
// To speed up formation of martensitic bands, random values between 0-0.01 is considered within two bands
// and 0-0.001 outside of them
if (((p[1]-p[0]<12) && (p[1]-p[0]> 8)) ||
((p[1]-p[0]<-8) && (p[1]-p[0]> -12)))
return (rand()%100+0.0)*0.0001;
else
return (rand()%100+0.0)*0.00001;
}
//////////////////////////////////////////////////
// constructor
template <int dim>
Solid<dim>::Solid(const std::string &input_file)
:
mpi_communicator (MPI_COMM_WORLD),
triangulation (mpi_communicator,
typename Triangulation<dim>::MeshSmoothing
(Triangulation<dim>::smoothing_on_refinement |
Triangulation<dim>::smoothing_on_coarsening)),
parameters(input_file),
time(parameters.end_time, parameters.delta_t),
pcout (std::cout,
(Utilities::MPI::this_mpi_process(mpi_communicator)
== 0)),
degree(parameters.poly_degree),
fe(FE_Q<dim>(parameters.poly_degree), dim),
dof_handler(triangulation),
dofs_per_cell (fe.dofs_per_cell),
u_fe(0),
qf_cell(parameters.quad_order),
qf_face(parameters.quad_order),
n_q_points (qf_cell.size()),
n_q_points_f (qf_face.size()),
degree_eta(parameters.poly_degree),
fe_eta (parameters.poly_degree),
dof_handler_eta (triangulation),
dofs_per_cell_eta(fe_eta.dofs_per_cell),
qf_cell_eta(parameters.quad_order),
n_q_points_eta (qf_cell_eta.size()),