-
Notifications
You must be signed in to change notification settings - Fork 471
/
maxwell.cpp
795 lines (700 loc) · 25.2 KB
/
maxwell.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
// Copyright (c) 2010-2024, Lawrence Livermore National Security, LLC. Produced
// at the Lawrence Livermore National Laboratory. All Rights reserved. See files
// LICENSE and NOTICE for details. LLNL-CODE-806117.
//
// This file is part of the MFEM library. For more information and source code
// availability visit https://mfem.org.
//
// MFEM is free software; you can redistribute it and/or modify it under the
// terms of the BSD-3 license. We welcome feedback and contributions, see file
// CONTRIBUTING.md for details.
//
// MFEM Ultraweak DPG example for Maxwell
//
// Compile with: make maxwell
//
// Sample runs
// maxwell -m ../../data/inline-tri.mesh -ref 4 -o 1 -rnum 1.0
// maxwell -m ../../data/amr-quad.mesh -ref 3 -o 2 -rnum 1.6 -sc
// maxwell -m ../../data/inline-quad.mesh -ref 2 -o 3 -rnum 4.2 -sc
// maxwell -m ../../data/inline-hex.mesh -ref 1 -o 2 -sc -rnum 1.0
// Description:
// This example code demonstrates the use of MFEM to define and solve
// the "ultraweak" (UW) DPG formulation for the (indefinite) Maxwell problem
// ∇×(1/μ ∇×E) - ω² ϵ E = Ĵ , in Ω
// E×n = E₀, on ∂Ω
// It solves a problem with a manufactured solution E_exact being a plane wave
// in the x-component and zero in y (and z) component.
// This example computes and prints out convergence rates for the L² error.
// The DPG UW deals with the First Order System
// i ω μ H + ∇ × E = 0, in Ω
// -i ω ϵ E + ∇ × H = J, in Ω (1)
// E × n = E₀, on ∂Ω
// Note: Ĵ = -iωJ
// in 2D
// E is vector valued and H is scalar.
// (∇ × E, F) = (E, ∇ × F) + < n × E , F>
// or (∇ ⋅ AE , F) = (AE, ∇ F) + < AE ⋅ n, F>
// where A = [0 1; -1 0];
// E ∈ (L²(Ω))² , H ∈ L²(Ω)
// Ê ∈ H^-1/2(Ω)(Γₕ), Ĥ ∈ H^1/2(Γₕ)
// i ω μ (H,F) + (E, ∇ × F) + < AÊ, F > = 0, ∀ F ∈ H¹
// -i ω ϵ (E,G) + (H,∇ × G) + < Ĥ, G × n > = (J,G) ∀ G ∈ H(curl,Ω)
// Ê = E₀ on ∂Ω
// -------------------------------------------------------------------------
// | | E | H | Ê | Ĥ | RHS |
// -------------------------------------------------------------------------
// | F | (E,∇ × F) | i ω μ (H,F) | < Ê, F > | | |
// | | | | | | |
// | G | -i ω ϵ (E,G) | (H,∇ × G) | | < Ĥ, G × n > | (J,G) |
// where (F,G) ∈ H¹ × H(curl,Ω)
// in 3D
// E,H ∈ ((L²(Ω)))³
// Ê ∈ H_0^1/2(Ω)(curl, Γₕ), Ĥ ∈ H^-1/2(curl, Γₕ)
// i ω μ (H,F) + (E,∇ × F) + < Ê, F × n > = 0, ∀ F ∈ H(curl,Ω)
// -i ω ϵ (E,G) + (H,∇ × G) + < Ĥ, G × n > = (J,G) ∀ G ∈ H(curl,Ω)
// Ê × n = E_0 on ∂Ω
// -------------------------------------------------------------------------
// | | E | H | Ê | Ĥ | RHS |
// -------------------------------------------------------------------------
// | F | (E,∇ × F) | i ω μ (H,F) | < n × Ê, F > | | |
// | | | | | | |
// | G | -i ω ϵ (E,G) | (H,∇ × G) | | < n × Ĥ, G > | (J,G) |
// where (F,G) ∈ H(curl,Ω) × H(curl,Ω)
// Here we use the "Adjoint Graph" norm on the test space i.e.,
// ||(F,G)||²ᵥ = ||A^*(F,G)||² + ||(F,G)||² where A is the
// Maxwell operator defined by (1)
// For more information see https://doi.org/10.1016/j.camwa.2021.01.017
#include "mfem.hpp"
#include "util/complexweakform.hpp"
#include "../common/mfem-common.hpp"
#include <fstream>
#include <iostream>
using namespace std;
using namespace mfem;
using namespace mfem::common;
void maxwell_solution(const Vector & X, std::vector<complex<real_t>> &E);
void maxwell_solution_curl(const Vector & X,
std::vector<complex<real_t>> &curlE);
void maxwell_solution_curlcurl(const Vector & X,
std::vector<complex<real_t>> &curlcurlE);
void E_exact_r(const Vector &x, Vector & E_r);
void E_exact_i(const Vector &x, Vector & E_i);
void H_exact_r(const Vector &x, Vector & H_r);
void H_exact_i(const Vector &x, Vector & H_i);
void curlE_exact_r(const Vector &x, Vector &curlE_r);
void curlE_exact_i(const Vector &x, Vector &curlE_i);
void curlH_exact_r(const Vector &x,Vector &curlH_r);
void curlH_exact_i(const Vector &x,Vector &curlH_i);
void curlcurlE_exact_r(const Vector &x, Vector & curlcurlE_r);
void curlcurlE_exact_i(const Vector &x, Vector & curlcurlE_i);
void hatE_exact_r(const Vector & X, Vector & hatE_r);
void hatE_exact_i(const Vector & X, Vector & hatE_i);
void hatH_exact_r(const Vector & X, Vector & hatH_r);
void hatH_exact_i(const Vector & X, Vector & hatH_i);
real_t hatH_exact_scalar_r(const Vector & X);
real_t hatH_exact_scalar_i(const Vector & X);
void rhs_func_r(const Vector &x, Vector & J_r);
void rhs_func_i(const Vector &x, Vector & J_i);
int dim;
int dimc;
real_t omega;
real_t mu = 1.0;
real_t epsilon = 1.0;
int main(int argc, char *argv[])
{
const char *mesh_file = "../../data/inline-quad.mesh";
int order = 1;
int delta_order = 1;
bool visualization = true;
real_t rnum=1.0;
int ref = 0;
bool static_cond = false;
OptionsParser args(argc, argv);
args.AddOption(&mesh_file, "-m", "--mesh",
"Mesh file to use.");
args.AddOption(&order, "-o", "--order",
"Finite element order (polynomial degree)");
args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
"--no-visualization",
"Enable or disable GLVis visualization.");
args.AddOption(&rnum, "-rnum", "--number-of-wavelengths",
"Number of wavelengths");
args.AddOption(&mu, "-mu", "--permeability",
"Permeability of free space (or 1/(spring constant)).");
args.AddOption(&epsilon, "-eps", "--permittivity",
"Permittivity of free space (or mass constant).");
args.AddOption(&delta_order, "-do", "--delta-order",
"Order enrichment for DPG test space.");
args.AddOption(&ref, "-ref", "--serial-ref",
"Number of serial refinements.");
args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
"--no-static-condensation", "Enable static condensation.");
args.Parse();
if (!args.Good())
{
args.PrintUsage(cout);
return 1;
}
args.PrintOptions(cout);
omega = 2.*M_PI*rnum;
Mesh mesh(mesh_file, 1, 1);
dim = mesh.Dimension();
MFEM_VERIFY(dim > 1, "Dimension = 1 is not supported in this example");
dimc = (dim == 3) ? 3 : 1;
int test_order = order+delta_order;
// Define spaces
enum TrialSpace
{
E_space = 0,
H_space = 1,
hatE_space = 2,
hatH_space = 3
};
enum TestSpace
{
F_space = 0,
G_space = 1
};
// L2 space for E
FiniteElementCollection *E_fec = new L2_FECollection(order-1,dim);
FiniteElementSpace *E_fes = new FiniteElementSpace(&mesh,E_fec,dim);
// Vector L2 space for H
FiniteElementCollection *H_fec = new L2_FECollection(order-1,dim);
FiniteElementSpace *H_fes = new FiniteElementSpace(&mesh,H_fec, dimc);
// H^-1/2 (curl) space for Ê
FiniteElementCollection * hatE_fec = nullptr;
FiniteElementCollection * hatH_fec = nullptr;
FiniteElementCollection * F_fec = nullptr;
if (dim == 3)
{
hatE_fec = new ND_Trace_FECollection(order,dim);
hatH_fec = new ND_Trace_FECollection(order,dim);
F_fec = new ND_FECollection(test_order, dim);
}
else
{
hatE_fec = new RT_Trace_FECollection(order-1,dim);
hatH_fec = new H1_Trace_FECollection(order,dim);
F_fec = new H1_FECollection(test_order, dim);
}
FiniteElementSpace *hatE_fes = new FiniteElementSpace(&mesh,hatE_fec);
FiniteElementSpace *hatH_fes = new FiniteElementSpace(&mesh,hatH_fec);
FiniteElementCollection * G_fec = new ND_FECollection(test_order, dim);
// Coefficients
ConstantCoefficient one(1.0);
ConstantCoefficient eps2omeg2(epsilon*epsilon*omega*omega);
ConstantCoefficient mu2omeg2(mu*mu*omega*omega);
ConstantCoefficient muomeg(mu*omega);
ConstantCoefficient negepsomeg(-epsilon*omega);
ConstantCoefficient epsomeg(epsilon*omega);
ConstantCoefficient negmuomeg(-mu*omega);
DenseMatrix rot_mat(2);
rot_mat(0,0) = 0.; rot_mat(0,1) = 1.;
rot_mat(1,0) = -1.; rot_mat(1,1) = 0.;
MatrixConstantCoefficient rot(rot_mat);
ScalarMatrixProductCoefficient epsrot(epsomeg,rot);
ScalarMatrixProductCoefficient negepsrot(negepsomeg,rot);
Array<FiniteElementSpace * > trial_fes;
Array<FiniteElementCollection * > test_fec;
trial_fes.Append(E_fes);
trial_fes.Append(H_fes);
trial_fes.Append(hatE_fes);
trial_fes.Append(hatH_fes);
test_fec.Append(F_fec);
test_fec.Append(G_fec);
ComplexDPGWeakForm * a = new ComplexDPGWeakForm(trial_fes,test_fec);
a->StoreMatrices();
// (E,∇ × F)
a->AddTrialIntegrator(new TransposeIntegrator(new MixedCurlIntegrator(one)),
nullptr,TrialSpace::E_space,TestSpace::F_space);
// -i ω ϵ (E , G)
a->AddTrialIntegrator(nullptr,
new TransposeIntegrator(new VectorFEMassIntegrator(negepsomeg)),
TrialSpace::E_space,TestSpace::G_space);
// (H,∇ × G)
a->AddTrialIntegrator(new TransposeIntegrator(new MixedCurlIntegrator(one)),
nullptr, TrialSpace::H_space,TestSpace::G_space);
if (dim == 3)
{
// i ω μ (H, F)
a->AddTrialIntegrator(nullptr,
new TransposeIntegrator(new VectorFEMassIntegrator(muomeg)),
TrialSpace::H_space,TestSpace::F_space);
// < n×Ê,F>
a->AddTrialIntegrator(new TangentTraceIntegrator,nullptr,
TrialSpace::hatE_space,TestSpace::F_space);
}
else
{
// i ω μ (H, F)
a->AddTrialIntegrator(nullptr,new MixedScalarMassIntegrator(muomeg),
TrialSpace::H_space,TestSpace::F_space);
// <Ê,F>
a->AddTrialIntegrator(new TraceIntegrator,nullptr,
TrialSpace::hatE_space, TestSpace::F_space);
}
// < n×Ĥ ,G>
a->AddTrialIntegrator(new TangentTraceIntegrator,nullptr,
TrialSpace::hatH_space, TestSpace::G_space);
// test integrators for the adjoint graph norm on the test space
// (∇×G ,∇× δG)
a->AddTestIntegrator(new CurlCurlIntegrator(one),nullptr,
TestSpace::G_space,TestSpace::G_space);
// (G,δG)
a->AddTestIntegrator(new VectorFEMassIntegrator(one),nullptr,
TestSpace::G_space, TestSpace::G_space);
if (dim == 3)
{
// (∇×F,∇×δF)
a->AddTestIntegrator(new CurlCurlIntegrator(one),nullptr,
TestSpace::F_space, TestSpace::F_space);
// (F,δF)
a->AddTestIntegrator(new VectorFEMassIntegrator(one),nullptr,
TestSpace::F_space,TestSpace::F_space);
// μ^2 ω^2 (F,δF)
a->AddTestIntegrator(new VectorFEMassIntegrator(mu2omeg2),nullptr,
TestSpace::F_space, TestSpace::F_space);
// -i ω μ (F,∇ × δG) = (F, ω μ ∇ × δ G)
a->AddTestIntegrator(nullptr,new MixedVectorWeakCurlIntegrator(negmuomeg),
TestSpace::F_space, TestSpace::G_space);
// -i ω ϵ (∇ × F, δG)
a->AddTestIntegrator(nullptr,new MixedVectorCurlIntegrator(negepsomeg),
TestSpace::F_space, TestSpace::G_space);
// i ω μ (∇ × G,δF)
a->AddTestIntegrator(nullptr,new MixedVectorCurlIntegrator(epsomeg),
TestSpace::G_space, TestSpace::F_space);
// i ω ϵ (G, ∇ × δF )
a->AddTestIntegrator(nullptr,new MixedVectorWeakCurlIntegrator(muomeg),
TestSpace::G_space, TestSpace::F_space);
// ϵ^2 ω^2 (G,δG)
a->AddTestIntegrator(new VectorFEMassIntegrator(eps2omeg2),nullptr,
TestSpace::G_space, TestSpace::G_space);
}
else
{
// (∇F,∇δF)
a->AddTestIntegrator(new DiffusionIntegrator(one),nullptr,
TestSpace::F_space, TestSpace::F_space);
// (F,δF)
a->AddTestIntegrator(new MassIntegrator(one),nullptr,
TestSpace::F_space, TestSpace::F_space);
// μ^2 ω^2 (F,δF)
a->AddTestIntegrator(new MassIntegrator(mu2omeg2),nullptr,
TestSpace::F_space, TestSpace::F_space);
// -i ω μ (F,∇ × δG) = i (F, -ω μ ∇ × δ G)
a->AddTestIntegrator(nullptr,
new TransposeIntegrator(new MixedCurlIntegrator(negmuomeg)),
TestSpace::F_space, TestSpace::G_space);
// -i ω ϵ (∇ × F, δG) = i (- ω ϵ A ∇ F,δG), A = [0 1; -1; 0]
a->AddTestIntegrator(nullptr,new MixedVectorGradientIntegrator(negepsrot),
TestSpace::F_space, TestSpace::G_space);
// i ω μ (∇ × G,δF) = i (ω μ ∇ × G, δF )
a->AddTestIntegrator(nullptr,new MixedCurlIntegrator(muomeg),
TestSpace::G_space, TestSpace::F_space);
// i ω ϵ (G, ∇ × δF ) = i (ω ϵ G, A ∇ δF) = i ( G , ω ϵ A ∇ δF)
a->AddTestIntegrator(nullptr,
new TransposeIntegrator(new MixedVectorGradientIntegrator(epsrot)),
TestSpace::G_space, TestSpace::F_space);
// ϵ^2 ω^2 (G,δG)
a->AddTestIntegrator(new VectorFEMassIntegrator(eps2omeg2),nullptr,
TestSpace::G_space, TestSpace::G_space);
}
// RHS
VectorFunctionCoefficient f_rhs_r(dim,rhs_func_r);
VectorFunctionCoefficient f_rhs_i(dim,rhs_func_i);
a->AddDomainLFIntegrator(new VectorFEDomainLFIntegrator(f_rhs_r),
new VectorFEDomainLFIntegrator(f_rhs_i),
TestSpace::G_space);
VectorFunctionCoefficient hatEex_r(dim,hatE_exact_r);
VectorFunctionCoefficient hatEex_i(dim,hatE_exact_i);
socketstream E_out_r;
socketstream E_out_i;
real_t err0 = 0.;
int dof0 = 0; // init to suppress gcc warning
std::cout << "\n Ref |"
<< " Dofs |"
<< " ω |"
<< " L2 Error |"
<< " Rate |"
<< " PCG it |" << endl;
std::cout << std::string(60,'-')
<< endl;
for (int it = 0; it<=ref; it++)
{
if (static_cond) { a->EnableStaticCondensation(); }
a->Assemble();
Array<int> ess_tdof_list;
Array<int> ess_bdr;
if (mesh.bdr_attributes.Size())
{
ess_bdr.SetSize(mesh.bdr_attributes.Max());
ess_bdr = 1;
hatE_fes->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
}
// shift the ess_tdofs
for (int j = 0; j < ess_tdof_list.Size(); j++)
{
ess_tdof_list[j] += E_fes->GetTrueVSize() + H_fes->GetTrueVSize();
}
Array<int> offsets(5);
offsets[0] = 0;
offsets[1] = E_fes->GetVSize();
offsets[2] = H_fes->GetVSize();
offsets[3] = hatE_fes->GetVSize();
offsets[4] = hatH_fes->GetVSize();
offsets.PartialSum();
Vector x(2*offsets.Last());
x = 0.;
GridFunction hatE_gf_r(hatE_fes, x, offsets[2]);
GridFunction hatE_gf_i(hatE_fes, x, offsets.Last() + offsets[2]);
if (dim == 3)
{
hatE_gf_r.ProjectBdrCoefficientTangent(hatEex_r, ess_bdr);
hatE_gf_i.ProjectBdrCoefficientTangent(hatEex_i, ess_bdr);
}
else
{
hatE_gf_r.ProjectBdrCoefficientNormal(hatEex_r, ess_bdr);
hatE_gf_i.ProjectBdrCoefficientNormal(hatEex_i, ess_bdr);
}
OperatorPtr Ah;
Vector X,B;
a->FormLinearSystem(ess_tdof_list,x,Ah, X,B);
ComplexOperator * Ahc = Ah.As<ComplexOperator>();
BlockMatrix * A_r = dynamic_cast<BlockMatrix *>(&Ahc->real());
BlockMatrix * A_i = dynamic_cast<BlockMatrix *>(&Ahc->imag());
int num_blocks = A_r->NumRowBlocks();
Array<int> tdof_offsets(2*num_blocks+1);
tdof_offsets[0] = 0;
int k = (static_cond) ? 2 : 0;
for (int i=0; i<num_blocks; i++)
{
tdof_offsets[i+1] = trial_fes[i+k]->GetTrueVSize();
tdof_offsets[num_blocks+i+1] = trial_fes[i+k]->GetTrueVSize();
}
tdof_offsets.PartialSum();
BlockOperator A(tdof_offsets);
for (int i = 0; i<num_blocks; i++)
{
for (int j = 0; j<num_blocks; j++)
{
A.SetBlock(i,j,&A_r->GetBlock(i,j));
A.SetBlock(i,j+num_blocks,&A_i->GetBlock(i,j), -1.0);
A.SetBlock(i+num_blocks,j+num_blocks,&A_r->GetBlock(i,j));
A.SetBlock(i+num_blocks,j,&A_i->GetBlock(i,j));
}
}
BlockDiagonalPreconditioner M(tdof_offsets);
M.owns_blocks = 1;
for (int i = 0; i<num_blocks; i++)
{
M.SetDiagonalBlock(i, new GSSmoother((SparseMatrix&)A_r->GetBlock(i,i)));
M.SetDiagonalBlock(num_blocks+i, new GSSmoother((SparseMatrix&)A_r->GetBlock(i,
i)));
}
CGSolver cg;
cg.SetRelTol(1e-10);
cg.SetMaxIter(2000);
cg.SetPrintLevel(0);
cg.SetPreconditioner(M);
cg.SetOperator(A);
cg.Mult(B, X);
a->RecoverFEMSolution(X,x);
GridFunction E_r(E_fes, x, 0);
GridFunction E_i(E_fes, x, offsets.Last());
VectorFunctionCoefficient E_ex_r(dim,E_exact_r);
VectorFunctionCoefficient E_ex_i(dim,E_exact_i);
GridFunction H_r(H_fes, x, offsets[1]);
GridFunction H_i(H_fes, x, offsets.Last() + offsets[1]);
VectorFunctionCoefficient H_ex_r(dimc,H_exact_r);
VectorFunctionCoefficient H_ex_i(dimc,H_exact_i);
int dofs = 0;
for (int i = 0; i<trial_fes.Size(); i++)
{
dofs += trial_fes[i]->GetTrueVSize();
}
real_t E_err_r = E_r.ComputeL2Error(E_ex_r);
real_t E_err_i = E_i.ComputeL2Error(E_ex_i);
real_t H_err_r = H_r.ComputeL2Error(H_ex_r);
real_t H_err_i = H_i.ComputeL2Error(H_ex_i);
real_t L2Error = sqrt(E_err_r*E_err_r + E_err_i*E_err_i
+ H_err_r*H_err_r + H_err_i*H_err_i);
real_t rate_err = (it) ? dim*log(err0/L2Error)/log((real_t)dof0/dofs) : 0.0;
err0 = L2Error;
dof0 = dofs;
std::ios oldState(nullptr);
oldState.copyfmt(std::cout);
std::cout << std::right << std::setw(5) << it << " | "
<< std::setw(10) << dof0 << " | "
<< std::setprecision(1) << std::fixed
<< std::setw(4) << 2*rnum << " π | "
<< std::setprecision(3)
<< std::setw(10) << std::scientific << err0 << " | "
<< std::setprecision(2)
<< std::setw(6) << std::fixed << rate_err << " | "
<< std::setw(6) << std::fixed << cg.GetNumIterations() << " | "
<< std::endl;
std::cout.copyfmt(oldState);
if (visualization)
{
const char * keys = (it == 0 && dim == 2) ? "jRcml\n" : nullptr;
char vishost[] = "localhost";
int visport = 19916;
VisualizeField(E_out_r,vishost, visport, E_r,
"Numerical Electric field (real part)", 0, 0, 500, 500, keys);
VisualizeField(E_out_i,vishost, visport, E_i,
"Numerical Electric field (imaginary part)", 501, 0, 500, 500, keys);
}
if (it == ref)
{
break;
}
mesh.UniformRefinement();
for (int i =0; i<trial_fes.Size(); i++)
{
trial_fes[i]->Update(false);
}
a->Update();
}
delete a;
delete F_fec;
delete G_fec;
delete hatH_fes;
delete hatH_fec;
delete hatE_fes;
delete hatE_fec;
delete H_fec;
delete E_fec;
delete H_fes;
delete E_fes;
return 0;
}
void E_exact_r(const Vector &x, Vector & E_r)
{
std::vector<std::complex<real_t>> E;
maxwell_solution(x, E);
E_r.SetSize(E.size());
for (unsigned i = 0; i < E.size(); i++)
{
E_r[i]= E[i].real();
}
}
void E_exact_i(const Vector &x, Vector & E_i)
{
std::vector<std::complex<real_t>> E;
maxwell_solution(x, E);
E_i.SetSize(E.size());
for (unsigned i = 0; i < E.size(); i++)
{
E_i[i]= E[i].imag();
}
}
void curlE_exact_r(const Vector &x, Vector &curlE_r)
{
std::vector<std::complex<real_t>> curlE;
maxwell_solution_curl(x, curlE);
curlE_r.SetSize(curlE.size());
for (unsigned i = 0; i < curlE.size(); i++)
{
curlE_r[i]= curlE[i].real();
}
}
void curlE_exact_i(const Vector &x, Vector &curlE_i)
{
std::vector<std::complex<real_t>> curlE;
maxwell_solution_curl(x, curlE);
curlE_i.SetSize(curlE.size());
for (unsigned i = 0; i < curlE.size(); i++)
{
curlE_i[i]= curlE[i].imag();
}
}
void curlcurlE_exact_r(const Vector &x, Vector & curlcurlE_r)
{
std::vector<std::complex<real_t>> curlcurlE;
maxwell_solution_curlcurl(x, curlcurlE);
curlcurlE_r.SetSize(curlcurlE.size());
for (unsigned i = 0; i < curlcurlE.size(); i++)
{
curlcurlE_r[i]= curlcurlE[i].real();
}
}
void curlcurlE_exact_i(const Vector &x, Vector & curlcurlE_i)
{
std::vector<std::complex<real_t>> curlcurlE;
maxwell_solution_curlcurl(x, curlcurlE);
curlcurlE_i.SetSize(curlcurlE.size());
for (unsigned i = 0; i < curlcurlE.size(); i++)
{
curlcurlE_i[i]= curlcurlE[i].imag();
}
}
void H_exact_r(const Vector &x, Vector & H_r)
{
// H = i ∇ × E / ω μ
// H_r = - ∇ × E_i / ω μ
Vector curlE_i;
curlE_exact_i(x,curlE_i);
H_r.SetSize(dimc);
for (int i = 0; i<dimc; i++)
{
H_r(i) = - curlE_i(i) / (omega * mu);
}
}
void H_exact_i(const Vector &x, Vector & H_i)
{
// H = i ∇ × E / ω μ
// H_i = ∇ × E_r / ω μ
Vector curlE_r;
curlE_exact_r(x,curlE_r);
H_i.SetSize(dimc);
for (int i = 0; i<dimc; i++)
{
H_i(i) = curlE_r(i) / (omega * mu);
}
}
void curlH_exact_r(const Vector &x,Vector &curlH_r)
{
// ∇ × H_r = - ∇ × ∇ × E_i / ω μ
Vector curlcurlE_i;
curlcurlE_exact_i(x,curlcurlE_i);
curlH_r.SetSize(dim);
for (int i = 0; i<dim; i++)
{
curlH_r(i) = -curlcurlE_i(i) / (omega * mu);
}
}
void curlH_exact_i(const Vector &x,Vector &curlH_i)
{
// ∇ × H_i = ∇ × ∇ × E_r / ω μ
Vector curlcurlE_r;
curlcurlE_exact_r(x,curlcurlE_r);
curlH_i.SetSize(dim);
for (int i = 0; i<dim; i++)
{
curlH_i(i) = curlcurlE_r(i) / (omega * mu);
}
}
void hatE_exact_r(const Vector & x, Vector & hatE_r)
{
if (dim == 3)
{
E_exact_r(x,hatE_r);
}
else
{
Vector E_r;
E_exact_r(x,E_r);
hatE_r.SetSize(hatE_r.Size());
// rotate E_hat
hatE_r[0] = E_r[1];
hatE_r[1] = -E_r[0];
}
}
void hatE_exact_i(const Vector & x, Vector & hatE_i)
{
if (dim == 3)
{
E_exact_i(x,hatE_i);
}
else
{
Vector E_i;
E_exact_i(x,E_i);
hatE_i.SetSize(hatE_i.Size());
// rotate E_hat
hatE_i[0] = E_i[1];
hatE_i[1] = -E_i[0];
}
}
void hatH_exact_r(const Vector & x, Vector & hatH_r)
{
H_exact_r(x,hatH_r);
}
void hatH_exact_i(const Vector & x, Vector & hatH_i)
{
H_exact_i(x,hatH_i);
}
real_t hatH_exact_scalar_r(const Vector & x)
{
Vector hatH_r;
H_exact_r(x,hatH_r);
return hatH_r[0];
}
real_t hatH_exact_scalar_i(const Vector & x)
{
Vector hatH_i;
H_exact_i(x,hatH_i);
return hatH_i[0];
}
// J = -i ω ϵ E + ∇ × H
// J_r + iJ_i = -i ω ϵ (E_r + i E_i) + ∇ × (H_r + i H_i)
void rhs_func_r(const Vector &x, Vector & J_r)
{
// J_r = ω ϵ E_i + ∇ × H_r
Vector E_i, curlH_r;
E_exact_i(x,E_i);
curlH_exact_r(x,curlH_r);
J_r.SetSize(dim);
for (int i = 0; i<dim; i++)
{
J_r(i) = omega * epsilon * E_i(i) + curlH_r(i);
}
}
void rhs_func_i(const Vector &x, Vector & J_i)
{
// J_i = - ω ϵ E_r + ∇ × H_i
Vector E_r, curlH_i;
E_exact_r(x,E_r);
curlH_exact_i(x,curlH_i);
J_i.SetSize(dim);
for (int i = 0; i<dim; i++)
{
J_i(i) = -omega * epsilon * E_r(i) + curlH_i(i);
}
}
void maxwell_solution(const Vector & X, std::vector<complex<real_t>> &E)
{
E.resize(dim);
std::complex<real_t> zi(0,1);
std::complex<real_t> pw = exp(-zi * omega * (X.Sum()));
E[0] = pw;
E[1] = 0.0;
if (dim == 3) { E[2] = 0.0; }
}
void maxwell_solution_curl(const Vector & X,
std::vector<complex<real_t>> &curlE)
{
curlE.resize(dimc);
std::complex<real_t> zi(0,1);
std::complex<real_t> pw = exp(-zi * omega * (X.Sum()));
if (dim == 3)
{
curlE[0] = 0.0;
curlE[1] = -zi * omega * pw;
curlE[2] = zi * omega * pw;
}
else
{
curlE[0] = zi * omega * pw;
}
}
void maxwell_solution_curlcurl(const Vector & X,
std::vector<complex<real_t>> &curlcurlE)
{
curlcurlE.resize(dim);
std::complex<real_t> zi(0,1);
std::complex<real_t> pw = exp(-zi * omega * (X.Sum()));
if (dim == 3)
{
curlcurlE[0] = 2_r * omega * omega * pw;
curlcurlE[1] = - omega * omega * pw;
curlcurlE[2] = - omega * omega * pw;
}
else
{
curlcurlE[0] = omega * omega * pw;
curlcurlE[1] = - omega * omega * pw ;
}
}