Choral
Data Types | Functions/Subroutines
geotsf_mod Module Reference

DERIVED TYPE geoTsf: geometrical transformation of reference cells More...

Data Types

interface  assemble
 1- Define the transformation T: K_ref –> K i.e. compute the ytansformation 'T' coefficients More...
 
interface  clear
 Destructor. More...
 
interface  geotsf
 Geometrical transformation \(T~: K_{\rm ref} \mapsto K \). More...
 

Functions/Subroutines

subroutine geotsf_clear (g)
 Destructor for geoTsf type. More...
 
type(geotsf) function geotsf_create (dim, nn, y)
  Constructor for the type geoTsf More...
 
subroutine geotsf_assemble (g, type, X, nNd)
 Define the transformation T: K_ref –> K. More...
 
subroutine, public compute_j (g)
 Compute Jy(i) at the nodes y(:;i) More...
 
subroutine, public compute_jc (g)
 Computes Jy(i) and Cy(:,:,i) at the nodes y(:;i) More...
 
subroutine, public compute_dtjc (g)
 Compute DTy(:,:,i), Jy(i) and Cy(:,:,i) at the nodes y(:,i) More...
 
subroutine, public compute_dtj (g)
 Compute DTy(:,:,i) and Jy(i) at the nodes y(:,i) More...
 
logical function, public belongstocell (x, Y, n, ct)
 Tests if the node x belongs to the cell with node coordinates Y and of type ct. More...
 
subroutine, public geotsf_t_inv (z, X, Y, n, ct)
 x = 3D node 3D cell of type ct with nodes = A, B, C, ... given by Y(:,1)=A, Y(:,2)=B, ..., Y(:,n) More...
 
subroutine, public cell_ms_itf_n (ms, n, itf_ms, nbItf, g)
 INPUT :: g = geoTsf, must have been set nbItf = number of interfaces OUTPUT : ms = measure of T(K_ref) n = unit normal to the interfaces of T(K_ref) itf_ms = interface measures. More...
 
subroutine compute_dt_cell_edg_2 (g)
 
subroutine compute_dt_cell_trg_2 (g)
 
subroutine compute_dt_cell_tet_2 (g)
 

Detailed Description

DERIVED TYPE geoTsf: geometrical transformation of reference cells

The derived type geoTsf defines geometrical transformations

\(~~~~~~~~~~ T~: K_{\rm ref} \mapsto K = T(K) ~~~~~~~~~~~~~~~~~~~~~~(1)\)

\(~~~~~~~~~~~~~~~~~~~~~~~ y \mapsto x = T(y) \)

A transformation \( T \) as in (1) is defined by

THE GOAL

Given a set \( y = \{y_1,\dots , y_n\}\) of \( n \) points in \( K_{\rm ref} \) (basically quadrature nodes),
we want to know at each point \( y_i \):

THE METHOD

To take into account both

we define an object 'g' of type geoTsf in two steps:

  1. g = geotsf(dim, n, y)
    • dim: the dimension of \( K_{\rm ref}\),
    • n: the number of points in \( y\)
    • y: the coordinates of the points in \( y\).
  2. call assemble(g, type, X, nNd)
    • type: the geometrical type of \( K \),
    • nNd: the number of nodes of \(K\),
    • X: the cell \(K\) node coordinates.

Which two steps compute:

Finally,

Complementary matrix definition:

\(~~~~~~~~~~ \nabla u(x) = C(y) \nabla u_{\rm ref}(y) / J(y) \)

Author
Charles Pierre

Function/Subroutine Documentation

◆ belongstocell()

logical function, public geotsf_mod::belongstocell ( real(rp), dimension(3), intent(in)  x,
real(rp), dimension(3,n), intent(inout)  Y,
integer, intent(in)  n,
integer, intent(in)  ct 
)

Tests if the node x belongs to the cell with node coordinates Y and of type ct.

WARNING : Y is modified (see geoTsf_T_Inv)

Definition at line 603 of file geoTsf_mod.F90.

Here is the call graph for this function:

◆ cell_ms_itf_n()

subroutine, public geotsf_mod::cell_ms_itf_n ( real(rp), intent(out)  ms,
real(rp), dimension(3,nbitf), intent(out)  n,
real(rp), dimension(nbitf), intent(out)  itf_ms,
integer, intent(in)  nbItf,
type(geotsf), intent(in)  g 
)

INPUT :: g = geoTsf, must have been set nbItf = number of interfaces OUTPUT : ms = measure of T(K_ref) n = unit normal to the interfaces of T(K_ref) itf_ms = interface measures.

