[cig-commits] r12461 - seismo/2D/SPECFEM2D/trunk
nlegoff at geodynamics.org
nlegoff at geodynamics.org
Mon Jul 21 14:14:12 PDT 2008
Author: nlegoff
Date: 2008-07-21 14:14:12 -0700 (Mon, 21 Jul 2008)
New Revision: 12461
Modified:
seismo/2D/SPECFEM2D/trunk/constants.h
seismo/2D/SPECFEM2D/trunk/specfem2D.F90
Log:
applied Cuthill-McKee on whole domain. ACTUALLY_IMPLEMENT_PERM_WHOLE should be set to true (in constants.h) in order to activate it.
Modified: seismo/2D/SPECFEM2D/trunk/constants.h
===================================================================
--- seismo/2D/SPECFEM2D/trunk/constants.h 2008-07-21 20:07:38 UTC (rev 12460)
+++ seismo/2D/SPECFEM2D/trunk/constants.h 2008-07-21 21:14:12 UTC (rev 12461)
@@ -31,8 +31,9 @@
integer, parameter :: LIMIT_MULTI_CUTHILL = 50
! implement Cuthill-McKee or replace with identity permutation
- logical, parameter :: ACTUALLY_IMPLEMENT_PERM_OUT = .true.
- logical, parameter :: ACTUALLY_IMPLEMENT_PERM_INN = .true.
+ logical, parameter :: ACTUALLY_IMPLEMENT_PERM_OUT = .false.
+ logical, parameter :: ACTUALLY_IMPLEMENT_PERM_INN = .false.
+ logical, parameter :: ACTUALLY_IMPLEMENT_PERM_WHOLE = .true.
! add MPI barriers and suppress seismograms if we generate traces of the run for analysis with "ParaVer"
logical, parameter :: GENERATE_PARAVER_TRACES = .false.
Modified: seismo/2D/SPECFEM2D/trunk/specfem2D.F90
===================================================================
--- seismo/2D/SPECFEM2D/trunk/specfem2D.F90 2008-07-21 20:07:38 UTC (rev 12460)
+++ seismo/2D/SPECFEM2D/trunk/specfem2D.F90 2008-07-21 21:14:12 UTC (rev 12461)
@@ -1470,6 +1470,28 @@
enddo
endif
+ if(ACTUALLY_IMPLEMENT_PERM_WHOLE) then
+
+ allocate(check_perm(nspec))
+ call get_perm(ibool,perm(1:nspec),LIMIT_MULTI_CUTHILL,nspec,npoin)
+! check that the permutation obtained is bijective
+ check_perm(:) = -1
+ do ispec = 1,nspec
+ check_perm(perm(ispec)) = ispec
+ enddo
+ if(minval(check_perm) /= 1) stop 'minval check_perm is incorrect for whole'
+ if(maxval(check_perm) /= nspec) stop 'maxval check_perm is incorrect for whole'
+ deallocate(check_perm)
+! add the right offset
+ perm(1:nspec) = perm(1:nspec) + 0
+
+ else
+! use identity permutation if flag is off
+ do ispec = 1,nspec
+ perm(ispec) = ispec
+ enddo
+ endif
+
endif
enddo ! end of further reduction of cache misses inner/outer in two passes
More information about the cig-commits
mailing list