reordered PETSc initialization (now first part done in interface, setting of parameters done in numerics), removed unnecessary includes for PETSc

still not running with gfortran, use at own risk!
This commit is contained in:
Martin Diehl 2012-11-06 16:00:51 +00:00
parent a86d528a4a
commit 3ada4897fb
5 changed files with 138 additions and 128 deletions

View File

@ -30,6 +30,10 @@ module DAMASK_interface
implicit none implicit none
private private
#ifdef PETSc
#include <finclude/petscsys.h>
#endif
logical, protected, public :: & logical, protected, public :: &
appendToOutFile = .false. !< Append to existing spectralOut file (in case of restart, not in case of regridding) appendToOutFile = .false. !< Append to existing spectralOut file (in case of restart, not in case of regridding)
integer(pInt), protected, public :: & integer(pInt), protected, public :: &
@ -61,6 +65,7 @@ subroutine DAMASK_interface_init(loadCaseParameterIn,geometryParameterIn)
use, intrinsic :: iso_fortran_env ! to get compiler_version and compiler_options (at least for gfortran 4.6 at the moment) use, intrinsic :: iso_fortran_env ! to get compiler_version and compiler_options (at least for gfortran 4.6 at the moment)
implicit none implicit none
character(len=1024), optional, intent(in) :: & character(len=1024), optional, intent(in) :: &
loadCaseParameterIn, & loadCaseParameterIn, &
geometryParameterIn geometryParameterIn
@ -83,6 +88,14 @@ subroutine DAMASK_interface_init(loadCaseParameterIn,geometryParameterIn)
gotLoadCase = .false., & gotLoadCase = .false., &
gotGeometry = .false. gotGeometry = .false.
#ifdef PETSc
PetscErrorCode :: ierr
!--------------------------------------------------------------------------------------------------
! PETSc Init
call PetscInitialize(PETSC_NULL_CHARACTER,ierr) ! according to PETSc manual, that should be the first line in the code
CHKERRQ(ierr) ! this is a macro definition, it is case sensitive
#endif
open(6, encoding='UTF-8') ! modern fortran compilers (gfortran >4.4, ifort >11 support it) open(6, encoding='UTF-8') ! modern fortran compilers (gfortran >4.4, ifort >11 support it)
write(6,'(a)') '' write(6,'(a)') ''
write(6,'(a)') ' <<<+- DAMASK_spectral_interface init -+>>>' write(6,'(a)') ' <<<+- DAMASK_spectral_interface init -+>>>'

View File

