MFEM v4.8.0
Finite element discretization library
Loading...
Searching...
No Matches
pgridfunc.cpp
Go to the documentation of this file.
1// Copyright (c) 2010-2025, 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 "../config/config.hpp"
13
14#ifdef MFEM_USE_MPI
15
16#include "fem.hpp"
17#include <iostream>
18#include <limits>
19#include "../general/forall.hpp"
20using namespace std;
21
22namespace mfem
23{
24
30
36
38 const int *partitioning)
39{
40 const FiniteElementSpace *glob_fes = gf->FESpace();
41 // duplicate the FiniteElementCollection from 'gf'
43 // create a local ParFiniteElementSpace from the global one:
44 fes = pfes = new ParFiniteElementSpace(pmesh, glob_fes, partitioning,
45 fec_owned);
47
48 if (partitioning)
49 {
50 // Assumption: the map "local element id" -> "global element id" is
51 // increasing, i.e. the local numbering preserves the element order from
52 // the global numbering.
53 Array<int> gvdofs, lvdofs;
54 Vector lnodes;
55 int element_counter = 0;
56 const int MyRank = pfes->GetMyRank();
57 const int glob_ne = glob_fes->GetNE();
58 for (int i = 0; i < glob_ne; i++)
59 {
60 if (partitioning[i] == MyRank)
61 {
62 const DofTransformation* const ltrans = pfes->GetElementVDofs(element_counter,
63 lvdofs);
64 const DofTransformation* const gtrans = glob_fes->GetElementVDofs(i, gvdofs);
65 gf->GetSubVector(gvdofs, lnodes);
66 if (gtrans)
67 {
68 gtrans->InvTransformPrimal(lnodes);
69 }
70 if (ltrans)
71 {
72 ltrans->TransformPrimal(lnodes);
73 }
74 SetSubVector(lvdofs, lnodes);
75 element_counter++;
76 }
77 }
78 }
79}
80
81ParGridFunction::ParGridFunction(ParMesh *pmesh, std::istream &input)
82 : GridFunction(pmesh, input)
83{
84 // Convert the FiniteElementSpace, fes, to a ParFiniteElementSpace:
86 fes->GetOrdering());
87 delete fes;
88 fes = pfes;
89}
90
96
98{
101 pfes = dynamic_cast<ParFiniteElementSpace*>(f);
102 MFEM_ASSERT(pfes != NULL, "not a ParFiniteElementSpace");
103}
104
111
113{
116 pfes = dynamic_cast<ParFiniteElementSpace*>(f);
117 MFEM_ASSERT(pfes != NULL, "not a ParFiniteElementSpace");
118}
119
126
128{
130 GridFunction::MakeRef(f, v, v_offset);
131 pfes = dynamic_cast<ParFiniteElementSpace*>(f);
132 MFEM_ASSERT(pfes != NULL, "not a ParFiniteElementSpace");
133}
134
136{
138 GridFunction::MakeRef(f, v, v_offset);
139 pfes = f;
140}
141
143{
144 const Operator *prolong = pfes->GetProlongationMatrix();
145 prolong->Mult(*tv, *this);
146}
147
149{
150 pfes->Dof_TrueDof_Matrix()->Mult(a, *tv, 1.0, *this);
151}
152
154{
156 GetTrueDofs(*tv);
157 return tv;
158}
159
161{
162 MFEM_VERIFY(pfes->Conforming(), "not implemented for NC meshes");
165}
166
168{
169 MFEM_VERIFY(pfes->Conforming(), "not implemented for NC meshes");
172}
173
180
182{
183 pfes->GetRestrictionMatrix()->Mult(*this, tv);
184}
185
190
197
202
207
214
216{
218
219 if (pfes->GetFaceNbrVSize() <= 0)
220 {
221 return;
222 }
223
224 ParMesh *pmesh = pfes->GetParMesh();
225
228
229 int *send_offset = pfes->send_face_nbr_ldof.GetI();
230 const int *d_send_ldof = mfem::Read(pfes->send_face_nbr_ldof.GetJMemory(),
231 send_data.Size());
232 int *recv_offset = pfes->face_nbr_ldof.GetI();
233 MPI_Comm MyComm = pfes->GetComm();
234
235 const int num_face_nbrs = pmesh->GetNFaceNeighbors();
236 MPI_Request *requests = new MPI_Request[2*num_face_nbrs];
237 MPI_Request *send_requests = requests;
238 MPI_Request *recv_requests = requests + num_face_nbrs;
239 MPI_Status *statuses = new MPI_Status[num_face_nbrs];
240
241 auto d_data = this->Read();
242 auto d_send_data = send_data.Write();
243 mfem::forall(send_data.Size(), [=] MFEM_HOST_DEVICE (int i)
244 {
245 const int ldof = d_send_ldof[i];
246 d_send_data[i] = d_data[ldof >= 0 ? ldof : -1-ldof];
247 });
248
249 const bool mpi_gpu_aware = Device::GetGPUAwareMPI();
250 auto send_data_ptr = mpi_gpu_aware ? send_data.Read() : send_data.HostRead();
251 auto face_nbr_data_ptr = mpi_gpu_aware ? face_nbr_data.Write() :
253 // Wait for the kernel to be done since it updates what's sent and it may be async
254 if (mpi_gpu_aware) { MFEM_STREAM_SYNC; }
255 for (int fn = 0; fn < num_face_nbrs; fn++)
256 {
257 int nbr_rank = pmesh->GetFaceNbrRank(fn);
258 int tag = 0;
259
260 MPI_Isend(&send_data_ptr[send_offset[fn]],
261 send_offset[fn+1] - send_offset[fn],
262 MPITypeMap<real_t>::mpi_type, nbr_rank, tag, MyComm, &send_requests[fn]);
263
264 MPI_Irecv(&face_nbr_data_ptr[recv_offset[fn]],
265 recv_offset[fn+1] - recv_offset[fn],
266 MPITypeMap<real_t>::mpi_type, nbr_rank, tag, MyComm, &recv_requests[fn]);
267 }
268
269 MPI_Waitall(num_face_nbrs, send_requests, statuses);
270 MPI_Waitall(num_face_nbrs, recv_requests, statuses);
271
272 delete [] statuses;
273 delete [] requests;
274}
275
277const
278{
279 Array<int> dofs;
280 Vector DofVal, LocVec;
281 const int nbr_el_no = i - pfes->GetParMesh()->GetNE();
282 if (nbr_el_no >= 0)
283 {
284 int fes_vdim = pfes->GetVDim();
285 const DofTransformation* const doftrans = pfes->GetFaceNbrElementVDofs(
286 nbr_el_no, dofs);
287 // Choose fe to be of the order whose number of DOFs matches dofs.Size(),
288 // in the variable order case.
289 const int ndofs = pfes->IsVariableOrder() ? dofs.Size() : 0;
290 const FiniteElement *fe = pfes->GetFaceNbrFE(nbr_el_no, fes_vdim * ndofs);
291
292 if (fes_vdim > 1)
293 {
294 int s = dofs.Size()/fes_vdim;
295 Array<int> dofs_(&dofs[(vdim-1)*s], s);
296 face_nbr_data.GetSubVector(dofs_, LocVec);
297
298 DofVal.SetSize(s);
299 }
300 else
301 {
302 face_nbr_data.GetSubVector(dofs, LocVec);
303 DofVal.SetSize(dofs.Size());
304 }
305 if (doftrans)
306 {
307 doftrans->InvTransformPrimal(LocVec);
308 }
309
310 if (fe->GetMapType() == FiniteElement::VALUE)
311 {
312 fe->CalcShape(ip, DofVal);
313 }
314 else
315 {
318 Tr->SetIntPoint(&ip);
319 fe->CalcPhysShape(*Tr, DofVal);
320 }
321 }
322 else
323 {
324 const DofTransformation* const doftrans = fes->GetElementDofs(i, dofs);
325 fes->DofsToVDofs(vdim-1, dofs);
326 DofVal.SetSize(dofs.Size());
327 const FiniteElement *fe = fes->GetFE(i);
328 if (fe->GetMapType() == FiniteElement::VALUE)
329 {
330 fe->CalcShape(ip, DofVal);
331 }
332 else
333 {
335 Tr->SetIntPoint(&ip);
336 fe->CalcPhysShape(*Tr, DofVal);
337 }
338 GetSubVector(dofs, LocVec);
339 if (doftrans)
340 {
341 doftrans->InvTransformPrimal(LocVec);
342 }
343 }
344
345 return (DofVal * LocVec);
346}
347
349 Vector &val) const
350{
351 const int nbr_el_no = i - pfes->GetParMesh()->GetNE();
352 if (nbr_el_no >= 0)
353 {
354 Array<int> dofs;
355 const DofTransformation* const doftrans = pfes->GetFaceNbrElementVDofs(
356 nbr_el_no,
357 dofs);
358 Vector loc_data;
359 face_nbr_data.GetSubVector(dofs, loc_data);
360 if (doftrans)
361 {
362 doftrans->InvTransformPrimal(loc_data);
363 }
364 const FiniteElement *FElem = pfes->GetFaceNbrFE(nbr_el_no);
365 int dof = FElem->GetDof();
366 if (FElem->GetRangeType() == FiniteElement::SCALAR)
367 {
368 Vector shape(dof);
369 if (FElem->GetMapType() == FiniteElement::VALUE)
370 {
371 FElem->CalcShape(ip, shape);
372 }
373 else
374 {
377 Tr->SetIntPoint(&ip);
378 FElem->CalcPhysShape(*Tr, shape);
379 }
380 int vdim = fes->GetVDim();
381 val.SetSize(vdim);
382 for (int k = 0; k < vdim; k++)
383 {
384 val(k) = shape * (&loc_data[dof * k]);
385 }
386 }
387 else
388 {
389 int spaceDim = fes->GetMesh()->SpaceDimension();
390 DenseMatrix vshape(dof, spaceDim);
393 Tr->SetIntPoint(&ip);
394 FElem->CalcVShape(*Tr, vshape);
395 val.SetSize(spaceDim);
396 vshape.MultTranspose(loc_data, val);
397 }
398 }
399 else
400 {
402 }
403}
404
406 const IntegrationPoint &ip,
407 int comp, Vector *tr) const
408{
409 // We can assume faces and edges are local
411 {
412 return GridFunction::GetValue(T, ip, comp, tr);
413 }
414
415 // Check for evaluation in a local element
416 const int nbr_el_no = T.ElementNo - pfes->GetParMesh()->GetNE();
417 if (nbr_el_no < 0)
418 {
419 return GridFunction::GetValue(T, ip, comp, tr);
420 }
421
422 // Evaluate using DoFs from a neighboring element
423 if (tr)
424 {
425 T.SetIntPoint(&ip);
426 T.Transform(ip, *tr);
427 }
428
429 Array<int> dofs;
430 const FiniteElement * fe = pfes->GetFaceNbrFE(nbr_el_no);
431 const DofTransformation* const doftrans = pfes->GetFaceNbrElementVDofs(
432 nbr_el_no, dofs);
433
434 pfes->DofsToVDofs(comp-1, dofs);
435 Vector DofVal(dofs.Size()), LocVec;
436 if (fe->GetMapType() == FiniteElement::VALUE)
437 {
438 fe->CalcShape(ip, DofVal);
439 }
440 else
441 {
442 fe->CalcPhysShape(T, DofVal);
443 }
444 face_nbr_data.GetSubVector(dofs, LocVec);
445 if (doftrans)
446 {
447 doftrans->InvTransformPrimal(LocVec);
448 }
449
450
451 return (DofVal * LocVec);
452}
453
455 const IntegrationPoint &ip,
456 Vector &val, Vector *tr) const
457{
458 // We can assume faces and edges are local
460 {
461 return GridFunction::GetVectorValue(T, ip, val, tr);
462 }
463
464 // Check for evaluation in a local element
465 const int nbr_el_no = T.ElementNo - pfes->GetParMesh()->GetNE();
466 if (nbr_el_no < 0)
467 {
468 return GridFunction::GetVectorValue(T, ip, val, tr);
469 }
470
471 // Evaluate using DoFs from a neighboring element
472 if (tr)
473 {
474 T.SetIntPoint(&ip);
475 T.Transform(ip, *tr);
476 }
477
478 Array<int> vdofs;
479 DofTransformation * doftrans = pfes->GetFaceNbrElementVDofs(nbr_el_no, vdofs);
480 Vector loc_data;
481 face_nbr_data.GetSubVector(vdofs, loc_data);
482 if (doftrans)
483 {
484 doftrans->InvTransformPrimal(loc_data);
485 }
486
487 const FiniteElement *fe = pfes->GetFaceNbrFE(nbr_el_no);
488 const int dof = fe->GetDof();
490 {
491 Vector shape(dof);
492 if (fe->GetMapType() == FiniteElement::VALUE)
493 {
494 fe->CalcShape(ip, shape);
495 }
496 else
497 {
498 fe->CalcPhysShape(T, shape);
499 }
500 int vdim = pfes->GetVDim();
501 val.SetSize(vdim);
502 for (int k = 0; k < vdim; k++)
503 {
504 val(k) = shape * (&loc_data[dof * k]);
505 }
506 }
507 else
508 {
509 int spaceDim = pfes->GetMesh()->SpaceDimension();
510 int vdim = std::max(spaceDim, fe->GetRangeDim());
511 DenseMatrix vshape(dof, vdim);
512 fe->CalcVShape(T, vshape);
513 val.SetSize(vdim);
514 vshape.MultTranspose(loc_data, val);
515 }
516}
517
519{
521 // Count the zones globally.
522 GroupCommunicator &gcomm = this->ParFESpace()->GroupComm();
523 gcomm.Reduce<int>(elem_per_vdof, GroupCommunicator::Sum);
524 gcomm.Bcast(elem_per_vdof);
525}
526
527void ParGridFunction::GetDerivative(int comp, int der_comp,
528 ParGridFunction &der) const
529{
530 Array<int> overlap;
531 AccumulateAndCountDerivativeValues(comp, der_comp, der, overlap);
532
533 // Count the zones globally.
534 GroupCommunicator &gcomm = der.ParFESpace()->GroupComm();
535 gcomm.Reduce<int>(overlap, GroupCommunicator::Sum);
536 gcomm.Bcast(overlap);
537
538 // Accumulate for all dofs.
540 gcomm.Bcast<real_t>(der.HostReadWrite());
541
542 for (int i = 0; i < overlap.Size(); i++)
543 {
544 der(i) /= overlap[i];
545 }
546}
547
548void ParGridFunction::GetElementDofValues(int el, Vector &dof_vals) const
549{
550 int ne = fes->GetNE();
551 if (el >= ne)
552 {
553 MFEM_ASSERT(face_nbr_data.Size() > 0,
554 "ParGridFunction::GetElementDofValues: ExchangeFaceNbrData "
555 "must be called before accessing face neighbor elements.");
556 // Face neighbor element
557 Array<int> dof_idx;
558 pfes->GetFaceNbrElementVDofs(el - ne, dof_idx);
559 face_nbr_data.GetSubVector(dof_idx, dof_vals);
560 }
561 else
562 {
564 }
565}
566
568{
569 DeltaCoefficient *delta_c = dynamic_cast<DeltaCoefficient *>(&coeff);
570
571 if (delta_c == NULL)
572 {
574 }
575 else
576 {
577 real_t loc_integral, glob_integral;
578
579 ProjectDeltaCoefficient(*delta_c, loc_integral);
580
581 MPI_Allreduce(&loc_integral, &glob_integral, 1, MPITypeMap<real_t>::mpi_type,
582 MPI_SUM,
583 pfes->GetComm());
584
585 (*this) *= (delta_c->Scale() / glob_integral);
586 }
587}
588
590{
591 // local maximal element attribute for each dof
592 Array<int> ldof_attr;
593
594 // local projection
595 GridFunction::ProjectDiscCoefficient(coeff, ldof_attr);
596
597 // global maximal element attribute for each dof
598 Array<int> gdof_attr;
599 ldof_attr.Copy(gdof_attr);
600 GroupCommunicator &gcomm = pfes->GroupComm();
601 gcomm.Reduce<int>(gdof_attr, GroupCommunicator::Max);
602 gcomm.Bcast(gdof_attr);
603
604 // set local value to zero if global maximal element attribute is larger than
605 // the local one, and mark (in gdof_attr) if we have the correct value
606 for (int i = 0; i < pfes->GetVSize(); i++)
607 {
608 if (gdof_attr[i] > ldof_attr[i])
609 {
610 (*this)(i) = 0.0;
611 gdof_attr[i] = 0;
612 }
613 else
614 {
615 gdof_attr[i] = 1;
616 }
617 }
618
619 // parallel averaging plus interpolation to determine final values
621 gcomm.Reduce<int>(gdof_attr, GroupCommunicator::Sum);
622 gcomm.Bcast(gdof_attr);
623 for (int i = 0; i < fes->GetVSize(); i++)
624 {
625 (*this)(i) /= gdof_attr[i];
626 }
627 this->ParallelAssemble(*tv);
628 this->Distribute(tv);
629 delete tv;
630}
631
632
634{
635 // Harmonic (x1 ... xn) = [ (1/x1 + ... + 1/xn) / n ]^-1.
636 // Arithmetic(x1 ... xn) = (x1 + ... + xn) / n.
637
638 // Number of zones that contain a given dof.
639 Array<int> zones_per_vdof;
640 AccumulateAndCountZones(coeff, type, zones_per_vdof);
641
642 // Count the zones globally.
643 GroupCommunicator &gcomm = pfes->GroupComm();
644 gcomm.Reduce<int>(zones_per_vdof, GroupCommunicator::Sum);
645 gcomm.Bcast(zones_per_vdof);
646
647 // Accumulate for all vdofs.
649 gcomm.Bcast<real_t>(data);
650
651 ComputeMeans(type, zones_per_vdof);
652}
653
655 AvgType type)
656{
657 // Harmonic (x1 ... xn) = [ (1/x1 + ... + 1/xn) / n ]^-1.
658 // Arithmetic(x1 ... xn) = (x1 + ... + xn) / n.
659
660 // Number of zones that contain a given dof.
661 Array<int> zones_per_vdof;
662 AccumulateAndCountZones(vcoeff, type, zones_per_vdof);
663
664 // Count the zones globally.
665 GroupCommunicator &gcomm = pfes->GroupComm();
666 gcomm.Reduce<int>(zones_per_vdof, GroupCommunicator::Sum);
667 gcomm.Bcast(zones_per_vdof);
668
669 // Accumulate for all vdofs.
671 gcomm.Bcast<real_t>(data);
672
673 ComputeMeans(type, zones_per_vdof);
674}
675
677 Coefficient *coeff[], VectorCoefficient *vcoeff, const Array<int> &attr)
678{
679 Array<int> values_counter;
680 AccumulateAndCountBdrValues(coeff, vcoeff, attr, values_counter);
681
682 Vector values(Size());
683 for (int i = 0; i < values.Size(); i++)
684 {
685 values(i) = values_counter[i] ? (*this)(i) : 0.0;
686 }
687
688 // Count the values globally.
689 GroupCommunicator &gcomm = pfes->GroupComm();
690 gcomm.Reduce<int>(values_counter.HostReadWrite(), GroupCommunicator::Sum);
691 // Accumulate the values globally.
693
694 for (int i = 0; i < values.Size(); i++)
695 {
696 if (values_counter[i])
697 {
698 (*this)(i) = values(i)/values_counter[i];
699 }
700 }
701 // Broadcast values to other processors to have a consistent GridFunction
702 gcomm.Bcast<real_t>((*this).HostReadWrite());
703
704#ifdef MFEM_DEBUG
705 Array<int> ess_vdofs_marker;
706 if (vcoeff) { pfes->GetEssentialVDofs(attr, ess_vdofs_marker); }
707 else
708 {
709 ess_vdofs_marker.SetSize(Size());
710 ess_vdofs_marker = 0;
711 for (int i = 0; i < fes->GetVDim(); i++)
712 {
713 if (!coeff[i]) { continue; }
714 Array<int> component_dof_marker;
715 pfes->GetEssentialVDofs(attr, component_dof_marker,i);
716 for (int j = 0; j<Size(); j++)
717 {
718 ess_vdofs_marker[j] = bool(ess_vdofs_marker[j]) ||
719 bool(component_dof_marker[j]);
720 }
721 }
722 }
723 gcomm.Bcast<int>(values_counter.HostReadWrite());
724 for (int i = 0; i < values_counter.Size(); i++)
725 {
726 MFEM_ASSERT(bool(values_counter[i]) == bool(ess_vdofs_marker[i]),
727 "internal error");
728 }
729#endif
730}
731
733 const Array<int> &bdr_attr)
734{
735 Array<int> values_counter;
736 AccumulateAndCountBdrTangentValues(vcoeff, bdr_attr, values_counter);
737
738 Vector values(Size());
739 for (int i = 0; i < values.Size(); i++)
740 {
741 values(i) = values_counter[i] ? (*this)(i) : 0.0;
742 }
743
744 // Count the values globally.
745 GroupCommunicator &gcomm = pfes->GroupComm();
746 gcomm.Reduce<int>(values_counter.HostReadWrite(), GroupCommunicator::Sum);
747 // Accumulate the values globally.
749
750 for (int i = 0; i < values.Size(); i++)
751 {
752 if (values_counter[i])
753 {
754 (*this)(i) = values(i)/values_counter[i];
755 }
756 }
757 // Broadcast values to other processors to have a consistent GridFunction
758 gcomm.Bcast<real_t>((*this).HostReadWrite());
759
760#ifdef MFEM_DEBUG
761 Array<int> ess_vdofs_marker;
762 pfes->GetEssentialVDofs(bdr_attr, ess_vdofs_marker);
763 gcomm.Bcast<int>(values_counter.HostReadWrite());
764 for (int i = 0; i < values_counter.Size(); i++)
765 {
766 MFEM_ASSERT(bool(values_counter[i]) == bool(ess_vdofs_marker[i]),
767 "internal error: " << pfes->GetLocalTDofNumber(i) << ' ' << bool(
768 values_counter[i]));
769 }
770#endif
771}
772
774 Coefficient *ell_coeff,
775 JumpScaling jump_scaling,
776 const IntegrationRule *irs[]) const
777{
778 const_cast<ParGridFunction *>(this)->ExchangeFaceNbrData();
779
780 int fdof, intorder, k;
781 ElementTransformation *transf;
782 Vector shape, el_dofs, err_val, ell_coeff_val;
783 Array<int> vdofs;
785 real_t error = 0.0;
786
787 ParMesh *mesh = pfes->GetParMesh();
788
789 std::map<int,int> local_to_shared;
790 for (int i = 0; i < mesh->GetNSharedFaces(); ++i)
791 {
792 int i_local = mesh->GetSharedFace(i);
793 local_to_shared[i_local] = i;
794 }
795
796 for (int i = 0; i < mesh->GetNumFaces(); i++)
797 {
798 real_t shared_face_factor = 1.0;
799 bool shared_face = false;
800 int iel1, iel2, info1, info2;
801 mesh->GetFaceElements(i, &iel1, &iel2);
802 mesh->GetFaceInfos(i, &info1, &info2);
803
804 real_t h = mesh->GetElementSize(iel1);
805 intorder = fes->GetFE(iel1)->GetOrder();
806
807 FaceElementTransformations *face_elem_transf;
808 const FiniteElement *fe1, *fe2;
809 if (info2 >= 0 && iel2 < 0)
810 {
811 int ishared = local_to_shared[i];
812 face_elem_transf = mesh->GetSharedFaceTransformations(ishared);
813 iel2 = face_elem_transf->Elem2No - mesh->GetNE();
814 fe2 = pfes->GetFaceNbrFE(iel2);
815 if ( (k = fe2->GetOrder()) > intorder )
816 {
817 intorder = k;
818 }
819 shared_face = true;
820 shared_face_factor = 0.5;
821 h = std::min(h, mesh->GetFaceNbrElementSize(iel2));
822 }
823 else
824 {
825 if (iel2 >= 0)
826 {
827 fe2 = pfes->GetFE(iel2);
828 if ( (k = fe2->GetOrder()) > intorder )
829 {
830 intorder = k;
831 }
832 h = std::min(h, mesh->GetElementSize(iel2));
833 }
834 else
835 {
836 fe2 = NULL;
837 }
838 face_elem_transf = mesh->GetFaceElementTransformations(i);
839 }
840 int p = intorder;
841
842 intorder = 2 * intorder; // <-------------
843 const IntegrationRule *ir;
844 if (irs)
845 {
846 ir = irs[face_elem_transf->GetGeometryType()];
847 }
848 else
849 {
850 ir = &(IntRules.Get(face_elem_transf->GetGeometryType(), intorder));
851 }
852 err_val.SetSize(ir->GetNPoints());
853 ell_coeff_val.SetSize(ir->GetNPoints());
854 // side 1
855 transf = face_elem_transf->Elem1;
856 fe1 = fes->GetFE(iel1);
857 fdof = fe1->GetDof();
858 fes->GetElementVDofs(iel1, vdofs);
859 shape.SetSize(fdof);
860 el_dofs.SetSize(fdof);
861 for (k = 0; k < fdof; k++)
862 if (vdofs[k] >= 0)
863 {
864 el_dofs(k) = (*this)(vdofs[k]);
865 }
866 else
867 {
868 el_dofs(k) = - (*this)(-1-vdofs[k]);
869 }
870 for (int j = 0; j < ir->GetNPoints(); j++)
871 {
872 face_elem_transf->Loc1.Transform(ir->IntPoint(j), eip);
873 fe1->CalcShape(eip, shape);
874 transf->SetIntPoint(&eip);
875 ell_coeff_val(j) = ell_coeff->Eval(*transf, eip);
876 err_val(j) = exsol->Eval(*transf, eip) - (shape * el_dofs);
877 }
878 if (fe2 != NULL)
879 {
880 // side 2
881 transf = face_elem_transf->Elem2;
882 fdof = fe2->GetDof();
883 shape.SetSize(fdof);
884 el_dofs.SetSize(fdof);
885 if (shared_face)
886 {
887 pfes->GetFaceNbrElementVDofs(iel2, vdofs);
888 for (k = 0; k < fdof; k++)
889 if (vdofs[k] >= 0)
890 {
891 el_dofs(k) = face_nbr_data[vdofs[k]];
892 }
893 else
894 {
895 el_dofs(k) = - face_nbr_data[-1-vdofs[k]];
896 }
897 }
898 else
899 {
900 pfes->GetElementVDofs(iel2, vdofs);
901 for (k = 0; k < fdof; k++)
902 if (vdofs[k] >= 0)
903 {
904 el_dofs(k) = (*this)(vdofs[k]);
905 }
906 else
907 {
908 el_dofs(k) = - (*this)(-1 - vdofs[k]);
909 }
910 }
911 for (int j = 0; j < ir->GetNPoints(); j++)
912 {
913 face_elem_transf->Loc2.Transform(ir->IntPoint(j), eip);
914 fe2->CalcShape(eip, shape);
915 transf->SetIntPoint(&eip);
916 ell_coeff_val(j) += ell_coeff->Eval(*transf, eip);
917 ell_coeff_val(j) *= 0.5;
918 err_val(j) -= (exsol->Eval(*transf, eip) - (shape * el_dofs));
919 }
920 }
921 real_t face_error = 0.0;
922 transf = face_elem_transf;
923 for (int j = 0; j < ir->GetNPoints(); j++)
924 {
925 const IntegrationPoint &ip = ir->IntPoint(j);
926 transf->SetIntPoint(&ip);
927 real_t nu = jump_scaling.Eval(h, p);
928 face_error += shared_face_factor*(ip.weight * nu * ell_coeff_val(j) *
929 transf->Weight() *
930 err_val(j) * err_val(j));
931 }
932 // negative quadrature weights may cause the error to be negative
933 error += fabs(face_error);
934 }
935
936 error = sqrt(error);
937 return GlobalLpNorm(2.0, error, pfes->GetComm());
938}
939
940void ParGridFunction::Save(std::ostream &os) const
941{
942 real_t *data_ = const_cast<real_t*>(HostRead());
943 for (int i = 0; i < size; i++)
944 {
945 if (pfes->GetDofSign(i) < 0) { data_[i] = -data_[i]; }
946 }
947
949
950 for (int i = 0; i < size; i++)
951 {
952 if (pfes->GetDofSign(i) < 0) { data_[i] = -data_[i]; }
953 }
954}
955
956void ParGridFunction::Save(const char *fname, int precision) const
957{
958 int rank = pfes->GetMyRank();
959 ostringstream fname_with_suffix;
960 fname_with_suffix << fname << "." << setfill('0') << setw(6) << rank;
961 ofstream ofs(fname_with_suffix.str().c_str());
962 ofs.precision(precision);
963 Save(ofs);
964}
965
966void ParGridFunction::SaveAsOne(const char *fname, int precision) const
967{
968 ofstream ofs;
969 int rank = pfes->GetMyRank();
970 if (rank == 0)
971 {
972 ofs.open(fname);
973 ofs.precision(precision);
974 }
975 SaveAsOne(ofs);
976}
977
978void ParGridFunction::SaveAsSerial(const char *fname, int precision,
979 int save_rank) const
980{
981 ParMesh *pmesh = ParFESpace()->GetParMesh();
982 Mesh serial_mesh = pmesh->GetSerialMesh(save_rank);
983 GridFunction serialgf = GetSerialGridFunction(save_rank, serial_mesh);
984
985 if (pmesh->GetMyRank() == save_rank)
986 {
987 serialgf.Save(fname, precision);
988 }
989 MPI_Barrier(pmesh->GetComm());
990}
991
993 int save_rank, FiniteElementSpace &serial_fes) const
994{
995 ParFiniteElementSpace *pfespace = ParFESpace();
996 ParMesh *pmesh = pfespace->GetParMesh();
997
998 GridFunction serial_gf(&serial_fes);
999
1000 Array<real_t> vals;
1001 Array<int> dofs;
1002 MPI_Status status;
1003
1004 const int vdim = pfespace->GetVDim();
1005
1006 const int my_rank = pmesh->GetMyRank();
1007 const int nranks = pmesh->GetNRanks();
1008 MPI_Comm comm = pmesh->GetComm();
1009
1010 if (my_rank == save_rank)
1011 {
1012 int elem_count = 0; // To keep track of element count in serial mesh
1013
1014 Vector nodeval;
1015 for (int e = 0; e < pmesh->GetNE(); e++)
1016 {
1017 GetElementDofValues(e, nodeval);
1018 serial_fes.GetElementVDofs(elem_count++, dofs);
1019 serial_gf.SetSubVector(dofs, nodeval);
1020 }
1021
1022 for (int p = 0; p < nranks; p++)
1023 {
1024 if (p == save_rank) { continue; }
1025 int n_send_recv;
1026 MPI_Recv(&n_send_recv, 1, MPI_INT, p, 448, comm, &status);
1027 vals.SetSize(n_send_recv);
1028 if (n_send_recv)
1029 {
1030 MPI_Recv(&vals[0], n_send_recv, MPITypeMap<real_t>::mpi_type, p, 449, comm,
1031 &status);
1032 }
1033 for (int i = 0; i < n_send_recv; )
1034 {
1035 serial_fes.GetElementVDofs(elem_count++, dofs);
1036 serial_gf.SetSubVector(dofs, &vals[i]);
1037 i += dofs.Size();
1038 }
1039 }
1040 } // my_rank == save_rank
1041 else
1042 {
1043 int n_send_recv = 0;
1044 Vector nodeval;
1045 for (int e = 0; e < pmesh->GetNE(); e++)
1046 {
1047 const FiniteElement *fe = pfespace->GetFE(e);
1048 n_send_recv += vdim*fe->GetDof();
1049 }
1050 MPI_Send(&n_send_recv, 1, MPI_INT, save_rank, 448, comm);
1051 vals.Reserve(n_send_recv);
1052 vals.SetSize(0);
1053 for (int e = 0; e < pmesh->GetNE(); e++)
1054 {
1055 GetElementDofValues(e, nodeval);
1056 for (int j = 0; j < nodeval.Size(); j++)
1057 {
1058 vals.Append(nodeval(j));
1059 }
1060 }
1061 if (n_send_recv)
1062 {
1063 MPI_Send(&vals[0], n_send_recv, MPITypeMap<real_t>::mpi_type, save_rank, 449,
1064 comm);
1065 }
1066 }
1067
1068 return serial_gf;
1069}
1070
1072 Mesh &serial_mesh) const
1073{
1074 auto *serial_fec = pfes->FEColl()->Clone(pfes->FEColl()->GetOrder());
1075 auto *serial_fes = new FiniteElementSpace(&serial_mesh,
1076 serial_fec,
1077 pfes->GetVDim(),
1078 pfes->GetOrdering());
1079 GridFunction serial_gf = GetSerialGridFunction(save_rank, *serial_fes);
1080 serial_gf.MakeOwner(serial_fec); // Also assumes ownership of serial_fes
1081 return serial_gf;
1082}
1083
1084#ifdef MFEM_USE_ADIOS2
1086 const std::string& variable_name,
1087 const adios2stream::data_type type) const
1088{
1089 real_t *data_ = const_cast<real_t*>(HostRead());
1090 for (int i = 0; i < size; i++)
1091 {
1092 if (pfes->GetDofSign(i) < 0) { data_[i] = -data_[i]; }
1093 }
1094
1095 GridFunction::Save(os, variable_name, type);
1096
1097 for (int i = 0; i < size; i++)
1098 {
1099 if (pfes->GetDofSign(i) < 0) { data_[i] = -data_[i]; }
1100 }
1101}
1102#endif
1103
1104void ParGridFunction::SaveAsOne(std::ostream &os) const
1105{
1106 int i, p;
1107
1108 MPI_Comm MyComm;
1109 MPI_Status status;
1110 int MyRank, NRanks;
1111
1112 MyComm = pfes -> GetComm();
1113
1114 MPI_Comm_size(MyComm, &NRanks);
1115 MPI_Comm_rank(MyComm, &MyRank);
1116
1117 real_t **values = new real_t*[NRanks];
1118 int *nv = new int[NRanks];
1119 int *nvdofs = new int[NRanks];
1120 int *nedofs = new int[NRanks];
1121 int *nfdofs = new int[NRanks];
1122 int *nrdofs = new int[NRanks];
1123
1124 real_t * h_data = const_cast<real_t *>(this->HostRead());
1125
1126 values[0] = h_data;
1127 nv[0] = pfes -> GetVSize();
1128 nvdofs[0] = pfes -> GetNVDofs();
1129 nedofs[0] = pfes -> GetNEDofs();
1130 nfdofs[0] = pfes -> GetNFDofs();
1131
1132 if (MyRank == 0)
1133 {
1134 pfes -> Save(os);
1135 os << '\n';
1136
1137 for (p = 1; p < NRanks; p++)
1138 {
1139 MPI_Recv(&nv[p], 1, MPI_INT, p, 455, MyComm, &status);
1140 MPI_Recv(&nvdofs[p], 1, MPI_INT, p, 456, MyComm, &status);
1141 MPI_Recv(&nedofs[p], 1, MPI_INT, p, 457, MyComm, &status);
1142 MPI_Recv(&nfdofs[p], 1, MPI_INT, p, 458, MyComm, &status);
1143 values[p] = new real_t[nv[p]];
1144 MPI_Recv(values[p], nv[p], MPITypeMap<real_t>::mpi_type, p, 460, MyComm,
1145 &status);
1146 }
1147
1148 int vdim = pfes -> GetVDim();
1149
1150 for (p = 0; p < NRanks; p++)
1151 {
1152 nrdofs[p] = nv[p]/vdim - nvdofs[p] - nedofs[p] - nfdofs[p];
1153 }
1154
1156 {
1157 for (int d = 0; d < vdim; d++)
1158 {
1159 for (p = 0; p < NRanks; p++)
1160 for (i = 0; i < nvdofs[p]; i++)
1161 {
1162 os << *values[p]++ << '\n';
1163 }
1164
1165 for (p = 0; p < NRanks; p++)
1166 for (i = 0; i < nedofs[p]; i++)
1167 {
1168 os << *values[p]++ << '\n';
1169 }
1170
1171 for (p = 0; p < NRanks; p++)
1172 for (i = 0; i < nfdofs[p]; i++)
1173 {
1174 os << *values[p]++ << '\n';
1175 }
1176
1177 for (p = 0; p < NRanks; p++)
1178 for (i = 0; i < nrdofs[p]; i++)
1179 {
1180 os << *values[p]++ << '\n';
1181 }
1182 }
1183 }
1184 else
1185 {
1186 for (p = 0; p < NRanks; p++)
1187 for (i = 0; i < nvdofs[p]; i++)
1188 for (int d = 0; d < vdim; d++)
1189 {
1190 os << *values[p]++ << '\n';
1191 }
1192
1193 for (p = 0; p < NRanks; p++)
1194 for (i = 0; i < nedofs[p]; i++)
1195 for (int d = 0; d < vdim; d++)
1196 {
1197 os << *values[p]++ << '\n';
1198 }
1199
1200 for (p = 0; p < NRanks; p++)
1201 for (i = 0; i < nfdofs[p]; i++)
1202 for (int d = 0; d < vdim; d++)
1203 {
1204 os << *values[p]++ << '\n';
1205 }
1206
1207 for (p = 0; p < NRanks; p++)
1208 for (i = 0; i < nrdofs[p]; i++)
1209 for (int d = 0; d < vdim; d++)
1210 {
1211 os << *values[p]++ << '\n';
1212 }
1213 }
1214
1215 for (p = 1; p < NRanks; p++)
1216 {
1217 values[p] -= nv[p];
1218 delete [] values[p];
1219 }
1220 os.flush();
1221 }
1222 else
1223 {
1224 MPI_Send(&nv[0], 1, MPI_INT, 0, 455, MyComm);
1225 MPI_Send(&nvdofs[0], 1, MPI_INT, 0, 456, MyComm);
1226 MPI_Send(&nedofs[0], 1, MPI_INT, 0, 457, MyComm);
1227 MPI_Send(&nfdofs[0], 1, MPI_INT, 0, 458, MyComm);
1228 MPI_Send(h_data, nv[0], MPITypeMap<real_t>::mpi_type, 0, 460, MyComm);
1229 }
1230
1231 delete [] values;
1232 delete [] nv;
1233 delete [] nvdofs;
1234 delete [] nedofs;
1235 delete [] nfdofs;
1236 delete [] nrdofs;
1237}
1238
1239real_t GlobalLpNorm(const real_t p, real_t loc_norm, MPI_Comm comm)
1240{
1241 real_t glob_norm;
1242
1243 // negative quadrature weights may cause the local norm to be negative
1244 loc_norm = fabs(loc_norm);
1245
1246 if (p < infinity())
1247 {
1248 loc_norm = pow(loc_norm, p);
1249
1250 MPI_Allreduce(&loc_norm, &glob_norm, 1, MPITypeMap<real_t>::mpi_type,
1251 MPI_SUM, comm);
1252
1253 glob_norm = pow(fabs(glob_norm), 1.0/p);
1254 }
1255 else
1256 {
1257 MPI_Allreduce(&loc_norm, &glob_norm, 1, MPITypeMap<real_t>::mpi_type,
1258 MPI_MAX, comm);
1259 }
1260
1261 return glob_norm;
1262}
1263
1266 GridFunction &flux, bool wcoef, int subdomain)
1267{
1268 ParFiniteElementSpace *ffes =
1269 dynamic_cast<ParFiniteElementSpace*>(flux.FESpace());
1270 MFEM_VERIFY(ffes, "the flux FE space must be ParFiniteElementSpace");
1271
1272 Array<int> count(flux.Size());
1273 SumFluxAndCount(blfi, flux, count, wcoef, subdomain);
1274
1275 // Accumulate flux and counts in parallel
1277 ffes->GroupComm().Bcast<real_t>(flux.HostReadWrite());
1278
1280 ffes->GroupComm().Bcast<int>(count.HostReadWrite());
1281
1282 // complete averaging
1283 for (int i = 0; i < count.Size(); i++)
1284 {
1285 if (count[i] != 0) { flux(i) /= count[i]; }
1286 }
1287
1288 if (ffes->Nonconforming())
1289 {
1290 // On a partially conforming flux space, project on the conforming space.
1291 // Using this code may lead to worse refinements in ex6, so we do not use
1292 // it by default.
1293
1294 // Vector conf_flux;
1295 // flux.ConformingProject(conf_flux);
1296 // flux.ConformingProlongate(conf_flux);
1297 }
1298}
1299
1300std::unique_ptr<ParGridFunction> ParGridFunction::ProlongateToMaxOrder() const
1301{
1302 ParMesh *mesh = pfes->GetParMesh();
1303 const FiniteElementCollection *pfesc = pfes->FEColl();
1304 const int vdim = pfes->GetVDim();
1305
1306 // Find the max order in the space
1307 const int maxOrder = pfes->GetMaxElementOrder();
1308
1309 // Create a visualization space of max order for all elements
1310 FiniteElementCollection *fecMax = pfesc->Clone(maxOrder);
1311 ParFiniteElementSpace *pfesMax = new ParFiniteElementSpace(mesh, fecMax, vdim,
1312 pfes->GetOrdering());
1313
1314 ParGridFunction *xMax = new ParGridFunction(pfesMax);
1315
1316 // Interpolate in the maximum-order space
1317 PRefinementTransferOperator P(*pfes, *pfesMax);
1318 P.Mult(*this, *xMax);
1319
1320 xMax->MakeOwner(fecMax);
1321 return std::unique_ptr<ParGridFunction>(xMax);
1322}
1323
1325 const ParGridFunction &x,
1326 ParFiniteElementSpace &smooth_flux_fes,
1327 ParFiniteElementSpace &flux_fes,
1328 Vector &errors,
1329 int norm_p, real_t solver_tol, int solver_max_it)
1330{
1331 // Compute fluxes in discontinuous space
1332 GridFunction flux(&flux_fes);
1333 flux = 0.0;
1334
1336 Array<int> xdofs, fdofs;
1337 Vector el_x, el_f;
1338
1339 for (int i = 0; i < xfes->GetNE(); i++)
1340 {
1341 const DofTransformation* const xtrans = xfes->GetElementVDofs(i, xdofs);
1342 x.GetSubVector(xdofs, el_x);
1343 if (xtrans)
1344 {
1345 xtrans->InvTransformPrimal(el_x);
1346 }
1347
1349 flux_integrator.ComputeElementFlux(*xfes->GetFE(i), *Transf, el_x,
1350 *flux_fes.GetFE(i), el_f, false);
1351
1352 const DofTransformation* const ftrans = flux_fes.GetElementVDofs(i, fdofs);
1353 if (ftrans)
1354 {
1355 ftrans->TransformPrimal(el_f);
1356 }
1357 flux.SetSubVector(fdofs, el_f);
1358 }
1359
1360 // Assemble the linear system for L2 projection into the "smooth" space
1361 ParBilinearForm *a = new ParBilinearForm(&smooth_flux_fes);
1362 ParLinearForm *b = new ParLinearForm(&smooth_flux_fes);
1364
1365 const FiniteElement *smooth_flux_fe = smooth_flux_fes.GetTypicalFE();
1366
1367 if (smooth_flux_fe->GetRangeType() == FiniteElement::SCALAR)
1368 {
1370 vmass->SetVDim(smooth_flux_fes.GetVDim());
1371 a->AddDomainIntegrator(vmass);
1372 b->AddDomainIntegrator(new VectorDomainLFIntegrator(f));
1373 }
1374 else
1375 {
1376 a->AddDomainIntegrator(new VectorFEMassIntegrator);
1377 b->AddDomainIntegrator(new VectorFEDomainLFIntegrator(f));
1378 }
1379
1380 b->Assemble();
1381 a->Assemble();
1382 a->Finalize();
1383
1384 // The destination of the projected discontinuous flux
1385 ParGridFunction smooth_flux(&smooth_flux_fes);
1386 smooth_flux = 0.0;
1387
1388 HypreParMatrix* A = a->ParallelAssemble();
1389 HypreParVector* B = b->ParallelAssemble();
1390 HypreParVector* X = smooth_flux.ParallelProject();
1391
1392 delete a;
1393 delete b;
1394
1395 // Define and apply a parallel PCG solver for AX=B with the BoomerAMG
1396 // preconditioner from hypre.
1397 HypreBoomerAMG *amg = new HypreBoomerAMG(*A);
1398 amg->SetPrintLevel(0);
1399 HyprePCG *pcg = new HyprePCG(*A);
1400 pcg->SetTol(solver_tol);
1401 pcg->SetMaxIter(solver_max_it);
1402 pcg->SetPrintLevel(0);
1403 pcg->SetPreconditioner(*amg);
1404 pcg->Mult(*B, *X);
1405
1406 // Extract the parallel grid function corresponding to the finite element
1407 // approximation X. This is the local solution on each processor.
1408 smooth_flux = *X;
1409
1410 delete A;
1411 delete B;
1412 delete X;
1413 delete amg;
1414 delete pcg;
1415
1416 // Proceed through the elements one by one, and find the Lp norm differences
1417 // between the flux as computed per element and the flux projected onto the
1418 // smooth_flux_fes space.
1419 real_t total_error = 0.0;
1420 errors.SetSize(xfes->GetNE());
1421 for (int i = 0; i < xfes->GetNE(); i++)
1422 {
1423 errors(i) = ComputeElementLpDistance(norm_p, i, smooth_flux, flux);
1424 total_error += pow(errors(i), norm_p);
1425 }
1426
1427 real_t glob_error;
1428 MPI_Allreduce(&total_error, &glob_error, 1, MPITypeMap<real_t>::mpi_type,
1429 MPI_SUM,
1430 xfes->GetComm());
1431
1432 return pow(glob_error, 1.0/norm_p);
1433}
1434
1435} // namespace mfem
1436
1437#endif // MFEM_USE_MPI
void Reserve(int capacity)
Ensures that the allocated size is at least the given size.
Definition array.hpp:165
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition array.hpp:758
int Size() const
Return the logical size of the array.
Definition array.hpp:147
int Append(const T &el)
Append element 'el' to array, resize if necessary.
Definition array.hpp:830
void Copy(Array &copy) const
Create a copy of the internal array to the provided copy.
Definition array.hpp:935
T * HostReadWrite()
Shortcut for mfem::ReadWrite(a.GetMemory(), a.Size(), false).
Definition array.hpp:357
Abstract base class BilinearFormIntegrator.
virtual void ComputeElementFlux(const FiniteElement &el, ElementTransformation &Trans, Vector &u, const FiniteElement &fluxelem, Vector &flux, bool with_coef=true, const IntegrationRule *ir=NULL)
Virtual method required for Zienkiewicz-Zhu type error estimators.
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the coefficient in the element described by T at the point ip.
Delta function coefficient optionally multiplied by a weight coefficient and a scaled time dependent ...
real_t Scale()
Return the scale factor times the optional time dependent function. Returns with when not set by th...
Data type dense matrix using column-major storage.
Definition densemat.hpp:24
void MultTranspose(const real_t *x, real_t *y) const
Multiply a vector with the transpose matrix.
Definition densemat.cpp:172
static bool GetGPUAwareMPI()
Definition device.hpp:291
void InvTransformPrimal(real_t *v) const
Definition doftrans.cpp:49
void TransformPrimal(real_t *v) const
Definition doftrans.cpp:17
Geometry::Type GetGeometryType() const
Return the Geometry::Type of the reference element.
Definition eltrans.hpp:175
real_t Weight()
Return the weight of the Jacobian matrix of the transformation at the currently set IntegrationPoint....
Definition eltrans.hpp:144
void SetIntPoint(const IntegrationPoint *ip)
Set the integration point ip that weights and Jacobians will be evaluated at.
Definition eltrans.hpp:106
virtual void Transform(const IntegrationPoint &, Vector &)=0
Transform integration point from reference coordinates to physical coordinates and store them in the ...
A specialized ElementTransformation class representing a face and its two neighboring elements.
Definition eltrans.hpp:750
ElementTransformation * Elem2
Definition eltrans.hpp:791
ElementTransformation * Elem1
Definition eltrans.hpp:791
IntegrationPointTransformation Loc1
Definition eltrans.hpp:793
IntegrationPointTransformation Loc2
Definition eltrans.hpp:793
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition fe_coll.hpp:27
static FiniteElementCollection * New(const char *name)
Factory method: return a newly allocated FiniteElementCollection according to the given name.
Definition fe_coll.cpp:124
int GetOrder() const
Return the order (polynomial degree) of the FE collection, corresponding to the order/degree returned...
Definition fe_coll.hpp:240
virtual FiniteElementCollection * Clone(int p) const
Instantiate a new collection of the same type with a different order.
Definition fe_coll.cpp:411
virtual const char * Name() const
Definition fe_coll.hpp:79
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:244
bool IsVariableOrder() const
Returns true if the space contains elements of varying polynomial orders.
Definition fespace.hpp:709
void DofsToVDofs(Array< int > &dofs, int ndofs=-1) const
Compute the full set of vdofs corresponding to each entry in dofs.
Definition fespace.cpp:258
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:3513
ElementTransformation * GetElementTransformation(int i) const
Returns ElementTransformation for the i-th element.
Definition fespace.hpp:924
DofTransformation * GetElementVDofs(int i, Array< int > &vdofs) const
Returns indices of degrees of freedom for the i'th element. The returned indices are offsets into an ...
Definition fespace.cpp:332
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:3835
Ordering::Type GetOrdering() const
Return the ordering method.
Definition fespace.hpp:876
int GetNE() const
Returns number of elements in the mesh.
Definition fespace.hpp:891
const FiniteElementCollection * FEColl() const
Definition fespace.hpp:878
Mesh * GetMesh() const
Returns the mesh.
Definition fespace.hpp:679
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition fespace.hpp:848
int GetVDim() const
Returns the vector dimension of the finite element space.
Definition fespace.hpp:841
const FiniteElement * GetTypicalFE() const
Return GetFE(0) if the local mesh is not empty; otherwise return a typical FE based on the Geometry t...
Definition fespace.cpp:3871
Abstract class for all finite elements.
Definition fe_base.hpp:244
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Evaluate the values of all shape functions of a vector finite element in reference space at the given...
Definition fe_base.cpp:40
int GetRangeDim() const
Returns the vector dimension for vector-valued finite elements, which is also the dimension of the in...
Definition fe_base.hpp:325
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Definition fe_base.hpp:338
int GetMapType() const
Returns the FiniteElement::MapType of the element describing how reference functions are mapped to ph...
Definition fe_base.hpp:360
int GetRangeType() const
Returns the FiniteElement::RangeType of the element, one of {SCALAR, VECTOR}.
Definition fe_base.hpp:351
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:334
void CalcPhysShape(ElementTransformation &Trans, Vector &shape) const
Evaluate the values of all shape functions of a scalar finite element in physical space at the point ...
Definition fe_base.cpp:182
Class for grid function - Vector with associated FE space.
Definition gridfunc.hpp:31
void ComputeMeans(AvgType type, Array< int > &zones_per_vdof)
void AccumulateAndCountBdrTangentValues(VectorCoefficient &vcoeff, const Array< int > &bdr_attr, Array< int > &values_counter)
virtual void CountElementsPerVDof(Array< int > &elem_per_vdof) const
For each vdof, counts how many elements contain the vdof, as containment is determined by FiniteEleme...
virtual real_t GetValue(int i, const IntegrationPoint &ip, int vdim=1) const
Definition gridfunc.cpp:446
void AccumulateAndCountBdrValues(Coefficient *coeff[], VectorCoefficient *vcoeff, const Array< int > &attr, Array< int > &values_counter)
virtual void Update()
Transform by the Space UpdateMatrix (e.g., on Mesh change).
Definition gridfunc.cpp:167
virtual void MakeRef(FiniteElementSpace *f, real_t *v)
Make the GridFunction reference external data on a new FiniteElementSpace.
Definition gridfunc.cpp:233
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
void MakeOwner(FiniteElementCollection *fec_)
Make the GridFunction the owner of fec_owned and fes.
Definition gridfunc.hpp:122
virtual void GetElementDofValues(int el, Vector &dof_vals) const
FiniteElementSpace * FESpace()
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
void ProjectDeltaCoefficient(DeltaCoefficient &delta_coeff, real_t &integral)
void ProjectDiscCoefficient(VectorCoefficient &coeff, Array< int > &dof_attr)
FiniteElementSpace * fes
FE space on which the grid function lives. Owned if fec_owned is not NULL.
Definition gridfunc.hpp:34
FiniteElementCollection * fec_owned
Used when the grid function is read from a file. It can also be set explicitly, see MakeOwner().
Definition gridfunc.hpp:40
void SumFluxAndCount(BilinearFormIntegrator &blfi, GridFunction &flux, Array< int > &counts, bool wcoef, int subdomain)
Definition gridfunc.cpp:281
void AccumulateAndCountDerivativeValues(int comp, int der_comp, GridFunction &der, Array< int > &zones_per_dof) const
Used for the serial and parallel implementations of the GetDerivative() method; see its documentation...
virtual void GetVectorValue(int i, const IntegrationPoint &ip, Vector &val) const
Definition gridfunc.cpp:473
void AccumulateAndCountZones(Coefficient &coeff, AvgType type, Array< int > &zones_per_vdof)
Accumulates (depending on type) the values of coeff at all shared vdofs and counts in how many zones ...
virtual void SetSpace(FiniteElementSpace *f)
Associate a new FiniteElementSpace with the GridFunction.
Definition gridfunc.cpp:225
Communicator performing operations within groups defined by a GroupTopology with arbitrary-size data ...
void Reduce(T *ldata, void(*Op)(OpData< T >)) const
Reduce within each group where the master is the root.
static void Sum(OpData< T >)
Reduce operation Sum, instantiated for int and double.
void Bcast(T *ldata, int layout) const
Broadcast within each group where the master is the root.
static void Max(OpData< T >)
Reduce operation Max, instantiated for int and double.
The BoomerAMG solver in hypre.
Definition hypre.hpp:1717
void SetPrintLevel(int print_level)
Definition hypre.hpp:1797
PCG solver in hypre.
Definition hypre.hpp:1301
void Mult(const HypreParVector &b, HypreParVector &x) const override
Solve Ax=b with hypre's PCG.
Definition hypre.cpp:4242
void SetPrintLevel(int print_lvl)
Definition hypre.cpp:4214
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
Definition hypre.cpp:4219
void SetMaxIter(int max_iter)
Definition hypre.cpp:4204
void SetTol(real_t tol)
Definition hypre.cpp:4194
Wrapper for hypre's ParCSR matrix class.
Definition hypre.hpp:408
HYPRE_Int Mult(HypreParVector &x, HypreParVector &y, real_t alpha=1.0, real_t beta=0.0) const
Computes y = alpha * A * x + beta * y.
Definition hypre.cpp:1861
Wrapper for hypre's parallel vector class.
Definition hypre.hpp:219
void Transform(const IntegrationPoint &, IntegrationPoint &)
Definition eltrans.cpp:587
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.
real_t Eval(real_t h, int p) const
Mesh data type.
Definition mesh.hpp:64
void GetFaceInfos(int Face, int *Inf1, int *Inf2) const
Definition mesh.cpp:1463
int GetNumFaces() const
Return the number of faces (3D), edges (2D) or vertices (1D).
Definition mesh.cpp:6531
int GetNE() const
Returns number of elements.
Definition mesh.hpp:1282
real_t GetElementSize(int i, int type=0)
Get the size of the i-th element relative to the perfect reference element.
Definition mesh.cpp:107
void GetFaceElements(int Face, int *Elem1, int *Elem2) const
Definition mesh.cpp:1457
int SpaceDimension() const
Dimension of the physical space containing the mesh.
Definition mesh.hpp:1219
Abstract operator.
Definition operator.hpp:25
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
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:567
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...
Class for parallel bilinear form.
Abstract parallel finite element space.
Definition pfespace.hpp:29
MPI_Comm GetComm() const
Definition pfespace.hpp:334
int GetMaxElementOrder() const override
Returns the maximum polynomial order over all elements globally.
int GetLocalTDofNumber(int ldof) const
void DivideByGroupSize(real_t *vec)
Scale a vector of true dofs.
HypreParVector * NewTrueDofVector()
Definition pfespace.hpp:398
const FiniteElement * GetFaceNbrFE(int i, int ndofs=0) const
GroupCommunicator & GroupComm()
Return a reference to the internal GroupCommunicator (on VDofs)
Definition pfespace.hpp:405
const Operator * GetProlongationMatrix() const override
void GetEssentialVDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_dofs, int component=-1) const override
Determine the boundary degrees of freedom.
HypreParMatrix * Dof_TrueDof_Matrix() const
The true dof-to-dof interpolation matrix.
Definition pfespace.hpp:388
void GetFaceNbrElementVDofs(int i, Array< int > &vdofs, DofTransformation &doftrans) const
const SparseMatrix * GetRestrictionMatrix() const override
Get the R matrix which restricts a local dof vector to true dof vector.
Definition pfespace.hpp:467
ParMesh * GetParMesh() const
Definition pfespace.hpp:338
const FiniteElement * GetFE(int i) const override
Definition pfespace.cpp:627
ElementTransformation * GetFaceNbrElementTransformation(int i) const
Definition pfespace.hpp:482
Class for parallel grid function.
Definition pgridfunc.hpp:50
void CountElementsPerVDof(Array< int > &elem_per_vdof) const override
For each vdof, counts how many elements contain the vdof, as containment is determined by FiniteEleme...
void GetDerivative(int comp, int der_comp, ParGridFunction &der) const
Parallel version of GridFunction::GetDerivative(); see its documentation.
HypreParVector * ParallelAverage() const
Returns a new vector averaged on the true dofs.
real_t ComputeDGFaceJumpError(Coefficient *exsol, Coefficient *ell_coeff, JumpScaling jump_scaling, const IntegrationRule *irs[]=NULL) const override
Returns the Face Jumps error for L2 elements.
void Save(std::ostream &out) const override
void ProjectDiscCoefficient(VectorCoefficient &coeff) override
Project a discontinuous vector coefficient as a grid function on a continuous finite element space....
real_t GetValue(int i, const IntegrationPoint &ip, int vdim=1) const override
HypreParVector * GetTrueDofs() const
Returns the true dofs in a new HypreParVector.
void ComputeFlux(BilinearFormIntegrator &blfi, GridFunction &flux, bool wcoef=true, int subdomain=-1) override
void ProjectBdrCoefficient(Coefficient *coeff[], VectorCoefficient *vcoeff, const Array< int > &attr)
ParFiniteElementSpace * pfes
Points to the same object as fes.
Definition pgridfunc.hpp:52
Vector send_data
Vector used as an MPI buffer to send face-neighbor data in ExchangeFaceNbrData() to neighboring proce...
Definition pgridfunc.hpp:61
HypreParVector * ParallelProject() const
Returns a new vector restricted to the true dofs.
HypreParVector * ParallelAssemble() const
Returns a new vector assembled on the true dofs.
ParFiniteElementSpace * ParFESpace() const
void AddDistribute(real_t a, const Vector *tv)
void MakeRef(FiniteElementSpace *f, real_t *v) override
Make the ParGridFunction reference external data on a new FiniteElementSpace.
void GetVectorValue(int i, const IntegrationPoint &ip, Vector &val) const override
void SaveAsSerial(const char *fname, int precision=16, int save_rank=0) const
Vector face_nbr_data
Vector used to store data from face-neighbor processors, initialized by ExchangeFaceNbrData().
Definition pgridfunc.hpp:56
void ProjectCoefficient(Coefficient &coeff) override
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
void SaveAsOne(const char *fname, int precision=16) const
GridFunction GetSerialGridFunction(int save_rank, Mesh &serial_mesh) const
Returns a GridFunction on MPI rank save_rank that does not have any duplication of vertices/nodes at ...
void ParallelProject(Vector &tv) const
Returns the vector restricted to the true dofs.
void Update() override
Transform by the Space UpdateMatrix (e.g., on Mesh change).
Definition pgridfunc.cpp:91
void SetSpace(FiniteElementSpace *f) override
Associate a new FiniteElementSpace with the ParGridFunction.
Definition pgridfunc.cpp:97
void Distribute(const Vector *tv)
std::unique_ptr< ParGridFunction > ProlongateToMaxOrder() const
Return a GridFunction with the values of this, prolongated to the maximum order of all elements in th...
void GetElementDofValues(int el, Vector &dof_vals) const override
void ProjectBdrCoefficientTangent(VectorCoefficient &vcoeff, const Array< int > &bdr_attr) override
Project the tangential components of the given VectorCoefficient on the boundary. Only boundary attri...
Class for parallel linear form.
Class for parallel meshes.
Definition pmesh.hpp:34
Mesh GetSerialMesh(int save_rank) const
Definition pmesh.cpp:5291
ElementTransformation * GetFaceNbrElementTransformation(int FaceNo)
Returns a pointer to the transformation defining the i-th face neighbor.
Definition pmesh.cpp:3087
MPI_Comm GetComm() const
Definition pmesh.hpp:402
int GetMyRank() const
Definition pmesh.hpp:404
int GetNSharedFaces() const
Return the number of shared faces (3D), edges (2D), vertices (1D)
Definition pmesh.cpp:3152
int GetNRanks() const
Definition pmesh.hpp:403
FaceElementTransformations * GetFaceElementTransformations(int FaceNo, int mask=31) override
Definition pmesh.cpp:2896
int GetSharedFace(int sface) const
Return the local face index for the given shared face.
Definition pmesh.cpp:3171
int GetNFaceNeighbors() const
Definition pmesh.hpp:573
real_t GetFaceNbrElementSize(int i, int type=0)
Definition pmesh.cpp:3147
int GetFaceNbrRank(int fn) const
Definition pmesh.cpp:2795
FaceElementTransformations * GetSharedFaceTransformations(int sf, bool fill2=true)
Get the FaceElementTransformations for the given shared face (edge 2D) using the shared face index sf...
Definition pmesh.cpp:2922
void Mult(const Vector &x, Vector &y) const override
Matrix vector multiplication.
int Size_of_connections() const
Definition table.hpp:106
Memory< int > & GetJMemory()
Definition table.hpp:127
int * GetI()
Definition table.hpp:121
Base class for vector Coefficients that optionally depend on time and space.
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Definition lininteg.hpp:344
Vector coefficient defined by a vector GridFunction.
Vector data type.
Definition vector.hpp:82
virtual const real_t * HostRead() const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), false).
Definition vector.hpp:498
virtual const real_t * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Definition vector.hpp:494
void SetDataAndSize(real_t *d, int s)
Set the Vector data and size.
Definition vector.hpp:183
void SetSubVector(const Array< int > &dofs, const real_t value)
Set the entries listed in dofs to the given value.
Definition vector.cpp:679
Memory< real_t > data
Definition vector.hpp:85
void Destroy()
Destroy a vector.
Definition vector.hpp:635
int Size() const
Returns the size of the vector.
Definition vector.hpp:226
void SetSize(int s)
Resize the vector to size s.
Definition vector.hpp:558
virtual real_t * HostWrite()
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), false).
Definition vector.hpp:506
real_t * GetData() const
Return a pointer to the beginning of the Vector data.
Definition vector.hpp:235
virtual real_t * HostReadWrite()
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), false).
Definition vector.hpp:514
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Extract entries listed in dofs to the output Vector elemvect.
Definition vector.cpp:653
virtual real_t * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
Definition vector.hpp:502
real_t b
Definition lissajous.cpp:42
real_t a
Definition lissajous.cpp:41
const T * Read(const Memory< T > &mem, int size, bool on_dev=true)
Get a pointer for read access to mem with the mfem::Device's DeviceMemoryClass, if on_dev = true,...
Definition device.hpp:341
real_t infinity()
Define a shortcut for std::numeric_limits<double>::infinity()
Definition vector.hpp:47
real_t L2ZZErrorEstimator(BilinearFormIntegrator &flux_integrator, const ParGridFunction &x, ParFiniteElementSpace &smooth_flux_fes, ParFiniteElementSpace &flux_fes, Vector &errors, int norm_p, real_t solver_tol, int solver_max_it)
real_t GlobalLpNorm(const real_t p, real_t loc_norm, MPI_Comm comm)
Compute a global Lp norm from the local Lp norms computed by each processor.
real_t ComputeElementLpDistance(real_t p, int i, GridFunction &gf1, GridFunction &gf2)
Compute the Lp distance between two grid functions on the given element.
float real_t
Definition config.hpp:43
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:753
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition intrules.hpp:492
real_t p(const Vector &x, real_t t)
Helper struct to convert a C++ type to an MPI type.