]> Gitweb @ Texas Instruments - Open Source Git Repositories - git.TI.com/gitweb - dense-linear-algebra-libraries/linalg.git/commitdiff
Pass NUM_DSP_CORES to ticblas build since it uses blis.h.
authorJianzhong Xu <xuj@ti.com>
Fri, 8 Apr 2016 19:13:42 +0000 (19:13 +0000)
committerJianzhong Xu <xuj@ti.com>
Fri, 8 Apr 2016 19:13:42 +0000 (19:13 +0000)
src/ti/linalg/Makefile
src/ti/linalg/ticblas/src/Makefile
src/ti/linalg/ticblas/src/ticblas.c

index b3af5bdee2a682408c565ad8a17232a54d21b33a..870866b1d008659a7f9a04d162ccfbd3a17b44e7 100644 (file)
@@ -53,7 +53,7 @@ DSPlibs:
        cd $(LINALG_CBLAS_DIR); make arch=C66 alllib; \
        cd ../$(LINALG_BLIS_DIR); ./configure -p install/$(BLIS_CFG) c66x; \
        make -j8 MEM_MODEL=$(MEM_MODEL) TARGET=$(TARGET) LIBOS=$(LIBOS) NUM_DSP_CORES=$(NUM_DSP_CORES); make install; \
-       cd ../$(LINALG_TICBLAS_DIR)/src; make MEM_MODEL=$(MEM_MODEL) TARGET=$(TARGET) LIBOS=$(LIBOS); cd ../lib; \
+       cd ../$(LINALG_TICBLAS_DIR)/src; make MEM_MODEL=$(MEM_MODEL) TARGET=$(TARGET) LIBOS=$(LIBOS) NUM_DSP_CORES=$(NUM_DSP_CORES); cd ../lib; \
        echo "combining BLIS, CBLAS, and TICBLAS libraries into one: libcblas.ae66"; \
        mkdir -p objs; cd objs; rm -f *; ar x ../../../blis/install/$(BLIS_CFG)/lib/libblis.ae66; mmv 'cblas*.o' 'blis_cblas#1.o'; \
        ar -x ../../../cblas/lib/C66/libcblas.ae66; ar -x ../libticblas.ae66; chmod +rw *;cd ../../..; \
index 6ac4ad54ff7fc1f295bdc2636e9741635692897a..1e45ca8697cf2c54a7201bf02bede69fcc98ebb0 100644 (file)
@@ -32,7 +32,7 @@ INCDIR += -I$(BLIS_INC)
 
 INCS = -I. -I$(strip $(subst ;, -I,$(subst $(space),$(space),$(INCDIR))))
 
-CL6X_FLAGS  = $(INCS) --openmp --use_g2 -D$(TARGET) -D$(LIBOS) -D$(BLIS_MEM_MODEL)
+CL6X_FLAGS  = $(INCS) --openmp --use_g2 -D$(TARGET) -D$(LIBOS) -D$(BLIS_MEM_MODEL) -DBLIS_MAX_NUM_THREADS=${NUM_DSP_CORES}
 
 DSP_LIB_DIR = ../lib
 DSP_LIB = $(DSP_LIB_DIR)/libticblas.ae66
index 6beabc2fb41c0c3389be4f72a2d7f4d9eac584c1..3af4e4a6716493a81f566ada50c18218e574879d 100644 (file)
@@ -86,10 +86,10 @@ void * blasGetMemHandle()
 void tiCblasGetSizes(size_t *smem_size_vfast,  size_t *smem_size_fast, \r
                      size_t *smem_size_medium, size_t *smem_size_slow)\r
 {\r
-    *smem_size_vfast  = BLAS_MEM_SIZE_VFAST;  // very fast scratch memory\r
-    *smem_size_fast   = BLAS_MEM_SIZE_FAST;   // fast scratch memory\r
-    *smem_size_medium = BLAS_MEM_SIZE_MEDIUM; // medium speed scratch memory\r
-    *smem_size_slow   = BLAS_MEM_SIZE_SLOW;   // slow scratch memory\r
+    *smem_size_vfast  = BLAS_MEM_SIZE_VFAST;  /* very fast scratch memory     */\r
+    *smem_size_fast   = BLAS_MEM_SIZE_FAST;   /* fast scratch memory          */\r
+    *smem_size_medium = BLAS_MEM_SIZE_MEDIUM; /* medium speed scratch memory  */ \r
+    *smem_size_slow   = BLAS_MEM_SIZE_SLOW;   /* slow scratch memory          */\r
 /*\r
     printf("BLIS_MK_POOL_SIZE_L1 is %d.\n", BLIS_MK_POOL_SIZE_L1);\r
     printf("BLIS_KN_POOL_SIZE_L1 is %d.\n", BLIS_KN_POOL_SIZE_L1);\r