new solver is now compiling without a PETSc installation, however only the plain basic solver is available then.

also removed reporting of PETSc related variables in the case it is not installed/found
This commit is contained in:
Martin Diehl 2012-08-28 19:19:47 +00:00
parent 0f64289d75
commit 5e9e8497e6
8 changed files with 63 additions and 31 deletions

View File

@ -61,8 +61,10 @@ program DAMASK_spectral_Driver
debugGeneral
use DAMASK_spectral_SolverBasic
#ifdef PETSc
use DAMASK_spectral_SolverBasicPETSC
use DAMASK_spectral_SolverAL
#endif
implicit none
@ -276,13 +278,17 @@ program DAMASK_spectral_Driver
case (DAMASK_spectral_SolverBasic_label)
call basic_init()
#ifdef PETSc
case (DAMASK_spectral_SolverBasicPETSC_label)
call BasicPETSC_init()
case (DAMASK_spectral_SolverAL_label)
call AL_init()
#endif
case default
call IO_error(error_ID = 891, ext_msg = trim(myspectralsolver))
end select
!--------------------------------------------------------------------------------------------------
! write header of output file
if (appendToOutFile) then
@ -352,7 +358,7 @@ program DAMASK_spectral_Driver
write(6,'(a)') '##################################################################'
write(6,'(A,I5.5,A,es12.5)') 'Increment ', totalIncsCounter, ' Time ',time
select case (myspectralsolver)
select case(myspectralsolver)
case (DAMASK_spectral_SolverBasic_label)
solres = basic_solution (&
@ -361,7 +367,7 @@ program DAMASK_spectral_Driver
F_BC = loadCases(currentLoadCase)%deformation, &
temperature_bc = loadCases(currentLoadCase)%temperature, &
rotation_BC = loadCases(currentLoadCase)%rotation)
#ifdef PETSc
case (DAMASK_spectral_SolverBasicPETSC_label)
solres = BasicPETSC_solution (&
guessmode,timeinc,timeinc_old, &
@ -377,7 +383,7 @@ program DAMASK_spectral_Driver
F_BC = loadCases(currentLoadCase)%deformation, &
temperature_bc = loadCases(currentLoadCase)%temperature, &
rotation_BC = loadCases(currentLoadCase)%rotation)
#endif
end select
write(6,'(a)') ''
@ -406,13 +412,13 @@ program DAMASK_spectral_Driver
case (DAMASK_spectral_SolverBasic_label)
call basic_destroy()
#ifdef PETSc
case (DAMASK_spectral_SolverBasicPETSC_label)
call BasicPETSC_destroy()
case (DAMASK_spectral_SolverAL_label)
call AL_destroy()
#endif
end select
write(6,'(a)') ''

View File

@ -18,7 +18,6 @@ module DAMASK_spectral_SolverAL
solutionState
implicit none
#ifdef PETSC
#include <finclude/petscsys.h>
#include <finclude/petscvec.h>
#include <finclude/petscdmda.h>
@ -30,7 +29,6 @@ module DAMASK_spectral_SolverAL
#include <finclude/petscvec.h90>
#include <finclude/petscdmda.h90>
#include <finclude/petscsnes.h90>
#endif
character (len=*), parameter, public :: &
DAMASK_spectral_SolverAL_label = 'AL'
@ -80,7 +78,7 @@ module DAMASK_spectral_SolverAL
!--------------------------------------------------------------------------------------------------
subroutine AL_init()
use, intrinsic :: iso_fortran_env ! to get compiler_version and compiler_options (at least for gfortran >4.6 at the moment)
use IO, only: &
IO_read_JobBinaryFile, &
IO_write_JobBinaryFile

View File

