xref: /petsc/src/mat/utils/hpl/hplcreate.c (revision 9baa572492f2d5c9a9469100eb79348cae8ec756)
1*3444e80cSBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h"  I*/
2*3444e80cSBarry Smith #include "hpl.h"
3*3444e80cSBarry Smith 
4*3444e80cSBarry Smith /*@
5*3444e80cSBarry Smith   MatSetHPL - fills a `MATSEQDENSE` matrix using the HPL 2.3 random matrix generation routine
6*3444e80cSBarry Smith 
7*3444e80cSBarry Smith   Collective
8*3444e80cSBarry Smith 
9*3444e80cSBarry Smith   Input Parameters:
10*3444e80cSBarry Smith + A     - the matrix
11*3444e80cSBarry Smith - iseed - the random number seed
12*3444e80cSBarry Smith 
13*3444e80cSBarry Smith   Level: intermediate
14*3444e80cSBarry Smith 
15*3444e80cSBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreate()`
16*3444e80cSBarry Smith @*/
MatSetHPL(Mat A,int iseed)17*3444e80cSBarry Smith PetscErrorCode MatSetHPL(Mat A, int iseed)
18*3444e80cSBarry Smith {
19*3444e80cSBarry Smith   PetscBool    isDense;
20*3444e80cSBarry Smith   PetscInt     M, N, LDA;
21*3444e80cSBarry Smith   PetscBLASInt bM, bN, bLDA;
22*3444e80cSBarry Smith   PetscScalar *values;
23*3444e80cSBarry Smith 
24*3444e80cSBarry Smith   PetscFunctionBegin;
25*3444e80cSBarry Smith   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
26*3444e80cSBarry Smith   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQDENSE, &isDense));
27*3444e80cSBarry Smith   PetscCheck(isDense, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Only supports sequential dense matrices");
28*3444e80cSBarry Smith   PetscCall(MatGetSize(A, &M, &N));
29*3444e80cSBarry Smith   PetscCall(PetscBLASIntCast(M, &bM));
30*3444e80cSBarry Smith   PetscCall(PetscBLASIntCast(N, &bN));
31*3444e80cSBarry Smith   PetscCall(MatDenseGetLDA(A, &LDA));
32*3444e80cSBarry Smith   PetscCall(PetscBLASIntCast(LDA, &bLDA));
33*3444e80cSBarry Smith   PetscCall(MatDenseGetArrayWrite(A, &values));
34*3444e80cSBarry Smith   PetscStackCallExternalVoid("HPL_dmatgen", HPL_dmatgen(bM, bN, values, bLDA, iseed));
35*3444e80cSBarry Smith   PetscCall(MatDenseRestoreArrayWrite(A, &values));
36*3444e80cSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
37*3444e80cSBarry Smith }
38*3444e80cSBarry Smith 
39*3444e80cSBarry Smith #include <petsc/private/bmimpl.h>
40*3444e80cSBarry Smith 
41*3444e80cSBarry Smith typedef struct {
42*3444e80cSBarry Smith   Mat A, F;
43*3444e80cSBarry Smith } PetscBench_HPL;
44*3444e80cSBarry Smith 
PetscBenchSetUp_HPL(PetscBench bm)45*3444e80cSBarry Smith static PetscErrorCode PetscBenchSetUp_HPL(PetscBench bm)
46*3444e80cSBarry Smith {
47*3444e80cSBarry Smith   PetscBench_HPL *hp = (PetscBench_HPL *)bm->data;
48*3444e80cSBarry Smith   PetscMPIInt     rank;
49*3444e80cSBarry Smith 
50*3444e80cSBarry Smith   PetscFunctionBegin;
51*3444e80cSBarry Smith   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)bm), &rank));
52*3444e80cSBarry Smith   if (rank > 0) PetscFunctionReturn(PETSC_SUCCESS);
53*3444e80cSBarry Smith   if (bm->size == PETSC_DECIDE) bm->size = 2500;
54*3444e80cSBarry Smith   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, bm->size, bm->size, NULL, &hp->A));
55*3444e80cSBarry Smith   PetscCall(MatSetHPL(hp->A, 0));
56*3444e80cSBarry Smith   PetscCall(MatGetFactor(hp->A, MATSOLVERPETSC, MAT_FACTOR_LU, &hp->F));
57*3444e80cSBarry Smith   PetscCall(MatLUFactorSymbolic(hp->F, hp->A, NULL, NULL, 0));
58*3444e80cSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
59*3444e80cSBarry Smith }
60*3444e80cSBarry Smith 
PetscBenchRun_HPL(PetscBench bm)61*3444e80cSBarry Smith static PetscErrorCode PetscBenchRun_HPL(PetscBench bm)
62*3444e80cSBarry Smith {
63*3444e80cSBarry Smith   PetscBench_HPL *hp = (PetscBench_HPL *)bm->data;
64*3444e80cSBarry Smith   PetscMPIInt     rank;
65*3444e80cSBarry Smith 
66*3444e80cSBarry Smith   PetscFunctionBegin;
67*3444e80cSBarry Smith   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)bm), &rank));
68*3444e80cSBarry Smith   if (rank > 0) PetscFunctionReturn(PETSC_SUCCESS);
69*3444e80cSBarry Smith   PetscCall(MatLUFactorNumeric(hp->F, hp->A, 0));
70*3444e80cSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
71*3444e80cSBarry Smith }
72*3444e80cSBarry Smith 
PetscBenchView_HPL(PetscBench bm,PetscViewer viewer)73*3444e80cSBarry Smith static PetscErrorCode PetscBenchView_HPL(PetscBench bm, PetscViewer viewer)
74*3444e80cSBarry Smith {
75*3444e80cSBarry Smith   PetscEventPerfInfo *info;
76*3444e80cSBarry Smith   PetscInt            numThreads;
77*3444e80cSBarry Smith 
78*3444e80cSBarry Smith   PetscFunctionBegin;
79*3444e80cSBarry Smith   PetscCall(PetscLogHandlerGetEventPerfInfo(bm->lhdlr, 0, MAT_LUFactor, &info)); // because symbolic is trivial it does not log numeric!
80*3444e80cSBarry Smith   PetscCall(PetscBLASGetNumThreads(&numThreads));
81*3444e80cSBarry Smith   PetscCall(PetscViewerASCIIPrintf(viewer, "HPL Benchmark, number of threads %d, matrix size %d, flop rate %g (megaflops)\n", (int)numThreads, (int)bm->size, (double)(info->flops / (1.e7 * info->time))));
82*3444e80cSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
83*3444e80cSBarry Smith }
84*3444e80cSBarry Smith 
PetscBenchReset_HPL(PetscBench bm)85*3444e80cSBarry Smith static PetscErrorCode PetscBenchReset_HPL(PetscBench bm)
86*3444e80cSBarry Smith {
87*3444e80cSBarry Smith   PetscBench_HPL *hp = (PetscBench_HPL *)bm->data;
88*3444e80cSBarry Smith   PetscMPIInt     rank;
89*3444e80cSBarry Smith 
90*3444e80cSBarry Smith   PetscFunctionBegin;
91*3444e80cSBarry Smith   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)bm), &rank));
92*3444e80cSBarry Smith   if (rank > 0) PetscFunctionReturn(PETSC_SUCCESS);
93*3444e80cSBarry Smith   PetscCall(MatDestroy(&hp->F));
94*3444e80cSBarry Smith   PetscCall(MatDestroy(&hp->A));
95*3444e80cSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
96*3444e80cSBarry Smith }
97*3444e80cSBarry Smith 
PetscBenchDestroy_HPL(PetscBench bm)98*3444e80cSBarry Smith static PetscErrorCode PetscBenchDestroy_HPL(PetscBench bm)
99*3444e80cSBarry Smith {
100*3444e80cSBarry Smith   PetscFunctionBegin;
101*3444e80cSBarry Smith   PetscCall(PetscFree(bm->data));
102*3444e80cSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
103*3444e80cSBarry Smith }
104*3444e80cSBarry Smith 
PetscBenchCreate_HPL(PetscBench bm)105*3444e80cSBarry Smith PETSC_INTERN PetscErrorCode PetscBenchCreate_HPL(PetscBench bm)
106*3444e80cSBarry Smith {
107*3444e80cSBarry Smith   PetscBench_HPL *hp;
108*3444e80cSBarry Smith 
109*3444e80cSBarry Smith   PetscFunctionBegin;
110*3444e80cSBarry Smith   PetscCall(PetscNew(&hp));
111*3444e80cSBarry Smith   bm->data         = hp;
112*3444e80cSBarry Smith   bm->ops->setup   = PetscBenchSetUp_HPL;
113*3444e80cSBarry Smith   bm->ops->run     = PetscBenchRun_HPL;
114*3444e80cSBarry Smith   bm->ops->view    = PetscBenchView_HPL;
115*3444e80cSBarry Smith   bm->ops->reset   = PetscBenchReset_HPL;
116*3444e80cSBarry Smith   bm->ops->destroy = PetscBenchDestroy_HPL;
117*3444e80cSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
118*3444e80cSBarry Smith }
119