Definition at line 749 of file geoTsf_mod.F90.

Here is the call graph for this function:

◆ compute_dt_cell_edg_2()

subroutine geotsf_mod::compute_dt_cell_edg_2 ( type(geotsf), intent(inout)  g)
private

Definition at line 806 of file geoTsf_mod.F90.

◆ compute_dt_cell_tet_2()

subroutine geotsf_mod::compute_dt_cell_tet_2 ( type(geotsf), intent(inout)  g)
private

Definition at line 833 of file geoTsf_mod.F90.

◆ compute_dt_cell_trg_2()

subroutine geotsf_mod::compute_dt_cell_trg_2 ( type(geotsf), intent(inout)  g)
private

Definition at line 818 of file geoTsf_mod.F90.

◆ compute_dtj()

subroutine, public geotsf_mod::compute_dtj ( type(geotsf), intent(inout)  g)

Compute DTy(:,:,i) and Jy(i) at the nodes y(:,i)

Definition at line 539 of file geoTsf_mod.F90.

Here is the call graph for this function:

◆ compute_dtjc()

subroutine, public geotsf_mod::compute_dtjc ( type(geotsf), intent(inout)  g)

Compute DTy(:,:,i), Jy(i) and Cy(:,:,i) at the nodes y(:,i)

Definition at line 490 of file geoTsf_mod.F90.

Here is the call graph for this function:

◆ compute_j()

subroutine, public geotsf_mod::compute_j ( type(geotsf), intent(inout)  g)

Compute Jy(i) at the nodes y(:;i)

Definition at line 383 of file geoTsf_mod.F90.

Here is the call graph for this function:

◆ compute_jc()

subroutine, public geotsf_mod::compute_jc ( type(geotsf), intent(inout)  g)

Computes Jy(i) and Cy(:,:,i) at the nodes y(:;i)

Definition at line 428 of file geoTsf_mod.F90.

Here is the call graph for this function:

◆ geotsf_assemble()

subroutine geotsf_mod::geotsf_assemble ( type(geotsf), intent(inout)  g,
integer, intent(in)  type,
real(rp), dimension(3,nnd), intent(in)  X,
integer, intent(in)  nNd 
)
private

Define the transformation T: K_ref –> K.

1- Define the transformation T: K_ref –> K i.e. compute the ttansformation 'T' coefficients depending on the type of the cell K and on its node coordinates 'X'

  1. computate the coordinates of the points 'Ty' in K

Definition at line 250 of file geoTsf_mod.F90.

Here is the call graph for this function:

◆ geotsf_clear()

subroutine geotsf_mod::geotsf_clear ( type(geotsf), intent(inout)  g)
private

Destructor for geoTsf type.

Definition at line 188 of file geoTsf_mod.F90.

◆ geotsf_create()

type(geotsf) function geotsf_mod::geotsf_create ( integer, intent(in)  dim,
integer, intent(in)  nn,
real(rp), dimension(dim,nn), intent(in)  y 
)
private

Constructor for the type geoTsf

  1. Allocate memory for g%y, g%Ty, g%Jy and g%Cy
  2. set g%y: the coordinates for the points \(y\) in \( K_{\rm ref}\)
  • OUTPUT:
    • g : geometrical transformation
  • INPUT:
    • dim : the dimension of \( K_{\rm ref}\),
    • nn : the number of points in \( y\)
    • y : the coordinates of the points in \( y\).

Definition at line 216 of file geoTsf_mod.F90.

Here is the call graph for this function:

◆ geotsf_t_inv()

subroutine, public geotsf_mod::geotsf_t_inv ( real(rp), dimension(3), intent(out)  z,
real(rp), dimension(3), intent(in)  X,
real(rp), dimension(3,n), intent(inout)  Y,
integer, intent(in)  n,
integer, intent(in)  ct 
)

x = 3D node 3D cell of type ct with nodes = A, B, C, ... given by Y(:,1)=A, Y(:,2)=B, ..., Y(:,n)

Output : z R3

EDG : AX = z(1)AB + n, n orthogonal to AB TRG : AX = z(1)AB + z(2)AC + n, n orthogonal to ABC TET : AX = AX = z(1)AB + z(2)AC + z(2)AD

Output : Y is modified : Y(:,1) = AX Y(:,2) = AB Y(:,3) = AC ...

Definition at line 683 of file geoTsf_mod.F90.

Here is the call graph for this function: