-
Notifications
You must be signed in to change notification settings - Fork 587
/
pympb.cpp
2677 lines (2225 loc) · 94.4 KB
/
pympb.cpp
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
#include <algorithm>
#include <complex>
#include <cstdlib>
#include <cstddef>
#include <iostream>
#include "config.h"
#include "pympb.hpp"
#include "meep/mympi.hpp"
// If the MPB lib is not new enough to have the mpb_verbosity global then make
// one here to give the swig wrapper something to link with.
#if MPB_VERSION_MAJOR > 1 || (MPB_VERSION_MAJOR == 1 && MPB_VERSION_MINOR >= 11)
// do nothing, libmpb should have the mpb_verbosity symbol
#else
extern "C" int mpb_verbosity = 2;
#endif
// xyz_loop.h
#ifndef HAVE_MPI
#define LOOP_XYZ(md) \
{ \
int n1 = md->nx, n2 = md->ny, n3 = md->nz, i1, i2, i3; \
for (i1 = 0; i1 < n1; ++i1) \
for (i2 = 0; i2 < n2; ++i2) \
for (i3 = 0; i3 < n3; ++i3) { \
int xyz_index = ((i1 * n2 + i2) * n3 + i3);
#else /* HAVE_MPI */
/* first two dimensions are transposed in MPI output: */
#define LOOP_XYZ(md) \
{ \
int n1 = md->nx, n2 = md->ny, n3 = md->nz, i1, i2_, i3; \
int local_n2 = md->local_ny, local_y_start = md->local_y_start; \
for (i2_ = 0; i2_ < local_n2; ++i2_) \
for (i1 = 0; i1 < n1; ++i1) \
for (i3 = 0; i3 < n3; ++i3) { \
int i2 = i2_ + local_y_start; \
int xyz_index = ((i2_ * n1 + i1) * n3 + i3);
#endif /* HAVE_MPI */
typedef mpb_real real; // needed for the CASSIGN macros below
// TODO: Support MPI
#define mpi_allreduce(sb, rb, n, ctype, t, op, comm) \
{ \
CHECK((sb) != (rb), "MPI_Allreduce doesn't work for sendbuf == recvbuf"); \
memcpy((rb), (sb), (n) * sizeof(ctype)); \
}
/* "in-place" Allreduce wrapper for reducing a single value */
#define mpi_allreduce_1(b, ctype, t, op, comm) \
{ \
ctype bbbb = *(b); \
mpi_allreduce(&bbbb, (b), 1, ctype, t, op, comm); \
}
#ifdef CHECK_DISABLE
#define CHECK(cond, s) // Do nothing
#else
#define CHECK(cond, s) \
if (!(cond)) { meep::abort(s "\n"); }
#endif
namespace py_mpb {
// TODO: Placeholder
int mpb_comm;
const double inf = 1.0e20;
// This is the function passed to `set_maxwell_dielectric`
static void dielectric_function(symmetric_matrix *eps, symmetric_matrix *eps_inv,
const mpb_real r[3], void *epsilon_data) {
mode_solver *ms = static_cast<mode_solver *>(epsilon_data);
meep_geom::material_type mat;
vector3 p;
// p needs to be in the lattice *unit* vector basis, while r is in the lattice
// vector basis. Also, shift origin to the center of the grid.
p.x = (r[0] - 0.5) * geometry_lattice.size.x;
p.y = (r[1] - 0.5) * geometry_lattice.size.y;
p.z = (r[2] - 0.5) * geometry_lattice.size.z;
// p = shift_to_unit_cell(p);
ms->get_material_pt(mat, p);
ms->material_epsmu(mat, eps, eps_inv, eps);
}
static int mean_epsilon_func(symmetric_matrix *meps, symmetric_matrix *meps_inv, mpb_real n[3],
mpb_real d1, mpb_real d2, mpb_real d3, mpb_real tol,
const mpb_real r[3], void *edata) {
mode_solver *ms = static_cast<mode_solver *>(edata);
return ms->mean_epsilon(meps, meps_inv, n, d1, d2, d3, tol, r);
}
/****** utils ******/
/* a couple of utilities to convert libctl data types to the data
types of the eigensolver & maxwell routines: */
void vector3_to_arr(mpb_real arr[3], vector3 v) {
arr[0] = v.x;
arr[1] = v.y;
arr[2] = v.z;
}
void matrix3x3_to_arr(mpb_real arr[3][3], matrix3x3 m) {
vector3_to_arr(arr[0], m.c0);
vector3_to_arr(arr[1], m.c1);
vector3_to_arr(arr[2], m.c2);
}
cnumber cscalar2cnumber(scalar_complex cs) { return make_cnumber(CSCALAR_RE(cs), CSCALAR_IM(cs)); }
// Return a string describing the current parity, used for frequency and filename
// prefixes
const char *parity_string(maxwell_data *d) {
static char s[128];
strcpy(s, "");
if (d->parity & EVEN_Z_PARITY) { strcat(s, (d->nz == 1) ? "te" : "zeven"); }
else if (d->parity & ODD_Z_PARITY) {
strcat(s, (d->nz == 1) ? "tm" : "zodd");
}
if (d->parity & EVEN_Y_PARITY) { strcat(s, "yeven"); }
else if (d->parity & ODD_Y_PARITY) {
strcat(s, "yodd");
}
return s;
}
/* Extract the mean epsilon from the effective inverse dielectric tensor,
which contains two eigenvalues that correspond to the mean epsilon,
and one which corresponds to the harmonic mean. */
mpb_real mean_medium_from_matrix(const symmetric_matrix *eps_inv) {
mpb_real eps_eigs[3];
maxwell_sym_matrix_eigs(eps_eigs, eps_inv);
/* the harmonic mean should be the largest eigenvalue (smallest
epsilon), so we'll ignore it and average the other two: */
return 2.0 / (eps_eigs[0] + eps_eigs[1]);
}
/* When we are solving for a few bands at a time, we solve for the
upper bands by "deflation"--by continually orthogonalizing them
against the already-computed lower bands. (This constraint
commutes with the eigen-operator, of course, so all is well.) */
typedef struct {
evectmatrix Y; /* the vectors to orthogonalize against; Y must
itself be normalized (Yt B Y = 1) */
evectmatrix BY; /* B * Y */
int p; /* the number of columns of Y to orthogonalize against */
scalar *S; /* a matrix for storing the dot products; should have
at least p * X.p elements (see below for X) */
scalar *S2; /* a scratch matrix the same size as S */
} deflation_data;
extern "C" {
void blasglue_gemm(char transa, char transb, int m, int n, int k, mpb_real a, scalar *A, int fdA,
scalar *B, int fdB, mpb_real b, scalar *C, int fdC);
}
static void deflation_constraint(evectmatrix X, void *data) {
deflation_data *d = (deflation_data *)data;
CHECK(X.n == d->BY.n && d->BY.p >= d->p && d->Y.p >= d->p, "invalid dimensions");
/* compute (1 - Y (BY)t) X = (1 - Y Yt B) X
= projection of X so that Yt B X = 0 */
/* (Sigh...call the BLAS functions directly since we are not
using all the columns of BY...evectmatrix is not set up for
this case.) */
/* compute S = Xt BY (i.e. all the dot products): */
blasglue_gemm('C', 'N', X.p, d->p, X.n, 1.0, X.data, X.p, d->BY.data, d->BY.p, 0.0, d->S2, d->p);
// TODO
// #if HAVE_MPI
// MPI_Allreduce(d->S2, d->S, d->p * X.p * SCALAR_NUMVALS, SCALAR_MPI_TYPE,
// MPI_SUM, mpb_comm);
// #else
memcpy(d->S, d->S2, sizeof(mpb_real) * d->p * X.p * SCALAR_NUMVALS);
// #endif
/* compute X = X - Y*St = (1 - BY Yt B) X */
blasglue_gemm('N', 'C', X.n, X.p, d->p, -1.0, d->Y.data, d->Y.p, d->S, d->p, 1.0, X.data, X.p);
}
/******* mode_solver *******/
mode_solver::mode_solver(int num_bands, double resolution[3], lattice lat, double tolerance,
int mesh_size, meep_geom::material_data *_default_material,
bool deterministic, double target_freq, int dims, bool verbose,
bool periodicity, double flops, bool negative_epsilon_ok,
std::string epsilon_input_file, std::string mu_input_file, bool force_mu,
bool use_simple_preconditioner, vector3 grid_size, int eigensolver_nwork,
int eigensolver_block_size)
: num_bands(num_bands), resolution{resolution[0], resolution[1], resolution[2]},
target_freq(target_freq), tolerance(tolerance), mesh_size(mesh_size),
negative_epsilon_ok(negative_epsilon_ok), epsilon_input_file(epsilon_input_file),
mu_input_file(mu_input_file), force_mu(force_mu),
use_simple_preconditioner(use_simple_preconditioner), grid_size(grid_size), nwork_alloc(0),
eigensolver_nwork(eigensolver_nwork), eigensolver_block_size(eigensolver_block_size),
last_parity(-2), iterations(0), eigensolver_flops(flops), geometry_list{},
geometry_tree(NULL), vol(0), R{}, G{}, mdata(NULL), mtdata(NULL),
curfield_band(0), H{}, Hblock{}, muinvH{}, W{}, freqs(num_bands), verbose(verbose),
deterministic(deterministic), kpoint_index(0), curfield(NULL), curfield_type('-'), eps(true) {
// See geom-ctl-io-defaults.c in libctl
geometry_lattice = lat;
dimensions = dims;
ensure_periodicity = periodicity;
#ifndef WITH_HERMITIAN_EPSILON
meep_geom::medium_struct *m;
if (meep_geom::is_medium(_default_material, &m)) { m->check_offdiag_im_zero_or_abort(); }
#else
(void)_default_material;
#endif
}
mode_solver::~mode_solver() {
destroy_maxwell_data(mdata);
destroy_maxwell_target_data(mtdata);
destroy_geom_box_tree(geometry_tree);
clear_geometry_list();
destroy_evectmatrix(H);
for (int i = 0; i < nwork_alloc; ++i) {
destroy_evectmatrix(W[i]);
}
if (Hblock.data != H.data) { destroy_evectmatrix(Hblock); }
if (muinvH.data != H.data) { destroy_evectmatrix(muinvH); }
}
int mode_solver::mean_epsilon(symmetric_matrix *meps, symmetric_matrix *meps_inv, mpb_real n[3],
mpb_real d1, mpb_real d2, mpb_real d3, mpb_real tol,
const mpb_real r[3]) {
const geometric_object *o1 = 0;
const geometric_object *o2 = 0;
geom_box pixel;
double fill;
meep_geom::material_type mat1;
meep_geom::material_type mat2;
int id1 = -1;
int id2 = -1;
const int num_neighbors[3] = {3, 5, 9};
const int neighbors[3][9][3] = {{{0, 0, 0},
{-1, 0, 0},
{1, 0, 0},
{0, 0, 0},
{0, 0, 0},
{0, 0, 0},
{0, 0, 0},
{0, 0, 0},
{0, 0, 0}},
{{0, 0, 0},
{-1, -1, 0},
{1, 1, 0},
{-1, 1, 0},
{1, -1, 0},
{0, 0, 0},
{0, 0, 0},
{0, 0, 0},
{0, 0, 0}},
{{0, 0, 0},
{1, 1, 1},
{1, 1, -1},
{1, -1, 1},
{1, -1, -1},
{-1, 1, 1},
{-1, 1, -1},
{-1, -1, 1},
{-1, -1, -1}}};
/* p needs to be in the lattice *unit* vector basis, while r is
in the lattice vector basis. Also, shift origin to the center
of the grid. */
vector3 p = {(r[0] - 0.5) * geometry_lattice.size.x, (r[1] - 0.5) * geometry_lattice.size.y,
(r[2] - 0.5) * geometry_lattice.size.z};
d1 *= geometry_lattice.size.x * 0.5;
d2 *= geometry_lattice.size.y * 0.5;
d3 *= geometry_lattice.size.z * 0.5;
vector3 shiftby1;
vector3 shiftby2;
vector3 normal;
for (int i = 0; i < num_neighbors[dimensions - 1]; ++i) {
const geometric_object *o;
vector3 q, z, shiftby;
int id;
q.x = p.x + neighbors[dimensions - 1][i][0] * d1;
q.y = p.y + neighbors[dimensions - 1][i][1] * d2;
q.z = p.z + neighbors[dimensions - 1][i][2] * d3;
geometry_lattice.size.x = geometry_lattice.size.x == 0 ? 1e-20 : geometry_lattice.size.x;
geometry_lattice.size.y = geometry_lattice.size.y == 0 ? 1e-20 : geometry_lattice.size.y;
geometry_lattice.size.z = geometry_lattice.size.z == 0 ? 1e-20 : geometry_lattice.size.z;
z = shift_to_unit_cell(q);
geometry_lattice.size.x = geometry_lattice.size.x == 1e-20 ? 0 : geometry_lattice.size.x;
geometry_lattice.size.y = geometry_lattice.size.y == 1e-20 ? 0 : geometry_lattice.size.y;
geometry_lattice.size.z = geometry_lattice.size.z == 1e-20 ? 0 : geometry_lattice.size.z;
o = object_of_point_in_tree(z, geometry_tree, &shiftby, &id);
shiftby = vector3_plus(shiftby, vector3_minus(q, z));
if ((id == id1 && vector3_equal(shiftby, shiftby1)) ||
(id == id2 && vector3_equal(shiftby, shiftby2))) {
continue;
}
meep_geom::material_type mat = (meep_geom::material_type)default_material;
if (o) {
meep_geom::material_data *md = (meep_geom::material_data *)o->material;
if (md->which_subclass != meep_geom::material_data::MATERIAL_FILE) { mat = md; }
}
if (id1 == -1) {
o1 = o;
shiftby1 = shiftby;
id1 = id;
mat1 = mat;
}
else if (id2 == -1 || ((id >= id1 && id >= id2) &&
(id1 == id2 || meep_geom::material_type_equal(mat1, mat2)))) {
o2 = o;
shiftby2 = shiftby;
id2 = id;
mat2 = mat;
}
else if (!(id1 < id2 && (id1 == id || meep_geom::material_type_equal(mat1, mat))) &&
!(id2 < id1 && (id2 == id || meep_geom::material_type_equal(mat2, mat)))) {
return 0; /* too many nearby objects for analysis */
}
}
CHECK(id1 > -1, "bug in object_of_point_in_tree?");
if (id2 == -1) { /* only one nearby object/material */
id2 = id1;
o2 = o1;
mat2 = mat1;
shiftby2 = shiftby1;
}
bool o1_is_var = o1 && meep_geom::is_variable(o1->material);
bool o2_is_var = o2 && meep_geom::is_variable(o2->material);
bool default_is_var = meep_geom::is_variable(default_material);
if (o1_is_var || o2_is_var || (default_is_var && (!o1 || !o2))) {
return 0; /* arbitrary material functions are non-analyzable */
}
material_epsmu(mat1, meps, meps_inv, eps);
/* check for trivial case of only one object/material */
if (id1 == id2 || meep_geom::material_type_equal(mat1, mat2)) {
n[0] = n[1] = n[2] = 0;
return 1;
}
if (id1 > id2) { normal = normal_to_fixed_object(vector3_minus(p, shiftby1), *o1); }
else {
normal = normal_to_fixed_object(vector3_minus(p, shiftby2), *o2);
}
n[0] = normal.x / (geometry_lattice.size.x == 0 ? 1e-20 : geometry_lattice.size.x);
n[1] = normal.y / (geometry_lattice.size.y == 0 ? 1e-20 : geometry_lattice.size.y);
n[2] = normal.z / (geometry_lattice.size.z == 0 ? 1e-20 : geometry_lattice.size.z);
pixel.low.x = p.x - d1;
pixel.high.x = p.x + d1;
pixel.low.y = p.y - d2;
pixel.high.y = p.y + d2;
pixel.low.z = p.z - d3;
pixel.high.z = p.z + d3;
tol = tol > 0.01 ? 0.01 : tol;
if (id1 > id2) {
pixel.low = vector3_minus(pixel.low, shiftby1);
pixel.high = vector3_minus(pixel.high, shiftby1);
fill = box_overlap_with_object(pixel, *o1, tol, 100 / tol);
}
else {
pixel.low = vector3_minus(pixel.low, shiftby2);
pixel.high = vector3_minus(pixel.high, shiftby2);
fill = 1 - box_overlap_with_object(pixel, *o2, tol, 100 / tol);
}
{
symmetric_matrix eps2, epsinv2;
symmetric_matrix eps1, delta;
double Rot[3][3], norm, n0, n1, n2;
material_epsmu(mat2, &eps2, &epsinv2, eps);
eps1 = *meps;
/* make Cartesian orthonormal frame relative to interface */
n0 = R[0][0] * n[0] + R[1][0] * n[1] + R[2][0] * n[2];
n1 = R[0][1] * n[0] + R[1][1] * n[1] + R[2][1] * n[2];
n2 = R[0][2] * n[0] + R[1][2] * n[1] + R[2][2] * n[2];
norm = sqrt(n0 * n0 + n1 * n1 + n2 * n2);
if (norm == 0.0) { return 0; }
norm = 1.0 / norm;
Rot[0][0] = n0 = n0 * norm;
Rot[1][0] = n1 = n1 * norm;
Rot[2][0] = n2 = n2 * norm;
if (fabs(n0) > 1e-2 || fabs(n1) > 1e-2) { /* (z x n) */
Rot[0][2] = n1;
Rot[1][2] = -n0;
Rot[2][2] = 0;
}
else { /* n is ~ parallel to z direction, use (x x n) instead */
Rot[0][2] = 0;
Rot[1][2] = -n2;
Rot[2][2] = n1;
}
{ /* normalize second column */
double s = Rot[0][2] * Rot[0][2] + Rot[1][2] * Rot[1][2] + Rot[2][2] * Rot[2][2];
s = 1.0 / sqrt(s);
Rot[0][2] *= s;
Rot[1][2] *= s;
Rot[2][2] *= s;
}
/* 1st column is 2nd column x 0th column */
Rot[0][1] = Rot[1][2] * Rot[2][0] - Rot[2][2] * Rot[1][0];
Rot[1][1] = Rot[2][2] * Rot[0][0] - Rot[0][2] * Rot[2][0];
Rot[2][1] = Rot[0][2] * Rot[1][0] - Rot[1][2] * Rot[0][0];
/* rotate epsilon tensors to surface parallel/perpendicular axes */
maxwell_sym_matrix_rotate(&eps1, &eps1, Rot);
maxwell_sym_matrix_rotate(&eps2, &eps2, Rot);
#define AVG (fill * (EXPR(eps1)) + (1 - fill) * (EXPR(eps2)))
#define EXPR(eps) (-1 / eps.m00)
delta.m00 = AVG;
#undef EXPR
#define EXPR(eps) (eps.m11 - ESCALAR_NORMSQR(eps.m01) / eps.m00)
delta.m11 = AVG;
#undef EXPR
#define EXPR(eps) (eps.m22 - ESCALAR_NORMSQR(eps.m02) / eps.m00)
delta.m22 = AVG;
#undef EXPR
#define EXPR(eps) (ESCALAR_RE(eps.m01) / eps.m00)
ESCALAR_RE(delta.m01) = AVG;
#undef EXPR
#define EXPR(eps) (ESCALAR_RE(eps.m02) / eps.m00)
ESCALAR_RE(delta.m02) = AVG;
#undef EXPR
#define EXPR(eps) (ESCALAR_RE(eps.m12) - ESCALAR_MULT_CONJ_RE(eps.m02, eps.m01) / eps.m00)
ESCALAR_RE(delta.m12) = AVG;
#undef EXPR
#ifdef WITH_HERMITIAN_EPSILON
#define EXPR(eps) (ESCALAR_IM(eps.m01) / eps.m00)
ESCALAR_IM(delta.m01) = AVG;
#undef EXPR
#define EXPR(eps) (ESCALAR_IM(eps.m02) / eps.m00)
ESCALAR_IM(delta.m02) = AVG;
#undef EXPR
#define EXPR(eps) (ESCALAR_IM(eps.m12) - ESCALAR_MULT_CONJ_IM(eps.m02, eps.m01) / eps.m00)
ESCALAR_IM(delta.m12) = AVG;
#undef EXPR
#endif /* WITH_HERMITIAN_EPSILON */
meps->m00 = -1 / delta.m00;
meps->m11 = delta.m11 - ESCALAR_NORMSQR(delta.m01) / delta.m00;
meps->m22 = delta.m22 - ESCALAR_NORMSQR(delta.m02) / delta.m00;
ASSIGN_ESCALAR(meps->m01, -ESCALAR_RE(delta.m01) / delta.m00,
-ESCALAR_IM(delta.m01) / delta.m00);
ASSIGN_ESCALAR(meps->m02, -ESCALAR_RE(delta.m02) / delta.m00,
-ESCALAR_IM(delta.m02) / delta.m00);
ASSIGN_ESCALAR(meps->m12,
ESCALAR_RE(delta.m12) - ESCALAR_MULT_CONJ_RE(delta.m02, delta.m01) / delta.m00,
ESCALAR_IM(delta.m12) - ESCALAR_MULT_CONJ_IM(delta.m02, delta.m01) / delta.m00);
#define SWAP(a, b) \
{ \
double xxx = a; \
a = b; \
b = xxx; \
}
/* invert rotation matrix = transpose */
SWAP(Rot[0][1], Rot[1][0]);
SWAP(Rot[0][2], Rot[2][0]);
SWAP(Rot[2][1], Rot[1][2]);
maxwell_sym_matrix_rotate(meps, meps, Rot); /* rotate back */
#undef SWAP
#ifdef DEBUG
CHECK(negative_epsilon_ok || maxwell_sym_matrix_positive_definite(meps),
"negative mean epsilon from Kottke algorithm");
#endif
}
return 1;
}
void mode_solver::material_epsmu(meep_geom::material_type material, symmetric_matrix *epsmu,
symmetric_matrix *epsmu_inv, bool eps) {
meep_geom::material_data *md = material;
#ifndef WITH_HERMITIAN_EPSILON
if (md->which_subclass == meep_geom::material_data::MATERIAL_USER ||
md->which_subclass == meep_geom::material_data::MATERIAL_FILE) {
md->medium.check_offdiag_im_zero_or_abort();
}
#endif
if (eps) {
switch (md->which_subclass) {
case meep_geom::material_data::MEDIUM:
case meep_geom::material_data::MATERIAL_FILE:
case meep_geom::material_data::MATERIAL_USER:
epsmu->m00 = md->medium.epsilon_diag.x;
epsmu->m11 = md->medium.epsilon_diag.y;
epsmu->m22 = md->medium.epsilon_diag.z;
#ifdef WITH_HERMITIAN_EPSILON
epsmu->m01.re = md->medium.epsilon_offdiag.x.re;
epsmu->m01.im = md->medium.epsilon_offdiag.x.im;
epsmu->m02.re = md->medium.epsilon_offdiag.y.re;
epsmu->m02.im = md->medium.epsilon_offdiag.y.im;
epsmu->m12.re = md->medium.epsilon_offdiag.z.re;
epsmu->m12.im = md->medium.epsilon_offdiag.z.im;
#else
epsmu->m01 = md->medium.epsilon_offdiag.x.re;
epsmu->m02 = md->medium.epsilon_offdiag.y.re;
epsmu->m12 = md->medium.epsilon_offdiag.z.re;
#endif
maxwell_sym_matrix_invert(epsmu_inv, epsmu);
break;
case meep_geom::material_data::PERFECT_METAL:
epsmu->m00 = -inf;
epsmu->m11 = -inf;
epsmu->m22 = -inf;
#ifdef WITH_HERMITIAN_EPSILON
epsmu->m01.re = 0.0;
epsmu->m01.im = 0.0;
epsmu->m02.re = 0.0;
epsmu->m02.im = 0.0;
epsmu->m12.re = 0.0;
epsmu->m12.im = 0.0;
epsmu_inv->m01.re = 0.0;
epsmu_inv->m01.im = 0.0;
epsmu_inv->m02.re = 0.0;
epsmu_inv->m02.im = 0.0;
epsmu_inv->m12.re = 0.0;
epsmu_inv->m12.im = 0.0;
#else
epsmu->m01 = 0.0;
epsmu->m02 = 0.0;
epsmu->m12 = 0.0;
epsmu_inv->m01 = 0.0;
epsmu_inv->m02 = 0.0;
epsmu_inv->m12 = 0.0;
#endif
epsmu_inv->m00 = -0.0;
epsmu_inv->m11 = -0.0;
epsmu_inv->m22 = -0.0;
break;
default: meep::abort("Unknown material type");
}
}
else {
switch (md->which_subclass) {
case meep_geom::material_data::MEDIUM:
case meep_geom::material_data::MATERIAL_FILE:
case meep_geom::material_data::MATERIAL_USER:
epsmu->m00 = md->medium.mu_diag.x;
epsmu->m11 = md->medium.mu_diag.y;
epsmu->m22 = md->medium.mu_diag.z;
#ifdef WITH_HERMITIAN_EPSILON
epsmu->m01.re = md->medium.mu_offdiag.x.re;
epsmu->m01.im = md->medium.mu_offdiag.x.im;
epsmu->m02.re = md->medium.mu_offdiag.y.re;
epsmu->m02.im = md->medium.mu_offdiag.y.im;
epsmu->m12.re = md->medium.mu_offdiag.z.re;
epsmu->m12.im = md->medium.mu_offdiag.z.im;
#else
epsmu->m01 = md->medium.mu_offdiag.x.re;
epsmu->m02 = md->medium.mu_offdiag.y.re;
epsmu->m12 = md->medium.mu_offdiag.z.re;
#endif
maxwell_sym_matrix_invert(epsmu_inv, epsmu);
break;
case meep_geom::material_data::PERFECT_METAL:
epsmu->m00 = 1.0;
epsmu->m11 = 1.0;
epsmu->m22 = 1.0;
epsmu_inv->m00 = 1.0;
epsmu_inv->m11 = 1.0;
epsmu_inv->m22 = 1.0;
#ifdef WITH_HERMITIAN_EPSILON
epsmu->m01.re = 0.0;
epsmu->m01.im = 0.0;
epsmu->m02.re = 0.0;
epsmu->m02.im = 0.0;
epsmu->m12.re = 0.0;
epsmu->m12.im = 0.0;
epsmu_inv->m01.re = 0.0;
epsmu_inv->m01.im = 0.0;
epsmu_inv->m02.re = 0.0;
epsmu_inv->m02.im = 0.0;
epsmu_inv->m12.re = 0.0;
epsmu_inv->m12.im = 0.0;
#else
epsmu->m01 = 0.0;
epsmu->m02 = 0.0;
epsmu->m12 = 0.0;
epsmu_inv->m01 = 0.0;
epsmu_inv->m02 = 0.0;
epsmu_inv->m12 = 0.0;
#endif
break;
default: meep::abort("unknown material type");
}
}
}
void mode_solver::get_material_pt(meep_geom::material_type &material, vector3 p) {
boolean inobject;
material = (meep_geom::material_type)material_of_unshifted_point_in_tree_inobject(
p, geometry_tree, &inobject);
meep_geom::material_data *md = material;
switch (md->which_subclass) {
// material read from file: interpolate to get properties at r
case meep_geom::material_data::MATERIAL_FILE:
if (md->epsilon_data) { meep_geom::epsilon_file_material(md, p); }
else {
material = (meep_geom::material_type)default_material;
}
return;
// material specified by user-supplied function: call user
// function to get properties at r.
// Note that we initialize the medium to vacuum, so that
// the user's function only needs to fill in whatever is
// different from vacuum.
case meep_geom::material_data::MATERIAL_USER:
md->medium = meep_geom::medium_struct();
md->user_func(p, md->user_data, &(md->medium));
return;
// position-independent material or metal: there is nothing to do
case meep_geom::material_data::MEDIUM:
case meep_geom::material_data::PERFECT_METAL: return;
default: meep::abort("unknown material type");
}
}
bool mode_solver::using_mu() { return mdata && mdata->mu_inv != NULL; }
void mode_solver::init(int p, bool reset_fields, geometric_object_list *geometry,
meep_geom::material_data *_default_material) {
int have_old_fields = 0;
default_material = _default_material;
n[0] = grid_size.x;
n[1] = grid_size.y;
n[2] = grid_size.z;
if (target_freq != 0.0 && mpb_verbosity >= 1) { meep::master_printf("Target frequency is %g\n", target_freq); }
int true_rank = n[2] > 1 ? 3 : (n[1] > 1 ? 2 : 1);
if (true_rank < dimensions) { dimensions = true_rank; }
else if (true_rank > dimensions) {
if (mpb_verbosity >= 1)
meep::master_printf("WARNING: rank of grid is > dimensions.\n"
" setting extra grid dims. to 1.\n");
// force extra dims to be 1
if (dimensions <= 2) { n[2] = 1; }
if (dimensions <= 1) { n[1] = 1; }
}
if (mpb_verbosity >= 1) {
meep::master_printf("Working in %d dimensions.\n", dimensions);
meep::master_printf("Grid size is %d x %d x %d.\n", n[0], n[1], n[2]);
}
int block_size;
if (eigensolver_block_size != 0 && eigensolver_block_size < num_bands) {
block_size = eigensolver_block_size;
if (block_size < 0) {
// Guess a block_size near -block_size, chosen so that all blocks are nearly equal in size
block_size = (num_bands - block_size - 1) / (-block_size);
block_size = (num_bands + block_size - 1) / block_size;
}
if (mpb_verbosity >= 1)
meep::master_printf("Solving for %d bands at a time.\n", block_size);
}
else {
block_size = num_bands;
}
if (mdata) {
if (n[0] == mdata->nx && n[1] == mdata->ny && n[2] == mdata->nz &&
block_size == Hblock.alloc_p && num_bands == H.p &&
eigensolver_nwork + (mdata->mu_inv != NULL) == nwork_alloc) {
have_old_fields = 1;
}
else {
destroy_evectmatrix(H);
for (int i = 0; i < nwork_alloc; ++i) {
destroy_evectmatrix(W[i]);
}
if (Hblock.data != H.data) { destroy_evectmatrix(Hblock); }
if (muinvH.data != H.data) { destroy_evectmatrix(muinvH); }
}
destroy_maxwell_target_data(mtdata);
mtdata = NULL;
destroy_maxwell_data(mdata);
mdata = NULL;
curfield_reset();
}
else {
srand(time(NULL));
}
if (deterministic) {
// seed should be the same for each run, although
// it should be different for each process.
// TODO: MPI
// int rank = meep::my_rank();
srand(314159); // * (rank + 1));
}
if (mpb_verbosity >= 1)
meep::master_printf("Creating Maxwell data...\n");
mdata = create_maxwell_data(n[0], n[1], n[2], &local_N, &N_start, &alloc_N, block_size,
NUM_FFT_BANDS);
if (target_freq != 0.0) { mtdata = create_maxwell_target_data(mdata, target_freq); }
init_epsilon(geometry);
if (check_maxwell_dielectric(mdata, 0)) { meep::abort("invalid dielectric function for MPB"); }
if (!have_old_fields) {
if (mpb_verbosity >= 1)
meep::master_printf("Allocating fields...\n");
int N = n[0] * n[1] * n[2];
int c = 2;
H = create_evectmatrix(N, c, num_bands, local_N, N_start, alloc_N);
nwork_alloc = eigensolver_nwork + (mdata->mu_inv != NULL);
for (int i = 0; i < nwork_alloc; ++i) {
W[i] = create_evectmatrix(N, c, block_size, local_N, N_start, alloc_N);
}
if (block_size < num_bands) {
Hblock = create_evectmatrix(N, c, block_size, local_N, N_start, alloc_N);
}
else {
Hblock = H;
}
if (using_mu() && block_size < num_bands) {
muinvH = create_evectmatrix(N, c, num_bands, local_N, N_start, alloc_N);
}
else {
muinvH = H;
}
}
set_parity(p);
if (!have_old_fields || reset_fields) { randomize_fields(); }
evectmatrix_flops = eigensolver_flops;
}
void mode_solver::init_epsilon(geometric_object_list *geometry_in) {
// Persist geometry data and move it out of the input argument.
clear_geometry_list();
if (geometry_in->num_items && geometry_in->items) {
geometry_list.items = geometry_in->items;
geometry_list.num_items = geometry_in->num_items;
geometry_in->items = NULL;
geometry_in->num_items = 0;
}
mpb_real no_size_x = geometry_lattice.size.x == 0 ? 1 : geometry_lattice.size.x;
mpb_real no_size_y = geometry_lattice.size.y == 0 ? 1 : geometry_lattice.size.y;
mpb_real no_size_z = geometry_lattice.size.z == 0 ? 1 : geometry_lattice.size.z;
if (mpb_verbosity >= 1)
meep::master_printf("Mesh size is %d.\n", mesh_size);
Rm.c0 = vector3_scale(no_size_x, geometry_lattice.basis.c0);
Rm.c1 = vector3_scale(no_size_y, geometry_lattice.basis.c1);
Rm.c2 = vector3_scale(no_size_z, geometry_lattice.basis.c2);
if (mpb_verbosity >= 1) {
meep::master_printf("Lattice vectors:\n");
meep::master_printf(" (%g, %g, %g)\n", Rm.c0.x, Rm.c0.y, Rm.c0.z);
meep::master_printf(" (%g, %g, %g)\n", Rm.c1.x, Rm.c1.y, Rm.c1.z);
meep::master_printf(" (%g, %g, %g)\n", Rm.c2.x, Rm.c2.y, Rm.c2.z);
}
vol = fabs(matrix3x3_determinant(Rm));
if (mpb_verbosity >= 1)
meep::master_printf("Cell volume = %g\n", vol);
Gm = matrix3x3_inverse(matrix3x3_transpose(Rm));
if (mpb_verbosity >= 1) {
meep::master_printf("Reciprocal lattice vectors (/ 2 pi):\n");
meep::master_printf(" (%g, %g, %g)\n", Gm.c0.x, Gm.c0.y, Gm.c0.z);
meep::master_printf(" (%g, %g, %g)\n", Gm.c1.x, Gm.c1.y, Gm.c1.z);
meep::master_printf(" (%g, %g, %g)\n", Gm.c2.x, Gm.c2.y, Gm.c2.z);
}
matrix3x3_to_arr(R, Rm);
matrix3x3_to_arr(G, Gm);
geom_fix_object_list(geometry_list);
if (mpb_verbosity >= 1)
meep::master_printf("Geometric objects:\n");
if (meep::am_master()) {
for (int i = 0; i < geometry_list.num_items; ++i) {
#ifndef WITH_HERMITIAN_EPSILON
meep_geom::medium_struct *mm;
if (meep_geom::is_medium(geometry_list.items[i].material, &mm)) { mm->check_offdiag_im_zero_or_abort(); }
#endif
display_geometric_object_info(5, geometry_list.items[i]);
// meep_geom::medium_struct *mm;
// if (meep_geom::is_medium(geometry.items[i].material, &mm)) {
// meep::master_printf("%*sdielectric constant epsilon diagonal = (%g,%g,%g)\n", 5 + 5, "",
// mm->epsilon_diag.x, mm->epsilon_diag.y, mm->epsilon_diag.z);
// }
}
}
{
// Replace 0 with 1e-20 for no size
vector3 tmp_size;
tmp_size.x = geometry_lattice.size.x == 0 ? 1e-20 : geometry_lattice.size.x;
tmp_size.y = geometry_lattice.size.y == 0 ? 1e-20 : geometry_lattice.size.y;
tmp_size.z = geometry_lattice.size.z == 0 ? 1e-20 : geometry_lattice.size.z;
geom_box b0;
b0.low = vector3_plus(geometry_center, vector3_scale(-0.5, tmp_size));
b0.high = vector3_plus(geometry_center, vector3_scale(0.5, tmp_size));
/* pad tree boundaries to allow for sub-pixel averaging */
b0.low.x -= tmp_size.x / mdata->nx;
b0.low.y -= tmp_size.y / mdata->ny;
b0.low.z -= tmp_size.z / mdata->nz;
b0.high.x += tmp_size.x / mdata->nx;
b0.high.y += tmp_size.y / mdata->ny;
b0.high.z += tmp_size.z / mdata->nz;
geometry_tree = create_geom_box_tree0(geometry_list, b0);
}
if (verbose && meep::am_master()) {
if (mpb_verbosity >= 1)
meep::master_printf("Geometry object bounding box tree:\n");
display_geom_box_tree(5, geometry_tree);
}
int tree_depth;
int tree_nobjects;
geom_box_tree_stats(geometry_tree, &tree_depth, &tree_nobjects);
if (mpb_verbosity >= 1)
meep::master_printf("Geometric object tree has depth %d and %d object nodes"
" (vs. %d actual objects)\n",
tree_depth, tree_nobjects, geometry_list.num_items);
// restricted_tree = geometry_tree;
reset_epsilon(&geometry_list);
}
void mode_solver::reset_epsilon(geometric_object_list *geometry) {
int mesh[3] = {
mesh_size,
(dimensions > 1) ? mesh_size : 1,
(dimensions > 2) ? mesh_size : 1,
};
if (!epsilon_input_file.empty()) {
default_material = meep_geom::make_file_material(epsilon_input_file.c_str());
}
// TODO: support mu_input_file
// if (!mu_input_file.empty()) {
// }
if (mpb_verbosity >= 1)
meep::master_printf("Initializing epsilon function...\n");
set_maxwell_dielectric(mdata, mesh, R, G, dielectric_function, mean_epsilon_func,
static_cast<void *>(this));
if (has_mu(geometry)) {
if (mpb_verbosity >= 1)
meep::master_printf("Initializing mu function...\n");
eps = false;
set_maxwell_mu(mdata, mesh, R, G, dielectric_function, mean_epsilon_func,
static_cast<void *>(this));
eps = true;
}
}
bool mode_solver::has_mu(geometric_object_list *geometry) {
// TODO: mu_file_func
if (material_has_mu(default_material) || force_mu) { return true; }
for (int i = 0; i < geometry->num_items; ++i) {
if (material_has_mu(geometry->items[i].material)) { return true; }
}
return false;
}
bool mode_solver::material_has_mu(void *mt) {
meep_geom::material_type mat = (meep_geom::material_type)mt;
meep_geom::medium_struct *m = &mat->medium;
if (mat->which_subclass != meep_geom::material_data::PERFECT_METAL) {
bool has_nonzero_mu_offdiag = false;
#ifdef WITH_HERMITIAN_EPSILON
if (m->mu_offdiag.x.re != 0 || m->mu_offdiag.x.im != 0 || m->mu_offdiag.y.re != 0 ||
m->mu_offdiag.y.im != 0 || m->mu_offdiag.z.re != 0 || m->mu_offdiag.z.im != 0) {
has_nonzero_mu_offdiag = true;
}
#else
if (m->mu_offdiag.x.re != 0 || m->mu_offdiag.y.re != 0 || m->mu_offdiag.z.re != 0) {
has_nonzero_mu_offdiag = true;
}
#endif
if (m->mu_diag.x != 1 || m->mu_diag.y != 1 || m->mu_diag.z != 1 || has_nonzero_mu_offdiag) {
return true;
}
}
return false;
}
void mode_solver::set_parity(integer p) {
if (!mdata) {
meep::master_fprintf(stderr, "init must be called before set-parity!\n");
return;
}
if (p == -1) { p = last_parity < 0 ? NO_PARITY : last_parity; }
set_maxwell_data_parity(mdata, p);
if (mdata->parity != p) {
meep::master_fprintf(stderr, "k vector incompatible with parity\n");
exit(EXIT_FAILURE);
}
if (mpb_verbosity >= 1)
meep::master_printf("Solving for band polarization: %s.\n", parity_string(mdata));
last_parity = p;
set_kpoint_index(0); /* reset index */
}
void mode_solver::set_num_bands(int nb) {
num_bands = nb;
freqs.resize(nb);
}
int mode_solver::get_kpoint_index() { return kpoint_index; }
void mode_solver::set_kpoint_index(int i) { kpoint_index = i; }
void mode_solver::randomize_fields() {
if (!mdata) { return; }
if (mpb_verbosity >= 1)
meep::master_printf("Initializing fields to random numbers...\n");
for (int i = 0; i < H.n * H.p; ++i) {