-
Notifications
You must be signed in to change notification settings - Fork 9
/
remhos_tools.cpp
1533 lines (1375 loc) · 46.2 KB
/
remhos_tools.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
// Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at
// the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights
// reserved. See files LICENSE and NOTICE for details.
//
// This file is part of CEED, a collection of benchmarks, miniapps, software
// libraries and APIs for efficient high-order finite element and spectral
// element discretizations for exascale applications. For more information and
// source code availability see http://github.com/ceed.
//
// The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
// a collaborative effort of two U.S. Department of Energy organizations (Office
// of Science and the National Nuclear Security Administration) responsible for
// the planning and preparation of a capable exascale ecosystem, including
// software, applications, hardware, advanced system engineering and early
// testbed platforms, in support of the nation's exascale computing imperative.
#include "remhos_tools.hpp"
using namespace std;
namespace mfem
{
SmoothnessIndicator::SmoothnessIndicator(int type_id,
ParMesh &subcell_mesh,
ParFiniteElementSpace &pfes_DG_,
ParGridFunction &u,
DofInfo &dof_info)
: type(type_id), param(type == 1 ? 5.0 : 3.0),
fec_sub(1, pfes_DG_.GetMesh()->Dimension(), BasisType::Positive),
pfes_CG_sub(&subcell_mesh, &fec_sub),
pfes_DG(pfes_DG_)
{
// TODO assemble SI matrices every RK stage for remap.
MFEM_VERIFY(type_id == 1 || type_id == 2, "Bad smoothness indicator id!");
BilinearForm massH1(&pfes_CG_sub);
massH1.AddDomainIntegrator(new MassIntegrator);
massH1.Assemble();
massH1.Finalize();
Mmat = massH1.SpMat();
ConstantCoefficient neg_one(-1.0);
BilinearForm lap(&pfes_CG_sub);
lap.AddDomainIntegrator(new DiffusionIntegrator(neg_one));
lap.AddBdrFaceIntegrator(new DGDiffusionIntegrator(neg_one, 0., 0.));
lap.Assemble();
lap.Finalize();
LaplaceOp = lap.SpMat();
MassMixed = new SparseMatrix(pfes_CG_sub.GetVSize(), pfes_DG.GetVSize());
MassInt = new MassIntegrator;
BilinearForm mlH1(&pfes_CG_sub);
mlH1.AddDomainIntegrator(new LumpedIntegrator(new MassIntegrator));
mlH1.Assemble();
mlH1.Finalize();
mlH1.SpMat().GetDiag(lumpedMH1);
Vector lumped_hv(pfes_CG_sub.GetTrueVSize());
pfes_CG_sub.Dof_TrueDof_Matrix()->MultTranspose(lumpedMH1, lumped_hv);
pfes_CG_sub.GetProlongationMatrix()->Mult(lumped_hv, lumpedMH1);
// Stores the index for the dof of H1-conforming for each node.
// If the node is on the boundary, the entry is -1.
const int dim = pfes_DG.GetMesh()->Dimension(),
nd = pfes_DG.GetFE(0)->GetDof(),
ne = pfes_DG.GetMesh()->GetNE();
DG2CG.SetSize(ne*nd);
Array<int> vdofs;
int e, i, j, e_id;
for (e = 0; e < ne; e++)
{
for (i = 0; i < dof_info.numSubcells; i++)
{
e_id = e*dof_info.numSubcells + i;
pfes_CG_sub.GetElementVDofs(e_id, vdofs);
// Switchero - numbering CG vs DG.
DG2CG(e*nd + dof_info.Sub2Ind(i, 0)) = vdofs[0];
DG2CG(e*nd + dof_info.Sub2Ind(i, 1)) = vdofs[1];
if (dim > 1)
{
DG2CG(e*nd + dof_info.Sub2Ind(i, 2)) = vdofs[3];
DG2CG(e*nd + dof_info.Sub2Ind(i, 3)) = vdofs[2];
}
if (dim == 3)
{
DG2CG(e*nd + dof_info.Sub2Ind(i, 4)) = vdofs[4];
DG2CG(e*nd + dof_info.Sub2Ind(i, 5)) = vdofs[5];
DG2CG(e*nd + dof_info.Sub2Ind(i, 6)) = vdofs[7];
DG2CG(e*nd + dof_info.Sub2Ind(i, 7)) = vdofs[6];
}
}
// Domain boundaries.
for (i = 0; i < dof_info.numBdrs; i++)
{
for (j = 0; j < dof_info.numFaceDofs; j++)
{
if (dof_info.NbrDof(e,i,j) < 0)
{
DG2CG(e*nd+dof_info.BdrDofs(j,i)) = -1;
}
}
}
}
ComputeVariationalMatrix(dof_info);
IntegrationRule *ir;
IntegrationRule irX;
QuadratureFunctions1D qf;
qf.ClosedUniform(pfes_DG.GetFE(0)->GetOrder() + 1, &irX);
Vector shape(nd);
if (dim == 1) { ir = &irX; }
if (dim == 2) { ir = new IntegrationRule(irX, irX); }
if (dim == 3) { ir = new IntegrationRule(irX, irX, irX); }
ShapeEval.SetSize(nd, nd); // nd equals ir->GetNPoints().
for (i = 0; i < ir->GetNPoints(); i++)
{
const IntegrationPoint &ip = ir->IntPoint(i);
pfes_DG.GetFE(0)->CalcShape(ip, shape);
ShapeEval.SetCol(i, shape);
}
ShapeEval.Transpose();
// Print the values of the smoothness indicator.
ParGridFunction si_val;
ComputeSmoothnessIndicator(u, si_val);
{
ofstream smth("si_init.gf");
smth.precision(8);
si_val.SaveAsOne(smth);
}
if (dim > 1) { delete ir; }
}
SmoothnessIndicator::~SmoothnessIndicator()
{
delete MassInt;
delete MassMixed;
}
void SmoothnessIndicator::ComputeSmoothnessIndicator(const Vector &u,
ParGridFunction &si_vals_u)
{
si_vals_u.SetSpace(&pfes_CG_sub);
ParGridFunction g(&pfes_CG_sub);
const int N = g.Size();
Vector gmin(N), gmax(N);
ApproximateLaplacian(u, g);
ComputeFromSparsity(Mmat, g, gmin, gmax);
if (type == 1)
{
const double eps = 1.0e-50;
for (int e = 0; e < N; e++)
{
si_vals_u(e) = 1.0 - pow( (abs(gmin(e) - gmax(e)) + eps) /
(abs(gmin(e)) + abs(gmax(e)) + eps), param);
}
}
else if (type == 2)
{
const double eps = 1.0e-15;
for (int e = 0; e < N; e++)
{
si_vals_u(e) = min(1.0, param *
max(0., gmin(e)*gmax(e)) /
(max(gmin(e)*gmin(e),gmax(e)*gmax(e)) + eps) );
}
}
}
void SmoothnessIndicator::UpdateBounds(int dof_id, double u_HO,
const ParGridFunction &si_vals,
double &u_min, double &u_max)
{
const double tmp = (DG2CG(dof_id) < 0.0) ? 1.0 : si_vals(DG2CG(dof_id));
u_min = max(0.0, tmp * u_HO + (1.0 - tmp) * u_min);
u_max = min(1.0, tmp * u_HO + (1.0 - tmp) * u_max);
}
void SmoothnessIndicator::ComputeVariationalMatrix(DofInfo &dof_info)
{
Mesh *subcell_mesh = pfes_CG_sub.GetMesh();
const int dim = subcell_mesh->Dimension();
const int nd = pfes_DG.GetFE(0)->GetDof(), ne = pfes_DG.GetMesh()->GetNE();
int k, m, e_id;
DenseMatrix elmat1;
Array <int> te_vdofs, tr_vdofs;
tr_vdofs.SetSize(dof_info.numDofsSubcell);
for (k = 0; k < ne; k++)
{
for (m = 0; m < dof_info.numSubcells; m++)
{
e_id = k*dof_info.numSubcells + m;
const FiniteElement *el = pfes_CG_sub.GetFE(e_id);
ElementTransformation *tr = subcell_mesh->GetElementTransformation(e_id);
MassInt->AssembleElementMatrix(*el, *tr, elmat1);
pfes_CG_sub.GetElementVDofs(e_id, te_vdofs);
// Switchero - numbering CG vs DG.
tr_vdofs[0] = k*nd + dof_info.Sub2Ind(m, 0);
tr_vdofs[1] = k*nd + dof_info.Sub2Ind(m, 1);
if (dim > 1)
{
tr_vdofs[2] = k*nd + dof_info.Sub2Ind(m, 3);
tr_vdofs[3] = k*nd + dof_info.Sub2Ind(m, 2);
}
if (dim == 3)
{
tr_vdofs[4] = k*nd + dof_info.Sub2Ind(m, 4);
tr_vdofs[5] = k*nd + dof_info.Sub2Ind(m, 5);
tr_vdofs[6] = k*nd + dof_info.Sub2Ind(m, 7);
tr_vdofs[7] = k*nd + dof_info.Sub2Ind(m, 6);
}
MassMixed->AddSubMatrix(te_vdofs, tr_vdofs, elmat1);
}
}
MassMixed->Finalize();
}
void SmoothnessIndicator::ApproximateLaplacian(const Vector &x,
ParGridFunction &y)
{
const int nd = pfes_DG.GetFE(0)->GetDof(), ne = pfes_DG.GetMesh()->GetNE();
int k, i, j, N = lumpedMH1.Size();
Array<int> eldofs;
Vector xDofs(nd), tmp(nd), xEval(ne*nd);
Vector rhs_tv(pfes_CG_sub.GetTrueVSize()), z_tv(pfes_CG_sub.GetTrueVSize());
eldofs.SetSize(nd);
y.SetSize(N);
Vector rhs(N), z(N);
// Approximate inversion, corresponding to Neumann series truncated after
// first two summands.
int iter, max_iter = 2;
const double abs_tol = 1.e-10;
double resid;
for (k = 0; k < ne; k++)
{
for (j = 0; j < nd; j++) { eldofs[j] = k*nd + j; }
x.GetSubVector(eldofs, xDofs);
ShapeEval.Mult(xDofs, tmp);
xEval.SetSubVector(eldofs, tmp);
}
MassMixed->Mult(xEval, rhs);
pfes_CG_sub.Dof_TrueDof_Matrix()->MultTranspose(rhs, rhs_tv);
pfes_CG_sub.GetProlongationMatrix()->Mult(rhs_tv, rhs);
y = 0.;
// Project x to a CG space (result is in y).
for (iter = 1; iter <= max_iter; iter++)
{
Mmat.Mult(y, z);
pfes_CG_sub.Dof_TrueDof_Matrix()->MultTranspose(z, z_tv);
z_tv -= rhs_tv;
pfes_CG_sub.GetProlongationMatrix()->Mult(z_tv, z);
double loc_res = z_tv.Norml2();
loc_res *= loc_res;
MPI_Allreduce(&loc_res, &resid, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
resid = sqrt(resid);
if (resid <= abs_tol) { break; }
y.HostReadWrite();
z.HostReadWrite();
const double *h_lumpedMH1 = lumpedMH1.HostRead();
for (i = 0; i < N; i++)
{
y(i) -= z(i) / h_lumpedMH1[i];
}
}
LaplaceOp.Mult(y, rhs);
pfes_CG_sub.Dof_TrueDof_Matrix()->MultTranspose(rhs, rhs_tv);
pfes_CG_sub.GetProlongationMatrix()->Mult(rhs_tv, rhs);
y = 0.;
for (iter = 1; iter <= max_iter; iter++)
{
Mmat.Mult(y, z);
pfes_CG_sub.Dof_TrueDof_Matrix()->MultTranspose(z, z_tv);
z_tv -= rhs_tv;
pfes_CG_sub.GetProlongationMatrix()->Mult(z_tv, z);
double loc_res = z_tv.Norml2();
loc_res *= loc_res;
MPI_Allreduce(&loc_res, &resid, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
resid = sqrt(resid);
if (resid <= abs_tol) { break; }
y.HostReadWrite();
z.HostReadWrite();
const double * h_lumpedMH1 = lumpedMH1.HostRead();
for (i = 0; i < N; i++)
{
y(i) -= z(i) / h_lumpedMH1[i];
}
}
}
void SmoothnessIndicator::ComputeFromSparsity(const SparseMatrix &K,
const ParGridFunction &x,
Vector &x_min, Vector &x_max)
{
const int *I = K.GetI(), *J = K.GetJ(), loc_size = K.Size();
int end;
for (int i = 0, k = 0; i < loc_size; i++)
{
x_min(i) = numeric_limits<double>::infinity();
x_max(i) = -x_min(i);
for (end = I[i+1]; k < end; k++)
{
const double x_j = x(J[k]);
x_max(i) = max(x_max(i), x_j);
x_min(i) = min(x_min(i), x_j);
}
}
GroupCommunicator &gcomm = x.ParFESpace()->GroupComm();
Array<double> minvals(x_min.GetData(), x_min.Size()),
maxvals(x_max.GetData(), x_max.Size());
gcomm.Reduce<double>(minvals, GroupCommunicator::Min);
gcomm.Bcast(minvals);
gcomm.Reduce<double>(maxvals, GroupCommunicator::Max);
gcomm.Bcast(maxvals);
}
DofInfo::DofInfo(ParFiniteElementSpace &pfes_sltn, int btype)
: bounds_type(btype),
pmesh(pfes_sltn.GetParMesh()), pfes(pfes_sltn),
fec_bounds(std::max(pfes.GetOrder(0), 1),
pmesh->Dimension(), BasisType::GaussLobatto),
pfes_bounds(pmesh, &fec_bounds),
x_min(&pfes_bounds), x_max(&pfes_bounds)
{
int n = pfes.GetVSize();
int ne = pmesh->GetNE();
xi_min.SetSize(n);
xi_max.SetSize(n);
xe_min.SetSize(ne);
xe_max.SetSize(ne);
ExtractBdrDofs(pfes.GetOrder(0),
pfes.GetFE(0)->GetGeomType(), BdrDofs);
numFaceDofs = BdrDofs.Height();
numBdrs = BdrDofs.Width();
FillNeighborDofs(); // Fill NbrDof.
FillSubcell2CellDof(); // Fill Sub2Ind.
}
void DofInfo::ComputeMatrixSparsityBounds(const Vector &el_min,
const Vector &el_max,
Vector &dof_min, Vector &dof_max,
Array<bool> *active_el)
{
ParMesh *pmesh = pfes.GetParMesh();
L2_FECollection fec_bounds(0, pmesh->Dimension());
ParFiniteElementSpace pfes_bounds(pmesh, &fec_bounds);
ParGridFunction x_min(&pfes_bounds), x_max(&pfes_bounds);
const int NE = pmesh->GetNE();
const int ndofs = dof_min.Size() / NE;
x_min.HostReadWrite();
x_max.HostReadWrite();
x_min = el_min;
x_max = el_max;
x_min.ExchangeFaceNbrData(); x_max.ExchangeFaceNbrData();
const Vector &minv = x_min.FaceNbrData(), &maxv = x_max.FaceNbrData();
const Table &el_to_el = pmesh->ElementToElementTable();
Array<int> face_nbr_el;
for (int i = 0; i < NE; i++)
{
double el_min = x_min(i), el_max = x_max(i);
el_to_el.GetRow(i, face_nbr_el);
for (int n = 0; n < face_nbr_el.Size(); n++)
{
if (face_nbr_el[n] < NE)
{
// Local neighbor.
el_min = std::min(el_min, x_min(face_nbr_el[n]));
el_max = std::max(el_max, x_max(face_nbr_el[n]));
}
else
{
// MPI face neighbor.
el_min = std::min(el_min, minv(face_nbr_el[n] - NE));
el_max = std::max(el_max, maxv(face_nbr_el[n] - NE));
}
}
for (int j = 0; j < ndofs; j++)
{
dof_min(i*ndofs + j) = el_min;
dof_max(i*ndofs + j) = el_max;
}
}
}
void DofInfo::ComputeOverlapBounds(const Vector &el_min,
const Vector &el_max,
Vector &dof_min, Vector &dof_max,
Array<bool> *active_el)
{
GroupCommunicator &gcomm = pfes_bounds.GroupComm();
Array<int> dofsCG;
const int NE = pfes.GetNE();
// Form min/max at each CG dof, considering element overlaps.
x_min = std::numeric_limits<double>::infinity();
x_max = - std::numeric_limits<double>::infinity();
for (int i = 0; i < NE; i++)
{
// Inactive elements don't affect the bounds.
if (active_el && (*active_el)[i] == false) { continue; }
x_min.HostReadWrite();
x_max.HostReadWrite();
pfes_bounds.GetElementDofs(i, dofsCG);
for (int j = 0; j < dofsCG.Size(); j++)
{
x_min(dofsCG[j]) = std::min(x_min(dofsCG[j]), el_min(i));
x_max(dofsCG[j]) = std::max(x_max(dofsCG[j]), el_max(i));
}
}
Array<double> minvals(x_min.GetData(), x_min.Size());
Array<double> maxvals(x_max.GetData(), x_max.Size());
gcomm.Reduce<double>(minvals, GroupCommunicator::Min);
gcomm.Bcast(minvals);
gcomm.Reduce<double>(maxvals, GroupCommunicator::Max);
gcomm.Bcast(maxvals);
// Use (x_min, x_max) to fill (dof_min, dof_max) for each DG dof.
const TensorBasisElement *fe_cg =
dynamic_cast<const TensorBasisElement *>(pfes_bounds.GetFE(0));
const Array<int> &dof_map = fe_cg->GetDofMap();
const int ndofs = dof_map.Size();
for (int i = 0; i < NE; i++)
{
// Comment about the case when active_el != null, i.e., when this function
// is used to compute the bounds of s:
//
// Note that this loop goes over all elements, including inactive ones.
// The following happens in an inactive element:
// - If a DOF is on the boundary with an active element, it will get the
// value that's propagated by the continuous functions x_min and x_max.
// - Otherwise, the DOF would get the inf values.
// This is the mechanism that allows new elements, that switch from
// inactive to active, to get some valid bounds. More specifically, this
// function is called on the old state, but the result from it is used
// to limit the new state, which has different active elements.
pfes_bounds.GetElementDofs(i, dofsCG);
for (int j = 0; j < dofsCG.Size(); j++)
{
dof_min(i*ndofs + j) = x_min(dofsCG[dof_map[j]]);
dof_max(i*ndofs + j) = x_max(dofsCG[dof_map[j]]);
}
}
}
void DofInfo::ComputeElementsMinMax(const Vector &u,
Vector &u_min, Vector &u_max,
Array<bool> *active_el,
Array<bool> *active_dof) const
{
const int NE = pfes.GetNE(), ndof = pfes.GetFE(0)->GetDof();
int dof_id;
u.HostRead(); u_min.HostReadWrite(); u_max.HostReadWrite();
for (int k = 0; k < NE; k++)
{
u_min(k) = numeric_limits<double>::infinity();
u_max(k) = -numeric_limits<double>::infinity();
// Inactive elements don't affect the bounds.
if (active_el && (*active_el)[k] == false) { continue; }
for (int i = 0; i < ndof; i++)
{
dof_id = k*ndof + i;
// Inactive dofs don't affect the bounds.
if (active_dof && (*active_dof)[dof_id] == false) { continue; }
u_min(k) = min(u_min(k), u(dof_id));
u_max(k) = max(u_max(k), u(dof_id));
}
}
}
void DofInfo::FillNeighborDofs()
{
// Use the first mesh element as indicator.
const FiniteElement &dummy = *pfes.GetFE(0);
const int dim = pmesh->Dimension();
int i, j, k, ne = pmesh->GetNE();
int nd = dummy.GetDof(), p = dummy.GetOrder();
Array <int> bdrs, orientation;
pmesh->ExchangeFaceNbrData();
Table *face_to_el = pmesh->GetFaceToAllElementTable();
NbrDof.SetSize(ne, numBdrs, numFaceDofs);
// Permutations of BdrDofs, taking into account all possible orientations.
// Assumes BdrDofs are ordered in xyz order, which is true for 3D hexes,
// but it isn't true for 2D quads.
// TODO: check other FEs, function ExtractBoundaryDofs().
int orient_cnt = 1;
if (dim == 2) { orient_cnt = 2; }
if (dim == 3) { orient_cnt = 8; }
const int dof1D_cnt = p+1;
DenseTensor fdof_ids(numFaceDofs, numBdrs, orient_cnt);
for (int ori = 0; ori < orient_cnt; ori++)
{
for (int face_id = 0; face_id < numBdrs; face_id++)
{
for (int fdof_id = 0; fdof_id < numFaceDofs; fdof_id++)
{
// Index of fdof_id in the current orientation.
const int ori_fdof_id = GetLocalFaceDofIndex(dim, face_id, ori,
fdof_id, dof1D_cnt);
fdof_ids(ori)(ori_fdof_id, face_id) = BdrDofs(fdof_id, face_id);
}
}
}
for (k = 0; k < ne; k++)
{
if (dim == 1)
{
pmesh->GetElementVertices(k, bdrs);
for (i = 0; i < numBdrs; i++)
{
const int nbr_cnt = face_to_el->RowSize(bdrs[i]);
if (nbr_cnt == 1)
{
// No neighbor element.
NbrDof(k, i, 0) = -1;
continue;
}
int el1_id, el2_id, nbr_id;
pmesh->GetFaceElements(bdrs[i], &el1_id, &el2_id);
if (el2_id < 0)
{
// This element is in a different mpi task.
el2_id = -1 - el2_id + ne;
}
nbr_id = (el1_id == k) ? el2_id : el1_id;
NbrDof(k,i,0) = nbr_id*nd + BdrDofs(0, (i+1) % 2);
}
}
else if (dim==2)
{
pmesh->GetElementEdges(k, bdrs, orientation);
for (i = 0; i < numBdrs; i++)
{
const int nbr_cnt = face_to_el->RowSize(bdrs[i]);
if (nbr_cnt == 1)
{
// No neighbor element.
for (j = 0; j < numFaceDofs; j++) { NbrDof(k, i, j) = -1; }
continue;
}
int el1_id, el2_id, nbr_id;
pmesh->GetFaceElements(bdrs[i], &el1_id, &el2_id);
if (el2_id < 0)
{
// This element is in a different mpi task.
el2_id = -1 - el2_id + ne;
}
nbr_id = (el1_id == k) ? el2_id : el1_id;
int el1_info, el2_info;
pmesh->GetFaceInfos(bdrs[i], &el1_info, &el2_info);
const int face_id_nbr = (nbr_id == el1_id) ? el1_info / 64
: el2_info / 64;
for (j = 0; j < numFaceDofs; j++)
{
// Here it is utilized that the orientations of the face for
// the two elements are opposite of each other.
NbrDof(k,i,j) = nbr_id*nd + BdrDofs(numFaceDofs - 1 - j,
face_id_nbr);
}
}
}
else if (dim==3)
{
pmesh->GetElementFaces(k, bdrs, orientation);
for (int f = 0; f < numBdrs; f++)
{
const int nbr_cnt = face_to_el->RowSize(bdrs[f]);
if (nbr_cnt == 1)
{
// No neighbor element.
for (j = 0; j < numFaceDofs; j++) { NbrDof(k, f, j) = -1; }
continue;
}
int el1_id, el2_id, nbr_id;
pmesh->GetFaceElements(bdrs[f], &el1_id, &el2_id);
if (el2_id < 0)
{
// This element is in a different mpi task.
el2_id = -1 - el2_id + ne;
}
nbr_id = (el1_id == k) ? el2_id : el1_id;
// Local index and orientation of the face, when considered in
// the neighbor element.
int el1_info, el2_info;
pmesh->GetFaceInfos(bdrs[f], &el1_info, &el2_info);
const int face_id_nbr = (nbr_id == el1_id) ? el1_info / 64
: el2_info / 64;
const int face_or_nbr = (nbr_id == el1_id) ? el1_info % 64
: el2_info % 64;
for (j = 0; j < numFaceDofs; j++)
{
// What is the index of the j-th dof on the face, given its
// orientation.
const int loc_face_dof_id =
GetLocalFaceDofIndex(dim, face_id_nbr, face_or_nbr,
j, dof1D_cnt);
// What is the corresponding local dof id on the element,
// given the face orientation.
const int nbr_dof_id =
fdof_ids(face_or_nbr)(loc_face_dof_id, face_id_nbr);
NbrDof(k, f, j) = nbr_id*nd + nbr_dof_id;
}
}
}
}
delete face_to_el;
}
void DofInfo::FillSubcell2CellDof()
{
const int dim = pmesh->Dimension(), p = pfes.GetFE(0)->GetOrder();
if (dim==1)
{
numSubcells = p;
numDofsSubcell = 2;
}
else if (dim==2)
{
numSubcells = p*p;
numDofsSubcell = 4;
}
else if (dim==3)
{
numSubcells = p*p*p;
numDofsSubcell = 8;
}
Sub2Ind.SetSize(numSubcells, numDofsSubcell);
int aux;
for (int m = 0; m < numSubcells; m++)
{
for (int j = 0; j < numDofsSubcell; j++)
{
if (dim == 1) { Sub2Ind(m,j) = m + j; }
else if (dim == 2)
{
aux = m + m/p;
switch (j)
{
case 0: Sub2Ind(m,j) = aux; break;
case 1: Sub2Ind(m,j) = aux + 1; break;
case 2: Sub2Ind(m,j) = aux + p+1; break;
case 3: Sub2Ind(m,j) = aux + p+2; break;
}
}
else if (dim == 3)
{
aux = m + m/p + (p+1)*(m/(p*p));
switch (j)
{
case 0: Sub2Ind(m,j) = aux; break;
case 1: Sub2Ind(m,j) = aux + 1; break;
case 2: Sub2Ind(m,j) = aux + p+1; break;
case 3: Sub2Ind(m,j) = aux + p+2; break;
case 4: Sub2Ind(m,j) = aux + (p+1)*(p+1); break;
case 5: Sub2Ind(m,j) = aux + (p+1)*(p+1)+1; break;
case 6: Sub2Ind(m,j) = aux + (p+1)*(p+1)+p+1; break;
case 7: Sub2Ind(m,j) = aux + (p+1)*(p+1)+p+2; break;
}
}
}
}
}
Assembly::Assembly(DofInfo &_dofs, LowOrderMethod &_lom,
const GridFunction &inflow,
ParFiniteElementSpace &pfes, ParMesh *submesh, int mode)
: exec_mode(mode), inflow_gf(inflow), x_gf(&pfes),
VolumeTerms(NULL),
fes(&pfes), SubFes0(NULL), SubFes1(NULL),
subcell_mesh(submesh), dofs(_dofs), lom(_lom)
{
Mesh *mesh = fes->GetMesh();
int k, i, m, dim = mesh->Dimension(), ne = fes->GetNE();
Array <int> bdrs, orientation;
FaceElementTransformations *Trans;
bdrInt.SetSize(ne, dofs.numBdrs, dofs.numFaceDofs*dofs.numFaceDofs);
bdrInt = 0.;
if (lom.subcell_scheme)
{
VolumeTerms = lom.VolumeTerms;
SubcellWeights.SetSize(dofs.numSubcells, dofs.numDofsSubcell, ne);
SubFes0 = lom.SubFes0;
SubFes1 = lom.SubFes1;
}
// Initialization for transport mode.
if (exec_mode == 0)
{
for (k = 0; k < ne; k++)
{
if (dim==1) { mesh->GetElementVertices(k, bdrs); }
else if (dim==2) { mesh->GetElementEdges(k, bdrs, orientation); }
else if (dim==3) { mesh->GetElementFaces(k, bdrs, orientation); }
for (i = 0; i < dofs.numBdrs; i++)
{
Trans = mesh->GetFaceElementTransformations(bdrs[i]);
ComputeFluxTerms(k, i, Trans, lom);
}
if (lom.subcell_scheme)
{
for (m = 0; m < dofs.numSubcells; m++)
{
ComputeSubcellWeights(k, m);
}
}
}
}
}
void Assembly::ComputeFluxTerms(const int e_id, const int BdrID,
FaceElementTransformations *Trans,
LowOrderMethod &lom)
{
Mesh *mesh = fes->GetMesh();
int i, j, l, dim = mesh->Dimension();
double aux, vn;
const FiniteElement &el = *fes->GetFE(e_id);
Vector vval, nor(dim), shape(el.GetDof());
for (l = 0; l < lom.irF->GetNPoints(); l++)
{
const IntegrationPoint &ip = lom.irF->IntPoint(l);
IntegrationPoint eip1;
Trans->Face->SetIntPoint(&ip);
if (dim == 1)
{
Trans->Loc1.Transform(ip, eip1);
nor(0) = 2.*eip1.x - 1.0;
}
else
{
CalcOrtho(Trans->Face->Jacobian(), nor);
}
if (Trans->Elem1No != e_id)
{
Trans->Loc2.Transform(ip, eip1);
el.CalcShape(eip1, shape);
Trans->Elem2->SetIntPoint(&eip1);
lom.coef->Eval(vval, *Trans->Elem2, eip1);
nor *= -1.;
}
else
{
Trans->Loc1.Transform(ip, eip1);
el.CalcShape(eip1, shape);
Trans->Elem1->SetIntPoint(&eip1);
lom.coef->Eval(vval, *Trans->Elem1, eip1);
}
nor /= nor.Norml2();
if (exec_mode == 0)
{
// Transport.
vn = std::min(0., vval * nor);
}
else
{
// Remap.
vn = std::max(0., vval * nor);
vn *= -1.0;
}
const double w = ip.weight * Trans->Face->Weight();
for (i = 0; i < dofs.numFaceDofs; i++)
{
aux = w * shape(dofs.BdrDofs(i,BdrID)) * vn;
for (j = 0; j < dofs.numFaceDofs; j++)
{
bdrInt(e_id, BdrID, i*dofs.numFaceDofs+j) -=
aux * shape(dofs.BdrDofs(j,BdrID));
}
}
}
}
void Assembly::ComputeSubcellWeights(const int k, const int m)
{
DenseMatrix elmat; // These are essentially the same.
const int e_id = k*dofs.numSubcells + m;
const FiniteElement *el0 = SubFes0->GetFE(e_id);
const FiniteElement *el1 = SubFes1->GetFE(e_id);
ElementTransformation *tr = subcell_mesh->GetElementTransformation(e_id);
VolumeTerms->AssembleElementMatrix2(*el1, *el0, *tr, elmat);
for (int j = 0; j < elmat.Width(); j++)
{
// Using the fact that elmat has just one row.
SubcellWeights(k)(m,j) = elmat(0,j);
}
}
void Assembly::LinearFluxLumping(const int k, const int nd, const int BdrID,
const Vector &x, Vector &y, const Vector &x_nd,
const Vector &alpha) const
{
int i, j, dofInd;
double xNeighbor;
Vector xDiff(dofs.numFaceDofs);
const int size_x = x.Size();
for (j = 0; j < dofs.numFaceDofs; j++)
{
dofInd = k*nd+dofs.BdrDofs(j,BdrID);
const int nbr_dof_id = dofs.NbrDof(k, BdrID, j);
// Note that if the boundary is outflow, we have bdrInt = 0 by definition,
// s.t. this value will not matter.
if (nbr_dof_id < 0) { xNeighbor = inflow_gf(dofInd); }
else
{
xNeighbor = (nbr_dof_id < size_x) ? x(nbr_dof_id)
: x_nd(nbr_dof_id - size_x);
}
xDiff(j) = xNeighbor - x(dofInd);
}
for (i = 0; i < dofs.numFaceDofs; i++)
{
dofInd = k*nd+dofs.BdrDofs(i,BdrID);
for (j = 0; j < dofs.numFaceDofs; j++)
{
// alpha=0 is the low order solution, alpha=1, the Galerkin solution.
// 0 < alpha < 1 can be used for limiting within the low order method.
y(dofInd) += bdrInt(k, BdrID, i*dofs.numFaceDofs + j) *
(xDiff(i) + (xDiff(j)-xDiff(i)) *
alpha(dofs.BdrDofs(i,BdrID)) *
alpha(dofs.BdrDofs(j,BdrID)));
}
}
}
void Assembly::NonlinFluxLumping(const int k, const int nd,
const int BdrID, const Vector &x,
Vector &y, const Vector &x_nd,
const Vector &alpha) const
{
int i, j, dofInd;
double xNeighbor, SumCorrP = 0., SumCorrN = 0., eps = 1.E-15;
const int size_x = x.Size();
Vector xDiff(dofs.numFaceDofs), BdrTermCorr(dofs.numFaceDofs);
BdrTermCorr = 0.;
for (j = 0; j < dofs.numFaceDofs; j++)
{
dofInd = k*nd+dofs.BdrDofs(j,BdrID);
const int nbr_dof_id = dofs.NbrDof(k, BdrID, j);
// Note that if the boundary is outflow, we have bdrInt = 0 by definition,
// s.t. this value will not matter.
if (nbr_dof_id < 0) { xNeighbor = inflow_gf(dofInd); }
else
{
xNeighbor = (nbr_dof_id < size_x) ? x(nbr_dof_id)
: x_nd(nbr_dof_id - size_x);
}
xDiff(j) = xNeighbor - x(dofInd);
}
y.HostReadWrite();
bdrInt.HostRead();
xDiff.HostReadWrite();
for (i = 0; i < dofs.numFaceDofs; i++)
{
dofInd = k*nd+dofs.BdrDofs(i,BdrID);
for (j = 0; j < dofs.numFaceDofs; j++)
{
y(dofInd) += bdrInt(k, BdrID, i*dofs.numFaceDofs + j) * xDiff(i);
BdrTermCorr(i) += bdrInt(k, BdrID,
i*dofs.numFaceDofs + j) * (xDiff(j)-xDiff(i));
}
BdrTermCorr(i) *= alpha(dofs.BdrDofs(i,BdrID));
SumCorrP += max(0., BdrTermCorr(i));
SumCorrN += min(0., BdrTermCorr(i));
}
for (i = 0; i < dofs.numFaceDofs; i++)
{
dofInd = k*nd+dofs.BdrDofs(i,BdrID);
if (SumCorrP + SumCorrN > eps)
{
BdrTermCorr(i) = min(0., BdrTermCorr(i)) -
max(0., BdrTermCorr(i)) * SumCorrN / SumCorrP;
}
else if (SumCorrP + SumCorrN < -eps)
{
BdrTermCorr(i) = max(0., BdrTermCorr(i)) -
min(0., BdrTermCorr(i)) * SumCorrP / SumCorrN;
}
y(dofInd) += BdrTermCorr(i);
}
}
void PrecondConvectionIntegrator::AssembleElementMatrix(
const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
{
int i, nd = el.GetDof(), dim = el.GetDim();
#ifdef MFEM_THREAD_SAFE
DenseMatrix dshape, adjJ, Q_ir;
Vector shape, vec2, BdFidxT;
#endif
elmat.SetSize(nd);
dshape.SetSize(nd,dim);
adjJ.SetSize(dim);
shape.SetSize(nd);
vec2.SetSize(dim);
BdFidxT.SetSize(nd);
double w;
Vector vec1;
DenseMatrix mass(nd,nd), conv(nd,nd), lumpedM(nd,nd), tmp(nd,nd);
const IntegrationRule *ir = IntRule;
if (ir == NULL)
{
int order = Trans.OrderGrad(&el) + Trans.Order() + el.GetOrder();
order = max(order, 2 * el.GetOrder() + Trans.OrderW());
ir = &IntRules.Get(el.GetGeomType(), order);
}
Q.Eval(Q_ir, Trans, *ir);
conv = mass = 0.0;
for (i = 0; i < ir->GetNPoints(); i++)