2018-08-17 14:53:24 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @author Pratheek Shanthraj, Max-Planck-Institut für Eisenforschung GmbH
|
|
|
|
!> @author Martin Diehl, Max-Planck-Institut für Eisenforschung GmbH
|
|
|
|
!> @author Philip Eisenlohr, Max-Planck-Institut für Eisenforschung GmbH
|
|
|
|
!> @author Franz Roters, Max-Planck-Institut für Eisenforschung GmbH
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2020-03-20 19:38:07 +05:30
|
|
|
module discretization_mesh
|
2018-08-30 15:25:13 +05:30
|
|
|
#include <petsc/finclude/petscdmplex.h>
|
2018-08-17 14:53:24 +05:30
|
|
|
#include <petsc/finclude/petscis.h>
|
|
|
|
#include <petsc/finclude/petscdmda.h>
|
2020-06-16 16:57:21 +05:30
|
|
|
use PETScdmplex
|
|
|
|
use PETScdmda
|
|
|
|
use PETScis
|
|
|
|
|
|
|
|
use DAMASK_interface
|
|
|
|
use IO
|
|
|
|
use debug
|
|
|
|
use discretization
|
|
|
|
use numerics
|
|
|
|
use FEsolving
|
|
|
|
use FEM_quadrature
|
2020-06-25 18:45:48 +05:30
|
|
|
use YAML_types
|
2020-06-16 16:57:21 +05:30
|
|
|
use prec
|
|
|
|
|
|
|
|
implicit none
|
|
|
|
private
|
|
|
|
|
|
|
|
integer, public, protected :: &
|
|
|
|
mesh_Nboundaries, &
|
|
|
|
mesh_NcpElemsGlobal
|
|
|
|
|
|
|
|
integer :: &
|
|
|
|
mesh_NcpElems !< total number of CP elements in mesh
|
2020-01-21 04:09:24 +05:30
|
|
|
|
2018-09-23 21:36:18 +05:30
|
|
|
!!!! BEGIN DEPRECATED !!!!!
|
2020-06-16 16:57:21 +05:30
|
|
|
integer, public, protected :: &
|
|
|
|
mesh_maxNips !< max number of IPs in any CP element
|
2018-09-23 21:36:18 +05:30
|
|
|
!!!! BEGIN DEPRECATED !!!!!
|
2018-08-17 14:53:24 +05:30
|
|
|
|
2020-06-16 16:57:21 +05:30
|
|
|
real(pReal), dimension(:,:), allocatable :: &
|
|
|
|
mesh_ipVolume, & !< volume associated with IP (initially!)
|
|
|
|
mesh_node0 !< node x,y,z coordinates (initially!)
|
2018-08-17 14:53:24 +05:30
|
|
|
|
2020-06-16 16:57:21 +05:30
|
|
|
real(pReal), dimension(:,:,:), allocatable :: &
|
|
|
|
mesh_ipCoordinates !< IP x,y,z coordinates (after deformation!)
|
2018-08-17 14:53:24 +05:30
|
|
|
|
2020-06-16 16:57:21 +05:30
|
|
|
DM, public :: geomMesh
|
2018-08-17 14:53:24 +05:30
|
|
|
|
2020-06-16 16:57:21 +05:30
|
|
|
PetscInt, dimension(:), allocatable, public, protected :: &
|
|
|
|
mesh_boundaries
|
2018-08-17 14:53:24 +05:30
|
|
|
|
2020-06-16 16:57:21 +05:30
|
|
|
public :: &
|
|
|
|
discretization_mesh_init, &
|
|
|
|
mesh_FEM_build_ipVolumes, &
|
|
|
|
mesh_FEM_build_ipCoordinates
|
2018-08-17 14:53:24 +05:30
|
|
|
|
|
|
|
contains
|
|
|
|
|
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief initializes the mesh by calling all necessary private routines the mesh module
|
|
|
|
!! Order and routines strongly depend on type of solver
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2020-05-06 21:20:59 +05:30
|
|
|
subroutine discretization_mesh_init(restart)
|
|
|
|
|
2020-05-07 01:45:09 +05:30
|
|
|
logical, intent(in) :: restart
|
2019-06-09 19:57:05 +05:30
|
|
|
|
2020-01-21 01:06:21 +05:30
|
|
|
integer, allocatable, dimension(:) :: chunkPos
|
|
|
|
integer :: dimPlex, &
|
2020-02-05 04:15:34 +05:30
|
|
|
mesh_Nnodes, & !< total number of nodes in mesh
|
2020-06-18 21:22:25 +05:30
|
|
|
j, l, &
|
|
|
|
debug_e, debug_i
|
2020-01-21 01:06:21 +05:30
|
|
|
PetscSF :: sf
|
|
|
|
DM :: globalMesh
|
2020-01-21 01:36:08 +05:30
|
|
|
PetscInt :: nFaceSets
|
2020-01-27 00:26:30 +05:30
|
|
|
PetscInt, pointer, dimension(:) :: pFaceSets
|
|
|
|
character(len=pStringLen), dimension(:), allocatable :: fileContent
|
2020-06-16 16:57:21 +05:30
|
|
|
IS :: faceSetIS
|
2020-01-21 01:06:21 +05:30
|
|
|
PetscErrorCode :: ierr
|
2020-06-16 17:29:59 +05:30
|
|
|
integer, dimension(:), allocatable :: &
|
|
|
|
homogenizationAt, &
|
|
|
|
microstructureAt
|
2020-06-16 22:45:01 +05:30
|
|
|
class(tNode), pointer :: &
|
2020-06-18 00:16:03 +05:30
|
|
|
num_mesh
|
2020-06-17 20:49:21 +05:30
|
|
|
integer :: integrationOrder !< order of quadrature rule required
|
2020-06-16 22:45:01 +05:30
|
|
|
|
2020-06-18 21:22:25 +05:30
|
|
|
write(6,'(/,a)') ' <<<+- mesh init -+>>>'
|
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------
|
|
|
|
! read numerics parameter
|
2020-06-18 00:16:03 +05:30
|
|
|
num_mesh => numerics_root%get('mesh',defaultVal=emptyDict)
|
|
|
|
integrationOrder = num_mesh%get_asInt('integrationorder',defaultVal = 2)
|
2020-06-16 22:45:01 +05:30
|
|
|
|
2020-06-18 21:22:25 +05:30
|
|
|
!---------------------------------------------------------------------------------
|
|
|
|
! read debug parameters
|
|
|
|
debug_e = debug_root%get_asInt('element',defaultVal=1)
|
|
|
|
debug_i = debug_root%get_asInt('integrationpoint',defaultVal=1)
|
2020-01-21 01:06:21 +05:30
|
|
|
|
2020-06-18 21:22:25 +05:30
|
|
|
|
2020-01-21 01:06:21 +05:30
|
|
|
call DMPlexCreateFromFile(PETSC_COMM_WORLD,geometryFile,PETSC_TRUE,globalMesh,ierr)
|
|
|
|
CHKERRQ(ierr)
|
|
|
|
call DMGetDimension(globalMesh,dimPlex,ierr)
|
|
|
|
CHKERRQ(ierr)
|
|
|
|
call DMGetStratumSize(globalMesh,'depth',dimPlex,mesh_NcpElemsGlobal,ierr)
|
|
|
|
CHKERRQ(ierr)
|
|
|
|
! get number of IDs in face sets (for boundary conditions?)
|
|
|
|
call DMGetLabelSize(globalMesh,'Face Sets',mesh_Nboundaries,ierr)
|
|
|
|
CHKERRQ(ierr)
|
|
|
|
call MPI_Bcast(mesh_Nboundaries,1,MPI_INTEGER,0,PETSC_COMM_WORLD,ierr)
|
|
|
|
call MPI_Bcast(mesh_NcpElemsGlobal,1,MPI_INTEGER,0,PETSC_COMM_WORLD,ierr)
|
|
|
|
call MPI_Bcast(dimPlex,1,MPI_INTEGER,0,PETSC_COMM_WORLD,ierr)
|
|
|
|
|
|
|
|
allocate(mesh_boundaries(mesh_Nboundaries), source = 0)
|
|
|
|
call DMGetLabelSize(globalMesh,'Face Sets',nFaceSets,ierr)
|
|
|
|
CHKERRQ(ierr)
|
|
|
|
call DMGetLabelIdIS(globalMesh,'Face Sets',faceSetIS,ierr)
|
|
|
|
CHKERRQ(ierr)
|
2020-01-22 00:03:22 +05:30
|
|
|
if (nFaceSets > 0) then
|
|
|
|
call ISGetIndicesF90(faceSetIS,pFaceSets,ierr)
|
|
|
|
CHKERRQ(ierr)
|
|
|
|
mesh_boundaries(1:nFaceSets) = pFaceSets
|
|
|
|
CHKERRQ(ierr)
|
|
|
|
call ISRestoreIndicesF90(faceSetIS,pFaceSets,ierr)
|
|
|
|
endif
|
2020-01-21 01:06:21 +05:30
|
|
|
call MPI_Bcast(mesh_boundaries,mesh_Nboundaries,MPI_INTEGER,0,PETSC_COMM_WORLD,ierr)
|
|
|
|
|
|
|
|
if (worldrank == 0) then
|
2020-06-05 18:14:31 +05:30
|
|
|
fileContent = IO_readlines(geometryFile)
|
2020-01-27 00:26:30 +05:30
|
|
|
l = 0
|
2020-01-21 01:06:21 +05:30
|
|
|
do
|
2020-01-27 00:26:30 +05:30
|
|
|
l = l + 1
|
2020-01-29 04:40:05 +05:30
|
|
|
if (IO_isBlank(fileContent(l))) cycle ! need also to ignore C and C++ style comments?
|
|
|
|
if (trim(fileContent(l)) == '$Elements') then
|
|
|
|
j = 0
|
|
|
|
l = l + 1
|
|
|
|
do
|
|
|
|
l = l + 1
|
|
|
|
if (trim(fileContent(l)) == '$EndElements') exit
|
|
|
|
chunkPos = IO_stringPos(fileContent(l))
|
|
|
|
if (chunkPos(1) == 3+IO_intValue(fileContent(l),chunkPos,3)+dimPlex+1) then
|
|
|
|
call DMSetLabelValue(globalMesh,'material',j,IO_intValue(fileContent(l),chunkPos,4),ierr)
|
|
|
|
CHKERRQ(ierr)
|
|
|
|
j = j + 1
|
|
|
|
endif
|
|
|
|
enddo
|
|
|
|
exit
|
2020-01-21 01:06:21 +05:30
|
|
|
endif
|
|
|
|
enddo
|
|
|
|
call DMClone(globalMesh,geomMesh,ierr)
|
|
|
|
CHKERRQ(ierr)
|
2020-06-16 16:57:21 +05:30
|
|
|
else
|
2020-01-21 01:06:21 +05:30
|
|
|
call DMPlexDistribute(globalMesh,0,sf,geomMesh,ierr)
|
|
|
|
CHKERRQ(ierr)
|
2020-06-16 16:57:21 +05:30
|
|
|
endif
|
2020-01-21 01:06:21 +05:30
|
|
|
|
|
|
|
call DMDestroy(globalMesh,ierr); CHKERRQ(ierr)
|
2020-06-16 16:57:21 +05:30
|
|
|
|
2020-01-21 01:06:21 +05:30
|
|
|
call DMGetStratumSize(geomMesh,'depth',dimPlex,mesh_NcpElems,ierr)
|
|
|
|
CHKERRQ(ierr)
|
|
|
|
call DMGetStratumSize(geomMesh,'depth',0,mesh_Nnodes,ierr)
|
|
|
|
CHKERRQ(ierr)
|
|
|
|
|
2020-06-16 16:57:21 +05:30
|
|
|
mesh_maxNips = FEM_nQuadrature(dimPlex,integrationOrder)
|
|
|
|
|
2020-03-20 19:25:10 +05:30
|
|
|
call mesh_FEM_build_ipCoordinates(dimPlex,FEM_quadrature_points(dimPlex,integrationOrder)%p)
|
2020-01-21 01:06:21 +05:30
|
|
|
call mesh_FEM_build_ipVolumes(dimPlex)
|
2020-06-16 16:57:21 +05:30
|
|
|
|
2020-06-16 17:29:59 +05:30
|
|
|
allocate(microstructureAt(mesh_NcpElems))
|
|
|
|
allocate(homogenizationAt(mesh_NcpElems),source=1)
|
2020-01-21 01:06:21 +05:30
|
|
|
do j = 1, mesh_NcpElems
|
2020-06-16 17:29:59 +05:30
|
|
|
call DMGetLabelValue(geomMesh,'material',j-1,microstructureAt(j),ierr)
|
2020-01-21 01:06:21 +05:30
|
|
|
CHKERRQ(ierr)
|
2020-06-16 16:57:21 +05:30
|
|
|
end do
|
|
|
|
|
|
|
|
if (debug_e < 1 .or. debug_e > mesh_NcpElems) call IO_error(602,ext_msg='element')
|
|
|
|
if (debug_i < 1 .or. debug_i > mesh_maxNips) call IO_error(602,ext_msg='IP')
|
|
|
|
|
|
|
|
FEsolving_execElem = [1,mesh_NcpElems] ! parallel loop bounds set to comprise all DAMASK elements
|
|
|
|
FEsolving_execIP = [1,mesh_maxNips]
|
2020-01-21 01:06:21 +05:30
|
|
|
|
|
|
|
allocate(mesh_node0(3,mesh_Nnodes),source=0.0_pReal)
|
|
|
|
|
2020-06-16 17:29:59 +05:30
|
|
|
call discretization_init(microstructureAt,homogenizationAt,&
|
2019-06-07 09:37:00 +05:30
|
|
|
reshape(mesh_ipCoordinates,[3,mesh_maxNips*mesh_NcpElems]), &
|
|
|
|
mesh_node0)
|
2020-06-16 16:57:21 +05:30
|
|
|
|
2020-03-20 19:38:07 +05:30
|
|
|
end subroutine discretization_mesh_init
|
2018-08-17 14:53:24 +05:30
|
|
|
|
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief Calculates IP volume. Allocates global array 'mesh_ipVolume'
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
subroutine mesh_FEM_build_ipVolumes(dimPlex)
|
2020-06-16 16:57:21 +05:30
|
|
|
|
2020-01-21 04:16:48 +05:30
|
|
|
PetscInt,intent(in):: dimPlex
|
2019-06-16 00:19:48 +05:30
|
|
|
PetscReal :: vol
|
2020-01-21 04:16:48 +05:30
|
|
|
PetscReal, pointer,dimension(:) :: pCent, pNorm
|
2019-06-16 00:19:48 +05:30
|
|
|
PetscInt :: cellStart, cellEnd, cell
|
|
|
|
PetscErrorCode :: ierr
|
2020-06-16 16:57:21 +05:30
|
|
|
|
2020-01-21 04:16:48 +05:30
|
|
|
allocate(mesh_ipVolume(mesh_maxNips,mesh_NcpElems),source=0.0_pReal)
|
2020-06-16 16:57:21 +05:30
|
|
|
|
2019-06-16 00:19:48 +05:30
|
|
|
call DMPlexGetHeightStratum(geomMesh,0,cellStart,cellEnd,ierr); CHKERRQ(ierr)
|
2020-01-21 04:16:48 +05:30
|
|
|
allocate(pCent(dimPlex))
|
|
|
|
allocate(pNorm(dimPlex))
|
2019-06-16 00:19:48 +05:30
|
|
|
do cell = cellStart, cellEnd-1
|
2020-06-16 16:57:21 +05:30
|
|
|
call DMPlexComputeCellGeometryFVM(geomMesh,cell,vol,pCent,pNorm,ierr)
|
2019-06-16 00:19:48 +05:30
|
|
|
CHKERRQ(ierr)
|
|
|
|
mesh_ipVolume(:,cell+1) = vol/real(mesh_maxNips,pReal)
|
2020-06-16 16:57:21 +05:30
|
|
|
enddo
|
2018-08-17 14:53:24 +05:30
|
|
|
|
|
|
|
end subroutine mesh_FEM_build_ipVolumes
|
|
|
|
|
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief Calculates IP Coordinates. Allocates global array 'mesh_ipCoordinates'
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
subroutine mesh_FEM_build_ipCoordinates(dimPlex,qPoints)
|
|
|
|
|
2019-06-16 00:19:48 +05:30
|
|
|
PetscInt, intent(in) :: dimPlex
|
|
|
|
PetscReal, intent(in) :: qPoints(mesh_maxNips*dimPlex)
|
2020-06-16 16:57:21 +05:30
|
|
|
|
2020-01-21 04:16:48 +05:30
|
|
|
PetscReal, pointer,dimension(:) :: pV0, pCellJ, pInvcellJ
|
2019-06-16 00:19:48 +05:30
|
|
|
PetscReal :: detJ
|
|
|
|
PetscInt :: cellStart, cellEnd, cell, qPt, dirI, dirJ, qOffset
|
|
|
|
PetscErrorCode :: ierr
|
2020-06-16 16:57:21 +05:30
|
|
|
|
|
|
|
|
2019-06-16 00:19:48 +05:30
|
|
|
allocate(mesh_ipCoordinates(3,mesh_maxNips,mesh_NcpElems),source=0.0_pReal)
|
2020-06-16 16:57:21 +05:30
|
|
|
|
2020-01-21 04:16:48 +05:30
|
|
|
allocate(pV0(dimPlex))
|
|
|
|
allocatE(pCellJ(dimPlex**2))
|
|
|
|
allocatE(pinvCellJ(dimPlex**2))
|
2019-06-16 00:19:48 +05:30
|
|
|
call DMPlexGetHeightStratum(geomMesh,0,cellStart,cellEnd,ierr); CHKERRQ(ierr)
|
2020-06-16 16:57:21 +05:30
|
|
|
do cell = cellStart, cellEnd-1 !< loop over all elements
|
2019-06-16 00:19:48 +05:30
|
|
|
call DMPlexComputeCellGeometryAffineFEM(geomMesh,cell,pV0,pCellJ,pInvcellJ,detJ,ierr)
|
|
|
|
CHKERRQ(ierr)
|
|
|
|
qOffset = 0
|
|
|
|
do qPt = 1, mesh_maxNips
|
|
|
|
do dirI = 1, dimPlex
|
|
|
|
mesh_ipCoordinates(dirI,qPt,cell+1) = pV0(dirI)
|
|
|
|
do dirJ = 1, dimPlex
|
|
|
|
mesh_ipCoordinates(dirI,qPt,cell+1) = mesh_ipCoordinates(dirI,qPt,cell+1) + &
|
|
|
|
pCellJ((dirI-1)*dimPlex+dirJ)*(qPoints(qOffset+dirJ) + 1.0)
|
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
qOffset = qOffset + dimPlex
|
|
|
|
enddo
|
2020-06-16 16:57:21 +05:30
|
|
|
enddo
|
2018-08-17 14:53:24 +05:30
|
|
|
|
|
|
|
end subroutine mesh_FEM_build_ipCoordinates
|
|
|
|
|
2020-03-20 19:38:07 +05:30
|
|
|
end module discretization_mesh
|