@ -19,14 +19,8 @@ module DAMASK_spectral_SolverAL
implicit none implicit none
#include <finclude/petscsys.h> #include <finclude/petscsys.h>
#include <finclude/petscvec.h>
#include <finclude/petscdmda.h> #include <finclude/petscdmda.h>
#include <finclude/petscis.h>
#include <finclude/petscmat.h>
#include <finclude/petscksp.h>
#include <finclude/petscpc.h>
#include <finclude/petscsnes.h> #include <finclude/petscsnes.h>
#include <finclude/petscvec.h90>
#include <finclude/petscdmda.h90> #include <finclude/petscdmda.h90>
#include <finclude/petscsnes.h90> #include <finclude/petscsnes.h90>
@ -44,9 +38,9 @@ module DAMASK_spectral_SolverAL
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! PETSc data ! PETSc data
DM, private :: da DM, private :: da
SNES, private :: snes SNES, private :: snes
Vec, private :: solution_vec Vec, private :: solution_vec
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! common pointwise data ! common pointwise data
@ -112,10 +106,9 @@ subroutine AL_init()
integer(pInt) :: i,j,k integer(pInt) :: i,j,k
real(pReal), dimension(:,:,:,:,:), allocatable :: P real(pReal), dimension(:,:,:,:,:), allocatable :: P
PetscErrorCode :: ierr_psc PetscErrorCode :: ierr
PetscObject :: dummy PetscObject :: dummy
PetscMPIInt :: rank PetscScalar, pointer, dimension(:,:,:,:) :: xx_psc, F, F_lambda
PetscScalar, pointer :: xx_psc(:,:,:,:), F(:,:,:,:), F_lambda(:,:,:,:)
call Utilities_init() call Utilities_init()
@ -135,24 +128,20 @@ subroutine AL_init()
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! PETSc Init ! PETSc Init
call PetscInitialize(PETSC_NULL_CHARACTER,ierr_psc) call SNESCreate(PETSC_COMM_WORLD,snes,ierr)
call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr_psc)
call SNESCreate(PETSC_COMM_WORLD,snes,ierr_psc)
call DMDACreate3d(PETSC_COMM_WORLD, & call DMDACreate3d(PETSC_COMM_WORLD, &
DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_NONE, & DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_NONE, &
DMDA_STENCIL_BOX,res(1),res(2),res(3),PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE, & DMDA_STENCIL_BOX,res(1),res(2),res(3),PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE, &
18,1,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,da,ierr_psc) 18,1,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,da,ierr)
call DMCreateGlobalVector(da,solution_vec,ierr)
call DMCreateGlobalVector(da,solution_vec,ierr_psc) call DMDASetLocalFunction(da,AL_formResidual,ierr)
call DMDASetLocalFunction(da,AL_formResidual,ierr_psc) call SNESSetDM(snes,da,ierr)
call SNESSetDM(snes,da,ierr_psc) call SNESSetConvergenceTest(snes,AL_converged,dummy,PETSC_NULL_FUNCTION,ierr)
call SNESSetConvergenceTest(snes,AL_converged,dummy,PETSC_NULL_FUNCTION,ierr_psc) call SNESSetFromOptions(snes,ierr)
call PetscOptionsInsertString(petsc_options,ierr_psc)
call SNESSetFromOptions(snes,ierr_psc)
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! init fields ! init fields
call DMDAVecGetArrayF90(da,solution_vec,xx_psc,ierr_psc) call DMDAVecGetArrayF90(da,solution_vec,xx_psc,ierr)
F => xx_psc(0:8,:,:,:) F => xx_psc(0:8,:,:,:)
F_lambda => xx_psc(9:17,:,:,:) F_lambda => xx_psc(9:17,:,:,:)
if (restartInc == 1_pInt) then ! no deformation (no restart) if (restartInc == 1_pInt) then ! no deformation (no restart)
@ -192,9 +181,9 @@ subroutine AL_init()
coordinates = 0.0 ! change it later!!! coordinates = 0.0 ! change it later!!!
endif endif
!print*, shape(F)
call Utilities_constitutiveResponse(coordinates,F,F,temperature,0.0_pReal,P,C,P_av,.false.,math_I3) call Utilities_constitutiveResponse(coordinates,F,F,temperature,0.0_pReal,P,C,P_av,.false.,math_I3)
call DMDAVecRestoreArrayF90(da,solution_vec,xx_psc,ierr_psc) call DMDAVecRestoreArrayF90(da,solution_vec,xx_psc,ierr)
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! reference stiffness ! reference stiffness
@ -258,9 +247,9 @@ subroutine AL_init()
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! !
PetscScalar, pointer :: xx_psc(:,:,:,:), F(:,:,:,:), F_lambda(:,:,:,:) PetscScalar, dimension(:,:,:,:), pointer :: xx_psc, F, F_lambda
PetscErrorCode ierr_psc PetscErrorCode :: ierr
SNESConvergedReason reason SNESConvergedReason ::reason
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! restart information for spectral solver ! restart information for spectral solver
@ -275,7 +264,7 @@ subroutine AL_init()
close(777) close(777)
endif endif
AL_solution%converged =.false. AL_solution%converged =.false.
call DMDAVecGetArrayF90(da,solution_vec,xx_psc,ierr_psc) call DMDAVecGetArrayF90(da,solution_vec,xx_psc,ierr)
F => xx_psc(0:8,:,:,:) F => xx_psc(0:8,:,:,:)
F_lambda => xx_psc(9:17,:,:,:) F_lambda => xx_psc(9:17,:,:,:)
@ -323,7 +312,7 @@ else
F = reshape(Utilities_forwardField(timeinc,F_aim,F_lastInc,Fdot),[9,res(1),res(2),res(3)]) F = reshape(Utilities_forwardField(timeinc,F_aim,F_lastInc,Fdot),[9,res(1),res(2),res(3)])
F_lambda = reshape(Utilities_forwardField(timeinc,F_aim,F_lambda_lastInc,F_lambdadot),[9,res(1),res(2),res(3)]) F_lambda = reshape(Utilities_forwardField(timeinc,F_aim,F_lambda_lastInc,F_lambdadot),[9,res(1),res(2),res(3)])
call DMDAVecRestoreArrayF90(da,solution_vec,xx_psc,ierr_psc) call DMDAVecRestoreArrayF90(da,solution_vec,xx_psc,ierr)
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! update stiffness (and gamma operator) ! update stiffness (and gamma operator)
@ -336,8 +325,8 @@ else
params%rotation_BC = rotation_BC params%rotation_BC = rotation_BC
params%timeinc = timeinc params%timeinc = timeinc
call SNESSolve(snes,PETSC_NULL_OBJECT,solution_vec,ierr_psc) call SNESSolve(snes,PETSC_NULL_OBJECT,solution_vec,ierr)
call SNESGetConvergedReason(snes,reason,ierr_psc) call SNESGetConvergedReason(snes,reason,ierr)
AL_solution%termIll = terminallyIll AL_solution%termIll = terminallyIll
terminallyIll = .false. terminallyIll = .false.
@ -352,7 +341,7 @@ else
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief forms the AL residual vector !> @brief forms the AL residual vector
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
subroutine AL_formResidual(in,x_scal,f_scal,dummy,ierr_psc) subroutine AL_formResidual(in,x_scal,f_scal,dummy,ierr)
use numerics, only: & use numerics, only: &
itmax, & itmax, &
@ -386,15 +375,15 @@ else
PetscScalar, pointer :: residual_F(:,:,:,:,:), residual_F_lambda(:,:,:,:,:) PetscScalar, pointer :: residual_F(:,:,:,:,:), residual_F_lambda(:,:,:,:,:)
PetscInt :: iter, nfuncs PetscInt :: iter, nfuncs
PetscObject :: dummy PetscObject :: dummy
PetscErrorCode :: ierr_psc PetscErrorCode :: ierr
F => x_scal(:,:,1,:,:,:) F => x_scal(:,:,1,:,:,:)
F_lambda => x_scal(:,:,2,:,:,:) F_lambda => x_scal(:,:,2,:,:,:)
residual_F => f_scal(:,:,1,:,:,:) residual_F => f_scal(:,:,1,:,:,:)
residual_F_lambda => f_scal(:,:,2,:,:,:) residual_F_lambda => f_scal(:,:,2,:,:,:)
call SNESGetNumberFunctionEvals(snes,nfuncs,ierr_psc) call SNESGetNumberFunctionEvals(snes,nfuncs,ierr)
call SNESGetIterationNumber(snes,iter,ierr_psc) call SNESGetIterationNumber(snes,iter,ierr)
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! report begin of new iteration ! report begin of new iteration
@ -454,7 +443,7 @@ else
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief convergence check !> @brief convergence check
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
subroutine AL_converged(snes_local,it,xnorm,snorm,fnorm,reason,dummy,ierr_psc) subroutine AL_converged(snes_local,it,xnorm,snorm,fnorm,reason,dummy,ierr)
use numerics, only: & use numerics, only: &
itmax, & itmax, &
@ -466,12 +455,12 @@ else
implicit none implicit none
SNES snes_local SNES :: snes_local
PetscInt it PetscInt :: it
PetscReal xnorm, snorm, fnorm PetscReal :: xnorm, snorm, fnorm
SNESConvergedReason reason SNESConvergedReason :: reason
PetscObject dummy PetscObject :: dummy
PetscErrorCode ierr_psc PetscErrorCode ::ierr
logical :: Converged logical :: Converged
Converged = (it > itmin .and. & Converged = (it > itmin .and. &
@ -501,12 +490,12 @@ else
use DAMASK_spectral_Utilities, only: & use DAMASK_spectral_Utilities, only: &
Utilities_destroy Utilities_destroy
implicit none implicit none
PetscErrorCode ierr_psc PetscErrorCode :: ierr
call VecDestroy(solution_vec,ierr_psc) call VecDestroy(solution_vec,ierr)
call SNESDestroy(snes,ierr_psc) call SNESDestroy(snes,ierr)
call DMDestroy(da,ierr_psc) call DMDestroy(da,ierr)
call PetscFinalize(ierr_psc) call PetscFinalize(ierr)
call Utilities_destroy() call Utilities_destroy()
end subroutine AL_destroy end subroutine AL_destroy

View File

@ -6,7 +6,7 @@
!> @author Philip Eisenlohr, Max-Planck-Institut für Eisenforschung GmbH !> @author Philip Eisenlohr, Max-Planck-Institut für Eisenforschung GmbH
!> @brief Basic scheme PETSc solver !> @brief Basic scheme PETSc solver
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
module DAMASK_spectral_SolverBasicPETSC module DAMASK_spectral_SolverBasicPETSc
use prec, only: & use prec, only: &
pInt, & pInt, &
pReal pReal
@ -19,14 +19,8 @@ module DAMASK_spectral_SolverBasicPETSC
implicit none implicit none
#include <finclude/petscsys.h> #include <finclude/petscsys.h>
#include <finclude/petscvec.h>
#include <finclude/petscdmda.h> #include <finclude/petscdmda.h>
#include <finclude/petscis.h>
#include <finclude/petscmat.h>
#include <finclude/petscksp.h>
#include <finclude/petscpc.h>
#include <finclude/petscsnes.h> #include <finclude/petscsnes.h>
#include <finclude/petscvec.h90>
#include <finclude/petscdmda.h90> #include <finclude/petscdmda.h90>
#include <finclude/petscsnes.h90> #include <finclude/petscsnes.h90>
@ -44,9 +38,9 @@ module DAMASK_spectral_SolverBasicPETSC
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! PETSc data ! PETSc data
DM, private :: da DM, private :: da
SNES, private :: snes SNES, private :: snes
Vec, private :: solution_vec Vec, private :: solution_vec
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! common pointwise data ! common pointwise data
@ -107,47 +101,46 @@ subroutine BasicPETSC_init()
implicit none implicit none
integer(pInt) :: i,j,k integer(pInt) :: i,j,k
real(pReal), dimension(:,:,:,:,:), allocatable :: P real(pReal), dimension(:,:,:,:,:), allocatable :: P
PetscScalar, dimension(:,:,:,:), pointer :: F
PetscErrorCode :: ierr
PetscObject :: dummy
PetscErrorCode :: ierr_psc call Utilities_init()
PetscObject :: dummy
PetscMPIInt :: rank
PetscScalar, pointer :: xx_psc(:,:,:,:), F(:,:,:,:)
call Utilities_init() write(6,'(/,a)') ' <<<+- DAMASK_spectral_solverBasicPETSc init -+>>>'
write(6,'(a)') ' $Id: DAMASK_spectral_SolverBasicPETSC.f90 1654 2012-08-03 09:25:48Z MPIE\m.diehl $'
write(6,'(/,a)') ' <<<+- DAMASK_spectral_solverBasicPETSc init -+>>>'
write(6,'(a)') ' $Id: DAMASK_spectral_SolverBasicPETSC.f90 1654 2012-08-03 09:25:48Z MPIE\m.diehl $'
#include "compilation_info.f90" #include "compilation_info.f90"
write(6,'(a)') '' write(6,'(a)') ''
allocate (F_lastInc (3,3, res(1), res(2),res(3)), source = 0.0_pReal) allocate (F_lastInc (3,3, res(1), res(2),res(3)), source = 0.0_pReal)
allocate (Fdot (3,3, res(1), res(2),res(3)), source = 0.0_pReal) allocate (Fdot (3,3, res(1), res(2),res(3)), source = 0.0_pReal)
allocate (P (3,3, res(1), res(2),res(3)), source = 0.0_pReal) allocate (P (3,3, res(1), res(2),res(3)), source = 0.0_pReal)
allocate (coordinates( res(1), res(2),res(3),3), source = 0.0_pReal) allocate (coordinates( res(1), res(2),res(3),3), source = 0.0_pReal)
allocate (temperature( res(1), res(2),res(3)), source = 0.0_pReal) allocate (temperature( res(1), res(2),res(3)), source = 0.0_pReal)
!-------------------------------------------------------------------------------------------------- call SNESCreate(PETSC_COMM_WORLD,snes,ierr)
! PETSc Init CHKERRQ(ierr)
call PetscInitialize(PETSC_NULL_CHARACTER,ierr_psc) call DMDACreate3d(PETSC_COMM_WORLD, &
call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr_psc)
call SNESCreate(PETSC_COMM_WORLD,snes,ierr_psc)
call DMDACreate3d(PETSC_COMM_WORLD, &
DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_NONE, & DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_NONE, &
DMDA_STENCIL_BOX,res(1),res(2),res(3),PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE, & DMDA_STENCIL_BOX,res(1),res(2),res(3),PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE, &
9,1,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,da,ierr_psc) 9,1,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,da,ierr)
CHKERRQ(ierr)
call DMCreateGlobalVector(da,solution_vec,ierr_psc) call DMCreateGlobalVector(da,solution_vec,ierr)
call DMDASetLocalFunction(da,BasicPETSC_formResidual,ierr_psc) CHKERRQ(ierr)
call SNESSetDM(snes,da,ierr_psc) call DMDASetLocalFunction(da,BasicPETSC_formResidual,ierr)
call SNESSetConvergenceTest(snes,BasicPETSC_converged,dummy,PETSC_NULL_FUNCTION,ierr_psc) CHKERRQ(ierr)
call PetscOptionsInsertString(petsc_options,ierr_psc) call SNESSetDM(snes,da,ierr)
call SNESSetFromOptions(snes,ierr_psc) CHKERRQ(ierr)
call SNESSetConvergenceTest(snes,BasicPETSC_converged,dummy,PETSC_NULL_FUNCTION,ierr)
CHKERRQ(ierr)
call SNESSetFromOptions(snes,ierr)
CHKERRQ(ierr)
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! init fields ! init fields
call DMDAVecGetArrayF90(da,solution_vec,xx_psc,ierr_psc) call DMDAVecGetArrayF90(da,solution_vec,F,ierr)
F => xx_psc(0:8,:,:,:) CHKERRQ(ierr)
if (restartInc == 1_pInt) then ! no deformation (no restart) if (restartInc == 1_pInt) then ! no deformation (no restart)
F_lastInc = spread(spread(spread(math_I3,3,res(1)),4,res(2)),5,res(3)) ! initialize to identity F_lastInc = spread(spread(spread(math_I3,3,res(1)),4,res(2)),5,res(3)) ! initialize to identity
F = reshape(F_lastInc,[9,res(1),res(2),res(3)]) F = reshape(F_lastInc,[9,res(1),res(2),res(3)])
@ -158,8 +151,9 @@ subroutine BasicPETSC_init()
elseif (restartInc > 1_pInt) then ! using old values from file elseif (restartInc > 1_pInt) then ! using old values from file
if (debugRestart) write(6,'(a,i6,a)') 'Reading values of increment ',& if (debugRestart) write(6,'(a,i6,a)') 'Reading values of increment ',&
restartInc - 1_pInt,' from file' restartInc - 1_pInt,' from file'
flush(6)
call IO_read_jobBinaryFile(777,'convergedSpectralDefgrad',& call IO_read_jobBinaryFile(777,'convergedSpectralDefgrad',&
trim(getSolverJobName()),size(F_lastInc)) trim(getSolverJobName()),size(F))
read (777,rec=1) F read (777,rec=1) F
close (777) close (777)
call IO_read_jobBinaryFile(777,'convergedSpectralDefgrad_lastInc',& call IO_read_jobBinaryFile(777,'convergedSpectralDefgrad_lastInc',&
@ -175,9 +169,18 @@ subroutine BasicPETSC_init()
coordinates = 0.0 ! change it later!!! coordinates = 0.0 ! change it later!!!
endif endif
print*, 'F', shape(F)
call Utilities_constitutiveResponse(coordinates,F,F,temperature,0.0_pReal,P,C,P_av,.false.,math_I3) print*, 'F_lastInc', shape(F_lastInc)
call DMDAVecRestoreArrayF90(da,solution_vec,xx_psc,ierr_psc) print*, 'gfortran runs till here'
flush(6)
call Utilities_constitutiveResponse(coordinates,&
reshape(F(0:8,0:res(1)-1_pInt,0:res(2)-1_pInt,0:res(3)-1_pInt),[3,3,res(1),res(2),res(3)]),&
reshape(F(0:8,0:res(1)-1_pInt,0:res(2)-1_pInt,0:res(3)-1_pInt),[3,3,res(1),res(2),res(3)]),&
temperature,0.0_pReal,P,C,P_av,.false.,math_I3)
print*, 'gfortran does not reach this point'
flush(6)
call DMDAVecRestoreArrayF90(da,solution_vec,F,ierr)
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! reference stiffness ! reference stiffness
@ -237,8 +240,8 @@ subroutine BasicPETSC_init()
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! !
PetscScalar, pointer :: xx_psc(:,:,:,:), F(:,:,:,:) PetscScalar, pointer :: F(:,:,:,:)
PetscErrorCode :: ierr_psc PetscErrorCode :: ierr
SNESConvergedReason :: reason SNESConvergedReason :: reason
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
@ -253,12 +256,10 @@ subroutine BasicPETSC_init()
write (777,rec=1) C write (777,rec=1) C
close(777) close(777)
endif endif
call DMDAVecGetArrayF90(da,solution_vec,xx_psc,ierr_psc) call DMDAVecGetArrayF90(da,solution_vec,F,ierr)
F => xx_psc(0:8,:,:,:)
if ( cutBack) then if ( cutBack) then
F_aim = F_aim_lastInc F_aim = F_aim_lastInc
F = reshape(F_lastInc,[9,res(1),res(2),res(3)]) F = reshape(F_lastInc,[9,res(1),res(2),res(3)])
C = C_lastInc C = C_lastInc
else else
@ -288,7 +289,7 @@ else
F = reshape(Utilities_forwardField(timeinc,F_aim,F_lastInc,Fdot),[9,res(1),res(2),res(3)]) F = reshape(Utilities_forwardField(timeinc,F_aim,F_lastInc,Fdot),[9,res(1),res(2),res(3)])
call DMDAVecRestoreArrayF90(da,solution_vec,xx_psc,ierr_psc) call DMDAVecRestoreArrayF90(da,solution_vec,F,ierr)
call deformed_fft(res,geomdim,math_rotate_backward33(F_aim,rotation_BC),1.0_pReal,F_lastInc,coordinates) call deformed_fft(res,geomdim,math_rotate_backward33(F_aim,rotation_BC),1.0_pReal,F_lastInc,coordinates)
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
@ -302,8 +303,8 @@ else
params%rotation_BC = rotation_BC params%rotation_BC = rotation_BC
params%timeinc = timeinc params%timeinc = timeinc
call SNESSolve(snes,PETSC_NULL_OBJECT,solution_vec,ierr_psc) call SNESSolve(snes,PETSC_NULL_OBJECT,solution_vec,ierr)
call SNESGetConvergedReason(snes,reason,ierr_psc) call SNESGetConvergedReason(snes,reason,ierr)
basicPETSc_solution%termIll = terminallyIll basicPETSc_solution%termIll = terminallyIll
terminallyIll = .false. terminallyIll = .false.
BasicPETSC_solution%converged =.false. BasicPETSC_solution%converged =.false.
@ -316,7 +317,7 @@ else
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief forms the AL residual vector !> @brief forms the AL residual vector
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
subroutine BasicPETSC_formResidual(in,x_scal,f_scal,dummy,ierr_psc) subroutine BasicPETSC_formResidual(in,x_scal,f_scal,dummy,ierr)
use numerics, only: & use numerics, only: &
itmax, & itmax, &
@ -340,19 +341,15 @@ else
real(pReal), dimension(3,3) :: F_aim_lab_lastIter, F_aim_lab real(pReal), dimension(3,3) :: F_aim_lab_lastIter, F_aim_lab
DMDALocalInfo :: in(DMDA_LOCAL_INFO_SIZE) DMDALocalInfo, dimension(DMDA_LOCAL_INFO_SIZE) :: in
PetscScalar, target :: x_scal(3,3,XG_RANGE,YG_RANGE,ZG_RANGE) PetscScalar, dimension(3,3,XG_RANGE,YG_RANGE,ZG_RANGE) :: x_scal
PetscScalar, target :: f_scal(3,3,X_RANGE,Y_RANGE,Z_RANGE) PetscScalar, dimension(3,3,X_RANGE,Y_RANGE,Z_RANGE):: f_scal
PetscScalar, pointer :: F(:,:,:,:,:)
PetscScalar, pointer :: residual_F(:,:,:,:,:)
PetscInt :: iter, nfuncs PetscInt :: iter, nfuncs
PetscObject :: dummy PetscObject :: dummy
PetscErrorCode :: ierr_psc PetscErrorCode :: ierr
F => x_scal(:,:,:,:,:)
residual_F => f_scal(:,:,:,:,:)
call SNESGetNumberFunctionEvals(snes,nfuncs,ierr_psc) call SNESGetNumberFunctionEvals(snes,nfuncs,ierr)
call SNESGetIterationNumber(snes,iter,ierr_psc) call SNESGetIterationNumber(snes,iter,ierr)
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! report begin of new iteration ! report begin of new iteration
@ -364,8 +361,8 @@ else
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! evaluate constitutive response ! evaluate constitutive response
call Utilities_constitutiveResponse(coordinates,F_lastInc,F,temperature,params%timeinc, & call Utilities_constitutiveResponse(coordinates,F_lastInc,x_scal,temperature,params%timeinc, &
residual_F,C,P_av,ForwardData,params%rotation_BC) f_scal,C,P_av,ForwardData,params%rotation_BC)
ForwardData = .false. ForwardData = .false.
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
@ -377,7 +374,7 @@ else
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! updated deformation gradient using fix point algorithm of basic scheme ! updated deformation gradient using fix point algorithm of basic scheme
field_real = 0.0_pReal field_real = 0.0_pReal
field_real(1:res(1),1:res(2),1:res(3),1:3,1:3) = reshape(residual_F,[res(1),res(2),res(3),3,3],& field_real(1:res(1),1:res(2),1:res(3),1:3,1:3) = reshape(f_scal,[res(1),res(2),res(3),3,3],&
order=[4,5,1,2,3]) ! field real has a different order order=[4,5,1,2,3]) ! field real has a different order
call Utilities_FFTforward() call Utilities_FFTforward()
err_div = Utilities_divergenceRMS() err_div = Utilities_divergenceRMS()
@ -386,14 +383,14 @@ else
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! constructing residual ! constructing residual
residual_F = reshape(field_real(1:res(1),1:res(2),1:res(3),1:3,1:3),shape(F),order=[3,4,5,1,2]) f_scal = reshape(field_real(1:res(1),1:res(2),1:res(3),1:3,1:3),shape(x_scal),order=[3,4,5,1,2])
write(6,'(/,a)') '==========================================================================' write(6,'(/,a)') '=========================================================================='
end subroutine BasicPETSC_formResidual end subroutine BasicPETSC_formResidual
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief convergence check !> @brief convergence check
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
subroutine BasicPETSC_converged(snes_local,it,xnorm,snorm,fnorm,reason,dummy,ierr_psc) subroutine BasicPETSC_converged(snes_local,it,xnorm,snorm,fnorm,reason,dummy,ierr)
use numerics, only: & use numerics, only: &
itmax, & itmax, &
@ -409,12 +406,12 @@ else
implicit none implicit none
SNES snes_local SNES :: snes_local
PetscInt it PetscInt :: it
PetscReal xnorm, snorm, fnorm PetscReal :: xnorm, snorm, fnorm
SNESConvergedReason reason SNESConvergedReason :: reason
PetscObject dummy PetscObject :: dummy
PetscErrorCode ierr_psc PetscErrorCode :: ierr
logical :: Converged logical :: Converged
real(pReal) :: pAvgDivL2 real(pReal) :: pAvgDivL2
@ -446,12 +443,12 @@ else
use DAMASK_spectral_Utilities, only: & use DAMASK_spectral_Utilities, only: &
Utilities_destroy Utilities_destroy
implicit none implicit none
PetscErrorCode ierr_psc PetscErrorCode :: ierr
call VecDestroy(solution_vec,ierr_psc) call VecDestroy(solution_vec,ierr)
call SNESDestroy(snes,ierr_psc) call SNESDestroy(snes,ierr)
call DMDestroy(da,ierr_psc) call DMDestroy(da,ierr)
call PetscFinalize(ierr_psc) call PetscFinalize(ierr)
call Utilities_destroy() call Utilities_destroy()
end subroutine BasicPETSC_destroy end subroutine BasicPETSC_destroy

View File

@ -657,6 +657,7 @@ subroutine utilities_constitutiveResponse(coordinates,F_lastInc,F,temperature,ti
endif endif
if (DebugGeneral) write(6,*) 'collect mode: ', collectMode,' calc mode: ', calcMode if (DebugGeneral) write(6,*) 'collect mode: ', collectMode,' calc mode: ', calcMode
flush(6)
ielem = 0_pInt ielem = 0_pInt
do k = 1_pInt, res(3); do j = 1_pInt, res(2); do i = 1_pInt, res(1) do k = 1_pInt, res(3); do j = 1_pInt, res(2); do i = 1_pInt, res(1)

View File

@ -25,6 +25,9 @@ module numerics
use prec, only: pInt, pReal use prec, only: pInt, pReal
implicit none implicit none
#ifdef PETSc
#include <finclude/petscsys.h>
#endif
character(len=64), parameter, private :: & character(len=64), parameter, private :: &
numerics_configFile = 'numerics.config' !< name of configuration file numerics_configFile = 'numerics.config' !< name of configuration file
@ -132,6 +135,9 @@ subroutine numerics_init
implicit none implicit none
#ifdef Marc ! use the non F90 standard include file because some versions of Marc and Abaqus crash when using the module #ifdef Marc ! use the non F90 standard include file because some versions of Marc and Abaqus crash when using the module
!$ include "omp_lib.h" !$ include "omp_lib.h"
#endif
#ifdef PETSc
PetscErrorCode :: ierr
#endif #endif
integer(pInt), parameter :: fileunit = 300_pInt ,& integer(pInt), parameter :: fileunit = 300_pInt ,&
maxNchunks = 2_pInt maxNchunks = 2_pInt
@ -326,7 +332,11 @@ subroutine numerics_init
fftw_planner_flag = 32_pInt fftw_planner_flag = 32_pInt
end select end select
#endif #endif
#ifdef PETSc
call PetscOptionsInsertString(petsc_options,ierr)
write(6,'(a)') ' Initializing PETSc'
CHKERRQ(ierr)
#endif
!* writing parameters to output file !* writing parameters to output file
write(6,'(a24,1x,es8.1)') ' relevantStrain: ',relevantStrain write(6,'(a24,1x,es8.1)') ' relevantStrain: ',relevantStrain