xref: /petsc/src/mat/impls/h2opus/cuda/math2opus.cu (revision 4dfa11a44d5adf2389f1d3acbc8f3c1116dc6c3a)
1f17f7ee2SSatish Balay #include <h2opusconf.h>
2f17f7ee2SSatish Balay /* skip compilation of this .cu file if H2OPUS is CPU only while PETSc has GPU support */
3f17f7ee2SSatish Balay #if !defined(__CUDACC__) || defined(H2OPUS_USE_GPU)
4f17f7ee2SSatish Balay #include <h2opus.h>
5f17f7ee2SSatish Balay #if defined(H2OPUS_USE_MPI)
6f17f7ee2SSatish Balay #include <h2opus/distributed/distributed_h2opus_handle.h>
7f17f7ee2SSatish Balay #include <h2opus/distributed/distributed_geometric_construction.h>
8f17f7ee2SSatish Balay #include <h2opus/distributed/distributed_hgemv.h>
9f17f7ee2SSatish Balay #include <h2opus/distributed/distributed_horthog.h>
10f17f7ee2SSatish Balay #include <h2opus/distributed/distributed_hcompress.h>
11f17f7ee2SSatish Balay #endif
12f17f7ee2SSatish Balay #include <h2opus/util/boxentrygen.h>
13f17f7ee2SSatish Balay #include <petsc/private/matimpl.h>
14f17f7ee2SSatish Balay #include <petsc/private/vecimpl.h>
15f17f7ee2SSatish Balay #include <petsc/private/deviceimpl.h>
16f17f7ee2SSatish Balay #include <petscsf.h>
17f17f7ee2SSatish Balay 
18f17f7ee2SSatish Balay /* math2opusutils */
19f17f7ee2SSatish Balay PETSC_INTERN PetscErrorCode PetscSFGetVectorSF(PetscSF, PetscInt, PetscInt, PetscInt, PetscSF *);
20f17f7ee2SSatish Balay PETSC_INTERN PetscErrorCode MatDenseGetH2OpusVectorSF(Mat, PetscSF, PetscSF *);
21f17f7ee2SSatish Balay PETSC_INTERN PetscErrorCode VecSign(Vec, Vec);
22f17f7ee2SSatish Balay PETSC_INTERN PetscErrorCode VecSetDelta(Vec, PetscInt);
23f17f7ee2SSatish Balay PETSC_INTERN PetscErrorCode MatApproximateNorm_Private(Mat, NormType, PetscInt, PetscReal *);
24f17f7ee2SSatish Balay 
25f17f7ee2SSatish Balay #define MatH2OpusGetThrustPointer(v) thrust::raw_pointer_cast((v).data())
26f17f7ee2SSatish Balay 
27f17f7ee2SSatish Balay /* Use GPU only if H2OPUS is configured for GPU */
28f17f7ee2SSatish Balay #if defined(PETSC_HAVE_CUDA) && defined(H2OPUS_USE_GPU)
29f17f7ee2SSatish Balay #define PETSC_H2OPUS_USE_GPU
30f17f7ee2SSatish Balay #endif
31f17f7ee2SSatish Balay #if defined(PETSC_H2OPUS_USE_GPU)
32f17f7ee2SSatish Balay #define MatH2OpusUpdateIfNeeded(A, B) MatBindToCPU(A, (PetscBool)((A)->boundtocpu || (B)))
33f17f7ee2SSatish Balay #else
34f17f7ee2SSatish Balay #define MatH2OpusUpdateIfNeeded(A, B) 0
35f17f7ee2SSatish Balay #endif
36f17f7ee2SSatish Balay 
37f17f7ee2SSatish Balay // TODO H2OPUS:
38f17f7ee2SSatish Balay // DistributedHMatrix
39f17f7ee2SSatish Balay //   unsymmetric ?
40f17f7ee2SSatish Balay //   transpose for distributed_hgemv?
41f17f7ee2SSatish Balay //   clearData()
42f17f7ee2SSatish Balay // Unify interface for sequential and parallel?
43f17f7ee2SSatish Balay // Reuse geometric construction (almost possible, only the unsymmetric case is explicitly handled)
44f17f7ee2SSatish Balay //
459371c9d4SSatish Balay template <class T>
469371c9d4SSatish Balay class PetscPointCloud : public H2OpusDataSet<T> {
47f17f7ee2SSatish Balay private:
48f17f7ee2SSatish Balay   int            dimension;
49f17f7ee2SSatish Balay   size_t         num_points;
50f17f7ee2SSatish Balay   std::vector<T> pts;
51f17f7ee2SSatish Balay 
52f17f7ee2SSatish Balay public:
539371c9d4SSatish Balay   PetscPointCloud(int dim, size_t num_pts, const T coords[]) {
54f17f7ee2SSatish Balay     dim              = dim > 0 ? dim : 1;
55f17f7ee2SSatish Balay     this->dimension  = dim;
56f17f7ee2SSatish Balay     this->num_points = num_pts;
57f17f7ee2SSatish Balay 
58f17f7ee2SSatish Balay     pts.resize(num_pts * dim);
59f17f7ee2SSatish Balay     if (coords) {
60036e4bdcSSatish Balay       for (size_t n = 0; n < num_pts; n++)
619371c9d4SSatish Balay         for (int i = 0; i < dim; i++) pts[n * dim + i] = coords[n * dim + i];
62f17f7ee2SSatish Balay     } else {
63036e4bdcSSatish Balay       PetscReal h = 1.0; //num_pts > 1 ? 1./(num_pts - 1) : 0.0;
64036e4bdcSSatish Balay       for (size_t n = 0; n < num_pts; n++) {
65036e4bdcSSatish Balay         pts[n * dim] = n * h;
669371c9d4SSatish Balay         for (int i = 1; i < dim; i++) pts[n * dim + i] = 0.0;
67036e4bdcSSatish Balay       }
68f17f7ee2SSatish Balay     }
69f17f7ee2SSatish Balay   }
70f17f7ee2SSatish Balay 
719371c9d4SSatish Balay   PetscPointCloud(const PetscPointCloud<T> &other) {
72f17f7ee2SSatish Balay     size_t N         = other.dimension * other.num_points;
73f17f7ee2SSatish Balay     this->dimension  = other.dimension;
74f17f7ee2SSatish Balay     this->num_points = other.num_points;
75f17f7ee2SSatish Balay     this->pts.resize(N);
769371c9d4SSatish Balay     for (size_t i = 0; i < N; i++) this->pts[i] = other.pts[i];
77f17f7ee2SSatish Balay   }
78f17f7ee2SSatish Balay 
799371c9d4SSatish Balay   int getDimension() const { return dimension; }
80f17f7ee2SSatish Balay 
819371c9d4SSatish Balay   size_t getDataSetSize() const { return num_points; }
82f17f7ee2SSatish Balay 
839371c9d4SSatish Balay   T getDataPoint(size_t idx, int dim) const {
84f17f7ee2SSatish Balay     assert(dim < dimension && idx < num_points);
85f17f7ee2SSatish Balay     return pts[idx * dimension + dim];
86f17f7ee2SSatish Balay   }
87f17f7ee2SSatish Balay 
889371c9d4SSatish Balay   void Print(std::ostream &out = std::cout) {
89f17f7ee2SSatish Balay     out << "Dimension: " << dimension << std::endl;
90f17f7ee2SSatish Balay     out << "NumPoints: " << num_points << std::endl;
91f17f7ee2SSatish Balay     for (size_t n = 0; n < num_points; n++) {
929371c9d4SSatish Balay       for (int d = 0; d < dimension; d++) out << pts[n * dimension + d] << " ";
93f17f7ee2SSatish Balay       out << std::endl;
94f17f7ee2SSatish Balay     }
95f17f7ee2SSatish Balay   }
96f17f7ee2SSatish Balay };
97f17f7ee2SSatish Balay 
989371c9d4SSatish Balay template <class T>
999371c9d4SSatish Balay class PetscFunctionGenerator {
100f17f7ee2SSatish Balay private:
101f17f7ee2SSatish Balay   MatH2OpusKernel k;
102f17f7ee2SSatish Balay   int             dim;
103f17f7ee2SSatish Balay   void           *ctx;
104f17f7ee2SSatish Balay 
105f17f7ee2SSatish Balay public:
1069371c9d4SSatish Balay   PetscFunctionGenerator(MatH2OpusKernel k, int dim, void *ctx) {
1079371c9d4SSatish Balay     this->k   = k;
1089371c9d4SSatish Balay     this->dim = dim;
1099371c9d4SSatish Balay     this->ctx = ctx;
110f17f7ee2SSatish Balay   }
1119371c9d4SSatish Balay   PetscFunctionGenerator(PetscFunctionGenerator &other) {
1129371c9d4SSatish Balay     this->k   = other.k;
1139371c9d4SSatish Balay     this->dim = other.dim;
1149371c9d4SSatish Balay     this->ctx = other.ctx;
1159371c9d4SSatish Balay   }
1169371c9d4SSatish Balay   T operator()(PetscReal *pt1, PetscReal *pt2) { return (T)((*this->k)(this->dim, pt1, pt2, this->ctx)); }
117f17f7ee2SSatish Balay };
118f17f7ee2SSatish Balay 
119f17f7ee2SSatish Balay #include <../src/mat/impls/h2opus/math2opussampler.hpp>
120f17f7ee2SSatish Balay 
121f17f7ee2SSatish Balay /* just to not clutter the code */
122f17f7ee2SSatish Balay #if !defined(H2OPUS_USE_GPU)
123f17f7ee2SSatish Balay typedef HMatrix HMatrix_GPU;
124f17f7ee2SSatish Balay #if defined(H2OPUS_USE_MPI)
125f17f7ee2SSatish Balay typedef DistributedHMatrix DistributedHMatrix_GPU;
126f17f7ee2SSatish Balay #endif
127f17f7ee2SSatish Balay #endif
128f17f7ee2SSatish Balay 
129f17f7ee2SSatish Balay typedef struct {
130f17f7ee2SSatish Balay #if defined(H2OPUS_USE_MPI)
131f17f7ee2SSatish Balay   distributedH2OpusHandle_t handle;
132f17f7ee2SSatish Balay #else
133f17f7ee2SSatish Balay   h2opusHandle_t handle;
134f17f7ee2SSatish Balay #endif
135f17f7ee2SSatish Balay   /* Sequential and parallel matrices are two different classes at the moment */
136f17f7ee2SSatish Balay   HMatrix *hmatrix;
137f17f7ee2SSatish Balay #if defined(H2OPUS_USE_MPI)
138f17f7ee2SSatish Balay   DistributedHMatrix *dist_hmatrix;
139f17f7ee2SSatish Balay #else
140f17f7ee2SSatish Balay   HMatrix       *dist_hmatrix;     /* just to not clutter the code */
141f17f7ee2SSatish Balay #endif
142f17f7ee2SSatish Balay   /* May use permutations */
143f17f7ee2SSatish Balay   PetscSF                           sf;
144f17f7ee2SSatish Balay   PetscLayout                       h2opus_rmap, h2opus_cmap;
145f17f7ee2SSatish Balay   IS                                h2opus_indexmap;
146f17f7ee2SSatish Balay   thrust::host_vector<PetscScalar> *xx, *yy;
147f17f7ee2SSatish Balay   PetscInt                          xxs, yys;
148f17f7ee2SSatish Balay   PetscBool                         multsetup;
149f17f7ee2SSatish Balay 
150f17f7ee2SSatish Balay   /* GPU */
151f17f7ee2SSatish Balay   HMatrix_GPU *hmatrix_gpu;
152f17f7ee2SSatish Balay #if defined(H2OPUS_USE_MPI)
153f17f7ee2SSatish Balay   DistributedHMatrix_GPU *dist_hmatrix_gpu;
154f17f7ee2SSatish Balay #else
155f17f7ee2SSatish Balay   HMatrix_GPU   *dist_hmatrix_gpu; /* just to not clutter the code */
156f17f7ee2SSatish Balay #endif
157f17f7ee2SSatish Balay #if defined(PETSC_H2OPUS_USE_GPU)
158f17f7ee2SSatish Balay   thrust::device_vector<PetscScalar> *xx_gpu, *yy_gpu;
159f17f7ee2SSatish Balay   PetscInt                            xxs_gpu, yys_gpu;
160f17f7ee2SSatish Balay #endif
161f17f7ee2SSatish Balay 
162f17f7ee2SSatish Balay   /* construction from matvecs */
163f17f7ee2SSatish Balay   PetscMatrixSampler *sampler;
164f17f7ee2SSatish Balay   PetscBool           nativemult;
165f17f7ee2SSatish Balay 
166f17f7ee2SSatish Balay   /* Admissibility */
167f17f7ee2SSatish Balay   PetscReal eta;
168f17f7ee2SSatish Balay   PetscInt  leafsize;
169f17f7ee2SSatish Balay 
170f17f7ee2SSatish Balay   /* for dof reordering */
171f17f7ee2SSatish Balay   PetscPointCloud<PetscReal> *ptcloud;
172f17f7ee2SSatish Balay 
173f17f7ee2SSatish Balay   /* kernel for generating matrix entries */
174f17f7ee2SSatish Balay   PetscFunctionGenerator<PetscScalar> *kernel;
175f17f7ee2SSatish Balay 
176f17f7ee2SSatish Balay   /* basis orthogonalized? */
177f17f7ee2SSatish Balay   PetscBool orthogonal;
178f17f7ee2SSatish Balay 
179f17f7ee2SSatish Balay   /* customization */
180f17f7ee2SSatish Balay   PetscInt  basisord;
181f17f7ee2SSatish Balay   PetscInt  max_rank;
182f17f7ee2SSatish Balay   PetscInt  bs;
183f17f7ee2SSatish Balay   PetscReal rtol;
184f17f7ee2SSatish Balay   PetscInt  norm_max_samples;
185f17f7ee2SSatish Balay   PetscBool check_construction;
186f17f7ee2SSatish Balay   PetscBool hara_verbose;
187036e4bdcSSatish Balay   PetscBool resize;
188f17f7ee2SSatish Balay 
189f17f7ee2SSatish Balay   /* keeps track of MatScale values */
190f17f7ee2SSatish Balay   PetscScalar s;
191f17f7ee2SSatish Balay } Mat_H2OPUS;
192f17f7ee2SSatish Balay 
1939371c9d4SSatish Balay static PetscErrorCode MatDestroy_H2OPUS(Mat A) {
194f17f7ee2SSatish Balay   Mat_H2OPUS *a = (Mat_H2OPUS *)A->data;
195f17f7ee2SSatish Balay 
196f17f7ee2SSatish Balay   PetscFunctionBegin;
197f17f7ee2SSatish Balay #if defined(H2OPUS_USE_MPI)
198f17f7ee2SSatish Balay   h2opusDestroyDistributedHandle(a->handle);
199f17f7ee2SSatish Balay #else
200f17f7ee2SSatish Balay   h2opusDestroyHandle(a->handle);
201f17f7ee2SSatish Balay #endif
202f17f7ee2SSatish Balay   delete a->dist_hmatrix;
203f17f7ee2SSatish Balay   delete a->hmatrix;
204f17f7ee2SSatish Balay   PetscCall(PetscSFDestroy(&a->sf));
205f17f7ee2SSatish Balay   PetscCall(PetscLayoutDestroy(&a->h2opus_rmap));
206f17f7ee2SSatish Balay   PetscCall(PetscLayoutDestroy(&a->h2opus_cmap));
207f17f7ee2SSatish Balay   PetscCall(ISDestroy(&a->h2opus_indexmap));
208f17f7ee2SSatish Balay   delete a->xx;
209f17f7ee2SSatish Balay   delete a->yy;
210f17f7ee2SSatish Balay   delete a->hmatrix_gpu;
211f17f7ee2SSatish Balay   delete a->dist_hmatrix_gpu;
212f17f7ee2SSatish Balay #if defined(PETSC_H2OPUS_USE_GPU)
213f17f7ee2SSatish Balay   delete a->xx_gpu;
214f17f7ee2SSatish Balay   delete a->yy_gpu;
215f17f7ee2SSatish Balay #endif
216f17f7ee2SSatish Balay   delete a->sampler;
217f17f7ee2SSatish Balay   delete a->ptcloud;
218f17f7ee2SSatish Balay   delete a->kernel;
219f17f7ee2SSatish Balay   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_h2opus_seqdense_C", NULL));
220f17f7ee2SSatish Balay   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_h2opus_seqdensecuda_C", NULL));
221f17f7ee2SSatish Balay   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_h2opus_mpidense_C", NULL));
222f17f7ee2SSatish Balay   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_h2opus_mpidensecuda_C", NULL));
223f17f7ee2SSatish Balay   PetscCall(PetscObjectChangeTypeName((PetscObject)A, NULL));
224f17f7ee2SSatish Balay   PetscCall(PetscFree(A->data));
225f17f7ee2SSatish Balay   PetscFunctionReturn(0);
226f17f7ee2SSatish Balay }
227f17f7ee2SSatish Balay 
2289371c9d4SSatish Balay PetscErrorCode MatH2OpusSetNativeMult(Mat A, PetscBool nm) {
229f17f7ee2SSatish Balay   Mat_H2OPUS *a = (Mat_H2OPUS *)A->data;
230f17f7ee2SSatish Balay   PetscBool   ish2opus;
231f17f7ee2SSatish Balay 
232f17f7ee2SSatish Balay   PetscFunctionBegin;
233f17f7ee2SSatish Balay   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
234f17f7ee2SSatish Balay   PetscValidLogicalCollectiveBool(A, nm, 2);
235f17f7ee2SSatish Balay   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATH2OPUS, &ish2opus));
236f17f7ee2SSatish Balay   if (ish2opus) {
237f17f7ee2SSatish Balay     if (a->h2opus_rmap) { /* need to swap layouts for vector creation */
238f17f7ee2SSatish Balay       if ((!a->nativemult && nm) || (a->nativemult && !nm)) {
239f17f7ee2SSatish Balay         PetscLayout t;
240f17f7ee2SSatish Balay         t              = A->rmap;
241f17f7ee2SSatish Balay         A->rmap        = a->h2opus_rmap;
242f17f7ee2SSatish Balay         a->h2opus_rmap = t;
243f17f7ee2SSatish Balay         t              = A->cmap;
244f17f7ee2SSatish Balay         A->cmap        = a->h2opus_cmap;
245f17f7ee2SSatish Balay         a->h2opus_cmap = t;
246f17f7ee2SSatish Balay       }
247f17f7ee2SSatish Balay     }
248f17f7ee2SSatish Balay     a->nativemult = nm;
249f17f7ee2SSatish Balay   }
250f17f7ee2SSatish Balay   PetscFunctionReturn(0);
251f17f7ee2SSatish Balay }
252f17f7ee2SSatish Balay 
2539371c9d4SSatish Balay PetscErrorCode MatH2OpusGetNativeMult(Mat A, PetscBool *nm) {
254f17f7ee2SSatish Balay   Mat_H2OPUS *a = (Mat_H2OPUS *)A->data;
255f17f7ee2SSatish Balay   PetscBool   ish2opus;
256f17f7ee2SSatish Balay 
257f17f7ee2SSatish Balay   PetscFunctionBegin;
258f17f7ee2SSatish Balay   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
259f17f7ee2SSatish Balay   PetscValidPointer(nm, 2);
260f17f7ee2SSatish Balay   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATH2OPUS, &ish2opus));
261f17f7ee2SSatish Balay   PetscCheck(ish2opus, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Not for type %s", ((PetscObject)A)->type_name);
262f17f7ee2SSatish Balay   *nm = a->nativemult;
263f17f7ee2SSatish Balay   PetscFunctionReturn(0);
264f17f7ee2SSatish Balay }
265f17f7ee2SSatish Balay 
2669371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode MatNorm_H2OPUS(Mat A, NormType normtype, PetscReal *n) {
267f17f7ee2SSatish Balay   PetscBool   ish2opus;
268f17f7ee2SSatish Balay   PetscInt    nmax = PETSC_DECIDE;
269f17f7ee2SSatish Balay   Mat_H2OPUS *a    = NULL;
270f17f7ee2SSatish Balay   PetscBool   mult = PETSC_FALSE;
271f17f7ee2SSatish Balay 
272f17f7ee2SSatish Balay   PetscFunctionBegin;
273f17f7ee2SSatish Balay   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATH2OPUS, &ish2opus));
274f17f7ee2SSatish Balay   if (ish2opus) { /* set userdefine number of samples and fastpath for mult (norms are order independent) */
275f17f7ee2SSatish Balay     a = (Mat_H2OPUS *)A->data;
276f17f7ee2SSatish Balay 
277f17f7ee2SSatish Balay     nmax = a->norm_max_samples;
278f17f7ee2SSatish Balay     mult = a->nativemult;
279f17f7ee2SSatish Balay     PetscCall(MatH2OpusSetNativeMult(A, PETSC_TRUE));
280f17f7ee2SSatish Balay   } else {
281f17f7ee2SSatish Balay     PetscCall(PetscOptionsGetInt(((PetscObject)A)->options, ((PetscObject)A)->prefix, "-mat_approximate_norm_samples", &nmax, NULL));
282f17f7ee2SSatish Balay   }
283f17f7ee2SSatish Balay   PetscCall(MatApproximateNorm_Private(A, normtype, nmax, n));
284f17f7ee2SSatish Balay   if (a) PetscCall(MatH2OpusSetNativeMult(A, mult));
285f17f7ee2SSatish Balay   PetscFunctionReturn(0);
286f17f7ee2SSatish Balay }
287f17f7ee2SSatish Balay 
2889371c9d4SSatish Balay static PetscErrorCode MatH2OpusResizeBuffers_Private(Mat A, PetscInt xN, PetscInt yN) {
289f17f7ee2SSatish Balay   Mat_H2OPUS *h2opus = (Mat_H2OPUS *)A->data;
290f17f7ee2SSatish Balay   PetscInt    n;
291f17f7ee2SSatish Balay   PetscBool   boundtocpu = PETSC_TRUE;
292f17f7ee2SSatish Balay 
293f17f7ee2SSatish Balay   PetscFunctionBegin;
294f17f7ee2SSatish Balay #if defined(PETSC_H2OPUS_USE_GPU)
295f17f7ee2SSatish Balay   boundtocpu = A->boundtocpu;
296f17f7ee2SSatish Balay #endif
297f17f7ee2SSatish Balay   PetscCall(PetscSFGetGraph(h2opus->sf, NULL, &n, NULL, NULL));
298f17f7ee2SSatish Balay   if (boundtocpu) {
2999371c9d4SSatish Balay     if (h2opus->xxs < xN) {
3009371c9d4SSatish Balay       h2opus->xx->resize(n * xN);
3019371c9d4SSatish Balay       h2opus->xxs = xN;
3029371c9d4SSatish Balay     }
3039371c9d4SSatish Balay     if (h2opus->yys < yN) {
3049371c9d4SSatish Balay       h2opus->yy->resize(n * yN);
3059371c9d4SSatish Balay       h2opus->yys = yN;
3069371c9d4SSatish Balay     }
307f17f7ee2SSatish Balay   }
308f17f7ee2SSatish Balay #if defined(PETSC_H2OPUS_USE_GPU)
309f17f7ee2SSatish Balay   if (!boundtocpu) {
3109371c9d4SSatish Balay     if (h2opus->xxs_gpu < xN) {
3119371c9d4SSatish Balay       h2opus->xx_gpu->resize(n * xN);
3129371c9d4SSatish Balay       h2opus->xxs_gpu = xN;
3139371c9d4SSatish Balay     }
3149371c9d4SSatish Balay     if (h2opus->yys_gpu < yN) {
3159371c9d4SSatish Balay       h2opus->yy_gpu->resize(n * yN);
3169371c9d4SSatish Balay       h2opus->yys_gpu = yN;
3179371c9d4SSatish Balay     }
318f17f7ee2SSatish Balay   }
319f17f7ee2SSatish Balay #endif
320f17f7ee2SSatish Balay   PetscFunctionReturn(0);
321f17f7ee2SSatish Balay }
322f17f7ee2SSatish Balay 
3239371c9d4SSatish Balay static PetscErrorCode MatMultNKernel_H2OPUS(Mat A, PetscBool transA, Mat B, Mat C) {
324f17f7ee2SSatish Balay   Mat_H2OPUS *h2opus = (Mat_H2OPUS *)A->data;
325f17f7ee2SSatish Balay #if defined(H2OPUS_USE_MPI)
326f17f7ee2SSatish Balay   h2opusHandle_t handle = h2opus->handle->handle;
327f17f7ee2SSatish Balay #else
328f17f7ee2SSatish Balay   h2opusHandle_t handle = h2opus->handle;
329f17f7ee2SSatish Balay #endif
330f17f7ee2SSatish Balay   PetscBool    boundtocpu = PETSC_TRUE;
331f17f7ee2SSatish Balay   PetscScalar *xx, *yy, *uxx, *uyy;
332f17f7ee2SSatish Balay   PetscInt     blda, clda;
333f17f7ee2SSatish Balay   PetscMPIInt  size;
334f17f7ee2SSatish Balay   PetscSF      bsf, csf;
335f17f7ee2SSatish Balay   PetscBool    usesf = (PetscBool)(h2opus->sf && !h2opus->nativemult);
336f17f7ee2SSatish Balay 
337f17f7ee2SSatish Balay   PetscFunctionBegin;
338f17f7ee2SSatish Balay   HLibProfile::clear();
339f17f7ee2SSatish Balay #if defined(PETSC_H2OPUS_USE_GPU)
340f17f7ee2SSatish Balay   boundtocpu = A->boundtocpu;
341f17f7ee2SSatish Balay #endif
342f17f7ee2SSatish Balay   PetscCall(MatDenseGetLDA(B, &blda));
343f17f7ee2SSatish Balay   PetscCall(MatDenseGetLDA(C, &clda));
344f17f7ee2SSatish Balay   if (usesf) {
345f17f7ee2SSatish Balay     PetscInt n;
346f17f7ee2SSatish Balay 
347f17f7ee2SSatish Balay     PetscCall(MatDenseGetH2OpusVectorSF(B, h2opus->sf, &bsf));
348f17f7ee2SSatish Balay     PetscCall(MatDenseGetH2OpusVectorSF(C, h2opus->sf, &csf));
349f17f7ee2SSatish Balay 
350f17f7ee2SSatish Balay     PetscCall(MatH2OpusResizeBuffers_Private(A, B->cmap->N, C->cmap->N));
351f17f7ee2SSatish Balay     PetscCall(PetscSFGetGraph(h2opus->sf, NULL, &n, NULL, NULL));
352f17f7ee2SSatish Balay     blda = n;
353f17f7ee2SSatish Balay     clda = n;
354f17f7ee2SSatish Balay   }
355f17f7ee2SSatish Balay   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
356f17f7ee2SSatish Balay   if (boundtocpu) {
357f17f7ee2SSatish Balay     PetscCall(MatDenseGetArrayRead(B, (const PetscScalar **)&xx));
358f17f7ee2SSatish Balay     PetscCall(MatDenseGetArrayWrite(C, &yy));
359f17f7ee2SSatish Balay     if (usesf) {
360f17f7ee2SSatish Balay       uxx = MatH2OpusGetThrustPointer(*h2opus->xx);
361f17f7ee2SSatish Balay       uyy = MatH2OpusGetThrustPointer(*h2opus->yy);
362f17f7ee2SSatish Balay       PetscCall(PetscSFBcastBegin(bsf, MPIU_SCALAR, xx, uxx, MPI_REPLACE));
363f17f7ee2SSatish Balay       PetscCall(PetscSFBcastEnd(bsf, MPIU_SCALAR, xx, uxx, MPI_REPLACE));
364f17f7ee2SSatish Balay     } else {
365f17f7ee2SSatish Balay       uxx = xx;
366f17f7ee2SSatish Balay       uyy = yy;
367f17f7ee2SSatish Balay     }
368f17f7ee2SSatish Balay     if (size > 1) {
369f17f7ee2SSatish Balay       PetscCheck(h2opus->dist_hmatrix, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing distributed CPU matrix");
370036e4bdcSSatish Balay       PetscCheck(!transA || A->symmetric, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "MatMultTranspose not yet coded in parallel");
371f17f7ee2SSatish Balay #if defined(H2OPUS_USE_MPI)
372f17f7ee2SSatish Balay       distributed_hgemv(/* transA ? H2Opus_Trans : H2Opus_NoTrans, */ h2opus->s, *h2opus->dist_hmatrix, uxx, blda, 0.0, uyy, clda, B->cmap->N, h2opus->handle);
373f17f7ee2SSatish Balay #endif
374f17f7ee2SSatish Balay     } else {
375f17f7ee2SSatish Balay       PetscCheck(h2opus->hmatrix, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing CPU matrix");
376f17f7ee2SSatish Balay       hgemv(transA ? H2Opus_Trans : H2Opus_NoTrans, h2opus->s, *h2opus->hmatrix, uxx, blda, 0.0, uyy, clda, B->cmap->N, handle);
377f17f7ee2SSatish Balay     }
378f17f7ee2SSatish Balay     PetscCall(MatDenseRestoreArrayRead(B, (const PetscScalar **)&xx));
379f17f7ee2SSatish Balay     if (usesf) {
380f17f7ee2SSatish Balay       PetscCall(PetscSFReduceBegin(csf, MPIU_SCALAR, uyy, yy, MPI_REPLACE));
381f17f7ee2SSatish Balay       PetscCall(PetscSFReduceEnd(csf, MPIU_SCALAR, uyy, yy, MPI_REPLACE));
382f17f7ee2SSatish Balay     }
383f17f7ee2SSatish Balay     PetscCall(MatDenseRestoreArrayWrite(C, &yy));
384f17f7ee2SSatish Balay #if defined(PETSC_H2OPUS_USE_GPU)
385f17f7ee2SSatish Balay   } else {
386f17f7ee2SSatish Balay     PetscBool ciscuda, biscuda;
387f17f7ee2SSatish Balay 
388f17f7ee2SSatish Balay     /* If not of type seqdensecuda, convert on the fly (i.e. allocate GPU memory) */
389f17f7ee2SSatish Balay     PetscCall(PetscObjectTypeCompareAny((PetscObject)B, &biscuda, MATSEQDENSECUDA, MATMPIDENSECUDA, ""));
39048a46eb9SPierre Jolivet     if (!biscuda) PetscCall(MatConvert(B, MATDENSECUDA, MAT_INPLACE_MATRIX, &B));
391f17f7ee2SSatish Balay     PetscCall(PetscObjectTypeCompareAny((PetscObject)C, &ciscuda, MATSEQDENSECUDA, MATMPIDENSECUDA, ""));
392f17f7ee2SSatish Balay     if (!ciscuda) {
393f17f7ee2SSatish Balay       C->assembled = PETSC_TRUE;
394f17f7ee2SSatish Balay       PetscCall(MatConvert(C, MATDENSECUDA, MAT_INPLACE_MATRIX, &C));
395f17f7ee2SSatish Balay     }
396f17f7ee2SSatish Balay     PetscCall(MatDenseCUDAGetArrayRead(B, (const PetscScalar **)&xx));
397f17f7ee2SSatish Balay     PetscCall(MatDenseCUDAGetArrayWrite(C, &yy));
398f17f7ee2SSatish Balay     if (usesf) {
399f17f7ee2SSatish Balay       uxx = MatH2OpusGetThrustPointer(*h2opus->xx_gpu);
400f17f7ee2SSatish Balay       uyy = MatH2OpusGetThrustPointer(*h2opus->yy_gpu);
401f17f7ee2SSatish Balay       PetscCall(PetscSFBcastBegin(bsf, MPIU_SCALAR, xx, uxx, MPI_REPLACE));
402f17f7ee2SSatish Balay       PetscCall(PetscSFBcastEnd(bsf, MPIU_SCALAR, xx, uxx, MPI_REPLACE));
403f17f7ee2SSatish Balay     } else {
404f17f7ee2SSatish Balay       uxx = xx;
405f17f7ee2SSatish Balay       uyy = yy;
406f17f7ee2SSatish Balay     }
407f17f7ee2SSatish Balay     PetscCall(PetscLogGpuTimeBegin());
408f17f7ee2SSatish Balay     if (size > 1) {
409f17f7ee2SSatish Balay       PetscCheck(h2opus->dist_hmatrix_gpu, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing distributed GPU matrix");
410036e4bdcSSatish Balay       PetscCheck(!transA || A->symmetric, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "MatMultTranspose not yet coded in parallel");
411f17f7ee2SSatish Balay #if defined(H2OPUS_USE_MPI)
412f17f7ee2SSatish Balay       distributed_hgemv(/* transA ? H2Opus_Trans : H2Opus_NoTrans, */ h2opus->s, *h2opus->dist_hmatrix_gpu, uxx, blda, 0.0, uyy, clda, B->cmap->N, h2opus->handle);
413f17f7ee2SSatish Balay #endif
414f17f7ee2SSatish Balay     } else {
415f17f7ee2SSatish Balay       PetscCheck(h2opus->hmatrix_gpu, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing GPU matrix");
416f17f7ee2SSatish Balay       hgemv(transA ? H2Opus_Trans : H2Opus_NoTrans, h2opus->s, *h2opus->hmatrix_gpu, uxx, blda, 0.0, uyy, clda, B->cmap->N, handle);
417f17f7ee2SSatish Balay     }
418f17f7ee2SSatish Balay     PetscCall(PetscLogGpuTimeEnd());
419f17f7ee2SSatish Balay     PetscCall(MatDenseCUDARestoreArrayRead(B, (const PetscScalar **)&xx));
420f17f7ee2SSatish Balay     if (usesf) {
421f17f7ee2SSatish Balay       PetscCall(PetscSFReduceBegin(csf, MPIU_SCALAR, uyy, yy, MPI_REPLACE));
422f17f7ee2SSatish Balay       PetscCall(PetscSFReduceEnd(csf, MPIU_SCALAR, uyy, yy, MPI_REPLACE));
423f17f7ee2SSatish Balay     }
424f17f7ee2SSatish Balay     PetscCall(MatDenseCUDARestoreArrayWrite(C, &yy));
42548a46eb9SPierre Jolivet     if (!biscuda) PetscCall(MatConvert(B, MATDENSE, MAT_INPLACE_MATRIX, &B));
42648a46eb9SPierre Jolivet     if (!ciscuda) PetscCall(MatConvert(C, MATDENSE, MAT_INPLACE_MATRIX, &C));
427f17f7ee2SSatish Balay #endif
428f17f7ee2SSatish Balay   }
429f17f7ee2SSatish Balay   { /* log flops */
430f17f7ee2SSatish Balay     double gops, time, perf, dev;
431f17f7ee2SSatish Balay     HLibProfile::getHgemvPerf(gops, time, perf, dev);
432f17f7ee2SSatish Balay #if defined(PETSC_H2OPUS_USE_GPU)
433f17f7ee2SSatish Balay     if (boundtocpu) {
434f17f7ee2SSatish Balay       PetscCall(PetscLogFlops(1e9 * gops));
435f17f7ee2SSatish Balay     } else {
436f17f7ee2SSatish Balay       PetscCall(PetscLogGpuFlops(1e9 * gops));
437f17f7ee2SSatish Balay     }
438f17f7ee2SSatish Balay #else
439f17f7ee2SSatish Balay     PetscCall(PetscLogFlops(1e9 * gops));
440f17f7ee2SSatish Balay #endif
441f17f7ee2SSatish Balay   }
442f17f7ee2SSatish Balay   PetscFunctionReturn(0);
443f17f7ee2SSatish Balay }
444f17f7ee2SSatish Balay 
4459371c9d4SSatish Balay static PetscErrorCode MatProductNumeric_H2OPUS(Mat C) {
446f17f7ee2SSatish Balay   Mat_Product *product = C->product;
447f17f7ee2SSatish Balay 
448f17f7ee2SSatish Balay   PetscFunctionBegin;
449f17f7ee2SSatish Balay   MatCheckProduct(C, 1);
450f17f7ee2SSatish Balay   switch (product->type) {
4519371c9d4SSatish Balay   case MATPRODUCT_AB: PetscCall(MatMultNKernel_H2OPUS(product->A, PETSC_FALSE, product->B, C)); break;
4529371c9d4SSatish Balay   case MATPRODUCT_AtB: PetscCall(MatMultNKernel_H2OPUS(product->A, PETSC_TRUE, product->B, C)); break;
4539371c9d4SSatish Balay   default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "MatProduct type %s is not supported", MatProductTypes[product->type]);
454f17f7ee2SSatish Balay   }
455f17f7ee2SSatish Balay   PetscFunctionReturn(0);
456f17f7ee2SSatish Balay }
457f17f7ee2SSatish Balay 
4589371c9d4SSatish Balay static PetscErrorCode MatProductSymbolic_H2OPUS(Mat C) {
459f17f7ee2SSatish Balay   Mat_Product *product = C->product;
460f17f7ee2SSatish Balay   PetscBool    cisdense;
461f17f7ee2SSatish Balay   Mat          A, B;
462f17f7ee2SSatish Balay 
463f17f7ee2SSatish Balay   PetscFunctionBegin;
464f17f7ee2SSatish Balay   MatCheckProduct(C, 1);
465f17f7ee2SSatish Balay   A = product->A;
466f17f7ee2SSatish Balay   B = product->B;
467f17f7ee2SSatish Balay   switch (product->type) {
468f17f7ee2SSatish Balay   case MATPRODUCT_AB:
469f17f7ee2SSatish Balay     PetscCall(MatSetSizes(C, A->rmap->n, B->cmap->n, A->rmap->N, B->cmap->N));
470f17f7ee2SSatish Balay     PetscCall(MatSetBlockSizesFromMats(C, product->A, product->B));
471f17f7ee2SSatish Balay     PetscCall(PetscObjectTypeCompareAny((PetscObject)C, &cisdense, MATSEQDENSE, MATMPIDENSE, MATSEQDENSECUDA, MATMPIDENSECUDA, ""));
472f17f7ee2SSatish Balay     if (!cisdense) PetscCall(MatSetType(C, ((PetscObject)product->B)->type_name));
473f17f7ee2SSatish Balay     PetscCall(MatSetUp(C));
474f17f7ee2SSatish Balay     break;
475f17f7ee2SSatish Balay   case MATPRODUCT_AtB:
476f17f7ee2SSatish Balay     PetscCall(MatSetSizes(C, A->cmap->n, B->cmap->n, A->cmap->N, B->cmap->N));
477f17f7ee2SSatish Balay     PetscCall(MatSetBlockSizesFromMats(C, product->A, product->B));
478f17f7ee2SSatish Balay     PetscCall(PetscObjectTypeCompareAny((PetscObject)C, &cisdense, MATSEQDENSE, MATMPIDENSE, MATSEQDENSECUDA, MATMPIDENSECUDA, ""));
479f17f7ee2SSatish Balay     if (!cisdense) PetscCall(MatSetType(C, ((PetscObject)product->B)->type_name));
480f17f7ee2SSatish Balay     PetscCall(MatSetUp(C));
481f17f7ee2SSatish Balay     break;
4829371c9d4SSatish Balay   default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "MatProduct type %s is not supported", MatProductTypes[product->type]);
483f17f7ee2SSatish Balay   }
484f17f7ee2SSatish Balay   C->ops->productsymbolic = NULL;
485f17f7ee2SSatish Balay   C->ops->productnumeric  = MatProductNumeric_H2OPUS;
486f17f7ee2SSatish Balay   PetscFunctionReturn(0);
487f17f7ee2SSatish Balay }
488f17f7ee2SSatish Balay 
4899371c9d4SSatish Balay static PetscErrorCode MatProductSetFromOptions_H2OPUS(Mat C) {
490f17f7ee2SSatish Balay   PetscFunctionBegin;
491f17f7ee2SSatish Balay   MatCheckProduct(C, 1);
492ad540459SPierre Jolivet   if (C->product->type == MATPRODUCT_AB || C->product->type == MATPRODUCT_AtB) C->ops->productsymbolic = MatProductSymbolic_H2OPUS;
493f17f7ee2SSatish Balay   PetscFunctionReturn(0);
494f17f7ee2SSatish Balay }
495f17f7ee2SSatish Balay 
4969371c9d4SSatish Balay static PetscErrorCode MatMultKernel_H2OPUS(Mat A, Vec x, PetscScalar sy, Vec y, PetscBool trans) {
497f17f7ee2SSatish Balay   Mat_H2OPUS *h2opus = (Mat_H2OPUS *)A->data;
498f17f7ee2SSatish Balay #if defined(H2OPUS_USE_MPI)
499f17f7ee2SSatish Balay   h2opusHandle_t handle = h2opus->handle->handle;
500f17f7ee2SSatish Balay #else
501f17f7ee2SSatish Balay   h2opusHandle_t handle = h2opus->handle;
502f17f7ee2SSatish Balay #endif
503f17f7ee2SSatish Balay   PetscBool    boundtocpu = PETSC_TRUE;
504f17f7ee2SSatish Balay   PetscInt     n;
505f17f7ee2SSatish Balay   PetscScalar *xx, *yy, *uxx, *uyy;
506f17f7ee2SSatish Balay   PetscMPIInt  size;
507f17f7ee2SSatish Balay   PetscBool    usesf = (PetscBool)(h2opus->sf && !h2opus->nativemult);
508f17f7ee2SSatish Balay 
509f17f7ee2SSatish Balay   PetscFunctionBegin;
510f17f7ee2SSatish Balay   HLibProfile::clear();
511f17f7ee2SSatish Balay   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
512f17f7ee2SSatish Balay #if defined(PETSC_H2OPUS_USE_GPU)
513f17f7ee2SSatish Balay   boundtocpu = A->boundtocpu;
514f17f7ee2SSatish Balay #endif
515f17f7ee2SSatish Balay   if (usesf) {
516f17f7ee2SSatish Balay     PetscCall(PetscSFGetGraph(h2opus->sf, NULL, &n, NULL, NULL));
517f17f7ee2SSatish Balay   } else n = A->rmap->n;
518f17f7ee2SSatish Balay   if (boundtocpu) {
519f17f7ee2SSatish Balay     PetscCall(VecGetArrayRead(x, (const PetscScalar **)&xx));
520f17f7ee2SSatish Balay     if (sy == 0.0) {
521f17f7ee2SSatish Balay       PetscCall(VecGetArrayWrite(y, &yy));
522f17f7ee2SSatish Balay     } else {
523f17f7ee2SSatish Balay       PetscCall(VecGetArray(y, &yy));
524f17f7ee2SSatish Balay     }
525f17f7ee2SSatish Balay     if (usesf) {
526f17f7ee2SSatish Balay       uxx = MatH2OpusGetThrustPointer(*h2opus->xx);
527f17f7ee2SSatish Balay       uyy = MatH2OpusGetThrustPointer(*h2opus->yy);
528f17f7ee2SSatish Balay 
529f17f7ee2SSatish Balay       PetscCall(PetscSFBcastBegin(h2opus->sf, MPIU_SCALAR, xx, uxx, MPI_REPLACE));
530f17f7ee2SSatish Balay       PetscCall(PetscSFBcastEnd(h2opus->sf, MPIU_SCALAR, xx, uxx, MPI_REPLACE));
531f17f7ee2SSatish Balay       if (sy != 0.0) {
532f17f7ee2SSatish Balay         PetscCall(PetscSFBcastBegin(h2opus->sf, MPIU_SCALAR, yy, uyy, MPI_REPLACE));
533f17f7ee2SSatish Balay         PetscCall(PetscSFBcastEnd(h2opus->sf, MPIU_SCALAR, yy, uyy, MPI_REPLACE));
534f17f7ee2SSatish Balay       }
535f17f7ee2SSatish Balay     } else {
536f17f7ee2SSatish Balay       uxx = xx;
537f17f7ee2SSatish Balay       uyy = yy;
538f17f7ee2SSatish Balay     }
539f17f7ee2SSatish Balay     if (size > 1) {
540f17f7ee2SSatish Balay       PetscCheck(h2opus->dist_hmatrix, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing distributed CPU matrix");
541036e4bdcSSatish Balay       PetscCheck(!trans || A->symmetric, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "MatMultTranspose not yet coded in parallel");
542f17f7ee2SSatish Balay #if defined(H2OPUS_USE_MPI)
543f17f7ee2SSatish Balay       distributed_hgemv(/*trans ? H2Opus_Trans : H2Opus_NoTrans, */ h2opus->s, *h2opus->dist_hmatrix, uxx, n, sy, uyy, n, 1, h2opus->handle);
544f17f7ee2SSatish Balay #endif
545f17f7ee2SSatish Balay     } else {
546f17f7ee2SSatish Balay       PetscCheck(h2opus->hmatrix, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing CPU matrix");
547f17f7ee2SSatish Balay       hgemv(trans ? H2Opus_Trans : H2Opus_NoTrans, h2opus->s, *h2opus->hmatrix, uxx, n, sy, uyy, n, 1, handle);
548f17f7ee2SSatish Balay     }
549f17f7ee2SSatish Balay     PetscCall(VecRestoreArrayRead(x, (const PetscScalar **)&xx));
550f17f7ee2SSatish Balay     if (usesf) {
551f17f7ee2SSatish Balay       PetscCall(PetscSFReduceBegin(h2opus->sf, MPIU_SCALAR, uyy, yy, MPI_REPLACE));
552f17f7ee2SSatish Balay       PetscCall(PetscSFReduceEnd(h2opus->sf, MPIU_SCALAR, uyy, yy, MPI_REPLACE));
553f17f7ee2SSatish Balay     }
554f17f7ee2SSatish Balay     if (sy == 0.0) {
555f17f7ee2SSatish Balay       PetscCall(VecRestoreArrayWrite(y, &yy));
556f17f7ee2SSatish Balay     } else {
557f17f7ee2SSatish Balay       PetscCall(VecRestoreArray(y, &yy));
558f17f7ee2SSatish Balay     }
559f17f7ee2SSatish Balay #if defined(PETSC_H2OPUS_USE_GPU)
560f17f7ee2SSatish Balay   } else {
561f17f7ee2SSatish Balay     PetscCall(VecCUDAGetArrayRead(x, (const PetscScalar **)&xx));
562f17f7ee2SSatish Balay     if (sy == 0.0) {
563f17f7ee2SSatish Balay       PetscCall(VecCUDAGetArrayWrite(y, &yy));
564f17f7ee2SSatish Balay     } else {
565f17f7ee2SSatish Balay       PetscCall(VecCUDAGetArray(y, &yy));
566f17f7ee2SSatish Balay     }
567f17f7ee2SSatish Balay     if (usesf) {
568f17f7ee2SSatish Balay       uxx = MatH2OpusGetThrustPointer(*h2opus->xx_gpu);
569f17f7ee2SSatish Balay       uyy = MatH2OpusGetThrustPointer(*h2opus->yy_gpu);
570f17f7ee2SSatish Balay 
571f17f7ee2SSatish Balay       PetscCall(PetscSFBcastBegin(h2opus->sf, MPIU_SCALAR, xx, uxx, MPI_REPLACE));
572f17f7ee2SSatish Balay       PetscCall(PetscSFBcastEnd(h2opus->sf, MPIU_SCALAR, xx, uxx, MPI_REPLACE));
573f17f7ee2SSatish Balay       if (sy != 0.0) {
574f17f7ee2SSatish Balay         PetscCall(PetscSFBcastBegin(h2opus->sf, MPIU_SCALAR, yy, uyy, MPI_REPLACE));
575f17f7ee2SSatish Balay         PetscCall(PetscSFBcastEnd(h2opus->sf, MPIU_SCALAR, yy, uyy, MPI_REPLACE));
576f17f7ee2SSatish Balay       }
577f17f7ee2SSatish Balay     } else {
578f17f7ee2SSatish Balay       uxx = xx;
579f17f7ee2SSatish Balay       uyy = yy;
580f17f7ee2SSatish Balay     }
581f17f7ee2SSatish Balay     PetscCall(PetscLogGpuTimeBegin());
582f17f7ee2SSatish Balay     if (size > 1) {
583f17f7ee2SSatish Balay       PetscCheck(h2opus->dist_hmatrix_gpu, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing distributed GPU matrix");
584036e4bdcSSatish Balay       PetscCheck(!trans || A->symmetric, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "MatMultTranspose not yet coded in parallel");
585f17f7ee2SSatish Balay #if defined(H2OPUS_USE_MPI)
586f17f7ee2SSatish Balay       distributed_hgemv(/*trans ? H2Opus_Trans : H2Opus_NoTrans, */ h2opus->s, *h2opus->dist_hmatrix_gpu, uxx, n, sy, uyy, n, 1, h2opus->handle);
587f17f7ee2SSatish Balay #endif
588f17f7ee2SSatish Balay     } else {
589f17f7ee2SSatish Balay       PetscCheck(h2opus->hmatrix_gpu, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing GPU matrix");
590f17f7ee2SSatish Balay       hgemv(trans ? H2Opus_Trans : H2Opus_NoTrans, h2opus->s, *h2opus->hmatrix_gpu, uxx, n, sy, uyy, n, 1, handle);
591f17f7ee2SSatish Balay     }
592f17f7ee2SSatish Balay     PetscCall(PetscLogGpuTimeEnd());
593f17f7ee2SSatish Balay     PetscCall(VecCUDARestoreArrayRead(x, (const PetscScalar **)&xx));
594f17f7ee2SSatish Balay     if (usesf) {
595f17f7ee2SSatish Balay       PetscCall(PetscSFReduceBegin(h2opus->sf, MPIU_SCALAR, uyy, yy, MPI_REPLACE));
596f17f7ee2SSatish Balay       PetscCall(PetscSFReduceEnd(h2opus->sf, MPIU_SCALAR, uyy, yy, MPI_REPLACE));
597f17f7ee2SSatish Balay     }
598f17f7ee2SSatish Balay     if (sy == 0.0) {
599f17f7ee2SSatish Balay       PetscCall(VecCUDARestoreArrayWrite(y, &yy));
600f17f7ee2SSatish Balay     } else {
601f17f7ee2SSatish Balay       PetscCall(VecCUDARestoreArray(y, &yy));
602f17f7ee2SSatish Balay     }
603f17f7ee2SSatish Balay #endif
604f17f7ee2SSatish Balay   }
605f17f7ee2SSatish Balay   { /* log flops */
606f17f7ee2SSatish Balay     double gops, time, perf, dev;
607f17f7ee2SSatish Balay     HLibProfile::getHgemvPerf(gops, time, perf, dev);
608f17f7ee2SSatish Balay #if defined(PETSC_H2OPUS_USE_GPU)
609f17f7ee2SSatish Balay     if (boundtocpu) {
610f17f7ee2SSatish Balay       PetscCall(PetscLogFlops(1e9 * gops));
611f17f7ee2SSatish Balay     } else {
612f17f7ee2SSatish Balay       PetscCall(PetscLogGpuFlops(1e9 * gops));
613f17f7ee2SSatish Balay     }
614f17f7ee2SSatish Balay #else
615f17f7ee2SSatish Balay     PetscCall(PetscLogFlops(1e9 * gops));
616f17f7ee2SSatish Balay #endif
617f17f7ee2SSatish Balay   }
618f17f7ee2SSatish Balay   PetscFunctionReturn(0);
619f17f7ee2SSatish Balay }
620f17f7ee2SSatish Balay 
6219371c9d4SSatish Balay static PetscErrorCode MatMultTranspose_H2OPUS(Mat A, Vec x, Vec y) {
622f17f7ee2SSatish Balay   PetscBool xiscuda, yiscuda;
623f17f7ee2SSatish Balay 
624f17f7ee2SSatish Balay   PetscFunctionBegin;
625f17f7ee2SSatish Balay   PetscCall(PetscObjectTypeCompareAny((PetscObject)x, &xiscuda, VECSEQCUDA, VECMPICUDA, ""));
626f17f7ee2SSatish Balay   PetscCall(PetscObjectTypeCompareAny((PetscObject)y, &yiscuda, VECSEQCUDA, VECMPICUDA, ""));
627f17f7ee2SSatish Balay   PetscCall(MatH2OpusUpdateIfNeeded(A, !xiscuda || !yiscuda));
628f17f7ee2SSatish Balay   PetscCall(MatMultKernel_H2OPUS(A, x, 0.0, y, PETSC_TRUE));
629f17f7ee2SSatish Balay   PetscFunctionReturn(0);
630f17f7ee2SSatish Balay }
631f17f7ee2SSatish Balay 
6329371c9d4SSatish Balay static PetscErrorCode MatMult_H2OPUS(Mat A, Vec x, Vec y) {
633f17f7ee2SSatish Balay   PetscBool xiscuda, yiscuda;
634f17f7ee2SSatish Balay 
635f17f7ee2SSatish Balay   PetscFunctionBegin;
636f17f7ee2SSatish Balay   PetscCall(PetscObjectTypeCompareAny((PetscObject)x, &xiscuda, VECSEQCUDA, VECMPICUDA, ""));
637f17f7ee2SSatish Balay   PetscCall(PetscObjectTypeCompareAny((PetscObject)y, &yiscuda, VECSEQCUDA, VECMPICUDA, ""));
638f17f7ee2SSatish Balay   PetscCall(MatH2OpusUpdateIfNeeded(A, !xiscuda || !yiscuda));
639f17f7ee2SSatish Balay   PetscCall(MatMultKernel_H2OPUS(A, x, 0.0, y, PETSC_FALSE));
640f17f7ee2SSatish Balay   PetscFunctionReturn(0);
641f17f7ee2SSatish Balay }
642f17f7ee2SSatish Balay 
6439371c9d4SSatish Balay static PetscErrorCode MatMultTransposeAdd_H2OPUS(Mat A, Vec x, Vec y, Vec z) {
644f17f7ee2SSatish Balay   PetscBool xiscuda, ziscuda;
645f17f7ee2SSatish Balay 
646f17f7ee2SSatish Balay   PetscFunctionBegin;
647f17f7ee2SSatish Balay   PetscCall(VecCopy(y, z));
648f17f7ee2SSatish Balay   PetscCall(PetscObjectTypeCompareAny((PetscObject)x, &xiscuda, VECSEQCUDA, VECMPICUDA, ""));
649f17f7ee2SSatish Balay   PetscCall(PetscObjectTypeCompareAny((PetscObject)z, &ziscuda, VECSEQCUDA, VECMPICUDA, ""));
650f17f7ee2SSatish Balay   PetscCall(MatH2OpusUpdateIfNeeded(A, !xiscuda || !ziscuda));
651f17f7ee2SSatish Balay   PetscCall(MatMultKernel_H2OPUS(A, x, 1.0, z, PETSC_TRUE));
652f17f7ee2SSatish Balay   PetscFunctionReturn(0);
653f17f7ee2SSatish Balay }
654f17f7ee2SSatish Balay 
6559371c9d4SSatish Balay static PetscErrorCode MatMultAdd_H2OPUS(Mat A, Vec x, Vec y, Vec z) {
656f17f7ee2SSatish Balay   PetscBool xiscuda, ziscuda;
657f17f7ee2SSatish Balay 
658f17f7ee2SSatish Balay   PetscFunctionBegin;
659f17f7ee2SSatish Balay   PetscCall(VecCopy(y, z));
660f17f7ee2SSatish Balay   PetscCall(PetscObjectTypeCompareAny((PetscObject)x, &xiscuda, VECSEQCUDA, VECMPICUDA, ""));
661f17f7ee2SSatish Balay   PetscCall(PetscObjectTypeCompareAny((PetscObject)z, &ziscuda, VECSEQCUDA, VECMPICUDA, ""));
662f17f7ee2SSatish Balay   PetscCall(MatH2OpusUpdateIfNeeded(A, !xiscuda || !ziscuda));
663f17f7ee2SSatish Balay   PetscCall(MatMultKernel_H2OPUS(A, x, 1.0, z, PETSC_FALSE));
664f17f7ee2SSatish Balay   PetscFunctionReturn(0);
665f17f7ee2SSatish Balay }
666f17f7ee2SSatish Balay 
6679371c9d4SSatish Balay static PetscErrorCode MatScale_H2OPUS(Mat A, PetscScalar s) {
668f17f7ee2SSatish Balay   Mat_H2OPUS *a = (Mat_H2OPUS *)A->data;
669f17f7ee2SSatish Balay 
670f17f7ee2SSatish Balay   PetscFunctionBegin;
671f17f7ee2SSatish Balay   a->s *= s;
672f17f7ee2SSatish Balay   PetscFunctionReturn(0);
673f17f7ee2SSatish Balay }
674f17f7ee2SSatish Balay 
6759371c9d4SSatish Balay static PetscErrorCode MatSetFromOptions_H2OPUS(Mat A, PetscOptionItems *PetscOptionsObject) {
676f17f7ee2SSatish Balay   Mat_H2OPUS *a = (Mat_H2OPUS *)A->data;
677f17f7ee2SSatish Balay 
678f17f7ee2SSatish Balay   PetscFunctionBegin;
679036e4bdcSSatish Balay   PetscOptionsHeadBegin(PetscOptionsObject, "H2OPUS options");
680f17f7ee2SSatish Balay   PetscCall(PetscOptionsInt("-mat_h2opus_leafsize", "Leaf size of cluster tree", NULL, a->leafsize, &a->leafsize, NULL));
681f17f7ee2SSatish Balay   PetscCall(PetscOptionsReal("-mat_h2opus_eta", "Admissibility condition tolerance", NULL, a->eta, &a->eta, NULL));
682f17f7ee2SSatish Balay   PetscCall(PetscOptionsInt("-mat_h2opus_order", "Basis order for off-diagonal sampling when constructed from kernel", NULL, a->basisord, &a->basisord, NULL));
683f17f7ee2SSatish Balay   PetscCall(PetscOptionsInt("-mat_h2opus_maxrank", "Maximum rank when constructed from matvecs", NULL, a->max_rank, &a->max_rank, NULL));
684f17f7ee2SSatish Balay   PetscCall(PetscOptionsInt("-mat_h2opus_samples", "Maximum number of samples to be taken concurrently when constructing from matvecs", NULL, a->bs, &a->bs, NULL));
685f17f7ee2SSatish Balay   PetscCall(PetscOptionsInt("-mat_h2opus_normsamples", "Maximum bumber of samples to be when estimating norms", NULL, a->norm_max_samples, &a->norm_max_samples, NULL));
686f17f7ee2SSatish Balay   PetscCall(PetscOptionsReal("-mat_h2opus_rtol", "Relative tolerance for construction from sampling", NULL, a->rtol, &a->rtol, NULL));
687f17f7ee2SSatish Balay   PetscCall(PetscOptionsBool("-mat_h2opus_check", "Check error when constructing from sampling during MatAssemblyEnd()", NULL, a->check_construction, &a->check_construction, NULL));
688f17f7ee2SSatish Balay   PetscCall(PetscOptionsBool("-mat_h2opus_hara_verbose", "Verbose output from hara construction", NULL, a->hara_verbose, &a->hara_verbose, NULL));
689036e4bdcSSatish Balay   PetscCall(PetscOptionsBool("-mat_h2opus_resize", "Resize after compression", NULL, a->resize, &a->resize, NULL));
690036e4bdcSSatish Balay   PetscOptionsHeadEnd();
691f17f7ee2SSatish Balay   PetscFunctionReturn(0);
692f17f7ee2SSatish Balay }
693f17f7ee2SSatish Balay 
694f17f7ee2SSatish Balay static PetscErrorCode MatH2OpusSetCoords_H2OPUS(Mat, PetscInt, const PetscReal[], PetscBool, MatH2OpusKernel, void *);
695f17f7ee2SSatish Balay 
6969371c9d4SSatish Balay static PetscErrorCode MatH2OpusInferCoordinates_Private(Mat A) {
697f17f7ee2SSatish Balay   Mat_H2OPUS        *a = (Mat_H2OPUS *)A->data;
698f17f7ee2SSatish Balay   Vec                c;
699f17f7ee2SSatish Balay   PetscInt           spacedim;
700f17f7ee2SSatish Balay   const PetscScalar *coords;
701f17f7ee2SSatish Balay 
702f17f7ee2SSatish Balay   PetscFunctionBegin;
703f17f7ee2SSatish Balay   if (a->ptcloud) PetscFunctionReturn(0);
704f17f7ee2SSatish Balay   PetscCall(PetscObjectQuery((PetscObject)A, "__math2opus_coords", (PetscObject *)&c));
705f17f7ee2SSatish Balay   if (!c && a->sampler) {
706f17f7ee2SSatish Balay     Mat S = a->sampler->GetSamplingMat();
707f17f7ee2SSatish Balay 
708f17f7ee2SSatish Balay     PetscCall(PetscObjectQuery((PetscObject)S, "__math2opus_coords", (PetscObject *)&c));
709f17f7ee2SSatish Balay   }
710f17f7ee2SSatish Balay   if (!c) {
711f17f7ee2SSatish Balay     PetscCall(MatH2OpusSetCoords_H2OPUS(A, -1, NULL, PETSC_FALSE, NULL, NULL));
712f17f7ee2SSatish Balay   } else {
713f17f7ee2SSatish Balay     PetscCall(VecGetArrayRead(c, &coords));
714f17f7ee2SSatish Balay     PetscCall(VecGetBlockSize(c, &spacedim));
715f17f7ee2SSatish Balay     PetscCall(MatH2OpusSetCoords_H2OPUS(A, spacedim, coords, PETSC_FALSE, NULL, NULL));
716f17f7ee2SSatish Balay     PetscCall(VecRestoreArrayRead(c, &coords));
717f17f7ee2SSatish Balay   }
718f17f7ee2SSatish Balay   PetscFunctionReturn(0);
719f17f7ee2SSatish Balay }
720f17f7ee2SSatish Balay 
7219371c9d4SSatish Balay static PetscErrorCode MatSetUpMultiply_H2OPUS(Mat A) {
722f17f7ee2SSatish Balay   MPI_Comm      comm;
723f17f7ee2SSatish Balay   PetscMPIInt   size;
724f17f7ee2SSatish Balay   Mat_H2OPUS   *a = (Mat_H2OPUS *)A->data;
725f17f7ee2SSatish Balay   PetscInt      n = 0, *idx = NULL;
726f17f7ee2SSatish Balay   int          *iidx = NULL;
727f17f7ee2SSatish Balay   PetscCopyMode own;
728f17f7ee2SSatish Balay   PetscBool     rid;
729f17f7ee2SSatish Balay 
730f17f7ee2SSatish Balay   PetscFunctionBegin;
731f17f7ee2SSatish Balay   if (a->multsetup) PetscFunctionReturn(0);
732f17f7ee2SSatish Balay   if (a->sf) { /* MatDuplicate_H2OPUS takes reference to the SF */
733f17f7ee2SSatish Balay     PetscCall(PetscSFGetGraph(a->sf, NULL, &n, NULL, NULL));
734f17f7ee2SSatish Balay #if defined(PETSC_H2OPUS_USE_GPU)
735f17f7ee2SSatish Balay     a->xx_gpu  = new thrust::device_vector<PetscScalar>(n);
736f17f7ee2SSatish Balay     a->yy_gpu  = new thrust::device_vector<PetscScalar>(n);
737f17f7ee2SSatish Balay     a->xxs_gpu = 1;
738f17f7ee2SSatish Balay     a->yys_gpu = 1;
739f17f7ee2SSatish Balay #endif
740f17f7ee2SSatish Balay     a->xx  = new thrust::host_vector<PetscScalar>(n);
741f17f7ee2SSatish Balay     a->yy  = new thrust::host_vector<PetscScalar>(n);
742f17f7ee2SSatish Balay     a->xxs = 1;
743f17f7ee2SSatish Balay     a->yys = 1;
744f17f7ee2SSatish Balay   } else {
745f17f7ee2SSatish Balay     IS is;
746f17f7ee2SSatish Balay     PetscCall(PetscObjectGetComm((PetscObject)A, &comm));
747f17f7ee2SSatish Balay     PetscCallMPI(MPI_Comm_size(comm, &size));
748f17f7ee2SSatish Balay     if (!a->h2opus_indexmap) {
749f17f7ee2SSatish Balay       if (size > 1) {
750f17f7ee2SSatish Balay         PetscCheck(a->dist_hmatrix, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing distributed CPU matrix");
751f17f7ee2SSatish Balay #if defined(H2OPUS_USE_MPI)
752f17f7ee2SSatish Balay         iidx = MatH2OpusGetThrustPointer(a->dist_hmatrix->basis_tree.basis_branch.index_map);
753f17f7ee2SSatish Balay         n    = a->dist_hmatrix->basis_tree.basis_branch.index_map.size();
754f17f7ee2SSatish Balay #endif
755f17f7ee2SSatish Balay       } else {
756f17f7ee2SSatish Balay         iidx = MatH2OpusGetThrustPointer(a->hmatrix->u_basis_tree.index_map);
757f17f7ee2SSatish Balay         n    = a->hmatrix->u_basis_tree.index_map.size();
758f17f7ee2SSatish Balay       }
759f17f7ee2SSatish Balay 
760f17f7ee2SSatish Balay       if (PetscDefined(USE_64BIT_INDICES)) {
761f17f7ee2SSatish Balay         PetscInt i;
762f17f7ee2SSatish Balay 
763f17f7ee2SSatish Balay         own = PETSC_OWN_POINTER;
764f17f7ee2SSatish Balay         PetscCall(PetscMalloc1(n, &idx));
765f17f7ee2SSatish Balay         for (i = 0; i < n; i++) idx[i] = iidx[i];
766f17f7ee2SSatish Balay       } else {
767f17f7ee2SSatish Balay         own = PETSC_COPY_VALUES;
768f17f7ee2SSatish Balay         idx = (PetscInt *)iidx;
769f17f7ee2SSatish Balay       }
770f17f7ee2SSatish Balay       PetscCall(ISCreateGeneral(comm, n, idx, own, &is));
771f17f7ee2SSatish Balay       PetscCall(ISSetPermutation(is));
772f17f7ee2SSatish Balay       PetscCall(ISViewFromOptions(is, (PetscObject)A, "-mat_h2opus_indexmap_view"));
773f17f7ee2SSatish Balay       a->h2opus_indexmap = is;
774f17f7ee2SSatish Balay     }
775f17f7ee2SSatish Balay     PetscCall(ISGetLocalSize(a->h2opus_indexmap, &n));
776f17f7ee2SSatish Balay     PetscCall(ISGetIndices(a->h2opus_indexmap, (const PetscInt **)&idx));
777f17f7ee2SSatish Balay     rid = (PetscBool)(n == A->rmap->n);
778f17f7ee2SSatish Balay     PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &rid, 1, MPIU_BOOL, MPI_LAND, comm));
77948a46eb9SPierre Jolivet     if (rid) PetscCall(ISIdentity(a->h2opus_indexmap, &rid));
780f17f7ee2SSatish Balay     if (!rid) {
781f17f7ee2SSatish Balay       if (size > 1) { /* Parallel distribution may be different, save it here for fast path in MatMult (see MatH2OpusSetNativeMult) */
782f17f7ee2SSatish Balay         PetscCall(PetscLayoutCreate(comm, &a->h2opus_rmap));
783f17f7ee2SSatish Balay         PetscCall(PetscLayoutSetLocalSize(a->h2opus_rmap, n));
784f17f7ee2SSatish Balay         PetscCall(PetscLayoutSetUp(a->h2opus_rmap));
785f17f7ee2SSatish Balay         PetscCall(PetscLayoutReference(a->h2opus_rmap, &a->h2opus_cmap));
786f17f7ee2SSatish Balay       }
787f17f7ee2SSatish Balay       PetscCall(PetscSFCreate(comm, &a->sf));
788f17f7ee2SSatish Balay       PetscCall(PetscSFSetGraphLayout(a->sf, A->rmap, n, NULL, PETSC_OWN_POINTER, idx));
789036e4bdcSSatish Balay       PetscCall(PetscSFSetUp(a->sf));
790f17f7ee2SSatish Balay       PetscCall(PetscSFViewFromOptions(a->sf, (PetscObject)A, "-mat_h2opus_sf_view"));
791f17f7ee2SSatish Balay #if defined(PETSC_H2OPUS_USE_GPU)
792f17f7ee2SSatish Balay       a->xx_gpu  = new thrust::device_vector<PetscScalar>(n);
793f17f7ee2SSatish Balay       a->yy_gpu  = new thrust::device_vector<PetscScalar>(n);
794f17f7ee2SSatish Balay       a->xxs_gpu = 1;
795f17f7ee2SSatish Balay       a->yys_gpu = 1;
796f17f7ee2SSatish Balay #endif
797f17f7ee2SSatish Balay       a->xx  = new thrust::host_vector<PetscScalar>(n);
798f17f7ee2SSatish Balay       a->yy  = new thrust::host_vector<PetscScalar>(n);
799f17f7ee2SSatish Balay       a->xxs = 1;
800f17f7ee2SSatish Balay       a->yys = 1;
801f17f7ee2SSatish Balay     }
802f17f7ee2SSatish Balay     PetscCall(ISRestoreIndices(a->h2opus_indexmap, (const PetscInt **)&idx));
803f17f7ee2SSatish Balay   }
804f17f7ee2SSatish Balay   a->multsetup = PETSC_TRUE;
805f17f7ee2SSatish Balay   PetscFunctionReturn(0);
806f17f7ee2SSatish Balay }
807f17f7ee2SSatish Balay 
8089371c9d4SSatish Balay static PetscErrorCode MatAssemblyEnd_H2OPUS(Mat A, MatAssemblyType assemblytype) {
809f17f7ee2SSatish Balay   Mat_H2OPUS *a = (Mat_H2OPUS *)A->data;
810f17f7ee2SSatish Balay #if defined(H2OPUS_USE_MPI)
811f17f7ee2SSatish Balay   h2opusHandle_t handle = a->handle->handle;
812f17f7ee2SSatish Balay #else
813f17f7ee2SSatish Balay   h2opusHandle_t handle = a->handle;
814f17f7ee2SSatish Balay #endif
815f17f7ee2SSatish Balay   PetscBool   kernel       = PETSC_FALSE;
816f17f7ee2SSatish Balay   PetscBool   boundtocpu   = PETSC_TRUE;
817f17f7ee2SSatish Balay   PetscBool   samplingdone = PETSC_FALSE;
818f17f7ee2SSatish Balay   MPI_Comm    comm;
819f17f7ee2SSatish Balay   PetscMPIInt size;
820f17f7ee2SSatish Balay 
821f17f7ee2SSatish Balay   PetscFunctionBegin;
822f17f7ee2SSatish Balay   PetscCall(PetscObjectGetComm((PetscObject)A, &comm));
823036e4bdcSSatish Balay   PetscCheck(A->rmap->n == A->cmap->n, PETSC_COMM_SELF, PETSC_ERR_SUP, "Different row and column local sizes are not supported");
824036e4bdcSSatish Balay   PetscCheck(A->rmap->N == A->cmap->N, comm, PETSC_ERR_SUP, "Rectangular matrices are not supported");
825f17f7ee2SSatish Balay 
826f17f7ee2SSatish Balay   /* XXX */
827f17f7ee2SSatish Balay   a->leafsize = PetscMin(a->leafsize, PetscMin(A->rmap->N, A->cmap->N));
828f17f7ee2SSatish Balay 
829f17f7ee2SSatish Balay   PetscCallMPI(MPI_Comm_size(comm, &size));
830f17f7ee2SSatish Balay   /* TODO REUSABILITY of geometric construction */
831f17f7ee2SSatish Balay   delete a->hmatrix;
832f17f7ee2SSatish Balay   delete a->dist_hmatrix;
833f17f7ee2SSatish Balay #if defined(PETSC_H2OPUS_USE_GPU)
834f17f7ee2SSatish Balay   delete a->hmatrix_gpu;
835f17f7ee2SSatish Balay   delete a->dist_hmatrix_gpu;
836f17f7ee2SSatish Balay #endif
837f17f7ee2SSatish Balay   a->orthogonal = PETSC_FALSE;
838f17f7ee2SSatish Balay 
839f17f7ee2SSatish Balay   /* TODO: other? */
840f17f7ee2SSatish Balay   H2OpusBoxCenterAdmissibility adm(a->eta);
841f17f7ee2SSatish Balay 
842f17f7ee2SSatish Balay   PetscCall(PetscLogEventBegin(MAT_H2Opus_Build, A, 0, 0, 0));
843f17f7ee2SSatish Balay   if (size > 1) {
844f17f7ee2SSatish Balay #if defined(H2OPUS_USE_MPI)
845f17f7ee2SSatish Balay     a->dist_hmatrix = new DistributedHMatrix(A->rmap->n /* ,A->symmetric */);
846f17f7ee2SSatish Balay #else
847f17f7ee2SSatish Balay     a->dist_hmatrix = NULL;
848f17f7ee2SSatish Balay #endif
849b94d7dedSBarry Smith   } else a->hmatrix = new HMatrix(A->rmap->n, A->symmetric == PETSC_BOOL3_TRUE);
850f17f7ee2SSatish Balay   PetscCall(MatH2OpusInferCoordinates_Private(A));
851f17f7ee2SSatish Balay   PetscCheck(a->ptcloud, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Missing pointcloud");
852f17f7ee2SSatish Balay   if (a->kernel) {
853f17f7ee2SSatish Balay     BoxEntryGen<PetscScalar, H2OPUS_HWTYPE_CPU, PetscFunctionGenerator<PetscScalar>> entry_gen(*a->kernel);
854f17f7ee2SSatish Balay     if (size > 1) {
855f17f7ee2SSatish Balay       PetscCheck(a->dist_hmatrix, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing distributed CPU matrix");
856f17f7ee2SSatish Balay #if defined(H2OPUS_USE_MPI)
857f17f7ee2SSatish Balay       buildDistributedHMatrix(*a->dist_hmatrix, a->ptcloud, adm, entry_gen, a->leafsize, a->basisord, a->handle);
858f17f7ee2SSatish Balay #endif
859f17f7ee2SSatish Balay     } else {
860f17f7ee2SSatish Balay       buildHMatrix(*a->hmatrix, a->ptcloud, adm, entry_gen, a->leafsize, a->basisord);
861f17f7ee2SSatish Balay     }
862f17f7ee2SSatish Balay     kernel = PETSC_TRUE;
863f17f7ee2SSatish Balay   } else {
864036e4bdcSSatish Balay     PetscCheck(size <= 1, comm, PETSC_ERR_SUP, "Construction from sampling not supported in parallel");
865f17f7ee2SSatish Balay     buildHMatrixStructure(*a->hmatrix, a->ptcloud, a->leafsize, adm);
866f17f7ee2SSatish Balay   }
867f17f7ee2SSatish Balay   PetscCall(MatSetUpMultiply_H2OPUS(A));
868f17f7ee2SSatish Balay 
869f17f7ee2SSatish Balay #if defined(PETSC_H2OPUS_USE_GPU)
870f17f7ee2SSatish Balay   boundtocpu = A->boundtocpu;
871f17f7ee2SSatish Balay   if (!boundtocpu) {
872f17f7ee2SSatish Balay     if (size > 1) {
873f17f7ee2SSatish Balay       PetscCheck(a->dist_hmatrix, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing distributed CPU matrix");
874f17f7ee2SSatish Balay #if defined(H2OPUS_USE_MPI)
875f17f7ee2SSatish Balay       a->dist_hmatrix_gpu = new DistributedHMatrix_GPU(*a->dist_hmatrix);
876f17f7ee2SSatish Balay #endif
877f17f7ee2SSatish Balay     } else {
878f17f7ee2SSatish Balay       a->hmatrix_gpu = new HMatrix_GPU(*a->hmatrix);
879f17f7ee2SSatish Balay     }
880f17f7ee2SSatish Balay   }
881f17f7ee2SSatish Balay #endif
882f17f7ee2SSatish Balay   if (size == 1) {
883f17f7ee2SSatish Balay     if (!kernel && a->sampler && a->sampler->GetSamplingMat()) {
884f17f7ee2SSatish Balay       PetscReal Anorm;
885f17f7ee2SSatish Balay       bool      verbose;
886f17f7ee2SSatish Balay 
887f17f7ee2SSatish Balay       PetscCall(PetscOptionsGetBool(((PetscObject)A)->options, ((PetscObject)A)->prefix, "-mat_h2opus_hara_verbose", &a->hara_verbose, NULL));
888f17f7ee2SSatish Balay       verbose = a->hara_verbose;
889f17f7ee2SSatish Balay       PetscCall(MatApproximateNorm_Private(a->sampler->GetSamplingMat(), NORM_2, a->norm_max_samples, &Anorm));
890f17f7ee2SSatish Balay       if (a->hara_verbose) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Sampling uses max rank %d, tol %g (%g*%g), %s samples %d\n", a->max_rank, a->rtol * Anorm, a->rtol, Anorm, boundtocpu ? "CPU" : "GPU", a->bs));
891ad540459SPierre Jolivet       if (a->sf && !a->nativemult) a->sampler->SetIndexMap(a->hmatrix->u_basis_tree.index_map.size(), a->hmatrix->u_basis_tree.index_map.data());
892f17f7ee2SSatish Balay       a->sampler->SetStream(handle->getMainStream());
893f17f7ee2SSatish Balay       if (boundtocpu) {
894f17f7ee2SSatish Balay         a->sampler->SetGPUSampling(false);
895f17f7ee2SSatish Balay         hara(a->sampler, *a->hmatrix, a->max_rank, 10 /* TODO */, a->rtol * Anorm, a->bs, handle, verbose);
896f17f7ee2SSatish Balay #if defined(PETSC_H2OPUS_USE_GPU)
897f17f7ee2SSatish Balay       } else {
898f17f7ee2SSatish Balay         a->sampler->SetGPUSampling(true);
899f17f7ee2SSatish Balay         hara(a->sampler, *a->hmatrix_gpu, a->max_rank, 10 /* TODO */, a->rtol * Anorm, a->bs, handle, verbose);
900f17f7ee2SSatish Balay #endif
901f17f7ee2SSatish Balay       }
902f17f7ee2SSatish Balay       samplingdone = PETSC_TRUE;
903f17f7ee2SSatish Balay     }
904f17f7ee2SSatish Balay   }
905f17f7ee2SSatish Balay #if defined(PETSC_H2OPUS_USE_GPU)
906f17f7ee2SSatish Balay   if (!boundtocpu) {
907f17f7ee2SSatish Balay     delete a->hmatrix;
908f17f7ee2SSatish Balay     delete a->dist_hmatrix;
909f17f7ee2SSatish Balay     a->hmatrix      = NULL;
910f17f7ee2SSatish Balay     a->dist_hmatrix = NULL;
911f17f7ee2SSatish Balay   }
912f17f7ee2SSatish Balay   A->offloadmask = boundtocpu ? PETSC_OFFLOAD_CPU : PETSC_OFFLOAD_GPU;
913f17f7ee2SSatish Balay #endif
914f17f7ee2SSatish Balay   PetscCall(PetscLogEventEnd(MAT_H2Opus_Build, A, 0, 0, 0));
915f17f7ee2SSatish Balay 
916f17f7ee2SSatish Balay   if (!a->s) a->s = 1.0;
917f17f7ee2SSatish Balay   A->assembled = PETSC_TRUE;
918f17f7ee2SSatish Balay 
919f17f7ee2SSatish Balay   if (samplingdone) {
920f17f7ee2SSatish Balay     PetscBool check  = a->check_construction;
921f17f7ee2SSatish Balay     PetscBool checke = PETSC_FALSE;
922f17f7ee2SSatish Balay 
923f17f7ee2SSatish Balay     PetscCall(PetscOptionsGetBool(((PetscObject)A)->options, ((PetscObject)A)->prefix, "-mat_h2opus_check", &check, NULL));
924f17f7ee2SSatish Balay     PetscCall(PetscOptionsGetBool(((PetscObject)A)->options, ((PetscObject)A)->prefix, "-mat_h2opus_check_explicit", &checke, NULL));
925f17f7ee2SSatish Balay     if (check) {
926f17f7ee2SSatish Balay       Mat       E, Ae;
927f17f7ee2SSatish Balay       PetscReal n1, ni, n2;
928f17f7ee2SSatish Balay       PetscReal n1A, niA, n2A;
929f17f7ee2SSatish Balay       void (*normfunc)(void);
930f17f7ee2SSatish Balay 
931f17f7ee2SSatish Balay       Ae = a->sampler->GetSamplingMat();
932f17f7ee2SSatish Balay       PetscCall(MatConvert(A, MATSHELL, MAT_INITIAL_MATRIX, &E));
933f17f7ee2SSatish Balay       PetscCall(MatShellSetOperation(E, MATOP_NORM, (void (*)(void))MatNorm_H2OPUS));
934f17f7ee2SSatish Balay       PetscCall(MatAXPY(E, -1.0, Ae, DIFFERENT_NONZERO_PATTERN));
935f17f7ee2SSatish Balay       PetscCall(MatNorm(E, NORM_1, &n1));
936f17f7ee2SSatish Balay       PetscCall(MatNorm(E, NORM_INFINITY, &ni));
937f17f7ee2SSatish Balay       PetscCall(MatNorm(E, NORM_2, &n2));
938f17f7ee2SSatish Balay       if (checke) {
939f17f7ee2SSatish Balay         Mat eA, eE, eAe;
940f17f7ee2SSatish Balay 
941f17f7ee2SSatish Balay         PetscCall(MatComputeOperator(A, MATAIJ, &eA));
942f17f7ee2SSatish Balay         PetscCall(MatComputeOperator(E, MATAIJ, &eE));
943f17f7ee2SSatish Balay         PetscCall(MatComputeOperator(Ae, MATAIJ, &eAe));
944f17f7ee2SSatish Balay         PetscCall(MatChop(eA, PETSC_SMALL));
945f17f7ee2SSatish Balay         PetscCall(MatChop(eE, PETSC_SMALL));
946f17f7ee2SSatish Balay         PetscCall(MatChop(eAe, PETSC_SMALL));
947f17f7ee2SSatish Balay         PetscCall(PetscObjectSetName((PetscObject)eA, "H2Mat"));
948f17f7ee2SSatish Balay         PetscCall(MatView(eA, NULL));
949f17f7ee2SSatish Balay         PetscCall(PetscObjectSetName((PetscObject)eAe, "S"));
950f17f7ee2SSatish Balay         PetscCall(MatView(eAe, NULL));
951f17f7ee2SSatish Balay         PetscCall(PetscObjectSetName((PetscObject)eE, "H2Mat - S"));
952f17f7ee2SSatish Balay         PetscCall(MatView(eE, NULL));
953f17f7ee2SSatish Balay         PetscCall(MatDestroy(&eA));
954f17f7ee2SSatish Balay         PetscCall(MatDestroy(&eE));
955f17f7ee2SSatish Balay         PetscCall(MatDestroy(&eAe));
956f17f7ee2SSatish Balay       }
957f17f7ee2SSatish Balay 
958f17f7ee2SSatish Balay       PetscCall(MatGetOperation(Ae, MATOP_NORM, &normfunc));
959f17f7ee2SSatish Balay       PetscCall(MatSetOperation(Ae, MATOP_NORM, (void (*)(void))MatNorm_H2OPUS));
960f17f7ee2SSatish Balay       PetscCall(MatNorm(Ae, NORM_1, &n1A));
961f17f7ee2SSatish Balay       PetscCall(MatNorm(Ae, NORM_INFINITY, &niA));
962f17f7ee2SSatish Balay       PetscCall(MatNorm(Ae, NORM_2, &n2A));
963f17f7ee2SSatish Balay       n1A = PetscMax(n1A, PETSC_SMALL);
964f17f7ee2SSatish Balay       n2A = PetscMax(n2A, PETSC_SMALL);
965f17f7ee2SSatish Balay       niA = PetscMax(niA, PETSC_SMALL);
966f17f7ee2SSatish Balay       PetscCall(MatSetOperation(Ae, MATOP_NORM, normfunc));
967f17f7ee2SSatish Balay       PetscCall(PetscPrintf(PetscObjectComm((PetscObject)A), "MATH2OPUS construction errors: NORM_1 %g, NORM_INFINITY %g, NORM_2 %g (%g %g %g)\n", (double)n1, (double)ni, (double)n2, (double)(n1 / n1A), (double)(ni / niA), (double)(n2 / n2A)));
968f17f7ee2SSatish Balay       PetscCall(MatDestroy(&E));
969f17f7ee2SSatish Balay     }
970f17f7ee2SSatish Balay     a->sampler->SetSamplingMat(NULL);
971f17f7ee2SSatish Balay   }
972f17f7ee2SSatish Balay   PetscFunctionReturn(0);
973f17f7ee2SSatish Balay }
974f17f7ee2SSatish Balay 
9759371c9d4SSatish Balay static PetscErrorCode MatZeroEntries_H2OPUS(Mat A) {
976f17f7ee2SSatish Balay   PetscMPIInt size;
977f17f7ee2SSatish Balay   Mat_H2OPUS *a = (Mat_H2OPUS *)A->data;
978f17f7ee2SSatish Balay 
979f17f7ee2SSatish Balay   PetscFunctionBegin;
980f17f7ee2SSatish Balay   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
981036e4bdcSSatish Balay   PetscCheck(size <= 1, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Not yet supported");
982f17f7ee2SSatish Balay   a->hmatrix->clearData();
983f17f7ee2SSatish Balay #if defined(PETSC_H2OPUS_USE_GPU)
984f17f7ee2SSatish Balay   if (a->hmatrix_gpu) a->hmatrix_gpu->clearData();
985f17f7ee2SSatish Balay #endif
986f17f7ee2SSatish Balay   PetscFunctionReturn(0);
987f17f7ee2SSatish Balay }
988f17f7ee2SSatish Balay 
9899371c9d4SSatish Balay static PetscErrorCode MatDuplicate_H2OPUS(Mat B, MatDuplicateOption op, Mat *nA) {
990f17f7ee2SSatish Balay   Mat         A;
991f17f7ee2SSatish Balay   Mat_H2OPUS *a, *b = (Mat_H2OPUS *)B->data;
992f17f7ee2SSatish Balay #if defined(PETSC_H2OPUS_USE_GPU)
993f17f7ee2SSatish Balay   PetscBool iscpu = PETSC_FALSE;
994f17f7ee2SSatish Balay #else
995f17f7ee2SSatish Balay   PetscBool iscpu = PETSC_TRUE;
996f17f7ee2SSatish Balay #endif
997f17f7ee2SSatish Balay   MPI_Comm comm;
998f17f7ee2SSatish Balay 
999f17f7ee2SSatish Balay   PetscFunctionBegin;
1000f17f7ee2SSatish Balay   PetscCall(PetscObjectGetComm((PetscObject)B, &comm));
1001f17f7ee2SSatish Balay   PetscCall(MatCreate(comm, &A));
1002f17f7ee2SSatish Balay   PetscCall(MatSetSizes(A, B->rmap->n, B->cmap->n, B->rmap->N, B->cmap->N));
1003f17f7ee2SSatish Balay   PetscCall(MatSetType(A, MATH2OPUS));
1004f17f7ee2SSatish Balay   PetscCall(MatPropagateSymmetryOptions(B, A));
1005f17f7ee2SSatish Balay   a = (Mat_H2OPUS *)A->data;
1006f17f7ee2SSatish Balay 
1007f17f7ee2SSatish Balay   a->eta              = b->eta;
1008f17f7ee2SSatish Balay   a->leafsize         = b->leafsize;
1009f17f7ee2SSatish Balay   a->basisord         = b->basisord;
1010f17f7ee2SSatish Balay   a->max_rank         = b->max_rank;
1011f17f7ee2SSatish Balay   a->bs               = b->bs;
1012f17f7ee2SSatish Balay   a->rtol             = b->rtol;
1013f17f7ee2SSatish Balay   a->norm_max_samples = b->norm_max_samples;
1014f17f7ee2SSatish Balay   if (op == MAT_COPY_VALUES) a->s = b->s;
1015f17f7ee2SSatish Balay 
1016f17f7ee2SSatish Balay   a->ptcloud = new PetscPointCloud<PetscReal>(*b->ptcloud);
1017f17f7ee2SSatish Balay   if (op == MAT_COPY_VALUES && b->kernel) a->kernel = new PetscFunctionGenerator<PetscScalar>(*b->kernel);
1018f17f7ee2SSatish Balay 
1019f17f7ee2SSatish Balay #if defined(H2OPUS_USE_MPI)
1020ad540459SPierre Jolivet   if (b->dist_hmatrix) a->dist_hmatrix = new DistributedHMatrix(*b->dist_hmatrix);
1021f17f7ee2SSatish Balay #if defined(PETSC_H2OPUS_USE_GPU)
1022ad540459SPierre Jolivet   if (b->dist_hmatrix_gpu) a->dist_hmatrix_gpu = new DistributedHMatrix_GPU(*b->dist_hmatrix_gpu);
1023f17f7ee2SSatish Balay #endif
1024f17f7ee2SSatish Balay #endif
1025f17f7ee2SSatish Balay   if (b->hmatrix) {
1026f17f7ee2SSatish Balay     a->hmatrix = new HMatrix(*b->hmatrix);
1027f17f7ee2SSatish Balay     if (op == MAT_DO_NOT_COPY_VALUES) a->hmatrix->clearData();
1028f17f7ee2SSatish Balay   }
1029f17f7ee2SSatish Balay #if defined(PETSC_H2OPUS_USE_GPU)
1030f17f7ee2SSatish Balay   if (b->hmatrix_gpu) {
1031f17f7ee2SSatish Balay     a->hmatrix_gpu = new HMatrix_GPU(*b->hmatrix_gpu);
1032f17f7ee2SSatish Balay     if (op == MAT_DO_NOT_COPY_VALUES) a->hmatrix_gpu->clearData();
1033f17f7ee2SSatish Balay   }
1034f17f7ee2SSatish Balay #endif
1035f17f7ee2SSatish Balay   if (b->sf) {
1036f17f7ee2SSatish Balay     PetscCall(PetscObjectReference((PetscObject)b->sf));
1037f17f7ee2SSatish Balay     a->sf = b->sf;
1038f17f7ee2SSatish Balay   }
1039f17f7ee2SSatish Balay   if (b->h2opus_indexmap) {
1040f17f7ee2SSatish Balay     PetscCall(PetscObjectReference((PetscObject)b->h2opus_indexmap));
1041f17f7ee2SSatish Balay     a->h2opus_indexmap = b->h2opus_indexmap;
1042f17f7ee2SSatish Balay   }
1043f17f7ee2SSatish Balay 
1044f17f7ee2SSatish Balay   PetscCall(MatSetUp(A));
1045f17f7ee2SSatish Balay   PetscCall(MatSetUpMultiply_H2OPUS(A));
1046f17f7ee2SSatish Balay   if (op == MAT_COPY_VALUES) {
1047f17f7ee2SSatish Balay     A->assembled  = PETSC_TRUE;
1048f17f7ee2SSatish Balay     a->orthogonal = b->orthogonal;
1049f17f7ee2SSatish Balay #if defined(PETSC_H2OPUS_USE_GPU)
1050f17f7ee2SSatish Balay     A->offloadmask = B->offloadmask;
1051f17f7ee2SSatish Balay #endif
1052f17f7ee2SSatish Balay   }
1053f17f7ee2SSatish Balay #if defined(PETSC_H2OPUS_USE_GPU)
1054f17f7ee2SSatish Balay   iscpu = B->boundtocpu;
1055f17f7ee2SSatish Balay #endif
1056f17f7ee2SSatish Balay   PetscCall(MatBindToCPU(A, iscpu));
1057f17f7ee2SSatish Balay 
1058f17f7ee2SSatish Balay   *nA = A;
1059f17f7ee2SSatish Balay   PetscFunctionReturn(0);
1060f17f7ee2SSatish Balay }
1061f17f7ee2SSatish Balay 
10629371c9d4SSatish Balay static PetscErrorCode MatView_H2OPUS(Mat A, PetscViewer view) {
1063f17f7ee2SSatish Balay   Mat_H2OPUS       *h2opus = (Mat_H2OPUS *)A->data;
1064036e4bdcSSatish Balay   PetscBool         isascii, vieweps;
1065f17f7ee2SSatish Balay   PetscMPIInt       size;
1066f17f7ee2SSatish Balay   PetscViewerFormat format;
1067f17f7ee2SSatish Balay 
1068f17f7ee2SSatish Balay   PetscFunctionBegin;
1069f17f7ee2SSatish Balay   PetscCall(PetscObjectTypeCompare((PetscObject)view, PETSCVIEWERASCII, &isascii));
1070f17f7ee2SSatish Balay   PetscCall(PetscViewerGetFormat(view, &format));
1071f17f7ee2SSatish Balay   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
1072f17f7ee2SSatish Balay   if (isascii) {
1073f17f7ee2SSatish Balay     if (format == PETSC_VIEWER_ASCII_MATLAB) {
1074f17f7ee2SSatish Balay       if (size == 1) {
1075f17f7ee2SSatish Balay         FILE *fp;
1076f17f7ee2SSatish Balay         PetscCall(PetscViewerASCIIGetPointer(view, &fp));
1077f17f7ee2SSatish Balay         dumpHMatrix(*h2opus->hmatrix, 6, fp);
1078f17f7ee2SSatish Balay       }
1079f17f7ee2SSatish Balay     } else {
1080f17f7ee2SSatish Balay       PetscCall(PetscViewerASCIIPrintf(view, "  H-Matrix constructed from %s\n", h2opus->kernel ? "Kernel" : "Mat"));
1081f17f7ee2SSatish Balay       PetscCall(PetscViewerASCIIPrintf(view, "  PointCloud dim %" PetscInt_FMT "\n", h2opus->ptcloud ? h2opus->ptcloud->getDimension() : 0));
1082f17f7ee2SSatish Balay       PetscCall(PetscViewerASCIIPrintf(view, "  Admissibility parameters: leaf size %" PetscInt_FMT ", eta %g\n", h2opus->leafsize, (double)h2opus->eta));
1083f17f7ee2SSatish Balay       if (!h2opus->kernel) {
1084f17f7ee2SSatish Balay         PetscCall(PetscViewerASCIIPrintf(view, "  Sampling parameters: max_rank %" PetscInt_FMT ", samples %" PetscInt_FMT ", tolerance %g\n", h2opus->max_rank, h2opus->bs, (double)h2opus->rtol));
1085f17f7ee2SSatish Balay       } else {
1086f17f7ee2SSatish Balay         PetscCall(PetscViewerASCIIPrintf(view, "  Offdiagonal blocks approximation order %" PetscInt_FMT "\n", h2opus->basisord));
1087f17f7ee2SSatish Balay       }
1088f17f7ee2SSatish Balay       PetscCall(PetscViewerASCIIPrintf(view, "  Number of samples for norms %" PetscInt_FMT "\n", h2opus->norm_max_samples));
1089f17f7ee2SSatish Balay       if (size == 1) {
1090f17f7ee2SSatish Balay         double dense_mem_cpu = h2opus->hmatrix ? h2opus->hmatrix->getDenseMemoryUsage() : 0;
1091f17f7ee2SSatish Balay         double low_rank_cpu  = h2opus->hmatrix ? h2opus->hmatrix->getLowRankMemoryUsage() : 0;
1092f17f7ee2SSatish Balay #if defined(PETSC_HAVE_CUDA)
1093f17f7ee2SSatish Balay         double dense_mem_gpu = h2opus->hmatrix_gpu ? h2opus->hmatrix_gpu->getDenseMemoryUsage() : 0;
1094f17f7ee2SSatish Balay         double low_rank_gpu  = h2opus->hmatrix_gpu ? h2opus->hmatrix_gpu->getLowRankMemoryUsage() : 0;
1095f17f7ee2SSatish Balay #endif
1096f17f7ee2SSatish Balay         PetscCall(PetscViewerASCIIPrintf(view, "  Memory consumption GB (CPU): %g (dense) %g (low rank) %g (total)\n", dense_mem_cpu, low_rank_cpu, low_rank_cpu + dense_mem_cpu));
1097f17f7ee2SSatish Balay #if defined(PETSC_HAVE_CUDA)
1098f17f7ee2SSatish Balay         PetscCall(PetscViewerASCIIPrintf(view, "  Memory consumption GB (GPU): %g (dense) %g (low rank) %g (total)\n", dense_mem_gpu, low_rank_gpu, low_rank_gpu + dense_mem_gpu));
1099f17f7ee2SSatish Balay #endif
1100f17f7ee2SSatish Balay       } else {
1101f17f7ee2SSatish Balay #if defined(PETSC_HAVE_CUDA)
1102f17f7ee2SSatish Balay         double      matrix_mem[4] = {0., 0., 0., 0.};
1103f17f7ee2SSatish Balay         PetscMPIInt rsize         = 4;
1104f17f7ee2SSatish Balay #else
1105f17f7ee2SSatish Balay         double      matrix_mem[2] = {0., 0.};
1106f17f7ee2SSatish Balay         PetscMPIInt rsize         = 2;
1107f17f7ee2SSatish Balay #endif
1108f17f7ee2SSatish Balay #if defined(H2OPUS_USE_MPI)
1109f17f7ee2SSatish Balay         matrix_mem[0] = h2opus->dist_hmatrix ? h2opus->dist_hmatrix->getLocalDenseMemoryUsage() : 0;
1110f17f7ee2SSatish Balay         matrix_mem[1] = h2opus->dist_hmatrix ? h2opus->dist_hmatrix->getLocalLowRankMemoryUsage() : 0;
1111f17f7ee2SSatish Balay #if defined(PETSC_HAVE_CUDA)
1112f17f7ee2SSatish Balay         matrix_mem[2] = h2opus->dist_hmatrix_gpu ? h2opus->dist_hmatrix_gpu->getLocalDenseMemoryUsage() : 0;
1113f17f7ee2SSatish Balay         matrix_mem[3] = h2opus->dist_hmatrix_gpu ? h2opus->dist_hmatrix_gpu->getLocalLowRankMemoryUsage() : 0;
1114f17f7ee2SSatish Balay #endif
1115f17f7ee2SSatish Balay #endif
1116f17f7ee2SSatish Balay         PetscCall(MPIU_Allreduce(MPI_IN_PLACE, matrix_mem, rsize, MPI_DOUBLE_PRECISION, MPI_SUM, PetscObjectComm((PetscObject)A)));
1117f17f7ee2SSatish Balay         PetscCall(PetscViewerASCIIPrintf(view, "  Memory consumption GB (CPU): %g (dense) %g (low rank) %g (total)\n", matrix_mem[0], matrix_mem[1], matrix_mem[0] + matrix_mem[1]));
1118f17f7ee2SSatish Balay #if defined(PETSC_HAVE_CUDA)
1119f17f7ee2SSatish Balay         PetscCall(PetscViewerASCIIPrintf(view, "  Memory consumption GB (GPU): %g (dense) %g (low rank) %g (total)\n", matrix_mem[2], matrix_mem[3], matrix_mem[2] + matrix_mem[3]));
1120f17f7ee2SSatish Balay #endif
1121f17f7ee2SSatish Balay       }
1122f17f7ee2SSatish Balay     }
1123f17f7ee2SSatish Balay   }
1124036e4bdcSSatish Balay   vieweps = PETSC_FALSE;
1125036e4bdcSSatish Balay   PetscCall(PetscOptionsGetBool(((PetscObject)A)->options, ((PetscObject)A)->prefix, "-mat_h2opus_vieweps", &vieweps, NULL));
1126036e4bdcSSatish Balay   if (vieweps) {
1127f17f7ee2SSatish Balay     char        filename[256];
1128f17f7ee2SSatish Balay     const char *name;
1129f17f7ee2SSatish Balay 
1130f17f7ee2SSatish Balay     PetscCall(PetscObjectGetName((PetscObject)A, &name));
1131f17f7ee2SSatish Balay     PetscCall(PetscSNPrintf(filename, sizeof(filename), "%s_structure.eps", name));
1132036e4bdcSSatish Balay     PetscCall(PetscOptionsGetString(((PetscObject)A)->options, ((PetscObject)A)->prefix, "-mat_h2opus_vieweps_filename", filename, sizeof(filename), NULL));
1133f17f7ee2SSatish Balay     outputEps(*h2opus->hmatrix, filename);
1134f17f7ee2SSatish Balay   }
1135f17f7ee2SSatish Balay   PetscFunctionReturn(0);
1136f17f7ee2SSatish Balay }
1137f17f7ee2SSatish Balay 
11389371c9d4SSatish Balay static PetscErrorCode MatH2OpusSetCoords_H2OPUS(Mat A, PetscInt spacedim, const PetscReal coords[], PetscBool cdist, MatH2OpusKernel kernel, void *kernelctx) {
1139f17f7ee2SSatish Balay   Mat_H2OPUS *h2opus = (Mat_H2OPUS *)A->data;
1140f17f7ee2SSatish Balay   PetscReal  *gcoords;
1141f17f7ee2SSatish Balay   PetscInt    N;
1142f17f7ee2SSatish Balay   MPI_Comm    comm;
1143f17f7ee2SSatish Balay   PetscMPIInt size;
1144f17f7ee2SSatish Balay   PetscBool   cong;
1145f17f7ee2SSatish Balay 
1146f17f7ee2SSatish Balay   PetscFunctionBegin;
1147f17f7ee2SSatish Balay   PetscCall(PetscLayoutSetUp(A->rmap));
1148f17f7ee2SSatish Balay   PetscCall(PetscLayoutSetUp(A->cmap));
1149f17f7ee2SSatish Balay   PetscCall(PetscObjectGetComm((PetscObject)A, &comm));
1150f17f7ee2SSatish Balay   PetscCall(MatHasCongruentLayouts(A, &cong));
1151f17f7ee2SSatish Balay   PetscCheck(cong, comm, PETSC_ERR_SUP, "Only for square matrices with congruent layouts");
1152f17f7ee2SSatish Balay   N = A->rmap->N;
1153f17f7ee2SSatish Balay   PetscCallMPI(MPI_Comm_size(comm, &size));
1154f17f7ee2SSatish Balay   if (spacedim > 0 && size > 1 && cdist) {
1155f17f7ee2SSatish Balay     PetscSF      sf;
1156f17f7ee2SSatish Balay     MPI_Datatype dtype;
1157f17f7ee2SSatish Balay 
1158f17f7ee2SSatish Balay     PetscCallMPI(MPI_Type_contiguous(spacedim, MPIU_REAL, &dtype));
1159f17f7ee2SSatish Balay     PetscCallMPI(MPI_Type_commit(&dtype));
1160f17f7ee2SSatish Balay 
1161f17f7ee2SSatish Balay     PetscCall(PetscSFCreate(comm, &sf));
1162f17f7ee2SSatish Balay     PetscCall(PetscSFSetGraphWithPattern(sf, A->rmap, PETSCSF_PATTERN_ALLGATHER));
1163f17f7ee2SSatish Balay     PetscCall(PetscMalloc1(spacedim * N, &gcoords));
1164f17f7ee2SSatish Balay     PetscCall(PetscSFBcastBegin(sf, dtype, coords, gcoords, MPI_REPLACE));
1165f17f7ee2SSatish Balay     PetscCall(PetscSFBcastEnd(sf, dtype, coords, gcoords, MPI_REPLACE));
1166f17f7ee2SSatish Balay     PetscCall(PetscSFDestroy(&sf));
1167f17f7ee2SSatish Balay     PetscCallMPI(MPI_Type_free(&dtype));
1168f17f7ee2SSatish Balay   } else gcoords = (PetscReal *)coords;
1169f17f7ee2SSatish Balay 
1170f17f7ee2SSatish Balay   delete h2opus->ptcloud;
1171f17f7ee2SSatish Balay   delete h2opus->kernel;
1172f17f7ee2SSatish Balay   h2opus->ptcloud = new PetscPointCloud<PetscReal>(spacedim, N, gcoords);
1173f17f7ee2SSatish Balay   if (kernel) h2opus->kernel = new PetscFunctionGenerator<PetscScalar>(kernel, spacedim, kernelctx);
1174f17f7ee2SSatish Balay   if (gcoords != coords) PetscCall(PetscFree(gcoords));
1175f17f7ee2SSatish Balay   A->preallocated = PETSC_TRUE;
1176f17f7ee2SSatish Balay   PetscFunctionReturn(0);
1177f17f7ee2SSatish Balay }
1178f17f7ee2SSatish Balay 
1179f17f7ee2SSatish Balay #if defined(PETSC_H2OPUS_USE_GPU)
11809371c9d4SSatish Balay static PetscErrorCode MatBindToCPU_H2OPUS(Mat A, PetscBool flg) {
1181f17f7ee2SSatish Balay   PetscMPIInt size;
1182f17f7ee2SSatish Balay   Mat_H2OPUS *a = (Mat_H2OPUS *)A->data;
1183f17f7ee2SSatish Balay 
1184f17f7ee2SSatish Balay   PetscFunctionBegin;
1185f17f7ee2SSatish Balay   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
1186f17f7ee2SSatish Balay   if (flg && A->offloadmask == PETSC_OFFLOAD_GPU) {
1187f17f7ee2SSatish Balay     if (size > 1) {
1188f17f7ee2SSatish Balay       PetscCheck(a->dist_hmatrix_gpu, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing GPU matrix");
1189f17f7ee2SSatish Balay #if defined(H2OPUS_USE_MPI)
1190f17f7ee2SSatish Balay       if (!a->dist_hmatrix) a->dist_hmatrix = new DistributedHMatrix(*a->dist_hmatrix_gpu);
1191f17f7ee2SSatish Balay       else *a->dist_hmatrix = *a->dist_hmatrix_gpu;
1192f17f7ee2SSatish Balay #endif
1193f17f7ee2SSatish Balay     } else {
1194f17f7ee2SSatish Balay       PetscCheck(a->hmatrix_gpu, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing GPU matrix");
1195f17f7ee2SSatish Balay       if (!a->hmatrix) a->hmatrix = new HMatrix(*a->hmatrix_gpu);
1196f17f7ee2SSatish Balay       else *a->hmatrix = *a->hmatrix_gpu;
1197f17f7ee2SSatish Balay     }
1198f17f7ee2SSatish Balay     delete a->hmatrix_gpu;
1199f17f7ee2SSatish Balay     delete a->dist_hmatrix_gpu;
1200f17f7ee2SSatish Balay     a->hmatrix_gpu      = NULL;
1201f17f7ee2SSatish Balay     a->dist_hmatrix_gpu = NULL;
1202f17f7ee2SSatish Balay     A->offloadmask      = PETSC_OFFLOAD_CPU;
1203f17f7ee2SSatish Balay   } else if (!flg && A->offloadmask == PETSC_OFFLOAD_CPU) {
1204f17f7ee2SSatish Balay     if (size > 1) {
1205f17f7ee2SSatish Balay       PetscCheck(a->dist_hmatrix, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing CPU matrix");
1206f17f7ee2SSatish Balay #if defined(H2OPUS_USE_MPI)
1207f17f7ee2SSatish Balay       if (!a->dist_hmatrix_gpu) a->dist_hmatrix_gpu = new DistributedHMatrix_GPU(*a->dist_hmatrix);
1208f17f7ee2SSatish Balay       else *a->dist_hmatrix_gpu = *a->dist_hmatrix;
1209f17f7ee2SSatish Balay #endif
1210f17f7ee2SSatish Balay     } else {
1211f17f7ee2SSatish Balay       PetscCheck(a->hmatrix, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing CPU matrix");
1212f17f7ee2SSatish Balay       if (!a->hmatrix_gpu) a->hmatrix_gpu = new HMatrix_GPU(*a->hmatrix);
1213f17f7ee2SSatish Balay       else *a->hmatrix_gpu = *a->hmatrix;
1214f17f7ee2SSatish Balay     }
1215f17f7ee2SSatish Balay     delete a->hmatrix;
1216f17f7ee2SSatish Balay     delete a->dist_hmatrix;
1217f17f7ee2SSatish Balay     a->hmatrix      = NULL;
1218f17f7ee2SSatish Balay     a->dist_hmatrix = NULL;
1219f17f7ee2SSatish Balay     A->offloadmask  = PETSC_OFFLOAD_GPU;
1220f17f7ee2SSatish Balay   }
1221f17f7ee2SSatish Balay   PetscCall(PetscFree(A->defaultvectype));
1222f17f7ee2SSatish Balay   if (!flg) {
1223f17f7ee2SSatish Balay     PetscCall(PetscStrallocpy(VECCUDA, &A->defaultvectype));
1224f17f7ee2SSatish Balay   } else {
1225f17f7ee2SSatish Balay     PetscCall(PetscStrallocpy(VECSTANDARD, &A->defaultvectype));
1226f17f7ee2SSatish Balay   }
1227f17f7ee2SSatish Balay   A->boundtocpu = flg;
1228f17f7ee2SSatish Balay   PetscFunctionReturn(0);
1229f17f7ee2SSatish Balay }
1230f17f7ee2SSatish Balay #endif
1231f17f7ee2SSatish Balay 
1232f17f7ee2SSatish Balay /*MC
1233f17f7ee2SSatish Balay      MATH2OPUS = "h2opus" - A matrix type for hierarchical matrices using the H2Opus package.
1234f17f7ee2SSatish Balay 
1235f17f7ee2SSatish Balay    Options Database Keys:
1236f17f7ee2SSatish Balay .     -mat_type h2opus - matrix type to "h2opus" during a call to MatSetFromOptions()
1237f17f7ee2SSatish Balay 
1238f17f7ee2SSatish Balay    Notes:
123911a5261eSBarry Smith      H2Opus implements hierarchical matrices in the H^2 flavour. It supports CPU or NVIDIA GPUs.
124011a5261eSBarry Smith 
1241f17f7ee2SSatish Balay      For CPU only builds, use ./configure --download-h2opus --download-thrust to install PETSc to use H2Opus.
1242f17f7ee2SSatish Balay      In order to run on NVIDIA GPUs, use ./configure --download-h2opus --download-magma --download-kblas.
1243f17f7ee2SSatish Balay 
1244f17f7ee2SSatish Balay    Level: beginner
1245f17f7ee2SSatish Balay 
124611a5261eSBarry Smith    Reference:
124711a5261eSBarry Smith .  * -  "H2Opus: A distributed-memory multi-GPU software package for non-local operators", https://arxiv.org/abs/2109.05451
124811a5261eSBarry Smith 
124911a5261eSBarry Smith .seealso: `MATH2OPUS`, `MATHTOOL`, `MATDENSE`, `MatCreateH2OpusFromKernel()`, `MatCreateH2OpusFromMat()`
1250f17f7ee2SSatish Balay M*/
12519371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode MatCreate_H2OPUS(Mat A) {
1252f17f7ee2SSatish Balay   Mat_H2OPUS *a;
1253f17f7ee2SSatish Balay   PetscMPIInt size;
1254f17f7ee2SSatish Balay 
1255f17f7ee2SSatish Balay   PetscFunctionBegin;
1256f17f7ee2SSatish Balay #if defined(PETSC_H2OPUS_USE_GPU)
1257f17f7ee2SSatish Balay   PetscCall(PetscDeviceInitialize(PETSC_DEVICE_CUDA));
1258f17f7ee2SSatish Balay #endif
1259*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&a));
1260f17f7ee2SSatish Balay   A->data = (void *)a;
1261f17f7ee2SSatish Balay 
1262f17f7ee2SSatish Balay   a->eta              = 0.9;
1263f17f7ee2SSatish Balay   a->leafsize         = 32;
1264f17f7ee2SSatish Balay   a->basisord         = 4;
1265f17f7ee2SSatish Balay   a->max_rank         = 64;
1266f17f7ee2SSatish Balay   a->bs               = 32;
1267f17f7ee2SSatish Balay   a->rtol             = 1.e-4;
1268f17f7ee2SSatish Balay   a->s                = 1.0;
1269f17f7ee2SSatish Balay   a->norm_max_samples = 10;
1270036e4bdcSSatish Balay   a->resize           = PETSC_TRUE; /* reallocate after compression */
1271f17f7ee2SSatish Balay #if defined(H2OPUS_USE_MPI)
1272f17f7ee2SSatish Balay   h2opusCreateDistributedHandleComm(&a->handle, PetscObjectComm((PetscObject)A));
1273f17f7ee2SSatish Balay #else
1274f17f7ee2SSatish Balay   h2opusCreateHandle(&a->handle);
1275f17f7ee2SSatish Balay #endif
1276f17f7ee2SSatish Balay   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
1277f17f7ee2SSatish Balay   PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATH2OPUS));
1278f17f7ee2SSatish Balay   PetscCall(PetscMemzero(A->ops, sizeof(struct _MatOps)));
1279f17f7ee2SSatish Balay 
1280f17f7ee2SSatish Balay   A->ops->destroy          = MatDestroy_H2OPUS;
1281f17f7ee2SSatish Balay   A->ops->view             = MatView_H2OPUS;
1282f17f7ee2SSatish Balay   A->ops->assemblyend      = MatAssemblyEnd_H2OPUS;
1283f17f7ee2SSatish Balay   A->ops->mult             = MatMult_H2OPUS;
1284f17f7ee2SSatish Balay   A->ops->multtranspose    = MatMultTranspose_H2OPUS;
1285f17f7ee2SSatish Balay   A->ops->multadd          = MatMultAdd_H2OPUS;
1286f17f7ee2SSatish Balay   A->ops->multtransposeadd = MatMultTransposeAdd_H2OPUS;
1287f17f7ee2SSatish Balay   A->ops->scale            = MatScale_H2OPUS;
1288f17f7ee2SSatish Balay   A->ops->duplicate        = MatDuplicate_H2OPUS;
1289f17f7ee2SSatish Balay   A->ops->setfromoptions   = MatSetFromOptions_H2OPUS;
1290f17f7ee2SSatish Balay   A->ops->norm             = MatNorm_H2OPUS;
1291f17f7ee2SSatish Balay   A->ops->zeroentries      = MatZeroEntries_H2OPUS;
1292f17f7ee2SSatish Balay #if defined(PETSC_H2OPUS_USE_GPU)
1293f17f7ee2SSatish Balay   A->ops->bindtocpu = MatBindToCPU_H2OPUS;
1294f17f7ee2SSatish Balay #endif
1295f17f7ee2SSatish Balay 
1296f17f7ee2SSatish Balay   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_h2opus_seqdense_C", MatProductSetFromOptions_H2OPUS));
1297f17f7ee2SSatish Balay   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_h2opus_seqdensecuda_C", MatProductSetFromOptions_H2OPUS));
1298f17f7ee2SSatish Balay   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_h2opus_mpidense_C", MatProductSetFromOptions_H2OPUS));
1299f17f7ee2SSatish Balay   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_h2opus_mpidensecuda_C", MatProductSetFromOptions_H2OPUS));
1300f17f7ee2SSatish Balay #if defined(PETSC_H2OPUS_USE_GPU)
1301f17f7ee2SSatish Balay   PetscCall(PetscFree(A->defaultvectype));
1302f17f7ee2SSatish Balay   PetscCall(PetscStrallocpy(VECCUDA, &A->defaultvectype));
1303f17f7ee2SSatish Balay #endif
1304f17f7ee2SSatish Balay   PetscFunctionReturn(0);
1305f17f7ee2SSatish Balay }
1306f17f7ee2SSatish Balay 
1307f17f7ee2SSatish Balay /*@C
1308f17f7ee2SSatish Balay      MatH2OpusOrthogonalize - Orthogonalize the basis tree of a hierarchical matrix.
1309f17f7ee2SSatish Balay 
1310f17f7ee2SSatish Balay    Input Parameter:
1311f17f7ee2SSatish Balay .     A - the matrix
1312f17f7ee2SSatish Balay 
1313f17f7ee2SSatish Balay    Level: intermediate
1314f17f7ee2SSatish Balay 
1315036e4bdcSSatish Balay .seealso: `MatCreate()`, `MATH2OPUS`, `MatCreateH2OpusFromMat()`, `MatCreateH2OpusFromKernel()`, `MatH2OpusCompress()`
1316f17f7ee2SSatish Balay @*/
13179371c9d4SSatish Balay PetscErrorCode MatH2OpusOrthogonalize(Mat A) {
1318f17f7ee2SSatish Balay   PetscBool   ish2opus;
1319f17f7ee2SSatish Balay   Mat_H2OPUS *a = (Mat_H2OPUS *)A->data;
1320f17f7ee2SSatish Balay   PetscMPIInt size;
1321f17f7ee2SSatish Balay   PetscBool   boundtocpu = PETSC_TRUE;
1322f17f7ee2SSatish Balay 
1323f17f7ee2SSatish Balay   PetscFunctionBegin;
1324f17f7ee2SSatish Balay   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
1325f17f7ee2SSatish Balay   PetscValidType(A, 1);
1326f17f7ee2SSatish Balay   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATH2OPUS, &ish2opus));
1327f17f7ee2SSatish Balay   if (!ish2opus) PetscFunctionReturn(0);
1328f17f7ee2SSatish Balay   if (a->orthogonal) PetscFunctionReturn(0);
1329f17f7ee2SSatish Balay   HLibProfile::clear();
1330f17f7ee2SSatish Balay   PetscCall(PetscLogEventBegin(MAT_H2Opus_Orthog, A, 0, 0, 0));
1331f17f7ee2SSatish Balay #if defined(PETSC_H2OPUS_USE_GPU)
1332f17f7ee2SSatish Balay   boundtocpu = A->boundtocpu;
1333f17f7ee2SSatish Balay #endif
1334f17f7ee2SSatish Balay   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
1335f17f7ee2SSatish Balay   if (size > 1) {
1336f17f7ee2SSatish Balay     if (boundtocpu) {
1337f17f7ee2SSatish Balay       PetscCheck(a->dist_hmatrix, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing CPU matrix");
1338f17f7ee2SSatish Balay #if defined(H2OPUS_USE_MPI)
1339f17f7ee2SSatish Balay       distributed_horthog(*a->dist_hmatrix, a->handle);
1340f17f7ee2SSatish Balay #endif
1341f17f7ee2SSatish Balay #if defined(PETSC_H2OPUS_USE_GPU)
1342f17f7ee2SSatish Balay       A->offloadmask = PETSC_OFFLOAD_CPU;
1343f17f7ee2SSatish Balay     } else {
1344f17f7ee2SSatish Balay       PetscCheck(a->dist_hmatrix_gpu, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing GPU matrix");
1345f17f7ee2SSatish Balay       PetscCall(PetscLogGpuTimeBegin());
1346f17f7ee2SSatish Balay #if defined(H2OPUS_USE_MPI)
1347f17f7ee2SSatish Balay       distributed_horthog(*a->dist_hmatrix_gpu, a->handle);
1348f17f7ee2SSatish Balay #endif
1349f17f7ee2SSatish Balay       PetscCall(PetscLogGpuTimeEnd());
1350f17f7ee2SSatish Balay #endif
1351f17f7ee2SSatish Balay     }
1352f17f7ee2SSatish Balay   } else {
1353f17f7ee2SSatish Balay #if defined(H2OPUS_USE_MPI)
1354f17f7ee2SSatish Balay     h2opusHandle_t handle = a->handle->handle;
1355f17f7ee2SSatish Balay #else
1356f17f7ee2SSatish Balay     h2opusHandle_t handle = a->handle;
1357f17f7ee2SSatish Balay #endif
1358f17f7ee2SSatish Balay     if (boundtocpu) {
1359f17f7ee2SSatish Balay       PetscCheck(a->hmatrix, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing CPU matrix");
1360f17f7ee2SSatish Balay       horthog(*a->hmatrix, handle);
1361f17f7ee2SSatish Balay #if defined(PETSC_H2OPUS_USE_GPU)
1362f17f7ee2SSatish Balay       A->offloadmask = PETSC_OFFLOAD_CPU;
1363f17f7ee2SSatish Balay     } else {
1364f17f7ee2SSatish Balay       PetscCheck(a->hmatrix_gpu, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing GPU matrix");
1365f17f7ee2SSatish Balay       PetscCall(PetscLogGpuTimeBegin());
1366f17f7ee2SSatish Balay       horthog(*a->hmatrix_gpu, handle);
1367f17f7ee2SSatish Balay       PetscCall(PetscLogGpuTimeEnd());
1368f17f7ee2SSatish Balay #endif
1369f17f7ee2SSatish Balay     }
1370f17f7ee2SSatish Balay   }
1371f17f7ee2SSatish Balay   a->orthogonal = PETSC_TRUE;
1372f17f7ee2SSatish Balay   { /* log flops */
1373f17f7ee2SSatish Balay     double gops, time, perf, dev;
1374f17f7ee2SSatish Balay     HLibProfile::getHorthogPerf(gops, time, perf, dev);
1375f17f7ee2SSatish Balay #if defined(PETSC_H2OPUS_USE_GPU)
1376f17f7ee2SSatish Balay     if (boundtocpu) {
1377f17f7ee2SSatish Balay       PetscCall(PetscLogFlops(1e9 * gops));
1378f17f7ee2SSatish Balay     } else {
1379f17f7ee2SSatish Balay       PetscCall(PetscLogGpuFlops(1e9 * gops));
1380f17f7ee2SSatish Balay     }
1381f17f7ee2SSatish Balay #else
1382f17f7ee2SSatish Balay     PetscCall(PetscLogFlops(1e9 * gops));
1383f17f7ee2SSatish Balay #endif
1384f17f7ee2SSatish Balay   }
1385f17f7ee2SSatish Balay   PetscCall(PetscLogEventEnd(MAT_H2Opus_Orthog, A, 0, 0, 0));
1386f17f7ee2SSatish Balay   PetscFunctionReturn(0);
1387f17f7ee2SSatish Balay }
1388f17f7ee2SSatish Balay 
1389f17f7ee2SSatish Balay /*@C
1390f17f7ee2SSatish Balay      MatH2OpusCompress - Compress a hierarchical matrix.
1391f17f7ee2SSatish Balay 
1392f17f7ee2SSatish Balay    Input Parameters:
1393f17f7ee2SSatish Balay +     A - the matrix
1394f17f7ee2SSatish Balay -     tol - the absolute truncation threshold
1395f17f7ee2SSatish Balay 
1396f17f7ee2SSatish Balay    Level: intermediate
1397f17f7ee2SSatish Balay 
1398036e4bdcSSatish Balay .seealso: `MatCreate()`, `MATH2OPUS`, `MatCreateH2OpusFromMat()`, `MatCreateH2OpusFromKernel()`, `MatH2OpusOrthogonalize()`
1399f17f7ee2SSatish Balay @*/
14009371c9d4SSatish Balay PetscErrorCode MatH2OpusCompress(Mat A, PetscReal tol) {
1401f17f7ee2SSatish Balay   PetscBool   ish2opus;
1402f17f7ee2SSatish Balay   Mat_H2OPUS *a = (Mat_H2OPUS *)A->data;
1403f17f7ee2SSatish Balay   PetscMPIInt size;
1404f17f7ee2SSatish Balay   PetscBool   boundtocpu = PETSC_TRUE;
1405f17f7ee2SSatish Balay 
1406f17f7ee2SSatish Balay   PetscFunctionBegin;
1407f17f7ee2SSatish Balay   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
1408f17f7ee2SSatish Balay   PetscValidType(A, 1);
1409f17f7ee2SSatish Balay   PetscValidLogicalCollectiveReal(A, tol, 2);
1410f17f7ee2SSatish Balay   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATH2OPUS, &ish2opus));
1411f17f7ee2SSatish Balay   if (!ish2opus || tol <= 0.0) PetscFunctionReturn(0);
1412f17f7ee2SSatish Balay   PetscCall(MatH2OpusOrthogonalize(A));
1413f17f7ee2SSatish Balay   HLibProfile::clear();
1414f17f7ee2SSatish Balay   PetscCall(PetscLogEventBegin(MAT_H2Opus_Compress, A, 0, 0, 0));
1415f17f7ee2SSatish Balay #if defined(PETSC_H2OPUS_USE_GPU)
1416f17f7ee2SSatish Balay   boundtocpu = A->boundtocpu;
1417f17f7ee2SSatish Balay #endif
1418f17f7ee2SSatish Balay   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
1419f17f7ee2SSatish Balay   if (size > 1) {
1420f17f7ee2SSatish Balay     if (boundtocpu) {
1421f17f7ee2SSatish Balay       PetscCheck(a->dist_hmatrix, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing CPU matrix");
1422f17f7ee2SSatish Balay #if defined(H2OPUS_USE_MPI)
1423f17f7ee2SSatish Balay       distributed_hcompress(*a->dist_hmatrix, tol, a->handle);
1424036e4bdcSSatish Balay       if (a->resize) {
1425036e4bdcSSatish Balay         DistributedHMatrix *dist_hmatrix = new DistributedHMatrix(*a->dist_hmatrix);
1426036e4bdcSSatish Balay         delete a->dist_hmatrix;
1427036e4bdcSSatish Balay         a->dist_hmatrix = dist_hmatrix;
1428036e4bdcSSatish Balay       }
1429f17f7ee2SSatish Balay #endif
1430f17f7ee2SSatish Balay #if defined(PETSC_H2OPUS_USE_GPU)
1431f17f7ee2SSatish Balay       A->offloadmask = PETSC_OFFLOAD_CPU;
1432f17f7ee2SSatish Balay     } else {
1433f17f7ee2SSatish Balay       PetscCheck(a->dist_hmatrix_gpu, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing GPU matrix");
1434f17f7ee2SSatish Balay       PetscCall(PetscLogGpuTimeBegin());
1435f17f7ee2SSatish Balay #if defined(H2OPUS_USE_MPI)
1436f17f7ee2SSatish Balay       distributed_hcompress(*a->dist_hmatrix_gpu, tol, a->handle);
1437036e4bdcSSatish Balay 
1438036e4bdcSSatish Balay       if (a->resize) {
1439036e4bdcSSatish Balay         DistributedHMatrix_GPU *dist_hmatrix_gpu = new DistributedHMatrix_GPU(*a->dist_hmatrix_gpu);
1440036e4bdcSSatish Balay         delete a->dist_hmatrix_gpu;
1441036e4bdcSSatish Balay         a->dist_hmatrix_gpu = dist_hmatrix_gpu;
1442036e4bdcSSatish Balay       }
1443f17f7ee2SSatish Balay #endif
1444f17f7ee2SSatish Balay       PetscCall(PetscLogGpuTimeEnd());
1445f17f7ee2SSatish Balay #endif
1446f17f7ee2SSatish Balay     }
1447f17f7ee2SSatish Balay   } else {
1448f17f7ee2SSatish Balay #if defined(H2OPUS_USE_MPI)
1449f17f7ee2SSatish Balay     h2opusHandle_t handle = a->handle->handle;
1450f17f7ee2SSatish Balay #else
1451f17f7ee2SSatish Balay     h2opusHandle_t handle = a->handle;
1452f17f7ee2SSatish Balay #endif
1453f17f7ee2SSatish Balay     if (boundtocpu) {
1454f17f7ee2SSatish Balay       PetscCheck(a->hmatrix, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing CPU matrix");
1455f17f7ee2SSatish Balay       hcompress(*a->hmatrix, tol, handle);
1456036e4bdcSSatish Balay 
1457036e4bdcSSatish Balay       if (a->resize) {
1458036e4bdcSSatish Balay         HMatrix *hmatrix = new HMatrix(*a->hmatrix);
1459036e4bdcSSatish Balay         delete a->hmatrix;
1460036e4bdcSSatish Balay         a->hmatrix = hmatrix;
1461036e4bdcSSatish Balay       }
1462f17f7ee2SSatish Balay #if defined(PETSC_H2OPUS_USE_GPU)
1463f17f7ee2SSatish Balay       A->offloadmask = PETSC_OFFLOAD_CPU;
1464f17f7ee2SSatish Balay     } else {
1465f17f7ee2SSatish Balay       PetscCheck(a->hmatrix_gpu, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing GPU matrix");
1466f17f7ee2SSatish Balay       PetscCall(PetscLogGpuTimeBegin());
1467f17f7ee2SSatish Balay       hcompress(*a->hmatrix_gpu, tol, handle);
1468f17f7ee2SSatish Balay       PetscCall(PetscLogGpuTimeEnd());
1469036e4bdcSSatish Balay 
1470036e4bdcSSatish Balay       if (a->resize) {
1471036e4bdcSSatish Balay         HMatrix_GPU *hmatrix_gpu = new HMatrix_GPU(*a->hmatrix_gpu);
1472036e4bdcSSatish Balay         delete a->hmatrix_gpu;
1473036e4bdcSSatish Balay         a->hmatrix_gpu = hmatrix_gpu;
1474036e4bdcSSatish Balay       }
1475f17f7ee2SSatish Balay #endif
1476f17f7ee2SSatish Balay     }
1477f17f7ee2SSatish Balay   }
1478f17f7ee2SSatish Balay   { /* log flops */
1479f17f7ee2SSatish Balay     double gops, time, perf, dev;
1480f17f7ee2SSatish Balay     HLibProfile::getHcompressPerf(gops, time, perf, dev);
1481f17f7ee2SSatish Balay #if defined(PETSC_H2OPUS_USE_GPU)
1482f17f7ee2SSatish Balay     if (boundtocpu) {
1483f17f7ee2SSatish Balay       PetscCall(PetscLogFlops(1e9 * gops));
1484f17f7ee2SSatish Balay     } else {
1485f17f7ee2SSatish Balay       PetscCall(PetscLogGpuFlops(1e9 * gops));
1486f17f7ee2SSatish Balay     }
1487f17f7ee2SSatish Balay #else
1488f17f7ee2SSatish Balay     PetscCall(PetscLogFlops(1e9 * gops));
1489f17f7ee2SSatish Balay #endif
1490f17f7ee2SSatish Balay   }
1491f17f7ee2SSatish Balay   PetscCall(PetscLogEventEnd(MAT_H2Opus_Compress, A, 0, 0, 0));
1492f17f7ee2SSatish Balay   PetscFunctionReturn(0);
1493f17f7ee2SSatish Balay }
1494f17f7ee2SSatish Balay 
1495f17f7ee2SSatish Balay /*@C
1496f17f7ee2SSatish Balay      MatH2OpusSetSamplingMat - Set a matrix to be sampled from matrix vector product to construct a hierarchical matrix.
1497f17f7ee2SSatish Balay 
1498f17f7ee2SSatish Balay    Input Parameters:
1499f17f7ee2SSatish Balay +     A - the hierarchical matrix
1500f17f7ee2SSatish Balay .     B - the matrix to be sampled
1501f17f7ee2SSatish Balay .     bs - maximum number of samples to be taken concurrently
1502f17f7ee2SSatish Balay -     tol - relative tolerance for construction
1503f17f7ee2SSatish Balay 
150411a5261eSBarry Smith    Notes:
150511a5261eSBarry Smith    You need to call `MatAssemblyBegin()` and `MatAssemblyEnd()` to update the hierarchical matrix.
1506f17f7ee2SSatish Balay 
1507f17f7ee2SSatish Balay    Level: intermediate
1508f17f7ee2SSatish Balay 
1509036e4bdcSSatish Balay .seealso: `MatCreate()`, `MATH2OPUS`, `MatCreateH2OpusFromMat()`, `MatCreateH2OpusFromKernel()`, `MatH2OpusCompress()`, `MatH2OpusOrthogonalize()`
1510f17f7ee2SSatish Balay @*/
15119371c9d4SSatish Balay PetscErrorCode MatH2OpusSetSamplingMat(Mat A, Mat B, PetscInt bs, PetscReal tol) {
1512f17f7ee2SSatish Balay   PetscBool ish2opus;
1513f17f7ee2SSatish Balay 
1514f17f7ee2SSatish Balay   PetscFunctionBegin;
1515f17f7ee2SSatish Balay   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
1516f17f7ee2SSatish Balay   PetscValidType(A, 1);
1517f17f7ee2SSatish Balay   if (B) PetscValidHeaderSpecific(B, MAT_CLASSID, 2);
1518f17f7ee2SSatish Balay   PetscValidLogicalCollectiveInt(A, bs, 3);
1519f17f7ee2SSatish Balay   PetscValidLogicalCollectiveReal(A, tol, 3);
1520f17f7ee2SSatish Balay   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATH2OPUS, &ish2opus));
1521f17f7ee2SSatish Balay   if (ish2opus) {
1522f17f7ee2SSatish Balay     Mat_H2OPUS *a = (Mat_H2OPUS *)A->data;
1523f17f7ee2SSatish Balay 
1524f17f7ee2SSatish Balay     if (!a->sampler) a->sampler = new PetscMatrixSampler();
1525f17f7ee2SSatish Balay     a->sampler->SetSamplingMat(B);
1526f17f7ee2SSatish Balay     if (bs > 0) a->bs = bs;
1527f17f7ee2SSatish Balay     if (tol > 0.) a->rtol = tol;
1528f17f7ee2SSatish Balay     delete a->kernel;
1529f17f7ee2SSatish Balay   }
1530f17f7ee2SSatish Balay   PetscFunctionReturn(0);
1531f17f7ee2SSatish Balay }
1532f17f7ee2SSatish Balay 
1533f17f7ee2SSatish Balay /*@C
153411a5261eSBarry Smith      MatCreateH2OpusFromKernel - Creates a `MATH2OPUS` from a user-supplied kernel.
1535f17f7ee2SSatish Balay 
1536f17f7ee2SSatish Balay    Input Parameters:
1537f17f7ee2SSatish Balay +     comm - MPI communicator
153811a5261eSBarry Smith .     m - number of local rows (or `PETSC_DECIDE` to have calculated if M is given)
153911a5261eSBarry Smith .     n - number of local columns (or `PETSC_DECIDE` to have calculated if N is given)
154011a5261eSBarry Smith .     M - number of global rows (or `PETSC_DETERMINE` to have calculated if m is given)
154111a5261eSBarry Smith .     N - number of global columns (or `PETSC_DETERMINE` to have calculated if n is given)
1542f17f7ee2SSatish Balay .     spacedim - dimension of the space coordinates
1543f17f7ee2SSatish Balay .     coords - coordinates of the points
1544f17f7ee2SSatish Balay .     cdist - whether or not coordinates are distributed
1545f17f7ee2SSatish Balay .     kernel - computational kernel (or NULL)
1546f17f7ee2SSatish Balay .     kernelctx - kernel context
1547f17f7ee2SSatish Balay .     eta - admissibility condition tolerance
1548f17f7ee2SSatish Balay .     leafsize - leaf size in cluster tree
1549f17f7ee2SSatish Balay -     basisord - approximation order for Chebychev interpolation of low-rank blocks
1550f17f7ee2SSatish Balay 
1551f17f7ee2SSatish Balay    Output Parameter:
1552f17f7ee2SSatish Balay .     nA - matrix
1553f17f7ee2SSatish Balay 
1554f17f7ee2SSatish Balay    Options Database Keys:
155511a5261eSBarry Smith +     -mat_h2opus_leafsize <`PetscInt`> - Leaf size of cluster tree
155611a5261eSBarry Smith .     -mat_h2opus_eta <`PetscReal`> - Admissibility condition tolerance
155711a5261eSBarry Smith .     -mat_h2opus_order <`PetscInt`> - Chebychev approximation order
155811a5261eSBarry Smith -     -mat_h2opus_normsamples <`PetscInt`> - Maximum number of samples to be used when estimating norms
1559f17f7ee2SSatish Balay 
1560f17f7ee2SSatish Balay    Level: intermediate
1561f17f7ee2SSatish Balay 
1562036e4bdcSSatish Balay .seealso: `MatCreate()`, `MATH2OPUS`, `MatCreateH2OpusFromMat()`
1563f17f7ee2SSatish Balay @*/
15649371c9d4SSatish Balay PetscErrorCode MatCreateH2OpusFromKernel(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt M, PetscInt N, PetscInt spacedim, const PetscReal coords[], PetscBool cdist, MatH2OpusKernel kernel, void *kernelctx, PetscReal eta, PetscInt leafsize, PetscInt basisord, Mat *nA) {
1565f17f7ee2SSatish Balay   Mat         A;
1566f17f7ee2SSatish Balay   Mat_H2OPUS *h2opus;
1567f17f7ee2SSatish Balay #if defined(PETSC_H2OPUS_USE_GPU)
1568f17f7ee2SSatish Balay   PetscBool iscpu = PETSC_FALSE;
1569f17f7ee2SSatish Balay #else
1570f17f7ee2SSatish Balay   PetscBool iscpu = PETSC_TRUE;
1571f17f7ee2SSatish Balay #endif
1572f17f7ee2SSatish Balay 
1573f17f7ee2SSatish Balay   PetscFunctionBegin;
1574036e4bdcSSatish Balay   PetscCheck(m == n, PETSC_COMM_SELF, PETSC_ERR_SUP, "Different row and column local sizes are not supported");
1575f17f7ee2SSatish Balay   PetscCall(MatCreate(comm, &A));
1576f17f7ee2SSatish Balay   PetscCall(MatSetSizes(A, m, n, M, N));
1577036e4bdcSSatish Balay   PetscCheck(M == N, comm, PETSC_ERR_SUP, "Rectangular matrices are not supported");
1578f17f7ee2SSatish Balay   PetscCall(MatSetType(A, MATH2OPUS));
1579f17f7ee2SSatish Balay   PetscCall(MatBindToCPU(A, iscpu));
1580f17f7ee2SSatish Balay   PetscCall(MatH2OpusSetCoords_H2OPUS(A, spacedim, coords, cdist, kernel, kernelctx));
1581f17f7ee2SSatish Balay 
1582f17f7ee2SSatish Balay   h2opus = (Mat_H2OPUS *)A->data;
1583f17f7ee2SSatish Balay   if (eta > 0.) h2opus->eta = eta;
1584f17f7ee2SSatish Balay   if (leafsize > 0) h2opus->leafsize = leafsize;
1585f17f7ee2SSatish Balay   if (basisord > 0) h2opus->basisord = basisord;
1586f17f7ee2SSatish Balay 
1587f17f7ee2SSatish Balay   *nA = A;
1588f17f7ee2SSatish Balay   PetscFunctionReturn(0);
1589f17f7ee2SSatish Balay }
1590f17f7ee2SSatish Balay 
1591f17f7ee2SSatish Balay /*@C
159211a5261eSBarry Smith      MatCreateH2OpusFromMat - Creates a `MATH2OPUS` sampling from a user-supplied operator.
1593f17f7ee2SSatish Balay 
1594f17f7ee2SSatish Balay    Input Parameters:
1595f17f7ee2SSatish Balay +     B - the matrix to be sampled
1596f17f7ee2SSatish Balay .     spacedim - dimension of the space coordinates
1597f17f7ee2SSatish Balay .     coords - coordinates of the points
1598f17f7ee2SSatish Balay .     cdist - whether or not coordinates are distributed
1599f17f7ee2SSatish Balay .     eta - admissibility condition tolerance
1600f17f7ee2SSatish Balay .     leafsize - leaf size in cluster tree
1601f17f7ee2SSatish Balay .     maxrank - maximum rank allowed
1602f17f7ee2SSatish Balay .     bs - maximum number of samples to be taken concurrently
1603f17f7ee2SSatish Balay -     rtol - relative tolerance for construction
1604f17f7ee2SSatish Balay 
1605f17f7ee2SSatish Balay    Output Parameter:
1606f17f7ee2SSatish Balay .     nA - matrix
1607f17f7ee2SSatish Balay 
1608f17f7ee2SSatish Balay    Options Database Keys:
160911a5261eSBarry Smith +     -mat_h2opus_leafsize <`PetscInt`> - Leaf size of cluster tree
161011a5261eSBarry Smith .     -mat_h2opus_eta <`PetscReal`> - Admissibility condition tolerance
161111a5261eSBarry Smith .     -mat_h2opus_maxrank <`PetscInt`> - Maximum rank when constructed from matvecs
161211a5261eSBarry Smith .     -mat_h2opus_samples <`PetscInt`> - Maximum number of samples to be taken concurrently when constructing from matvecs
161311a5261eSBarry Smith .     -mat_h2opus_rtol <`PetscReal`> - Relative tolerance for construction from sampling
161411a5261eSBarry Smith .     -mat_h2opus_check <`PetscBool`> - Check error when constructing from sampling during MatAssemblyEnd()
161511a5261eSBarry Smith .     -mat_h2opus_hara_verbose <`PetscBool`> - Verbose output from hara construction
161611a5261eSBarry Smith -     -mat_h2opus_normsamples <`PetscInt`> - Maximum bumber of samples to be when estimating norms
1617f17f7ee2SSatish Balay 
161811a5261eSBarry Smith    Note:
161911a5261eSBarry Smith    Not available in parallel
1620f17f7ee2SSatish Balay 
1621f17f7ee2SSatish Balay    Level: intermediate
1622f17f7ee2SSatish Balay 
1623036e4bdcSSatish Balay .seealso: `MatCreate()`, `MATH2OPUS`, `MatCreateH2OpusFromKernel()`
1624f17f7ee2SSatish Balay @*/
16259371c9d4SSatish Balay PetscErrorCode MatCreateH2OpusFromMat(Mat B, PetscInt spacedim, const PetscReal coords[], PetscBool cdist, PetscReal eta, PetscInt leafsize, PetscInt maxrank, PetscInt bs, PetscReal rtol, Mat *nA) {
1626f17f7ee2SSatish Balay   Mat         A;
1627f17f7ee2SSatish Balay   Mat_H2OPUS *h2opus;
1628f17f7ee2SSatish Balay   MPI_Comm    comm;
1629f17f7ee2SSatish Balay   PetscBool   boundtocpu = PETSC_TRUE;
1630f17f7ee2SSatish Balay 
1631f17f7ee2SSatish Balay   PetscFunctionBegin;
1632f17f7ee2SSatish Balay   PetscValidHeaderSpecific(B, MAT_CLASSID, 1);
1633f17f7ee2SSatish Balay   PetscValidLogicalCollectiveInt(B, spacedim, 2);
1634f17f7ee2SSatish Balay   PetscValidLogicalCollectiveBool(B, cdist, 4);
1635f17f7ee2SSatish Balay   PetscValidLogicalCollectiveReal(B, eta, 5);
1636f17f7ee2SSatish Balay   PetscValidLogicalCollectiveInt(B, leafsize, 6);
1637f17f7ee2SSatish Balay   PetscValidLogicalCollectiveInt(B, maxrank, 7);
1638f17f7ee2SSatish Balay   PetscValidLogicalCollectiveInt(B, bs, 8);
1639f17f7ee2SSatish Balay   PetscValidLogicalCollectiveReal(B, rtol, 9);
1640f17f7ee2SSatish Balay   PetscValidPointer(nA, 10);
1641f17f7ee2SSatish Balay   PetscCall(PetscObjectGetComm((PetscObject)B, &comm));
1642036e4bdcSSatish Balay   PetscCheck(B->rmap->n == B->cmap->n, PETSC_COMM_SELF, PETSC_ERR_SUP, "Different row and column local sizes are not supported");
1643036e4bdcSSatish Balay   PetscCheck(B->rmap->N == B->cmap->N, comm, PETSC_ERR_SUP, "Rectangular matrices are not supported");
1644f17f7ee2SSatish Balay   PetscCall(MatCreate(comm, &A));
1645f17f7ee2SSatish Balay   PetscCall(MatSetSizes(A, B->rmap->n, B->cmap->n, B->rmap->N, B->cmap->N));
1646f17f7ee2SSatish Balay #if defined(PETSC_H2OPUS_USE_GPU)
1647f17f7ee2SSatish Balay   {
1648f17f7ee2SSatish Balay     PetscBool iscuda;
1649f17f7ee2SSatish Balay     VecType   vtype;
1650f17f7ee2SSatish Balay 
1651f17f7ee2SSatish Balay     PetscCall(MatGetVecType(B, &vtype));
1652f17f7ee2SSatish Balay     PetscCall(PetscStrcmp(vtype, VECCUDA, &iscuda));
1653f17f7ee2SSatish Balay     if (!iscuda) {
1654f17f7ee2SSatish Balay       PetscCall(PetscStrcmp(vtype, VECSEQCUDA, &iscuda));
165548a46eb9SPierre Jolivet       if (!iscuda) PetscCall(PetscStrcmp(vtype, VECMPICUDA, &iscuda));
1656f17f7ee2SSatish Balay     }
1657f17f7ee2SSatish Balay     if (iscuda && !B->boundtocpu) boundtocpu = PETSC_FALSE;
1658f17f7ee2SSatish Balay   }
1659f17f7ee2SSatish Balay #endif
1660f17f7ee2SSatish Balay   PetscCall(MatSetType(A, MATH2OPUS));
1661f17f7ee2SSatish Balay   PetscCall(MatBindToCPU(A, boundtocpu));
166248a46eb9SPierre Jolivet   if (spacedim) PetscCall(MatH2OpusSetCoords_H2OPUS(A, spacedim, coords, cdist, NULL, NULL));
1663f17f7ee2SSatish Balay   PetscCall(MatPropagateSymmetryOptions(B, A));
1664f17f7ee2SSatish Balay   /* PetscCheck(A->symmetric,comm,PETSC_ERR_SUP,"Unsymmetric sampling does not work"); */
1665f17f7ee2SSatish Balay 
1666f17f7ee2SSatish Balay   h2opus          = (Mat_H2OPUS *)A->data;
1667f17f7ee2SSatish Balay   h2opus->sampler = new PetscMatrixSampler(B);
1668f17f7ee2SSatish Balay   if (eta > 0.) h2opus->eta = eta;
1669f17f7ee2SSatish Balay   if (leafsize > 0) h2opus->leafsize = leafsize;
1670f17f7ee2SSatish Balay   if (maxrank > 0) h2opus->max_rank = maxrank;
1671f17f7ee2SSatish Balay   if (bs > 0) h2opus->bs = bs;
1672f17f7ee2SSatish Balay   if (rtol > 0.) h2opus->rtol = rtol;
1673f17f7ee2SSatish Balay   *nA             = A;
1674f17f7ee2SSatish Balay   A->preallocated = PETSC_TRUE;
1675f17f7ee2SSatish Balay   PetscFunctionReturn(0);
1676f17f7ee2SSatish Balay }
1677f17f7ee2SSatish Balay 
1678f17f7ee2SSatish Balay /*@C
1679f17f7ee2SSatish Balay      MatH2OpusGetIndexMap - Access reordering index set.
1680f17f7ee2SSatish Balay 
1681f17f7ee2SSatish Balay    Input Parameters:
1682f17f7ee2SSatish Balay .     A - the matrix
1683f17f7ee2SSatish Balay 
1684f17f7ee2SSatish Balay    Output Parameter:
1685f17f7ee2SSatish Balay .     indexmap - the index set for the reordering
1686f17f7ee2SSatish Balay 
1687f17f7ee2SSatish Balay    Level: intermediate
1688f17f7ee2SSatish Balay 
1689036e4bdcSSatish Balay .seealso: `MatCreate()`, `MATH2OPUS`, `MatCreateH2OpusFromMat()`, `MatCreateH2OpusFromKernel()`
1690f17f7ee2SSatish Balay @*/
16919371c9d4SSatish Balay PetscErrorCode MatH2OpusGetIndexMap(Mat A, IS *indexmap) {
1692f17f7ee2SSatish Balay   PetscBool   ish2opus;
1693f17f7ee2SSatish Balay   Mat_H2OPUS *a = (Mat_H2OPUS *)A->data;
1694f17f7ee2SSatish Balay 
1695f17f7ee2SSatish Balay   PetscFunctionBegin;
1696f17f7ee2SSatish Balay   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
1697f17f7ee2SSatish Balay   PetscValidType(A, 1);
1698f17f7ee2SSatish Balay   PetscValidPointer(indexmap, 2);
1699f17f7ee2SSatish Balay   PetscCheck(A->assembled, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
1700f17f7ee2SSatish Balay   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATH2OPUS, &ish2opus));
1701f17f7ee2SSatish Balay   PetscCheck(ish2opus, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Not for type %s", ((PetscObject)A)->type_name);
1702f17f7ee2SSatish Balay   *indexmap = a->h2opus_indexmap;
1703f17f7ee2SSatish Balay   PetscFunctionReturn(0);
1704f17f7ee2SSatish Balay }
1705f17f7ee2SSatish Balay 
1706f17f7ee2SSatish Balay /*@C
1707f17f7ee2SSatish Balay      MatH2OpusMapVec - Maps a vector between PETSc and H2Opus ordering
1708f17f7ee2SSatish Balay 
1709f17f7ee2SSatish Balay    Input Parameters:
1710f17f7ee2SSatish Balay +     A - the matrix
1711f17f7ee2SSatish Balay .     nativetopetsc - if true, maps from H2Opus ordering to PETSc ordering. If false, applies the reverse map
1712f17f7ee2SSatish Balay -     in - the vector to be mapped
1713f17f7ee2SSatish Balay 
1714f17f7ee2SSatish Balay    Output Parameter:
1715f17f7ee2SSatish Balay .     out - the newly created mapped vector
1716f17f7ee2SSatish Balay 
1717f17f7ee2SSatish Balay    Level: intermediate
1718f17f7ee2SSatish Balay 
1719036e4bdcSSatish Balay .seealso: `MatCreate()`, `MATH2OPUS`, `MatCreateH2OpusFromMat()`, `MatCreateH2OpusFromKernel()`
1720f17f7ee2SSatish Balay @*/
17219371c9d4SSatish Balay PetscErrorCode MatH2OpusMapVec(Mat A, PetscBool nativetopetsc, Vec in, Vec *out) {
1722f17f7ee2SSatish Balay   PetscBool    ish2opus;
1723f17f7ee2SSatish Balay   Mat_H2OPUS  *a = (Mat_H2OPUS *)A->data;
1724f17f7ee2SSatish Balay   PetscScalar *xin, *xout;
1725f17f7ee2SSatish Balay   PetscBool    nm;
1726f17f7ee2SSatish Balay 
1727f17f7ee2SSatish Balay   PetscFunctionBegin;
1728f17f7ee2SSatish Balay   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
1729f17f7ee2SSatish Balay   PetscValidType(A, 1);
1730f17f7ee2SSatish Balay   PetscValidLogicalCollectiveBool(A, nativetopetsc, 2);
1731f17f7ee2SSatish Balay   PetscValidHeaderSpecific(in, VEC_CLASSID, 3);
1732f17f7ee2SSatish Balay   PetscValidPointer(out, 4);
1733f17f7ee2SSatish Balay   PetscCheck(A->assembled, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
1734f17f7ee2SSatish Balay   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATH2OPUS, &ish2opus));
1735f17f7ee2SSatish Balay   PetscCheck(ish2opus, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Not for type %s", ((PetscObject)A)->type_name);
1736f17f7ee2SSatish Balay   nm = a->nativemult;
1737f17f7ee2SSatish Balay   PetscCall(MatH2OpusSetNativeMult(A, (PetscBool)!nativetopetsc));
1738f17f7ee2SSatish Balay   PetscCall(MatCreateVecs(A, out, NULL));
1739f17f7ee2SSatish Balay   PetscCall(MatH2OpusSetNativeMult(A, nm));
1740f17f7ee2SSatish Balay   if (!a->sf) { /* same ordering */
1741f17f7ee2SSatish Balay     PetscCall(VecCopy(in, *out));
1742f17f7ee2SSatish Balay     PetscFunctionReturn(0);
1743f17f7ee2SSatish Balay   }
1744f17f7ee2SSatish Balay   PetscCall(VecGetArrayRead(in, (const PetscScalar **)&xin));
1745f17f7ee2SSatish Balay   PetscCall(VecGetArrayWrite(*out, &xout));
1746f17f7ee2SSatish Balay   if (nativetopetsc) {
1747f17f7ee2SSatish Balay     PetscCall(PetscSFReduceBegin(a->sf, MPIU_SCALAR, xin, xout, MPI_REPLACE));
1748f17f7ee2SSatish Balay     PetscCall(PetscSFReduceEnd(a->sf, MPIU_SCALAR, xin, xout, MPI_REPLACE));
1749f17f7ee2SSatish Balay   } else {
1750f17f7ee2SSatish Balay     PetscCall(PetscSFBcastBegin(a->sf, MPIU_SCALAR, xin, xout, MPI_REPLACE));
1751f17f7ee2SSatish Balay     PetscCall(PetscSFBcastEnd(a->sf, MPIU_SCALAR, xin, xout, MPI_REPLACE));
1752f17f7ee2SSatish Balay   }
1753f17f7ee2SSatish Balay   PetscCall(VecRestoreArrayRead(in, (const PetscScalar **)&xin));
1754f17f7ee2SSatish Balay   PetscCall(VecRestoreArrayWrite(*out, &xout));
1755f17f7ee2SSatish Balay   PetscFunctionReturn(0);
1756f17f7ee2SSatish Balay }
1757f17f7ee2SSatish Balay 
1758f17f7ee2SSatish Balay /*@C
1759f17f7ee2SSatish Balay      MatH2OpusLowRankUpdate - Perform a low-rank update of the form A = A + s * U * V^T
1760f17f7ee2SSatish Balay 
1761f17f7ee2SSatish Balay    Input Parameters:
176211a5261eSBarry Smith +     A - the hierarchical `MATH2OPUS` matrix
1763f17f7ee2SSatish Balay .     s - the scaling factor
1764f17f7ee2SSatish Balay .     U - the dense low-rank update matrix
1765f17f7ee2SSatish Balay -     V - (optional) the dense low-rank update matrix (if NULL, then V = U is assumed)
1766f17f7ee2SSatish Balay 
176711a5261eSBarry Smith    Note:
176811a5261eSBarry Smith    The U and V matrices must be in dense format
1769f17f7ee2SSatish Balay 
1770f17f7ee2SSatish Balay    Level: intermediate
1771f17f7ee2SSatish Balay 
1772036e4bdcSSatish Balay .seealso: `MatCreate()`, `MATH2OPUS`, `MatCreateH2OpusFromMat()`, `MatCreateH2OpusFromKernel()`, `MatH2OpusCompress()`, `MatH2OpusOrthogonalize()`, `MATDENSE`
1773f17f7ee2SSatish Balay @*/
17749371c9d4SSatish Balay PetscErrorCode MatH2OpusLowRankUpdate(Mat A, Mat U, Mat V, PetscScalar s) {
1775f17f7ee2SSatish Balay   PetscBool flg;
1776f17f7ee2SSatish Balay 
1777f17f7ee2SSatish Balay   PetscFunctionBegin;
1778f17f7ee2SSatish Balay   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
1779f17f7ee2SSatish Balay   PetscValidType(A, 1);
1780f17f7ee2SSatish Balay   PetscCheck(A->assembled, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
1781f17f7ee2SSatish Balay   PetscValidHeaderSpecific(U, MAT_CLASSID, 2);
1782f17f7ee2SSatish Balay   PetscCheckSameComm(A, 1, U, 2);
1783f17f7ee2SSatish Balay   if (V) {
1784f17f7ee2SSatish Balay     PetscValidHeaderSpecific(V, MAT_CLASSID, 3);
1785f17f7ee2SSatish Balay     PetscCheckSameComm(A, 1, V, 3);
1786f17f7ee2SSatish Balay   }
1787f17f7ee2SSatish Balay   PetscValidLogicalCollectiveScalar(A, s, 4);
1788f17f7ee2SSatish Balay 
1789f17f7ee2SSatish Balay   if (!V) V = U;
1790036e4bdcSSatish Balay   PetscCheck(U->cmap->N == V->cmap->N, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Non matching rank update %" PetscInt_FMT " != %" PetscInt_FMT, U->cmap->N, V->cmap->N);
1791f17f7ee2SSatish Balay   if (!U->cmap->N) PetscFunctionReturn(0);
1792f17f7ee2SSatish Balay   PetscCall(PetscLayoutCompare(U->rmap, A->rmap, &flg));
1793f17f7ee2SSatish Balay   PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "A and U must have the same row layout");
1794f17f7ee2SSatish Balay   PetscCall(PetscLayoutCompare(V->rmap, A->cmap, &flg));
1795f17f7ee2SSatish Balay   PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "A column layout must match V row column layout");
1796f17f7ee2SSatish Balay   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATH2OPUS, &flg));
1797f17f7ee2SSatish Balay   if (flg) {
1798f17f7ee2SSatish Balay     Mat_H2OPUS        *a = (Mat_H2OPUS *)A->data;
1799f17f7ee2SSatish Balay     const PetscScalar *u, *v, *uu, *vv;
1800f17f7ee2SSatish Balay     PetscInt           ldu, ldv;
1801f17f7ee2SSatish Balay     PetscMPIInt        size;
1802f17f7ee2SSatish Balay #if defined(H2OPUS_USE_MPI)
1803f17f7ee2SSatish Balay     h2opusHandle_t handle = a->handle->handle;
1804f17f7ee2SSatish Balay #else
1805f17f7ee2SSatish Balay     h2opusHandle_t handle = a->handle;
1806f17f7ee2SSatish Balay #endif
1807f17f7ee2SSatish Balay     PetscBool usesf = (PetscBool)(a->sf && !a->nativemult);
1808f17f7ee2SSatish Balay     PetscSF   usf, vsf;
1809f17f7ee2SSatish Balay 
1810f17f7ee2SSatish Balay     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
1811036e4bdcSSatish Balay     PetscCheck(size <= 1, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Not yet implemented in parallel");
1812f17f7ee2SSatish Balay     PetscCall(PetscLogEventBegin(MAT_H2Opus_LR, A, 0, 0, 0));
1813f17f7ee2SSatish Balay     PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)U, &flg, MATSEQDENSE, MATMPIDENSE, ""));
1814f17f7ee2SSatish Balay     PetscCheck(flg, PetscObjectComm((PetscObject)U), PETSC_ERR_SUP, "Not for U of type %s", ((PetscObject)U)->type_name);
1815f17f7ee2SSatish Balay     PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)V, &flg, MATSEQDENSE, MATMPIDENSE, ""));
1816f17f7ee2SSatish Balay     PetscCheck(flg, PetscObjectComm((PetscObject)V), PETSC_ERR_SUP, "Not for V of type %s", ((PetscObject)V)->type_name);
1817f17f7ee2SSatish Balay     PetscCall(MatDenseGetLDA(U, &ldu));
1818f17f7ee2SSatish Balay     PetscCall(MatDenseGetLDA(V, &ldv));
1819f17f7ee2SSatish Balay     PetscCall(MatBoundToCPU(A, &flg));
1820f17f7ee2SSatish Balay     if (usesf) {
1821f17f7ee2SSatish Balay       PetscInt n;
1822f17f7ee2SSatish Balay 
1823f17f7ee2SSatish Balay       PetscCall(MatDenseGetH2OpusVectorSF(U, a->sf, &usf));
1824f17f7ee2SSatish Balay       PetscCall(MatDenseGetH2OpusVectorSF(V, a->sf, &vsf));
1825f17f7ee2SSatish Balay       PetscCall(MatH2OpusResizeBuffers_Private(A, U->cmap->N, V->cmap->N));
1826f17f7ee2SSatish Balay       PetscCall(PetscSFGetGraph(a->sf, NULL, &n, NULL, NULL));
1827f17f7ee2SSatish Balay       ldu = n;
1828f17f7ee2SSatish Balay       ldv = n;
1829f17f7ee2SSatish Balay     }
1830f17f7ee2SSatish Balay     if (flg) {
1831f17f7ee2SSatish Balay       PetscCheck(a->hmatrix, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing CPU matrix");
1832f17f7ee2SSatish Balay       PetscCall(MatDenseGetArrayRead(U, &u));
1833f17f7ee2SSatish Balay       PetscCall(MatDenseGetArrayRead(V, &v));
1834f17f7ee2SSatish Balay       if (usesf) {
1835f17f7ee2SSatish Balay         vv = MatH2OpusGetThrustPointer(*a->yy);
1836f17f7ee2SSatish Balay         PetscCall(PetscSFBcastBegin(vsf, MPIU_SCALAR, v, (PetscScalar *)vv, MPI_REPLACE));
1837f17f7ee2SSatish Balay         PetscCall(PetscSFBcastEnd(vsf, MPIU_SCALAR, v, (PetscScalar *)vv, MPI_REPLACE));
1838f17f7ee2SSatish Balay         if (U != V) {
1839f17f7ee2SSatish Balay           uu = MatH2OpusGetThrustPointer(*a->xx);
1840f17f7ee2SSatish Balay           PetscCall(PetscSFBcastBegin(usf, MPIU_SCALAR, u, (PetscScalar *)uu, MPI_REPLACE));
1841f17f7ee2SSatish Balay           PetscCall(PetscSFBcastEnd(usf, MPIU_SCALAR, u, (PetscScalar *)uu, MPI_REPLACE));
1842f17f7ee2SSatish Balay         } else uu = vv;
18439371c9d4SSatish Balay       } else {
18449371c9d4SSatish Balay         uu = u;
18459371c9d4SSatish Balay         vv = v;
18469371c9d4SSatish Balay       }
1847f17f7ee2SSatish Balay       hlru_global(*a->hmatrix, uu, ldu, vv, ldv, U->cmap->N, s, handle);
1848f17f7ee2SSatish Balay       PetscCall(MatDenseRestoreArrayRead(U, &u));
1849f17f7ee2SSatish Balay       PetscCall(MatDenseRestoreArrayRead(V, &v));
1850f17f7ee2SSatish Balay     } else {
1851f17f7ee2SSatish Balay #if defined(PETSC_H2OPUS_USE_GPU)
1852f17f7ee2SSatish Balay       PetscBool flgU, flgV;
1853f17f7ee2SSatish Balay 
1854f17f7ee2SSatish Balay       PetscCheck(a->hmatrix_gpu, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing GPU matrix");
1855f17f7ee2SSatish Balay       PetscCall(PetscObjectTypeCompareAny((PetscObject)U, &flgU, MATSEQDENSE, MATMPIDENSE, ""));
1856f17f7ee2SSatish Balay       if (flgU) PetscCall(MatConvert(U, MATDENSECUDA, MAT_INPLACE_MATRIX, &U));
1857f17f7ee2SSatish Balay       PetscCall(PetscObjectTypeCompareAny((PetscObject)V, &flgV, MATSEQDENSE, MATMPIDENSE, ""));
1858f17f7ee2SSatish Balay       if (flgV) PetscCall(MatConvert(V, MATDENSECUDA, MAT_INPLACE_MATRIX, &V));
1859f17f7ee2SSatish Balay       PetscCall(MatDenseCUDAGetArrayRead(U, &u));
1860f17f7ee2SSatish Balay       PetscCall(MatDenseCUDAGetArrayRead(V, &v));
1861f17f7ee2SSatish Balay       if (usesf) {
1862f17f7ee2SSatish Balay         vv = MatH2OpusGetThrustPointer(*a->yy_gpu);
1863f17f7ee2SSatish Balay         PetscCall(PetscSFBcastBegin(vsf, MPIU_SCALAR, v, (PetscScalar *)vv, MPI_REPLACE));
1864f17f7ee2SSatish Balay         PetscCall(PetscSFBcastEnd(vsf, MPIU_SCALAR, v, (PetscScalar *)vv, MPI_REPLACE));
1865f17f7ee2SSatish Balay         if (U != V) {
1866f17f7ee2SSatish Balay           uu = MatH2OpusGetThrustPointer(*a->xx_gpu);
1867f17f7ee2SSatish Balay           PetscCall(PetscSFBcastBegin(usf, MPIU_SCALAR, u, (PetscScalar *)uu, MPI_REPLACE));
1868f17f7ee2SSatish Balay           PetscCall(PetscSFBcastEnd(usf, MPIU_SCALAR, u, (PetscScalar *)uu, MPI_REPLACE));
1869f17f7ee2SSatish Balay         } else uu = vv;
18709371c9d4SSatish Balay       } else {
18719371c9d4SSatish Balay         uu = u;
18729371c9d4SSatish Balay         vv = v;
18739371c9d4SSatish Balay       }
1874f17f7ee2SSatish Balay #else
1875f17f7ee2SSatish Balay       SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "This should not happen");
1876f17f7ee2SSatish Balay #endif
1877f17f7ee2SSatish Balay       hlru_global(*a->hmatrix_gpu, uu, ldu, vv, ldv, U->cmap->N, s, handle);
1878f17f7ee2SSatish Balay #if defined(PETSC_H2OPUS_USE_GPU)
1879f17f7ee2SSatish Balay       PetscCall(MatDenseCUDARestoreArrayRead(U, &u));
1880f17f7ee2SSatish Balay       PetscCall(MatDenseCUDARestoreArrayRead(V, &v));
1881f17f7ee2SSatish Balay       if (flgU) PetscCall(MatConvert(U, MATDENSE, MAT_INPLACE_MATRIX, &U));
1882f17f7ee2SSatish Balay       if (flgV) PetscCall(MatConvert(V, MATDENSE, MAT_INPLACE_MATRIX, &V));
1883f17f7ee2SSatish Balay #endif
1884f17f7ee2SSatish Balay     }
1885f17f7ee2SSatish Balay     PetscCall(PetscLogEventEnd(MAT_H2Opus_LR, A, 0, 0, 0));
1886f17f7ee2SSatish Balay     a->orthogonal = PETSC_FALSE;
1887f17f7ee2SSatish Balay   }
1888f17f7ee2SSatish Balay   PetscFunctionReturn(0);
1889f17f7ee2SSatish Balay }
1890f17f7ee2SSatish Balay #endif
1891