forked from FranckLab/STAQ-DIC
-
Notifications
You must be signed in to change notification settings - Fork 2
/
main_ALDIC.m
701 lines (603 loc) · 38.1 KB
/
main_ALDIC.m
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
% ---------------------------------------------
% Augmented Lagrangian Digital Image Correlation (AL-DIC)
% working with an adaptive quadtree mesh
%
% Author: Jin Yang, PhD @Caltech
% Contact and support: [email protected] -or- [email protected]
% Date: 2015.04,06,07; 2016.03,04; 2020.11
% ---------------------------------------------
%% Section 1: Clear MATLAB environment & mex set up Spline interpolation
close all; clear; clc; clearvars -global
fprintf('------------ Section 1 Start ------------ \n')
setenv('MW_MINGW64_LOC','C:\TDM-GCC-64');
try mex -O ba_interp2.cpp; catch; end % mex set up ba_interp2.cpp script
% [Comment]: If this line reports error but it works before,
% Change line 14 to: "try mex -O ba_interp2.cpp; catch; end"
addpath("./func",'./src','./plotFiles','./func_quadtree','./func_quadtree/refinement','./plotFiles/export_fig-d966721/');
% TODO: addpath("./YOUR IMAGE FOLDER");
fprintf('------------ Section 1 Done ------------ \n \n')
%% Section 2: Load DIC parameters and set up DIC parameters
fprintf('------------ Section 2 Start ------------ \n')
% ====== Read images ======
[file_name,Img,DICpara] = ReadImage; % Load DIC raw images
% %%%%%% Uncomment the line below to change the DIC computing region (ROI) manually %%%%%%
% DICpara.gridxROIRange = [gridxROIRange1,gridxROIRange2]; DICpara.gridyROIRange = [Val1, Val2];
DICpara.gridxyROIRange.gridx = [292,764]; DICpara.gridxyROIRange.gridy = [44,1144];
% E.g., gridxROIRange = [224,918]; gridyROIRange = [787,1162];
% ====== Normalize images: fNormalized = (f-f_avg)/(f_std) ======
[ImgNormalized,DICpara.gridxyROIRange] = funNormalizeImg(Img,DICpara.gridxyROIRange);
fNormalized = ImgNormalized{1}; % Load the first referece image
% ====== Compute image gradients ======
Df = funImgGradient(fNormalized,fNormalized); % Finite difference to compute image grayscale gradients;
% ====== Initialize variable storage ======
ResultDisp = cell(length(ImgNormalized)-1,1); ResultDefGrad = cell(length(ImgNormalized)-1,1);
ResultStrain = cell(length(ImgNormalized)-1,1); ResultStress = cell(length(ImgNormalized)-1,1);
ResultFEMeshEachFrame = cell(length(ImgNormalized)-1,1); % Needs future improvment: to store FE-mesh for each frame
ResultFEMesh = cell(ceil((length(ImgNormalized)-1)/DICpara.ImgSeqIncUnit),1); % For incremental DIC mode
fprintf('------------ Section 2 Done ------------ \n \n')
%% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
% To solve each frame in an image sequence
% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
for ImgSeqNum = 2 : length(ImgNormalized)
disp(['Current image frame #: ', num2str(ImgSeqNum),'/',num2str(length(ImgNormalized))]);
gNormalized = ImgNormalized{ImgSeqNum}; % Load current deformed image frame
%% Section 3: Compute an initial guess of the unknown displacement field
fprintf('\n'); fprintf('------------ Section 3 Start ------------ \n')
% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
% This section is to find or update an initial guess of the unknown displacements.
% The key idea is to either to use a new FFT-based cross correlation peak fitting,
% or use the results from the last frame as the new initial guess for the next frame;
% Particularly in the incremental mode DIC, the reference image can also be updated, e.g.,
% " fNormalized = ImgNormalized{ImgSeqNum-mod(ImgSeqNum-1,ImgSeqIncUnit)}; "
%
% DICpara.NewFFTSearch = 0; % If you want to apply the FFT-based cross correlation to
% compute the initial guess for each frame, please make sure that "DICpara.NewFFTSearch = 0".
% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
if ImgSeqNum == 2 || DICpara.NewFFTSearch == 1 % Apply FFT-based cross correlation to compute the initial guess
% ====== Integer Search ======
[DICpara,x0temp,y0temp,u,v,cc]= IntegerSearch(fNormalized,gNormalized,file_name,DICpara);
% ====== FEM mesh set up ======
[DICmesh] = MeshSetUp(x0temp,y0temp,DICpara); clear x0temp y0temp;
% ====== Initial Value ======
U0 = Init(u,v,cc.max,DICmesh.x0,DICmesh.y0,0); %PlotuvInit; [x0temp,y0temp,u,v,cc]= IntegerSearchMg(fNormalized,gNormalized,file_name,DICpara);
% ====== Deal with incremental mode ======
fNormalizedNewIndex = ImgSeqNum-mod(ImgSeqNum-2,DICpara.ImgSeqIncUnit)-1;
if DICpara.ImgSeqIncUnit == 1, fNormalizedNewIndex = fNormalizedNewIndex-1; end
ResultFEMesh{1+floor(fNormalizedNewIndex/DICpara.ImgSeqIncUnit)} = ... % To save first mesh info
struct( 'coordinatesFEM',DICmesh.coordinatesFEM,'elementsFEM',DICmesh.elementsFEM, ...
'winsize',DICpara.winsize,'winstepsize',DICpara.winstepsize,'gridxyROIRange',DICpara.gridxyROIRange );
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
elseif mod(ImgSeqNum-2,DICpara.ImgSeqIncUnit) == 0 % To update ref image in incremental mode
fNormalizedNewIndex = ImgSeqNum-mod(ImgSeqNum-2,DICpara.ImgSeqIncUnit)-1;
if DICpara.ImgSeqIncUnit == 1, fNormalizedNewIndex = fNormalizedNewIndex-1; end
fNormalized = ImgNormalized{fNormalizedNewIndex}; % Update reference
[DICpara,DICmesh] = ReadImageRefUpdate(file_name,ImgSeqNum,ResultDisp{ImgSeqNum-2}.U,DICpara,DICmesh); % Update reference image if needed;
U0 = zeros(2*size(DICmesh.coordinatesFEM,1),1); % PlotuvInit;
ResultFEMesh{1+floor(fNormalizedNewIndex/DICpara.ImgSeqIncUnit)} = ... % To save first mesh info
struct( 'coordinatesFEM',DICmesh.coordinatesFEM,'elementsFEM',DICmesh.elementsFEM, ...
'winsize',DICpara.winsize,'winstepsize',DICpara.winstepsize,'gridxyROIRange',DICpara.gridxyROIRange );
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
else % Use the solved results from the last frame as the new initial guess
if ImgSeqNum < 7 % Import previous U for ImgSeqNum [2,6]
U0 = ResultDisp{ImgSeqNum-2}.U;
else % When ImgSeqNum > 6: POD predicts next disp U0 from previous results of (ImgSeqNum+[-5:1:-1])
nTime = 5; np = length(ResultDisp{ImgSeqNum-2}.U)/2; % "nTime" value 5 is an empirical value, can be changed.
T_data_u = zeros(nTime,np); T_data_v = zeros(nTime,np);
for tempi = 1:nTime
T_data_u(tempi,:) = ResultDisp{ImgSeqNum-(2+nTime)+tempi, 1}.U(1:2:np*2)';
T_data_v(tempi,:) = ResultDisp{ImgSeqNum-(2+nTime)+tempi, 1}.U(2:2:np*2)';
end
nB = 3; t_train = [ImgSeqNum-1-nTime:ImgSeqNum-2]'; t_pre = [ImgSeqNum-1]';
[u_pred,~,~,~] = funPOR_GPR(T_data_u,t_train,t_pre,nB);
[v_pred,~,~,~] = funPOR_GPR(T_data_v,t_train,t_pre,nB);
tempu = u_pred(1,:); tempv = v_pred(1,:);
U0 = [tempu(:),tempv(:)]'; U0 = U0(:);
% %%%%% After running the new ImgSeqNum, you can uncomment these
% %%%%% lines to compare how the initial guess has been improved.
% Plotdisp_show(U0-ResultDisp{ImgSeqNum-1}.U,DICmesh.coordinatesFEMWorld,DICmesh.elementsFEM(:,1:4));
% Plotdisp_show(ResultDisp{ImgSeqNum-2}.U-ResultDisp{ImgSeqNum-1}.U,DICmesh.coordinatesFEMWorld,DICmesh.elementsFEM(:,1:4));
end
end
% ====== Compute f(X)-g(x+u) ======
% PlotImgDiff(x0,y0,u,v,fNormalized,gNormalized);
ResultFEMeshEachFrame{ImgSeqNum-1} = struct( 'coordinatesFEM',DICmesh.coordinatesFEM,'elementsFEM',DICmesh.elementsFEM);
fprintf('------------ Section 3 Done ------------ \n \n')
%% Section 4: Subproblem 1 -or- Local ICGN Subset DIC
fprintf('------------ Section 4 Start ------------ \n')
% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
% This section is to solve the first local step in ALDIC: Subproblem 1
% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
% ====== ALStep 1 Subproblem1: Local Subset DIC ======
mu=0; beta=0; tol=1e-2; ALSolveStep=1; ALSub1Time=zeros(6,1); ALSub2Time=zeros(6,1);
ConvItPerEle=zeros(size(DICmesh.coordinatesFEM,1),6); ALSub1BadPtNum=zeros(6,1);
disp(['***** Start step',num2str(ALSolveStep),' Subproblem1 *****'])
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
% ------ Start Local DIC IC-GN iteration ------
[USubpb1,FSubpb1,HtempPar,ALSub1Timetemp,ConvItPerEletemp,LocalICGNBadPtNumtemp] = ...
LocalICGN(U0,DICmesh.coordinatesFEM,Df,fNormalized,gNormalized,DICpara,'GaussNewton',tol);
ALSub1Time(ALSolveStep) = ALSub1Timetemp; ConvItPerEle(:,ALSolveStep) = ConvItPerEletemp; ALSub1BadPtNum(ALSolveStep) = LocalICGNBadPtNumtemp; toc
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
% ------ Manually find some bad points from Local Subset ICGN step ------
% Comment these lines below if you don't need to manually remove local bad points
% %%%%% Comment START %%%%%
% close all; USubpb1World = USubpb1; USubpb1World(2:2:end) = -USubpb1(2:2:end);
% Plotuv(USubpb1World,DICmesh.x0,DICmesh.y0World);
% disp('--- Start to manually remove bad points ---')
% u = reshape(USubpb1(1:2:end),size(DICmesh.x0,1),size(DICmesh.x0,2));
% v = reshape(USubpb1(2:2:end),size(DICmesh.x0,1),size(DICmesh.x0,2));
% [u,v,~,Local_BadptRow,Local_BadptCol,RemoveOutliersList] = funRemoveOutliers(u',v',[],0.5,100); u=u';v=v';
% USubpb1(1:2:end) = reshape(u,size(DICmesh.coordinatesFEM,1),1); USubpb1(2:2:end) = reshape(v,size(DICmesh.coordinatesFEM,1),1);
% f11 = reshape(FSubpb1(1:4:end),size(DICmesh.x0,1),size(DICmesh.x0,2));
% f21 = reshape(FSubpb1(2:4:end),size(DICmesh.x0,1),size(DICmesh.x0,2));
% f12 = reshape(FSubpb1(3:4:end),size(DICmesh.x0,1),size(DICmesh.x0,2));
% f22 = reshape(FSubpb1(4:4:end),size(DICmesh.x0,1),size(DICmesh.x0,2));
% f11=f11'; f11(RemoveOutliersList) = NaN; f11 = inpaint_nans(f11,4); f11=f11';
% f21=f21'; f21(RemoveOutliersList) = NaN; f21 = inpaint_nans(f21,4); f21=f21';
% f12=f12'; f12(RemoveOutliersList) = NaN; f12 = inpaint_nans(f12,4); f12=f12';
% f22=f22'; f22(RemoveOutliersList) = NaN; f22 = inpaint_nans(f22,4); f22=f22';
% FSubpb1(1:4:end) = reshape(f11,size(DICmesh.coordinatesFEM,1),1);
% FSubpb1(2:4:end) = reshape(f21,size(DICmesh.coordinatesFEM,1),1);
% FSubpb1(3:4:end) = reshape(f12,size(DICmesh.coordinatesFEM,1),1);
% FSubpb1(4:4:end) = reshape(f22,size(DICmesh.coordinatesFEM,1),1);
% disp('--- Remove bad points done ---')
% %%%%% Comment END %%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
% ------ Plot ------
USubpb1World = USubpb1; USubpb1World(2:2:end) = -USubpb1(2:2:end); FSubpb1World = FSubpb1;
close all; Plotuv(USubpb1World,DICmesh.x0,DICmesh.y0World);
Plotdisp_show(USubpb1World,DICmesh.coordinatesFEMWorld,DICmesh.elementsFEM);
Plotstrain_show(FSubpb1World,DICmesh.coordinatesFEMWorld,DICmesh.elementsFEM);
save(['Subpb1_step',num2str(ALSolveStep)],'USubpb1','FSubpb1');
fprintf('------------ Section 4 Done ------------ \n \n')
%% Section 5: Subproblem 2 -- solve the global compatible displacement field
fprintf('------------ Section 5 Start ------------ \n'); tic;
% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
% This section is to solve the global step in ALDIC Subproblem 2
% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
% ======= ALStep 1 Subproblem 2: Global constraint =======
% ------ Smooth displacements for a better F ------
DICpara.DispFilterSize=0; DICpara.DispFilterStd=0; DICpara.StrainFilterSize=0; DICpara.StrainFilterStd=0; LevelNo=1;
FSubpb1 = funSmoothStrain(FSubpb1,DICmesh,DICpara);
% ====== Define penalty parameter ======
mu = 1e-3; udual = 0*FSubpb1; vdual = 0*USubpb1;
betaList = [1e-3,sqrt(1e-5),1e-2,sqrt(1e-3),1e-1,sqrt(1e-1)]*mean(DICpara.winstepsize).^2.*mu; % Tune beta in the betaList
Err1 = zeros(length(betaList),1); Err2 = Err1;
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
disp(['***** Start step',num2str(ALSolveStep),' Subproblem2 *****'])
% ====== Check to use FD or FE methods to solve Subpb2 step ======
if DICpara.Subpb2FDOrFEM == 1 % Using FD method
% ====== Build sparse finite difference operator ======
disp('Assemble finite difference operator D');
M = size(DICmesh.x0,1); N = size(DICmesh.x0,2);
tic; Rad = 1; D = funDerivativeOp((M-2*Rad),(N-2*Rad),DICpara.winstepsize); % D = sparse(4*(M-2*Rad)*(N-2*Rad), 2*(M-2*Rad)*(N-2*Rad));
D2 = funDerivativeOp(M,N,DICpara.winstepsize); toc
disp('Finish assembling finite difference operator D');
% ===== Solver using finite difference approximation ======
tic; a = FSubpb1-udual; b = USubpb1-vdual;
Rad = 1; [temp3,temp4] = funFDNeumannBCInd(size(DICmesh.coordinatesFEM,1),M,N,Rad); % Find coordinatesFEM that belong to (x(Rad+1:M-Rad,Rad+1:N-Rad),y(Rad+1:M-Rad,Rad+1:N-Rad))
atemp = a(temp3); btemp = b(temp4); hbar = waitbar(0,'Please wait for Subproblem 2 global step!');
for tempk = 1:length(betaList)
beta = betaList(tempk);
tempAMatrixSub2 = (beta*(D')*D) + mu*speye(2*(M-2*Rad)*(N-2*Rad));
USubpb2temp = (tempAMatrixSub2) \ (beta*D'*atemp + mu*btemp ) ;
USubpb2 = USubpb1; USubpb2(temp4) = USubpb2temp;
FSubpb2 = D2*USubpb2;
Err1(tempk) = norm(USubpb1-USubpb2,2);
Err2(tempk) = norm(FSubpb1-FSubpb2,2);
waitbar(tempk/(length(betaList)+1));
end
ErrSum = Err1+Err2*mean(DICpara.winstepsize)^2; [~,indexOfbeta] = min(ErrSum);
try % Tune the best beta by a quadratic polynomial fitting
[fitobj] = fit(log10(betaList(indexOfbeta-1:1:indexOfbeta+1))',ErrSum(indexOfbeta-1:1:indexOfbeta+1),'poly2');
p = coeffvalues(fitobj); beta = 10^(-p(2)/2/p(1));
catch, beta = betaList(indexOfbeta);
end
% Using the optimal beta to solve the ALDIC Subproblem 2 again
tempAMatrixSub2 = (beta*(D')*D) + mu*speye(2*(M-2*Rad)*(N-2*Rad));
USubpb2temp = (tempAMatrixSub2) \ (beta*D'*atemp + mu*btemp) ;
USubpb2 = USubpb1; USubpb2(temp4) = USubpb2temp;
waitbar(1); close(hbar);
%%%%%%%%%%%%%% End of using finite difference approximation %%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
else %Subpb2FDOrFEM: Using FE method
M = size(DICmesh.x0,1); N = size(DICmesh.x0,2); GaussPtOrder = 2; alpha = 0;
close all; hbar = waitbar(0,'Please wait for Subproblem 2 global step!');
% ====== Solver using finite element method ======
for tempk = 1:length(betaList)
beta = betaList(tempk);
GaussPtOrder = 2; alpha = 0; [USubpb2] = Subpb2(DICmesh,beta,mu,USubpb1,FSubpb1,udual,vdual,alpha,GaussPtOrder);
[FSubpb2] = funGlobal_NodalStrainAvg(DICmesh.coordinatesFEM,DICmesh.elementsFEM,USubpb2,GaussPtOrder);
Err1(tempk) = norm(USubpb1-USubpb2,2);
Err2(tempk) = norm(FSubpb1-FSubpb2,2);
waitbar(tempk/(length(betaList)+1));
end
Err1Norm = (Err1-mean(Err1))/std(Err1); figure, plot(Err1Norm);
Err2Norm = (Err2-mean(Err2))/std(Err2); figure, plot(Err2Norm);
ErrSum = Err1Norm+Err2Norm; figure,plot(ErrSum); [~,indexOfbeta] = min(ErrSum);
try % Tune the best beta by a quadratic polynomial fitting
[fitobj] = fit(log10(betaList(indexOfbeta-1:1:indexOfbeta+1))',ErrSum(indexOfbeta-1:1:indexOfbeta+1),'poly2');
p = coeffvalues(fitobj); beta = 10^(-p(2)/2/p(1));
catch, beta = betaList(indexOfbeta);
end
% Using the optimal beta to solve the ALDIC Subproblem 2 again
[USubpb2] = Subpb2(DICmesh,beta,mu,USubpb1,FSubpb1,udual,vdual,alpha,GaussPtOrder);
USubpb2 = full(USubpb2);
waitbar(1); close(hbar);
end
ALSub2Time(ALSolveStep) = toc; toc
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
% ------- Before computing strain, we smooth the displacement field -------
%USubpb2 = funSmoothDisp(USubpb2,coordinatesFEM,elementsFEM,x0,y0,winstepsize,DispFilterSize,DispFilterStd);
% ------- Compute strain field --------
if DICpara.Subpb2FDOrFEM == 1 %FD
FSubpb2 = D2*USubpb2; % D2 = funDerivativeOp(M,N,winstepsize);
else %FEM
[FSubpb2] = funGlobal_NodalStrainAvg(DICmesh.coordinatesFEM,DICmesh.elementsFEM,USubpb2,GaussPtOrder);
end
% ------- Smooth strain field --------
FSubpb2 = funSmoothStrain(FSubpb2,DICmesh,DICpara);
% ------- Save data ------
save(['Subpb2_step',num2str(ALSolveStep)],'USubpb2','FSubpb2');
% ------ Plot ------
USubpb2World = USubpb2; USubpb2World(2:2:end) = -USubpb2(2:2:end); FSubpb2World = FSubpb2;
close all; Plotuv(USubpb2World,DICmesh.x0,DICmesh.y0World); Plotdisp_show(USubpb2World,DICmesh.coordinatesFEMWorld,DICmesh.elementsFEM);
%Plotstrain_show(FSubpb2World,DICmesh.coordinatesFEMWorld,DICmesh.elementsFEM);
% ======= Update dual variables =======
if DICpara.Subpb2FDOrFEM == 1 %FD
udualtemp1 = (FSubpb2 - FSubpb1); udualtemp2 = udualtemp1(temp3);
vdualtemp1 = (USubpb2 - USubpb1); vdualtemp2 = vdualtemp1(temp4);
udual = zeros(4*M*N,1); vdual = zeros(2*M*N,1);
udual(temp3) = udualtemp2; vdual(temp4) = vdualtemp2;
else % FEM or other methods
udual = FSubpb2 - FSubpb1; vdual = USubpb2 - USubpb1;
end
save(['uvdual_step',num2str(ALSolveStep)],'udual','vdual');
fprintf('------------ Section 5 Done ------------ \n \n')
%% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%% Section 6: ADMM iterations
fprintf('------------ Section 6 Start ------------ \n')
% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
% This section is to run ADMM iteration: Subproblem 1 & 2
% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
% ==================== ADMM AL Loop ==========================
ALSolveStep = 1; tol2 = 1e-4; UpdateY = 1e4; CrackOrNot = 0; CrackPath1 = [0,0]; CrackPath2 = [0,0]; CrackTip = [0,0];
HPar = cell(21,1); for tempj = 1:21, HPar{tempj} = HtempPar(:,tempj); end
while (ALSolveStep < 5)
ALSolveStep = ALSolveStep + 1; % Update using the last step
%%%%%%%%%%%%%%%%%%%%%%% Subproblem 1 %%%%%%%%%%%%%%%%%%%%%%%%%
disp(['***** Start step',num2str(ALSolveStep),' Subproblem1 *****']);
tic;[USubpb1,~,ALSub1Timetemp,ConvItPerEletemp,LocalICGNBadPtNumtemp] = Subpb1( ...
USubpb2,FSubpb2,udual,vdual,DICmesh.coordinatesFEM,...
Df,fNormalized,gNormalized,mu,beta,HPar,ALSolveStep,DICpara,'GaussNewton',tol);
FSubpb1 = FSubpb2; toc
ALSub1Time(ALSolveStep) = ALSub1Timetemp; ConvItPerEle(:,ALSolveStep) = ConvItPerEletemp; ALSub1BadPtNum(ALSolveStep) = LocalICGNBadPtNumtemp;
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
% ------ Manually find some bad points from Local Subset ICGN step ------
% disp('--- Start to manually remove bad points --- \n')
% disp(' Comment codes here if you do not have bad local points \n')
% %%%%% Comment START %%%%%
% close all; Plotuv(USubpb1,DICmesh.x0,DICmesh.y0);
% u = reshape(USubpb1(1:2:end),M,N); v = reshape(USubpb1(2:2:end),M,N);
% [u,v,~,Subpb1_BadptRow,Subpb1_BadptCol] = funRemoveOutliers(u',v',[],0.5,100,Local_BadptRow,Local_BadptCol); u=u';v=v';
% disp('--- Remove bad points done ---')
% USubpb1(1:2:end) = reshape(u,size(DICmesh.coordinatesFEM,1),1);
% USubpb1(2:2:end) = reshape(v,size(DICmesh.coordinatesFEM,1),1);
% close all; Plotuv(USubpb1,DICmesh.x0,DICmesh.y0); Plotdisp_show(USubpb1,DICmesh.coordinatesFEM,DICmesh.elementsFEM);
% %%%%% Comment END %%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
save(['Subpb1_step',num2str(ALSolveStep)],'USubpb1','FSubpb1');
USubpb1 = funSmoothDisp(USubpb1,DICmesh,DICpara);
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
% ============== Subproblem 2 ==============
disp(['***** Start step',num2str(ALSolveStep),' Subproblem2 *****'])
if DICpara.Subpb2FDOrFEM == 1 %FD
% ------- using finite difference approximation --------
tic; a = FSubpb1-udual; b = USubpb1-vdual; atemp = a(temp3); btemp = b(temp4);
USubpb2temp = (tempAMatrixSub2) \ (beta*D'*atemp + mu*btemp ) ;
USubpb2 = USubpb1; USubpb2(temp4) = USubpb2temp; %toc
% ------- End of using finite difference approximation --------
else % FEM
tic; [USubpb2] = Subpb2(DICmesh,beta,mu,USubpb1,FSubpb1,udual,vdual,alpha,GaussPtOrder);
USubpb2 = full(USubpb2);
end
ALSub2Time(ALSolveStep) = toc; toc
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
% ------- Before computing strain, we smooth the displacement field -------
USubpb2 = funSmoothDisp(USubpb2,DICmesh,DICpara);
% ------- Compute strain field --------
if DICpara.Subpb2FDOrFEM == 1 %FD
FSubpb2 = D2*USubpb2; % D2 = funDerivativeOp(M,N,winstepsize);
else %FEM
GaussPtOrder = 2; [FSubpb2] = funGlobal_NodalStrainAvg(DICmesh.coordinatesFEM,DICmesh.elementsFEM,USubpb2,GaussPtOrder);
end
% ------- Smooth strain field --------
FSubpb2 = funSmoothStrain(FSubpb2,DICmesh,DICpara);
save(['Subpb2_step',num2str(ALSolveStep)],'USubpb2','FSubpb2');
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
% Compute norm of UpdateY
USubpb2_Old = load(['Subpb2_step',num2str(ALSolveStep-1)],'USubpb2');
USubpb2_New = load(['Subpb2_step',num2str(ALSolveStep)],'USubpb2');
USubpb1_Old = load(['Subpb1_step',num2str(ALSolveStep-1)],'USubpb1');
USubpb1_New = load(['Subpb1_step',num2str(ALSolveStep)],'USubpb1');
if (mod(ImgSeqNum-2,DICpara.ImgSeqIncUnit) ~= 0 && (ImgSeqNum>2)) || (ImgSeqNum < DICpara.ImgSeqIncUnit)
UpdateY = norm((USubpb2_Old.USubpb2 - USubpb2_New.USubpb2), 2)/sqrt(size(USubpb2_Old.USubpb2,1));
try
UpdateY2 = norm((USubpb1_Old.USubpb1 - USubpb1_New.USubpb1), 2)/sqrt(size(USubpb1_Old.USubpb1,1));
catch
end
end
try
disp(['Update local step = ',num2str(UpdateY2)]);
disp(['Update global step = ',num2str(UpdateY)]);
catch
end
fprintf('*********************************** \n \n');
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
% Update dual variables------------------------------
if DICpara.Subpb2FDOrFEM == 1 %FD
udualtemp1 = (FSubpb2 - FSubpb1); udualtemp2 = udualtemp1(temp3);
vdualtemp1 = (USubpb2 - USubpb1); vdualtemp2 = vdualtemp1(temp4);
udual(temp3) = udual(temp3)+udualtemp2;
vdual(temp4) = vdual(temp4)+vdualtemp2;
else %FEM
udual = FSubpb2 - FSubpb1; vdual = USubpb2 - USubpb1;
end
save(['uvdual_step',num2str(ALSolveStep)],'udual','vdual');
try
if UpdateY < tol2 || UpdateY2 < tol2
break
end
catch
end
end
fprintf('------------ Section 6 Done ------------ \n \n')
% Save data
ResultDisp{ImgSeqNum-1}.U = full(USubpb2);
ResultDisp{ImgSeqNum-1}.ALSub1BadPtNum = ALSub1BadPtNum;
ResultDefGrad{ImgSeqNum-1}.F = full(FSubpb2); % tempFoamAL;
end
% ------ Plot ------
USubpb2World = USubpb2; USubpb2World(2:2:end) = -USubpb2(2:2:end); FSubpb2World = FSubpb2;
close all; Plotuv(USubpb2World,DICmesh.x0,DICmesh.y0World);
Plotdisp_show(USubpb2World,DICmesh.coordinatesFEMWorld,DICmesh.elementsFEM);
Plotstrain_show(FSubpb2World,DICmesh.coordinatesFEMWorld,DICmesh.elementsFEM);
% ------ Save results ------
% Find img name and save all the results
[~,imgname,imgext] = fileparts(file_name{1,end});
results_name = ['results_',imgname,'_ws',num2str(DICpara.winsize),'_st',num2str(DICpara.winstepsize),'.mat'];
save(results_name, 'file_name','DICpara','DICmesh','ResultDisp','ResultDefGrad','ResultFEMesh','ResultFEMeshEachFrame','ALSub1Time','ALSub2Time','ALSolveStep');
%% Section 7: Check convergence
% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
% This section is to check convergence of ADMM
% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
fprintf('------------ Section 7 Start ------------ \n')
% ====== Check convergence ======
fprintf('***** Check convergence ***** \n');
ALSolveStep1 = min(6,ALSolveStep);
disp('==== uhat^(k) - u^(k) ====');
for ALSolveStep = 1:ALSolveStep1
USubpb2 = load(['Subpb2_step',num2str(ALSolveStep )],'USubpb2');
USubpb1 = load(['Subpb1_step',num2str(ALSolveStep )],'USubpb1');
UpdateY = norm((USubpb2.USubpb2 - USubpb1.USubpb1), 2)/sqrt(length(USubpb2.USubpb2));
disp(num2str(UpdateY));
end
disp('==== Fhat^(k) - F^(k) ====');
for ALSolveStep = 1:ALSolveStep1
FSubpb1 = load(['Subpb1_step',num2str(ALSolveStep )],'FSubpb1');
FSubpb2 = load(['Subpb2_step',num2str(ALSolveStep )],'FSubpb2');
UpdateF = norm((FSubpb1.FSubpb1 - FSubpb2.FSubpb2), 2)/sqrt(length(FSubpb1.FSubpb1));
disp(num2str(UpdateF));
end
disp('==== uhat^(k) - uhat^(k-1) ====');
for ALSolveStep = 2:ALSolveStep1
USubpb2_Old = load(['Subpb2_step',num2str(ALSolveStep-1)],'USubpb2');
USubpb2_New = load(['Subpb2_step',num2str(ALSolveStep)],'USubpb2');
UpdateY = norm((USubpb2_Old.USubpb2 - USubpb2_New.USubpb2), 2)/sqrt(length(USubpb2.USubpb2));
disp(num2str(UpdateY));
end
disp('==== udual^(k) - udual^(k-1) ====');
for ALSolveStep = 2:ALSolveStep1
uvdual_Old = load(['uvdual_step',num2str(ALSolveStep-1)],'udual');
uvdual_New = load(['uvdual_step',num2str(ALSolveStep)],'udual');
UpdateW = norm((uvdual_Old.udual - uvdual_New.udual), 2)/sqrt(length(uvdual_Old.udual));
disp(num2str(UpdateW));
end
disp('==== vdual^(k) - vdual^(k-1) ====');
for ALSolveStep = 2:ALSolveStep1
uvdual_Old = load(['uvdual_step',num2str(ALSolveStep-1)],'vdual');
uvdual_New = load(['uvdual_step',num2str(ALSolveStep)],'vdual');
Updatev = norm((uvdual_Old.vdual - uvdual_New.vdual), 2)/sqrt(length(uvdual_Old.vdual));
disp(num2str(Updatev));
end
fprintf('------------ Section 7 Done ------------ \n \n')
% ------ Delete temp files ------
%%%%% Comment START %%%%%
% Uncomment these lines to delete temporary files
% for tempi = 1:ALSolveStep
% file_name_Subpb1 = ['Subpb1_step',num2str(tempi),'.mat'];
% file_name_Subpb2 = ['Subpb2_step',num2str(tempi),'.mat'];
% file_name_dual = ['uvdual_step',num2str(tempi),'.mat'];
% delete(file_name_Subpb1); delete(file_name_Subpb2); delete(file_name_dual);
% end
%%%%% Comment END %%%%%
% ------ clear temp variables ------
clear a ALSub1BadPtNum ALSub1Timetemp atemp b btemp cc ConvItPerEletemp hbar Hbar
%% Section 8: Compute strains
fprintf('------------ Section 8 Start ------------ \n')
% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
% This section is to compute strain fields and plot disp and strain results
% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
% ------ Convert units from pixels to the physical world ------
DICpara.um2px = funParaInput('ConvertUnit');
% ------ Smooth displacements ------
DICpara.DoYouWantToSmoothOnceMore = funParaInput('SmoothDispOrNot');
% ------ Choose strain computation method ------
DICpara.MethodToComputeStrain = funParaInput('StrainMethodOp');
% ------ Choose strain type (infinitesimal, Eulerian, Green-Lagrangian) ------
DICpara.StrainType = funParaInput('StrainType');
% ------ Choose image to plot (first only, second and next images) ------
if length(ImgNormalized)==2, DICpara.Image2PlotResults = funParaInput('Image2PlotResults');
else DICpara.Image2PlotResults = 1; % Plot over current, deformed image by default
end
% ------ Save fig format ------
DICpara.MethodToSaveFig = funParaInput('SaveFigFormat');
% ------ Choose overlay image transparency ------
DICpara.OrigDICImgTransparency = 1;
if DICpara.MethodToSaveFig == 1
DICpara.OrigDICImgTransparency = funParaInput('OrigDICImgTransparency');
end
%% ------ Start main part ------
for ImgSeqNum = 2 % : length(ImgNormalized)
disp(['Current image frame #: ', num2str(ImgSeqNum),'/',num2str(length(ImgNormalized))]);
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
fNormalizedNewIndex = ImgSeqNum-mod(ImgSeqNum-2,DICpara.ImgSeqIncUnit)-1;
if DICpara.ImgSeqIncUnit > 1
FEMeshIndLast = floor(fNormalizedNewIndex/DICpara.ImgSeqIncUnit);
elseif DICpara.ImgSeqIncUnit == 1
FEMeshIndLast = floor(fNormalizedNewIndex/DICpara.ImgSeqIncUnit)-1;
end
FEMeshInd = FEMeshIndLast + 1;
if FEMeshInd == 1
USubpb2 = ResultDisp{ImgSeqNum-1}.U; %+ ResultDisp{10}.U + ResultDisp{20}.U;
coordinatesFEM = ResultFEMesh{1}.coordinatesFEM;
elementsFEM = ResultFEMesh{1}.elementsFEM;
if (ImgSeqNum-1 == 1) || (DICpara.ImgSeqIncROIUpdateOrNot==1), UFEMesh = 0*USubpb2; end
else
USubpb2 = ResultDisp{ImgSeqNum-1}.U;
if mod(ImgSeqNum-2,DICpara.ImgSeqIncUnit) == 0
coordinatesFEM = ResultFEMesh{FEMeshInd}.coordinatesFEM;
elementsFEM = ResultFEMesh{FEMeshInd}.elementsFEM;
coordinatesFEMLast = ResultFEMesh{FEMeshIndLast}.coordinatesFEM;
UFEMeshLast = ResultDisp{ImgSeqNum-2}.U + UFEMesh;
xq = coordinatesFEM(:,1); yq = coordinatesFEM(:,2);
UFEMesh = 0*USubpb2;
UFEMesh(1:2:end) = griddata(coordinatesFEMLast(:,1),coordinatesFEMLast(:,2),UFEMeshLast(1:2:end),xq,yq,'v4');
UFEMesh(2:2:end) = griddata(coordinatesFEMLast(:,1),coordinatesFEMLast(:,2),UFEMeshLast(2:2:end),xq,yq,'v4');
end
USubpb2 = USubpb2 + UFEMesh;
end
FSubpb2 = ResultDefGrad{ImgSeqNum-1}.F;
coordinatesFEM = ResultFEMeshEachFrame{ImgSeqNum-1}.coordinatesFEM;
elementsFEM = ResultFEMeshEachFrame{ImgSeqNum-1}.elementsFEM;
xList = min(coordinatesFEM(:,1)):DICpara.winstepsize:max(coordinatesFEM(:,1)); M = length(xList);
yList = min(coordinatesFEM(:,2)):DICpara.winstepsize:max(coordinatesFEM(:,2)); N = length(yList);
[x0,y0] = ndgrid(xList,yList);
x0 = x0-reshape(UFEMesh(1:2:end),size(x0,1),size(x0,2));
y0 = y0-reshape(UFEMesh(2:2:end),size(y0,1),size(y0,2));
x0World = DICpara.um2px*x0;
y0World = DICpara.um2px*y0; % Ignore this: (size(ImgNormalized{1},2)+1-y0);
coordinatesFEMWorld = DICpara.um2px*[coordinatesFEM(:,1),size(ImgNormalized{1},2)+1-coordinatesFEM(:,2)];
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
% ------ Plotting and Compute Strain-------
if size(USubpb2,1) == 1
ULocal = USubpb2_New.USubpb2; FLocal = FSubpb2.FSubpb2;
else
ULocal = USubpb2; FLocal = FSubpb2;
end
UWorld = DICpara.um2px*ULocal; UWorld(2:2:end) = -UWorld(2:2:end); FWorld = FLocal; % close all; Plotuv(UWorld,x0,y0World);
% ------ Smooth displacements ------
%prompt = 'Do you want to smooth displacement? (0-yes; 1-no)';
%DoYouWantToSmoothOnceMore = input(prompt);
SmoothTimes = 0;
try
while DICpara.DoYouWantToSmoothOnceMore == 0 && SmoothTimes < 3
ULocal = funSmoothDisp(ULocal,DICmesh,DICpara);
%%DICpara.DoYouWantToSmoothOnceMore = input(prompt);
SmoothTimes = SmoothTimes + 1;
end
catch
end
% ----- Compute strain field ------
ComputeStrain; % Compute strain
% %%%%% Add filter and plot strain field %%%%%
% %%%%% Plotstrain_Fij; %%%%%
% ------ Plot disp and strain ------
close all; Plotuv(ULocal,x0,y0);
if DICpara.OrigDICImgTransparency == 1
EdgeColorOrNot = 'NoEdgeColor'; % {'EdgeColor','NoEdgeColor'}
Plotdisp_show(UWorld,coordinatesFEMWorld,elementsFEM,DICpara,EdgeColorOrNot);
[strainxCoord,strainyCoord,dispu,dispv,dudx,dvdx,dudy,dvdy,strain_exx,strain_exy,strain_eyy,strain_principal_max, ...
strain_principal_min,strain_maxshear,strain_vonMises] = Plotstrain0( ...
UWorld,FStraintemp,Rad,x0World,y0World,size(ImgNormalized{1}),DICpara);
else % Plot over raw DIC images
if DICpara.Image2PlotResults == 0 % Plot over the first image; "file_name{1,1}" corresponds to the first image
Plotdisp(UWorld,x0World,y0World,size(ImgNormalized{1}),file_name{1,1},DICpara);
[strainxCoord,strainyCoord,dispu,dispv,dudx,dvdx,dudy,dvdy,strain_exx,strain_exy,strain_eyy,strain_principal_max, ...
strain_principal_min,strain_maxshear,strain_vonMises] = Plotstrain( ...
UWorld,FStraintemp,Rad,x0World,y0World,size(ImgNormalized{1}),file_name{1,1},DICpara);
else % Plot over second or next deformed images
Plotdisp(UWorld,x0World,y0World,size(ImgNormalized{1}),file_name{1,ImgSeqNum},DICpara);
[strainxCoord,strainyCoord,dispu,dispv,dudx,dvdx,dudy,dvdy,strain_exx,strain_exy,strain_eyy,strain_principal_max, ...
strain_principal_min,strain_maxshear,strain_vonMises] = Plotstrain( ...
UWorld,FStraintemp,Rad,x0World,y0World,size(ImgNormalized{1}),file_name{1,ImgSeqNum},DICpara);
end
end
% ----- Save strain results ------
ResultStrain{ImgSeqNum-1} = struct('strainxCoord',strainxCoord,'strainyCoord',strainyCoord, ...
'dispu',dispu,'dispv',dispv,'dudx',dudx,'dvdx',dvdx,'dudy',dudy,'dvdy',dvdy, ...
'strain_exx',strain_exx,'strain_exy',strain_exy,'strain_eyy',strain_eyy, ...
'strain_principal_max',strain_principal_max,'strain_principal_min',strain_principal_min, ...
'strain_maxshear',strain_maxshear,'strain_vonMises',strain_vonMises);
% ------ Save figures for tracked displacement and strain fields ------
SaveFigFilesDispAndStrain;
end
% ------ END of for-loop {ImgSeqNum = 2:length(ImgNormalized)} ------
fprintf('------------ Section 8 Done ------------ \n \n')
% ------ Save data again including solved strain fields ------
results_name = ['results_',imgname,'_ws',num2str(DICpara.winsize),'_st',num2str(DICpara.winstepsize),'.mat'];
save(results_name, 'file_name','DICpara','DICmesh','ResultDisp','ResultDefGrad','ResultFEMesh','ResultFEMeshEachFrame',...
'ALSub1Time','ALSub2Time','ALSolveStep','ResultStrain');
%% Section 9: Compute stress
fprintf('------------ Section 9 Start ------------ \n')
% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
% This section is to compute stress fields and plot stress fields
% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
% ------ Choose material model ------
DICpara.MaterialModel = funParaInput('MaterialModel');
% ------ Define parameters in material models ------
if (DICpara.MaterialModel == 1) || (DICpara.MaterialModel == 2) % Linear elasticity
fprintf('Define Linear elasticity parameters \n')
fprintf("Young's modulus (unit: Pa): \n"); prompt = 'Input here (e.g., 69e9): ';
DICpara.MaterialModelPara.YoungsModulus = input(prompt);
fprintf("Poisson's ratio: \n"); prompt = 'Input here (e.g., 0.3): ';
DICpara.MaterialModelPara.PoissonsRatio = input(prompt);
fprintf('------------------------------------- \n');
end
% ------ Start main part ------
for ImgSeqNum = 2 : length(ImgNormalized)
disp(['Current image frame #: ', num2str(ImgSeqNum),'/',num2str(length(ImgNormalized))]); close all;
% ------ Plot stress ------
if DICpara.OrigDICImgTransparency == 1
[stress_sxx,stress_sxy,stress_syy, stress_principal_max_xyplane, ...
stress_principal_min_xyplane, stress_maxshear_xyplane, ...
stress_maxshear_xyz3d, stress_vonMises] = Plotstress0( ...
DICpara,ResultStrain{ImgSeqNum-1},size(ImgNormalized{1}));
else % Plot over raw DIC images
if DICpara.Image2PlotResults == 0 % Plot over the first image; "file_name{1,1}" corresponds to the first image
[stress_sxx,stress_sxy,stress_syy, stress_principal_max_xyplane, ...
stress_principal_min_xyplane, stress_maxshear_xyplane, ...
stress_maxshear_xyz3d, stress_vonMises] = Plotstress( ...
DICpara,ResultStrain{ImgSeqNum-1},size(ImgNormalized{1}),file_name{1,1});
else % Plot over second or next deformed images
[stress_sxx,stress_sxy,stress_syy, stress_principal_max_xyplane, ...
stress_principal_min_xyplane, stress_maxshear_xyplane, ...
stress_maxshear_xyz3d, stress_vonMises] = Plotstress( ...
DICpara,ResultStrain{ImgSeqNum-1},size(ImgNormalized{1}),file_name{1,ImgSeqNum});
end
end
% ------ Save figures for computed stress fields ------
SaveFigFilesStress;
% ----- Save strain results ------
ResultStress{ImgSeqNum-1} = struct('stressxCoord',ResultStrain{ImgSeqNum-1}.strainxCoord,'stressyCoord',ResultStrain{ImgSeqNum-1}.strainyCoord, ...
'stress_sxx',stress_sxx,'stress_sxy',stress_sxy,'stress_syy',stress_syy, ...
'stress_principal_max_xyplane',stress_principal_max_xyplane, 'stress_principal_min_xyplane',stress_principal_min_xyplane, ...
'stress_maxshear_xyplane',stress_maxshear_xyplane,'stress_maxshear_xyz3d',stress_maxshear_xyz3d, ...
'stress_vonMises',stress_vonMises);
end
% ------ END of for-loop {ImgSeqNum = 2:length(ImgNormalized)} ------
fprintf('------------ Section 9 Done ------------ \n \n')
% ------ Save data again including solved stress fields ------
results_name = ['results_',imgname,'_ws',num2str(DICpara.winsize),'_st',num2str(DICpara.winstepsize),'.mat'];
save(results_name, 'file_name','DICpara','DICmesh','ResultDisp','ResultDefGrad','ResultFEMesh','ResultFEMeshEachFrame',...
'ALSub1Time','ALSub2Time','ALSolveStep','ResultStrain','ResultStress');