MFEM  v4.5.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
Public Member Functions | List of all members
mfem::QVectorFuncAutoDiff< TFunctor, vector_size, state_size, param_size > Class Template Reference

#include <admfem.hpp>

Public Member Functions

void VectorFunc (const mfem::Vector &vparam, mfem::Vector &uu, mfem::Vector &rr)
 
void Jacobian (mfem::Vector &vparam, mfem::Vector &uu, mfem::DenseMatrix &jac)
 
void VectorFunc (const Vector &vparam, Vector &uu, Vector &rr)
 
void Jacobian (mfem::Vector &vparam, mfem::Vector &uu, mfem::DenseMatrix &jac)
 

Detailed Description

template<template< typename, typename, typename, int, int, int > class TFunctor, int vector_size = 1, int state_size = 1, int param_size = 0>
class mfem::QVectorFuncAutoDiff< TFunctor, vector_size, state_size, param_size >

The class provides an evaluation of the Jacobian of a templated vector function provided as a functor TFunctor. The Jacobian is evaluated with the help of automatic differentiation (AD). The template parameters specify the size of the return vector (vector_size), the size of the input vector (state_size), and the size of the parameters supplied to the function. The TFunctor functor is a template class with parameters [Float data type], [Vector type for the additional parameters], [Vector type for the state vector and the return residual]. The integer template parameters are the same ones passed to QVectorFuncAutoDiff.

The class provides an evaluation of the Jacobian of a templated vector
function provided as a functor TFunctor. The Jacobian is evaluated with the
help of automatic differentiation (AD). The template parameters specify the
size of the return vector (vector_size), the size of the input vector
(state_size), and the size of the parameters supplied to the function. The
TFunctor functor is a template class with parameters [Float data type],
[Vector type for the additional parameters], [Vector type for the state
vector and the return residual].
The integer template parameters are the same ones
passed to QVectorFuncAutoDiff. \n
Example: f={sin(a*x*y), cos(b*x*y*z), x*x+y*x} \n
The vector function has vector_size=3, and state_size=3, i.e., it has
three arguments [x,y,z]. The parameters [a,b] size is 2.
The functor class will have the following form
template<typename TDataType, typename TParamVector, typename TStateVector,
int residual_size, int state_size, int param_size>
class MyVectorFunction{
public:
TDataType operator() (TParamVector& vparam, TStateVector& uu, TStateVector& rr)
{
auto a=vparam[0];
auto b=vparam[1];
rr[0]=sin(a*uu[0]*uu[1]);
rr[1]=cos(b*uu[0]*uu[1]*uu[2]);
rr[2]=uu[0]*uu[0]+uu[0]*uu[1];
}
//
};

Definition at line 146 of file admfem.hpp.

Member Function Documentation

template<template< typename, typename, typename, int, int, int > class TFunctor, int vector_size = 1, int state_size = 1, int param_size = 0>
void mfem::QVectorFuncAutoDiff< TFunctor, vector_size, state_size, param_size >::Jacobian ( mfem::Vector vparam,
mfem::Vector uu,
mfem::DenseMatrix jac 
)
inline

Returns the gradient of TFunctor(...) in the dense matrix jac. The dimensions of jac are vector_size x state_size, where state_size is the length of vector uu.

Definition at line 159 of file admfem.hpp.

template<template< typename, typename, typename, int, int, int > class TFunctor, int vector_size = 1, int state_size = 1, int param_size = 0>
void mfem::QVectorFuncAutoDiff< TFunctor, vector_size, state_size, param_size >::Jacobian ( mfem::Vector vparam,
mfem::Vector uu,
mfem::DenseMatrix jac 
)
inline

Returns the gradient of TFunctor(...) residual in the dense matrix jac. The dimensions of jac are vector_size x state_size, where state_size is the length of vector uu.

Definition at line 553 of file admfem.hpp.

template<template< typename, typename, typename, int, int, int > class TFunctor, int vector_size = 1, int state_size = 1, int param_size = 0>
void mfem::QVectorFuncAutoDiff< TFunctor, vector_size, state_size, param_size >::VectorFunc ( const mfem::Vector vparam,
mfem::Vector uu,
mfem::Vector rr 
)
inline

Evaluates the vector function for given set of parameters and state values in vector uu. The result is returned in vector rr.

Definition at line 151 of file admfem.hpp.

template<template< typename, typename, typename, int, int, int > class TFunctor, int vector_size = 1, int state_size = 1, int param_size = 0>
void mfem::QVectorFuncAutoDiff< TFunctor, vector_size, state_size, param_size >::VectorFunc ( const Vector vparam,
Vector uu,
Vector rr 
)
inline

Returns a vector valued function rr for supplied passive arguments vparam and active arguments uu. The evaluation is based on the user supplied TFunctor template class.

Definition at line 545 of file admfem.hpp.


The documentation for this class was generated from the following file: