specialisation of NiHu::double_integral for the collocational formalism
More...
#include <double_integral.hpp>
|
template<class dual_iterator_t > |
static result_t & | eval_on_accelerator (result_t &result, kernel_base< Kernel > const &kernel, field_base< TestField > const &test_field, field_base< TrialField > const &trial_field, dual_iterator_t first, dual_iterator_t last) |
| evaluate regular collocational integral with selected trial field accelerator More...
|
|
static result_t & | eval (WITHOUT_SINGULARITY_CHECK, result_t &result, kernel_base< Kernel > const &kernel, field_base< TestField > const &test_field, field_base< TrialField > const &trial_field) |
| evaluate single integral of a kernel on specific fields without singularity check More...
|
|
static result_t & | eval (WITH_SINGULARITY_CHECK, result_t &result, kernel_base< Kernel > const &kernel, field_base< TestField > const &test_field, field_base< TrialField > const &trial_field) |
| evaluate single integral of a kernel on specific fields with singularity check More...
|
|
template<class Kernel, class TestField, class TrialField>
class NiHu::double_integral< Kernel, TestField, TrialField, formalism::collocational >
specialisation of NiHu::double_integral for the collocational formalism
- Template Parameters
-
Kernel | type of the kernel to integrate |
TestField | type of the test field |
TrialField | type of the trial field |
Definition at line 395 of file double_integral.hpp.
◆ anonymous enum
template<class Kernel , class TestField , class TrialField >
Enumerator |
---|
kernel_rows | number of rows of the kernel result
|
kernel_cols | number of columns of the kernel result
|
Definition at line 427 of file double_integral.hpp.
◆ eval() [1/3]
template<class Kernel , class TestField , class TrialField >
template<class OnSameMesh >
evaluate collocational integral on given fields
- Template Parameters
-
OnSameMesh | indicates that the two fields are defined over the same mesh and singularity check may be needed |
- Parameters
-
[in] | kernel | the kernel to integrate |
[in] | test_field | the test field defining the collocation points |
[in] | trial_field | the trial field to integrate on |
- Returns
- the integration result by value
Definition at line 660 of file double_integral.hpp.
◆ eval() [2/3]
template<class Kernel , class TestField , class TrialField >
evaluate single integral of a kernel on specific fields with singularity check
- Parameters
-
[out] | result | reference to the integration result matrix |
[in] | kernel | the kernel to integrate |
[in] | test_field | reference to the test field defining the collocation points |
[in] | trial_field | reference to the trial field to integrate on |
- Returns
- reference to the stored result
Definition at line 608 of file double_integral.hpp.
◆ eval() [3/3]
template<class Kernel , class TestField , class TrialField >
evaluate single integral of a kernel on specific fields without singularity check
- Parameters
-
[out] | result | reference to the integration result matrix |
[in] | kernel | the kernel to integrate |
[in] | test_field | reference to the test field defining the collocation points |
[in] | trial_field | reference to the trial field to integrate on |
- Returns
- reference to the stored result
Definition at line 552 of file double_integral.hpp.
◆ eval_on_accelerator()
template<class Kernel , class TestField , class TrialField >
template<class dual_iterator_t >
evaluate regular collocational integral with selected trial field accelerator
- Template Parameters
-
- Parameters
-
[out] | result | reference to the integration result matrix |
[in] | kernel | the kernel to integrate |
[in] | test_field | defining the collocation points |
[in] | trial_field | the trial field to integrate on |
[in] | first | the begin iterator of the accelerator |
[in] | last | the end iterator of the accelerator |
- Returns
- reference to the integration result
Definition at line 446 of file double_integral.hpp.
The documentation for this class was generated from the following file: