Merge branch '231-no-petscscalar-and-_range' into 'development'
Resolve "no PetscScalar and ?_RANGE" Closes #231 See merge request damask/DAMASK!668
This commit is contained in:
commit
99b2a509d9
|
@ -75,7 +75,7 @@ subroutine grid_damage_spectral_init()
|
|||
PetscInt, dimension(0:worldsize-1) :: localK
|
||||
integer :: i, j, k, ce
|
||||
DM :: damage_grid
|
||||
PetscScalar, dimension(:,:,:), pointer :: phi_PETSc
|
||||
real(pReal), dimension(:,:,:), pointer :: phi_PETSc
|
||||
Vec :: uBound, lBound
|
||||
integer(MPI_INTEGER_KIND) :: err_MPI
|
||||
PetscErrorCode :: err_PETSc
|
||||
|
@ -273,7 +273,7 @@ subroutine grid_damage_spectral_forward(cutBack)
|
|||
|
||||
integer :: i, j, k, ce
|
||||
DM :: dm_local
|
||||
PetscScalar, dimension(:,:,:), pointer :: phi_PETSc
|
||||
real(pReal), dimension(:,:,:), pointer :: phi_PETSc
|
||||
PetscErrorCode :: err_PETSc
|
||||
|
||||
|
||||
|
@ -335,15 +335,13 @@ end subroutine grid_damage_spectral_restartWrite
|
|||
!--------------------------------------------------------------------------------------------------
|
||||
!> @brief Construct the residual vector.
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
subroutine formResidual(in,x_scal,r,dummy,err_PETSc)
|
||||
subroutine formResidual(residual_subdomain,x_scal,r,dummy,err_PETSc)
|
||||
|
||||
DMDALocalInfo, dimension(DMDA_LOCAL_INFO_SIZE) :: &
|
||||
in
|
||||
PetscScalar, dimension( &
|
||||
XG_RANGE,YG_RANGE,ZG_RANGE), intent(in) :: &
|
||||
residual_subdomain
|
||||
real(pReal), dimension(cells(1),cells(2),cells3), intent(in) :: &
|
||||
x_scal
|
||||
PetscScalar, dimension( &
|
||||
X_RANGE,Y_RANGE,Z_RANGE), intent(out) :: &
|
||||
real(pReal), dimension(cells(1),cells(2),cells3), intent(out) :: &
|
||||
r !< residual
|
||||
PetscObject :: dummy
|
||||
PetscErrorCode, intent(out) :: err_PETSc
|
||||
|
|
|
@ -539,8 +539,8 @@ subroutine formResidual(da_local,x_local, &
|
|||
|
||||
DM :: da_local
|
||||
Vec :: x_local, f_local
|
||||
PetscScalar, pointer,dimension(:,:,:,:) :: x_scal, r
|
||||
PetscScalar, dimension(8,3) :: x_elem, f_elem
|
||||
real(pReal), pointer,dimension(:,:,:,:) :: x_scal, r
|
||||
real(pReal), dimension(8,3) :: x_elem, f_elem
|
||||
PetscInt :: i, ii, j, jj, k, kk, ctr, ele
|
||||
PetscInt :: &
|
||||
PETScIter, &
|
||||
|
@ -659,12 +659,12 @@ subroutine formJacobian(da_local,x_local,Jac_pre,Jac,dummy,err_PETSc)
|
|||
Vec :: x_local, coordinates
|
||||
Mat :: Jac_pre, Jac
|
||||
MatStencil,dimension(4,24) :: row, col
|
||||
PetscScalar,pointer,dimension(:,:,:,:) :: x_scal
|
||||
PetscScalar,dimension(24,24) :: K_ele
|
||||
PetscScalar,dimension(9,24) :: BMatFull
|
||||
real(pReal),pointer,dimension(:,:,:,:) :: x_scal
|
||||
real(pReal),dimension(24,24) :: K_ele
|
||||
real(pReal),dimension(9,24) :: BMatFull
|
||||
PetscInt :: i, ii, j, jj, k, kk, ctr, ce
|
||||
PetscInt,dimension(3),parameter :: rows = [0, 1, 2]
|
||||
PetscScalar :: diag
|
||||
real(pReal) :: diag
|
||||
PetscObject :: dummy
|
||||
MatNullSpace :: matnull
|
||||
PetscErrorCode :: err_PETSc
|
||||
|
|
|
@ -109,7 +109,7 @@ subroutine grid_mechanical_spectral_basic_init()
|
|||
real(pReal), dimension(3,3,cells(1),cells(2),cells3) :: P
|
||||
PetscErrorCode :: err_PETSc
|
||||
integer(MPI_INTEGER_KIND) :: err_MPI
|
||||
PetscScalar, pointer, dimension(:,:,:,:) :: &
|
||||
real(pReal), pointer, dimension(:,:,:,:) :: &
|
||||
F ! pointer to solution data
|
||||
PetscInt, dimension(0:worldsize-1) :: localK
|
||||
real(pReal), dimension(3,3,product(cells(1:2))*cells3) :: temp33n
|
||||
|
@ -332,7 +332,7 @@ subroutine grid_mechanical_spectral_basic_forward(cutBack,guess,Delta_t,Delta_t_
|
|||
type(tRotation), intent(in) :: &
|
||||
rotation_BC
|
||||
PetscErrorCode :: err_PETSc
|
||||
PetscScalar, pointer, dimension(:,:,:,:) :: F
|
||||
real(pReal), pointer, dimension(:,:,:,:) :: F
|
||||
|
||||
|
||||
call DMDAVecGetArrayF90(da,solution_vec,F,err_PETSc)
|
||||
|
@ -397,7 +397,7 @@ end subroutine grid_mechanical_spectral_basic_forward
|
|||
subroutine grid_mechanical_spectral_basic_updateCoords
|
||||
|
||||
PetscErrorCode :: err_PETSc
|
||||
PetscScalar, dimension(:,:,:,:), pointer :: F
|
||||
real(pReal), dimension(:,:,:,:), pointer :: F
|
||||
|
||||
call DMDAVecGetArrayF90(da,solution_vec,F,err_PETSc)
|
||||
CHKERRQ(err_PETSc)
|
||||
|
@ -415,7 +415,7 @@ subroutine grid_mechanical_spectral_basic_restartWrite
|
|||
|
||||
PetscErrorCode :: err_PETSc
|
||||
integer(HID_T) :: fileHandle, groupHandle
|
||||
PetscScalar, dimension(:,:,:,:), pointer :: F
|
||||
real(pReal), dimension(:,:,:,:), pointer :: F
|
||||
|
||||
call DMDAVecGetArrayF90(da,solution_vec,F,err_PETSc)
|
||||
CHKERRQ(err_PETSc)
|
||||
|
@ -496,14 +496,15 @@ end subroutine converged
|
|||
!--------------------------------------------------------------------------------------------------
|
||||
!> @brief Construct the residual vector.
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
subroutine formResidual(in, F, &
|
||||
subroutine formResidual(residual_subdomain, F, &
|
||||
r, dummy, err_PETSc)
|
||||
|
||||
DMDALocalInfo, dimension(DMDA_LOCAL_INFO_SIZE) :: in !< DMDA info (needs to be named "in" for macros like XRANGE to work)
|
||||
PetscScalar, dimension(3,3,XG_RANGE,YG_RANGE,ZG_RANGE), &
|
||||
intent(in) :: F !< deformation gradient field
|
||||
PetscScalar, dimension(3,3,X_RANGE,Y_RANGE,Z_RANGE), &
|
||||
intent(out) :: r !< residuum field
|
||||
DMDALocalInfo, dimension(DMDA_LOCAL_INFO_SIZE) :: &
|
||||
residual_subdomain !< DMDA info (needs to be named "in" for macros like XRANGE to work)
|
||||
real(pReal), dimension(3,3,cells(1),cells(2),cells3), intent(in) :: &
|
||||
F !< deformation gradient field
|
||||
real(pReal), dimension(3,3,cells(1),cells(2),cells3), intent(out) :: &
|
||||
r !< residuum field
|
||||
PetscObject :: dummy
|
||||
PetscErrorCode :: err_PETSc
|
||||
|
||||
|
|
|
@ -120,7 +120,7 @@ subroutine grid_mechanical_spectral_polarisation_init()
|
|||
real(pReal), dimension(3,3,cells(1),cells(2),cells3) :: P
|
||||
PetscErrorCode :: err_PETSc
|
||||
integer(MPI_INTEGER_KIND) :: err_MPI
|
||||
PetscScalar, pointer, dimension(:,:,:,:) :: &
|
||||
real(pReal), pointer, dimension(:,:,:,:) :: &
|
||||
FandF_tau, & ! overall pointer to solution data
|
||||
F, & ! specific (sub)pointer
|
||||
F_tau ! specific (sub)pointer
|
||||
|
@ -365,7 +365,7 @@ subroutine grid_mechanical_spectral_polarisation_forward(cutBack,guess,Delta_t,D
|
|||
type(tRotation), intent(in) :: &
|
||||
rotation_BC
|
||||
PetscErrorCode :: err_PETSc
|
||||
PetscScalar, pointer, dimension(:,:,:,:) :: FandF_tau, F, F_tau
|
||||
real(pReal), pointer, dimension(:,:,:,:) :: FandF_tau, F, F_tau
|
||||
integer :: i, j, k
|
||||
real(pReal), dimension(3,3) :: F_lambda33
|
||||
|
||||
|
@ -452,7 +452,7 @@ end subroutine grid_mechanical_spectral_polarisation_forward
|
|||
subroutine grid_mechanical_spectral_polarisation_updateCoords
|
||||
|
||||
PetscErrorCode :: err_PETSc
|
||||
PetscScalar, dimension(:,:,:,:), pointer :: FandF_tau
|
||||
real(pReal), dimension(:,:,:,:), pointer :: FandF_tau
|
||||
|
||||
call DMDAVecGetArrayF90(da,solution_vec,FandF_tau,err_PETSc)
|
||||
CHKERRQ(err_PETSc)
|
||||
|
@ -470,7 +470,7 @@ subroutine grid_mechanical_spectral_polarisation_restartWrite
|
|||
|
||||
PetscErrorCode :: err_PETSc
|
||||
integer(HID_T) :: fileHandle, groupHandle
|
||||
PetscScalar, dimension(:,:,:,:), pointer :: FandF_tau, F, F_tau
|
||||
real(pReal), dimension(:,:,:,:), pointer :: FandF_tau, F, F_tau
|
||||
|
||||
call DMDAVecGetArrayF90(da,solution_vec,FandF_tau,err_PETSc)
|
||||
CHKERRQ(err_PETSc)
|
||||
|
@ -558,18 +558,18 @@ end subroutine converged
|
|||
!--------------------------------------------------------------------------------------------------
|
||||
!> @brief Construct the residual vector.
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
subroutine formResidual(in, FandF_tau, &
|
||||
subroutine formResidual(residual_subdomain, FandF_tau, &
|
||||
r, dummy,err_PETSc)
|
||||
|
||||
DMDALocalInfo, dimension(DMDA_LOCAL_INFO_SIZE) :: in !< DMDA info (needs to be named "in" for macros like XRANGE to work)
|
||||
PetscScalar, dimension(3,3,2,XG_RANGE,YG_RANGE,ZG_RANGE), &
|
||||
target, intent(in) :: FandF_tau
|
||||
PetscScalar, dimension(3,3,2,X_RANGE,Y_RANGE,Z_RANGE),&
|
||||
target, intent(out) :: r !< residuum field
|
||||
DMDALocalInfo, dimension(DMDA_LOCAL_INFO_SIZE) :: residual_subdomain !< DMDA info (needs to be named "in" for macros like XRANGE to work)
|
||||
real(pReal), dimension(3,3,2,cells(1),cells(2),cells3), target, intent(in) :: &
|
||||
FandF_tau !< deformation gradient field
|
||||
real(pReal), dimension(3,3,2,cells(1),cells(2),cells3), target, intent(out) :: &
|
||||
r !< residuum field
|
||||
PetscObject :: dummy
|
||||
PetscErrorCode :: err_PETSc
|
||||
|
||||
PetscScalar, pointer, dimension(:,:,:,:,:) :: &
|
||||
real(pReal), pointer, dimension(:,:,:,:,:) :: &
|
||||
F, &
|
||||
F_tau, &
|
||||
r_F, &
|
||||
|
@ -582,14 +582,10 @@ subroutine formResidual(in, FandF_tau, &
|
|||
i, j, k, e
|
||||
|
||||
|
||||
F => FandF_tau(1:3,1:3,1,&
|
||||
XG_RANGE,YG_RANGE,ZG_RANGE)
|
||||
F_tau => FandF_tau(1:3,1:3,2,&
|
||||
XG_RANGE,YG_RANGE,ZG_RANGE)
|
||||
r_F => r(1:3,1:3,1,&
|
||||
X_RANGE, Y_RANGE, Z_RANGE)
|
||||
r_F_tau => r(1:3,1:3,2,&
|
||||
X_RANGE, Y_RANGE, Z_RANGE)
|
||||
F => FandF_tau(1:3,1:3,1,1:cells(1),1:cells(2),1:cells3)
|
||||
F_tau => FandF_tau(1:3,1:3,2,1:cells(1),1:cells(2),1:cells3)
|
||||
r_F => r(1:3,1:3,1,1:cells(1),1:cells(2),1:cells3)
|
||||
r_F_tau => r(1:3,1:3,2,1:cells(1),1:cells(2),1:cells3)
|
||||
|
||||
F_av = sum(sum(sum(F,dim=5),dim=4),dim=3) * wgt
|
||||
call MPI_Allreduce(MPI_IN_PLACE,F_av,9_MPI_INTEGER_KIND,MPI_DOUBLE,MPI_SUM,MPI_COMM_WORLD,err_MPI)
|
||||
|
|
|
@ -74,7 +74,7 @@ subroutine grid_thermal_spectral_init()
|
|||
PetscInt, dimension(0:worldsize-1) :: localK
|
||||
integer :: i, j, k, ce
|
||||
DM :: thermal_grid
|
||||
PetscScalar, dimension(:,:,:), pointer :: T_PETSc
|
||||
real(pReal), dimension(:,:,:), pointer :: T_PETSc
|
||||
integer(MPI_INTEGER_KIND) :: err_MPI
|
||||
PetscErrorCode :: err_PETSc
|
||||
integer(HID_T) :: fileHandle, groupHandle
|
||||
|
@ -249,7 +249,7 @@ subroutine grid_thermal_spectral_forward(cutBack)
|
|||
|
||||
integer :: i, j, k, ce
|
||||
DM :: dm_local
|
||||
PetscScalar, dimension(:,:,:), pointer :: T_PETSc
|
||||
real(pReal), dimension(:,:,:), pointer :: T_PETSc
|
||||
PetscErrorCode :: err_PETSc
|
||||
|
||||
|
||||
|
@ -288,7 +288,7 @@ subroutine grid_thermal_spectral_restartWrite
|
|||
PetscErrorCode :: err_PETSc
|
||||
DM :: dm_local
|
||||
integer(HID_T) :: fileHandle, groupHandle
|
||||
PetscScalar, dimension(:,:,:), pointer :: T
|
||||
real(pReal), dimension(:,:,:), pointer :: T
|
||||
|
||||
call SNESGetDM(SNES_thermal,dm_local,err_PETSc);
|
||||
CHKERRQ(err_PETSc)
|
||||
|
@ -315,15 +315,13 @@ end subroutine grid_thermal_spectral_restartWrite
|
|||
!--------------------------------------------------------------------------------------------------
|
||||
!> @brief Construct the residual vector.
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
subroutine formResidual(in,x_scal,r,dummy,err_PETSc)
|
||||
subroutine formResidual(residual_subdomain,x_scal,r,dummy,err_PETSc)
|
||||
|
||||
DMDALocalInfo, dimension(DMDA_LOCAL_INFO_SIZE) :: &
|
||||
in
|
||||
PetscScalar, dimension( &
|
||||
XG_RANGE,YG_RANGE,ZG_RANGE), intent(in) :: &
|
||||
residual_subdomain
|
||||
real(pReal), dimension(cells(1),cells(2),cells3), intent(in) :: &
|
||||
x_scal
|
||||
PetscScalar, dimension( &
|
||||
X_RANGE,Y_RANGE,Z_RANGE), intent(out) :: &
|
||||
real(pReal), dimension(cells(1),cells(2),cells3), intent(out) :: &
|
||||
r !< residual
|
||||
PetscObject :: dummy
|
||||
PetscErrorCode, intent(out) :: err_PETSc
|
||||
|
|
|
@ -182,8 +182,8 @@ subroutine utilities_projectBCValues(localVec,section,field,comp,bcPointsIS,BCVa
|
|||
PetscSection :: section
|
||||
IS :: bcPointsIS
|
||||
PetscInt, pointer :: bcPoints(:)
|
||||
PetscScalar, pointer :: localArray(:)
|
||||
PetscScalar :: BCValue,BCDotValue,timeinc
|
||||
real(pReal), pointer :: localArray(:)
|
||||
real(pReal) :: BCValue,BCDotValue,timeinc
|
||||
PetscErrorCode :: err_PETSc
|
||||
|
||||
|
||||
|
|
|
@ -120,8 +120,8 @@ subroutine FEM_mechanical_init(fieldBC)
|
|||
PetscReal :: detJ
|
||||
PetscReal, allocatable, target :: cellJMat(:,:)
|
||||
|
||||
PetscScalar, pointer, dimension(:) :: px_scal
|
||||
PetscScalar, allocatable, target, dimension(:) :: x_scal
|
||||
real(pReal), pointer, dimension(:) :: px_scal
|
||||
real(pReal), allocatable, target, dimension(:) :: x_scal
|
||||
|
||||
character(len=*), parameter :: prefix = 'mechFE_'
|
||||
PetscErrorCode :: err_PETSc
|
||||
|
@ -369,8 +369,8 @@ subroutine FEM_mechanical_formResidual(dm_local,xx_local,f_local,dummy,err_PETSc
|
|||
PetscDS :: prob
|
||||
Vec :: x_local, f_local, xx_local
|
||||
PetscSection :: section
|
||||
PetscScalar, dimension(:), pointer :: x_scal, pf_scal
|
||||
PetscScalar, dimension(cellDof), target :: f_scal
|
||||
real(pReal), dimension(:), pointer :: x_scal, pf_scal
|
||||
real(pReal), dimension(cellDof), target :: f_scal
|
||||
PetscReal :: IcellJMat(dimPlex,dimPlex)
|
||||
PetscReal, dimension(:),pointer :: pV0, pCellJ, pInvcellJ, basisField, basisFieldDer
|
||||
PetscInt :: cellStart, cellEnd, cell, field, face, &
|
||||
|
@ -517,10 +517,10 @@ subroutine FEM_mechanical_formJacobian(dm_local,xx_local,Jac_pre,Jac,dummy,err_P
|
|||
PetscReal, dimension(:), pointer :: basisField, basisFieldDer, &
|
||||
pV0, pCellJ, pInvcellJ
|
||||
|
||||
PetscScalar, dimension(:), pointer :: pK_e, x_scal
|
||||
real(pReal), dimension(:), pointer :: pK_e, x_scal
|
||||
|
||||
PetscScalar,dimension(cellDOF,cellDOF), target :: K_e
|
||||
PetscScalar,dimension(cellDOF,cellDOF) :: K_eA, K_eB
|
||||
real(pReal),dimension(cellDOF,cellDOF), target :: K_e
|
||||
real(pReal),dimension(cellDOF,cellDOF) :: K_eA, K_eB
|
||||
|
||||
PetscInt :: cellStart, cellEnd, cell, field, face, &
|
||||
qPt, basis, comp, cidx,bcSize, m, i
|
||||
|
@ -777,7 +777,7 @@ subroutine FEM_mechanical_updateCoords()
|
|||
PetscQuadrature :: mechQuad
|
||||
PetscReal, dimension(:), pointer :: basisField, basisFieldDer, &
|
||||
nodeCoords_linear !< nodal coordinates (dimPlex*Nnodes)
|
||||
PetscScalar, dimension(:), pointer :: x_scal
|
||||
real(pReal), dimension(:), pointer :: x_scal
|
||||
|
||||
call SNESGetDM(mechanical_snes,dm_local,err_PETSc)
|
||||
CHKERRQ(err_PETSc)
|
||||
|
|
|
@ -23,8 +23,10 @@ module prec
|
|||
integer, parameter :: pI32 = selected_int_kind(9) !< number with at least up to +-1e9 (typically 32 bit)
|
||||
integer, parameter :: pI64 = selected_int_kind(18) !< number with at least up to +-1e18 (typically 64 bit)
|
||||
#ifdef PETSC
|
||||
PetscInt, private :: dummy
|
||||
integer, parameter :: pPETSCINT = kind(dummy)
|
||||
PetscInt, private :: dummy_int
|
||||
integer, parameter :: pPETSCINT = kind(dummy_int)
|
||||
PetscScalar, private :: dummy_scalar
|
||||
real(pReal), parameter, private :: pPETSCSCALAR = kind(dummy_scalar)
|
||||
#endif
|
||||
integer, parameter :: pSTRINGLEN = 256 !< default string length
|
||||
integer, parameter :: pPATHLEN = 4096 !< maximum length of a path name on linux
|
||||
|
@ -253,6 +255,9 @@ subroutine selfTest()
|
|||
real(pReal), dimension(2) :: r
|
||||
|
||||
|
||||
#ifdef PETSC
|
||||
if (pReal /= pPETSCSCALAR) error stop 'PETSc and DAMASK scalar datatypes do not match'
|
||||
#endif
|
||||
realloc_lhs_test = [1,2]
|
||||
if (any(realloc_lhs_test/=[1,2])) error stop 'LHS allocation'
|
||||
|
||||
|
|
Loading…
Reference in New Issue