|
| type(mformspace) function, dimension(0:3) | derhamsequence3c_1d (bspline_x, domain) |
| | Create a 1D DeRham sequence of m-forms from a tensor product B-spline basis with 3 components for vector-valued forms.
|
| |
| type(mformspace) function, dimension(0:3) | derhamsequence3c_2d (bspline_x, bspline_y, domain) |
| | Create a 2D DeRham sequence of m-forms from a tensor product B-spline basis with 3 components for vector-valued forms.
|
| |
| type(mformspace) function, dimension(0:3) | derhamsequence_3d (bspline_x, bspline_y, bspline_z, domain) |
| | Create a DeRham sequence of m-forms from a tensor product B-spline basis.
|
| |
◆ derhamsequence3c_1d()
| type(mformspace) function, dimension(0:3) m_mform_derham::derhamsequence3c::derhamsequence3c_1d |
( |
type(bsplinespace), intent(in) | bspline_x, |
|
|
type(domaindecomp), intent(in), optional | domain ) |
Create a 1D DeRham sequence of m-forms from a tensor product B-spline basis with 3 components for vector-valued forms.
- Parameters
-
| [in] | bspline_x | B-spline basis in x direction |
| [in] | domain | _(optional)_ Domain decomposition |
- Returns
- DeRham sequence (space of 0-forms, 1-forms, 2-forms, 3-forms)
◆ derhamsequence3c_2d()
| type(mformspace) function, dimension(0:3) m_mform_derham::derhamsequence3c::derhamsequence3c_2d |
( |
type(bsplinespace), intent(in) | bspline_x, |
|
|
type(bsplinespace), intent(in) | bspline_y, |
|
|
type(domaindecomp), intent(in), optional | domain ) |
Create a 2D DeRham sequence of m-forms from a tensor product B-spline basis with 3 components for vector-valued forms.
- Parameters
-
| [in] | bspline_x | B-spline basis in x direction |
| [in] | bspline_y | B-spline basis in y direction |
| [in] | domain | _(optional)_ Domain decomposition |
- Returns
- DeRham sequence (space of 0-forms, 1-forms, 2-forms, 3-forms)
◆ derhamsequence_3d()
| type(mformspace) function, dimension(0:3) m_mform_derham::derhamsequence3c::derhamsequence_3d |
( |
type(bsplinespace), intent(in) | bspline_x, |
|
|
type(bsplinespace), intent(in) | bspline_y, |
|
|
type(bsplinespace), intent(in) | bspline_z, |
|
|
type(domaindecomp), intent(in), optional | domain ) |
Create a DeRham sequence of m-forms from a tensor product B-spline basis.
- Parameters
-
| [in] | bspline_x | B-spline basis in x direction |
| [in] | bspline_y | B-spline basis in y direction |
| [in] | bspline_z | B-spline basis in z direction |
| [in] | domain | _(optional)_ Domain decomposition |
- Returns
- DeRham sequence (space of 0-forms, 1-forms, 2-forms, 3-forms)
The documentation for this interface was generated from the following file:
- /builds/rwr/bspline_feec/src/mform/m_mform_derham.f90