MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
transfer.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 "transfer.hpp"
13#include "bilinearform.hpp"
14#include "../general/forall.hpp"
15
16namespace mfem
17{
18
20 FiniteElementSpace &ran_fes_)
21 : dom_fes(dom_fes_), ran_fes(ran_fes_),
22 oper_type(Operator::ANY_TYPE),
23 fw_t_oper(), bw_t_oper()
24{
25#ifdef MFEM_USE_MPI
26 const bool par_dom = dynamic_cast<ParFiniteElementSpace*>(&dom_fes);
27 const bool par_ran = dynamic_cast<ParFiniteElementSpace*>(&ran_fes);
28 MFEM_VERIFY(par_dom == par_ran, "the domain and range FE spaces must both"
29 " be either serial or parallel");
30 parallel = par_dom;
31#endif
32}
33
35 FiniteElementSpace &fes_in, FiniteElementSpace &fes_out,
36 const Operator &oper, OperatorHandle &t_oper)
37{
38 if (t_oper.Ptr())
39 {
40 return *t_oper.Ptr();
41 }
42
43 if (!Parallel())
44 {
45 const SparseMatrix *in_cP = fes_in.GetConformingProlongation();
46 const SparseMatrix *out_cR = fes_out.GetConformingRestriction();
48 {
49 const SparseMatrix *mat = dynamic_cast<const SparseMatrix *>(&oper);
50 MFEM_VERIFY(mat != NULL, "Operator is not a SparseMatrix");
51 if (!out_cR)
52 {
53 t_oper.Reset(const_cast<SparseMatrix*>(mat), false);
54 }
55 else
56 {
57 t_oper.Reset(mfem::Mult(*out_cR, *mat));
58 }
59 if (in_cP)
60 {
61 t_oper.Reset(mfem::Mult(*t_oper.As<SparseMatrix>(), *in_cP));
62 }
63 }
64 else if (oper_type == Operator::ANY_TYPE)
65 {
66 const int RP_case = bool(out_cR) + 2*bool(in_cP);
67 switch (RP_case)
68 {
69 case 0:
70 t_oper.Reset(const_cast<Operator*>(&oper), false);
71 break;
72 case 1:
73 t_oper.Reset(
74 new ProductOperator(out_cR, &oper, false, false));
75 break;
76 case 2:
77 t_oper.Reset(
78 new ProductOperator(&oper, in_cP, false, false));
79 break;
80 case 3:
81 t_oper.Reset(
83 out_cR, &oper, in_cP, false, false, false));
84 break;
85 }
86 }
87 else
88 {
89 MFEM_ABORT("Operator::Type is not supported: " << oper_type);
90 }
91 }
92 else // Parallel() == true
93 {
94#ifdef MFEM_USE_MPI
96 {
97 const SparseMatrix *out_R = fes_out.GetRestrictionMatrix();
98 const ParFiniteElementSpace *pfes_in =
99 dynamic_cast<const ParFiniteElementSpace *>(&fes_in);
100 const ParFiniteElementSpace *pfes_out =
101 dynamic_cast<const ParFiniteElementSpace *>(&fes_out);
102 const SparseMatrix *sp_mat = dynamic_cast<const SparseMatrix *>(&oper);
103 const HypreParMatrix *hy_mat;
104 if (sp_mat)
105 {
106 SparseMatrix *RA = mfem::Mult(*out_R, *sp_mat);
107 t_oper.Reset(pfes_in->Dof_TrueDof_Matrix()->
108 LeftDiagMult(*RA, pfes_out->GetTrueDofOffsets()));
109 delete RA;
110 }
111 else if ((hy_mat = dynamic_cast<const HypreParMatrix *>(&oper)))
112 {
113 HypreParMatrix *RA =
114 hy_mat->LeftDiagMult(*out_R, pfes_out->GetTrueDofOffsets());
115 t_oper.Reset(mfem::ParMult(RA, pfes_in->Dof_TrueDof_Matrix()));
116 delete RA;
117 }
118 else
119 {
120 MFEM_ABORT("unknown Operator type");
121 }
122 }
123 else if (oper_type == Operator::ANY_TYPE)
124 {
125 const Operator *out_R = fes_out.GetRestrictionOperator();
126 t_oper.Reset(new TripleProductOperator(
127 out_R, &oper, fes_in.GetProlongationMatrix(),
128 false, false, false));
129 }
130 else
131 {
132 MFEM_ABORT("Operator::Type is not supported: " << oper_type);
133 }
134#endif
135 }
136
137 return *t_oper.Ptr();
138}
139
140
145
147 BilinearFormIntegrator *mass_integ_, bool own_mass_integ_)
148{
149 if (own_mass_integ) { delete mass_integ; }
150
151 mass_integ = mass_integ_;
152 own_mass_integ = own_mass_integ_;
153}
154
156{
157 if (F.Ptr())
158 {
159 return *F.Ptr();
160 }
161
162 // Construct F
164 {
166 }
168 {
169 Mesh::GeometryList elem_geoms(*ran_fes.GetMesh());
170
172 for (int i = 0; i < elem_geoms.Size(); i++)
173 {
175 localP[elem_geoms[i]]);
176 }
180 }
181 else
182 {
183 MFEM_ABORT("Operator::Type is not supported: " << oper_type);
184 }
185
186 return *F.Ptr();
187}
188
190{
191 if (B.Ptr())
192 {
193 return *B.Ptr();
194 }
195
196 // Construct B, if not set, define a suitable mass_integ
197 if (!mass_integ && ran_fes.GetNE() > 0)
198 {
199 const FiniteElement *f_fe_0 = ran_fes.GetFE(0);
200 const int map_type = f_fe_0->GetMapType();
201 if (map_type == FiniteElement::VALUE ||
202 map_type == FiniteElement::INTEGRAL)
203 {
205 }
206 else if (map_type == FiniteElement::H_DIV ||
207 map_type == FiniteElement::H_CURL)
208 {
210 }
211 else
212 {
213 MFEM_ABORT("unknown type of FE space");
214 }
215 own_mass_integ = true;
216 }
218 {
221 }
222 else
223 {
224 MFEM_ABORT("Operator::Type is not supported: " << oper_type);
225 }
226
227 return *B.Ptr();
228}
229
230
232 const FiniteElementSpace &fes_ho_, const FiniteElementSpace &fes_lor_)
233 : Operator(fes_lor_.GetVSize(), fes_ho_.GetVSize()),
234 fes_ho(fes_ho_),
235 fes_lor(fes_lor_)
236{ }
237
239 int nel_ho, int nel_lor, const CoarseFineTransformations& cf_tr)
240{
241 // Construct the mapping from HO to LOR
242 // ho2lor.GetRow(iho) will give all the LOR elements contained in iho
243 ho2lor.MakeI(nel_ho);
244 for (int ilor = 0; ilor < nel_lor; ++ilor)
245 {
246 int iho = cf_tr.embeddings[ilor].parent;
247 ho2lor.AddAColumnInRow(iho);
248 }
249 ho2lor.MakeJ();
250 for (int ilor = 0; ilor < nel_lor; ++ilor)
251 {
252 int iho = cf_tr.embeddings[ilor].parent;
253 ho2lor.AddConnection(iho, ilor);
254 }
255 ho2lor.ShiftUpI();
256}
257
259 Geometry::Type geom, const FiniteElement& fe_ho,
260 const FiniteElement& fe_lor, ElementTransformation* el_tr,
262 DenseMatrix& M_mixed_el) const
263{
264 int order = fe_lor.GetOrder() + fe_ho.GetOrder() + el_tr->OrderW();
265 const IntegrationRule* ir = &IntRules.Get(geom, order);
266 M_mixed_el = 0.0;
267 for (int i = 0; i < ir->GetNPoints(); i++)
268 {
269 const IntegrationPoint& ip_lor = ir->IntPoint(i);
270 IntegrationPoint ip_ho;
271 ip_tr.Transform(ip_lor, ip_ho);
272 Vector shape_lor(fe_lor.GetDof());
273 fe_lor.CalcShape(ip_lor, shape_lor);
274 Vector shape_ho(fe_ho.GetDof());
275 fe_ho.CalcShape(ip_ho, shape_ho);
276 el_tr->SetIntPoint(&ip_lor);
277 // For now we use the geometry information from the LOR space, which means
278 // we won't be mass conservative if the mesh is curved
279 real_t w = el_tr->Weight() * ip_lor.weight;
280 shape_lor *= w;
281 AddMultVWt(shape_lor, shape_ho, M_mixed_el);
282 }
283}
284
286 const FiniteElementSpace &fes_ho_, const FiniteElementSpace &fes_lor_)
287 : L2Projection(fes_ho_, fes_lor_)
288{
289 Mesh *mesh_ho = fes_ho.GetMesh();
290 Mesh *mesh_lor = fes_lor.GetMesh();
291 int nel_ho = mesh_ho->GetNE();
292 int nel_lor = mesh_lor->GetNE();
293
294 // The prolongation operation is only well-defined when the LOR space has at
295 // least as many DOFs as the high-order space.
296 const bool build_P = fes_lor.GetTrueVSize() >= fes_ho.GetTrueVSize();
297
298 // If the local mesh is empty, skip all computations
299 if (nel_ho == 0) { return; }
300
301 const CoarseFineTransformations &cf_tr = mesh_lor->GetRefinementTransforms();
302
303 int nref_max = 0;
305 mesh_ho->GetGeometries(mesh_ho->Dimension(), geoms);
306 for (int ig = 0; ig < geoms.Size(); ++ig)
307 {
308 Geometry::Type geom = geoms[ig];
309 nref_max = std::max(nref_max, cf_tr.point_matrices[geom].SizeK());
310 }
311
312 BuildHo2Lor(nel_ho, nel_lor, cf_tr);
313
314 offsets.SetSize(nel_ho+1);
315 offsets[0] = 0;
316 for (int iho = 0; iho < nel_ho; ++iho)
317 {
318 int nref = ho2lor.RowSize(iho);
319 const FiniteElement &fe_ho = *fes_ho.GetFE(iho);
320 const FiniteElement &fe_lor = *fes_lor.GetFE(ho2lor.GetRow(iho)[0]);
321 offsets[iho+1] = offsets[iho] + fe_ho.GetDof()*fe_lor.GetDof()*nref;
322 }
323 // R will contain the restriction (L^2 projection operator) defined on each
324 // coarse HO element (and corresponding patch of LOR elements)
325 R.SetSize(offsets[nel_ho]);
326 if (build_P)
327 {
328 // P will contain the corresponding prolongation operator
329 P.SetSize(offsets[nel_ho]);
330 }
331
333 IsoparametricTransformation &emb_tr = ip_tr.Transf;
334
335 for (int iho = 0; iho < nel_ho; ++iho)
336 {
337 Array<int> lor_els;
338 ho2lor.GetRow(iho, lor_els);
339 int nref = ho2lor.RowSize(iho);
340
341 Geometry::Type geom = mesh_ho->GetElementBaseGeometry(iho);
342 const FiniteElement &fe_ho = *fes_ho.GetFE(iho);
343 const FiniteElement &fe_lor = *fes_lor.GetFE(lor_els[0]);
344 int ndof_ho = fe_ho.GetDof();
345 int ndof_lor = fe_lor.GetDof();
346
347 emb_tr.SetIdentityTransformation(geom);
348 const DenseTensor &pmats = cf_tr.point_matrices[geom];
349
350 DenseMatrix R_iho(&R[offsets[iho]], ndof_lor*nref, ndof_ho);
351
352 DenseMatrix Minv_lor(ndof_lor*nref, ndof_lor*nref);
353 DenseMatrix M_mixed(ndof_lor*nref, ndof_ho);
354
356 DenseMatrix M_lor_el(ndof_lor, ndof_lor);
357 DenseMatrixInverse Minv_lor_el(&M_lor_el);
358 DenseMatrix M_lor(ndof_lor*nref, ndof_lor*nref);
359 DenseMatrix M_mixed_el(ndof_lor, ndof_ho);
360
361 Minv_lor = 0.0;
362 M_lor = 0.0;
363
364 DenseMatrix RtMlor(ndof_ho, ndof_lor*nref);
365 DenseMatrix RtMlorR(ndof_ho, ndof_ho);
366 DenseMatrixInverse RtMlorR_inv(&RtMlorR);
367
368 for (int iref = 0; iref < nref; ++iref)
369 {
370 // Assemble the low-order refined mass matrix and invert locally
371 int ilor = lor_els[iref];
373 mi.AssembleElementMatrix(fe_lor, *el_tr, M_lor_el);
374 M_lor.CopyMN(M_lor_el, iref*ndof_lor, iref*ndof_lor);
375 Minv_lor_el.Factor();
376 Minv_lor_el.GetInverseMatrix(M_lor_el);
377 // Insert into the diagonal of the patch LOR mass matrix
378 Minv_lor.CopyMN(M_lor_el, iref*ndof_lor, iref*ndof_lor);
379
380 // Now assemble the block-row of the mixed mass matrix associated
381 // with integrating HO functions against LOR functions on the LOR
382 // sub-element.
383
384 // Create the transformation that embeds the fine low-order element
385 // within the coarse high-order element in reference space
386 emb_tr.SetPointMat(pmats(cf_tr.embeddings[ilor].matrix));
387
388 ElemMixedMass(geom, fe_ho, fe_lor, el_tr, ip_tr, M_mixed_el);
389
390 M_mixed.CopyMN(M_mixed_el, iref*ndof_lor, 0);
391 }
392 mfem::Mult(Minv_lor, M_mixed, R_iho);
393
394 if (build_P)
395 {
396 DenseMatrix P_iho(&P[offsets[iho]], ndof_ho, ndof_lor*nref);
397
398 mfem::MultAtB(R_iho, M_lor, RtMlor);
399 mfem::Mult(RtMlor, R_iho, RtMlorR);
400 RtMlorR_inv.Factor();
401 RtMlorR_inv.Mult(RtMlor, P_iho);
402 }
403 }
404}
405
407 const Vector &x, Vector &y) const
408{
409 int vdim = fes_ho.GetVDim();
410 Array<int> vdofs;
411 DenseMatrix xel_mat, yel_mat;
412 for (int iho = 0; iho < fes_ho.GetNE(); ++iho)
413 {
414 int nref = ho2lor.RowSize(iho);
415 int ndof_ho = fes_ho.GetFE(iho)->GetDof();
416 int ndof_lor = fes_lor.GetFE(ho2lor.GetRow(iho)[0])->GetDof();
417 xel_mat.SetSize(ndof_ho, vdim);
418 yel_mat.SetSize(ndof_lor*nref, vdim);
419 DenseMatrix R_iho(&R[offsets[iho]], ndof_lor*nref, ndof_ho);
420
421 fes_ho.GetElementVDofs(iho, vdofs);
422 x.GetSubVector(vdofs, xel_mat.GetData());
423 mfem::Mult(R_iho, xel_mat, yel_mat);
424 // Place result correctly into the low-order vector
425 for (int iref = 0; iref < nref; ++iref)
426 {
427 int ilor = ho2lor.GetRow(iho)[iref];
428 for (int vd=0; vd<vdim; ++vd)
429 {
430 fes_lor.GetElementDofs(ilor, vdofs);
431 fes_lor.DofsToVDofs(vd, vdofs);
432 y.SetSubVector(vdofs, &yel_mat(iref*ndof_lor,vd));
433 }
434 }
435 }
436}
437
439 const Vector &x, Vector &y) const
440{
441 int vdim = fes_ho.GetVDim();
442 Array<int> vdofs;
443 DenseMatrix xel_mat, yel_mat;
444 y = 0.0;
445 for (int iho = 0; iho < fes_ho.GetNE(); ++iho)
446 {
447 int nref = ho2lor.RowSize(iho);
448 int ndof_ho = fes_ho.GetFE(iho)->GetDof();
449 int ndof_lor = fes_lor.GetFE(ho2lor.GetRow(iho)[0])->GetDof();
450 xel_mat.SetSize(ndof_lor*nref, vdim);
451 yel_mat.SetSize(ndof_ho, vdim);
452 DenseMatrix R_iho(&R[offsets[iho]], ndof_lor*nref, ndof_ho);
453
454 // Extract the LOR DOFs
455 for (int iref=0; iref<nref; ++iref)
456 {
457 int ilor = ho2lor.GetRow(iho)[iref];
458 for (int vd=0; vd<vdim; ++vd)
459 {
460 fes_lor.GetElementDofs(ilor, vdofs);
461 fes_lor.DofsToVDofs(vd, vdofs);
462 x.GetSubVector(vdofs, &xel_mat(iref*ndof_lor, vd));
463 }
464 }
465 // Multiply locally by the transpose
466 mfem::MultAtB(R_iho, xel_mat, yel_mat);
467 // Place the result in the HO vector
468 fes_ho.GetElementVDofs(iho, vdofs);
469 y.AddElementVector(vdofs, yel_mat.GetData());
470 }
471}
472
474 const Vector &x, Vector &y) const
475{
476 if (fes_ho.GetNE() == 0) { return; }
477 MFEM_VERIFY(P.Size() > 0, "Prolongation not supported for these spaces.")
478 int vdim = fes_ho.GetVDim();
479 Array<int> vdofs;
480 DenseMatrix xel_mat,yel_mat;
481 y = 0.0;
482 for (int iho = 0; iho < fes_ho.GetNE(); ++iho)
483 {
484 int nref = ho2lor.RowSize(iho);
485 int ndof_ho = fes_ho.GetFE(iho)->GetDof();
486 int ndof_lor = fes_lor.GetFE(ho2lor.GetRow(iho)[0])->GetDof();
487 xel_mat.SetSize(ndof_lor*nref, vdim);
488 yel_mat.SetSize(ndof_ho, vdim);
489 DenseMatrix P_iho(&P[offsets[iho]], ndof_ho, ndof_lor*nref);
490
491 // Extract the LOR DOFs
492 for (int iref = 0; iref < nref; ++iref)
493 {
494 int ilor = ho2lor.GetRow(iho)[iref];
495 for (int vd = 0; vd < vdim; ++vd)
496 {
497 fes_lor.GetElementDofs(ilor, vdofs);
498 fes_lor.DofsToVDofs(vd, vdofs);
499 x.GetSubVector(vdofs, &xel_mat(iref*ndof_lor, vd));
500 }
501 }
502 // Locally prolongate
503 mfem::Mult(P_iho, xel_mat, yel_mat);
504 // Place the result in the HO vector
505 fes_ho.GetElementVDofs(iho, vdofs);
506 y.AddElementVector(vdofs, yel_mat.GetData());
507 }
508}
509
511 const Vector &x, Vector &y) const
512{
513 if (fes_ho.GetNE() == 0) { return; }
514 MFEM_VERIFY(P.Size() > 0, "Prolongation not supported for these spaces.")
515 int vdim = fes_ho.GetVDim();
516 Array<int> vdofs;
517 DenseMatrix xel_mat,yel_mat;
518 for (int iho = 0; iho < fes_ho.GetNE(); ++iho)
519 {
520 int nref = ho2lor.RowSize(iho);
521 int ndof_ho = fes_ho.GetFE(iho)->GetDof();
522 int ndof_lor = fes_lor.GetFE(ho2lor.GetRow(iho)[0])->GetDof();
523 xel_mat.SetSize(ndof_ho, vdim);
524 yel_mat.SetSize(ndof_lor*nref, vdim);
525 DenseMatrix P_iho(&P[offsets[iho]], ndof_ho, ndof_lor*nref);
526
527 fes_ho.GetElementVDofs(iho, vdofs);
528 x.GetSubVector(vdofs, xel_mat.GetData());
529 mfem::MultAtB(P_iho, xel_mat, yel_mat);
530
531 // Place result correctly into the low-order vector
532 for (int iref = 0; iref < nref; ++iref)
533 {
534 int ilor = ho2lor.GetRow(iho)[iref];
535 for (int vd=0; vd<vdim; ++vd)
536 {
537 fes_lor.GetElementDofs(ilor, vdofs);
538 fes_lor.DofsToVDofs(vd, vdofs);
539 y.SetSubVector(vdofs, &yel_mat(iref*ndof_lor,vd));
540 }
541 }
542 }
543}
544
546 const FiniteElementSpace& fes_ho_, const FiniteElementSpace& fes_lor_)
547 : L2Projection(fes_ho_, fes_lor_)
548{
549 std::unique_ptr<SparseMatrix> R_mat, M_LH_mat;
550 std::tie(R_mat, M_LH_mat) = ComputeSparseRAndM_LH();
551
552 FiniteElementSpace fes_ho_scalar(fes_ho.GetMesh(), fes_ho.FEColl(), 1);
553 FiniteElementSpace fes_lor_scalar(fes_lor.GetMesh(), fes_lor.FEColl(), 1);
554
555 const SparseMatrix *P_ho = fes_ho_scalar.GetConformingProlongation();
556 const SparseMatrix *P_lor = fes_lor_scalar.GetConformingProlongation();
557
558 if (P_ho || P_lor)
559 {
560 if (P_ho && P_lor)
561 {
562 R_mat.reset(RAP(*P_lor, *R_mat, *P_ho));
563 M_LH_mat.reset(RAP(*P_lor, *M_LH_mat, *P_ho));
564 }
565 else if (P_ho)
566 {
567 R_mat.reset(mfem::Mult(*R_mat, *P_ho));
568 M_LH_mat.reset(mfem::Mult(*M_LH_mat, *P_ho));
569 }
570 else // P_lor != nullptr
571 {
572 R_mat.reset(mfem::Mult(*P_lor, *R_mat));
573 M_LH_mat.reset(mfem::Mult(*P_lor, *M_LH_mat));
574 }
575 }
576
577 SparseMatrix *RTxM_LH_mat = TransposeMult(*R_mat, *M_LH_mat);
578 precon.reset(new DSmoother(*RTxM_LH_mat));
579
580 // Set ownership
581 RTxM_LH.reset(RTxM_LH_mat);
582 R = std::move(R_mat);
583 M_LH = std::move(M_LH_mat);
584
585 SetupPCG();
586}
587
588#ifdef MFEM_USE_MPI
589
591 const ParFiniteElementSpace& pfes_ho, const ParFiniteElementSpace& pfes_lor)
592 : L2Projection(pfes_ho, pfes_lor),
593 pcg(pfes_ho.GetComm())
594{
595 std::tie(R, M_LH) = ComputeSparseRAndM_LH();
596
597 ParFiniteElementSpace pfes_ho_scalar(pfes_ho.GetParMesh(),
598 pfes_ho.FEColl(), 1);
599 ParFiniteElementSpace pfes_lor_scalar(pfes_lor.GetParMesh(),
600 pfes_lor.FEColl(), 1);
601
602 HypreParMatrix R_local = HypreParMatrix(pfes_ho.GetComm(),
603 pfes_lor_scalar.GlobalVSize(),
604 pfes_ho_scalar.GlobalVSize(),
605 pfes_lor_scalar.GetDofOffsets(),
606 pfes_ho_scalar.GetDofOffsets(),
607 static_cast<SparseMatrix*>(R.get()));
608 HypreParMatrix M_LH_local = HypreParMatrix(pfes_ho.GetComm(),
609 pfes_lor_scalar.GlobalVSize(),
610 pfes_ho_scalar.GlobalVSize(),
611 pfes_lor_scalar.GetDofOffsets(),
612 pfes_ho_scalar.GetDofOffsets(),
613 static_cast<SparseMatrix*>(M_LH.get()));
614
615 HypreParMatrix *R_mat = RAP(pfes_lor_scalar.Dof_TrueDof_Matrix(),
616 &R_local, pfes_ho_scalar.Dof_TrueDof_Matrix());
617 HypreParMatrix *M_LH_mat = RAP(pfes_lor_scalar.Dof_TrueDof_Matrix(),
618 &M_LH_local, pfes_ho_scalar.Dof_TrueDof_Matrix());
619
620 std::unique_ptr<HypreParMatrix> R_T(R_mat->Transpose());
621 HypreParMatrix *RTxM_LH_mat = ParMult(R_T.get(), M_LH_mat, true);
622
623 HypreBoomerAMG *amg = new HypreBoomerAMG(*RTxM_LH_mat);
624 amg->SetPrintLevel(0);
625
626 R.reset(R_mat);
627 M_LH.reset(M_LH_mat);
628 RTxM_LH.reset(RTxM_LH_mat);
629 precon.reset(amg);
630
631 SetupPCG();
634}
635
636#endif
637
639{
640 // Basic PCG solver setup
641 pcg.SetPrintLevel(0);
642 // pcg.SetPrintLevel(IterativeSolver::PrintLevel().Summary());
643 pcg.SetMaxIter(1000);
644 // initial values for relative and absolute tolerance
645 pcg.SetRelTol(1e-13);
646 pcg.SetAbsTol(1e-13);
647 pcg.SetPreconditioner(*precon);
648 pcg.SetOperator(*RTxM_LH);
649}
650
652 const Vector& x, Vector& y) const
653{
654 Vector X(fes_ho.GetTrueVSize());
655 Vector X_dim(R->Width());
656
657 Vector Y_dim(R->Height());
658 Vector Y(fes_lor.GetTrueVSize());
659
660 Array<int> vdofs_list;
661
662 GetTDofs(fes_ho, x, X);
663
664 for (int d = 0; d < fes_ho.GetVDim(); ++d)
665 {
666 TDofsListByVDim(fes_ho, d, vdofs_list);
667 X.GetSubVector(vdofs_list, X_dim);
668 R->Mult(X_dim, Y_dim);
669 TDofsListByVDim(fes_lor, d, vdofs_list);
670 Y.SetSubVector(vdofs_list, Y_dim);
671 }
672
673 SetFromTDofs(fes_lor, Y, y);
674}
675
677 const Vector& x, Vector& y) const
678{
679 Vector X(fes_lor.GetTrueVSize());
680 Vector X_dim(R->Height());
681
682 Vector Y_dim(R->Width());
683 Vector Y(fes_ho.GetTrueVSize());
684
685 Array<int> vdofs_list;
686
687 GetTDofsTranspose(fes_lor, x, X);
688
689 for (int d = 0; d < fes_ho.GetVDim(); ++d)
690 {
691 TDofsListByVDim(fes_lor, d, vdofs_list);
692 X.GetSubVector(vdofs_list, X_dim);
693 R->MultTranspose(X_dim, Y_dim);
694 TDofsListByVDim(fes_ho, d, vdofs_list);
695 Y.SetSubVector(vdofs_list, Y_dim);
696 }
697
698 SetFromTDofsTranspose(fes_ho, Y, y);
699}
700
702 const Vector& x, Vector& y) const
703{
704 Vector X(fes_lor.GetTrueVSize());
705 Vector X_dim(M_LH->Height());
706 Vector Xbar(pcg.Width());
707
708 Vector Y_dim(pcg.Height());
709 Vector Y(fes_ho.GetTrueVSize());
710
711 Array<int> vdofs_list;
712
713 GetTDofs(fes_lor, x, X);
714
715 for (int d = 0; d < fes_ho.GetVDim(); ++d)
716 {
717 TDofsListByVDim(fes_lor, d, vdofs_list);
718 X.GetSubVector(vdofs_list, X_dim);
719 // Compute y = P x = (R^T M_LH)^(-1) M_LH^T X = (R^T M_LH)^(-1) Xbar
720 M_LH->MultTranspose(X_dim, Xbar);
721 Y_dim = 0.0;
722 pcg.Mult(Xbar, Y_dim);
723 TDofsListByVDim(fes_ho, d, vdofs_list);
724 Y.SetSubVector(vdofs_list, Y_dim);
725 }
726
727 SetFromTDofs(fes_ho, Y, y);
728}
729
731 const Vector& x, Vector& y) const
732{
733 Vector X(fes_ho.GetTrueVSize());
734 Vector X_dim(pcg.Width());
735 Vector Xbar(pcg.Height());
736
737 Vector Y_dim(M_LH->Height());
738 Vector Y(fes_lor.GetTrueVSize());
739
740 Array<int> vdofs_list;
741
742 GetTDofsTranspose(fes_ho, x, X);
743
744 for (int d = 0; d < fes_ho.GetVDim(); ++d)
745 {
746 TDofsListByVDim(fes_ho, d, vdofs_list);
747 X.GetSubVector(vdofs_list, X_dim);
748 // Compute y = P^T x = M_LH (R^T M_LH)^(-1) X = M_LH Xbar
749 Xbar = 0.0;
750 pcg.Mult(X_dim, Xbar);
751 M_LH->Mult(Xbar, Y_dim);
752 TDofsListByVDim(fes_lor, d, vdofs_list);
753 Y.SetSubVector(vdofs_list, Y_dim);
754 }
755
756 SetFromTDofsTranspose(fes_lor, Y, y);
757}
758
760{
761 pcg.SetRelTol(p_rtol_);
762}
763
765{
766 pcg.SetAbsTol(p_atol_);
767}
768
769std::pair<
770std::unique_ptr<SparseMatrix>,
771std::unique_ptr<SparseMatrix>>
773{
774 std::pair<std::unique_ptr<SparseMatrix>,
775 std::unique_ptr<SparseMatrix>> r_and_mlh;
776
777 Mesh* mesh_ho = fes_ho.GetMesh();
778 Mesh* mesh_lor = fes_lor.GetMesh();
779 int nel_ho = mesh_ho->GetNE();
780 int nel_lor = mesh_lor->GetNE();
781 int ndof_lor = fes_lor.GetNDofs();
782
783 // If the local mesh is empty, skip all computations
784 if (nel_ho == 0)
785 {
786 return std::make_pair(
787 std::unique_ptr<SparseMatrix>(new SparseMatrix),
788 std::unique_ptr<SparseMatrix>(new SparseMatrix)
789 );
790 }
791
792 const CoarseFineTransformations& cf_tr = mesh_lor->GetRefinementTransforms();
793
794 int nref_max = 0;
796 mesh_ho->GetGeometries(mesh_ho->Dimension(), geoms);
797 for (int ig = 0; ig < geoms.Size(); ++ig)
798 {
799 Geometry::Type geom = geoms[ig];
800 nref_max = std::max(nref_max, cf_tr.point_matrices[geom].SizeK());
801 }
802
803 BuildHo2Lor(nel_ho, nel_lor, cf_tr);
804
805 // ML_inv contains the inverse lumped (row sum) mass matrix. Note that the
806 // method will also work with a full (consistent) mass matrix, though this is
807 // not implemented here. L refers to the low-order refined mesh
808 Vector ML_inv(ndof_lor);
809 ML_inv = 0.0;
810
811 // Compute ML_inv
812 for (int iho = 0; iho < nel_ho; ++iho)
813 {
814 Array<int> lor_els;
815 ho2lor.GetRow(iho, lor_els);
816 int nref = ho2lor.RowSize(iho);
817
818 Geometry::Type geom = mesh_ho->GetElementBaseGeometry(iho);
819 const FiniteElement& fe_lor = *fes_lor.GetFE(lor_els[0]);
820 int nedof_lor = fe_lor.GetDof();
821
822 // Instead of using a MassIntegrator, manually loop over integration
823 // points so we can row sum and store the diagonal as a Vector.
824 Vector ML_el(nedof_lor);
825 Vector shape_lor(nedof_lor);
826 Array<int> dofs_lor(nedof_lor);
827
828 for (int iref = 0; iref < nref; ++iref)
829 {
830 int ilor = lor_els[iref];
831 ElementTransformation* el_tr = fes_lor.GetElementTransformation(ilor);
832
833 int order = 2 * fe_lor.GetOrder() + el_tr->OrderW();
834 const IntegrationRule* ir = &IntRules.Get(geom, order);
835 ML_el = 0.0;
836 for (int i = 0; i < ir->GetNPoints(); ++i)
837 {
838 const IntegrationPoint& ip_lor = ir->IntPoint(i);
839 fe_lor.CalcShape(ip_lor, shape_lor);
840 el_tr->SetIntPoint(&ip_lor);
841 ML_el += (shape_lor *= (el_tr->Weight() * ip_lor.weight));
842 }
843 fes_lor.GetElementDofs(ilor, dofs_lor);
844 ML_inv.AddElementVector(dofs_lor, ML_el);
845 }
846 }
847 // DOF by DOF inverse of non-zero entries
848 LumpedMassInverse(ML_inv);
849
850 // Compute sparsity pattern for R = M_L^(-1) M_LH and allocate
851 r_and_mlh.first = AllocR();
852 // Allocate M_LH (same sparsity pattern as R)
853 // L refers to the low-order refined mesh (DOFs correspond to rows)
854 // H refers to the higher-order mesh (DOFs correspond to columns)
855 Memory<int> I(r_and_mlh.first->Height() + 1);
856 for (int icol = 0; icol < r_and_mlh.first->Height() + 1; ++icol)
857 {
858 I[icol] = r_and_mlh.first->GetI()[icol];
859 }
860 Memory<int> J(r_and_mlh.first->NumNonZeroElems());
861 for (int jcol = 0; jcol < r_and_mlh.first->NumNonZeroElems(); ++jcol)
862 {
863 J[jcol] = r_and_mlh.first->GetJ()[jcol];
864 }
865 r_and_mlh.second = std::unique_ptr<SparseMatrix>(
866 new SparseMatrix(I, J, NULL, r_and_mlh.first->Height(),
867 r_and_mlh.first->Width(), true, true, true));
868
870 IsoparametricTransformation& emb_tr = ip_tr.Transf;
871
872 // Compute M_LH and R
873 for (int iho = 0; iho < nel_ho; ++iho)
874 {
875 Array<int> lor_els;
876 ho2lor.GetRow(iho, lor_els);
877 int nref = ho2lor.RowSize(iho);
878
879 Geometry::Type geom = mesh_ho->GetElementBaseGeometry(iho);
880 const FiniteElement& fe_ho = *fes_ho.GetFE(iho);
881 const FiniteElement& fe_lor = *fes_lor.GetFE(lor_els[0]);
882
883 emb_tr.SetIdentityTransformation(geom);
884 const DenseTensor& pmats = cf_tr.point_matrices[geom];
885
886 int nedof_ho = fe_ho.GetDof();
887 int nedof_lor = fe_lor.GetDof();
888 DenseMatrix M_LH_el(nedof_lor, nedof_ho);
889 DenseMatrix R_el(nedof_lor, nedof_ho);
890
891 for (int iref = 0; iref < nref; ++iref)
892 {
893 int ilor = lor_els[iref];
894 ElementTransformation* el_tr = fes_lor.GetElementTransformation(ilor);
895
896 // Create the transformation that embeds the fine low-order element
897 // within the coarse high-order element in reference space
898 emb_tr.SetPointMat(pmats(cf_tr.embeddings[ilor].matrix));
899
900 ElemMixedMass(geom, fe_ho, fe_lor, el_tr, ip_tr, M_LH_el);
901
902 Array<int> dofs_lor(nedof_lor);
903 fes_lor.GetElementDofs(ilor, dofs_lor);
904 Vector R_row;
905 for (int i = 0; i < nedof_lor; ++i)
906 {
907 M_LH_el.GetRow(i, R_row);
908 R_el.SetRow(i, R_row.Set(ML_inv[dofs_lor[i]], R_row));
909 }
910 Array<int> dofs_ho(nedof_ho);
911 fes_ho.GetElementDofs(iho, dofs_ho);
912 r_and_mlh.second->AddSubMatrix(dofs_lor, dofs_ho, M_LH_el);
913 r_and_mlh.first->AddSubMatrix(dofs_lor, dofs_ho, R_el);
914 }
915 }
916
917 return r_and_mlh;
918}
919
921 const FiniteElementSpace& fes, const Vector& x, Vector& X) const
922{
923 const Operator* res = fes.GetRestrictionOperator();
924 if (res)
925 {
926 res->Mult(x, X);
927 }
928 else
929 {
930 X = x;
931 }
932}
933
935 const FiniteElementSpace& fes, const Vector &X, Vector& x) const
936{
937 const Operator* P = fes.GetProlongationMatrix();
938 if (P)
939 {
940 P->Mult(X, x);
941 }
942 else
943 {
944 x = X;
945 }
946}
947
949 const FiniteElementSpace& fes, const Vector& x, Vector& X) const
950{
951 const Operator* P = fes.GetProlongationMatrix();
952 if (P)
953 {
954 P->MultTranspose(x, X);
955 }
956 else
957 {
958 X = x;
959 }
960}
961
963 const FiniteElementSpace& fes, const Vector &X, Vector& x) const
964{
965 const Operator *R_op = fes.GetRestrictionOperator();
966 if (R_op)
967 {
968 R_op->MultTranspose(X, x);
969 }
970 else
971 {
972 x = X;
973 }
974}
975
977 const FiniteElementSpace& fes, int vdim, Array<int>& vdofs_list) const
978{
979 const SparseMatrix *R_mat = fes.GetRestrictionMatrix();
980 if (R_mat)
981 {
982 Array<int> x_vdofs_list(fes.GetNDofs());
983 Array<int> x_vdofs_marker(fes.GetVSize());
984 Array<int> X_vdofs_marker(fes.GetTrueVSize());
985 fes.GetVDofs(vdim, x_vdofs_list);
986 FiniteElementSpace::ListToMarker(x_vdofs_list, fes.GetVSize(), x_vdofs_marker);
987 R_mat->BooleanMult(x_vdofs_marker, X_vdofs_marker);
988 FiniteElementSpace::MarkerToList(X_vdofs_marker, vdofs_list);
989 }
990 else
991 {
992 vdofs_list.SetSize(fes.GetNDofs());
993 fes.GetVDofs(vdim, vdofs_list);
994 }
995}
996
998 Vector& ML_inv) const
999{
1000 Vector ML_inv_full(fes_lor.GetVSize());
1001 // set ML_inv on dofs for vdim = 0
1002 Array<int> vdofs_list(fes_lor.GetNDofs());
1003 fes_lor.GetVDofs(0, vdofs_list);
1004 ML_inv_full.SetSubVector(vdofs_list, ML_inv);
1005
1006 Vector ML_inv_true(fes_lor.GetTrueVSize());
1007 const Operator *P = fes_lor.GetProlongationMatrix();
1008 if (P) { P->MultTranspose(ML_inv_full, ML_inv_true); }
1009 else { ML_inv_true = ML_inv_full; }
1010
1011 for (int i = 0; i < ML_inv_true.Size(); ++i)
1012 {
1013 ML_inv_true[i] = 1.0 / ML_inv_true[i];
1014 }
1015
1016 if (P) { P->Mult(ML_inv_true, ML_inv_full); }
1017 else { ML_inv_full = ML_inv_true; }
1018
1019 ML_inv_full.GetSubVector(vdofs_list, ML_inv);
1020}
1021
1022std::unique_ptr<SparseMatrix>
1024{
1025 const Table& elem_dof_ho = fes_ho.GetElementToDofTable();
1026 const Table& elem_dof_lor = fes_lor.GetElementToDofTable();
1027 const int ndof_ho = fes_ho.GetNDofs();
1028 const int ndof_lor = fes_lor.GetNDofs();
1029
1030 Table dof_elem_lor;
1031 Transpose(elem_dof_lor, dof_elem_lor, ndof_lor);
1032
1033 Mesh* mesh_lor = fes_lor.GetMesh();
1034 const CoarseFineTransformations& cf_tr = mesh_lor->GetRefinementTransforms();
1035
1036 // mfem::Mult but uses ho2lor to map HO elements to LOR elements
1037 const int* elem_dof_hoI = elem_dof_ho.GetI();
1038 const int* elem_dof_hoJ = elem_dof_ho.GetJ();
1039 const int* dof_elem_lorI = dof_elem_lor.GetI();
1040 const int* dof_elem_lorJ = dof_elem_lor.GetJ();
1041
1042 Array<int> I(ndof_lor + 1);
1043
1044 // figure out the size of J
1045 Array<int> dof_used_ho;
1046 dof_used_ho.SetSize(ndof_ho, -1);
1047
1048 int sizeJ = 0;
1049 for (int ilor = 0; ilor < ndof_lor; ++ilor)
1050 {
1051 for (int jlor = dof_elem_lorI[ilor]; jlor < dof_elem_lorI[ilor + 1]; ++jlor)
1052 {
1053 int el_lor = dof_elem_lorJ[jlor];
1054 int iho = cf_tr.embeddings[el_lor].parent;
1055 for (int jho = elem_dof_hoI[iho]; jho < elem_dof_hoI[iho + 1]; ++jho)
1056 {
1057 int dof_ho = elem_dof_hoJ[jho];
1058 if (dof_used_ho[dof_ho] != ilor)
1059 {
1060 dof_used_ho[dof_ho] = ilor;
1061 ++sizeJ;
1062 }
1063 }
1064 }
1065 }
1066
1067 // initialize dof_ho_dof_lor
1068 Table dof_lor_dof_ho;
1069 dof_lor_dof_ho.SetDims(ndof_lor, sizeJ);
1070
1071 for (int i = 0; i < ndof_ho; ++i)
1072 {
1073 dof_used_ho[i] = -1;
1074 }
1075
1076 // set values of J
1077 int* dof_dofI = dof_lor_dof_ho.GetI();
1078 int* dof_dofJ = dof_lor_dof_ho.GetJ();
1079 sizeJ = 0;
1080 for (int ilor = 0; ilor < ndof_lor; ++ilor)
1081 {
1082 dof_dofI[ilor] = sizeJ;
1083 for (int jlor = dof_elem_lorI[ilor]; jlor < dof_elem_lorI[ilor + 1]; ++jlor)
1084 {
1085 int el_lor = dof_elem_lorJ[jlor];
1086 int iho = cf_tr.embeddings[el_lor].parent;
1087 for (int jho = elem_dof_hoI[iho]; jho < elem_dof_hoI[iho + 1]; ++jho)
1088 {
1089 int dof_ho = elem_dof_hoJ[jho];
1090 if (dof_used_ho[dof_ho] != ilor)
1091 {
1092 dof_used_ho[dof_ho] = ilor;
1093 dof_dofJ[sizeJ] = dof_ho;
1094 ++sizeJ;
1095 }
1096 }
1097 }
1098 }
1099
1100 dof_lor_dof_ho.SortRows();
1101 real_t* data = Memory<real_t>(dof_dofI[ndof_lor]);
1102
1103 std::unique_ptr<SparseMatrix> R_local(new SparseMatrix(
1104 dof_dofI, dof_dofJ, data, ndof_lor,
1105 ndof_ho, true, true, true));
1106 (*R_local) = 0.0;
1107
1108 dof_lor_dof_ho.LoseData();
1109
1110 return R_local;
1111}
1112
1114{
1115 delete F;
1116 delete B;
1117}
1118
1120{
1121 if (!F) { BuildF(); }
1122 return *F;
1123}
1124
1126{
1127 if (!B)
1128 {
1129 if (!F) { BuildF(); }
1130 B = new L2Prolongation(*F);
1131 }
1132 return *B;
1133}
1134
1135void L2ProjectionGridTransfer::BuildF()
1136{
1137 if (!force_l2_space &&
1139 {
1140 if (!Parallel())
1141 {
1142 F = new L2ProjectionH1Space(dom_fes, ran_fes);
1143 }
1144 else
1145 {
1146#ifdef MFEM_USE_MPI
1147 const mfem::ParFiniteElementSpace& dom_pfes =
1148 static_cast<mfem::ParFiniteElementSpace&>(dom_fes);
1149 const mfem::ParFiniteElementSpace& ran_pfes =
1150 static_cast<mfem::ParFiniteElementSpace&>(ran_fes);
1151 F = new L2ProjectionH1Space(dom_pfes, ran_pfes);
1152#endif
1153 }
1154 }
1155 else
1156 {
1157 F = new L2ProjectionL2Space(dom_fes, ran_fes);
1158 }
1159}
1160
1165
1166
1168 const FiniteElementSpace& hFESpace_)
1169 : Operator(hFESpace_.GetVSize(), lFESpace_.GetVSize())
1170{
1171 bool isvar_order = lFESpace_.IsVariableOrder() || hFESpace_.IsVariableOrder();
1172 if (lFESpace_.FEColl() == hFESpace_.FEColl() && !isvar_order)
1173 {
1175 hFESpace_.GetTransferOperator(lFESpace_, P);
1176 P.SetOperatorOwner(false);
1177 opr = P.Ptr();
1178 }
1179 else if (lFESpace_.GetMesh()->GetNE() > 0
1180 && hFESpace_.GetMesh()->GetNE() > 0
1181 && lFESpace_.GetVDim() == 1
1182 && hFESpace_.GetVDim() == 1
1183 && dynamic_cast<const TensorBasisElement*>(lFESpace_.GetFE(0))
1184 && dynamic_cast<const TensorBasisElement*>(hFESpace_.GetFE(0))
1185 && !isvar_order
1186 && (hFESpace_.FEColl()->GetContType() ==
1188 hFESpace_.FEColl()->GetContType() ==
1190 {
1191 opr = new TensorProductPRefinementTransferOperator(lFESpace_, hFESpace_);
1192 }
1193 else
1194 {
1195 opr = new PRefinementTransferOperator(lFESpace_, hFESpace_);
1196 }
1197}
1198
1200
1201void TransferOperator::Mult(const Vector& x, Vector& y) const
1202{
1203 opr->Mult(x, y);
1204}
1205
1207{
1208 opr->MultTranspose(x, y);
1209}
1210
1211
1213 const FiniteElementSpace& lFESpace_, const FiniteElementSpace& hFESpace_)
1214 : Operator(hFESpace_.GetVSize(), lFESpace_.GetVSize()), lFESpace(lFESpace_),
1215 hFESpace(hFESpace_)
1216{
1217 isvar_order = lFESpace_.IsVariableOrder() || hFESpace_.IsVariableOrder();
1218}
1219
1221
1223{
1224 Mesh* mesh = hFESpace.GetMesh();
1225 Array<int> l_dofs, h_dofs, l_vdofs, h_vdofs;
1226 DenseMatrix loc_prol;
1227 Vector subY, subX;
1228
1229 Geometry::Type cached_geom = Geometry::INVALID;
1230 const FiniteElement* h_fe = NULL;
1231 const FiniteElement* l_fe = NULL;
1233
1234 int vdim = lFESpace.GetVDim();
1235
1236 for (int i = 0; i < mesh->GetNE(); i++)
1237 {
1238 DofTransformation * doftrans_h = hFESpace.GetElementDofs(i, h_dofs);
1239 DofTransformation * doftrans_l = lFESpace.GetElementDofs(i, l_dofs);
1240
1241 const Geometry::Type geom = mesh->GetElementBaseGeometry(i);
1242 if (geom != cached_geom || isvar_order)
1243 {
1244 h_fe = hFESpace.GetFE(i);
1245 l_fe = lFESpace.GetFE(i);
1247 h_fe->GetTransferMatrix(*l_fe, T, loc_prol);
1248 subY.SetSize(loc_prol.Height());
1249 cached_geom = geom;
1250 }
1251
1252 for (int vd = 0; vd < vdim; vd++)
1253 {
1254 l_dofs.Copy(l_vdofs);
1255 lFESpace.DofsToVDofs(vd, l_vdofs);
1256 h_dofs.Copy(h_vdofs);
1257 hFESpace.DofsToVDofs(vd, h_vdofs);
1258 x.GetSubVector(l_vdofs, subX);
1259 if (doftrans_l)
1260 {
1261 doftrans_l->InvTransformPrimal(subX);
1262 }
1263 loc_prol.Mult(subX, subY);
1264 if (doftrans_h)
1265 {
1266 doftrans_h->TransformPrimal(subY);
1267 }
1268 y.SetSubVector(h_vdofs, subY);
1269 }
1270 }
1271}
1272
1274 Vector& y) const
1275{
1276 y = 0.0;
1277
1278 Mesh* mesh = hFESpace.GetMesh();
1279 Array<int> l_dofs, h_dofs, l_vdofs, h_vdofs;
1280 DenseMatrix loc_prol;
1281 Vector subY, subX;
1282
1283 Array<char> processed(hFESpace.GetVSize());
1284 processed = 0;
1285
1286 Geometry::Type cached_geom = Geometry::INVALID;
1287 const FiniteElement* h_fe = NULL;
1288 const FiniteElement* l_fe = NULL;
1290
1291 int vdim = lFESpace.GetVDim();
1292
1293 for (int i = 0; i < mesh->GetNE(); i++)
1294 {
1295 DofTransformation * doftrans_h = hFESpace.GetElementDofs(i, h_dofs);
1296 DofTransformation * doftrans_l = lFESpace.GetElementDofs(i, l_dofs);
1297
1298 const Geometry::Type geom = mesh->GetElementBaseGeometry(i);
1299 if (geom != cached_geom || isvar_order)
1300 {
1301 h_fe = hFESpace.GetFE(i);
1302 l_fe = lFESpace.GetFE(i);
1304 h_fe->GetTransferMatrix(*l_fe, T, loc_prol);
1305 loc_prol.Transpose();
1306 subY.SetSize(loc_prol.Height());
1307 cached_geom = geom;
1308 }
1309
1310 for (int vd = 0; vd < vdim; vd++)
1311 {
1312 l_dofs.Copy(l_vdofs);
1313 lFESpace.DofsToVDofs(vd, l_vdofs);
1314 h_dofs.Copy(h_vdofs);
1315 hFESpace.DofsToVDofs(vd, h_vdofs);
1316
1317 x.GetSubVector(h_vdofs, subX);
1318 if (doftrans_h)
1319 {
1320 doftrans_h->InvTransformDual(subX);
1321 }
1322 for (int p = 0; p < h_dofs.Size(); ++p)
1323 {
1324 if (processed[lFESpace.DecodeDof(h_dofs[p])])
1325 {
1326 subX[p] = 0.0;
1327 }
1328 }
1329
1330 loc_prol.Mult(subX, subY);
1331 if (doftrans_l)
1332 {
1333 doftrans_l->TransformDual(subY);
1334 }
1335 y.AddElementVector(l_vdofs, subY);
1336 }
1337
1338 for (int p = 0; p < h_dofs.Size(); ++p)
1339 {
1340 processed[lFESpace.DecodeDof(h_dofs[p])] = 1;
1341 }
1342 }
1343}
1344
1345
1348 const FiniteElementSpace& lFESpace_,
1349 const FiniteElementSpace& hFESpace_)
1350 : Operator(hFESpace_.GetVSize(), lFESpace_.GetVSize()), lFESpace(lFESpace_),
1351 hFESpace(hFESpace_)
1352{
1353 // Assuming the same element type
1354 Mesh* mesh = lFESpace.GetMesh();
1355 dim = mesh->Dimension();
1356 if (mesh->GetNE() == 0)
1357 {
1358 return;
1359 }
1360 const FiniteElement& el = *lFESpace.GetFE(0);
1361
1362 const TensorBasisElement* ltel =
1363 dynamic_cast<const TensorBasisElement*>(&el);
1364 MFEM_VERIFY(ltel, "Low order FE space must be tensor product space");
1365
1366 const TensorBasisElement* htel =
1367 dynamic_cast<const TensorBasisElement*>(hFESpace.GetFE(0));
1368 MFEM_VERIFY(htel, "High order FE space must be tensor product space");
1369 const Array<int>& hdofmap = htel->GetDofMap();
1370
1371 const IntegrationRule& ir = hFESpace.GetFE(0)->GetNodes();
1372 IntegrationRule irLex = ir;
1373
1374 // The quadrature points, or equivalently, the dofs of the high order space
1375 // must be sorted in lexicographical order
1376 for (int i = 0; i < ir.GetNPoints(); ++i)
1377 {
1378 int j = hdofmap[i] >=0 ? hdofmap[i] : -1 - hdofmap[i];
1379 irLex.IntPoint(i) = ir.IntPoint(j);
1380 }
1381
1382 NE = lFESpace.GetNE();
1383 const DofToQuad& maps = el.GetDofToQuad(irLex, DofToQuad::TENSOR);
1384
1385 D1D = maps.ndof;
1386 Q1D = maps.nqpt;
1387 B = maps.B;
1388 Bt = maps.Bt;
1389
1390 elem_restrict_lex_l =
1392
1393 MFEM_VERIFY(elem_restrict_lex_l,
1394 "Low order ElementRestriction not available");
1395
1396 elem_restrict_lex_h =
1398
1399 MFEM_VERIFY(elem_restrict_lex_h,
1400 "High order ElementRestriction not available");
1401
1402 localL.SetSize(elem_restrict_lex_l->Height(), Device::GetMemoryType());
1403 localH.SetSize(elem_restrict_lex_h->Height(), Device::GetMemoryType());
1404 localL.UseDevice(true);
1405 localH.UseDevice(true);
1406
1407 MFEM_VERIFY(dynamic_cast<const ElementRestriction*>(elem_restrict_lex_h),
1408 "High order element restriction is of unsupported type");
1409
1410 mask.SetSize(localH.Size(), Device::GetMemoryType());
1411 static_cast<const ElementRestriction*>(elem_restrict_lex_h)
1412 ->BooleanMask(mask);
1413 mask.UseDevice(true);
1414}
1415
1416namespace TransferKernels
1417{
1418void Prolongation2D(const int NE, const int D1D, const int Q1D,
1419 const Vector& localL, Vector& localH,
1420 const Array<real_t>& B, const Vector& mask)
1421{
1422 auto x_ = Reshape(localL.Read(), D1D, D1D, NE);
1423 auto y_ = Reshape(localH.Write(), Q1D, Q1D, NE);
1424 auto B_ = Reshape(B.Read(), Q1D, D1D);
1425 auto m_ = Reshape(mask.Read(), Q1D, Q1D, NE);
1426
1427 mfem::forall(NE, [=] MFEM_HOST_DEVICE (int e)
1428 {
1429 for (int qy = 0; qy < Q1D; ++qy)
1430 {
1431 for (int qx = 0; qx < Q1D; ++qx)
1432 {
1433 y_(qx, qy, e) = 0.0;
1434 }
1435 }
1436
1437 for (int dy = 0; dy < D1D; ++dy)
1438 {
1439 real_t sol_x[DofQuadLimits::MAX_Q1D];
1440 for (int qy = 0; qy < Q1D; ++qy)
1441 {
1442 sol_x[qy] = 0.0;
1443 }
1444 for (int dx = 0; dx < D1D; ++dx)
1445 {
1446 const real_t s = x_(dx, dy, e);
1447 for (int qx = 0; qx < Q1D; ++qx)
1448 {
1449 sol_x[qx] += B_(qx, dx) * s;
1450 }
1451 }
1452 for (int qy = 0; qy < Q1D; ++qy)
1453 {
1454 const real_t d2q = B_(qy, dy);
1455 for (int qx = 0; qx < Q1D; ++qx)
1456 {
1457 y_(qx, qy, e) += d2q * sol_x[qx];
1458 }
1459 }
1460 }
1461 for (int qy = 0; qy < Q1D; ++qy)
1462 {
1463 for (int qx = 0; qx < Q1D; ++qx)
1464 {
1465 y_(qx, qy, e) *= m_(qx, qy, e);
1466 }
1467 }
1468 });
1469}
1470
1471void Prolongation3D(const int NE, const int D1D, const int Q1D,
1472 const Vector& localL, Vector& localH,
1473 const Array<real_t>& B, const Vector& mask)
1474{
1475 auto x_ = Reshape(localL.Read(), D1D, D1D, D1D, NE);
1476 auto y_ = Reshape(localH.Write(), Q1D, Q1D, Q1D, NE);
1477 auto B_ = Reshape(B.Read(), Q1D, D1D);
1478 auto m_ = Reshape(mask.Read(), Q1D, Q1D, Q1D, NE);
1479
1480 mfem::forall(NE, [=] MFEM_HOST_DEVICE (int e)
1481 {
1482 for (int qz = 0; qz < Q1D; ++qz)
1483 {
1484 for (int qy = 0; qy < Q1D; ++qy)
1485 {
1486 for (int qx = 0; qx < Q1D; ++qx)
1487 {
1488 y_(qx, qy, qz, e) = 0.0;
1489 }
1490 }
1491 }
1492
1493 for (int dz = 0; dz < D1D; ++dz)
1494 {
1495 real_t sol_xy[DofQuadLimits::MAX_Q1D][DofQuadLimits::MAX_Q1D];
1496 for (int qy = 0; qy < Q1D; ++qy)
1497 {
1498 for (int qx = 0; qx < Q1D; ++qx)
1499 {
1500 sol_xy[qy][qx] = 0.0;
1501 }
1502 }
1503 for (int dy = 0; dy < D1D; ++dy)
1504 {
1505 real_t sol_x[DofQuadLimits::MAX_Q1D];
1506 for (int qx = 0; qx < Q1D; ++qx)
1507 {
1508 sol_x[qx] = 0;
1509 }
1510 for (int dx = 0; dx < D1D; ++dx)
1511 {
1512 const real_t s = x_(dx, dy, dz, e);
1513 for (int qx = 0; qx < Q1D; ++qx)
1514 {
1515 sol_x[qx] += B_(qx, dx) * s;
1516 }
1517 }
1518 for (int qy = 0; qy < Q1D; ++qy)
1519 {
1520 const real_t wy = B_(qy, dy);
1521 for (int qx = 0; qx < Q1D; ++qx)
1522 {
1523 sol_xy[qy][qx] += wy * sol_x[qx];
1524 }
1525 }
1526 }
1527 for (int qz = 0; qz < Q1D; ++qz)
1528 {
1529 const real_t wz = B_(qz, dz);
1530 for (int qy = 0; qy < Q1D; ++qy)
1531 {
1532 for (int qx = 0; qx < Q1D; ++qx)
1533 {
1534 y_(qx, qy, qz, e) += wz * sol_xy[qy][qx];
1535 }
1536 }
1537 }
1538 }
1539 for (int qz = 0; qz < Q1D; ++qz)
1540 {
1541 for (int qy = 0; qy < Q1D; ++qy)
1542 {
1543 for (int qx = 0; qx < Q1D; ++qx)
1544 {
1545 y_(qx, qy, qz, e) *= m_(qx, qy, qz, e);
1546 }
1547 }
1548 }
1549 });
1550}
1551
1552void Restriction2D(const int NE, const int D1D, const int Q1D,
1553 const Vector& localH, Vector& localL,
1554 const Array<real_t>& Bt, const Vector& mask)
1555{
1556 auto x_ = Reshape(localH.Read(), Q1D, Q1D, NE);
1557 auto y_ = Reshape(localL.Write(), D1D, D1D, NE);
1558 auto Bt_ = Reshape(Bt.Read(), D1D, Q1D);
1559 auto m_ = Reshape(mask.Read(), Q1D, Q1D, NE);
1560
1561 mfem::forall(NE, [=] MFEM_HOST_DEVICE (int e)
1562 {
1563 for (int dy = 0; dy < D1D; ++dy)
1564 {
1565 for (int dx = 0; dx < D1D; ++dx)
1566 {
1567 y_(dx, dy, e) = 0.0;
1568 }
1569 }
1570
1571 for (int qy = 0; qy < Q1D; ++qy)
1572 {
1573 real_t sol_x[DofQuadLimits::MAX_D1D];
1574 for (int dx = 0; dx < D1D; ++dx)
1575 {
1576 sol_x[dx] = 0.0;
1577 }
1578 for (int qx = 0; qx < Q1D; ++qx)
1579 {
1580 const real_t s = m_(qx, qy, e) * x_(qx, qy, e);
1581 for (int dx = 0; dx < D1D; ++dx)
1582 {
1583 sol_x[dx] += Bt_(dx, qx) * s;
1584 }
1585 }
1586 for (int dy = 0; dy < D1D; ++dy)
1587 {
1588 const real_t q2d = Bt_(dy, qy);
1589 for (int dx = 0; dx < D1D; ++dx)
1590 {
1591 y_(dx, dy, e) += q2d * sol_x[dx];
1592 }
1593 }
1594 }
1595 });
1596}
1597void Restriction3D(const int NE, const int D1D, const int Q1D,
1598 const Vector& localH, Vector& localL,
1599 const Array<real_t>& Bt, const Vector& mask)
1600{
1601 auto x_ = Reshape(localH.Read(), Q1D, Q1D, Q1D, NE);
1602 auto y_ = Reshape(localL.Write(), D1D, D1D, D1D, NE);
1603 auto Bt_ = Reshape(Bt.Read(), D1D, Q1D);
1604 auto m_ = Reshape(mask.Read(), Q1D, Q1D, Q1D, NE);
1605
1606 mfem::forall(NE, [=] MFEM_HOST_DEVICE (int e)
1607 {
1608 for (int dz = 0; dz < D1D; ++dz)
1609 {
1610 for (int dy = 0; dy < D1D; ++dy)
1611 {
1612 for (int dx = 0; dx < D1D; ++dx)
1613 {
1614 y_(dx, dy, dz, e) = 0.0;
1615 }
1616 }
1617 }
1618
1619 for (int qz = 0; qz < Q1D; ++qz)
1620 {
1621 real_t sol_xy[DofQuadLimits::MAX_D1D][DofQuadLimits::MAX_D1D];
1622 for (int dy = 0; dy < D1D; ++dy)
1623 {
1624 for (int dx = 0; dx < D1D; ++dx)
1625 {
1626 sol_xy[dy][dx] = 0;
1627 }
1628 }
1629 for (int qy = 0; qy < Q1D; ++qy)
1630 {
1631 real_t sol_x[DofQuadLimits::MAX_D1D];
1632 for (int dx = 0; dx < D1D; ++dx)
1633 {
1634 sol_x[dx] = 0;
1635 }
1636 for (int qx = 0; qx < Q1D; ++qx)
1637 {
1638 const real_t s = m_(qx, qy, qz, e) * x_(qx, qy, qz, e);
1639 for (int dx = 0; dx < D1D; ++dx)
1640 {
1641 sol_x[dx] += Bt_(dx, qx) * s;
1642 }
1643 }
1644 for (int dy = 0; dy < D1D; ++dy)
1645 {
1646 const real_t wy = Bt_(dy, qy);
1647 for (int dx = 0; dx < D1D; ++dx)
1648 {
1649 sol_xy[dy][dx] += wy * sol_x[dx];
1650 }
1651 }
1652 }
1653 for (int dz = 0; dz < D1D; ++dz)
1654 {
1655 const real_t wz = Bt_(dz, qz);
1656 for (int dy = 0; dy < D1D; ++dy)
1657 {
1658 for (int dx = 0; dx < D1D; ++dx)
1659 {
1660 y_(dx, dy, dz, e) += wz * sol_xy[dy][dx];
1661 }
1662 }
1663 }
1664 }
1665 });
1666}
1667} // namespace TransferKernels
1668
1669
1674
1676 Vector& y) const
1677{
1678 if (lFESpace.GetMesh()->GetNE() == 0)
1679 {
1680 return;
1681 }
1682
1683 elem_restrict_lex_l->Mult(x, localL);
1684 if (dim == 2)
1685 {
1686 TransferKernels::Prolongation2D(NE, D1D, Q1D, localL, localH, B, mask);
1687 }
1688 else if (dim == 3)
1689 {
1690 TransferKernels::Prolongation3D(NE, D1D, Q1D, localL, localH, B, mask);
1691 }
1692 else
1693 {
1694 MFEM_ABORT("TensorProductPRefinementTransferOperator::Mult not "
1695 "implemented for dim = "
1696 << dim);
1697 }
1698 elem_restrict_lex_h->MultTranspose(localH, y);
1699}
1700
1702 Vector& y) const
1703{
1704 if (lFESpace.GetMesh()->GetNE() == 0)
1705 {
1706 return;
1707 }
1708
1709 elem_restrict_lex_h->Mult(x, localH);
1710 if (dim == 2)
1711 {
1712 TransferKernels::Restriction2D(NE, D1D, Q1D, localH, localL, Bt, mask);
1713 }
1714 else if (dim == 3)
1715 {
1716 TransferKernels::Restriction3D(NE, D1D, Q1D, localH, localL, Bt, mask);
1717 }
1718 else
1719 {
1720 MFEM_ABORT("TensorProductPRefinementTransferOperator::MultTranspose not "
1721 "implemented for dim = "
1722 << dim);
1723 }
1724 elem_restrict_lex_l->MultTranspose(localL, y);
1725}
1726
1727
1729 const FiniteElementSpace& hFESpace_)
1730 : Operator(hFESpace_.GetTrueVSize(), lFESpace_.GetTrueVSize()),
1731 lFESpace(lFESpace_),
1732 hFESpace(hFESpace_)
1733{
1734 localTransferOperator = new TransferOperator(lFESpace_, hFESpace_);
1735
1736 P = lFESpace.GetProlongationMatrix();
1737 R = hFESpace.IsVariableOrder() ? hFESpace.GetHpRestrictionMatrix() :
1738 hFESpace.GetRestrictionMatrix();
1739
1740 // P and R can be both null
1741 // P can be null and R not null
1742 // If P is not null it is assumed that R is not null as well
1743 if (P) { MFEM_VERIFY(R, "Both P and R have to be not NULL") }
1744
1745 if (P)
1746 {
1747 tmpL.SetSize(lFESpace_.GetVSize());
1748 tmpH.SetSize(hFESpace_.GetVSize());
1749 }
1750 // P can be null and R not null
1751 else if (R)
1752 {
1753 tmpH.SetSize(hFESpace_.GetVSize());
1754 }
1755}
1756
1758{
1759 delete localTransferOperator;
1760}
1761
1763{
1764 if (P)
1765 {
1766 P->Mult(x, tmpL);
1767 localTransferOperator->Mult(tmpL, tmpH);
1768 R->Mult(tmpH, y);
1769 }
1770 else if (R)
1771 {
1772 localTransferOperator->Mult(x, tmpH);
1773 R->Mult(tmpH, y);
1774 }
1775 else
1776 {
1777 localTransferOperator->Mult(x, y);
1778 }
1779}
1780
1782{
1783 if (P)
1784 {
1785 R->MultTranspose(x, tmpH);
1786 localTransferOperator->MultTranspose(tmpH, tmpL);
1787 P->MultTranspose(tmpL, y);
1788 }
1789 else if (R)
1790 {
1791 R->MultTranspose(x, tmpH);
1792 localTransferOperator->MultTranspose(tmpH, y);
1793 }
1794 else
1795 {
1796 localTransferOperator->MultTranspose(x, y);
1797 }
1798}
1799
1800} // namespace mfem
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 Copy(Array &copy) const
Create a copy of the internal array to the provided copy.
Definition array.hpp:874
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Definition array.hpp:317
Abstract base class BilinearFormIntegrator.
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition solvers.hpp:526
Data type for scaled Jacobi-type smoother of sparse matrix.
void Factor()
Factor the current DenseMatrix, *a.
void Mult(const real_t *x, real_t *y) const
Matrix vector multiplication with the inverse of dense matrix.
void GetInverseMatrix(DenseMatrix &Ainv) const
Compute and return the inverse matrix in Ainv.
Data type dense matrix using column-major storage.
Definition densemat.hpp:24
void Mult(const real_t *x, real_t *y) const
Matrix vector multiplication.
Definition densemat.cpp:220
void Transpose()
(*this) = (*this)^t
void SetRow(int r, const real_t *row)
real_t * GetData() const
Returns the matrix data array.
Definition densemat.hpp:115
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition densemat.hpp:105
void CopyMN(const DenseMatrix &A, int m, int n, int Aro, int Aco)
Copy the m x n submatrix of A at row/col offsets Aro/Aco to *this.
void GetRow(int r, Vector &row) const
Rank 3 tensor (array of matrices)
int SizeK() const
static MemoryType GetMemoryType()
(DEPRECATED) Equivalent to GetDeviceMemoryType().
Definition device.hpp:278
Structure representing the matrices/tensors needed to evaluate (in reference space) the values,...
Definition fe_base.hpp:137
@ TENSOR
Tensor product representation using 1D matrices/tensors with dimensions using 1D number of quadrature...
Definition fe_base.hpp:161
Array< real_t > B
Basis functions evaluated at quadrature points.
Definition fe_base.hpp:189
int ndof
Number of degrees of freedom = number of basis functions. When mode is TENSOR, this is the 1D number.
Definition fe_base.hpp:174
int nqpt
Number of quadrature points. When mode is TENSOR, this is the 1D number.
Definition fe_base.hpp:178
Array< real_t > Bt
Transpose of B.
Definition fe_base.hpp:195
void TransformDual(real_t *v) const
Definition doftrans.cpp:81
void InvTransformPrimal(real_t *v) const
Definition doftrans.cpp:49
void InvTransformDual(real_t *v) const
Definition doftrans.cpp:113
void TransformPrimal(real_t *v) const
Definition doftrans.cpp:17
Operator that converts FiniteElementSpace L-vectors to E-vectors.
real_t Weight()
Return the weight of the Jacobian matrix of the transformation at the currently set IntegrationPoint....
Definition eltrans.hpp:131
virtual int OrderW() const =0
Return the order of the determinant of the Jacobian (weight) of the transformation.
void SetIntPoint(const IntegrationPoint *ip)
Set the integration point ip that weights and Jacobians will be evaluated at.
Definition eltrans.hpp:93
@ DISCONTINUOUS
Field is discontinuous across element interfaces.
Definition fe_coll.hpp:48
@ CONTINUOUS
Field is continuous across element interfaces.
Definition fe_coll.hpp:45
virtual int GetContType() const =0
Derefinement operator, used by the friend class InterpolationGridTransfer.
Definition fespace.hpp:443
GridFunction interpolation operator applicable after mesh refinement.
Definition fespace.hpp:417
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:220
void GetVDofs(int vd, Array< int > &dofs, int ndofs=-1) const
Returns the indices of all of the VDofs for the specified dimension 'vd'.
Definition fespace.cpp:193
const SparseMatrix * GetConformingRestriction() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition fespace.cpp:1309
bool IsVariableOrder() const
Returns true if the space contains elements of varying polynomial orders.
Definition fespace.hpp:581
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
void DofsToVDofs(Array< int > &dofs, int ndofs=-1) const
Compute the full set of vdofs corresponding to each entry in dofs.
Definition fespace.cpp:213
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition fespace.hpp:716
DofTransformation * GetElementDofs(int elem, Array< int > &dofs) const
Returns indices of degrees of freedom of element 'elem'. The returned indices are offsets into an ldo...
Definition fespace.cpp:2846
ElementTransformation * GetElementTransformation(int i) const
Returns ElementTransformation for the i-th element.
Definition fespace.hpp:773
virtual const SparseMatrix * GetRestrictionMatrix() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition fespace.hpp:620
static void ListToMarker(const Array< int > &list, int marker_size, Array< int > &marker, int mark_val=-1)
Convert an array of indices (list) to a Boolean marker array where all indices in the list are marked...
Definition fespace.cpp:667
int GetNDofs() const
Returns number of degrees of freedom. This is the number of Local Degrees of Freedom.
Definition fespace.hpp:710
void GetLocalRefinementMatrices(Geometry::Type geom, DenseTensor &localP) const
Definition fespace.cpp:1528
void GetTransferOperator(const FiniteElementSpace &coarse_fes, OperatorHandle &T) const
Construct and return an Operator that can be used to transfer GridFunction data from coarse_fes,...
Definition fespace.cpp:3349
virtual const Operator * GetProlongationMatrix() const
The returned Operator is owned by the FiniteElementSpace.
Definition fespace.hpp:597
virtual const Operator * GetRestrictionOperator() const
An abstract operator that performs the same action as GetRestrictionMatrix.
Definition fespace.hpp:616
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
int GetNE() const
Returns number of elements in the mesh.
Definition fespace.hpp:740
const ElementRestrictionOperator * GetElementRestriction(ElementDofOrdering e_ordering) const
Return an Operator that converts L-vectors to E-vectors.
Definition fespace.cpp:1336
SparseMatrix * RefinementMatrix_main(const int coarse_ndofs, const Table &coarse_elem_dof, const Table *coarse_elem_fos, const DenseTensor localP[]) const
Definition fespace.cpp:1464
static void MarkerToList(const Array< int > &marker, Array< int > &list)
Convert a Boolean marker array to a list containing all marked indices.
Definition fespace.cpp:648
const SparseMatrix * GetConformingProlongation() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition fespace.cpp:1302
const FiniteElementCollection * FEColl() const
Definition fespace.hpp:727
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
const Table * GetElementToFaceOrientationTable() const
Definition fespace.hpp:1132
int GetVDim() const
Returns vector dimension.
Definition fespace.hpp:706
virtual const SparseMatrix * GetHpRestrictionMatrix() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition fespace.hpp:624
static int DecodeDof(int dof)
Helper to return the DOF associated with a sign encoded DOF.
Definition fespace.hpp:1017
Abstract class for all finite elements.
Definition fe_base.hpp:239
virtual const DofToQuad & GetDofToQuad(const IntegrationRule &ir, DofToQuad::Mode mode) const
Return a DofToQuad structure corresponding to the given IntegrationRule using the given DofToQuad::Mo...
Definition fe_base.cpp:365
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 GetTransferMatrix(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &I) const
Return interpolation matrix, I, which maps dofs from a coarse element, fe, to the fine dofs on this f...
Definition fe_base.cpp:119
int GetMapType() const
Returns the FiniteElement::MapType of the element describing how reference functions are mapped to ph...
Definition fe_base.hpp:355
const IntegrationRule & GetNodes() const
Get a const reference to the nodes of the element.
Definition fe_base.hpp:395
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
Definition fe_base.hpp:326
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const =0
Evaluate the values of all shape functions of a scalar finite element in reference space at the given...
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition fe_base.hpp:329
static const int NumGeom
Definition geom.hpp:42
bool Parallel() const
Definition transfer.hpp:46
FiniteElementSpace & dom_fes
Domain FE space.
Definition transfer.hpp:32
FiniteElementSpace & ran_fes
Range FE space.
Definition transfer.hpp:33
Operator::Type oper_type
Desired Operator::Type for the construction of all operators defined by the underlying transfer algor...
Definition transfer.hpp:38
const Operator & MakeTrueOperator(FiniteElementSpace &fes_in, FiniteElementSpace &fes_out, const Operator &oper, OperatorHandle &t_oper)
Definition transfer.cpp:34
GridTransfer(FiniteElementSpace &dom_fes_, FiniteElementSpace &ran_fes_)
Definition transfer.cpp:19
The BoomerAMG solver in hypre.
Definition hypre.hpp:1691
void SetPrintLevel(int print_level)
Definition hypre.hpp:1771
Wrapper for hypre's ParCSR matrix class.
Definition hypre.hpp:388
HypreParMatrix * LeftDiagMult(const SparseMatrix &D, HYPRE_BigInt *row_starts=NULL) const
Multiply the HypreParMatrix on the left by a block-diagonal parallel matrix D and return the result a...
Definition hypre.cpp:1994
HypreParMatrix * Transpose() const
Returns the transpose of *this.
Definition hypre.cpp:1684
IsoparametricTransformation Transf
Definition eltrans.hpp:467
void Transform(const IntegrationPoint &, IntegrationPoint &)
Definition eltrans.cpp:543
Class for integration point with weight.
Definition intrules.hpp:35
Class for an integration rule - an Array of IntegrationPoint.
Definition intrules.hpp:100
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition intrules.hpp:256
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition intrules.hpp:259
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
BilinearFormIntegrator * mass_integ
Ownership depends on own_mass_integ.
Definition transfer.hpp:126
OperatorHandle F
Forward, coarse-to-fine, operator.
Definition transfer.hpp:129
bool own_mass_integ
Ownership flag for mass_integ.
Definition transfer.hpp:127
virtual const Operator & BackwardOperator()
Return an Operator that transfers GridFunctions from the range FE space back to GridFunctions in the ...
Definition transfer.cpp:189
OperatorHandle B
Backward, fine-to-coarse, operator.
Definition transfer.hpp:130
void SetMassIntegrator(BilinearFormIntegrator *mass_integ_, bool own_mass_integ_=true)
Assign a mass integrator to be used in the construction of the backward, fine-to-coarse,...
Definition transfer.cpp:146
virtual const Operator & ForwardOperator()
Return an Operator that transfers GridFunctions from the domain FE space to GridFunctions in the rang...
Definition transfer.cpp:155
A standard isoparametric element transformation.
Definition eltrans.hpp:363
void SetPointMat(const DenseMatrix &pm)
Set the underlying point matrix describing the transformation.
Definition eltrans.hpp:402
void SetIdentityTransformation(Geometry::Type GeomType)
Set the FiniteElement Geometry for the reference elements being used.
Definition eltrans.cpp:379
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition solvers.cpp:173
void SetFromTDofsTranspose(const FiniteElementSpace &fes, const Vector &X, Vector &x) const
Sets dual field coefficients given a vector of dual field coefficients on the tdofs and a finite elem...
Definition transfer.cpp:962
virtual void Prolongate(const Vector &x, Vector &y) const
Definition transfer.cpp:701
void TDofsListByVDim(const FiniteElementSpace &fes, int vdim, Array< int > &vdofs_list) const
Fills the vdofs_list array with a list of vdofs for a given vdim and a given finite element space.
Definition transfer.cpp:976
void SetupPCG()
Sets up the PCG solver (sets parameters, operator, and preconditioner)
Definition transfer.cpp:638
virtual void ProlongateTranspose(const Vector &x, Vector &y) const
Definition transfer.cpp:730
void GetTDofs(const FiniteElementSpace &fes, const Vector &x, Vector &X) const
Recovers vector of tdofs given a vector of dofs and a finite element space.
Definition transfer.cpp:920
void SetFromTDofs(const FiniteElementSpace &fes, const Vector &X, Vector &x) const
Sets dof values given a vector of tdofs and a finite element space.
Definition transfer.cpp:934
virtual void SetAbsTol(real_t p_atol_)
Sets absolute tolerance in preconditioned conjugate gradient solver.
Definition transfer.cpp:764
L2ProjectionH1Space(const FiniteElementSpace &fes_ho_, const FiniteElementSpace &fes_lor_)
Definition transfer.cpp:545
void GetTDofsTranspose(const FiniteElementSpace &fes, const Vector &x, Vector &X) const
Recovers a vector of dual field coefficients on the tdofs given a vector of dual coefficients and a f...
Definition transfer.cpp:948
virtual void Mult(const Vector &x, Vector &y) const
Definition transfer.cpp:651
std::unique_ptr< SparseMatrix > AllocR()
Computes sparsity pattern and initializes R matrix.
void LumpedMassInverse(Vector &ML_inv) const
Returns the inverse of an on-rank lumped mass matrix.
Definition transfer.cpp:997
std::pair< std::unique_ptr< SparseMatrix >, std::unique_ptr< SparseMatrix > > ComputeSparseRAndM_LH()
Computes on-rank R and M_LH matrices.
Definition transfer.cpp:772
virtual void MultTranspose(const Vector &x, Vector &y) const
Definition transfer.cpp:676
virtual void SetRelTol(real_t p_rtol_)
Sets relative tolerance in preconditioned conjugate gradient solver.
Definition transfer.cpp:759
L2ProjectionL2Space(const FiniteElementSpace &fes_ho_, const FiniteElementSpace &fes_lor_)
Definition transfer.cpp:285
virtual void Prolongate(const Vector &x, Vector &y) const
Definition transfer.cpp:473
virtual void MultTranspose(const Vector &x, Vector &y) const
Definition transfer.cpp:438
virtual void Mult(const Vector &x, Vector &y) const
Definition transfer.cpp:406
virtual void ProlongateTranspose(const Vector &x, Vector &y) const
Definition transfer.cpp:510
void BuildHo2Lor(int nel_ho, int nel_lor, const CoarseFineTransformations &cf_tr)
Definition transfer.cpp:238
void ElemMixedMass(Geometry::Type geom, const FiniteElement &fe_ho, const FiniteElement &fe_lor, ElementTransformation *el_tr, IntegrationPointTransformation &ip_tr, DenseMatrix &M_mixed_el) const
Definition transfer.cpp:258
L2Projection(const FiniteElementSpace &fes_ho_, const FiniteElementSpace &fes_lor_)
Definition transfer.cpp:231
L2Projection * F
Forward, coarse-to-fine, operator.
Definition transfer.hpp:374
virtual bool SupportsBackwardsOperator() const
virtual const Operator & ForwardOperator()
Return an Operator that transfers GridFunctions from the domain FE space to GridFunctions in the rang...
L2Prolongation * B
Backward, fine-to-coarse, operator.
Definition transfer.hpp:375
virtual const Operator & BackwardOperator()
Return an Operator that transfers GridFunctions from the range FE space back to GridFunctions in the ...
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Class used by MFEM to store pointers to host and/or device memory.
List of mesh geometries stored as Array<Geometry::Type>.
Definition mesh.hpp:1419
Mesh data type.
Definition mesh.hpp:56
void GetGeometries(int dim, Array< Geometry::Type > &el_geoms) const
Return all element geometries of the given dimension present in the mesh.
Definition mesh.cpp:6972
const CoarseFineTransformations & GetRefinementTransforms() const
Definition mesh.cpp:11082
int GetNE() const
Returns number of elements.
Definition mesh.hpp:1226
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1160
Geometry::Type GetElementBaseGeometry(int i) const
Definition mesh.hpp:1385
Pointer to an Operator of a specified type.
Definition handle.hpp:34
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get().
Definition handle.hpp:104
void SetOperatorOwner(bool own=true)
Set the ownership flag for the held Operator.
Definition handle.hpp:120
Operator * Ptr() const
Access the underlying Operator pointer.
Definition handle.hpp:87
void Reset(OpType *A, bool own_A=true)
Reset the OperatorHandle to the given OpType pointer, A.
Definition handle.hpp:145
Abstract operator.
Definition operator.hpp:25
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition operator.hpp:66
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
@ ANY_TYPE
ID for the base class Operator, i.e. any type.
Definition operator.hpp:285
@ MFEM_SPARSEMAT
ID for class SparseMatrix.
Definition operator.hpp:286
@ Hypre_ParCSR
ID for class HypreParMatrix.
Definition operator.hpp:287
virtual void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
Definition operator.hpp:93
Matrix-free transfer operator between finite element spaces on the same mesh.
Definition transfer.hpp:429
PRefinementTransferOperator(const FiniteElementSpace &lFESpace_, const FiniteElementSpace &hFESpace_)
Constructs a transfer operator from lFESpace to hFESpace which have different FE collections.
virtual ~PRefinementTransferOperator()
Destructor.
virtual void MultTranspose(const Vector &x, Vector &y) const override
Restriction by applying the transpose of the Mult method.
virtual void Mult(const Vector &x, Vector &y) const override
Interpolation or prolongation of a vector x corresponding to the coarse space to the vector y corresp...
Abstract parallel finite element space.
Definition pfespace.hpp:29
MPI_Comm GetComm() const
Definition pfespace.hpp:273
HYPRE_BigInt * GetTrueDofOffsets() const
Definition pfespace.hpp:282
HYPRE_BigInt GlobalVSize() const
Definition pfespace.hpp:283
HYPRE_BigInt * GetDofOffsets() const
Definition pfespace.hpp:281
HypreParMatrix * Dof_TrueDof_Matrix() const
The true dof-to-dof interpolation matrix.
Definition pfespace.hpp:327
ParMesh * GetParMesh() const
Definition pfespace.hpp:277
General product operator: x -> (A*B)(x) = A(B(x)).
Definition operator.hpp:797
Data type sparse matrix.
Definition sparsemat.hpp:51
virtual void MultTranspose(const Vector &x, Vector &y) const
Multiply a vector with the transposed matrix. y = At * x.
void BooleanMult(const Array< int > &x, Array< int > &y) const
y = A * x, treating all entries as booleans (zero=false, nonzero=true).
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
void LoseData()
Call this if data has been stolen.
Definition table.hpp:180
int * GetJ()
Definition table.hpp:114
int RowSize(int i) const
Definition table.hpp:108
void GetRow(int i, Array< int > &row) const
Return row i in array row (the Table must be finalized)
Definition table.cpp:187
int * GetI()
Definition table.hpp:113
void SortRows()
Sort the column (TYPE II) indices in each row.
Definition table.cpp:199
void SetDims(int rows, int nnz)
Definition table.cpp:140
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
Matrix-free transfer operator between finite element spaces on the same mesh exploiting the tensor pr...
Definition transfer.hpp:460
virtual void Mult(const Vector &x, Vector &y) const override
Interpolation or prolongation of a vector x corresponding to the coarse space to the vector y corresp...
TensorProductPRefinementTransferOperator(const FiniteElementSpace &lFESpace_, const FiniteElementSpace &hFESpace_)
Constructs a transfer operator from lFESpace to hFESpace which have different FE collections.
virtual ~TensorProductPRefinementTransferOperator()
Destructor.
virtual void MultTranspose(const Vector &x, Vector &y) const override
Restriction by applying the transpose of the Mult method.
Matrix-free transfer operator between finite element spaces.
Definition transfer.hpp:398
virtual ~TransferOperator()
Destructor.
virtual void Mult(const Vector &x, Vector &y) const override
Interpolation or prolongation of a vector x corresponding to the coarse space to the vector y corresp...
virtual void MultTranspose(const Vector &x, Vector &y) const override
Restriction by applying the transpose of the Mult method.
TransferOperator(const FiniteElementSpace &lFESpace, const FiniteElementSpace &hFESpace)
Constructs a transfer operator from lFESpace to hFESpace.
General triple product operator x -> A*B*C*x, with ownership of the factors.
Definition operator.hpp:861
~TrueTransferOperator()
Destructor.
virtual void MultTranspose(const Vector &x, Vector &y) const override
Restriction by applying the transpose of the Mult method.
TrueTransferOperator(const FiniteElementSpace &lFESpace_, const FiniteElementSpace &hFESpace_)
Constructs a transfer operator working on true degrees of freedom from lFESpace to hFESpace.
virtual void Mult(const Vector &x, Vector &y) const override
Interpolation or prolongation of a true dof vector x to a true dof vector y.
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
void SetSubVector(const Array< int > &dofs, const real_t value)
Set the entries listed in dofs to the given value.
Definition vector.cpp:604
void AddElementVector(const Array< int > &dofs, const Vector &elemvect)
Add elements of the elemvect Vector to the entries listed in dofs. Negative dof values cause the -dof...
Definition vector.cpp:670
Vector & Set(const real_t a, const Vector &x)
(*this) = a * x
Definition vector.cpp:262
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
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Extract entries listed in dofs to the output Vector elemvect.
Definition vector.cpp:578
virtual real_t * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
Definition vector.hpp:482
void Prolongation2D(const int NE, const int D1D, const int Q1D, const Vector &localL, Vector &localH, const Array< real_t > &B, const Vector &mask)
void Restriction3D(const int NE, const int D1D, const int Q1D, const Vector &localH, Vector &localL, const Array< real_t > &Bt, const Vector &mask)
void Restriction2D(const int NE, const int D1D, const int Q1D, const Vector &localH, Vector &localL, const Array< real_t > &Bt, const Vector &mask)
void Prolongation3D(const int NE, const int D1D, const int Q1D, const Vector &localL, Vector &localH, const Array< real_t > &B, const Vector &mask)
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition table.cpp:476
void Transpose(const Table &A, Table &At, int ncols_A_)
Transpose a Table.
Definition table.cpp:414
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
void RAP(const DenseMatrix &A, const DenseMatrix &P, DenseMatrix &RAP)
void AddMultVWt(const Vector &v, const Vector &w, DenseMatrix &VWt)
VWt += v w^t.
HypreParMatrix * ParMult(const HypreParMatrix *A, const HypreParMatrix *B, bool own_matrix)
Definition hypre.cpp:2949
ComplexDenseMatrix * MultAtB(const ComplexDenseMatrix &A, const ComplexDenseMatrix &B)
Multiply the complex conjugate transpose of a matrix A with a matrix B. A^H*B.
float real_t
Definition config.hpp:43
@ LEXICOGRAPHIC
Lexicographic ordering for tensor-product FiniteElements.
SparseMatrix * TransposeMult(const SparseMatrix &A, const SparseMatrix &B)
C = A^T B.
void forall(int N, lambda &&body)
Definition forall.hpp:754
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition intrules.hpp:486
real_t p(const Vector &x, real_t t)
RefCoord s[3]
Defines the coarse-fine transformations of all fine elements.
Definition ncmesh.hpp:72
Array< Embedding > embeddings
Fine element positions in their parents.
Definition ncmesh.hpp:74
DenseTensor point_matrices[Geometry::NumGeom]
Definition ncmesh.hpp:78