! Copyright 2012 Max-Planck-Institut für Eisenforschung GmbH ! ! This file is part of DAMASK, ! the Düsseldorf Advanced Material Simulation Kit. ! ! DAMASK is free software: you can redistribute it and/or modify ! it under the terms of the GNU General Public License as published by ! the Free Software Foundation, either version 3 of the License, or ! (at your option) any later version. ! ! DAMASK is distributed in the hope that it will be useful, ! but WITHOUT ANY WARRANTY; without even the implied warranty of ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the ! GNU General Public License for more details. ! ! You should have received a copy of the GNU General Public License ! along with DAMASK. If not, see . ! !################################################################################################## !* $Id$ !################################################################################################## ! Material subroutine for BVP solution using spectral method ! ! Run 'DAMASK_spectral.exe --help' to get usage hints ! ! written by P. Eisenlohr, ! F. Roters, ! L. Hantcherli, ! W.A. Counts, ! D.D. Tjahjanto, ! C. Kords, ! M. Diehl, ! R. Lebensohn ! ! MPI fuer Eisenforschung, Duesseldorf #include "spectral_quit.f90" program DAMASK_spectral use, intrinsic :: iso_fortran_env ! to get compiler_version and compiler_options (at least for gfortran >4.6 at the moment) use DAMASK_interface, only: & DAMASK_interface_init, & loadCaseFile, & geometryFile, & getSolverWorkingDirectoryName, & getSolverJobName, & appendToOutFile use prec, only: & pInt, & pReal, & DAMASK_NaN use IO, only: & IO_isBlank, & IO_open_file, & IO_stringPos, & IO_stringValue, & IO_floatValue, & IO_intValue, & IO_error, & IO_lc, & IO_read_jobBinaryFile, & IO_write_jobBinaryFile use debug, only: & debug_level, & debug_spectral, & debug_levelBasic, & debug_spectralDivergence, & debug_spectralRestart, & debug_spectralFFTW, & debug_reset, & debug_info use math use mesh, only : & mesh_spectral_getResolution, & mesh_spectral_getDimension, & mesh_spectral_getHomogenization use CPFEM, only: & CPFEM_general, & CPFEM_initAll use FEsolving, only: & restartWrite, & restartInc use numerics, only: & err_div_tol, & err_stress_tolrel, & err_stress_tolabs, & rotation_tol, & itmax,& itmin, & memory_efficient, & divergence_correction, & DAMASK_NumThreadsInt, & fftw_planner_flag, & fftw_timelimit use homogenization, only: & materialpoint_sizeResults, & materialpoint_results implicit none ! field in real an fourier space real(pReal), dimension(:,:,:,:,:), pointer :: P_real, deltaF_real ! field in real space (pointer) complex(pReal), dimension(:,:,:,:,:), pointer :: P_fourier,deltaF_fourier ! field in fourier space (pointer) !-------------------------------------------------------------------------------------------------- ! loop variables, convergence etc. real(pReal) :: time = 0.0_pReal, time0 = 0.0_pReal, timeinc = 1.0_pReal, timeinc_old = 0.0_pReal ! elapsed time, begin of interval, time interval real(pReal) :: guessmode, err_div, err_stress, err_stress_tol real(pReal), dimension(3,3), parameter :: ones = 1.0_pReal, zeroes = 0.0_pReal complex(pReal), dimension(3) :: temp3_Complex complex(pReal), dimension(3,3) :: temp33_Complex real(pReal), dimension(3,3) :: temp33_Real integer(pInt) :: i, j, k, l, m, n, p, errorID integer(pInt) :: N_Loadcases, loadcase = 0_pInt, inc, iter, ielem, CPFEM_mode=1_pInt, & ierr, totalIncsCounter = 0_pInt,& notConvergedCounter = 0_pInt, convergedCounter = 0_pInt logical :: errmatinv real(pReal) :: defgradDet character(len=6) :: loadcase_string !-------------------------------------------------------------------------------------------------- !variables controlling debugging logical :: debugGeneral, debugDivergence, debugRestart, debugFFTW !-------------------------------------------------------------------------------------------------- !variables for additional output due to general debugging real(pReal) :: defgradDetMax, defgradDetMin, maxCorrectionSym, maxCorrectionSkew !-------------------------------------------------------------------------------------------------- ! variables for additional output of divergence calculations type(C_PTR) :: divergence, plan_divergence real(pReal), dimension(:,:,:,:), pointer :: divergence_real complex(pReal), dimension(:,:,:,:), pointer :: divergence_fourier real(pReal), dimension(:,:,:,:), allocatable :: divergence_post real(pReal) :: pstress_av_L2, err_div_RMS, err_real_div_RMS, err_post_div_RMS,& err_div_max, err_real_div_max !-------------------------------------------------------------------------------------------------- ! variables for debugging fft using a scalar field type(C_PTR) :: scalarField_realC, scalarField_fourierC,& plan_scalarField_forth, plan_scalarField_back complex(pReal), dimension(:,:,:), pointer :: scalarField_real complex(pReal), dimension(:,:,:), pointer :: scalarField_fourier integer(pInt) :: row, column !################################################################################################## ! reading of information from load case file and geometry file !################################################################################################## subroutine init #ifdef PETSC integer :: ierr_psc call PetscInitialize(PETSC_NULL_CHARACTER, ierr_psc) #endif call DAMASK_interface_init write(6,'(a)') '' write(6,'(a)') ' <<<+- DAMASK_spectral init -+>>>' write(6,'(a)') ' $Id$' #include "compilation_info.f90" write(6,'(a)') '' !-------------------------------------------------------------------------------------------------- ! debugging parameters debugGeneral = iand(debug_level(debug_spectral),debug_levelBasic) /= 0 debugDivergence = iand(debug_level(debug_spectral),debug_spectralDivergence) /= 0 debugRestart = iand(debug_level(debug_spectral),debug_spectralRestart) /= 0 debugFFTW = iand(debug_level(debug_spectral),debug_spectralFFTW) /= 0 !################################################################################################## ! initialization !################################################################################################## call c_f_pointer(tensorField, P_real, [ res(1)+2_pInt,res(2),res(3),3,3]) ! place a pointer for a real representation on tensorField call c_f_pointer(tensorField, deltaF_real, [ res(1)+2_pInt,res(2),res(3),3,3]) ! place a pointer for a real representation on tensorField call c_f_pointer(tensorField, P_fourier, [ res1_red, res(2),res(3),3,3]) ! place a pointer for a complex representation on tensorField call c_f_pointer(tensorField, deltaF_fourier, [ res1_red, res(2),res(3),3,3]) ! place a pointer for a complex representation on tensorField !-------------------------------------------------------------------------------------------------- ! creating plans plan_stress = fftw_plan_many_dft_r2c(3,[ res(3),res(2) ,res(1)],9,& ! dimensions , length in each dimension in reversed order P_real,[ res(3),res(2) ,res(1)+2_pInt],& ! input data , physical length in each dimension in reversed order 1, res(3)*res(2)*(res(1)+2_pInt),& ! striding , product of physical lenght in the 3 dimensions P_fourier,[ res(3),res(2) ,res1_red],& 1, res(3)*res(2)* res1_red,fftw_planner_flag) plan_correction =fftw_plan_many_dft_c2r(3,[ res(3),res(2) ,res(1)],9,& deltaF_fourier,[ res(3),res(2) ,res1_red],& 1, res(3)*res(2)* res1_red,& deltaF_real,[ res(3),res(2) ,res(1)+2_pInt],& 1, res(3)*res(2)*(res(1)+2_pInt),fftw_planner_flag) !-------------------------------------------------------------------------------------------------- ! in case of no restart get reference material stiffness and init fields to no deformation if (restartInc == 1_pInt) then ielem = 0_pInt do k = 1_pInt, res(3); do j = 1_pInt, res(2); do i = 1_pInt, res(1) ielem = ielem + 1_pInt F(i,j,k,1:3,1:3) = math_I3 F_lastInc(i,j,k,1:3,1:3) = math_I3 coordinates(i,j,k,1:3) = geomdim/real(res,pReal)*real([i,j,k],pReal) - geomdim/real(2_pInt*res,pReal) call CPFEM_general(2_pInt,coordinates(i,j,k,1:3),math_I3,math_I3,temperature(i,j,k),& 0.0_pReal,ielem,1_pInt,sigma,dsde,P_real(i,j,k,1:3,1:3),dPdF) C = C + dPdF enddo; enddo; enddo C = C * wgt C_ref = C call IO_write_jobBinaryFile(777,'C_ref',size(C_ref)) write (777,rec=1) C_ref close(777) !-------------------------------------------------------------------------------------------------- ! restore deformation gradient and stiffness from saved state elseif (restartInc > 1_pInt) then ! using old values from file if (debugRestart) write(6,'(a,i6,a)') 'Reading values of increment ',& restartInc - 1_pInt,' from file' call IO_read_jobBinaryFile(777,'convergedSpectralDefgrad',& trim(getSolverJobName()),size(F)) read (777,rec=1) F close (777) F_lastInc = F F_aim = 0.0_pReal do k = 1_pInt, res(3); do j = 1_pInt, res(2); do i = 1_pInt, res(1) F_aim = F_aim + F(i,j,k,1:3,1:3) ! calculating old average deformation enddo; enddo; enddo F_aim = F_aim * wgt F_aim_lastInc = F_aim coordinates = 0.0 ! change it later!!! call IO_read_jobBinaryFile(777,'C_ref',trim(getSolverJobName()),size(C_ref)) read (777,rec=1) C_ref close (777) call IO_read_jobBinaryFile(777,'C',trim(getSolverJobName()),size(C)) read (777,rec=1) C close (777) CPFEM_mode = 2_pInt endif end subroutine init subroutine solution(guessmode, F_aim,F_aimLastInc, BC_stress, mask_stress) !-------------------------------------------------------------------------------------------------- ! update local deformation gradient and coordinates deltaF_aim = math_rotate_backward33(deltaF_aim,bc(loadcase)%rotation) do k = 1_pInt, res(3); do j = 1_pInt, res(2); do i = 1_pInt, res(1) temp33_Real = F(i,j,k,1:3,1:3) F(i,j,k,1:3,1:3) = F(i,j,k,1:3,1:3) & ! decide if guessing along former trajectory or apply homogeneous addon + guessmode * (F(i,j,k,1:3,1:3) - F_lastInc(i,j,k,1:3,1:3))& ! guessing... *timeinc/timeinc_old & + (1.0_pReal-guessmode) * deltaF_aim ! if not guessing, use prescribed average deformation where applicable F_lastInc(i,j,k,1:3,1:3) = temp33_Real enddo; enddo; enddo call deformed_fft(res,geomdim,math_rotate_backward33(F_aim,bc(loadcase)%rotation),& ! calculate current coordinates 1.0_pReal,F_lastInc,coordinates) guessmode = 1.0_pReal ! keep guessing along former trajectory during same loadcase iter = 0_pInt err_div = huge(err_div_tol) ! go into loop !################################################################################################## ! convergence loop (looping over iterations) !################################################################################################## do while((iter < itmax .and. (err_div > err_div_tol .or. err_stress > err_stress_tol))& .or. iter < itmin) iter = iter + 1_pInt !-------------------------------------------------------------------------------------------------- ! report begin of new iteration write(6,'(a)') '' write(6,'(a)') '==================================================================' write(6,'(6(a,i6.6))') 'Loadcase ',loadcase,' Inc. ',inc,'/',bc(loadcase)%incs,& ' @ Iter. ',itmin,' < ',iter,' < ',itmax write(6,'(a,/,3(3(f12.7,1x)/))',advance='no') 'deformation gradient aim =',& math_transpose33(F_aim) write(6,'(a)') '' write(6,'(a)') '... update stress field P(F) .....................................' if (restartWrite) write(6,'(a)') 'writing restart info for last increment' F_aim_lab_lastIter = math_rotate_backward33(F_aim,bc(loadcase)%rotation) !-------------------------------------------------------------------------------------------------- ! evaluate constitutive response ielem = 0_pInt do k = 1_pInt, res(3); do j = 1_pInt, res(2); do i = 1_pInt, res(1) ielem = ielem + 1_pInt call CPFEM_general(3_pInt,& ! collect cycle coordinates(i,j,k,1:3), F_lastInc(i,j,k,1:3,1:3),F(i,j,k,1:3,1:3), & temperature(i,j,k),timeinc,ielem,1_pInt,sigma,dsde,& P_real(i,j,k,1:3,1:3),dPdF) enddo; enddo; enddo P_real = 0.0_pReal ! needed because of the padding for FFTW C = 0.0_pReal ielem = 0_pInt call debug_reset() do k = 1_pInt, res(3); do j = 1_pInt, res(2); do i = 1_pInt, res(1) ielem = ielem + 1_pInt call CPFEM_general(CPFEM_mode,& ! first element in first iteration retains CPFEM_mode 1, coordinates(i,j,k,1:3),F_lastInc(i,j,k,1:3,1:3), F(i,j,k,1:3,1:3), & ! others get 2 (saves winding forward effort) temperature(i,j,k),timeinc,ielem,1_pInt,sigma,dsde, & P_real(i,j,k,1:3,1:3),dPdF) CPFEM_mode = 2_pInt C = C + dPdF enddo; enddo; enddo call debug_info() ! for test of regridding ! if( bc(loadcase)%restartFrequency > 0_pInt .and. & ! mod(inc-1,bc(loadcase)%restartFrequency) == 0_pInt .and. & ! restartInc/=inc) call quit(-1*(restartInc+1)) ! trigger exit to regrid !-------------------------------------------------------------------------------------------------- ! copy one component of the stress field to to a single FT and check for mismatch if (debugFFTW) then row = (mod(totalIncsCounter+iter-2_pInt,9_pInt))/3_pInt + 1_pInt ! go through the elements of the tensors, controlled by totalIncsCounter and iter, starting at 1 column = (mod(totalIncsCounter+iter-2_pInt,3_pInt)) + 1_pInt scalarField_real(1:res(1),1:res(2),1:res(3)) =& ! store the selected component cmplx(P_real(1:res(1),1:res(2),1:res(3),row,column),0.0_pReal,pReal) endif !-------------------------------------------------------------------------------------------------- ! call function to calculate divergence from math (for post processing) to check results if (debugDivergence) & call divergence_fft(res,virt_dim,3_pInt,& P_real(1:res(1),1:res(2),1:res(3),1:3,1:3),divergence_post) ! padding !-------------------------------------------------------------------------------------------------- ! doing the FT because it simplifies calculation of average stress in real space also call fftw_execute_dft_r2c(plan_stress,P_real,P_fourier) P_av_lab = real(P_fourier(1,1,1,1:3,1:3),pReal)*wgt P_av = math_rotate_forward33(P_av_lab,bc(loadcase)%rotation) write (6,'(a,/,3(3(f12.7,1x)/))',advance='no') 'Piola-Kirchhoff stress / MPa =',& math_transpose33(P_av)/1.e6_pReal !-------------------------------------------------------------------------------------------------- ! comparing 1 and 3x3 FT results if (debugFFTW) then call fftw_execute_dft(plan_scalarField_forth,scalarField_real,scalarField_fourier) write(6,'(a,i1,1x,i1)') 'checking FT results of compontent ', row, column write(6,'(a,2(es11.4,1x))') 'max FT relative error = ',& maxval( real((scalarField_fourier(1:res1_red,1:res(2),1:res(3))-& P_fourier(1:res1_red,1:res(2),1:res(3),row,column))/& scalarField_fourier(1:res1_red,1:res(2),1:res(3)))), & maxval(aimag((scalarField_fourier(1:res1_red,1:res(2),1:res(3))-& P_fourier(1:res1_red,1:res(2),1:res(3),row,column))/& scalarField_fourier(1:res1_red,1:res(2),1:res(3)))) endif !-------------------------------------------------------------------------------------------------- ! removing highest frequencies P_fourier ( res1_red,1:res(2) , 1:res(3) ,1:3,1:3)& = cmplx(0.0_pReal,0.0_pReal,pReal) P_fourier (1:res1_red, res(2)/2_pInt+1_pInt,1:res(3) ,1:3,1:3)& = cmplx(0.0_pReal,0.0_pReal,pReal) if(res(3)>1_pInt) & P_fourier (1:res1_red,1:res(2), res(3)/2_pInt+1_pInt,1:3,1:3)& = cmplx(0.0_pReal,0.0_pReal,pReal) !-------------------------------------------------------------------------------------------------- ! stress BC handling if(size_reduced > 0_pInt) then ! calculate stress BC if applied err_stress = maxval(abs(mask_stress * (P_av - bc(loadcase)%stress))) ! maximum deviaton (tensor norm not applicable) err_stress_tol = min(maxval(abs(P_av)) * err_stress_tolrel,err_stress_tolabs) ! don't use any tensor norm for the relative criterion because the comparison should be coherent write(6,'(a)') '' write(6,'(a)') '... correcting deformation gradient to fulfill BCs ...............' write(6,'(a,f6.2,a,es11.4,a)') 'error stress = ', err_stress/err_stress_tol, & ' (',err_stress,' Pa)' F_aim = F_aim - math_mul3333xx33(S_lastInc, ((P_av - bc(loadcase)%stress))) ! residual on given stress components write(6,'(a,1x,es11.4)')'determinant of new deformation = ',math_det33(F_aim) else err_stress_tol = +huge(1.0_pReal) endif F_aim_lab = math_rotate_backward33(F_aim,bc(loadcase)%rotation) ! boundary conditions from load frame into lab (Fourier) frame !-------------------------------------------------------------------------------------------------- ! actual spectral method write(6,'(a)') '' write(6,'(a)') '... calculating equilibrium with spectral method .................' !-------------------------------------------------------------------------------------------------- ! calculating RMS divergence criterion in Fourier space pstress_av_L2 = sqrt(maxval(math_eigenvalues33(math_mul33x33(P_av_lab,& ! L_2 norm of average stress (http://mathworld.wolfram.com/SpectralNorm.html) math_transpose33(P_av_lab))))) err_div_RMS = 0.0_pReal do k = 1_pInt, res(3); do j = 1_pInt, res(2) do i = 2_pInt, res1_red -1_pInt ! Has somewhere a conj. complex counterpart. Therefore count it twice. err_div_RMS = err_div_RMS & + 2.0_pReal*(sum (real(math_mul33x3_complex(P_fourier(i,j,k,1:3,1:3),& ! (sqrt(real(a)**2 + aimag(a)**2))**2 = real(a)**2 + aimag(a)**2. do not take square root and square again xi(1:3,i,j,k))*TWOPIIMG)**2.0_pReal)& ! --> sum squared L_2 norm of vector +sum(aimag(math_mul33x3_complex(P_fourier(i,j,k,1:3,1:3),& xi(1:3,i,j,k))*TWOPIIMG)**2.0_pReal)) enddo err_div_RMS = err_div_RMS & ! Those two layers (DC and Nyquist) do not have a conjugate complex counterpart + sum( real(math_mul33x3_complex(P_fourier(1 ,j,k,1:3,1:3),& xi(1:3,1 ,j,k))*TWOPIIMG)**2.0_pReal)& + sum(aimag(math_mul33x3_complex(P_fourier(1 ,j,k,1:3,1:3),& xi(1:3,1 ,j,k))*TWOPIIMG)**2.0_pReal)& + sum( real(math_mul33x3_complex(P_fourier(res1_red,j,k,1:3,1:3),& xi(1:3,res1_red,j,k))*TWOPIIMG)**2.0_pReal)& + sum(aimag(math_mul33x3_complex(P_fourier(res1_red,j,k,1:3,1:3),& xi(1:3,res1_red,j,k))*TWOPIIMG)**2.0_pReal) enddo; enddo err_div_RMS = sqrt(err_div_RMS)*wgt ! RMS in real space calculated with Parsevals theorem from Fourier space if (err_div_RMS/pstress_av_L2 > err_div & .and. err_stress < err_stress_tol & .and. iter >= itmin ) then write(6,'(a)') 'Increasing divergence, stopping iterations' iter = itmax endif err_div = err_div_RMS/pstress_av_L2 ! criterion to stop iterations !-------------------------------------------------------------------------------------------------- ! calculate additional divergence criteria and report if (debugDivergence) then ! calculate divergence again err_div_max = 0.0_pReal do k = 1_pInt, res(3); do j = 1_pInt, res(2); do i = 1_pInt, res1_red temp3_Complex = math_mul33x3_complex(P_fourier(i,j,k,1:3,1:3)*wgt,& ! weighting P_fourier xi(1:3,i,j,k))*TWOPIIMG err_div_max = max(err_div_max,sum(abs(temp3_Complex)**2.0_pReal)) divergence_fourier(i,j,k,1:3) = temp3_Complex ! need divergence NOT squared enddo; enddo; enddo call fftw_execute_dft_c2r(plan_divergence,divergence_fourier,divergence_real) ! already weighted err_real_div_RMS = 0.0_pReal err_post_div_RMS = 0.0_pReal err_real_div_max = 0.0_pReal do k = 1_pInt, res(3); do j = 1_pInt, res(2); do i = 1_pInt, res(1) err_real_div_RMS = err_real_div_RMS + sum(divergence_real(i,j,k,1:3)**2.0_pReal) ! avg of squared L_2 norm of div(stress) in real space err_post_div_RMS = err_post_div_RMS + sum(divergence_post(i,j,k,1:3)**2.0_pReal) ! avg of squared L_2 norm of div(stress) in real space err_real_div_max = max(err_real_div_max,sum(divergence_real(i,j,k,1:3)**2.0_pReal)) ! max of squared L_2 norm of div(stress) in real space enddo; enddo; enddo err_real_div_RMS = sqrt(wgt*err_real_div_RMS) ! RMS in real space err_post_div_RMS = sqrt(wgt*err_post_div_RMS) ! RMS in real space err_real_div_max = sqrt( err_real_div_max) ! max in real space err_div_max = sqrt( err_div_max) ! max in Fourier space write(6,'(a,es11.4)') 'error divergence FT RMS = ',err_div_RMS write(6,'(a,es11.4)') 'error divergence Real RMS = ',err_real_div_RMS write(6,'(a,es11.4)') 'error divergence post RMS = ',err_post_div_RMS write(6,'(a,es11.4)') 'error divergence FT max = ',err_div_max write(6,'(a,es11.4)') 'error divergence Real max = ',err_real_div_max endif write(6,'(a,f6.2,a,es11.4,a)') 'error divergence = ', err_div/err_div_tol,& ' (',err_div,' N/m³)' !-------------------------------------------------------------------------------------------------- ! to the actual spectral method calculation (mechanical equilibrium) if(memory_efficient) then ! memory saving version, on-the-fly calculation of gamma_hat do k = 1_pInt, res(3); do j = 1_pInt, res(2) ;do i = 1_pInt, res1_red if(any([i,j,k] /= 1_pInt)) then ! singular point at xi=(0.0,0.0,0.0) i.e. i=j=k=1 forall(l = 1_pInt:3_pInt, m = 1_pInt:3_pInt) & xiDyad(l,m) = xi(l, i,j,k)*xi(m, i,j,k) forall(l = 1_pInt:3_pInt, m = 1_pInt:3_pInt) & temp33_Real(l,m) = sum(C_ref(l,m,1:3,1:3)*xiDyad) temp33_Real = math_inv33(temp33_Real) forall(l=1_pInt:3_pInt, m=1_pInt:3_pInt, n=1_pInt:3_pInt, p=1_pInt:3_pInt)& gamma_hat(1,1,1, l,m,n,p) = temp33_Real(l,n)*xiDyad(m,p) forall(l = 1_pInt:3_pInt, m = 1_pInt:3_pInt) & temp33_Complex(l,m) = sum(gamma_hat(1,1,1, l,m, 1:3,1:3) *& P_fourier(i,j,k,1:3,1:3)) deltaF_fourier(i,j,k,1:3,1:3) = temp33_Complex endif enddo; enddo; enddo else ! use precalculated gamma-operator do k = 1_pInt, res(3); do j = 1_pInt, res(2); do i = 1_pInt,res1_red forall( m = 1_pInt:3_pInt, n = 1_pInt:3_pInt) & temp33_Complex(m,n) = sum(gamma_hat(i,j,k, m,n, 1:3,1:3) *& P_fourier(i,j,k,1:3,1:3)) deltaF_fourier(i,j,k, 1:3,1:3) = temp33_Complex enddo; enddo; enddo endif deltaF_fourier(1,1,1,1:3,1:3) = cmplx((F_aim_lab_lastIter - F_aim_lab) & ! assign (negative) average deformation gradient change to zero frequency (real part) * real(Npoints,pReal),0.0_pReal,pReal) ! singular point at xi=(0.0,0.0,0.0) i.e. i=j=k=1 !-------------------------------------------------------------------------------------------------- ! comparing 1 and 3x3 inverse FT results if (debugFFTW) then do k = 1_pInt, res(3); do j = 1_pInt, res(2); do i = 1_pInt, res1_red scalarField_fourier(i,j,k) = deltaF_fourier(i,j,k,row,column) enddo; enddo; enddo do i = 0_pInt, res(1)/2_pInt-2_pInt ! unpack fft data for conj complex symmetric part m = 1_pInt do k = 1_pInt, res(3) n = 1_pInt do j = 1_pInt, res(2) scalarField_fourier(res(1)-i,j,k) = conjg(scalarField_fourier(2+i,n,m)) if(n == 1_pInt) n = res(2) + 1_pInt n = n-1_pInt enddo if(m == 1_pInt) m = res(3) + 1_pInt m = m -1_pInt enddo; enddo endif !-------------------------------------------------------------------------------------------------- ! doing the inverse FT call fftw_execute_dft_c2r(plan_correction,deltaF_fourier,deltaF_real) ! back transform of fluct deformation gradient !-------------------------------------------------------------------------------------------------- ! comparing 1 and 3x3 inverse FT results if (debugFFTW) then write(6,'(a,i1,1x,i1)') 'checking iFT results of compontent ', row, column call fftw_execute_dft(plan_scalarField_back,scalarField_fourier,scalarField_real) write(6,'(a,es11.4)') 'max iFT relative error = ',& maxval((real(scalarField_real(1:res(1),1:res(2),1:res(3)))-& deltaF_real(1:res(1),1:res(2),1:res(3),row,column))/& real(scalarField_real(1:res(1),1:res(2),1:res(3)))) endif !-------------------------------------------------------------------------------------------------- ! calculate some additional output if(debugGeneral) then maxCorrectionSkew = 0.0_pReal maxCorrectionSym = 0.0_pReal temp33_Real = 0.0_pReal do k = 1_pInt, res(3); do j = 1_pInt, res(2); do i = 1_pInt, res(1) maxCorrectionSym = max(maxCorrectionSym,& maxval(math_symmetric33(deltaF_real(i,j,k,1:3,1:3)))) maxCorrectionSkew = max(maxCorrectionSkew,& maxval(math_skew33(deltaF_real(i,j,k,1:3,1:3)))) temp33_Real = temp33_Real + deltaF_real(i,j,k,1:3,1:3) enddo; enddo; enddo write(6,'(a,1x,es11.4)') 'max symmetric correction of deformation =',& maxCorrectionSym*wgt write(6,'(a,1x,es11.4)') 'max skew correction of deformation =',& maxCorrectionSkew*wgt write(6,'(a,1x,es11.4)') 'max sym/skew of avg correction = ',& maxval(math_symmetric33(temp33_real))/& maxval(math_skew33(temp33_real)) endif !-------------------------------------------------------------------------------------------------- ! updated deformation gradient do k = 1_pInt, res(3); do j = 1_pInt, res(2); do i = 1_pInt, res(1) F(i,j,k,1:3,1:3) = F(i,j,k,1:3,1:3) - deltaF_real(i,j,k,1:3,1:3)*wgt ! F(x)^(n+1) = F(x)^(n) + correction; *wgt: correcting for missing normalization enddo; enddo; enddo !-------------------------------------------------------------------------------------------------- ! calculate bounds of det(F) and report if(debugGeneral) then defgradDetMax = -huge(1.0_pReal) defgradDetMin = +huge(1.0_pReal) do k = 1_pInt, res(3); do j = 1_pInt, res(2); do i = 1_pInt, res(1) defgradDet = math_det33(F(i,j,k,1:3,1:3)) defgradDetMax = max(defgradDetMax,defgradDet) defgradDetMin = min(defgradDetMin,defgradDet) enddo; enddo; enddo write(6,'(a,1x,es11.4)') 'max determinant of deformation =', defgradDetMax write(6,'(a,1x,es11.4)') 'min determinant of deformation =', defgradDetMin endif enddo ! end looping when convergency is achieved CPFEM_mode = 1_pInt ! winding forward C = C * wgt write(6,'(a)') '' write(6,'(a)') '==================================================================' if(err_div > err_div_tol .or. err_stress > err_stress_tol) then write(6,'(A,I5.5,A)') 'increment ', totalIncsCounter, ' NOT converged' notConvergedCounter = notConvergedCounter + 1_pInt else convergedCounter = convergedCounter + 1_pInt write(6,'(A,I5.5,A)') 'increment ', totalIncsCounter, ' converged' endif if (mod(inc,bc(loadcase)%outputFrequency) == 0_pInt) then ! at output frequency write(6,'(a)') '' write(6,'(a)') '... writing results to file ......................................' write(538) materialpoint_results(1_pInt:materialpoint_sizeResults,1,1_pInt:Npoints) ! write result to file flush(538) endif if( bc(loadcase)%restartFrequency > 0_pInt .and. & mod(inc,bc(loadcase)%restartFrequency) == 0_pInt) then ! at frequency of writing restart information set restart parameter for FEsolving (first call to CPFEM_general will write ToDo: true?) restartInc=totalIncsCounter restartWrite = .true. write(6,'(a)') 'writing converged results for restart' call IO_write_jobBinaryFile(777,'convergedSpectralDefgrad',size(F)) ! writing deformation gradient field to file write (777,rec=1) F close (777) call IO_write_jobBinaryFile(777,'C',size(C)) write (777,rec=1) C close(777) endif endif ! end calculation/forwarding enddo ! end looping over incs in current loadcase