From 5a85af050adb9c3971d5cff37ee76d65674f0e5a Mon Sep 17 00:00:00 2001 From: Salvatore Filippone Date: Thu, 31 Jul 2008 14:15:58 +0000 Subject: [PATCH] psblas2-typext: Beginning of experiments with type, extends. --- Make.inc.in | 3 +- base/modules/psb_spmat_type.f90 | 51 +++++++++++++++------------------ configure.ac | 16 +++++++++++ 3 files changed, 41 insertions(+), 29 deletions(-) diff --git a/Make.inc.in b/Make.inc.in index 62851d15..540af565 100755 --- a/Make.inc.in +++ b/Make.inc.in @@ -6,11 +6,12 @@ # These lines are quite portable. .mod=@MODEXT@ .fh=.fh -.SUFFIXES: .f90 $(.mod) .F90 .F +.SUFFIXES: .f90 $(.mod) .F90 .F .f03 .F03 # The following ones are the variables used by the PSBLAS make scripts. +F03=@F90@ F90=@F90@ FC=@FC@ CC=@CC@ diff --git a/base/modules/psb_spmat_type.f90 b/base/modules/psb_spmat_type.f90 index 54ea1c59..641d79ca 100644 --- a/base/modules/psb_spmat_type.f90 +++ b/base/modules/psb_spmat_type.f90 @@ -120,45 +120,30 @@ module psb_spmat_type character(len=5) :: psb_fidef_='CSR' - - type psb_sspmat_type - integer :: m, k - character(len=5) :: fida + type psb_base_spmat_type + integer :: m, k + character(len=5) :: fida character(len=11) :: descra - integer :: infoa(psb_ifasize_) - real(psb_spk_), allocatable :: aspk(:) + integer :: infoa(psb_ifasize_) integer, allocatable :: ia1(:), ia2(:) integer, allocatable :: pl(:), pr(:) + end type psb_base_spmat_type + + + type, extends(psb_base_spmat_type) :: psb_sspmat_type + real(psb_spk_), allocatable :: aspk(:) end type psb_sspmat_type - type psb_cspmat_type - integer :: m, k - character(len=5) :: fida - character(len=11) :: descra - integer :: infoa(psb_ifasize_) + type, extends(psb_base_spmat_type) :: psb_cspmat_type complex(psb_spk_), allocatable :: aspk(:) - integer, allocatable :: ia1(:), ia2(:) - integer, allocatable :: pl(:), pr(:) end type psb_cspmat_type - type psb_dspmat_type - integer :: m, k - character(len=5) :: fida - character(len=11) :: descra - integer :: infoa(psb_ifasize_) + type, extends(psb_base_spmat_type) :: psb_dspmat_type real(psb_dpk_), allocatable :: aspk(:) - integer, allocatable :: ia1(:), ia2(:) - integer, allocatable :: pl(:), pr(:) end type psb_dspmat_type - type psb_zspmat_type - integer :: m, k - character(len=5) :: fida - character(len=11) :: descra - integer :: infoa(psb_ifasize_) + type, extends(psb_base_spmat_type) :: psb_zspmat_type complex(psb_dpk_), allocatable :: aspk(:) - integer, allocatable :: ia1(:), ia2(:) - integer, allocatable :: pl(:), pr(:) end type psb_zspmat_type interface psb_nullify_sp @@ -464,7 +449,17 @@ contains psb_get_zsp_nnz_row = 0 end if end function psb_get_zsp_nnz_row - +!!$ +!!$ subroutine psb_nullify_base_sp(mat) +!!$ implicit none +!!$ class(psb_base_spmat_type), intent(inout) :: mat +!!$ mat%infoa(:)=0 +!!$ mat%m=0 +!!$ mat%k=0 +!!$ mat%fida='' +!!$ mat%descra='' +!!$ +!!$ end subroutine psb_nullify_base_sp subroutine psb_nullify_ssp(mat) implicit none diff --git a/configure.ac b/configure.ac index 5baec23b..0315c441 100755 --- a/configure.ac +++ b/configure.ac @@ -671,6 +671,10 @@ $(.mod).o: $(F90) $(F90COPT) $(FINCLUDES) -c $< %$(.mod): %.f90 $(F90) $(F90COPT) $(FINCLUDES) -c $< +%.o: %.f03 + $(F90) $(F90COPT) $(FINCLUDES) -c $< +%$(.mod): %.f03 + $(F90) $(F90COPT) $(FINCLUDES) -c $< %.o: %.F $(FC) $(FCOPT) $(FINCLUDES) $(FDEFINES) -c $< %$(.mod): %.F @@ -679,6 +683,10 @@ $(.mod).o: $(F90) $(F90COPT) $(FINCLUDES) $(FDEFINES) -c $< %$(.mod): %.F90 $(F90) $(F90COPT) $(FINCLUDES) $(FDEFINES) -c $<' +%.o: %.F03 + $(F03) $(F90COPT) $(FINCLUDES) $(FDEFINES) -c $< +%$(.mod): %.F03 + $(F03) $(F90COPT) $(FINCLUDES) $(FDEFINES) -c $<' else @@ -701,6 +709,10 @@ $(.mod).o: $(F90) $(F90COPT) $(FINCLUDES) -c $< .f90.o: $(F90) $(F90COPT) $(FINCLUDES) -c $< +.f03$(.mod): + $(F03) $(F90COPT) $(FINCLUDES) -c $< +.f03.o: + $(F03) $(F90COPT) $(FINCLUDES) -c $< .F.o: $(FC) $(FCOPT) $(FINCLUDES) $(FDEFINES) -c $< .F$(.mod): @@ -709,6 +721,10 @@ $(.mod).o: $(F90) $(F90COPT) $(FINCLUDES) $(FDEFINES) -c $< .F90$(.mod): $(F90) $(F90COPT) $(FINCLUDES) $(FDEFINES) -c $<' +.F03.o: + $(F03) $(F90COPT) $(FINCLUDES) $(FDEFINES) -c $< +.F03$(.mod): + $(F03) $(F90COPT) $(FINCLUDES) $(FDEFINES) -c $<' fi AC_SUBST(PSBLASRULES)