MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
restriction.cpp
Go to the documentation of this file.
1// Copyright (c) 2010-2024, Lawrence Livermore National Security, LLC. Produced
2// at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3// LICENSE and NOTICE for details. LLNL-CODE-806117.
4//
5// This file is part of the MFEM library. For more information and source code
6// availability visit https://mfem.org.
7//
8// MFEM is free software; you can redistribute it and/or modify it under the
9// terms of the BSD-3 license. We welcome feedback and contributions, see file
10// CONTRIBUTING.md for details.
11
12#include "restriction.hpp"
14#include "gridfunc.hpp"
15#include "fespace.hpp"
16#include "pgridfunc.hpp"
17#include "qspace.hpp"
18#include "fe/face_map_utils.hpp"
19#include "../general/forall.hpp"
20
21#include <climits>
22
23namespace mfem
24{
25
27 ElementDofOrdering e_ordering)
28 : fes(f),
29 ne(fes.GetNE()),
30 vdim(fes.GetVDim()),
31 byvdim(fes.GetOrdering() == Ordering::byVDIM),
32 ndofs(fes.GetNDofs()),
33 dof(ne > 0 ? fes.GetFE(0)->GetDof() : 0),
34 nedofs(ne*dof),
35 offsets(ndofs+1),
36 indices(ne*dof),
37 gather_map(ne*dof)
38{
39 // Assuming all finite elements are the same.
40 height = vdim*ne*dof;
41 width = fes.GetVSize();
42 const bool dof_reorder = (e_ordering == ElementDofOrdering::LEXICOGRAPHIC);
43 const int *dof_map = NULL;
44 if (dof_reorder && ne > 0)
45 {
46 for (int e = 0; e < ne; ++e)
47 {
48 const FiniteElement *fe = fes.GetFE(e);
49 const TensorBasisElement* el =
50 dynamic_cast<const TensorBasisElement*>(fe);
51 if (el) { continue; }
52 MFEM_ABORT("Finite element not suitable for lexicographic ordering");
53 }
54 const FiniteElement *fe = fes.GetFE(0);
55 const TensorBasisElement* el =
56 dynamic_cast<const TensorBasisElement*>(fe);
57 const Array<int> &fe_dof_map = el->GetDofMap();
58 MFEM_VERIFY(fe_dof_map.Size() > 0, "invalid dof map");
59 dof_map = fe_dof_map.GetData();
60 }
61 const Table& e2dTable = fes.GetElementToDofTable();
62 const int* element_map = e2dTable.GetJ();
63 // We will be keeping a count of how many local nodes point to its global dof
64 for (int i = 0; i <= ndofs; ++i)
65 {
66 offsets[i] = 0;
67 }
68 for (int e = 0; e < ne; ++e)
69 {
70 for (int d = 0; d < dof; ++d)
71 {
72 const int sgid = element_map[dof*e + d]; // signed
73 const int gid = (sgid >= 0) ? sgid : -1 - sgid;
74 ++offsets[gid + 1];
75 }
76 }
77 // Aggregate to find offsets for each global dof
78 for (int i = 1; i <= ndofs; ++i)
79 {
80 offsets[i] += offsets[i - 1];
81 }
82 // For each global dof, fill in all local nodes that point to it
83 for (int e = 0; e < ne; ++e)
84 {
85 for (int d = 0; d < dof; ++d)
86 {
87 const int sdid = dof_reorder ? dof_map[d] : 0; // signed
88 const int did = (!dof_reorder)?d:(sdid >= 0 ? sdid : -1-sdid);
89 const int sgid = element_map[dof*e + did]; // signed
90 const int gid = (sgid >= 0) ? sgid : -1-sgid;
91 const int lid = dof*e + d;
92 const bool plus = (sgid >= 0 && sdid >= 0) || (sgid < 0 && sdid < 0);
93 gather_map[lid] = plus ? gid : -1-gid;
94 indices[offsets[gid]++] = plus ? lid : -1-lid;
95 }
96 }
97 // We shifted the offsets vector by 1 by using it as a counter.
98 // Now we shift it back.
99 for (int i = ndofs; i > 0; --i)
100 {
101 offsets[i] = offsets[i - 1];
102 }
103 offsets[0] = 0;
104}
105
106void ElementRestriction::Mult(const Vector& x, Vector& y) const
107{
108 // Assumes all elements have the same number of dofs
109 const int nd = dof;
110 const int vd = vdim;
111 const bool t = byvdim;
112 auto d_x = Reshape(x.Read(), t?vd:ndofs, t?ndofs:vd);
113 auto d_y = Reshape(y.Write(), nd, vd, ne);
114 auto d_gather_map = gather_map.Read();
115 mfem::forall(dof*ne, [=] MFEM_HOST_DEVICE (int i)
116 {
117 const int gid = d_gather_map[i];
118 const bool plus = gid >= 0;
119 const int j = plus ? gid : -1-gid;
120 for (int c = 0; c < vd; ++c)
121 {
122 const real_t dof_value = d_x(t?c:j, t?j:c);
123 d_y(i % nd, c, i / nd) = plus ? dof_value : -dof_value;
124 }
125 });
126}
127
129{
130 // Assumes all elements have the same number of dofs
131 const int nd = dof;
132 const int vd = vdim;
133 const bool t = byvdim;
134 auto d_x = Reshape(x.Read(), t?vd:ndofs, t?ndofs:vd);
135 auto d_y = Reshape(y.Write(), nd, vd, ne);
136 auto d_gather_map = gather_map.Read();
137
138 mfem::forall(dof*ne, [=] MFEM_HOST_DEVICE (int i)
139 {
140 const int gid = d_gather_map[i];
141 const int j = gid >= 0 ? gid : -1-gid;
142 for (int c = 0; c < vd; ++c)
143 {
144 d_y(i % nd, c, i / nd) = d_x(t?c:j, t?j:c);
145 }
146 });
147}
148
149template <bool ADD>
150void ElementRestriction::TAddMultTranspose(const Vector& x, Vector& y) const
151{
152 // Assumes all elements have the same number of dofs
153 const int nd = dof;
154 const int vd = vdim;
155 const bool t = byvdim;
156 auto d_offsets = offsets.Read();
157 auto d_indices = indices.Read();
158 auto d_x = Reshape(x.Read(), nd, vd, ne);
159 auto d_y = Reshape(ADD ? y.ReadWrite() : y.Write(), t?vd:ndofs, t?ndofs:vd);
160 mfem::forall(ndofs, [=] MFEM_HOST_DEVICE (int i)
161 {
162 const int offset = d_offsets[i];
163 const int next_offset = d_offsets[i + 1];
164 for (int c = 0; c < vd; ++c)
165 {
166 real_t dof_value = 0;
167 for (int j = offset; j < next_offset; ++j)
168 {
169 const int idx_j = (d_indices[j] >= 0) ? d_indices[j] : -1 - d_indices[j];
170 dof_value += ((d_indices[j] >= 0) ? d_x(idx_j % nd, c, idx_j / nd) :
171 -d_x(idx_j % nd, c, idx_j / nd));
172 }
173 if (ADD) { d_y(t?c:i,t?i:c) += dof_value; }
174 else { d_y(t?c:i,t?i:c) = dof_value; }
175 }
176 });
177}
178
180{
181 constexpr bool ADD = false;
182 TAddMultTranspose<ADD>(x, y);
183}
184
186 const real_t a) const
187{
188 MFEM_VERIFY(a == 1.0, "General coefficient case is not yet supported!");
189 constexpr bool ADD = true;
190 TAddMultTranspose<ADD>(x, y);
191}
192
194{
195 // Assumes all elements have the same number of dofs
196 const int nd = dof;
197 const int vd = vdim;
198 const bool t = byvdim;
199 auto d_offsets = offsets.Read();
200 auto d_indices = indices.Read();
201 auto d_x = Reshape(x.Read(), nd, vd, ne);
202 auto d_y = Reshape(y.Write(), t?vd:ndofs, t?ndofs:vd);
203 mfem::forall(ndofs, [=] MFEM_HOST_DEVICE (int i)
204 {
205 const int offset = d_offsets[i];
206 const int next_offset = d_offsets[i + 1];
207 for (int c = 0; c < vd; ++c)
208 {
209 real_t dof_value = 0;
210 for (int j = offset; j < next_offset; ++j)
211 {
212 const int idx_j = (d_indices[j] >= 0) ? d_indices[j] : -1 - d_indices[j];
213 dof_value += d_x(idx_j % nd, c, idx_j / nd);
214 }
215 d_y(t?c:i,t?i:c) = dof_value;
216 }
217 });
218}
219
221{
222 // Assumes all elements have the same number of dofs
223 const int nd = dof;
224 const int vd = vdim;
225 const bool t = byvdim;
226 auto d_offsets = offsets.Read();
227 auto d_indices = indices.Read();
228 auto d_x = Reshape(x.Read(), nd, vd, ne);
229 auto d_y = Reshape(y.Write(), t?vd:ndofs, t?ndofs:vd);
230 mfem::forall(ndofs, [=] MFEM_HOST_DEVICE (int i)
231 {
232 const int next_offset = d_offsets[i + 1];
233 for (int c = 0; c < vd; ++c)
234 {
235 real_t dof_value = 0;
236 const int j = next_offset - 1;
237 const int idx_j = (d_indices[j] >= 0) ? d_indices[j] : -1 - d_indices[j];
238 dof_value = (d_indices[j] >= 0) ? d_x(idx_j % nd, c, idx_j / nd) :
239 -d_x(idx_j % nd, c, idx_j / nd);
240 d_y(t?c:i,t?i:c) = dof_value;
241 }
242 });
243}
244
246{
247 // Assumes all elements have the same number of dofs
248 const int nd = dof;
249 const int vd = vdim;
250 const bool t = byvdim;
251
252 Array<char> processed(vd * ndofs);
253 processed = 0;
254
255 auto d_offsets = offsets.HostRead();
256 auto d_indices = indices.HostRead();
257 auto d_x = Reshape(processed.HostReadWrite(), t?vd:ndofs, t?ndofs:vd);
258 auto d_y = Reshape(y.HostWrite(), nd, vd, ne);
259 for (int i = 0; i < ndofs; ++i)
260 {
261 const int offset = d_offsets[i];
262 const int next_offset = d_offsets[i+1];
263 for (int c = 0; c < vd; ++c)
264 {
265 for (int j = offset; j < next_offset; ++j)
266 {
267 const int idx_j = d_indices[j];
268 if (d_x(t?c:i,t?i:c))
269 {
270 d_y(idx_j % nd, c, idx_j / nd) = 0.0;
271 }
272 else
273 {
274 d_y(idx_j % nd, c, idx_j / nd) = 1.0;
275 d_x(t?c:i,t?i:c) = 1;
276 }
277 }
278 }
279 }
280}
281
283 SparseMatrix &mat) const
284{
285 mat.GetMemoryI().New(mat.Height()+1, mat.GetMemoryI().GetMemoryType());
286 const int nnz = FillI(mat);
287 mat.GetMemoryJ().New(nnz, mat.GetMemoryJ().GetMemoryType());
288 mat.GetMemoryData().New(nnz, mat.GetMemoryData().GetMemoryType());
289 FillJAndData(mat_ea, mat);
290}
291
292static MFEM_HOST_DEVICE int GetMinElt(const int *my_elts, const int nbElts,
293 const int *nbr_elts, const int nbrNbElts)
294{
295 // Find the minimal element index found in both my_elts[] and nbr_elts[]
296 int min_el = INT_MAX;
297 for (int i = 0; i < nbElts; i++)
298 {
299 const int e_i = my_elts[i];
300 if (e_i >= min_el) { continue; }
301 for (int j = 0; j < nbrNbElts; j++)
302 {
303 if (e_i==nbr_elts[j])
304 {
305 min_el = e_i; // we already know e_i < min_el
306 break;
307 }
308 }
309 }
310 return min_el;
311}
312
313/** Returns the index where a non-zero entry should be added and increment the
314 number of non-zeros for the row i_L. */
315static MFEM_HOST_DEVICE int GetAndIncrementNnzIndex(const int i_L, int* I)
316{
317 int ind = AtomicAdd(I[i_L],1);
318 return ind;
319}
320
322{
323 static constexpr int Max = MaxNbNbr;
324 const int all_dofs = ndofs;
325 const int vd = vdim;
326 const int elt_dofs = dof;
327 auto I = mat.ReadWriteI();
328 auto d_offsets = offsets.Read();
329 auto d_indices = indices.Read();
330 auto d_gather_map = gather_map.Read();
331 mfem::forall(vd*all_dofs+1, [=] MFEM_HOST_DEVICE (int i_L)
332 {
333 I[i_L] = 0;
334 });
335 mfem::forall(ne*elt_dofs, [=] MFEM_HOST_DEVICE (int l_dof)
336 {
337 const int e = l_dof/elt_dofs;
338 const int i = l_dof%elt_dofs;
339
340 int i_elts[Max];
341 const int i_gm = e*elt_dofs + i;
342 const int i_L = d_gather_map[i_gm];
343 const int i_offset = d_offsets[i_L];
344 const int i_next_offset = d_offsets[i_L+1];
345 const int i_nbElts = i_next_offset - i_offset;
346 MFEM_ASSERT_KERNEL(
347 i_nbElts <= Max,
348 "The connectivity of this mesh is beyond the max, increase the "
349 "MaxNbNbr variable to comply with your mesh.");
350 for (int e_i = 0; e_i < i_nbElts; ++e_i)
351 {
352 const int i_E = d_indices[i_offset+e_i];
353 i_elts[e_i] = i_E/elt_dofs;
354 }
355 for (int j = 0; j < elt_dofs; j++)
356 {
357 const int j_gm = e*elt_dofs + j;
358 const int j_L = d_gather_map[j_gm];
359 const int j_offset = d_offsets[j_L];
360 const int j_next_offset = d_offsets[j_L+1];
361 const int j_nbElts = j_next_offset - j_offset;
362 MFEM_ASSERT_KERNEL(
363 j_nbElts <= Max,
364 "The connectivity of this mesh is beyond the max, increase the "
365 "MaxNbNbr variable to comply with your mesh.");
366 if (i_nbElts == 1 || j_nbElts == 1) // no assembly required
367 {
368 GetAndIncrementNnzIndex(i_L, I);
369 }
370 else // assembly required
371 {
372 int j_elts[Max];
373 for (int e_j = 0; e_j < j_nbElts; ++e_j)
374 {
375 const int j_E = d_indices[j_offset+e_j];
376 const int elt = j_E/elt_dofs;
377 j_elts[e_j] = elt;
378 }
379 int min_e = GetMinElt(i_elts, i_nbElts, j_elts, j_nbElts);
380 if (e == min_e) // add the nnz only once
381 {
382 GetAndIncrementNnzIndex(i_L, I);
383 }
384 }
385 }
386 });
387 // We need to sum the entries of I, we do it on CPU as it is very sequential.
388 auto h_I = mat.HostReadWriteI();
389 const int nTdofs = vd*all_dofs;
390 int sum = 0;
391 for (int i = 0; i < nTdofs; i++)
392 {
393 const int nnz = h_I[i];
394 h_I[i] = sum;
395 sum+=nnz;
396 }
397 h_I[nTdofs] = sum;
398 // We return the number of nnz
399 return h_I[nTdofs];
400}
401
403 SparseMatrix &mat) const
404{
405 static constexpr int Max = MaxNbNbr;
406 const int all_dofs = ndofs;
407 const int vd = vdim;
408 const int elt_dofs = dof;
409 auto I = mat.ReadWriteI();
410 auto J = mat.WriteJ();
411 auto Data = mat.WriteData();
412 auto d_offsets = offsets.Read();
413 auto d_indices = indices.Read();
414 auto d_gather_map = gather_map.Read();
415 auto mat_ea = Reshape(ea_data.Read(), elt_dofs, elt_dofs, ne);
416 mfem::forall(ne*elt_dofs, [=] MFEM_HOST_DEVICE (int l_dof)
417 {
418 const int e = l_dof/elt_dofs;
419 const int i = l_dof%elt_dofs;
420
421 int i_elts[Max];
422 int i_B[Max];
423 const int i_gm = e*elt_dofs + i;
424 const int i_L = d_gather_map[i_gm];
425 const int i_offset = d_offsets[i_L];
426 const int i_next_offset = d_offsets[i_L+1];
427 const int i_nbElts = i_next_offset - i_offset;
428 MFEM_ASSERT_KERNEL(
429 i_nbElts <= Max,
430 "The connectivity of this mesh is beyond the max, increase the "
431 "MaxNbNbr variable to comply with your mesh.");
432 for (int e_i = 0; e_i < i_nbElts; ++e_i)
433 {
434 const int i_E = d_indices[i_offset+e_i];
435 i_elts[e_i] = i_E/elt_dofs;
436 i_B[e_i] = i_E%elt_dofs;
437 }
438 for (int j = 0; j < elt_dofs; j++)
439 {
440 const int j_gm = e*elt_dofs + j;
441 const int j_L = d_gather_map[j_gm];
442 const int j_offset = d_offsets[j_L];
443 const int j_next_offset = d_offsets[j_L+1];
444 const int j_nbElts = j_next_offset - j_offset;
445 if (i_nbElts == 1 || j_nbElts == 1) // no assembly required
446 {
447 const int nnz = GetAndIncrementNnzIndex(i_L, I);
448 J[nnz] = j_L;
449 Data[nnz] = mat_ea(j,i,e);
450 }
451 else // assembly required
452 {
453 int j_elts[Max];
454 int j_B[Max];
455 for (int e_j = 0; e_j < j_nbElts; ++e_j)
456 {
457 const int j_E = d_indices[j_offset+e_j];
458 const int elt = j_E/elt_dofs;
459 j_elts[e_j] = elt;
460 j_B[e_j] = j_E%elt_dofs;
461 }
462 int min_e = GetMinElt(i_elts, i_nbElts, j_elts, j_nbElts);
463 if (e == min_e) // add the nnz only once
464 {
465 real_t val = 0.0;
466 for (int k = 0; k < i_nbElts; k++)
467 {
468 const int e_i = i_elts[k];
469 const int i_Bloc = i_B[k];
470 for (int l = 0; l < j_nbElts; l++)
471 {
472 const int e_j = j_elts[l];
473 const int j_Bloc = j_B[l];
474 if (e_i == e_j)
475 {
476 val += mat_ea(j_Bloc, i_Bloc, e_i);
477 }
478 }
479 }
480 const int nnz = GetAndIncrementNnzIndex(i_L, I);
481 J[nnz] = j_L;
482 Data[nnz] = val;
483 }
484 }
485 }
486 });
487 // We need to shift again the entries of I, we do it on CPU as it is very
488 // sequential.
489 auto h_I = mat.HostReadWriteI();
490 const int size = vd*all_dofs;
491 for (int i = 0; i < size; i++)
492 {
493 h_I[size-i] = h_I[size-(i+1)];
494 }
495 h_I[0] = 0;
496}
497
499 : ne(fes.GetNE()),
500 vdim(fes.GetVDim()),
501 byvdim(fes.GetOrdering() == Ordering::byVDIM),
502 ndof(ne > 0 ? fes.GetFE(0)->GetDof() : 0),
503 ndofs(fes.GetNDofs())
504{
505 height = vdim*ne*ndof;
506 width = vdim*ne*ndof;
507}
508
510{
511 const int nd = ndof;
512 const int vd = vdim;
513 const bool t = byvdim;
514 auto d_x = Reshape(x.Read(), t?vd:ndofs, t?ndofs:vd);
515 auto d_y = Reshape(y.Write(), nd, vd, ne);
516 mfem::forall(ndofs, [=] MFEM_HOST_DEVICE (int i)
517 {
518 const int idx = i;
519 const int dof = idx % nd;
520 const int e = idx / nd;
521 for (int c = 0; c < vd; ++c)
522 {
523 d_y(dof, c, e) = d_x(t?c:idx, t?idx:c);
524 }
525 });
526}
527
528template <bool ADD>
529void L2ElementRestriction::TAddMultTranspose(const Vector &x, Vector &y) const
530{
531 const int nd = ndof;
532 const int vd = vdim;
533 const bool t = byvdim;
534 auto d_x = Reshape(x.Read(), nd, vd, ne);
535 auto d_y = Reshape(ADD ? y.ReadWrite() : y.Write(), t?vd:ndofs, t?ndofs:vd);
536 mfem::forall(ndofs, [=] MFEM_HOST_DEVICE (int i)
537 {
538 const int idx = i;
539 const int dof = idx % nd;
540 const int e = idx / nd;
541 for (int c = 0; c < vd; ++c)
542 {
543 if (ADD) { d_y(t?c:idx,t?idx:c) += d_x(dof, c, e); }
544 else { d_y(t?c:idx,t?idx:c) = d_x(dof, c, e); }
545 }
546 });
547}
548
550{
551 constexpr bool ADD = false;
552 TAddMultTranspose<ADD>(x, y);
553}
554
556 const real_t a) const
557{
558 MFEM_VERIFY(a == 1.0, "General coefficient case is not yet supported!");
559 constexpr bool ADD = true;
560 TAddMultTranspose<ADD>(x, y);
561}
562
564{
565 const int elem_dofs = ndof;
566 const int vd = vdim;
567 auto I = mat.WriteI();
568 const int isize = mat.Height() + 1;
569 const int interior_dofs = ne*elem_dofs*vd;
570 mfem::forall(isize, [=] MFEM_HOST_DEVICE (int dof)
571 {
572 I[dof] = dof<interior_dofs ? elem_dofs : 0;
573 });
574}
575
576static MFEM_HOST_DEVICE int AddNnz(const int iE, int *I, const int dofs)
577{
578 int val = AtomicAdd(I[iE],dofs);
579 return val;
580}
581
583 SparseMatrix &mat) const
584{
585 const int elem_dofs = ndof;
586 const int vd = vdim;
587 auto I = mat.ReadWriteI();
588 auto J = mat.WriteJ();
589 auto Data = mat.WriteData();
590 auto mat_ea = Reshape(ea_data.Read(), elem_dofs, elem_dofs, ne);
591 mfem::forall(ne*elem_dofs*vd, [=] MFEM_HOST_DEVICE (int iE)
592 {
593 const int offset = AddNnz(iE,I,elem_dofs);
594 const int e = iE/elem_dofs;
595 const int i = iE%elem_dofs;
596 for (int j = 0; j < elem_dofs; j++)
597 {
598 J[offset+j] = e*elem_dofs+j;
599 Data[offset+j] = mat_ea(j,i,e);
600 }
601 });
602}
603
605 const FiniteElementSpace &fes,
606 const ElementDofOrdering f_ordering,
607 const FaceType type,
608 bool build)
609 : fes(fes),
610 nf(fes.GetNFbyType(type)),
611 vdim(fes.GetVDim()),
612 byvdim(fes.GetOrdering() == Ordering::byVDIM),
613 face_dofs(nf > 0 ? fes.GetFaceElement(0)->GetDof() : 0),
614 elem_dofs(fes.GetFE(0)->GetDof()),
615 nfdofs(nf*face_dofs),
616 ndofs(fes.GetNDofs()),
617 scatter_indices(nf*face_dofs),
618 gather_offsets(ndofs+1),
619 gather_indices(nf*face_dofs),
620 face_map(face_dofs)
621{
623 width = fes.GetVSize();
624 if (nf==0) { return; }
625
626 CheckFESpace(f_ordering);
627
628 // Get the mapping from lexicographic DOF ordering to native ordering.
629 const TensorBasisElement* el =
630 dynamic_cast<const TensorBasisElement*>(fes.GetFE(0));
631 const Array<int> &dof_map_ = el->GetDofMap();
632 if (dof_map_.Size() > 0)
633 {
634 vol_dof_map.MakeRef(dof_map_);
635 }
636 else
637 {
638 // For certain types of elements dof_map_ is empty. In this case, that
639 // means the element is already ordered lexicographically, so the
640 // permutation is the identity.
642 for (int i = 0; i < elem_dofs; ++i) { vol_dof_map[i] = i; }
643 }
644
645 if (!build) { return; }
646 ComputeScatterIndicesAndOffsets(f_ordering, type);
647 ComputeGatherIndices(f_ordering,type);
648}
649
651 const FiniteElementSpace &fes,
652 const ElementDofOrdering f_ordering,
653 const FaceType type)
654 : ConformingFaceRestriction(fes, f_ordering, type, true)
655{ }
656
658{
659 if (nf==0) { return; }
660 // Assumes all elements have the same number of dofs
661 const int nface_dofs = face_dofs;
662 const int vd = vdim;
663 const bool t = byvdim;
664 auto d_indices = scatter_indices.Read();
665 auto d_x = Reshape(x.Read(), t?vd:ndofs, t?ndofs:vd);
666 auto d_y = Reshape(y.Write(), nface_dofs, vd, nf);
667 mfem::forall(nfdofs, [=] MFEM_HOST_DEVICE (int i)
668 {
669 const int s_idx = d_indices[i];
670 const int sgn = (s_idx >= 0) ? 1 : -1;
671 const int idx = (s_idx >= 0) ? s_idx : -1 - s_idx;
672 const int dof = i % nface_dofs;
673 const int face = i / nface_dofs;
674 for (int c = 0; c < vd; ++c)
675 {
676 d_y(dof, c, face) = sgn*d_x(t?c:idx, t?idx:c);
677 }
678 });
679}
680
681static void ConformingFaceRestriction_AddMultTranspose(
682 const int ndofs,
683 const int face_dofs,
684 const int nf,
685 const int vdim,
686 const bool by_vdim,
687 const Array<int> &gather_offsets,
688 const Array<int> &gather_indices,
689 const Vector &x,
690 Vector &y,
691 bool use_signs,
692 const real_t a)
693{
694 MFEM_VERIFY(a == 1.0, "General coefficient case is not yet supported!");
695 if (nf==0) { return; }
696 // Assumes all elements have the same number of dofs
697 auto d_offsets = gather_offsets.Read();
698 auto d_indices = gather_indices.Read();
699 auto d_x = Reshape(x.Read(), face_dofs, vdim, nf);
700 auto d_y = Reshape(y.ReadWrite(), by_vdim?vdim:ndofs, by_vdim?ndofs:vdim);
701 mfem::forall(ndofs, [=] MFEM_HOST_DEVICE (int i)
702 {
703 const int offset = d_offsets[i];
704 const int next_offset = d_offsets[i + 1];
705 for (int c = 0; c < vdim; ++c)
706 {
707 real_t dof_value = 0;
708 for (int j = offset; j < next_offset; ++j)
709 {
710 const int s_idx_j = d_indices[j];
711 const real_t sgn = (s_idx_j >= 0 || !use_signs) ? 1.0 : -1.0;
712 const int idx_j = (s_idx_j >= 0) ? s_idx_j : -1 - s_idx_j;
713 dof_value += sgn*d_x(idx_j % face_dofs, c, idx_j / face_dofs);
714 }
715 d_y(by_vdim?c:i,by_vdim?i:c) += dof_value;
716 }
717 });
718}
719
721 const Vector& x, Vector& y, const real_t a) const
722{
723 ConformingFaceRestriction_AddMultTranspose(
725 true, a);
726}
727
729 const Vector& x, Vector& y, const real_t a) const
730{
731 ConformingFaceRestriction_AddMultTranspose(
733 false, a);
734}
735
737 f_ordering)
738{
739#ifdef MFEM_USE_MPI
740
741 // If the underlying finite element space is parallel, ensure the face
742 // neighbor information is generated.
743 if (const ParFiniteElementSpace *pfes
744 = dynamic_cast<const ParFiniteElementSpace*>(&fes))
745 {
746 pfes->GetParMesh()->ExchangeFaceNbrData();
747 }
748
749#endif
750
751#ifdef MFEM_DEBUG
752 const FiniteElement *fe0 = fes.GetFE(0);
753 const TensorBasisElement *tfe = dynamic_cast<const TensorBasisElement*>(fe0);
754 MFEM_VERIFY(tfe != NULL &&
757 "Only Gauss-Lobatto and Bernstein basis are supported in "
758 "ConformingFaceRestriction.");
759
760 // Assuming all finite elements are using Gauss-Lobatto.
761 const bool dof_reorder = (f_ordering == ElementDofOrdering::LEXICOGRAPHIC);
762 if (dof_reorder && nf > 0)
763 {
764 for (int f = 0; f < fes.GetNF(); ++f)
765 {
766 const FiniteElement *fe = fes.GetFaceElement(f);
767 const TensorBasisElement* el =
768 dynamic_cast<const TensorBasisElement*>(fe);
769 if (el) { continue; }
770 MFEM_ABORT("Finite element not suitable for lexicographic ordering");
771 }
772 }
773#endif
774}
775
776void ConformingFaceRestriction::ComputeScatterIndicesAndOffsets(
777 const ElementDofOrdering f_ordering,
778 const FaceType type)
779{
780 Mesh &mesh = *fes.GetMesh();
781
782 // Initialization of the offsets
783 for (int i = 0; i <= ndofs; ++i)
784 {
785 gather_offsets[i] = 0;
786 }
787
788 // Computation of scatter indices and offsets
789 int f_ind = 0;
790 for (int f = 0; f < fes.GetNF(); ++f)
791 {
792 Mesh::FaceInformation face = mesh.GetFaceInformation(f);
793 if ( face.IsNonconformingCoarse() )
794 {
795 // We skip nonconforming coarse faces as they are treated
796 // by the corresponding nonconforming fine faces.
797 continue;
798 }
799 else if ( face.IsOfFaceType(type) )
800 {
801 SetFaceDofsScatterIndices(face, f_ind, f_ordering);
802 f_ind++;
803 }
804 }
805 MFEM_VERIFY(f_ind==nf, "Unexpected number of faces.");
806
807 // Summation of the offsets
808 for (int i = 1; i <= ndofs; ++i)
809 {
810 gather_offsets[i] += gather_offsets[i - 1];
811 }
812}
813
814void ConformingFaceRestriction::ComputeGatherIndices(
815 const ElementDofOrdering f_ordering,
816 const FaceType type)
817{
818 Mesh &mesh = *fes.GetMesh();
819
820 // Computation of gather_indices
821 int f_ind = 0;
822 for (int f = 0; f < fes.GetNF(); ++f)
823 {
824 Mesh::FaceInformation face = mesh.GetFaceInformation(f);
825 if ( face.IsNonconformingCoarse() )
826 {
827 // We skip nonconforming coarse faces as they are treated
828 // by the corresponding nonconforming fine faces.
829 continue;
830 }
831 else if ( face.IsOfFaceType(type) )
832 {
833 SetFaceDofsGatherIndices(face, f_ind, f_ordering);
834 f_ind++;
835 }
836 }
837 MFEM_VERIFY(f_ind==nf, "Unexpected number of faces.");
838
839 // Reset offsets to their initial value
840 for (int i = ndofs; i > 0; --i)
841 {
842 gather_offsets[i] = gather_offsets[i - 1];
843 }
844 gather_offsets[0] = 0;
845}
846
847static inline int absdof(int i) { return i < 0 ? -1-i : i; }
848
850 const Mesh::FaceInformation &face,
851 const int face_index,
852 const ElementDofOrdering f_ordering)
853{
854 MFEM_ASSERT(!(face.IsNonconformingCoarse()),
855 "This method should not be used on nonconforming coarse faces.");
856 MFEM_ASSERT(face.element[0].orientation==0,
857 "FaceRestriction used on degenerated mesh.");
858 MFEM_CONTRACT_VAR(f_ordering); // not supported yet
859
861
862 const Table& e2dTable = fes.GetElementToDofTable();
863 const int* elem_map = e2dTable.GetJ();
864 const int elem_index = face.element[0].index;
865
866 for (int face_dof = 0; face_dof < face_dofs; ++face_dof)
867 {
868 const int lex_volume_dof = face_map[face_dof];
869 const int s_volume_dof = AsConst(vol_dof_map)[lex_volume_dof]; // signed
870 const int volume_dof = absdof(s_volume_dof);
871 const int s_global_dof = elem_map[elem_index*elem_dofs + volume_dof];
872 const int global_dof = absdof(s_global_dof);
873 const int restriction_dof = face_dofs*face_index + face_dof;
874 scatter_indices[restriction_dof] = s_global_dof;
875 ++gather_offsets[global_dof + 1];
876 }
877}
878
880 const Mesh::FaceInformation &face,
881 const int face_index,
882 const ElementDofOrdering f_ordering)
883{
884 MFEM_ASSERT(!(face.IsNonconformingCoarse()),
885 "This method should not be used on nonconforming coarse faces.");
886 MFEM_CONTRACT_VAR(f_ordering); // not supported yet
887
889
890 const Table& e2dTable = fes.GetElementToDofTable();
891 const int* elem_map = e2dTable.GetJ();
892 const int elem_index = face.element[0].index;
893
894 for (int face_dof = 0; face_dof < face_dofs; ++face_dof)
895 {
896 const int lex_volume_dof = face_map[face_dof];
897 const int s_volume_dof = AsConst(vol_dof_map)[lex_volume_dof];
898 const int volume_dof = absdof(s_volume_dof);
899 const int s_global_dof = elem_map[elem_index*elem_dofs + volume_dof];
900 const int sgn = (s_global_dof >= 0) ? 1 : -1;
901 const int global_dof = absdof(s_global_dof);
902 const int restriction_dof = face_dofs*face_index + face_dof;
903 const int s_restriction_dof = (sgn >= 0) ? restriction_dof : -1 -
904 restriction_dof;
905 gather_indices[gather_offsets[global_dof]++] = s_restriction_dof;
906 }
907}
908
909// Permute dofs or quads on a face for e2 to match with the ordering of e1
910int PermuteFaceL2(const int dim, const int face_id1,
911 const int face_id2, const int orientation,
912 const int size1d, const int index)
913{
914 switch (dim)
915 {
916 case 1:
917 return 0;
918 case 2:
919 return internal::PermuteFace2D(face_id1, face_id2, orientation, size1d, index);
920 case 3:
921 return internal::PermuteFace3D(face_id1, face_id2, orientation, size1d, index);
922 default:
923 MFEM_ABORT("Unsupported dimension.");
924 return 0;
925 }
926}
927
929 const ElementDofOrdering f_ordering,
930 const FaceType type,
931 const L2FaceValues m,
932 bool build)
933 : fes(fes),
934 ordering(f_ordering),
935 nf(fes.GetNFbyType(type)),
936 ne(fes.GetNE()),
937 vdim(fes.GetVDim()),
938 byvdim(fes.GetOrdering() == Ordering::byVDIM),
939 face_dofs(nf > 0 ?
940 fes.GetTraceElement(0, fes.GetMesh()->GetFaceGeometry(0))->GetDof()
941 : 0),
942 elem_dofs(fes.GetFE(0)->GetDof()),
943 nfdofs(nf*face_dofs),
944 ndofs(fes.GetNDofs()),
945 type(type),
946 m(m),
947 scatter_indices1(nf*face_dofs),
948 scatter_indices2(m==L2FaceValues::DoubleValued?nf*face_dofs:0),
949 gather_offsets(ndofs+1),
950 gather_indices((m==L2FaceValues::DoubleValued? 2 : 1)*nf*face_dofs),
951 face_map(face_dofs)
952{
954 width = fes.GetVSize();
955 if (!build) { return; }
956
957 CheckFESpace();
958 ComputeScatterIndicesAndOffsets();
959 ComputeGatherIndices();
960}
961
963 const ElementDofOrdering f_ordering,
964 const FaceType type,
965 const L2FaceValues m)
966 : L2FaceRestriction(fes, f_ordering, type, m, true)
967{ }
968
970 Vector& y) const
971{
972 if (nf == 0) { return; }
973 MFEM_ASSERT(
975 "This method should be called when m == L2FaceValues::SingleValued.");
976 // Assumes all elements have the same number of dofs
977 const int nface_dofs = face_dofs;
978 const int vd = vdim;
979 const bool t = byvdim;
980 auto d_indices1 = scatter_indices1.Read();
981 auto d_x = Reshape(x.Read(), t?vd:ndofs, t?ndofs:vd);
982 auto d_y = Reshape(y.Write(), nface_dofs, vd, nf);
983 mfem::forall(nfdofs, [=] MFEM_HOST_DEVICE (int i)
984 {
985 const int dof = i % nface_dofs;
986 const int face = i / nface_dofs;
987 const int idx1 = d_indices1[i];
988 for (int c = 0; c < vd; ++c)
989 {
990 d_y(dof, c, face) = d_x(t?c:idx1, t?idx1:c);
991 }
992 });
993}
994
996 Vector& y) const
997{
998 MFEM_ASSERT(
1000 "This method should be called when m == L2FaceValues::DoubleValued.");
1001 // Assumes all elements have the same number of dofs
1002 const int nface_dofs = face_dofs;
1003 const int vd = vdim;
1004 const bool t = byvdim;
1005 auto d_indices1 = scatter_indices1.Read();
1006 auto d_indices2 = scatter_indices2.Read();
1007 auto d_x = Reshape(x.Read(), t?vd:ndofs, t?ndofs:vd);
1008 auto d_y = Reshape(y.Write(), nface_dofs, vd, 2, nf);
1009 mfem::forall(nfdofs, [=] MFEM_HOST_DEVICE (int i)
1010 {
1011 const int dof = i % nface_dofs;
1012 const int face = i / nface_dofs;
1013 const int idx1 = d_indices1[i];
1014 for (int c = 0; c < vd; ++c)
1015 {
1016 d_y(dof, c, 0, face) = d_x(t?c:idx1, t?idx1:c);
1017 }
1018 const int idx2 = d_indices2[i];
1019 for (int c = 0; c < vd; ++c)
1020 {
1021 d_y(dof, c, 1, face) = idx2==-1 ? 0.0 : d_x(t?c:idx2, t?idx2:c);
1022 }
1023 });
1024}
1025
1026void L2FaceRestriction::Mult(const Vector& x, Vector& y) const
1027{
1028 if (nf==0) { return; }
1030 {
1032 }
1033 else
1034 {
1036 }
1037}
1038
1040 const Vector& x, Vector& y) const
1041{
1042 // Assumes all elements have the same number of dofs
1043 const int nface_dofs = face_dofs;
1044 const int vd = vdim;
1045 const bool t = byvdim;
1046 auto d_offsets = gather_offsets.Read();
1047 auto d_indices = gather_indices.Read();
1048 auto d_x = Reshape(x.Read(), nface_dofs, vd, nf);
1049 auto d_y = Reshape(y.ReadWrite(), t?vd:ndofs, t?ndofs:vd);
1050 mfem::forall(ndofs, [=] MFEM_HOST_DEVICE (int i)
1051 {
1052 const int offset = d_offsets[i];
1053 const int next_offset = d_offsets[i + 1];
1054 for (int c = 0; c < vd; ++c)
1055 {
1056 real_t dof_value = 0;
1057 for (int j = offset; j < next_offset; ++j)
1058 {
1059 int idx_j = d_indices[j];
1060 dof_value += d_x(idx_j % nface_dofs, c, idx_j / nface_dofs);
1061 }
1062 d_y(t?c:i,t?i:c) += dof_value;
1063 }
1064 });
1065}
1066
1068 const Vector& x, Vector& y) const
1069{
1070 // Assumes all elements have the same number of dofs
1071 const int nface_dofs = face_dofs;
1072 const int vd = vdim;
1073 const bool t = byvdim;
1074 const int dofs = nfdofs;
1075 auto d_offsets = gather_offsets.Read();
1076 auto d_indices = gather_indices.Read();
1077 auto d_x = Reshape(x.Read(), nface_dofs, vd, 2, nf);
1078 auto d_y = Reshape(y.ReadWrite(), t?vd:ndofs, t?ndofs:vd);
1079 mfem::forall(ndofs, [=] MFEM_HOST_DEVICE (int i)
1080 {
1081 const int offset = d_offsets[i];
1082 const int next_offset = d_offsets[i + 1];
1083 for (int c = 0; c < vd; ++c)
1084 {
1085 real_t dof_value = 0;
1086 for (int j = offset; j < next_offset; ++j)
1087 {
1088 int idx_j = d_indices[j];
1089 bool isE1 = idx_j < dofs;
1090 idx_j = isE1 ? idx_j : idx_j - dofs;
1091 dof_value += isE1 ?
1092 d_x(idx_j % nface_dofs, c, 0, idx_j / nface_dofs)
1093 :d_x(idx_j % nface_dofs, c, 1, idx_j / nface_dofs);
1094 }
1095 d_y(t?c:i,t?i:c) += dof_value;
1096 }
1097 });
1098}
1099
1101 const real_t a) const
1102{
1103 MFEM_VERIFY(a == 1.0, "General coefficient case is not yet supported!");
1104 if (nf==0) { return; }
1106 {
1108 }
1109 else
1110 {
1112 }
1113}
1114
1116 const bool keep_nbr_block) const
1117{
1118 const int nface_dofs = face_dofs;
1119 auto d_indices1 = scatter_indices1.Read();
1120 auto d_indices2 = scatter_indices2.Read();
1121 auto I = mat.ReadWriteI();
1122 mfem::forall(nf*nface_dofs, [=] MFEM_HOST_DEVICE (int fdof)
1123 {
1124 const int iE1 = d_indices1[fdof];
1125 const int iE2 = d_indices2[fdof];
1126 AddNnz(iE1,I,nface_dofs);
1127 AddNnz(iE2,I,nface_dofs);
1128 });
1129}
1130
1132 SparseMatrix &mat,
1133 const bool keep_nbr_block) const
1134{
1135 const int nface_dofs = face_dofs;
1136 auto d_indices1 = scatter_indices1.Read();
1137 auto d_indices2 = scatter_indices2.Read();
1138 auto I = mat.ReadWriteI();
1139 auto mat_fea = Reshape(fea_data.Read(), nface_dofs, nface_dofs, 2, nf);
1140 auto J = mat.WriteJ();
1141 auto Data = mat.WriteData();
1142 mfem::forall(nf*nface_dofs, [=] MFEM_HOST_DEVICE (int fdof)
1143 {
1144 const int f = fdof/nface_dofs;
1145 const int iF = fdof%nface_dofs;
1146 const int iE1 = d_indices1[f*nface_dofs+iF];
1147 const int iE2 = d_indices2[f*nface_dofs+iF];
1148 const int offset1 = AddNnz(iE1,I,nface_dofs);
1149 const int offset2 = AddNnz(iE2,I,nface_dofs);
1150 for (int jF = 0; jF < nface_dofs; jF++)
1151 {
1152 const int jE1 = d_indices1[f*nface_dofs+jF];
1153 const int jE2 = d_indices2[f*nface_dofs+jF];
1154 J[offset2+jF] = jE1;
1155 J[offset1+jF] = jE2;
1156 Data[offset2+jF] = mat_fea(jF,iF,0,f);
1157 Data[offset1+jF] = mat_fea(jF,iF,1,f);
1158 }
1159 });
1160}
1161
1163 Vector &ea_data) const
1164{
1165 const int nface_dofs = face_dofs;
1166 const int nelem_dofs = elem_dofs;
1167 const int NE = ne;
1169 {
1170 auto d_indices1 = scatter_indices1.Read();
1171 auto d_indices2 = scatter_indices2.Read();
1172 auto mat_fea = Reshape(fea_data.Read(), nface_dofs, nface_dofs, 2, nf);
1173 auto mat_ea = Reshape(ea_data.ReadWrite(), nelem_dofs, nelem_dofs, ne);
1174 mfem::forall(nf, [=] MFEM_HOST_DEVICE (int f)
1175 {
1176 const int e1 = d_indices1[f*nface_dofs]/nelem_dofs;
1177 const int e2 = d_indices2[f*nface_dofs]/nelem_dofs;
1178 for (int j = 0; j < nface_dofs; j++)
1179 {
1180 const int jB1 = d_indices1[f*nface_dofs+j]%nelem_dofs;
1181 for (int i = 0; i < nface_dofs; i++)
1182 {
1183 const int iB1 = d_indices1[f*nface_dofs+i]%nelem_dofs;
1184 AtomicAdd(mat_ea(iB1,jB1,e1), mat_fea(i,j,0,f));
1185 }
1186 }
1187 if (e2 < NE)
1188 {
1189 for (int j = 0; j < nface_dofs; j++)
1190 {
1191 const int jB2 = d_indices2[f*nface_dofs+j]%nelem_dofs;
1192 for (int i = 0; i < nface_dofs; i++)
1193 {
1194 const int iB2 = d_indices2[f*nface_dofs+i]%nelem_dofs;
1195 AtomicAdd(mat_ea(iB2,jB2,e2), mat_fea(i,j,1,f));
1196 }
1197 }
1198 }
1199 });
1200 }
1201 else
1202 {
1203 auto d_indices = scatter_indices1.Read();
1204 auto mat_fea = Reshape(fea_data.Read(), nface_dofs, nface_dofs, nf);
1205 auto mat_ea = Reshape(ea_data.ReadWrite(), nelem_dofs, nelem_dofs, ne);
1206 mfem::forall(nf, [=] MFEM_HOST_DEVICE (int f)
1207 {
1208 const int e = d_indices[f*nface_dofs]/nelem_dofs;
1209 for (int j = 0; j < nface_dofs; j++)
1210 {
1211 const int jE = d_indices[f*nface_dofs+j]%nelem_dofs;
1212 for (int i = 0; i < nface_dofs; i++)
1213 {
1214 const int iE = d_indices[f*nface_dofs+i]%nelem_dofs;
1215 AtomicAdd(mat_ea(iE,jE,e), mat_fea(i,j,f));
1216 }
1217 }
1218 });
1219 }
1220}
1221
1223{
1224#ifdef MFEM_USE_MPI
1225
1226 // If the underlying finite element space is parallel, ensure the face
1227 // neighbor information is generated.
1228 if (const ParFiniteElementSpace *pfes
1229 = dynamic_cast<const ParFiniteElementSpace*>(&fes))
1230 {
1231 pfes->GetParMesh()->ExchangeFaceNbrData();
1232 }
1233
1234#endif
1235
1236#ifdef MFEM_DEBUG
1237 // If fespace == L2
1238 const FiniteElement *fe0 = fes.GetFE(0);
1239 const TensorBasisElement *tfe = dynamic_cast<const TensorBasisElement*>(fe0);
1240 MFEM_VERIFY(tfe != NULL &&
1243 "Only Gauss-Lobatto and Bernstein basis are supported in "
1244 "L2FaceRestriction.");
1245 if (nf==0) { return; }
1246 const bool dof_reorder = (ordering == ElementDofOrdering::LEXICOGRAPHIC);
1247 if (!dof_reorder)
1248 {
1249 MFEM_ABORT("Non-Tensor L2FaceRestriction not yet implemented.");
1250 }
1251 if (dof_reorder && nf > 0)
1252 {
1253 for (int f = 0; f < fes.GetNF(); ++f)
1254 {
1255 const FiniteElement *fe =
1257 const TensorBasisElement* el =
1258 dynamic_cast<const TensorBasisElement*>(fe);
1259 if (el) { continue; }
1260 MFEM_ABORT("Finite element not suitable for lexicographic ordering");
1261 }
1262 }
1263#endif
1264}
1265
1266void L2FaceRestriction::ComputeScatterIndicesAndOffsets()
1267{
1268 Mesh &mesh = *fes.GetMesh();
1269 // Initialization of the offsets
1270 for (int i = 0; i <= ndofs; ++i)
1271 {
1272 gather_offsets[i] = 0;
1273 }
1274
1275 // Computation of scatter indices and offsets
1276 int f_ind=0;
1277 for (int f = 0; f < fes.GetNF(); ++f)
1278 {
1279 Mesh::FaceInformation face = mesh.GetFaceInformation(f);
1280 MFEM_ASSERT(!face.IsShared(),
1281 "Unexpected shared face in L2FaceRestriction.");
1282 if ( face.IsOfFaceType(type) )
1283 {
1284 SetFaceDofsScatterIndices1(face,f_ind);
1286 {
1287 if ( type==FaceType::Interior && face.IsInterior() )
1288 {
1290 }
1291 else if ( type==FaceType::Boundary && face.IsBoundary() )
1292 {
1294 }
1295 }
1296 f_ind++;
1297 }
1298 }
1299 MFEM_VERIFY(f_ind==nf, "Unexpected number of faces.");
1300
1301 // Summation of the offsets
1302 for (int i = 1; i <= ndofs; ++i)
1303 {
1304 gather_offsets[i] += gather_offsets[i - 1];
1305 }
1306}
1307
1308void L2FaceRestriction::ComputeGatherIndices()
1309{
1310 Mesh &mesh = *fes.GetMesh();
1311 // Computation of gather_indices
1312 int f_ind = 0;
1313 for (int f = 0; f < fes.GetNF(); ++f)
1314 {
1315 Mesh::FaceInformation face = mesh.GetFaceInformation(f);
1316 MFEM_ASSERT(!face.IsShared(),
1317 "Unexpected shared face in L2FaceRestriction.");
1318 if ( face.IsOfFaceType(type) )
1319 {
1320 SetFaceDofsGatherIndices1(face,f_ind);
1323 face.IsLocal())
1324 {
1326 }
1327 f_ind++;
1328 }
1329 }
1330 MFEM_VERIFY(f_ind==nf, "Unexpected number of faces.");
1331
1332 // Reset offsets to their correct value
1333 for (int i = ndofs; i > 0; --i)
1334 {
1335 gather_offsets[i] = gather_offsets[i - 1];
1336 }
1337 gather_offsets[0] = 0;
1338}
1339
1341 const Mesh::FaceInformation &face,
1342 const int face_index)
1343{
1344 MFEM_ASSERT(!(face.IsNonconformingCoarse()),
1345 "This method should not be used on nonconforming coarse faces.");
1346 const Table& e2dTable = fes.GetElementToDofTable();
1347 const int* elem_map = e2dTable.GetJ();
1348 const int face_id1 = face.element[0].local_face_id;
1349 const int elem_index = face.element[0].index;
1350 fes.GetFE(0)->GetFaceMap(face_id1, face_map);
1351
1352 for (int face_dof_elem1 = 0; face_dof_elem1 < face_dofs; ++face_dof_elem1)
1353 {
1354 const int volume_dof_elem1 = face_map[face_dof_elem1];
1355 const int global_dof_elem1 = elem_map[elem_index*elem_dofs + volume_dof_elem1];
1356 const int restriction_dof_elem1 = face_dofs*face_index + face_dof_elem1;
1357 scatter_indices1[restriction_dof_elem1] = global_dof_elem1;
1358 ++gather_offsets[global_dof_elem1 + 1];
1359 }
1360}
1361
1363 const Mesh::FaceInformation &face,
1364 const int face_index)
1365{
1366 MFEM_ASSERT(face.IsLocal(),
1367 "This method should only be used on local faces.");
1368 const Table& e2dTable = fes.GetElementToDofTable();
1369 const int* elem_map = e2dTable.GetJ();
1370 const int elem_index = face.element[1].index;
1371 const int face_id1 = face.element[0].local_face_id;
1372 const int face_id2 = face.element[1].local_face_id;
1373 const int orientation = face.element[1].orientation;
1374 const int dim = fes.GetMesh()->Dimension();
1375 const int dof1d = fes.GetFE(0)->GetOrder()+1;
1376 fes.GetFE(0)->GetFaceMap(face_id2, face_map);
1377
1378 for (int face_dof_elem1 = 0; face_dof_elem1 < face_dofs; ++face_dof_elem1)
1379 {
1380 const int face_dof_elem2 = PermuteFaceL2(dim, face_id1, face_id2,
1381 orientation, dof1d,
1382 face_dof_elem1);
1383 const int volume_dof_elem2 = face_map[face_dof_elem2];
1384 const int global_dof_elem2 = elem_map[elem_index*elem_dofs + volume_dof_elem2];
1385 const int restriction_dof_elem2 = face_dofs*face_index + face_dof_elem1;
1386 scatter_indices2[restriction_dof_elem2] = global_dof_elem2;
1387 ++gather_offsets[global_dof_elem2 + 1];
1388 }
1389}
1390
1392 const Mesh::FaceInformation &face,
1393 const int face_index)
1394{
1395#ifdef MFEM_USE_MPI
1396 MFEM_ASSERT(face.IsShared(),
1397 "This method should only be used on shared faces.");
1398 const int elem_index = face.element[1].index;
1399 const int face_id1 = face.element[0].local_face_id;
1400 const int face_id2 = face.element[1].local_face_id;
1401 const int orientation = face.element[1].orientation;
1402 const int dim = fes.GetMesh()->Dimension();
1403 const int dof1d = fes.GetFE(0)->GetOrder()+1;
1404 fes.GetFE(0)->GetFaceMap(face_id2, face_map);
1405 Array<int> face_nbr_dofs;
1406 const ParFiniteElementSpace &pfes =
1407 static_cast<const ParFiniteElementSpace&>(this->fes);
1408 pfes.GetFaceNbrElementVDofs(elem_index, face_nbr_dofs);
1409
1410 for (int face_dof_elem1 = 0; face_dof_elem1 < face_dofs; ++face_dof_elem1)
1411 {
1412 const int face_dof_elem2 = PermuteFaceL2(dim, face_id1, face_id2,
1413 orientation, dof1d, face_dof_elem1);
1414 const int volume_dof_elem2 = face_map[face_dof_elem2];
1415 const int global_dof_elem2 = face_nbr_dofs[volume_dof_elem2];
1416 const int restriction_dof_elem2 = face_dofs*face_index + face_dof_elem1;
1417 // Trick to differentiate dof location inter/shared
1418 scatter_indices2[restriction_dof_elem2] = ndofs+global_dof_elem2;
1419 }
1420#endif
1421}
1422
1424 const Mesh::FaceInformation &face,
1425 const int face_index)
1426{
1427 MFEM_ASSERT(face.IsBoundary(),
1428 "This method should only be used on boundary faces.");
1429
1430 for (int d = 0; d < face_dofs; ++d)
1431 {
1432 const int restriction_dof_elem2 = face_dofs*face_index + d;
1433 scatter_indices2[restriction_dof_elem2] = -1;
1434 }
1435}
1436
1438 const Mesh::FaceInformation &face,
1439 const int face_index)
1440{
1441 MFEM_ASSERT(!(face.IsNonconformingCoarse()),
1442 "This method should not be used on nonconforming coarse faces.");
1443 const Table& e2dTable = fes.GetElementToDofTable();
1444 const int* elem_map = e2dTable.GetJ();
1445 const int face_id1 = face.element[0].local_face_id;
1446 const int elem_index = face.element[0].index;
1447 fes.GetFE(0)->GetFaceMap(face_id1, face_map);
1448
1449 for (int face_dof_elem1 = 0; face_dof_elem1 < face_dofs; ++face_dof_elem1)
1450 {
1451 const int volume_dof_elem1 = face_map[face_dof_elem1];
1452 const int global_dof_elem1 = elem_map[elem_index*elem_dofs + volume_dof_elem1];
1453 const int restriction_dof_elem1 = face_dofs*face_index + face_dof_elem1;
1454 // We don't shift restriction_dof_elem1 to express that it's elem1 of the face
1455 gather_indices[gather_offsets[global_dof_elem1]++] = restriction_dof_elem1;
1456 }
1457}
1458
1460 const Mesh::FaceInformation &face,
1461 const int face_index)
1462{
1463 MFEM_ASSERT(face.IsLocal(),
1464 "This method should only be used on local faces.");
1465 const Table& e2dTable = fes.GetElementToDofTable();
1466 const int* elem_map = e2dTable.GetJ();
1467 const int elem_index = face.element[1].index;
1468 const int face_id1 = face.element[0].local_face_id;
1469 const int face_id2 = face.element[1].local_face_id;
1470 const int orientation = face.element[1].orientation;
1471 const int dim = fes.GetMesh()->Dimension();
1472 const int dof1d = fes.GetFE(0)->GetOrder()+1;
1473 fes.GetFE(0)->GetFaceMap(face_id2, face_map);
1474
1475 for (int face_dof_elem1 = 0; face_dof_elem1 < face_dofs; ++face_dof_elem1)
1476 {
1477 const int face_dof_elem2 = PermuteFaceL2(dim, face_id1, face_id2,
1478 orientation, dof1d,
1479 face_dof_elem1);
1480 const int volume_dof_elem2 = face_map[face_dof_elem2];
1481 const int global_dof_elem2 = elem_map[elem_index*elem_dofs + volume_dof_elem2];
1482 const int restriction_dof_elem2 = face_dofs*face_index + face_dof_elem1;
1483 // We shift restriction_dof_elem2 to express that it's elem2 of the face
1484 gather_indices[gather_offsets[global_dof_elem2]++] = nfdofs +
1485 restriction_dof_elem2;
1486 }
1487}
1488
1490{
1491 EnsureNormalDerivativeRestriction();
1492 normal_deriv_restr->Mult(x, y);
1493}
1494
1496 Vector &y) const
1497{
1498 EnsureNormalDerivativeRestriction();
1499 normal_deriv_restr->AddMultTranspose(x, y);
1500}
1501
1502void L2FaceRestriction::EnsureNormalDerivativeRestriction() const
1503{
1504 if (!normal_deriv_restr)
1505 {
1506 normal_deriv_restr.reset(
1508 }
1509}
1510
1512 ElementDofOrdering ordering,
1513 FaceType type)
1514 : fes(fes),
1515 ordering(ordering),
1516 interp_config( fes.GetNFbyType(type) ),
1517 nc_cpt(0)
1518{ }
1519
1521 const Mesh::FaceInformation &face,
1522 int face_index)
1523{
1524 interp_config[face_index] = InterpConfig();
1525}
1526
1528 const Mesh::FaceInformation &face,
1529 int face_index)
1530{
1531 MFEM_ASSERT(!face.IsConforming(),
1532 "Registering face as nonconforming even though it is not.");
1533 const DenseMatrix* ptMat = face.point_matrix;
1534 // In the case of nonconforming slave shared face the master face is elem1.
1535 const int master_side =
1537 const int face_key = (master_side == 0 ? 1000 : 0) +
1538 face.element[0].local_face_id +
1539 6*face.element[1].local_face_id +
1540 36*face.element[1].orientation ;
1541 // Unfortunately we can't trust unicity of the ptMat to identify the transformation.
1542 Key key(ptMat, face_key);
1543 auto itr = interp_map.find(key);
1544 if ( itr == interp_map.end() )
1545 {
1546 const DenseMatrix* interpolator =
1547 GetCoarseToFineInterpolation(face,ptMat);
1548 interp_map[key] = {nc_cpt, interpolator};
1549 interp_config[face_index] = {master_side, nc_cpt};
1550 nc_cpt++;
1551 }
1552 else
1553 {
1554 interp_config[face_index] = {master_side, itr->second.first};
1555 }
1556}
1557
1558const DenseMatrix* InterpolationManager::GetCoarseToFineInterpolation(
1559 const Mesh::FaceInformation &face,
1560 const DenseMatrix* ptMat)
1561{
1563 "The following interpolation operator is only implemented for"
1564 "lexicographic ordering.");
1565 MFEM_VERIFY(!face.IsConforming(),
1566 "This method should not be called on conforming faces.")
1567 const int face_id1 = face.element[0].local_face_id;
1568 const int face_id2 = face.element[1].local_face_id;
1569
1570 const bool is_ghost_slave =
1571 face.element[0].conformity == Mesh::ElementConformity::Superset;
1572 const int master_face_id = is_ghost_slave ? face_id1 : face_id2;
1573
1574 // Computation of the interpolation matrix from master
1575 // (coarse) face to slave (fine) face.
1576 // Assumes all trace elements are the same.
1577 const FiniteElement *trace_fe =
1578 fes.GetTraceElement(0, fes.GetMesh()->GetFaceGeometry(0));
1579 const int face_dofs = trace_fe->GetDof();
1580 const TensorBasisElement* el =
1581 dynamic_cast<const TensorBasisElement*>(trace_fe);
1582 const auto dof_map = el->GetDofMap();
1583 DenseMatrix* interpolator = new DenseMatrix(face_dofs,face_dofs);
1584 Vector shape(face_dofs);
1585
1587 isotr.SetIdentityTransformation(trace_fe->GetGeomType());
1588 isotr.SetPointMat(*ptMat);
1589 DenseMatrix& trans_pt_mat = isotr.GetPointMat();
1590 // PointMatrix needs to be flipped in 2D
1591 if ( trace_fe->GetGeomType()==Geometry::SEGMENT && !is_ghost_slave )
1592 {
1593 std::swap(trans_pt_mat(0,0),trans_pt_mat(0,1));
1594 }
1595 DenseMatrix native_interpolator(face_dofs,face_dofs);
1596 trace_fe->GetLocalInterpolation(isotr, native_interpolator);
1597 const int dim = trace_fe->GetDim()+1;
1598 const int dof1d = trace_fe->GetOrder()+1;
1599 const int orientation = face.element[1].orientation;
1600 for (int i = 0; i < face_dofs; i++)
1601 {
1602 const int ni = (dof_map.Size()==0) ? i : dof_map[i];
1603 int li = ToLexOrdering(dim, master_face_id, dof1d, i);
1604 if ( !is_ghost_slave )
1605 {
1606 // master side is elem 2, so we permute to order dofs as elem 1.
1607 li = PermuteFaceL2(dim, face_id2, face_id1,
1608 orientation, dof1d, li);
1609 }
1610 for (int j = 0; j < face_dofs; j++)
1611 {
1612 int lj = ToLexOrdering(dim, master_face_id, dof1d, j);
1613 if ( !is_ghost_slave )
1614 {
1615 // master side is elem 2, so we permute to order dofs as elem 1.
1616 lj = PermuteFaceL2(dim, face_id2, face_id1,
1617 orientation, dof1d, lj);
1618 }
1619 const int nj = (dof_map.Size()==0) ? j : dof_map[j];
1620 (*interpolator)(li,lj) = native_interpolator(ni,nj);
1621 }
1622 }
1623 return interpolator;
1624}
1625
1627{
1628 // Assumes all trace elements are the same.
1629 const FiniteElement *trace_fe =
1631 const int face_dofs = trace_fe->GetDof();
1632 const int nc_size = interp_map.size();
1633 MFEM_VERIFY(nc_cpt==nc_size, "Unexpected number of interpolators.");
1634 interpolators.SetSize(face_dofs*face_dofs*nc_size);
1635 auto d_interp = Reshape(interpolators.HostWrite(),face_dofs,face_dofs,nc_size);
1636 for (auto val : interp_map)
1637 {
1638 const int idx = val.second.first;
1639 const DenseMatrix &interpolator = *val.second.second;
1640 for (int i = 0; i < face_dofs; i++)
1641 {
1642 for (int j = 0; j < face_dofs; j++)
1643 {
1644 d_interp(i,j,idx) = interpolator(i,j);
1645 }
1646 }
1647 delete val.second.second;
1648 }
1649 interp_map.clear();
1650}
1651
1653{
1654 // Count nonconforming faces
1655 int num_nc_faces = 0;
1656 for (int i = 0; i < interp_config.Size(); i++)
1657 {
1658 if ( interp_config[i].is_non_conforming )
1659 {
1660 num_nc_faces++;
1661 }
1662 }
1663 // Set nc_interp_config
1664 nc_interp_config.SetSize(num_nc_faces);
1665 int nc_index = 0;
1666 for (int i = 0; i < interp_config.Size(); i++)
1667 {
1668 auto & config = interp_config[i];
1669 if ( config.is_non_conforming )
1670 {
1671 nc_interp_config[nc_index] = NCInterpConfig(i, config);
1672 nc_index++;
1673 }
1674 }
1675}
1676
1678 const ElementDofOrdering f_ordering,
1679 const FaceType type,
1680 const L2FaceValues m,
1681 bool build)
1682 : L2FaceRestriction(fes, f_ordering, type, m, false),
1683 interpolations(fes, f_ordering, type)
1684{
1685 if (!build) { return; }
1686 x_interp.UseDevice(true);
1687
1688 CheckFESpace();
1689
1690 ComputeScatterIndicesAndOffsets();
1691
1692 ComputeGatherIndices();
1693}
1694
1696 const ElementDofOrdering f_ordering,
1697 const FaceType type,
1698 const L2FaceValues m)
1699 : NCL2FaceRestriction(fes, f_ordering, type, m, true)
1700{ }
1701
1708
1710 Vector& y) const
1711{
1712 if (nf == 0) { return; }
1713 // Assumes all elements have the same number of dofs
1714 const int nface_dofs = face_dofs;
1715 const int vd = vdim;
1716 auto d_y = Reshape(y.ReadWrite(), nface_dofs, vd, 2, nf);
1717 auto &nc_interp_config = interpolations.GetNCFaceInterpConfig();
1718 const int num_nc_faces = nc_interp_config.Size();
1719 if ( num_nc_faces == 0 ) { return; }
1720 auto interp_config_ptr = nc_interp_config.Read();
1721 const int nc_size = interpolations.GetNumInterpolators();
1722 auto d_interp = Reshape(interpolations.GetInterpolators().Read(),
1723 nface_dofs, nface_dofs, nc_size);
1724 static constexpr int max_nd = 16*16;
1725 MFEM_VERIFY(nface_dofs<=max_nd, "Too many degrees of freedom.");
1726 mfem::forall_2D(num_nc_faces, nface_dofs, 1, [=] MFEM_HOST_DEVICE (int nc_face)
1727 {
1728 MFEM_SHARED real_t dof_values[max_nd];
1729 const NCInterpConfig conf = interp_config_ptr[nc_face];
1730 if ( conf.is_non_conforming )
1731 {
1732 const int master_side = conf.master_side;
1733 const int interp_index = conf.index;
1734 const int face = conf.face_index;
1735 for (int c = 0; c < vd; ++c)
1736 {
1737 MFEM_FOREACH_THREAD(dof,x,nface_dofs)
1738 {
1739 dof_values[dof] = d_y(dof, c, master_side, face);
1740 }
1741 MFEM_SYNC_THREAD;
1742 MFEM_FOREACH_THREAD(dof_out,x,nface_dofs)
1743 {
1744 real_t res = 0.0;
1745 for (int dof_in = 0; dof_in<nface_dofs; dof_in++)
1746 {
1747 res += d_interp(dof_out, dof_in, interp_index)*dof_values[dof_in];
1748 }
1749 d_y(dof_out, c, master_side, face) = res;
1750 }
1751 MFEM_SYNC_THREAD;
1752 }
1753 }
1754 });
1755}
1756
1757void NCL2FaceRestriction::Mult(const Vector& x, Vector& y) const
1758{
1759 if ( type==FaceType::Interior && m==L2FaceValues::DoubleValued )
1760 {
1761 DoubleValuedNonconformingMult(x, y);
1762 }
1763 else if ( type==FaceType::Boundary && m==L2FaceValues::DoubleValued )
1764 {
1765 DoubleValuedConformingMult(x, y);
1766 }
1767 else // Single valued (assumes no nonconforming master on elem1)
1768 {
1769 SingleValuedConformingMult(x, y);
1770 }
1771}
1772
1773void NCL2FaceRestriction::SingleValuedNonconformingTransposeInterpolation(
1774 const Vector& x) const
1775{
1776 MFEM_ASSERT(
1777 m == L2FaceValues::SingleValued,
1778 "This method should be called when m == L2FaceValues::SingleValued.");
1779 if (x_interp.Size()==0)
1780 {
1781 x_interp.SetSize(x.Size());
1782 }
1783 x_interp = x;
1784 SingleValuedNonconformingTransposeInterpolationInPlace(x_interp);
1785}
1786
1787
1788void NCL2FaceRestriction::SingleValuedNonconformingTransposeInterpolationInPlace(
1789 Vector& x) const
1790{
1791 // Assumes all elements have the same number of dofs
1792 const int nface_dofs = face_dofs;
1793 const int vd = vdim;
1794 // Interpolation
1795 auto d_x = Reshape(x_interp.ReadWrite(), nface_dofs, vd, nf);
1796 auto &nc_interp_config = interpolations.GetNCFaceInterpConfig();
1797 const int num_nc_faces = nc_interp_config.Size();
1798 if ( num_nc_faces == 0 ) { return; }
1799 auto interp_config_ptr = nc_interp_config.Read();
1800 auto interpolators = interpolations.GetInterpolators().Read();
1801 const int nc_size = interpolations.GetNumInterpolators();
1802 auto d_interp = Reshape(interpolators, nface_dofs, nface_dofs, nc_size);
1803 static constexpr int max_nd = 16*16;
1804 MFEM_VERIFY(nface_dofs<=max_nd, "Too many degrees of freedom.");
1805 mfem::forall_2D(num_nc_faces, nface_dofs, 1, [=] MFEM_HOST_DEVICE (int nc_face)
1806 {
1807 MFEM_SHARED real_t dof_values[max_nd];
1808 const NCInterpConfig conf = interp_config_ptr[nc_face];
1809 const int master_side = conf.master_side;
1810 const int interp_index = conf.index;
1811 const int face = conf.face_index;
1812 if ( conf.is_non_conforming && master_side==0 )
1813 {
1814 // Interpolation from fine to coarse
1815 for (int c = 0; c < vd; ++c)
1816 {
1817 MFEM_FOREACH_THREAD(dof,x,nface_dofs)
1818 {
1819 dof_values[dof] = d_x(dof, c, face);
1820 }
1821 MFEM_SYNC_THREAD;
1822 MFEM_FOREACH_THREAD(dof_out,x,nface_dofs)
1823 {
1824 real_t res = 0.0;
1825 for (int dof_in = 0; dof_in<nface_dofs; dof_in++)
1826 {
1827 res += d_interp(dof_in, dof_out, interp_index)*dof_values[dof_in];
1828 }
1829 d_x(dof_out, c, face) = res;
1830 }
1831 MFEM_SYNC_THREAD;
1832 }
1833 }
1834 });
1835}
1836
1837void NCL2FaceRestriction::DoubleValuedNonconformingTransposeInterpolation(
1838 const Vector& x) const
1839{
1840 MFEM_ASSERT(
1841 m == L2FaceValues::DoubleValued,
1842 "This method should be called when m == L2FaceValues::DoubleValued.");
1843 if (x_interp.Size()==0)
1844 {
1845 x_interp.SetSize(x.Size());
1846 }
1847 x_interp = x;
1848 DoubleValuedNonconformingTransposeInterpolationInPlace(x_interp);
1849}
1850
1851void NCL2FaceRestriction::DoubleValuedNonconformingTransposeInterpolationInPlace(
1852 Vector& x) const
1853{
1854 // Assumes all elements have the same number of dofs
1855 const int nface_dofs = face_dofs;
1856 const int vd = vdim;
1857 // Interpolation
1858 auto d_x = Reshape(x.ReadWrite(), nface_dofs, vd, 2, nf);
1859 auto &nc_interp_config = interpolations.GetNCFaceInterpConfig();
1860 const int num_nc_faces = nc_interp_config.Size();
1861 if ( num_nc_faces == 0 ) { return; }
1862 auto interp_config_ptr = nc_interp_config.Read();
1863 auto interpolators = interpolations.GetInterpolators().Read();
1864 const int nc_size = interpolations.GetNumInterpolators();
1865 auto d_interp = Reshape(interpolators, nface_dofs, nface_dofs, nc_size);
1866 static constexpr int max_nd = 16*16;
1867 MFEM_VERIFY(nface_dofs<=max_nd, "Too many degrees of freedom.");
1868 mfem::forall_2D(num_nc_faces, nface_dofs, 1, [=] MFEM_HOST_DEVICE (int nc_face)
1869 {
1870 MFEM_SHARED real_t dof_values[max_nd];
1871 const NCInterpConfig conf = interp_config_ptr[nc_face];
1872 const int master_side = conf.master_side;
1873 const int interp_index = conf.index;
1874 const int face = conf.face_index;
1875 if ( conf.is_non_conforming )
1876 {
1877 // Interpolation from fine to coarse
1878 for (int c = 0; c < vd; ++c)
1879 {
1880 MFEM_FOREACH_THREAD(dof,x,nface_dofs)
1881 {
1882 dof_values[dof] = d_x(dof, c, master_side, face);
1883 }
1884 MFEM_SYNC_THREAD;
1885 MFEM_FOREACH_THREAD(dof_out,x,nface_dofs)
1886 {
1887 real_t res = 0.0;
1888 for (int dof_in = 0; dof_in<nface_dofs; dof_in++)
1889 {
1890 res += d_interp(dof_in, dof_out, interp_index)*dof_values[dof_in];
1891 }
1892 d_x(dof_out, c, master_side, face) = res;
1893 }
1894 MFEM_SYNC_THREAD;
1895 }
1896 }
1897 });
1898}
1899
1900void NCL2FaceRestriction::AddMultTranspose(const Vector& x, Vector& y,
1901 const real_t a) const
1902{
1903 MFEM_VERIFY(a == 1.0, "General coefficient case is not yet supported!");
1904 if (nf==0) { return; }
1905 if (type==FaceType::Interior)
1906 {
1907 if ( m==L2FaceValues::DoubleValued )
1908 {
1909 DoubleValuedNonconformingTransposeInterpolation(x);
1910 DoubleValuedConformingAddMultTranspose(x_interp, y);
1911 }
1912 else if ( m==L2FaceValues::SingleValued )
1913 {
1914 SingleValuedNonconformingTransposeInterpolation(x);
1915 SingleValuedConformingAddMultTranspose(x_interp, y);
1916 }
1917 }
1918 else
1919 {
1920 if ( m==L2FaceValues::DoubleValued )
1921 {
1922 DoubleValuedConformingAddMultTranspose(x, y);
1923 }
1924 else if ( m==L2FaceValues::SingleValued )
1925 {
1926 SingleValuedConformingAddMultTranspose(x, y);
1927 }
1928 }
1929}
1930
1931void NCL2FaceRestriction::AddMultTransposeInPlace(Vector& x, Vector& y) const
1932{
1933 if (nf==0) { return; }
1934 if (type==FaceType::Interior)
1935 {
1936 if ( m==L2FaceValues::DoubleValued )
1937 {
1938 DoubleValuedNonconformingTransposeInterpolationInPlace(x);
1939 DoubleValuedConformingAddMultTranspose(x, y);
1940 }
1941 else if ( m==L2FaceValues::SingleValued )
1942 {
1943 SingleValuedNonconformingTransposeInterpolationInPlace(x);
1944 SingleValuedConformingAddMultTranspose(x, y);
1945 }
1946 }
1947 else
1948 {
1949 if ( m==L2FaceValues::DoubleValued )
1950 {
1951 DoubleValuedConformingAddMultTranspose(x, y);
1952 }
1953 else if ( m==L2FaceValues::SingleValued )
1954 {
1955 SingleValuedConformingAddMultTranspose(x, y);
1956 }
1957 }
1958}
1959
1960void NCL2FaceRestriction::FillI(SparseMatrix &mat,
1961 const bool keep_nbr_block) const
1962{
1963 const int nface_dofs = face_dofs;
1964 auto d_indices1 = scatter_indices1.Read();
1965 auto d_indices2 = scatter_indices2.Read();
1966 auto I = mat.ReadWriteI();
1967 mfem::forall(nf*nface_dofs, [=] MFEM_HOST_DEVICE (int fdof)
1968 {
1969 const int iE1 = d_indices1[fdof];
1970 const int iE2 = d_indices2[fdof];
1971 AddNnz(iE1,I,nface_dofs);
1972 AddNnz(iE2,I,nface_dofs);
1973 });
1974}
1975
1976void NCL2FaceRestriction::FillJAndData(const Vector &fea_data,
1977 SparseMatrix &mat,
1978 const bool keep_nbr_block) const
1979{
1980 const int nface_dofs = face_dofs;
1981 auto d_indices1 = scatter_indices1.Read();
1982 auto d_indices2 = scatter_indices2.Read();
1983 auto I = mat.ReadWriteI();
1984 auto mat_fea = Reshape(fea_data.Read(), nface_dofs, nface_dofs, 2, nf);
1985 auto J = mat.WriteJ();
1986 auto Data = mat.WriteData();
1987 auto interp_config_ptr = interpolations.GetFaceInterpConfig().Read();
1988 auto interpolators = interpolations.GetInterpolators().Read();
1989 const int nc_size = interpolations.GetNumInterpolators();
1990 auto d_interp = Reshape(interpolators, nface_dofs, nface_dofs, nc_size);
1991 mfem::forall(nf*nface_dofs, [=] MFEM_HOST_DEVICE (int fdof)
1992 {
1993 const int f = fdof/nface_dofs;
1994 const InterpConfig conf = interp_config_ptr[f];
1995 const int master_side = conf.master_side;
1996 const int interp_index = conf.index;
1997 const int iF = fdof%nface_dofs;
1998 const int iE1 = d_indices1[f*nface_dofs+iF];
1999 const int iE2 = d_indices2[f*nface_dofs+iF];
2000 const int offset1 = AddNnz(iE1,I,nface_dofs);
2001 const int offset2 = AddNnz(iE2,I,nface_dofs);
2002 for (int jF = 0; jF < nface_dofs; jF++)
2003 {
2004 const int jE1 = d_indices1[f*nface_dofs+jF];
2005 const int jE2 = d_indices2[f*nface_dofs+jF];
2006 J[offset2+jF] = jE1;
2007 J[offset1+jF] = jE2;
2008 real_t val1 = 0.0;
2009 real_t val2 = 0.0;
2010 if ( conf.is_non_conforming && master_side==0 )
2011 {
2012 for (int kF = 0; kF < nface_dofs; kF++)
2013 {
2014 val1 += mat_fea(kF,iF,0,f) * d_interp(kF, jF, interp_index);
2015 val2 += d_interp(kF, iF, interp_index) * mat_fea(jF,kF,1,f);
2016 }
2017 }
2018 else if ( conf.is_non_conforming && master_side==1 )
2019 {
2020 for (int kF = 0; kF < nface_dofs; kF++)
2021 {
2022 val1 += d_interp(kF, iF, interp_index) * mat_fea(jF,kF,0,f);
2023 val2 += mat_fea(kF,iF,1,f) * d_interp(kF, jF, interp_index);
2024 }
2025 }
2026 else
2027 {
2028 val1 = mat_fea(jF,iF,0,f);
2029 val2 = mat_fea(jF,iF,1,f);
2030 }
2031 Data[offset2+jF] = val1;
2032 Data[offset1+jF] = val2;
2033 }
2034 });
2035}
2036
2037void NCL2FaceRestriction::AddFaceMatricesToElementMatrices(
2038 const Vector &fea_data,
2039 Vector &ea_data)
2040const
2041{
2042 const int nface_dofs = face_dofs;
2043 const int nelem_dofs = elem_dofs;
2044 const int NE = ne;
2045 if (m==L2FaceValues::DoubleValued)
2046 {
2047 auto d_indices1 = scatter_indices1.Read();
2048 auto d_indices2 = scatter_indices2.Read();
2049 auto mat_fea = Reshape(fea_data.Read(), nface_dofs, nface_dofs, 2, nf);
2050 auto mat_ea = Reshape(ea_data.ReadWrite(), nelem_dofs, nelem_dofs, ne);
2051 auto interp_config_ptr = interpolations.GetFaceInterpConfig().Read();
2052 auto interpolators = interpolations.GetInterpolators().Read();
2053 const int nc_size = interpolations.GetNumInterpolators();
2054 auto d_interp = Reshape(interpolators, nface_dofs, nface_dofs, nc_size);
2055 mfem::forall(nf, [=] MFEM_HOST_DEVICE (int f)
2056 {
2057 const InterpConfig conf = interp_config_ptr[f];
2058 const int master_side = conf.master_side;
2059 const int interp_index = conf.index;
2060 const int e1 = d_indices1[f*nface_dofs]/nelem_dofs;
2061 const int e2 = d_indices2[f*nface_dofs]/nelem_dofs;
2062 for (int j = 0; j < nface_dofs; j++)
2063 {
2064 const int jB1 = d_indices1[f*nface_dofs+j]%nelem_dofs;
2065 for (int i = 0; i < nface_dofs; i++)
2066 {
2067 const int iB1 = d_indices1[f*nface_dofs+i]%nelem_dofs;
2068 real_t val = 0.0;
2069 if ( conf.is_non_conforming && master_side==0 )
2070 {
2071 for (int k = 0; k < nface_dofs; k++)
2072 {
2073 for (int l = 0; l < nface_dofs; l++)
2074 {
2075 val += d_interp(l, j, interp_index)
2076 * mat_fea(k,l,0,f)
2077 * d_interp(k, i, interp_index);
2078 }
2079 }
2080 }
2081 else
2082 {
2083 val = mat_fea(i,j,0,f);
2084 }
2085 AtomicAdd(mat_ea(iB1,jB1,e1), val);
2086 }
2087 }
2088 if (e2 < NE)
2089 {
2090 for (int j = 0; j < nface_dofs; j++)
2091 {
2092 const int jB2 = d_indices2[f*nface_dofs+j]%nelem_dofs;
2093 for (int i = 0; i < nface_dofs; i++)
2094 {
2095 const int iB2 = d_indices2[f*nface_dofs+i]%nelem_dofs;
2096 real_t val = 0.0;
2097 if ( conf.is_non_conforming && master_side==1 )
2098 {
2099 for (int k = 0; k < nface_dofs; k++)
2100 {
2101 for (int l = 0; l < nface_dofs; l++)
2102 {
2103 val += d_interp(l, j, interp_index)
2104 * mat_fea(k,l,1,f)
2105 * d_interp(k, i, interp_index);
2106 }
2107 }
2108 }
2109 else
2110 {
2111 val = mat_fea(i,j,1,f);
2112 }
2113 AtomicAdd(mat_ea(iB2,jB2,e2), val);
2114 }
2115 }
2116 }
2117 });
2118 }
2119 else
2120 {
2121 auto d_indices = scatter_indices1.Read();
2122 auto mat_fea = Reshape(fea_data.Read(), nface_dofs, nface_dofs, nf);
2123 auto mat_ea = Reshape(ea_data.ReadWrite(), nelem_dofs, nelem_dofs, ne);
2124 auto interp_config_ptr = interpolations.GetFaceInterpConfig().Read();
2125 auto interpolators = interpolations.GetInterpolators().Read();
2126 const int nc_size = interpolations.GetNumInterpolators();
2127 auto d_interp = Reshape(interpolators, nface_dofs, nface_dofs, nc_size);
2128 mfem::forall(nf, [=] MFEM_HOST_DEVICE (int f)
2129 {
2130 const InterpConfig conf = interp_config_ptr[f];
2131 const int master_side = conf.master_side;
2132 const int interp_index = conf.index;
2133 const int e = d_indices[f*nface_dofs]/nelem_dofs;
2134 for (int j = 0; j < nface_dofs; j++)
2135 {
2136 const int jE = d_indices[f*nface_dofs+j]%nelem_dofs;
2137 for (int i = 0; i < nface_dofs; i++)
2138 {
2139 const int iE = d_indices[f*nface_dofs+i]%nelem_dofs;
2140 real_t val = 0.0;
2141 if ( conf.is_non_conforming && master_side==0 )
2142 {
2143 for (int k = 0; k < nface_dofs; k++)
2144 {
2145 for (int l = 0; l < nface_dofs; l++)
2146 {
2147 val += d_interp(l, j, interp_index)
2148 * mat_fea(k,l,f)
2149 * d_interp(k, i, interp_index);
2150 }
2151 }
2152 }
2153 else
2154 {
2155 val = mat_fea(i,j,f);
2156 }
2157 AtomicAdd(mat_ea(iE,jE,e), val);
2158 }
2159 }
2160 });
2161 }
2162}
2163
2164int ToLexOrdering(const int dim, const int face_id, const int size1d,
2165 const int index)
2166{
2167 switch (dim)
2168 {
2169 case 1:
2170 return 0;
2171 case 2:
2172 return internal::ToLexOrdering2D(face_id, size1d, index);
2173 case 3:
2174 return internal::ToLexOrdering3D(face_id, size1d, index%size1d, index/size1d);
2175 default:
2176 MFEM_ABORT("Unsupported dimension.");
2177 return 0;
2178 }
2179}
2180
2181void NCL2FaceRestriction::ComputeScatterIndicesAndOffsets()
2182{
2183 Mesh &mesh = *fes.GetMesh();
2184
2185 // Initialization of the offsets
2186 for (int i = 0; i <= ndofs; ++i)
2187 {
2188 gather_offsets[i] = 0;
2189 }
2190
2191 // Computation of scatter and offsets indices
2192 int f_ind=0;
2193 for (int f = 0; f < fes.GetNF(); ++f)
2194 {
2195 Mesh::FaceInformation face = mesh.GetFaceInformation(f);
2196 if ( face.IsNonconformingCoarse() )
2197 {
2198 // We skip nonconforming coarse faces as they are treated
2199 // by the corresponding nonconforming fine faces.
2200 continue;
2201 }
2202 else if ( type==FaceType::Interior && face.IsInterior() )
2203 {
2204 SetFaceDofsScatterIndices1(face,f_ind);
2205 if ( m==L2FaceValues::DoubleValued )
2206 {
2207 PermuteAndSetFaceDofsScatterIndices2(face,f_ind);
2208 }
2209 if ( face.IsConforming() )
2210 {
2211 interpolations.RegisterFaceConformingInterpolation(face,f_ind);
2212 }
2213 else // Non-conforming face
2214 {
2215 interpolations.RegisterFaceCoarseToFineInterpolation(face,f_ind);
2216 }
2217 f_ind++;
2218 }
2219 else if ( type==FaceType::Boundary && face.IsBoundary() )
2220 {
2221 SetFaceDofsScatterIndices1(face,f_ind);
2222 if ( m==L2FaceValues::DoubleValued )
2223 {
2224 SetBoundaryDofsScatterIndices2(face,f_ind);
2225 }
2226 interpolations.RegisterFaceConformingInterpolation(face,f_ind);
2227 f_ind++;
2228 }
2229 }
2230 MFEM_VERIFY(f_ind==nf, "Unexpected number of " <<
2231 (type==FaceType::Interior? "interior" : "boundary") <<
2232 " faces: " << f_ind << " vs " << nf );
2233
2234 // Summation of the offsets
2235 for (int i = 1; i <= ndofs; ++i)
2236 {
2237 gather_offsets[i] += gather_offsets[i - 1];
2238 }
2239
2240 // Transform the interpolation matrix map into a contiguous memory structure.
2241 interpolations.LinearizeInterpolatorMapIntoVector();
2242 interpolations.InitializeNCInterpConfig();
2243}
2244
2245void NCL2FaceRestriction::ComputeGatherIndices()
2246{
2247 Mesh &mesh = *fes.GetMesh();
2248 // Computation of gather_indices
2249 int f_ind = 0;
2250 for (int f = 0; f < fes.GetNF(); ++f)
2251 {
2252 Mesh::FaceInformation face = mesh.GetFaceInformation(f);
2253 MFEM_ASSERT(!face.IsShared(),
2254 "Unexpected shared face in NCL2FaceRestriction.");
2255 if ( face.IsNonconformingCoarse() )
2256 {
2257 // We skip nonconforming coarse faces as they are treated
2258 // by the corresponding nonconforming fine faces.
2259 continue;
2260 }
2261 else if ( face.IsOfFaceType(type) )
2262 {
2263 SetFaceDofsGatherIndices1(face,f_ind);
2264 if ( m==L2FaceValues::DoubleValued &&
2265 type==FaceType::Interior &&
2266 face.IsInterior() )
2267 {
2268 PermuteAndSetFaceDofsGatherIndices2(face,f_ind);
2269 }
2270 f_ind++;
2271 }
2272 }
2273 MFEM_VERIFY(f_ind==nf, "Unexpected number of " <<
2274 (type==FaceType::Interior? "interior" : "boundary") <<
2275 " faces: " << f_ind << " vs " << nf );
2276
2277 // Switch back offsets to their correct value
2278 for (int i = ndofs; i > 0; --i)
2279 {
2280 gather_offsets[i] = gather_offsets[i - 1];
2281 }
2282 gather_offsets[0] = 0;
2283}
2284
2286 const FiniteElementSpace &fes, const Vector &x, FaceType ftype)
2287{
2288#ifdef MFEM_USE_MPI
2289 if (ftype == FaceType::Interior)
2290 {
2291 if (auto *pfes = const_cast<ParFiniteElementSpace*>
2292 (dynamic_cast<const ParFiniteElementSpace*>(&fes)))
2293 {
2294 if (auto *x_gf = const_cast<ParGridFunction*>
2295 (dynamic_cast<const ParGridFunction*>(&x)))
2296 {
2297 Vector &gf_face_nbr = x_gf->FaceNbrData();
2298 if (gf_face_nbr.Size() == 0) { x_gf->ExchangeFaceNbrData(); }
2299 gf_face_nbr.Read();
2300 return Vector(gf_face_nbr, 0, gf_face_nbr.Size());
2301 }
2302 else
2303 {
2304 ParGridFunction gf(pfes, const_cast<Vector&>(x));
2306 return std::move(gf.FaceNbrData());
2307 }
2308 }
2309 }
2310#endif
2311 return Vector();
2312}
2313
2314} // namespace mfem
MFEM_HOST_DEVICE T AtomicAdd(T &add, const T val)
Definition backends.hpp:92
const T * HostRead() const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), false).
Definition array.hpp:321
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition array.hpp:697
int Size() const
Return the logical size of the array.
Definition array.hpp:144
void MakeRef(T *data_, int size_, bool own_data=false)
Make this Array a reference to a pointer.
Definition array.hpp:882
T * GetData()
Returns the data.
Definition array.hpp:118
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Definition array.hpp:317
T * HostReadWrite()
Shortcut for mfem::ReadWrite(a.GetMemory(), a.Size(), false).
Definition array.hpp:337
@ GaussLobatto
Closed type.
Definition fe_base.hpp:32
@ Positive
Bernstein polynomials.
Definition fe_base.hpp:33
Operator that extracts face degrees of freedom for H1, ND, or RT FiniteElementSpaces.
ConformingFaceRestriction(const FiniteElementSpace &fes, const ElementDofOrdering f_ordering, const FaceType type, bool build)
Construct a ConformingFaceRestriction.
void CheckFESpace(const ElementDofOrdering f_ordering)
Verify that ConformingFaceRestriction is built from a supported finite element space.
void SetFaceDofsGatherIndices(const Mesh::FaceInformation &face, const int face_index, const ElementDofOrdering f_ordering)
Set the gathering indices of elem1 for the interior face described by the face.
void AddMultTransposeUnsigned(const Vector &x, Vector &y, const real_t a=1.0) const override
Gather the degrees of freedom, i.e. goes from face E-Vector to L-Vector not taking into account signs...
const FiniteElementSpace & fes
void AddMultTranspose(const Vector &x, Vector &y, const real_t a=1.0) const override
Gather the degrees of freedom, i.e. goes from face E-Vector to L-Vector.
void SetFaceDofsScatterIndices(const Mesh::FaceInformation &face, const int face_index, const ElementDofOrdering f_ordering)
Set the scattering indices of elem1, and increment the offsets for the face described by the face.
void Mult(const Vector &x, Vector &y) const override
Scatter the degrees of freedom, i.e. goes from L-Vector to face E-Vector.
Data type dense matrix using column-major storage.
Definition densemat.hpp:24
ElementRestriction(const FiniteElementSpace &, ElementDofOrdering)
const FiniteElementSpace & fes
void Mult(const Vector &x, Vector &y) const override
Operator application: y=A(x).
void FillSparseMatrix(const Vector &mat_ea, SparseMatrix &mat) const
Fill a Sparse Matrix with Element Matrices.
void AddMultTranspose(const Vector &x, Vector &y, const real_t a=1.0) const override
Add the E-vector degrees of freedom x to the L-vector degrees of freedom y.
void FillJAndData(const Vector &ea_data, SparseMatrix &mat) const
void MultTranspose(const Vector &x, Vector &y) const override
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
void MultTransposeUnsigned(const Vector &x, Vector &y) const
Compute MultTranspose without applying signs based on DOF orientations.
void MultLeftInverse(const Vector &x, Vector &y) const
void BooleanMask(Vector &y) const
Fills the E-vector y with boolean values 0.0 and 1.0 such that each each entry of the L-vector is uni...
void MultUnsigned(const Vector &x, Vector &y) const
Compute Mult without applying signs based on DOF orientations.
int FillI(SparseMatrix &mat) const
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:220
const Table & GetElementToDofTable() const
Return a reference to the internal Table that stores the lists of scalar dofs, for each mesh element,...
Definition fespace.hpp:1136
int GetNF() const
Returns number of faces (i.e. co-dimension 1 entities) in the mesh.
Definition fespace.hpp:746
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
Definition fespace.cpp:3168
const FiniteElement * GetTraceElement(int i, Geometry::Type geom_type) const
Return the trace element from element 'i' to the given 'geom_type'.
Definition fespace.cpp:3276
const FiniteElement * GetFaceElement(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th face in the ...
Definition fespace.cpp:3237
Mesh * GetMesh() const
Returns the mesh.
Definition fespace.hpp:559
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition fespace.hpp:713
Abstract class for all finite elements.
Definition fe_base.hpp:239
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Definition fe_base.hpp:333
virtual void GetFaceMap(const int face_id, Array< int > &face_map) const
Return the mapping from lexicographic face DOFs to lexicographic element DOFs for the given local fac...
Definition fe_base.cpp:507
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition fe_base.hpp:329
std::pair< const DenseMatrix *, int > Key
const Array< NCInterpConfig > & GetNCFaceInterpConfig() const
Return an array containing the interpolation configuration for each face registered with RegisterFace...
void RegisterFaceCoarseToFineInterpolation(const Mesh::FaceInformation &face, int face_index)
Register the face with face and index face_index as a conforming face for the interpolation of the de...
int GetNumInterpolators() const
Return the total number of interpolators.
Array< NCInterpConfig > nc_interp_config
const Vector & GetInterpolators() const
Return an mfem::Vector containing the interpolators in the following format: face_dofs x face_dofs x ...
void LinearizeInterpolatorMapIntoVector()
Transform the interpolation matrix map into a contiguous memory structure.
const FiniteElementSpace & fes
const ElementDofOrdering ordering
Array< InterpConfig > interp_config
void RegisterFaceConformingInterpolation(const Mesh::FaceInformation &face, int face_index)
Register the face with face and index face_index as a conforming face for the interpolation of the de...
A standard isoparametric element transformation.
Definition eltrans.hpp:363
L2ElementRestriction(const FiniteElementSpace &)
void AddMultTranspose(const Vector &x, Vector &y, const real_t a=1.0) const override
Add the E-vector degrees of freedom x to the L-vector degrees of freedom y.
void FillJAndData(const Vector &ea_data, SparseMatrix &mat) const
void MultTranspose(const Vector &x, Vector &y) const override
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
void FillI(SparseMatrix &mat) const
void Mult(const Vector &x, Vector &y) const override
Operator application: y=A(x).
Operator that extracts Face degrees of freedom for L2 spaces.
void PermuteAndSetFaceDofsGatherIndices2(const Mesh::FaceInformation &face, const int face_index)
Permute and set the gathering indices of elem2 for the interior face described by the face....
virtual void FillI(SparseMatrix &mat, const bool keep_nbr_block=false) const
Fill the I array of SparseMatrix corresponding to the sparsity pattern given by this L2FaceRestrictio...
Array< int > scatter_indices2
void SingleValuedConformingMult(const Vector &x, Vector &y) const
Scatter the degrees of freedom, i.e. goes from L-Vector to face E-Vector. Should only be used with co...
L2FaceRestriction(const FiniteElementSpace &fes, const ElementDofOrdering f_ordering, const FaceType type, const L2FaceValues m, bool build)
Constructs an L2FaceRestriction.
void NormalDerivativeMult(const Vector &x, Vector &y) const override
Scatter the degrees of freedom, i.e. goes from L-Vector to face E-Vector.
void AddMultTranspose(const Vector &x, Vector &y, const real_t a=1.0) const override
Gather the degrees of freedom, i.e. goes from face E-Vector to L-Vector.
void NormalDerivativeAddMultTranspose(const Vector &x, Vector &y) const override
Add the face reference-normal derivative degrees of freedom in x to the element degrees of freedom in...
void Mult(const Vector &x, Vector &y) const override
Scatter the degrees of freedom, i.e. goes from L-Vector to face E-Vector.
void SingleValuedConformingAddMultTranspose(const Vector &x, Vector &y) const
Gather the degrees of freedom, i.e. goes from face E-Vector to L-Vector. Should only be used with con...
std::unique_ptr< L2NormalDerivativeFaceRestriction > normal_deriv_restr
void PermuteAndSetSharedFaceDofsScatterIndices2(const Mesh::FaceInformation &face, const int face_index)
Permute and set the scattering indices of elem2 for the shared face described by the face....
void CheckFESpace()
Verify that L2FaceRestriction is built from an L2 FESpace.
virtual void DoubleValuedConformingMult(const Vector &x, Vector &y) const
Scatter the degrees of freedom, i.e. goes from L-Vector to face E-Vector. Should only be used with co...
void PermuteAndSetFaceDofsScatterIndices2(const Mesh::FaceInformation &face, const int face_index)
Permute and set the scattering indices of elem2, and increment the offsets for the face described by ...
Array< int > scatter_indices1
const L2FaceValues m
const FiniteElementSpace & fes
void SetBoundaryDofsScatterIndices2(const Mesh::FaceInformation &face, const int face_index)
Set the scattering indices of elem2 for the boundary face described by the face.
void SetFaceDofsScatterIndices1(const Mesh::FaceInformation &face, const int face_index)
Set the scattering indices of elem1, and increment the offsets for the face described by the face....
void SetFaceDofsGatherIndices1(const Mesh::FaceInformation &face, const int face_index)
Set the gathering indices of elem1 for the interior face described by the face.
const ElementDofOrdering ordering
void DoubleValuedConformingAddMultTranspose(const Vector &x, Vector &y) const
Gather the degrees of freedom, i.e. goes from face E-Vector to L-Vector. Should only be used with con...
virtual void FillJAndData(const Vector &fea_data, SparseMatrix &mat, const bool keep_nbr_block=false) const
Fill the J and Data arrays of the SparseMatrix corresponding to the sparsity pattern given by this L2...
virtual void AddFaceMatricesToElementMatrices(const Vector &fea_data, Vector &ea_data) const
This methods adds the DG face matrices to the element matrices.
Class to compute face normal derivatives (in reference coordinate) of an L2 grid function (used inter...
MemoryType GetMemoryType() const
Return a MemoryType that is currently valid. If both the host and the device pointers are currently v...
void New(int size)
Allocate host memory for size entries with the current host memory type returned by MemoryManager::Ge...
Mesh data type.
Definition mesh.hpp:56
Geometry::Type GetFaceGeometry(int i) const
Return the Geometry::Type associated with face i.
Definition mesh.cpp:1452
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1160
FaceInformation GetFaceInformation(int f) const
Definition mesh.cpp:1169
Operator that extracts face degrees of freedom for L2 nonconforming spaces.
InterpolationManager interpolations
virtual void DoubleValuedNonconformingMult(const Vector &x, Vector &y) const
Scatter the degrees of freedom, i.e. goes from L-Vector to face E-Vector. Should only be used with no...
void DoubleValuedNonconformingInterpolation(Vector &x) const
Apply a change of basis from coarse element basis to fine element basis for the coarse face dofs.
NCL2FaceRestriction(const FiniteElementSpace &fes, const ElementDofOrdering f_ordering, const FaceType type, const L2FaceValues m, bool build)
Constructs an NCL2FaceRestriction, this is a specialization of a L2FaceRestriction for nonconforming ...
int width
Dimension of the input / number of columns in the matrix.
Definition operator.hpp:28
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition operator.hpp:66
int height
Dimension of the output / number of rows in the matrix.
Definition operator.hpp:27
The ordering method used when the number of unknowns per mesh node (vector dimension) is bigger than ...
Definition fespace.hpp:30
Abstract parallel finite element space.
Definition pfespace.hpp:29
void GetFaceNbrElementVDofs(int i, Array< int > &vdofs, DofTransformation &doftrans) const
Class for parallel grid function.
Definition pgridfunc.hpp:33
Data type sparse matrix.
Definition sparsemat.hpp:51
int * ReadWriteI(bool on_dev=true)
int * WriteJ(bool on_dev=true)
Memory< int > & GetMemoryI()
int * WriteI(bool on_dev=true)
Memory< int > & GetMemoryJ()
Memory< real_t > & GetMemoryData()
real_t * WriteData(bool on_dev=true)
int * GetJ()
Definition table.hpp:114
const Array< int > & GetDofMap() const
Get an Array<int> that maps lexicographically ordered indices to the indices of the respective nodes/...
Definition fe_base.hpp:1222
Vector data type.
Definition vector.hpp:80
virtual const real_t * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Definition vector.hpp:474
virtual real_t * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
Definition vector.hpp:490
int Size() const
Returns the size of the vector.
Definition vector.hpp:218
virtual void UseDevice(bool use_dev) const
Enable execution of Vector operations using the mfem::Device.
Definition vector.hpp:136
void SetSize(int s)
Resize the vector to size s.
Definition vector.hpp:538
virtual real_t * HostWrite()
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), false).
Definition vector.hpp:486
virtual real_t * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
Definition vector.hpp:482
int dim
Definition ex24.cpp:53
Mesh * GetMesh(int type)
Definition ex29.cpp:218
int index(int i, int j, int nx, int ny)
Definition life.cpp:235
real_t a
Definition lissajous.cpp:41
real_t f(const Vector &p)
Vector GetLVectorFaceNbrData(const FiniteElementSpace &fes, const Vector &x, FaceType ftype)
Return the face-neighbor data given the L-vector x.
T * Write(Memory< T > &mem, int size, bool on_dev=true)
Get a pointer for write access to mem with the mfem::Device's DeviceMemoryClass, if on_dev = true,...
Definition device.hpp:337
MFEM_HOST_DEVICE DeviceTensor< sizeof...(Dims), T > Reshape(T *ptr, Dims... dims)
Wrap a pointer as a DeviceTensor with automatically deduced template parameters.
Definition dtensor.hpp:131
int ToLexOrdering(const int dim, const int face_id, const int size1d, const int index)
Convert a dof face index from Native ordering to lexicographic ordering for quads and hexes.
void forall_2D(int N, int X, int Y, lambda &&body)
Definition forall.hpp:763
@ byVDIM
VDIM x NQPT x NE (values) / VDIM x DIM x NQPT x NE (grads)
const T & AsConst(const T &a)
Utility function similar to std::as_const in c++17.
Definition array.hpp:360
float real_t
Definition config.hpp:43
int PermuteFaceL2(const int dim, const int face_id1, const int face_id2, const int orientation, const int size1d, const int index)
Compute the dof face index of elem2 corresponding to the given dof face index.
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
Definition fespace.hpp:75
@ LEXICOGRAPHIC
Lexicographic ordering for tensor-product FiniteElements.
std::function< real_t(const Vector &)> f(real_t mass_coeff)
Definition lor_mms.hpp:30
void forall(int N, lambda &&body)
Definition forall.hpp:754
FaceType
Definition mesh.hpp:47
RefCoord t[3]
This structure is used as a human readable output format that deciphers the information contained in ...
Definition mesh.hpp:1894
bool IsBoundary() const
Return true if the face is a boundary face.
Definition mesh.hpp:1935
bool IsNonconformingCoarse() const
Return true if the face is a nonconforming coarse face.
Definition mesh.hpp:1972
struct mfem::Mesh::FaceInformation::@13 element[2]
bool IsLocal() const
Return true if the face is a local interior face which is NOT a master nonconforming face.
Definition mesh.hpp:1912
bool IsConforming() const
Return true if the face is a conforming face.
Definition mesh.hpp:1955
ElementConformity conformity
Definition mesh.hpp:1900
bool IsShared() const
Return true if the face is a shared interior face which is NOT a master nonconforming face.
Definition mesh.hpp:1919
const DenseMatrix * point_matrix
Definition mesh.hpp:1908