planc
Parallel Lowrank Approximation with Non-negativity Constraints
dimtrees.h
Go to the documentation of this file.
1 /* Copyright Koby Hayashi 2018 */
2 
3 #ifndef DIMTREE_DIMTREES_H_
4 #define DIMTREE_DIMTREES_H_
5 
6 #include "dimtree/ddttensor.h"
7 
16 typedef enum { RowMajor, ColMajor } Output_Layout;
17 
22 void partial_MTTKRP(Output_Layout OL, long int s, direction D, tensor *T,
23  double *A, long int r, double *C, long int num_threads);
24 void partial_MTTKRP_with_KRP(Output_Layout OL, long int s, direction D,
25  ktensor *Y, tensor *T, double *C,
26  long int num_threads);
27 void partial_MTTKRP_with_KRP_output_FM(direction D, ktensor *Y, tensor *T,
28  long int num_threads);
29 void partial_MTTKRP_with_KRP_output_T(long int s, direction D,
30  ktensor *input_ktensor,
31  tensor *input_tensor,
32  tensor *output_tensor,
33  long int num_threads);
34 
35 void multi_TTV(Output_Layout OL, long int s, direction D, tensor *T, double *A,
36  long int r, double *C, long int num_threads);
37 void multi_TTV_with_KRP(Output_Layout OL, long int s, direction D, tensor *T,
38  ktensor *Y, double *C, long int num_threads);
39 void multi_TTV_with_KRP_output_FM(direction D, tensor *input_tensor,
40  ktensor *input_ktensor, long int num_threads);
41 void multi_TTV_with_KRP_output_T(long int s, direction D, tensor *input_tensor,
42  ktensor *input_ktensor, tensor *output_tensor,
43  long int num_threads);
44 #endif // DIMTREE_DIMTREES_H_
void multi_TTV_with_KRP_output_FM(direction D, tensor *input_tensor, ktensor *input_ktensor, long int num_threads)
multi_TTV_with_KRP_output_FM(); Wrapper function for performing a multi_TTV with a KRP and outputin...
Definition: dimtrees.hpp:464
void partial_MTTKRP(Output_Layout OL, long int s, direction D, tensor *T, double *A, long int r, double *C, long int num_threads)
Partial MTTKRP functions.
Definition: dimtrees.hpp:31
void multi_TTV_with_KRP_output_T(long int s, direction D, tensor *input_tensor, ktensor *input_ktensor, tensor *output_tensor, long int num_threads)
Definition: dimtrees.hpp:481
void multi_TTV_with_KRP(Output_Layout OL, long int s, direction D, tensor *T, ktensor *Y, double *C, long int num_threads)
multi_TTV_with_KRP(); KRP wrapper for the general multi_TTV function.
Definition: dimtrees.hpp:418
void partial_MTTKRP_with_KRP(Output_Layout OL, long int s, direction D, ktensor *Y, tensor *T, double *C, long int num_threads)
Wrapper function for partial_MTTKRP() Forms a desired KRP, manages all the memory for the KRP Passes ...
Definition: dimtrees.hpp:176
void partial_MTTKRP_with_KRP_output_FM(direction D, ktensor *Y, tensor *T, long int num_threads)
partial_MTTKRP_with_KRP_output_FM(); Wrapper function for computing a partial_MTTKRP with KRP and o...
Definition: dimtrees.hpp:250
void partial_MTTKRP_with_KRP_output_T(long int s, direction D, ktensor *input_ktensor, tensor *input_tensor, tensor *output_tensor, long int num_threads)
partial_MTTKRP_with_KRP_output_T() Wrapper function for computing partial_MTTKRP with KRP and outpu...
Definition: dimtrees.hpp:281
void multi_TTV(Output_Layout OL, long int s, direction D, tensor *T, double *A, long int r, double *C, long int num_threads)
multi_TTV(); performs a multi_ttv with between subtensors of T->data and the column of the matrix ...
Definition: dimtrees.hpp:315
Output_Layout
This .h file contains all the functions need to construct dimension trees for fast_CP_ALS algorithms...
Definition: dimtrees.h:16