made analytic tangent default for spectral solver

removed possibility to compile without PETSc (a lot of effort for little/no use)
This commit is contained in:
Martin Diehl 2014-10-01 12:29:12 +00:00
parent 1996105c7f
commit baa08d8155
4 changed files with 17 additions and 70 deletions

View File

@ -60,11 +60,9 @@ program DAMASK_spectral_Driver
tSolutionState, &
cutBack
use DAMASK_spectral_SolverBasic
#ifdef PETSc
use DAMASK_spectral_SolverBasicPETSC
use DAMASK_spectral_SolverAL
use DAMASK_spectral_SolverPolarisation
#endif
implicit none
type tLoadCase
@ -316,7 +314,6 @@ program DAMASK_spectral_Driver
select case (spectral_solver)
case (DAMASK_spectral_SolverBasic_label)
call basic_init(loadCases(1)%temperature)
#ifdef PETSc
case (DAMASK_spectral_SolverBasicPETSc_label)
call basicPETSc_init(loadCases(1)%temperature)
case (DAMASK_spectral_SolverAL_label)
@ -327,7 +324,6 @@ program DAMASK_spectral_Driver
if(iand(debug_level(debug_spectral),debug_levelBasic)/= 0) &
call IO_warning(42_pInt, ext_msg='debug Divergence')
call Polarisation_init(loadCases(1)%temperature)
#endif
case default
call IO_error(error_ID = 891, ext_msg = trim(spectral_solver))
end select
@ -410,7 +406,6 @@ program DAMASK_spectral_Driver
! forward solution
select case(spectral_solver)
case (DAMASK_spectral_SolverBasic_label)
#ifdef PETSc
case (DAMASK_spectral_SolverBasicPETSC_label)
call BasicPETSC_forward (&
guess,timeinc,timeIncOld,remainingLoadCaseTime, &
@ -431,7 +426,6 @@ program DAMASK_spectral_Driver
P_BC = loadCases(currentLoadCase)%P, &
F_BC = loadCases(currentLoadCase)%deformation, &
rotation_BC = loadCases(currentLoadCase)%rotation)
#endif
end select
!--------------------------------------------------------------------------------------------------
@ -461,7 +455,6 @@ 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 (&
incInfo,guess,timeinc,timeIncOld,remainingLoadCaseTime, &
@ -486,7 +479,6 @@ program DAMASK_spectral_Driver
temperature_bc = loadCases(currentLoadCase)%temperature, &
rotation_BC = loadCases(currentLoadCase)%rotation, &
density = loadCases(currentLoadCase)%density)
#endif
end select
!--------------------------------------------------------------------------------------------------
@ -549,14 +541,12 @@ program DAMASK_spectral_Driver
select case (spectral_solver)
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()
case (DAMASK_spectral_SolverPolarisation_label)
call Polarisation_destroy()
#endif
end select
!--------------------------------------------------------------------------------------------------

View File

