38 const int *partitioning)
55 int element_counter = 0;
57 const int glob_ne = glob_fes->
GetNE();
58 for (
int i = 0; i < glob_ne; i++)
60 if (partitioning[i] == MyRank)
102 MFEM_ASSERT(
pfes != NULL,
"not a ParFiniteElementSpace");
117 MFEM_ASSERT(
pfes != NULL,
"not a ParFiniteElementSpace");
132 MFEM_ASSERT(
pfes != NULL,
"not a ParFiniteElementSpace");
145 prolong->
Mult(*tv, *
this);
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];
241 auto d_data = this->
Read();
245 const int ldof = d_send_ldof[i];
246 d_send_data[i] = d_data[ldof >= 0 ? ldof : -1-ldof];
254 if (mpi_gpu_aware) { MFEM_STREAM_SYNC; }
255 for (
int fn = 0; fn < num_face_nbrs; fn++)
260 MPI_Isend(&send_data_ptr[send_offset[fn]],
261 send_offset[fn+1] - send_offset[fn],
264 MPI_Irecv(&face_nbr_data_ptr[recv_offset[fn]],
265 recv_offset[fn+1] - recv_offset[fn],
269 MPI_Waitall(num_face_nbrs, send_requests, statuses);
270 MPI_Waitall(num_face_nbrs, recv_requests, statuses);
294 int s = dofs.
Size()/fes_vdim;
345 return (DofVal * LocVec);
365 int dof = FElem->
GetDof();
382 for (
int k = 0; k < vdim; k++)
384 val(k) = shape * (&loc_data[dof * k]);
407 int comp,
Vector *tr)
const
451 return (DofVal * LocVec);
488 const int dof = fe->
GetDof();
502 for (
int k = 0; k < vdim; k++)
504 val(k) = shape * (&loc_data[dof * k]);
524 gcomm.
Bcast(elem_per_vdof);
536 gcomm.
Bcast(overlap);
542 for (
int i = 0; i < overlap.
Size(); i++)
544 der(i) /= overlap[i];
554 "ParGridFunction::GetElementDofValues: ExchangeFaceNbrData "
555 "must be called before accessing face neighbor elements.");
577 real_t loc_integral, glob_integral;
585 (*this) *= (delta_c->
Scale() / glob_integral);
599 ldof_attr.
Copy(gdof_attr);
602 gcomm.
Bcast(gdof_attr);
608 if (gdof_attr[i] > ldof_attr[i])
622 gcomm.
Bcast(gdof_attr);
625 (*this)(i) /= gdof_attr[i];
645 gcomm.
Bcast(zones_per_vdof);
667 gcomm.
Bcast(zones_per_vdof);
683 for (
int i = 0; i < values.
Size(); i++)
685 values(i) = values_counter[i] ? (*this)(i) : 0.0;
694 for (
int i = 0; i < values.
Size(); i++)
696 if (values_counter[i])
698 (*this)(i) = values(i)/values_counter[i];
709 ess_vdofs_marker.SetSize(
Size());
710 ess_vdofs_marker = 0;
713 if (!coeff[i]) {
continue; }
716 for (
int j = 0; j<
Size(); j++)
718 ess_vdofs_marker[j] = bool(ess_vdofs_marker[j]) ||
719 bool(component_dof_marker[j]);
724 for (
int i = 0; i < values_counter.
Size(); i++)
726 MFEM_ASSERT(
bool(values_counter[i]) ==
bool(ess_vdofs_marker[i]),
739 for (
int i = 0; i < values.
Size(); i++)
741 values(i) = values_counter[i] ? (*this)(i) : 0.0;
750 for (
int i = 0; i < values.
Size(); i++)
752 if (values_counter[i])
754 (*this)(i) = values(i)/values_counter[i];
764 for (
int i = 0; i < values_counter.
Size(); i++)
766 MFEM_ASSERT(
bool(values_counter[i]) ==
bool(ess_vdofs_marker[i]),
780 int fdof, intorder, k;
782 Vector shape, el_dofs, err_val, ell_coeff_val;
789 std::map<int,int> local_to_shared;
793 local_to_shared[i_local] = i;
798 real_t shared_face_factor = 1.0;
799 bool shared_face =
false;
800 int iel1, iel2, info1, info2;
809 if (info2 >= 0 && iel2 < 0)
811 int ishared = local_to_shared[i];
815 if ( (k = fe2->
GetOrder()) > intorder )
820 shared_face_factor = 0.5;
828 if ( (k = fe2->
GetOrder()) > intorder )
842 intorder = 2 * intorder;
855 transf = face_elem_transf->
Elem1;
861 for (k = 0; k < fdof; k++)
864 el_dofs(k) = (*this)(vdofs[k]);
868 el_dofs(k) = - (*this)(-1-vdofs[k]);
875 ell_coeff_val(j) = ell_coeff->
Eval(*transf, eip);
876 err_val(j) = exsol->
Eval(*transf, eip) - (shape * el_dofs);
881 transf = face_elem_transf->
Elem2;
888 for (k = 0; k < fdof; k++)
901 for (k = 0; k < fdof; k++)
904 el_dofs(k) = (*this)(vdofs[k]);
908 el_dofs(k) = - (*this)(-1 - vdofs[k]);
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));
922 transf = face_elem_transf;
928 face_error += shared_face_factor*(ip.
weight * nu * ell_coeff_val(j) *
930 err_val(j) * err_val(j));
933 error += fabs(face_error);
943 for (
int i = 0; i <
size; i++)
950 for (
int i = 0; i <
size; i++)
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);
973 ofs.precision(precision);
987 serialgf.
Save(fname, precision);
1004 const int vdim = pfespace->
GetVDim();
1008 MPI_Comm comm = pmesh->
GetComm();
1010 if (my_rank == save_rank)
1015 for (
int e = 0; e < pmesh->
GetNE(); e++)
1022 for (
int p = 0;
p < nranks;
p++)
1024 if (
p == save_rank) {
continue; }
1026 MPI_Recv(&n_send_recv, 1, MPI_INT,
p, 448, comm, &status);
1033 for (
int i = 0; i < n_send_recv; )
1043 int n_send_recv = 0;
1045 for (
int e = 0; e < pmesh->
GetNE(); e++)
1048 n_send_recv += vdim*fe->
GetDof();
1050 MPI_Send(&n_send_recv, 1, MPI_INT, save_rank, 448, comm);
1053 for (
int e = 0; e < pmesh->
GetNE(); e++)
1056 for (
int j = 0; j < nodeval.
Size(); j++)
1072 Mesh &serial_mesh)
const
1084#ifdef MFEM_USE_ADIOS2
1086 const std::string& variable_name,
1090 for (
int i = 0; i <
size; i++)
1097 for (
int i = 0; i <
size; i++)
1112 MyComm =
pfes -> GetComm();
1114 MPI_Comm_size(MyComm, &NRanks);
1115 MPI_Comm_rank(MyComm, &MyRank);
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];
1127 nv[0] =
pfes -> GetVSize();
1128 nvdofs[0] =
pfes -> GetNVDofs();
1129 nedofs[0] =
pfes -> GetNEDofs();
1130 nfdofs[0] =
pfes -> GetNFDofs();
1137 for (
p = 1;
p < NRanks;
p++)
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);
1148 int vdim =
pfes -> GetVDim();
1150 for (
p = 0;
p < NRanks;
p++)
1152 nrdofs[
p] = nv[
p]/vdim - nvdofs[
p] - nedofs[
p] - nfdofs[
p];
1157 for (
int d = 0; d < vdim; d++)
1159 for (
p = 0;
p < NRanks;
p++)
1160 for (i = 0; i < nvdofs[
p]; i++)
1162 os << *values[
p]++ <<
'\n';
1165 for (
p = 0;
p < NRanks;
p++)
1166 for (i = 0; i < nedofs[
p]; i++)
1168 os << *values[
p]++ <<
'\n';
1171 for (
p = 0;
p < NRanks;
p++)
1172 for (i = 0; i < nfdofs[
p]; i++)
1174 os << *values[
p]++ <<
'\n';
1177 for (
p = 0;
p < NRanks;
p++)
1178 for (i = 0; i < nrdofs[
p]; i++)
1180 os << *values[
p]++ <<
'\n';
1186 for (
p = 0;
p < NRanks;
p++)
1187 for (i = 0; i < nvdofs[
p]; i++)
1188 for (
int d = 0; d < vdim; d++)
1190 os << *values[
p]++ <<
'\n';
1193 for (
p = 0;
p < NRanks;
p++)
1194 for (i = 0; i < nedofs[
p]; i++)
1195 for (
int d = 0; d < vdim; d++)
1197 os << *values[
p]++ <<
'\n';
1200 for (
p = 0;
p < NRanks;
p++)
1201 for (i = 0; i < nfdofs[
p]; i++)
1202 for (
int d = 0; d < vdim; d++)
1204 os << *values[
p]++ <<
'\n';
1207 for (
p = 0;
p < NRanks;
p++)
1208 for (i = 0; i < nrdofs[
p]; i++)
1209 for (
int d = 0; d < vdim; d++)
1211 os << *values[
p]++ <<
'\n';
1215 for (
p = 1;
p < NRanks;
p++)
1218 delete [] values[
p];
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);
1244 loc_norm = fabs(loc_norm);
1248 loc_norm = pow(loc_norm,
p);
1253 glob_norm = pow(fabs(glob_norm), 1.0/
p);
1270 MFEM_VERIFY(ffes,
"the flux FE space must be ParFiniteElementSpace");
1283 for (
int i = 0; i < count.
Size(); i++)
1285 if (count[i] != 0) { flux(i) /= count[i]; }
1318 P.
Mult(*
this, *xMax);
1321 return std::unique_ptr<ParGridFunction>(xMax);
1329 int norm_p,
real_t solver_tol,
int solver_max_it)
1339 for (
int i = 0; i < xfes->
GetNE(); i++)
1350 *flux_fes.
GetFE(i), el_f,
false);
1371 a->AddDomainIntegrator(vmass);
1419 real_t total_error = 0.0;
1421 for (
int i = 0; i < xfes->
GetNE(); i++)
1424 total_error += pow(errors(i), norm_p);
1432 return pow(glob_error, 1.0/norm_p);
void Reserve(int capacity)
Ensures that the allocated size is at least the given size.
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
int Size() const
Return the logical size of the array.
int Append(const T &el)
Append element 'el' to array, resize if necessary.
void Copy(Array ©) const
Create a copy of the internal array to the provided copy.
T * HostReadWrite()
Shortcut for mfem::ReadWrite(a.GetMemory(), a.Size(), false).
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.
void MultTranspose(const real_t *x, real_t *y) const
Multiply a vector with the transpose matrix.
static bool GetGPUAwareMPI()
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
static FiniteElementCollection * New(const char *name)
Factory method: return a newly allocated FiniteElementCollection according to the given name.
int GetOrder() const
Return the order (polynomial degree) of the FE collection, corresponding to the order/degree returned...
virtual FiniteElementCollection * Clone(int p) const
Instantiate a new collection of the same type with a different order.
virtual const char * Name() const
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
bool IsVariableOrder() const
Returns true if the space contains elements of varying polynomial orders.
void DofsToVDofs(Array< int > &dofs, int ndofs=-1) const
Compute the full set of vdofs corresponding to each entry in dofs.
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...
ElementTransformation * GetElementTransformation(int i) const
Returns ElementTransformation for the i-th element.
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 ...
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
Ordering::Type GetOrdering() const
Return the ordering method.
int GetNE() const
Returns number of elements in the mesh.
const FiniteElementCollection * FEColl() const
Mesh * GetMesh() const
Returns the mesh.
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
int GetVDim() const
Returns the vector dimension of the finite element space.
const FiniteElement * GetTypicalFE() const
Return GetFE(0) if the local mesh is not empty; otherwise return a typical FE based on the Geometry t...
Abstract class for all finite elements.
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...
int GetRangeDim() const
Returns the vector dimension for vector-valued finite elements, which is also the dimension of the in...
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
int GetMapType() const
Returns the FiniteElement::MapType of the element describing how reference functions are mapped to ph...
int GetRangeType() const
Returns the FiniteElement::RangeType of the element, one of {SCALAR, VECTOR}.
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.
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 ...
Class for grid function - Vector with associated FE space.
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
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).
virtual void MakeRef(FiniteElementSpace *f, real_t *v)
Make the GridFunction reference external data on a new FiniteElementSpace.
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.
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.
FiniteElementCollection * fec_owned
Used when the grid function is read from a file. It can also be set explicitly, see MakeOwner().
void SumFluxAndCount(BilinearFormIntegrator &blfi, GridFunction &flux, Array< int > &counts, bool wcoef, int subdomain)
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
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.
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.
void SetPrintLevel(int print_level)
void Mult(const HypreParVector &b, HypreParVector &x) const override
Solve Ax=b with hypre's PCG.
void SetPrintLevel(int print_lvl)
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
void SetMaxIter(int max_iter)
Wrapper for hypre's ParCSR matrix class.
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.
Wrapper for hypre's parallel vector class.
Class for integration point with weight.
Class for an integration rule - an Array of IntegrationPoint.
int GetNPoints() const
Returns the number of the points in the integration rule.
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
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
void GetFaceInfos(int Face, int *Inf1, int *Inf2) const
int GetNumFaces() const
Return the number of faces (3D), edges (2D) or vertices (1D).
int GetNE() const
Returns number of elements.
real_t GetElementSize(int i, int type=0)
Get the size of the i-th element relative to the perfect reference element.
void GetFaceElements(int Face, int *Elem1, int *Elem2) const
int SpaceDimension() const
Dimension of the physical space containing the mesh.
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 ...
Matrix-free transfer operator between finite element spaces on the same mesh.
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.
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.
void ExchangeFaceNbrData()
HypreParVector * NewTrueDofVector()
const FiniteElement * GetFaceNbrFE(int i, int ndofs=0) const
GroupCommunicator & GroupComm()
Return a reference to the internal GroupCommunicator (on VDofs)
bool Nonconforming() const
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.
void GetFaceNbrElementVDofs(int i, Array< int > &vdofs, DofTransformation &doftrans) const
int GetFaceNbrVSize() const
const SparseMatrix * GetRestrictionMatrix() const override
Get the R matrix which restricts a local dof vector to true dof vector.
ParMesh * GetParMesh() const
const FiniteElement * GetFE(int i) const override
ElementTransformation * GetFaceNbrElementTransformation(int i) const
Class for parallel grid function.
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 ExchangeFaceNbrData()
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.
Vector send_data
Vector used as an MPI buffer to send face-neighbor data in ExchangeFaceNbrData() to neighboring proce...
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().
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).
void SetSpace(FiniteElementSpace *f) override
Associate a new FiniteElementSpace with the ParGridFunction.
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 meshes.
Mesh GetSerialMesh(int save_rank) const
ElementTransformation * GetFaceNbrElementTransformation(int FaceNo)
Returns a pointer to the transformation defining the i-th face neighbor.
int GetNSharedFaces() const
Return the number of shared faces (3D), edges (2D), vertices (1D)
FaceElementTransformations * GetFaceElementTransformations(int FaceNo, int mask=31) override
int GetSharedFace(int sface) const
Return the local face index for the given shared face.
int GetNFaceNeighbors() const
real_t GetFaceNbrElementSize(int i, int type=0)
int GetFaceNbrRank(int fn) const
FaceElementTransformations * GetSharedFaceTransformations(int sf, bool fill2=true)
Get the FaceElementTransformations for the given shared face (edge 2D) using the shared face index sf...
void Mult(const Vector &x, Vector &y) const override
Matrix vector multiplication.
int Size_of_connections() const
Memory< int > & GetJMemory()
Base class for vector Coefficients that optionally depend on time and space.
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Vector coefficient defined by a vector GridFunction.
virtual const real_t * HostRead() const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), false).
virtual const real_t * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
void SetDataAndSize(real_t *d, int s)
Set the Vector data and size.
void SetSubVector(const Array< int > &dofs, const real_t value)
Set the entries listed in dofs to the given value.
void Destroy()
Destroy a vector.
int Size() const
Returns the size of the vector.
void SetSize(int s)
Resize the vector to size s.
virtual real_t * HostWrite()
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), false).
real_t * GetData() const
Return a pointer to the beginning of the Vector data.
virtual real_t * HostReadWrite()
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), false).
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Extract entries listed in dofs to the output Vector elemvect.
virtual real_t * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
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,...
real_t infinity()
Define a shortcut for std::numeric_limits<double>::infinity()
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.
std::function< real_t(const Vector &)> f(real_t mass_coeff)
void forall(int N, lambda &&body)
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
real_t p(const Vector &x, real_t t)
Helper struct to convert a C++ type to an MPI type.