-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathFFT_NLIE.m
337 lines (335 loc) · 11.5 KB
/
FFT_NLIE.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
%% FFT-Non_Linear_Inductance_Extractor:
close all
clear global
clear
clc
restoredefaultpath
format short
warning off
delete log_Green_K.txt
warning on
%%
%% BEGIN USER SETTINGS
%%
%% Directory
name_dir='test';
%% Selections
plot_vectorsJ_flag = 1; %quiver plot of real and imag of J
plot_potential_flag = 1; %color plot of phi real and imag
paraview_export_flag = 1; % export to paraviw
refine.flag = 0; refine.x=1; refine.y=1; refine.z=1; % refine
Integration_flag = 'NumAn'; %'NumAn'; 'NumNum' (Integration: NumericalNumerical or AnalyticalNumerical)
%% Solver parameters
tol = 1e-6;
inner_it = 40;
outer_it = 5;
NR_tol = 1e-3;
NR_iter=50;
%%
%% END USER SETTINGS
%%
%% Add Path
dad = pwd;
cd('fun'); addpath(genpath(pwd)); cd(dad)
cd('fortran'); addpath(pwd); cd(dad)
cd('data'); cd(name_dir); load('data.mat');
fileList = dir('*.stl');
figure
hold on
xmin=[];xmax=[];ymin=[];ymax=[];zmin=[];zmax=[];ccolor=distinguishable_colors(size(fileList,1));
for ii = 1:size(fileList,1)
[stlcoords] = READ_stl(fileList(ii).name);
xco = squeeze( stlcoords(:,1,:) )';
yco = squeeze( stlcoords(:,2,:) )';
zco = squeeze( stlcoords(:,3,:) )';
[hpat] = patch(xco,yco,zco,ccolor(ii,:));
axis equal
xlabel('x');
ylabel('y');
zlabel('z');
view(3)
title('stl (original, not scaled)')
drawnow
end
cd(dad)
modelname = name_dir;
%% refine
if refine.flag
warning('refine on')
mymod=1;
for ii = 1:refine.x
[Ind,L,M,N,xyz,smeshx,smeshy,smeshz,Nmat,nVoxel] = fun_refine(Ind,xyz,smeshx,smeshy,smeshz,Nmat,L,M,N,1,mymod);
end
for ii = 1:refine.y
[Ind,L,M,N,xyz,smeshx,smeshy,smeshz,Nmat,nVoxel] = fun_refine(Ind,xyz,smeshx,smeshy,smeshz,Nmat,L,M,N,2,mymod);
end
for ii = 1:refine.z
[Ind,L,M,N,xyz,smeshx,smeshy,smeshz,Nmat,nVoxel] = fun_refine(Ind,xyz,smeshx,smeshy,smeshz,Nmat,L,M,N,3,mymod);
end
end
%% EM constants
mu = 4*pi*1e-7;
co = 299792458;
eo = 1/co^2/mu;
%% extract data information
rhoVoxel=zeros(nVoxel,1);
idxVe=[]; rhomin=Inf; ind_c=[]; val_c=[]; k=1;
idxVm=[];
for ii = 1:Nmat
Ind(ii).ind=reshape(Ind(ii).ind,length(Ind(ii).ind),1);
if strcmp(Ind(ii).tag,'air') || strcmp(Ind(ii).tag,'diel')
% nothing to do here (?)
elseif strcmp(Ind(ii).tag,'cond')
idxVe=[idxVe;Ind(ii).ind];
rhoVoxel(Ind(ii).ind,1)=Ind(ii).rho;
rhomin=min([rhomin,Ind(ii).rho]);
elseif strcmp(Ind(ii).tag,'mag')
idxVm=[idxVm;Ind(ii).ind];
BH_tab=Ind(ii).BH_tab;
elseif strcmp(Ind(ii).tag,'terminal')
idxVe=[idxVe;Ind(ii).ind];
rhoVoxel(Ind(ii).ind,1)=Ind(ii).rho;
rhomin=min([rhomin,Ind(ii).rho]);
ind_c(k)=Ind(ii).ind(1);
val_c(k)=Ind(ii).cur;
k=k+1;
end
end
idxVe=unique(idxVe);
idxVm=unique(idxVm);
%% Grid Definition
disp('----DOMAIN--------------------------------')
%%Grid resolution
disp([' Number of voxels in x direction: ', num2str(L)])
disp([' Number of voxels in y direction: ', num2str(M)])
disp([' Number of voxels in z direction: ', num2str(N)])
disp(' Resolution:')
dx = smeshx; dy = smeshy; dz = smeshz;
disp([' dx = ',num2str(dx),' m']); disp([' dy = ',num2str(dy),' m']); disp([' dz = ',num2str(dz),' m'])
d = [dx dy dz];
Kt = nVoxel; %total number of voxels
K = length(idxVe); %number of non-empty voxels
Km = length(idxVm); %number of non-empty voxels
%% Set Material Properties
rho_eV=reshape(rhoVoxel,L,M,N); %
clear rhoVoxel
%%
disp([' Total number of voxels: ', num2str(Kt)])
disp([' Number of non-empty voxels: ', num2str(K)])
disp(' ')
%% Incidence Matix A
disp('----COMPUTING INCIDENCE--------------------------------')
mytic=tic;
[Ae,Aee,idxF,idxFx,idxFy,idxFz,Ae1x,Ae1y,Ae1z] = ...
incidence_matrix3(Kt,[L M N],idxVe);
ind=zeros(L*M*N,1);
ind(ind_c)=val_c;
ind=ind(idxVe);
ind=setdiff(1:length(idxVe),find(ind));
Aee(max(ind),:)=[];
idxVe(max(ind))=[];
K = length(idxVe); %number of non-empty voxels
[~,Amm,idxFm,idxFxm,idxFym,idxFzm,~,~,~] = ...
incidence_matrix3(Kt,[L M N],idxVm);
disp([' Number of DoFs (cond): ', num2str(size(Aee,1)+size(Aee,2))])
disp([' Number of DoFs (mag): ', num2str(size(Amm,1)+size(Amm,2))])
disp([' Time for computing incidence ::: ' ,num2str(toc(mytic))]);
disp(' ')
%% Forcing Term
Vx=zeros(L*M*N,1);
Vy=zeros(L*M*N,1);
Vz=zeros(L*M*N,1);
%% Matrices Z_real and Z_imag
rho_eF=0.5*(abs(Ae(:,:)).'*rho_eV(:)); clear rho_eV
z_realF=rho_eF;
indFneq=setdiff([1:3*Kt].',idxF);
z_realF(indFneq,:)=0;
z_realx=zeros(L,M,N);
z_realx(idxFx)=z_realF(idxFx);
z_realy=zeros(L,M,N);
z_realy(idxFy)=z_realF(Kt+idxFy);
z_realz=zeros(L,M,N);
z_realz(idxFz)=z_realF(2*Kt+idxFz);
%
%% Compute Green Tensor
disp('----COMPUTING GREEN TENSOR--------------------------------')
mytic_G=tic;
[Gmn] = computeGREEN(d,L,M,N,Integration_flag);
[Kmn] = computeGREEN_GK(d,L,M,N,Integration_flag,'mex');
[opCirculantL_all,st_sparse_preconL] = computeCIRCULANT(Gmn,d,'L');
opCirculantL_all = (mu)*opCirculantL_all;
disp([' Time for getting Green tensor ::: ' ,num2str(toc(mytic_G))]);
disp(' ')
%% Compute Circulant Tensors
disp('----COMPUTING CIRCULANT TENSOR--------------------------------')
disp(' Circulant Tensors related to P,L matrices')
mytic_cir=tic;
[opCirculantP_all,st_sparse_preconP] = computeCIRCULANT(Gmn,d,'P');
opCirculantP_all=opCirculantP_all/mu;
st_sparse_preconP=st_sparse_preconP/mu;
[opCirculantK_all] = computeCIRCULANT_K(Kmn);
%%Add constants to Circulants
disp([' Time for getting circulant tensors ::: ' ,num2str(toc(mytic_cir))])
clear Gmn %Green tensor is not used anymore
disp(' ')
%% Generating RHS vector
num_node = size(Aee,1); %all potential nodes in non-empty voxels
num_curr = size(Aee,2); %all currens in non-empty voxels
num_nodem = size(Amm,1); %all potential nodes in non-empty voxels
num_currm = size(Amm,2); %all currens in non-empty voxels
%%Define RHS: (injected currents)
iinj=zeros(L*M*N,1);
iinj(ind_c)=val_c;
iinj=iinj(idxVe);
rhs_vect = [Vx(idxFx);Vy(idxFy);Vz(idxFz);-iinj];
clear Vx Vy Vz
%% Computing Preconditioner
disp('----COMPUTING PRECONDITIONER--------------------------------')
mytic_prec=tic;
[Y_inv,LL,UU,PP,QQ,RR] = preparePREC_NEW(d,z_realF,idxFx,...
idxFy,idxFz,Aee,Kt);
fPMV = @(JOut_full_in)multiplyPREC_NEW(JOut_full_in,Aee,Y_inv,LL,UU,PP,QQ,RR);
disp([' Time for computing preconditioner ::: ' ,num2str(toc(mytic_prec))]);
disp(' ')
%% Solution of Linear System
disp('----SOLVING LINEAR SYSTEM-------------------------------')
zloc=[(dx/(dy*dz))*z_realx(idxFx);...
(dy/(dx*dz))*z_realy(idxFy);...
(dz/(dx*dy))*z_realz(idxFz)];
fMVM = @(J) multiplyMATVECT(J,zloc,Aee);
mytic_solver=tic;
[vsol] = pgmres_mod(@(J)fMVM(J),rhs_vect, inner_it, tol, outer_it, @(JOut_full_in)fPMV(JOut_full_in) );
disp([' Time for solving system with gmres ::: ' ,num2str(toc(mytic_solver))]);
disp(' ')
%% extract solution
Jout = zeros(L,M,N,3);
Jout(idxF) = vsol(1:num_curr) ; % return to global variables
%%
[A,a_dof] = post_processing_A(Jout,xyz,d,opCirculantL_all,Ae1x,Ae1y,Ae1z);
Wm=0.5*0.5*sum((vsol(1:num_curr).*conj(a_dof(idxF))));
L1=2*2*Wm/(abs(val_c(1))^2);
disp('-------------------------------------------------------------------')
disp(' ')
disp([' Coil Inductance ' ,num2str(L1),' H [without magnetic core]']);
disp(' ')
disp('-------------------------------------------------------------------')
%%
[H] = post_processing_H(Jout,xyz,d,Ae1x,Ae1y,Ae1z,opCirculantK_all);
Hx = reshape(H(:,1),L,M,N);
Hy = reshape(H(:,2),L,M,N);
Hz = reshape(H(:,3),L,M,N);
Gram = dx*dy*dz; %volume of cubic element
Vmx = (Gram.*Hx)./(dy*dz);
Vmy = (Gram.*Hy)./(dx*dz);
Vmz = (Gram.*Hz)./(dx*dy);
clear Hx Hy Hz
rhs_vectm = [Vmx(idxFxm);Vmy(idxFym);Vmz(idxFzm);zeros(num_nodem,1)];
%%
[fun_J_to_coeff,fun_J2_to_dcoeff] = fun_from_BH_to_funs3(BH_tab);
fun2.F=@(x) fun_Ax_minus_b ...
(x,d,idxVm,idxFm,idxFxm,idxFym,...
idxFzm,Kt,Ae,Ae1x,Ae1y,Ae1z,L,M,N,...
fun_J_to_coeff,...
opCirculantP_all,Amm,Amm,idxVm,rhs_vectm);
fun2.solvJ=@(x,FF) fun_resolve_Jac(x,d,idxVm,idxFm,idxFxm,idxFym,...
idxFzm,Kt,Ae,Ae1x,Ae1y,Ae1z,L,M,N,...
fun_J_to_coeff,...
opCirculantP_all,Amm,Amm,idxVm,FF,...
fun_J2_to_dcoeff,Km,st_sparse_preconP,'schur_invert','gmres'); % funzione che risolve il problema con jacobiano Jac(x)\b
fun2.JTFF_T = @(x,FF) fun_JacTdotF_T(x,d,idxVm,idxFm,idxFxm,idxFym,...
idxFzm,Kt,Ae,Ae1x,Ae1y,Ae1z,L,M,N,...
fun_J_to_coeff,...
opCirculantP_all,Amm,Amm,idxVm,FF,...
fun_J2_to_dcoeff,Km);
global NLcounter
NLcounter=1;
options = optimset('TolFun',NR_tol,'MaxIter',NR_iter); % set TolX
mytic_sol_NR=tic;
[vsolm, resnorm, f, exitflag, ...
output, jacob,resnormstore,xstore] ...
= newtonraphson_mod2(fun2, zeros(num_currm+num_nodem,1), options);
disp([' Total time for solving NON LINEAR SYSTEM ::: ' ,num2str(toc(mytic_sol_NR))])
%% extract solution
Joutm = zeros(L,M,N,3);
Joutm(idxFm) = vsolm(1:num_currm) ; % return to global variables
%%
[~,a_dof2] = post_processing_H(Joutm,xyz,d,Ae1x,Ae1y,Ae1z,opCirculantK_all);
Wm=0.5*0.5*sum(((vsol(1:num_curr)).*conj(a_dof(idxF)+a_dof2(idxF))));
L2=2*2*Wm/(abs(val_c(1))^2);
disp(' ')
disp('-------------------------------------------------------------------')
disp(' ')
disp([' Coil Inductance ' ,num2str(L2),' H [with magnetic core]']);
disp(' ')
disp('-------------------------------------------------------------------')
disp(' ')
%%
%% POST PROCESSING
%%
%% Post Processing J
disp('----POST PROCESSING J------------------------------')
mytic_prec=tic;
[J,XYZ] = fun_my_postRT2(Jout,Kt,Ae1x,Ae1y,Ae1z,xyz,L,M,N,d);
[MM,XYZ] = fun_my_postRT2(Joutm,Kt,Ae1x,Ae1y,Ae1z,xyz,L,M,N,d);
potval=zeros(Kt,1);
potval(idxVe)=vsol(num_curr+1:end);
disp([' Total time for post processing J and mu0M ::: ' ,num2str(toc(mytic_prec))]);
disp(' ')
%% Plot Vectors J
if plot_vectorsJ_flag
jjR = (J);
figure
normJR=sqrt(jjR(:,1).^2+jjR(:,2).^2+jjR(:,3).^2);
quiver3_c_scal(XYZ(:,1),XYZ(:,2),XYZ(:,3),jjR(:,1),jjR(:,2),jjR(:,3),...
normJR,4);
axis equal
c1=colorbar;
caxis([min(normJR) max(normJR)]);
xlabel('x')
ylabel('y')
zlabel('z')
title('Current Density Vector ')
c1.Location = 'southoutside';
xlim([min(XYZ(:,1))-dx max(XYZ(:,1))+dx])
ylim([min(XYZ(:,2))-dy max(XYZ(:,2))+dy])
zlim([min(XYZ(:,3))-dz max(XYZ(:,3))+dz])
%% Plot Vectors H
hhR = H;%(M);
figure
normHR=sqrt(hhR(:,1).^2+hhR(:,2).^2+hhR(:,3).^2);
quiver3_c_scal(XYZ(:,1),XYZ(:,2),XYZ(:,3),hhR(:,1),hhR(:,2),hhR(:,3),...
normHR,4);
axis equal
c1=colorbar;
caxis([min(normHR) max(normHR)]);
xlabel('x')
ylabel('y')
zlabel('z')
title('H produced by coil only')
c1.Location = 'southoutside';
xlim([min(XYZ(:,1))-dx max(XYZ(:,1))+dx])
ylim([min(XYZ(:,2))-dy max(XYZ(:,2))+dy])
zlim([min(XYZ(:,3))-dz max(XYZ(:,3))+dz])
%
%% Plot Vectors H
mmR = MM;
figure
normMR=sqrt(mmR(:,1).^2+mmR(:,2).^2+mmR(:,3).^2);
quiver3_c_scal(XYZ(:,1),XYZ(:,2),XYZ(:,3),mmR(:,1),mmR(:,2),mmR(:,3),...
normMR,4);
axis equal
c1=colorbar;
caxis([min(normMR) max(normMR)]);
xlabel('x')
ylabel('y')
zlabel('z')
title('\mu_0M Vector')
c1.Location = 'southoutside';
xlim([min(XYZ(:,1))-dx max(XYZ(:,1))+dx])
ylim([min(XYZ(:,2))-dy max(XYZ(:,2))+dy])
zlim([min(XYZ(:,3))-dz max(XYZ(:,3))+dz])
%
end
warning off
delete log_Green_K.txt
warning on