@ -1717,8 +1717,6 @@ subroutine IO_warning(warning_ID,el,ip,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 (42_pInt)
msg = 'parameter has no effect'
case (43_pInt)

View File

@ -26,18 +26,11 @@ else
include /etc/damask.conf
endif
ifdef PETSC_DIR
include $(PETSC_DIR)/conf/variables
INCLUDE_DIRS :=$(PETSC_FC_INCLUDES) -DPETSc -I../lib
LIBRARIES :=$(PETSC_WITH_EXTERNAL_LIB)
COMPILERNAME ?= $(FC)
LINKERNAME ?= $(FLINKER)
else
COMPILERNAME ?= $(F90)
LINKERNAME ?= $(F90)
INCLUDE_DIRS :=-I../lib
LIBRARIES :=-lfftw3
endif
LIB_DIRS :=-L$(FFTW_ROOT)/lib64 -L$(FFTW_ROOT)/lib
RUN_PATH :=-Wl,-rpath,$(FFTW_ROOT)/lib64,-rpath,$(FFTW_ROOT)/lib
@ -90,27 +83,6 @@ LIBRARIES +=-liomp5
endif
endif
ifndef PETSC_DIR #petsc provides linking options agains selected blas/lapack already
ifneq "x$(IMKL_ROOT)" "x"
LIB_DIRS +=-L$(IMKL_ROOT)/lib/intel64
RUN_PATH :=$(RUN_PATH),-rpath,$(IMKL_ROOT)/lib/intel64
INCLUDE_DIRS +=-I$(IMKL_ROOT)/include
LIBRARIES +=-lmkl_$(IMKL_COMPILER_$(F90))_lp64 -lmkl_core -lmkl_sequential -lm
else
ifneq "x$(ACML_ROOT)" "x"
LIB_DIRS +=-L$(ACML_ROOT)/$(F90)64/lib
RUN_PATH :=$(RUN_PATH),-rpath,$(ACML_ROOT)/$(F90)64/lib
LIBRARIES +=-lacml
else
ifneq "x$(LAPACK_ROOT)" "x"
LIB_DIRS +=-L$(LAPACK_ROOT)/lib64 -L$(LAPACK_ROOT)/lib
RUN_PATH :=$(RUN_PATH),-rpath,$(LAPACK_ROOT)/lib64,-rpath,$(LAPACK_ROOT)/lib
LIBRARIES +=-llapack
endif
endif
endif
endif
#hdf5
ifeq "$(HDF5)" "ON"
LIBRARIES +=-lhdf5hl_fortran -lhdf5_hl -lhdf5_fortran -lhdf5
@ -380,20 +352,18 @@ SPECTRAL_FILES = prec.o DAMASK_interface.o IO.o libs.o numerics.o debug.o math.o
FEsolving.o mesh.o material.o lattice.o \
$(DAMAGE_FILES) $(THERMAL_FILES) $(CONSTITUTIVE_FILES) \
crystallite.o $(HOMOGENIZATION_FILES) CPFEM.o \
DAMASK_spectral_utilities.o DAMASK_spectral_solverBasic.o
ifdef PETSC_DIR
PETSC_FILES = DAMASK_spectral_solverAL.o \
DAMASK_spectral_solverBasicPETSc.o \
DAMASK_spectral_solverPolarisation.o
SPECTRAL_FILES += $(PETSC_FILES)
endif
DAMASK_spectral_utilities.o DAMASK_spectral_solverBasic.o \
DAMASK_spectral_solverAL.o DAMASK_spectral_solverBasicPETSc.o DAMASK_spectral_solverPolarisation.o
DAMASK_spectral.exe: DAMASK_spectral_driver.o
$(PREFIX) $(LINKERNAME) $(OPENMP_FLAG_$(F90)) $(LINK_OPTIONS_$(F90)) $(STANDARD_CHECK_$(F90)) $(OPTIMIZATION_$(MAXOPTI)_$(F90)) \
-o DAMASK_spectral.exe DAMASK_spectral_driver.o \
$(SPECTRAL_FILES) $(LIBRARIES) $(LIB_DIRS) $(RUN_PATH) $(SUFFIX)
DAMASK_spectral_driver.o: DAMASK_spectral_driver.f90 DAMASK_spectral_solverBasic.o $(PETSC_FILES)
DAMASK_spectral_driver.o: DAMASK_spectral_driver.f90 DAMASK_spectral_solverBasic.o \
DAMASK_spectral_solverAL.o \
DAMASK_spectral_solverBasicPETSc.o \
DAMASK_spectral_solverPolarisation.o
$(PREFIX) $(COMPILERNAME) $(COMPILE_MAXOPTI) -c DAMASK_spectral_driver.f90 $(SUFFIX)
DAMASK_spectral_solverAL.o: DAMASK_spectral_solverAL.f90 \

View File

@ -61,11 +61,15 @@ module numerics
maxVolDiscr_RGC = 1.0e-5_pReal, & !< threshold of maximum volume discrepancy allowed
volDiscrMod_RGC = 1.0e+12_pReal, & !< stiffness of RGC volume discrepancy (zero = without volume discrepancy constraint)
volDiscrPow_RGC = 5.0_pReal, & !< powerlaw penalty for volume discrepancy
charLength = 1.0_pReal !< characteristic length scale for gradient problems
charLength = 1.0_pReal !< characteristic length scale for gradient problems
logical, protected, public :: &
analyticJaco = .false., & !< use analytic Jacobian or perturbation, Default .false.: calculate Jacobian using perturbations
#ifdef Spectral
analyticJaco = .true., & !< use analytic Jacobian or perturbation, Default for Spectral solver .true.:
#else
analyticJaco = .false., & !< use analytic Jacobian or perturbation, Default .false.: calculate Jacobian using perturbations
#endif
usePingPong = .true., &
numerics_timeSyncing = .false. !< flag indicating if time synchronization in crystallite is used for nonlocal plasticity
numerics_timeSyncing = .false. !< flag indicating if time synchronization in crystallite is used for nonlocal plasticity
!--------------------------------------------------------------------------------------------------
! spectral parameters:
@ -182,13 +186,11 @@ subroutine numerics_init
#if defined(Spectral) || defined(FEM)
!$ use OMP_LIB, only: omp_set_num_threads ! Use the standard conforming module file for omp if using the spectral solver
#endif
implicit none
#ifndef Spectral
#ifndef FEM
#else
implicit none
!$ include "omp_lib.h" ! use the not F90 standard conforming include file to prevent crashes with some versions of MSC.Marc
#endif
#endif
integer(pInt), parameter :: FILEUNIT = 300_pInt ,&
maxNchunks = 2_pInt
!$ integer :: gotDAMASK_NUM_THREADS = 1
@ -360,7 +362,6 @@ subroutine numerics_init
divergence_correction = IO_intValue(line,positions,2_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 ('spectralsolver','myspectralsolver')
@ -373,14 +374,7 @@ subroutine numerics_init
polarAlpha = IO_floatValue(line,positions,2_pInt)
case ('polarbeta')
polarBeta = IO_floatValue(line,positions,2_pInt)
#endif
#ifndef PETSc
case ('spectralsolver', 'myspectralsolver', 'petsc_options', &
'err_curl_tolabs','err_curl_tolrel','polaralpha','polarBeta') ! found PETSc parameter, but compiled without PETSc
call IO_warning(41_pInt,ext_msg=tag)
#endif
#endif
#ifndef Spectral
#else
case ('err_div_tolabs','err_div_tolrel','err_stress_tolrel','err_stress_tolabs',& ! found spectral parameter for FEM build
'itmax', 'itmin','memory_efficient','fftw_timelimit','fftw_plan_mode', &
'divergence_correction','update_gamma','spectralfilter','myfilter', &
@ -418,8 +412,7 @@ subroutine numerics_init
damageorder = IO_intValue(line,positions,2_pInt)
case ('petsc_optionsfem')
petsc_optionsFEM = trim(line(positions(4):))
#endif
#ifndef FEM
#else
case ('err_struct_tolabs','err_struct_tolrel','err_thermal_tol','err_damage_tol','residualstiffness',& ! found spectral parameter for FEM build
'itmaxfem', 'itminfem','maxcutbackfem','integrationorder','structorder','thermalorder', &
'damageorder','petsc_optionsfem')
@ -529,7 +522,6 @@ subroutine numerics_init
write(6,'(a24,1x,es8.1)') ' err_stress_tolRel: ',err_stress_tolRel
write(6,'(a24,1x,es8.1)') ' err_div_tolAbs: ',err_div_tolAbs
write(6,'(a24,1x,es8.1)') ' err_div_tolRel: ',err_div_tolRel
#ifdef PETSc
write(6,'(a24,1x,es8.1)') ' err_curl_tolAbs: ',err_curl_tolAbs
write(6,'(a24,1x,es8.1)') ' err_curl_tolRel: ',err_curl_tolRel
write(6,'(a24,1x,es8.1)') ' polarAlpha: ',polarAlpha
@ -537,7 +529,6 @@ subroutine numerics_init
write(6,'(a24,1x,a)') ' spectral solver: ',trim(spectral_solver)
write(6,'(a24,1x,a)') ' PETSc_options: ',trim(petsc_options)
#endif
#endif
!--------------------------------------------------------------------------------------------------
! spectral parameters
@ -610,7 +601,6 @@ subroutine numerics_init
if (err_stress_tolabs <= 0.0_pReal) call IO_error(301_pInt,ext_msg='err_stress_tolAbs')
if (err_div_tolRel <= 0.0_pReal) call IO_error(301_pInt,ext_msg='err_div_tolRel')
if (err_div_tolAbs <= 0.0_pReal) call IO_error(301_pInt,ext_msg='err_div_tolAbs')
#ifdef PETSc
if (err_curl_tolRel <= 0.0_pReal) call IO_error(301_pInt,ext_msg='err_curl_tolRel')
if (err_curl_tolAbs <= 0.0_pReal) call IO_error(301_pInt,ext_msg='err_curl_tolAbs')
if (polarAlpha <= 0.0_pReal .or. &
@ -618,7 +608,6 @@ subroutine numerics_init
if (polarBeta < 0.0_pReal .or. &
polarBeta > 2.0_pReal) call IO_error(301_pInt,ext_msg='polarBeta')
#endif
#endif
#ifdef FEM
if (itmaxFEM <= 1_pInt) call IO_error(301_pInt,ext_msg='itmaxFEM')
if (itminFEM > itmaxFEM .or. &