12 #ifndef MFEM_COMPLEX_FEM
13 #define MFEM_COMPLEX_FEM
15 #include "../linalg/complex_operator.hpp"
38 void Destroy() {
delete gfr;
delete gfi; }
49 { *gfr = value.real(); *gfi = value.imag();
return *
this; }
125 convention) { conv = convention; }
241 convention) { conv = convention; }
316 int copy_interior = 0);
361 { *pgfr = value.real(); *pgfi = value.imag();
return *
this; }
409 return sqrt(err_r * err_r + err_i * err_i);
419 return sqrt(err_r * err_r + err_i * err_i);
471 convention) { conv = convention; }
589 convention) { conv = convention; }
670 int copy_interior = 0);
685 #endif // MFEM_USE_MPI
689 #endif // MFEM_COMPLEX_FEM
AssemblyLevel
Enumeration defining the assembly level for bilinear and nonlinear form classes derived from Operator...
virtual void ProjectBdrCoefficientTangent(VectorCoefficient &real_coeff, VectorCoefficient &imag_coeff, Array< int > &attr)
Class for an integration rule - an Array of IntegrationPoint.
Class for grid function - Vector with associated FE space.
const GridFunction & imag() const
Base class for vector Coefficients that optionally depend on time and space.
virtual void ProjectBdrCoefficient(Coefficient &real_coeff, Coefficient &imag_coeff, Array< int > &attr)
virtual void ProjectCoefficient(Coefficient &real_coeff, Coefficient &imag_coeff)
Pointer to an Operator of a specified type.
FiniteElementSpace * FESpace()
virtual double ComputeL2Error(Coefficient &exsolr, Coefficient &exsoli, const IntegrationRule *irs[]=NULL) const
Abstract parallel finite element space.
const GridFunction & real() const
void SyncAliasMemory(const Vector &v) const
Update the alias memory location of the vector to match v.
const ParFiniteElementSpace * ParFESpace() const
virtual void ProjectBdrCoefficientNormal(VectorCoefficient &real_coeff, VectorCoefficient &imag_coeff, Array< int > &attr)
ParFiniteElementSpace * ParFESpace()
ParComplexGridFunction & operator=(const std::complex< double > &value)
Assign constant values to the ParComplexGridFunction data.
double f(const Vector &xvec)
void Distribute(const Vector &tv)
ComplexGridFunction(FiniteElementSpace *f)
Construct a ComplexGridFunction associated with the FiniteElementSpace *f.
virtual void ProjectBdrCoefficient(Coefficient &real_coeff, Coefficient &imag_coeff, Array< int > &attr)
virtual double ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const
virtual double ComputeL2Error(VectorCoefficient &exsolr, VectorCoefficient &exsoli, const IntegrationRule *irs[]=NULL, Array< int > *elems=NULL) const
virtual ~ComplexGridFunction()
Destroys the grid function.
Set the diagonal value to one.
FiniteElementSpace * FESpace()
void SyncMemory(const Vector &v) const
Update the memory location of the vector to match v.
Wrapper for hypre's parallel vector class.
void Distribute(const Vector *tv)
virtual void ProjectBdrCoefficientNormal(VectorCoefficient &real_coeff, VectorCoefficient &imag_coeff, Array< int > &attr)
Specialization of the ComplexOperator built from a pair of Sparse Matrices.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
virtual ~ParComplexGridFunction()
Destroys grid function.
ParComplexGridFunction(ParFiniteElementSpace *pf)
Construct a ParComplexGridFunction associated with the ParFiniteElementSpace *pf. ...
FiniteElementSpace * FESpace()
virtual void ProjectBdrCoefficientTangent(VectorCoefficient &real_coeff, VectorCoefficient &imag_coeff, Array< int > &attr)
const ParGridFunction & imag() const
DiagonalPolicy
Defines operator diagonal policy upon elimination of rows and/or columns.
Specialization of the ComplexOperator built from a pair of HypreParMatrices.
Native convention for Hermitian operators.
const FiniteElementSpace * FESpace() const
const FiniteElementSpace * FESpace() const
const ParGridFunction & real() const
ComplexGridFunction & operator=(const std::complex< double > &value)
Assign constant values to the ComplexGridFunction data.
Class for parallel grid function.
void ParallelProject(Vector &tv) const
Returns the vector restricted to the true dofs.
virtual void ProjectCoefficient(Coefficient &real_coeff, Coefficient &imag_coeff)
ParFiniteElementSpace * ParFESpace() const