Skip to main content

ConvectiveMatrix example 2

!!! note "" This example shows how to USE the SUBROUTINE called ConvectiveMatrix to create a convective matrix in space domain for Line2 element.

Here, we want to DO the following.

M(I,J)=∫ΩNIckβˆ‚NJβˆ‚xkdΞ©M(I,J) = \int_{\Omega} N^I c_k \frac{\partial N^J}{\partial x_k} d{\Omega} M(I,J)=βˆ‚uiI∫Ωckβˆ‚NIβˆ‚xkNJdΞ©M(I,J) = \partial u_{iI} \int_{\Omega} c_k \frac{\partial N^I}{\partial x_k} N^J d{\Omega}

!!! warning "" cc is convective velocity, which can be a constant, or a FUNCTION of spatial coordinates, or some nonlinear FUNCTION.

In this example, convective matrix is formed for

  • [[ReferenceLine_]] element
  • [[QuadraturePoint_]] GaussLegendre
  • velocity is given by nodal values

Modules and classes​

  • [[ElemshapeData_]]
  • [[QuadraturePoint_]]
  • [[ReferenceLine_]]
  • [[FEVariable_]]

Usage​

PROGRAM main
USE easifemBase
IMPLICIT NONE
TYPE(ElemshapeData_) :: test, trial
TYPE(QuadraturePoint_) :: quad
TYPE(ReferenceLine_) :: refelem
TYPE(FEVariable_) :: cvar
REAL(DFP), PARAMETER :: c(1,2) = RESHAPE([1.0,1.0],[1,2])
REAL(DFP), ALLOCATABLE :: mat(:, :)
REAL(DFP), ALLOCATABLE :: XiJ(:, :)
INTEGER( I4B ), PARAMETER :: order = 1

!!! note "" Let us now create the physical coordinate of the line element.

    XiJ = RESHAPE([-1, 1], [1, 2])

!!! note "" Now we create an instance of [[ReferenceLine_]].

    refelem = referenceline(nsd=1)

!!! note "" Here, we create the quadrature points.

    CALL initiate( obj=quad, refelem=refelem, order=order+order-1, &
& quadratureType='GaussLegendre' )

!!! note "" Initiate an instance of [[ElemshapeData_]]. You can learn more about it from [[ElemshapeData_test]]

    CALL initiate(obj=test, &
& quad=quad, &
& refelem=refelem, &
& ContinuityType=typeH1, &
& InterpolType=typeLagrangeInterpolation)
CALL Set( obj=test, val=xij, N=test%N, dNdXi=test%dNdXi)

!!! note "" Let us now create an instance of [[FEVariable_]] to wrape the nodal values of cc. You can learn more about this at [[FEVariable_test]]

    cvar = NodalVariable(c, typeFEVariableVector, typeFEVariableSpace)

!!! note "" Let us now create the convective matrix.

    mat=ConvectiveMatrix(test=test, trial=test, term1=1, term2=0, c=cvar)
CALL Display(mat, "mat:")

??? example "Results"

        mat:
--------------------
-0.500000 -0.500000
0.500000 0.500000

!!! note "" Let us now create the convective matrix.

    mat=ConvectiveMatrix(test=test, trial=test, term1=0, term2=1, c=cvar)
CALL Display(mat, "mat:")

??? example "Results"

        mat:
-------------------
-0.500000 0.500000
-0.500000 0.500000

!!! settings "Cleanup"

END PROGRAM main