@ -5,6 +5,8 @@
!> @author Martin Diehl, Max-Planck-Institut für Eisenforschung GmbH
!> @author Philip Eisenlohr, Max-Planck-Institut für Eisenforschung GmbH
!> @brief Basic scheme solver
!> @details this solver follows closely the original large strain formulation presented by
!> Suquet. The iterative procedure is solved using a fix-point iteration
!--------------------------------------------------------------------------------------------------
module DAMASK_spectral_SolverBasic
use prec, only: &
@ -35,7 +37,7 @@ module DAMASK_spectral_SolverBasic
real(pReal), private,dimension(3,3,3,3) :: &
C = 0.0_pReal
contains
contains
!--------------------------------------------------------------------------------------------------
!> @brief allocates all neccessary fields and fills them with data, potentially from restart info
@ -76,7 +78,7 @@ subroutine basic_init()
#include "compilation_info.f90"
write(6,'(a)') ''
allocate (F ( 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 (P ( 3,3,res(1), res(2),res(3)), source = 0.0_pReal)

View File

@ -18,7 +18,6 @@ module DAMASK_spectral_SolverBasicPETSC
solutionState
implicit none
#ifdef PETSC
#include <finclude/petscsys.h>
#include <finclude/petscvec.h>
#include <finclude/petscdmda.h>
@ -30,7 +29,6 @@ module DAMASK_spectral_SolverBasicPETSC
#include <finclude/petscvec.h90>
#include <finclude/petscdmda.h90>
#include <finclude/petscsnes.h90>
#endif
character (len=*), parameter, public :: &
DAMASK_spectral_SolverBasicPETSC_label = 'basicpetsc'

View File

@ -627,6 +627,9 @@ real(pReal) function Utilities_getFilter(k)
*(1.0_pReal + cos(pi*k(2)/res(2))) &
*(1.0_pReal + cos(pi*k(1)/res(1)))/8.0_pReal
case default
call IO_error(error_ID = 892_pInt, ext_msg = trim(myfilter))
end select
end function Utilities_getFilter

View File

@ -1420,8 +1420,11 @@ subroutine IO_error(error_ID,e,i,g,ext_msg)
msg = 'mismatch of microstructure count and a*b*c in geom file'
case (890_pInt)
msg = 'invalid input for regridding'
case (891_pInt)
msg = 'unknown solver type selected'
case (892_pInt)
msg = 'unknown filter type selected'
!* Error messages related to parsing of Abaqus input file
case (900_pInt)
@ -1501,6 +1504,8 @@ subroutine IO_warning(warning_ID,e,i,g,ext_msg)
msg = 'could not get $DAMASK_NUM_THREADS'
case (40_pInt)
msg = 'Found Spectral solver parameter '
case (41_pInt)
msg = 'Found PETSc solver parameter '
case (47_pInt)
msg = 'No valid parameter for FFTW given, using FFTW_PATIENT'
case (101_pInt)

View File

@ -99,7 +99,7 @@ LIB_DIRS +=-L$(FFTWROOT)/lib
ifeq "$(SOLVER)" "NEW"
ifdef PETSC_DIR
include ${PETSC_DIR}/conf/variables
INCLUDE_DIRS +=${PETSC_FC_INCLUDES} -DPETSC
INCLUDE_DIRS +=${PETSC_FC_INCLUDES} -DPETSc
LIBRARIES +=${PETSC_WITH_EXTERNAL_LIB}
endif
endif
@ -277,15 +277,20 @@ COMPILED_FILES = prec.o DAMASK_spectral_interface.o IO.o numerics.o debug.o math
constitutive_titanmod.o constitutive_nonlocal.o constitutive_none.o constitutive.o crystallite.o \
homogenization_RGC.o homogenization_isostrain.o homogenization.o CPFEM.o
ifeq "$(SOLVER)" "NEW"
COMPILED_FILES += DAMASK_spectral_Utilities.o DAMASK_spectral_SolverBasic.o DAMASK_spectral_SolverAL.o DAMASK_spectral_SolverBasicPETSC.o
ifdef PETSC_DIR
PETSC_FILES = DAMASK_spectral_SolverAL.o DAMASK_spectral_SolverBasicPETSC.o
endif
COMPILED_FILES += DAMASK_spectral_Utilities.o DAMASK_spectral_SolverBasic.o $(PETSC_FILES)
DAMASK_spectral.exe: DAMASK_spectral_Driver.o
$(PREFIX) $(COMPILERNAME) $(OPENMP_FLAG_$(F90)) $(OPTIMIZATION_$(MAXOPTI)_$(F90)) $(STANDARD_CHECK_$(F90)) \
-o DAMASK_spectral.exe DAMASK_spectral_Driver.o \
$(COMPILED_FILES) $(LIB_DIRS) $(LIBRARIES) $(SUFFIX)
DAMASK_spectral_Driver.o: DAMASK_spectral_Driver.f90 DAMASK_spectral_SolverBasic.o DAMASK_spectral_SolverAL.o DAMASK_spectral_SolverBasicPETSC.o
DAMASK_spectral_Driver.o: DAMASK_spectral_Driver.f90 DAMASK_spectral_SolverBasic.o $(PETSC_FILES)
$(PREFIX) $(COMPILERNAME) $(COMPILE_MAXOPTI) -c DAMASK_spectral_Driver.f90 $(SUFFIX)
DAMASK_spectral_SolverAL.o: DAMASK_spectral_SolverAL.f90\

View File

@ -242,10 +242,6 @@ subroutine numerics_init
err_stress_tolrel = IO_floatValue(line,positions,2_pInt)
case ('err_stress_tolabs')
err_stress_tolabs = IO_floatValue(line,positions,2_pInt)
case ('err_f_tol')
err_f_tol = IO_floatValue(line,positions,2_pInt)
case ('err_p_tol')
err_p_tol = IO_floatValue(line,positions,2_pInt)
case ('itmax')
itmax = IO_intValue(line,positions,2_pInt)
case ('itmin')
@ -256,23 +252,34 @@ subroutine numerics_init
fftw_timelimit = IO_floatValue(line,positions,2_pInt)
case ('fftw_plan_mode')
fftw_plan_mode = IO_stringValue(line,positions,2_pInt)
case ('myspectralsolver')
myspectralsolver = IO_stringValue(line,positions,2_pInt)
case ('myfilter')
myfilter = IO_stringValue(line,positions,2_pInt)
case ('petsc_options')
petsc_options = trim(line(positions(4):))
case ('rotation_tol')
rotation_tol = IO_floatValue(line,positions,2_pInt)
case ('divergence_correction')
divergence_correction = IO_intValue(line,positions,2_pInt) > 0_pInt
case ('update_gamma')
update_gamma = IO_intValue(line,positions,2_pInt) > 0_pInt
#ifdef PETSc
case ('petsc_options')
petsc_options = trim(line(positions(4):))
case ('myspectralsolver')
myspectralsolver = IO_stringValue(line,positions,2_pInt)
case ('err_f_tol')
err_f_tol = IO_floatValue(line,positions,2_pInt)
case ('err_p_tol')
err_p_tol = IO_floatValue(line,positions,2_pInt)
#endif
#ifndef PETSc
case ('myspectralsolver', 'petsc_options','err_f_tol', 'err_p_tol')
call IO_warning(41_pInt,ext_msg=tag)
#endif
#endif
#ifndef Spectral
case ('err_div_tol','err_stress_tolrel','err_stress_tolabs',&
'itmax', 'itmin','memory_efficient','fftw_timelimit','fftw_plan_mode','myspectralsolver', &
'rotation_tol','divergence_correction','update_gamma','petsc_options','myfilter')
'rotation_tol','divergence_correction','update_gamma','petsc_options','myfilter', &
'err_f_tol', 'err_p_tol')
call IO_warning(40_pInt,ext_msg=tag)
#endif
case default
@ -354,8 +361,7 @@ subroutine numerics_init
write(6,'(a24,1x,es8.1)') ' err_div_tol: ',err_div_tol
write(6,'(a24,1x,es8.1)') ' err_stress_tolrel: ',err_stress_tolrel
write(6,'(a24,1x,es8.1)') ' err_stress_tolabs: ',err_stress_tolabs
write(6,'(a24,1x,es8.1)') ' err_f_tol: ',err_f_tol
write(6,'(a24,1x,es8.1)') ' err_p_tol: ',err_p_tol
write(6,'(a24,1x,i8)') ' itmax: ',itmax
write(6,'(a24,1x,i8)') ' itmin: ',itmin
write(6,'(a24,1x,L8)') ' memory_efficient: ',memory_efficient
@ -365,13 +371,18 @@ subroutine numerics_init
write(6,'(a24,1x,es8.1)') ' fftw_timelimit: ',fftw_timelimit
endif
write(6,'(a24,1x,a)') ' fftw_plan_mode: ',trim(fftw_plan_mode)
write(6,'(a24,1x,a)') ' myspectralsolver: ',trim(myspectralsolver)
write(6,'(a24,1x,a)') ' myfilter: ',trim(myfilter)
write(6,'(a24,1x,a)') ' PetSc_options: ',trim(petsc_options)
write(6,'(a24,1x,i8)') ' fftw_planner_flag: ',fftw_planner_flag
write(6,'(a24,1x,es8.1)') ' rotation_tol: ',rotation_tol
write(6,'(a24,1x,L8,/)') ' divergence_correction: ',divergence_correction
write(6,'(a24,1x,L8,/)') ' update_gamma: ',update_gamma
#ifdef PETSc
write(6,'(a24,1x,es8.1)') ' err_f_tol: ',err_f_tol
write(6,'(a24,1x,es8.1)') ' err_p_tol: ',err_p_tol
write(6,'(a24,1x,a)') ' myspectralsolver: ',trim(myspectralsolver)
write(6,'(a24,1x,a)') ' PETSc_options: ',trim(petsc_options)
#endif
#endif
!* sanity check
@ -425,6 +436,10 @@ subroutine numerics_init
if (itmin > itmax .or. itmin < 1_pInt) call IO_error(301_pInt,ext_msg='itmin')
if (update_gamma .and. &
.not. memory_efficient) call IO_error(error_ID = 847_pInt)
#ifdef PETSc
if (err_f_tol <= 0.0_pReal) call IO_error(301_pInt,ext_msg='err_f_tol')
if (err_p_tol <= 0.0_pReal) call IO_error(301_pInt,ext_msg='err_p_tol')
#endif
#endif
if (fixedSeed <= 0_pInt) then
write(6,'(a,/)') ' Random is random!'