Reduction in pushing

Reduction in 10-20% of time spend in pushing in 2DCyl thanks to
rewriting fPsi and dPsi.
This commit is contained in:
Jorge Gonzalez 2023-01-05 16:47:13 +01:00
commit 2486ef6316
18 changed files with 1289 additions and 1280 deletions

View file

@ -37,26 +37,19 @@ MODULE moduleMesh2DCart
END TYPE meshEdge2DCart
TYPE, PUBLIC, ABSTRACT, EXTENDS(meshVol):: meshVol2DCart
TYPE, PUBLIC, ABSTRACT, EXTENDS(meshCell):: meshCell2DCart
CONTAINS
PROCEDURE, PASS:: detJac => detJ2DCart
PROCEDURE, PASS:: invJac => invJ2DCart
PROCEDURE(dPsi_interface), DEFERRED, NOPASS:: dPsi
PROCEDURE(partialDer_interface), DEFERRED, PASS:: partialDer
PROCEDURE(partialDer_interface), DEFERRED, PASS, PRIVATE:: partialDer
END TYPE meshVol2DCart
END TYPE meshCell2DCart
ABSTRACT INTERFACE
PURE FUNCTION dPsi_interface(Xi) RESULT(dPsi)
REAL(8), INTENT(in):: Xi(1:3)
REAL(8), ALLOCATABLE:: dPsi(:,:)
END FUNCTION dPsi_interface
PURE SUBROUTINE partialDer_interface(self, dPsi, dx, dy)
IMPORT meshVol2DCart
CLASS(meshVol2DCart), INTENT(in):: self
REAL(8), INTENT(in):: dPsi(1:,1:)
IMPORT meshCell2DCart
CLASS(meshCell2DCart), INTENT(in):: self
REAL(8), INTENT(in):: dPsi(1:3,1:self%nNodes)
REAL(8), INTENT(out), DIMENSION(1:2):: dx, dy
END SUBROUTINE partialDer_interface
@ -64,7 +57,7 @@ MODULE moduleMesh2DCart
END INTERFACE
!Quadrilateral volume element
TYPE, PUBLIC, EXTENDS(meshVol2DCart):: meshVol2DCartQuad
TYPE, PUBLIC, EXTENDS(meshCell2DCart):: meshCell2DCartQuad
!Element coordinates
REAL(8):: x(1:4) = 0.D0, y(1:4) = 0.D0
!Connectivity to nodes
@ -73,27 +66,27 @@ MODULE moduleMesh2DCart
CLASS(meshElement), POINTER:: e1 => NULL(), e2 => NULL(), e3 => NULL(), e4 => NULL()
REAL(8):: arNodes(1:4) = 0.D0
CONTAINS
PROCEDURE, PASS:: init => initVolQuad2DCart
PROCEDURE, PASS:: randPos => randPosVolQuad
PROCEDURE, PASS:: init => initCellQuad2DCart
PROCEDURE, PASS:: randPos => randPosCellQuad
PROCEDURE, PASS:: area => areaQuad
PROCEDURE, NOPASS:: fPsi => fPsiQuad
PROCEDURE, NOPASS:: dPsi => dPsiQuad
PROCEDURE, NOPASS:: dPsiXi1 => dPsiQuadXi1
PROCEDURE, NOPASS:: dPsiXi2 => dPsiQuadXi2
PROCEDURE, PASS:: partialDer => partialDerQuad
PROCEDURE, PASS:: fPsi => fPsiQuad
PROCEDURE, PASS:: dPsi => dPsiQuad
PROCEDURE, NOPASS, PRIVATE:: dPsiXi1 => dPsiQuadXi1
PROCEDURE, NOPASS, PRIVATE:: dPsiXi2 => dPsiQuadXi2
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:: gatherEF => gatherEFQuad
PROCEDURE, PASS:: gatherMF => gatherMFQuad
PROCEDURE, PASS:: getNodes => getNodesQuad
PROCEDURE, PASS:: phy2log => phy2logQuad
PROCEDURE, PASS:: nextElement => nextElementQuad
END TYPE meshVol2DCartQuad
END TYPE meshCell2DCartQuad
!Triangular volume element
TYPE, PUBLIC, EXTENDS(meshVol2DCart):: meshVol2DCartTria
TYPE, PUBLIC, EXTENDS(meshCell2DCart):: meshCell2DCartTria
!Element coordinates
REAL(8):: x(1:3) = 0.D0, y(1:3) = 0.D0
!Connectivity to nodes
@ -103,24 +96,24 @@ MODULE moduleMesh2DCart
REAL(8):: arNodes(1:3) = 0.D0
CONTAINS
PROCEDURE, PASS:: init => initVolTria2DCart
PROCEDURE, PASS:: randPos => randPosVolTria
PROCEDURE, PASS:: init => initCellTria2DCart
PROCEDURE, PASS:: randPos => randPosCellTria
PROCEDURE, PASS:: area => areaTria
PROCEDURE, NOPASS:: fPsi => fPsiTria
PROCEDURE, NOPASS:: dPsi => dPsiTria
PROCEDURE, PASS:: fPsi => fPsiTria
PROCEDURE, PASS:: dPsi => dPsiTria
PROCEDURE, NOPASS:: dPsiXi1 => dPsiTriaXi1
PROCEDURE, NOPASS:: dPsiXi2 => dPsiTriaXi2
PROCEDURE, PASS:: partialDer => partialDerTria
PROCEDURE, PASS:: elemK => elemKTria
PROCEDURE, PASS:: elemF => elemFTria
PROCEDURE, PASS:: gatherElectricField => gatherEFTria
PROCEDURE, PASS:: gatherMagneticField => gatherMFTria
PROCEDURE, NOPASS:: inside => insideTria
PROCEDURE, PASS:: gatherEF => gatherEFTria
PROCEDURE, PASS:: gatherMF => gatherMFTria
PROCEDURE, PASS:: getNodes => getNodesTria
PROCEDURE, PASS:: phy2log => phy2logTria
PROCEDURE, PASS:: nextElement => nextElementTria
END TYPE meshVol2DCartTria
END TYPE meshCell2DCartTria
CONTAINS
!NODE FUNCTIONS
@ -204,26 +197,26 @@ MODULE moduleMesh2DCart
END SUBROUTINE initEdge2DCart
!Random position in quadrilateral volume
FUNCTION randPosVolQuad(self) RESULT(r)
FUNCTION randPosCellQuad(self) RESULT(r)
USE moduleRandom
IMPLICIT NONE
CLASS(meshVol2DCartQuad), INTENT(in):: self
CLASS(meshCell2DCartQuad), INTENT(in):: self
REAL(8):: r(1:3)
REAL(8):: Xi(1:3)
REAL(8), ALLOCATABLE:: fPsi(:)
REAL(8):: fPsi(1:4)
Xi(1) = random(-1.D0, 1.D0)
Xi(2) = random(-1.D0, 1.D0)
Xi(3) = 0.D0
CALL self%fPsi(Xi, fPsi)
fPsi = self%fPsi(Xi)
r(1) = DOT_PRODUCT(fPsi, self%x)
r(2) = DOT_PRODUCT(fPsi, self%y)
r(3) = 0.D0
END FUNCTION randposVolQuad
END FUNCTION randposCellQuad
!Get nodes from edge
PURE FUNCTION getNodes2DCart(self) RESULT(n)
@ -232,7 +225,6 @@ MODULE moduleMesh2DCart
CLASS(meshEdge2DCart), INTENT(in):: self
INTEGER, ALLOCATABLE:: n(:)
ALLOCATE(n(1:2))
n = (/self%n1%n, self%n2%n /)
END FUNCTION getNodes2DCart
@ -277,17 +269,18 @@ MODULE moduleMesh2DCart
!VOLUME FUNCTIONS
!QUAD FUNCTIONS
!Inits quadrilateral element
SUBROUTINE initVolQuad2DCart(self, n, p, nodes)
SUBROUTINE initCellQuad2DCart(self, n, p, nodes)
USE moduleRefParam
IMPLICIT NONE
CLASS(meshVol2DCartQuad), INTENT(out):: self
CLASS(meshCell2DCartQuad), 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
@ -312,13 +305,13 @@ MODULE moduleMesh2DCart
ALLOCATE(self%listPart_in(1:nSpecies))
ALLOCATE(self%totalWeight(1:nSpecies))
END SUBROUTINE initVolQuad2DCart
END SUBROUTINE initCellQuad2DCart
!Computes element area
PURE SUBROUTINE areaQuad(self)
IMPLICIT NONE
CLASS(meshVol2DCartQuad), INTENT(inout):: self
CLASS(meshCell2DCartQuad), INTENT(inout):: self
REAL(8):: Xi(1:3)
REAL(8):: detJ
REAL(8):: fPsi(1:4)
@ -328,18 +321,19 @@ MODULE moduleMesh2DCart
!2D 1 point Gauss Quad Integral
Xi = 0.D0
detJ = self%detJac(Xi)*4.D0 !4
CALL self%fPsi(Xi, fPsi)
fPsi = self%fPsi(Xi)
self%volume = detJ
self%arNodes = fPsi*detJ
END SUBROUTINE areaQuad
!Computes element functions in point Xi
PURE SUBROUTINE fPsiQuad(Xi, fPsi)
PURE FUNCTION fPsiQuad(self, Xi) RESULT(fPsi)
IMPLICIT NONE
CLASS(meshCell2DCartQuad), INTENT(in):: self
REAL(8), INTENT(in):: Xi(1:3)
REAL(8), INTENT(out):: fPsi(:)
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))
@ -348,16 +342,17 @@ MODULE moduleMesh2DCart
fPsi = fPsi*0.25D0
END SUBROUTINE fPsiQuad
END FUNCTION fPsiQuad
!Derivative element function at coordinates Xi
PURE FUNCTION dPsiQuad(Xi) RESULT(dPsi)
PURE FUNCTION dPsiQuad(self, Xi) RESULT(dPsi)
IMPLICIT NONE
CLASS(meshCell2DCartQuad), INTENT(in):: self
REAL(8), INTENT(in):: Xi(1:3)
REAL(8), ALLOCATABLE:: dPsi(:,:)
REAL(8):: dPsi(1:3,1:self%nNodes)
ALLOCATE(dPsi(1:2,1:4))
dPsi = 0.D0
dPsi(1,:) = dPsiQuadXi1(Xi(2))
dPsi(2,:) = dPsiQuadXi2(Xi(1))
@ -397,8 +392,8 @@ MODULE moduleMesh2DCart
PURE SUBROUTINE partialDerQuad(self, dPsi, dx, dy)
IMPLICIT NONE
CLASS(meshVol2DCartQuad), INTENT(in):: self
REAL(8), INTENT(in):: dPsi(1:,1:)
CLASS(meshCell2DCartQuad), INTENT(in):: self
REAL(8), INTENT(in):: dPsi(1:3,1:self%nNodes)
REAL(8), INTENT(out), DIMENSION(1:2):: dx, dy
dx(1) = DOT_PRODUCT(dPsi(1,:),self%x)
@ -412,14 +407,13 @@ MODULE moduleMesh2DCart
PURE FUNCTION elemKQuad(self) RESULT(localK)
IMPLICIT NONE
CLASS(meshVol2DCartQuad), INTENT(in):: self
REAL(8), ALLOCATABLE:: localK(:,:)
CLASS(meshCell2DCartQuad), INTENT(in):: self
REAL(8):: localK(1:self%nNodes,1:self%nNodes)
REAL(8):: Xi(1:3)
REAL(8):: fPsi(1:4), dPsi(1:2,1:4)
REAL(8):: invJ(1:2,1:2), detJ
REAL(8):: fPsi(1:4), dPsi(1:3,1:4)
REAL(8):: invJ(1:3,1:3), detJ
INTEGER:: l, m
ALLOCATE(localK(1:4, 1:4))
localK=0.D0
Xi=0.D0
!Start 2D Gauss Quad Integral
@ -429,7 +423,7 @@ MODULE moduleMesh2DCart
DO m = 1, 3
Xi(1) = corQuad(m)
dPsi(2,:) = self%dPsiXi2(Xi(1))
CALL self%fPsi(Xi, fPsi)
fPsi = self%fPsi(Xi)
detJ = self%detJac(Xi,dPsi)
invJ = self%invJac(Xi,dPsi)
localK = localK + MATMUL(TRANSPOSE(MATMUL(invJ,dPsi)),MATMUL(invJ,dPsi))*wQuad(l)*wQuad(m)/detJ
@ -443,15 +437,14 @@ MODULE moduleMesh2DCart
PURE FUNCTION elemFQuad(self, source) RESULT(localF)
IMPLICIT NONE
CLASS(meshVol2DCartQuad), INTENT(in):: self
REAL(8), INTENT(in):: source(1:)
REAL(8), ALLOCATABLE:: localF(:)
CLASS(meshCell2DCartQuad), 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):: detJ, f
INTEGER:: l, m
ALLOCATE(localF(1:4))
localF = 0.D0
Xi = 0.D0
DO l=1, 3
@ -459,7 +452,7 @@ MODULE moduleMesh2DCart
DO m = 1, 3
Xi(2) = corQuad(m)
detJ = self%detJac(Xi)
CALL self%fPsi(Xi, fPsi)
fPsi = self%fPsi(Xi)
f = DOT_PRODUCT(fPsi,source)
localF = localF + f*fPsi*wQuad(l)*wQuad(m)*detJ
@ -468,6 +461,48 @@ MODULE moduleMesh2DCart
END FUNCTION elemFQuad
PURE FUNCTION gatherEFQuad(self, Xi) RESULT(array)
IMPLICIT NONE
CLASS(meshCell2DCartQuad), 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(meshCell2DCartQuad), 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
@ -480,97 +515,42 @@ MODULE moduleMesh2DCart
END FUNCTION insideQuad
!Gathers the electric field at position Xi
PURE FUNCTION gatherEFQuad(self,Xi) RESULT(EF)
IMPLICIT NONE
CLASS(meshVol2DCartQuad), INTENT(in):: self
REAL(8), INTENT(in):: Xi(1:3)
REAL(8):: dPsi(1:2,1:4)
REAL(8):: dPsiR(1:2,1:4)!Derivative of shpae functions in global coordinates
REAL(8):: invJ(1:2,1:2), detJ
REAL(8):: phi(1:4)
REAL(8):: EF(1:3)
phi = (/self%n1%emData%phi, &
self%n2%emData%phi, &
self%n3%emData%phi, &
self%n4%emData%phi /)
dPsi = self%dPsi(Xi)
detJ = self%detJac(Xi,dPsi)
invJ = self%invJac(Xi,dPsi)
dPsiR = MATMUL(invJ, dPsi)/detJ
EF(1) = -DOT_PRODUCT(dPsiR(1,:), phi)
EF(2) = -DOT_PRODUCT(dPsiR(2,:), phi)
EF(3) = 0.D0
END FUNCTION gatherEFQuad
PURE FUNCTION gatherMFQuad(self,Xi) RESULT(MF)
IMPLICIT NONE
CLASS(meshVol2DCartQuad), INTENT(in):: self
REAL(8), INTENT(in):: Xi(1:3)
REAL(8):: fPsi(1:4)
REAL(8):: MF_Nodes(1:4,1:3)
REAL(8):: MF(1:3)
MF_Nodes(1:4,1) = (/self%n1%emData%B(1), &
self%n2%emData%B(1), &
self%n3%emData%B(1), &
self%n4%emData%B(1) /)
MF_Nodes(1:4,2) = (/self%n1%emData%B(2), &
self%n2%emData%B(2), &
self%n3%emData%B(2), &
self%n4%emData%B(2) /)
MF_Nodes(1:4,3) = (/self%n1%emData%B(3), &
self%n2%emData%B(3), &
self%n3%emData%B(3), &
self%n4%emData%B(3) /)
CALL self%fPsi(Xi, fPsi)
MF = MATMUL(fPsi(:), MF_Nodes)
END FUNCTION gatherMFQuad
!Gets nodes from quadrilateral element
PURE FUNCTION getNodesQuad(self) RESULT(n)
IMPLICIT NONE
CLASS(meshVol2DCartQuad), INTENT(in):: self
INTEGER, ALLOCATABLE:: n(:)
CLASS(meshCell2DCartQuad), INTENT(in):: self
INTEGER:: n(1:self%nNodes)
ALLOCATE(n(1:4))
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(XiN)
PURE FUNCTION phy2logQuad(self,r) RESULT(Xi)
IMPLICIT NONE
CLASS(meshVol2DCartQuad), INTENT(in):: self
CLASS(meshCell2DCartQuad), INTENT(in):: self
REAL(8), INTENT(in):: r(1:3)
REAL(8):: XiN(1:3)
REAL(8):: Xi(1:3)
REAL(8):: XiO(1:3), detJ, invJ(1:2,1:2), f(1:2)
REAL(8):: dPsi(1:2,1:4), fPsi(1:4)
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-4)
DO WHILE(conv>1.D-3)
dPsi = self%dPsi(XiO)
invJ = self%invJac(XiO, dPsi)
CALL self%fPsi(XiO, fPsi)
fPsi = self%fPsi(XiO)
f(1) = DOT_PRODUCT(fPsi,self%x)-r(1)
f(2) = DOT_PRODUCT(fPsi,self%y)-r(2)
detJ = self%detJac(XiO,dPsi)
XiN(1:2)=XiO(1:2) - MATMUL(invJ, f)/detJ
conv=MAXVAL(DABS(XiN-XiO),1)
XiO=XiN
Xi(1:2)=XiO(1:2) - MATMUL(invJ, f)/detJ
conv=MAXVAL(DABS(Xi-XiO),1)
XiO=Xi
END DO
@ -580,7 +560,7 @@ MODULE moduleMesh2DCart
SUBROUTINE nextElementQuad(self, Xi, nextElement)
IMPLICIT NONE
CLASS(meshVol2DCartQuad), INTENT(in):: self
CLASS(meshCell2DCartQuad), INTENT(in):: self
REAL(8), INTENT(in):: Xi(1:3)
CLASS(meshElement), POINTER, INTENT(out):: nextElement
REAL(8):: XiArray(1:4)
@ -605,11 +585,11 @@ MODULE moduleMesh2DCart
!TRIA ELEMENT
!Init tria element
SUBROUTINE initVolTria2DCart(self, n, p, nodes)
SUBROUTINE initCellTria2DCart(self, n, p, nodes)
USE moduleRefParam
IMPLICIT NONE
CLASS(meshVol2DCartTria), INTENT(out):: self
CLASS(meshCell2DCartTria), INTENT(out):: self
INTEGER, INTENT(in):: n
INTEGER, INTENT(in):: p(:)
TYPE(meshNodeCont), INTENT(in), TARGET:: nodes(:)
@ -618,6 +598,9 @@ MODULE moduleMesh2DCart
!Assign node index
self%n = n
!Assign number of nodes of cell
self%nNodes = SIZE(p)
!Assign nodes to element
self%n1 => nodes(p(1))%obj
self%n2 => nodes(p(2))%obj
@ -639,14 +622,14 @@ MODULE moduleMesh2DCart
ALLOCATE(self%listPart_in(1:nSpecies))
ALLOCATE(self%totalWeight(1:nSpecies))
END SUBROUTINE initVolTria2DCart
END SUBROUTINE initCellTria2DCart
!Random position in quadrilateral volume
FUNCTION randPosVolTria(self) RESULT(r)
FUNCTION randPosCellTria(self) RESULT(r)
USE moduleRandom
IMPLICIT NONE
CLASS(meshVol2DCartTria), INTENT(in):: self
CLASS(meshCell2DCartTria), INTENT(in):: self
REAL(8):: r(1:3)
REAL(8):: Xi(1:3)
REAL(8):: fPsi(1:3)
@ -655,19 +638,19 @@ MODULE moduleMesh2DCart
Xi(2) = random( 0.D0, 1.D0 - Xi(1))
Xi(3) = 0.D0
CALL self%fPsi(Xi, fPsi)
fPsi = self%fPsi(Xi)
r(1) = DOT_PRODUCT(fPsi, self%x)
r(2) = DOT_PRODUCT(fPsi, self%y)
r(3) = 0.D0
END FUNCTION randposVolTria
END FUNCTION randposCellTria
!Calculates area for triangular element
PURE SUBROUTINE areaTria(self)
IMPLICIT NONE
CLASS(meshVol2DCartTria), INTENT(inout):: self
CLASS(meshCell2DCartTria), INTENT(inout):: self
REAL(8):: Xi(1:3)
REAL(8):: detJ
REAL(8):: fPsi(1:3)
@ -677,33 +660,35 @@ MODULE moduleMesh2DCart
!2D 1 point Gauss Quad Integral
Xi = (/1.D0/3.D0, 1.D0/3.D0, 0.D0 /)
detJ = self%detJac(Xi)/2.D0
CALL self%fPsi(Xi, fPsi)
fPsi = self%fPsi(Xi)
self%volume = detJ
self%arNodes = fPsi*detJ
END SUBROUTINE areaTria
!Shape functions for triangular element
PURE SUBROUTINE fPsiTria(Xi, fPsi)
PURE FUNCTION fPsiTria(self, Xi) RESULT(fPsi)
IMPLICIT NONE
CLASS(meshCell2DCartTria), INTENT(in):: self
REAL(8), INTENT(in):: Xi(1:3)
REAL(8), INTENT(out):: fPsi(:)
REAL(8):: fPsi(1:self%nNodes)
fPsi(1) = 1.D0 - Xi(1) - Xi(2)
fPsi(2) = Xi(1)
fPsi(3) = Xi(2)
END SUBROUTINE fPsiTria
END FUNCTION fPsiTria
!Derivative element function at coordinates Xi
PURE FUNCTION dPsiTria(Xi) RESULT(dPsi)
PURE FUNCTION dPsiTria(self, Xi) RESULT(dPsi)
IMPLICIT NONE
CLASS(meshCell2DCartTria), INTENT(in):: self
REAL(8), INTENT(in):: Xi(1:3)
REAL(8), ALLOCATABLE:: dPsi(:,:)
REAL(8):: dPsi(1:3,1:self%nNodes)
ALLOCATE(dPsi(1:2,1:3))
dPsi = 0.D0
dPsi(1,:) = dPsiTriaXi1(Xi(2))
dPsi(2,:) = dPsiTriaXi2(Xi(1))
@ -739,8 +724,8 @@ MODULE moduleMesh2DCart
PURE SUBROUTINE partialDerTria(self, dPsi, dx, dy)
IMPLICIT NONE
CLASS(meshVol2DCartTria), INTENT(in):: self
REAL(8), INTENT(in):: dPsi(1:,1:)
CLASS(meshCell2DCartTria), INTENT(in):: self
REAL(8), INTENT(in):: dPsi(1:3,1:self%nNodes)
REAL(8), INTENT(out), DIMENSION(1:2):: dx, dy
dx(1) = DOT_PRODUCT(dPsi(1,:),self%x)
@ -754,14 +739,13 @@ MODULE moduleMesh2DCart
PURE FUNCTION elemKTria(self) RESULT(localK)
IMPLICIT NONE
CLASS(meshVol2DCartTria), INTENT(in):: self
REAL(8), ALLOCATABLE:: localK(:,:)
CLASS(meshCell2DCartTria), INTENT(in):: self
REAL(8):: localK(1:self%nNodes,1:self%nNodes)
REAL(8):: Xi(1:3)
REAL(8):: fPsi(1:3), dPsi(1:2,1:3)
REAL(8):: invJ(1:2,1:2), detJ
REAL(8):: fPsi(1:3), dPsi(1:3,1:3)
REAL(8):: invJ(1:3,1:3), detJ
INTEGER:: l
ALLOCATE(localK(1:4, 1:4))
localK=0.D0
Xi=0.D0
!Start 2D Gauss Quad Integral
@ -771,7 +755,7 @@ MODULE moduleMesh2DCart
dPsi = self%dPsi(Xi)
detJ = self%detJac(Xi,dPsi)
invJ = self%invJac(Xi,dPsi)
CALL self%fPsi(Xi, fPsi)
fPsi = self%fPsi(Xi)
localK = localK + MATMUL(TRANSPOSE(MATMUL(invJ,dPsi)),MATMUL(invJ,dPsi))*wTria(l)/detJ
END DO
@ -782,15 +766,14 @@ MODULE moduleMesh2DCart
PURE FUNCTION elemFTria(self, source) RESULT(localF)
IMPLICIT NONE
CLASS(meshVol2DCartTria), INTENT(in):: self
REAL(8), INTENT(in):: source(1:)
REAL(8), ALLOCATABLE:: localF(:)
CLASS(meshCell2DCartTria), 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):: detJ, f
INTEGER:: l
ALLOCATE(localF(1:3))
localF = 0.D0
Xi = 0.D0
!Start 2D Gauss Quad Integral
@ -798,7 +781,7 @@ MODULE moduleMesh2DCart
Xi(1) = Xi1Tria(l)
Xi(2) = Xi2Tria(l)
detJ = self%detJac(Xi)
CALL self%fPsi(Xi, fPsi)
fPsi = self%fPsi(Xi)
f = DOT_PRODUCT(fPsi,source)
localF = localF + f*fPsi*wTria(l)*detJ
@ -806,6 +789,44 @@ MODULE moduleMesh2DCart
END FUNCTION elemFTria
PURE FUNCTION gatherEFTria(self, Xi) RESULT(array)
IMPLICIT NONE
CLASS(meshCell2DCartTria), 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(meshCell2DCartTria), 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
@ -818,64 +839,13 @@ MODULE moduleMesh2DCart
END FUNCTION insideTria
!Gathers the electric field at position Xi
PURE FUNCTION gatherEFTria(self,Xi) RESULT(EF)
IMPLICIT NONE
CLASS(meshVol2DCartTria), INTENT(in):: self
REAL(8), INTENT(in):: Xi(1:3)
REAL(8):: dPsi(1:2,1:3)
REAL(8):: dPsiR(1:2,1:3)!Derivative of shpae functions in global coordinates
REAL(8):: invJ(1:2,1:2), detJ
REAL(8):: phi(1:3)
REAL(8):: EF(1:3)
phi = (/self%n1%emData%phi, &
self%n2%emData%phi, &
self%n3%emData%phi /)
dPsi = self%dPsi(Xi)
detJ = self%detJac(Xi,dPsi)
invJ = self%invJac(Xi,dPsi)
dPsiR = MATMUL(invJ, dPsi)/detJ
EF(1) = -DOT_PRODUCT(dPsiR(1,:), phi)
EF(2) = -DOT_PRODUCT(dPsiR(2,:), phi)
EF(3) = 0.D0
END FUNCTION gatherEFTria
PURE FUNCTION gatherMFTria(self,Xi) RESULT(MF)
IMPLICIT NONE
CLASS(meshVol2DCartTria), INTENT(in):: self
REAL(8), INTENT(in):: Xi(1:3)
REAL(8):: fPsi(1:3)
REAL(8):: MF_Nodes(1:3,1:3)
REAL(8):: MF(1:3)
MF_Nodes(1:3,1) = (/self%n1%emData%B(1), &
self%n2%emData%B(1), &
self%n3%emData%B(1) /)
MF_Nodes(1:3,2) = (/self%n1%emData%B(2), &
self%n2%emData%B(2), &
self%n3%emData%B(2) /)
MF_Nodes(1:3,3) = (/self%n1%emData%B(3), &
self%n2%emData%B(3), &
self%n3%emData%B(3) /)
CALL self%fPsi(Xi, fPsi)
MF = MATMUL(fPsi, MF_Nodes)
END FUNCTION gatherMFTria
!Gets node indexes from triangular element
PURE FUNCTION getNodesTria(self) RESULT(n)
IMPLICIT NONE
CLASS(meshVol2DCartTria), INTENT(in):: self
INTEGER, ALLOCATABLE:: n(:)
CLASS(meshCell2DCartTria), INTENT(in):: self
INTEGER:: n(1:self%nNodes)
ALLOCATE(n(1:3))
n = (/self%n1%n, self%n2%n, self%n3%n /)
END FUNCTION getNodesTria
@ -884,27 +854,27 @@ MODULE moduleMesh2DCart
PURE FUNCTION phy2logTria(self,r) RESULT(Xi)
IMPLICIT NONE
CLASS(meshVol2DCartTria), INTENT(in):: self
CLASS(meshCell2DCartTria), INTENT(in):: self
REAL(8), INTENT(in):: r(1:3)
REAL(8):: Xi(1:3)
REAL(8):: invJ(1:2,1:2), detJ
REAL(8):: deltaR(1:2)
REAL(8):: dPsi(1:2,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 !Irrelevant, required for input
deltaR = (/ r(1) - self%x(1), r(2) - self%y(1) /)
dPsi = self%dPsi(Xi)
invJ = self%invJac(Xi, dPsi)
detJ = self%detJac(Xi, dPsi)
Xi(1:2) = MATMUL(invJ,deltaR)/detJ
Xi = 0.D0
deltaR = (/ r(1) - self%x(1), r(2) - self%y(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(meshVol2DCartTria), INTENT(in):: self
CLASS(meshCell2DCartTria), INTENT(in):: self
REAL(8), INTENT(in):: Xi(1:3)
CLASS(meshElement), POINTER, INTENT(out):: nextElement
REAL(8):: XiArray(1:3)
@ -929,10 +899,10 @@ MODULE moduleMesh2DCart
PURE FUNCTION detJ2DCart(self, Xi, dPsi_in) RESULT(dJ)
IMPLICIT NONE
CLASS(meshVol2DCart), INTENT(in):: self
CLASS(meshCell2DCart), INTENT(in):: self
REAL(8), INTENT(in):: Xi(1:3)
REAL(8), INTENT(in), OPTIONAL:: dPsi_in(1:,1:)
REAL(8), ALLOCATABLE:: dPsi(:,:)
REAL(8), INTENT(in), OPTIONAL:: dPsi_in(1:3,1:self%nNodes)
REAL(8):: dPsi(1:3,1:self%nNodes)
REAL(8):: dJ
REAL(8):: dx(1:2), dy(1:2)
@ -953,12 +923,12 @@ MODULE moduleMesh2DCart
PURE FUNCTION invJ2DCart(self,Xi,dPsi_in) RESULT(invJ)
IMPLICIT NONE
CLASS(meshVol2DCart), INTENT(in):: self
CLASS(meshCell2DCart), INTENT(in):: self
REAL(8), INTENT(in):: Xi(1:3)
REAL(8), INTENT(in), OPTIONAL:: dPsi_in(1:,1:)
REAL(8), ALLOCATABLE:: dPsi(:,:)
REAL(8), INTENT(in), OPTIONAL:: dPsi_in(1:3,1:self%nNodes)
REAL(8):: dPsi(1:3,1:self%nNodes)
REAL(8):: dx(1:2), dy(1:2)
REAL(8):: invJ(1:2,1:2)
REAL(8):: invJ(1:3,1:3)
IF(PRESENT(dPsi_in)) THEN
dPsi=dPsi_in
@ -968,9 +938,12 @@ MODULE moduleMesh2DCart
END IF
invJ = 0.D0
CALL self%partialDer(dPsi, dx, dy)
invJ(1,:) = (/ dy(2), -dx(2) /)
invJ(2,:) = (/ -dy(1), dx(1) /)
invJ(1,1:2) = (/ dy(2), -dx(2) /)
invJ(2,1:2) = (/ -dy(1), dx(1) /)
END FUNCTION invJ2DCart
@ -980,11 +953,11 @@ MODULE moduleMesh2DCart
CLASS(meshGeneric), INTENT(inout):: self
INTEGER:: e, et
DO e = 1, self%numVols
!Connect Vol-Vol
DO et = 1, self%numVols
DO e = 1, self%numCells
!Connect Cell-Cell
DO et = 1, self%numCells
IF (e /= et) THEN
CALL connectVolVol(self%vols(e)%obj, self%vols(et)%obj)
CALL connectCellCell(self%cells(e)%obj, self%cells(et)%obj)
END IF
@ -992,9 +965,9 @@ MODULE moduleMesh2DCart
SELECT TYPE(self)
TYPE IS(meshParticles)
!Connect Vol-Edge
!Connect Cell-Edge
DO et = 1, self%numEdges
CALL connectVolEdge(self%vols(e)%obj, self%edges(et)%obj)
CALL connectCellEdge(self%cells(e)%obj, self%edges(et)%obj)
END DO
@ -1005,34 +978,34 @@ MODULE moduleMesh2DCart
END SUBROUTINE connectMesh2DCart
!Selects type of elements to build connection
SUBROUTINE connectVolVol(elemA, elemB)
SUBROUTINE connectCellCell(elemA, elemB)
IMPLICIT NONE
CLASS(meshVol), INTENT(inout):: elemA
CLASS(meshVol), INTENT(inout):: elemB
CLASS(meshCell), INTENT(inout):: elemA
CLASS(meshCell), INTENT(inout):: elemB
SELECT TYPE(elemA)
TYPE IS(meshVol2DCartQuad)
TYPE IS(meshCell2DCartQuad)
!Element A is a quadrilateral
SELECT TYPE(elemB)
TYPE IS(meshVol2DCartQuad)
TYPE IS(meshCell2DCartQuad)
!Element B is a quadrilateral
CALL connectQuadQuad(elemA, elemB)
TYPE IS(meshVol2DCartTria)
TYPE IS(meshCell2DCartTria)
!Element B is a triangle
CALL connectQuadTria(elemA, elemB)
END SELECT
TYPE IS(meshVol2DCartTria)
TYPE IS(meshCell2DCartTria)
!Element A is a Triangle
SELECT TYPE(elemB)
TYPE IS(meshVol2DCartQuad)
TYPE IS(meshCell2DCartQuad)
!Element B is a quadrilateral
CALL connectQuadTria(elemB, elemA)
TYPE IS(meshVol2DCartTria)
TYPE IS(meshCell2DCartTria)
!Element B is a triangle
CALL connectTriaTria(elemA, elemB)
@ -1040,22 +1013,22 @@ MODULE moduleMesh2DCart
END SELECT
END SUBROUTINE connectVolVol
END SUBROUTINE connectCellCell
SUBROUTINE connectVolEdge(elemA, elemB)
SUBROUTINE connectCellEdge(elemA, elemB)
IMPLICIT NONE
CLASS(meshVol), INTENT(inout):: elemA
CLASS(meshCell), INTENT(inout):: elemA
CLASS(meshEdge), INTENT(inout):: elemB
SELECT TYPE(elemB)
CLASS IS(meshEdge2DCart)
SELECT TYPE(elemA)
TYPE IS(meshVol2DCartQuad)
TYPE IS(meshCell2DCartQuad)
!Element A is a quadrilateral
CALL connectQuadEdge(elemA, elemB)
TYPE IS(meshVol2DCartTria)
TYPE IS(meshCell2DCartTria)
!Element A is a triangle
CALL connectTriaEdge(elemA, elemB)
@ -1063,13 +1036,13 @@ MODULE moduleMesh2DCart
END SELECT
END SUBROUTINE connectVolEdge
END SUBROUTINE connectCellEdge
SUBROUTINE connectQuadQuad(elemA, elemB)
IMPLICIT NONE
CLASS(meshVol2DCartQuad), INTENT(inout), TARGET:: elemA
CLASS(meshVol2DCartQuad), INTENT(inout), TARGET:: elemB
CLASS(meshCell2DCartQuad), INTENT(inout), TARGET:: elemA
CLASS(meshCell2DCartQuad), INTENT(inout), TARGET:: elemB
!Check direction 1
IF (.NOT. ASSOCIATED(elemA%e1) .AND. &
@ -1112,8 +1085,8 @@ MODULE moduleMesh2DCart
SUBROUTINE connectQuadTria(elemA, elemB)
IMPLICIT NONE
CLASS(meshVol2DCartQuad), INTENT(inout), TARGET:: elemA
CLASS(meshVol2DCartTria), INTENT(inout), TARGET:: elemB
CLASS(meshCell2DCartQuad), INTENT(inout), TARGET:: elemA
CLASS(meshCell2DCartTria), INTENT(inout), TARGET:: elemB
!Check direction 1
IF (.NOT. ASSOCIATED(elemA%e1)) THEN
@ -1204,8 +1177,8 @@ MODULE moduleMesh2DCart
SUBROUTINE connectTriaTria(elemA, elemB)
IMPLICIT NONE
CLASS(meshVol2DCartTria), INTENT(inout), TARGET:: elemA
CLASS(meshVol2DCartTria), INTENT(inout), TARGET:: elemB
CLASS(meshCell2DCartTria), INTENT(inout), TARGET:: elemA
CLASS(meshCell2DCartTria), INTENT(inout), TARGET:: elemB
!Check direction 1
IF (.NOT. ASSOCIATED(elemA%e1)) THEN
@ -1277,7 +1250,7 @@ MODULE moduleMesh2DCart
SUBROUTINE connectQuadEdge(elemA, elemB)
IMPLICIT NONE
CLASS(meshVol2DCartQuad), INTENT(inout), TARGET:: elemA
CLASS(meshCell2DCartQuad), INTENT(inout), TARGET:: elemA
CLASS(meshEdge2DCart), INTENT(inout), TARGET:: elemB
!Check direction 1
@ -1361,7 +1334,7 @@ MODULE moduleMesh2DCart
SUBROUTINE connectTriaEdge(elemA, elemB)
IMPLICIT NONE
CLASS(meshVol2DCartTria), INTENT(inout), TARGET:: elemA
CLASS(meshCell2DCartTria), INTENT(inout), TARGET:: elemA
CLASS(meshEdge2DCart), INTENT(inout), TARGET:: elemB
!Check direction 1