28 template<
typename TDataType,
typename TParamVector,
typename TStateVector,
29 int residual_size,
int state_size,
int param_size>
37 void operator()(TParamVector &vparam, TStateVector &uu, TStateVector &rr)
39 MFEM_ASSERT(residual_size==4,
40 "PLaplacianResidual residual_size should be equal to 4!");
41 double pp = vparam[0];
42 double ee = vparam[1];
43 double ff = vparam[2];
50 TDataType norm2 = uu[0] * uu[0] + uu[1] * uu[1] + uu[2] * uu[2];
51 TDataType tvar = pow(ee * ee + norm2, (pp - 2.0) / 2.0);
70 template<
typename TDataType,
typename TParamVector,
typename TStateVector
71 ,
int state_size,
int param_size>
77 TDataType
operator()(TParamVector &vparam, TStateVector &uu)
79 MFEM_ASSERT(state_size==4,
"MyEnergyFunctor state_size should be equal to 4!");
80 MFEM_ASSERT(param_size==3,
"MyEnergyFunctor param_size should be equal to 3!");
81 double pp = vparam[0];
82 double ee = vparam[1];
83 double ff = vparam[2];
86 TDataType norm2 = uu[0] * uu[0] + uu[1] * uu[1] + uu[2] * uu[2];
88 TDataType rez = pow(ee * ee + norm2, pp / 2.0) / pp - ff *
u;
103 template<
class CQVectAutoDiff>
151 const int ndof = el.
GetDof();
152 const int ndim = el.
GetDim();
154 bool square = (ndim == spaceDim);
180 detJ = (square ? w : w * w);
195 vparam[0] =
pp->
Eval(trans, ip);
199 if (
coeff !=
nullptr)
209 for (
int jj = 0; jj < spaceDim; jj++)
211 uu[jj] = grad[jj] / detJ;
213 uu[3] = shapef * elfun;
215 energy = energy + w * qfunc(vparam,uu);
225 MFEM_PERF_BEGIN(
"AssembleElementVector");
226 const int ndof = el.
GetDof();
227 const int ndim = el.
GetDim();
258 for (
int jj = 0; jj < spaceDim; jj++)
267 vparam[0] =
pp->
Eval(trans, ip);
270 if (
coeff !=
nullptr)
283 rdf.VectorFunc(vparam,uu,du);
287 MFEM_PERF_END(
"AssembleElementVector");
295 MFEM_PERF_BEGIN(
"AssembleElementGrad");
296 const int ndof = el.
GetDof();
297 const int ndim = el.
GetDim();
329 for (
int jj = 0; jj < spaceDim; jj++)
338 vparam[0] =
pp->
Eval(trans, ip);
341 if (
coeff !=
nullptr)
354 rdf.Jacobian(vparam,uu,duu);
359 MFEM_PERF_END(
"AssembleElementGrad");
397 const int ndof = el.
GetDof();
398 const int ndim = el.
GetDim();
400 bool square = (ndim == spaceDim);
420 detJ = (square ? w : w * w);
431 nrgrad2 = grad * grad / (detJ * detJ);
436 ppp =
pp->
Eval(trans, ip);
440 if (
coeff !=
nullptr)
445 energy = energy + w * std::pow(nrgrad2 + eee * eee, ppp / 2.0) / ppp;
450 energy = energy - w * (shapef * elfun) *
load->
Eval(trans, ip);
461 MFEM_PERF_BEGIN(
"AssembleElementVector");
462 const int ndof = el.
GetDof();
463 const int ndim = el.
GetDim();
465 bool square = (ndim == spaceDim);
489 detJ = (square ? w : w * w);
501 nrgrad = grad.
Norml2() / detJ;
507 ppp =
pp->
Eval(trans, ip);
511 if (
coeff !=
nullptr)
516 aa = nrgrad * nrgrad + eee * eee;
517 aa = std::pow(aa, (ppp - 2.0) / 2.0);
518 dshape_xyz.
Mult(grad, lvec);
519 elvect.
Add(w * aa / (detJ * detJ), lvec);
527 MFEM_PERF_END(
"AssembleElementVector");
535 MFEM_PERF_BEGIN(
"AssembleElementGrad");
536 const int ndof = el.
GetDof();
537 const int ndim = el.
GetDim();
539 bool square = (ndim == spaceDim);
564 detJ = (square ? w : w * w);
577 ppp =
pp->
Eval(trans, ip);
580 if (
coeff !=
nullptr)
587 nrgrad = grad.
Norml2() / detJ;
589 aa0 = nrgrad * nrgrad + eee * eee;
590 aa1 = std::pow(aa0, (ppp - 2.0) / 2.0);
591 aa0 = (ppp - 2.0) * std::pow(aa0, (ppp - 4.0) / 2.0);
592 dshape_xyz.
Mult(grad, lvec);
593 w = w / (detJ * detJ);
598 MFEM_PERF_END(
"AssembleElementGrad");
608 template<
int sizeres=10>
637 const int ndof = el.
GetDof();
638 const int ndim = el.
GetDim();
640 bool square = (ndim == spaceDim);
660 detJ = (square ? w : w * w);
671 nrgrad2 = grad * grad / (detJ * detJ);
676 ppp =
pp->
Eval(trans, ip);
680 if (
coeff !=
nullptr)
685 energy = energy + w * std::pow(nrgrad2 + eee * eee, ppp / 2.0) / ppp;
690 energy = energy - w * (shapef * elfun) *
load->
Eval(trans, ip);
701 MFEM_PERF_BEGIN(
"AssembleElementVector");
702 const int ndof = el.
GetDof();
703 const int ndim = el.
GetDim();
705 bool square = (ndim == spaceDim);
729 detJ = (square ? w : w * w);
741 nrgrad = grad.
Norml2() / detJ;
747 ppp =
pp->
Eval(trans, ip);
751 if (
coeff !=
nullptr)
756 aa = nrgrad * nrgrad + eee * eee;
757 aa = std::pow(aa, (ppp - 2.0) / 2.0);
758 dshape_xyz.
Mult(grad, lvec);
759 elvect.
Add(w * aa / (detJ * detJ), lvec);
767 MFEM_PERF_END(
"AssembleElementVector");
775 MFEM_PERF_BEGIN(
"AssembleElementGrad");
776 const int ndof = el.
GetDof();
777 const int ndim = el.
GetDim();
779 bool square = (ndim == spaceDim);
801 vres.
SetSize(uu.Size()); vres=0.0;
814 detJ = (square ? w : w * w);
827 ppp =
pp->
Eval(trans, ip);
830 if (
coeff !=
nullptr)
837 for (
int i=0; i<spaceDim; i++)
839 for (
int j=0; j<ndof; j++)
841 grad[i]= grad[i]+ dshape_xyz(j,i)*uu[j];
846 nrgrad= (grad*grad)/(detJ*detJ);
848 aa = nrgrad + eee * eee;
849 aa = pow(aa, (ppp - 2.0) / 2.0);
851 for (
int i=0; i<spaceDim; i++)
853 for (
int j=0; j<ndof; j++)
855 lvec[j] = lvec[j] + dshape_xyz(j,i) * grad[i];
859 for (
int j=0; j<ndof; j++)
861 vres[j]=vres[j] + lvec[j] * (w*aa/(detJ*detJ));
870 MFEM_PERF_END(
"AssembleElementGrad");
int GetNPoints() const
Returns the number of the points in the integration rule.
Abstract class for all finite elements.
int GetDim() const
Returns the reference space dimension for the finite element.
void trans(const Vector &u, Vector &x)
Class for an integration rule - an Array of IntegrationPoint.
void SetCol(int c, const double *col)
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
virtual double GetElementEnergy(const FiniteElement &el, ElementTransformation &trans, const Vector &elfun)
Compute the local energy.
void SetSize(int s)
Resize the vector to size s.
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
pLaplaceSL(Coefficient &pp_)
double Norml2() const
Returns the l2 norm of the vector.
virtual void AssembleElementGrad(const FiniteElement &el, ElementTransformation &trans, const Vector &elfun, DenseMatrix &elmat)
Assemble the local gradient matrix.
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
virtual void AssembleElementGrad(const FiniteElement &el, ElementTransformation &trans, const Vector &elfun, DenseMatrix &elmat)
Assemble the local gradient matrix.
Data type dense matrix using column-major storage.
void AddMult_a_ABt(double a, const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &ABt)
ABt += a * A * B^t.
virtual double GetElementEnergy(const FiniteElement &el, ElementTransformation &trans, const Vector &elfun)
Compute the local energy.
virtual void AssembleElementVector(const FiniteElement &el, ElementTransformation &trans, const Vector &elfun, Vector &elvect)
Perform the local action of the NonlinearFormIntegrator.
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
virtual double GetElementEnergy(const FiniteElement &el, ElementTransformation &trans, const Vector &elfun)
Compute the local energy.
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
void MultTranspose(const double *x, double *y) const
Multiply a vector with the transpose matrix.
void operator()(TParamVector &vparam, TStateVector &uu, TStateVector &rr)
codi::RealForward ADFloatType
Forward AD type declaration.
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...
pLaplace(Coefficient &pp_, Coefficient &q, Coefficient &ld_)
void AddMult_a_VVt(const double a, const Vector &v, DenseMatrix &VVt)
VVt += a * v v^t.
void GetColumn(int c, Vector &col) const
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
virtual void AssembleElementVector(const FiniteElement &el, ElementTransformation &trans, const Vector &elfun, Vector &elvect)
Perform the local action of the NonlinearFormIntegrator.
TDataType operator()(TParamVector &vparam, TStateVector &uu)
virtual void AssembleElementGrad(const FiniteElement &el, ElementTransformation &trans, const Vector &elfun, DenseMatrix &elmat)
Assemble the local gradient matrix.
Vector & Add(const double a, const Vector &Va)
(*this) += a * Va
Class for integration point with weight.
pLaplaceAD(Coefficient &pp_)
pLaplace(Coefficient &pp_)
pLaplaceSL(Coefficient &pp_, Coefficient &q, Coefficient &ld_)
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the coefficient in the element described by T at the point ip.
Templated vector data type.
void Mult(const double *x, double *y) const
Matrix vector multiplication.
pLaplaceAD(Coefficient &pp_, Coefficient &q, Coefficient &ld_)
double u(const Vector &xvec)
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
void Jacobian(mfem::Vector &vparam, mfem::Vector &vstate, mfem::DenseMatrix &jac)
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const =0
Evaluate the gradients of all shape functions of a scalar finite element in reference space at the gi...
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
void AddMult_a_AAt(double a, const DenseMatrix &A, DenseMatrix &AAt)
AAt += a * A * A^t.
virtual void AssembleElementVector(const FiniteElement &el, ElementTransformation &trans, const Vector &elfun, Vector &elvect)
Perform the local action of the NonlinearFormIntegrator.