-
Notifications
You must be signed in to change notification settings - Fork 641
/
Copy pathstructure.cpp
1541 lines (1369 loc) · 56.4 KB
/
structure.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 "meep-ctl.hpp"
#include <ctlgeom.h>
#include <string.h>
using namespace ctlio;
#define master_printf meep::master_printf
#define MTS material_type_struct
typedef struct {
double m00, m01, m02, m11, m12, m22;
} symmetric_matrix;
/* rotate A by a unitary (real) rotation matrix R:
RAR = transpose(R) * A * R
*/
void sym_matrix_rotate(symmetric_matrix *RAR, const symmetric_matrix *A_, const double R[3][3]) {
int i, j;
double A[3][3], AR[3][3];
A[0][0] = A_->m00;
A[1][1] = A_->m11;
A[2][2] = A_->m22;
A[0][1] = A[1][0] = A_->m01;
A[0][2] = A[2][0] = A_->m02;
A[1][2] = A[2][1] = A_->m12;
for (i = 0; i < 3; ++i)
for (j = 0; j < 3; ++j)
AR[i][j] = A[i][0] * R[0][j] + A[i][1] * R[1][j] + A[i][2] * R[2][j];
for (i = 0; i < 3; ++i)
for (j = i; j < 3; ++j)
A[i][j] = R[0][i] * AR[0][j] + R[1][i] * AR[1][j] + R[2][i] * AR[2][j];
RAR->m00 = A[0][0];
RAR->m11 = A[1][1];
RAR->m22 = A[2][2];
RAR->m01 = A[0][1];
RAR->m02 = A[0][2];
RAR->m12 = A[1][2];
}
/* Set Vinv = inverse of V, where both V and Vinv are real-symmetric matrices.*/
void sym_matrix_invert(symmetric_matrix *Vinv, const symmetric_matrix *V) {
double m00 = V->m00, m11 = V->m11, m22 = V->m22;
double m01 = V->m01, m02 = V->m02, m12 = V->m12;
if (m01 == 0.0 && m02 == 0.0 && m12 == 0.0) {
/* optimize common case of a diagonal matrix: */
Vinv->m00 = 1.0 / m00;
Vinv->m11 = 1.0 / m11;
Vinv->m22 = 1.0 / m22;
Vinv->m01 = Vinv->m02 = Vinv->m12 = 0.0;
} else {
double detinv;
/* compute the determinant: */
detinv = m00 * m11 * m22 - m02 * m11 * m02 + 2.0 * m01 * m12 * m02 - m01 * m01 * m22 -
m12 * m12 * m00;
if (detinv == 0.0) meep::abort("singular 3x3 matrix");
detinv = 1.0 / detinv;
Vinv->m00 = detinv * (m11 * m22 - m12 * m12);
Vinv->m11 = detinv * (m00 * m22 - m02 * m02);
Vinv->m22 = detinv * (m11 * m00 - m01 * m01);
Vinv->m02 = detinv * (m01 * m12 - m11 * m02);
Vinv->m01 = detinv * (m12 * m02 - m01 * m22);
Vinv->m12 = detinv * (m01 * m02 - m00 * m12);
}
}
/* Returns whether or not V is positive-definite. */
int sym_matrix_positive_definite(symmetric_matrix *V) {
double det2, det3;
double m00 = V->m00, m11 = V->m11, m22 = V->m22;
#if defined(WITH_HERMITIAN_EPSILON)
scalar_complex m01 = V->m01, m02 = V->m02, m12 = V->m12;
det2 = m00 * m11 - CSCALAR_NORMSQR(m01);
det3 = det2 * m22 - m11 * CSCALAR_NORMSQR(m02) - CSCALAR_NORMSQR(m12) * m00 +
2.0 * ((m01.re * m12.re - m01.im * m12.im) * m02.re +
(m01.re * m12.im + m01.im * m12.re) * m02.im);
#else /* real matrix */
double m01 = V->m01, m02 = V->m02, m12 = V->m12;
det2 = m00 * m11 - m01 * m01;
det3 = det2 * m22 - m02 * m11 * m02 + 2.0 * m01 * m12 * m02 - m12 * m12 * m00;
#endif /* real matrix */
return (m00 > 0.0 && det2 > 0.0 && det3 > 0.0);
}
static meep::ndim dim = meep::D3;
/***********************************************************************/
void set_dimensions(int dims) {
if (dims == CYLINDRICAL) {
dimensions = 2;
dim = meep::Dcyl;
} else {
dimensions = dims;
dim = meep::ndim(dims - 1);
}
}
vector3 vec_to_vector3(const meep::vec &pt) {
vector3 v3;
switch (pt.dim) {
case meep::D1:
v3.x = 0;
v3.y = 0;
v3.z = pt.z();
break;
case meep::D2:
v3.x = pt.x();
v3.y = pt.y();
v3.z = 0;
break;
case meep::D3:
v3.x = pt.x();
v3.y = pt.y();
v3.z = pt.z();
break;
case meep::Dcyl:
v3.x = pt.r();
v3.y = 0;
v3.z = pt.z();
break;
}
return v3;
}
meep::vec vector3_to_vec(const vector3 v3) {
switch (dim) {
case meep::D1: return meep::vec(v3.z);
case meep::D2: return meep::vec(v3.x, v3.y);
case meep::D3: return meep::vec(v3.x, v3.y, v3.z);
case meep::Dcyl: return meep::veccyl(v3.x, v3.z);
default: meep::abort("unknown dimensionality in vector3_to_vec");
}
}
static meep::vec geometry_edge; // geometry_lattice.size / 2
static geom_box gv2box(const meep::volume &v) {
geom_box box;
box.low = vec_to_vector3(v.get_min_corner());
box.high = vec_to_vector3(v.get_max_corner());
return box;
}
/***********************************************************************/
static meep::realnum *epsilon_data = NULL;
static size_t epsilon_dims[3] = {0, 0, 0};
static void read_epsilon_file(const char *eps_input_file) {
delete[] epsilon_data;
epsilon_data = NULL;
epsilon_dims[0] = epsilon_dims[1] = epsilon_dims[2] = 1;
if (eps_input_file && eps_input_file[0]) { // file specified
char *fname = new char[strlen(eps_input_file) + 1];
strcpy(fname, eps_input_file);
// parse epsilon-input-file as "fname.h5:dataname"
char *dataname = strrchr(fname, ':');
if (dataname) *(dataname++) = 0;
meep::h5file eps_file(fname, meep::h5file::READONLY, false);
int rank; // ignored since rank < 3 is equivalent to singleton dims
epsilon_data = eps_file.read(dataname, &rank, epsilon_dims, 3);
master_printf("read in %zdx%zdx%zd epsilon-input-file \"%s\"\n", epsilon_dims[0],
epsilon_dims[1], epsilon_dims[2], eps_input_file);
delete[] fname;
}
}
// return material of the point p from the file (assumed already read)
static void epsilon_file_material(material_type &m, vector3 p) {
material_type_copy(&default_material, &m);
if (!epsilon_data) return;
if (m.which_subclass != MTS::MEDIUM)
meep::abort("epsilon-input-file only works with a type=medium default-material");
medium *mm = m.subclass.medium_data;
double rx =
geometry_lattice.size.x == 0 ? 0 : 0.5 + (p.x - geometry_center.x) / geometry_lattice.size.x;
double ry =
geometry_lattice.size.y == 0 ? 0 : 0.5 + (p.y - geometry_center.y) / geometry_lattice.size.y;
double rz =
geometry_lattice.size.z == 0 ? 0 : 0.5 + (p.z - geometry_center.z) / geometry_lattice.size.z;
mm->epsilon_diag.x = mm->epsilon_diag.y = mm->epsilon_diag.z = meep::linear_interpolate(
rx, ry, rz, epsilon_data, epsilon_dims[0], epsilon_dims[1], epsilon_dims[2], 1);
mm->epsilon_offdiag.x = mm->epsilon_offdiag.y = mm->epsilon_offdiag.z = 0;
}
/***********************************************************************/
struct pol {
susceptibility user_s;
struct pol *next;
};
// structure to hold a conductivity profile (for scalar absorbing layers)
struct cond_profile {
double L; // thickness
int N; // number of points prof[n] from 0..N corresponding to 0..L
double *prof; // (NULL if none)
};
class geom_epsilon : public meep::material_function {
geometric_object_list geometry;
geom_box_tree geometry_tree;
geom_box_tree restricted_tree;
cond_profile cond[5][2]; // [direction][side]
public:
geom_epsilon(geometric_object_list g, material_type_list mlist, const meep::volume &v);
virtual ~geom_epsilon();
virtual void set_cond_profile(meep::direction, meep::boundary_side, double L, double dx,
double (*prof)(int, double *, void *), void *, double R);
virtual void set_volume(const meep::volume &v);
virtual void unset_volume(void);
virtual bool has_chi3(meep::component c);
virtual double chi3(meep::component c, const meep::vec &r);
virtual bool has_chi2(meep::component c);
virtual double chi2(meep::component c, const meep::vec &r);
virtual bool has_mu();
virtual bool has_conductivity(meep::component c);
virtual double conductivity(meep::component c, const meep::vec &r);
virtual double chi1p1(meep::field_type ft, const meep::vec &r);
virtual void eff_chi1inv_row(meep::component c, double chi1inv_row[3], const meep::volume &v,
double tol, int maxeval);
void fallback_chi1inv_row(meep::component c, double chi1inv_row[3], const meep::volume &v,
double tol, int maxeval);
virtual void sigma_row(meep::component c, double sigrow[3], const meep::vec &r);
void add_susceptibilities(meep::structure *s);
void add_susceptibilities(meep::field_type ft, meep::structure *s);
private:
bool get_material_pt(material_type &material, const meep::vec &r);
material_type_list extra_materials;
pol *current_pol;
};
geom_epsilon::geom_epsilon(geometric_object_list g, material_type_list mlist,
const meep::volume &v) {
geometry = g; // don't bother making a copy, only used in one place
extra_materials = mlist;
current_pol = NULL;
FOR_DIRECTIONS(d) FOR_SIDES(b) { cond[d][b].prof = NULL; }
if (meep::am_master()) {
for (int i = 0; i < geometry.num_items; ++i) {
display_geometric_object_info(5, geometry.items[i]);
if (geometry.items[i].material.which_subclass == MTS::MEDIUM)
printf("%*sdielectric constant epsilon diagonal = (%g,%g,%g)\n", 5 + 5, "",
geometry.items[i].material.subclass.medium_data->epsilon_diag.x,
geometry.items[i].material.subclass.medium_data->epsilon_diag.y,
geometry.items[i].material.subclass.medium_data->epsilon_diag.z);
}
}
geom_fix_object_list(geometry);
geom_box box = gv2box(v);
geometry_tree = create_geom_box_tree0(geometry, box);
if (verbose && meep::am_master()) {
printf("Geometric-object bounding-box tree:\n");
display_geom_box_tree(5, geometry_tree);
int tree_depth, tree_nobjects;
geom_box_tree_stats(geometry_tree, &tree_depth, &tree_nobjects);
master_printf("Geometric object tree has depth %d "
"and %d object nodes (vs. %d actual objects)\n",
tree_depth, tree_nobjects, geometry.num_items);
}
restricted_tree = geometry_tree;
}
geom_epsilon::~geom_epsilon() {
unset_volume();
destroy_geom_box_tree(geometry_tree);
FOR_DIRECTIONS(d) FOR_SIDES(b) {
if (cond[d][b].prof) delete[] cond[d][b].prof;
}
}
void geom_epsilon::set_cond_profile(meep::direction dir, meep::boundary_side side, double L,
double dx, double (*P)(int, double *, void *), void *data,
double R) {
if (cond[dir][side].prof) delete[] cond[dir][side].prof;
int N = int(L / dx + 0.5);
cond[dir][side].L = L;
cond[dir][side].N = N;
double *prof = cond[dir][side].prof = new double[N + 1];
double umin = 0, umax = 1, esterr;
int errflag;
double prof_int =
adaptive_integration(P, &umin, &umax, 1, data, 1e-9, 1e-4, 50000, &esterr, &errflag);
double prefac = (-log(R)) / (4 * L * prof_int);
for (int i = 0; i <= N; ++i) {
double u = double(i) / N;
prof[i] = prefac * P(1, &u, data);
}
}
void geom_epsilon::unset_volume(void) {
if (restricted_tree != geometry_tree) {
destroy_geom_box_tree(restricted_tree);
restricted_tree = geometry_tree;
}
}
void geom_epsilon::set_volume(const meep::volume &v) {
unset_volume();
geom_box box = gv2box(v);
restricted_tree = create_geom_box_tree0(geometry, box);
}
static material_type eval_material_func(function material_func, vector3 p) {
SCM pscm = ctl_convert_vector3_to_scm(p);
material_type material;
SCM mo;
mo = gh_call1(material_func, pscm);
material_type_input(mo, &material);
while (material.which_subclass == MTS::MATERIAL_FUNCTION) {
material_type m;
mo = gh_call1(material.subclass.material_function_data->material_func, pscm);
material_type_input(mo, &m);
material_type_destroy(material);
material = m;
}
if (material.which_subclass == MTS::MATERIAL_TYPE_SELF) { epsilon_file_material(material, p); }
CK(material.which_subclass != MTS::MATERIAL_FUNCTION, "infinite loop in material functions");
return material;
}
static int variable_material(int which_subclass) {
return (which_subclass == MTS::MATERIAL_FUNCTION);
}
static bool is_metal(meep::field_type ft, const material_type *material) {
if (ft == meep::E_stuff) switch (material->which_subclass) {
case MTS::MEDIUM:
return (material->subclass.medium_data->epsilon_diag.x < 0 ||
material->subclass.medium_data->epsilon_diag.y < 0 ||
material->subclass.medium_data->epsilon_diag.z < 0);
case MTS::PERFECT_METAL: return true;
default: meep::abort("unknown material type");
}
else
switch (material->which_subclass) {
case MTS::MEDIUM:
return (material->subclass.medium_data->mu_diag.x < 0 ||
material->subclass.medium_data->mu_diag.y < 0 ||
material->subclass.medium_data->mu_diag.z < 0);
case MTS::PERFECT_METAL:
return false; // is an electric conductor, but not a magnetic conductor
default: meep::abort("unknown material type");
}
}
static void material_epsmu(meep::field_type ft, material_type material, symmetric_matrix *epsmu,
symmetric_matrix *epsmu_inv) {
if (ft == meep::E_stuff) switch (material.which_subclass) {
case MTS::MEDIUM: {
epsmu->m00 = material.subclass.medium_data->epsilon_diag.x;
epsmu->m11 = material.subclass.medium_data->epsilon_diag.y;
epsmu->m22 = material.subclass.medium_data->epsilon_diag.z;
epsmu->m01 = material.subclass.medium_data->epsilon_offdiag.x;
epsmu->m02 = material.subclass.medium_data->epsilon_offdiag.y;
epsmu->m12 = material.subclass.medium_data->epsilon_offdiag.z;
sym_matrix_invert(epsmu_inv, epsmu);
break;
}
case MTS::PERFECT_METAL: {
epsmu->m00 = -meep::infinity;
epsmu->m11 = -meep::infinity;
epsmu->m22 = -meep::infinity;
epsmu_inv->m00 = -0.0;
epsmu_inv->m11 = -0.0;
epsmu_inv->m22 = -0.0;
epsmu->m01 = epsmu->m02 = epsmu->m12 = 0.0;
epsmu_inv->m01 = epsmu_inv->m02 = epsmu_inv->m12 = 0.0;
break;
}
default: meep::abort("unknown material type");
}
else
switch (material.which_subclass) {
case MTS::MEDIUM: {
epsmu->m00 = material.subclass.medium_data->mu_diag.x;
epsmu->m11 = material.subclass.medium_data->mu_diag.y;
epsmu->m22 = material.subclass.medium_data->mu_diag.z;
epsmu->m01 = material.subclass.medium_data->mu_offdiag.x;
epsmu->m02 = material.subclass.medium_data->mu_offdiag.y;
epsmu->m12 = material.subclass.medium_data->mu_offdiag.z;
sym_matrix_invert(epsmu_inv, epsmu);
break;
}
case MTS::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;
epsmu->m01 = epsmu->m02 = epsmu->m12 = 0.0;
epsmu_inv->m01 = epsmu_inv->m02 = epsmu_inv->m12 = 0.0;
break;
}
default: meep::abort("unknown material type");
}
}
bool geom_epsilon::get_material_pt(material_type &material, const meep::vec &r) {
vector3 p = vec_to_vector3(r);
boolean inobject;
material = material_of_unshifted_point_in_tree_inobject(p, restricted_tree, &inobject);
bool destroy_material = false;
if (!inobject && epsilon_data) {
epsilon_file_material(material, p);
destroy_material = true;
} else if (material.which_subclass == MTS::MATERIAL_TYPE_SELF) {
if (epsilon_data) {
epsilon_file_material(material, p);
destroy_material = true;
} else
material = default_material;
} else if (material.which_subclass == MTS::MATERIAL_FUNCTION) {
material = eval_material_func(material.subclass.material_function_data->material_func, p);
destroy_material = true;
}
return destroy_material;
}
// returns trace of the tensor diagonal
double geom_epsilon::chi1p1(meep::field_type ft, const meep::vec &r) {
symmetric_matrix chi1p1, chi1p1_inv;
#ifdef DEBUG
vector3 p = vec_to_vector3(r);
if (p.x < restricted_tree->b.low.x || p.y < restricted_tree->b.low.y ||
p.z < restricted_tree->b.low.z || p.x > restricted_tree->b.high.x ||
p.y > restricted_tree->b.high.y || p.z > restricted_tree->b.high.z)
meep::abort("invalid point (%g,%g,%g)\n", p.x, p.y, p.z);
#endif
material_type material;
bool destroy_material = get_material_pt(material, r);
material_epsmu(ft, material, &chi1p1, &chi1p1_inv);
if (destroy_material) material_type_destroy(material);
return (chi1p1.m00 + chi1p1.m11 + chi1p1.m22) / 3;
}
/* Find frontmost object in v, along with the constant material behind it.
Returns false if material behind the object is not constant.
Requires moderately horrifying logic to figure things out properly,
stolen from MPB. */
static bool get_front_object(const meep::volume &v, geom_box_tree geometry_tree, vector3 &pcenter,
const geometric_object **o_front, vector3 &shiftby_front,
material_type &mat_front, material_type &mat_behind) {
vector3 p;
const geometric_object *o1 = 0, *o2 = 0;
vector3 shiftby1 = {0, 0, 0}, shiftby2 = {0, 0, 0};
geom_box pixel;
material_type mat1, mat2;
int id1 = -1, id2 = -1;
const int num_neighbors[3] = {3, 5, 9};
const int neighbors[3][9][3] = {{{0, 0, 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},
{-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}}};
pixel = gv2box(v);
pcenter = p = vec_to_vector3(v.center());
double d1, d2, d3;
d1 = (pixel.high.x - pixel.low.x) * 0.5;
d2 = (pixel.high.y - pixel.low.y) * 0.5;
d3 = (pixel.high.z - pixel.low.z) * 0.5;
for (int i = 0; i < num_neighbors[dimensions - 1]; ++i) {
const geometric_object *o;
material_type mat;
vector3 q, 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;
o = object_of_point_in_tree(q, geometry_tree, &shiftby, &id);
if ((id == id1 && vector3_equal(shiftby, shiftby1)) ||
(id == id2 && vector3_equal(shiftby, shiftby2)))
continue;
mat = (o && o->material.which_subclass != MTS::MATERIAL_TYPE_SELF) ? o->material
: default_material;
if (id1 == -1) {
o1 = o;
shiftby1 = shiftby;
id1 = id;
mat1 = mat;
} else if (id2 == -1 ||
((id >= id1 && id >= id2) && (id1 == id2 || material_type_equal(&mat1, &mat2)))) {
o2 = o;
shiftby2 = shiftby;
id2 = id;
mat2 = mat;
} else if (!(id1 < id2 && (id1 == id || material_type_equal(&mat1, &mat))) &&
!(id2 < id1 && (id2 == id || material_type_equal(&mat2, &mat))))
return false;
}
// 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;
}
if ((o1 && variable_material(o1->material.which_subclass)) ||
(o2 && variable_material(o2->material.which_subclass)) ||
((variable_material(default_material.which_subclass) || epsilon_data) &&
(!o1 || !o2 || o1->material.which_subclass == MTS::MATERIAL_TYPE_SELF ||
o2->material.which_subclass == MTS::MATERIAL_TYPE_SELF)))
return false;
if (id1 >= id2) {
*o_front = o1;
shiftby_front = shiftby1;
mat_front = mat1;
if (id1 == id2)
mat_behind = mat1;
else
mat_behind = mat2;
}
if (id2 > id1) {
*o_front = o2;
shiftby_front = shiftby2;
mat_front = mat2;
mat_behind = mat1;
}
return true;
}
void geom_epsilon::eff_chi1inv_row(meep::component c, double chi1inv_row[3], const meep::volume &v,
double tol, int maxeval) {
const geometric_object *o;
material_type mat, mat_behind;
symmetric_matrix meps, meps_inv;
vector3 p, shiftby, normal;
bool destroy_material = false;
if (maxeval == 0 || !get_front_object(v, geometry_tree, p, &o, shiftby, mat, mat_behind)) {
noavg:
destroy_material = get_material_pt(mat, v.center());
trivial:
material_epsmu(meep::type(c), mat, &meps, &meps_inv);
switch (component_direction(c)) {
case meep::X:
case meep::R:
chi1inv_row[0] = meps_inv.m00;
chi1inv_row[1] = meps_inv.m01;
chi1inv_row[2] = meps_inv.m02;
break;
case meep::Y:
case meep::P:
chi1inv_row[0] = meps_inv.m01;
chi1inv_row[1] = meps_inv.m11;
chi1inv_row[2] = meps_inv.m12;
break;
case meep::Z:
chi1inv_row[0] = meps_inv.m02;
chi1inv_row[1] = meps_inv.m12;
chi1inv_row[2] = meps_inv.m22;
break;
case meep::NO_DIRECTION: chi1inv_row[0] = chi1inv_row[1] = chi1inv_row[2] = 0;
}
if (destroy_material) material_type_destroy(mat);
return;
}
// FIXME: reimplement support for fallback integration, without
// messing up anisotropic support
// if (!get_front_object(v, geometry_tree,
// p, &o, shiftby, mat, mat_behind)) {
// fallback_chi1inv_row(c, chi1inv_row, v, tol, maxeval);
// return;
// }
/* check for trivial case of only one object/material */
if (material_type_equal(&mat, &mat_behind)) goto trivial;
// it doesn't make sense to average metals (electric or magnetic)
if (is_metal(meep::type(c), &mat) || is_metal(meep::type(c), &mat_behind)) goto noavg;
normal = unit_vector3(normal_to_fixed_object(vector3_minus(p, shiftby), *o));
if (normal.x == 0 && normal.y == 0 && normal.z == 0)
goto noavg; // couldn't get normal vector for this point, punt
geom_box pixel = gv2box(v);
pixel.low = vector3_minus(pixel.low, shiftby);
pixel.high = vector3_minus(pixel.high, shiftby);
double fill = box_overlap_with_object(pixel, *o, tol, maxeval);
material_epsmu(meep::type(c), mat, &meps, &meps_inv);
symmetric_matrix eps2, epsinv2;
symmetric_matrix eps1, delta;
double Rot[3][3];
material_epsmu(meep::type(c), mat_behind, &eps2, &epsinv2);
eps1 = meps;
Rot[0][0] = normal.x;
Rot[1][0] = normal.y;
Rot[2][0] = normal.z;
if (fabs(normal.x) > 1e-2 || fabs(normal.y) > 1e-2) {
Rot[0][2] = normal.y;
Rot[1][2] = -normal.x;
Rot[2][2] = 0;
} else { /* n is not parallel to z direction, use (x x n) instead */
Rot[0][2] = 0;
Rot[1][2] = -normal.z;
Rot[2][2] = normal.y;
}
{ /* 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 */
sym_matrix_rotate(&eps1, &eps1, Rot);
sym_matrix_rotate(&eps2, &eps2, Rot);
#define AVG (fill * (EXPR(eps1)) + (1 - fill) * (EXPR(eps2)))
#define SQR(x) ((x) * (x))
#define EXPR(eps) (-1 / eps.m00)
delta.m00 = AVG;
#undef EXPR
#define EXPR(eps) (eps.m11 - SQR(eps.m01) / eps.m00)
delta.m11 = AVG;
#undef EXPR
#define EXPR(eps) (eps.m22 - SQR(eps.m02) / eps.m00)
delta.m22 = AVG;
#undef EXPR
#define EXPR(eps) (eps.m01 / eps.m00)
delta.m01 = AVG;
#undef EXPR
#define EXPR(eps) (eps.m02 / eps.m00)
delta.m02 = AVG;
#undef EXPR
#define EXPR(eps) (eps.m12 - eps.m02 * eps.m01 / eps.m00)
delta.m12 = AVG;
#undef EXPR
meps.m00 = -1 / delta.m00;
meps.m11 = delta.m11 - SQR(delta.m01) / delta.m00;
meps.m22 = delta.m22 - SQR(delta.m02) / delta.m00;
meps.m01 = -delta.m01 / delta.m00;
meps.m02 = -delta.m02 / delta.m00;
meps.m12 = delta.m12 - (delta.m02 * delta.m01) / delta.m00;
#undef SQR
#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]);
sym_matrix_rotate(&meps, &meps, Rot); /* rotate back */
#undef SWAP
#ifdef DEBUG
if (!sym_matrix_positive_definite(&meps))
meep::abort("negative mean epsilon from Kottke algorithm");
#endif
sym_matrix_invert(&meps_inv, &meps);
switch (component_direction(c)) {
case meep::X:
case meep::R:
chi1inv_row[0] = meps_inv.m00;
chi1inv_row[1] = meps_inv.m01;
chi1inv_row[2] = meps_inv.m02;
break;
case meep::Y:
case meep::P:
chi1inv_row[0] = meps_inv.m01;
chi1inv_row[1] = meps_inv.m11;
chi1inv_row[2] = meps_inv.m12;
break;
case meep::Z:
chi1inv_row[0] = meps_inv.m02;
chi1inv_row[1] = meps_inv.m12;
chi1inv_row[2] = meps_inv.m22;
break;
case meep::NO_DIRECTION: chi1inv_row[0] = chi1inv_row[1] = chi1inv_row[2] = 0;
}
}
static int eps_ever_negative = 0;
static meep::field_type func_ft = meep::E_stuff;
#ifdef CTL_HAS_COMPLEX_INTEGRATION
static cnumber ceps_func(int n, number *x, void *geomeps_) {
geom_epsilon *geomeps = (geom_epsilon *)geomeps_;
vector3 p = {0, 0, 0};
p.x = x[0];
p.y = n > 1 ? x[1] : 0;
p.z = n > 2 ? x[2] : 0;
double s = 1;
if (dim == meep::Dcyl) {
double py = p.y;
p.y = p.z;
p.z = py;
s = p.x;
}
cnumber ret;
double ep = geomeps->chi1p1(func_ft, vector3_to_vec(p));
if (ep < 0) eps_ever_negative = 1;
ret.re = ep * s;
ret.im = s / ep;
return ret;
}
#else
static number eps_func(int n, number *x, void *geomeps_) {
geom_epsilon *geomeps = (geom_epsilon *)geomeps_;
vector3 p = {0, 0, 0};
double s = 1;
p.x = x[0];
p.y = n > 1 ? x[1] : 0;
p.z = n > 2 ? x[2] : 0;
if (dim == meep::Dcyl) {
double py = p.y;
p.y = p.z;
p.z = py;
s = p.x;
}
double ep = geomeps->chi1p1(func_ft, vector3_to_vec(p));
if (ep < 0) eps_ever_negative = 1;
return ep * s;
}
static number inveps_func(int n, number *x, void *geomeps_) {
geom_epsilon *geomeps = (geom_epsilon *)geomeps_;
vector3 p = {0, 0, 0};
double s = 1;
p.x = x[0];
p.y = n > 1 ? x[1] : 0;
p.z = n > 2 ? x[2] : 0;
if (dim == meep::Dcyl) {
double py = p.y;
p.y = p.z;
p.z = py;
s = p.x;
}
double ep = geomeps->chi1p1(func_ft, vector3_to_vec(p));
if (ep < 0) eps_ever_negative = 1;
return s / ep;
}
#endif
// fallback meaneps using libctl's adaptive cubature routine
void geom_epsilon::fallback_chi1inv_row(meep::component c, double chi1inv_row[3],
const meep::volume &v, double tol, int maxeval) {
symmetric_matrix chi1p1, chi1p1_inv;
material_type material;
bool destroy_material = get_material_pt(material, v.center());
material_epsmu(meep::type(c), material, &chi1p1, &chi1p1_inv);
if (destroy_material) material_type_destroy(material);
if (chi1p1.m01 != 0 || chi1p1.m02 != 0 || chi1p1.m12 != 0 || chi1p1.m00 != chi1p1.m11 ||
chi1p1.m11 != chi1p1.m22 || chi1p1.m00 != chi1p1.m22) {
int rownum = meep::component_direction(c) % 3;
if (rownum == 0) {
chi1inv_row[0] = chi1p1.m00;
chi1inv_row[1] = chi1p1.m01;
chi1inv_row[2] = chi1p1.m02;
} else if (rownum == 1) {
chi1inv_row[0] = chi1p1.m01;
chi1inv_row[1] = chi1p1.m11;
chi1inv_row[2] = chi1p1.m12;
} else {
chi1inv_row[0] = chi1p1.m02;
chi1inv_row[1] = chi1p1.m12;
chi1inv_row[2] = chi1p1.m22;
}
return;
}
number esterr;
integer errflag, n;
number xmin[3], xmax[3];
vector3 gvmin, gvmax;
gvmin = vec_to_vector3(v.get_min_corner());
gvmax = vec_to_vector3(v.get_max_corner());
xmin[0] = gvmin.x;
xmax[0] = gvmax.x;
if (dim == meep::Dcyl) {
xmin[1] = gvmin.z;
xmin[2] = gvmin.y;
xmax[1] = gvmax.z;
xmax[2] = gvmax.y;
} else {
xmin[1] = gvmin.y;
xmin[2] = gvmin.z;
xmax[1] = gvmax.y;
xmax[2] = gvmax.z;
}
if (xmin[2] == xmax[2])
n = xmin[1] == xmax[1] ? 1 : 2;
else
n = 3;
double vol = 1;
for (int i = 0; i < n; ++i)
vol *= xmax[i] - xmin[i];
if (dim == meep::Dcyl) vol *= (xmin[0] + xmax[0]) * 0.5;
eps_ever_negative = 0;
func_ft = meep::type(c);
double meps, minveps;
#ifdef CTL_HAS_COMPLEX_INTEGRATION
cnumber ret = cadaptive_integration(ceps_func, xmin, xmax, n, (void *)this, 0, tol, maxeval,
&esterr, &errflag);
meps = ret.re / vol;
minveps = ret.im / vol;
#else
meps = adaptive_integration(eps_func, xmin, xmax, n, (void *)this, 0, tol, maxeval, &esterr,
&errflag) /
vol;
minveps = adaptive_integration(inveps_func, xmin, xmax, n, (void *)this, 0, tol, maxeval, &esterr,
&errflag) /
vol;
#endif
if (eps_ever_negative) // averaging negative eps causes instability
minveps = 1.0 / (meps = eps(v.center()));
{
meep::vec gradient(normal_vector(meep::type(c), v));
double n[3] = {0, 0, 0};
double nabsinv = 1.0 / meep::abs(gradient);
LOOP_OVER_DIRECTIONS(gradient.dim, k) { n[k % 3] = gradient.in_direction(k) * nabsinv; }
int rownum = meep::component_direction(c) % 3;
for (int i = 0; i < 3; ++i)
chi1inv_row[i] = n[rownum] * n[i] * (minveps - 1 / meps);
chi1inv_row[rownum] += 1 / meps;
}
}
static double get_chi3(meep::component c, const medium *m) {
switch (c) {
case meep::Er:
case meep::Ex: return m->E_chi3_diag.x;
case meep::Ep:
case meep::Ey: return m->E_chi3_diag.y;
case meep::Ez: return m->E_chi3_diag.z;
case meep::Hr:
case meep::Hx: return m->H_chi3_diag.x;
case meep::Hp:
case meep::Hy: return m->H_chi3_diag.y;
case meep::Hz: return m->H_chi3_diag.z;
default: return 0;
}
}
bool geom_epsilon::has_chi3(meep::component c) {
for (int i = 0; i < geometry.num_items; ++i) {
if (geometry.items[i].material.which_subclass == MTS::MEDIUM) {
if (get_chi3(c, geometry.items[i].material.subclass.medium_data) != 0) return true;
}
}
for (int i = 0; i < extra_materials.num_items; ++i)
if (extra_materials.items[i].which_subclass == MTS::MEDIUM)
if (get_chi3(c, extra_materials.items[i].subclass.medium_data) != 0) return true;
return (default_material.which_subclass == MTS::MEDIUM &&
get_chi3(c, default_material.subclass.medium_data) != 0);
}
double geom_epsilon::chi3(meep::component c, const meep::vec &r) {
material_type material;
bool destroy_material = get_material_pt(material, r);
double chi3_val;
switch (material.which_subclass) {
case MTS::MEDIUM: chi3_val = get_chi3(c, material.subclass.medium_data); break;
default: chi3_val = 0;
}
if (destroy_material) material_type_destroy(material);
return chi3_val;
}
static double get_chi2(meep::component c, const medium *m) {
switch (c) {
case meep::Er:
case meep::Ex: return m->E_chi2_diag.x;
case meep::Ep:
case meep::Ey: return m->E_chi2_diag.y;
case meep::Ez: return m->E_chi2_diag.z;
case meep::Hr:
case meep::Hx: return m->H_chi2_diag.x;
case meep::Hp:
case meep::Hy: return m->H_chi2_diag.y;
case meep::Hz: return m->H_chi2_diag.z;
default: return 0;
}
}
bool geom_epsilon::has_chi2(meep::component c) {
for (int i = 0; i < geometry.num_items; ++i) {
if (geometry.items[i].material.which_subclass == MTS::MEDIUM) {
if (get_chi2(c, geometry.items[i].material.subclass.medium_data) != 0) return true;
}
}
for (int i = 0; i < extra_materials.num_items; ++i)
if (extra_materials.items[i].which_subclass == MTS::MEDIUM)
if (get_chi2(c, extra_materials.items[i].subclass.medium_data) != 0) return true;
return (default_material.which_subclass == MTS::MEDIUM &&
get_chi2(c, default_material.subclass.medium_data) != 0);
}
double geom_epsilon::chi2(meep::component c, const meep::vec &r) {
material_type material;
bool destroy_material = get_material_pt(material, r);
double chi2_val;
switch (material.which_subclass) {
case MTS::MEDIUM: chi2_val = get_chi2(c, material.subclass.medium_data); break;
default: chi2_val = 0;
}
if (destroy_material) material_type_destroy(material);
return chi2_val;
}
static bool mu_not_1(material_type &m) {
return (m.which_subclass == MTS::MEDIUM &&
(m.subclass.medium_data->mu_diag.x != 1 || m.subclass.medium_data->mu_diag.y != 1 ||
m.subclass.medium_data->mu_diag.z != 1 || m.subclass.medium_data->mu_offdiag.x != 0 ||
m.subclass.medium_data->mu_offdiag.y != 0 || m.subclass.medium_data->mu_offdiag.z != 0));