fpakc/src/modules/mesh/2DCyl/moduleMesh2DCyl.f90
JGonzalez 6f24b5f1f6 Small changes before trying something big
I think that creating arrays with self%nNodes takes a lot of time.
I'm trying now to pass the number of nodes as argument.
2023-01-05 20:32:45 +01:00

1400 lines
38 KiB
Fortran

!moduleMesh2DCyl: 2D aXial symmetric extension of generic mesh from GMSH format.
! x == z
! y == r
! z == theta (unused)
MODULE moduleMesh2DCyl
USE moduleMesh
USE moduleMeshBoundary
IMPLICIT NONE
!Values for Gauss integral
REAL(8), PARAMETER:: corQuad(1:3) = (/ -DSQRT(3.D0/5.D0), 0.D0, DSQRT(3.D0/5.D0) /)
REAL(8), PARAMETER:: wQuad(1:3) = (/ 5.D0/9.D0, 8.D0/9.D0, 5.D0/9.D0 /)
REAL(8), PARAMETER:: Xi1Tria(1:4) = (/ 1.D0/3.D0, 1.D0/5.D0, 3.D0/5.D0, 1.D0/5.D0 /)
REAL(8), PARAMETER:: Xi2Tria(1:4) = (/ 1.D0/3.D0, 1.D0/5.D0, 1.D0/5.D0, 3.D0/5.D0 /)
REAL(8), PARAMETER:: wTria(1:4) = (/ -27.D0/96.D0, 25.D0/96.D0, 25.D0/96.D0, 25.D0/96.D0 /)
TYPE, PUBLIC, EXTENDS(meshNode):: meshNode2DCyl
!Element coordinates
REAL(8):: r = 0.D0, z = 0.D0
CONTAINS
PROCEDURE, PASS:: init => initNode2DCyl
PROCEDURE, PASS:: getCoordinates => getCoord2DCyl
END TYPE meshNode2DCyl
TYPE, PUBLIC, EXTENDS(meshEdge):: meshEdge2DCyl
!Element coordinates
REAL(8):: r(1:2) = 0.D0, z(1:2) = 0.D0
!Connectivity to nodes
CLASS(meshNode), POINTER:: n1 => NULL(), n2 => NULL()
CONTAINS
PROCEDURE, PASS:: init => initEdge2DCyl
PROCEDURE, PASS:: getNodes => getNodes2DCyl
PROCEDURE, PASS:: intersection => intersection2DCylEdge
PROCEDURE, PASS:: randPos => randPosEdge
END TYPE meshEdge2DCyl
TYPE, PUBLIC, ABSTRACT, EXTENDS(meshCell):: meshCell2DCyl
CONTAINS
PROCEDURE, PASS:: detJac => detJ2DCyl
PROCEDURE, PASS:: invJac => invJ2DCyl
PROCEDURE(partialDer_interface), DEFERRED, PASS, PRIVATE:: partialDer
END TYPE meshCell2DCyl
ABSTRACT INTERFACE
PURE SUBROUTINE partialDer_interface(self, dPsi, dz, dr)
IMPORT meshCell2DCyl
CLASS(meshCell2DCyl), INTENT(in):: self
REAL(8), INTENT(in):: dPsi(1:3,1:self%nNodes)
REAL(8), INTENT(out), DIMENSION(1:2):: dz, dr
END SUBROUTINE partialDer_interface
END INTERFACE
!Quadrilateral volume element
TYPE, PUBLIC, EXTENDS(meshCell2DCyl):: meshCell2DCylQuad
!Element coordinates
REAL(8):: r(1:4) = 0.D0, z(1:4) = 0.D0
!Connectivity to nodes
CLASS(meshNode), POINTER:: n1 => NULL(), n2 => NULL(), n3 => NULL(), n4 => NULL()
!Connectivity to adjacent elements
CLASS(meshElement), POINTER:: e1 => NULL(), e2 => NULL(), e3 => NULL(), e4 => NULL()
REAL(8):: arNodes(1:4) = 0.D0
CONTAINS
PROCEDURE, PASS:: init => initCellQuad2DCyl
PROCEDURE, PASS:: randPos => randPosCellQuad
PROCEDURE, PASS:: area => areaQuad
PROCEDURE, PASS:: fPsi => fPsiQuad
PROCEDURE, PASS:: dPsi => dPsiQuad
PROCEDURE, PASS, PRIVATE:: partialDer => partialDerQuad
PROCEDURE, PASS:: elemK => elemKQuad
PROCEDURE, PASS:: elemF => elemFQuad
PROCEDURE, PASS:: gatherElectricField => gatherEFQuad
PROCEDURE, PASS:: gatherMagneticField => gatherMFQuad
PROCEDURE, NOPASS:: inside => insideQuad
PROCEDURE, PASS:: getNodes => getNodesQuad
PROCEDURE, PASS:: phy2log => phy2logQuad
PROCEDURE, PASS:: nextElement => nextElementQuad
END TYPE meshCell2DCylQuad
!Triangular volume element
TYPE, PUBLIC, EXTENDS(meshCell2DCyl):: meshCell2DCylTria
!Element coordinates
REAL(8):: r(1:3) = 0.D0, z(1:3) = 0.D0
!Connectivity to nodes
CLASS(meshNode), POINTER:: n1 => NULL(), n2 => NULL(), n3 => NULL()
!Connectivity to adjacent elements
CLASS(meshElement), POINTER:: e1 => NULL(), e2 => NULL(), e3 => NULL()
REAL(8):: arNodes(1:3) = 0.D0
CONTAINS
PROCEDURE, PASS:: init => initCellTria2DCyl
PROCEDURE, PASS:: randPos => randPosCellTria
PROCEDURE, PASS:: area => areaTria
PROCEDURE, PASS:: fPsi => fPsiTria
PROCEDURE, PASS:: dPsi => dPsiTria
PROCEDURE, PASS, PRIVATE:: partialDer => partialDerTria
PROCEDURE, PASS:: elemK => elemKTria
PROCEDURE, PASS:: elemF => elemFTria
PROCEDURE, PASS:: gatherElectricField => gatherEFTria
PROCEDURE, PASS:: gatherMagneticField => gatherMFTria
PROCEDURE, NOPASS:: inside => insideTria
PROCEDURE, PASS:: getNodes => getNodesTria
PROCEDURE, PASS:: phy2log => phy2logTria
PROCEDURE, PASS:: nextElement => nextElementTria
END TYPE meshCell2DCylTria
CONTAINS
!NODE FUNCTIONS
!Inits node element
SUBROUTINE initNode2DCyl(self, n, r)
USE moduleSpecies
USE moduleRefParam
USE OMP_LIB
IMPLICIT NONE
CLASS(meshNode2DCyl), INTENT(out):: self
INTEGER, INTENT(in):: n
REAL(8), INTENT(in):: r(1:3)
self%n = n
self%z = r(1)/L_ref
self%r = r(2)/L_ref
!Node volume, to be determined in mesh
self%v = 0.D0
!Allocates output:
ALLOCATE(self%output(1:nSpecies))
CALL OMP_INIT_LOCK(self%lock)
END SUBROUTINE initNode2DCyl
!Get coordinates from node
PURE FUNCTION getCoord2DCyl(self) RESULT(r)
IMPLICIT NONE
CLASS(meshNode2DCyl), INTENT(in):: self
REAL(8):: r(1:3)
r = (/self%z, self%r, 0.D0/)
END FUNCTION getCoord2DCyl
!EDGE FUNCTIONS
!Inits edge element
SUBROUTINE initEdge2DCyl(self, n, p, bt, physicalSurface)
USE moduleSpecies
USE moduleBoundary
USE moduleErrors
IMPLICIT NONE
CLASS(meshEdge2DCyl), INTENT(out):: self
INTEGER, INTENT(in):: n
INTEGER, INTENT(in):: p(:)
INTEGER, INTENT(in):: bt
INTEGER, INTENT(in):: physicalSurface
REAL(8), DIMENSION(1:3):: r1, r2
INTEGER:: s
self%n = n
self%n1 => mesh%nodes(p(1))%obj
self%n2 => mesh%nodes(p(2))%obj
!Get element coordinates
r1 = self%n1%getCoordinates()
r2 = self%n2%getCoordinates()
self%z = (/r1(1), r2(1)/)
self%r = (/r1(2), r2(2)/)
self%weight = SUM(self%r)*5.D-1
!Normal vector
self%normal = (/ -(self%r(2)-self%r(1)), &
self%z(2)-self%z(1) , &
0.D0 /)
self%normal = self%normal/NORM2(self%normal)
!Boundary index
self%boundary => boundary(bt)
ALLOCATE(self%fboundary(1:nSpecies))
!Assign functions to boundary
DO s = 1, nSpecies
CALL pointBoundaryFunction(self, s)
END DO
!Physical surface
self%physicalSurface = physicalSurface
END SUBROUTINE initEdge2DCyl
!Get nodes from edge
PURE FUNCTION getNodes2DCyl(self) RESULT(n)
IMPLICIT NONE
CLASS(meshEdge2DCyl), INTENT(in):: self
INTEGER, ALLOCATABLE:: n(:)
ALLOCATE(n(1:2))
n = (/self%n1%n, self%n2%n /)
END FUNCTION getNodes2DCyl
PURE FUNCTION intersection2DCylEdge(self, r0) RESULT(r)
IMPLICIT NONE
CLASS(meshEdge2DCyl), INTENT(in):: self
REAL(8), DIMENSION(1:3), INTENT(in):: r0
REAL(8), DIMENSION(1:3):: r
REAL(8), DIMENSION(1:3):: edge0, edgeV
REAL(8):: tI
edge0 = (/self%z(1), self%r(1), 0.D0 /)
edgeV = (/self%z(2), self%r(2), 0.D0 /) - edge0
tI = DOT_PRODUCT(r0 - edge0, edgeV)/DOT_PRODUCT(edgeV, edgeV)
r = edge0 + tI*edgeV
END FUNCTION intersection2DCylEdge
!Calculates a random position in edge
FUNCTION randPosEdge(self) RESULT(r)
USE moduleRandom
IMPLICIT NONE
CLASS(meshEdge2DCyl), INTENT(in):: self
REAL(8):: rnd
REAL(8):: dr, dz
REAL(8):: r(1:3)
rnd = random()
dr = self%r(2) - self%r(1)
dz = self%z(2) - self%z(1)
IF (dr /= 0.D0) THEN
r(2) = dr * DSQRT(rnd) + self%r(1)
r(1) = dz * (r(2) - self%r(1))/dr + self%z(1)
ELSE
r(2) = self%r(1)
r(1) = dz * rnd + self%z(1)
END IF
r(3) = 0.D0
END FUNCTION randPosEdge
!VOLUME FUNCTIONS
!QUAD FUNCTIONS
!Inits quadrilateral element
SUBROUTINE initCellQuad2DCyl(self, n, p, nodes)
USE moduleRefParam
IMPLICIT NONE
CLASS(meshCell2DCylQuad), INTENT(out):: self
INTEGER, INTENT(in):: n
INTEGER, INTENT(in):: p(:)
TYPE(meshNodeCont), INTENT(in), TARGET:: nodes(:)
REAL(8), DIMENSION(1:3):: r1, r2, r3, r4
self%n = n
self%nNodes = SIZE(p)
self%n1 => nodes(p(1))%obj
self%n2 => nodes(p(2))%obj
self%n3 => nodes(p(3))%obj
self%n4 => nodes(p(4))%obj
!Get element coordinates
r1 = self%n1%getCoordinates()
r2 = self%n2%getCoordinates()
r3 = self%n3%getCoordinates()
r4 = self%n4%getCoordinates()
self%z = (/r1(1), r2(1), r3(1), r4(1)/)
self%r = (/r1(2), r2(2), r3(2), r4(2)/)
!Assign node volume
CALL self%area()
self%n1%v = self%n1%v + self%arNodes(1)
self%n2%v = self%n2%v + self%arNodes(2)
self%n3%v = self%n3%v + self%arNodes(3)
self%n4%v = self%n4%v + self%arNodes(4)
CALL OMP_INIT_LOCK(self%lock)
ALLOCATE(self%listPart_in(1:nSpecies))
ALLOCATE(self%totalWeight(1:nSpecies))
END SUBROUTINE initCellQuad2DCyl
!Computes element area
PURE SUBROUTINE areaQuad(self)
USE moduleConstParam, ONLY: PI8
IMPLICIT NONE
CLASS(meshCell2DCylQuad), INTENT(inout):: self
REAL(8):: r, Xi(1:3)
REAL(8):: detJ
REAL(8):: fPsi(1:4), fPsi_node(1:4)
self%volume = 0.D0
self%arNodes = 0.D0
!2D 1 point Gauss Quad Integral
Xi = 0.D0
detJ = self%detJac(Xi)*PI8 !4*2*pi
fPsi = self%fPsi(Xi)
!Computes total volume of the cell
r = DOT_PRODUCT(fPsi,self%r)
self%volume = r*detJ
!Computes volume per node
Xi = (/-5.D-1, -5.D-1, 0.D0/)
r = self%gatherF(Xi, self%r)
self%arNodes(1) = fPsi(1)*r*detJ
Xi = (/ 5.D-1, -5.D-1, 0.D0/)
r = self%gatherF(Xi, self%r)
self%arNodes(2) = fPsi(2)*r*detJ
Xi = (/ 5.D-1, 5.D-1, 0.D0/)
r = self%gatherF(Xi, self%r)
self%arNodes(3) = fPsi(3)*r*detJ
Xi = (/-5.D-1, 5.D-1, 0.D0/)
r = self%gatherF(Xi, self%r)
self%arNodes(4) = fPsi(4)*r*detJ
END SUBROUTINE areaQuad
!Computes element functions in point Xi
PURE FUNCTION fPsiQuad(self, Xi) RESULT(fPsi)
IMPLICIT NONE
CLASS(meshCell2DCylQuad), INTENT(in):: self
REAL(8), INTENT(in):: Xi(1:3)
REAL(8):: fPsi(1:self%nNodes)
fPsi(1) = (1.D0-Xi(1)) * (1.D0-Xi(2))
fPsi(2) = (1.D0+Xi(1)) * (1.D0-Xi(2))
fPsi(3) = (1.D0+Xi(1)) * (1.D0+Xi(2))
fPsi(4) = (1.D0-Xi(1)) * (1.D0+Xi(2))
fPsi = fPsi*0.25D0
END FUNCTION fPsiQuad
!Derivative element function at coordinates Xi
PURE FUNCTION dPsiQuad(self, Xi) RESULT(dPsi)
IMPLICIT NONE
CLASS(meshCell2DCylQuad), INTENT(in):: self
REAL(8), INTENT(in):: Xi(1:3)
REAL(8):: dPsi(1:3,1:self%nNodes)
dPsi = 0.D0
dPsi(1,:) = (/ -(1.D0 - Xi(2)), &
(1.D0 - Xi(2)), &
(1.D0 + Xi(2)), &
-(1.D0 + Xi(2)) /)
dPsi(2,:) = (/ -(1.D0 - Xi(1)), &
-(1.D0 + Xi(1)), &
(1.D0 + Xi(1)), &
(1.D0 - Xi(1)) /)
dPsi = dPsi * 0.25D0
END FUNCTION dPsiQuad
!Partial derivative in global coordinates
PURE SUBROUTINE partialDerQuad(self, dPsi, dz, dr)
IMPLICIT NONE
CLASS(meshCell2DCylQuad), INTENT(in):: self
REAL(8), INTENT(in):: dPsi(1:3,1:self%nNodes)
REAL(8), INTENT(out), DIMENSION(1:2):: dz, dr
dz = (/ DOT_PRODUCT(dPsi(1,1:4),self%z(1:4)), &
DOT_PRODUCT(dPsi(2,1:4),self%z(1:4)) /)
dr = (/ DOT_PRODUCT(dPsi(1,1:4),self%r(1:4)), &
DOT_PRODUCT(dPsi(2,1:4),self%r(1:4)) /)
END SUBROUTINE partialDerQuad
!Random position in quadrilateral volume
FUNCTION randPosCellQuad(self) RESULT(r)
USE moduleRandom
IMPLICIT NONE
CLASS(meshCell2DCylQuad), INTENT(in):: self
REAL(8):: r(1:3)
REAL(8):: Xi(1:3)
REAL(8):: fPsi(1:4)
Xi(1) = random(-1.D0, 1.D0)
Xi(2) = random(-1.D0, 1.D0)
Xi(3) = 0.D0
fPsi = self%fPsi(Xi)
r(1) = DOT_PRODUCT(fPsi, self%z)
r(2) = DOT_PRODUCT(fPsi, self%r)
r(3) = 0.D0
END FUNCTION randPosCellQuad
!Computes element local stiffness matrix
PURE FUNCTION elemKQuad(self) RESULT(localK)
USE moduleConstParam, ONLY: PI2
IMPLICIT NONE
CLASS(meshCell2DCylQuad), INTENT(in):: self
REAL(8):: localK(1:self%nNodes,1:self%nNodes)
REAL(8):: Xi(1:3)
REAL(8):: fPsi(1:4), dPsi(1:3,1:4)
REAL(8):: r
REAL(8):: invJ(1:3,1:3), detJ
INTEGER:: l, m
localK=0.D0
Xi=0.D0
!Start 2D Gauss Quad Integral
DO l=1, 3
Xi(2) = corQuad(l)
DO m = 1, 3
Xi(1) = corQuad(m)
fPsi = self%fPsi(Xi)
dPsi = self%dPsi(Xi)
detJ = self%detJac(Xi,dPsi)
invJ = self%invJac(Xi,dPsi)
r = DOT_PRODUCT(fPsi,self%r)
localK = localK + MATMUL(TRANSPOSE(MATMUL(invJ,dPsi)), &
MATMUL(invJ,dPsi))* &
r*wQuad(l)*wQuad(m)/detJ
END DO
END DO
localK = localK*PI2
END FUNCTION elemKQuad
!Computes the local source vector for a force f
PURE FUNCTION elemFQuad(self, source) RESULT(localF)
USE moduleConstParam, ONLY: PI2
IMPLICIT NONE
CLASS(meshCell2DCylQuad), INTENT(in):: self
REAL(8), INTENT(in):: source(1:self%nNodes)
REAL(8):: localF(1:self%nNodes)
REAL(8):: Xi(1:3)
REAL(8):: fPsi(1:4)
REAL(8):: r
REAL(8):: detJ, f
INTEGER:: l, m
localF = 0.D0
Xi = 0.D0
DO l=1, 3
Xi(1) = corQuad(l)
DO m = 1, 3
Xi(2) = corQuad(m)
detJ = self%detJac(Xi)
fPsi = self%fPsi(Xi)
r = DOT_PRODUCT(fPsi,self%r)
f = DOT_PRODUCT(fPsi,source)
localF = localF + r*f*fPsi*wQuad(l)*wQuad(m)*detJ
END DO
END DO
localF = localF*PI2
END FUNCTION elemFQuad
PURE FUNCTION gatherEFQuad(self, Xi) RESULT(array)
IMPLICIT NONE
CLASS(meshCell2DCylQuad), INTENT(in):: self
REAL(8), INTENT(in):: Xi(1:3)
REAL(8):: array(1:3)
REAL(8):: phi(1:4)
phi = (/ self%n1%emData%phi, &
self%n2%emData%phi, &
self%n3%emData%phi, &
self%n4%emData%phi /)
array = -self%gatherDF(Xi, phi)
END FUNCTION gatherEFQuad
PURE FUNCTION gatherMFQuad(self, Xi) RESULT(array)
IMPLICIT NONE
CLASS(meshCell2DCylQuad), INTENT(in):: self
REAL(8), INTENT(in):: Xi(1:3)
REAL(8):: array(1:3)
REAL(8):: B(1:4,1:3)
B(:,1) = (/ self%n1%emData%B(1), &
self%n2%emData%B(1), &
self%n3%emData%B(1), &
self%n4%emData%B(1) /)
B(:,2) = (/ self%n1%emData%B(2), &
self%n2%emData%B(2), &
self%n3%emData%B(2), &
self%n4%emData%B(2) /)
B(:,3) = (/ self%n1%emData%B(3), &
self%n2%emData%B(3), &
self%n3%emData%B(3), &
self%n4%emData%B(3) /)
array = self%gatherF(Xi, 3, B)
END FUNCTION gatherMFQuad
!Checks if a particle is inside a quad element
PURE FUNCTION insideQuad(Xi) RESULT(ins)
IMPLICIT NONE
REAL(8), INTENT(in):: Xi(1:3)
LOGICAL:: ins
ins = (Xi(1) >= -1.D0 .AND. Xi(1) <= 1.D0) .AND. &
(Xi(2) >= -1.D0 .AND. Xi(2) <= 1.D0)
END FUNCTION insideQuad
!Gets nodes from quadrilateral element
PURE FUNCTION getNodesQuad(self) RESULT(n)
IMPLICIT NONE
CLASS(meshCell2DCylQuad), INTENT(in):: self
INTEGER:: n(1:self%nNodes)
n = (/self%n1%n, self%n2%n, self%n3%n, self%n4%n /)
END FUNCTION getNodesQuad
!Transforms physical coordinates to element coordinates
PURE FUNCTION phy2logQuad(self,r) RESULT(Xi)
IMPLICIT NONE
CLASS(meshCell2DCylQuad), INTENT(in):: self
REAL(8), INTENT(in):: r(1:3)
REAL(8):: Xi(1:3)
REAL(8):: XiO(1:3), detJ, invJ(1:3,1:3), f(1:3)
REAL(8):: dPsi(1:3,1:4), fPsi(1:4)
REAL(8):: conv
!Iterative newton method to transform coordinates
conv = 1.D0
XiO = 0.D0
DO WHILE(conv > 1.D-2)
dPsi = self%dPsi(XiO)
invJ = self%invJac(XiO, dPsi)
detJ = self%detJac(XiO, dPsi)
fPsi = self%fPsi(XiO)
f = (/ DOT_PRODUCT(fPsi,self%z), &
DOT_PRODUCT(fPsi,self%r), &
0.D0 /)
f = f - r
Xi = XiO - MATMUL(invJ, f)/detJ
conv = MAXVAL(DABS(Xi-XiO),1)
XiO = Xi
END DO
END FUNCTION phy2logQuad
!Gets the next element for a logical position Xi
SUBROUTINE nextElementQuad(self, Xi, nextElement)
IMPLICIT NONE
CLASS(meshCell2DCylQuad), INTENT(in):: self
REAL(8), INTENT(in):: Xi(1:3)
CLASS(meshElement), POINTER, INTENT(out):: nextElement
REAL(8):: XiArray(1:4)
INTEGER:: nextInt
XiArray = (/ -Xi(2), Xi(1), Xi(2), -Xi(1) /)
nextInt = MAXLOC(XiArray,1)
!Selects the higher value of directions and searches in that direction
NULLIFY(nextElement)
SELECT CASE (nextInt)
CASE (1)
nextElement => self%e1
CASE (2)
nextElement => self%e2
CASE (3)
nextElement => self%e3
CASE (4)
nextElement => self%e4
END SELECT
END SUBROUTINE nextElementQuad
!TRIA ELEMENT
!Init tria element
SUBROUTINE initCellTria2DCyl(self, n, p, nodes)
USE moduleRefParam
IMPLICIT NONE
CLASS(meshCell2DCylTria), INTENT(out):: self
INTEGER, INTENT(in):: n
INTEGER, INTENT(in):: p(:)
TYPE(meshNodeCont), INTENT(in), TARGET:: nodes(:)
REAL(8), DIMENSION(1:3):: r1, r2, r3
!Assign node index
self%n = n
!Assign nomber of nodes to cell
self%nNodes = SIZE(p)
!Assign nodes to element
self%n1 => nodes(p(1))%obj
self%n2 => nodes(p(2))%obj
self%n3 => nodes(p(3))%obj
!Get element coordinates
r1 = self%n1%getCoordinates()
r2 = self%n2%getCoordinates()
r3 = self%n3%getCoordinates()
self%z = (/r1(1), r2(1), r3(1)/)
self%r = (/r1(2), r2(2), r3(2)/)
!Assign node volume
CALL self%area()
self%n1%v = self%n1%v + self%arNodes(1)
self%n2%v = self%n2%v + self%arNodes(2)
self%n3%v = self%n3%v + self%arNodes(3)
CALL OMP_INIT_LOCK(self%lock)
ALLOCATE(self%listPart_in(1:nSpecies))
ALLOCATE(self%totalWeight(1:nSpecies))
END SUBROUTINE initCellTria2DCyl
!Random position in quadrilateral volume
FUNCTION randPosCellTria(self) RESULT(r)
USE moduleRandom
IMPLICIT NONE
CLASS(meshCell2DCylTria), INTENT(in):: self
REAL(8):: r(1:3)
REAL(8):: Xi(1:3)
REAL(8):: fPsi(1:3)
Xi(1) = random( 0.D0, 1.D0)
Xi(2) = random( 0.D0, 1.D0 - Xi(1))
Xi(3) = 0.D0
fPsi = self%fPsi(Xi)
r(1) = DOT_PRODUCT(fPsi, self%z)
r(2) = DOT_PRODUCT(fPsi, self%r)
r(3) = 0.D0
END FUNCTION randPosCellTria
!Calculates area for triangular element
PURE SUBROUTINE areaTria(self)
USE moduleConstParam, ONLY: PI
IMPLICIT NONE
CLASS(meshCell2DCylTria), INTENT(inout):: self
REAL(8):: Xi(1:3)
REAL(8):: r
REAL(8):: detJ
REAL(8):: fPsi(1:3)
self%volume = 0.D0
self%arNodes = 0.D0
!2D 1 point Gauss Quad Integral
Xi = (/1.D0/3.D0, 1.D0/3.D0, 0.D0 /)
detJ = self%detJac(Xi)*PI !2PI*1/2
fPsi = self%fPsi(Xi)
!Computes total volume of the cell
r = DOT_PRODUCT(fPsi,self%r)
self%volume = r*detJ
!Computes volume per node
self%arNodes = fPsi*r*detJ
END SUBROUTINE areaTria
!Shape functions for triangular element
PURE FUNCTION fPsiTria(self, Xi) RESULT(fPsi)
IMPLICIT NONE
CLASS(meshCell2DCylTria), INTENT(in):: self
REAL(8), INTENT(in):: Xi(1:3)
REAL(8):: fPsi(1:self%nNodes)
fPsi(1) = 1.D0 - Xi(1) - Xi(2)
fPsi(2) = Xi(1)
fPsi(3) = Xi(2)
END FUNCTION fPsiTria
!Derivative element function at coordinates Xi
PURE FUNCTION dPsiTria(self, Xi) RESULT(dPsi)
IMPLICIT NONE
CLASS(meshCell2DCylTria), INTENT(in):: self
REAL(8), INTENT(in):: Xi(1:3)
REAL(8):: dPsi(1:3,1:self%nNodes)
dPsi = 0.D0
dPsi(1,:) = (/ -1.D0, 1.D0, 0.D0 /)
dPsi(2,:) = (/ -1.D0, 0.D0, 1.D0 /)
END FUNCTION dPsiTria
PURE SUBROUTINE partialDerTria(self, dPsi, dz, dr)
IMPLICIT NONE
CLASS(meshCell2DCylTria), INTENT(in):: self
REAL(8), INTENT(in):: dPsi(1:3,1:self%nNodes)
REAL(8), INTENT(out), DIMENSION(1:2):: dz, dr
dz = (/ DOT_PRODUCT(dPsi(1,:),self%z), &
DOT_PRODUCT(dPsi(2,:),self%z) /)
dr = (/ DOT_PRODUCT(dPsi(1,:),self%r), &
DOT_PRODUCT(dPsi(2,:),self%r) /)
END SUBROUTINE partialDerTria
!Computes element local stiffness matrix
PURE FUNCTION elemKTria(self) RESULT(localK)
USE moduleConstParam, ONLY: PI2
IMPLICIT NONE
CLASS(meshCell2DCylTria), INTENT(in):: self
REAL(8):: localK(1:self%nNodes,1:self%nNodes)
REAL(8):: Xi(1:3)
REAL(8):: r
REAL(8):: fPsi(1:3), dPsi(1:3,1:3)
REAL(8):: invJ(1:3,1:3), detJ
INTEGER:: l
localK=0.D0
Xi=0.D0
!Start 2D Gauss Quad Integral
DO l=1, 4
Xi(1) = Xi1Tria(l)
Xi(2) = Xi2Tria(l)
dPsi = self%dPsi(Xi)
detJ = self%detJac(Xi,dPsi)
invJ = self%invJac(Xi,dPsi)
fPsi = self%fPsi(Xi)
r = DOT_PRODUCT(fPsi,self%r)
localK = localK + MATMUL(TRANSPOSE(MATMUL(invJ,dPsi)),MATMUL(invJ,dPsi))*r*wTria(l)/detJ
END DO
localK = localK*PI2
END FUNCTION elemKTria
!Computes element local source vector
PURE FUNCTION elemFTria(self, source) RESULT(localF)
USE moduleConstParam, ONLY: PI2
IMPLICIT NONE
CLASS(meshCell2DCylTria), INTENT(in):: self
REAL(8), INTENT(in):: source(1:self%nNodes)
REAL(8):: localF(1:self%nNodes)
REAL(8):: fPsi(1:3)
REAL(8):: Xi(1:3)
REAL(8):: r
REAL(8):: detJ, f
INTEGER:: l
localF = 0.D0
Xi = 0.D0
!Start 2D Gauss Quad Integral
DO l=1, 4
Xi(1) = Xi1Tria(l)
Xi(2) = Xi2Tria(l)
detJ = self%detJac(Xi)
fPsi = self%fPsi(Xi)
r = DOT_PRODUCT(fPsi,self%r)
f = DOT_PRODUCT(fPsi,source)
localF = localF + r*f*fPsi*wTria(l)*detJ
END DO
localF = localF*PI2
END FUNCTION elemFTria
PURE FUNCTION gatherEFTria(self, Xi) RESULT(array)
IMPLICIT NONE
CLASS(meshCell2DCylTria), INTENT(in):: self
REAL(8), INTENT(in):: Xi(1:3)
REAL(8):: array(1:3)
REAL(8):: phi(1:3)
phi = (/ self%n1%emData%phi, &
self%n2%emData%phi, &
self%n3%emData%phi /)
array = -self%gatherDF(Xi, phi)
END FUNCTION gatherEFTria
PURE FUNCTION gatherMFTria(self, Xi) RESULT(array)
IMPLICIT NONE
CLASS(meshCell2DCylTria), INTENT(in):: self
REAL(8), INTENT(in):: Xi(1:3)
REAL(8):: array(1:3)
REAL(8):: B(1:3,1:3)
B(:,1) = (/ self%n1%emData%B(1), &
self%n2%emData%B(1), &
self%n3%emData%B(1) /)
B(:,2) = (/ self%n1%emData%B(2), &
self%n2%emData%B(2), &
self%n3%emData%B(2) /)
B(:,3) = (/ self%n1%emData%B(3), &
self%n2%emData%B(3), &
self%n3%emData%B(3) /)
array = self%gatherF(Xi, 3, B)
END FUNCTION gatherMFTria
PURE FUNCTION insideTria(Xi) RESULT(ins)
IMPLICIT NONE
REAL(8), INTENT(in):: Xi(1:3)
LOGICAL:: ins
ins = Xi(1) >= 0.D0 .AND. &
Xi(2) >= 0.D0 .AND. &
1.D0 - Xi(1) - Xi(2) >= 0.D0
END FUNCTION insideTria
!Gets node indexes from triangular element
PURE FUNCTION getNodesTria(self) RESULT(n)
IMPLICIT NONE
CLASS(meshCell2DCylTria), INTENT(in):: self
INTEGER:: n(1:self%nNodes)
n = (/self%n1%n, self%n2%n, self%n3%n /)
END FUNCTION getNodesTria
!Transforms physical coordinates to element coordinates
PURE FUNCTION phy2logTria(self,r) RESULT(Xi)
IMPLICIT NONE
CLASS(meshCell2DCylTria), INTENT(in):: self
REAL(8), INTENT(in):: r(1:3)
REAL(8):: Xi(1:3)
REAL(8):: invJ(1:3,1:3), detJ
REAL(8):: deltaR(1:3)
REAL(8):: dPsi(1:3,1:3)
!Direct method to convert coordinates
Xi = 0.D0
deltaR = (/ r(1) - self%z(1), r(2) - self%r(1), 0.D0 /)
dPsi = self%dPsi(Xi)
invJ = self%invJac(Xi, dPsi)
detJ = self%detJac(Xi, dPsi)
Xi = MATMUL(invJ,deltaR)/detJ
END FUNCTION phy2logTria
SUBROUTINE nextElementTria(self, Xi, nextElement)
IMPLICIT NONE
CLASS(meshCell2DCylTria), INTENT(in):: self
REAL(8), INTENT(in):: Xi(1:3)
CLASS(meshElement), POINTER, INTENT(out):: nextElement
REAL(8):: XiArray(1:3)
INTEGER:: nextInt
XiArray = (/ Xi(2), 1.D0-Xi(1)-Xi(2), Xi(1) /)
nextInt = MINLOC(XiArray,1)
NULLIFY(nextElement)
SELECT CASE (nextInt)
CASE (1)
nextElement => self%e1
CASE (2)
nextElement => self%e2
CASE (3)
nextElement => self%e3
END SELECT
END SUBROUTINE nextElementTria
!COMMON FUNCTIONS FOR CYLINDRICAL VOLUME ELEMENTS
!Computes element Jacobian determinant
PURE FUNCTION detJ2DCyl(self, Xi, dPsi_in) RESULT(dJ)
IMPLICIT NONE
CLASS(meshCell2DCyl), INTENT(in):: self
REAL(8), INTENT(in):: Xi(1:3)
REAL(8), INTENT(in), OPTIONAL:: dPsi_in(1:3,1:self%nNodes)
REAL(8):: dJ
REAL(8):: dPsi(1:3,1:self%nNodes)
REAL(8):: dz(1:2), dr(1:2)
IF(PRESENT(dPsi_in)) THEN
dPsi = dPsi_in
ELSE
dPsi = self%dPsi(Xi)
END IF
CALL self%partialDer(dPsi, dz, dr)
dJ = dz(1)*dr(2)-dz(2)*dr(1)
END FUNCTION detJ2DCyl
!Computes element Jacobian inverse matrix (without determinant)
PURE FUNCTION invJ2DCyl(self,Xi,dPsi_in) RESULT(invJ)
IMPLICIT NONE
CLASS(meshCell2DCyl), INTENT(in):: self
REAL(8), INTENT(in):: Xi(1:3)
REAL(8), INTENT(in), OPTIONAL:: dPsi_in(1:3,1:self%nNodes)
REAL(8):: invJ(1:3,1:3)
REAL(8):: dPsi(1:3,1:self%nNodes)
REAL(8):: dz(1:2), dr(1:2)
IF(PRESENT(dPsi_in)) THEN
dPsi=dPsi_in
ELSE
dPsi = self%dPsi(Xi)
END IF
invJ = 0.D0
CALL self%partialDer(dPsi, dz, dr)
invJ(1,1:2) = (/ dr(2), -dz(2) /)
invJ(2,1:2) = (/ -dr(1), dz(1) /)
END FUNCTION invJ2DCyl
SUBROUTINE connectMesh2DCyl(self)
IMPLICIT NONE
CLASS(meshGeneric), INTENT(inout):: self
INTEGER:: e, et
DO e = 1, self%numCells
!Connect Cell-Cell
DO et = 1, self%numCells
IF (e /= et) THEN
CALL connectCellCell(self%cells(e)%obj, self%cells(et)%obj)
END IF
END DO
SELECT TYPE(self)
TYPE IS(meshParticles)
!Connect Cell-Edge
DO et = 1, self%numEdges
CALL connectCellEdge(self%cells(e)%obj, self%edges(et)%obj)
END DO
END SELECT
END DO
END SUBROUTINE connectMesh2DCyl
!Selects type of elements to build connection
SUBROUTINE connectCellCell(elemA, elemB)
IMPLICIT NONE
CLASS(meshCell), INTENT(inout):: elemA
CLASS(meshCell), INTENT(inout):: elemB
SELECT TYPE(elemA)
TYPE IS(meshCell2DCylQuad)
!Element A is a quadrilateral
SELECT TYPE(elemB)
TYPE IS(meshCell2DCylQuad)
!Element B is a quadrilateral
CALL connectQuadQuad(elemA, elemB)
TYPE IS(meshCell2DCylTria)
!Element B is a triangle
CALL connectQuadTria(elemA, elemB)
END SELECT
TYPE IS(meshCell2DCylTria)
!Element A is a Triangle
SELECT TYPE(elemB)
TYPE IS(meshCell2DCylQuad)
!Element B is a quadrilateral
CALL connectQuadTria(elemB, elemA)
TYPE IS(meshCell2DCylTria)
!Element B is a triangle
CALL connectTriaTria(elemA, elemB)
END SELECT
END SELECT
END SUBROUTINE connectCellCell
SUBROUTINE connectCellEdge(elemA, elemB)
IMPLICIT NONE
CLASS(meshCell), INTENT(inout):: elemA
CLASS(meshEdge), INTENT(inout):: elemB
SELECT TYPE(elemB)
CLASS IS(meshEdge2DCyl)
SELECT TYPE(elemA)
TYPE IS(meshCell2DCylQuad)
!Element A is a quadrilateral
CALL connectQuadEdge(elemA, elemB)
TYPE IS(meshCell2DCylTria)
!Element A is a triangle
CALL connectTriaEdge(elemA, elemB)
END SELECT
END SELECT
END SUBROUTINE connectCellEdge
SUBROUTINE connectQuadQuad(elemA, elemB)
IMPLICIT NONE
CLASS(meshCell2DCylQuad), INTENT(inout), TARGET:: elemA
CLASS(meshCell2DCylQuad), INTENT(inout), TARGET:: elemB
!Check direction 1
IF (.NOT. ASSOCIATED(elemA%e1) .AND. &
elemA%n1%n == elemB%n4%n .AND. &
elemA%n2%n == elemB%n3%n) THEN
elemA%e1 => elemB
elemB%e3 => elemA
END IF
!Check direction 2
IF (.NOT. ASSOCIATED(elemA%e2) .AND. &
elemA%n2%n == elemB%n1%n .AND. &
elemA%n3%n == elemB%n4%n) THEN
elemA%e2 => elemB
elemB%e4 => elemA
END IF
!Check direction 3
IF (.NOT. ASSOCIATED(elemA%e3) .AND. &
elemA%n3%n == elemB%n2%n .AND. &
elemA%n4%n == elemB%n1%n) THEN
elemA%e3 => elemB
elemB%e1 => elemA
END IF
!Check direction 4
IF (.NOT. ASSOCIATED(elemA%e4) .AND. &
elemA%n4%n == elemB%n3%n .AND. &
elemA%n1%n == elemB%n2%n) THEN
elemA%e4 => elemB
elemB%e2 => elemA
END IF
END SUBROUTINE connectQuadQuad
SUBROUTINE connectQuadTria(elemA, elemB)
IMPLICIT NONE
CLASS(meshCell2DCylQuad), INTENT(inout), TARGET:: elemA
CLASS(meshCell2DCylTria), INTENT(inout), TARGET:: elemB
!Check direction 1
IF (.NOT. ASSOCIATED(elemA%e1)) THEN
IF (elemA%n1%n == elemB%n1%n .AND. &
elemA%n2%n == elemB%n3%n) THEN
elemA%e1 => elemB
elemB%e3 => elemA
ELSEIF (elemA%n1%n == elemB%n3%n .AND. &
elemA%n2%n == elemB%n2%n) THEN
elemA%e1 => elemB
elemB%e2 => elemA
ELSEIF (elemA%n1%n == elemB%n2%n .AND. &
elemA%n2%n == elemB%n1%n) THEN
elemA%e1 => elemB
elemB%e1 => elemA
END IF
END IF
!Check direction 2
IF (.NOT. ASSOCIATED(elemA%e2)) THEN
IF (elemA%n2%n == elemB%n1%n .AND. &
elemA%n3%n == elemB%n3%n) THEN
elemA%e2 => elemB
elemB%e3 => elemA
ELSEIF (elemA%n2%n == elemB%n3%n .AND. &
elemA%n3%n == elemB%n2%n) THEN
elemA%e2 => elemB
elemB%e2 => elemA
ELSEIF (elemA%n2%n == elemB%n2%n .AND. &
elemA%n3%n == elemB%n1%n) THEN
elemA%e2 => elemB
elemB%e1 => elemA
END IF
END IF
!Check direction 3
IF (.NOT. ASSOCIATED(elemA%e3)) THEN
IF (elemA%n3%n == elemB%n1%n .AND. &
elemA%n4%n == elemB%n3%n) THEN
elemA%e3 => elemB
elemB%e3 => elemA
ELSEIF (elemA%n3%n == elemB%n3%n .AND. &
elemA%n4%n == elemB%n2%n) THEN
elemA%e3 => elemB
elemB%e2 => elemA
ELSEIF (elemA%n3%n == elemB%n2%n .AND. &
elemA%n4%n == elemB%n1%n) THEN
elemA%e3 => elemB
elemB%e1 => elemA
END IF
END IF
!Check direction 4
IF (.NOT. ASSOCIATED(elemA%e4)) THEN
IF (elemA%n4%n == elemB%n1%n .AND. &
elemA%n1%n == elemB%n3%n) THEN
elemA%e4 => elemB
elemB%e3 => elemA
ELSEIF (elemA%n4%n == elemB%n3%n .AND. &
elemA%n1%n == elemB%n2%n) THEN
elemA%e4 => elemB
elemB%e2 => elemA
ELSEIF (elemA%n4%n == elemB%n2%n .AND. &
elemA%n1%n == elemB%n1%n) THEN
elemA%e4 => elemB
elemB%e1 => elemA
END IF
END IF
END SUBROUTINE connectQuadTria
SUBROUTINE connectTriaTria(elemA, elemB)
IMPLICIT NONE
CLASS(meshCell2DCylTria), INTENT(inout), TARGET:: elemA
CLASS(meshCell2DCylTria), INTENT(inout), TARGET:: elemB
!Check direction 1
IF (.NOT. ASSOCIATED(elemA%e1)) THEN
IF (elemA%n1%n == elemB%n1%n .AND. &
elemA%n2%n == elemB%n3%n) THEN
elemA%e1 => elemB
elemB%e3 => elemA
ELSEIF (elemA%n1%n == elemB%n2%n .AND. &
elemA%n2%n == elemB%n1%n) THEN
elemA%e1 => elemB
elemB%e1 => elemA
ELSEIF (elemA%n1%n == elemB%n3%n .AND. &
elemA%n2%n == elemB%n2%n) THEN
elemA%e1 => elemB
elemB%e2 => elemA
END IF
END IF
!Check direction 2
IF (.NOT. ASSOCIATED(elemA%e2)) THEN
IF (elemA%n2%n == elemB%n1%n .AND. &
elemA%n3%n == elemB%n3%n) THEN
elemA%e2 => elemB
elemB%e3 => elemA
ELSEIF (elemA%n2%n == elemB%n2%n .AND. &
elemA%n3%n == elemB%n1%n) THEN
elemA%e2 => elemB
elemB%e1 => elemA
ELSEIF (elemA%n2%n == elemB%n3%n .AND. &
elemA%n3%n == elemB%n2%n) THEN
elemA%e2 => elemB
elemB%e2 => elemA
END IF
END IF
!Check direction 3
IF (.NOT. ASSOCIATED(elemA%e3)) THEN
IF (elemA%n3%n == elemB%n1%n .AND. &
elemA%n1%n == elemB%n3%n) THEN
elemA%e3 => elemB
elemB%e3 => elemA
ELSEIF (elemA%n3%n == elemB%n2%n .AND. &
elemA%n1%n == elemB%n1%n) THEN
elemA%e3 => elemB
elemB%e1 => elemA
ELSEIF (elemA%n3%n == elemB%n3%n .AND. &
elemA%n1%n == elemB%n2%n) THEN
elemA%e3 => elemB
elemB%e2 => elemA
END IF
END IF
END SUBROUTINE connectTriaTria
SUBROUTINE connectQuadEdge(elemA, elemB)
IMPLICIT NONE
CLASS(meshCell2DCylQuad), INTENT(inout), TARGET:: elemA
CLASS(meshEdge2DCyl), INTENT(inout), TARGET:: elemB
!Check direction 1
IF (.NOT. ASSOCIATED(elemA%e1)) THEN
IF (elemA%n1%n == elemB%n1%n .AND. &
elemA%n2%n == elemB%n2%n) THEN
elemA%e1 => elemB
elemB%e1 => elemA
ELSEIF (elemA%n1%n == elemB%n2%n .AND. &
elemA%n2%n == elemB%n1%n) THEN
elemA%e1 => elemB
elemB%e2 => elemA
!Revers the normal to point inside the domain
elemB%normal = - elemB%normal
END IF
END IF
!Check direction 2
IF (.NOT. ASSOCIATED(elemA%e2)) THEN
IF (elemA%n2%n == elemB%n1%n .AND. &
elemA%n3%n == elemB%n2%n) THEN
elemA%e2 => elemB
elemB%e1 => elemA
ELSEIF (elemA%n2%n == elemB%n2%n .AND. &
elemA%n3%n == elemB%n1%n) THEN
elemA%e2 => elemB
elemB%e2 => elemA
!Revers the normal to point inside the domain
elemB%normal = - elemB%normal
END IF
END IF
!Check direction 3
IF (.NOT. ASSOCIATED(elemA%e3)) THEN
IF (elemA%n3%n == elemB%n1%n .AND. &
elemA%n4%n == elemB%n2%n) THEN
elemA%e3 => elemB
elemB%e1 => elemA
ELSEIF (elemA%n3%n == elemB%n2%n .AND. &
elemA%n4%n == elemB%n1%n) THEN
elemA%e3 => elemB
elemB%e2 => elemA
!Revers the normal to point inside the domain
elemB%normal = - elemB%normal
END IF
END IF
!Check direction 4
IF (.NOT. ASSOCIATED(elemA%e4)) THEN
IF (elemA%n4%n == elemB%n1%n .AND. &
elemA%n1%n == elemB%n2%n) THEN
elemA%e4 => elemB
elemB%e1 => elemA
ELSEIF (elemA%n4%n == elemB%n2%n .AND. &
elemA%n1%n == elemB%n1%n) THEN
elemA%e4 => elemB
elemB%e2 => elemA
!Revers the normal to point inside the domain
elemB%normal = - elemB%normal
END IF
END IF
END SUBROUTINE connectQuadEdge
SUBROUTINE connectTriaEdge(elemA, elemB)
IMPLICIT NONE
CLASS(meshCell2DCylTria), INTENT(inout), TARGET:: elemA
CLASS(meshEdge2DCyl), INTENT(inout), TARGET:: elemB
!Check direction 1
IF (.NOT. ASSOCIATED(elemA%e1)) THEN
IF (elemA%n1%n == elemB%n1%n .AND. &
elemA%n2%n == elemB%n2%n) THEN
elemA%e1 => elemB
elemB%e1 => elemA
ELSEIF (elemA%n1%n == elemB%n2%n .AND. &
elemA%n2%n == elemB%n1%n) THEN
elemA%e1 => elemB
elemB%e2 => elemA
!Revers the normal to point inside the domain
elemB%normal = - elemB%normal
END IF
END IF
!Check direction 2
IF (.NOT. ASSOCIATED(elemA%e2)) THEN
IF (elemA%n2%n == elemB%n1%n .AND. &
elemA%n3%n == elemB%n2%n) THEN
elemA%e2 => elemB
elemB%e1 => elemA
ELSEIF (elemA%n2%n == elemB%n2%n .AND. &
elemA%n3%n == elemB%n1%n) THEN
elemA%e2 => elemB
elemB%e2 => elemA
!Revers the normal to point inside the domain
elemB%normal = - elemB%normal
END IF
END IF
!Check direction 3
IF (.NOT. ASSOCIATED(elemA%e3)) THEN
IF (elemA%n3%n == elemB%n1%n .AND. &
elemA%n1%n == elemB%n2%n) THEN
elemA%e3 => elemB
elemB%e1 => elemA
ELSEIF (elemA%n3%n == elemB%n2%n .AND. &
elemA%n1%n == elemB%n1%n) THEN
elemA%e3 => elemB
elemB%e2 => elemA
!Revers the normal to point inside the domain
elemB%normal = - elemB%normal
END IF
END IF
END SUBROUTINE connectTriaEdge
END MODULE moduleMesh2DCyl