[cig-commits] r8550 - seismo/2D/SPECFEM2D/trunk

walter at geodynamics.org walter at geodynamics.org
Fri Dec 7 15:55:39 PST 2007


Author: walter
Date: 2007-12-07 15:55:38 -0800 (Fri, 07 Dec 2007)
New Revision: 8550

Modified:
   seismo/2D/SPECFEM2D/trunk/specfem2D.F90
Log:
added sorting of ibool to reduce cache misses


Modified: seismo/2D/SPECFEM2D/trunk/specfem2D.F90
===================================================================
--- seismo/2D/SPECFEM2D/trunk/specfem2D.F90	2007-07-03 23:39:30 UTC (rev 8549)
+++ seismo/2D/SPECFEM2D/trunk/specfem2D.F90	2007-12-07 23:55:38 UTC (rev 8550)
@@ -276,8 +276,6 @@
   integer, dimension(:,:), allocatable  :: acoustic_surface
   integer, dimension(:,:), allocatable  :: acoustic_edges
 
-
-
   integer  :: ixmin, ixmax, izmin, izmax
 
   integer  :: ie, num_interface
@@ -285,6 +283,10 @@
   integer  :: nrecloc, irecloc
   integer, dimension(:), allocatable :: recloc, which_proc_receiver
 
+! mask to sort ibool
+  integer, dimension(:), allocatable :: mask_ibool
+  integer, dimension(:,:,:), allocatable :: copy_ibool_ori
+  integer :: inumber
 
 !***********************************************************************
 !
@@ -658,7 +660,33 @@
     call createnum_slow(knods,ibool,npoin,nspec,ngnod)
   endif
 
-!---- compute shape functions and their derivatives for regular !interpolated display grid
+! create a new indirect addressing array instead, to reduce cache misses
+! in memory access in the solver
+  allocate(copy_ibool_ori(NGLLX,NGLLZ,nspec))
+  allocate(mask_ibool(npoin))
+  mask_ibool(:) = -1
+  copy_ibool_ori(:,:,:) = ibool(:,:,:)
+
+  inumber = 0
+  do ispec=1,nspec
+    do j=1,NGLLZ
+      do i=1,NGLLX
+        if(mask_ibool(copy_ibool_ori(i,j,ispec)) == -1) then
+! create a new point
+          inumber = inumber + 1
+          ibool(i,j,ispec) = inumber
+          mask_ibool(copy_ibool_ori(i,j,ispec)) = inumber
+        else
+! use an existing point created previously
+          ibool(i,j,ispec) = mask_ibool(copy_ibool_ori(i,j,ispec))
+        endif
+      enddo
+    enddo
+  enddo
+  deallocate(copy_ibool_ori)
+  deallocate(mask_ibool)
+
+!---- compute shape functions and their derivatives for regular interpolated display grid
   do j = 1,pointsdisp
     do i = 1,pointsdisp
       xirec  = 2.d0*dble(i-1)/dble(pointsdisp-1) - 1.d0



More information about the cig-commits mailing list